linearFEM::Laplace2d Class Referenceabstract
Discrete equivalent of the Laplacian in 2D for linear FEM. More...
#include <bilinearForm2D.hh>
Public Member Functions | |
virtual Laplace2d * | clone () const |
Virtual constructor. More... | |
virtual BilinearForm * | clone () const=0 |
Virtual constructor. More... | |
virtual void | operator() (const concepts::Element< Real > &elmX, const concepts::Element< Real > &elmY, concepts::ElementMatrix< Real > &em) const |
virtual void | operator() (const Element< typename Realtype< Real >::type > &elmX, const Element< typename Realtype< Real >::type > &elmY, ElementMatrix< Real > &em) const=0 |
Evaluates the bilinear form for all shape functions on elmX and elmY and stores the result in the matrix em . More... | |
virtual void | operator() (const Element< typename Realtype< Real >::type > &elmX, const Element< typename Realtype< Real >::type > &elmY, ElementMatrix< Real > &em, const ElementPair< typename Realtype< Real >::type > &ep) const |
Evaluates the bilinear form for all shape functions on elmX and elmY and stores the result in the matrix em . More... | |
void | operator() (const linearFEM::Quad &elmX, const linearFEM::Quad &elmY, concepts::ElementMatrix< Real > &em) const |
Computes the element stiffness matrix for a quadrilateral. More... | |
void | operator() (const linearFEM::Triangle &elmX, const linearFEM::Triangle &elmY, concepts::ElementMatrix< Real > &em) const |
Computes the element stiffness matrix for a triangle. More... | |
Protected Member Functions | |
virtual std::ostream & | info (std::ostream &os) const |
Returns information in an output stream. More... | |
Private Attributes | |
concepts::Array< concepts::Real2d > | shpfct_ |
Storage area for values of derivatives of shape functions. More... | |
Detailed Description
Discrete equivalent of the Laplacian in 2D for linear FEM.
This bilinear form computes the stiffness matrix resulting from the discretization of the Laplacian (with integration by parts):
for the element shape functions .
Definition at line 30 of file bilinearForm2D.hh.
Member Function Documentation
◆ clone() [1/2]
|
inlinevirtual |
Virtual constructor.
Returns a pointer to a copy of itself. The caller is responsible to destroy this copy.
Implements concepts::Cloneable.
Definition at line 48 of file bilinearForm2D.hh.
◆ clone() [2/2]
|
pure virtualinherited |
Virtual constructor.
Returns a pointer to a copy of itself. The caller is responsible to destroy this copy.
◆ info()
|
protectedvirtualinherited |
Returns information in an output stream.
Reimplemented from concepts::OutputOperator.
Reimplemented in hp3D::RotRot, hp3D::Hook, hp3D::DivDiv< Weight >, hp3D::Identity< F >, hp3D::Laplace< F >, hp3D::BilinearFormTwoPartDeriv< F >, hp2Dedge::EdgeIdentity, hp2Dedge::Rotuv, hp2Dedge::RotRot< F >, hp2Dedge::IdentityMatrix< F >, hp2Dedge::Identity< F >, hp2Dedge::GraduvMatrix< F >, hp2Dedge::Graduv< F >, hp2D::RotRot, hp2D::DivDiv< Weight >, hp2D::BilinearFormTwoPartDeriv< F >, hp2D::BilinearFormOnePartDeriv< F >, hp2D::LaplaceMatrix< F >, and hp2D::Laplace< F >.
◆ operator()() [1/5]
|
virtual |
◆ operator()() [2/5]
|
pure virtualinherited |
Evaluates the bilinear form for all shape functions on elmX
and elmY
and stores the result in the matrix em
.
- Postcondition
- The returned matrix
em
has the correct size.
- Parameters
-
elmX Left element (test functions) elmY Right element (trial functions) em Return element matrix
◆ operator()() [3/5]
|
inlinevirtualinherited |
Evaluates the bilinear form for all shape functions on elmX
and elmY
and stores the result in the matrix em
.
If this method is not reimplemented in a derived class, the default behaviour is to call the application operator without ep
.
- Postcondition
- The returned matrix
em
has the correct size.
- Parameters
-
elmX Left element elmY Right element em Return element matrix ep Element pair holding more information on the pair elmX
andelmY
Definition at line 57 of file bilinearForm.hh.
◆ operator()() [4/5]
void linearFEM::Laplace2d::operator() | ( | const linearFEM::Quad & | elmX, |
const linearFEM::Quad & | elmY, | ||
concepts::ElementMatrix< Real > & | em | ||
) | const |
Computes the element stiffness matrix for a quadrilateral.
The stiffness matrix has to be integrated numericaly since the Jacobian of the element map is not constant.
◆ operator()() [5/5]
void linearFEM::Laplace2d::operator() | ( | const linearFEM::Triangle & | elmX, |
const linearFEM::Triangle & | elmY, | ||
concepts::ElementMatrix< Real > & | em | ||
) | const |
Computes the element stiffness matrix for a triangle.
The stiffness matrix is precomputed since the element map is constant.
Member Data Documentation
◆ shpfct_
|
mutableprivate |
Storage area for values of derivatives of shape functions.
Definition at line 51 of file bilinearForm2D.hh.
The documentation for this class was generated from the following file:
- linearFEM/bilinearForm2D.hh