Discontinuous Galerkin Library
#include "dg/algorithm.h"
Loading...
Searching...
No Matches
dg::aRealMPIGeometry2d< real_type > Struct Template Referenceabstract

This is the abstract interface class for a two-dimensional Geometry. More...

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

Public Member Functions

SparseTensor< MPI_Vector< thrust::host_vector< real_type > > > jacobian () const
 The Jacobian of the coordinate transformation from physical to computational space.
 
SparseTensor< MPI_Vector< thrust::host_vector< real_type > > > metric () const
 The (inverse) metric tensor of the coordinate system.
 
std::vector< MPI_Vector< thrust::host_vector< real_type > > > map () const
 The coordinate map from computational to physical space.
 
virtual aRealMPIGeometry2dclone () const =0
 Geometries are cloneable.
 
virtual aRealGeometry2d< real_type > * global_geometry () const =0
 Construct the global non-MPI geometry.
 
virtual ~aRealMPIGeometry2d ()=default
 allow deletion through base class pointer
 
- Public Member Functions inherited from dg::aRealMPITopology< real_type, Nd >
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.
 

Protected Member Functions

 aRealMPIGeometry2d (const aRealMPIGeometry2d &src)=default
 
aRealMPIGeometry2doperator= (const aRealMPIGeometry2d &src)=default
 
- 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
 
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.
 

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.
 

Detailed Description

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

This is the abstract interface class for a two-dimensional Geometry.

Constructor & Destructor Documentation

◆ ~aRealMPIGeometry2d()

template<class real_type >
virtual dg::aRealMPIGeometry2d< real_type >::~aRealMPIGeometry2d ( )
virtualdefault

allow deletion through base class pointer

◆ aRealMPIGeometry2d()

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

explicit copy constructor (default)

Parameters
srcsource

Member Function Documentation

◆ clone()

template<class real_type >
virtual aRealMPIGeometry2d * dg::aRealMPIGeometry2d< real_type >::clone ( ) const
pure virtual

Geometries are cloneable.

Implemented in dg::RealCartesianMPIGrid2d< real_type >.

◆ global_geometry()

template<class real_type >
virtual aRealGeometry2d< real_type > * dg::aRealMPIGeometry2d< real_type >::global_geometry ( ) const
pure virtual

Construct the global non-MPI geometry.

Implemented in dg::RealCartesianMPIGrid2d< real_type >.

◆ jacobian()

template<class real_type >
SparseTensor< MPI_Vector< thrust::host_vector< real_type > > > dg::aRealMPIGeometry2d< real_type >::jacobian ( ) const
inline

The Jacobian of the coordinate transformation from physical to computational space.

The elements of the Tensor are (if x,y are the coordinates in computational space and R,Z are the physical space coordinates)

\[ J = \begin{pmatrix} x_R(x,y) & x_Z(x,y) \\ y_R(x,y) & y_Z(x,y) \end{pmatrix} \]

Returns
Jacobian
Note
per default this will be the identity tensor

◆ map()

template<class real_type >
std::vector< MPI_Vector< thrust::host_vector< real_type > > > dg::aRealMPIGeometry2d< real_type >::map ( ) const
inline

The coordinate map from computational to physical space.

The elements of the map are (if x,y are the coordinates in computational space and R,Z are the physical space coordinates)

\[ R(x,y) \\ Z(x,y) \]

Returns
a vector of size 2
Note
per default this will be the identity map

◆ metric()

template<class real_type >
SparseTensor< MPI_Vector< thrust::host_vector< real_type > > > dg::aRealMPIGeometry2d< real_type >::metric ( ) const
inline

The (inverse) metric tensor of the coordinate system.

The elements of the inverse metric tensor are the contravariant elements of the metric \(g\). If x,y are the coordinates, then

\[ g^{-1} = \begin{pmatrix} g^{xx}(x,y) & g^{xy}(x,y) \\ & g^{yy}(x,y) \end{pmatrix} \]

Returns
symmetric tensor
Note
use the dg::tensor::volume2d function to compute the volume element from here
per default this will be the identity tensor

◆ operator=()

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

explicit assignment operator (default)

Parameters
srcsource

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