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

Affine quadrilateral = parallelogram. More...

#include <lf/geometry/quad_o1.h>

Inheritance diagram for lf::geometry::Parallelogram:
lf::geometry::Geometry

Public Member Functions

 Parallelogram (Eigen::Matrix< double, Eigen::Dynamic, 4 > coords)
 Constructor building a parallelogram from four vertex coordinates.
 
 Parallelogram (const Eigen::VectorXd &p0, const Eigen::VectorXd &p1, const Eigen::VectorXd &p2)
 Constructor building a parallelogram from four vertex 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
 
bool isAffine () const override
 element shape by affine mapping from reference element
 
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 ~Geometry ()=default
 Virtual destructor.
 

Private Member Functions

void init ()
 performs initialization of data members
 

Private Attributes

Eigen::Matrix< double, Eigen::Dynamic, 4 > coords_
 Coordinates of the a four vertices, stored in matrix columns.
 
Eigen::Matrix< double, Eigen::Dynamic, 2 > jacobian_
 Matrix of affine mapping generating the parallelogram, agrees with constant Jacobian.
 
Eigen::Matrix< double, Eigen::Dynamic, 2 > jacobian_inverse_gramian_
 Constant matrix ( \( J (J^T J)^{-1}\))
 
double integrationElement_
 constant Gramian determinant
 

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

Affine quadrilateral = parallelogram.

Obtained as affine image of the unit square.

Definition at line 97 of file quad_o1.h.

Constructor & Destructor Documentation

◆ Parallelogram() [1/2]

lf::geometry::Parallelogram::Parallelogram ( Eigen::Matrix< double, Eigen::Dynamic, 4 > coords)
explicit

Constructor building a parallelogram from four vertex coordinates.

Parameters
coordsw x 4 matrix, w = world dimension, whose columns contain the world coordinates of the vertices.
Note
third vertex position is ignored, because it is already fixed by the three others in the case of a parallelogram

Definition at line 258 of file quad_o1.cc.

References lf::geometry::assertNonDegenerateQuad(), coords_, and init().

◆ Parallelogram() [2/2]

lf::geometry::Parallelogram::Parallelogram ( const Eigen::VectorXd & p0,
const Eigen::VectorXd & p1,
const Eigen::VectorXd & p2 )
explicit

Constructor building a parallelogram from four vertex coordinates.

Parameters
p0coordinate vector for vertex 0
p1coordinate vector for vertex 1
p2coordinate vector for vertex 2

Constructs a parallelogram from the coordinate vectors of three spanning vertices.

Definition at line 273 of file quad_o1.cc.

References lf::geometry::assertNonDegenerateQuad(), coords_, and init().

Member Function Documentation

◆ ChildGeometry()

std::vector< std::unique_ptr< Geometry > > lf::geometry::Parallelogram::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

For a detailed description of the indexing of the vertices of child entities see Refinement.xoj.

Implements lf::geometry::Geometry.

Definition at line 351 of file quad_o1.cc.

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

◆ DimGlobal()

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

Dimension of the image of this mapping.

Implements lf::geometry::Geometry.

Definition at line 123 of file quad_o1.h.

References coords_.

Referenced by SubGeometry().

◆ DimLocal()

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

Dimension of the domain of this mapping.

Implements lf::geometry::Geometry.

Definition at line 122 of file quad_o1.h.

◆ Global()

Eigen::MatrixXd lf::geometry::Parallelogram::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 306 of file quad_o1.cc.

References coords_.

Referenced by ChildGeometry().

◆ init()

void lf::geometry::Parallelogram::init ( )
private

performs initialization of data members

Definition at line 290 of file quad_o1.cc.

References coords_, integrationElement_, jacobian_, and jacobian_inverse_gramian_.

Referenced by Parallelogram(), and Parallelogram().

◆ IntegrationElement()

Eigen::VectorXd lf::geometry::Parallelogram::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 321 of file quad_o1.cc.

References integrationElement_.

◆ isAffine()

bool lf::geometry::Parallelogram::isAffine ( ) const
inlineoverridevirtual

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.
In contrast to a general quadrilateral a parallelogram is is the affine image of a square.

Reimplemented from lf::geometry::Geometry.

Definition at line 146 of file quad_o1.h.

◆ Jacobian()

Eigen::MatrixXd lf::geometry::Parallelogram::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 312 of file quad_o1.cc.

References jacobian_.

◆ JacobianInverseGramian()

Eigen::MatrixXd lf::geometry::Parallelogram::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 quad_o1.cc:316

More explanations in Paragraph 2.8.3.14.

Implements lf::geometry::Geometry.

Definition at line 316 of file quad_o1.cc.

References jacobian_inverse_gramian_.

◆ RefEl()

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

The Reference element that defines the domain of this mapping.

Implements lf::geometry::Geometry.

Definition at line 124 of file quad_o1.h.

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

Referenced by SubGeometry().

◆ SubGeometry()

std::unique_ptr< Geometry > lf::geometry::Parallelogram::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 327 of file quad_o1.cc.

References coords_, DimGlobal(), and RefEl().

Member Data Documentation

◆ coords_

Eigen::Matrix<double, Eigen::Dynamic, 4> lf::geometry::Parallelogram::coords_
private

Coordinates of the a four vertices, stored in matrix columns.

In fact, three vertex coordinates of the spanning vertices would be sufficient to define the shape of a parallelogram. For convenience all vertices are stored as in QuadO1.

Definition at line 167 of file quad_o1.h.

Referenced by DimGlobal(), Global(), init(), Parallelogram(), Parallelogram(), and SubGeometry().

◆ integrationElement_

double lf::geometry::Parallelogram::integrationElement_
private

constant Gramian determinant

Definition at line 174 of file quad_o1.h.

Referenced by init(), and IntegrationElement().

◆ jacobian_

Eigen::Matrix<double, Eigen::Dynamic, 2> lf::geometry::Parallelogram::jacobian_
private

Matrix of affine mapping generating the parallelogram, agrees with constant Jacobian.

Definition at line 170 of file quad_o1.h.

Referenced by init(), and Jacobian().

◆ jacobian_inverse_gramian_

Eigen::Matrix<double, Eigen::Dynamic, 2> lf::geometry::Parallelogram::jacobian_inverse_gramian_
private

Constant matrix ( \( J (J^T J)^{-1}\))

Definition at line 172 of file quad_o1.h.

Referenced by init(), and JacobianInverseGramian().


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