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-2020,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 "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{
43}
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 << " determining closedness for surface with "
105 << triSurface::size() << " triangles" << endl;
106 }
107
108 const pointField& pts = triSurface::points();
109
110/*
111 // Construct pointFaces. Let's hope surface has compact point
112 // numbering ...
113 labelListList pointFaces;
114 invertManyToMany(pts.size(), *this, pointFaces);
115
116 // Loop over all faces surrounding point. Count edges emanating from point.
117 // Every edge should be used by two faces exactly.
118 // To prevent doing work twice per edge only look at edges to higher
119 // point
120 EdgeMap<label> facesPerEdge(128);
121 forAll(pointFaces, pointi)
122 {
123 const labelList& pFaces = pointFaces[pointi];
124
125 facesPerEdge.clear();
126 for (const label facei : pFaces)
127 {
128 const triSurface::face_type& f = triSurface::operator[](facei);
129 const label fp = f.find(pointi);
130
131 // Something weird: if I expand the code of addFaceToEdge in both
132 // below instances it gives a segmentation violation on some
133 // surfaces. Compiler (4.3.2) problem?
134
135
136 // Forward edge
137 const label nextPointi = f[f.fcIndex(fp)];
138
139 if (nextPointi > pointi)
140 {
141 bool okFace = addFaceToEdge
142 (
143 edge(pointi, nextPointi),
144 facesPerEdge
145 );
146
147 if (!okFace)
148 {
149 if (debug)
150 {
151 Pout<< "triSurfaceMesh::isSurfaceClosed :"
152 << " surface is open" << endl;
153 }
154 return false;
155 }
156 }
157
158 // Reverse edge
159 const label prevPointi = f[f.rcIndex(fp)];
160
161 if (prevPointi > pointi)
162 {
163 bool okFace = addFaceToEdge
164 (
165 edge(pointi, prevPointi),
166 facesPerEdge
167 );
168
169 if (!okFace)
170 {
171 if (debug)
172 {
173 Pout<< "triSurfaceMesh::isSurfaceClosed :"
174 << " surface is open" << endl;
175 }
176 return false;
177 }
178 }
179 }
180
181 // Check for any edges used only once.
182 forAllConstIters(facesPerEdge, iter)
183 {
184 if (iter.val() == -1)
185 {
186 // Special value for incorrect orientation
187 if (debug)
188 {
189 Pout<< "triSurfaceMesh::isSurfaceClosed :"
190 << " surface is closed but has inconsistent"
191 << " face orientation" << endl;
192 }
193 WarningInFunction << "Surface " << searchableSurface::name()
194 << " is closed but has inconsistent face orientation"
195 << " at edge " << pts[iter.key().first()]
196 << pts[iter.key().second()] << endl;
197 return false;
198 }
199 else if (iter.val() != 2)
200 {
201 if (debug)
202 {
203 Pout<< "triSurfaceMesh::isSurfaceClosed :"
204 << " surface is open" << endl;
205 }
206 return false;
207 }
208 }
209 }
210*/
211
212 const triSurface& ts = *this;
213 EdgeMap<label> facesPerEdge(2*ts.size());
214 for (const auto& f : ts)
215 {
216 forAll(f, fp)
217 {
218 // Count number of occurences of edge between fp and fp+1
219 const bool okFace = addFaceToEdge(f.edge(fp), facesPerEdge);
220
221 if (!okFace)
222 {
223 if (debug)
224 {
225 Pout<< "triSurfaceMesh::isSurfaceClosed :"
226 << " surface is non-manifold" << endl;
227 }
228 return false;
229 }
230 }
231 }
232
233
234 // Check for any edges used only once.
235 bool haveWarned = false;
236 forAllConstIters(facesPerEdge, iter)
237 {
238 if (iter.val() != 2 && iter.val() != -1)
239 {
240 if (debug)
241 {
242 Pout<< "triSurfaceMesh::isSurfaceClosed :"
243 << " surface is open" << endl;
244 }
245 return false;
246 }
247 }
248
249 // Check for any edges with inconsistent normal orientation.
250 forAllConstIters(facesPerEdge, iter)
251 {
252 if (iter.val() == -1)
253 {
254 // Special value for incorrect orientation
255 if (!haveWarned)
256 {
258 << "Surface " << searchableSurface::name()
259 << " is closed but has inconsistent face orientation"
260 << nl
261 << " at edge " << pts[iter.key().first()]
262 << pts[iter.key().second()]
263 << nl
264 << " This means it probably cannot be used"
265 << " for inside/outside queries."
266 << " Suppressing further messages." << endl;
267 haveWarned = true;
268 }
269 //- Return open?
270 //return false;
271 }
272 }
273
274 if (debug)
275 {
276 Pout<< "triSurfaceMesh::isSurfaceClosed :"
277 << " surface is closed" << endl;
278 }
279 return true;
280}
281
282
283// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
284
286:
289 (
291 (
292 io.name(),
293 io.instance(),
294 io.local(),
295 io.db(),
296 io.readOpt(),
297 io.writeOpt(),
298 false // searchableSurface already registered under name
299 )
300 ),
301 triSurface(s),
302 triSurfaceRegionSearch(static_cast<const triSurface&>(*this)),
303 minQuality_(-1),
304 surfaceClosed_(-1),
305 outsideVolType_(volumeType::UNKNOWN)
306{
307 const pointField& pts = triSurface::points();
308
309 bounds() = boundBox(pts, false);
310}
311
312
314:
315 // Find instance for triSurfaceMesh
317 // Reused found instance in objectRegistry
319 (
321 (
322 io.name(),
323 searchableSurface::instance(),
324 io.local(),
325 io.db(),
326 io.readOpt(),
327 io.writeOpt(),
328 false // searchableSurface already registered under name
329 )
330 ),
331 triSurface(static_cast<const searchableSurface&>(*this), dictionary::null),
332 triSurfaceRegionSearch(static_cast<const triSurface&>(*this)),
333 minQuality_(-1),
334 surfaceClosed_(-1),
335 outsideVolType_(volumeType::UNKNOWN)
336{
337 const pointField& pts = triSurface::points();
338
339 bounds() = boundBox(pts, false);
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
372 {
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 {
387 << " : using scale " << scaleFactor << endl;
388 }
389
390 const pointField& pts = triSurface::points();
391
392 bounds() = boundBox(pts, false);
393
394 // Have optional minimum quality for normal calculation
395 if (dict.readIfPresent("minQuality", minQuality_) && minQuality_ > 0)
396 {
398 << " : ignoring triangles with quality < "
399 << minQuality_ << " for normals calculation." << endl;
400 }
401}
402
403
405:
406 // Find instance for triSurfaceMesh
408 // Reused found instance in objectRegistry
410 (
412 (
413 io.name(),
414 searchableSurface::instance(),
415 io.local(),
416 io.db(),
417 io.readOpt(),
418 io.writeOpt(),
419 false // searchableSurface already registered under name
420 )
421 ),
422 triSurface(),
423 triSurfaceRegionSearch(static_cast<const triSurface&>(*this)),
424 minQuality_(-1),
425 surfaceClosed_(-1),
426 outsideVolType_(volumeType::UNKNOWN)
427{
428 // Check IO flags
430 {
431 const bool searchGlobal(r == localOrGlobal || r == masterOnly);
432
433 const fileName actualFile
434 (
435 searchGlobal
436 ? io.globalFilePath(typeName)
437 : io.localFilePath(typeName)
438 );
439
440 if (debug)
441 {
442 Pout<< "triSurfaceMesh(const IOobject& io) :"
443 << " loading surface " << io.objectPath()
444 << " local filePath:" << io.localFilePath(typeName)
445 << " from:" << actualFile << endl;
446 }
447
448 if (searchGlobal && Pstream::parRun())
449 {
450 // Check where surface was found
451 const fileName localFile(io.localFilePath(typeName));
452
453 if (r == masterOnly && (actualFile != localFile))
454 {
455 // Found undecomposed surface. Load on master only
456 if (Pstream::master())
457 {
458 triSurface s2(actualFile);
460 }
462 if (debug)
463 {
464 Pout<< "triSurfaceMesh(const IOobject& io) :"
465 << " loaded triangles:" << triSurface::size() << endl;
466 }
467 }
468 else
469 {
470 // Read on all processors
471 triSurface s2(actualFile);
473 if (debug)
474 {
475 Pout<< "triSurfaceMesh(const IOobject& io) :"
476 << " loaded triangles:" << triSurface::size() << endl;
477 }
478 }
479 }
480 else
481 {
482 // Read on all processors
483 triSurface s2(actualFile);
485 if (debug)
486 {
487 Pout<< "triSurfaceMesh(const IOobject& io) :"
488 << " loaded triangles:" << triSurface::size() << endl;
489 }
490 }
491 }
492
493 const pointField& pts = triSurface::points();
494 bounds() = boundBox(pts, false);
495}
496
497
499(
500 const IOobject& io,
501 const dictionary& dict,
502 const readAction r
503)
504:
506 // Reused found instance in objectRegistry
508 (
510 (
511 io.name(),
512 searchableSurface::instance(),
513 io.local(),
514 io.db(),
515 io.readOpt(),
516 io.writeOpt(),
517 false // searchableSurface already registered under name
518 )
519 ),
520 triSurface(),
521 triSurfaceRegionSearch(static_cast<const triSurface&>(*this), dict),
522 minQuality_(-1),
523 surfaceClosed_(-1),
524 outsideVolType_(volumeType::UNKNOWN)
525{
527 {
528 // Surface type (optional)
529 const word surfType(dict.getOrDefault<word>("fileType", word::null));
530
531 // Scale factor (optional)
532 const scalar scaleFactor(dict.getOrDefault<scalar>("scale", 0));
533
534 const bool searchGlobal(r == localOrGlobal || r == masterOnly);
535
536 fileName actualFile
537 (
538 searchGlobal
539 ? io.globalFilePath(typeName)
540 : io.localFilePath(typeName)
541 );
542
543 // Reading from supplied file name instead of objectPath/filePath
545 {
547 (
548 static_cast<const searchableSurface&>(*this),
549 fName_,
550 searchGlobal
551 );
552 actualFile = fName_;
553 }
554
555 if (debug)
556 {
557 Pout<< "triSurfaceMesh(const IOobject& io, const dictionary&) :"
558 << " loading surface " << io.objectPath()
559 << " local filePath:" << io.localFilePath(typeName)
560 << " from:" << actualFile << endl;
561 }
562
563
564 if (searchGlobal && Pstream::parRun())
565 {
566 // Check where surface was found
567 const fileName localFile(io.localFilePath(typeName));
568
569 if (r == masterOnly && (actualFile != localFile))
570 {
571 // Surface not loaded from processor directories -> undecomposed
572 // surface. Load on master only
573 if (Pstream::master())
574 {
575 triSurface s2(actualFile, surfType, scaleFactor);
577 }
579 if (debug)
580 {
581 Pout<< "triSurfaceMesh(const IOobject& io) :"
582 << " loaded triangles:" << triSurface::size() << endl;
583 }
584 }
585 else
586 {
587 // Read on all processors
588 triSurface s2(actualFile, surfType, scaleFactor);
590 if (debug)
591 {
592 Pout<< "triSurfaceMesh(const IOobject& io) :"
593 << " loaded triangles:" << triSurface::size() << endl;
594 }
595 }
596 }
597 else
598 {
599 // Read on all processors
600 triSurface s2(actualFile, surfType, scaleFactor);
602 if (debug)
603 {
604 Pout<< "triSurfaceMesh(const IOobject& io) :"
605 << " loaded triangles:" << triSurface::size() << endl;
606 }
607 }
608
609
610 // Report optional scale factor (eg, mm to m),
611 // which was already applied during triSurface reading
612 if (scaleFactor > 0)
613 {
615 << " : using scale " << scaleFactor << endl;
616 }
617 }
618
619
620 const pointField& pts = triSurface::points();
621 bounds() = boundBox(pts, false);
622
623 // Have optional minimum quality for normal calculation
624 if (dict.readIfPresent("minQuality", minQuality_) && minQuality_ > 0)
625 {
627 << " : ignoring triangles with quality < "
628 << minQuality_ << " for normals calculation." << endl;
629 }
630}
631
632
633// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
634
636{
637 clearOut();
638}
639
640
642{
643 // Do not clear closedness status
645 edgeTree_.clear();
647}
648
649
650// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
651
653{
654 auto tpts = tmp<pointField>::New();
655 auto& pts = tpts.ref();
656
658 {
659 // Can reuse existing values
661 }
662 else
663 {
665
666 // Calculate face centres from a copy to avoid incurring
667 // additional storage
669 (
672 ).faceCentres();
673 }
674
675 return tpts;
676}
677
678
680(
681 pointField& centres,
682 scalarField& radiusSqr
683) const
684{
685 centres = coordinates();
686 radiusSqr.setSize(size());
687 radiusSqr = 0.0;
688
689 const pointField& pts = triSurface::points();
690
691 forAll(*this, facei)
692 {
693 const labelledTri& f = triSurface::operator[](facei);
694 const point& fc = centres[facei];
695 for (const label pointi : f)
696 {
697 const point& pt = pts[pointi];
698 radiusSqr[facei] = max(radiusSqr[facei], Foam::magSqr(fc-pt));
699 }
700 }
701
702 // Add a bit to make sure all points are tested inside
703 radiusSqr += Foam::sqr(SMALL);
704}
705
706
708{
709 return triSurface::points();
710}
711
712
714{
715 const indexedOctree<treeDataTriSurface>& octree = tree();
716
717 labelList indices = octree.findBox(treeBoundBox(bb));
718
719 return !indices.empty();
720}
721
722
724{
725 if (debug)
726 {
727 Pout<< "triSurfaceMesh::movePoints :"
728 << " moving at time " << objectRegistry::time().timeName()
729 << endl;
730 }
731
732 // Preserve topological point status (surfaceClosed_, outsideVolType_)
733
734 // Update local information (instance, event number)
737
738 const label event = getEvent();
741
742 // Clear additional addressing
744 edgeTree_.clear();
745 triSurface::movePoints(newPoints);
746
747 bounds() = boundBox(triSurface::points(), false);
748 if (debug)
749 {
750 Pout<< "triSurfaceMesh::movePoints: finished moving points" << endl;
751 }
752}
753
754
757{
758 if (!edgeTree_)
759 {
760 if (debug)
761 {
762 Pout<< "triSurfaceMesh::edgeTree :"
763 << " constructing tree for " << nEdges() - nInternalEdges()
764 << " boundary edges" << endl;
765 }
766
767 // Boundary edges
768 labelList bEdges
769 (
770 identity(nEdges() - nInternalEdges(), nInternalEdges())
771 );
772
774
775 if (bEdges.size())
776 {
777 label nPoints;
778 PatchTools::calcBounds
779 (
780 *this,
781 bb,
782 nPoints
783 );
784
785 // Random number generator. Bit dodgy since not exactly random ;-)
786 Random rndGen(65431);
787
788 // Slightly extended bb. Slightly off-centred just so on symmetric
789 // geometry there are less face/edge aligned items.
790
791 bb = bb.extend(rndGen, 1e-4);
792 bb.min() -= point::uniform(ROOTVSMALL);
793 bb.max() += point::uniform(ROOTVSMALL);
794 }
795
796
797 if (debug)
798 {
799 Pout<< "triSurfaceMesh::edgeTree : "
800 << "calculating edge tree for bb:" << bb << endl;
801 }
802
805
806 edgeTree_.reset
807 (
809 (
811 (
812 false, // cachebb
813 edges(), // edges
814 localPoints(), // points
815 bEdges // selected edges
816 ),
817 bb, // bb
818 maxTreeDepth(), // maxLevel
819 10, // leafsize
820 3.0 // duplicity
821 )
822 );
823
825
826 if (debug)
827 {
828 Pout<< "triSurfaceMesh::edgeTree :"
829 << " finished constructing tree for "
830 << nEdges() - nInternalEdges()
831 << " boundary edges" << endl;
832 }
833 }
834
835 return *edgeTree_;
836}
837
838
840{
841 if (regions_.empty())
842 {
843 regions_.setSize(patches().size());
844 forAll(regions_, regioni)
845 {
846 regions_[regioni] = patches()[regioni].name();
847 }
848 }
849 return regions_;
850}
851
852
854{
855 if (surfaceClosed_ == -1)
856 {
857 if (isSurfaceClosed())
858 {
859 surfaceClosed_ = 1;
860 }
861 else
862 {
863 surfaceClosed_ = 0;
864 }
865 }
866
867 return surfaceClosed_ == 1;
868}
869
870
872{
873 if (outsideVolType_ == volumeType::UNKNOWN)
874 {
875 // Get point outside bounds()
876 const point outsidePt(bounds().max() + 0.5*bounds().span());
877
878 if (debug)
879 {
880 Pout<< "triSurfaceMesh::outsideVolumeType :"
881 << " triggering outsidePoint:" << outsidePt
882 << " orientation" << endl;
883 }
884
885 //outsideVolType_ = tree().shapes().getVolumeType(tree(), outsidePt);
886 // Note: do not use tree directly so e.g. distributedTriSurfaceMesh
887 // has opportunity to intercept
888 List<volumeType> outsideVolTypes;
889 getVolumeType(pointField(1, outsidePt), outsideVolTypes);
890 outsideVolType_ = outsideVolTypes[0];
891
892 if (debug)
893 {
894 Pout<< "triSurfaceMesh::outsideVolumeType :"
895 << " finished outsidePoint:" << outsidePt
896 << " orientation:" << volumeType::names[outsideVolType_]
897 << endl;
898 }
899 }
900
901 return outsideVolType_;
902}
903
904
906(
907 const pointField& samples,
908 const scalarField& nearestDistSqr,
910) const
911{
912 if (debug)
913 {
914 Pout<< "triSurfaceMesh::findNearest :"
915 << " trying to find nearest for " << samples.size()
916 << " samples with max sphere "
917 << (samples.size() ? Foam::sqrt(max(nearestDistSqr)) : Zero)
918 << endl;
919 }
920 triSurfaceSearch::findNearest(samples, nearestDistSqr, info);
921 if (debug)
922 {
923 Pout<< "triSurfaceMesh::findNearest :"
924 << " finished trying to find nearest for " << samples.size()
925 << " samples" << endl;
926 }
927}
928
929
931(
932 const pointField& samples,
933 const scalarField& nearestDistSqr,
934 const labelList& regionIndices,
936) const
937{
938 if (debug)
939 {
940 Pout<< "triSurfaceMesh::findNearest :"
941 << " trying to find nearest and region for " << samples.size()
942 << " samples with max sphere "
943 << (samples.size() ? Foam::sqrt(max(nearestDistSqr)) : Zero)
944 << endl;
945 }
947 (
948 samples,
949 nearestDistSqr,
950 regionIndices,
951 info
952 );
953 if (debug)
954 {
955 Pout<< "triSurfaceMesh::findNearest :"
956 << " finished trying to find nearest and region for "
957 << samples.size() << " samples" << endl;
958 }
959}
960
961
963(
964 const pointField& start,
965 const pointField& end,
967) const
968{
969 if (debug)
970 {
971 Pout<< "triSurfaceMesh::findLine :"
972 << " intersecting with "
973 << start.size() << " rays" << endl;
974 }
975 triSurfaceSearch::findLine(start, end, info);
976 if (debug)
977 {
978 Pout<< "triSurfaceMesh::findLine :"
979 << " finished intersecting with "
980 << start.size() << " rays" << endl;
981 }
982}
983
984
986(
987 const pointField& start,
988 const pointField& end,
990) const
991{
992 if (debug)
993 {
994 Pout<< "triSurfaceMesh::findLineAny :"
995 << " intersecting with "
996 << start.size() << " rays" << endl;
997 }
998 triSurfaceSearch::findLineAny(start, end, info);
999 if (debug)
1000 {
1001 Pout<< "triSurfaceMesh::findLineAny :"
1002 << " finished intersecting with "
1003 << start.size() << " rays" << endl;
1004 }
1005}
1006
1007
1009(
1010 const pointField& start,
1011 const pointField& end,
1013) const
1014{
1015 if (debug)
1016 {
1017 Pout<< "triSurfaceMesh::findLineAll :"
1018 << " intersecting with "
1019 << start.size() << " rays" << endl;
1020 }
1021 triSurfaceSearch::findLineAll(start, end, info);
1022 if (debug)
1023 {
1024 Pout<< "triSurfaceMesh::findLineAll :"
1025 << " finished intersecting with "
1026 << start.size() << " rays" << endl;
1027 }
1028}
1029
1030
1032(
1033 const List<pointIndexHit>& info,
1034 labelList& region
1035) const
1036{
1037 if (debug)
1038 {
1039 Pout<< "triSurfaceMesh::getRegion :"
1040 << " getting region for "
1041 << info.size() << " triangles" << endl;
1042 }
1043 region.setSize(info.size());
1044 forAll(info, i)
1045 {
1046 if (info[i].hit())
1047 {
1048 region[i] = triSurface::operator[](info[i].index()).region();
1049 }
1050 else
1051 {
1052 region[i] = -1;
1053 }
1054 }
1055 if (debug)
1056 {
1057 Pout<< "triSurfaceMesh::getRegion :"
1058 << " finished getting region for "
1059 << info.size() << " triangles" << endl;
1060 }
1061}
1062
1063
1065(
1066 const List<pointIndexHit>& info,
1067 vectorField& normal
1068) const
1069{
1070 if (debug)
1071 {
1072 Pout<< "triSurfaceMesh::getNormal :"
1073 << " getting normal for "
1074 << info.size() << " triangles" << endl;
1075 }
1076
1077 const triSurface& s = *this;
1078 const pointField& pts = s.points();
1079
1080 normal.setSize(info.size());
1081
1082 if (minQuality_ >= 0)
1083 {
1084 // Make sure we don't use triangles with low quality since
1085 // normal is not reliable.
1086
1087 const labelListList& faceFaces = s.faceFaces();
1088
1089 forAll(info, i)
1090 {
1091 if (info[i].hit())
1092 {
1093 const label facei = info[i].index();
1094 normal[i] = s[facei].unitNormal(pts);
1095
1096 scalar qual = s[facei].tri(pts).quality();
1097
1098 if (qual < minQuality_)
1099 {
1100 // Search neighbouring triangles
1101 const labelList& fFaces = faceFaces[facei];
1102
1103 for (const label nbri : fFaces)
1104 {
1105 scalar nbrQual = s[nbri].tri(pts).quality();
1106 if (nbrQual > qual)
1107 {
1108 qual = nbrQual;
1109 normal[i] = s[nbri].unitNormal(pts);
1110 }
1111 }
1112 }
1113 }
1114 else
1115 {
1116 // Set to what?
1117 normal[i] = Zero;
1118 }
1119 }
1120 }
1121 else
1122 {
1123 forAll(info, i)
1124 {
1125 if (info[i].hit())
1126 {
1127 const label facei = info[i].index();
1128
1129 // Uncached
1130 normal[i] = s[facei].unitNormal(pts);
1131 }
1132 else
1133 {
1134 // Set to what?
1135 normal[i] = Zero;
1136 }
1137 }
1138 }
1139 if (debug)
1140 {
1141 Pout<< "triSurfaceMesh::getNormal :"
1142 << " finished getting normal for "
1143 << info.size() << " triangles" << endl;
1144 }
1145}
1146
1147
1149{
1150 auto* fldPtr = getObjectPtr<triSurfaceLabelField>("values");
1151
1152 if (fldPtr)
1153 {
1154 (*fldPtr).field() = values;
1155 }
1156 else
1157 {
1158 fldPtr = new triSurfaceLabelField
1159 (
1160 IOobject
1161 (
1162 "values",
1163 objectRegistry::time().timeName(), // instance
1164 meshSubDir, // local
1165 *this,
1168 ),
1169 *this,
1170 dimless,
1171 labelField(values)
1172 );
1173
1174 // Store field on triMesh
1175 fldPtr->store();
1176 }
1177 if (debug)
1178 {
1179 Pout<< "triSurfaceMesh::setField :"
1180 << " finished setting field for "
1181 << values.size() << " triangles" << endl;
1182 }
1183}
1184
1185
1187(
1188 const List<pointIndexHit>& info,
1189 labelList& values
1190) const
1191{
1192 const auto* fldPtr = getObjectPtr<triSurfaceLabelField>("values");
1193
1194 if (fldPtr)
1195 {
1196 const auto& fld = *fldPtr;
1197
1198 values.setSize(info.size());
1199
1200 forAll(info, i)
1201 {
1202 if (info[i].hit())
1203 {
1204 values[i] = fld[info[i].index()];
1205 }
1206 }
1207 }
1208 if (debug)
1209 {
1210 Pout<< "triSurfaceMesh::setField :"
1211 << " finished getting field for "
1212 << info.size() << " triangles" << endl;
1213 }
1214}
1215
1216
1218(
1219 const pointField& points,
1220 List<volumeType>& volType
1221) const
1222{
1223 const scalar oldTol = indexedOctree<treeDataTriSurface>::perturbTol();
1225
1226 if (debug)
1227 {
1228 Pout<< "triSurfaceMesh::getVolumeType :"
1229 << " finding orientation for " << points.size()
1230 << " samples" << endl;
1231 }
1232
1233 volType.setSize(points.size());
1234
1235 forAll(points, pointi)
1236 {
1237 const point& pt = points[pointi];
1238
1239 if (tree().bb().contains(pt))
1240 {
1241 // Use cached volume type per each tree node
1242 volType[pointi] = tree().getVolumeType(pt);
1243 }
1244 else if (hasVolumeType())
1245 {
1246 // Precalculate and cache value for this outside point
1247 if (outsideVolType_ == volumeType::UNKNOWN)
1248 {
1249 outsideVolType_ = tree().shapes().getVolumeType(tree(), pt);
1250 }
1251 volType[pointi] = outsideVolType_;
1252 }
1253 else
1254 {
1255 // Have to calculate directly as outside the octree
1256 volType[pointi] = tree().shapes().getVolumeType(tree(), pt);
1257 }
1258 }
1259
1261 if (debug)
1262 {
1263 Pout<< "triSurfaceMesh::getVolumeType :"
1264 << " finished finding orientation for " << points.size()
1265 << " samples" << endl;
1266 }
1267}
1268
1269
1271(
1273 const bool valid
1274) const
1275{
1277 const fileName& instance = searchableSurface::instance();
1278
1279 if
1280 (
1281 instance != runTime.timeName()
1282 && instance != runTime.system()
1283 && instance != runTime.caseSystem()
1284 && instance != runTime.constant()
1285 && instance != runTime.caseConstant()
1286 )
1287 {
1288 const_cast<triSurfaceMesh&>(*this).searchableSurface::instance() =
1289 runTime.timeName();
1290 const_cast<triSurfaceMesh&>(*this).objectRegistry::instance() =
1291 runTime.timeName();
1292 }
1293
1294 fileName fullPath;
1295 if (fName_.size())
1296 {
1297 // Override file name
1298
1299 fullPath = fName_;
1300
1301 fullPath.expand();
1302 if (!fullPath.isAbsolute())
1303 {
1304 // Add directory from regIOobject
1305 fullPath = searchableSurface::objectPath().path()/fullPath;
1306 }
1307 }
1308 else
1309 {
1310 fullPath = searchableSurface::objectPath();
1311 }
1312
1313 if (!mkDir(fullPath.path()))
1314 {
1315 return false;
1316 }
1317
1318 triSurface::write(fullPath);
1319
1320 if (!isFile(fullPath))
1321 {
1322 return false;
1323 }
1324
1325 return true;
1326}
1327
1328
1329// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Info<< nl<< "Wrote faMesh in vtk format: "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.mesh().time().globalPath()/"finiteArea-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. For easier forward declaration it is currently i...
Definition: EdgeMap.H:54
bool insert(const Key &key, const T &obj)
Copy insert a new entry, not overwriting existing entries.
Definition: HashTableI.H:180
iterator find(const Key &key)
Find and return an iterator set at the hashed entry.
Definition: HashTableI.H:114
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:65
readOption readOpt() const noexcept
The read option.
Definition: IOobjectI.H:164
fileName globalFilePath(const word &typeName, const bool search=true) const
Helper for filePath that searches up if in parallel.
Definition: IOobject.C:593
fileName localFilePath(const word &typeName, const bool search=true) const
Helper for filePath that searches locally.
Definition: IOobject.C:582
const fileName & instance() const noexcept
Read access to instance path component.
Definition: IOobjectI.H:196
fileName objectPath() const
The complete path + object name.
Definition: IOobjectI.H:214
The IOstreamOption is a simple container for options an IOstream can normally have.
void setSize(const label n)
Alias for resize()
Definition: List.H:218
int overlaps
Flag to control which overlap calculations are performed.
Definition: PDRparams.H:97
A list of faces which address into the list of points.
const Field< point_type > & faceCentres() const
Return face centres for patch.
static void scatter(const List< commsStruct > &comms, T &value, const int tag, const label comm)
Broadcast data: Distribute without modification.
Random number generator.
Definition: Random.H:60
A List obtained as a section of another List.
Definition: SubList.H:70
const word & system() const
Return system name.
Definition: TimePathsI.H:102
fileName caseConstant() const
Definition: TimePathsI.H:108
fileName caseSystem() const
Definition: TimePathsI.H:119
const word & constant() const
Return constant name.
Definition: TimePathsI.H:96
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:80
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
bool empty() const noexcept
True if the UList is empty (ie, size() is zero)
Definition: UListI.H:427
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:433
A bounding box defined in terms of min/max extrema points.
Definition: boundBox.H:64
const point & min() const
Minimum describing the bounding box.
Definition: boundBoxI.H:91
const point & max() const
Maximum describing the bounding box.
Definition: boundBoxI.H:97
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Particle-size distribution model wherein random samples are drawn from the doubly-truncated uniform p...
Definition: uniform.H:164
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:66
static int compare(const edge &a, const edge &b)
Compare edges.
Definition: edgeI.H:33
label operator[](const label i) const
Processor-local element id from linear-list of addresses.
Definition: ensightPart.H:189
A class for handling file names.
Definition: fileName.H:76
static std::string path(const std::string &str)
Return directory path name (part before last /)
Definition: fileNameI.H:176
static bool isAbsolute(const std::string &str)
Definition: fileNameI.H:136
virtual bool write()
Write the output fields.
scalar time
Injection time - set at collection [s].
Non-pointer based hierarchical recursive searching.
Definition: indexedOctree.H:74
static scalar & perturbTol()
Get the perturbation tolerance.
@ LITERAL
String literal.
Definition: keyType.H:81
A triFace with additional (region) index.
Definition: labelledTri.H:60
void movePoints()
Update for new mesh geometry.
Registry of regIOobjects.
const Time & time() const noexcept
Return time registry.
label count(const char *clsName) const
The number of objects of the given class name.
const vectorField & faceCentres() const
bool hasFaceCentres() const noexcept
label eventNo() const noexcept
Event number at last update.
Definition: regIOobjectI.H:191
Base class of (analytical or triangulated) surface. Encapsulates all the search routines....
virtual const boundBox & bounds() const
Return const reference to boundBox.
splitCell * master() const
Definition: splitCell.H:113
string & expand(const bool allowEmpty=false)
Definition: string.C:173
A class for managing temporary objects.
Definition: tmp.H:65
Standard boundBox with extra functionality for use in octree.
Definition: treeBoundBox.H:89
treeBoundBox extend(Random &rndGen, const scalar s) const
Return slightly wider bounding box.
Holds data for octree to work on an edges subset.
Definition: treeDataEdge.H:57
IOoject and searching on triSurface.
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.
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.
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.
virtual bool writeObject(IOstreamOption streamOpt, const bool valid) const
Write using stream options.
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 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.
void findLine(const pointField &start, const pointField &end, List< pointIndexHit > &info) const
void findLineAny(const pointField &start, const pointField &end, List< pointIndexHit > &info) const
void findNearest(const pointField &samples, const scalarField &nearestDistSqr, List< pointIndexHit > &info) const
Triangulated surface description with patch information.
Definition: triSurface.H:79
const geometricSurfacePatchList & patches() const noexcept
Definition: triSurface.H:399
static fileName relativeFilePath(const IOobject &io, const fileName &f, const bool isGlobal=true)
Return fileName.
Definition: triSurfaceIO.C:111
An enumeration wrapper for classification of a location as being inside/outside of a volume.
Definition: volumeType.H:61
static const Enum< volumeType::type > names
Names for the classification enumeration.
Definition: volumeType.H:76
@ UNKNOWN
Unknown state.
Definition: volumeType.H:67
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
const polyBoundaryMesh & patches
engineTime & runTime
PtrList< coordinateSystem > coordinates(solidRegions.size())
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
const labelList nEdges(UPstream::listGatherValues< label >(aMesh.nEdges()))
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))
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, false)
word timeName
Definition: getTimeIndex.H:3
#define WarningInFunction
Report a warning using Foam::Warning.
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:47
Foam::DimensionedField< label, triSurfaceGeoMesh > triSurfaceLabelField
const dimensionSet dimless
Dimensionless.
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i)
Definition: labelList.C:38
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: MSwindows.C:515
messageStream Info
Information stream (stdout output on master, null elsewhere)
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensionedScalar sqrt(const dimensionedScalar &ds)
Field< label > labelField
Specialisation of Field<T> for label.
Definition: labelField.H:52
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
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: MSwindows.C:666
error FatalError
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
labelList f(nPoints)
dictionary dict
volScalarField & e
Definition: createFields.H:11
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
Definition: stdFoam.H:278
Fields for triSurface.
scalarField samples(nIntervals, Zero)
Random rndGen
Definition: createFields.H:23