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

3D MPI Grid class More...

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

Public Types

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

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...
 

Protected Member Functions

 ~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::aRealMPITopology3d< real_type >

3D MPI Grid class

Represents the global grid coordinates and the process topology. It just divides the given (global) box into nonoverlapping (local) subboxes that are attributed to each process

Note
a single cell is never divided across processes.
although it is abstract, objects are not meant to be hold on the heap via a base class pointer ( we protected the destructor)
Attention
The access functions nx() Nx() ,... all return the global parameters. If you want to have the local ones call the local() function.

Member Typedef Documentation

◆ host_grid

template<class real_type >
using dg::aRealMPITopology3d< real_type >::host_grid = RealMPIGrid3d<real_type>

◆ host_vector

template<class real_type >
using dg::aRealMPITopology3d< real_type >::host_vector = MPI_Vector<thrust::host_vector<real_type> >

The host vector type used by host functions like evaluate.

◆ value_type

template<class real_type >
using dg::aRealMPITopology3d< real_type >::value_type = real_type

Constructor & Destructor Documentation

◆ ~aRealMPITopology3d()

template<class real_type >
dg::aRealMPITopology3d< real_type >::~aRealMPITopology3d ( )
protecteddefault

disallow deletion through base class pointer

◆ aRealMPITopology3d() [1/2]

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

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

◆ aRealMPITopology3d() [2/2]

template<class real_type >
dg::aRealMPITopology3d< real_type >::aRealMPITopology3d ( const aRealMPITopology3d< real_type > &  src)
protecteddefault

explicit copy constructor (default)

Parameters
srcsource

Member Function Documentation

◆ bcx()

template<class real_type >
bc dg::aRealMPITopology3d< real_type >::bcx ( ) const
inline

global x boundary

Returns
boundary condition

◆ bcy()

template<class real_type >
bc dg::aRealMPITopology3d< real_type >::bcy ( ) const
inline

global y boundary

Returns
boundary condition

◆ bcz()

template<class real_type >
bc dg::aRealMPITopology3d< real_type >::bcz ( ) const
inline

global z boundary

Returns
boundary condition

◆ communicator()

template<class real_type >
MPI_Comm dg::aRealMPITopology3d< real_type >::communicator ( ) const
inline

Return mpi cartesian communicator that is used in this grid.

Returns
Communicator

◆ display()

template<class real_type >
void dg::aRealMPITopology3d< real_type >::display ( std::ostream &  os = std::cout) const
inline

Display global and local grid paramters.

Parameters
osoutput stream

◆ dlt()

template<class real_type >
const DLT< real_type > & dg::aRealMPITopology3d< real_type >::dlt ( ) const
inline

The Discrete Legendre Transformation.

Returns
DLT corresponding to n given in the constructor

◆ dltx()

template<class real_type >
const DLT< real_type > & dg::aRealMPITopology3d< real_type >::dltx ( ) const
inline

◆ dlty()

template<class real_type >
const DLT< real_type > & dg::aRealMPITopology3d< real_type >::dlty ( ) const
inline

◆ dltz()

template<class real_type >
const DLT< real_type > & dg::aRealMPITopology3d< real_type >::dltz ( ) const
inline

◆ do_set()

template<class real_type >
virtual void dg::aRealMPITopology3d< real_type >::do_set ( unsigned  new_nx,
unsigned  new_Nx,
unsigned  new_ny,
unsigned  new_Ny,
unsigned  new_nz,
unsigned  new_Nz 
)
protectedpure virtual

◆ get_perp_comm()

template<class real_type >
MPI_Comm dg::aRealMPITopology3d< real_type >::get_perp_comm ( ) const
inline

MPI Cartesian communicator in the first two dimensions (x and y)

Returns
2d Cartesian Communicator

◆ global()

template<class real_type >
const RealGrid3d< real_type > & dg::aRealMPITopology3d< real_type >::global ( ) const
inline

Return the global non-MPI grid.

The global grid contains the global boundaries and cell numbers. This is the grid that we would have to use in a non-MPI implementation. The global grid returns the same values for x0(), x1(), ..., Nx(), Ny(), ... as the grid class itself

Returns
non-MPI Grid object

◆ global2localIdx()

template<class real_type >
bool dg::aRealMPITopology3d< real_type >::global2localIdx ( int  globalIdx,
int &  localIdx,
int &  PID 
) const
inline

Map a global vector index to a local vector Index and the corresponding PID.

Parameters
globalIdxa global vector Index
localIdxcontains local vector index on output
PIDcontains corresponding PID in the communicator on output
Returns
true if successful, false if globalIdx is not part of the grid

◆ hx()

template<class real_type >
real_type dg::aRealMPITopology3d< real_type >::hx ( ) const
inline

Return global hx.

Returns
global grid constant

◆ hy()

template<class real_type >
real_type dg::aRealMPITopology3d< real_type >::hy ( ) const
inline

Return global hy.

Returns
global grid constant

◆ hz()

template<class real_type >
real_type dg::aRealMPITopology3d< real_type >::hz ( ) const
inline

Return global hz.

Returns
global grid constant

◆ local()

template<class real_type >
const RealGrid3d< real_type > & dg::aRealMPITopology3d< real_type >::local ( ) const
inline

Return a non-MPI grid local for the calling process.

The local grid contains the boundaries and cell numbers the calling process sees and is in charge of.

Returns
Grid object
Note
the boundary conditions in the local grid are not well defined since there might not actually be any boundaries

◆ local2globalIdx()

template<class real_type >
bool dg::aRealMPITopology3d< real_type >::local2globalIdx ( int  localIdx,
int  PID,
int &  globalIdx 
) const
inline

Map a local index plus the PID to a global vector index.

Parameters
localIdxa local vector index
PIDa PID in the communicator
globalIdxthe corresponding global vector Index (contains result on output)
Returns
true if successful, false if localIdx or PID is not part of the grid

◆ local_size()

template<class real_type >
unsigned dg::aRealMPITopology3d< real_type >::local_size ( ) const
inline

The total local number of points.

Returns
equivalent to local.size()

◆ lx()

template<class real_type >
real_type dg::aRealMPITopology3d< real_type >::lx ( ) const
inline

Return global lx.

Returns
global length

◆ ly()

template<class real_type >
real_type dg::aRealMPITopology3d< real_type >::ly ( ) const
inline

Return global ly.

Returns
global length

◆ lz()

template<class real_type >
real_type dg::aRealMPITopology3d< real_type >::lz ( ) const
inline

Return global lz.

Returns
global length

◆ multiplyCellNumbers()

template<class real_type >
void dg::aRealMPITopology3d< real_type >::multiplyCellNumbers ( real_type  fx,
real_type  fy 
)
inline

Multiply the number of cells with a given factor.

With this function you can resize the grid ignorantly of its current size

Parameters
fxnew global number of cells is fx*global().Nx()
fynew global number of cells is fy*global().Ny()

◆ n()

template<class real_type >
unsigned dg::aRealMPITopology3d< real_type >::n ( ) const
inline

Return n.

Returns
number of polynomial coefficients

◆ nx()

template<class real_type >
unsigned dg::aRealMPITopology3d< real_type >::nx ( ) const
inline

number of polynomial coefficients in x

◆ Nx()

template<class real_type >
unsigned dg::aRealMPITopology3d< real_type >::Nx ( ) const
inline

Return the global number of cells.

Returns
number of cells

◆ ny()

template<class real_type >
unsigned dg::aRealMPITopology3d< real_type >::ny ( ) const
inline

number of polynomial coefficients in y

◆ Ny()

template<class real_type >
unsigned dg::aRealMPITopology3d< real_type >::Ny ( ) const
inline

Return the global number of cells.

Returns
number of cells

◆ nz()

template<class real_type >
unsigned dg::aRealMPITopology3d< real_type >::nz ( ) const
inline

number of polynomial coefficients in z

◆ Nz()

template<class real_type >
unsigned dg::aRealMPITopology3d< real_type >::Nz ( ) const
inline

Return the global number of cells.

Returns
number of cells

◆ operator=()

template<class real_type >
aRealMPITopology3d & dg::aRealMPITopology3d< real_type >::operator= ( const aRealMPITopology3d< real_type > &  src)
protecteddefault

explicit assignment operator (default)

Parameters
srcsource

◆ pidOf()

template<class real_type >
int dg::aRealMPITopology3d< real_type >::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.

Parameters
xX-coord
yY-coord
zZ-coord
Returns
pid of a process, or -1 if non of the grids matches

◆ set() [1/2]

template<class real_type >
void dg::aRealMPITopology3d< real_type >::set ( unsigned  new_n,
unsigned  new_Nx,
unsigned  new_Ny,
unsigned  new_Nz 
)
inline

Set the number of polynomials and cells.

Set nz to 1 Same as set(new_n,new_Nx,new_n,new_Ny,1,new_Nz)

Parameters
new_nnew number of Gaussian nodes in x and y
new_Nxnew number of cells in x
new_Nynew number of cells in y
new_Nznew number of cells in z

◆ set() [2/2]

template<class real_type >
void dg::aRealMPITopology3d< real_type >::set ( unsigned  new_nx,
unsigned  new_Nx,
unsigned  new_ny,
unsigned  new_Ny,
unsigned  new_nz,
unsigned  new_Nz 
)
inline

Set the number of polynomials and cells.

Parameters
new_nxnew number of Gaussian nodes in x
new_Nxnew number of cells in x
new_nynew number of Gaussian nodes in y
new_Nynew number of cells in y
new_nznew number of Gaussian nodes in z
new_Nznew number of cells in z

◆ size()

template<class real_type >
unsigned dg::aRealMPITopology3d< real_type >::size ( ) const
inline

The total global number of points.

Returns
equivalent to nx()*ny()*nz()*Nx()*Ny()*Nz()

◆ x0()

template<class real_type >
real_type dg::aRealMPITopology3d< real_type >::x0 ( ) const
inline

Return global x0.

Returns
global left boundary

◆ x1()

template<class real_type >
real_type dg::aRealMPITopology3d< real_type >::x1 ( ) const
inline

Return global x1.

Returns
global right boundary

◆ y0()

template<class real_type >
real_type dg::aRealMPITopology3d< real_type >::y0 ( ) const
inline

Return global y0.

Returns
global left boundary

◆ y1()

template<class real_type >
real_type dg::aRealMPITopology3d< real_type >::y1 ( ) const
inline

Return global y1.

Returns
global right boundary

◆ z0()

template<class real_type >
real_type dg::aRealMPITopology3d< real_type >::z0 ( ) const
inline

Return global z0.

Returns
global left boundary

◆ z1()

template<class real_type >
real_type dg::aRealMPITopology3d< real_type >::z1 ( ) const
inline

Return global z1.

Returns
global right boundary

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