Loading...
Searching...
No Matches
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-2022 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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"
34#include "surfaceFields.H"
36#include "OBJstream.H"
37
38// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39
40namespace Foam
41{
43}
44
45// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
46
47void Foam::faceShading::calculate()
48{
49 const auto& pbm = mesh_.boundaryMesh();
50
51 const bitSet isOpaqueFace
52 (
53 selectOpaqueFaces
54 (
56 patchIDs_,
57 zoneIDs_
58 )
59 );
60
61 // Find faces potentially hit by solar rays
62 // - correct normal
63 // - transmissivity 0
64 labelList hitFacesIds;
65 bitSet hitFacesFlips;
66 selectFaces
67 (
68 true, // use normal to do first filtering
69 isOpaqueFace,
70 patchIDs_,
71 zoneIDs_,
72
73 hitFacesIds,
74 hitFacesFlips
75 );
76
77 Info<< "Number of 'potential' direct hits : "
78 << returnReduce(hitFacesIds.size(), sumOp<label>()) << endl;
79
80
81 // * * * * * * * * * * * * * * *
82 // Create distributedTriSurfaceMesh
83 Random rndGen(653213);
84
85 // Find potential obstructions. Include all faces that might potentially
86 // block (so ignore normal)
87 labelList blockingFacesIds;
88 bitSet blockingFacesFlips;
89 selectFaces
90 (
91 false, // use normal to do first filtering
92 isOpaqueFace,
93 patchIDs_,
94 zoneIDs_,
95
96 blockingFacesIds,
97 blockingFacesFlips
98 );
99
100 const triSurface localSurface = triangulate
101 (
102 blockingFacesIds,
103 blockingFacesFlips
104 );
105
106 // Determine mesh bounding boxes:
107 List<treeBoundBox> meshBb
108 (
109 1,
110 treeBoundBox(mesh_.points()).extend(rndGen, 1e-3)
111 );
112
113 // Dummy bounds dictionary
115 dict.add("bounds", meshBb);
116 dict.add
117 (
118 "distributionType",
120 [
122 ]
123 );
124 dict.add("mergeDistance", SMALL);
125
126 distributedTriSurfaceMesh surfacesMesh
127 (
128 IOobject
129 (
130 "opaqueSurface.stl",
131 mesh_.time().constant(), // directory
132 "triSurface", // instance
133 mesh_.time(), // registry
136 ),
138 dict
139 );
140
141 if (debug)
142 {
143 surfacesMesh.searchableSurface::write();
144 }
145
146 const scalar maxBounding =
147 returnReduce(5.0*mesh_.bounds().mag(), maxOp<scalar>());
148
149 // Calculate index of faces which have a direct hit (local)
150
151 // Shoot Rays
152 // * * * * * * * * * * * * * * * *
153 {
154
155 DynamicField<point> start(hitFacesIds.size());
156 DynamicField<point> end(start.size());
157 DynamicList<label> startIndex(start.size());
158
159 const pointField& faceCentres = mesh_.faceCentres();
160
161 const vector d(direction_*maxBounding);
162
163 forAll(hitFacesIds, i)
164 {
165 const label facei = hitFacesIds[i];
166 const point& fc = faceCentres[facei];
167
168 start.append(fc - 0.001*d);
169 startIndex.append(facei);
170 end.append(fc - d);
171 }
172
173 List<pointIndexHit> hitInfo(startIndex.size());
174 surfacesMesh.findLine(start, end, hitInfo);
175
176 // Collect the rays which has 'only one not wall' obstacle between
177 // start and end.
178 // If the ray hit itself get stored in dRayIs
179
180 label nVisible = 0;
181 forAll(hitInfo, rayI)
182 {
183 if (!hitInfo[rayI].hit())
184 {
185 nVisible++;
186 }
187 }
188 rayStartFaces_.setSize(nVisible);
189 nVisible = 0;
190
191 forAll(hitInfo, rayI)
192 {
193 if (!hitInfo[rayI].hit())
194 {
195 rayStartFaces_[nVisible++] = startIndex[rayI];
196 }
197 }
198
199 // Plot all rays between visible faces.
200 if (debug)
201 {
202 writeRays
203 (
204 mesh_.time().path()/"allVisibleFaces.obj",
205 end,
206 start
207 );
208 }
209
210 start.clear();
211 startIndex.clear();
212 end.clear();
213 }
214
215 if (debug)
216 {
217 auto thitFaces = surfaceScalarField::New
218 (
219 "hitFaces",
221 mesh_,
223 );
224 auto& hitFaces = thitFaces.ref();
225
226 surfaceScalarField::Boundary& hitFacesBf = hitFaces.boundaryFieldRef();
227
228 hitFacesBf = 0.0;
229 for (const label facei : rayStartFaces_)
230 {
231 const label patchID = pbm.whichPatch(facei);
232 if (patchID == -1)
233 {
234 hitFaces[facei] = 1.0;
235 }
236 else
237 {
238 const polyPatch& pp = pbm[patchID];
239 hitFacesBf[patchID][facei - pp.start()] = 1.0;
240 }
241 }
242 hitFaces.write();
243 }
244
245 Info<< "Total number of hit faces : "
246 << returnReduce(rayStartFaces_.size(), sumOp<label>()) << endl;
247}
248
249
250Foam::triSurface Foam::faceShading::triangulate
251(
252 const labelUList& faceIDs,
253 const bitSet& flipMap
254) const
255{
256 if (faceIDs.size() != flipMap.size())
257 {
258 FatalErrorInFunction << "Size problem :"
259 << "faceIDs:" << faceIDs.size()
260 << "flipMap:" << flipMap.size()
261 << exit(FatalError);
262 }
263
264 const auto& points = mesh_.points();
265 const auto& faces = mesh_.faces();
266 const auto& bMesh = mesh_.boundaryMesh();
267 const auto& fzs = mesh_.faceZones();
268
269 // Patching of surface:
270 // - non-processor patches
271 // - faceZones
272 // Note: check for faceZones on boundary? Who has priority?
273 geometricSurfacePatchList surfPatches(bMesh.nNonProcessor()+fzs.size());
274 labelList patchID(mesh_.nFaces(), -1);
275 {
276 label newPatchi = 0;
277 for (label patchi = 0; patchi < bMesh.nNonProcessor(); ++patchi)
278 {
279 const auto& pp = bMesh[patchi];
280
281 surfPatches[newPatchi] = geometricSurfacePatch
282 (
283 pp.name(),
284 newPatchi,
285 pp.type()
286 );
288 (
289 patchID,
290 pp.size(),
291 pp.start()
292 ) = newPatchi;
293
294 newPatchi++;
295 }
296 for (const auto& fz : fzs)
297 {
298 surfPatches[newPatchi] = geometricSurfacePatch
299 (
300 fz.name(),
301 newPatchi,
302 fz.type()
303 );
304 UIndirectList<label>(patchID, fz) = newPatchi;
305
306 newPatchi++;
307 }
308 }
309
310
311 // Storage for surfaceMesh. Size estimate.
312 DynamicList<labelledTri> triangles(2*faceIDs.size());
313
314 // Work array
315 faceList triFaces;
316
317 forAll(faceIDs, i)
318 {
319 const label facei = faceIDs[i];
320 const bool flip = flipMap[i];
321 const label patchi = patchID[facei];
322 const face& f = faces[facei];
323
324 // Triangulate face
325 triFaces.setSize(f.nTriangles(points));
326 label nTri = 0;
327 f.triangles(points, nTri, triFaces);
328
329 for (const face& f : triFaces)
330 {
331 if (!flip)
332 {
333 triangles.append(labelledTri(f[0], f[1], f[2], patchi));
334 }
335 else
336 {
337 triangles.append(labelledTri(f[0], f[2], f[1], patchi));
338 }
339 }
340 }
341
342 triangles.shrink();
343
344 // Create globally numbered tri surface
345 triSurface rawSurface(triangles, mesh_.points());
346
347 // Create locally numbered tri surface
349 (
350 rawSurface.localFaces(),
351 rawSurface.localPoints()
352 );
353
354 // Add patch names to surface
355 surface.patches().transfer(surfPatches);
356
357 return surface;
358}
359
360
361Foam::bitSet Foam::faceShading::selectOpaqueFaces
362(
363 const radiation::boundaryRadiationProperties& boundaryRadiation,
364 const labelUList& patchIDs,
365 const labelUList& zoneIDs
366) const
367{
368 const auto& pbm = mesh_.boundaryMesh();
369
370 bitSet isOpaqueFace(mesh_.nFaces(), false);
371
372 // Check selected patches
373 for (const label patchi : patchIDs)
374 {
375 const auto& pp = pbm[patchi];
376 tmp<scalarField> tt = boundaryRadiation.transmissivity(patchi);
377 const scalarField& t = tt.cref();
378
379 forAll(t, i)
380 {
381 isOpaqueFace[i + pp.start()] = (t[i] == 0.0);
382 }
383 }
384
385 // Check selected faceZones
386 const auto& fzs = mesh_.faceZones();
387
388 for (const label zonei : zoneIDs)
389 {
390 const auto& fz = fzs[zonei];
391
392 //- Note: slice mesh face centres preferentially
393 tmp<scalarField> tt = boundaryRadiation.zoneTransmissivity
394 (
395 zonei,
396 fz
397 );
398 const scalarField& t = tt.cref();
399
400 forAll(t, i)
401 {
402 isOpaqueFace[fz[i]] = (t[i] == 0.0);
403 }
404 }
405
406 return isOpaqueFace;
407}
408
409
410void Foam::faceShading::selectFaces
411(
412 const bool useNormal,
413 const bitSet& isCandidateFace,
414 const labelUList& patchIDs,
415 const labelUList& zoneIDs,
416
417 labelList& faceIDs,
418 bitSet& flipMap
419) const
420{
421 const auto& pbm = mesh_.boundaryMesh();
422
423 bitSet isSelected(mesh_.nFaces());
424 DynamicList<label> dynFaces(mesh_.nBoundaryFaces());
425 bitSet isFaceFlipped(mesh_.nFaces());
426
427 // Add patches
428 for (const label patchi : patchIDs)
429 {
430 const auto& pp = pbm[patchi];
431 const vectorField& n = pp.faceNormals();
432
433 forAll(n, i)
434 {
435 const label meshFacei = i + pp.start();
436 if
437 (
438 isCandidateFace[meshFacei]
439 && (
440 !useNormal
441 || ((direction_ & n[i]) > 0)
442 )
443 )
444 {
445 isSelected.set(meshFacei);
446 isFaceFlipped[meshFacei] = false;
447 dynFaces.append(meshFacei);
448 }
449 }
450 }
451
452
453 // Add faceZones
454 const auto& fzs = mesh_.faceZones();
455
456 for (const label zonei : zoneIDs)
457 {
458 const auto& fz = fzs[zonei];
459 const primitiveFacePatch& pp = fz();
460 const vectorField& n = pp.faceNormals();
461
462 forAll(n, i)
463 {
464 const label meshFacei = fz[i];
465
466 if
467 (
468 !isSelected[meshFacei]
469 && isCandidateFace[meshFacei]
470 && (
471 !useNormal
472 || ((direction_ & n[i]) > 0)
473 )
474 )
475 {
476 isSelected.set(meshFacei);
477 dynFaces.append(meshFacei);
478 isFaceFlipped[meshFacei] = fz.flipMap()[i];
479 }
480 }
481 }
482 faceIDs = std::move(dynFaces);
483 flipMap = bitSet(isFaceFlipped, faceIDs);
484}
485
486
487void Foam::faceShading::writeRays
488(
489 const fileName& fName,
490 const DynamicField<point>& endCf,
491 const pointField& myFc
492)
493{
494 OBJstream os(fName);
495
496 Pout<< "Dumping rays to " << os.name() << endl;
497
498 forAll(myFc, facei)
499 {
500 os.writeLine(myFc[facei], endCf[facei]);
501 }
503
504
505// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
506
507Foam::faceShading::faceShading
508(
509 const fvMesh& mesh,
510 const vector& dir
511)
512:
513 mesh_(mesh),
514 patchIDs_(nonCoupledPatches(mesh)),
515 zoneIDs_(0),
516 direction_(dir),
517 rayStartFaces_(0)
519 calculate();
520}
521
522
523Foam::faceShading::faceShading
524(
525 const fvMesh& mesh,
526 const labelList& patchIDs,
527 const labelList& zoneIDs,
528 const vector& dir
529)
530:
531 mesh_(mesh),
532 patchIDs_(patchIDs),
533 zoneIDs_(zoneIDs),
534 direction_(dir),
535 rayStartFaces_(0)
536{
537 calculate();
539
540
541// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
542
544{
545 const auto& pbm = mesh.boundaryMesh();
546
547 DynamicList<label> ncPatches;
548 forAll(pbm, patchi)
549 {
550 const polyPatch& pp = pbm[patchi];
551 if (!pp.coupled() && !isA<cyclicAMIPolyPatch>(pp))
552 {
553 ncPatches.append(patchi);
554 }
556 return ncPatches;
557}
558
559
561{
562 calculate();
563}
564
565
566// ************************************************************************* //
label n
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
labelList patchIDs
const polyBoundaryMesh & pbm
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.
static tmp< GeometricField< scalar, fvsPatchField, surfaceMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvsPatchField< scalar >::calculatedType())
GeometricBoundaryField< scalar, fvsPatchField, surfaceMesh > Boundary
@ NO_REGISTER
Do not request registration (bool: false).
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
label size() const noexcept
The number of elements in the list.
void transfer(List< T > &list)
Transfer the contents of the argument List into this list and annul the argument list.
Definition List.C:347
void setSize(label n)
Alias for resize().
Definition List.H:536
static FOAM_NO_DANGLING_REFERENCE const boundaryRadiationProperties & New(const fvMesh &mesh, Args &&... args)
An OFstream that keeps track of vertices and provides convenience output methods for OBJ files.
Definition OBJstream.H:58
Ostream & writeLine(const point &p0, const point &p1)
Write line joining two points.
Definition OBJstream.C:234
virtual const fileName & name() const override
Read/write access to the name of the stream.
Definition OSstream.H:134
const Field< point_type > & faceNormals() const
Return face unit normals for patch.
A non-owning sub-view of a List (allocated or unallocated storage).
Definition SubList.H:61
const word & constant() const noexcept
Return constant name.
Definition TimePathsI.H:131
A List with indirect addressing. Like IndirectList but does not store addressing.
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition bitSet.H:61
static const Enum< distributionType > distributionTypeNames_
Helper class to calculate visible faces for global, sun-like illumination.
Definition faceShading.H:58
void correct()
Recalculate rayStartFaces using direction vector.
static labelList nonCoupledPatches(const polyMesh &mesh)
Helper: return all uncoupled patches.
A face is a list of labels corresponding to mesh vertices.
Definition face.H:71
A class for handling file names.
Definition fileName.H:75
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
Identifies a surface patch/zone by name and index, with geometric type.
A triFace with additional (region) index.
Definition labelledTri.H:56
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition polyMesh.H:609
virtual const pointField & points() const
Return raw points.
Definition polyMesh.C:1063
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
A class for managing temporary objects.
Definition tmp.H:75
Triangulated surface description with patch information.
Definition triSurface.H:74
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
const pointField & points
const labelIOList & zoneIDs
Definition correctPhi.H:59
const char * end
Definition SVGTools.H:223
const wordList surface
Standard surface field types (scalar, vector, tensor, etc).
Namespace for OpenFOAM.
const dimensionSet dimless
Dimensionless.
List< label > labelList
A List of labels.
Definition List.H:62
messageStream Info
Information stream (stdout output on master, null elsewhere).
List< face > faceList
List of faces.
Definition faceListFwd.H:41
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.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
PrimitivePatch< List< face >, const pointField & > primitiveFacePatch
A PrimitivePatch with List storage for the faces, const reference for the point field.
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
Field< vector > vectorField
Specialisation of Field<T> for vector.
vector point
Point is a vector.
Definition point.H:37
List< geometricSurfacePatch > geometricSurfacePatchList
List of geometricSurfacePatch.
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...
PrimitivePatch< List< face >, const pointField > bMesh
Holder of faceList and points. (v.s. e.g. primitivePatch which references points).
Definition bMesh.H:39
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
vectorField pointField
pointField is a vectorField.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Vector< scalar > vector
Definition vector.H:57
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
UList< label > labelUList
A UList of labels.
Definition UList.H:75
labelList f(nPoints)
dictionary dict
const triSurface localSurface
List< treeBoundBox > meshBb(1, treeBoundBox(coarseMesh.points()).extend(rndGen, 1e-3))
distributedTriSurfaceMesh surfacesMesh(IOobject("wallSurface.stl", runTime.constant(), "triSurface", runTime, IOobject::NO_READ, IOobject::NO_WRITE), localSurface, dict)
std::vector< Triangle > triangles
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Foam::surfaceFields.
Random rndGen