28#ifndef EWOMS_SIMULATOR_HH
29#define EWOMS_SIMULATOR_HH
32#define RESERVOIR_COUPLING_ENABLED
34#ifdef RESERVOIR_COUPLING_ENABLED
35#include <opm/simulators/flow/ReservoirCouplingMaster.hpp>
36#include <opm/simulators/flow/ReservoirCouplingSlave.hpp>
39#include <dune/common/parallel/mpihelper.hh>
48#include <opm/models/utils/simulatorutils.hpp>
52#include <opm/simulators/utils/DeferredLoggingErrorHelpers.hpp>
66 static constexpr T constexpr_max(T
a, T
b) {
67 return (
a >
b) ?
a :
b;
82template <
class TypeTag>
91 using MPIComm =
typename Dune::MPIHelper::MPICommunicator;
92 using Communication = Dune::Communication<MPIComm>;
95 static constexpr Scalar eps =
96 constexpr_max(std::numeric_limits<Scalar>::epsilon(),
static_cast<Scalar
>(1.0e-9));
113 verbose_ =
verbose && comm.rank() == 0;
118 endTime_ = Parameters::Get<Parameters::EndTime<Scalar>>();
119 timeStepSize_ = Parameters::Get<Parameters::InitialTimeStepSize<Scalar>>();
120 assert(timeStepSize_ > 0);
122 Parameters::Get<Parameters::PredeterminedTimeStepsFile>();
128 episodeStartTime_ = 0;
129 episodeLength_ = std::numeric_limits<Scalar>::max();
134 std::cout <<
"Allocating the simulation vanguard\n" << std::flush;
138 OPM_BEGIN_PARALLEL_TRY_CATCH();
139 vanguard_ = std::make_unique<Vanguard>(*
this);
140 OPM_END_PARALLEL_TRY_CATCH(
"Allocating the simulation vanguard failed: ", comm);
144 std::cout <<
"Distributing the vanguard's data\n" << std::flush;
148 OPM_BEGIN_PARALLEL_TRY_CATCH();
149 vanguard_->loadBalance();
150 OPM_END_PARALLEL_TRY_CATCH(
"Could not distribute the vanguard data: ", comm);
155 std::cout <<
"Adding LGRs, if any, in serial run\n" << std::flush;
159 OPM_BEGIN_PARALLEL_TRY_CATCH();
160 vanguard_->addLgrs();
161 OPM_END_PARALLEL_TRY_CATCH(
"Adding LGRs to the simulation vanguard in serial run failed: ", comm);
165 std::cout <<
"Allocating the model\n" << std::flush;
169 OPM_BEGIN_PARALLEL_TRY_CATCH();
170 model_ = std::make_unique<Model>(*
this);
171 OPM_END_PARALLEL_TRY_CATCH(
"Could not allocate model: ", comm);
175 std::cout <<
"Allocating the problem\n" << std::flush;
179 OPM_BEGIN_PARALLEL_TRY_CATCH();
180 problem_ = std::make_unique<Problem>(*
this);
181 OPM_END_PARALLEL_TRY_CATCH(
"Could not allocate the problem: ", comm);
185 std::cout <<
"Initializing the model\n" << std::flush;
189 OPM_BEGIN_PARALLEL_TRY_CATCH();
190 model_->finishInit();
191 OPM_END_PARALLEL_TRY_CATCH(
"Could not initialize the model: ", comm);
195 std::cout <<
"Initializing the problem\n" << std::flush;
199 OPM_BEGIN_PARALLEL_TRY_CATCH();
200 problem_->finishInit();
201 OPM_END_PARALLEL_TRY_CATCH(
"Could not initialize the problem: ", comm);
207 std::cout <<
"Simulator successfully set up\n" << std::flush;
216 Parameters::Register<Parameters::EndTime<Scalar>>
217 (
"The simulation time at which the simulation is finished [s]");
218 Parameters::Register<Parameters::InitialTimeStepSize<Scalar>>
219 (
"The size of the initial time step [s]");
220 Parameters::Register<Parameters::RestartTime<Scalar>>
221 (
"The simulation time at which a restart should be attempted [s]");
222 Parameters::Register<Parameters::PredeterminedTimeStepsFile>
223 (
"A file with a list of predetermined time step sizes "
224 "(one time step per line)");
226 Vanguard::registerParameters();
227 Model::registerParameters();
228 Problem::registerParameters();
235 {
return *vanguard_; }
241 {
return *vanguard_; }
247 {
return vanguard_->gridView(); }
266 {
return *problem_; }
273 {
return *problem_; }
287 {
return startTime_; }
340 {
return setupTimer_; }
347 {
return executionTimer_; }
350 {
return executionTimer_; }
357 {
return prePostProcessTimer_; }
364 {
return linearizeTimer_; }
371 {
return solveTimer_; }
378 {
return updateTimer_; }
385 {
return writeTimer_; }
398 { timeStepSize_ = value; }
406 { timeStepIdx_ = value; }
414 {
return timeStepSize_; }
421 {
return timeStepIdx_; }
431 { finished_ =
yesno; }
441 assert(timeStepSize_ >= 0.0);
442 return finished_ || (this->
time() * (1.0 + eps) >=
endTime());
451 return finished_ || (this->
time() + timeStepSize_) * (1.0 + eps) >=
endTime();
491 episodeStartTime_ = startTime_ + time_;
492 episodeLength_ =
len;
509 {
return episodeIdx_; }
516 {
return episodeStartTime_; }
524 { episodeLength_ =
dt; }
531 {
return episodeLength_; }
539 return this->
time() <= (episodeStartTime_ -
startTime()) * (1 + eps);
577 return std::max<Scalar>(0.0,
601 const Scalar
restartTime = Parameters::Get<Parameters::RestartTime<Scalar>>();
606 OPM_BEGIN_PARALLEL_TRY_CATCH();
611 std::cout <<
"Deserialize from file '" <<
res.fileName() <<
"'\n" << std::flush;
615 problem_->deserialize(
res);
616 model_->deserialize(
res);
617 res.deserializeEnd();
618 OPM_END_PARALLEL_TRY_CATCH(
"Deserialization failed: ",
619 Dune::MPIHelper::getCommunication());
621 std::cout <<
"Deserialization done."
625 <<
"\n" << std::flush;
631 std::cout <<
"Applying the initial solution of the \"" << problem_->name()
632 <<
"\" problem\n" << std::flush;
641 OPM_BEGIN_PARALLEL_TRY_CATCH();
642 model_->applyInitialSolution();
643 OPM_END_PARALLEL_TRY_CATCH(
"Apply initial solution failed: ",
644 Dune::MPIHelper::getCommunication());
648 if (problem_->shouldWriteOutput()) {
649 OPM_BEGIN_PARALLEL_TRY_CATCH();
650 problem_->writeOutput(
true);
651 OPM_END_PARALLEL_TRY_CATCH(
"Write output failed: ",
652 Dune::MPIHelper::getCommunication());
660 executionTimer_.
start();
664 prePostProcessTimer_.
start();
669 OPM_BEGIN_PARALLEL_TRY_CATCH();
670 problem_->beginEpisode();
671 OPM_END_PARALLEL_TRY_CATCH(
"Begin episode failed: ",
672 Dune::MPIHelper::getCommunication());
678 OPM_BEGIN_PARALLEL_TRY_CATCH();
679 problem_->endEpisode();
680 OPM_END_PARALLEL_TRY_CATCH(
"End episode failed: ",
681 Dune::MPIHelper::getCommunication());
682 prePostProcessTimer_.
stop();
690 std::cout <<
"Begin time step " <<
timeStepIndex() + 1 <<
". "
698 OPM_BEGIN_PARALLEL_TRY_CATCH();
699 problem_->beginTimeStep();
700 OPM_END_PARALLEL_TRY_CATCH(
"Begin timestep failed: ",
701 Dune::MPIHelper::getCommunication());
707 OPM_BEGIN_PARALLEL_TRY_CATCH();
708 problem_->endTimeStep();
709 problem_->endEpisode();
710 OPM_END_PARALLEL_TRY_CATCH(
"Finish failed: ",
711 Dune::MPIHelper::getCommunication());
712 prePostProcessTimer_.
stop();
716 prePostProcessTimer_.
stop();
720 problem_->timeIntegration();
725 const auto&
pmodel = problem_->model();
726 prePostProcessTimer_ +=
pmodel.prePostProcessTimer();
727 linearizeTimer_ +=
pmodel.linearizeTimer();
728 solveTimer_ +=
pmodel.solveTimer();
729 updateTimer_ +=
pmodel.updateTimer();
734 const auto&
pmodel = problem_->model();
735 prePostProcessTimer_ +=
pmodel.prePostProcessTimer();
736 linearizeTimer_ +=
pmodel.linearizeTimer();
737 solveTimer_ +=
pmodel.solveTimer();
738 updateTimer_ +=
pmodel.updateTimer();
741 prePostProcessTimer_.
start();
743 OPM_BEGIN_PARALLEL_TRY_CATCH();
744 problem_->endTimeStep();
745 OPM_END_PARALLEL_TRY_CATCH(
"End timestep failed: ",
746 Dune::MPIHelper::getCommunication());
748 prePostProcessTimer_.
stop();
752 if (problem_->shouldWriteOutput()) {
753 OPM_BEGIN_PARALLEL_TRY_CATCH();
754 problem_->writeOutput(
true);
755 OPM_END_PARALLEL_TRY_CATCH(
"Write output failed: ",
756 Dune::MPIHelper::getCommunication());
763 OPM_BEGIN_PARALLEL_TRY_CATCH();
764 problem_->advanceTimeLevel();
765 OPM_END_PARALLEL_TRY_CATCH(
"Advance time level failed: ",
766 Dune::MPIHelper::getCommunication());
770 std::cout <<
"Time step " <<
timeStepIndex() + 1 <<
" done. "
773 <<
", end time: " << this->
time() + oldDt <<
" seconds"
776 <<
"\n" << std::flush;
783 prePostProcessTimer_.
start();
787 OPM_BEGIN_PARALLEL_TRY_CATCH();
788 problem_->endEpisode();
789 OPM_END_PARALLEL_TRY_CATCH(
"End episode failed: ",
790 Dune::MPIHelper::getCommunication());
795 if (timeStepIdx_ <
static_cast<int>(forcedTimeSteps_.size())) {
797 dt = forcedTimeSteps_[timeStepIdx_];
806 prePostProcessTimer_.
stop();
810 if (problem_->shouldWriteRestartFile()) {
811 OPM_BEGIN_PARALLEL_TRY_CATCH();
813 OPM_END_PARALLEL_TRY_CATCH(
"Serialize failed: ",
814 Dune::MPIHelper::getCommunication());
818 executionTimer_.
stop();
821 OPM_BEGIN_PARALLEL_TRY_CATCH();
822 problem_->finalize();
823 OPM_END_PARALLEL_TRY_CATCH(
"Finalize failed: ",
824 Dune::MPIHelper::getCommunication());
828#ifdef RESERVOIR_COUPLING_ENABLED
866 if (
gridView().comm().rank() == 0) {
867 std::cout <<
"Serialize to file '" <<
res.fileName() <<
"'"
869 <<
"\n" << std::flush;
873 problem_->serialize(
res);
874 model_->serialize(
res);
885 template <
class Restarter>
888 restarter.serializeSectionBegin(
"Simulator");
890 << episodeIdx_ <<
" "
891 << episodeStartTime_ <<
" "
892 << episodeLength_ <<
" "
895 << timeStepIdx_ <<
" ";
906 template <
class Restarter>
909 restarter.deserializeSectionBegin(
"Simulator");
920 template<
class Serializer>
935 std::unique_ptr<Vanguard> vanguard_;
936 std::unique_ptr<Model> model_;
937 std::unique_ptr<Problem> problem_;
940 Scalar episodeStartTime_;
941 Scalar episodeLength_;
944 Timer executionTimer_;
945 Timer prePostProcessTimer_;
946 Timer linearizeTimer_;
951 std::vector<Scalar> forcedTimeSteps_;
956 Scalar timeStepSize_;
962#ifdef RESERVOIR_COUPLING_ENABLED
969namespace Properties {
970template<
class TypeTag>
971struct Simulator<TypeTag, TTag::NumericModel>
Defines a type tags and some fundamental properties all models.
Definition ReservoirCouplingMaster.hpp:35
Load or save a state of a problem to/from the harddisk.
Definition restart.hpp:45
void serializeBegin(Simulator &simulator)
Write the current state of the model to disk.
Definition restart.hpp:92
void deserializeBegin(Simulator &simulator, Scalar t)
Start reading a restart file at a certain simulated time.
Definition restart.hpp:147
Manages the initializing and running of time dependent problems.
Definition simulator.hh:84
const Timer & writeTimer() const
Returns a reference to the timer object which measures the time needed to write the visualization out...
Definition simulator.hh:384
Scalar timeStepSize() const
Returns the time step length so that we don't miss the beginning of the next episode or cross the en...
Definition simulator.hh:413
Scalar startTime() const
Return the time of the start of the simulation.
Definition simulator.hh:286
int timeStepIndex() const
Returns number of time steps which have been executed since the beginning of the simulation.
Definition simulator.hh:420
void serialize()
This method writes the complete state of the simulation to the harddisk.
Definition simulator.hh:861
Scalar episodeLength() const
Returns the length of the current episode in simulated time .
Definition simulator.hh:530
const Timer & prePostProcessTimer() const
Returns a reference to the timer object which measures the time needed for pre- and postprocessing of...
Definition simulator.hh:356
const Timer & solveTimer() const
Returns a reference to the timer object which measures the time needed by the solver.
Definition simulator.hh:370
void startNextEpisode(Scalar len=std::numeric_limits< Scalar >::max())
Start the next episode, but don't change the episode identifier.
Definition simulator.hh:488
const Vanguard & vanguard() const
Return a reference to the grid manager of simulation.
Definition simulator.hh:240
void serialize(Restarter &restarter)
Write the time manager's state to a restart file.
Definition simulator.hh:886
const Timer & updateTimer() const
Returns a reference to the timer object which measures the time needed to the solutions of the non-li...
Definition simulator.hh:377
const Timer & executionTimer() const
Returns a reference to the timer object which measures the time needed to run the simulation.
Definition simulator.hh:346
void startNextEpisode(Scalar episodeStartTime, Scalar episodeLength)
Change the current episode of the simulation.
Definition simulator.hh:474
void setEndTime(Scalar t)
Set the time of simulated seconds at which the simulation runs.
Definition simulator.hh:325
const Timer & linearizeTimer() const
Returns a reference to the timer object which measures the time needed for linarizing the solutions.
Definition simulator.hh:363
void setStartTime(Scalar t)
Set the time of the start of the simulation.
Definition simulator.hh:280
void deserialize(Restarter &restarter)
Read the time manager's state from a restart file.
Definition simulator.hh:907
void setTimeStepSize(Scalar value)
Set the current time step size to a given value.
Definition simulator.hh:397
bool episodeStarts() const
Returns true if the current episode has just been started at the current time.
Definition simulator.hh:537
void run()
Runs the simulation using a given problem class.
Definition simulator.hh:592
Vanguard & vanguard()
Return a reference to the grid manager of simulation.
Definition simulator.hh:234
int episodeIndex() const
Returns the index of the current episode.
Definition simulator.hh:508
void setTimeStepIndex(unsigned value)
Set the current time step index to a given value.
Definition simulator.hh:405
void setFinished(bool yesno=true)
Specify whether the simulation is finished.
Definition simulator.hh:430
Problem & problem()
Return the object which specifies the pysical setup of the simulation.
Definition simulator.hh:265
static void registerParameters()
Registers all runtime parameters used by the simulation.
Definition simulator.hh:214
void setTime(Scalar t)
Set the current simulated time, don't change the current time step index.
Definition simulator.hh:295
bool willBeFinished() const
Returns true if the simulation is finished after the time level is incremented by the current time st...
Definition simulator.hh:449
bool finished() const
Returns true if the simulation is finished.
Definition simulator.hh:439
Scalar maxTimeStepSize() const
Aligns the time step size to the episode boundary and to the end time of the simulation.
Definition simulator.hh:458
const Model & model() const
Return the physical model used in the simulation.
Definition simulator.hh:258
void setTime(Scalar t, unsigned stepIdx)
Set the current simulated time and the time step index.
Definition simulator.hh:304
bool episodeIsOver() const
Returns true if the current episode is finished at the current time.
Definition simulator.hh:546
Scalar endTime() const
Returns the number of (simulated) seconds which the simulation runs.
Definition simulator.hh:332
bool episodeWillBeOver() const
Returns true if the current episode will be finished after the current time step.
Definition simulator.hh:555
const GridView & gridView() const
Return the grid view for which the simulation is done.
Definition simulator.hh:246
void setEpisodeIndex(int episodeIdx)
Sets the index of the current episode.
Definition simulator.hh:500
Scalar time() const
Return the number of seconds of simulated time which have elapsed since the start time.
Definition simulator.hh:317
void setEpisodeLength(Scalar dt)
Sets the length in seconds of the current episode.
Definition simulator.hh:523
const Problem & problem() const
Return the object which specifies the pysical setup of the simulation.
Definition simulator.hh:272
Model & model()
Return the physical model used in the simulation.
Definition simulator.hh:252
Scalar episodeMaxTimeStepSize() const
Aligns the time step size to the episode boundary if the current time step crosses the boundary of th...
Definition simulator.hh:565
Scalar episodeStartTime() const
Returns the absolute time when the current episode started .
Definition simulator.hh:515
const Timer & setupTimer() const
Returns a reference to the timer object which measures the time needed to set up and initialize the s...
Definition simulator.hh:339
A simple class which makes sure that a timer gets stopped if an exception is thrown.
Definition timerguard.hh:42
Provides an encapsulation to measure the system time.
Definition timer.hpp:46
void start()
Start counting the time resources used by the simulation.
Definition timer.cpp:46
double realTimeElapsed() const
Return the real time [s] elapsed during the periods the timer was active since the last reset.
Definition timer.cpp:90
double stop()
Stop counting the time resources.
Definition timer.cpp:52
Declare the properties used by the infrastructure code of the finite volume discretizations.
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition blackoilbioeffectsmodules.hh:43
std::string humanReadableTime(double timeInSeconds, bool isAmendment)
Given a time step size in seconds, return it in a format which is more easily parsable by humans.
Definition simulatorutils.cpp:45
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:240
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition propertysystem.hh:233
This file provides the infrastructure to retrieve run-time parameters.
The Opm property system, traits with inheritance.
Load or save a state of a problem to/from the harddisk.
Provides an encapsulation to measure the system time.
A simple class which makes sure that a timer gets stopped if an exception is thrown.