// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
// vi: set et ts=4 sw=4 sts=4:
/*****************************************************************************
* See the file COPYING for full copying permissions. *
* *
* This program is free software: you can redistribute it and/or modify *
* it under the terms of the GNU General Public License as published by *
* the Free Software Foundation, either version 3 of the License, or *
* (at your option) any later version. *
* *
* This program is distributed in the hope that it will be useful, *
* but WITHOUT ANY WARRANTY; without even the implied warranty of *
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
* GNU General Public License for more details. *
* *
* You should have received a copy of the GNU General Public License *
* along with this program. If not, see . *
*****************************************************************************/
/*!
* \file
* \ingroup CO2Tests
* \brief Test for the two-phase two-component CC model.
*/
#include
#include
#include
#include
#include
#include
#include
#include
#include
#include
#include
#include
#include
#include
#include
#include
#include
#include
#include
#include
// the problem definitions
#include "problem.hh"
int main(int argc, char** argv) try
{
using namespace Dumux;
// define the type tag for this problem
using TypeTag = Properties::TTag::TYPETAG;
// 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);
// try to create a grid (from the given grid file or the input file)
GridManager> gridManager;
gridManager.init();
////////////////////////////////////////////////////////////
// run instationary non-linear problem on this grid
////////////////////////////////////////////////////////////
// we compute on the leaf grid view
const auto& leafGridView = gridManager.grid().leafGridView();
// create the finite volume grid geometry
using GridGeometry = GetPropType;
auto gridGeometry = std::make_shared(leafGridView);
gridGeometry->update();
// the spatial parameters
using SpatialParams = GetPropType;
auto spatialParams = std::make_shared(gridGeometry, gridManager.getGridData());
// the problem (initial and boundary conditions)
using Problem = GetPropType;
auto problem = std::make_shared(gridGeometry, spatialParams);
// the solution vector
using SolutionVector = GetPropType;
SolutionVector x(gridGeometry->numDofs());
problem->applyInitialSolution(x);
auto xOld = x;
// the grid variables
using GridVariables = GetPropType;
auto gridVariables = std::make_shared(problem, gridGeometry);
gridVariables->init(x);
// get some time loop parameters
using Scalar = GetPropType;
const auto tEnd = getParam("TimeLoop.TEnd");
const auto maxDt = getParam("TimeLoop.MaxTimeStepSize");
auto dt = getParam("TimeLoop.DtInitial");
// intialize the vtk output module
using IOFields = GetPropType;
VtkOutputModule vtkWriter(*gridVariables, x, problem->name());
using VelocityOutput = GetPropType;
vtkWriter.addVelocityOutput(std::make_shared(*gridVariables));
IOFields::initOutputModule(vtkWriter); // Add model specific output fields
problem->addFieldsToWriter(vtkWriter); //!< Add some more problem dependent fields
vtkWriter.write(0.0);
// instantiate time loop
auto timeLoop = std::make_shared>(0, dt, tEnd);
timeLoop->setMaxTimeStepSize(maxDt);
// the assembler with time loop for instationary problem
using Assembler = FVAssembler;
auto assembler = std::make_shared(problem, gridGeometry, gridVariables, timeLoop, xOld);
// the linear solver
using LinearSolver = AMGBiCGSTABBackend>;
auto linearSolver = std::make_shared(leafGridView, gridGeometry->dofMapper());
// the non-linear solver
using NewtonSolver = NewtonSolver;
NewtonSolver nonLinearSolver(assembler, linearSolver);
// time loop
timeLoop->start(); do
{
// solve the non-linear system with time step control
nonLinearSolver.solve(x, *timeLoop);
// make the new solution the old solution
xOld = x;
gridVariables->advanceTimeStep();
// advance to the time loop to the next step
timeLoop->advanceTimeStep();
// report statistics of this time step
timeLoop->reportTimeStep();
// double the timestep each time (it is still limited by maximum time step size)
timeLoop->setTimeStepSize(timeLoop->timeStepSize()*2.0);
// write vtk output
vtkWriter.write(timeLoop->time());
} 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;
} // end main
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 (...)
{
std::cerr << "Unknown exception thrown! ---> Abort!" << std::endl;
return 4;
}