dune-fem  2.4.1-rc
explicit.hh
Go to the documentation of this file.
1 #ifndef DUNE_FEM_SOLVER_RUNGEKUTTA_EXPLICIT_HH
2 #define DUNE_FEM_SOLVER_RUNGEKUTTA_EXPLICIT_HH
3 
4 //- system includes
5 #include <iostream>
6 #include <cmath>
7 #include <vector>
8 #include <cassert>
9 
10 //- Dune includes
11 #include <dune/fem/io/parameter.hh>
15 
16 namespace DuneODE
17 {
18 
19  using namespace Dune;
20  using namespace Fem;
21  using namespace std;
22 
61  template<class DestinationImp>
63  : public OdeSolverInterface<DestinationImp>
64  {
65  public:
66  typedef DestinationImp DestinationType;
68  typedef typename DestinationType :: DiscreteFunctionSpaceType SpaceType;
69 
71 
72  using OdeSolverInterface<DestinationImp> :: solve ;
73  protected:
74  std::vector< std::vector<double> > a;
75  std::vector<double> b;
76  std::vector<double> c;
77  std::vector<DestinationType*> Upd;
78  const int ord_;
79 
80  public:
87  ExplicitRungeKuttaSolver(OperatorType& op,
88  TimeProviderBase& tp,
89  const int pord,
90  bool verbose ) :
91  a(0),b(0),c(0), Upd(0),
92  ord_(pord),
93  op_(op),
94  tp_(tp),
95  initialized_(false)
96  {
97  assert(ord_>0);
98  a.resize(ord_);
99  for (int i=0; i<ord_; ++i)
100  {
101  a[i].resize(ord_);
102  }
103  b.resize(ord_);
104  c.resize(ord_);
105 
106  switch (ord_)
107  {
108  case 4 :
109  a[0][0]=0.; a[0][1]=0.; a[0][2]=0.; a[0][3]=0.;
110  a[1][0]=1.0; a[1][1]=0.; a[1][2]=0.; a[1][3]=0.;
111  a[2][0]=0.25; a[2][1]=0.25; a[2][2]=0.; a[2][3]=0.;
112  a[3][0]=1./6.; a[3][1]=1./6.; a[3][2]=2./3.; a[3][3]=0.;
113  b[0]=1./6.; b[1]=1./6.; b[2]=2./3.; b[3]=0.;
114  c[0]=0.; c[1]=1.0; c[2]=0.5; c[3]=1.0;
115  break;
116  case 3 :
117  a[0][0]=0.; a[0][1]=0.; a[0][2]=0.;
118  a[1][0]=1.0; a[1][1]=0.; a[1][2]=0.;
119  a[2][0]=0.25; a[2][1]=0.25; a[2][2]=0.;
120  b[0]=1./6.; b[1]=1./6.; b[2]=2./3.;
121  c[0]=0.; c[1]=1; c[2]=0.5;
122  break;
123  case 2 :
124  a[0][0]=0.; a[0][1]=0.;
125  a[1][0]=1.0; a[1][1]=0.;
126  b[0]=0.5; b[1]=0.5;
127  c[0]=0; c[1]=1;
128  break;
129  case 1:
130  a[0][0]=0.;
131  b[0]=1.;
132  c[0]=0.;
133  break;
134  default : std::cerr << "Runge-Kutta method of this order not implemented"
135  << std::endl;
136  abort();
137  }
138 
139  // create update memory
140  for (int i=0; i<ord_; ++i)
141  {
142  Upd.push_back(new DestinationType("URK",op_.space()) );
143  }
144  Upd.push_back(new DestinationType("Ustep",op_.space()) );
145  }
146 
147 
148  ExplicitRungeKuttaSolver(OperatorType& op,
149  TimeProviderBase& tp,
150  const int pord,
152  : ExplicitRungeKuttaSolver( op, tp, pord, true )
153  {}
154 
157  {
158  for(size_t i=0; i<Upd.size(); ++i)
159  delete Upd[i];
160  }
161 
163  void initialize(const DestinationType& U0)
164  {
165  if( ! initialized_ )
166  {
167  // Compute Steps
168  op_(U0, *(Upd[0]));
169  initialized_ = true;
170 
171  // provide operators time step estimate
172  tp_.provideTimeStepEstimate( op_.timeStepEstimate() );
173  }
174  }
175 
177  void solve(DestinationType& U0, MonitorType& monitor )
178  {
179  // no information here
180  monitor.reset();
181 
182  // initialize
183  if( ! initialized_ )
184  {
185  DUNE_THROW(InvalidStateException,"ExplicitRungeKuttaSolver wasn't initialized before first call!");
186  }
187 
188  // get cfl * timeStepEstimate
189  const double dt = tp_.deltaT();
190  // get time
191  const double t = tp_.time();
192 
193  // set new time
194  op_.setTime( t );
195 
196  // Compute Steps
197  op_(U0, *(Upd[0]));
198 
199  // provide operators time step estimate
200  tp_.provideTimeStepEstimate( op_.timeStepEstimate() );
201 
202  for (int i=1; i<ord_; ++i)
203  {
204  (Upd[ord_])->assign(U0);
205  for (int j=0; j<i ; ++j)
206  {
207  (Upd[ord_])->axpy((a[i][j]*dt), *(Upd[j]));
208  }
209 
210  // set new time
211  op_.setTime( t + c[i]*dt );
212 
213  // apply operator
214  op_( *(Upd[ord_]), *(Upd[i]) );
215 
216  // provide operators time step estimate
217  tp_.provideTimeStepEstimate( op_.timeStepEstimate() );
218  }
219 
220  // Perform Update
221  for (int j=0; j<ord_; ++j)
222  {
223  U0.axpy((b[j]*dt), *(Upd[j]));
224  }
225  }
226 
227  void description(std::ostream& out) const
228  {
229  out << "ExplRungeKutta, steps: " << ord_
230  //<< ", cfl: " << this->tp_.factor()
231  << "\\\\" <<std::endl;
232  }
233 
234  protected:
235  // operator to solve for
236  OperatorType& op_;
237  // time provider
239  // init flag
241  };
242 
245 } // namespace DuneODE
246 
247 #endif // #ifndef DUNE_FEM_SOLVER_RUNGEKUTTA_EXPLICIT_HH
std::vector< std::vector< double > > a
Definition: explicit.hh:74
TimeProviderBase & tp_
Definition: explicit.hh:238
void description(std::ostream &out) const
print description of ODE solver to out stream
Definition: explicit.hh:227
Interface class for ODE Solver.
Definition: odesolverinterface.hh:11
DestinationImp DestinationType
Definition: explicit.hh:66
void initialize(const DestinationType &U0)
apply operator once to get dt estimate
Definition: explicit.hh:163
~ExplicitRungeKuttaSolver()
destructor
Definition: explicit.hh:156
void solve(DestinationType &U0, MonitorType &monitor)
solve the system
Definition: explicit.hh:177
Definition: multistep.hh:16
void reset()
Definition: odesolverinterface.hh:34
void axpy(const T &a, const T &x, T &y)
Definition: space/basisfunctionset/functor.hh:37
DestinationType::DiscreteFunctionSpaceType SpaceType
Definition: explicit.hh:68
general base for time providers
Definition: timeprovider.hh:35
const int ord_
Definition: explicit.hh:78
Definition: coordinate.hh:4
ExplicitRungeKuttaSolver(OperatorType &op, TimeProviderBase &tp, const int pord, const Dune::Fem::ParameterReader &parameter=Dune::Fem::Parameter::container())
Definition: explicit.hh:148
OdeSolverInterface< DestinationImp >::MonitorType MonitorType
Definition: explicit.hh:70
static ParameterContainer & container()
Definition: io/parameter.hh:190
SpaceOperatorInterface< DestinationImp > OperatorType
Definition: explicit.hh:67
STL namespace.
bool initialized_
Definition: explicit.hh:240
std::vector< double > b
Definition: explicit.hh:75
std::vector< DestinationType * > Upd
Definition: explicit.hh:77
ExplicitRungeKuttaSolver(OperatorType &op, TimeProviderBase &tp, const int pord, bool verbose)
constructor
Definition: explicit.hh:87
Definition: odesolverinterface.hh:17
Exlicit RungeKutta ODE solver.
Definition: explicit.hh:62
std::vector< double > c
Definition: explicit.hh:76
OperatorType & op_
Definition: explicit.hh:236