Dune Core Modules (2.10.0)
Implementations of the GlobalBasis concept. More...
Classes | |
class | Dune::Functions::BSplineLocalBasis< GV, R > |
LocalBasis class in the sense of dune-localfunctions, presenting the restriction of a B-spline patch to a knot span. More... | |
class | Dune::Functions::BSplineLocalCoefficients< dim > |
Attaches a shape function to an entity. More... | |
class | Dune::Functions::BSplineLocalInterpolation< dim, LB > |
Local interpolation in the sense of dune-localfunctions, for the B-spline basis on tensor-product grids. More... | |
class | Dune::Functions::BSplineLocalFiniteElement< GV, R > |
LocalFiniteElement in the sense of dune-localfunctions, for the B-spline basis on tensor-product grids. More... | |
class | Dune::Functions::BSplinePreBasis< GV > |
Pre-basis for B-spline basis. More... | |
class | Dune::Functions::CompositePreBasis< IMS, SPB > |
A pre-basis for composite bases. More... | |
class | Dune::Functions::HierarchicalLagrangePreBasis< GV, k, R > |
A pre-basis for a hierarchical Lagrange basis. More... | |
class | Dune::Functions::LagrangePreBasis< GV, k, R > |
A pre-basis for a PQ-lagrange bases with given order. More... | |
class | Dune::Functions::LFEPreBasisMixin< GV, LFE > |
A pre-basis mixin class parametrized with a local finite-element and a DOF layout. More... | |
class | Dune::Functions::RannacherTurekPreBasis< GV > |
Pre-basis for a Rannacher-Turek basis. More... | |
class | Dune::Functions::RefinedLagrangePreBasis< GV, k, R > |
A pre-basis for a refined Lagrange bases. More... | |
class | Dune::Functions::TaylorHoodPreBasis< GV, HI > |
Pre-basis for lowest order Taylor-Hood basis. More... | |
class | Dune::Functions::CubicHermitePreBasis< GV, reduced > |
Pre-basis for a CubicHermite basis. More... | |
class | Dune::Functions::Experimental::TransformedIndexPreBasis< RPB, T > |
A pre-basis transforming multi-indices. More... | |
Typedefs | |
template<typename GV > | |
using | Dune::Functions::BSplineBasis = DefaultGlobalBasis< BSplinePreBasis< GV > > |
A global B-spline basis. More... | |
template<typename GV , int k, typename R = double> | |
using | Dune::Functions::HierarchicalLagrangeBasis = DefaultGlobalBasis< HierarchicalLagrangePreBasis< GV, k, R > > |
Basis of a scalar Hierarchical Lagrange finite element space. More... | |
template<typename GV , int k = -1, typename R = double> | |
using | Dune::Functions::LagrangeBasis = DefaultGlobalBasis< LagrangePreBasis< GV, k, R > > |
Nodal basis of a scalar k-th-order Lagrangean finite element space. More... | |
template<typename GV , int k> | |
using | Dune::Functions::LagrangeDGBasis = DefaultGlobalBasis< LagrangeDGPreBasis< GV, k > > |
Basis of a scalar k-th-order Lagrangean-DG finite element space. More... | |
template<typename GV > | |
using | Dune::Functions::RannacherTurekBasis = DefaultGlobalBasis< RannacherTurekPreBasis< GV > > |
Rannacher-Turek basis. More... | |
template<typename GV , int k, typename R = double> | |
using | Dune::Functions::RefinedLagrangeBasis = DefaultGlobalBasis< RefinedLagrangePreBasis< GV, k, R > > |
Nodal basis of a continuous Lagrange finite-element space on a uniformly refined simplex element. More... | |
template<typename GV > | |
using | Dune::Functions::TaylorHoodBasis = DefaultGlobalBasis< TaylorHoodPreBasis< GV > > |
Nodal basis for a lowest order Taylor-Hood Lagrangean finite element space. More... | |
template<typename GV > | |
using | Dune::Functions::CubicHermiteBasis = DefaultGlobalBasis< CubicHermitePreBasis< GV > > |
CubicHermite basis. More... | |
template<typename GV > | |
using | Dune::Functions::ReducedCubicHermiteBasis = DefaultGlobalBasis< CubicHermitePreBasis< GV, true > > |
CubicHermite basis. More... | |
Functions | |
template<std::size_t k> | |
auto | Dune::Functions::BasisFactory::brezziDouglasMarini () |
Create a pre-basis factory that can create a Brezzi-Douglas-Marini pre-basis. More... | |
auto | Dune::Functions::BasisFactory::bSpline (const std::vector< double > &knotVector, unsigned int order, bool makeOpen=true) |
Create a pre-basis factory that can create a B-spline pre-basis. | |
template<typename... Args, std::enable_if_t< Concept::isIndexMergingStrategy< typename LastType< Args... >::type >(), int > = 0> | |
auto | Dune::Functions::BasisFactory::composite (Args &&... args) |
Create a factory builder that can build a CompositePreBasis. More... | |
template<class ChildPreBasisFactory , class IndexMergingStrategy > | |
auto | Dune::Functions::BasisFactory::power (ChildPreBasisFactory &&childPreBasisFactory, std::size_t k, const IndexMergingStrategy &) |
Create a pre-basis factory that can build a PowerPreBasis. More... | |
template<class ChildPreBasisFactory > | |
auto | Dune::Functions::BasisFactory::power (ChildPreBasisFactory &&childPreBasisFactory, std::size_t k) |
Create a factory builder that can build a PowerPreBasis. More... | |
template<int k, typename R = double> | |
auto | Dune::Functions::BasisFactory::hierarchicalLagrange () |
A factory that can create a HierarchicalLagrange pre-basis. More... | |
template<std::size_t k, typename R = double> | |
auto | Dune::Functions::BasisFactory::lagrange () |
Create a pre-basis factory that can create a Lagrange pre-basis. More... | |
template<typename R = double> | |
auto | Dune::Functions::BasisFactory::lagrange (int order) |
Create a pre-basis factory that can create a Lagrange pre-basis with a run-time order. More... | |
template<std::size_t k> | |
auto | Dune::Functions::BasisFactory::lagrangeDG () |
Create a pre-basis factory that can create a LagrangeDG pre-basis. More... | |
template<std::size_t kind, std::size_t order, typename Range = double> | |
auto | Dune::Functions::BasisFactory::nedelec () |
Create a pre-basis factory that can create a Nédélec pre-basis. More... | |
template<class RawPreBasisIndicator , class PIS > | |
auto | Dune::Functions::BasisFactory::Experimental::periodic (RawPreBasisIndicator &&rawPreBasisIndicator, PIS &&periodicIndexSet) |
Create a pre-basis factory that can create a periodic pre-basis. More... | |
template<std::size_t k, class ChildPreBasisFactory , class IndexMergingStrategy > | |
auto | Dune::Functions::BasisFactory::power (ChildPreBasisFactory &&childPreBasisFactory, const IndexMergingStrategy &) |
Create a pre-basis factory that can build a PowerPreBasis. More... | |
template<std::size_t k, class ChildPreBasisFactory > | |
auto | Dune::Functions::BasisFactory::power (ChildPreBasisFactory &&childPreBasisFactory) |
Create a factory builder that can build a PowerPreBasis. More... | |
template<class Dummy = void> | |
auto | Dune::Functions::BasisFactory::rannacherTurek () |
Create a pre-basis factory that can create a Rannacher-Turek pre-basis. | |
template<std::size_t k> | |
auto | Dune::Functions::BasisFactory::raviartThomas () |
Create a pre-basis factory that can create a Raviart-Thomas pre-basis. More... | |
template<int k, typename R = double> | |
auto | Dune::Functions::BasisFactory::refinedLagrange () |
Create a pre-basis factory that can create a RefinedLagrange pre-basis. More... | |
auto | Dune::Functions::BasisFactory::taylorHood () |
Create a pre-basis factory that can create a Taylor-Hood pre-basis. | |
template<class Dummy = void> | |
auto | Dune::Functions::BasisFactory::cubicHermite () |
Create a pre-basis factory that can create a CubicHermite pre-basis. | |
template<class Dummy = void> | |
auto | Dune::Functions::BasisFactory::reducedCubicHermite () |
Create a pre-basis factory that can create a CubicHermite pre-basis. | |
Detailed Description
Implementations of the GlobalBasis concept.
Typedef Documentation
◆ BSplineBasis
using Dune::Functions::BSplineBasis = typedef DefaultGlobalBasis<BSplinePreBasis<GV> > |
A global B-spline basis.
- Template Parameters
-
GV The GridView that the space is defined on
◆ CubicHermiteBasis
using Dune::Functions::CubicHermiteBasis = typedef DefaultGlobalBasis<CubicHermitePreBasis<GV> > |
CubicHermite basis.
- Template Parameters
-
GV The GridView that the space is defined on
◆ HierarchicalLagrangeBasis
using Dune::Functions::HierarchicalLagrangeBasis = typedef DefaultGlobalBasis<HierarchicalLagrangePreBasis<GV, k, R> > |
Basis of a scalar Hierarchical Lagrange finite element space.
- Template Parameters
-
GV The GridView that the space is defined on k The order of the basis (0 < k <= 2) R The range type of the local basis
- Note
- currently only supports simplex grids
◆ LagrangeBasis
using Dune::Functions::LagrangeBasis = typedef DefaultGlobalBasis<LagrangePreBasis<GV, k, R> > |
Nodal basis of a scalar k-th-order Lagrangean finite element space.
- 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
All arguments passed to the constructor will be forwarded to the constructor of LagrangePreBasis.
- Warning
- The implementation of the basis with run-time order order uses the LagrangeFiniteElement implementation of dune-localfunctions, which is known to violate strict-aliasing rules (see https://gitlab.dune-project.org/core/dune-localfunctions/issues/14) Keep this in mind if ever you experience difficult-to-explain crashes or wrong results.
- Template Parameters
-
GV The GridView that the space is defined on k The order of the basis; -1 means 'order determined at run-time' R The range type of the local basis
◆ LagrangeDGBasis
using Dune::Functions::LagrangeDGBasis = typedef DefaultGlobalBasis<LagrangeDGPreBasis<GV, k> > |
Basis of a scalar k-th-order Lagrangean-DG finite element space.
- Template Parameters
-
GV The GridView that the space is defined on k The order of the basis
◆ RannacherTurekBasis
using Dune::Functions::RannacherTurekBasis = typedef DefaultGlobalBasis<RannacherTurekPreBasis<GV> > |
Rannacher-Turek basis.
These are Crouzeix-Raviart-elements for quadrilateral elements. See Rolf Rannacher and Stefan Turek. Simple nonconforming quadrilateral Stokes element. Numerical Methods for Partial Differential Equations, 8:97–111, 1992.
- Template Parameters
-
GV The GridView that the space is defined on
◆ ReducedCubicHermiteBasis
using Dune::Functions::ReducedCubicHermiteBasis = typedef DefaultGlobalBasis<CubicHermitePreBasis<GV, true> > |
CubicHermite basis.
- Template Parameters
-
GV The GridView that the space is defined on
◆ RefinedLagrangeBasis
using Dune::Functions::RefinedLagrangeBasis = typedef DefaultGlobalBasis<RefinedLagrangePreBasis<GV,k,R> > |
Nodal basis of a continuous Lagrange finite-element space on a uniformly refined simplex element.
- Note
- This only works for simplex grids.
All arguments passed to the constructor will be forwarded to the constructor of RefinedLagrangeBasis.
- Template Parameters
-
GV The GridView that the space is defined on k The order of the basis R The range type of the local basis
◆ TaylorHoodBasis
using Dune::Functions::TaylorHoodBasis = typedef DefaultGlobalBasis<TaylorHoodPreBasis<GV> > |
Nodal basis for a lowest order Taylor-Hood Lagrangean finite element space.
- Template Parameters
-
GV The GridView that the space is defined on.
- Note
- This mainly serves as an example, since you can construct a basis with the same functionality manually using static const int k = 1;auto taylorHoodBasis = makeBasis(gridView,power<dim>(lagrange<k+1>(),flatInterleaved()),lagrange<k>()));auto composite(Args &&... args)Create a factory builder that can build a CompositePreBasis.Definition: compositebasis.hh:397constexpr FlatInterleaved flatInterleaved()Creates an interleaved merging of direct children without blocking.Definition: basistags.hh:202
Function Documentation
◆ brezziDouglasMarini()
auto Dune::Functions::BasisFactory::brezziDouglasMarini | ( | ) |
Create a pre-basis factory that can create a Brezzi-Douglas-Marini pre-basis.
- Template Parameters
-
k Order of the Brezzi-Douglas-Marini element
◆ composite()
auto Dune::Functions::BasisFactory::composite | ( | Args &&... | args | ) |
Create a factory builder that can build a CompositePreBasis.
- Template Parameters
-
Args Types of child factory builders and IndexMergingStrategy type
- Parameters
-
args Child factory builder objects and an IndexMergingStrategy
This is the overload used if the last argument is an IndexMergingStrategy.
- Template Parameters
-
Args Types of child factory builders
- Parameters
-
args Child factory builder objects
This is the overload used if no IndexMergingStrategy is supplied. In this case the BasisFactory::BlockedLexicographic strategy is used.
References Dune::applyPartial().
◆ hierarchicalLagrange()
auto Dune::Functions::BasisFactory::hierarchicalLagrange | ( | ) |
A factory that can create a HierarchicalLagrange pre-basis.
- Template Parameters
-
k The polynomial order of the ansatz functions (0 < k <= 2) R The range field-type of the local basis
◆ lagrange() [1/2]
auto Dune::Functions::BasisFactory::lagrange | ( | ) |
Create a pre-basis factory that can create a Lagrange pre-basis.
- Template Parameters
-
k The polynomial order of the ansatz functions; -1 means 'order determined at run-time' R The range type of the local basis
◆ lagrange() [2/2]
auto Dune::Functions::BasisFactory::lagrange | ( | int | order | ) |
Create a pre-basis factory that can create a Lagrange pre-basis with a run-time order.
- Template Parameters
-
R The range type of the local basis
◆ lagrangeDG()
auto Dune::Functions::BasisFactory::lagrangeDG | ( | ) |
Create a pre-basis factory that can create a LagrangeDG pre-basis.
- Template Parameters
-
k The polynomial order of the ansatz functions
◆ nedelec()
auto Dune::Functions::BasisFactory::nedelec | ( | ) |
Create a pre-basis factory that can create a Nédélec pre-basis.
- Template Parameters
-
kind Kind of the Nédélec element (1 or 2) order Order of the Nédélec element (lowest order is '1') Range Number type used for shape function values
◆ periodic()
auto Dune::Functions::BasisFactory::Experimental::periodic | ( | RawPreBasisIndicator && | rawPreBasisIndicator, |
PIS && | periodicIndexSet | ||
) |
Create a pre-basis factory that can create a periodic pre-basis.
- Parameters
-
rawPreBasisIndicator Object encoding the raw non-periodic basis periodicIndexSet A PeriodicIndexSet containing the indices to be identified
The rawPreBasisIndicator can either be a PreBasisFactory, a PreBasis, or a GlobalBasis. In the latter two cases the multi index type of those bases and the periodic basis to be constructed, has to coincide. Both arguments will be copied. Currently only wrapped bases with flat indices are supported.
Referenced by Dune::Torus< Communication, d >::is_neighbor(), Dune::YaspGrid< dim, Coordinates >::makelevel(), Dune::BackupRestoreFacility< Dune::YaspGrid< dim, Coordinates > >::restore(), Dune::BackupRestoreFacility< YaspGrid< dim, TensorProductCoordinates< ctype, dim > > >::restore(), and Dune::YaspGrid< dim, Coordinates >::YaspGrid().
◆ power() [1/4]
auto Dune::Functions::BasisFactory::power | ( | ChildPreBasisFactory && | childPreBasisFactory | ) |
Create a factory builder that can build a PowerPreBasis.
- Template Parameters
-
ChildPreBasisFactory Types of child pre-basis factory
- Parameters
-
childPreBasisFactory Child pre-basis factory
This overload will select the BasisFactory::BlockedInterleaved strategy.
◆ power() [2/4]
auto Dune::Functions::BasisFactory::power | ( | ChildPreBasisFactory && | childPreBasisFactory, |
const IndexMergingStrategy & | |||
) |
Create a pre-basis factory that can build a PowerPreBasis.
- Template Parameters
-
ChildPreBasisFactory Types of child pre-basis factory IndexMergingStrategy An IndexMergingStrategy type
- Parameters
-
childPreBasisFactory Child pre-basis factory ims IndexMergingStrategy to be used
This overload can be used to explicitly supply an IndexMergingStrategy.
◆ power() [3/4]
auto Dune::Functions::BasisFactory::power | ( | ChildPreBasisFactory && | childPreBasisFactory, |
std::size_t | k | ||
) |
Create a factory builder that can build a PowerPreBasis.
- Template Parameters
-
ChildPreBasisFactory Types of child pre-basis factory
- Parameters
-
childPreBasisFactory Child pre-basis factory
This overload will select the BasisFactory::BlockedInterleaved strategy.
◆ power() [4/4]
auto Dune::Functions::BasisFactory::power | ( | ChildPreBasisFactory && | childPreBasisFactory, |
std::size_t | k, | ||
const IndexMergingStrategy & | |||
) |
Create a pre-basis factory that can build a PowerPreBasis.
- Template Parameters
-
ChildPreBasisFactory Types of child pre-basis factory IndexMergingStrategy An IndexMergingStrategy type
- Parameters
-
childPreBasisFactory Child pre-basis factory ims IndexMergingStrategy to be used
This overload can be used to explicitly supply an IndexMergingStrategy.
◆ raviartThomas()
auto Dune::Functions::BasisFactory::raviartThomas | ( | ) |
Create a pre-basis factory that can create a Raviart-Thomas pre-basis.
- Template Parameters
-
k Order of the Raviart-Thomas element
◆ refinedLagrange()
auto Dune::Functions::BasisFactory::refinedLagrange | ( | ) |
Create a pre-basis factory that can create a RefinedLagrange pre-basis.
- Template Parameters
-
R The range type of the local basis