DUNE-FEM (unstable)

Dune::Fem::DefaultBasisFunctionSet< Entity, ShapeFunctionSet > Class Template Reference

#include <dune/fem/space/basisfunctionset/default.hh>

Public Types

typedef BaseType::EntityType EntityType
 entity type
 
typedef BaseType::Geometry Geometry
 geometry
 
typedef Geometry::ctype ctype
 type of coordinate field
 
typedef BaseType::ShapeFunctionSetType ShapeFunctionSetType
 shape function set type
 
typedef ToNewDimDomainFunctionSpace< LocalFunctionSpaceType, Geometry::coorddimension >::Type FunctionSpaceType
 type of function space
 
typedef FunctionSpaceType::RangeType RangeType
 range type
 
typedef FunctionSpaceType::DomainType DomainType
 domain type
 
typedef FunctionSpaceType::JacobianRangeType JacobianRangeType
 jacobian range type
 
typedef FunctionSpaceType::HessianRangeType HessianRangeType
 hessian range type
 
typedef std::decay_t< decltype(Dune::ReferenceElements< ctype, Geometry::coorddimension >::general(std::declval< const Dune::GeometryType & >())) > ReferenceElementType
 type of reference element
 

Public Member Functions

 DefaultBasisFunctionSet ()
 constructor
 
 DefaultBasisFunctionSet (const EntityType &entity, const ShapeFunctionSet &shapeFunctionSet=ShapeFunctionSet())
 constructor
 
template<class QuadratureType , class Vector , class DofVector >
void axpy (const QuadratureType &quad, const Vector &values, DofVector &dofs) const
 evaluate all basis function and multiply with given values and add to dofs
 
template<class QuadratureType , class VectorA , class VectorB , class DofVector >
void axpy (const QuadratureType &quad, const VectorA &valuesA, const VectorB &valuesB, DofVector &dofs) const
 evaluate all basis function and multiply with given values and add to dofs More...
 
template<class Point , class DofVector >
void axpy (const Point &x, const RangeType &valueFactor, DofVector &dofs) const
 evaluate all basis function and multiply with given values and add to dofs
 
template<class Point , class DofVector >
void axpy (const Point &x, const JacobianRangeType &jacobianFactor, DofVector &dofs) const
 evaluate all derivatives of all basis function and multiply with given values and add to dofs
 
template<class Point , class DofVector >
void axpy (const Point &x, const RangeType &valueFactor, const JacobianRangeType &jacobianFactor, DofVector &dofs) const
 evaluate all basis function and derivatives and multiply with given values and add to dofs
 
template<class Point , class DofVector >
void axpy (const Point &x, const HessianRangeType &hessianFactor, DofVector &dofs) const
 Add H:D^2phi to each dof.
 
template<class QuadratureType , class DofVector , class RangeArray >
void evaluateAll (const QuadratureType &quad, const DofVector &dofs, RangeArray &ranges) const
 
template<class Point , class DofVector >
void evaluateAll (const Point &x, const DofVector &dofs, RangeType &value) const
 
template<class Point , class RangeArray >
void evaluateAll (const Point &x, RangeArray &values) const
 
template<class QuadratureType , class DofVector , class JacobianArray >
void jacobianAll (const QuadratureType &quad, const DofVector &dofs, JacobianArray &jacobians) const
 
template<class Point , class DofVector >
void jacobianAll (const Point &x, const DofVector &dofs, JacobianRangeType &jacobian) const
 
template<class Point , class JacobianRangeArray >
void jacobianAll (const Point &x, JacobianRangeArray &jacobians) const
 
template<class QuadratureType , class DofVector , class HessianArray >
void hessianAll (const QuadratureType &quad, const DofVector &dofs, HessianArray &hessians) const
 
template<class Point , class DofVector >
void hessianAll (const Point &x, const DofVector &dofs, HessianRangeType &hessian) const
 
template<class Point , class HessianRangeArray >
void hessianAll (const Point &x, HessianRangeArray &hessians) const
 
const ShapeFunctionSetTypeshapeFunctionSet () const
 return shape function set
 
int order () const
 return order of basis function set
 
std::size_t size () const
 return size of basis function set
 
const Entityentity () const
 return entity
 
bool valid () const
 return true if entity pointer is set
 
const Geometrygeometry () const
 return geometry
 
Dune::GeometryType type () const
 return geometry type
 
const ReferenceElementTypereferenceElement () const
 return reference element
 
void bind (const EntityType &entity)
 set new entity object and geometry if enabled
 
void unbind ()
 release entity and geometry object
 

Protected Member Functions

template<class QuadratureType , class RangeArray , class DofVector >
void axpyImpl (const QuadratureType &quad, const RangeArray &rangeFactors, DofVector &dofs, const RangeType &) const
 evaluate all basis function and multiply with given values and add to dofs
 
template<class QuadratureType , class JacobianArray , class DofVector >
void axpyImpl (const QuadratureType &quad, const JacobianArray &jacobianFactors, DofVector &dofs, const JacobianRangeType &) const
 evaluate all basis function and multiply with given values and add to dofs
 
template<class QuadratureType , class HessianArray , class DofVector >
void axpyImpl (const QuadratureType &quad, const HessianArray &hessianFactors, DofVector &dofs, const HessianRangeType &) const
 evaluate all basis function and multiply with given values and add to dofs
 

Detailed Description

template<class Entity, class ShapeFunctionSet>
class Dune::Fem::DefaultBasisFunctionSet< Entity, ShapeFunctionSet >
   \brief implementation of a basis function set for given entity

   \tparam  Entity            entity type
   \tparam  ShapeFunctionSet  shape function set

   \note ShapeFunctionSet must be a copyable object. For most
         non-trivial implementations, you may want to use a
         proxy, see file
<dune/fem/space/shapefunctionset/proxy.hh>

Member Function Documentation

◆ axpy()

template<class Entity , class ShapeFunctionSet >
template<class QuadratureType , class VectorA , class VectorB , class DofVector >
void Dune::Fem::DefaultBasisFunctionSet< Entity, ShapeFunctionSet >::axpy ( const QuadratureType &  quad,
const VectorA &  valuesA,
const VectorB &  valuesB,
DofVector &  dofs 
) const
inline

evaluate all basis function and multiply with given values and add to dofs

Note
valuesA and valuesB can be vectors of RangeType or JacobianRangeType

References Dune::Fem::DefaultBasisFunctionSet< Entity, ShapeFunctionSet >::axpyImpl().

◆ evaluateAll() [1/3]

template<class Entity , class ShapeFunctionSet >
template<class Point , class DofVector >
void Dune::Fem::DefaultBasisFunctionSet< Entity, ShapeFunctionSet >::evaluateAll ( const Point &  x,
const DofVector &  dofs,
RangeType value 
) const
inline
Todo:
please doc me

◆ evaluateAll() [2/3]

template<class Entity , class ShapeFunctionSet >
template<class Point , class RangeArray >
void Dune::Fem::DefaultBasisFunctionSet< Entity, ShapeFunctionSet >::evaluateAll ( const Point &  x,
RangeArray &  values 
) const
inline
Todo:
please doc me

◆ evaluateAll() [3/3]

template<class Entity , class ShapeFunctionSet >
template<class QuadratureType , class DofVector , class RangeArray >
void Dune::Fem::DefaultBasisFunctionSet< Entity, ShapeFunctionSet >::evaluateAll ( const QuadratureType &  quad,
const DofVector &  dofs,
RangeArray &  ranges 
) const
inline

◆ hessianAll() [1/3]

template<class Entity , class ShapeFunctionSet >
template<class Point , class DofVector >
void Dune::Fem::DefaultBasisFunctionSet< Entity, ShapeFunctionSet >::hessianAll ( const Point &  x,
const DofVector &  dofs,
HessianRangeType hessian 
) const
inline

◆ hessianAll() [2/3]

template<class Entity , class ShapeFunctionSet >
template<class Point , class HessianRangeArray >
void Dune::Fem::DefaultBasisFunctionSet< Entity, ShapeFunctionSet >::hessianAll ( const Point &  x,
HessianRangeArray &  hessians 
) const
inline

◆ hessianAll() [3/3]

template<class Entity , class ShapeFunctionSet >
template<class QuadratureType , class DofVector , class HessianArray >
void Dune::Fem::DefaultBasisFunctionSet< Entity, ShapeFunctionSet >::hessianAll ( const QuadratureType &  quad,
const DofVector &  dofs,
HessianArray &  hessians 
) const
inline

◆ jacobianAll() [1/3]

template<class Entity , class ShapeFunctionSet >
template<class Point , class DofVector >
void Dune::Fem::DefaultBasisFunctionSet< Entity, ShapeFunctionSet >::jacobianAll ( const Point &  x,
const DofVector &  dofs,
JacobianRangeType jacobian 
) const
inline

◆ jacobianAll() [2/3]

template<class Entity , class ShapeFunctionSet >
template<class Point , class JacobianRangeArray >
void Dune::Fem::DefaultBasisFunctionSet< Entity, ShapeFunctionSet >::jacobianAll ( const Point &  x,
JacobianRangeArray &  jacobians 
) const
inline

◆ jacobianAll() [3/3]

template<class Entity , class ShapeFunctionSet >
template<class QuadratureType , class DofVector , class JacobianArray >
void Dune::Fem::DefaultBasisFunctionSet< Entity, ShapeFunctionSet >::jacobianAll ( const QuadratureType &  quad,
const DofVector &  dofs,
JacobianArray &  jacobians 
) const
inline

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 (Nov 12, 23:30, 2024)