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

Computing the element matrix for the (negative) Laplacian and linear finite elements. More...

#include <lf/uscalfe/uscalfe.h>

Public Types

using ElemMat = Eigen::Matrix<double, 4, 4>
 

Public Member Functions

 LinearFELaplaceElementMatrix ()=default
 Idle constructor.
 
virtual bool isActive (const lf::mesh::Entity &) const
 All cells are considered active in the default implementation.
 
ElemMat Eval (const lf::mesh::Entity &cell) const
 main routine for the computation of element matrices
 

Static Public Member Functions

static std::shared_ptr< spdlog::logger > & Logger ()
 Used by LinearFELaplaceElementMatrix to log additional (debug) info.
 

Static Private Member Functions

static Eigen::Matrix< double, 4, 2 > DervRefShapFncts (const Eigen::Vector2d &xh)
 

Private Attributes

const double kSqrt3 = 1.0 / std::sqrt(3.0)
 quadrature points on reference quadrilateral
 
const double kZeta_0 = 0.5 - 0.5 * kSqrt3
 
const double kZeta_1 = 0.5 + 0.5 * kSqrt3
 
const std::array< Eigen::Vector2d, 4 > kQuadPoints
 
const Eigen::Matrix< double, 2, 1 > kTriabc
 
const Eigen::Matrix< double, 2, 1 > kQuadpt
 

Detailed Description

Computing the element matrix for the (negative) Laplacian and linear finite elements.

The main purpose of this class is to compute the element matrix for the Laplacian on affine triangles or bilinearly mapped quadrilaterals. These element matrices are provided by the Eval() method.

Note
the Eval() method will always return a reference to a 4x4 matrix also for triangles. In this case the last row and column must be ignored.

This class complies with the requirements for the type ENTITY_MATRIX_PROVIDER given as a template parameter to define an incarnation of the function Cell-Oriented Assembly of Galerkin Matrices .

Definition at line 47 of file lin_fe.h.

Member Typedef Documentation

◆ ElemMat

using lf::uscalfe::LinearFELaplaceElementMatrix::ElemMat = Eigen::Matrix<double, 4, 4>

Definition at line 49 of file lin_fe.h.

Constructor & Destructor Documentation

◆ LinearFELaplaceElementMatrix()

lf::uscalfe::LinearFELaplaceElementMatrix::LinearFELaplaceElementMatrix ( )
default

Idle constructor.

Member Function Documentation

◆ DervRefShapFncts()

Eigen::Matrix< double, 4, 2 > lf::uscalfe::LinearFELaplaceElementMatrix::DervRefShapFncts ( const Eigen::Vector2d & xh)
staticprivate

Definition at line 27 of file lin_fe.cc.

Referenced by Eval().

◆ Eval()

LinearFELaplaceElementMatrix::ElemMat lf::uscalfe::LinearFELaplaceElementMatrix::Eval ( const lf::mesh::Entity & cell) const

main routine for the computation of element matrices

Parameters
cellreference to the (triangular or quadrilateral) cell for which the element matrix should be computed.
Returns
a 4x4 matrix, containing the element matrix. The bottom row/column is not used in the case of a triangle.

Definition at line 39 of file lin_fe.cc.

References DervRefShapFncts(), lf::geometry::Geometry::DimGlobal(), lf::geometry::Geometry::DimLocal(), lf::mesh::Entity::Geometry(), lf::geometry::Geometry::Global(), lf::geometry::Geometry::IntegrationElement(), lf::geometry::Geometry::JacobianInverseGramian(), lf::base::RefEl::kQuad(), kQuadPoints, kQuadpt, lf::base::RefEl::kTria(), kTriabc, Logger(), and lf::mesh::Entity::RefEl().

◆ isActive()

virtual bool lf::uscalfe::LinearFELaplaceElementMatrix::isActive ( const lf::mesh::Entity & ) const
inlinevirtual

All cells are considered active in the default implementation.

Definition at line 58 of file lin_fe.h.

◆ Logger()

std::shared_ptr< spdlog::logger > & lf::uscalfe::LinearFELaplaceElementMatrix::Logger ( )
static

Used by LinearFELaplaceElementMatrix to log additional (debug) info.

Definition at line 21 of file lin_fe.cc.

References lf::base::InitLogger().

Referenced by Eval().

Member Data Documentation

◆ kQuadPoints

const std::array<Eigen::Vector2d, 4> lf::uscalfe::LinearFELaplaceElementMatrix::kQuadPoints
private
Initial value:
{
Eigen::Vector2d{kZeta_0, kZeta_0}, Eigen::Vector2d{kZeta_0, kZeta_1},
Eigen::Vector2d{kZeta_1, kZeta_0}, Eigen::Vector2d{kZeta_1, kZeta_1}}

Definition at line 76 of file lin_fe.h.

Referenced by Eval().

◆ kQuadpt

const Eigen::Matrix<double, 2, 1> lf::uscalfe::LinearFELaplaceElementMatrix::kQuadpt
private
Initial value:
{
(Eigen::Matrix<double, 2, 1>() << 0.7, 0.3).finished()}

Definition at line 86 of file lin_fe.h.

Referenced by Eval().

◆ kSqrt3

const double lf::uscalfe::LinearFELaplaceElementMatrix::kSqrt3 = 1.0 / std::sqrt(3.0)
private

quadrature points on reference quadrilateral

Definition at line 73 of file lin_fe.h.

◆ kTriabc

const Eigen::Matrix<double, 2, 1> lf::uscalfe::LinearFELaplaceElementMatrix::kTriabc
private
Initial value:
{
(Eigen::Matrix<double, 2, 1>() << (1.0 / 3), (1.0 / 3)).finished()}

Definition at line 83 of file lin_fe.h.

Referenced by Eval().

◆ kZeta_0

const double lf::uscalfe::LinearFELaplaceElementMatrix::kZeta_0 = 0.5 - 0.5 * kSqrt3
private

Definition at line 74 of file lin_fe.h.

◆ kZeta_1

const double lf::uscalfe::LinearFELaplaceElementMatrix::kZeta_1 = 0.5 + 0.5 * kSqrt3
private

Definition at line 75 of file lin_fe.h.


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