LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
Loading...
Searching...
No Matches
Public Types | Public Member Functions | Protected Member Functions | Related Symbols | List of all members
lf::geometry::Geometry Class Referenceabstract

Interface class for shape information on a mesh cell in the spirit of parametric finite element methods More...

#include <lf/geometry/geometry_interface.h>

Inheritance diagram for lf::geometry::Geometry:
lf::geometry::Parallelogram lf::geometry::Point lf::geometry::QuadO1 lf::geometry::QuadO2 lf::geometry::SegmentO1 lf::geometry::SegmentO2 lf::geometry::TriaO1 lf::geometry::TriaO2

Public Types

using dim_t = base::RefEl::dim_t
 

Public Member Functions

virtual dim_t DimLocal () const =0
 Dimension of the domain of this mapping.
 
virtual dim_t DimGlobal () const =0
 Dimension of the image of this mapping.
 
virtual base::RefEl RefEl () const =0
 The Reference element that defines the domain of this mapping.
 
virtual Eigen::MatrixXd Global (const Eigen::MatrixXd &local) const =0
 Map a number of points in local coordinates into the global coordinate system.
 
virtual Eigen::MatrixXd Jacobian (const Eigen::MatrixXd &local) const =0
 Evaluate the jacobian of the mapping simultaneously at numPoints points.
 
virtual Eigen::MatrixXd JacobianInverseGramian (const Eigen::MatrixXd &local) const =0
 Evaluate the Jacobian * Inverse Gramian ( \( J (J^T J)^{-1}\)) simultaneously at numPoints.
 
virtual Eigen::VectorXd IntegrationElement (const Eigen::MatrixXd &local) const =0
 The integration element (factor appearing in integral transformation formula, see below) at number of evaluation points (specified in local coordinates).
 
virtual std::unique_ptr< GeometrySubGeometry (dim_t codim, dim_t i) const =0
 Construct a new Geometry() object that describes the geometry of the i-th sub-entity with codimension=codim
 
virtual std::vector< std::unique_ptr< Geometry > > ChildGeometry (const RefinementPattern &ref_pat, lf::base::dim_t codim) const =0
 Generate geometry objects for child entities created in the course of refinement.
 
virtual bool isAffine () const
 element shape by affine mapping from reference element
 
virtual ~Geometry ()=default
 Virtual destructor.
 

Protected Member Functions

 Geometry ()=default
 
 Geometry (const Geometry &)=default
 
 Geometry (Geometry &&)=default
 
Geometryoperator= (const Geometry &)=default
 
Geometryoperator= (Geometry &&)=default
 

Related Symbols

(Note that these are not member symbols.)

void PrintInfo (std::ostream &o, const Geometry &geom, int output_ctrl=0)
 Diagnostic output operator for Geometry.
 
std::ostream & operator<< (std::ostream &stream, const Geometry &geom)
 Operator overload to print a Geometry to a stream, such as std::cout
 

Detailed Description

Interface class for shape information on a mesh cell in the spirit of parametric finite element methods

This abstract base class is supposed to provide complete information about the mapping \(\Phi:\widehat{K}\to K\) that takes the reference cell to a particular cell. The reference cell is defined through the topological type of a cell, see lf::base::RefEl.

Definition at line 21 of file geometry_interface.h.

Member Typedef Documentation

◆ dim_t

Definition at line 30 of file geometry_interface.h.

Constructor & Destructor Documentation

◆ Geometry() [1/3]

lf::geometry::Geometry::Geometry ( )
protecteddefault

◆ Geometry() [2/3]

lf::geometry::Geometry::Geometry ( const Geometry & )
protecteddefault

◆ Geometry() [3/3]

lf::geometry::Geometry::Geometry ( Geometry && )
protecteddefault

◆ ~Geometry()

virtual lf::geometry::Geometry::~Geometry ( )
virtualdefault

Virtual destructor.

Member Function Documentation

◆ ChildGeometry()

virtual std::vector< std::unique_ptr< Geometry > > lf::geometry::Geometry::ChildGeometry ( const RefinementPattern & ref_pat,
lf::base::dim_t codim ) const
pure virtual

Generate geometry objects for child entities created in the course of refinement.

Parameters
ref_patAn object encoding the details of refinement
codimrelative codimension of child entities whose shape is requested
Returns
an array of unique pointers pointers to geometry objects for child entities of the specified co-dimension and the given refinement pattern. The numbering of child entities is a convention.

Implementation of local mesh refinement entails splitting of elements into parts ("children"), whose shape will depend on the refinement pattern. This method creates the geometry objects describing the shape of children. The details of subdivisions corresponding to particular refinement patterns are fixed by the method RefinementPattern::ChildPolygons() and should be documented there.

See also
class lf::geometry::RefinementPattern

Implemented in lf::geometry::SegmentO1, lf::geometry::SegmentO2, lf::geometry::QuadO1, lf::geometry::Parallelogram, lf::geometry::QuadO2, lf::geometry::TriaO1, lf::geometry::TriaO2, and lf::geometry::Point.

Referenced by lf::geometry::test_utils::checkChildGeometry(), lf::geometry::test_utils::checkChildGeometryVolume(), and lf::refinement::MeshHierarchy::PerformRefinement().

◆ DimGlobal()

virtual dim_t lf::geometry::Geometry::DimGlobal ( ) const
pure virtual

◆ DimLocal()

virtual dim_t lf::geometry::Geometry::DimLocal ( ) const
pure virtual

◆ Global()

virtual Eigen::MatrixXd lf::geometry::Geometry::Global ( const Eigen::MatrixXd & local) const
pure virtual

Map a number of points in local coordinates into the global coordinate system.

Parameters
localA Matrix of size DimLocal() x numPoints that contains in its columns the coordinates of the points at which the mapping function should be evaluated.
Returns
A Matrix of size DimGlobal() x numPoints that contains the mapped points as column vectors. Here numPoints is the number of columns of the matrix passed in the local argument.

\[ \mathtt{Global}\left(\left[\widehat{x}^1,\ldots,\widehat{x}^k\right]\right) = \left[ \mathbf{\Phi}_K(\widehat{x}^1),\ldots,\mathbf{\Phi}_K(\widehat{x}^k)\right]\;,\quad \widehat{x}^{\ell}\in\widehat{K}\;, \]

where \(\mathbf{\Phi}\) is the mapping taking the reference element to the current element \(K\).

This method provides a complete description of the shape of an entity through a parameterization over the corresponding reference element = parameter domain. The method takes as arguments a number of coordinate vectors of points in the reference element. For the sake of efficiency, these coordinate vectors are passed as the columns of a dynamic matrix type as supplied by Eigen.

For instance, this method is used in lf::geometry::Corners()

inline Eigen::MatrixXd Corners(const Geometry& geo) {
return geo.Global(geo.RefEl().NodeCoords()); }
const Eigen::MatrixXd & NodeCoords() const
Get the coordinates of the nodes of this reference element.
Definition ref_el.h:241
Interface class for shape information on a mesh cell in the spirit of parametric finite element metho...
virtual base::RefEl RefEl() const =0
The Reference element that defines the domain of this mapping.
virtual Eigen::MatrixXd Global(const Eigen::MatrixXd &local) const =0
Map a number of points in local coordinates into the global coordinate system.
Eigen::MatrixXd Corners(const Geometry &geo)
The corners of a shape with piecewise smooth boundary.

Additional explanations in Lecture Document Paragraph 2.7.5.17.

Implemented in lf::geometry::Point, lf::geometry::QuadO1, lf::geometry::Parallelogram, lf::geometry::QuadO2, lf::geometry::SegmentO1, lf::geometry::SegmentO2, lf::geometry::TriaO1, and lf::geometry::TriaO2.

Referenced by lf::geometry::test_utils::checkChildGeometry(), lf::mesh::test_utils::checkGeometryOrientation(), lf::geometry::test_utils::checkJacobian(), lf::geometry::test_utils::checkSubGeometry(), lf::geometry::Corners(), lf::uscalfe::ReactionDiffusionElementMatrixProvider< SCALAR, DIFF_COEFF, REACTION_COEFF >::Eval(), lf::uscalfe::ScalarLoadElementVectorProvider< SCALAR, MESH_FUNCTION >::Eval(), lf::fe::DiffusionElementMatrixProvider< SCALAR, DIFF_COEFF >::Eval(), lf::fe::MassElementMatrixProvider< SCALAR, REACTION_COEFF >::Eval(), lf::fe::ScalarLoadElementVectorProvider< SCALAR, MESH_FUNCTION >::Eval(), lf::uscalfe::LinearFELaplaceElementMatrix::Eval(), lf::uscalfe::LinearFELocalLoadVector< SCALAR, FUNCTOR >::Eval(), lf::uscalfe::ScalarLoadEdgeVectorProvider< SCALAR, FUNCTOR, EDGESELECTOR >::Eval(), lf::refinement::EntityCenterPositionSelector< POSPRED >::operator()(), lf::mesh::utils::MeshFunctionGlobal< F >::operator()(), PrintInfo(), lf::mesh::Mesh::PrintInfo(), lf::mesh::Entity::PrintInfo(), lf::io::writeMatlab(), lf::io::writeMatplotlib(), and lf::io::writeTikZ().

◆ IntegrationElement()

virtual Eigen::VectorXd lf::geometry::Geometry::IntegrationElement ( const Eigen::MatrixXd & local) const
pure virtual

The integration element (factor appearing in integral transformation formula, see below) at number of evaluation points (specified in local coordinates).

Parameters
localA Matrix of size DimLocal() x numPoints that contains the evaluation points (in local = reference coordinates) as column vectors.
Returns
A Vector of size numPoints x 1 that contains the integration elements at every evaluation point.

For a transformation \( \Phi : K \mapsto R^{\text{DimGlobal}}\) with Jacobian \( D\Phi : K \mapsto R^{\text{DimGlobal} \times \text{DimLocal}} \) the integration element \( g \) at point \( \xi \in K \) is defined by

\[ g(\xi) := \sqrt{\mathrm{det}\left|D\Phi^T(\xi) D\Phi(\xi) \right|} \]

More information also related to the use of lovcal quadrature rules is given in Lecture Document Paragraph 2.7.5.24.

Implemented in lf::geometry::Point, lf::geometry::QuadO1, lf::geometry::Parallelogram, lf::geometry::QuadO2, lf::geometry::SegmentO1, lf::geometry::SegmentO2, lf::geometry::TriaO1, and lf::geometry::TriaO2.

Referenced by lf::geometry::test_utils::checkChildGeometryVolume(), lf::geometry::test_utils::checkIntegrationElement(), lf::uscalfe::ReactionDiffusionElementMatrixProvider< SCALAR, DIFF_COEFF, REACTION_COEFF >::Eval(), lf::uscalfe::ScalarLoadElementVectorProvider< SCALAR, MESH_FUNCTION >::Eval(), lf::fe::DiffusionElementMatrixProvider< SCALAR, DIFF_COEFF >::Eval(), lf::fe::MassElementMatrixProvider< SCALAR, REACTION_COEFF >::Eval(), lf::fe::ScalarLoadElementVectorProvider< SCALAR, MESH_FUNCTION >::Eval(), lf::uscalfe::LinearFELaplaceElementMatrix::Eval(), lf::uscalfe::MassEdgeMatrixProvider< SCALAR, COEFF, EDGESELECTOR >::Eval(), lf::uscalfe::ScalarLoadEdgeVectorProvider< SCALAR, FUNCTOR, EDGESELECTOR >::Eval(), lf::fe::MassEdgeMatrixProvider< SCALAR, COEFF, EDGESELECTOR >::Eval(), lf::fe::ScalarLoadEdgeVectorProvider< SCALAR, FUNCTOR, EDGESELECTOR >::Eval(), and lf::geometry::Volume().

◆ isAffine()

virtual bool lf::geometry::Geometry::isAffine ( ) const
inlinevirtual

element shape by affine mapping from reference element

Returns
true, if the element is the affine image of a reference element

An affine map is a linear mapping plus a translation. For a 2D mesh an entity has an affine geometry, if

  • a segement is a straight line
  • a triangle is flat
  • a quadrilateral is a flat parallelogram.
Note
The Jacobian/Gramian for an affine element are constant.

Reimplemented in lf::geometry::QuadO1, and lf::geometry::Parallelogram.

Definition at line 220 of file geometry_interface.h.

◆ Jacobian()

virtual Eigen::MatrixXd lf::geometry::Geometry::Jacobian ( const Eigen::MatrixXd & local) const
pure virtual

Evaluate the jacobian of the mapping simultaneously at numPoints points.

Parameters
localA Matrix of size DimLocal x numPoints that contains the evaluation points as column vectors
Returns
A Matrix of size DimGlobal() x (DimLocal() * numPoints) that contains the jacobians at the evaluation points.

This method allows access to the derivative of the parametrization mapping in a number of points, passed as the columns of a dynamic matrix. The derivative of the parametrization in a point is a Jacobian matrix of size ‘DimGlobal() x DimLocal()’. For the sake of efficiency, these matrices are stacked horizontally and returned as one big dynamic matrix. Use Eigen's ‘block()’ method of Eigen::MatrixXd to extract the individual Jacobians from the returned matrix.

Implemented in lf::geometry::Point, lf::geometry::QuadO1, lf::geometry::Parallelogram, lf::geometry::QuadO2, lf::geometry::SegmentO1, lf::geometry::SegmentO2, lf::geometry::TriaO1, and lf::geometry::TriaO2.

Referenced by lf::geometry::test_utils::checkIntegrationElement(), lf::geometry::test_utils::checkJacobian(), and lf::geometry::test_utils::checkJacobianInverseGramian().

◆ JacobianInverseGramian()

virtual Eigen::MatrixXd lf::geometry::Geometry::JacobianInverseGramian ( const Eigen::MatrixXd & local) const
pure virtual

Evaluate the Jacobian * Inverse Gramian ( \( J (J^T J)^{-1}\)) simultaneously at numPoints.

Parameters
localA Matrix of size DimLocal() x numPoints that contains the evaluation points as column vectors.
Returns
A Matrix of size DimGlobal() x (DimLocal() * numPoints) that contains the Jacobian multiplied with the Inverse Gramian ( \( J (J^T J)^{-1}\)) at every evaluation point.
Note
If DimLocal() == DimGlobal() then \( J (J^T J)^{-1} = J^{-T} \), i.e. this method returns the inverse of the transposed jacobian.

Example for recovering a single transposed inverse Jacobian.

If both dimensions agree and have the value D, then the method returns the transposed of the inverse Jacobians of the transformation at the passed points. These are square DxD matrices.

To retrieve the j-th inverse transposed Jacobian from the returned matrix, use the block methdod of Eigen (case (D = DimLocal()) == DimGlobal())

JacobianInverseGramian(local).block(0,i*D,D,D)
virtual Eigen::MatrixXd JacobianInverseGramian(const Eigen::MatrixXd &local) const =0
Evaluate the Jacobian * Inverse Gramian ( ) simultaneously at numPoints.

More explanations in Paragraph 2.8.3.14.

Implemented in lf::geometry::Point, lf::geometry::QuadO1, lf::geometry::Parallelogram, lf::geometry::QuadO2, lf::geometry::SegmentO1, lf::geometry::SegmentO2, lf::geometry::TriaO1, and lf::geometry::TriaO2.

Referenced by lf::geometry::test_utils::checkJacobianInverseGramian(), lf::uscalfe::ReactionDiffusionElementMatrixProvider< SCALAR, DIFF_COEFF, REACTION_COEFF >::Eval(), lf::fe::DiffusionElementMatrixProvider< SCALAR, DIFF_COEFF >::Eval(), lf::uscalfe::LinearFELaplaceElementMatrix::Eval(), and lf::fe::MeshFunctionGradFE< SCALAR_FE, SCALAR_COEFF >::operator()().

◆ operator=() [1/2]

Geometry & lf::geometry::Geometry::operator= ( const Geometry & )
protecteddefault

◆ operator=() [2/2]

Geometry & lf::geometry::Geometry::operator= ( Geometry && )
protecteddefault

◆ RefEl()

virtual base::RefEl lf::geometry::Geometry::RefEl ( ) const
pure virtual

◆ SubGeometry()

virtual std::unique_ptr< Geometry > lf::geometry::Geometry::SubGeometry ( dim_t codim,
dim_t i ) const
pure virtual

Construct a new Geometry() object that describes the geometry of the i-th sub-entity with codimension=codim

Parameters
codimThe codimension of the sub-entity (w.r.t. DimLocal())
iThe zero-based index of the sub-entity.
Returns
A new Geometry object that describes the geometry of the specified sub-entity.

Let \( \mathbf{\Phi} : K \mapsto \mathbb{R}^\text{DimGlobal} \) be the mapping of this Geometry object and let \( \mathbf{\xi} : \mathbb{R}^{\text{DimLocal}-codim} \mapsto K\) be the first-order mapping that maps the reference element RefEl().SubType(codim,i) to the i-th sub-entity of RefEl(). I.e. for every node \( \mathbf{n_j} \) of RefEl().SubType(codim,i) it holds that \( \mathbf{\xi}(\mathbf{n_j}) = \) RefEl().NodeCoords(RefEl().SubSubEntity2SubEntity(codim, i, DimLocal()-codim, j)).

Then the geometry element returned by this method describes exactly the mapping \( \mathbf{\Phi} \circ \mathbf{\xi} \)

Implemented in lf::geometry::Point, lf::geometry::QuadO1, lf::geometry::Parallelogram, lf::geometry::QuadO2, lf::geometry::SegmentO1, lf::geometry::SegmentO2, lf::geometry::TriaO1, and lf::geometry::TriaO2.

Referenced by lf::geometry::test_utils::checkSubGeometry().

Friends And Related Symbol Documentation

◆ operator<<()

std::ostream & operator<< ( std::ostream & stream,
const Geometry & geom )
related

Operator overload to print a Geometry to a stream, such as std::cout

Parameters
streamThe stream to which this function should output
geomThe geometry to write to stream.
Returns
The stream itself.
Note
At the moment this will just write the type of the reference element to the stream.

Definition at line 26 of file print_info.cc.

References RefEl().

◆ PrintInfo()

void PrintInfo ( std::ostream & o,
const Geometry & geom,
int output_ctrl = 0 )
related

Diagnostic output operator for Geometry.

Parameters
geomThe geometry to print info about
oThe stream to which this function should output
output_ctrlControls the level of detail of the written output (see below)

Output levels

  • output_ctrl == 0: Reference element type of geometry is printed
  • output_ctrl > 10: The above and global and local dimension and derived type of reference element is printed.
  • output_ctrl > 90: The above and coordinates og the geometry is printed.

Definition at line 5 of file print_info.cc.

References DimGlobal(), DimLocal(), Global(), lf::base::RefEl::NodeCoords(), and RefEl().


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