#include <bformAdapt.hh>
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 AdaptLaplaceDL01 * | clone () 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. | |
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 . For an element on level the number of integration points in each dimension is given by .
F | Field (Real or Cmplx) |
Definition at line 153 of file bformAdapt.hh.
|
inline |
3-dim quadrature (3-d Gauss)
gaussX | number of Gauss points in the integration formula for the outer ( ) integration. |
gaussY | number of Gauss points in the integration formula for the inner ( ) integration and the case where the two panels intersect. |
deltaG | Additional number of Gauss points per level |
dist | distance 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) |
eta | ratio between element size and distance for recursive subdivision ( ) |
Definition at line 232 of file bformAdapt.hh.
|
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.
|
protectedvirtualinherited |
Returns information in an output stream.
Reimplemented from concepts::OutputOperator.
Reimplemented in constraints::ConstraintsList< F >, hp1D::Laplace< F >, hp1D::Identity< F >, hp1D::IdentityParallel< F >, hp1D::BiLaplace< F >, hp1D::Jump1Jump1< F >, hp1D::Mean2Jump1< F >, hp1D::Jump1Mean2< F >, hp2D::Advection< F >, hp2D::Identity< F >, hp2D::Laplace< F >, hp2D::LaplaceMatrix< F >, hp2D::BilinearFormOnePartDeriv< F >, hp2D::BilinearFormTwoPartDeriv< F >, hp2D::DivDiv< Weight >, hp2D::RotRot, hp2Dedge::Graduv< F >, hp2Dedge::GraduvMatrix< F >, hp2Dedge::Identity< F >, hp2Dedge::IdentityMatrix< F >, hp2Dedge::RotRot< F >, hp2Dedge::Rotuv, hp2Dedge::EdgeIdentity, hp3D::LinearElasticity< F >, hp3D::BilinearFormTwoPartDeriv< F >, hp3D::Laplace< F >, hp3D::Identity< F >, hp3D::Advection< F >, hp3D::DivDiv< Weight >, hp3D::Hook, hp3D::RotRot, concepts::BilinearFormLiCo< F, G >, concepts::BilinearFormContainer< F, G >, concepts::BilinearFormContainer< F, typename Realtype< F >::type >, concepts::BilinearFormContainer< Real, Real >, concepts::BilinearF_Sum< F, H, J, G >, concepts::BilinearF_W< F, H, J, G >, and vectorial::BilinearForm< F, G >.
void bem::AdaptLaplaceDL01< F >::operator() | ( | const concepts::Element< F > & | elmX, |
const concepts::Element< F > & | elmY, | ||
concepts::ElementMatrix< F > & | em | ||
) |
Application operator.
MissingFeature |
elmX | Element |
elmY | Element |
em | Element matrix for the two given elements. |
|
pure virtualinherited |
Evaluates the bilinear form for all shape functions on elmX
and elmY
and stores the result in the matrix em
.
em
has the correct size. elmX | Left element (test functions) |
elmY | Right element (trial functions) |
em | Return 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 >.
|
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
.
em
has the correct size. elmX | Left element |
elmY | Right element |
em | Return element matrix |
ep | Element pair holding more information on the pair elmX and elmY |
Reimplemented in vectorial::BilinearForm< F, G >.
Definition at line 57 of file bilinearForm.hh.