bilinearForm.hh

Go to the documentation of this file.
1 
6 #ifndef BilinearForm_hh
7 #define BilinearForm_hh
8 
9 #include <memory>
10 #include "basics/outputOperator.hh"
11 #include "basics/cloneable.hh"
12 #include "space/elementPairs.hh"
14 #include "basics/exceptions.hh"
15 
16 namespace concepts {
17 
18  // forward declaration
19  template<class F>
20  class Element;
21 
22  template<class F>
23  class ElementMatrix;
24 
25  // ********************************************************** BilinearForm **
26 
32  template<class F, class G = typename Realtype<F>::type>
33  class BilinearForm : public Cloneable, public virtual OutputOperator {
34  public:
43  virtual void operator()(const Element<G>& elmX, const Element<G>& elmY,
44  ElementMatrix<F>& em) const = 0;
57  virtual void operator()(const Element<G>& elmX, const Element<G>& elmY,
58  ElementMatrix<F>& em, const ElementPair<G>& ep) const{
59  operator()(elmX, elmY, em);
60  }
61 
62  virtual BilinearForm* clone() const = 0;
63 
64  protected:
65  virtual std::ostream& info(std::ostream& os) const;
66  };
67 
68 
69 
70 
71 
72 
73 
74  // ****************************************************** BilinearFormLiCo **
75 
82  template<class F, class G = typename Realtype<F>::type>
83  class BilinearFormLiCo : public BilinearForm<F,G> {
84  public:
92  BilinearForm<F,G>& bfB,
93  const F cA = 1.0, const F cB = 1.0) :
94 
95  bfAptr_(nullptr), bfBptr_(nullptr),
96  bfA_(&bfA), bfB_(&bfB),
97 
98  cA_(cA), cB_(cB), em_() {}
101  virtual ~BilinearFormLiCo() {}
102  virtual void operator()(const Element<G>& elmX,
103  const Element<G>& elmY,
104  ElementMatrix<F>& em) const;
105  virtual BilinearFormLiCo<F,G>* clone() const;
106  protected:
107  virtual std::ostream& info(std::ostream& os) const;
108  private:
109  std::unique_ptr<BilinearForm<F,G> > bfAptr_, bfBptr_;
112 
113 
115 
117  const F cA_;
119  const F cB_;
122  };
123 
124 
125 
126 
127 
128 
129 
130 
131 
132  // ************************************************* BilinearFormContainer **
133  // @autor Robert Gruhlke, 2016
134 
135  template<class F, typename G = typename Realtype<F>::type>
136  class BilinearFormContainer : public BilinearForm<F,G> {
137  public:
141 
145 
147 
149 
152 
153 
154  virtual void operator()(const Element<G>& elmX,
155  const Element<G>& elmY,
156  ElementMatrix<F>& em) const;
157 
160  std::cout << "in clone()" << std::endl;
161  return new BilinearFormContainer<F,G>(*this);
162  }
163 
165  bool empty() const { return bf_.get() == 0; }
166 
167 
169  protected:
170  virtual std::ostream& info(std::ostream& os) const;
171  private:
172 
175 
176 
177  };
178 
179 
180  // ******************************************* PointerToEmptyBilinearForm **
181 
187  public:
190  {
191  errorMessage_ = "Empty RCP<BilinearForm>";
192  }
193  virtual ~PointerToEmptyBilinearForm() throw() {}
194  protected:
195  virtual std::ostream& info(std::ostream& os) const throw();
196  private:
197  std::string errorMessage_;
198  };
199 
200 
201 
202  // ********************************************************* BilinearF_Sum **
209 template<class F, class H = F, class J = F, class G = typename Realtype<F>::type>
210  class BilinearF_Sum : public BilinearForm<F,G> {
211  public:
213  const BilinearFormContainer<J,G> bf2)
214  : bf1_(bf1), bf2_(bf2)
215  {
216  if (bf1_.empty()) throw(PointerToEmptyBilinearForm());
217  if (bf2_.empty()) throw(PointerToEmptyBilinearForm());
218  }
219 
220  virtual ~BilinearF_Sum(){};
221 
222  void operator()(const Element<G>& elmX,
223  const Element<G>& elmY,
224  ElementMatrix<F>& em) const
225  {
226  //intermediate Matrix
227  ElementMatrix<H> em_H(elmX.T().m(), elmY.T().m());
228  ElementMatrix<J> em_J(elmX.T().m(), elmY.T().m());
229  em.resize(elmX.T().m(), elmY.T().m());
230  em_H.zeros();
231  em_J.zeros();
232  em.zeros();
233 
234  bf1_(elmX, elmY, em_H);
235  bf2_(elmX, elmY, em_J);
236  F* emData = (F*)em;
237 
238  H* em_H_Data = (H*)em_H;
239  J* em_J_Data = (J*)em_J;
240 
241  for (uint i = 0; i < elmX.T().m()*elmY.T().m(); ++i, emData++, em_H_Data++, em_J_Data++) {
242  *emData = *em_H_Data + *em_J_Data;
243  }
244  }
245 
246 
247  virtual BilinearF_Sum<F,H,J,G>* clone() const {
248  return new BilinearF_Sum<F,H,J,G>(bf1_, bf2_);
249  }
250 
251  protected:
252  virtual std::ostream& info(std::ostream& os) const {
253  return os << typeOf(*this) << "(" << bf1_
254  << " , " << bf2_ << ")";
255  }
256  private:
261  };
262 
263 
264 
265 
266 // ********************************************************* BilinearF_W **
273  template<class F, class H = F, class J = F, class G = typename Realtype<F>::type>
274  class BilinearF_W : public BilinearForm<F,G> {
275  public:
277  const J w)
278  : bf1_(bf1), w_(w)
279  {
280  if (bf1_.empty()) throw(PointerToEmptyBilinearForm());
281  }
282 
283  virtual ~BilinearF_W(){};
284 
285  void operator()(const Element<G>& elmX,
286  const Element<G>& elmY,
287  ElementMatrix<F>& em) const
288  {
289  concepts::ElementMatrix<H> tmp(elmX.T().m(), elmY.T().m());
290  tmp.zeros();
291 
292  bf1_(elmX, elmY, tmp);
293 
294  F* emData = (F*)em;
295  H* tmpData = (H*)tmp;
296 
297  for (uint i = 0; i < elmX.T().m()*elmY.T().m(); ++i, emData++, tmpData++) {
298  *emData = *tmpData * w_;
299  }
300  }
301 
302 
303  virtual BilinearF_W<F,H,J,G>* clone() const {
304  return new BilinearF_W<F,H,J,G>(bf1_, w_);
305  }
306 
307  protected:
308  virtual std::ostream& info(std::ostream& os) const {
309  return os << typeOf(*this) << "(" << bf1_
310  << " , " << w_ << ")";
311  }
312  private:
316  J w_;
317  };
318 
319 
320 
321 
322 
323 
328  const BilinearFormContainer<Real> frm2);
329 
332  const BilinearFormContainer<Real> frm2);
333 
336  const BilinearFormContainer<Cmplx> frm2);
337 
340  const BilinearFormContainer<Cmplx> frm2);
341 
342 
343 
348  const Real w);
349 
352  const Cmplx w);
353 
356  const Real w);
357 
360  const Cmplx w);
361 
362  //Multiplikation from left
364  operator*(const Real w,
365  const BilinearFormContainer<Real> frm1);
366 
368  operator*(const Cmplx w,
369  const BilinearFormContainer<Real> frm1);
370 
372  operator*(const Real w,
373  const BilinearFormContainer<Cmplx> frm1);
374 
376  operator*(const Cmplx w,
377  const BilinearFormContainer<Cmplx> frm1);
378 
379 
380 
385  const BilinearFormContainer<Real> frm2);
386 
389  const BilinearFormContainer<Real> frm2);
390 
393  const BilinearFormContainer<Cmplx> frm2);
394 
397  const BilinearFormContainer<Cmplx> frm2);
398 
399 
400 
401 
402 } // namespace concepts
403 
404 #endif // BilinearForm_hh
Product of scalar and a bilinear form with possible different field type F.
virtual BilinearFormContainer< F, G > * clone() const
Virtual copy constructor.
BilinearFormContainer(const BilinearForm< F, G > &bf)
Constructor for an element formula, takes a clone.
virtual std::ostream & info(std::ostream &os) const
Returns information in an output stream.
A linear combination of bilinear forms.
Definition: bilinearForm.hh:83
Gives access to a pair of elements.
Definition: elementPairs.hh:25
J w_
the scalar value
BilinearFormContainer< J, G > bf2_
Second element formula.
const F cA_
Coefficient of bfA.
BilinearFormContainer(const BilinearFormContainer< F, G > &frm)
Copy constructor.
Base class for exceptions.
Definition: exceptions.hh:86
Exception class to express that the RCP pointer points to 0.
bool empty() const
Returns true if no formula is stored.
ElementFormulaContainer< Real > operator-(const ElementFormulaContainer< Real > frm, const Real a)
Simple subtracting of a element formulas and a constant via –operator.
GenericElement< KarniadakisMixin< F > > Element
template aliases for backwards compatibility
Definition: element.hh:270
Point< typename Combtype< F, Real >::type, dim > operator*(const Real x, const Point< F, dim > &y)
virtual void operator()(const Element< G > &elmX, const Element< G > &elmY, ElementMatrix< F > &em) const
virtual const TMatrixBase< F > & T() const =0
Returns the T matrix of the element.
virtual std::ostream & info(std::ostream &os) const
Returns information in an output stream.
BilinearFormContainer< H, G > bf1_
First element formula.
concepts::RCP< BilinearForm< F, G > > bA
virtual std::ostream & info(std::ostream &os) const
Returns information in an output stream.
Abstract function class to evaluate a bilinear form.
Definition: bilinearForm.hh:33
std::unique_ptr< BilinearForm< F, G > > bfAptr_
virtual BilinearForm * clone() const =0
Virtual constructor.
BilinearF_Sum(const BilinearFormContainer< H, G > bf1, const BilinearFormContainer< J, G > bf2)
Cloneable interface.
Definition: cloneable.hh:16
virtual std::ostream & info(std::ostream &os) const
Returns information in an output stream.
BilinearFormContainer()
Constructor with no argument.
virtual BilinearFormLiCo< F, G > * clone() const
Virtual constructor.
BilinearForm< F, G > * bfB_
void operator()(const Element< G > &elmX, const Element< G > &elmY, ElementMatrix< F > &em) const
const F cB_
Coefficient of bfB.
BilinearFormLiCo(const BilinearFormLiCo &b)
Copy constructor. This copy constructor implements a deep copy.
std::complex< Real > Cmplx
Type for a complex number. It also depends on the setting of Real.
Definition: typedefs.hh:39
BilinearF_W(const BilinearFormContainer< H, G > bf1, const J w)
virtual void operator()(const Element< G > &elmX, const Element< G > &elmY, ElementMatrix< F > &em) const
virtual BilinearF_Sum< F, H, J, G > * clone() const
Virtual constructor.
concepts::RCP< BilinearForm< F, G > > bB
Reference-counting pointer.
Definition: bf_iddiv.hh:15
BilinearForm< F, G > * bfA_
Bilinear form.
virtual std::ostream & info(std::ostream &os) const
Returns information in an output stream.
virtual void resize(uint m, uint n)
Sets a new size.
Element matrix.
Definition: linearForm.hh:18
virtual void operator()(const Element< G > &elmX, const Element< G > &elmY, ElementMatrix< F > &em) const =0
Evaluates the bilinear form for all shape functions on elmX and elmY and stores the result in the mat...
virtual std::ostream & info(std::ostream &os) const
Returns information in an output stream.
BilinearFormContainer< H, G > bf1_
First element formula.
BilinearFormContainer(const BilinearForm< F, G > *frm)
void operator()(const Element< G > &elmX, const Element< G > &elmY, ElementMatrix< F > &em) const
BilinearFormContainer(const RCP< const BilinearForm< F, G > > frm)
void zeros()
Fills the matrix with zeros.
Definition: element.hh:295
Frm_Sum< Real > operator+(const Formula< Real > &frm, const Real a)
Simple adding two formulas by +-operator.
std::unique_ptr< BilinearForm< F, G > > bfBptr_
ElementMatrix< F > em_
Local element matrix used as temporary storage.
std::string typeOf(const T &t)
Return the typeid name of a class object.
Definition: output.hh:43
RCP< const BilinearForm< F, G > > bf_
bilinear form is stored as a RCP
virtual BilinearF_W< F, H, J, G > * clone() const
Virtual constructor.
Class providing an output operator.
virtual void operator()(const Element< G > &elmX, const Element< G > &elmY, ElementMatrix< F > &em, const ElementPair< G > &ep) const
Evaluates the bilinear form for all shape functions on elmX and elmY and stores the result in the mat...
Definition: bilinearForm.hh:57
Sum of two bilinear forms of possible different underlying field type F.
BilinearFormLiCo(BilinearForm< F, G > &bfA, BilinearForm< F, G > &bfB, const F cA=1.0, const F cB=1.0)
Constructor.
Definition: bilinearForm.hh:91
double Real
Type normally used for a floating point number.
Definition: typedefs.hh:36
Basic namespace for Concepts-2.
Definition: pml_formula.h:16
Page URL: http://wiki.math.ethz.ch/bin/view/Concepts/WebHome
21 August 2020
© 2020 Eidgenössische Technische Hochschule Zürich