Dune Core Modules (2.6.0)
Dune::QkLocalBasis< D, R, k, d > Class Template Reference
Lagrange shape functions of order k on the reference cube. More...
#include <dune/localfunctions/lagrange/qk/qklocalbasis.hh>
Public Member Functions | |
unsigned int | size () const |
number of shape functions | |
void | evaluateFunction (const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const |
Evaluate all shape functions. | |
void | evaluateJacobian (const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const |
Evaluate Jacobian of all shape functions. More... | |
void | partial (const std::array< unsigned int, d > &order, const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const |
Evaluate partial derivatives of any order of all shape functions. More... | |
unsigned int | order () const |
Polynomial order of the shape functions. | |
Detailed Description
template<class D, class R, int k, int d>
class Dune::QkLocalBasis< D, R, k, d >
class Dune::QkLocalBasis< D, R, k, d >
Lagrange shape functions of order k on the reference cube.
Also known as \(Q^k\).
- Template Parameters
-
D Type to represent the field in the domain. R Type to represent the field in the range. k Polynomial degree d Dimension of the cube
Member Function Documentation
◆ evaluateJacobian()
template<class D , class R , int k, int d>
|
inline |
Evaluate Jacobian of all shape functions.
- Parameters
-
in position where to evaluate out The return value
References Dune::QkLocalBasis< D, R, k, d >::size().
◆ partial()
template<class D , class R , int k, int d>
|
inline |
Evaluate partial derivatives of any order of all shape functions.
- Parameters
-
order Order of the partial derivatives, in the classic multi-index notation in Position where to evaluate the derivatives [out] out Return value: the desired partial derivatives
References DUNE_THROW, Dune::QkLocalBasis< D, R, k, d >::order(), and Dune::QkLocalBasis< D, R, k, d >::size().
The documentation for this class was generated from the following file:
- dune/localfunctions/lagrange/qk/qklocalbasis.hh
|
Legal Statements / Impressum |
Hosted by TU Dresden |
generated with Hugo v0.111.3
(Nov 24, 23:30, 2024)