Class documentation of Concepts

Loading...
Searching...
No Matches
bem::LaplaceSL< F > Class Template Referenceabstract

#include <bform.hh>

Inheritance diagram for bem::LaplaceSL< F >:
concepts::BilinearForm< F, G > concepts::Cloneable concepts::OutputOperator

Public Member Functions

 LaplaceSL (uint stroud=0, uint gauss=0, concepts::Real dist=0.0)
 
void operator() (const concepts::Element< F > &elmX, const concepts::Element< F > &elmY, concepts::ElementMatrix< F > &em)
 
void operator() (const Dirac3d000< F > &elmX, const Linear3d000< F > &elmY, concepts::ElementMatrix< F > &em)
 
void operator() (const Constant3d000< F > &elmX, const Constant3d000< F > &elmY, concepts::ElementMatrix< F > &em)
 
void operator() (const Constant3d001< F > &elmX, const Constant3d001< F > &elmY, concepts::ElementMatrix< F > &em)
 
void operator() (const Constant3d002< F > &elmX, const Constant3d002< F > &elmY, concepts::ElementMatrix< F > &em)
 
void operator() (const Linear3d000< F > &elmX, const Linear3d000< F > &elmY, concepts::ElementMatrix< F > &em)
 
virtual LaplaceSLclone () const
 
virtual void operator() (const Element< G > &elmX, const Element< G > &elmY, ElementMatrix< F > &em) const =0
 
virtual void operator() (const Element< G > &elmX, const Element< G > &elmY, ElementMatrix< F > &em, const ElementPair< G > &ep) const
 

Protected Member Functions

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

Detailed Description

template<class F = concepts::Real>
class bem::LaplaceSL< F >

Bilinear form to compute the Laplace single layer potential.

Parameters
FField (Real or Cmplx)

Definition at line 21 of file bform.hh.

Constructor & Destructor Documentation

◆ LaplaceSL()

template<class F >
bem::LaplaceSL< F >::LaplaceSL ( uint  stroud = 0,
uint  gauss = 0,
concepts::Real  dist = 0.0 
)
inline

Constructor.

Parameters
stroudNumber of integration points
gaussNumber of integration points
distMinimal distance to apply one point quadrature rule.

Definition at line 75 of file bform.hh.

Member Function Documentation

◆ clone()

template<class F = concepts::Real>
virtual LaplaceSL * bem::LaplaceSL< F >::clone ( ) const
inlinevirtual

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

Implements concepts::BilinearForm< F, G >.

Definition at line 70 of file bform.hh.

◆ info()

◆ operator()() [1/8]

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

Application operator.

Exceptions
MissingFeature
Parameters
elmXElement
elmYElement
emElement matrix for the two given elements.

◆ operator()() [2/8]

template<class F >
void bem::LaplaceSL< F >::operator() ( const Constant3d000< F > &  elmX,
const Constant3d000< F > &  elmY,
concepts::ElementMatrix< F > &  em 
)
inline

Definition at line 96 of file bform.hh.

◆ operator()() [3/8]

template<class F >
void bem::LaplaceSL< F >::operator() ( const Constant3d001< F > &  elmX,
const Constant3d001< F > &  elmY,
concepts::ElementMatrix< F > &  em 
)
inline

Definition at line 107 of file bform.hh.

◆ operator()() [4/8]

template<class F >
void bem::LaplaceSL< F >::operator() ( const Constant3d002< F > &  elmX,
const Constant3d002< F > &  elmY,
concepts::ElementMatrix< F > &  em 
)
inline

Definition at line 118 of file bform.hh.

◆ operator()() [5/8]

template<class F >
void bem::LaplaceSL< F >::operator() ( const Dirac3d000< F > &  elmX,
const Linear3d000< F > &  elmY,
concepts::ElementMatrix< F > &  em 
)
inline

Definition at line 80 of file bform.hh.

◆ operator()() [6/8]

template<class F , class G = typename Realtype<F>::type>
virtual void concepts::BilinearForm< F, G >::operator() ( const Element< G > &  elmX,
const Element< G > &  elmY,
ElementMatrix< F > &  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

Implemented in vectorial::BilinearForm< F, G >, concepts::BilinearFormLiCo< F, G >, concepts::BilinearFormContainer< F, G >, concepts::BilinearF_Sum< F, H, J, G >, and concepts::BilinearF_W< F, H, J, G >.

◆ operator()() [7/8]

template<class F , class G = typename Realtype<F>::type>
virtual void concepts::BilinearForm< F, G >::operator() ( const Element< G > &  elmX,
const Element< G > &  elmY,
ElementMatrix< F > &  em,
const ElementPair< G > &  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

Reimplemented in vectorial::BilinearForm< F, G >.

Definition at line 57 of file bilinearForm.hh.

◆ operator()() [8/8]

template<class F >
void bem::LaplaceSL< F >::operator() ( const Linear3d000< F > &  elmX,
const Linear3d000< F > &  elmY,
concepts::ElementMatrix< F > &  em 
)
inline

Definition at line 129 of file bform.hh.


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