20#ifndef DUMUX_ONE_MODEL_PROBLEM_HH
21#define DUMUX_ONE_MODEL_PROBLEM_HH
23#include <dune/common/shared_ptr.hh>
44template<
class TypeTag>
50 using Grid =
typename GridView::Grid;
63 using VertexMapper =
typename SolutionTypes::VertexMapper;
64 using ElementMapper =
typename SolutionTypes::ElementMapper;
68 dim = GridView::dimension,
69 dimWorld = GridView::dimensionworld
73 wetting = 0, nonwetting = 1
76 using GlobalPosition = Dune::FieldVector<Scalar, dimWorld>;
77 using Element =
typename GridView::Traits::template Codim<0>::Entity;
78 using Intersection =
typename GridView::Intersection;
80 using PrimaryVariables =
typename SolutionTypes::PrimaryVariables;
95 : gridView_(grid.leafGridView()),
96 bBoxMin_(std::numeric_limits<double>::max()),
97 bBoxMax_(-std::numeric_limits<double>::max()),
98 variables_(grid.leafGridView()),
100 outputTimeInterval_(0)
105 for (
const auto& vertex : vertices(grid.leafGridView())) {
106 for (
int i=0; i<dim; i++) {
107 bBoxMin_[i] = min(bBoxMin_[i], vertex.geometry().center()[i]);
108 bBoxMax_[i] = max(bBoxMax_[i], vertex.geometry().center()[i]);
112 timeManager_ = std::make_shared<TimeManager>(verbose);
114 model_ = std::make_shared<Model>(asImp_()) ;
115 maxTimeStepSize_ = getParam<Scalar>(
"TimeManager.MaxTimeStepSize", std::numeric_limits<Scalar>::max());
124 : gridView_(grid.leafGridView()),
125 bBoxMin_(std::numeric_limits<double>::max()),
126 bBoxMax_(-std::numeric_limits<double>::max()),
127 variables_(grid.leafGridView()),
129 outputTimeInterval_(0)
134 for (
const auto& vertex : vertices(grid.leafGridView())) {
135 for (
int i=0; i<dim; i++) {
136 bBoxMin_[i] = min(bBoxMin_[i], vertex.geometry().center()[i]);
137 bBoxMax_[i] = max(bBoxMax_[i], vertex.geometry().center()[i]);
141 timeManager_ = Dune::stackobject_to_shared_ptr<TimeManager>(
timeManager);
143 model_ = std::make_shared<Model>(asImp_()) ;
144 maxTimeStepSize_ = getParam<Scalar>(
"TimeManager.MaxTimeStepSize", std::numeric_limits<Scalar>::max());
155 const Intersection &intersection)
const
158 asImp_().boundaryTypesAtPos(bcTypes, intersection.geometry().center());
169 const GlobalPosition &globalPos)
const
173 DUNE_THROW(Dune::InvalidStateException,
174 "The problem does not provide "
175 "a boundaryTypesAtPos() method.");
188 const Intersection &intersection)
const
191 asImp_().dirichletAtPos(values, intersection.geometry().center());
204 const GlobalPosition &globalPos)
const
208 DUNE_THROW(Dune::InvalidStateException,
209 "The problem specifies that some boundary "
210 "segments are dirichlet, but does not provide "
211 "a dirichletAtPos() method.");
225 const Intersection &intersection)
const
228 asImp_().neumannAtPos(values, intersection.geometry().center());
242 const GlobalPosition &globalPos)
const
246 DUNE_THROW(Dune::InvalidStateException,
247 "The problem specifies that some boundary "
248 "segments are neumann, but does not provide "
249 "a neumannAtPos() method.");
263 const Element &element)
const
266 asImp_().sourceAtPos(values, element.geometry().center());
283 const GlobalPosition &globalPos)
const
285 DUNE_THROW(Dune::InvalidStateException,
286 "The problem does not provide "
287 "a sourceAtPos() method.");
300 const Element &element)
const
303 asImp_().initialAtPos(values, element.geometry().center());
317 const GlobalPosition &globalPos)
const
330 variables_.initialize();
331 model().initialize();
373 {
return maxTimeStepSize_; }
417 outputTimeInterval_ = (timeInterval > 0.0) ? timeInterval : 1e100;
418 timeManager().startNextEpisode(outputTimeInterval_);
427 { outputInterval_ = interval; }
440 if (outputInterval_ > 0)
442 if (
timeManager().timeStepIndex() % outputInterval_ == 0
463 if (verbose &&
gridView().comm().rank() == 0)
464 std::cout <<
"Writing result file for current time step\n";
466 resultWriter_ = std::make_shared<VtkMultiWriter>(
gridView(), asImp_().
name());
468 model().addOutputVtkFields(*resultWriter_);
469 asImp_().addOutputVtkFields();
470 resultWriter_->endWrite();
478 if (outputTimeInterval_ > 0.0 && !
timeManager().finished())
480 timeManager().startNextEpisode(outputTimeInterval_);
484 std::cerr <<
"The end of an episode is reached, but the problem "
485 <<
"does not override the episodeEnd() method. "
486 <<
"Doing nothing!\n";
499 const std::string&
name()
const
520 {
return gridView_; }
526 {
return variables_.vertexMapper(); }
532 {
return variables_.elementMapper(); }
552 {
return *timeManager_; }
558 {
return *timeManager_; }
564 {
return variables_; }
570 {
return variables_; }
606 std::cout <<
"Serialize to file " << res.fileName() <<
"\n";
610 res.template deserializeEntities<0> (
model(), gridView_);
627 std::cout <<
"Deserialize from file " << res.fileName() <<
"\n";
631 res.template deserializeEntities<0> (
model(), gridView_);
633 res.deserializeEnd();
642 resultWriter_ = std::make_shared<VtkMultiWriter>(gridView_, asImp_().
name());
643 return *resultWriter_;
649 resultWriter_ = std::make_shared<VtkMultiWriter>(gridView_, asImp_().
name());
650 return *resultWriter_;
655 Implementation &asImp_()
656 {
return *
static_cast<Implementation *
>(
this); }
659 const Implementation &asImp_()
const
660 {
return *
static_cast<const Implementation *
>(
this); }
662 std::string simname_;
665 const GridView gridView_;
667 GlobalPosition bBoxMin_;
668 GlobalPosition bBoxMax_;
670 std::shared_ptr<TimeManager> timeManager_;
671 Scalar maxTimeStepSize_;
673 Variables variables_;
675 std::shared_ptr<Model> model_;
677 std::shared_ptr<VtkMultiWriter> resultWriter_;
679 Scalar outputTimeInterval_;
Generic class to read/write restart files.
Simplifies writing multi-file VTK datasets.
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type GetProp
get the type of a property (equivalent to old macro GET_PROP(...))
Definition: propertysystem.hh:140
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:149
Load or save a state of a model to/from the harddisk.
Definition: restart.hh:44
void deserializeBegin(Problem &problem, double t)
Start reading a restart file at a certain simulated time.
Definition: restart.hh:168
void serializeBegin(Problem &problem)
Write the current state of the model to disk.
Definition: restart.hh:94
Simplifies writing multi-file VTK datasets.
Definition: vtkmultiwriter.hh:61
void serialize(Restarter &res)
Write the multi-writer's state to a restart file.
Definition: vtkmultiwriter.hh:289
void deserialize(Restarter &res)
Read the multi-writer's state from a restart file.
Definition: vtkmultiwriter.hh:323
Base class for definition of an sequential diffusion (pressure) or transport problem.
Definition: onemodelproblem.hh:46
Scalar timeStepSize() const
Returns the current time step size [seconds].
Definition: onemodelproblem.hh:378
const VertexMapper & vertexMapper() const
Returns the mapper for vertices to indices.
Definition: onemodelproblem.hh:525
void initialAtPos(PrimaryVariables &values, const GlobalPosition &globalPos) const
Evaluate the initial value for a control volume.
Definition: onemodelproblem.hh:316
void timeIntegration()
Called by TimeManager in order to do a time integration on the model.
Definition: onemodelproblem.hh:345
void setName(const std::string &newName)
Set the problem name.
Definition: onemodelproblem.hh:511
void setTimeStepSize(Scalar dt)
Sets the current time step size [seconds].
Definition: onemodelproblem.hh:384
void dirichlet(PrimaryVariables &values, const Intersection &intersection) const
Evaluate the boundary conditions for a dirichlet control volume.
Definition: onemodelproblem.hh:187
const ElementMapper & elementMapper() const
Returns the mapper for elements to indices.
Definition: onemodelproblem.hh:531
TimeManager & timeManager()
Returns TimeManager object used by the simulation.
Definition: onemodelproblem.hh:551
const TimeManager & timeManager() const
Returns TimeManager object used by the simulation.
Definition: onemodelproblem.hh:557
Scalar maxTimeStepSize() const
Returns the user specified maximum time step size.
Definition: onemodelproblem.hh:372
bool shouldWriteRestartFile() const
Returns true if a restart file should be written to disk.
Definition: onemodelproblem.hh:403
VtkMultiWriter & resultWriter()
Definition: onemodelproblem.hh:639
void restart(double tRestart)
This method restores the complete state of the problem from disk.
Definition: onemodelproblem.hh:621
void postTimeStep()
Called by TimeManager whenever a solution for a timestep has been computed and the simulation time ha...
Definition: onemodelproblem.hh:356
void writeOutput(bool verbose=true)
Write the fields current solution into an VTK output file.
Definition: onemodelproblem.hh:461
Model & model()
Returns numerical model used for the problem.
Definition: onemodelproblem.hh:575
void advanceTimeLevel()
Called by the time manager after everything which can be done about the current time step is finished...
Definition: onemodelproblem.hh:364
void neumann(PrimaryVariables &values, const Intersection &intersection) const
Evaluate the boundary conditions for a neumann boundary segment.
Definition: onemodelproblem.hh:224
void setOutputInterval(int interval)
Sets the interval for Output.
Definition: onemodelproblem.hh:426
bool shouldWriteOutput() const
Returns true if the current solution should be written to disk (i.e. as a VTK file)
Definition: onemodelproblem.hh:438
VtkMultiWriter & resultWriter() const
Definition: onemodelproblem.hh:646
void episodeEnd()
Called when the end of an simulation episode is reached.
Definition: onemodelproblem.hh:476
OneModelProblem(TimeManager &timeManager, Grid &grid)
Constructs an object of type OneModelProblemProblem.
Definition: onemodelproblem.hh:123
const Variables & variables() const
Returns variables object.
Definition: onemodelproblem.hh:569
void setOutputTimeInterval(const Scalar timeInterval)
Sets a time interval for Output.
Definition: onemodelproblem.hh:415
void init()
Called by the TimeManager in order to initialize the problem.
Definition: onemodelproblem.hh:327
void dirichletAtPos(PrimaryVariables &values, const GlobalPosition &globalPos) const
Evaluate the boundary conditions for a dirichlet control volume.
Definition: onemodelproblem.hh:203
Scalar nextTimeStepSize(Scalar dt)
Called by TimeManager whenever a solution for a timestep has been computed and the simulation time ha...
Definition: onemodelproblem.hh:392
void addOutputVtkFields()
Definition: onemodelproblem.hh:457
void source(PrimaryVariables &values, const Element &element) const
Evaluate the source term.
Definition: onemodelproblem.hh:262
void boundaryTypesAtPos(BoundaryTypes &bcTypes, const GlobalPosition &globalPos) const
Specifies which kind of boundary condition should be used for which equation on a given boundary segm...
Definition: onemodelproblem.hh:168
const GlobalPosition & bBoxMin() const
The coordinate of the corner of the GridView's bounding box with the smallest values.
Definition: onemodelproblem.hh:538
const Model & model() const
Returns numerical model used for the problem.
Definition: onemodelproblem.hh:581
void boundaryTypes(BoundaryTypes &bcTypes, const Intersection &intersection) const
Specifies which kind of boundary condition should be used for which equation on a given boundary segm...
Definition: onemodelproblem.hh:154
void initial(PrimaryVariables &values, const Element &element) const
Evaluate the initial value for a control volume.
Definition: onemodelproblem.hh:299
const GlobalPosition & bBoxMax() const
The coordinate of the corner of the GridView's bounding box with the largest values.
Definition: onemodelproblem.hh:545
OneModelProblem(Grid &grid, bool verbose=true)
Constructs an object of type OneModelProblemProblem.
Definition: onemodelproblem.hh:94
const GridView & gridView() const
The GridView which used by the problem.
Definition: onemodelproblem.hh:519
const std::string & name() const
The problem name.
Definition: onemodelproblem.hh:499
void preTimeStep()
Called by TimeManager just before the time integration.
Definition: onemodelproblem.hh:338
void sourceAtPos(PrimaryVariables &values, const GlobalPosition &globalPos) const
Evaluate the source term for all phases within a given sub-control-volume.
Definition: onemodelproblem.hh:282
void serialize()
This method writes the complete state of the problem to the harddisk.
Definition: onemodelproblem.hh:600
void neumannAtPos(PrimaryVariables &values, const GlobalPosition &globalPos) const
Evaluate the boundary conditions for a neumann boundary segment.
Definition: onemodelproblem.hh:241
Variables & variables()
Returns variables object.
Definition: onemodelproblem.hh:563
Declares all properties used in Dumux.
Base file for properties related to sequential models.