24#ifndef DUMUX_CAKE_GRID_MANAGER_HH
25#define DUMUX_CAKE_GRID_MANAGER_HH
32#include <dune/common/dynvector.hh>
33#include <dune/common/float_cmp.hh>
34#include <dune/common/math.hh>
35#include <dune/grid/common/gridfactory.hh>
49 using Scalar =
typename Grid::ctype;
51 using GridFactory = Dune::GridFactory<Grid>;
52 using GridPointer = std::shared_ptr<Grid>;
54 enum { dim = Grid::dimension,
55 dimWorld = Grid::dimensionworld };
61 void init(
const std::string& modelParamGroup =
"")
63 static_assert(dim == 2 || dim == 3,
"The CakeGridManager is only implemented for 2D and 3D.");
65 const bool verbose = getParamFromGroup<bool>(modelParamGroup,
"Grid.Verbosity",
false);
67 std::array<std::vector<Scalar>, dim> polarCoordinates;
69 Dune::FieldVector<int, dim> indices(-1);
70 createVectors(polarCoordinates, indices, modelParamGroup, verbose);
96 static void createVectors(std::array<std::vector<Scalar>, dim> &polarCoordinates,
97 Dune::FieldVector<int, dim> &indices,
98 const std::string& modelParamGroup,
102 std::array<std::vector<Scalar>, dim> positions;
104 for (
int i = 0; i < dim; ++i)
106 const bool hasRadial =
hasParamInGroup(modelParamGroup,
"Grid.Radial" + std::to_string(i));
107 const bool hasAngular =
hasParamInGroup(modelParamGroup,
"Grid.Angular" + std::to_string(i));
108 const bool hasAxial = (dim == 3) &&
hasParamInGroup(modelParamGroup,
"Grid.Axial" + std::to_string(i));
109 if (
static_cast<int>(hasRadial) +
static_cast<int>(hasAngular) +
static_cast<int>(hasAxial) != 1)
110 DUNE_THROW(Dune::RangeError,
"Multiple or no position vectors (radial, angular, axial) specified in coord direction: " << i << std::endl);
114 positions[i] = getParamFromGroup<std::vector<Scalar>>(modelParamGroup,
"Grid.Radial" + std::to_string(i));
119 positions[i] = getParamFromGroup<std::vector<Scalar>>(modelParamGroup,
"Grid.Angular" + std::to_string(i));
124 positions[i] = getParamFromGroup<std::vector<Scalar>>(modelParamGroup,
"Grid.Axial" + std::to_string(i));
128 if (!std::is_sorted(positions[i].begin(), positions[i].end()))
129 DUNE_THROW(Dune::GridError,
"Make sure to specify a monotone increasing \"Positions\" array");
131 if(positions[i].size() < 2)
132 DUNE_THROW(Dune::GridError,
"Make sure to specify position arrays with at least two entries (min and max value).");
136 for (
int i = 0; i < dim; ++i)
138 DUNE_THROW(Dune::RangeError,
"Please specify Positions Angular and Radial and Axial correctly and unambiguously!" << std::endl);
141 std::array<std::vector<int>, dim> cells;
142 for (
int i = 0; i < dim; ++i)
144 cells[i].resize(positions[i].size()-1, 1.0);
145 cells[i] = getParamFromGroup<std::vector<int>>(modelParamGroup,
"Grid.Cells" + std::to_string(i), cells[i]);
146 if (cells[i].size() + 1 != positions[i].size())
147 DUNE_THROW(Dune::RangeError,
"Make sure to specify the correct length of \"Cells\" and \"Positions\" arrays");
151 std::array<std::vector<Scalar>, dim> grading;
152 for (
int i = 0; i < dim; ++i)
154 grading[i].resize(positions[i].size()-1, 1.0);
155 grading[i] = getParamFromGroup<std::vector<Scalar>>(modelParamGroup,
"Grid.Grading" + std::to_string(i), grading[i]);
156 if (grading[i].size() + 1 != positions[i].size())
157 DUNE_THROW(Dune::RangeError,
"Make sure to specify the correct length of \"Grading\" and \"Positions\" arrays");
161 std::array<std::vector<Scalar>, dim> globalPositions;
164 for (
int dimIdx = 0; dimIdx < dim; dimIdx++)
167 std::size_t numGlobalPositions = 1;
168 for (
int zoneIdx = 0; zoneIdx < cells[dimIdx].size(); ++zoneIdx)
169 numGlobalPositions += cells[dimIdx][zoneIdx];
171 globalPositions[dimIdx].resize(numGlobalPositions);
172 std::size_t posIdx = 0;
173 for (
int zoneIdx = 0; zoneIdx < cells[dimIdx].size(); ++zoneIdx)
175 const Scalar lower = positions[dimIdx][zoneIdx];
176 const Scalar upper = positions[dimIdx][zoneIdx+1];
177 const int numCells = cells[dimIdx][zoneIdx];
178 Scalar gradingFactor = grading[dimIdx][zoneIdx];
179 const Scalar length = upper - lower;
181 bool increasingCellSize =
false;
185 std::cout <<
"dim " << dimIdx
186 <<
" lower " << lower
187 <<
" upper " << upper
188 <<
" numCells " << numCells
189 <<
" grading " << gradingFactor;
192 if (gradingFactor > 1.0)
193 increasingCellSize =
true;
197 if (gradingFactor < 0.0)
200 gradingFactor = abs(gradingFactor);
202 if (gradingFactor < 1.0)
203 increasingCellSize =
true;
206 const bool useGrading = Dune::FloatCmp::eq(gradingFactor, 1.0) ? false :
true;
211 height = 1.0 / numCells;
213 std::cout <<
" -> h " << height * length << std::endl;
218 height = (1.0 - gradingFactor) / (1.0 - power(gradingFactor, numCells));
222 std::cout <<
" -> grading_eff " << gradingFactor
223 <<
" h_min " << height * power(gradingFactor, 0) * length
224 <<
" h_max " << height * power(gradingFactor, numCells-1) * length
230 Dune::DynamicVector<Scalar> localPositions(numCells, 0.0);
231 for (
int i = 1; i < numCells; i++)
236 if (increasingCellSize)
237 hI *= power(gradingFactor, i-1);
240 hI *= power(gradingFactor, numCells-i);
242 localPositions[i] = localPositions[i-1] + hI;
245 localPositions *= length;
246 localPositions += lower;
248 for (
int i = 0; i < numCells; ++i)
249 globalPositions[dimIdx][posIdx++] = localPositions[i];
252 globalPositions[dimIdx][posIdx] = positions[dimIdx].back();
255 polarCoordinates[0] = globalPositions[indices[0]];
256 polarCoordinates[1] = globalPositions[indices[1]];
258 polarCoordinates[2] = globalPositions[dim - indices[0] - indices[1]];
261 std::transform(polarCoordinates[1].begin(), polarCoordinates[1].end(),
262 polarCoordinates[1].begin(),
263 [](Scalar s){
return s*M_PI/180; });
274 std::unique_ptr<Grid>
createCakeGrid(std::array<std::vector<Scalar>, dim> &polarCoordinates,
275 Dune::FieldVector<int, dim> &indices,
276 const std::string& modelParamGroup,
277 bool verbose =
false)
279 std::vector<Scalar> dR = polarCoordinates[0];
280 std::vector<Scalar> dA = polarCoordinates[1];
286 int maxdA = dA.size() - 1;
287 if (Dune::FloatCmp::eq(dA[dA.size()-1], 2*M_PI))
289 maxdA = dA.size() - 2;
292 GridFactory gridFactory;
293 constexpr auto type = Dune::GeometryTypes::cube(dim);
296 if(dR[0] < 1.0e-8*dR.back())
302 constexpr auto prismType = Dune::GeometryTypes::prism;
303 std::vector<Scalar> dZ = polarCoordinates[2];
304 for (
int j = 0; j <= maxdA; ++j)
306 for (
int l = 0; l <= dZ.size() - 1; ++l)
308 for (
int i = hasHole ? 0 : 1; i <= dR.size()- 1; ++i)
313 Dune::FieldVector <double, dim> v(0.0);
314 v[indices[2]] = dZ[l];
315 v[indices[0]] = cos(dA[j])*dR[i];
316 v[indices[1]] = sin(dA[j])*dR[i];
319 gridFactory.insertVertex(v);
328 for (
int l = 0; l <= dZ.size() - 1; ++l)
330 Dune::FieldVector <double, dim> v(0.0);
331 v[indices[2]] = dZ[l];
337 gridFactory.insertVertex(v);
342 std::cout <<
"Filled node vector" << std::endl;
347 unsigned int rSize = hasHole ? dR.size() : dR.size()-1;
348 unsigned int zSize = dZ.size();
349 for (
int j = 0; j < dA.size() - 1; ++j)
351 for (
int l = 0; l < dZ.size() - 1; ++l)
355 for (
int i = 0; i < dR.size() - 1; ++i)
359 std::vector<unsigned int> vid({rSize*zSize*(maxdA+1) + l, z, z+rSize*zSize,
360 rSize*zSize*(maxdA+1) + l+1, z+rSize, z+rSize*zSize+rSize});
365 gridFactory.insertElement(prismType, vid);
369 std::vector<unsigned int> vid({z, z+1,
370 z+rSize*zSize, z+rSize*zSize+1,
372 z+rSize*zSize+rSize, z+rSize*zSize+rSize+1});
377 gridFactory.insertElement(type, vid);
387 for (
int i = 0; i < dR.size() - 1; ++i)
391 std::vector<unsigned int> vid({rSize*zSize*(maxdA+1) + l, z, t,
392 rSize*zSize*(maxdA+1) + l+1, z+rSize, t+rSize});
397 gridFactory.insertElement(prismType, vid);
401 std::vector<unsigned int> vid({z, z+1,
404 t+rSize, t+rSize+1});
409 gridFactory.insertElement(type, vid);
418 std::cout <<
"assign nodes 360° ends..." << std::endl;
427 constexpr auto triangleType = Dune::GeometryTypes::simplex(dim);
428 for (
int j = 0; j <= maxdA; ++j)
430 for (
int i = hasHole ? 0 : 1; i <= dR.size()- 1; ++i)
433 Dune::FieldVector <double, dim> v(0.0);
435 v[indices[0]] = cos(dA[j])*dR[i];
436 v[indices[1]] = sin(dA[j])*dR[i];
439 gridFactory.insertVertex(v);
447 Dune::FieldVector <double, dim> v(0.0);
453 gridFactory.insertVertex(v);
456 std::cout <<
"Filled node vector" << std::endl;
461 unsigned int rSize = hasHole ? dR.size() : dR.size()-1;
462 for (
int j = 0; j < dA.size() - 1; ++j)
466 for (
int i = 0; i < dR.size() - 1; ++i)
470 std::vector<unsigned int> vid({rSize*(maxdA+1), z, z+rSize});
475 gridFactory.insertElement(triangleType, vid);
479 std::vector<unsigned int> vid({z, z+1, z+rSize, z+rSize+1});
484 gridFactory.insertElement(type, vid);
493 for (
int i = 0; i < dR.size() - 1; ++i)
497 std::vector<unsigned int> vid({rSize*(maxdA+1), z, t});
502 gridFactory.insertElement(triangleType, vid);
506 std::vector<unsigned int> vid({z, z+1, t, t+1});
511 gridFactory.insertElement(type, vid);
520 std::cout <<
"assign nodes 360 ends..." << std::endl;
525 return std::unique_ptr<Grid>(gridFactory.createGrid());
548 std::cout <<
"Coordinates of : ";
549 for (
int k = 0; k < v.size(); ++k)
550 std::cout << v[k] <<
" ";
551 std::cout << std::endl;
556 std::cout <<
"element vertex indices: ";
557 for (
int k = 0; k < vid.size(); ++k)
558 std::cout << vid[k] <<
" ";
559 std::cout << std::endl;
571 GridPointer cakeGrid_;
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
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:390
Provides a grid manager with a method for creating creating vectors with polar Coordinates and one fo...
Definition: cakegridmanager.hh:48
void init(const std::string &modelParamGroup="")
Make the grid.
Definition: cakegridmanager.hh:61
void loadBalance()
Distributes the grid on all processes of a parallel computation.
Definition: cakegridmanager.hh:540
GridPointer & gridPtr()
Returns a reference to the shared pointer to the grid.
Definition: cakegridmanager.hh:565
static void printCoordinate(const Dune::FieldVector< double, dim > &v)
Definition: cakegridmanager.hh:546
static void printIndices(const std::vector< unsigned int > &vid)
Definition: cakegridmanager.hh:554
Grid & grid()
Returns a reference to the grid.
Definition: cakegridmanager.hh:531
static void createVectors(std::array< std::vector< Scalar >, dim > &polarCoordinates, Dune::FieldVector< int, dim > &indices, const std::string &modelParamGroup, bool verbose=false)
Create vectors containing polar coordinates of all points.
Definition: cakegridmanager.hh:96
std::unique_ptr< Grid > createCakeGrid(std::array< std::vector< Scalar >, dim > &polarCoordinates, Dune::FieldVector< int, dim > &indices, const std::string &modelParamGroup, bool verbose=false)
Creates Cartesian grid from polar coordinates.
Definition: cakegridmanager.hh:274