bilinearForm2D.hh

Go to the documentation of this file.
1 
6 #ifndef bilinearLinFEM_hh
7 #define bilinearLinFEM_hh
8 
9 #include "basics/typedefs.hh"
11 #include "operator/bilinearForm.hh"
12 #include "toolbox/array.hh"
13 
14 namespace linearFEM {
15 
16  using concepts::Real;
17 
18  class Quad;
19  class Triangle;
20 
21  // ************************************************************* Laplace2d **
22 
30  class Laplace2d : public concepts::BilinearForm<Real> {
31  public:
32  virtual void operator()(const concepts::Element<Real>& elmX,
33  const concepts::Element<Real>& elmY,
39  void operator()(const linearFEM::Quad& elmX,
40  const linearFEM::Quad& elmY,
45  void operator()(const linearFEM::Triangle& elmX,
46  const linearFEM::Triangle& elmY,
48  virtual Laplace2d* clone() const { return new Laplace2d(); }
49  private:
52  };
53 
54  // ************************************************************ Identity2d **
55 
63  class Identity2d : public concepts::BilinearForm<Real> {
64  public:
65  virtual void operator()(const concepts::Element<Real>& elmX,
66  const concepts::Element<Real>& elmY,
72  void operator()(const linearFEM::Quad& elmX,
73  const linearFEM::Quad& elmY,
78  void operator()(const linearFEM::Triangle& elmX,
79  const linearFEM::Triangle& elmY,
81  virtual Identity2d* clone() const { return new Identity2d(); }
82  private:
85  };
86 
87 } // namespace linearFEM
88 
89 #endif // bilinearLinFEM_hh
void operator()(const linearFEM::Triangle &elmX, const linearFEM::Triangle &elmY, concepts::ElementMatrix< Real > &em) const
Computes the element stiffness matrix for a triangle.
virtual void operator()(const concepts::Element< Real > &elmX, const concepts::Element< Real > &elmY, concepts::ElementMatrix< Real > &em) const
virtual void operator()(const concepts::Element< Real > &elmX, const concepts::Element< Real > &elmY, concepts::ElementMatrix< Real > &em) const
virtual Laplace2d * clone() const
Virtual constructor.
Discrete equivalent of the Laplacian in 2D for linear FEM.
Abstract function class to evaluate a bilinear form.
Definition: bilinearForm.hh:33
An array of objects.
Definition: bilinearForm.hh:23
Linear FEM in 1D, 2D and 3D.
Definition: spaceTraits.hh:19
concepts::Array< Real > shpfct_
Storage area for values of shape functions.
void operator()(const linearFEM::Triangle &elmX, const linearFEM::Triangle &elmY, concepts::ElementMatrix< Real > &em) const
Computes the element mass matrix for a triangle.
virtual Identity2d * clone() const
Virtual constructor.
Quadrilateral element with bilinear shape functions in 2D.
Definition: element2D.hh:79
Discrete equivalent for a reaction term in 2D for linear FEM.
Triangular element with linear shape functions in 2D.
Definition: element2D.hh:29
concepts::Array< concepts::Real2d > shpfct_
Storage area for values of derivatives of shape functions.
double Real
Type normally used for a floating point number.
Definition: typedefs.hh:36
void operator()(const linearFEM::Quad &elmX, const linearFEM::Quad &elmY, concepts::ElementMatrix< Real > &em) const
Computes the element stiffness matrix for a quadrilateral.
void operator()(const linearFEM::Quad &elmX, const linearFEM::Quad &elmY, concepts::ElementMatrix< Real > &em) const
Computes the element mass matrix for a quadrilateral.
Page URL: http://wiki.math.ethz.ch/bin/view/Concepts/WebHome
21 August 2020
© 2020 Eidgenössische Technische Hochschule Zürich