cakegridmanager.hh 20 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
// vi: set et ts=4 sw=4 sts=4:
/*****************************************************************************
 *   See the file COPYING for full copying permissions.                      *
 *                                                                           *
 *   This program is free software: you can redistribute it and/or modify    *
 *   it under the terms of the GNU General Public License as published by    *
 *   the Free Software Foundation, either version 3 of the License, or       *
 *   (at your option) any later version.                                     *
 *                                                                           *
 *   This program is distributed in the hope that it will be useful,         *
 *   but WITHOUT ANY WARRANTY; without even the implied warranty of          *
 *   MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the            *
 *   GNU General Public License for more details.                            *
 *                                                                           *
 *   You should have received a copy of the GNU General Public License       *
 *   along with this program.  If not, see <http://www.gnu.org/licenses/>.   *
 *****************************************************************************/
/*!
 * \file
 * \ingroup InputOutput
 * \brief Provides a grid manager for a piece of cake grid
 */
#ifndef DUMUX_CAKE_GRID_MANAGER_HH
#define DUMUX_CAKE_GRID_MANAGER_HH

#include <vector>
#include <iostream>
#include <cmath>
#include <algorithm>

#include <dune/common/dynvector.hh>
#include <dune/common/float_cmp.hh>
#include <dune/grid/common/gridfactory.hh>
#include <dumux/common/parameters.hh>

namespace Dumux {

/*!
 * \ingroup InputOutput
 * \brief Provides a grid manager with a method for creating creating vectors
 *        with polar Coordinates and one for creating a cartesian grid from
 *        these polar coordinates.
 */
template <class Grid>
class CakeGridManager
{
    using Scalar = typename Grid::ctype;

    using GridFactory = Dune::GridFactory<Grid>;
    using GridPointer = std::shared_ptr<Grid>;

    enum { dim = Grid::dimension,
           dimWorld = Grid::dimensionworld };

public:
    /*!
     * \brief Make the grid.
     */
    void init(const std::string& modelParamGroup = "")
    {
        static_assert(dim == 2 || dim == 3, "The CakeGridManager is only implemented for 2D and 3D.");

        const bool verbose = getParamFromGroup<bool>(modelParamGroup, "Grid.Verbosity", false);

        std::array<std::vector<Scalar>, dim> polarCoordinates;
        // Indices specifing in which direction the piece of cake is oriented
        Dune::FieldVector<int, dim> indices(-1);
        createVectors(polarCoordinates, indices, modelParamGroup, verbose);

        gridPtr() = createCakeGrid(polarCoordinates, indices, modelParamGroup, verbose);
        loadBalance();
    }

    /*!
     * \brief Create vectors containing polar coordinates of all points.
     *
     * All keys are expected to be in group GridParameterGroup.
     * The following keys are recognized:
     * - Radial : min/max value for radial coordinate
     * - Angular : min/max value for angular coordinate
     * - Axial : min/max value for axial coordinate
     *   Adding 0, 1 (or 2 in 3D) specifies in which direction (x, y and z, respectively)
     *   the radial, angular and axial direction are oriented
     * - Cells : number of cells array for x-coordinate (Again, an added 0, 1 or 3 specifies x, y or z
     * - Grading : grading factor array for x-coordinate (Same here)
     * - Verbosity : whether the grid construction should output to standard out
     *
     * The grading factor \f$ g \f$ specifies the ratio between the next and the current cell size:
     * \f$ g = \frac{h_{i+1}}{h_i} \f$.
     * Negative grading factors are converted to
     * \f$ g = -\frac{1}{g_\textrm{negative}} \f$
     * to avoid issues with imprecise fraction numbers.
     */
    static void createVectors(std::array<std::vector<Scalar>, dim> &polarCoordinates,
                              Dune::FieldVector<int, dim> &indices,
                              const std::string& modelParamGroup,
                              bool verbose = false)
    {
        // The positions
        std::array<std::vector<Scalar>, dim> positions;

        for (int i = 0; i < dim; ++i)
        {
            const bool hasRadial = hasParamInGroup(modelParamGroup, "Grid.Radial" + std::to_string(i));
            const bool hasAngular = hasParamInGroup(modelParamGroup, "Grid.Angular" + std::to_string(i));
            const bool hasAxial = (dim == 3) && hasParamInGroup(modelParamGroup, "Grid.Axial" + std::to_string(i));
            if (static_cast<int>(hasRadial) + static_cast<int>(hasAngular) + static_cast<int>(hasAxial) != 1)
                DUNE_THROW(Dune::RangeError, "Multiple or no position vectors (radial, angular, axial) specified in coord direction: " << i << std::endl);

            if (hasRadial)
            {
                positions[i] = getParamFromGroup<std::vector<Scalar>>(modelParamGroup, "Grid.Radial" + std::to_string(i));
                indices[0] = i; // Index specifing radial direction
            }

            else if (hasAngular)
            {
                positions[i] = getParamFromGroup<std::vector<Scalar>>(modelParamGroup, "Grid.Angular" + std::to_string(i));
                indices[1] = i; // Index specifing angular direction
            }

            else // hasAxial
            {
                positions[i] = getParamFromGroup<std::vector<Scalar>>(modelParamGroup, "Grid.Axial" + std::to_string(i));
                indices[2] = i; // Index specifing axial direction
            }

            if (!std::is_sorted(positions[i].begin(), positions[i].end()))
                DUNE_THROW(Dune::GridError, "Make sure to specify a monotone increasing \"Positions\" array");
        }

        // check that all indices are specified i.e. > 0
        for (int i = 0; i < dim; ++i)
            if (indices[i] < 0)
                DUNE_THROW(Dune::RangeError, "Please specify Positions Angular and Radial and Axial correctly and unambiguously!" << std::endl);

        // the number of cells (has a default)
        std::array<std::vector<int>, dim> cells;
        for (int i = 0; i < dim; ++i)
        {
            cells[i].resize(positions[i].size()-1, 1.0);
            cells[i] = getParamFromGroup<std::vector<int>>(modelParamGroup, "Grid.Cells" +  std::to_string(i), cells[i]);
            if (cells[i].size() + 1 != positions[i].size())
                DUNE_THROW(Dune::RangeError, "Make sure to specify the correct length of \"Cells\" and \"Positions\" arrays");
        }

        // grading factor (has a default)
        std::array<std::vector<Scalar>, dim> grading;
        for (int i = 0; i < dim; ++i)
        {
            grading[i].resize(positions[i].size()-1, 1.0);
            grading[i] = getParamFromGroup<std::vector<Scalar>>(modelParamGroup, "Grid.Grading" +  std::to_string(i), grading[i]);
            if (grading[i].size() + 1 != positions[i].size())
                DUNE_THROW(Dune::RangeError, "Make sure to specify the correct length of \"Grading\" and \"Positions\" arrays");
        }

        // make the grid
        std::array<std::vector<Scalar>, dim> globalPositions;
        using std::pow;
        for (int dimIdx = 0; dimIdx < dim; dimIdx++)
        {
163
164
            // Each grid direction is subdivided into (numCells + 1) points
            std::size_t numGlobalPositions = 1;
165
            for (int zoneIdx = 0; zoneIdx < cells[dimIdx].size(); ++zoneIdx)
166
                numGlobalPositions += cells[dimIdx][zoneIdx];
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498

            globalPositions[dimIdx].resize(numGlobalPositions);
            std::size_t posIdx = 0;
            for (int zoneIdx = 0; zoneIdx < cells[dimIdx].size(); ++zoneIdx)
            {
                const Scalar lower = positions[dimIdx][zoneIdx];
                const Scalar upper = positions[dimIdx][zoneIdx+1];
                const int numCells = cells[dimIdx][zoneIdx];
                Scalar gradingFactor = grading[dimIdx][zoneIdx];
                const Scalar length = upper - lower;
                Scalar height = 1.0;
                bool increasingCellSize = false;

                if (verbose)
                {
                    std::cout << "dim " << dimIdx
                              << " lower "  << lower
                              << " upper "  << upper
                              << " numCells "  << numCells
                              << " grading "  << gradingFactor;
                }

                if (gradingFactor > 1.0)
                    increasingCellSize = true;

                // take absolute values and reverse cell size increment to achieve
                // reverse behavior for negative values
                if (gradingFactor < 0.0)
                {
                    using std::abs;
                    gradingFactor = abs(gradingFactor);

                    if (gradingFactor < 1.0)
                        increasingCellSize = true;
                }

                // if the grading factor is exactly 1.0 do equal spacing
                if (gradingFactor > 1.0 - 1e-7 && gradingFactor < 1.0 + 1e-7)
                {
                    height = 1.0 / numCells;
                    if (verbose)
                        std::cout << " -> h "  << height * length << std::endl;
                }

                // if grading factor is not 1.0, do power law spacing
                else
                {
                    height = (1.0 - gradingFactor) / (1.0 - pow(gradingFactor, numCells));

                    if (verbose)
                    {
                        std::cout << " -> grading_eff "  << gradingFactor
                                  << " h_min "  << height * pow(gradingFactor, 0) * length
                                  << " h_max "  << height * pow(gradingFactor, numCells-1) * length
                                  << std::endl;
                    }
                }

                // the positions inside a global segment
                Dune::DynamicVector<Scalar> localPositions(numCells, 0.0);
                for (int i = 1; i < numCells; i++)
                {
                    Scalar hI = height;
                    if (!(gradingFactor < 1.0 + 1e-7 && gradingFactor > 1.0 - 1e-7))
                    {
                        if (increasingCellSize)
                            hI *= pow(gradingFactor, i-1);

                        else
                            hI *= pow(gradingFactor, numCells-i);
                    }
                    localPositions[i] = localPositions[i-1] + hI;
                }

                localPositions *= length;
                localPositions += lower;

                for (int i = 0; i < numCells; ++i)
                    globalPositions[dimIdx][posIdx++] = localPositions[i];
            }

            globalPositions[dimIdx][posIdx] = positions[dimIdx].back();
        }

        polarCoordinates[0] = globalPositions[indices[0]];
        polarCoordinates[1] = globalPositions[indices[1]];
        if (dim == 3)
            polarCoordinates[2] = globalPositions[dim - indices[0] - indices[1]];

        // convert angular coordinates into radians
        std::transform(polarCoordinates[1].begin(), polarCoordinates[1].end(),
                       polarCoordinates[1].begin(),
                       [](Scalar s){ return s*M_PI/180; });
    }

    /*!
     * \brief Creates cartesian grid from polar coordinates.
     *
     * \param polarCoordinates Vector containing radial, angular and axial coordinates (in this order)
     * \param indices Indices specifing the radial, angular and axial direction (in this order)
     * \param modelParamGroup name of the model parameter group
     * \param verbose if the output should be verbose
     */
    std::unique_ptr<Grid> createCakeGrid(std::array<std::vector<Scalar>, dim> &polarCoordinates,
                                         Dune::FieldVector<int, dim> &indices,
                                         const std::string& modelParamGroup,
                                         bool verbose = false)
    {
        std::vector<Scalar> dR = polarCoordinates[0];
        std::vector<Scalar> dA = polarCoordinates[1];

        // For the special case of 36o°, the last element is connected with the first one.
        // Thus, one needs to distinguish here between all other cases, were the normal procedure
        // is applied for all elements until the last one  (= dA.size() - 1) and the 360° case,
        // where it is only applied until the second to last element (dA.size() - 2).
        int maxdA = dA.size() - 1;
        if (Dune::FloatCmp::eq(polarCoordinates[1][polarCoordinates[1].size()-1], 360.0))
        {
            maxdA = dA.size() - 2;
        }

        GridFactory gridFactory;
        constexpr auto type = Dune::GeometryTypes::cube(dim);

        // create nodes
        if (dim == 3)
        {
            std::vector<Scalar> dZ = polarCoordinates[2];
            for (int j = 0; j <= dA.size() - 1; ++j)
            {
                for (int l = 0; l <= dZ.size() - 1; ++l)
                {
                    for (int i = 0; i <= dR.size()- 1; ++i)
                    {
                        // Get radius for the well (= a hole) in the center
                        const auto wellRadius = getParamFromGroup<Scalar>(modelParamGroup, "Grid.WellRadius");

                        // transform into cartesian Coordinates
                        using std::cos;
                        using std::sin;
                        Dune::FieldVector <double, dim> v(0.0);
                        v[indices[2]] = dZ[l];
                        v[indices[0]] = cos(dA[j])*wellRadius + cos(dA[j])*dR[i];
                        v[indices[1]] = sin(dA[j])*wellRadius + sin(dA[j])*dR[i];
                        if (verbose)
                        {
                            std::cout << "Coordinates of : "
                                      << v[0] << " " << v[1] << " " << v[2] << std::endl;
                        }
                        gridFactory.insertVertex(v);
                    }
                }
            }

            if (verbose)
                std::cout << "Filled node vector" << std::endl;

            // assign nodes
            unsigned int z = 0;
            unsigned int t = 0;
            for (int j = 0; j < dA.size() - 1; ++j)
            {
                for (int l = 0; l < dZ.size() - 1; ++l)
                {
                    if (j < maxdA)
                    {
                        for (int i = 0; i < dR.size() - 1; ++i)
                        {
                            unsigned int rSize = dR.size();
                            unsigned int zSize = dZ.size();
                            std::vector<unsigned int> vid({z, z+1, z+rSize*zSize,
                                                           z+rSize*zSize+1, z+rSize, z+rSize+1,
                                                           z+rSize*zSize+rSize, z+rSize*zSize+rSize+1});

                            if (verbose)
                            {
                                std::cout << "element vertex ids: ";
                                for (int k = 0; k < vid.size(); ++k)
                                    std::cout << vid[k] << " ";
                                std::cout << std::endl;
                            }

                            gridFactory.insertElement(type, vid);

                            z = z+1;
                        }
                        z = z+1;
                    }
                    else
                    {
                        // assign nodes for 360°-cake
                        for (int i = 0; i < dR.size() - 1; ++i)
                        {
                            // z = z + 1;
                            unsigned int rSize = dR.size();
                            std::vector<unsigned int> vid({z, z+1, t,
                                                           t+1, z+rSize, z+rSize+1,
                                                           t+rSize, t+rSize+1});

                            if (verbose)
                            {
                                std::cout << "element vertex ids: ";
                                for (int k = 0; k < vid.size(); ++k)
                                    std::cout << vid[k] << " ";
                                std::cout << std::endl;
                            }

                            gridFactory.insertElement(type, vid);
                            t = t + 1;
                            z = z+1;
                        }
                        t = t + 1;
                        z = z+1;
                    }
                    if (verbose)
                        std::cout << "assign nodes 360° ends..." << std::endl;
                }

                z = z + dR.size();
            }
        }

        // for dim = 2
        else
        {
            for (int j = 0; j <= dA.size() - 1; ++j)
            {
                for (int i = 0; i <= dR.size()- 1; ++i)
                {
                    // Get radius for the well (= a hole) in the center
                    const Scalar wellRadius = getParamFromGroup<Scalar>(modelParamGroup, "Grid.WellRadius");

                    // transform into cartesian Coordinates
                    Dune::FieldVector <double, dim> v(0.0);

                    v[indices[0]] = cos(dA[j])*wellRadius + cos(dA[j])*dR[i];
                    v[indices[1]] = sin(dA[j])*wellRadius + sin(dA[j])*dR[i];
                    if(verbose) std::cout << "Coordinates of : " << v[0] << " " << v[1] << std::endl;
                    gridFactory.insertVertex(v);
                }
            }
            std::cout << "Filled node vector" << std::endl;

            // assign nodes
            unsigned int z = 0;
            unsigned int t = 0;
            for (int j = 0; j < dA.size() - 1; ++j)
            {
                if (j < maxdA)
                {
                    for (int i = 0; i < dR.size() - 1; ++i)
                    {
                        unsigned int rSize = dR.size();
                        std::vector<unsigned int> vid({z, z+1, z+rSize, z+rSize+1});

                        if (verbose)
                        {
                            std::cout << "element vertex ids: ";
                            for (int k = 0; k < vid.size(); ++k)
                                std::cout << vid[k] << " ";
                            std::cout << std::endl;
                        }

                        gridFactory.insertElement(type, vid);
                        z = z+1;
                    }
                    z = z+1;
                }
                else
                {
                    // assign nodes for 360°-cake
                    for (int i = 0; i < dR.size() - 1; ++i)
                    {
                        // z = z + 1;
                        std::vector<unsigned int> vid({z, z+1, t, t+1});

                        if (verbose)
                        {
                            std::cout << "element vertex ids: ";
                            for (int k = 0; k < vid.size(); ++k)
                                std::cout << vid[k] << " ";
                            std::cout << std::endl;
                        }

                        gridFactory.insertElement(type, vid);
                        t = t + 1;
                        z = z+1;
                    }
                    t = t + 1;
                    z = z+1;
                }
                std::cout << "assign nodes 360 ends..." << std::endl;
            }
        }
        // return the grid pointer
        return std::unique_ptr<Grid>(gridFactory.createGrid());
    }

    /*!
     * \brief Returns a reference to the grid.
     */
    Grid& grid()
    {
        return *gridPtr();
    }

    /*!
     * \brief Distributes the grid on all processes of a parallel
     *        computation.
     */
    void loadBalance()
    {
        gridPtr()->loadBalance();
    }

protected:

    /*!
     * \brief Returns a reference to the shared pointer to the grid.
     */
    GridPointer& gridPtr()
    {
        return cakeGrid_;
    }

private:
    GridPointer cakeGrid_;
};

} // end namespace Dumux

#endif