LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
Loading...
Searching...
No Matches
Public Types | Public Member Functions | Private Attributes | List of all members
lf::fe::DiffusionElementMatrixProvider< SCALAR, DIFF_COEFF > Class Template Referencefinal

Class for computing element matrices for general scalar-valued finite elements and homogeneous 2nd-order elliptic bilinear forms. More...

#include <lf/fe/fe.h>

Public Types

using Scalar
 type of returned element matrix
 
using ElemMat = Eigen::Matrix<Scalar, Eigen::Dynamic, Eigen::Dynamic>
 

Public Member Functions

 DiffusionElementMatrixProvider (std::shared_ptr< const ScalarFESpace< SCALAR > > fe_space, DIFF_COEFF alpha)
 Constructor.
 
bool isActive (const lf::mesh::Entity &) const
 All cells are considered active.
 
ElemMat Eval (const lf::mesh::Entity &cell) const
 Routine for the computation of element matrices.
 
 ~DiffusionElementMatrixProvider ()=default
 
 DiffusionElementMatrixProvider (const DiffusionElementMatrixProvider &)=delete
 standard constructors
 
 DiffusionElementMatrixProvider (DiffusionElementMatrixProvider &&) noexcept=default
 
DiffusionElementMatrixProvideroperator= (const DiffusionElementMatrixProvider &)=delete
 
DiffusionElementMatrixProvideroperator= (DiffusionElementMatrixProvider &&)=delete
 

Private Attributes

std::shared_ptr< const ScalarFESpace< SCALAR > > fe_space_
 
quad::QuadRuleCache qr_cache_
 
functors providing coefficient functions
DIFF_COEFF alpha_
 

Detailed Description

template<base::Scalar SCALAR, mesh::utils::MeshFunction DIFF_COEFF>
class lf::fe::DiffusionElementMatrixProvider< SCALAR, DIFF_COEFF >

Class for computing element matrices for general scalar-valued finite elements and homogeneous 2nd-order elliptic bilinear forms.

Template Parameters
SCALARscalar type of the ScalarFESpace Must be a field type such as double or std::complex<double>
DIFF_COEFFa MeshFunction that defines the diffusion coefficient \( \mathbf{\alpha} \). It should be either scalar- or matrix-valued.
Note
This class complies with the type requirements for the template argument ENTITY_MATRIX_PROVIDER of the function lf::assemble::AssembleMatrixLocally().

The element matrix corresponds to the (local) bilinear form

\[ (u,v) \mapsto\int\limits_{K}\boldsymbol{\alpha}(\mathbf{x})\mathbf{grad}\,u \cdot\mathbf{grad}\,v\,\mathrm{d}\mathbf{x} \;, \]

with diffusion coefficient \(\mathbf{\alpha}\), see also Example 2.8.3.29

Template parameter requirement

Note
If you intend to use this matrix provider with FE Spaces from lf::uscalfe please consider switching to lf::uscalfe::ReactionDiffusionElementMatrixProvider provided there, as it is specifically optimized for uniform FE Spaces.

Logger

This class logs additional information to DiffusionElementMatrixProviderLogger(). See Loggers and Debug output for more information.

Example usage

The following code snippet computes the solution of the BVP

\begin{align} - \Delta u &= 1 && \text{on }\Omega := [0,1]^2 \\ u &= 0 && \text{on }\partial \Omega \end{align}

// Get a mesh
// Create HierarchicalFESpace with uniform degree 5.
const unsigned degree = 5;
const auto fe_space =
std::make_shared<lf::fe::HierarchicScalarFESpace<double>>(mesh_p, degree);
// define diffusion coefficient
// define rhs load
// Assemble the system matrix and right hand side
Eigen::VectorXd rhs = Eigen::VectorXd::Zero(fe_space->LocGlobMap().NumDofs());
const assemble::DofHandler &dofh = fe_space->LocGlobMap();
assemble::COOMatrix<double> A_COO(dofh.NumDofs(), dofh.NumDofs());
DiffusionElementMatrixProvider element_matrix_provider(fe_space, mf_alpha);
AssembleMatrixLocally(0, dofh, dofh, element_matrix_provider, A_COO);
ScalarLoadElementVectorProvider element_vector_provider(fe_space, mf_load);
AssembleVectorLocally(0, dofh, element_vector_provider, rhs);
// Enforce zero dirichlet boundary conditions
const auto boundary = lf::mesh::utils::flagEntitiesOnBoundary(mesh_p);
const auto selector = [&](unsigned int idx) -> std::pair<bool, double> {
const auto &entity = dofh.Entity(idx);
return {entity.Codim() > 0 && boundary(entity), 0};
};
FixFlaggedSolutionComponents(selector, A_COO, rhs);
// Solve the LSE using the cholesky decomposition
Eigen::SparseMatrix<double> A = A_COO.makeSparse();
Eigen::SimplicialLDLT<Eigen::SparseMatrix<double>> solver(A);
const Eigen::VectorXd solution = solver.solve(rhs);
// visualize the solution:
io::VtkWriter vtk(mesh_p, "solution.vtk", 0, 5);
vtk.WritePointData("solution", MeshFunctionFE(fe_space, solution));

Definition at line 76 of file loc_comp_ellbvp.h.

Member Typedef Documentation

◆ ElemMat

template<base::Scalar SCALAR, mesh::utils::MeshFunction DIFF_COEFF>
using lf::fe::DiffusionElementMatrixProvider< SCALAR, DIFF_COEFF >::ElemMat = Eigen::Matrix<Scalar, Eigen::Dynamic, Eigen::Dynamic>

Definition at line 84 of file loc_comp_ellbvp.h.

◆ Scalar

template<base::Scalar SCALAR, mesh::utils::MeshFunction DIFF_COEFF>
using lf::fe::DiffusionElementMatrixProvider< SCALAR, DIFF_COEFF >::Scalar
Initial value:
typename decltype(mesh::utils::MeshFunctionReturnType<DIFF_COEFF>() *
Eigen::Matrix<SCALAR, Eigen::Dynamic, 1>())::Scalar
typename decltype(mesh::utils::MeshFunctionReturnType< DIFF_COEFF >() * Eigen::Matrix< SCALAR, Eigen::Dynamic, 1 >())::Scalar Scalar
type of returned element matrix

type of returned element matrix

Definition at line 81 of file loc_comp_ellbvp.h.

Constructor & Destructor Documentation

◆ DiffusionElementMatrixProvider() [1/3]

template<base::Scalar SCALAR, mesh::utils::MeshFunction DIFF_COEFF>
lf::fe::DiffusionElementMatrixProvider< SCALAR, DIFF_COEFF >::DiffusionElementMatrixProvider ( const DiffusionElementMatrixProvider< SCALAR, DIFF_COEFF > & )
delete

standard constructors

◆ DiffusionElementMatrixProvider() [2/3]

template<base::Scalar SCALAR, mesh::utils::MeshFunction DIFF_COEFF>
lf::fe::DiffusionElementMatrixProvider< SCALAR, DIFF_COEFF >::DiffusionElementMatrixProvider ( DiffusionElementMatrixProvider< SCALAR, DIFF_COEFF > && )
defaultnoexcept

◆ DiffusionElementMatrixProvider() [3/3]

template<base::Scalar SCALAR, mesh::utils::MeshFunction DIFF_COEFF>
lf::fe::DiffusionElementMatrixProvider< SCALAR, DIFF_COEFF >::DiffusionElementMatrixProvider ( std::shared_ptr< const ScalarFESpace< SCALAR > > fe_space,
DIFF_COEFF alpha )

Constructor.

Parameters
fe_spacecollection of specifications for scalar-valued parametric reference elements
alphamesh function for the (possibly matrix-valued) diffusion coefficient
Note
If your fe_space is from lf::uscalfe, please consider using the lf::uscalfe::DiffusionElementMatrixProvider provided there, as that is specifically optimized for uniform FE spaces.

Definition at line 161 of file loc_comp_ellbvp.h.

◆ ~DiffusionElementMatrixProvider()

template<base::Scalar SCALAR, mesh::utils::MeshFunction DIFF_COEFF>
lf::fe::DiffusionElementMatrixProvider< SCALAR, DIFF_COEFF >::~DiffusionElementMatrixProvider ( )
default

destructor

Member Function Documentation

◆ Eval()

template<base::Scalar SCALAR, mesh::utils::MeshFunction DIFF_COEFF>
lf::fe::DiffusionElementMatrixProvider< SCALAR, DIFF_COEFF >::ElemMat lf::fe::DiffusionElementMatrixProvider< SCALAR, DIFF_COEFF >::Eval ( const lf::mesh::Entity & cell) const

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 small dense element matrix.

Evaluation of the integral over the basis functions is done using numerical quadrature with double the degree of exactness as the polynomial degree of the basis functions on the provided cell.

Throws an assertion in case the finite element specification is missing for the type of the cell.

Definition at line 171 of file loc_comp_ellbvp.h.

References lf::fe::DiffusionElementMatrixProviderLogger(), 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::NodeCoords(), lf::quad::QuadRule::NumPoints(), lf::quad::QuadRule::Points(), lf::mesh::Entity::RefEl(), and lf::quad::QuadRule::Weights().

◆ isActive()

template<base::Scalar SCALAR, mesh::utils::MeshFunction DIFF_COEFF>
bool lf::fe::DiffusionElementMatrixProvider< SCALAR, DIFF_COEFF >::isActive ( const lf::mesh::Entity & ) const
inline

All cells are considered active.

Definition at line 116 of file loc_comp_ellbvp.h.

◆ operator=() [1/2]

template<base::Scalar SCALAR, mesh::utils::MeshFunction DIFF_COEFF>
DiffusionElementMatrixProvider & lf::fe::DiffusionElementMatrixProvider< SCALAR, DIFF_COEFF >::operator= ( const DiffusionElementMatrixProvider< SCALAR, DIFF_COEFF > & )
delete

◆ operator=() [2/2]

template<base::Scalar SCALAR, mesh::utils::MeshFunction DIFF_COEFF>
DiffusionElementMatrixProvider & lf::fe::DiffusionElementMatrixProvider< SCALAR, DIFF_COEFF >::operator= ( DiffusionElementMatrixProvider< SCALAR, DIFF_COEFF > && )
delete

Member Data Documentation

◆ alpha_

template<base::Scalar SCALAR, mesh::utils::MeshFunction DIFF_COEFF>
DIFF_COEFF lf::fe::DiffusionElementMatrixProvider< SCALAR, DIFF_COEFF >::alpha_
private

Diffusion coefficient

Definition at line 140 of file loc_comp_ellbvp.h.

◆ fe_space_

template<base::Scalar SCALAR, mesh::utils::MeshFunction DIFF_COEFF>
std::shared_ptr<const ScalarFESpace<SCALAR> > lf::fe::DiffusionElementMatrixProvider< SCALAR, DIFF_COEFF >::fe_space_
private

FE Space

Definition at line 144 of file loc_comp_ellbvp.h.

◆ qr_cache_

template<base::Scalar SCALAR, mesh::utils::MeshFunction DIFF_COEFF>
quad::QuadRuleCache lf::fe::DiffusionElementMatrixProvider< SCALAR, DIFF_COEFF >::qr_cache_
private

Definition at line 146 of file loc_comp_ellbvp.h.


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