surfaceInflate.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) 2015 OpenFOAM Foundation
9  Copyright (C) 2020-2023 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 Application
28  surfaceInflate
29 
30 Group
31  grpSurfaceUtilities
32 
33 Description
34  Inflates surface. WIP. Checks for overlaps and locally lowers
35  inflation distance
36 
37 Usage
38  - surfaceInflate [OPTION]
39 
40  \param -checkSelfIntersection \n
41  Includes checks for self-intersection.
42 
43  \param -nSmooth
44  Specify number of smoothing iterations
45 
46  \param -featureAngle
47  Specify a feature angle
48 
49 
50  E.g. inflate surface by 20mm with 1.5 safety factor:
51  surfaceInflate DTC-scaled.obj 0.02 1.5 -featureAngle 45 -nSmooth 2
52 
53 \*---------------------------------------------------------------------------*/
54 
55 #include "argList.H"
56 #include "Time.H"
57 #include "triSurfaceFields.H"
58 #include "triSurfaceMesh.H"
59 #include "triSurfaceGeoMesh.H"
60 #include "bitSet.H"
61 #include "OBJstream.H"
62 #include "surfaceFeatures.H"
63 
64 using namespace Foam;
65 
66 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
67 
68 scalar calcVertexNormalWeight
69 (
70  const triFace& f,
71  const label pI,
72  const vector& fN,
73  const pointField& points
74 )
75 {
76  label index = f.find(pI);
77 
78  if (index == -1)
79  {
81  << "Point not in face" << abort(FatalError);
82  }
83 
84  const vector e1 = points[f[index]] - points[f[f.fcIndex(index)]];
85  const vector e2 = points[f[index]] - points[f[f.rcIndex(index)]];
86 
87  return mag(fN)/(magSqr(e1)*magSqr(e2) + VSMALL);
88 }
89 
90 
91 tmp<vectorField> calcVertexNormals(const triSurface& surf)
92 {
93  // Weighted average of normals of faces attached to the vertex
94  // Weight = fA / (mag(e0)^2 * mag(e1)^2);
95  auto tpointNormals = tmp<vectorField>::New(surf.nPoints(), Zero);
96  auto& pointNormals = tpointNormals.ref();
97 
98  const pointField& points = surf.points();
99  const labelListList& pointFaces = surf.pointFaces();
100  const labelList& meshPoints = surf.meshPoints();
101 
102  forAll(pointFaces, pI)
103  {
104  const labelList& pFaces = pointFaces[pI];
105 
106  forAll(pFaces, fI)
107  {
108  const label faceI = pFaces[fI];
109  const triFace& f = surf[faceI];
110 
111  vector areaNorm = f.areaNormal(points);
112 
113  scalar weight = calcVertexNormalWeight
114  (
115  f,
116  meshPoints[pI],
117  areaNorm,
118  points
119  );
120 
121  pointNormals[pI] += weight * areaNorm;
122  }
123 
124  pointNormals[pI].normalise();
125  }
126 
127  return tpointNormals;
128 }
129 
130 
131 tmp<vectorField> calcPointNormals
132 (
133  const triSurface& s,
134  const bitSet& isFeaturePoint,
135  const List<surfaceFeatures::edgeStatus>& edgeStat
136 )
137 {
138  //const pointField pointNormals(s.pointNormals());
139  tmp<vectorField> tpointNormals(calcVertexNormals(s));
140  vectorField& pointNormals = tpointNormals.ref();
141 
142 
143  // feature edges: create edge normals from edgeFaces only.
144  {
145  const labelListList& edgeFaces = s.edgeFaces();
146 
147  labelList nNormals(s.nPoints(), Zero);
148  forAll(edgeStat, edgeI)
149  {
150  if (edgeStat[edgeI] != surfaceFeatures::NONE)
151  {
152  const edge& e = s.edges()[edgeI];
153  forAll(e, i)
154  {
155  if (!isFeaturePoint.test(e[i]))
156  {
157  pointNormals[e[i]] = Zero;
158  }
159  }
160  }
161  }
162 
163  forAll(edgeStat, edgeI)
164  {
165  if (edgeStat[edgeI] != surfaceFeatures::NONE)
166  {
167  const labelList& eFaces = edgeFaces[edgeI];
168 
169  // Get average edge normal
170  vector n = Zero;
171  for (const label facei : eFaces)
172  {
173  n += s.faceNormals()[facei];
174  }
175  n /= eFaces.size();
176 
177 
178  // Sum to points
179  const edge& e = s.edges()[edgeI];
180  forAll(e, i)
181  {
182  if (!isFeaturePoint.test(e[i]))
183  {
184  pointNormals[e[i]] += n;
185  nNormals[e[i]]++;
186  }
187  }
188  }
189  }
190 
191  forAll(nNormals, pointI)
192  {
193  if (nNormals[pointI] > 0)
194  {
195  pointNormals[pointI].normalise();
196  }
197  }
198  }
199 
200 
201  forAll(pointNormals, pointI)
202  {
203  if (mag(mag(pointNormals[pointI])-1) > SMALL)
204  {
206  << "unitialised"
207  << exit(FatalError);
208  }
209  }
210 
211  return tpointNormals;
212 }
213 
214 
215 void detectSelfIntersections
216 (
217  const triSurfaceMesh& s,
218  bitSet& isEdgeIntersecting
219 )
220 {
221  const edgeList& edges = s.edges();
222  const indexedOctree<treeDataTriSurface>& tree = s.tree();
223  const labelList& meshPoints = s.meshPoints();
224  const tmp<pointField> tpoints(s.points());
225  const pointField& points = tpoints();
226 
227  isEdgeIntersecting.resize_nocopy(edges.size());
228  isEdgeIntersecting = false;
229 
230  forAll(edges, edgeI)
231  {
232  const edge& e = edges[edgeI];
233 
234  pointIndexHit hitInfo
235  (
236  tree.findLine
237  (
238  points[meshPoints[e[0]]],
239  points[meshPoints[e[1]]],
241  )
242  );
243 
244  if (hitInfo.hit())
245  {
246  isEdgeIntersecting.set(edgeI);
247  }
248  }
249 }
250 
251 
252 label detectIntersectionPoints
253 (
254  const scalar tolerance,
255  const triSurfaceMesh& s,
256 
257  const scalar extendFactor,
258  const pointField& initialPoints,
259  const vectorField& displacement,
260 
261  const bool checkSelfIntersect,
262  const bitSet& initialIsEdgeIntersecting,
263 
264  bitSet& isPointOnHitEdge,
265  scalarField& scale
266 )
267 {
268  const pointField initialLocalPoints(initialPoints, s.meshPoints());
269  const List<labelledTri>& localFaces = s.localFaces();
270 
271 
272  label nHits = 0;
273  isPointOnHitEdge.setSize(s.nPoints());
274  isPointOnHitEdge = false;
275 
276 
277  // 1. Extrusion offset vectors intersecting new surface location
278  {
279  scalar tol = max(tolerance, 10*s.tolerance());
280 
281  // Collect all the edge vectors. Slightly shorten the edges to prevent
282  // finding lots of intersections. The fast triangle intersection routine
283  // has problems with rays passing through a point of the
284  // triangle.
285  pointField start(initialLocalPoints+tol*displacement);
286  pointField end(initialLocalPoints+extendFactor*displacement);
287 
288  List<pointIndexHit> hits;
289  s.findLineAny(start, end, hits);
290 
291  forAll(hits, pointI)
292  {
293  if
294  (
295  hits[pointI].hit()
296  && !localFaces[hits[pointI].index()].found(pointI)
297  )
298  {
299  scale[pointI] = max(0.0, scale[pointI]-0.2);
300 
301  isPointOnHitEdge.set(pointI);
302  nHits++;
303  }
304  }
305  }
306 
307 
308  // 2. (new) surface self intersections
309  if (checkSelfIntersect)
310  {
311  bitSet isEdgeIntersecting;
312  detectSelfIntersections(s, isEdgeIntersecting);
313 
314  const edgeList& edges = s.edges();
315  const tmp<pointField> tpoints(s.points());
316  const pointField& points = tpoints();
317 
318  forAll(edges, edgeI)
319  {
320  const edge& e = edges[edgeI];
321 
322  if (isEdgeIntersecting[edgeI] && !initialIsEdgeIntersecting[edgeI])
323  {
324  if (isPointOnHitEdge.set(e[0]))
325  {
326  label start = s.meshPoints()[e[0]];
327  const point& pt = points[start];
328 
329  Pout<< "Additional self intersection at "
330  << pt
331  << endl;
332 
333  scale[e[0]] = max(0.0, scale[e[0]]-0.2);
334  nHits++;
335  }
336  if (isPointOnHitEdge.set(e[1]))
337  {
338  label end = s.meshPoints()[e[1]];
339  const point& pt = points[end];
340 
341  Pout<< "Additional self intersection at "
342  << pt
343  << endl;
344 
345  scale[e[1]] = max(0.0, scale[e[1]]-0.2);
346  nHits++;
347  }
348  }
349  }
350  }
351 
352  return nHits;
353 }
354 
355 
357 (
358  const triSurface& s,
359  const scalarField& fld,
360  const scalarField& edgeWeights
361 )
362 {
363  auto tres = tmp<scalarField>::New(s.nPoints(), Zero);
364  auto& res = tres.ref();
365 
366  scalarField sumWeight(s.nPoints(), Zero);
367 
368  const edgeList& edges = s.edges();
369 
370  forAll(edges, edgeI)
371  {
372  const edge& e = edges[edgeI];
373  const scalar w = edgeWeights[edgeI];
374 
375  res[e[0]] += w*fld[e[1]];
376  sumWeight[e[0]] += w;
377 
378  res[e[1]] += w*fld[e[0]];
379  sumWeight[e[1]] += w;
380  }
381 
382  // Average
383  // ~~~~~~~
384 
385  forAll(res, pointI)
386  {
387  if (mag(sumWeight[pointI]) < VSMALL)
388  {
389  // Unconnected point. Take over original value
390  res[pointI] = fld[pointI];
391  }
392  else
393  {
394  res[pointI] /= sumWeight[pointI];
395  }
396  }
397 
398  return tres;
399 }
400 
401 
402 void minSmooth
403 (
404  const triSurface& s,
405  const bitSet& isAffectedPoint,
406  const scalarField& fld,
407  scalarField& newFld
408 )
409 {
410 
411  const edgeList& edges = s.edges();
412  const pointField& points = s.points();
413  const labelList& mp = s.meshPoints();
414 
415  scalarField edgeWeights(edges.size());
416  forAll(edges, edgeI)
417  {
418  const edge& e = edges[edgeI];
419  scalar w = mag(points[mp[e[0]]]-points[mp[e[1]]]);
420 
421  edgeWeights[edgeI] = 1.0/(max(w, SMALL));
422  }
423 
424  tmp<scalarField> tavgFld = avg(s, fld, edgeWeights);
425 
426  const scalarField& avgFld = tavgFld();
427 
428  forAll(fld, pointI)
429  {
430  if (isAffectedPoint.test(pointI))
431  {
432  newFld[pointI] = min
433  (
434  fld[pointI],
435  0.5*fld[pointI] + 0.5*avgFld[pointI]
436  //avgFld[pointI]
437  );
438  }
439  }
440 }
441 
442 
443 void lloydsSmoothing
444 (
445  const label nSmooth,
446  triSurface& s,
447  const bitSet& isFeaturePoint,
448  const List<surfaceFeatures::edgeStatus>& edgeStat,
449  const bitSet& isAffectedPoint
450 )
451 {
452  const labelList& meshPoints = s.meshPoints();
453  const edgeList& edges = s.edges();
454 
455 
456  bitSet isSmoothPoint(isAffectedPoint);
457  // Extend isSmoothPoint
458  {
459  bitSet newIsSmoothPoint(isSmoothPoint);
460  forAll(edges, edgeI)
461  {
462  const edge& e = edges[edgeI];
463  if (isSmoothPoint.test(e[0]))
464  {
465  newIsSmoothPoint.set(e[1]);
466  }
467  if (isSmoothPoint.test(e[1]))
468  {
469  newIsSmoothPoint.set(e[0]);
470  }
471  }
472  Info<< "From points-to-smooth " << isSmoothPoint.count()
473  << " to " << newIsSmoothPoint.count()
474  << endl;
475  isSmoothPoint.transfer(newIsSmoothPoint);
476  }
477 
478  // Do some smoothing (Lloyds algorithm) around problematic points
479  for (label i = 0; i < nSmooth; i++)
480  {
481  const labelListList& pointFaces = s.pointFaces();
482  const pointField& faceCentres = s.faceCentres();
483 
484  pointField newPoints(s.points());
485  forAll(isSmoothPoint, pointI)
486  {
487  if (isSmoothPoint[pointI] && !isFeaturePoint[pointI])
488  {
489  const labelList& pFaces = pointFaces[pointI];
490 
491  point avg(Zero);
492  forAll(pFaces, pFaceI)
493  {
494  avg += faceCentres[pFaces[pFaceI]];
495  }
496  avg /= pFaces.size();
497  newPoints[meshPoints[pointI]] = avg;
498  }
499  }
500 
501 
502  // Move points on feature edges only according to feature edges.
503 
504  const pointField& points = s.points();
505 
506  vectorField pointSum(s.nPoints(), Zero);
507  labelList nPointSum(s.nPoints(), Zero);
508 
509  forAll(edges, edgeI)
510  {
511  if (edgeStat[edgeI] != surfaceFeatures::NONE)
512  {
513  const edge& e = edges[edgeI];
514  const point& start = points[meshPoints[e[0]]];
515  const point& end = points[meshPoints[e[1]]];
516  point eMid = 0.5*(start+end);
517  pointSum[e[0]] += eMid;
518  nPointSum[e[0]]++;
519  pointSum[e[1]] += eMid;
520  nPointSum[e[1]]++;
521  }
522  }
523 
524  forAll(pointSum, pointI)
525  {
526  if
527  (
528  isSmoothPoint[pointI]
529  && isFeaturePoint[pointI]
530  && nPointSum[pointI] >= 2
531  )
532  {
533  newPoints[meshPoints[pointI]] =
534  pointSum[pointI]/nPointSum[pointI];
535  }
536  }
537 
538 
539  s.movePoints(newPoints);
540 
541 
542  // Extend isSmoothPoint
543  {
544  bitSet newIsSmoothPoint(isSmoothPoint);
545  forAll(edges, edgeI)
546  {
547  const edge& e = edges[edgeI];
548  if (isSmoothPoint[e[0]])
549  {
550  newIsSmoothPoint.set(e[1]);
551  }
552  if (isSmoothPoint[e[1]])
553  {
554  newIsSmoothPoint.set(e[0]);
555  }
556  }
557  Info<< "From points-to-smooth " << isSmoothPoint.count()
558  << " to " << newIsSmoothPoint.count()
559  << endl;
560  isSmoothPoint.transfer(newIsSmoothPoint);
561  }
562  }
563 }
564 
565 
566 
567 // Main program:
568 
569 int main(int argc, char *argv[])
570 {
572  (
573  "Inflates surface according to point normals."
574  );
575 
578  (
579  "Creates inflated version of surface using point normals. "
580  "Takes surface, distance to inflate and additional safety factor"
581  );
583  (
584  "checkSelfIntersection",
585  "Also check for self-intersection"
586  );
588  (
589  "nSmooth",
590  "integer",
591  "Number of smoothing iterations (default 20)"
592  );
594  (
595  "featureAngle",
596  "scalar",
597  "Feature angle"
598  );
600  (
601  "debug",
602  "Switch on additional debug information"
603  );
604 
605  argList::addArgument("input", "The input surface file");
606  argList::addArgument("distance", "The inflate distance");
607  argList::addArgument("factor", "The extend safety factor [1,10]");
608 
609  argList::noFunctionObjects(); // Never use function objects
610 
611  #include "setRootCase.H"
612  #include "createTime.H"
613 
614  const auto inputName = args.get<word>(1);
615  const auto distance = args.get<scalar>(2);
616  const auto extendFactor = args.get<scalar>(3);
617  const bool checkSelfIntersect = args.found("checkSelfIntersection");
618  const auto nSmooth = args.getOrDefault<label>("nSmooth", 10);
619  const auto featureAngle = args.getOrDefault<scalar>("featureAngle", 180);
620  const bool debug = args.found("debug");
621 
622 
623  Info<< "Inflating surface " << inputName
624  << " according to point normals" << nl
625  << " distance : " << distance << nl
626  << " safety factor : " << extendFactor << nl
627  << " self intersection test : " << checkSelfIntersect << nl
628  << " smoothing iterations : " << nSmooth << nl
629  << " feature angle : " << featureAngle << nl
630  << endl;
631 
632 
633  if (extendFactor < 1 || extendFactor > 10)
634  {
636  << "Illegal safety factor " << extendFactor
637  << ". It is usually 1..2"
638  << exit(FatalError);
639  }
640 
641 
642 
643  // Load triSurfaceMesh
645  (
646  IOobject
647  (
648  inputName, // name
649  runTime.constant(), // instance
650  "triSurface", // local
651  runTime, // registry
654  )
655  );
656 
657  // Mark features
658  const surfaceFeatures features(s, 180.0-featureAngle);
659 
660  Info<< "Detected features:" << nl
661  << " feature points : " << features.featurePoints().size()
662  << " out of " << s.nPoints() << nl
663  << " feature edges : " << features.featureEdges().size()
664  << " out of " << s.nEdges() << nl
665  << endl;
666 
667  bitSet isFeaturePoint(s.nPoints(), features.featurePoints());
668 
669  const List<surfaceFeatures::edgeStatus> edgeStat(features.toStatus());
670 
671 
672  const pointField initialPoints(s.points());
673 
674 
675  // Construct scale
676  Info<< "Constructing field scale\n" << endl;
678  (
679  IOobject
680  (
681  "scale", // name
682  runTime.timeName(), // instance
683  s, // registry
686  ),
687  s,
688  scalar(1),
689  dimLength
690  );
691 
692 
693  // Construct unit normals
694 
695  Info<< "Calculating vertex normals\n" << endl;
696  const pointField pointNormals
697  (
698  calcPointNormals
699  (
700  s,
701  isFeaturePoint,
702  edgeStat
703  )
704  );
705 
706 
707  // Construct pointDisplacement
708  Info<< "Constructing field pointDisplacement\n" << endl;
709  triSurfacePointVectorField pointDisplacement
710  (
711  IOobject
712  (
713  "pointDisplacement", // name
714  runTime.timeName(), // instance
715  s, // registry
718  ),
719  s,
720  dimLength,
721  // - calculate with primitive fields (#2758)
722  (distance*scale.field()*pointNormals)
723  );
724 
725 
726  const labelList& meshPoints = s.meshPoints();
727 
728 
729  // Any point on any intersected edge in any of the iterations
730  bitSet isScaledPoint(s.nPoints());
731 
732 
733  // Detect any self intersections on initial mesh
734  bitSet initialIsEdgeIntersecting;
735  if (checkSelfIntersect)
736  {
737  detectSelfIntersections(s, initialIsEdgeIntersecting);
738  }
739  else
740  {
741  // Mark all edges as already self intersecting so avoid detecting any
742  // new ones
743  initialIsEdgeIntersecting.setSize(s.nEdges(), true);
744  }
745 
746 
747  // Inflate
748  while (runTime.loop())
749  {
750  Info<< "Time = " << runTime.timeName() << nl << endl;
751 
752  // Move to new location
753  pointField newPoints(initialPoints);
754  forAll(meshPoints, i)
755  {
756  newPoints[meshPoints[i]] += pointDisplacement[i];
757  }
758  s.movePoints(newPoints);
759  Info<< "Bounding box : " << s.bounds() << endl;
760 
761 
762  // Work out scale from pointDisplacement
763  forAll(scale, pointI)
764  {
765  if (s.pointFaces()[pointI].size())
766  {
767  scale[pointI] = mag(pointDisplacement[pointI])/distance;
768  }
769  else
770  {
771  scale[pointI] = 1.0;
772  }
773  }
774 
775 
776  Info<< "Scale : " << gAverage(scale) << endl;
777  Info<< "Displacement : " << gAverage(pointDisplacement) << endl;
778 
779 
780 
781  // Detect any intersections and scale back
782  bitSet isAffectedPoint;
783  label nIntersections = detectIntersectionPoints
784  (
785  1e-9, // intersection tolerance
786  s,
787  extendFactor,
788  initialPoints,
789  pointDisplacement,
790 
791  checkSelfIntersect,
792  initialIsEdgeIntersecting,
793 
794  isAffectedPoint,
795  scale
796  );
797  Info<< "Detected " << nIntersections << " intersections" << nl
798  << endl;
799 
800  if (nIntersections == 0)
801  {
802  runTime.write();
803  break;
804  }
805 
806 
807  // Accumulate all affected points
808  forAll(isAffectedPoint, pointI)
809  {
810  if (isAffectedPoint.test(pointI))
811  {
812  isScaledPoint.set(pointI);
813  }
814  }
815 
816  // Smear out lowering of scale so any edges not found are
817  // still included
818  for (label i = 0; i < nSmooth; i++)
819  {
820  triSurfacePointScalarField oldScale(scale);
821  oldScale.rename("oldScale");
822  minSmooth
823  (
824  s,
825  bitSet(s.nPoints(), true),
826  oldScale,
827  scale
828  );
829  }
830 
831 
832  // From scale update the pointDisplacement
833  // - calculate with primitive fields (#2758)
834  pointDisplacement.normalise();
835  pointDisplacement.field() *= distance*scale.field();
836 
837 
838  // Do some smoothing (Lloyds algorithm)
839  lloydsSmoothing(nSmooth, s, isFeaturePoint, edgeStat, isAffectedPoint);
840 
841  // Update pointDisplacement
842  const tmp<pointField> tpoints(s.points());
843  const pointField& pts = tpoints();
844 
845  forAll(meshPoints, i)
846  {
847  label meshPointI = meshPoints[i];
848  pointDisplacement[i] = pts[meshPointI]-initialPoints[meshPointI];
849  }
850 
851 
852  // Write
853  runTime.write();
854 
856  }
857 
858 
859  Info<< "Detected overall intersecting points " << isScaledPoint.count()
860  << " out of " << s.nPoints() << nl << endl;
861 
862 
863  if (debug)
864  {
865  OBJstream str(runTime.path()/"isScaledPoint.obj");
866 
867  for (const label pointi : isScaledPoint)
868  {
869  str.write(initialPoints[meshPoints[pointi]]);
870  }
871  }
872 
873 
874  Info<< "End\n" << endl;
875 
876  return 0;
877 }
878 
879 
880 // ************************************************************************* //
label nPoints() const
Number of points supporting patch faces.
static void noFunctionObjects(bool addWithOption=false)
Remove &#39;-noFunctionObjects&#39; option and ignore any occurrences.
Definition: argList.C:547
static void addNote(const string &note)
Add extra notes for the usage information.
Definition: argList.C:462
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
fileName path() const
Return path = rootPath/caseName. Same as TimePaths::path()
Definition: Time.H:503
void set(const bitSet &bitset)
Set specified bits from another bitset.
Definition: bitSetI.H:498
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
virtual Ostream & write(const char c) override
Write character.
Definition: OBJstream.C:69
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:608
Fields for triSurface.
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition: tmpI.H:235
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
std::enable_if< std::is_same< bool, TypeT >::value, bool >::type set(const label i, bool val=true)
A bitSet::set() method for a list of bool.
Definition: List.H:493
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
void setSize(const label n, unsigned int val=0u)
Alias for resize()
Definition: PackedList.H:807
virtual bool loop()
Return true if run should continue and if so increment time.
Definition: Time.C:856
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
static void addBoolOption(const word &optName, const string &usage="", bool advanced=false)
Add a bool option to validOptions with usage information.
Definition: argList.C:374
static void noParallel()
Remove the parallel options.
Definition: argList.C:584
scalar distance(const vector &p1, const vector &p2)
Definition: curveTools.C:12
This class describes the interaction of an object (often a face) and a point. It carries the info of ...
Definition: pointIndexHit.H:44
T getOrDefault(const word &optName, const T &deflt) const
Get a value from the named option if present, or return default.
Definition: argListI.H:300
const labelList & meshPoints() const
Return labelList of mesh points in patch.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
IOoject and searching on triSurface.
label fcIndex(const label i) const noexcept
The forward circular index. The next index in the list which returns to the first at the end of the l...
Definition: UListI.H:90
Not a classified feature edge.
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
const pointField & points
An edge is a list of two vertex labels. This can correspond to a directed graph edge or an edge on a ...
Definition: edge.H:59
A class for handling words, derived from Foam::string.
Definition: word.H:63
A triangular face using a FixedList of labels corresponding to mesh vertices.
Definition: triFace.H:65
Tree tree(triangles.begin(), triangles.end())
Reading is optional [identical to LAZY_READ].
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
static void addOption(const word &optName, const string &param="", const string &usage="", bool advanced=false)
Add an option to validOptions with usage information.
Definition: argList.C:385
const Field< point_type > & points() const noexcept
Return reference to global points.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
errorManip< error > abort(error &err)
Definition: errorManip.H:139
bool test(const label pos) const
Test for True value at specified position, never auto-vivify entries.
Definition: bitSet.H:329
label find(const T &val) const
Find index of the first occurrence of the value.
Definition: UList.C:173
const labelListList & pointFaces() const
Return point-face addressing.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
An OFstream that keeps track of vertices and provides convenience output methods for OBJ files...
Definition: OBJstream.H:55
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition: Time.C:714
const word & constant() const noexcept
Return constant name.
Definition: TimePathsI.H:112
int debug
Static debugging option.
constexpr auto end(C &c) -> decltype(c.end())
Return iterator to the end of the container c.
Definition: stdFoam.H:201
Ostream & printExecutionTime(OSstream &os) const
Print the elapsed ExecutionTime (cpu-time), ClockTime.
Definition: TimeIO.C:607
labelList f(nPoints)
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
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))
T get(const label index) const
Get a value from the argument at index.
Definition: argListI.H:271
Non-pointer based hierarchical recursive searching.
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:59
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
Type gAverage(const FieldField< Field, Type > &f)
void resize_nocopy(const label numElem)
Currently identical to resize. Subject to future change (Oct-2021)
Definition: PackedListI.H:445
Nothing to be read.
Automatically write from objectRegistry::writeObject()
static void addArgument(const string &argName, const string &usage="")
Append a (mandatory) argument to validArgs.
Definition: argList.C:351
messageStream Info
Information stream (stdout output on master, null elsewhere)
label n
A class for managing temporary objects.
Definition: HashPtrTable.H:50
void normalise()
Inplace normalise this field. Generally a no-op except for vector fields.
Definition: Field.C:601
Triangulated surface description with patch information.
Definition: triSurface.H:71
label rcIndex(const label i) const noexcept
The reverse circular index. The previous index in the list which returns to the last at the beginning...
Definition: UListI.H:97
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))
Foam::argList args(argc, argv)
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:180
bool found
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
const dimensionedScalar mp
Proton mass.
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
bool found(const word &optName) const
Return true if the named option is found.
Definition: argListI.H:171
Holds feature edges/points of surface.
Namespace for OpenFOAM.
const pointField & pts
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127