3 #ifndef DUNE_PDELAB_TIMESTEPPINGPARAMETERINTERFACE_HH
4 #define DUNE_PDELAB_TIMESTEPPINGPARAMETERINTERFACE_HH
8 #include <dune/common/fvector.hh>
9 #include <dune/common/fmatrix.hh>
30 virtual unsigned s ()
const = 0;
35 virtual R
a (
int r,
int i)
const = 0;
40 virtual R
b (
int r,
int i)
const = 0;
45 virtual R
d (
int r)
const = 0;
49 virtual std::string
name ()
const = 0;
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;
81 virtual unsigned s ()
const
89 virtual R
a (
int r,
int i)
const
97 virtual R
b (
int r,
int i)
const
105 virtual R
d (
int i)
const
112 virtual std::string
name ()
const
114 return std::string(
"implicit Euler");
118 Dune::FieldVector<R,2> D;
119 Dune::FieldMatrix<R,1,2> A;
120 Dune::FieldMatrix<R,1,2> B;
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