1#ifndef DUNE_FEM_SPACE_SHAPEFUNCTIONSET_ORTHONORMAL_HH
2#define DUNE_FEM_SPACE_SHAPEFUNCTIONSET_ORTHONORMAL_HH
10#include <dune/fem/common/coordinate.hh>
11#include <dune/fem/quadrature/quadrature.hh>
13#include <dune/fem/space/common/functionspace.hh>
14#include <dune/fem/space/shapefunctionset/orthonormal/orthonormalbase_1d.hh>
15#include <dune/fem/space/shapefunctionset/orthonormal/orthonormalbase_2d.hh>
16#include <dune/fem/space/shapefunctionset/orthonormal/orthonormalbase_3d.hh>
18#if HAVE_DUNE_LOCALFUNCTIONS
19#include <dune/localfunctions/orthonormal.hh>
33 template<
int dimension >
34 struct OrthonormalShapeFunctions;
37 struct OrthonormalShapeFunctions< 1 >
39 static constexpr std::size_t
size (
int order )
41 return static_cast< std::size_t
>( order+1 );
46 struct OrthonormalShapeFunctions< 2 >
48 static constexpr std::size_t
size (
int order )
50 return static_cast< std::size_t
>( (order+2)*(order+1)/2 );
55 struct OrthonormalShapeFunctions< 3 >
57 static constexpr std::size_t
size (
int order )
59 return static_cast< std::size_t
>( ((order+1)*(order+2)*(2*order+3)/6+(order+1)*(order+2)/2)/2 );
70 template<
class FunctionSpace >
71 class OrthonormalShapeFunctionSet
73 typedef OrthonormalShapeFunctionSet< FunctionSpace > ThisType;
75 static_assert(
FunctionSpace::dimDomain <= 3,
"OrthonormalShapeFunctionSet only implemented up to domain dimension 3" );
76 static_assert(
FunctionSpace::dimRange == 1,
"OrthonormalShapeFunctionSet only implemented for scalar function spaces" );
86 typedef RangeFieldType (*Evaluate) (
const int,
const DomainFieldType * );
87 typedef void (*Jacobian) (
const int i,
const DomainFieldType *, RangeFieldType * );
90 typedef RangeFieldType Array[ 3 ];
91 typedef void (*Hessian) (
const int i,
const DomainFieldType *, Array & );
93 typedef OrthonormalBase_1D< DomainFieldType, RangeFieldType > OrthonormalBase1d;
94 typedef OrthonormalBase_2D< DomainFieldType, RangeFieldType > OrthonormalBase2d;
95 typedef OrthonormalBase_3D< DomainFieldType, RangeFieldType > OrthonormalBase3d;
98 Evaluate &evaluate, Jacobian &jacobian )
102 evaluate = OrthonormalBase1d::eval_line;
103 jacobian = OrthonormalBase1d::grad_line;
108 evaluate = OrthonormalBase2d::eval_quadrilateral_2d;
109 jacobian = OrthonormalBase2d::grad_quadrilateral_2d;
114 evaluate = OrthonormalBase2d::eval_triangle_2d;
115 jacobian = OrthonormalBase2d::grad_triangle_2d;
120 evaluate = OrthonormalBase3d::eval_hexahedron_3d;
121 jacobian = OrthonormalBase3d::grad_hexahedron_3d;
126 evaluate = OrthonormalBase3d::eval_tetrahedron_3d;
127 jacobian = OrthonormalBase3d::grad_tetrahedron_3d;
132 evaluate = OrthonormalBase3d::eval_prism_3d;
133 jacobian = OrthonormalBase3d::grad_prism_3d;
138 evaluate = OrthonormalBase3d::eval_pyramid_3d;
139 jacobian = OrthonormalBase3d::grad_pyramid_3d;
143 DUNE_THROW(InvalidStateException,
"Invalid geometry type " << geomType );
159#if HAVE_DUNE_LOCALFUNCTIONS
160 typedef OrthonormalLocalFiniteElement< dimension,
161 typename DomainType::field_type,
typename RangeType::field_type> OrthonormalLocalFiniteElementType;
162 static const bool haveLocalFunctions_ =
true;
164 typedef int OrthonormalLocalFiniteElementType;
165 static const bool haveLocalFunctions_ =
false;
171 OrthonormalShapeFunctionSet () =
default;
173 OrthonormalShapeFunctionSet ( GeometryType type,
int order )
175 evaluate_( nullptr ),
176 jacobian_( nullptr ),
186 setFunctionPointers( type, evaluate_, jacobian_ );
190 if constexpr ( dimension == 2 )
192 if( type.isTriangle() )
193 hessian_ = OrthonormalBase2d::hess_triangle_2d;
194 else if( type.isQuadrilateral() )
195 hessian_ = OrthonormalBase2d::hess_quadrilateral_2d;
197 else if constexpr ( dimension == 3 && haveLocalFunctions_ )
198 onbElem_.reset(
new OrthonormalLocalFiniteElementType( type, order ) );
202 OrthonormalShapeFunctionSet (
const ThisType & ) =
default;
204 OrthonormalShapeFunctionSet ( ThisType && ) =
default;
212 OrthonormalShapeFunctionSet &operator= (
const ThisType & ) =
default;
214 OrthonormalShapeFunctionSet &operator= ( ThisType && ) =
default;
223 int order ()
const {
return order_; }
226 std::size_t
constexpr size ()
const {
return size( order() ); }
229 static std::size_t
constexpr size (
int order )
231 return OrthonormalShapeFunctions< dimension >::size( order );
235 template<
class Po
int,
class Functor >
236 void evaluateEach (
const Point &x, Functor functor )
const
238 const DomainType y = Dune::Fem::coordinate( x );
240 const std::size_t
size = this->
size();
241 for( std::size_t i = 0; i <
size; ++i )
243 evaluate( i, y, value );
249 template<
class Po
int,
class Functor >
250 void jacobianEach (
const Point &x, Functor functor )
const
252 const DomainType y = Dune::Fem::coordinate( x );
253 JacobianRangeType jacobian;
254 const std::size_t
size = this->
size();
255 for( std::size_t i = 0; i <
size; ++i )
257 this->jacobian( i, y, jacobian );
258 functor( i, jacobian );
263 template<
class Po
int,
class Functor >
264 void hessianEach (
const Point &x, Functor functor )
const
266 const DomainType y = Dune::Fem::coordinate( x );
267 if constexpr ( dimension == 2 )
269 HessianRangeType hessian;
270 const std::size_t
size = this->
size();
271 for( std::size_t i = 0; i <
size; ++i )
273 this->hessian( i, y, hessian );
274 functor( i, hessian );
277 else if constexpr( dimension == 3 && haveLocalFunctions_ )
279 const std::size_t
size = this->
size();
283 HessianRangeType hessian(0);
284 for( std::size_t i = 0; i <
size; ++i )
288 functor( i, hessian );
293 typedef typename OrthonormalLocalFiniteElementType::Traits::LocalBasisType::HessianType HessianType;
294 std::vector< HessianType > hessians;
295 onbElem_->localBasis().evaluateHessian( y, hessians );
296 HessianRangeType hessian;
297 for( std::size_t i = 0; i <
size; ++i )
301 hessian = hessians[ i ];
302 functor( i, hessian );
314 void evaluate ( std::size_t i,
const DomainType &x, RangeType &value )
const
316 value[ 0 ] = evaluate_( i, &x[ 0 ] );
319 void jacobian ( std::size_t i,
const DomainType &x, JacobianRangeType &jacobian )
const
321 jacobian_( i, &x[ 0 ], &jacobian[ 0 ][ 0 ] );
324 void hessian ( std::size_t i,
const DomainType &x, HessianRangeType &hessian )
const
327 if constexpr (dimension == 2 )
329 RangeFieldType values[] = { 0, 0, 0 };
330 hessian_( i , &x[ 0 ], values );
333 hessian[ 0 ][ j ][ k ] = values[ j + k ];
346 mutable std::shared_ptr< OrthonormalLocalFiniteElementType > onbElem_;
Definition: explicitfieldvector.hh:75
FunctionSpaceTraits::DomainFieldType DomainFieldType
Intrinsic type used for values in the domain field (usually a double)
Definition: functionspaceinterface.hh:60
FunctionSpaceTraits::RangeType RangeType
Type of range vector (using type of range field) has a Dune::FieldVector type interface.
Definition: functionspaceinterface.hh:71
@ dimDomain
dimension of domain vector space
Definition: functionspaceinterface.hh:46
@ dimRange
dimension of range vector space
Definition: functionspaceinterface.hh:48
FunctionSpaceTraits::LinearMappingType JacobianRangeType
Intrinsic type used for the jacobian values has a Dune::FieldMatrix type interface.
Definition: functionspaceinterface.hh:75
FunctionSpaceTraits::DomainType DomainType
Type of domain vector (using type of domain field) has a Dune::FieldVector type interface.
Definition: functionspaceinterface.hh:67
FunctionSpaceTraits::RangeFieldType RangeFieldType
Intrinsic type used for values in the range field (usually a double)
Definition: functionspaceinterface.hh:63
A vector valued function space.
Definition: functionspace.hh:60
Unique label for each type of entities that can occur in DUNE grids.
Definition: type.hh:114
constexpr bool isPyramid() const
Return true if entity is a pyramid.
Definition: type.hh:304
constexpr bool isTetrahedron() const
Return true if entity is a tetrahedron.
Definition: type.hh:299
constexpr bool isPrism() const
Return true if entity is a prism.
Definition: type.hh:309
constexpr bool isTriangle() const
Return true if entity is a triangle.
Definition: type.hh:289
constexpr bool isLine() const
Return true if entity is a line segment.
Definition: type.hh:284
constexpr bool isQuadrilateral() const
Return true if entity is a quadrilateral.
Definition: type.hh:294
constexpr bool isHexahedron() const
Return true if entity is a hexahedron.
Definition: type.hh:314
#define DUNE_THROW(E,...)
Definition: exceptions.hh:314
constexpr GeometryType cube(unsigned int dim)
Returns a GeometryType representing a hypercube of dimension dim.
Definition: type.hh:462
Dune namespace.
Definition: alignedallocator.hh:13
constexpr std::integral_constant< std::size_t, sizeof...(II)> size(std::integer_sequence< T, II... >)
Return the size of the sequence.
Definition: integersequence.hh:75
A unique label for each type of element that can occur in a grid.