3.3.0
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
gridmanager_yasp.hh
Go to the documentation of this file.
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_IO_GRID_MANAGER_YASP_HH
25#define DUMUX_IO_GRID_MANAGER_YASP_HH
26
27#include <dune/common/math.hh>
28
29#include <dune/grid/yaspgrid.hh>
30#include <dune/grid/io/file/dgfparser/dgfyasp.hh>
31
32#ifndef DUMUX_IO_GRID_MANAGER_BASE_HH
34#endif
35
36namespace Dumux {
37
57template<class Coordinates, int dim>
58class GridManager<Dune::YaspGrid<dim, Coordinates>>
59: public GridManagerBase<Dune::YaspGrid<dim, Coordinates>>
60{
61 using ct = typename Dune::YaspGrid<dim, Coordinates>::ctype;
62 using GlobalPosition = Dune::FieldVector<ct, dim>;
63public:
64 using Grid = typename Dune::YaspGrid<dim, Coordinates>;
66
70 void init(const std::string& modelParamGroup = "")
71 {
72 // First try to create it from a DGF file in GridParameterGroup.File
73 if (hasParamInGroup(modelParamGroup, "Grid.File"))
74 {
75 ParentType::makeGridFromDgfFile(getParamFromGroup<std::string>(modelParamGroup, "Grid.File"));
76 postProcessing_(modelParamGroup);
77 return;
78 }
79
80 // Then look for the necessary keys to construct from the input file
81 else if (hasParamInGroup(modelParamGroup, "Grid.UpperRight"))
82 {
83 // get the upper right corner coordinates
84 const auto upperRight = getParamFromGroup<GlobalPosition>(modelParamGroup, "Grid.UpperRight");
85
86 // number of cells in each direction
87 std::array<int, dim> cells; cells.fill(1);
88 cells = getParamFromGroup<std::array<int, dim>>(modelParamGroup, "Grid.Cells", cells);
89
90 // \todo TODO periodic boundaries with yasp (the periodicity concept of yasp grid is currently not supported, use dune-spgrid)
91 // const auto periodic = getParamFromGroup<std::bitset<dim>>(modelParamGroup, "Grid.Periodic", std::bitset<dim>());
92 const std::bitset<dim> periodic;
93
94 // get the overlap
95 const int overlap = getParamFromGroup<int>(modelParamGroup, "Grid.Overlap", 1);
96
97 // make the grid
98 ParentType::gridPtr() = createGrid_(modelParamGroup, upperRight, cells, periodic, overlap, Coordinates{});
99
100 postProcessing_(modelParamGroup);
101 }
102
103 // Didn't find a way to construct the grid
104 else
105 {
106 const auto prefix = modelParamGroup.empty() ? modelParamGroup : modelParamGroup + ".";
107 DUNE_THROW(ParameterException, "Please supply one of the parameters "
108 << prefix + "Grid.UpperRight"
109 << ", or a grid file in " << prefix + "Grid.File");
110
111 }
112 }
113
114private:
118 std::unique_ptr<Grid> createGrid_(const std::string& modelParamGroup,
119 const GlobalPosition& upperRight,
120 const std::array<int, dim>& cells,
121 const std::bitset<dim>& periodic,
122 const int overlap,
123 Dune::EquidistantCoordinates<ct, dim>) const
124 {
125 if (!hasParamInGroup(modelParamGroup, "Grid.Partitioning"))
126 {
127 // construct using default load balancing
128 return std::make_unique<Grid>(upperRight, cells, periodic, overlap);
129 }
130 else
131 {
132 // construct using user defined partitioning
133 const auto partitioning = getParamFromGroup<std::array<int, dim>>(modelParamGroup, "Grid.Partitioning");
134 Dune::YaspFixedSizePartitioner<dim> lb(partitioning);
135 return std::make_unique<Grid>(upperRight, cells, periodic, overlap, typename Grid::CollectiveCommunicationType(), &lb);
136 }
137 }
138
142 std::unique_ptr<Grid> createGrid_(const std::string& modelParamGroup,
143 const GlobalPosition& upperRight,
144 const std::array<int, dim>& cells,
145 const std::bitset<dim>& periodic,
146 const int overlap,
147 Dune::EquidistantOffsetCoordinates<ct, dim>) const
148 {
149 const auto lowerLeft = getParamFromGroup<GlobalPosition>(modelParamGroup, "Grid.LowerLeft", GlobalPosition(0.0));
150
151 if (!hasParamInGroup(modelParamGroup, "Grid.Partitioning"))
152 {
153 // construct using default load balancing
154 return std::make_unique<Grid>(lowerLeft, upperRight, cells, periodic, overlap);
155 }
156 else
157 {
158 // construct using user defined partitioning
159 const auto partitioning = getParamFromGroup<std::array<int, dim>>(modelParamGroup, "Grid.Partitioning");
160 Dune::YaspFixedSizePartitioner<dim> lb(partitioning);
161 return std::make_unique<Grid>(lowerLeft, upperRight, cells, periodic, overlap, typename Grid::CollectiveCommunicationType(), &lb);
162 }
163 }
164
168 void postProcessing_(const std::string& modelParamGroup)
169 {
170 // Check if should refine the grid
171 const bool keepPhysicalOverlap = getParamFromGroup<bool>(modelParamGroup, "Grid.KeepPhysicalOverlap", true);
172 ParentType::grid().refineOptions(keepPhysicalOverlap);
173 ParentType::maybeRefineGrid(modelParamGroup);
174 ParentType::loadBalance();
175 }
176};
177
207template<class ctype, int dim>
208class GridManager<Dune::YaspGrid<dim, Dune::TensorProductCoordinates<ctype, dim> >>
209: public GridManagerBase<Dune::YaspGrid<dim, Dune::TensorProductCoordinates<ctype, dim> > >
210{
211public:
212 using Grid = typename Dune::YaspGrid<dim, Dune::TensorProductCoordinates<ctype, dim> >;
214
218 void init(const std::string& modelParamGroup = "")
219 {
220 // Only construction from the input file is possible
221 // Look for the necessary keys to construct from the input file
222 // The positions
223 std::array<std::vector<ctype>, dim> positions;
224 for (int i = 0; i < dim; ++i)
225 positions[i] = getParamFromGroup<std::vector<ctype>>(modelParamGroup, "Grid.Positions" + std::to_string(i));
226
227 // the number of cells (has a default)
228 std::array<std::vector<int>, dim> cells;
229 for (int i = 0; i < dim; ++i)
230 {
231 cells[i].resize(positions[i].size()-1, 1.0);
232 cells[i] = getParamFromGroup<std::vector<int>>(modelParamGroup, "Grid.Cells" + std::to_string(i), cells[i]);
233 }
234
235 // grading factor (has a default)
236 std::array<std::vector<ctype>, dim> grading;
237 for (int i = 0; i < dim; ++i)
238 {
239 grading[i].resize(positions[i].size()-1, 1.0);
240 grading[i] = getParamFromGroup<std::vector<ctype>>(modelParamGroup, "Grid.Grading" + std::to_string(i), grading[i]);
241 }
242
243 // call the generic function
244 init(positions, cells, grading, modelParamGroup);
245 }
246
250 void init(const std::array<std::vector<ctype>, dim>& positions,
251 const std::array<std::vector<int>, dim>& cells,
252 const std::array<std::vector<ctype>, dim>& grading,
253 const std::string& modelParamGroup = "")
254 {
255
256
257 // Additional arameters (they have a default)
258 const int overlap = getParamFromGroup<int>(modelParamGroup, "Grid.Overlap", 1);
259 const bool verbose = getParamFromGroup<bool>(modelParamGroup, "Grid.Verbosity", false);
260 // \todo TODO periodic boundaries with yasp (the periodicity concept of yasp grid is currently not supported, use dune-spgrid)
261 // const auto periodic = getParamFromGroup<std::bitset<dim>>(modelParamGroup, "Grid.Periodic", std::bitset<dim>());
262 const std::bitset<dim> periodic;
263
264 // Some sanity checks
265 for (unsigned int dimIdx = 0; dimIdx < dim; ++dimIdx)
266 {
267 if (cells[dimIdx].size() + 1 != positions[dimIdx].size())
268 {
269 DUNE_THROW(Dune::RangeError, "Make sure to specify correct \"Cells\" and \"Positions\" arrays");
270 }
271 if (grading[dimIdx].size() + 1 != positions[dimIdx].size())
272 {
273 DUNE_THROW(Dune::RangeError, "Make sure to specify correct \"Grading\" and \"Positions\" arrays");
274 }
275 ctype temp = std::numeric_limits<ctype>::lowest();
276 for (unsigned int posIdx = 0; posIdx < positions[dimIdx].size(); ++posIdx)
277 {
278 if (temp > positions[dimIdx][posIdx])
279 {
280 DUNE_THROW(Dune::RangeError, "Make sure to specify a monotone increasing \"Positions\" array");
281 }
282 temp = positions[dimIdx][posIdx];
283 }
284 }
285
286 const auto globalPositions = computeGlobalPositions_(positions, cells, grading, verbose);
287
288 // make the grid
289 if (!hasParamInGroup(modelParamGroup, "Grid.Partitioning"))
290 {
291 // construct using default load balancing
292 ParentType::gridPtr() = std::make_shared<Grid>(globalPositions, periodic, overlap);
293 }
294 else
295 {
296 // construct using user defined partitioning
297 const auto partitioning = getParamFromGroup<std::array<int, dim>>(modelParamGroup, "Grid.Partitioning");
298 Dune::YaspFixedSizePartitioner<dim> lb(partitioning);
299 ParentType::gridPtr() = std::make_shared<Grid>(globalPositions, periodic, overlap, typename Grid::CollectiveCommunicationType(), &lb);
300 }
301
302 postProcessing_(modelParamGroup);
303 }
304
305private:
309 void postProcessing_(const std::string& modelParamGroup)
310 {
311 // Check if should refine the grid
312 const bool keepPhysicalOverlap = getParamFromGroup<bool>(modelParamGroup, "Grid.KeepPhysicalOverlap", true);
313 ParentType::grid().refineOptions(keepPhysicalOverlap);
314 ParentType::maybeRefineGrid(modelParamGroup);
315 ParentType::loadBalance();
316 }
317
319 std::array<std::vector<ctype>, dim>
320 computeGlobalPositions_(const std::array<std::vector<ctype>, dim>& positions,
321 const std::array<std::vector<int>, dim>& cells,
322 const std::array<std::vector<ctype>, dim>& grading,
323 bool verbose = false)
324 {
325 std::array<std::vector<ctype>, dim> globalPositions;
326 using std::pow;
327 using Dune::power;
328 for (int dimIdx = 0; dimIdx < dim; dimIdx++)
329 {
330 for (int zoneIdx = 0; zoneIdx < cells[dimIdx].size(); ++zoneIdx)
331 {
332 ctype lower = positions[dimIdx][zoneIdx];
333 ctype upper = positions[dimIdx][zoneIdx+1];
334 int numCells = cells[dimIdx][zoneIdx];
335 ctype gradingFactor = grading[dimIdx][zoneIdx];
336 ctype length = upper - lower;
337 ctype height = 1.0;
338 bool increasingCellSize = false;
339
340 if (verbose)
341 {
342 std::cout << "dim " << dimIdx
343 << " lower " << lower
344 << " upper " << upper
345 << " numCells " << numCells
346 << " grading " << gradingFactor;
347 }
348
349 if (gradingFactor > 1.0)
350 {
351 increasingCellSize = true;
352 }
353
354 // take absolute values and reverse cell size increment to achieve
355 // reverse behavior for negative values
356 if (gradingFactor < 0.0)
357 {
358 using std::abs;
359 gradingFactor = abs(gradingFactor);
360 if (gradingFactor < 1.0)
361 {
362 increasingCellSize = true;
363 }
364 }
365
366 // if the grading factor is exactly 1.0 do equal spacing
367 if (gradingFactor > 1.0 - 1e-7 && gradingFactor < 1.0 + 1e-7)
368 {
369 height = 1.0 / numCells;
370 if (verbose)
371 {
372 std::cout << " -> h " << height * length << std::endl;
373 }
374 }
375 // if grading factor is not 1.0, do power law spacing
376 else
377 {
378 height = (1.0 - gradingFactor) / (1.0 - power(gradingFactor, numCells));
379
380 if (verbose)
381 {
382 std::cout << " -> grading_eff " << gradingFactor
383 << " h_min " << height * power(gradingFactor, 0) * length
384 << " h_max " << height * power(gradingFactor, numCells-1) * length
385 << std::endl;
386 }
387 }
388
389 std::vector<ctype> localPositions;
390 localPositions.push_back(0);
391 for (int i = 0; i < numCells-1; i++)
392 {
393 ctype hI = height;
394 if (!(gradingFactor < 1.0 + 1e-7 && gradingFactor > 1.0 - 1e-7))
395 {
396 if (increasingCellSize)
397 {
398 hI *= power(gradingFactor, i);
399 }
400 else
401 {
402 hI *= power(gradingFactor, numCells-i-1);
403 }
404 }
405 localPositions.push_back(localPositions[i] + hI);
406 }
407
408 for (int i = 0; i < localPositions.size(); i++)
409 {
410 localPositions[i] *= length;
411 localPositions[i] += lower;
412 }
413
414 for (unsigned int i = 0; i < localPositions.size(); ++i)
415 {
416 globalPositions[dimIdx].push_back(localPositions[i]);
417 }
418 }
419 globalPositions[dimIdx].push_back(positions[dimIdx].back());
420 }
421
422 return globalPositions;
423 }
424};
425
426} // end namespace Dumux
427
428#endif
Provides a grid manager for all supported grid managers with input file interfaces....
T getParamFromGroup(Args &&... args)
A free function to get a parameter from the parameter tree singleton with a model group.
Definition: parameters.hh:374
bool hasParamInGroup(const std::string &paramGroup, const std::string &param)
Check whether a key exists in the parameter tree with a model group prefix.
Definition: parameters.hh:390
Definition: adapt.hh:29
Definition: common/pdesolver.hh:35
Exception thrown if a run-time parameter is not specified correctly.
Definition: exceptions.hh:60
The grid manager (this is the class used by the user) for all supported grid managers that constructs...
Definition: gridmanager_base.hh:313
The grid manager base interface (public) and methods common to most grid manager specializations (pro...
Definition: gridmanager_base.hh:66
void init(const std::string &modelParamGroup="")
Make the grid. Implement this method in the specialization of this class for a grid type.
Definition: gridmanager_base.hh:74
typename Dune::YaspGrid< dim, Coordinates > Grid
Definition: gridmanager_yasp.hh:64
void init(const std::string &modelParamGroup="")
Make the grid. This is implemented by specializations of this method.
Definition: gridmanager_yasp.hh:70
typename Dune::YaspGrid< dim, Dune::TensorProductCoordinates< ctype, dim > > Grid
Definition: gridmanager_yasp.hh:212
void init(const std::array< std::vector< ctype >, dim > &positions, const std::array< std::vector< int >, dim > &cells, const std::array< std::vector< ctype >, dim > &grading, const std::string &modelParamGroup="")
Make the grid using input data not read from the input file.
Definition: gridmanager_yasp.hh:250
void init(const std::string &modelParamGroup="")
Make the grid. This is implemented by specializations of this method.
Definition: gridmanager_yasp.hh:218