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 | |
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. | |
unsigned int | order () const |
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 local coordinates in the subtriangle. | |
static void | getSubElement (const FieldVector< D, 2 > &global, int &subElement, FieldVector< D, 2 > &local) |
Get local coordinates in the subtriangle. |
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)
D | Type to represent the field in the domain. | |
R | Type to represent the field in the range. |
typedef LocalBasisTraits<D,2,Dune::FieldVector<D,2>,R,1,Dune::FieldVector<R,1>, Dune::FieldMatrix<R,1,2> > Dune::RefinedP1LocalBasis< D, R, 2 >::Traits |
export type traits for function signature
void Dune::RefinedP1LocalBasis< D, R, 2 >::evaluateFunction | ( | const typename Traits::DomainType & | in, | |
std::vector< typename Traits::RangeType > & | out | |||
) | const [inline] |
Evaluate all shape functions.
void Dune::RefinedP1LocalBasis< D, R, 2 >::evaluateJacobian | ( | const typename Traits::DomainType & | in, | |
std::vector< typename Traits::JacobianType > & | out | |||
) | const [inline] |
Evaluate Jacobian of all shape functions.
static void Dune::RefinedSimplexLocalBasis< D, 2 >::getSubElement | ( | const FieldVector< D, 2 > & | global, | |
int & | subElement, | |||
FieldVector< D, 2 > & | local | |||
) | [inline, static, protected, inherited] |
Get local coordinates in the subtriangle.
[in] | global | Coordinates in the reference triangle |
[out] | subElement | Which of the four subtriangles is global in? |
[out] | local | The local coordinates in the subtriangle |
static int Dune::RefinedSimplexLocalBasis< D, 2 >::getSubElement | ( | const FieldVector< D, 2 > & | global | ) | [inline, static, protected, inherited] |
Get local coordinates in the subtriangle.
The triangles are ordered according to
|\ |2\ |--\ ||\ |0\|1\ ------
[in] | global | Coordinates in the reference triangle |
unsigned int Dune::RefinedP1LocalBasis< D, R, 2 >::order | ( | ) | const [inline] |
Polynomial order of the shape functions Doesn't really apply: these shape functions are only piecewise linear.
unsigned int Dune::RefinedP1LocalBasis< D, R, 2 >::size | ( | ) | const [inline] |
number of shape functions