BEM++  2.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Groups Pages
geometry.hpp
1 // Copyright (C) 2011-2012 by the BEM++ Authors
2 //
3 // Permission is hereby granted, free of charge, to any person obtaining a copy
4 // of this software and associated documentation files (the "Software"), to deal
5 // in the Software without restriction, including without limitation the rights
6 // to use, copy, modify, merge, publish, distribute, sublicense, and/or sell
7 // copies of the Software, and to permit persons to whom the Software is
8 // furnished to do so, subject to the following conditions:
9 //
10 // The above copyright notice and this permission notice shall be included in
11 // all copies or substantial portions of the Software.
12 //
13 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR
14 // IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
15 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE
16 // AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
17 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING FROM,
18 // OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN
19 // THE SOFTWARE.
20 
21 #ifndef bempp_geometry_hpp
22 #define bempp_geometry_hpp
23 
24 #include "../common/common.hpp"
25 
26 #include "geometry_type.hpp"
27 
28 #include "../common/armadillo_fwd.hpp"
29 #include "../fiber/_3d_array.hpp"
30 
31 namespace Fiber
32 {
33 
35 template <typename CoordinateType> struct GeometricalData;
38 } // namespace Fiber
39 
40 namespace Bempp
41 {
42 
45 class Geometry
46 {
47 public:
49  virtual ~Geometry() {}
50 
52  virtual int dim() const = 0;
53 
55  virtual int dimWorld() const = 0;
56 
66  void setup(const arma::Mat<double>& corners,
67  const arma::Col<char>& auxData);
69  void setup(const arma::Mat<float>& corners,
70  const arma::Col<char>& auxData);
71 
76  virtual GeometryType type() const = 0;
77 
79  virtual bool affine() const = 0;
80 
88  virtual int cornerCount() const = 0;
89 
96  void getCorners(arma::Mat<double>& c) const;
98  void getCorners(arma::Mat<float>& c) const;
99 
107  void local2global(const arma::Mat<double>& local,
108  arma::Mat<double>& global) const;
110  void local2global(const arma::Mat<float>& local,
111  arma::Mat<float>& global) const;
112 
124  void global2local(const arma::Mat<double>& global,
125  arma::Mat<double>& local) const;
127  void global2local(const arma::Mat<float>& global,
128  arma::Mat<float>& local) const;
129 
156  void getIntegrationElements(const arma::Mat<double>& local,
157  arma::Row<double>& int_element) const;
159  void getIntegrationElements(const arma::Mat<float>& local,
160  arma::Row<float>& int_element) const;
161 
163  virtual double volume() const = 0;
164 
179  void getCenter(arma::Col<double>& c) const;
181  void getCenter(arma::Col<float>& c) const;
182 
196  const arma::Mat<double>& local,
197  arma::Cube<double>& jacobian_t) const;
200  const arma::Mat<float>& local,
201  arma::Cube<float>& jacobian_t) const;
204  const arma::Mat<double>& local,
205  Fiber::_3dArray<double>& jacobian_t) const;
208  const arma::Mat<float>& local,
209  Fiber::_3dArray<float>& jacobian_t) const;
210 
236  const arma::Mat<double>& local,
237  arma::Cube<double>& jacobian_inv_t) const;
240  const arma::Mat<float>& local,
241  arma::Cube<float>& jacobian_inv_t) const;
244  const arma::Mat<double>& local,
245  Fiber::_3dArray<double>& jacobian_inv_t) const;
248  const arma::Mat<float>& local,
249  Fiber::_3dArray<float>& jacobian_inv_t) const;
250 
262  void getNormals(const arma::Mat<double>& local,
263  arma::Mat<double>& normal) const;
265  void getNormals(const arma::Mat<float>& local,
266  arma::Mat<float>& normal) const;
267 
268  void getData(size_t what, const arma::Mat<double>& local,
269  Fiber::GeometricalData<double>& data) const;
271  void getData(size_t what, const arma::Mat<float>& local,
272  Fiber::GeometricalData<float>& data) const;
273 
274 private:
275  // Virtual functions to be implemented in subclasses
276  virtual void setupImpl(
277  const arma::Mat<double>& corners,
278  const arma::Col<char>& auxData) = 0;
279  virtual void getCornersImpl(arma::Mat<double>& c) const = 0;
280  virtual void local2globalImpl(
281  const arma::Mat<double>& local,
282  arma::Mat<double>& global) const = 0;
283  virtual void global2localImpl(
284  const arma::Mat<double>& global,
285  arma::Mat<double>& local) const = 0;
286  virtual void getIntegrationElementsImpl(
287  const arma::Mat<double>& local,
288  arma::Row<double>& int_element) const = 0;
289  virtual void getCenterImpl(arma::Col<double>& c) const = 0;
290  virtual void getJacobiansTransposedImpl(
291  const arma::Mat<double>& local,
292  Fiber::_3dArray<double>& jacobian_t) const = 0;
293  virtual void getJacobianInversesTransposedImpl(
294  const arma::Mat<double>& local,
295  Fiber::_3dArray<double>& jacobian_inv_t) const = 0;
296  virtual void getNormalsImpl(
297  const arma::Mat<double>& local,
298  arma::Mat<double>& normal) const = 0;
299  virtual void getDataImpl(
300  size_t what, const arma::Mat<double>& local,
301  Fiber::GeometricalData<double>& data) const = 0;
302 
303  // Helper functions for implementation
304  template <typename T1, typename T2>
305  void convertMat(const arma::Mat<T1>& in, arma::Mat<T2>& out) const;
306  template <typename T1, typename T2>
307  void convertCube(const arma::Cube<T1>& in, arma::Cube<T2>& out) const;
308  template <typename T1, typename T2>
309  void convertCube(const Fiber::_3dArray<T1>& in, Fiber::_3dArray<T2>& out) const;
310  template <typename T1, typename T2>
311  void convertCube(const Fiber::_3dArray<T1>& in, arma::Cube<T2>& out) const;
312 };
313 
314 } // namespace Bempp
315 
316 #include "geometry_imp.hpp"
317 
318 #endif
void local2global(const arma::Mat< double > &local, arma::Mat< double > &global) const
Convert local (logical) to global (physical) coordinates.
Definition: geometry_imp.hpp:60
virtual int dim() const =0
Dimension of the geometry.
void setup(const arma::Mat< double > &corners, const arma::Col< char > &auxData)
Set up geometry of an entity.
Definition: geometry_imp.hpp:34
Dune::GeometryType GeometryType
Identifier of geometry type.
Definition: geometry_type.hpp:34
void getJacobiansTransposed(const arma::Mat< double > &local, arma::Cube< double > &jacobian_t) const
Get transposed Jacobian matrices at specified points.
Definition: geometry_imp.hpp:118
virtual bool affine() const =0
True if the geometry mapping is affine and false otherwise.
Storage of geometrical data.
Definition: geometrical_data.hpp:54
virtual int cornerCount() const =0
Number of corners of the reference element.
virtual GeometryType type() const =0
Type of the reference element.
virtual ~Geometry()
Destructor.
Definition: geometry.hpp:49
void global2local(const arma::Mat< double > &global, arma::Mat< double > &local) const
Convert global (physical) to local (logical) coordinates.
Definition: geometry_imp.hpp:75
void getIntegrationElements(const arma::Mat< double > &local, arma::Row< double > &int_element) const
Get the factor appearing in the integral transformation formula at specified points.
Definition: geometry_imp.hpp:90
void getCenter(arma::Col< double > &c) const
Get center of geometry.
Definition: geometry_imp.hpp:106
void getJacobianInversesTransposed(const arma::Mat< double > &local, arma::Cube< double > &jacobian_inv_t) const
Get inverses of transposed Jacobian matrices at specified points.
Definition: geometry_imp.hpp:161
virtual double volume() const =0
Volume of geometry.
virtual int dimWorld() const =0
Dimension of the space containing the geometry.
void getNormals(const arma::Mat< double > &local, arma::Mat< double > &normal) const
Get unit vectors normal to the entity at specified points.
Definition: geometry_imp.hpp:204
Simple implementation of a 3D Fortran-ordered array.
Definition: _3d_array.hpp:48
void getCorners(arma::Mat< double > &c) const
Get the positions of the geometry corners.
Definition: geometry_imp.hpp:48
Abstract wrapper of a geometry.
Definition: geometry.hpp:45