bilinearFormHelper.hh

Go to the documentation of this file.
1 
6 #ifndef hp2dbilinearformhelper_hh
7 #define hp2dbilinearformhelper_hh
8 
9 #include "linearFormHelper.hh"
11 
12 namespace concepts {
13 
14  // forward declarations
15  class SharedJacobianDet;
16 
17  template<uint dim>
18  class SharedJacobianAdj;
19 
20  template<class F>
21  class Array;
22 
23  template<class F, class G>
24  class ElementFormula;
25 
26  template<class F, uint dim>
27  class Point;
28 }
29 
30 namespace hp2D {
31 
32  // forward declarations
33  template<class F>
34  class BaseQuad;
35 
36  template<class F>
37  class ArrayElementFormula;
38 
39  using concepts::Real;
40 
41  // *********************************************** BilinearFormHelper_0_0 **
42 
49  template<class F>
51  public:
55  : LinearFormHelper_0<F>(frm) {}
56 
58  return this->frm_; }
59  };
60 
61  // ******************************************** BilinearFormHelper_0_1_Part **
62 
77  template<class F>
79  public:
81 
83 
88  protected:
96  void computeIntermediate_(const BaseQuad<concepts::Real>& elm, const int i) const;
104 
107  };
108 
109  // ************************************************ BilinearFormHelper_1_1 **
110 
139  template<class F, class G = typename concepts::Realtype<F>::type>
141  public:
147 
150 
155  protected:
170  const int i = -1, const int j = -1) const;
196  private:
202  concepts::Array<Real>& detJ_inv) const;
208  (concepts::Array<concepts::Mapping<Real,2> >& intermediateMatrix,
209  const int i, const int j) const;
210 
213  };
214 
215  // ************************************************ BilinearFormHelper_2_1 **
216 
228  template<class F>
230  public:
233 
236 
238  protected:
240  void computeIntermediate_(const BaseQuad<Real>& elm) const;
246 
249  };
250 
251  // ************************************************ BilinearFormHelper_2_2 **
252 
261  template<class F>
263  public:
267  protected:
269  void computeIntermediate_(const BaseQuad<Real>& elm) const;
277  };
278 
279 } // namespace hp2D
280 
281 #endif // hp2dbilinearformhelper_hh
Basic class for a 2D or 3D map.
BilinearFormHelper_0_0(const concepts::ElementFormulaContainer< F > frm=concepts::ElementFormulaContainer< F >())
BilinearFormHelper_2_2(const concepts::ElementFormulaContainer< F > frm=concepts::ElementFormulaContainer< F >())
concepts::ElementFormulaContainer< concepts::Mapping< G, 2 > > frmM_
Matrix element formula.
concepts::Array< F > intermediateValue_
Intermediate value.
void computeJacobianMatrix_(const BaseQuad< Real > &elm, concepts::Array< concepts::Mapping< Real, 2 > > &J, concepts::Array< Real > &detJ_inv) const
Compute the Jacobian matrix and the inverse of its determinant on each quadrature point.
concepts::RCP< concepts::SharedJacobianAdj< 2 > > data() const
Gets the pointer to the shared data.
void data(const concepts::RCP< concepts::SharedJacobianAdj< 2 > > d)
Set the pointer to the shared data.
concepts::ElementFormulaContainer< F > frm_
Element formula.
concepts::RCP< concepts::SharedJacobianAdj< 2 > > sharedData_
Shared data for vectorial bilinear forms.
concepts::ElementFormulaContainer< F > frm_
ElementFormula.
concepts::ElementFormulaContainer< F > frm_
ElementFormula.
concepts::Array< concepts::Mapping< G, 2 > > intermediateMatrix_
Intermediate matrix.
const concepts::ElementFormulaContainer< F > formula() const
Helper class for bilinearforms a(u,v), where u and v are 2-forms.
Helper class for linearforms l(v), where v is a 0-form.
virtual ~BilinearFormHelper_1_1()
Destructor.
concepts::RCP< concepts::SharedJacobianAdj< 2 > > data() const
Gets the pointer to the shared data.
Array of formula values on quadrature points.
BilinearFormHelper_0_1_Part(const concepts::ElementFormulaContainer< F > frm)
void computeIntermediate_(const BaseQuad< Real > &elm, const int i=-1, const int j=-1) const
Compute the intermediate data for element matrix computation.
BilinearFormHelper_2_1(const concepts::ElementFormulaContainer< F > frm1, const concepts::ElementFormulaContainer< F > frm2)
concepts::Array< F > intermediateValue_
Intermediate value.
Helper class for bilinearforms a(u,v), where u and v are 1-forms, which computes intermediate data fo...
ArrayElementFormula< concepts::Point< F, 2 > > intermediateVector_
Intermediate vector (on each quadrature point)
void computeIntermediate_(const BaseQuad< Real > &elm) const
Compute the intermediate data for element matrix computation.
concepts::ElementFormulaContainer< concepts::Point< F, 2 > > frm_
ElementFormula.
2D hp-FEM for H1-conforming elements.
Helper class for bilinearforms , where is a 2-form and a 1-form.
BilinearFormHelper_1_1(const concepts::ElementFormulaContainer< concepts::Mapping< G, 2 > > frmM)
Constructor for matrix valued function.
An array of objects.
Definition: bilinearForm.hh:23
concepts::RCP< concepts::SharedJacobianAdj< 2 > > sharedData_
Shared data for vectorial bilinear forms.
Helper class for bilinearforms a(u,v), where u and v are 0-forms.
BilinearFormHelper_1_1(const concepts::ElementFormulaContainer< F > frm)
Constructor for scalar function.
void computeIntermediate_(const BaseQuad< concepts::Real > &elm, const int i) const
Compute the intermediate data for element matrix computation.
void data(const concepts::RCP< concepts::SharedJacobianAdj< 2 > > d)
Set the pointer to the shared data.
Reference-counting pointer.
Definition: bf_iddiv.hh:15
A base of a 2D quad FEM element for different basis functions.
void computeIntermediate_(const BaseQuad< Real > &elm) const
Compute the intermediate data for element matrix computation.
void computeadjJ_adjJT_rank1_(concepts::Array< concepts::Mapping< Real, 2 > > &intermediateMatrix, const int i, const int j) const
Computes the either adj(J)*adj(J)^T or in the case of partial derivatives (i > NONE,...
concepts::ElementFormulaContainer< F > frm_
ElementFormula.
Helper class for bilinear forms a(u,v) where u is a 0-form and v is a 1-form.
ArrayElementFormula< concepts::Point< F, 2 > > intermediateVector_
Intermediate vector (on each quadrature point)
BilinearFormHelper_2_1(const concepts::ElementFormulaContainer< concepts::Point< F, 2 > > frm)
double Real
Type normally used for a floating point number.
Definition: typedefs.hh:36
Shared data for bilinear forms on vectorial spaces, like hp2D::RotRot and hp2D::DivDiv.
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