Dune Core Modules (2.4.1)

mappings.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_ALU3DGRIDMAPPINGS_HH
4#define DUNE_ALU3DGRIDMAPPINGS_HH
5
6// System includes
7#include <limits>
8#include <cmath>
9
10// Dune includes
14
15// Local includes
16#include "alu3dinclude.hh"
17
18namespace Dune {
19
20 static const alu3d_ctype ALUnumericEpsilon = 10.0 * std::numeric_limits< alu3d_ctype >::epsilon();
21
22 template<int mydim, int coorddim, class GridImp>
23 class ALU3dGridGeometry;
24
25 template< ALU3dGridElementType, class >
26 class ALU3dGrid;
27
31 {
32 public:
33 typedef alu3d_ctype double_t[3];
36 private:
37 static const double _epsilon ;
38
39 // the internal mapping
40 alu3d_ctype a [8][3] ;
41 mat_t Df;
42 mat_t Dfi;
43 mat_t invTransposed_;
44 alu3d_ctype DetDf ;
45
46 bool calcedDet_;
47 bool calcedLinear_;
48 bool calcedInv_;
49 bool affine_;
50
51 void linear (const alu3d_ctype, const alu3d_ctype, const alu3d_ctype) ;
52 void linear (const coord_t&) ;
53 void inverse (const coord_t&) ;
54 public:
55 TrilinearMapping (const coord_t&, const coord_t&,
56 const coord_t&, const coord_t&,
57 const coord_t&, const coord_t&,
58 const coord_t&, const coord_t&);
59
60 // only to call from geometry class
62
64
66 alu3d_ctype det (const coord_t&) ;
67 const mat_t& jacobianInverseTransposed(const coord_t&);
68 const mat_t& jacobianTransposed(const coord_t&);
69 void map2world (const coord_t&, coord_t&) const ;
70 void map2world (const alu3d_ctype , const alu3d_ctype , const alu3d_ctype ,
71 coord_t&) const ;
72 void world2map (const coord_t&, coord_t&) ;
73
74 template <class vector_t>
75 void buildMapping(const vector_t&, const vector_t&,
76 const vector_t&, const vector_t&,
77 const vector_t&, const vector_t&,
78 const vector_t&, const vector_t&);
79
80 // returns true if mapping is affine
81 inline bool affine () const { return affine_; }
82 };
83
85 // NOTE: this class is different to the BilinearSurfaceMapping in
86 // ALUGrid, for example the reference elements differ
87 // here we have [0,1]^2 and in ALUGrid its [-1,1]^2
88 // also the point numbering is different
90 {
91 public:
92 // our coordinate types
95
96 // type of coordinate vectors from elements
97 typedef alu3d_ctype double3_t[3];
98 protected:
99
100 alu3d_ctype _n [3][3] ;
101
102 static const double _epsilon ;
103
104 bool _affine;
105
106 public:
109
112
113 // returns true if mapping is affine
114 inline bool affine () const { return _affine ; }
115
116 // calcuates normal
117 void normal(const coord2_t&, coord3_t&) const ;
118 void normal(const alu3d_ctype, const alu3d_ctype, coord3_t&) const;
119
120 void negativeNormal(const coord2_t&, coord3_t&) const ;
121 void negativeNormal(const alu3d_ctype, const alu3d_ctype, coord3_t&) const;
122
123 public:
124 // builds _b and _n, called from the constructors
125 // public because also used in faceutility
126 template <class vector_t>
127 void buildMapping (const vector_t & , const vector_t & ,
128 const vector_t & , const vector_t & );
129 protected:
130 // builds _b and _n, called from the constructors
131 // public because also used in faceutility
132 template <class vector_t>
133 void buildMapping (const vector_t & , const vector_t & ,
134 const vector_t & , const vector_t & ,
135 alu3d_ctype (&_b)[4][3] );
136 } ;
137
138
140 // NOTE: this class is different to the BilinearSurfaceMapping in
141 // ALUGrid, for example the reference elements differ
142 // here we have [0,1]^2 and in ALUGrid its [-1,1]^2
143 // also the point numbering is different
145 {
146 protected:
148
149 using BaseType :: _n;
150 static const double _epsilon;
151
152 // our coordinate types
155
156 // type of coordinate vectors from elements
157 typedef alu3d_ctype double3_t[3];
158
159 // type for helper matrices
161
162 // type for inverse matrices
164
165 // type for inverse matrices
167
168 alu3d_ctype _b [4][3] ;
169
170 mutable mat3_t Df,Dfi;
171 mutable inv_t invTransposed_;
172 mutable matrix_t matrix_;
173 mutable alu3d_ctype DetDf;
174
175 mutable coord3_t normal_;
176 mutable coord3_t tmp_;
177
178 mutable bool _calcedInv;
179 mutable bool _calcedTransposed;
180 mutable bool _calcedMatrix;
181
182 public:
185
188 const coord3_t&, const coord3_t&) ;
190 BilinearSurfaceMapping (const double3_t &, const double3_t &,
191 const double3_t &, const double3_t &) ;
194
195 void inverse (const coord3_t&) const;
196 const inv_t& jacobianInverseTransposed(const coord2_t&) const ;
197
198 const matrix_t& jacobianTransposed(const coord2_t&) const ;
199
200 // calculates determinant of face mapping using the normal
201 alu3d_ctype det(const coord2_t&) const;
202
203 // maps from local coordinates to global coordinates
204 void world2map(const coord3_t &, coord2_t & ) const;
205
206 // maps form global coordinates to local (within reference element)
207 // coordinates
208 void map2world(const coord2_t&, coord3_t&) const ;
209 void map2world(const alu3d_ctype ,const alu3d_ctype , coord3_t&) const ;
210
211 private:
212 void map2worldnormal(const alu3d_ctype, const alu3d_ctype, const alu3d_ctype , coord3_t&) const;
213 void map2worldlinear(const alu3d_ctype, const alu3d_ctype, const alu3d_ctype ) const;
214
215 public:
216 // builds _b and _n, called from the constructors
217 // public because also used in faceutility
218 template <class vector_t>
219 void buildMapping (const vector_t & , const vector_t & ,
220 const vector_t & , const vector_t & );
221 } ;
222
223
224
226 template< int cdim >
228 {
229 public:
230 typedef alu3d_ctype ctype;
231
234
237
238 protected:
239 ctype _b [4][cdim];
240
241 mutable ctype det_;
242 mutable matrix_t matrix_;
243 mutable inv_t invTransposed_;
244
245 mutable bool affine_;
246 mutable bool calcedMatrix_;
247 mutable bool calcedDet_;
248 mutable bool calcedInv_;
249
250 public:
252 BilinearMapping ( const world_t &p0, const world_t &p1,
253 const world_t &p2, const world_t &p3 );
254 BilinearMapping ( const ctype (&p0)[ cdim ], const ctype (&p1)[ cdim ],
255 const ctype (&p2)[ cdim ], const ctype (&p3)[ cdim ] );
256
257 bool affine () const;
258
259 void world2map ( const world_t &, map_t & ) const;
260 void map2world ( const ctype x, const ctype y, world_t &w ) const;
261 void map2world ( const map_t &, world_t & ) const;
262
263 ctype det ( const map_t & ) const;
264
265 const matrix_t &jacobianTransposed ( const map_t & ) const;
266 const inv_t &jacobianInverseTransposed ( const map_t & ) const;
267
268 template< class vector_t >
269 void buildMapping ( const vector_t &, const vector_t &,
270 const vector_t &, const vector_t & );
271
272 protected:
273 static void multTransposedMatrix ( const matrix_t &, FieldMatrix< ctype, 2, 2 > & );
274 static void multMatrix ( const matrix_t &, const FieldMatrix< ctype, 2, 2 > &, inv_t & );
275
276 void map2worldlinear ( const ctype, const ctype ) const;
277 void inverse ( const map_t & ) const;
278 };
279
280
281
283 template< int cdim, int mydim >
285 {
286 public:
287 typedef alu3d_ctype ctype;
288
289 typedef ctype double_t[ cdim ];
290
293
296
297 protected:
300 world_t _p0;
301
303 mutable ctype _det;
304
306 mutable bool _calcedInv;
307
309 mutable bool _calcedDet;
310
311 public:
314
317
318 // returns true if mapping is affine (which is always true)
319 inline bool affine () const { return true ; }
320
321 // return reference to transposed jacobian
322 const matrix_t& jacobianTransposed(const map_t &) const ;
323
324 // return reference to transposed jacobian inverse
325 const inv_t& jacobianInverseTransposed(const map_t &) const ;
326
327 // calculates determinant of mapping
328 ctype det(const map_t&) const;
329
330 // maps from local coordinates to global coordinates
331 void world2map(const world_t &, map_t &) const;
332
333 // maps form global coordinates to local (within reference element)
334 // coordinates
335 void map2world(const map_t &, world_t &) const ;
336
337 protected:
338 // calculate inverse
339 void inverse (const map_t&) const;
340
341 // calculate inverse one codim one entity
342 void inverseCodimOne (const map_t&) const;
343
344 // calculate determinant
345 void calculateDeterminant (const map_t&) const;
346
347 void multTransposedMatrix(const matrix_t& matrix,
349
350 void multMatrix ( const matrix_t& A,
352 inv_t& ret ) const ;
353
354 public:
355 // builds _b and _n, called from the constructors
356 // public because also used in faceutility
357 template <class vector_t>
358 void buildMapping (const vector_t & , const vector_t & ,
359 const vector_t & , const vector_t & );
360
361 // builds _b and _n, called from the constructors
362 // public because also used in faceutility
363 template <class vector_t>
364 void buildMapping (const vector_t & , const vector_t & ,
365 const vector_t & );
366
367 // builds _b and _n, called from the constructors
368 // public because also used in faceutility
369 template <class vector_t>
370 void buildMapping (const vector_t & , const vector_t & );
371
372 template <class vector_t>
373 void buildMapping (const vector_t & );
374 } ;
375
376
378 //
379 // NonConforming Mappings
380 //
382
383
386 template< ALU3dGridElementType type, class Comm >
388
390 template< class Comm >
391 class NonConformingFaceMapping< tetra, Comm >
392 {
393 public:
395 typedef typename ALU3dImplTraits< tetra, Comm >::HfaceRuleType RefinementRuleType;
396
397 NonConformingFaceMapping ( RefinementRuleType rule, int nChild )
398 : rule_( rule ), nChild_( nChild )
399 {}
400
401 void child2parent ( const CoordinateType &childCoordinates,
402 CoordinateType &parentCoordinates) const;
403
404 CoordinateType child2parent ( const FieldVector< alu3d_ctype, 2 > &childCoordinates ) const;
405
406 private:
407 void child2parentNosplit(const CoordinateType& childCoordinates,
408 CoordinateType& parentCoordinates) const;
409 void child2parentE01(const CoordinateType& childCoordinates,
410 CoordinateType& parentCoordinates) const;
411 void child2parentE12(const CoordinateType& childCoordinates,
412 CoordinateType& parentCoordinates) const;
413 void child2parentE20(const CoordinateType& childCoordinates,
414 CoordinateType& parentCoordinates) const;
415 void child2parentIso4(const CoordinateType& childCoordinates,
416 CoordinateType& parentCoordinates) const;
417
418 RefinementRuleType rule_;
419 int nChild_;
420 };
421
423 template< class Comm >
424 class NonConformingFaceMapping< hexa, Comm >
425 {
426 public:
428 typedef typename ALU3dImplTraits< hexa, Comm >::HfaceRuleType RefinementRuleType;
429
430 NonConformingFaceMapping ( RefinementRuleType rule, int nChild )
431 : rule_( rule ), nChild_( nChild )
432 {}
433
434 void child2parent ( const CoordinateType &childCoordinates,
435 CoordinateType &parentCoordinates) const;
436
437 CoordinateType child2parent ( const FieldVector< alu3d_ctype, 2 > &childCoordinates ) const;
438
439 private:
440 void child2parentNosplit(const CoordinateType& childCoordinates,
441 CoordinateType& parentCoordinates) const;
442 void child2parentIso4(const CoordinateType& childCoordinates,
443 CoordinateType& parentCoordinates) const;
444
445 RefinementRuleType rule_;
446 int nChild_;
447 };
448
449} // end namespace Dune
450
451#if COMPILE_ALUGRID_INLINE
452 #include "mappings_imp.cc"
453#endif
454#endif
A bilinear mapping.
Definition: mappings.hh:228
A bilinear surface mapping.
Definition: mappings.hh:145
BilinearSurfaceMapping(const coord3_t &, const coord3_t &, const coord3_t &, const coord3_t &)
Constructor getting FieldVectors.
BilinearSurfaceMapping()
Constructor creating empty mapping with double , i.e. zero.
BilinearSurfaceMapping(const double3_t &, const double3_t &, const double3_t &, const double3_t &)
Constructor for double[3].
A linear mapping.
Definition: mappings.hh:285
inv_t _invTransposed
storage for inverse of jacobian (transposed)
Definition: mappings.hh:299
LinearMapping(const LinearMapping &)
copy constructor
LinearMapping()
Constructor creating empty mapping with double , i.e. zero.
bool _calcedInv
true if inverse has been calculated
Definition: mappings.hh:306
bool _calcedDet
true if determinant has been calculated
Definition: mappings.hh:309
matrix_t _matrix
transformation matrix (transposed)
Definition: mappings.hh:298
ctype _det
P[0].
Definition: mappings.hh:303
Definition: mappings.hh:387
A bilinear surface mapping.
Definition: mappings.hh:90
SurfaceNormalCalculator()
Constructor creating empty mapping with double , i.e. zero.
Definition: mappings.hh:31
A few common exception classes.
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.
Dune namespace.
Definition: alignment.hh:10
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Dec 22, 23:30, 2024)