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

Collects data structures and algorithms designed for scalar finite element methods primarily meant for second-order elliptic boundary value problems. More...

Classes

class  FeLagrangeO1Quad
 Linear Lagrange finite element on the quadrilateral reference element. More...
 
class  FeLagrangeO1Segment
 Linear Lagrange finite element on a line segment. More...
 
class  FeLagrangeO1Tria
 Linear Lagrange finite element on triangular reference element. More...
 
class  FeLagrangeO2Quad
 Quadratic Lagrangian finite element on a quadrilateral reference element. More...
 
class  FeLagrangeO2Segment
 Quadratic Lagrangian finite element on a line segment. More...
 
class  FeLagrangeO2Tria
 Quadratic Lagrangian finite element on a triangular reference element. More...
 
class  FeLagrangeO3Quad
 Cubic Lagrangian finite element on a quadrilateral reference element. More...
 
class  FeLagrangeO3Segment
 Cubic Lagrangian finite element on a line segment. More...
 
class  FeLagrangeO3Tria
 Cubic Lagrangian finite elment on a triangular reference element. More...
 
class  FeSpaceLagrangeO1
 (Bi)Linear Lagrangian Finite Element space More...
 
class  FeSpaceLagrangeO2
 Quadratic Lagrangian Finite Element space. More...
 
class  FeSpaceLagrangeO3
 Cubic Lagrangian Finite Element space. More...
 
class  LinearFELaplaceElementMatrix
 Computing the element matrix for the (negative) Laplacian and linear finite elements. More...
 
class  LinearFELocalLoadVector
 Class for computation of local load vector for linear finite elements. More...
 
class  MassEdgeMatrixProvider
 Quadrature-based computation of local mass matrix for an edge. More...
 
class  PrecomputedScalarReferenceFiniteElement
 Helper class which stores a ScalarReferenceFiniteElement with precomputed values at the nodes of a quadrature rule. More...
 
class  ReactionDiffusionElementMatrixProvider
 Class for local quadrature based computations for Lagrangian finite elements and second-order scalar elliptic BVPs. More...
 
class  ScalarLoadEdgeVectorProvider
 Local edge contributions to element vector. More...
 
class  ScalarLoadElementVectorProvider
 Local computation of general element (load) vector for scalar finite elements; volume contributions only. More...
 
class  UniformScalarFESpace
 Space of scalar valued finite element functions on a hybrid 2D mesh More...
 

Typedefs

using gdof_idx_t = lf::assemble::gdof_idx_t
 
using ldof_idx_t = lf::assemble::ldof_idx_t
 
using size_type = lf::assemble::size_type
 
using dim_t = lf::assemble::dim_t
 
using glb_idx_t = lf::assemble::glb_idx_t
 
using sub_idx_t = lf::base::sub_idx_t
 
using quad_rule_collection_t = std::map<lf::base::RefEl, lf::quad::QuadRule>
 Auxiliary data structure for passing collections of quadrature rules.
 

Functions

std::shared_ptr< spdlog::logger > & LinearFeLocalLoadVectorLogger ()
 Used by LinearFELocalLoadVector to log additional (debug) information during vector assembly.
 
std::shared_ptr< spdlog::logger > & ReactionDiffusionElementMatrixProviderLogger ()
 logger for ReactionDiffusionElementMatrixProvider
 
std::shared_ptr< spdlog::logger > & MassEdgeMatrixProviderLogger ()
 logger for MassEdgeMatrixProvider
 
std::shared_ptr< spdlog::logger > & ScalarLoadElementVectorProviderLogger ()
 logger used by ScalarLoadElementVectorProvider
 
std::shared_ptr< spdlog::logger > & ScalarLoadEdgeVectorProviderLogger ()
 logger for ScalarLoadEdgeVectorProvider class template.
 
template<class PTR , class DIFF_COEFF , class REACTION_COEFF >
 ReactionDiffusionElementMatrixProvider (PTR fe_space, DIFF_COEFF alpha, REACTION_COEFF gamma) -> ReactionDiffusionElementMatrixProvider< typename PTR::element_type::Scalar, DIFF_COEFF, REACTION_COEFF >
 
template<class PTR , class DIFF_COEFF , class REACTION_COEFF >
 ReactionDiffusionElementMatrixProvider (PTR fe_space, DIFF_COEFF alpha, REACTION_COEFF gamma, std::map< lf::base::RefEl, lf::quad::QuadRule >) -> ReactionDiffusionElementMatrixProvider< typename PTR::element_type::Scalar, DIFF_COEFF, REACTION_COEFF >
 
template<class PTR , mesh::utils::MeshFunction COEFF, class EDGESELECTOR = base::PredicateTrue>
 MassEdgeMatrixProvider (PTR, COEFF coeff, EDGESELECTOR edge_predicate=base::PredicateTrue{}) -> MassEdgeMatrixProvider< typename PTR::element_type::Scalar, COEFF, EDGESELECTOR >
 
template<class PTR , class MESH_FUNCTION >
 ScalarLoadElementVectorProvider (PTR fe_space, MESH_FUNCTION mf) -> ScalarLoadElementVectorProvider< typename PTR::element_type::Scalar, MESH_FUNCTION >
 
template<class PTR , class FUNCTOR , class EDGESELECTOR = base::PredicateTrue>
 ScalarLoadEdgeVectorProvider (PTR, FUNCTOR, EDGESELECTOR=base::PredicateTrue{}) -> ScalarLoadEdgeVectorProvider< typename PTR::element_type::Scalar, FUNCTOR, EDGESELECTOR >
 
template<typename SCALAR >
void PrintInfo (std::ostream &o, const UniformScalarFESpace< SCALAR > &fes, unsigned ctrl)
 
template<typename SCALAR >
std::ostream & operator<< (std::ostream &o, const UniformScalarFESpace< SCALAR > &fes)
 

Variables

const unsigned int kUniformScalarFESpaceOutMesh = 1
 mesh information will be printed
 
const unsigned int kUniformScalarFESpaceOutDofh = 2
 information about the dof handler will be printed.
 
const unsigned int kUniformScalarFESpaceOutRsfs = 4
 information about the reference shape functions will be printed
 

Detailed Description

Collects data structures and algorithms designed for scalar finite element methods primarily meant for second-order elliptic boundary value problems.

This namespace contains a number of classes/functions which can be used to solve boundary value problems with uniform, scalar Finite elements:

Examples of approximation spaces that the methods/classes in this namespace can represent/handle are:

Here are examples of use cases not supported by lf::uscalfe:

Typedef Documentation

◆ dim_t

Type for (co-)dimensions

Definition at line 32 of file lagr_fe.h.

◆ gdof_idx_t

Type for indices into global matrices/vectors

Definition at line 26 of file lagr_fe.h.

◆ glb_idx_t

Type for global index of entities

Definition at line 34 of file lagr_fe.h.

◆ ldof_idx_t

Type for indices referring to entity matrices/vectors

Definition at line 28 of file lagr_fe.h.

◆ quad_rule_collection_t

Auxiliary data structure for passing collections of quadrature rules.

This type can be used to pass several quadrature rules to a function, when different quadrature rules for different types of entities are required.

Definition at line 32 of file loc_comp_ellbvp.h.

◆ size_type

Type for vector length/matrix sizes

Definition at line 30 of file lagr_fe.h.

◆ sub_idx_t

Type for indexing sub-entities

Definition at line 36 of file lagr_fe.h.

Function Documentation

◆ LinearFeLocalLoadVectorLogger()

std::shared_ptr< spdlog::logger > & lf::uscalfe::LinearFeLocalLoadVectorLogger ( )

Used by LinearFELocalLoadVector to log additional (debug) information during vector assembly.

Definition at line 167 of file lin_fe.cc.

References lf::base::InitLogger().

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

◆ MassEdgeMatrixProvider()

template<class PTR , mesh::utils::MeshFunction COEFF, class EDGESELECTOR = base::PredicateTrue>
lf::uscalfe::MassEdgeMatrixProvider ( PTR ,
COEFF coeff,
EDGESELECTOR edge_predicate = base::PredicateTrue{} ) -> MassEdgeMatrixProvider< typename PTR::element_type::Scalar, COEFF, EDGESELECTOR >

◆ MassEdgeMatrixProviderLogger()

std::shared_ptr< spdlog::logger > & lf::uscalfe::MassEdgeMatrixProviderLogger ( )

logger for MassEdgeMatrixProvider

Definition at line 27 of file loc_comp_ellbvp.cc.

References lf::base::InitLogger().

◆ operator<<()

template<typename SCALAR >
std::ostream & lf::uscalfe::operator<< ( std::ostream & o,
const UniformScalarFESpace< SCALAR > & fes )
related

output operator for scalar parametric finite element space

Definition at line 445 of file uniform_scalar_fe_space.h.

◆ PrintInfo()

template<typename SCALAR >
void lf::uscalfe::PrintInfo ( std::ostream & o,
const UniformScalarFESpace< SCALAR > & fes,
unsigned ctrl )

◆ ReactionDiffusionElementMatrixProvider() [1/2]

template<class PTR , class DIFF_COEFF , class REACTION_COEFF >
lf::uscalfe::ReactionDiffusionElementMatrixProvider ( PTR fe_space,
DIFF_COEFF alpha,
REACTION_COEFF gamma ) -> ReactionDiffusionElementMatrixProvider< typename PTR::element_type::Scalar, DIFF_COEFF, REACTION_COEFF >

◆ ReactionDiffusionElementMatrixProvider() [2/2]

template<class PTR , class DIFF_COEFF , class REACTION_COEFF >
lf::uscalfe::ReactionDiffusionElementMatrixProvider ( PTR fe_space,
DIFF_COEFF alpha,
REACTION_COEFF gamma,
std::map< lf::base::RefEl, lf::quad::QuadRule >  ) -> ReactionDiffusionElementMatrixProvider< typename PTR::element_type::Scalar, DIFF_COEFF, REACTION_COEFF >

◆ ReactionDiffusionElementMatrixProviderLogger()

std::shared_ptr< spdlog::logger > & lf::uscalfe::ReactionDiffusionElementMatrixProviderLogger ( )

◆ ScalarLoadEdgeVectorProvider()

template<class PTR , class FUNCTOR , class EDGESELECTOR = base::PredicateTrue>
lf::uscalfe::ScalarLoadEdgeVectorProvider ( PTR ,
FUNCTOR ,
EDGESELECTOR = base::PredicateTrue{} ) -> ScalarLoadEdgeVectorProvider< typename PTR::element_type::Scalar, FUNCTOR, EDGESELECTOR >

◆ ScalarLoadEdgeVectorProviderLogger()

std::shared_ptr< spdlog::logger > & lf::uscalfe::ScalarLoadEdgeVectorProviderLogger ( )

logger for ScalarLoadEdgeVectorProvider class template.

Definition at line 39 of file loc_comp_ellbvp.cc.

References lf::base::InitLogger().

◆ ScalarLoadElementVectorProvider()

template<class PTR , class MESH_FUNCTION >
lf::uscalfe::ScalarLoadElementVectorProvider ( PTR fe_space,
MESH_FUNCTION mf ) -> ScalarLoadElementVectorProvider< typename PTR::element_type::Scalar, MESH_FUNCTION >

◆ ScalarLoadElementVectorProviderLogger()

std::shared_ptr< spdlog::logger > & lf::uscalfe::ScalarLoadElementVectorProviderLogger ( )

Variable Documentation

◆ kUniformScalarFESpaceOutDofh

const unsigned int lf::uscalfe::kUniformScalarFESpaceOutDofh = 2

information about the dof handler will be printed.

To be used with PrintInfo(std::ostream&, const UniformScalarFESpace<SCALAR>&, unsigned int)

Definition at line 199 of file uniform_scalar_fe_space.h.

Referenced by PrintInfo().

◆ kUniformScalarFESpaceOutMesh

const unsigned int lf::uscalfe::kUniformScalarFESpaceOutMesh = 1

mesh information will be printed

Output control variables for PrintInfo:

To be used with PrintInfo(std::ostream&, const UniformScalarFESpace<SCALAR>&, unsigned int)

Definition at line 191 of file uniform_scalar_fe_space.h.

Referenced by PrintInfo().

◆ kUniformScalarFESpaceOutRsfs

const unsigned int lf::uscalfe::kUniformScalarFESpaceOutRsfs = 4

information about the reference shape functions will be printed

To be used with PrintInfo(std::ostream&, const UniformScalarFESpace<SCALAR>&, unsigned int)

Definition at line 207 of file uniform_scalar_fe_space.h.

Referenced by PrintInfo().