ug4
ug::BaseLSFS< TImpl, TDim, TShape, TGrad > Class Template Reference

static interface for trial spaces More...

#include <local_shape_function_set.h>

Public Types

typedef TGrad grad_type
 Gradient type. More...
 
typedef TImpl ImplType
 type of implementation More...
 
typedef TShape shape_type
 Shape type. More...
 

Public Member Functions

void grad (grad_type &g, size_t i, const MathVector< dim > &x) const
 evaluates the gradient of the shape function More...
 
void grads (grad_type *vGrad, const MathVector< dim > &x) const
 returns all gradients evaluated at a point More...
 
void grads (std::vector< grad_type > &vGrad, const MathVector< dim > &x) const
 returns all gradients evaluated at a point More...
 
void grads (std::vector< std::vector< grad_type > > &vvGrad, const std::vector< MathVector< dim > > &vLocPos) const
 returns all gradients evaluated at a point More...
 
void shape (shape_type &sh, size_t i, const MathVector< dim > &x) const
 evaluates the shape function More...
 
shape_type shape (size_t i, const MathVector< dim > &x) const
 evaluates the shape function More...
 
void shapes (shape_type *vShape, const MathVector< dim > &x) const
 returns all shape functions evaluated at a point More...
 
void shapes (std::vector< shape_type > &vShape, const MathVector< dim > &x) const
 returns all shape functions evaluated at a point More...
 
void shapes (std::vector< std::vector< shape_type > > &vvShape, const std::vector< MathVector< dim > > &vLocPos) const
 returns all shape functions evaluated at a point More...
 

Static Public Attributes

static const int dim = TDim
 dimension of reference element More...
 

Protected Member Functions

ImplTypegetImpl ()
 access to implementation More...
 
const ImplTypegetImpl () const
 const access to implementation More...
 

Detailed Description

template<typename TImpl, int TDim, typename TShape = number, typename TGrad = MathVector<TDim>>
class ug::BaseLSFS< TImpl, TDim, TShape, TGrad >

static interface for trial spaces

Member Typedef Documentation

◆ grad_type

template<typename TImpl , int TDim, typename TShape = number, typename TGrad = MathVector<TDim>>
typedef TGrad ug::BaseLSFS< TImpl, TDim, TShape, TGrad >::grad_type

Gradient type.

◆ ImplType

template<typename TImpl , int TDim, typename TShape = number, typename TGrad = MathVector<TDim>>
typedef TImpl ug::BaseLSFS< TImpl, TDim, TShape, TGrad >::ImplType

type of implementation

◆ shape_type

template<typename TImpl , int TDim, typename TShape = number, typename TGrad = MathVector<TDim>>
typedef TShape ug::BaseLSFS< TImpl, TDim, TShape, TGrad >::shape_type

Shape type.

Member Function Documentation

◆ getImpl() [1/2]

template<typename TImpl , int TDim, typename TShape = number, typename TGrad = MathVector<TDim>>
ImplType& ug::BaseLSFS< TImpl, TDim, TShape, TGrad >::getImpl ( )
inlineprotected

◆ getImpl() [2/2]

template<typename TImpl , int TDim, typename TShape = number, typename TGrad = MathVector<TDim>>
const ImplType& ug::BaseLSFS< TImpl, TDim, TShape, TGrad >::getImpl ( ) const
inlineprotected

const access to implementation

◆ grad()

template<typename TImpl , int TDim, typename TShape = number, typename TGrad = MathVector<TDim>>
void ug::BaseLSFS< TImpl, TDim, TShape, TGrad >::grad ( grad_type g,
size_t  i,
const MathVector< dim > &  x 
) const
inline

evaluates the gradient of the shape function

This function returns the gradient of Shape Function i at an element-local evaluation point.

Parameters
[in]inumber of DoF
[in]xPosition on reference element (evaluation point)
Returns
gradient at point

References ug::BaseLSFS< TImpl, TDim, TShape, TGrad >::getImpl().

Referenced by ug::BaseLSFS< TImpl, TDim, TShape, TGrad >::grads().

◆ grads() [1/3]

template<typename TImpl , int TDim, typename TShape = number, typename TGrad = MathVector<TDim>>
void ug::BaseLSFS< TImpl, TDim, TShape, TGrad >::grads ( grad_type vGrad,
const MathVector< dim > &  x 
) const
inline

returns all gradients evaluated at a point

This function returns the gradients of all Shape Functions at an element-local evaluation point in an array.

Parameters
[out]vGradVector of gradients
[in]xPosition on reference element (evaluation point)

References ug::BaseLSFS< TImpl, TDim, TShape, TGrad >::getImpl(), and ug::BaseLSFS< TImpl, TDim, TShape, TGrad >::grad().

Referenced by ug::BaseLSFS< TImpl, TDim, TShape, TGrad >::grads().

◆ grads() [2/3]

template<typename TImpl , int TDim, typename TShape = number, typename TGrad = MathVector<TDim>>
void ug::BaseLSFS< TImpl, TDim, TShape, TGrad >::grads ( std::vector< grad_type > &  vGrad,
const MathVector< dim > &  x 
) const
inline

returns all gradients evaluated at a point

This function returns the gradients of all Shape Functions at an element-local evaluation point in an array.

Parameters
[out]vGradVector of gradients
[in]xPosition on reference element (evaluation point)

References ug::BaseLSFS< TImpl, TDim, TShape, TGrad >::getImpl(), and ug::BaseLSFS< TImpl, TDim, TShape, TGrad >::grads().

◆ grads() [3/3]

template<typename TImpl , int TDim, typename TShape = number, typename TGrad = MathVector<TDim>>
void ug::BaseLSFS< TImpl, TDim, TShape, TGrad >::grads ( std::vector< std::vector< grad_type > > &  vvGrad,
const std::vector< MathVector< dim > > &  vLocPos 
) const
inline

returns all gradients evaluated at a point

This function returns the gradients of all Shape Functions at an element-local evaluation point in an array.

Parameters
[out]vGradVector of gradients
[in]xPosition on reference element (evaluation point)

References ug::BaseLSFS< TImpl, TDim, TShape, TGrad >::grads().

◆ shape() [1/2]

template<typename TImpl , int TDim, typename TShape = number, typename TGrad = MathVector<TDim>>
void ug::BaseLSFS< TImpl, TDim, TShape, TGrad >::shape ( shape_type sh,
size_t  i,
const MathVector< dim > &  x 
) const
inline

evaluates the shape function

This function returns the value of Shape Function i at an element-local evaluation point.

Parameters
[in]inumber of DoF
[in]xPosition on reference element (evaluation point)
Returns
shape function value at point

References ug::BaseLSFS< TImpl, TDim, TShape, TGrad >::shape().

◆ shape() [2/2]

template<typename TImpl , int TDim, typename TShape = number, typename TGrad = MathVector<TDim>>
shape_type ug::BaseLSFS< TImpl, TDim, TShape, TGrad >::shape ( size_t  i,
const MathVector< dim > &  x 
) const
inline

evaluates the shape function

This function returns the value of Shape Function i at an element-local evaluation point.

Parameters
[in]inumber of DoF
[in]xPosition on reference element (evaluation point)
Returns
shape function value at point

References ug::BaseLSFS< TImpl, TDim, TShape, TGrad >::getImpl().

Referenced by ug::BaseLSFS< TImpl, TDim, TShape, TGrad >::shape(), and ug::BaseLSFS< TImpl, TDim, TShape, TGrad >::shapes().

◆ shapes() [1/3]

template<typename TImpl , int TDim, typename TShape = number, typename TGrad = MathVector<TDim>>
void ug::BaseLSFS< TImpl, TDim, TShape, TGrad >::shapes ( shape_type vShape,
const MathVector< dim > &  x 
) const
inline

returns all shape functions evaluated at a point

This function returns the values of all Shape Functions at an element-local evaluation point in an array.

Parameters
[out]vShapeVector of Shapes
[in]xPosition on reference element (evaluation point)

References ug::BaseLSFS< TImpl, TDim, TShape, TGrad >::getImpl(), and ug::BaseLSFS< TImpl, TDim, TShape, TGrad >::shape().

Referenced by ug::FV1InnerBoundaryElemDisc< TImpl, TDomain >::prep_err_est_elem(), ug::FV1InnerBoundaryElemDisc< TImpl, TDomain >::prep_err_est_elem_loop(), and ug::BaseLSFS< TImpl, TDim, TShape, TGrad >::shapes().

◆ shapes() [2/3]

template<typename TImpl , int TDim, typename TShape = number, typename TGrad = MathVector<TDim>>
void ug::BaseLSFS< TImpl, TDim, TShape, TGrad >::shapes ( std::vector< shape_type > &  vShape,
const MathVector< dim > &  x 
) const
inline

returns all shape functions evaluated at a point

This function returns the values of all Shape Functions at an element-local evaluation point in an array.

Parameters
[out]vShapeVector of Shapes
[in]xPosition on reference element (evaluation point)

References ug::BaseLSFS< TImpl, TDim, TShape, TGrad >::getImpl(), and ug::BaseLSFS< TImpl, TDim, TShape, TGrad >::shapes().

◆ shapes() [3/3]

template<typename TImpl , int TDim, typename TShape = number, typename TGrad = MathVector<TDim>>
void ug::BaseLSFS< TImpl, TDim, TShape, TGrad >::shapes ( std::vector< std::vector< shape_type > > &  vvShape,
const std::vector< MathVector< dim > > &  vLocPos 
) const
inline

returns all shape functions evaluated at a point

This function returns the values of all Shape Functions at an element-local evaluation point in an array.

Parameters
[out]vShapeVector of Shapes
[in]xPosition on reference element (evaluation point)

References ug::BaseLSFS< TImpl, TDim, TShape, TGrad >::shapes().

Member Data Documentation

◆ dim

template<typename TImpl , int TDim, typename TShape = number, typename TGrad = MathVector<TDim>>
const int ug::BaseLSFS< TImpl, TDim, TShape, TGrad >::dim = TDim
static

dimension of reference element


The documentation for this class was generated from the following file: