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