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 | Copyright (C) 2015 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 Application
25  surfaceInflate
26 
27 Description
28  Inflates surface. WIP. Checks for overlaps and locally lowers
29  inflation distance
30 
31 Usage
32  - surfaceInflate [OPTION]
33 
34  \param -checkSelfIntersection \n
35  Includes checks for self-intersection.
36 
37  \param -nSmooth
38  Specify number of smoothing iterations
39 
40  \param -featureAngle
41  Specify a feature angle
42 
43 
44  E.g. inflate surface by 2cm with 1.5 safety factor:
45  surfaceInflate DTC-scaled.obj 0.02 1.5 -featureAngle 45 -nSmooth 2
46 
47 \*---------------------------------------------------------------------------*/
48 
49 #include "argList.H"
50 #include "Time.H"
51 #include "triSurfaceFields.H"
52 #include "triSurfaceMesh.H"
53 #include "triSurfaceGeoMesh.H"
54 #include "PackedBoolList.H"
55 #include "OBJstream.H"
56 #include "surfaceFeatures.H"
57 
58 using namespace Foam;
59 
60 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
61 
63 (
64  const triFace& f,
65  const label pI,
66  const vector& fN,
67  const pointField& points
68 )
69 {
70  label index = findIndex(f, pI);
71 
72  if (index == -1)
73  {
75  << "Point not in face" << abort(FatalError);
76  }
77 
78  const vector e1 = points[f[index]] - points[f[f.fcIndex(index)]];
79  const vector e2 = points[f[index]] - points[f[f.rcIndex(index)]];
80 
81  return mag(fN)/(magSqr(e1)*magSqr(e2) + VSMALL);
82 }
83 
84 
86 {
87  // Weighted average of normals of faces attached to the vertex
88  // Weight = fA / (mag(e0)^2 * mag(e1)^2);
89  tmp<vectorField> tpointNormals
90  (
91  new pointField(surf.nPoints(), vector::zero)
92  );
93  vectorField& pointNormals = tpointNormals();
94 
95  const pointField& points = surf.points();
96  const labelListList& pointFaces = surf.pointFaces();
97  const labelList& meshPoints = surf.meshPoints();
98 
99  forAll(pointFaces, pI)
100  {
101  const labelList& pFaces = pointFaces[pI];
102 
103  forAll(pFaces, fI)
104  {
105  const label faceI = pFaces[fI];
106  const triFace& f = surf[faceI];
107 
108  vector fN = f.normal(points);
109 
110  scalar weight = calcVertexNormalWeight
111  (
112  f,
113  meshPoints[pI],
114  fN,
115  points
116  );
117 
118  pointNormals[pI] += weight*fN;
119  }
120 
121  pointNormals[pI] /= mag(pointNormals[pI]) + VSMALL;
122  }
123 
124  return tpointNormals;
125 }
126 
127 
128 tmp<vectorField> calcPointNormals
129 (
130  const triSurface& s,
131  const PackedBoolList& isFeaturePoint,
132  const List<surfaceFeatures::edgeStatus>& edgeStat
133 )
134 {
135  //const pointField pointNormals(s.pointNormals());
136  tmp<vectorField> tpointNormals(calcVertexNormals(s));
137  vectorField& pointNormals = tpointNormals();
138 
139 
140  // feature edges: create edge normals from edgeFaces only.
141  {
142  const labelListList& edgeFaces = s.edgeFaces();
143 
144  labelList nNormals(s.nPoints(), 0);
145  forAll(edgeStat, edgeI)
146  {
147  if (edgeStat[edgeI] != surfaceFeatures::NONE)
148  {
149  const edge& e = s.edges()[edgeI];
150  forAll(e, i)
151  {
152  if (!isFeaturePoint[e[i]])
153  {
154  pointNormals[e[i]] = vector::zero;
155  }
156  }
157  }
158  }
159 
160  forAll(edgeStat, edgeI)
161  {
162  if (edgeStat[edgeI] != surfaceFeatures::NONE)
163  {
164  const labelList& eFaces = edgeFaces[edgeI];
165 
166  // Get average edge normal
168  forAll(eFaces, i)
169  {
170  n += s.faceNormals()[eFaces[i]];
171  }
172  n /= eFaces.size();
173 
174 
175  // Sum to points
176  const edge& e = s.edges()[edgeI];
177  forAll(e, i)
178  {
179  if (!isFeaturePoint[e[i]])
180  {
181  pointNormals[e[i]] += n;
182  nNormals[e[i]]++;
183  }
184  }
185  }
186  }
187 
188  forAll(nNormals, pointI)
189  {
190  if (nNormals[pointI] > 0)
191  {
192  pointNormals[pointI] /= mag(pointNormals[pointI]);
193  }
194  }
195  }
196 
197 
198  forAll(pointNormals, pointI)
199  {
200  if (mag(mag(pointNormals[pointI])-1) > SMALL)
201  {
203  << "unitialised"
204  << exit(FatalError);
205  }
206  }
207 
208  return tpointNormals;
209 }
210 
211 
212 void detectSelfIntersections
213 (
214  const triSurfaceMesh& s,
215  PackedBoolList& isEdgeIntersecting
216 )
217 {
218  const edgeList& edges = s.edges();
219  const indexedOctree<treeDataTriSurface>& tree = s.tree();
220  const labelList& meshPoints = s.meshPoints();
221  const pointField& points = s.points();
222 
223  isEdgeIntersecting.setSize(edges.size());
224  isEdgeIntersecting = false;
225 
226  forAll(edges, edgeI)
227  {
228  const edge& e = edges[edgeI];
229 
230  pointIndexHit hitInfo
231  (
232  tree.findLine
233  (
234  points[meshPoints[e[0]]],
235  points[meshPoints[e[1]]],
237  )
238  );
239 
240  if (hitInfo.hit())
241  {
242  isEdgeIntersecting[edgeI] = true;
243  }
244  }
245 }
246 
247 
248 label detectIntersectionPoints
249 (
250  const scalar tolerance,
251  const triSurfaceMesh& s,
252 
253  const scalar extendFactor,
254  const pointField& initialPoints,
255  const vectorField& displacement,
256 
257  const bool checkSelfIntersect,
258  const PackedBoolList& initialIsEdgeIntersecting,
259 
260  PackedBoolList& isPointOnHitEdge,
261  scalarField& scale
262 )
263 {
264  const pointField initialLocalPoints(initialPoints, s.meshPoints());
265  const List<labelledTri>& localFaces = s.localFaces();
266 
267 
268  label nHits = 0;
269  isPointOnHitEdge.setSize(s.nPoints());
270  isPointOnHitEdge = false;
271 
272 
273  // 1. Extrusion offset vectors intersecting new surface location
274  {
275  scalar tol = max(tolerance, 10*s.tolerance());
276 
277  // Collect all the edge vectors. Slightly shorten the edges to prevent
278  // finding lots of intersections. The fast triangle intersection routine
279  // has problems with rays passing through a point of the
280  // triangle.
281  pointField start(initialLocalPoints+tol*displacement);
282  pointField end(initialLocalPoints+extendFactor*displacement);
283 
284  List<pointIndexHit> hits;
285  s.findLineAny(start, end, hits);
286 
287  forAll(hits, pointI)
288  {
289  if
290  (
291  hits[pointI].hit()
292  && findIndex(localFaces[hits[pointI].index()], pointI) == -1
293  )
294  {
295  scale[pointI] = max(0.0, scale[pointI]-0.2);
296 
297  isPointOnHitEdge[pointI] = true;
298  nHits++;
299  }
300  }
301  }
302 
303 
304  // 2. (new) surface self intersections
305  if (checkSelfIntersect)
306  {
307  PackedBoolList isEdgeIntersecting;
308  detectSelfIntersections(s, isEdgeIntersecting);
309 
310  const edgeList& edges = s.edges();
311  const pointField& points = s.points();
312 
313  forAll(edges, edgeI)
314  {
315  const edge& e = edges[edgeI];
316 
317  if (isEdgeIntersecting[edgeI] && !initialIsEdgeIntersecting[edgeI])
318  {
319  if (isPointOnHitEdge.set(e[0]))
320  {
321  label start = s.meshPoints()[e[0]];
322  const point& pt = points[start];
323 
324  Pout<< "Additional self intersection at "
325  << pt
326  << endl;
327 
328  scale[e[0]] = max(0.0, scale[e[0]]-0.2);
329  nHits++;
330  }
331  if (isPointOnHitEdge.set(e[1]))
332  {
333  label end = s.meshPoints()[e[1]];
334  const point& pt = points[end];
335 
336  Pout<< "Additional self intersection at "
337  << pt
338  << endl;
339 
340  scale[e[1]] = max(0.0, scale[e[1]]-0.2);
341  nHits++;
342  }
343  }
344  }
345  }
346 
347  return nHits;
348 }
349 
350 
352 (
353  const triSurface& s,
354  const scalarField& fld,
355  const scalarField& edgeWeights
356 )
357 {
358  tmp<scalarField> tres(new scalarField(s.nPoints(), 0.0));
359  scalarField& res = tres();
360 
361  scalarField sumWeight(s.nPoints(), 0.0);
362 
363  const edgeList& edges = s.edges();
364 
365  forAll(edges, edgeI)
366  {
367  const edge& e = edges[edgeI];
368  const scalar w = edgeWeights[edgeI];
369 
370  res[e[0]] += w*fld[e[1]];
371  sumWeight[e[0]] += w;
372 
373  res[e[1]] += w*fld[e[0]];
374  sumWeight[e[1]] += w;
375  }
376 
377  // Average
378  // ~~~~~~~
379 
380  forAll(res, pointI)
381  {
382  if (mag(sumWeight[pointI]) < VSMALL)
383  {
384  // Unconnected point. Take over original value
385  res[pointI] = fld[pointI];
386  }
387  else
388  {
389  res[pointI] /= sumWeight[pointI];
390  }
391  }
392 
393  return tres;
394 }
395 
396 
397 void minSmooth
398 (
399  const triSurface& s,
400  const PackedBoolList& isAffectedPoint,
401  const scalarField& fld,
402  scalarField& newFld
403 )
404 {
405 
406  const edgeList& edges = s.edges();
407  const pointField& points = s.points();
408  const labelList& mp = s.meshPoints();
409 
410  scalarField edgeWeights(edges.size());
411  forAll(edges, edgeI)
412  {
413  const edge& e = edges[edgeI];
414  scalar w = mag(points[mp[e[0]]]-points[mp[e[1]]]);
415 
416  edgeWeights[edgeI] = 1.0/(max(w, SMALL));
417  }
418 
419  tmp<scalarField> tavgFld = avg(s, fld, edgeWeights);
420 
421  const scalarField& avgFld = tavgFld();
422 
423  forAll(fld, pointI)
424  {
425  if (isAffectedPoint.get(pointI) == 1)
426  {
427  newFld[pointI] = min
428  (
429  fld[pointI],
430  0.5*fld[pointI] + 0.5*avgFld[pointI]
431  //avgFld[pointI]
432  );
433  }
434  }
435 }
436 
437 
438 void lloydsSmoothing
439 (
440  const label nSmooth,
441  triSurface& s,
442  const PackedBoolList& isFeaturePoint,
443  const List<surfaceFeatures::edgeStatus>& edgeStat,
444  const PackedBoolList& isAffectedPoint
445 )
446 {
447  const labelList& meshPoints = s.meshPoints();
448  const edgeList& edges = s.edges();
449 
450 
451  PackedBoolList isSmoothPoint(isAffectedPoint);
452  // Extend isSmoothPoint
453  {
454  PackedBoolList newIsSmoothPoint(isSmoothPoint);
455  forAll(edges, edgeI)
456  {
457  const edge& e = edges[edgeI];
458  if (isSmoothPoint[e[0]])
459  {
460  newIsSmoothPoint[e[1]] = true;
461  }
462  if (isSmoothPoint[e[1]])
463  {
464  newIsSmoothPoint[e[0]] = true;
465  }
466  }
467  Info<< "From points-to-smooth " << isSmoothPoint.count()
468  << " to " << newIsSmoothPoint.count()
469  << endl;
470  isSmoothPoint.transfer(newIsSmoothPoint);
471  }
472 
473  // Do some smoothing (Lloyds algorithm) around problematic points
474  for (label i = 0; i < nSmooth; i++)
475  {
476  const labelListList& pointFaces = s.pointFaces();
477  const pointField& faceCentres = s.faceCentres();
478 
479  pointField newPoints(s.points());
480  forAll(isSmoothPoint, pointI)
481  {
482  if (isSmoothPoint[pointI] && !isFeaturePoint[pointI])
483  {
484  const labelList& pFaces = pointFaces[pointI];
485 
486  point avg = point::zero;
487  forAll(pFaces, pFaceI)
488  {
489  avg += faceCentres[pFaces[pFaceI]];
490  }
491  avg /= pFaces.size();
492  newPoints[meshPoints[pointI]] = avg;
493  }
494  }
495 
496 
497  // Move points on feature edges only according to feature edges.
498 
499  const pointField& points = s.points();
500 
501  vectorField pointSum(s.nPoints(), vector::zero);
502  labelList nPointSum(s.nPoints(), 0);
503 
504  forAll(edges, edgeI)
505  {
506  if (edgeStat[edgeI] != surfaceFeatures::NONE)
507  {
508  const edge& e = edges[edgeI];
509  const point& start = points[meshPoints[e[0]]];
510  const point& end = points[meshPoints[e[1]]];
511  point eMid = 0.5*(start+end);
512  pointSum[e[0]] += eMid;
513  nPointSum[e[0]]++;
514  pointSum[e[1]] += eMid;
515  nPointSum[e[1]]++;
516  }
517  }
518 
519  forAll(pointSum, pointI)
520  {
521  if
522  (
523  isSmoothPoint[pointI]
524  && isFeaturePoint[pointI]
525  && nPointSum[pointI] >= 2
526  )
527  {
528  newPoints[meshPoints[pointI]] =
529  pointSum[pointI]/nPointSum[pointI];
530  }
531  }
532 
533 
534  s.movePoints(newPoints);
535 
536 
537  // Extend isSmoothPoint
538  {
539  PackedBoolList newIsSmoothPoint(isSmoothPoint);
540  forAll(edges, edgeI)
541  {
542  const edge& e = edges[edgeI];
543  if (isSmoothPoint[e[0]])
544  {
545  newIsSmoothPoint[e[1]] = true;
546  }
547  if (isSmoothPoint[e[1]])
548  {
549  newIsSmoothPoint[e[0]] = true;
550  }
551  }
552  Info<< "From points-to-smooth " << isSmoothPoint.count()
553  << " to " << newIsSmoothPoint.count()
554  << endl;
555  isSmoothPoint.transfer(newIsSmoothPoint);
556  }
557  }
558 }
559 
560 
561 
562 // Main program:
563 
564 int main(int argc, char *argv[])
565 {
566  argList::addNote("Inflates surface according to point normals.");
567 
570  (
571  "Creates inflated version of surface using point normals."
572  " Takes surface, distance to inflate and additional safety factor"
573  );
575  (
576  "checkSelfIntersection",
577  "also check for self-intersection"
578  );
580  (
581  "nSmooth",
582  "integer",
583  "number of smoothing iterations (default 20)"
584  );
586  (
587  "featureAngle",
588  "scalar",
589  "feature angle"
590  );
592  (
593  "debug",
594  "switch on additional debug information"
595  );
596 
597  argList::validArgs.append("inputFile");
598  argList::validArgs.append("distance");
599  argList::validArgs.append("safety factor [1..]");
600 
601  #include "setRootCase.H"
602  #include "createTime.H"
603  runTime.functionObjects().off();
604 
605  const word inputName(args.args()[1]);
606  const scalar distance(args.argRead<scalar>(2));
607  const scalar extendFactor(args.argRead<scalar>(3));
608  const bool checkSelfIntersect = args.optionFound("checkSelfIntersection");
609  const label nSmooth = args.optionLookupOrDefault("nSmooth", 10);
610  const scalar featureAngle = args.optionLookupOrDefault<scalar>
611  (
612  "featureAngle",
613  180
614  );
615  const bool debug = args.optionFound("debug");
616 
617 
618  Info<< "Inflating surface " << inputName
619  << " according to point normals" << nl
620  << " distance : " << distance << nl
621  << " safety factor : " << extendFactor << nl
622  << " self intersection test : " << checkSelfIntersect << nl
623  << " smoothing iterations : " << nSmooth << nl
624  << " feature angle : " << featureAngle << nl
625  << endl;
626 
627 
628  if (extendFactor < 1 || extendFactor > 10)
629  {
631  << "Illegal safety factor " << extendFactor
632  << ". It is usually 1..2"
633  << exit(FatalError);
634  }
635 
636 
637 
638  // Load triSurfaceMesh
640  (
641  IOobject
642  (
643  inputName, // name
644  runTime.constant(), // instance
645  "triSurface", // local
646  runTime, // registry
649  )
650  );
651 
652  // Mark features
653  const surfaceFeatures features(s, 180.0-featureAngle);
654 
655  Info<< "Detected features:" << nl
656  << " feature points : " << features.featurePoints().size()
657  << " out of " << s.nPoints() << nl
658  << " feature edges : " << features.featureEdges().size()
659  << " out of " << s.nEdges() << nl
660  << endl;
661 
662  PackedBoolList isFeaturePoint(s.nPoints());
663  forAll(features.featurePoints(), i)
664  {
665  label pointI = features.featurePoints()[i];
666  isFeaturePoint[pointI] = true;
667  }
668  const List<surfaceFeatures::edgeStatus> edgeStat(features.toStatus());
669 
670 
671 
672 
673  const pointField initialPoints(s.points());
674 
675 
676  // Construct scale
677  Info<< "Constructing field scale\n" << endl;
679  (
680  IOobject
681  (
682  "scale", // name
683  runTime.timeName(), // instance
684  s, // registry
687  ),
688  s,
689  dimensionedScalar("scale", dimLength, 1.0)
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  distance*scale*pointNormals
722  );
723 
724 
725  const labelList& meshPoints = s.meshPoints();
726 
727 
728  // Any point on any intersected edge in any of the iterations
729  PackedBoolList isScaledPoint(s.nPoints());
730 
731 
732  // Detect any self intersections on initial mesh
733  PackedBoolList initialIsEdgeIntersecting;
734  if (checkSelfIntersect)
735  {
736  detectSelfIntersections(s, initialIsEdgeIntersecting);
737  }
738  else
739  {
740  // Mark all edges as already self intersecting so avoid detecting any
741  // new ones
742  initialIsEdgeIntersecting.setSize(s.nEdges(), true);
743  }
744 
745 
746  // Inflate
747  while (runTime.loop())
748  {
749  Info<< "Time = " << runTime.timeName() << nl << endl;
750 
751  // Move to new location
752  pointField newPoints(initialPoints);
753  forAll(meshPoints, i)
754  {
755  newPoints[meshPoints[i]] += pointDisplacement[i];
756  }
757  s.movePoints(newPoints);
758  Info<< "Bounding box : " << s.bounds() << endl;
759 
760 
761  // Work out scale from pointDisplacement
762  forAll(scale, pointI)
763  {
764  if (s.pointFaces()[pointI].size())
765  {
766  scale[pointI] = mag(pointDisplacement[pointI])/distance;
767  }
768  else
769  {
770  scale[pointI] = 1.0;
771  }
772  }
773 
774 
775  Info<< "Scale : " << gAverage(scale) << endl;
776  Info<< "Displacement : " << gAverage(pointDisplacement) << endl;
777 
778 
779 
780  // Detect any intersections and scale back
781  PackedBoolList isAffectedPoint;
782  label nIntersections = detectIntersectionPoints
783  (
784  1e-9, // intersection tolerance
785  s,
786  extendFactor,
787  initialPoints,
788  pointDisplacement,
789 
790  checkSelfIntersect,
791  initialIsEdgeIntersecting,
792 
793  isAffectedPoint,
794  scale
795  );
796  Info<< "Detected " << nIntersections << " intersections" << nl
797  << endl;
798 
799  if (nIntersections == 0)
800  {
801  runTime.write();
802  break;
803  }
804 
805 
806  // Accumulate all affected points
807  forAll(isAffectedPoint, pointI)
808  {
809  if (isAffectedPoint[pointI])
810  {
811  isScaledPoint[pointI] = true;
812  }
813  }
814 
815  // Smear out lowering of scale so any edges not found are
816  // still included
817  for (label i = 0; i < nSmooth; i++)
818  {
819  triSurfacePointScalarField oldScale(scale);
820  oldScale.rename("oldScale");
821  minSmooth
822  (
823  s,
824  PackedBoolList(s.nPoints(), true),
825  oldScale,
826  scale
827  );
828  }
829 
830 
831  // From scale update the pointDisplacement
832  pointDisplacement *= distance*scale/mag(pointDisplacement);
833 
834 
835  // Do some smoothing (Lloyds algorithm)
836  lloydsSmoothing(nSmooth, s, isFeaturePoint, edgeStat, isAffectedPoint);
837 
838 
839  // Update pointDisplacement
840  const pointField& pts = s.points();
841  forAll(meshPoints, i)
842  {
843  label meshPointI = meshPoints[i];
844  pointDisplacement[i] = pts[meshPointI]-initialPoints[meshPointI];
845  }
846 
847 
848  // Write
849  runTime.write();
850 
851  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
852  << " ClockTime = " << runTime.elapsedClockTime() << " s"
853  << nl << endl;
854  }
855 
856 
857  Info<< "Detected overall intersecting points " << isScaledPoint.count()
858  << " out of " << s.nPoints() << nl << endl;
859 
860 
861  if (debug)
862  {
863  OBJstream str(runTime.path()/"isScaledPoint.obj");
864  forAll(isScaledPoint, pointI)
865  {
866  if (isScaledPoint[pointI])
867  {
868  str.write(initialPoints[meshPoints[pointI]]);
869  }
870  }
871  }
872 
873 
874  Info << "End\n" << endl;
875 
876  return 0;
877 }
878 
879 
880 // ************************************************************************* //
Foam::argList::validArgs
static SLList< string > validArgs
A list of valid (mandatory) arguments.
Definition: argList.H:143
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:48
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::constant::atomic::mp
const dimensionedScalar mp
Proton mass.
Foam::PackedBoolList
A bit-packed bool list.
Definition: PackedBoolList.H:63
Foam::PrimitivePatch::pointFaces
const labelListList & pointFaces() const
Return point-face addressing.
Definition: PrimitivePatchTemplate.C:352
w
volScalarField w(IOobject("w", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE), mesh, dimensionedScalar("w", dimensionSet(0, 0, 0, 0, 0, 0, 0), 0.0))
Foam::PrimitivePatch::points
const Field< PointType > & points() const
Return reference to global points.
Definition: PrimitivePatchTemplate.H:282
Foam::Vector< scalar >::zero
static const Vector zero
Definition: Vector.H:80
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::argList::args
const stringList & args() const
Return arguments.
Definition: argListI.H:66
Foam::OBJstream
OFstream which keeps track of vertices.
Definition: OBJstream.H:53
Foam::IOobject::AUTO_WRITE
@ AUTO_WRITE
Definition: IOobject.H:117
Foam::argList::addOption
static void addOption(const word &opt, const string &param="", const string &usage="")
Add to an option to validOptions with usage information.
Definition: argList.C:108
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
Foam::argList::addNote
static void addNote(const string &)
Add extra notes for the usage information.
Definition: argList.C:139
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::findIndex
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurence of given element and return index,.
Foam::gAverage
Type gAverage(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:571
surfaceFeatures.H
Foam::PackedBoolList::set
void set(const PackedList< 1 > &)
Set specified bits.
Definition: PackedBoolList.C:169
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
Foam::argList::addBoolOption
static void addBoolOption(const word &opt, const string &usage="")
Add to a bool option to validOptions with usage information.
Definition: argList.C:98
triSurfaceFields.H
Fields for triSurface.
Foam::IOobject::MUST_READ
@ MUST_READ
Definition: IOobject.H:108
Foam::triSurfaceMesh
IOoject and searching on triSurface.
Definition: triSurfaceMesh.H:63
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
triSurfaceGeoMesh.H
calcVertexNormals
vectorField calcVertexNormals(const triSurface &surf)
Definition: surfaceFeatureExtract.C:189
Foam::PackedList::setSize
void setSize(const label, const unsigned int &val=0u)
Alias for resize()
Definition: PackedListI.H:823
pFaces
Info<< "Finished reading KIVA file"<< endl;cellShapeList cellShapes(nPoints);labelList cellZoning(nPoints, -1);const cellModel &hex=*(cellModeller::lookup("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]=cellShape(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
calcVertexNormalWeight
scalar calcVertexNormalWeight(const triFace &f, const label pI, const vector &fN, const pointField &points)
Definition: surfaceFeatureExtract.C:62
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::surfaceFeatures::NONE
@ NONE
Definition: surfaceFeatures.H:74
Foam::PointIndexHit
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:53
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::triSurface
Triangulated surface description with patch information.
Definition: triSurface.H:57
Foam::indexedOctree::findLine
pointIndexHit findLine(const bool findAny, const point &treeStart, const point &treeEnd, const label startNodeI, const direction startOctantI, const FindIntersectOp &fiOp, const bool verbose=false) const
Find any or nearest intersection.
Foam::Ostream::write
virtual Ostream & write(const token &)=0
Write next token to stream.
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
PackedBoolList.H
argList.H
Foam::surfaceFeatures
Holds feature edges/points of surface.
Definition: surfaceFeatures.H:68
Foam::indexedOctree
Non-pointer based hierarchical recursive searching.
Definition: treeDataTriSurface.H:47
Foam::treeDataPrimitivePatch::findSelfIntersectOp
Definition: treeDataPrimitivePatch.H:171
main
int main(int argc, char *argv[])
Definition: postCalc.C:54
Foam::argList::optionLookupOrDefault
T optionLookupOrDefault(const word &opt, const T &deflt) const
Read a value from the named option if present.
Definition: argListI.H:237
Foam::PrimitivePatch::nPoints
label nPoints() const
Return number of points supporting patch faces.
Definition: PrimitivePatchTemplate.H:293
Foam::FatalError
error FatalError
fld
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))
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
Foam::PackedList::get
unsigned int get(const label) const
Get value at index I.
Definition: PackedListI.H:964
Foam::e
const double e
Elementary charge.
Definition: doubleFloat.H:94
s
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Foam::distance
scalar distance(const vector &p1, const vector &p2)
Definition: curveTools.C:12
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
setRootCase.H
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::triFace
A triangular face using a FixedList of labels corresponding to mesh vertices.
Definition: triFace.H:68
Foam::Pout
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
f
labelList f(nPoints)
Foam::Vector< scalar >
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
points
const pointField & points
Definition: gmvOutputHeader.H:1
createTime.H
Foam::argList::optionFound
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:108
Foam::argList::argRead
T argRead(const label index) const
Read a value from the argument at index.
Definition: argListI.H:177
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::argList::noParallel
static void noParallel()
Remove the parallel options.
Definition: argList.C:161
Foam::PrimitivePatch::meshPoints
const labelList & meshPoints() const
Return labelList of mesh points in patch.
Definition: PrimitivePatchTemplate.C:392
args
Foam::argList args(argc, argv)
Foam::IOobject::READ_IF_PRESENT
@ READ_IF_PRESENT
Definition: IOobject.H:110
triSurfaceMesh.H
Foam::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::magSqr
dimensioned< scalar > magSqr(const dimensioned< Type > &)
OBJstream.H
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:51