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