3.2-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
impetproblem.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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program 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 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program 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 this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
19#ifndef DUMUX_IMPETPROBLEM_HH
20#define DUMUX_IMPETPROBLEM_HH
21
22#include <dune/common/float_cmp.hh>
23#include "impetproperties.hh"
25#include <dumux/io/restart.hh>
27
29
35namespace Dumux
36{
44template<class TypeTag>
46{
47private:
48 using Implementation = GetPropType<TypeTag, Properties::Problem>;
51
53
55
57
59 using VertexMapper = typename SolutionTypes::VertexMapper;
60 using ElementMapper = typename SolutionTypes::ElementMapper;
61
66
68
69 enum
70 {
71 dim = GridView::dimension,
72 dimWorld = GridView::dimensionworld
73 };
74 enum
75 {
76 wetting = 0, nonwetting = 1,
77 adaptiveGrid = getPropValue<TypeTag, Properties::AdaptiveGrid>()
78 };
79
80 using GlobalPosition = Dune::FieldVector<Scalar, dimWorld>;
81 using Element = typename GridView::Traits::template Codim<0>::Entity;
82 using Intersection = typename GridView::Intersection;
83 // The module to adapt grid. If adaptiveGrid is false, this model does nothing.
85
86 using PrimaryVariables = typename SolutionTypes::PrimaryVariables;
88
89 //private!! copy constructor
91 {}
92
93public:
98 IMPETProblem(TimeManager &timeManager, Grid& grid)
99 : IMPETProblem(timeManager, grid, grid.leafGridView())
100 {}
101
102 IMPETProblem(TimeManager &timeManager, Grid& grid, const GridView& gridView)
103 : grid_(&grid),
104 gridView_(gridView),
105 bBoxMin_(std::numeric_limits<double>::max()),
106 bBoxMax_(-std::numeric_limits<double>::max()),
107 timeManager_(&timeManager),
108 outputInterval_(1),
109 outputTimeInterval_(0.0),
110 vtkOutputLevel_(-1)
111 {
112 // calculate the bounding box of the grid view
113 using std::min;
114 using std::max;
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]);
119 }
120 }
121
122 // communicate to get the bounding box of the whole domain
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]);
127 }
128
129 variables_ = std::make_shared<Variables>(this->gridView());
130 pressModel_ = std::make_shared<PressureModel>(asImp_());
131
132 transportModel_ = std::make_shared<TransportModel>(asImp_());
133 model_ = std::make_shared<IMPETModel>(asImp_()) ;
134
135 // create an Object to handle adaptive grids
136 if (adaptiveGrid)
137 gridAdapt_ = std::make_shared<GridAdaptModel>(asImp_());
138
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());
142 }
143
151 void boundaryTypes(BoundaryTypes &bcTypes,
152 const Intersection &intersection) const
153 {
154 // forward it to the method which only takes the global coordinate
155 asImp_().boundaryTypesAtPos(bcTypes, intersection.geometry().center());
156 }
157
165 void boundaryTypesAtPos(BoundaryTypes &bcTypes,
166 const GlobalPosition &globalPos) const
167 {
168 // Throw an exception (there is no reasonable default value
169 // for Dirichlet conditions)
170 DUNE_THROW(Dune::InvalidStateException,
171 "The problem does not provide "
172 "a boundaryTypesAtPos() method.");
173 }
174
184 void dirichlet(PrimaryVariables &values,
185 const Intersection &intersection) const
186 {
187 // forward it to the method which only takes the global coordinate
188 asImp_().dirichletAtPos(values, intersection.geometry().center());
189 }
190
200 void dirichletAtPos(PrimaryVariables &values,
201 const GlobalPosition &globalPos) const
202 {
203 // Throw an exception (there is no reasonable default value
204 // for Dirichlet conditions)
205 DUNE_THROW(Dune::InvalidStateException,
206 "The problem specifies that some boundary "
207 "segments are dirichlet, but does not provide "
208 "a dirichletAtPos() method.");
209 }
210
221 void neumann(PrimaryVariables &values,
222 const Intersection &intersection) const
223 {
224 // forward it to the interface with only the global position
225 asImp_().neumannAtPos(values, intersection.geometry().center());
226 }
227
238 void neumannAtPos(PrimaryVariables &values,
239 const GlobalPosition &globalPos) const
240 {
241 // Throw an exception (there is no reasonable default value
242 // for Neumann conditions)
243 DUNE_THROW(Dune::InvalidStateException,
244 "The problem specifies that some boundary "
245 "segments are neumann, but does not provide "
246 "a neumannAtPos() method.");
247 }
248
259 void source(PrimaryVariables &values,
260 const Element &element) const
261 {
262 // forward to generic interface
263 asImp_().sourceAtPos(values, element.geometry().center());
264 }
265
279 void sourceAtPos(PrimaryVariables &values,
280 const GlobalPosition &globalPos) const
281 { // Throw an exception (there is no initial condition)
282 DUNE_THROW(Dune::InvalidStateException,
283 "The problem does not provide "
284 "a sourceAtPos() method.");
285 }
286
296 void initial(PrimaryVariables &values,
297 const Element &element) const
298 {
299 // forward to generic interface
300 asImp_().initialAtPos(values, element.geometry().center());
301 }
302
313 void initialAtPos(PrimaryVariables &values,
314 const GlobalPosition &globalPos) const
315 {
316 // initialize with 0 by default
317 values = 0;
318 }
319
324 void init()
325 {
326 // set the initial condition of the model
327 variables_->initialize();
328 model().initialize();
329 if (adaptiveGrid)
330 gridAdapt().init();
331 }
332
338 {
339 // if adaptivity is used, this method adapts the grid.
340 // if it is not used, this method does nothing.
341 if (adaptiveGrid && timeManager().timeStepIndex() > 0)
342 this->gridAdapt().adaptGrid();
343 asImp_().pressureModel().updateMaterialLaws();
344 }
345
358 {
359 // allocate temporary vectors for the updates
360 TransportSolutionType updateVector;//(this->transportModel().solutionSize());
361
362 Scalar t = timeManager().time();
363 Scalar dt = std::numeric_limits<Scalar>::max();
364
365 // obtain the first update and the time step size
366 model().update(t, dt, updateVector);
367
368 using std::max;
369 using std::min;
370 if (t > 0 && timeManager().timeStepSize()> 0.)
371 {
372 Scalar oldDt = timeManager().timeStepSize();
373 Scalar minDt = max(oldDt - dtVariationRestrictionFactor_ * oldDt, 0.0);
374 Scalar maxDt = oldDt + dtVariationRestrictionFactor_ * oldDt;
375 dt = max(min(dt, maxDt), minDt);
376 }
377
378 //make sure t_old + dt is not larger than tend
379 dt = min(dt, timeManager().episodeMaxTimeStepSize());
380
381 // check if we are in first TS and an initialDt was assigned
382 if (Dune::FloatCmp::eq<Scalar, Dune::FloatCmp::absolute>(t, 0.0, 1.0e-30)
383 && Dune::FloatCmp::ne<Scalar, Dune::FloatCmp::absolute>(timeManager().timeStepSize(), 0.0, 1.0e-30))
384 {
385 if (gridView().comm().size() > 1)
386 dt = gridView().comm().min(dt);
387
388 // check if assigned initialDt is in accordance with dt from first transport step
389 if (timeManager().timeStepSize() > dt
390 && gridView().comm().rank() == 0)
391 Dune::dwarn << "Initial timestep of size " << timeManager().timeStepSize()
392 << " is larger then dt=" << dt <<" from transport" << std::endl;
393 // internally assign next timestep size
394 dt = min(dt, timeManager().timeStepSize());
395 }
396
397 //make sure the right time-step is used by all processes in the parallel case
398 if (gridView().comm().size() > 1)
399 dt = gridView().comm().min(dt);
400
401 // check maximum allowed time step size
402 timeManager().setTimeStepSize(dt);
403
404 // explicit Euler: Sat <- Sat + dt*N(Sat)
405 transportModel().updateTransportedQuantity(updateVector);
406 }
407
417 {}
418
425 {}
426
430 Scalar timeStepSize() const
431 { return timeManager().timeStepSize(); }
432
436 void setTimeStepSize(const Scalar dt)
437 { timeManager().setTimeStepSize(dt); }
438
444 Scalar nextTimeStepSize(const Scalar dt)
445 { return timeManager().timeStepSize();}
446
452 Scalar maxTimeStepSize() const
453 { return maxTimeStepSize_; }
454
464 {
465 if (outputInterval_ > 0)
466 {
467 return
468 timeManager().timeStepIndex() > 0 &&
469 (timeManager().timeStepIndex() % int(100*outputInterval_) == 0);
470 }
471 else
472 return
474 }
475
481 void setOutputTimeInterval(const Scalar timeInterval)
482 {
483 outputTimeInterval_ = (timeInterval > 0.0) ? timeInterval : 1e100;
484 timeManager().startNextEpisode(outputTimeInterval_);
485 }
486
492 void setOutputInterval(const int interval)
493 {
494 using std::max;
495 outputInterval_ = max(interval, 0);
496 }
497
506 bool shouldWriteOutput() const
507 {
508 if (outputInterval_ > 0)
509 {
510 if (timeManager().timeStepIndex() % outputInterval_ == 0
511 || timeManager().willBeFinished()
512 || timeManager().episodeWillBeFinished())
513 {
514 return true;
515 }
516 }
517 else if (timeManager().willBeFinished()
518 || timeManager().episodeWillBeFinished() || timeManager().timeStepIndex() == 0)
519 {
520 return true;
521 }
522 return false;
523 }
524
529 {
530 if (outputTimeInterval_ > 0.0 && !timeManager().finished())
531 {
532 timeManager().startNextEpisode(outputTimeInterval_);
533 }
534 else if (!timeManager().finished())
535 {
536 std::cerr << "The end of an episode is reached, but the problem "
537 << "does not override the episodeEnd() method. "
538 << "Doing nothing!\n";
539 }
540 }
541
542 // \}
543
551 const std::string& name() const
552 {
553 return simname_;
554 }
555
565 void setName(std::string newName)
566 {
567 simname_ = newName;
568 }
569
573 const GridView& gridView() const
574 {
575 return gridView_;
576 }
577
581 Grid &grid()
582 {
583 if (!grid_)
584 {
585 DUNE_THROW(Dune::InvalidStateException, "Grid was called in problemclass, "
586 << "although it is not specified. Do so by using setGrid() method!");
587 }
588 return *grid_;
589 }
594 void setGrid(Grid &grid)
595 {
596 grid_ = &grid;
597 }
598
602 GridAdaptModel& gridAdapt()
603 {
604 return *gridAdapt_;
605 }
606
614 void preAdapt()
615 {
616 }
617
626 {
627 // write out new grid
628 // Dune::VTKWriter<LeafGridView> vtkwriter(problem_.gridView());
629 // vtkwriter.write("latestgrid",Dune::VTK::appendedraw);
630 }
631
635 const VertexMapper &vertexMapper() const
636 { return variables_->vertexMapper(); }
637
641 const ElementMapper &elementMapper() const
642 { return variables_->elementMapper(); }
643
648 const GlobalPosition &bBoxMin() const
649 { return bBoxMin_; }
650
655 const GlobalPosition &bBoxMax() const
656 { return bBoxMax_; }
657
659
660
663 TimeManager &timeManager()
664 { return *timeManager_; }
665
667 const TimeManager &timeManager() const
668 { return *timeManager_; }
669
676 Variables& variables ()
677 { return *variables_; }
678
680 const Variables& variables () const
681 { return *variables_; }
682
686 IMPETModel &model()
687 { return *model_; }
688
690 const IMPETModel &model() const
691 { return *model_; }
692
696 PressureModel &pressureModel()
697 { return *pressModel_; }
698
700 const PressureModel &pressureModel() const
701 { return *pressModel_; }
702
706 TransportModel &transportModel()
707 { return *transportModel_; }
708
710 const TransportModel &transportModel() const
711 { return *transportModel_; }
712 // \}
713
717 // \{
718
729 {
730 using Restarter = Restart;
731
732 Restarter res;
733 res.serializeBegin(asImp_());
734 std::cout << "Serialize to file " << res.fileName() << "\n";
735
736 timeManager().serialize(res);
737 resultWriter().serialize(res);
738
739 // do the actual serialization process: write primary variables
740 res.template serializeEntities<0> (*pressModel_, gridView());
741 res.template serializeEntities<0> (*transportModel_, gridView());
742
743 res.serializeEnd();
744
745 if (adaptiveGrid)
746 {
748 }
749 }
750
758 void restart(const double tRestart)
759 {
760 if (adaptiveGrid)
761 {
763 variables().initialize();
764 model().initialize();
765 }
766
767 using Restarter = Restart;
768
769 Restarter res;
770 res.deserializeBegin(asImp_(), tRestart);
771 std::cout << "Deserialize from file " << res.fileName() << "\n";
772
773 timeManager().deserialize(res);
774
776
777 // do the actual serialization process: get primary variables
778 res.template deserializeEntities<0> (*pressModel_, gridView());
779 res.template deserializeEntities<0> (*transportModel_, gridView());
780 pressureModel().updateMaterialLaws();
781
782 res.deserializeEnd();
783 }
784
786 {
787 }
793 int vtkOutputLevel() const
794 {
795 return vtkOutputLevel_;
796 }
797
799 void writeOutput(bool verbose = true)
800 {
801 if (verbose && gridView().comm().rank() == 0)
802 std::cout << "Writing result file for current time step\n";
803
804 if (!resultWriter_)
805 resultWriter_ = std::make_shared<VtkMultiWriter>(gridView(), asImp_().name());
806 if (adaptiveGrid)
807 resultWriter_->gridChanged();
808 resultWriter_->beginWrite(timeManager().time() + timeManager().timeStepSize());
809 model().addOutputVtkFields(*resultWriter_);
810 asImp_().addOutputVtkFields();
811 resultWriter_->endWrite();
812 }
813
814 // \}
815
816protected:
819 {
820 if (!resultWriter_)
821 resultWriter_ = std::make_shared<VtkMultiWriter>(gridView(), asImp_().name());
822 return *resultWriter_;
823 }
826 {
827 if (!resultWriter_)
828 resultWriter_ = std::make_shared<VtkMultiWriter>(gridView(), asImp_().name());
829 return *resultWriter_;
830 }
831
832private:
834 Implementation &asImp_()
835 { return *static_cast<Implementation *>(this); }
836
838 const Implementation &asImp_() const
839 { return *static_cast<const Implementation *>(this); }
840
841 std::string simname_; // a string for the name of the current simulation,
842 // which could be set by means of a program argument, for example.
843
844 // pointer to a possibly adaptive grid.
845 Grid *grid_;
846 GridView gridView_;
847
848 GlobalPosition bBoxMin_;
849 GlobalPosition bBoxMax_;
850
851 TimeManager *timeManager_;
852 Scalar maxTimeStepSize_;
853
854 std::shared_ptr<Variables> variables_;
855
856 std::shared_ptr<PressureModel> pressModel_;
857 std::shared_ptr<TransportModel> transportModel_;
858 std::shared_ptr<IMPETModel> model_;
859
860 std::shared_ptr<VtkMultiWriter> resultWriter_;
861 int outputInterval_;
862 Scalar outputTimeInterval_;
863 int vtkOutputLevel_;
864 std::shared_ptr<GridAdaptModel> gridAdapt_;
865
866 Scalar dtVariationRestrictionFactor_;
867};
868}
869#endif
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.
Definition: adapt.hh:29
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