version 3.11-dev
Loading...
Searching...
No Matches
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-FileCopyrightText: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
13#ifndef DUMUX_IO_GRID_MANAGER_YASP_HH
14#define DUMUX_IO_GRID_MANAGER_YASP_HH
15
16#include <dune/common/math.hh>
17
18#include <dune/grid/yaspgrid.hh>
19#include <dune/grid/io/file/dgfparser/dgfyasp.hh>
20
21#ifndef DUMUX_IO_GRID_MANAGER_BASE_HH
23#endif
24
26
27namespace Dumux {
28
49template<class Coordinates, int dim>
50class GridManager<Dune::YaspGrid<dim, Coordinates>>
51: public GridManagerBase<Dune::YaspGrid<dim, Coordinates>>
52{
54 using GlobalPosition = Dune::FieldVector<ct, dim>;
55public:
58
62 void init(const std::string& modelParamGroup = "")
63 {
64 // First, try to create it from a DGF or vti file in GridParameterGroup.File
65 if (hasParamInGroup(modelParamGroup, "Grid.File"))
66 {
67 const auto filename = getParamFromGroup<std::string>(modelParamGroup, "Grid.File");
68 const std::string extension = ParentType::getFileExtension(filename);
69
70 if (extension != "dgf" && extension != "vti")
71 DUNE_THROW(Dune::IOError, "Grid type " << Dune::className<Grid>() << " doesn't support grid files with extension: *."<< extension);
72
73 // VTK file formats for image grids
74 if (extension == "vti")
75 {
76#if DUMUX_HAVE_GRIDFORMAT
77 VTKReader vtkReader(filename);
78 VTKReader::Data cellData, pointData;
79
80 const bool verbose = getParamFromGroup<bool>(modelParamGroup, "Grid.Verbosity", false);
81 std::shared_ptr<Detail::VTKReader::ImageGrid<ct, dim>> gridInfo
82 = vtkReader.readStructuredGrid<ct, dim>(cellData, pointData, verbose);
83 const std::bitset<dim> periodic;
84 const int overlap = getParamFromGroup<int>(modelParamGroup, "Grid.Overlap", 1);
85 const auto upperRight = gridInfo->upperRight();
86 const auto cells = gridInfo->cells();
87 const auto lowerLeft = gridInfo->lowerLeft();
88
89 if constexpr (std::is_same_v<Dune::EquidistantCoordinates<ct, dim>, Coordinates>)
90 init_(upperRight, cells, modelParamGroup, overlap, periodic);
91 else
92 init_(lowerLeft, upperRight, cells, modelParamGroup, overlap, periodic);
93
94 ParentType::gridData_ = std::make_shared<GridData<Grid>>(
95 ParentType::gridPtr(), gridInfo, std::move(cellData), std::move(pointData)
96 );
97
99#else
100 DUNE_THROW(Dune::IOError, "Grid type " << Dune::className<Grid>() << " doesn't support grid files with extension: *."<< extension);
101#endif
102 }
103 else
104 {
106 getParamFromGroup<std::string>(modelParamGroup, "Grid.File")
107 );
108 }
109
110 postProcessing_(modelParamGroup);
111 return;
112 }
113
114 // Then look for the necessary keys to construct from the input file
115 else if (hasParamInGroup(modelParamGroup, "Grid.UpperRight"))
116 {
117 // get the upper right corner coordinates
118 const auto upperRight = getParamFromGroup<GlobalPosition>(modelParamGroup, "Grid.UpperRight");
119
120 // number of cells in each direction
121 std::array<int, dim> cells; cells.fill(1);
122 cells = getParamFromGroup<std::array<int, dim>>(modelParamGroup, "Grid.Cells", cells);
123
124 // \todo TODO periodic boundaries with yasp (the periodicity concept of yasp grid is currently not supported, use dune-spgrid)
125 // const auto periodic = getParamFromGroup<std::bitset<dim>>(modelParamGroup, "Grid.Periodic", std::bitset<dim>{});
126 const std::bitset<dim> periodic;
127
128 // get the overlap
129 const int overlap = getParamFromGroup<int>(modelParamGroup, "Grid.Overlap", 1);
130
131 if constexpr (std::is_same_v<Dune::EquidistantCoordinates<ct, dim>, Coordinates>)
132 init(upperRight, cells, modelParamGroup, overlap, periodic);
133 else
134 {
135 const auto lowerLeft = getParamFromGroup<GlobalPosition>(modelParamGroup, "Grid.LowerLeft", GlobalPosition(0.0));
136 init(lowerLeft, upperRight, cells, modelParamGroup, overlap, periodic);
137 }
138 }
139
140 // Didn't find a way to construct the grid
141 else
142 {
143 const auto prefix = modelParamGroup.empty() ? modelParamGroup : modelParamGroup + ".";
144 DUNE_THROW(ParameterException, "Please supply one of the parameters "
145 << prefix + "Grid.UpperRight"
146 << ", or a grid file in " << prefix + "Grid.File");
147
148 }
149 }
150
155 void init(const GlobalPosition& upperRight,
156 const std::array<int, dim>& cells,
157 const std::string& modelParamGroup = "",
158 const int overlap = 1,
159 const std::bitset<dim> periodic = std::bitset<dim>{})
160 {
161 init_(upperRight, cells, modelParamGroup, overlap, periodic);
162 postProcessing_(modelParamGroup);
163 }
164
169 void init(const GlobalPosition& lowerLeft,
170 const GlobalPosition& upperRight,
171 const std::array<int, dim>& cells,
172 const std::string& modelParamGroup = "",
173 const int overlap = 1,
174 const std::bitset<dim> periodic = std::bitset<dim>{})
175 {
176 init_(lowerLeft, upperRight, cells, modelParamGroup, overlap, periodic);
177 postProcessing_(modelParamGroup);
178 }
179
180private:
185 void init_(const GlobalPosition& upperRight,
186 const std::array<int, dim>& cells,
187 const std::string& modelParamGroup = "",
188 const int overlap = 1,
189 const std::bitset<dim> periodic = std::bitset<dim>{})
190 {
191 static_assert(std::is_same_v<Dune::EquidistantCoordinates<ct, dim>, Coordinates>,
192 "Use init function taking lowerLeft as argument when working with EquidistantOffsetCoordinates");
193
194 if (!hasParamInGroup(modelParamGroup, "Grid.Partitioning"))
195 {
196 // construct using default load balancing
197 ParentType::gridPtr() = std::make_unique<Grid>(upperRight, cells, periodic, overlap);
198 }
199 else
200 {
201 // construct using user defined partitioning
202 const auto partitioning = getParamFromGroup<std::array<int, dim>>(modelParamGroup, "Grid.Partitioning");
203 Dune::Yasp::FixedSizePartitioning<dim> lb(partitioning);
204 ParentType::gridPtr() = std::make_unique<Grid>(upperRight, cells, periodic, overlap, typename Grid::Communication(), &lb);
205 }
206 }
207
212 void init_(const GlobalPosition& lowerLeft,
213 const GlobalPosition& upperRight,
214 const std::array<int, dim>& cells,
215 const std::string& modelParamGroup = "",
216 const int overlap = 1,
217 const std::bitset<dim> periodic = std::bitset<dim>{})
218 {
219 static_assert(std::is_same_v<Dune::EquidistantOffsetCoordinates<ct, dim>, Coordinates>,
220 "LowerLeft can only be specified with EquidistantOffsetCoordinates");
221
222 if (!hasParamInGroup(modelParamGroup, "Grid.Partitioning"))
223 {
224 // construct using default load balancing
225 ParentType::gridPtr() = std::make_unique<Grid>(lowerLeft, upperRight, cells, periodic, overlap);
226 }
227 else
228 {
229 // construct using user defined partitioning
230 const auto partitioning = getParamFromGroup<std::array<int, dim>>(modelParamGroup, "Grid.Partitioning");
231 Dune::Yasp::FixedSizePartitioning<dim> lb(partitioning);
232 ParentType::gridPtr() = std::make_unique<Grid>(lowerLeft, upperRight, cells, periodic, overlap, typename Grid::Communication(), &lb);
233 }
234 }
235
239 void postProcessing_(const std::string& modelParamGroup)
240 {
241 // Check if should refine the grid
242 const bool keepPhysicalOverlap = getParamFromGroup<bool>(modelParamGroup, "Grid.KeepPhysicalOverlap", true);
243 ParentType::grid().refineOptions(keepPhysicalOverlap);
244 ParentType::maybeRefineGrid(modelParamGroup);
245 ParentType::loadBalance();
246 }
247};
248
278template<class ctype, int dim>
279class GridManager<Dune::YaspGrid<dim, Dune::TensorProductCoordinates<ctype, dim> >>
280: public GridManagerBase<Dune::YaspGrid<dim, Dune::TensorProductCoordinates<ctype, dim> > >
281{
282public:
285
289 void init(const std::string& modelParamGroup = "")
290 {
291 // Only construction from the input file is possible
292 // Look for the necessary keys to construct from the input file
293 // The positions
294 std::array<std::vector<ctype>, dim> positions;
295 for (int i = 0; i < dim; ++i)
296 positions[i] = getParamFromGroup<std::vector<ctype>>(modelParamGroup, "Grid.Positions" + std::to_string(i));
297
298 // the number of cells (has a default)
299 std::array<std::vector<int>, dim> cells;
300 for (int i = 0; i < dim; ++i)
301 {
302 cells[i].resize(positions[i].size()-1, 1.0);
303 cells[i] = getParamFromGroup<std::vector<int>>(modelParamGroup, "Grid.Cells" + std::to_string(i), cells[i]);
304 }
305
306 // grading factor (has a default)
307 std::array<std::vector<ctype>, dim> grading;
308 for (int i = 0; i < dim; ++i)
309 {
310 grading[i].resize(positions[i].size()-1, 1.0);
311 grading[i] = getParamFromGroup<std::vector<ctype>>(modelParamGroup, "Grid.Grading" + std::to_string(i), grading[i]);
312 }
313
314 // call the generic function
315 init(positions, cells, grading, modelParamGroup);
316 }
317
321 void init(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 const std::string& modelParamGroup = "")
325 {
326 // Additional parameters (they have a default)
327 const int overlap = getParamFromGroup<int>(modelParamGroup, "Grid.Overlap", 1);
328 const bool verbose = getParamFromGroup<bool>(modelParamGroup, "Grid.Verbosity", false);
329 // \todo TODO periodic boundaries with yasp (the periodicity concept of yasp grid is currently not supported, use dune-spgrid)
330 // const auto periodic = getParamFromGroup<std::bitset<dim>>(modelParamGroup, "Grid.Periodic", std::bitset<dim>{});
331 const std::bitset<dim> periodic;
332
333 // Some sanity checks
334 for (unsigned int dimIdx = 0; dimIdx < dim; ++dimIdx)
335 {
336 if (cells[dimIdx].size() + 1 != positions[dimIdx].size())
337 {
338 DUNE_THROW(Dune::RangeError, "Make sure to specify correct \"Cells\" and \"Positions\" arrays");
339 }
340 if (grading[dimIdx].size() + 1 != positions[dimIdx].size())
341 {
342 DUNE_THROW(Dune::RangeError, "Make sure to specify correct \"Grading\" and \"Positions\" arrays");
343 }
344 ctype temp = std::numeric_limits<ctype>::lowest();
345 for (unsigned int posIdx = 0; posIdx < positions[dimIdx].size(); ++posIdx)
346 {
347 if (temp > positions[dimIdx][posIdx])
348 {
349 DUNE_THROW(Dune::RangeError, "Make sure to specify a monotone increasing \"Positions\" array");
350 }
351 temp = positions[dimIdx][posIdx];
352 }
353 }
354
355 const auto globalPositions = computeGlobalPositions_(positions, cells, grading, verbose);
356
357 // make the grid
358 if (!hasParamInGroup(modelParamGroup, "Grid.Partitioning"))
359 {
360 // construct using default load balancing
361 ParentType::gridPtr() = std::make_shared<Grid>(globalPositions, periodic, overlap);
362 }
363 else
364 {
365 // construct using user defined partitioning
366 const auto partitioning = getParamFromGroup<std::array<int, dim>>(modelParamGroup, "Grid.Partitioning");
367 Dune::Yasp::FixedSizePartitioning<dim> lb(partitioning);
368 ParentType::gridPtr() = std::make_shared<Grid>(globalPositions, periodic, overlap, typename Grid::Communication(), &lb);
369 }
370
371 postProcessing_(modelParamGroup);
372 }
373
374private:
378 void postProcessing_(const std::string& modelParamGroup)
379 {
380 // Check if should refine the grid
381 const bool keepPhysicalOverlap = getParamFromGroup<bool>(modelParamGroup, "Grid.KeepPhysicalOverlap", true);
382 ParentType::grid().refineOptions(keepPhysicalOverlap);
383 ParentType::maybeRefineGrid(modelParamGroup);
384 ParentType::loadBalance();
385 }
386
388 std::array<std::vector<ctype>, dim>
389 computeGlobalPositions_(const std::array<std::vector<ctype>, dim>& positions,
390 const std::array<std::vector<int>, dim>& cells,
391 const std::array<std::vector<ctype>, dim>& grading,
392 bool verbose = false)
393 {
394 std::array<std::vector<ctype>, dim> globalPositions;
395 using std::pow;
396 using Dune::power;
397 for (int dimIdx = 0; dimIdx < dim; dimIdx++)
398 {
399 for (int zoneIdx = 0; zoneIdx < cells[dimIdx].size(); ++zoneIdx)
400 {
401 ctype lower = positions[dimIdx][zoneIdx];
402 ctype upper = positions[dimIdx][zoneIdx+1];
403 int numCells = cells[dimIdx][zoneIdx];
404 ctype gradingFactor = grading[dimIdx][zoneIdx];
405 ctype length = upper - lower;
406 ctype height = 1.0;
407 bool increasingCellSize = false;
408
409 if (verbose)
410 {
411 std::cout << "dim " << dimIdx
412 << " lower " << lower
413 << " upper " << upper
414 << " numCells " << numCells
415 << " grading " << gradingFactor;
416 }
417
418 if (gradingFactor > 1.0)
419 {
420 increasingCellSize = true;
421 }
422
423 // take absolute values and reverse cell size increment to achieve
424 // reverse behavior for negative values
425 if (gradingFactor < 0.0)
426 {
427 using std::abs;
428 gradingFactor = abs(gradingFactor);
429 if (gradingFactor < 1.0)
430 {
431 increasingCellSize = true;
432 }
433 }
434
435 // if the grading factor is exactly 1.0 do equal spacing
436 if (gradingFactor > 1.0 - 1e-7 && gradingFactor < 1.0 + 1e-7)
437 {
438 height = 1.0 / numCells;
439 if (verbose)
440 {
441 std::cout << " -> h " << height * length << std::endl;
442 }
443 }
444 // if grading factor is not 1.0, do power law spacing
445 else
446 {
447 height = (1.0 - gradingFactor) / (1.0 - power(gradingFactor, numCells));
448
449 if (verbose)
450 {
451 std::cout << " -> grading_eff " << gradingFactor
452 << " h_min " << height * power(gradingFactor, 0) * length
453 << " h_max " << height * power(gradingFactor, numCells-1) * length
454 << std::endl;
455 }
456 }
457
458 std::vector<ctype> localPositions;
459 localPositions.push_back(0);
460 for (int i = 0; i < numCells-1; i++)
461 {
462 ctype hI = height;
463 if (!(gradingFactor < 1.0 + 1e-7 && gradingFactor > 1.0 - 1e-7))
464 {
465 if (increasingCellSize)
466 {
467 hI *= power(gradingFactor, i);
468 }
469 else
470 {
471 hI *= power(gradingFactor, numCells-i-1);
472 }
473 }
474 localPositions.push_back(localPositions[i] + hI);
475 }
476
477 for (int i = 0; i < localPositions.size(); i++)
478 {
479 localPositions[i] *= length;
480 localPositions[i] += lower;
481 }
482
483 for (unsigned int i = 0; i < localPositions.size(); ++i)
484 {
485 globalPositions[dimIdx].push_back(localPositions[i]);
486 }
487 }
488 globalPositions[dimIdx].push_back(positions[dimIdx].back());
489 }
490
491 return globalPositions;
492 }
493};
494
495namespace Grid::Capabilities {
496
497// To the best of our knowledge YaspGrid is view thread-safe
498// This specialization can be removed after we depend on Dune release 2.9 and this is guaranteed by YaspGrid itself
499template<class Coordinates, int dim>
500struct MultithreadingSupported<Dune::YaspGrid<dim, Coordinates>>
501{
502 template<class GV>
503 static bool eval(const GV& gv) // default is independent of the grid view
504 { return true; }
505};
506
507} // end namespace Grid::Capabilities
508
509} // end namespace Dumux
510
511#endif
typename Dune::YaspGrid< dim, Coordinates > Grid
Definition gridmanager_yasp.hh:56
void init(const std::string &modelParamGroup="")
Make the grid. This is implemented by specializations of this method.
Definition gridmanager_yasp.hh:62
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:155
GridManagerBase< Grid > ParentType
Definition gridmanager_yasp.hh:57
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:169
GridManagerBase< Grid > ParentType
Definition gridmanager_yasp.hh:284
typename Dune::YaspGrid< dim, Dune::TensorProductCoordinates< ctype, dim > > Grid
Definition gridmanager_yasp.hh:283
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:321
void init(const std::string &modelParamGroup="")
Make the grid. This is implemented by specializations of this method.
Definition gridmanager_yasp.hh:289
The grid manager base interface (public) and methods common to most grid manager specializations (pro...
Definition gridmanager_base.hh:56
bool enableVtkData_
Definition gridmanager_base.hh:335
std::shared_ptr< GridData > gridData_
Definition gridmanager_base.hh:340
void makeGridFromDgfFile(const std::string &fileName)
Definition gridmanager_base.hh:265
std::shared_ptr< Grid > & gridPtr()
Definition gridmanager_base.hh:163
void init(const std::string &modelParamGroup="")
Definition gridmanager_base.hh:64
std::string getFileExtension(const std::string &fileName) const
Definition gridmanager_base.hh:185
The grid manager (this is the class used by the user) for all supported grid managers that constructs...
Definition gridmanager_base.hh:344
Exception thrown if a run-time parameter is not specified correctly.
Definition exceptions.hh:48
A vtk file reader using tinyxml2 as xml backend.
Definition vtkreader.hh:360
std::unordered_map< std::string, std::vector< double > > Data
the cell / point data type for point data read from a grid file
Definition vtkreader.hh:368
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 gridcapabilities.hh:17
Definition adapt.hh:17
Definition common/pdesolver.hh:24
static bool eval(const GV &gv)
Definition gridmanager_yasp.hh:503