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-2023 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 {
43 
45  (
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 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(identityMap(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 =
362  coeffDict.lookup<label>("nSmoothSurfaceNormals");
363 
364  //- When is medial axis
365  word angleKey = "minMedialAxisAngle";
366  if (!coeffDict.found(angleKey))
367  {
368  // Backwards compatibility
369  angleKey = "minMedianAxisAngle";
370  }
371  scalar minMedialAxisAngleCos = Foam::cos
372  (
373  degToRad(coeffDict.lookup<scalar>(angleKey))
374  );
375 
376  //- Feature angle when to stop adding layers
377  const scalar featureAngle = coeffDict.lookup<scalar>("featureAngle");
378 
379  //- When to slip along wall
380  const scalar slipFeatureAngle =
381  (
382  coeffDict.found("slipFeatureAngle")
383  ? coeffDict.lookup<scalar>("slipFeatureAngle")
384  : 0.5*featureAngle
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(adaptPatchIDs_);
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  << slipFeatureAngle << " degrees." << endl;
690 
691  scalar slipFeatureAngleCos = Foam::cos
692  (
693  degToRad(slipFeatureAngle)
694  );
695  pointField pointNormals
696  (
697  PatchTools::pointNormals(mesh(), pp)
698  );
699 
700  forAll(meshPoints, i)
701  {
702  label pointi = meshPoints[i];
703 
704  if
705  (
706  pointWallDist[pointi].valid(dummyTrackData)
707  && !pointMedialDist[pointi].valid(dummyTrackData)
708  )
709  {
710  // Check if angle not too large.
711  scalar cosAngle =
712  (
713  -pointWallDist[pointi].v()
714  & pointNormals[i]
715  );
716  if (cosAngle > slipFeatureAngleCos)
717  {
718  // Extrusion direction practically perpendicular
719  // to the patch. Disable movement at the patch.
720 
721  maxPoints.append(pointi);
722  maxInfo.append
723  (
724  pointData
725  (
726  points[pointi],
727  0.0,
728  pointi, // passive data
729  Zero // passive data
730  )
731  );
732  pointMedialDist[pointi] = maxInfo.last();
733  }
734  else
735  {
736  // Extrusion direction makes angle with patch
737  // so allow sliding - don't insert zero points
738  }
739  }
740  }
741  }
742  }
743  }
744 
745  maxInfo.shrink();
746  maxPoints.shrink();
747 
748  // Do all calculations
749  PointEdgeWave<pointData> medialDistCalc
750  (
751  mesh(),
752  maxPoints,
753  maxInfo,
754 
755  pointMedialDist,
756  edgeMedialDist,
757  0, // max iterations
758  dummyTrackData
759  );
760  medialDistCalc.iterate(2*nMedialAxisIter);
761 
762 
763  // Extract medial axis distance as pointScalarField
764  forAll(pointMedialDist, pointi)
765  {
766  if (pointMedialDist[pointi].valid(dummyTrackData))
767  {
768  medialDist_[pointi] = Foam::sqrt
769  (
770  pointMedialDist[pointi].distSqr()
771  );
772  medialVec_[pointi] = pointMedialDist[pointi].origin();
773  }
774  else
775  {
776  // Unvisited. Do as if on medial axis so unmoving
777  medialDist_[pointi] = 0.0;
778  medialVec_[pointi] = point(1, 0, 0);
779  }
780  }
781  }
782 
783  // Extract transported surface normals as pointVectorField
784  forAll(dispVec_, i)
785  {
786  if (!pointWallDist[i].valid(dummyTrackData))
787  {
788  dispVec_[i] = vector(1, 0, 0);
789  }
790  else
791  {
792  dispVec_[i] = pointWallDist[i].v();
793  }
794  }
795 
796  // Smooth normal vectors. Do not change normals on pp.meshPoints
797  smoothNormals
798  (
799  nSmoothNormals,
800  isMeshMasterPoint,
801  isMeshMasterEdge,
802  meshPoints,
803  dispVec_
804  );
805 
806  // Calculate ratio point medial distance to point wall distance
807  forAll(medialRatio_, pointi)
808  {
809  if (!pointWallDist[pointi].valid(dummyTrackData))
810  {
811  medialRatio_[pointi] = 0.0;
812  }
813  else
814  {
815  scalar wDist2 = pointWallDist[pointi].distSqr();
816  scalar mDist = medialDist_[pointi];
817 
818  if (wDist2 < sqr(small) && mDist < small)
819  //- Note: maybe less strict:
820  //(
821  // wDist2 < sqr(meshRefiner_.mergeDistance())
822  // && mDist < meshRefiner_.mergeDistance()
823  //)
824  {
825  medialRatio_[pointi] = 0.0;
826  }
827  else
828  {
829  medialRatio_[pointi] = mDist / (Foam::sqrt(wDist2) + mDist);
830  }
831  }
832  }
833 
834 
835  if (debug)
836  {
837  Info<< typeName
838  << " : Writing medial axis fields:" << nl
839  << incrIndent
840  << "ratio of medial distance to wall distance : "
841  << medialRatio_.name() << nl
842  << "distance to nearest medial axis : "
843  << medialDist_.name() << nl
844  << "nearest medial axis location : "
845  << medialVec_.name() << nl
846  << "normal at nearest wall : "
847  << dispVec_.name() << nl
848  << decrIndent << nl
849  << endl;
850 
851  dispVec_.write();
852  medialRatio_.write();
853  medialDist_.write();
854  medialVec_.write();
855  }
856 }
857 
858 
859 bool Foam::medialAxisMeshMover::unmarkExtrusion
860 (
861  const label patchPointi,
862  pointField& patchDisp,
863  List<snappyLayerDriver::extrudeMode>& extrudeStatus
864 )
865 {
866  if (extrudeStatus[patchPointi] == snappyLayerDriver::EXTRUDE)
867  {
868  extrudeStatus[patchPointi] = snappyLayerDriver::NOEXTRUDE;
869  patchDisp[patchPointi] = Zero;
870  return true;
871  }
872  else if (extrudeStatus[patchPointi] == snappyLayerDriver::EXTRUDEREMOVE)
873  {
874  extrudeStatus[patchPointi] = snappyLayerDriver::NOEXTRUDE;
875  patchDisp[patchPointi] = Zero;
876  return true;
877  }
878  else
879  {
880  return false;
881  }
882 }
883 
884 
885 void Foam::medialAxisMeshMover::syncPatchDisplacement
886 (
887  const scalarField& minThickness,
888  pointField& patchDisp,
889  List<snappyLayerDriver::extrudeMode>& extrudeStatus
890 ) const
891 {
892  const indirectPrimitivePatch& pp = adaptPatchPtr_();
893  const labelList& meshPoints = pp.meshPoints();
894 
895  while (true)
896  {
897  label nChanged = 0;
898 
899  // Sync displacement (by taking min)
901  (
902  mesh(),
903  meshPoints,
904  patchDisp,
905  minMagSqrEqOp<vector>(),
906  point::rootMax // null value
907  );
908 
909  // Unmark if displacement too small
910  forAll(patchDisp, i)
911  {
912  if (mag(patchDisp[i]) < minThickness[i])
913  {
914  if (unmarkExtrusion(i, patchDisp, extrudeStatus))
915  {
916  nChanged++;
917  }
918  }
919  }
920 
921  if (!returnReduce(nChanged, sumOp<label>()))
922  {
923  break;
924  }
925  }
926 }
927 
928 
929 void Foam::medialAxisMeshMover::minSmoothField
930 (
931  const label nSmoothDisp,
932  const PackedBoolList& isPatchMasterPoint,
933  const PackedBoolList& isPatchMasterEdge,
934  const scalarField& fieldMin,
935  scalarField& field
936 ) const
937 {
938  const indirectPrimitivePatch& pp = adaptPatchPtr_();
939  const edgeList& edges = pp.edges();
940  const labelList& meshPoints = pp.meshPoints();
941 
942  scalarField edgeWeights(edges.size());
943  scalarField invSumWeight(meshPoints.size());
945  (
946  mesh(),
947  isPatchMasterEdge,
948  meshPoints,
949  edges,
950  edgeWeights,
951  invSumWeight
952  );
953 
954  // Get smoothly varying patch field.
955  Info<< typeName << " : Smoothing field ..." << endl;
956 
957  for (label iter = 0; iter < nSmoothDisp; iter++)
958  {
959  scalarField average(pp.nPoints());
961  (
962  mesh(),
963  isPatchMasterEdge,
964  meshPoints,
965  edges,
966  edgeWeights,
967  field,
968  average
969  );
970  average *= invSumWeight;
971 
972  // Transfer to field
973  forAll(field, pointi)
974  {
975  // full smoothing neighbours + point value
976  average[pointi] = 0.5*(field[pointi]+average[pointi]);
977 
978  // perform monotonic smoothing
979  if
980  (
981  average[pointi] < field[pointi]
982  && average[pointi] >= fieldMin[pointi]
983  )
984  {
985  field[pointi] = average[pointi];
986  }
987  }
988 
989  // Do residual calculation every so often.
990  if ((iter % 10) == 0)
991  {
992  scalar resid = meshRefinement::gAverage
993  (
994  isPatchMasterPoint,
995  mag(field-average)()
996  );
997  Info<< " Iteration " << iter << " residual " << resid << endl;
998  }
999  }
1000 }
1001 
1002 // Stop layer growth where mesh wraps around edge with a
1003 // large feature angle
1004 void Foam::medialAxisMeshMover::
1005 handleFeatureAngleLayerTerminations
1006 (
1007  const scalar minCos,
1008  const PackedBoolList& isPatchMasterPoint,
1009  const labelList& meshEdges,
1010  List<snappyLayerDriver::extrudeMode>& extrudeStatus,
1011  pointField& patchDisp,
1012  label& nPointCounter
1013 ) const
1014 {
1015  const indirectPrimitivePatch& pp = adaptPatchPtr_();
1016 
1017  // Mark faces that have all points extruded
1018  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1019 
1020  boolList extrudedFaces(pp.size(), true);
1021 
1022  forAll(pp.localFaces(), facei)
1023  {
1024  const face& f = pp.localFaces()[facei];
1025 
1026  forAll(f, fp)
1027  {
1028  if (extrudeStatus[f[fp]] == snappyLayerDriver::NOEXTRUDE)
1029  {
1030  extrudedFaces[facei] = false;
1031  break;
1032  }
1033  }
1034  }
1035 
1036 
1037 
1038  // label nOldPointCounter = nPointCounter;
1039 
1040  // Detect situation where two featureedge-neighbouring faces are partly or
1041  // not extruded and the edge itself is extruded. In this case unmark the
1042  // edge for extrusion.
1043 
1044 
1045  List<List<point>> edgeFaceNormals(pp.nEdges());
1046  List<List<bool>> edgeFaceExtrude(pp.nEdges());
1047 
1048  const labelListList& edgeFaces = pp.edgeFaces();
1049  const vectorField& faceNormals = pp.faceNormals();
1050 
1051  forAll(edgeFaces, edgei)
1052  {
1053  const labelList& eFaces = edgeFaces[edgei];
1054 
1055  edgeFaceNormals[edgei].setSize(eFaces.size());
1056  edgeFaceExtrude[edgei].setSize(eFaces.size());
1057  forAll(eFaces, i)
1058  {
1059  label facei = eFaces[i];
1060  edgeFaceNormals[edgei][i] = faceNormals[facei];
1061  edgeFaceExtrude[edgei][i] = extrudedFaces[facei];
1062  }
1063  }
1064 
1066  (
1067  mesh(),
1068  meshEdges,
1069  edgeFaceNormals,
1070  globalMeshData::ListPlusEqOp<List<point>>(), // combine operator
1071  List<point>() // null value
1072  );
1073 
1075  (
1076  mesh(),
1077  meshEdges,
1078  edgeFaceExtrude,
1079  globalMeshData::ListPlusEqOp<List<bool>>(), // combine operator
1080  List<bool>() // null value
1081  );
1082 
1083 
1084  forAll(edgeFaceNormals, edgei)
1085  {
1086  const List<point>& eFaceNormals = edgeFaceNormals[edgei];
1087  const List<bool>& eFaceExtrude = edgeFaceExtrude[edgei];
1088 
1089  if (eFaceNormals.size() == 2)
1090  {
1091  const edge& e = pp.edges()[edgei];
1092  label v0 = e[0];
1093  label v1 = e[1];
1094 
1095  if
1096  (
1097  extrudeStatus[v0] != snappyLayerDriver::NOEXTRUDE
1098  || extrudeStatus[v1] != snappyLayerDriver::NOEXTRUDE
1099  )
1100  {
1101  if (!eFaceExtrude[0] || !eFaceExtrude[1])
1102  {
1103  const vector& n0 = eFaceNormals[0];
1104  const vector& n1 = eFaceNormals[1];
1105 
1106  if ((n0 & n1) < minCos)
1107  {
1108  if (unmarkExtrusion(v0, patchDisp, extrudeStatus))
1109  {
1110  if (isPatchMasterPoint[v0])
1111  {
1112  nPointCounter++;
1113  }
1114  }
1115  if (unmarkExtrusion(v1, patchDisp, extrudeStatus))
1116  {
1117  if (isPatchMasterPoint[v1])
1118  {
1119  nPointCounter++;
1120  }
1121  }
1122  }
1123  }
1124  }
1125  }
1126  }
1127 
1128  // Info<< "Added "
1129  // << returnReduce(nPointCounter-nOldPointCounter, sumOp<label>())
1130  // << " point not to extrude." << endl;
1131 }
1132 
1133 
1134 // Find isolated islands (points, edges and faces and layer terminations)
1135 // in the layer mesh and stop any layer growth at these points.
1136 void Foam::medialAxisMeshMover::findIsolatedRegions
1137 (
1138  const scalar minCosLayerTermination,
1139  const bool detectExtrusionIsland,
1140  const PackedBoolList& isPatchMasterPoint,
1141  const PackedBoolList& isPatchMasterEdge,
1142  const labelList& meshEdges,
1143  const scalarField& minThickness,
1144  List<snappyLayerDriver::extrudeMode>& extrudeStatus,
1145  pointField& patchDisp
1146 ) const
1147 {
1148  const indirectPrimitivePatch& pp = adaptPatchPtr_();
1149  const labelListList& pointFaces = pp.pointFaces();
1150  const labelList& meshPoints = pp.meshPoints();
1151 
1152  Info<< typeName << " : Removing isolated regions ..." << endl;
1153 
1154  // Keep count of number of points unextruded
1155  label nPointCounter = 0;
1156 
1157 
1158  autoPtr<OBJstream> str;
1159  if (debug)
1160  {
1161  str.reset
1162  (
1163  new OBJstream
1164  (
1165  mesh().time().path()
1166  / "islandExcludePoints_"
1167  + mesh().time().name()
1168  + ".obj"
1169  )
1170  );
1171  Info<< typeName
1172  << " : Writing points surrounded by non-extruded points to "
1173  << str().name() << endl;
1174  }
1175 
1176  while (true)
1177  {
1178  // Stop layer growth where mesh wraps around edge with a
1179  // large feature angle
1180  handleFeatureAngleLayerTerminations
1181  (
1182  minCosLayerTermination,
1183  isPatchMasterPoint,
1184  meshEdges,
1185 
1186  extrudeStatus,
1187  patchDisp,
1188  nPointCounter
1189  );
1190 
1191  syncPatchDisplacement(minThickness, patchDisp, extrudeStatus);
1192 
1193 
1194 
1195  // Detect either:
1196  // - point where all surrounding points are not extruded
1197  // (detectExtrusionIsland)
1198  // or
1199  // - point where all the faces surrounding it are not fully
1200  // extruded
1201 
1202  boolList keptPoints(pp.nPoints(), false);
1203 
1204  if (detectExtrusionIsland)
1205  {
1206  // Do not extrude from point where all neighbouring
1207  // points are not grown
1208  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1209 
1210  labelList islandPoint(pp.size(), -1);
1211  forAll(pp, facei)
1212  {
1213  const face& f = pp.localFaces()[facei];
1214 
1215  forAll(f, fp)
1216  {
1217  if (extrudeStatus[f[fp]] != snappyLayerDriver::NOEXTRUDE)
1218  {
1219  if (islandPoint[facei] == -1)
1220  {
1221  // First point to extrude
1222  islandPoint[facei] = f[fp];
1223  }
1224  else if (islandPoint[facei] != -2)
1225  {
1226  // Second or more point to extrude
1227  islandPoint[facei] = -2;
1228  }
1229  }
1230  }
1231  }
1232 
1233  // islandPoint:
1234  // -1 : no point extruded on face
1235  // -2 : >= 2 points extruded on face
1236  // >=0: label of point extruded
1237 
1238  // Check all surrounding faces that I am the islandPoint
1239  forAll(pointFaces, patchPointi)
1240  {
1241  if (extrudeStatus[patchPointi] != snappyLayerDriver::NOEXTRUDE)
1242  {
1243  const labelList& pFaces = pointFaces[patchPointi];
1244 
1245  forAll(pFaces, i)
1246  {
1247  label facei = pFaces[i];
1248  if (islandPoint[facei] != patchPointi)
1249  {
1250  keptPoints[patchPointi] = true;
1251  break;
1252  }
1253  }
1254  }
1255  }
1256  }
1257  else
1258  {
1259  // Do not extrude from point where all neighbouring
1260  // faces are not grown
1261  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1262 
1263  boolList extrudedFaces(pp.size(), true);
1264  forAll(pp.localFaces(), facei)
1265  {
1266  const face& f = pp.localFaces()[facei];
1267  forAll(f, fp)
1268  {
1269  if (extrudeStatus[f[fp]] == snappyLayerDriver::NOEXTRUDE)
1270  {
1271  extrudedFaces[facei] = false;
1272  break;
1273  }
1274  }
1275  }
1276 
1277  const labelListList& pointFaces = pp.pointFaces();
1278 
1279  forAll(keptPoints, patchPointi)
1280  {
1281  const labelList& pFaces = pointFaces[patchPointi];
1282 
1283  forAll(pFaces, i)
1284  {
1285  label facei = pFaces[i];
1286  if (extrudedFaces[facei])
1287  {
1288  keptPoints[patchPointi] = true;
1289  break;
1290  }
1291  }
1292  }
1293  }
1294 
1295 
1297  (
1298  mesh(),
1299  meshPoints,
1300  keptPoints,
1301  orEqOp<bool>(),
1302  false // null value
1303  );
1304 
1305  label nChanged = 0;
1306 
1307  forAll(keptPoints, patchPointi)
1308  {
1309  if (!keptPoints[patchPointi])
1310  {
1311  if (unmarkExtrusion(patchPointi, patchDisp, extrudeStatus))
1312  {
1313  nPointCounter++;
1314  nChanged++;
1315 
1316  if (str.valid())
1317  {
1318  str().write(pp.points()[meshPoints[patchPointi]]);
1319  }
1320  }
1321  }
1322  }
1323 
1324 
1325  if (returnReduce(nChanged, sumOp<label>()) == 0)
1326  {
1327  break;
1328  }
1329  }
1330 
1331  const edgeList& edges = pp.edges();
1332 
1333 
1334  // Count number of mesh edges using a point
1335  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1336 
1337  labelList isolatedPoint(pp.nPoints(),0);
1338 
1339  forAll(edges, edgei)
1340  {
1341  if (isPatchMasterEdge[edgei])
1342  {
1343  const edge& e = edges[edgei];
1344 
1345  label v0 = e[0];
1346  label v1 = e[1];
1347 
1348  if (extrudeStatus[v1] != snappyLayerDriver::NOEXTRUDE)
1349  {
1350  isolatedPoint[v0] += 1;
1351  }
1352  if (extrudeStatus[v0] != snappyLayerDriver::NOEXTRUDE)
1353  {
1354  isolatedPoint[v1] += 1;
1355  }
1356  }
1357  }
1358 
1360  (
1361  mesh(),
1362  meshPoints,
1363  isolatedPoint,
1364  plusEqOp<label>(),
1365  label(0) // null value
1366  );
1367 
1368  // stop layer growth on isolated faces
1369  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1370  forAll(pp, facei)
1371  {
1372  const face& f = pp.localFaces()[facei];
1373  bool failed = false;
1374  forAll(f, fp)
1375  {
1376  if (isolatedPoint[f[fp]] > 2)
1377  {
1378  failed = true;
1379  break;
1380  }
1381  }
1382  bool allPointsExtruded = true;
1383  if (!failed)
1384  {
1385  forAll(f, fp)
1386  {
1387  if (extrudeStatus[f[fp]] == snappyLayerDriver::NOEXTRUDE)
1388  {
1389  allPointsExtruded = false;
1390  break;
1391  }
1392  }
1393 
1394  if (allPointsExtruded)
1395  {
1396  forAll(f, fp)
1397  {
1398  if
1399  (
1400  unmarkExtrusion
1401  (
1402  f[fp],
1403  patchDisp,
1404  extrudeStatus
1405  )
1406  )
1407  {
1408  nPointCounter++;
1409 
1410  if (str.valid())
1411  {
1412  str().write(pp.points()[meshPoints[f[fp]]]);
1413  }
1414  }
1415  }
1416  }
1417  }
1418  }
1419 
1420  reduce(nPointCounter, sumOp<label>());
1421  Info<< typeName
1422  << " : Number of isolated points extrusion stopped : "<< nPointCounter
1423  << endl;
1424 }
1425 
1426 
1427 void Foam::medialAxisMeshMover::smoothLambdaMuDisplacement
1428 (
1429  const label nSmoothDisp,
1430  const PackedBoolList& isMeshMasterPoint,
1431  const PackedBoolList& isMeshMasterEdge,
1432  vectorField& displacement
1433 ) const
1434 {
1435  const edgeList& edges = mesh().edges();
1436 
1437  // Correspondence between local edges/points and mesh edges/points
1438  const labelList meshPoints(identityMap(mesh().nPoints()));
1439 
1440  // Calculate inverse sum of weights
1441  scalarField edgeWeights(mesh().nEdges());
1442  scalarField invSumWeight(meshPoints.size());
1444  (
1445  mesh(),
1446  isMeshMasterEdge,
1447  meshPoints,
1448  edges,
1449  edgeWeights,
1450  invSumWeight
1451  );
1452 
1453  // Get smoothly varying patch field.
1454  Info<< typeName << " : Smoothing displacement ..." << endl;
1455 
1456  const scalar lambda = 0.33;
1457  const scalar mu = -0.34;
1458 
1460 
1461  for (label iter = 0; iter < nSmoothDisp; iter++)
1462  {
1464  (
1465  mesh(),
1466  isMeshMasterEdge,
1467  meshPoints,
1468  edges,
1469  edgeWeights,
1470  displacement,
1471  average
1472  );
1473  average *= invSumWeight;
1474 
1475  forAll(displacement, i)
1476  {
1477  if (medialRatio_[i] > small && medialRatio_[i] < 1-small)
1478  {
1479  displacement[i] = (1-lambda)*displacement[i]+lambda*average[i];
1480  }
1481  }
1482 
1484  (
1485  mesh(),
1486  isMeshMasterEdge,
1487  meshPoints,
1488  edges,
1489  edgeWeights,
1490  displacement,
1491  average
1492  );
1493  average *= invSumWeight;
1494 
1495 
1496  forAll(displacement, i)
1497  {
1498  if (medialRatio_[i] > small && medialRatio_[i] < 1-small)
1499  {
1500  displacement[i] = (1-mu)*displacement[i]+mu*average[i];
1501  }
1502  }
1503 
1504 
1505  // Do residual calculation every so often.
1506  if ((iter % 10) == 0)
1507  {
1508  scalar resid = meshRefinement::gAverage
1509  (
1510  isMeshMasterPoint,
1511  mag(displacement-average)()
1512  );
1513  Info<< " Iteration " << iter << " residual " << resid << endl;
1514  }
1515  }
1516 }
1517 
1518 
1519 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
1520 
1522 (
1523  const dictionary& dict,
1524  const List<labelPair>& baffles,
1525  pointVectorField& pointDisplacement
1526 )
1527 :
1528  externalDisplacementMeshMover(dict, baffles, pointDisplacement),
1529  adaptPatchIDs_(getFixedValueBCs(pointDisplacement)),
1530  adaptPatchPtr_(getPatch(mesh(), adaptPatchIDs_)),
1531  scale_
1532  (
1533  IOobject
1534  (
1535  "scale",
1536  pointDisplacement.time().name(),
1537  pointDisplacement.db(),
1538  IOobject::NO_READ,
1539  IOobject::AUTO_WRITE
1540  ),
1541  pMesh(),
1543  ),
1544  oldPoints_(mesh().points()),
1545  meshMover_
1546  (
1547  const_cast<polyMesh&>(mesh()),
1548  const_cast<pointMesh&>(pMesh()),
1549  adaptPatchPtr_(),
1550  pointDisplacement,
1551  scale_,
1552  oldPoints_,
1553  adaptPatchIDs_,
1554  dict
1555  ),
1556  dispVec_
1557  (
1558  IOobject
1559  (
1560  "dispVec",
1561  pointDisplacement.time().name(),
1562  pointDisplacement.db(),
1563  IOobject::NO_READ,
1564  IOobject::NO_WRITE,
1565  false
1566  ),
1567  pMesh(),
1569  ),
1570  medialRatio_
1571  (
1572  IOobject
1573  (
1574  "medialRatio",
1575  pointDisplacement.time().name(),
1576  pointDisplacement.db(),
1577  IOobject::NO_READ,
1578  IOobject::NO_WRITE,
1579  false
1580  ),
1581  pMesh(),
1583  ),
1584  medialDist_
1585  (
1586  IOobject
1587  (
1588  "pointMedialDist",
1589  pointDisplacement.time().name(),
1590  pointDisplacement.db(),
1591  IOobject::NO_READ,
1592  IOobject::NO_WRITE,
1593  false
1594  ),
1595  pMesh(),
1597  ),
1598  medialVec_
1599  (
1600  IOobject
1601  (
1602  "medialVec",
1603  pointDisplacement.time().name(),
1604  pointDisplacement.db(),
1605  IOobject::NO_READ,
1606  IOobject::NO_WRITE,
1607  false
1608  ),
1609  pMesh(),
1611  )
1612 {
1613  update(dict);
1614 }
1615 
1616 
1617 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
1618 
1620 {}
1621 
1622 
1623 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1624 
1625 void Foam::medialAxisMeshMover::calculateDisplacement
1626 (
1627  const dictionary& coeffDict,
1628  const scalarField& minThickness,
1629  List<snappyLayerDriver::extrudeMode>& extrudeStatus,
1630  pointField& patchDisp
1631 )
1632 {
1633  Info<< typeName << " : Smoothing using Medial Axis ..." << endl;
1634 
1635  const indirectPrimitivePatch& pp = adaptPatchPtr_;
1636  const labelList& meshPoints = pp.meshPoints();
1637 
1638 
1639  // Read settings
1640  // ~~~~~~~~~~~~~
1641 
1642  //- (lambda-mu) smoothing of internal displacement
1643  const label nSmoothDisplacement = coeffDict.lookupOrDefault
1644  (
1645  "nSmoothDisplacement",
1646  0
1647  );
1648 
1649  //- Layer thickness too big
1650  const scalar maxThicknessToMedialRatio =
1651  coeffDict.lookup<scalar>("maxThicknessToMedialRatio");
1652 
1653  //- Feature angle when to stop adding layers
1654  const scalar featureAngle = coeffDict.lookup<scalar>("featureAngle");
1655 
1656  //- Stop layer growth where mesh wraps around sharp edge
1657  const scalar minCosLayerTermination = Foam::cos
1658  (
1659  degToRad(0.5*featureAngle)
1660  );
1661 
1662  //- Smoothing wanted patch thickness
1663  const label nSmoothPatchThickness =
1664  coeffDict.lookup<label>("nSmoothThickness");
1665 
1666  //- Number of edges walking out
1667  const label nMedialAxisIter = coeffDict.lookupOrDefault<label>
1668  (
1669  "nMedialAxisIter",
1670  mesh().globalData().nTotalPoints()
1671  );
1672 
1673  //- Use strict extrusionIsland detection
1674  const Switch detectExtrusionIsland = coeffDict.lookupOrDefault<Switch>
1675  (
1676  "detectExtrusionIsland",
1677  false
1678  );
1679 
1680 
1681  // Precalculate master points/edge (only relevant for shared points/edges)
1682  const PackedBoolList isMeshMasterPoint(syncTools::getMasterPoints(mesh()));
1683  const PackedBoolList isMeshMasterEdge(syncTools::getMasterEdges(mesh()));
1684  // Precalculate meshEdge per pp edge
1685  const labelList meshEdges
1686  (
1687  pp.meshEdges
1688  (
1689  mesh().edges(),
1690  mesh().pointEdges()
1691  )
1692  );
1693 
1694  // Precalculate (patch) master point/edge
1695  const PackedBoolList isPatchMasterPoint
1696  (
1698  (
1699  mesh(),
1700  meshPoints
1701  )
1702  );
1703  const PackedBoolList isPatchMasterEdge
1704  (
1706  (
1707  mesh(),
1708  meshEdges
1709  )
1710  );
1711 
1712 
1713  scalarField thickness(patchDisp.size());
1714 
1715  thickness = mag(patchDisp);
1716 
1717  forAll(thickness, patchPointi)
1718  {
1719  if (extrudeStatus[patchPointi] == snappyLayerDriver::NOEXTRUDE)
1720  {
1721  thickness[patchPointi] = 0.0;
1722  }
1723  }
1724 
1725  label numThicknessRatioExclude = 0;
1726 
1727  // reduce thickness where thickness/medial axis distance large
1728  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1729 
1730  autoPtr<OBJstream> str;
1731  if (debug)
1732  {
1733  str.reset
1734  (
1735  new OBJstream
1736  (
1737  mesh().time().path()
1738  / "thicknessRatioExcludePoints_"
1739  + mesh().time().name()
1740  + ".obj"
1741  )
1742  );
1743  Info<< typeName
1744  << " : Writing points with too large an extrusion distance to "
1745  << str().name() << endl;
1746  }
1747 
1748  autoPtr<OBJstream> medialVecStr;
1749  if (debug)
1750  {
1751  medialVecStr.reset
1752  (
1753  new OBJstream
1754  (
1755  mesh().time().path()
1756  / "thicknessRatioExcludeMedialVec_"
1757  + mesh().time().name()
1758  + ".obj"
1759  )
1760  );
1761  Info<< typeName
1762  << " : Writing medial axis vectors on points with too large"
1763  << " an extrusion distance to " << medialVecStr().name() << endl;
1764  }
1765 
1766  forAll(meshPoints, patchPointi)
1767  {
1768  if (extrudeStatus[patchPointi] != snappyLayerDriver::NOEXTRUDE)
1769  {
1770  label pointi = meshPoints[patchPointi];
1771 
1772  //- Option 1: look only at extrusion thickness v.s. distance
1773  // to nearest (medial axis or static) point.
1774  scalar mDist = medialDist_[pointi];
1775  scalar thicknessRatio = thickness[patchPointi]/(mDist+vSmall);
1776 
1777  //- Option 2: Look at component in the direction
1778  // of nearest (medial axis or static) point.
1779  vector n =
1780  patchDisp[patchPointi]
1781  / (mag(patchDisp[patchPointi]) + vSmall);
1782  vector mVec = mesh().points()[pointi]-medialVec_[pointi];
1783  mVec /= mag(mVec)+vSmall;
1784  thicknessRatio *= (n&mVec);
1785 
1786  if (thicknessRatio > maxThicknessToMedialRatio)
1787  {
1788  // Truncate thickness.
1789  if (debug&2)
1790  {
1791  Pout<< "truncating displacement at "
1792  << mesh().points()[pointi]
1793  << " from " << thickness[patchPointi]
1794  << " to "
1795  << 0.5
1796  *(
1797  minThickness[patchPointi]
1798  +thickness[patchPointi]
1799  )
1800  << " medial direction:" << mVec
1801  << " extrusion direction:" << n
1802  << " with thicknessRatio:" << thicknessRatio
1803  << endl;
1804  }
1805 
1806  thickness[patchPointi] =
1807  0.5*(minThickness[patchPointi]+thickness[patchPointi]);
1808 
1809  patchDisp[patchPointi] = thickness[patchPointi]*n;
1810 
1811  if (isPatchMasterPoint[patchPointi])
1812  {
1813  numThicknessRatioExclude++;
1814  }
1815 
1816  if (str.valid())
1817  {
1818  const point& pt = mesh().points()[pointi];
1819  str().write(linePointRef(pt, pt+patchDisp[patchPointi]));
1820  }
1821  if (medialVecStr.valid())
1822  {
1823  const point& pt = mesh().points()[pointi];
1824  medialVecStr().write
1825  (
1826  linePointRef
1827  (
1828  pt,
1829  medialVec_[pointi]
1830  )
1831  );
1832  }
1833  }
1834  }
1835  }
1836 
1837  reduce(numThicknessRatioExclude, sumOp<label>());
1838  Info<< typeName << " : Reducing layer thickness at "
1839  << numThicknessRatioExclude
1840  << " nodes where thickness to medial axis distance is large " << endl;
1841 
1842 
1843  // find points where layer growth isolated to a lone point, edge or face
1844 
1845  findIsolatedRegions
1846  (
1847  minCosLayerTermination,
1848  detectExtrusionIsland,
1849 
1850  isPatchMasterPoint,
1851  isPatchMasterEdge,
1852  meshEdges,
1853  minThickness,
1854 
1855  extrudeStatus,
1856  patchDisp
1857  );
1858 
1859  // Update thickness for changed extrusion
1860  forAll(thickness, patchPointi)
1861  {
1862  if (extrudeStatus[patchPointi] == snappyLayerDriver::NOEXTRUDE)
1863  {
1864  thickness[patchPointi] = 0.0;
1865  }
1866  }
1867 
1868 
1869  // smooth layer thickness on moving patch
1870  minSmoothField
1871  (
1872  nSmoothPatchThickness,
1873  isPatchMasterPoint,
1874  isPatchMasterEdge,
1875  minThickness,
1876 
1877  thickness
1878  );
1879 
1880 
1881  // Dummy additional info for PointEdgeWave
1882  int dummyTrackData = 0;
1883 
1884  List<pointData> pointWallDist(mesh().nPoints());
1885 
1886  const pointField& points = mesh().points();
1887  // 1. Calculate distance to points where displacement is specified.
1888  // This wave is used to transport layer thickness
1889  {
1890  // Distance to wall and medial axis on edges.
1891  List<pointData> edgeWallDist(mesh().nEdges());
1892  labelList wallPoints(meshPoints.size());
1893 
1894  // Seed data.
1895  List<pointData> wallInfo(meshPoints.size());
1896 
1897  forAll(meshPoints, patchPointi)
1898  {
1899  label pointi = meshPoints[patchPointi];
1900  wallPoints[patchPointi] = pointi;
1901  wallInfo[patchPointi] = pointData
1902  (
1903  points[pointi],
1904  0.0,
1905  thickness[patchPointi], // transport layer thickness
1906  Zero // passive vector
1907  );
1908  }
1909 
1910  // Do all calculations
1911  PointEdgeWave<pointData> wallDistCalc
1912  (
1913  mesh(),
1914  wallPoints,
1915  wallInfo,
1916  pointWallDist,
1917  edgeWallDist,
1918  0, // max iterations
1919  dummyTrackData
1920  );
1921  wallDistCalc.iterate(nMedialAxisIter);
1922  }
1923 
1924 
1925  // Calculate scaled displacement vector
1926  pointField& displacement = pointDisplacement_;
1927 
1928  forAll(displacement, pointi)
1929  {
1930  if (!pointWallDist[pointi].valid(dummyTrackData))
1931  {
1932  displacement[pointi] = Zero;
1933  }
1934  else
1935  {
1936  // 1) displacement on nearest wall point, scaled by medialRatio
1937  // (wall distance / medial distance)
1938  // 2) pointWallDist[pointi].s() is layer thickness transported
1939  // from closest wall point.
1940  // 3) shrink in opposite direction of addedPoints
1941  displacement[pointi] =
1942  -medialRatio_[pointi]
1943  *pointWallDist[pointi].s()
1944  *dispVec_[pointi];
1945  }
1946  }
1947 
1948 
1949  // Smear displacement away from fixed values (medialRatio=0 or 1)
1950  if (nSmoothDisplacement > 0)
1951  {
1952  smoothLambdaMuDisplacement
1953  (
1954  nSmoothDisplacement,
1955  isMeshMasterPoint,
1956  isMeshMasterEdge,
1957  displacement
1958  );
1959  }
1960 }
1961 
1962 
1963 bool Foam::medialAxisMeshMover::shrinkMesh
1964 (
1965  const dictionary& meshQualityDict,
1966  const label nAllowableErrors,
1967  const labelList& checkFaces
1968 )
1969 {
1970  //- Number of attempts shrinking the mesh
1971  const label nSnap = meshQualityDict.lookup<label>("nRelaxIter");
1972 
1973  // Make sure displacement boundary conditions is up-to-date with
1974  // internal field
1975  meshMover_.setDisplacementPatchFields();
1976 
1977  Info<< typeName << " : Moving mesh ..." << endl;
1978  scalar oldErrorReduction = -1;
1979 
1980  bool meshOk = false;
1981 
1982  for (label iter = 0; iter < 2*nSnap ; iter++)
1983  {
1984  Info<< typeName
1985  << " : Iteration " << iter << endl;
1986  if (iter == nSnap)
1987  {
1988  Info<< typeName
1989  << " : Displacement scaling for error reduction set to 0."
1990  << endl;
1991  oldErrorReduction = meshMover_.setErrorReduction(0.0);
1992  }
1993 
1994  if
1995  (
1996  meshMover_.scaleMesh
1997  (
1998  checkFaces,
1999  baffles_,
2000  meshMover_.paramDict(),
2001  meshQualityDict,
2002  true,
2003  nAllowableErrors
2004  )
2005  )
2006  {
2007  Info<< typeName << " : Successfully moved mesh" << endl;
2008  meshOk = true;
2009  break;
2010  }
2011  }
2012 
2013  if (oldErrorReduction >= 0)
2014  {
2015  meshMover_.setErrorReduction(oldErrorReduction);
2016  }
2017 
2018  Info<< typeName << " : Finished moving mesh ..." << endl;
2019 
2020  return meshOk;
2021 }
2022 
2023 
2025 (
2026  const dictionary& moveDict,
2027  const label nAllowableErrors,
2028  const labelList& checkFaces
2029 )
2030 {
2031  // Read a few settings
2032  // ~~~~~~~~~~~~~~~~~~~
2033 
2034  //- Name of field specifying min thickness
2035  const word minThicknessName = word(moveDict.lookup("minThicknessName"));
2036 
2037 
2038  // The points have moved so before calculation update
2039  // the mesh and motionSolver accordingly
2040  movePoints(mesh().points());
2041  //
2043  // pointDisplacement_.boundaryField().updateCoeffs();
2044 
2045 
2046  // Extract out patch-wise displacement
2047  const indirectPrimitivePatch& pp = adaptPatchPtr_();
2048 
2049  scalarField zeroMinThickness;
2050  if (minThicknessName == "none")
2051  {
2052  zeroMinThickness = scalarField(pp.nPoints(), 0.0);
2053  }
2054  const scalarField& minThickness =
2055  (
2056  (minThicknessName == "none")
2057  ? zeroMinThickness
2058  : mesh().lookupObject<scalarField>(minThicknessName)
2059  );
2060 
2061 
2062  pointField patchDisp
2063  (
2064  pointDisplacement_.primitiveField(),
2065  pp.meshPoints()
2066  );
2067 
2069  (
2070  pp.nPoints(),
2072  );
2073  forAll(extrudeStatus, pointi)
2074  {
2075  if (mag(patchDisp[pointi]) <= minThickness[pointi]+small)
2076  {
2077  extrudeStatus[pointi] = snappyLayerDriver::NOEXTRUDE;
2078  }
2079  }
2080 
2081 
2082  // Solve displacement
2083  calculateDisplacement(moveDict, minThickness, extrudeStatus, patchDisp);
2084 
2085  //- Move mesh according to calculated displacement
2086  return shrinkMesh
2087  (
2088  moveDict, // meshQualityDict,
2089  nAllowableErrors, // nAllowableErrors
2090  checkFaces
2091  );
2092 }
2093 
2094 
2096 {
2098 
2099  // Update local data for new geometry
2100  adaptPatchPtr_().clearGeom();
2101 
2102  // Update motionSmoother for new geometry
2103  meshMover_.movePoints();
2104 
2105  // Assume current mesh location is correct
2106  meshMover_.correct();
2107 }
2108 
2109 
2110 // ************************************************************************* //
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:117
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:160
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:860
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:52
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:105
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:235
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
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
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:228
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:139
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:211
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:260
dimensioned< scalar > magSqr(const dimensioned< Type > &)
dimensionedScalar cos(const dimensionedScalar &ds)
scalar degToRad(const scalar deg)
Conversion from degrees to radians.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
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
Unit conversion functions.