triSurfaceMesh.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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2015-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "triSurfaceMesh.H"
30 #include "Random.H"
32 #include "edgeHashes.H"
33 #include "triSurfaceFields.H"
34 #include "Time.H"
35 #include "PatchTools.H"
36 
37 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
41  defineTypeNameAndDebug(triSurfaceMesh, 0);
42  addToRunTimeSelectionTable(searchableSurface, triSurfaceMesh, dict);
43 }
44 
46 
47 
48 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
49 
51 (
52  const edge& e,
53  EdgeMap<label>& facesPerEdge
54 )
55 {
56  label& count = facesPerEdge(e, 0); // lookup or new entry
57  if (count == 2)
58  {
59  return false;
60  }
61 
62  ++count;
63  return true;
64 }
65 
66 
68 {
69  if (debug)
70  {
71  Pout<< "triSurfaceMesh::isSurfaceClosed:"
72  << " determining closedness for surface with "
73  << triSurface::size() << " triangles" << endl;
74  }
75 
76  const pointField& pts = triSurface::points();
77 
78  // Construct pointFaces. Let's hope surface has compact point
79  // numbering ...
81  invertManyToMany(pts.size(), *this, pointFaces);
82 
83  // Loop over all faces surrounding point. Count edges emanating from point.
84  // Every edge should be used by two faces exactly.
85  // To prevent doing work twice per edge only look at edges to higher
86  // point
87  EdgeMap<label> facesPerEdge(128);
88  forAll(pointFaces, pointi)
89  {
90  const labelList& pFaces = pointFaces[pointi];
91 
92  facesPerEdge.clear();
93  for (const label facei : pFaces)
94  {
95  const triSurface::face_type& f = triSurface::operator[](facei);
96  const label fp = f.find(pointi);
97 
98  // Something weird: if I expand the code of addFaceToEdge in both
99  // below instances it gives a segmentation violation on some
100  // surfaces. Compiler (4.3.2) problem?
101 
102 
103  // Forward edge
104  const label nextPointi = f[f.fcIndex(fp)];
105 
106  if (nextPointi > pointi)
107  {
108  bool okFace = addFaceToEdge
109  (
110  edge(pointi, nextPointi),
111  facesPerEdge
112  );
113 
114  if (!okFace)
115  {
116  if (debug)
117  {
118  Pout<< "triSurfaceMesh::isSurfaceClosed :"
119  << " surface is open" << endl;
120  }
121  return false;
122  }
123  }
124 
125  // Reverse edge
126  const label prevPointi = f[f.rcIndex(fp)];
127 
128  if (prevPointi > pointi)
129  {
130  bool okFace = addFaceToEdge
131  (
132  edge(pointi, prevPointi),
133  facesPerEdge
134  );
135 
136  if (!okFace)
137  {
138  if (debug)
139  {
140  Pout<< "triSurfaceMesh::isSurfaceClosed :"
141  << " surface is open" << endl;
142  }
143  return false;
144  }
145  }
146  }
147 
148  // Check for any edges used only once.
149  forAllConstIters(facesPerEdge, iter)
150  {
151  if (iter.val() != 2)
152  {
153  if (debug)
154  {
155  Pout<< "triSurfaceMesh::isSurfaceClosed :"
156  << " surface is open" << endl;
157  }
158  return false;
159  }
160  }
161  }
162 
163  if (debug)
164  {
165  Pout<< "triSurfaceMesh::isSurfaceClosed :"
166  << " surface is closed" << endl;
167  }
168  return true;
169 }
170 
171 
172 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
173 
175 :
176  searchableSurface(io),
178  (
179  IOobject
180  (
181  io.name(),
182  io.instance(),
183  io.local(),
184  io.db(),
185  io.readOpt(),
186  io.writeOpt(),
187  false // searchableSurface already registered under name
188  )
189  ),
190  triSurface(s),
191  triSurfaceRegionSearch(static_cast<const triSurface&>(*this)),
192  minQuality_(-1),
193  surfaceClosed_(-1),
194  outsideVolType_(volumeType::UNKNOWN)
195 {
196  const pointField& pts = triSurface::points();
197 
198  bounds() = boundBox(pts, false);
199 }
200 
201 
203 :
204  // Find instance for triSurfaceMesh
205  searchableSurface(io),
206  // Reused found instance in objectRegistry
208  (
209  IOobject
210  (
211  io.name(),
212  searchableSurface::instance(),
213  io.local(),
214  io.db(),
215  io.readOpt(),
216  io.writeOpt(),
217  false // searchableSurface already registered under name
218  )
219  ),
220  triSurface(static_cast<const searchableSurface&>(*this), dictionary::null),
221  triSurfaceRegionSearch(static_cast<const triSurface&>(*this)),
222  minQuality_(-1),
223  surfaceClosed_(-1),
224  outsideVolType_(volumeType::UNKNOWN)
225 {
226  const pointField& pts = triSurface::points();
227 
228  bounds() = boundBox(pts, false);
229 }
230 
231 
233 (
234  const IOobject& io,
235  const dictionary& dict
236 )
237 :
238  searchableSurface(io),
239  // Reused found instance in objectRegistry
241  (
242  IOobject
243  (
244  io.name(),
246  io.local(),
247  io.db(),
248  io.readOpt(),
249  io.writeOpt(),
250  false // searchableSurface already registered under name
251  )
252  ),
253  triSurface(static_cast<const searchableSurface&>(*this), dict),
254  triSurfaceRegionSearch(static_cast<const triSurface&>(*this), dict),
255  minQuality_(-1),
256  surfaceClosed_(-1),
257  outsideVolType_(volumeType::UNKNOWN)
258 {
259  // Read with supplied file name instead of objectPath/filePath
260  if (dict.readIfPresent("file", fName_, keyType::LITERAL))
261  {
263  (
264  static_cast<const searchableSurface&>(*this),
265  fName_,
266  true
267  );
268  }
269 
270  // Report optional scale factor (eg, mm to m),
271  // which was already applied during triSurface construction
272  scalar scaleFactor{0};
273  if (dict.getOrDefault("scale", scaleFactor) && scaleFactor > 0)
274  {
276  << " : using scale " << scaleFactor << endl;
277  }
278 
279  const pointField& pts = triSurface::points();
280 
281  bounds() = boundBox(pts, false);
282 
283  // Have optional minimum quality for normal calculation
284  if (dict.readIfPresent("minQuality", minQuality_) && minQuality_ > 0)
285  {
287  << " : ignoring triangles with quality < "
288  << minQuality_ << " for normals calculation." << endl;
289  }
290 }
291 
292 
294 :
295  // Find instance for triSurfaceMesh
296  searchableSurface(io),
297  // Reused found instance in objectRegistry
299  (
300  IOobject
301  (
302  io.name(),
303  searchableSurface::instance(),
304  io.local(),
305  io.db(),
306  io.readOpt(),
307  io.writeOpt(),
308  false // searchableSurface already registered under name
309  )
310  ),
311  triSurface(),
312  triSurfaceRegionSearch(static_cast<const triSurface&>(*this)),
313  minQuality_(-1),
314  surfaceClosed_(-1),
315  outsideVolType_(volumeType::UNKNOWN)
316 {
317  // Check IO flags
318  if (io.readOpt() != IOobject::NO_READ)
319  {
320  const bool searchGlobal(r == localOrGlobal || r == masterOnly);
321 
322  const fileName actualFile
323  (
324  searchGlobal
325  ? io.globalFilePath(typeName)
326  : io.localFilePath(typeName)
327  );
328 
329  if (debug)
330  {
331  Pout<< "triSurfaceMesh(const IOobject& io) :"
332  << " loading surface " << io.objectPath()
333  << " local filePath:" << io.localFilePath(typeName)
334  << " from:" << actualFile << endl;
335  }
336 
337  if (searchGlobal && Pstream::parRun())
338  {
339  // Check where surface was found
340  const fileName localFile(io.localFilePath(typeName));
341 
342  if (r == masterOnly && (actualFile != localFile))
343  {
344  // Found undecomposed surface. Load on master only
345  if (Pstream::master())
346  {
347  triSurface s2(actualFile);
349  }
351  if (debug)
352  {
353  Pout<< "triSurfaceMesh(const IOobject& io) :"
354  << " loaded triangles:" << triSurface::size() << endl;
355  }
356  }
357  else
358  {
359  // Read on all processors
360  triSurface s2(actualFile);
362  if (debug)
363  {
364  Pout<< "triSurfaceMesh(const IOobject& io) :"
365  << " loaded triangles:" << triSurface::size() << endl;
366  }
367  }
368  }
369  else
370  {
371  // Read on all processors
372  triSurface s2(actualFile);
374  if (debug)
375  {
376  Pout<< "triSurfaceMesh(const IOobject& io) :"
377  << " loaded triangles:" << triSurface::size() << endl;
378  }
379  }
380  }
381 
382  const pointField& pts = triSurface::points();
383  bounds() = boundBox(pts, false);
384 }
385 
386 
388 (
389  const IOobject& io,
390  const dictionary& dict,
391  const readAction r
392 )
393 :
394  searchableSurface(io),
395  // Reused found instance in objectRegistry
397  (
398  IOobject
399  (
400  io.name(),
402  io.local(),
403  io.db(),
404  io.readOpt(),
405  io.writeOpt(),
406  false // searchableSurface already registered under name
407  )
408  ),
409  triSurface(),
410  triSurfaceRegionSearch(static_cast<const triSurface&>(*this), dict),
411  minQuality_(-1),
412  surfaceClosed_(-1),
413  outsideVolType_(volumeType::UNKNOWN)
414 {
415  if (io.readOpt() != IOobject::NO_READ)
416  {
417  // Surface type (optional)
418  const word surfType(dict.getOrDefault<word>("fileType", word::null));
419 
420  // Scale factor (optional)
421  const scalar scaleFactor(dict.getOrDefault<scalar>("scale", 0));
422 
423  const bool searchGlobal(r == localOrGlobal || r == masterOnly);
424 
425  fileName actualFile
426  (
427  searchGlobal
428  ? io.globalFilePath(typeName)
429  : io.localFilePath(typeName)
430  );
431 
432  // Reading from supplied file name instead of objectPath/filePath
433  if (dict.readIfPresent("file", fName_, keyType::LITERAL))
434  {
435  fName_ = relativeFilePath
436  (
437  static_cast<const searchableSurface&>(*this),
438  fName_,
439  searchGlobal
440  );
441  actualFile = fName_;
442  }
443 
444  if (debug)
445  {
446  Pout<< "triSurfaceMesh(const IOobject& io, const dictionary&) :"
447  << " loading surface " << io.objectPath()
448  << " local filePath:" << io.localFilePath(typeName)
449  << " from:" << actualFile << endl;
450  }
451 
452 
453  if (searchGlobal && Pstream::parRun())
454  {
455  // Check where surface was found
456  const fileName localFile(io.localFilePath(typeName));
457 
458  if (r == masterOnly && (actualFile != localFile))
459  {
460  // Surface not loaded from processor directories -> undecomposed
461  // surface. Load on master only
462  if (Pstream::master())
463  {
464  triSurface s2(actualFile, surfType, scaleFactor);
466  }
468  if (debug)
469  {
470  Pout<< "triSurfaceMesh(const IOobject& io) :"
471  << " loaded triangles:" << triSurface::size() << endl;
472  }
473  }
474  else
475  {
476  // Read on all processors
477  triSurface s2(actualFile, surfType, scaleFactor);
479  if (debug)
480  {
481  Pout<< "triSurfaceMesh(const IOobject& io) :"
482  << " loaded triangles:" << triSurface::size() << endl;
483  }
484  }
485  }
486  else
487  {
488  // Read on all processors
489  triSurface s2(actualFile, surfType, scaleFactor);
491  if (debug)
492  {
493  Pout<< "triSurfaceMesh(const IOobject& io) :"
494  << " loaded triangles:" << triSurface::size() << endl;
495  }
496  }
497 
498 
499  // Report optional scale factor (eg, mm to m),
500  // which was already applied during triSurface reading
501  if (scaleFactor > 0)
502  {
504  << " : using scale " << scaleFactor << endl;
505  }
506  }
507 
508 
509  const pointField& pts = triSurface::points();
510  bounds() = boundBox(pts, false);
511 
512  // Have optional minimum quality for normal calculation
513  if (dict.readIfPresent("minQuality", minQuality_) && minQuality_ > 0)
514  {
516  << " : ignoring triangles with quality < "
517  << minQuality_ << " for normals calculation." << endl;
518  }
519 }
520 
521 
522 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
523 
525 {
526  clearOut();
527 }
528 
529 
531 {
532  // Do not clear closedness status
534  edgeTree_.clear();
536 }
537 
538 
539 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
540 
542 {
543  auto tpts = tmp<pointField>::New();
544  auto& pts = tpts.ref();
545 
547  {
548  // Can reuse existing values
549  pts = triSurface::faceCentres();
550  }
551  else
552  {
554 
555  // Calculate face centres from a copy to avoid incurring
556  // additional storage
558  (
559  FaceListType(*this, triSurface::size()),
561  ).faceCentres();
562  }
563 
564  return tpts;
565 }
566 
567 
569 (
570  pointField& centres,
571  scalarField& radiusSqr
572 ) const
573 {
574  centres = coordinates();
575  radiusSqr.setSize(size());
576  radiusSqr = 0.0;
577 
578  const pointField& pts = triSurface::points();
579 
580  forAll(*this, facei)
581  {
582  const labelledTri& f = triSurface::operator[](facei);
583  const point& fc = centres[facei];
584  for (const label pointi : f)
585  {
586  const point& pt = pts[pointi];
587  radiusSqr[facei] = max(radiusSqr[facei], Foam::magSqr(fc-pt));
588  }
589  }
590 
591  // Add a bit to make sure all points are tested inside
592  radiusSqr += Foam::sqr(SMALL);
593 }
594 
595 
597 {
598  return triSurface::points();
599 }
600 
601 
603 {
604  const indexedOctree<treeDataTriSurface>& octree = tree();
605 
606  labelList indices = octree.findBox(treeBoundBox(bb));
607 
608  return !indices.empty();
609 }
610 
611 
613 {
614  if (debug)
615  {
616  Pout<< "triSurfaceMesh::movePoints :"
617  << " moving at time " << objectRegistry::time().timeName()
618  << endl;
619  }
620 
621  // Preserve topological point status (surfaceClosed_, outsideVolType_)
622 
623  // Update local information (instance, event number)
626 
627  const label event = getEvent();
628  searchableSurface::eventNo() = event;
630 
631  // Clear additional addressing
633  edgeTree_.clear();
634  triSurface::movePoints(newPoints);
635 
636  bounds() = boundBox(triSurface::points(), false);
637  if (debug)
638  {
639  Pout<< "triSurfaceMesh::movePoints: finished moving points" << endl;
640  }
641 }
642 
643 
646 {
647  if (!edgeTree_)
648  {
649  if (debug)
650  {
651  Pout<< "triSurfaceMesh::edgeTree :"
652  << " constructing tree for " << nEdges() - nInternalEdges()
653  << " boundary edges" << endl;
654  }
655 
656  // Boundary edges
657  labelList bEdges
658  (
659  identity(nEdges() - nInternalEdges(), nInternalEdges())
660  );
661 
662  treeBoundBox bb(Zero, Zero);
663 
664  if (bEdges.size())
665  {
666  label nPoints;
668  (
669  *this,
670  bb,
671  nPoints
672  );
673 
674  // Random number generator. Bit dodgy since not exactly random ;-)
675  Random rndGen(65431);
676 
677  // Slightly extended bb. Slightly off-centred just so on symmetric
678  // geometry there are less face/edge aligned items.
679 
680  bb = bb.extend(rndGen, 1e-4);
681  bb.min() -= point::uniform(ROOTVSMALL);
682  bb.max() += point::uniform(ROOTVSMALL);
683  }
684 
685 
686  if (debug)
687  {
688  Pout<< "triSurfaceMesh::edgeTree : "
689  << "calculating edge tree for bb:" << bb << endl;
690  }
691 
692  scalar oldTol = indexedOctree<treeDataEdge>::perturbTol();
694 
695  edgeTree_.reset
696  (
698  (
700  (
701  false, // cachebb
702  edges(), // edges
703  localPoints(), // points
704  bEdges // selected edges
705  ),
706  bb, // bb
707  maxTreeDepth(), // maxLevel
708  10, // leafsize
709  3.0 // duplicity
710  )
711  );
712 
714 
715  if (debug)
716  {
717  Pout<< "triSurfaceMesh::edgeTree :"
718  << " finished constructing tree for "
719  << nEdges() - nInternalEdges()
720  << " boundary edges" << endl;
721  }
722  }
723 
724  return *edgeTree_;
725 }
726 
727 
729 {
730  if (regions_.empty())
731  {
732  regions_.setSize(patches().size());
733  forAll(regions_, regioni)
734  {
735  regions_[regioni] = patches()[regioni].name();
736  }
737  }
738  return regions_;
739 }
740 
741 
743 {
744  if (surfaceClosed_ == -1)
745  {
746  if (isSurfaceClosed())
747  {
748  surfaceClosed_ = 1;
749  }
750  else
751  {
752  surfaceClosed_ = 0;
753  }
754  }
755 
756  return surfaceClosed_ == 1;
757 }
758 
759 
761 {
762  if (outsideVolType_ == volumeType::UNKNOWN)
763  {
764  // Get point outside bounds()
765  const point outsidePt(bounds().max() + 0.5*bounds().span());
766 
767  if (debug)
768  {
769  Pout<< "triSurfaceMesh::outsideVolumeType :"
770  << " triggering outsidePoint:" << outsidePt
771  << " orientation" << endl;
772  }
773 
774  //outsideVolType_ = tree().shapes().getVolumeType(tree(), outsidePt);
775  // Note: do not use tree directly so e.g. distributedTriSurfaceMesh
776  // has opportunity to intercept
777  List<volumeType> outsideVolTypes;
778  getVolumeType(pointField(1, outsidePt), outsideVolTypes);
779  outsideVolType_ = outsideVolTypes[0];
780 
781  if (debug)
782  {
783  Pout<< "triSurfaceMesh::outsideVolumeType :"
784  << " finished outsidePoint:" << outsidePt
785  << " orientation:" << volumeType::names[outsideVolType_]
786  << endl;
787  }
788  }
789 
790  return outsideVolType_;
791 }
792 
793 
795 (
796  const pointField& samples,
797  const scalarField& nearestDistSqr,
798  List<pointIndexHit>& info
799 ) const
800 {
801  if (debug)
802  {
803  Pout<< "triSurfaceMesh::findNearest :"
804  << " trying to find nearest for " << samples.size()
805  << " samples with max sphere "
806  << (samples.size() ? Foam::sqrt(max(nearestDistSqr)) : Zero)
807  << endl;
808  }
809  triSurfaceSearch::findNearest(samples, nearestDistSqr, info);
810  if (debug)
811  {
812  Pout<< "triSurfaceMesh::findNearest :"
813  << " finished trying to find nearest for " << samples.size()
814  << " samples" << endl;
815  }
816 }
817 
818 
820 (
821  const pointField& samples,
822  const scalarField& nearestDistSqr,
823  const labelList& regionIndices,
824  List<pointIndexHit>& info
825 ) const
826 {
827  if (debug)
828  {
829  Pout<< "triSurfaceMesh::findNearest :"
830  << " trying to find nearest and region for " << samples.size()
831  << " samples with max sphere "
832  << (samples.size() ? Foam::sqrt(max(nearestDistSqr)) : Zero)
833  << endl;
834  }
836  (
837  samples,
838  nearestDistSqr,
839  regionIndices,
840  info
841  );
842  if (debug)
843  {
844  Pout<< "triSurfaceMesh::findNearest :"
845  << " finished trying to find nearest and region for "
846  << samples.size() << " samples" << endl;
847  }
848 }
849 
850 
852 (
853  const pointField& start,
854  const pointField& end,
855  List<pointIndexHit>& info
856 ) const
857 {
858  if (debug)
859  {
860  Pout<< "triSurfaceMesh::findLine :"
861  << " intersecting with "
862  << start.size() << " rays" << endl;
863  }
864  triSurfaceSearch::findLine(start, end, info);
865  if (debug)
866  {
867  Pout<< "triSurfaceMesh::findLine :"
868  << " finished intersecting with "
869  << start.size() << " rays" << endl;
870  }
871 }
872 
873 
875 (
876  const pointField& start,
877  const pointField& end,
878  List<pointIndexHit>& info
879 ) const
880 {
881  if (debug)
882  {
883  Pout<< "triSurfaceMesh::findLineAny :"
884  << " intersecting with "
885  << start.size() << " rays" << endl;
886  }
887  triSurfaceSearch::findLineAny(start, end, info);
888  if (debug)
889  {
890  Pout<< "triSurfaceMesh::findLineAny :"
891  << " finished intersecting with "
892  << start.size() << " rays" << endl;
893  }
894 }
895 
896 
898 (
899  const pointField& start,
900  const pointField& end,
902 ) const
903 {
904  if (debug)
905  {
906  Pout<< "triSurfaceMesh::findLineAll :"
907  << " intersecting with "
908  << start.size() << " rays" << endl;
909  }
910  triSurfaceSearch::findLineAll(start, end, info);
911  if (debug)
912  {
913  Pout<< "triSurfaceMesh::findLineAll :"
914  << " finished intersecting with "
915  << start.size() << " rays" << endl;
916  }
917 }
918 
919 
921 (
922  const List<pointIndexHit>& info,
923  labelList& region
924 ) const
925 {
926  if (debug)
927  {
928  Pout<< "triSurfaceMesh::getRegion :"
929  << " getting region for "
930  << info.size() << " triangles" << endl;
931  }
932  region.setSize(info.size());
933  forAll(info, i)
934  {
935  if (info[i].hit())
936  {
937  region[i] = triSurface::operator[](info[i].index()).region();
938  }
939  else
940  {
941  region[i] = -1;
942  }
943  }
944  if (debug)
945  {
946  Pout<< "triSurfaceMesh::getRegion :"
947  << " finished getting region for "
948  << info.size() << " triangles" << endl;
949  }
950 }
951 
952 
954 (
955  const List<pointIndexHit>& info,
956  vectorField& normal
957 ) const
958 {
959  if (debug)
960  {
961  Pout<< "triSurfaceMesh::getNormal :"
962  << " getting normal for "
963  << info.size() << " triangles" << endl;
964  }
965 
966  const triSurface& s = *this;
967  const pointField& pts = s.points();
968 
969  normal.setSize(info.size());
970 
971  if (minQuality_ >= 0)
972  {
973  // Make sure we don't use triangles with low quality since
974  // normal is not reliable.
975 
976  const labelListList& faceFaces = s.faceFaces();
977 
978  forAll(info, i)
979  {
980  if (info[i].hit())
981  {
982  const label facei = info[i].index();
983  normal[i] = s[facei].unitNormal(pts);
984 
985  scalar qual = s[facei].tri(pts).quality();
986 
987  if (qual < minQuality_)
988  {
989  // Search neighbouring triangles
990  const labelList& fFaces = faceFaces[facei];
991 
992  for (const label nbri : fFaces)
993  {
994  scalar nbrQual = s[nbri].tri(pts).quality();
995  if (nbrQual > qual)
996  {
997  qual = nbrQual;
998  normal[i] = s[nbri].unitNormal(pts);
999  }
1000  }
1001  }
1002  }
1003  else
1004  {
1005  // Set to what?
1006  normal[i] = Zero;
1007  }
1008  }
1009  }
1010  else
1011  {
1012  forAll(info, i)
1013  {
1014  if (info[i].hit())
1015  {
1016  const label facei = info[i].index();
1017 
1018  // Uncached
1019  normal[i] = s[facei].unitNormal(pts);
1020  }
1021  else
1022  {
1023  // Set to what?
1024  normal[i] = Zero;
1025  }
1026  }
1027  }
1028  if (debug)
1029  {
1030  Pout<< "triSurfaceMesh::getNormal :"
1031  << " finished getting normal for "
1032  << info.size() << " triangles" << endl;
1033  }
1034 }
1035 
1036 
1038 {
1039  auto* fldPtr = getObjectPtr<triSurfaceLabelField>("values");
1040 
1041  if (fldPtr)
1042  {
1043  (*fldPtr).field() = values;
1044  }
1045  else
1046  {
1047  fldPtr = new triSurfaceLabelField
1048  (
1049  IOobject
1050  (
1051  "values",
1052  objectRegistry::time().timeName(), // instance
1053  meshSubDir, // local
1054  *this,
1057  ),
1058  *this,
1059  dimless,
1061  );
1062 
1063  // Store field on triMesh
1064  fldPtr->store();
1065  }
1066  if (debug)
1067  {
1068  Pout<< "triSurfaceMesh::setField :"
1069  << " finished setting field for "
1070  << values.size() << " triangles" << endl;
1071  }
1072 }
1073 
1074 
1077  const List<pointIndexHit>& info,
1078  labelList& values
1079 ) const
1080 {
1081  const auto* fldPtr = getObjectPtr<triSurfaceLabelField>("values");
1082 
1083  if (fldPtr)
1084  {
1085  const auto& fld = *fldPtr;
1086 
1087  values.setSize(info.size());
1088 
1089  forAll(info, i)
1090  {
1091  if (info[i].hit())
1092  {
1093  values[i] = fld[info[i].index()];
1094  }
1095  }
1096  }
1097  if (debug)
1098  {
1099  Pout<< "triSurfaceMesh::setField :"
1100  << " finished getting field for "
1101  << info.size() << " triangles" << endl;
1102  }
1103 }
1104 
1105 
1108  const pointField& points,
1109  List<volumeType>& volType
1110 ) const
1111 {
1112  const scalar oldTol = indexedOctree<treeDataTriSurface>::perturbTol();
1114 
1115  if (debug)
1116  {
1117  Pout<< "triSurfaceMesh::getVolumeType :"
1118  << " finding orientation for " << points.size()
1119  << " samples" << endl;
1120  }
1121 
1122  volType.setSize(points.size());
1123 
1124  forAll(points, pointi)
1125  {
1126  const point& pt = points[pointi];
1127 
1128  if (tree().bb().contains(pt))
1129  {
1130  // Use cached volume type per each tree node
1131  volType[pointi] = tree().getVolumeType(pt);
1132  }
1133  else if (hasVolumeType())
1134  {
1135  // Precalculate and cache value for this outside point
1136  if (outsideVolType_ == volumeType::UNKNOWN)
1137  {
1138  outsideVolType_ = tree().shapes().getVolumeType(tree(), pt);
1139  }
1140  volType[pointi] = outsideVolType_;
1141  }
1142  else
1143  {
1144  // Have to calculate directly as outside the octree
1145  volType[pointi] = tree().shapes().getVolumeType(tree(), pt);
1146  }
1147  }
1148 
1150  if (debug)
1151  {
1152  Pout<< "triSurfaceMesh::getVolumeType :"
1153  << " finished finding orientation for " << points.size()
1154  << " samples" << endl;
1155  }
1156 }
1157 
1158 
1162  const bool valid
1163 ) const
1164 {
1166  const fileName& instance = searchableSurface::instance();
1167 
1168  if
1169  (
1170  instance != runTime.timeName()
1171  && instance != runTime.system()
1172  && instance != runTime.caseSystem()
1173  && instance != runTime.constant()
1174  && instance != runTime.caseConstant()
1175  )
1176  {
1177  const_cast<triSurfaceMesh&>(*this).searchableSurface::instance() =
1178  runTime.timeName();
1179  const_cast<triSurfaceMesh&>(*this).objectRegistry::instance() =
1180  runTime.timeName();
1181  }
1182 
1183  fileName fullPath;
1184  if (fName_.size())
1185  {
1186  // Override file name
1187 
1188  fullPath = fName_;
1189 
1190  fullPath.expand();
1191  if (!fullPath.isAbsolute())
1192  {
1193  // Add directory from regIOobject
1194  fullPath = searchableSurface::objectPath().path()/fullPath;
1195  }
1196  }
1197  else
1198  {
1199  fullPath = searchableSurface::objectPath();
1200  }
1201 
1202  if (!mkDir(fullPath.path()))
1203  {
1204  return false;
1205  }
1206 
1207  triSurface::write(fullPath);
1208 
1209  if (!isFile(fullPath))
1210  {
1211  return false;
1212  }
1213 
1214  return true;
1215 }
1216 
1217 
1218 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::PrimitivePatch<::Foam::List< labelledTri >, pointField >::points
const Field< point_type > & points() const noexcept
Return reference to global points.
Definition: PrimitivePatch.H:299
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
runTime
engineTime & runTime
Definition: createEngineTime.H:13
Foam::Random
Random number generator.
Definition: Random.H:59
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::treeBoundBox::extend
treeBoundBox extend(Random &rndGen, const scalar s) const
Return slightly wider bounding box.
Definition: treeBoundBoxI.H:325
Foam::PrimitivePatch<::Foam::List< labelledTri >, pointField >::pointFaces
const labelListList & pointFaces() const
Return point-face addressing.
Definition: PrimitivePatch.C:301
Foam::IOobject::AUTO_WRITE
Definition: IOobject.H:194
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:73
Foam::triSurface::clearOut
void clearOut()
Definition: triSurface.C:566
Foam::IOobject::localFilePath
fileName localFilePath(const word &typeName, const bool search=true) const
Helper for filePath that searches locally.
Definition: IOobject.C:569
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::VectorSpace< Vector< Cmpt >, Cmpt, 3 >::uniform
static Vector< Cmpt > uniform(const Cmpt &s)
Return a VectorSpace with all elements = s.
Definition: VectorSpaceI.H:164
Foam::fileName
A class for handling file names.
Definition: fileName.H:73
s
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.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))
Definition: gmvOutputSpray.H:25
Foam::IOobject::instance
const fileName & instance() const noexcept
Definition: IOobjectI.H:196
Foam::triSurfaceMesh::localOrGlobal
Definition: triSurfaceMesh.H:200
Foam::fileName::path
static std::string path(const std::string &str)
Return directory path name (part before last /)
Definition: fileNameI.H:176
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
PatchTools.H
Foam::triSurface::relativeFilePath
static fileName relativeFilePath(const IOobject &io, const fileName &f, const bool isGlobal=true)
Return fileName.
Definition: triSurfaceIO.C:111
Foam::triSurfaceMesh::getNormal
virtual void getNormal(const List< pointIndexHit > &, vectorField &normal) const
From a set of points and indices get the normal.
Definition: triSurfaceMesh.C:954
Foam::triSurfaceRegionSearch::findNearest
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.
Definition: triSurfaceRegionSearch.C:191
Foam::HashTableOps::values
List< T > values(const HashTable< T, Key, Hash > &tbl, const bool doSort=false)
List of values from HashTable, optionally sorted.
Definition: HashOps.H:149
Foam::triSurfaceMesh::getRegion
virtual void getRegion(const List< pointIndexHit > &, labelList &region) const
From a set of points and indices get the region.
Definition: triSurfaceMesh.C:921
Foam::triSurfaceSearch::findLineAll
void findLineAll(const pointField &start, const pointField &end, List< List< pointIndexHit >> &info) const
Calculate all intersections from start to end.
Definition: triSurfaceSearch.C:379
Foam::SubList
A List obtained as a section of another List.
Definition: SubList.H:54
Foam::treeBoundBox
Standard boundBox with extra functionality for use in octree.
Definition: treeBoundBox.H:86
Foam::PrimitivePatch<::Foam::List< labelledTri >, pointField >::hasFaceCentres
bool hasFaceCentres() const
Definition: PrimitivePatch.H:465
Foam::triSurfaceMesh::edgeTree
const indexedOctree< treeDataEdge > & edgeTree() const
Demand driven construction of octree for boundary edges.
Definition: triSurfaceMesh.C:645
Foam::Time::timeName
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:63
Foam::triSurfaceRegionSearch::clearOut
void clearOut()
Clear storage.
Definition: triSurfaceRegionSearch.C:64
Foam::TimePaths::caseSystem
fileName caseSystem() const
Definition: TimePathsI.H:119
Foam::triSurfaceSearch::findNearest
void findNearest(const pointField &samples, const scalarField &nearestDistSqr, List< pointIndexHit > &info) const
Definition: triSurfaceSearch.C:290
Foam::PatchTools::calcBounds
static void calcBounds(const PrimitivePatch< FaceList, PointField > &p, boundBox &bb, label &nPoints)
Definition: PatchToolsSearch.C:178
Foam::UPstream::master
static bool master(const label communicator=worldComm)
Am I the master process.
Definition: UPstream.H:457
Foam::Pstream::scatter
static void scatter(const List< commsStruct > &comms, T &Value, const int tag, const label comm)
Scatter data. Distribute without modification. Reverse of gather.
Definition: gatherScatter.C:150
triSurfaceFields.H
Fields for triSurface.
Foam::isFile
bool isFile(const fileName &name, const bool checkGzip=true, const bool followLink=true)
Does the name exist as a FILE in the file system?
Definition: MSwindows.C:658
Foam::triSurfaceMesh
IOoject and searching on triSurface.
Definition: triSurfaceMesh.H:106
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::boundBox::max
const point & max() const
Maximum describing the bounding box.
Definition: boundBoxI.H:97
Foam::regIOobject::eventNo
label eventNo() const
Event number at last update.
Definition: regIOobjectI.H:185
Foam::triSurfaceRegionSearch
Helper class to search on triSurface. Creates an octree for each region of the surface and only searc...
Definition: triSurfaceRegionSearch.H:58
Foam::Pout
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Foam::triSurfaceMesh::readAction
readAction
Definition: triSurfaceMesh.H:197
Foam::IOobject::writeOpt
writeOption writeOpt() const noexcept
The write option.
Definition: IOobjectI.H:179
Foam::IOobject::time
const Time & time() const
Return Time associated with the objectRegistry.
Definition: IOobject.C:493
Foam::triSurface::patches
const geometricSurfacePatchList & patches() const noexcept
Definition: triSurface.H:399
Foam::triSurfaceMesh::outsideVolumeType
virtual volumeType outsideVolumeType() const
If surface is closed, what is type of points outside bounds.
Definition: triSurfaceMesh.C:760
Foam::triSurfaceMesh::findLineAny
virtual void findLineAny(const pointField &start, const pointField &end, List< pointIndexHit > &) const
Return any intersection on segment from start to end.
Definition: triSurfaceMesh.C:875
Foam::triSurfaceMesh::triSurfaceMesh
triSurfaceMesh(const triSurfaceMesh &)=delete
No copy construct.
Foam::boundBox::min
const point & min() const
Minimum describing the bounding box.
Definition: boundBoxI.H:91
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::indexedOctree::perturbTol
static scalar & perturbTol()
Get the perturbation tolerance.
Definition: indexedOctree.C:2383
Foam::magSqr
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Foam::invertManyToMany
void invertManyToMany(const label len, const UList< InputIntListType > &input, List< OutputIntListType > &output)
Invert many-to-many.
Definition: ListOpsTemplates.C:720
Foam::triSurfaceMesh::movePoints
virtual void movePoints(const pointField &)
Move points.
Definition: triSurfaceMesh.C:612
nPoints
label nPoints
Definition: gmvOutputHeader.H:2
Foam::objectRegistry
Registry of regIOobjects.
Definition: objectRegistry.H:60
Foam::volumeType
An enumeration wrapper for classification of a location as being inside/outside of a volume.
Definition: volumeType.H:60
Foam::triSurface::movePoints
virtual void movePoints(const pointField &pts)
Move points.
Definition: triSurface.C:612
Foam::Field< vector >
Foam::labelField
Field< label > labelField
Specialisation of Field<T> for label.
Definition: labelField.H:52
Foam::triSurface::write
void write(Ostream &os) const
Write to Ostream in simple OpenFOAM format.
Definition: triSurfaceIO.C:336
Foam::triSurface
Triangulated surface description with patch information.
Definition: triSurface.H:76
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::IOobject::local
const fileName & local() const noexcept
Definition: IOobjectI.H:208
Foam::volumeType::UNKNOWN
Unknown state.
Definition: volumeType.H:67
Foam::volumeType::names
static const Enum< volumeType::type > names
Names for the classification enumeration.
Definition: volumeType.H:76
Foam::List::setSize
void setSize(const label n)
Alias for resize()
Definition: List.H:222
Foam::searchableSurface
Base class of (analytical or triangulated) surface. Encapsulates all the search routines....
Definition: searchableSurface.H:69
coordinates
PtrList< coordinateSystem > coordinates(solidRegions.size())
Foam::IOstreamOption
The IOstreamOption is a simple container for options an IOstream can normally have.
Definition: IOstreamOption.H:63
Foam::triSurfaceMesh::findLine
virtual void findLine(const pointField &start, const pointField &end, List< pointIndexHit > &) const
Find first intersection on segment from start to end.
Definition: triSurfaceMesh.C:852
Foam::indexedOctree
Non-pointer based hierarchical recursive searching.
Definition: treeDataEdge.H:50
Foam::triSurface::transfer
void transfer(triSurface &surf)
Alter contents by transferring (triangles, points) components.
Definition: triSurface.C:961
Foam::treeDataEdge
Holds data for octree to work on an edges subset.
Definition: treeDataEdge.H:56
samples
scalarField samples(nIntervals, Zero)
Foam::triSurfaceMesh::getField
virtual void getField(const List< pointIndexHit > &, labelList &) const
WIP. From a set of hits (points and.
Definition: triSurfaceMesh.C:1076
Foam::triSurfaceMesh::regions
virtual const wordList & regions() const
Names of regions.
Definition: triSurfaceMesh.C:728
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
fld
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< ' ';}gmvFile<< nl;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Definition: gmvOutputLagrangian.H:23
Foam::triSurfaceMesh::boundingSpheres
virtual void boundingSpheres(pointField &centres, scalarField &radiusSqr) const
Get bounding spheres (centre and radius squared). Any point.
Definition: triSurfaceMesh.C:569
timeName
word timeName
Definition: getTimeIndex.H:3
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::triSurfaceMesh::clearOut
void clearOut()
Clear storage.
Definition: triSurfaceMesh.C:530
Foam::IOobject::globalFilePath
fileName globalFilePath(const word &typeName, const bool search=true) const
Helper for filePath that searches up if in parallel.
Definition: IOobject.C:580
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
stdFoam::end
constexpr auto end(C &c) -> decltype(c.end())
Return iterator to the end of the container c.
Definition: stdFoam.H:121
Foam::triSurfaceMesh::points
virtual tmp< pointField > points() const
Get the points that define the surface.
Definition: triSurfaceMesh.C:596
Foam::triSurfaceMesh::coordinates
virtual tmp< pointField > coordinates() const
Get representative set of element coordinates.
Definition: triSurfaceMesh.C:541
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::triSurfaceMesh::writeObject
virtual bool writeObject(IOstreamOption streamOpt, const bool valid) const
Write using stream options.
Definition: triSurfaceMesh.C:1160
Foam::triSurfaceMesh::setField
virtual void setField(const labelList &values)
WIP. Store element-wise field.
Definition: triSurfaceMesh.C:1037
Random.H
Foam::IOobject::readOpt
readOption readOpt() const noexcept
The read option.
Definition: IOobjectI.H:164
Foam::IOobject::name
const word & name() const noexcept
Return name.
Definition: IOobjectI.H:65
Time.H
Foam::EdgeMap< label >
Foam::TimePaths::system
const word & system() const
Return system name.
Definition: TimePathsI.H:102
Foam::triSurfaceSearch::findLine
void findLine(const pointField &start, const pointField &end, List< pointIndexHit > &info) const
Definition: triSurfaceSearch.C:333
edgeHashes.H
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::triSurfaceMesh::~triSurfaceMesh
virtual ~triSurfaceMesh()
Destructor.
Definition: triSurfaceMesh.C:524
forAllConstIters
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
f
labelList f(nPoints)
Foam::BitOps::count
unsigned int count(const UList< bool > &bools, const bool val=true)
Count number of 'true' entries.
Definition: BitOps.H:77
Foam::triSurfaceMesh::findNearest
virtual void findNearest(const pointField &sample, const scalarField &nearestDistSqr, List< pointIndexHit > &) const
Definition: triSurfaceMesh.C:795
Foam::Vector< scalar >
Foam::UPstream::parRun
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:433
Foam::List< labelList >
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::TimePaths::caseConstant
fileName caseConstant() const
Definition: TimePathsI.H:108
Foam::triSurfaceMesh::addFaceToEdge
static bool addFaceToEdge(const edge &, EdgeMap< label > &)
Helper function for isSurfaceClosed.
Definition: triSurfaceMesh.C:51
Foam::string::expand
string & expand(const bool allowEmpty=false)
Definition: string.C:173
Foam::triSurfaceMesh::meshSubDir
static word meshSubDir
Return the mesh sub-directory name (usually "triSurface")
Definition: triSurfaceMesh.H:174
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::triSurfaceMesh::overlaps
virtual bool overlaps(const boundBox &bb) const
Does any part of the surface overlap the supplied bound box?
Definition: triSurfaceMesh.C:602
Foam::identity
labelList identity(const label len, label start=0)
Create identity map of the given length with (map[i] == i)
Definition: labelList.C:38
Foam::List::clear
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:116
Foam::tmp::New
static tmp< T > New(Args &&... args)
Construct tmp of T with forwarding arguments.
Foam::word::null
static const word null
An empty word.
Definition: word.H:80
patches
const polyBoundaryMesh & patches
Definition: convertProcessorPatches.H:65
Foam::triSurfaceLabelField
Foam::DimensionedField< label, triSurfaceGeoMesh > triSurfaceLabelField
Definition: triSurfaceFieldsFwd.H:44
Foam::labelledTri
A triFace with additional (region) index.
Definition: labelledTri.H:57
Foam::triSurfaceMesh::isSurfaceClosed
bool isSurfaceClosed() const
Check whether surface is closed without calculating any permanent.
Definition: triSurfaceMesh.C:67
Foam::boundBox
A bounding box defined in terms of min/max extrema points.
Definition: boundBox.H:63
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
rndGen
Random rndGen
Definition: createFields.H:23
Foam::triSurfaceSearch::findLineAny
void findLineAny(const pointField &start, const pointField &end, List< pointIndexHit > &info) const
Definition: triSurfaceSearch.C:356
Foam::searchableSurface::bounds
virtual const boundBox & bounds() const
Return const reference to boundBox.
Definition: searchableSurface.H:179
Foam::keyType::LITERAL
String literal.
Definition: keyType.H:81
Foam::triSurfaceMesh::getVolumeType
virtual void getVolumeType(const pointField &points, List< volumeType > &volType) const
Determine type (inside/outside/mixed) for point.
Definition: triSurfaceMesh.C:1107
Foam::IOobject::objectPath
fileName objectPath() const
The complete path + object name.
Definition: IOobjectI.H:214
Foam::triSurfaceMesh::masterOnly
Definition: triSurfaceMesh.H:201
Foam::TimePaths::constant
const word & constant() const
Return constant name.
Definition: TimePathsI.H:96
Foam::dictionary::getOrDefault
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:148
Foam::IOobject::NO_READ
Definition: IOobject.H:188
Foam::mkDir
bool mkDir(const fileName &pathName, mode_t mode=0777)
Make a directory and return an error if it could not be created.
Definition: MSwindows.C:507
triSurfaceMesh.H
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::triSurfaceMesh::hasVolumeType
virtual bool hasVolumeType() const
Whether supports volume type (below) - i.e. whether is closed.
Definition: triSurfaceMesh.C:742
Foam::PrimitivePatch<::Foam::List< labelledTri >, pointField >::faceCentres
const Field< point_type > & faceCentres() const
Return face centres for patch.
Definition: PrimitivePatch.C:400
Foam::fileName::isAbsolute
static bool isAbsolute(const std::string &str)
Definition: fileNameI.H:136
Foam::objectRegistry::time
const Time & time() const noexcept
Return time registry.
Definition: objectRegistry.H:178
pFaces
Info<< "Finished reading KIVA file"<< endl;cellShapeList cellShapes(nPoints);labelList cellZoning(nPoints, -1);const cellModel &hex=cellModel::ref(cellModel::HEX);labelList hexLabels(8);label activeCells=0;labelList pointMap(nPoints);forAll(pointMap, i){ pointMap[i]=i;}for(label i=0;i< nPoints;i++){ if(f[i] > 0.0) { hexLabels[0]=i;hexLabels[1]=i1tab[i];hexLabels[2]=i3tab[i1tab[i]];hexLabels[3]=i3tab[i];hexLabels[4]=i8tab[i];hexLabels[5]=i1tab[i8tab[i]];hexLabels[6]=i3tab[i1tab[i8tab[i]]];hexLabels[7]=i3tab[i8tab[i]];cellShapes[activeCells].reset(hex, hexLabels);edgeList edges=cellShapes[activeCells].edges();forAll(edges, ei) { if(edges[ei].mag(points)< SMALL) { label start=pointMap[edges[ei].start()];while(start !=pointMap[start]) { start=pointMap[start];} label end=pointMap[edges[ei].end()];while(end !=pointMap[end]) { end=pointMap[end];} label minLabel=min(start, end);pointMap[start]=pointMap[end]=minLabel;} } cellZoning[activeCells]=idreg[i];activeCells++;}}cellShapes.setSize(activeCells);cellZoning.setSize(activeCells);forAll(cellShapes, celli){ cellShape &cs=cellShapes[celli];forAll(cs, i) { cs[i]=pointMap[cs[i]];} cs.collapse();}label bcIDs[11]={-1, 0, 2, 4, -1, 5, -1, 6, 7, 8, 9};const label nBCs=12;const word *kivaPatchTypes[nBCs]={ &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &symmetryPolyPatch::typeName, &wedgePolyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &symmetryPolyPatch::typeName, &oldCyclicPolyPatch::typeName};enum patchTypeNames{ PISTON, VALVE, LINER, CYLINDERHEAD, AXIS, WEDGE, INFLOW, OUTFLOW, PRESIN, PRESOUT, SYMMETRYPLANE, CYCLIC};const char *kivaPatchNames[nBCs]={ "piston", "valve", "liner", "cylinderHead", "axis", "wedge", "inflow", "outflow", "presin", "presout", "symmetryPlane", "cyclic"};List< SLList< face > > pFaces[nBCs]
Definition: readKivaGrid.H:235
Foam::triSurfaceMesh::findLineAll
virtual void findLineAll(const pointField &start, const pointField &end, List< List< pointIndexHit >> &) const
Get all intersections in order from start to end.
Definition: triSurfaceMesh.C:898
Foam::dictionary::readIfPresent
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:405
Foam::dimless
const dimensionSet dimless
Dimensionless.
Definition: dimensionSets.C:189
Foam::PrimitivePatch
A list of faces which address into the list of points.
Definition: PrimitivePatch.H:79
Foam::IOobject::db
const objectRegistry & db() const noexcept
Return the local objectRegistry.
Definition: IOobject.C:487