3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
Modules | Files | Classes | Typedefs | Functions
Discretization schemes

The discretization schemes available in DuMux More...

Description

The discretization schemes available in DuMux

Modules

 Box FV scheme
 The box method is a collocated finite volume scheme with control volumes centered at grid nodes.
 
 Cell-centered FV scheme
 Finite volume schemes with degrees of freedom located at grid cell centers.
 
 CVFE
 Control-volume finite element schemes (e.g. box method) Control-volume finite element schemes are based on finite element basis functions for interpolation but define control volumes to construct a finite volume scheme. They can be interpreted both as finite volume or as (Petrov-Galerkin) finite element scheme.
 
 Face-centered staggered FV scheme
 A staggered finite volume scheme with degrees of freedom at cell-centers and facets. In this implementation, momentum control volumes exist.
 
 Staggered FV scheme
 A staggered finite volume scheme with degrees of freedom at cell-centers and facets. In this implementation, momentum control volumes do not explicitly exist, but the implementation uses workarounds.
 
 Finite element method
 The finite element method.
 
 Pore network model discretization
 The pore-network model discretization.
 
 PQ1 bubble scheme
 Control-volume finite element scheme based on P1/Q1 basis function with enrichment by a bubble function.
 

Files

file  basegridgeometry.hh
 Base class for grid geometries.
 
file  basicgridgeometry.hh
 A basic implementation of a grid geometry with some common interfaces.
 
file  box.hh
 Defines a type tag and some properties for models using the box scheme.
 
file  ccmpfa.hh
 Properties for all models using cell-centered finite volume scheme with mpfa.
 
file  cctpfa.hh
 Properties for all models using cell-centered finite volume scheme with TPFA.
 
file  checkoverlapsize.hh
 Check the overlap size for different discretization methods.
 
file  elementsolution.hh
 Element solution classes and factory functions.
 
file  evalgradients.hh
 free functions for the evaluation of primary variable gradients inside elements.
 
file  evalsolution.hh
 free functions for the evaluation of primary variables inside elements.
 
file  extrusion.hh
 Helper classes to compute the integration elements.
 
file  fcdiamond.hh
 Defines a type tag and some properties for models using the diamond scheme. This scheme features degrees of freedom at the elements' centers and intersections (faces).
 
file  fcstaggered.hh
 Defines a type tag and some properties for models using the staggered scheme. This scheme features degrees of freedom at the elements' centers and intersections (faces). TODO: detailed documentation and figures.
 
file  fluxstencil.hh
 The flux stencil specialized for different discretization schemes.
 
file  functionspacebasis.hh
 Provides helper aliases and functionality to obtain the types and instances of Dune::Functions function space bases that underlie different discretization schemes.
 
file  discretization/fvgridvariables.hh
 The grid variable class for finite volume schemes, storing variables on scv and scvf (volume and flux variables)
 
file  fvproperties.hh
 Declares properties required for finite-volume models models.
 
file  localview.hh
 Free function to get the local view of a grid cache object.
 
file  method.hh
 The available discretization methods in Dumux.
 
file  nonconformingfecache.hh
 A finite element cache for the non-conforming FE spaces RT and CR.
 
file  pq1bubble.hh
 Defines a type tag and some properties for models using the pq1bubble scheme.
 
file  projector.hh
 Contains functionality for L2-projections from one function space into another, which can live both on the same or different grids of potentially different dimensionality.
 
file  scvandscvfiterators.hh
 Class providing iterators over sub control volumes and sub control volume faces of an element.
 
file  staggered.hh
 Defines a type tag and some properties for models using the staggered scheme. This scheme features degrees of freedom at the elements' centers and intersections (faces). TODO: detailed documentation and figures.
 
file  subcontrolvolumebase.hh
 Base class for a sub control volume.
 
file  subcontrolvolumefacebase.hh
 Base class for a sub control volume face.
 
file  walldistance.hh
 
file  experimental/discretization/fvgridvariables.hh
 The grid variable class for finite volume schemes, storing variables on scv and scvf (volume and flux variables)
 
file  experimental/discretization/gridvariables.hh
 Base class for grid variables.
 

Classes

class  Dumux::BaseGridGeometry< GV, Traits >
 Base class for all grid geometries. More...
 
class  Dumux::BasicGridGeometry< GV, EM, VM >
 An implementation of a grid geometry with some basic features. More...
 
struct  Dumux::CheckOverlapSize< DiscretizationMethod >
 Check if the overlap size is valid for a given discretization method. More...
 
struct  Dumux::NoExtrusion
 Default implementation that performs no extrusion (extrusion with identity) More...
 
struct  Dumux::RotationalExtrusion< radAx >
 Rotation symmetric extrusion policy for rotating about an external axis. More...
 
struct  Dumux::SphericalExtrusion
 Rotation symmetric extrusion policy for spherical rotation. More...
 
class  Dumux::FluxStencil< FVElementGeometry, DiscretizationMethod >
 The flux stencil specialized for different discretization schemes. More...
 
class  Dumux::FVGridVariables< GG, GVV, GFVC >
 The grid variable class for finite volume schemes storing variables on scv and scvf (volume and flux variables) More...
 
class  Dumux::Projector< ScalarType >
 Does an L2-projection from one discrete function space into another. The convenience functions makeProjectorPair or makeProjector can be used to create such a projection. More...
 
class  Dumux::ScvIterator< SubControlVolume, Vector, FVElementGeometry >
 Iterators over sub control volumes. More...
 
class  Dumux::ScvfIterator< SubControlVolumeFace, Vector, FVElementGeometry >
 Iterators over sub control volume faces of an fv geometry. More...
 
class  Dumux::SkippingScvfIterator< SubControlVolumeFace, Vector, FVElementGeometry >
 Iterators over sub control volume faces of an fv geometry and a given sub control volume. More...
 
class  Dumux::SubControlVolumeBase< Imp, ScvGeometryTraits >
 Base class for a sub control volume, i.e a part of the control volume we are making the balance for. Defines the general interface. More...
 
class  Dumux::SubControlVolumeFaceBase< Imp, ScvfGeometryTraits >
 Base class for a sub control volume face, i.e a part of the boundary of a sub control volume we computing a flux on. More...
 
class  Dumux::WallDistance< GridGeometry, DistanceField >
 Class to calculate the wall distance at every element or vertex of a grid. More...
 
class  Dumux::Experimental::FVGridVariablesLocalView< GV >
 Finite volume-specific local view on grid variables. More...
 
class  Dumux::Experimental::FVGridVariables< GVV, GFVC, X >
 The grid variable class for finite volume schemes, storing variables on scv and scvf (volume and flux variables). More...
 
class  Dumux::Experimental::GridVariables< GG, X >
 Base class for grid variables. More...
 

Typedefs

template<class GV , class T >
using Dumux::BasicGridGeometry_t = Dune::Std::detected_or_t< Dumux::BasicGridGeometry< GV, typename T::ElementMapper, typename T::VertexMapper >, Detail::SpecifiesBaseGridGeometry, T >
 Type of the basic grid geometry implementation used as backend. More...
 

Functions

 Dumux::BaseGridGeometry< GV, Traits >::BaseGridGeometry (std::shared_ptr< BaseImplementation > impl)
 Constructor from a BaseImplementation. More...
 
 Dumux::BaseGridGeometry< GV, Traits >::BaseGridGeometry (const GridView &gridView)
 Constructor from a grid view. More...
 
 Dumux::BasicGridGeometry< GV, EM, VM >::BasicGridGeometry (const GridView &gridView)
 Constructor computes the bounding box of the entire domain, for e.g. setting boundary conditions. More...
 
template<class Element , class GridGeometry , class CVFEElemSol >
auto Dumux::Detail::evalCVFEGradients (const Element &element, const typename Element::Geometry &geometry, const GridGeometry &gridGeometry, const CVFEElemSol &elemSol, const typename Element::Geometry::GlobalCoordinate &globalPos, bool ignoreState=false)
 Evaluates the gradient of a control-volume finite element solution to a given global position. More...
 
template<class Element , class FVElementGeometry , class PrimaryVariables >
auto Dumux::evalGradients (const Element &element, const typename Element::Geometry &geometry, const typename FVElementGeometry::GridGeometry &gridGeometry, const BoxElementSolution< FVElementGeometry, PrimaryVariables > &elemSol, const typename Element::Geometry::GlobalCoordinate &globalPos, bool ignoreState=false)
 Evaluates the gradient of a given box element solution to a given global position. More...
 
template<class Element , class FVElementGeometry , class PrimaryVariables >
auto Dumux::evalGradients (const Element &element, const typename Element::Geometry &geometry, const typename FVElementGeometry::GridGeometry &gridGeometry, const FaceCenteredDiamondElementSolution< FVElementGeometry, PrimaryVariables > &elemSol, const typename Element::Geometry::GlobalCoordinate &globalPos, bool ignoreState=false)
 Evaluates the gradient of a given diamond scheme element solution to a given global position. More...
 
template<class Element , class FVElementGeometry , class PrimaryVariables >
auto Dumux::evalGradients (const Element &element, const typename Element::Geometry &geometry, const typename FVElementGeometry::GridGeometry &gridGeometry, const PQ1BubbleElementSolution< FVElementGeometry, PrimaryVariables > &elemSol, const typename Element::Geometry::GlobalCoordinate &globalPos, bool ignoreState=false)
 Evaluates the gradient of a given pq1bubble scheme element solution to a given global position. More...
 
template<class Element , class FVElementGeometry , class PrimaryVariables >
Dune::FieldVector< typename Element::Geometry::GlobalCoordinate, PrimaryVariables::dimension > Dumux::evalGradients (const Element &element, const typename Element::Geometry &geometry, const typename FVElementGeometry::GridGeometry &gridGeometry, const CCElementSolution< FVElementGeometry, PrimaryVariables > &elemSol, const typename Element::Geometry::GlobalCoordinate &globalPos)
 Evaluates the gradient of a given CCElementSolution to a given global position. This function is only here for (compilation) compatibility reasons with the box scheme. The solution within the control volumes is constant and thus gradients are zero. One can compute gradients towards the sub-control volume faces after reconstructing the solution on the faces. More...
 
template<class Element , class GridGeometry , class CVFEElemSol >
CVFEElemSol::PrimaryVariables Dumux::Detail::evalCVFESolution (const Element &element, const typename Element::Geometry &geometry, const GridGeometry &gridGeometry, const CVFEElemSol &elemSol, const typename Element::Geometry::GlobalCoordinate &globalPos, bool ignoreState=false)
 Interpolates a given control-volume finite element solution at a given global position. Uses the finite element cache of the grid geometry. More...
 
template<class Element , class FVElementGeometry , class PrimaryVariables >
PrimaryVariables Dumux::evalSolution (const Element &element, const typename Element::Geometry &geometry, const typename FVElementGeometry::GridGeometry &gridGeometry, const BoxElementSolution< FVElementGeometry, PrimaryVariables > &elemSol, const typename Element::Geometry::GlobalCoordinate &globalPos, bool ignoreState=false)
 Interpolates a given box element solution at a given global position. Uses the finite element cache of the grid geometry. More...
 
template<class Element , class FVElementGeometry , class PrimaryVariables >
PrimaryVariables Dumux::evalSolution (const Element &element, const typename Element::Geometry &geometry, const BoxElementSolution< FVElementGeometry, PrimaryVariables > &elemSol, const typename Element::Geometry::GlobalCoordinate &globalPos, bool ignoreState=false)
 Interpolates a given box element solution at a given global position. More...
 
template<class Element , class FVElementGeometry , class PrimaryVariables >
PrimaryVariables Dumux::evalSolution (const Element &element, const typename Element::Geometry &geometry, const typename FVElementGeometry::GridGeometry &gridGeometry, const CCElementSolution< FVElementGeometry, PrimaryVariables > &elemSol, const typename Element::Geometry::GlobalCoordinate &globalPos, bool ignoreState=false)
 Interpolates a given cell-centered element solution at a given global position. More...
 
template<class Element , class FVElementGeometry , class PrimaryVariables >
PrimaryVariables Dumux::evalSolution (const Element &element, const typename Element::Geometry &geometry, const CCElementSolution< FVElementGeometry, PrimaryVariables > &elemSol, const typename Element::Geometry::GlobalCoordinate &globalPos, bool ignoreState=false)
 Interpolates a given cell-centered element solution at a given global position. Overload of the above evalSolution() function without a given gridGeometry. For compatibility reasons with the box scheme. More...
 
template<class Element , class FVElementGeometry , class PrimaryVariables >
PrimaryVariables Dumux::evalSolution (const Element &element, const typename Element::Geometry &geometry, const typename FVElementGeometry::GridGeometry &gridGeometry, const FaceCenteredDiamondElementSolution< FVElementGeometry, PrimaryVariables > &elemSol, const typename Element::Geometry::GlobalCoordinate &globalPos, bool ignoreState=false)
 Interpolates a given diamond scheme element solution at a given global position. Uses the finite element cache of the grid geometry. More...
 
template<class Element , class FVElementGeometry , class PrimaryVariables >
PrimaryVariables Dumux::evalSolution (const Element &element, const typename Element::Geometry &geometry, const typename FVElementGeometry::GridGeometry &gridGeometry, const PQ1BubbleElementSolution< FVElementGeometry, PrimaryVariables > &elemSol, const typename Element::Geometry::GlobalCoordinate &globalPos, bool ignoreState=false)
 Interpolates a given pq1bubble scheme element solution at a given global position. Uses the finite element cache of the grid geometry. More...
 
template<class GridCache >
GridCache::LocalView Dumux::localView (const GridCache &gridCache)
 Free function to get the local view of a grid cache object. More...
 
template<class FEBasisDomain , class FEBasisTarget , class GlueType >
auto Dumux::makeProjectorPair (const FEBasisDomain &feBasisDomain, const FEBasisTarget &feBasisTarget, GlueType glue)
 Creates a pair of projectors between the space with basis feBasisDomain to the space with basis feBasisTarget. More...
 
template<class FEBasisDomain , class FEBasisTarget , class GlueType >
auto Dumux::makeProjector (const FEBasisDomain &feBasisDomain, const FEBasisTarget &feBasisTarget, GlueType glue)
 Creates a forward projector from the space feBasisDomain to the space with basis feBasisTarget. More...
 

Typedef Documentation

◆ BasicGridGeometry_t

template<class GV , class T >
using Dumux::BasicGridGeometry_t = typedef Dune::Std::detected_or_t< Dumux::BasicGridGeometry<GV, typename T::ElementMapper, typename T::VertexMapper>, Detail::SpecifiesBaseGridGeometry, T >

Type of the basic grid geometry implementation used as backend.

Function Documentation

◆ BaseGridGeometry() [1/2]

template<class GV , class Traits >
Dumux::BaseGridGeometry< GV, Traits >::BaseGridGeometry ( const GridView gridView)
inline

Constructor from a grid view.

Parameters
gridViewthe grid view on which to construct the grid geometry

◆ BaseGridGeometry() [2/2]

template<class GV , class Traits >
Dumux::BaseGridGeometry< GV, Traits >::BaseGridGeometry ( std::shared_ptr< BaseImplementation >  impl)
inline

Constructor from a BaseImplementation.

◆ BasicGridGeometry()

template<class GV , class EM , class VM >
Dumux::BasicGridGeometry< GV, EM, VM >::BasicGridGeometry ( const GridView gridView)
inline

Constructor computes the bounding box of the entire domain, for e.g. setting boundary conditions.

Parameters
gridViewthe grid view on which to construct the grid geometry

◆ evalCVFEGradients()

template<class Element , class GridGeometry , class CVFEElemSol >
auto Dumux::Detail::evalCVFEGradients ( const Element &  element,
const typename Element::Geometry &  geometry,
const GridGeometry &  gridGeometry,
const CVFEElemSol &  elemSol,
const typename Element::Geometry::GlobalCoordinate &  globalPos,
bool  ignoreState = false 
)

Evaluates the gradient of a control-volume finite element solution to a given global position.

Parameters
elementThe element
geometryThe element geometry
gridGeometryThe finite volume grid geometry
elemSolThe primary variables at the dofs of the element
globalPosThe global position
ignoreStateIf true, the state of primary variables is ignored
Returns
Dune::FieldVector with as many entries as dimension of the PrimaryVariables object (i.e. numEq). Each entry is a GlobalCoordinate object holding the priVar gradient.

◆ evalCVFESolution()

template<class Element , class GridGeometry , class CVFEElemSol >
CVFEElemSol::PrimaryVariables Dumux::Detail::evalCVFESolution ( const Element &  element,
const typename Element::Geometry &  geometry,
const GridGeometry &  gridGeometry,
const CVFEElemSol &  elemSol,
const typename Element::Geometry::GlobalCoordinate &  globalPos,
bool  ignoreState = false 
)

Interpolates a given control-volume finite element solution at a given global position. Uses the finite element cache of the grid geometry.

Returns
the interpolated primary variables
Parameters
elementThe element
geometryThe element geometry
gridGeometryThe finite volume grid geometry
elemSolThe primary variables at the dofs of the element
globalPosThe global position
ignoreStateIf true, the state of primary variables is ignored

◆ evalGradients() [1/4]

template<class Element , class FVElementGeometry , class PrimaryVariables >
auto Dumux::evalGradients ( const Element &  element,
const typename Element::Geometry &  geometry,
const typename FVElementGeometry::GridGeometry &  gridGeometry,
const BoxElementSolution< FVElementGeometry, PrimaryVariables > &  elemSol,
const typename Element::Geometry::GlobalCoordinate &  globalPos,
bool  ignoreState = false 
)

Evaluates the gradient of a given box element solution to a given global position.

Parameters
elementThe element
geometryThe element geometry
gridGeometryThe finite volume grid geometry
elemSolThe primary variables at the dofs of the element
globalPosThe global position
ignoreStateIf true, the state of primary variables is ignored

◆ evalGradients() [2/4]

template<class Element , class FVElementGeometry , class PrimaryVariables >
Dune::FieldVector< typename Element::Geometry::GlobalCoordinate, PrimaryVariables::dimension > Dumux::evalGradients ( const Element &  element,
const typename Element::Geometry &  geometry,
const typename FVElementGeometry::GridGeometry &  gridGeometry,
const CCElementSolution< FVElementGeometry, PrimaryVariables > &  elemSol,
const typename Element::Geometry::GlobalCoordinate &  globalPos 
)

Evaluates the gradient of a given CCElementSolution to a given global position. This function is only here for (compilation) compatibility reasons with the box scheme. The solution within the control volumes is constant and thus gradients are zero. One can compute gradients towards the sub-control volume faces after reconstructing the solution on the faces.

Parameters
elementThe element
geometryThe element geometry
gridGeometryThe finite volume grid geometry
elemSolThe primary variables at the dofs of the element
globalPosThe global position
Exceptions
Dune::NotImplemented
Returns
Dune::FieldVector with as many entries as dimension of the PrimaryVariables object (i.e. numEq). Each entry is a GlobalCoordinate object holding the priVar gradient.

◆ evalGradients() [3/4]

template<class Element , class FVElementGeometry , class PrimaryVariables >
auto Dumux::evalGradients ( const Element &  element,
const typename Element::Geometry &  geometry,
const typename FVElementGeometry::GridGeometry &  gridGeometry,
const FaceCenteredDiamondElementSolution< FVElementGeometry, PrimaryVariables > &  elemSol,
const typename Element::Geometry::GlobalCoordinate &  globalPos,
bool  ignoreState = false 
)

Evaluates the gradient of a given diamond scheme element solution to a given global position.

Parameters
elementThe element
geometryThe element geometry
gridGeometryThe finite volume grid geometry
elemSolThe primary variables at the dofs of the element
globalPosThe global position
ignoreStateIf true, the state of primary variables is ignored

◆ evalGradients() [4/4]

template<class Element , class FVElementGeometry , class PrimaryVariables >
auto Dumux::evalGradients ( const Element &  element,
const typename Element::Geometry &  geometry,
const typename FVElementGeometry::GridGeometry &  gridGeometry,
const PQ1BubbleElementSolution< FVElementGeometry, PrimaryVariables > &  elemSol,
const typename Element::Geometry::GlobalCoordinate &  globalPos,
bool  ignoreState = false 
)

Evaluates the gradient of a given pq1bubble scheme element solution to a given global position.

Parameters
elementThe element
geometryThe element geometry
gridGeometryThe finite volume grid geometry
elemSolThe primary variables at the dofs of the element
globalPosThe global position
ignoreStateIf true, the state of primary variables is ignored

◆ evalSolution() [1/6]

template<class Element , class FVElementGeometry , class PrimaryVariables >
PrimaryVariables Dumux::evalSolution ( const Element &  element,
const typename Element::Geometry &  geometry,
const BoxElementSolution< FVElementGeometry, PrimaryVariables > &  elemSol,
const typename Element::Geometry::GlobalCoordinate &  globalPos,
bool  ignoreState = false 
)

Interpolates a given box element solution at a given global position.

Overload of the above evalSolution() function without a given gridGeometry. The local basis is computed on the fly.

Returns
the interpolated primary variables
Parameters
elementThe element
geometryThe element geometry
elemSolThe primary variables at the dofs of the element
globalPosThe global position
ignoreStateIf true, the state of primary variables is ignored

◆ evalSolution() [2/6]

template<class Element , class FVElementGeometry , class PrimaryVariables >
PrimaryVariables Dumux::evalSolution ( const Element &  element,
const typename Element::Geometry &  geometry,
const CCElementSolution< FVElementGeometry, PrimaryVariables > &  elemSol,
const typename Element::Geometry::GlobalCoordinate &  globalPos,
bool  ignoreState = false 
)

Interpolates a given cell-centered element solution at a given global position. Overload of the above evalSolution() function without a given gridGeometry. For compatibility reasons with the box scheme.

Returns
the primary variables (constant over the element)
Parameters
elementThe element
geometryThe element geometry
elemSolThe primary variables at the dofs of the element
globalPosThe global position
ignoreStateIf true, the state of primary variables is ignored

◆ evalSolution() [3/6]

template<class Element , class FVElementGeometry , class PrimaryVariables >
PrimaryVariables Dumux::evalSolution ( const Element &  element,
const typename Element::Geometry &  geometry,
const typename FVElementGeometry::GridGeometry &  gridGeometry,
const BoxElementSolution< FVElementGeometry, PrimaryVariables > &  elemSol,
const typename Element::Geometry::GlobalCoordinate &  globalPos,
bool  ignoreState = false 
)

Interpolates a given box element solution at a given global position. Uses the finite element cache of the grid geometry.

Returns
the interpolated primary variables
Parameters
elementThe element
geometryThe element geometry
gridGeometryThe finite volume grid geometry
elemSolThe primary variables at the dofs of the element
globalPosThe global position
ignoreStateIf true, the state of primary variables is ignored

◆ evalSolution() [4/6]

template<class Element , class FVElementGeometry , class PrimaryVariables >
PrimaryVariables Dumux::evalSolution ( const Element &  element,
const typename Element::Geometry &  geometry,
const typename FVElementGeometry::GridGeometry &  gridGeometry,
const CCElementSolution< FVElementGeometry, PrimaryVariables > &  elemSol,
const typename Element::Geometry::GlobalCoordinate &  globalPos,
bool  ignoreState = false 
)

Interpolates a given cell-centered element solution at a given global position.

Returns
the primary variables (constant over the element)
Parameters
elementThe element
geometryThe element geometry
gridGeometryThe finite volume grid geometry
elemSolThe primary variables at the dofs of the element
globalPosThe global position
ignoreStateIf true, the state of primary variables is ignored

◆ evalSolution() [5/6]

template<class Element , class FVElementGeometry , class PrimaryVariables >
PrimaryVariables Dumux::evalSolution ( const Element &  element,
const typename Element::Geometry &  geometry,
const typename FVElementGeometry::GridGeometry &  gridGeometry,
const FaceCenteredDiamondElementSolution< FVElementGeometry, PrimaryVariables > &  elemSol,
const typename Element::Geometry::GlobalCoordinate &  globalPos,
bool  ignoreState = false 
)

Interpolates a given diamond scheme element solution at a given global position. Uses the finite element cache of the grid geometry.

Returns
the interpolated primary variables
Parameters
elementThe element
geometryThe element geometry
gridGeometryThe finite volume grid geometry
elemSolThe primary variables at the dofs of the element
globalPosThe global position
ignoreStateIf true, the state of primary variables is ignored

◆ evalSolution() [6/6]

template<class Element , class FVElementGeometry , class PrimaryVariables >
PrimaryVariables Dumux::evalSolution ( const Element &  element,
const typename Element::Geometry &  geometry,
const typename FVElementGeometry::GridGeometry &  gridGeometry,
const PQ1BubbleElementSolution< FVElementGeometry, PrimaryVariables > &  elemSol,
const typename Element::Geometry::GlobalCoordinate &  globalPos,
bool  ignoreState = false 
)

Interpolates a given pq1bubble scheme element solution at a given global position. Uses the finite element cache of the grid geometry.

Returns
the interpolated primary variables
Parameters
elementThe element
geometryThe element geometry
gridGeometryThe finite volume grid geometry
elemSolThe primary variables at the dofs of the element
globalPosThe global position
ignoreStateIf true, the state of primary variables is ignored

◆ localView()

template<class GridCache >
GridCache::LocalView Dumux::localView ( const GridCache &  gridCache)
inline

Free function to get the local view of a grid cache object.

Note
A local object is only functional after calling its bind/bindElement method.
Template Parameters
GridCachethe grid caching type (such as GridGeometry)
Parameters
gridCachethe grid caching object we want to localView from

◆ makeProjector()

template<class FEBasisDomain , class FEBasisTarget , class GlueType >
auto Dumux::makeProjector ( const FEBasisDomain &  feBasisDomain,
const FEBasisTarget &  feBasisTarget,
GlueType  glue 
)

Creates a forward projector from the space feBasisDomain to the space with basis feBasisTarget.

Parameters
feBasisDomainThe domain finite element space basis
feBasisTargetThe target finite element space basis
glueThe glue object containing the intersections between the grids.
Returns
The forward projector from the space with basis feBasisDomain to the space with basis feBasisTarget.

◆ makeProjectorPair()

template<class FEBasisDomain , class FEBasisTarget , class GlueType >
auto Dumux::makeProjectorPair ( const FEBasisDomain &  feBasisDomain,
const FEBasisTarget &  feBasisTarget,
GlueType  glue 
)

Creates a pair of projectors between the space with basis feBasisDomain to the space with basis feBasisTarget.

Parameters
feBasisDomainThe domain finite element space basis
feBasisTargetThe target finite element space basis
glueThe glue object containing the intersections between the grids.
Returns
An std::pair of projectors where the first is the forward projector from the space with basis feBasisDomain to the space with basis feBasisTarget and the second does the backward projection.