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>
37#include <opm/simulators/flow/NewtonIterationContext.hpp>
42#include <opm/models/utils/simulatorutils.hpp>
51namespace Opm::Properties {
53template <
class TypeTag,
class MyTypeTag>
69template<
class TypeTag>
94 dim = GridView::dimension,
95 dimWorld = GridView::dimensionworld
98 using Element =
typename GridView::template Codim<0>::Entity;
99 using Vertex =
typename GridView::template Codim<dim>::Entity;
100 using VertexIterator =
typename GridView::template Codim<dim>::Iterator;
102 using CoordScalar =
typename GridView::Grid::ctype;
103 using GlobalPosition = Dune::FieldVector<CoordScalar, dimWorld>;
111 FvBaseProblem(
const FvBaseProblem&) =
delete;
122 : nextTimeStepSize_(0.0)
124 , elementMapper_(gridView_, Dune::mcmgElementLayout())
125 , vertexMapper_(gridView_, Dune::mcmgVertexLayout())
126 , boundingBoxMin_(std::numeric_limits<double>::max())
127 , boundingBoxMax_(-std::numeric_limits<double>::max())
131 for (
const auto& vertex : vertices(gridView_)) {
132 for (
unsigned i = 0; i < dim; ++i) {
133 boundingBoxMin_[i] = std::min(boundingBoxMin_[i], vertex.geometry().corner(0)[i]);
134 boundingBoxMax_[i] = std::max(boundingBoxMax_[i], vertex.geometry().corner(0)[i]);
139 for (
unsigned i = 0; i < dim; ++i) {
140 boundingBoxMin_[i] = gridView_.comm().min(boundingBoxMin_[i]);
141 boundingBoxMax_[i] = gridView_.comm().max(boundingBoxMax_[i]);
144 if (enableVtkOutput_()) {
145 const bool asyncVtkOutput =
146 simulator_.gridView().comm().size() == 1 &&
153 if (asyncVtkOutput && enableGridAdaptation) {
154 throw std::runtime_error(
"Asynchronous VTK output currently cannot be used "
155 "at the same time as grid adaptivity");
158 defaultVtkWriter_ = std::make_unique<VtkMultiWriter>(asyncVtkOutput,
171 Model::registerParameters();
173 (
"The maximum size to which all time steps are limited to [s]");
175 (
"The minimum size to which all time steps are limited to [s]");
177 (
"The maximum number of divisions by two of the timestep size "
178 "before the simulation bails out");
180 (
"Dispatch a separate thread to write the VTK output");
182 (
"Continue with a non-converged solution instead of giving up "
183 "if we encounter a time step size smaller than the minimum time "
233 std::string desc = Implementation::briefDescription();
238 return "Usage: " + std::string(argv[0]) +
" [OPTIONS]\n" + desc;
269 const std::string&)>,
270 std::set<std::string>&,
271 std::string& errorMsg,
277 errorMsg = std::string(
"Illegal parameter \"") + argv[paramIdx] +
"\".";
303 elementMapper_.update(gridView_);
304 vertexMapper_.update(gridView_);
306 if (enableVtkOutput_()) {
307 defaultVtkWriter_->gridChanged();
320 template <
class Context>
325 {
throw std::logic_error(
"Problem does not provide a boundary() method"); }
337 template <
class Context>
342 {
throw std::logic_error(
"Problem does not provide a constraints() method"); }
356 template <
class Context>
361 {
throw std::logic_error(
"Problem does not provide a source() method"); }
373 template <
class Context>
378 {
throw std::logic_error(
"Problem does not provide a initial() method"); }
395 template <
class Context>
399 {
return asImp_().extrusionFactor(); }
451 std::cerr <<
"The end of episode " <<
simulator().episodeIndex() + 1 <<
" has been "
452 <<
"reached, but the problem does not override the endEpisode() method. "
453 <<
"Doing nothing!\n";
461 const auto& executionTimer =
simulator().executionTimer();
463 const Scalar executionTime = executionTimer.realTimeElapsed();
464 const Scalar setupTime =
simulator().setupTimer().realTimeElapsed();
465 const Scalar prePostProcessTime =
simulator().prePostProcessTimer().realTimeElapsed();
466 const Scalar localCpuTime = executionTimer.cpuTimeElapsed();
467 const Scalar globalCpuTime = executionTimer.globalCpuTimeElapsed();
468 const Scalar writeTime =
simulator().writeTimer().realTimeElapsed();
469 const Scalar linearizeTime =
simulator().linearizeTimer().realTimeElapsed();
470 const Scalar solveTime =
simulator().solveTimer().realTimeElapsed();
471 const Scalar updateTime =
simulator().updateTimer().realTimeElapsed();
472 const unsigned numProcesses =
static_cast<unsigned>(this->
gridView().comm().size());
474 if (
gridView().comm().rank() == 0) {
475 std::cout << std::setprecision(3)
476 <<
"Simulation of problem '" << asImp_().name() <<
"' finished.\n"
478 <<
"------------------------ Timing ------------------------\n"
479 <<
"Setup time: " << setupTime <<
" seconds"
481 <<
", " << setupTime / (executionTime + setupTime) * 100 <<
"%\n"
482 <<
"Simulation time: " << executionTime <<
" seconds"
484 <<
", " << executionTime / (executionTime + setupTime) * 100 <<
"%\n"
485 <<
" Linearization time: " << linearizeTime <<
" seconds"
487 <<
", " << linearizeTime / executionTime * 100 <<
"%\n"
488 <<
" Linear solve time: " << solveTime <<
" seconds"
490 <<
", " << solveTime / executionTime * 100 <<
"%\n"
491 <<
" Newton update time: " << updateTime <<
" seconds"
493 <<
", " << updateTime / executionTime * 100 <<
"%\n"
494 <<
" Pre/postprocess time: " << prePostProcessTime <<
" seconds"
496 <<
", " << prePostProcessTime / executionTime * 100 <<
"%\n"
497 <<
" Output write time: " << writeTime <<
" seconds"
499 <<
", " << writeTime / executionTime * 100 <<
"%\n"
500 <<
"First process' simulation CPU time: " << localCpuTime <<
" seconds"
502 <<
"Number of processes: " << numProcesses <<
"\n"
503 <<
"Threads per processes: " << threadsPerProcess <<
"\n"
504 <<
"Total CPU time: " << globalCpuTime <<
" seconds"
507 <<
"----------------------------------------------------------------\n"
518 const unsigned maxFails = asImp_().maxTimeIntegrationFailures();
519 Scalar minTimeStep = asImp_().minTimeStepSize();
521 std::string errorMessage;
522 for (
unsigned i = 0; i < maxFails; ++i) {
523 bool converged =
model().update();
528 const Scalar dt =
simulator().timeStepSize();
529 Scalar nextDt = dt / 2.0;
530 if (dt < minTimeStep * (1.0 + 1e-9)) {
532 if (
gridView().comm().rank() == 0) {
533 std::cout <<
"Newton solver did not converge with minimum time step of "
534 << dt <<
" seconds. Continuing with unconverged solution!\n"
540 errorMessage =
"Time integration did not succeed with the minumum time step size of " +
541 std::to_string(
double(minTimeStep)) +
" seconds. Giving up!";
545 else if (nextDt < minTimeStep) {
546 nextDt = minTimeStep;
551 if (
gridView().comm().rank() == 0) {
552 std::cout <<
"Newton solver did not converge with "
553 <<
"dt=" << dt <<
" seconds. Retrying with time step of "
554 << nextDt <<
" seconds\n" << std::flush;
558 if (errorMessage.empty()) {
559 errorMessage =
"Newton solver didn't converge after " +
560 std::to_string(maxFails) +
" time-step divisions. dt=" +
561 std::to_string(
double(
simulator().timeStepSize()));
563 throw std::runtime_error(errorMessage);
591 { nextTimeStepSize_ = dt; }
600 if (nextTimeStepSize_ > 0.0) {
601 return nextTimeStepSize_;
607 if (dtNext <
simulator().maxTimeStepSize() &&
608 simulator().maxTimeStepSize() < dtNext * 2)
610 dtNext =
simulator().maxTimeStepSize() / 2 * 1.01;
626 return simulator().timeStepIndex() > 0 &&
647 {
model().advanceTimeLevel(); }
663 {
return gridView_; }
670 {
return boundingBoxMin_; }
677 {
return boundingBoxMax_; }
683 {
return vertexMapper_; }
689 {
return elementMapper_; }
695 {
return simulator_; }
701 {
return simulator_; }
707 {
return simulator_.model(); }
713 {
return simulator_.model(); }
719 {
return model().newtonMethod(); }
725 {
return model().newtonMethod(); }
734 return RestrictProlongOperator();
762 template <
class Restarter>
765 if (enableVtkOutput_()) {
766 defaultVtkWriter_->serialize(res);
780 template <
class Restarter>
783 if (enableVtkOutput_()) {
784 defaultVtkWriter_->deserialize(res);
796 if (!enableVtkOutput_()) {
800 if (verbose &&
gridView().comm().rank() == 0) {
801 std::cout <<
"Writing visualization results for the current time step.\n"
808 defaultVtkWriter_->beginWrite(t);
809 model().prepareOutputFields();
810 model().appendOutputFields(*defaultVtkWriter_);
811 defaultVtkWriter_->endWrite(
false);
819 {
return *defaultVtkWriter_; }
825 {
return iterationContext_; }
831 { iterationContext_.resetForNewTimestep(); }
837 { iterationContext_.advanceIteration(); }
843 { iterationContext_.markTimestepInitialized(); }
854 {
return iterationContext_; }
859 Scalar nextTimeStepSize_;
860 NewtonIterationContext iterationContext_;
862 bool enableVtkOutput_()
const
863 {
return Parameters::Get<Parameters::EnableVtkOutput>(); }
867 Implementation& asImp_()
868 {
return *
static_cast<Implementation*
>(
this); }
871 const Implementation& asImp_()
const
872 {
return *
static_cast<const Implementation*
>(
this); }
875 const GridView gridView_;
876 ElementMapper elementMapper_;
877 VertexMapper vertexMapper_;
878 GlobalPosition boundingBoxMin_;
879 GlobalPosition boundingBoxMax_;
882 Simulator& simulator_;
883 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:268
std::string name() const
The problem name.
Definition fvbaseproblem.hh:656
const Simulator & simulator() const
Returns Simulator object used by the simulation.
Definition fvbaseproblem.hh:700
void writeOutput(bool verbose=true)
Write the relevant secondary variables of the current solution into an VTK output file.
Definition fvbaseproblem.hh:794
void boundary(BoundaryRateVector &, const Context &, unsigned, unsigned) const
Evaluate the boundary conditions for a boundary segment.
Definition fvbaseproblem.hh:321
void finalize()
Called after the simulation has been run sucessfully.
Definition fvbaseproblem.hh:459
Simulator & simulator()
Returns Simulator object used by the simulation.
Definition fvbaseproblem.hh:694
const Model & model() const
Returns numerical model used for the problem.
Definition fvbaseproblem.hh:712
static void registerParameters()
Registers all available parameters for the problem and the model.
Definition fvbaseproblem.hh:169
void source(RateVector &, const Context &, unsigned, unsigned) const
Evaluate the source term for all phases within a given sub-control-volume.
Definition fvbaseproblem.hh:357
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:194
unsigned maxTimeIntegrationFailures() const
Returns the maximum number of subsequent failures for the time integration before giving up.
Definition fvbaseproblem.hh:576
unsigned markForGridAdaptation()
Mark grid cells for refinement or coarsening.
Definition fvbaseproblem.hh:744
void advanceIteration()
Advance the iteration counter.
Definition fvbaseproblem.hh:836
void serialize(Restarter &res)
This method writes the complete state of the problem to the harddisk.
Definition fvbaseproblem.hh:763
void beginIteration()
Called by the simulator before each Newton-Raphson iteration.
Definition fvbaseproblem.hh:426
void markTimestepInitialized()
Mark timestep initialization as complete.
Definition fvbaseproblem.hh:842
unsigned intensiveQuantityHistorySize() const
Returns the required history size for intensive quantities cache.
Definition fvbaseproblem.hh:204
void endTimeStep()
Called by the simulator after each time integration.
Definition fvbaseproblem.hh:441
FvBaseProblem(Simulator &simulator)
Definition fvbaseproblem.hh:121
const GlobalPosition & boundingBoxMax() const
The coordinate of the corner of the GridView's bounding box with the largest values.
Definition fvbaseproblem.hh:676
RestrictProlongOperator restrictProlongOperator()
return restriction and prolongation operator
Definition fvbaseproblem.hh:732
const NewtonIterationContext & iterationContext() const
Returns the iteration context for iteration-dependent decisions.
Definition fvbaseproblem.hh:824
void timeIntegration()
Called by Opm::Simulator in order to do a time integration on the model.
Definition fvbaseproblem.hh:516
Model & model()
Returns numerical model used for the problem.
Definition fvbaseproblem.hh:706
void endEpisode()
Called when the end of an simulation episode is reached.
Definition fvbaseproblem.hh:449
void beginTimeStep()
Called by the simulator before each time integration.
Definition fvbaseproblem.hh:420
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:231
void deserialize(Restarter &res)
This method restores the complete state of the problem from disk.
Definition fvbaseproblem.hh:781
void prefetch(const Element &) const
Allows to improve the performance by prefetching all data which is associated with a given element.
Definition fvbaseproblem.hh:293
Scalar extrusionFactor(const Context &, unsigned, unsigned) const
Return how much the domain is extruded at a given sub-control volume.
Definition fvbaseproblem.hh:396
bool shouldWriteRestartFile() const
Returns true if a restart file should be written to disk.
Definition fvbaseproblem.hh:624
Scalar nextTimeStepSize() const
Called by Opm::Simulator whenever a solution for a time step has been computed and the simulation tim...
Definition fvbaseproblem.hh:598
Scalar minTimeStepSize() const
Returns the minimum allowable size of a time step.
Definition fvbaseproblem.hh:569
void resetIterationForNewTimestep()
Reset the iteration context for a new timestep.
Definition fvbaseproblem.hh:830
void beginEpisode()
Called at the beginning of an simulation episode.
Definition fvbaseproblem.hh:414
void initial(PrimaryVariables &, const Context &, unsigned, unsigned) const
Evaluate the initial value for a control volume.
Definition fvbaseproblem.hh:374
std::string outputDir() const
Determine the directory for simulation output.
Definition fvbaseproblem.hh:220
const VertexMapper & vertexMapper() const
Returns the mapper for vertices to indices.
Definition fvbaseproblem.hh:682
void gridChanged()
Handle changes of the grid.
Definition fvbaseproblem.hh:301
const GlobalPosition & boundingBoxMin() const
The coordinate of the corner of the GridView's bounding box with the smallest values.
Definition fvbaseproblem.hh:669
static std::string briefDescription()
Returns a human readable description of the problem for the help message.
Definition fvbaseproblem.hh:247
void endIteration()
Called by the simulator after each Newton-Raphson update.
Definition fvbaseproblem.hh:432
void setNextTimeStepSize(Scalar dt)
Impose the next time step size to be used externally.
Definition fvbaseproblem.hh:590
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:584
void constraints(Constraints &, const Context &, unsigned, unsigned) const
Evaluate the constraints for a control volume.
Definition fvbaseproblem.hh:338
void finishInit()
Called by the Opm::Simulator in order to initialize the problem.
Definition fvbaseproblem.hh:286
const NewtonMethod & newtonMethod() const
Returns object which implements the Newton method.
Definition fvbaseproblem.hh:724
const GridView & gridView() const
The GridView which used by the problem.
Definition fvbaseproblem.hh:662
NewtonMethod & newtonMethod()
Returns object which implements the Newton method.
Definition fvbaseproblem.hh:718
bool shouldWriteOutput() const
Returns true if the current solution should be written to disk (i.e.
Definition fvbaseproblem.hh:638
void initialSolutionApplied()
Callback used by the model to indicate that the initial solution has been determined for all degrees ...
Definition fvbaseproblem.hh:408
void advanceTimeLevel()
Called by the simulator after everything which can be done about the current time step is finished an...
Definition fvbaseproblem.hh:646
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:818
const ElementMapper & elementMapper() const
Returns the mapper for elements to indices.
Definition fvbaseproblem.hh:688
RAII guard for NLDD domain-local iteration context.
Definition NewtonIterationContext.hpp:152
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:45
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.
Context for iteration-dependent decisions in the Newton solver.
Definition NewtonIterationContext.hpp:43
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.