DUNE PDELab (2.8)

monomialbasis.hh
1// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2// vi: set et ts=4 sw=2 sts=2:
3#ifndef DUNE_MONOMIALBASIS_HH
4#define DUNE_MONOMIALBASIS_HH
5
6#include <vector>
7
10
11#include <dune/geometry/type.hh>
12#include <dune/geometry/topologyfactory.hh>
13
14#include <dune/localfunctions/utility/field.hh>
15#include <dune/localfunctions/utility/multiindex.hh>
16#include <dune/localfunctions/utility/tensor.hh>
17
18namespace Dune
19{
20 /************************************************
21 * Classes for evaluating ''Monomials'' on any order
22 * for all reference element type.
23 * For a simplex topology these are the normal
24 * monomials for cube topologies the bimonomials.
25 * The construction follows the construction of the
26 * generic geometries using tensor products for
27 * prism generation and duffy transform for pyramid
28 * construction.
29 * A derivative argument can be applied, in which case
30 * all derivatives up to the desired order are
31 * evaluated. Note that for higher order derivatives
32 * only the ''lower'' part of the symmetric tensor
33 * is evaluated, e.g., passing derivative equal to 2
34 * to the class will provide the vector
35 * (d/dxdx p, d/dxydx p, d/dydy p,
36 * d/dx p, d/dy p, p)
37 * Important:
38 * So far the computation of the derivatives has not
39 * been fully implemented for general pyramid
40 * construction, i.e., in the case where a pyramid is
41 * build over a non simplex base geometry.
42 *
43 * Central classes:
44 * 1) template< GeometryType::Id geometryId, class F >
45 * class MonomialBasisImpl;
46 * Implementation of the monomial evaluation for
47 * a given topology and field type.
48 * The method evaluate fills a F* vector
49 * 2) template< GeometryType::Id geometryId, class F >
50 * class MonomialBasis
51 * The base class for the static monomial evaluation
52 * providing addiional evaluate methods including
53 * one taking std::vector<F>.
54 * 3) template< int dim, class F >
55 * class VirtualMonomialBasis
56 * Virtualization of the MonomialBasis.
57 * 4) template< int dim, class F >
58 * struct MonomialBasisFactory;
59 * A factory class for the VirtualMonomialBasis
60 * 5) template< int dim, class F >
61 * struct MonomialBasisProvider
62 * A singleton container for the virtual monomial
63 * basis
64 ************************************************/
65
66 // Internal Forward Declarations
67 // -----------------------------
68
69 template< GeometryType::Id geometryId >
70 class MonomialBasisSize;
71
72 template< GeometryType::Id geometryId, class F >
73 class MonomialBasis;
74
75
76
77 // MonomialBasisSize
78 // -----------------
79
80 template< GeometryType::Id geometryId >
81 class MonomialBasisSize
82 {
83 typedef MonomialBasisSize< geometryId > This;
84
85 public:
86 static This &instance ()
87 {
88 static This _instance;
89 return _instance;
90 }
91
92 unsigned int maxOrder_;
93
94 // sizes_[ k ]: number of basis functions of exactly order k
95 mutable unsigned int *sizes_;
96
97 // numBaseFunctions_[ k ] = sizes_[ 0 ] + ... + sizes_[ k ]
98 mutable unsigned int *numBaseFunctions_;
99
100 MonomialBasisSize ()
101 : maxOrder_( 0 ),
102 sizes_( 0 ),
103 numBaseFunctions_( 0 )
104 {
105 computeSizes( 2 );
106 }
107
108 ~MonomialBasisSize ()
109 {
110 delete[] sizes_;
111 delete[] numBaseFunctions_;
112 }
113
114 unsigned int operator() ( const unsigned int order ) const
115 {
116 return numBaseFunctions_[ order ];
117 }
118
119 unsigned int maxOrder() const
120 {
121 return maxOrder_;
122 }
123
124 void computeSizes ( unsigned int order )
125 {
126 if (order <= maxOrder_)
127 return;
128
129 maxOrder_ = order;
130
131 delete[] sizes_;
132 delete[] numBaseFunctions_;
133 sizes_ = new unsigned int[ order+1 ];
134 numBaseFunctions_ = new unsigned int[ order+1 ];
135
136 constexpr GeometryType geometry = geometryId;
137 constexpr auto dim = geometry.dim();
138
139 sizes_[ 0 ] = 1;
140 for( unsigned int k = 1; k <= order; ++k )
141 sizes_[ k ] = 0;
142
143 std::fill(numBaseFunctions_, numBaseFunctions_+order+1, 1);
144
145 for( int codim=dim-1; codim>=0; codim--)
146 {
147 if (Impl::isPrism(geometry.id(),dim,codim))
148 {
149 for( unsigned int k = 1; k <= order; ++k )
150 {
151 sizes_[ k ] = numBaseFunctions_[ k ] + k*sizes_[ k ];
152 numBaseFunctions_[ k ] = numBaseFunctions_[ k-1 ] + sizes_[ k ];
153 }
154 }
155 else
156 {
157 for( unsigned int k = 1; k <= order; ++k )
158 {
159 sizes_[ k ] = numBaseFunctions_[ k ];
160 numBaseFunctions_[ k ] = numBaseFunctions_[ k-1 ] + sizes_[ k ];
161 }
162 }
163 }
164 }
165 };
166
167
168
169 // MonomialBasisHelper
170 // -------------------
171
172
173 template< int mydim, int dim, class F >
174 struct MonomialBasisHelper
175 {
176 typedef MonomialBasisSize< GeometryTypes::simplex(mydim).toId() > MySize;
177 typedef MonomialBasisSize< GeometryTypes::simplex(dim).toId() > Size;
178
179 static void copy ( const unsigned int deriv, F *&wit, F *&rit,
180 const unsigned int numBaseFunctions, const F &z )
181 {
182 // n(d,k) = size<k>[d];
183 MySize &mySize = MySize::instance();
184 Size &size = Size::instance();
185
186 const F *const rend = rit + size( deriv )*numBaseFunctions;
187 for( ; rit != rend; )
188 {
189 F *prit = rit;
190
191 *wit = z * *rit;
192 ++rit, ++wit;
193
194 for( unsigned d = 1; d <= deriv; ++d )
195 {
196 #ifndef NDEBUG
197 const F *const derivEnd = rit + mySize.sizes_[ d ];
198 #endif
199
200 {
201 const F *const drend = rit + mySize.sizes_[ d ] - mySize.sizes_[ d-1 ];
202 for( ; rit != drend ; ++rit, ++wit )
203 *wit = z * *rit;
204 }
205
206 for (unsigned int j=1; j<d; ++j)
207 {
208 const F *const drend = rit + mySize.sizes_[ d-j ] - mySize.sizes_[ d-j-1 ];
209 for( ; rit != drend ; ++prit, ++rit, ++wit )
210 *wit = F(j) * *prit + z * *rit;
211 }
212 *wit = F(d) * *prit + z * *rit;
213 ++prit, ++rit, ++wit;
214 assert(derivEnd == rit);
215 rit += size.sizes_[d] - mySize.sizes_[d];
216 prit += size.sizes_[d-1] - mySize.sizes_[d-1];
217 const F *const emptyWitEnd = wit + size.sizes_[d] - mySize.sizes_[d];
218 for ( ; wit != emptyWitEnd; ++wit )
219 *wit = Zero<F>();
220 }
221 }
222 }
223 };
224
225
226
227 // MonomialBasisImpl
228 // -----------------
229
230 template< GeometryType::Id geometryId, class F>
231 class MonomialBasisImpl
232 {
233 public:
234 typedef F Field;
235
236 static constexpr GeometryType geometry = geometryId;
237
238 static const unsigned int dimDomain = geometry.dim();
239
240 typedef FieldVector< Field, dimDomain > DomainVector;
241
242 private:
243 friend class MonomialBasis< geometryId, Field >;
244
245 MonomialBasisImpl ()
246 {}
247
248 template< int dimD >
249 void evaluate ( const unsigned int deriv, const unsigned int order,
250 const FieldVector< Field, dimD > &x,
251 const unsigned int block, const unsigned int *const offsets,
252 Field *const values ) const
253 {
254 //start with vertex
255 *values = Unity< F >();
256 F *const end = values + block;
257 for( Field *it = values+1 ; it != end; ++it )
258 *it = Zero< F >();
259
261
262 if constexpr ( geometry == gt)
263 return;
264 else
265 evaluate<gt,dimD>(deriv, order, x, block, offsets, values );
266 }
267
268 template<GeometryType::Id baseGeometryId, int dimD >
269 void evaluate ( const unsigned int deriv, const unsigned int order,
270 const FieldVector< Field, dimD > &x,
271 const unsigned int block, const unsigned int *const offsets,
272 Field *const values ) const
273 {
274
275 static constexpr GeometryType baseGeometry = baseGeometryId;
276
277 auto constexpr isPrismatic = geometry.isPrismatic(baseGeometry.dim());
278
279 // compute
280 typedef MonomialBasisHelper< baseGeometry.dim() + 1, dimD, Field > Helper;
281 typedef MonomialBasisSize<baseGeometryId> BaseSize;
282
283 const BaseSize &size = BaseSize::instance();
284 const_cast<BaseSize&>(size).computeSizes(order);
285
286 const Field &z = x[ baseGeometry.dim() ];
287
288 Field *row0 = values;
289 for( unsigned int k = 1; k <= order; ++k )
290 {
291 Field *row1 = values + block*offsets[ k-1 ];
292 Field *wit = row1 + block*size.sizes_[ k ];
293 if constexpr ( isPrismatic )
294 Helper::copy( deriv, wit, row1, k*size.sizes_[ k ], z );
295 Helper::copy( deriv, wit, row0, size( k-1 ), z );
296 row0 = row1;
297 }
298
299 // stop if desired dimension is reached
300 if constexpr( baseGeometry.dim() == dimDomain-1)
301 return;
302 else
303 {
304 constexpr GeometryType nextGeometry = isPrismatic ? GeometryTypes::prismaticExtension(baseGeometry)
305 : GeometryTypes::conicalExtension(baseGeometry);
306
307 evaluate<nextGeometry.toId(),dimD>(deriv, order, x, block, offsets, values );
308 }
309 }
310
311 void integrate ( const unsigned int order,
312 const unsigned int *const offsets,
313 Field *const values ) const
314 {
315 //start with vertex
316 values[ 0 ] = Unity< Field >();
317 static constexpr GeometryType gt = GeometryTypes::vertex;
318
319 if constexpr ( geometry == gt)
320 return;
321 else
322 integrate<gt>(order, offsets, values);
323 }
324
325 template<GeometryType::Id baseGeometryId>
326 void integrate ( const unsigned int order,
327 const unsigned int *const offsets,
328 Field *const values) const
329 {
330 static constexpr GeometryType baseGeometry = baseGeometryId;
331
332 auto constexpr isPrismatic = geometry.isPrismatic(baseGeometry.dim());
333
334 // decide which kind of integration should be performed
335 if constexpr ( isPrismatic )
336 integratePrismatic<baseGeometry>(order, offsets, values);
337 else
338 integrateConical<baseGeometry>(order, offsets, values);
339
340 // stop if the desired dimension is reached
341 if constexpr( baseGeometry.dim() == dimDomain-1)
342 return;
343 else
344 {
345 static constexpr GeometryType nextGeometry = (isPrismatic ? GeometryTypes::prismaticExtension(baseGeometry)
346 : GeometryTypes::conicalExtension(baseGeometry));
347
348 integrate<nextGeometry.toId()>(order, offsets, values);
349 }
350
351 }
352
353 template<GeometryType::Id baseGeometryId>
354 void integratePrismatic ( const unsigned int order,
355 const unsigned int *const offsets,
356 Field *const values ) const
357 {
358 typedef MonomialBasisSize<baseGeometryId> BaseSize;
359 static const BaseSize &size = BaseSize::instance();
360 const unsigned int *const baseSizes = size.sizes_;
361
362 static constexpr GeometryType baseGeometry = baseGeometryId;
363 static constexpr GeometryType nextGeometry = GeometryTypes::prismaticExtension(baseGeometry);
364
365 typedef MonomialBasisSize<nextGeometry.toId()> Size;
366 static const Size &mySize = Size::instance();
367
368 Field *row0 = values;
369 for( unsigned int k = 1; k <= order; ++k )
370 {
371 Field *const row1begin = values + offsets[ k-1 ];
372 Field *const row1End = row1begin + mySize.sizes_[ k ];
373 assert( (unsigned int)(row1End - values) <= offsets[ k ] );
374
375 Field *row1 = row1begin;
376 Field *it = row1begin + baseSizes[ k ];
377 for( unsigned int j = 1; j <= k; ++j )
378 {
379 Field *const end = it + baseSizes[ k ];
380 assert( (unsigned int)(end - values) <= offsets[ k ] );
381 for( ; it != end; ++row1, ++it )
382 *it = (Field( j ) / Field( j+1 )) * (*row1);
383 }
384 for( ; it != row1End; ++row0, ++it )
385 *it = (Field( k ) / Field( k+1 )) * (*row0);
386 row0 = row1;
387 }
388 }
389
390
391 template<GeometryType::Id baseGeometryId>
392 void integrateConical ( const unsigned int order,
393 const unsigned int *const offsets,
394 Field *const values) const
395 {
396 typedef MonomialBasisSize<baseGeometryId> BaseSize;
397 static const BaseSize &size = BaseSize::instance();
398 const unsigned int *const baseSizes = size.sizes_;
399
400 static constexpr GeometryType baseGeometry = baseGeometryId;
401
402 {
403 Field *const col0End = values + baseSizes[ 0 ];
404 for( Field *it = values; it != col0End; ++it )
405 *it *= Field( 1 ) / Field( int(baseGeometry.dim()+1) );
406 }
407
408 Field *row0 = values;
409 for( unsigned int k = 1; k <= order; ++k )
410 {
411 const Field factor = (Field( 1 ) / Field( k + baseGeometry.dim()+1));
412
413 Field *const row1 = values+offsets[ k-1 ];
414 Field *const col0End = row1 + baseSizes[ k ];
415 Field *it = row1;
416 for( ; it != col0End; ++it )
417 *it *= factor;
418 for( unsigned int i = 1; i <= k; ++i )
419 {
420 Field *const end = it + baseSizes[ k-i ];
421 assert( (unsigned int)(end - values) <= offsets[ k ] );
422 for( ; it != end; ++row0, ++it )
423 *it = (*row0) * (Field( i ) * factor);
424 }
425 row0 = row1;
426 }
427 }
428
429 };
430
431
432 // MonomialBasis
433 // -------------
434
435 template< GeometryType::Id geometryId, class F >
436 class MonomialBasis
437 : public MonomialBasisImpl< geometryId, F >
438 {
439 static constexpr GeometryType geometry = geometryId;
440 typedef MonomialBasis< geometryId, F > This;
441 typedef MonomialBasisImpl< geometryId, F > Base;
442
443 public:
444 static const unsigned int dimension = Base::dimDomain;
445 static const unsigned int dimRange = 1;
446
447 typedef typename Base::Field Field;
448
449 typedef typename Base::DomainVector DomainVector;
450
451 typedef Dune::FieldVector<Field,dimRange> RangeVector;
452
453 typedef MonomialBasisSize<geometryId> Size;
454
455 MonomialBasis (unsigned int order)
456 : Base(),
457 order_(order),
458 size_(Size::instance())
459 {
460 assert(order<=1024); // avoid wrapping of unsigned int (0-1) order=1024 is quite hight...)
461 }
462
463 const unsigned int *sizes ( unsigned int order ) const
464 {
465 size_.computeSizes( order );
466 return size_.numBaseFunctions_;
467 }
468
469 const unsigned int *sizes () const
470 {
471 return sizes( order_ );
472 }
473
474 unsigned int size () const
475 {
476 size_.computeSizes( order_ );
477 return size_( order_ );
478 }
479
480 unsigned int derivSize ( const unsigned int deriv ) const
481 {
482 MonomialBasisSize< GeometryTypes::simplex(dimension).toId() >::instance().computeSizes( deriv );
483 return MonomialBasisSize< GeometryTypes::simplex(dimension).toId() >::instance() ( deriv );
484 }
485
486 unsigned int order () const
487 {
488 return order_ ;
489 }
490
491 unsigned int topologyId ( ) const
492 {
493 return geometry.id();
494 }
495
496 void evaluate ( const unsigned int deriv, const DomainVector &x,
497 Field *const values ) const
498 {
499 Base::evaluate( deriv, order_, x, derivSize( deriv ), sizes( order_ ), values );
500 }
501
502 template <unsigned int deriv>
503 void evaluate ( const DomainVector &x,
504 Field *const values ) const
505 {
506 evaluate( deriv, x, values );
507 }
508
509 template<unsigned int deriv, class Vector >
510 void evaluate ( const DomainVector &x,
511 Vector &values ) const
512 {
513 evaluate<deriv>(x,&(values[0]));
514 }
515 template<unsigned int deriv, DerivativeLayoutNS::DerivativeLayout layout >
516 void evaluate ( const DomainVector &x,
517 Derivatives<Field,dimension,1,deriv,layout> *values ) const
518 {
519 evaluate<deriv>(x,&(values->block()));
520 }
521 template< unsigned int deriv >
522 void evaluate ( const DomainVector &x,
523 FieldVector<Field,Derivatives<Field,dimension,1,deriv,DerivativeLayoutNS::value>::size> *values ) const
524 {
525 evaluate(0,x,&(values[0][0]));
526 }
527
528 template<class Vector >
529 void evaluate ( const DomainVector &x,
530 Vector &values ) const
531 {
532 evaluate<0>(x,&(values[0]));
533 }
534
535 template< class DVector, class RVector >
536 void evaluate ( const DVector &x, RVector &values ) const
537 {
538 assert( DVector::dimension == dimension);
539 DomainVector bx;
540 for( int d = 0; d < dimension; ++d )
541 field_cast( x[ d ], bx[ d ] );
542 evaluate<0>( bx, values );
543 }
544
545 void integrate ( Field *const values ) const
546 {
547 Base::integrate( order_, sizes( order_ ), values );
548 }
549 template <class Vector>
550 void integrate ( Vector &values ) const
551 {
552 integrate( &(values[ 0 ]) );
553 }
554 private:
555 MonomialBasis(const This&);
556 This& operator=(const This&);
557 unsigned int order_;
558 Size &size_;
559 };
560
561
562
563 // StdMonomialBasis
564 // ----------------
565
566 template< int dim,class F >
567 class StandardMonomialBasis
568 : public MonomialBasis< GeometryTypes::simplex(dim).toId() , F >
569 {
570 typedef StandardMonomialBasis< dim, F > This;
571 typedef MonomialBasis< GeometryTypes::simplex(dim).toId(), F > Base;
572
573 public:
574 static constexpr GeometryType geometry = GeometryTypes::simplex(dim);
575 static const int dimension = dim;
576
577 StandardMonomialBasis ( unsigned int order )
578 : Base( order )
579 {}
580 };
581
582
583
584 // StandardBiMonomialBasis
585 // -----------------------
586
587 template< int dim, class F >
588 class StandardBiMonomialBasis
589 : public MonomialBasis< GeometryTypes::cube(dim).toId() , F >
590 {
591 typedef StandardBiMonomialBasis< dim, F > This;
592 typedef MonomialBasis< GeometryTypes::cube(dim).toId() , F > Base;
593
594 public:
595 static constexpr GeometryType geometry = GeometryTypes::cube(dim);
596 static const int dimension = dim;
597
598 StandardBiMonomialBasis ( unsigned int order )
599 : Base( order )
600 {}
601 };
602
603 // -----------------------------------------------------------
604 // -----------------------------------------------------------
605 // VirtualMonomialBasis
606 // -------------------
607
608 template< int dim, class F >
609 class VirtualMonomialBasis
610 {
611 typedef VirtualMonomialBasis< dim, F > This;
612
613 public:
614 typedef F Field;
615 typedef F StorageField;
616 static const int dimension = dim;
617 static const unsigned int dimRange = 1;
618
619 typedef FieldVector<Field,dimension> DomainVector;
620 typedef FieldVector<Field,dimRange> RangeVector;
621
622 [[deprecated("Use VirtualMonomialBasis(GeometryType gt, unsigned int order) instead.")]]
623 explicit VirtualMonomialBasis(unsigned int topologyId,
624 unsigned int order)
625 : order_(order), geometry_(GeometryType(topologyId,dim)) {}
626
627 explicit VirtualMonomialBasis(const GeometryType& gt,
628 unsigned int order)
629 : order_(order), geometry_(gt) {}
630
631 virtual ~VirtualMonomialBasis() {}
632
633 virtual const unsigned int *sizes ( ) const = 0;
634
635 unsigned int size ( ) const
636 {
637 return sizes( )[ order_ ];
638 }
639
640 unsigned int order () const
641 {
642 return order_;
643 }
644
645 [[deprecated("Use type().id() instead.")]]
646 unsigned int topologyId ( ) const
647 {
648 return type().id();
649 }
650
651 GeometryType type() const
652 {
653 return geometry_;
654 }
655
656 virtual void evaluate ( const unsigned int deriv, const DomainVector &x,
657 Field *const values ) const = 0;
658 template < unsigned int deriv >
659 void evaluate ( const DomainVector &x,
660 Field *const values ) const
661 {
662 evaluate( deriv, x, values );
663 }
664 template < unsigned int deriv, int size >
665 void evaluate ( const DomainVector &x,
666 Dune::FieldVector<Field,size> *const values ) const
667 {
668 evaluate( deriv, x, &(values[0][0]) );
669 }
670 template<unsigned int deriv, DerivativeLayoutNS::DerivativeLayout layout >
671 void evaluate ( const DomainVector &x,
672 Derivatives<Field,dimension,1,deriv,layout> *values ) const
673 {
674 evaluate<deriv>(x,&(values->block()));
675 }
676 template <unsigned int deriv, class Vector>
677 void evaluate ( const DomainVector &x,
678 Vector &values ) const
679 {
680 evaluate<deriv>( x, &(values[ 0 ]) );
681 }
682 template< class Vector >
683 void evaluate ( const DomainVector &x,
684 Vector &values ) const
685 {
686 evaluate<0>(x,values);
687 }
688 template< class DVector, class RVector >
689 void evaluate ( const DVector &x, RVector &values ) const
690 {
691 assert( DVector::dimension == dimension);
692 DomainVector bx;
693 for( int d = 0; d < dimension; ++d )
694 field_cast( x[ d ], bx[ d ] );
695 evaluate<0>( bx, values );
696 }
697 template< unsigned int deriv, class DVector, class RVector >
698 void evaluate ( const DVector &x, RVector &values ) const
699 {
700 assert( DVector::dimension == dimension);
701 DomainVector bx;
702 for( int d = 0; d < dimension; ++d )
703 field_cast( x[ d ], bx[ d ] );
704 evaluate<deriv>( bx, values );
705 }
706
707 virtual void integrate ( Field *const values ) const = 0;
708 template <class Vector>
709 void integrate ( Vector &values ) const
710 {
711 integrate( &(values[ 0 ]) );
712 }
713 protected:
714 unsigned int order_;
715 GeometryType geometry_;
716 };
717
718 template< GeometryType::Id geometryId, class F >
719 class VirtualMonomialBasisImpl
720 : public VirtualMonomialBasis< static_cast<GeometryType>(geometryId).dim(), F >
721 {
722 static constexpr GeometryType geometry = geometryId;
723 typedef VirtualMonomialBasis< geometry.dim(), F > Base;
724 typedef VirtualMonomialBasisImpl< geometryId, F > This;
725
726 public:
727 typedef typename Base::Field Field;
728 typedef typename Base::DomainVector DomainVector;
729
730 VirtualMonomialBasisImpl(unsigned int order)
731 : Base(geometry,order), basis_(order)
732 {}
733
734 const unsigned int *sizes ( ) const
735 {
736 return basis_.sizes(order_);
737 }
738
739 void evaluate ( const unsigned int deriv, const DomainVector &x,
740 Field *const values ) const
741 {
742 basis_.evaluate(deriv,x,values);
743 }
744
745 void integrate ( Field *const values ) const
746 {
747 basis_.integrate(values);
748 }
749
750 private:
751 MonomialBasis<geometryId,Field> basis_;
752 using Base::order_;
753 };
754
755 // MonomialBasisFactory
756 // --------------------
757
758 template< int dim, class F >
759 struct MonomialBasisFactory
760 {
761 static const unsigned int dimension = dim;
762 typedef F StorageField;
763
764 typedef unsigned int Key;
765 typedef const VirtualMonomialBasis< dimension, F > Object;
766
767 template < int dd, class FF >
768 struct EvaluationBasisFactory
769 {
770 typedef MonomialBasisFactory<dd,FF> Type;
771 };
772
773 template< GeometryType::Id geometryId >
774 static Object* create ( const Key &order )
775 {
776 return new VirtualMonomialBasisImpl< geometryId, StorageField >( order );
777 }
778 static void release( Object *object ) { delete object; }
779 };
780
781
782
783 // MonomialBasisProvider
784 // ---------------------
785
786 template< int dim, class SF >
787 struct MonomialBasisProvider
788 : public TopologySingletonFactory< MonomialBasisFactory< dim, SF > >
789 {
790 static const unsigned int dimension = dim;
791 typedef SF StorageField;
792 template < int dd, class FF >
793 struct EvaluationBasisFactory
794 {
795 typedef MonomialBasisProvider<dd,FF> Type;
796 };
797 };
798
799}
800
801#endif
vector space out of a tensor product of fields.
Definition: fvector.hh:95
constexpr Id toId() const
Create an Id representation of this GeometryType.
Definition: type.hh:219
GeometryType
Type representing VTK's entity geometry types.
Definition: common.hh:130
Implements a matrix constructed from a given type representing a field and compile-time given number ...
Implements a vector constructed from a given type representing a field and a compile-time given size.
bool gt(const T &first, const T &second, typename EpsilonType< T >::Type epsilon)
test if first greater than second
Definition: float_cmp.cc:156
constexpr GeometryType cube(unsigned int dim)
Returns a GeometryType representing a hypercube of dimension dim.
Definition: type.hh:470
constexpr GeometryType simplex(unsigned int dim)
Returns a GeometryType representing a simplex of dimension dim.
Definition: type.hh:461
constexpr GeometryType vertex
GeometryType representing a vertex.
Definition: type.hh:504
constexpr GeometryType prismaticExtension(const GeometryType &gt)
Return GeometryType of a prismatic construction with gt as base
Definition: type.hh:491
constexpr GeometryType conicalExtension(const GeometryType &gt)
Return GeometryType of a conical construction with gt as base
Definition: type.hh:485
Dune namespace.
Definition: alignedallocator.hh:11
void field_cast(const F1 &f1, F2 &f2)
a helper class to cast from one field to another
Definition: field.hh:157
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)