dune-fem  2.4.1-rc
space/shapefunctionset/orthonormal.hh
Go to the documentation of this file.
1 #ifndef DUNE_FEM_SPACE_SHAPEFUNCTIONSET_ORTHONORMAL_HH
2 #define DUNE_FEM_SPACE_SHAPEFUNCTIONSET_ORTHONORMAL_HH
3 
4 #include <cassert>
5 #include <cstdlib>
6 
7 #include <dune/common/exceptions.hh>
8 
9 #include <dune/geometry/genericgeometry/topologytypes.hh>
10 #include <dune/geometry/type.hh>
11 
16 
24 namespace Dune
25 {
26 
27  namespace Fem
28  {
29 
30  // Forward declaration
31  // -------------------
32 
33  template< class FunctionSpace, int polOrder >
35 
36 
37 
38  // OrthonormalShapeFunctionSetSize
39  // -------------------------------
40 
41  template< class FunctionSpace, int polOrder >
43  {
44  static_assert( (FunctionSpace::dimDomain <= 3), "Shape function set only implemented up to dimension 3." );
45 
46  template< int order, int dimension >
47  struct ShapeFunctionSetSize;
48 
49  template< int order >
50  struct ShapeFunctionSetSize< order, 1 >
51  {
52  static const std::size_t v = order + 1;
53  };
54 
55  template< int order >
56  struct ShapeFunctionSetSize< order, 2 >
57  {
58  static const std::size_t v = (order + 2) * (order + 1) / 2;
59  };
60 
61  template< int order >
62  struct ShapeFunctionSetSize< order, 3 >
63  {
64  static const size_t v = ((order+1)*(order+2)*(2*order+3)/6
65  + (order+1)*(order+2)/2)/2;
66  };
67 
68  public:
69  static const std::size_t v = ShapeFunctionSetSize< polOrder, FunctionSpace::dimDomain >::v;
70  };
71 
72 
73 
74  // OrthonormalShapeFunctionHelper
75  // ------------------------------
76 
77  template< class FunctionSpace, int polOrder >
79  {
80  static_assert( (FunctionSpace::dimRange == 1),
81  "FunctionSpace must be scalar (i.e., dimRange = 1)." );
82 
84 
87 
92 
93  protected:
94  // line
95  typedef Dune::GenericGeometry::Prism< Dune::GenericGeometry::Point > Line;
96  // quadrilateral
97  typedef Dune::GenericGeometry::Prism< Line > Quadrilateral;
98  // triangle
99  typedef Dune::GenericGeometry::Pyramid< Line > Triangle;
100  // pyramid
101  typedef Dune::GenericGeometry::Prism< Quadrilateral > Hexahedron;
102  // hexahedron
103  typedef Dune::GenericGeometry::Pyramid< Quadrilateral > Pyramid;
104  // prism
105  typedef Dune::GenericGeometry::Prism< Triangle > Prism;
106  // tetrahedron
107  typedef Dune::GenericGeometry::Pyramid< Triangle > Tetrahedron;
108 
109  // mapping from topology to basic geometry type (see list above)
110  template< class Topology >
112  {
113  template< unsigned int topologyId >
114  struct UniqueId
115  {
116  static const unsigned int v = topologyId | (unsigned int)Dune::GenericGeometry::prismConstruction;
117  };
118 
119  public:
120  typedef typename Dune::GenericGeometry::Topology< UniqueId< Topology::id >::v, Topology::dimension >::type Type;
121  };
122 
123  // instantiate a basic geometry type
124  template< class Topology >
126  {
127  return typename BasicGeometryType< Topology >::Type();
128  }
129 
130  public:
131  template< class Topology >
132  struct EvaluateEach;
133 
134  template< class Topology >
135  struct JacobianEach;
136 
137  template< class Topology >
138  struct HessianEach;
139  };
140 
141 
142 
143  // Implementation of OrthonormalShapeFunctionHelper::EvaluateEach
144  // --------------------------------------------------------------
145 
146  template< class FunctionSpace, int polOrder >
147  template< class Topology >
149  {
150  template< class Functor >
151  static void apply ( const DomainType &x, Functor functor )
152  {
154  for( std::size_t i = 0; i < size; ++i )
155  {
156  assert( x.size() == Topology::dimension );
157  functor( i, evaluate( basicGeometryType< Topology>(), i, x ) );
158  }
159  }
160 
161 
162  protected:
166 
167  static RangeType evaluate ( const Line &, std::size_t i, const DomainType &x )
168  {
169  return OrthonormalBase1d::eval_line( i, &x[ 0 ] );
170  }
171 
172  static RangeType evaluate ( const Quadrilateral &, std::size_t i, const DomainType &x )
173  {
174  return OrthonormalBase2d::eval_quadrilateral_2d( i, &x[ 0 ] );
175  }
176 
177  static RangeType evaluate ( const Triangle &, std::size_t i, const DomainType &x )
178  {
179  return OrthonormalBase2d::eval_triangle_2d( i, &x[ 0 ] );
180  }
181 
182  static RangeType evaluate ( const Pyramid &, std::size_t i, const DomainType &x )
183  {
184  return OrthonormalBase3d::eval_pyramid_3d( i, &x[ 0 ] );
185  }
186 
187  static RangeType evaluate ( const Hexahedron &, std::size_t i, const DomainType &x )
188  {
189  return OrthonormalBase3d::eval_hexahedron_3d( i, &x[ 0 ] );
190  }
191 
192  static RangeType evaluate ( const Prism &, std::size_t i, const DomainType &x )
193  {
194  return OrthonormalBase3d::eval_prism_3d( i, &x[ 0 ] );
195  }
196 
197  static RangeType evaluate ( const Tetrahedron &, std::size_t i, const DomainType &x )
198  {
199  return OrthonormalBase3d::eval_tetrahedron_3d( i, &x[ 0 ] );
200  }
201  };
202 
203 
204 
205  // Implementation of OrthonormalShapeFunctionHelper::JacobianEach
206  // --------------------------------------------------------------
207 
208  template< class FunctionSpace, int polOrder >
209  template< class Topology >
210  struct OrthonormalShapeFunctionHelper< FunctionSpace, polOrder >::JacobianEach
211  {
212  template< class Functor >
213  static void apply ( const DomainType &x, Functor functor )
214  {
215  JacobianRangeType jacobian;
217  for( std::size_t i = 0; i < size; ++i )
218  {
219  assert( x.size() == Topology::dimension );
220  evaluate( basicGeometryType< Topology >(), i, x, jacobian );
221  functor( i, jacobian );
222  }
223  }
224 
225  protected:
229 
230  static void evaluate ( const Line &, std::size_t i, const DomainType &x, JacobianRangeType &jacobian )
231  {
232  OrthonormalBase1d::grad_line( i , &x[ 0 ], &jacobian[ 0 ][ 0 ] );
233  }
234 
235  static void evaluate ( const Quadrilateral &, std::size_t i, const DomainType &x, JacobianRangeType &jacobian )
236  {
237  OrthonormalBase2d::grad_quadrilateral_2d( i , &x[ 0 ], &jacobian[ 0 ][ 0 ] );
238  }
239 
240  static void evaluate ( const Triangle &, std::size_t i, const DomainType &x, JacobianRangeType &jacobian )
241  {
242  OrthonormalBase2d::grad_triangle_2d( i , &x[ 0 ], &jacobian[ 0 ][ 0 ] );
243  }
244 
245  static void evaluate ( const Pyramid &, std::size_t i, const DomainType &x, JacobianRangeType &jacobian )
246  {
247  OrthonormalBase3d::grad_pyramid_3d( i , &x[ 0 ], &jacobian[ 0 ][ 0 ] );
248  }
249 
250  static void evaluate ( const Hexahedron &, std::size_t i, const DomainType &x, JacobianRangeType &jacobian )
251  {
252  OrthonormalBase3d::grad_hexahedron_3d( i , &x[ 0 ], &jacobian[ 0 ][ 0 ] );
253  }
254 
255  static void evaluate ( const Prism &, std::size_t i, const DomainType &x, JacobianRangeType &jacobian )
256  {
257  OrthonormalBase3d::grad_prism_3d( i , &x[ 0 ], &jacobian[ 0 ][ 0 ] );
258  }
259 
260  static void evaluate ( const Tetrahedron &, std::size_t i, const DomainType &x, JacobianRangeType &jacobian )
261  {
262  OrthonormalBase3d::grad_tetrahedron_3d( i , &x[ 0 ], &jacobian[ 0 ][ 0 ] );
263  }
264  };
265 
266 
267 
268  // Implementation of OrthonormalShapeFunctionHelper::HessianEach
269  // -------------------------------------------------------------
270 
271  template< class FunctionSpace, int polOrder >
272  template< class Topology >
273  struct OrthonormalShapeFunctionHelper< FunctionSpace, polOrder >::HessianEach
274  {
275  template< class Functor >
276  static void apply ( const DomainType &x, Functor functor )
277  {
278  HessianRangeType hessian;
280  for( std::size_t i = 0; i < size; ++i )
281  {
282  assert( x.size() == Topology::dimension );
283  evaluate( basicGeometryType< Topology >(), i, x, hessian );
284  functor( i, hessian );
285  }
286  }
287 
288  protected:
289  static void evaluate ( const Quadrilateral &, std::size_t i, const DomainType &x, HessianRangeType &hessian )
290  {
292  RangeFieldType values[] = { 0, 0, 0 };
293  OrthonormalBase2d::hess_quadrilateral_2d( i , &x[ 0 ], values );
294 
295  for( unsigned int j = 0; j < FunctionSpace::dimDomain; ++j )
296  for( unsigned int k = 0; k < FunctionSpace::dimDomain; ++k )
297  hessian[ 0 ][ j ][ k ] = values[ j + k ];
298  }
299 
300  static void evaluate ( const Triangle &, std::size_t i, const DomainType &x, HessianRangeType &hessian )
301  {
303  RangeFieldType values[] = { 0, 0, 0 };
304  OrthonormalBase2d::hess_triangle_2d( i , &x[ 0 ], values );
305 
306  for( unsigned int j = 0; j < FunctionSpace::dimDomain; ++j )
307  for( unsigned int k = 0; k < FunctionSpace::dimDomain; ++k )
308  hessian[ 0 ][ j ][ k ] = values[ j + k ];
309  }
310 
311  template< class Other >
312  static void evaluate ( const Other &, std::size_t i, const DomainType &x, HessianRangeType &hessian )
313  {
314  DUNE_THROW( NotImplemented, "On orthonormal shape function set HessianAll() is only implemented for triangles" );
315  }
316  };
317 
318 
319 
320  // OrthonormalShapeFunctionSet
321  // ---------------------------
322 
323  template< class FunctionSpace, int polOrder >
325  {
326  static_assert( (FunctionSpace::dimRange == 1),
327  "FunctionSpace must be scalar (i.e., dimRange = 1)." );
328 
329  // this type
331 
332  // helper class
333  typedef OrthonormalShapeFunctionHelper< FunctionSpace, polOrder > ShapeFunctionSetHelperType;
334 
335  public:
337  typedef FunctionSpace FunctionSpaceType;
338 
340  static const int dimension = FunctionSpaceType::dimDomain;
341 
350 
351  explicit OrthonormalShapeFunctionSet ( const GeometryType &type )
352  : topologyId_( type.id() )
353  {}
354 
356  int order () const { return polOrder; }
357 
358 
360  std::size_t size () const
361  {
363  }
364 
366  template< class Point, class Functor >
367  void evaluateEach ( const Point &x, Functor functor ) const
368  {
369  const DomainType y = coordinate( x );
370  Dune::GenericGeometry::IfTopology< ShapeFunctionSetHelperType::template EvaluateEach, dimension >
371  ::apply( topologyId_, y, functor );
372  }
373 
375  template< class Point, class Functor >
376  void jacobianEach ( const Point &x, Functor functor ) const
377  {
378  const DomainType y = coordinate( x );
379  Dune::GenericGeometry::IfTopology< ShapeFunctionSetHelperType::template JacobianEach, dimension >
380  ::apply( topologyId_, y, functor );
381  }
382 
384  template< class Point, class Functor >
385  void hessianEach ( const Point &x, Functor functor ) const
386  {
387  const DomainType y = coordinate( x );
388  Dune::GenericGeometry::IfTopology< ShapeFunctionSetHelperType::template HessianEach, dimension >
389  ::apply( topologyId_, y, functor );
390  }
391 
392  private:
393  unsigned int topologyId_;
394  };
395 
396  } // namespace Fem
397 
398 } // namespace Dune
399 
400 #endif // #ifndef DUNE_FEM_SPACE_SHAPEFUNCTIONSET_ORTHONORMAL_HH
static RangeType evaluate(const Hexahedron &, std::size_t i, const DomainType &x)
Definition: space/shapefunctionset/orthonormal.hh:187
void jacobianEach(const Point &x, Functor functor) const
evalute jacobian of each shape function
Definition: space/shapefunctionset/orthonormal.hh:376
OrthonormalBase_1D< DomainFieldType, RangeFieldType > OrthonormalBase1d
Definition: space/shapefunctionset/orthonormal.hh:226
static BasicGeometryType< Topology >::Type basicGeometryType()
Definition: space/shapefunctionset/orthonormal.hh:125
VectorSpaceTraits< DomainField, RangeField, dimD, dimR >::RangeFieldType RangeFieldType
Intrinsic type used for values in the range field (usually a double)
Definition: functionspaceinterface.hh:62
Definition: space/shapefunctionset/orthonormal.hh:138
static void apply(const DomainType &x, Functor functor)
Definition: space/shapefunctionset/orthonormal.hh:151
OrthonormalBase_2D< DomainFieldType, RangeFieldType > OrthonormalBase2d
Definition: space/shapefunctionset/orthonormal.hh:227
A vector valued function space.
Definition: functionspace.hh:16
Dune::GenericGeometry::Prism< Dune::GenericGeometry::Point > Line
Definition: space/shapefunctionset/orthonormal.hh:95
Dune::GenericGeometry::Prism< Line > Quadrilateral
Definition: space/shapefunctionset/orthonormal.hh:97
Definition: space/shapefunctionset/orthonormal.hh:42
static RangeType evaluate(const Tetrahedron &, std::size_t i, const DomainType &x)
Definition: space/shapefunctionset/orthonormal.hh:197
static RangeType evaluate(const Line &, std::size_t i, const DomainType &x)
Definition: space/shapefunctionset/orthonormal.hh:167
VectorSpaceTraits< DomainField, RangeField, dimD, dimR >::LinearMappingType JacobianRangeType
Intrinsic type used for the jacobian values has a Dune::FieldMatrix type interface.
Definition: functionspaceinterface.hh:74
Definition: space/shapefunctionset/orthonormal.hh:78
FieldVector< FieldMatrix< RangeFieldType, dimDomain, dimDomain >, dimRange > HessianRangeType
Intrinsic type used for the hessian values has a Dune::FieldMatrix type interface.
Definition: functionspaceinterface.hh:78
std::size_t size() const
return number of shape functions
Definition: space/shapefunctionset/orthonormal.hh:360
OrthonormalBase_1D< DomainFieldType, RangeFieldType > OrthonormalBase1d
Definition: space/shapefunctionset/orthonormal.hh:163
FunctionSpaceType::JacobianRangeType JacobianRangeType
Definition: space/shapefunctionset/orthonormal.hh:90
FunctionSpaceType::RangeFieldType RangeFieldType
Definition: space/shapefunctionset/orthonormal.hh:86
FunctionSpaceType::RangeType RangeType
range type
Definition: space/shapefunctionset/orthonormal.hh:345
Definition: orthonormalbase_2d.hh:16
dimension of domain vector space
Definition: functionspaceinterface.hh:45
FunctionSpaceType::DomainType DomainType
domain type
Definition: space/shapefunctionset/orthonormal.hh:343
static RangeType evaluate(const Quadrilateral &, std::size_t i, const DomainType &x)
Definition: space/shapefunctionset/orthonormal.hh:172
Dune::GenericGeometry::Pyramid< Line > Triangle
Definition: space/shapefunctionset/orthonormal.hh:99
dimension of range vector space
Definition: functionspaceinterface.hh:47
void evaluateEach(const Point &x, Functor functor) const
evalute each shape function
Definition: space/shapefunctionset/orthonormal.hh:367
VectorSpaceTraits< DomainField, RangeField, dimD, dimR >::DomainFieldType DomainFieldType
Intrinsic type used for values in the domain field (usually a double)
Definition: functionspaceinterface.hh:59
static void apply(const DomainType &x, Functor functor)
Definition: space/shapefunctionset/orthonormal.hh:276
static const std::size_t v
Definition: space/shapefunctionset/orthonormal.hh:69
Dune::GenericGeometry::Topology< UniqueId< Topology::id >::v, Topology::dimension >::type Type
Definition: space/shapefunctionset/orthonormal.hh:120
static void evaluate(const Pyramid &, std::size_t i, const DomainType &x, JacobianRangeType &jacobian)
Definition: space/shapefunctionset/orthonormal.hh:245
Definition: orthonormalbase_1d.hh:14
FunctionSpaceType::DomainFieldType DomainFieldType
Definition: space/shapefunctionset/orthonormal.hh:85
int order() const
return order of shape functions
Definition: space/shapefunctionset/orthonormal.hh:356
Definition: coordinate.hh:4
FunctionSpace FunctionSpaceType
Definition: space/shapefunctionset/orthonormal.hh:81
OrthonormalBase_2D< DomainFieldType, RangeFieldType > OrthonormalBase2d
Definition: space/shapefunctionset/orthonormal.hh:164
FunctionSpaceType::HessianRangeType HessianRangeType
Definition: space/shapefunctionset/orthonormal.hh:91
VectorSpaceTraits< DomainField, RangeField, dimD, dimR >::DomainType DomainType
Type of domain vector (using type of domain field) has a Dune::FieldVector type interface.
Definition: functionspaceinterface.hh:66
FunctionSpaceType::JacobianRangeType JacobianRangeType
jacobian range type
Definition: space/shapefunctionset/orthonormal.hh:347
Definition: space/shapefunctionset/orthonormal.hh:111
OrthonormalBase_3D< DomainFieldType, RangeFieldType > OrthonormalBase3d
Definition: space/shapefunctionset/orthonormal.hh:165
static RangeType evaluate(const Prism &, std::size_t i, const DomainType &x)
Definition: space/shapefunctionset/orthonormal.hh:192
FunctionSpaceType::RangeType RangeType
Definition: space/shapefunctionset/orthonormal.hh:89
static void evaluate(const Line &, std::size_t i, const DomainType &x, JacobianRangeType &jacobian)
Definition: space/shapefunctionset/orthonormal.hh:230
static void evaluate(const Triangle &, std::size_t i, const DomainType &x, HessianRangeType &hessian)
Definition: space/shapefunctionset/orthonormal.hh:300
OrthonormalShapeFunctionSet(const GeometryType &type)
Definition: space/shapefunctionset/orthonormal.hh:351
static RangeType evaluate(const Triangle &, std::size_t i, const DomainType &x)
Definition: space/shapefunctionset/orthonormal.hh:177
Dune::GenericGeometry::Prism< Quadrilateral > Hexahedron
Definition: space/shapefunctionset/orthonormal.hh:101
Definition: space/shapefunctionset/orthonormal.hh:132
Dune::GenericGeometry::Pyramid< Triangle > Tetrahedron
Definition: space/shapefunctionset/orthonormal.hh:107
static void apply(const DomainType &x, Functor functor)
Definition: space/shapefunctionset/orthonormal.hh:213
static RangeType evaluate(const Pyramid &, std::size_t i, const DomainType &x)
Definition: space/shapefunctionset/orthonormal.hh:182
VectorSpaceTraits< DomainField, RangeField, dimD, dimR >::RangeType RangeType
Type of range vector (using type of range field) has a Dune::FieldVector type interface.
Definition: functionspaceinterface.hh:70
FunctionSpaceType::HessianRangeType HessianRangeType
hessian range type
Definition: space/shapefunctionset/orthonormal.hh:349
Dune::GenericGeometry::Pyramid< Quadrilateral > Pyramid
Definition: space/shapefunctionset/orthonormal.hh:103
FunctionSpace FunctionSpaceType
function space type
Definition: space/shapefunctionset/orthonormal.hh:337
FunctionSpaceType::DomainType DomainType
Definition: space/shapefunctionset/orthonormal.hh:88
static void evaluate(const Triangle &, std::size_t i, const DomainType &x, JacobianRangeType &jacobian)
Definition: space/shapefunctionset/orthonormal.hh:240
static void evaluate(const Other &, std::size_t i, const DomainType &x, HessianRangeType &hessian)
Definition: space/shapefunctionset/orthonormal.hh:312
Definition: space/shapefunctionset/orthonormal.hh:135
static void evaluate(const Tetrahedron &, std::size_t i, const DomainType &x, JacobianRangeType &jacobian)
Definition: space/shapefunctionset/orthonormal.hh:260
static void evaluate(const Quadrilateral &, std::size_t i, const DomainType &x, HessianRangeType &hessian)
Definition: space/shapefunctionset/orthonormal.hh:289
Definition: orthonormalbase_3d.hh:16
static void evaluate(const Quadrilateral &, std::size_t i, const DomainType &x, JacobianRangeType &jacobian)
Definition: space/shapefunctionset/orthonormal.hh:235
Dune::GenericGeometry::Prism< Triangle > Prism
Definition: space/shapefunctionset/orthonormal.hh:105
static const Point & coordinate(const Point &x)
Definition: coordinate.hh:11
static void evaluate(const Hexahedron &, std::size_t i, const DomainType &x, JacobianRangeType &jacobian)
Definition: space/shapefunctionset/orthonormal.hh:250
OrthonormalBase_3D< DomainFieldType, RangeFieldType > OrthonormalBase3d
Definition: space/shapefunctionset/orthonormal.hh:228
void hessianEach(const Point &x, Functor functor) const
evalute hessian of each shape function
Definition: space/shapefunctionset/orthonormal.hh:385
static void evaluate(const Prism &, std::size_t i, const DomainType &x, JacobianRangeType &jacobian)
Definition: space/shapefunctionset/orthonormal.hh:255
Definition: space/shapefunctionset/orthonormal.hh:34