3.2-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_CAKE_GRID_MANAGER_HH
25#define DUMUX_CAKE_GRID_MANAGER_HH
26
27#include <vector>
28#include <iostream>
29#include <cmath>
30#include <algorithm>
31
32#include <dune/common/dynvector.hh>
33#include <dune/common/float_cmp.hh>
34#include <dune/grid/common/gridfactory.hh>
36
37namespace Dumux {
38
45template <class Grid>
47{
48 using Scalar = typename Grid::ctype;
49
50 using GridFactory = Dune::GridFactory<Grid>;
51 using GridPointer = std::shared_ptr<Grid>;
52
53 enum { dim = Grid::dimension,
54 dimWorld = Grid::dimensionworld };
55
56public:
60 void init(const std::string& modelParamGroup = "")
61 {
62 static_assert(dim == 2 || dim == 3, "The CakeGridManager is only implemented for 2D and 3D.");
63
64 const bool verbose = getParamFromGroup<bool>(modelParamGroup, "Grid.Verbosity", false);
65
66 std::array<std::vector<Scalar>, dim> polarCoordinates;
67 // Indices specifying in which direction the piece of cake is oriented
68 Dune::FieldVector<int, dim> indices(-1);
69 createVectors(polarCoordinates, indices, modelParamGroup, verbose);
70
71 gridPtr() = createCakeGrid(polarCoordinates, indices, modelParamGroup, verbose);
73 }
74
95 static void createVectors(std::array<std::vector<Scalar>, dim> &polarCoordinates,
96 Dune::FieldVector<int, dim> &indices,
97 const std::string& modelParamGroup,
98 bool verbose = false)
99 {
100 // The positions
101 std::array<std::vector<Scalar>, dim> positions;
102
103 for (int i = 0; i < dim; ++i)
104 {
105 const bool hasRadial = hasParamInGroup(modelParamGroup, "Grid.Radial" + std::to_string(i));
106 const bool hasAngular = hasParamInGroup(modelParamGroup, "Grid.Angular" + std::to_string(i));
107 const bool hasAxial = (dim == 3) && hasParamInGroup(modelParamGroup, "Grid.Axial" + std::to_string(i));
108 if (static_cast<int>(hasRadial) + static_cast<int>(hasAngular) + static_cast<int>(hasAxial) != 1)
109 DUNE_THROW(Dune::RangeError, "Multiple or no position vectors (radial, angular, axial) specified in coord direction: " << i << std::endl);
110
111 if (hasRadial)
112 {
113 positions[i] = getParamFromGroup<std::vector<Scalar>>(modelParamGroup, "Grid.Radial" + std::to_string(i));
114 indices[0] = i; // Index specifying radial direction
115 }
116 else if (hasAngular)
117 {
118 positions[i] = getParamFromGroup<std::vector<Scalar>>(modelParamGroup, "Grid.Angular" + std::to_string(i));
119 indices[1] = i; // Index specifying angular direction
120 }
121 else // hasAxial
122 {
123 positions[i] = getParamFromGroup<std::vector<Scalar>>(modelParamGroup, "Grid.Axial" + std::to_string(i));
124 indices[2] = i; // Index specifying axial direction
125 }
126
127 if (!std::is_sorted(positions[i].begin(), positions[i].end()))
128 DUNE_THROW(Dune::GridError, "Make sure to specify a monotone increasing \"Positions\" array");
129
130 if(positions[i].size() < 2)
131 DUNE_THROW(Dune::GridError, "Make sure to specify position arrays with at least two entries (min and max value).");
132 }
133
134 // check that all indices are specified i.e. > 0
135 for (int i = 0; i < dim; ++i)
136 if (indices[i] < 0)
137 DUNE_THROW(Dune::RangeError, "Please specify Positions Angular and Radial and Axial correctly and unambiguously!" << std::endl);
138
139 if(hasParam("Grid.WellRadius"))
140 {
141 std::cerr << "Deprecation warning: parameter Grid.WellRadius is deprecated. "
142 << "Specify the WellRadius as the first radial coordinate." << std::endl;
143
144 positions[indices[0]][0] = getParamFromGroup<Scalar>(modelParamGroup, "Grid.WellRadius");
145 }
146
147 // the number of cells (has a default)
148 std::array<std::vector<int>, dim> cells;
149 for (int i = 0; i < dim; ++i)
150 {
151 cells[i].resize(positions[i].size()-1, 1.0);
152 cells[i] = getParamFromGroup<std::vector<int>>(modelParamGroup, "Grid.Cells" + std::to_string(i), cells[i]);
153 if (cells[i].size() + 1 != positions[i].size())
154 DUNE_THROW(Dune::RangeError, "Make sure to specify the correct length of \"Cells\" and \"Positions\" arrays");
155 }
156
157 // grading factor (has a default)
158 std::array<std::vector<Scalar>, dim> grading;
159 for (int i = 0; i < dim; ++i)
160 {
161 grading[i].resize(positions[i].size()-1, 1.0);
162 grading[i] = getParamFromGroup<std::vector<Scalar>>(modelParamGroup, "Grid.Grading" + std::to_string(i), grading[i]);
163 if (grading[i].size() + 1 != positions[i].size())
164 DUNE_THROW(Dune::RangeError, "Make sure to specify the correct length of \"Grading\" and \"Positions\" arrays");
165 }
166
167 // make the grid
168 std::array<std::vector<Scalar>, dim> globalPositions;
169 using std::pow;
170 for (int dimIdx = 0; dimIdx < dim; dimIdx++)
171 {
172 // Each grid direction is subdivided into (numCells + 1) points
173 std::size_t numGlobalPositions = 1;
174 for (int zoneIdx = 0; zoneIdx < cells[dimIdx].size(); ++zoneIdx)
175 numGlobalPositions += cells[dimIdx][zoneIdx];
176
177 globalPositions[dimIdx].resize(numGlobalPositions);
178 std::size_t posIdx = 0;
179 for (int zoneIdx = 0; zoneIdx < cells[dimIdx].size(); ++zoneIdx)
180 {
181 const Scalar lower = positions[dimIdx][zoneIdx];
182 const Scalar upper = positions[dimIdx][zoneIdx+1];
183 const int numCells = cells[dimIdx][zoneIdx];
184 Scalar gradingFactor = grading[dimIdx][zoneIdx];
185 const Scalar length = upper - lower;
186 Scalar height = 1.0;
187 bool increasingCellSize = false;
188
189 if (verbose)
190 {
191 std::cout << "dim " << dimIdx
192 << " lower " << lower
193 << " upper " << upper
194 << " numCells " << numCells
195 << " grading " << gradingFactor;
196 }
197
198 if (gradingFactor > 1.0)
199 increasingCellSize = true;
200
201 // take absolute values and reverse cell size increment to achieve
202 // reverse behavior for negative values
203 if (gradingFactor < 0.0)
204 {
205 using std::abs;
206 gradingFactor = abs(gradingFactor);
207
208 if (gradingFactor < 1.0)
209 increasingCellSize = true;
210 }
211
212 const bool useGrading = Dune::FloatCmp::eq(gradingFactor, 1.0) ? false : true;
213
214 // if the grading factor is exactly 1.0 do equal spacing
215 if (!useGrading)
216 {
217 height = 1.0 / numCells;
218 if (verbose)
219 std::cout << " -> h " << height * length << std::endl;
220 }
221 // if grading factor is not 1.0, do power law spacing
222 else
223 {
224 height = (1.0 - gradingFactor) / (1.0 - pow(gradingFactor, numCells));
225
226 if (verbose)
227 {
228 std::cout << " -> grading_eff " << gradingFactor
229 << " h_min " << height * pow(gradingFactor, 0) * length
230 << " h_max " << height * pow(gradingFactor, numCells-1) * length
231 << std::endl;
232 }
233 }
234
235 // the positions inside a global segment
236 Dune::DynamicVector<Scalar> localPositions(numCells, 0.0);
237 for (int i = 1; i < numCells; i++)
238 {
239 Scalar hI = height;
240 if (useGrading)
241 {
242 if (increasingCellSize)
243 hI *= pow(gradingFactor, i-1);
244
245 else
246 hI *= pow(gradingFactor, numCells-i);
247 }
248 localPositions[i] = localPositions[i-1] + hI;
249 }
250
251 localPositions *= length;
252 localPositions += lower;
253
254 for (int i = 0; i < numCells; ++i)
255 globalPositions[dimIdx][posIdx++] = localPositions[i];
256 }
257
258 globalPositions[dimIdx][posIdx] = positions[dimIdx].back();
259 }
260
261 polarCoordinates[0] = globalPositions[indices[0]];
262 polarCoordinates[1] = globalPositions[indices[1]];
263 if (dim == 3)
264 polarCoordinates[2] = globalPositions[dim - indices[0] - indices[1]];
265
266 // convert angular coordinates into radians
267 std::transform(polarCoordinates[1].begin(), polarCoordinates[1].end(),
268 polarCoordinates[1].begin(),
269 [](Scalar s){ return s*M_PI/180; });
270 }
271
280 std::unique_ptr<Grid> createCakeGrid(std::array<std::vector<Scalar>, dim> &polarCoordinates,
281 Dune::FieldVector<int, dim> &indices,
282 const std::string& modelParamGroup,
283 bool verbose = false)
284 {
285 std::vector<Scalar> dR = polarCoordinates[0];
286 std::vector<Scalar> dA = polarCoordinates[1];
287
288 // For the special case of 360°, the last element is connected with the first one.
289 // Thus, one needs to distinguish here between all other cases, were the normal procedure
290 // is applied for all elements until the last one (= dA.size() - 1) and the 360° case,
291 // where it is only applied until the second to last element (dA.size() - 2).
292 int maxdA = dA.size() - 1;
293 if (Dune::FloatCmp::eq(dA[dA.size()-1], 2*M_PI))
294 {
295 maxdA = dA.size() - 2;
296 }
297
298 GridFactory gridFactory;
299 constexpr auto type = Dune::GeometryTypes::cube(dim);
300
301 bool hasHole = true;
302 if(dR[0] < 1.0e-8*dR.back())
303 hasHole = false;
304
305 // create nodes
306 if (dim == 3)
307 {
308 constexpr auto prismType = Dune::GeometryTypes::prism;
309 std::vector<Scalar> dZ = polarCoordinates[2];
310 for (int j = 0; j <= maxdA; ++j)
311 {
312 for (int l = 0; l <= dZ.size() - 1; ++l)
313 {
314 for (int i = hasHole ? 0 : 1; i <= dR.size()- 1; ++i)
315 {
316 // transform into Cartesian coordinates
317 using std::cos;
318 using std::sin;
319 Dune::FieldVector <double, dim> v(0.0);
320 v[indices[2]] = dZ[l];
321 v[indices[0]] = cos(dA[j])*dR[i];
322 v[indices[1]] = sin(dA[j])*dR[i];
323 if (verbose)
325 gridFactory.insertVertex(v);
326 }
327 }
328 }
329
330 // if the well is not considered, we have to add the points lying on the center-line
331 // this has to be done separately, otherwise these points would occur multiple times
332 if(!hasHole)
333 {
334 for (int l = 0; l <= dZ.size() - 1; ++l)
335 {
336 Dune::FieldVector <double, dim> v(0.0);
337 v[indices[2]] = dZ[l];
338 v[indices[0]] = 0.0;
339 v[indices[1]] = 0.0;
340
341 if (verbose)
343 gridFactory.insertVertex(v);
344 }
345 }
346
347 if (verbose)
348 std::cout << "Filled node vector" << std::endl;
349
350 // assign nodes
351 unsigned int z = 0;
352 unsigned int t = 0;
353 unsigned int rSize = hasHole ? dR.size() : dR.size()-1;
354 unsigned int zSize = dZ.size();
355 for (int j = 0; j < dA.size() - 1; ++j)
356 {
357 for (int l = 0; l < dZ.size() - 1; ++l)
358 {
359 if (j < maxdA)
360 {
361 for (int i = 0; i < dR.size() - 1; ++i)
362 {
363 if(!hasHole && i==0)
364 {
365 std::vector<unsigned int> vid({rSize*zSize*(maxdA+1) + l, z, z+rSize*zSize,
366 rSize*zSize*(maxdA+1) + l+1, z+rSize, z+rSize*zSize+rSize});
367
368 if (verbose)
369 printIndices(vid);
370
371 gridFactory.insertElement(prismType, vid);
372 }
373 else
374 {
375 std::vector<unsigned int> vid({z, z+1,
376 z+rSize*zSize, z+rSize*zSize+1,
377 z+rSize, z+rSize+1,
378 z+rSize*zSize+rSize, z+rSize*zSize+rSize+1});
379
380 if (verbose)
381 printIndices(vid);
382
383 gridFactory.insertElement(type, vid);
384
385 z++;
386 }
387 }
388 z++;
389 }
390 else
391 {
392 // assign nodes for 360°-cake
393 for (int i = 0; i < dR.size() - 1; ++i)
394 {
395 if(!hasHole && i==0)
396 {
397 std::vector<unsigned int> vid({rSize*zSize*(maxdA+1) + l, z, t,
398 rSize*zSize*(maxdA+1) + l+1, z+rSize, t+rSize});
399
400 if (verbose)
401 printIndices(vid);
402
403 gridFactory.insertElement(prismType, vid);
404 }
405 else
406 {
407 std::vector<unsigned int> vid({z, z+1,
408 t, t+1,
409 z+rSize, z+rSize+1,
410 t+rSize, t+rSize+1});
411
412 if (verbose)
413 printIndices(vid);
414
415 gridFactory.insertElement(type, vid);
416 t++;
417 z++;
418 }
419 }
420 t++;
421 z++;
422
423 if (verbose)
424 std::cout << "assign nodes 360° ends..." << std::endl;
425 }
426 }
427 z += rSize;
428 }
429 }
430 // for dim = 2
431 else
432 {
433 constexpr auto triangleType = Dune::GeometryTypes::simplex(dim);
434 for (int j = 0; j <= maxdA; ++j)
435 {
436 for (int i = hasHole ? 0 : 1; i <= dR.size()- 1; ++i)
437 {
438 // transform into Cartesian coordinates
439 Dune::FieldVector <double, dim> v(0.0);
440
441 v[indices[0]] = cos(dA[j])*dR[i];
442 v[indices[1]] = sin(dA[j])*dR[i];
443 if(verbose)
445 gridFactory.insertVertex(v);
446 }
447 }
448
449 // if the well is not considered, we have to add the point located at the center
450 // this has to be done separately, otherwise this point would occur multiple times
451 if(!hasHole)
452 {
453 Dune::FieldVector <double, dim> v(0.0);
454 v[indices[0]] = 0.0;
455 v[indices[1]] = 0.0;
456
457 if(verbose)
459 gridFactory.insertVertex(v);
460 }
461
462 std::cout << "Filled node vector" << std::endl;
463
464 // assign nodes
465 unsigned int z = 0;
466 unsigned int t = 0;
467 unsigned int rSize = hasHole ? dR.size() : dR.size()-1;
468 for (int j = 0; j < dA.size() - 1; ++j)
469 {
470 if (j < maxdA)
471 {
472 for (int i = 0; i < dR.size() - 1; ++i)
473 {
474 if(!hasHole && i==0)
475 {
476 std::vector<unsigned int> vid({rSize*(maxdA+1), z, z+rSize});
477
478 if (verbose)
479 printIndices(vid);
480
481 gridFactory.insertElement(triangleType, vid);
482 }
483 else
484 {
485 std::vector<unsigned int> vid({z, z+1, z+rSize, z+rSize+1});
486
487 if (verbose)
488 printIndices(vid);
489
490 gridFactory.insertElement(type, vid);
491 z++;
492 }
493 }
494 z++;
495 }
496 else
497 {
498 // assign nodes for 360°-cake
499 for (int i = 0; i < dR.size() - 1; ++i)
500 {
501 if(!hasHole && i==0)
502 {
503 std::vector<unsigned int> vid({rSize*(maxdA+1), z, t});
504
505 if (verbose)
506 printIndices(vid);
507
508 gridFactory.insertElement(triangleType, vid);
509 }
510 else
511 {
512 std::vector<unsigned int> vid({z, z+1, t, t+1});
513
514 if (verbose)
515 printIndices(vid);
516
517 gridFactory.insertElement(type, vid);
518 t++;
519 z++;
520 }
521 }
522 t++;
523 z++;
524
525 if (verbose)
526 std::cout << "assign nodes 360 ends..." << std::endl;
527 }
528 }
529 }
530 // return the grid pointer
531 return std::unique_ptr<Grid>(gridFactory.createGrid());
532 }
533
537 Grid& grid()
538 {
539 return *gridPtr();
540 }
541
547 {
548 gridPtr()->loadBalance();
549 }
550
551protected:
552 static void printCoordinate(const Dune::FieldVector <double, dim>& v)
553 {
554 std::cout << "Coordinates of : ";
555 for (int k = 0; k < v.size(); ++k)
556 std::cout << v[k] << " ";
557 std::cout << std::endl;
558 }
559
560 static void printIndices(const std::vector<unsigned int>& vid)
561 {
562 std::cout << "element vertex indices: ";
563 for (int k = 0; k < vid.size(); ++k)
564 std::cout << vid[k] << " ";
565 std::cout << std::endl;
566 }
567
571 GridPointer& gridPtr()
572 {
573 return cakeGrid_;
574 }
575
576private:
577 GridPointer cakeGrid_;
578};
579
580} // end namespace Dumux
581
582#endif
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
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:391
bool hasParam(const std::string &param)
Check whether a key exists in the parameter tree.
Definition: parameters.hh:383
Definition: adapt.hh:29
Provides a grid manager with a method for creating creating vectors with polar Coordinates and one fo...
Definition: cakegridmanager.hh:47
void init(const std::string &modelParamGroup="")
Make the grid.
Definition: cakegridmanager.hh:60
void loadBalance()
Distributes the grid on all processes of a parallel computation.
Definition: cakegridmanager.hh:546
GridPointer & gridPtr()
Returns a reference to the shared pointer to the grid.
Definition: cakegridmanager.hh:571
static void printCoordinate(const Dune::FieldVector< double, dim > &v)
Definition: cakegridmanager.hh:552
static void printIndices(const std::vector< unsigned int > &vid)
Definition: cakegridmanager.hh:560
Grid & grid()
Returns a reference to the grid.
Definition: cakegridmanager.hh:537
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:95
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:280