DUNE PDELab (git)

referenceelementimplementation.hh
1// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2// vi: set et ts=4 sw=2 sts=2:
3// SPDX-FileCopyrightInfo: Copyright © DUNE Project contributors, see file LICENSE.md in module root
4// SPDX-License-Identifier: LicenseRef-GPL-2.0-only-with-DUNE-exception
5#ifndef DUNE_GEOMETRY_REFERENCEELEMENTIMPLEMENTATION_HH
6#define DUNE_GEOMETRY_REFERENCEELEMENTIMPLEMENTATION_HH
7
8#include <cassert>
9
10#include <algorithm>
11#include <limits>
12#include <tuple>
13#include <utility>
14#include <vector>
15#include <array>
16#include <bitset>
17
20#include <dune/common/hybridutilities.hh>
22#include <dune/common/iteratorrange.hh>
23#include <dune/common/math.hh>
24
25#include <dune/geometry/referenceelement.hh>
27#include <dune/geometry/type.hh>
28
29namespace Dune
30{
31
32 namespace Geo
33 {
34
35#ifndef DOXYGEN
36
37 // Internal Forward Declarations
38 // -----------------------------
39
40 namespace Impl
41 {
42 template< class ctype, int dim >
43 class ReferenceElementContainer;
44 }
45
46 template< class ctype, int dim >
47 struct ReferenceElements;
48
49
50
51 namespace Impl
52 {
53
54 using Dune::Impl::isPrism;
55 using Dune::Impl::isPyramid;
56 using Dune::Impl::baseTopologyId;
57 using Dune::Impl::prismConstruction;
58 using Dune::Impl::pyramidConstruction;
59 using Dune::Impl::numTopologies;
60
62 unsigned int size ( unsigned int topologyId, int dim, int codim );
63
64
65
73 unsigned int subTopologyId ( unsigned int topologyId, int dim, int codim, unsigned int i );
74
75
76
77 // subTopologyNumbering
78 // --------------------
79
80 void subTopologyNumbering ( unsigned int topologyId, int dim, int codim, unsigned int i, int subcodim,
81 unsigned int *beginOut, unsigned int *endOut );
82
83
84
85
86 // checkInside
87 // -----------
88
89 template< class ct, int cdim >
90 inline bool
91 checkInside ( unsigned int topologyId, int dim, const FieldVector< ct, cdim > &x, ct tolerance, ct factor = ct( 1 ) )
92 {
93 assert( (dim >= 0) && (dim <= cdim) );
94 assert( topologyId < numTopologies( dim ) );
95
96 if( dim > 0 )
97 {
98 const ct baseFactor = (isPrism( topologyId, dim ) ? factor : factor - x[ dim-1 ]);
99 if( (x[ dim-1 ] > -tolerance) && (factor - x[ dim-1 ] > -tolerance) )
100 return checkInside< ct, cdim >( baseTopologyId( topologyId, dim ), dim-1, x, tolerance, baseFactor );
101 else
102 return false;
103 }
104 else
105 return true;
106 }
107
108
109
110 // referenceCorners
111 // ----------------
112
113 template< class ct, int cdim >
114 inline unsigned int
115 referenceCorners ( unsigned int topologyId, int dim, FieldVector< ct, cdim > *corners )
116 {
117 assert( (dim >= 0) && (dim <= cdim) );
118 assert( topologyId < numTopologies( dim ) );
119
120 if( dim > 0 )
121 {
122 const unsigned int nBaseCorners
123 = referenceCorners( baseTopologyId( topologyId, dim ), dim-1, corners );
124 assert( nBaseCorners == size( baseTopologyId( topologyId, dim ), dim-1, dim-1 ) );
125 if( isPrism( topologyId, dim ) )
126 {
127 std::copy( corners, corners + nBaseCorners, corners + nBaseCorners );
128 for( unsigned int i = 0; i < nBaseCorners; ++i )
129 corners[ i+nBaseCorners ][ dim-1 ] = ct( 1 );
130 return 2*nBaseCorners;
131 }
132 else
133 {
134 corners[ nBaseCorners ] = FieldVector< ct, cdim >( ct( 0 ) );
135 corners[ nBaseCorners ][ dim-1 ] = ct( 1 );
136 return nBaseCorners+1;
137 }
138 }
139 else
140 {
141 *corners = FieldVector< ct, cdim >( ct( 0 ) );
142 return 1;
143 }
144 }
145
146
147
148 // referenceVolume
149 // ---------------
150
151 unsigned long referenceVolumeInverse ( unsigned int topologyId, int dim );
152
153 template< class ct >
154 inline ct referenceVolume ( unsigned int topologyId, int dim )
155 {
156 return ct( 1 ) / ct( referenceVolumeInverse( topologyId, dim ) );
157 }
158
159
160
161 // referenceOrigins
162 // ----------------
163
164 template< class ct, int cdim >
165 inline unsigned int
166 referenceOrigins ( unsigned int topologyId, int dim, int codim, FieldVector< ct, cdim > *origins )
167 {
168 assert( (dim >= 0) && (dim <= cdim) );
169 assert( topologyId < numTopologies( dim ) );
170 assert( (codim >= 0) && (codim <= dim) );
171
172 if( codim > 0 )
173 {
174 const unsigned int baseId = baseTopologyId( topologyId, dim );
175 if( isPrism( topologyId, dim ) )
176 {
177 const unsigned int n = (codim < dim ? referenceOrigins( baseId, dim-1, codim, origins ) : 0);
178 const unsigned int m = referenceOrigins( baseId, dim-1, codim-1, origins+n );
179 for( unsigned int i = 0; i < m; ++i )
180 {
181 origins[ n+m+i ] = origins[ n+i ];
182 origins[ n+m+i ][ dim-1 ] = ct( 1 );
183 }
184 return n+2*m;
185 }
186 else
187 {
188 const unsigned int m = referenceOrigins( baseId, dim-1, codim-1, origins );
189 if( codim == dim )
190 {
191 origins[ m ] = FieldVector< ct, cdim >( ct( 0 ) );
192 origins[ m ][ dim-1 ] = ct( 1 );
193 return m+1;
194 }
195 else
196 return m+referenceOrigins( baseId, dim-1, codim, origins+m );
197 }
198 }
199 else
200 {
201 origins[ 0 ] = FieldVector< ct, cdim >( ct( 0 ) );
202 return 1;
203 }
204 }
205
206
207
208 // referenceEmbeddings
209 // -------------------
210
211 template< class ct, int cdim, int mydim >
212 inline unsigned int
213 referenceEmbeddings ( unsigned int topologyId, int dim, int codim,
214 FieldVector< ct, cdim > *origins,
215 FieldMatrix< ct, mydim, cdim > *jacobianTransposeds )
216 {
217 assert( (0 <= codim) && (codim <= dim) && (dim <= cdim) );
218 assert( (dim - codim <= mydim) && (mydim <= cdim) );
219 assert( topologyId < numTopologies( dim ) );
220
221 if( (0 < codim) && (codim <= dim) )
222 {
223 const unsigned int baseId = baseTopologyId( topologyId, dim );
224 if( isPrism( topologyId, dim ) )
225 {
226 const unsigned int n = (codim < dim ? referenceEmbeddings( baseId, dim-1, codim, origins, jacobianTransposeds ) : 0);
227 for( unsigned int i = 0; i < n; ++i )
228 jacobianTransposeds[ i ][ dim-codim-1 ][ dim-1 ] = ct( 1 );
229
230 const unsigned int m = referenceEmbeddings( baseId, dim-1, codim-1, origins+n, jacobianTransposeds+n );
231 std::copy( origins+n, origins+n+m, origins+n+m );
232 std::copy( jacobianTransposeds+n, jacobianTransposeds+n+m, jacobianTransposeds+n+m );
233 for( unsigned int i = 0; i < m; ++i )
234 origins[ n+m+i ][ dim-1 ] = ct( 1 );
235
236 return n+2*m;
237 }
238 else // !isPrism
239 {
240 const unsigned int m = referenceEmbeddings( baseId, dim-1, codim-1, origins, jacobianTransposeds );
241 if( codim == dim )
242 {
243 origins[ m ] = FieldVector< ct, cdim >( ct( 0 ) );
244 origins[ m ][ dim-1 ] = ct( 1 );
245 jacobianTransposeds[ m ] = FieldMatrix< ct, mydim, cdim >( ct( 0 ) );
246 return m+1;
247 }
248 else if( codim < dim )
249 {
250 const unsigned int n = referenceEmbeddings( baseId, dim-1, codim, origins+m, jacobianTransposeds+m );
251 for( unsigned int i = 0; i < n; ++i )
252 {
253 for( int k = 0; k < dim-1; ++k )
254 jacobianTransposeds[ m+i ][ dim-codim-1 ][ k ] = -origins[ m+i ][ k ];
255 jacobianTransposeds[ m+i ][ dim-codim-1 ][ dim-1 ] = ct( 1 );
256 }
257 return m+n;
258 }
259 }
260 }
261 else if( codim == 0 )
262 {
263 origins[ 0 ] = FieldVector< ct, cdim >( ct( 0 ) );
264 jacobianTransposeds[ 0 ] = FieldMatrix< ct, mydim, cdim >( ct( 0 ) );
265 for( int k = 0; k < dim; ++k )
266 jacobianTransposeds[ 0 ][ k ][ k ] = ct( 1 );
267 return 1;
268 }
269
270 // this point should not be reached since all cases are handled before.
271 std::abort();
272 return 0;
273 }
274
275
276
277 // referenceIntegrationOuterNormals
278 // --------------------------------
279
280 template< class ct, int cdim >
281 inline unsigned int
282 referenceIntegrationOuterNormals ( unsigned int topologyId, int dim,
283 const FieldVector< ct, cdim > *origins,
284 FieldVector< ct, cdim > *normals )
285 {
286 assert( (dim > 0) && (dim <= cdim) );
287 assert( topologyId < numTopologies( dim ) );
288
289 if( dim > 1 )
290 {
291 const unsigned int baseId = baseTopologyId( topologyId, dim );
292 if( isPrism( topologyId, dim ) )
293 {
294 const unsigned int numBaseFaces
295 = referenceIntegrationOuterNormals( baseId, dim-1, origins, normals );
296
297 for( unsigned int i = 0; i < 2; ++i )
298 {
299 normals[ numBaseFaces+i ] = FieldVector< ct, cdim >( ct( 0 ) );
300 normals[ numBaseFaces+i ][ dim-1 ] = ct( 2*int( i )-1 );
301 }
302
303 return numBaseFaces+2;
304 }
305 else
306 {
307 normals[ 0 ] = FieldVector< ct, cdim >( ct( 0 ) );
308 normals[ 0 ][ dim-1 ] = ct( -1 );
309
310 const unsigned int numBaseFaces
311 = referenceIntegrationOuterNormals( baseId, dim-1, origins+1, normals+1 );
312 for( unsigned int i = 1; i <= numBaseFaces; ++i )
313 normals[ i ][ dim-1 ] = normals[ i ]*origins[ i ];
314
315 return numBaseFaces+1;
316 }
317 }
318 else
319 {
320 for( unsigned int i = 0; i < 2; ++i )
321 {
322 normals[ i ] = FieldVector< ct, cdim >( ct( 0 ) );
323 normals[ i ][ 0 ] = ct( 2*int( i )-1 );
324 }
325
326 return 2;
327 }
328 }
329
330 template< class ct, int cdim >
331 inline unsigned int
332 referenceIntegrationOuterNormals ( unsigned int topologyId, int dim,
333 FieldVector< ct, cdim > *normals )
334 {
335 assert( (dim > 0) && (dim <= cdim) );
336
337 FieldVector< ct, cdim > *origins
338 = new FieldVector< ct, cdim >[ size( topologyId, dim, 1 ) ];
339 referenceOrigins( topologyId, dim, 1, origins );
340
341 const unsigned int numFaces
342 = referenceIntegrationOuterNormals( topologyId, dim, origins, normals );
343 assert( numFaces == size( topologyId, dim, 1 ) );
344
345 delete[] origins;
346
347 return numFaces;
348 }
349
350 } // namespace Impl
351
352
353
354 // ReferenceElement
355 // ----------------
356
375 template< class ctype_, int dim >
376 class ReferenceElementImplementation
377 {
378
379 public:
380
382 using ctype = ctype_;
383
385 using CoordinateField = ctype;
386
388 using Coordinate = Dune::FieldVector<ctype,dim>;
389
391 static constexpr int dimension = dim;
392
394 typedef ctype Volume;
395
396 private:
397
398 friend class Impl::ReferenceElementContainer< ctype, dim >;
399
400 struct SubEntityInfo;
401
402 template< int codim > struct CreateGeometries;
403
404 public:
406 template< int codim >
407 struct Codim
408 {
410 typedef AffineGeometry< ctype, dim-codim, dim > Geometry;
411 };
412
413 // ReferenceElement cannot be copied.
414 ReferenceElementImplementation ( const ReferenceElementImplementation& ) = delete;
415
416 // ReferenceElementImplementation cannot be copied.
417 ReferenceElementImplementation& operator= ( const ReferenceElementImplementation& ) = delete;
418
419 // ReferenceElementImplementation is default-constructible (required for storage in std::array)
420 ReferenceElementImplementation () = default;
421
426 int size ( int c ) const
427 {
428 assert( (c >= 0) && (c <= dim) );
429 return info_[ c ].size();
430 }
431
443 int size ( int i, int c, int cc ) const
444 {
445 assert( (i >= 0) && (i < size( c )) );
446 return info_[ c ][ i ].size( cc );
447 }
448
462 int subEntity ( int i, int c, int ii, int cc ) const
463 {
464 assert( (i >= 0) && (i < size( c )) );
465 return info_[ c ][ i ].number( ii, cc );
466 }
467
483 auto subEntities ( int i, int c, int cc ) const
484 {
485 assert( (i >= 0) && (i < size( c )) );
486 return info_[ c ][ i ].numbers( cc );
487 }
488
497 const GeometryType &type ( int i, int c ) const
498 {
499 assert( (i >= 0) && (i < size( c )) );
500 return info_[ c ][ i ].type();
501 }
502
504 const GeometryType &type () const { return type( 0, 0 ); }
505
515 const Coordinate &position( int i, int c ) const
516 {
517 assert( (c >= 0) && (c <= dim) );
518 return baryCenters_[ c ][ i ];
519 }
520
528 bool checkInside ( const Coordinate &local ) const
529 {
530 const ctype tolerance = ctype( 64 ) * std::numeric_limits< ctype >::epsilon();
531 return Impl::template checkInside< ctype, dim >( type().id(), dim, local, tolerance );
532 }
533
545 template< int codim >
546 typename Codim< codim >::Geometry geometry ( int i ) const
547 {
548 return std::get< codim >( geometries_ )[ i ];
549 }
550
552 Volume volume () const
553 {
554 return volume_;
555 }
556
564 const Coordinate &integrationOuterNormal ( int face ) const
565 {
566 assert( (face >= 0) && (face < int( integrationNormals_.size() )) );
567 return integrationNormals_[ face ];
568 }
569
570 private:
571 void initialize ( unsigned int topologyId )
572 {
573 assert( topologyId < Impl::numTopologies( dim ) );
574
575 // set up subentities
576 for( int codim = 0; codim <= dim; ++codim )
577 {
578 const unsigned int size = Impl::size( topologyId, dim, codim );
579 info_[ codim ].resize( size );
580 for( unsigned int i = 0; i < size; ++i )
581 info_[ codim ][ i ].initialize( topologyId, codim, i );
582 }
583
584 // compute corners
585 const unsigned int numVertices = size( dim );
586 baryCenters_[ dim ].resize( numVertices );
587 Impl::referenceCorners( topologyId, dim, &(baryCenters_[ dim ][ 0 ]) );
588
589 // compute barycenters
590 for( int codim = 0; codim < dim; ++codim )
591 {
592 baryCenters_[ codim ].resize( size(codim) );
593 for( int i = 0; i < size( codim ); ++i )
594 {
595 baryCenters_[ codim ][ i ] = Coordinate( ctype( 0 ) );
596 const unsigned int numCorners = size( i, codim, dim );
597 for( unsigned int j = 0; j < numCorners; ++j )
598 baryCenters_[ codim ][ i ] += baryCenters_[ dim ][ subEntity( i, codim, j, dim ) ];
599 baryCenters_[ codim ][ i ] *= ctype( 1 ) / ctype( numCorners );
600 }
601 }
602
603 // compute reference element volume
604 volume_ = Impl::template referenceVolume< ctype >( topologyId, dim );
605
606 // compute integration outer normals
607 if( dim > 0 )
608 {
609 integrationNormals_.resize( size( 1 ) );
610 Impl::referenceIntegrationOuterNormals( topologyId, dim, &(integrationNormals_[ 0 ]) );
611 }
612
613 // set up geometries
614 Hybrid::forEach( std::make_index_sequence< dim+1 >{}, [ & ]( auto i ){ CreateGeometries< i >::apply( *this, geometries_ ); } );
615 }
616
617 template< int... codim >
618 static std::tuple< std::vector< typename Codim< codim >::Geometry >... >
619 makeGeometryTable ( std::integer_sequence< int, codim... > );
620
622 typedef decltype( makeGeometryTable( std::make_integer_sequence< int, dim+1 >() ) ) GeometryTable;
623
625 ctype volume_;
626
627 std::vector< Coordinate > baryCenters_[ dim+1 ];
628 std::vector< Coordinate > integrationNormals_;
629
631 GeometryTable geometries_;
632
633 std::vector< SubEntityInfo > info_[ dim+1 ];
634 };
635
637 template< class ctype, int dim >
638 struct ReferenceElementImplementation< ctype, dim >::SubEntityInfo
639 {
640 // Compute upper bound for the number of subsentities.
641 // If someone knows an explicit formal feel free to
642 // implement it here.
643 static constexpr std::size_t maxSubEntityCount()
644 {
645 std::size_t maxCount=0;
646 for(std::size_t codim=0; codim<=dim; ++codim)
647 maxCount = std::max(maxCount, binomial(std::size_t(dim),codim)*(1 << codim));
648 return maxCount;
649 }
650
651 using SubEntityFlags = std::bitset<maxSubEntityCount()>;
652
653 class SubEntityRange
654 : public Dune::IteratorRange<const unsigned int*>
655 {
656 using Base = typename Dune::IteratorRange<const unsigned int*>;
657
658 public:
659
660 using iterator = Base::iterator;
661 using const_iterator = Base::const_iterator;
662
663 SubEntityRange(const iterator& begin, const iterator& end, const SubEntityFlags& contains) :
664 Base(begin, end),
665 containsPtr_(&contains),
666 size_(end-begin)
667 {}
668
669 SubEntityRange() :
670 Base(),
671 containsPtr_(nullptr),
672 size_(0)
673 {}
674
675 std::size_t size() const
676 {
677 return size_;
678 }
679
680 bool contains(std::size_t i) const
681 {
682 return (*containsPtr_)[i];
683 }
684
685 private:
686 const SubEntityFlags* containsPtr_;
687 std::size_t size_;
688 std::size_t offset_;
689 };
690
691 using NumberRange = typename Dune::IteratorRange<const unsigned int*>;
692
693 SubEntityInfo ()
694 : numbering_( nullptr )
695 {
696 std::fill( offset_.begin(), offset_.end(), 0 );
697 }
698
699 SubEntityInfo ( const SubEntityInfo &other )
700 : offset_( other.offset_ ),
701 type_( other.type_ ),
702 containsSubentity_( other.containsSubentity_ )
703 {
704 numbering_ = allocate();
705 std::copy( other.numbering_, other.numbering_ + capacity(), numbering_ );
706 }
707
708 ~SubEntityInfo () { deallocate( numbering_ ); }
709
710 const SubEntityInfo &operator= ( const SubEntityInfo &other )
711 {
712 type_ = other.type_;
713 offset_ = other.offset_;
714
715 deallocate( numbering_ );
716 numbering_ = allocate();
717 std::copy( other.numbering_, other.numbering_ + capacity(), numbering_ );
718
719 containsSubentity_ = other.containsSubentity_;
720
721 return *this;
722 }
723
724 int size ( int cc ) const
725 {
726 assert( (cc >= 0) && (cc <= dim) );
727 return (offset_[ cc+1 ] - offset_[ cc ]);
728 }
729
730 int number ( int ii, int cc ) const
731 {
732 assert( (ii >= 0) && (ii < size( cc )) );
733 return numbering_[ offset_[ cc ] + ii ];
734 }
735
736 auto numbers ( int cc ) const
737 {
738 return SubEntityRange( numbering_ + offset_[ cc ], numbering_ + offset_[ cc+1 ], containsSubentity_[cc]);
739 }
740
741 const GeometryType &type () const { return type_; }
742
743 void initialize ( unsigned int topologyId, int codim, unsigned int i )
744 {
745 const unsigned int subId = Impl::subTopologyId( topologyId, dim, codim, i );
746 type_ = GeometryType( subId, dim-codim );
747
748 // compute offsets
749 for( int cc = 0; cc <= codim; ++cc )
750 offset_[ cc ] = 0;
751 for( int cc = codim; cc <= dim; ++cc )
752 offset_[ cc+1 ] = offset_[ cc ] + Impl::size( subId, dim-codim, cc-codim );
753
754 // compute subnumbering
755 deallocate( numbering_ );
756 numbering_ = allocate();
757 for( int cc = codim; cc <= dim; ++cc )
758 Impl::subTopologyNumbering( topologyId, dim, codim, i, cc-codim, numbering_+offset_[ cc ], numbering_+offset_[ cc+1 ] );
759
760 // initialize containsSubentity lookup-table
761 for(std::size_t cc=0; cc<= dim; ++cc)
762 {
763 containsSubentity_[cc].reset();
764 for(std::size_t idx=0; idx<std::size_t(size(cc)); ++idx)
765 containsSubentity_[cc][number(idx,cc)] = true;
766 }
767 }
768
769 protected:
770 int codim () const { return dim - type().dim(); }
771
772 unsigned int *allocate () { return (capacity() != 0 ? new unsigned int[ capacity() ] : nullptr); }
773 void deallocate ( unsigned int *ptr ) { delete[] ptr; }
774 unsigned int capacity () const { return offset_[ dim+1 ]; }
775
776 private:
777 unsigned int *numbering_;
778 std::array< unsigned int, dim+2 > offset_;
779 GeometryType type_;
780 std::array< SubEntityFlags, dim+1> containsSubentity_;
781 };
782
783
784 template< class ctype, int dim >
785 template< int codim >
786 struct ReferenceElementImplementation< ctype, dim >::CreateGeometries
787 {
788 template< int cc >
789 static typename ReferenceElements< ctype, dim-cc >::ReferenceElement
790 subRefElement( const ReferenceElementImplementation< ctype, dim > &refElement, int i, std::integral_constant< int, cc > )
791 {
792 return ReferenceElements< ctype, dim-cc >::general( refElement.type( i, cc ) );
793 }
794
796 subRefElement(const ReferenceElementImplementation< ctype, dim > &refElement,
797 [[maybe_unused]] int i, std::integral_constant<int, 0>)
798 {
799 return refElement;
800 }
801
802 static void apply ( const ReferenceElementImplementation< ctype, dim > &refElement, GeometryTable &geometries )
803 {
804 const int size = refElement.size( codim );
805 std::vector< FieldVector< ctype, dim > > origins( size );
806 std::vector< FieldMatrix< ctype, dim - codim, dim > > jacobianTransposeds( size );
807 Impl::referenceEmbeddings( refElement.type().id(), dim, codim, &(origins[ 0 ]), &(jacobianTransposeds[ 0 ]) );
808
809 std::get< codim >( geometries ).reserve( size );
810 for( int i = 0; i < size; ++i )
811 {
812 typename Codim< codim >::Geometry geometry( subRefElement( refElement, i, std::integral_constant< int, codim >() ), origins[ i ], jacobianTransposeds[ i ] );
813 std::get< codim >( geometries ).push_back( geometry );
814 }
815 }
816 };
817
818#endif // DOXYGEN
819
820 } // namespace Geo
821
822} // namespace Dune
823
824#endif // #ifndef DUNE_GEOMETRY_REFERENCEELEMENTIMPLEMENTATION_HH
An implementation of the Geometry interface for affine geometries.
Simple range between a begin and an end iterator.
Definition: iteratorrange.hh:22
Traits for type conversions and type information.
GeometryType
Type representing VTK's entity geometry types.
Definition: common.hh:132
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.
constexpr void forEach(Range &&range, F &&f)
Range based for loop.
Definition: hybridutilities.hh:256
constexpr auto max
Function object that returns the greater of the given values.
Definition: hybridutilities.hh:484
Some useful basic math stuff.
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
static constexpr T binomial(const T &n, const T &k) noexcept
calculate the binomial coefficient n over k as a constexpr
Definition: math.hh:131
constexpr std::bool_constant<((II==value)||...)> contains(std::integer_sequence< T, II... >, std::integral_constant< T, value >)
Checks whether or not a given sequence contains a value.
Definition: integersequence.hh:137
typename Container::ReferenceElement ReferenceElement
The reference element type.
Definition: referenceelements.hh:146
static const ReferenceElement & general(const GeometryType &type)
get general reference elements
Definition: referenceelements.hh:156
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 (Nov 24, 23:30, 2024)