Dune Core Modules (2.9.0)

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 (C) 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( codim > 0 )
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
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
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
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
271
272
273 // referenceIntegrationOuterNormals
274 // --------------------------------
275
276 template< class ct, int cdim >
277 inline unsigned int
278 referenceIntegrationOuterNormals ( unsigned int topologyId, int dim,
279 const FieldVector< ct, cdim > *origins,
280 FieldVector< ct, cdim > *normals )
281 {
282 assert( (dim > 0) && (dim <= cdim) );
283 assert( topologyId < numTopologies( dim ) );
284
285 if( dim > 1 )
286 {
287 const unsigned int baseId = baseTopologyId( topologyId, dim );
288 if( isPrism( topologyId, dim ) )
289 {
290 const unsigned int numBaseFaces
291 = referenceIntegrationOuterNormals( baseId, dim-1, origins, normals );
292
293 for( unsigned int i = 0; i < 2; ++i )
294 {
295 normals[ numBaseFaces+i ] = FieldVector< ct, cdim >( ct( 0 ) );
296 normals[ numBaseFaces+i ][ dim-1 ] = ct( 2*int( i )-1 );
297 }
298
299 return numBaseFaces+2;
300 }
301 else
302 {
303 normals[ 0 ] = FieldVector< ct, cdim >( ct( 0 ) );
304 normals[ 0 ][ dim-1 ] = ct( -1 );
305
306 const unsigned int numBaseFaces
307 = referenceIntegrationOuterNormals( baseId, dim-1, origins+1, normals+1 );
308 for( unsigned int i = 1; i <= numBaseFaces; ++i )
309 normals[ i ][ dim-1 ] = normals[ i ]*origins[ i ];
310
311 return numBaseFaces+1;
312 }
313 }
314 else
315 {
316 for( unsigned int i = 0; i < 2; ++i )
317 {
318 normals[ i ] = FieldVector< ct, cdim >( ct( 0 ) );
319 normals[ i ][ 0 ] = ct( 2*int( i )-1 );
320 }
321
322 return 2;
323 }
324 }
325
326 template< class ct, int cdim >
327 inline unsigned int
328 referenceIntegrationOuterNormals ( unsigned int topologyId, int dim,
329 FieldVector< ct, cdim > *normals )
330 {
331 assert( (dim > 0) && (dim <= cdim) );
332
333 FieldVector< ct, cdim > *origins
334 = new FieldVector< ct, cdim >[ size( topologyId, dim, 1 ) ];
335 referenceOrigins( topologyId, dim, 1, origins );
336
337 const unsigned int numFaces
338 = referenceIntegrationOuterNormals( topologyId, dim, origins, normals );
339 assert( numFaces == size( topologyId, dim, 1 ) );
340
341 delete[] origins;
342
343 return numFaces;
344 }
345
346 } // namespace Impl
347
348
349
350 // ReferenceElement
351 // ----------------
352
371 template< class ctype_, int dim >
372 class ReferenceElementImplementation
373 {
374
375 public:
376
378 using ctype = ctype_;
379
381 using CoordinateField = ctype;
382
384 using Coordinate = Dune::FieldVector<ctype,dim>;
385
387 static constexpr int dimension = dim;
388
390 typedef ctype Volume;
391
392 private:
393
394 friend class Impl::ReferenceElementContainer< ctype, dim >;
395
396 struct SubEntityInfo;
397
398 template< int codim > struct CreateGeometries;
399
400 public:
402 template< int codim >
403 struct Codim
404 {
406 typedef AffineGeometry< ctype, dim-codim, dim > Geometry;
407 };
408
409 // ReferenceElement cannot be copied.
410 ReferenceElementImplementation ( const ReferenceElementImplementation& ) = delete;
411
412 // ReferenceElementImplementation cannot be copied.
413 ReferenceElementImplementation& operator= ( const ReferenceElementImplementation& ) = delete;
414
415 // ReferenceElementImplementation is default-constructible (required for storage in std::array)
416 ReferenceElementImplementation () = default;
417
422 int size ( int c ) const
423 {
424 assert( (c >= 0) && (c <= dim) );
425 return info_[ c ].size();
426 }
427
439 int size ( int i, int c, int cc ) const
440 {
441 assert( (i >= 0) && (i < size( c )) );
442 return info_[ c ][ i ].size( cc );
443 }
444
458 int subEntity ( int i, int c, int ii, int cc ) const
459 {
460 assert( (i >= 0) && (i < size( c )) );
461 return info_[ c ][ i ].number( ii, cc );
462 }
463
479 auto subEntities ( int i, int c, int cc ) const
480 {
481 assert( (i >= 0) && (i < size( c )) );
482 return info_[ c ][ i ].numbers( cc );
483 }
484
493 const GeometryType &type ( int i, int c ) const
494 {
495 assert( (i >= 0) && (i < size( c )) );
496 return info_[ c ][ i ].type();
497 }
498
500 const GeometryType &type () const { return type( 0, 0 ); }
501
511 const Coordinate &position( int i, int c ) const
512 {
513 assert( (c >= 0) && (c <= dim) );
514 return baryCenters_[ c ][ i ];
515 }
516
524 bool checkInside ( const Coordinate &local ) const
525 {
526 const ctype tolerance = ctype( 64 ) * std::numeric_limits< ctype >::epsilon();
527 return Impl::template checkInside< ctype, dim >( type().id(), dim, local, tolerance );
528 }
529
541 template< int codim >
542 typename Codim< codim >::Geometry geometry ( int i ) const
543 {
544 return std::get< codim >( geometries_ )[ i ];
545 }
546
548 Volume volume () const
549 {
550 return volume_;
551 }
552
560 const Coordinate &integrationOuterNormal ( int face ) const
561 {
562 assert( (face >= 0) && (face < int( integrationNormals_.size() )) );
563 return integrationNormals_[ face ];
564 }
565
566 private:
567 void initialize ( unsigned int topologyId )
568 {
569 assert( topologyId < Impl::numTopologies( dim ) );
570
571 // set up subentities
572 for( int codim = 0; codim <= dim; ++codim )
573 {
574 const unsigned int size = Impl::size( topologyId, dim, codim );
575 info_[ codim ].resize( size );
576 for( unsigned int i = 0; i < size; ++i )
577 info_[ codim ][ i ].initialize( topologyId, codim, i );
578 }
579
580 // compute corners
581 const unsigned int numVertices = size( dim );
582 baryCenters_[ dim ].resize( numVertices );
583 Impl::referenceCorners( topologyId, dim, &(baryCenters_[ dim ][ 0 ]) );
584
585 // compute barycenters
586 for( int codim = 0; codim < dim; ++codim )
587 {
588 baryCenters_[ codim ].resize( size(codim) );
589 for( int i = 0; i < size( codim ); ++i )
590 {
591 baryCenters_[ codim ][ i ] = Coordinate( ctype( 0 ) );
592 const unsigned int numCorners = size( i, codim, dim );
593 for( unsigned int j = 0; j < numCorners; ++j )
594 baryCenters_[ codim ][ i ] += baryCenters_[ dim ][ subEntity( i, codim, j, dim ) ];
595 baryCenters_[ codim ][ i ] *= ctype( 1 ) / ctype( numCorners );
596 }
597 }
598
599 // compute reference element volume
600 volume_ = Impl::template referenceVolume< ctype >( topologyId, dim );
601
602 // compute integration outer normals
603 if( dim > 0 )
604 {
605 integrationNormals_.resize( size( 1 ) );
606 Impl::referenceIntegrationOuterNormals( topologyId, dim, &(integrationNormals_[ 0 ]) );
607 }
608
609 // set up geometries
610 Hybrid::forEach( std::make_index_sequence< dim+1 >{}, [ & ]( auto i ){ CreateGeometries< i >::apply( *this, geometries_ ); } );
611 }
612
613 template< int... codim >
614 static std::tuple< std::vector< typename Codim< codim >::Geometry >... >
615 makeGeometryTable ( std::integer_sequence< int, codim... > );
616
618 typedef decltype( makeGeometryTable( std::make_integer_sequence< int, dim+1 >() ) ) GeometryTable;
619
621 ctype volume_;
622
623 std::vector< Coordinate > baryCenters_[ dim+1 ];
624 std::vector< Coordinate > integrationNormals_;
625
627 GeometryTable geometries_;
628
629 std::vector< SubEntityInfo > info_[ dim+1 ];
630 };
631
633 template< class ctype, int dim >
634 struct ReferenceElementImplementation< ctype, dim >::SubEntityInfo
635 {
636 // Compute upper bound for the number of subsentities.
637 // If someone knows an explicit formal feel free to
638 // implement it here.
639 static constexpr std::size_t maxSubEntityCount()
640 {
641 std::size_t maxCount=0;
642 for(std::size_t codim=0; codim<=dim; ++codim)
643 maxCount = std::max(maxCount, binomial(std::size_t(dim),codim)*(1 << codim));
644 return maxCount;
645 }
646
647 using SubEntityFlags = std::bitset<maxSubEntityCount()>;
648
649 class SubEntityRange
650 : public Dune::IteratorRange<const unsigned int*>
651 {
652 using Base = typename Dune::IteratorRange<const unsigned int*>;
653
654 public:
655
656 using iterator = Base::iterator;
657 using const_iterator = Base::const_iterator;
658
659 SubEntityRange(const iterator& begin, const iterator& end, const SubEntityFlags& contains) :
660 Base(begin, end),
661 containsPtr_(&contains),
662 size_(end-begin)
663 {}
664
665 SubEntityRange() :
666 Base(),
667 containsPtr_(nullptr),
668 size_(0)
669 {}
670
671 std::size_t size() const
672 {
673 return size_;
674 }
675
676 bool contains(std::size_t i) const
677 {
678 return (*containsPtr_)[i];
679 }
680
681 private:
682 const SubEntityFlags* containsPtr_;
683 std::size_t size_;
684 std::size_t offset_;
685 };
686
687 using NumberRange = typename Dune::IteratorRange<const unsigned int*>;
688
689 SubEntityInfo ()
690 : numbering_( nullptr )
691 {
692 std::fill( offset_.begin(), offset_.end(), 0 );
693 }
694
695 SubEntityInfo ( const SubEntityInfo &other )
696 : offset_( other.offset_ ),
697 type_( other.type_ ),
698 containsSubentity_( other.containsSubentity_ )
699 {
700 numbering_ = allocate();
701 std::copy( other.numbering_, other.numbering_ + capacity(), numbering_ );
702 }
703
704 ~SubEntityInfo () { deallocate( numbering_ ); }
705
706 const SubEntityInfo &operator= ( const SubEntityInfo &other )
707 {
708 type_ = other.type_;
709 offset_ = other.offset_;
710
711 deallocate( numbering_ );
712 numbering_ = allocate();
713 std::copy( other.numbering_, other.numbering_ + capacity(), numbering_ );
714
715 containsSubentity_ = other.containsSubentity_;
716
717 return *this;
718 }
719
720 int size ( int cc ) const
721 {
722 assert( (cc >= 0) && (cc <= dim) );
723 return (offset_[ cc+1 ] - offset_[ cc ]);
724 }
725
726 int number ( int ii, int cc ) const
727 {
728 assert( (ii >= 0) && (ii < size( cc )) );
729 return numbering_[ offset_[ cc ] + ii ];
730 }
731
732 auto numbers ( int cc ) const
733 {
734 return SubEntityRange( numbering_ + offset_[ cc ], numbering_ + offset_[ cc+1 ], containsSubentity_[cc]);
735 }
736
737 const GeometryType &type () const { return type_; }
738
739 void initialize ( unsigned int topologyId, int codim, unsigned int i )
740 {
741 const unsigned int subId = Impl::subTopologyId( topologyId, dim, codim, i );
742 type_ = GeometryType( subId, dim-codim );
743
744 // compute offsets
745 for( int cc = 0; cc <= codim; ++cc )
746 offset_[ cc ] = 0;
747 for( int cc = codim; cc <= dim; ++cc )
748 offset_[ cc+1 ] = offset_[ cc ] + Impl::size( subId, dim-codim, cc-codim );
749
750 // compute subnumbering
751 deallocate( numbering_ );
752 numbering_ = allocate();
753 for( int cc = codim; cc <= dim; ++cc )
754 Impl::subTopologyNumbering( topologyId, dim, codim, i, cc-codim, numbering_+offset_[ cc ], numbering_+offset_[ cc+1 ] );
755
756 // initialize containsSubentity lookup-table
757 for(std::size_t cc=0; cc<= dim; ++cc)
758 {
759 containsSubentity_[cc].reset();
760 for(std::size_t idx=0; idx<std::size_t(size(cc)); ++idx)
761 containsSubentity_[cc][number(idx,cc)] = true;
762 }
763 }
764
765 protected:
766 int codim () const { return dim - type().dim(); }
767
768 unsigned int *allocate () { return (capacity() != 0 ? new unsigned int[ capacity() ] : nullptr); }
769 void deallocate ( unsigned int *ptr ) { delete[] ptr; }
770 unsigned int capacity () const { return offset_[ dim+1 ]; }
771
772 private:
773 unsigned int *numbering_;
774 std::array< unsigned int, dim+2 > offset_;
775 GeometryType type_;
776 std::array< SubEntityFlags, dim+1> containsSubentity_;
777 };
778
779
780 template< class ctype, int dim >
781 template< int codim >
782 struct ReferenceElementImplementation< ctype, dim >::CreateGeometries
783 {
784 template< int cc >
785 static typename ReferenceElements< ctype, dim-cc >::ReferenceElement
786 subRefElement( const ReferenceElementImplementation< ctype, dim > &refElement, int i, std::integral_constant< int, cc > )
787 {
788 return ReferenceElements< ctype, dim-cc >::general( refElement.type( i, cc ) );
789 }
790
792 subRefElement(const ReferenceElementImplementation< ctype, dim > &refElement,
793 [[maybe_unused]] int i, std::integral_constant<int, 0>)
794 {
795 return refElement;
796 }
797
798 static void apply ( const ReferenceElementImplementation< ctype, dim > &refElement, GeometryTable &geometries )
799 {
800 const int size = refElement.size( codim );
801 std::vector< FieldVector< ctype, dim > > origins( size );
802 std::vector< FieldMatrix< ctype, dim - codim, dim > > jacobianTransposeds( size );
803 Impl::referenceEmbeddings( refElement.type().id(), dim, codim, &(origins[ 0 ]), &(jacobianTransposeds[ 0 ]) );
804
805 std::get< codim >( geometries ).reserve( size );
806 for( int i = 0; i < size; ++i )
807 {
808 typename Codim< codim >::Geometry geometry( subRefElement( refElement, i, std::integral_constant< int, codim >() ), origins[ i ], jacobianTransposeds[ i ] );
809 std::get< codim >( geometries ).push_back( geometry );
810 }
811 }
812 };
813
814#endif // DOXYGEN
815
816 } // namespace Geo
817
818} // namespace Dune
819
820#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.
unspecified-type ReferenceElement
Returns the type of reference element for the argument type T.
Definition: referenceelements.hh:497
constexpr void forEach(Range &&range, F &&f)
Range based for loop.
Definition: hybridutilities.hh:268
Some useful basic math stuff.
Dune namespace.
Definition: alignedallocator.hh:13
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
typename Container::ReferenceElement ReferenceElement
The reference element type.
Definition: referenceelements.hh:188
static const ReferenceElement & general(const GeometryType &type)
get general reference elements
Definition: referenceelements.hh:198
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 (Dec 21, 23:30, 2024)