medialAxisMeshMover.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) 2014-2015 OpenFOAM Foundation
6  \\/ M anipulation | Copyright (C) 2015 OpenCFD Ltd.
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "medialAxisMeshMover.H"
28 #include "pointFields.H"
29 #include "valuePointPatchFields.H"
30 #include "PointEdgeWave.H"
31 #include "meshRefinement.H"
32 #include "unitConversion.H"
33 #include "PatchTools.H"
34 #include "OBJstream.H"
35 #include "PointData.H"
37 
38 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39 
40 namespace Foam
41 {
42  defineTypeNameAndDebug(medialAxisMeshMover, 0);
43 
45  (
46  externalDisplacementMeshMover,
47  medialAxisMeshMover,
48  dictionary
49  );
50 }
51 
52 
53 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
54 
55 // Tries and find a medial axis point. Done by comparing vectors to nearest
56 // wall point for both vertices of edge.
58 (
59  const List<PointData<vector> >& pointWallDist,
60  const label edgeI,
61  const scalar minCos
62 ) const
63 {
64  const pointField& points = mesh().points();
65 
66  // Do not mark edges with one side on moving wall.
67 
68  const edge& e = mesh().edges()[edgeI];
69 
70  vector v0(points[e[0]] - pointWallDist[e[0]].origin());
71  scalar magV0(mag(v0));
72 
73  if (magV0 < SMALL)
74  {
75  return false;
76  }
77 
78  vector v1(points[e[1]] - pointWallDist[e[1]].origin());
79  scalar magV1(mag(v1));
80 
81  if (magV1 < SMALL)
82  {
83  return false;
84  }
85 
86 
87  //- Detect based on vector to nearest point differing for both endpoints
88  //v0 /= magV0;
89  //v1 /= magV1;
90  //
92  //if ((v0 & v1) < minCos)
93  //{
94  // return true;
95  //}
96  //else
97  //{
98  // return false;
99  //}
100 
101  //- Detect based on extrusion vector differing for both endpoints
102  // the idea is that e.g. a sawtooth wall can still be extruded
103  // successfully as long as it is done all to the same direction.
104  if ((pointWallDist[e[0]].data() & pointWallDist[e[1]].data()) < minCos)
105  {
106  return true;
107  }
108  else
109  {
110  return false;
111  }
112 }
113 
114 
116 {
117  Info<< typeName
118  << " : Calculating distance to Medial Axis ..." << endl;
119 
120  const pointField& points = mesh().points();
121 
123  const labelList& meshPoints = pp.meshPoints();
124 
125 
126  // Read a few parameters
127  // ~~~~~~~~~~~~~~~~~~~~~
128 
129  //- Smooth surface normals
130  const label nSmoothSurfaceNormals = readLabel
131  (
132  coeffDict.lookup("nSmoothSurfaceNormals")
133  );
134 
135  //- When is medial axis
136  word angleKey = "minMedialAxisAngle";
137  if (!coeffDict.found(angleKey))
138  {
139  // Backwards compatibility
140  angleKey = "minMedianAxisAngle";
141  }
142  scalar minMedialAxisAngleCos = Foam::cos
143  (
144  degToRad(readScalar(coeffDict.lookup(angleKey)))
145  );
146 
147  //- Feature angle when to stop adding layers
148  const scalar featureAngle = readScalar(coeffDict.lookup("featureAngle"));
149 
150  //- When to slip along wall
151  const scalar slipFeatureAngle =
152  (
153  coeffDict.found("slipFeatureAngle")
154  ? readScalar(coeffDict.lookup("slipFeatureAngle"))
155  : 0.5*featureAngle
156  );
157 
158  //- Smooth internal normals
159  const label nSmoothNormals = readLabel
160  (
161  coeffDict.lookup("nSmoothNormals")
162  );
163 
164  //- Number of edges walking out
165  const label nMedialAxisIter = coeffDict.lookupOrDefault<label>
166  (
167  "nMedialAxisIter",
169  );
170 
171 
172  // Predetermine mesh edges
173  // ~~~~~~~~~~~~~~~~~~~~~~~
174 
175  // Precalulate (mesh) master point/edge (only relevant for shared pts/edges)
176  const PackedBoolList isMeshMasterPoint(syncTools::getMasterPoints(mesh()));
177  const PackedBoolList isMeshMasterEdge(syncTools::getMasterEdges(mesh()));
178  // Precalculate meshEdge per pp edge
179  const labelList meshEdges
180  (
181  pp.meshEdges
182  (
183  mesh().edges(),
184  mesh().pointEdges()
185  )
186  );
187 
188  // Precalulate (patch) master point/edge
189  const PackedBoolList isPatchMasterPoint
190  (
192  (
193  mesh(),
194  meshPoints
195  )
196  );
197  const PackedBoolList isPatchMasterEdge
198  (
200  (
201  mesh(),
202  meshEdges
203  )
204  );
205 
206  // Determine pointNormal
207  // ~~~~~~~~~~~~~~~~~~~~~
208 
209  pointField pointNormals(PatchTools::pointNormals(mesh(), pp));
210 
211  // Smooth patch normal vectors
213  (
214  nSmoothSurfaceNormals,
215  isPatchMasterPoint,
216  isPatchMasterEdge,
217  pp,
218  pointNormals
219  );
220 
221 
222  // Calculate distance to pp points
223  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
224 
225  // Distance to wall
226  List<PointData<vector> > pointWallDist(mesh().nPoints());
227 
228  // Dummy additional info for PointEdgeWave
229  int dummyTrackData = 0;
230 
231 
232  // 1. Calculate distance to points where displacement is specified.
233  {
234  // Seed data.
235  List<PointData<vector> > wallInfo(meshPoints.size());
236 
237  forAll(meshPoints, patchPointI)
238  {
239  label pointI = meshPoints[patchPointI];
240  wallInfo[patchPointI] = PointData<vector>
241  (
242  points[pointI],
243  0.0,
244  pointNormals[patchPointI] // surface normals
245  );
246  }
247 
248  // Do all calculations
249  List<PointData<vector> > edgeWallDist(mesh().nEdges());
250  PointEdgeWave<PointData<vector> > wallDistCalc
251  (
252  mesh(),
253  meshPoints,
254  wallInfo,
255  pointWallDist,
256  edgeWallDist,
257  0, // max iterations
258  dummyTrackData
259  );
260  wallDistCalc.iterate(nMedialAxisIter);
261 
262  label nUnvisit = returnReduce
263  (
264  wallDistCalc.getUnsetPoints(),
265  sumOp<label>()
266  );
267 
268  if (nUnvisit > 0)
269  {
270  if (nMedialAxisIter > 0)
271  {
272  Info<< typeName
273  << " : Limited walk to " << nMedialAxisIter
274  << " steps. Not visited " << nUnvisit
275  << " out of " << mesh().globalData().nTotalPoints()
276  << " points" << endl;
277  }
278  else
279  {
281  << "Walking did not visit all points." << nl
282  << " Did not visit " << nUnvisit
283  << " out of " << mesh().globalData().nTotalPoints()
284  << " points. This is not necessarily a problem" << nl
285  << " and might be due to faceZones splitting of part"
286  << " of the domain." << nl << endl;
287  }
288  }
289  }
290 
291 
292  // 2. Find points with max distance and transport information back to
293  // wall.
294  {
295  List<pointEdgePoint> pointMedialDist(mesh().nPoints());
296  List<pointEdgePoint> edgeMedialDist(mesh().nEdges());
297 
298  // Seed point data.
299  DynamicList<pointEdgePoint> maxInfo(meshPoints.size());
300  DynamicList<label> maxPoints(meshPoints.size());
301 
302  // 1. Medial axis points
303 
304  const edgeList& edges = mesh().edges();
305 
306  forAll(edges, edgeI)
307  {
308  const edge& e = edges[edgeI];
309 
310  if
311  (
312  !pointWallDist[e[0]].valid(dummyTrackData)
313  || !pointWallDist[e[1]].valid(dummyTrackData)
314  )
315  {
316  // Unvisited point. See above about nUnvisit warning
317  forAll(e, ep)
318  {
319  label pointI = e[ep];
320 
321  if (!pointMedialDist[pointI].valid(dummyTrackData))
322  {
323  maxPoints.append(pointI);
324  maxInfo.append
325  (
327  (
328  points[pointI],
329  0.0
330  )
331  );
332  pointMedialDist[pointI] = maxInfo.last();
333  }
334  }
335 
336  }
337  else if (isMaxEdge(pointWallDist, edgeI, minMedialAxisAngleCos))
338  {
339  // Both end points of edge have very different nearest wall
340  // point. Mark both points as medial axis points.
341 
342  // Approximate medial axis location on edge.
343  //const point medialAxisPt = e.centre(points);
344  vector eVec = e.vec(points);
345  scalar eMag = mag(eVec);
346  if (eMag > VSMALL)
347  {
348  eVec /= eMag;
349 
350  // Calculate distance along edge
351  const point& p0 = points[e[0]];
352  const point& p1 = points[e[1]];
353  scalar dist0 = (p0-pointWallDist[e[0]].origin()) & eVec;
354  scalar dist1 = (pointWallDist[e[1]].origin()-p1) & eVec;
355  scalar s = 0.5*(dist1+eMag+dist0);
356 
357  point medialAxisPt;
358  if (s <= dist0)
359  {
360  medialAxisPt = p0;
361  }
362  else if (s >= dist0+eMag)
363  {
364  medialAxisPt = p1;
365  }
366  else
367  {
368  medialAxisPt = p0+(s-dist0)*eVec;
369  }
370 
371  forAll(e, ep)
372  {
373  label pointI = e[ep];
374 
375  if (!pointMedialDist[pointI].valid(dummyTrackData))
376  {
377  maxPoints.append(pointI);
378  maxInfo.append
379  (
381  (
382  medialAxisPt, //points[pointI],
383  magSqr(points[pointI]-medialAxisPt)//0.0
384  )
385  );
386  pointMedialDist[pointI] = maxInfo.last();
387  }
388  }
389  }
390  }
391  }
392 
393 
394  // 2. Seed non-adapt patches
396 
397  labelHashSet adaptPatches(adaptPatchIDs_);
398 
399 
400  forAll(patches, patchI)
401  {
402  const polyPatch& pp = patches[patchI];
403  const pointPatchVectorField& pvf =
404  pointDisplacement().boundaryField()[patchI];
405 
406  if
407  (
408  !pp.coupled()
409  && !isA<emptyPolyPatch>(pp)
410  && !adaptPatches.found(patchI)
411  )
412  {
413  const labelList& meshPoints = pp.meshPoints();
414 
415  // Check the type of the patchField. The types are
416  // - fixedValue (0 or more layers) but the >0 layers have
417  // already been handled in the adaptPatches loop
418  // - constraint (but not coupled) types, e.g. symmetryPlane,
419  // slip.
420  if (pvf.fixesValue())
421  {
422  // Disable all movement on fixedValue patchFields
423  Info<< typeName
424  << " : Inserting all points on patch " << pp.name()
425  << endl;
426 
427  forAll(meshPoints, i)
428  {
429  label pointI = meshPoints[i];
430  if (!pointMedialDist[pointI].valid(dummyTrackData))
431  {
432  maxPoints.append(pointI);
433  maxInfo.append
434  (
436  (
437  points[pointI],
438  0.0
439  )
440  );
441  pointMedialDist[pointI] = maxInfo.last();
442  }
443  }
444  }
445  else
446  {
447  // Based on geometry: analyse angle w.r.t. nearest moving
448  // point. In the pointWallDist we transported the
449  // normal as the passive vector. Note that this points
450  // out of the originating wall so inside of the domain
451  // on this patch.
452  Info<< typeName
453  << " : Inserting points on patch " << pp.name()
454  << " if angle to nearest layer patch > "
455  << slipFeatureAngle << " degrees." << endl;
456 
457  scalar slipFeatureAngleCos = Foam::cos
458  (
459  degToRad(slipFeatureAngle)
460  );
461  pointField pointNormals
462  (
464  );
465 
466  forAll(meshPoints, i)
467  {
468  label pointI = meshPoints[i];
469 
470  if
471  (
472  pointWallDist[pointI].valid(dummyTrackData)
473  && !pointMedialDist[pointI].valid(dummyTrackData)
474  )
475  {
476  // Check if angle not too large.
477  scalar cosAngle =
478  (
479  -pointWallDist[pointI].data()
480  & pointNormals[i]
481  );
482  if (cosAngle > slipFeatureAngleCos)
483  {
484  // Extrusion direction practically perpendicular
485  // to the patch. Disable movement at the patch.
486 
487  maxPoints.append(pointI);
488  maxInfo.append
489  (
491  (
492  points[pointI],
493  0.0
494  )
495  );
496  pointMedialDist[pointI] = maxInfo.last();
497  }
498  else
499  {
500  // Extrusion direction makes angle with patch
501  // so allow sliding - don't insert zero points
502  }
503  }
504  }
505  }
506  }
507  }
508 
509  maxInfo.shrink();
510  maxPoints.shrink();
511 
512  // Do all calculations
513  PointEdgeWave<pointEdgePoint> medialDistCalc
514  (
515  mesh(),
516  maxPoints,
517  maxInfo,
518 
519  pointMedialDist,
520  edgeMedialDist,
521  0, // max iterations
522  dummyTrackData
523  );
524  medialDistCalc.iterate(2*nMedialAxisIter);
525 
526 
527  // Extract medial axis distance as pointScalarField
528  forAll(pointMedialDist, pointI)
529  {
530  if (pointMedialDist[pointI].valid(dummyTrackData))
531  {
532  medialDist_[pointI] = Foam::sqrt
533  (
534  pointMedialDist[pointI].distSqr()
535  );
536  medialVec_[pointI] = pointMedialDist[pointI].origin();
537  }
538  else
539  {
540  // Unvisited. Do as if on medial axis so unmoving
541  medialDist_[pointI] = 0.0;
542  medialVec_[pointI] = point(1, 0, 0);
543  }
544  }
545  }
546 
547  // Extract transported surface normals as pointVectorField
548  forAll(dispVec_, i)
549  {
550  if (!pointWallDist[i].valid(dummyTrackData))
551  {
552  dispVec_[i] = vector(1, 0, 0);
553  }
554  else
555  {
556  dispVec_[i] = pointWallDist[i].data();
557  }
558  }
559 
560  // Smooth normal vectors. Do not change normals on pp.meshPoints
562  (
563  nSmoothNormals,
564  isMeshMasterPoint,
565  isMeshMasterEdge,
566  meshPoints,
567  dispVec_
568  );
569 
570  // Calculate ratio point medial distance to point wall distance
571  forAll(medialRatio_, pointI)
572  {
573  if (!pointWallDist[pointI].valid(dummyTrackData))
574  {
575  medialRatio_[pointI] = 0.0;
576  }
577  else
578  {
579  scalar wDist2 = pointWallDist[pointI].distSqr();
580  scalar mDist = medialDist_[pointI];
581 
582  if (wDist2 < sqr(SMALL) && mDist < SMALL)
583  //- Note: maybe less strict:
584  //(
585  // wDist2 < sqr(meshRefiner_.mergeDistance())
586  // && mDist < meshRefiner_.mergeDistance()
587  //)
588  {
589  medialRatio_[pointI] = 0.0;
590  }
591  else
592  {
593  medialRatio_[pointI] = mDist / (Foam::sqrt(wDist2) + mDist);
594  }
595  }
596  }
597 
598 
599  if (debug)
600  {
601  Info<< typeName
602  << " : Writing medial axis fields:" << nl
603  << incrIndent
604  << "ratio of medial distance to wall distance : "
605  << medialRatio_.name() << nl
606  << "distance to nearest medial axis : "
607  << medialDist_.name() << nl
608  << "nearest medial axis location : "
609  << medialVec_.name() << nl
610  << "normal at nearest wall : "
611  << dispVec_.name() << nl
612  << decrIndent << nl
613  << endl;
614 
615  dispVec_.write();
616  medialRatio_.write();
617  medialDist_.write();
618  medialVec_.write();
619  }
620 }
621 
622 
624 (
625  const label patchPointI,
626  pointField& patchDisp,
628 )
629 {
630  if (extrudeStatus[patchPointI] == autoLayerDriver::EXTRUDE)
631  {
632  extrudeStatus[patchPointI] = autoLayerDriver::NOEXTRUDE;
633  patchDisp[patchPointI] = vector::zero;
634  return true;
635  }
636  else if (extrudeStatus[patchPointI] == autoLayerDriver::EXTRUDEREMOVE)
637  {
638  extrudeStatus[patchPointI] = autoLayerDriver::NOEXTRUDE;
639  patchDisp[patchPointI] = vector::zero;
640  return true;
641  }
642  else
643  {
644  return false;
645  }
646 }
647 
648 
650 (
651  const scalarField& minThickness,
652  pointField& patchDisp,
654 ) const
655 {
656  const indirectPrimitivePatch& pp = adaptPatchPtr_();
657  const labelList& meshPoints = pp.meshPoints();
658 
659  label nChangedTotal = 0;
660 
661  while (true)
662  {
663  label nChanged = 0;
664 
665  // Sync displacement (by taking min)
667  (
668  mesh(),
669  meshPoints,
670  patchDisp,
672  point::rootMax // null value
673  );
674 
675  // Unmark if displacement too small
676  forAll(patchDisp, i)
677  {
678  if (mag(patchDisp[i]) < minThickness[i])
679  {
680  if (unmarkExtrusion(i, patchDisp, extrudeStatus))
681  {
682  nChanged++;
683  }
684  }
685  }
686 
687  //labelList syncPatchNLayers(patchNLayers);
688  //
689  //syncTools::syncPointList
690  //(
691  // mesh(),
692  // meshPoints,
693  // syncPatchNLayers,
694  // minEqOp<label>(),
695  // labelMax // null value
696  //);
697  //
700  //forAll(syncPatchNLayers, i)
701  //{
702  // if (syncPatchNLayers[i] != patchNLayers[i])
703  // {
704  // if
705  // (
706  // unmarkExtrusion
707  // (
708  // i,
709  // patchDisp,
710  // patchNLayers,
711  // extrudeStatus
712  // )
713  // )
714  // {
715  // nChanged++;
716  // }
717  // }
718  //}
719  //
720  //syncTools::syncPointList
721  //(
722  // mesh(),
723  // meshPoints,
724  // syncPatchNLayers,
725  // maxEqOp<label>(),
726  // labelMin // null value
727  //);
728  //
731  //forAll(syncPatchNLayers, i)
732  //{
733  // if (syncPatchNLayers[i] != patchNLayers[i])
734  // {
735  // if
736  // (
737  // unmarkExtrusion
738  // (
739  // i,
740  // patchDisp,
741  // patchNLayers,
742  // extrudeStatus
743  // )
744  // )
745  // {
746  // nChanged++;
747  // }
748  // }
749  //}
750 
751  nChangedTotal += nChanged;
752 
753  if (!returnReduce(nChanged, sumOp<label>()))
754  {
755  break;
756  }
757  }
758 
759  //Info<< "Prevented extrusion on "
760  // << returnReduce(nChangedTotal, sumOp<label>())
761  // << " coupled patch points during syncPatchDisplacement." << endl;
762 }
763 
764 
765 // Stop layer growth where mesh wraps around edge with a
766 // large feature angle
769 (
770  const scalar minCos,
771  const PackedBoolList& isPatchMasterPoint,
772  const labelList& meshEdges,
773  List<autoLayerDriver::extrudeMode>& extrudeStatus,
774  pointField& patchDisp,
775  label& nPointCounter
776 ) const
777 {
778  const indirectPrimitivePatch& pp = adaptPatchPtr_();
779 
780  // Mark faces that have all points extruded
781  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
782 
783  boolList extrudedFaces(pp.size(), true);
784 
785  forAll(pp.localFaces(), faceI)
786  {
787  const face& f = pp.localFaces()[faceI];
788 
789  forAll(f, fp)
790  {
791  if (extrudeStatus[f[fp]] == autoLayerDriver::NOEXTRUDE)
792  {
793  extrudedFaces[faceI] = false;
794  break;
795  }
796  }
797  }
798 
799 
800 
801  //label nOldPointCounter = nPointCounter;
802 
803  // Detect situation where two featureedge-neighbouring faces are partly or
804  // not extruded and the edge itself is extruded. In this case unmark the
805  // edge for extrusion.
806 
807 
808  List<List<point> > edgeFaceNormals(pp.nEdges());
809  List<List<bool> > edgeFaceExtrude(pp.nEdges());
810 
811  const labelListList& edgeFaces = pp.edgeFaces();
812  const vectorField& faceNormals = pp.faceNormals();
813 
814  forAll(edgeFaces, edgeI)
815  {
816  const labelList& eFaces = edgeFaces[edgeI];
817 
818  edgeFaceNormals[edgeI].setSize(eFaces.size());
819  edgeFaceExtrude[edgeI].setSize(eFaces.size());
820  forAll(eFaces, i)
821  {
822  label faceI = eFaces[i];
823  edgeFaceNormals[edgeI][i] = faceNormals[faceI];
824  edgeFaceExtrude[edgeI][i] = extrudedFaces[faceI];
825  }
826  }
827 
829  (
830  mesh(),
831  meshEdges,
832  edgeFaceNormals,
833  globalMeshData::ListPlusEqOp<List<point> >(), // combine operator
834  List<point>() // null value
835  );
836 
838  (
839  mesh(),
840  meshEdges,
841  edgeFaceExtrude,
842  globalMeshData::ListPlusEqOp<List<bool> >(), // combine operator
843  List<bool>() // null value
844  );
845 
846 
847  forAll(edgeFaceNormals, edgeI)
848  {
849  const List<point>& eFaceNormals = edgeFaceNormals[edgeI];
850  const List<bool>& eFaceExtrude = edgeFaceExtrude[edgeI];
851 
852  if (eFaceNormals.size() == 2)
853  {
854  const edge& e = pp.edges()[edgeI];
855  label v0 = e[0];
856  label v1 = e[1];
857 
858  if
859  (
860  extrudeStatus[v0] != autoLayerDriver::NOEXTRUDE
861  || extrudeStatus[v1] != autoLayerDriver::NOEXTRUDE
862  )
863  {
864  if (!eFaceExtrude[0] || !eFaceExtrude[1])
865  {
866  const vector& n0 = eFaceNormals[0];
867  const vector& n1 = eFaceNormals[1];
868 
869  if ((n0 & n1) < minCos)
870  {
871  if (unmarkExtrusion(v0, patchDisp, extrudeStatus))
872  {
873  if (isPatchMasterPoint[v0])
874  {
875  nPointCounter++;
876  }
877  }
878  if (unmarkExtrusion(v1, patchDisp, extrudeStatus))
879  {
880  if (isPatchMasterPoint[v1])
881  {
882  nPointCounter++;
883  }
884  }
885  }
886  }
887  }
888  }
889  }
890 
891  //Info<< "Added "
892  // << returnReduce(nPointCounter-nOldPointCounter, sumOp<label>())
893  // << " point not to extrude due to minCos "
894  // << minCos << endl;
895 }
896 
897 
898 // Find isolated islands (points, edges and faces and layer terminations)
899 // in the layer mesh and stop any layer growth at these points.
901 (
902  const scalar minCosLayerTermination,
903  const bool detectExtrusionIsland,
904  const PackedBoolList& isPatchMasterPoint,
905  const PackedBoolList& isPatchMasterEdge,
906  const labelList& meshEdges,
907  const scalarField& minThickness,
908  List<autoLayerDriver::extrudeMode>& extrudeStatus,
909  pointField& patchDisp
910 ) const
911 {
912  const indirectPrimitivePatch& pp = adaptPatchPtr_();
913  const labelListList& pointFaces = pp.pointFaces();
914  const labelList& meshPoints = pp.meshPoints();
915 
916 
917  Info<< typeName << " : Removing isolated regions ..." << nl
918  << indent << "- if partially extruded faces make angle < "
919  << Foam::radToDeg(Foam::acos(minCosLayerTermination)) << nl;
920  if (detectExtrusionIsland)
921  {
922  Info<< indent << "- if exclusively surrounded by non-extruded points"
923  << nl;
924  }
925  else
926  {
927  Info<< indent << "- if exclusively surrounded by non-extruded faces"
928  << nl;
929  }
930 
931  // Keep count of number of points unextruded
932  label nPointCounter = 0;
933 
934  while (true)
935  {
936  // Stop layer growth where mesh wraps around edge with a
937  // large feature angle
938  if (minCosLayerTermination > -1)
939  {
940  handleFeatureAngleLayerTerminations
941  (
942  minCosLayerTermination,
943  isPatchMasterPoint,
944  meshEdges,
945 
946  extrudeStatus,
947  patchDisp,
948  nPointCounter
949  );
950 
951  syncPatchDisplacement(minThickness, patchDisp, extrudeStatus);
952  }
953 
954 
955  // Detect either:
956  // - point where all surrounding points are not extruded
957  // (detectExtrusionIsland)
958  // or
959  // - point where all the faces surrounding it are not fully
960  // extruded
961 
962  boolList keptPoints(pp.nPoints(), false);
963 
964  if (detectExtrusionIsland)
965  {
966  // Do not extrude from point where all neighbouring
967  // points are not grown
968  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
969 
970  labelList islandPoint(pp.size(), -1);
971  forAll(pp, faceI)
972  {
973  const face& f = pp.localFaces()[faceI];
974 
975  forAll(f, fp)
976  {
977  if (extrudeStatus[f[fp]] != autoLayerDriver::NOEXTRUDE)
978  {
979  if (islandPoint[faceI] == -1)
980  {
981  // First point to extrude
982  islandPoint[faceI] = f[fp];
983  }
984  else if (islandPoint[faceI] != -2)
985  {
986  // Second or more point to extrude
987  islandPoint[faceI] = -2;
988  }
989  }
990  }
991  }
992 
993  // islandPoint:
994  // -1 : no point extruded on face
995  // -2 : >= 2 points extruded on face
996  // >=0: label of point extruded
997 
998  // Check all surrounding faces that I am the islandPoint
999  forAll(pointFaces, patchPointI)
1000  {
1001  if (extrudeStatus[patchPointI] != autoLayerDriver::NOEXTRUDE)
1002  {
1003  const labelList& pFaces = pointFaces[patchPointI];
1004 
1005  forAll(pFaces, i)
1006  {
1007  label faceI = pFaces[i];
1008  if (islandPoint[faceI] != patchPointI)
1009  {
1010  keptPoints[patchPointI] = true;
1011  break;
1012  }
1013  }
1014  }
1015  }
1016  }
1017  else
1018  {
1019  // Do not extrude from point where all neighbouring
1020  // faces are not grown
1021  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1022 
1023  boolList extrudedFaces(pp.size(), true);
1024  forAll(pp.localFaces(), faceI)
1025  {
1026  const face& f = pp.localFaces()[faceI];
1027  forAll(f, fp)
1028  {
1029  if (extrudeStatus[f[fp]] == autoLayerDriver::NOEXTRUDE)
1030  {
1031  extrudedFaces[faceI] = false;
1032  break;
1033  }
1034  }
1035  }
1036 
1037  const labelListList& pointFaces = pp.pointFaces();
1038 
1039  forAll(keptPoints, patchPointI)
1040  {
1041  const labelList& pFaces = pointFaces[patchPointI];
1042 
1043  forAll(pFaces, i)
1044  {
1045  label faceI = pFaces[i];
1046  if (extrudedFaces[faceI])
1047  {
1048  keptPoints[patchPointI] = true;
1049  break;
1050  }
1051  }
1052  }
1053  }
1054 
1055 
1057  (
1058  mesh(),
1059  meshPoints,
1060  keptPoints,
1061  orEqOp<bool>(),
1062  false // null value
1063  );
1064 
1065  label nChanged = 0;
1066 
1067  forAll(keptPoints, patchPointI)
1068  {
1069  if (!keptPoints[patchPointI])
1070  {
1071  if (unmarkExtrusion(patchPointI, patchDisp, extrudeStatus))
1072  {
1073  nPointCounter++;
1074  nChanged++;
1075  }
1076  }
1077  }
1078 
1079 
1080  if (returnReduce(nChanged, sumOp<label>()) == 0)
1081  {
1082  break;
1083  }
1084  }
1085 
1086  const edgeList& edges = pp.edges();
1087 
1088 
1089  // Count number of mesh edges using a point
1090  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1091 
1092  labelList isolatedPoint(pp.nPoints(),0);
1093 
1094  forAll(edges, edgeI)
1095  {
1096  if (isPatchMasterEdge[edgeI])
1097  {
1098  const edge& e = edges[edgeI];
1099 
1100  label v0 = e[0];
1101  label v1 = e[1];
1102 
1103  if (extrudeStatus[v1] != autoLayerDriver::NOEXTRUDE)
1104  {
1105  isolatedPoint[v0] += 1;
1106  }
1107  if (extrudeStatus[v0] != autoLayerDriver::NOEXTRUDE)
1108  {
1109  isolatedPoint[v1] += 1;
1110  }
1111  }
1112  }
1113 
1115  (
1116  mesh(),
1117  meshPoints,
1118  isolatedPoint,
1119  plusEqOp<label>(),
1120  label(0) // null value
1121  );
1122 
1123  // stop layer growth on isolated faces
1124  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1125  forAll(pp, faceI)
1126  {
1127  const face& f = pp.localFaces()[faceI];
1128  bool failed = false;
1129  forAll(f, fp)
1130  {
1131  if (isolatedPoint[f[fp]] > 2)
1132  {
1133  failed = true;
1134  break;
1135  }
1136  }
1137  bool allPointsExtruded = true;
1138  if (!failed)
1139  {
1140  forAll(f, fp)
1141  {
1142  if (extrudeStatus[f[fp]] == autoLayerDriver::NOEXTRUDE)
1143  {
1144  allPointsExtruded = false;
1145  break;
1146  }
1147  }
1148 
1149  if (allPointsExtruded)
1150  {
1151  forAll(f, fp)
1152  {
1153  if
1154  (
1155  unmarkExtrusion
1156  (
1157  f[fp],
1158  patchDisp,
1159  extrudeStatus
1160  )
1161  )
1162  {
1163  nPointCounter++;
1164  }
1165  }
1166  }
1167  }
1168  }
1169 
1170  reduce(nPointCounter, sumOp<label>());
1171  Info<< typeName
1172  << " : Number of isolated points extrusion stopped : "<< nPointCounter
1173  << endl;
1174 }
1175 
1176 
1177 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
1178 
1181  const dictionary& dict,
1182  const List<labelPair>& baffles,
1183  pointVectorField& pointDisplacement
1184 )
1185 :
1186  externalDisplacementMeshMover(dict, baffles, pointDisplacement),
1187  adaptPatchIDs_(getFixedValueBCs(pointDisplacement)),
1188  adaptPatchPtr_(getPatch(mesh(), adaptPatchIDs_)),
1189  scale_
1190  (
1191  IOobject
1192  (
1193  "scale",
1194  pointDisplacement.time().timeName(),
1195  pointDisplacement.db(),
1198  ),
1199  pMesh(),
1200  dimensionedScalar("scale", dimless, 1.0)
1201  ),
1202  oldPoints_(mesh().points()),
1203  meshMover_
1204  (
1205  const_cast<polyMesh&>(mesh()),
1206  const_cast<pointMesh&>(pMesh()),
1207  adaptPatchPtr_(),
1208  pointDisplacement,
1209  scale_,
1210  oldPoints_,
1211  adaptPatchIDs_,
1212  dict
1213  ),
1214  fieldSmoother_(mesh()),
1215  dispVec_
1216  (
1217  IOobject
1218  (
1219  "dispVec",
1220  pointDisplacement.time().timeName(),
1221  pointDisplacement.db(),
1224  false
1225  ),
1226  pMesh(),
1228  ),
1229  medialRatio_
1230  (
1231  IOobject
1232  (
1233  "medialRatio",
1234  pointDisplacement.time().timeName(),
1235  pointDisplacement.db(),
1238  false
1239  ),
1240  pMesh(),
1241  dimensionedScalar("medialRatio", dimless, 0.0)
1242  ),
1243  medialDist_
1244  (
1245  IOobject
1246  (
1247  "pointMedialDist",
1248  pointDisplacement.time().timeName(),
1249  pointDisplacement.db(),
1252  false
1253  ),
1254  pMesh(),
1255  dimensionedScalar("pointMedialDist", dimLength, 0.0)
1256  ),
1257  medialVec_
1258  (
1259  IOobject
1260  (
1261  "medialVec",
1262  pointDisplacement.time().timeName(),
1263  pointDisplacement.db(),
1266  false
1267  ),
1268  pMesh(),
1270  )
1271 {
1272  update(dict);
1273 }
1274 
1275 
1276 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
1277 
1279 {}
1280 
1281 
1282 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1283 
1286  const dictionary& coeffDict,
1287  const scalarField& minThickness,
1288  List<autoLayerDriver::extrudeMode>& extrudeStatus,
1289  pointField& patchDisp
1290 )
1291 {
1292  Info<< typeName << " : Smoothing using Medial Axis ..." << endl;
1293 
1294  const indirectPrimitivePatch& pp = adaptPatchPtr_;
1295  const labelList& meshPoints = pp.meshPoints();
1296 
1297 
1298  // Read settings
1299  // ~~~~~~~~~~~~~
1300 
1301  //- (lambda-mu) smoothing of internal displacement
1302  const label nSmoothDisplacement = coeffDict.lookupOrDefault
1303  (
1304  "nSmoothDisplacement",
1305  0
1306  );
1307 
1308  //- Layer thickness too big
1309  const scalar maxThicknessToMedialRatio = readScalar
1310  (
1311  coeffDict.lookup("maxThicknessToMedialRatio")
1312  );
1313 
1314  //- Feature angle when to stop adding layers
1315  const scalar featureAngle = readScalar(coeffDict.lookup("featureAngle"));
1316 
1317  //- Stop layer growth where mesh wraps around sharp edge
1318  scalar layerTerminationAngle = coeffDict.lookupOrDefault<scalar>
1319  (
1320  "layerTerminationAngle",
1321  0.5*featureAngle
1322  );
1323  scalar minCosLayerTermination = Foam::cos(degToRad(layerTerminationAngle));
1324 
1325  //- Smoothing wanted patch thickness
1326  const label nSmoothPatchThickness = readLabel
1327  (
1328  coeffDict.lookup("nSmoothThickness")
1329  );
1330 
1331  //- Number of edges walking out
1332  const label nMedialAxisIter = coeffDict.lookupOrDefault<label>
1333  (
1334  "nMedialAxisIter",
1336  );
1337 
1338  //- Use strict extrusionIsland detection
1339  const Switch detectExtrusionIsland = coeffDict.lookupOrDefault<Switch>
1340  (
1341  "detectExtrusionIsland",
1342  false
1343  );
1344 
1345 
1346  // Precalulate master points/edge (only relevant for shared points/edges)
1347  const PackedBoolList isMeshMasterPoint(syncTools::getMasterPoints(mesh()));
1348  const PackedBoolList isMeshMasterEdge(syncTools::getMasterEdges(mesh()));
1349  // Precalculate meshEdge per pp edge
1350  const labelList meshEdges
1351  (
1352  pp.meshEdges
1353  (
1354  mesh().edges(),
1355  mesh().pointEdges()
1356  )
1357  );
1358 
1359  // Precalulate (patch) master point/edge
1360  const PackedBoolList isPatchMasterPoint
1361  (
1363  (
1364  mesh(),
1365  meshPoints
1366  )
1367  );
1368  const PackedBoolList isPatchMasterEdge
1369  (
1371  (
1372  mesh(),
1373  meshEdges
1374  )
1375  );
1376 
1377 
1378  scalarField thickness(mag(patchDisp));
1379 
1380  forAll(thickness, patchPointI)
1381  {
1382  if (extrudeStatus[patchPointI] == autoLayerDriver::NOEXTRUDE)
1383  {
1384  thickness[patchPointI] = 0.0;
1385  }
1386  }
1387 
1388  label numThicknessRatioExclude = 0;
1389 
1390  // reduce thickness where thickness/medial axis distance large
1391  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1392 
1393  autoPtr<OBJstream> str;
1394  if (debug)
1395  {
1396  str.reset
1397  (
1398  new OBJstream
1399  (
1400  mesh().time().path()
1401  / "thicknessRatioExcludePoints_"
1402  + mesh().time().timeName()
1403  + ".obj"
1404  )
1405  );
1406  Info<< typeName
1407  << " : Writing points with too large an extrusion distance to "
1408  << str().name() << endl;
1409  }
1410 
1411  autoPtr<OBJstream> medialVecStr;
1412  if (debug)
1413  {
1414  medialVecStr.reset
1415  (
1416  new OBJstream
1417  (
1418  mesh().time().path()
1419  / "thicknessRatioExcludeMedialVec_"
1420  + mesh().time().timeName()
1421  + ".obj"
1422  )
1423  );
1424  Info<< typeName
1425  << " : Writing medial axis vectors on points with too large"
1426  << " an extrusion distance to " << medialVecStr().name() << endl;
1427  }
1428 
1429  forAll(meshPoints, patchPointI)
1430  {
1431  if (extrudeStatus[patchPointI] != autoLayerDriver::NOEXTRUDE)
1432  {
1433  label pointI = meshPoints[patchPointI];
1434 
1435  //- Option 1: look only at extrusion thickness v.s. distance
1436  // to nearest (medial axis or static) point.
1437  scalar mDist = medialDist_[pointI];
1438  scalar thicknessRatio = thickness[patchPointI]/(mDist+VSMALL);
1439 
1440  //- Option 2: Look at component in the direction
1441  // of nearest (medial axis or static) point.
1442  vector n =
1443  patchDisp[patchPointI]
1444  / (mag(patchDisp[patchPointI]) + VSMALL);
1445  vector mVec = medialVec_[pointI]-mesh().points()[pointI];
1446  mVec /= mag(mVec)+VSMALL;
1447  thicknessRatio *= (n&mVec);
1448 
1449  if (thicknessRatio > maxThicknessToMedialRatio)
1450  {
1451  // Truncate thickness.
1452  if (debug&2)
1453  {
1454  Pout<< "truncating displacement at "
1455  << mesh().points()[pointI]
1456  << " from " << thickness[patchPointI]
1457  << " to "
1458  << 0.5
1459  *(
1460  minThickness[patchPointI]
1461  +thickness[patchPointI]
1462  )
1463  << " medial direction:" << mVec
1464  << " extrusion direction:" << n
1465  << " with thicknessRatio:" << thicknessRatio
1466  << endl;
1467  }
1468 
1469  thickness[patchPointI] =
1470  0.5*(minThickness[patchPointI]+thickness[patchPointI]);
1471 
1472  patchDisp[patchPointI] = thickness[patchPointI]*n;
1473 
1474  if (isPatchMasterPoint[patchPointI])
1475  {
1476  numThicknessRatioExclude++;
1477  }
1478 
1479  if (str.valid())
1480  {
1481  const point& pt = mesh().points()[pointI];
1482  str().write(linePointRef(pt, pt+patchDisp[patchPointI]));
1483  }
1484  if (medialVecStr.valid())
1485  {
1486  const point& pt = mesh().points()[pointI];
1487  medialVecStr().write
1488  (
1489  linePointRef
1490  (
1491  pt,
1492  medialVec_[pointI]
1493  )
1494  );
1495  }
1496  }
1497  }
1498  }
1499 
1500  reduce(numThicknessRatioExclude, sumOp<label>());
1501  Info<< typeName << " : Reducing layer thickness at "
1502  << numThicknessRatioExclude
1503  << " nodes where thickness to medial axis distance is large " << endl;
1504 
1505 
1506  // find points where layer growth isolated to a lone point, edge or face
1507 
1508  findIsolatedRegions
1509  (
1510  minCosLayerTermination,
1511  detectExtrusionIsland,
1512 
1513  isPatchMasterPoint,
1514  isPatchMasterEdge,
1515  meshEdges,
1516  minThickness,
1517 
1518  extrudeStatus,
1519  patchDisp
1520  );
1521 
1522  // Update thickess for changed extrusion
1523  forAll(thickness, patchPointI)
1524  {
1525  if (extrudeStatus[patchPointI] == autoLayerDriver::NOEXTRUDE)
1526  {
1527  thickness[patchPointI] = 0.0;
1528  }
1529  }
1530 
1531 
1532  // Smooth layer thickness on moving patch. Since some locations will have
1533  // disabled the extrusion this might cause big jumps in wanted displacement
1534  // for neighbouring patch points. So smooth the wanted displacement
1535  // before actually trying to move the mesh.
1536  fieldSmoother_.minSmoothField
1537  (
1538  nSmoothPatchThickness,
1539  isPatchMasterPoint,
1540  isPatchMasterEdge,
1541  pp,
1542  minThickness,
1543  thickness
1544  );
1545 
1546 
1547  // Dummy additional info for PointEdgeWave
1548  int dummyTrackData = 0;
1549 
1550  List<PointData<scalar> > pointWallDist(mesh().nPoints());
1551 
1552  const pointField& points = mesh().points();
1553  // 1. Calculate distance to points where displacement is specified.
1554  // This wave is used to transport layer thickness
1555  {
1556  // Distance to wall and medial axis on edges.
1557  List<PointData<scalar> > edgeWallDist(mesh().nEdges());
1558  labelList wallPoints(meshPoints.size());
1559 
1560  // Seed data.
1561  List<PointData<scalar> > wallInfo(meshPoints.size());
1562 
1563  forAll(meshPoints, patchPointI)
1564  {
1565  label pointI = meshPoints[patchPointI];
1566  wallPoints[patchPointI] = pointI;
1567  wallInfo[patchPointI] = PointData<scalar>
1568  (
1569  points[pointI],
1570  0.0,
1571  thickness[patchPointI] // transport layer thickness
1572  );
1573  }
1574 
1575  // Do all calculations
1576  PointEdgeWave<PointData<scalar> > wallDistCalc
1577  (
1578  mesh(),
1579  wallPoints,
1580  wallInfo,
1581  pointWallDist,
1582  edgeWallDist,
1583  0, // max iterations
1584  dummyTrackData
1585  );
1586  wallDistCalc.iterate(nMedialAxisIter);
1587  }
1588 
1589 
1590  // Calculate scaled displacement vector
1591  pointField& displacement = pointDisplacement_;
1592 
1593  forAll(displacement, pointI)
1594  {
1595  if (!pointWallDist[pointI].valid(dummyTrackData))
1596  {
1597  displacement[pointI] = vector::zero;
1598  }
1599  else
1600  {
1601  // 1) displacement on nearest wall point, scaled by medialRatio
1602  // (wall distance / medial distance)
1603  // 2) pointWallDist[pointI].data() is layer thickness transported
1604  // from closest wall point.
1605  // 3) shrink in opposite direction of addedPoints
1606  displacement[pointI] =
1607  -medialRatio_[pointI]
1608  *pointWallDist[pointI].data()
1609  *dispVec_[pointI];
1610  }
1611  }
1612 
1613 
1614  // Smear displacement away from fixed values (medialRatio=0 or 1)
1615  if (nSmoothDisplacement > 0)
1616  {
1617  PackedBoolList isToBeSmoothed(displacement.size(), false);
1618 
1619  forAll(displacement, i)
1620  {
1621  isToBeSmoothed[i] =
1622  medialRatio_[i] > SMALL && medialRatio_[i] < 1-SMALL;
1623  }
1624 
1625  fieldSmoother_.smoothLambdaMuDisplacement
1626  (
1627  nSmoothDisplacement,
1628  isMeshMasterPoint,
1629  isMeshMasterEdge,
1630  isToBeSmoothed,
1631  displacement
1632  );
1633  }
1634 }
1635 
1636 
1639  const dictionary& meshQualityDict,
1640  const label nAllowableErrors,
1641  labelList& checkFaces
1642 )
1643 {
1644  //- Number of attempts shrinking the mesh
1645  const label nSnap = readLabel(meshQualityDict.lookup("nRelaxIter"));
1646 
1647 
1648  // Make sure displacement boundary conditions is uptodate with
1649  // internal field
1650  meshMover_.setDisplacementPatchFields();
1651 
1652  Info<< typeName << " : Moving mesh ..." << endl;
1653  scalar oldErrorReduction = -1;
1654 
1655  bool meshOk = false;
1656 
1657  for (label iter = 0; iter < 2*nSnap ; iter++)
1658  {
1659  Info<< typeName
1660  << " : Iteration " << iter << endl;
1661  if (iter == nSnap)
1662  {
1663  Info<< typeName
1664  << " : Displacement scaling for error reduction set to 0."
1665  << endl;
1666  oldErrorReduction = meshMover_.setErrorReduction(0.0);
1667  }
1668 
1669  if
1670  (
1671  meshMover_.scaleMesh
1672  (
1673  checkFaces,
1674  baffles_,
1675  meshMover_.paramDict(),
1676  meshQualityDict,
1677  true,
1678  nAllowableErrors
1679  )
1680  )
1681  {
1682  Info<< typeName << " : Successfully moved mesh" << endl;
1683  meshOk = true;
1684  break;
1685  }
1686  }
1687 
1688  if (oldErrorReduction >= 0)
1689  {
1690  meshMover_.setErrorReduction(oldErrorReduction);
1691  }
1692 
1693  Info<< typeName << " : Finished moving mesh ..." << endl;
1694 
1695  return meshOk;
1696 }
1697 
1698 
1701  const dictionary& moveDict,
1702  const label nAllowableErrors,
1703  labelList& checkFaces
1704 )
1705 {
1706  // Read a few settings
1707  // ~~~~~~~~~~~~~~~~~~~
1708 
1709  //- Name of field specifying min thickness
1710  const word minThicknessName = word(moveDict.lookup("minThicknessName"));
1711 
1712 
1713  // Extract out patch-wise displacement
1714  const indirectPrimitivePatch& pp = adaptPatchPtr_();
1715 
1716  scalarField zeroMinThickness;
1717  if (minThicknessName == "none")
1718  {
1719  zeroMinThickness = scalarField(pp.nPoints(), 0.0);
1720  }
1721  const scalarField& minThickness =
1722  (
1723  (minThicknessName == "none")
1724  ? zeroMinThickness
1725  : mesh().lookupObject<scalarField>(minThicknessName)
1726  );
1727 
1728 
1729  pointField patchDisp
1730  (
1731  pointDisplacement_.internalField(),
1732  pp.meshPoints()
1733  );
1734 
1736  (
1737  pp.nPoints(),
1739  );
1740  forAll(extrudeStatus, pointI)
1741  {
1742  if (mag(patchDisp[pointI]) <= minThickness[pointI]+SMALL)
1743  {
1744  extrudeStatus[pointI] = autoLayerDriver::NOEXTRUDE;
1745  }
1746  }
1747 
1748 
1749  // Solve displacement
1750  calculateDisplacement(moveDict, minThickness, extrudeStatus, patchDisp);
1751 
1752  //- Move mesh according to calculated displacement
1753  return shrinkMesh
1754  (
1755  moveDict, // meshQualityDict,
1756  nAllowableErrors, // nAllowableErrors
1757  checkFaces
1758  );
1759 }
1760 
1761 
1763 {
1765 
1766  // Update motionSmoother for new geometry (moves adaptPatchPtr_)
1767  meshMover_.movePoints();
1768 
1769  // Assume corrent mesh location is correct (reset oldPoints, scale)
1770  meshMover_.correct();
1771 }
1772 
1773 
1774 // ************************************************************************* //
Foam::medialAxisMeshMover::handleFeatureAngleLayerTerminations
void handleFeatureAngleLayerTerminations(const scalar minCos, const PackedBoolList &isMasterPoint, const labelList &meshEdges, List< autoLayerDriver::extrudeMode > &extrudeStatus, pointField &patchDisp, label &nPointCounter) const
Stop layer growth at feature edges.
Definition: medialAxisMeshMover.C:769
Foam::medialAxisMeshMover::syncPatchDisplacement
void syncPatchDisplacement(const scalarField &minThickness, pointField &patchDisp, List< autoLayerDriver::extrudeMode > &extrudeStatus) const
Synchronise extrusion.
Definition: medialAxisMeshMover.C:650
Foam::polyMesh::points
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:979
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::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
Foam::Switch
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:60
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
Foam::PrimitivePatch::edgeFaces
const labelListList & edgeFaces() const
Return edge-face addressing.
Definition: PrimitivePatchTemplate.C:292
Foam::medialAxisMeshMover::medialRatio_
pointScalarField medialRatio_
Ratio of medial distance to wall distance.
Definition: medialAxisMeshMover.H:94
Foam::Vector< scalar >::zero
static const Vector zero
Definition: Vector.H:80
p
p
Definition: pEqn.H:62
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:86
Foam::OBJstream
OFstream which keeps track of vertices.
Definition: OBJstream.H:53
Foam::polyBoundaryMesh
Foam::polyBoundaryMesh.
Definition: polyBoundaryMesh.H:60
Foam::PrimitivePatch::edges
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
Definition: PrimitivePatchTemplate.C:212
Foam::IOobject::AUTO_WRITE
@ AUTO_WRITE
Definition: IOobject.H:117
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::DynamicList
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:56
Foam::autoLayerDriver::EXTRUDEREMOVE
@ EXTRUDEREMOVE
Definition: autoLayerDriver.H:67
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
Foam::PointEdgeWave
Wave propagation of information through grid. Every iteration information goes through one layer of e...
Definition: PointEdgeWave.H:85
Foam::PrimitivePatch::nEdges
label nEdges() const
Return number of edges in patch.
Definition: PrimitivePatchTemplate.H:299
Foam::medialAxisMeshMover::~medialAxisMeshMover
virtual ~medialAxisMeshMover()
Definition: medialAxisMeshMover.C:1278
Foam::globalMeshData::nTotalPoints
label nTotalPoints() const
Return total number of points in decomposed mesh. Not.
Definition: globalMeshData.H:381
Foam::polyPatch::coupled
virtual bool coupled() const
Return true if this patch is geometrically coupled (i.e. faces and.
Definition: polyPatch.H:322
zeroFixedValuePointPatchFields.H
Foam::dictionary::lookup
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:449
Foam::PrimitivePatch::meshEdges
labelList meshEdges(const edgeList &allEdges, const labelListList &cellEdges, const labelList &faceCells) const
Return labels of patch edges in the global edge list using.
Definition: PrimitivePatchMeshEdges.C:41
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::medialAxisMeshMover::calculateDisplacement
void calculateDisplacement(const dictionary &, const scalarField &minThickness, List< autoLayerDriver::extrudeMode > &extrudeStatus, pointField &patchDisp)
Calculate desired displacement. Modifies per-patch displacement.
Definition: medialAxisMeshMover.C:1285
Foam::dictionary::lookupOrDefault
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
Definition: dictionaryTemplates.C:33
unitConversion.H
Unit conversion functions.
Foam::GeometricField::boundaryField
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
Definition: GeometricField.C:735
valuePointPatchFields.H
medialAxisMeshMover.H
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::primitiveMesh::edges
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
Definition: primitiveMeshEdges.C:353
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::HashSet< label, Hash< label > >
Foam::fieldSmoother::smoothNormals
void smoothNormals(const label nIter, const PackedBoolList &isMeshMasterPoint, const PackedBoolList &isMeshMasterEdge, const labelList &fixedPoints, pointVectorField &normals) const
Smooth interior normals.
Definition: fieldSmoother.C:53
Foam::incrIndent
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:228
Foam::PointData
Variant of pointEdgePoint with some transported additional data. Templated on the transported data ty...
Definition: medialAxisMeshMover.H:57
Foam::pointPatchField
Abstract base class for point-mesh patch fields.
Definition: pointMVCWeight.H:58
Foam::globalMeshData::ListPlusEqOp
Definition: globalMeshData.H:320
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
Foam::medialAxisMeshMover::update
void update(const dictionary &)
Initialise medial axis. Uses pointDisplacement field only.
Definition: medialAxisMeshMover.C:115
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
Foam::externalDisplacementMeshMover::pointDisplacement
pointVectorField & pointDisplacement()
Return reference to the point motion displacement field.
Definition: externalDisplacementMeshMover.H:151
Foam::dimensionedVector
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Definition: dimensionedVector.H:48
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
nPoints
label nPoints
Definition: gmvOutputHeader.H:2
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::medialAxisMeshMover::fieldSmoother_
fieldSmoother fieldSmoother_
Field smoothing.
Definition: medialAxisMeshMover.H:83
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:43
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::dictionary::found
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:304
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::medialAxisMeshMover::isMaxEdge
bool isMaxEdge(const List< PointData< vector > > &pointWallDist, const label edgeI, const scalar minCos) const
Is mesh edge on a cusp of displacement.
Definition: medialAxisMeshMover.C:58
Foam::PointEdgeWave::iterate
label iterate(const label maxIter)
Iterate until no changes or maxIter reached. Returns actual.
Definition: PointEdgeWave.C:914
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
PointEdgeWave.H
Foam::externalDisplacementMeshMover::movePoints
virtual void movePoints(const pointField &)
Update local data for geometry changes.
Definition: externalDisplacementMeshMover.C:172
Foam::externalDisplacementMeshMover
Virtual base class for mesh movers with externally provided displacement field giving the boundary co...
Definition: externalDisplacementMeshMover.H:56
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:41
Foam::DynamicList::shrink
DynamicList< T, SizeInc, SizeMult, SizeDiv > & shrink()
Shrink the allocated space to the number of elements used.
Definition: DynamicListI.H:258
Foam::orEqOp
Definition: ops.H:82
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::PrimitivePatch::nPoints
label nPoints() const
Return number of points supporting patch faces.
Definition: PrimitivePatchTemplate.H:293
Foam::medialAxisMeshMover::movePoints
virtual void movePoints(const pointField &)
Update local data for geometry changes.
Definition: medialAxisMeshMover.C:1762
Foam::medialAxisMeshMover::medialDist_
pointScalarField medialDist_
Distance to nearest medial axis point.
Definition: medialAxisMeshMover.H:97
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::PointEdgeWave::getUnsetPoints
label getUnsetPoints() const
Definition: PointEdgeWave.C:717
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
meshRefinement.H
Foam::plusEqOp
Definition: ops.H:71
Foam::HashTable::found
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:109
Foam::autoLayerDriver::NOEXTRUDE
@ NOEXTRUDE
Definition: autoLayerDriver.H:65
Foam::autoLayerDriver::EXTRUDE
@ EXTRUDE
Definition: autoLayerDriver.H:66
Foam::medialAxisMeshMover::dispVec_
pointVectorField dispVec_
Normal of nearest wall. Where this normal changes direction.
Definition: medialAxisMeshMover.H:90
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::pointMesh
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:48
Foam::decrIndent
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:235
Foam::fieldSmoother::smoothPatchNormals
void smoothPatchNormals(const label nIter, const PackedBoolList &isPatchMasterPoint, const PackedBoolList &isPatchMasterEdge, const indirectPrimitivePatch &adaptPatch, pointField &normals) const
Smooth patch normals.
Definition: fieldSmoother.C:141
Foam::e
const double e
Elementary charge.
Definition: doubleFloat.H:94
Foam::indent
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:221
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
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::meshRefinement::getMasterEdges
static PackedBoolList getMasterEdges(const polyMesh &mesh, const labelList &meshEdges)
Determine master edge for subset of edges. If coupled.
Definition: meshRefinement.C:2784
Foam::pointEdgePoint
Holds information regarding nearest wall point. Used in PointEdgeWave. (so not standard FaceCellWave)...
Definition: pointEdgePoint.H:59
Foam::medialAxisMeshMover::medialAxisMeshMover
medialAxisMeshMover(const medialAxisMeshMover &)
Disallow default bitwise copy construct.
Foam::DynamicList::append
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Foam::syncTools::getMasterEdges
static PackedBoolList getMasterEdges(const polyMesh &)
Get per edge whether it is uncoupled or a master of a.
Definition: syncTools.C:109
Foam::List::setSize
void setSize(const label)
Reset size of List.
Foam::autoPtr
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:117
Foam::linePointRef
line< point, const point & > linePointRef
Line using referred points.
Definition: linePointRef.H:45
Foam::meshRefinement::getMasterPoints
static PackedBoolList getMasterPoints(const polyMesh &mesh, const labelList &meshPoints)
Determine master point for subset of points. If coupled.
Definition: meshRefinement.C:2747
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:49
Foam::syncTools::syncEdgeList
static void syncEdgeList(const polyMesh &, List< T > &, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronize values on all mesh edges.
Definition: syncToolsTemplates.C:1141
Foam::sumOp
Definition: ops.H:162
Foam::readScalar
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
Foam::Pout
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
Foam::PrimitivePatch::faceNormals
const Field< PointType > & faceNormals() const
Return face normals for patch.
Definition: PrimitivePatchTemplate.C:520
f
labelList f(nPoints)
PointData.H
Foam::PatchTools::pointNormals
static tmp< pointField > pointNormals(const polyMesh &, const PrimitivePatch< Face, FaceList, PointField, PointType > &)
Return parallel consistent point normals for patches using mesh points.
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
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:142
Foam::acos
dimensionedScalar acos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:259
Foam::autoPtr::valid
bool valid() const
Return true if the autoPtr valid (ie, the pointer is set).
Definition: autoPtrI.H:83
Foam::Vector< scalar >::rootMax
static const Vector rootMax
Definition: Vector.H:84
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: doubleFloat.H:94
Foam::medialAxisMeshMover::unmarkExtrusion
static bool unmarkExtrusion(const label patchPointI, pointField &patchDisp, List< autoLayerDriver::extrudeMode > &extrudeStatus)
Unmark extrusion at given point.
Definition: medialAxisMeshMover.C:624
path
fileName path(UMean.rootPath()/UMean.caseName()/"graphs"/UMean.instance())
timeName
word timeName
Definition: getTimeIndex.H:3
Foam::syncTools::getMasterPoints
static PackedBoolList getMasterPoints(const polyMesh &)
Get per point whether it is uncoupled or a master of a.
Definition: syncTools.C:65
Foam::readLabel
label readLabel(Istream &is)
Definition: label.H:64
Foam::line
A line primitive.
Definition: line.H:56
Foam::radToDeg
scalar radToDeg(const scalar rad)
Conversion from radians to degrees.
Definition: unitConversion.H:51
Foam::medialAxisMeshMover::adaptPatchPtr_
autoPtr< indirectPrimitivePatch > adaptPatchPtr_
Definition: medialAxisMeshMover.H:71
Foam::medialAxisMeshMover::move
virtual bool move(const dictionary &, const label nAllowableErrors, labelList &checkFaces)
Move mesh using current pointDisplacement boundary values.
Definition: medialAxisMeshMover.C:1700
Foam::autoPtr::reset
void reset(T *=0)
If object pointer already set, delete object and set to given.
Definition: autoPtrI.H:114
Foam::PrimitivePatch::localFaces
const List< Face > & localFaces() const
Return patch faces addressing into local point list.
Definition: PrimitivePatchTemplate.C:372
patches
patches[0]
Definition: createSingleCellMesh.H:36
Foam::pointPatchField::fixesValue
virtual bool fixesValue() const
Return true if this patch field fixes a value.
Definition: pointPatchField.H:302
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
Foam::minMagSqrEqOp
Definition: ops.H:79
Foam::objectRegistry::lookupObject
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
Definition: objectRegistryTemplates.C:165
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::medialAxisMeshMover::medialVec_
pointVectorField medialVec_
Location on nearest medial axis point.
Definition: medialAxisMeshMover.H:100
Foam::polyMesh::globalData
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1143
Foam::medialAxisMeshMover::adaptPatchIDs_
const labelList adaptPatchIDs_
Definition: medialAxisMeshMover.H:69
Foam::point
vector point
Point is a vector.
Definition: point.H:41
Foam::PrimitivePatch::meshPoints
const labelList & meshPoints() const
Return labelList of mesh points in patch.
Definition: PrimitivePatchTemplate.C:392
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::patchIdentifier::name
const word & name() const
Return name.
Definition: patchIdentifier.H:109
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:259
Foam::data
Database for solution data, solver performance and other reduced data.
Definition: data.H:52
Foam::magSqr
dimensioned< scalar > magSqr(const dimensioned< Type > &)
pointFields.H
OBJstream.H
Foam::medialAxisMeshMover::findIsolatedRegions
void findIsolatedRegions(const scalar minCosLayerTermination, const bool detectExtrusionIsland, const PackedBoolList &isMasterPoint, const PackedBoolList &isMasterEdge, const labelList &meshEdges, const scalarField &minThickness, List< autoLayerDriver::extrudeMode > &extrudeStatus, pointField &patchDisp) const
Find isolated islands (points, edges and faces and layer.
Definition: medialAxisMeshMover.C:901
Foam::syncTools::syncPointList
static void syncPointList(const polyMesh &, List< T > &, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronize values on all mesh points.
Definition: syncToolsTemplates.C:984
Foam::externalDisplacementMeshMover::mesh
const polyMesh & mesh() const
Definition: externalDisplacementMeshMover.H:167
Foam::medialAxisMeshMover::shrinkMesh
bool shrinkMesh(const dictionary &meshQualityDict, const label nAllowableErrors, labelList &checkFaces)
Move mesh according to calculated displacement.
Definition: medialAxisMeshMover.C:1638
Foam::PrimitivePatch
A list of faces which address into the list of points.
Definition: PrimitivePatchTemplate.H:88
Foam::degToRad
scalar degToRad(const scalar deg)
Conversion from degrees to radians.
Definition: unitConversion.H:45
Foam::cos
dimensionedScalar cos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:256