snappyLayerDriverShrink.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) 2011-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 Description
25  Shrinking mesh (part of adding cell layers)
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "snappyLayerDriver.H"
30 #include "fvMesh.H"
31 #include "Time.H"
32 #include "pointFields.H"
33 #include "motionSmoother.H"
34 #include "pointData.H"
35 #include "PointEdgeWave.H"
36 #include "OBJstream.H"
37 #include "meshTools.H"
38 #include "PatchTools.H"
39 
40 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
41 
42 // Calculate inverse sum of edge weights (currently always 1.0)
43 void Foam::snappyLayerDriver::sumWeights
44 (
45  const PackedBoolList& isMasterEdge,
46  const labelList& meshEdges,
47  const labelList& meshPoints,
48  const edgeList& edges,
49  scalarField& invSumWeight
50 ) const
51 {
52  const pointField& pts = meshRefiner_.mesh().points();
53 
54  invSumWeight = 0;
55 
56  forAll(edges, edgei)
57  {
58  if (isMasterEdge.get(meshEdges[edgei]) == 1)
59  {
60  const edge& e = edges[edgei];
61  // scalar eWeight = edgeWeights[edgei];
62  // scalar eWeight = 1.0;
63 
64  scalar eMag = max
65  (
66  vSmall,
67  mag
68  (
69  pts[meshPoints[e[1]]]
70  - pts[meshPoints[e[0]]]
71  )
72  );
73  scalar eWeight = 1.0/eMag;
74 
75  invSumWeight[e[0]] += eWeight;
76  invSumWeight[e[1]] += eWeight;
77  }
78  }
79 
81  (
82  meshRefiner_.mesh(),
83  meshPoints,
84  invSumWeight,
85  plusEqOp<scalar>(),
86  scalar(0) // null value
87  );
88 
89  forAll(invSumWeight, pointi)
90  {
91  scalar w = invSumWeight[pointi];
92 
93  if (w > 0.0)
94  {
95  invSumWeight[pointi] = 1.0/w;
96  }
97  }
98 }
99 
100 
101 // Smooth field on moving patch
102 void Foam::snappyLayerDriver::smoothField
103 (
104  const motionSmoother& meshMover,
105  const PackedBoolList& isMasterPoint,
106  const PackedBoolList& isMasterEdge,
107  const labelList& meshEdges,
108  const scalarField& fieldMin,
109  const label nSmoothDisp,
110  scalarField& field
111 ) const
112 {
113  const indirectPrimitivePatch& pp = meshMover.patch();
114  const edgeList& edges = pp.edges();
115  const labelList& meshPoints = pp.meshPoints();
116 
117  scalarField invSumWeight(pp.nPoints());
118  sumWeights
119  (
120  isMasterEdge,
121  meshEdges,
122  meshPoints,
123  edges,
124  invSumWeight
125  );
126 
127  // Get smoothly varying patch field.
128  Info<< "shrinkMeshDistance : Smoothing field ..." << endl;
129 
130  for (label iter = 0; iter < nSmoothDisp; iter++)
131  {
132  scalarField average(pp.nPoints());
133  averageNeighbours
134  (
135  meshMover.mesh(),
136  isMasterEdge,
137  meshEdges,
138  meshPoints,
139  edges,
140  invSumWeight,
141  field,
142  average
143  );
144 
145  // Transfer to field
146  forAll(field, pointi)
147  {
148  // full smoothing neighbours + point value
149  average[pointi] = 0.5*(field[pointi]+average[pointi]);
150 
151  // perform monotonic smoothing
152  if
153  (
154  average[pointi] < field[pointi]
155  && average[pointi] >= fieldMin[pointi]
156  )
157  {
158  field[pointi] = average[pointi];
159  }
160  }
161 
162  // Do residual calculation every so often.
163  if ((iter % 10) == 0)
164  {
165  scalar resid = meshRefinement::gAverage
166  (
167  meshMover.mesh(),
168  isMasterPoint,
169  meshPoints,
170  mag(field-average)()
171  );
172  Info<< " Iteration " << iter << " residual " << resid << endl;
173  }
174  }
175 }
176 //XXXXXXXXX
177 //void Foam::snappyLayerDriver::smoothField
178 //(
179 // const motionSmoother& meshMover,
180 // const PackedBoolList& isMasterEdge,
181 // const labelList& meshEdges,
182 // const scalarField& fieldMin,
183 // const label nSmoothDisp,
184 // scalarField& field
185 //) const
186 //{
187 // const indirectPrimitivePatch& pp = meshMover.patch();
188 // const edgeList& edges = pp.edges();
189 // const labelList& meshPoints = pp.meshPoints();
190 //
191 // scalarField invSumWeight(pp.nPoints());
192 // sumWeights
193 // (
194 // isMasterEdge,
195 // meshEdges,
196 // meshPoints,
197 // edges,
198 // invSumWeight
199 // );
200 //
201 // // Get smoothly varying patch field.
202 // Info<< "shrinkMeshDistance : (lambda-mu) Smoothing field ..." << endl;
203 //
204 //
205 // const scalar lambda = 0.33;
206 // const scalar mu = -0.34;
207 //
208 // for (label iter = 0; iter < 90; iter++)
209 // {
210 // scalarField average(pp.nPoints());
211 //
212 // // Calculate average of field
213 // averageNeighbours
214 // (
215 // meshMover.mesh(),
216 // isMasterEdge,
217 // meshEdges,
218 // meshPoints,
219 // pp.edges(),
220 // invSumWeight,
221 // field,
222 // average
223 // );
224 //
225 // forAll(field, i)
226 // {
227 // if (field[i] >= fieldMin[i])
228 // {
229 // field[i] = (1-lambda)*field[i]+lambda*average[i];
230 // }
231 // }
232 //
233 //
234 // // Calculate average of field
235 // averageNeighbours
236 // (
237 // meshMover.mesh(),
238 // isMasterEdge,
239 // meshEdges,
240 // meshPoints,
241 // pp.edges(),
242 // invSumWeight,
243 // field,
244 // average
245 // );
246 //
247 // forAll(field, i)
248 // {
249 // if (field[i] >= fieldMin[i])
250 // {
251 // field[i] = (1-mu)*field[i]+mu*average[i];
252 // }
253 // }
254 //
255 //
256 // // Do residual calculation every so often.
257 // if ((iter % 10) == 0)
258 // {
259 // Info<< " Iteration " << iter << " residual "
260 // << gSum(mag(field-average))
261 // /returnReduce(average.size(), sumOp<label>())
262 // << endl;
263 // }
264 // }
265 //}
266 //XXXXXXXXX
267 
268 // Smooth normals on moving patch.
269 void Foam::snappyLayerDriver::smoothPatchNormals
270 (
271  const motionSmoother& meshMover,
272  const PackedBoolList& isMasterPoint,
273  const PackedBoolList& isMasterEdge,
274  const labelList& meshEdges,
275  const label nSmoothDisp,
276  pointField& normals
277 ) const
278 {
279  const indirectPrimitivePatch& pp = meshMover.patch();
280  const edgeList& edges = pp.edges();
281  const labelList& meshPoints = pp.meshPoints();
282 
283  // Calculate inverse sum of weights
284 
285  scalarField invSumWeight(pp.nPoints());
286  sumWeights
287  (
288  isMasterEdge,
289  meshEdges,
290  meshPoints,
291  edges,
292  invSumWeight
293  );
294 
295  // Get smoothly varying internal normals field.
296  Info<< "shrinkMeshDistance : Smoothing normals ..." << endl;
297 
298  for (label iter = 0; iter < nSmoothDisp; iter++)
299  {
300  vectorField average(pp.nPoints());
301  averageNeighbours
302  (
303  meshMover.mesh(),
304  isMasterEdge,
305  meshEdges,
306  meshPoints,
307  edges,
308  invSumWeight,
309  normals,
310  average
311  );
312 
313  // Do residual calculation every so often.
314  if ((iter % 10) == 0)
315  {
316  scalar resid = meshRefinement::gAverage
317  (
318  meshMover.mesh(),
319  isMasterPoint,
320  meshPoints,
321  mag(normals-average)()
322  );
323  Info<< " Iteration " << iter << " residual " << resid << endl;
324  }
325 
326  // Transfer to normals vector field
327  forAll(average, pointi)
328  {
329  // full smoothing neighbours + point value
330  average[pointi] = 0.5*(normals[pointi]+average[pointi]);
331  normals[pointi] = average[pointi];
332  normals[pointi] /= mag(normals[pointi]) + vSmall;
333  }
334  }
335 }
336 
337 
338 // Smooth normals in interior.
339 void Foam::snappyLayerDriver::smoothNormals
340 (
341  const label nSmoothDisp,
342  const PackedBoolList& isMasterPoint,
343  const PackedBoolList& isMasterEdge,
344  const labelList& fixedPoints,
345  pointVectorField& normals
346 ) const
347 {
348  // Get smoothly varying internal normals field.
349  Info<< "shrinkMeshDistance : Smoothing normals in interior ..." << endl;
350 
351  const fvMesh& mesh = meshRefiner_.mesh();
352  const edgeList& edges = mesh.edges();
353 
354  // Points that do not change.
355  PackedBoolList isFixedPoint(mesh.nPoints());
356 
357  // Internal points that are fixed
358  forAll(fixedPoints, i)
359  {
360  label meshPointi = fixedPoints[i];
361  isFixedPoint.set(meshPointi, 1);
362  }
363 
364  // Make sure that points that are coupled to meshPoints but not on a patch
365  // are fixed as well
366  syncTools::syncPointList(mesh, isFixedPoint, maxEqOp<unsigned int>(), 0);
367 
368 
369  // Correspondence between local edges/points and mesh edges/points
370  const labelList meshEdges(identityMap(mesh.nEdges()));
371  const labelList meshPoints(identityMap(mesh.nPoints()));
372 
373  // Calculate inverse sum of weights
374 
375  scalarField invSumWeight(mesh.nPoints(), 0);
376  sumWeights
377  (
378  isMasterEdge,
379  meshEdges,
380  meshPoints,
381  edges,
382  invSumWeight
383  );
384 
385  for (label iter = 0; iter < nSmoothDisp; iter++)
386  {
387  vectorField average(mesh.nPoints());
388  averageNeighbours
389  (
390  mesh,
391  isMasterEdge,
392  meshEdges,
393  meshPoints,
394  edges,
395  invSumWeight,
396  normals,
397  average
398  );
399 
400  // Do residual calculation every so often.
401  if ((iter % 10) == 0)
402  {
403  scalar resid = meshRefinement::gAverage
404  (
405  mesh,
406  isMasterPoint,
407  mag(normals-average)()
408  );
409  Info<< " Iteration " << iter << " residual " << resid << endl;
410  }
411 
412 
413  // Transfer to normals vector field
414  forAll(average, pointi)
415  {
416  if (isFixedPoint.get(pointi) == 0)
417  {
418  // full smoothing neighbours + point value
419  average[pointi] = 0.5*(normals[pointi]+average[pointi]);
420  normals[pointi] = average[pointi];
421  normals[pointi] /= mag(normals[pointi]) + vSmall;
422  }
423  }
424  }
425 }
426 
427 
428 // Tries and find a medial axis point. Done by comparing vectors to nearest
429 // wall point for both vertices of edge.
430 bool Foam::snappyLayerDriver::isMaxEdge
431 (
432  const List<pointData>& pointWallDist,
433  const label edgei,
434  const scalar minCos
435 ) const
436 {
437  const fvMesh& mesh = meshRefiner_.mesh();
438  const pointField& points = mesh.points();
439 
440  // Do not mark edges with one side on moving wall.
441 
442  const edge& e = mesh.edges()[edgei];
443 
444  vector v0(points[e[0]] - pointWallDist[e[0]].origin());
445  scalar magV0(mag(v0));
446 
447  if (magV0 < small)
448  {
449  return false;
450  }
451 
452  vector v1(points[e[1]] - pointWallDist[e[1]].origin());
453  scalar magV1(mag(v1));
454 
455  if (magV1 < small)
456  {
457  return false;
458  }
459 
460 
461  //- Detect based on vector to nearest point differing for both endpoints
462  // v0 /= magV0;
463  // v1 /= magV1;
464  //
466  // if ((v0 & v1) < minCos)
467  //{
468  // return true;
469  //}
470  // else
471  //{
472  // return false;
473  //}
474 
475  //- Detect based on extrusion vector differing for both endpoints
476  // the idea is that e.g. a sawtooth wall can still be extruded
477  // successfully as long as it is done all to the same direction.
478  if ((pointWallDist[e[0]].v() & pointWallDist[e[1]].v()) < minCos)
479  {
480  return true;
481  }
482  else
483  {
484  return false;
485  }
486 }
487 
488 
489 // Stop layer growth where mesh wraps around edge with a
490 // large feature angle
491 void Foam::snappyLayerDriver::handleFeatureAngleLayerTerminations
492 (
493  const scalar minCos,
494  const PackedBoolList& isMasterPoint,
495  const indirectPrimitivePatch& pp,
496  const labelList& meshEdges,
497  List<extrudeMode>& extrudeStatus,
498  pointField& patchDisp,
499  labelList& patchNLayers,
500  label& nPointCounter
501 ) const
502 {
503  const fvMesh& mesh = meshRefiner_.mesh();
504 
505  // Mark faces that have all points extruded
506  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
507 
508  boolList extrudedFaces(pp.size(), true);
509 
510  forAll(pp.localFaces(), facei)
511  {
512  const face& f = pp.localFaces()[facei];
513 
514  forAll(f, fp)
515  {
516  if (extrudeStatus[f[fp]] == NOEXTRUDE)
517  {
518  extrudedFaces[facei] = false;
519  break;
520  }
521  }
522  }
523 
524 
525 
526  // label nOldPointCounter = nPointCounter;
527 
528  // Detect situation where two featureedge-neighbouring faces are partly or
529  // not extruded and the edge itself is extruded. In this case unmark the
530  // edge for extrusion.
531 
532 
533  List<List<point>> edgeFaceNormals(pp.nEdges());
534  List<List<bool>> edgeFaceExtrude(pp.nEdges());
535 
536  const labelListList& edgeFaces = pp.edgeFaces();
537  const vectorField& faceNormals = pp.faceNormals();
538  const labelList& meshPoints = pp.meshPoints();
539 
540  forAll(edgeFaces, edgei)
541  {
542  const labelList& eFaces = edgeFaces[edgei];
543 
544  edgeFaceNormals[edgei].setSize(eFaces.size());
545  edgeFaceExtrude[edgei].setSize(eFaces.size());
546  forAll(eFaces, i)
547  {
548  label facei = eFaces[i];
549  edgeFaceNormals[edgei][i] = faceNormals[facei];
550  edgeFaceExtrude[edgei][i] = extrudedFaces[facei];
551  }
552  }
553 
555  (
556  mesh,
557  meshEdges,
558  edgeFaceNormals,
559  globalMeshData::ListPlusEqOp<List<point>>(), // combine operator
560  List<point>() // null value
561  );
562 
564  (
565  mesh,
566  meshEdges,
567  edgeFaceExtrude,
568  globalMeshData::ListPlusEqOp<List<bool>>(), // combine operator
569  List<bool>() // null value
570  );
571 
572 
573  forAll(edgeFaceNormals, edgei)
574  {
575  const List<point>& eFaceNormals = edgeFaceNormals[edgei];
576  const List<bool>& eFaceExtrude = edgeFaceExtrude[edgei];
577 
578  if (eFaceNormals.size() == 2)
579  {
580  const edge& e = pp.edges()[edgei];
581  label v0 = e[0];
582  label v1 = e[1];
583 
584  if
585  (
586  extrudeStatus[v0] != NOEXTRUDE
587  || extrudeStatus[v1] != NOEXTRUDE
588  )
589  {
590  if (!eFaceExtrude[0] || !eFaceExtrude[1])
591  {
592  const vector& n0 = eFaceNormals[0];
593  const vector& n1 = eFaceNormals[1];
594 
595  if ((n0 & n1) < minCos)
596  {
597  if
598  (
599  unmarkExtrusion
600  (
601  v0,
602  patchDisp,
603  patchNLayers,
604  extrudeStatus
605  )
606  )
607  {
608  if (isMasterPoint[meshPoints[v0]])
609  {
610  nPointCounter++;
611  }
612  }
613  if
614  (
615  unmarkExtrusion
616  (
617  v1,
618  patchDisp,
619  patchNLayers,
620  extrudeStatus
621  )
622  )
623  {
624  if (isMasterPoint[meshPoints[v1]])
625  {
626  nPointCounter++;
627  }
628  }
629  }
630  }
631  }
632  }
633  }
634 
635  // Info<< "Added "
636  // << returnReduce(nPointCounter-nOldPointCounter, sumOp<label>())
637  // << " point not to extrude." << endl;
638 }
639 
640 
641 // Find isolated islands (points, edges and faces and layer terminations)
642 // in the layer mesh and stop any layer growth at these points.
643 void Foam::snappyLayerDriver::findIsolatedRegions
644 (
645  const scalar minCosLayerTermination,
646  const PackedBoolList& isMasterPoint,
647  const PackedBoolList& isMasterEdge,
648  const indirectPrimitivePatch& pp,
649  const labelList& meshEdges,
650  const scalarField& minThickness,
651  List<extrudeMode>& extrudeStatus,
652  pointField& patchDisp,
653  labelList& patchNLayers
654 ) const
655 {
656  const fvMesh& mesh = meshRefiner_.mesh();
657 
658  Info<< "shrinkMeshDistance : Removing isolated regions ..." << endl;
659 
660  // Keep count of number of points unextruded
661  label nPointCounter = 0;
662 
663  while (true)
664  {
665  // Stop layer growth where mesh wraps around edge with a
666  // large feature angle
667  handleFeatureAngleLayerTerminations
668  (
669  minCosLayerTermination,
670  isMasterPoint,
671  pp,
672  meshEdges,
673 
674  extrudeStatus,
675  patchDisp,
676  patchNLayers,
677  nPointCounter
678  );
679 
680  syncPatchDisplacement
681  (
682  pp,
683  minThickness,
684  patchDisp,
685  patchNLayers,
686  extrudeStatus
687  );
688 
689 
690  // Do not extrude from point where all neighbouring
691  // faces are not grown
692  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
693  boolList extrudedFaces(pp.size(), true);
694  forAll(pp.localFaces(), facei)
695  {
696  const face& f = pp.localFaces()[facei];
697  forAll(f, fp)
698  {
699  if (extrudeStatus[f[fp]] == NOEXTRUDE)
700  {
701  extrudedFaces[facei] = false;
702  break;
703  }
704  }
705  }
706 
707  const labelListList& pointFaces = pp.pointFaces();
708 
709  boolList keptPoints(pp.nPoints(), false);
710  forAll(keptPoints, patchPointi)
711  {
712  const labelList& pFaces = pointFaces[patchPointi];
713 
714  forAll(pFaces, i)
715  {
716  label facei = pFaces[i];
717  if (extrudedFaces[facei])
718  {
719  keptPoints[patchPointi] = true;
720  break;
721  }
722  }
723  }
724 
726  (
727  mesh,
728  pp.meshPoints(),
729  keptPoints,
730  orEqOp<bool>(),
731  false // null value
732  );
733 
734  label nChanged = 0;
735 
736  forAll(keptPoints, patchPointi)
737  {
738  if (!keptPoints[patchPointi])
739  {
740  if
741  (
742  unmarkExtrusion
743  (
744  patchPointi,
745  patchDisp,
746  patchNLayers,
747  extrudeStatus
748  )
749  )
750  {
751  nPointCounter++;
752  nChanged++;
753  }
754  }
755  }
756 
757 
758  if (returnReduce(nChanged, sumOp<label>()) == 0)
759  {
760  break;
761  }
762  }
763 
764  const edgeList& edges = pp.edges();
765 
766 
767  // Count number of mesh edges using a point
768  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
769 
770  labelList isolatedPoint(pp.nPoints(),0);
771 
772  forAll(edges, edgei)
773  {
774  if (isMasterEdge.get(meshEdges[edgei]) == 1)
775  {
776  const edge& e = edges[edgei];
777 
778  label v0 = e[0];
779  label v1 = e[1];
780 
781  if (extrudeStatus[v1] != NOEXTRUDE)
782  {
783  isolatedPoint[v0] += 1;
784  }
785  if (extrudeStatus[v0] != NOEXTRUDE)
786  {
787  isolatedPoint[v1] += 1;
788  }
789  }
790  }
791 
793  (
794  mesh,
795  pp.meshPoints(),
796  isolatedPoint,
797  plusEqOp<label>(),
798  label(0) // null value
799  );
800 
801  // stop layer growth on isolated faces
802  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
803  forAll(pp, facei)
804  {
805  const face& f = pp.localFaces()[facei];
806  bool failed = false;
807  forAll(f, fp)
808  {
809  if (isolatedPoint[f[fp]] > 2)
810  {
811  failed = true;
812  break;
813  }
814  }
815  bool allPointsExtruded = true;
816  if (!failed)
817  {
818  forAll(f, fp)
819  {
820  if (extrudeStatus[f[fp]] == NOEXTRUDE)
821  {
822  allPointsExtruded = false;
823  break;
824  }
825  }
826 
827  if (allPointsExtruded)
828  {
829  forAll(f, fp)
830  {
831  if
832  (
833  unmarkExtrusion
834  (
835  f[fp],
836  patchDisp,
837  patchNLayers,
838  extrudeStatus
839  )
840  )
841  {
842  nPointCounter++;
843  }
844  }
845  }
846  }
847  }
848 
849  reduce(nPointCounter, sumOp<label>());
850  Info<< "Number isolated points extrusion stopped : "<< nPointCounter
851  << endl;
852 }
853 
854 
855 // Calculates medial axis fields:
856 // dispVec : normal of nearest wall. Where this normal changes direction
857 // defines the medial axis
858 // medialDist : distance to medial axis
859 // medialRatio : ratio of medial distance to wall distance.
860 // (1 at wall, 0 at medial axis)
861 void Foam::snappyLayerDriver::medialAxisSmoothingInfo
862 (
863  const motionSmoother& meshMover,
864  const label nSmoothNormals,
865  const label nSmoothSurfaceNormals,
866  const scalar minMedialAxisAngleCos,
867  const scalar featureAngle,
868 
869  pointVectorField& dispVec,
870  pointScalarField& medialRatio,
871  pointScalarField& medialDist,
872  pointVectorField& medialVec
873 ) const
874 {
875 
876  Info<< "medialAxisSmoothingInfo :"
877  << " Calculate distance to Medial Axis ..." << endl;
878 
879  const polyMesh& mesh = meshMover.mesh();
880  const pointField& points = mesh.points();
881 
882  const indirectPrimitivePatch& pp = meshMover.patch();
883  const labelList& meshPoints = pp.meshPoints();
884 
885  // Predetermine mesh edges
886  // ~~~~~~~~~~~~~~~~~~~~~~~
887 
888  // Precalculate master point/edge (only relevant for shared points/edges)
889  const PackedBoolList isMasterPoint(syncTools::getMasterPoints(mesh));
890  const PackedBoolList isMasterEdge(syncTools::getMasterEdges(mesh));
891  // Precalculate meshEdge per pp edge
892  const labelList meshEdges
893  (
894  pp.meshEdges
895  (
896  mesh.edges(),
897  mesh.pointEdges()
898  )
899  );
900 
901 
902  // Determine pointNormal
903  // ~~~~~~~~~~~~~~~~~~~~~
904 
905  pointField pointNormals(PatchTools::pointNormals(mesh, pp));
906 
907  // pointNormals
909  {
910  pointField meshPointNormals(mesh.nPoints(), point(1, 0, 0));
911  UIndirectList<point>(meshPointNormals, meshPoints) = pointNormals;
913  (
914  "pointNormals",
915  mesh,
916  meshPointNormals
917  );
918  }
919 
920  // Smooth patch normal vectors
921  smoothPatchNormals
922  (
923  meshMover,
924  isMasterPoint,
925  isMasterEdge,
926  meshEdges,
927  nSmoothSurfaceNormals,
928  pointNormals
929  );
930 
931  // smoothed pointNormals
933  {
934  pointField meshPointNormals(mesh.nPoints(), point(1, 0, 0));
935  UIndirectList<point>(meshPointNormals, meshPoints) = pointNormals;
937  (
938  "smoothed pointNormals",
939  mesh,
940  meshPointNormals
941  );
942  }
943 
944  // Calculate distance to pp points
945  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
946 
947  // Distance to wall
948  List<pointData> pointWallDist(mesh.nPoints());
949 
950  // Dummy additional info for PointEdgeWave
951  int dummyTrackData = 0;
952 
953 
954  // 1. Calculate distance to points where displacement is specified.
955  {
956  // Seed data.
957  List<pointData> wallInfo(meshPoints.size());
958 
959  forAll(meshPoints, patchPointi)
960  {
961  label pointi = meshPoints[patchPointi];
962  wallInfo[patchPointi] = pointData
963  (
964  points[pointi],
965  0.0,
966  pointi, // passive scalar
967  pointNormals[patchPointi] // surface normals
968  );
969  }
970 
971  // Do all calculations
972  List<pointData> edgeWallDist(mesh.nEdges());
973  PointEdgeWave<pointData> wallDistCalc
974  (
975  mesh,
976  meshPoints,
977  wallInfo,
978  pointWallDist,
979  edgeWallDist,
980  mesh.globalData().nTotalPoints(), // max iterations
981  dummyTrackData
982  );
983 
984 
985  label nUnvisit = returnReduce
986  (
987  wallDistCalc.getUnsetPoints(),
988  sumOp<label>()
989  );
990 
991  if (nUnvisit > 0)
992  {
994  << "Walking did not visit all points." << nl
995  << " Did not visit " << nUnvisit
996  << " out of " << mesh.globalData().nTotalPoints()
997  << " points. This is not necessarily a problem" << nl
998  << " and might be due to faceZones splitting of part"
999  << " of the domain." << nl << endl;
1000  }
1001  }
1002 
1003 
1004  // Check sync of wall distance
1006  {
1007  pointField origin(pointWallDist.size());
1008  scalarField distSqr(pointWallDist.size());
1009  // NA scalarField passiveS(pointWallDist.size());
1010  pointField passiveV(pointWallDist.size());
1011  forAll(pointWallDist, pointi)
1012  {
1013  origin[pointi] = pointWallDist[pointi].origin();
1014  distSqr[pointi] = pointWallDist[pointi].distSqr();
1015  // passiveS[pointi] = pointWallDist[pointi].s();
1016  passiveV[pointi] = pointWallDist[pointi].v();
1017  }
1018  meshRefinement::testSyncPointList("origin", mesh, origin);
1019  meshRefinement::testSyncPointList("distSqr", mesh, distSqr);
1020  // meshRefinement::testSyncPointList("passiveS", mesh, passiveS);
1021  meshRefinement::testSyncPointList("passiveV", mesh, passiveV);
1022  }
1023 
1024 
1025  // 2. Find points with max distance and transport information back to
1026  // wall.
1027  {
1028  List<pointData> pointMedialDist(mesh.nPoints());
1029  List<pointData> edgeMedialDist(mesh.nEdges());
1030 
1031  // Seed point data.
1032  DynamicList<pointData> maxInfo(meshPoints.size());
1033  DynamicList<label> maxPoints(meshPoints.size());
1034 
1035  // 1. Medial axis points
1036 
1037  const edgeList& edges = mesh.edges();
1038 
1039  forAll(edges, edgei)
1040  {
1041  const edge& e = edges[edgei];
1042 
1043  if
1044  (
1045  !pointWallDist[e[0]].valid(dummyTrackData)
1046  || !pointWallDist[e[1]].valid(dummyTrackData)
1047  )
1048  {
1049  // Unvisited point. See above about nUnvisit warning
1050  }
1051  else if (isMaxEdge(pointWallDist, edgei, minMedialAxisAngleCos))
1052  {
1053  // Both end points of edge have very different nearest wall
1054  // point. Mark both points as medial axis points.
1055 
1056  // Approximate medial axis location on edge.
1057  // const point medialAxisPt = e.centre(points);
1058  vector eVec = e.vec(points);
1059  scalar eMag = mag(eVec);
1060  if (eMag > vSmall)
1061  {
1062  eVec /= eMag;
1063 
1064  // Calculate distance along edge
1065  const point& p0 = points[e[0]];
1066  const point& p1 = points[e[1]];
1067  scalar dist0 = (p0-pointWallDist[e[0]].origin()) & eVec;
1068  scalar dist1 = (pointWallDist[e[1]].origin()-p1) & eVec;
1069  scalar s = 0.5*(dist1+eMag+dist0);
1070 
1071  point medialAxisPt;
1072  if (s <= dist0)
1073  {
1074  medialAxisPt = p0;
1075  }
1076  else if (s >= dist0+eMag)
1077  {
1078  medialAxisPt = p1;
1079  }
1080  else
1081  {
1082  medialAxisPt = p0+(s-dist0)*eVec;
1083  }
1084 
1085  forAll(e, ep)
1086  {
1087  label pointi = e[ep];
1088 
1089  if (!pointMedialDist[pointi].valid(dummyTrackData))
1090  {
1091  maxPoints.append(pointi);
1092  maxInfo.append
1093  (
1094  pointData
1095  (
1096  medialAxisPt, // points[pointi],
1097  magSqr(points[pointi]-medialAxisPt),//0.0,
1098  pointi, // passive data
1099  Zero // passive data
1100  )
1101  );
1102  pointMedialDist[pointi] = maxInfo.last();
1103  }
1104  }
1105  }
1106  }
1107  }
1108 
1109 
1110  // 2. Seed non-adapt patches
1111  const polyBoundaryMesh& patches = mesh.boundaryMesh();
1112 
1113  labelHashSet adaptPatches(meshMover.adaptPatchIndices());
1114 
1115 
1117  {
1118  const polyPatch& pp = patches[patchi];
1119  const pointPatchVectorField& pvf =
1120  meshMover.displacement().boundaryField()[patchi];
1121 
1122  if
1123  (
1124  !pp.coupled()
1125  && !isA<emptyPolyPatch>(pp)
1126  && !adaptPatches.found(patchi)
1127  )
1128  {
1129  const labelList& meshPoints = pp.meshPoints();
1130 
1131  // Check the type of the patchField. The types are
1132  // - fixedValue (0 or more layers) but the >0 layers have
1133  // already been handled in the adaptPatches loop
1134  // - constraint (but not coupled) types, e.g. symmetryPlane,
1135  // slip.
1136  if (pvf.fixesValue())
1137  {
1138  // Disable all movement on fixedValue patchFields
1139  Info<< "Inserting all points on patch " << pp.name()
1140  << endl;
1141 
1142  forAll(meshPoints, i)
1143  {
1144  label pointi = meshPoints[i];
1145  if (!pointMedialDist[pointi].valid(dummyTrackData))
1146  {
1147  maxPoints.append(pointi);
1148  maxInfo.append
1149  (
1150  pointData
1151  (
1152  points[pointi],
1153  0.0,
1154  pointi, // passive data
1155  Zero // passive data
1156  )
1157  );
1158  pointMedialDist[pointi] = maxInfo.last();
1159  }
1160  }
1161  }
1162  else
1163  {
1164  // Based on geometry: analyse angle w.r.t. nearest moving
1165  // point. In the pointWallDist we transported the
1166  // normal as the passive vector. Note that this points
1167  // out of the originating wall so inside of the domain
1168  // on this patch.
1169  Info<< "Inserting points on patch " << pp.name()
1170  << " if angle to nearest layer patch > "
1171  << radToDeg(featureAngle) << " degrees." << endl;
1172 
1173  scalar featureAngleCos = Foam::cos(featureAngle);
1174  pointField pointNormals(PatchTools::pointNormals(mesh, pp));
1175 
1176  forAll(meshPoints, i)
1177  {
1178  label pointi = meshPoints[i];
1179 
1180  if
1181  (
1182  pointWallDist[pointi].valid(dummyTrackData)
1183  && !pointMedialDist[pointi].valid(dummyTrackData)
1184  )
1185  {
1186  // Check if angle not too large.
1187  scalar cosAngle =
1188  (
1189  -pointWallDist[pointi].v()
1190  & pointNormals[i]
1191  );
1192  if (cosAngle > featureAngleCos)
1193  {
1194  // Extrusion direction practically perpendicular
1195  // to the patch. Disable movement at the patch.
1196 
1197  maxPoints.append(pointi);
1198  maxInfo.append
1199  (
1200  pointData
1201  (
1202  points[pointi],
1203  0.0,
1204  pointi, // passive data
1205  Zero // passive data
1206  )
1207  );
1208  pointMedialDist[pointi] = maxInfo.last();
1209  }
1210  else
1211  {
1212  // Extrusion direction makes angle with patch
1213  // so allow sliding - don't insert zero points
1214  }
1215  }
1216  }
1217  }
1218  }
1219  }
1220 
1221  maxInfo.shrink();
1222  maxPoints.shrink();
1223 
1224  // Do all calculations
1225  PointEdgeWave<pointData> medialDistCalc
1226  (
1227  mesh,
1228  maxPoints,
1229  maxInfo,
1230 
1231  pointMedialDist,
1232  edgeMedialDist,
1233  mesh.globalData().nTotalPoints(), // max iterations
1234  dummyTrackData
1235  );
1236 
1237  // Extract medial axis distance as pointScalarField
1238  forAll(pointMedialDist, pointi)
1239  {
1240  medialDist[pointi] = Foam::sqrt(pointMedialDist[pointi].distSqr());
1241  medialVec[pointi] = pointMedialDist[pointi].origin();
1242  }
1243 
1244  // Check
1246  {
1247  meshRefinement::testSyncPointList("medialDist", mesh, medialDist);
1248  meshRefinement::testSyncPointList("medialVec", mesh, medialVec);
1249  }
1250  }
1251 
1252  // Extract transported surface normals as pointVectorField
1253  forAll(dispVec, i)
1254  {
1255  if (!pointWallDist[i].valid(dummyTrackData))
1256  {
1257  dispVec[i] = vector(1, 0, 0);
1258  }
1259  else
1260  {
1261  dispVec[i] = pointWallDist[i].v();
1262  }
1263  }
1264 
1265  // Smooth normal vectors. Do not change normals on pp.meshPoints
1266  smoothNormals
1267  (
1268  nSmoothNormals,
1269  isMasterPoint,
1270  isMasterEdge,
1271  meshPoints,
1272  dispVec
1273  );
1274 
1276  {
1277  meshRefinement::testSyncPointList("smoothed dispVec", mesh, dispVec);
1278  }
1279 
1280  // Calculate ratio point medial distance to point wall distance
1281  forAll(medialRatio, pointi)
1282  {
1283  if (!pointWallDist[pointi].valid(dummyTrackData))
1284  {
1285  medialRatio[pointi] = 0.0;
1286  }
1287  else
1288  {
1289  scalar wDist2 = pointWallDist[pointi].distSqr();
1290  scalar mDist = medialDist[pointi];
1291 
1292  if (wDist2 < sqr(small) && mDist < small)
1293  //- Note: maybe less strict:
1294  //(
1295  // wDist2 < sqr(meshRefiner_.mergeDistance())
1296  // && mDist < meshRefiner_.mergeDistance()
1297  //)
1298  {
1299  medialRatio[pointi] = 0.0;
1300  }
1301  else
1302  {
1303  medialRatio[pointi] = mDist / (Foam::sqrt(wDist2) + mDist);
1304  }
1305  }
1306  }
1307 
1308 
1310  {
1311  // medialRatio
1312  meshRefinement::testSyncPointList("medialRatio", mesh, medialRatio);
1313  }
1314 
1315 
1317  {
1318  Info<< "medialAxisSmoothingInfo :"
1319  << " Writing:" << nl
1320  << " " << dispVec.name()
1321  << " : normalised direction of nearest displacement" << nl
1322  << " " << medialDist.name()
1323  << " : distance to medial axis" << nl
1324  << " " << medialVec.name()
1325  << " : nearest point on medial axis" << nl
1326  << " " << medialRatio.name()
1327  << " : ratio of medial distance to wall distance" << nl
1328  << endl;
1329  meshRefiner_.mesh().setInstance(meshRefiner_.name());
1330  meshRefiner_.write
1331  (
1334  (
1337  ),
1338  mesh.time().path()/meshRefiner_.name()
1339  );
1340  dispVec.write();
1341  medialDist.write();
1342  medialVec.write();
1343  medialRatio.write();
1344  }
1345 }
1346 
1347 
1348 void Foam::snappyLayerDriver::shrinkMeshMedialDistance
1349 (
1350  motionSmoother& meshMover,
1351  const dictionary& meshQualityDict,
1352  const List<labelPair>& baffles,
1353  const label nSmoothPatchThickness,
1354  const label nSmoothDisplacement,
1355  const scalar maxThicknessToMedialRatio,
1356  const label nAllowableErrors,
1357  const label nSnap,
1358  const scalar minCosLayerTermination,
1359 
1360  const scalarField& layerThickness,
1361  const scalarField& minThickness,
1362 
1363  const pointVectorField& dispVec,
1364  const pointScalarField& medialRatio,
1365  const pointScalarField& medialDist,
1366  const pointVectorField& medialVec,
1367 
1368  List<extrudeMode>& extrudeStatus,
1369  pointField& patchDisp,
1370  labelList& patchNLayers
1371 ) const
1372 {
1373  Info<< "shrinkMeshMedialDistance : Smoothing using Medial Axis ..." << endl;
1374 
1375  const polyMesh& mesh = meshMover.mesh();
1376 
1377  const indirectPrimitivePatch& pp = meshMover.patch();
1378  const labelList& meshPoints = pp.meshPoints();
1379 
1380  // Precalculate master points/edge (only relevant for shared points/edges)
1381  const PackedBoolList isMasterPoint(syncTools::getMasterPoints(mesh));
1382  const PackedBoolList isMasterEdge(syncTools::getMasterEdges(mesh));
1383  // Precalculate meshEdge per pp edge
1384  const labelList meshEdges
1385  (
1386  pp.meshEdges
1387  (
1388  mesh.edges(),
1389  mesh.pointEdges()
1390  )
1391  );
1392 
1393 
1394  scalarField thickness(layerThickness.size());
1395 
1396  thickness = mag(patchDisp);
1397 
1398  forAll(thickness, patchPointi)
1399  {
1400  if (extrudeStatus[patchPointi] == NOEXTRUDE)
1401  {
1402  thickness[patchPointi] = 0.0;
1403  }
1404  }
1405 
1406  label numThicknessRatioExclude = 0;
1407 
1408  // reduce thickness where thickness/medial axis distance large
1409  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1410 
1411  autoPtr<OBJstream> str;
1412  if (debug)
1413  {
1414  str.reset
1415  (
1416  new OBJstream
1417  (
1418  mesh.time().path()
1419  / "thicknessRatioExcludePoints_"
1420  + meshRefiner_.name()
1421  + ".obj"
1422  )
1423  );
1424  Info<< "Writing points with too large an extrusion distance to "
1425  << str().name() << endl;
1426  }
1427 
1428  autoPtr<OBJstream> medialVecStr;
1429  if (debug)
1430  {
1431  medialVecStr.reset
1432  (
1433  new OBJstream
1434  (
1435  mesh.time().path()
1436  / "thicknessRatioExcludeMedialVec_"
1437  + meshRefiner_.name()
1438  + ".obj"
1439  )
1440  );
1441  Info<< "Writing points with too large an extrusion distance to "
1442  << medialVecStr().name() << endl;
1443  }
1444 
1445  forAll(meshPoints, patchPointi)
1446  {
1447  if (extrudeStatus[patchPointi] != NOEXTRUDE)
1448  {
1449  label pointi = meshPoints[patchPointi];
1450 
1451  //- Option 1: look only at extrusion thickness v.s. distance
1452  // to nearest (medial axis or static) point.
1453  scalar mDist = medialDist[pointi];
1454  scalar thicknessRatio = thickness[patchPointi]/(mDist+vSmall);
1455 
1456  //- Option 2: Look at component in the direction
1457  // of nearest (medial axis or static) point.
1458  vector n =
1459  patchDisp[patchPointi]
1460  / (mag(patchDisp[patchPointi]) + vSmall);
1461  vector mVec = mesh.points()[pointi]-medialVec[pointi];
1462  mVec /= mag(mVec)+vSmall;
1463  thicknessRatio *= (n&mVec);
1464 
1465  if (thicknessRatio > maxThicknessToMedialRatio)
1466  {
1467  // Truncate thickness.
1468  if (debug)
1469  {
1470  Pout<< "truncating displacement at "
1471  << mesh.points()[pointi]
1472  << " from " << thickness[patchPointi]
1473  << " to "
1474  << 0.5
1475  *(
1476  minThickness[patchPointi]
1477  +thickness[patchPointi]
1478  )
1479  << " medial direction:" << mVec
1480  << " extrusion direction:" << n
1481  << " with thicknessRatio:" << thicknessRatio
1482  << endl;
1483  }
1484 
1485  thickness[patchPointi] =
1486  0.5*(minThickness[patchPointi]+thickness[patchPointi]);
1487 
1488  patchDisp[patchPointi] = thickness[patchPointi]*n;
1489 
1490  if (isMasterPoint[pointi])
1491  {
1492  numThicknessRatioExclude++;
1493  }
1494 
1495  if (str.valid())
1496  {
1497  const point& pt = mesh.points()[pointi];
1498  str().write(linePointRef(pt, pt+patchDisp[patchPointi]));
1499  }
1500  if (medialVecStr.valid())
1501  {
1502  const point& pt = mesh.points()[pointi];
1503  medialVecStr().write
1504  (
1505  linePointRef
1506  (
1507  pt,
1508  medialVec[pointi]
1509  )
1510  );
1511  }
1512  }
1513  }
1514  }
1515 
1516  reduce(numThicknessRatioExclude, sumOp<label>());
1517  Info<< "shrinkMeshMedialDistance : Reduce layer thickness at "
1518  << numThicknessRatioExclude
1519  << " nodes where thickness to medial axis distance is large " << endl;
1520 
1521 
1522  // find points where layer growth isolated to a lone point, edge or face
1523 
1524  findIsolatedRegions
1525  (
1526  minCosLayerTermination,
1527  isMasterPoint,
1528  isMasterEdge,
1529  pp,
1530  meshEdges,
1531  minThickness,
1532 
1533  extrudeStatus,
1534  patchDisp,
1535  patchNLayers
1536  );
1537 
1538  // Update thickness for changed extrusion
1539  forAll(thickness, patchPointi)
1540  {
1541  if (extrudeStatus[patchPointi] == NOEXTRUDE)
1542  {
1543  thickness[patchPointi] = 0.0;
1544  }
1545  }
1546 
1547 
1548  // smooth layer thickness on moving patch
1549  smoothField
1550  (
1551  meshMover,
1552  isMasterPoint,
1553  isMasterEdge,
1554  meshEdges,
1555  minThickness,
1556  nSmoothPatchThickness,
1557 
1558  thickness
1559  );
1560 
1561  // Dummy additional info for PointEdgeWave
1562  int dummyTrackData = 0;
1563 
1564  List<pointData> pointWallDist(mesh.nPoints());
1565 
1566  const pointField& points = mesh.points();
1567  // 1. Calculate distance to points where displacement is specified.
1568  // This wave is used to transport layer thickness
1569  {
1570  // Distance to wall and medial axis on edges.
1571  List<pointData> edgeWallDist(mesh.nEdges());
1572  labelList wallPoints(meshPoints.size());
1573 
1574  // Seed data.
1575  List<pointData> wallInfo(meshPoints.size());
1576 
1577  forAll(meshPoints, patchPointi)
1578  {
1579  label pointi = meshPoints[patchPointi];
1580  wallPoints[patchPointi] = pointi;
1581  wallInfo[patchPointi] = pointData
1582  (
1583  points[pointi],
1584  0.0,
1585  thickness[patchPointi], // transport layer thickness
1586  Zero // passive vector
1587  );
1588  }
1589 
1590  // Do all calculations
1591  PointEdgeWave<pointData> wallDistCalc
1592  (
1593  mesh,
1594  wallPoints,
1595  wallInfo,
1596  pointWallDist,
1597  edgeWallDist,
1598  mesh.globalData().nTotalPoints(), // max iterations
1599  dummyTrackData
1600  );
1601  }
1602 
1603  // Calculate scaled displacement vector
1604  pointVectorField& displacement = meshMover.displacement();
1605 
1606  forAll(displacement, pointi)
1607  {
1608  if (!pointWallDist[pointi].valid(dummyTrackData))
1609  {
1610  displacement[pointi] = Zero;
1611  }
1612  else
1613  {
1614  // 1) displacement on nearest wall point, scaled by medialRatio
1615  // (wall distance / medial distance)
1616  // 2) pointWallDist[pointi].s() is layer thickness transported
1617  // from closest wall point.
1618  // 3) shrink in opposite direction of addedPoints
1619  displacement[pointi] =
1620  -medialRatio[pointi]
1621  *pointWallDist[pointi].s()
1622  *dispVec[pointi];
1623  }
1624  }
1625 
1626 
1627 
1628  // Smear displacement away from fixed values (medialRatio=0 or 1)
1629  if (nSmoothDisplacement > 0)
1630  {
1631  scalarField invSumWeight(mesh.nPoints());
1632  sumWeights
1633  (
1634  isMasterEdge,
1635  identityMap(mesh.nEdges()),
1636  identityMap(mesh.nPoints()),
1637  mesh.edges(),
1638  invSumWeight
1639  );
1640 
1641 
1642  // Get smoothly varying patch field.
1643  Info<< "shrinkMeshDistance : Smoothing displacement ..." << endl;
1644 
1645  const scalar lambda = 0.33;
1646  const scalar mu = -0.34;
1647 
1648  pointField average(mesh.nPoints());
1649  for (label iter = 0; iter < nSmoothDisplacement; iter++)
1650  {
1651  // Calculate average of field
1652  averageNeighbours
1653  (
1654  mesh,
1655  isMasterEdge,
1656  identityMap(mesh.nEdges()), // meshEdges,
1657  identityMap(mesh.nPoints()), // meshPoints,
1658  mesh.edges(), // edges,
1659  invSumWeight,
1660  displacement,
1661  average
1662  );
1663 
1664  forAll(displacement, i)
1665  {
1666  if (medialRatio[i] > small && medialRatio[i] < 1-small)
1667  {
1668  displacement[i] =
1669  (1-lambda)*displacement[i]
1670  +lambda*average[i];
1671  }
1672  }
1673 
1674 
1675  // Calculate average of field
1676  averageNeighbours
1677  (
1678  mesh,
1679  isMasterEdge,
1680  identityMap(mesh.nEdges()), // meshEdges,
1681  identityMap(mesh.nPoints()), // meshPoints,
1682  mesh.edges(), // edges,
1683  invSumWeight,
1684  displacement,
1685  average
1686  );
1687 
1688  forAll(displacement, i)
1689  {
1690  if (medialRatio[i] > small && medialRatio[i] < 1-small)
1691  {
1692  displacement[i] = (1-mu)*displacement[i]+mu*average[i];
1693  }
1694  }
1695 
1696 
1697  // Do residual calculation every so often.
1698  if ((iter % 10) == 0)
1699  {
1700  scalar resid = meshRefinement::gAverage
1701  (
1702  mesh,
1704  mag(displacement-average)()
1705  );
1706  Info<< " Iteration " << iter << " residual " << resid
1707  << endl;
1708  }
1709  }
1710  }
1711 
1712  // Make sure displacement boundary conditions is up-to-date with
1713  // internal field
1714  meshMover.setDisplacementPatchFields();
1715 
1716 
1717  // Check a bit the sync of displacements
1719  {
1720  // initial mesh
1721  meshRefinement::testSyncPointList("mesh.points()", mesh, mesh.points());
1722 
1723  // pointWallDist
1724  scalarField pWallDist(pointWallDist.size());
1725  forAll(pointWallDist, pointi)
1726  {
1727  pWallDist[pointi] = pointWallDist[pointi].s();
1728  }
1729  meshRefinement::testSyncPointList("pointWallDist", mesh, pWallDist);
1730 
1731  // dispVec
1732  meshRefinement::testSyncPointList("dispVec", mesh, dispVec);
1733 
1734  // displacement before and after correction
1736  (
1737  "displacement BEFORE",
1738  mesh,
1739  displacement
1740  );
1741 
1742  meshMover.correctBoundaryConditions(displacement);
1744  (
1745  "displacement AFTER",
1746  mesh,
1747  displacement
1748  );
1749  }
1750 
1751 
1753  {
1754  const_cast<Time&>(mesh.time())++;
1755  Info<< "Writing wanted-displacement mesh (possibly illegal) to "
1756  << meshRefiner_.name() << endl;
1757  pointField oldPoints(mesh.points());
1758 
1759  meshRefiner_.mesh().movePoints(meshMover.curPoints());
1760  // Warn meshMover for changed geometry
1761  meshMover.movePoints();
1762 
1763  // Above move will have changed the instance only on the points (which
1764  // is correct).
1765  // However the previous mesh written will be the mesh with layers
1766  // (see snappyLayerDriver.C) so we now have to force writing all files
1767  // so we can easily step through time steps. Note that if you
1768  // don't write the mesh with layers this is not necessary.
1769  meshRefiner_.mesh().setInstance(meshRefiner_.name());
1770 
1771  meshRefiner_.write
1772  (
1775  (
1778  ),
1779  mesh.time().path()/meshRefiner_.name()
1780  );
1781  dispVec.write();
1782  medialDist.write();
1783  medialRatio.write();
1784 
1785  // Move mesh back
1786  meshRefiner_.mesh().movePoints(oldPoints);
1787  // Warn meshMover for changed geometry
1788  meshMover.movePoints();
1789  }
1790 
1791 
1792  // Current faces to check. Gets modified in meshMover.scaleMesh
1793  labelList checkFaces(identityMap(mesh.nFaces()));
1794 
1795  Info<< "shrinkMeshMedialDistance : Moving mesh ..." << endl;
1796  scalar oldErrorReduction = -1;
1797 
1798  for (label iter = 0; iter < 2*nSnap ; iter++)
1799  {
1800  Info<< "Iteration " << iter << endl;
1801  if (iter == nSnap)
1802  {
1803  Info<< "Displacement scaling for error reduction set to 0." << endl;
1804  oldErrorReduction = meshMover.setErrorReduction(0.0);
1805  }
1806 
1807  if
1808  (
1809  meshMover.scaleMesh
1810  (
1811  checkFaces,
1812  baffles,
1813  meshMover.paramDict(),
1814  meshQualityDict,
1815  true,
1816  nAllowableErrors
1817  )
1818  )
1819  {
1820  Info<< "shrinkMeshMedialDistance : Successfully moved mesh" << endl;
1821  break;
1822  }
1823  }
1824 
1825  if (oldErrorReduction >= 0)
1826  {
1827  meshMover.setErrorReduction(oldErrorReduction);
1828  }
1829 
1830  Info<< "shrinkMeshMedialDistance : Finished moving mesh ..." << endl;
1831 }
1832 
1833 
1834 // ************************************************************************* //
label n
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
void setSize(const label)
Reset size of List.
Definition: List.C:281
virtual Ostream & write(const char)=0
Write character.
static tmp< pointField > pointNormals(const polyMesh &, const PrimitivePatch< FaceList, PointField > &)
Return parallel consistent point normals for patches using mesh points.
static void testSyncPointList(const string &msg, const polyMesh &mesh, const List< scalar > &fld)
const fvMesh & mesh() const
Reference to mesh.
static T gAverage(const PackedBoolList &isMasterElem, const UList< T > &values)
Helper: calculate average.
static writeType writeLevel()
Get/set write level.
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1313
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
volScalarField scalarField(fieldObject, mesh)
label patchi
const pointField & points
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.
static const zero Zero
Definition: zero.H:97
const doubleScalar e
Definition: doubleScalar.H:106
scalar radToDeg(const scalar rad)
Convert radians to degrees.
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:257
line< point, const point & > linePointRef
Line using referred points.
Definition: linePointRef.H:45
messageStream Info
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &df)
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.
PointField< scalar > pointScalarField
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
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 > &)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
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)
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