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