DUNE PDELab (git)

A pre-basis for a refined Lagrange bases. More...

#include <dune/functions/functionspacebases/refinedlagrangebasis.hh>

Public Types

using GridView = GV
 The grid view that the FE basis is defined on.
 
using Node = RefinedLagrangeNode< GV, k, R >
 Type of the refined Lagrange tree node.
 
using size_type = std::size_t
 Type used for index digits.
 

Public Member Functions

 RefinedLagrangePreBasis (const GridView &gv)
 Constructor for a given grid view object. More...
 
Node makeNode () const
 Create tree node.
 
void initializeIndices ()
 Initialize the global index information.
 
const GridViewgridView () const
 Export the stored GridView.
 
void update (const GridView &gv)
 Update the stored GridView.
 
size_type dimension () const
 Return total number of basis functions.
 
size_type maxNodeSize () const
 Return maximal number of basis functions per element.
 
template<class Node , class It >
It indices (const Node &node, It it) const
 Fill cache with global indices of DOFs associated to the given bound node.
 
size_type size (const SizePrefix &prefix) const
 Return number of possible values for next position in multi index.
 
size_type size () const
 Get the total dimension of the space spanned by this basis.
 
auto containerDescriptor () const
 Return a flat container-descriptor.
 

Static Public Member Functions

static constexpr unsigned int order ()
 Polynomial order used in the local Lagrange finite-elements. More...
 

Static Public Attributes

static constexpr size_type maxMultiIndexSize
 Maximal length of global multi-indices.
 
static constexpr size_type minMultiIndexSize
 Minimal length of global multi-indices.
 
static constexpr size_type multiIndexBufferSize
 Size required temporarily when constructing global multi-indices.
 

Detailed Description

template<typename GV, int k, typename R = double>
class Dune::Functions::RefinedLagrangePreBasis< GV, k, R >

A pre-basis for a refined Lagrange bases.

Template Parameters
GVThe grid view that the FE basis is defined on
kThe polynomial order of ansatz functions
RRange field-type used for shape function values
Note
This only works for simplex grids.

Constructor & Destructor Documentation

◆ RefinedLagrangePreBasis()

template<typename GV , int k, typename R = double>
Dune::Functions::RefinedLagrangePreBasis< GV, k, R >::RefinedLagrangePreBasis ( const GridView gv)
inline

Constructor for a given grid view object.

Parameters
gvThe GridView the basis is defined on.
Exceptions
Dune::NotImplementedIf an element of type !simplex is found.

References DUNE_THROW, and Dune::FloatCmp::gt().

Member Function Documentation

◆ order()

template<typename GV , int k, typename R = double>
static constexpr unsigned int Dune::Functions::RefinedLagrangePreBasis< GV, k, R >::order ( )
inlinestaticconstexpr

Polynomial order used in the local Lagrange finite-elements.

Note
The local function is of order k only in subdomains of the element. It might be necessary to use a subdivided quadrature rule for integration.

The documentation for this class was generated from the following file:
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 12, 23:30, 2024)