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-2022 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
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 
33 
34 using namespace Foam::constant::mathematical;
35 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40  defineTypeNameAndDebug(faceReflecting, 0);
41 }
42 
43 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
44 
45 void Foam::faceReflecting::initialise(const dictionary& coeffs)
46 {
47 
48  forAll(qreflective_, bandI)
49  {
50  qreflective_.set
51  (
52  bandI,
53  new volScalarField
54  (
55  IOobject
56  (
57  "qreflective_" + Foam::name(bandI) ,
58  mesh_.time().timeName(),
59  mesh_,
62  ),
63  mesh_,
65  )
66  );
67  }
68 
69  label rayI = 0;
70  if (mesh_.nSolutionD() == 3)
71  {
72  nRay_ = 4*nPhi_*nTheta_;
73  refDiscAngles_.resize(nRay_);
74  const scalar deltaPhi = pi/(2.0*nPhi_);
75  const scalar deltaTheta = pi/nTheta_;
76 
77  for (label n = 1; n <= nTheta_; n++)
78  {
79  for (label m = 1; m <= 4*nPhi_; m++)
80  {
81  const scalar thetai = (2*n - 1)*deltaTheta/2.0;
82  const scalar phii = (2*m - 1)*deltaPhi/2.0;
83 
84  scalar sinTheta = Foam::sin(thetai);
85  scalar cosTheta = Foam::cos(thetai);
86  scalar sinPhi = Foam::sin(phii);
87  scalar cosPhi = Foam::cos(phii);
88  refDiscAngles_[rayI++] =
89  vector(sinTheta*sinPhi, sinTheta*cosPhi, cosTheta);
90 
91  }
92  }
93 
94  }
95  else if (mesh_.nSolutionD() == 2)
96  {
97  nRay_ = 4*nPhi_;
98  refDiscAngles_.resize(nRay_);
99  const scalar thetai = piByTwo;
100  //const scalar deltaTheta = pi;
101  const scalar deltaPhi = pi/(2.0*nPhi_);
102  for (label m = 1; m <= 4*nPhi_; m++)
103  {
104  const scalar phii = (2*m - 1)*deltaPhi/2.0;
105 
106  scalar sinTheta = Foam::sin(thetai);
107  scalar cosTheta = Foam::cos(thetai);
108  scalar sinPhi = Foam::sin(phii);
109  scalar cosPhi = Foam::cos(phii);
110 
111  refDiscAngles_[rayI++] =
112  vector(sinTheta*sinPhi, sinTheta*cosPhi, cosTheta);
113  }
114  }
115  else
116  {
118  << "The reflected rays are available in 2D or 3D "
119  << abort(FatalError);
120  }
121 
122  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
123 
124  const radiation::boundaryRadiationProperties& boundaryRadiation =
126 
127  // global face index
128  globalIndex globalNumbering(mesh_.nFaces());
129 
130 
131  // Collect faces with t = 0, r = 0 and a > 0 to shoot rays
132  // and patches to construct the triSurface
133  DynamicList<point> dynCf;
134  DynamicList<vector> dynNf;
135  DynamicList<label> dynFacesI;
136  forAll(patches, patchI)
137  {
138  const polyPatch& pp = patches[patchI];
139  const vectorField::subField cf = pp.faceCentres();
140 
141  if (!pp.coupled() && !isA<cyclicAMIPolyPatch>(pp))
142  {
143  const tmp<scalarField> tt =
144  boundaryRadiation.transmissivity(patchI);
145 
146  const tmp<scalarField> tr =
147  boundaryRadiation.specReflectivity(patchI);
148 
149  const tmp<scalarField> ta =
150  boundaryRadiation.absorptivity(patchI);
151 
152  const scalarField& t = tt();
153  const scalarField& r = tr();
154  const scalarField& a = ta();
155 
156  const vectorField& n = pp.faceNormals();
157 
158  forAll(pp, faceI)
159  {
160  //const vector nf(n[faceI]);
161  // Opaque, non-reflective, absortived faces to shoot
162  if
163  (
164  t[faceI] == 0
165  && r[faceI] == 0
166  && a[faceI] > 0
167  )
168  {
169  dynFacesI.append(faceI + pp.start());
170  dynCf.append(cf[faceI]);
171  dynNf.append(n[faceI]);
172  }
173 
174  // relfective opaque patches to build reflective surface
175  // plus opaque non-reflective
176  if
177  (
178  (r[faceI] > 0 && t[faceI] == 0) ||
179  (t[faceI] == 0 && a[faceI] > 0 && r[faceI] == 0)
180  )
181  {
182  includePatches_.insert(patchI);
183  }
184  }
185  }
186  }
187 
188  shootFacesIds_.reset(new labelList(dynFacesI));
189  Cfs_.reset(new pointField(dynCf));
190  Nfs_.reset(new vectorField(dynNf));
191 
192  // * * * * * * * * * * * * * * *
193  // Create distributedTriSurfaceMesh
194  Random rndGen(653213);
195 
196  // Determine mesh bounding boxes:
197  List<treeBoundBox> meshBb
198  (
199  1,
200  treeBoundBox(mesh_.points()).extend(rndGen, 1e-3)
201  );
202 
203  // Dummy bounds dictionary
205  dict.add("bounds", meshBb);
206  dict.add
207  (
208  "distributionType",
210  [
212  ]
213  );
214  dict.add("mergeDistance", SMALL);
215 
216 
218  (
219  mesh_.boundaryMesh(),
220  includePatches_,
221  mapTriToGlobal_
222  );
223 
224  surfacesMesh_.reset
225  (
226  new distributedTriSurfaceMesh
227  (
228  IOobject
229  (
230  "reflectiveSurface.stl",
231  mesh_.time().constant(),
232  "triSurface",
233  mesh_.time(),
236  ),
237  localSurface,
238  dict
239  )
240  );
241 
242  if (debug)
243  {
244  surfacesMesh_->searchableSurface::write();
245  }
246 }
247 
248 
249 void Foam::faceReflecting::calculate()
250 {
251  const radiation::boundaryRadiationProperties& boundaryRadiation =
253 
254  label nFaces = 0;
255 
256  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
257 
258  const fvBoundaryMesh& fvPatches = mesh_.boundary();
259 
260  label nBands = spectralDistribution_.size();
261 
262  // Collect reflected directions from reflecting surfaces on direct hit
263  // faces
264  const vector sunDir = directHitFaces_.direction();
265  const labelList& directHits = directHitFaces_.rayStartFaces();
266 
267  globalIndex globalNumbering(mesh_.nFaces());
268 
269  Map<label> refFacesDirIndex;
270  labelList refDisDirsIndex(nRay_, -1);
271 
272  forAll(patches, patchI)
273  {
274  const polyPatch& pp = patches[patchI];
275 
276  if (!pp.coupled() && !isA<cyclicAMIPolyPatch>(pp))
277  {
278  const tmp<scalarField> tr =
279  boundaryRadiation.specReflectivity(patchI);
280 
281  const scalarField& r = tr();
282  const vectorField n(fvPatches[patchI].nf());
283 
284  forAll(pp, faceI)
285  {
286  label globalID = faceI + pp.start();
287 
288  if (r[faceI] > 0.0 && directHits.found(globalID))
289  {
290  vector refDir =
291  sunDir + 2.0*(-sunDir & n[faceI]) * n[faceI];
292 
293  // Look for the discrete direction
294  scalar dev(-GREAT);
295  label rayIndex = -1;
296  forAll(refDiscAngles_, iDisc)
297  {
298  scalar dotProd = refDir & refDiscAngles_[iDisc];
299  if (dev < dotProd)
300  {
301  dev = dotProd;
302  rayIndex = iDisc;
303  }
304  }
305 
306  if (rayIndex >= 0)
307  {
308  if (refDisDirsIndex[rayIndex] == -1)
309  {
310  refDisDirsIndex[rayIndex] = 1;
311  }
312  }
313 
314  refFacesDirIndex.insert
315  (
316  globalNumbering.toGlobal(globalID),
317  rayIndex
318  );
319 
320  nFaces++;
321  }
322  }
323  }
324  }
325 
326  // Distribute ray indexes to all proc's
327  // Make sure all the processors have the same information
328 
329  Pstream::listCombineReduce(refDisDirsIndex, maxEqOp<label>());
330  Pstream::mapCombineReduce(refFacesDirIndex, minEqOp<label>());
331 
332  const scalar maxBounding =
333  returnReduce(5.0*mesh_.bounds().mag(), maxOp<scalar>());
334 
335  // Shoot Rays
336  // From faces t = 0, r = 0 and a > 0 to all 'used' discrete reflected
337  // directions
338 
339  DynamicField<point> start(nFaces);
340  DynamicField<point> end(start.size());
341  DynamicList<label> startIndex(start.size());
342  DynamicField<label> dirStartIndex(start.size());
343 
344  label i = 0;
345  do
346  {
347  for (; i < Cfs_->size(); i++)
348  {
349  const point& fc = Cfs_()[i];
350 
351  const vector nf = Nfs_()[i];
352 
353  const label myFaceId = shootFacesIds_()[i];
354 
355  forAll(refDisDirsIndex, dirIndex)
356  {
357  if (refDisDirsIndex[dirIndex] > -1)
358  {
359  if ((nf & refDiscAngles_[dirIndex]) > 0)
360  {
361  const vector direction = -refDiscAngles_[dirIndex];
362 
363  start.append(fc + 0.001*direction);
364 
365  startIndex.append(myFaceId);
366  dirStartIndex.append(dirIndex);
367 
368  end.append(fc + maxBounding*direction);
369  }
370  }
371  }
372  }
373 
374  } while (returnReduceOr(i < Cfs_->size()));
375 
376  List<pointIndexHit> hitInfo(startIndex.size());
377 
378  surfacesMesh_->findLine(start, end, hitInfo);
379 
380  // Query the local trigId on hit faces
381  labelList triangleIndex;
382  autoPtr<mapDistribute> mapPtr
383  (
384  surfacesMesh_->localQueries
385  (
386  hitInfo,
387  triangleIndex
388  )
389  );
390  const mapDistribute& map = mapPtr();
391 
392  PtrList<List<scalarField>> patchr(patches.size());
393  PtrList<List<scalarField>> patcha(patches.size());
394  forAll(patchr, patchi)
395  {
396  patchr.set
397  (
398  patchi,
399  new List<scalarField>(nBands)
400  );
401 
402  patcha.set
403  (
404  patchi,
405  new List<scalarField>(nBands)
406  );
407  }
408 
409  // Fill patchr
410  forAll(patchr, patchi)
411  {
412  const polyPatch& pp = patches[patchi];
413 
414  if (!pp.coupled() && !isA<cyclicAMIPolyPatch>(pp))
415  {
416  for (label bandI = 0; bandI < nBands; bandI++)
417  {
418  patchr[patchi][bandI] =
419  boundaryRadiation.specReflectivity
420  (
421  patchi,
422  bandI,
423  new vectorField(patches[patchi].size(), sunDir)
424  );
425 
426  patcha[patchi][bandI] =
427  boundaryRadiation.absorptivity
428  (
429  patchi,
430  bandI,
431  new vectorField(patches[patchi].size(), sunDir)
432  );
433  }
434  }
435  }
436 
437  List<scalarField> r(nBands);
438  for (label bandI = 0; bandI < nBands; bandI++)
439  {
440  r[bandI].setSize(triangleIndex.size());
441  }
442  labelList refDirIndex(triangleIndex.size(), -1);
443  labelList refIndex(triangleIndex.size(), -1);
444  // triangleIndex includes hits on non-reflecting and reflecting faces
445  forAll(triangleIndex, i)
446  {
447  label trii = triangleIndex[i];
448  label facei = mapTriToGlobal_[trii];
449  label patchI = patches.whichPatch(facei);
450  const polyPatch& pp = patches[patchI];
451  label localFaceI = pp.whichFace(facei);
452 
453  label globalFace = globalNumbering.toGlobal(Pstream::myProcNo(), facei);
454  if (refFacesDirIndex.found(globalFace))
455  {
456  refDirIndex[i] = refFacesDirIndex.find(globalFace)();
457  refIndex[i] = globalFace;
458  }
459  for (label bandI = 0; bandI < nBands; bandI++)
460  {
461  r[bandI][i] = patchr[patchI][bandI][localFaceI];
462  }
463  }
464  map.reverseDistribute(hitInfo.size(), refDirIndex);
465  map.reverseDistribute(hitInfo.size(), refIndex);
466  for (label bandI = 0; bandI < nBands; bandI++)
467  {
468  map.reverseDistribute(hitInfo.size(), r[bandI]);
469  }
470 
471  for (label bandI = 0; bandI < nBands; bandI++)
472  {
473  volScalarField::Boundary& qrefBf =
474  qreflective_[bandI].boundaryFieldRef();
475  qrefBf = 0.0;
476  }
477 
478  const vector qPrim(solarCalc_.directSolarRad()*solarCalc_.direction());
479 
480  // Collect rays with a hit (hitting reflecting surfaces)
481  // and whose reflected direction are equal to the shot ray
482  forAll(hitInfo, rayI)
483  {
484  if (hitInfo[rayI].hit())
485  {
486  if
487  (
488  dirStartIndex[rayI] == refDirIndex[rayI]
489  && refFacesDirIndex.found(refIndex[rayI])
490  )
491  {
492  for (label bandI = 0; bandI < nBands; bandI++)
493  {
494  volScalarField::Boundary& qrefBf =
495  qreflective_[bandI].boundaryFieldRef();
496 
497  label startFaceId = startIndex[rayI];
498  label startPatchI = patches.whichPatch(startFaceId);
499 
500  const polyPatch& ppStart = patches[startPatchI];
501  label localStartFaceI = ppStart.whichFace(startFaceId);
502 
503  scalar a = patcha[startPatchI][bandI][localStartFaceI];
504 
505  const vectorField& nStart = ppStart.faceNormals();
506 
507  vector rayIn = refDiscAngles_[dirStartIndex[rayI]];
508 
509  rayIn /= mag(rayIn);
510 
511  qrefBf[startPatchI][localStartFaceI] +=
512  (
513  (
514  mag(qPrim)
515  *r[bandI][rayI]
516  *spectralDistribution_[bandI]
517  *a
518  *rayIn
519  )
520  & nStart[localStartFaceI]
521  );
522  }
523  }
524  }
525  }
526 
527  start.clear();
528  startIndex.clear();
529  end.clear();
530  dirStartIndex.clear();
531 }
532 
533 
534 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
535 
536 Foam::faceReflecting::faceReflecting
537 (
538  const fvMesh& mesh,
539  const faceShading& directHiyFaces,
540  const solarCalculator& solar,
541  const scalarList& spectralDistribution,
542  const dictionary& dict
543 )
544 :
545  mesh_(mesh),
546  nTheta_(dict.subDict("reflecting").getOrDefault<label>("nTheta", 10)),
547  nPhi_(dict.subDict("reflecting").getOrDefault<label>("nPhi", 10)),
548  nRay_(0),
549  refDiscAngles_(0),
550  spectralDistribution_(spectralDistribution),
551  qreflective_(spectralDistribution_.size()),
552  directHitFaces_(directHiyFaces),
553  surfacesMesh_(),
554  shootFacesIds_(),
555  Cfs_(),
556  Nfs_(),
557  solarCalc_(solar),
558  includePatches_(),
559  mapTriToGlobal_()
560 {
561  initialise(dict);
562 }
563 
564 
565 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
566 
568 {
569  calculate();
570 }
571 
572 
573 // ************************************************************************* //
static const Enum< distributionType > distributionTypeNames_
dictionary dict
const triSurface localSurface
A solar calculator model providing models for the solar direction and solar loads.
uint8_t direction
Definition: direction.H:46
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
static const boundaryRadiationProperties & New(const fvMesh &mesh, Args &&... args)
Get existing or create a new MeshObject. Registered with typeName.
Definition: MeshObject.C:53
static triSurface triangulate(const polyBoundaryMesh &mBesh, const labelHashSet &includePatches, labelList &faceMap, const bool verbose=false)
Simple triangulation of (selected patches of) boundaryMesh. Needs.
Random rndGen
Definition: createFields.H:23
Ignore writing from objectRegistry::writeObject()
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Can be negative if the process i...
Definition: UPstream.H:1074
SubField is a Field obtained as a section of another Field, without its own allocation. SubField is derived from a SubList rather than a List.
Definition: Field.H:63
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
T cosPhi(const Vector< T > &a, const Vector< T > &b, const T &tolerance=SMALL)
Calculate angle between a and b in radians.
Definition: vectorTools.H:115
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
dynamicFvMesh & mesh
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
dimensionedScalar cos(const dimensionedScalar &ds)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
Mathematical constants.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
constexpr scalar pi(M_PI)
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
Vector< scalar > vector
Definition: vector.H:57
constexpr scalar piByTwo(0.5 *M_PI)
errorManip< error > abort(error &err)
Definition: errorManip.H:139
Helper class to calculate visible faces for global, sun-like illumination.
Definition: faceShading.H:57
dimensionedScalar sin(const dimensionedScalar &ds)
int debug
Static debugging option.
constexpr auto end(C &c) -> decltype(c.end())
Return iterator to the end of the container c.
Definition: stdFoam.H:201
defineTypeNameAndDebug(combustionModel, 0)
void correct()
Correct reflected flux.
dimensionedScalar pow3(const dimensionedScalar &ds)
vector point
Point is a vector.
Definition: point.H:37
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
const polyBoundaryMesh & patches
Nothing to be read.
Automatically write from objectRegistry::writeObject()
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
void clear()
Clear the PtrList. Delete allocated entries and set size to zero.
Definition: PtrListI.H:81
label n
Field< vector > vectorField
Specialisation of Field<T> for vector.
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
List< label > labelList
A List of labels.
Definition: List.H:62
constexpr scalar e(M_E)
static void mapCombineReduce(Container &values, const CombineOp &cop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce) applying cop to inplace combine map values from different processo...
bool returnReduceOr(const bool value, const label comm=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
List< treeBoundBox > meshBb(1, treeBoundBox(coarseMesh.points()).extend(rndGen, 1e-3))
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Namespace for OpenFOAM.
static void listCombineReduce(List< T > &values, const CombineOp &cop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
After completion all processors have the same data.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127