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