Dune Core Modules (2.7.1)

multilineargeometry.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_MULTILINEARGEOMETRY_HH
4 #define DUNE_GEOMETRY_MULTILINEARGEOMETRY_HH
5 
6 #include <cassert>
7 #include <functional>
8 #include <iterator>
9 #include <limits>
10 #include <vector>
11 
12 #include <dune/common/fmatrix.hh>
13 #include <dune/common/fvector.hh>
15 
17 #include <dune/geometry/referenceelements.hh>
18 #include <dune/geometry/type.hh>
19 
20 namespace Dune
21 {
22 
23  // MultiLinearGeometryTraits
24  // -------------------------
25 
35  template< class ct >
37  {
56  typedef Impl::FieldMatrixHelper< ct > MatrixHelper;
57 
59  static ct tolerance () { return ct( 16 ) * std::numeric_limits< ct >::epsilon(); }
60 
125  template< int mydim, int cdim >
127  {
128  typedef std::vector< FieldVector< ct, cdim > > Type;
129  };
130 
144  template< int dim >
146  {
147  static const bool v = false;
148  static const unsigned int topologyId = ~0u;
149  };
150  };
151 
152 
153 
154  // MultiLinearGeometry
155  // -------------------
156 
177  template< class ct, int mydim, int cdim, class Traits = MultiLinearGeometryTraits< ct > >
179  {
181 
182  public:
184  typedef ct ctype;
185 
187  static const int mydimension= mydim;
189  static const int coorddimension = cdim;
190 
196  typedef ctype Volume;
197 
200 
202  class JacobianInverseTransposed;
203 
204  protected:
205 
207 
208  public:
209 
212 
213  private:
214  static const bool hasSingleGeometryType = Traits::template hasSingleGeometryType< mydimension >::v;
215 
216  protected:
217  typedef typename Traits::MatrixHelper MatrixHelper;
218  typedef typename std::conditional< hasSingleGeometryType, std::integral_constant< unsigned int, Traits::template hasSingleGeometryType< mydimension >::topologyId >, unsigned int >::type TopologyId;
219 
220  public:
230  template< class Corners >
232  const Corners &corners )
233  : refElement_( refElement ),
234  corners_( corners )
235  {}
236 
246  template< class Corners >
248  const Corners &corners )
249  : refElement_( ReferenceElements::general( gt ) ),
250  corners_( corners )
251  {}
252 
254  bool affine () const
255  {
257  return affine( jt );
258  }
259 
261  Dune::GeometryType type () const { return GeometryType( toUnsignedInt(topologyId()), mydimension ); }
262 
264  int corners () const { return refElement().size( mydimension ); }
265 
267  GlobalCoordinate corner ( int i ) const
268  {
269  assert( (i >= 0) && (i < corners()) );
270  return std::cref(corners_).get()[ i ];
271  }
272 
274  GlobalCoordinate center () const { return global( refElement().position( 0, 0 ) ); }
275 
282  GlobalCoordinate global ( const LocalCoordinate &local ) const
283  {
284  using std::begin;
285 
286  auto cit = begin(std::cref(corners_).get());
288  global< false >( topologyId(), std::integral_constant< int, mydimension >(), cit, ctype( 1 ), local, ctype( 1 ), y );
289  return y;
290  }
291 
304  LocalCoordinate local ( const GlobalCoordinate &globalCoord ) const
305  {
306  const ctype tolerance = Traits::tolerance();
307  LocalCoordinate x = refElement().position( 0, 0 );
308  LocalCoordinate dx;
309  const bool affineMapping = this->affine();
310  do
311  {
312  // Newton's method: DF^n dx^n = F^n, x^{n+1} -= dx^n
313  const GlobalCoordinate dglobal = (*this).global( x ) - globalCoord;
314  const bool invertible =
315  MatrixHelper::template xTRightInvA< mydimension, coorddimension >( jacobianTransposed( x ), dglobal, dx );
316  if( ! invertible )
318 
319  // update x with correction
320  x -= dx;
321 
322  // for affine mappings only one iteration is needed
323  if ( affineMapping ) break;
324  } while( dx.two_norm2() > tolerance );
325  return x;
326  }
327 
342  ctype integrationElement ( const LocalCoordinate &local ) const
343  {
344  return MatrixHelper::template sqrtDetAAT< mydimension, coorddimension >( jacobianTransposed( local ) );
345  }
346 
355  Volume volume () const
356  {
357  return integrationElement( refElement().position( 0, 0 ) ) * refElement().volume();
358  }
359 
370  {
371  using std::begin;
372 
374  auto cit = begin(std::cref(corners_).get());
375  jacobianTransposed< false >( topologyId(), std::integral_constant< int, mydimension >(), cit, ctype( 1 ), local, ctype( 1 ), jt );
376  return jt;
377  }
378 
385  JacobianInverseTransposed jacobianInverseTransposed ( const LocalCoordinate &local ) const;
386 
387  friend ReferenceElement referenceElement ( const MultiLinearGeometry &geometry )
388  {
389  return geometry.refElement();
390  }
391 
392  protected:
393 
394  ReferenceElement refElement () const
395  {
396  return refElement_;
397  }
398 
399  TopologyId topologyId () const
400  {
401  return topologyId( std::integral_constant< bool, hasSingleGeometryType >() );
402  }
403 
404  template< bool add, int dim, class CornerIterator >
405  static void global ( TopologyId topologyId, std::integral_constant< int, dim >,
406  CornerIterator &cit, const ctype &df, const LocalCoordinate &x,
407  const ctype &rf, GlobalCoordinate &y );
408  template< bool add, class CornerIterator >
409  static void global ( TopologyId topologyId, std::integral_constant< int, 0 >,
410  CornerIterator &cit, const ctype &df, const LocalCoordinate &x,
411  const ctype &rf, GlobalCoordinate &y );
412 
413  template< bool add, int rows, int dim, class CornerIterator >
414  static void jacobianTransposed ( TopologyId topologyId, std::integral_constant< int, dim >,
415  CornerIterator &cit, const ctype &df, const LocalCoordinate &x,
416  const ctype &rf, FieldMatrix< ctype, rows, cdim > &jt );
417  template< bool add, int rows, class CornerIterator >
418  static void jacobianTransposed ( TopologyId topologyId, std::integral_constant< int, 0 >,
419  CornerIterator &cit, const ctype &df, const LocalCoordinate &x,
420  const ctype &rf, FieldMatrix< ctype, rows, cdim > &jt );
421 
422  template< int dim, class CornerIterator >
423  static bool affine ( TopologyId topologyId, std::integral_constant< int, dim >, CornerIterator &cit, JacobianTransposed &jt );
424  template< class CornerIterator >
425  static bool affine ( TopologyId topologyId, std::integral_constant< int, 0 >, CornerIterator &cit, JacobianTransposed &jt );
426 
427  bool affine ( JacobianTransposed &jacobianT ) const
428  {
429  using std::begin;
430 
431  auto cit = begin(std::cref(corners_).get());
432  return affine( topologyId(), std::integral_constant< int, mydimension >(), cit, jacobianT );
433  }
434 
435  private:
436  // The following methods are needed to convert the return type of topologyId to
437  // unsigned int with g++-4.4. It has problems casting integral_constant to the
438  // integral type.
439  static unsigned int toUnsignedInt(unsigned int i) { return i; }
440  template<unsigned int v>
441  static unsigned int toUnsignedInt(std::integral_constant<unsigned int,v> i) { return v; }
442  TopologyId topologyId ( std::integral_constant< bool, true > ) const { return TopologyId(); }
443  unsigned int topologyId ( std::integral_constant< bool, false > ) const { return refElement().type().id(); }
444 
445  ReferenceElement refElement_;
446  typename Traits::template CornerStorage< mydimension, coorddimension >::Type corners_;
447  };
448 
449 
450 
451  // MultiLinearGeometry::JacobianInverseTransposed
452  // ----------------------------------------------
453 
454  template< class ct, int mydim, int cdim, class Traits >
455  class MultiLinearGeometry< ct, mydim, cdim, Traits >::JacobianInverseTransposed
456  : public FieldMatrix< ctype, coorddimension, mydimension >
457  {
459 
460  public:
461  void setup ( const JacobianTransposed &jt )
462  {
463  detInv_ = MatrixHelper::template rightInvA< mydimension, coorddimension >( jt, static_cast< Base & >( *this ) );
464  }
465 
466  void setupDeterminant ( const JacobianTransposed &jt )
467  {
468  detInv_ = MatrixHelper::template sqrtDetAAT< mydimension, coorddimension >( jt );
469  }
470 
471  ctype det () const { return ctype( 1 ) / detInv_; }
472  ctype detInv () const { return detInv_; }
473 
474  private:
475  ctype detInv_;
476  };
477 
478 
479 
492  template< class ct, int mydim, int cdim, class Traits = MultiLinearGeometryTraits< ct > >
494  : public MultiLinearGeometry< ct, mydim, cdim, Traits >
495  {
498 
499  protected:
500  typedef typename Base::MatrixHelper MatrixHelper;
501 
502  public:
503  typedef typename Base::ReferenceElement ReferenceElement;
504 
505  typedef typename Base::ctype ctype;
506 
507  using Base::mydimension;
508  using Base::coorddimension;
509 
510  typedef typename Base::LocalCoordinate LocalCoordinate;
511  typedef typename Base::GlobalCoordinate GlobalCoordinate;
512  typedef typename Base::Volume Volume;
513 
515  typedef typename Base::JacobianInverseTransposed JacobianInverseTransposed;
516 
517  template< class CornerStorage >
518  CachedMultiLinearGeometry ( const ReferenceElement &referenceElement, const CornerStorage &cornerStorage )
519  : Base( referenceElement, cornerStorage ),
520  affine_( Base::affine( jacobianTransposed_ ) ),
521  jacobianInverseTransposedComputed_( false ),
522  integrationElementComputed_( false )
523  {}
524 
525  template< class CornerStorage >
526  CachedMultiLinearGeometry ( Dune::GeometryType gt, const CornerStorage &cornerStorage )
527  : Base( gt, cornerStorage ),
528  affine_( Base::affine( jacobianTransposed_ ) ),
529  jacobianInverseTransposedComputed_( false ),
530  integrationElementComputed_( false )
531  {}
532 
534  bool affine () const { return affine_; }
535 
536  using Base::corner;
537 
539  GlobalCoordinate center () const { return global( refElement().position( 0, 0 ) ); }
540 
547  GlobalCoordinate global ( const LocalCoordinate &local ) const
548  {
549  if( affine() )
550  {
552  jacobianTransposed_.umtv( local, global );
553  return global;
554  }
555  else
556  return Base::global( local );
557  }
558 
571  LocalCoordinate local ( const GlobalCoordinate &global ) const
572  {
573  if( affine() )
574  {
575  LocalCoordinate local;
576  if( jacobianInverseTransposedComputed_ )
577  jacobianInverseTransposed_.mtv( global - corner( 0 ), local );
578  else
579  MatrixHelper::template xTRightInvA< mydimension, coorddimension >( jacobianTransposed_, global - corner( 0 ), local );
580  return local;
581  }
582  else
583  return Base::local( global );
584  }
585 
600  ctype integrationElement ( const LocalCoordinate &local ) const
601  {
602  if( affine() )
603  {
604  if( !integrationElementComputed_ )
605  {
606  jacobianInverseTransposed_.setupDeterminant( jacobianTransposed_ );
607  integrationElementComputed_ = true;
608  }
609  return jacobianInverseTransposed_.detInv();
610  }
611  else
612  return Base::integrationElement( local );
613  }
614 
616  Volume volume () const
617  {
618  if( affine() )
619  return integrationElement( refElement().position( 0, 0 ) ) * refElement().volume();
620  else
621  return Base::volume();
622  }
623 
634  {
635  if( affine() )
636  return jacobianTransposed_;
637  else
638  return Base::jacobianTransposed( local );
639  }
640 
647  JacobianInverseTransposed jacobianInverseTransposed ( const LocalCoordinate &local ) const
648  {
649  if( affine() )
650  {
651  if( !jacobianInverseTransposedComputed_ )
652  {
653  jacobianInverseTransposed_.setup( jacobianTransposed_ );
654  jacobianInverseTransposedComputed_ = true;
655  integrationElementComputed_ = true;
656  }
657  return jacobianInverseTransposed_;
658  }
659  else
660  return Base::jacobianInverseTransposed( local );
661  }
662 
663  protected:
664  using Base::refElement;
665 
666  private:
667  mutable JacobianTransposed jacobianTransposed_;
668  mutable JacobianInverseTransposed jacobianInverseTransposed_;
669 
670  mutable bool affine_ : 1;
671 
672  mutable bool jacobianInverseTransposedComputed_ : 1;
673  mutable bool integrationElementComputed_ : 1;
674  };
675 
676 
677 
678  // Implementation of MultiLinearGeometry
679  // -------------------------------------
680 
681  template< class ct, int mydim, int cdim, class Traits >
682  inline typename MultiLinearGeometry< ct, mydim, cdim, Traits >::JacobianInverseTransposed
684  {
685  JacobianInverseTransposed jit;
686  jit.setup( jacobianTransposed( local ) );
687  return jit;
688  }
689 
690 
691  template< class ct, int mydim, int cdim, class Traits >
692  template< bool add, int dim, class CornerIterator >
694  ::global ( TopologyId topologyId, std::integral_constant< int, dim >,
695  CornerIterator &cit, const ctype &df, const LocalCoordinate &x,
696  const ctype &rf, GlobalCoordinate &y )
697  {
698  const ctype xn = df*x[ dim-1 ];
699  const ctype cxn = ctype( 1 ) - xn;
700 
701  if( Impl::isPrism( toUnsignedInt(topologyId), mydimension, mydimension-dim ) )
702  {
703  // apply (1-xn) times mapping for bottom
704  global< add >( topologyId, std::integral_constant< int, dim-1 >(), cit, df, x, rf*cxn, y );
705  // apply xn times mapping for top
706  global< true >( topologyId, std::integral_constant< int, dim-1 >(), cit, df, x, rf*xn, y );
707  }
708  else
709  {
710  assert( Impl::isPyramid( toUnsignedInt(topologyId), mydimension, mydimension-dim ) );
711  // apply (1-xn) times mapping for bottom (with argument x/(1-xn))
712  if( cxn > Traits::tolerance() || cxn < -Traits::tolerance() )
713  global< add >( topologyId, std::integral_constant< int, dim-1 >(), cit, df/cxn, x, rf*cxn, y );
714  else
715  global< add >( topologyId, std::integral_constant< int, dim-1 >(), cit, df, x, ctype( 0 ), y );
716  // apply xn times the tip
717  y.axpy( rf*xn, *cit );
718  ++cit;
719  }
720  }
721 
722  template< class ct, int mydim, int cdim, class Traits >
723  template< bool add, class CornerIterator >
725  ::global ( TopologyId topologyId, std::integral_constant< int, 0 >,
726  CornerIterator &cit, const ctype &df, const LocalCoordinate &x,
727  const ctype &rf, GlobalCoordinate &y )
728  {
729  const GlobalCoordinate &origin = *cit;
730  ++cit;
731  for( int i = 0; i < coorddimension; ++i )
732  y[ i ] = (add ? y[ i ] + rf*origin[ i ] : rf*origin[ i ]);
733  }
734 
735 
736  template< class ct, int mydim, int cdim, class Traits >
737  template< bool add, int rows, int dim, class CornerIterator >
739  ::jacobianTransposed ( TopologyId topologyId, std::integral_constant< int, dim >,
740  CornerIterator &cit, const ctype &df, const LocalCoordinate &x,
741  const ctype &rf, FieldMatrix< ctype, rows, cdim > &jt )
742  {
743  assert( rows >= dim );
744 
745  const ctype xn = df*x[ dim-1 ];
746  const ctype cxn = ctype( 1 ) - xn;
747 
748  auto cit2( cit );
749  if( Impl::isPrism( toUnsignedInt(topologyId), mydimension, mydimension-dim ) )
750  {
751  // apply (1-xn) times Jacobian for bottom
752  jacobianTransposed< add >( topologyId, std::integral_constant< int, dim-1 >(), cit2, df, x, rf*cxn, jt );
753  // apply xn times Jacobian for top
754  jacobianTransposed< true >( topologyId, std::integral_constant< int, dim-1 >(), cit2, df, x, rf*xn, jt );
755  // compute last row as difference between top value and bottom value
756  global< add >( topologyId, std::integral_constant< int, dim-1 >(), cit, df, x, -rf, jt[ dim-1 ] );
757  global< true >( topologyId, std::integral_constant< int, dim-1 >(), cit, df, x, rf, jt[ dim-1 ] );
758  }
759  else
760  {
761  assert( Impl::isPyramid( toUnsignedInt(topologyId), mydimension, mydimension-dim ) );
762  /*
763  * In the pyramid case, we need a transformation Tb: B -> R^n for the
764  * base B \subset R^{n-1}. The pyramid transformation is then defined as
765  * T: P \subset R^n -> R^n
766  * (x, xn) |-> (1-xn) Tb(x*) + xn t (x \in R^{n-1}, xn \in R)
767  * with the tip of the pyramid mapped to t and x* = x/(1-xn)
768  * the projection of (x,xn) onto the base.
769  *
770  * For the Jacobi matrix DT we get
771  * DT = ( A | b )
772  * with A = DTb(x*) (n x n-1 matrix)
773  * and b = dT/dxn (n-dim column vector).
774  * Furthermore
775  * b = -Tb(x*) + t + \sum_i dTb/dx_i(x^*) x_i/(1-xn)
776  *
777  * Note that both A and b are not defined in the pyramid tip (x=0, xn=1)!
778  * Indeed for B the unit square, Tb mapping B to the quadrilateral given
779  * by the vertices (0,0,0), (2,0,0), (0,1,0), (1,1,0) and t=(0,0,1), we get
780  *
781  * T(x,y,xn) = ( x(2-y/(1-xn)), y, xn )
782  * / 2-y/(1-xn) -x 0 \
783  * DT(x,y,xn) = | 0 1 0 |
784  * \ 0 0 1 /
785  * which is not continuous for xn -> 1, choose for example
786  * x=0, y=1-xn, xn -> 1 --> DT -> diag(1,1,1)
787  * x=1-xn, y=0, xn -> 1 --> DT -> diag(2,1,1)
788  *
789  * However, for Tb affine-linear, Tb(y) = My + y0, DTb = M:
790  * A = M
791  * b = -M x* - y0 + t + \sum_i M_i x_i/(1-xn)
792  * = -M x* - y0 + t + M x*
793  * = -y0 + t
794  * which is continuous for xn -> 1. Note that this b is also given by
795  * b = -Tb(0) + t + \sum_i dTb/dx_i(0) x_i/1
796  * that is replacing x* by 1 and 1-xn by 1 in the formular above.
797  *
798  * For xn -> 1, we can thus set x*=0, "1-xn"=1 (or anything != 0) and get
799  * the right result in case Tb is affine-linear.
800  */
801 
802  /* The second case effectively results in x* = 0 */
803  ctype dfcxn = (cxn > Traits::tolerance() || cxn < -Traits::tolerance()) ? ctype(df / cxn) : ctype(0);
804 
805  // initialize last row
806  // b = -Tb(x*)
807  // (b = -Tb(0) = -y0 in case xn -> 1 and Tb affine-linear)
808  global< add >( topologyId, std::integral_constant< int, dim-1 >(), cit, dfcxn, x, -rf, jt[ dim-1 ] );
809  // b += t
810  jt[ dim-1 ].axpy( rf, *cit );
811  ++cit;
812  // apply Jacobian for bottom (with argument x/(1-xn)) and correct last row
813  if( add )
814  {
815  FieldMatrix< ctype, dim-1, coorddimension > jt2;
816  // jt2 = dTb/dx_i(x*)
817  jacobianTransposed< false >( topologyId, std::integral_constant< int, dim-1 >(), cit2, dfcxn, x, rf, jt2 );
818  // A = dTb/dx_i(x*) (jt[j], j=0..dim-1)
819  // b += \sum_i dTb/dx_i(x*) x_i/(1-xn) (jt[dim-1])
820  // (b += 0 in case xn -> 1)
821  for( int j = 0; j < dim-1; ++j )
822  {
823  jt[ j ] += jt2[ j ];
824  jt[ dim-1 ].axpy( dfcxn*x[ j ], jt2[ j ] );
825  }
826  }
827  else
828  {
829  // jt = dTb/dx_i(x*)
830  jacobianTransposed< false >( topologyId, std::integral_constant< int, dim-1 >(), cit2, dfcxn, x, rf, jt );
831  // b += \sum_i dTb/dx_i(x*) x_i/(1-xn)
832  for( int j = 0; j < dim-1; ++j )
833  jt[ dim-1 ].axpy( dfcxn*x[ j ], jt[ j ] );
834  }
835  }
836  }
837 
838  template< class ct, int mydim, int cdim, class Traits >
839  template< bool add, int rows, class CornerIterator >
841  ::jacobianTransposed ( TopologyId topologyId, std::integral_constant< int, 0 >,
842  CornerIterator &cit, const ctype &df, const LocalCoordinate &x,
843  const ctype &rf, FieldMatrix< ctype, rows, cdim > &jt )
844  {
845  ++cit;
846  }
847 
848 
849 
850  template< class ct, int mydim, int cdim, class Traits >
851  template< int dim, class CornerIterator >
853  ::affine ( TopologyId topologyId, std::integral_constant< int, dim >, CornerIterator &cit, JacobianTransposed &jt )
854  {
855  const GlobalCoordinate &orgBottom = *cit;
856  if( !affine( topologyId, std::integral_constant< int, dim-1 >(), cit, jt ) )
857  return false;
858  const GlobalCoordinate &orgTop = *cit;
859 
860  if( Impl::isPrism( toUnsignedInt(topologyId), mydimension, mydimension-dim ) )
861  {
862  JacobianTransposed jtTop;
863  if( !affine( topologyId, std::integral_constant< int, dim-1 >(), cit, jtTop ) )
864  return false;
865 
866  // check whether both jacobians are identical
867  ctype norm( 0 );
868  for( int i = 0; i < dim-1; ++i )
869  norm += (jtTop[ i ] - jt[ i ]).two_norm2();
870  if( norm >= Traits::tolerance() )
871  return false;
872  }
873  else
874  ++cit;
875  jt[ dim-1 ] = orgTop - orgBottom;
876  return true;
877  }
878 
879  template< class ct, int mydim, int cdim, class Traits >
880  template< class CornerIterator >
882  ::affine ( TopologyId topologyId, std::integral_constant< int, 0 >, CornerIterator &cit, JacobianTransposed &jt )
883  {
884  ++cit;
885  return true;
886  }
887 
888 } // namespace Dune
889 
890 #endif // #ifndef DUNE_GEOMETRY_MULTILINEARGEOMETRY_HH
An implementation of the Geometry interface for affine geometries.
Implement a MultiLinearGeometry with additional caching.
Definition: multilineargeometry.hh:495
GlobalCoordinate global(const LocalCoordinate &local) const
evaluate the mapping
Definition: multilineargeometry.hh:547
bool affine() const
is this mapping affine?
Definition: multilineargeometry.hh:534
JacobianTransposed jacobianTransposed(const LocalCoordinate &local) const
obtain the transposed of the Jacobian
Definition: multilineargeometry.hh:633
GlobalCoordinate corner(int i) const
obtain coordinates of the i-th corner
Definition: multilineargeometry.hh:267
Volume volume() const
obtain the volume of the mapping's image
Definition: multilineargeometry.hh:616
ctype integrationElement(const LocalCoordinate &local) const
obtain the integration element
Definition: multilineargeometry.hh:600
GlobalCoordinate center() const
obtain the centroid of the mapping's image
Definition: multilineargeometry.hh:539
JacobianInverseTransposed jacobianInverseTransposed(const LocalCoordinate &local) const
obtain the transposed of the Jacobian's inverse
Definition: multilineargeometry.hh:647
void umtv(const X &x, Y &y) const
y += A^T x
Definition: densematrix.hh:446
FieldTraits< value_type >::real_type two_norm2() const
square of two norm (sum over squared values of entries), need for block recursion
Definition: densevector.hh:651
vector space out of a tensor product of fields.
Definition: fvector.hh:96
Unique label for each type of entities that can occur in DUNE grids.
Definition: type.hh:279
generic geometry implementation based on corner coordinates
Definition: multilineargeometry.hh:179
static const int mydimension
geometry dimension
Definition: multilineargeometry.hh:187
Dune::GeometryType type() const
obtain the name of the reference element
Definition: multilineargeometry.hh:261
FieldVector< ctype, coorddimension > GlobalCoordinate
type of global coordinates
Definition: multilineargeometry.hh:194
JacobianTransposed jacobianTransposed(const LocalCoordinate &local) const
obtain the transposed of the Jacobian
Definition: multilineargeometry.hh:369
GlobalCoordinate global(const LocalCoordinate &local) const
evaluate the mapping
Definition: multilineargeometry.hh:282
GlobalCoordinate center() const
obtain the centroid of the mapping's image
Definition: multilineargeometry.hh:274
GlobalCoordinate corner(int i) const
obtain coordinates of the i-th corner
Definition: multilineargeometry.hh:267
ct ctype
coordinate type
Definition: multilineargeometry.hh:184
static const int coorddimension
coordinate dimension
Definition: multilineargeometry.hh:189
int corners() const
obtain number of corners of the corresponding reference element
Definition: multilineargeometry.hh:264
Volume volume() const
obtain the volume of the mapping's image
Definition: multilineargeometry.hh:355
FieldVector< ctype, mydimension > LocalCoordinate
type of local coordinates
Definition: multilineargeometry.hh:192
MultiLinearGeometry(const ReferenceElement &refElement, const Corners &corners)
constructor
Definition: multilineargeometry.hh:231
ctype Volume
type of volume
Definition: multilineargeometry.hh:196
MultiLinearGeometry(Dune::GeometryType gt, const Corners &corners)
constructor
Definition: multilineargeometry.hh:247
ctype integrationElement(const LocalCoordinate &local) const
obtain the integration element
Definition: multilineargeometry.hh:342
FieldMatrix< ctype, mydimension, coorddimension > JacobianTransposed
type of jacobian transposed
Definition: multilineargeometry.hh:199
ReferenceElements::ReferenceElement ReferenceElement
type of reference element
Definition: multilineargeometry.hh:211
JacobianInverseTransposed jacobianInverseTransposed(const LocalCoordinate &local) const
obtain the transposed of the Jacobian's inverse
Definition: multilineargeometry.hh:683
bool affine() const
is this mapping affine?
Definition: multilineargeometry.hh:254
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:156
unspecified value type referenceElement(T &&... t)
Returns a reference element for the objects t....
unspecified-type ReferenceElement
Returns the type of reference element for the argument type T.
Definition: referenceelements.hh:495
auto max(ADLTag< 0 >, const V &v1, const V &v2)
implements binary Simd::max()
Definition: defaults.hh:79
Dune namespace.
Definition: alignedallocator.hh:14
Class providing access to the singletons of the reference elements.
Definition: referenceelements.hh:168
typename Container::ReferenceElement ReferenceElement
The reference element type.
Definition: referenceelements.hh:186
template specifying the storage for the corners
Definition: multilineargeometry.hh:127
will there be only one geometry type for a dimension?
Definition: multilineargeometry.hh:146
default traits class for MultiLinearGeometry
Definition: multilineargeometry.hh:37
Impl::FieldMatrixHelper< ct > MatrixHelper
helper structure containing some matrix routines
Definition: multilineargeometry.hh:56
static ct tolerance()
tolerance to numerical algorithms
Definition: multilineargeometry.hh:59
A unique label for each type of element that can occur in a grid.
Traits for type conversions and type information.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.80.0 (May 16, 22:29, 2024)