Loading...
Searching...
No Matches
triSurfaceMesh.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) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2015-2024 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 "triSurfaceMesh.H"
30#include "Random.H"
32#include "edgeHashes.H"
33#include "triSurfaceFields.H"
34#include "Time.H"
35#include "PatchTools.H"
36
37// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38
39namespace Foam
40{
44
46
47
48// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
49
51(
52 const edge& e,
53 EdgeMap<label>& facesPerEdge
54)
55{
56 //label& count = facesPerEdge(e, 0); // lookup or new entry
57 //if (count == 2)
58 //{
59 // return false;
60 //}
61 //++count;
62 //return true;
63
64 auto fnd = facesPerEdge.find(e);
65 if (fnd)
66 {
67 label& count = fnd.val();
68 const int dir = edge::compare(e, fnd.key());
69 if (dir == 1)
70 {
71 // Incorrect order. Mark with special value
72 count = -1;
73 }
74 else if (dir == 0)
75 {
77 << "Incorrect matched edge " << fnd.key()
78 << " to edge " << e << exit(FatalError);
79 return false;
80 }
81 else if (count != -1)
82 {
83 if (count == 2)
84 {
85 return false;
86 }
87 ++count;
88 }
89 }
90 else
91 {
92 facesPerEdge.insert(e, 1);
93 }
94
95 return true;
96}
97
98
100{
101 if (debug)
102 {
103 Pout<< "triSurfaceMesh::isSurfaceClosed:"
104 << " surface:" << searchableSurface::name()
105 << " determining closedness for surface with "
106 << triSurface::size() << " triangles" << endl;
107 }
108
110
111/*
112 // Construct pointFaces. Let's hope surface has compact point
113 // numbering ...
114 labelListList pointFaces;
115 invertManyToMany(pts.size(), *this, pointFaces);
116
117 // Loop over all faces surrounding point. Count edges emanating from point.
118 // Every edge should be used by two faces exactly.
119 // To prevent doing work twice per edge only look at edges to higher
120 // point
121 EdgeMap<label> facesPerEdge;
122 forAll(pointFaces, pointi)
123 {
124 const labelList& pFaces = pointFaces[pointi];
125
126 facesPerEdge.clear();
127 for (const label facei : pFaces)
128 {
129 const triSurface::face_type& f = triSurface::operator[](facei);
130 const label fp = f.find(pointi);
131
132 // Something weird: if I expand the code of addFaceToEdge in both
133 // below instances it gives a segmentation violation on some
134 // surfaces. Compiler (4.3.2) problem?
135
136
137 // Forward edge
138 const label nextPointi = f[f.fcIndex(fp)];
139
140 if (nextPointi > pointi)
141 {
142 bool okFace = addFaceToEdge
143 (
144 edge(pointi, nextPointi),
145 facesPerEdge
146 );
147
148 if (!okFace)
149 {
150 if (debug)
151 {
152 Pout<< "triSurfaceMesh::isSurfaceClosed :"
153 << " surface is open" << endl;
154 }
155 return false;
156 }
157 }
158
159 // Reverse edge
160 const label prevPointi = f[f.rcIndex(fp)];
161
162 if (prevPointi > pointi)
163 {
164 bool okFace = addFaceToEdge
165 (
166 edge(pointi, prevPointi),
167 facesPerEdge
168 );
169
170 if (!okFace)
171 {
172 if (debug)
173 {
174 Pout<< "triSurfaceMesh::isSurfaceClosed :"
175 << " surface is open" << endl;
176 }
177 return false;
178 }
179 }
180 }
181
182 // Check for any edges used only once.
183 forAllConstIters(facesPerEdge, iter)
184 {
185 if (iter.val() == -1)
186 {
187 // Special value for incorrect orientation
188 if (debug)
189 {
190 Pout<< "triSurfaceMesh::isSurfaceClosed :"
191 << " surface is closed but has inconsistent"
192 << " face orientation" << endl;
193 }
194 WarningInFunction << "Surface " << searchableSurface::name()
195 << " is closed but has inconsistent face orientation"
196 << " at edge " << pts[iter.key().first()]
197 << pts[iter.key().second()] << endl;
198 return false;
199 }
200 else if (iter.val() != 2)
201 {
202 if (debug)
203 {
204 Pout<< "triSurfaceMesh::isSurfaceClosed :"
205 << " surface is open" << endl;
206 }
207 return false;
208 }
209 }
210 }
211*/
212
213 const triSurface& ts = *this;
214 EdgeMap<label> facesPerEdge(2*ts.size());
215 for (const auto& f : ts)
216 {
217 forAll(f, fp)
218 {
219 // Count number of occurences of edge between fp and fp+1
220 const bool okFace = addFaceToEdge(f.edge(fp), facesPerEdge);
221
222 if (!okFace)
223 {
224 if (debug)
225 {
226 Pout<< "triSurfaceMesh::isSurfaceClosed :"
227 << " surface:" << searchableSurface::name()
228 << " surface is non-manifold" << endl;
229 }
230 return false;
231 }
232 }
233 }
234
235
236 // Check for any edges used only once.
237 bool haveWarned = false;
238 forAllConstIters(facesPerEdge, iter)
239 {
240 if (iter.val() != 2 && iter.val() != -1)
241 {
242 if (debug)
243 {
244 Pout<< "triSurfaceMesh::isSurfaceClosed :"
245 << " surface:" << searchableSurface::name()
246 << " surface is open" << endl;
247 }
248 return false;
249 }
250 }
251
252 // Check for any edges with inconsistent normal orientation.
253 forAllConstIters(facesPerEdge, iter)
254 {
255 if (iter.val() == -1)
256 {
257 // Special value for incorrect orientation
258 if (!haveWarned)
259 {
261 << "Surface " << searchableSurface::name()
262 << " is closed but has inconsistent face orientation"
263 << nl
264 << " at edge " << pts[iter.key().first()]
265 << pts[iter.key().second()]
266 << nl
267 << " This means it probably cannot be used"
268 << " for inside/outside queries."
269 << " Suppressing further messages." << endl;
270 haveWarned = true;
271 }
272 //- Return open?
273 //return false;
274 }
275 }
276
277 if (debug)
278 {
279 Pout<< "triSurfaceMesh::isSurfaceClosed :"
280 << " surface:" << searchableSurface::name()
281 << " surface is closed" << endl;
282 }
283 return true;
284}
286
287// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
288
290:
291 searchableSurface(io),
292 objectRegistry
293 (
295 (
296 io.name(),
297 io.instance(),
298 io.local(),
299 io.db(),
300 io.readOpt(),
301 io.writeOpt(),
302 false // searchableSurface already registered under name
303 )
304 ),
305 triSurface(s),
306 triSurfaceRegionSearch(static_cast<const triSurface&>(*this)),
307 minQuality_(-1),
308 surfaceClosed_(-1),
310{
312}
313
314
316:
317 // Find instance for triSurfaceMesh
319 // Reused found instance in objectRegistry
321 (
323 (
324 io.name(),
325 searchableSurface::instance(),
326 io.local(),
327 io.db(),
328 io.readOpt(),
329 io.writeOpt(),
330 false // searchableSurface already registered under name
331 )
332 ),
333 triSurface(static_cast<const searchableSurface&>(*this), dictionary::null),
334 triSurfaceRegionSearch(static_cast<const triSurface&>(*this)),
335 minQuality_(-1),
336 surfaceClosed_(-1),
337 outsideVolType_(volumeType::UNKNOWN)
338{
340}
341
342
344(
345 const IOobject& io,
346 const dictionary& dict
347)
348:
350 // Reused found instance in objectRegistry
352 (
354 (
355 io.name(),
356 searchableSurface::instance(),
357 io.local(),
358 io.db(),
359 io.readOpt(),
360 io.writeOpt(),
361 false // searchableSurface already registered under name
362 )
363 ),
364 triSurface(static_cast<const searchableSurface&>(*this), dict),
365 triSurfaceRegionSearch(static_cast<const triSurface&>(*this), dict),
366 minQuality_(-1),
367 surfaceClosed_(-1),
368 outsideVolType_(volumeType::UNKNOWN)
369{
370 // Read with supplied file name instead of objectPath/filePath
371 if (dict.readIfPresent("file", fName_, keyType::LITERAL))
372 {
373 fName_ = triSurface::relativeFilePath
374 (
375 static_cast<const searchableSurface&>(*this),
376 fName_,
377 true
378 );
379 }
380
381 // Report optional scale factor (eg, mm to m),
382 // which was already applied during triSurface construction
383 scalar scaleFactor{0};
384 if (dict.getOrDefault("scale", scaleFactor) && scaleFactor > 0)
385 {
386 Info<< searchableSurface::name()
387 << " : using scale " << scaleFactor << endl;
388 }
389
390 bounds() = boundBox(triSurface::points(), false);
391
392 // Have optional minimum quality for normal calculation
393 if (dict.readIfPresent("minQuality", minQuality_) && minQuality_ > 0)
394 {
395 Info<< searchableSurface::name()
396 << " : ignoring triangles with quality < "
397 << minQuality_ << " for normals calculation." << endl;
399}
400
401
402Foam::triSurfaceMesh::triSurfaceMesh(const IOobject& io, const readAction r)
403:
404 // Find instance for triSurfaceMesh
406 // Reused found instance in objectRegistry
408 (
410 (
411 io.name(),
412 searchableSurface::instance(),
413 io.local(),
414 io.db(),
415 io.readOpt(),
416 io.writeOpt(),
417 false // searchableSurface already registered under name
418 )
419 ),
420 triSurface(),
421 triSurfaceRegionSearch(static_cast<const triSurface&>(*this)),
422 minQuality_(-1),
423 surfaceClosed_(-1),
424 outsideVolType_(volumeType::UNKNOWN)
425{
426 // Check IO flags
427 if (io.isAnyRead())
428 {
429 const bool searchGlobal(r == localOrGlobal || r == masterOnly);
430
431 const fileName actualFile
432 (
433 searchGlobal
434 ? io.globalFilePath(typeName)
435 : io.localFilePath(typeName)
436 );
437
438 // const fileName actualFile
439 // (
440 // io.filePath(searchGlobal, typeName)
441 // );
442
443 if (debug)
444 {
445 Pout<< "triSurfaceMesh(const IOobject& io) :"
446 << " loading surface " << io.objectPath()
447 << " local filePath:" << io.localFilePath(typeName)
448 << " from:" << actualFile << endl;
449 }
450
451 if (searchGlobal && Pstream::parRun())
452 {
453 // Check where surface was found
454 const fileName localFile(io.localFilePath(typeName));
455
456 if
457 (
458 r == masterOnly
459 && ((actualFile.empty() || actualFile != localFile))
460 )
461 {
462 // Found undecomposed surface. Load on master only
463 if (Pstream::master())
464 {
465 triSurface s2(actualFile);
467 }
469 if (debug)
470 {
471 Pout<< "triSurfaceMesh(const IOobject& io) :"
472 << " loaded triangles:" << triSurface::size() << endl;
473 }
474 }
475 else
476 {
477 // Read on all processors
478 triSurface s2(actualFile);
480 if (debug)
481 {
482 Pout<< "triSurfaceMesh(const IOobject& io) :"
483 << " loaded triangles:" << triSurface::size() << endl;
484 }
485 }
486 }
487 else
488 {
489 // Read on all processors
490 triSurface s2(actualFile);
492 if (debug)
493 {
494 Pout<< "triSurfaceMesh(const IOobject& io) :"
495 << " loaded triangles:" << triSurface::size() << endl;
496 }
497 }
498 }
499
501}
502
503
505(
506 const IOobject& io,
507 const dictionary& dict,
508 const readAction r
509)
510:
512 // Reused found instance in objectRegistry
514 (
516 (
517 io.name(),
518 searchableSurface::instance(),
519 io.local(),
520 io.db(),
521 io.readOpt(),
522 io.writeOpt(),
523 false // searchableSurface already registered under name
524 )
525 ),
526 triSurface(),
527 triSurfaceRegionSearch(static_cast<const triSurface&>(*this), dict),
528 minQuality_(-1),
529 surfaceClosed_(-1),
530 outsideVolType_(volumeType::UNKNOWN)
531{
532 if (io.isAnyRead())
533 {
534 // Surface type (optional)
535 const word surfType(dict.getOrDefault<word>("fileType", word::null));
536
537 // Scale factor (optional)
538 const scalar scaleFactor(dict.getOrDefault<scalar>("scale", 0));
539
540 const bool searchGlobal(r == localOrGlobal || r == masterOnly);
541
542 fileName actualFile
543 (
544 searchGlobal
545 ? io.globalFilePath(typeName)
546 : io.localFilePath(typeName)
547 );
548
549 // const fileName actualFile
550 // (
551 // io.filePath(searchGlobal, typeName)
552 // );
553
554 // Reading from supplied file name instead of objectPath/filePath
555 if (dict.readIfPresent("file", fName_, keyType::LITERAL))
556 {
557 fName_ = relativeFilePath
558 (
559 static_cast<const searchableSurface&>(*this),
560 fName_,
561 searchGlobal
562 );
563 actualFile = fName_;
564 }
565
566 if (debug)
567 {
568 Pout<< "triSurfaceMesh(const IOobject& io, const dictionary&) :"
569 << " loading surface " << io.objectPath()
570 << " local filePath:" << io.localFilePath(typeName)
571 << " from:" << actualFile << endl;
572 }
573
574
575 if (searchGlobal && Pstream::parRun())
576 {
577 // Check where surface was found. Bit tricky:
578 // - master will have actualFile (in parent directory)
579 // different from localFilePath (in processor0/)
580 // - slave might have actualFile empty and localFile empty
581
582 const fileName localFile(io.localFilePath(typeName));
583
584 if
585 (
586 r == masterOnly
587 && ((actualFile.empty() || actualFile != localFile))
588 )
589 {
590 // Surface not loaded from processor directories -> undecomposed
591 // surface. Load on master only
592 if (Pstream::master())
593 {
594 triSurface s2(actualFile, surfType, scaleFactor);
596 }
598 if (debug)
599 {
600 Pout<< "triSurfaceMesh(const IOobject& io) :"
601 << " loaded triangles:" << triSurface::size() << endl;
602 }
603 }
604 else
605 {
606 // Read on all processors
607 triSurface s2(actualFile, surfType, scaleFactor);
609 if (debug)
610 {
611 Pout<< "triSurfaceMesh(const IOobject& io) :"
612 << " loaded triangles:" << triSurface::size() << endl;
613 }
614 }
615 }
616 else
617 {
618 // Read on all processors
619 triSurface s2(actualFile, surfType, scaleFactor);
621 if (debug)
622 {
623 Pout<< "triSurfaceMesh(const IOobject& io) :"
624 << " loaded triangles:" << triSurface::size() << endl;
625 }
626 }
627
628
629 // Report optional scale factor (eg, mm to m),
630 // which was already applied during triSurface reading
631 if (scaleFactor > 0)
632 {
634 << " : using scale " << scaleFactor << endl;
635 }
636 }
637
638
639 bounds() = boundBox(triSurface::points(), false);
640
641 // Have optional minimum quality for normal calculation
642 if (dict.readIfPresent("minQuality", minQuality_) && minQuality_ > 0)
643 {
645 << " : ignoring triangles with quality < "
646 << minQuality_ << " for normals calculation." << endl;
647 }
648}
650
651// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
652
654{
656}
657
658
660{
661 // Do not clear closedness status
663 edgeTree_.clear();
665}
667
668// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
669
671{
672 auto tpts = tmp<pointField>::New();
673 auto& pts = tpts.ref();
674
676 {
677 // Can reuse existing values
679 }
680 else
681 {
683
684 // Calculate face centres from a copy to avoid incurring
685 // additional storage
687 (
690 ).faceCentres();
691 }
692
693 return tpts;
694}
695
696
698(
699 pointField& centres,
700 scalarField& radiusSqr
701) const
702{
703 centres = coordinates();
704 radiusSqr.setSize(size());
705 radiusSqr = 0.0;
706
708
709 forAll(*this, facei)
710 {
711 const labelledTri& f = triSurface::operator[](facei);
712 const point& fc = centres[facei];
713 for (const label pointi : f)
714 {
715 radiusSqr[facei] = max(radiusSqr[facei], fc.distSqr(pts[pointi]));
716 }
717 }
718
719 // Add a bit to make sure all points are tested inside
720 radiusSqr += Foam::sqr(SMALL);
721}
722
723
725{
727}
728
729
730bool Foam::triSurfaceMesh::overlaps(const boundBox& bb) const
731{
732 const indexedOctree<treeDataTriSurface>& octree = tree();
733
734 labelList indices = octree.findBox(treeBoundBox(bb));
735
736 return !indices.empty();
737}
738
739
740void Foam::triSurfaceMesh::movePoints(const pointField& newPoints)
741{
742 if (debug)
743 {
744 Pout<< "triSurfaceMesh::movePoints :"
745 << " surface:" << searchableSurface::name()
746 << " moving at time " << objectRegistry::time().timeName()
747 << endl;
748 }
749
750 // Preserve topological point status (surfaceClosed_, outsideVolType_)
751
752 // Update local information (instance, event number)
755
756 const label event = getEvent();
759
760 // Clear additional addressing
762 edgeTree_.clear();
763 triSurface::movePoints(newPoints);
764
765 bounds() = boundBox(triSurface::points(), false);
766 if (debug)
767 {
768 Pout<< "triSurfaceMesh::movePoints: finished moving points" << endl;
769 }
771
772
775{
776 if (!edgeTree_)
777 {
778 if (debug)
779 {
780 Pout<< "triSurfaceMesh::edgeTree :"
781 << " surface:" << searchableSurface::name()
782 << " constructing tree for " << nEdges() - nInternalEdges()
783 << " boundary edges" << endl;
784 }
785
786 // Boundary edges
787 const labelRange bEdges(nInternalEdges(), nBoundaryEdges());
788
790
791 if (bEdges.size())
792 {
793 label nPoints;
795 (
796 *this,
797 bb,
798 nPoints
799 );
800
801 // Random number generator. Bit dodgy since not exactly random ;-)
802 Random rndGen(65431);
803
804 // Slightly extended bb. Slightly off-centred just so on symmetric
805 // geometry there are less face/edge aligned items.
806
807 bb.inflate(rndGen, 1e-4, ROOTVSMALL);
808 }
809
810
811 if (debug)
812 {
813 Pout<< "triSurfaceMesh::edgeTree : "
814 << "calculating edge tree for bb:" << bb << endl;
815 }
816
817 const scalar oldTol =
819
820 edgeTree_.reset
821 (
823 (
824 // Boundary edges
825 treeDataEdge(edges(), localPoints(), bEdges),
826
827 bb, // bb
828 maxTreeDepth(), // maxLevel
829 10, // leafsize
830 3.0 // duplicity
831 )
832 );
833
835
836 if (debug)
837 {
838 Pout<< "triSurfaceMesh::edgeTree :"
839 << " finished constructing tree for "
840 << nEdges() - nInternalEdges()
841 << " boundary edges" << endl;
842 }
843 }
844
845 return *edgeTree_;
846}
847
848
850{
851 if (regions_.empty())
852 {
853 regions_.setSize(patches().size());
854 forAll(regions_, regioni)
855 {
856 regions_[regioni] = patches()[regioni].name();
857 }
858 }
859 return regions_;
860}
861
862
864{
865 if (surfaceClosed_ == -1)
866 {
867 if (isSurfaceClosed())
868 {
869 surfaceClosed_ = 1;
870 }
871 else
872 {
873 surfaceClosed_ = 0;
874 }
875 }
876
877 return surfaceClosed_ == 1;
878}
879
880
882{
883 if (outsideVolType_ == volumeType::UNKNOWN)
884 {
885 // Get point outside bounds()
886 const point outsidePt(bounds().max() + 0.5*bounds().span());
887
888 if (debug)
889 {
890 Pout<< "triSurfaceMesh::outsideVolumeType :"
891 << " surface:" << searchableSurface::name()
892 << " triggering outsidePoint:" << outsidePt
893 << " orientation" << endl;
894 }
895
896 //outsideVolType_ = tree().shapes().getVolumeType(tree(), outsidePt);
897 // Note: do not use tree directly so e.g. distributedTriSurfaceMesh
898 // has opportunity to intercept
899 List<volumeType> outsideVolTypes;
900 getVolumeType(pointField(1, outsidePt), outsideVolTypes);
901 outsideVolType_ = outsideVolTypes[0];
902
903 if (debug)
904 {
905 Pout<< "triSurfaceMesh::outsideVolumeType :"
906 << " finished outsidePoint:" << outsidePt
907 << " orientation:" << volumeType::names[outsideVolType_]
908 << endl;
909 }
910 }
911
913}
914
915
917{
918 if (debug)
919 {
920 Pout<< "triSurfaceMesh::flip :"
921 << " surface:" << searchableSurface::name()
922 << " with current orientation "
923 << volumeType::names[outsideVolType_]
924 << endl;
925 }
926
927 // Don't bother getting nearest etc. Just flip the triangles.
928
929 // triSurface
930 {
931 triSurface& s = *this;
932 s.clearOut();
933 for (auto& tri : s)
934 {
935 tri.flip();
936 }
937 }
938
939 // triSurfaceRegionSearch (if cached volume type)
941
942 // edge tree not relevant
943
944 if (hasVolumeType())
945 {
946 // outsideVolType_
947 if (outsideVolType_ == volumeType::INSIDE)
948 {
949 outsideVolType_ = volumeType::OUTSIDE;
950 }
951 else if (outsideVolType_ == volumeType::OUTSIDE)
952 {
953 outsideVolType_ = volumeType::INSIDE;
954 }
956}
957
958
960(
961 const pointField& samples,
962 const scalarField& nearestDistSqr,
963 List<pointIndexHit>& info
964) const
965{
966 if (debug)
967 {
968 Pout<< "triSurfaceMesh::findNearest :"
969 << " surface:" << searchableSurface::name()
970 << " trying to find nearest for " << samples.size()
971 << " samples with max sphere "
972 << (samples.size() ? Foam::sqrt(max(nearestDistSqr)) : Zero)
973 << endl;
974 }
975 triSurfaceSearch::findNearest(samples, nearestDistSqr, info);
976 if (debug)
977 {
978 Pout<< "triSurfaceMesh::findNearest :"
979 << " finished trying to find nearest for " << samples.size()
980 << " samples" << endl;
982}
983
984
986(
987 const pointField& samples,
988 const scalarField& nearestDistSqr,
989 const labelList& regionIndices,
990 List<pointIndexHit>& info
991) const
992{
993 if (debug)
994 {
995 Pout<< "triSurfaceMesh::findNearest :"
996 << " surface:" << searchableSurface::name()
997 << " trying to find nearest and region for " << samples.size()
998 << " samples with max sphere "
999 << (samples.size() ? Foam::sqrt(max(nearestDistSqr)) : Zero)
1000 << endl;
1001 }
1003 (
1004 samples,
1005 nearestDistSqr,
1006 regionIndices,
1007 info
1008 );
1009 if (debug)
1010 {
1011 Pout<< "triSurfaceMesh::findNearest :"
1012 << " finished trying to find nearest and region for "
1013 << samples.size() << " samples" << endl;
1015}
1016
1017
1019(
1020 const pointField& start,
1021 const pointField& end,
1022 List<pointIndexHit>& info
1023) const
1024{
1025 if (debug)
1026 {
1027 Pout<< "triSurfaceMesh::findLine :"
1028 << " surface:" << searchableSurface::name()
1029 << " intersecting with "
1030 << start.size() << " rays" << endl;
1031 }
1032 triSurfaceSearch::findLine(start, end, info);
1033 if (debug)
1034 {
1035 Pout<< "triSurfaceMesh::findLine :"
1036 << " finished intersecting with "
1037 << start.size() << " rays" << endl;
1039}
1040
1041
1043(
1044 const pointField& start,
1045 const pointField& end,
1046 List<pointIndexHit>& info
1047) const
1048{
1049 if (debug)
1050 {
1051 Pout<< "triSurfaceMesh::findLineAny :"
1052 << " surface:" << searchableSurface::name()
1053 << " intersecting with "
1054 << start.size() << " rays" << endl;
1055 }
1056 triSurfaceSearch::findLineAny(start, end, info);
1057 if (debug)
1058 {
1059 Pout<< "triSurfaceMesh::findLineAny :"
1060 << " finished intersecting with "
1061 << start.size() << " rays" << endl;
1063}
1064
1065
1067(
1068 const pointField& start,
1069 const pointField& end,
1070 List<List<pointIndexHit>>& info
1071) const
1072{
1073 if (debug)
1074 {
1075 Pout<< "triSurfaceMesh::findLineAll :"
1076 << " surface:" << searchableSurface::name()
1077 << " intersecting with "
1078 << start.size() << " rays" << endl;
1079 }
1080 triSurfaceSearch::findLineAll(start, end, info);
1081 if (debug)
1082 {
1083 Pout<< "triSurfaceMesh::findLineAll :"
1084 << " finished intersecting with "
1085 << start.size() << " rays" << endl;
1087}
1088
1089
1091(
1092 const List<pointIndexHit>& info,
1093 labelList& region
1094) const
1095{
1096 if (debug)
1097 {
1098 Pout<< "triSurfaceMesh::getRegion :"
1099 << " surface:" << searchableSurface::name()
1100 << " getting region for "
1101 << info.size() << " triangles" << endl;
1102 }
1103 region.setSize(info.size());
1104 forAll(info, i)
1105 {
1106 if (info[i].hit())
1107 {
1108 region[i] = triSurface::operator[](info[i].index()).region();
1109 }
1110 else
1111 {
1112 region[i] = -1;
1113 }
1114 }
1115 if (debug)
1116 {
1117 Pout<< "triSurfaceMesh::getRegion :"
1118 << " finished getting region for "
1119 << info.size() << " triangles" << endl;
1121}
1122
1123
1125(
1126 const List<pointIndexHit>& info,
1127 vectorField& normal
1128) const
1129{
1130 if (debug)
1131 {
1132 Pout<< "triSurfaceMesh::getNormal :"
1133 << " surface:" << searchableSurface::name()
1134 << " getting normal for "
1135 << info.size() << " triangles" << endl;
1136 }
1137
1138 const triSurface& s = *this;
1139 const pointField& pts = s.points();
1140
1141 normal.setSize(info.size());
1142
1143 if (minQuality_ >= 0)
1144 {
1145 // Make sure we don't use triangles with low quality since
1146 // normal is not reliable.
1147
1148 const labelListList& faceFaces = s.faceFaces();
1149
1150 forAll(info, i)
1151 {
1152 if (info[i].hit())
1153 {
1154 const label facei = info[i].index();
1155 normal[i] = s[facei].unitNormal(pts);
1156
1157 scalar qual = s[facei].tri(pts).quality();
1158
1159 if (qual < minQuality_)
1160 {
1161 // Search neighbouring triangles
1162 const labelList& fFaces = faceFaces[facei];
1163
1164 for (const label nbri : fFaces)
1165 {
1166 scalar nbrQual = s[nbri].tri(pts).quality();
1167 if (nbrQual > qual)
1168 {
1169 qual = nbrQual;
1170 normal[i] = s[nbri].unitNormal(pts);
1171 }
1172 }
1173 }
1174 }
1175 else
1176 {
1177 // Set to what?
1178 normal[i] = Zero;
1179 }
1180 }
1181 }
1182 else
1183 {
1184 forAll(info, i)
1185 {
1186 if (info[i].hit())
1187 {
1188 const label facei = info[i].index();
1189
1190 // Uncached
1191 normal[i] = s[facei].unitNormal(pts);
1192 }
1193 else
1194 {
1195 // Set to what?
1196 normal[i] = Zero;
1197 }
1198 }
1199 }
1200 if (debug)
1201 {
1202 Pout<< "triSurfaceMesh::getNormal :"
1203 << " finished getting normal for "
1204 << info.size() << " triangles" << endl;
1206}
1207
1208
1209void Foam::triSurfaceMesh::setField(const labelList& values)
1210{
1211 auto* fldPtr = getObjectPtr<triSurfaceLabelField>("values");
1212
1213 if (fldPtr)
1214 {
1215 (*fldPtr).field() = values;
1216 }
1217 else
1218 {
1219 fldPtr = new triSurfaceLabelField
1220 (
1221 IOobject
1222 (
1223 "values",
1224 objectRegistry::time().timeName(), // instance
1225 meshSubDir, // local
1226 *this,
1230 ),
1231 *this,
1232 dimless,
1233 labelField(values)
1234 );
1235
1236 // Store field on triMesh
1237 fldPtr->store();
1238 }
1239 if (debug)
1240 {
1241 Pout<< "triSurfaceMesh::setField :"
1242 << " surface:" << searchableSurface::name()
1243 << " finished setting field for "
1244 << values.size() << " triangles" << endl;
1246}
1247
1248
1250(
1251 const List<pointIndexHit>& info,
1252 labelList& values
1253) const
1254{
1255 const auto* fldPtr = getObjectPtr<triSurfaceLabelField>("values");
1256
1257 if (fldPtr)
1258 {
1259 const auto& fld = *fldPtr;
1260
1261 values.setSize(info.size());
1262
1263 forAll(info, i)
1264 {
1265 if (info[i].hit())
1266 {
1267 values[i] = fld[info[i].index()];
1268 }
1269 }
1270 }
1271 if (debug)
1272 {
1273 Pout<< "triSurfaceMesh::setField :"
1274 << " surface:" << searchableSurface::name()
1275 << " finished getting field for "
1276 << info.size() << " triangles" << endl;
1278}
1279
1280
1282(
1283 const pointField& points,
1284 List<volumeType>& volType
1285) const
1286{
1287 const scalar oldTol =
1289
1290 if (debug)
1291 {
1292 Pout<< "triSurfaceMesh::getVolumeType :"
1293 << " surface:" << searchableSurface::name()
1294 << " finding orientation for " << points.size()
1295 << " samples" << endl;
1296 }
1297
1298 volType.setSize(points.size());
1299
1300 forAll(points, pointi)
1301 {
1302 const point& pt = points[pointi];
1303
1304 if (tree().bb().contains(pt))
1305 {
1306 // Use cached volume type per each tree node
1307 volType[pointi] = tree().getVolumeType(pt);
1308 }
1309 else if (hasVolumeType())
1310 {
1311 // Precalculate and cache value for this outside point
1312 if (outsideVolType_ == volumeType::UNKNOWN)
1313 {
1314 outsideVolType_ = tree().shapes().getVolumeType(tree(), pt);
1315 }
1316 volType[pointi] = outsideVolType_;
1317 }
1318 else
1319 {
1320 // Have to calculate directly as outside the octree
1321 volType[pointi] = tree().shapes().getVolumeType(tree(), pt);
1322 }
1323 }
1324
1326
1327 if (debug)
1328 {
1329 Pout<< "triSurfaceMesh::getVolumeType :"
1330 << " finished finding orientation for " << points.size()
1331 << " samples" << endl;
1333}
1334
1335
1337(
1338 IOstreamOption,
1339 const bool writeOnProc
1340) const
1341{
1342 const Time& runTime = searchableSurface::time();
1343 const fileName& instance = searchableSurface::instance();
1344
1345 if
1346 (
1347 instance != runTime.timeName()
1348 && instance != runTime.system()
1349 && instance != runTime.caseSystem()
1350 && instance != runTime.constant()
1351 && instance != runTime.caseConstant()
1352 )
1353 {
1354 const_cast<triSurfaceMesh&>(*this).searchableSurface::instance() =
1355 runTime.timeName();
1356 const_cast<triSurfaceMesh&>(*this).objectRegistry::instance() =
1357 runTime.timeName();
1358 }
1359
1360 fileName fullPath;
1361 if (fName_.size())
1362 {
1363 // Override file name
1364
1365 fullPath = fName_;
1366
1367 fullPath.expand();
1368 if (!fullPath.isAbsolute())
1369 {
1370 // Add directory from regIOobject
1371 fullPath = searchableSurface::objectPath().path()/fullPath;
1372 }
1373 }
1374 else
1375 {
1376 fullPath = searchableSurface::objectPath();
1377 }
1378
1379 if (!mkDir(fullPath.path()))
1380 {
1381 return false;
1382 }
1383
1384 triSurface::write(fullPath);
1385
1386 if (!isFile(fullPath))
1387 {
1388 return false;
1389 }
1390
1391 return true;
1392}
1393
1394
1395// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
if(patchID !=-1)
Info<< nl;Info<< "Write faMesh in vtk format:"<< nl;{ vtk::uindirectPatchWriter writer(aMesh.patch(), fileName(aMesh.time().globalPath()/vtkBaseFileName));writer.writeGeometry();globalIndex procAddr(aMesh.nFaces());labelList cellIDs;if(UPstream::master()) { cellIDs.resize(procAddr.totalSize());for(const labelRange &range :procAddr.ranges()) { auto slice=cellIDs.slice(range);slice=identity(range);} } writer.beginCellData(4);writer.writeProcIDs();writer.write("cellID", cellIDs);writer.write("area", aMesh.S().field());writer.write("normal", aMesh.faceAreaNormals());writer.beginPointData(1);writer.write("normal", aMesh.pointAreaNormals());Info<< " "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edges")));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
Map from edge (expressed as its endpoints) to value. Hashing (and ==) on an edge is symmetric.
Definition edgeHashes.H:59
bool insert(const Key &key, const T &obj)
Copy insert a new entry, not overwriting existing entries.
Definition HashTableI.H:152
iterator find(const Key &key)
Find and return an iterator set at the hashed entry.
Definition HashTableI.H:86
@ REGISTER
Request registration (bool: true).
readOption readOpt() const noexcept
Get the read option.
writeOption writeOpt() const noexcept
Get the write option.
@ NO_READ
Nothing to be read.
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const Time & time() const noexcept
Return Time associated with the objectRegistry.
Definition IOobject.C:456
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
const objectRegistry & db() const noexcept
Return the local objectRegistry.
Definition IOobject.C:450
InfoProxy< IOobject > info() const noexcept
Return info proxy, for printing information to a stream.
Definition IOobject.H:1041
const fileName & local() const noexcept
Read access to local path component.
Definition IOobjectI.H:301
const fileName & instance() const noexcept
Read access to instance path component.
Definition IOobjectI.H:289
fileName objectPath() const
The complete path + object name.
Definition IOobjectI.H:313
A simple container for options an IOstream can normally have.
IntType size() const noexcept
The size of the range.
Definition IntRange.H:208
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
void setSize(label n)
Alias for resize().
Definition List.H:536
static const List< T > & null() noexcept
Return a null List (reference to a nullObject). Behaves like an empty List.
Definition List.H:138
static void calcBounds(const PrimitivePatch< FaceList, PointField > &p, boundBox &bb, label &nPoints)
A list of faces which address into the list of points.
const Field< point_type > & points() const noexcept
Random number generator.
Definition Random.H:56
A non-owning sub-view of a List (allocated or unallocated storage).
Definition SubList.H:61
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition Time.C:714
bool empty() const noexcept
True if List is empty (ie, size() is zero).
Definition UList.H:701
label size() const noexcept
The number of elements in the container.
Definition UList.H:706
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
T & operator[](const label i)
Return element of UList.
Definition UListI.H:363
static bool parRun(const bool on) noexcept
Set as parallel run on/off.
Definition UPstream.H:1669
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition UPstream.H:1714
@ broadcast
broadcast [MPI]
Definition UPstream.H:189
static bool & parRun() noexcept
Test if this a parallel run.
Definition UPstream.H:1681
static const Form zero
scalar distSqr(const Vector< Cmpt > &v2) const
The L2-norm distance squared from another vector. The magSqr() of the difference.
Definition VectorI.H:95
A bounding box defined in terms of min/max extrema points.
Definition boundBox.H:71
void inflate(const scalar factor)
Expand box by factor*mag(span) in all dimensions.
Definition boundBoxI.H:381
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry if present, and assign to T val. FatalIOError if it is found and the number of tokens i...
An edge is a list of two vertex labels. This can correspond to a directed graph edge or an edge on a ...
Definition edge.H:62
static int compare(const edge &a, const edge &b)
Compare edges.
Definition edgeI.H:26
A class for handling file names.
Definition fileName.H:75
static std::string path(const std::string &str)
Return directory path name (part before last /).
Definition fileNameI.H:169
static bool isAbsolute(const std::string &str)
Return true if filename starts with a '/' or '\' or (windows-only) with a filesystem-root.
Definition fileNameI.H:129
static scalar & perturbTol() noexcept
Get the perturbation tolerance.
Non-pointer based hierarchical recursive searching.
@ LITERAL
String literal.
Definition keyType.H:82
A range or interval of labels defined by a start and a size.
Definition labelRange.H:66
A triFace with additional (region) index.
Definition labelledTri.H:56
Registry of regIOobjects.
bool contains(const word &name, const bool recursive=false) const
Does the registry contain the regIOobject object (by name).
const Time & time() const noexcept
Return time registry.
label count(const char *clsName) const
The number of objects of the given class name.
label getEvent() const
Return new event number.
Type * getObjectPtr(const word &name, const bool recursive=false) const
Return non-const pointer to the object of the given Type, using a const-cast to have it behave like a...
label eventNo() const noexcept
Event number at last update.
Base class of (analytical or triangulated) surface. Encapsulates all the search routines....
virtual const boundBox & bounds() const
Return const reference to boundBox.
string & expand(const bool allowEmpty=false)
Inplace expand initial tags, tildes, and all occurrences of environment variables as per stringOps::e...
Definition string.C:166
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
Standard boundBox with extra functionality for use in octree.
Holds data for octree to work on an edges subset.
IOoject and searching on triSurface.
label surfaceClosed_
Is surface closed.
virtual label size() const
Range of local indices that can be returned.
virtual bool writeObject(IOstreamOption streamOpt, const bool writeOnProc) const
Write using stream options.
virtual void flip()
Flip triangles, outsideVolumeType and all cached inside/outside.
wordList regions_
Names of regions.
volumeType outsideVolType_
If surface is closed, what is type of outside points.
virtual ~triSurfaceMesh()
Destructor.
virtual void findLine(const pointField &start, const pointField &end, List< pointIndexHit > &) const
Find first intersection on segment from start to end.
virtual void findLineAll(const pointField &start, const pointField &end, List< List< pointIndexHit > > &) const
Get all intersections in order from start to end.
virtual bool overlaps(const boundBox &bb) const
Does any part of the surface overlap the supplied bound box?
const indexedOctree< treeDataEdge > & edgeTree() const
Demand driven construction of octree for boundary edges.
virtual void findLineAny(const pointField &start, const pointField &end, List< pointIndexHit > &) const
Return any intersection on segment from start to end.
virtual void boundingSpheres(pointField &centres, scalarField &radiusSqr) const
Get bounding spheres (centre and radius squared). Any point.
virtual void getVolumeType(const pointField &points, List< volumeType > &volType) const
Determine type (inside/outside/mixed) for point.
scalar minQuality_
Optional min triangle quality. Triangles below this get.
autoPtr< indexedOctree< treeDataEdge > > edgeTree_
Search tree for boundary edges.
virtual void setField(const labelList &values)
WIP. Store element-wise field.
virtual void getField(const List< pointIndexHit > &, labelList &) const
WIP. From a set of hits (points and.
triSurfaceMesh(const triSurfaceMesh &)=delete
No copy construct.
virtual void getNormal(const List< pointIndexHit > &, vectorField &normal) const
From a set of points and indices get the normal.
fileName fName_
Supplied fileName override.
virtual void getRegion(const List< pointIndexHit > &, labelList &region) const
From a set of points and indices get the region.
virtual void movePoints(const pointField &)
Move points.
virtual const wordList & regions() const
Names of regions.
virtual bool hasVolumeType() const
Whether supports volume type (below) - i.e. whether is closed.
bool isSurfaceClosed() const
Check whether surface is closed without calculating any permanent.
virtual void findNearest(const pointField &sample, const scalarField &nearestDistSqr, List< pointIndexHit > &) const
virtual tmp< pointField > coordinates() const
Get representative set of element coordinates.
static word meshSubDir
Return the mesh sub-directory name (usually "triSurface").
virtual volumeType outsideVolumeType() const
If surface is closed, what is type of points outside bounds.
static bool addFaceToEdge(const edge &, EdgeMap< label > &)
Helper function for isSurfaceClosed.
void clearOut()
Clear storage.
virtual tmp< pointField > points() const
Get the points that define the surface.
Helper class to search on triSurface. Creates an octree for each region of the surface and only searc...
void findNearest(const pointField &samples, const scalarField &nearestDistSqr, const labelList &regionIndices, List< pointIndexHit > &info) const
Find the nearest point on the surface out of the regions.
void findLineAll(const pointField &start, const pointField &end, List< List< pointIndexHit > > &info) const
Calculate all intersections from start to end.
scalar tolerance() const
Return tolerance to use in searches.
void findLine(const pointField &start, const pointField &end, List< pointIndexHit > &info) const
void findLineAny(const pointField &start, const pointField &end, List< pointIndexHit > &info) const
label maxTreeDepth() const
Return max tree depth of octree.
void findNearest(const pointField &samples, const scalarField &nearestDistSqr, List< pointIndexHit > &info) const
Triangulated surface description with patch information.
Definition triSurface.H:74
triSurface()
Default construct.
Definition triSurface.C:426
void transfer(triSurface &surf)
Alter contents by transferring (triangles, points) components.
Definition triSurface.C:955
void write(Ostream &os) const
Write to Ostream in simple OpenFOAM format.
const geometricSurfacePatchList & patches() const noexcept
Definition triSurface.H:509
virtual void movePoints(const pointField &pts)
Move points.
Definition triSurface.C:606
An enumeration wrapper for classification of a location as being inside/outside of a volume.
Definition volumeType.H:56
static const Enum< volumeType::type > names
Names for the classification enumeration.
Definition volumeType.H:75
@ OUTSIDE
A location outside the volume.
Definition volumeType.H:66
@ UNKNOWN
Unknown state.
Definition volumeType.H:64
@ INSIDE
A location inside the volume.
Definition volumeType.H:65
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
bool local
Definition EEqn.H:20
const polyBoundaryMesh & patches
engineTime & runTime
PtrList< coordinateSystem > coordinates(solidRegions.size())
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
const auto & io
const pointField & points
label nPoints
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
word timeName
Definition getTimeIndex.H:3
#define WarningInFunction
Report a warning using Foam::Warning.
List< T > values(const HashTable< T, Key, Hash > &tbl, const bool doSort=false)
List of values from HashTable, optionally sorted.
Definition HashOps.H:164
Namespace for bounding specifications. At the moment, mostly for tables.
Namespace for handling debugging switches.
Definition debug.C:45
Namespace for OpenFOAM.
List< word > wordList
List of word.
Definition fileName.H:60
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
DimensionedField< label, triSurfaceGeoMesh > triSurfaceLabelField
const dimensionSet dimless
Dimensionless.
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
List< label > labelList
A List of labels.
Definition List.H:62
dimensionedSymmTensor sqr(const dimensionedVector &dv)
bool mkDir(const fileName &pathName, mode_t mode=0777)
Make a directory and return an error if it could not be created.
Definition POSIX.C:616
messageStream Info
Information stream (stdout output on master, null elsewhere).
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensionedScalar sqrt(const dimensionedScalar &ds)
Field< vector > vectorField
Specialisation of Field<T> for vector.
Field< label > labelField
Specialisation of Field<T> for label.
Definition labelField.H:48
vector point
Point is a vector.
Definition point.H:37
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
bool isFile(const fileName &name, const bool checkGzip=true, const bool followLink=true)
Does the name exist as a FILE in the file system?
Definition POSIX.C:879
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.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
labelList f(nPoints)
dictionary dict
const pointField & pts
Tree tree(triangles.begin(), triangles.end())
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
Definition stdFoam.H:235
Fields for triSurface.
scalarField samples(nIntervals, Zero)
Random rndGen