EclProblem: clean up the time management code

- when an episode/report step is over, the next is started by endEpisode()
- the problem does not deal with updating the simulation time anymore
- rename `episodeIdx` in to `reportStepIdx` the 'EclWriter' because
  this variable is -- and always has been -- the report step number
  used by some parts of `opm-output`'s ECL writing code (the report
  step number is equivalent to the episode index plus 1). IMO, the
  output and parser code should be made more consistent in regard of
  whether it expects 0-based or 1-based indices, but this is a story
  for another day.
This commit is contained in:
Andreas Lauser 2019-04-03 17:26:57 +02:00
parent 17a4092c82
commit d329547463
8 changed files with 109 additions and 105 deletions

View File

@ -610,6 +610,19 @@ public:
ParentType::finishInit();
auto& simulator = this->simulator();
const auto& eclState = simulator.vanguard().eclState();
const auto& schedule = simulator.vanguard().schedule();
const auto& timeMap = schedule.getTimeMap();
// Set the start time of the simulation
simulator.setStartTime(timeMap.getStartTime(/*reportStepIdx=*/0));
// We want the episode index to be the same as the report step index to make
// things simpler, so we have to set the episode index to -1 because it is
// incremented by endEpisode(). The size of the initial time step and
// length of the initial episode is set to zero for the same reason.
simulator.setEpisodeIndex(-1);
simulator.setEpisodeLength(0.0);
// the "NOGRAV" keyword from Frontsim or setting the EnableGravity to false
// disables gravity, else the standard value of the gravity constant at sea level
@ -634,7 +647,6 @@ public:
}
// deal with DRSDT
const auto& eclState = simulator.vanguard().eclState();
unsigned ntpvt = eclState.runspec().tabdims().getNumPVTTables();
maxDRs_.resize(ntpvt, 1e30);
dRsDtOnlyFreeGas_.resize(ntpvt, false);
@ -651,19 +663,7 @@ public:
readThermalParameters_();
transmissibilities_.finishInit();
const auto& timeMap = simulator.vanguard().schedule().getTimeMap();
readInitialCondition_();
// Set the start time of the simulation
simulator.setStartTime(timeMap.getStartTime(/*timeStepIdx=*/0));
// We want the episode index to be the same as the report step index to make
// things simpler, so we have to set the episode index to -1 because it is
// incremented inside beginEpisode(). The size of the initial time step and
// length of the initial episode is set to zero for the same reason.
simulator.setEpisodeIndex(-1);
simulator.setEpisodeLength(0.0);
simulator.setTimeStepSize(0.0);
updatePffDofData_();
@ -691,6 +691,11 @@ public:
eclWriter_->writeInit();
simulator.vanguard().releaseGlobalTransmissibilities();
// after finishing the initialization and writing the initial solution, we move
// to the first "real" episode/report step
simulator.startNextEpisode(timeMap.getTimeStepLength(0));
simulator.setEpisodeIndex(0);
}
void prefetch(const Element& elem) const
@ -739,37 +744,21 @@ public:
void beginEpisode(bool isOnRestart = false)
{
// Proceed to the next report step
Simulator& simulator = this->simulator();
auto& simulator = this->simulator();
auto& eclState = simulator.vanguard().eclState();
const auto& schedule = simulator.vanguard().schedule();
const auto& events = schedule.getEvents();
const auto& timeMap = schedule.getTimeMap();
int episodeIdx = simulator.episodeIndex();
// The first thing to do in the morning of an episode is update update the
// eclState and the deck if they need to be changed.
int nextEpisodeIdx = simulator.episodeIndex();
if (enableExperiments && this->gridView().comm().rank() == 0) {
boost::posix_time::ptime curDateTime =
boost::posix_time::from_time_t(timeMap.getStartTime(nextEpisodeIdx+1));
std::cout << "Report step " << nextEpisodeIdx + 2
<< "/" << timeMap.numTimesteps()
<< " at day " << timeMap.getTimePassedUntil(nextEpisodeIdx+1)/(24*3600)
<< "/" << timeMap.getTotalTime()/(24*3600)
<< ", date = " << curDateTime.date()
<< "\n ";
}
if (nextEpisodeIdx > 0 &&
events.hasEvent(Opm::ScheduleEvents::GEO_MODIFIER, nextEpisodeIdx))
{
if (episodeIdx >= 0 && events.hasEvent(Opm::ScheduleEvents::GEO_MODIFIER, episodeIdx)) {
// bring the contents of the keywords to the current state of the SCHEDULE
// section
// section.
//
// TODO (?): make grid topology changes possible (depending on what exactly
// has changed, the grid may need be re-created which has some serious
// implications on e.g., the solution of the simulation.)
const auto& miniDeck = schedule.getModifierDeck(nextEpisodeIdx);
const auto& miniDeck = schedule.getModifierDeck(episodeIdx);
eclState.applyModifierDeck(miniDeck);
// re-compute all quantities which may possibly be affected.
@ -778,44 +767,32 @@ public:
updatePffDofData_();
}
// Opm::TimeMap deals with points in time, so the number of time intervals (i.e.,
// report steps) is one less!
int numReportSteps = timeMap.size() - 1;
// start the next episode if there are additional report steps, else finish the
// simulation
while (nextEpisodeIdx < numReportSteps &&
simulator.time() >= timeMap.getTimePassedUntil(nextEpisodeIdx + 1)*(1 - 1e-10))
{
++ nextEpisodeIdx;
if (enableExperiments && this->gridView().comm().rank() == 0 && episodeIdx >= 0) {
// print some useful information in experimental mode. (the production
// simulator does this externally.)
boost::posix_time::ptime curDateTime =
boost::posix_time::from_time_t(timeMap.getStartTime(episodeIdx));
std::cout << "Report step " << episodeIdx + 1
<< "/" << timeMap.numTimesteps()
<< " at day " << timeMap.getTimePassedUntil(episodeIdx)/(24*3600)
<< "/" << timeMap.getTotalTime()/(24*3600)
<< ", date = " << curDateTime.date()
<< "\n ";
}
// react to TUNING changes
bool tuningEvent = false;
if (nextEpisodeIdx > 0
&& enableTuning_
&& events.hasEvent(Opm::ScheduleEvents::TUNING_CHANGE, nextEpisodeIdx))
if (episodeIdx > 0 && enableTuning_ && events.hasEvent(Opm::ScheduleEvents::TUNING_CHANGE, episodeIdx))
{
const auto& tuning = schedule.getTuning();
initialTimeStepSize_ = tuning.getTSINIT(nextEpisodeIdx);
maxTimeStepAfterWellEvent_ = tuning.getTMAXWC(nextEpisodeIdx);
maxTimeStepSize_ = tuning.getTSMAXZ(nextEpisodeIdx);
restartShrinkFactor_ = 1./tuning.getTSFCNV(nextEpisodeIdx);
minTimeStepSize_ = tuning.getTSMINZ(nextEpisodeIdx);
initialTimeStepSize_ = tuning.getTSINIT(episodeIdx);
maxTimeStepAfterWellEvent_ = tuning.getTMAXWC(episodeIdx);
maxTimeStepSize_ = tuning.getTSMAXZ(episodeIdx);
restartShrinkFactor_ = 1./tuning.getTSFCNV(episodeIdx);
minTimeStepSize_ = tuning.getTSMINZ(episodeIdx);
tuningEvent = true;
}
Scalar episodeLength = timeMap.getTimeStepLength(nextEpisodeIdx);
simulator.startNextEpisode(episodeLength);
Scalar dt = limitNextTimeStepSize_(episodeLength);
if (nextEpisodeIdx == 0 || tuningEvent)
// allow the size of the initial time step to be set via an external parameter
// if TUNING is enabled, also limit the time step size after a tuning event to TSINIT
dt = std::min(dt, initialTimeStepSize_);
simulator.setTimeStepSize(dt);
const bool invalidateFromHyst = updateHysteresis_();
const bool invalidateFromMaxOilSat = updateMaxOilSaturation_();
const bool doInvalidate = invalidateFromHyst || invalidateFromMaxOilSat;
@ -824,10 +801,19 @@ public:
updateMaxPolymerAdsorption_();
// set up the wells for the next episode.
wellModel_.beginEpisode(isOnRestart);
wellModel_.beginEpisode();
// set up the aquifers for the next episode.
aquiferModel_.beginEpisode();
// set the size of the initial time step of the episode
Scalar dt = limitNextTimeStepSize_(simulator.episodeLength());
if (episodeIdx == 0 || tuningEvent)
// allow the size of the initial time step to be set via an external parameter
// if TUNING is enabled, also limit the time step size after a tuning event to TSINIT
dt = std::min(dt, initialTimeStepSize_);
simulator.setTimeStepSize(dt);
if (doInvalidate)
this->model().invalidateIntensiveQuantitiesCache(/*timeIdx=*/0);
}
@ -927,15 +913,18 @@ public:
{
auto& simulator = this->simulator();
const auto& schedule = simulator.vanguard().schedule();
const auto& timeMap = schedule.getTimeMap();
int episodeIdx = simulator.episodeIndex();
const auto& timeMap = schedule.getTimeMap();
int numReportSteps = timeMap.size() - 1;
if (episodeIdx + 1 >= numReportSteps) {
// check if we're finished ...
if (episodeIdx + 1 >= static_cast<int>(timeMap.numTimesteps())) {
simulator.setFinished(true);
return;
}
// .. if we're not yet done, start the next episode (report step)
simulator.startNextEpisode(timeMap.getTimeStepLength(episodeIdx + 1));
}
/*!
@ -1442,7 +1431,8 @@ public:
void initialSolutionApplied()
{
const auto& simulator = this->simulator();
const auto& eclState = simulator.vanguard().eclState();
const auto& vanguard = simulator.vanguard();
const auto& eclState = vanguard.eclState();
// initialize the wells. Note that this needs to be done after initializing the
// intrinsic permeabilities and the after applying the initial solution because
@ -2109,9 +2099,13 @@ private:
simulator.setStartTime(timeMap.getStartTime(/*timeStepIdx=*/0));
simulator.setTime(timeMap.getTimePassedUntil(episodeIdx));
simulator.startNextEpisode(simulator.startTime() + simulator.time(),
timeMap.getTimeStepLength(episodeIdx));
simulator.setEpisodeIndex(episodeIdx);
simulator.setEpisodeLength(timeMap.getTimeStepLength(episodeIdx));
simulator.setTimeStepSize(eclWriter_->restartTimeStepSize());
Scalar dt = std::min(eclWriter_->restartTimeStepSize(), simulator.episodeLength());
simulator.setTimeStepSize(dt);
eclWriter_->beginRestart();
@ -2195,7 +2189,7 @@ private:
initialFluidStates_ = tmpInitialFs;
// make sure that the stuff which needs to be done at the beginning of an episode
// is run.
this->beginEpisode(/*isOnRestart=*/true);
this->beginEpisode();
eclWriter_->endRestart();
}

View File

@ -592,7 +592,7 @@ public:
void deserialize(Restarter& res OPM_UNUSED)
{
// initialize the wells for the current episode
beginEpisode(simulator_.vanguard().eclState(), simulator_.vanguard().schedule(), /*wasRestarted=*/true);
beginEpisode(/*wasRestarted=*/true);
}
/*!

View File

@ -160,11 +160,11 @@ public:
// output using eclWriter if enabled
Opm::data::Wells localWellData = simulator_.problem().wellModel().wellData();
int episodeIdx = simulator_.episodeIndex() + 1;
int reportStepNum = simulator_.episodeIndex() + 1;
const auto& gridView = simulator_.vanguard().gridView();
int numElements = gridView.size(/*codim=*/0);
bool log = collectToIORank_.isIORank();
eclOutputModule_.allocBuffers(numElements, episodeIdx, isSubStep, log);
eclOutputModule_.allocBuffers(numElements, reportStepNum, isSubStep, log);
ElementContext elemCtx(simulator_);
ElementIterator elemIt = gridView.template begin</*codim=*/0>();
@ -184,7 +184,7 @@ public:
// add cell data to perforations for Rft output
if (!isSubStep)
eclOutputModule_.addRftDataToWells(localWellData, episodeIdx);
eclOutputModule_.addRftDataToWells(localWellData, reportStepNum);
if (collectToIORank_.isParallel())
collectToIORank_.collect(localCellData, eclOutputModule_.getBlockData(), localWellData);
@ -221,7 +221,7 @@ public:
// first, create a tasklet to write the data for the current time step to disk
auto eclWriteTasklet = std::make_shared<EclWriteTasklet>(*eclIO_,
episodeIdx,
reportStepNum,
isSubStep,
curTime,
restartValue,
@ -437,7 +437,7 @@ private:
: public TaskletInterface
{
Opm::EclipseIO& eclIO_;
int episodeIdx_;
int reportStepNum_;
bool isSubStep_;
double secondsElapsed_;
Opm::RestartValue restartValue_;
@ -447,7 +447,7 @@ private:
bool writeDoublePrecision_;
explicit EclWriteTasklet(Opm::EclipseIO& eclIO,
int episodeIdx,
int reportStepNum,
bool isSubStep,
double secondsElapsed,
Opm::RestartValue restartValue,
@ -456,7 +456,7 @@ private:
const std::map<std::pair<std::string, int>, double>& blockSummaryValues,
bool writeDoublePrecision)
: eclIO_(eclIO)
, episodeIdx_(episodeIdx)
, reportStepNum_(reportStepNum)
, isSubStep_(isSubStep)
, secondsElapsed_(secondsElapsed)
, restartValue_(restartValue)
@ -469,7 +469,7 @@ private:
// callback to eclIO serial writeTimeStep method
void run()
{
eclIO_.writeTimeStep(episodeIdx_,
eclIO_.writeTimeStep(reportStepNum_,
isSubStep_,
secondsElapsed_,
restartValue_,

View File

@ -182,7 +182,6 @@ namespace Opm {
/// \param[in] timer simulation timer
void prepareStep(const SimulatorTimerInterface& timer)
{
// update the solution variables in ebos
if ( timer.lastStepFailed() ) {
ebosSimulator_.model().updateFailed();
@ -194,13 +193,8 @@ namespace Opm {
// know the report step/episode index because of timing dependend data
// despide the fact that flow uses its own time stepper. (The length of the
// episode does not matter, though.)
Scalar t = timer.simulationTimeElapsed();
ebosSimulator_.startNextEpisode(/*episodeStartTime=*/t, /*episodeLength=*/1e30);
ebosSimulator_.setEpisodeIndex(timer.reportStepNum());
ebosSimulator_.setTime(t);
ebosSimulator_.setTime(timer.simulationTimeElapsed());
ebosSimulator_.setTimeStepSize(timer.currentStepLength());
ebosSimulator_.setTimeStepIndex(ebosSimulator_.timeStepIndex() + 1);
ebosSimulator_.problem().beginTimeStep();
unsigned numDof = ebosSimulator_.model().numGridDof();

View File

@ -153,16 +153,9 @@ namespace Opm {
// TODO (?)
}
void beginEpisode(bool isRestart)
void beginEpisode()
{
size_t episodeIdx = ebosSimulator_.episodeIndex();
// beginEpisode in eclProblem advances the episode index
// we don't want this when we are at the beginning of an
// restart.
if (isRestart)
episodeIdx -= 1;
beginReportStep(episodeIdx);
beginReportStep(ebosSimulator_.episodeIndex());
}
void beginTimeStep();

View File

@ -41,7 +41,6 @@ namespace Opm {
init()
{
const Opm::EclipseState& eclState = ebosSimulator_.vanguard().eclState();
const Opm::Schedule& schedule = ebosSimulator_.vanguard().schedule();
gravity_ = ebosSimulator_.problem().gravity()[2];

View File

@ -133,6 +133,8 @@ public:
{
failureReport_ = SimulatorReport();
ebosSimulator_.setEpisodeIndex(-1);
// handle restarts
std::unique_ptr<RestartValue> restartValues;
if (isRestart()) {
@ -164,6 +166,21 @@ public:
double suggestedStepSize = -1.0;
if (isRestart()) {
// Set the start time of the simulation
const auto& schedule = ebosSimulator_.vanguard().schedule();
const auto& eclState = ebosSimulator_.vanguard().eclState();
const auto& timeMap = schedule.getTimeMap();
const auto& initconfig = eclState.getInitConfig();
int episodeIdx = initconfig.getRestartStep() - 1;
ebosSimulator_.setStartTime(timeMap.getStartTime(/*timeStepIdx=*/0));
ebosSimulator_.setTime(timeMap.getTimePassedUntil(episodeIdx));
ebosSimulator_.startNextEpisode(ebosSimulator_.startTime() + ebosSimulator_.time(),
timeMap.getTimeStepLength(episodeIdx));
ebosSimulator_.setEpisodeIndex(episodeIdx);
// This is a restart, determine the time step size from the restart data
if (restartValues->hasExtra("OPMEXTRA")) {
std::vector<double> opmextra = restartValues->getExtra("OPMEXTRA");
@ -206,12 +223,18 @@ public:
Dune::Timer perfTimer;
perfTimer.start();
ebosSimulator_.setEpisodeIndex(-1);
ebosSimulator_.setEpisodeLength(0.0);
ebosSimulator_.setTimeStepSize(0.0);
wellModel_().beginReportStep(timer.currentStepNum());
ebosSimulator_.problem().writeOutput(false);
report.output_write_time += perfTimer.stop();
}
ebosSimulator_.setEpisodeIndex(timer.currentStepNum());
// Run a multiple steps of the solver depending on the time step control.
solverTimer.start();
@ -265,6 +288,14 @@ public:
}
}
// write simulation state at the report stage
Dune::Timer perfTimer;
perfTimer.start();
const double nextstep = adaptiveTimeStepping ? adaptiveTimeStepping->suggestedNextStep() : -1.0;
ebosSimulator_.problem().setNextTimeStepSize(nextstep);
ebosSimulator_.problem().writeOutput(false);
report.output_write_time += perfTimer.stop();
solver->model().endReportStep();
// take time that was used to solve system for this reportStep
@ -284,14 +315,6 @@ public:
}
}
// write simulation state at the report stage
Dune::Timer perfTimer;
perfTimer.start();
const double nextstep = adaptiveTimeStepping ? adaptiveTimeStepping->suggestedNextStep() : -1.0;
ebosSimulator_.problem().setNextTimeStepSize(nextstep);
ebosSimulator_.problem().writeOutput(false);
report.output_write_time += perfTimer.stop();
if (terminalOutput_) {
std::string msg =
"Time step took " + std::to_string(solverTimer.secsSinceStart()) + " seconds; "

View File

@ -131,6 +131,7 @@ void test_summary()
simulator->model().applyInitialSolution();
Opm::data::Wells dw;
bool substep = false;
simulator->startNextEpisode(0.0, 1e30);
simulator->setEpisodeIndex(0);
eclWriter->writeOutput(substep);
simulator->setEpisodeIndex(1);