![]() |
LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
|
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) |
| template<class SCALAR> | |
| void | PrintInfo (std::ostream &o, const UniformScalarFESpace< SCALAR > &fes, unsigned int ctrl=0) |
| Print information about a UniformScalarFESpace to the given stream object. | |
| template<typename SCALAR> | |
| std::ostream & | operator<< (std::ostream &o, const UniformScalarFESpace< SCALAR > &fes) |
| output operator for scalar parametric finite element space | |
| template<MeshFunction A, MeshFunction B> | |
| auto | operator* (const A &a, const B &b) -> MeshFunctionBinary< internal::OperatorMultiplication, A, B > |
| Multiply two mesh functions with each other. | |
| template<MeshFunction A, MeshFunction B> | |
| auto | operator+ (const A &a, const B &b) -> MeshFunctionBinary< internal::OperatorAddition, A, B > |
| Add's two mesh functions. | |
| template<MeshFunction A, MeshFunction B> | |
| auto | operator- (const A &a, const B &b) -> MeshFunctionBinary< internal::OperatorSubtraction, A, B > |
| Subtracts two mesh functions. | |
| template<MeshFunction A> | |
| auto | adjoint (const A &a) -> MeshFunctionUnary< internal::UnaryOpAdjoint, A > |
Pointwise adjoint of an Eigen::Matrix, i.e. the conjugate transposed of the matrix. | |
| template<MeshFunction A> | |
| auto | conjugate (const A &a) -> MeshFunctionUnary< internal::UnaryOpConjugate, A > |
Pointwise conjuagte of an Eigen::Matrix, Eigen::Array or scalar valued mesh function. | |
| template<MeshFunction A> | |
| auto | squaredNorm (const A &a) -> MeshFunctionUnary< internal::UnaryOpSquaredNorm, A > |
| Pointwise squared norm of another mesh function. | |
| template<MeshFunction A> | |
| auto | transpose (const A &a) -> MeshFunctionUnary< internal::UnaryOpTranspose, A > |
Pointwise transpose of an Eigen::Matrix or Eigen::Array | |
Variables | |
| template<typename SCALAR> | |
| const FeLagrangeO2Segment< SCALAR > | FeLagrangeO2Quad< SCALAR >::krsf_segment_ |
| template<typename SCALAR> | |
| const Eigen::Matrix< int, 9, 2 > | FeLagrangeO2Quad< SCALAR >::ksegment_to_quad_mapping_ |
| template<typename SCALAR> | |
| const FeLagrangeO3Segment< SCALAR > | FeLagrangeO3Quad< SCALAR >::krsf_segment_ |
| template<typename SCALAR> | |
| const Eigen::Matrix< int, 16, 2 > | FeLagrangeO3Quad< SCALAR >::ksegment_to_quad_mapping_ |
| 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 | |
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:
| using lf::uscalfe::quad_rule_collection_t = std::map<lf::base::RefEl, lf::quad::QuadRule> |
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.
|
Pointwise adjoint of an Eigen::Matrix, i.e. the conjugate transposed of the matrix.
| A | The type of the original mesh function. |
| a | The original mesh function whose adjoint should be taken. |
a. Definition at line 308 of file mesh_function_unary.h.
|
Pointwise conjuagte of an Eigen::Matrix, Eigen::Array or scalar valued mesh function.
| A | The type of the original mesh function. |
| a | The original mesh function that should be conjugated. |
a. Definition at line 325 of file mesh_function_unary.h.
| 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().
| lf::uscalfe::MassEdgeMatrixProvider | ( | PTR | , |
| COEFF | coeff, | ||
| EDGESELECTOR | edge_predicate = base::PredicateTrue{} ) -> MassEdgeMatrixProvider< typename PTR::element_type::Scalar, COEFF, EDGESELECTOR > |
| 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().
| auto operator* | ( | const A & | a, |
| const B & | b ) -> MeshFunctionBinary<internal::OperatorMultiplication, A, B> |
Multiply two mesh functions with each other.
| A | The type of the lhs MeshFunction |
| B | The type of the rhs MeshFunction |
| a | the lhs MeshFunction |
| b | the rhs MeshFunction |
a*b, i.e. a new mesh function which represents the pointwise product of a and b.A and B are Eigen::Matrix valued, the resulting MeshFunction will also be Matrix/Vector valued and will represent the Matrix product of a and b.lf::uscalfe namespace, it will not be found by Argument Dependent Lookup (ADL). You can get around this by explicitly importing the operator overload: using lf::uscalfe::operator*;Definition at line 722 of file mesh_function_binary.h.
| auto operator+ | ( | const A & | a, |
| const B & | b ) -> MeshFunctionBinary<internal::OperatorAddition, A, B> |
Add's two mesh functions.
| A | Type of the lhs MeshFunction |
| B | Type of the rhs MeshFunction |
| a | the lhs MeshFunction |
| b | the rhs MeshFunction |
a + b, i.e. a new mesh function which represents the pointwise addition of a and ba and b should produce the same type of values, e.g. both should be scalar valued or both matrix/vector/array valued.lf::uscalfe namespace, it will not be found by Argument Dependent Lookup. You can get around this by explictly importing the operator overload: using lf::uscalfe::operator+;Definition at line 664 of file mesh_function_binary.h.
| auto operator- | ( | const A & | a, |
| const B & | b ) -> MeshFunctionBinary<internal::OperatorSubtraction, A, B> |
Subtracts two mesh functions.
| A | Type of the lhs MeshFunction |
| B | Type of the rhs MeshFunction |
| a | the lhs MeshFunction |
| b | the rhs MeshFunction |
a - b, i.e. a new mesh function which represents the pointwise difference of a minus ba and b should produce the same type of values, e.g. both should be scalar valued or both matrix/vector valued.lf::uscalfe namespace, it will not be found by Argument Dependent Lookup (ADL). You can get around this by explicitly importing the operator overload: using lf::uscalfe::operator-;Definition at line 693 of file mesh_function_binary.h.
|
output operator for scalar parametric finite element space
Definition at line 445 of file uniform_scalar_fe_space.h.
| void lf::uscalfe::PrintInfo | ( | std::ostream & | o, |
| const UniformScalarFESpace< SCALAR > & | fes, | ||
| unsigned | ctrl ) |
Definition at line 420 of file uniform_scalar_fe_space.h.
References lf::base::RefEl::kPoint(), lf::base::RefEl::kQuad(), lf::base::RefEl::kSegment(), lf::base::RefEl::kTria(), kUniformScalarFESpaceOutDofh, kUniformScalarFESpaceOutMesh, kUniformScalarFESpaceOutRsfs, lf::uscalfe::UniformScalarFESpace< SCALAR >::LocGlobMap(), lf::uscalfe::UniformScalarFESpace< SCALAR >::Mesh(), lf::assemble::DofHandler::NumDofs(), and lf::uscalfe::UniformScalarFESpace< SCALAR >::NumRefShapeFunctions().
| lf::uscalfe::ReactionDiffusionElementMatrixProvider | ( | PTR | fe_space, |
| DIFF_COEFF | alpha, | ||
| REACTION_COEFF | gamma ) -> ReactionDiffusionElementMatrixProvider< typename PTR::element_type::Scalar, DIFF_COEFF, 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 > |
| std::shared_ptr< spdlog::logger > & lf::uscalfe::ReactionDiffusionElementMatrixProviderLogger | ( | ) |
logger for ReactionDiffusionElementMatrixProvider
Definition at line 21 of file loc_comp_ellbvp.cc.
References lf::base::InitLogger().
Referenced by lf::uscalfe::ReactionDiffusionElementMatrixProvider< SCALAR, DIFF_COEFF, REACTION_COEFF >::Eval().
| lf::uscalfe::ScalarLoadEdgeVectorProvider | ( | PTR | , |
| FUNCTOR | , | ||
| EDGESELECTOR | = base::PredicateTrue{} ) -> ScalarLoadEdgeVectorProvider< typename PTR::element_type::Scalar, FUNCTOR, EDGESELECTOR > |
| 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().
| lf::uscalfe::ScalarLoadElementVectorProvider | ( | PTR | fe_space, |
| MESH_FUNCTION | mf ) -> ScalarLoadElementVectorProvider< typename PTR::element_type::Scalar, MESH_FUNCTION > |
| std::shared_ptr< spdlog::logger > & lf::uscalfe::ScalarLoadElementVectorProviderLogger | ( | ) |
logger used by ScalarLoadElementVectorProvider
Definition at line 33 of file loc_comp_ellbvp.cc.
References lf::base::InitLogger().
Referenced by lf::uscalfe::ScalarLoadElementVectorProvider< SCALAR, MESH_FUNCTION >::Eval().
|
Pointwise squared norm of another mesh function.
| A | The type of the wrapped mesh function. |
| a | The original mesh function |
|a|^2 (pointwise)a to be either scalar or (eigen-) vector valued. Definition at line 278 of file mesh_function_unary.h.
|
Pointwise transpose of an Eigen::Matrix or Eigen::Array
| A | The type of the wrapped mesh function. |
| a | The original MeshFunction |
a. Definition at line 291 of file mesh_function_unary.h.
Referenced by lf::uscalfe::LinearFELaplaceElementMatrix::Eval().
| const FeLagrangeO2Segment<SCALAR> lf::uscalfe::FeLagrangeO2Quad< SCALAR >::krsf_segment_ |
| const Eigen::Matrix<int, 9, 2> lf::uscalfe::FeLagrangeO2Quad< SCALAR >::ksegment_to_quad_mapping_ |
| const FeLagrangeO3Segment<SCALAR> lf::uscalfe::FeLagrangeO3Quad< SCALAR >::krsf_segment_ |
| const Eigen::Matrix<int, 16, 2> lf::uscalfe::FeLagrangeO3Quad< SCALAR >::ksegment_to_quad_mapping_ |
| 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().
| 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().
| 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().