triSurfaceMesh.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2021 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 "triSurfaceMesh.H"
27 #include "Random.H"
29 #include "EdgeMap.H"
30 #include "triSurfaceFields.H"
31 #include "Time.H"
32 #include "PatchTools.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38  defineTypeNameAndDebug(triSurfaceMesh, 0);
39  addToRunTimeSelectionTable(searchableSurface, triSurfaceMesh, dict);
40 }
41 
42 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
43 
44 Foam::fileName Foam::triSurfaceMesh::checkFile
45 (
46  const regIOobject& io,
47  const bool isGlobal
48 )
49 {
50  const fileName fName(io.filePath(typeName, isGlobal));
51 
52  if (fName.empty())
53  {
55  << "Cannot find triSurfaceMesh starting from "
56  << io.objectPath() << exit(FatalError);
57  }
58 
59  return fName;
60 }
61 
62 
63 Foam::fileName Foam::triSurfaceMesh::relativeFilePath
64 (
65  const regIOobject& io,
66  const fileName& f,
67  const bool isGlobal
68 )
69 {
70  fileName fName(f);
71  fName.expand();
72  if (!fName.isAbsolute())
73  {
74  // Is the specified file:
75  // - local to the cwd?
76  // - local to the case dir?
77  // - or just another name?
78  fName = fileHandler().filePath
79  (
80  isGlobal,
81  IOobject(io, fName),
83  );
84  }
85  return fName;
86 }
87 
88 
89 Foam::fileName Foam::triSurfaceMesh::checkFile
90 (
91  const regIOobject& io,
92  const dictionary& dict,
93  const bool isGlobal
94 )
95 {
96  fileName dictFName, fName;
97 
98  if (dict.readIfPresent("file", dictFName, false, false))
99  {
100  fName = relativeFilePath(io, dictFName, isGlobal);
101 
102  if (!exists(fName))
103  {
105  << "Cannot find triSurfaceMesh at " << io.path()/dictFName
106  << exit(FatalError);
107  }
108  }
109  else
110  {
111  fName = io.filePath(typeName, isGlobal);
112 
113  if (!exists(fName))
114  {
116  << "Cannot find triSurfaceMesh starting from "
117  << io.objectPath() << exit(FatalError);
118  }
119  }
120 
121  return fName;
122 }
123 
124 
125 bool Foam::triSurfaceMesh::addFaceToEdge
126 (
127  const edge& e,
128  EdgeMap<label>& facesPerEdge
129 )
130 {
131  EdgeMap<label>::iterator eFnd = facesPerEdge.find(e);
132  if (eFnd != facesPerEdge.end())
133  {
134  if (eFnd() == 2)
135  {
136  return false;
137  }
138  eFnd()++;
139  }
140  else
141  {
142  facesPerEdge.insert(e, 1);
143  }
144  return true;
145 }
146 
147 
148 bool Foam::triSurfaceMesh::isSurfaceClosed() const
149 {
150  const pointField& pts = triSurface::points();
151 
152  // Construct pointFaces. Let's hope surface has compact point
153  // numbering ...
154  labelListList pointFaces;
155  invertManyToMany(pts.size(), *this, pointFaces);
156 
157  // Loop over all faces surrounding point. Count edges emanating from point.
158  // Every edge should be used by two faces exactly.
159  // To prevent doing work twice per edge only look at edges to higher
160  // point
161  EdgeMap<label> facesPerEdge(100);
162  forAll(pointFaces, pointi)
163  {
164  const labelList& pFaces = pointFaces[pointi];
165 
166  facesPerEdge.clear();
167  forAll(pFaces, i)
168  {
169  const triSurface::FaceType& f = triSurface::operator[](pFaces[i]);
170  label fp = findIndex(f, pointi);
171 
172  // Something weird: if I expand the code of addFaceToEdge in both
173  // below instances it gives a segmentation violation on some
174  // surfaces. Compiler (4.3.2) problem?
175 
176 
177  // Forward edge
178  label nextPointi = f[f.fcIndex(fp)];
179 
180  if (nextPointi > pointi)
181  {
182  bool okFace = addFaceToEdge
183  (
184  edge(pointi, nextPointi),
185  facesPerEdge
186  );
187 
188  if (!okFace)
189  {
190  return false;
191  }
192  }
193  // Reverse edge
194  label prevPointi = f[f.rcIndex(fp)];
195 
196  if (prevPointi > pointi)
197  {
198  bool okFace = addFaceToEdge
199  (
200  edge(pointi, prevPointi),
201  facesPerEdge
202  );
203 
204  if (!okFace)
205  {
206  return false;
207  }
208  }
209  }
210 
211  // Check for any edges used only once.
212  forAllConstIter(EdgeMap<label>, facesPerEdge, iter)
213  {
214  if (iter() != 2)
215  {
216  return false;
217  }
218  }
219  }
220 
221  return true;
222 }
223 
224 
225 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
226 
228 :
229  searchableSurface(io),
231  (
232  IOobject
233  (
234  io.name(),
235  io.instance(),
236  io.local(),
237  io.db(),
238  io.readOpt(),
239  io.writeOpt(),
240  false // searchableSurface already registered under name
241  )
242  ),
243  triSurface(s),
245  minQuality_(-1),
246  surfaceClosed_(-1)
247 {
248  const pointField& pts = triSurface::points();
249 
250  bounds() = boundBox(pts);
251 }
252 
253 
255 :
256  // Find instance for triSurfaceMesh
257  searchableSurface(io),
258  // Reused found instance in objectRegistry
260  (
261  IOobject
262  (
263  io.name(),
265  io.local(),
266  io.db(),
267  io.readOpt(),
268  io.writeOpt(),
269  false // searchableSurface already registered under name
270  )
271  ),
272  triSurface(checkFile(static_cast<const searchableSurface&>(*this), true)),
273  triSurfaceRegionSearch(static_cast<const triSurface&>(*this)),
274  minQuality_(-1),
275  surfaceClosed_(-1)
276 {
277  const pointField& pts = triSurface::points();
278 
279  bounds() = boundBox(pts);
280 }
281 
282 
284 (
285  const IOobject& io,
286  const dictionary& dict
287 )
288 :
289  searchableSurface(io),
290  // Reused found instance in objectRegistry
292  (
293  IOobject
294  (
295  io.name(),
297  io.local(),
298  io.db(),
299  io.readOpt(),
300  io.writeOpt(),
301  false // searchableSurface already registered under name
302  )
303  ),
304  triSurface
305  (
306  checkFile(static_cast<const searchableSurface&>(*this), dict, true)
307  ),
308  triSurfaceRegionSearch(static_cast<const triSurface&>(*this), dict),
309  minQuality_(-1),
310  surfaceClosed_(-1)
311 {
312  // Reading from supplied file name instead of objectPath/filePath
313  if (dict.readIfPresent("file", fName_, false, false))
314  {
315  fName_ = relativeFilePath
316  (
317  static_cast<const searchableSurface&>(*this),
318  fName_,
319  true
320  );
321  }
322 
323  scalar scaleFactor = 0;
324 
325  // Allow rescaling of the surface points
326  // eg, CAD geometries are often done in millimeters
327  if (dict.readIfPresent("scale", scaleFactor) && scaleFactor > 0)
328  {
329  Info<< searchableSurface::name() << " : using scale " << scaleFactor
330  << endl;
331  triSurface::scalePoints(scaleFactor);
332  }
333 
334  const pointField& pts = triSurface::points();
335 
336  bounds() = boundBox(pts);
337 
338  // Have optional minimum quality for normal calculation
339  if (dict.readIfPresent("minQuality", minQuality_) && minQuality_ > 0)
340  {
342  << " : ignoring triangles with quality < "
343  << minQuality_ << " for normals calculation." << endl;
344  }
345 }
346 
347 
348 Foam::triSurfaceMesh::triSurfaceMesh(const IOobject& io, const bool isGlobal)
349 :
350  // Find instance for triSurfaceMesh
351  searchableSurface(io),
352  // Reused found instance in objectRegistry
354  (
355  IOobject
356  (
357  io.name(),
359  io.local(),
360  io.db(),
361  io.readOpt(),
362  io.writeOpt(),
363  false // searchableSurface already registered under name
364  )
365  ),
366  triSurface
367  (
368  checkFile(static_cast<const searchableSurface&>(*this), isGlobal)
369  ),
370  triSurfaceRegionSearch(static_cast<const triSurface&>(*this)),
371  minQuality_(-1),
372  surfaceClosed_(-1)
373 {
374  const pointField& pts = triSurface::points();
375 
376  bounds() = boundBox(pts);
377 }
378 
379 
381 (
382  const IOobject& io,
383  const dictionary& dict,
384  const bool isGlobal
385 )
386 :
387  searchableSurface(io),
388  // Reused found instance in objectRegistry
390  (
391  IOobject
392  (
393  io.name(),
395  io.local(),
396  io.db(),
397  io.readOpt(),
398  io.writeOpt(),
399  false // searchableSurface already registered under name
400  )
401  ),
402  triSurface
403  (
404  checkFile(static_cast<const searchableSurface&>(*this), dict, isGlobal)
405  ),
406  triSurfaceRegionSearch(static_cast<const triSurface&>(*this), dict),
407  minQuality_(-1),
408  surfaceClosed_(-1)
409 {
410  // Reading from supplied file name instead of objectPath/filePath
411  if (dict.readIfPresent("file", fName_, false, false))
412  {
413  fName_ = relativeFilePath
414  (
415  static_cast<const searchableSurface&>(*this),
416  fName_,
417  isGlobal
418  );
419  }
420 
421  scalar scaleFactor = 0;
422 
423  // Allow rescaling of the surface points
424  // eg, CAD geometries are often done in millimeters
425  if (dict.readIfPresent("scale", scaleFactor) && scaleFactor > 0)
426  {
427  Info<< searchableSurface::name() << " : using scale " << scaleFactor
428  << endl;
429  triSurface::scalePoints(scaleFactor);
430  }
431 
432  const pointField& pts = triSurface::points();
433 
434  bounds() = boundBox(pts);
435 
436  // Have optional minimum quality for normal calculation
437  if (dict.readIfPresent("minQuality", minQuality_) && minQuality_ > 0)
438  {
440  << " : ignoring triangles with quality < "
441  << minQuality_ << " for normals calculation." << endl;
442  }
443 }
444 
445 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
446 
448 {
449  clearOut();
450 }
451 
452 
454 {
456  edgeTree_.clear();
458 }
459 
460 
461 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
462 
464 {
465  tmp<pointField> tPts(new pointField(8));
466  pointField& pt = tPts.ref();
467 
468  // Use copy to calculate face centres so they don't get stored
470  (
473  ).faceCentres();
474 
475  return tPts;
476 }
477 
478 
480 (
481  pointField& centres,
482  scalarField& radiusSqr
483 ) const
484 {
485  centres = coordinates();
486  radiusSqr.setSize(size());
487  radiusSqr = 0.0;
488 
489  const pointField& pts = triSurface::points();
490 
491  forAll(*this, facei)
492  {
493  const labelledTri& f = triSurface::operator[](facei);
494  const point& fc = centres[facei];
495  forAll(f, fp)
496  {
497  const point& pt = pts[f[fp]];
498  radiusSqr[facei] = max(radiusSqr[facei], Foam::magSqr(fc-pt));
499  }
500  }
501 
502  // Add a bit to make sure all points are tested inside
503  radiusSqr += Foam::sqr(small);
504 }
505 
506 
508 {
509  return triSurface::points();
510 }
511 
512 
514 {
515  const indexedOctree<treeDataTriSurface>& octree = tree();
516 
517  labelList indices = octree.findBox(treeBoundBox(bb));
518 
519  return !indices.empty();
520 }
521 
522 
524 {
526  edgeTree_.clear();
527  triSurface::movePoints(newPoints);
528 }
529 
530 
533 {
534  if (edgeTree_.empty())
535  {
536  // Boundary edges
537  labelList bEdges
538  (
539  identity
540  (
541  nEdges()
542  -nInternalEdges()
543  )
544  + nInternalEdges()
545  );
546 
547  treeBoundBox bb(Zero, Zero);
548 
549  if (bEdges.size())
550  {
551  label nPoints;
553  (
554  *this,
555  bb,
556  nPoints
557  );
558 
559  // Slightly extended bb. Slightly off-centred just so on symmetric
560  // geometry there are less face/edge aligned items.
561  bb = bb.extend(1e-4);
562  }
563 
564  scalar oldTol = indexedOctree<treeDataEdge>::perturbTol();
566 
567  edgeTree_.reset
568  (
570  (
572  (
573  false, // cachebb
574  edges(), // edges
575  localPoints(), // points
576  bEdges // selected edges
577  ),
578  bb, // bb
579  maxTreeDepth(), // maxLevel
580  10, // leafsize
581  3.0 // duplicity
582  )
583  );
584 
586  }
587  return edgeTree_();
588 }
589 
590 
592 {
593  if (regions_.empty())
594  {
595  regions_.setSize(patches().size());
596  forAll(regions_, regionI)
597  {
598  regions_[regionI] = patches()[regionI].name();
599  }
600  }
601  return regions_;
602 }
603 
604 
606 {
607  if (surfaceClosed_ == -1)
608  {
609  if (isSurfaceClosed())
610  {
611  surfaceClosed_ = 1;
612  }
613  else
614  {
615  surfaceClosed_ = 0;
616  }
617  }
618 
619  return surfaceClosed_ == 1;
620 }
621 
622 
624 (
625  const pointField& samples,
626  const scalarField& nearestDistSqr,
627  List<pointIndexHit>& info
628 ) const
629 {
630  triSurfaceSearch::findNearest(samples, nearestDistSqr, info);
631 }
632 
633 
635 (
636  const pointField& samples,
637  const scalarField& nearestDistSqr,
638  const labelList& regionIndices,
639  List<pointIndexHit>& info
640 ) const
641 {
643  (
644  samples,
645  nearestDistSqr,
646  regionIndices,
647  info
648  );
649 }
650 
651 
653 (
654  const pointField& start,
655  const pointField& end,
656  List<pointIndexHit>& info
657 ) const
658 {
659  triSurfaceSearch::findLine(start, end, info);
660 }
661 
662 
664 (
665  const pointField& start,
666  const pointField& end,
667  List<pointIndexHit>& info
668 ) const
669 {
670  triSurfaceSearch::findLineAny(start, end, info);
671 }
672 
673 
675 (
676  const pointField& start,
677  const pointField& end,
679 ) const
680 {
681  triSurfaceSearch::findLineAll(start, end, info);
682 }
683 
684 
686 (
687  const List<pointIndexHit>& info,
688  labelList& region
689 ) const
690 {
691  region.setSize(info.size());
692  forAll(info, i)
693  {
694  if (info[i].hit())
695  {
696  region[i] = triSurface::operator[](info[i].index()).region();
697  }
698  else
699  {
700  region[i] = -1;
701  }
702  }
703 }
704 
705 
707 (
708  const List<pointIndexHit>& info,
709  vectorField& normal
710 ) const
711 {
712  const triSurface& s = *this;
713  const pointField& pts = s.points();
714 
715  normal.setSize(info.size());
716 
717  if (minQuality_ >= 0)
718  {
719  // Make sure we don't use triangles with low quality since
720  // normal is not reliable.
721 
722  const labelListList& faceFaces = s.faceFaces();
723 
724  forAll(info, i)
725  {
726  if (info[i].hit())
727  {
728  label facei = info[i].index();
729  normal[i] = s[facei].area(pts);
730 
731  scalar qual = s[facei].tri(pts).quality();
732 
733  if (qual < minQuality_)
734  {
735  // Search neighbouring triangles
736  const labelList& fFaces = faceFaces[facei];
737 
738  forAll(fFaces, j)
739  {
740  label nbrI = fFaces[j];
741  scalar nbrQual = s[nbrI].tri(pts).quality();
742  if (nbrQual > qual)
743  {
744  qual = nbrQual;
745  normal[i] = s[nbrI].area(pts);
746  }
747  }
748  }
749 
750  normal[i] /= mag(normal[i]) + vSmall;
751  }
752  else
753  {
754  // Set to what?
755  normal[i] = Zero;
756  }
757  }
758  }
759  else
760  {
761  forAll(info, i)
762  {
763  if (info[i].hit())
764  {
765  label facei = info[i].index();
766  // Cached:
767  // normal[i] = faceNormals()[facei];
768 
769  // Uncached
770  normal[i] = s[facei].normal(pts);
771  }
772  else
773  {
774  // Set to what?
775  normal[i] = Zero;
776  }
777  }
778  }
779 }
780 
781 
783 (
784  const pointField& points,
785  List<volumeType>& volType
786 ) const
787 {
788  volType.setSize(points.size());
789 
792 
793  forAll(points, pointi)
794  {
795  const point& pt = points[pointi];
796 
797  if (!tree().bb().contains(pt))
798  {
799  // Have to calculate directly as outside the octree
800  volType[pointi] = tree().shapes().getVolumeType(tree(), pt);
801  }
802  else
803  {
804  // - use cached volume type per each tree node
805  volType[pointi] = tree().getVolumeType(pt);
806  }
807  }
808 
810 }
811 
812 
814 {
816  (
818  (
819  IOobject
820  (
821  "values",
824  *this,
827  ),
828  *this,
829  dimless,
830  labelField(values)
831  )
832  );
833 
834  // Store field on triMesh
835  fldPtr.ptr()->store();
836 }
837 
838 
840 (
841  const List<pointIndexHit>& info,
842  labelList& values
843 ) const
844 {
845  if (foundObject<triSurfaceLabelField>("values"))
846  {
847  values.setSize(info.size());
848 
849  const triSurfaceLabelField& fld = lookupObject<triSurfaceLabelField>
850  (
851  "values"
852  );
853 
854  forAll(info, i)
855  {
856  if (info[i].hit())
857  {
858  values[i] = fld[info[i].index()];
859  }
860  }
861  }
862 }
863 
864 
866 (
870  const bool write
871 ) const
872 {
873  fileName fullPath;
874  if (fName_.size())
875  {
876  // Override file name
877 
878  fullPath = fName_;
879 
880  fullPath.expand();
881  if (!fullPath.isAbsolute())
882  {
883  // Add directory from regIOobject
884  fullPath = searchableSurface::objectPath().path()/fullPath;
885  }
886  }
887  else
888  {
889  fullPath = searchableSurface::objectPath();
890  }
891 
892  if (!mkDir(fullPath.path()))
893  {
894  return false;
895  }
896 
897  triSurface::write(fullPath);
898 
899  if (!isFile(fullPath))
900  {
901  return false;
902  }
903 
904  return true;
905 }
906 
907 
908 // ************************************************************************* //
virtual void findLine(const pointField &start, const pointField &end, List< pointIndexHit > &) const
Find first intersection on segment from start to end.
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
label nPoints() const
Return number of points supporting patch faces.
label maxTreeDepth() const
Return max tree depth of octree.
Field< label > labelField
Specialisation of Field<T> for label.
Definition: labelField.H:49
bool exists(const fileName &, const bool checkVariants=true, const bool followLink=true)
Does the name exist (as directory or file) in the file system?
Definition: POSIX.C:520
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
virtual tmp< pointField > coordinates() const
Get representative set of element coordinates.
bool empty() const
Return true if the UList is empty (ie, size() is zero)
Definition: UListI.H:325
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
const word & name() const
Return name.
Definition: IOobject.H:315
A class for handling file names.
Definition: fileName.H:79
triSurfaceMesh(const IOobject &, const triSurface &)
Construct from triSurface.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
void clearOut()
Clear storage.
void findNearest(const pointField &samples, const scalarField &nearestDistSqr, const labelList &regionIndices, List< pointIndexHit > &info) const
Find the nearest point on the surface out of the regions.
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
virtual void boundingSpheres(pointField &centres, scalarField &radiusSqr) const
Get bounding spheres (centre and radius squared). Any point.
T & operator[](const label)
Return element of UList.
Definition: UListI.H:167
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
const Field< PointType > & faceCentres() const
Return face centres for patch.
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
triSurfaceRegionSearch(const triSurface &)
Construct from surface. Holds reference to surface!
Fields for triSurface.
virtual tmp< pointField > points() const
Get the points that define the surface.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
virtual ~triSurfaceMesh()
Destructor.
bool isFile(const fileName &, const bool checkVariants=true, const bool followLink=true)
Does the name exist as a file in the file system?
Definition: POSIX.C:555
dimensionedSymmTensor sqr(const dimensionedVector &dv)
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
virtual void findNearest(const pointField &sample, const scalarField &nearestDistSqr, List< pointIndexHit > &) const
virtual void setField(const labelList &values)
WIP. Store element-wise field.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
static scalar & perturbTol()
Get the perturbation tolerance.
fileName objectPath() const
Return complete path + object name.
Definition: regIOobject.H:159
Helper class to search on triSurface. Creates an octree for each region of the surface and only searc...
virtual bool overlaps(const boundBox &bb) const
Does any part of the surface overlap the supplied bound box?
Holds data for octree to work on an edges subset.
Definition: treeDataEdge.H:53
label nInternalEdges() const
Number of internal edges.
labelList identity(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
A bounding box defined in terms of the points at its extremities.
Definition: boundBox.H:58
Base class of (analytical or triangulated) surface. Encapsulates all the search routines. WIP.
const dimensionSet dimless
virtual void scalePoints(const scalar)
Scale points. A non-positive factor is ignored.
Definition: triSurface.C:749
T * ptr()
Return object pointer for reuse.
Definition: autoPtrI.H:90
const Field< PointType > & localPoints() const
Return pointField of points in patch.
Macros for easy insertion into run-time selection tables.
void findNearest(const pointField &samples, const scalarField &nearestDistSqr, List< pointIndexHit > &info) const
virtual bool hasVolumeType() const
Whether supports volume type below. I.e. whether is closed.
virtual void movePoints(const pointField &)
Move points.
Definition: triSurface.C:736
virtual void getVolumeType(const pointField &, List< volumeType > &) const
Determine type (inside/outside/mixed) for point. unknown if.
objectRegistry(const Time &db, const label nIoObjects=128)
Construct the time objectRegistry given an initial estimate.
A list of faces which address into the list of points.
A List obtained as a section of another List.
Definition: SubList.H:53
const labelListList & faceFaces() const
Return face-face addressing.
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))
virtual bool writeObject(IOstream::streamFormat fmt, IOstream::versionNumber ver, IOstream::compressionType cmp, const bool write=true) const
Write using given format, version and compression.
static void calcBounds(const PrimitivePatch< FaceList, PointField > &p, boundBox &bb, label &nPoints)
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
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){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
void findLine(const pointField &start, const pointField &end, List< pointIndexHit > &info) const
void findLineAll(const pointField &start, const pointField &end, List< List< pointIndexHit >> &info) const
Calculate all intersections from start to end.
bool isAbsolute() const
Return true if file name is absolute.
Definition: fileName.C:73
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:125
static const word & geometryDir()
Return the geometry directory name.
virtual fileName filePath(const bool globalFile, const IOobject &, const word &typeName) const =0
Search for an object. globalFile : also check undecomposed case.
void findLineAny(const pointField &start, const pointField &end, List< pointIndexHit > &info) const
const fileName & local() const
Definition: IOobject.H:409
virtual void findLineAll(const pointField &start, const pointField &end, List< List< pointIndexHit >> &) const
Get all intersections in order from start to end.
bool readIfPresent(const word &, T &, bool recursive=false, bool patternMatch=true) const
Find an entry if present, and assign to T.
streamFormat
Enumeration for the format of data in the stream.
Definition: IOstream.H:86
const Field< PointType > & points() const
Return reference to global points.
static const word null
An empty word.
Definition: word.H:77
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
word timeName
Definition: getTimeIndex.H:3
List< label > labelList
A List of labels.
Definition: labelList.H:56
static const zero Zero
Definition: zero.H:97
Triangle with additional region number.
Definition: labelledTri.H:57
const fileOperation & fileHandler()
Get current file handler.
searchableSurface(const IOobject &io)
compressionType
Enumeration for the format of data in the stream.
Definition: IOstream.H:193
dimensioned< scalar > magSqr(const dimensioned< Type > &)
label nEdges() const
Return number of edges in patch.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
const indexedOctree< treeDataEdge > & edgeTree() const
Demand driven construction of octree for boundary edges.
const Time & time() const
Return time.
defineTypeNameAndDebug(combustionModel, 0)
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
bool mkDir(const fileName &, mode_t=0777)
Make a directory and return an error if it could not be created.
Definition: POSIX.C:290
virtual void getRegion(const List< pointIndexHit > &, labelList &region) const
From a set of points and indices get the region.
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
virtual void getNormal(const List< pointIndexHit > &, vectorField &normal) const
From a set of points and indices get the normal.
void setSize(const label)
Reset size of List.
Definition: List.C:281
const geometricSurfacePatchList & patches() const
Definition: triSurface.H:322
writeOption writeOpt() const
Definition: IOobject.H:375
triSurface()
Construct null.
Definition: triSurface.C:534
std::remove_reference< ::Foam::List< labelledTri > >::type::value_type FaceType
Non-pointer based hierarchical recursive searching.
Definition: treeDataEdge.H:47
string & expand(const bool allowEmpty=false)
Expand initial tildes and all occurrences of environment variables.
Definition: string.C:125
virtual label size() const
Range of local indices that can be returned.
virtual void findLineAny(const pointField &start, const pointField &end, List< pointIndexHit > &) const
Return any intersection on segment from start to end.
virtual void movePoints(const pointField &)
Move points.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
scalar tolerance() const
Return tolerance to use in searches.
void invertManyToMany(const label len, const UList< InList > &, List< OutList > &)
Invert many-to-many.
Version number type.
Definition: IOstream.H:96
Standard boundBox + extra functionality for use in octree.
Definition: treeBoundBox.H:87
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
fileName path() const
Return directory path name (part before last /)
Definition: fileName.C:265
const boundBox & bounds() const
Return const reference to boundBox.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
fileName & instance() const
Return the instance directory, constant, system, <time> etc.
Definition: IOobject.C:355
A class for managing temporary objects.
Definition: PtrList.H:53
Registry of regIOobjects.
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:312
friend class iterator
Declare friendship with the iterator.
Definition: HashTable.H:194
Triangulated surface description with patch information.
Definition: triSurface.H:66
readOption readOpt() const
Definition: IOobject.H:365
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:98
label size() const
Return the number of elements in the UList.
Definition: ListI.H:171
const indexedOctree< treeDataTriSurface > & tree() const
Demand driven construction of the octree.
virtual const wordList & regions() const
Names of regions.
treeBoundBox extend(const scalar s) const
Return asymetrically extended bounding box, with guaranteed.
Namespace for OpenFOAM.
virtual void getField(const List< pointIndexHit > &, labelList &) const
WIP. From a set of hits (points and.