24#ifndef DUMUX_IO_GRID_MANAGER_YASP_HH
25#define DUMUX_IO_GRID_MANAGER_YASP_HH
27#include <dune/common/math.hh>
29#include <dune/grid/yaspgrid.hh>
30#include <dune/grid/io/file/dgfparser/dgfyasp.hh>
32#ifndef DUMUX_IO_GRID_MANAGER_BASE_HH
57template<
class Coordinates,
int dim>
62 using GlobalPosition = Dune::FieldVector<ct, dim>;
70 void init(
const std::string& modelParamGroup =
"")
75 ParentType::makeGridFromDgfFile(getParamFromGroup<std::string>(modelParamGroup,
"Grid.File"));
76 postProcessing_(modelParamGroup);
84 const auto upperRight = getParamFromGroup<GlobalPosition>(modelParamGroup,
"Grid.UpperRight");
87 std::array<int, dim> cells; cells.fill(1);
88 cells = getParamFromGroup<std::array<int, dim>>(modelParamGroup,
"Grid.Cells", cells);
92 const std::bitset<dim> periodic;
95 const int overlap = getParamFromGroup<int>(modelParamGroup,
"Grid.Overlap", 1);
97 if constexpr (std::is_same_v<Dune::EquidistantCoordinates<ct, dim>, Coordinates>)
98 init(upperRight, cells, modelParamGroup, overlap, periodic);
101 const auto lowerLeft = getParamFromGroup<GlobalPosition>(modelParamGroup,
"Grid.LowerLeft", GlobalPosition(0.0));
102 init(lowerLeft, upperRight, cells, modelParamGroup, overlap, periodic);
109 const auto prefix = modelParamGroup.empty() ? modelParamGroup : modelParamGroup +
".";
111 << prefix +
"Grid.UpperRight"
112 <<
", or a grid file in " << prefix +
"Grid.File");
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>{})
127 static_assert(std::is_same_v<Dune::EquidistantCoordinates<ct, dim>, Coordinates>,
128 "Use init function taking lowerLeft as argument when working with EquidistantOffsetCoordinates");
133 ParentType::gridPtr() = std::make_unique<Grid>(upperRight, cells, periodic, overlap);
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);
143 postProcessing_(modelParamGroup);
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>{})
157 static_assert(std::is_same_v<Dune::EquidistantOffsetCoordinates<ct, dim>, Coordinates>,
158 "LowerLeft can only be specified with EquidistantOffsetCoordinates");
163 ParentType::gridPtr() = std::make_unique<Grid>(lowerLeft, upperRight, cells, periodic, overlap);
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);
173 postProcessing_(modelParamGroup);
181 void postProcessing_(
const std::string& modelParamGroup)
184 const bool keepPhysicalOverlap = getParamFromGroup<bool>(modelParamGroup,
"Grid.KeepPhysicalOverlap",
true);
185 ParentType::grid().refineOptions(keepPhysicalOverlap);
186 ParentType::maybeRefineGrid(modelParamGroup);
187 ParentType::loadBalance();
220template<
class ctype,
int dim>
222:
public GridManagerBase<Dune::YaspGrid<dim, Dune::TensorProductCoordinates<ctype, dim> > >
231 void init(
const std::string& modelParamGroup =
"")
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));
241 std::array<std::vector<int>, dim> cells;
242 for (
int i = 0; i < dim; ++i)
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]);
249 std::array<std::vector<ctype>, dim> grading;
250 for (
int i = 0; i < dim; ++i)
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]);
257 init(positions, cells, grading, modelParamGroup);
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 =
"")
269 const int overlap = getParamFromGroup<int>(modelParamGroup,
"Grid.Overlap", 1);
270 const bool verbose = getParamFromGroup<bool>(modelParamGroup,
"Grid.Verbosity",
false);
273 const std::bitset<dim> periodic;
276 for (
unsigned int dimIdx = 0; dimIdx < dim; ++dimIdx)
278 if (cells[dimIdx].size() + 1 != positions[dimIdx].size())
280 DUNE_THROW(Dune::RangeError,
"Make sure to specify correct \"Cells\" and \"Positions\" arrays");
282 if (grading[dimIdx].size() + 1 != positions[dimIdx].size())
284 DUNE_THROW(Dune::RangeError,
"Make sure to specify correct \"Grading\" and \"Positions\" arrays");
286 ctype temp = std::numeric_limits<ctype>::lowest();
287 for (
unsigned int posIdx = 0; posIdx < positions[dimIdx].size(); ++posIdx)
289 if (temp > positions[dimIdx][posIdx])
291 DUNE_THROW(Dune::RangeError,
"Make sure to specify a monotone increasing \"Positions\" array");
293 temp = positions[dimIdx][posIdx];
297 const auto globalPositions = computeGlobalPositions_(positions, cells, grading, verbose);
303 ParentType::gridPtr() = std::make_shared<Grid>(globalPositions, periodic, overlap);
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);
313 postProcessing_(modelParamGroup);
320 void postProcessing_(
const std::string& modelParamGroup)
323 const bool keepPhysicalOverlap = getParamFromGroup<bool>(modelParamGroup,
"Grid.KeepPhysicalOverlap",
true);
324 ParentType::grid().refineOptions(keepPhysicalOverlap);
325 ParentType::maybeRefineGrid(modelParamGroup);
326 ParentType::loadBalance();
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)
336 std::array<std::vector<ctype>, dim> globalPositions;
339 for (
int dimIdx = 0; dimIdx < dim; dimIdx++)
341 for (
int zoneIdx = 0; zoneIdx < cells[dimIdx].size(); ++zoneIdx)
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;
349 bool increasingCellSize =
false;
353 std::cout <<
"dim " << dimIdx
354 <<
" lower " << lower
355 <<
" upper " << upper
356 <<
" numCells " << numCells
357 <<
" grading " << gradingFactor;
360 if (gradingFactor > 1.0)
362 increasingCellSize =
true;
367 if (gradingFactor < 0.0)
370 gradingFactor = abs(gradingFactor);
371 if (gradingFactor < 1.0)
373 increasingCellSize =
true;
378 if (gradingFactor > 1.0 - 1e-7 && gradingFactor < 1.0 + 1e-7)
380 height = 1.0 / numCells;
383 std::cout <<
" -> h " << height * length << std::endl;
389 height = (1.0 - gradingFactor) / (1.0 - power(gradingFactor, numCells));
393 std::cout <<
" -> grading_eff " << gradingFactor
394 <<
" h_min " << height * power(gradingFactor, 0) * length
395 <<
" h_max " << height * power(gradingFactor, numCells-1) * length
400 std::vector<ctype> localPositions;
401 localPositions.push_back(0);
402 for (
int i = 0; i < numCells-1; i++)
405 if (!(gradingFactor < 1.0 + 1e-7 && gradingFactor > 1.0 - 1e-7))
407 if (increasingCellSize)
409 hI *= power(gradingFactor, i);
413 hI *= power(gradingFactor, numCells-i-1);
416 localPositions.push_back(localPositions[i] + hI);
419 for (
int i = 0; i < localPositions.size(); i++)
421 localPositions[i] *= length;
422 localPositions[i] += lower;
425 for (
unsigned int i = 0; i < localPositions.size(); ++i)
427 globalPositions[dimIdx].push_back(localPositions[i]);
430 globalPositions[dimIdx].push_back(positions[dimIdx].back());
433 return globalPositions;
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:161
bool hasParamInGroup(const std::string ¶mGroup, const std::string ¶m)
Check whether a key exists in the parameter tree with a model group prefix.
Definition: parameters.hh:177
Definition: common/pdesolver.hh:36
Exception thrown if a run-time parameter is not specified correctly.
Definition: exceptions.hh:60
Definition: consistentlyorientedgrid.hh:32
The grid manager (this is the class used by the user) for all supported grid managers that constructs...
Definition: gridmanager_base.hh:324
The grid manager base interface (public) and methods common to most grid manager specializations (pro...
Definition: gridmanager_base.hh:67
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:75
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