bilinearForm.hh

Go to the documentation of this file.
1 
6 #ifndef hp3dbilinearform_hh
7 #define hp3dbilinearform_hh
8 
9 #include <memory>
10 #include "basics/typedefs.hh"
12 #include "formula/formula.hh"
14 #include "operator/bilinearForm.hh"
15 #include "linearFormHelper.hh"
16 
17 
18 namespace concepts {
19  // forward declarations
20  template<class F>
21  class Element;
22 
23  template<class F>
24  class ElementMatrix;
25 
26  template<class F>
27  class Array;
28 
29  class InOutParameters;
30 }
31 
32 namespace vectorial {
33  //forward declaration
34  template<class F, class G>
35  class BilinearForm;
36 }
37 
38 namespace hp3D {
39  // forward declarations
40  class Hexahedron;
41 
42  using concepts::Real;
43 
44  // *************************************************************** Laplace **
45 
50  template <class F = Real>
51  class Laplace : public concepts::BilinearForm<F> {
52  public:
55  virtual ~Laplace();
56 
57  virtual void operator()(const concepts::Element<Real>& elmX,
58  const concepts::Element<Real>& elmY,
59  concepts::ElementMatrix<F>& em) const;
60  void operator()(const Hexahedron& elmX, const Hexahedron& elmY,
61  concepts::ElementMatrix<F>& em) const;
62  virtual Laplace<F>* clone() const { return new Laplace(); }
63  protected:
64  virtual std::ostream& info(std::ostream& os) const;
65  private:
69 
70  mutable const Hexahedron* oldElm_;
71  };
72 
73  // ************************************************************** Identity **
74 
83  template <class F = Real>
84  class Identity : public concepts::BilinearForm<F> {
85  public:
89  virtual ~Identity();
90 
91  virtual void operator()(const concepts::Element<Real>& elmX,
92  const concepts::Element<Real>& elmY,
93  concepts::ElementMatrix<F>& em) const;
94  void operator()(const Hexahedron& elmX, const Hexahedron& elmY,
95  concepts::ElementMatrix<F>& em) const;
96 
127  static bool timings();
129  virtual Identity* clone() const { return new Identity(frm_); }
130  protected:
131  virtual std::ostream& info(std::ostream& os) const;
132  private:
135 
136  mutable const Hexahedron* oldElm_;
137 
140 
144  static uint timeCntr_;
145  };
146 
147  // *********************************************************** setupIdentity **
148 
153  template<class F>
156 
163  template<class F>
167  bool transposed = false);
168 
169  // ************************************************************* Advection **
170 
185  template<class F = Real>
186  class Advection : public concepts::BilinearForm<F,Real>,
187  public hp3D::LinearFormHelper_1<F>
188  {
189  public:
193  bool all = false)
194  : hp3D::LinearFormHelper_1<F>(frm1, frm2, frm3), all_(all)
195  {}
196 
199  bool all = false)
200  : hp3D::LinearFormHelper_1<F>(frm), all_(all)
201  {}
202 
203  virtual ~Advection() {}
204 
205  virtual Advection<F>* clone() const {
206  return new Advection<F>(this->frm_, all_);
207  }
208 
209  virtual void operator()(const concepts::Element<Real>& elmX,
210  const concepts::Element<Real>& elmY,
211  concepts::ElementMatrix<F>& em) const;
212 
213  protected:
214  virtual std::ostream& info(std::ostream& os) const;
215  private:
216  bool all_;
217 
219  bool operator()(const Hexahedron* elmX,
220  const Hexahedron* elmY,
221  concepts::ElementMatrix<F>& em) const;
222  };
223 
224  // ********************************************************** setupAdvection **
225 
246  template<class F>
250 
251 } // namespace hp3D
252 
253 #endif // hp3dbilinearform_hh
Basic class for a 2D or 3D map.
concepts::Array< concepts::MapReal3d > jacobianInv_
Definition: bilinearForm.hh:68
const Hexahedron * oldElm_
Intermediate data for element matrix computation.
A 3D FEM element: a hexahedron.
Definition: hexahedron.hh:37
Helper class for linearforms l(v), where v is a one form.
A function class to calculate element matrices of the bilinear form.
void operator()(const Hexahedron &elmX, const Hexahedron &elmY, concepts::ElementMatrix< F > &em) const
virtual Advection< F > * clone() const
Virtual constructor.
void operator()(const Hexahedron &elmX, const Hexahedron &elmY, concepts::ElementMatrix< F > &em) const
virtual Identity * clone() const
Intermediate data for element matrix computation.
Vector valued problems.
Definition: spaceTraits.hh:24
virtual void operator()(const concepts::Element< Real > &elmX, const concepts::Element< Real > &elmY, concepts::ElementMatrix< F > &em) const
Vector valued bilinear form.
Definition: bf_advection.hh:38
Holds parameters in hashes.
Definition: inputOutput.hh:75
GenericElement< KarniadakisMixin< F > > Element
template aliases for backwards compatibility
Definition: element.hh:270
virtual Laplace< F > * clone() const
Virtual constructor.
Definition: bilinearForm.hh:62
virtual std::ostream & info(std::ostream &os) const
Intermediate data for element matrix computation.
Identity(const concepts::ElementFormulaContainer< F > frm=concepts::ElementFormulaContainer< F >())
Constructor.
Abstract function class to evaluate a bilinear form.
Definition: bilinearForm.hh:33
const Hexahedron * oldElm_
Definition: bilinearForm.hh:70
Laplace()
Constructor.
virtual void operator()(const concepts::Element< Real > &elmX, const concepts::Element< Real > &elmY, concepts::ElementMatrix< F > &em) const
void setupIdentity(vectorial::BilinearForm< F, typename concepts::Realtype< F >::type > &bf)
Function to setup a bilinear form related to the vector Identity, namely.
An array of objects.
Definition: bilinearForm.hh:23
static void setTimings(concepts::InOutParameters *timings)
Sets the class to store the timing values in.
A function class to calculate element matrices for the Laplacian.
Definition: bilinearForm.hh:51
bool operator()(const Hexahedron *elmX, const Hexahedron *elmY, concepts::ElementMatrix< F > &em) const
Assembling for hp3D::Hexahedron.
Advection(const concepts::ElementFormulaContainer< concepts::Point< F, 3 > > frm, bool all=false)
virtual std::ostream & info(std::ostream &os) const
Returns information in an output stream.
concepts::ElementFormulaContainer< F > frm_
Element formula.
virtual ~Laplace()
Advection(const concepts::ElementFormulaContainer< F > frm1, const concepts::ElementFormulaContainer< F > frm2, const concepts::ElementFormulaContainer< F > frm3, bool all=false)
concepts::Array< F > jacobian_
Intermediate data for element matrix computation.
virtual ~Advection()
static concepts::InOutParameters * timings_
Place to store timing values.
Element matrix.
Definition: linearForm.hh:18
virtual std::ostream & info(std::ostream &os) const
Returns information in an output stream.
static uint timeCntr_
Counter for timing table.
3D hp-FEM for H1-conforming elements.
Definition: meshDX.hh:23
concepts::ElementFormulaContainer< concepts::Point< Real, 3 > > frm_
ElementFormula.
static bool timings()
Returns true if the class is able to do timings.
A function class to calculate element matrices for the mass matrix.
Definition: bilinearForm.hh:84
virtual ~Identity()
concepts::Array< F > jacobian_
Intermediate data for element matrix computation.
Definition: bilinearForm.hh:67
double Real
Type normally used for a floating point number.
Definition: typedefs.hh:36
void setupAdvection(vectorial::BilinearForm< F, typename concepts::Realtype< F >::type > &bf, const concepts::ElementFormulaContainer< concepts::Point< F, 3 > > frm)
Function to setup a bilinear form related to the vector Advection, namely.
Basic namespace for Concepts-2.
Definition: pml_formula.h:16
virtual void operator()(const concepts::Element< Real > &elmX, const concepts::Element< Real > &elmY, concepts::ElementMatrix< F > &em) const
Page URL: http://wiki.math.ethz.ch/bin/view/Concepts/WebHome
21 August 2020
© 2020 Eidgenössische Technische Hochschule Zürich