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_,
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  ),
238  localSurface,
239  dict
240  )
241  );
242 
243  if (debug)
244  {
245  surfacesMesh_->searchableSurface::write();
246  }
247 }
248 
249 
250 void 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::listCombineReduce(refDisDirsIndex, maxEqOp<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;
383  autoPtr<mapDistribute> mapPtr
384  (
385  surfacesMesh_->localQueries
386  (
387  hitInfo,
388  triangleIndex
389  )
390  );
391  const mapDistribute& map = mapPtr();
392 
393  PtrList<List<scalarField>> patchr(patches.size());
394  PtrList<List<scalarField>> patcha(patches.size());
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  {
465  volScalarField::Boundary& qrefBf =
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  {
486  volScalarField::Boundary& qrefBf =
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 
528 Foam::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_()
552 {
553  initialise(dict);
554 }
555 
556 
557 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
558 
560 {
561  calculate();
562 }
563 
564 
565 // ************************************************************************* //
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:608
static const boundaryRadiationProperties & New(const fvMesh &mesh, Args &&... args)
Get existing or create 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:1086
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:72
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.
Request registration (bool: true)
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)
Combines List elements. After completion all processors have the same data.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127