1 #ifndef DUNE_FEM_SOLVER_RUNGEKUTTA_BUTCHERTABLE_HH 2 #define DUNE_FEM_SOLVER_RUNGEKUTTA_BUTCHERTABLE_HH 5 #include <dune/common/dynmatrix.hh> 6 #include <dune/common/dynvector.hh> 14 template<
class Field >
35 static Dune::DynamicMatrix< FieldType >
makeMatrix (
int m,
int n,
const FieldType *data )
37 Dune::DynamicMatrix< FieldType >
A( m, n );
38 for(
int i = 0; i < m; ++i )
39 std::copy( data + i*n, data + (i+1)*n, A[ i ].begin() );
43 static Dune::DynamicVector< FieldType >
makeVector (
int n,
const FieldType *data )
45 Dune::DynamicVector< FieldType > v( n );
46 std::copy( data, data + n, v.begin() );
80 template<
class Field >
90 : Base(stages,order,a,b,c)
94 Dune::DynamicMatrix< FieldType >
B ()
const {
return Base::makeMatrix(
stages_,
stages_, a2_ ); }
104 #endif // #ifndef DUNE_FEM_SOLVER_RUNGEKUTTA_BUTCHERTABLE_HH static Dune::DynamicMatrix< FieldType > makeMatrix(int m, int n, const FieldType *data)
Definition: butchertable.hh:35
Definition: butchertable.hh:15
int order() const
Definition: butchertable.hh:31
static Dune::DynamicVector< FieldType > makeVector(int n, const FieldType *data)
Definition: butchertable.hh:43
SimpleButcherTable< double > implicitEulerButcherTable()
Definition: butchertable.cc:58
SimpleButcherTable< double > semiImplicitEulerButcherTable(bool expl)
Definition: butchertable.cc:93
Definition: multistep.hh:16
SimpleButcherTable< double > semiImplicitSSP222ButcherTable(bool expl)
Definition: butchertable.cc:189
SimpleButcherTable< double > gauss2ButcherTable()
Definition: butchertable.cc:73
ROWSimpleButcherTable< double > row3ButcherTable()
Definition: butchertable.cc:362
const FieldType * c_
Definition: butchertable.hh:51
ROWSimpleButcherTable(int stages, int order, const FieldType *a, const FieldType *b, const FieldType *c, const FieldType *a2)
Definition: butchertable.hh:89
SimpleButcherTable< double > semiImplicitARK46ButcherTable(bool expl)
SimpleButcherTable< double > semiImplicit33ButcherTable(bool expl)
Definition: butchertable.cc:159
Dune::DynamicMatrix< FieldType > B() const
Definition: butchertable.hh:94
Dune::DynamicMatrix< FieldType > A() const
Definition: butchertable.hh:27
SimpleButcherTable< double > implicit3ButcherTable()
Definition: butchertable.cc:44
Field FieldType
Definition: butchertable.hh:87
int stages() const
Definition: butchertable.hh:32
Definition: butchertable.hh:81
SimpleButcherTable< double > semiImplicitARK34ButcherTable(bool expl)
Field FieldType
Definition: butchertable.hh:20
SimpleButcherTable< double > implicit34ButcherTable()
Definition: butchertable.cc:24
SimpleButcherTable(int stages, int order, const FieldType *a, const FieldType *b, const FieldType *c)
Definition: butchertable.hh:22
const FieldType * b_
Definition: butchertable.hh:51
int order_
Definition: butchertable.hh:50
const FieldType * a_
Definition: butchertable.hh:51
Dune::DynamicVector< FieldType > c() const
Definition: butchertable.hh:29
int stages_
Definition: butchertable.hh:50
Dune::DynamicVector< FieldType > b() const
Definition: butchertable.hh:28
SimpleButcherTable< double > semiImplicit23ButcherTable(bool expl)
Definition: butchertable.cc:124
SimpleButcherTable< double > semiImplicitIERK45ButcherTable(bool expl)
Definition: butchertable.cc:326
ROWSimpleButcherTable< double > row2ButcherTable()
Definition: butchertable.cc:341