fvbaseproblem.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*
4 This file is part of the Open Porous Media project (OPM).
5
6 OPM is free software: you can redistribute it and/or modify
7 it under the terms of the GNU General Public License as published by
8 the Free Software Foundation, either version 2 of the License, or
9 (at your option) any later version.
10
11 OPM is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU General Public License for more details.
15
16 You should have received a copy of the GNU General Public License
17 along with OPM. If not, see <http://www.gnu.org/licenses/>.
18
19 Consult the COPYING file in the top-level source directory of this
20 module for the precise wording of the license and the list of
21 copyright holders.
22*/
28#ifndef EWOMS_FV_BASE_PROBLEM_HH
29#define EWOMS_FV_BASE_PROBLEM_HH
30
31#include <dune/common/fvector.hh>
32
36
38
41
43
44#include <functional>
45#include <iomanip>
46#include <iostream>
47#include <limits>
48#include <memory>
49#include <string>
50
51namespace Opm::Properties {
52
53template <class TypeTag, class MyTypeTag>
54struct NewtonMethod;
55
56} // namespace Opm::Properties
57
58namespace Opm {
59
69template<class TypeTag>
71{
72private:
73 using Implementation = GetPropType<TypeTag, Properties::Problem>;
75
76 static constexpr auto vtkOutputFormat = getPropValue<TypeTag, Properties::VtkOutputFormat>();
78
84
87
92
93 enum {
94 dim = GridView::dimension,
95 dimWorld = GridView::dimensionworld
96 };
97
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;
101
102 using CoordScalar = typename GridView::Grid::ctype;
103 using GlobalPosition = Dune::FieldVector<CoordScalar, dimWorld>;
104
105public:
106 // the default restriction and prolongation for adaptation is simply an empty one
108
109private:
110 // copying a problem is not a good idea
111 FvBaseProblem(const FvBaseProblem&) = delete;
112
113public:
121 explicit FvBaseProblem(Simulator& simulator)
122 : nextTimeStepSize_(0.0)
123 , gridView_(simulator.gridView())
124 , elementMapper_(gridView_, Dune::mcmgElementLayout())
125 , vertexMapper_(gridView_, Dune::mcmgVertexLayout())
126 , boundingBoxMin_(std::numeric_limits<double>::max())
127 , boundingBoxMax_(-std::numeric_limits<double>::max())
128 , simulator_(simulator)
129 {
130 // calculate the bounding box of the local partition of the grid view
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]);
135 }
136 }
137
138 // communicate to get the bounding box of the whole domain
139 for (unsigned i = 0; i < dim; ++i) {
140 boundingBoxMin_[i] = gridView_.comm().min(boundingBoxMin_[i]);
141 boundingBoxMax_[i] = gridView_.comm().max(boundingBoxMax_[i]);
142 }
143
144 if (enableVtkOutput_()) {
145 const bool asyncVtkOutput =
146 simulator_.gridView().comm().size() == 1 &&
147 Parameters::Get<Parameters::EnableAsyncVtkOutput>();
148
149 // asynchonous VTK output currently does not work in conjunction with grid
150 // adaptivity because the async-IO code assumes that the grid stays
151 // constant. complain about that case.
152 const bool enableGridAdaptation = Parameters::Get<Parameters::EnableGridAdaptation>();
153 if (asyncVtkOutput && enableGridAdaptation) {
154 throw std::runtime_error("Asynchronous VTK output currently cannot be used "
155 "at the same time as grid adaptivity");
156 }
157
158 defaultVtkWriter_ = std::make_unique<VtkMultiWriter>(asyncVtkOutput,
159 gridView_,
160 asImp_().outputDir(),
161 asImp_().name());
162 }
163 }
164
169 static void registerParameters()
170 {
171 Model::registerParameters();
172 Parameters::Register<Parameters::MaxTimeStepSize<Scalar>>
173 ("The maximum size to which all time steps are limited to [s]");
174 Parameters::Register<Parameters::MinTimeStepSize<Scalar>>
175 ("The minimum size to which all time steps are limited to [s]");
176 Parameters::Register<Parameters::MaxTimeStepDivisions>
177 ("The maximum number of divisions by two of the timestep size "
178 "before the simulation bails out");
179 Parameters::Register<Parameters::EnableAsyncVtkOutput>
180 ("Dispatch a separate thread to write the VTK output");
181 Parameters::Register<Parameters::ContinueOnConvergenceError>
182 ("Continue with a non-converged solution instead of giving up "
183 "if we encounter a time step size smaller than the minimum time "
184 "step size.");
185 }
186
195 { return true; }
196
205 {
206 if (Parameters::Get<Parameters::EnableStorageCache>() && asImp_().recycleFirstIterationStorage()) {
207 return 1;
208 }
209 return 2;
210 }
211
220 std::string outputDir() const
221 {
222 return simulatorOutputDir();
223 }
224
231 static std::string helpPreamble(int, const char **argv)
232 {
233 std::string desc = Implementation::briefDescription();
234 if (!desc.empty()) {
235 desc = desc + "\n";
236 }
237
238 return "Usage: " + std::string(argv[0]) + " [OPTIONS]\n" + desc;
239 }
240
247 static std::string briefDescription()
248 { return ""; }
249
250 // TODO (?): detailedDescription()
251
268 static int handlePositionalParameter(std::function<void(const std::string&,
269 const std::string&)>,
270 std::set<std::string>&,
271 std::string& errorMsg,
272 int,
273 const char** argv,
274 int paramIdx,
275 int)
276 {
277 errorMsg = std::string("Illegal parameter \"") + argv[paramIdx] + "\".";
278 return 0;
279 }
280
287 {}
288
293 void prefetch(const Element&) const
294 {
295 // do nothing by default
296 }
297
302 {
303 elementMapper_.update(gridView_);
304 vertexMapper_.update(gridView_);
305
306 if (enableVtkOutput_()) {
307 defaultVtkWriter_->gridChanged();
308 }
309 }
310
320 template <class Context>
321 void boundary(BoundaryRateVector&,
322 const Context&,
323 unsigned,
324 unsigned) const
325 { throw std::logic_error("Problem does not provide a boundary() method"); }
326
337 template <class Context>
338 void constraints(Constraints&,
339 const Context&,
340 unsigned,
341 unsigned) const
342 { throw std::logic_error("Problem does not provide a constraints() method"); }
343
356 template <class Context>
357 void source(RateVector&,
358 const Context&,
359 unsigned,
360 unsigned) const
361 { throw std::logic_error("Problem does not provide a source() method"); }
362
373 template <class Context>
374 void initial(PrimaryVariables&,
375 const Context&,
376 unsigned,
377 unsigned) const
378 { throw std::logic_error("Problem does not provide a initial() method"); }
379
395 template <class Context>
396 Scalar extrusionFactor(const Context&,
397 unsigned,
398 unsigned) const
399 { return asImp_().extrusionFactor(); }
400
401 Scalar extrusionFactor() const
402 { return 1.0; }
403
409 {}
410
415 {}
416
421 {}
422
427 {}
428
433 {}
434
442 {}
443
450 {
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";
454 }
455
459 void finalize()
460 {
461 const auto& executionTimer = simulator().executionTimer();
462
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());
473 const unsigned threadsPerProcess = ThreadManager::maxThreads();
474 if (gridView().comm().rank() == 0) {
475 std::cout << std::setprecision(3)
476 << "Simulation of problem '" << asImp_().name() << "' finished.\n"
477 << "\n"
478 << "------------------------ Timing ------------------------\n"
479 << "Setup time: " << setupTime << " seconds"
480 << humanReadableTime(setupTime)
481 << ", " << setupTime / (executionTime + setupTime) * 100 << "%\n"
482 << "Simulation time: " << executionTime << " seconds"
483 << humanReadableTime(executionTime)
484 << ", " << executionTime / (executionTime + setupTime) * 100 << "%\n"
485 << " Linearization time: " << linearizeTime << " seconds"
486 << humanReadableTime(linearizeTime)
487 << ", " << linearizeTime / executionTime * 100 << "%\n"
488 << " Linear solve time: " << solveTime << " seconds"
489 << humanReadableTime(solveTime)
490 << ", " << solveTime / executionTime * 100 << "%\n"
491 << " Newton update time: " << updateTime << " seconds"
492 << humanReadableTime(updateTime)
493 << ", " << updateTime / executionTime * 100 << "%\n"
494 << " Pre/postprocess time: " << prePostProcessTime << " seconds"
495 << humanReadableTime(prePostProcessTime)
496 << ", " << prePostProcessTime / executionTime * 100 << "%\n"
497 << " Output write time: " << writeTime << " seconds"
498 << humanReadableTime(writeTime)
499 << ", " << writeTime / executionTime * 100 << "%\n"
500 << "First process' simulation CPU time: " << localCpuTime << " seconds"
501 << humanReadableTime(localCpuTime) << "\n"
502 << "Number of processes: " << numProcesses << "\n"
503 << "Threads per processes: " << threadsPerProcess << "\n"
504 << "Total CPU time: " << globalCpuTime << " seconds"
505 << humanReadableTime(globalCpuTime) << "\n"
506 << "\n"
507 << "----------------------------------------------------------------\n"
508 << std::endl;
509 }
510 }
511
517 {
518 const unsigned maxFails = asImp_().maxTimeIntegrationFailures();
519 Scalar minTimeStep = asImp_().minTimeStepSize();
520
521 std::string errorMessage;
522 for (unsigned i = 0; i < maxFails; ++i) {
523 bool converged = model().update();
524 if (converged) {
525 return;
526 }
527
528 const Scalar dt = simulator().timeStepSize();
529 Scalar nextDt = dt / 2.0;
530 if (dt < minTimeStep * (1.0 + 1e-9)) {
531 if (asImp_().continueOnConvergenceError()) {
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"
535 << std::flush;
536 }
537 return;
538 }
539 else {
540 errorMessage = "Time integration did not succeed with the minumum time step size of " +
541 std::to_string(double(minTimeStep)) + " seconds. Giving up!";
542 break; // give up: we can't make the time step smaller anymore!
543 }
544 }
545 else if (nextDt < minTimeStep) {
546 nextDt = minTimeStep;
547 }
548 simulator().setTimeStepSize(nextDt);
549
550 // update failed
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;
555 }
556 }
557
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()));
562 }
563 throw std::runtime_error(errorMessage);
564 }
565
569 Scalar minTimeStepSize() const
570 { return Parameters::Get<Parameters::MinTimeStepSize<Scalar>>(); }
571
577 { return Parameters::Get<Parameters::MaxTimeStepDivisions>(); }
578
585 { return Parameters::Get<Parameters::ContinueOnConvergenceError>(); }
586
590 void setNextTimeStepSize(Scalar dt)
591 { nextTimeStepSize_ = dt; }
592
598 Scalar nextTimeStepSize() const
599 {
600 if (nextTimeStepSize_ > 0.0) {
601 return nextTimeStepSize_;
602 }
603
604 Scalar dtNext = std::min(Parameters::Get<Parameters::MaxTimeStepSize<Scalar>>(),
605 newtonMethod().suggestTimeStepSize(simulator().timeStepSize()));
606
607 if (dtNext < simulator().maxTimeStepSize() &&
608 simulator().maxTimeStepSize() < dtNext * 2)
609 {
610 dtNext = simulator().maxTimeStepSize() / 2 * 1.01;
611 }
612
613 return dtNext;
614 }
615
625 {
626 return simulator().timeStepIndex() > 0 &&
627 (simulator().timeStepIndex() % 10 == 0);
628 }
629
638 bool shouldWriteOutput() const
639 { return true; }
640
647 { model().advanceTimeLevel(); }
648
656 std::string name() const
657 { return "sim"; }
658
662 const GridView& gridView() const
663 { return gridView_; }
664
669 const GlobalPosition& boundingBoxMin() const
670 { return boundingBoxMin_; }
671
676 const GlobalPosition& boundingBoxMax() const
677 { return boundingBoxMax_; }
678
682 const VertexMapper& vertexMapper() const
683 { return vertexMapper_; }
684
688 const ElementMapper& elementMapper() const
689 { return elementMapper_; }
690
694 Simulator& simulator()
695 { return simulator_; }
696
700 const Simulator& simulator() const
701 { return simulator_; }
702
706 Model& model()
707 { return simulator_.model(); }
708
712 const Model& model() const
713 { return simulator_.model(); }
714
718 NewtonMethod& newtonMethod()
719 { return model().newtonMethod(); }
720
724 const NewtonMethod& newtonMethod() const
725 { return model().newtonMethod(); }
726 // \}
727
733 {
735 }
736
745 {
746 return 0;
747 }
748
762 template <class Restarter>
763 void serialize(Restarter& res)
764 {
765 if (enableVtkOutput_()) {
766 defaultVtkWriter_->serialize(res);
767 }
768 }
769
780 template <class Restarter>
781 void deserialize(Restarter& res)
782 {
783 if (enableVtkOutput_()) {
784 defaultVtkWriter_->deserialize(res);
785 }
786 }
787
794 void writeOutput(bool verbose = true)
795 {
796 if (!enableVtkOutput_()) {
797 return;
798 }
799
800 if (verbose && gridView().comm().rank() == 0) {
801 std::cout << "Writing visualization results for the current time step.\n"
802 << std::flush;
803 }
804
805 // calculate the time _after_ the time was updated
806 const Scalar t = simulator().time() + simulator().timeStepSize();
807
808 defaultVtkWriter_->beginWrite(t);
809 model().prepareOutputFields();
810 model().appendOutputFields(*defaultVtkWriter_);
811 defaultVtkWriter_->endWrite(false);
812 }
813
819 { return *defaultVtkWriter_; }
820
825 { return iterationContext_; }
826
832
838
844
845private:
846 // LocalContextGuard needs mutable access to save/restore the context
847 // during NLDD domain-local solves.
848 template<class P> friend class LocalContextGuard;
849
853 NewtonIterationContext& mutableIterationContext()
854 { return iterationContext_; }
855
856public:
857
858protected:
861
862 bool enableVtkOutput_() const
863 { return Parameters::Get<Parameters::EnableVtkOutput>(); }
864
865private:
867 Implementation& asImp_()
868 { return *static_cast<Implementation*>(this); }
869
871 const Implementation& asImp_() const
872 { return *static_cast<const Implementation*>(this); }
873
874 // Grid management stuff
875 const GridView gridView_;
876 ElementMapper elementMapper_;
877 VertexMapper vertexMapper_;
878 GlobalPosition boundingBoxMin_;
879 GlobalPosition boundingBoxMax_;
880
881 // Attributes required for the actual simulation
882 Simulator& simulator_;
883 std::unique_ptr<VtkMultiWriter> defaultVtkWriter_{};
884};
885
886} // namespace Opm
887
888#endif
Definition: restrictprolong.hh:142
Base class for all problems which use a finite volume spatial discretization.
Definition: fvbaseproblem.hh:71
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
NewtonIterationContext iterationContext_
Definition: fvbaseproblem.hh:860
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
Scalar extrusionFactor() const
Definition: fvbaseproblem.hh:401
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
Scalar nextTimeStepSize_
Definition: fvbaseproblem.hh:859
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
bool enableVtkOutput_() const
Definition: fvbaseproblem.hh:862
EmptyRestrictProlong RestrictProlongOperator
Definition: fvbaseproblem.hh:107
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. as a VTK file)
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
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.
Definition: fvbaseprimaryvariables.hh:161
auto Get(bool errorIfNotRegistered=true)
Retrieve a runtime parameter.
Definition: parametersystem.hpp:187
Definition: blackoilmodel.hh:80
Definition: blackoilbioeffectsmodules.hh:45
std::string humanReadableTime(double timeInSeconds, bool isAmendment=true)
Given a time step size in seconds, return it in a format which is more easily parsable by humans.
std::string simulatorOutputDir()
Determine and check the configured directory for simulation output.
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 to_string(const ConvergenceReport::ReservoirFailure::Type t)
Context for iteration-dependent decisions in the Newton solver.
Definition: NewtonIterationContext.hpp:43
void markTimestepInitialized()
State Mutations.
Definition: NewtonIterationContext.hpp:104
void advanceIteration()
Definition: NewtonIterationContext.hpp:112
void resetForNewTimestep()
Reset all state for a new timestep.
Definition: NewtonIterationContext.hpp:122
Specify the maximum size of a time integration [s].
Definition: fvbaseparameters.hh:106