faceShading.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) 2015 OpenFOAM Foundation
9  Copyright (C) 2017-2019 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "faceShading.H"
30 #include "fvMesh.H"
32 #include "cyclicAMIPolyPatch.H"
33 #include "volFields.H"
35 #include "OBJstream.H"
36 
37 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
41  defineTypeNameAndDebug(faceShading, 0);
42 }
43 
44 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
45 
46 void Foam::faceShading::writeRays
47 (
48  const fileName& fName,
49  const DynamicField<point>& endCf,
50  const pointField& myFc
51 )
52 {
53  OBJstream os(fName);
54 
55  Pout<< "Dumping rays to " << os.name() << endl;
56 
57  forAll(myFc, faceI)
58  {
59  os.write(linePointRef(myFc[faceI], endCf[faceI]));
60  }
61 }
62 
63 
64 Foam::triSurface Foam::faceShading::triangulate
65 (
66  const labelHashSet& includePatches,
67  const List<labelHashSet>& includeAllFacesPerPatch
68 )
69 {
70  const polyBoundaryMesh& bMesh = mesh_.boundaryMesh();
71 
72  // Storage for surfaceMesh. Size estimate.
73  DynamicList<labelledTri> triangles(mesh_.nBoundaryFaces());
74 
75  label newPatchI = 0;
76 
77  for (const label patchI : includePatches)
78  {
79  const polyPatch& patch = bMesh[patchI];
80  const pointField& points = patch.points();
81 
82  label nTriTotal = 0;
83 
84  if (includeAllFacesPerPatch[patchI].size())
85  {
86  for (const label patchFaceI : includeAllFacesPerPatch[patchI])
87  {
88  const face& f = patch[patchFaceI];
89 
90  faceList triFaces(f.nTriangles(points));
91 
92  label nTri = 0;
93 
94  f.triangles(points, nTri, triFaces);
95 
96  for (const face& f : triFaces)
97  {
98  triangles.append
99  (
100  labelledTri(f[0], f[1], f[2], newPatchI)
101  );
102  nTriTotal++;
103  }
104  }
105  newPatchI++;
106  }
107  }
108 
109  triangles.shrink();
110 
111  // Create globally numbered tri surface
112  triSurface rawSurface(triangles, mesh_.points());
113 
114  // Create locally numbered tri surface
115  triSurface surface
116  (
117  rawSurface.localFaces(),
118  rawSurface.localPoints()
119  );
120 
121  // Add patch names to surface
122  surface.patches().setSize(newPatchI);
123 
124  newPatchI = 0;
125 
126  for (const label patchI : includePatches)
127  {
128  const polyPatch& patch = bMesh[patchI];
129 
130  if (includeAllFacesPerPatch[patchI].size())
131  {
132  surface.patches()[newPatchI].name() = patch.name();
133  surface.patches()[newPatchI].geometricType() = patch.type();
134 
135  newPatchI++;
136  }
137  }
138 
139  return surface;
140 }
141 
142 
143 void Foam::faceShading::calculate()
144 {
145  const radiation::boundaryRadiationProperties& boundaryRadiation =
147 
148  label nFaces = 0; //total number of direct hit faces
149 
150  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
151 
152  DynamicList<point> dynCf(nFaces);
153  DynamicList<label> dynFacesI;
154 
155  forAll(patches, patchI)
156  {
157  const polyPatch& pp = patches[patchI];
158  const pointField& cf = pp.faceCentres();
159 
160  if (!pp.coupled() && !isA<cyclicAMIPolyPatch>(pp))
161  {
162  const tmp<scalarField> tt =
163  boundaryRadiation.transmissivity(patchI);
164  const scalarField& t = tt();
165  const vectorField& n = pp.faceNormals();
166 
167  forAll(n, faceI)
168  {
169  const vector nf(n[faceI]);
170  if (((direction_ & nf) > 0) && (t[faceI] == 0.0))
171  {
172  dynFacesI.append(faceI + pp.start());
173  dynCf.append(cf[faceI]);
174  nFaces++;
175  }
176  }
177  }
178  }
179 
180  label numberPotentialHits = nFaces;
181 
182  reduce(numberPotentialHits, sumOp<label>());
183 
184  Info<< "Number of 'potential' direct hits : "
185  << numberPotentialHits << endl;
186 
187  labelList hitFacesIds(nFaces);
188  hitFacesIds.transfer(dynFacesI);
189 
190  pointField Cfs(hitFacesIds.size());
191  Cfs.transfer(dynCf);
192 
193  // * * * * * * * * * * * * * * *
194  // Create distributedTriSurfaceMesh
195  Random rndGen(653213);
196 
197  // Determine mesh bounding boxes:
198  List<treeBoundBox> meshBb
199  (
200  1,
201  treeBoundBox
202  (
203  boundBox(mesh_.points(), false)
204  ).extend(rndGen, 1e-3)
205  );
206 
207  // Dummy bounds dictionary
208  dictionary dict;
209  dict.add("bounds", meshBb);
210  dict.add
211  (
212  "distributionType",
214  [
216  ]
217  );
218  dict.add("mergeDistance", SMALL);
219 
220  labelHashSet includePatches;
221  List<labelHashSet> includeAllFacesPerPatch(patches.size());
222 
223  forAll(patches, patchI)
224  {
225  const polyPatch& pp = patches[patchI];
226  if (!pp.coupled() && !isA<cyclicAMIPolyPatch>(pp))
227  {
228  includePatches.insert(patchI);
229 
230  const tmp<scalarField> tt =
231  boundaryRadiation.transmissivity(patchI);
232  const scalarField& tau = tt();
233 
234  forAll(pp, faceI)
235  {
236  if (tau[faceI] == 0.0)
237  {
238  includeAllFacesPerPatch[patchI].insert(faceI);
239  }
240  }
241  }
242  }
243 
244  triSurface localSurface = triangulate
245  (
246  includePatches,
247  includeAllFacesPerPatch
248  );
249 
250  distributedTriSurfaceMesh surfacesMesh
251  (
252  IOobject
253  (
254  "opaqueSurface.stl",
255  mesh_.time().constant(), // directory
256  "triSurface", // instance
257  mesh_.time(), // registry
260  ),
261  localSurface,
262  dict
263  );
264 
265  if (debug)
266  {
267  surfacesMesh.searchableSurface::write();
268  }
269 
270  scalar maxBounding = 5.0*mag(mesh_.bounds().max() - mesh_.bounds().min());
271 
272  reduce(maxBounding, maxOp<scalar>());
273 
274  // Calculate index of faces which have a direct hit (local)
275  DynamicList<label> rayStartFace(nFaces + 0.01*nFaces);
276 
277  // Shoot Rays
278  // * * * * * * * * * * * * * * * *
279  {
280 
281  DynamicField<point> start(nFaces);
282  DynamicField<point> end(start.size());
283  DynamicList<label> startIndex(start.size());
284 
285  label i = 0;
286  do
287  {
288  for (; i < Cfs.size(); i++)
289  {
290  const point& fc = Cfs[i];
291 
292  const label myFaceId = hitFacesIds[i];
293 
294  const vector d(direction_*maxBounding);
295 
296  start.append(fc - 0.001*d);
297 
298  startIndex.append(myFaceId);
299 
300  end.append(fc - d);
301 
302  }
303 
304  }while (returnReduce(i < Cfs.size(), orOp<bool>()));
305 
306  List<pointIndexHit> hitInfo(startIndex.size());
307  surfacesMesh.findLine(start, end, hitInfo);
308 
309  // Collect the rays which has 'only one not wall' obstacle between
310  // start and end.
311  // If the ray hit itself get stored in dRayIs
312  forAll(hitInfo, rayI)
313  {
314  if (!hitInfo[rayI].hit())
315  {
316  rayStartFace.append(startIndex[rayI]);
317  }
318  }
319 
320  // Plot all rays between visible faces.
321  if (debug)
322  {
323  writeRays
324  (
325  mesh_.time().path()/"allVisibleFaces.obj",
326  end,
327  start
328  );
329  }
330 
331  start.clear();
332  startIndex.clear();
333  end.clear();
334  }
335 
336  rayStartFaces_.transfer(rayStartFace);
337 
338  if (debug)
339  {
340  tmp<volScalarField> thitFaces
341  (
342  new volScalarField
343  (
344  IOobject
345  (
346  "hitFaces",
347  mesh_.time().timeName(),
348  mesh_,
351  ),
352  mesh_,
354  )
355  );
356 
357  volScalarField& hitFaces = thitFaces.ref();
358  volScalarField::Boundary& hitFacesBf = hitFaces.boundaryFieldRef();
359 
360  hitFacesBf = 0.0;
361  for (const label faceI : rayStartFaces_)
362  {
363  label patchID = patches.whichPatch(faceI);
364  const polyPatch& pp = patches[patchID];
365  hitFacesBf[patchID][faceI - pp.start()] = 1.0;
366  }
367  hitFaces.write();
368  }
369 
370  label totalHitFaces = rayStartFaces_.size();
371 
372  reduce(totalHitFaces, sumOp<label>());
373 
374  Info<< "Total number of hit faces : " << totalHitFaces << endl;
375 }
376 
377 
378 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
379 
380 Foam::faceShading::faceShading
381 (
382  const fvMesh& mesh,
383  const vector dir,
384  const labelList& hitFaceList
385 )
386 :
387  mesh_(mesh),
388  direction_(dir),
389  rayStartFaces_(hitFaceList)
390 {}
391 
392 
393 Foam::faceShading::faceShading
394 (
395  const fvMesh& mesh,
396  const vector dir
397 )
398 :
399  mesh_(mesh),
400  direction_(dir),
401  rayStartFaces_(0)
402 {
403  calculate();
404 }
405 
406 
407 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
408 
410 {
411  calculate();
412 }
413 
414 
415 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
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
distributedTriSurfaceMesh.H
Foam::distributedTriSurfaceMesh::distributionTypeNames_
static const Enum< distributionType > distributionTypeNames_
Definition: distributedTriSurfaceMesh.H:101
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
surfacesMesh
distributedTriSurfaceMesh surfacesMesh(IOobject("wallSurface.stl", runTime.constant(), "triSurface", runTime, IOobject::NO_READ, IOobject::NO_WRITE), localSurface, dict)
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::List::append
void append(const T &val)
Append an element at the end of the list.
Definition: ListI.H:175
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
meshBb
List< treeBoundBox > meshBb(1, treeBoundBox(boundBox(coarseMesh.points(), false)).extend(rndGen, 1e-3))
Foam::Pout
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
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::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:51
Foam::regIOobject::write
virtual bool write(const bool valid=true) const
Write using setting from DB.
Definition: regIOobjectWrite.C:132
Foam::triSurface
Triangulated surface description with patch information.
Definition: triSurface.H:76
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::bMesh
PrimitivePatch< List< face >, const pointField > bMesh
Holder of faceList and points. (v.s. e.g. primitivePatch which references points)
Definition: bMesh.H:48
patchID
label patchID
Definition: boundaryProcessorFaPatchPoints.H:5
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
Foam::polyBoundaryMesh::whichPatch
label whichPatch(const label faceIndex) const
Return patch index for a given face label.
Definition: polyBoundaryMesh.C:812
faceShading.H
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
dict
dictionary dict
Definition: searchingEngine.H:14
os
OBJstream os(runTime.globalPath()/outputName)
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
fvMesh.H
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
cyclicAMIPolyPatch.H
Foam::linePointRef
line< point, const point & > linePointRef
A line using referred points.
Definition: linePointRef.H:47
f
labelList f(nPoints)
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::faceList
List< face > faceList
A List of faces.
Definition: faceListFwd.H:47
Foam::GeometricField::ref
Internal & ref(const bool updateAccessTime=true)
Return a reference to the dimensioned internal field.
Definition: GeometricField.C:749
Foam::Vector< scalar >
Foam::List< label >
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::faceShading::correct
void correct()
Recalculate rayStartFaces using direction vector.
Definition: faceShading.C:409
patches
const polyBoundaryMesh & patches
Definition: convertProcessorPatches.H:65
rndGen
Random rndGen
Definition: createFields.H:23
Foam::dictionary::add
entry * add(entry *entryPtr, bool mergeEntry=false)
Add a new entry.
Definition: dictionary.C:640
Foam::point
vector point
Point is a vector.
Definition: point.H:43
Foam::labelHashSet
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition: HashSet.H:85
Foam::IOobject::NO_READ
Definition: IOobject.H:188
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
OBJstream.H
boundaryRadiationProperties.H
Foam::dimless
const dimensionSet dimless
Dimensionless.
Definition: dimensionSets.C:189