Class documentation of Concepts

Loading...
Searching...
No Matches
linearFEM::Laplace2d Class Referenceabstract

#include <bilinearForm2D.hh>

Inheritance diagram for linearFEM::Laplace2d:
concepts::BilinearForm< Real > concepts::Cloneable concepts::OutputOperator

Public Member Functions

virtual void operator() (const concepts::Element< Real > &elmX, const concepts::Element< Real > &elmY, concepts::ElementMatrix< Real > &em) const
 
void operator() (const linearFEM::Quad &elmX, const linearFEM::Quad &elmY, concepts::ElementMatrix< Real > &em) const
 
void operator() (const linearFEM::Triangle &elmX, const linearFEM::Triangle &elmY, concepts::ElementMatrix< Real > &em) const
 
virtual Laplace2dclone () const
 
virtual void operator() (const Element< typename Realtype< Real >::type > &elmX, const Element< typename Realtype< Real >::type > &elmY, ElementMatrix< Real > &em) const=0
 
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
 

Protected Member Functions

virtual std::ostream & info (std::ostream &os) const
 Returns information in an output stream.
 

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):

\[ \int_K \nabla \phi_i \cdot \nabla \phi_j \, dx \]

for the element shape functions $\phi_i$.

Author
Philipp Frauenfelder, 2002

Definition at line 30 of file bilinearForm2D.hh.

Member Function Documentation

◆ clone()

virtual Laplace2d * linearFEM::Laplace2d::clone ( ) const
inlinevirtual

Virtual constructor. Returns a pointer to a copy of itself. The caller is responsible to destroy this copy.

Implements concepts::BilinearForm< Real >.

Definition at line 48 of file bilinearForm2D.hh.

◆ info()

virtual std::ostream & concepts::BilinearForm< Real , typename Realtype<Real >::type >::info ( std::ostream &  os) const
protectedvirtualinherited

Returns information in an output stream.

Reimplemented from concepts::OutputOperator.

Reimplemented in hp2D::DivDiv< Weight >, hp2D::RotRot, hp2Dedge::Rotuv, hp2Dedge::EdgeIdentity, hp3D::DivDiv< Weight >, hp3D::Hook, and hp3D::RotRot.

◆ operator()() [1/4]

virtual void concepts::BilinearForm< Real , typename Realtype<Real >::type >::operator() ( const Element< typename Realtype<Real >::type > &  elmX,
const Element< typename Realtype<Real >::type > &  elmY,
ElementMatrix< Real > &  em 
) const
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
elmXLeft element (test functions)
elmYRight element (trial functions)
emReturn element matrix

◆ operator()() [2/4]

virtual void concepts::BilinearForm< Real , typename Realtype<Real >::type >::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
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
elmXLeft element
elmYRight element
emReturn element matrix
epElement pair holding more information on the pair elmX and elmY

Definition at line 57 of file bilinearForm.hh.

◆ operator()() [3/4]

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()() [4/4]

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.


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