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 "fvbaseproperties.hh"
32
36
37#include <dune/common/fvector.hh>
38
39#include <iostream>
40#include <limits>
41#include <string>
42
43#include <sys/stat.h>
44
45namespace Opm::Properties {
46
47template <class TypeTag, class MyTypeTag>
48struct NewtonMethod;
49
50} // namespace Opm::Properties
51
52namespace Opm {
53
63template<class TypeTag>
65{
66private:
67 using Implementation = GetPropType<TypeTag, Properties::Problem>;
69
70 static const int vtkOutputFormat = getPropValue<TypeTag, Properties::VtkOutputFormat>();
72
78
81
86
87 enum {
88 dim = GridView::dimension,
89 dimWorld = GridView::dimensionworld
90 };
91
92 using Element = typename GridView::template Codim<0>::Entity;
93 using Vertex = typename GridView::template Codim<dim>::Entity;
94 using VertexIterator = typename GridView::template Codim<dim>::Iterator;
95
96 using CoordScalar = typename GridView::Grid::ctype;
97 using GlobalPosition = Dune::FieldVector<CoordScalar, dimWorld>;
98
99public:
100 // the default restriction and prolongation for adaptation is simply an empty one
102
103private:
104 // copying a problem is not a good idea
105 FvBaseProblem(const FvBaseProblem& ) = delete;
106
107public:
116 : nextTimeStepSize_(0.0)
117 , gridView_(simulator.gridView())
118 , elementMapper_(gridView_, Dune::mcmgElementLayout())
119 , vertexMapper_(gridView_, Dune::mcmgVertexLayout())
120 , boundingBoxMin_(std::numeric_limits<double>::max())
121 , boundingBoxMax_(-std::numeric_limits<double>::max())
122 , simulator_(simulator)
123 , defaultVtkWriter_(0)
124 {
125 // calculate the bounding box of the local partition of the grid view
126 VertexIterator vIt = gridView_.template begin<dim>();
127 const VertexIterator vEndIt = gridView_.template end<dim>();
128 for (; vIt!=vEndIt; ++vIt) {
129 for (unsigned i=0; i<dim; i++) {
130 boundingBoxMin_[i] = std::min(boundingBoxMin_[i], vIt->geometry().corner(0)[i]);
131 boundingBoxMax_[i] = std::max(boundingBoxMax_[i], vIt->geometry().corner(0)[i]);
132 }
133 }
134
135 // communicate to get the bounding box of the whole domain
136 for (unsigned i = 0; i < dim; ++i) {
137 boundingBoxMin_[i] = gridView_.comm().min(boundingBoxMin_[i]);
138 boundingBoxMax_[i] = gridView_.comm().max(boundingBoxMax_[i]);
139 }
140
141 if (enableVtkOutput_()) {
142 bool asyncVtkOutput =
143 simulator_.gridView().comm().size() == 1 &&
144 Parameters::get<TypeTag, Properties::EnableAsyncVtkOutput>();
145
146 // asynchonous VTK output currently does not work in conjunction with grid
147 // adaptivity because the async-IO code assumes that the grid stays
148 // constant. complain about that case.
149 bool enableGridAdaptation = Parameters::get<TypeTag, Properties::EnableGridAdaptation>();
150 if (asyncVtkOutput && enableGridAdaptation)
151 throw std::runtime_error("Asynchronous VTK output currently cannot be used "
152 "at the same time as grid adaptivity");
153
154 std::string outputDir = asImp_().outputDir();
155
156 defaultVtkWriter_ =
157 new VtkMultiWriter(asyncVtkOutput, gridView_, outputDir, asImp_().name());
158 }
159 }
160
162 { delete defaultVtkWriter_; }
163
168 static void registerParameters()
169 {
170 Model::registerParameters();
171 Parameters::registerParam<TypeTag, Properties::MaxTimeStepSize>
172 ("The maximum size to which all time steps are limited to [s]");
173 Parameters::registerParam<TypeTag, Properties::MinTimeStepSize>
174 ("The minimum size to which all time steps are limited to [s]");
175 Parameters::registerParam<TypeTag, Properties::MaxTimeStepDivisions>
176 ("The maximum number of divisions by two of the timestep size "
177 "before the simulation bails out");
178 Parameters::registerParam<TypeTag, Properties::EnableAsyncVtkOutput>
179 ("Dispatch a separate thread to write the VTK output");
180 Parameters::registerParam<TypeTag, Properties::ContinueOnConvergenceError>
181 ("Continue with a non-converged solution instead of giving up "
182 "if we encounter a time step size smaller than the minimum time "
183 "step size.");
184 }
185
194 { return true; }
195
204 std::string outputDir() const
205 {
206 std::string outputDir = Parameters::get<TypeTag, Properties::OutputDir>();
207
208 if (outputDir.empty())
209 outputDir = ".";
210
211 // TODO: replace this by std::filesystem once we require c++-2017
212 struct stat st;
213 if (::stat(outputDir.c_str(), &st) != 0)
214 throw std::runtime_error("Could not access output directory '"+outputDir+"':"
215 +strerror(errno));
216 if (!S_ISDIR(st.st_mode))
217 throw std::runtime_error("Path to output directory '"+outputDir+"' exists but is not a directory");
218
219 if (access(outputDir.c_str(), W_OK) != 0)
220 throw std::runtime_error("Output directory '"+outputDir+"' exists but is not writeable");
221
222 return outputDir;
223 }
224
231 static std::string helpPreamble(int,
232 const char **argv)
233 {
234 std::string desc = Implementation::briefDescription();
235 if (!desc.empty())
236 desc = desc + "\n";
237
238 return
239 "Usage: "+std::string(argv[0]) + " [OPTIONS]\n"
240 + desc;
241 }
242
249 static std::string briefDescription()
250 { return ""; }
251
252 // TODO (?): detailedDescription()
253
270 static int handlePositionalParameter(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#if DUNE_VERSION_NEWER(DUNE_GRID, 2, 8)
304 elementMapper_.update(gridView_);
305 vertexMapper_.update(gridView_);
306#else
307 elementMapper_.update();
308 vertexMapper_.update();
309#endif
310
311 if (enableVtkOutput_())
312 defaultVtkWriter_->gridChanged();
313 }
314
324 template <class Context>
325 void boundary(BoundaryRateVector&,
326 const Context&,
327 unsigned,
328 unsigned) const
329 { throw std::logic_error("Problem does not provide a boundary() method"); }
330
341 template <class Context>
342 void constraints(Constraints&,
343 const Context&,
344 unsigned,
345 unsigned) const
346 { throw std::logic_error("Problem does not provide a constraints() method"); }
347
360 template <class Context>
361 void source(RateVector&,
362 const Context&,
363 unsigned,
364 unsigned) const
365 { throw std::logic_error("Problem does not provide a source() method"); }
366
377 template <class Context>
378 void initial(PrimaryVariables&,
379 const Context&,
380 unsigned,
381 unsigned) const
382 { throw std::logic_error("Problem does not provide a initial() method"); }
383
399 template <class Context>
400 Scalar extrusionFactor(const Context&,
401 unsigned,
402 unsigned) const
403 { return asImp_().extrusionFactor(); }
404
405 Scalar extrusionFactor() const
406 { return 1.0; }
407
413 {}
414
419 { }
420
425 { }
426
431 { }
432
437 { }
438
446 { }
447
454 {
455 std::cerr << "The end of episode " << simulator().episodeIndex() + 1 << " has been "
456 << "reached, but the problem does not override the endEpisode() method. "
457 << "Doing nothing!\n";
458 }
459
463 void finalize()
464 {
465 const auto& executionTimer = simulator().executionTimer();
466
467 Scalar executionTime = executionTimer.realTimeElapsed();
468 Scalar setupTime = simulator().setupTimer().realTimeElapsed();
469 Scalar prePostProcessTime = simulator().prePostProcessTimer().realTimeElapsed();
470 Scalar localCpuTime = executionTimer.cpuTimeElapsed();
471 Scalar globalCpuTime = executionTimer.globalCpuTimeElapsed();
472 Scalar writeTime = simulator().writeTimer().realTimeElapsed();
473 Scalar linearizeTime = simulator().linearizeTimer().realTimeElapsed();
474 Scalar solveTime = simulator().solveTimer().realTimeElapsed();
475 Scalar updateTime = simulator().updateTimer().realTimeElapsed();
476 unsigned numProcesses = static_cast<unsigned>(this->gridView().comm().size());
477 unsigned threadsPerProcess = ThreadManager::maxThreads();
478 if (gridView().comm().rank() == 0) {
479 std::cout << std::setprecision(3)
480 << "Simulation of problem '" << asImp_().name() << "' finished.\n"
481 << "\n"
482 << "------------------------ Timing ------------------------\n"
483 << "Setup time: " << setupTime << " seconds" << Simulator::humanReadableTime(setupTime)
484 << ", " << setupTime/(executionTime + setupTime)*100 << "%\n"
485 << "Simulation time: " << executionTime << " seconds" << Simulator::humanReadableTime(executionTime)
486 << ", " << executionTime/(executionTime + setupTime)*100 << "%\n"
487 << " Linearization time: " << linearizeTime << " seconds" << Simulator::humanReadableTime(linearizeTime)
488 << ", " << linearizeTime/executionTime*100 << "%\n"
489 << " Linear solve time: " << solveTime << " seconds" << Simulator::humanReadableTime(solveTime)
490 << ", " << solveTime/executionTime*100 << "%\n"
491 << " Newton update time: " << updateTime << " seconds" << Simulator::humanReadableTime(updateTime)
492 << ", " << updateTime/executionTime*100 << "%\n"
493 << " Pre/postprocess time: " << prePostProcessTime << " seconds" << Simulator::humanReadableTime(prePostProcessTime)
494 << ", " << prePostProcessTime/executionTime*100 << "%\n"
495 << " Output write time: " << writeTime << " seconds" << Simulator::humanReadableTime(writeTime)
496 << ", " << writeTime/executionTime*100 << "%\n"
497 << "First process' simulation CPU time: " << localCpuTime << " seconds" << Simulator::humanReadableTime(localCpuTime) << "\n"
498 << "Number of processes: " << numProcesses << "\n"
499 << "Threads per processes: " << threadsPerProcess << "\n"
500 << "Total CPU time: " << globalCpuTime << " seconds" << Simulator::humanReadableTime(globalCpuTime) << "\n"
501 << "\n"
502 << "----------------------------------------------------------------\n"
503 << std::endl;
504 }
505 }
506
512 {
513 unsigned maxFails = asImp_().maxTimeIntegrationFailures();
514 Scalar minTimeStepSize = asImp_().minTimeStepSize();
515
516 std::string errorMessage;
517 for (unsigned i = 0; i < maxFails; ++i) {
518 bool converged = model().update();
519 if (converged)
520 return;
521
522 Scalar dt = simulator().timeStepSize();
523 Scalar nextDt = dt / 2.0;
524 if (dt < minTimeStepSize*(1 + 1e-9)) {
525 if (asImp_().continueOnConvergenceError()) {
526 if (gridView().comm().rank() == 0)
527 std::cout << "Newton solver did not converge with minimum time step of "
528 << dt << " seconds. Continuing with unconverged solution!\n"
529 << std::flush;
530 return;
531 }
532 else {
533 errorMessage =
534 "Time integration did not succeed with the minumum time step size of "
535 + std::to_string(double(minTimeStepSize)) + " seconds. Giving up!";
536 break; // give up: we can't make the time step smaller anymore!
537 }
538 }
539 else if (nextDt < minTimeStepSize)
540 nextDt = minTimeStepSize;
541 simulator().setTimeStepSize(nextDt);
542
543 // update failed
544 if (gridView().comm().rank() == 0)
545 std::cout << "Newton solver did not converge with "
546 << "dt=" << dt << " seconds. Retrying with time step of "
547 << nextDt << " seconds\n" << std::flush;
548 }
549
550 if (errorMessage.empty())
551 errorMessage =
552 "Newton solver didn't converge after "
553 +std::to_string(maxFails)+" time-step divisions. dt="
554 +std::to_string(double(simulator().timeStepSize()));
555 throw std::runtime_error(errorMessage);
556 }
557
561 Scalar minTimeStepSize() const
562 { return Parameters::get<TypeTag, Properties::MinTimeStepSize>(); }
563
569 { return Parameters::get<TypeTag, Properties::MaxTimeStepDivisions>(); }
570
577 { return Parameters::get<TypeTag, Properties::ContinueOnConvergenceError>(); }
578
582 void setNextTimeStepSize(Scalar dt)
583 { nextTimeStepSize_ = dt; }
584
590 Scalar nextTimeStepSize() const
591 {
592 if (nextTimeStepSize_ > 0.0)
593 return nextTimeStepSize_;
594
595 Scalar dtNext = std::min(Parameters::get<TypeTag, Properties::MaxTimeStepSize>(),
596 newtonMethod().suggestTimeStepSize(simulator().timeStepSize()));
597
598 if (dtNext < simulator().maxTimeStepSize()
599 && simulator().maxTimeStepSize() < dtNext*2)
600 {
601 dtNext = simulator().maxTimeStepSize()/2 * 1.01;
602 }
603
604 return dtNext;
605 }
606
616 {
617 return simulator().timeStepIndex() > 0 &&
618 (simulator().timeStepIndex() % 10 == 0);
619 }
620
629 bool shouldWriteOutput() const
630 { return true; }
631
638 { model().advanceTimeLevel(); }
639
647 std::string name() const
648 { return "sim"; }
649
653 const GridView& gridView() const
654 { return gridView_; }
655
660 const GlobalPosition& boundingBoxMin() const
661 { return boundingBoxMin_; }
662
667 const GlobalPosition& boundingBoxMax() const
668 { return boundingBoxMax_; }
669
673 const VertexMapper& vertexMapper() const
674 { return vertexMapper_; }
675
679 const ElementMapper& elementMapper() const
680 { return elementMapper_; }
681
685 Simulator& simulator()
686 { return simulator_; }
687
691 const Simulator& simulator() const
692 { return simulator_; }
693
697 Model& model()
698 { return simulator_.model(); }
699
703 const Model& model() const
704 { return simulator_.model(); }
705
709 NewtonMethod& newtonMethod()
710 { return model().newtonMethod(); }
711
715 const NewtonMethod& newtonMethod() const
716 { return model().newtonMethod(); }
717 // \}
718
724 {
726 }
727
736 {
737 return 0;
738 }
739
753 template <class Restarter>
754 void serialize(Restarter& res)
755 {
756 if (enableVtkOutput_())
757 defaultVtkWriter_->serialize(res);
758 }
759
770 template <class Restarter>
771 void deserialize(Restarter& res)
772 {
773 if (enableVtkOutput_())
774 defaultVtkWriter_->deserialize(res);
775 }
776
783 void writeOutput(bool verbose = true)
784 {
785 if (!enableVtkOutput_())
786 return;
787
788 if (verbose && gridView().comm().rank() == 0)
789 std::cout << "Writing visualization results for the current time step.\n"
790 << std::flush;
791
792 // calculate the time _after_ the time was updated
793 Scalar t = simulator().time() + simulator().timeStepSize();
794
795 defaultVtkWriter_->beginWrite(t);
796 model().prepareOutputFields();
797 model().appendOutputFields(*defaultVtkWriter_);
798 defaultVtkWriter_->endWrite();
799
800 }
801
807 { return defaultVtkWriter_; }
808
809protected:
811
812private:
813 bool enableVtkOutput_() const
814 { return Parameters::get<TypeTag, Properties::EnableVtkOutput>(); }
815
817 Implementation& asImp_()
818 { return *static_cast<Implementation *>(this); }
819
821 const Implementation& asImp_() const
822 { return *static_cast<const Implementation *>(this); }
823
824 // Grid management stuff
825 const GridView gridView_;
826 ElementMapper elementMapper_;
827 VertexMapper vertexMapper_;
828 GlobalPosition boundingBoxMin_;
829 GlobalPosition boundingBoxMax_;
830
831 // Attributes required for the actual simulation
832 Simulator& simulator_;
833 mutable VtkMultiWriter *defaultVtkWriter_;
834};
835
836} // namespace Opm
837
838#endif
Definition: restrictprolong.hh:142
Base class for all problems which use a finite volume spatial discretization.
Definition: fvbaseproblem.hh:65
std::string name() const
The problem name.
Definition: fvbaseproblem.hh:647
const Simulator & simulator() const
Returns Simulator object used by the simulation.
Definition: fvbaseproblem.hh:691
void writeOutput(bool verbose=true)
Write the relevant secondary variables of the current solution into an VTK output file.
Definition: fvbaseproblem.hh:783
void boundary(BoundaryRateVector &, const Context &, unsigned, unsigned) const
Evaluate the boundary conditions for a boundary segment.
Definition: fvbaseproblem.hh:325
void finalize()
Called after the simulation has been run sucessfully.
Definition: fvbaseproblem.hh:463
Simulator & simulator()
Returns Simulator object used by the simulation.
Definition: fvbaseproblem.hh:685
const Model & model() const
Returns numerical model used for the problem.
Definition: fvbaseproblem.hh:703
static void registerParameters()
Registers all available parameters for the problem and the model.
Definition: fvbaseproblem.hh:168
void source(RateVector &, const Context &, unsigned, unsigned) const
Evaluate the source term for all phases within a given sub-control-volume.
Definition: fvbaseproblem.hh:361
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:193
unsigned maxTimeIntegrationFailures() const
Returns the maximum number of subsequent failures for the time integration before giving up.
Definition: fvbaseproblem.hh:568
unsigned markForGridAdaptation()
Mark grid cells for refinement or coarsening.
Definition: fvbaseproblem.hh:735
void serialize(Restarter &res)
This method writes the complete state of the problem to the harddisk.
Definition: fvbaseproblem.hh:754
void beginIteration()
Called by the simulator before each Newton-Raphson iteration.
Definition: fvbaseproblem.hh:430
~FvBaseProblem()
Definition: fvbaseproblem.hh:161
void endTimeStep()
Called by the simulator after each time integration.
Definition: fvbaseproblem.hh:445
FvBaseProblem(Simulator &simulator)
Definition: fvbaseproblem.hh:115
const GlobalPosition & boundingBoxMax() const
The coordinate of the corner of the GridView's bounding box with the largest values.
Definition: fvbaseproblem.hh:667
RestrictProlongOperator restrictProlongOperator()
return restriction and prolongation operator
Definition: fvbaseproblem.hh:723
void timeIntegration()
Called by Opm::Simulator in order to do a time integration on the model.
Definition: fvbaseproblem.hh:511
Model & model()
Returns numerical model used for the problem.
Definition: fvbaseproblem.hh:697
Scalar extrusionFactor() const
Definition: fvbaseproblem.hh:405
void endEpisode()
Called when the end of an simulation episode is reached.
Definition: fvbaseproblem.hh:453
void beginTimeStep()
Called by the simulator before each time integration.
Definition: fvbaseproblem.hh:424
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:771
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:400
bool shouldWriteRestartFile() const
Returns true if a restart file should be written to disk.
Definition: fvbaseproblem.hh:615
Scalar nextTimeStepSize() const
Called by Opm::Simulator whenever a solution for a time step has been computed and the simulation tim...
Definition: fvbaseproblem.hh:590
Scalar minTimeStepSize() const
Returns the minimum allowable size of a time step.
Definition: fvbaseproblem.hh:561
static int handlePositionalParameter(std::set< std::string > &, std::string &errorMsg, int, const char **argv, int paramIdx, int)
Handles positional command line parameters.
Definition: fvbaseproblem.hh:270
void beginEpisode()
Called at the beginning of an simulation episode.
Definition: fvbaseproblem.hh:418
void initial(PrimaryVariables &, const Context &, unsigned, unsigned) const
Evaluate the initial value for a control volume.
Definition: fvbaseproblem.hh:378
std::string outputDir() const
Determine the directory for simulation output.
Definition: fvbaseproblem.hh:204
Scalar nextTimeStepSize_
Definition: fvbaseproblem.hh:810
const VertexMapper & vertexMapper() const
Returns the mapper for vertices to indices.
Definition: fvbaseproblem.hh:673
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:660
static std::string briefDescription()
Returns a human readable description of the problem for the help message.
Definition: fvbaseproblem.hh:249
void endIteration()
Called by the simulator after each Newton-Raphson update.
Definition: fvbaseproblem.hh:436
void setNextTimeStepSize(Scalar dt)
Impose the next time step size to be used externally.
Definition: fvbaseproblem.hh:582
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:576
void constraints(Constraints &, const Context &, unsigned, unsigned) const
Evaluate the constraints for a control volume.
Definition: fvbaseproblem.hh:342
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:715
const GridView & gridView() const
The GridView which used by the problem.
Definition: fvbaseproblem.hh:653
EmptyRestrictProlong RestrictProlongOperator
Definition: fvbaseproblem.hh:101
NewtonMethod & newtonMethod()
Returns object which implements the Newton method.
Definition: fvbaseproblem.hh:709
bool shouldWriteOutput() const
Returns true if the current solution should be written to disk (i.e. as a VTK file)
Definition: fvbaseproblem.hh:629
void initialSolutionApplied()
Callback used by the model to indicate that the initial solution has been determined for all degrees ...
Definition: fvbaseproblem.hh:412
void advanceTimeLevel()
Called by the simulator after everything which can be done about the current time step is finished an...
Definition: fvbaseproblem.hh:637
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:806
const ElementMapper & elementMapper() const
Returns the mapper for elements to indices.
Definition: fvbaseproblem.hh:679
static std::string humanReadableTime(Scalar timeInSeconds, bool isAmendment=true)
Given a time step size in seconds, return it in a format which is more easily parsable by humans.
Definition: simulator.hh:822
static unsigned maxThreads()
Return the maximum number of threads of the current process.
Definition: threadmanager.hh:118
Simplifies writing multi-file VTK datasets.
Definition: vtkmultiwriter.hh:66
void endWrite(bool onlyDiscard=false)
Finalizes the current writer.
Definition: vtkmultiwriter.hh:389
void deserialize(Restarter &res)
Read the multi-writer's state from a restart file.
Definition: vtkmultiwriter.hh:442
void gridChanged()
Updates the internal data structures after mesh refinement.
Definition: vtkmultiwriter.hh:165
void beginWrite(double t)
Called whenever a new time step must be written.
Definition: vtkmultiwriter.hh:179
void serialize(Restarter &res)
Write the multi-writer's state to a restart file.
Definition: vtkmultiwriter.hh:408
Declare the properties used by the infrastructure code of the finite volume discretizations.
Definition: fvbaseprimaryvariables.hh:141
Definition: blackoilmodel.hh:72
Definition: blackoilboundaryratevector.hh:37
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:242