PDRblock.C 14.8 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
/*---------------------------------------------------------------------------*\
  =========                 |
  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
   \\    /   O peration     |
    \\  /    A nd           | Copyright (C) 2019 OpenCFD Ltd.
     \\/     M anipulation  |
-------------------------------------------------------------------------------
License
    This file is part of OpenFOAM.

    OpenFOAM 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.

    OpenFOAM 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 OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.

\*---------------------------------------------------------------------------*/

#include "PDRblock.H"
27
#include "ListOps.H"
28
29
30
31
32
33
34
35
36
37
38
39
40
41
#include "emptyPolyPatch.H"

// * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //

namespace Foam
{
    //- Calculate geometric expansion factor from expansion ratio
    inline scalar calcGexp(const scalar expRatio, const label nDiv)
    {
        return nDiv > 1 ? pow(expRatio, 1.0/(nDiv - 1)) : 0.0;
    }
}


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
// * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //

bool Foam::PDRblock::checkMonotonic
(
    const direction cmpt,
    const UList<scalar>& pts
)
{
    const label len = pts.size();

    if (!len)
    {
        return false;
    }

    const scalar& minVal = pts[0];

    for (label i=1; i < len; ++i)
    {
        if (pts[i] <= minVal)
        {
            FatalErrorInFunction
                << "Points in " << vector::componentNames[cmpt]
                << " direction do not increase monotonically" << nl
                << flatOutput(pts) << nl << nl
                << exit(FatalError);
        }
    }

    return true;
}


75
76
77
78
79
80
const Foam::PDRblock& Foam::PDRblock::null()
{
    return NullObjectRef<PDRblock>();
}


81
82
// * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //

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
void Foam::PDRblock::adjustSizes()
{
    // Adjust i-j-k addressing
    sizes().x() = grid_.x().nCells();
    sizes().y() = grid_.y().nCells();
    sizes().z() = grid_.z().nCells();

    if (sizes().x() <= 0 || sizes().y() <= 0 || sizes().z() <= 0)
    {
        // Sanity check. Silently disallow bad sizing
        ijkMesh::clear();

        grid_.x().clear();
        grid_.y().clear();
        grid_.z().clear();

        bounds_ = boundBox::invertedBox;
        minEdgeLen_ = Zero;
        return;
    }

    // Adjust boundBox
    bounds_.min().x() = grid_.x().first();
    bounds_.min().y() = grid_.y().first();
    bounds_.min().z() = grid_.z().first();

    bounds_.max().x() = grid_.x().last();
    bounds_.max().y() = grid_.y().last();
    bounds_.max().z() = grid_.z().last();

    // Min edge length
    minEdgeLen_ = GREAT;

    for (direction cmpt=0; cmpt < vector::nComponents; ++cmpt)
    {
        const label nEdge = grid_[cmpt].nCells();

        for (label edgei=0; edgei < nEdge; ++edgei)
        {
            const scalar len = grid_[cmpt].width(edgei);
            minEdgeLen_ = min(minEdgeLen_, len);
        }
    }
}


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
163
164
165
166
167
168
169
170
171
172
173
174
175
176
void Foam::PDRblock::readGridControl
(
    const direction cmpt,
    const dictionary& dict,
    const scalar scaleFactor
)
{
    //- The begin/end nodes for each segment
    scalarList knots;

    // The number of division per segment
    labelList divisions;

    // The expansion ratio per segment
    scalarList expansion;

    if (verbose_)
    {
        Info<< "Reading grid control for "
            << vector::componentNames[cmpt] << " direction" << nl;
    }

    // Points
    // ~~~~~~

    dict.readEntry("points", knots);

    const label nSegments = (knots.size()-1);

    if (nSegments < 1)
    {
        FatalIOErrorInFunction(dict)
            << "Must define at least two control points:"
            << " in " << vector::componentNames[cmpt]
            << " direction" << nl
            << exit(FatalIOError);
    }

    // Apply scaling
    if (scaleFactor > 0)
    {
        for (scalar& pt : knots)
        {
            pt *= scaleFactor;
        }
    }

    // Do points increase monotonically?
177
    checkMonotonic(cmpt, knots);
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

    if (verbose_)
    {
        Info<< "    points : " << flatOutput(knots) << nl;
    }


    // Divisions
    // ~~~~~~~~~

    dict.readEntry("nCells", divisions);

    label nTotalDivs = 0;
    for (const label ndiv : divisions)
    {
        nTotalDivs += ndiv;

        if (ndiv < 1)
        {
            FatalIOErrorInFunction(dict)
                << "Negative or zero divisions:"
                << " in " << vector::componentNames[cmpt]
                << " direction" << nl
                << exit(FatalIOError);
        }
    }

    if (divisions.size() != nSegments)
    {
        FatalIOErrorInFunction(dict)
            << "Mismatch in number of divisions and number of points:"
            << " in " << vector::componentNames[cmpt]
            << " direction" << nl
            << exit(FatalIOError);
    }
    else if (!nTotalDivs)
    {
        FatalIOErrorInFunction(dict)
            << "No grid divisions at all:"
            << " in " << vector::componentNames[cmpt]
            << " direction" << nl
            << exit(FatalIOError);
    }


    if (verbose_)
    {
        Info<< "    nCells : " << flatOutput(divisions) << nl;
    }


    // Expansion ratios
    // ~~~~~~~~~~~~~~~~

    dict.readIfPresent("ratios", expansion);

    // Correct expansion ratios - negative is the same as inverse.
    for (scalar& expRatio : expansion)
    {
        if (expRatio < 0)
        {
            expRatio = 1.0/(-expRatio);
        }
    }

    if (expansion.size() && expansion.size() != nSegments)
    {
        FatalIOErrorInFunction(dict)
            << "Mismatch in number of expansion ratios and divisions:"
            << " in " << vector::componentNames[cmpt]
            << " direction" << nl
            << exit(FatalIOError);
    }

    if (expansion.empty())
    {
        expansion.resize(nSegments, scalar(1));

        if (verbose_)
        {
            Info<< "Warning: 'ratios' not specified, using constant spacing"
                << nl;
        }
    }
    else
    {
        if (verbose_)
        {
            Info<< "    ratios : " << flatOutput(expansion) << nl;
        }
    }



    // Define the grid points

    scalarList& gridPoint = grid_[cmpt];
    gridPoint.resize(nTotalDivs+1);

    label start = 0;

    for (label segmenti=0; segmenti < nSegments; ++segmenti)
    {
        const label nDiv = divisions[segmenti];
        const scalar expRatio = expansion[segmenti];

        SubList<scalar> subPoint(gridPoint, nDiv+1, start);
        start += nDiv;

        subPoint[0] = knots[segmenti];
        subPoint[nDiv] = knots[segmenti+1];

        const scalar dist = (subPoint.last() - subPoint.first());

        if (equal(expRatio, 1.0))
        {
            for (label i=1; i < nDiv; ++i)
            {
                subPoint[i] = (subPoint[0] + (dist * i)/nDiv);
            }
        }
        else
        {
            // Geometric expansion factor from the expansion ratio
            const scalar expFact = calcGexp(expRatio, nDiv);

            for (label i=1; i < nDiv; ++i)
            {
                subPoint[i] =
                (
                    subPoint[0]
                  + dist * (1.0 - pow(expFact, i))/(1.0 - pow(expFact, nDiv))
                );
            }
        }
    }
}


void Foam::PDRblock::readBoundary(const dictionary& dict)
{
    if (verbose_)
    {
        Info<< "Reading boundary entries" << nl;
    }

    PtrList<entry> patchEntries;

    if (dict.found("boundary"))
    {
        PtrList<entry> newEntries(dict.lookup("boundary"));
        patchEntries.transfer(newEntries);
    }

    patches_.clear();
    patches_.resize(patchEntries.size());

    // Hex cell has 6 sides
    const labelRange validFaceId(0, 6);

    // Track which sides have been handled
    labelHashSet handled;

    forAll(patchEntries, patchi)
    {
        if (!patchEntries.set(patchi))
        {
            continue;
        }

        const entry& e = patchEntries[patchi];

        if (!e.isDict())
        {
            FatalIOErrorInFunction(dict)
                << "Entry " << e
                << " in boundary section is not a valid dictionary."
                << exit(FatalIOError);
        }

        const dictionary& patchDict = e.dict();

        const word& patchName = e.keyword();

        const word patchType(patchDict.get<word>("type"));

        labelList faceLabels(patchDict.get<labelList>("faces"));

        labelHashSet patchFaces(faceLabels);

        if (faceLabels.size() != patchFaces.size())
        {
            WarningInFunction
                << "Patch: " << patchName
                << ": Ignoring duplicate face ids" << nl;
        }

        label nErased = patchFaces.filterKeys(validFaceId);
        if (nErased)
        {
            WarningInFunction
                << "Patch: " << patchName << ": Ignoring " << nErased
                << " faces with invalid identifiers" << nl;
        }

        nErased = patchFaces.erase(handled);
        if (nErased)
        {
            WarningInFunction
                << "Patch: " << patchName << ": Ignoring " << nErased
                << " faces, which were already handled" << nl;
        }

        // Mark these faces as having been handled
        handled += patchFaces;


        // Save information for later access during mesh creation.

        patches_.set(patchi, new boundaryEntry());

        boundaryEntry& bentry = patches_[patchi];

        bentry.name_ = patchName;
        bentry.type_ = patchType;
        bentry.size_ = 0;
        bentry.faces_ = patchFaces.sortedToc();

        // Count patch faces
        for (const label shapeFacei : bentry.faces_)
        {
            bentry.size_ += nBoundaryFaces(shapeFacei);
        }
    }


    // Deal with unhandled faces

    labelHashSet missed(identity(6));
    missed.erase(handled);

    if (missed.size())
    {
        patches_.append(new boundaryEntry());

        boundaryEntry& bentry = patches_.last();

        bentry.name_ = "defaultFaces";
        bentry.type_ = emptyPolyPatch::typeName;
        bentry.size_ = 0;
        bentry.faces_ = missed.sortedToc();

        // Count patch faces
        for (const label shapeFacei : bentry.faces_)
        {
            bentry.size_ += nBoundaryFaces(shapeFacei);
        }


        // Use name/type from "defaultPatch" entry if available
        // - be generous with error handling

        const dictionary* defaultEntry = dict.findDict("defaultPatch");
        if (defaultEntry)
        {
            defaultEntry->readIfPresent("name", bentry.name_);
            defaultEntry->readIfPresent("type", bentry.type_);
        }
    }

    if (verbose_)
    {
        label patchi = 0;
        for (const boundaryEntry& bentry : patches_)
        {
            Info<< "    patch: " << patchi
                << " (size: " << bentry.size_
                << " type: " << bentry.type_
                << ") name: " << bentry.name_
                << " faces: " << flatOutput(bentry.faces_) << nl;

            ++patchi;
        }
    }
}


// * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //

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
Foam::PDRblock::PDRblock
(
    const UList<scalar>& xgrid,
    const UList<scalar>& ygrid,
    const UList<scalar>& zgrid
)
:
    PDRblock()
{
    // Default boundaries with patchi == shapeFacei
    patches_.resize(6);
    for (label patchi=0; patchi < 6; ++patchi)
    {
        patches_.set(patchi, new boundaryEntry());

        boundaryEntry& bentry = patches_[patchi];

        bentry.name_ = "patch" + Foam::name(patchi);
        bentry.type_ = "patch";
        bentry.size_ = 0;
        bentry.faces_ = labelList(one(), patchi);
    }

    reset(xgrid, ygrid, zgrid);
}


494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
Foam::PDRblock::PDRblock(const dictionary& dict, bool verbose)
:
    PDRblock()
{
    verbose_ = verbose;

    read(dict);
}


// * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //

bool Foam::PDRblock::read(const dictionary& dict)
{
    // Mark no scaling with invalid value
    const scalar scaleFactor = dict.lookupOrDefault<scalar>("scale", -1);

    readGridControl(0, dict.subDict("x"), scaleFactor);
    readGridControl(1, dict.subDict("y"), scaleFactor);
    readGridControl(2, dict.subDict("z"), scaleFactor);

515
    adjustSizes();
516

517
    readBoundary(dict);
518

519
520
    return true;
}
521
522


523
524
525
526
527
528
529
530
531
532
void Foam::PDRblock::reset
(
    const UList<scalar>& xgrid,
    const UList<scalar>& ygrid,
    const UList<scalar>& zgrid
)
{
    static_cast<scalarList&>(grid_.x()) = xgrid;
    static_cast<scalarList&>(grid_.y()) = ygrid;
    static_cast<scalarList&>(grid_.z()) = zgrid;
533

534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
    #ifdef FULLDEBUG
    for (direction cmpt=0; cmpt < vector::nComponents; ++cmpt)
    {
        checkMonotonic(cmpt, grid_[cmpt]);
    }
    #endif

    adjustSizes();

    // Adjust boundaries
    for (boundaryEntry& bentry : patches_)
    {
        bentry.size_ = 0;

        // Count patch faces
        for (const label shapeFacei : bentry.faces_)
        {
            bentry.size_ += nBoundaryFaces(shapeFacei);
        }
    }
554
555
556
}


557
bool Foam::PDRblock::findCell(const point& pt, labelVector& pos) const
558
559
560
561
562
{
    // Out-of-bounds is handled explicitly, for efficiency and consistency,
    // but principally to ensure that findLower() returns a valid
    // result when the point is to the right of the bounds.

563
564
565
566
    // Since findLower returns the lower index, it corresponds to the
    // cell in which the point is found

    if (!bounds_.contains(pt))
567
    {
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
        return false;
    }

    for (direction cmpt=0; cmpt < labelVector::nComponents; ++cmpt)
    {
        // Binary search
        pos[cmpt] = findLower(grid_[cmpt], pt[cmpt]);
    }

    return true;
}


bool Foam::PDRblock::gridIndex
(
    const point& pt,
    labelVector& pos,
    const scalar relTol
) const
{
    const scalar tol = relTol * minEdgeLen_;

    for (direction cmpt=0; cmpt < labelVector::nComponents; ++cmpt)
    {
        // Linear search
        pos[cmpt] = grid_[cmpt].findIndex(pt[cmpt], tol);

        if (pos[cmpt] < 0) return false;
    }

    return true;
}


Foam::labelVector Foam::PDRblock::findCell(const point& pt) const
{
    labelVector pos;

    if (findCell(pt, pos))
    {
        return pos;
    }

    return labelVector(-1,-1,-1);
}


Foam::labelVector Foam::PDRblock::gridIndex
(
    const point& pt,
    const scalar relTol
) const
{
    labelVector pos;

    if (gridIndex(pt, pos, relTol))
    {
        return pos;
626
627
    }

628
    return labelVector(-1,-1,-1);
629
630
631
}


632
// ************************************************************************* //