faceReflecting.C 16.5 KB
Newer Older
sergio's avatar
ENH:    
sergio committed
1
2
3
4
/*---------------------------------------------------------------------------*\
  =========                 |
  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
   \\    /   O peration     |
5
    \\  /    A nd           | Copyright (C) 2018-2019 OpenCFD Ltd.
sergio's avatar
ENH:    
sergio committed
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
     \\/     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 "faceReflecting.H"
#include "boundaryRadiationProperties.H"
#include "cyclicAMIPolyPatch.H"
#include "volFields.H"


using namespace Foam::constant::mathematical;

// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //

namespace Foam
{
    defineTypeNameAndDebug(faceReflecting, 0);
}

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

void Foam::faceReflecting::initialise(const dictionary& coeffs)
{

    forAll(qreflective_, bandI)
    {
        qreflective_.set
        (
            bandI,
            new volScalarField
            (
                IOobject
                (
                    "qreflective_" + Foam::name(bandI) ,
                    mesh_.time().timeName(),
                    mesh_,
                    IOobject::NO_READ,
                    IOobject::AUTO_WRITE
                ),
                mesh_,
                dimensionedScalar(dimMass/pow3(dimTime), Zero)
            )
        );
    }

    label rayI = 0;
    if (mesh_.nSolutionD() == 3)
    {
        nRay_ = 4*nPhi_*nTheta_;
        refDiscAngles_.resize(nRay_);
        const scalar deltaPhi = pi/(2.0*nPhi_);
        const scalar deltaTheta = pi/nTheta_;

        for (label n = 1; n <= nTheta_; n++)
        {
            for (label m = 1; m <= 4*nPhi_; m++)
            {
                const scalar thetai = (2*n - 1)*deltaTheta/2.0;
                const scalar phii = (2*m - 1)*deltaPhi/2.0;

                scalar sinTheta = Foam::sin(thetai);
                scalar cosTheta = Foam::cos(thetai);
                scalar sinPhi = Foam::sin(phii);
                scalar cosPhi = Foam::cos(phii);
                refDiscAngles_[rayI++] =
                    vector(sinTheta*sinPhi, sinTheta*cosPhi, cosTheta);

            }
        }

    }
    else if (mesh_.nSolutionD() == 2)
    {
        nRay_ = 4*nPhi_;
        refDiscAngles_.resize(nRay_);
        const scalar thetai = piByTwo;
        //const scalar deltaTheta = pi;
        const scalar deltaPhi = pi/(2.0*nPhi_);
        for (label m = 1; m <= 4*nPhi_; m++)
        {
            const scalar phii = (2*m - 1)*deltaPhi/2.0;

            scalar sinTheta = Foam::sin(thetai);
            scalar cosTheta = Foam::cos(thetai);
            scalar sinPhi = Foam::sin(phii);
            scalar cosPhi = Foam::cos(phii);

            refDiscAngles_[rayI++] =
                vector(sinTheta*sinPhi, sinTheta*cosPhi, cosTheta);
        }
    }
    else
    {
        FatalErrorInFunction
            << "The reflected rays are available in 2D or 3D "
            << abort(FatalError);
    }

    const polyBoundaryMesh& patches = mesh_.boundaryMesh();

    const radiation::boundaryRadiationProperties& boundaryRadiation =
        radiation::boundaryRadiationProperties::New(mesh_);

    // global face index
    globalIndex globalNumbering(mesh_.nFaces());

128

sergio's avatar
ENH:    
sergio committed
129
130
131
    // Collect faces with t = 0, r = 0 and a > 0 to shoot rays
    // and patches to construct the triSurface
    DynamicList<point> dynCf;
132
    DynamicList<vector> dynNf;
sergio's avatar
ENH:    
sergio committed
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
    DynamicList<label> dynFacesI;
    forAll(patches, patchI)
    {
        const polyPatch& pp = patches[patchI];
        const pointField& cf = pp.faceCentres();

        if (!pp.coupled() && !isA<cyclicAMIPolyPatch>(pp))
        {
            const tmp<scalarField> tt =
                boundaryRadiation.transmissivity(patchI);

            const tmp<scalarField> tr =
                boundaryRadiation.specReflectivity(patchI);

            const tmp<scalarField> ta =
                boundaryRadiation.absorptivity(patchI);

            const scalarField& t = tt();
            const scalarField& r = tr();
            const scalarField& a = ta();

154
155
            const vectorField& n = pp.faceNormals();

sergio's avatar
ENH:    
sergio committed
156
157
            forAll(pp, faceI)
            {
158
                //const vector nf(n[faceI]);
sergio's avatar
ENH:    
sergio committed
159
                // Opaque, non-reflective, absortived faces to shoot
160
161
162
163
164
165
                if
                (
                    t[faceI] == 0
                 && r[faceI] == 0
                 && a[faceI] > 0
                )
sergio's avatar
ENH:    
sergio committed
166
167
168
                {
                    dynFacesI.append(faceI + pp.start());
                    dynCf.append(cf[faceI]);
169
                    dynNf.append(n[faceI]);
sergio's avatar
ENH:    
sergio committed
170
171
172
                }

                // relfective opaque patches to build reflective surface
173
174
175
176
177
178
                // plus opaque non-reflective
                if
                (
                    (r[faceI] > 0 && t[faceI] == 0) ||
                    (t[faceI] == 0 && a[faceI] > 0 && r[faceI] == 0)
                )
sergio's avatar
ENH:    
sergio committed
179
180
181
182
183
184
185
186
187
                {
                    includePatches_.insert(patchI);
                }
            }
        }
    }

    shootFacesIds_.reset(new labelList(dynFacesI));
    Cfs_.reset(new pointField(dynCf));
188
    Nfs_.reset(new vectorField(dynNf));
sergio's avatar
ENH:    
sergio committed
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

    // * * * * * * * * * * * * * * *
    // Create distributedTriSurfaceMesh
    Random rndGen(653213);

    // Determine mesh bounding boxes:
    List<treeBoundBox> meshBb
    (
        1,
        treeBoundBox
        (
            boundBox(mesh_.points(), false)
        ).extend(rndGen, 1e-3)
    );

    // Dummy bounds dictionary
    dictionary dict;
    dict.add("bounds", meshBb);
    dict.add
    (
        "distributionType",
        distributedTriSurfaceMesh::distributionTypeNames_
        [
            distributedTriSurfaceMesh::FROZEN
        ]
    );
    dict.add("mergeDistance", SMALL);


    triSurface localSurface = triSurfaceTools::triangulate
    (
        mesh_.boundaryMesh(),
        includePatches_,
        mapTriToGlobal_
    );

    surfacesMesh_.reset
    (
        new distributedTriSurfaceMesh
        (
            IOobject
            (
                "reflectiveSurface.stl",
                mesh_.time().constant(),
                "triSurface",
                mesh_.time(),
                IOobject::NO_READ,
                IOobject::NO_WRITE
            ),
            localSurface,
            dict
        )
    );

    if (debug)
    {
        surfacesMesh_->searchableSurface::write();
    }
}

249

sergio's avatar
ENH:    
sergio committed
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
void Foam::faceReflecting::calculate()
{
    const radiation::boundaryRadiationProperties& boundaryRadiation =
        radiation::boundaryRadiationProperties::New(mesh_);

    label nFaces = 0;

    const polyBoundaryMesh& patches = mesh_.boundaryMesh();

    const fvBoundaryMesh& fvPatches = mesh_.boundary();

    label nBands = spectralDistribution_.size();

    // Collect reflected directions from reflecting surfaces on direct hit
    // faces
    const vector sunDir = directHitFaces_.direction();
    const labelList& directHits = directHitFaces_.rayStartFaces();

    globalIndex globalNumbering(mesh_.nFaces());

    Map<label> refFacesDirIndex;
    labelList refDisDirsIndex(nRay_, -1);

    forAll(patches, patchI)
    {
        const polyPatch& pp = patches[patchI];

        if (!pp.coupled() && !isA<cyclicAMIPolyPatch>(pp))
        {
            const tmp<scalarField> tr =
                boundaryRadiation.specReflectivity(patchI);

            const scalarField& r = tr();
            const vectorField n(fvPatches[patchI].nf());

            forAll(pp, faceI)
            {
                label globalID = faceI + pp.start();

                if (r[faceI] > 0.0 && directHits.found(globalID))
                {
                    vector refDir =
                        sunDir + 2.0*(-sunDir & n[faceI]) * n[faceI];

                    // Look for the discrete direction
                    scalar dev(-GREAT);
                    label rayIndex = -1;
                    forAll(refDiscAngles_, iDisc)
                    {
                        scalar dotProd = refDir & refDiscAngles_[iDisc];
                        if (dev < dotProd)
                        {
                            dev = dotProd;
                            rayIndex = iDisc;
                        }
                    }

                    if (rayIndex >= 0)
                    {
                        if (refDisDirsIndex[rayIndex] == -1)
                        {
                            refDisDirsIndex[rayIndex] = 1;
                        }
                    }

                    refFacesDirIndex.insert
                    (
                        globalNumbering.toGlobal(globalID),
                        rayIndex
                    );

                    nFaces++;
                }
            }
        }
    }

    // Distribute ray indexes to all proc's
    Pstream::listCombineGather(refDisDirsIndex, maxEqOp<label>());
    Pstream::listCombineScatter(refDisDirsIndex);

    // Make sure all the processors have the same map
    Pstream::mapCombineGather(refFacesDirIndex, minEqOp<label>());
    Pstream::mapCombineScatter(refFacesDirIndex);

    scalar maxBounding = 5.0*mag(mesh_.bounds().max() - mesh_.bounds().min());

    reduce(maxBounding, maxOp<scalar>());

    // Shoot Rays
    // From faces t = 0, r = 0 and a > 0 to all 'used' discrete reflected
    // directions

    DynamicField<point> start(nFaces);
    DynamicField<point> end(start.size());
    DynamicList<label> startIndex(start.size());
    DynamicField<label> dirStartIndex(start.size());

    label i = 0;
    do
    {
        for (; i < Cfs_->size(); i++)
        {
            const point& fc = Cfs_()[i];

355
356
            const vector nf = Nfs_()[i];

sergio's avatar
ENH:    
sergio committed
357
358
            const label myFaceId = shootFacesIds_()[i];

359
            forAll(refDisDirsIndex, dirIndex)
sergio's avatar
ENH:    
sergio committed
360
361
362
            {
                if (refDisDirsIndex[dirIndex] > -1)
                {
363
                    if ((nf & refDiscAngles_[dirIndex]) > 0)
364
365
                    {
                        const vector direction = -refDiscAngles_[dirIndex];
sergio's avatar
ENH:    
sergio committed
366

367
                        start.append(fc + 0.001*direction);
sergio's avatar
ENH:    
sergio committed
368

369
370
                        startIndex.append(myFaceId);
                        dirStartIndex.append(dirIndex);
sergio's avatar
ENH:    
sergio committed
371

372
373
                        end.append(fc + maxBounding*direction);
                    }
sergio's avatar
ENH:    
sergio committed
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
                }
            }
        }

    }while (returnReduce(i < Cfs_->size(), orOp<bool>()));

    List<pointIndexHit> hitInfo(startIndex.size());

    surfacesMesh_->findLine(start, end, hitInfo);

    // Query the local trigId on hit faces
    labelList triangleIndex;
    autoPtr<mapDistribute> mapPtr
    (
        surfacesMesh_->localQueries
        (
            hitInfo,
            triangleIndex
        )
    );
    const mapDistribute& map = mapPtr();

396
397
398
    PtrList<List<scalarField>> patchr(patches.size());
    PtrList<List<scalarField>> patcha(patches.size());
    forAll(patchr, patchi)
sergio's avatar
ENH:    
sergio committed
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
    {
        patchr.set
        (
            patchi,
            new List<scalarField>(nBands)
        );

        patcha.set
        (
            patchi,
            new List<scalarField>(nBands)
        );
    }

    // Fill patchr
414
    forAll(patchr, patchi)
sergio's avatar
ENH:    
sergio committed
415
    {
416
417
418
        const polyPatch& pp = patches[patchi];

        if (!pp.coupled() && !isA<cyclicAMIPolyPatch>(pp))
sergio's avatar
ENH:    
sergio committed
419
        {
420
421
422
423
424
425
426
427
428
            for (label bandI = 0; bandI < nBands; bandI++)
            {
                patchr[patchi][bandI] =
                    boundaryRadiation.specReflectivity
                    (
                        patchi,
                        bandI,
                        new vectorField(patches[patchi].size(), sunDir)
                    );
sergio's avatar
ENH:    
sergio committed
429

430
431
432
433
434
435
436
437
                patcha[patchi][bandI] =
                    boundaryRadiation.absorptivity
                    (
                        patchi,
                        bandI,
                        new vectorField(patches[patchi].size(), sunDir)
                    );
            }
sergio's avatar
ENH:    
sergio committed
438
439
440
        }
    }

441
    List<scalarField> r(nBands);
442
    for (label bandI = 0; bandI < nBands; bandI++)
443
444
445
446
447
448
    {
        r[bandI].setSize(triangleIndex.size());
    }
    labelList refDirIndex(triangleIndex.size());
    labelList refIndex(triangleIndex.size());
    // triangleIndex includes hits on non-reflecting and reflecting faces
sergio's avatar
ENH:    
sergio committed
449
450
451
452
453
454
455
456
457
458
459
460
    forAll(triangleIndex, i)
    {
        label trii = triangleIndex[i];
        label facei = mapTriToGlobal_[trii];
        label patchI = patches.whichPatch(facei);
        const polyPatch& pp = patches[patchI];
        label localFaceI = pp.whichFace(facei);

        label globalFace = globalNumbering.toGlobal(Pstream::myProcNo(), facei);
        if (refFacesDirIndex.found(globalFace))
        {
            refDirIndex[i] = refFacesDirIndex.find(globalFace)();
461
            refIndex[i] = globalFace;
sergio's avatar
ENH:    
sergio committed
462
463
464
465
466
467
468
        }
        for (label bandI = 0; bandI < nBands; bandI++)
        {
            r[bandI][i] = patchr[patchI][bandI][localFaceI];
        }
    }
    map.reverseDistribute(hitInfo.size(), refDirIndex);
469
    map.reverseDistribute(hitInfo.size(), refIndex);
sergio's avatar
ENH:    
sergio committed
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
    for (label bandI = 0; bandI < nBands; bandI++)
    {
        map.reverseDistribute(hitInfo.size(), r[bandI]);
    }

    for (label bandI = 0; bandI < nBands; bandI++)
    {
        volScalarField::Boundary& qrefBf =
            qreflective_[bandI].boundaryFieldRef();
        qrefBf = 0.0;
    }

    const vector qPrim(solarCalc_.directSolarRad()*solarCalc_.direction());

    // Collect rays with a hit (hitting reflecting surfaces)
    // and whose reflected direction are equal to the shot ray
    forAll(hitInfo, rayI)
    {
        if (hitInfo[rayI].hit())
        {
490
491
            if
            (
492
                dirStartIndex[rayI] == refDirIndex[rayI]
493
494
             && refFacesDirIndex.found(refIndex[rayI])
            )
sergio's avatar
ENH:    
sergio committed
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
            {
                for (label bandI = 0; bandI < nBands; bandI++)
                {
                    volScalarField::Boundary& qrefBf =
                        qreflective_[bandI].boundaryFieldRef();

                    label startFaceId = startIndex[rayI];
                    label startPatchI = patches.whichPatch(startFaceId);

                    const polyPatch& ppStart = patches[startPatchI];
                    label localStartFaceI = ppStart.whichFace(startFaceId);

                    scalar a = patcha[startPatchI][bandI][localStartFaceI];

                    const vectorField& nStart = ppStart.faceNormals();

511
512
513
514
                    vector rayIn = refDiscAngles_[dirStartIndex[rayI]];

                    rayIn /= mag(rayIn);

sergio's avatar
ENH:    
sergio committed
515
                    qrefBf[startPatchI][localStartFaceI] +=
516
517
                    (
                        (
518
519
520
521
522
                            mag(qPrim)
                           *r[bandI][rayI]
                           *spectralDistribution_[bandI]
                           *a
                           *rayIn
523
524
525
                        )
                        & nStart[localStartFaceI]
                    );
sergio's avatar
ENH:    
sergio committed
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
                }
            }
        }
    }

    start.clear();
    startIndex.clear();
    end.clear();
    dirStartIndex.clear();
}


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

Foam::faceReflecting::faceReflecting
(
    const fvMesh& mesh,
    const faceShading& directHiyFaces,
    const solarCalculator& solar,
    const scalarList& spectralDistribution,
    const dictionary& dict
)
:
    mesh_(mesh),
    nTheta_(dict.subDict("reflecting").lookupOrDefault<label>("nTheta", 10)),
    nPhi_(dict.subDict("reflecting").lookupOrDefault<label>("nPhi", 10)),
    nRay_(0),
    refDiscAngles_(0),
    spectralDistribution_(spectralDistribution),
    qreflective_(spectralDistribution_.size()),
    directHitFaces_(directHiyFaces),
    surfacesMesh_(),
    shootFacesIds_(),
    Cfs_(),
560
    Nfs_(),
sergio's avatar
ENH:    
sergio committed
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
    solarCalc_(solar),
    includePatches_(),
    mapTriToGlobal_()
{
    initialise(dict);
}


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

void Foam::faceReflecting::correct()
{
    calculate();
}


// ************************************************************************* //