searchableSurfaces.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2023 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "searchableSurfaces.H"
28 #include "ListOps.H"
29 #include "Time.H"
30 #include "DynamicField.H"
31 #include "PatchTools.H"
32 #include "triSurfaceMesh.H"
33 #include "vtkWritePolyData.H"
34 
35 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
40 }
41 
42 
43 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
44 
45 bool Foam::searchableSurfaces::connected
46 (
47  const triSurface& s,
48  const label edgeI,
49  const pointIndexHit& hit
50 )
51 {
52  const triFace& localFace = s.localFaces()[hit.index()];
53  const edge& e = s.edges()[edgeI];
54 
55  forAll(localFace, i)
56  {
57  if (e.otherVertex(localFace[i]) != -1)
58  {
59  return true;
60  }
61  }
62 
63  return false;
64 }
65 
66 
67 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
68 
69 // Construct with length.
71 :
73  regionNames_(size),
74  allSurfaces_(identityMap(size))
75 {}
76 
77 
78 //Foam::searchableSurfaces::searchableSurfaces
79 //(
80 // const IOobject& io,
81 // const PtrList<dictionary>& dicts
82 //)
83 //:
84 // PtrList<searchableSurface>(dicts.size()),
85 // regionNames_(dicts.size()),
86 // allSurfaces_(identityMap(dicts.size()))
87 //{
88 // forAll(dicts, surfI)
89 // {
90 // const dictionary& dict = dicts[surfI];
91 //
92 // // Make IOobject with correct name
93 // autoPtr<IOobject> namedIO(io.clone());
94 // namedIO().rename(dict.lookup("name"));
95 //
96 // // Create and hook surface
97 // set
98 // (
99 // surfI,
100 // searchableSurface::New
101 // (
102 // dict.lookup("type"),
103 // namedIO(),
104 // dict
105 // )
106 // );
107 // const searchableSurface& s = operator[](surfI);
108 //
109 // // Construct default region names by prepending surface name
110 // // to region name.
111 // const wordList& localNames = s.regions();
112 //
113 // wordList globalNames(localNames.size());
114 // forAll(localNames, regionI)
115 // {
116 // globalNames[regionI] = s.name() + '_' + localNames[regionI];
117 // }
118 //
119 // // See if dictionary provides any global region names.
120 // if (dict.found("regions"))
121 // {
122 // const dictionary& regionsDict = dict.subDict("regions");
123 //
124 // forAllConstIter(dictionary, regionsDict, iter)
125 // {
126 // const word& key = iter().keyword();
127 //
128 // if (regionsDict.isDict(key))
129 // {
130 // // Get the dictionary for region iter.key()
131 // const dictionary& regionDict = regionsDict.subDict(key);
132 //
133 // label index = findIndex(localNames, key);
134 //
135 // if (index == -1)
136 // {
137 // FatalErrorInFunction
138 // << "Unknown region name " << key
139 // << " for surface " << s.name() << endl
140 // << "Valid region names are " << localNames
141 // << exit(FatalError);
142 // }
143 //
144 // globalNames[index] = word(regionDict.lookup("name"));
145 // }
146 // }
147 // }
148 //
149 // // Now globalNames contains the names of the regions.
150 // Info<< "Surface:" << s.name() << " has regions:"
151 // << endl;
152 // forAll(globalNames, regionI)
153 // {
154 // Info<< " " << globalNames[regionI] << endl;
155 // }
156 //
157 // // Create reverse lookup
158 // forAll(globalNames, regionI)
159 // {
160 // regionNames_.insert
161 // (
162 // globalNames[regionI],
163 // labelPair(surfI, regionI)
164 // );
165 // }
166 // }
167 //}
168 
169 
171 (
172  const IOobject& io,
173  const dictionary& topDict,
174  const bool singleRegionName
175 )
176 :
177  PtrList<searchableSurface>(topDict.size()),
178  names_(topDict.size()),
179  regionNames_(topDict.size()),
180  allSurfaces_(identityMap(topDict.size()))
181 {
182  label surfI = 0;
183  forAllConstIter(dictionary, topDict, iter)
184  {
185  const word& key = iter().keyword();
186 
187  if (!topDict.isDict(key))
188  {
190  << "Found non-dictionary entry " << iter()
191  << " in top-level dictionary " << topDict
192  << exit(FatalError);
193  }
194 
195  const dictionary& dict = topDict.subDict(key);
196 
197  names_[surfI] = key;
198  dict.readIfPresent("name", names_[surfI]);
199 
200  // Make IOobject with correct name
201  autoPtr<IOobject> namedIO(io.clone());
202  // Note: we would like to e.g. register triSurface 'sphere.stl' as
203  // 'sphere'. Unfortunately
204  // no support for having object read from different location than
205  // their object name. Maybe have stlTriSurfaceMesh which appends .stl
206  // when reading/writing?
207  namedIO().rename(key); // names_[surfI]
208 
209  // Create and hook surface
210  set
211  (
212  surfI,
214  (
215  dict.lookup("type"),
216  namedIO(),
217  dict
218  )
219  );
220  const searchableSurface& s = operator[](surfI);
221 
222  // Construct default region names by prepending surface name
223  // to region name.
224  const wordList& localNames = s.regions();
225 
226  wordList& rNames = regionNames_[surfI];
227  rNames.setSize(localNames.size());
228 
229  if (singleRegionName && localNames.size() == 1)
230  {
231  rNames[0] = names_[surfI];
232  }
233  else
234  {
235  forAll(localNames, regionI)
236  {
237  rNames[regionI] = names_[surfI] + '_' + localNames[regionI];
238  }
239  }
240 
241  // See if dictionary provides any global region names.
242  if (dict.found("regions"))
243  {
244  const dictionary& regionsDict = dict.subDict("regions");
245 
246  forAllConstIter(dictionary, regionsDict, iter)
247  {
248  const word& key = iter().keyword();
249 
250  if (regionsDict.isDict(key))
251  {
252  // Get the dictionary for region iter.keyword()
253  const dictionary& regionDict = regionsDict.subDict(key);
254 
255  label index = findIndex(localNames, key);
256 
257  if (index == -1)
258  {
260  << "Unknown region name " << key
261  << " for surface " << s.name() << endl
262  << "Valid region names are " << localNames
263  << exit(FatalError);
264  }
265 
266  rNames[index] = word(regionDict.lookup("name"));
267  }
268  }
269  }
270 
271  surfI++;
272  }
273 
274  // Trim (not really necessary since we don't allow non-dictionary entries)
276  names_.setSize(surfI);
277  regionNames_.setSize(surfI);
278  allSurfaces_.setSize(surfI);
279 }
280 
281 
282 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
283 
285 (
286  const word& wantedName
287 ) const
288 {
289  return findIndex(names_, wantedName);
290 }
291 
292 
294 (
295  const word& surfaceName,
296  const word& regionName
297 ) const
298 {
299  label surfaceIndex = findSurfaceID(surfaceName);
300 
301  return findIndex(this->operator[](surfaceIndex).regions(), regionName);
302 }
303 
304 
305 // Find any intersection
307 (
308  const pointField& start,
309  const pointField& end,
310  labelList& hitSurfaces,
311  List<pointIndexHit>& hitInfo
312 ) const
313 {
315  (
316  *this,
317  allSurfaces_,
318  start,
319  end,
320  hitSurfaces,
321  hitInfo
322  );
323 }
324 
325 
327 (
328  const pointField& start,
329  const pointField& end,
330  labelListList& hitSurfaces,
331  List<List<pointIndexHit>>& hitInfo
332 ) const
333 {
335  (
336  *this,
337  allSurfaces_,
338  start,
339  end,
340  hitSurfaces,
341  hitInfo
342  );
343 }
344 
345 
346 //Find intersections of edge nearest to both endpoints.
348 (
349  const pointField& start,
350  const pointField& end,
351  labelList& surface1,
352  List<pointIndexHit>& hit1,
353  labelList& surface2,
354  List<pointIndexHit>& hit2
355 ) const
356 {
358  (
359  *this,
360  allSurfaces_,
361  start,
362  end,
363  surface1,
364  hit1,
365  surface2,
366  hit2
367  );
368 }
369 
370 
372 (
373  const pointField& samples,
374  const scalarField& nearestDistSqr,
375  labelList& nearestSurfaces,
376  List<pointIndexHit>& nearestInfo
377 ) const
378 {
380  (
381  *this,
382  allSurfaces_,
383  samples,
384  nearestDistSqr,
385  nearestSurfaces,
386  nearestInfo
387  );
388 }
389 
390 
392 (
393  const pointField& samples,
394  const scalarField& nearestDistSqr,
395  const labelList& regionIndices,
396  labelList& nearestSurfaces,
397  List<pointIndexHit>& nearestInfo
398 ) const
399 {
401  (
402  *this,
403  allSurfaces_,
404  samples,
405  nearestDistSqr,
406  regionIndices,
407  nearestSurfaces,
408  nearestInfo
409  );
410 }
411 
412 
414 {
416  (
417  *this,
418  allSurfaces_
419  );
420 }
421 
422 
423 bool Foam::searchableSurfaces::checkClosed(const bool report) const
424 {
425  if (report)
426  {
427  Info<< "Checking for closedness." << endl;
428  }
429 
430  bool hasError = false;
431 
432  forAll(*this, surfI)
433  {
434  if (!operator[](surfI).hasVolumeType())
435  {
436  hasError = true;
437 
438  if (report)
439  {
440  Info<< " " << names()[surfI]
441  << " : not closed" << endl;
442  }
443 
444  if (isA<triSurface>(operator[](surfI)))
445  {
446  const triSurface& s = dynamic_cast<const triSurface&>
447  (
448  operator[](surfI)
449  );
450  const labelListList& edgeFaces = s.edgeFaces();
451 
452  label nSingleEdges = 0;
453  forAll(edgeFaces, edgeI)
454  {
455  if (edgeFaces[edgeI].size() == 1)
456  {
457  nSingleEdges++;
458  }
459  }
460 
461  label nMultEdges = 0;
462  forAll(edgeFaces, edgeI)
463  {
464  if (edgeFaces[edgeI].size() > 2)
465  {
466  nMultEdges++;
467  }
468  }
469 
470  if (report && (nSingleEdges != 0 || nMultEdges != 0))
471  {
472  Info<< " connected to one face : "
473  << nSingleEdges << nl
474  << " connected to >2 faces : "
475  << nMultEdges << endl;
476  }
477  }
478  }
479  }
480 
481  if (report)
482  {
483  Info<< endl;
484  }
485 
486  return returnReduce(hasError, orOp<bool>());
487 }
488 
489 
491 {
492  if (report)
493  {
494  Info<< "Checking for normal orientation." << endl;
495  }
496 
497  bool hasError = false;
498 
499  forAll(*this, surfI)
500  {
501  if (isA<triSurface>(operator[](surfI)))
502  {
503  const triSurface& s = dynamic_cast<const triSurface&>
504  (
505  operator[](surfI)
506  );
507 
508  labelHashSet borderEdge(s.size()/1000);
509  PatchTools::checkOrientation(s, false, &borderEdge);
510  // Colour all faces into zones using borderEdge
511  labelList normalZone;
512  label nZones = PatchTools::markZones(s, borderEdge, normalZone);
513  if (nZones > 1)
514  {
515  hasError = true;
516 
517  if (report)
518  {
519  Info<< " " << names()[surfI]
520  << " : has multiple orientation zones ("
521  << nZones << ")" << endl;
522  }
523  }
524  }
525  }
526 
527  if (report)
528  {
529  Info<< endl;
530  }
531 
532  return returnReduce(hasError, orOp<bool>());
533 }
534 
535 
537 (
538  const scalar maxRatio,
539  const bool report
540 ) const
541 {
542  if (report)
543  {
544  Info<< "Checking for size." << endl;
545  }
546 
547  bool hasError = false;
548 
549  forAll(*this, i)
550  {
551  const boundBox& bb = operator[](i).bounds();
552 
553  for (label j = i+1; j < size(); j++)
554  {
555  scalar ratio = bb.mag()/operator[](j).bounds().mag();
556 
557  if (ratio > maxRatio || ratio < 1.0/maxRatio)
558  {
559  hasError = true;
560 
561  if (report)
562  {
563  Info<< " " << names()[i]
564  << " bounds differ from " << names()[j]
565  << " by more than a factor 100:" << nl
566  << " bounding box : " << bb << nl
567  << " bounding box : " << operator[](j).bounds()
568  << endl;
569  }
570  break;
571  }
572  }
573  }
574 
575  if (report)
576  {
577  Info<< endl;
578  }
579 
580  return returnReduce(hasError, orOp<bool>());
581 }
582 
583 
585 (
586  const scalar tolerance,
587  const bool report
588 ) const
589 {
590  if (report)
591  {
592  Info<< "Checking for intersection." << endl;
593  }
594 
595  // cpuTime timer;
596 
597  bool hasError = false;
598 
599  forAll(*this, i)
600  {
601  if (isA<triSurfaceMesh>(operator[](i)))
602  {
603  const triSurfaceMesh& s0 = dynamic_cast<const triSurfaceMesh&>
604  (
605  operator[](i)
606  );
607  const edgeList& edges0 = s0.edges();
608  const pointField& localPoints0 = s0.localPoints();
609 
610  // Collect all the edge vectors
611  pointField start(edges0.size());
612  pointField end(edges0.size());
613  forAll(edges0, edgeI)
614  {
615  const edge& e = edges0[edgeI];
616  start[edgeI] = localPoints0[e[0]];
617  end[edgeI] = localPoints0[e[1]];
618  }
619 
620  // Test all other surfaces for intersection
621  forAll(*this, j)
622  {
623  List<pointIndexHit> hits;
624 
625  if (i == j)
626  {
627  // Slightly shorten the edges to prevent finding lots of
628  // intersections. The fast triangle intersection routine
629  // has problems with rays passing through a point of the
630  // triangle.
631  vectorField d
632  (
633  max(tolerance, 10*s0.tolerance())
634  *(end-start)
635  );
636  start += d;
637  end -= d;
638  }
639 
640  operator[](j).findLineAny(start, end, hits);
641 
642  label nHits = 0;
643  DynamicField<point> intersections(edges0.size()/100);
644  DynamicField<scalar> intersectionEdge(intersections.capacity());
645 
646  forAll(hits, edgeI)
647  {
648  if
649  (
650  hits[edgeI].hit()
651  && (i != j || !connected(s0, edgeI, hits[edgeI]))
652  )
653  {
654  intersections.append(hits[edgeI].hitPoint());
655  intersectionEdge.append(1.0*edgeI);
656  nHits++;
657  }
658  }
659 
660  if (nHits > 0)
661  {
662  if (report)
663  {
664  Info<< " " << names()[i]
665  << " intersects " << names()[j]
666  << " at " << nHits
667  << " locations."
668  << endl;
669 
670  const fileName fName
671  (
672  names()[i] + '_' + names()[j] + "_edgeIndex.vtk"
673  );
674 
675  Info<< " Writing intersection locations to "
676  << fName << endl;
677 
679  (
680  fName,
681  names()[i] + '_' + names()[j],
682  false,
683  intersections,
684  identityMap(intersections.size()),
685  edgeList(),
686  faceList(),
687  "edgeIndex",
688  true,
689  intersectionEdge
690  );
691  }
692 
693  hasError = true;
694  break;
695  }
696  }
697  }
698  }
699 
700  if (report)
701  {
702  Info<< endl;
703  }
704 
705  return returnReduce(hasError, orOp<bool>());
706 }
707 
708 
710 (
711  const scalar minQuality,
712  const bool report
713 ) const
714 {
715  if (report)
716  {
717  Info<< "Checking for triangle quality." << endl;
718  }
719 
720  bool hasError = false;
721 
722  forAll(*this, surfI)
723  {
724  if (isA<triSurface>(operator[](surfI)))
725  {
726  const triSurface& s = dynamic_cast<const triSurface&>
727  (
728  operator[](surfI)
729  );
730 
731  label nBadTris = 0;
732  forAll(s, facei)
733  {
734  const labelledTri& f = s[facei];
735 
736  scalar q = triPointRef
737  (
738  s.points()[f[0]],
739  s.points()[f[1]],
740  s.points()[f[2]]
741  ).quality();
742 
743  if (q < minQuality)
744  {
745  nBadTris++;
746  }
747  }
748 
749  if (nBadTris > 0)
750  {
751  hasError = true;
752 
753  if (report)
754  {
755  Info<< " " << names()[surfI]
756  << " : has " << nBadTris << " bad quality triangles "
757  << " (quality < " << minQuality << ")" << endl;
758  }
759  }
760  }
761  }
762 
763  if (report)
764  {
765  Info<< endl;
766  }
767 
768  return returnReduce(hasError, orOp<bool>());
769 
770 }
771 
772 
773 // Check all surfaces. Return number of failures.
775 (
776  const bool report
777 ) const
778 {
779  label noFailedChecks = 0;
780 
781  if (checkClosed(report))
782  {
783  noFailedChecks++;
784  }
785 
786  if (checkNormalOrientation(report))
787  {
788  noFailedChecks++;
789  }
790  return noFailedChecks;
791 }
792 
793 
795 (
796  const scalar maxRatio,
797  const scalar tol,
798  const scalar minQuality,
799  const bool report
800 ) const
801 {
802  label noFailedChecks = 0;
803 
804  if (maxRatio > 0 && checkSizes(maxRatio, report))
805  {
806  noFailedChecks++;
807  }
808 
809  if (checkIntersection(tol, report))
810  {
811  noFailedChecks++;
812  }
813 
814  if (checkQuality(minQuality, report))
815  {
816  noFailedChecks++;
817  }
818 
819  return noFailedChecks;
820 }
821 
822 
824 (
825  const List<wordList>& patchTypes,
826  Ostream& os
827 ) const
828 {
829  Info<< "Statistics." << endl;
830  forAll(*this, surfI)
831  {
832  Info<< " " << names()[surfI] << ':' << endl;
833 
834  const searchableSurface& s = operator[](surfI);
835 
836  Info<< " type : " << s.type() << nl
837  << " size : " << s.globalSize() << nl;
838  if (isA<triSurfaceMesh>(s))
839  {
840  const triSurfaceMesh& ts = dynamic_cast<const triSurfaceMesh&>(s);
841  Info<< " edges : " << ts.nEdges() << nl
842  << " points : " << ts.points()().size() << nl;
843  }
844  Info<< " bounds : " << s.bounds() << nl
845  << " closed : " << Switch(s.hasVolumeType()) << endl;
846 
847  if (patchTypes.size() && patchTypes[surfI].size() >= 1)
848  {
849  wordList unique(HashSet<word>(patchTypes[surfI]).sortedToc());
850  Info<< " patches : ";
851  forAll(unique, i)
852  {
853  Info<< unique[i];
854  if (i < unique.size()-1)
855  {
856  Info<< ',';
857  }
858  }
859  Info<< endl;
860  }
861  }
862  Info<< endl;
863 }
864 
865 
866 // * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * * //
867 
869 (
870  const word& surfName
871 ) const
872 {
873  const label surfI = findSurfaceID(surfName);
874 
875  if (surfI < 0)
876  {
878  << "Surface named " << surfName << " not found." << nl
879  << "Available surface names: " << names_ << endl
880  << abort(FatalError);
881  }
882 
883  return operator[](surfI);
884 }
885 
886 
888 (
889  const word& surfName
890 )
891 {
892  const label surfI = findSurfaceID(surfName);
893 
894  if (surfI < 0)
895  {
897  << "Surface named " << surfName << " not found." << nl
898  << "Available surface names: " << names_ << endl
899  << abort(FatalError);
900  }
901 
902  return operator[](surfI);
903 }
904 
905 
906 // ************************************************************************* //
Various functions to operate on Lists.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
Dynamically sized Field.
Definition: DynamicField.H:72
DynamicField< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
label capacity() const
Size of the underlying storage.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
autoPtr< IOobject > clone() const
Clone.
Definition: IOobject.H:283
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
void setSize(const label)
Reset size of List.
Definition: List.C:281
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
static bool checkOrientation(const PrimitivePatch< FaceList, PointField > &, const bool report=false, labelHashSet *marked=0)
Check for orientation issues.
static label markZones(const PrimitivePatch< FaceList, PointField > &, const BoolListType &borderEdge, labelList &faceZone)
Size and fills faceZone with zone of face.
label nEdges() const
Return number of edges in patch.
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
const Field< PointType > & localPoints() const
Return pointField of points in patch.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
bool set(const label) const
Is element set.
Definition: PtrListI.H:65
void setSize(const label)
Reset size of PtrList. If extending the PtrList, new entries are.
Definition: PtrList.C:131
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:61
friend Ostream & operator(Ostream &, const UPtrList< T > &)
Write UPtrList to Ostream.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
A bounding box defined in terms of the points at its extremities.
Definition: boundBox.H:59
scalar mag() const
The magnitude of the bounding box span.
Definition: boundBoxI.H:90
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:860
bool isDict(const word &) const
Check if entry is a sub-dictionary.
Definition: dictionary.C:952
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:998
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:61
A class for handling file names.
Definition: fileName.H:82
Triangle with additional region number.
Definition: labelledTri.H:60
Base class of (analytical or triangulated) surface. Encapsulates all the search routines....
static autoPtr< searchableSurface > New(const word &surfaceType, const IOobject &io, const dictionary &dict)
Return a reference to the selected searchableSurface.
static boundBox bounds(const PtrList< searchableSurface > &allSurfaces, const labelList &surfacesToTest)
Find the boundBox of the selected surfaces.
static void findNearestIntersection(const PtrList< searchableSurface > &allSurfaces, const labelList &surfacesToTest, const pointField &start, const pointField &end, labelList &surface1, List< pointIndexHit > &hit1, labelList &surface2, List< pointIndexHit > &hit2)
Find intersections of edge nearest to both endpoints.
static void findAnyIntersection(const PtrList< searchableSurface > &, const labelList &surfacesToTest, const pointField &start, const pointField &end, labelList &surfaces, List< pointIndexHit > &)
Find any intersection. Return hit point information and.
static void findAllIntersections(const PtrList< searchableSurface > &, const labelList &surfacesToTest, const pointField &start, const pointField &end, labelListList &surfaces, List< List< pointIndexHit >> &surfaceHits)
Find all intersections in order from start to end. Returns for.
static void findNearest(const PtrList< searchableSurface > &, const labelList &surfacesToTest, const pointField &, const scalarField &nearestDistSqr, labelList &surfaces, List< pointIndexHit > &)
Find nearest. Return -1 (and a miss()) or surface and nearest.
Container for searchableSurfaces.
bool checkQuality(const scalar minQuality, const bool report) const
Check triangle quality.
label checkGeometry(const scalar maxRatio, const scalar tolerance, const scalar minQuality, const bool report) const
All geometric checks. Return number of failed checks.
boundBox bounds() const
Calculate bounding box.
label findSurfaceRegionID(const word &surfaceName, const word &regionName) const
bool checkClosed(const bool report) const
Are all surfaces closed and manifold.
bool checkSizes(const scalar maxRatio, const bool report) const
Are all bounding boxes of similar size.
const searchableSurface & operator[](const word &) const
Return const reference to searchableSurface by name.
void findAllIntersections(const pointField &start, const pointField &end, labelListList &surfaces, List< List< pointIndexHit >> &) const
Find all intersections in order from start to end. Returns for.
label checkTopology(const bool report) const
All topological checks. Return number of failed checks.
void findAnyIntersection(const pointField &start, const pointField &end, labelList &surfaces, List< pointIndexHit > &) const
Find any intersection. Return hit point information and.
void writeStats(const List< wordList > &, Ostream &) const
Write some stats.
bool checkNormalOrientation(const bool report) const
Are all (triangulated) surfaces consistent normal orientation.
bool checkIntersection(const scalar tol, const bool report) const
Do surfaces self-intersect or intersect others.
searchableSurfaces(const label)
Construct with length specified. Fill later.
void findNearest(const pointField &, const scalarField &nearestDistSqr, labelList &surfaces, List< pointIndexHit > &) const
Find nearest. Return -1 (and a miss()) or surface and nearest.
label findSurfaceID(const word &name) const
Find index of surface. Return -1 if not found.
void findNearestIntersection(const pointField &start, const pointField &end, labelList &surface1, List< pointIndexHit > &hit1, labelList &surface2, List< pointIndexHit > &hit2) const
A surface geometry formed of discrete facets, e.g. triangles and/or quadrilaterals,...
virtual tmp< pointField > points() const
Get the points that define the surface.
scalar tolerance() const
Return tolerance to use in searches.
Triangulated surface description with patch information.
Definition: triSurface.H:69
A class for handling words, derived from string.
Definition: word.H:62
Foam::word regionName
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.name(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
void write(const fileName &file, const word &title, const bool binary, const PointField &points, const VertexList &vertices, const LineList &lines, const FaceList &faces, const wordList &fieldNames, const boolList &fieldIsPointValues, const UPtrList< const Field< label >> &fieldLabelValues #define FieldTypeValuesConstArg(Type, nullArg))
Write VTK polygonal data to a file. Takes a PtrList of fields of labels and.
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
const doubleScalar e
Definition: doubleScalar.H:105
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
PointIndexHit< point > pointIndexHit
Definition: pointIndexHit.H:42
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
errorManip< error > abort(error &err)
Definition: errorManip.H:131
messageStream Info
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
defineTypeNameAndDebug(combustionModel, 0)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
triangle< point, const point & > triPointRef
Definition: triPointRef.H:44
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
error FatalError
labelList identityMap(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
List< face > faceList
Definition: faceListFwd.H:43
static const char nl
Definition: Ostream.H:260
wordList patchTypes(nPatches)
face triFace(3)
labelList f(nPoints)
dictionary dict
scalarField samples(nIntervals, 0)