3.6-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/common/math.hh>
35#include <dune/grid/common/gridfactory.hh>
37
38namespace Dumux {
39
46template <class Grid>
48{
49 using Scalar = typename Grid::ctype;
50
51 using GridFactory = Dune::GridFactory<Grid>;
52 using GridPointer = std::shared_ptr<Grid>;
53
54 enum { dim = Grid::dimension,
55 dimWorld = Grid::dimensionworld };
56
57public:
61 void init(const std::string& modelParamGroup = "")
62 {
63 static_assert(dim == 2 || dim == 3, "The CakeGridManager is only implemented for 2D and 3D.");
64
65 const bool verbose = getParamFromGroup<bool>(modelParamGroup, "Grid.Verbosity", false);
66
67 std::array<std::vector<Scalar>, dim> polarCoordinates;
68 // Indices specifying in which direction the piece of cake is oriented
69 Dune::FieldVector<int, dim> indices(-1);
70 createVectors(polarCoordinates, indices, modelParamGroup, verbose);
71
72 gridPtr() = createCakeGrid(polarCoordinates, indices, modelParamGroup, verbose);
74 }
75
96 static void createVectors(std::array<std::vector<Scalar>, dim> &polarCoordinates,
97 Dune::FieldVector<int, dim> &indices,
98 const std::string& modelParamGroup,
99 bool verbose = false)
100 {
101 // The positions
102 std::array<std::vector<Scalar>, dim> positions;
103
104 for (int i = 0; i < dim; ++i)
105 {
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);
111
112 if (hasRadial)
113 {
114 positions[i] = getParamFromGroup<std::vector<Scalar>>(modelParamGroup, "Grid.Radial" + std::to_string(i));
115 indices[0] = i; // Index specifying radial direction
116 }
117 else if (hasAngular)
118 {
119 positions[i] = getParamFromGroup<std::vector<Scalar>>(modelParamGroup, "Grid.Angular" + std::to_string(i));
120 indices[1] = i; // Index specifying angular direction
121 }
122 else // hasAxial
123 {
124 positions[i] = getParamFromGroup<std::vector<Scalar>>(modelParamGroup, "Grid.Axial" + std::to_string(i));
125 indices[2] = i; // Index specifying axial direction
126 }
127
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");
130
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).");
133 }
134
135 // check that all indices are specified i.e. > 0
136 for (int i = 0; i < dim; ++i)
137 if (indices[i] < 0)
138 DUNE_THROW(Dune::RangeError, "Please specify Positions Angular and Radial and Axial correctly and unambiguously!" << std::endl);
139
140 // the number of cells (has a default)
141 std::array<std::vector<int>, dim> cells;
142 for (int i = 0; i < dim; ++i)
143 {
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");
148 }
149
150 // grading factor (has a default)
151 std::array<std::vector<Scalar>, dim> grading;
152 for (int i = 0; i < dim; ++i)
153 {
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");
158 }
159
160 // make the grid
161 std::array<std::vector<Scalar>, dim> globalPositions;
162 using std::pow;
163 using Dune::power;
164 for (int dimIdx = 0; dimIdx < dim; dimIdx++)
165 {
166 // Each grid direction is subdivided into (numCells + 1) points
167 std::size_t numGlobalPositions = 1;
168 for (int zoneIdx = 0; zoneIdx < cells[dimIdx].size(); ++zoneIdx)
169 numGlobalPositions += cells[dimIdx][zoneIdx];
170
171 globalPositions[dimIdx].resize(numGlobalPositions);
172 std::size_t posIdx = 0;
173 for (int zoneIdx = 0; zoneIdx < cells[dimIdx].size(); ++zoneIdx)
174 {
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;
180 Scalar height = 1.0;
181 bool increasingCellSize = false;
182
183 if (verbose)
184 {
185 std::cout << "dim " << dimIdx
186 << " lower " << lower
187 << " upper " << upper
188 << " numCells " << numCells
189 << " grading " << gradingFactor;
190 }
191
192 if (gradingFactor > 1.0)
193 increasingCellSize = true;
194
195 // take absolute values and reverse cell size increment to achieve
196 // reverse behavior for negative values
197 if (gradingFactor < 0.0)
198 {
199 using std::abs;
200 gradingFactor = abs(gradingFactor);
201
202 if (gradingFactor < 1.0)
203 increasingCellSize = true;
204 }
205
206 const bool useGrading = Dune::FloatCmp::eq(gradingFactor, 1.0) ? false : true;
207
208 // if the grading factor is exactly 1.0 do equal spacing
209 if (!useGrading)
210 {
211 height = 1.0 / numCells;
212 if (verbose)
213 std::cout << " -> h " << height * length << std::endl;
214 }
215 // if grading factor is not 1.0, do power law spacing
216 else
217 {
218 height = (1.0 - gradingFactor) / (1.0 - power(gradingFactor, numCells));
219
220 if (verbose)
221 {
222 std::cout << " -> grading_eff " << gradingFactor
223 << " h_min " << height * power(gradingFactor, 0) * length
224 << " h_max " << height * power(gradingFactor, numCells-1) * length
225 << std::endl;
226 }
227 }
228
229 // the positions inside a global segment
230 Dune::DynamicVector<Scalar> localPositions(numCells, 0.0);
231 for (int i = 1; i < numCells; i++)
232 {
233 Scalar hI = height;
234 if (useGrading)
235 {
236 if (increasingCellSize)
237 hI *= power(gradingFactor, i-1);
238
239 else
240 hI *= power(gradingFactor, numCells-i);
241 }
242 localPositions[i] = localPositions[i-1] + hI;
243 }
244
245 localPositions *= length;
246 localPositions += lower;
247
248 for (int i = 0; i < numCells; ++i)
249 globalPositions[dimIdx][posIdx++] = localPositions[i];
250 }
251
252 globalPositions[dimIdx][posIdx] = positions[dimIdx].back();
253 }
254
255 polarCoordinates[0] = globalPositions[indices[0]];
256 polarCoordinates[1] = globalPositions[indices[1]];
257 if (dim == 3)
258 polarCoordinates[2] = globalPositions[dim - indices[0] - indices[1]];
259
260 // convert angular coordinates into radians
261 std::transform(polarCoordinates[1].begin(), polarCoordinates[1].end(),
262 polarCoordinates[1].begin(),
263 [](Scalar s){ return s*M_PI/180; });
264 }
265
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)
278 {
279 GridFactory gridFactory;
280
281 // create grid on rank 0, return empty grid on all other ranks
282 if (gridFactory.comm().rank() != 0)
283 return std::unique_ptr<Grid>(gridFactory.createGrid());
284
285 const auto& dR = polarCoordinates[0];
286 const auto& 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 const int dAMax = Dune::FloatCmp::eq(dA[dA.size()-1], 2*M_PI) ? dA.size() - 2 : dA.size() - 1;
293 constexpr auto type = Dune::GeometryTypes::cube(dim);
294 const bool hasHole = dR[0] >= 1.0e-8*dR.back();
295
296 // create nodes
297 if constexpr (dim == 3)
298 {
299 constexpr auto prismType = Dune::GeometryTypes::prism;
300 std::vector<Scalar> dZ = polarCoordinates[2];
301 for (int j = 0; j <= dAMax; ++j)
302 {
303 for (int l = 0; l <= dZ.size() - 1; ++l)
304 {
305 for (int i = hasHole ? 0 : 1; i <= dR.size()- 1; ++i)
306 {
307 // transform into Cartesian coordinates
308 using std::cos;
309 using std::sin;
310 Dune::FieldVector <double, dim> v(0.0);
311 v[indices[2]] = dZ[l];
312 v[indices[0]] = cos(dA[j])*dR[i];
313 v[indices[1]] = sin(dA[j])*dR[i];
314 if (verbose)
316 gridFactory.insertVertex(v);
317 }
318 }
319 }
320
321 // if the well is not considered, we have to add the points lying on the center-line
322 // this has to be done separately, otherwise these points would occur multiple times
323 if(!hasHole)
324 {
325 for (int l = 0; l <= dZ.size() - 1; ++l)
326 {
327 Dune::FieldVector <double, dim> v(0.0);
328 v[indices[2]] = dZ[l];
329 v[indices[0]] = 0.0;
330 v[indices[1]] = 0.0;
331
332 if (verbose)
334 gridFactory.insertVertex(v);
335 }
336 }
337
338 if (verbose)
339 std::cout << "Filled node vector" << std::endl;
340
341 // assign nodes
342 unsigned int z = 0;
343 unsigned int t = 0;
344 unsigned int rSize = hasHole ? dR.size() : dR.size()-1;
345 unsigned int zSize = dZ.size();
346 for (int j = 0; j < dA.size() - 1; ++j)
347 {
348 for (int l = 0; l < dZ.size() - 1; ++l)
349 {
350 if (j < dAMax)
351 {
352 for (int i = 0; i < dR.size() - 1; ++i)
353 {
354 if(!hasHole && i==0)
355 {
356 std::vector<unsigned int> vid({rSize*zSize*(dAMax+1) + l, z, z+rSize*zSize,
357 rSize*zSize*(dAMax+1) + l+1, z+rSize, z+rSize*zSize+rSize});
358
359 if (verbose)
360 printIndices(vid);
361
362 gridFactory.insertElement(prismType, vid);
363 }
364 else
365 {
366 std::vector<unsigned int> vid({z, z+1,
367 z+rSize*zSize, z+rSize*zSize+1,
368 z+rSize, z+rSize+1,
369 z+rSize*zSize+rSize, z+rSize*zSize+rSize+1});
370
371 if (verbose)
372 printIndices(vid);
373
374 gridFactory.insertElement(type, vid);
375
376 z++;
377 }
378 }
379 z++;
380 }
381 else
382 {
383 // assign nodes for 360°-cake
384 for (int i = 0; i < dR.size() - 1; ++i)
385 {
386 if(!hasHole && i==0)
387 {
388 std::vector<unsigned int> vid({rSize*zSize*(dAMax+1) + l, z, t,
389 rSize*zSize*(dAMax+1) + l+1, z+rSize, t+rSize});
390
391 if (verbose)
392 printIndices(vid);
393
394 gridFactory.insertElement(prismType, vid);
395 }
396 else
397 {
398 std::vector<unsigned int> vid({z, z+1,
399 t, t+1,
400 z+rSize, z+rSize+1,
401 t+rSize, t+rSize+1});
402
403 if (verbose)
404 printIndices(vid);
405
406 gridFactory.insertElement(type, vid);
407 t++;
408 z++;
409 }
410 }
411 t++;
412 z++;
413
414 if (verbose)
415 std::cout << "assign nodes 360° ends..." << std::endl;
416 }
417 }
418 z += rSize;
419 }
420 }
421
422 // for dim = 2
423 else
424 {
425 constexpr auto triangleType = Dune::GeometryTypes::simplex(dim);
426 for (int j = 0; j <= dAMax; ++j)
427 {
428 for (int i = hasHole ? 0 : 1; i <= dR.size()- 1; ++i)
429 {
430 // transform into Cartesian coordinates
431 Dune::FieldVector <double, dim> v(0.0);
432
433 v[indices[0]] = cos(dA[j])*dR[i];
434 v[indices[1]] = sin(dA[j])*dR[i];
435 if(verbose)
437 gridFactory.insertVertex(v);
438 }
439 }
440
441 // if the well is not considered, we have to add the point located at the center
442 // this has to be done separately, otherwise this point would occur multiple times
443 if(!hasHole)
444 {
445 Dune::FieldVector <double, dim> v(0.0);
446 v[indices[0]] = 0.0;
447 v[indices[1]] = 0.0;
448
449 if(verbose)
451 gridFactory.insertVertex(v);
452 }
453
454 std::cout << "Filled node vector" << std::endl;
455
456 // assign nodes
457 unsigned int z = 0;
458 unsigned int t = 0;
459 unsigned int rSize = hasHole ? dR.size() : dR.size()-1;
460 for (int j = 0; j < dA.size() - 1; ++j)
461 {
462 if (j < dAMax)
463 {
464 for (int i = 0; i < dR.size() - 1; ++i)
465 {
466 if(!hasHole && i==0)
467 {
468 std::vector<unsigned int> vid({rSize*(dAMax+1), z, z+rSize});
469
470 if (verbose)
471 printIndices(vid);
472
473 gridFactory.insertElement(triangleType, vid);
474 }
475 else
476 {
477 std::vector<unsigned int> vid({z, z+1, z+rSize, z+rSize+1});
478
479 if (verbose)
480 printIndices(vid);
481
482 gridFactory.insertElement(type, vid);
483 z++;
484 }
485 }
486 z++;
487 }
488 else
489 {
490 // assign nodes for 360°-cake
491 for (int i = 0; i < dR.size() - 1; ++i)
492 {
493 if(!hasHole && i==0)
494 {
495 std::vector<unsigned int> vid({rSize*(dAMax+1), z, t});
496
497 if (verbose)
498 printIndices(vid);
499
500 gridFactory.insertElement(triangleType, vid);
501 }
502 else
503 {
504 std::vector<unsigned int> vid({z, z+1, t, t+1});
505
506 if (verbose)
507 printIndices(vid);
508
509 gridFactory.insertElement(type, vid);
510 t++;
511 z++;
512 }
513 }
514 t++;
515 z++;
516
517 if (verbose)
518 std::cout << "assign nodes 360 ends..." << std::endl;
519 }
520 }
521 }
522
523 // create the grid
524 return std::unique_ptr<Grid>(gridFactory.createGrid());
525 }
526
530 Grid& grid()
531 {
532 return *gridPtr();
533 }
534
540 {
541 gridPtr()->loadBalance();
542 }
543
544protected:
545 static void printCoordinate(const Dune::FieldVector <double, dim>& v)
546 { std::cout << "Coordinates of: " << v << std::endl; }
547
548 static void printIndices(const std::vector<unsigned int>& vid)
549 {
550 std::cout << "element vertex indices: ";
551 for (int k = 0; k < vid.size(); ++k)
552 std::cout << vid[k] << " ";
553 std::cout << std::endl;
554 }
555
559 GridPointer& gridPtr()
560 {
561 return cakeGrid_;
562 }
563
564private:
565 GridPointer cakeGrid_;
566};
567
568} // end namespace Dumux
569
570#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:177
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
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:539
GridPointer & gridPtr()
Returns a reference to the shared pointer to the grid.
Definition: cakegridmanager.hh:559
static void printCoordinate(const Dune::FieldVector< double, dim > &v)
Definition: cakegridmanager.hh:545
static void printIndices(const std::vector< unsigned int > &vid)
Definition: cakegridmanager.hh:548
Grid & grid()
Returns a reference to the grid.
Definition: cakegridmanager.hh:530
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