Loading...
Searching...
No Matches
faceReflecting.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2018-2025 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28#include "faceReflecting.H"
30#include "cyclicAMIPolyPatch.H"
31#include "volFields.H"
32
34using namespace Foam::constant::mathematical;
35
36// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37
38namespace Foam
39{
41}
42
43// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
44
45void Foam::faceReflecting::initialise(const dictionary& coeffs)
46{
47
48 forAll(qreflective_, bandI)
49 {
50 qreflective_.set
51 (
52 bandI,
54 (
56 (
57 "qreflective_" + Foam::name(bandI) ,
58 mesh_.time().timeName(),
59 mesh_,
63 ),
64 mesh_,
66 )
67 );
68 }
69
70 label rayI = 0;
71 if (mesh_.nSolutionD() == 3)
72 {
73 nRay_ = 4*nPhi_*nTheta_;
74 refDiscAngles_.resize(nRay_);
75 const scalar deltaPhi = pi/(2.0*nPhi_);
76 const scalar deltaTheta = pi/nTheta_;
77
78 for (label n = 1; n <= nTheta_; n++)
79 {
80 for (label m = 1; m <= 4*nPhi_; m++)
81 {
82 const scalar thetai = (2*n - 1)*deltaTheta/2.0;
83 const scalar phii = (2*m - 1)*deltaPhi/2.0;
84
85 scalar sinTheta = Foam::sin(thetai);
86 scalar cosTheta = Foam::cos(thetai);
87 scalar sinPhi = Foam::sin(phii);
88 scalar cosPhi = Foam::cos(phii);
89 refDiscAngles_[rayI++] =
90 vector(sinTheta*sinPhi, sinTheta*cosPhi, cosTheta);
91
92 }
93 }
94
95 }
96 else if (mesh_.nSolutionD() == 2)
97 {
98 nRay_ = 4*nPhi_;
99 refDiscAngles_.resize(nRay_);
100 const scalar thetai = piByTwo;
101 //const scalar deltaTheta = pi;
102 const scalar deltaPhi = pi/(2.0*nPhi_);
103 for (label m = 1; m <= 4*nPhi_; m++)
104 {
105 const scalar phii = (2*m - 1)*deltaPhi/2.0;
106
107 scalar sinTheta = Foam::sin(thetai);
108 scalar cosTheta = Foam::cos(thetai);
109 scalar sinPhi = Foam::sin(phii);
110 scalar cosPhi = Foam::cos(phii);
111
112 refDiscAngles_[rayI++] =
113 vector(sinTheta*sinPhi, sinTheta*cosPhi, cosTheta);
114 }
115 }
116 else
117 {
119 << "The reflected rays are available in 2D or 3D "
120 << abort(FatalError);
121 }
122
123 const polyBoundaryMesh& patches = mesh_.boundaryMesh();
124
125 const radiation::boundaryRadiationProperties& boundaryRadiation =
127
128 // global face index
129 globalIndex globalNumbering(mesh_.nFaces());
130
131
132 // Collect faces with t = 0, r = 0 and a > 0 to shoot rays
133 // and patches to construct the triSurface
134 DynamicList<point> dynCf;
135 DynamicList<vector> dynNf;
136 DynamicList<label> dynFacesI;
137 forAll(patches, patchI)
138 {
139 const polyPatch& pp = patches[patchI];
140 const vectorField::subField cf = pp.faceCentres();
141
142 if (!pp.coupled() && !isA<cyclicAMIPolyPatch>(pp))
143 {
144 const tmp<scalarField> tt =
145 boundaryRadiation.transmissivity(patchI);
146
147 const tmp<scalarField> tr =
148 boundaryRadiation.specReflectivity(patchI);
149
150 const tmp<scalarField> ta =
151 boundaryRadiation.absorptivity(patchI);
152
153 const scalarField& t = tt();
154 const scalarField& r = tr();
155 const scalarField& a = ta();
156
157 const vectorField& n = pp.faceNormals();
158
159 forAll(pp, faceI)
160 {
161 //const vector nf(n[faceI]);
162 // Opaque, non-reflective, absortived faces to shoot
163 if
164 (
165 t[faceI] == 0
166 && r[faceI] == 0
167 && a[faceI] > 0
168 )
169 {
170 dynFacesI.append(faceI + pp.start());
171 dynCf.append(cf[faceI]);
172 dynNf.append(n[faceI]);
173 }
174
175 // relfective opaque patches to build reflective surface
176 // plus opaque non-reflective
177 if
178 (
179 (r[faceI] > 0 && t[faceI] == 0) ||
180 (t[faceI] == 0 && a[faceI] > 0 && r[faceI] == 0)
181 )
182 {
183 includePatches_.insert(patchI);
184 }
185 }
186 }
187 }
188
189 shootFacesIds_.reset(new labelList(dynFacesI));
190 Cfs_.reset(new pointField(dynCf));
191 Nfs_.reset(new vectorField(dynNf));
192
193 // * * * * * * * * * * * * * * *
194 // Create distributedTriSurfaceMesh
195 Random rndGen(653213);
196
197 // Determine mesh bounding boxes:
198 List<treeBoundBox> meshBb
199 (
200 1,
201 treeBoundBox(mesh_.points()).extend(rndGen, 1e-3)
202 );
203
204 // Dummy bounds dictionary
206 dict.add("bounds", meshBb);
207 dict.add
208 (
209 "distributionType",
211 [
213 ]
214 );
215 dict.add("mergeDistance", SMALL);
216
217
219 (
220 mesh_.boundaryMesh(),
221 includePatches_,
222 mapTriToGlobal_
223 );
224
225 surfacesMesh_.reset
226 (
227 new distributedTriSurfaceMesh
228 (
229 IOobject
230 (
231 "reflectiveSurface.stl",
232 mesh_.time().constant(),
233 "triSurface",
234 mesh_.time(),
237 ),
239 dict
240 )
241 );
242
243 if (debug)
244 {
245 surfacesMesh_->searchableSurface::write();
246 }
247}
248
249
250void Foam::faceReflecting::calculate()
251{
252 const radiation::boundaryRadiationProperties& boundaryRadiation =
254
255 label nFaces = 0;
256
257 const polyBoundaryMesh& patches = mesh_.boundaryMesh();
258
259 const fvBoundaryMesh& fvPatches = mesh_.boundary();
260
261 label nBands = spectralDistribution_.size();
262
263 // Collect reflected directions from reflecting surfaces on direct hit
264 // faces
265 const vector sunDir = directHitFaces_.direction();
266 const labelList& directHits = directHitFaces_.rayStartFaces();
267
268 globalIndex globalNumbering(mesh_.nFaces());
269
270 Map<label> refFacesDirIndex;
271 labelList refDisDirsIndex(nRay_, -1);
272
273 forAll(patches, patchI)
274 {
275 const polyPatch& pp = patches[patchI];
276
277 if (!pp.coupled() && !isA<cyclicAMIPolyPatch>(pp))
278 {
279 const tmp<scalarField> tr =
280 boundaryRadiation.specReflectivity(patchI);
281
282 const scalarField& r = tr();
283 const vectorField n(fvPatches[patchI].nf());
284
285 forAll(pp, faceI)
286 {
287 label globalID = faceI + pp.start();
288
289 if (r[faceI] > 0.0 && directHits.found(globalID))
290 {
291 vector refDir =
292 sunDir + 2.0*(-sunDir & n[faceI]) * n[faceI];
293
294 // Look for the discrete direction
295 scalar dev(-GREAT);
296 label rayIndex = -1;
297 forAll(refDiscAngles_, iDisc)
298 {
299 scalar dotProd = refDir & refDiscAngles_[iDisc];
300 if (dev < dotProd)
301 {
302 dev = dotProd;
303 rayIndex = iDisc;
304 }
305 }
306
307 if (rayIndex >= 0)
308 {
309 if (refDisDirsIndex[rayIndex] == -1)
310 {
311 refDisDirsIndex[rayIndex] = 1;
312 }
313 }
314
315 refFacesDirIndex.insert
316 (
317 globalNumbering.toGlobal(globalID),
318 rayIndex
319 );
320
321 nFaces++;
322 }
323 }
324 }
325 }
326
327 // Distribute ray indexes to all proc's
328 // Make sure all the processors have the same information
329
330 Pstream::listReduce(refDisDirsIndex, maxOp<label>());
331 Pstream::mapCombineReduce(refFacesDirIndex, minEqOp<label>());
332
333 const scalar maxBounding =
334 returnReduce(5.0*mesh_.bounds().mag(), maxOp<scalar>());
335
336 // Shoot Rays
337 // From faces t = 0, r = 0 and a > 0 to all 'used' discrete reflected
338 // directions
339
340 DynamicField<point> start(nFaces);
341 DynamicField<point> end(start.size());
342 DynamicList<label> startIndex(start.size());
343 DynamicField<label> dirStartIndex(start.size());
344
345 label i = 0;
346 do
347 {
348 for (; i < Cfs_->size(); i++)
349 {
350 const point& fc = Cfs_()[i];
351
352 const vector nf = Nfs_()[i];
353
354 const label myFaceId = shootFacesIds_()[i];
355
356 forAll(refDisDirsIndex, dirIndex)
357 {
358 if (refDisDirsIndex[dirIndex] > -1)
359 {
360 if ((nf & refDiscAngles_[dirIndex]) > 0)
361 {
362 const vector direction = -refDiscAngles_[dirIndex];
363
364 start.append(fc + 0.001*direction);
365
366 startIndex.append(myFaceId);
367 dirStartIndex.append(dirIndex);
368
369 end.append(fc + maxBounding*direction);
370 }
371 }
372 }
373 }
374
375 } while (returnReduceOr(i < Cfs_->size()));
376
377 List<pointIndexHit> hitInfo(startIndex.size());
378
379 surfacesMesh_->findLine(start, end, hitInfo);
380
381 // Query the local trigId on hit faces
382 labelList triangleIndex;
384 (
385 surfacesMesh_->localQueries
386 (
387 hitInfo,
388 triangleIndex
389 )
390 );
391 const mapDistribute& map = mapPtr();
392
395 forAll(patchr, patchi)
396 {
397 patchr.emplace_set(patchi, nBands);
398 patcha.emplace_set(patchi, nBands);
399 }
400
401 // Fill patchr
402 forAll(patchr, patchi)
403 {
404 const polyPatch& pp = patches[patchi];
405
406 if (!pp.coupled() && !isA<cyclicAMIPolyPatch>(pp))
407 {
408 for (label bandI = 0; bandI < nBands; bandI++)
409 {
410 patchr[patchi][bandI] =
411 boundaryRadiation.specReflectivity
412 (
413 patchi,
414 bandI,
415 new vectorField(patches[patchi].size(), sunDir)
416 );
417
418 patcha[patchi][bandI] =
419 boundaryRadiation.absorptivity
420 (
421 patchi,
422 bandI,
423 new vectorField(patches[patchi].size(), sunDir)
424 );
425 }
426 }
427 }
428
429 List<scalarField> r(nBands);
430 for (label bandI = 0; bandI < nBands; bandI++)
431 {
432 r[bandI].setSize(triangleIndex.size());
433 }
434 labelList refDirIndex(triangleIndex.size(), -1);
435 labelList refIndex(triangleIndex.size(), -1);
436 // triangleIndex includes hits on non-reflecting and reflecting faces
437 forAll(triangleIndex, i)
438 {
439 label trii = triangleIndex[i];
440 label facei = mapTriToGlobal_[trii];
441 label patchI = patches.whichPatch(facei);
442 const polyPatch& pp = patches[patchI];
443 label localFaceI = pp.whichFace(facei);
444
445 label globalFace = globalNumbering.toGlobal(Pstream::myProcNo(), facei);
446 if (refFacesDirIndex.found(globalFace))
447 {
448 refDirIndex[i] = refFacesDirIndex.find(globalFace)();
449 refIndex[i] = globalFace;
450 }
451 for (label bandI = 0; bandI < nBands; bandI++)
452 {
453 r[bandI][i] = patchr[patchI][bandI][localFaceI];
454 }
455 }
456 map.reverseDistribute(hitInfo.size(), refDirIndex);
457 map.reverseDistribute(hitInfo.size(), refIndex);
458 for (label bandI = 0; bandI < nBands; bandI++)
459 {
460 map.reverseDistribute(hitInfo.size(), r[bandI]);
461 }
462
463 for (label bandI = 0; bandI < nBands; bandI++)
464 {
466 qreflective_[bandI].boundaryFieldRef();
467 qrefBf = 0.0;
468 }
469
470 const vector qPrim(solarCalc_.directSolarRad()*solarCalc_.direction());
471
472 // Collect rays with a hit (hitting reflecting surfaces)
473 // and whose reflected direction are equal to the shot ray
474 forAll(hitInfo, rayI)
475 {
476 if (hitInfo[rayI].hit())
477 {
478 if
479 (
480 dirStartIndex[rayI] == refDirIndex[rayI]
481 && refFacesDirIndex.found(refIndex[rayI])
482 )
483 {
484 for (label bandI = 0; bandI < nBands; bandI++)
485 {
487 qreflective_[bandI].boundaryFieldRef();
488
489 label startFaceId = startIndex[rayI];
490 label startPatchI = patches.whichPatch(startFaceId);
491
492 const polyPatch& ppStart = patches[startPatchI];
493 label localStartFaceI = ppStart.whichFace(startFaceId);
494
495 scalar a = patcha[startPatchI][bandI][localStartFaceI];
496
497 const vectorField& nStart = ppStart.faceNormals();
498
499 vector rayIn = refDiscAngles_[dirStartIndex[rayI]];
500
501 rayIn /= mag(rayIn);
502
503 qrefBf[startPatchI][localStartFaceI] +=
504 (
505 (
506 mag(qPrim)
507 *r[bandI][rayI]
508 *spectralDistribution_[bandI]
509 *a
510 *rayIn
511 )
512 & nStart[localStartFaceI]
513 );
514 }
515 }
516 }
517 }
518
519 start.clear();
520 startIndex.clear();
521 end.clear();
522 dirStartIndex.clear();
523}
524
525
526// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
527
528Foam::faceReflecting::faceReflecting
529(
530 const fvMesh& mesh,
531 const faceShading& directHiyFaces,
532 const solarCalculator& solar,
533 const scalarList& spectralDistribution,
534 const dictionary& dict
535)
536:
537 mesh_(mesh),
538 nTheta_(dict.subDict("reflecting").getOrDefault<label>("nTheta", 10)),
539 nPhi_(dict.subDict("reflecting").getOrDefault<label>("nPhi", 10)),
540 nRay_(0),
541 refDiscAngles_(0),
542 spectralDistribution_(spectralDistribution),
543 qreflective_(spectralDistribution_.size()),
544 directHitFaces_(directHiyFaces),
545 surfacesMesh_(),
546 shootFacesIds_(),
547 Cfs_(),
548 Nfs_(),
549 solarCalc_(solar),
550 includePatches_(),
551 mapTriToGlobal_()
553 initialise(dict);
554}
555
556
557// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
558
560{
561 calculate();
562}
563
564
565// ************************************************************************* //
label n
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Dynamically sized Field. Similar to DynamicList, but inheriting from a Field instead of a List.
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
void append(const T &val)
Copy append an element to the end of this list.
SubField< vector > subField
Definition Field.H:183
GeometricBoundaryField< scalar, fvPatchField, volMesh > Boundary
bool found(const Key &key) const
Same as contains().
Definition HashTable.H:1370
bool insert(const Key &key, const T &obj)
Copy insert a new entry, not overwriting existing entries.
Definition HashTableI.H:152
iterator find(const Key &key)
Find and return an iterator set at the hashed entry.
Definition HashTableI.H:86
@ REGISTER
Request registration (bool: true).
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
void setSize(label n)
Alias for resize().
Definition List.H:536
A HashTable to objects of type <T> with a label key.
Definition Map.H:54
static FOAM_NO_DANGLING_REFERENCE const boundaryRadiationProperties & New(const fvMesh &mesh, Args &&... args)
const Field< point_type > & faceNormals() const
Return face unit normals for patch.
static void mapCombineReduce(Container &values, CombineOp cop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Forwards to Pstream::mapReduce with an in-place cop.
static void listReduce(UList< T > &values, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce list elements (list must be equal size on all ranks), applying bop to each list element.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
UList< T > & reset(std::nullptr_t) noexcept
Reset to zero-sized and nullptr.
Definition SubListI.H:109
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition Time.C:714
bool found(const T &val, label pos=0) const
Same as contains().
Definition UList.H:983
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Negative if the process is not a...
Definition UPstream.H:1706
label size() const noexcept
The number of entries in the list.
Definition UPtrListI.H:106
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
static const Enum< distributionType > distributionTypeNames_
Calculates the reflecting faces from specular surfaces. It only takes into account the first reflecti...
void correct()
Correct reflected flux.
Helper class to calculate visible faces for global, sun-like illumination.
Definition faceShading.H:58
A fvBoundaryMesh is a fvPatch list with a reference to the associated fvMesh, with additional search ...
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
const Time & time() const
Return the top-level database.
Definition fvMesh.H:360
Calculates a non-overlapping list of offsets based on an input size (eg, number of cells) from differ...
Definition globalIndex.H:77
Class containing processor-to-processor mapping information.
void reverseDistribute(const label constructSize, List< T > &fld, const bool dummyTransform=true, const int tag=UPstream::msgType()) const
Reverse distribute data using default commsType.
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh,...
label whichPatch(const label meshFacei) const
Return patch index for a given mesh face index. Uses binary search.
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
label whichFace(const label facei) const noexcept
Return label of face in patch from global face label.
Definition polyPatch.H:608
tmp< scalarField > transmissivity(const label patchI, const label bandI=0, const vectorField *incomingDirection=nullptr, const scalarField *T=nullptr) const
Access boundary transmissivity on patch.
tmp< scalarField > absorptivity(const label patchI, const label bandI=0, const vectorField *incomingDirection=nullptr, const scalarField *T=nullptr) const
Access boundary absorptivity on patch.
tmp< scalarField > specReflectivity(const label patchI, const label bandI=0, const vectorField *incomingDirection=nullptr, const scalarField *T=nullptr) const
Access boundary specular reflectivity on patch.
A solar calculator model providing models for the solar direction and solar loads.
A class for managing temporary objects.
Definition tmp.H:75
static triSurface triangulate(const polyBoundaryMesh &mBesh, const labelHashSet &includePatches, labelList &faceMap, const bool verbose=false)
Simple triangulation of (selected patches of) boundaryMesh. Needs.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
const polyBoundaryMesh & patches
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
const char * end
Definition SVGTools.H:223
Mathematical constants.
constexpr scalar pi(M_PI)
constexpr scalar piByTwo(0.5 *M_PI)
T cosPhi(const Vector< T > &a, const Vector< T > &b, const T &tolerance=SMALL)
Calculate angle between a and b in radians.
Namespace for OpenFOAM.
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
bool returnReduceOr(const bool value, const int communicator=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
List< label > labelList
A List of labels.
Definition List.H:62
dimensionedScalar pow3(const dimensionedScalar &ds)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
dimensionedScalar sin(const dimensionedScalar &ds)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
T returnReduce(const T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
errorManip< error > abort(error &err)
Definition errorManip.H:139
Field< vector > vectorField
Specialisation of Field<T> for vector.
uint8_t direction
Definition direction.H:49
vector point
Point is a vector.
Definition point.H:37
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
vectorField pointField
pointField is a vectorField.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Vector< scalar > vector
Definition vector.H:57
List< scalar > scalarList
List of scalar.
Definition scalarList.H:32
dimensionedScalar cos(const dimensionedScalar &ds)
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
dictionary dict
const triSurface localSurface
List< treeBoundBox > meshBb(1, treeBoundBox(coarseMesh.points()).extend(rndGen, 1e-3))
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Random rndGen