hp1D::BiLaplace< F > Class Template Referenceabstract

A function class to calculate element matrices for the Bi-Laplacian. More...

#include <bilinearForm.hh>

Inheritance diagram for hp1D::BiLaplace< F >:
concepts::BilinearForm< Real, Real > hp1D::BilinearFormHelper< 2, 2, Real > concepts::Cloneable concepts::OutputOperator

Public Member Functions

 BiLaplace (const concepts::ElementFormulaContainer< F > frm=concepts::ElementFormulaContainer< F >())
 Constructor. More...
 
virtual BiLaplace< F > * clone () const
 Virtual constructor. More...
 
virtual BilinearForm * clone () const=0
 Virtual constructor. More...
 
void operator() (const concepts::Element< Real > &elmX, const concepts::Element< Real > &elmY, concepts::ElementMatrix< F > &em) const
 
virtual void operator() (const Element< Real > &elmX, const Element< Real > &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< Real > &elmX, const Element< Real > &elmY, ElementMatrix< Real > &em, const ElementPair< Real > &ep) const
 Evaluates the bilinear form for all shape functions on elmX and elmY and stores the result in the matrix em. More...
 
virtual ~BiLaplace ()
 Destructor. More...
 

Protected Member Functions

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

Detailed Description

template<class F = Real>
class hp1D::BiLaplace< F >

A function class to calculate element matrices for the Bi-Laplacian.

Author
Christian Heier, 2012

Definition at line 133 of file bilinearForm.hh.

Constructor & Destructor Documentation

◆ BiLaplace()

template<class F = Real>
hp1D::BiLaplace< F >::BiLaplace ( const concepts::ElementFormulaContainer< F >  frm = concepts::ElementFormulaContainer<F>())
inline

Constructor.

The formula frm is evaluated in each quadrature point.

Definition at line 139 of file bilinearForm.hh.

◆ ~BiLaplace()

template<class F = Real>
virtual hp1D::BiLaplace< F >::~BiLaplace ( )
virtual

Destructor.

Member Function Documentation

◆ clone() [1/2]

template<class F = Real>
virtual BiLaplace<F>* hp1D::BiLaplace< F >::clone ( ) const
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 150 of file bilinearForm.hh.

◆ clone() [2/2]

virtual BilinearForm* concepts::BilinearForm< Real , Real >::clone
pure virtualinherited

Virtual constructor.

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

◆ info()

template<class F = Real>
virtual std::ostream& hp1D::BiLaplace< F >::info ( std::ostream &  os) const
protectedvirtual

Returns information in an output stream.

Reimplemented from concepts::BilinearForm< Real, Real >.

◆ operator()() [1/3]

template<class F = Real>
void hp1D::BiLaplace< F >::operator() ( const concepts::Element< Real > &  elmX,
const concepts::Element< Real > &  elmY,
concepts::ElementMatrix< F > &  em 
) const

◆ operator()() [2/3]

virtual void concepts::BilinearForm< Real , Real >::operator() ( const Element< Real > &  elmX,
const Element< Real > &  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()() [3/3]

virtual void concepts::BilinearForm< Real , Real >::operator() ( const Element< Real > &  elmX,
const Element< Real > &  elmY,
ElementMatrix< Real > &  em,
const ElementPair< Real > &  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.


The documentation for this class was generated from the following file:
Page URL: http://wiki.math.ethz.ch/bin/view/Concepts/WebHome
21 August 2020
© 2020 Eidgenössische Technische Hochschule Zürich