Dune Core Modules (2.9.0)

Uniformly refined linear Lagrange shape functions on the triangle. More...

#include <dune/localfunctions/refined/refinedp1/refinedp1localbasis.hh>

Public Types

typedef LocalBasisTraits< D, 2, Dune::FieldVector< D, 2 >, R, 1, Dune::FieldVector< R, 1 >, Dune::FieldMatrix< R, 1, 2 > > Traits
 export type traits for function signature
 

Public Member 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.
 
void partial (const std::array< unsigned int, 2 > &order, const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
 Evaluate partial derivatives of all shape functions.
 

Static Public Member Functions

static constexpr unsigned int size ()
 number of shape functions
 
static constexpr unsigned int order ()
 Polynomial order of the shape functions Doesn't really apply: these shape functions are only piecewise linear.
 

Static Protected Member Functions

static int getSubElement (const FieldVector< D, 2 > &global)
 Get the number of the subtriangle containing a given point. More...
 
static void getSubElement (const FieldVector< D, 2 > &global, int &subElement, FieldVector< D, 2 > &local)
 Get local coordinates in the subtriangle. More...
 

Detailed Description

template<class D, class R>
class Dune::RefinedP1LocalBasis< D, R, 2 >

Uniformly refined linear Lagrange shape functions on the triangle.

This shape function set mimicks the P1 shape functions that you would get on a uniformly refined grid. Hence these shape functions are only piecewise linear! The data layout is identical to P2 shape functions.

Shape functions like these are necessary for hierarchical error estimators for certain nonlinear problems.

The functions are associated to points by:

f_0 ~ (0.0, 0.0) f_1 ~ (0.5, 0.0) f_2 ~ (1.0, 0.0) f_3 ~ (0.0, 0.5) f_4 ~ (0.5, 0.5) f_5 ~ (0.0, 1.0)

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

Member Function Documentation

◆ getSubElement() [1/2]

template<class D >
static int Dune::RefinedSimplexLocalBasis< D, 2 >::getSubElement ( const FieldVector< D, 2 > &  global)
inlinestaticprotectedinherited

Get the number of the subtriangle containing a given point.

The triangles are ordered according to

|\
|2\
|--\
|\3|\
|0\|1\
------
Parameters
[in]globalCoordinates in the reference triangle
Returns
Number of the subtriangle containing global

◆ getSubElement() [2/2]

template<class D >
static void Dune::RefinedSimplexLocalBasis< D, 2 >::getSubElement ( const FieldVector< D, 2 > &  global,
int &  subElement,
FieldVector< D, 2 > &  local 
)
inlinestaticprotectedinherited

Get local coordinates in the subtriangle.

Parameters
[in]globalCoordinates in the reference triangle
[out]subElementNumber of the subtriangle containing global
[out]localThe local coordinates in the subtriangle

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.80.0 (Apr 19, 22:31, 2024)