AdaptiveTimeStepping_impl.hpp
Go to the documentation of this file.
1 /*
2  Copyright 2014 IRIS AS
3 
4  This file is part of the Open Porous Media project (OPM).
5 
6  OPM is free software: you can redistribute it and/or modify
7  it under the terms of the GNU General Public License as published by
8  the Free Software Foundation, either version 3 of the License, or
9  (at your option) any later version.
10 
11  OPM is distributed in the hope that it will be useful,
12  but WITHOUT ANY WARRANTY; without even the implied warranty of
13  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14  GNU General Public License for more details.
15 
16  You should have received a copy of the GNU General Public License
17  along with OPM. If not, see <http://www.gnu.org/licenses/>.
18 */
19 #ifndef OPM_ADAPTIVETIMESTEPPING_IMPL_HEADER_INCLUDED
20 #define OPM_ADAPTIVETIMESTEPPING_IMPL_HEADER_INCLUDED
21 
22 #include <iostream>
23 #include <string>
24 #include <utility>
25 
29 #include <dune/istl/istlexception.hh>
30 #include <dune/istl/ilu.hh> // For MatrixBlockException
31 
32 namespace Opm {
33 
34  // AdaptiveTimeStepping
35  //---------------------
36 
38  const boost::any& parallel_information,
39  const bool terminal_output )
40  : timeStepControl_()
41  , restart_factor_( param.getDefault("solver.restartfactor", double(0.33) ) )
42  , growth_factor_( param.getDefault("solver.growthfactor", double(2) ) )
43  , max_growth_( param.getDefault("timestep.control.maxgrowth", double(3.0) ) )
44  // default is 1 year, convert to seconds
45  , max_time_step_( unit::convert::from(param.getDefault("timestep.max_timestep_in_days", 365.0 ), unit::day) )
46  , solver_restart_max_( param.getDefault("solver.restart", int(10) ) )
47  , solver_verbose_( param.getDefault("solver.verbose", bool(true) ) && terminal_output )
48  , timestep_verbose_( param.getDefault("timestep.verbose", bool(true) ) && terminal_output )
49  , suggested_next_timestep_( -1.0 )
50  , full_timestep_initially_( param.getDefault("full_timestep_initially", bool(false) ) )
51  {
52  // valid are "pid" and "pid+iteration"
53  std::string control = param.getDefault("timestep.control", std::string("pid") );
54  // iterations is the accumulation of all linear iterations over all newton steops per time step
55  const int defaultTargetIterations = 30;
56 
57  const double tol = param.getDefault("timestep.control.tol", double(1e-1) );
58  if( control == "pid" ) {
59  timeStepControl_ = TimeStepControlType( new PIDTimeStepControl( tol, parallel_information ) );
60  }
61  else if ( control == "pid+iteration" )
62  {
63  const int iterations = param.getDefault("timestep.control.targetiteration", defaultTargetIterations );
64  timeStepControl_ = TimeStepControlType( new PIDAndIterationCountTimeStepControl( iterations, tol, parallel_information ) );
65  }
66  else if ( control == "iterationcount" )
67  {
68  const int iterations = param.getDefault("timestep.control.targetiteration", defaultTargetIterations );
69  const double decayrate = param.getDefault("timestep.control.decayrate", double(0.75) );
70  const double growthrate = param.getDefault("timestep.control.growthrate", double(1.25) );
71  timeStepControl_ = TimeStepControlType( new SimpleIterationCountTimeStepControl( iterations, decayrate, growthrate ) );
72  }
73  else
74  OPM_THROW(std::runtime_error,"Unsupported time step control selected "<< control );
75 
76  // make sure growth factor is something reasonable
77  assert( growth_factor_ >= 1.0 );
78  }
79 
80 
81  template <class Solver, class State, class WellState>
83  step( const SimulatorTimer& simulatorTimer, Solver& solver, State& state, WellState& well_state )
84  {
85  stepImpl( simulatorTimer, solver, state, well_state );
86  }
87 
88  template <class Solver, class State, class WellState>
90  step( const SimulatorTimer& simulatorTimer, Solver& solver, State& state, WellState& well_state,
91  OutputWriter& outputWriter )
92  {
93  stepImpl( simulatorTimer, solver, state, well_state, &outputWriter );
94  }
95 
96  // implementation of the step method
97  template <class Solver, class State, class WState>
99  stepImpl( const SimulatorTimer& simulatorTimer,
100  Solver& solver, State& state, WState& well_state,
101  OutputWriter* outputWriter )
102  {
103  const double timestep = simulatorTimer.currentStepLength();
104 
105  // init last time step as a fraction of the given time step
106  if( suggested_next_timestep_ < 0 ) {
108  }
109 
111  suggested_next_timestep_ = timestep;
112  }
113 
114  // TODO
115  // take change in well state into account
116 
117  // create adaptive step timer with previously used sub step size
118  AdaptiveSimulatorTimer substepTimer( simulatorTimer, suggested_next_timestep_, max_time_step_ );
119 
120  // copy states in case solver has to be restarted (to be revised)
121  State last_state( state );
122  WState last_well_state( well_state );
123 
124  // counter for solver restarts
125  int restarts = 0;
126 
127  // sub step time loop
128  while( ! substepTimer.done() )
129  {
130  // get current delta t
131  const double dt = substepTimer.currentStepLength() ;
132 
133  // initialize time step control in case current state is needed later
134  timeStepControl_->initialize( state );
135 
136  if( timestep_verbose_ )
137  {
138  std::cout <<"Substep( " << substepTimer.currentStepNum() << " ), try with stepsize "
139  << unit::convert::to(substepTimer.currentStepLength(), unit::day) << " (days)." << std::endl;
140  }
141 
142  int linearIterations = -1;
143  try {
144  // (linearIterations < 0 means on convergence in solver)
145  linearIterations = solver.step( dt, state, well_state);
146 
147  if( solver_verbose_ ) {
148  // report number of linear iterations
149  std::cout << "Overall linear iterations used: " << linearIterations << std::endl;
150  }
151  }
152  catch (const Opm::NumericalProblem& e) {
153  std::cerr << e.what() << std::endl;
154  // since linearIterations is < 0 this will restart the solver
155  }
156  catch (const std::runtime_error& e) {
157  std::cerr << e.what() << std::endl;
158  // also catch linear solver not converged
159  }
160  catch (const Dune::ISTLError& e) {
161  std::cerr << e.what() << std::endl;
162  // also catch errors in ISTL AMG that occur when time step is too large
163  }
164  catch (const Dune::MatrixBlockError& e) {
165  std::cerr << e.what() << std::endl;
166  // this can be thrown by ISTL's ILU0 in block mode, yet is not an ISTLError
167  }
168 
169  // (linearIterations < 0 means no convergence in solver)
170  if( linearIterations >= 0 )
171  {
172  // advance by current dt
173  ++substepTimer;
174 
175  // compute new time step estimate
176  double dtEstimate =
177  timeStepControl_->computeTimeStepSize( dt, linearIterations, state );
178 
179  // limit the growth of the timestep size by the growth factor
180  dtEstimate = std::min( dtEstimate, double(max_growth_ * dt) );
181 
182  // further restrict time step size growth after convergence problems
183  if( restarts > 0 ) {
184  dtEstimate = std::min( growth_factor_ * dt, dtEstimate );
185  // solver converged, reset restarts counter
186  restarts = 0;
187  }
188 
189  if( timestep_verbose_ )
190  {
191  std::cout << "Substep( " << substepTimer.currentStepNum()-1 // it was already advanced by ++
192  << " ) finished at time " << unit::convert::to(substepTimer.simulationTimeElapsed(),unit::day) << " (days)." << std::endl << std::endl;
193  }
194 
195  // write data if outputWriter was provided
196  if( outputWriter ) {
197  bool substep = true;
198  outputWriter->writeTimeStep( substepTimer, state, well_state, substep);
199  }
200 
201  // set new time step length
202  substepTimer.provideTimeStepEstimate( dtEstimate );
203 
204  // update states
205  last_state = state ;
206  last_well_state = well_state;
207 
208  }
209  else // in case of no convergence (linearIterations < 0)
210  {
211  // increase restart counter
212  if( restarts >= solver_restart_max_ ) {
213  OPM_THROW(Opm::NumericalProblem,"Solver failed to converge after " << restarts << " restarts.");
214  }
215 
216  const double newTimeStep = restart_factor_ * dt;
217  // we need to revise this
218  substepTimer.provideTimeStepEstimate( newTimeStep );
219  if( solver_verbose_ )
220  std::cerr << "Solver convergence failed, restarting solver with new time step ("
221  << unit::convert::to( newTimeStep, unit::day ) <<" days)." << std::endl;
222 
223  // reset states
224  state = last_state;
225  well_state = last_well_state;
226 
227  ++restarts;
228  }
229  }
230 
231 
232  // store estimated time step for next reportStep
234  if( timestep_verbose_ )
235  {
236  substepTimer.report( std::cout );
237  std::cout << "Suggested next step size = " << unit::convert::to( suggested_next_timestep_, unit::day ) << " (days)" << std::endl;
238  }
239 
240  if( ! std::isfinite( suggested_next_timestep_ ) ) { // check for NaN
241  suggested_next_timestep_ = timestep;
242  }
243  }
244 }
245 
246 #endif
const double day
Definition: Units.hpp:96
const double growth_factor_
factor to multiply time step when solver recovered from failed convergence
Definition: AdaptiveTimeStepping.hpp:85
Definition: TimeStepControl.hpp:75
Definition: AnisotropicEikonal.hpp:43
virtual void writeTimeStep(const SimulatorTimerInterface &timer, const SimulatorState &reservoirState, const WellState &wellState, bool isSubstep)=0
Write a blackoil reservoir state to disk for later inspection with visualization tools like ResInsigh...
Definition: ParameterGroup.hpp:109
double currentStepLength() const
double simulationTimeElapsed() const
void report(std::ostream &os) const
report start and end time as well as used steps so far
double suggested_next_timestep_
suggested size of next timestep
Definition: AdaptiveTimeStepping.hpp:91
void provideTimeStepEstimate(const double dt_estimate)
provide and estimate for new time step size
const bool solver_verbose_
solver verbosity
Definition: AdaptiveTimeStepping.hpp:89
const double max_time_step_
maximal allowed time step size
Definition: AdaptiveTimeStepping.hpp:87
T getDefault(const std::string &name, const T &default_value) const
This method is used to read a parameter from the parameter group.
Definition: ParameterGroup_impl.hpp:226
Simulation timer for adaptive time stepping.
Definition: AdaptiveSimulatorTimer.hpp:39
double to(const double q, const double unit)
Definition: Units.hpp:214
TimeStepControlType timeStepControl_
time step control object
Definition: AdaptiveTimeStepping.hpp:83
const double max_growth_
factor that limits the maximum growth of a time step
Definition: AdaptiveTimeStepping.hpp:86
double currentStepLength() const
Definition: OutputWriter.hpp:61
const bool timestep_verbose_
timestep verbosity
Definition: AdaptiveTimeStepping.hpp:90
std::unique_ptr< TimeStepControlInterface > TimeStepControlType
Definition: AdaptiveTimeStepping.hpp:81
The state of a set of wells.
Definition: WellState.hpp:36
const int solver_restart_max_
how many restart of solver are allowed
Definition: AdaptiveTimeStepping.hpp:88
Definition: TimeStepControl.hpp:147
void stepImpl(const SimulatorTimer &timer, Solver &solver, State &state, WellState &well_state, OutputWriter *outputWriter)
const double restart_factor_
factor to multiply time step with when solver fails to converge
Definition: AdaptiveTimeStepping.hpp:84
Definition: SimulatorTimer.hpp:34
Definition: TimeStepControl.hpp:38
double from(const double q, const double unit)
Definition: Units.hpp:191
AdaptiveTimeStepping(const parameter::ParameterGroup &param, const boost::any &pinfo=boost::any(), const bool terminal_output=true)
contructor taking parameter object
Definition: AdaptiveTimeStepping_impl.hpp:37
bool full_timestep_initially_
beginning with the size of the time step from data file
Definition: AdaptiveTimeStepping.hpp:92
void step(const SimulatorTimer &timer, Solver &solver, State &state, WellState &well_state)
step method that acts like the solver::step method in a sub cycle of time steps
Definition: AdaptiveTimeStepping_impl.hpp:83