Skip to content
Snippets Groups Projects
Commit 56ec173c authored by Theresa Schollenberger's avatar Theresa Schollenberger
Browse files

[execise-basic] update exercise and readme to dumux 3.0

parent 55bd53dc
No related branches found
No related tags found
1 merge request!54update dumux-course to dumux 3.0
......@@ -60,7 +60,7 @@ paraview injection-2p2c.pvd
* Copy the main file `exercise_basic_2p.cc` and rename it to `exercise_basic_2pni.cc`
* In `exercise_basic_2pni.cc`, include the header `injection2pniproblem.hh` instead of `injection2pproblem.hh`.
* In `exercise_basic_2pni.cc`, change `Injection2pCCTypeTag` to `Injection2pNICCTypeTag` in the line `using TypeTag = TTAG(Injection2pCCTypeTag);`
* In `exercise_basic_2pni.cc`, change `Injection2pCCTypeTag` to `Injection2pNICCTypeTag` in the line `using TypeTag = Properties::TTag::Injection2pNICCTypeTag;`
* Add a new executable in `CMakeLists.txt` by adding the lines
```cmake
......
......@@ -41,10 +41,11 @@
#include <dumux/assembly/fvassembler.hh>
#include <dumux/assembly/diffmethod.hh>
#include <dumux/discretization/methods.hh>
#include <dumux/discretization/method.hh>
#include <dumux/io/vtkoutputmodule.hh>
#include <dumux/io/grid/gridmanager.hh>
#include <dumux/io/loadsolution.hh>
// The problem file, where setup-specific boundary and initial conditions are defined.
#include "injection2pproblem.hh"
......@@ -57,7 +58,7 @@ int main(int argc, char** argv) try
using namespace Dumux;
// define the type tag for this problem
using TypeTag = TTAG(Injection2pCCTypeTag);
using TypeTag = Properties::TTag::Injection2pCCTypeTag;
// initialize MPI, finalize is done automatically on exit
const auto& mpiHelper = Dune::MPIHelper::instance(argc, argv);
......@@ -70,7 +71,7 @@ int main(int argc, char** argv) try
Parameters::init(argc, argv);
// try to create a grid (from the given grid file or the input file)
GridManager<typename GET_PROP_TYPE(TypeTag, Grid)> gridManager;
GridManager<GetPropType<TypeTag, Properties::Grid>> gridManager;
gridManager.init();
////////////////////////////////////////////////////////////
......@@ -81,35 +82,57 @@ int main(int argc, char** argv) try
const auto& leafGridView = gridManager.grid().leafGridView();
// create the finite volume grid geometry
using FVGridGeometry = typename GET_PROP_TYPE(TypeTag, FVGridGeometry);
using FVGridGeometry = GetPropType<TypeTag, Properties::FVGridGeometry>;
auto fvGridGeometry = std::make_shared<FVGridGeometry>(leafGridView);
fvGridGeometry->update();
// the problem (initial and boundary conditions)
using Problem = typename GET_PROP_TYPE(TypeTag, Problem);
using Problem = GetPropType<TypeTag, Properties::Problem>;
auto problem = std::make_shared<Problem>(fvGridGeometry);
// check if we are about to restart a previously interrupted simulation
using Scalar = GetPropType<TypeTag, Properties::Scalar>;
Scalar restartTime = getParam<Scalar>("Restart.Time", 0);
// the solution vector
using SolutionVector = typename GET_PROP_TYPE(TypeTag, SolutionVector);
using SolutionVector = GetPropType<TypeTag, Properties::SolutionVector>;
SolutionVector x(fvGridGeometry->numDofs());
problem->applyInitialSolution(x);
if (restartTime > 0)
{
using IOFields = GetPropType<TypeTag, Properties::IOFields>;
using PrimaryVariables = GetPropType<TypeTag, Properties::PrimaryVariables>;
using ModelTraits = GetPropType<TypeTag, Properties::ModelTraits>;
using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
const auto fileName = getParam<std::string>("Restart.File");
const auto pvName = createPVNameFunction<IOFields, PrimaryVariables, ModelTraits, FluidSystem>();
loadSolution(x, fileName, pvName, *fvGridGeometry);
}
else
problem->applyInitialSolution(x);
auto xOld = x;
// the grid variables
using GridVariables = typename GET_PROP_TYPE(TypeTag, GridVariables);
using GridVariables = GetPropType<TypeTag, Properties::GridVariables>;
auto gridVariables = std::make_shared<GridVariables>(problem, fvGridGeometry);
gridVariables->init(x, xOld);
gridVariables->init(x);
// get some time loop parameters
using Scalar = typename GET_PROP_TYPE(TypeTag, Scalar);
using Scalar = GetPropType<TypeTag, Properties::Scalar>;
const auto tEnd = getParam<Scalar>("TimeLoop.TEnd");
const auto maxDt = getParam<Scalar>("TimeLoop.MaxTimeStepSize");
auto dt = getParam<Scalar>("TimeLoop.DtInitial");
// intialize the vtk output module
using VtkOutputFields = typename GET_PROP_TYPE(TypeTag, VtkOutputFields);
VtkOutputModule<TypeTag> vtkWriter(*problem, *fvGridGeometry, *gridVariables, x, problem->name());
VtkOutputFields::init(vtkWriter); //! Add model specific output fields
using IOFields = GetPropType<TypeTag, Properties::IOFields>;
// use non-conforming output for the test with interface solver
const auto ncOutput = getParam<bool>("Problem.UseNonConformingOutput", false);
VtkOutputModule<GridVariables, SolutionVector> vtkWriter(*gridVariables, x, problem->name(), "",
ncOutput ? Dune::VTK::nonconforming : Dune::VTK::conforming);
using VelocityOutput = GetPropType<TypeTag, Properties::VelocityOutput>;
vtkWriter.addVelocityOutput(std::make_shared<VelocityOutput>(*gridVariables));
IOFields::initOutputModule(vtkWriter); //!< Add model specific output fields
vtkWriter.write(restartTime);
// instantiate time loop
auto timeLoop = std::make_shared<TimeLoop<Scalar>>(0.0, dt, tEnd);
......
......@@ -36,15 +36,16 @@
#include <dumux/common/dumuxmessage.hh>
#include <dumux/linear/amgbackend.hh>
#include <dumux/nonlinear/privarswitchnewtonsolver.hh>
#include <dumux/nonlinear/newtonsolver.hh>
#include <dumux/assembly/fvassembler.hh>
#include <dumux/assembly/diffmethod.hh>
#include <dumux/discretization/methods.hh>
#include <dumux/discretization/method.hh>
#include <dumux/io/vtkoutputmodule.hh>
#include <dumux/io/grid/gridmanager.hh>
#include <dumux/io/loadsolution.hh>
// The problem file, where setup-specific boundary and initial conditions are defined.
#include "injection2p2cproblem.hh"
......@@ -57,7 +58,7 @@ int main(int argc, char** argv) try
using namespace Dumux;
// define the type tag for this problem
using TypeTag = TTAG(Injection2p2cCCTypeTag);
using TypeTag = Properties::TTag::Injection2p2cCCTypeTag;
// initialize MPI, finalize is done automatically on exit
const auto& mpiHelper = Dune::MPIHelper::instance(argc, argv);
......@@ -70,7 +71,7 @@ int main(int argc, char** argv) try
Parameters::init(argc, argv);
// try to create a grid (from the given grid file or the input file)
GridManager<typename GET_PROP_TYPE(TypeTag, Grid)> gridManager;
GridManager<GetPropType<TypeTag, Properties::Grid>> gridManager;
gridManager.init();
////////////////////////////////////////////////////////////
......@@ -81,38 +82,62 @@ int main(int argc, char** argv) try
const auto& leafGridView = gridManager.grid().leafGridView();
// create the finite volume grid geometry
using FVGridGeometry = typename GET_PROP_TYPE(TypeTag, FVGridGeometry);
using FVGridGeometry = GetPropType<TypeTag, Properties::FVGridGeometry>;
auto fvGridGeometry = std::make_shared<FVGridGeometry>(leafGridView);
fvGridGeometry->update();
// the problem (initial and boundary conditions)
using Problem = typename GET_PROP_TYPE(TypeTag, Problem);
using Problem = GetPropType<TypeTag, Properties::Problem>;
auto problem = std::make_shared<Problem>(fvGridGeometry);
// check if we are about to restart a previously interrupted simulation
using Scalar = GetPropType<TypeTag, Properties::Scalar>;
Scalar restartTime = getParam<Scalar>("Restart.Time", 0);
// the solution vector
using SolutionVector = typename GET_PROP_TYPE(TypeTag, SolutionVector);
using SolutionVector = GetPropType<TypeTag, Properties::SolutionVector>;
SolutionVector x(fvGridGeometry->numDofs());
problem->applyInitialSolution(x);
// problem->applyInitialSolution(x);
if (restartTime > 0)
{
using IOFields = GetPropType<TypeTag, Properties::IOFields>;
using PrimaryVariables = GetPropType<TypeTag, Properties::PrimaryVariables>;
using ModelTraits = GetPropType<TypeTag, Properties::ModelTraits>;
using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
const auto fileName = getParam<std::string>("Restart.File");
const auto pvName = createPVNameFunction<IOFields, PrimaryVariables, ModelTraits, FluidSystem>();
loadSolution(x, fileName, pvName, *fvGridGeometry);
}
else
problem->applyInitialSolution(x);
auto xOld = x;
// the grid variables
using GridVariables = typename GET_PROP_TYPE(TypeTag, GridVariables);
using GridVariables = GetPropType<TypeTag, Properties::GridVariables>;
auto gridVariables = std::make_shared<GridVariables>(problem, fvGridGeometry);
gridVariables->init(x, xOld);
gridVariables->init(x);
// get some time loop parameters
using Scalar = typename GET_PROP_TYPE(TypeTag, Scalar);
using Scalar = GetPropType<TypeTag, Properties::Scalar>;
const auto tEnd = getParam<Scalar>("TimeLoop.TEnd");
const auto maxDt = getParam<Scalar>("TimeLoop.MaxTimeStepSize");
auto dt = getParam<Scalar>("TimeLoop.DtInitial");
// intialize the vtk output module
using VtkOutputFields = typename GET_PROP_TYPE(TypeTag, VtkOutputFields);
VtkOutputModule<TypeTag> vtkWriter(*problem, *fvGridGeometry, *gridVariables, x, problem->name());
VtkOutputFields::init(vtkWriter); //! Add model specific output fields
// using VtkOutputFields = GetPropType<TypeTag, Properties::VtkOutputFields>;
// VtkOutputModule<TypeTag> vtkWriter(*problem, *fvGridGeometry, *gridVariables, x, problem->name());
// VtkOutputFields::init(vtkWriter); //! Add model specific output fields
// intialize the vtk output module
using VtkOutputFields = GetPropType<TypeTag, Properties::VtkOutputFields>;
VtkOutputModule<GridVariables, SolutionVector> vtkWriter(*gridVariables, x, problem->name());
using VelocityOutput = GetPropType<TypeTag, Properties::VelocityOutput>;
vtkWriter.addVelocityOutput(std::make_shared<VelocityOutput>(*gridVariables));
VtkOutputFields::initOutputModule(vtkWriter); //!< Add model specific output fields
vtkWriter.write(restartTime);
// instantiate time loop
auto timeLoop = std::make_shared<TimeLoop<Scalar>>(0.0, dt, tEnd);
auto timeLoop = std::make_shared<TimeLoop<Scalar>>(restartTime, dt, tEnd);
timeLoop->setMaxTimeStepSize(maxDt);
// the assembler with time loop for instationary problem
......@@ -124,8 +149,8 @@ int main(int argc, char** argv) try
auto linearSolver = std::make_shared<LinearSolver>(leafGridView, fvGridGeometry->dofMapper());
// the non-linear solver
using PrimaryVariableSwitch = typename GET_PROP_TYPE(TypeTag, PrimaryVariableSwitch);
using NewtonSolver = Dumux::PriVarSwitchNewtonSolver<Assembler, LinearSolver, PrimaryVariableSwitch>;
// using PrimaryVariableSwitch = GetPropType<TypeTag, Properties::PrimaryVariableSwitch>;
using NewtonSolver = NewtonSolver<Assembler, LinearSolver>;
NewtonSolver nonLinearSolver(assembler, linearSolver);
// time loop
......
......@@ -38,25 +38,32 @@ template <class TypeTag>
class Injection2p2cProblem;
namespace Properties {
NEW_TYPE_TAG(Injection2p2cTypeTag, INHERITS_FROM(TwoPTwoC));
NEW_TYPE_TAG(Injection2p2cCCTypeTag, INHERITS_FROM(CCTpfaModel, Injection2p2cTypeTag));
// Create new type tags
namespace TTag {
struct Injection2p2cTypeTag { using InheritsFrom = std::tuple<TwoPTwoC>; };
struct Injection2p2cCCTypeTag { using InheritsFrom = std::tuple<CCTpfaModel, Injection2p2cTypeTag>; };
} // end namespace TTag
// Set the grid type
SET_TYPE_PROP(Injection2p2cTypeTag, Grid, Dune::YaspGrid<2>);
template<class TypeTag>
struct Grid<TypeTag, TTag::Injection2p2cTypeTag> { using type = Dune::YaspGrid<2>; };
// Set the problem property
SET_TYPE_PROP(Injection2p2cTypeTag, Problem, Injection2p2cProblem<TypeTag>);
template<class TypeTag>
struct Problem<TypeTag, TTag::Injection2p2cTypeTag> { using type = Injection2p2cProblem<TypeTag>; };
// Set the spatial parameters
SET_TYPE_PROP(Injection2p2cTypeTag, SpatialParams,
InjectionSpatialParams<typename GET_PROP_TYPE(TypeTag, FVGridGeometry),
typename GET_PROP_TYPE(TypeTag, Scalar)>);
InjectionSpatialParams<GetPropType<TypeTag, Properties::FVGridGeometry>,
GetPropType<TypeTag, Properties::Scalar>>);
// Set fluid configuration
SET_TYPE_PROP(Injection2p2cTypeTag, FluidSystem, FluidSystems::H2ON2<typename GET_PROP_TYPE(TypeTag, Scalar), FluidSystems::H2ON2DefaultPolicy</*fastButSimplifiedRelations=*/ true>>);
template<class TypeTag>
struct FluidSystem<TypeTag, TTag::Injection2p2cTypeTag> { using type = FluidSystems::H2ON2<GetPropType<TypeTag, Properties::Scalar>, FluidSystems::H2ON2DefaultPolicy</*fastButSimplifiedRelations=*/ true>>; };
// Define whether mole(true) or mass (false) fractions are used
SET_BOOL_PROP(Injection2p2cTypeTag, UseMoles, true);
template<class TypeTag>
struct UseMoles<TypeTag, TTag::Injection2p2cTypeTag> { static constexpr bool value = true; };
} // end namespace Properties
/*!
......@@ -84,15 +91,15 @@ template<class TypeTag>
class Injection2p2cProblem : public PorousMediumFlowProblem<TypeTag>
{
using ParentType = PorousMediumFlowProblem<TypeTag>;
using GridView = typename GET_PROP_TYPE(TypeTag, GridView);
using Scalar = typename GET_PROP_TYPE(TypeTag, Scalar);
using Indices = typename GET_PROP_TYPE(TypeTag, ModelTraits)::Indices;
using PrimaryVariables = typename GET_PROP_TYPE(TypeTag, PrimaryVariables);
using BoundaryTypes = typename GET_PROP_TYPE(TypeTag, BoundaryTypes);
using FVGridGeometry = typename GET_PROP_TYPE(TypeTag, FVGridGeometry);
using FVElementGeometry = typename GET_PROP_TYPE(TypeTag, FVGridGeometry)::LocalView;
using FluidSystem = typename GET_PROP_TYPE(TypeTag, FluidSystem);
using NumEqVector = typename GET_PROP_TYPE(TypeTag, NumEqVector);
using GridView = GetPropType<TypeTag, Properties::GridView>;
using Scalar = GetPropType<TypeTag, Properties::Scalar>;
using Indices = typename GetPropType<TypeTag, Properties::ModelTraits>::Indices;
using PrimaryVariables = GetPropType<TypeTag, Properties::PrimaryVariables>;
using BoundaryTypes = GetPropType<TypeTag, Properties::BoundaryTypes>;
using FVGridGeometry = GetPropType<TypeTag, Properties::FVGridGeometry>;
using FVElementGeometry = typename GetPropType<TypeTag, Properties::FVGridGeometry>::LocalView;
using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
using NumEqVector = GetPropType<TypeTag, Properties::NumEqVector>;
enum { dimWorld = GridView::dimensionworld };
using Element = typename GridView::template Codim<0>::Entity;
......
......@@ -44,22 +44,28 @@ namespace Properties
* TODO:dumux-course-task:
* inherit from the TwoPNI model instead of TwoP here
*/
NEW_TYPE_TAG(Injection2pNITypeTag, INHERITS_FROM(TwoP));
NEW_TYPE_TAG(Injection2pNICCTypeTag, INHERITS_FROM(CCTpfaModel, Injection2pNITypeTag));
// Create new type tags
namespace TTag {
struct Injection2pNITypeTag { using InheritsFrom = std::tuple<TwoP>; };
struct Injection2pNICCTypeTag { using InheritsFrom = std::tuple<Injection2pNITypeTag, CCTpfaModel>; };
} // end namespace TTag
// Set the grid type
SET_TYPE_PROP(Injection2pNITypeTag, Grid, Dune::YaspGrid<2>);
template<class TypeTag>
struct Grid<TypeTag, TTag::Injection2pNITypeTag> { using type = Dune::YaspGrid<2>; };
// Set the problem property
SET_TYPE_PROP(Injection2pNITypeTag, Problem, InjectionProblem2PNI<TypeTag>);
template<class TypeTag>
struct Problem<TypeTag, TTag::Injection2pNITypeTag> { using type = InjectionProblem2PNI<TypeTag>; };
// Set the spatial parameters
SET_TYPE_PROP(Injection2pNITypeTag, SpatialParams,
InjectionSpatialParams<typename GET_PROP_TYPE(TypeTag, FVGridGeometry),
typename GET_PROP_TYPE(TypeTag, Scalar)>);
InjectionSpatialParams<GetPropType<TypeTag, Properties::FVGridGeometry>,
GetPropType<TypeTag, Properties::Scalar>>);
// Set fluid configuration
SET_TYPE_PROP(Injection2pNITypeTag, FluidSystem, FluidSystems::H2ON2<typename GET_PROP_TYPE(TypeTag, Scalar), FluidSystems::H2ON2DefaultPolicy</*fastButSimplifiedRelations=*/ true>>);
template<class TypeTag>
struct FluidSystem<TypeTag, TTag::Injection2pNITypeTag> { using type = FluidSystems::H2ON2<GetPropType<TypeTag, Properties::Scalar>, FluidSystems::H2ON2DefaultPolicy</*fastButSimplifiedRelations=*/ true>>; };
} // end namespace Properties
/*!
......@@ -87,15 +93,15 @@ template<class TypeTag>
class InjectionProblem2PNI : public PorousMediumFlowProblem<TypeTag>
{
using ParentType = PorousMediumFlowProblem<TypeTag>;
using GridView = typename GET_PROP_TYPE(TypeTag, GridView);
using Scalar = typename GET_PROP_TYPE(TypeTag, Scalar);
using Indices = typename GET_PROP_TYPE(TypeTag, ModelTraits)::Indices;
using PrimaryVariables = typename GET_PROP_TYPE(TypeTag, PrimaryVariables);
using BoundaryTypes = typename GET_PROP_TYPE(TypeTag, BoundaryTypes);
using FVGridGeometry = typename GET_PROP_TYPE(TypeTag, FVGridGeometry);
using FVElementGeometry = typename GET_PROP_TYPE(TypeTag, FVGridGeometry)::LocalView;
using FluidSystem = typename GET_PROP_TYPE(TypeTag, FluidSystem);
using NumEqVector = typename GET_PROP_TYPE(TypeTag, NumEqVector);
using GridView = GetPropType<TypeTag, Properties::GridView>;
using Scalar = GetPropType<TypeTag, Properties::Scalar>;
using Indices = typename GetPropType<TypeTag, Properties::ModelTraits>::Indices;
using PrimaryVariables = GetPropType<TypeTag, Properties::PrimaryVariables>;
using BoundaryTypes = GetPropType<TypeTag, Properties::BoundaryTypes>;
using FVGridGeometry = GetPropType<TypeTag, Properties::FVGridGeometry>;
using FVElementGeometry = typename GetPropType<TypeTag, Properties::FVGridGeometry>::LocalView;
using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
using NumEqVector = GetPropType<TypeTag, Properties::NumEqVector>;
enum { dimWorld = GridView::dimensionworld };
using Element = typename GridView::template Codim<0>::Entity;
......
......@@ -40,22 +40,28 @@ class InjectionProblem2P;
namespace Properties {
// define the TypeTag for this problem with a cell-centered two-point flux approximation spatial discretization.
NEW_TYPE_TAG(Injection2pTypeTag, INHERITS_FROM(TwoP));
NEW_TYPE_TAG(Injection2pCCTypeTag, INHERITS_FROM(CCTpfaModel, Injection2pTypeTag));
// Create new type tags
namespace TTag {
struct Injection2pTypeTag { using InheritsFrom = std::tuple<TwoP>; };
struct Injection2pCCTypeTag { using InheritsFrom = std::tuple<Injection2pTypeTag, CCTpfaModel>; };
} // end namespace TTag
// Set the grid type
SET_TYPE_PROP(Injection2pTypeTag, Grid, Dune::YaspGrid<2>);
template<class TypeTag>
struct Grid<TypeTag, TTag::Injection2pTypeTag> { using type = Dune::YaspGrid<2>; };
// Set the problem property
SET_TYPE_PROP(Injection2pTypeTag, Problem, InjectionProblem2P<TypeTag>);
template<class TypeTag>
struct Problem<TypeTag, TTag::Injection2pTypeTag> { using type = InjectionProblem2P<TypeTag>; };
// Set the spatial parameters
SET_TYPE_PROP(Injection2pTypeTag, SpatialParams,
InjectionSpatialParams<typename GET_PROP_TYPE(TypeTag, FVGridGeometry),
typename GET_PROP_TYPE(TypeTag, Scalar)>);
InjectionSpatialParams<GetPropType<TypeTag, Properties::FVGridGeometry>,
GetPropType<TypeTag, Properties::Scalar>>);
// Set fluid configuration
SET_TYPE_PROP(Injection2pTypeTag, FluidSystem, FluidSystems::H2ON2<typename GET_PROP_TYPE(TypeTag, Scalar), FluidSystems::H2ON2DefaultPolicy</*fastButSimplifiedRelations=*/ true>>);
template<class TypeTag>
struct FluidSystem<TypeTag, TTag::Injection2pTypeTag> { using type = FluidSystems::H2ON2<GetPropType<TypeTag, Properties::Scalar>, FluidSystems::H2ON2DefaultPolicy</*fastButSimplifiedRelations=*/ true>>; };
} // end namespace Properties
/*!
......@@ -83,15 +89,15 @@ template<class TypeTag>
class InjectionProblem2P : public PorousMediumFlowProblem<TypeTag>
{
using ParentType = PorousMediumFlowProblem<TypeTag>;
using GridView = typename GET_PROP_TYPE(TypeTag, GridView);
using Scalar = typename GET_PROP_TYPE(TypeTag, Scalar);
using Indices = typename GET_PROP_TYPE(TypeTag, ModelTraits)::Indices;
using PrimaryVariables = typename GET_PROP_TYPE(TypeTag, PrimaryVariables);
using BoundaryTypes = typename GET_PROP_TYPE(TypeTag, BoundaryTypes);
using FVGridGeometry = typename GET_PROP_TYPE(TypeTag, FVGridGeometry);
using FVElementGeometry = typename GET_PROP_TYPE(TypeTag, FVGridGeometry)::LocalView;
using FluidSystem = typename GET_PROP_TYPE(TypeTag, FluidSystem);
using NumEqVector = typename GET_PROP_TYPE(TypeTag, NumEqVector);
using GridView = GetPropType<TypeTag, Properties::GridView>;
using Scalar = GetPropType<TypeTag, Properties::Scalar>;
using Indices = typename GetPropType<TypeTag, Properties::ModelTraits>::Indices;
using PrimaryVariables = GetPropType<TypeTag, Properties::PrimaryVariables>;
using BoundaryTypes = GetPropType<TypeTag, Properties::BoundaryTypes>;
using FVGridGeometry = GetPropType<TypeTag, Properties::FVGridGeometry>;
using FVElementGeometry = typename GetPropType<TypeTag, Properties::FVGridGeometry>::LocalView;
using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
using NumEqVector = GetPropType<TypeTag, Properties::NumEqVector>;
enum { dimWorld = GridView::dimensionworld };
using Element = typename GridView::template Codim<0>::Entity;
......
......@@ -41,10 +41,11 @@
#include <dumux/assembly/fvassembler.hh>
#include <dumux/assembly/diffmethod.hh>
#include <dumux/discretization/methods.hh>
#include <dumux/discretization/method.hh>
#include <dumux/io/vtkoutputmodule.hh>
#include <dumux/io/grid/gridmanager.hh>
#include <dumux/io/loadsolution.hh>
// The problem file, where setup-specific boundary and initial conditions are defined.
#include "injection2pniproblem.hh"
......@@ -57,7 +58,7 @@ int main(int argc, char** argv) try
using namespace Dumux;
// define the type tag for this problem
using TypeTag = TTAG(Injection2pNICCTypeTag);
using TypeTag = Properties::TTag::Injection2pNICCTypeTag;
// initialize MPI, finalize is done automatically on exit
const auto& mpiHelper = Dune::MPIHelper::instance(argc, argv);
......@@ -70,7 +71,7 @@ int main(int argc, char** argv) try
Parameters::init(argc, argv);
// try to create a grid (from the given grid file or the input file)
GridManager<typename GET_PROP_TYPE(TypeTag, Grid)> gridManager;
GridManager<GetPropType<TypeTag, Properties::Grid>> gridManager;
gridManager.init();
////////////////////////////////////////////////////////////
......@@ -81,35 +82,53 @@ int main(int argc, char** argv) try
const auto& leafGridView = gridManager.grid().leafGridView();
// create the finite volume grid geometry
using FVGridGeometry = typename GET_PROP_TYPE(TypeTag, FVGridGeometry);
using FVGridGeometry = GetPropType<TypeTag, Properties::FVGridGeometry>;
auto fvGridGeometry = std::make_shared<FVGridGeometry>(leafGridView);
fvGridGeometry->update();
// the problem (initial and boundary conditions)
using Problem = typename GET_PROP_TYPE(TypeTag, Problem);
using Problem = GetPropType<TypeTag, Properties::Problem>;
auto problem = std::make_shared<Problem>(fvGridGeometry);
// check if we are about to restart a previously interrupted simulation
using Scalar = GetPropType<TypeTag, Properties::Scalar>;
Scalar restartTime = getParam<Scalar>("Restart.Time", 0);
// the solution vector
using SolutionVector = typename GET_PROP_TYPE(TypeTag, SolutionVector);
using SolutionVector = GetPropType<TypeTag, Properties::SolutionVector>;
SolutionVector x(fvGridGeometry->numDofs());
problem->applyInitialSolution(x);
if (restartTime > 0)
{
using IOFields = GetPropType<TypeTag, Properties::IOFields>;
using PrimaryVariables = GetPropType<TypeTag, Properties::PrimaryVariables>;
using ModelTraits = GetPropType<TypeTag, Properties::ModelTraits>;
using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
const auto fileName = getParam<std::string>("Restart.File");
const auto pvName = createPVNameFunction<IOFields, PrimaryVariables, ModelTraits, FluidSystem>();
loadSolution(x, fileName, pvName, *fvGridGeometry);
}
else
problem->applyInitialSolution(x);
auto xOld = x;
// the grid variables
using GridVariables = typename GET_PROP_TYPE(TypeTag, GridVariables);
using GridVariables = GetPropType<TypeTag, Properties::GridVariables>;
auto gridVariables = std::make_shared<GridVariables>(problem, fvGridGeometry);
gridVariables->init(x, xOld);
gridVariables->init(x);
// get some time loop parameters
using Scalar = typename GET_PROP_TYPE(TypeTag, Scalar);
using Scalar = GetPropType<TypeTag, Properties::Scalar>;
const auto tEnd = getParam<Scalar>("TimeLoop.TEnd");
const auto maxDt = getParam<Scalar>("TimeLoop.MaxTimeStepSize");
auto dt = getParam<Scalar>("TimeLoop.DtInitial");
// intialize the vtk output module
using VtkOutputFields = typename GET_PROP_TYPE(TypeTag, VtkOutputFields);
VtkOutputModule<TypeTag> vtkWriter(*problem, *fvGridGeometry, *gridVariables, x, problem->name());
VtkOutputFields::init(vtkWriter); //! Add model specific output fields
using VtkOutputFields = GetPropType<TypeTag, Properties::VtkOutputFields>;
VtkOutputModule<GridVariables, SolutionVector> vtkWriter(*gridVariables, x, problem->name());
using VelocityOutput = GetPropType<TypeTag, Properties::VelocityOutput>;
vtkWriter.addVelocityOutput(std::make_shared<VelocityOutput>(*gridVariables));
VtkOutputFields::initOutputModule(vtkWriter); //!< Add model specific output fields
vtkWriter.write(restartTime);
// instantiate time loop
auto timeLoop = std::make_shared<TimeLoop<Scalar>>(0.0, dt, tEnd);
......
......@@ -39,23 +39,28 @@ template <class TypeTag>
class InjectionProblem2PNI;
namespace Properties {
NEW_TYPE_TAG(Injection2pNITypeTag, INHERITS_FROM(TwoPNI));
NEW_TYPE_TAG(Injection2pNICCTypeTag, INHERITS_FROM(CCTpfaModel, Injection2pNITypeTag));
// Create new type tags
namespace TTag {
struct Injection2pNITypeTag { using InheritsFrom = std::tuple<TwoPNI>; };
struct Injection2pNICCTypeTag { using InheritsFrom = std::tuple<Injection2pNITypeTag, CCTpfaModel>; };
} // end namespace TTag
// Set the grid type
SET_TYPE_PROP(Injection2pNITypeTag, Grid, Dune::YaspGrid<2>);
template<class TypeTag>
struct Grid<TypeTag, TTag::Injection2pNITypeTag> { using type = Dune::YaspGrid<2>; };
// Set the problem property
SET_TYPE_PROP(Injection2pNITypeTag, Problem, InjectionProblem2PNI<TypeTag>);
template<class TypeTag>
struct Problem<TypeTag, TTag::Injection2pNITypeTag> { using type = InjectionProblem2PNI<TypeTag>; };
// Set the spatial parameters
SET_TYPE_PROP(Injection2pNITypeTag, SpatialParams,
InjectionSpatialParams<typename GET_PROP_TYPE(TypeTag, FVGridGeometry),
typename GET_PROP_TYPE(TypeTag, Scalar)>);
InjectionSpatialParams<GetPropType<TypeTag, Properties::FVGridGeometry>,
GetPropType<TypeTag, Properties::Scalar>>);
// Set fluid configuration
SET_TYPE_PROP(Injection2pNITypeTag, FluidSystem,
FluidSystems::H2ON2<typename GET_PROP_TYPE(TypeTag, Scalar), FluidSystems::H2ON2DefaultPolicy</*fastButSimplifiedRelations=*/ true>>);
FluidSystems::H2ON2<GetPropType<TypeTag, Properties::Scalar>, FluidSystems::H2ON2DefaultPolicy</*fastButSimplifiedRelations=*/ true>>);
} // end namespace Properties
/*!
......@@ -83,15 +88,15 @@ template<class TypeTag>
class InjectionProblem2PNI : public PorousMediumFlowProblem<TypeTag>
{
using ParentType = PorousMediumFlowProblem<TypeTag>;
using GridView = typename GET_PROP_TYPE(TypeTag, GridView);
using Scalar = typename GET_PROP_TYPE(TypeTag, Scalar);
using Indices = typename GET_PROP_TYPE(TypeTag, ModelTraits)::Indices;
using PrimaryVariables = typename GET_PROP_TYPE(TypeTag, PrimaryVariables);
using BoundaryTypes = typename GET_PROP_TYPE(TypeTag, BoundaryTypes);
using FVGridGeometry = typename GET_PROP_TYPE(TypeTag, FVGridGeometry);
using FVElementGeometry = typename GET_PROP_TYPE(TypeTag, FVGridGeometry)::LocalView;
using FluidSystem = typename GET_PROP_TYPE(TypeTag, FluidSystem);
using NumEqVector = typename GET_PROP_TYPE(TypeTag, NumEqVector);
using GridView = GetPropType<TypeTag, Properties::GridView>;
using Scalar = GetPropType<TypeTag, Properties::Scalar>;
using Indices = typename GetPropType<TypeTag, Properties::ModelTraits>::Indices;
using PrimaryVariables = GetPropType<TypeTag, Properties::PrimaryVariables>;
using BoundaryTypes = GetPropType<TypeTag, Properties::BoundaryTypes>;
using FVGridGeometry = GetPropType<TypeTag, Properties::FVGridGeometry>;
using FVElementGeometry = typename GetPropType<TypeTag, Properties::FVGridGeometry>::LocalView;
using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
using NumEqVector = GetPropType<TypeTag, Properties::NumEqVector>;
enum { dimWorld = GridView::dimensionworld };
using Element = typename GridView::template Codim<0>::Entity;
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment