Discontinuous Galerkin Library
#include "dg/algorithm.h"
dg::NestedGrids< Geometry, Matrix, Container > Struct Template Reference

Hold nested grids and provide dg fast interpolation and projection matrices. More...

Public Types

using geometry_type = Geometry
 
using matrix_type = Matrix
 
using container_type = Container
 
using value_type = get_value_type< Container >
 

Public Member Functions

 NestedGrids ()
 Allocate nothing, Call construct method before usage. More...
 
template<class ... ContainerParams>
 NestedGrids (const Geometry &grid, const unsigned stages, ContainerParams &&...ps)
 Construct the grids and the interpolation/projection operators. More...
 
template<class ... Params>
void construct (Params &&...ps)
 Perfect forward parameters to one of the constructors. More...
 
template<class ContainerType0 >
void project (const ContainerType0 &src, std::vector< ContainerType0 > &out) const
 Project vector to all involved grids. More...
 
template<class ContainerType0 >
std::vector< ContainerType0 > project (const ContainerType0 &src) const
 Project vector to all involved grids (allocate memory version) More...
 
const Container & copyable () const
 Return an object of same size as the object used for construction on the finest grid. More...
 
unsigned stages () const
 
unsigned num_stages () const
 
const Geometry & grid (unsigned stage) const
 return the grid at given stage More...
 
const MultiMatrix< Matrix, Container > & interpolation (unsigned stage) const
 return the interpolation matrix at given stage More...
 
const MultiMatrix< Matrix, Container > & projection (unsigned stage) const
 return the projection matrix at given stage More...
 
Container & x (unsigned stage)
 
const Container & x (unsigned stage) const
 
Container & r (unsigned stage)
 
const Container & r (unsigned stage) const
 
Container & b (unsigned stage)
 
const Container & b (unsigned stage) const
 
Container & w (unsigned stage)
 
const Container & w (unsigned stage) const
 

Detailed Description

template<class Geometry, class Matrix, class Container>
struct dg::NestedGrids< Geometry, Matrix, Container >

Hold nested grids and provide dg fast interpolation and projection matrices.

Template Parameters
GeometryA type that is or derives from one of the abstract geometry base classes ( aGeometry2d, aGeometry3d, aMPIGeometry2d, ...). Geometry determines which Matrix and Container types can be used:
MatrixA class for which the dg::blas2::symv functions are callable in connection with the Container class and to which the return type of dg::create::dx() can be converted using dg::blas2::transfer. The Matrix type can be one of:
ContainerA data container class for which the blas1 functionality is overloaded and to which the return type of blas1::subroutine() can be converted using dg::assign. We assume that Container is copyable/assignable and has a swap member function. In connection with Geometry this is one of

Member Typedef Documentation

◆ container_type

template<class Geometry , class Matrix , class Container >
using dg::NestedGrids< Geometry, Matrix, Container >::container_type = Container

◆ geometry_type

template<class Geometry , class Matrix , class Container >
using dg::NestedGrids< Geometry, Matrix, Container >::geometry_type = Geometry

◆ matrix_type

template<class Geometry , class Matrix , class Container >
using dg::NestedGrids< Geometry, Matrix, Container >::matrix_type = Matrix

◆ value_type

template<class Geometry , class Matrix , class Container >
using dg::NestedGrids< Geometry, Matrix, Container >::value_type = get_value_type<Container>

Constructor & Destructor Documentation

◆ NestedGrids() [1/2]

template<class Geometry , class Matrix , class Container >
dg::NestedGrids< Geometry, Matrix, Container >::NestedGrids ( )
inline

Allocate nothing, Call construct method before usage.

◆ NestedGrids() [2/2]

template<class Geometry , class Matrix , class Container >
template<class ... ContainerParams>
dg::NestedGrids< Geometry, Matrix, Container >::NestedGrids ( const Geometry &  grid,
const unsigned  stages,
ContainerParams &&...  ps 
)
inline

Construct the grids and the interpolation/projection operators.

Parameters
gridthe original grid (Nx() and Ny() must be evenly divisable by pow(2, stages-1)
stagesnumber of grids in total (The second grid contains half the points of the original grids, The third grid contains half of the second grid ...). Must be >= 1
psparameters necessary for dg::construct to construct a Container from a dg::HVec

Member Function Documentation

◆ b() [1/2]

template<class Geometry , class Matrix , class Container >
Container & dg::NestedGrids< Geometry, Matrix, Container >::b ( unsigned  stage)
inline

◆ b() [2/2]

template<class Geometry , class Matrix , class Container >
const Container & dg::NestedGrids< Geometry, Matrix, Container >::b ( unsigned  stage) const
inline

◆ construct()

template<class Geometry , class Matrix , class Container >
template<class ... Params>
void dg::NestedGrids< Geometry, Matrix, Container >::construct ( Params &&...  ps)
inline

Perfect forward parameters to one of the constructors.

Template Parameters
Paramsdeduced by the compiler
Parameters
psparameters forwarded to constructors

◆ copyable()

template<class Geometry , class Matrix , class Container >
const Container & dg::NestedGrids< Geometry, Matrix, Container >::copyable ( ) const
inline

Return an object of same size as the object used for construction on the finest grid.

Returns
A copyable object; what it contains is undefined, its size is important

◆ grid()

template<class Geometry , class Matrix , class Container >
const Geometry & dg::NestedGrids< Geometry, Matrix, Container >::grid ( unsigned  stage) const
inline

return the grid at given stage

Parameters
stagemust fulfill 0 <= stage < stages()

◆ interpolation()

template<class Geometry , class Matrix , class Container >
const MultiMatrix< Matrix, Container > & dg::NestedGrids< Geometry, Matrix, Container >::interpolation ( unsigned  stage) const
inline

return the interpolation matrix at given stage

Parameters
stagemust fulfill 0 <= stage < stages()-1

◆ num_stages()

template<class Geometry , class Matrix , class Container >
unsigned dg::NestedGrids< Geometry, Matrix, Container >::num_stages ( ) const
inline
Returns
number of stages (same as stages)

◆ project() [1/2]

template<class Geometry , class Matrix , class Container >
template<class ContainerType0 >
std::vector< ContainerType0 > dg::NestedGrids< Geometry, Matrix, Container >::project ( const ContainerType0 &  src) const
inline

Project vector to all involved grids (allocate memory version)

Parameters
srcthe input vector
Returns
the input vector projected to all grids ( index 0 contains a copy of src, 1 is the projetion to the first coarse grid, 2 is the next coarser grid, ...)

◆ project() [2/2]

template<class Geometry , class Matrix , class Container >
template<class ContainerType0 >
void dg::NestedGrids< Geometry, Matrix, Container >::project ( const ContainerType0 &  src,
std::vector< ContainerType0 > &  out 
) const
inline

Project vector to all involved grids.

Parameters
srcthe input vector (may alias first element of out)
outthe input vector projected to all grids ( index 0 contains a copy of src, 1 is the projetion to the first coarse grid, 2 is the next coarser grid, ...)
Note
out is not resized

◆ projection()

template<class Geometry , class Matrix , class Container >
const MultiMatrix< Matrix, Container > & dg::NestedGrids< Geometry, Matrix, Container >::projection ( unsigned  stage) const
inline

return the projection matrix at given stage

Parameters
stagemust fulfill 0 <= stage < stages()-1

◆ r() [1/2]

template<class Geometry , class Matrix , class Container >
Container & dg::NestedGrids< Geometry, Matrix, Container >::r ( unsigned  stage)
inline

◆ r() [2/2]

template<class Geometry , class Matrix , class Container >
const Container & dg::NestedGrids< Geometry, Matrix, Container >::r ( unsigned  stage) const
inline

◆ stages()

template<class Geometry , class Matrix , class Container >
unsigned dg::NestedGrids< Geometry, Matrix, Container >::stages ( ) const
inline
Returns
number of stages (same as num_stages)

◆ w() [1/2]

template<class Geometry , class Matrix , class Container >
Container & dg::NestedGrids< Geometry, Matrix, Container >::w ( unsigned  stage)
inline

◆ w() [2/2]

template<class Geometry , class Matrix , class Container >
const Container & dg::NestedGrids< Geometry, Matrix, Container >::w ( unsigned  stage) const
inline

◆ x() [1/2]

template<class Geometry , class Matrix , class Container >
Container & dg::NestedGrids< Geometry, Matrix, Container >::x ( unsigned  stage)
inline

◆ x() [2/2]

template<class Geometry , class Matrix , class Container >
const Container & dg::NestedGrids< Geometry, Matrix, Container >::x ( unsigned  stage) const
inline

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