ug4
reference_mapping_provider.h
Go to the documentation of this file.
1 /*
2  * Copyright (c) 2012-2015: G-CSC, Goethe University Frankfurt
3  * Author: Andreas Vogel
4  *
5  * This file is part of UG4.
6  *
7  * UG4 is free software: you can redistribute it and/or modify it under the
8  * terms of the GNU Lesser General Public License version 3 (as published by the
9  * Free Software Foundation) with the following additional attribution
10  * requirements (according to LGPL/GPL v3 §7):
11  *
12  * (1) The following notice must be displayed in the Appropriate Legal Notices
13  * of covered and combined works: "Based on UG4 (www.ug4.org/license)".
14  *
15  * (2) The following notice must be displayed at a prominent place in the
16  * terminal output of covered works: "Based on UG4 (www.ug4.org/license)".
17  *
18  * (3) The following bibliography is recommended for citation and must be
19  * preserved in all covered files:
20  * "Reiter, S., Vogel, A., Heppner, I., Rupp, M., and Wittum, G. A massively
21  * parallel geometric multigrid solver on hierarchically distributed grids.
22  * Computing and visualization in science 16, 4 (2013), 151-164"
23  * "Vogel, A., Reiter, S., Rupp, M., Nägel, A., and Wittum, G. UG4 -- a novel
24  * flexible software system for simulating pde based models on high performance
25  * computers. Computing and visualization in science 16, 4 (2013), 165-179"
26  *
27  * This program is distributed in the hope that it will be useful,
28  * but WITHOUT ANY WARRANTY; without even the implied warranty of
29  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
30  * GNU Lesser General Public License for more details.
31  */
32 
33 #ifndef __H__UG__LIB_DISC__REFERENCE_ELEMENT__REFERENCE_MAPPING_PROVIDER__
34 #define __H__UG__LIB_DISC__REFERENCE_ELEMENT__REFERENCE_MAPPING_PROVIDER__
35 
36 
37 #include "common/common.h"
38 #include "common/math/ugmath.h"
40 
41 namespace ug{
42 
44 
51 template <int TDim, int TWorldDim = TDim>
53 {
54  public:
56  static const int worldDim = TWorldDim;
57 
59  static const int dim = TDim;
60 
61  public:
63  virtual bool is_linear() const = 0;
64 
66  virtual void update(const MathVector<worldDim>* vCornerCoord) = 0;
67 
69  virtual void update(const std::vector<MathVector<worldDim> >& vCornerCoord) = 0;
70 
72  virtual void local_to_global(MathVector<worldDim>& globPos,
73  const MathVector<dim>& locPos) const = 0;
74 
76  virtual void local_to_global(MathVector<worldDim>* vGlobPos,
77  const MathVector<dim>* vLocPos, size_t n) const = 0;
78 
80  virtual void local_to_global(std::vector<MathVector<worldDim> >& vGlobPos,
81  const std::vector<MathVector<dim> >& vLocPos) const = 0;
82 
84  virtual void global_to_local(MathVector<dim>& locPos,
85  const MathVector<worldDim>& globPos,
86  const size_t maxIter = 1000,
87  const number tol = 1e-10) const = 0;
88 
90  virtual void global_to_local(MathVector<dim>* vLocPos,
91  const MathVector<worldDim>* vGlobPos, size_t n,
92  const size_t maxIter = 1000,
93  const number tol = 1e-10) const = 0;
94 
96  virtual void global_to_local(std::vector<MathVector<dim> >& vLocPos,
97  const std::vector<MathVector<worldDim> >& vGlobPos,
98  const size_t maxIter = 1000,
99  const number tol = 1e-10) const = 0;
100 
103  const MathVector<dim>& locPos) const = 0;
104 
107  const MathVector<dim>* vLocPos, size_t n) const = 0;
108 
110  virtual void jacobian(std::vector<MathMatrix<worldDim, dim> >& vJ,
111  const std::vector<MathVector<dim> >& vLocPos) const = 0;
112 
115  const MathVector<dim>& locPos) const = 0;
116 
119  const MathVector<dim>* vLocPos, size_t n) const = 0;
120 
122  virtual void jacobian_transposed(std::vector<MathMatrix<dim, worldDim> >& vJT,
123  const std::vector<MathVector<dim> >& vLocPos) const = 0;
124 
127  const MathVector<dim>& locPos) const = 0;
128 
131  const MathVector<dim>* vLocPos, size_t n) const = 0;
132 
135  number* vDet,
136  const MathVector<dim>* vLocPos, size_t n) const = 0;
137 
139  virtual void jacobian_transposed_inverse(std::vector<MathMatrix<worldDim, dim> >& vJTInv,
140  const std::vector<MathVector<dim> >& vLocPos) const = 0;
141 
143  virtual void jacobian_transposed_inverse(std::vector<MathMatrix<worldDim, dim> >& vJTInv,
144  std::vector<number>& vDet,
145  const std::vector<MathVector<dim> >& vLocPos) const = 0;
146 
148  virtual number sqrt_gram_det(const MathVector<dim>& locPos) const = 0;
149 
151  virtual void sqrt_gram_det(number* vDet, const MathVector<dim>* vLocPos, size_t n) const = 0;
152 
154  virtual void sqrt_gram_det(std::vector<number>& vDet,
155  const std::vector<MathVector<dim> >& vLocPos) const = 0;
156 
158  virtual ~DimReferenceMapping() {}
159 };
160 
163 
165 
169  private:
170  // disallow private constructor
172 
173  // disallow copy and assignment (intentionally left unimplemented)
176 
177  // private destructor
179 
180  // Singleton provider
182  {
183  static ReferenceMappingProvider myInst;
184  return myInst;
185  };
186 
187  // This is very dirty implementation, since casting to void. But, it is
188  // efficient, easy and typesafe. Maybe it should be changed to something
189  // inherent typesafe (not using casts) some day
190  // holding all mappings (worldDim x dim x roid)
192 
193  // casts void to map
194  template <int TDim, int TWorldDim>
196  {
197  UG_STATIC_ASSERT(TDim <= 3, only_implemented_for_ref_dim_smaller_equal_3);
198  UG_STATIC_ASSERT(TWorldDim <= 3, only_implemented_for_ref_dim_smaller_equal_3);
199  UG_ASSERT(roid < NUM_REFERENCE_OBJECTS, "Roid specified incorrectly.");
200  UG_ASSERT(roid >= 0, "Roid specified incorrectly.");
201  return reinterpret_cast<DimReferenceMapping<TDim, TWorldDim>*>
202  (m_vvvMapping[TDim][TWorldDim][roid]);
203  }
204 
205  // casts map to void
206  template <int TDim, int TWorldDim>
208  {
209  m_vvvMapping[TDim][TWorldDim][roid] = reinterpret_cast<void*>(&map);
210  }
211 
212  public:
214 
224  template <int TDim, int TWorldDim>
226  {
227  DimReferenceMapping<TDim, TWorldDim>* pMap = inst().get_mapping<TDim, TWorldDim>(roid);
228  if(!pMap){
229  UG_THROW("ReferenceMappingProvider: ReferenceMapping not found for "
230  <<roid<<" from R^"<<TDim<<" to R^"<<TWorldDim);
231  }
232  else return *pMap;
233  }
234 
236 
247  template <int TDim, int TWorldDim>
249  const std::vector<MathVector<TWorldDim> >& vCornerCoord)
250  {
251  DimReferenceMapping<TDim, TWorldDim>& map = get<TDim, TWorldDim>(roid);
252  map.update(vCornerCoord);
253  return map;
254  }
255 
257 
268  template <int TDim, int TWorldDim>
270  const MathVector<TWorldDim>* vCornerCoord)
271  {
272  DimReferenceMapping<TDim, TWorldDim>& map = get<TDim, TWorldDim>(roid);
273  map.update(vCornerCoord);
274  return map;
275  }
276 
277 };
278 
279 } // end namespace ug
280 
281 #endif /* __H__UG__LIB_DISC__REFERENCE_ELEMENT__REFERENCE_MAPPING_PROVIDER__ */
virtual base class for reference mappings
Definition: reference_mapping_provider.h:53
virtual void jacobian_transposed(MathMatrix< dim, worldDim > *vJT, const MathVector< dim > *vLocPos, size_t n) const =0
returns transposed of jacobian for n local positions
virtual void jacobian_transposed_inverse(MathMatrix< worldDim, dim > *vJTInv, number *vDet, const MathVector< dim > *vLocPos, size_t n) const =0
returns transposed of the inverse of the jacobian for n local positions
virtual void jacobian_transposed(MathMatrix< dim, worldDim > &JT, const MathVector< dim > &locPos) const =0
returns transposed of jacobian
virtual void jacobian(std::vector< MathMatrix< worldDim, dim > > &vJ, const std::vector< MathVector< dim > > &vLocPos) const =0
returns jacobian for a vector of local positions
virtual void local_to_global(std::vector< MathVector< worldDim > > &vGlobPos, const std::vector< MathVector< dim > > &vLocPos) const =0
map local coordinate to global coordinate for a vector of local positions
virtual void jacobian_transposed_inverse(std::vector< MathMatrix< worldDim, dim > > &vJTInv, const std::vector< MathVector< dim > > &vLocPos) const =0
returns transposed of the inverse of the jacobian for a vector of positions
virtual void global_to_local(std::vector< MathVector< dim > > &vLocPos, const std::vector< MathVector< worldDim > > &vGlobPos, const size_t maxIter=1000, const number tol=1e-10) const =0
map global coordinate to local coordinate for a vector of local positions
static const int dim
reference dimension (domain space dimension)
Definition: reference_mapping_provider.h:59
virtual number sqrt_gram_det(const MathVector< dim > &locPos) const =0
returns the determinate of the jacobian
static const int worldDim
world dimension (range space dimension)
Definition: reference_mapping_provider.h:56
virtual void sqrt_gram_det(number *vDet, const MathVector< dim > *vLocPos, size_t n) const =0
returns the determinate of the jacobian for n local positions
virtual void global_to_local(MathVector< dim > *vLocPos, const MathVector< worldDim > *vGlobPos, size_t n, const size_t maxIter=1000, const number tol=1e-10) const =0
map global coordinate to local coordinate for n local positions
virtual void sqrt_gram_det(std::vector< number > &vDet, const std::vector< MathVector< dim > > &vLocPos) const =0
returns the determinate of the jacobian for a vector of local positions
virtual ~DimReferenceMapping()
virtual destructor
Definition: reference_mapping_provider.h:158
virtual number jacobian_transposed_inverse(MathMatrix< worldDim, dim > &JTInv, const MathVector< dim > &locPos) const =0
returns transposed of the inverse of the jacobian and returns sqrt of gram determinante
virtual void jacobian_transposed_inverse(MathMatrix< worldDim, dim > *vJTInv, const MathVector< dim > *vLocPos, size_t n) const =0
returns transposed of the inverse of the jacobian for n local positions
virtual void jacobian_transposed(std::vector< MathMatrix< dim, worldDim > > &vJT, const std::vector< MathVector< dim > > &vLocPos) const =0
returns transposed of jacobian for a vector of positions
virtual void local_to_global(MathVector< worldDim > &globPos, const MathVector< dim > &locPos) const =0
map local coordinate to global coordinate
virtual void local_to_global(MathVector< worldDim > *vGlobPos, const MathVector< dim > *vLocPos, size_t n) const =0
map n local coordinate to global coordinate
virtual void update(const MathVector< worldDim > *vCornerCoord)=0
refresh mapping for new set of corners
virtual void update(const std::vector< MathVector< worldDim > > &vCornerCoord)=0
refresh mapping for new set of corners
virtual void jacobian_transposed_inverse(std::vector< MathMatrix< worldDim, dim > > &vJTInv, std::vector< number > &vDet, const std::vector< MathVector< dim > > &vLocPos) const =0
returns transposed of the inverse of the jacobian for a vector of positions
virtual void global_to_local(MathVector< dim > &locPos, const MathVector< worldDim > &globPos, const size_t maxIter=1000, const number tol=1e-10) const =0
map global coordinate to local coordinate
virtual bool is_linear() const =0
returns if mapping is affine
virtual void jacobian(MathMatrix< worldDim, dim > &J, const MathVector< dim > &locPos) const =0
returns jacobian
virtual void jacobian(MathMatrix< worldDim, dim > *vJ, const MathVector< dim > *vLocPos, size_t n) const =0
returns jacobian for n local positions
class to provide reference mappings
Definition: reference_mapping_provider.h:168
ReferenceMappingProvider & operator=(const ReferenceMappingProvider &)
void * m_vvvMapping[4][4][NUM_REFERENCE_OBJECTS]
Definition: reference_mapping_provider.h:185
void set_mapping(ReferenceObjectID roid, DimReferenceMapping< TDim, TWorldDim > &map)
Definition: reference_mapping_provider.h:207
DimReferenceMapping< TDim, TWorldDim > * get_mapping(ReferenceObjectID roid)
Definition: reference_mapping_provider.h:195
ReferenceMappingProvider(const ReferenceMappingProvider &)
static DimReferenceMapping< TDim, TWorldDim > & get(ReferenceObjectID roid, const MathVector< TWorldDim > *vCornerCoord)
returns a reference to a DimReferenceMapping with updated element corners
Definition: reference_mapping_provider.h:269
ReferenceMappingProvider()
Definition: reference_mapping_provider.cpp:222
static DimReferenceMapping< TDim, TWorldDim > & get(ReferenceObjectID roid, const std::vector< MathVector< TWorldDim > > &vCornerCoord)
returns a reference to a DimReferenceMapping with updated element corners
Definition: reference_mapping_provider.h:248
static DimReferenceMapping< TDim, TWorldDim > & get(ReferenceObjectID roid)
returns a reference to a DimReferenceMapping
Definition: reference_mapping_provider.h:225
~ReferenceMappingProvider()
Definition: reference_mapping_provider.h:178
static ReferenceMappingProvider & inst()
Definition: reference_mapping_provider.h:181
#define UG_ASSERT(expr, msg)
Definition: assert.h:70
#define UG_THROW(msg)
Definition: error.h:57
#define UG_STATIC_ASSERT(expr, msg)
Checks an expression at compile-time and raises a compile-error if the expression equals 0.
Definition: static_assert.h:63
double number
Definition: types.h:124
the ug namespace
ReferenceObjectID
these ids are used to identify the shape of a geometric object.
Definition: grid_base_objects.h:74
@ NUM_REFERENCE_OBJECTS
Definition: grid_base_objects.h:85