dune-pdelab  2.0.0
timesteppingparameterinterface.hh
Go to the documentation of this file.
1 // -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2 // vi: set et ts=4 sw=2 sts=2:
3 #ifndef DUNE_PDELAB_TIMESTEPPINGPARAMETERINTERFACE_HH
4 #define DUNE_PDELAB_TIMESTEPPINGPARAMETERINTERFACE_HH
5 
6 #include <string>
7 
8 #include <dune/common/fvector.hh>
9 #include <dune/common/fmatrix.hh>
10 
11 namespace Dune {
12  namespace PDELab {
13 
15 
18  template<class R>
20  {
21  public:
22  typedef R RealType;
23 
26  virtual bool implicit () const = 0;
27 
30  virtual unsigned s () const = 0;
31 
35  virtual R a (int r, int i) const = 0;
36 
40  virtual R b (int r, int i) const = 0;
41 
45  virtual R d (int r) const = 0;
46 
49  virtual std::string name () const = 0;
50 
53  };
54 
55 
57 
60  template<class R>
62  {
63  public:
64 
66  {
67  D[0] = 0.0; D[1] = 1.0;
68  A[0][0] = -1.0; A[0][1] = 1.0;
69  B[0][0] = 0.0; B[0][1] = 1.0;
70  }
71 
74  virtual bool implicit () const
75  {
76  return true;
77  }
78 
81  virtual unsigned s () const
82  {
83  return 1;
84  }
85 
89  virtual R a (int r, int i) const
90  {
91  return A[r-1][i];
92  }
93 
97  virtual R b (int r, int i) const
98  {
99  return B[r-1][i];
100  }
101 
105  virtual R d (int i) const
106  {
107  return D[i];
108  }
109 
112  virtual std::string name () const
113  {
114  return std::string("implicit Euler");
115  }
116 
117  private:
118  Dune::FieldVector<R,2> D;
119  Dune::FieldMatrix<R,1,2> A;
120  Dune::FieldMatrix<R,1,2> B;
121  };
122 
124  } // namespace PDELab
125 } // namespace Dune
126 
127 #endif
virtual std::string name() const =0
Return name of the scheme.
Base parameter class for time stepping scheme parameters.
Definition: timesteppingparameterinterface.hh:19
virtual R b(int r, int i) const =0
Return entries of the B matrix.
virtual ~TimeSteppingParameterInterface()
every abstract base class has a virtual destructor
Definition: timesteppingparameterinterface.hh:52
virtual R b(int r, int i) const
Return entries of the B matrix.
Definition: timesteppingparameterinterface.hh:97
virtual unsigned s() const =0
Return number of stages of the method.
virtual R a(int r, int i) const =0
Return entries of the A matrix.
virtual bool implicit() const =0
Return true if method is implicit.
virtual R a(int r, int i) const
Return entries of the A matrix.
Definition: timesteppingparameterinterface.hh:89
virtual R d(int r) const =0
Return entries of the d Vector.
ImplicitEulerParameter()
Definition: timesteppingparameterinterface.hh:65
Parameters specifying implicit euler.
Definition: timesteppingparameterinterface.hh:61
virtual R d(int i) const
Return entries of the d Vector.
Definition: timesteppingparameterinterface.hh:105
virtual bool implicit() const
Return true if method is implicit.
Definition: timesteppingparameterinterface.hh:74
virtual std::string name() const
Return name of the scheme.
Definition: timesteppingparameterinterface.hh:112
R RealType
Definition: timesteppingparameterinterface.hh:22
virtual unsigned s() const
Return number of stages s of the method.
Definition: timesteppingparameterinterface.hh:81