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(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 }
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.readOpt() != IOobject::NO_READ)
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  if (debug)
435  {
436  Pout<< "triSurfaceMesh(const IOobject& io) :"
437  << " loading surface " << io.objectPath()
438  << " local filePath:" << io.localFilePath(typeName)
439  << " from:" << actualFile << endl;
440  }
441 
442  if (searchGlobal && Pstream::parRun())
443  {
444  // Check where surface was found
445  const fileName localFile(io.localFilePath(typeName));
446 
447  if
448  (
449  r == masterOnly
450  && ((actualFile.empty() || actualFile != localFile))
451  )
452  {
453  // Found undecomposed surface. Load on master only
454  if (Pstream::master())
455  {
456  triSurface s2(actualFile);
458  }
460  if (debug)
461  {
462  Pout<< "triSurfaceMesh(const IOobject& io) :"
463  << " loaded triangles:" << triSurface::size() << endl;
464  }
465  }
466  else
467  {
468  // Read on all processors
469  triSurface s2(actualFile);
471  if (debug)
472  {
473  Pout<< "triSurfaceMesh(const IOobject& io) :"
474  << " loaded triangles:" << triSurface::size() << endl;
475  }
476  }
477  }
478  else
479  {
480  // Read on all processors
481  triSurface s2(actualFile);
483  if (debug)
484  {
485  Pout<< "triSurfaceMesh(const IOobject& io) :"
486  << " loaded triangles:" << triSurface::size() << endl;
487  }
488  }
489  }
490 
491  bounds() = boundBox(triSurface::points(), false);
492 }
493 
494 
496 (
497  const IOobject& io,
498  const dictionary& dict,
499  const readAction r
500 )
501 :
503  // Reused found instance in objectRegistry
505  (
506  IOobject
507  (
508  io.name(),
509  searchableSurface::instance(),
510  io.local(),
511  io.db(),
512  io.readOpt(),
513  io.writeOpt(),
514  false // searchableSurface already registered under name
515  )
516  ),
517  triSurface(),
518  triSurfaceRegionSearch(static_cast<const triSurface&>(*this), dict),
519  minQuality_(-1),
520  surfaceClosed_(-1),
521  outsideVolType_(volumeType::UNKNOWN)
522 {
523  if (io.readOpt() != IOobject::NO_READ)
524  {
525  // Surface type (optional)
526  const word surfType(dict.getOrDefault<word>("fileType", word::null));
527 
528  // Scale factor (optional)
529  const scalar scaleFactor(dict.getOrDefault<scalar>("scale", 0));
530 
531  const bool searchGlobal(r == localOrGlobal || r == masterOnly);
532 
533  fileName actualFile
534  (
535  searchGlobal
536  ? io.globalFilePath(typeName)
537  : io.localFilePath(typeName)
538  );
539 
540  // Reading from supplied file name instead of objectPath/filePath
541  if (dict.readIfPresent("file", fName_, keyType::LITERAL))
542  {
544  (
545  static_cast<const searchableSurface&>(*this),
546  fName_,
547  searchGlobal
548  );
549  actualFile = fName_;
550  }
551 
552  if (debug)
553  {
554  Pout<< "triSurfaceMesh(const IOobject& io, const dictionary&) :"
555  << " loading surface " << io.objectPath()
556  << " local filePath:" << io.localFilePath(typeName)
557  << " from:" << actualFile << endl;
558  }
559 
560 
561  if (searchGlobal && Pstream::parRun())
562  {
563  // Check where surface was found. Bit tricky:
564  // - master will have actualFile (in parent directory)
565  // different from localFilePath (in processor0/)
566  // - slave might have actualFile empty and localFile empty
567 
568  const fileName localFile(io.localFilePath(typeName));
569 
570  if
571  (
572  r == masterOnly
573  && ((actualFile.empty() || actualFile != localFile))
574  )
575  {
576  // Surface not loaded from processor directories -> undecomposed
577  // surface. Load on master only
578  if (Pstream::master())
579  {
580  triSurface s2(actualFile, surfType, scaleFactor);
582  }
584  if (debug)
585  {
586  Pout<< "triSurfaceMesh(const IOobject& io) :"
587  << " loaded triangles:" << triSurface::size() << endl;
588  }
589  }
590  else
591  {
592  // Read on all processors
593  triSurface s2(actualFile, surfType, scaleFactor);
595  if (debug)
596  {
597  Pout<< "triSurfaceMesh(const IOobject& io) :"
598  << " loaded triangles:" << triSurface::size() << endl;
599  }
600  }
601  }
602  else
603  {
604  // Read on all processors
605  triSurface s2(actualFile, surfType, scaleFactor);
607  if (debug)
608  {
609  Pout<< "triSurfaceMesh(const IOobject& io) :"
610  << " loaded triangles:" << triSurface::size() << endl;
611  }
612  }
613 
614 
615  // Report optional scale factor (eg, mm to m),
616  // which was already applied during triSurface reading
617  if (scaleFactor > 0)
618  {
620  << " : using scale " << scaleFactor << endl;
621  }
622  }
623 
624 
625  bounds() = boundBox(triSurface::points(), false);
626 
627  // Have optional minimum quality for normal calculation
628  if (dict.readIfPresent("minQuality", minQuality_) && minQuality_ > 0)
629  {
631  << " : ignoring triangles with quality < "
632  << minQuality_ << " for normals calculation." << endl;
633  }
634 }
636 
637 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
638 
640 {
641  clearOut();
642 }
643 
644 
646 {
647  // Do not clear closedness status
649  edgeTree_.clear();
651 }
653 
654 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
655 
657 {
658  auto tpts = tmp<pointField>::New();
659  auto& pts = tpts.ref();
660 
662  {
663  // Can reuse existing values
665  }
666  else
667  {
669 
670  // Calculate face centres from a copy to avoid incurring
671  // additional storage
673  (
674  FaceListType(*this, triSurface::size()),
676  ).faceCentres();
677  }
678 
679  return tpts;
680 }
681 
682 
684 (
685  pointField& centres,
686  scalarField& radiusSqr
687 ) const
688 {
689  centres = coordinates();
690  radiusSqr.setSize(size());
691  radiusSqr = 0.0;
692 
693  const pointField& pts = triSurface::points();
694 
695  forAll(*this, facei)
696  {
697  const labelledTri& f = triSurface::operator[](facei);
698  const point& fc = centres[facei];
699  for (const label pointi : f)
700  {
701  radiusSqr[facei] = max(radiusSqr[facei], fc.distSqr(pts[pointi]));
702  }
703  }
704 
705  // Add a bit to make sure all points are tested inside
706  radiusSqr += Foam::sqr(SMALL);
707 }
708 
709 
711 {
713 }
714 
715 
716 bool Foam::triSurfaceMesh::overlaps(const boundBox& bb) const
717 {
718  const indexedOctree<treeDataTriSurface>& octree = tree();
719 
720  labelList indices = octree.findBox(treeBoundBox(bb));
721 
722  return !indices.empty();
723 }
724 
725 
726 void Foam::triSurfaceMesh::movePoints(const pointField& newPoints)
727 {
728  if (debug)
729  {
730  Pout<< "triSurfaceMesh::movePoints :"
731  << " moving at time " << objectRegistry::time().timeName()
732  << endl;
733  }
734 
735  // Preserve topological point status (surfaceClosed_, outsideVolType_)
736 
737  // Update local information (instance, event number)
740 
741  const label event = getEvent();
742  searchableSurface::eventNo() = event;
744 
745  // Clear additional addressing
747  edgeTree_.clear();
748  triSurface::movePoints(newPoints);
749 
750  bounds() = boundBox(triSurface::points(), false);
751  if (debug)
752  {
753  Pout<< "triSurfaceMesh::movePoints: finished moving points" << endl;
754  }
755 }
756 
757 
760 {
761  if (!edgeTree_)
762  {
763  if (debug)
764  {
765  Pout<< "triSurfaceMesh::edgeTree :"
766  << " constructing tree for " << nEdges() - nInternalEdges()
767  << " boundary edges" << endl;
768  }
769 
770  // Boundary edges
771  const labelRange bEdges(nInternalEdges(), nBoundaryEdges());
772 
773  treeBoundBox bb(point::zero);
774 
775  if (bEdges.size())
776  {
777  label nPoints;
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.inflate(rndGen, 1e-4, ROOTVSMALL);
792  }
793 
794 
795  if (debug)
796  {
797  Pout<< "triSurfaceMesh::edgeTree : "
798  << "calculating edge tree for bb:" << bb << endl;
799  }
800 
801  const scalar oldTol =
803 
804  edgeTree_.reset
805  (
807  (
808  // Boundary edges
809  treeDataEdge(edges(), localPoints(), bEdges),
810 
811  bb, // bb
812  maxTreeDepth(), // maxLevel
813  10, // leafsize
814  3.0 // duplicity
815  )
816  );
817 
819 
820  if (debug)
821  {
822  Pout<< "triSurfaceMesh::edgeTree :"
823  << " finished constructing tree for "
824  << nEdges() - nInternalEdges()
825  << " boundary edges" << endl;
826  }
827  }
828 
829  return *edgeTree_;
830 }
831 
832 
834 {
835  if (regions_.empty())
836  {
837  regions_.setSize(patches().size());
838  forAll(regions_, regioni)
839  {
840  regions_[regioni] = patches()[regioni].name();
841  }
842  }
843  return regions_;
844 }
845 
846 
848 {
849  if (surfaceClosed_ == -1)
850  {
851  if (isSurfaceClosed())
852  {
853  surfaceClosed_ = 1;
854  }
855  else
856  {
857  surfaceClosed_ = 0;
858  }
859  }
860 
861  return surfaceClosed_ == 1;
862 }
863 
864 
866 {
867  if (outsideVolType_ == volumeType::UNKNOWN)
868  {
869  // Get point outside bounds()
870  const point outsidePt(bounds().max() + 0.5*bounds().span());
871 
872  if (debug)
873  {
874  Pout<< "triSurfaceMesh::outsideVolumeType :"
875  << " triggering outsidePoint:" << outsidePt
876  << " orientation" << endl;
877  }
878 
879  //outsideVolType_ = tree().shapes().getVolumeType(tree(), outsidePt);
880  // Note: do not use tree directly so e.g. distributedTriSurfaceMesh
881  // has opportunity to intercept
882  List<volumeType> outsideVolTypes;
883  getVolumeType(pointField(1, outsidePt), outsideVolTypes);
884  outsideVolType_ = outsideVolTypes[0];
885 
886  if (debug)
887  {
888  Pout<< "triSurfaceMesh::outsideVolumeType :"
889  << " finished outsidePoint:" << outsidePt
890  << " orientation:" << volumeType::names[outsideVolType_]
891  << endl;
892  }
893  }
894 
895  return outsideVolType_;
896 }
897 
898 
900 (
901  const pointField& samples,
902  const scalarField& nearestDistSqr,
903  List<pointIndexHit>& info
904 ) const
905 {
906  if (debug)
907  {
908  Pout<< "triSurfaceMesh::findNearest :"
909  << " trying to find nearest for " << samples.size()
910  << " samples with max sphere "
911  << (samples.size() ? Foam::sqrt(max(nearestDistSqr)) : Zero)
912  << endl;
913  }
914  triSurfaceSearch::findNearest(samples, nearestDistSqr, info);
915  if (debug)
916  {
917  Pout<< "triSurfaceMesh::findNearest :"
918  << " finished trying to find nearest for " << samples.size()
919  << " samples" << endl;
920  }
921 }
922 
923 
925 (
926  const pointField& samples,
927  const scalarField& nearestDistSqr,
928  const labelList& regionIndices,
929  List<pointIndexHit>& info
930 ) const
931 {
932  if (debug)
933  {
934  Pout<< "triSurfaceMesh::findNearest :"
935  << " trying to find nearest and region for " << samples.size()
936  << " samples with max sphere "
937  << (samples.size() ? Foam::sqrt(max(nearestDistSqr)) : Zero)
938  << endl;
939  }
941  (
942  samples,
943  nearestDistSqr,
944  regionIndices,
945  info
946  );
947  if (debug)
948  {
949  Pout<< "triSurfaceMesh::findNearest :"
950  << " finished trying to find nearest and region for "
951  << samples.size() << " samples" << endl;
952  }
953 }
954 
955 
957 (
958  const pointField& start,
959  const pointField& end,
960  List<pointIndexHit>& info
961 ) const
962 {
963  if (debug)
964  {
965  Pout<< "triSurfaceMesh::findLine :"
966  << " intersecting with "
967  << start.size() << " rays" << endl;
968  }
969  triSurfaceSearch::findLine(start, end, info);
970  if (debug)
971  {
972  Pout<< "triSurfaceMesh::findLine :"
973  << " finished intersecting with "
974  << start.size() << " rays" << endl;
975  }
976 }
977 
978 
980 (
981  const pointField& start,
982  const pointField& end,
983  List<pointIndexHit>& info
984 ) const
985 {
986  if (debug)
987  {
988  Pout<< "triSurfaceMesh::findLineAny :"
989  << " intersecting with "
990  << start.size() << " rays" << endl;
991  }
992  triSurfaceSearch::findLineAny(start, end, info);
993  if (debug)
994  {
995  Pout<< "triSurfaceMesh::findLineAny :"
996  << " finished intersecting with "
997  << start.size() << " rays" << endl;
998  }
999 }
1000 
1001 
1003 (
1004  const pointField& start,
1005  const pointField& end,
1006  List<List<pointIndexHit>>& info
1007 ) const
1008 {
1009  if (debug)
1010  {
1011  Pout<< "triSurfaceMesh::findLineAll :"
1012  << " intersecting with "
1013  << start.size() << " rays" << endl;
1014  }
1015  triSurfaceSearch::findLineAll(start, end, info);
1016  if (debug)
1017  {
1018  Pout<< "triSurfaceMesh::findLineAll :"
1019  << " finished intersecting with "
1020  << start.size() << " rays" << endl;
1021  }
1023 
1024 
1026 (
1027  const List<pointIndexHit>& info,
1028  labelList& region
1029 ) const
1030 {
1031  if (debug)
1032  {
1033  Pout<< "triSurfaceMesh::getRegion :"
1034  << " getting region for "
1035  << info.size() << " triangles" << endl;
1036  }
1037  region.setSize(info.size());
1038  forAll(info, i)
1039  {
1040  if (info[i].hit())
1041  {
1042  region[i] = triSurface::operator[](info[i].index()).region();
1043  }
1044  else
1045  {
1046  region[i] = -1;
1047  }
1048  }
1049  if (debug)
1050  {
1051  Pout<< "triSurfaceMesh::getRegion :"
1052  << " finished getting region for "
1053  << info.size() << " triangles" << endl;
1054  }
1056 
1057 
1059 (
1060  const List<pointIndexHit>& info,
1061  vectorField& normal
1062 ) const
1063 {
1064  if (debug)
1065  {
1066  Pout<< "triSurfaceMesh::getNormal :"
1067  << " getting normal for "
1068  << info.size() << " triangles" << endl;
1069  }
1070 
1071  const triSurface& s = *this;
1072  const pointField& pts = s.points();
1073 
1074  normal.setSize(info.size());
1075 
1076  if (minQuality_ >= 0)
1077  {
1078  // Make sure we don't use triangles with low quality since
1079  // normal is not reliable.
1080 
1081  const labelListList& faceFaces = s.faceFaces();
1082 
1083  forAll(info, i)
1084  {
1085  if (info[i].hit())
1086  {
1087  const label facei = info[i].index();
1088  normal[i] = s[facei].unitNormal(pts);
1089 
1090  scalar qual = s[facei].tri(pts).quality();
1091 
1092  if (qual < minQuality_)
1093  {
1094  // Search neighbouring triangles
1095  const labelList& fFaces = faceFaces[facei];
1096 
1097  for (const label nbri : fFaces)
1098  {
1099  scalar nbrQual = s[nbri].tri(pts).quality();
1100  if (nbrQual > qual)
1101  {
1102  qual = nbrQual;
1103  normal[i] = s[nbri].unitNormal(pts);
1104  }
1105  }
1106  }
1107  }
1108  else
1109  {
1110  // Set to what?
1111  normal[i] = Zero;
1112  }
1113  }
1114  }
1115  else
1116  {
1117  forAll(info, i)
1118  {
1119  if (info[i].hit())
1120  {
1121  const label facei = info[i].index();
1122 
1123  // Uncached
1124  normal[i] = s[facei].unitNormal(pts);
1125  }
1126  else
1127  {
1128  // Set to what?
1129  normal[i] = Zero;
1130  }
1131  }
1132  }
1133  if (debug)
1134  {
1135  Pout<< "triSurfaceMesh::getNormal :"
1136  << " finished getting normal for "
1137  << info.size() << " triangles" << endl;
1138  }
1139 }
1140 
1141 
1143 {
1144  auto* fldPtr = getObjectPtr<triSurfaceLabelField>("values");
1145 
1146  if (fldPtr)
1147  {
1148  (*fldPtr).field() = values;
1149  }
1150  else
1151  {
1152  fldPtr = new triSurfaceLabelField
1153  (
1154  IOobject
1155  (
1156  "values",
1157  objectRegistry::time().timeName(), // instance
1158  meshSubDir, // local
1159  *this,
1162  ),
1163  *this,
1164  dimless,
1166  );
1167 
1168  // Store field on triMesh
1169  fldPtr->store();
1170  }
1171  if (debug)
1172  {
1173  Pout<< "triSurfaceMesh::setField :"
1174  << " finished setting field for "
1175  << values.size() << " triangles" << endl;
1176  }
1178 
1179 
1181 (
1182  const List<pointIndexHit>& info,
1183  labelList& values
1184 ) const
1185 {
1186  const auto* fldPtr = getObjectPtr<triSurfaceLabelField>("values");
1187 
1188  if (fldPtr)
1189  {
1190  const auto& fld = *fldPtr;
1191 
1192  values.setSize(info.size());
1193 
1194  forAll(info, i)
1195  {
1196  if (info[i].hit())
1197  {
1198  values[i] = fld[info[i].index()];
1199  }
1200  }
1201  }
1202  if (debug)
1203  {
1204  Pout<< "triSurfaceMesh::setField :"
1205  << " finished getting field for "
1206  << info.size() << " triangles" << endl;
1207  }
1209 
1210 
1212 (
1213  const pointField& points,
1214  List<volumeType>& volType
1215 ) const
1216 {
1217  const scalar oldTol =
1219 
1220  if (debug)
1221  {
1222  Pout<< "triSurfaceMesh::getVolumeType :"
1223  << " finding orientation for " << points.size()
1224  << " samples" << endl;
1225  }
1226 
1227  volType.setSize(points.size());
1228 
1229  forAll(points, pointi)
1230  {
1231  const point& pt = points[pointi];
1232 
1233  if (tree().bb().contains(pt))
1234  {
1235  // Use cached volume type per each tree node
1236  volType[pointi] = tree().getVolumeType(pt);
1237  }
1238  else if (hasVolumeType())
1239  {
1240  // Precalculate and cache value for this outside point
1241  if (outsideVolType_ == volumeType::UNKNOWN)
1242  {
1243  outsideVolType_ = tree().shapes().getVolumeType(tree(), pt);
1244  }
1245  volType[pointi] = outsideVolType_;
1246  }
1247  else
1248  {
1249  // Have to calculate directly as outside the octree
1250  volType[pointi] = tree().shapes().getVolumeType(tree(), pt);
1251  }
1252  }
1253 
1255 
1256  if (debug)
1257  {
1258  Pout<< "triSurfaceMesh::getVolumeType :"
1259  << " finished finding orientation for " << points.size()
1260  << " samples" << endl;
1261  }
1263 
1264 
1266 (
1268  const bool valid
1269 ) const
1270 {
1272  const fileName& instance = searchableSurface::instance();
1273 
1274  if
1275  (
1276  instance != runTime.timeName()
1277  && instance != runTime.system()
1278  && instance != runTime.caseSystem()
1279  && instance != runTime.constant()
1280  && instance != runTime.caseConstant()
1281  )
1282  {
1283  const_cast<triSurfaceMesh&>(*this).searchableSurface::instance() =
1284  runTime.timeName();
1285  const_cast<triSurfaceMesh&>(*this).objectRegistry::instance() =
1286  runTime.timeName();
1287  }
1288 
1289  fileName fullPath;
1290  if (fName_.size())
1291  {
1292  // Override file name
1293 
1294  fullPath = fName_;
1295 
1296  fullPath.expand();
1297  if (!fullPath.isAbsolute())
1298  {
1299  // Add directory from regIOobject
1300  fullPath = searchableSurface::objectPath().path()/fullPath;
1301  }
1302  }
1303  else
1304  {
1305  fullPath = searchableSurface::objectPath();
1306  }
1307 
1308  if (!mkDir(fullPath.path()))
1309  {
1310  return false;
1311  }
1312 
1313  triSurface::write(fullPath);
1314 
1315  if (!isFile(fullPath))
1316  {
1317  return false;
1318  }
1319 
1320  return true;
1321 }
1322 
1323 
1324 // ************************************************************************* //
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.
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:51
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:118
virtual bool writeObject(IOstreamOption streamOpt, const bool valid) const
Write using stream options.
virtual tmp< pointField > coordinates() const
Get representative set of element coordinates.
A class for handling file names.
Definition: fileName.H:71
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:120
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:578
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:150
void inflate(const scalar factor)
Expand box by factor*mag(span) in all dimensions.
Definition: boundBoxI.H:367
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:51
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:49
label eventNo() const noexcept
Event number at last update.
Definition: regIOobjectI.H:184
bool empty() const noexcept
True if the UList is empty (ie, size() is zero)
Definition: UListI.H:420
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:487
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:639
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.
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:239
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 processes in communicator.
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:413
IOoject and searching on triSurface.
T & operator[](const label i)
Return element of UList.
Definition: UListI.H:292
fileName caseSystem() const
Return the system name for the case, which is ../system() for parallel runs.
Definition: TimePathsI.H:112
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:173
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:510
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:289
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, false)
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:567
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 INVALID.
Definition: exprTraits.C:52
const pointField & points
An edge is a list of two point labels. The functionality it provides supports the discretisation 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:107
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:95
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:203
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
fileName caseConstant() const
Return the constant name for the case, which is ../constant() for parallel runs.
Definition: TimePathsI.H:101
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 time name of given scalar time formatted with the given precision.
Definition: Time.C:760
const indexedOctree< treeDataEdge > & edgeTree() const
Demand driven construction of octree for boundary edges.
const word & constant() const noexcept
Return constant name.
Definition: TimePathsI.H:89
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:193
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:221
scalar minQuality_
Optional min triangle quality. Triangles below this get.
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
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:830
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
const Time & time() const
Return Time associated with the objectRegistry.
Definition: IOobject.C:447
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)
Am I the master rank.
Definition: UPstream.H:672
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)
label size() const noexcept
The number of elements in the UList.
Definition: UListI.H:413
List< label > labelList
A List of labels.
Definition: List.H:62
IntType size() const noexcept
The size of the range.
Definition: IntRange.H:189
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:72
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:166
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.
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:157