LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
|
Interface class for shape information on a mesh cell in the spirit of parametric finite element methods More...
#include <lf/geometry/geometry_interface.h>
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< Geometry > | SubGeometry (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 | |
Geometry & | operator= (const Geometry &)=default |
Geometry & | operator= (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 | |
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.
Definition at line 30 of file geometry_interface.h.
|
protecteddefault |
Referenced by lf::mesh::Mesh::PrintInfo(), lf::io::writeMatlab(), lf::io::writeMatplotlib(), and lf::io::writeTikZ().
|
protecteddefault |
|
protecteddefault |
|
virtualdefault |
Virtual destructor.
|
pure virtual |
Generate geometry objects for child entities created in the course of refinement.
ref_pat | An object encoding the details of refinement |
codim | relative codimension of child entities whose shape is requested |
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.
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().
|
pure virtual |
Dimension of the image of this mapping.
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::geometry::test_utils::checkIntegrationElement(), lf::geometry::test_utils::checkJacobian(), lf::geometry::test_utils::checkJacobianInverseGramian(), lf::uscalfe::ReactionDiffusionElementMatrixProvider< SCALAR, DIFF_COEFF, REACTION_COEFF >::Eval(), lf::fe::DiffusionElementMatrixProvider< SCALAR, DIFF_COEFF >::Eval(), lf::fe::MassElementMatrixProvider< SCALAR, REACTION_COEFF >::Eval(), lf::uscalfe::LinearFELaplaceElementMatrix::Eval(), lf::uscalfe::LinearFELocalLoadVector< SCALAR, FUNCTOR >::Eval(), and PrintInfo().
|
pure virtual |
Dimension of the domain of this mapping.
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::geometry::test_utils::checkIntegrationElement(), lf::geometry::test_utils::checkJacobian(), lf::geometry::test_utils::checkJacobianInverseGramian(), lf::geometry::test_utils::checkSubGeometry(), 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::fe::ScalarLoadEdgeVectorProvider< SCALAR, FUNCTOR, EDGESELECTOR >::Eval(), PrintInfo(), lf::mesh::Mesh::PrintInfo(), and lf::geometry::Volume().
|
pure virtual |
Map a number of points in local coordinates into the global coordinate system.
local | A Matrix of size DimLocal() x numPoints that contains in its columns the coordinates of the points at which the mapping function should be evaluated. |
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()
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().
|
pure virtual |
The integration element (factor appearing in integral transformation formula, see below) at number of evaluation points (specified in local coordinates).
local | A Matrix of size DimLocal() x numPoints that contains the evaluation points (in local = reference coordinates) as column vectors. |
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().
|
inlinevirtual |
element shape by affine mapping from reference element
An affine map is a linear mapping plus a translation. For a 2D mesh an entity has an affine geometry, if
Reimplemented in lf::geometry::QuadO1, and lf::geometry::Parallelogram.
Definition at line 220 of file geometry_interface.h.
|
pure virtual |
Evaluate the jacobian of the mapping simultaneously at numPoints
points.
local | A Matrix of size DimLocal x numPoints that contains the evaluation points as column vectors |
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().
|
pure virtual |
Evaluate the Jacobian * Inverse Gramian ( \( J (J^T J)^{-1}\)) simultaneously at numPoints
.
local | A Matrix of size DimLocal() x numPoints that contains the evaluation points as column vectors. |
DimGlobal() x (DimLocal() * numPoints)
that contains the Jacobian multiplied with the Inverse Gramian ( \( J (J^T J)^{-1}\)) at every evaluation point.DimLocal() == DimGlobal()
then \( J (J^T J)^{-1} = J^{-T} \), i.e. this method returns the inverse of the transposed 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()
)
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()().
|
pure virtual |
The Reference element that defines the domain of this mapping.
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::geometry::test_utils::checkChildGeometryVolume(), lf::geometry::test_utils::checkSubGeometry(), lf::geometry::Corners(), operator<<(), PrintInfo(), lf::mesh::Mesh::PrintInfo(), and lf::geometry::Volume().
|
pure virtual |
Construct a new Geometry() object that describes the geometry of the i
-th sub-entity with codimension=codim
codim | The codimension of the sub-entity (w.r.t. DimLocal() ) |
i | The zero-based index of the 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().
|
related |
Operator overload to print a Geometry
to a stream, such as std::cout
stream | The stream to which this function should output |
geom | The geometry to write to stream . |
Definition at line 26 of file print_info.cc.
References RefEl().
|
related |
Diagnostic output operator for Geometry.
geom | The geometry to print info about |
o | The stream to which this function should output |
output_ctrl | Controls the level of detail of the written output (see below) |
Definition at line 5 of file print_info.cc.
References DimGlobal(), DimLocal(), Global(), lf::base::RefEl::NodeCoords(), and RefEl().