Line | Branch | Exec | Source |
---|---|---|---|
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 | // SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder | ||
5 | // SPDX-License-Identifier: GPL-3.0-or-later | ||
6 | // | ||
7 | // ## The main file (`main.cc`) | ||
8 | // [[content]] | ||
9 | // | ||
10 | // ### Included header files | ||
11 | // [[details]] includes | ||
12 | // [[exclude]] | ||
13 | // Some generic includes. | ||
14 | #include <config.h> | ||
15 | #include <ctime> | ||
16 | #include <iostream> | ||
17 | // [[/exclude]] | ||
18 | |||
19 | // These time measurements and parallel backend initialization | ||
20 | #include <dune/common/timer.hh> | ||
21 | #include <dumux/common/initialize.hh> | ||
22 | |||
23 | // The following headers include functionality related to property definition or retrieval, as well as | ||
24 | // the retrieval of input parameters specified in the input file or via the command line. | ||
25 | #include <dumux/common/properties.hh> | ||
26 | #include <dumux/common/parameters.hh> | ||
27 | |||
28 | // The following files contains the available linear solver backends and the assembler for the linear | ||
29 | // systems arising from finite volume discretizations (box-scheme, tpfa-approximation, mpfa-approximation). | ||
30 | #include <dumux/linear/istlsolvers.hh> | ||
31 | #include <dumux/linear/linearsolvertraits.hh> | ||
32 | #include <dumux/linear/linearalgebratraits.hh> | ||
33 | #include <dumux/assembly/fvassembler.hh> | ||
34 | #include <dumux/assembly/diffmethod.hh> // analytic or numeric differentiation | ||
35 | |||
36 | // The following class provides a convenient way of writing DuMu<sup>x</sup> simulation results to VTK format. | ||
37 | #include <dumux/io/vtkoutputmodule.hh> | ||
38 | // The gridmanager constructs a grid from the information in the input or grid file. | ||
39 | // Many different Dune grid implementations are supported, of which a list can be found | ||
40 | // in `gridmanager.hh`. | ||
41 | #include <dumux/io/grid/gridmanager.hh> | ||
42 | |||
43 | // For both the single-phase and the tracer problem, `TypeTags` are defined, which collect | ||
44 | // the properties that are required for the simulation. These type tags are defined | ||
45 | // in the headers that we include here. For detailed information, please have a look | ||
46 | // at the documentation provided therein. | ||
47 | #include "properties_1p.hh" | ||
48 | #include "properties_tracer.hh" | ||
49 | // [[/details]] | ||
50 | // | ||
51 | |||
52 | // ### The main function | ||
53 | // We will now discuss the main program flow implemented within the `main` function. | ||
54 | // At the beginning of each program using Dune, an instance of `Dune::MPIHelper` has to | ||
55 | // be created. Moreover, we parse the run-time arguments from the command line and the | ||
56 | // input file: | ||
57 | // [[codeblock]] | ||
58 | 1 | int main(int argc, char** argv) try | |
59 | { | ||
60 | 1 | using namespace Dumux; | |
61 | |||
62 | // maybe initialize MPI and/or multithreading backend | ||
63 |
1/2✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
|
1 | Dumux::initialize(argc, argv); |
64 | |||
65 | // parse command line arguments and input file | ||
66 |
9/28✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 1 times.
✗ Branch 8 not taken.
✓ Branch 10 taken 1 times.
✗ Branch 11 not taken.
✓ Branch 13 taken 1 times.
✗ Branch 14 not taken.
✓ Branch 15 taken 1 times.
✗ Branch 16 not taken.
✗ Branch 17 not taken.
✓ Branch 18 taken 1 times.
✓ Branch 19 taken 1 times.
✗ Branch 20 not taken.
✓ Branch 21 taken 1 times.
✗ Branch 22 not taken.
✗ Branch 23 not taken.
✗ Branch 24 not taken.
✗ Branch 25 not taken.
✗ Branch 26 not taken.
✗ Branch 27 not taken.
✗ Branch 28 not taken.
✗ Branch 29 not taken.
✗ Branch 30 not taken.
✗ Branch 31 not taken.
✗ Branch 32 not taken.
|
4 | Parameters::init(argc, argv); |
67 | // [[/codeblock]] | ||
68 | |||
69 | // We define convenience aliases for the type tags of the two problems. The type | ||
70 | // tags contain all the properties that are needed to define the model and the problem | ||
71 | // setup. Throughout the main file, we will obtain types defined for these type tags | ||
72 | // using the property system, i.e. with `GetPropType`. A more detailed documentation | ||
73 | // for the type tags and properties can be found in the documentation of the single-phase | ||
74 | // and the tracer transport setups, respectively. | ||
75 | 1 | using OnePTypeTag = Properties::TTag::IncompressibleTest; | |
76 | 1 | using TracerTypeTag = Properties::TTag::TracerTest; | |
77 | |||
78 | // #### Step 1: Create the grid | ||
79 | // The `GridManager` class creates the grid from information given in the input file. | ||
80 | // This can either be a grid file, or in the case of structured grids, one can specify the coordinates | ||
81 | // of the corners of the grid and the number of cells to be used to discretize each spatial direction. | ||
82 | // Here, we solve both the single-phase and the tracer problem on the same grid, and thus, | ||
83 | // the grid is only created once using the grid type defined by the `OnePTypeTag` of the single-phase problem. | ||
84 | // [[codeblock]] | ||
85 |
1/2✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
|
2 | GridManager<GetPropType<OnePTypeTag, Properties::Grid>> gridManager; |
86 |
5/12✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 1 times.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✓ Branch 10 taken 1 times.
✓ Branch 12 taken 1 times.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
2 | gridManager.init(); |
87 | |||
88 | // We compute on the leaf grid view. | ||
89 |
2/4✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
|
1 | const auto& leafGridView = gridManager.grid().leafGridView(); |
90 | // [[/codeblock]] | ||
91 | |||
92 | // #### Step 2: Solving the single-phase problem | ||
93 | // First, a finite volume grid geometry is constructed from the grid that was created above. | ||
94 | // This builds the sub-control volumes (scv) and sub-control volume faces (scvf) for each element | ||
95 | // of the grid partition. | ||
96 | 1 | using GridGeometry = GetPropType<OnePTypeTag, Properties::GridGeometry>; | |
97 |
1/4✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
|
2 | auto gridGeometry = std::make_shared<GridGeometry>(leafGridView); |
98 | |||
99 | // We now instantiate the problem, in which we define the boundary and initial conditions. | ||
100 | 1 | using OnePProblem = GetPropType<OnePTypeTag, Properties::Problem>; | |
101 |
2/6✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 1 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
|
2 | auto problemOneP = std::make_shared<OnePProblem>(gridGeometry); |
102 | |||
103 | // The jacobian matrix (`A`), the solution vector (`p`) and the residual (`r`) make up the linear system. | ||
104 | 1 | using JacobianMatrix = GetPropType<OnePTypeTag, Properties::JacobianMatrix>; | |
105 | 1 | using SolutionVector = GetPropType<OnePTypeTag, Properties::SolutionVector>; | |
106 |
3/8✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 1 times.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
|
2 | SolutionVector p(leafGridView.size(0)); |
107 | |||
108 |
2/6✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 1 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
|
2 | auto A = std::make_shared<JacobianMatrix>(); |
109 |
2/6✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 1 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
|
2 | auto r = std::make_shared<SolutionVector>(); |
110 | |||
111 | // The grid variables are used store variables (primary and secondary variables) on sub-control volumes and faces (volume and flux variables). | ||
112 | 1 | using OnePGridVariables = GetPropType<OnePTypeTag, Properties::GridVariables>; | |
113 |
2/6✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 1 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
|
2 | auto onePGridVariables = std::make_shared<OnePGridVariables>(problemOneP, gridGeometry); |
114 |
2/4✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
|
2 | onePGridVariables->init(p); |
115 | |||
116 | // We now instantiate the assembler class, assemble the linear system and solve it with the linear | ||
117 | // solver AMGCGIstlSolver (a conjugate gradient solver preconditioned by an algebraic multigrid). | ||
118 | // Besides that, the time needed for assembly and solve is measured and printed. | ||
119 | 1 | using OnePAssembler = FVAssembler<OnePTypeTag, DiffMethod::analytic>; | |
120 |
2/4✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 1 times.
✗ Branch 4 not taken.
|
2 | auto assemblerOneP = std::make_shared<OnePAssembler>(problemOneP, gridGeometry, onePGridVariables); |
121 |
3/12✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 1 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 1 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
3 | assemblerOneP->setLinearSystem(A, r); // tell assembler to use our previously defined system |
122 | |||
123 | 1 | Dune::Timer timer; | |
124 |
1/2✓ Branch 3 taken 1 times.
✗ Branch 4 not taken.
|
2 | Dune::Timer assemblyTimer; std::cout << "Assembling linear system ..." << std::flush; |
125 |
2/4✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
|
2 | assemblerOneP->assembleJacobianAndResidual(p); // assemble linear system around current solution |
126 |
2/4✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✓ Branch 8 taken 1 times.
✗ Branch 9 not taken.
|
3 | assemblyTimer.stop(); std::cout << " took " << assemblyTimer.elapsed() << " seconds." << std::endl; |
127 | |||
128 | 2 | (*r) *= -1.0; // We want to solve `Ax = -r`. | |
129 | |||
130 | 1 | using OnePLinearSolver = AMGCGIstlSolver<LinearSolverTraits<GridGeometry>, LinearAlgebraTraitsFromAssembler<OnePAssembler>>; | |
131 |
1/2✓ Branch 3 taken 1 times.
✗ Branch 4 not taken.
|
2 | Dune::Timer solverTimer; std::cout << "Solving linear system ..." << std::flush; |
132 |
6/14✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 1 times.
✗ Branch 8 not taken.
✓ Branch 10 taken 1 times.
✗ Branch 11 not taken.
✓ Branch 13 taken 1 times.
✗ Branch 14 not taken.
✓ Branch 15 taken 1 times.
✗ Branch 16 not taken.
✗ Branch 17 not taken.
✗ Branch 18 not taken.
|
6 | auto onePLinearSolver = std::make_shared<OnePLinearSolver>(gridGeometry->gridView(), gridGeometry->dofMapper()); |
133 |
4/8✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 1 times.
✗ Branch 8 not taken.
✓ Branch 10 taken 1 times.
✗ Branch 11 not taken.
|
4 | onePLinearSolver->solve(*A, p, *r); |
134 |
2/4✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✓ Branch 8 taken 1 times.
✗ Branch 9 not taken.
|
3 | solverTimer.stop(); std::cout << " took " << solverTimer.elapsed() << " seconds." << std::endl; |
135 | |||
136 |
1/2✓ Branch 3 taken 1 times.
✗ Branch 4 not taken.
|
2 | Dune::Timer updateTimer; std::cout << "Updating variables ..." << std::flush; |
137 |
2/4✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
|
2 | onePGridVariables->update(p); // update grid variables to new pressure distribution |
138 |
2/4✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 1 times.
✗ Branch 8 not taken.
|
2 | updateTimer.elapsed(); std::cout << " took " << updateTimer.elapsed() << std::endl; |
139 | |||
140 | // The solution vector `p` now contains the pressure field, i.e.the solution to the single-phase | ||
141 | // problem defined in `problem_1p.hh`. Let us now write this solution to a VTK file using the Dune | ||
142 | // `VTKWriter`. Moreover, we add the permeability distribution to the writer. | ||
143 | // [[codeblock]] | ||
144 | 1 | using GridView = typename GridGeometry::GridView; | |
145 |
2/4✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 1 times.
✗ Branch 4 not taken.
|
2 | Dune::VTKWriter<GridView> onepWriter(leafGridView); |
146 |
5/12✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 1 times.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✓ Branch 10 taken 1 times.
✓ Branch 12 taken 1 times.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
2 | onepWriter.addCellData(p, "p"); |
147 | |||
148 |
3/6✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 1 times.
✗ Branch 8 not taken.
|
3 | const auto& k = problemOneP->spatialParams().getKField(); // defined in spatialparams_1p.hh |
149 |
5/12✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 1 times.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✓ Branch 10 taken 1 times.
✓ Branch 12 taken 1 times.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
2 | onepWriter.addCellData(k, "permeability"); // add permeability to writer |
150 |
5/12✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 1 times.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✓ Branch 10 taken 1 times.
✗ Branch 11 not taken.
✓ Branch 12 taken 1 times.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
2 | onepWriter.write("1p"); // write the file "1p.vtk" |
151 | |||
152 | // print overall CPU time required for assembling and solving the 1p problem. | ||
153 | 1 | timer.stop(); | |
154 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 1 times.
|
1 | const auto& comm = leafGridView.comm(); |
155 |
1/2✓ Branch 3 taken 1 times.
✗ Branch 4 not taken.
|
2 | std::cout << "Simulation took " << timer.elapsed() << " seconds on " |
156 | << comm.size() << " processes.\n" | ||
157 |
5/10✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✓ Branch 10 taken 1 times.
✗ Branch 11 not taken.
✓ Branch 13 taken 1 times.
✗ Branch 14 not taken.
✓ Branch 16 taken 1 times.
✗ Branch 17 not taken.
|
4 | << "The cumulative CPU time was " << timer.elapsed()*comm.size() << " seconds.\n"; |
158 | // [[/codeblock]] | ||
159 | |||
160 | // #### Step 3: Computation of the volume fluxes | ||
161 | // We use the results of the 1p problem to calculate the volume fluxes across all sub-control volume | ||
162 | // faces of the discretization and store them in the vector `volumeFlux`. In order to do so, we iterate | ||
163 | // over all elements of the grid, and in each element compute the volume fluxes for all sub-control volume | ||
164 | // faces embedded in that element. | ||
165 | // [[codeblock]] | ||
166 | 1 | using Scalar = GetPropType<OnePTypeTag, Properties::Scalar>; // type for scalar values | |
167 |
5/12✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 1 times.
✗ Branch 8 not taken.
✓ Branch 10 taken 1 times.
✗ Branch 11 not taken.
✓ Branch 13 taken 1 times.
✗ Branch 14 not taken.
✗ Branch 16 not taken.
✗ Branch 17 not taken.
|
5 | std::vector<Scalar> volumeFlux(gridGeometry->numScvf(), 0.0); |
168 | |||
169 | 1 | using FluxVariables = GetPropType<OnePTypeTag, Properties::FluxVariables>; | |
170 | ✗ | auto upwindTerm = [](const auto& volVars) { return volVars.mobility(0); }; | |
171 | |||
172 | // We iterate over all elements | ||
173 |
2/4✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 2501 times.
✗ Branch 4 not taken.
|
5001 | for (const auto& element : elements(leafGridView)) |
174 | { | ||
175 | // Compute the element-local views on geometry, primary and secondary variables | ||
176 | // as well as variables needed for flux computations. | ||
177 | |||
178 | // This creates instances of the local views | ||
179 | 5000 | auto fvGeometry = localView(*gridGeometry); | |
180 | 10000 | auto elemVolVars = localView(onePGridVariables->curGridVolVars()); | |
181 | 7500 | auto elemFluxVars = localView(onePGridVariables->gridFluxVarsCache()); | |
182 | |||
183 | // we now have to bind the views to the current element | ||
184 | 2500 | fvGeometry.bind(element); | |
185 |
1/2✓ Branch 1 taken 2500 times.
✗ Branch 2 not taken.
|
2500 | elemVolVars.bind(element, fvGeometry, p); |
186 | 2500 | elemFluxVars.bind(element, fvGeometry, elemVolVars); | |
187 | |||
188 | // We calculate the volume fluxes for all sub-control volume faces except for Neumann boundary faces | ||
189 |
6/6✓ Branch 0 taken 10000 times.
✓ Branch 1 taken 2500 times.
✓ Branch 2 taken 10000 times.
✓ Branch 3 taken 2500 times.
✓ Branch 4 taken 200 times.
✓ Branch 5 taken 9800 times.
|
15000 | for (const auto& scvf : scvfs(fvGeometry)) |
190 | { | ||
191 | // skip Neumann boundary faces | ||
192 |
10/10✓ Branch 0 taken 200 times.
✓ Branch 1 taken 9800 times.
✓ Branch 2 taken 150 times.
✓ Branch 3 taken 50 times.
✓ Branch 4 taken 150 times.
✓ Branch 5 taken 50 times.
✓ Branch 6 taken 100 times.
✓ Branch 7 taken 100 times.
✓ Branch 8 taken 100 times.
✓ Branch 9 taken 100 times.
|
10300 | if (scvf.boundary() && problemOneP->boundaryTypes(element, scvf).hasNeumann()) |
193 | 100 | continue; | |
194 | |||
195 | // let the FluxVariables class do the flux computation. | ||
196 |
1/2✓ Branch 1 taken 9900 times.
✗ Branch 2 not taken.
|
9900 | FluxVariables fluxVars; |
197 |
2/4✓ Branch 1 taken 9900 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 9900 times.
✗ Branch 5 not taken.
|
19800 | fluxVars.init(*problemOneP, element, fvGeometry, elemVolVars, scvf, elemFluxVars); |
198 |
1/2✓ Branch 1 taken 9900 times.
✗ Branch 2 not taken.
|
9900 | volumeFlux[scvf.index()] = fluxVars.advectiveFlux(0, upwindTerm); |
199 | } | ||
200 | } | ||
201 | // [[/codeblock]] | ||
202 | |||
203 | // #### Step 4: Solving the tracer transport problem | ||
204 | // First, we instantiate the tracer problem containing initial and boundary conditions, | ||
205 | // and pass to it the previously computed volume fluxes (see the documentation of the | ||
206 | // file `spatialparams_tracer.hh` for more details). | ||
207 | 1 | using TracerProblem = GetPropType<TracerTypeTag, Properties::Problem>; | |
208 |
2/6✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 1 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
|
2 | auto tracerProblem = std::make_shared<TracerProblem>(gridGeometry); |
209 |
3/6✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 1 times.
✗ Branch 8 not taken.
|
3 | tracerProblem->spatialParams().setVolumeFlux(volumeFlux); |
210 | |||
211 | // We create and initialize the solution vector. In contrast to the steady-state single-phase problem, | ||
212 | // the tracer problem is transient, and the initial solution defined in the problem is applied to the | ||
213 | // solution vector. On the basis of this solution, we initialize then the grid variables. | ||
214 |
3/8✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 1 times.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
|
2 | SolutionVector x(leafGridView.size(0)); |
215 |
2/4✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
|
2 | tracerProblem->applyInitialSolution(x); |
216 |
2/6✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 1 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
|
2 | auto xOld = x; |
217 | |||
218 | 1 | using GridVariables = GetPropType<TracerTypeTag, Properties::GridVariables>; | |
219 |
2/4✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 1 times.
✗ Branch 4 not taken.
|
2 | auto gridVariables = std::make_shared<GridVariables>(tracerProblem, gridGeometry); |
220 |
2/4✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
|
2 | gridVariables->init(x); |
221 | |||
222 | // Let us now instantiate the time loop. Therefore, we read in some time loop parameters from the input file. | ||
223 | // The parameter `tEnd` defines the duration of the simulation, `dt` the initial time step size and `maxDt` the maximal time step size. | ||
224 | // Moreover, we define 10 check points in the time loop at which we will write the solution to vtk files. | ||
225 | // [[codeblock]] | ||
226 |
1/2✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
|
1 | const auto tEnd = getParam<Scalar>("TimeLoop.TEnd"); |
227 |
1/2✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
|
1 | auto dt = getParam<Scalar>("TimeLoop.DtInitial"); |
228 |
1/2✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
|
1 | const auto maxDt = getParam<Scalar>("TimeLoop.MaxTimeStepSize"); |
229 | |||
230 | // We instantiate the time loop. | ||
231 |
2/8✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 1 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
|
2 | auto timeLoop = std::make_shared<CheckPointTimeLoop<Scalar>>(0.0, dt, tEnd); |
232 |
2/4✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
|
2 | timeLoop->setMaxTimeStepSize(maxDt); |
233 |
2/4✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
|
2 | timeLoop->setPeriodicCheckPoint(tEnd/10.0); |
234 | // [[/codeblock]] | ||
235 | |||
236 | // We create and initialize the assembler with a time loop for the transient problem. | ||
237 | // Within the time loop, we will use this assembler in each time step to assemble the linear system. | ||
238 | // As solver we use the ILUBiCGSTABIstlSolver (a bi-conjugate gradient solver preconditioned by an incomplete LU-factorization). | ||
239 | 1 | using TracerAssembler = FVAssembler<TracerTypeTag, DiffMethod::analytic, /*implicit=*/false>; | |
240 |
2/4✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 1 times.
✗ Branch 4 not taken.
|
2 | auto assembler = std::make_shared<TracerAssembler>(tracerProblem, gridGeometry, gridVariables, timeLoop, xOld); |
241 |
3/12✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 1 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 1 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
3 | assembler->setLinearSystem(A, r); |
242 | |||
243 | 1 | using TracerLinearSolver = ILUBiCGSTABIstlSolver<LinearSolverTraits<GridGeometry>, LinearAlgebraTraitsFromAssembler<TracerAssembler>>; | |
244 |
6/14✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 1 times.
✗ Branch 8 not taken.
✓ Branch 10 taken 1 times.
✗ Branch 11 not taken.
✓ Branch 13 taken 1 times.
✗ Branch 14 not taken.
✓ Branch 15 taken 1 times.
✗ Branch 16 not taken.
✗ Branch 17 not taken.
✗ Branch 18 not taken.
|
6 | auto tracerLinearSolver = std::make_shared<TracerLinearSolver>(gridGeometry->gridView(), gridGeometry->dofMapper()); |
245 | |||
246 | // The following lines of code initialize the vtk output module, add the velocity output facility | ||
247 | // and write out the initial solution. At each checkpoint, we will use the output module to write | ||
248 | // the solution of a time step into a corresponding vtk file. | ||
249 | // [[codeblock]] | ||
250 |
8/16✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 1 times.
✗ Branch 8 not taken.
✓ Branch 10 taken 1 times.
✗ Branch 11 not taken.
✓ Branch 13 taken 1 times.
✗ Branch 14 not taken.
✓ Branch 16 taken 1 times.
✗ Branch 17 not taken.
✓ Branch 19 taken 1 times.
✗ Branch 20 not taken.
✓ Branch 21 taken 1 times.
✗ Branch 22 not taken.
|
4 | VtkOutputModule<GridVariables, SolutionVector> vtkWriter(*gridVariables, x, tracerProblem->name()); |
251 | |||
252 | // add model-specific output fields to the writer | ||
253 | 1 | using IOFields = GetPropType<TracerTypeTag, Properties::IOFields>; | |
254 |
1/2✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
|
1 | IOFields::initOutputModule(vtkWriter); |
255 | |||
256 | // add velocity output facility | ||
257 | 1 | using VelocityOutput = GetPropType<TracerTypeTag, Properties::VelocityOutput>; | |
258 |
4/8✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✓ Branch 8 taken 1 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 1 times.
|
3 | vtkWriter.addVelocityOutput(std::make_shared<VelocityOutput>(*gridVariables)); |
259 | |||
260 | // write initial solution | ||
261 |
1/2✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
|
1 | vtkWriter.write(0.0); |
262 | // [[/codeblock]] | ||
263 | |||
264 | // ##### The time loop | ||
265 | // We start the time loop and solve a new time step as long as `tEnd` is not reached. In every time step, | ||
266 | // the problem is assembled and solved, the solution is updated, and when a checkpoint is reached the solution | ||
267 | // is written to a new vtk file. In addition, statistics related to CPU time, the current simulation time | ||
268 | // and the time step sizes used is printed to the terminal. | ||
269 | // [[codeblock]] | ||
270 |
2/4✓ Branch 0 taken 1 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
|
2 | timeLoop->start(); do |
271 | { | ||
272 | // Then the linear system is assembled. | ||
273 | 500 | Dune::Timer assembleTimer; | |
274 |
2/4✓ Branch 1 taken 500 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 500 times.
✗ Branch 5 not taken.
|
1000 | assembler->assembleJacobianAndResidual(x); |
275 | 500 | assembleTimer.stop(); | |
276 | |||
277 | // We solve the linear system `A(xOld-xNew) = r`. | ||
278 | 500 | Dune::Timer solveTimer; | |
279 |
2/6✓ Branch 1 taken 500 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 500 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
|
1000 | SolutionVector xDelta(x); |
280 |
4/8✓ Branch 1 taken 500 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 500 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 500 times.
✗ Branch 8 not taken.
✓ Branch 10 taken 500 times.
✗ Branch 11 not taken.
|
2000 | tracerLinearSolver->solve(*A, xDelta, *r); |
281 | 500 | solveTimer.stop(); | |
282 | |||
283 | // update the solution vector and the grid variables. | ||
284 | 500 | updateTimer.reset(); | |
285 | 500 | x -= xDelta; | |
286 |
2/4✓ Branch 1 taken 500 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 500 times.
✗ Branch 5 not taken.
|
1000 | gridVariables->update(x); |
287 | 500 | updateTimer.stop(); | |
288 | |||
289 | // display the statistics of the actual time step. | ||
290 | 500 | const auto elapsedTot = assembleTimer.elapsed() + solveTimer.elapsed() + updateTimer.elapsed(); | |
291 | 500 | std::cout << "Assemble/solve/update time: " | |
292 |
2/4✓ Branch 2 taken 500 times.
✗ Branch 3 not taken.
✓ Branch 7 taken 500 times.
✗ Branch 8 not taken.
|
1000 | << assembleTimer.elapsed() << "(" << 100*assembleTimer.elapsed()/elapsedTot << "%)/" |
293 |
2/4✓ Branch 2 taken 500 times.
✗ Branch 3 not taken.
✓ Branch 7 taken 500 times.
✗ Branch 8 not taken.
|
1000 | << solveTimer.elapsed() << "(" << 100*solveTimer.elapsed()/elapsedTot << "%)/" |
294 |
2/4✓ Branch 2 taken 500 times.
✗ Branch 3 not taken.
✓ Branch 7 taken 500 times.
✗ Branch 8 not taken.
|
1000 | << updateTimer.elapsed() << "(" << 100*updateTimer.elapsed()/elapsedTot << "%)" |
295 |
1/2✓ Branch 1 taken 500 times.
✗ Branch 2 not taken.
|
500 | << std::endl; |
296 | |||
297 | // Update the old solution with the one computed in this time step and move to the next one | ||
298 |
1/2✓ Branch 1 taken 500 times.
✗ Branch 2 not taken.
|
500 | xOld = x; |
299 |
2/4✓ Branch 1 taken 500 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 500 times.
✗ Branch 5 not taken.
|
1000 | gridVariables->advanceTimeStep(); |
300 |
2/4✓ Branch 1 taken 500 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 500 times.
✗ Branch 5 not taken.
|
1000 | timeLoop->advanceTimeStep(); |
301 | |||
302 | // Write the Vtk output on check points. | ||
303 |
4/4✓ Branch 0 taken 10 times.
✓ Branch 1 taken 490 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 490 times.
|
1000 | if (timeLoop->isCheckPoint()) |
304 |
3/6✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 10 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 10 times.
✗ Branch 8 not taken.
|
30 | vtkWriter.write(timeLoop->time()); |
305 | |||
306 | // report statistics of this time step | ||
307 |
2/4✓ Branch 1 taken 500 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 500 times.
✗ Branch 5 not taken.
|
1000 | timeLoop->reportTimeStep(); |
308 | |||
309 | // set the time step size for the next time step | ||
310 |
2/4✓ Branch 1 taken 500 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 500 times.
✗ Branch 5 not taken.
|
1000 | timeLoop->setTimeStepSize(dt); |
311 | |||
312 |
4/6✓ Branch 1 taken 500 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 500 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 499 times.
✓ Branch 7 taken 1 times.
|
1000 | } while (!timeLoop->finished()); |
313 | // [[/codeblock]] | ||
314 | |||
315 | // The following piece of code prints a final status report of the time loop | ||
316 | // before the program is terminated. | ||
317 |
3/6✗ Branch 0 not taken.
✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 1 times.
✓ Branch 5 taken 1 times.
✗ Branch 6 not taken.
|
2 | timeLoop->finalize(leafGridView.comm()); |
318 | |||
319 | 1 | return 0; | |
320 | } | ||
321 | // #### Exception handling | ||
322 | // In this part of the main file we catch and print possible exceptions that could | ||
323 | // occur during the simulation. | ||
324 | // [[details]] exception handler | ||
325 | // [[codeblock]] | ||
326 | // errors related to run-time parameters | ||
327 | ✗ | catch (Dumux::ParameterException &e) | |
328 | { | ||
329 | ✗ | std::cerr << std::endl << e << " ---> Abort!" << std::endl; | |
330 | ✗ | return 1; | |
331 | } | ||
332 | // errors related to the parsing of Dune grid files | ||
333 | ✗ | catch (Dune::DGFException & e) | |
334 | { | ||
335 | ✗ | std::cerr << "DGF exception thrown (" << e << | |
336 | "). Most likely, the DGF file name is wrong " | ||
337 | "or the DGF file is corrupted, " | ||
338 | "e.g. missing hash at end of file or wrong number (dimensions) of entries." | ||
339 | ✗ | << " ---> Abort!" << std::endl; | |
340 | ✗ | return 2; | |
341 | } | ||
342 | // generic error handling with Dune::Exception | ||
343 | ✗ | catch (Dune::Exception &e) | |
344 | { | ||
345 | ✗ | std::cerr << "Dune reported error: " << e << " ---> Abort!" << std::endl; | |
346 | ✗ | return 3; | |
347 | } | ||
348 | // other exceptions | ||
349 | ✗ | catch (...) | |
350 | { | ||
351 | ✗ | std::cerr << "Unknown exception thrown! ---> Abort!" << std::endl; | |
352 | ✗ | return 4; | |
353 | } | ||
354 | // [[/codeblock]] | ||
355 | // [[/details]] | ||
356 | // [[/content]] | ||
357 |