Discontinuous Galerkin Library
#include "dg/algorithm.h"
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
dg::RealMPIGrid< real_type, Nd > Struct Template Reference

The simplest implementation of aRealMPITopology3d. More...

Inheritance diagram for dg::RealMPIGrid< real_type, Nd >:
[legend]

Public Member Functions

 RealMPIGrid ()=default
 construct an empty grid this leaves the access functions undefined
 
template<size_t Md = Nd>
 RealMPIGrid (real_type x0, real_type x1, unsigned n, unsigned N, MPI_Comm comm)
 
template<size_t Md = Nd>
 RealMPIGrid (real_type x0, real_type x1, unsigned n, unsigned Nx, dg::bc bcx, MPI_Comm comm)
 1D grid
 
template<size_t Md = Nd>
 RealMPIGrid (real_type x0, real_type x1, real_type y0, real_type y1, unsigned n, unsigned Nx, unsigned Ny, MPI_Comm comm)
 Construct with equal polynomial coefficients.
 
template<size_t Md = Nd>
 RealMPIGrid (real_type x0, real_type x1, real_type y0, real_type y1, unsigned n, unsigned Nx, unsigned Ny, dg::bc bcx, dg::bc bcy, MPI_Comm comm)
 Construct with equal polynomial coefficients.
 
template<size_t Md = Nd>
 RealMPIGrid (real_type x0, real_type x1, real_type y0, real_type y1, real_type z0, real_type z1, unsigned n, unsigned Nx, unsigned Ny, unsigned Nz, MPI_Comm comm)
 Construct with equal polynomial coefficients.
 
template<size_t Md = Nd>
 RealMPIGrid (real_type x0, real_type x1, real_type y0, real_type y1, real_type z0, real_type z1, unsigned n, unsigned Nx, unsigned Ny, unsigned Nz, dg::bc bcx, dg::bc bcy, dg::bc bcz, MPI_Comm comm)
 Construct with equal polynomial coefficients.
 
 RealMPIGrid (const std::array< RealMPIGrid< real_type, 1 >, Nd > &axes)
 Construct a topology as the product of 1d axes grids.
 
template<class ... Grid1ds>
 RealMPIGrid (const RealMPIGrid< real_type, 1 > &g0, const Grid1ds &...gs)
 Construct from given 1d grids Equivalent to RealMPIGrid( std::array{g0,gs...})
 
 RealMPIGrid (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.
 
 RealMPIGrid (const aRealMPITopology< real_type, Nd > &src)
 
- Public Member Functions inherited from dg::aRealMPITopology< real_type, Nd >
unsigned shape (unsigned u=0) const
 nuNu 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
 nuNu 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 p.
 
std::array< real_type, Nd > get_q () const
 Get right boundary point q.
 
std::array< real_type, Nd > get_l () const
 Get grid length lu=qupu for all axes.
 
std::array< real_type, Nd > get_h () const
 Get grid constant hu=qupuNu for all axes.
 
std::array< unsigned, Nd > get_N () const
 Get number of cells Nu for all axes.
 
std::array< unsigned, Nd > get_n () const
 Get number of polynomial coefficients nu for all axes.
 
std::array< dg::bc, Nd > get_bc () const
 Get boundary condition bu for all axes.
 
std::array< MPI_Comm, Nd > get_comms () const
 Get 1d Cartesian communicator cu for all axes.
 
real_type p (unsigned u=0) const
 Get left boundary point pu for axis u.
 
real_type q (unsigned u=0) const
 Get right boundary point qu for axis u.
 
real_type h (unsigned u=0) const
 Get grid constant hu=qupuNu for axis u.
 
real_type l (unsigned u=0) const
 Get grid length lu=qupu for axis u.
 
unsigned n (unsigned u=0) const
 Get number of polynomial coefficients nu for axis u.
 
unsigned N (unsigned u=0) const
 Get number of cells Nu for axis u.
 
dg::bc bc (unsigned u=0) const
 Get boundary condition bu for axis u.
 
MPI_Comm comm (unsigned u) const
 Get 1d Cartesian communicator cu 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.
 

Additional Inherited Members

- Public Types inherited from dg::aRealMPITopology< real_type, Nd >
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>
 
- Static Public Member Functions inherited from dg::aRealMPITopology< real_type, Nd >
static constexpr unsigned ndim ()
 Dimensionality == Nd.
 
- Protected Member Functions inherited from dg::aRealMPITopology< real_type, Nd >
 ~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
 
aRealMPITopologyoperator= (const aRealMPITopology &src)=default
 

Detailed Description

template<class real_type, size_t Nd>
struct dg::RealMPIGrid< real_type, Nd >

The simplest implementation of aRealMPITopology3d.

// This code snippet demonstrates how to discretize and compute the norm of a
// function on both shared and distributed memory systems
#ifdef WITH_MPI
// In an MPI environment define a 2d Cartesian communicator
MPI_Comm comm2d = dg::mpi_cart_create( MPI_COMM_WORLD, {0,0}, {1,1});
#endif
// create a grid of the domain [0,2]x[0,2] with 20 cells in x and y and
// 3 polynomial coefficients
dg::x::Grid2d g2d( 0, 2, 0, 2, 3, 20, 20
#ifdef WITH_MPI
, comm2d // in MPI distribute among processes
#endif
);
// define the function to integrate
const double amp = 2.;
auto function = [=] (double x, double y){
return amp*exp(x)*exp(y);
};
// create the Gaussian weights (volume form) for the integration
const dg::x::DVec w2d = dg::create::weights( g2d);
// discretize the function on the grid
const dg::x::DVec vec = dg::evaluate( function, g2d);
// now compute the scalar product (the L2 norm)
double square_norm = dg::blas2::dot( vec, w2d, vec);
// norm is now (or at least converges to): (e^4-1)^2
CHECK( fabs( square_norm - (exp(4)-1)*(exp(4)-1)) < 1e-6);

Constructor & Destructor Documentation

◆ RealMPIGrid() [1/11]

template<class real_type , size_t Nd>
dg::RealMPIGrid< real_type, Nd >::RealMPIGrid ( )
default

construct an empty grid this leaves the access functions undefined

◆ RealMPIGrid() [2/11]

template<class real_type , size_t Nd>
template<size_t Md = Nd>
dg::RealMPIGrid< real_type, Nd >::RealMPIGrid ( real_type x0,
real_type x1,
unsigned n,
unsigned N,
MPI_Comm comm )
inline

◆ RealMPIGrid() [3/11]

template<class real_type , size_t Nd>
template<size_t Md = Nd>
dg::RealMPIGrid< real_type, Nd >::RealMPIGrid ( real_type x0,
real_type x1,
unsigned n,
unsigned Nx,
dg::bc bcx,
MPI_Comm comm )
inline

1D grid

Parameters
x0left boundary
x1right boundary
n# of polynomial coefficients (1<=n<=20)
Nx# of cells
bcxboundary conditions
Parameters
comma one-dimensional Cartesian communicator
Note
the parameters given in the constructor are global parameters

◆ RealMPIGrid() [4/11]

template<class real_type , size_t Nd>
template<size_t Md = Nd>
dg::RealMPIGrid< real_type, Nd >::RealMPIGrid ( real_type x0,
real_type x1,
real_type y0,
real_type y1,
unsigned n,
unsigned Nx,
unsigned Ny,
MPI_Comm comm )
inline

Construct with equal polynomial coefficients.

Parameters
x0left boundary in x
x1right boundary in x
y0lower boundary in y
y1upper boundary in y
n# of polynomial coefficients for both x and y dimension (1<=n<=20)
Nx# of points in x
Ny# of points in y
Parameters
comma two-dimensional Cartesian communicator
Note
the parameters given in the constructor are global parameters

◆ RealMPIGrid() [5/11]

template<class real_type , size_t Nd>
template<size_t Md = Nd>
dg::RealMPIGrid< real_type, Nd >::RealMPIGrid ( real_type x0,
real_type x1,
real_type y0,
real_type y1,
unsigned n,
unsigned Nx,
unsigned Ny,
dg::bc bcx,
dg::bc bcy,
MPI_Comm comm )
inline

Construct with equal polynomial coefficients.

Parameters
x0left boundary in x
x1right boundary in x
y0lower boundary in y
y1upper boundary in y
n# of polynomial coefficients for both x and y dimension (1<=n<=20)
Nx# of points in x
Ny# of points in y
Parameters
bcxboundary condition in x
bcyboundary condition in y
Parameters
comma two-dimensional Cartesian communicator
Note
the parameters given in the constructor are global parameters

◆ RealMPIGrid() [6/11]

template<class real_type , size_t Nd>
template<size_t Md = Nd>
dg::RealMPIGrid< real_type, Nd >::RealMPIGrid ( real_type x0,
real_type x1,
real_type y0,
real_type y1,
real_type z0,
real_type z1,
unsigned n,
unsigned Nx,
unsigned Ny,
unsigned Nz,
MPI_Comm comm )
inline

Construct with equal polynomial coefficients.

Parameters
x0left boundary in x
x1right boundary in x
y0lower boundary in y
y1upper boundary in y
z0lower boundary in z
z1upper boundary in z
n# of polynomial coefficients for x and y dimension ( z-direction is set to 1) (1<=n<=20)
Nx# of points in x
Ny# of points in y
Nz# of points in z
Parameters
comma three-dimensional Cartesian communicator
Note
the parameters given in the constructor are global parameters

◆ RealMPIGrid() [7/11]

template<class real_type , size_t Nd>
template<size_t Md = Nd>
dg::RealMPIGrid< real_type, Nd >::RealMPIGrid ( real_type x0,
real_type x1,
real_type y0,
real_type y1,
real_type z0,
real_type z1,
unsigned n,
unsigned Nx,
unsigned Ny,
unsigned Nz,
dg::bc bcx,
dg::bc bcy,
dg::bc bcz,
MPI_Comm comm )
inline

Construct with equal polynomial coefficients.

Parameters
x0left boundary in x
x1right boundary in x
y0lower boundary in y
y1upper boundary in y
z0lower boundary in z
z1upper boundary in z
n# of polynomial coefficients for x and y dimension ( z-direction is set to 1) (1<=n<=20)
Nx# of points in x
Ny# of points in y
Nz# of points in z
Parameters
bcxboundary condition in x
bcyboundary condition in y
bczboundary condition in z
Parameters
comma three-dimensional Cartesian communicator
Note
the parameters given in the constructor are global parameters

◆ RealMPIGrid() [8/11]

template<class real_type , size_t Nd>
dg::RealMPIGrid< real_type, Nd >::RealMPIGrid ( const std::array< RealMPIGrid< real_type, 1 >, Nd > & axes)
inline

Construct a topology as the product of 1d axes grids.

Parameters
axesOne-dimensional grids for each dimension

◆ RealMPIGrid() [9/11]

template<class real_type , size_t Nd>
template<class ... Grid1ds>
dg::RealMPIGrid< real_type, Nd >::RealMPIGrid ( const RealMPIGrid< real_type, 1 > & g0,
const Grid1ds &... gs )
inline

Construct from given 1d grids Equivalent to RealMPIGrid( std::array{g0,gs...})

Parameters
g0Axis 0 grid
gsmore axes

◆ RealMPIGrid() [10/11]

template<class real_type , size_t Nd>
dg::RealMPIGrid< real_type, Nd >::RealMPIGrid ( 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 )
inline

Construct a topology directly from points and dimensions.

Parameters
pleft boundary point
qright boundary point
nnumber of polynomial coefficients for each axis
Nnumber of cells for each axis
bcsboundary condition for each axis
commsThe 1d MPI Cartesian communicators for each axis that make up an Nd dimensional Cartesian communicator (see dg::mpi_cart_split, dg::mpi_cart_kron)

◆ RealMPIGrid() [11/11]

template<class real_type , size_t Nd>
dg::RealMPIGrid< real_type, Nd >::RealMPIGrid ( const aRealMPITopology< real_type, Nd > & src)
inlineexplicit

allow explicit type conversion from any other topology

Parameters
srcsource

The documentation for this struct was generated from the following file: