Dune Core Modules (2.6.0)

Dune::Pk1DLocalBasis< D, R, k > Class Template Reference

Lagrange shape functions of arbitrary order on the 1D reference triangle. More...

#include <dune/localfunctions/lagrange/pk1d/pk1dlocalbasis.hh>

Public Types

enum  
 Export the number of degrees of freedom.
 
enum  
 Export the element order.
 

Public Member Functions

 Pk1DLocalBasis ()
 Standard constructor.
 
unsigned int size () const
 number of shape functions
 
void evaluateFunction (const typename Traits::DomainType &x, std::vector< typename Traits::RangeType > &out) const
 Evaluate all shape functions.
 
void evaluateJacobian (const typename Traits::DomainType &x, std::vector< typename Traits::JacobianType > &out) const
 Evaluate Jacobian of all shape functions.
 
void partial (const std::array< unsigned int, 1 > &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, unsigned int k>
class Dune::Pk1DLocalBasis< D, R, k >

Lagrange shape functions of arbitrary order on the 1D reference triangle.

Lagrange shape functions of arbitrary order have the property that \(\hat\phi^i(x_j) = \delta_{i,j}\) for certain points \(x_j\).

Template Parameters
DType to represent the field in the domain.
RType to represent the field in the range.
kPolynomial order.

Member Function Documentation

◆ partial()

template<class D , class R , unsigned int k>
void Dune::Pk1DLocalBasis< D, R, k >::partial ( const std::array< unsigned int, 1 > &  order,
const typename Traits::DomainType in,
std::vector< typename Traits::RangeType > &  out 
) const
inline

Evaluate partial derivatives of any order of all shape functions.

Parameters
orderOrder of the partial derivatives, in the classic multi-index notation
inPosition where to evaluate the derivatives
[out]outReturn value: the desired partial derivatives

References DUNE_THROW, Dune::Pk1DLocalBasis< D, R, k >::evaluateFunction(), and Dune::Pk1DLocalBasis< D, R, k >::order().


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 (Dec 26, 23:30, 2024)