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-2020 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 pointField& 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
201  (
202  boundBox(mesh_.points(), false)
203  ).extend(rndGen, 1e-3)
204  );
205 
206  // Dummy bounds dictionary
207  dictionary dict;
208  dict.add("bounds", meshBb);
209  dict.add
210  (
211  "distributionType",
213  [
215  ]
216  );
217  dict.add("mergeDistance", SMALL);
218 
219 
221  (
222  mesh_.boundaryMesh(),
223  includePatches_,
224  mapTriToGlobal_
225  );
226 
227  surfacesMesh_.reset
228  (
229  new distributedTriSurfaceMesh
230  (
231  IOobject
232  (
233  "reflectiveSurface.stl",
234  mesh_.time().constant(),
235  "triSurface",
236  mesh_.time(),
239  ),
240  localSurface,
241  dict
242  )
243  );
244 
245  if (debug)
246  {
247  surfacesMesh_->searchableSurface::write();
248  }
249 }
250 
251 
252 void Foam::faceReflecting::calculate()
253 {
254  const radiation::boundaryRadiationProperties& boundaryRadiation =
256 
257  label nFaces = 0;
258 
259  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
260 
261  const fvBoundaryMesh& fvPatches = mesh_.boundary();
262 
263  label nBands = spectralDistribution_.size();
264 
265  // Collect reflected directions from reflecting surfaces on direct hit
266  // faces
267  const vector sunDir = directHitFaces_.direction();
268  const labelList& directHits = directHitFaces_.rayStartFaces();
269 
270  globalIndex globalNumbering(mesh_.nFaces());
271 
272  Map<label> refFacesDirIndex;
273  labelList refDisDirsIndex(nRay_, -1);
274 
275  forAll(patches, patchI)
276  {
277  const polyPatch& pp = patches[patchI];
278 
279  if (!pp.coupled() && !isA<cyclicAMIPolyPatch>(pp))
280  {
281  const tmp<scalarField> tr =
282  boundaryRadiation.specReflectivity(patchI);
283 
284  const scalarField& r = tr();
285  const vectorField n(fvPatches[patchI].nf());
286 
287  forAll(pp, faceI)
288  {
289  label globalID = faceI + pp.start();
290 
291  if (r[faceI] > 0.0 && directHits.found(globalID))
292  {
293  vector refDir =
294  sunDir + 2.0*(-sunDir & n[faceI]) * n[faceI];
295 
296  // Look for the discrete direction
297  scalar dev(-GREAT);
298  label rayIndex = -1;
299  forAll(refDiscAngles_, iDisc)
300  {
301  scalar dotProd = refDir & refDiscAngles_[iDisc];
302  if (dev < dotProd)
303  {
304  dev = dotProd;
305  rayIndex = iDisc;
306  }
307  }
308 
309  if (rayIndex >= 0)
310  {
311  if (refDisDirsIndex[rayIndex] == -1)
312  {
313  refDisDirsIndex[rayIndex] = 1;
314  }
315  }
316 
317  refFacesDirIndex.insert
318  (
319  globalNumbering.toGlobal(globalID),
320  rayIndex
321  );
322 
323  nFaces++;
324  }
325  }
326  }
327  }
328 
329  // Distribute ray indexes to all proc's
330  Pstream::listCombineGather(refDisDirsIndex, maxEqOp<label>());
331  Pstream::listCombineScatter(refDisDirsIndex);
332 
333  // Make sure all the processors have the same map
334  Pstream::mapCombineGather(refFacesDirIndex, minEqOp<label>());
335  Pstream::mapCombineScatter(refFacesDirIndex);
336 
337  scalar maxBounding = 5.0*mag(mesh_.bounds().max() - mesh_.bounds().min());
338 
339  reduce(maxBounding, maxOp<scalar>());
340 
341  // Shoot Rays
342  // From faces t = 0, r = 0 and a > 0 to all 'used' discrete reflected
343  // directions
344 
345  DynamicField<point> start(nFaces);
346  DynamicField<point> end(start.size());
347  DynamicList<label> startIndex(start.size());
348  DynamicField<label> dirStartIndex(start.size());
349 
350  label i = 0;
351  do
352  {
353  for (; i < Cfs_->size(); i++)
354  {
355  const point& fc = Cfs_()[i];
356 
357  const vector nf = Nfs_()[i];
358 
359  const label myFaceId = shootFacesIds_()[i];
360 
361  forAll(refDisDirsIndex, dirIndex)
362  {
363  if (refDisDirsIndex[dirIndex] > -1)
364  {
365  if ((nf & refDiscAngles_[dirIndex]) > 0)
366  {
367  const vector direction = -refDiscAngles_[dirIndex];
368 
369  start.append(fc + 0.001*direction);
370 
371  startIndex.append(myFaceId);
372  dirStartIndex.append(dirIndex);
373 
374  end.append(fc + maxBounding*direction);
375  }
376  }
377  }
378  }
379 
380  }while (returnReduce(i < Cfs_->size(), orOp<bool>()));
381 
382  List<pointIndexHit> hitInfo(startIndex.size());
383 
384  surfacesMesh_->findLine(start, end, hitInfo);
385 
386  // Query the local trigId on hit faces
387  labelList triangleIndex;
388  autoPtr<mapDistribute> mapPtr
389  (
390  surfacesMesh_->localQueries
391  (
392  hitInfo,
393  triangleIndex
394  )
395  );
396  const mapDistribute& map = mapPtr();
397 
398  PtrList<List<scalarField>> patchr(patches.size());
399  PtrList<List<scalarField>> patcha(patches.size());
400  forAll(patchr, patchi)
401  {
402  patchr.set
403  (
404  patchi,
405  new List<scalarField>(nBands)
406  );
407 
408  patcha.set
409  (
410  patchi,
411  new List<scalarField>(nBands)
412  );
413  }
414 
415  // Fill patchr
416  forAll(patchr, patchi)
417  {
418  const polyPatch& pp = patches[patchi];
419 
420  if (!pp.coupled() && !isA<cyclicAMIPolyPatch>(pp))
421  {
422  for (label bandI = 0; bandI < nBands; bandI++)
423  {
424  patchr[patchi][bandI] =
425  boundaryRadiation.specReflectivity
426  (
427  patchi,
428  bandI,
429  new vectorField(patches[patchi].size(), sunDir)
430  );
431 
432  patcha[patchi][bandI] =
433  boundaryRadiation.absorptivity
434  (
435  patchi,
436  bandI,
437  new vectorField(patches[patchi].size(), sunDir)
438  );
439  }
440  }
441  }
442 
443  List<scalarField> r(nBands);
444  for (label bandI = 0; bandI < nBands; bandI++)
445  {
446  r[bandI].setSize(triangleIndex.size());
447  }
448  labelList refDirIndex(triangleIndex.size());
449  labelList refIndex(triangleIndex.size());
450  // triangleIndex includes hits on non-reflecting and reflecting faces
451  forAll(triangleIndex, i)
452  {
453  label trii = triangleIndex[i];
454  label facei = mapTriToGlobal_[trii];
455  label patchI = patches.whichPatch(facei);
456  const polyPatch& pp = patches[patchI];
457  label localFaceI = pp.whichFace(facei);
458 
459  label globalFace = globalNumbering.toGlobal(Pstream::myProcNo(), facei);
460  if (refFacesDirIndex.found(globalFace))
461  {
462  refDirIndex[i] = refFacesDirIndex.find(globalFace)();
463  refIndex[i] = globalFace;
464  }
465  for (label bandI = 0; bandI < nBands; bandI++)
466  {
467  r[bandI][i] = patchr[patchI][bandI][localFaceI];
468  }
469  }
470  map.reverseDistribute(hitInfo.size(), refDirIndex);
471  map.reverseDistribute(hitInfo.size(), refIndex);
472  for (label bandI = 0; bandI < nBands; bandI++)
473  {
474  map.reverseDistribute(hitInfo.size(), r[bandI]);
475  }
476 
477  for (label bandI = 0; bandI < nBands; bandI++)
478  {
479  volScalarField::Boundary& qrefBf =
480  qreflective_[bandI].boundaryFieldRef();
481  qrefBf = 0.0;
482  }
483 
484  const vector qPrim(solarCalc_.directSolarRad()*solarCalc_.direction());
485 
486  // Collect rays with a hit (hitting reflecting surfaces)
487  // and whose reflected direction are equal to the shot ray
488  forAll(hitInfo, rayI)
489  {
490  if (hitInfo[rayI].hit())
491  {
492  if
493  (
494  dirStartIndex[rayI] == refDirIndex[rayI]
495  && refFacesDirIndex.found(refIndex[rayI])
496  )
497  {
498  for (label bandI = 0; bandI < nBands; bandI++)
499  {
500  volScalarField::Boundary& qrefBf =
501  qreflective_[bandI].boundaryFieldRef();
502 
503  label startFaceId = startIndex[rayI];
504  label startPatchI = patches.whichPatch(startFaceId);
505 
506  const polyPatch& ppStart = patches[startPatchI];
507  label localStartFaceI = ppStart.whichFace(startFaceId);
508 
509  scalar a = patcha[startPatchI][bandI][localStartFaceI];
510 
511  const vectorField& nStart = ppStart.faceNormals();
512 
513  vector rayIn = refDiscAngles_[dirStartIndex[rayI]];
514 
515  rayIn /= mag(rayIn);
516 
517  qrefBf[startPatchI][localStartFaceI] +=
518  (
519  (
520  mag(qPrim)
521  *r[bandI][rayI]
522  *spectralDistribution_[bandI]
523  *a
524  *rayIn
525  )
526  & nStart[localStartFaceI]
527  );
528  }
529  }
530  }
531  }
532 
533  start.clear();
534  startIndex.clear();
535  end.clear();
536  dirStartIndex.clear();
537 }
538 
539 
540 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
541 
542 Foam::faceReflecting::faceReflecting
543 (
544  const fvMesh& mesh,
545  const faceShading& directHiyFaces,
546  const solarCalculator& solar,
547  const scalarList& spectralDistribution,
548  const dictionary& dict
549 )
550 :
551  mesh_(mesh),
552  nTheta_(dict.subDict("reflecting").getOrDefault<label>("nTheta", 10)),
553  nPhi_(dict.subDict("reflecting").getOrDefault<label>("nPhi", 10)),
554  nRay_(0),
555  refDiscAngles_(0),
556  spectralDistribution_(spectralDistribution),
557  qreflective_(spectralDistribution_.size()),
558  directHitFaces_(directHiyFaces),
559  surfacesMesh_(),
560  shootFacesIds_(),
561  Cfs_(),
562  Nfs_(),
563  solarCalc_(solar),
564  includePatches_(),
565  mapTriToGlobal_()
566 {
567  initialise(dict);
568 }
569 
570 
571 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
572 
574 {
575  calculate();
576 }
577 
578 
579 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::Pstream::mapCombineGather
static void mapCombineGather(const List< commsStruct > &comms, Container &Values, const CombineOp &cop, const int tag, const label comm)
Definition: combineGatherScatter.C:551
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
volFields.H
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
localSurface
const triSurface localSurface
Definition: searchingEngine.H:28
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::IOobject::AUTO_WRITE
Definition: IOobject.H:194
Foam::distributedTriSurfaceMesh::distributionTypeNames_
static const Enum< distributionType > distributionTypeNames_
Definition: distributedTriSurfaceMesh.H:101
faceReflecting.H
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:94
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::sin
dimensionedScalar sin(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:264
Foam::MeshObject< fvMesh, Foam::GeometricMeshObject, boundaryRadiationProperties >::New
static const boundaryRadiationProperties & New(const fvMesh &mesh, Args &&... args)
Get existing or create a new MeshObject.
Definition: MeshObject.C:48
Foam::vectorTools::cosPhi
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:107
Foam::constant::mathematical::e
constexpr scalar e(M_E)
meshBb
List< treeBoundBox > meshBb(1, treeBoundBox(boundBox(coarseMesh.points(), false)).extend(rndGen, 1e-3))
Foam::triSurfaceTools::triangulate
static triSurface triangulate(const polyBoundaryMesh &mBesh, const labelHashSet &includePatches, labelList &faceMap, const bool verbose=false)
Simple triangulation of (selected patches of) boundaryMesh. Needs.
Definition: triSurfaceTools.C:2192
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::distributedTriSurfaceMesh::FROZEN
Definition: distributedTriSurfaceMesh.H:98
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:54
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::solarCalculator
A solar calculator model providing models for the solar direction and solar loads.
Definition: solarCalculator.H:444
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:51
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:89
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
Foam::faceShading
faceShading uses the transmissivity value in the boundaryRadiationProperties in order to evaluate whi...
Definition: faceShading.H:61
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
stdFoam::end
constexpr auto end(C &c) -> decltype(c.end())
Return iterator to the end of the container c.
Definition: stdFoam.H:121
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::dimMass
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:51
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
cyclicAMIPolyPatch.H
Foam::Pstream::listCombineGather
static void listCombineGather(const List< commsStruct > &comms, List< T > &Value, const CombineOp &cop, const int tag, const label comm)
Definition: combineGatherScatter.C:290
Foam::constant::mathematical::piByTwo
constexpr scalar piByTwo(0.5 *M_PI)
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::constant::mathematical
Mathematical constants.
Foam::UPstream::myProcNo
static int myProcNo(const label communicator=worldComm)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:463
Foam::constant::mathematical::pi
constexpr scalar pi(M_PI)
Foam::Pstream::mapCombineScatter
static void mapCombineScatter(const List< commsStruct > &comms, Container &Values, const int tag, const label comm)
Scatter data. Reverse of combineGather.
Definition: combineGatherScatter.C:666
Foam::List< scalar >
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::direction
uint8_t direction
Definition: direction.H:52
patches
const polyBoundaryMesh & patches
Definition: convertProcessorPatches.H:65
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
rndGen
Random rndGen
Definition: createFields.H:23
Foam::tr
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
Definition: dimensionedSphericalTensor.C:51
Foam::Pstream::listCombineScatter
static void listCombineScatter(const List< commsStruct > &comms, List< T > &Value, const int tag, const label comm)
Scatter data. Reverse of combineGather.
Definition: combineGatherScatter.C:432
Foam::point
vector point
Point is a vector.
Definition: point.H:43
Foam::IOobject::NO_READ
Definition: IOobject.H:188
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
boundaryRadiationProperties.H
Foam::faceReflecting::correct
void correct()
Correct reflected flux.
Definition: faceReflecting.C:573
Foam::dev
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:106
Foam::cos
dimensionedScalar cos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:265