28#ifndef EWOMS_FV_BASE_PROBLEM_HH
29#define EWOMS_FV_BASE_PROBLEM_HH
31#include <dune/common/fvector.hh>
35#include <opm/models/discretization/common/restrictprolong.hh>
40#include <opm/models/utils/simulatorutils.hpp>
49namespace Opm::Properties {
51template <
class TypeTag,
class MyTypeTag>
67template<
class TypeTag>
92 dim = GridView::dimension,
93 dimWorld = GridView::dimensionworld
96 using Element =
typename GridView::template Codim<0>::Entity;
97 using Vertex =
typename GridView::template Codim<dim>::Entity;
98 using VertexIterator =
typename GridView::template Codim<dim>::Iterator;
100 using CoordScalar =
typename GridView::Grid::ctype;
101 using GlobalPosition = Dune::FieldVector<CoordScalar, dimWorld>;
109 FvBaseProblem(
const FvBaseProblem&) =
delete;
120 : nextTimeStepSize_(0.0)
122 , elementMapper_(gridView_, Dune::mcmgElementLayout())
123 , vertexMapper_(gridView_, Dune::mcmgVertexLayout())
124 , boundingBoxMin_(std::numeric_limits<double>::max())
125 , boundingBoxMax_(-std::numeric_limits<double>::max())
129 for (
const auto& vertex : vertices(gridView_)) {
130 for (
unsigned i = 0; i < dim; ++i) {
131 boundingBoxMin_[i] = std::min(boundingBoxMin_[i], vertex.geometry().corner(0)[i]);
132 boundingBoxMax_[i] = std::max(boundingBoxMax_[i], vertex.geometry().corner(0)[i]);
137 for (
unsigned i = 0; i < dim; ++i) {
138 boundingBoxMin_[i] = gridView_.comm().min(boundingBoxMin_[i]);
139 boundingBoxMax_[i] = gridView_.comm().max(boundingBoxMax_[i]);
142 if (enableVtkOutput_()) {
143 const bool asyncVtkOutput =
144 simulator_.gridView().comm().size() == 1 &&
151 if (asyncVtkOutput && enableGridAdaptation) {
152 throw std::runtime_error(
"Asynchronous VTK output currently cannot be used "
153 "at the same time as grid adaptivity");
156 defaultVtkWriter_ = std::make_unique<VtkMultiWriter>(asyncVtkOutput,
169 Model::registerParameters();
171 (
"The maximum size to which all time steps are limited to [s]");
173 (
"The minimum size to which all time steps are limited to [s]");
175 (
"The maximum number of divisions by two of the timestep size "
176 "before the simulation bails out");
178 (
"Dispatch a separate thread to write the VTK output");
180 (
"Continue with a non-converged solution instead of giving up "
181 "if we encounter a time step size smaller than the minimum time "
231 std::string desc = Implementation::briefDescription();
236 return "Usage: " + std::string(argv[0]) +
" [OPTIONS]\n" + desc;
267 const std::string&)>,
268 std::set<std::string>&,
269 std::string& errorMsg,
275 errorMsg = std::string(
"Illegal parameter \"") + argv[paramIdx] +
"\".";
301 elementMapper_.update(gridView_);
302 vertexMapper_.update(gridView_);
304 if (enableVtkOutput_()) {
305 defaultVtkWriter_->gridChanged();
318 template <
class Context>
323 {
throw std::logic_error(
"Problem does not provide a boundary() method"); }
335 template <
class Context>
340 {
throw std::logic_error(
"Problem does not provide a constraints() method"); }
354 template <
class Context>
359 {
throw std::logic_error(
"Problem does not provide a source() method"); }
371 template <
class Context>
376 {
throw std::logic_error(
"Problem does not provide a initial() method"); }
393 template <
class Context>
397 {
return asImp_().extrusionFactor(); }
449 std::cerr <<
"The end of episode " <<
simulator().episodeIndex() + 1 <<
" has been "
450 <<
"reached, but the problem does not override the endEpisode() method. "
451 <<
"Doing nothing!\n";
459 const auto& executionTimer =
simulator().executionTimer();
461 const Scalar executionTime = executionTimer.realTimeElapsed();
462 const Scalar setupTime =
simulator().setupTimer().realTimeElapsed();
463 const Scalar prePostProcessTime =
simulator().prePostProcessTimer().realTimeElapsed();
464 const Scalar localCpuTime = executionTimer.cpuTimeElapsed();
465 const Scalar globalCpuTime = executionTimer.globalCpuTimeElapsed();
466 const Scalar writeTime =
simulator().writeTimer().realTimeElapsed();
467 const Scalar linearizeTime =
simulator().linearizeTimer().realTimeElapsed();
468 const Scalar solveTime =
simulator().solveTimer().realTimeElapsed();
469 const Scalar updateTime =
simulator().updateTimer().realTimeElapsed();
470 const unsigned numProcesses =
static_cast<unsigned>(this->
gridView().comm().size());
472 if (
gridView().comm().rank() == 0) {
473 std::cout << std::setprecision(3)
474 <<
"Simulation of problem '" << asImp_().name() <<
"' finished.\n"
476 <<
"------------------------ Timing ------------------------\n"
477 <<
"Setup time: " << setupTime <<
" seconds"
479 <<
", " << setupTime / (executionTime + setupTime) * 100 <<
"%\n"
480 <<
"Simulation time: " << executionTime <<
" seconds"
482 <<
", " << executionTime / (executionTime + setupTime) * 100 <<
"%\n"
483 <<
" Linearization time: " << linearizeTime <<
" seconds"
485 <<
", " << linearizeTime / executionTime * 100 <<
"%\n"
486 <<
" Linear solve time: " << solveTime <<
" seconds"
488 <<
", " << solveTime / executionTime * 100 <<
"%\n"
489 <<
" Newton update time: " << updateTime <<
" seconds"
491 <<
", " << updateTime / executionTime * 100 <<
"%\n"
492 <<
" Pre/postprocess time: " << prePostProcessTime <<
" seconds"
494 <<
", " << prePostProcessTime / executionTime * 100 <<
"%\n"
495 <<
" Output write time: " << writeTime <<
" seconds"
497 <<
", " << writeTime / executionTime * 100 <<
"%\n"
498 <<
"First process' simulation CPU time: " << localCpuTime <<
" seconds"
500 <<
"Number of processes: " << numProcesses <<
"\n"
501 <<
"Threads per processes: " << threadsPerProcess <<
"\n"
502 <<
"Total CPU time: " << globalCpuTime <<
" seconds"
505 <<
"----------------------------------------------------------------\n"
516 const unsigned maxFails = asImp_().maxTimeIntegrationFailures();
517 Scalar minTimeStep = asImp_().minTimeStepSize();
519 std::string errorMessage;
520 for (
unsigned i = 0; i < maxFails; ++i) {
521 bool converged =
model().update();
526 const Scalar dt =
simulator().timeStepSize();
527 Scalar nextDt = dt / 2.0;
528 if (dt < minTimeStep * (1.0 + 1e-9)) {
530 if (
gridView().comm().rank() == 0) {
531 std::cout <<
"Newton solver did not converge with minimum time step of "
532 << dt <<
" seconds. Continuing with unconverged solution!\n"
538 errorMessage =
"Time integration did not succeed with the minumum time step size of " +
539 std::to_string(
double(minTimeStep)) +
" seconds. Giving up!";
543 else if (nextDt < minTimeStep) {
544 nextDt = minTimeStep;
549 if (
gridView().comm().rank() == 0) {
550 std::cout <<
"Newton solver did not converge with "
551 <<
"dt=" << dt <<
" seconds. Retrying with time step of "
552 << nextDt <<
" seconds\n" << std::flush;
556 if (errorMessage.empty()) {
557 errorMessage =
"Newton solver didn't converge after " +
558 std::to_string(maxFails) +
" time-step divisions. dt=" +
559 std::to_string(
double(
simulator().timeStepSize()));
561 throw std::runtime_error(errorMessage);
589 { nextTimeStepSize_ = dt; }
598 if (nextTimeStepSize_ > 0.0) {
599 return nextTimeStepSize_;
605 if (dtNext <
simulator().maxTimeStepSize() &&
606 simulator().maxTimeStepSize() < dtNext * 2)
608 dtNext =
simulator().maxTimeStepSize() / 2 * 1.01;
624 return simulator().timeStepIndex() > 0 &&
645 {
model().advanceTimeLevel(); }
661 {
return gridView_; }
668 {
return boundingBoxMin_; }
675 {
return boundingBoxMax_; }
681 {
return vertexMapper_; }
687 {
return elementMapper_; }
693 {
return simulator_; }
699 {
return simulator_; }
705 {
return simulator_.model(); }
711 {
return simulator_.model(); }
717 {
return model().newtonMethod(); }
723 {
return model().newtonMethod(); }
732 return RestrictProlongOperator();
760 template <
class Restarter>
763 if (enableVtkOutput_()) {
764 defaultVtkWriter_->serialize(res);
778 template <
class Restarter>
781 if (enableVtkOutput_()) {
782 defaultVtkWriter_->deserialize(res);
794 if (!enableVtkOutput_()) {
798 if (verbose &&
gridView().comm().rank() == 0) {
799 std::cout <<
"Writing visualization results for the current time step.\n"
806 defaultVtkWriter_->beginWrite(t);
807 model().prepareOutputFields();
808 model().appendOutputFields(*defaultVtkWriter_);
809 defaultVtkWriter_->endWrite(
false);
817 {
return *defaultVtkWriter_; }
820 Scalar nextTimeStepSize_;
822 bool enableVtkOutput_()
const
827 Implementation& asImp_()
828 {
return *
static_cast<Implementation*
>(
this); }
831 const Implementation& asImp_()
const
832 {
return *
static_cast<const Implementation*
>(
this); }
835 const GridView gridView_;
836 ElementMapper elementMapper_;
837 VertexMapper vertexMapper_;
838 GlobalPosition boundingBoxMin_;
839 GlobalPosition boundingBoxMax_;
842 Simulator& simulator_;
843 std::unique_ptr<VtkMultiWriter> defaultVtkWriter_{};
Definition restrictprolong.hh:142
static int handlePositionalParameter(std::function< void(const std::string &, const std::string &)>, std::set< std::string > &, std::string &errorMsg, int, const char **argv, int paramIdx, int)
Handles positional command line parameters.
Definition fvbaseproblem.hh:266
std::string name() const
The problem name.
Definition fvbaseproblem.hh:654
const Simulator & simulator() const
Returns Simulator object used by the simulation.
Definition fvbaseproblem.hh:698
void writeOutput(bool verbose=true)
Write the relevant secondary variables of the current solution into an VTK output file.
Definition fvbaseproblem.hh:792
void boundary(BoundaryRateVector &, const Context &, unsigned, unsigned) const
Evaluate the boundary conditions for a boundary segment.
Definition fvbaseproblem.hh:319
void finalize()
Called after the simulation has been run sucessfully.
Definition fvbaseproblem.hh:457
Simulator & simulator()
Returns Simulator object used by the simulation.
Definition fvbaseproblem.hh:692
const Model & model() const
Returns numerical model used for the problem.
Definition fvbaseproblem.hh:710
static void registerParameters()
Registers all available parameters for the problem and the model.
Definition fvbaseproblem.hh:167
void source(RateVector &, const Context &, unsigned, unsigned) const
Evaluate the source term for all phases within a given sub-control-volume.
Definition fvbaseproblem.hh:355
bool recycleFirstIterationStorage() const
Return if the storage term of the first iteration is identical to the storage term for the solution o...
Definition fvbaseproblem.hh:192
unsigned maxTimeIntegrationFailures() const
Returns the maximum number of subsequent failures for the time integration before giving up.
Definition fvbaseproblem.hh:574
unsigned markForGridAdaptation()
Mark grid cells for refinement or coarsening.
Definition fvbaseproblem.hh:742
void serialize(Restarter &res)
This method writes the complete state of the problem to the harddisk.
Definition fvbaseproblem.hh:761
void beginIteration()
Called by the simulator before each Newton-Raphson iteration.
Definition fvbaseproblem.hh:424
unsigned intensiveQuantityHistorySize() const
Returns the required history size for intensive quantities cache.
Definition fvbaseproblem.hh:202
void endTimeStep()
Called by the simulator after each time integration.
Definition fvbaseproblem.hh:439
FvBaseProblem(Simulator &simulator)
Definition fvbaseproblem.hh:119
const GlobalPosition & boundingBoxMax() const
The coordinate of the corner of the GridView's bounding box with the largest values.
Definition fvbaseproblem.hh:674
RestrictProlongOperator restrictProlongOperator()
return restriction and prolongation operator
Definition fvbaseproblem.hh:730
void timeIntegration()
Called by Opm::Simulator in order to do a time integration on the model.
Definition fvbaseproblem.hh:514
Model & model()
Returns numerical model used for the problem.
Definition fvbaseproblem.hh:704
void endEpisode()
Called when the end of an simulation episode is reached.
Definition fvbaseproblem.hh:447
void beginTimeStep()
Called by the simulator before each time integration.
Definition fvbaseproblem.hh:418
static std::string helpPreamble(int, const char **argv)
Returns the string that is printed before the list of command line parameters in the help message.
Definition fvbaseproblem.hh:229
void deserialize(Restarter &res)
This method restores the complete state of the problem from disk.
Definition fvbaseproblem.hh:779
void prefetch(const Element &) const
Allows to improve the performance by prefetching all data which is associated with a given element.
Definition fvbaseproblem.hh:291
Scalar extrusionFactor(const Context &, unsigned, unsigned) const
Return how much the domain is extruded at a given sub-control volume.
Definition fvbaseproblem.hh:394
bool shouldWriteRestartFile() const
Returns true if a restart file should be written to disk.
Definition fvbaseproblem.hh:622
Scalar nextTimeStepSize() const
Called by Opm::Simulator whenever a solution for a time step has been computed and the simulation tim...
Definition fvbaseproblem.hh:596
Scalar minTimeStepSize() const
Returns the minimum allowable size of a time step.
Definition fvbaseproblem.hh:567
void beginEpisode()
Called at the beginning of an simulation episode.
Definition fvbaseproblem.hh:412
void initial(PrimaryVariables &, const Context &, unsigned, unsigned) const
Evaluate the initial value for a control volume.
Definition fvbaseproblem.hh:372
std::string outputDir() const
Determine the directory for simulation output.
Definition fvbaseproblem.hh:218
const VertexMapper & vertexMapper() const
Returns the mapper for vertices to indices.
Definition fvbaseproblem.hh:680
void gridChanged()
Handle changes of the grid.
Definition fvbaseproblem.hh:299
const GlobalPosition & boundingBoxMin() const
The coordinate of the corner of the GridView's bounding box with the smallest values.
Definition fvbaseproblem.hh:667
static std::string briefDescription()
Returns a human readable description of the problem for the help message.
Definition fvbaseproblem.hh:245
void endIteration()
Called by the simulator after each Newton-Raphson update.
Definition fvbaseproblem.hh:430
void setNextTimeStepSize(Scalar dt)
Impose the next time step size to be used externally.
Definition fvbaseproblem.hh:588
bool continueOnConvergenceError() const
Returns if we should continue with a non-converged solution instead of giving up if we encounter a ti...
Definition fvbaseproblem.hh:582
void constraints(Constraints &, const Context &, unsigned, unsigned) const
Evaluate the constraints for a control volume.
Definition fvbaseproblem.hh:336
void finishInit()
Called by the Opm::Simulator in order to initialize the problem.
Definition fvbaseproblem.hh:284
const NewtonMethod & newtonMethod() const
Returns object which implements the Newton method.
Definition fvbaseproblem.hh:722
const GridView & gridView() const
The GridView which used by the problem.
Definition fvbaseproblem.hh:660
NewtonMethod & newtonMethod()
Returns object which implements the Newton method.
Definition fvbaseproblem.hh:716
bool shouldWriteOutput() const
Returns true if the current solution should be written to disk (i.e.
Definition fvbaseproblem.hh:636
void initialSolutionApplied()
Callback used by the model to indicate that the initial solution has been determined for all degrees ...
Definition fvbaseproblem.hh:406
void advanceTimeLevel()
Called by the simulator after everything which can be done about the current time step is finished an...
Definition fvbaseproblem.hh:644
VtkMultiWriter & defaultVtkWriter() const
Method to retrieve the VTK writer which should be used to write the default ouput after each time ste...
Definition fvbaseproblem.hh:816
const ElementMapper & elementMapper() const
Returns the mapper for elements to indices.
Definition fvbaseproblem.hh:686
static unsigned maxThreads()
Return the maximum number of threads of the current process.
Definition threadmanager.hpp:66
Simplifies writing multi-file VTK datasets.
Definition vtkmultiwriter.hh:65
Declare the properties used by the infrastructure code of the finite volume discretizations.
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
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
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
std::string simulatorOutputDir()
Determine and check the configured directory for simulation output.
Definition simulatorutils.cpp:119
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:240
void Register(const char *usageString)
Register a run-time parameter.
Definition parametersystem.hpp:292
auto Get(bool errorIfNotRegistered=true)
Retrieve a runtime parameter.
Definition parametersystem.hpp:187
Load or save a state of a problem to/from the harddisk.
Specify the maximum size of a time integration [s].
Definition fvbaseparameters.hh:106
Specifies the type of the actual Newton method.
Definition newtonmethodproperties.hh:32
Simplifies writing multi-file VTK datasets.