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 | /*! | ||
8 | * \file | ||
9 | * \ingroup InputOutput | ||
10 | * \brief Provides a grid manager for all supported grid managers with | ||
11 | * input file interfaces. Manages data via the grid data member. | ||
12 | */ | ||
13 | #ifndef DUMUX_IO_GRID_MANAGER_BASE_HH | ||
14 | #define DUMUX_IO_GRID_MANAGER_BASE_HH | ||
15 | |||
16 | #include <array> | ||
17 | #include <bitset> | ||
18 | #include <memory> | ||
19 | #include <sstream> | ||
20 | |||
21 | #include <dune/common/exceptions.hh> | ||
22 | #include <dune/common/classname.hh> | ||
23 | #include <dune/common/parallel/communication.hh> | ||
24 | #include <dune/common/parallel/mpihelper.hh> | ||
25 | #include <dune/grid/io/file/dgfparser/dgfparser.hh> | ||
26 | #include <dune/grid/io/file/gmshreader.hh> | ||
27 | #include <dune/grid/common/gridfactory.hh> | ||
28 | #include <dune/grid/utility/structuredgridfactory.hh> | ||
29 | |||
30 | #include <dumux/common/typetraits/typetraits.hh> | ||
31 | #include <dumux/common/parameters.hh> | ||
32 | #include <dumux/discretization/method.hh> | ||
33 | #include <dumux/io/vtk/vtkreader.hh> | ||
34 | |||
35 | #include "griddata.hh" | ||
36 | |||
37 | namespace Dumux { | ||
38 | |||
39 | /*! | ||
40 | * \ingroup InputOutput | ||
41 | * \brief The grid manager (this is the class used by the user) for all supported grid managers that constructs a grid | ||
42 | * from information in the input file and handles the data. | ||
43 | * \note This class is specialised below for all supported grid managers. It inherits the functionality of the base class Dumux::GridManagerBase. | ||
44 | */ | ||
45 | template <class Grid> | ||
46 | class GridManager; | ||
47 | |||
48 | /*! | ||
49 | * \ingroup InputOutput | ||
50 | * \brief The grid manager base interface (public) and methods common | ||
51 | * to most grid manager specializations (protected). | ||
52 | */ | ||
53 | template <class GridType> | ||
54 | class GridManagerBase | ||
55 | { | ||
56 | public: | ||
57 | using Grid = GridType; | ||
58 | using GridData = Dumux::GridData<Grid>; | ||
59 | |||
60 | /*! | ||
61 | * \brief Make the grid. Implement this method in the specialization of this class for a grid type. | ||
62 | */ | ||
63 | void init(const std::string& modelParamGroup = "") | ||
64 | { | ||
65 | static_assert(AlwaysFalse<GridType>::value, | ||
66 | "The header with the GridManager specialization for your grid type is not included " | ||
67 | "or no specialization has been implemented!" | ||
68 | " In case of the latter, consider providing your own GridManager." | ||
69 | ); | ||
70 | } | ||
71 | |||
72 | /*! | ||
73 | * \brief Returns a reference to the grid. | ||
74 | */ | ||
75 | 81788 | Grid& grid() | |
76 | { | ||
77 |
2/2✓ Branch 0 taken 147 times.
✓ Branch 1 taken 41299 times.
|
81788 | if(enableDgfGridPointer_) |
78 | 159 | return *dgfGridPtr(); | |
79 | else | ||
80 | 81629 | return *gridPtr(); | |
81 | } | ||
82 | |||
83 | /*! | ||
84 | * \brief Returns a const reference to the grid. | ||
85 | */ | ||
86 | const Grid& grid() const | ||
87 | { | ||
88 | if(enableDgfGridPointer_) | ||
89 | return *dgfGridPtr(); | ||
90 | else | ||
91 | return *gridPtr(); | ||
92 | } | ||
93 | |||
94 | /*! | ||
95 | * \brief Call loadBalance() function of the grid. | ||
96 | */ | ||
97 | 650 | void loadBalance() | |
98 | { | ||
99 |
4/4✓ Branch 1 taken 54 times.
✓ Branch 2 taken 465 times.
✓ Branch 3 taken 54 times.
✓ Branch 4 taken 465 times.
|
650 | if (Dune::MPIHelper::instance().size() > 1) |
100 | { | ||
101 | // if we may have dgf parameters use load balancing of the dgf pointer | ||
102 |
2/2✓ Branch 0 taken 14 times.
✓ Branch 1 taken 40 times.
|
56 | if(enableDgfGridPointer_) |
103 | { | ||
104 | 14 | dgfGridPtr().loadBalance(); | |
105 | // update the grid data | ||
106 |
2/4✓ Branch 2 taken 14 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 14 times.
|
28 | gridData_ = std::make_shared<GridData>(dgfGridPtr()); |
107 | } | ||
108 | |||
109 | // if we have gmsh parameters we have to manually load balance the data | ||
110 |
2/2✓ Branch 0 taken 4 times.
✓ Branch 1 taken 36 times.
|
42 | else if (enableGmshDomainMarkers_) |
111 | { | ||
112 | // element and face markers are communicated during load balance | ||
113 | 8 | auto dh = gridData_->createGmshDataHandle(); | |
114 |
3/6✓ Branch 1 taken 4 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 4 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 4 times.
✗ Branch 8 not taken.
|
4 | gridPtr()->loadBalance(dh.interface()); |
115 |
3/6✓ Branch 1 taken 4 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 4 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 4 times.
✗ Branch 8 not taken.
|
4 | gridPtr()->communicate(dh.interface(), Dune::InteriorBorder_All_Interface, Dune::ForwardCommunication); |
116 | } | ||
117 | |||
118 | // if we have VTK parameters we have to manually load balance the data | ||
119 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 36 times.
|
38 | else if (enableVtkData_) |
120 | { | ||
121 | // cell and point data is communicated during load balance | ||
122 | ✗ | auto dh = gridData_->createVtkDataHandle(); | |
123 | ✗ | gridPtr()->loadBalance(dh.interface()); | |
124 | ✗ | gridPtr()->communicate(dh.interface(), Dune::InteriorBorder_All_Interface, Dune::ForwardCommunication); | |
125 | } | ||
126 | |||
127 | else | ||
128 | 64 | gridPtr()->loadBalance(); | |
129 | } | ||
130 | 650 | } | |
131 | |||
132 | /*! | ||
133 | * \brief Get an owning pointer to grid data associated with the grid | ||
134 | * \note Throws if no grid data is available | ||
135 | */ | ||
136 | 49 | std::shared_ptr<GridData> getGridData() const | |
137 | { | ||
138 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 49 times.
|
49 | if (!gridData_) |
139 | ✗ | DUNE_THROW(Dune::IOError, "No grid data available"); | |
140 | |||
141 | 49 | return gridData_; | |
142 | } | ||
143 | |||
144 | /*! | ||
145 | * \brief Check whether there is data associated with the grid | ||
146 | */ | ||
147 | bool hasGridData() const | ||
148 | { return static_cast<bool>(gridData_); } | ||
149 | |||
150 | protected: | ||
151 | |||
152 | /*! | ||
153 | * \brief Returns a reference to the grid pointer (std::shared_ptr<Grid>) | ||
154 | */ | ||
155 | 99904 | std::shared_ptr<Grid>& gridPtr() | |
156 | { | ||
157 |
1/2✓ Branch 0 taken 59389 times.
✗ Branch 1 not taken.
|
99904 | if(!enableDgfGridPointer_) |
158 | 99904 | return gridPtr_; | |
159 | else | ||
160 | ✗ | DUNE_THROW(Dune::InvalidStateException, "You are using DGF. To get the grid pointer use method dgfGridPtr()!"); | |
161 | } | ||
162 | |||
163 | /*! | ||
164 | * \brief Returns a reference to the DGF grid pointer (Dune::GridPtr<Grid>). | ||
165 | */ | ||
166 | 284 | Dune::GridPtr<Grid>& dgfGridPtr() | |
167 | { | ||
168 |
1/2✓ Branch 0 taken 266 times.
✗ Branch 1 not taken.
|
284 | if(enableDgfGridPointer_) |
169 | 284 | return dgfGridPtr_; | |
170 | else | ||
171 | ✗ | DUNE_THROW(Dune::InvalidStateException, "The DGF grid pointer is only available if the grid was constructed with a DGF file!"); | |
172 | } | ||
173 | |||
174 | /*! | ||
175 | * \brief Returns the filename extension of a given filename | ||
176 | */ | ||
177 | 119 | std::string getFileExtension(const std::string& fileName) const | |
178 | { | ||
179 |
1/2✓ Branch 0 taken 107 times.
✗ Branch 1 not taken.
|
119 | std::size_t i = fileName.rfind('.', fileName.length()); |
180 |
1/2✓ Branch 0 taken 107 times.
✗ Branch 1 not taken.
|
119 | if (i != std::string::npos) |
181 | { | ||
182 | 119 | return(fileName.substr(i+1, fileName.length() - i)); | |
183 | } | ||
184 | else | ||
185 | { | ||
186 | ✗ | DUNE_THROW(Dune::IOError, "Please provide and extension for your grid file ('"<< fileName << "')!"); | |
187 | } | ||
188 | return ""; | ||
189 | } | ||
190 | |||
191 | /*! | ||
192 | * \brief Makes a grid from a file. We currently support | ||
193 | * - dgf (Dune Grid Format) | ||
194 | * - msh (Gmsh mesh format) | ||
195 | * - vtp/vtu (VTK file formats) | ||
196 | */ | ||
197 | 64 | void makeGridFromFile(const std::string& fileName, | |
198 | const std::string& modelParamGroup) | ||
199 | { | ||
200 | // We found a file in the input file...does it have a supported extension? | ||
201 |
1/2✗ Branch 1 not taken.
✓ Branch 2 taken 63 times.
|
128 | const std::string extension = getFileExtension(fileName); |
202 |
3/8✓ Branch 1 taken 38 times.
✓ Branch 2 taken 25 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 38 times.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
|
64 | if (extension != "dgf" && extension != "msh" && extension != "vtu" && extension != "vtp") |
203 | ✗ | DUNE_THROW(Dune::IOError, "Grid type " << Dune::className<Grid>() << " doesn't support grid files with extension: *."<< extension); | |
204 | |||
205 | // Dune Grid Format (DGF) files | ||
206 |
2/2✓ Branch 0 taken 25 times.
✓ Branch 1 taken 38 times.
|
64 | if (extension == "dgf") |
207 | { | ||
208 | 26 | enableDgfGridPointer_ = true; | |
209 |
10/20✓ Branch 1 taken 25 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 25 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 25 times.
✗ Branch 8 not taken.
✓ Branch 10 taken 25 times.
✗ Branch 11 not taken.
✓ Branch 13 taken 25 times.
✗ Branch 14 not taken.
✓ Branch 16 taken 25 times.
✗ Branch 17 not taken.
✓ Branch 19 taken 25 times.
✗ Branch 20 not taken.
✓ Branch 22 taken 23 times.
✓ Branch 23 taken 2 times.
✓ Branch 25 taken 25 times.
✗ Branch 26 not taken.
✗ Branch 28 not taken.
✗ Branch 29 not taken.
|
128 | dgfGridPtr() = Dune::GridPtr<Grid>(fileName.c_str(), Dune::MPIHelper::getCommunicator()); |
210 |
3/6✓ Branch 1 taken 25 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✓ Branch 4 taken 25 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 25 times.
|
26 | gridData_ = std::make_shared<GridData>(dgfGridPtr_); |
211 | } | ||
212 | |||
213 | // Gmsh mesh format | ||
214 |
1/2✓ Branch 1 taken 38 times.
✗ Branch 2 not taken.
|
38 | else if (extension == "msh") |
215 | { | ||
216 | // get some optional parameters | ||
217 |
1/2✓ Branch 1 taken 38 times.
✗ Branch 2 not taken.
|
38 | const bool verbose = getParamFromGroup<bool>(modelParamGroup, "Grid.Verbosity", false); |
218 |
1/2✓ Branch 1 taken 38 times.
✗ Branch 2 not taken.
|
38 | const bool boundarySegments = getParamFromGroup<bool>(modelParamGroup, "Grid.BoundarySegments", false); |
219 |
1/4✓ Branch 1 taken 38 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
|
38 | const bool domainMarkers = getParamFromGroup<bool>(modelParamGroup, "Grid.DomainMarkers", false); |
220 | |||
221 |
2/2✓ Branch 0 taken 9 times.
✓ Branch 1 taken 29 times.
|
38 | if (domainMarkers) |
222 | { | ||
223 | 9 | enableGmshDomainMarkers_ = true; | |
224 |
4/10✓ Branch 1 taken 9 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 9 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 5 times.
✓ Branch 7 taken 4 times.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
|
32 | std::vector<int> boundaryMarkers, elementMarkers; |
225 |
3/4✓ Branch 1 taken 9 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 7 times.
✓ Branch 4 taken 2 times.
|
18 | auto gridFactory = std::make_unique<Dune::GridFactory<Grid>>(); |
226 | // the gmsh reader reads the grid on rank 0 | ||
227 |
2/4✓ Branch 1 taken 9 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 9 times.
✗ Branch 5 not taken.
|
18 | Dune::GmshReader<Grid>::read(*gridFactory, fileName, boundaryMarkers, elementMarkers, verbose, boundarySegments); |
228 |
7/20✓ Branch 1 taken 9 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 9 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 9 times.
✗ Branch 8 not taken.
✓ Branch 10 taken 9 times.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✓ Branch 13 taken 9 times.
✗ Branch 14 not taken.
✓ Branch 15 taken 9 times.
✗ Branch 16 not taken.
✓ Branch 17 taken 9 times.
✗ Branch 18 not taken.
✗ Branch 19 not taken.
✗ Branch 20 not taken.
✗ Branch 21 not taken.
✗ Branch 22 not taken.
✗ Branch 23 not taken.
|
18 | gridPtr() = std::shared_ptr<Grid>(gridFactory->createGrid()); |
229 |
4/8✓ Branch 1 taken 9 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✓ Branch 4 taken 9 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 9 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 9 times.
|
9 | gridData_ = std::make_shared<GridData>(gridPtr_, std::move(gridFactory), std::move(elementMarkers), std::move(boundaryMarkers)); |
230 | } | ||
231 | else | ||
232 | { | ||
233 |
1/2✓ Branch 1 taken 29 times.
✗ Branch 2 not taken.
|
58 | auto gridFactory = std::make_unique<Dune::GridFactory<Grid>>(); |
234 |
2/4✓ Branch 1 taken 29 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 29 times.
✗ Branch 5 not taken.
|
58 | Dune::GmshReader<Grid>::read(*gridFactory, fileName, verbose, boundarySegments); |
235 |
8/22✓ Branch 1 taken 29 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 29 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 29 times.
✗ Branch 8 not taken.
✓ Branch 10 taken 29 times.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✓ Branch 13 taken 29 times.
✗ Branch 14 not taken.
✓ Branch 15 taken 29 times.
✗ Branch 16 not taken.
✓ Branch 17 taken 29 times.
✓ Branch 18 taken 29 times.
✗ Branch 19 not taken.
✗ Branch 20 not taken.
✗ Branch 21 not taken.
✗ Branch 22 not taken.
✗ Branch 23 not taken.
✗ Branch 24 not taken.
✗ Branch 25 not taken.
|
58 | gridPtr() = std::shared_ptr<Grid>(gridFactory->createGrid()); |
236 | } | ||
237 | } | ||
238 | |||
239 | // VTK file formats for unstructured grids | ||
240 | ✗ | else if (extension == "vtu" || extension == "vtp") | |
241 | { | ||
242 | ✗ | if (Dune::MPIHelper::getCommunication().size() > 1) | |
243 | ✗ | DUNE_THROW(Dune::NotImplemented, "Reading grids in parallel from VTK file formats is currently not supported!"); | |
244 | |||
245 | ✗ | VTKReader vtkReader(fileName); | |
246 | ✗ | VTKReader::Data cellData, pointData; | |
247 | ✗ | auto gridFactory = std::make_unique<Dune::GridFactory<Grid>>(); | |
248 | ✗ | const bool verbose = getParamFromGroup<bool>(modelParamGroup, "Grid.Verbosity", false); | |
249 | ✗ | gridPtr() = vtkReader.readGrid(*gridFactory, cellData, pointData, verbose); | |
250 | ✗ | gridData_ = std::make_shared<GridData>(gridPtr_, std::move(gridFactory), std::move(cellData), std::move(pointData)); | |
251 | } | ||
252 | 64 | } | |
253 | |||
254 | /*! | ||
255 | * \brief Makes a grid from a DGF file. This is used by grid managers that only support DGF. | ||
256 | */ | ||
257 | 6 | void makeGridFromDgfFile(const std::string& fileName) | |
258 | { | ||
259 | // We found a file in the input file...does it have a supported extension? | ||
260 |
0/2✗ Branch 1 not taken.
✗ Branch 2 not taken.
|
6 | const std::string extension = getFileExtension(fileName); |
261 |
1/2✗ Branch 1 not taken.
✓ Branch 2 taken 6 times.
|
6 | if(extension != "dgf") |
262 | ✗ | DUNE_THROW(Dune::IOError, "Grid type " << Dune::className<Grid>() << " only supports DGF (*.dgf) but the specified filename has extension: *."<< extension); | |
263 | |||
264 | 6 | enableDgfGridPointer_ = true; | |
265 |
9/20✓ Branch 1 taken 6 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 6 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 6 times.
✗ Branch 8 not taken.
✓ Branch 10 taken 6 times.
✗ Branch 11 not taken.
✓ Branch 13 taken 6 times.
✗ Branch 14 not taken.
✓ Branch 16 taken 6 times.
✗ Branch 17 not taken.
✓ Branch 19 taken 6 times.
✗ Branch 20 not taken.
✗ Branch 22 not taken.
✓ Branch 23 taken 6 times.
✓ Branch 25 taken 6 times.
✗ Branch 26 not taken.
✗ Branch 28 not taken.
✗ Branch 29 not taken.
|
24 | dgfGridPtr() = Dune::GridPtr<Grid>(fileName.c_str(), Dune::MPIHelper::getCommunicator()); |
266 |
4/8✓ Branch 1 taken 6 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✓ Branch 4 taken 6 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 6 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 6 times.
|
6 | gridData_ = std::make_shared<GridData>(dgfGridPtr_); |
267 | 6 | } | |
268 | |||
269 | /*! | ||
270 | * \brief The cell types for structured grids | ||
271 | */ | ||
272 | enum CellType {Simplex, Cube}; | ||
273 | |||
274 | /*! | ||
275 | * \brief Makes a structured cube grid using the structured grid factory | ||
276 | */ | ||
277 | template <int dim, int dimworld> | ||
278 | 71 | void makeStructuredGrid(CellType cellType, | |
279 | const std::string& modelParamGroup) | ||
280 | { | ||
281 | using GlobalPosition = Dune::FieldVector<typename Grid::ctype, dimworld>; | ||
282 | 71 | const auto upperRight = getParamFromGroup<GlobalPosition>(modelParamGroup, "Grid.UpperRight"); | |
283 | 142 | const auto lowerLeft = getParamFromGroup<GlobalPosition>(modelParamGroup, "Grid.LowerLeft", GlobalPosition(0.0)); | |
284 | |||
285 | using CellArray = std::array<unsigned int, dim>; | ||
286 | 71 | CellArray cells; cells.fill(1); | |
287 | 71 | cells = getParamFromGroup<CellArray>(modelParamGroup, "Grid.Cells", cells); | |
288 | |||
289 | // make the grid | ||
290 |
2/2✓ Branch 0 taken 60 times.
✓ Branch 1 taken 11 times.
|
71 | if (cellType == CellType::Cube) |
291 | { | ||
292 |
6/10✓ Branch 2 taken 42 times.
✓ Branch 3 taken 18 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 42 times.
✓ Branch 6 taken 18 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 42 times.
✓ Branch 9 taken 18 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
|
78 | gridPtr() = Dune::StructuredGridFactory<Grid>::createCubeGrid(lowerLeft, upperRight, cells); |
293 | } | ||
294 |
1/2✓ Branch 0 taken 11 times.
✗ Branch 1 not taken.
|
11 | else if (cellType == CellType::Simplex) |
295 | { | ||
296 |
6/10✓ Branch 2 taken 2 times.
✓ Branch 3 taken 9 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✓ Branch 6 taken 9 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✓ Branch 9 taken 9 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
|
20 | gridPtr() = Dune::StructuredGridFactory<Grid>::createSimplexGrid(lowerLeft, upperRight, cells); |
297 | } | ||
298 | else | ||
299 | { | ||
300 | ✗ | DUNE_THROW(Dune::GridError, "Unknown cell type for making structured grid! Choose Cube or Simplex."); | |
301 | } | ||
302 | 71 | } | |
303 | |||
304 | /*! | ||
305 | * \brief Refines a grid after construction if GridParameterGroup.Refinement is set in the input file | ||
306 | */ | ||
307 | 730 | void maybeRefineGrid(const std::string& modelParamGroup) | |
308 | { | ||
309 |
8/14✓ Branch 1 taken 597 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 597 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 597 times.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✓ Branch 10 taken 597 times.
✓ Branch 11 taken 134 times.
✓ Branch 12 taken 463 times.
✓ Branch 13 taken 134 times.
✓ Branch 14 taken 463 times.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
|
2190 | if (hasParamInGroup(modelParamGroup, "Grid.Refinement")) |
310 | 182 | grid().globalRefine(getParamFromGroup<int>(modelParamGroup, "Grid.Refinement")); | |
311 | 730 | } | |
312 | |||
313 | /*! | ||
314 | * \brief A state variable if the DGF Dune::GridPtr has been enabled. | ||
315 | * It is always enabled if a DGF grid file was used to create the grid. | ||
316 | */ | ||
317 | bool enableDgfGridPointer_ = false; | ||
318 | |||
319 | /*! | ||
320 | * \brief A state variable if domain markers have been read from a Gmsh file. | ||
321 | */ | ||
322 | bool enableGmshDomainMarkers_ = false; | ||
323 | |||
324 | /*! | ||
325 | * \brief A state variable if cell or point data have been read from a VTK file. | ||
326 | */ | ||
327 | bool enableVtkData_ = false; | ||
328 | |||
329 | std::shared_ptr<Grid> gridPtr_; | ||
330 | Dune::GridPtr<Grid> dgfGridPtr_; | ||
331 | |||
332 | std::shared_ptr<GridData> gridData_; | ||
333 | }; | ||
334 | |||
335 | template <class Grid> | ||
336 | ✗ | class GridManager : public GridManagerBase<Grid> {}; | |
337 | |||
338 | } // end namespace Dumux | ||
339 | |||
340 | #endif | ||
341 |