DUNE PDELab (git)
Dune::Functions::BSplineLocalFiniteElement< GV, R > Class Template Reference
LocalFiniteElement in the sense of dune-localfunctions, for the B-spline basis on tensor-product grids. More...
#include <dune/functions/functionspacebases/bsplinebasis.hh>
Public Types | |
typedef LocalFiniteElementTraits< BSplineLocalBasis< GV, R >, BSplineLocalCoefficients< dim >, BSplineLocalInterpolation< dim, BSplineLocalBasis< GV, R > > > | Traits |
Export various types related to this LocalFiniteElement. | |
Public Member Functions | |
BSplineLocalFiniteElement (const BSplinePreBasis< GV > &preBasis) | |
Constructor with a given B-spline basis. | |
BSplineLocalFiniteElement (const BSplineLocalFiniteElement &other) | |
Copy constructor. | |
void | bind (const std::array< unsigned, dim > &elementIdx) |
Bind LocalFiniteElement to a specific knot span of the spline patch. More... | |
const BSplineLocalBasis< GV, R > & | localBasis () const |
Hand out a LocalBasis object. | |
const BSplineLocalCoefficients< dim > & | localCoefficients () const |
Hand out a LocalCoefficients object. | |
const BSplineLocalInterpolation< dim, BSplineLocalBasis< GV, R > > & | localInterpolation () const |
Hand out a LocalInterpolation object. | |
unsigned | size () const |
Number of shape functions in this finite element. | |
GeometryType | type () const |
Return the reference element that the local finite element is defined on (here, a hypercube) | |
unsigned int | size (int i) const |
Number of degrees of freedom for one coordinate direction. | |
Detailed Description
template<class GV, class R>
class Dune::Functions::BSplineLocalFiniteElement< GV, R >
class Dune::Functions::BSplineLocalFiniteElement< GV, R >
LocalFiniteElement in the sense of dune-localfunctions, for the B-spline basis on tensor-product grids.
This class ties together the implementation classes BSplineLocalBasis, BSplineLocalCoefficients, and BSplineLocalInterpolation
- Template Parameters
-
D Number type used for domain coordinates R Number type used for spline function values
Member Function Documentation
◆ bind()
template<class GV , class R >
|
inline |
Bind LocalFiniteElement to a specific knot span of the spline patch.
Elements are the non-empty knot spans, here we do the renumbering
- Parameters
-
ijk Integer coordinates in the tensor product patch
References Dune::Functions::BSplineLocalFiniteElement< GV, R >::size().
The documentation for this class was generated from the following file:
- dune/functions/functionspacebases/bsplinebasis.hh
|
Legal Statements / Impressum |
Hosted by TU Dresden |
generated with Hugo v0.111.3
(Nov 12, 23:30, 2024)