dune-geometry  2.3beta2
multilineargeometry.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_MULTILINEARGEOMETRY_HH
4 #define DUNE_GEOMETRY_MULTILINEARGEOMETRY_HH
5 
6 #include <cassert>
7 #include <limits>
8 #include <vector>
9 
10 #include <dune/common/fmatrix.hh>
11 #include <dune/common/fvector.hh>
12 #include <dune/common/typetraits.hh>
13 
15 #include <dune/geometry/type.hh>
18 
19 namespace Dune
20 {
21 
22  // External Forward Declarations
23  // -----------------------------
24 
25  template< class ctype, int dim >
26  class ReferenceElement;
27 
28  template< class ctype, int dim >
29  struct ReferenceElements;
30 
31 
32 
33  // MultiLinearGeometryTraits
34  // -------------------------
35 
45  template< class ct >
47  {
67 
69  static ct tolerance () { return ct( 16 ) * std::numeric_limits< ct >::epsilon(); }
70 
94  template< int mydim, int cdim >
96  {
97  typedef std::vector< FieldVector< ct, cdim > > Type;
98  };
99 
113  template< int dim >
115  {
116  static const bool v = false;
117  static const unsigned int topologyId = ~0u;
118  };
119  };
120 
121 
122 
123  // MultiLinearGeometry
124  // -------------------
125 
149  template< class ct, int mydim, int cdim, class Traits = MultiLinearGeometryTraits< ct > >
151  {
153 
154  public:
156  typedef ct ctype;
157 
159  static const int mydimension= mydim;
161  static const int coorddimension = cdim;
162 
167 
168  typedef FieldVector< ctype, mydimension > LocalCoordinate;
170  typedef FieldVector< ctype, coorddimension > GlobalCoordinate;
171 
173  typedef FieldMatrix< ctype, mydimension, coorddimension > JacobianTransposed;
174 
177 
182 
185 
186  private:
187  static const bool hasSingleGeometryType = Traits::template hasSingleGeometryType< mydimension >::v;
188 
189  protected:
190  typedef typename Traits::MatrixHelper MatrixHelper;
191  typedef typename conditional< hasSingleGeometryType, integral_constant< unsigned int, Traits::template hasSingleGeometryType< mydimension >::topologyId >, unsigned int >::type TopologyId;
192 
194 
195  private:
196  typedef typename Traits::template CornerStorage< mydimension, coorddimension >::Type::const_iterator CornerIterator;
197 
198  public:
208  template< class Corners >
210  const Corners &corners )
211  : refElement_( &refElement ),
212  corners_( corners )
213  {}
214 
224  template< class Corners >
226  const Corners &corners )
227  : refElement_( &ReferenceElements::general( gt ) ),
228  corners_( corners )
229  {}
230 
232  bool affine () const
233  {
234  CornerIterator cit = corners_.begin();
235  return affine( topologyId(), integral_constant< int, mydimension >(), cit, jacobianTransposed_ );
236  }
237 
240 
242  int corners () const { return refElement().size( mydimension ); }
243 
245  GlobalCoordinate corner ( int i ) const
246  {
247  assert( (i >= 0) && (i < corners()) );
248  return corners_[ i ];
249  }
250 
252  GlobalCoordinate center () const { return global( refElement().position( 0, 0 ) ); }
253 
261  {
262  CornerIterator cit = corners_.begin();
264  global< false >( topologyId(), integral_constant< int, mydimension >(), cit, ctype( 1 ), local, ctype( 1 ), y );
265  return y;
266  }
267 
280  {
281  const ctype tolerance = Traits::tolerance();
282  LocalCoordinate x = refElement().position( 0, 0 );
283  LocalCoordinate dx;
284  do
285  {
286  // Newton's method: DF^n dx^n = F^n, x^{n+1} -= dx^n
287  const GlobalCoordinate dglobal = (*this).global( x ) - global;
288  MatrixHelper::template xTRightInvA< mydimension, coorddimension >( jacobianTransposed( x ), dglobal, dx );
289  x -= dx;
290  assert( refElement().checkInside( x ) );
291  } while( dx.two_norm2() > tolerance );
292  return x;
293  }
294 
310  {
311  return MatrixHelper::template sqrtDetAAT< mydimension, coorddimension >( jacobianTransposed( local ) );
312  }
313 
322  ctype volume () const
323  {
324  return integrationElement( refElement().position( 0, 0 ) ) * refElement().volume();
325  }
326 
337  {
338  CornerIterator cit = corners_.begin();
339  jacobianTransposed< false >( topologyId(), integral_constant< int, mydimension >(), cit, ctype( 1 ), local, ctype( 1 ), jacobianTransposed_ );
340  return jacobianTransposed_;
341  }
342 
349  const JacobianInverseTransposed &jacobianInverseTransposed ( const LocalCoordinate &local ) const;
350 
351  protected:
352  const ReferenceElement &refElement () const { return *refElement_; }
353 
355  {
356  return topologyId( integral_constant< bool, hasSingleGeometryType >() );
357  }
358 
359  template< bool add, int dim >
360  static void global ( TopologyId topologyId, integral_constant< int, dim >,
361  CornerIterator &cit, const ctype &df, const LocalCoordinate &x,
362  const ctype &rf, GlobalCoordinate &y );
363  template< bool add >
364  static void global ( TopologyId topologyId, integral_constant< int, 0 >,
365  CornerIterator &cit, const ctype &df, const LocalCoordinate &x,
366  const ctype &rf, GlobalCoordinate &y );
367 
368  template< bool add, int rows, int dim >
369  static void jacobianTransposed ( TopologyId topologyId, integral_constant< int, dim >,
370  CornerIterator &cit, const ctype &df, const LocalCoordinate &x,
371  const ctype &rf, FieldMatrix< ctype, rows, cdim > &jt );
372  template< bool add, int rows >
373  static void jacobianTransposed ( TopologyId topologyId, integral_constant< int, 0 >,
374  CornerIterator &cit, const ctype &df, const LocalCoordinate &x,
375  const ctype &rf, FieldMatrix< ctype, rows, cdim > &jt );
376 
377  template< int dim >
378  static bool affine ( TopologyId topologyId, integral_constant< int, dim >, CornerIterator &cit, JacobianTransposed &jt );
379  static bool affine ( TopologyId topologyId, integral_constant< int, 0 >, CornerIterator &cit, JacobianTransposed &jt );
380 
381  protected:
382  TopologyId topologyId ( integral_constant< bool, true > ) const { return TopologyId(); }
383  unsigned int topologyId ( integral_constant< bool, false > ) const { return refElement().type().id(); }
384 
387 
388  private:
389  const ReferenceElement *refElement_;
390  typename Traits::template CornerStorage< mydimension, coorddimension >::Type corners_;
391  };
392 
393 
394 
395  // MultiLinearGeometry::JacobianInverseTransposed
396  // ----------------------------------------------
397 
398  template< class ct, int mydim, int cdim, class Traits >
399  class MultiLinearGeometry< ct, mydim, cdim, Traits >::JacobianInverseTransposed
400  : public FieldMatrix< ctype, coorddimension, mydimension >
401  {
402  typedef FieldMatrix< ctype, coorddimension, mydimension > Base;
403 
404  public:
405  void setup ( const JacobianTransposed &jt )
406  {
407  detInv_ = MatrixHelper::template rightInvA< mydimension, coorddimension >( jt, static_cast< Base & >( *this ) );
408  }
409 
410  void setupDeterminant ( const JacobianTransposed &jt )
411  {
412  detInv_ = MatrixHelper::template sqrtDetAAT< mydimension, coorddimension >( jt );
413  }
414 
415  ctype det () const { return ctype( 1 ) / detInv_; }
416  ctype detInv () const { return detInv_; }
417 
418  private:
419  ctype detInv_;
420  };
421 
422 
423 
436  template< class ct, int mydim, int cdim, class Traits = MultiLinearGeometryTraits< ct > >
438  : public MultiLinearGeometry< ct, mydim, cdim, Traits >
439  {
442 
443  protected:
445 
446  public:
448 
449  typedef typename Base::ctype ctype;
450 
451  using Base::mydimension;
452  using Base::coorddimension;
453 
456 
459 
460  template< class CornerStorage >
461  CachedMultiLinearGeometry ( const ReferenceElement &refElement, const CornerStorage &cornerStorage )
462  : Base( refElement, cornerStorage ),
463  affine_( Base::affine() ),
464  jacobianInverseTransposedComputed_( false ),
465  integrationElementComputed_( false )
466  {}
467 
468  template< class CornerStorage >
469  CachedMultiLinearGeometry ( Dune::GeometryType gt, const CornerStorage &cornerStorage )
470  : Base( gt, cornerStorage ),
471  affine_( Base::affine() ),
472  jacobianInverseTransposedComputed_( false ),
473  integrationElementComputed_( false )
474  {}
475 
477  bool affine () const { return affine_; }
478 
479  using Base::corner;
480 
482  GlobalCoordinate center () const { return global( refElement().position( 0, 0 ) ); }
483 
491  {
492  if( affine() )
493  {
495  jacobianTransposed_.umtv( local, global );
496  return global;
497  }
498  else
499  return Base::global( local );
500  }
501 
514  {
515  if( affine() )
516  {
518  if( jacobianInverseTransposedComputed_ )
519  jacobianInverseTransposed_.mtv( global - corner( 0 ), local );
520  else
521  MatrixHelper::template xTRightInvA< mydimension, coorddimension >( jacobianTransposed_, global - corner( 0 ), local );
522  return local;
523  }
524  else
525  return Base::local( global );
526  }
527 
543  {
544  if( affine() )
545  {
546  if( !integrationElementComputed_ )
547  {
549  integrationElementComputed_ = true;
550  }
552  }
553  else
554  return Base::integrationElement( local );
555  }
556 
558  ctype volume () const
559  {
560  if( affine() )
561  return integrationElement( refElement().position( 0, 0 ) ) * refElement().volume();
562  else
563  return Base::volume();
564  }
565 
576  {
577  if( affine() )
578  return jacobianTransposed_;
579  else
580  return Base::jacobianTransposed( local );
581  }
582 
590  {
591  if( affine() )
592  {
593  if( !jacobianInverseTransposedComputed_ )
594  {
596  jacobianInverseTransposedComputed_ = true;
597  integrationElementComputed_ = true;
598  }
600  }
601  else
602  return Base::jacobianInverseTransposed( local );
603  }
604 
605  protected:
606  using Base::refElement;
607 
610 
611  private:
612  mutable bool affine_ : 1;
613 
614  mutable bool jacobianInverseTransposedComputed_ : 1;
615  mutable bool integrationElementComputed_ : 1;
616  };
617 
618 
619 
620  // Implementation of MultiLinearGeometry
621  // -------------------------------------
622 
623  template< class ct, int mydim, int cdim, class Traits >
624  inline const typename MultiLinearGeometry< ct, mydim, cdim, Traits >::JacobianInverseTransposed &
627  {
628  jacobianInverseTransposed_.setup( jacobianTransposed( local ) );
629  return jacobianInverseTransposed_;
630  }
631 
632 
633  template< class ct, int mydim, int cdim, class Traits >
634  template< bool add, int dim >
636  ::global ( TopologyId topologyId, integral_constant< int, dim >,
637  CornerIterator &cit, const ctype &df, const LocalCoordinate &x,
638  const ctype &rf, GlobalCoordinate &y )
639  {
640  const ctype xn = df*x[ dim-1 ];
641  const ctype cxn = ctype( 1 ) - xn;
642  assert( (xn > -Traits::tolerance()) && (cxn > -Traits::tolerance()) );
643 
644  if( GenericGeometry::isPrism( topologyId, mydimension, mydimension-dim ) )
645  {
646  // apply (1-xn) times mapping for bottom
647  global< add >( topologyId, integral_constant< int, dim-1 >(), cit, df, x, rf*cxn, y );
648  // apply xn times mapping for top
649  global< true >( topologyId, integral_constant< int, dim-1 >(), cit, df, x, rf*xn, y );
650  }
651  else
652  {
653  assert( GenericGeometry::isPyramid( topologyId, mydimension, mydimension-dim ) );
654  // apply (1-xn) times mapping for bottom (with argument x/(1-xn))
655  if( cxn > Traits::tolerance() )
656  global< add >( topologyId, integral_constant< int, dim-1 >(), cit, df/cxn, x, rf*cxn, y );
657  else
658  global< add >( topologyId, integral_constant< int, dim-1 >(), cit, df, x, ctype( 0 ), y );
659  // apply xn times the tip
660  y.axpy( rf*xn, *cit );
661  ++cit;
662  }
663  }
664 
665  template< class ct, int mydim, int cdim, class Traits >
666  template< bool add >
668  ::global ( TopologyId topologyId, integral_constant< int, 0 >,
669  CornerIterator &cit, const ctype &df, const LocalCoordinate &x,
670  const ctype &rf, GlobalCoordinate &y )
671  {
672  const GlobalCoordinate &origin = *cit;
673  ++cit;
674  for( int i = 0; i < coorddimension; ++i )
675  y[ i ] = (add ? y[ i ] + rf*origin[ i ] : rf*origin[ i ]);
676  }
677 
678 
679  template< class ct, int mydim, int cdim, class Traits >
680  template< bool add, int rows, int dim >
682  ::jacobianTransposed ( TopologyId topologyId, integral_constant< int, dim >,
683  CornerIterator &cit, const ctype &df, const LocalCoordinate &x,
684  const ctype &rf, FieldMatrix< ctype, rows, cdim > &jt )
685  {
686  assert( rows >= dim );
687 
688  const ctype xn = df*x[ dim-1 ];
689  const ctype cxn = ctype( 1 ) - xn;
690  assert( (xn > -Traits::tolerance()) && (cxn > -Traits::tolerance()) );
691 
692  CornerIterator cit2( cit );
693  if( GenericGeometry::isPrism( topologyId, mydimension, mydimension-dim ) )
694  {
695  // apply (1-xn) times Jacobian for bottom
696  jacobianTransposed< add >( topologyId, integral_constant< int, dim-1 >(), cit2, df, x, rf*cxn, jt );
697  // apply xn times Jacobian for top
698  jacobianTransposed< true >( topologyId, integral_constant< int, dim-1 >(), cit2, df, x, rf*xn, jt );
699  // compute last row as difference between top value and bottom value
700  global< add >( topologyId, integral_constant< int, dim-1 >(), cit, df, x, -rf, jt[ dim-1 ] );
701  global< true >( topologyId, integral_constant< int, dim-1 >(), cit, df, x, rf, jt[ dim-1 ] );
702  }
703  else
704  {
705  assert( GenericGeometry::isPyramid( topologyId, mydimension, mydimension-dim ) );
706  // initialize last row
707  global< add >( topologyId, integral_constant< int, dim-1 >(), cit, df/cxn, x, -rf, jt[ dim-1 ] );
708  jt[ dim-1 ].axpy( rf, *cit );
709  ++cit;
710  // apply Jacobian for bottom (with argument x/(1-xn)) and correct last row
711  if( add )
712  {
713  FieldMatrix< ctype, dim-1, coorddimension > jt2;
714  jacobianTransposed< false >( topologyId, integral_constant< int, dim-1 >(), cit2, df/cxn, x, rf, jt2 );
715  for( int j = 0; j < dim-1; ++j )
716  {
717  jt[ j ] += jt2[ j ];
718  jt[ dim-1 ].axpy( (df/cxn)*x[ j ], jt2[ j ] );
719  }
720  }
721  else
722  {
723  jacobianTransposed< false >( topologyId, integral_constant< int, dim-1 >(), cit2, df/cxn, x, rf, jt );
724  for( int j = 0; j < dim-1; ++j )
725  jt[ dim-1 ].axpy( (df/cxn)*x[ j ], jt[ j ] );
726  }
727  }
728  }
729 
730  template< class ct, int mydim, int cdim, class Traits >
731  template< bool add, int rows >
733  ::jacobianTransposed ( TopologyId topologyId, integral_constant< int, 0 >,
734  CornerIterator &cit, const ctype &df, const LocalCoordinate &x,
735  const ctype &rf, FieldMatrix< ctype, rows, cdim > &jt )
736  {
737  ++cit;
738  }
739 
740 
741 
742  template< class ct, int mydim, int cdim, class Traits >
743  template< int dim >
745  ::affine ( TopologyId topologyId, integral_constant< int, dim >, CornerIterator &cit, JacobianTransposed &jt )
746  {
747  const GlobalCoordinate &orgBottom = *cit;
748  if( !affine( topologyId, integral_constant< int, dim-1 >(), cit, jt ) )
749  return false;
750  const GlobalCoordinate &orgTop = *cit;
751 
752  if( GenericGeometry::isPrism( topologyId, mydimension, mydimension-dim ) )
753  {
754  JacobianTransposed jtTop;
755  if( !affine( topologyId, integral_constant< int, dim-1 >(), cit, jtTop ) )
756  return false;
757 
758  // check whether both jacobians are identical
759  ctype norm( 0 );
760  for( int i = 0; i < dim-1; ++i )
761  norm += (jtTop[ i ] - jt[ i ]).two_norm2();
762  if( norm >= Traits::tolerance() )
763  return false;
764  }
765  else
766  ++cit;
767  jt[ dim-1 ] = orgTop - orgBottom;
768  return true;
769  }
770 
771  template< class ct, int mydim, int cdim, class Traits >
773  ::affine ( TopologyId topologyId, integral_constant< int, 0 >, CornerIterator &cit, JacobianTransposed &jt )
774  {
775  ++cit;
776  return true;
777  }
778 
779 } // namespace Dune
780 
781 #endif // #ifndef DUNE_GEOMETRY_MULTILINEARGEOMETRY_HH