Loading...
Searching...
No Matches
voxelRaySearchEngine.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) 2023-2024 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
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
29#include "processorPolyPatch.H"
31
32// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33
34namespace Foam
35{
36namespace VF
37{
40}
41}
42
43// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
44
45void Foam::VF::voxel::setTriangulation(const fvMesh& mesh)
46{
47 Info<< "\nCreating triangulated surface" << endl;
48
49 // Storage for surfaceMesh. Size estimate.
50 DynamicList<labelledTri> triangles(mesh.nBoundaryFaces());
51 DynamicList<label> globalFaces(mesh.nBoundaryFaces());
52 label nFace = 0;
53
54 const auto& pbm = mesh.boundaryMesh();
55
57 {
58 const label patchi = patchIDs_[i];
59 const polyPatch& patch = pbm[patchi];
60 const pointField& points = patch.points();
61
62 for (const face& f : patch)
63 {
64 label nTri = 0;
65 faceList triFaces(f.nTriangles(points));
66 f.triangles(points, nTri, triFaces);
67
68 const label globalFacei =
69 globalNumbering_.toGlobal(Pstream::myProcNo(), nFace++);
70
71 for (const face& f : triFaces)
72 {
73 triangles.push_back(labelledTri(f[0], f[1], f[2], i));
74 globalFaces.push_back(globalFacei);
75 }
76 }
77 }
78
79 triToGlobalFace_.transfer(globalFaces);
80
81 Info<< " Total number of triangles: "
82 << returnReduce(triangles.size(), sumOp<label>())
83 << endl;
84
85 triangles.shrink();
86 surface_ = triSurface(triangles, mesh.points());
87 surface_.compactPoints();
88}
89
90
91Foam::labelList Foam::VF::voxel::setFaceVertexHits
92(
93 const fvMesh& mesh,
94 const labelList& patchIDs
95)
96{
97 labelList vertHits(mesh.points().size(), Zero);
98
99 if (mesh.nSolutionD() == 3)
100 {
101 const auto& pbm = mesh.boundaryMesh();
102
103 // Create a new triangulation based on the surface agglomeration
104 for (const label patchI : patchIDs)
105 {
106 const polyPatch& patch = pbm[patchI];
107 for (const face& f : patch)
108 {
109 for (const label fpi : f)
110 {
111 ++vertHits[fpi];
112 }
113 }
114 }
115
116 for (const auto& pp : pbm)
117 {
118 const labelList& meshPts = pp.meshPoints();
119
120 if (pp.size())
121 {
123 {
124 // Add all processor patch points
125 for (const label pi : meshPts)
126 {
127 ++vertHits[pi];
128 }
129 }
130 else
131 {
132 // Add boundary points
133
134 const auto& bndyPts = pp.boundaryPoints();
135
136 for (const label pi : bndyPts)
137 {
138 ++vertHits[meshPts[pi]];
139 }
140 }
141 }
142 }
143 }
144
145 return vertHits;
146}
147
148
149void Foam::VF::voxel::setCoarseTriangulation(const fvMesh& mesh)
150{
151 Info<< "\nCreating triangulated surface" << endl;
152
153
154 // Filter out fine mesh points along coarse mesh faces
155
156
157 // Storage for surfaceMesh. Size estimate.
158 DynamicList<labelledTri> triangles(mesh.nBoundaryFaces());
159 DynamicList<label> globalFaces(mesh.nBoundaryFaces());
160 labelList vertHits = setFaceVertexHits(mesh, patchIDs_);
161
162
163 // Only simplifying edges for 3D
164 const label nVertMin = mesh.nSolutionD() == 3 ? 2 : 0;
165
166 label nInvalid = 0;
167 label nFace = 0;
168
169 const auto& pbm = mesh.boundaryMesh();
170
171 for (const label patchi : patchIDs_)
172 {
173 const polyPatch& patch = pbm[patchi];
174 const pointField& points = patch.points();
175
176 for (const face& f : patch)
177 {
178 DynamicList<label> faceVerts;
179 for (const label fpi : f)
180 {
181 if (vertHits[fpi] > nVertMin)
182 {
183 faceVerts.push_back(fpi);
184 }
185 }
186
187 if (faceVerts.size() < 3)
188 {
189 ++nInvalid;
190 continue;
191 }
192
193 label nTri = 0;
194 const face reducedFace(faceVerts);
195 faceList triFaces(reducedFace.nTriangles(points));
196 reducedFace.triangles(points, nTri, triFaces);
197
198 const label globalFacei =
199 globalNumbering_.toGlobal(Pstream::myProcNo(), nFace++);
200
201 for (const face& f : triFaces)
202 {
203 triangles.push_back(labelledTri(f[0], f[1], f[2], patchi));
204 globalFaces.push_back(globalFacei);
205 }
206 }
207 }
208
209 triToGlobalFace_.transfer(globalFaces);
210
211 Info<< " Total number of triangles: "
212 << returnReduce(triangles.size(), sumOp<label>())
213 << "\n Number of invalid (removed) triangles: "
214 << returnReduce(nInvalid, sumOp<label>())
215 << endl;
216
217 triangles.shrink();
218 surface_ = triSurface(triangles, mesh.points());
219 surface_.compactPoints();
220}
221
222
223void Foam::VF::voxel::broadcast()
224{
225 // Every processor has whole surface
226 const globalIndex globalFaceIdx(globalIndex::gatherOnly{}, surface_.size());
227 const globalIndex globalPointIdx
228 (
229 globalIndex::gatherOnly{},
230 surface_.points().size()
231 );
232
233 List<labelledTri> globalSurfaceTris(globalFaceIdx.gather(surface_));
234 pointField globalSurfacePoints(globalPointIdx.gather(surface_.points()));
235 List<label> globalTriToGlobalFace(globalFaceIdx.gather(triToGlobalFace_));
236
237
238 for (const label proci : globalPointIdx.allProcs())
239 {
240 const label offset = globalPointIdx.localStart(proci);
241
242 if (offset)
243 {
244 for
245 (
246 labelledTri& tri
247 : globalSurfaceTris.slice(globalFaceIdx.range(proci))
248 )
249 {
250 tri[0] += offset;
251 tri[1] += offset;
252 tri[2] += offset;
253 }
254 }
255 }
256
257 surface_ =
258 triSurface
259 (
260 std::move(globalSurfaceTris),
261 std::move(globalSurfacePoints)
262 );
263
264 Pstream::broadcast(surface_);
265
266 triToGlobalFace_ = std::move(globalTriToGlobalFace);
267 Pstream::broadcast(triToGlobalFace_);
268}
269
270
271Foam::pointHit Foam::VF::voxel::rayTriIntersect
272(
273 const label trii,
274 const point& origin,
275 const vector& direction
276) const
277{
278 // Note: origin was made relative to voxel mesh on entry to hit function
279 // - need to convert back into global position to be consistent with
280 // triangles for intersection tests
281
282 const auto& ind = surface_[trii];
283 const auto& pts = surface_.points();
284
285 // Note: flipped orientation of triangle (into domain) so that we can use
286 // visibility check when performing ray-triangle intersections
287 const triPointRef tri(pts[ind[0]], pts[ind[2]], pts[ind[1]]);
288
289 static scalar tolerance = 1e-6;
290
291 return
292 tri.intersection
293 (
294 globalPosition(origin),
295 direction,
297 tolerance
298 );
299}
300
301
302void Foam::VF::voxel::writeBox
303(
304 OBJstream& os,
305 bool lines,
306 const boundBox& bb
307) const
308{
309 os.write(treeBoundBox(bb), lines);
310}
311
312
313void Foam::VF::voxel::writeVoxels(const word& fName) const
314{
315 if (!UPstream::master()) return;
316
317 OBJstream os(fName);
318 Info<< "Writing voxels to " << os.name() << endl;
319
320 boundBox bb;
321 const bool lines = true;
322 for (label i = 0; i < nijk_[0]; ++i)
323 {
324 for (label j = 0; j < nijk_[1]; ++j)
325 {
326 for (label k = 0; k < nijk_[2]; ++k)
327 {
328 bb.min() = voxelMin(i, j, k);
329 bb.max() = voxelMax(i, j, k);
330 writeBox(os, lines, bb);
331 }
332 }
333 }
334
335 Info<< "- done" << endl;
336}
337
338
339void Foam::VF::voxel::writeTriBoundBoxes(const word& fName) const
340{
341 if (!UPstream::master()) return;
342
343 OBJstream os(fName);
344 Info<< "Writing triangle boundBoxes to " << os.name() << endl;
345
346 const bool lines = true;
347 for (const auto& voxeli : objects_)
348 {
349 for (const label trii : voxeli)
350 {
351 writeBox(os, lines, objectBbs_[trii]);
352 }
353 }
354
355 Info<< "- done" << endl;
356}
357
358
359void Foam::VF::voxel::writeHitObjects
360(
361 const label voxeli,
362 const point& origin,
363 const vector& dir
364) const
365{
366 OBJstream os("voxel_" + Foam::name(voxeli) + ".obj");
367
368 // Write voxel
369 labelVector ijk = this->ijk(voxeli);
370
371 boundBox voxelBb;
372 voxelBb.min() = voxelMin(ijk[0], ijk[1], ijk[2]);
373 voxelBb.max() = voxelMax(ijk[0], ijk[1], ijk[2]);
374
375 writeBox(os, true, voxelBb);
376
377 for (const auto& trii : objects_[voxeli])
378 {
379 writeBox(os, true, objectBbs_[trii]);
380
381 const auto& ind = surface_[trii];
382 const auto& pts = surface_.points();
383 const triPointRef tri(pts[ind[0]], pts[ind[1]], pts[ind[2]]);
384 os.write(tri, false);
385 }
386}
387
388
389Foam::pointIndexHit Foam::VF::voxel::hitObject
390(
391 const label voxeli,
392 const point& origin,
393 const vector& dir,
394 const vector& t,
395 const scalar minDistance
396) const
397{
398 if (objects_[voxeli].empty()) return pointIndexHit();
399
400 // boundBox rayBb;
401 // rayBb.add(origin);
402 // rayBb.add(origin + dir*(dir & t));
403
404 label triHiti = -1;
405 //rayBb.add(origin + dir);
406 //rayBb.inflate(0.01);
407
408 if (debug > 2)
409 {
410 writeHitObjects(voxeli, origin, dir);
411 }
412
413 // Determine all triangles that intersect with ray
414 // - only keep nearest hit
415
416 pointHit closestHit;
417 for (const auto& trii : objects_[voxeli])
418 {
419 // Only perform ray/tri intersection if bound boxes intersect
420 //if (objectBbs_[trii].overlaps(rayBb))
421 {
422 pointHit pi = rayTriIntersect(trii, origin, dir);
423
424 if
425 (
426 pi.hit()
427 && (
428 pi.distance() > minDistance
429 && pi.distance() < closestHit.distance()
430 )
431 )
432 {
433 triHiti = trii;
434 closestHit = pi;
435 }
436 }
437 }
438
439 return pointIndexHit(closestHit, triHiti);
440}
441
442
443// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
444
445Foam::VF::voxel::voxel(const fvMesh& mesh, const dictionary& dict)
446:
448 bb0_(),
449 span0_(Zero),
450 nijk_(Zero),
451 dxyz_(Zero),
452 nRayPerFace_(dict.get<label>("nRayPerFace")),
453 nTriPerVoxelMax_(dict.getOrDefault<label>("nTriPerVoxelMax", 50)),
454 depthMax_(dict.getOrDefault<label>("depthMax", 5)),
455 objects_(),
456 objectBbs_()
457{
458 if (agglomMeshPtr_)
459 {
460 setCoarseTriangulation(*agglomMeshPtr_);
461 }
462 else
463 {
464 setTriangulation(mesh);
465 }
466
467 broadcast();
468
469 objectBbs_.resize_nocopy(surface_.size());
470
471 bb0_.add(surface_.points());
472 bb0_.inflate(0.01);
473 span0_ = bb0_.span();
474
475 {
476 scalar maxDx = span0_.x();
477 scalar maxDy = span0_.y();
478 scalar maxDz = span0_.z();
479
480 const label refDn = 8;
481 scalar maxDim = max(maxDx, max(maxDy, maxDz));
482
483 setVoxelDims
484 (
485 refDn*maxDx/maxDim,
486 refDn*maxDy/maxDim,
487 refDn*maxDz/maxDim
488 );
489
490 objects_.resize_nocopy(nVoxel());
491 }
492
493 label depth = 0;
494 label trii = 0;
495 voxelise(objects_, trii, depth);
496
497 Info<< "\nCreated voxel mesh: " << nijk_ << endl;
498
499 if ((debug > 3) && UPstream::master())
500 {
501 writeVoxels("voxels.obj");
502 writeTriBoundBoxes("triBoundBoxes.obj");
503 }
504}
505
506
507// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
508
509void Foam::VF::voxel::refineObjects
510(
511 List<DynamicList<label>>& objects,
512 const label triMax
513)
514{
515 refineVoxelDims();
516
517 if (debug > 2) Pout<< "Refining voxels: n=" << nijk_ << endl;
518
519 List<DynamicList<label>> objectsNew(objects.size()*8);
520
521 for (label trii = 0; trii <= triMax; ++trii)
522 {
523 addBbToVoxels(objectBbs_[trii], trii, objectsNew);
524 }
525
526 objects.transfer(objectsNew);
527}
528
529
530Foam::label Foam::VF::voxel::addBbToVoxels
531(
532 const boundBox& bb,
533 const label trii,
534 List<DynamicList<label>>& objects
535) const
536{
537 //Pout<< "addBbToVoxels: trii=" << trii << " n=" << nijk_ << endl;
538
539 const point minbb(localPosition(bb.min()));
540 const label i0 = max(0, floor(minbb.x()/dxyz_[0]));
541 const label j0 = max(0, floor(minbb.y()/dxyz_[1]));
542 const label k0 = max(0, floor(minbb.z()/dxyz_[2]));
543
544 const point maxbb(localPosition(bb.max()));
545 const label i1 = min(nijk_[0], ceil(maxbb.x()/dxyz_[0]));
546 const label j1 = min(nijk_[1], ceil(maxbb.y()/dxyz_[1]));
547 const label k1 = min(nijk_[2], ceil(maxbb.z()/dxyz_[2]));
548
549 label nTriMax = 0;
550
551 for (label ii = i0; ii < i1; ++ii)
552 {
553 for (label jj = j0; jj < j1; ++jj)
554 {
555 for (label kk = k0; kk < k1; ++kk)
556 {
557 const label voxeli = this->voxeli(ii, jj, kk);
558
559 objects[voxeli].push_back(trii);
560 nTriMax = max(nTriMax, objects[voxeli].size());
561 }
562 }
563 }
564
565 return nTriMax;
566}
567
568
569void Foam::VF::voxel::voxelise
570(
571 List<DynamicList<label>>& objects,
572 const label trii0,
573 const label depth
574)
575{
576 if (debug > 2)
577 {
578 Pout<< "voxelise - start at tri=" << trii0
579 << " depth=" << depth
580 << endl;
581 }
582
583 const auto& points = surface_.points();
584
585 for (label trii = trii0; trii < surface_.size(); ++trii)
586 {
587 // Triangle bounding box
588 boundBox bb(points, surface_[trii]);
589 bb.inflate(0.01);
590 objectBbs_[trii] = bb;
591
592 const label nVoxelMax = addBbToVoxels(bb, trii, objects);
593
594 // Number of triangles per voxel - if exceed limit, refine voxels...
595 if (nVoxelMax > nTriPerVoxelMax_ && depth < depthMax_)
596 {
597 refineObjects(objects, trii);
598 voxelise(objects, trii + 1, depth + 1);
599 break;
600 }
601 }
602}
603
604
606(
607 const label tri0,
608 const vector& direction0
609) const
610{
611 if (tri0 > surface_.size()-1)
612 {
614 << "Index tri0 out of bounds: " << tri0
615 << ". Surface size: " << surface_.size()
616 << abort(FatalError);
617 }
618
619 return hit(surface_.faceCentres()[tri0], direction0);
620}
621
622
624(
625 const point& origin0,
626 const vector& direction0
627) const
628{
629 // Initialise return value
630 pointIndexHit hitInfo;
631
632 const point origin(localPosition(origin0));
633
634 if (cmptMin(origin) < 0)
635 {
637 << "Point located outside voxel mesh"
638 << " - possible coarsening problem?"
639 << abort(FatalError);
640 }
641
642 if (magSqr(direction0) < ROOTVSMALL)
643 {
645 << "Supplied direction has zero size"
646 << endl;
647
648 return hitInfo;
649 }
650
651 const vector dir(normalised(direction0));
652
653 labelVector ijk(Zero);
654 labelVector step(Zero);
655 vector tDelta(vector::max);
656 vector tMax(vector::max);
657
658 for (direction d = 0; d < 3; ++d)
659 {
660 ijk[d] = floor(origin[d]/dxyz_[d]);
661 step[d] = sign0(dir[d]);
662 if (step[d] != 0)
663 {
664 tDelta[d] = mag(dxyz_[d]/dir[d]);
665
666 scalar voxelMax = (1 + ijk[d] - neg(dir[d]))*dxyz_[d];
667 tMax[d] = (voxelMax - origin[d])/dir[d];
668 }
669 }
670
671 if (debug > 2)
672 {
673 Info<< "surfBb:" << boundBox(surface_.points())
674 << " bb:" << bb0_
675 << " origin" << origin0
676 << " voxel_origin:" << origin
677 << " ijk:" << ijk
678 << " step:" << step
679 << " dxyz:" << dxyz_
680 << " tDelta:" << tDelta
681 << " tMax:" << tMax
682 << endl;
683 }
684
685 auto traverse = [&](const label i)
686 {
687 ijk[i] += step[i];
688 if (outOfBounds(ijk, i)) return false;
689 tMax[i] += tDelta[i];
690 return true;
691 };
692
693
694 while (true)
695 {
696 const label voxeli = this->voxeli(ijk);
697
698 if (debug > 2)
699 {
700 Info<< "ijk:" << ijk
701 << " voxeli:" << voxeli
702 << " t:" << tMax
703 << " objs:" << objects_[voxeli].size()
704 << endl;
705 }
706
707 hitInfo = hitObject(voxeli, origin, dir, tMax);
708
709 if (hitInfo.hit())
710 {
711 // Valid hit
712 break;
713 }
714 else
715 {
716 if (tMax[0] < tMax[1] && tMax[0] < tMax[2])
717 {
718 if (!traverse(0)) break;
719 }
720 else if (tMax[1] < tMax[2])
721 {
722 if (!traverse(1)) break;
723 }
724 else
725 {
726 if (!traverse(2)) break;
727 }
728 }
729 }
730
731 return hitInfo;
732}
733
734
736(
737 labelList& rayStartFaceOut,
738 labelList& rayEndFaceOut
739) const
740{
741 (void)mesh_.time().cpuTimeIncrement();
742
743 const pointField& myFc = allCf_[Pstream::myProcNo()];
744 const vectorField& myArea = allSf_[Pstream::myProcNo()];
745
746 const label nTotalRays = myFc.size()*nRayPerFace_;
747 if (nTotalRays > maxDynListLength)
748 {
750 << "Dynamic list needs more capacity to support " << nRayPerFace_
751 << " rays per face (" << nTotalRays << "). "
752 << "Limit set by maxDynListLength: " << maxDynListLength
753 << abort(FatalError);
754 }
755
756 DynamicList<label> rayStartFace(nTotalRays);
757 DynamicList<label> rayEndFace(rayStartFace.size());
758
759 DynamicList<label> endFacei(nTotalRays);
760 DynamicList<label> startFacei(nTotalRays);
761
762 const pointField hemi(createHemiPoints(nRayPerFace_));
763
764 Info<< "\nShooting rays" << endl;
765
766 const scalar nDiv = 10;
767 const label nStep = floor(myFc.size()/nDiv);
768 labelHashSet localFaceHits;
769
770 for (label stepi=0; stepi<nDiv; ++stepi)
771 {
772 const label step0 = stepi*nStep;
773 const label step1 = stepi == nDiv - 1 ? myFc.size() : (stepi+1)*nStep;
774
775 for (label i = step0; i < step1; ++i)
776 {
777 // Info<< "i/N = " << i << "/" << (myFc.size()-1)
778 // << " step0:" << step0 << " step1:" << step1
779 // << endl;
780
781 localFaceHits.clear();
782
783 const point& origin = myFc[i];
784 const vector dir(-normalised(myArea[i]));
785
786 const coordSystem::cartesian cs = createCoordSystem(origin, dir);
787
788 const vectorField pts(cs.transformPoint(hemi));
789
790 for (const auto& p : pts)
791 {
792 const pointIndexHit hitInfo = hit(origin, p-origin);
793
794 if (hitInfo.hit())
795 {
796 label hitFacei = triToGlobalFace_[hitInfo.index()];
797
798 if (localFaceHits.insert(hitFacei))
799 {
800 endFacei.push_back(hitFacei);
801 startFacei.push_back(i);
802 }
803 }
804 else
805 {
806 // Miss
807 }
808 }
809 }
810
811 const label globalStepi = returnReduce(stepi, minOp<label>()) + 1;
812
813 Info<< " " << globalStepi/nDiv*100 << "% complete" << endl;
814 }
815
816
817 // Ensure all rays are unique/filter out duplicates
818 // - add symmetric connections for non-self-intersecting rays
819
820 if (UPstream::parRun())
821 {
822 edgeHashSet uniqueRays;
823 List<DynamicList<label>> remoteStartFace(Pstream::nProcs());
824 List<DynamicList<label>> remoteEndFace(Pstream::nProcs());
825
826 const labelList globalStartFaceIDs
827 (
828 globalNumbering_.toGlobal(startFacei)
829 );
830
831 forAll(globalStartFaceIDs, rayi)
832 {
833 label start = globalStartFaceIDs[rayi];
834 label end = endFacei[rayi];
835
836 const label endProci = globalNumbering_.whichProcID(end);
837
838 if (start > end) Swap(start, end);
839
840 if (uniqueRays.insert(edge(start, end)))
841 {
842 if (endProci != UPstream::myProcNo())
843 {
844 // Convert local face into global face and vice-versa
845 remoteStartFace[endProci].push_back(start);
846 remoteEndFace[endProci].push_back(end);
847 }
848 }
849 }
850
851 PstreamBuffers pBufs(Pstream::commsTypes::nonBlocking);
852
853 // Send remote data
854 for (const int domain : Pstream::allProcs())
855 {
856 if (domain != Pstream::myProcNo())
857 {
858 UOPstream str(domain, pBufs);
859 str << remoteStartFace[domain]
860 << remoteEndFace[domain];
861 }
862 }
863
864 pBufs.finishedSends();
865
866 // Consume
867 for (const int domain : Pstream::allProcs())
868 {
869 if (domain != Pstream::myProcNo())
870 {
871 UIPstream is(domain, pBufs);
872 is >> remoteStartFace[domain]
873 >> remoteEndFace[domain];
874
875 forAll(remoteStartFace[domain], i)
876 {
877 const label start = remoteStartFace[domain][i];
878 const label end = remoteEndFace[domain][i];
879 uniqueRays.insert(edge(start, end));
880 }
881 }
882 }
883
884 // Populate ray addressing based on uniqueRays
885 for (const edge& e : uniqueRays)
886 {
887 const label start = e.first();
888 const label end = e.second();
889
890 bool sameFace = start == end;
891
892 if (globalNumbering_.isLocal(start))
893 {
894 // Ray originates from this processor
895 const label localStart = globalNumbering_.toLocal(start);
896
897 rayStartFace.append(localStart);
898 rayEndFace.append(end);
899 }
900
901 if (!sameFace && globalNumbering_.isLocal(end))
902 {
903 // Ray hits this processor
904 // - add symmetric ray from end->start
905 const label localEnd = globalNumbering_.toLocal(end);
906
907 rayStartFace.append(localEnd);
908 rayEndFace.append(start);
909 }
910 }
911 }
912 else
913 {
914 // Populate ray addressing based on uniqueRays
915
916 edgeHashSet uniqueRays;
917
918 forAll(startFacei, rayi)
919 {
920 label start = startFacei[rayi];
921 label end = endFacei[rayi];
922
923 if (start > end) Swap(start, end);
924
925 if (uniqueRays.insert(edge(start, end)))
926 {
927 rayStartFace.append(start);
928 rayEndFace.append(end);
929
930 if (start != end)
931 {
932 // Add symmetric ray from end->start
933 rayStartFace.append(end);
934 rayEndFace.append(start);
935 }
936 }
937 }
938 }
939
940 rayStartFaceOut.transfer(rayStartFace);
941 rayEndFaceOut.transfer(rayEndFace);
942
943 Info<< " Time taken: " << mesh_.time().cpuTimeIncrement() << " s"
944 << endl;
945}
946
947
948// ************************************************************************* //
label k
constexpr scalar pi(M_PI)
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
labelList patchIDs
const polyBoundaryMesh & pbm
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Negative if the process is not a...
Definition UPstream.H:1706
@ nonBlocking
"nonBlocking" (immediate) : (MPI_Isend, MPI_Irecv)
Definition UPstream.H:84
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition UPstream.H:1714
static label nProcs(const label communicator=worldComm)
Number of ranks in parallel run (for given communicator). It is 1 for serial run.
Definition UPstream.H:1697
static rangeType allProcs(const label communicator=worldComm)
Range of process indices for all processes.
Definition UPstream.H:1857
@ broadcast
broadcast [MPI]
Definition UPstream.H:189
static bool & parRun() noexcept
Test if this a parallel run.
Definition UPstream.H:1681
Base class for ray search engines.
globalIndex globalNumbering_
Global numbering.
const fvMesh & mesh() const noexcept
Reference to the mesh.
labelList patchIDs_
List of participating patch IDs.
Ray search engine based on discretising space into uniform voxels.
pointIndexHit hit(const point &origin, const vector &dir) const
voxel(const fvMesh &mesh, const dictionary &dict)
Constructor.
virtual void shootRays(labelList &rayStartFaceOut, labelList &rayEndFaceOut) const
Shoot rays; returns lists of ray start and end faces.
static const Form max
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
volScalarField & p
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
const pointField & points
#define WarningInFunction
Report a warning using Foam::Warning.
bool broadcast(Type *values, int count, MPI_Datatype datatype, const int communicator, const int root=0)
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
const char * end
Definition SVGTools.H:223
A namespace for various viewFactor model implementations.
const std::string patch
OpenFOAM patch number as a std::string.
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
PointHit< point > pointHit
A PointHit with a 3D point.
Definition pointHit.H:51
HashSet< edge, Hash< edge > > edgeHashSet
A HashSet with edge for its key. Hashing (and ==) on an edge is symmetric.
Definition edgeHashes.H:48
List< label > labelList
A List of labels.
Definition List.H:62
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
Vector< label > labelVector
Vector of labels.
Definition labelVector.H:47
dimensionedScalar j1(const dimensionedScalar &ds)
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition HashSet.H:85
void cmptMin(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
messageStream Info
Information stream (stdout output on master, null elsewhere).
List< face > faceList
List of faces.
Definition faceListFwd.H:41
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.
void Swap(DynamicList< T, SizeMinA > &a, DynamicList< T, SizeMinB > &b)
Exchange contents of lists - see DynamicList::swap().
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
errorManip< error > abort(error &err)
Definition errorManip.H:139
Field< vector > vectorField
Specialisation of Field<T> for vector.
uint8_t direction
Definition direction.H:49
vector point
Point is a vector.
Definition point.H:37
dimensionedScalar neg(const dimensionedScalar &ds)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
triangle< point, const point & > triPointRef
A triangle using referred points.
Definition triangleFwd.H:46
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
vectorField pointField
pointField is a vectorField.
Vector< scalar > vector
Definition vector.H:57
PointIndexHit< point > pointIndexHit
A PointIndexHit with a 3D point.
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
dimensionedScalar j0(const dimensionedScalar &ds)
labelList f(nPoints)
dictionary dict
const pointField & pts
std::vector< Triangle > triangles
const label maxDynListLength(viewFactorDict.getOrDefault< label >("maxDynListLength", 1000000))
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299