Dune Core Modules (2.7.1)

affinegeometry.hh
Go to the documentation of this file.
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_AFFINEGEOMETRY_HH
4#define DUNE_GEOMETRY_AFFINEGEOMETRY_HH
5
11#include <cmath>
12
15
16#include <dune/geometry/type.hh>
17
18namespace Dune
19{
20
21 // External Forward Declarations
22 // -----------------------------
23
24 namespace Geo
25 {
26
27 template< typename Implementation >
28 class ReferenceElement;
29
30 template< class ctype, int dim >
31 class ReferenceElementImplementation;
32
33 template< class ctype, int dim >
34 struct ReferenceElements;
35
36 }
37
38
39 namespace Impl
40 {
41
42 // FieldMatrixHelper
43 // -----------------
44
45 template< class ct >
46 struct FieldMatrixHelper
47 {
48 typedef ct ctype;
49
50 template< int m, int n >
51 static void Ax ( const FieldMatrix< ctype, m, n > &A, const FieldVector< ctype, n > &x, FieldVector< ctype, m > &ret )
52 {
53 for( int i = 0; i < m; ++i )
54 {
55 ret[ i ] = ctype( 0 );
56 for( int j = 0; j < n; ++j )
57 ret[ i ] += A[ i ][ j ] * x[ j ];
58 }
59 }
60
61 template< int m, int n >
62 static void ATx ( const FieldMatrix< ctype, m, n > &A, const FieldVector< ctype, m > &x, FieldVector< ctype, n > &ret )
63 {
64 for( int i = 0; i < n; ++i )
65 {
66 ret[ i ] = ctype( 0 );
67 for( int j = 0; j < m; ++j )
68 ret[ i ] += A[ j ][ i ] * x[ j ];
69 }
70 }
71
72 template< int m, int n, int p >
73 static void AB ( const FieldMatrix< ctype, m, n > &A, const FieldMatrix< ctype, n, p > &B, FieldMatrix< ctype, m, p > &ret )
74 {
75 for( int i = 0; i < m; ++i )
76 {
77 for( int j = 0; j < p; ++j )
78 {
79 ret[ i ][ j ] = ctype( 0 );
80 for( int k = 0; k < n; ++k )
81 ret[ i ][ j ] += A[ i ][ k ] * B[ k ][ j ];
82 }
83 }
84 }
85
86 template< int m, int n, int p >
87 static void ATBT ( const FieldMatrix< ctype, m, n > &A, const FieldMatrix< ctype, p, m > &B, FieldMatrix< ctype, n, p > &ret )
88 {
89 for( int i = 0; i < n; ++i )
90 {
91 for( int j = 0; j < p; ++j )
92 {
93 ret[ i ][ j ] = ctype( 0 );
94 for( int k = 0; k < m; ++k )
95 ret[ i ][ j ] += A[ k ][ i ] * B[ j ][ k ];
96 }
97 }
98 }
99
100 template< int m, int n >
101 static void ATA_L ( const FieldMatrix< ctype, m, n > &A, FieldMatrix< ctype, n, n > &ret )
102 {
103 for( int i = 0; i < n; ++i )
104 {
105 for( int j = 0; j <= i; ++j )
106 {
107 ret[ i ][ j ] = ctype( 0 );
108 for( int k = 0; k < m; ++k )
109 ret[ i ][ j ] += A[ k ][ i ] * A[ k ][ j ];
110 }
111 }
112 }
113
114 template< int m, int n >
115 static void ATA ( const FieldMatrix< ctype, m, n > &A, FieldMatrix< ctype, n, n > &ret )
116 {
117 for( int i = 0; i < n; ++i )
118 {
119 for( int j = 0; j <= i; ++j )
120 {
121 ret[ i ][ j ] = ctype( 0 );
122 for( int k = 0; k < m; ++k )
123 ret[ i ][ j ] += A[ k ][ i ] * A[ k ][ j ];
124 ret[ j ][ i ] = ret[ i ][ j ];
125 }
126
127 ret[ i ][ i ] = ctype( 0 );
128 for( int k = 0; k < m; ++k )
129 ret[ i ][ i ] += A[ k ][ i ] * A[ k ][ i ];
130 }
131 }
132
133 template< int m, int n >
134 static void AAT_L ( const FieldMatrix< ctype, m, n > &A, FieldMatrix< ctype, m, m > &ret )
135 {
136 /*
137 if (m==2) {
138 ret[0][0] = A[0]*A[0];
139 ret[1][1] = A[1]*A[1];
140 ret[1][0] = A[0]*A[1];
141 }
142 else
143 */
144 for( int i = 0; i < m; ++i )
145 {
146 for( int j = 0; j <= i; ++j )
147 {
148 ctype &retij = ret[ i ][ j ];
149 retij = A[ i ][ 0 ] * A[ j ][ 0 ];
150 for( int k = 1; k < n; ++k )
151 retij += A[ i ][ k ] * A[ j ][ k ];
152 }
153 }
154 }
155
156 template< int m, int n >
157 static void AAT ( const FieldMatrix< ctype, m, n > &A, FieldMatrix< ctype, m, m > &ret )
158 {
159 for( int i = 0; i < m; ++i )
160 {
161 for( int j = 0; j < i; ++j )
162 {
163 ret[ i ][ j ] = ctype( 0 );
164 for( int k = 0; k < n; ++k )
165 ret[ i ][ j ] += A[ i ][ k ] * A[ j ][ k ];
166 ret[ j ][ i ] = ret[ i ][ j ];
167 }
168 ret[ i ][ i ] = ctype( 0 );
169 for( int k = 0; k < n; ++k )
170 ret[ i ][ i ] += A[ i ][ k ] * A[ i ][ k ];
171 }
172 }
173
174 template< int n >
175 static void Lx ( const FieldMatrix< ctype, n, n > &L, const FieldVector< ctype, n > &x, FieldVector< ctype, n > &ret )
176 {
177 for( int i = 0; i < n; ++i )
178 {
179 ret[ i ] = ctype( 0 );
180 for( int j = 0; j <= i; ++j )
181 ret[ i ] += L[ i ][ j ] * x[ j ];
182 }
183 }
184
185 template< int n >
186 static void LTx ( const FieldMatrix< ctype, n, n > &L, const FieldVector< ctype, n > &x, FieldVector< ctype, n > &ret )
187 {
188 for( int i = 0; i < n; ++i )
189 {
190 ret[ i ] = ctype( 0 );
191 for( int j = i; j < n; ++j )
192 ret[ i ] += L[ j ][ i ] * x[ j ];
193 }
194 }
195
196 template< int n >
197 static void LTL ( const FieldMatrix< ctype, n, n > &L, FieldMatrix< ctype, n, n > &ret )
198 {
199 for( int i = 0; i < n; ++i )
200 {
201 for( int j = 0; j < i; ++j )
202 {
203 ret[ i ][ j ] = ctype( 0 );
204 for( int k = i; k < n; ++k )
205 ret[ i ][ j ] += L[ k ][ i ] * L[ k ][ j ];
206 ret[ j ][ i ] = ret[ i ][ j ];
207 }
208 ret[ i ][ i ] = ctype( 0 );
209 for( int k = i; k < n; ++k )
210 ret[ i ][ i ] += L[ k ][ i ] * L[ k ][ i ];
211 }
212 }
213
214 template< int n >
215 static void LLT ( const FieldMatrix< ctype, n, n > &L, FieldMatrix< ctype, n, n > &ret )
216 {
217 for( int i = 0; i < n; ++i )
218 {
219 for( int j = 0; j < i; ++j )
220 {
221 ret[ i ][ j ] = ctype( 0 );
222 for( int k = 0; k <= j; ++k )
223 ret[ i ][ j ] += L[ i ][ k ] * L[ j ][ k ];
224 ret[ j ][ i ] = ret[ i ][ j ];
225 }
226 ret[ i ][ i ] = ctype( 0 );
227 for( int k = 0; k <= i; ++k )
228 ret[ i ][ i ] += L[ i ][ k ] * L[ i ][ k ];
229 }
230 }
231
232 template< int n >
233 static bool cholesky_L ( const FieldMatrix< ctype, n, n > &A, FieldMatrix< ctype, n, n > &ret, const bool checkSingular = false )
234 {
235 for( int i = 0; i < n; ++i )
236 {
237 ctype &rii = ret[ i ][ i ];
238
239 ctype xDiag = A[ i ][ i ];
240 for( int j = 0; j < i; ++j )
241 xDiag -= ret[ i ][ j ] * ret[ i ][ j ];
242
243 // in some cases A can be singular, e.g. when checking local for
244 // outside points during checkInside
245 if( checkSingular && ! ( xDiag > ctype( 0 )) )
246 return false ;
247
248 // otherwise this should be true always
249 assert( xDiag > ctype( 0 ) );
250 rii = sqrt( xDiag );
251
252 ctype invrii = ctype( 1 ) / rii;
253 for( int k = i+1; k < n; ++k )
254 {
255 ctype x = A[ k ][ i ];
256 for( int j = 0; j < i; ++j )
257 x -= ret[ i ][ j ] * ret[ k ][ j ];
258 ret[ k ][ i ] = invrii * x;
259 }
260 }
261
262 // return true for meaning A is non-singular
263 return true;
264 }
265
266 template< int n >
267 static ctype detL ( const FieldMatrix< ctype, n, n > &L )
268 {
269 ctype det( 1 );
270 for( int i = 0; i < n; ++i )
271 det *= L[ i ][ i ];
272 return det;
273 }
274
275 template< int n >
276 static ctype invL ( FieldMatrix< ctype, n, n > &L )
277 {
278 ctype det( 1 );
279 for( int i = 0; i < n; ++i )
280 {
281 ctype &lii = L[ i ][ i ];
282 det *= lii;
283 lii = ctype( 1 ) / lii;
284 for( int j = 0; j < i; ++j )
285 {
286 ctype &lij = L[ i ][ j ];
287 ctype x = lij * L[ j ][ j ];
288 for( int k = j+1; k < i; ++k )
289 x += L[ i ][ k ] * L[ k ][ j ];
290 lij = (-lii) * x;
291 }
292 }
293 return det;
294 }
295
296 // calculates x := L^{-1} x
297 template< int n >
298 static void invLx ( FieldMatrix< ctype, n, n > &L, FieldVector< ctype, n > &x )
299 {
300 for( int i = 0; i < n; ++i )
301 {
302 for( int j = 0; j < i; ++j )
303 x[ i ] -= L[ i ][ j ] * x[ j ];
304 x[ i ] /= L[ i ][ i ];
305 }
306 }
307
308 // calculates x := L^{-T} x
309 template< int n >
310 static void invLTx ( FieldMatrix< ctype, n, n > &L, FieldVector< ctype, n > &x )
311 {
312 for( int i = n; i > 0; --i )
313 {
314 for( int j = i; j < n; ++j )
315 x[ i-1 ] -= L[ j ][ i-1 ] * x[ j ];
316 x[ i-1 ] /= L[ i-1 ][ i-1 ];
317 }
318 }
319
320 template< int n >
321 static ctype spdDetA ( const FieldMatrix< ctype, n, n > &A )
322 {
323 // return A[0][0]*A[1][1]-A[1][0]*A[1][0];
324 FieldMatrix< ctype, n, n > L;
325 cholesky_L( A, L );
326 return detL( L );
327 }
328
329 template< int n >
330 static ctype spdInvA ( FieldMatrix< ctype, n, n > &A )
331 {
332 FieldMatrix< ctype, n, n > L;
333 cholesky_L( A, L );
334 const ctype det = invL( L );
335 LTL( L, A );
336 return det;
337 }
338
339 // calculate x := A^{-1} x
340 template< int n >
341 static bool spdInvAx ( FieldMatrix< ctype, n, n > &A, FieldVector< ctype, n > &x, const bool checkSingular = false )
342 {
343 FieldMatrix< ctype, n, n > L;
344 const bool invertible = cholesky_L( A, L, checkSingular );
345 if( ! invertible ) return invertible ;
346 invLx( L, x );
347 invLTx( L, x );
348 return invertible;
349 }
350
351 template< int m, int n >
352 static ctype detATA ( const FieldMatrix< ctype, m, n > &A )
353 {
354 if( m >= n )
355 {
356 FieldMatrix< ctype, n, n > ata;
357 ATA_L( A, ata );
358 return spdDetA( ata );
359 }
360 else
361 return ctype( 0 );
362 }
363
369 template< int m, int n >
370 static ctype sqrtDetAAT ( const FieldMatrix< ctype, m, n > &A )
371 {
372 using std::abs;
373 using std::sqrt;
374 // These special cases are here not only for speed reasons:
375 // The general implementation aborts if the matrix is almost singular,
376 // and the special implementation provide a stable way to handle that case.
377 if( (n == 2) && (m == 2) )
378 {
379 // Special implementation for 2x2 matrices: faster and more stable
380 return abs( A[ 0 ][ 0 ]*A[ 1 ][ 1 ] - A[ 1 ][ 0 ]*A[ 0 ][ 1 ] );
381 }
382 else if( (n == 3) && (m == 3) )
383 {
384 // Special implementation for 3x3 matrices
385 const ctype v0 = A[ 0 ][ 1 ] * A[ 1 ][ 2 ] - A[ 1 ][ 1 ] * A[ 0 ][ 2 ];
386 const ctype v1 = A[ 0 ][ 2 ] * A[ 1 ][ 0 ] - A[ 1 ][ 2 ] * A[ 0 ][ 0 ];
387 const ctype v2 = A[ 0 ][ 0 ] * A[ 1 ][ 1 ] - A[ 1 ][ 0 ] * A[ 0 ][ 1 ];
388 return abs( v0 * A[ 2 ][ 0 ] + v1 * A[ 2 ][ 1 ] + v2 * A[ 2 ][ 2 ] );
389 }
390 else if ( (n == 3) && (m == 2) )
391 {
392 // Special implementation for 2x3 matrices
393 const ctype v0 = A[ 0 ][ 0 ] * A[ 1 ][ 1 ] - A[ 0 ][ 1 ] * A[ 1 ][ 0 ];
394 const ctype v1 = A[ 0 ][ 0 ] * A[ 1 ][ 2 ] - A[ 1 ][ 0 ] * A[ 0 ][ 2 ];
395 const ctype v2 = A[ 0 ][ 1 ] * A[ 1 ][ 2 ] - A[ 0 ][ 2 ] * A[ 1 ][ 1 ];
396 return sqrt( v0*v0 + v1*v1 + v2*v2);
397 }
398 else if( n >= m )
399 {
400 // General case
401 FieldMatrix< ctype, m, m > aat;
402 AAT_L( A, aat );
403 return spdDetA( aat );
404 }
405 else
406 return ctype( 0 );
407 }
408
409 // A^{-1}_L = (A^T A)^{-1} A^T
410 // => A^{-1}_L A = I
411 template< int m, int n >
412 static ctype leftInvA ( const FieldMatrix< ctype, m, n > &A, FieldMatrix< ctype, n, m > &ret )
413 {
414 static_assert((m >= n), "Matrix has no left inverse.");
415 FieldMatrix< ctype, n, n > ata;
416 ATA_L( A, ata );
417 const ctype det = spdInvA( ata );
418 ATBT( ata, A, ret );
419 return det;
420 }
421
422 template< int m, int n >
423 static void leftInvAx ( const FieldMatrix< ctype, m, n > &A, const FieldVector< ctype, m > &x, FieldVector< ctype, n > &y )
424 {
425 static_assert((m >= n), "Matrix has no left inverse.");
426 FieldMatrix< ctype, n, n > ata;
427 ATx( A, x, y );
428 ATA_L( A, ata );
429 spdInvAx( ata, y );
430 }
431
433 template< int m, int n >
434 static ctype rightInvA ( const FieldMatrix< ctype, m, n > &A, FieldMatrix< ctype, n, m > &ret )
435 {
436 static_assert((n >= m), "Matrix has no right inverse.");
437 using std::abs;
438 if( (n == 2) && (m == 2) )
439 {
440 const ctype det = (A[ 0 ][ 0 ]*A[ 1 ][ 1 ] - A[ 1 ][ 0 ]*A[ 0 ][ 1 ]);
441 const ctype detInv = ctype( 1 ) / det;
442 ret[ 0 ][ 0 ] = A[ 1 ][ 1 ] * detInv;
443 ret[ 1 ][ 1 ] = A[ 0 ][ 0 ] * detInv;
444 ret[ 1 ][ 0 ] = -A[ 1 ][ 0 ] * detInv;
445 ret[ 0 ][ 1 ] = -A[ 0 ][ 1 ] * detInv;
446 return abs( det );
447 }
448 else
449 {
450 FieldMatrix< ctype, m , m > aat;
451 AAT_L( A, aat );
452 const ctype det = spdInvA( aat );
453 ATBT( A , aat , ret );
454 return det;
455 }
456 }
457
458 template< int m, int n >
459 static bool xTRightInvA ( const FieldMatrix< ctype, m, n > &A, const FieldVector< ctype, n > &x, FieldVector< ctype, m > &y )
460 {
461 static_assert((n >= m), "Matrix has no right inverse.");
462 FieldMatrix< ctype, m, m > aat;
463 Ax( A, x, y );
464 AAT_L( A, aat );
465 // check whether aat is singular and return true if non-singular
466 return spdInvAx( aat, y, true );
467 }
468 };
469
470 } // namespace Impl
471
472
473
479 template< class ct, int mydim, int cdim>
481 {
482 public:
483
485 typedef ct ctype;
486
488 static const int mydimension= mydim;
489
491 static const int coorddimension = cdim;
492
495
498
500 typedef ctype Volume;
501
504
507
508 private:
511
513
514 // Helper class to compute a matrix pseudo inverse
515 typedef Impl::FieldMatrixHelper< ct > MatrixHelper;
516
517 public:
519 AffineGeometry ( const ReferenceElement &refElement, const GlobalCoordinate &origin,
520 const JacobianTransposed &jt )
521 : refElement_(refElement), origin_(origin), jacobianTransposed_(jt)
522 {
523 integrationElement_ = MatrixHelper::template rightInvA< mydimension, coorddimension >( jacobianTransposed_, jacobianInverseTransposed_ );
524 }
525
528 const JacobianTransposed &jt )
529 : AffineGeometry(ReferenceElements::general( gt ), origin, jt)
530 { }
531
533 template< class CoordVector >
534 AffineGeometry ( const ReferenceElement &refElement, const CoordVector &coordVector )
535 : refElement_(refElement), origin_(coordVector[0])
536 {
537 for( int i = 0; i < mydimension; ++i )
538 jacobianTransposed_[ i ] = coordVector[ i+1 ] - origin_;
539 integrationElement_ = MatrixHelper::template rightInvA< mydimension, coorddimension >( jacobianTransposed_, jacobianInverseTransposed_ );
540 }
541
543 template< class CoordVector >
544 AffineGeometry ( Dune::GeometryType gt, const CoordVector &coordVector )
545 : AffineGeometry(ReferenceElements::general( gt ), coordVector)
546 { }
547
549 bool affine () const { return true; }
550
552 Dune::GeometryType type () const { return refElement_.type(); }
553
555 int corners () const { return refElement_.size( mydimension ); }
556
558 GlobalCoordinate corner ( int i ) const
559 {
560 return global( refElement_.position( i, mydimension ) );
561 }
562
564 GlobalCoordinate center () const { return global( refElement_.position( 0, 0 ) ); }
565
573 {
574 GlobalCoordinate global( origin_ );
575 jacobianTransposed_.umtv( local, global );
576 return global;
577 }
578
593 {
594 LocalCoordinate local;
595 jacobianInverseTransposed_.mtv( global - origin_, local );
596 return local;
597 }
598
610 {
612 return integrationElement_;
613 }
614
616 Volume volume () const
617 {
618 return integrationElement_ * refElement_.volume();
619 }
620
628 {
630 return jacobianTransposed_;
631 }
632
640 {
642 return jacobianInverseTransposed_;
643 }
644
645 friend ReferenceElement referenceElement ( const AffineGeometry &geometry )
646 {
647 return geometry.refElement_;
648 }
649
650 private:
651 ReferenceElement refElement_;
652 GlobalCoordinate origin_;
653 JacobianTransposed jacobianTransposed_;
654 JacobianInverseTransposed jacobianInverseTransposed_;
655 ctype integrationElement_;
656 };
657
658} // namespace Dune
659
660#endif // #ifndef DUNE_GEOMETRY_AFFINEGEOMETRY_HH
Implementation of the Geometry interface for affine geometries.
Definition: affinegeometry.hh:481
AffineGeometry(const ReferenceElement &refElement, const CoordVector &coordVector)
Create affine geometry from reference element and a vector of vertex coordinates.
Definition: affinegeometry.hh:534
AffineGeometry(Dune::GeometryType gt, const GlobalCoordinate &origin, const JacobianTransposed &jt)
Create affine geometry from GeometryType, one vertex, and the Jacobian matrix.
Definition: affinegeometry.hh:527
FieldVector< ctype, mydimension > LocalCoordinate
Type for local coordinate vector.
Definition: affinegeometry.hh:494
Dune::GeometryType type() const
Obtain the type of the reference element.
Definition: affinegeometry.hh:552
static const int mydimension
Dimension of the geometry.
Definition: affinegeometry.hh:488
AffineGeometry(const ReferenceElement &refElement, const GlobalCoordinate &origin, const JacobianTransposed &jt)
Create affine geometry from reference element, one vertex, and the Jacobian matrix.
Definition: affinegeometry.hh:519
ctype Volume
Type used for volume.
Definition: affinegeometry.hh:500
AffineGeometry(Dune::GeometryType gt, const CoordVector &coordVector)
Create affine geometry from GeometryType and a vector of vertex coordinates.
Definition: affinegeometry.hh:544
ctype integrationElement(const LocalCoordinate &local) const
Obtain the integration element.
Definition: affinegeometry.hh:609
const JacobianInverseTransposed & jacobianInverseTransposed(const LocalCoordinate &local) const
Obtain the transposed of the Jacobian's inverse.
Definition: affinegeometry.hh:639
FieldMatrix< ctype, mydimension, coorddimension > JacobianTransposed
Type for the transposed Jacobian matrix.
Definition: affinegeometry.hh:503
GlobalCoordinate corner(int i) const
Obtain coordinates of the i-th corner.
Definition: affinegeometry.hh:558
int corners() const
Obtain number of corners of the corresponding reference element.
Definition: affinegeometry.hh:555
FieldMatrix< ctype, coorddimension, mydimension > JacobianInverseTransposed
Type for the transposed inverse Jacobian matrix.
Definition: affinegeometry.hh:506
static const int coorddimension
Dimension of the world space.
Definition: affinegeometry.hh:491
GlobalCoordinate global(const LocalCoordinate &local) const
Evaluate the mapping.
Definition: affinegeometry.hh:572
GlobalCoordinate center() const
Obtain the centroid of the mapping's image.
Definition: affinegeometry.hh:564
ct ctype
Type used for coordinates.
Definition: affinegeometry.hh:485
FieldVector< ctype, coorddimension > GlobalCoordinate
Type for coordinate vector in world space.
Definition: affinegeometry.hh:497
bool affine() const
Always true: this is an affine geometry.
Definition: affinegeometry.hh:549
const JacobianTransposed & jacobianTransposed(const LocalCoordinate &local) const
Obtain the transposed of the Jacobian.
Definition: affinegeometry.hh:627
Volume volume() const
Obtain the volume of the element.
Definition: affinegeometry.hh:616
void mtv(const X &x, Y &y) const
y = A^T x
Definition: densematrix.hh:414
void umtv(const X &x, Y &y) const
y += A^T x
Definition: densematrix.hh:446
vector space out of a tensor product of fields.
Definition: fvector.hh:96
This class provides access to geometric and topological properties of a reference element.
Definition: referenceelement.hh:56
CoordinateField volume() const
obtain the volume of the reference element
Definition: referenceelement.hh:245
decltype(auto) type(int i, int c) const
obtain the type of subentity (i,c)
Definition: referenceelement.hh:175
int size(int c) const
number of subentities of codimension c
Definition: referenceelement.hh:98
decltype(auto) position(int i, int c) const
position of the barycenter of entity (i,c)
Definition: referenceelement.hh:207
Unique label for each type of entities that can occur in DUNE grids.
Definition: type.hh:279
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.
#define DUNE_UNUSED_PARAMETER(parm)
A macro to mark intentionally unused function parameters with.
Definition: unused.hh:25
bool gt(const T &first, const T &second, typename EpsilonType< T >::Type epsilon)
test if first greater than second
Definition: float_cmp.cc:156
unspecified-type ReferenceElement
Returns the type of reference element for the argument type T.
Definition: referenceelements.hh:495
Dune namespace.
Definition: alignedallocator.hh:14
Class providing access to the singletons of the reference elements.
Definition: referenceelements.hh:168
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 12, 23:30, 2024)