Dune Core Modules (2.5.2)

referenceelements.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_GEOMETRY_REFERENCEELEMENTS_HH
4#define DUNE_GEOMETRY_REFERENCEELEMENTS_HH
5
6#include <cassert>
7
8#include <algorithm>
9#include <limits>
10#include <tuple>
11#include <utility>
12#include <vector>
13
14#include <dune/common/array.hh>
17#include <dune/common/hybridutilities.hh>
18#include <dune/common/std/utility.hh>
21#include <dune/common/unused.hh>
22
24#include <dune/geometry/type.hh>
25
26namespace Dune
27{
28
29 // Internal Forward Declarations
30 // -----------------------------
31
32 template< class ctype, int dim >
33 class ReferenceElementContainer;
34
35 template< class ctype, int dim >
36 struct ReferenceElements;
37
38
39
40 namespace Impl
41 {
42
44 unsigned int size ( unsigned int topologyId, int dim, int codim );
45
46
47
55 unsigned int subTopologyId ( unsigned int topologyId, int dim, int codim, unsigned int i );
56
57
58
59 // subTopologyNumbering
60 // --------------------
61
62 void subTopologyNumbering ( unsigned int topologyId, int dim, int codim, unsigned int i, int subcodim,
63 unsigned int *beginOut, unsigned int *endOut );
64
65
66
67
68 // checkInside
69 // -----------
70
71 template< class ct, int cdim >
72 inline bool
73 checkInside ( unsigned int topologyId, int dim, const FieldVector< ct, cdim > &x, ct tolerance, ct factor = ct( 1 ) )
74 {
75 assert( (dim >= 0) && (dim <= cdim) );
76 assert( topologyId < numTopologies( dim ) );
77
78 if( dim > 0 )
79 {
80 const ct baseFactor = (isPrism( topologyId, dim ) ? factor : factor - x[ dim-1 ]);
81 if( (x[ dim-1 ] > -tolerance) && (factor - x[ dim-1 ] > -tolerance) )
82 return checkInside< ct, cdim >( baseTopologyId( topologyId, dim ), dim-1, x, tolerance, baseFactor );
83 else
84 return false;
85 }
86 else
87 return true;
88 }
89
90
91
92 // referenceCorners
93 // ----------------
94
95 template< class ct, int cdim >
96 inline unsigned int
97 referenceCorners ( unsigned int topologyId, int dim, FieldVector< ct, cdim > *corners )
98 {
99 assert( (dim >= 0) && (dim <= cdim) );
100 assert( topologyId < numTopologies( dim ) );
101
102 if( dim > 0 )
103 {
104 const unsigned int nBaseCorners
105 = referenceCorners( baseTopologyId( topologyId, dim ), dim-1, corners );
106 assert( nBaseCorners == size( baseTopologyId( topologyId, dim ), dim-1, dim-1 ) );
107 if( isPrism( topologyId, dim ) )
108 {
109 std::copy( corners, corners + nBaseCorners, corners + nBaseCorners );
110 for( unsigned int i = 0; i < nBaseCorners; ++i )
111 corners[ i+nBaseCorners ][ dim-1 ] = ct( 1 );
112 return 2*nBaseCorners;
113 }
114 else
115 {
116 corners[ nBaseCorners ] = FieldVector< ct, cdim >( ct( 0 ) );
117 corners[ nBaseCorners ][ dim-1 ] = ct( 1 );
118 return nBaseCorners+1;
119 }
120 }
121 else
122 {
123 *corners = FieldVector< ct, cdim >( ct( 0 ) );
124 return 1;
125 }
126 }
127
128
129
130 // referenceVolume
131 // ---------------
132
133 unsigned long referenceVolumeInverse ( unsigned int topologyId, int dim );
134
135 template< class ct >
136 inline ct referenceVolume ( unsigned int topologyId, int dim )
137 {
138 return ct( 1 ) / ct( referenceVolumeInverse( topologyId, dim ) );
139 }
140
141
142
143 // referenceOrigins
144 // ----------------
145
146 template< class ct, int cdim >
147 inline unsigned int
148 referenceOrigins ( unsigned int topologyId, int dim, int codim, FieldVector< ct, cdim > *origins )
149 {
150 assert( (dim >= 0) && (dim <= cdim) );
151 assert( topologyId < numTopologies( dim ) );
152 assert( (codim >= 0) && (codim <= dim) );
153
154 if( codim > 0 )
155 {
156 const unsigned int baseId = baseTopologyId( topologyId, dim );
157 if( isPrism( topologyId, dim ) )
158 {
159 const unsigned int n = (codim < dim ? referenceOrigins( baseId, dim-1, codim, origins ) : 0);
160 const unsigned int m = referenceOrigins( baseId, dim-1, codim-1, origins+n );
161 for( unsigned int i = 0; i < m; ++i )
162 {
163 origins[ n+m+i ] = origins[ n+i ];
164 origins[ n+m+i ][ dim-1 ] = ct( 1 );
165 }
166 return n+2*m;
167 }
168 else
169 {
170 const unsigned int m = referenceOrigins( baseId, dim-1, codim-1, origins );
171 if( codim == dim )
172 {
173 origins[ m ] = FieldVector< ct, cdim >( ct( 0 ) );
174 origins[ m ][ dim-1 ] = ct( 1 );
175 return m+1;
176 }
177 else
178 return m+referenceOrigins( baseId, dim-1, codim, origins+m );
179 }
180 }
181 else
182 {
183 origins[ 0 ] = FieldVector< ct, cdim >( ct( 0 ) );
184 return 1;
185 }
186 }
187
188
189
190 // referenceEmbeddings
191 // -------------------
192
193 template< class ct, int cdim, int mydim >
194 inline unsigned int
195 referenceEmbeddings ( unsigned int topologyId, int dim, int codim,
196 FieldVector< ct, cdim > *origins,
197 FieldMatrix< ct, mydim, cdim > *jacobianTransposeds )
198 {
199 assert( (0 <= codim) && (codim <= dim) && (dim <= cdim) );
200 assert( (dim - codim <= mydim) && (mydim <= cdim) );
201 assert( topologyId < numTopologies( dim ) );
202
203 if( codim > 0 )
204 {
205 const unsigned int baseId = baseTopologyId( topologyId, dim );
206 if( isPrism( topologyId, dim ) )
207 {
208 const unsigned int n = (codim < dim ? referenceEmbeddings( baseId, dim-1, codim, origins, jacobianTransposeds ) : 0);
209 for( unsigned int i = 0; i < n; ++i )
210 jacobianTransposeds[ i ][ dim-codim-1 ][ dim-1 ] = ct( 1 );
211
212 const unsigned int m = referenceEmbeddings( baseId, dim-1, codim-1, origins+n, jacobianTransposeds+n );
213 std::copy( origins+n, origins+n+m, origins+n+m );
214 std::copy( jacobianTransposeds+n, jacobianTransposeds+n+m, jacobianTransposeds+n+m );
215 for( unsigned int i = 0; i < m; ++i )
216 origins[ n+m+i ][ dim-1 ] = ct( 1 );
217
218 return n+2*m;
219 }
220 else
221 {
222 const unsigned int m = referenceEmbeddings( baseId, dim-1, codim-1, origins, jacobianTransposeds );
223 if( codim == dim )
224 {
225 origins[ m ] = FieldVector< ct, cdim >( ct( 0 ) );
226 origins[ m ][ dim-1 ] = ct( 1 );
227 jacobianTransposeds[ m ] = FieldMatrix< ct, mydim, cdim >( ct( 0 ) );
228 return m+1;
229 }
230 else
231 {
232 const unsigned int n = referenceEmbeddings( baseId, dim-1, codim, origins+m, jacobianTransposeds+m );
233 for( unsigned int i = 0; i < n; ++i )
234 {
235 for( int k = 0; k < dim-1; ++k )
236 jacobianTransposeds[ m+i ][ dim-codim-1 ][ k ] = -origins[ m+i ][ k ];
237 jacobianTransposeds[ m+i ][ dim-codim-1 ][ dim-1 ] = ct( 1 );
238 }
239 return m+n;
240 }
241 }
242 }
243 else
244 {
245 origins[ 0 ] = FieldVector< ct, cdim >( ct( 0 ) );
246 jacobianTransposeds[ 0 ] = FieldMatrix< ct, mydim, cdim >( ct( 0 ) );
247 for( int k = 0; k < dim; ++k )
248 jacobianTransposeds[ 0 ][ k ][ k ] = ct( 1 );
249 return 1;
250 }
251 }
252
253
254
255 // referenceIntegrationOuterNormals
256 // --------------------------------
257
258 template< class ct, int cdim >
259 inline unsigned int
260 referenceIntegrationOuterNormals ( unsigned int topologyId, int dim,
261 const FieldVector< ct, cdim > *origins,
262 FieldVector< ct, cdim > *normals )
263 {
264 assert( (dim > 0) && (dim <= cdim) );
265 assert( topologyId < numTopologies( dim ) );
266
267 if( dim > 1 )
268 {
269 const unsigned int baseId = baseTopologyId( topologyId, dim );
270 if( isPrism( topologyId, dim ) )
271 {
272 const unsigned int numBaseFaces
273 = referenceIntegrationOuterNormals( baseId, dim-1, origins, normals );
274
275 for( unsigned int i = 0; i < 2; ++i )
276 {
277 normals[ numBaseFaces+i ] = FieldVector< ct, cdim >( ct( 0 ) );
278 normals[ numBaseFaces+i ][ dim-1 ] = ct( 2*int( i )-1 );
279 }
280
281 return numBaseFaces+2;
282 }
283 else
284 {
285 normals[ 0 ] = FieldVector< ct, cdim >( ct( 0 ) );
286 normals[ 0 ][ dim-1 ] = ct( -1 );
287
288 const unsigned int numBaseFaces
289 = referenceIntegrationOuterNormals( baseId, dim-1, origins+1, normals+1 );
290 for( unsigned int i = 1; i <= numBaseFaces; ++i )
291 normals[ i ][ dim-1 ] = normals[ i ]*origins[ i ];
292
293 return numBaseFaces+1;
294 }
295 }
296 else
297 {
298 for( unsigned int i = 0; i < 2; ++i )
299 {
300 normals[ i ] = FieldVector< ct, cdim >( ct( 0 ) );
301 normals[ i ][ 0 ] = ct( 2*int( i )-1 );
302 }
303
304 return 2;
305 }
306 }
307
308 template< class ct, int cdim >
309 inline unsigned int
310 referenceIntegrationOuterNormals ( unsigned int topologyId, int dim,
311 FieldVector< ct, cdim > *normals )
312 {
313 assert( (dim > 0) && (dim <= cdim) );
314
315 FieldVector< ct, cdim > *origins
316 = new FieldVector< ct, cdim >[ size( topologyId, dim, 1 ) ];
317 referenceOrigins( topologyId, dim, 1, origins );
318
319 const unsigned int numFaces
320 = referenceIntegrationOuterNormals( topologyId, dim, origins, normals );
321 assert( numFaces == size( topologyId, dim, 1 ) );
322
323 delete[] origins;
324
325 return numFaces;
326 }
327
328 } // namespace Impl
329
330
331
332 // ReferenceElement
333 // ----------------
334
353 template< class ctype, int dim >
355 {
357
358 friend class ReferenceElementContainer< ctype, dim >;
359
360 struct SubEntityInfo;
361
362 // make copy constructor private
363 ReferenceElement ( const This & );
364
365 ReferenceElement () {}
366
367 template< int codim > struct CreateGeometries;
368
369 public:
371 template< int codim >
372 struct Codim
373 {
375 typedef AffineGeometry< ctype, dim-codim, dim > Geometry;
376 };
377
382 int size ( int c ) const
383 {
384 assert( (c >= 0) && (c <= dim) );
385 return info_[ c ].size();
386 }
387
399 int size ( int i, int c, int cc ) const
400 {
401 assert( (i >= 0) && (i < size( c )) );
402 return info_[ c ][ i ].size( cc );
403 }
404
418 int subEntity ( int i, int c, int ii, int cc ) const
419 {
420 assert( (i >= 0) && (i < size( c )) );
421 return info_[ c ][ i ].number( ii, cc );
422 }
423
432 const GeometryType &type ( int i, int c ) const
433 {
434 assert( (i >= 0) && (i < size( c )) );
435 return info_[ c ][ i ].type();
436 }
437
439 const GeometryType &type () const { return type( 0, 0 ); }
440
450 const FieldVector< ctype, dim > &position( int i, int c ) const
451 {
452 assert( (c >= 0) && (c <= dim) );
453 return baryCenters_[ c ][ i ];
454 }
455
463 bool checkInside ( const FieldVector< ctype, dim > &local ) const
464 {
465 const ctype tolerance = ctype( 64 ) * std::numeric_limits< ctype >::epsilon();
466 return Impl::template checkInside< ctype, dim >( type().id(), dim, local, tolerance );
467 }
468
480 template< int codim >
481 typename Codim< codim >::Geometry geometry ( int i ) const
482 {
483 return std::get< codim >( geometries_ )[ i ];
484 }
485
487 ctype volume () const
488 {
489 return volume_;
490 }
491
500 {
501 assert( (face >= 0) && (face < int( integrationNormals_.size() )) );
502 return integrationNormals_[ face ];
503 }
504
505 private:
506 void initialize ( unsigned int topologyId )
507 {
508 assert( topologyId < Impl::numTopologies( dim ) );
509
510 // set up subentities
511 for( int codim = 0; codim <= dim; ++codim )
512 {
513 const unsigned int size = Impl::size( topologyId, dim, codim );
514 info_[ codim ].resize( size );
515 for( unsigned int i = 0; i < size; ++i )
516 info_[ codim ][ i ].initialize( topologyId, codim, i );
517 }
518
519 // compute corners
520 const unsigned int numVertices = size( dim );
521 baryCenters_[ dim ].resize( numVertices );
522 Impl::referenceCorners( topologyId, dim, &(baryCenters_[ dim ][ 0 ]) );
523
524 // compute barycenters
525 for( int codim = 0; codim < dim; ++codim )
526 {
527 baryCenters_[ codim ].resize( size(codim) );
528 for( int i = 0; i < size( codim ); ++i )
529 {
530 baryCenters_[ codim ][ i ] = FieldVector< ctype, dim >( ctype( 0 ) );
531 const unsigned int numCorners = size( i, codim, dim );
532 for( unsigned int j = 0; j < numCorners; ++j )
533 baryCenters_[ codim ][ i ] += baryCenters_[ dim ][ subEntity( i, codim, j, dim ) ];
534 baryCenters_[ codim ][ i ] *= ctype( 1 ) / ctype( numCorners );
535 }
536 }
537
538 // compute reference element volume
539 volume_ = Impl::template referenceVolume< ctype >( topologyId, dim );
540
541 // compute integration outer normals
542 if( dim > 0 )
543 {
544 integrationNormals_.resize( size( 1 ) );
545 Impl::referenceIntegrationOuterNormals( topologyId, dim, &(integrationNormals_[ 0 ]) );
546 }
547
548 // set up geometries
549 Hybrid::forEach( Std::make_index_sequence< dim+1 >{}, [ & ]( auto i ){ CreateGeometries< i >::apply( *this, geometries_ ); } );
550 }
551
552 template< int... codim >
553 static std::tuple< std::vector< typename Codim< codim >::Geometry >... >
554 makeGeometryTable ( std::integer_sequence< int, codim... > );
555
557 typedef decltype( makeGeometryTable( std::make_integer_sequence< int, dim+1 >() ) ) GeometryTable;
558
560 ctype volume_;
561
562 std::vector< FieldVector< ctype, dim > > baryCenters_[ dim+1 ];
563 std::vector< FieldVector< ctype, dim > > integrationNormals_;
564
566 GeometryTable geometries_;
567
568 std::vector< SubEntityInfo > info_[ dim+1 ];
569 };
570
572 template< class ctype, int dim >
573 struct ReferenceElement< ctype, dim >::SubEntityInfo
574 {
576 : numbering_( nullptr )
577 {
578 std::fill( offset_.begin(), offset_.end(), 0 );
579 }
580
581 SubEntityInfo ( const SubEntityInfo &other )
582 : offset_( other.offset_ ),
583 type_( other.type_ )
584 {
585 numbering_ = allocate();
586 std::copy( other.numbering_, other.numbering_ + capacity(), numbering_ );
587 }
588
589 ~SubEntityInfo () { deallocate( numbering_ ); }
590
591 const SubEntityInfo &operator= ( const SubEntityInfo &other )
592 {
593 type_ = other.type_;
594 offset_ = other.offset_;
595
596 deallocate( numbering_ );
597 numbering_ = allocate();
598 std::copy( other.numbering_, other.numbering_ + capacity(), numbering_ );
599
600 return *this;
601 }
602
603 int size ( int cc ) const
604 {
605 assert( (cc >= codim()) && (cc <= dim) );
606 return (offset_[ cc+1 ] - offset_[ cc ]);
607 }
608
609 int number ( int ii, int cc ) const
610 {
611 assert( (ii >= 0) && (ii < size( cc )) );
612 return numbering_[ offset_[ cc ] + ii ];
613 }
614
615 const GeometryType &type () const { return type_; }
616
617 void initialize ( unsigned int topologyId, int codim, unsigned int i )
618 {
619 const unsigned int subId = Impl::subTopologyId( topologyId, dim, codim, i );
620 type_ = GeometryType( subId, dim-codim );
621
622 // compute offsets
623 for( int cc = 0; cc <= codim; ++cc )
624 offset_[ cc ] = 0;
625 for( int cc = codim; cc <= dim; ++cc )
626 offset_[ cc+1 ] = offset_[ cc ] + Impl::size( subId, dim-codim, cc-codim );
627
628 // compute subnumbering
629 deallocate( numbering_ );
630 numbering_ = allocate();
631 for( int cc = codim; cc <= dim; ++cc )
632 Impl::subTopologyNumbering( topologyId, dim, codim, i, cc-codim, numbering_+offset_[ cc ], numbering_+offset_[ cc+1 ] );
633 }
634
635 protected:
636 int codim () const { return dim - type().dim(); }
637
638 unsigned int *allocate () { return (capacity() != 0 ? new unsigned int[ capacity() ] : nullptr); }
639 void deallocate ( unsigned int *ptr ) { delete[] ptr; }
640 unsigned int capacity () const { return offset_[ dim+1 ]; }
641
642 private:
643 unsigned int *numbering_;
644 std::array< unsigned int, dim+2 > offset_;
645 GeometryType type_;
646 };
647
648
649 template< class ctype, int dim >
650 template< int codim >
651 struct ReferenceElement< ctype, dim >::CreateGeometries
652 {
653 template< int cc >
654 static const ReferenceElement< ctype, dim-cc > &
655 subRefElement( const ReferenceElement< ctype, dim > &refElement, int i, std::integral_constant< int, cc > )
656 {
657 return ReferenceElements< ctype, dim-cc >::general( refElement.type( i, cc ) );
658 }
659
660 static const ReferenceElement< ctype, dim > &
661 subRefElement( const ReferenceElement< ctype, dim > &refElement, int i, std::integral_constant< int, 0 > )
662 {
664 return refElement;
665 }
666
667 static void apply ( const ReferenceElement< ctype, dim > &refElement, GeometryTable &geometries )
668 {
669 const int size = refElement.size( codim );
670 std::vector< FieldVector< ctype, dim > > origins( size );
671 std::vector< FieldMatrix< ctype, dim - codim, dim > > jacobianTransposeds( size );
672 Impl::referenceEmbeddings( refElement.type().id(), dim, codim, &(origins[ 0 ]), &(jacobianTransposeds[ 0 ]) );
673
674 std::get< codim >( geometries ).reserve( size );
675 for( int i = 0; i < size; ++i )
676 {
677 typename Codim< codim >::Geometry geometry( subRefElement( refElement, i, std::integral_constant< int, codim >() ), origins[ i ], jacobianTransposeds[ i ] );
678 std::get< codim >( geometries ).push_back( geometry );
679 }
680 }
681 };
682
683
684
685 // ReferenceElementContainer
686 // -------------------------
687
688 template< class ctype, int dim >
689 class ReferenceElementContainer
690 {
691 static const unsigned int numTopologies = (1u << dim);
692
693 public:
694 typedef ReferenceElement< ctype, dim > value_type;
695 typedef const value_type *const_iterator;
696
697 ReferenceElementContainer ()
698 {
699 for( unsigned int topologyId = 0; topologyId < numTopologies; ++topologyId )
700 values_[ topologyId ].initialize( topologyId );
701 }
702
703 const value_type &operator() ( const GeometryType &type ) const
704 {
705 assert( type.dim() == dim );
706 return values_[ type.id() ];
707 }
708
709 const value_type &simplex () const
710 {
711 return values_[ Impl::SimplexTopology< dim >::type::id ];
712 }
713
714 const value_type &cube () const
715 {
716 return values_[ Impl::CubeTopology< dim >::type::id ];
717 }
718
719 const value_type &pyramid () const
720 {
721 return values_[ Impl::PyramidTopology< dim >::type::id ];
722 }
723
724 const value_type &prism () const
725 {
726 return values_[ Impl::PrismTopology< dim >::type::id ];
727 }
728
729 const_iterator begin () const { return values_; }
730 const_iterator end () const { return values_ + numTopologies; }
731
732 private:
733 value_type values_[ numTopologies ];
734 };
735
736
737
738 // ReferenceElements
739 // ------------------------
740
751 template< class ctype, int dim >
753 {
755
757 static const ReferenceElement< ctype, dim > &
758 general ( const GeometryType &type )
759 {
760 return container() ( type );
761 }
762
765 {
766 return container().simplex();
767 }
768
771 {
772 return container().cube();
773 }
774
775 static Iterator begin () { return container().begin(); }
776 static Iterator end () { return container().end(); }
777
778 private:
779 DUNE_EXPORT static const ReferenceElementContainer< ctype, dim > &container ()
780 {
781 static ReferenceElementContainer< ctype, dim > container;
782 return container;
783 }
784 };
785
786} // namespace Dune
787
788#endif // #ifndef DUNE_GEOMETRY_REFERENCEELEMENTS_HH
An implementation of the Geometry interface for affine geometries.
Fallback implementation of the std::array class (a static array)
Implementation of the Geometry interface for affine geometries.
Definition: affinegeometry.hh:461
Unique label for each type of entities that can occur in DUNE grids.
Definition: type.hh:268
unsigned int dim() const
Return dimension of the type.
Definition: type.hh:565
This class provides access to geometric and topological properties of a reference element.
Definition: referenceelements.hh:355
int size(int c) const
number of subentities of codimension c
Definition: referenceelements.hh:382
int size(int i, int c, int cc) const
number of subentities of codimension cc of subentity (i,c)
Definition: referenceelements.hh:399
const FieldVector< ctype, dim > & position(int i, int c) const
position of the barycenter of entity (i,c)
Definition: referenceelements.hh:450
ctype volume() const
obtain the volume of the reference element
Definition: referenceelements.hh:487
int subEntity(int i, int c, int ii, int cc) const
obtain number of ii-th subentity with codim cc of (i,c)
Definition: referenceelements.hh:418
bool checkInside(const FieldVector< ctype, dim > &local) const
check if a coordinate is in the reference element
Definition: referenceelements.hh:463
const GeometryType & type(int i, int c) const
obtain the type of subentity (i,c)
Definition: referenceelements.hh:432
Codim< codim >::Geometry geometry(int i) const
obtain the embedding of subentity (i,codim) into the reference element
Definition: referenceelements.hh:481
const FieldVector< ctype, dim > & integrationOuterNormal(int face) const
obtain the integration outer normal of the reference element
Definition: referenceelements.hh:499
const GeometryType & type() const
obtain the type of this reference element
Definition: referenceelements.hh:439
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:314
Dune namespace.
Definition: alignment.hh:11
Collection of types depending on the codimension.
Definition: referenceelements.hh:373
AffineGeometry< ctype, dim-codim, dim > Geometry
type of geometry embedding a subentity into the reference element
Definition: referenceelements.hh:375
topological information about the subentities of a reference element
Definition: referenceelements.hh:574
Class providing access to the singletons of the reference elements.
Definition: referenceelements.hh:753
static const ReferenceElement< ctype, dim > & simplex()
get simplex reference elements
Definition: referenceelements.hh:764
static const ReferenceElement< ctype, dim > & cube()
get hypercube reference elements
Definition: referenceelements.hh:770
static const ReferenceElement< ctype, dim > & general(const GeometryType &type)
get general reference elements
Definition: referenceelements.hh:758
A unique label for each type of element that can occur in a grid.
Traits for type conversions and type information.
Definition of the DUNE_UNUSED macro for the case that config.h is not available.
#define DUNE_UNUSED_PARAMETER(parm)
A macro to mark intentionally unused function parameters with.
Definition: unused.hh:18
Definition of macros controlling symbol visibility at the ABI level.
#define DUNE_EXPORT
Export a symbol as part of the public ABI.
Definition: visibility.hh:18
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Sep 26, 22:29, 2024)