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