the adaptive time stepping utility classes.
This commit is contained in:
parent
879dddc9b6
commit
5af49ed90b
@ -57,12 +57,13 @@ namespace Opm
|
||||
}
|
||||
|
||||
/// \brief advance time by currentStepLength
|
||||
void advance()
|
||||
AdaptiveSimulatorTimer& operator++ ()
|
||||
{
|
||||
++current_step_;
|
||||
current_time_ += dt_;
|
||||
// store used time step sizes
|
||||
steps_.push_back( dt_ );
|
||||
return *this;
|
||||
}
|
||||
|
||||
/// \brief provide and estimate for new time step size
|
||||
|
48
opm/core/simulator/AdaptiveTimeStepping.hpp
Normal file
48
opm/core/simulator/AdaptiveTimeStepping.hpp
Normal file
@ -0,0 +1,48 @@
|
||||
#ifndef OPM_SUBSTEPPING_HEADER_INCLUDED
|
||||
#define OPM_SUBSTEPPING_HEADER_INCLUDED
|
||||
|
||||
#include <iostream>
|
||||
#include <utility>
|
||||
|
||||
#include <opm/core/utility/parameters/ParameterGroup.hpp>
|
||||
#include <opm/core/utility/ErrorMacros.hpp>
|
||||
#include <opm/core/simulator/TimeStepControlInterface.hpp>
|
||||
|
||||
namespace Opm {
|
||||
|
||||
// AdaptiveTimeStepping
|
||||
//---------------------
|
||||
|
||||
class AdaptiveTimeStepping
|
||||
{
|
||||
public:
|
||||
//! \brief contructor taking parameter object
|
||||
AdaptiveTimeStepping( const parameter::ParameterGroup& param );
|
||||
|
||||
/** \brief step method that acts like the solver::step method
|
||||
in a sub cycle of time steps
|
||||
|
||||
\param solver solver object that must implement a method step( dt, state, well_state )
|
||||
\param state current state of the solution variables
|
||||
\param well_state additional well state object
|
||||
\param time current simulation time
|
||||
\param timestep current time step length that is to be sub cycled
|
||||
*/
|
||||
template <class Solver, class State, class WellState>
|
||||
void step( Solver& solver, State& state, WellState& well_state,
|
||||
const double time, const double timestep );
|
||||
|
||||
protected:
|
||||
typedef std::unique_ptr< TimeStepControlInterface > TimeStepControlType;
|
||||
|
||||
TimeStepControlType timeStepControl_; //!< time step control object
|
||||
const double restart_factor_; //!< factor to multiply time step with when solver fails to converge
|
||||
const int solver_restart_max_; //!< how many restart of solver are allowed
|
||||
const bool solver_verbose_; //!< solver verbosity
|
||||
const bool timestep_verbose_; //!< timestep verbosity
|
||||
double last_timestep_; //!< size of last timestep
|
||||
};
|
||||
}
|
||||
|
||||
#include <opm/core/simulator/AdaptiveTimeStepping_impl.hpp>
|
||||
#endif
|
130
opm/core/simulator/AdaptiveTimeStepping_impl.hpp
Normal file
130
opm/core/simulator/AdaptiveTimeStepping_impl.hpp
Normal file
@ -0,0 +1,130 @@
|
||||
#ifndef OPM_ADAPTIVETIMESTEPPING_IMPL_HEADER_INCLUDED
|
||||
#define OPM_ADAPTIVETIMESTEPPING_IMPL_HEADER_INCLUDED
|
||||
|
||||
#include <iostream>
|
||||
#include <string>
|
||||
#include <utility>
|
||||
|
||||
#include <opm/core/simulator/AdaptiveSimulatorTimer.hpp>
|
||||
#include <opm/core/simulator/PIDTimeStepControl.hpp>
|
||||
|
||||
namespace Opm {
|
||||
|
||||
// AdaptiveTimeStepping
|
||||
//---------------------
|
||||
|
||||
AdaptiveTimeStepping::AdaptiveTimeStepping( const parameter::ParameterGroup& param )
|
||||
: timeStepControl_()
|
||||
, restart_factor_( param.getDefault("solver.restartfactor", double(0.1) ) )
|
||||
, solver_restart_max_( param.getDefault("solver.restart", int(3) ) )
|
||||
, solver_verbose_( param.getDefault("solver.verbose", bool(false) ) )
|
||||
, timestep_verbose_( param.getDefault("timestep.verbose", bool(false) ) )
|
||||
, last_timestep_( std::numeric_limits< double >::max() )
|
||||
{
|
||||
// valid are "pid" and "pid+iteration"
|
||||
std::string control = param.getDefault("timestep.control", std::string("pid") );
|
||||
|
||||
const double tol = param.getDefault("timestep.control.tol", double(1e-3) );
|
||||
if( control == "pid" )
|
||||
timeStepControl_ = TimeStepControlType( new PIDTimeStepControl( tol ) );
|
||||
else if ( control == "pid+iteration" )
|
||||
{
|
||||
const int iterations = param.getDefault("timestep.control.targetiteration", int(25) );
|
||||
timeStepControl_ = TimeStepControlType( new PIDAndIterationCountTimeStepControl( iterations, tol ) );
|
||||
}
|
||||
else
|
||||
OPM_THROW(std::runtime_error,"Unsupported time step control selected "<< control );
|
||||
}
|
||||
|
||||
|
||||
template <class Solver, class State, class WellState>
|
||||
void AdaptiveTimeStepping::
|
||||
step( Solver& solver, State& state, WellState& well_state,
|
||||
const double time, const double timestep )
|
||||
{
|
||||
// create adaptive step timer with previously used sub step size
|
||||
AdaptiveSimulatorTimer timer( time, time+timestep, last_timestep_ );
|
||||
|
||||
// copy states in case solver has to be restarted (to be revised)
|
||||
State last_state( state );
|
||||
WellState last_well_state( well_state );
|
||||
|
||||
// counter for solver restarts
|
||||
int restarts = 0;
|
||||
|
||||
// sub step time loop
|
||||
while( ! timer.done() )
|
||||
{
|
||||
// initialize time step control in case current state is needed later
|
||||
timeStepControl_->initialize( state );
|
||||
|
||||
int linearIterations = -1;
|
||||
try {
|
||||
// (linearIterations < 0 means on convergence in solver)
|
||||
linearIterations = solver.step(timer.currentStepLength(), state, well_state);
|
||||
|
||||
if( solver_verbose_ ) {
|
||||
// report number of linear iterations
|
||||
std::cout << "Overall linear iterations used: " << linearIterations << std::endl;
|
||||
}
|
||||
}
|
||||
catch (Opm::NumericalProblem)
|
||||
{
|
||||
// since linearIterations is < 0 this will restart the solver
|
||||
}
|
||||
|
||||
// (linearIterations < 0 means on convergence in solver)
|
||||
if( linearIterations >= 0 )
|
||||
{
|
||||
// advance by current dt
|
||||
++timer;
|
||||
|
||||
// compute new time step estimate
|
||||
const double dtEstimate =
|
||||
timeStepControl_->computeTimeStepSize( timer.currentStepLength(), linearIterations, state );
|
||||
if( timestep_verbose_ )
|
||||
std::cout << "Suggested time step size = " << dtEstimate/86400.0 << " (days)" << std::endl;
|
||||
|
||||
// set new time step length
|
||||
timer.provideTimeStepEstimate( dtEstimate );
|
||||
|
||||
// update states
|
||||
last_state = state ;
|
||||
last_well_state = well_state;
|
||||
}
|
||||
else // in case of no convergence
|
||||
{
|
||||
// increase restart counter
|
||||
if( restarts >= solver_restart_max_ ) {
|
||||
OPM_THROW(Opm::NumericalProblem,"Solver failed to converge after " << restarts << " restarts.");
|
||||
}
|
||||
|
||||
const double newTimeStep = restart_factor_ * timer.currentStepLength();
|
||||
// we need to revise this
|
||||
timer.provideTimeStepEstimate( newTimeStep );
|
||||
if( solver_verbose_ )
|
||||
std::cerr << "Solver convergence failed, restarting solver with new time step ("<< newTimeStep <<" days)." << std::endl;
|
||||
|
||||
// reset states
|
||||
state = last_state;
|
||||
well_state = last_well_state;
|
||||
|
||||
++restarts;
|
||||
}
|
||||
}
|
||||
|
||||
|
||||
// store last small time step for next reportStep
|
||||
last_timestep_ = timer.suggestedAverage();
|
||||
if( timestep_verbose_ )
|
||||
{
|
||||
timer.report( std::cout );
|
||||
std::cout << "Last suggested step size = " << last_timestep_/86400.0 << " (days)" << std::endl;
|
||||
}
|
||||
|
||||
if( ! std::isfinite( last_timestep_ ) ) // check for NaN
|
||||
last_timestep_ = timestep;
|
||||
}
|
||||
}
|
||||
|
||||
#endif
|
121
opm/core/simulator/PIDTimeStepControl.cpp
Normal file
121
opm/core/simulator/PIDTimeStepControl.cpp
Normal file
@ -0,0 +1,121 @@
|
||||
/*
|
||||
Copyright 2014 IRIS AS
|
||||
|
||||
This file is part of the Open Porous Media project (OPM).
|
||||
|
||||
OPM is free software: you can redistribute it and/or modify
|
||||
it under the terms of the GNU General Public License as published by
|
||||
the Free Software Foundation, either version 3 of the License, or
|
||||
(at your option) any later version.
|
||||
|
||||
OPM is distributed in the hope that it will be useful,
|
||||
but WITHOUT ANY WARRANTY; without even the implied warranty of
|
||||
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
|
||||
GNU General Public License for more details.
|
||||
|
||||
You should have received a copy of the GNU General Public License
|
||||
along with OPM. If not, see <http://www.gnu.org/licenses/>.
|
||||
*/
|
||||
|
||||
#include <opm/core/simulator/PIDTimeStepControl.hpp>
|
||||
|
||||
namespace Opm
|
||||
{
|
||||
/// \brief constructor
|
||||
/// \param tol tolerance for the relative changes of the numerical solution to be accepted
|
||||
/// in one time step (default is 1e-3)
|
||||
PIDTimeStepControl::PIDTimeStepControl( const double, const bool verbose )
|
||||
: p0_()
|
||||
, sat0_()
|
||||
, tol_( tol )
|
||||
, errors_( 3, tol_ )
|
||||
, verbose_( verbose )
|
||||
{}
|
||||
|
||||
/// \brief \copydoc TimeStepControlInterface::initialize
|
||||
void PIDTimeStepControl::initialize( const SimulatorState& state )
|
||||
{
|
||||
// store current state for later time step computation
|
||||
p0_ = state.pressure();
|
||||
sat0_ = state.saturation();
|
||||
}
|
||||
|
||||
/// \brief \copydoc TimeStepControlInterface::computeTimeStepSize
|
||||
double PIDTimeStepControl::computeTimeStepSize( const double dt, const int /* iterations */, const SimulatorState& state ) const
|
||||
{
|
||||
const size_t size = p0_.size();
|
||||
assert( state.pressure().size() == size );
|
||||
assert( state.saturation().size() == size );
|
||||
assert( sat0_.size() == size );
|
||||
|
||||
// compute u^n - u^n+1
|
||||
for( size_t i=0; i<size; ++i )
|
||||
{
|
||||
p0_[ i ] -= state.pressure()[ i ];
|
||||
sat0_[ i ] -= state.saturation()[ i ];
|
||||
}
|
||||
|
||||
// compute || u^n - u^n+1 ||
|
||||
const double stateOld = inner_product( p0_.begin(), p0_.end() ) +
|
||||
inner_product( sat0_.begin(), sat0_.end() );
|
||||
|
||||
// compute || u^n+1 ||
|
||||
const double stateNew = inner_product( state.pressure().begin(), state.pressure().end() ) +
|
||||
inner_product( state.saturation().begin(), state.saturation().end() );
|
||||
|
||||
// shift errors
|
||||
for( int i=0; i<2; ++i )
|
||||
errors_[ i ] = errors_[i+1];
|
||||
|
||||
// store new error
|
||||
const double error = stateOld / stateNew;
|
||||
errors_[ 2 ] = error ;
|
||||
|
||||
if( error > tol_ )
|
||||
{
|
||||
// adjust dt by given tolerance
|
||||
if( verbose_ )
|
||||
std::cout << "Computed step size (tol): " << (dt * tol_ / error )/86400.0 << " (days)" << std::endl;
|
||||
return (dt * tol_ / error );
|
||||
}
|
||||
else
|
||||
{
|
||||
// values taking from turek time stepping paper
|
||||
const double kP = 0.075 ;
|
||||
const double kI = 0.175 ;
|
||||
const double kD = 0.01 ;
|
||||
double newDt = (dt * std::pow( errors_[ 1 ] / errors_[ 2 ], kP ) *
|
||||
std::pow( tol_ / errors_[ 2 ], kI ) *
|
||||
std::pow( errors_[0]*errors_[0]/errors_[ 1 ]/errors_[ 2 ], kD ));
|
||||
if( verbose_ )
|
||||
std::cout << "Computed step size (pow): " << newDt/86400.0 << " (days)" << std::endl;
|
||||
return newDt;
|
||||
}
|
||||
}
|
||||
|
||||
|
||||
PIDAndIterationCountTimeStepControl::
|
||||
PIDAndIterationCountTimeStepControl( const int target_iterations,
|
||||
const double tol,
|
||||
const bool verbose)
|
||||
: BaseType( tol, verbose )
|
||||
, target_iterations_( target_iterations )
|
||||
{}
|
||||
|
||||
double PIDAndIterationCountTimeStepControl::
|
||||
computeTimeStepSize( const double dt, const int iterations, const SimulatorState& state ) const
|
||||
{
|
||||
double dtEstimate = BaseType :: computeTimeStepSize( dt, iterations, state );
|
||||
|
||||
// further reduce step size if to many iterations were used
|
||||
if( iterations > target_iterations_ )
|
||||
{
|
||||
// if iterations was the same or dts were the same, do some magic
|
||||
dtEstimate *= double( target_iterations_ ) / double(iterations);
|
||||
}
|
||||
|
||||
return dtEstimate;
|
||||
}
|
||||
|
||||
} // end namespace Opm
|
||||
#endif
|
107
opm/core/simulator/PIDTimeStepControl.hpp
Normal file
107
opm/core/simulator/PIDTimeStepControl.hpp
Normal file
@ -0,0 +1,107 @@
|
||||
/*
|
||||
Copyright 2014 IRIS AS
|
||||
|
||||
This file is part of the Open Porous Media project (OPM).
|
||||
|
||||
OPM is free software: you can redistribute it and/or modify
|
||||
it under the terms of the GNU General Public License as published by
|
||||
the Free Software Foundation, either version 3 of the License, or
|
||||
(at your option) any later version.
|
||||
|
||||
OPM is distributed in the hope that it will be useful,
|
||||
but WITHOUT ANY WARRANTY; without even the implied warranty of
|
||||
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
|
||||
GNU General Public License for more details.
|
||||
|
||||
You should have received a copy of the GNU General Public License
|
||||
along with OPM. If not, see <http://www.gnu.org/licenses/>.
|
||||
*/
|
||||
#ifndef OPM_PIDTIMESTEPCONTROL_HEADER_INCLUDED
|
||||
#define OPM_PIDTIMESTEPCONTROL_HEADER_INCLUDED
|
||||
|
||||
#include <vector>
|
||||
|
||||
#include <opm/core/simulator/TimeStepControlInterface.hpp>
|
||||
|
||||
namespace Opm
|
||||
{
|
||||
///////////////////////////////////////////////////////////////////////////////////////////////////////////////
|
||||
///
|
||||
/// PID controller based adaptive time step control as suggested in:
|
||||
/// Turek and Kuzmin. Algebraic Flux Correction III. Incompressible Flow Problems. Uni Dortmund.
|
||||
///
|
||||
/// See also:
|
||||
/// D. Kuzmin and S.Turek. Numerical simulation of turbulent bubbly flows. Techreport Uni Dortmund. 2004
|
||||
///
|
||||
/// and the original article:
|
||||
/// Valli, Coutinho, and Carey. Adaptive Control for Time Step Selection in Finite Element
|
||||
/// Simulation of Coupled Viscous Flow and Heat Transfer. Proc of the 10th
|
||||
/// International Conference on Numerical Methods in Fluids. 1998.
|
||||
///
|
||||
///////////////////////////////////////////////////////////////////////////////////////////////////////////////
|
||||
class PIDTimeStepControl : public TimeStepControlInterface
|
||||
{
|
||||
public:
|
||||
/// \brief constructor
|
||||
/// \param tol tolerance for the relative changes of the numerical solution to be accepted
|
||||
/// in one time step (default is 1e-3)
|
||||
/// \param verbose if true get some output (default = false)
|
||||
PIDTimeStepControl( const double tol = 1e-3, const bool verbose = false );
|
||||
|
||||
/// \brief \copydoc TimeStepControlInterface::initialize
|
||||
void initialize( const SimulatorState& state );
|
||||
|
||||
/// \brief \copydoc TimeStepControlInterface::computeTimeStepSize
|
||||
double computeTimeStepSize( const double dt, const int /* iterations */, const SimulatorState& state ) const;
|
||||
|
||||
protected:
|
||||
// return inner product for given container, here std::vector
|
||||
template <class Iterator>
|
||||
double inner_product( Iterator it, const Iterator end ) const
|
||||
{
|
||||
double product = 0.0 ;
|
||||
for( ; it != end; ++it )
|
||||
product += ( *it * *it );
|
||||
return product;
|
||||
}
|
||||
|
||||
protected:
|
||||
mutable std::vector<double> p0_;
|
||||
mutable std::vector<double> sat0_;
|
||||
|
||||
const double tol_;
|
||||
mutable std::vector< double > errors_;
|
||||
|
||||
const bool verbose_;
|
||||
};
|
||||
|
||||
///////////////////////////////////////////////////////////////////////////////////////////////////////////////
|
||||
///
|
||||
/// PID controller based adaptive time step control as above that also takes
|
||||
/// an target iteration into account.
|
||||
//
|
||||
///////////////////////////////////////////////////////////////////////////////////////////////////////////////
|
||||
class PIDAndIterationCountTimeStepControl : public PIDTimeStepControl
|
||||
{
|
||||
typedef PIDTimeStepControl BaseType;
|
||||
public:
|
||||
/// \brief constructor
|
||||
/// \param target_iterations number of desired iterations per time step
|
||||
/// \param tol tolerance for the relative changes of the numerical solution to be accepted
|
||||
/// in one time step (default is 1e-3)
|
||||
/// \param verbose if true get some output (default = false)
|
||||
PIDAndIterationCountTimeStepControl( const int target_iterations = 20,
|
||||
const double tol = 1e-3,
|
||||
const bool verbose = false);
|
||||
|
||||
/// \brief \copydoc TimeStepControlInterface::computeTimeStepSize
|
||||
double computeTimeStepSize( const double dt, const int iterations, const SimulatorState& state ) const;
|
||||
|
||||
protected:
|
||||
const int target_iterations_;
|
||||
};
|
||||
|
||||
|
||||
} // end namespace Opm
|
||||
#endif
|
||||
|
51
opm/core/simulator/TimeStepControlInterface.hpp
Normal file
51
opm/core/simulator/TimeStepControlInterface.hpp
Normal file
@ -0,0 +1,51 @@
|
||||
/*
|
||||
Copyright 2014 IRIS AS
|
||||
|
||||
This file is part of the Open Porous Media project (OPM).
|
||||
|
||||
OPM is free software: you can redistribute it and/or modify
|
||||
it under the terms of the GNU General Public License as published by
|
||||
the Free Software Foundation, either version 3 of the License, or
|
||||
(at your option) any later version.
|
||||
|
||||
OPM is distributed in the hope that it will be useful,
|
||||
but WITHOUT ANY WARRANTY; without even the implied warranty of
|
||||
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
|
||||
GNU General Public License for more details.
|
||||
|
||||
You should have received a copy of the GNU General Public License
|
||||
along with OPM. If not, see <http://www.gnu.org/licenses/>.
|
||||
*/
|
||||
#ifndef OPM_TIMESTEPCONTROLINTERFACE_HEADER_INCLUDED
|
||||
#define OPM_TIMESTEPCONTROLINTERFACE_HEADER_INCLUDED
|
||||
|
||||
namespace Opm
|
||||
{
|
||||
|
||||
///////////////////////////////////////////////////////////////////
|
||||
///
|
||||
/// TimeStepControlInterface
|
||||
///
|
||||
///////////////////////////////////////////////////////////////////
|
||||
class TimeStepControlInterface
|
||||
{
|
||||
protected:
|
||||
TimeStepControlInterface() {}
|
||||
public:
|
||||
/// \param state simulation state before computing update in the solver (default is empty)
|
||||
virtual void initialize( const SimulatorState& state ) {}
|
||||
|
||||
/// compute new time step size suggestions based on the PID controller
|
||||
/// \param dt time step size used in the current step
|
||||
/// \param iterations number of iterations used (linear/nonlinear)
|
||||
/// \param state new solution state
|
||||
///
|
||||
/// \return suggested time step size for the next step
|
||||
virtual double computeTimeStepSize( const double dt, const int iterations, const SimulatorState& ) const = 0;
|
||||
|
||||
/// virtual destructor (empty)
|
||||
virtual ~TimeStepControlInterface () {}
|
||||
};
|
||||
|
||||
}
|
||||
#endif
|
Loading…
Reference in New Issue
Block a user