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::ScalarLoadElementVectorProvider< SCALAR, MESH_FUNCTION > Class Template Referencefinal

Local computation of general element (load) vector for scalar finite elements; volume contributions only. More...

#include <lf/fe/fe.h>

Public Types

using scalar_t
 
using ElemVec = Eigen::Matrix<scalar_t, Eigen::Dynamic, 1>
 

Public Member Functions

 ScalarLoadElementVectorProvider (std::shared_ptr< const ScalarFESpace< SCALAR > > fe_space, MESH_FUNCTION f)
 Constructor, performs precomputations.
 
bool isActive (const lf::mesh::Entity &) const
 all cells are active
 
ElemVec Eval (const lf::mesh::Entity &cell) const
 
 ~ScalarLoadElementVectorProvider ()=default
 
standard constructors
 ScalarLoadElementVectorProvider (const ScalarLoadElementVectorProvider &)=delete
 
 ScalarLoadElementVectorProvider (ScalarLoadElementVectorProvider &&) noexcept=default
 
ScalarLoadElementVectorProvideroperator= (const ScalarLoadElementVectorProvider &)=delete
 
ScalarLoadElementVectorProvideroperator= (ScalarLoadElementVectorProvider &&)=delete
 

Private Attributes

MESH_FUNCTION f_
 An object providing the source function.
 
std::shared_ptr< const ScalarFESpace< SCALAR > > fe_space_
 
quad::QuadRuleCache qr_cache_
 

Detailed Description

template<base::Scalar SCALAR, mesh::utils::MeshFunction MESH_FUNCTION>
class lf::fe::ScalarLoadElementVectorProvider< SCALAR, MESH_FUNCTION >

Local computation of general element (load) vector for scalar finite elements; volume contributions only.

Template Parameters
SCALARunderlying scalar type of the ScalarFESpace, usually double or complex<double>
MESH_FUNCTIONMeshFunction which defines the source function \( f \)

The underlying local linear form is

\[ v \mapsto \int_K f(\mathbf{x})\,\overline{v(\mathbf{x})}\,\mathrm{d}\mathbf{x}\;, \]

where \(f\) is supposed to be a locally continuous source function.

Computation is based on a quadrature rules supplied by the LehrFEM++ lf::quad::QuadRule module.

This class complies with the requirements for the template parameter ELEM_VEC_COMP of the function assemble::AssembleVectorLocally().

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 569 of file loc_comp_ellbvp.h.

Member Typedef Documentation

◆ ElemVec

template<base::Scalar SCALAR, mesh::utils::MeshFunction MESH_FUNCTION>
using lf::fe::ScalarLoadElementVectorProvider< SCALAR, MESH_FUNCTION >::ElemVec = Eigen::Matrix<scalar_t, Eigen::Dynamic, 1>

Definition at line 575 of file loc_comp_ellbvp.h.

◆ scalar_t

template<base::Scalar SCALAR, mesh::utils::MeshFunction MESH_FUNCTION>
using lf::fe::ScalarLoadElementVectorProvider< SCALAR, MESH_FUNCTION >::scalar_t
Initial value:
decltype(static_cast<SCALAR>(0) *
static_cast<mesh::utils::MeshFunctionReturnType<MESH_FUNCTION>>(
0))

Definition at line 571 of file loc_comp_ellbvp.h.

Constructor & Destructor Documentation

◆ ScalarLoadElementVectorProvider() [1/3]

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

◆ ScalarLoadElementVectorProvider() [2/3]

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

◆ ScalarLoadElementVectorProvider() [3/3]

template<base::Scalar SCALAR, mesh::utils::MeshFunction MESH_FUNCTION>
lf::fe::ScalarLoadElementVectorProvider< SCALAR, MESH_FUNCTION >::ScalarLoadElementVectorProvider ( std::shared_ptr< const ScalarFESpace< SCALAR > > fe_space,
MESH_FUNCTION f )

Constructor, performs precomputations.

Parameters
fe_spacespecification of local shape functions
ffunctor object for source function

Uses quadrature rule of double the degree of exactness compared to the degree of the finite element space.

Definition at line 634 of file loc_comp_ellbvp.h.

◆ ~ScalarLoadElementVectorProvider()

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

Member Function Documentation

◆ Eval()

template<base::Scalar SCALAR, mesh::utils::MeshFunction MESH_FUNCTION>
ScalarLoadElementVectorProvider< SCALAR, MESH_FUNCTION >::ElemVec lf::fe::ScalarLoadElementVectorProvider< SCALAR, MESH_FUNCTION >::Eval ( const lf::mesh::Entity & cell) const

◆ isActive()

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

all cells are active

Definition at line 600 of file loc_comp_ellbvp.h.

◆ operator=() [1/2]

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

◆ operator=() [2/2]

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

Member Data Documentation

◆ f_

template<base::Scalar SCALAR, mesh::utils::MeshFunction MESH_FUNCTION>
MESH_FUNCTION lf::fe::ScalarLoadElementVectorProvider< SCALAR, MESH_FUNCTION >::f_
private

An object providing the source function.

Definition at line 614 of file loc_comp_ellbvp.h.

◆ fe_space_

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

Definition at line 616 of file loc_comp_ellbvp.h.

◆ qr_cache_

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

Definition at line 618 of file loc_comp_ellbvp.h.


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