DUNE-FEM (unstable)

vectorial.hh
1#ifndef DUNE_FEM_BASISFUNCTIONSET_VECTORIAL_HH
2#define DUNE_FEM_BASISFUNCTIONSET_VECTORIAL_HH
3
4#include <cstddef>
5#include <utility>
6#include <vector>
7
8#include <dune/geometry/referenceelements.hh>
10
11#include <dune/fem/space/basisfunctionset/functor.hh>
12#include <dune/fem/space/common/functionspace.hh>
13
14namespace Dune
15{
16
17 namespace Fem
18 {
19
20 // DofAlignment
21 // ------------
22
32 template< class Implementation >
34 {
35 public:
37 typedef std::size_t GlobalDofType;
41 typedef std::pair< int, std::size_t > LocalDofType;
42
43 protected:
44 DofAlignment () = default;
45
46 public:
56 {
57 return impl().globalDof( localDof );
58 }
59
69 {
70 return impl().localDof( globalDof );
71 }
72
73 protected:
74 const Implementation &impl () const
75 {
76 return static_cast< const Implementation & >( *this );
77 }
78 };
79
80
81
82 // HorizontalDofAlignment
83 // ----------------------
84
93 template< class ScalarBasisFunctionSet, class Range >
95 : public DofAlignment< HorizontalDofAlignment< ScalarBasisFunctionSet, Range > >
96 {
99
100 public:
101 typedef typename BaseType::GlobalDofType GlobalDofType;
102 typedef typename BaseType::LocalDofType LocalDofType;
103
104 HorizontalDofAlignment () = default;
105
106 explicit HorizontalDofAlignment ( const ScalarBasisFunctionSet &scalarBasisFunctionSet )
107 : scalarSize_( scalarBasisFunctionSet.size() )
108 {}
109
111 GlobalDofType globalDof ( const LocalDofType &localDof ) const
112 {
113 return GlobalDofType( localDof.first*scalarSize_ + localDof.second );
114 }
115
117 LocalDofType localDof ( const GlobalDofType &globalDof ) const
118 {
119 return LocalDofType( globalDof / scalarSize_, globalDof % scalarSize_ );
120 }
121
122 private:
123 std::size_t scalarSize_;
124 };
125
126
127
128 // VerticalDofAlignment
129 // --------------------
130
138 template< class ScalarBasisFunctionSet, class Range >
140 : public DofAlignment< VerticalDofAlignment< ScalarBasisFunctionSet, Range > >
141 {
144
145 static const int dimRange = Range::dimension;
146
147 public:
148 typedef typename BaseType::GlobalDofType GlobalDofType;
149 typedef typename BaseType::LocalDofType LocalDofType;
150
151 VerticalDofAlignment () = default;
152
153 explicit VerticalDofAlignment ( const ScalarBasisFunctionSet & ) {}
154
156 GlobalDofType globalDof ( const LocalDofType &localDof ) const
157 {
158 return GlobalDofType( localDof.first + localDof.second*dimRange );
159 }
160
162 LocalDofType localDof ( const GlobalDofType &globalDof ) const
163 {
164 return LocalDofType( globalDof % dimRange, globalDof / dimRange );
165 }
166 };
167
168
169
170 // SubDofVector
171 // ------------
172
178 template< class DofVector, class DofAlignment >
180
181 // note: DofVector can be also const DofVector
182 template< class DofVector, class ScalarBasisFunctionSet, class Range >
183 class SubDofVector< DofVector, HorizontalDofAlignment< ScalarBasisFunctionSet, Range > >
184 {
186
187 typedef typename Range::value_type RangeFieldType;
188
190 typedef typename DofAlignmentType::LocalDofType LocalDofType;
191
192 // extract correct RangeFieldType for const or non-const version
193 typedef typename std::conditional<
194 std::is_const< DofVector > :: value,
195 const RangeFieldType,
196 RangeFieldType > :: type DofType;
197
198 public:
199 typedef RangeFieldType value_type;
200
201 SubDofVector( DofVector &dofs, int coordinate, const DofAlignmentType &dofAlignment )
202 : dofs_( &(dofs[ dofAlignment.globalDof( LocalDofType( coordinate, 0 ) ) ] ) )
203 {}
204
205 DofType &operator[] ( std::size_t i )
206 {
207 return dofs_[ i ];
208 }
209
210 // const RFT& leads to problem with returning reference to temporary
211 // with dofs_ = PetscDF::LocalFunction::DofVector
212 RangeFieldType operator[] ( std::size_t i ) const
213 {
214 return dofs_[ i ];
215 }
216
217 private:
218 DofType *dofs_;
219 };
220
221 // note: DofVector can be also const DofVector
222 template< class DofVector, class ScalarBasisFunctionSet, class Range >
223 class SubDofVector< DofVector, VerticalDofAlignment< ScalarBasisFunctionSet, Range > >
224 {
225 typedef SubDofVector< DofVector, VerticalDofAlignment< ScalarBasisFunctionSet, Range > > ThisType;
226
227 typedef typename Range::value_type RangeFieldType;
228
229 typedef VerticalDofAlignment< ScalarBasisFunctionSet, Range > DofAlignmentType;
230 typedef typename DofAlignmentType::LocalDofType LocalDofType;
231
232 // extract correct RangeFieldType for const or non-const version
233 typedef typename std::conditional<
234 std::is_const< DofVector > :: value,
235 const RangeFieldType,
236 RangeFieldType > :: type DofType;
237
238 public:
239 typedef RangeFieldType value_type;
240
241 SubDofVector( DofVector &dofs, int coordinate, const DofAlignmentType &dofAlignment )
242 : dofs_( dofs ),
243 coordinate_( coordinate ),
244 dofAlignment_( dofAlignment )
245 {}
246
247 DofType &operator[] ( std::size_t i )
248 {
249 return dofs_[ dofAlignment_.globalDof( LocalDofType( coordinate_, i ) ) ];
250 }
251
252 // const RFT& leads to problem with returning reference to temporary
253 // with dofs_ = PetscDF::LocalFunction::DofVector
254 RangeFieldType operator[] ( std::size_t i ) const
255 {
256 return dofs_[ dofAlignment_.globalDof( LocalDofType( coordinate_, i ) ) ];
257 }
258
259 private:
260 DofVector &dofs_;
261 int coordinate_;
262 DofAlignmentType dofAlignment_;
263 };
264
265
266
267 // VectorialBasisFunctionSet
268 // -------------------------
269
277 template< class ScalarBasisFunctionSet, class Range, template< class, class > class DofAlignment = VerticalDofAlignment >
279 {
281
282 public:
283 typedef ScalarBasisFunctionSet ScalarBasisFunctionSetType;
284
285 typedef typename ScalarBasisFunctionSetType::EntityType EntityType;
286 typedef typename ScalarBasisFunctionSetType::ReferenceElementType ReferenceElementType;
287
288 private:
289 typedef typename ScalarBasisFunctionSetType::FunctionSpaceType ScalarFunctionSpaceType;
290 static const int dimRange = Range::dimension;
291
292 public:
294
295 typedef typename FunctionSpaceType::DomainType DomainType;
296 typedef typename FunctionSpaceType::RangeType RangeType;
297 typedef typename FunctionSpaceType::JacobianRangeType JacobianRangeType;
299
301
302 typedef typename DofAlignmentType::GlobalDofType GlobalDofType;
303 typedef typename DofAlignmentType::LocalDofType LocalDofType;
304
305 typedef typename FunctionSpaceType::DomainFieldType DomainFieldType;
306 typedef typename FunctionSpaceType::RangeFieldType RangeFieldType;
307
308 private:
309 struct EvaluateAll;
310 struct JacobianAll;
311 struct HessianAll;
312
313 public:
315
316 explicit VectorialBasisFunctionSet ( const ScalarBasisFunctionSetType &scalarBasisFunctionSet )
317 : scalarBasisFunctionSet_( scalarBasisFunctionSet ),
318 dofAlignment_( scalarBasisFunctionSet_ )
319 {}
320
321 int order () const { return scalarBasisFunctionSet().order(); }
322
323 std::size_t size () const { return dimRange*scalarBasisFunctionSet().size(); }
324
325 const ReferenceElementType &referenceElement () const { return scalarBasisFunctionSet().referenceElement(); }
326
327 template< class Point, class DofVector >
328 void axpy ( const Point &x, const RangeType &valueFactor, DofVector &dofs ) const
329 {
330 axpy< EvaluateAll >( x, valueFactor, dofs );
331 }
332
333 template< class Point, class DofVector >
334 void axpy ( const Point &x, const JacobianRangeType &jacobianFactor, DofVector &dofs ) const
335 {
336 axpy< JacobianAll >( x, jacobianFactor, dofs );
337 }
338
339 template< class Point, class DofVector >
340 void axpy ( const Point &x, const HessianRangeType &hessianFactor, DofVector &dofs ) const
341 {
342 axpyH( x, hessianFactor, dofs );
343 }
344
345 template< class Point, class DofVector >
346 void axpy ( const Point &x, const RangeType &valueFactor, const JacobianRangeType &jacobianFactor,
347 DofVector &dofs ) const
348 {
349 axpy( x, valueFactor, dofs );
350 axpy( x, jacobianFactor, dofs );
351 }
352
353 template< class Quadrature, class Vector, class DofVector >
354 void axpy ( const Quadrature &quad, const Vector &values, DofVector & dofs ) const
355 {
356 const unsigned int nop = quad.nop();
357 for( unsigned int qp = 0; qp < nop ; ++qp )
358 axpy( quad[ qp ], values[ qp ], dofs );
359 }
360
361 template< class Quadrature, class VectorA, class VectorB, class DofVector >
362 void axpy ( const Quadrature &quad, const VectorA &valuesA, const VectorB &valuesB, DofVector & dofs ) const
363 {
364 const unsigned int nop = quad.nop();
365 for( unsigned int qp = 0; qp < nop ; ++qp )
366 {
367 axpy( quad[ qp ], valuesA[ qp ], dofs );
368 axpy( quad[ qp ], valuesB[ qp ], dofs );
369 }
370 }
371
372 template< class Point, class DofVector >
373 void evaluateAll ( const Point &x, const DofVector &dofs, RangeType &value ) const
374 {
375 evaluateAll< EvaluateAll >( x, dofs, value );
376 }
377
378 template< class Point, class RangeArray >
379 void evaluateAll ( const Point &x, RangeArray &values ) const
380 {
381 evaluateAll< EvaluateAll >( x, values );
382 }
383
384 template< class Quadrature, class DofVector, class RangeArray >
385 void evaluateAll ( const Quadrature &quad, const DofVector &dofs, RangeArray &ranges ) const
386 {
387 const unsigned int nop = quad.nop();
388 for( unsigned int qp = 0; qp < nop ; ++qp )
389 evaluateAll( quad[ qp ], dofs, ranges[ qp ] );
390 }
391
392 template< class Point, class DofVector >
393 void jacobianAll ( const Point &x, const DofVector &dofs, JacobianRangeType &jacobian ) const
394 {
395 evaluateAll< JacobianAll >( x, dofs, jacobian );
396 }
397
398 template< class Point, class JacobianRangeArray >
399 void jacobianAll ( const Point &x, JacobianRangeArray &jacobians ) const
400 {
401 evaluateAll< JacobianAll >( x, jacobians );
402 }
403
404 template< class Quadrature, class DofVector, class JacobianArray >
405 void jacobianAll ( const Quadrature &quad, const DofVector &dofs, JacobianArray &jacobians ) const
406 {
407 const unsigned int nop = quad.nop();
408 for( unsigned int qp = 0; qp < nop ; ++qp )
409 jacobianAll( quad[ qp ], dofs, jacobians[ qp ] );
410 }
411
412 template< class Point, class DofVector >
413 void hessianAll ( const Point &x, const DofVector &dofs, HessianRangeType &hessian ) const
414 {
415 evaluateAll< HessianAll >( x, dofs, hessian );
416 }
417
418 template< class Point, class HessianRangeArray >
419 void hessianAll ( const Point &x, HessianRangeArray &hessians ) const
420 {
421 evaluateAll< HessianAll >( x, hessians );
422 }
423
424 template< class Quadrature, class DofVector, class HessianArray >
425 void hessianAll ( const Quadrature &quad, const DofVector &dofs, HessianArray &hessians ) const
426 {
427 const unsigned int nop = quad.nop();
428 for( unsigned int qp = 0; qp < nop ; ++qp )
429 hessianAll( quad[ qp ], dofs, hessians[ qp ] );
430 }
431
432 const EntityType &entity () const { return scalarBasisFunctionSet().entity(); }
433 bool valid () const { return scalarBasisFunctionSet().valid(); }
434
435 DofAlignmentType dofAlignment () const { return dofAlignment_ ; }
436
437 const ScalarBasisFunctionSetType &scalarBasisFunctionSet () const
438 {
439 return scalarBasisFunctionSet_;
440 }
441
442 private:
443 template< class Evaluate, class Point, class DofVector >
444 void axpy ( const Point &x, const typename Evaluate::Vector &factor, DofVector &dofs ) const
445 {
446 const std::size_t size = scalarBasisFunctionSet().size();
447 std::vector< typename Evaluate::Scalar > scalars( size );
448 Evaluate::apply( scalarBasisFunctionSet(), x, scalars );
449
450 for( int r = 0; r < dimRange; ++r )
451 {
452 for( std::size_t i = 0; i < size; ++i )
453 {
454 const GlobalDofType globalDof = dofAlignment_.globalDof( LocalDofType( r, i ) );
455 dofs[ globalDof ] += factor[ r ] * scalars[ i ][ 0 ];
456 }
457 }
458 }
459 template< class Point, class DofVector >
460 void axpyH ( const Point &x, const HessianRangeType &factor, DofVector &dofs ) const
461 {
462 const std::size_t size = scalarBasisFunctionSet().size();
463 std::vector< typename HessianAll::Scalar > scalars( size );
464 HessianAll::apply( scalarBasisFunctionSet(), x, scalars );
465
466 const int xSize = DomainType::dimension;
467 for( int r = 0; r < dimRange; ++r )
468 {
469 for( std::size_t i = 0; i < size; ++i )
470 {
471 const GlobalDofType globalDof = dofAlignment_.globalDof( LocalDofType( r, i ) );
472 for ( int j = 0; j < xSize; ++j )
473 dofs[ globalDof ] += factor[ r ][ j ] * scalars[ i ][ 0 ][ j ];
474 }
475 }
476 }
477
478 template< class Evaluate, class Point, class DofVector >
479 void evaluateAll ( const Point &x, const DofVector &dofs, typename Evaluate::Vector &vector ) const
480 {
481 typename Evaluate::Scalar scalar;
482 for( int r = 0; r < dimRange; ++r )
483 {
484 SubDofVector< const DofVector, DofAlignmentType > subDofs( dofs, r, dofAlignment_ );
485 Evaluate::apply( scalarBasisFunctionSet(), x, subDofs, scalar );
486 vector[ r ] = scalar[ 0 ];
487 }
488 }
489
490 template< class Evaluate, class Point, class VectorArray >
491 void evaluateAll ( const Point &x, VectorArray &vectorials ) const
492 {
493 const std::size_t size = scalarBasisFunctionSet().size();
494 std::vector< typename Evaluate::Scalar > scalars( size );
495 Evaluate::apply( scalarBasisFunctionSet(), x, scalars );
496
497 typedef typename Evaluate::Vector Vector;
498
499 for( int r = 0; r < dimRange; ++r )
500 {
501 for( std::size_t i = 0; i < size; ++i )
502 {
503 const GlobalDofType globalDof = dofAlignment_.globalDof( LocalDofType( r, i ) );
504 Vector &vector = vectorials[ globalDof ];
505 vector = Vector( typename Vector::value_type( 0 ) );
506 vector[ r ] = scalars[ i ][ 0 ];
507 }
508 }
509 }
510
511 ScalarBasisFunctionSetType scalarBasisFunctionSet_;
512 DofAlignmentType dofAlignment_;
513 };
514
515
516
517 // VectorialBasisFunctionSet< ScalarBasisFunctionSet, Range >::EvaluateAll
518 // -----------------------------------------------------------------------
519
520 template< class ScalarBasisFunctionSet, class Range, template< class, class > class DofAlignment >
521 struct VectorialBasisFunctionSet< ScalarBasisFunctionSet, Range, DofAlignment >::EvaluateAll
522 {
523 typedef typename ScalarFunctionSpaceType::RangeType Scalar;
524 typedef RangeType Vector;
525
526 template< class Point, class SubDofVector >
527 static void apply ( const ScalarBasisFunctionSetType &scalarBasisFunctionSet,
528 const Point &x, const SubDofVector &dofs, Scalar &scalar )
529 {
530 scalarBasisFunctionSet.evaluateAll( x, dofs, scalar );
531 }
532
533 template< class Point, class ScalarArray >
534 static void apply ( const ScalarBasisFunctionSetType &scalarBasisFunctionSet,
535 const Point &x, ScalarArray &scalars )
536 {
537 scalarBasisFunctionSet.evaluateAll( x, scalars );
538 }
539 };
540
541
542
543 // VectorialBasisFunctionSet< ScalarBasisFunctionSet, Range >::JacobianAll
544 // -----------------------------------------------------------------------
545
546 template< class ScalarBasisFunctionSet, class Range, template< class, class > class DofAlignment >
547 struct VectorialBasisFunctionSet< ScalarBasisFunctionSet, Range, DofAlignment >::JacobianAll
548 {
549 typedef typename ScalarFunctionSpaceType::JacobianRangeType Scalar;
550 typedef JacobianRangeType Vector;
551
552 template< class Point, class SubDofVector >
553 static void apply ( const ScalarBasisFunctionSetType &scalarBasisFunctionSet,
554 const Point &x, const SubDofVector &dofs, Scalar &scalar )
555 {
556 scalarBasisFunctionSet.jacobianAll( x, dofs, scalar );
557 }
558
559 template< class Point, class ScalarArray >
560 static void apply ( const ScalarBasisFunctionSetType &scalarBasisFunctionSet,
561 const Point &x, ScalarArray &scalars )
562 {
563 scalarBasisFunctionSet.jacobianAll( x, scalars );
564 }
565 };
566
567
568
569 // VectorialBasisFunctionSet< ScalarBasisFunctionSet, Range >::HessianAll
570 // ----------------------------------------------------------------------
571
572 template< class ScalarBasisFunctionSet, class Range, template< class, class > class DofAlignment >
573 struct VectorialBasisFunctionSet< ScalarBasisFunctionSet, Range, DofAlignment >::HessianAll
574 {
575 typedef typename ScalarFunctionSpaceType::HessianRangeType Scalar;
576 typedef HessianRangeType Vector;
577
578 template< class Point, class SubDofVector >
579 static void apply ( const ScalarBasisFunctionSetType &scalarBasisFunctionSet,
580 const Point &x, const SubDofVector &dofs, Scalar &scalar )
581 {
582 scalarBasisFunctionSet.hessianAll( x, dofs, scalar );
583 }
584
585 template< class Point, class ScalarArray >
586 static void apply ( const ScalarBasisFunctionSetType &scalarBasisFunctionSet,
587 const Point &x, ScalarArray &scalars )
588 {
589 scalarBasisFunctionSet.hessianAll( x, scalars );
590 }
591 };
592
593 } // namespace Fem
594
595} // namespace Dune
596
597#endif // #ifndef DUNE_FEM_BASISFUNCTIONSET_VECTORIAL_HH
Interface documentation for Dof alignment classes used in VectorialBasisFunctionSet.
Definition: vectorial.hh:34
std::size_t GlobalDofType
global Dof type
Definition: vectorial.hh:37
LocalDofType localDof(const GlobalDofType &globalDof) const
map global to local Dof
Definition: vectorial.hh:68
std::pair< int, std::size_t > LocalDofType
local Dof type consists of coordinate number and Dof number in scalar basis function set
Definition: vectorial.hh:41
GlobalDofType globalDof(const LocalDofType &localDof) const
map local to global Dof
Definition: vectorial.hh:55
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
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
Implementation of DofAlignment.
Definition: vectorial.hh:96
GlobalDofType globalDof(const LocalDofType &localDof) const
map local to global Dof
Definition: vectorial.hh:111
LocalDofType localDof(const GlobalDofType &globalDof) const
map global to local Dof
Definition: vectorial.hh:117
int nop() const
obtain the number of integration points
Definition: quadrature.hh:291
actual interface class for quadratures
Definition: quadrature.hh:401
Extract Sub dof vector for single coordinate.
Definition: vectorial.hh:179
Builds a vectorial basis function set from given scalar basis function set.
Definition: vectorial.hh:279
Implementation of DofAlignment.
Definition: vectorial.hh:141
LocalDofType localDof(const GlobalDofType &globalDof) const
map global to local Dof
Definition: vectorial.hh:162
GlobalDofType globalDof(const LocalDofType &localDof) const
map local to global Dof
Definition: vectorial.hh:156
unspecified value type referenceElement(T &&... t)
Returns a reference element for the objects t....
typename Overloads::ScalarType< std::decay_t< V > >::type Scalar
Element type of some SIMD type.
Definition: interface.hh:235
Dune namespace.
Definition: alignedallocator.hh:13
convert functions space to space with new dim range
Definition: functionspace.hh:250
A unique label for each type of element that can occur in a grid.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Jul 27, 22:29, 2024)