DUNE PDELab (2.8)
Dune::Functions::LagrangePreBasis< GV, k, MI, R > Class Template Reference
A pre-basis for a PQ-lagrange bases with given order. More...
#include <dune/functions/functionspacebases/lagrangebasis.hh>
Public Types | |
using | GridView = GV |
The grid view that the FE basis is defined on. | |
using | size_type = std::size_t |
Type used for indices and size information. | |
using | Node = LagrangeNode< GV, k, R > |
Template mapping root tree path to type of created tree node. | |
using | IndexSet = Impl::DefaultNodeIndexSet< LagrangePreBasis > |
Type of created tree node index set. More... | |
using | MultiIndex = MI |
Type used for global numbering of the basis vectors. | |
using | SizePrefix = Dune::ReservedVector< size_type, 1 > |
Type used for prefixes handed to the size() method. | |
Public Member Functions | |
LagrangePreBasis (const GridView &gv) | |
Constructor for a given grid view object with compile-time order. | |
LagrangePreBasis (const GridView &gv, unsigned int order) | |
Constructor for a given grid view object and run-time order. | |
void | initializeIndices () |
Initialize the global indices. | |
const GridView & | gridView () const |
Obtain the grid view that the basis is defined on. | |
void | update (const GridView &gv) |
Update the stored grid view, to be called if the grid has changed. | |
Node | makeNode () const |
Create tree node. | |
IndexSet | makeIndexSet () const |
Create tree node index set. More... | |
size_type | size () const |
Same as size(prefix) with empty prefix. | |
size_type | size (const SizePrefix prefix) const |
Return number of possible values for next position in multi index. | |
size_type | dimension () const |
Get the total dimension of the space spanned by this basis. | |
size_type | maxNodeSize () const |
Get the maximal number of DOFs associated to node for any element. | |
Protected Member Functions | |
size_type | dofsPerSimplex (std::size_t simplexDim) const |
Number of degrees of freedom assigned to a simplex (without the ones assigned to its faces!) | |
size_type | dofsPerCube (std::size_t cubeDim) const |
Number of degrees of freedom assigned to a cube (without the ones assigned to its faces!) | |
size_type | computeDofsPerSimplex (std::size_t simplexDim) const |
Number of degrees of freedom assigned to a simplex (without the ones assigned to its faces!) | |
size_type | computeDofsPerCube (std::size_t cubeDim) const |
Number of degrees of freedom assigned to a cube (without the ones assigned to its faces!) | |
Detailed Description
template<typename GV, int k, class MI, typename R>
class Dune::Functions::LagrangePreBasis< GV, k, MI, R >
class Dune::Functions::LagrangePreBasis< GV, k, MI, R >
A pre-basis for a PQ-lagrange bases with given order.
- Template Parameters
-
GV The grid view that the FE basis is defined on k The polynomial order of ansatz functions; -1 means 'order determined at run-time' MI Type to be used for multi-indices R Range type used for shape function values
- Note
- This only works for certain grids. The following restrictions hold
- If k is no larger than 2, then the grids can have any dimension
- If k is larger than 3 then the grid must be two-dimensional
- If k is 3, then the grid can be 3d if it is a simplex grid
Member Typedef Documentation
◆ IndexSet
template<typename GV , int k, class MI , typename R >
using Dune::Functions::LagrangePreBasis< GV, k, MI, R >::IndexSet = Impl::DefaultNodeIndexSet<LagrangePreBasis> |
Type of created tree node index set.
Member Function Documentation
◆ makeIndexSet()
template<typename GV , int k, class MI , typename R >
|
inline |
Create tree node index set.
Create an index set suitable for the tree node obtained by makeNode().
The documentation for this class was generated from the following file:
- dune/functions/functionspacebases/lagrangebasis.hh
|
Legal Statements / Impressum |
Hosted by TU Dresden |
generated with Hugo v0.111.3
(Dec 21, 23:30, 2024)