Discontinuous Galerkin Library
#include "dg/algorithm.h"
dg::RealCartesianMPIGrid3d< real_type > Struct Template Reference

The mpi version of RealCartesianGrid3d. More...

Inheritance diagram for dg::RealCartesianMPIGrid3d< real_type >:
[legend]

Public Types

using perpendicular_grid = RealCartesianMPIGrid2d< real_type >
 
- Public Types inherited from dg::aRealMPITopology3d< real_type >
using value_type = real_type
 
using host_vector = MPI_Vector< thrust::host_vector< real_type > >
 The host vector type used by host functions like evaluate. More...
 
using host_grid = RealMPIGrid3d< real_type >
 

Public Member Functions

 RealCartesianMPIGrid3d (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)
 Equal polynomial coefficients. More...
 
 RealCartesianMPIGrid3d (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, bc bcx, bc bcy, bc bcz, MPI_Comm comm)
 Equal polynomial coefficients. More...
 
 RealCartesianMPIGrid3d (const dg::RealMPIGrid3d< real_type > &g)
 Implicit type conversion from RealMPIGrid3d. More...
 
virtual RealCartesianMPIGrid3dclone () const override final
 Geometries are cloneable. More...
 
 RealCartesianMPIGrid3d (RealGrid1d< real_type > gx, RealGrid1d< real_type > gy, RealGrid1d< real_type > gz, MPI_Comm comm)
 Construct a 3d topology as the product of three 1d grids. More...
 
virtual RealCartesianGrid3d< real_type > * global_geometry () const override final
 Construct the global non-MPI geometry. More...
 
- Public Member Functions inherited from dg::aRealProductMPIGeometry3d< real_type >
aRealMPIGeometry2d< real_type > * perp_grid () const
 The grid made up by the first two dimensions. More...
 
virtual ~aRealProductMPIGeometry3d ()=default
 allow deletion through base class pointer More...
 
virtual aRealProductMPIGeometry3dclone () const =0
 Geometries are cloneable. More...
 
- Public Member Functions inherited from dg::aRealMPIGeometry3d< real_type >
SparseTensor< MPI_Vector< thrust::host_vector< real_type > > > jacobian () const
 The Jacobian of the coordinate transformation from physical to computational space. More...
 
SparseTensor< MPI_Vector< thrust::host_vector< real_type > > > metric () const
 The (inverse) metric tensor of the coordinate system. More...
 
std::vector< MPI_Vector< thrust::host_vector< real_type > > > map () const
 The coordinate map from computational to physical space. More...
 
virtual aRealMPIGeometry3dclone () const =0
 Geometries are cloneable. More...
 
virtual aRealGeometry3d< real_type > * global_geometry () const =0
 Construct the global non-MPI geometry. More...
 
virtual ~aRealMPIGeometry3d ()=default
 allow deletion through base class pointer More...
 
- Public Member Functions inherited from dg::aRealMPITopology3d< real_type >
real_type x0 () const
 Return global x0. More...
 
real_type x1 () const
 Return global x1. More...
 
real_type y0 () const
 Return global y0. More...
 
real_type y1 () const
 Return global y1. More...
 
real_type z0 () const
 Return global z0. More...
 
real_type z1 () const
 Return global z1. More...
 
real_type lx () const
 Return global lx. More...
 
real_type ly () const
 Return global ly. More...
 
real_type lz () const
 Return global lz. More...
 
real_type hx () const
 Return global hx. More...
 
real_type hy () const
 Return global hy. More...
 
real_type hz () const
 Return global hz. More...
 
unsigned n () const
 Return n. More...
 
unsigned nx () const
 number of polynomial coefficients in x More...
 
unsigned ny () const
 number of polynomial coefficients in y More...
 
unsigned nz () const
 number of polynomial coefficients in z More...
 
unsigned Nx () const
 Return the global number of cells. More...
 
unsigned Ny () const
 Return the global number of cells. More...
 
unsigned Nz () const
 Return the global number of cells. More...
 
bc bcx () const
 global x boundary More...
 
bc bcy () const
 global y boundary More...
 
bc bcz () const
 global z boundary More...
 
MPI_Comm communicator () const
 Return mpi cartesian communicator that is used in this grid. More...
 
MPI_Comm get_perp_comm () const
 MPI Cartesian communicator in the first two dimensions (x and y) More...
 
const DLT< real_type > & dlt () const
 The Discrete Legendre Transformation. More...
 
const DLT< real_type > & dltx () const
 
const DLT< real_type > & dlty () const
 
const DLT< real_type > & dltz () const
 
unsigned size () const
 The total global number of points. More...
 
unsigned local_size () const
 The total local number of points. More...
 
void display (std::ostream &os=std::cout) const
 Display global and local grid paramters. More...
 
int pidOf (real_type x, real_type y, real_type z) const
 Returns the pid of the process that holds the local grid surrounding the given point. More...
 
void multiplyCellNumbers (real_type fx, real_type fy)
 Multiply the number of cells with a given factor. More...
 
void set (unsigned new_n, unsigned new_Nx, unsigned new_Ny, unsigned new_Nz)
 Set the number of polynomials and cells. More...
 
void set (unsigned new_nx, unsigned new_Nx, unsigned new_ny, unsigned new_Ny, unsigned new_nz, unsigned new_Nz)
 Set the number of polynomials and cells. More...
 
bool local2globalIdx (int localIdx, int PID, int &globalIdx) const
 Map a local index plus the PID to a global vector index. More...
 
bool global2localIdx (int globalIdx, int &localIdx, int &PID) const
 Map a global vector index to a local vector Index and the corresponding PID. More...
 
const RealGrid3d< real_type > & local () const
 Return a non-MPI grid local for the calling process. More...
 
const RealGrid3d< real_type > & global () const
 Return the global non-MPI grid. More...
 

Additional Inherited Members

- Protected Member Functions inherited from dg::aRealProductMPIGeometry3d< real_type >
 aRealProductMPIGeometry3d (const aRealProductMPIGeometry3d &src)=default
 
aRealProductMPIGeometry3doperator= (const aRealProductMPIGeometry3d &src)=default
 
- Protected Member Functions inherited from dg::aRealMPIGeometry3d< real_type >
 aRealMPIGeometry3d (const aRealMPIGeometry3d &src)=default
 
aRealMPIGeometry3doperator= (const aRealMPIGeometry3d &src)=default
 
- Protected Member Functions inherited from dg::aRealMPITopology3d< real_type >
 ~aRealMPITopology3d ()=default
 disallow deletion through base class pointer More...
 
 aRealMPITopology3d (RealGrid1d< real_type > gx, RealGrid1d< real_type > gy, RealGrid1d< real_type > gz, MPI_Comm comm)
 Construct a 3d topology as the product of three 1d grids. More...
 
 aRealMPITopology3d (const aRealMPITopology3d &src)=default
 
aRealMPITopology3doperator= (const aRealMPITopology3d &src)=default
 
virtual void do_set (unsigned new_nx, unsigned new_Nx, unsigned new_ny, unsigned new_Ny, unsigned new_nz, unsigned new_Nz)=0
 

Detailed Description

template<class real_type>
struct dg::RealCartesianMPIGrid3d< real_type >

The mpi version of RealCartesianGrid3d.

Member Typedef Documentation

◆ perpendicular_grid

template<class real_type >
using dg::RealCartesianMPIGrid3d< real_type >::perpendicular_grid = RealCartesianMPIGrid2d<real_type>

Constructor & Destructor Documentation

◆ RealCartesianMPIGrid3d() [1/4]

template<class real_type >
dg::RealCartesianMPIGrid3d< real_type >::RealCartesianMPIGrid3d ( 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

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
comma three-dimensional Cartesian communicator
Note
the paramateres given in the constructor are global parameters

◆ RealCartesianMPIGrid3d() [2/4]

template<class real_type >
dg::RealCartesianMPIGrid3d< real_type >::RealCartesianMPIGrid3d ( 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,
bc  bcx,
bc  bcy,
bc  bcz,
MPI_Comm  comm 
)
inline

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
bcxboundary condition in x
bcyboundary condition in y
bczboundary condition in z
comma three-dimensional Cartesian communicator
Note
the paramateres given in the constructor are global parameters

◆ RealCartesianMPIGrid3d() [3/4]

template<class real_type >
dg::RealCartesianMPIGrid3d< real_type >::RealCartesianMPIGrid3d ( const dg::RealMPIGrid3d< real_type > &  g)
inline

Implicit type conversion from RealMPIGrid3d.

Parameters
gexisting grid object

◆ RealCartesianMPIGrid3d() [4/4]

template<class real_type >
dg::RealCartesianMPIGrid3d< real_type >::RealCartesianMPIGrid3d ( RealGrid1d< real_type >  gx,
RealGrid1d< real_type >  gy,
RealGrid1d< real_type >  gz,
MPI_Comm  comm 
)
inline

Construct a 3d topology as the product of three 1d grids.

The simplest implementation of aRealTopology3d.
Definition: grid.h:844
bc bcz() const
global z boundary
Definition: mpi_grid.h:454
unsigned nz() const
number of polynomial coefficients in z
Definition: mpi_grid.h:418
unsigned Nz() const
Return the global number of cells.
Definition: mpi_grid.h:436
real_type x1() const
Return global x1.
Definition: mpi_grid.h:346
unsigned Ny() const
Return the global number of cells.
Definition: mpi_grid.h:430
real_type y1() const
Return global y1.
Definition: mpi_grid.h:358
bc bcy() const
global y boundary
Definition: mpi_grid.h:448
real_type z1() const
Return global z1.
Definition: mpi_grid.h:370
bc bcx() const
global x boundary
Definition: mpi_grid.h:442
unsigned ny() const
number of polynomial coefficients in y
Definition: mpi_grid.h:416
real_type z0() const
Return global z0.
Definition: mpi_grid.h:364
unsigned Nx() const
Return the global number of cells.
Definition: mpi_grid.h:424
unsigned nx() const
number of polynomial coefficients in x
Definition: mpi_grid.h:414
real_type y0() const
Return global y0.
Definition: mpi_grid.h:352
real_type x0() const
Return global x0.
Definition: mpi_grid.h:340
Parameters
gxa Grid1d in x - direction
gya Grid1d in y - direction
gza Grid1d in z - direction
comma three-dimensional Cartesian communicator
Note
the paramateres given in the constructor are global parameters

Member Function Documentation

◆ clone()

template<class real_type >
virtual RealCartesianMPIGrid3d * dg::RealCartesianMPIGrid3d< real_type >::clone ( ) const
inlinefinaloverridevirtual

Geometries are cloneable.

Implements dg::aRealProductMPIGeometry3d< real_type >.

◆ global_geometry()

template<class real_type >
virtual RealCartesianGrid3d< real_type > * dg::RealCartesianMPIGrid3d< real_type >::global_geometry ( ) const
inlinefinaloverridevirtual

Construct the global non-MPI geometry.

Implements dg::aRealMPIGeometry3d< real_type >.


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