2pmain.cc 7.09 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
// -*- 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 2 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 <http://www.gnu.org/licenses/>.   *
 *****************************************************************************/
/*!
 * \file
21
 * \brief The main file for the two-phase porousmediumflow problem of exercise-basic
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
 */
#include <config.h>

#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 <dumux/common/properties.hh>
#include <dumux/common/parameters.hh>
#include <dumux/common/dumuxmessage.hh>
#include <dumux/common/defaultusagemessage.hh>

#include <dumux/linear/amgbackend.hh>
39
#include <dumux/linear/linearsolvertraits.hh>
40
41
42
43
44
#include <dumux/nonlinear/newtonsolver.hh>

#include <dumux/assembly/fvassembler.hh>
#include <dumux/assembly/diffmethod.hh>

45
#include <dumux/discretization/method.hh>
46
47
48
49

#include <dumux/io/vtkoutputmodule.hh>
#include <dumux/io/grid/gridmanager.hh>

50
51
// The properties file, where the compile time options are defined
#include "properties2p.hh"
52
53
54
55

////////////////////////
// the main function
////////////////////////
56
int main(int argc, char** argv)
57
58
59
60
{
    using namespace Dumux;

    // define the type tag for this problem
Felix Weinhardt's avatar
Felix Weinhardt committed
61
    using TypeTag = Properties::TTag::Injection2pCC;
62
63
64
65
66
67
68
69
70
71
72
73

    // 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)
74
    GridManager<GetPropType<TypeTag, Properties::Grid>> gridManager;
75
76
77
78
79
80
81
82
83
84
    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
85
    using FVGridGeometry = GetPropType<TypeTag, Properties::GridGeometry>;
86
87
88
89
    auto fvGridGeometry = std::make_shared<FVGridGeometry>(leafGridView);
    fvGridGeometry->update();

    // the problem (initial and boundary conditions)
90
    using Problem = GetPropType<TypeTag, Properties::Problem>;
91
92
93
    auto problem = std::make_shared<Problem>(fvGridGeometry);

    // the solution vector
94
    using SolutionVector = GetPropType<TypeTag, Properties::SolutionVector>;
Timo Koch's avatar
Timo Koch committed
95
    SolutionVector x;
Timo Koch's avatar
Timo Koch committed
96
    problem->applyInitialSolution(x);
97
98
99
    auto xOld = x;

    // the grid variables
100
    using GridVariables = GetPropType<TypeTag, Properties::GridVariables>;
101
    auto gridVariables = std::make_shared<GridVariables>(problem, fvGridGeometry);
102
    gridVariables->init(x);
103
104

    // get some time loop parameters
105
    using Scalar = GetPropType<TypeTag, Properties::Scalar>;
106
107
108
109
110
    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
111
112
113
114
115
116
117
118
119
    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
Timo Koch's avatar
Timo Koch committed
120
    vtkWriter.write(0.0);
121
122
123
124
125
126
127

    // instantiate time loop
    auto timeLoop = std::make_shared<TimeLoop<Scalar>>(0.0, dt, tEnd);
    timeLoop->setMaxTimeStepSize(maxDt);

    // the assembler with time loop for instationary problem
    using Assembler = FVAssembler<TypeTag, DiffMethod::numeric>;
128
    auto assembler = std::make_shared<Assembler>(problem, fvGridGeometry, gridVariables, timeLoop, xOld);
129
130

    // the linear solver
131
    using LinearSolver = AMGBiCGSTABBackend<LinearSolverTraits<FVGridGeometry>>;
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
    auto linearSolver = std::make_shared<LinearSolver>(leafGridView, fvGridGeometry->dofMapper());

    // the non-linear solver
    using NewtonSolver = Dumux::NewtonSolver<Assembler, LinearSolver>;
    NewtonSolver nonLinearSolver(assembler, linearSolver);

    // time loop
    timeLoop->start();
    while (!timeLoop->finished())
    {
        // set previous solution for storage evaluations
        assembler->setPreviousSolution(xOld);

        //set time in problem (is used in time-dependent Neumann boundary condition)
        problem->setTime(timeLoop->time()+timeLoop->timeStepSize());

        // 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();

        // set new dt as suggested by the newton solver
        timeLoop->setTimeStepSize(nonLinearSolver.suggestTimeStepSize(timeLoop->timeStepSize()));

        // output to vtk
        vtkWriter.write(timeLoop->time());
    }

    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