bformAdapt.hh

Go to the documentation of this file.
1 
6 #ifndef bformAdapt_hh
7 #define bformAdapt_hh
8 
10 #include "bemInt.hh"
11 #include "bem/element.hh"
12 
13 // debugging bilinear forms
14 #define AdaptLaplaceDL00_D 0
15 
16 namespace bem {
17 
18  // ****************************************************** AdaptLaplaceDL00 **
19 
25  template <class F>
27  public:
37  inline AdaptLaplaceDL00(uint stroud = 2, uint gauss = 2,
38  concepts::Real dist = 0.0,
39  concepts::Real eta = 0.5);
40 #if AdaptLaplaceDL00_D
41  virtual ~AdaptLaplaceDL00() {
42  std::cout << "AdaptLaplaceDL00(samel=" << samel_ << ", far=" << far_;
43  std::cout << ", recX=" << recX_ << ", recY=" << recY_;
44  std::cout << ", rec=" << rec_ << ", farrec=" << farrec_ << ")\n";
45  };
46 #endif
47 
54  void operator()(const concepts::Element<F>& elmX,
55  const concepts::Element<F>& elmY,
57  void operator()(const Constant3d001<F>& elmX,
58  const Constant3d001<F>& elmY,
60 
61  virtual AdaptLaplaceDL00* clone() const {
62  return new AdaptLaplaceDL00(stroud_, gauss_, dist_, eta_); }
63  private:
66  concepts::Real* r) const;
68  concepts::Real* r) const;
75  inline bool admissible_(const concepts::Real3d* cX,
76  const concepts::Real* rX,
77  const concepts::Real3d* cY,
78  const concepts::Real* rY) const;
86  uint deltal, F* m, bool swp) const;
87 
89  uint stroud_;
91  uint gauss_;
98 
102 
103 #if AdaptLaplaceDL00_D
104  // Debug information
105  // # integration (on the same level)
106  uint samel_;
107  // # integration (elements far away from each other)
108  uint far_;
109  // # integration (elements far away in recursive call)
110  mutable uint farrec_;
111  // # integration (recursive with \f$l_X > l_Y\f$)
112  uint recX_;
113  // # integration (recursive with \f$l_Y > l_X\f$)
114  uint recY_;
115  // # recursive integration calls
116  mutable uint rec_;
117 #endif
118  };
119 
120  template <class F>
121  inline AdaptLaplaceDL00<F>::AdaptLaplaceDL00(uint stroud, uint gauss,
122  concepts::Real dist,
123  concepts::Real eta)
124  : stroud_(stroud), gauss_(gauss), dist_(dist), eta_(eta), c_(0), r_(0) {
125 #if AdaptLaplaceDL00_D
126  samel_ = 0; far_ = 0; recX_ = 0; recY_ = 0; rec_ = 0; farrec_ = 0;
127 #endif
128  }
129 
130  template <class F>
131  inline bool
133  const concepts::Real* rX,
134  const concepts::Real3d* cY,
135  const concepts::Real* rY) const {
136  return concepts::Real3d(*cX - *cY).l2() * eta_ > *rX + *rY;
137  }
138 
139  // ****************************************************** AdaptLaplaceDL01 **
140 
152  template <class F>
154  public:
168  inline AdaptLaplaceDL01(uint gaussX = 2, uint gaussY = 2,
169  concepts::Real deltaG = 1.0,
170  concepts::Real dist = 0.0,
171  concepts::Real eta = 0.5);
172 
180  const concepts::Element<F>& elmY,
182  void operator()(const Constant3d001<F>& elmX,
183  const Constant3d001<F>& elmY,
185 
186  virtual AdaptLaplaceDL01* clone() const {
187  return new AdaptLaplaceDL01(gaussX_, gaussY_, deltaG_, dist_, eta_); }
188  private:
191  concepts::Real* r) const;
193  concepts::Real* r) const;
200  inline bool admissible_(const concepts::Real3d* cX,
201  const concepts::Real* rX,
202  const concepts::Real3d* cY,
203  const concepts::Real* rY) const;
211  uint deltal, F* m, bool swp) const;
212 
214  uint gaussX_;
215  uint gaussY_;
222 
225 
229  };
230 
231  template <class F>
232  inline AdaptLaplaceDL01<F>::AdaptLaplaceDL01(uint gaussX, uint gaussY,
233  concepts::Real deltaG,
234  concepts::Real dist,
235  concepts::Real eta)
236  : gaussX_(gaussX), gaussY_(gaussY), deltaG_(deltaG), dist_(dist),
237  eta_(eta), c_(0), r_(0) {
238  }
239 
240  template <class F>
241  inline bool
243  const concepts::Real* rX,
244  const concepts::Real3d* cY,
245  const concepts::Real* rY) const {
246  return concepts::Real3d(*cX - *cY).l2() * eta_ > *rX + *rY;
247  }
248 
249  // ****************************************************** AdaptLaplaceSL01 **
250 
262  template <class F>
264  public:
278  inline AdaptLaplaceSL01(uint stroud = 0, uint gauss = 0,
279  concepts::Real deltaI = 1.0,
280  concepts::Real dist = 0.0,
281  concepts::Real eta = 0.5);
282 
290  const concepts::Element<F>& elmY,
292  void operator()(const Constant3d001<F>& elmX,
293  const Constant3d001<F>& elmY,
295 
296  virtual AdaptLaplaceSL01* clone() const {
297  return new AdaptLaplaceSL01(stroud_, gauss_, deltaI_, dist_, eta_); }
298  private:
301  concepts::Real* r) const;
303  concepts::Real* r) const;
310  inline bool admissible_(const concepts::Real3d* cX,
311  const concepts::Real* rX,
312  const concepts::Real3d* cY,
313  const concepts::Real* rY) const;
321  uint deltal, F* m, bool swp) const;
322 
324  uint stroud_;
326  uint gauss_;
333 
336 
340  };
341 
342  template <class F>
343  inline AdaptLaplaceSL01<F>::AdaptLaplaceSL01(uint stroud, uint gauss,
344  concepts::Real deltaI,
345  concepts::Real dist,
346  concepts::Real eta)
347  : stroud_(stroud), gauss_(gauss), deltaI_(deltaI), dist_(dist),
348  eta_(eta), c_(0), r_(0) {
349  }
350 
351  template <class F>
352  inline bool
354  const concepts::Real* rX,
355  const concepts::Real3d* cY,
356  const concepts::Real* rY) const {
357  return concepts::Real3d(*cX - *cY).l2() * eta_ > *rX + *rY;
358  }
359 
360 } // namespace bem
361 
362 #endif // bformAdapt_hh
uint gauss_
Number of Integration points (inner $y$ integration)
Definition: bformAdapt.hh:326
virtual AdaptLaplaceDL01 * clone() const
Virtual constructor.
Definition: bformAdapt.hh:186
void operator()(const Constant3d001< F > &elmX, const Constant3d001< F > &elmY, concepts::ElementMatrix< F > &em)
Evaluation of the Laplace single layer potential with constant test/trial functions.
Definition: lplGal018.hh:28
AdaptLaplaceDL00(uint stroud=2, uint gauss=2, concepts::Real dist=0.0, concepts::Real eta=0.5)
3-dim quadrature (Stroud + 1-d Gauss)
Definition: bformAdapt.hh:121
An abstract class for an element of a space.
Definition: exceptions.hh:15
concepts::Real dist_
Distance to distinguish between near/far field.
Definition: bformAdapt.hh:219
void operator()(const Constant3d001< F > &elmX, const Constant3d001< F > &elmY, concepts::ElementMatrix< F > &em)
virtual AdaptLaplaceDL00 * clone() const
Virtual constructor.
Definition: bformAdapt.hh:61
concepts::Real eta_
Ratio between element size and distance for recursive subdivision.
Definition: bformAdapt.hh:221
AdaptLaplaceSL01(uint stroud=0, uint gauss=0, concepts::Real deltaI=1.0, concepts::Real dist=0.0, concepts::Real eta=0.5)
3-dim quadrature (3-d Gauss)
Definition: bformAdapt.hh:343
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.
concepts::Real3d * c_
Space for some intermediate values.
Definition: bformAdapt.hh:227
uint stroud_
Number of Integration points (outer $x$ integration)
Definition: bformAdapt.hh:324
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.
Definition: bformAdapt.hh:132
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.
Definition: bformAdapt.hh:353
void ball_(const Constant3d001< F > &elm, concepts::Real3d *c, concepts::Real *r) const
Computes the center and radius of a panels.
Point< Real, 3 > Real3d
Bilinear form for the Laplace double layer potential with piecewise constant shape functions and hang...
Definition: bformAdapt.hh:26
LplGal018< F > qr_
Integration.
Definition: bformAdapt.hh:335
Abstract function class to evaluate a bilinear form.
Definition: bilinearForm.hh:33
LplGal014< F > qr_
Integration.
Definition: bformAdapt.hh:97
LplGal011< F > qr_
Integration.
Definition: bformAdapt.hh:224
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.
Definition: bformAdapt.hh:242
concepts::Real3d * c_
Space for some intermediate values.
Definition: bformAdapt.hh:100
Constant space element with a level dependent key.
Definition: element.hh:335
concepts::Real deltaI_
Additional number of Integration points per level.
Definition: bformAdapt.hh:328
Real l2() const
Returns the Euclidian norm of the vector.
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.
concepts::Real deltaG_
Additional number of Integration points per level.
Definition: bformAdapt.hh:217
void operator()(const concepts::Element< F > &elmX, const concepts::Element< F > &elmY, concepts::ElementMatrix< F > &em)
Application operator.
virtual AdaptLaplaceSL01 * clone() const
Virtual constructor.
Definition: bformAdapt.hh:296
void ball_(const concepts::Triangle3d &cell, concepts::Real3d *c, concepts::Real *r) const
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.
concepts::Real * r_
Definition: bformAdapt.hh:339
concepts::Real3d * c_
Space for some intermediate values.
Definition: bformAdapt.hh:338
concepts::Real eta_
Ratio between element size and distance for recursive subdivision.
Definition: bformAdapt.hh:332
Bilinear form for the Laplace double layer potential with piecewise constant shape functions and hang...
Definition: bformAdapt.hh:153
Used for the basic classes of the boundary element method.
Definition: bform.hh:13
void operator()(const Constant3d001< F > &elmX, const Constant3d001< F > &elmY, concepts::ElementMatrix< F > &em)
concepts::Real dist_
Distance to distinguish between near/far field.
Definition: bformAdapt.hh:93
concepts::Real eta_
Ratio between element size and distance for recursive subdivision.
Definition: bformAdapt.hh:95
concepts::Real dist_
Distance to distinguish between near/far field.
Definition: bformAdapt.hh:330
uint gauss_
Number of Gauss points.
Definition: bformAdapt.hh:91
void operator()(const concepts::Element< F > &elmX, const concepts::Element< F > &elmY, concepts::ElementMatrix< F > &em)
Application operator.
uint stroud_
Number of Stroud points.
Definition: bformAdapt.hh:89
void ball_(const concepts::Triangle3d &cell, concepts::Real3d *c, concepts::Real *r) const
Element matrix.
Definition: linearForm.hh:18
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.
void operator()(const concepts::Element< F > &elmX, const concepts::Element< F > &elmY, concepts::ElementMatrix< F > &em)
Application operator.
A 3D cell: triangle.
Definition: cell2D.hh:719
Evaluation of the Laplace double layer potential with constant test/trial functions.
Definition: lplGal014.hh:29
Bilinear form for the Laplace Single Layer potential with piecewise constant shape functions and hang...
Definition: bformAdapt.hh:263
uint gaussX_
Number of Gauss points.
Definition: bformAdapt.hh:214
concepts::Real * r_
Definition: bformAdapt.hh:228
concepts::Real * r_
Definition: bformAdapt.hh:101
Evaluation of the Laplace double layer potential with constant test/trial functions.
Definition: lplGal011.hh:29
double Real
Type normally used for a floating point number.
Definition: typedefs.hh:36
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)
Definition: bformAdapt.hh:232
void ball_(const Constant3d001< F > &elm, concepts::Real3d *c, concepts::Real *r) const
Computes the center and radius of a panels.
Page URL: http://wiki.math.ethz.ch/bin/view/Concepts/WebHome
21 August 2020
© 2020 Eidgenössische Technische Hochschule Zürich