Dune Core Modules (2.6.0)

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
15
17#include <dune/geometry/referenceelements.hh>
18#include <dune/geometry/type.hh>
19
20namespace 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
195
198
200 class JacobianInverseTransposed;
201
202 protected:
203
205
206 public:
207
210
211 private:
212 static const bool hasSingleGeometryType = Traits::template hasSingleGeometryType< mydimension >::v;
213
214 protected:
215 typedef typename Traits::MatrixHelper MatrixHelper;
216 typedef typename std::conditional< hasSingleGeometryType, std::integral_constant< unsigned int, Traits::template hasSingleGeometryType< mydimension >::topologyId >, unsigned int >::type TopologyId;
217
218 public:
228 template< class Corners >
230 const Corners &corners )
231 : refElement_( refElement ),
232 corners_( corners )
233 {}
234
244 template< class Corners >
246 const Corners &corners )
247 : refElement_( ReferenceElements::general( gt ) ),
248 corners_( corners )
249 {}
250
252 bool affine () const
253 {
255 return affine( jt );
256 }
257
259 Dune::GeometryType type () const { return GeometryType( toUnsignedInt(topologyId()), mydimension ); }
260
262 int corners () const { return refElement().size( mydimension ); }
263
265 GlobalCoordinate corner ( int i ) const
266 {
267 assert( (i >= 0) && (i < corners()) );
268 return std::cref(corners_).get()[ i ];
269 }
270
272 GlobalCoordinate center () const { return global( refElement().position( 0, 0 ) ); }
273
281 {
282 using std::begin;
283
284 auto cit = begin(std::cref(corners_).get());
286 global< false >( topologyId(), std::integral_constant< int, mydimension >(), cit, ctype( 1 ), local, ctype( 1 ), y );
287 return y;
288 }
289
302 LocalCoordinate local ( const GlobalCoordinate &globalCoord ) const
303 {
304 const ctype tolerance = Traits::tolerance();
305 LocalCoordinate x = refElement().position( 0, 0 );
307 do
308 {
309 // Newton's method: DF^n dx^n = F^n, x^{n+1} -= dx^n
310 const GlobalCoordinate dglobal = (*this).global( x ) - globalCoord;
311 MatrixHelper::template xTRightInvA< mydimension, coorddimension >( jacobianTransposed( x ), dglobal, dx );
312 x -= dx;
313 } while( dx.two_norm2() > tolerance );
314 return x;
315 }
316
332 {
333 return MatrixHelper::template sqrtDetAAT< mydimension, coorddimension >( jacobianTransposed( local ) );
334 }
335
344 ctype volume () const
345 {
346 return integrationElement( refElement().position( 0, 0 ) ) * refElement().volume();
347 }
348
359 {
360 using std::begin;
361
363 auto cit = begin(std::cref(corners_).get());
364 jacobianTransposed< false >( topologyId(), std::integral_constant< int, mydimension >(), cit, ctype( 1 ), local, ctype( 1 ), jt );
365 return jt;
366 }
367
374 JacobianInverseTransposed jacobianInverseTransposed ( const LocalCoordinate &local ) const;
375
376 friend ReferenceElement referenceElement ( const MultiLinearGeometry &geometry )
377 {
378 return geometry.refElement();
379 }
380
381 protected:
382
383 ReferenceElement refElement () const
384 {
385 return refElement_;
386 }
387
388 TopologyId topologyId () const
389 {
390 return topologyId( std::integral_constant< bool, hasSingleGeometryType >() );
391 }
392
393 template< bool add, int dim, class CornerIterator >
394 static void global ( TopologyId topologyId, std::integral_constant< int, dim >,
395 CornerIterator &cit, const ctype &df, const LocalCoordinate &x,
396 const ctype &rf, GlobalCoordinate &y );
397 template< bool add, class CornerIterator >
398 static void global ( TopologyId topologyId, std::integral_constant< int, 0 >,
399 CornerIterator &cit, const ctype &df, const LocalCoordinate &x,
400 const ctype &rf, GlobalCoordinate &y );
401
402 template< bool add, int rows, int dim, class CornerIterator >
403 static void jacobianTransposed ( TopologyId topologyId, std::integral_constant< int, dim >,
404 CornerIterator &cit, const ctype &df, const LocalCoordinate &x,
405 const ctype &rf, FieldMatrix< ctype, rows, cdim > &jt );
406 template< bool add, int rows, class CornerIterator >
407 static void jacobianTransposed ( TopologyId topologyId, std::integral_constant< int, 0 >,
408 CornerIterator &cit, const ctype &df, const LocalCoordinate &x,
409 const ctype &rf, FieldMatrix< ctype, rows, cdim > &jt );
410
411 template< int dim, class CornerIterator >
412 static bool affine ( TopologyId topologyId, std::integral_constant< int, dim >, CornerIterator &cit, JacobianTransposed &jt );
413 template< class CornerIterator >
414 static bool affine ( TopologyId topologyId, std::integral_constant< int, 0 >, CornerIterator &cit, JacobianTransposed &jt );
415
416 bool affine ( JacobianTransposed &jacobianT ) const
417 {
418 using std::begin;
419
420 auto cit = begin(std::cref(corners_).get());
421 return affine( topologyId(), std::integral_constant< int, mydimension >(), cit, jacobianT );
422 }
423
424 private:
425 // The following methods are needed to convert the return type of topologyId to
426 // unsigned int with g++-4.4. It has problems casting integral_constant to the
427 // integral type.
428 static unsigned int toUnsignedInt(unsigned int i) { return i; }
429 template<unsigned int v>
430 static unsigned int toUnsignedInt(std::integral_constant<unsigned int,v> i) { return v; }
431 TopologyId topologyId ( std::integral_constant< bool, true > ) const { return TopologyId(); }
432 unsigned int topologyId ( std::integral_constant< bool, false > ) const { return refElement().type().id(); }
433
434 ReferenceElement refElement_;
435 typename Traits::template CornerStorage< mydimension, coorddimension >::Type corners_;
436 };
437
438
439
440 // MultiLinearGeometry::JacobianInverseTransposed
441 // ----------------------------------------------
442
443 template< class ct, int mydim, int cdim, class Traits >
444 class MultiLinearGeometry< ct, mydim, cdim, Traits >::JacobianInverseTransposed
445 : public FieldMatrix< ctype, coorddimension, mydimension >
446 {
448
449 public:
450 void setup ( const JacobianTransposed &jt )
451 {
452 detInv_ = MatrixHelper::template rightInvA< mydimension, coorddimension >( jt, static_cast< Base & >( *this ) );
453 }
454
455 void setupDeterminant ( const JacobianTransposed &jt )
456 {
457 detInv_ = MatrixHelper::template sqrtDetAAT< mydimension, coorddimension >( jt );
458 }
459
460 ctype det () const { return ctype( 1 ) / detInv_; }
461 ctype detInv () const { return detInv_; }
462
463 private:
464 ctype detInv_;
465 };
466
467
468
481 template< class ct, int mydim, int cdim, class Traits = MultiLinearGeometryTraits< ct > >
483 : public MultiLinearGeometry< ct, mydim, cdim, Traits >
484 {
487
488 protected:
489 typedef typename Base::MatrixHelper MatrixHelper;
490
491 public:
492 typedef typename Base::ReferenceElement ReferenceElement;
493
494 typedef typename Base::ctype ctype;
495
496 using Base::mydimension;
498
499 typedef typename Base::LocalCoordinate LocalCoordinate;
501
503 typedef typename Base::JacobianInverseTransposed JacobianInverseTransposed;
504
505 template< class CornerStorage >
506 CachedMultiLinearGeometry ( const ReferenceElement &referenceElement, const CornerStorage &cornerStorage )
507 : Base( referenceElement, cornerStorage ),
508 affine_( Base::affine( jacobianTransposed_ ) ),
509 jacobianInverseTransposedComputed_( false ),
510 integrationElementComputed_( false )
511 {}
512
513 template< class CornerStorage >
514 CachedMultiLinearGeometry ( Dune::GeometryType gt, const CornerStorage &cornerStorage )
515 : Base( gt, cornerStorage ),
516 affine_( Base::affine( jacobianTransposed_ ) ),
517 jacobianInverseTransposedComputed_( false ),
518 integrationElementComputed_( false )
519 {}
520
522 bool affine () const { return affine_; }
523
524 using Base::corner;
525
527 GlobalCoordinate center () const { return global( refElement().position( 0, 0 ) ); }
528
536 {
537 if( affine() )
538 {
540 jacobianTransposed_.umtv( local, global );
541 return global;
542 }
543 else
544 return Base::global( local );
545 }
546
560 {
561 if( affine() )
562 {
563 LocalCoordinate local;
564 if( jacobianInverseTransposedComputed_ )
565 jacobianInverseTransposed_.mtv( global - corner( 0 ), local );
566 else
567 MatrixHelper::template xTRightInvA< mydimension, coorddimension >( jacobianTransposed_, global - corner( 0 ), local );
568 return local;
569 }
570 else
571 return Base::local( global );
572 }
573
588 ctype integrationElement ( const LocalCoordinate &local ) const
589 {
590 if( affine() )
591 {
592 if( !integrationElementComputed_ )
593 {
594 jacobianInverseTransposed_.setupDeterminant( jacobianTransposed_ );
595 integrationElementComputed_ = true;
596 }
597 return jacobianInverseTransposed_.detInv();
598 }
599 else
600 return Base::integrationElement( local );
601 }
602
604 ctype volume () const
605 {
606 if( affine() )
607 return integrationElement( refElement().position( 0, 0 ) ) * refElement().volume();
608 else
609 return Base::volume();
610 }
611
622 {
623 if( affine() )
624 return jacobianTransposed_;
625 else
626 return Base::jacobianTransposed( local );
627 }
628
635 JacobianInverseTransposed jacobianInverseTransposed ( const LocalCoordinate &local ) const
636 {
637 if( affine() )
638 {
639 if( !jacobianInverseTransposedComputed_ )
640 {
641 jacobianInverseTransposed_.setup( jacobianTransposed_ );
642 jacobianInverseTransposedComputed_ = true;
643 integrationElementComputed_ = true;
644 }
645 return jacobianInverseTransposed_;
646 }
647 else
648 return Base::jacobianInverseTransposed( local );
649 }
650
651 protected:
652 using Base::refElement;
653
654 private:
655 mutable JacobianTransposed jacobianTransposed_;
656 mutable JacobianInverseTransposed jacobianInverseTransposed_;
657
658 mutable bool affine_ : 1;
659
660 mutable bool jacobianInverseTransposedComputed_ : 1;
661 mutable bool integrationElementComputed_ : 1;
662 };
663
664
665
666 // Implementation of MultiLinearGeometry
667 // -------------------------------------
668
669 template< class ct, int mydim, int cdim, class Traits >
670 inline typename MultiLinearGeometry< ct, mydim, cdim, Traits >::JacobianInverseTransposed
672 {
673 JacobianInverseTransposed jit;
674 jit.setup( jacobianTransposed( local ) );
675 return jit;
676 }
677
678
679 template< class ct, int mydim, int cdim, class Traits >
680 template< bool add, int dim, class CornerIterator >
682 ::global ( TopologyId topologyId, std::integral_constant< int, dim >,
683 CornerIterator &cit, const ctype &df, const LocalCoordinate &x,
684 const ctype &rf, GlobalCoordinate &y )
685 {
686 const ctype xn = df*x[ dim-1 ];
687 const ctype cxn = ctype( 1 ) - xn;
688
689 if( Impl::isPrism( toUnsignedInt(topologyId), mydimension, mydimension-dim ) )
690 {
691 // apply (1-xn) times mapping for bottom
692 global< add >( topologyId, std::integral_constant< int, dim-1 >(), cit, df, x, rf*cxn, y );
693 // apply xn times mapping for top
694 global< true >( topologyId, std::integral_constant< int, dim-1 >(), cit, df, x, rf*xn, y );
695 }
696 else
697 {
698 assert( Impl::isPyramid( toUnsignedInt(topologyId), mydimension, mydimension-dim ) );
699 // apply (1-xn) times mapping for bottom (with argument x/(1-xn))
700 if( cxn > Traits::tolerance() || cxn < -Traits::tolerance() )
701 global< add >( topologyId, std::integral_constant< int, dim-1 >(), cit, df/cxn, x, rf*cxn, y );
702 else
703 global< add >( topologyId, std::integral_constant< int, dim-1 >(), cit, df, x, ctype( 0 ), y );
704 // apply xn times the tip
705 y.axpy( rf*xn, *cit );
706 ++cit;
707 }
708 }
709
710 template< class ct, int mydim, int cdim, class Traits >
711 template< bool add, class CornerIterator >
713 ::global ( TopologyId topologyId, std::integral_constant< int, 0 >,
714 CornerIterator &cit, const ctype &df, const LocalCoordinate &x,
715 const ctype &rf, GlobalCoordinate &y )
716 {
717 const GlobalCoordinate &origin = *cit;
718 ++cit;
719 for( int i = 0; i < coorddimension; ++i )
720 y[ i ] = (add ? y[ i ] + rf*origin[ i ] : rf*origin[ i ]);
721 }
722
723
724 template< class ct, int mydim, int cdim, class Traits >
725 template< bool add, int rows, int dim, class CornerIterator >
727 ::jacobianTransposed ( TopologyId topologyId, std::integral_constant< int, dim >,
728 CornerIterator &cit, const ctype &df, const LocalCoordinate &x,
729 const ctype &rf, FieldMatrix< ctype, rows, cdim > &jt )
730 {
731 assert( rows >= dim );
732
733 const ctype xn = df*x[ dim-1 ];
734 const ctype cxn = ctype( 1 ) - xn;
735
736 auto cit2( cit );
737 if( Impl::isPrism( toUnsignedInt(topologyId), mydimension, mydimension-dim ) )
738 {
739 // apply (1-xn) times Jacobian for bottom
740 jacobianTransposed< add >( topologyId, std::integral_constant< int, dim-1 >(), cit2, df, x, rf*cxn, jt );
741 // apply xn times Jacobian for top
742 jacobianTransposed< true >( topologyId, std::integral_constant< int, dim-1 >(), cit2, df, x, rf*xn, jt );
743 // compute last row as difference between top value and bottom value
744 global< add >( topologyId, std::integral_constant< int, dim-1 >(), cit, df, x, -rf, jt[ dim-1 ] );
745 global< true >( topologyId, std::integral_constant< int, dim-1 >(), cit, df, x, rf, jt[ dim-1 ] );
746 }
747 else
748 {
749 assert( Impl::isPyramid( toUnsignedInt(topologyId), mydimension, mydimension-dim ) );
750 /*
751 * In the pyramid case, we need a transformation Tb: B -> R^n for the
752 * base B \subset R^{n-1}. The pyramid transformation is then defined as
753 * T: P \subset R^n -> R^n
754 * (x, xn) |-> (1-xn) Tb(x*) + xn t (x \in R^{n-1}, xn \in R)
755 * with the tip of the pyramid mapped to t and x* = x/(1-xn)
756 * the projection of (x,xn) onto the base.
757 *
758 * For the Jacobi matrix DT we get
759 * DT = ( A | b )
760 * with A = DTb(x*) (n x n-1 matrix)
761 * and b = dT/dxn (n-dim column vector).
762 * Furthermore
763 * b = -Tb(x*) + t + \sum_i dTb/dx_i(x^*) x_i/(1-xn)
764 *
765 * Note that both A and b are not defined in the pyramid tip (x=0, xn=1)!
766 * Indeed for B the unit square, Tb mapping B to the quadrilateral given
767 * by the vertices (0,0,0), (2,0,0), (0,1,0), (1,1,0) and t=(0,0,1), we get
768 *
769 * T(x,y,xn) = ( x(2-y/(1-xn)), y, xn )
770 * / 2-y/(1-xn) -x 0 \
771 * DT(x,y,xn) = | 0 1 0 |
772 * \ 0 0 1 /
773 * which is not continuous for xn -> 1, choose for example
774 * x=0, y=1-xn, xn -> 1 --> DT -> diag(1,1,1)
775 * x=1-xn, y=0, xn -> 1 --> DT -> diag(2,1,1)
776 *
777 * However, for Tb affine-linear, Tb(y) = My + y0, DTb = M:
778 * A = M
779 * b = -M x* - y0 + t + \sum_i M_i x_i/(1-xn)
780 * = -M x* - y0 + t + M x*
781 * = -y0 + t
782 * which is continuous for xn -> 1. Note that this b is also given by
783 * b = -Tb(0) + t + \sum_i dTb/dx_i(0) x_i/1
784 * that is replacing x* by 1 and 1-xn by 1 in the formular above.
785 *
786 * For xn -> 1, we can thus set x*=0, "1-xn"=1 (or anything != 0) and get
787 * the right result in case Tb is affine-linear.
788 */
789
790 /* The second case effectively results in x* = 0 */
791 ctype dfcxn = (cxn > Traits::tolerance() || cxn < -Traits::tolerance()) ? ctype(df / cxn) : ctype(0);
792
793 // initialize last row
794 // b = -Tb(x*)
795 // (b = -Tb(0) = -y0 in case xn -> 1 and Tb affine-linear)
796 global< add >( topologyId, std::integral_constant< int, dim-1 >(), cit, dfcxn, x, -rf, jt[ dim-1 ] );
797 // b += t
798 jt[ dim-1 ].axpy( rf, *cit );
799 ++cit;
800 // apply Jacobian for bottom (with argument x/(1-xn)) and correct last row
801 if( add )
802 {
803 FieldMatrix< ctype, dim-1, coorddimension > jt2;
804 // jt2 = dTb/dx_i(x*)
805 jacobianTransposed< false >( topologyId, std::integral_constant< int, dim-1 >(), cit2, dfcxn, x, rf, jt2 );
806 // A = dTb/dx_i(x*) (jt[j], j=0..dim-1)
807 // b += \sum_i dTb/dx_i(x*) x_i/(1-xn) (jt[dim-1])
808 // (b += 0 in case xn -> 1)
809 for( int j = 0; j < dim-1; ++j )
810 {
811 jt[ j ] += jt2[ j ];
812 jt[ dim-1 ].axpy( dfcxn*x[ j ], jt2[ j ] );
813 }
814 }
815 else
816 {
817 // jt = dTb/dx_i(x*)
818 jacobianTransposed< false >( topologyId, std::integral_constant< int, dim-1 >(), cit2, dfcxn, x, rf, jt );
819 // b += \sum_i dTb/dx_i(x*) x_i/(1-xn)
820 for( int j = 0; j < dim-1; ++j )
821 jt[ dim-1 ].axpy( dfcxn*x[ j ], jt[ j ] );
822 }
823 }
824 }
825
826 template< class ct, int mydim, int cdim, class Traits >
827 template< bool add, int rows, class CornerIterator >
829 ::jacobianTransposed ( TopologyId topologyId, std::integral_constant< int, 0 >,
830 CornerIterator &cit, const ctype &df, const LocalCoordinate &x,
831 const ctype &rf, FieldMatrix< ctype, rows, cdim > &jt )
832 {
833 ++cit;
834 }
835
836
837
838 template< class ct, int mydim, int cdim, class Traits >
839 template< int dim, class CornerIterator >
841 ::affine ( TopologyId topologyId, std::integral_constant< int, dim >, CornerIterator &cit, JacobianTransposed &jt )
842 {
843 const GlobalCoordinate &orgBottom = *cit;
844 if( !affine( topologyId, std::integral_constant< int, dim-1 >(), cit, jt ) )
845 return false;
846 const GlobalCoordinate &orgTop = *cit;
847
848 if( Impl::isPrism( toUnsignedInt(topologyId), mydimension, mydimension-dim ) )
849 {
850 JacobianTransposed jtTop;
851 if( !affine( topologyId, std::integral_constant< int, dim-1 >(), cit, jtTop ) )
852 return false;
853
854 // check whether both jacobians are identical
855 ctype norm( 0 );
856 for( int i = 0; i < dim-1; ++i )
857 norm += (jtTop[ i ] - jt[ i ]).two_norm2();
858 if( norm >= Traits::tolerance() )
859 return false;
860 }
861 else
862 ++cit;
863 jt[ dim-1 ] = orgTop - orgBottom;
864 return true;
865 }
866
867 template< class ct, int mydim, int cdim, class Traits >
868 template< class CornerIterator >
870 ::affine ( TopologyId topologyId, std::integral_constant< int, 0 >, CornerIterator &cit, JacobianTransposed &jt )
871 {
872 ++cit;
873 return true;
874 }
875
876} // namespace Dune
877
878#endif // #ifndef DUNE_GEOMETRY_MULTILINEARGEOMETRY_HH
An implementation of the Geometry interface for affine geometries.
Implement a MultiLinearGeometry with additional caching.
Definition: multilineargeometry.hh:484
GlobalCoordinate global(const LocalCoordinate &local) const
evaluate the mapping
Definition: multilineargeometry.hh:535
bool affine() const
is this mapping affine?
Definition: multilineargeometry.hh:522
JacobianTransposed jacobianTransposed(const LocalCoordinate &local) const
obtain the transposed of the Jacobian
Definition: multilineargeometry.hh:621
GlobalCoordinate corner(int i) const
obtain coordinates of the i-th corner
Definition: multilineargeometry.hh:265
ctype volume() const
obtain the volume of the mapping's image
Definition: multilineargeometry.hh:604
ctype integrationElement(const LocalCoordinate &local) const
obtain the integration element
Definition: multilineargeometry.hh:588
GlobalCoordinate center() const
obtain the centroid of the mapping's image
Definition: multilineargeometry.hh:527
JacobianInverseTransposed jacobianInverseTransposed(const LocalCoordinate &local) const
obtain the transposed of the Jacobian's inverse
Definition: multilineargeometry.hh:635
void umtv(const X &x, Y &y) const
y += A^T x
Definition: densematrix.hh:423
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:599
vector space out of a tensor product of fields.
Definition: fvector.hh:93
Unique label for each type of entities that can occur in DUNE grids.
Definition: type.hh:277
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:259
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:358
GlobalCoordinate global(const LocalCoordinate &local) const
evaluate the mapping
Definition: multilineargeometry.hh:280
GlobalCoordinate center() const
obtain the centroid of the mapping's image
Definition: multilineargeometry.hh:272
GlobalCoordinate corner(int i) const
obtain coordinates of the i-th corner
Definition: multilineargeometry.hh:265
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:262
FieldVector< ctype, mydimension > LocalCoordinate
type of local coordinates
Definition: multilineargeometry.hh:192
ctype volume() const
obtain the volume of the mapping's image
Definition: multilineargeometry.hh:344
MultiLinearGeometry(const ReferenceElement &refElement, const Corners &corners)
constructor
Definition: multilineargeometry.hh:229
MultiLinearGeometry(Dune::GeometryType gt, const Corners &corners)
constructor
Definition: multilineargeometry.hh:245
ctype integrationElement(const LocalCoordinate &local) const
obtain the integration element
Definition: multilineargeometry.hh:331
FieldMatrix< ctype, mydimension, coorddimension > JacobianTransposed
type of jacobian transposed
Definition: multilineargeometry.hh:197
ReferenceElements::ReferenceElement ReferenceElement
type of reference element
Definition: multilineargeometry.hh:209
JacobianInverseTransposed jacobianInverseTransposed(const LocalCoordinate &local) const
obtain the transposed of the Jacobian's inverse
Definition: multilineargeometry.hh:671
bool affine() const
is this mapping affine?
Definition: multilineargeometry.hh:252
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:147
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
Dune namespace.
Definition: alignedallocator.hh:10
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.111.3 (Nov 24, 23:30, 2024)