AMDiS  0.1
The Adaptive Multi-Dimensional Simulation Toolbox
GridFunctionOperatorBase< Derived, LC, GF > Class Template Reference

The main implementation of an CRTP-base class for operators using a grid-function coefficient to be used in an Assembler. More...

#include <GridFunctionOperator.hpp>

Inherits LocalOperator< Derived, LC >.

Public Member Functions

template<class GridFct >
 GridFunctionOperatorBase (GridFct &&gridFct, int termOrder)
 Constructor. Stores a copy of gridFct. More...
 
template<class PQF >
void setQuadFactory (PQF &&pre)
 
- Public Member Functions inherited from LocalOperator< Derived, LC >
template<class GridView >
void init (GridView const &gridView)
 Initialize the local operator on the current gridView.
 
void bind (Element const &element, Geometry const &geometry)
 Binds operator to element and geometry. More...
 
void unbind ()
 Unbinds operator from element.
 
template<class CG , class RN , class CN , class Mat >
void calculateElementMatrix (CG const &contextGeo, RN const &rowNode, CN const &colNode, Mat &elementMatrix)
 Assemble a local element matrix on the element that is bound. More...
 
template<class CG , class Node , class Vec >
void calculateElementVector (CG const &contextGeo, Node const &node, Vec &elementVector)
 Assemble a local element vector on the element that is bound. More...
 
Derived & derived ()
 
Derived const & derived () const
 

Protected Member Functions

auto coefficient (LocalCoordinate const &local) const
 Return expression value at LocalCoordinates.
 
template<class... Nodes>
auto const & getQuadratureRule (Dune::GeometryType type, Nodes const &... nodes) const
 
- Protected Member Functions inherited from LocalOperator< Derived, LC >
template<class GridView >
void init_impl (GridView const &)
 
template<class Element , class Geometry >
void bind_impl (Element const &, Geometry const &)
 
void unbind_impl ()
 
template<class RN , class CN >
int getDegree (int derivOrder, int coeffDegree, RN const &rowNode, CN const &colNode) const
 Return the quadrature degree for a matrix operator. More...
 
template<class Node >
int getDegree (int derivOrder, int coeffDegree, Node const &node) const
 Return the quadrature degree for a vector operator. More...
 

Friends

template<class , class >
class LocalOperator
 

Additional Inherited Members

- Public Types inherited from LocalOperator< Derived, LC >
using LocalContext = LC
 The element or intersection the operator is assembled on.
 
using Element = typename ContextType::Entity
 The codim=0 grid entity.
 
using Geometry = typename Element::Geometry
 The geometry of the Element.
 
- Protected Attributes inherited from LocalOperator< Derived, LC >
bool isSimplex_ = false
 
bool isAffine_ = false
 
bool bound_ = false
 

Detailed Description

template<class Derived, class LC, class GF>
class AMDiS::GridFunctionOperatorBase< Derived, LC, GF >

The main implementation of an CRTP-base class for operators using a grid-function coefficient to be used in an Assembler.

An Operator that takes a GridFunction as coefficient. Provides quadrature rules and handles the evaluation of the GridFunction at local coordinates.

The class is specialized, by deriving from it, in GridFunctionOperator.

Template Parameters
DerivedThe class derived from GridFunctionOperatorBase
LCThe Element or Intersection type
GFThe GridFunction, a LocalFunction is created from, and that is evaluated at quadrature points.

Requirements:

  • LC models either Entity (of codim 0) or Intersection
  • GF models the Concepts::GridFunction

Constructor & Destructor Documentation

◆ GridFunctionOperatorBase()

GridFunctionOperatorBase ( GridFct &&  gridFct,
int  termOrder 
)
inline

Constructor. Stores a copy of gridFct.

A GridFunctionOperator takes a gridFunction and the differentiation order of the operator, to calculate the quadrature degree in getDegree.

Member Function Documentation

◆ getQuadratureRule()

auto const& getQuadratureRule ( Dune::GeometryType  type,
Nodes const &...  nodes 
) const
inlineprotected

Create a quadrature rule using the QuadratureFactory by calculating the quadrature order necessary to integrate the (bi)linear-form accurately.

◆ setQuadFactory()

void setQuadFactory ( PQF &&  pre)
inline

Create a quadrature factory from a PreQuadratureFactory, e.g. class derived from QuadratureFactory

Template Parameters
PQFA PreQuadratureFactory

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