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