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-2015 OpenFOAM Foundation
6  \\/ M anipulation | Copyright (C) 2015 OpenCFD Ltd.
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 //- Is edge connected to triangle
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_(identity(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_(identity(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_(identity(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 
326 //- Find all intersections in order from start to end. Returns for
327 // every hit the surface and the hit info.
329 (
330  const pointField& start,
331  const pointField& end,
332  labelListList& hitSurfaces,
333  List<List<pointIndexHit> >& hitInfo
334 ) const
335 {
337  (
338  *this,
339  allSurfaces_,
340  start,
341  end,
342  hitSurfaces,
343  hitInfo
344  );
345 }
346 
347 
348 //Find intersections of edge nearest to both endpoints.
350 (
351  const pointField& start,
352  const pointField& end,
353  labelList& surface1,
354  List<pointIndexHit>& hit1,
355  labelList& surface2,
356  List<pointIndexHit>& hit2
357 ) const
358 {
360  (
361  *this,
362  allSurfaces_,
363  start,
364  end,
365  surface1,
366  hit1,
367  surface2,
368  hit2
369  );
370 }
371 
372 
373 // Find nearest. Return -1 or nearest point
375 (
376  const pointField& samples,
377  const scalarField& nearestDistSqr,
378  labelList& nearestSurfaces,
379  List<pointIndexHit>& nearestInfo
380 ) const
381 {
383  (
384  *this,
385  allSurfaces_,
386  samples,
387  nearestDistSqr,
388  nearestSurfaces,
389  nearestInfo
390  );
391 }
392 
393 
394 // Find nearest. Return -1 or nearest point
396 (
397  const labelListList& regionIndices,
398  const pointField& samples,
399  const scalarField& nearestDistSqr,
400  labelList& nearestSurfaces,
401  List<pointIndexHit>& nearestInfo
402 ) const
403 {
405  (
406  *this,
407  allSurfaces_,
408  regionIndices,
409 
410  samples,
411  nearestDistSqr,
412 
413  nearestSurfaces,
414  nearestInfo
415  );
416 }
417 
418 //- Calculate bounding box
420 {
422  (
423  *this,
424  allSurfaces_
425  );
426 }
427 
428 
429 //- Calculate point which is on a set of surfaces.
431 (
432  const scalar initDistSqr,
433  const scalar convergenceDistSqr,
434  const point& start
435 ) const
436 {
438  (
439  *this,
440  allSurfaces_,
441  initDistSqr,
442  convergenceDistSqr,
443  start
444  );
445 }
446 
447 
448 bool Foam::searchableSurfaces::checkClosed(const bool report) const
449 {
450  if (report)
451  {
452  Info<< "Checking for closedness." << endl;
453  }
454 
455  bool hasError = false;
456 
457  forAll(*this, surfI)
458  {
459  if (!operator[](surfI).hasVolumeType())
460  {
461  hasError = true;
462 
463  if (report)
464  {
465  Info<< " " << names()[surfI]
466  << " : not closed" << endl;
467  }
468 
469  if (isA<triSurface>(operator[](surfI)))
470  {
471  const triSurface& s = dynamic_cast<const triSurface&>
472  (
473  operator[](surfI)
474  );
475  const labelListList& edgeFaces = s.edgeFaces();
476 
477  label nSingleEdges = 0;
478  forAll(edgeFaces, edgeI)
479  {
480  if (edgeFaces[edgeI].size() == 1)
481  {
482  nSingleEdges++;
483  }
484  }
485 
486  label nMultEdges = 0;
487  forAll(edgeFaces, edgeI)
488  {
489  if (edgeFaces[edgeI].size() > 2)
490  {
491  nMultEdges++;
492  }
493  }
494 
495  if (report && (nSingleEdges != 0 || nMultEdges != 0))
496  {
497  Info<< " connected to one face : "
498  << nSingleEdges << nl
499  << " connected to >2 faces : "
500  << nMultEdges << endl;
501  }
502  }
503  }
504  }
505 
506  if (report)
507  {
508  Info<< endl;
509  }
510 
511  return returnReduce(hasError, orOp<bool>());
512 }
513 
514 
516 {
517  if (report)
518  {
519  Info<< "Checking for normal orientation." << endl;
520  }
521 
522  bool hasError = false;
523 
524  forAll(*this, surfI)
525  {
526  if (isA<triSurface>(operator[](surfI)))
527  {
528  const triSurface& s = dynamic_cast<const triSurface&>
529  (
530  operator[](surfI)
531  );
532 
533  labelHashSet borderEdge(s.size()/1000);
534  PatchTools::checkOrientation(s, false, &borderEdge);
535  // Colour all faces into zones using borderEdge
536  labelList normalZone;
537  label nZones = PatchTools::markZones(s, borderEdge, normalZone);
538  if (nZones > 1)
539  {
540  hasError = true;
541 
542  if (report)
543  {
544  Info<< " " << names()[surfI]
545  << " : has multiple orientation zones ("
546  << nZones << ")" << endl;
547  }
548  }
549  }
550  }
551 
552  if (report)
553  {
554  Info<< endl;
555  }
556 
557  return returnReduce(hasError, orOp<bool>());
558 }
559 
560 
562 (
563  const scalar maxRatio,
564  const bool report
565 ) const
566 {
567  if (report)
568  {
569  Info<< "Checking for size." << endl;
570  }
571 
572  bool hasError = false;
573 
574  forAll(*this, i)
575  {
576  const boundBox& bb = operator[](i).bounds();
577 
578  for (label j = i+1; j < size(); j++)
579  {
580  scalar ratio = bb.mag()/operator[](j).bounds().mag();
581 
582  if (ratio > maxRatio || ratio < 1.0/maxRatio)
583  {
584  hasError = true;
585 
586  if (report)
587  {
588  Info<< " " << names()[i]
589  << " bounds differ from " << names()[j]
590  << " by more than a factor 100:" << nl
591  << " bounding box : " << bb << nl
592  << " bounding box : " << operator[](j).bounds()
593  << endl;
594  }
595  break;
596  }
597  }
598  }
599 
600  if (report)
601  {
602  Info<< endl;
603  }
604 
605  return returnReduce(hasError, orOp<bool>());
606 }
607 
608 
610 (
611  const scalar tolerance,
612  const autoPtr<writer<scalar> >& setWriter,
613  const bool report
614 ) const
615 {
616  if (report)
617  {
618  Info<< "Checking for intersection." << endl;
619  }
620 
621  //cpuTime timer;
622 
623  bool hasError = false;
624 
625  forAll(*this, i)
626  {
627  if (isA<triSurfaceMesh>(operator[](i)))
628  {
629  const triSurfaceMesh& s0 = dynamic_cast<const triSurfaceMesh&>
630  (
631  operator[](i)
632  );
633  const edgeList& edges0 = s0.edges();
634  const pointField& localPoints0 = s0.localPoints();
635 
636  // Collect all the edge vectors
637  pointField start(edges0.size());
638  pointField end(edges0.size());
639  forAll(edges0, edgeI)
640  {
641  const edge& e = edges0[edgeI];
642  start[edgeI] = localPoints0[e[0]];
643  end[edgeI] = localPoints0[e[1]];
644  }
645 
646  // Test all other surfaces for intersection
647  forAll(*this, j)
648  {
649  List<pointIndexHit> hits;
650 
651  if (i == j)
652  {
653  // Slightly shorten the edges to prevent finding lots of
654  // intersections. The fast triangle intersection routine
655  // has problems with rays passing through a point of the
656  // triangle.
657  vectorField d
658  (
659  max(tolerance, 10*s0.tolerance())
660  *(end-start)
661  );
662  start += d;
663  end -= d;
664  }
665 
666  operator[](j).findLineAny(start, end, hits);
667 
668  label nHits = 0;
669  DynamicField<point> intersections(edges0.size()/100);
670  DynamicField<scalar> intersectionEdge(intersections.capacity());
671 
672  forAll(hits, edgeI)
673  {
674  if
675  (
676  hits[edgeI].hit()
677  && (i != j || !connected(s0, edgeI, hits[edgeI]))
678  )
679  {
680  intersections.append(hits[edgeI].hitPoint());
681  intersectionEdge.append(1.0*edgeI);
682  nHits++;
683  }
684  }
685 
686  if (nHits > 0)
687  {
688  if (report)
689  {
690  Info<< " " << names()[i]
691  << " intersects " << names()[j]
692  << " at " << nHits
693  << " locations."
694  << endl;
695 
696  //vtkSetWriter<scalar> setWriter;
697  if (setWriter.valid())
698  {
699  scalarField dist(mag(intersections));
700  coordSet track
701  (
702  names()[i] + '_' + names()[j],
703  "xyz",
704  intersections.xfer(),
705  dist
706  );
707  wordList valueSetNames(1, "edgeIndex");
708  List<const scalarField*> valueSets
709  (
710  1,
711  &intersectionEdge
712  );
713 
714  fileName fName
715  (
716  setWriter().getFileName(track, valueSetNames)
717  );
718  Info<< " Writing intersection locations to "
719  << fName << endl;
720  OFstream os
721  (
722  s0.searchableSurface::time().path()
723  /fName
724  );
725  setWriter().write
726  (
727  track,
728  valueSetNames,
729  valueSets,
730  os
731  );
732  }
733  }
734 
735  hasError = true;
736  break;
737  }
738  }
739  }
740  }
741 
742  if (report)
743  {
744  Info<< endl;
745  }
746 
747  return returnReduce(hasError, orOp<bool>());
748 }
749 
750 
752 (
753  const scalar minQuality,
754  const bool report
755 ) const
756 {
757  if (report)
758  {
759  Info<< "Checking for triangle quality." << endl;
760  }
761 
762  bool hasError = false;
763 
764  forAll(*this, surfI)
765  {
766  if (isA<triSurface>(operator[](surfI)))
767  {
768  const triSurface& s = dynamic_cast<const triSurface&>
769  (
770  operator[](surfI)
771  );
772 
773  label nBadTris = 0;
774  forAll(s, faceI)
775  {
776  const labelledTri& f = s[faceI];
777 
778  scalar q = triPointRef
779  (
780  s.points()[f[0]],
781  s.points()[f[1]],
782  s.points()[f[2]]
783  ).quality();
784 
785  if (q < minQuality)
786  {
787  nBadTris++;
788  }
789  }
790 
791  if (nBadTris > 0)
792  {
793  hasError = true;
794 
795  if (report)
796  {
797  Info<< " " << names()[surfI]
798  << " : has " << nBadTris << " bad quality triangles "
799  << " (quality < " << minQuality << ")" << endl;
800  }
801  }
802  }
803  }
804 
805  if (report)
806  {
807  Info<< endl;
808  }
809 
810  return returnReduce(hasError, orOp<bool>());
811 
812 }
813 
814 
815 // Check all surfaces. Return number of failures.
817 (
818  const bool report
819 ) const
820 {
821  label noFailedChecks = 0;
822 
823  if (checkClosed(report))
824  {
825  noFailedChecks++;
826  }
827 
828  if (checkNormalOrientation(report))
829  {
830  noFailedChecks++;
831  }
832  return noFailedChecks;
833 }
834 
835 
837 (
838  const scalar maxRatio,
839  const scalar tol,
840  const autoPtr<writer<scalar> >& setWriter,
841  const scalar minQuality,
842  const bool report
843 ) const
844 {
845  label noFailedChecks = 0;
846 
847  if (maxRatio > 0 && checkSizes(maxRatio, report))
848  {
849  noFailedChecks++;
850  }
851 
852  if (checkIntersection(tol, setWriter, report))
853  {
854  noFailedChecks++;
855  }
856 
857  if (checkQuality(minQuality, report))
858  {
859  noFailedChecks++;
860  }
861 
862  return noFailedChecks;
863 }
864 
865 
867 (
868  const List<wordList>& patchTypes,
869  Ostream& os
870 ) const
871 {
872  Info<< "Statistics." << endl;
873  forAll(*this, surfI)
874  {
875  Info<< " " << names()[surfI] << ':' << endl;
876 
877  const searchableSurface& s = operator[](surfI);
878 
879  Info<< " type : " << s.type() << nl
880  << " size : " << s.globalSize() << nl;
881  if (isA<triSurfaceMesh>(s))
882  {
883  const triSurfaceMesh& ts = dynamic_cast<const triSurfaceMesh&>(s);
884  Info<< " edges : " << ts.nEdges() << nl
885  << " points : " << ts.points()().size() << nl;
886  }
887  Info<< " bounds : " << s.bounds() << nl
888  << " closed : " << Switch(s.hasVolumeType()) << endl;
889 
890  if (patchTypes.size() && patchTypes[surfI].size() >= 1)
891  {
892  wordList unique(HashSet<word>(patchTypes[surfI]).sortedToc());
893  Info<< " patches : ";
894  forAll(unique, i)
895  {
896  Info<< unique[i];
897  if (i < unique.size()-1)
898  {
899  Info<< ',';
900  }
901  }
902  Info<< endl;
903  }
904  }
905  Info<< endl;
906 }
907 
908 
909 // * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * * //
910 
912 (
913  const word& surfName
914 ) const
915 {
916  const label surfI = findSurfaceID(surfName);
917 
918  if (surfI < 0)
919  {
921  << "Surface named " << surfName << " not found." << nl
922  << "Available surface names: " << names_ << endl
923  << abort(FatalError);
924  }
925 
926  return operator[](surfI);
927 }
928 
929 
931 (
932  const word& surfName
933 )
934 {
935  const label surfI = findSurfaceID(surfName);
936 
937  if (surfI < 0)
938  {
940  << "Surface named " << surfName << " not found." << nl
941  << "Available surface names: " << names_ << endl
942  << abort(FatalError);
943  }
944 
945  return operator[](surfI);
946 }
947 
948 
949 // ************************************************************************* //
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::Switch
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:60
Foam::PointIndexHit::index
label index() const
Return index.
Definition: PointIndexHit.H:111
Foam::searchableSurfaces::checkIntersection
bool checkIntersection(const scalar tol, const autoPtr< writer< scalar > > &, const bool report) const
Do surfaces self-intersect or intersect others.
Definition: searchableSurfaces.C:610
Foam::boundBox::mag
scalar mag() const
The magnitude of the bounding box span.
Definition: boundBoxI.H:90
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::fileName
A class for handling file names.
Definition: fileName.H:69
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:86
Foam::PrimitivePatch::edges
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
Definition: PrimitivePatchTemplate.C:212
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::findIndex
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurence of given element and return index,.
Foam::dictionary::readIfPresent
bool readIfPresent(const word &, T &, bool recursive=false, bool patternMatch=true) const
Find an entry if present, and assign to T.
Definition: dictionaryTemplates.C:94
Foam::searchableSurfaces::checkQuality
bool checkQuality(const scalar minQuality, const bool report) const
Check triangle quality.
Definition: searchableSurfaces.C:752
Foam::PrimitivePatch::localPoints
const Field< PointType > & localPoints() const
Return pointField of points in patch.
Definition: PrimitivePatchTemplate.C:432
Foam::searchableSurfaces::findNearestIntersection
void findNearestIntersection(const pointField &start, const pointField &end, labelList &surface1, List< pointIndexHit > &hit1, labelList &surface2, List< pointIndexHit > &hit2) const
Definition: searchableSurfaces.C:350
Foam::searchableSurfaces::connected
static bool connected(const triSurface &s, const label edgeI, const pointIndexHit &hit)
Is edge on face.
Definition: searchableSurfaces.C:46
Foam::PrimitivePatch::nEdges
label nEdges() const
Return number of edges in patch.
Definition: PrimitivePatchTemplate.H:299
Foam::searchableSurfaces::checkTopology
label checkTopology(const bool report) const
All topological checks. Return number of failed checks.
Definition: searchableSurfaces.C:817
Foam::dictionary::lookup
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:449
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
Foam::searchableSurfacesQueries::findAnyIntersection
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.
Definition: searchableSurfacesQueries.C:396
Foam::triSurfaceMesh
IOoject and searching on triSurface.
Definition: triSurfaceMesh.H:63
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::searchableSurfaces::writeStats
void writeStats(const List< wordList > &, Ostream &) const
Write some stats.
Definition: searchableSurfaces.C:867
Foam::searchableSurfacesQueries::findAllIntersections
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.
Definition: searchableSurfacesQueries.C:457
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::searchableSurfaces::searchableSurfaces
searchableSurfaces(const searchableSurfaces &)
Disallow default bitwise copy construct.
Foam::HashSet< label, Hash< label > >
searchableSurfacesQueries.H
Foam::dictionary::isDict
bool isDict(const word &) const
Check if entry is a sub-dictionary.
Definition: dictionary.C:600
patchTypes
wordList patchTypes(nPatches)
forAllConstIter
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
Foam::DynamicField
Dynamically sized Field.
Definition: DynamicField.H:49
samples
scalarField samples(nIntervals, 0)
searchableSurfaces.H
Foam::searchableSurfaces::facesIntersection
pointIndexHit facesIntersection(const scalar initialDistSqr, const scalar convergenceDistSqr, const point &start) const
Calculate point which is on a set of surfaces.
Definition: searchableSurfaces.C:431
regionName
Foam::word regionName
Definition: createNamedDynamicFvMesh.H:1
Foam::searchableSurfaces::findAnyIntersection
void findAnyIntersection(const pointField &start, const pointField &end, labelList &surfaces, List< pointIndexHit > &) const
Find any intersection. Return hit point information and.
Definition: searchableSurfaces.C:307
Foam::triSurfaceSearch::tolerance
scalar tolerance() const
Return tolerance to use in searches.
Definition: triSurfaceSearch.H:131
Foam::searchableSurfaces::bounds
boundBox bounds() const
Calculate bounding box.
Definition: searchableSurfaces.C:419
Foam::PointIndexHit
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:53
Foam::searchableSurfacesQueries::findNearestIntersection
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)
Definition: searchableSurfacesQueries.C:537
Foam::label
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
Foam::dictionary::found
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:304
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::searchableSurfaces::checkClosed
bool checkClosed(const bool report) const
Are all surfaces closed and manifold.
Definition: searchableSurfaces.C:448
Foam::triSurface
Triangulated surface description with patch information.
Definition: triSurface.H:57
Foam::searchableSurfaces::findNearest
void findNearest(const pointField &, const scalarField &nearestDistSqr, labelList &surfaces, List< pointIndexHit > &) const
Find nearest. Return -1 (and a miss()) or surface and nearest.
Definition: searchableSurfaces.C:375
Foam::orOp
Definition: ops.H:178
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Foam::searchableSurface
Base class of (analytical or triangulated) surface. Encapsulates all the search routines....
Definition: searchableSurface.H:66
Foam::DynamicField::append
DynamicField< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Foam::searchableSurfacesQueries::bounds
static boundBox bounds(const PtrList< searchableSurface > &allSurfaces, const labelList &surfacesToTest)
Find the boundBox of the selected surfaces.
Definition: searchableSurfacesQueries.C:834
Foam::PtrList
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:61
Foam::searchableSurfacesQueries::facesIntersection
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.
Definition: searchableSurfacesQueries.C:856
Foam::searchableSurfaces::findSurfaceRegionID
label findSurfaceRegionID(const word &surfaceName, const word &regionName) const
Definition: searchableSurfaces.C:294
Foam::searchableSurfaces::findSurfaceID
label findSurfaceID(const word &name) const
Find index of surface. Return -1 if not found.
Definition: searchableSurfaces.C:285
Foam::identity
labelList identity(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::DynamicField::capacity
label capacity() const
Size of the underlying storage.
Definition: DynamicFieldI.H:135
Foam::writer
Base class for graphics format writing. Entry points are.
Definition: writer.H:78
Foam::coordSet
Holds list of sampling positions.
Definition: coordSet.H:49
Foam::triSurfaceMesh::points
virtual tmp< pointField > points() const
Get the points that define the surface.
Definition: triSurfaceMesh.C:433
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
Foam::e
const double e
Elementary charge.
Definition: doubleFloat.H:94
Foam::searchableSurface::New
static autoPtr< searchableSurface > New(const word &surfaceType, const IOobject &io, const dictionary &dict)
Return a reference to the selected searchableSurface.
Definition: searchableSurface.C:40
s
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))
Foam::searchableSurfaces::checkSizes
bool checkSizes(const scalar maxRatio, const bool report) const
Are all bounding boxes of similar size.
Definition: searchableSurfaces.C:562
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::DynamicField::xfer
Xfer< List< T > > xfer()
Transfer contents to the Xfer container as a plain List.
Definition: DynamicFieldI.H:312
Foam::OFstream
Output to file stream.
Definition: OFstream.H:81
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::List::setSize
void setSize(const label)
Reset size of List.
Foam::autoPtr
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:117
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
DynamicField.H
Foam::triFace
A triangular face using a FixedList of labels corresponding to mesh vertices.
Definition: triFace.H:68
f
labelList f(nPoints)
Foam::Vector< scalar >
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::IOobject::clone
Foam::autoPtr< IOobject > clone() const
Clone.
Definition: IOobject.H:252
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: doubleFloat.H:94
Foam::PatchTools::markZones
static label markZones(const PrimitivePatch< Face, FaceList, PointField, PointType > &, const BoolListType &borderEdge, labelList &faceZone)
Size and fills faceZone with zone of face.
Foam::PtrList::setSize
void setSize(const label)
Reset size of PtrList. If extending the PtrList, new entries are.
Definition: PtrList.C:142
Foam::searchableSurfaces::checkNormalOrientation
bool checkNormalOrientation(const bool report) const
Are all (triangulated) surfaces consistent normal orientation.
Definition: searchableSurfaces.C:515
Foam::PtrList< searchableSurface >::operator
friend Ostream & operator(Ostream &, const PtrList< T > &)
Foam::labelledTri
Triangle with additional region number.
Definition: labelledTri.H:49
Foam::boundBox
A bounding box defined in terms of the points at its extremities.
Definition: boundBox.H:55
Foam::searchableSurfaces::findAllIntersections
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.
Definition: searchableSurfaces.C:329
ListOps.H
Various functions to operate on Lists.
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
Foam::dictionary::subDict
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:631
Foam::IOobject::path
fileName path() const
Return complete path.
Definition: IOobject.C:293
triSurfaceMesh.H
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::searchableSurfaces::checkGeometry
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.
Definition: searchableSurfaces.C:837
Foam::searchableSurfacesQueries::findNearest
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.
Definition: searchableSurfacesQueries.C:625
Foam::PatchTools::checkOrientation
static bool checkOrientation(const PrimitivePatch< Face, FaceList, PointField, PointType > &, const bool report=false, labelHashSet *marked=0)
Check for orientation issues.
Definition: PatchToolsCheck.C:40
Foam::triPointRef
triangle< point, const point & > triPointRef
Definition: triangle.H:75