version 3.10-dev
cakegridmanager.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_CAKE_GRID_MANAGER_HH
13#define DUMUX_CAKE_GRID_MANAGER_HH
14
15#include <vector>
16#include <iostream>
17#include <cmath>
18#include <algorithm>
19
20#include <dune/common/dynvector.hh>
21#include <dune/common/float_cmp.hh>
22#include <dune/common/math.hh>
23#include <dune/grid/common/gridfactory.hh>
25
26namespace Dumux {
27
34template <class Grid>
36{
37 using Scalar = typename Grid::ctype;
38
39 using GridFactory = Dune::GridFactory<Grid>;
40 using GridPointer = std::shared_ptr<Grid>;
41
42 enum { dim = Grid::dimension,
43 dimWorld = Grid::dimensionworld };
44
45public:
49 void init(const std::string& modelParamGroup = "")
50 {
51 static_assert(dim == 2 || dim == 3, "The CakeGridManager is only implemented for 2D and 3D.");
52
53 const bool verbose = getParamFromGroup<bool>(modelParamGroup, "Grid.Verbosity", false);
54
55 std::array<std::vector<Scalar>, dim> polarCoordinates;
56 // Indices specifying in which direction the piece of cake is oriented
57 Dune::FieldVector<int, dim> indices(-1);
58 createVectors(polarCoordinates, indices, modelParamGroup, verbose);
59
60 gridPtr() = createCakeGrid(polarCoordinates, indices, modelParamGroup, verbose);
62 }
63
84 static void createVectors(std::array<std::vector<Scalar>, dim> &polarCoordinates,
85 Dune::FieldVector<int, dim> &indices,
86 const std::string& modelParamGroup,
87 bool verbose = false)
88 {
89 // The positions
90 std::array<std::vector<Scalar>, dim> positions;
91
92 for (int i = 0; i < dim; ++i)
93 {
94 const bool hasRadial = hasParamInGroup(modelParamGroup, "Grid.Radial" + std::to_string(i));
95 const bool hasAngular = hasParamInGroup(modelParamGroup, "Grid.Angular" + std::to_string(i));
96 const bool hasAxial = (dim == 3) && hasParamInGroup(modelParamGroup, "Grid.Axial" + std::to_string(i));
97 if (static_cast<int>(hasRadial) + static_cast<int>(hasAngular) + static_cast<int>(hasAxial) != 1)
98 DUNE_THROW(Dune::RangeError, "Multiple or no position vectors (radial, angular, axial) specified in coord direction: " << i << std::endl);
99
100 if (hasRadial)
101 {
102 positions[i] = getParamFromGroup<std::vector<Scalar>>(modelParamGroup, "Grid.Radial" + std::to_string(i));
103 indices[0] = i; // Index specifying radial direction
104 }
105 else if (hasAngular)
106 {
107 positions[i] = getParamFromGroup<std::vector<Scalar>>(modelParamGroup, "Grid.Angular" + std::to_string(i));
108 indices[1] = i; // Index specifying angular direction
109 }
110 else // hasAxial
111 {
112 positions[i] = getParamFromGroup<std::vector<Scalar>>(modelParamGroup, "Grid.Axial" + std::to_string(i));
113 indices[2] = i; // Index specifying axial direction
114 }
115
116 if (!std::is_sorted(positions[i].begin(), positions[i].end()))
117 DUNE_THROW(Dune::GridError, "Make sure to specify a monotone increasing \"Positions\" array");
118
119 if(positions[i].size() < 2)
120 DUNE_THROW(Dune::GridError, "Make sure to specify position arrays with at least two entries (min and max value).");
121 }
122
123 // check that all indices are specified i.e. > 0
124 for (int i = 0; i < dim; ++i)
125 if (indices[i] < 0)
126 DUNE_THROW(Dune::RangeError, "Please specify Positions Angular and Radial and Axial correctly and unambiguously!" << std::endl);
127
128 // the number of cells (has a default)
129 std::array<std::vector<int>, dim> cells;
130 for (int i = 0; i < dim; ++i)
131 {
132 cells[i].resize(positions[i].size()-1, 1.0);
133 cells[i] = getParamFromGroup<std::vector<int>>(modelParamGroup, "Grid.Cells" + std::to_string(i), cells[i]);
134 if (cells[i].size() + 1 != positions[i].size())
135 DUNE_THROW(Dune::RangeError, "Make sure to specify the correct length of \"Cells\" and \"Positions\" arrays");
136 }
137
138 // grading factor (has a default)
139 std::array<std::vector<Scalar>, dim> grading;
140 for (int i = 0; i < dim; ++i)
141 {
142 grading[i].resize(positions[i].size()-1, 1.0);
143 grading[i] = getParamFromGroup<std::vector<Scalar>>(modelParamGroup, "Grid.Grading" + std::to_string(i), grading[i]);
144 if (grading[i].size() + 1 != positions[i].size())
145 DUNE_THROW(Dune::RangeError, "Make sure to specify the correct length of \"Grading\" and \"Positions\" arrays");
146 }
147
148 // make the grid
149 std::array<std::vector<Scalar>, dim> globalPositions;
150 using std::pow;
151 using Dune::power;
152 for (int dimIdx = 0; dimIdx < dim; dimIdx++)
153 {
154 // Each grid direction is subdivided into (numCells + 1) points
155 std::size_t numGlobalPositions = 1;
156 for (int zoneIdx = 0; zoneIdx < cells[dimIdx].size(); ++zoneIdx)
157 numGlobalPositions += cells[dimIdx][zoneIdx];
158
159 globalPositions[dimIdx].resize(numGlobalPositions);
160 std::size_t posIdx = 0;
161 for (int zoneIdx = 0; zoneIdx < cells[dimIdx].size(); ++zoneIdx)
162 {
163 const Scalar lower = positions[dimIdx][zoneIdx];
164 const Scalar upper = positions[dimIdx][zoneIdx+1];
165 const int numCells = cells[dimIdx][zoneIdx];
166 Scalar gradingFactor = grading[dimIdx][zoneIdx];
167 const Scalar length = upper - lower;
168 Scalar height = 1.0;
169 bool increasingCellSize = false;
170
171 if (verbose)
172 {
173 std::cout << "dim " << dimIdx
174 << " lower " << lower
175 << " upper " << upper
176 << " numCells " << numCells
177 << " grading " << gradingFactor;
178 }
179
180 if (gradingFactor > 1.0)
181 increasingCellSize = true;
182
183 // take absolute values and reverse cell size increment to achieve
184 // reverse behavior for negative values
185 if (gradingFactor < 0.0)
186 {
187 using std::abs;
188 gradingFactor = abs(gradingFactor);
189
190 if (gradingFactor < 1.0)
191 increasingCellSize = true;
192 }
193
194 const bool useGrading = Dune::FloatCmp::eq(gradingFactor, 1.0) ? false : true;
195
196 // if the grading factor is exactly 1.0 do equal spacing
197 if (!useGrading)
198 {
199 height = 1.0 / numCells;
200 if (verbose)
201 std::cout << " -> h " << height * length << std::endl;
202 }
203 // if grading factor is not 1.0, do power law spacing
204 else
205 {
206 height = (1.0 - gradingFactor) / (1.0 - power(gradingFactor, numCells));
207
208 if (verbose)
209 {
210 std::cout << " -> grading_eff " << gradingFactor
211 << " h_min " << height * power(gradingFactor, 0) * length
212 << " h_max " << height * power(gradingFactor, numCells-1) * length
213 << std::endl;
214 }
215 }
216
217 // the positions inside a global segment
218 Dune::DynamicVector<Scalar> localPositions(numCells, 0.0);
219 for (int i = 1; i < numCells; i++)
220 {
221 Scalar hI = height;
222 if (useGrading)
223 {
224 if (increasingCellSize)
225 hI *= power(gradingFactor, i-1);
226
227 else
228 hI *= power(gradingFactor, numCells-i);
229 }
230 localPositions[i] = localPositions[i-1] + hI;
231 }
232
233 localPositions *= length;
234 localPositions += lower;
235
236 for (int i = 0; i < numCells; ++i)
237 globalPositions[dimIdx][posIdx++] = localPositions[i];
238 }
239
240 globalPositions[dimIdx][posIdx] = positions[dimIdx].back();
241 }
242
243 polarCoordinates[0] = globalPositions[indices[0]];
244 polarCoordinates[1] = globalPositions[indices[1]];
245 if (dim == 3)
246 polarCoordinates[2] = globalPositions[dim - indices[0] - indices[1]];
247
248 // convert angular coordinates into radians
249 std::transform(polarCoordinates[1].begin(), polarCoordinates[1].end(),
250 polarCoordinates[1].begin(),
251 [](Scalar s){ return s*M_PI/180; });
252 }
253
262 std::unique_ptr<Grid> createCakeGrid(std::array<std::vector<Scalar>, dim> &polarCoordinates,
263 Dune::FieldVector<int, dim> &indices,
264 const std::string& modelParamGroup,
265 bool verbose = false)
266 {
267 GridFactory gridFactory;
268
269 // create grid on rank 0, return empty grid on all other ranks
270 if (gridFactory.comm().rank() != 0)
271 return std::unique_ptr<Grid>(gridFactory.createGrid());
272
273 const auto& dR = polarCoordinates[0];
274 const auto& dA = polarCoordinates[1];
275
276 // For the special case of 360°, the last element is connected with the first one.
277 // Thus, one needs to distinguish here between all other cases, were the normal procedure
278 // is applied for all elements until the last one (= dA.size() - 1) and the 360° case,
279 // where it is only applied until the second to last element (dA.size() - 2).
280 const int dAMax = Dune::FloatCmp::eq(dA[dA.size()-1], 2*M_PI) ? dA.size() - 2 : dA.size() - 1;
281 constexpr auto type = Dune::GeometryTypes::cube(dim);
282 const bool hasHole = dR[0] >= 1.0e-8*dR.back();
283
284 // create nodes
285 if constexpr (dim == 3)
286 {
287 constexpr auto prismType = Dune::GeometryTypes::prism;
288 std::vector<Scalar> dZ = polarCoordinates[2];
289 for (int j = 0; j <= dAMax; ++j)
290 {
291 for (int l = 0; l <= dZ.size() - 1; ++l)
292 {
293 for (int i = hasHole ? 0 : 1; i <= dR.size()- 1; ++i)
294 {
295 // transform into Cartesian coordinates
296 using std::cos;
297 using std::sin;
298 Dune::FieldVector <double, dim> v(0.0);
299 v[indices[2]] = dZ[l];
300 v[indices[0]] = cos(dA[j])*dR[i];
301 v[indices[1]] = sin(dA[j])*dR[i];
302 if (verbose)
304 gridFactory.insertVertex(v);
305 }
306 }
307 }
308
309 // if the well is not considered, we have to add the points lying on the center-line
310 // this has to be done separately, otherwise these points would occur multiple times
311 if(!hasHole)
312 {
313 for (int l = 0; l <= dZ.size() - 1; ++l)
314 {
315 Dune::FieldVector <double, dim> v(0.0);
316 v[indices[2]] = dZ[l];
317 v[indices[0]] = 0.0;
318 v[indices[1]] = 0.0;
319
320 if (verbose)
322 gridFactory.insertVertex(v);
323 }
324 }
325
326 if (verbose)
327 std::cout << "Filled node vector" << std::endl;
328
329 // assign nodes
330 unsigned int z = 0;
331 unsigned int t = 0;
332 unsigned int rSize = hasHole ? dR.size() : dR.size()-1;
333 unsigned int zSize = dZ.size();
334 for (int j = 0; j < dA.size() - 1; ++j)
335 {
336 for (int l = 0; l < dZ.size() - 1; ++l)
337 {
338 if (j < dAMax)
339 {
340 for (int i = 0; i < dR.size() - 1; ++i)
341 {
342 if(!hasHole && i==0)
343 {
344 std::vector<unsigned int> vid({rSize*zSize*(dAMax+1) + l, z, z+rSize*zSize,
345 rSize*zSize*(dAMax+1) + l+1, z+rSize, z+rSize*zSize+rSize});
346
347 if (verbose)
348 printIndices(vid);
349
350 gridFactory.insertElement(prismType, vid);
351 }
352 else
353 {
354 std::vector<unsigned int> vid({z, z+1,
355 z+rSize*zSize, z+rSize*zSize+1,
356 z+rSize, z+rSize+1,
357 z+rSize*zSize+rSize, z+rSize*zSize+rSize+1});
358
359 if (verbose)
360 printIndices(vid);
361
362 gridFactory.insertElement(type, vid);
363
364 z++;
365 }
366 }
367 z++;
368 }
369 else
370 {
371 // assign nodes for 360°-cake
372 for (int i = 0; i < dR.size() - 1; ++i)
373 {
374 if(!hasHole && i==0)
375 {
376 std::vector<unsigned int> vid({rSize*zSize*(dAMax+1) + l, z, t,
377 rSize*zSize*(dAMax+1) + l+1, z+rSize, t+rSize});
378
379 if (verbose)
380 printIndices(vid);
381
382 gridFactory.insertElement(prismType, vid);
383 }
384 else
385 {
386 std::vector<unsigned int> vid({z, z+1,
387 t, t+1,
388 z+rSize, z+rSize+1,
389 t+rSize, t+rSize+1});
390
391 if (verbose)
392 printIndices(vid);
393
394 gridFactory.insertElement(type, vid);
395 t++;
396 z++;
397 }
398 }
399 t++;
400 z++;
401
402 if (verbose)
403 std::cout << "assign nodes 360° ends..." << std::endl;
404 }
405 }
406 z += rSize;
407 }
408 }
409
410 // for dim = 2
411 else
412 {
413 constexpr auto triangleType = Dune::GeometryTypes::simplex(dim);
414 for (int j = 0; j <= dAMax; ++j)
415 {
416 for (int i = hasHole ? 0 : 1; i <= dR.size()- 1; ++i)
417 {
418 // transform into Cartesian coordinates
419 Dune::FieldVector <double, dim> v(0.0);
420
421 v[indices[0]] = cos(dA[j])*dR[i];
422 v[indices[1]] = sin(dA[j])*dR[i];
423 if(verbose)
425 gridFactory.insertVertex(v);
426 }
427 }
428
429 // if the well is not considered, we have to add the point located at the center
430 // this has to be done separately, otherwise this point would occur multiple times
431 if(!hasHole)
432 {
433 Dune::FieldVector <double, dim> v(0.0);
434 v[indices[0]] = 0.0;
435 v[indices[1]] = 0.0;
436
437 if(verbose)
439 gridFactory.insertVertex(v);
440 }
441
442 std::cout << "Filled node vector" << std::endl;
443
444 // assign nodes
445 unsigned int z = 0;
446 unsigned int t = 0;
447 unsigned int rSize = hasHole ? dR.size() : dR.size()-1;
448 for (int j = 0; j < dA.size() - 1; ++j)
449 {
450 if (j < dAMax)
451 {
452 for (int i = 0; i < dR.size() - 1; ++i)
453 {
454 if(!hasHole && i==0)
455 {
456 std::vector<unsigned int> vid({rSize*(dAMax+1), z, z+rSize});
457
458 if (verbose)
459 printIndices(vid);
460
461 gridFactory.insertElement(triangleType, vid);
462 }
463 else
464 {
465 std::vector<unsigned int> vid({z, z+1, z+rSize, z+rSize+1});
466
467 if (verbose)
468 printIndices(vid);
469
470 gridFactory.insertElement(type, vid);
471 z++;
472 }
473 }
474 z++;
475 }
476 else
477 {
478 // assign nodes for 360°-cake
479 for (int i = 0; i < dR.size() - 1; ++i)
480 {
481 if(!hasHole && i==0)
482 {
483 std::vector<unsigned int> vid({rSize*(dAMax+1), z, t});
484
485 if (verbose)
486 printIndices(vid);
487
488 gridFactory.insertElement(triangleType, vid);
489 }
490 else
491 {
492 std::vector<unsigned int> vid({z, z+1, t, t+1});
493
494 if (verbose)
495 printIndices(vid);
496
497 gridFactory.insertElement(type, vid);
498 t++;
499 z++;
500 }
501 }
502 t++;
503 z++;
504
505 if (verbose)
506 std::cout << "assign nodes 360 ends..." << std::endl;
507 }
508 }
509 }
510
511 // create the grid
512 return std::unique_ptr<Grid>(gridFactory.createGrid());
513 }
514
518 Grid& grid()
519 {
520 return *gridPtr();
521 }
522
528 {
529 gridPtr()->loadBalance();
530 }
531
532protected:
533 static void printCoordinate(const Dune::FieldVector <double, dim>& v)
534 { std::cout << "Coordinates of: " << v << std::endl; }
535
536 static void printIndices(const std::vector<unsigned int>& vid)
537 {
538 std::cout << "element vertex indices: ";
539 for (int k = 0; k < vid.size(); ++k)
540 std::cout << vid[k] << " ";
541 std::cout << std::endl;
542 }
543
547 GridPointer& gridPtr()
548 {
549 return cakeGrid_;
550 }
551
552private:
553 GridPointer cakeGrid_;
554};
555
556} // end namespace Dumux
557
558#endif
Provides a grid manager with a method for creating creating vectors with polar Coordinates and one fo...
Definition: cakegridmanager.hh:36
void init(const std::string &modelParamGroup="")
Make the grid.
Definition: cakegridmanager.hh:49
void loadBalance()
Distributes the grid on all processes of a parallel computation.
Definition: cakegridmanager.hh:527
GridPointer & gridPtr()
Returns a reference to the shared pointer to the grid.
Definition: cakegridmanager.hh:547
static void printCoordinate(const Dune::FieldVector< double, dim > &v)
Definition: cakegridmanager.hh:533
static void printIndices(const std::vector< unsigned int > &vid)
Definition: cakegridmanager.hh:536
Grid & grid()
Returns a reference to the grid.
Definition: cakegridmanager.hh:518
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:84
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:262
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
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.