BEM++
2.0
|
Abstract wrapper of a geometry. More...
#include </home/wojtek/Projects/BEM/bempp-sven/bempp/lib/grid/geometry.hpp>
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 |
Abstract wrapper of a geometry.
|
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 >.
|
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.
[out] | c | Coordinates 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().
|
inline |
This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
|
inline |
Get the positions of the geometry corners.
[out] | c | Matrix whose ![]() ![]() |
Referenced by Bempp::Grid::areInside(), Bempp::UnitScalarSpace< BasisFunctionType >::getFlatLocalDofBoundingBoxes(), Bempp::PiecewiseConstantDualGridDiscontinuousScalarSpace< BasisFunctionType >::getFlatLocalDofBoundingBoxes(), Bempp::PiecewiseLinearContinuousScalarSpace< BasisFunctionType >::getFlatLocalDofBoundingBoxes(), Bempp::PiecewiseConstantDualGridScalarSpace< BasisFunctionType >::getFlatLocalDofBoundingBoxes(), Bempp::PiecewiseLinearDiscontinuousScalarSpaceBarycentric< BasisFunctionType >::getFlatLocalDofBoundingBoxes(), Bempp::PiecewiseLinearContinuousScalarSpaceBarycentric< BasisFunctionType >::getFlatLocalDofBoundingBoxes(), Bempp::RaviartThomas0VectorSpace< BasisFunctionType >::getFlatLocalDofBoundingBoxes(), Bempp::UnitScalarSpace< BasisFunctionType >::getGlobalDofBoundingBoxes(), and Bempp::PiecewiseConstantScalarSpace< BasisFunctionType >::getGlobalDofBoundingBoxes().
|
inline |
This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
|
inline |
This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
|
inline |
Get the factor appearing in the integral transformation formula at specified points.
Let denote the transformation described by the Geometry. Then the jacobian of the transformation is defined as the
matrix
Here we abbreviated and
for readability.
The integration element for any
is then defined as
[in] | local | Matrix whose ![]() ![]() |
[out] | int_element | Row vector whose ![]() ![]() |
|
inline |
This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
|
inline |
Get inverses of transposed Jacobian matrices at specified points.
The Jacobian matrix is defined in the documentation of getIntegrationElements().
[in] | local | Matrix whose ![]() ![]() |
[out] | jacobian_inv_t | 3D array whose ![]() ![]() ![]() |
The use of this function is to compute the gradient of some function at some position
, where
and
the transformation of the Geometry. When we set
and apply the chain rule we obtain
References dim(), and dimWorld().
Referenced by getJacobianInversesTransposed().
|
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().
|
inline |
This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
|
inline |
This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
|
inline |
Get transposed Jacobian matrices at specified points.
The Jacobian matrix is defined in the documentation of getIntegrationElements().
[in] | local | Matrix whose ![]() ![]() |
[out] | jacobian_t | 3D array whose ![]() ![]() ![]() |
References dim(), and dimWorld().
|
inline |
This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
|
inline |
This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
|
inline |
This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
|
inline |
Get unit vectors normal to the entity at specified points.
An exception is thrown if dim() != dimWorld() - 1.
[in] | local | Matrix whose ![]() ![]() |
[out] | normal | Matrix whose ![]() ![]() |
Referenced by Bempp::PiecewiseConstantDualGridDiscontinuousScalarSpace< BasisFunctionType >::getFlatLocalDofNormals(), Bempp::PiecewiseLinearContinuousScalarSpace< BasisFunctionType >::getFlatLocalDofNormals(), Bempp::PiecewiseConstantDualGridScalarSpace< BasisFunctionType >::getFlatLocalDofNormals(), Bempp::PiecewiseConstantScalarSpaceBarycentric< BasisFunctionType >::getFlatLocalDofNormals(), Bempp::PiecewiseConstantDiscontinuousScalarSpaceBarycentric< BasisFunctionType >::getFlatLocalDofNormals(), Bempp::PiecewiseLinearDiscontinuousScalarSpaceBarycentric< BasisFunctionType >::getFlatLocalDofNormals(), Bempp::PiecewiseLinearContinuousScalarSpaceBarycentric< BasisFunctionType >::getFlatLocalDofNormals(), Bempp::RaviartThomas0VectorSpace< BasisFunctionType >::getFlatLocalDofNormals(), Bempp::PiecewiseLinearDiscontinuousScalarSpaceBarycentric< BasisFunctionType >::getGlobalDofNormals(), Bempp::PiecewiseLinearContinuousScalarSpaceBarycentric< BasisFunctionType >::getGlobalDofNormals(), and Bempp::RaviartThomas0VectorSpace< BasisFunctionType >::getGlobalDofNormals().
|
inline |
This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
|
inline |
Convert global (physical) to local (logical) coordinates.
[in] | global | Matrix whose ![]() ![]() |
[out] | local | Matrix whose ![]() ![]() ![]() |
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.
|
inline |
This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
|
inline |
Convert local (logical) to global (physical) coordinates.
[in] | local | Matrix whose ![]() ![]() |
[out] | global | Matrix whose ![]() ![]() ![]() |
Referenced by Bempp::GridFunction< BasisFunctionType, ResultType >::exportToGmsh().
|
inline |
This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
|
inline |
Set up geometry of an entity.
[in] | corners | Coordinates of the entity's vertices, stored columnwise. |
[in] | auxData | Auxiliary 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. |
|
inline |
This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
|
pure virtual |
Type of the reference element.
The type can be used to access the Dune::GenericReferenceElement
.
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 >.