bem::LaplaceDL< F > Class Template Referenceabstract

Bilinear form to compute the Laplace double layer potential. More...

#include <bform.hh>

Inheritance diagram for bem::LaplaceDL< F >:
concepts::BilinearForm< concepts::Real > concepts::Cloneable concepts::OutputOperator

Public Member Functions

virtual LaplaceDLclone () const
 Virtual constructor. More...
 
virtual BilinearForm * clone () const=0
 Virtual constructor. More...
 
 LaplaceDL (uint stroud=0, uint gauss=0, concepts::Real dist=0.0)
 Constructor. More...
 
void operator() (const concepts::Element< F > &elmX, const concepts::Element< F > &elmY, concepts::ElementMatrix< F > &em)
 Application operator. More...
 
void operator() (const Constant3d000< F > &elmX, const Constant3d000< F > &elmY, concepts::ElementMatrix< F > &em)
 
void operator() (const Constant3d002< F > &elmX, const Constant3d002< F > &elmY, concepts::ElementMatrix< F > &em)
 
void operator() (const Dirac3d000< F > &elmX, const Linear3d000< F > &elmY, concepts::ElementMatrix< F > &em)
 
virtual void operator() (const Element< typename Realtype< concepts::Real >::type > &elmX, const Element< typename Realtype< concepts::Real >::type > &elmY, ElementMatrix< concepts::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< concepts::Real >::type > &elmX, const Element< typename Realtype< concepts::Real >::type > &elmY, ElementMatrix< concepts::Real > &em, const ElementPair< typename Realtype< concepts::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 Linear3d000< F > &elmX, const Linear3d000< F > &elmY, concepts::ElementMatrix< F > &em)
 

Protected Member Functions

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

Private Attributes

concepts::Real dist_
 Distance to distinguish between one point integration formula and integration using the given number of integration points. More...
 
uint gauss_
 
LplCol007< F > qrA_
 Classes used for integration. More...
 
LplGal007< F > qrB_
 
LplGal014< F > qrC_
 
uint stroud_
 Number of integration points. More...
 

Detailed Description

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

Bilinear form to compute the Laplace double layer potential.

Parameters
FField (Real or Cmplx)

Definition at line 147 of file bform.hh.

Constructor & Destructor Documentation

◆ LaplaceDL()

template<class F >
bem::LaplaceDL< F >::LaplaceDL ( 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 198 of file bform.hh.

Member Function Documentation

◆ clone() [1/2]

template<class F = concepts::Real>
virtual LaplaceDL* bem::LaplaceDL< 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 193 of file bform.hh.

◆ clone() [2/2]

virtual BilinearForm* concepts::BilinearForm< concepts::Real , typename Realtype<concepts::Real >::type >::clone ( ) const
pure virtualinherited

Virtual constructor.

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

◆ info()

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

Returns information in an output stream.

Reimplemented from concepts::OutputOperator.

◆ operator()() [1/7]

template<class F = concepts::Real>
void bem::LaplaceDL< 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/7]

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

Definition at line 219 of file bform.hh.

◆ operator()() [3/7]

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

Definition at line 230 of file bform.hh.

◆ operator()() [4/7]

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

Definition at line 203 of file bform.hh.

◆ operator()() [5/7]

virtual void concepts::BilinearForm< concepts::Real , typename Realtype<concepts::Real >::type >::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

◆ operator()() [6/7]

virtual void concepts::BilinearForm< concepts::Real , typename Realtype<concepts::Real >::type >::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

Definition at line 57 of file bilinearForm.hh.

◆ operator()() [7/7]

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

Definition at line 241 of file bform.hh.

Member Data Documentation

◆ dist_

template<class F = concepts::Real>
concepts::Real bem::LaplaceDL< F >::dist_
private

Distance to distinguish between one point integration formula and integration using the given number of integration points.

Definition at line 159 of file bform.hh.

◆ gauss_

template<class F = concepts::Real>
uint bem::LaplaceDL< F >::gauss_
private

Definition at line 155 of file bform.hh.

◆ qrA_

template<class F = concepts::Real>
LplCol007<F> bem::LaplaceDL< F >::qrA_
private

Classes used for integration.

Definition at line 150 of file bform.hh.

◆ qrB_

template<class F = concepts::Real>
LplGal007<F> bem::LaplaceDL< F >::qrB_
private

Definition at line 151 of file bform.hh.

◆ qrC_

template<class F = concepts::Real>
LplGal014<F> bem::LaplaceDL< F >::qrC_
private

Definition at line 152 of file bform.hh.

◆ stroud_

template<class F = concepts::Real>
uint bem::LaplaceDL< F >::stroud_
private

Number of integration points.

Definition at line 155 of file bform.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