Dune Core Modules (2.6.0)

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 void cholesky_L ( const FieldMatrix< ctype, n, n > &A, FieldMatrix< ctype, n, n > &ret )
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 assert( xDiag > ctype( 0 ) );
243 rii = sqrt( xDiag );
244
245 ctype invrii = ctype( 1 ) / rii;
246 for( int k = i+1; k < n; ++k )
247 {
248 ctype x = A[ k ][ i ];
249 for( int j = 0; j < i; ++j )
250 x -= ret[ i ][ j ] * ret[ k ][ j ];
251 ret[ k ][ i ] = invrii * x;
252 }
253 }
254 }
255
256 template< int n >
257 static ctype detL ( const FieldMatrix< ctype, n, n > &L )
258 {
259 ctype det( 1 );
260 for( int i = 0; i < n; ++i )
261 det *= L[ i ][ i ];
262 return det;
263 }
264
265 template< int n >
266 static ctype invL ( FieldMatrix< ctype, n, n > &L )
267 {
268 ctype det( 1 );
269 for( int i = 0; i < n; ++i )
270 {
271 ctype &lii = L[ i ][ i ];
272 det *= lii;
273 lii = ctype( 1 ) / lii;
274 for( int j = 0; j < i; ++j )
275 {
276 ctype &lij = L[ i ][ j ];
277 ctype x = lij * L[ j ][ j ];
278 for( int k = j+1; k < i; ++k )
279 x += L[ i ][ k ] * L[ k ][ j ];
280 lij = (-lii) * x;
281 }
282 }
283 return det;
284 }
285
286 // calculates x := L^{-1} x
287 template< int n >
288 static void invLx ( FieldMatrix< ctype, n, n > &L, FieldVector< ctype, n > &x )
289 {
290 for( int i = 0; i < n; ++i )
291 {
292 for( int j = 0; j < i; ++j )
293 x[ i ] -= L[ i ][ j ] * x[ j ];
294 x[ i ] /= L[ i ][ i ];
295 }
296 }
297
298 // calculates x := L^{-T} x
299 template< int n >
300 static void invLTx ( FieldMatrix< ctype, n, n > &L, FieldVector< ctype, n > &x )
301 {
302 for( int i = n; i > 0; --i )
303 {
304 for( int j = i; j < n; ++j )
305 x[ i-1 ] -= L[ j ][ i-1 ] * x[ j ];
306 x[ i-1 ] /= L[ i-1 ][ i-1 ];
307 }
308 }
309
310 template< int n >
311 static ctype spdDetA ( const FieldMatrix< ctype, n, n > &A )
312 {
313 // return A[0][0]*A[1][1]-A[1][0]*A[1][0];
314 FieldMatrix< ctype, n, n > L;
315 cholesky_L( A, L );
316 return detL( L );
317 }
318
319 template< int n >
320 static ctype spdInvA ( FieldMatrix< ctype, n, n > &A )
321 {
322 FieldMatrix< ctype, n, n > L;
323 cholesky_L( A, L );
324 const ctype det = invL( L );
325 LTL( L, A );
326 return det;
327 }
328
329 // calculate x := A^{-1} x
330 template< int n >
331 static void spdInvAx ( FieldMatrix< ctype, n, n > &A, FieldVector< ctype, n > &x )
332 {
333 FieldMatrix< ctype, n, n > L;
334 cholesky_L( A, L );
335 invLx( L, x );
336 invLTx( L, x );
337 }
338
339 template< int m, int n >
340 static ctype detATA ( const FieldMatrix< ctype, m, n > &A )
341 {
342 if( m >= n )
343 {
344 FieldMatrix< ctype, n, n > ata;
345 ATA_L( A, ata );
346 return spdDetA( ata );
347 }
348 else
349 return ctype( 0 );
350 }
351
357 template< int m, int n >
358 static ctype sqrtDetAAT ( const FieldMatrix< ctype, m, n > &A )
359 {
360 using std::abs;
361 using std::sqrt;
362 // These special cases are here not only for speed reasons:
363 // The general implementation aborts if the matrix is almost singular,
364 // and the special implementation provide a stable way to handle that case.
365 if( (n == 2) && (m == 2) )
366 {
367 // Special implementation for 2x2 matrices: faster and more stable
368 return abs( A[ 0 ][ 0 ]*A[ 1 ][ 1 ] - A[ 1 ][ 0 ]*A[ 0 ][ 1 ] );
369 }
370 else if( (n == 3) && (m == 3) )
371 {
372 // Special implementation for 3x3 matrices
373 const ctype v0 = A[ 0 ][ 1 ] * A[ 1 ][ 2 ] - A[ 1 ][ 1 ] * A[ 0 ][ 2 ];
374 const ctype v1 = A[ 0 ][ 2 ] * A[ 1 ][ 0 ] - A[ 1 ][ 2 ] * A[ 0 ][ 0 ];
375 const ctype v2 = A[ 0 ][ 0 ] * A[ 1 ][ 1 ] - A[ 1 ][ 0 ] * A[ 0 ][ 1 ];
376 return abs( v0 * A[ 2 ][ 0 ] + v1 * A[ 2 ][ 1 ] + v2 * A[ 2 ][ 2 ] );
377 }
378 else if ( (n == 3) && (m == 2) )
379 {
380 // Special implementation for 2x3 matrices
381 const ctype v0 = A[ 0 ][ 0 ] * A[ 1 ][ 1 ] - A[ 0 ][ 1 ] * A[ 1 ][ 0 ];
382 const ctype v1 = A[ 0 ][ 0 ] * A[ 1 ][ 2 ] - A[ 1 ][ 0 ] * A[ 0 ][ 2 ];
383 const ctype v2 = A[ 0 ][ 1 ] * A[ 1 ][ 2 ] - A[ 0 ][ 2 ] * A[ 1 ][ 1 ];
384 return sqrt( v0*v0 + v1*v1 + v2*v2);
385 }
386 else if( n >= m )
387 {
388 // General case
389 FieldMatrix< ctype, m, m > aat;
390 AAT_L( A, aat );
391 return spdDetA( aat );
392 }
393 else
394 return ctype( 0 );
395 }
396
397 // A^{-1}_L = (A^T A)^{-1} A^T
398 // => A^{-1}_L A = I
399 template< int m, int n >
400 static ctype leftInvA ( const FieldMatrix< ctype, m, n > &A, FieldMatrix< ctype, n, m > &ret )
401 {
402 static_assert((m >= n), "Matrix has no left inverse.");
403 FieldMatrix< ctype, n, n > ata;
404 ATA_L( A, ata );
405 const ctype det = spdInvA( ata );
406 ATBT( ata, A, ret );
407 return det;
408 }
409
410 template< int m, int n >
411 static void leftInvAx ( const FieldMatrix< ctype, m, n > &A, const FieldVector< ctype, m > &x, FieldVector< ctype, n > &y )
412 {
413 static_assert((m >= n), "Matrix has no left inverse.");
414 FieldMatrix< ctype, n, n > ata;
415 ATx( A, x, y );
416 ATA_L( A, ata );
417 spdInvAx( ata, y );
418 }
419
421 template< int m, int n >
422 static ctype rightInvA ( const FieldMatrix< ctype, m, n > &A, FieldMatrix< ctype, n, m > &ret )
423 {
424 static_assert((n >= m), "Matrix has no right inverse.");
425 using std::abs;
426 if( (n == 2) && (m == 2) )
427 {
428 const ctype det = (A[ 0 ][ 0 ]*A[ 1 ][ 1 ] - A[ 1 ][ 0 ]*A[ 0 ][ 1 ]);
429 const ctype detInv = ctype( 1 ) / det;
430 ret[ 0 ][ 0 ] = A[ 1 ][ 1 ] * detInv;
431 ret[ 1 ][ 1 ] = A[ 0 ][ 0 ] * detInv;
432 ret[ 1 ][ 0 ] = -A[ 1 ][ 0 ] * detInv;
433 ret[ 0 ][ 1 ] = -A[ 0 ][ 1 ] * detInv;
434 return abs( det );
435 }
436 else
437 {
438 FieldMatrix< ctype, m , m > aat;
439 AAT_L( A, aat );
440 const ctype det = spdInvA( aat );
441 ATBT( A , aat , ret );
442 return det;
443 }
444 }
445
446 template< int m, int n >
447 static void xTRightInvA ( const FieldMatrix< ctype, m, n > &A, const FieldVector< ctype, n > &x, FieldVector< ctype, m > &y )
448 {
449 static_assert((n >= m), "Matrix has no right inverse.");
450 FieldMatrix< ctype, m, m > aat;
451 Ax( A, x, y );
452 AAT_L( A, aat );
453 spdInvAx( aat, y );
454 }
455 };
456
457 } // namespace Impl
458
459
460
466 template< class ct, int mydim, int cdim>
468 {
469 public:
470
472 typedef ct ctype;
473
475 static const int mydimension= mydim;
476
478 static const int coorddimension = cdim;
479
482
485
488
491
492 private:
495
497
498 // Helper class to compute a matrix pseudo inverse
499 typedef Impl::FieldMatrixHelper< ct > MatrixHelper;
500
501 public:
503 AffineGeometry ( const ReferenceElement &refElement, const GlobalCoordinate &origin,
504 const JacobianTransposed &jt )
505 : refElement_(refElement), origin_(origin), jacobianTransposed_(jt)
506 {
507 integrationElement_ = MatrixHelper::template rightInvA< mydimension, coorddimension >( jacobianTransposed_, jacobianInverseTransposed_ );
508 }
509
512 const JacobianTransposed &jt )
513 : AffineGeometry(ReferenceElements::general( gt ), origin, jt)
514 { }
515
517 template< class CoordVector >
518 AffineGeometry ( const ReferenceElement &refElement, const CoordVector &coordVector )
519 : refElement_(refElement), origin_(coordVector[0])
520 {
521 for( int i = 0; i < mydimension; ++i )
522 jacobianTransposed_[ i ] = coordVector[ i+1 ] - origin_;
523 integrationElement_ = MatrixHelper::template rightInvA< mydimension, coorddimension >( jacobianTransposed_, jacobianInverseTransposed_ );
524 }
525
527 template< class CoordVector >
528 AffineGeometry ( Dune::GeometryType gt, const CoordVector &coordVector )
529 : AffineGeometry(ReferenceElements::general( gt ), coordVector)
530 { }
531
533 bool affine () const { return true; }
534
536 Dune::GeometryType type () const { return refElement_.type(); }
537
539 int corners () const { return refElement_.size( mydimension ); }
540
542 GlobalCoordinate corner ( int i ) const
543 {
544 return global( refElement_.position( i, mydimension ) );
545 }
546
548 GlobalCoordinate center () const { return global( refElement_.position( 0, 0 ) ); }
549
557 {
558 GlobalCoordinate global( origin_ );
559 jacobianTransposed_.umtv( local, global );
560 return global;
561 }
562
577 {
578 LocalCoordinate local;
579 jacobianInverseTransposed_.mtv( global - origin_, local );
580 return local;
581 }
582
594 {
596 return integrationElement_;
597 }
598
600 ctype volume () const
601 {
602 return integrationElement_ * refElement_.volume();
603 }
604
612 {
614 return jacobianTransposed_;
615 }
616
624 {
626 return jacobianInverseTransposed_;
627 }
628
629 friend ReferenceElement referenceElement ( const AffineGeometry &geometry )
630 {
631 return geometry.refElement_;
632 }
633
634 private:
635 ReferenceElement refElement_;
636 GlobalCoordinate origin_;
637 JacobianTransposed jacobianTransposed_;
638 JacobianInverseTransposed jacobianInverseTransposed_;
639 ctype integrationElement_;
640 };
641
642} // namespace Dune
643
644#endif // #ifndef DUNE_GEOMETRY_AFFINEGEOMETRY_HH
Implementation of the Geometry interface for affine geometries.
Definition: affinegeometry.hh:468
AffineGeometry(const ReferenceElement &refElement, const CoordVector &coordVector)
Create affine geometry from reference element and a vector of vertex coordinates.
Definition: affinegeometry.hh:518
AffineGeometry(Dune::GeometryType gt, const GlobalCoordinate &origin, const JacobianTransposed &jt)
Create affine geometry from GeometryType, one vertex, and the Jacobian matrix.
Definition: affinegeometry.hh:511
FieldVector< ctype, mydimension > LocalCoordinate
Type for local coordinate vector.
Definition: affinegeometry.hh:481
Dune::GeometryType type() const
Obtain the type of the reference element.
Definition: affinegeometry.hh:536
static const int mydimension
Dimension of the geometry.
Definition: affinegeometry.hh:475
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:503
ctype volume() const
Obtain the volume of the element.
Definition: affinegeometry.hh:600
AffineGeometry(Dune::GeometryType gt, const CoordVector &coordVector)
Create affine geometry from GeometryType and a vector of vertex coordinates.
Definition: affinegeometry.hh:528
ctype integrationElement(const LocalCoordinate &local) const
Obtain the integration element.
Definition: affinegeometry.hh:593
const JacobianInverseTransposed & jacobianInverseTransposed(const LocalCoordinate &local) const
Obtain the transposed of the Jacobian's inverse.
Definition: affinegeometry.hh:623
FieldMatrix< ctype, mydimension, coorddimension > JacobianTransposed
Type for the transposed Jacobian matrix.
Definition: affinegeometry.hh:487
GlobalCoordinate corner(int i) const
Obtain coordinates of the i-th corner.
Definition: affinegeometry.hh:542
int corners() const
Obtain number of corners of the corresponding reference element.
Definition: affinegeometry.hh:539
FieldMatrix< ctype, coorddimension, mydimension > JacobianInverseTransposed
Type for the transposed inverse Jacobian matrix.
Definition: affinegeometry.hh:490
static const int coorddimension
Dimension of the world space.
Definition: affinegeometry.hh:478
GlobalCoordinate global(const LocalCoordinate &local) const
Evaluate the mapping.
Definition: affinegeometry.hh:556
GlobalCoordinate center() const
Obtain the centroid of the mapping's image.
Definition: affinegeometry.hh:548
ct ctype
Type used for coordinates.
Definition: affinegeometry.hh:472
FieldVector< ctype, coorddimension > GlobalCoordinate
Type for coordinate vector in world space.
Definition: affinegeometry.hh:484
bool affine() const
Always true: this is an affine geometry.
Definition: affinegeometry.hh:533
const JacobianTransposed & jacobianTransposed(const LocalCoordinate &local) const
Obtain the transposed of the Jacobian.
Definition: affinegeometry.hh:611
void mtv(const X &x, Y &y) const
y = A^T x
Definition: densematrix.hh:395
void umtv(const X &x, Y &y) const
y += A^T x
Definition: densematrix.hh:423
vector space out of a tensor product of fields.
Definition: fvector.hh:93
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:219
decltype(auto) type(int i, int c) const
obtain the type of subentity (i,c)
Definition: referenceelement.hh:149
int size(int c) const
number of subentities of codimension c
Definition: referenceelement.hh:95
decltype(auto) position(int i, int c) const
position of the barycenter of entity (i,c)
Definition: referenceelement.hh:181
Unique label for each type of entities that can occur in DUNE grids.
Definition: type.hh:277
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:147
unspecified-type ReferenceElement
Returns the type of reference element for the argument type T.
Definition: referenceelements.hh:495
Dune namespace.
Definition: alignedallocator.hh:10
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 (Jul 15, 22:36, 2024)