bform.hh

Go to the documentation of this file.
1 
6 #ifndef aglowav2bform_hh
7 #define aglowav2bform_hh
8 
9 #ifdef __GUNG__
10 #pragma interface
11 #endif
12 
13 #include "operator/bilinearForm.hh"
14 #include "aglowav2/xfy.hh"
15 
16 namespace aglowav2 {
17 
18  // *********************************************************** AglowavBF00 **
19 
28  template<class F = concepts::Real>
29  class AglowavBF00 : public concepts::BilinearForm<F> {
30  private:
32  typedef struct Nfld {
34  struct Nfld* lnk;
36  const Haar3d000<F>* elm;
38  F* valXY;
40  F* valYX;
41  } Nfld;
42 
44  struct Cell {
45  F* mtsttrl;
46  F* mtrltst;
47 
48  const Haar3d000<F>* elm;
49 
50  Cell() : mtsttrl(0), mtrltst(0), elm(0) {}
51  };
52 
53  virtual AglowavBF00* clone() const { return new AglowavBF00(bf_, XFY_); }
54  public:
62 
63  void operator()(const concepts::Element<F>& elmX,
64  const concepts::Element<F>& elmY,
66  void operator()(const Haar3d000<F>& elmX, const Haar3d000<F>& elmY,
68 
69  uint memory() const {return mem_;}
70 
71  private:
81  void nearfield_(const Haar3d000<F>& elmX, const Haar3d000<F>& elmY,
82  F* valXY, uint ldXY, F* valYX, uint ldYX);
83 
93  void updateChld_(const Haar3d000<F>& tst, const Haar3d000<F>& trl,
94  F* valtsttrl, uint ldtsttrl, F* valtrltst, uint ldtrltst);
95 
105  void updateElm_(const Haar3d000<F>& tst, const Haar3d000<F>& trl,
106  F* valtsttrl, uint ldtsttrl, F* valtrltst, uint ldtrltst);
107 
115  uint mem_;
116 
121 
123  F* valXY_;
124  F* vXY_;
125  F* valYX_;
126  F* vYX_;
127  F* foo_;
128 
134  bool reuse_;
135  };
136 
137  // *********************************************************** WavIdentity **
138 
142  template<class F = concepts::Real>
144  public:
146  const concepts::Element<F>& elmY,
148  void operator()(const Haar3d000<F>& elmX, const Haar3d000<F>& elmY,
150  virtual WavIdentity* clone() const { return new WavIdentity(); }
151  };
152 
153 } // namespace aglowav2
154 
155 #endif // aglowav2bform_hh
Cell * buf_
Buffer to store wavelet entries "for a short time".
Definition: bform.hh:132
Class for the far field matrix F.
Definition: xfy.hh:153
An abstract class for an element of a space.
Definition: exceptions.hh:15
Nfld *** nfldLast_
Array to hold pointers to the last list entry of the near field (used only by nearfield_)
Definition: bform.hh:120
AglowavBF00(concepts::BilinearForm< F > &bf, F00< F > &XFY)
Constructor.
Structure to hold the one scale near field one each level.
Definition: bform.hh:32
Bilinear form for the stiffness matrix compression with aglomerated wavelets.
Definition: bform.hh:29
F * valYX
Near field entries stored row wise.
Definition: bform.hh:40
const Haar3d000< F > * elm
Pointer to element in the near field.
Definition: bform.hh:36
virtual WavIdentity * clone() const
Virtual constructor.
Definition: bform.hh:150
Abstract function class to evaluate a bilinear form.
Definition: bilinearForm.hh:33
void nearfield_(const Haar3d000< F > &elmX, const Haar3d000< F > &elmY, F *valXY, uint ldXY, F *valYX, uint ldYX)
Computation of the one scale near field on each level (used by the constructor).
struct aglowav2::AglowavBF00::Nfld Nfld
Structure to hold the one scale near field one each level.
concepts::BilinearForm< F > & bf_
Bilinear form to compute the one scale matrix entries.
Definition: bform.hh:109
Nfld ** nfld_
Array to hold the near field for each element.
Definition: bform.hh:113
Used for the aglowav2 classes for the boundary element method.
Definition: bform.hh:16
uint mem_
memory
Definition: bform.hh:115
Space element for the agglomerated wavelets.
Definition: element.hh:346
bool reuse_
Entry reused and thus saved in cache_ otherwise stored in buf_.
Definition: bform.hh:134
void operator()(const Haar3d000< F > &elmX, const Haar3d000< F > &elmY, concepts::ElementMatrix< F > &em)
const Haar3d000< F > * elm
Definition: bform.hh:48
void operator()(const concepts::Element< F > &elmX, const concepts::Element< F > &elmY, concepts::ElementMatrix< F > &em)
F * valXY
Near field entries stored row wise.
Definition: bform.hh:38
void operator()(const Haar3d000< F > &elmX, const Haar3d000< F > &elmY, concepts::ElementMatrix< F > &em)
uint memory() const
Definition: bform.hh:69
virtual AglowavBF00 * clone() const
Virtual constructor.
Definition: bform.hh:53
F * valXY_
Auxiliary vectors to hold intermediate values.
Definition: bform.hh:123
Element matrix.
Definition: linearForm.hh:18
F00< F > & XFY_
Cluster approximation of the far field.
Definition: bform.hh:111
void updateChld_(const Haar3d000< F > &tst, const Haar3d000< F > &trl, F *valtsttrl, uint ldtsttrl, F *valtrltst, uint ldtrltst)
Computation of one scale entries.
void operator()(const concepts::Element< F > &elmX, const concepts::Element< F > &elmY, concepts::ElementMatrix< F > &em)
Entry in the cache to compute wavelet entries.
Definition: bform.hh:44
Cell * cache_
Cache to store wavelet entries for later usage.
Definition: bform.hh:130
void updateElm_(const Haar3d000< F > &tst, const Haar3d000< F > &trl, F *valtsttrl, uint ldtsttrl, F *valtrltst, uint ldtrltst)
Computation of one scale entries.
Identity bilinear form for wavelets.
Definition: bform.hh:143
struct Nfld * lnk
Pointer to next near field entry.
Definition: bform.hh:34
Page URL: http://wiki.math.ethz.ch/bin/view/Concepts/WebHome
21 August 2020
© 2020 Eidgenössische Technische Hochschule Zürich