LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
Loading...
Searching...
No Matches
Public Member Functions | Private Attributes | List of all members
lf::geometry::TriaO2 Class Reference

A second-order triangle in the plane or in 3D space. More...

#include <lf/geometry/tria_o2.h>

Inheritance diagram for lf::geometry::TriaO2:
lf::geometry::Geometry

Public Member Functions

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

Private Attributes

Eigen::Matrix< double, Eigen::Dynamic, 6 > coords_
 Coordinates of the 6 vertices/midpoints, stored in matrix columns.
 
Eigen::Matrix< double, Eigen::Dynamic, 1 > alpha_
 
Eigen::Matrix< double, Eigen::Dynamic, 2 > beta_
 
Eigen::Matrix< double, Eigen::Dynamic, 2 > gamma_
 
Eigen::Matrix< double, Eigen::Dynamic, 1 > delta_
 
Eigen::Matrix< double, Eigen::Dynamic, 2 > gamma_x_2_
 

Additional Inherited Members

- Public Types inherited from lf::geometry::Geometry
using dim_t = base::RefEl::dim_t
 
- Protected Member Functions inherited from lf::geometry::Geometry
 Geometry ()=default
 
 Geometry (const Geometry &)=default
 
 Geometry (Geometry &&)=default
 
Geometryoperator= (const Geometry &)=default
 
Geometryoperator= (Geometry &&)=default
 

Detailed Description

A second-order triangle in the plane or in 3D space.

Coordinates \( coords = [A, B, C, D, E, F] \) are mapped to the reference element as follows:

(0.0, 1.0)                                            C
     |     \                                          | \
(0.0, 0.5)   (0.5, 0.5)                     ->        F   E
     |                  \                             |     \
(0.0, 0.0) - (0.5, 0.0) - (1.0, 0.0)                  A - D - B

Definition at line 29 of file tria_o2.h.

Constructor & Destructor Documentation

◆ TriaO2()

lf::geometry::TriaO2::TriaO2 ( Eigen::Matrix< double, Eigen::Dynamic, 6 > coords)
explicit

Constructor building triangle from vertex/midpoint coordinates.

Parameters
coordsw x 6 matrix, w = world dimension, whose columns contain the world coordinates of the vertices/midpoints.

The numbering convention for vertices and midpoints follows the usual rules for numbering local interpolation nodes in LehrFEM++.

Definition at line 18 of file tria_o2.cc.

References alpha_, beta_, coords_, delta_, gamma_, and gamma_x_2_.

Member Function Documentation

◆ ChildGeometry()

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

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

Implements lf::geometry::Geometry.

Definition at line 173 of file tria_o2.cc.

References lf::geometry::RefinementPattern::ChildPolygons(), Global(), lf::base::RefEl::kTria(), lf::geometry::RefinementPattern::LatticeConst(), lf::geometry::RefinementPattern::NumChildren(), lf::geometry::RefinementPattern::RefEl(), and lf::base::RefEl::ToString().

◆ DimGlobal()

dim_t lf::geometry::TriaO2::DimGlobal ( ) const
inlineoverridevirtual

Dimension of the image of this mapping.

Implements lf::geometry::Geometry.

Definition at line 42 of file tria_o2.h.

References coords_.

Referenced by IntegrationElement(), JacobianInverseGramian(), and SubGeometry().

◆ DimLocal()

dim_t lf::geometry::TriaO2::DimLocal ( ) const
inlineoverridevirtual

Dimension of the domain of this mapping.

Implements lf::geometry::Geometry.

Definition at line 41 of file tria_o2.h.

◆ Global()

Eigen::MatrixXd lf::geometry::TriaO2::Global ( const Eigen::MatrixXd & local) const
overridevirtual

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.

Implements lf::geometry::Geometry.

Definition at line 50 of file tria_o2.cc.

References alpha_, beta_, delta_, and gamma_.

Referenced by ChildGeometry().

◆ IntegrationElement()

Eigen::VectorXd lf::geometry::TriaO2::IntegrationElement ( const Eigen::MatrixXd & local) const
overridevirtual

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.

Implements lf::geometry::Geometry.

Definition at line 123 of file tria_o2.cc.

References DimGlobal(), and Jacobian().

◆ Jacobian()

Eigen::MatrixXd lf::geometry::TriaO2::Jacobian ( const Eigen::MatrixXd & local) const
overridevirtual

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.

Implements lf::geometry::Geometry.

Definition at line 64 of file tria_o2.cc.

References beta_, delta_, gamma_, and gamma_x_2_.

Referenced by IntegrationElement(), and JacobianInverseGramian().

◆ JacobianInverseGramian()

Eigen::MatrixXd lf::geometry::TriaO2::JacobianInverseGramian ( const Eigen::MatrixXd & local) const
overridevirtual

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)
Eigen::MatrixXd JacobianInverseGramian(const Eigen::MatrixXd &local) const override
Evaluate the Jacobian * Inverse Gramian ( ) simultaneously at numPoints.
Definition tria_o2.cc:87

More explanations in Paragraph 2.8.3.14.

Implements lf::geometry::Geometry.

Definition at line 87 of file tria_o2.cc.

References DimGlobal(), and Jacobian().

◆ RefEl()

base::RefEl lf::geometry::TriaO2::RefEl ( ) const
inlineoverridevirtual

The Reference element that defines the domain of this mapping.

Implements lf::geometry::Geometry.

Definition at line 43 of file tria_o2.h.

References lf::base::RefEl::kTria().

◆ SubGeometry()

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

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} \)

Implements lf::geometry::Geometry.

Definition at line 149 of file tria_o2.cc.

References coords_, and DimGlobal().

Member Data Documentation

◆ alpha_

Eigen::Matrix<double, Eigen::Dynamic, 1> lf::geometry::TriaO2::alpha_
private

Definition at line 72 of file tria_o2.h.

Referenced by Global(), and TriaO2().

◆ beta_

Eigen::Matrix<double, Eigen::Dynamic, 2> lf::geometry::TriaO2::beta_
private

Definition at line 73 of file tria_o2.h.

Referenced by Global(), Jacobian(), and TriaO2().

◆ coords_

Eigen::Matrix<double, Eigen::Dynamic, 6> lf::geometry::TriaO2::coords_
private

Coordinates of the 6 vertices/midpoints, stored in matrix columns.

Definition at line 66 of file tria_o2.h.

Referenced by DimGlobal(), SubGeometry(), and TriaO2().

◆ delta_

Eigen::Matrix<double, Eigen::Dynamic, 1> lf::geometry::TriaO2::delta_
private

Definition at line 75 of file tria_o2.h.

Referenced by Global(), Jacobian(), and TriaO2().

◆ gamma_

Eigen::Matrix<double, Eigen::Dynamic, 2> lf::geometry::TriaO2::gamma_
private

Definition at line 74 of file tria_o2.h.

Referenced by Global(), Jacobian(), and TriaO2().

◆ gamma_x_2_

Eigen::Matrix<double, Eigen::Dynamic, 2> lf::geometry::TriaO2::gamma_x_2_
private

Definition at line 78 of file tria_o2.h.

Referenced by Jacobian(), and TriaO2().


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