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 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "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 
39 namespace Foam
40 {
41  defineTypeNameAndDebug(triSurfaceMesh, 0);
42  addToRunTimeSelectionTable(searchableSurface, triSurfaceMesh, dict);
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;
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 }
282 
283 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
284 
286 :
289  (
290  IOobject
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 {
308 }
309 
310 
312 :
313  // Find instance for triSurfaceMesh
315  // Reused found instance in objectRegistry
317  (
318  IOobject
319  (
320  io.name(),
321  searchableSurface::instance(),
322  io.local(),
323  io.db(),
324  io.readOpt(),
325  io.writeOpt(),
326  false // searchableSurface already registered under name
327  )
328  ),
329  triSurface(static_cast<const searchableSurface&>(*this), dictionary::null),
330  triSurfaceRegionSearch(static_cast<const triSurface&>(*this)),
331  minQuality_(-1),
332  surfaceClosed_(-1),
333  outsideVolType_(volumeType::UNKNOWN)
334 {
335  bounds() = boundBox(triSurface::points(), false);
336 }
337 
338 
340 (
341  const IOobject& io,
342  const dictionary& dict
343 )
344 :
346  // Reused found instance in objectRegistry
348  (
349  IOobject
350  (
351  io.name(),
352  searchableSurface::instance(),
353  io.local(),
354  io.db(),
355  io.readOpt(),
356  io.writeOpt(),
357  false // searchableSurface already registered under name
358  )
359  ),
360  triSurface(static_cast<const searchableSurface&>(*this), dict),
361  triSurfaceRegionSearch(static_cast<const triSurface&>(*this), dict),
362  minQuality_(-1),
363  surfaceClosed_(-1),
364  outsideVolType_(volumeType::UNKNOWN)
365 {
366  // Read with supplied file name instead of objectPath/filePath
367  if (dict.readIfPresent("file", fName_, keyType::LITERAL))
368  {
370  (
371  static_cast<const searchableSurface&>(*this),
372  fName_,
373  true
374  );
375  }
376 
377  // Report optional scale factor (eg, mm to m),
378  // which was already applied during triSurface construction
379  scalar scaleFactor{0};
380  if (dict.getOrDefault("scale", scaleFactor) && scaleFactor > 0)
381  {
383  << " : using scale " << scaleFactor << endl;
384  }
385 
386  bounds() = boundBox(triSurface::points(), false);
387 
388  // Have optional minimum quality for normal calculation
389  if (dict.readIfPresent("minQuality", minQuality_) && minQuality_ > 0)
390  {
392  << " : ignoring triangles with quality < "
393  << minQuality_ << " for normals calculation." << endl;
394  }
395 }
396 
397 
398 Foam::triSurfaceMesh::triSurfaceMesh(const IOobject& io, const readAction r)
399 :
400  // Find instance for triSurfaceMesh
401  searchableSurface(io),
402  // Reused found instance in objectRegistry
403  objectRegistry
404  (
405  IOobject
406  (
407  io.name(),
408  searchableSurface::instance(),
409  io.local(),
410  io.db(),
411  io.readOpt(),
412  io.writeOpt(),
413  false // searchableSurface already registered under name
414  )
415  ),
416  triSurface(),
417  triSurfaceRegionSearch(static_cast<const triSurface&>(*this)),
418  minQuality_(-1),
419  surfaceClosed_(-1),
420  outsideVolType_(volumeType::UNKNOWN)
421 {
422  // Check IO flags
423  if (io.isAnyRead())
424  {
425  const bool searchGlobal(r == localOrGlobal || r == masterOnly);
426 
427  const fileName actualFile
428  (
429  searchGlobal
430  ? io.globalFilePath(typeName)
431  : io.localFilePath(typeName)
432  );
433 
434  // const fileName actualFile
435  // (
436  // io.filePath(searchGlobal, typeName)
437  // );
438 
439  if (debug)
440  {
441  Pout<< "triSurfaceMesh(const IOobject& io) :"
442  << " loading surface " << io.objectPath()
443  << " local filePath:" << io.localFilePath(typeName)
444  << " from:" << actualFile << endl;
445  }
446 
447  if (searchGlobal && Pstream::parRun())
448  {
449  // Check where surface was found
450  const fileName localFile(io.localFilePath(typeName));
451 
452  if
453  (
454  r == masterOnly
455  && ((actualFile.empty() || actualFile != localFile))
456  )
457  {
458  // Found undecomposed surface. Load on master only
459  if (Pstream::master())
460  {
461  triSurface s2(actualFile);
463  }
465  if (debug)
466  {
467  Pout<< "triSurfaceMesh(const IOobject& io) :"
468  << " loaded triangles:" << triSurface::size() << endl;
469  }
470  }
471  else
472  {
473  // Read on all processors
474  triSurface s2(actualFile);
476  if (debug)
477  {
478  Pout<< "triSurfaceMesh(const IOobject& io) :"
479  << " loaded triangles:" << triSurface::size() << endl;
480  }
481  }
482  }
483  else
484  {
485  // Read on all processors
486  triSurface s2(actualFile);
488  if (debug)
489  {
490  Pout<< "triSurfaceMesh(const IOobject& io) :"
491  << " loaded triangles:" << triSurface::size() << endl;
492  }
493  }
494  }
495 
496  bounds() = boundBox(triSurface::points(), false);
497 }
498 
499 
501 (
502  const IOobject& io,
503  const dictionary& dict,
504  const readAction r
505 )
506 :
508  // Reused found instance in objectRegistry
510  (
511  IOobject
512  (
513  io.name(),
514  searchableSurface::instance(),
515  io.local(),
516  io.db(),
517  io.readOpt(),
518  io.writeOpt(),
519  false // searchableSurface already registered under name
520  )
521  ),
522  triSurface(),
523  triSurfaceRegionSearch(static_cast<const triSurface&>(*this), dict),
524  minQuality_(-1),
525  surfaceClosed_(-1),
526  outsideVolType_(volumeType::UNKNOWN)
527 {
528  if (io.isAnyRead())
529  {
530  // Surface type (optional)
531  const word surfType(dict.getOrDefault<word>("fileType", word::null));
532 
533  // Scale factor (optional)
534  const scalar scaleFactor(dict.getOrDefault<scalar>("scale", 0));
535 
536  const bool searchGlobal(r == localOrGlobal || r == masterOnly);
537 
538  fileName actualFile
539  (
540  searchGlobal
541  ? io.globalFilePath(typeName)
542  : io.localFilePath(typeName)
543  );
544 
545  // const fileName actualFile
546  // (
547  // io.filePath(searchGlobal, typeName)
548  // );
549 
550  // Reading from supplied file name instead of objectPath/filePath
551  if (dict.readIfPresent("file", fName_, keyType::LITERAL))
552  {
554  (
555  static_cast<const searchableSurface&>(*this),
556  fName_,
557  searchGlobal
558  );
559  actualFile = fName_;
560  }
561 
562  if (debug)
563  {
564  Pout<< "triSurfaceMesh(const IOobject& io, const dictionary&) :"
565  << " loading surface " << io.objectPath()
566  << " local filePath:" << io.localFilePath(typeName)
567  << " from:" << actualFile << endl;
568  }
569 
570 
571  if (searchGlobal && Pstream::parRun())
572  {
573  // Check where surface was found. Bit tricky:
574  // - master will have actualFile (in parent directory)
575  // different from localFilePath (in processor0/)
576  // - slave might have actualFile empty and localFile empty
577 
578  const fileName localFile(io.localFilePath(typeName));
579 
580  if
581  (
582  r == masterOnly
583  && ((actualFile.empty() || actualFile != localFile))
584  )
585  {
586  // Surface not loaded from processor directories -> undecomposed
587  // surface. Load on master only
588  if (Pstream::master())
589  {
590  triSurface s2(actualFile, surfType, scaleFactor);
592  }
594  if (debug)
595  {
596  Pout<< "triSurfaceMesh(const IOobject& io) :"
597  << " loaded triangles:" << triSurface::size() << endl;
598  }
599  }
600  else
601  {
602  // Read on all processors
603  triSurface s2(actualFile, surfType, scaleFactor);
605  if (debug)
606  {
607  Pout<< "triSurfaceMesh(const IOobject& io) :"
608  << " loaded triangles:" << triSurface::size() << endl;
609  }
610  }
611  }
612  else
613  {
614  // Read on all processors
615  triSurface s2(actualFile, surfType, scaleFactor);
617  if (debug)
618  {
619  Pout<< "triSurfaceMesh(const IOobject& io) :"
620  << " loaded triangles:" << triSurface::size() << endl;
621  }
622  }
623 
624 
625  // Report optional scale factor (eg, mm to m),
626  // which was already applied during triSurface reading
627  if (scaleFactor > 0)
628  {
630  << " : using scale " << scaleFactor << endl;
631  }
632  }
633 
634 
635  bounds() = boundBox(triSurface::points(), false);
636 
637  // Have optional minimum quality for normal calculation
638  if (dict.readIfPresent("minQuality", minQuality_) && minQuality_ > 0)
639  {
641  << " : ignoring triangles with quality < "
642  << minQuality_ << " for normals calculation." << endl;
643  }
644 }
646 
647 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
648 
650 {
651  clearOut();
652 }
653 
654 
656 {
657  // Do not clear closedness status
659  edgeTree_.clear();
661 }
663 
664 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
665 
667 {
668  auto tpts = tmp<pointField>::New();
669  auto& pts = tpts.ref();
670 
672  {
673  // Can reuse existing values
675  }
676  else
677  {
679 
680  // Calculate face centres from a copy to avoid incurring
681  // additional storage
683  (
684  FaceListType(*this, triSurface::size()),
686  ).faceCentres();
687  }
688 
689  return tpts;
690 }
691 
692 
694 (
695  pointField& centres,
696  scalarField& radiusSqr
697 ) const
698 {
699  centres = coordinates();
700  radiusSqr.setSize(size());
701  radiusSqr = 0.0;
702 
703  const pointField& pts = triSurface::points();
704 
705  forAll(*this, facei)
706  {
707  const labelledTri& f = triSurface::operator[](facei);
708  const point& fc = centres[facei];
709  for (const label pointi : f)
710  {
711  radiusSqr[facei] = max(radiusSqr[facei], fc.distSqr(pts[pointi]));
712  }
713  }
714 
715  // Add a bit to make sure all points are tested inside
716  radiusSqr += Foam::sqr(SMALL);
717 }
718 
719 
721 {
723 }
724 
725 
726 bool Foam::triSurfaceMesh::overlaps(const boundBox& bb) const
727 {
728  const indexedOctree<treeDataTriSurface>& octree = tree();
729 
730  labelList indices = octree.findBox(treeBoundBox(bb));
731 
732  return !indices.empty();
733 }
734 
735 
736 void Foam::triSurfaceMesh::movePoints(const pointField& newPoints)
737 {
738  if (debug)
739  {
740  Pout<< "triSurfaceMesh::movePoints :"
741  << " moving at time " << objectRegistry::time().timeName()
742  << endl;
743  }
744 
745  // Preserve topological point status (surfaceClosed_, outsideVolType_)
746 
747  // Update local information (instance, event number)
750 
751  const label event = getEvent();
752  searchableSurface::eventNo() = event;
754 
755  // Clear additional addressing
757  edgeTree_.clear();
758  triSurface::movePoints(newPoints);
759 
760  bounds() = boundBox(triSurface::points(), false);
761  if (debug)
762  {
763  Pout<< "triSurfaceMesh::movePoints: finished moving points" << endl;
764  }
765 }
766 
767 
770 {
771  if (!edgeTree_)
772  {
773  if (debug)
774  {
775  Pout<< "triSurfaceMesh::edgeTree :"
776  << " constructing tree for " << nEdges() - nInternalEdges()
777  << " boundary edges" << endl;
778  }
779 
780  // Boundary edges
781  const labelRange bEdges(nInternalEdges(), nBoundaryEdges());
782 
783  treeBoundBox bb(point::zero);
784 
785  if (bEdges.size())
786  {
787  label nPoints;
789  (
790  *this,
791  bb,
792  nPoints
793  );
794 
795  // Random number generator. Bit dodgy since not exactly random ;-)
796  Random rndGen(65431);
797 
798  // Slightly extended bb. Slightly off-centred just so on symmetric
799  // geometry there are less face/edge aligned items.
800 
801  bb.inflate(rndGen, 1e-4, ROOTVSMALL);
802  }
803 
804 
805  if (debug)
806  {
807  Pout<< "triSurfaceMesh::edgeTree : "
808  << "calculating edge tree for bb:" << bb << endl;
809  }
810 
811  const scalar oldTol =
813 
814  edgeTree_.reset
815  (
817  (
818  // Boundary edges
819  treeDataEdge(edges(), localPoints(), bEdges),
820 
821  bb, // bb
822  maxTreeDepth(), // maxLevel
823  10, // leafsize
824  3.0 // duplicity
825  )
826  );
827 
829 
830  if (debug)
831  {
832  Pout<< "triSurfaceMesh::edgeTree :"
833  << " finished constructing tree for "
834  << nEdges() - nInternalEdges()
835  << " boundary edges" << endl;
836  }
837  }
838 
839  return *edgeTree_;
840 }
841 
842 
844 {
845  if (regions_.empty())
846  {
847  regions_.setSize(patches().size());
848  forAll(regions_, regioni)
849  {
850  regions_[regioni] = patches()[regioni].name();
851  }
852  }
853  return regions_;
854 }
855 
856 
858 {
859  if (surfaceClosed_ == -1)
860  {
861  if (isSurfaceClosed())
862  {
863  surfaceClosed_ = 1;
864  }
865  else
866  {
867  surfaceClosed_ = 0;
868  }
869  }
870 
871  return surfaceClosed_ == 1;
872 }
873 
874 
876 {
877  if (outsideVolType_ == volumeType::UNKNOWN)
878  {
879  // Get point outside bounds()
880  const point outsidePt(bounds().max() + 0.5*bounds().span());
881 
882  if (debug)
883  {
884  Pout<< "triSurfaceMesh::outsideVolumeType :"
885  << " triggering outsidePoint:" << outsidePt
886  << " orientation" << endl;
887  }
888 
889  //outsideVolType_ = tree().shapes().getVolumeType(tree(), outsidePt);
890  // Note: do not use tree directly so e.g. distributedTriSurfaceMesh
891  // has opportunity to intercept
892  List<volumeType> outsideVolTypes;
893  getVolumeType(pointField(1, outsidePt), outsideVolTypes);
894  outsideVolType_ = outsideVolTypes[0];
895 
896  if (debug)
897  {
898  Pout<< "triSurfaceMesh::outsideVolumeType :"
899  << " finished outsidePoint:" << outsidePt
900  << " orientation:" << volumeType::names[outsideVolType_]
901  << endl;
902  }
903  }
904 
905  return outsideVolType_;
906 }
907 
908 
910 (
911  const pointField& samples,
912  const scalarField& nearestDistSqr,
913  List<pointIndexHit>& info
914 ) const
915 {
916  if (debug)
917  {
918  Pout<< "triSurfaceMesh::findNearest :"
919  << " trying to find nearest for " << samples.size()
920  << " samples with max sphere "
921  << (samples.size() ? Foam::sqrt(max(nearestDistSqr)) : Zero)
922  << endl;
923  }
924  triSurfaceSearch::findNearest(samples, nearestDistSqr, info);
925  if (debug)
926  {
927  Pout<< "triSurfaceMesh::findNearest :"
928  << " finished trying to find nearest for " << samples.size()
929  << " samples" << endl;
930  }
931 }
932 
933 
935 (
936  const pointField& samples,
937  const scalarField& nearestDistSqr,
938  const labelList& regionIndices,
939  List<pointIndexHit>& info
940 ) const
941 {
942  if (debug)
943  {
944  Pout<< "triSurfaceMesh::findNearest :"
945  << " trying to find nearest and region for " << samples.size()
946  << " samples with max sphere "
947  << (samples.size() ? Foam::sqrt(max(nearestDistSqr)) : Zero)
948  << endl;
949  }
951  (
952  samples,
953  nearestDistSqr,
954  regionIndices,
955  info
956  );
957  if (debug)
958  {
959  Pout<< "triSurfaceMesh::findNearest :"
960  << " finished trying to find nearest and region for "
961  << samples.size() << " samples" << endl;
962  }
963 }
964 
965 
967 (
968  const pointField& start,
969  const pointField& end,
970  List<pointIndexHit>& info
971 ) const
972 {
973  if (debug)
974  {
975  Pout<< "triSurfaceMesh::findLine :"
976  << " intersecting with "
977  << start.size() << " rays" << endl;
978  }
979  triSurfaceSearch::findLine(start, end, info);
980  if (debug)
981  {
982  Pout<< "triSurfaceMesh::findLine :"
983  << " finished intersecting with "
984  << start.size() << " rays" << endl;
985  }
986 }
987 
988 
990 (
991  const pointField& start,
992  const pointField& end,
993  List<pointIndexHit>& info
994 ) const
995 {
996  if (debug)
997  {
998  Pout<< "triSurfaceMesh::findLineAny :"
999  << " intersecting with "
1000  << start.size() << " rays" << endl;
1001  }
1002  triSurfaceSearch::findLineAny(start, end, info);
1003  if (debug)
1004  {
1005  Pout<< "triSurfaceMesh::findLineAny :"
1006  << " finished intersecting with "
1007  << start.size() << " rays" << endl;
1008  }
1010 
1011 
1013 (
1014  const pointField& start,
1015  const pointField& end,
1016  List<List<pointIndexHit>>& info
1017 ) const
1018 {
1019  if (debug)
1020  {
1021  Pout<< "triSurfaceMesh::findLineAll :"
1022  << " intersecting with "
1023  << start.size() << " rays" << endl;
1024  }
1025  triSurfaceSearch::findLineAll(start, end, info);
1026  if (debug)
1027  {
1028  Pout<< "triSurfaceMesh::findLineAll :"
1029  << " finished intersecting with "
1030  << start.size() << " rays" << endl;
1031  }
1033 
1034 
1036 (
1037  const List<pointIndexHit>& info,
1038  labelList& region
1039 ) const
1040 {
1041  if (debug)
1042  {
1043  Pout<< "triSurfaceMesh::getRegion :"
1044  << " getting region for "
1045  << info.size() << " triangles" << endl;
1046  }
1047  region.setSize(info.size());
1048  forAll(info, i)
1049  {
1050  if (info[i].hit())
1051  {
1052  region[i] = triSurface::operator[](info[i].index()).region();
1053  }
1054  else
1055  {
1056  region[i] = -1;
1057  }
1058  }
1059  if (debug)
1060  {
1061  Pout<< "triSurfaceMesh::getRegion :"
1062  << " finished getting region for "
1063  << info.size() << " triangles" << endl;
1064  }
1066 
1067 
1069 (
1070  const List<pointIndexHit>& info,
1071  vectorField& normal
1072 ) const
1073 {
1074  if (debug)
1075  {
1076  Pout<< "triSurfaceMesh::getNormal :"
1077  << " getting normal for "
1078  << info.size() << " triangles" << endl;
1079  }
1080 
1081  const triSurface& s = *this;
1082  const pointField& pts = s.points();
1083 
1084  normal.setSize(info.size());
1085 
1086  if (minQuality_ >= 0)
1087  {
1088  // Make sure we don't use triangles with low quality since
1089  // normal is not reliable.
1090 
1091  const labelListList& faceFaces = s.faceFaces();
1092 
1093  forAll(info, i)
1094  {
1095  if (info[i].hit())
1096  {
1097  const label facei = info[i].index();
1098  normal[i] = s[facei].unitNormal(pts);
1099 
1100  scalar qual = s[facei].tri(pts).quality();
1101 
1102  if (qual < minQuality_)
1103  {
1104  // Search neighbouring triangles
1105  const labelList& fFaces = faceFaces[facei];
1106 
1107  for (const label nbri : fFaces)
1108  {
1109  scalar nbrQual = s[nbri].tri(pts).quality();
1110  if (nbrQual > qual)
1111  {
1112  qual = nbrQual;
1113  normal[i] = s[nbri].unitNormal(pts);
1114  }
1115  }
1116  }
1117  }
1118  else
1119  {
1120  // Set to what?
1121  normal[i] = Zero;
1122  }
1123  }
1124  }
1125  else
1126  {
1127  forAll(info, i)
1128  {
1129  if (info[i].hit())
1130  {
1131  const label facei = info[i].index();
1132 
1133  // Uncached
1134  normal[i] = s[facei].unitNormal(pts);
1135  }
1136  else
1137  {
1138  // Set to what?
1139  normal[i] = Zero;
1140  }
1141  }
1142  }
1143  if (debug)
1144  {
1145  Pout<< "triSurfaceMesh::getNormal :"
1146  << " finished getting normal for "
1147  << info.size() << " triangles" << endl;
1148  }
1149 }
1150 
1151 
1153 {
1154  auto* fldPtr = getObjectPtr<triSurfaceLabelField>("values");
1155 
1156  if (fldPtr)
1157  {
1158  (*fldPtr).field() = values;
1159  }
1160  else
1161  {
1162  fldPtr = new triSurfaceLabelField
1163  (
1164  IOobject
1165  (
1166  "values",
1167  objectRegistry::time().timeName(), // instance
1168  meshSubDir, // local
1169  *this,
1173  ),
1174  *this,
1175  dimless,
1177  );
1178 
1179  // Store field on triMesh
1180  fldPtr->store();
1181  }
1182  if (debug)
1183  {
1184  Pout<< "triSurfaceMesh::setField :"
1185  << " finished setting field for "
1186  << values.size() << " triangles" << endl;
1187  }
1189 
1190 
1192 (
1193  const List<pointIndexHit>& info,
1194  labelList& values
1195 ) const
1196 {
1197  const auto* fldPtr = getObjectPtr<triSurfaceLabelField>("values");
1198 
1199  if (fldPtr)
1200  {
1201  const auto& fld = *fldPtr;
1202 
1203  values.setSize(info.size());
1204 
1205  forAll(info, i)
1206  {
1207  if (info[i].hit())
1208  {
1209  values[i] = fld[info[i].index()];
1210  }
1211  }
1212  }
1213  if (debug)
1214  {
1215  Pout<< "triSurfaceMesh::setField :"
1216  << " finished getting field for "
1217  << info.size() << " triangles" << endl;
1218  }
1220 
1221 
1223 (
1224  const pointField& points,
1225  List<volumeType>& volType
1226 ) const
1227 {
1228  const scalar oldTol =
1230 
1231  if (debug)
1232  {
1233  Pout<< "triSurfaceMesh::getVolumeType :"
1234  << " finding orientation for " << points.size()
1235  << " samples" << endl;
1236  }
1237 
1238  volType.setSize(points.size());
1239 
1240  forAll(points, pointi)
1241  {
1242  const point& pt = points[pointi];
1243 
1244  if (tree().bb().contains(pt))
1245  {
1246  // Use cached volume type per each tree node
1247  volType[pointi] = tree().getVolumeType(pt);
1248  }
1249  else if (hasVolumeType())
1250  {
1251  // Precalculate and cache value for this outside point
1252  if (outsideVolType_ == volumeType::UNKNOWN)
1253  {
1254  outsideVolType_ = tree().shapes().getVolumeType(tree(), pt);
1255  }
1256  volType[pointi] = outsideVolType_;
1257  }
1258  else
1259  {
1260  // Have to calculate directly as outside the octree
1261  volType[pointi] = tree().shapes().getVolumeType(tree(), pt);
1262  }
1263  }
1264 
1266 
1267  if (debug)
1268  {
1269  Pout<< "triSurfaceMesh::getVolumeType :"
1270  << " finished finding orientation for " << points.size()
1271  << " samples" << endl;
1272  }
1274 
1275 
1277 (
1279  const bool writeOnProc
1280 ) const
1281 {
1283  const fileName& instance = searchableSurface::instance();
1284 
1285  if
1286  (
1287  instance != runTime.timeName()
1288  && instance != runTime.system()
1289  && instance != runTime.caseSystem()
1290  && instance != runTime.constant()
1291  && instance != runTime.caseConstant()
1292  )
1293  {
1294  const_cast<triSurfaceMesh&>(*this).searchableSurface::instance() =
1295  runTime.timeName();
1296  const_cast<triSurfaceMesh&>(*this).objectRegistry::instance() =
1297  runTime.timeName();
1298  }
1299 
1300  fileName fullPath;
1301  if (fName_.size())
1302  {
1303  // Override file name
1304 
1305  fullPath = fName_;
1306 
1307  fullPath.expand();
1308  if (!fullPath.isAbsolute())
1309  {
1310  // Add directory from regIOobject
1311  fullPath = searchableSurface::objectPath().path()/fullPath;
1312  }
1313  }
1314  else
1315  {
1316  fullPath = searchableSurface::objectPath();
1317  }
1318 
1319  if (!mkDir(fullPath.path()))
1320  {
1321  return false;
1322  }
1323 
1324  triSurface::write(fullPath);
1325 
1326  if (!isFile(fullPath))
1327  {
1328  return false;
1329  }
1330 
1331  return true;
1332 }
1333 
1334 
1335 // ************************************************************************* //
static scalar & perturbTol() noexcept
Get the perturbation tolerance.
virtual void findLine(const pointField &start, const pointField &end, List< pointIndexHit > &) const
Find first intersection on segment from start to end.
dictionary dict
Field< label > labelField
Specialisation of Field<T> for label.
Definition: labelField.H:48
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
virtual tmp< pointField > coordinates() const
Get representative set of element coordinates.
A class for handling file names.
Definition: fileName.H:72
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
void clearOut()
Clear storage.
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.
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
virtual void boundingSpheres(pointField &centres, scalarField &radiusSqr) const
Get bounding spheres (centre and radius squared). Any point.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:195
void inflate(const scalar factor)
Expand box by factor*mag(span) in all dimensions.
Definition: boundBoxI.H:381
Fields for triSurface.
virtual tmp< pointField > points() const
Get the points that define the surface.
virtual ~triSurfaceMesh()
Destructor.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
dimensionedSymmTensor sqr(const dimensionedVector &dv)
A range or interval of labels defined by a start and a size.
Definition: labelRange.H:52
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
label eventNo() const noexcept
Event number at last update.
Definition: regIOobjectI.H:186
bool empty() const noexcept
True if List is empty (ie, size() is zero)
Definition: UList.H:666
Random rndGen
Definition: createFields.H:23
virtual void findNearest(const pointField &sample, const scalarField &nearestDistSqr, List< pointIndexHit > &) const
virtual void movePoints(const pointField &pts)
Move points.
Definition: triSurface.C:605
scalarField samples(nIntervals, Zero)
virtual void setField(const labelList &values)
WIP. Store element-wise field.
engineTime & runTime
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:1049
static bool addFaceToEdge(const edge &, EdgeMap< label > &)
Helper function for isSurfaceClosed.
Foam::DimensionedField< label, triSurfaceGeoMesh > triSurfaceLabelField
void transfer(triSurface &surf)
Alter contents by transferring (triangles, points) components.
Definition: triSurface.C:954
Helper class to search on triSurface. Creates an octree for each region of the surface and only searc...
virtual bool overlaps(const boundBox &bb) const
Does any part of the surface overlap the supplied bound box?
Holds data for octree to work on an edges subset.
Definition: treeDataEdge.H:53
A bounding box defined in terms of min/max extrema points.
Definition: boundBox.H:63
A simple container for options an IOstream can normally have.
static std::string path(const std::string &str)
Return directory path name (part before last /)
Definition: fileNameI.H:169
Base class of (analytical or triangulated) surface. Encapsulates all the search routines. WIP.
const dimensionSet dimless
Dimensionless.
An enumeration wrapper for classification of a location as being inside/outside of a volume...
Definition: volumeType.H:55
bool isSurfaceClosed() const
Check whether surface is closed without calculating any permanent.
virtual const boundBox & bounds() const
Return const reference to boundBox.
List< labelList > labelListList
List of labelList.
Definition: labelList.H:38
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
void write(Ostream &os) const
Write to Ostream in simple OpenFOAM format.
Definition: triSurfaceIO.C:331
fileName objectPath() const
The complete path + object name.
Definition: IOobjectI.H:284
Macros for easy insertion into run-time selection tables.
Unknown state.
Definition: volumeType.H:64
static word meshSubDir
Return the mesh sub-directory name (usually "triSurface")
static void broadcast(Type &value, const label comm=UPstream::worldComm)
Broadcast content (contiguous or non-contiguous) to all communicator ranks. Does nothing in non-paral...
void findNearest(const pointField &samples, const scalarField &nearestDistSqr, List< pointIndexHit > &info) const
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
IOoject and searching on triSurface.
T & operator[](const label i)
Return element of UList.
Definition: UListI.H:361
fileName caseSystem() const
Return the system name for the case, which is ../system() for parallel runs.
Definition: TimePathsI.H:135
virtual bool hasVolumeType() const
Whether supports volume type (below) - i.e. whether is closed.
bool insert(const Key &key, const T &obj)
Copy insert a new entry, not overwriting existing entries.
Definition: HashTableI.H:152
List< T > values(const HashTable< T, Key, Hash > &tbl, const bool doSort=false)
List of values from HashTable, optionally sorted.
Definition: HashOps.H:164
word timeName
Definition: getTimeIndex.H:3
static fileName relativeFilePath(const IOobject &io, const fileName &f, const bool isGlobal=true)
Return fileName.
Definition: triSurfaceIO.C:106
A list of faces which address into the list of points.
A List obtained as a section of another List.
Definition: SubList.H:50
static void calcBounds(const PrimitivePatch< FaceList, PointField > &p, boundBox &bb, label &nPoints)
const geometricSurfacePatchList & patches() const noexcept
Definition: triSurface.H:509
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
void setSize(const label n)
Alias for resize()
Definition: List.H:316
virtual void getVolumeType(const pointField &points, List< volumeType > &volType) const
Determine type (inside/outside/mixed) for point.
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:614
void findLine(const pointField &start, const pointField &end, List< pointIndexHit > &info) const
label count(const char *clsName) const
The number of objects of the given class name.
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
const pointField & points
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:59
void findLineAll(const pointField &start, const pointField &end, List< List< pointIndexHit >> &info) const
Calculate all intersections from start to end.
A class for handling words, derived from Foam::string.
Definition: word.H:63
label nPoints
iterator find(const Key &key)
Find and return an iterator set at the hashed entry.
Definition: HashTableI.H:86
const Time & time() const noexcept
Return time registry.
void findLineAny(const pointField &start, const pointField &end, List< pointIndexHit > &info) const
const word & system() const noexcept
Return system name.
Definition: TimePathsI.H:118
const Field< point_type > & faceCentres() const
Return face centres for patch.
Tree tree(triangles.begin(), triangles.end())
virtual void findLineAll(const pointField &start, const pointField &end, List< List< pointIndexHit >> &) const
Get all intersections in order from start to end.
static const word null
An empty word.
Definition: word.H:84
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
bool local
Definition: EEqn.H:20
A triFace with additional (region) index.
Definition: labelledTri.H:53
const Field< point_type > & points() const noexcept
Return reference to global points.
String literal.
Definition: keyType.H:82
Random number generator.
Definition: Random.H:55
const Time & time() const noexcept
Return Time associated with the objectRegistry.
Definition: IOobject.C:456
fileName caseConstant() const
Return the constant name for the case, which is ../constant() for parallel runs.
Definition: TimePathsI.H:124
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
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
label size() const noexcept
The number of elements in the container.
Definition: UList.H:671
const indexedOctree< treeDataEdge > & edgeTree() const
Demand driven construction of octree for boundary edges.
const word & constant() const noexcept
Return constant name.
Definition: TimePathsI.H:112
int debug
Static debugging option.
constexpr auto end(C &c) -> decltype(c.end())
Return iterator to the end of the container c.
Definition: stdFoam.H:201
defineTypeNameAndDebug(combustionModel, 0)
labelList f(nPoints)
virtual void getRegion(const List< pointIndexHit > &, labelList &region) const
From a set of points and indices get the region.
const fileName & instance() const noexcept
Read access to instance path component.
Definition: IOobjectI.H:266
scalar minQuality_
Optional min triangle quality. Triangles below this get.
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;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
virtual void getNormal(const List< pointIndexHit > &, vectorField &normal) const
From a set of points and indices get the normal.
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:877
triSurface()
Default construct.
Definition: triSurface.C:425
PtrList< coordinateSystem > coordinates(solidRegions.size())
vector point
Point is a vector.
Definition: point.H:37
Non-pointer based hierarchical recursive searching.
#define WarningInFunction
Report a warning using Foam::Warning.
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
virtual bool writeObject(IOstreamOption streamOpt, const bool writeOnProc) const
Write using stream options.
virtual void findLineAny(const pointField &start, const pointField &end, List< pointIndexHit > &) const
Return any intersection on segment from start to end.
virtual void movePoints(const pointField &)
Move points.
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition: UPstream.H:1082
const polyBoundaryMesh & patches
Nothing to be read.
Automatically write from objectRegistry::writeObject()
Standard boundBox with extra functionality for use in octree.
Definition: treeBoundBox.H:90
messageStream Info
Information stream (stdout output on master, null elsewhere)
List< label > labelList
A List of labels.
Definition: List.H:62
IntType size() const noexcept
The size of the range.
Definition: IntRange.H:194
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, IOobject::NO_REGISTER)
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Registry of regIOobjects.
virtual volumeType outsideVolumeType() const
If surface is closed, what is type of points outside bounds.
triSurfaceMesh(const triSurfaceMesh &)=delete
No copy construct.
Triangulated surface description with patch information.
Definition: triSurface.H:71
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))
fileName fName_
Supplied fileName override.
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
Request registration (bool: true)
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
virtual const wordList & regions() const
Names of regions.
Namespace for OpenFOAM.
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
virtual void getField(const List< pointIndexHit > &, labelList &) const
WIP. From a set of hits (points and.
addToRunTimeSelectionTable(functionObject, pointHistory, dictionary)
static const Enum< volumeType::type > names
Names for the classification enumeration.
Definition: volumeType.H:75
static int compare(const edge &a, const edge &b)
Compare edges.
Definition: edgeI.H:26
const pointField & pts
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127