version 3.8
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// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_IO_GRID_MANAGER_YASP_HH
13#define DUMUX_IO_GRID_MANAGER_YASP_HH
14
15#include <dune/common/math.hh>
16
17#include <dune/grid/yaspgrid.hh>
18#include <dune/grid/io/file/dgfparser/dgfyasp.hh>
19
20#ifndef DUMUX_IO_GRID_MANAGER_BASE_HH
22#endif
23
25
26namespace Dumux {
27
47template<class Coordinates, int dim>
48class GridManager<Dune::YaspGrid<dim, Coordinates>>
49: public GridManagerBase<Dune::YaspGrid<dim, Coordinates>>
50{
52 using GlobalPosition = Dune::FieldVector<ct, dim>;
53public:
56
60 void init(const std::string& modelParamGroup = "")
61 {
62 // First try to create it from a DGF file in GridParameterGroup.File
63 if (hasParamInGroup(modelParamGroup, "Grid.File"))
64 {
65 ParentType::makeGridFromDgfFile(getParamFromGroup<std::string>(modelParamGroup, "Grid.File"));
66 postProcessing_(modelParamGroup);
67 return;
68 }
69
70 // Then look for the necessary keys to construct from the input file
71 else if (hasParamInGroup(modelParamGroup, "Grid.UpperRight"))
72 {
73 // get the upper right corner coordinates
74 const auto upperRight = getParamFromGroup<GlobalPosition>(modelParamGroup, "Grid.UpperRight");
75
76 // number of cells in each direction
77 std::array<int, dim> cells; cells.fill(1);
78 cells = getParamFromGroup<std::array<int, dim>>(modelParamGroup, "Grid.Cells", cells);
79
80 // \todo TODO periodic boundaries with yasp (the periodicity concept of yasp grid is currently not supported, use dune-spgrid)
81 // const auto periodic = getParamFromGroup<std::bitset<dim>>(modelParamGroup, "Grid.Periodic", std::bitset<dim>{});
82 const std::bitset<dim> periodic;
83
84 // get the overlap
85 const int overlap = getParamFromGroup<int>(modelParamGroup, "Grid.Overlap", 1);
86
87 if constexpr (std::is_same_v<Dune::EquidistantCoordinates<ct, dim>, Coordinates>)
88 init(upperRight, cells, modelParamGroup, overlap, periodic);
89 else
90 {
91 const auto lowerLeft = getParamFromGroup<GlobalPosition>(modelParamGroup, "Grid.LowerLeft", GlobalPosition(0.0));
92 init(lowerLeft, upperRight, cells, modelParamGroup, overlap, periodic);
93 }
94 }
95
96 // Didn't find a way to construct the grid
97 else
98 {
99 const auto prefix = modelParamGroup.empty() ? modelParamGroup : modelParamGroup + ".";
100 DUNE_THROW(ParameterException, "Please supply one of the parameters "
101 << prefix + "Grid.UpperRight"
102 << ", or a grid file in " << prefix + "Grid.File");
103
104 }
105 }
106
111 void init(const GlobalPosition& upperRight,
112 const std::array<int, dim>& cells,
113 const std::string& modelParamGroup = "",
114 const int overlap = 1,
115 const std::bitset<dim> periodic = std::bitset<dim>{})
116 {
117 static_assert(std::is_same_v<Dune::EquidistantCoordinates<ct, dim>, Coordinates>,
118 "Use init function taking lowerLeft as argument when working with EquidistantOffsetCoordinates");
119
120 if (!hasParamInGroup(modelParamGroup, "Grid.Partitioning"))
121 {
122 // construct using default load balancing
123 ParentType::gridPtr() = std::make_unique<Grid>(upperRight, cells, periodic, overlap);
124 }
125 else
126 {
127 // construct using user defined partitioning
128 const auto partitioning = getParamFromGroup<std::array<int, dim>>(modelParamGroup, "Grid.Partitioning");
129 Dune::Yasp::FixedSizePartitioning<dim> lb(partitioning);
130 ParentType::gridPtr() = std::make_unique<Grid>(upperRight, cells, periodic, overlap, typename Grid::CollectiveCommunication(), &lb);
131 }
132
133 postProcessing_(modelParamGroup);
134 }
135
140 void init(const GlobalPosition& lowerLeft,
141 const GlobalPosition& upperRight,
142 const std::array<int, dim>& cells,
143 const std::string& modelParamGroup = "",
144 const int overlap = 1,
145 const std::bitset<dim> periodic = std::bitset<dim>{})
146 {
147 static_assert(std::is_same_v<Dune::EquidistantOffsetCoordinates<ct, dim>, Coordinates>,
148 "LowerLeft can only be specified with EquidistantOffsetCoordinates");
149
150 if (!hasParamInGroup(modelParamGroup, "Grid.Partitioning"))
151 {
152 // construct using default load balancing
153 ParentType::gridPtr() = std::make_unique<Grid>(lowerLeft, upperRight, cells, periodic, overlap);
154 }
155 else
156 {
157 // construct using user defined partitioning
158 const auto partitioning = getParamFromGroup<std::array<int, dim>>(modelParamGroup, "Grid.Partitioning");
159 Dune::Yasp::FixedSizePartitioning<dim> lb(partitioning);
160 ParentType::gridPtr() = std::make_unique<Grid>(lowerLeft, upperRight, cells, periodic, overlap, typename Grid::CollectiveCommunication(), &lb);
161 }
162
163 postProcessing_(modelParamGroup);
164 }
165
166private:
167
171 void postProcessing_(const std::string& modelParamGroup)
172 {
173 // Check if should refine the grid
174 const bool keepPhysicalOverlap = getParamFromGroup<bool>(modelParamGroup, "Grid.KeepPhysicalOverlap", true);
175 ParentType::grid().refineOptions(keepPhysicalOverlap);
176 ParentType::maybeRefineGrid(modelParamGroup);
177 ParentType::loadBalance();
178 }
179};
180
210template<class ctype, int dim>
211class GridManager<Dune::YaspGrid<dim, Dune::TensorProductCoordinates<ctype, dim> >>
212: public GridManagerBase<Dune::YaspGrid<dim, Dune::TensorProductCoordinates<ctype, dim> > >
213{
214public:
217
221 void init(const std::string& modelParamGroup = "")
222 {
223 // Only construction from the input file is possible
224 // Look for the necessary keys to construct from the input file
225 // The positions
226 std::array<std::vector<ctype>, dim> positions;
227 for (int i = 0; i < dim; ++i)
228 positions[i] = getParamFromGroup<std::vector<ctype>>(modelParamGroup, "Grid.Positions" + std::to_string(i));
229
230 // the number of cells (has a default)
231 std::array<std::vector<int>, dim> cells;
232 for (int i = 0; i < dim; ++i)
233 {
234 cells[i].resize(positions[i].size()-1, 1.0);
235 cells[i] = getParamFromGroup<std::vector<int>>(modelParamGroup, "Grid.Cells" + std::to_string(i), cells[i]);
236 }
237
238 // grading factor (has a default)
239 std::array<std::vector<ctype>, dim> grading;
240 for (int i = 0; i < dim; ++i)
241 {
242 grading[i].resize(positions[i].size()-1, 1.0);
243 grading[i] = getParamFromGroup<std::vector<ctype>>(modelParamGroup, "Grid.Grading" + std::to_string(i), grading[i]);
244 }
245
246 // call the generic function
247 init(positions, cells, grading, modelParamGroup);
248 }
249
253 void init(const std::array<std::vector<ctype>, dim>& positions,
254 const std::array<std::vector<int>, dim>& cells,
255 const std::array<std::vector<ctype>, dim>& grading,
256 const std::string& modelParamGroup = "")
257 {
258 // Additional parameters (they have a default)
259 const int overlap = getParamFromGroup<int>(modelParamGroup, "Grid.Overlap", 1);
260 const bool verbose = getParamFromGroup<bool>(modelParamGroup, "Grid.Verbosity", false);
261 // \todo TODO periodic boundaries with yasp (the periodicity concept of yasp grid is currently not supported, use dune-spgrid)
262 // const auto periodic = getParamFromGroup<std::bitset<dim>>(modelParamGroup, "Grid.Periodic", std::bitset<dim>{});
263 const std::bitset<dim> periodic;
264
265 // Some sanity checks
266 for (unsigned int dimIdx = 0; dimIdx < dim; ++dimIdx)
267 {
268 if (cells[dimIdx].size() + 1 != positions[dimIdx].size())
269 {
270 DUNE_THROW(Dune::RangeError, "Make sure to specify correct \"Cells\" and \"Positions\" arrays");
271 }
272 if (grading[dimIdx].size() + 1 != positions[dimIdx].size())
273 {
274 DUNE_THROW(Dune::RangeError, "Make sure to specify correct \"Grading\" and \"Positions\" arrays");
275 }
276 ctype temp = std::numeric_limits<ctype>::lowest();
277 for (unsigned int posIdx = 0; posIdx < positions[dimIdx].size(); ++posIdx)
278 {
279 if (temp > positions[dimIdx][posIdx])
280 {
281 DUNE_THROW(Dune::RangeError, "Make sure to specify a monotone increasing \"Positions\" array");
282 }
283 temp = positions[dimIdx][posIdx];
284 }
285 }
286
287 const auto globalPositions = computeGlobalPositions_(positions, cells, grading, verbose);
288
289 // make the grid
290 if (!hasParamInGroup(modelParamGroup, "Grid.Partitioning"))
291 {
292 // construct using default load balancing
293 ParentType::gridPtr() = std::make_shared<Grid>(globalPositions, periodic, overlap);
294 }
295 else
296 {
297 // construct using user defined partitioning
298 const auto partitioning = getParamFromGroup<std::array<int, dim>>(modelParamGroup, "Grid.Partitioning");
299 Dune::Yasp::FixedSizePartitioning<dim> lb(partitioning);
300 ParentType::gridPtr() = std::make_shared<Grid>(globalPositions, periodic, overlap, typename Grid::CollectiveCommunication(), &lb);
301 }
302
303 postProcessing_(modelParamGroup);
304 }
305
306private:
310 void postProcessing_(const std::string& modelParamGroup)
311 {
312 // Check if should refine the grid
313 const bool keepPhysicalOverlap = getParamFromGroup<bool>(modelParamGroup, "Grid.KeepPhysicalOverlap", true);
314 ParentType::grid().refineOptions(keepPhysicalOverlap);
315 ParentType::maybeRefineGrid(modelParamGroup);
316 ParentType::loadBalance();
317 }
318
320 std::array<std::vector<ctype>, dim>
321 computeGlobalPositions_(const std::array<std::vector<ctype>, dim>& positions,
322 const std::array<std::vector<int>, dim>& cells,
323 const std::array<std::vector<ctype>, dim>& grading,
324 bool verbose = false)
325 {
326 std::array<std::vector<ctype>, dim> globalPositions;
327 using std::pow;
328 using Dune::power;
329 for (int dimIdx = 0; dimIdx < dim; dimIdx++)
330 {
331 for (int zoneIdx = 0; zoneIdx < cells[dimIdx].size(); ++zoneIdx)
332 {
333 ctype lower = positions[dimIdx][zoneIdx];
334 ctype upper = positions[dimIdx][zoneIdx+1];
335 int numCells = cells[dimIdx][zoneIdx];
336 ctype gradingFactor = grading[dimIdx][zoneIdx];
337 ctype length = upper - lower;
338 ctype height = 1.0;
339 bool increasingCellSize = false;
340
341 if (verbose)
342 {
343 std::cout << "dim " << dimIdx
344 << " lower " << lower
345 << " upper " << upper
346 << " numCells " << numCells
347 << " grading " << gradingFactor;
348 }
349
350 if (gradingFactor > 1.0)
351 {
352 increasingCellSize = true;
353 }
354
355 // take absolute values and reverse cell size increment to achieve
356 // reverse behavior for negative values
357 if (gradingFactor < 0.0)
358 {
359 using std::abs;
360 gradingFactor = abs(gradingFactor);
361 if (gradingFactor < 1.0)
362 {
363 increasingCellSize = true;
364 }
365 }
366
367 // if the grading factor is exactly 1.0 do equal spacing
368 if (gradingFactor > 1.0 - 1e-7 && gradingFactor < 1.0 + 1e-7)
369 {
370 height = 1.0 / numCells;
371 if (verbose)
372 {
373 std::cout << " -> h " << height * length << std::endl;
374 }
375 }
376 // if grading factor is not 1.0, do power law spacing
377 else
378 {
379 height = (1.0 - gradingFactor) / (1.0 - power(gradingFactor, numCells));
380
381 if (verbose)
382 {
383 std::cout << " -> grading_eff " << gradingFactor
384 << " h_min " << height * power(gradingFactor, 0) * length
385 << " h_max " << height * power(gradingFactor, numCells-1) * length
386 << std::endl;
387 }
388 }
389
390 std::vector<ctype> localPositions;
391 localPositions.push_back(0);
392 for (int i = 0; i < numCells-1; i++)
393 {
394 ctype hI = height;
395 if (!(gradingFactor < 1.0 + 1e-7 && gradingFactor > 1.0 - 1e-7))
396 {
397 if (increasingCellSize)
398 {
399 hI *= power(gradingFactor, i);
400 }
401 else
402 {
403 hI *= power(gradingFactor, numCells-i-1);
404 }
405 }
406 localPositions.push_back(localPositions[i] + hI);
407 }
408
409 for (int i = 0; i < localPositions.size(); i++)
410 {
411 localPositions[i] *= length;
412 localPositions[i] += lower;
413 }
414
415 for (unsigned int i = 0; i < localPositions.size(); ++i)
416 {
417 globalPositions[dimIdx].push_back(localPositions[i]);
418 }
419 }
420 globalPositions[dimIdx].push_back(positions[dimIdx].back());
421 }
422
423 return globalPositions;
424 }
425};
426
427namespace Grid::Capabilities {
428
429// To the best of our knowledge YaspGrid is view thread-safe
430// This specialization can be removed after we depend on Dune release 2.9 and this is guaranteed by YaspGrid itself
431template<class Coordinates, int dim>
432struct MultithreadingSupported<Dune::YaspGrid<dim, Coordinates>>
433{
434 template<class GV>
435 static bool eval(const GV& gv) // default is independent of the grid view
436 { return true; }
437};
438
439} // end namespace Grid::Capabilities
440
441} // end namespace Dumux
442
443#endif
typename Dune::YaspGrid< dim, Coordinates > Grid
Definition: gridmanager_yasp.hh:54
void init(const std::string &modelParamGroup="")
Make the grid. This is implemented by specializations of this method.
Definition: gridmanager_yasp.hh:60
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:111
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:140
typename Dune::YaspGrid< dim, Dune::TensorProductCoordinates< ctype, dim > > Grid
Definition: gridmanager_yasp.hh:215
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:253
void init(const std::string &modelParamGroup="")
Make the grid. This is implemented by specializations of this method.
Definition: gridmanager_yasp.hh:221
The grid manager base interface (public) and methods common to most grid manager specializations (pro...
Definition: gridmanager_base.hh:55
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:63
The grid manager (this is the class used by the user) for all supported grid managers that constructs...
Definition: gridmanager_base.hh:323
Exception thrown if a run-time parameter is not specified correctly.
Definition: exceptions.hh:48
Definition: consistentlyorientedgrid.hh:20
dune-grid capabilities compatibility layer
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:149
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:165
Definition: adapt.hh:17
Definition: common/pdesolver.hh:24
static bool eval(const GV &gv)
Definition: gridmanager_yasp.hh:435
Definition: gridcapabilities.hh:67