19#ifndef DUMUX_IMPETPROBLEM_HH
20#define DUMUX_IMPETPROBLEM_HH
22#include <dune/common/float_cmp.hh>
44template<
class TypeTag>
59 using VertexMapper =
typename SolutionTypes::VertexMapper;
60 using ElementMapper =
typename SolutionTypes::ElementMapper;
71 dim = GridView::dimension,
72 dimWorld = GridView::dimensionworld
76 wetting = 0, nonwetting = 1,
77 adaptiveGrid = getPropValue<TypeTag, Properties::AdaptiveGrid>()
80 using GlobalPosition = Dune::FieldVector<Scalar, dimWorld>;
81 using Element =
typename GridView::Traits::template Codim<0>::Entity;
82 using Intersection =
typename GridView::Intersection;
86 using PrimaryVariables =
typename SolutionTypes::PrimaryVariables;
105 bBoxMin_(std::numeric_limits<double>::max()),
106 bBoxMax_(-std::numeric_limits<double>::max()),
109 outputTimeInterval_(0.0),
115 for (
const auto& vertex : vertices(this->
gridView())) {
116 for (
int i=0; i<dim; i++) {
117 bBoxMin_[i] = min(bBoxMin_[i], vertex.geometry().center()[i]);
118 bBoxMax_[i] = max(bBoxMax_[i], vertex.geometry().center()[i]);
123 if (this->
gridView().comm().size() > 1)
124 for (
int i = 0; i < dim; ++i) {
125 bBoxMin_[i] = this->
gridView().comm().min(bBoxMin_[i]);
126 bBoxMax_[i] = this->
gridView().comm().max(bBoxMax_[i]);
129 variables_ = std::make_shared<Variables>(this->
gridView());
130 pressModel_ = std::make_shared<PressureModel>(asImp_());
132 transportModel_ = std::make_shared<TransportModel>(asImp_());
133 model_ = std::make_shared<IMPETModel>(asImp_()) ;
137 gridAdapt_ = std::make_shared<GridAdaptModel>(asImp_());
139 vtkOutputLevel_ = getParam<int>(
"Vtk.OutputLevel");
140 dtVariationRestrictionFactor_ = getParam<Scalar>(
"Impet.DtVariationRestrictionFactor", std::numeric_limits<Scalar>::max());
141 maxTimeStepSize_ = getParam<Scalar>(
"TimeManager.MaxTimeStepSize", std::numeric_limits<Scalar>::max());
152 const Intersection &intersection)
const
155 asImp_().boundaryTypesAtPos(bcTypes, intersection.geometry().center());
166 const GlobalPosition &globalPos)
const
170 DUNE_THROW(Dune::InvalidStateException,
171 "The problem does not provide "
172 "a boundaryTypesAtPos() method.");
185 const Intersection &intersection)
const
188 asImp_().dirichletAtPos(values, intersection.geometry().center());
201 const GlobalPosition &globalPos)
const
205 DUNE_THROW(Dune::InvalidStateException,
206 "The problem specifies that some boundary "
207 "segments are dirichlet, but does not provide "
208 "a dirichletAtPos() method.");
222 const Intersection &intersection)
const
225 asImp_().neumannAtPos(values, intersection.geometry().center());
239 const GlobalPosition &globalPos)
const
243 DUNE_THROW(Dune::InvalidStateException,
244 "The problem specifies that some boundary "
245 "segments are neumann, but does not provide "
246 "a neumannAtPos() method.");
260 const Element &element)
const
263 asImp_().sourceAtPos(values, element.geometry().center());
280 const GlobalPosition &globalPos)
const
282 DUNE_THROW(Dune::InvalidStateException,
283 "The problem does not provide "
284 "a sourceAtPos() method.");
297 const Element &element)
const
300 asImp_().initialAtPos(values, element.geometry().center());
314 const GlobalPosition &globalPos)
const
327 variables_->initialize();
328 model().initialize();
341 if (adaptiveGrid &&
timeManager().timeStepIndex() > 0)
343 asImp_().pressureModel().updateMaterialLaws();
360 TransportSolutionType updateVector;
363 Scalar dt = std::numeric_limits<Scalar>::max();
366 model().update(t, dt, updateVector);
373 Scalar minDt = max(oldDt - dtVariationRestrictionFactor_ * oldDt, 0.0);
374 Scalar maxDt = oldDt + dtVariationRestrictionFactor_ * oldDt;
375 dt = max(min(dt, maxDt), minDt);
379 dt = min(dt,
timeManager().episodeMaxTimeStepSize());
382 if (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(t, 0.0, 1.0e-30)
391 Dune::dwarn <<
"Initial timestep of size " <<
timeManager().timeStepSize()
392 <<
" is larger then dt=" << dt <<
" from transport" << std::endl;
453 {
return maxTimeStepSize_; }
465 if (outputInterval_ > 0)
469 (
timeManager().timeStepIndex() % int(100*outputInterval_) == 0);
483 outputTimeInterval_ = (timeInterval > 0.0) ? timeInterval : 1e100;
484 timeManager().startNextEpisode(outputTimeInterval_);
495 outputInterval_ = max(interval, 0);
508 if (outputInterval_ > 0)
510 if (
timeManager().timeStepIndex() % outputInterval_ == 0
530 if (outputTimeInterval_ > 0.0 && !
timeManager().finished())
532 timeManager().startNextEpisode(outputTimeInterval_);
536 std::cerr <<
"The end of an episode is reached, but the problem "
537 <<
"does not override the episodeEnd() method. "
538 <<
"Doing nothing!\n";
551 const std::string&
name()
const
585 DUNE_THROW(Dune::InvalidStateException,
"Grid was called in problemclass, "
586 <<
"although it is not specified. Do so by using setGrid() method!");
636 {
return variables_->vertexMapper(); }
642 {
return variables_->elementMapper(); }
664 {
return *timeManager_; }
668 {
return *timeManager_; }
677 {
return *variables_; }
681 {
return *variables_; }
697 {
return *pressModel_; }
701 {
return *pressModel_; }
707 {
return *transportModel_; }
711 {
return *transportModel_; }
734 std::cout <<
"Serialize to file " << res.fileName() <<
"\n";
740 res.template serializeEntities<0> (*pressModel_,
gridView());
741 res.template serializeEntities<0> (*transportModel_,
gridView());
764 model().initialize();
771 std::cout <<
"Deserialize from file " << res.fileName() <<
"\n";
778 res.template deserializeEntities<0> (*pressModel_,
gridView());
779 res.template deserializeEntities<0> (*transportModel_,
gridView());
782 res.deserializeEnd();
795 return vtkOutputLevel_;
801 if (verbose &&
gridView().comm().rank() == 0)
802 std::cout <<
"Writing result file for current time step\n";
805 resultWriter_ = std::make_shared<VtkMultiWriter>(
gridView(), asImp_().
name());
807 resultWriter_->gridChanged();
809 model().addOutputVtkFields(*resultWriter_);
810 asImp_().addOutputVtkFields();
811 resultWriter_->endWrite();
821 resultWriter_ = std::make_shared<VtkMultiWriter>(
gridView(), asImp_().
name());
822 return *resultWriter_;
828 resultWriter_ = std::make_shared<VtkMultiWriter>(
gridView(), asImp_().
name());
829 return *resultWriter_;
834 Implementation &asImp_()
835 {
return *
static_cast<Implementation *
>(
this); }
838 const Implementation &asImp_()
const
839 {
return *
static_cast<const Implementation *
>(
this); }
841 std::string simname_;
848 GlobalPosition bBoxMin_;
849 GlobalPosition bBoxMax_;
851 TimeManager *timeManager_;
852 Scalar maxTimeStepSize_;
854 std::shared_ptr<Variables> variables_;
856 std::shared_ptr<PressureModel> pressModel_;
857 std::shared_ptr<TransportModel> transportModel_;
858 std::shared_ptr<IMPETModel> model_;
860 std::shared_ptr<VtkMultiWriter> resultWriter_;
862 Scalar outputTimeInterval_;
864 std::shared_ptr<GridAdaptModel> gridAdapt_;
866 Scalar dtVariationRestrictionFactor_;
Provides a restart functionality for adaptive grids.
Generic class to read/write restart files.
Simplifies writing multi-file VTK datasets.
Base class for h-adaptive sequential models.
Base file for properties related to sequential IMPET algorithms.
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
static void restartGrid(Problem &problem)
Restart the grid from the file.
Definition: adaptivegridrestart.hh:83
static void serializeGrid(Problem &problem)
Write the grid to a file.
Definition: adaptivegridrestart.hh:73
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 problems using a sequential implicit-explicit strategy
Definition: impetproblem.hh:46
const GridView & gridView() const
The GridView which used by the problem.
Definition: impetproblem.hh:573
GridAdaptModel & gridAdapt()
Returns adaptivity model used for the problem.
Definition: impetproblem.hh:602
void initial(PrimaryVariables &values, const Element &element) const
Evaluate the initial value for a control volume.
Definition: impetproblem.hh:296
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: impetproblem.hh:165
TimeManager & timeManager()
Returns TimeManager object used by the simulation.
Definition: impetproblem.hh:663
void setTimeStepSize(const Scalar dt)
Sets the current time step size [seconds].
Definition: impetproblem.hh:436
void serialize()
This method writes the complete state of the problem to the harddisk.
Definition: impetproblem.hh:728
void preAdapt()
Capability to introduce problem-specific routines at the beginning of the grid adaptation.
Definition: impetproblem.hh:614
void neumannAtPos(PrimaryVariables &values, const GlobalPosition &globalPos) const
Evaluate the boundary conditions for a neumann boundary segment.
Definition: impetproblem.hh:238
const PressureModel & pressureModel() const
Returns the pressure model used for the problem.
Definition: impetproblem.hh:700
void setOutputTimeInterval(const Scalar timeInterval)
Sets a time interval for Output.
Definition: impetproblem.hh:481
void source(PrimaryVariables &values, const Element &element) const
Evaluate the source term.
Definition: impetproblem.hh:259
PressureModel & pressureModel()
Returns the pressure model used for the problem.
Definition: impetproblem.hh:696
void init()
Called by the TimeManager in order to initialize the problem.
Definition: impetproblem.hh:324
int vtkOutputLevel() const
Returns the vtk output verbosity level.
Definition: impetproblem.hh:793
Scalar timeStepSize() const
Returns the current time step size [seconds].
Definition: impetproblem.hh:430
const GlobalPosition & bBoxMin() const
The coordinate of the corner of the GridView's bounding box with the smallest values.
Definition: impetproblem.hh:648
void sourceAtPos(PrimaryVariables &values, const GlobalPosition &globalPos) const
Evaluate the source term for all phases within a given sub-control-volume.
Definition: impetproblem.hh:279
void setName(std::string newName)
Set the problem name.
Definition: impetproblem.hh:565
void postAdapt()
Capability to introduce problem-specific routines after grid adaptation.
Definition: impetproblem.hh:625
VtkMultiWriter & resultWriter() const
Returns the applied VTK-writer for the output.
Definition: impetproblem.hh:825
IMPETProblem(TimeManager &timeManager, Grid &grid)
Constructs an object of type IMPETProblemProblem.
Definition: impetproblem.hh:98
void setGrid(Grid &grid)
Specifies the grid from outside the problem.
Definition: impetproblem.hh:594
const TransportModel & transportModel() const
Returns transport model used for the problem.
Definition: impetproblem.hh:710
void addOutputVtkFields()
Definition: impetproblem.hh:785
void setOutputInterval(const int interval)
Sets the interval for Output.
Definition: impetproblem.hh:492
TransportModel & transportModel()
Returns transport model used for the problem.
Definition: impetproblem.hh:706
Scalar maxTimeStepSize() const
Returns the maximum allowed time step size [s].
Definition: impetproblem.hh:452
IMPETProblem(TimeManager &timeManager, Grid &grid, const GridView &gridView)
Definition: impetproblem.hh:102
VtkMultiWriter & resultWriter()
Returns the applied VTK-writer for the output.
Definition: impetproblem.hh:818
const TimeManager & timeManager() const
Returns TimeManager object used by the simulation.
Definition: impetproblem.hh:667
const Variables & variables() const
Returns variables container. ()
Definition: impetproblem.hh:680
void episodeEnd()
Called when the end of an simulation episode is reached.
Definition: impetproblem.hh:528
const VertexMapper & vertexMapper() const
Returns the mapper for vertices to indices.
Definition: impetproblem.hh:635
void preTimeStep()
Called by TimeManager just before the time integration.
Definition: impetproblem.hh:337
bool shouldWriteRestartFile() const
Returns true if a restart file should be written to disk.
Definition: impetproblem.hh:463
void writeOutput(bool verbose=true)
Write the fields current solution into an VTK output file.
Definition: impetproblem.hh:799
bool shouldWriteOutput() const
Returns true if the current solution should be written to disk (i.e. as a VTK file)
Definition: impetproblem.hh:506
void initialAtPos(PrimaryVariables &values, const GlobalPosition &globalPos) const
Evaluate the initial value for a control volume.
Definition: impetproblem.hh:313
const std::string & name() const
The problem name.
Definition: impetproblem.hh:551
Scalar nextTimeStepSize(const Scalar dt)
Called by TimeManager whenever a solution for a timestep has been computed and the simulation time ha...
Definition: impetproblem.hh:444
IMPETModel & model()
Returns numerical model used for the problem.
Definition: impetproblem.hh:686
void advanceTimeLevel()
Called by the time manager after everything which can be done about the current time step is finished...
Definition: impetproblem.hh:424
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: impetproblem.hh:151
const ElementMapper & elementMapper() const
Returns the mapper for elements to indices.
Definition: impetproblem.hh:641
const GlobalPosition & bBoxMax() const
The coordinate of the corner of the GridView's bounding box with the largest values.
Definition: impetproblem.hh:655
Variables & variables()
Returns variables container.
Definition: impetproblem.hh:676
void neumann(PrimaryVariables &values, const Intersection &intersection) const
Evaluate the boundary conditions for a neumann boundary segment.
Definition: impetproblem.hh:221
void postTimeStep()
Called by TimeManager whenever a solution for a timestep has been computed and the simulation time ha...
Definition: impetproblem.hh:416
void dirichletAtPos(PrimaryVariables &values, const GlobalPosition &globalPos) const
Evaluate the boundary conditions for a dirichlet control volume.
Definition: impetproblem.hh:200
void timeIntegration()
Called by TimeManager in order to do a time integration on the model.
Definition: impetproblem.hh:357
void restart(const double tRestart)
This method restores the complete state of the problem from disk.
Definition: impetproblem.hh:758
const IMPETModel & model() const
Returns numerical model used for the problem.
Definition: impetproblem.hh:690
Grid & grid()
Returns the current grid which used by the problem.
Definition: impetproblem.hh:581
void dirichlet(PrimaryVariables &values, const Intersection &intersection) const
Evaluate the boundary conditions for a dirichlet control volume.
Definition: impetproblem.hh:184