diff --git a/test/porousmediumflow/richardsnc/implicit/richardswelltracerproblem.hh b/test/porousmediumflow/richardsnc/implicit/richardswelltracerproblem.hh index 5252fc5fe21cf811d5851103b9f535e895a801a4..ab6a5e7f8398190f30bffd14da1ccdb279f03fad 100644 --- a/test/porousmediumflow/richardsnc/implicit/richardswelltracerproblem.hh +++ b/test/porousmediumflow/richardsnc/implicit/richardswelltracerproblem.hh @@ -26,8 +26,9 @@ #ifndef DUMUX_RICHARDS_NC_WELL_TRACER_PROBLEM_HH #define DUMUX_RICHARDS_NC_WELL_TRACER_PROBLEM_HH -#include <dumux/implicit/cellcentered/tpfa/properties.hh> -#include <dumux/porousmediumflow/implicit/problem.hh> +#include <dumux/discretization/cellcentered/tpfa/properties.hh> +#include <dumux/discretization/box/properties.hh> +#include <dumux/porousmediumflow/problem.hh> #include <dumux/porousmediumflow/richardsnc/implicit/model.hh> #include "richardswelltracerspatialparams.hh" @@ -55,8 +56,6 @@ SET_TYPE_PROP(RichardsWellTracerProblem, Problem, RichardsWellTracerProblem<Type // Set the physical problem to be solved SET_TYPE_PROP(RichardsWellTracerProblem, PointSource, SolDependentPointSource<TypeTag>); -// Enable gravity -SET_BOOL_PROP(RichardsWellTracerProblem, ProblemEnableGravity, true); } /*! @@ -86,9 +85,9 @@ SET_BOOL_PROP(RichardsWellTracerProblem, ProblemEnableGravity, true); * simulation time is 10,000,000 seconds (115.7 days) */ template <class TypeTag> -class RichardsWellTracerProblem : public ImplicitPorousMediaProblem<TypeTag> +class RichardsWellTracerProblem : public PorousMediumFlowProblem<TypeTag> { - using ParentType = ImplicitPorousMediaProblem<TypeTag>; + using ParentType = PorousMediumFlowProblem<TypeTag>; using Problem = typename GET_PROP_TYPE(TypeTag, Problem); using GridView = typename GET_PROP_TYPE(TypeTag, GridView); using PrimaryVariables = typename GET_PROP_TYPE(TypeTag, PrimaryVariables); @@ -97,11 +96,13 @@ class RichardsWellTracerProblem : public ImplicitPorousMediaProblem<TypeTag> using SubControlVolume = typename GET_PROP_TYPE(TypeTag, SubControlVolume); using MaterialLaw = typename GET_PROP_TYPE(TypeTag, MaterialLaw); using BoundaryTypes = typename GET_PROP_TYPE(TypeTag, BoundaryTypes); - using TimeManager = typename GET_PROP_TYPE(TypeTag, TimeManager); using PointSource = typename GET_PROP_TYPE(TypeTag, PointSource); using FluidSystem = typename GET_PROP_TYPE(TypeTag, FluidSystem); using Scalar = typename GET_PROP_TYPE(TypeTag, Scalar); using Indices = typename GET_PROP_TYPE(TypeTag, Indices); + using FVGridGeometry = typename GET_PROP_TYPE(TypeTag, FVGridGeometry); + using SolutionVector = typename GET_PROP_TYPE(TypeTag, SolutionVector); + using GridVariables = typename GET_PROP_TYPE(TypeTag, GridVariables); enum { // copy some indices for convenience pressureIdx = Indices::pressureIdx, @@ -122,36 +123,37 @@ public: * \param timeManager The Dumux TimeManager for simulation management. * \param gridView The grid view on the spatial domain of the problem */ - RichardsWellTracerProblem(TimeManager &timeManager, const GridView &gridView) - : ParentType(timeManager, gridView) + RichardsWellTracerProblem(std::shared_ptr<const FVGridGeometry> fvGridGeometry) + : ParentType(fvGridGeometry) { - name_ = GET_RUNTIME_PARAM_FROM_GROUP(TypeTag, std::string, Problem, Name); - contaminantMoleFraction_ = GET_RUNTIME_PARAM_FROM_GROUP(TypeTag, Scalar, Problem, ContaminantMoleFraction); - pumpRate_ = GET_RUNTIME_PARAM_FROM_GROUP(TypeTag, Scalar, Problem, PumpRate); // in kg/s + name_ = getParam<std::string>("Problem.Name"); + contaminantMoleFraction_ = getParam<Scalar>("Problem.ContaminantMoleFraction"); + pumpRate_ = getParam<Scalar>("Problem.PumpRate"); // in kg/s // for initial conditions const Scalar sw = 0.4; // start with 80% saturation on top - pcTop_ = MaterialLaw::pc(this->spatialParams().materialLawParamsAtPos(this->bBoxMax()), sw); + pcTop_ = MaterialLaw::pc(this->spatialParams().materialLawParamsAtPos(this->fvGridGeometry().bBoxMax()), sw); // for post time step mass balance accumulatedSource_ = 0.0; } - void postTimeStep() - { - ParentType::postTimeStep(); + void postTimeStep(const SolutionVector& curSol, + const GridVariables& gridVariables, + const Scalar timeStepSize) + { // compute the mass in the entire domain to make sure the tracer is conserved Scalar tracerMass = 0.0; // bulk elements - for (const auto& element : elements(this->gridView())) + for (const auto& element : elements(this->fvGridGeometry().gridView())) { - auto fvGeometry = localView(this->model().fvGridGeometry()); + auto fvGeometry = localView(this->fvGridGeometry()); fvGeometry.bindElement(element); - auto elemVolVars = localView(this->model().curGlobalVolVars()); - elemVolVars.bindElement(element, fvGeometry, this->model().curSol()); + auto elemVolVars = localView(gridVariables.curGridVolVars()); + elemVolVars.bindElement(element, fvGeometry, curSol); for (auto&& scv : scvs(fvGeometry)) { @@ -162,7 +164,7 @@ public: accumulatedSource_ += this->scvPointSources(element, fvGeometry, elemVolVars, scv)[compIdx] * scv.volume() * volVars.extrusionFactor() * FluidSystem::molarMass(compIdx) - * this->timeManager().timeStepSize(); + * timeStepSize; } } @@ -269,7 +271,7 @@ public: */ void addPointSources(std::vector<PointSource>& pointSources) const { - auto globalPos = this->bBoxMax()-this->bBoxMin(); + auto globalPos = this->fvGridGeometry().bBoxMax()-this->fvGridGeometry().bBoxMin(); globalPos *= 0.5; //! Add point source in middle of domain pointSources.emplace_back(globalPos, @@ -305,8 +307,8 @@ private: { const auto xTracer = [&,this]() { - const GlobalPosition contaminationPos({0.2*this->bBoxMax()[0], 0.5*this->bBoxMax()[1]}); - if ((globalPos - contaminationPos).two_norm() < 0.1*(this->bBoxMax()-this->bBoxMin()).two_norm() + eps_) + const GlobalPosition contaminationPos({0.2*this->fvGridGeometry().bBoxMax()[0], 0.5*this->fvGridGeometry().bBoxMax()[1]}); + if ((globalPos - contaminationPos).two_norm() < 0.1*(this->fvGridGeometry().bBoxMax()-this->fvGridGeometry().bBoxMin()).two_norm() + eps_) return contaminantMoleFraction_; else return 0.0; @@ -315,29 +317,29 @@ private: PrimaryVariables values(0.0); //! hydrostatic pressure profile values[pressureIdx] = (nonWettingReferencePressure() - pcTop_) - - 9.81*1000*(globalPos[dimWorld-1] - this->bBoxMax()[dimWorld-1]); + - 9.81*1000*(globalPos[dimWorld-1] - this->fvGridGeometry().bBoxMax()[dimWorld-1]); values[compIdx] = xTracer; return values; } bool onLeftBoundary_(const GlobalPosition &globalPos) const { - return globalPos[0] < this->bBoxMin()[0] + eps_; + return globalPos[0] < this->fvGridGeometry().bBoxMin()[0] + eps_; } bool onRightBoundary_(const GlobalPosition &globalPos) const { - return globalPos[0] > this->bBoxMax()[0] - eps_; + return globalPos[0] > this->fvGridGeometry().bBoxMax()[0] - eps_; } bool onLowerBoundary_(const GlobalPosition &globalPos) const { - return globalPos[1] < this->bBoxMin()[1] + eps_; + return globalPos[1] < this->fvGridGeometry().bBoxMin()[1] + eps_; } bool onUpperBoundary_(const GlobalPosition &globalPos) const { - return globalPos[1] > this->bBoxMax()[1] - eps_; + return globalPos[1] > this->fvGridGeometry().bBoxMax()[1] - eps_; } static constexpr Scalar eps_ = 1.5e-7; diff --git a/test/porousmediumflow/richardsnc/implicit/richardswelltracerspatialparams.hh b/test/porousmediumflow/richardsnc/implicit/richardswelltracerspatialparams.hh index 9fe7171f4b37d46989332f78383ddd344c55e6d2..3fdda1d37e543a6e8db151e3ac255f1154e67655 100644 --- a/test/porousmediumflow/richardsnc/implicit/richardswelltracerspatialparams.hh +++ b/test/porousmediumflow/richardsnc/implicit/richardswelltracerspatialparams.hh @@ -89,12 +89,12 @@ public: * \param gridView The DUNE GridView representing the spatial * domain of the problem. */ - RichardsWellTracerSpatialParams(const Problem& problem, const GridView& gridView) - : ParentType(problem, gridView) + RichardsWellTracerSpatialParams(const Problem& problem) + : ParentType(problem) { - lensLowerLeft_ = GET_RUNTIME_PARAM_FROM_GROUP(TypeTag, GlobalPosition, Problem, LensLowerLeft); - lensUpperRight_ = GET_RUNTIME_PARAM_FROM_GROUP(TypeTag, GlobalPosition, Problem, LensUpperRight); + lensLowerLeft_ = getParam<GlobalPosition>("Problem.LensLowerLeft"); + lensUpperRight_ = getParam<GlobalPosition>("Problem.LensUpperRight"); // residual saturations lensMaterialParams_.setSwr(0.18); diff --git a/test/porousmediumflow/richardsnc/implicit/test_boxrichardsnc.cc b/test/porousmediumflow/richardsnc/implicit/test_boxrichardsnc.cc index 1e5b947b19956d7fe0a7627452f1946522c38fa4..d1092acbca73477187180b79b1863a3439569a0b 100644 --- a/test/porousmediumflow/richardsnc/implicit/test_boxrichardsnc.cc +++ b/test/porousmediumflow/richardsnc/implicit/test_boxrichardsnc.cc @@ -19,12 +19,35 @@ /*! * \file * - * \brief Test for the Richards n-component Box model. + * \brief Test for the Richards CC model. */ #include <config.h> + #include "richardswelltracerproblem.hh" -#include <dumux/common/start.hh> +#include <ctime> +#include <iostream> + +#include <dune/common/parallel/mpihelper.hh> +#include <dune/common/timer.hh> +#include <dune/grid/io/file/dgfparser/dgfexception.hh> +#include <dune/grid/io/file/vtk.hh> +#include <dune/istl/io.hh> + +#include <dumux/common/propertysystem.hh> +#include <dumux/common/parameters.hh> +#include <dumux/common/valgrind.hh> +#include <dumux/common/dumuxmessage.hh> +#include <dumux/common/defaultusagemessage.hh> + +#include <dumux/linear/amgbackend.hh> +#include <dumux/nonlinear/newtonmethod.hh> +#include <dumux/nonlinear/newtoncontroller.hh> +#include <dumux/porousmediumflow/richards/implicit/newtoncontroller.hh> + +#include <dumux/assembly/fvassembler.hh> + +#include <dumux/io/vtkoutputmodule.hh> /*! * \brief Provides an interface for customizing error messages associated with * reading in parameters. @@ -43,16 +66,183 @@ void usage(const char *progName, const std::string &errorMsg) errorMessageOut += "\n\nThe list of mandatory options for this program is:\n" "\t-TimeManager.TEnd End of the simulation [s] \n" "\t-TimeManager.DtInitial Initial timestep size [s] \n" - "\t-Grid.UpperRight Upper right corner coordinates [m] \n" - "\t-Grid.Cells Number of cells in each coordinate direction [-] \n"; + "\t-Grid.File Name of the file containing the grid \n" + "\t definition in DGF format\n"; std::cout << errorMessageOut << "\n"; } } -int main(int argc, char** argv) +//////////////////////// +// the main function +//////////////////////// +int main(int argc, char** argv) try +{ + using namespace Dumux; + + // define the type tag for this problem + using TypeTag = TTAG(RichardsWellTracerBoxProblem); + + // initialize MPI, finalize is done automatically on exit + const auto& mpiHelper = Dune::MPIHelper::instance(argc, argv); + + // print dumux start message + if (mpiHelper.rank() == 0) + DumuxMessage::print(/*firstCall=*/true); + + // parse command line arguments and input file + Parameters::init(argc, argv, usage); + + // try to create a grid (from the given grid file or the input file) + using GridCreator = typename GET_PROP_TYPE(TypeTag, GridCreator); + GridCreator::makeGrid(); + GridCreator::loadBalance(); + + //////////////////////////////////////////////////////////// + // run instationary non-linear problem on this grid + //////////////////////////////////////////////////////////// + + // we compute on the leaf grid view + const auto& leafGridView = GridCreator::grid().leafGridView(); + + // create the finite volume grid geometry + using FVGridGeometry = typename GET_PROP_TYPE(TypeTag, FVGridGeometry); + auto fvGridGeometry = std::make_shared<FVGridGeometry>(leafGridView); + fvGridGeometry->update(); + + // the problem (initial and boundary conditions) + using Problem = typename GET_PROP_TYPE(TypeTag, Problem); + auto problem = std::make_shared<Problem>(fvGridGeometry); + + // the solution vector + using SolutionVector = typename GET_PROP_TYPE(TypeTag, SolutionVector); + SolutionVector x(fvGridGeometry->numDofs()); + problem->applyInitialSolution(x); + auto xOld = x; + + // the grid variables + using GridVariables = typename GET_PROP_TYPE(TypeTag, GridVariables); + auto gridVariables = std::make_shared<GridVariables>(problem, fvGridGeometry); + gridVariables->init(x, xOld); + + // get some time loop parameters + using Scalar = typename GET_PROP_TYPE(TypeTag, Scalar); + const auto tEnd = getParam<Scalar>("TimeLoop.TEnd"); + const auto maxDivisions = getParam<int>("TimeLoop.MaxTimeStepDivisions"); + const auto maxDt = getParam<Scalar>("TimeLoop.MaxTimeStepSize"); + auto dt = getParam<Scalar>("TimeLoop.DtInitial"); + + // check if we are about to restart a previously interrupted simulation + Scalar restartTime = 0; + if (Parameters::getTree().hasKey("Restart") || Parameters::getTree().hasKey("TimeLoop.Restart")) + restartTime = getParam<Scalar>("TimeLoop.Restart"); + + // 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 + vtkWriter.write(0.0); + + // instantiate time loop + auto timeLoop = std::make_shared<TimeLoop<Scalar>>(restartTime, dt, tEnd); + timeLoop->setMaxTimeStepSize(maxDt); + + // the assembler with time loop for instationary problem + using Assembler = FVAssembler<TypeTag, DiffMethod::numeric>; + auto assembler = std::make_shared<Assembler>(problem, fvGridGeometry, gridVariables, timeLoop); + + // the linear solver + using LinearSolver = Dumux::ILU0BiCGSTABBackend<TypeTag>; + auto linearSolver = std::make_shared<LinearSolver>(); + + // the non-linear solver + using NewtonController = Dumux::RichardsNewtonController<TypeTag>; + using NewtonMethod = Dumux::NewtonMethod<NewtonController, Assembler, LinearSolver>; + auto newtonController = std::make_shared<NewtonController>(leafGridView.comm(), timeLoop); + NewtonMethod nonLinearSolver(newtonController, assembler, linearSolver); + + // time loop + timeLoop->start(); do + { + // set previous solution for storage evaluations + assembler->setPreviousSolution(xOld); + + // try solving the non-linear system + for (int i = 0; i < maxDivisions; ++i) + { + // linearize & solve + auto converged = nonLinearSolver.solve(x); + + if (converged) + break; + + if (!converged && i == maxDivisions-1) + DUNE_THROW(Dune::MathError, + "Newton solver didn't converge after " + << maxDivisions + << " time-step divisions. dt=" + << timeLoop->timeStepSize() + << ".\nThe solutions of the current and the previous time steps " + << "have been saved to restart files."); + } + + // make the new solution the old solution + xOld = x; + problem->postTimeStep(x, *gridVariables, timeLoop->timeStepSize()); + gridVariables->advanceTimeStep(); + + // advance to the time loop to the next step + timeLoop->advanceTimeStep(); + + // write vtk output + vtkWriter.write(timeLoop->time()); + + // report statistics of this time step + timeLoop->reportTimeStep(); + + // set new dt as suggested by newton controller + timeLoop->setTimeStepSize(newtonController->suggestTimeStepSize(timeLoop->timeStepSize())); + + } while (!timeLoop->finished()); + + timeLoop->finalize(leafGridView.comm()); + + //////////////////////////////////////////////////////////// + // finalize, print dumux message to say goodbye + //////////////////////////////////////////////////////////// + + // print dumux end message + if (mpiHelper.rank() == 0) + { + Parameters::print(); + DumuxMessage::print(/*firstCall=*/false); + } + + return 0; +} + +catch (Dumux::ParameterException &e) +{ + std::cerr << std::endl << e << " ---> Abort!" << std::endl; + return 1; +} +catch (Dune::DGFException & e) +{ + std::cerr << "DGF exception thrown (" << e << + "). Most likely, the DGF file name is wrong " + "or the DGF file is corrupted, " + "e.g. missing hash at end of file or wrong number (dimensions) of entries." + << " ---> Abort!" << std::endl; + return 2; +} +catch (Dune::Exception &e) +{ + std::cerr << "Dune reported error: " << e << " ---> Abort!" << std::endl; + return 3; +} +catch (...) { - using ProblemTypeTag = TTAG(RichardsWellTracerBoxProblem); - return Dumux::start<ProblemTypeTag>(argc, argv, usage); + std::cerr << "Unknown exception thrown! ---> Abort!" << std::endl; + return 4; } diff --git a/test/porousmediumflow/richardsnc/implicit/test_boxrichardsnc.input b/test/porousmediumflow/richardsnc/implicit/test_boxrichardsnc.input index 24c6c481a5d27f969aee238d710885eb84a7e152..ad06d912a51f3e6bcc4fc6e8fe8627c7db71bde1 100644 --- a/test/porousmediumflow/richardsnc/implicit/test_boxrichardsnc.input +++ b/test/porousmediumflow/richardsnc/implicit/test_boxrichardsnc.input @@ -1,4 +1,4 @@ -[TimeManager] +[TimeLoop] DtInitial = 10 # [s] TEnd = 1e6 # [s] @@ -22,3 +22,6 @@ Name = "D2O" [Vtk] AddVelocity = true + +[Newton] +EnableChop = false # chop for better convergence diff --git a/test/porousmediumflow/richardsnc/implicit/test_ccrichardsnc.cc b/test/porousmediumflow/richardsnc/implicit/test_ccrichardsnc.cc index b779e6123ace5bd6c25cf88184403f7e3d1601ae..dce20959bce306bf2c61925a1f82655a7cfefc63 100644 --- a/test/porousmediumflow/richardsnc/implicit/test_ccrichardsnc.cc +++ b/test/porousmediumflow/richardsnc/implicit/test_ccrichardsnc.cc @@ -19,12 +19,35 @@ /*! * \file * - * \brief Test for the Richards n-component CC model. + * \brief Test for the Richards CC model. */ #include <config.h> + #include "richardswelltracerproblem.hh" -#include <dumux/common/start.hh> +#include <ctime> +#include <iostream> + +#include <dune/common/parallel/mpihelper.hh> +#include <dune/common/timer.hh> +#include <dune/grid/io/file/dgfparser/dgfexception.hh> +#include <dune/grid/io/file/vtk.hh> +#include <dune/istl/io.hh> + +#include <dumux/common/propertysystem.hh> +#include <dumux/common/parameters.hh> +#include <dumux/common/valgrind.hh> +#include <dumux/common/dumuxmessage.hh> +#include <dumux/common/defaultusagemessage.hh> + +#include <dumux/linear/amgbackend.hh> +#include <dumux/nonlinear/newtonmethod.hh> +#include <dumux/nonlinear/newtoncontroller.hh> +#include <dumux/porousmediumflow/richards/implicit/newtoncontroller.hh> + +#include <dumux/assembly/fvassembler.hh> + +#include <dumux/io/vtkoutputmodule.hh> /*! * \brief Provides an interface for customizing error messages associated with * reading in parameters. @@ -43,16 +66,183 @@ void usage(const char *progName, const std::string &errorMsg) errorMessageOut += "\n\nThe list of mandatory options for this program is:\n" "\t-TimeManager.TEnd End of the simulation [s] \n" "\t-TimeManager.DtInitial Initial timestep size [s] \n" - "\t-Grid.UpperRight Upper right corner coordinates [m] \n" - "\t-Grid.Cells Number of cells in each coordinate direction [-] \n"; + "\t-Grid.File Name of the file containing the grid \n" + "\t definition in DGF format\n"; std::cout << errorMessageOut << "\n"; } } -int main(int argc, char** argv) +//////////////////////// +// the main function +//////////////////////// +int main(int argc, char** argv) try +{ + using namespace Dumux; + + // define the type tag for this problem + using TypeTag = TTAG(RichardsWellTracerCCProblem); + + // initialize MPI, finalize is done automatically on exit + const auto& mpiHelper = Dune::MPIHelper::instance(argc, argv); + + // print dumux start message + if (mpiHelper.rank() == 0) + DumuxMessage::print(/*firstCall=*/true); + + // parse command line arguments and input file + Parameters::init(argc, argv, usage); + + // try to create a grid (from the given grid file or the input file) + using GridCreator = typename GET_PROP_TYPE(TypeTag, GridCreator); + GridCreator::makeGrid(); + GridCreator::loadBalance(); + + //////////////////////////////////////////////////////////// + // run instationary non-linear problem on this grid + //////////////////////////////////////////////////////////// + + // we compute on the leaf grid view + const auto& leafGridView = GridCreator::grid().leafGridView(); + + // create the finite volume grid geometry + using FVGridGeometry = typename GET_PROP_TYPE(TypeTag, FVGridGeometry); + auto fvGridGeometry = std::make_shared<FVGridGeometry>(leafGridView); + fvGridGeometry->update(); + + // the problem (initial and boundary conditions) + using Problem = typename GET_PROP_TYPE(TypeTag, Problem); + auto problem = std::make_shared<Problem>(fvGridGeometry); + + // the solution vector + using SolutionVector = typename GET_PROP_TYPE(TypeTag, SolutionVector); + SolutionVector x(fvGridGeometry->numDofs()); + problem->applyInitialSolution(x); + auto xOld = x; + + // the grid variables + using GridVariables = typename GET_PROP_TYPE(TypeTag, GridVariables); + auto gridVariables = std::make_shared<GridVariables>(problem, fvGridGeometry); + gridVariables->init(x, xOld); + + // get some time loop parameters + using Scalar = typename GET_PROP_TYPE(TypeTag, Scalar); + const auto tEnd = getParam<Scalar>("TimeLoop.TEnd"); + const auto maxDivisions = getParam<int>("TimeLoop.MaxTimeStepDivisions"); + const auto maxDt = getParam<Scalar>("TimeLoop.MaxTimeStepSize"); + auto dt = getParam<Scalar>("TimeLoop.DtInitial"); + + // check if we are about to restart a previously interrupted simulation + Scalar restartTime = 0; + if (Parameters::getTree().hasKey("Restart") || Parameters::getTree().hasKey("TimeLoop.Restart")) + restartTime = getParam<Scalar>("TimeLoop.Restart"); + + // 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 + vtkWriter.write(0.0); + + // instantiate time loop + auto timeLoop = std::make_shared<TimeLoop<Scalar>>(restartTime, dt, tEnd); + timeLoop->setMaxTimeStepSize(maxDt); + + // the assembler with time loop for instationary problem + using Assembler = FVAssembler<TypeTag, DiffMethod::numeric>; + auto assembler = std::make_shared<Assembler>(problem, fvGridGeometry, gridVariables, timeLoop); + + // the linear solver + using LinearSolver = Dumux::ILU0BiCGSTABBackend<TypeTag>; + auto linearSolver = std::make_shared<LinearSolver>(); + + // the non-linear solver + using NewtonController = Dumux::RichardsNewtonController<TypeTag>; + using NewtonMethod = Dumux::NewtonMethod<NewtonController, Assembler, LinearSolver>; + auto newtonController = std::make_shared<NewtonController>(leafGridView.comm(), timeLoop); + NewtonMethod nonLinearSolver(newtonController, assembler, linearSolver); + + // time loop + timeLoop->start(); do + { + // set previous solution for storage evaluations + assembler->setPreviousSolution(xOld); + + // try solving the non-linear system + for (int i = 0; i < maxDivisions; ++i) + { + // linearize & solve + auto converged = nonLinearSolver.solve(x); + + if (converged) + break; + + if (!converged && i == maxDivisions-1) + DUNE_THROW(Dune::MathError, + "Newton solver didn't converge after " + << maxDivisions + << " time-step divisions. dt=" + << timeLoop->timeStepSize() + << ".\nThe solutions of the current and the previous time steps " + << "have been saved to restart files."); + } + + // make the new solution the old solution + xOld = x; + problem->postTimeStep(x, *gridVariables, timeLoop->timeStepSize()); + gridVariables->advanceTimeStep(); + + // advance to the time loop to the next step + timeLoop->advanceTimeStep(); + + // write vtk output + vtkWriter.write(timeLoop->time()); + + // report statistics of this time step + timeLoop->reportTimeStep(); + + // set new dt as suggested by newton controller + timeLoop->setTimeStepSize(newtonController->suggestTimeStepSize(timeLoop->timeStepSize())); + + } while (!timeLoop->finished()); + + timeLoop->finalize(leafGridView.comm()); + + //////////////////////////////////////////////////////////// + // finalize, print dumux message to say goodbye + //////////////////////////////////////////////////////////// + + // print dumux end message + if (mpiHelper.rank() == 0) + { + Parameters::print(); + DumuxMessage::print(/*firstCall=*/false); + } + + return 0; +} + +catch (Dumux::ParameterException &e) +{ + std::cerr << std::endl << e << " ---> Abort!" << std::endl; + return 1; +} +catch (Dune::DGFException & e) +{ + std::cerr << "DGF exception thrown (" << e << + "). Most likely, the DGF file name is wrong " + "or the DGF file is corrupted, " + "e.g. missing hash at end of file or wrong number (dimensions) of entries." + << " ---> Abort!" << std::endl; + return 2; +} +catch (Dune::Exception &e) +{ + std::cerr << "Dune reported error: " << e << " ---> Abort!" << std::endl; + return 3; +} +catch (...) { - using ProblemTypeTag = TTAG(RichardsWellTracerCCProblem); - return Dumux::start<ProblemTypeTag>(argc, argv, usage); + std::cerr << "Unknown exception thrown! ---> Abort!" << std::endl; + return 4; } diff --git a/test/porousmediumflow/richardsnc/implicit/test_ccrichardsnc.input b/test/porousmediumflow/richardsnc/implicit/test_ccrichardsnc.input index e2e52d593f5114612a8c40449e2986a930ba3f61..71d8dc9342fc8d1a87e54eccf06679a0f6601ade 100644 --- a/test/porousmediumflow/richardsnc/implicit/test_ccrichardsnc.input +++ b/test/porousmediumflow/richardsnc/implicit/test_ccrichardsnc.input @@ -1,4 +1,4 @@ -[TimeManager] +[TimeLoop] DtInitial = 10 # [s] TEnd = 1e6 # [s] @@ -22,3 +22,6 @@ Name = "D2O" [Vtk] AddVelocity = true + +[Newton] +EnableChop = false # chop for better convergences