LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
Loading...
Searching...
No Matches
Namespaces | Classes | Functions
lf::geometry Namespace Reference

Defines the Geometry interface and provides a number of classes that implement this interface + additional geometry related helper routines. More...

Namespaces

namespace  test_utils
 Defines the Geometry::test_utils module and provides a number of test functions to check geometry objects.
 

Classes

class  Geometry
 Interface class for shape information on a mesh cell in the spirit of parametric finite element methods More...
 
class  Parallelogram
 Affine quadrilateral = parallelogram. More...
 
class  Point
 
class  QuadO1
 Bilinear quadrilateral element shape. More...
 
class  QuadO2
 A second-order quadrilateral in the plane or in 3D space. More...
 
class  RefinementPattern
 Abstract interface class for encoding topological local refinement More...
 
class  SegmentO1
 A straight edge defined by the location of its two endpoints. More...
 
class  SegmentO2
 A second-order segment in the plane or in 3D space. More...
 
class  TriaO1
 An affine triangle in the plane or in 3D space. More...
 
class  TriaO2
 A second-order triangle in the plane or in 3D space. More...
 

Functions

std::unique_ptr< GeometryCompose (std::unique_ptr< Geometry > &&a, std::unique_ptr< Geometry > &&b)
 Create a new geometry object that is the composition of \( \Phi_a \circ \Phi_b \).
 
double Volume (const Geometry &geo)
 Compute the (approximate) volume (area) of a shape.
 
Eigen::MatrixXd Corners (const Geometry &geo)
 The corners of a shape with piecewise smooth boundary.
 
void PrintInfo (std::ostream &o, const Geometry &geom, int output_ctrl)
 
std::ostream & operator<< (std::ostream &stream, const Geometry &geom)
 
bool assertNonDegenerateQuad (const Eigen::Matrix< double, Eigen::Dynamic, 4 > &coords, double tol=1.0E-8)
 Asserting a non-degenerate bilinear quadrilateral.
 
bool isParallelogram (const Eigen::Matrix< int, Eigen::Dynamic, Eigen::Dynamic > &polygon)
 Test whether a lattice polygon describes a logical parallelogram.
 
bool assertNonDegenerateTriangle (const Eigen::Matrix< double, Eigen::Dynamic, 3 > &coords, double tol=1.0E-8)
 Asserting non-degenerate shape of a flat triangle.
 

Detailed Description

Defines the Geometry interface and provides a number of classes that implement this interface + additional geometry related helper routines.

Function Documentation

◆ assertNonDegenerateQuad()

bool lf::geometry::assertNonDegenerateQuad ( const Eigen::Matrix< double, Eigen::Dynamic, 4 > & coords,
double tol = 1.0E-8 )

Asserting a non-degenerate bilinear quadrilateral.

Parameters
coordsw x 4 Eigen matrix whose columns contain the vertex coordinates of the quadrilateral, w = world dimension
tolrelative tolerance for numerical tests of equality with zero

Terminates execution in case degenerate shape is detected.

Definition at line 12 of file quad_o1.cc.

Referenced by lf::geometry::Parallelogram::Parallelogram(), lf::geometry::Parallelogram::Parallelogram(), and lf::geometry::QuadO1::QuadO1().

◆ assertNonDegenerateTriangle()

bool lf::geometry::assertNonDegenerateTriangle ( const Eigen::Matrix< double, Eigen::Dynamic, 3 > & coords,
double tol = 1.0E-8 )

Asserting non-degenerate shape of a flat triangle.

Parameters
coordsw x 3 Eigen matrix whose columns contain the vertex coordinates of the quadrilateral, w = world dimension
tolrelative tolerance for numerical tests of equality with zero

Terminates execution in case degenerate shape is detected.

Definition at line 10 of file tria_o1.cc.

Referenced by lf::geometry::TriaO1::TriaO1().

◆ Compose()

std::unique_ptr< Geometry > lf::geometry::Compose ( std::unique_ptr< Geometry > && a,
std::unique_ptr< Geometry > && b )

Create a new geometry object that is the composition of \( \Phi_a \circ \Phi_b \).

Parameters
aThe mapping \( \Phi_a \)
bThe mapping \( \Phi_b \)
Returns
The composition of \(\Phi_a \circ \Phi_b\)

◆ Corners()

Eigen::MatrixXd lf::geometry::Corners ( const Geometry & geo)
inline

The corners of a shape with piecewise smooth boundary.

Parameters
geoThe geometry object
Returns
the coordinates vectors for the corners of the shape represented by the geometry object packed into the columns of a dxn-matrix, where d is the dimension of physical space, and n stands for the number of corners.
Note
Only for 2D shapes with all straight edges the positions of the corners will completely define the shape. This is the case with shapes represented by the types lf::geometry::TriaO1 and lf::geometry::QuadO1

Demonstration code [Lecture Document]

(https://www.sam.math.ethz.ch/~grsam/NUMPDEFL/NUMPDE.pdf) Code 2.7.2.23

void PrintGeometryInfo(const lf::mesh::Mesh &mesh, lf::base::dim_t codim) {
// loop over all entities of the specified codimension
for (const lf::mesh::Entity *ent : mesh.Entities(codim)) {
// Number of nodes = number of corner points
const lf::base::size_type num_nodes = ent->RefEl().NumNodes();
// Obtain pointer to geometry object associated with entity
const lf::geometry::Geometry *geo_ptr = ent->Geometry();
// Fetch coordinates of corner points in packed format \cref{par:coords}
Eigen::MatrixXd corners = lf::geometry::Corners(*geo_ptr);
std::cout << ent->RefEl() << "(" << mesh.Index(*ent) << ") pts: ";
for (int l = 0; l < num_nodes; ++l) {
std::cout << l << " =[" << corners.col(l).transpose() << "], ";
}
std::cout << std::endl;
}
}

Additional explanations can be found in [Lecture Document] (https://www.sam.math.ethz.ch/~grsam/NUMPDEFL/NUMPDE.pdf) Paragraph 2.7.2.21.

Definition at line 261 of file geometry_interface.h.

References lf::geometry::Geometry::Global(), lf::base::RefEl::NodeCoords(), and lf::geometry::Geometry::RefEl().

◆ isParallelogram()

bool lf::geometry::isParallelogram ( const Eigen::Matrix< int, Eigen::Dynamic, Eigen::Dynamic > & polygon)

Test whether a lattice polygon describes a logical parallelogram.

Parameters
polygonan integer matrix whose column contain the lattice coordinates of the vertices of the polygon.

A polygon passes the parallelogram test, if

  1. it has four vertices
  2. the vectors \(x_1-x_0\) and \(x_2-x_3\) agree.

Definition at line 18 of file refinement_pattern.cc.

Referenced by lf::geometry::Parallelogram::ChildGeometry().

◆ operator<<()

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

Definition at line 26 of file print_info.cc.

◆ PrintInfo()

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

Definition at line 5 of file print_info.cc.

◆ Volume()

double lf::geometry::Volume ( const Geometry & geo)

Compute the (approximate) volume (area) of a shape.

Parameters
geoThe geometry object
Returns
approximate volume
Note
the volume can be computed exactly only for planar affine/bilinear (2D) shapes Otherwise this functions returns a one-point quadrature approximation

Definition at line 11 of file geometry_interface.cc.

References lf::geometry::Geometry::DimLocal(), lf::geometry::Geometry::IntegrationElement(), lf::base::RefEl::kPoint(), lf::base::RefEl::kQuad(), lf::base::RefEl::kSegment(), lf::base::RefEl::kTria(), lf::geometry::Geometry::RefEl(), and lf::base::RefEl::ToString().

Referenced by lf::uscalfe::LinearFELocalLoadVector< SCALAR, FUNCTOR >::Eval().