linearForm.hh

Go to the documentation of this file.
1 
6 #pragma once
7 
8 #include <memory>
9 #include "basics/typedefs.hh"
10 #include "function/linearForm.hh"
11 #include "formula.hh"
12 
13 namespace concepts {
14  // forward declarations
15  template<class F>
16  class Formula;
17 
18  class BoundaryConditions;
19 }
20 
21 namespace hp1D {
22 
23  using concepts::Real;
24 
25  // *************************************************************** Neumann **
26 
38  class Neumann : public concepts::Neumann<Real> {
39  public:
45  virtual ~Neumann();
46 
56  };
57 
58  // ***************************************************************** Riesz **
59 
66  template<class F = Real>
67  class Riesz : public concepts::LinearForm<F>, public LinearFormHelper<0,F>, public concepts::LinearFormChoice {
68  public:
73  const concepts::BoundaryConditions* bc = 0);
74  virtual ~Riesz();
75 
83  concepts::ElementMatrix<F>& em) const;
84  protected:
85  virtual std::ostream& info(std::ostream& os) const;
86  private:
89 
90  };
91 
92 
93  // ********************************************** LinearFormGradInterp_Grad **
94 
103  template<typename F = Real>
105  public LinearFormHelper<1,F> {
106  public:
113 
122  concepts::ElementMatrix<F>& em) const;
123  protected:
124  virtual std::ostream& info(std::ostream& os) const;
125  };
126 
127  // ********************************************** GradLinearForm **
128 
137  template<typename F = Real>
139  public LinearFormHelper<1,F> {
140  public:
146  virtual ~GradLinearForm();
147 
156  concepts::ElementMatrix<F>& em) const;
157  protected:
158  virtual std::ostream& info(std::ostream& os) const;
159  };
160 
161  // ********************************************** PointLinearForm **
162 
171  template<typename F = Real>
173  public:
181  : attrib_(attrib), frm_(frm) {};
190  concepts::ElementMatrix<F>& em) const;
191  protected:
192  virtual std::ostream& info(std::ostream& os) const;
193  private:
198 
199  };
200 
201 } // namespace hp1D
Helper class for Linearforms l(v), where v is i-form.
Definition: formula.hh:94
Linear form on edges in nD.
Definition: linearForm.hh:139
void operator()(const concepts::Element< Real > &elm, concepts::ElementMatrix< F > &em) const
Computes the element load vector.
Interface class for Linearform in that one can choice for the basis evaluation type,...
Definition: linearForm.hh:64
virtual ~GradLinearForm()
Class for a constant formula.
Definition: constFormula.hh:26
virtual std::ostream & info(std::ostream &os) const
Riesz(const concepts::ElementFormulaContainer< F > frm, const concepts::BoundaryConditions *bc=0)
Constructor.
const concepts::ElementFormulaContainer< F > frm_
The element formula.
Definition: linearForm.hh:197
virtual ~Riesz()
Abstract class for a linear form.
void operator()(const concepts::Element< Real > &elm, concepts::ElementMatrix< F > &em) const
Computes the element load vector.
Linear form on edges in nD.
Definition: linearForm.hh:67
Abstract class for the Neumann boundary term.
Definition: linearForm.hh:89
virtual ~Neumann()
Point evaluation in 1D.
Definition: linearForm.hh:172
LinearFormGradInterp_Grad(const concepts::ElementFormulaContainer< F > frm)
Constructor.
GradLinearForm(const concepts::ElementFormulaContainer< F > frm)
Constructor.
void operator()(const concepts::Element< Real > &elm, concepts::ElementMatrix< F > &em) const
Computes the element load vector.
void operator()(const concepts::Element< Real > &elm, concepts::ElementMatrix< Real > &em) const
Computes the element load vector.
const concepts::Set< concepts::Attribute > attrib_
Set of vertex attributes.
Definition: linearForm.hh:195
virtual std::ostream & info(std::ostream &os) const
Neumann neumann_
Reference to the linear form of the Neumann condition.
Definition: linearForm.hh:88
virtual std::ostream & info(std::ostream &os) const
PointEvaluation(const concepts::Set< concepts::Attribute > attrib, const concepts::ElementFormulaContainer< F > frm=concepts::ConstFormula< F >(1.0))
Constructor.
Definition: linearForm.hh:179
Neumann(const concepts::BoundaryConditions *bc)
Constructor.
void operator()(const concepts::Element< Real > &elm, concepts::ElementMatrix< F > &em) const
Computes the element load vector.
Linear form on edges in nD.
Definition: linearForm.hh:105
virtual std::ostream & info(std::ostream &os) const
Linear form on edges in nD.
Definition: linearForm.hh:38
double Real
Type normally used for a floating point number.
Definition: typedefs.hh:36
Basic namespace for Concepts-2.
Definition: pml_formula.h:16
1D hp-FEM
Page URL: http://wiki.math.ethz.ch/bin/view/Concepts/WebHome
21 August 2020
© 2020 Eidgenössische Technische Hochschule Zürich