bem::AdaptLaplaceDL00< F > Class Template Referenceabstract

Bilinear form for the Laplace double layer potential with piecewise constant shape functions and hanging nodes (=> recursive subdivision of the larger triangle). More...

#include <bformAdapt.hh>

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

Public Member Functions

 AdaptLaplaceDL00 (uint stroud=2, uint gauss=2, concepts::Real dist=0.0, concepts::Real eta=0.5)
 3-dim quadrature (Stroud + 1-d Gauss) More...
 
virtual AdaptLaplaceDL00clone () const
 Virtual constructor. More...
 
virtual BilinearForm * clone () const=0
 Virtual constructor. More...
 
void operator() (const concepts::Element< F > &elmX, const concepts::Element< F > &elmY, concepts::ElementMatrix< F > &em)
 Application operator. More...
 
void operator() (const Constant3d001< F > &elmX, const Constant3d001< F > &elmY, concepts::ElementMatrix< F > &em)
 
virtual void operator() (const Element< typename Realtype< F >::type > &elmX, const Element< typename Realtype< F >::type > &elmY, ElementMatrix< F > &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< F >::type > &elmX, const Element< typename Realtype< F >::type > &elmY, ElementMatrix< F > &em, const ElementPair< typename Realtype< F >::type > &ep) const
 Evaluates the bilinear form for all shape functions on elmX and elmY and stores the result in the matrix em. More...
 

Protected Member Functions

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

Private Member Functions

bool admissible_ (const concepts::Real3d *cX, const concepts::Real *rX, const concepts::Real3d *cY, const concepts::Real *rY) const
 Admissible condition of two panels for integration. More...
 
void ball_ (const concepts::Triangle3d &cell, concepts::Real3d *c, concepts::Real *r) const
 
void ball_ (const Constant3d001< F > &elm, concepts::Real3d *c, concepts::Real *r) const
 Computes the center and radius of a panels. More...
 
void update_ (concepts::Triangle3d &cell, const Constant3d001< F > &elm, uint deltal, F *m, bool swp) const
 Recursive subdivison of the panel on the higher level. More...
 

Private Attributes

concepts::Real3dc_
 Space for some intermediate values. More...
 
concepts::Real dist_
 Distance to distinguish between near/far field. More...
 
concepts::Real eta_
 Ratio between element size and distance for recursive subdivision. More...
 
uint gauss_
 Number of Gauss points. More...
 
LplGal014< F > qr_
 Integration. More...
 
concepts::Realr_
 
uint stroud_
 Number of Stroud points. More...
 

Detailed Description

template<class F>
class bem::AdaptLaplaceDL00< F >

Bilinear form for the Laplace double layer potential with piecewise constant shape functions and hanging nodes (=> recursive subdivision of the larger triangle).

Parameters
FField (Real or Cmplx)

Definition at line 26 of file bformAdapt.hh.

Constructor & Destructor Documentation

◆ AdaptLaplaceDL00()

template<class F >
bem::AdaptLaplaceDL00< F >::AdaptLaplaceDL00 ( uint  stroud = 2,
uint  gauss = 2,
concepts::Real  dist = 0.0,
concepts::Real  eta = 0.5 
)
inline

3-dim quadrature (Stroud + 1-d Gauss)

Parameters
stroudnumber of Stroud points in the integration formula
gaussnumber of Gauss points in the integration formula
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, 2/3]$)

Definition at line 121 of file bformAdapt.hh.

Member Function Documentation

◆ admissible_()

template<class F >
bool bem::AdaptLaplaceDL00< F >::admissible_ ( const concepts::Real3d cX,
const concepts::Real rX,
const concepts::Real3d cY,
const concepts::Real rY 
) const
inlineprivate

Admissible condition of two panels for integration.

Parameters
cXcenter of the first panel
rXradius of the first panel
cYcenter of the second panel
rYradius of the second panel

Definition at line 132 of file bformAdapt.hh.

◆ ball_() [1/2]

template<class F >
void bem::AdaptLaplaceDL00< F >::ball_ ( const concepts::Triangle3d cell,
concepts::Real3d c,
concepts::Real r 
) const
private

◆ ball_() [2/2]

template<class F >
void bem::AdaptLaplaceDL00< F >::ball_ ( const Constant3d001< F > &  elm,
concepts::Real3d c,
concepts::Real r 
) const
private

Computes the center and radius of a panels.

◆ clone() [1/2]

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

◆ clone() [2/2]

virtual BilinearForm* concepts::BilinearForm< F, typename Realtype<F>::type >::clone
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< F, typename Realtype<F>::type >::info ( std::ostream &  os) const
protectedvirtualinherited

Returns information in an output stream.

Reimplemented from concepts::OutputOperator.

Reimplemented in hp3D::LinearElasticity< F >, and constraints::ConstraintsList< F >.

◆ operator()() [1/4]

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

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

◆ operator()() [3/4]

virtual void concepts::BilinearForm< F, typename Realtype<F>::type >::operator() ( const Element< typename Realtype<F>::type > &  elmX,
const Element< typename Realtype<F>::type > &  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()() [4/4]

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

◆ update_()

template<class F >
void bem::AdaptLaplaceDL00< F >::update_ ( concepts::Triangle3d cell,
const Constant3d001< F > &  elm,
uint  deltal,
F *  m,
bool  swp 
) const
private

Recursive subdivison of the panel on the higher level.

Parameters
deltallevel difference between cell and elm
melement matrix (update_ just sums the all parts up => set m to 0.0 before calling update_ for the first time!)
swpswap cell and elm for integration

Member Data Documentation

◆ c_

template<class F >
concepts::Real3d* bem::AdaptLaplaceDL00< F >::c_
private

Space for some intermediate values.

Definition at line 100 of file bformAdapt.hh.

◆ dist_

template<class F >
concepts::Real bem::AdaptLaplaceDL00< F >::dist_
private

Distance to distinguish between near/far field.

Definition at line 93 of file bformAdapt.hh.

◆ eta_

template<class F >
concepts::Real bem::AdaptLaplaceDL00< F >::eta_
private

Ratio between element size and distance for recursive subdivision.

Definition at line 95 of file bformAdapt.hh.

◆ gauss_

template<class F >
uint bem::AdaptLaplaceDL00< F >::gauss_
private

Number of Gauss points.

Definition at line 91 of file bformAdapt.hh.

◆ qr_

template<class F >
LplGal014<F> bem::AdaptLaplaceDL00< F >::qr_
private

Integration.

Definition at line 97 of file bformAdapt.hh.

◆ r_

template<class F >
concepts::Real* bem::AdaptLaplaceDL00< F >::r_
private

Definition at line 101 of file bformAdapt.hh.

◆ stroud_

template<class F >
uint bem::AdaptLaplaceDL00< F >::stroud_
private

Number of Stroud points.

Definition at line 89 of file bformAdapt.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