BEM++  2.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Groups Pages
Public Member Functions | Private Member Functions | List of all members
Bempp::Geometry Class Referenceabstract

Abstract wrapper of a geometry. More...

#include </home/wojtek/Projects/BEM/bempp-sven/bempp/lib/grid/geometry.hpp>

Inheritance diagram for Bempp::Geometry:
Bempp::ConcreteGeometry< DuneGeometry > Bempp::ConcreteGeometry< typename DuneEntity::Bempp::Geometry > Bempp::ConcreteGeometry< typename DuneSubentity::Bempp::Geometry > Bempp::ConcreteGeometry< typename typename DuneEntityIt::Bempp::Entity::Bempp::Geometry >

Public Member Functions

virtual ~Geometry ()
 Destructor.
 
virtual int dim () const =0
 Dimension of the geometry.
 
virtual int dimWorld () const =0
 Dimension of the space containing the geometry.
 
void setup (const arma::Mat< double > &corners, const arma::Col< char > &auxData)
 Set up geometry of an entity. More...
 
void setup (const arma::Mat< float > &corners, const arma::Col< char > &auxData)
 
virtual GeometryType type () const =0
 Type of the reference element. More...
 
virtual bool affine () const =0
 True if the geometry mapping is affine and false otherwise.
 
virtual int cornerCount () const =0
 Number of corners of the reference element. More...
 
void getCorners (arma::Mat< double > &c) const
 Get the positions of the geometry corners. More...
 
void getCorners (arma::Mat< float > &c) const
 
void local2global (const arma::Mat< double > &local, arma::Mat< double > &global) const
 Convert local (logical) to global (physical) coordinates. More...
 
void local2global (const arma::Mat< float > &local, arma::Mat< float > &global) const
 
void global2local (const arma::Mat< double > &global, arma::Mat< double > &local) const
 Convert global (physical) to local (logical) coordinates. More...
 
void global2local (const arma::Mat< float > &global, arma::Mat< float > &local) const
 
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. More...
 
void getIntegrationElements (const arma::Mat< float > &local, arma::Row< float > &int_element) const
 
virtual double volume () const =0
 Volume of geometry.
 
void getCenter (arma::Col< double > &c) const
 Get center of geometry. More...
 
void getCenter (arma::Col< float > &c) const
 
void getJacobiansTransposed (const arma::Mat< double > &local, arma::Cube< double > &jacobian_t) const
 Get transposed Jacobian matrices at specified points. More...
 
void getJacobiansTransposed (const arma::Mat< float > &local, arma::Cube< float > &jacobian_t) const
 
void getJacobiansTransposed (const arma::Mat< double > &local, Fiber::_3dArray< double > &jacobian_t) const
 
void getJacobiansTransposed (const arma::Mat< float > &local, Fiber::_3dArray< float > &jacobian_t) const
 
void getJacobianInversesTransposed (const arma::Mat< double > &local, arma::Cube< double > &jacobian_inv_t) const
 Get inverses of transposed Jacobian matrices at specified points. More...
 
void getJacobianInversesTransposed (const arma::Mat< float > &local, arma::Cube< float > &jacobian_inv_t) const
 
void getJacobianInversesTransposed (const arma::Mat< double > &local, Fiber::_3dArray< double > &jacobian_inv_t) const
 
void getJacobianInversesTransposed (const arma::Mat< float > &local, Fiber::_3dArray< float > &jacobian_inv_t) const
 
void getNormals (const arma::Mat< double > &local, arma::Mat< double > &normal) const
 Get unit vectors normal to the entity at specified points. More...
 
void getNormals (const arma::Mat< float > &local, arma::Mat< float > &normal) const
 
void getData (size_t what, const arma::Mat< double > &local, Fiber::GeometricalData< double > &data) const
 
void getData (size_t what, const arma::Mat< float > &local, Fiber::GeometricalData< float > &data) const
 

Private Member Functions

virtual void setupImpl (const arma::Mat< double > &corners, const arma::Col< char > &auxData)=0
 
virtual void getCornersImpl (arma::Mat< double > &c) const =0
 
virtual void local2globalImpl (const arma::Mat< double > &local, arma::Mat< double > &global) const =0
 
virtual void global2localImpl (const arma::Mat< double > &global, arma::Mat< double > &local) const =0
 
virtual void getIntegrationElementsImpl (const arma::Mat< double > &local, arma::Row< double > &int_element) const =0
 
virtual void getCenterImpl (arma::Col< double > &c) const =0
 
virtual void getJacobiansTransposedImpl (const arma::Mat< double > &local, Fiber::_3dArray< double > &jacobian_t) const =0
 
virtual void getJacobianInversesTransposedImpl (const arma::Mat< double > &local, Fiber::_3dArray< double > &jacobian_inv_t) const =0
 
virtual void getNormalsImpl (const arma::Mat< double > &local, arma::Mat< double > &normal) const =0
 
virtual void getDataImpl (size_t what, const arma::Mat< double > &local, Fiber::GeometricalData< double > &data) const =0
 
template<typename T1 , typename T2 >
void convertMat (const arma::Mat< T1 > &in, arma::Mat< T2 > &out) const
 
template<typename T1 , typename T2 >
void convertCube (const arma::Cube< T1 > &in, arma::Cube< T2 > &out) const
 
template<typename T1 , typename T2 >
void convertCube (const Fiber::_3dArray< T1 > &in, Fiber::_3dArray< T2 > &out) const
 
template<typename T1 , typename T2 >
void convertCube (const Fiber::_3dArray< T1 > &in, arma::Cube< T2 > &out) const
 

Detailed Description

Abstract wrapper of a geometry.

Member Function Documentation

virtual int Bempp::Geometry::cornerCount ( ) const
pure virtual

Number of corners of the reference element.

Since a geometry is a convex polytope the number of corners is a well-defined concept. The method is redundant because this information is also available via the reference element. It is here for efficiency and ease of use.

Implemented in Bempp::ConcreteGeometry< DuneGeometry >, Bempp::ConcreteGeometry< typename typename DuneEntityIt::Bempp::Entity::Bempp::Geometry >, Bempp::ConcreteGeometry< typename DuneEntity::Bempp::Geometry >, and Bempp::ConcreteGeometry< typename DuneSubentity::Bempp::Geometry >.

void Bempp::Geometry::getCenter ( arma::Col< double > &  c) const
inline

Get center of geometry.

Note that this method is still subject to a change of name and semantics. At the moment, the center is not required to be the centroid of the geometry, or even the centroid of its corners. This makes acceptable the current default implementation, which maps the centroid of the reference element to the geometry.

We may change the name (and semantic) of the method to centroid() if Dune's developers find reasonably efficient ways to implement it properly.

Parameters
[out]cCoordinates of the center of geometry.

Referenced by Bempp::ElementaryPotentialOperator< BasisFunctionType_, KernelType_, ResultType_ >::evaluateOnGrid(), Bempp::UnitScalarSpace< BasisFunctionType >::getFlatLocalDofBoundingBoxes(), Bempp::UnitScalarSpace< BasisFunctionType >::getFlatLocalDofPositions(), Bempp::UnitScalarSpace< BasisFunctionType >::getGlobalDofBoundingBoxes(), Bempp::PiecewiseConstantScalarSpace< BasisFunctionType >::getGlobalDofBoundingBoxes(), and Bempp::UnitScalarSpace< BasisFunctionType >::getGlobalDofPositions().

void Bempp::Geometry::getCenter ( arma::Col< float > &  c) const
inline

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

void Bempp::Geometry::getCorners ( arma::Mat< double > &  c) const
inline
void Bempp::Geometry::getCorners ( arma::Mat< float > &  c) const
inline

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

void Bempp::Geometry::getData ( size_t  what,
const arma::Mat< float > &  local,
Fiber::GeometricalData< float > &  data 
) const
inline

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

void Bempp::Geometry::getIntegrationElements ( const arma::Mat< double > &  local,
arma::Row< double > &  int_element 
) const
inline

Get the factor appearing in the integral transformation formula at specified points.

Let $ g : D \to W$ denote the transformation described by the Geometry. Then the jacobian of the transformation is defined as the $\textrm{cdim}\times\textrm{mydim}$ matrix

\[ J_g(x) = \left( \begin{array}{ccc} \frac{\partial g_0}{\partial x_0} & \cdots & \frac{\partial g_0}{\partial x_{n-1}} \\ \vdots & \ddots & \vdots \\ \frac{\partial g_{m-1}}{\partial x_0} & \cdots & \frac{\partial g_{m-1}}{\partial x_{n-1}} \end{array} \right).\]

Here we abbreviated $m=\textrm{cdim}$ and $n=\textrm{mydim}$ for readability.

The integration element $\mu(x)$ for any $x\in D$ is then defined as

\[ \mu(x) = \sqrt{|\det J_g^T(x)J_g(x)|}.\]

Parameters
[in]localMatrix whose $i$th column contains the local coordinates of a point $x_i \in D$.
[out]int_elementRow vector whose $i$th entry contains the integration element $\mu(x_i)$.
Note
Each implementation computes the integration element with optimal efficiency. For example in an equidistant structured mesh it may be as simple as $h^\textrm{mydim}$.
void Bempp::Geometry::getIntegrationElements ( const arma::Mat< float > &  local,
arma::Row< float > &  int_element 
) const
inline

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

void Bempp::Geometry::getJacobianInversesTransposed ( const arma::Mat< double > &  local,
arma::Cube< double > &  jacobian_inv_t 
) const
inline

Get inverses of transposed Jacobian matrices at specified points.

The Jacobian matrix is defined in the documentation of getIntegrationElements().

Parameters
[in]localMatrix whose $i$th column contains the local coordinates of a point $x_i \in D$.
[out]jacobian_inv_t3D array whose $i$th slice (i.e. jacobian_inv_t(:,:,i)) contains the inverse of the transposed Jacobian matrix at $x_i$, i.e. $J_g^{-T}(x_i)$.

The use of this function is to compute the gradient of some function $f : W \to \textbf{R}$ at some position $y=g(x)$, where $x\in D$ and $g$ the transformation of the Geometry. When we set $\hat{f}(x) = f(g(x))$ and apply the chain rule we obtain

\[\nabla f(g(x)) = J_g^{-T}(x) \nabla \hat{f}(x).\]

Note
In the non-symmetric case $\textrm{cdim} \neq \textrm{mydim}$, the pseudoinverse of $J_g^T(x)$ is returned. This means that it is inverse for all tangential vectors in $g(x)$ while mapping all normal vectors to zero.

References dim(), and dimWorld().

Referenced by getJacobianInversesTransposed().

void Bempp::Geometry::getJacobianInversesTransposed ( const arma::Mat< float > &  local,
arma::Cube< float > &  jacobian_inv_t 
) const
inline

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

References getJacobianInversesTransposed().

void Bempp::Geometry::getJacobianInversesTransposed ( const arma::Mat< double > &  local,
Fiber::_3dArray< double > &  jacobian_inv_t 
) const
inline

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

void Bempp::Geometry::getJacobianInversesTransposed ( const arma::Mat< float > &  local,
Fiber::_3dArray< float > &  jacobian_inv_t 
) const
inline

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

void Bempp::Geometry::getJacobiansTransposed ( const arma::Mat< double > &  local,
arma::Cube< double > &  jacobian_t 
) const
inline

Get transposed Jacobian matrices at specified points.

The Jacobian matrix is defined in the documentation of getIntegrationElements().

Parameters
[in]localMatrix whose $i$th column contains the local coordinates of a point $x_i \in D$.
[out]jacobian_t3D array whose $i$th slice (i.e. jacobian_t(:,:,i)) contains the transposed Jacobian matrix at $x_i$, i.e. $J_g^T(x_i)$.

References dim(), and dimWorld().

void Bempp::Geometry::getJacobiansTransposed ( const arma::Mat< float > &  local,
arma::Cube< float > &  jacobian_t 
) const
inline

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

void Bempp::Geometry::getJacobiansTransposed ( const arma::Mat< double > &  local,
Fiber::_3dArray< double > &  jacobian_t 
) const
inline

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

void Bempp::Geometry::getJacobiansTransposed ( const arma::Mat< float > &  local,
Fiber::_3dArray< float > &  jacobian_t 
) const
inline

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

void Bempp::Geometry::getNormals ( const arma::Mat< double > &  local,
arma::Mat< double > &  normal 
) const
inline
void Bempp::Geometry::getNormals ( const arma::Mat< float > &  local,
arma::Mat< float > &  normal 
) const
inline

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

void Bempp::Geometry::global2local ( const arma::Mat< double > &  global,
arma::Mat< double > &  local 
) const
inline

Convert global (physical) to local (logical) coordinates.

Parameters
[in]globalMatrix whose $i$th column contains the global coordinates of a point $x_i \in W$.
[out]localMatrix whose $i$th column contains the local coordinates of $x_i$, i.e. $g^{-1}(x_i)$.

This is going to be tricky to implement for dimGrid < dimWorld. Maybe the docstring should say that we convert some sort of projection of global to local.

void Bempp::Geometry::global2local ( const arma::Mat< float > &  global,
arma::Mat< float > &  local 
) const
inline

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

void Bempp::Geometry::local2global ( const arma::Mat< double > &  local,
arma::Mat< double > &  global 
) const
inline

Convert local (logical) to global (physical) coordinates.

Parameters
[in]localMatrix whose $i$th column contains the local coordinates of a point $x_i \in D$.
[out]globalMatrix whose $i$th column contains the global coordinates of $x_i$, i.e. $g(x_i)$.

Referenced by Bempp::GridFunction< BasisFunctionType, ResultType >::exportToGmsh().

void Bempp::Geometry::local2global ( const arma::Mat< float > &  local,
arma::Mat< float > &  global 
) const
inline

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

void Bempp::Geometry::setup ( const arma::Mat< double > &  corners,
const arma::Col< char > &  auxData 
)
inline

Set up geometry of an entity.

Parameters
[in]cornersCoordinates of the entity's vertices, stored columnwise.
[in]auxDataAuxiliary data necessary for the description of the entity. Interpretation of these data in subclasses of Geometry may vary. They can be used for example to define a curvilinear element.
void Bempp::Geometry::setup ( const arma::Mat< float > &  corners,
const arma::Col< char > &  auxData 
)
inline

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

virtual GeometryType Bempp::Geometry::type ( ) const
pure virtual

The documentation for this class was generated from the following files: