3.4
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 if constexpr (std::is_same_v<Dune::EquidistantCoordinates<ct, dim>, Coordinates>)
98 init(upperRight, cells, modelParamGroup, overlap, periodic);
99 else
100 {
101 const auto lowerLeft = getParamFromGroup<GlobalPosition>(modelParamGroup, "Grid.LowerLeft", GlobalPosition(0.0));
102 init(lowerLeft, upperRight, cells, modelParamGroup, overlap, periodic);
103 }
104 }
105
106 // Didn't find a way to construct the grid
107 else
108 {
109 const auto prefix = modelParamGroup.empty() ? modelParamGroup : modelParamGroup + ".";
110 DUNE_THROW(ParameterException, "Please supply one of the parameters "
111 << prefix + "Grid.UpperRight"
112 << ", or a grid file in " << prefix + "Grid.File");
113
114 }
115 }
116
121 void init(const GlobalPosition& upperRight,
122 const std::array<int, dim>& cells,
123 const std::string& modelParamGroup = "",
124 const int overlap = 1,
125 const std::bitset<dim> periodic = std::bitset<dim>{})
126 {
127 static_assert(std::is_same_v<Dune::EquidistantCoordinates<ct, dim>, Coordinates>,
128 "Use init function taking lowerLeft as argument when working with EquidistantOffsetCoordinates");
129
130 if (!hasParamInGroup(modelParamGroup, "Grid.Partitioning"))
131 {
132 // construct using default load balancing
133 ParentType::gridPtr() = std::make_unique<Grid>(upperRight, cells, periodic, overlap);
134 }
135 else
136 {
137 // construct using user defined partitioning
138 const auto partitioning = getParamFromGroup<std::array<int, dim>>(modelParamGroup, "Grid.Partitioning");
139 Dune::YaspFixedSizePartitioner<dim> lb(partitioning);
140 ParentType::gridPtr() = std::make_unique<Grid>(upperRight, cells, periodic, overlap, typename Grid::CollectiveCommunicationType(), &lb);
141 }
142
143 postProcessing_(modelParamGroup);
144 }
145
150 void init(const GlobalPosition& lowerLeft,
151 const GlobalPosition& upperRight,
152 const std::array<int, dim>& cells,
153 const std::string& modelParamGroup = "",
154 const int overlap = 1,
155 const std::bitset<dim> periodic = std::bitset<dim>{})
156 {
157 static_assert(std::is_same_v<Dune::EquidistantOffsetCoordinates<ct, dim>, Coordinates>,
158 "LowerLeft can only be specified with EquidistantOffsetCoordinates");
159
160 if (!hasParamInGroup(modelParamGroup, "Grid.Partitioning"))
161 {
162 // construct using default load balancing
163 ParentType::gridPtr() = std::make_unique<Grid>(lowerLeft, upperRight, cells, periodic, overlap);
164 }
165 else
166 {
167 // construct using user defined partitioning
168 const auto partitioning = getParamFromGroup<std::array<int, dim>>(modelParamGroup, "Grid.Partitioning");
169 Dune::YaspFixedSizePartitioner<dim> lb(partitioning);
170 ParentType::gridPtr() = std::make_unique<Grid>(lowerLeft, upperRight, cells, periodic, overlap, typename Grid::CollectiveCommunicationType(), &lb);
171 }
172
173 postProcessing_(modelParamGroup);
174 }
175
176private:
177
181 void postProcessing_(const std::string& modelParamGroup)
182 {
183 // Check if should refine the grid
184 const bool keepPhysicalOverlap = getParamFromGroup<bool>(modelParamGroup, "Grid.KeepPhysicalOverlap", true);
185 ParentType::grid().refineOptions(keepPhysicalOverlap);
186 ParentType::maybeRefineGrid(modelParamGroup);
187 ParentType::loadBalance();
188 }
189};
190
220template<class ctype, int dim>
221class GridManager<Dune::YaspGrid<dim, Dune::TensorProductCoordinates<ctype, dim> >>
222: public GridManagerBase<Dune::YaspGrid<dim, Dune::TensorProductCoordinates<ctype, dim> > >
223{
224public:
225 using Grid = typename Dune::YaspGrid<dim, Dune::TensorProductCoordinates<ctype, dim> >;
227
231 void init(const std::string& modelParamGroup = "")
232 {
233 // Only construction from the input file is possible
234 // Look for the necessary keys to construct from the input file
235 // The positions
236 std::array<std::vector<ctype>, dim> positions;
237 for (int i = 0; i < dim; ++i)
238 positions[i] = getParamFromGroup<std::vector<ctype>>(modelParamGroup, "Grid.Positions" + std::to_string(i));
239
240 // the number of cells (has a default)
241 std::array<std::vector<int>, dim> cells;
242 for (int i = 0; i < dim; ++i)
243 {
244 cells[i].resize(positions[i].size()-1, 1.0);
245 cells[i] = getParamFromGroup<std::vector<int>>(modelParamGroup, "Grid.Cells" + std::to_string(i), cells[i]);
246 }
247
248 // grading factor (has a default)
249 std::array<std::vector<ctype>, dim> grading;
250 for (int i = 0; i < dim; ++i)
251 {
252 grading[i].resize(positions[i].size()-1, 1.0);
253 grading[i] = getParamFromGroup<std::vector<ctype>>(modelParamGroup, "Grid.Grading" + std::to_string(i), grading[i]);
254 }
255
256 // call the generic function
257 init(positions, cells, grading, modelParamGroup);
258 }
259
263 void init(const std::array<std::vector<ctype>, dim>& positions,
264 const std::array<std::vector<int>, dim>& cells,
265 const std::array<std::vector<ctype>, dim>& grading,
266 const std::string& modelParamGroup = "")
267 {
268 // Additional parameters (they have a default)
269 const int overlap = getParamFromGroup<int>(modelParamGroup, "Grid.Overlap", 1);
270 const bool verbose = getParamFromGroup<bool>(modelParamGroup, "Grid.Verbosity", false);
271 // \todo TODO periodic boundaries with yasp (the periodicity concept of yasp grid is currently not supported, use dune-spgrid)
272 // const auto periodic = getParamFromGroup<std::bitset<dim>>(modelParamGroup, "Grid.Periodic", std::bitset<dim>());
273 const std::bitset<dim> periodic;
274
275 // Some sanity checks
276 for (unsigned int dimIdx = 0; dimIdx < dim; ++dimIdx)
277 {
278 if (cells[dimIdx].size() + 1 != positions[dimIdx].size())
279 {
280 DUNE_THROW(Dune::RangeError, "Make sure to specify correct \"Cells\" and \"Positions\" arrays");
281 }
282 if (grading[dimIdx].size() + 1 != positions[dimIdx].size())
283 {
284 DUNE_THROW(Dune::RangeError, "Make sure to specify correct \"Grading\" and \"Positions\" arrays");
285 }
286 ctype temp = std::numeric_limits<ctype>::lowest();
287 for (unsigned int posIdx = 0; posIdx < positions[dimIdx].size(); ++posIdx)
288 {
289 if (temp > positions[dimIdx][posIdx])
290 {
291 DUNE_THROW(Dune::RangeError, "Make sure to specify a monotone increasing \"Positions\" array");
292 }
293 temp = positions[dimIdx][posIdx];
294 }
295 }
296
297 const auto globalPositions = computeGlobalPositions_(positions, cells, grading, verbose);
298
299 // make the grid
300 if (!hasParamInGroup(modelParamGroup, "Grid.Partitioning"))
301 {
302 // construct using default load balancing
303 ParentType::gridPtr() = std::make_shared<Grid>(globalPositions, periodic, overlap);
304 }
305 else
306 {
307 // construct using user defined partitioning
308 const auto partitioning = getParamFromGroup<std::array<int, dim>>(modelParamGroup, "Grid.Partitioning");
309 Dune::YaspFixedSizePartitioner<dim> lb(partitioning);
310 ParentType::gridPtr() = std::make_shared<Grid>(globalPositions, periodic, overlap, typename Grid::CollectiveCommunicationType(), &lb);
311 }
312
313 postProcessing_(modelParamGroup);
314 }
315
316private:
320 void postProcessing_(const std::string& modelParamGroup)
321 {
322 // Check if should refine the grid
323 const bool keepPhysicalOverlap = getParamFromGroup<bool>(modelParamGroup, "Grid.KeepPhysicalOverlap", true);
324 ParentType::grid().refineOptions(keepPhysicalOverlap);
325 ParentType::maybeRefineGrid(modelParamGroup);
326 ParentType::loadBalance();
327 }
328
330 std::array<std::vector<ctype>, dim>
331 computeGlobalPositions_(const std::array<std::vector<ctype>, dim>& positions,
332 const std::array<std::vector<int>, dim>& cells,
333 const std::array<std::vector<ctype>, dim>& grading,
334 bool verbose = false)
335 {
336 std::array<std::vector<ctype>, dim> globalPositions;
337 using std::pow;
338 using Dune::power;
339 for (int dimIdx = 0; dimIdx < dim; dimIdx++)
340 {
341 for (int zoneIdx = 0; zoneIdx < cells[dimIdx].size(); ++zoneIdx)
342 {
343 ctype lower = positions[dimIdx][zoneIdx];
344 ctype upper = positions[dimIdx][zoneIdx+1];
345 int numCells = cells[dimIdx][zoneIdx];
346 ctype gradingFactor = grading[dimIdx][zoneIdx];
347 ctype length = upper - lower;
348 ctype height = 1.0;
349 bool increasingCellSize = false;
350
351 if (verbose)
352 {
353 std::cout << "dim " << dimIdx
354 << " lower " << lower
355 << " upper " << upper
356 << " numCells " << numCells
357 << " grading " << gradingFactor;
358 }
359
360 if (gradingFactor > 1.0)
361 {
362 increasingCellSize = true;
363 }
364
365 // take absolute values and reverse cell size increment to achieve
366 // reverse behavior for negative values
367 if (gradingFactor < 0.0)
368 {
369 using std::abs;
370 gradingFactor = abs(gradingFactor);
371 if (gradingFactor < 1.0)
372 {
373 increasingCellSize = true;
374 }
375 }
376
377 // if the grading factor is exactly 1.0 do equal spacing
378 if (gradingFactor > 1.0 - 1e-7 && gradingFactor < 1.0 + 1e-7)
379 {
380 height = 1.0 / numCells;
381 if (verbose)
382 {
383 std::cout << " -> h " << height * length << std::endl;
384 }
385 }
386 // if grading factor is not 1.0, do power law spacing
387 else
388 {
389 height = (1.0 - gradingFactor) / (1.0 - power(gradingFactor, numCells));
390
391 if (verbose)
392 {
393 std::cout << " -> grading_eff " << gradingFactor
394 << " h_min " << height * power(gradingFactor, 0) * length
395 << " h_max " << height * power(gradingFactor, numCells-1) * length
396 << std::endl;
397 }
398 }
399
400 std::vector<ctype> localPositions;
401 localPositions.push_back(0);
402 for (int i = 0; i < numCells-1; i++)
403 {
404 ctype hI = height;
405 if (!(gradingFactor < 1.0 + 1e-7 && gradingFactor > 1.0 - 1e-7))
406 {
407 if (increasingCellSize)
408 {
409 hI *= power(gradingFactor, i);
410 }
411 else
412 {
413 hI *= power(gradingFactor, numCells-i-1);
414 }
415 }
416 localPositions.push_back(localPositions[i] + hI);
417 }
418
419 for (int i = 0; i < localPositions.size(); i++)
420 {
421 localPositions[i] *= length;
422 localPositions[i] += lower;
423 }
424
425 for (unsigned int i = 0; i < localPositions.size(); ++i)
426 {
427 globalPositions[dimIdx].push_back(localPositions[i]);
428 }
429 }
430 globalPositions[dimIdx].push_back(positions[dimIdx].back());
431 }
432
433 return globalPositions;
434 }
435};
436
437} // end namespace Dumux
438
439#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:358
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:374
Definition: adapt.hh:29
Definition: common/pdesolver.hh:36
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:322
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
void init(const GlobalPosition &upperRight, const std::array< int, dim > &cells, const std::string &modelParamGroup="", const int overlap=1, const std::bitset< dim > periodic=std::bitset< dim >{})
Make the grid using input data not read from the input file.
Definition: gridmanager_yasp.hh:121
void init(const GlobalPosition &lowerLeft, const GlobalPosition &upperRight, const std::array< int, dim > &cells, const std::string &modelParamGroup="", const int overlap=1, const std::bitset< dim > periodic=std::bitset< dim >{})
Make the grid using input data not read from the input file.
Definition: gridmanager_yasp.hh:150
typename Dune::YaspGrid< dim, Dune::TensorProductCoordinates< ctype, dim > > Grid
Definition: gridmanager_yasp.hh:225
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:263
void init(const std::string &modelParamGroup="")
Make the grid. This is implemented by specializations of this method.
Definition: gridmanager_yasp.hh:231