Class documentation of Concepts

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

#include <bformAdapt.hh>

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

Public Member Functions

 AdaptLaplaceDL01 (uint gaussX=2, uint gaussY=2, concepts::Real deltaG=1.0, concepts::Real dist=0.0, concepts::Real eta=0.5)
 
void operator() (const concepts::Element< F > &elmX, const concepts::Element< F > &elmY, concepts::ElementMatrix< F > &em)
 
void operator() (const Constant3d001< F > &elmX, const Constant3d001< F > &elmY, concepts::ElementMatrix< F > &em)
 
virtual AdaptLaplaceDL01clone () 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>
class bem::AdaptLaplaceDL01< F >

Bilinear form for the Laplace double layer potential with piecewise constant shape functions and hanging nodes (=> recursive subdivision of the larger triangle). Number of integration points depends on the level of the element. The number of integration points given in the constructor is used for the highest level $L$. For an element on level $l$ the number of integration points $n_I$ in each dimension is given by $n_I = {\rm gauss} + (L-l)*\Delta_G$.

Parameters
FField (Real or Cmplx)

Definition at line 153 of file bformAdapt.hh.

Constructor & Destructor Documentation

◆ AdaptLaplaceDL01()

template<class F >
bem::AdaptLaplaceDL01< F >::AdaptLaplaceDL01 ( uint  gaussX = 2,
uint  gaussY = 2,
concepts::Real  deltaG = 1.0,
concepts::Real  dist = 0.0,
concepts::Real  eta = 0.5 
)
inline

3-dim quadrature (3-d Gauss)

Parameters
gaussXnumber of Gauss points in the integration formula for the outer ( $x$) integration.
gaussYnumber of Gauss points in the integration formula for the inner ( $y$) integration and the case where the two panels intersect.
deltaGAdditional number of Gauss points per level
distdistance to distinguish between near/far field (far field if the distance between the two centers of the panels is larger than dist => 1-point formular)
etaratio between element size and distance for recursive subdivision ( $\eta \sim [1/2, 1]$)

Definition at line 232 of file bformAdapt.hh.

Member Function Documentation

◆ clone()

template<class F >
virtual AdaptLaplaceDL01 * bem::AdaptLaplaceDL01< 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 186 of file bformAdapt.hh.

◆ info()

◆ operator()() [1/3]

template<class F >
void bem::AdaptLaplaceDL01< 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/3]

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

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.


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