bem::AdaptLaplaceDL01< 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>
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) | |
3-dim quadrature (3-d Gauss) More... | |
virtual AdaptLaplaceDL01 * | clone () 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::Real3d * | c_ |
Space for some intermediate values. More... | |
concepts::Real | deltaG_ |
Additional number of Integration points per level. 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 | gaussX_ |
Number of Gauss points. More... | |
uint | gaussY_ |
LplGal011< F > | qr_ |
Integration. More... | |
concepts::Real * | r_ |
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 . For an element on level the number of integration points in each dimension is given by .
- Parameters
-
F Field (Real or Cmplx)
Definition at line 153 of file bformAdapt.hh.
Constructor & Destructor Documentation
◆ AdaptLaplaceDL01()
|
inline |
3-dim quadrature (3-d Gauss)
- Parameters
-
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.
Member Function Documentation
◆ admissible_()
|
inlineprivate |
Admissible condition of two panels for integration.
- Parameters
-
cX center of the first panel rX radius of the first panel cY center of the second panel rY radius of the second panel
Definition at line 242 of file bformAdapt.hh.
◆ ball_() [1/2]
|
private |
◆ ball_() [2/2]
|
private |
Computes the center and radius of a panels.
◆ clone() [1/2]
|
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 186 of file bformAdapt.hh.
◆ clone() [2/2]
|
pure virtualinherited |
Virtual constructor.
Returns a pointer to a copy of itself. The caller is responsible to destroy this copy.
◆ info()
|
protectedvirtualinherited |
Returns information in an output stream.
Reimplemented from concepts::OutputOperator.
Reimplemented in hp3D::LinearElasticity< F >, and constraints::ConstraintsList< F >.
◆ operator()() [1/4]
void bem::AdaptLaplaceDL01< F >::operator() | ( | const concepts::Element< F > & | elmX, |
const concepts::Element< F > & | elmY, | ||
concepts::ElementMatrix< F > & | em | ||
) |
Application operator.
- Exceptions
-
MissingFeature
- Parameters
-
elmX Element elmY Element em Element matrix for the two given elements.
◆ operator()() [2/4]
void bem::AdaptLaplaceDL01< F >::operator() | ( | const Constant3d001< F > & | elmX, |
const Constant3d001< F > & | elmY, | ||
concepts::ElementMatrix< F > & | em | ||
) |
◆ operator()() [3/4]
|
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
-
elmX Left element (test functions) elmY Right element (trial functions) em Return element matrix
◆ operator()() [4/4]
|
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
-
elmX Left element elmY Right element em Return element matrix ep Element pair holding more information on the pair elmX
andelmY
Definition at line 57 of file bilinearForm.hh.
◆ update_()
|
private |
Recursive subdivison of the panel on the higher level.
- Parameters
-
deltal level difference between cell and elm m element matrix (update_ just sums the all parts up => set m to 0.0 before calling update_ for the first time!) swp swap cell and elm for integration
Member Data Documentation
◆ c_
|
private |
Space for some intermediate values.
Definition at line 227 of file bformAdapt.hh.
◆ deltaG_
|
private |
Additional number of Integration points per level.
Definition at line 217 of file bformAdapt.hh.
◆ dist_
|
private |
Distance to distinguish between near/far field.
Definition at line 219 of file bformAdapt.hh.
◆ eta_
|
private |
Ratio between element size and distance for recursive subdivision.
Definition at line 221 of file bformAdapt.hh.
◆ gaussX_
|
private |
Number of Gauss points.
Definition at line 214 of file bformAdapt.hh.
◆ gaussY_
|
private |
Definition at line 215 of file bformAdapt.hh.
◆ qr_
|
private |
Integration.
Definition at line 224 of file bformAdapt.hh.
◆ r_
|
private |
Definition at line 228 of file bformAdapt.hh.
The documentation for this class was generated from the following file:
- bem/bformAdapt.hh