Discontinuous Galerkin Library
#include "dg/algorithm.h"
|
An abstract base class for MPI distributed Nd-dimensional dG grids. More...
Public Types | |
using | value_type = real_type |
value type of abscissas and weights | |
using | host_vector = MPI_Vector<thrust::host_vector<real_type>> |
using | host_grid = RealMPIGrid<real_type, Nd> |
Public Member Functions | |
unsigned | shape (unsigned u=0) const |
\( n_u N_u\) the total number of points of an axis | |
host_vector | abscissas (unsigned u=0) const |
Get the grid abscissas of the u axis. | |
host_vector | weights (unsigned u=0) const |
Get the weights of the u axis. | |
std::array< unsigned, Nd > | get_shape () const |
\( n_u N_u\) the total number of points of an axis | |
std::array< host_vector, Nd > | get_abscissas () const |
Construct abscissas for all axes. | |
std::array< host_vector, Nd > | get_weights () const |
Construct weights for all axes. | |
std::array< real_type, Nd > | get_p () const |
Get left boundary point \( \vec p\). | |
std::array< real_type, Nd > | get_q () const |
Get right boundary point \( \vec q\). | |
std::array< real_type, Nd > | get_l () const |
Get grid length \( l_u = q_u - p_u\) for all axes. | |
std::array< real_type, Nd > | get_h () const |
Get grid constant \( h_u = \frac{q_u - p_u}{N_u}\) for all axes. | |
std::array< unsigned, Nd > | get_N () const |
Get number of cells \( N_u\) for all axes. | |
std::array< unsigned, Nd > | get_n () const |
Get number of polynomial coefficients \( n_u\) for all axes. | |
std::array< dg::bc, Nd > | get_bc () const |
Get boundary condition \( b_u\) for all axes. | |
std::array< MPI_Comm, Nd > | get_comms () const |
Get 1d Cartesian communicator \( c_u\) for all axes. | |
real_type | p (unsigned u=0) const |
Get left boundary point \( p_u\) for axis u . | |
real_type | q (unsigned u=0) const |
Get right boundary point \( q_u\) for axis u . | |
real_type | h (unsigned u=0) const |
Get grid constant \( h_u = \frac{q_u - p_u}{N_u}\) for axis u . | |
real_type | l (unsigned u=0) const |
Get grid length \( l_u = q_u - p_u\) for axis u . | |
unsigned | n (unsigned u=0) const |
Get number of polynomial coefficients \( n_u\) for axis u . | |
unsigned | N (unsigned u=0) const |
Get number of cells \( N_u\) for axis u . | |
dg::bc | bc (unsigned u=0) const |
Get boundary condition \( b_u\) for axis u . | |
MPI_Comm | comm (unsigned u) const |
Get 1d Cartesian communicator \( c_u\) for axis u . | |
template<size_t Md = Nd> | |
std::enable_if_t< Md==1, MPI_Comm > | comm () const |
Equivalent to comm(0) | |
MPI_Comm | communicator () const |
Return Nd dimensional MPI cartesian communicator that is used in this grid. | |
template<size_t Md = Nd> | |
std::enable_if_t<(Md >=2), MPI_Comm > | get_perp_comm () const |
MPI Cartesian communicator in the first two dimensions (x and y) | |
RealMPIGrid< real_type, 1 > | grid (unsigned u) const |
Get axis u as a 1d grid. | |
RealMPIGrid< real_type, 1 > | axis (unsigned u) const |
An alias for "grid". | |
const RealGrid< real_type, Nd > & | local () const |
The local grid as a shared memory grid. | |
const RealGrid< real_type, Nd > & | global () const |
The global grid as a shared memory grid. | |
template<size_t Md = Nd> | |
real_type | x0 () const |
Equivalent to p(0) | |
template<size_t Md = Nd> | |
real_type | x1 () const |
Equivalent to p(1) | |
template<size_t Md = Nd> | |
real_type | y0 () const |
Equivalent to p(2) | |
template<size_t Md = Nd> | |
real_type | y1 () const |
Equivalent to q(0) | |
template<size_t Md = Nd> | |
real_type | z0 () const |
Equivalent to q(1) | |
template<size_t Md = Nd> | |
real_type | z1 () const |
Equivalent to q(2) | |
template<size_t Md = Nd> | |
real_type | lx () const |
Equivalent to l(0) | |
template<size_t Md = Nd> | |
real_type | ly () const |
Equivalent to l(1) | |
template<size_t Md = Nd> | |
real_type | lz () const |
Equivalent to l(2) | |
template<size_t Md = Nd> | |
real_type | hx () const |
Equivalent to h(0) | |
template<size_t Md = Nd> | |
real_type | hy () const |
Equivalent to h(1) | |
template<size_t Md = Nd> | |
real_type | hz () const |
Equivalent to h(2) | |
template<size_t Md = Nd> | |
unsigned | nx () const |
Equivalent to n(0) | |
template<size_t Md = Nd> | |
unsigned | ny () const |
Equivalent to n(1) | |
template<size_t Md = Nd> | |
unsigned | nz () const |
Equivalent to n(2) | |
template<size_t Md = Nd> | |
unsigned | Nx () const |
Equivalent to N(0) | |
template<size_t Md = Nd> | |
unsigned | Ny () const |
Equivalent to N(1) | |
template<size_t Md = Nd> | |
unsigned | Nz () const |
Equivalent to N(2) | |
template<size_t Md = Nd> | |
dg::bc | bcx () const |
Equivalent to bc(0) | |
template<size_t Md = Nd> | |
dg::bc | bcy () const |
Equivalent to bc(1) | |
template<size_t Md = Nd> | |
dg::bc | bcz () const |
Equivalent to bc(2) | |
template<size_t Md = Nd> | |
RealMPIGrid< real_type, 1 > | gx () const |
Equivalent to grid(0) | |
template<size_t Md = Nd> | |
RealMPIGrid< real_type, 1 > | gy () const |
Equivalent to grid(1) | |
template<size_t Md = Nd> | |
RealMPIGrid< real_type, 1 > | gz () const |
Equivalent to grid(2) | |
template<size_t Md = Nd> | |
std::enable_if_t<(Md >=2), void > | multiplyCellNumbers (real_type fx, real_type fy) |
Multiply the number of cells in the first two dimensions with a given factor. | |
template<size_t Md = Nd> | |
std::enable_if_t<(Md==1), void > | set (unsigned new_n, unsigned new_Nx) |
Set n and N in a 1-dimensional grid. | |
template<size_t Md = Nd> | |
std::enable_if_t<(Md==2), void > | set (unsigned new_n, unsigned new_Nx, unsigned new_Ny) |
Set n and N in a 2-dimensional grid. | |
template<size_t Md = Nd> | |
std::enable_if_t<(Md==3), void > | set (unsigned new_n, unsigned new_Nx, unsigned new_Ny, unsigned new_Nz) |
Set n and N in a 3-dimensional grid. | |
void | set (unsigned new_n, std::array< unsigned, Nd > new_N) |
Same as set( {new_n, new_n,...}, new_N); | |
void | set_axis (unsigned coord, unsigned new_n, unsigned new_N) |
Set n and N for axis coord . | |
void | set (std::array< unsigned, Nd > new_n, std::array< unsigned, Nd > new_N) |
Set the number of polynomials and cells. | |
void | set_pq (std::array< real_type, Nd > new_p, std::array< real_type, Nd > new_q) |
Reset the boundaries of the grid. | |
void | set_bcs (std::array< dg::bc, Nd > new_bcs) |
Reset the boundary conditions of the grid. | |
unsigned | size () const |
The total global number of points. | |
unsigned | local_size () const |
The total local number of points. | |
void | display (std::ostream &os=std::cout) const |
Display global and local grid paramters. | |
template<size_t Md = Nd> | |
std::enable_if_t<(Md==1), bool > | contains (real_type x) const |
Check if the grid contains a point. | |
template<class Vector > | |
bool | contains (const Vector &x) const |
Check if the grid contains a point. | |
bool | local2globalIdx (int localIdx, int rank, int &globalIdx) const |
Convert the index of a local vector and the rank of the containing process to a global index. | |
bool | global2localIdx (int globalIdx, int &localIdx, int &rank) const |
Convert the global index of a vector to a local index and the rank of the containing process. | |
std::array< unsigned, Nd > | start () const |
The global start coordinate in C-order of dg::file::MPINcHyperslab that the local grid represents. | |
std::array< unsigned, Nd > | count () const |
Count vector in C-order for dg::file::MPINcHyperslab . | |
Static Public Member Functions | |
static constexpr unsigned | ndim () |
Dimensionality == Nd. | |
Protected Member Functions | |
~aRealMPITopology ()=default | |
disallow deletion through base class pointer | |
aRealMPITopology ()=default | |
default constructor | |
aRealMPITopology (std::array< real_type, Nd > p, std::array< real_type, Nd > q, std::array< unsigned, Nd > n, std::array< unsigned, Nd > N, std::array< dg::bc, Nd > bcs, std::array< MPI_Comm, Nd > comms) | |
Construct a topology directly from points and dimensions. | |
aRealMPITopology (const std::array< RealMPIGrid< real_type, 1 >, Nd > &axes) | |
Construct a topology as the product of 1d axes grids. | |
aRealMPITopology (const aRealMPITopology &src)=default | |
aRealMPITopology & | operator= (const aRealMPITopology &src)=default |
virtual void | do_set (std::array< unsigned, Nd > new_n, std::array< unsigned, Nd > new_N)=0 |
Set the number of polynomials and cells. | |
virtual void | do_set_pq (std::array< real_type, Nd > new_p, std::array< real_type, Nd > new_q)=0 |
Reset the boundaries of the grid. | |
virtual void | do_set (std::array< dg::bc, Nd > new_bcs)=0 |
Reset the boundary conditions of the grid. | |
An abstract base class for MPI distributed Nd-dimensional dG grids.
In the dg library we have MPI grids as a separate type from shared memory grids. The design goal is that MPI grids generally can be used the same way as just a single global shared memory grid. A notable exception is the constructor, when MPI grids depend on a Cartesian MPI communicator (that has to be constructed beforehand with e.g. dg::mpi_cart_create
)
This grid defines a discretization of the \(N_d\) dimensional hypercube given by
\[ [\vec p, \vec q] = [p_0, p_1] \times [p_1,q_1] \times ... \times [p_{N_d-1}, q_{N_d-1}]\]
Each axis \( [p_i, q_i] \) is discretized using \( N_i\) equidistant cells. Each cells is then further discretized using \( n_i\) Gauss-Legendre nodes. The Gauss-Legendre nodes are tabulated by the dg::DLT
class for \( 1\leq n_i \leq 20\). Each axis further can have a boundary condition \( b_i \) that is given by dg::bc
For more information dG methods see
In MPI we want to distribute the above hypercube among processes in a Cartesian communicator of same dimensionality \( N_d\). This is done by evenly distributing the global number of cells \( N_u\) in each dimension among the processes in the corresponding dimension in the Cartesian communicator. Each process gets a local number of cells \( N_{ur} \approx N_u / s_u\), where \( s_u \) is the size of the Cartesian communicator in direction \( u\). The approximation becomes an equality if \( s_u\) evenly divides \( N_u\) otherwise the remainder is distributed among the participating processes (such that some have one cell more or less than others). We have
\[ N_u = \sum_{r=0}^{s-1} N_{ur}\]
The number of polynomial coefficients and the boundary condition is the same for all processes of an axis. Each axis among processes is thus
\[ [p_u,q_u] = [p_{u0}, q_{u0}],[p_{u1}, q_{u1}],...,[p_{us-1}, q_{us-1}] \]
The local boundaries are determined such that the local \( h_{ur} = \frac{q_{ur}-p_{ur}}{N_{ur}} = \frac{q_u-p_u}{N_u}\)
This class in essence provides a collection of getters and setters for the aforementioned parameters together with the abscissas
and weights
members that are necessary for dg::evaluate
and dg::create::weights
.
For code readability in many physical contexts the first 3 dimensions get special names \( x,\ y,\ z\) i.e. \( \vec p = (x_0, y_0, z_0)\) and \(
\vec q = (x_1, y_1, z_1)\). The class provides coresponding getters and setters like e.g. x0()
or Nx()
for p(0)
and N(0)
local()
function, which provides a grid with all local quantitiesLastly, we provide start
and count
members such that the grid can be used as a dg::file::MPINcHyperslab
in NetCDF output in dg::file
real_type | Determines value type of abscissas and weights |
Nd | The number of dimensions \( N_d\) |
using dg::aRealMPITopology< real_type, Nd >::host_grid = RealMPIGrid<real_type, Nd> |
using dg::aRealMPITopology< real_type, Nd >::host_vector = MPI_Vector<thrust::host_vector<real_type>> |
vector type of abscissas and weights; Can be used to recognize MPI grid via: dg::is_vector_v< typename Topology::host_vector, MPIVectorTag>
using dg::aRealMPITopology< real_type, Nd >::value_type = real_type |
value type of abscissas and weights
|
protecteddefault |
disallow deletion through base class pointer
|
protecteddefault |
default constructor
|
inlineprotected |
Construct a topology directly from points and dimensions.
p | left boundary point |
q | right boundary point |
n | number of polynomial coefficients for each axis |
N | number of cells for each axis |
bcs | boundary condition for each axis |
comms | The 1d MPI Cartesian communicators for each axis that make up an Nd dimensional Cartesian communicator (see dg::mpi_cart_split , dg::mpi_cart_kron ) |
|
inlineprotected |
Construct a topology as the product of 1d axes grids.
axes | One-dimensional grids for each dimension |
|
protecteddefault |
explicit copy constructor (default)
src | source |
|
inline |
Get the grid abscissas of the u
axis.
global()
grid. u | Axis number u<Nd |
|
inline |
An alias for "grid".
u | Axis number u<Nd |
|
inline |
Get boundary condition \( b_u\) for axis u
.
u | Axis number u<Nd |
u
|
inline |
Equivalent to bc(0)
|
inline |
Equivalent to bc(1)
|
inline |
Equivalent to bc(2)
|
inline |
Equivalent to comm(0)
|
inline |
Get 1d Cartesian communicator \( c_u\) for axis u
.
u | Axis number u<Nd |
|
inline |
Return Nd dimensional MPI cartesian communicator that is used in this grid.
|
inline |
Check if the grid contains a point.
Used for example in integrate_in_domain
method of dg::AdaptiveTimeloop
x | point to check |
|
inline |
Check if the grid contains a point.
Used for example in integrate_in_domain
method of dg::AdaptiveTimeloop
x | point to check |
|
inline |
Count vector in C-order for dg::file::MPINcHyperslab
.
Used to construct dg::file::MPINcHyperslab
together with start()
dg::evaluate
and dg::kronecker
make the 0 dimension/ 1st argument the fastest varying, so we return the reverse order of local()
.get_shape()
|
inline |
Display global and local grid paramters.
os | output stream |
|
protectedpure virtual |
Reset the boundary conditions of the grid.
new_bcs | new boundary condition |
|
protectedpure virtual |
Set the number of polynomials and cells.
new_n | new number of Gaussian nodes in each dimension |
new_N | new number of cells in each dimension |
|
protectedpure virtual |
Reset the boundaries of the grid.
new_p | new left boundary |
new_q | new right boundary ( > x0) |
|
inline |
Construct abscissas for all axes.
|
inline |
Get boundary condition \( b_u\) for all axes.
u
|
inline |
Get 1d Cartesian communicator \( c_u\) for all axes.
u
|
inline |
Get grid constant \( h_u = \frac{q_u - p_u}{N_u}\) for all axes.
u
|
inline |
Get grid length \( l_u = q_u - p_u\) for all axes.
u
|
inline |
Get number of cells \( N_u\) for all axes.
u
|
inline |
Get number of polynomial coefficients \( n_u\) for all axes.
u
|
inline |
Get left boundary point \( \vec p\).
|
inline |
MPI Cartesian communicator in the first two dimensions (x and y)
|
inline |
Get right boundary point \( \vec q\).
|
inline |
\( n_u N_u\) the total number of points of an axis
|
inline |
Construct weights for all axes.
|
inline |
The global grid as a shared memory grid.
|
inline |
Convert the global index of a vector to a local index and the rank of the containing process.
globalIdx | Global index of the element |
localIdx | (write only) Index in a local chunk of a vector |
rank | (write only) Rank of the process that contains the globalIdx |
|
inline |
Get axis u
as a 1d grid.
u | Axis number u<Nd |
|
inline |
Equivalent to grid(0)
|
inline |
Equivalent to grid(1)
|
inline |
Equivalent to grid(2)
|
inline |
Get grid constant \( h_u = \frac{q_u - p_u}{N_u}\) for axis u
.
u | Axis number u<Nd |
u
|
inline |
Equivalent to h(0)
|
inline |
Equivalent to h(1)
|
inline |
Equivalent to h(2)
|
inline |
Get grid length \( l_u = q_u - p_u\) for axis u
.
u | Axis number u<Nd |
u
|
inline |
The local grid as a shared memory grid.
|
inline |
Convert the index of a local vector and the rank of the containing process to a global index.
localIdx | Index in a local chunk of a vector |
rank | Rank of the process that contains the localIdx |
globalIdx | (write only) Global index of the element |
|
inline |
|
inline |
Equivalent to l(0)
|
inline |
Equivalent to l(1)
|
inline |
Equivalent to l(2)
|
inline |
Multiply the number of cells in the first two dimensions with a given factor.
With this function you can resize the grid ignorantly of its current size
fx | new global number of cells is fx*Nx() |
fy | new global number of cells is fy*Ny() The remaining dimensions are left unchanged |
|
inline |
Get number of cells \( N_u\) for axis u
.
u | Axis number u<Nd |
u
|
inline |
Get number of polynomial coefficients \( n_u\) for axis u
.
u | Axis number u<Nd |
u
|
inlinestaticconstexpr |
Dimensionality == Nd.
|
inline |
Equivalent to N(0)
|
inline |
Equivalent to n(0)
|
inline |
Equivalent to N(1)
|
inline |
Equivalent to n(1)
|
inline |
Equivalent to N(2)
|
inline |
Equivalent to n(2)
|
protecteddefault |
explicit assignment operator (default)
src | source |
|
inline |
Get left boundary point \( p_u\) for axis u
.
u | Axis number u<Nd |
u
|
inline |
Get right boundary point \( q_u\) for axis u
.
u | Axis number u<Nd |
u
|
inline |
Set the number of polynomials and cells.
new_n | new number of Gaussian nodes in each dimension |
new_N | new number of cells in each dimension |
|
inline |
Same as set( {new_n, new_n,...}, new_N);
|
inline |
Set n and N in a 1-dimensional grid.
new_n | new number of Gaussian nodes |
new_Nx | new number of cells |
|
inline |
Set n and N in a 2-dimensional grid.
new_n | new number of Gaussian nodes in x and y |
new_Nx | new number of cells |
new_Ny | new number of cells |
|
inline |
Set n and N in a 3-dimensional grid.
Same as set({new_n,new_n,1}, {new_Nx,new_Ny,new_Nz})
new_n | new number of Gaussian nodes in x and y |
nz
to 1 new_Nx | new number of cells in x |
new_Ny | new number of cells in y |
new_Nz | new number of cells in z |
|
inline |
Set n and N for axis coord
.
coord | Axis coord<Nd |
new_n | new number of Gaussian nodes of axis coord |
new_N | new number of cells |
|
inline |
Reset the boundary conditions of the grid.
new_bcs | new boundary condition |
|
inline |
Reset the boundaries of the grid.
new_p | new left boundary |
new_q | new right boundary ( > x0) |
|
inline |
\( n_u N_u\) the total number of points of an axis
u | Axis number u<Nd |
|
inline |
The total global number of points.
|
inline |
The global start coordinate in C-order of dg::file::MPINcHyperslab
that the local grid represents.
Used to construct dg::file::MPINcHyperslab
together with count()
dg::evaluate
and dg::kronecker
make the 0 dimension/ 1st argument the fastest varying.
|
inline |
Get the weights of the u
axis.
global()
grid. u | Axis number u<Nd |
|
inline |
Equivalent to p(0)
|
inline |
Equivalent to p(1)
|
inline |
Equivalent to p(2)
|
inline |
Equivalent to q(0)
|
inline |
Equivalent to q(1)
|
inline |
Equivalent to q(2)