21 #ifndef bempp_geometry_hpp
22 #define bempp_geometry_hpp
24 #include "../common/common.hpp"
28 #include "../common/armadillo_fwd.hpp"
29 #include "../fiber/_3d_array.hpp"
35 template <
typename CoordinateType>
struct GeometricalData;
52 virtual int dim()
const = 0;
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);
79 virtual bool affine()
const = 0;
108 arma::Mat<double>& global)
const;
111 arma::Mat<float>& global)
const;
125 arma::Mat<double>& local)
const;
128 arma::Mat<float>& local)
const;
157 arma::Row<double>& int_element)
const;
160 arma::Row<float>& int_element)
const;
163 virtual double volume()
const = 0;
179 void getCenter(arma::Col<double>& c)
const;
181 void getCenter(arma::Col<float>& c)
const;
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,
208 const arma::Mat<float>& local,
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,
248 const arma::Mat<float>& local,
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;
268 void getData(
size_t what,
const arma::Mat<double>& local,
271 void getData(
size_t what,
const arma::Mat<float>& local,
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,
293 virtual void getJacobianInversesTransposedImpl(
294 const arma::Mat<double>& local,
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,
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>
310 template <
typename T1,
typename T2>
316 #include "geometry_imp.hpp"
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