medialAxisMeshMover.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2014-2019 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "medialAxisMeshMover.H"
28 #include "pointFields.H"
29 #include "valuePointPatchFields.H"
30 #include "PointEdgeWave.H"
31 #include "meshRefinement.H"
32 #include "unitConversion.H"
33 #include "PatchTools.H"
34 #include "OBJstream.H"
35 #include "pointData.H"
37 
38 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39 
40 namespace Foam
41 {
42  defineTypeNameAndDebug(medialAxisMeshMover, 0);
43 
45  (
46  externalDisplacementMeshMover,
47  medialAxisMeshMover,
48  dictionary
49  );
50 }
51 
52 
53 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
54 
55 Foam::labelList Foam::medialAxisMeshMover::getFixedValueBCs
56 (
57  const pointVectorField& fld
58 )
59 {
60  DynamicList<label> adaptPatchIDs;
61  forAll(fld.boundaryField(), patchi)
62  {
63  const pointPatchField<vector>& patchFld =
64  fld.boundaryField()[patchi];
65 
66  if (isA<valuePointPatchField<vector>>(patchFld))
67  {
68  if (isA<zeroFixedValuePointPatchField<vector>>(patchFld))
69  {
70  // Special condition of fixed boundary condition. Does not
71  // get adapted
72  }
73  else
74  {
75  adaptPatchIDs.append(patchi);
76  }
77  }
78  }
79 
80  return Foam::move(adaptPatchIDs);
81 }
82 
83 
85 Foam::medialAxisMeshMover::getPatch
86 (
87  const polyMesh& mesh,
88  const labelList& patchIDs
89 )
90 {
91  const polyBoundaryMesh& patches = mesh.boundaryMesh();
92 
93  // Count faces.
94  label nFaces = 0;
95 
96  forAll(patchIDs, i)
97  {
98  const polyPatch& pp = patches[patchIDs[i]];
99 
100  nFaces += pp.size();
101  }
102 
103  // Collect faces.
104  labelList addressing(nFaces);
105  nFaces = 0;
106 
107  forAll(patchIDs, i)
108  {
109  const polyPatch& pp = patches[patchIDs[i]];
110 
111  label meshFacei = pp.start();
112 
113  forAll(pp, i)
114  {
115  addressing[nFaces++] = meshFacei++;
116  }
117  }
118 
119  return autoPtr<indirectPrimitivePatch>
120  (
122  (
123  IndirectList<face>(mesh.faces(), addressing),
124  mesh.points()
125  )
126  );
127 }
128 
129 
130 void Foam::medialAxisMeshMover::smoothPatchNormals
131 (
132  const label nSmoothDisp,
133  const PackedBoolList& isPatchMasterPoint,
134  const PackedBoolList& isPatchMasterEdge,
135  pointField& normals
136 ) const
137 {
138  const indirectPrimitivePatch& pp = adaptPatchPtr_();
139  const edgeList& edges = pp.edges();
140  const labelList& meshPoints = pp.meshPoints();
141 
142  // Get smoothly varying internal normals field.
143  Info<< typeName << " : Smoothing normals ..." << endl;
144 
145  scalarField edgeWeights(edges.size());
146  scalarField invSumWeight(meshPoints.size());
148  (
149  mesh(),
150  isPatchMasterEdge,
151  meshPoints,
152  edges,
153  edgeWeights,
154  invSumWeight
155  );
156 
157 
159  for (label iter = 0; iter < nSmoothDisp; iter++)
160  {
162  (
163  mesh(),
164  isPatchMasterEdge,
165  meshPoints,
166  edges,
167  edgeWeights,
168  normals,
169  average
170  );
171  average *= invSumWeight;
172 
173  // Do residual calculation every so often.
174  if ((iter % 10) == 0)
175  {
176  scalar resid = meshRefinement::gAverage
177  (
178  isPatchMasterPoint,
179  mag(normals-average)()
180  );
181  Info<< " Iteration " << iter << " residual " << resid << endl;
182  }
183 
184  // Transfer to normals vector field
185  forAll(average, pointi)
186  {
187  // full smoothing neighbours + point value
188  average[pointi] = 0.5*(normals[pointi]+average[pointi]);
189  normals[pointi] = average[pointi];
190  normals[pointi] /= mag(normals[pointi]) + vSmall;
191  }
192  }
193 }
194 
195 
196 // Smooth normals in interior.
197 void Foam::medialAxisMeshMover::smoothNormals
198 (
199  const label nSmoothDisp,
200  const PackedBoolList& isMeshMasterPoint,
201  const PackedBoolList& isMeshMasterEdge,
202  const labelList& fixedPoints,
203  pointVectorField& normals
204 ) const
205 {
206  // Get smoothly varying internal normals field.
207  Info<< typeName
208  << " : Smoothing normals in interior ..." << endl;
209 
210  const edgeList& edges = mesh().edges();
211 
212  // Points that do not change.
213  PackedBoolList isFixedPoint(mesh().nPoints());
214 
215  // Internal points that are fixed
216  forAll(fixedPoints, i)
217  {
218  label meshPointi = fixedPoints[i];
219  isFixedPoint.set(meshPointi, 1);
220  }
221 
222  // Make sure that points that are coupled to meshPoints but not on a patch
223  // are fixed as well
224  syncTools::syncPointList(mesh(), isFixedPoint, maxEqOp<unsigned int>(), 0);
225 
226 
227  // Correspondence between local edges/points and mesh edges/points
228  const labelList meshPoints(identity(mesh().nPoints()));
229 
230  // Calculate inverse sum of weights
231 
232  scalarField edgeWeights(mesh().nEdges());
233  scalarField invSumWeight(meshPoints.size());
235  (
236  mesh(),
237  isMeshMasterEdge,
238  meshPoints,
239  edges,
240  edgeWeights,
241  invSumWeight
242  );
243 
245  for (label iter = 0; iter < nSmoothDisp; iter++)
246  {
248  (
249  mesh(),
250  isMeshMasterEdge,
251  meshPoints,
252  edges,
253  edgeWeights,
254  normals,
255  average
256  );
257  average *= invSumWeight;
258 
259  // Do residual calculation every so often.
260  if ((iter % 10) == 0)
261  {
262  scalar resid = meshRefinement::gAverage
263  (
264  isMeshMasterPoint,
265  mag(normals-average)()
266  );
267  Info<< " Iteration " << iter << " residual " << resid << endl;
268  }
269 
270 
271  // Transfer to normals vector field
272  forAll(average, pointi)
273  {
274  if (isFixedPoint.get(pointi) == 0)
275  {
276  // full smoothing neighbours + point value
277  average[pointi] = 0.5*(normals[pointi]+average[pointi]);
278  normals[pointi] = average[pointi];
279  normals[pointi] /= mag(normals[pointi]) + vSmall;
280  }
281  }
282  }
283 }
284 
285 
286 // Tries and find a medial axis point. Done by comparing vectors to nearest
287 // wall point for both vertices of edge.
288 bool Foam::medialAxisMeshMover::isMaxEdge
289 (
290  const List<pointData>& pointWallDist,
291  const label edgeI,
292  const scalar minCos
293 ) const
294 {
295  const pointField& points = mesh().points();
296 
297  // Do not mark edges with one side on moving wall.
298 
299  const edge& e = mesh().edges()[edgeI];
300 
301  vector v0(points[e[0]] - pointWallDist[e[0]].origin());
302  scalar magV0(mag(v0));
303 
304  if (magV0 < small)
305  {
306  return false;
307  }
308 
309  vector v1(points[e[1]] - pointWallDist[e[1]].origin());
310  scalar magV1(mag(v1));
311 
312  if (magV1 < small)
313  {
314  return false;
315  }
316 
317 
318  //- Detect based on vector to nearest point differing for both endpoints
319  // v0 /= magV0;
320  // v1 /= magV1;
321  //
323  // if ((v0 & v1) < minCos)
324  //{
325  // return true;
326  //}
327  // else
328  //{
329  // return false;
330  //}
331 
332  //- Detect based on extrusion vector differing for both endpoints
333  // the idea is that e.g. a sawtooth wall can still be extruded
334  // successfully as long as it is done all to the same direction.
335  if ((pointWallDist[e[0]].v() & pointWallDist[e[1]].v()) < minCos)
336  {
337  return true;
338  }
339  else
340  {
341  return false;
342  }
343 }
344 
345 
346 void Foam::medialAxisMeshMover::update(const dictionary& coeffDict)
347 {
348  Info<< typeName
349  << " : Calculating distance to Medial Axis ..." << endl;
350 
351  const pointField& points = mesh().points();
352 
353  const indirectPrimitivePatch& pp = adaptPatchPtr_();
354  const labelList& meshPoints = pp.meshPoints();
355 
356 
357  // Read a few parameters
358  // ~~~~~~~~~~~~~~~~~~~~~
359 
360  //- Smooth surface normals
361  const label nSmoothSurfaceNormals = readLabel
362  (
363  coeffDict.lookup("nSmoothSurfaceNormals")
364  );
365 
366  //- When is medial axis
367  word angleKey = "minMedialAxisAngle";
368  if (!coeffDict.found(angleKey))
369  {
370  // Backwards compatibility
371  angleKey = "minMedianAxisAngle";
372  }
373  scalar minMedialAxisAngleCos = Foam::cos
374  (
375  degToRad(readScalar(coeffDict.lookup(angleKey)))
376  );
377 
378  //- Feature angle when to stop adding layers
379  const scalar featureAngle = readScalar(coeffDict.lookup("featureAngle"));
380 
381  //- When to slip along wall
382  const scalar slipFeatureAngle =
383  (
384  coeffDict.found("slipFeatureAngle")
385  ? readScalar(coeffDict.lookup("slipFeatureAngle"))
386  : 0.5*featureAngle
387  );
388 
389  //- Smooth internal normals
390  const label nSmoothNormals = readLabel
391  (
392  coeffDict.lookup("nSmoothNormals")
393  );
394 
395  //- Number of edges walking out
396  const label nMedialAxisIter = coeffDict.lookupOrDefault<label>
397  (
398  "nMedialAxisIter",
400  );
401 
402 
403  // Predetermine mesh edges
404  // ~~~~~~~~~~~~~~~~~~~~~~~
405 
406  // Precalulate (mesh) master point/edge (only relevant for shared pts/edges)
407  const PackedBoolList isMeshMasterPoint(syncTools::getMasterPoints(mesh()));
408  const PackedBoolList isMeshMasterEdge(syncTools::getMasterEdges(mesh()));
409  // Precalculate meshEdge per pp edge
410  const labelList meshEdges
411  (
412  pp.meshEdges
413  (
414  mesh().edges(),
415  mesh().pointEdges()
416  )
417  );
418 
419  // Precalulate (patch) master point/edge
420  const PackedBoolList isPatchMasterPoint
421  (
423  (
424  mesh(),
425  meshPoints
426  )
427  );
428  const PackedBoolList isPatchMasterEdge
429  (
431  (
432  mesh(),
433  meshEdges
434  )
435  );
436 
437  // Determine pointNormal
438  // ~~~~~~~~~~~~~~~~~~~~~
439 
440  pointField pointNormals(PatchTools::pointNormals(mesh(), pp));
441 
442  // Smooth patch normal vectors
443  smoothPatchNormals
444  (
445  nSmoothSurfaceNormals,
446  isPatchMasterPoint,
447  isPatchMasterEdge,
448  pointNormals
449  );
450 
451 
452  // Calculate distance to pp points
453  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
454 
455  // Distance to wall
456  List<pointData> pointWallDist(mesh().nPoints());
457 
458  // Dummy additional info for PointEdgeWave
459  int dummyTrackData = 0;
460 
461 
462  // 1. Calculate distance to points where displacement is specified.
463  {
464  // Seed data.
465  List<pointData> wallInfo(meshPoints.size());
466 
467  forAll(meshPoints, patchPointi)
468  {
469  label pointi = meshPoints[patchPointi];
470  wallInfo[patchPointi] = pointData
471  (
472  points[pointi],
473  0.0,
474  pointi, // passive scalar
475  pointNormals[patchPointi] // surface normals
476  );
477  }
478 
479  // Do all calculations
480  List<pointData> edgeWallDist(mesh().nEdges());
481  PointEdgeWave<pointData> wallDistCalc
482  (
483  mesh(),
484  meshPoints,
485  wallInfo,
486  pointWallDist,
487  edgeWallDist,
488  0, // max iterations
489  dummyTrackData
490  );
491  wallDistCalc.iterate(nMedialAxisIter);
492 
493  label nUnvisit = returnReduce
494  (
495  wallDistCalc.getUnsetPoints(),
496  sumOp<label>()
497  );
498 
499  if (nUnvisit > 0)
500  {
501  if (nMedialAxisIter > 0)
502  {
503  Info<< typeName
504  << " : Limited walk to " << nMedialAxisIter
505  << " steps. Not visited " << nUnvisit
506  << " out of " << mesh().globalData().nTotalPoints()
507  << " points" << endl;
508  }
509  else
510  {
512  << "Walking did not visit all points." << nl
513  << " Did not visit " << nUnvisit
514  << " out of " << mesh().globalData().nTotalPoints()
515  << " points. This is not necessarily a problem" << nl
516  << " and might be due to faceZones splitting of part"
517  << " of the domain." << nl << endl;
518  }
519  }
520  }
521 
522 
523  // 2. Find points with max distance and transport information back to
524  // wall.
525  {
526  List<pointData> pointMedialDist(mesh().nPoints());
527  List<pointData> edgeMedialDist(mesh().nEdges());
528 
529  // Seed point data.
530  DynamicList<pointData> maxInfo(meshPoints.size());
531  DynamicList<label> maxPoints(meshPoints.size());
532 
533  // 1. Medial axis points
534 
535  const edgeList& edges = mesh().edges();
536 
537  forAll(edges, edgeI)
538  {
539  const edge& e = edges[edgeI];
540 
541  if
542  (
543  !pointWallDist[e[0]].valid(dummyTrackData)
544  || !pointWallDist[e[1]].valid(dummyTrackData)
545  )
546  {
547  // Unvisited point. See above about nUnvisit warning
548  forAll(e, ep)
549  {
550  label pointi = e[ep];
551 
552  if (!pointMedialDist[pointi].valid(dummyTrackData))
553  {
554  maxPoints.append(pointi);
555  maxInfo.append
556  (
557  pointData
558  (
559  points[pointi],
560  0.0,
561  pointi, // passive data
562  Zero // passive data
563  )
564  );
565  pointMedialDist[pointi] = maxInfo.last();
566  }
567  }
568 
569  }
570  else if (isMaxEdge(pointWallDist, edgeI, minMedialAxisAngleCos))
571  {
572  // Both end points of edge have very different nearest wall
573  // point. Mark both points as medial axis points.
574 
575  // Approximate medial axis location on edge.
576  // const point medialAxisPt = e.centre(points);
577  vector eVec = e.vec(points);
578  scalar eMag = mag(eVec);
579  if (eMag > vSmall)
580  {
581  eVec /= eMag;
582 
583  // Calculate distance along edge
584  const point& p0 = points[e[0]];
585  const point& p1 = points[e[1]];
586  scalar dist0 = (p0-pointWallDist[e[0]].origin()) & eVec;
587  scalar dist1 = (pointWallDist[e[1]].origin()-p1) & eVec;
588  scalar s = 0.5*(dist1+eMag+dist0);
589 
590  point medialAxisPt;
591  if (s <= dist0)
592  {
593  medialAxisPt = p0;
594  }
595  else if (s >= dist0+eMag)
596  {
597  medialAxisPt = p1;
598  }
599  else
600  {
601  medialAxisPt = p0+(s-dist0)*eVec;
602  }
603 
604  forAll(e, ep)
605  {
606  label pointi = e[ep];
607 
608  if (!pointMedialDist[pointi].valid(dummyTrackData))
609  {
610  maxPoints.append(pointi);
611  maxInfo.append
612  (
613  pointData
614  (
615  medialAxisPt, // points[pointi],
616  magSqr(points[pointi]-medialAxisPt),//0.0,
617  pointi, // passive data
618  Zero // passive data
619  )
620  );
621  pointMedialDist[pointi] = maxInfo.last();
622  }
623  }
624  }
625  }
626  }
627 
628 
629  // 2. Seed non-adapt patches
630  const polyBoundaryMesh& patches = mesh().boundaryMesh();
631 
632  labelHashSet adaptPatches(adaptPatchIDs_);
633 
634 
635  forAll(patches, patchi)
636  {
637  const polyPatch& pp = patches[patchi];
638  const pointPatchVectorField& pvf =
639  pointDisplacement().boundaryField()[patchi];
640 
641  if
642  (
643  !pp.coupled()
644  && !isA<emptyPolyPatch>(pp)
645  && !adaptPatches.found(patchi)
646  )
647  {
648  const labelList& meshPoints = pp.meshPoints();
649 
650  // Check the type of the patchField. The types are
651  // - fixedValue (0 or more layers) but the >0 layers have
652  // already been handled in the adaptPatches loop
653  // - constraint (but not coupled) types, e.g. symmetryPlane,
654  // slip.
655  if (pvf.fixesValue())
656  {
657  // Disable all movement on fixedValue patchFields
658  Info<< typeName
659  << " : Inserting all points on patch " << pp.name()
660  << endl;
661 
662  forAll(meshPoints, i)
663  {
664  label pointi = meshPoints[i];
665  if (!pointMedialDist[pointi].valid(dummyTrackData))
666  {
667  maxPoints.append(pointi);
668  maxInfo.append
669  (
670  pointData
671  (
672  points[pointi],
673  0.0,
674  pointi, // passive data
675  Zero // passive data
676  )
677  );
678  pointMedialDist[pointi] = maxInfo.last();
679  }
680  }
681  }
682  else
683  {
684  // Based on geometry: analyse angle w.r.t. nearest moving
685  // point. In the pointWallDist we transported the
686  // normal as the passive vector. Note that this points
687  // out of the originating wall so inside of the domain
688  // on this patch.
689  Info<< typeName
690  << " : Inserting points on patch " << pp.name()
691  << " if angle to nearest layer patch > "
692  << slipFeatureAngle << " degrees." << endl;
693 
694  scalar slipFeatureAngleCos = Foam::cos
695  (
696  degToRad(slipFeatureAngle)
697  );
698  pointField pointNormals
699  (
701  );
702 
703  forAll(meshPoints, i)
704  {
705  label pointi = meshPoints[i];
706 
707  if
708  (
709  pointWallDist[pointi].valid(dummyTrackData)
710  && !pointMedialDist[pointi].valid(dummyTrackData)
711  )
712  {
713  // Check if angle not too large.
714  scalar cosAngle =
715  (
716  -pointWallDist[pointi].v()
717  & pointNormals[i]
718  );
719  if (cosAngle > slipFeatureAngleCos)
720  {
721  // Extrusion direction practically perpendicular
722  // to the patch. Disable movement at the patch.
723 
724  maxPoints.append(pointi);
725  maxInfo.append
726  (
727  pointData
728  (
729  points[pointi],
730  0.0,
731  pointi, // passive data
732  Zero // passive data
733  )
734  );
735  pointMedialDist[pointi] = maxInfo.last();
736  }
737  else
738  {
739  // Extrusion direction makes angle with patch
740  // so allow sliding - don't insert zero points
741  }
742  }
743  }
744  }
745  }
746  }
747 
748  maxInfo.shrink();
749  maxPoints.shrink();
750 
751  // Do all calculations
752  PointEdgeWave<pointData> medialDistCalc
753  (
754  mesh(),
755  maxPoints,
756  maxInfo,
757 
758  pointMedialDist,
759  edgeMedialDist,
760  0, // max iterations
761  dummyTrackData
762  );
763  medialDistCalc.iterate(2*nMedialAxisIter);
764 
765 
766  // Extract medial axis distance as pointScalarField
767  forAll(pointMedialDist, pointi)
768  {
769  if (pointMedialDist[pointi].valid(dummyTrackData))
770  {
771  medialDist_[pointi] = Foam::sqrt
772  (
773  pointMedialDist[pointi].distSqr()
774  );
775  medialVec_[pointi] = pointMedialDist[pointi].origin();
776  }
777  else
778  {
779  // Unvisited. Do as if on medial axis so unmoving
780  medialDist_[pointi] = 0.0;
781  medialVec_[pointi] = point(1, 0, 0);
782  }
783  }
784  }
785 
786  // Extract transported surface normals as pointVectorField
787  forAll(dispVec_, i)
788  {
789  if (!pointWallDist[i].valid(dummyTrackData))
790  {
791  dispVec_[i] = vector(1, 0, 0);
792  }
793  else
794  {
795  dispVec_[i] = pointWallDist[i].v();
796  }
797  }
798 
799  // Smooth normal vectors. Do not change normals on pp.meshPoints
800  smoothNormals
801  (
802  nSmoothNormals,
803  isMeshMasterPoint,
804  isMeshMasterEdge,
805  meshPoints,
806  dispVec_
807  );
808 
809  // Calculate ratio point medial distance to point wall distance
810  forAll(medialRatio_, pointi)
811  {
812  if (!pointWallDist[pointi].valid(dummyTrackData))
813  {
814  medialRatio_[pointi] = 0.0;
815  }
816  else
817  {
818  scalar wDist2 = pointWallDist[pointi].distSqr();
819  scalar mDist = medialDist_[pointi];
820 
821  if (wDist2 < sqr(small) && mDist < small)
822  //- Note: maybe less strict:
823  //(
824  // wDist2 < sqr(meshRefiner_.mergeDistance())
825  // && mDist < meshRefiner_.mergeDistance()
826  //)
827  {
828  medialRatio_[pointi] = 0.0;
829  }
830  else
831  {
832  medialRatio_[pointi] = mDist / (Foam::sqrt(wDist2) + mDist);
833  }
834  }
835  }
836 
837 
838  if (debug)
839  {
840  Info<< typeName
841  << " : Writing medial axis fields:" << nl
842  << incrIndent
843  << "ratio of medial distance to wall distance : "
844  << medialRatio_.name() << nl
845  << "distance to nearest medial axis : "
846  << medialDist_.name() << nl
847  << "nearest medial axis location : "
848  << medialVec_.name() << nl
849  << "normal at nearest wall : "
850  << dispVec_.name() << nl
851  << decrIndent << nl
852  << endl;
853 
854  dispVec_.write();
855  medialRatio_.write();
856  medialDist_.write();
857  medialVec_.write();
858  }
859 }
860 
861 
862 bool Foam::medialAxisMeshMover::unmarkExtrusion
863 (
864  const label patchPointi,
865  pointField& patchDisp,
866  List<snappyLayerDriver::extrudeMode>& extrudeStatus
867 )
868 {
869  if (extrudeStatus[patchPointi] == snappyLayerDriver::EXTRUDE)
870  {
871  extrudeStatus[patchPointi] = snappyLayerDriver::NOEXTRUDE;
872  patchDisp[patchPointi] = Zero;
873  return true;
874  }
875  else if (extrudeStatus[patchPointi] == snappyLayerDriver::EXTRUDEREMOVE)
876  {
877  extrudeStatus[patchPointi] = snappyLayerDriver::NOEXTRUDE;
878  patchDisp[patchPointi] = Zero;
879  return true;
880  }
881  else
882  {
883  return false;
884  }
885 }
886 
887 
888 void Foam::medialAxisMeshMover::syncPatchDisplacement
889 (
890  const scalarField& minThickness,
891  pointField& patchDisp,
892  List<snappyLayerDriver::extrudeMode>& extrudeStatus
893 ) const
894 {
895  const indirectPrimitivePatch& pp = adaptPatchPtr_();
896  const labelList& meshPoints = pp.meshPoints();
897 
898  label nChangedTotal = 0;
899 
900  while (true)
901  {
902  label nChanged = 0;
903 
904  // Sync displacement (by taking min)
906  (
907  mesh(),
908  meshPoints,
909  patchDisp,
910  minMagSqrEqOp<vector>(),
911  point::rootMax // null value
912  );
913 
914  // Unmark if displacement too small
915  forAll(patchDisp, i)
916  {
917  if (mag(patchDisp[i]) < minThickness[i])
918  {
919  if (unmarkExtrusion(i, patchDisp, extrudeStatus))
920  {
921  nChanged++;
922  }
923  }
924  }
925 
926  // labelList syncPatchNLayers(patchNLayers);
927  //
928  // syncTools::syncPointList
929  //(
930  // mesh(),
931  // meshPoints,
932  // syncPatchNLayers,
933  // minEqOp<label>(),
934  // labelMax // null value
935  //);
936  //
939  // forAll(syncPatchNLayers, i)
940  //{
941  // if (syncPatchNLayers[i] != patchNLayers[i])
942  // {
943  // if
944  // (
945  // unmarkExtrusion
946  // (
947  // i,
948  // patchDisp,
949  // patchNLayers,
950  // extrudeStatus
951  // )
952  // )
953  // {
954  // nChanged++;
955  // }
956  // }
957  //}
958  //
959  // syncTools::syncPointList
960  //(
961  // mesh(),
962  // meshPoints,
963  // syncPatchNLayers,
964  // maxEqOp<label>(),
965  // labelMin // null value
966  //);
967  //
970  // forAll(syncPatchNLayers, i)
971  //{
972  // if (syncPatchNLayers[i] != patchNLayers[i])
973  // {
974  // if
975  // (
976  // unmarkExtrusion
977  // (
978  // i,
979  // patchDisp,
980  // patchNLayers,
981  // extrudeStatus
982  // )
983  // )
984  // {
985  // nChanged++;
986  // }
987  // }
988  //}
989 
990  nChangedTotal += nChanged;
991 
992  if (!returnReduce(nChanged, sumOp<label>()))
993  {
994  break;
995  }
996  }
997 
998  // Info<< "Prevented extrusion on "
999  // << returnReduce(nChangedTotal, sumOp<label>())
1000  // << " coupled patch points during syncPatchDisplacement." << endl;
1001 }
1002 
1003 
1004 void Foam::medialAxisMeshMover::minSmoothField
1005 (
1006  const label nSmoothDisp,
1007  const PackedBoolList& isPatchMasterPoint,
1008  const PackedBoolList& isPatchMasterEdge,
1009  const scalarField& fieldMin,
1010  scalarField& field
1011 ) const
1012 {
1013  const indirectPrimitivePatch& pp = adaptPatchPtr_();
1014  const edgeList& edges = pp.edges();
1015  const labelList& meshPoints = pp.meshPoints();
1016 
1017  scalarField edgeWeights(edges.size());
1018  scalarField invSumWeight(meshPoints.size());
1020  (
1021  mesh(),
1022  isPatchMasterEdge,
1023  meshPoints,
1024  edges,
1025  edgeWeights,
1026  invSumWeight
1027  );
1028 
1029  // Get smoothly varying patch field.
1030  Info<< typeName << " : Smoothing field ..." << endl;
1031 
1032  for (label iter = 0; iter < nSmoothDisp; iter++)
1033  {
1034  scalarField average(pp.nPoints());
1036  (
1037  mesh(),
1038  isPatchMasterEdge,
1039  meshPoints,
1040  edges,
1041  edgeWeights,
1042  field,
1043  average
1044  );
1045  average *= invSumWeight;
1046 
1047  // Transfer to field
1048  forAll(field, pointi)
1049  {
1050  // full smoothing neighbours + point value
1051  average[pointi] = 0.5*(field[pointi]+average[pointi]);
1052 
1053  // perform monotonic smoothing
1054  if
1055  (
1056  average[pointi] < field[pointi]
1057  && average[pointi] >= fieldMin[pointi]
1058  )
1059  {
1060  field[pointi] = average[pointi];
1061  }
1062  }
1063 
1064  // Do residual calculation every so often.
1065  if ((iter % 10) == 0)
1066  {
1067  scalar resid = meshRefinement::gAverage
1068  (
1069  isPatchMasterPoint,
1070  mag(field-average)()
1071  );
1072  Info<< " Iteration " << iter << " residual " << resid << endl;
1073  }
1074  }
1075 }
1076 
1077 // Stop layer growth where mesh wraps around edge with a
1078 // large feature angle
1079 void Foam::medialAxisMeshMover::
1080 handleFeatureAngleLayerTerminations
1081 (
1082  const scalar minCos,
1083  const PackedBoolList& isPatchMasterPoint,
1084  const labelList& meshEdges,
1085  List<snappyLayerDriver::extrudeMode>& extrudeStatus,
1086  pointField& patchDisp,
1087  label& nPointCounter
1088 ) const
1089 {
1090  const indirectPrimitivePatch& pp = adaptPatchPtr_();
1091 
1092  // Mark faces that have all points extruded
1093  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1094 
1095  boolList extrudedFaces(pp.size(), true);
1096 
1097  forAll(pp.localFaces(), facei)
1098  {
1099  const face& f = pp.localFaces()[facei];
1100 
1101  forAll(f, fp)
1102  {
1103  if (extrudeStatus[f[fp]] == snappyLayerDriver::NOEXTRUDE)
1104  {
1105  extrudedFaces[facei] = false;
1106  break;
1107  }
1108  }
1109  }
1110 
1111 
1112 
1113  // label nOldPointCounter = nPointCounter;
1114 
1115  // Detect situation where two featureedge-neighbouring faces are partly or
1116  // not extruded and the edge itself is extruded. In this case unmark the
1117  // edge for extrusion.
1118 
1119 
1120  List<List<point>> edgeFaceNormals(pp.nEdges());
1121  List<List<bool>> edgeFaceExtrude(pp.nEdges());
1122 
1123  const labelListList& edgeFaces = pp.edgeFaces();
1124  const vectorField& faceNormals = pp.faceNormals();
1125 
1126  forAll(edgeFaces, edgeI)
1127  {
1128  const labelList& eFaces = edgeFaces[edgeI];
1129 
1130  edgeFaceNormals[edgeI].setSize(eFaces.size());
1131  edgeFaceExtrude[edgeI].setSize(eFaces.size());
1132  forAll(eFaces, i)
1133  {
1134  label facei = eFaces[i];
1135  edgeFaceNormals[edgeI][i] = faceNormals[facei];
1136  edgeFaceExtrude[edgeI][i] = extrudedFaces[facei];
1137  }
1138  }
1139 
1141  (
1142  mesh(),
1143  meshEdges,
1144  edgeFaceNormals,
1145  globalMeshData::ListPlusEqOp<List<point>>(), // combine operator
1146  List<point>() // null value
1147  );
1148 
1150  (
1151  mesh(),
1152  meshEdges,
1153  edgeFaceExtrude,
1154  globalMeshData::ListPlusEqOp<List<bool>>(), // combine operator
1155  List<bool>() // null value
1156  );
1157 
1158 
1159  forAll(edgeFaceNormals, edgeI)
1160  {
1161  const List<point>& eFaceNormals = edgeFaceNormals[edgeI];
1162  const List<bool>& eFaceExtrude = edgeFaceExtrude[edgeI];
1163 
1164  if (eFaceNormals.size() == 2)
1165  {
1166  const edge& e = pp.edges()[edgeI];
1167  label v0 = e[0];
1168  label v1 = e[1];
1169 
1170  if
1171  (
1172  extrudeStatus[v0] != snappyLayerDriver::NOEXTRUDE
1173  || extrudeStatus[v1] != snappyLayerDriver::NOEXTRUDE
1174  )
1175  {
1176  if (!eFaceExtrude[0] || !eFaceExtrude[1])
1177  {
1178  const vector& n0 = eFaceNormals[0];
1179  const vector& n1 = eFaceNormals[1];
1180 
1181  if ((n0 & n1) < minCos)
1182  {
1183  if (unmarkExtrusion(v0, patchDisp, extrudeStatus))
1184  {
1185  if (isPatchMasterPoint[v0])
1186  {
1187  nPointCounter++;
1188  }
1189  }
1190  if (unmarkExtrusion(v1, patchDisp, extrudeStatus))
1191  {
1192  if (isPatchMasterPoint[v1])
1193  {
1194  nPointCounter++;
1195  }
1196  }
1197  }
1198  }
1199  }
1200  }
1201  }
1202 
1203  // Info<< "Added "
1204  // << returnReduce(nPointCounter-nOldPointCounter, sumOp<label>())
1205  // << " point not to extrude." << endl;
1206 }
1207 
1208 
1209 // Find isolated islands (points, edges and faces and layer terminations)
1210 // in the layer mesh and stop any layer growth at these points.
1211 void Foam::medialAxisMeshMover::findIsolatedRegions
1212 (
1213  const scalar minCosLayerTermination,
1214  const bool detectExtrusionIsland,
1215  const PackedBoolList& isPatchMasterPoint,
1216  const PackedBoolList& isPatchMasterEdge,
1217  const labelList& meshEdges,
1218  const scalarField& minThickness,
1219  List<snappyLayerDriver::extrudeMode>& extrudeStatus,
1220  pointField& patchDisp
1221 ) const
1222 {
1223  const indirectPrimitivePatch& pp = adaptPatchPtr_();
1224  const labelListList& pointFaces = pp.pointFaces();
1225  const labelList& meshPoints = pp.meshPoints();
1226 
1227  Info<< typeName << " : Removing isolated regions ..." << endl;
1228 
1229  // Keep count of number of points unextruded
1230  label nPointCounter = 0;
1231 
1232 
1233  autoPtr<OBJstream> str;
1234  if (debug)
1235  {
1236  str.reset
1237  (
1238  new OBJstream
1239  (
1240  mesh().time().path()
1241  / "islandExcludePoints_"
1242  + mesh().time().timeName()
1243  + ".obj"
1244  )
1245  );
1246  Info<< typeName
1247  << " : Writing points surrounded by non-extruded points to "
1248  << str().name() << endl;
1249  }
1250 
1251  while (true)
1252  {
1253  // Stop layer growth where mesh wraps around edge with a
1254  // large feature angle
1255  handleFeatureAngleLayerTerminations
1256  (
1257  minCosLayerTermination,
1258  isPatchMasterPoint,
1259  meshEdges,
1260 
1261  extrudeStatus,
1262  patchDisp,
1263  nPointCounter
1264  );
1265 
1266  syncPatchDisplacement(minThickness, patchDisp, extrudeStatus);
1267 
1268 
1269 
1270  // Detect either:
1271  // - point where all surrounding points are not extruded
1272  // (detectExtrusionIsland)
1273  // or
1274  // - point where all the faces surrounding it are not fully
1275  // extruded
1276 
1277  boolList keptPoints(pp.nPoints(), false);
1278 
1279  if (detectExtrusionIsland)
1280  {
1281  // Do not extrude from point where all neighbouring
1282  // points are not grown
1283  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1284 
1285  labelList islandPoint(pp.size(), -1);
1286  forAll(pp, facei)
1287  {
1288  const face& f = pp.localFaces()[facei];
1289 
1290  forAll(f, fp)
1291  {
1292  if (extrudeStatus[f[fp]] != snappyLayerDriver::NOEXTRUDE)
1293  {
1294  if (islandPoint[facei] == -1)
1295  {
1296  // First point to extrude
1297  islandPoint[facei] = f[fp];
1298  }
1299  else if (islandPoint[facei] != -2)
1300  {
1301  // Second or more point to extrude
1302  islandPoint[facei] = -2;
1303  }
1304  }
1305  }
1306  }
1307 
1308  // islandPoint:
1309  // -1 : no point extruded on face
1310  // -2 : >= 2 points extruded on face
1311  // >=0: label of point extruded
1312 
1313  // Check all surrounding faces that I am the islandPoint
1314  forAll(pointFaces, patchPointi)
1315  {
1316  if (extrudeStatus[patchPointi] != snappyLayerDriver::NOEXTRUDE)
1317  {
1318  const labelList& pFaces = pointFaces[patchPointi];
1319 
1320  forAll(pFaces, i)
1321  {
1322  label facei = pFaces[i];
1323  if (islandPoint[facei] != patchPointi)
1324  {
1325  keptPoints[patchPointi] = true;
1326  break;
1327  }
1328  }
1329  }
1330  }
1331  }
1332  else
1333  {
1334  // Do not extrude from point where all neighbouring
1335  // faces are not grown
1336  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1337 
1338  boolList extrudedFaces(pp.size(), true);
1339  forAll(pp.localFaces(), facei)
1340  {
1341  const face& f = pp.localFaces()[facei];
1342  forAll(f, fp)
1343  {
1344  if (extrudeStatus[f[fp]] == snappyLayerDriver::NOEXTRUDE)
1345  {
1346  extrudedFaces[facei] = false;
1347  break;
1348  }
1349  }
1350  }
1351 
1352  const labelListList& pointFaces = pp.pointFaces();
1353 
1354  forAll(keptPoints, patchPointi)
1355  {
1356  const labelList& pFaces = pointFaces[patchPointi];
1357 
1358  forAll(pFaces, i)
1359  {
1360  label facei = pFaces[i];
1361  if (extrudedFaces[facei])
1362  {
1363  keptPoints[patchPointi] = true;
1364  break;
1365  }
1366  }
1367  }
1368  }
1369 
1370 
1372  (
1373  mesh(),
1374  meshPoints,
1375  keptPoints,
1376  orEqOp<bool>(),
1377  false // null value
1378  );
1379 
1380  label nChanged = 0;
1381 
1382  forAll(keptPoints, patchPointi)
1383  {
1384  if (!keptPoints[patchPointi])
1385  {
1386  if (unmarkExtrusion(patchPointi, patchDisp, extrudeStatus))
1387  {
1388  nPointCounter++;
1389  nChanged++;
1390 
1391  if (str.valid())
1392  {
1393  str().write(pp.points()[meshPoints[patchPointi]]);
1394  }
1395  }
1396  }
1397  }
1398 
1399 
1400  if (returnReduce(nChanged, sumOp<label>()) == 0)
1401  {
1402  break;
1403  }
1404  }
1405 
1406  const edgeList& edges = pp.edges();
1407 
1408 
1409  // Count number of mesh edges using a point
1410  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1411 
1412  labelList isolatedPoint(pp.nPoints(),0);
1413 
1414  forAll(edges, edgeI)
1415  {
1416  if (isPatchMasterEdge[edgeI])
1417  {
1418  const edge& e = edges[edgeI];
1419 
1420  label v0 = e[0];
1421  label v1 = e[1];
1422 
1423  if (extrudeStatus[v1] != snappyLayerDriver::NOEXTRUDE)
1424  {
1425  isolatedPoint[v0] += 1;
1426  }
1427  if (extrudeStatus[v0] != snappyLayerDriver::NOEXTRUDE)
1428  {
1429  isolatedPoint[v1] += 1;
1430  }
1431  }
1432  }
1433 
1435  (
1436  mesh(),
1437  meshPoints,
1438  isolatedPoint,
1439  plusEqOp<label>(),
1440  label(0) // null value
1441  );
1442 
1443  // stop layer growth on isolated faces
1444  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1445  forAll(pp, facei)
1446  {
1447  const face& f = pp.localFaces()[facei];
1448  bool failed = false;
1449  forAll(f, fp)
1450  {
1451  if (isolatedPoint[f[fp]] > 2)
1452  {
1453  failed = true;
1454  break;
1455  }
1456  }
1457  bool allPointsExtruded = true;
1458  if (!failed)
1459  {
1460  forAll(f, fp)
1461  {
1462  if (extrudeStatus[f[fp]] == snappyLayerDriver::NOEXTRUDE)
1463  {
1464  allPointsExtruded = false;
1465  break;
1466  }
1467  }
1468 
1469  if (allPointsExtruded)
1470  {
1471  forAll(f, fp)
1472  {
1473  if
1474  (
1475  unmarkExtrusion
1476  (
1477  f[fp],
1478  patchDisp,
1479  extrudeStatus
1480  )
1481  )
1482  {
1483  nPointCounter++;
1484 
1485  if (str.valid())
1486  {
1487  str().write(pp.points()[meshPoints[f[fp]]]);
1488  }
1489  }
1490  }
1491  }
1492  }
1493  }
1494 
1495  reduce(nPointCounter, sumOp<label>());
1496  Info<< typeName
1497  << " : Number of isolated points extrusion stopped : "<< nPointCounter
1498  << endl;
1499 }
1500 
1501 
1502 void Foam::medialAxisMeshMover::smoothLambdaMuDisplacement
1503 (
1504  const label nSmoothDisp,
1505  const PackedBoolList& isMeshMasterPoint,
1506  const PackedBoolList& isMeshMasterEdge,
1507  vectorField& displacement
1508 ) const
1509 {
1510  const edgeList& edges = mesh().edges();
1511 
1512  // Correspondence between local edges/points and mesh edges/points
1513  const labelList meshPoints(identity(mesh().nPoints()));
1514 
1515  // Calculate inverse sum of weights
1516  scalarField edgeWeights(mesh().nEdges());
1517  scalarField invSumWeight(meshPoints.size());
1519  (
1520  mesh(),
1521  isMeshMasterEdge,
1522  meshPoints,
1523  edges,
1524  edgeWeights,
1525  invSumWeight
1526  );
1527 
1528  // Get smoothly varying patch field.
1529  Info<< typeName << " : Smoothing displacement ..." << endl;
1530 
1531  const scalar lambda = 0.33;
1532  const scalar mu = -0.34;
1533 
1535 
1536  for (label iter = 0; iter < nSmoothDisp; iter++)
1537  {
1539  (
1540  mesh(),
1541  isMeshMasterEdge,
1542  meshPoints,
1543  edges,
1544  edgeWeights,
1545  displacement,
1546  average
1547  );
1548  average *= invSumWeight;
1549 
1550  forAll(displacement, i)
1551  {
1552  if (medialRatio_[i] > small && medialRatio_[i] < 1-small)
1553  {
1554  displacement[i] = (1-lambda)*displacement[i]+lambda*average[i];
1555  }
1556  }
1557 
1559  (
1560  mesh(),
1561  isMeshMasterEdge,
1562  meshPoints,
1563  edges,
1564  edgeWeights,
1565  displacement,
1566  average
1567  );
1568  average *= invSumWeight;
1569 
1570 
1571  forAll(displacement, i)
1572  {
1573  if (medialRatio_[i] > small && medialRatio_[i] < 1-small)
1574  {
1575  displacement[i] = (1-mu)*displacement[i]+mu*average[i];
1576  }
1577  }
1578 
1579 
1580  // Do residual calculation every so often.
1581  if ((iter % 10) == 0)
1582  {
1583  scalar resid = meshRefinement::gAverage
1584  (
1585  isMeshMasterPoint,
1586  mag(displacement-average)()
1587  );
1588  Info<< " Iteration " << iter << " residual " << resid << endl;
1589  }
1590  }
1591 }
1592 
1593 
1594 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
1595 
1598  const dictionary& dict,
1599  const List<labelPair>& baffles,
1600  pointVectorField& pointDisplacement
1601 )
1602 :
1603  externalDisplacementMeshMover(dict, baffles, pointDisplacement),
1604  adaptPatchIDs_(getFixedValueBCs(pointDisplacement)),
1605  adaptPatchPtr_(getPatch(mesh(), adaptPatchIDs_)),
1606  scale_
1607  (
1608  IOobject
1609  (
1610  "scale",
1611  pointDisplacement.time().timeName(),
1612  pointDisplacement.db(),
1615  ),
1616  pMesh(),
1618  ),
1619  oldPoints_(mesh().points()),
1620  meshMover_
1621  (
1622  const_cast<polyMesh&>(mesh()),
1623  const_cast<pointMesh&>(pMesh()),
1624  adaptPatchPtr_(),
1625  pointDisplacement,
1626  scale_,
1627  oldPoints_,
1628  adaptPatchIDs_,
1629  dict
1630  ),
1631  dispVec_
1632  (
1633  IOobject
1634  (
1635  "dispVec",
1636  pointDisplacement.time().timeName(),
1637  pointDisplacement.db(),
1640  false
1641  ),
1642  pMesh(),
1644  ),
1645  medialRatio_
1646  (
1647  IOobject
1648  (
1649  "medialRatio",
1650  pointDisplacement.time().timeName(),
1651  pointDisplacement.db(),
1654  false
1655  ),
1656  pMesh(),
1658  ),
1659  medialDist_
1660  (
1661  IOobject
1662  (
1663  "pointMedialDist",
1664  pointDisplacement.time().timeName(),
1665  pointDisplacement.db(),
1668  false
1669  ),
1670  pMesh(),
1672  ),
1673  medialVec_
1674  (
1675  IOobject
1676  (
1677  "medialVec",
1678  pointDisplacement.time().timeName(),
1679  pointDisplacement.db(),
1682  false
1683  ),
1684  pMesh(),
1686  )
1687 {
1688  update(dict);
1689 }
1690 
1691 
1692 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
1693 
1695 {}
1696 
1697 
1698 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1699 
1700 void Foam::medialAxisMeshMover::calculateDisplacement
1701 (
1702  const dictionary& coeffDict,
1703  const scalarField& minThickness,
1704  List<snappyLayerDriver::extrudeMode>& extrudeStatus,
1705  pointField& patchDisp
1706 )
1707 {
1708  Info<< typeName << " : Smoothing using Medial Axis ..." << endl;
1709 
1710  const indirectPrimitivePatch& pp = adaptPatchPtr_;
1711  const labelList& meshPoints = pp.meshPoints();
1712 
1713 
1714  // Read settings
1715  // ~~~~~~~~~~~~~
1716 
1717  //- (lambda-mu) smoothing of internal displacement
1718  const label nSmoothDisplacement = coeffDict.lookupOrDefault
1719  (
1720  "nSmoothDisplacement",
1721  0
1722  );
1723 
1724  //- Layer thickness too big
1725  const scalar maxThicknessToMedialRatio = readScalar
1726  (
1727  coeffDict.lookup("maxThicknessToMedialRatio")
1728  );
1729 
1730  //- Feature angle when to stop adding layers
1731  const scalar featureAngle = readScalar(coeffDict.lookup("featureAngle"));
1732 
1733  //- Stop layer growth where mesh wraps around sharp edge
1734  const scalar minCosLayerTermination = Foam::cos
1735  (
1736  degToRad(0.5*featureAngle)
1737  );
1738 
1739  //- Smoothing wanted patch thickness
1740  const label nSmoothPatchThickness = readLabel
1741  (
1742  coeffDict.lookup("nSmoothThickness")
1743  );
1744 
1745  //- Number of edges walking out
1746  const label nMedialAxisIter = coeffDict.lookupOrDefault<label>
1747  (
1748  "nMedialAxisIter",
1750  );
1751 
1752  //- Use strict extrusionIsland detection
1753  const Switch detectExtrusionIsland = coeffDict.lookupOrDefault<Switch>
1754  (
1755  "detectExtrusionIsland",
1756  false
1757  );
1758 
1759 
1760  // Precalulate master points/edge (only relevant for shared points/edges)
1761  const PackedBoolList isMeshMasterPoint(syncTools::getMasterPoints(mesh()));
1762  const PackedBoolList isMeshMasterEdge(syncTools::getMasterEdges(mesh()));
1763  // Precalculate meshEdge per pp edge
1764  const labelList meshEdges
1765  (
1766  pp.meshEdges
1767  (
1768  mesh().edges(),
1769  mesh().pointEdges()
1770  )
1771  );
1772 
1773  // Precalulate (patch) master point/edge
1774  const PackedBoolList isPatchMasterPoint
1775  (
1777  (
1778  mesh(),
1779  meshPoints
1780  )
1781  );
1782  const PackedBoolList isPatchMasterEdge
1783  (
1785  (
1786  mesh(),
1787  meshEdges
1788  )
1789  );
1790 
1791 
1792  scalarField thickness(patchDisp.size());
1793 
1794  thickness = mag(patchDisp);
1795 
1796  forAll(thickness, patchPointi)
1797  {
1798  if (extrudeStatus[patchPointi] == snappyLayerDriver::NOEXTRUDE)
1799  {
1800  thickness[patchPointi] = 0.0;
1801  }
1802  }
1803 
1804  label numThicknessRatioExclude = 0;
1805 
1806  // reduce thickness where thickness/medial axis distance large
1807  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1808 
1809  autoPtr<OBJstream> str;
1810  if (debug)
1811  {
1812  str.reset
1813  (
1814  new OBJstream
1815  (
1816  mesh().time().path()
1817  / "thicknessRatioExcludePoints_"
1818  + mesh().time().timeName()
1819  + ".obj"
1820  )
1821  );
1822  Info<< typeName
1823  << " : Writing points with too large an extrusion distance to "
1824  << str().name() << endl;
1825  }
1826 
1827  autoPtr<OBJstream> medialVecStr;
1828  if (debug)
1829  {
1830  medialVecStr.reset
1831  (
1832  new OBJstream
1833  (
1834  mesh().time().path()
1835  / "thicknessRatioExcludeMedialVec_"
1836  + mesh().time().timeName()
1837  + ".obj"
1838  )
1839  );
1840  Info<< typeName
1841  << " : Writing medial axis vectors on points with too large"
1842  << " an extrusion distance to " << medialVecStr().name() << endl;
1843  }
1844 
1845  forAll(meshPoints, patchPointi)
1846  {
1847  if (extrudeStatus[patchPointi] != snappyLayerDriver::NOEXTRUDE)
1848  {
1849  label pointi = meshPoints[patchPointi];
1850 
1851  //- Option 1: look only at extrusion thickness v.s. distance
1852  // to nearest (medial axis or static) point.
1853  scalar mDist = medialDist_[pointi];
1854  scalar thicknessRatio = thickness[patchPointi]/(mDist+vSmall);
1855 
1856  //- Option 2: Look at component in the direction
1857  // of nearest (medial axis or static) point.
1858  vector n =
1859  patchDisp[patchPointi]
1860  / (mag(patchDisp[patchPointi]) + vSmall);
1861  vector mVec = mesh().points()[pointi]-medialVec_[pointi];
1862  mVec /= mag(mVec)+vSmall;
1863  thicknessRatio *= (n&mVec);
1864 
1865  if (thicknessRatio > maxThicknessToMedialRatio)
1866  {
1867  // Truncate thickness.
1868  if (debug&2)
1869  {
1870  Pout<< "truncating displacement at "
1871  << mesh().points()[pointi]
1872  << " from " << thickness[patchPointi]
1873  << " to "
1874  << 0.5
1875  *(
1876  minThickness[patchPointi]
1877  +thickness[patchPointi]
1878  )
1879  << " medial direction:" << mVec
1880  << " extrusion direction:" << n
1881  << " with thicknessRatio:" << thicknessRatio
1882  << endl;
1883  }
1884 
1885  thickness[patchPointi] =
1886  0.5*(minThickness[patchPointi]+thickness[patchPointi]);
1887 
1888  patchDisp[patchPointi] = thickness[patchPointi]*n;
1889 
1890  if (isPatchMasterPoint[patchPointi])
1891  {
1892  numThicknessRatioExclude++;
1893  }
1894 
1895  if (str.valid())
1896  {
1897  const point& pt = mesh().points()[pointi];
1898  str().write(linePointRef(pt, pt+patchDisp[patchPointi]));
1899  }
1900  if (medialVecStr.valid())
1901  {
1902  const point& pt = mesh().points()[pointi];
1903  medialVecStr().write
1904  (
1905  linePointRef
1906  (
1907  pt,
1908  medialVec_[pointi]
1909  )
1910  );
1911  }
1912  }
1913  }
1914  }
1915 
1916  reduce(numThicknessRatioExclude, sumOp<label>());
1917  Info<< typeName << " : Reducing layer thickness at "
1918  << numThicknessRatioExclude
1919  << " nodes where thickness to medial axis distance is large " << endl;
1920 
1921 
1922  // find points where layer growth isolated to a lone point, edge or face
1923 
1924  findIsolatedRegions
1925  (
1926  minCosLayerTermination,
1927  detectExtrusionIsland,
1928 
1929  isPatchMasterPoint,
1930  isPatchMasterEdge,
1931  meshEdges,
1932  minThickness,
1933 
1934  extrudeStatus,
1935  patchDisp
1936  );
1937 
1938  // Update thickess for changed extrusion
1939  forAll(thickness, patchPointi)
1940  {
1941  if (extrudeStatus[patchPointi] == snappyLayerDriver::NOEXTRUDE)
1942  {
1943  thickness[patchPointi] = 0.0;
1944  }
1945  }
1946 
1947 
1948  // smooth layer thickness on moving patch
1949  minSmoothField
1950  (
1951  nSmoothPatchThickness,
1952  isPatchMasterPoint,
1953  isPatchMasterEdge,
1954  minThickness,
1955 
1956  thickness
1957  );
1958 
1959 
1960  // Dummy additional info for PointEdgeWave
1961  int dummyTrackData = 0;
1962 
1963  List<pointData> pointWallDist(mesh().nPoints());
1964 
1965  const pointField& points = mesh().points();
1966  // 1. Calculate distance to points where displacement is specified.
1967  // This wave is used to transport layer thickness
1968  {
1969  // Distance to wall and medial axis on edges.
1970  List<pointData> edgeWallDist(mesh().nEdges());
1971  labelList wallPoints(meshPoints.size());
1972 
1973  // Seed data.
1974  List<pointData> wallInfo(meshPoints.size());
1975 
1976  forAll(meshPoints, patchPointi)
1977  {
1978  label pointi = meshPoints[patchPointi];
1979  wallPoints[patchPointi] = pointi;
1980  wallInfo[patchPointi] = pointData
1981  (
1982  points[pointi],
1983  0.0,
1984  thickness[patchPointi], // transport layer thickness
1985  Zero // passive vector
1986  );
1987  }
1988 
1989  // Do all calculations
1990  PointEdgeWave<pointData> wallDistCalc
1991  (
1992  mesh(),
1993  wallPoints,
1994  wallInfo,
1995  pointWallDist,
1996  edgeWallDist,
1997  0, // max iterations
1998  dummyTrackData
1999  );
2000  wallDistCalc.iterate(nMedialAxisIter);
2001  }
2002 
2003 
2004  // Calculate scaled displacement vector
2005  pointField& displacement = pointDisplacement_;
2006 
2007  forAll(displacement, pointi)
2008  {
2009  if (!pointWallDist[pointi].valid(dummyTrackData))
2010  {
2011  displacement[pointi] = Zero;
2012  }
2013  else
2014  {
2015  // 1) displacement on nearest wall point, scaled by medialRatio
2016  // (wall distance / medial distance)
2017  // 2) pointWallDist[pointi].s() is layer thickness transported
2018  // from closest wall point.
2019  // 3) shrink in opposite direction of addedPoints
2020  displacement[pointi] =
2021  -medialRatio_[pointi]
2022  *pointWallDist[pointi].s()
2023  *dispVec_[pointi];
2024  }
2025  }
2026 
2027 
2028  // Smear displacement away from fixed values (medialRatio=0 or 1)
2029  if (nSmoothDisplacement > 0)
2030  {
2031  smoothLambdaMuDisplacement
2032  (
2033  nSmoothDisplacement,
2034  isMeshMasterPoint,
2035  isMeshMasterEdge,
2036  displacement
2037  );
2038  }
2039 }
2040 
2041 
2042 bool Foam::medialAxisMeshMover::shrinkMesh
2043 (
2044  const dictionary& meshQualityDict,
2045  const label nAllowableErrors,
2046  labelList& checkFaces
2047 )
2048 {
2049  //- Number of attempts shrinking the mesh
2050  const label nSnap = readLabel(meshQualityDict.lookup("nRelaxIter"));
2051 
2052 
2053 
2054 
2055  // Make sure displacement boundary conditions is uptodate with
2056  // internal field
2057  meshMover_.setDisplacementPatchFields();
2058 
2059  Info<< typeName << " : Moving mesh ..." << endl;
2060  scalar oldErrorReduction = -1;
2061 
2062  bool meshOk = false;
2063 
2064  for (label iter = 0; iter < 2*nSnap ; iter++)
2065  {
2066  Info<< typeName
2067  << " : Iteration " << iter << endl;
2068  if (iter == nSnap)
2069  {
2070  Info<< typeName
2071  << " : Displacement scaling for error reduction set to 0."
2072  << endl;
2073  oldErrorReduction = meshMover_.setErrorReduction(0.0);
2074  }
2075 
2076  if
2077  (
2078  meshMover_.scaleMesh
2079  (
2080  checkFaces,
2081  baffles_,
2082  meshMover_.paramDict(),
2083  meshQualityDict,
2084  true,
2085  nAllowableErrors
2086  )
2087  )
2088  {
2089  Info<< typeName << " : Successfully moved mesh" << endl;
2090  meshOk = true;
2091  break;
2092  }
2093  }
2094 
2095  if (oldErrorReduction >= 0)
2096  {
2097  meshMover_.setErrorReduction(oldErrorReduction);
2098  }
2099 
2100  Info<< typeName << " : Finished moving mesh ..." << endl;
2101 
2102  return meshOk;
2103 }
2104 
2105 
2108  const dictionary& moveDict,
2109  const label nAllowableErrors,
2110  labelList& checkFaces
2111 )
2112 {
2113  // Read a few settings
2114  // ~~~~~~~~~~~~~~~~~~~
2115 
2116  //- Name of field specifying min thickness
2117  const word minThicknessName = word(moveDict.lookup("minThicknessName"));
2118 
2119 
2120  // The points have moved so before calculation update
2121  // the mesh and motionSolver accordingly
2122  movePoints(mesh().points());
2123  //
2125  // pointDisplacement_.boundaryField().updateCoeffs();
2126 
2127 
2128  // Extract out patch-wise displacement
2129  const indirectPrimitivePatch& pp = adaptPatchPtr_();
2130 
2131  scalarField zeroMinThickness;
2132  if (minThicknessName == "none")
2133  {
2134  zeroMinThickness = scalarField(pp.nPoints(), 0.0);
2135  }
2136  const scalarField& minThickness =
2137  (
2138  (minThicknessName == "none")
2139  ? zeroMinThickness
2140  : mesh().lookupObject<scalarField>(minThicknessName)
2141  );
2142 
2143 
2144  pointField patchDisp
2145  (
2146  pointDisplacement_.primitiveField(),
2147  pp.meshPoints()
2148  );
2149 
2151  (
2152  pp.nPoints(),
2154  );
2155  forAll(extrudeStatus, pointi)
2156  {
2157  if (mag(patchDisp[pointi]) <= minThickness[pointi]+small)
2158  {
2159  extrudeStatus[pointi] = snappyLayerDriver::NOEXTRUDE;
2160  }
2161  }
2162 
2163 
2164  // Solve displacement
2165  calculateDisplacement(moveDict, minThickness, extrudeStatus, patchDisp);
2166 
2167  //- Move mesh according to calculated displacement
2168  return shrinkMesh
2169  (
2170  moveDict, // meshQualityDict,
2171  nAllowableErrors, // nAllowableErrors
2172  checkFaces
2173  );
2174 }
2175 
2176 
2178 {
2180 
2181  // Update local data for new geometry
2182  adaptPatchPtr_().movePoints(p);
2183 
2184  // Update motionSmoother for new geometry
2185  meshMover_.movePoints();
2186 
2187  // Assume current mesh location is correct
2188  meshMover_.correct();
2189 }
2190 
2191 
2192 // ************************************************************************* //
Variant of pointEdgePoint with some transported additional data. WIP - should be templated on data li...
Definition: pointData.H:61
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:434
bool isA(const Type &t)
Check if a dynamic_cast to typeid is possible.
Definition: typeInfo.H:134
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
label nPoints() const
Return number of points supporting patch faces.
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
virtual void movePoints(const pointField &)
Update local data for geometry changes.
A line primitive.
Definition: line.H:56
void reset(T *=nullptr)
If object pointer already set, delete object and set to given.
Definition: autoPtrI.H:114
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
const Boundary & boundaryField() const
Return const-reference to the boundary field.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Unit conversion functions.
label nTotalPoints() const
Return total number of points in decomposed mesh. Not.
label iterate(const label maxIter)
Iterate until no changes or maxIter reached. Returns actual.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
A simple wrapper around bool so that it can be read as a word: true/false, on/off, yes/no, y/n, t/f, or none/any.
Definition: Switch.H:60
labelList identity(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:626
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
PrimitivePatch< IndirectList< face >, const pointField & > indirectPrimitivePatch
Foam::indirectPrimitivePatch.
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
Macros for easy insertion into run-time selection tables.
scalar degToRad(const scalar deg)
Conversion from degrees to radians.
const labelList & meshPoints() const
Return labelList of mesh points in patch. They are constructed.
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1131
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:48
List< bool > boolList
Bool container classes.
Definition: boolList.H:50
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:208
medialAxisMeshMover(const dictionary &dict, const List< labelPair > &baffles, pointVectorField &pointDisplacement)
Construct from dictionary and displacement field.
static PackedBoolList getMasterEdges(const polyMesh &)
Get per edge whether it is uncoupled or a master of a.
Definition: syncTools.C:109
A list of faces which address into the list of points.
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
dynamicFvMesh & mesh
dimensionedScalar cos(const dimensionedScalar &ds)
line< point, const point & > linePointRef
Line using referred points.
Definition: linePointRef.H:45
const pointField & points
List< edge > edgeList
Definition: edgeList.H:38
static void weightedSum(const polyMesh &mesh, const PackedBoolList &isMasterEdge, const labelList &meshPoints, const edgeList &edges, const scalarField &edgeWeights, const Field< Type > &data, Field< Type > &sum)
Helper: weighted sum (over all subset of mesh points) by.
static void calculateEdgeWeights(const polyMesh &mesh, const PackedBoolList &isMasterEdge, const labelList &meshPoints, const edgeList &edges, scalarField &edgeWeights, scalarField &invSumWeight)
Helper: calculate edge weights (1/length)
static void syncEdgeList(const polyMesh &, List< T > &, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronize values on all mesh edges.
virtual bool move(const dictionary &, const label nAllowableErrors, labelList &checkFaces)
Move mesh using current pointDisplacement boundary values.
A class for handling words, derived from string.
Definition: word.H:59
label nPoints
Do not extrude. No layers added.
static PackedBoolList getMasterEdges(const polyMesh &mesh, const labelList &meshEdges)
Determine master edge for subset of edges. If coupled.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
virtual const fileName & name() const
Return the name of the stream.
Definition: IOstream.H:297
virtual void movePoints(const pointField &)
Update local data for geometry changes.
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1394
word timeName
Definition: getTimeIndex.H:3
List< label > labelList
A List of labels.
Definition: labelList.H:56
static const zero Zero
Definition: zero.H:97
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &df)
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if successful.
Definition: doubleScalar.H:68
bool valid() const
Return true if the autoPtr valid (ie, the pointer is set)
Definition: autoPtrI.H:83
label readLabel(Istream &is)
Definition: label.H:64
dimensioned< scalar > magSqr(const dimensioned< Type > &)
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
OFstream which keeps track of vertices.
Definition: OBJstream.H:53
static const char nl
Definition: Ostream.H:265
defineTypeNameAndDebug(combustionModel, 0)
static PackedBoolList getMasterPoints(const polyMesh &)
Get per point whether it is uncoupled or a master of a.
Definition: syncTools.C:65
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:240
pointPatchField< vector > pointPatchVectorField
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
static tmp< pointField > pointNormals(const polyMesh &, const PrimitivePatch< FaceList, PointField > &)
Return parallel consistent point normals for patches using mesh points.
const dimensionedScalar mu
Atomic mass unit.
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
Virtual base class for mesh movers with externally provided displacement field giving the boundary co...
void setSize(const label)
Reset size of List.
Definition: List.C:281
A bit-packed bool list.
static T gAverage(const PackedBoolList &isMasterElem, const UList< T > &values)
Helper: calculate average.
label patchi
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
vector point
Point is a vector.
Definition: point.H:41
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
#define WarningInFunction
Report a warning using Foam::Warning.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const Time & time() const
Return time.
Definition: IOobject.C:360
labelList meshEdges(const edgeList &allEdges, const labelListList &cellEdges, const labelList &faceCells) const
Return labels of patch edges in the global edge list using.
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
dimensionedScalar lambda(laminarTransport.lookup("lambda"))
static void syncPointList(const polyMesh &, List< T > &, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronize values on all mesh points.
GeometricField< vector, pointPatchField, pointMesh > pointVectorField
Definition: pointFields.H:50
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
label n
Field< vector > vectorField
Specialisation of Field<T> for vector.
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
virtual Ostream & write(const token &)=0
Write next token to stream.
volScalarField & p
static PackedBoolList getMasterPoints(const polyMesh &mesh, const labelList &meshPoints)
Determine master point for subset of points. If coupled.
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:354
Wave propagation of information through grid. Every iteration information goes through one layer of e...
Definition: PointEdgeWave.H:85
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:233
Info<< "Finished reading KIVA file"<< endl;cellShapeList cellShapes(nPoints);labelList cellZoning(nPoints, -1);const cellModel &hex=*(cellModeller::lookup("hex"));labelList hexLabels(8);label activeCells=0;labelList pointMap(nPoints);forAll(pointMap, i){ pointMap[i]=i;}for(label i=0;i< nPoints;i++){ if(f[i] > 0.0) { hexLabels[0]=i;hexLabels[1]=i1tab[i];hexLabels[2]=i3tab[i1tab[i]];hexLabels[3]=i3tab[i];hexLabels[4]=i8tab[i];hexLabels[5]=i1tab[i8tab[i]];hexLabels[6]=i3tab[i1tab[i8tab[i]]];hexLabels[7]=i3tab[i8tab[i]];cellShapes[activeCells]=cellShape(hex, hexLabels);edgeList edges=cellShapes[activeCells].edges();forAll(edges, ei) { if(edges[ei].mag(points)< small) { label start=pointMap[edges[ei].start()];while(start !=pointMap[start]) { start=pointMap[start];} label end=pointMap[edges[ei].end()];while(end !=pointMap[end]) { end=pointMap[end];} label minLabel=min(start, end);pointMap[start]=pointMap[end]=minLabel;} } cellZoning[activeCells]=idreg[i];activeCells++;}}cellShapes.setSize(activeCells);cellZoning.setSize(activeCells);forAll(cellShapes, celli){ cellShape &cs=cellShapes[celli];forAll(cs, i) { cs[i]=pointMap[cs[i]];} cs.collapse();}label bcIDs[11]={-1, 0, 2, 4, -1, 5, -1, 6, 7, 8, 9};const label nBCs=12;const word *kivaPatchTypes[nBCs]={ &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &symmetryPolyPatch::typeName, &wedgePolyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &symmetryPolyPatch::typeName, &oldCyclicPolyPatch::typeName};enum patchTypeNames{ PISTON, VALVE, LINER, CYLINDERHEAD, AXIS, WEDGE, INFLOW, OUTFLOW, PRESIN, PRESOUT, SYMMETRYPLANE, CYCLIC};const char *kivaPatchNames[nBCs]={ "piston", "valve", "liner", "cylinderHead", "axis", "wedge", "inflow", "outflow", "presin", "presout", "symmetryPlane", "cyclic"};List< SLList< face > > pFaces[nBCs]
Definition: readKivaGrid.H:235
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:583
fileName path(UMean.rootPath()/UMean.caseName()/functionObjects::writeFile::outputPrefix/"graphs"/UMean.instance())