snappySnapDriver.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2015 OpenFOAM Foundation
9  Copyright (C) 2015-2022 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Description
28  All to do with snapping to the surface
29 
30 \*----------------------------------------------------------------------------*/
31 
32 #include "snappySnapDriver.H"
33 #include "motionSmoother.H"
34 #include "polyTopoChange.H"
35 #include "syncTools.H"
36 #include "fvMesh.H"
37 #include "Time.H"
38 #include "OFstream.H"
39 #include "OBJstream.H"
40 #include "mapPolyMesh.H"
41 #include "pointEdgePoint.H"
42 #include "PointEdgeWave.H"
43 #include "mergePoints.H"
44 #include "snapParameters.H"
45 #include "refinementSurfaces.H"
46 #include "searchableSurfaces.H"
47 #include "unitConversion.H"
48 #include "localPointRegion.H"
49 #include "PatchTools.H"
50 #include "refinementFeatures.H"
51 #include "weightedPosition.H"
52 #include "profiling.H"
53 
54 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
55 
56 namespace Foam
57 {
58 
59 defineTypeNameAndDebug(snappySnapDriver, 0);
60 
61 } // End namespace Foam
62 
63 
64 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
65 
66 // Calculate geometrically collocated points, Requires bitSet to be
67 // sized and initialised!
68 Foam::label Foam::snappySnapDriver::getCollocatedPoints
69 (
70  const scalar tol,
71  const pointField& points,
72  bitSet& isCollocatedPoint
73 )
74 {
75  labelList pointMap;
76  label nUnique = Foam::mergePoints
77  (
78  points, // points
79  tol, // mergeTol
80  false, // verbose
81  pointMap
82  );
83  bool hasMerged = (nUnique < points.size());
84 
85  if (!returnReduceOr(hasMerged))
86  {
87  return 0;
88  }
89 
90  // Determine which merged points are referenced more than once
91  label nCollocated = 0;
92 
93  // Per old point the newPoint. Or -1 (not set yet) or -2 (already seen
94  // twice)
95  labelList firstOldPoint(nUnique, -1);
96  forAll(pointMap, oldPointi)
97  {
98  label newPointi = pointMap[oldPointi];
99 
100  if (firstOldPoint[newPointi] == -1)
101  {
102  // First use of oldPointi. Store.
103  firstOldPoint[newPointi] = oldPointi;
104  }
105  else if (firstOldPoint[newPointi] == -2)
106  {
107  // Third or more reference of oldPointi -> non-manifold
108  isCollocatedPoint.set(oldPointi);
109  nCollocated++;
110  }
111  else
112  {
113  // Second reference of oldPointi -> non-manifold
114  isCollocatedPoint.set(firstOldPoint[newPointi]);
115  nCollocated++;
116 
117  isCollocatedPoint.set(oldPointi);
118  nCollocated++;
119 
120  // Mark with special value to save checking next time round
121  firstOldPoint[newPointi] = -2;
122  }
123  }
124  return returnReduce(nCollocated, sumOp<label>());
125 }
126 
127 
128 Foam::tmp<Foam::pointField> Foam::snappySnapDriver::smoothInternalDisplacement
129 (
130  const meshRefinement& meshRefiner,
131  const motionSmoother& meshMover
132 )
133 {
134  const indirectPrimitivePatch& pp = meshMover.patch();
135  const polyMesh& mesh = meshMover.mesh();
136 
137  // Get neighbour refinement
138  const hexRef8& cutter = meshRefiner.meshCutter();
139  const labelList& cellLevel = cutter.cellLevel();
140 
141 
142  // Get the faces on the boundary
143  bitSet isFront(mesh.nFaces(), pp.addressing());
144 
145  // Walk out from the surface a bit. Poor man's FaceCellWave.
146  // Commented out for now - not sure if needed and if so how much
147  //for (label iter = 0; iter < 2; iter++)
148  //{
149  // bitSet newIsFront(mesh.nFaces());
150  //
151  // forAll(isFront, facei)
152  // {
153  // if (isFront.test(facei))
154  // {
155  // label own = mesh.faceOwner()[facei];
156  // const cell& ownFaces = mesh.cells()[own];
157  // newIsFront.set(ownFaces);
158  //
159  // if (mesh.isInternalFace(facei))
160  // {
161  // label nei = mesh.faceNeighbour()[facei];
162  // const cell& neiFaces = mesh.cells()[nei];
163  // newIsFront.set(neiFaces);
164  // }
165  // }
166  // }
167  //
168  // syncTools::syncFaceList
169  // (
170  // mesh,
171  // newIsFront,
172  // orEqOp<unsigned int>()
173  // );
174  //
175  // isFront = newIsFront;
176  //}
177 
178  // Mark all points on faces
179  // - not on the boundary
180  // - inbetween differing refinement levels
181  bitSet isMovingPoint(mesh.nPoints());
182 
183  label nInterface = 0;
184 
185  for (label facei = 0; facei < mesh.nInternalFaces(); facei++)
186  {
187  label ownLevel = cellLevel[mesh.faceOwner()[facei]];
188  label neiLevel = cellLevel[mesh.faceNeighbour()[facei]];
189 
190  if (!isFront.test(facei) && ownLevel != neiLevel)
191  {
192  const face& f = mesh.faces()[facei];
193  isMovingPoint.set(f);
194 
195  ++nInterface;
196  }
197  }
198 
199  labelList neiCellLevel;
200  syncTools::swapBoundaryCellList(mesh, cellLevel, neiCellLevel);
201 
202  for (label facei = mesh.nInternalFaces(); facei < mesh.nFaces(); facei++)
203  {
204  label ownLevel = cellLevel[mesh.faceOwner()[facei]];
205  label neiLevel = neiCellLevel[facei-mesh.nInternalFaces()];
206 
207  if (!isFront.test(facei) && ownLevel != neiLevel)
208  {
209  const face& f = mesh.faces()[facei];
210  isMovingPoint.set(f);
211 
212  ++nInterface;
213  }
214  }
215 
216  if (debug)
217  {
218  Info<< "Found " << returnReduce(nInterface, sumOp<label>())
219  << " faces out of " << mesh.globalData().nTotalFaces()
220  << " inbetween refinement regions." << endl;
221  }
222 
223  // Make sure that points that are coupled to a moving point are marked
224  // as well
225  syncTools::syncPointList(mesh, isMovingPoint, maxEqOp<unsigned int>(), 0);
226 
227  // Unmark any point on the boundary. If we're doing zero iterations of
228  // face-cell wave we might have coupled points not being unmarked.
229  isMovingPoint.unset(pp.meshPoints());
230 
231  // Make sure that points that are coupled to meshPoints but not on a patch
232  // are unmarked as well
233  syncTools::syncPointList(mesh, isMovingPoint, minEqOp<unsigned int>(), 1);
234 
235 
236  // Calculate average of connected cells
237  Field<weightedPosition> sumLocation
238  (
239  mesh.nPoints(),
241  );
242 
243  forAll(isMovingPoint, pointi)
244  {
245  if (isMovingPoint.test(pointi))
246  {
247  const labelList& pCells = mesh.pointCells(pointi);
248 
249  sumLocation[pointi].first() = pCells.size();
250  for (const label celli : pCells)
251  {
252  sumLocation[pointi].second() += mesh.cellCentres()[celli];
253  }
254  }
255  }
256 
257  // Add coupled contributions
258  weightedPosition::syncPoints(mesh, sumLocation);
259 
260  tmp<pointField> tdisplacement(new pointField(mesh.nPoints(), Zero));
261  pointField& displacement = tdisplacement.ref();
262 
263  label nAdapted = 0;
264 
265  forAll(displacement, pointi)
266  {
267  const weightedPosition& wp = sumLocation[pointi];
268  if (mag(wp.first()) > VSMALL)
269  {
270  displacement[pointi] =
271  wp.second()/wp.first()
272  - mesh.points()[pointi];
273  nAdapted++;
274  }
275  }
276 
277  Info<< "Smoothing " << returnReduce(nAdapted, sumOp<label>())
278  << " points inbetween refinement regions."
279  << endl;
280 
281  return tdisplacement;
282 }
283 
284 
285 // Calculate displacement as average of patch points.
286 Foam::tmp<Foam::pointField> Foam::snappySnapDriver::smoothPatchDisplacement
287 (
288  const motionSmoother& meshMover,
289  const List<labelPair>& baffles
290 )
291 {
292  const indirectPrimitivePatch& pp = meshMover.patch();
293 
294  // Calculate geometrically non-manifold points on the patch to be moved.
295  bitSet nonManifoldPoint(pp.nPoints());
296  label nNonManifoldPoints = getCollocatedPoints
297  (
298  SMALL,
299  pp.localPoints(),
300  nonManifoldPoint
301  );
302  Info<< "Found " << nNonManifoldPoints << " non-manifold point(s)."
303  << endl;
304 
305 
306  // Average points
307  // ~~~~~~~~~~~~~~
308 
309  // We determine three points:
310  // - average of (centres of) connected patch faces
311  // - average of (centres of) connected internal mesh faces
312  // - as fallback: centre of any connected cell
313  // so we can do something moderately sensible for non/manifold points.
314 
315  // Note: the averages are calculated properly parallel. This is
316  // necessary to get the points shared by processors correct.
317 
318 
319  const labelListList& pointFaces = pp.pointFaces();
320  const labelList& meshPoints = pp.meshPoints();
321  const pointField& points = pp.points();
322  const polyMesh& mesh = meshMover.mesh();
323 
324  // Get labels of faces to count (master of coupled faces and baffle pairs)
325  bitSet isMasterFace(syncTools::getMasterFaces(mesh));
326 
327  {
328  forAll(baffles, i)
329  {
330  label f0 = baffles[i].first();
331  label f1 = baffles[i].second();
332 
333  if (isMasterFace.test(f0))
334  {
335  // Make f1 a slave
336  isMasterFace.unset(f1);
337  }
338  else if (isMasterFace.test(f1))
339  {
340  isMasterFace.unset(f0);
341  }
342  else
343  {
345  << "Both sides of baffle consisting of faces " << f0
346  << " and " << f1 << " are already slave faces."
347  << abort(FatalError);
348  }
349  }
350  }
351 
352 
353  // Get average position of boundary face centres
354  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
355 
356  Field<weightedPosition> avgBoundary
357  (
358  pointFaces.size(),
360  );
361  {
362  forAll(pointFaces, patchPointi)
363  {
364  const labelList& pFaces = pointFaces[patchPointi];
365 
366  forAll(pFaces, pfi)
367  {
368  label facei = pFaces[pfi];
369 
370  if (isMasterFace.test(pp.addressing()[facei]))
371  {
372  avgBoundary[patchPointi].first() += 1.0;
373  avgBoundary[patchPointi].second() +=
374  pp[facei].centre(points);
375  }
376  }
377  }
378 
379  // Add coupled contributions
381 
382  // Normalise
383  forAll(avgBoundary, i)
384  {
385  // Note: what if there is no master boundary face?
386  if (mag(avgBoundary[i].first()) > VSMALL)
387  {
388  avgBoundary[i].second() /= avgBoundary[i].first();
389  }
390  }
391  }
392 
393 
394  // Get average position of internal face centres
395  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
396 
397  Field<weightedPosition> avgInternal;
398  {
399  Field<weightedPosition> globalSum
400  (
401  mesh.nPoints(),
403  );
404 
405  // Note: no use of pointFaces
406  const faceList& faces = mesh.faces();
407 
408  for (label facei = 0; facei < mesh.nInternalFaces(); facei++)
409  {
410  const face& f = faces[facei];
411  const point& fc = mesh.faceCentres()[facei];
412 
413  forAll(f, fp)
414  {
415  weightedPosition& wp = globalSum[f[fp]];
416  wp.first() += 1.0;
417  wp.second() += fc;
418  }
419  }
420 
421  // Count coupled faces as internal ones (but only once)
422  const polyBoundaryMesh& patches = mesh.boundaryMesh();
423 
424  forAll(patches, patchi)
425  {
426  if
427  (
428  patches[patchi].coupled()
429  && refCast<const coupledPolyPatch>(patches[patchi]).owner()
430  )
431  {
432  const coupledPolyPatch& pp =
433  refCast<const coupledPolyPatch>(patches[patchi]);
434 
435  const vectorField::subField faceCentres = pp.faceCentres();
436 
437  forAll(pp, i)
438  {
439  const face& f = pp[i];
440  const point& fc = faceCentres[i];
441 
442  forAll(f, fp)
443  {
444  weightedPosition& wp = globalSum[f[fp]];
445  wp.first() += 1.0;
446  wp.second() += fc;
447  }
448  }
449  }
450  }
451 
452  // Add coupled contributions
454 
455  avgInternal.setSize(meshPoints.size());
456 
457  forAll(avgInternal, patchPointi)
458  {
459  label meshPointi = meshPoints[patchPointi];
460  const weightedPosition& wp = globalSum[meshPointi];
461 
462  avgInternal[patchPointi].first() = wp.first();
463  if (mag(wp.first()) < VSMALL)
464  {
465  // Set to zero?
466  avgInternal[patchPointi].second() = wp.second();
467  }
468  else
469  {
470  avgInternal[patchPointi].second() = wp.second()/wp.first();
471  }
472  }
473  }
474 
475 
476  // Precalculate any cell using mesh point (replacement of pointCells()[])
477  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
478 
479  labelList anyCell(mesh.nPoints(), -1);
480  forAll(mesh.faceOwner(), facei)
481  {
482  label own = mesh.faceOwner()[facei];
483  const face& f = mesh.faces()[facei];
484 
485  forAll(f, fp)
486  {
487  anyCell[f[fp]] = own;
488  }
489  }
490 
491 
492  // Displacement to calculate.
493  tmp<pointField> tpatchDisp(new pointField(meshPoints.size(), Zero));
494  pointField& patchDisp = tpatchDisp.ref();
495 
496  forAll(pointFaces, i)
497  {
498  label meshPointi = meshPoints[i];
499  const point& currentPos = pp.points()[meshPointi];
500 
501  // Now we have the two average points and their counts:
502  // avgBoundary and avgInternal
503  // Do some blending between the two.
504  // Note: the following section has some reasoning behind it but the
505  // blending factors can be experimented with.
506 
507  const weightedPosition& internal = avgInternal[i];
508  const weightedPosition& boundary = avgBoundary[i];
509 
510  point newPos;
511 
512  if (!nonManifoldPoint.test(i))
513  {
514  // Points that are manifold. Weight the internal and boundary
515  // by their number of faces and blend with
516  scalar internalBlend = 0.1;
517  scalar blend = 0.1;
518 
519  point avgPos =
520  (
521  internalBlend*internal.first()*internal.second()
522  +(1-internalBlend)*boundary.first()*boundary.second()
523  )
524  / (
525  internalBlend*internal.first()
526  +(1-internalBlend)*boundary.first()
527  );
528 
529  newPos = (1-blend)*avgPos + blend*currentPos;
530  }
531  else if (internal.first() == 0)
532  {
533  // Non-manifold without internal faces. Use any connected cell
534  // as internal point instead. Use precalculated any cell to avoid
535  // e.g. pointCells()[meshPointi][0]
536 
537  const point& cc = mesh.cellCentres()[anyCell[meshPointi]];
538 
539  scalar cellCBlend = 0.8;
540  scalar blend = 0.1;
541 
542  point avgPos = (1-cellCBlend)*boundary.second() + cellCBlend*cc;
543 
544  newPos = (1-blend)*avgPos + blend*currentPos;
545  }
546  else
547  {
548  // Non-manifold point with internal faces connected to them
549  scalar internalBlend = 0.9;
550  scalar blend = 0.1;
551 
552  point avgPos =
553  internalBlend*internal.second()
554  + (1-internalBlend)*boundary.second();
555 
556  newPos = (1-blend)*avgPos + blend*currentPos;
557  }
558 
559  patchDisp[i] = newPos - currentPos;
560  }
561 
562  return tpatchDisp;
563 }
564 //XXXXXXX
565 //Foam::tmp<Foam::pointField> Foam::snappySnapDriver::avg
566 //(
567 // const indirectPrimitivePatch& pp,
568 // const pointField& localPoints
569 //)
570 //{
571 // const labelListList& pointEdges = pp.pointEdges();
572 // const edgeList& edges = pp.edges();
573 //
574 // tmp<pointField> tavg(new pointField(pointEdges.size(), Zero));
575 // pointField& avg = tavg();
576 //
577 // forAll(pointEdges, verti)
578 // {
579 // vector& avgPos = avg[verti];
580 //
581 // const labelList& pEdges = pointEdges[verti];
582 //
583 // forAll(pEdges, myEdgei)
584 // {
585 // const edge& e = edges[pEdges[myEdgei]];
586 //
587 // label otherVerti = e.otherVertex(verti);
588 //
589 // avgPos += localPoints[otherVerti];
590 // }
591 //
592 // avgPos /= pEdges.size();
593 // }
594 // return tavg;
595 //}
596 //Foam::tmp<Foam::pointField>
597 //Foam::snappySnapDriver::smoothLambdaMuPatchDisplacement
598 //(
599 // const motionSmoother& meshMover,
600 // const List<labelPair>& baffles
601 //)
602 //{
603 // const indirectPrimitivePatch& pp = meshMover.patch();
604 // pointField newLocalPoints(pp.localPoints());
605 //
606 // const label iters = 90;
607 // const scalar lambda = 0.33;
608 // const scalar mu = 0.34;
609 //
610 // for (label iter = 0; iter < iters; iter++)
611 // {
612 // // Lambda
613 // newLocalPoints =
614 // (1 - lambda)*newLocalPoints
615 // + lambda*avg(pp, newLocalPoints);
616 //
617 // // Mu
618 // newLocalPoints =
619 // (1 + mu)*newLocalPoints
620 // - mu*avg(pp, newLocalPoints);
621 // }
622 // return newLocalPoints-pp.localPoints();
623 //}
624 //XXXXXXX
625 
626 
627 Foam::tmp<Foam::scalarField> Foam::snappySnapDriver::edgePatchDist
628 (
629  const pointMesh& pMesh,
631 )
632 {
633  const polyMesh& mesh = pMesh();
634 
635  // Set initial changed points to all the patch points
636  List<pointEdgePoint> wallInfo(pp.nPoints());
637 
638  forAll(pp.localPoints(), ppi)
639  {
640  wallInfo[ppi] = pointEdgePoint(pp.localPoints()[ppi], 0.0);
641  }
642 
643  // Current info on points
644  List<pointEdgePoint> allPointInfo(mesh.nPoints());
645 
646  // Current info on edges
647  List<pointEdgePoint> allEdgeInfo(mesh.nEdges());
648 
649  PointEdgeWave<pointEdgePoint> wallCalc
650  (
651  mesh,
652  pp.meshPoints(),
653  wallInfo,
654 
655  allPointInfo,
656  allEdgeInfo,
657  mesh.globalData().nTotalPoints() // max iterations
658  );
659 
660  // Copy edge values into scalarField
661  tmp<scalarField> tedgeDist(new scalarField(mesh.nEdges()));
662  scalarField& edgeDist = tedgeDist.ref();
663 
664  forAll(allEdgeInfo, edgei)
665  {
666  edgeDist[edgei] = Foam::sqrt(allEdgeInfo[edgei].distSqr());
667  }
668 
669  return tedgeDist;
670 }
671 
672 
673 void Foam::snappySnapDriver::dumpMove
674 (
675  const fileName& fName,
676  const pointField& meshPts,
677  const pointField& surfPts
678 )
679 {
680  // Dump direction of growth into file
681  Info<< "Dumping move direction to " << fName << endl;
682 
683  OFstream nearestStream(fName);
684 
685  label verti = 0;
686 
687  forAll(meshPts, pti)
688  {
689  meshTools::writeOBJ(nearestStream, meshPts[pti]);
690  verti++;
691 
692  meshTools::writeOBJ(nearestStream, surfPts[pti]);
693  verti++;
694 
695  nearestStream<< "l " << verti-1 << ' ' << verti << nl;
696  }
697 }
698 
699 
700 // Check whether all displacement vectors point outwards of patch. Return true
701 // if so.
702 bool Foam::snappySnapDriver::outwardsDisplacement
703 (
704  const indirectPrimitivePatch& pp,
705  const vectorField& patchDisp
706 )
707 {
708  const vectorField& faceNormals = pp.faceNormals();
709  const labelListList& pointFaces = pp.pointFaces();
710 
711  forAll(pointFaces, pointi)
712  {
713  const labelList& pFaces = pointFaces[pointi];
714 
715  vector disp(patchDisp[pointi]);
716 
717  scalar magDisp = mag(disp);
718 
719  if (magDisp > SMALL)
720  {
721  disp /= magDisp;
722 
723  bool outwards = meshTools::visNormal(disp, faceNormals, pFaces);
724 
725  if (!outwards)
726  {
727  Warning<< "Displacement " << patchDisp[pointi]
728  << " at mesh point " << pp.meshPoints()[pointi]
729  << " coord " << pp.points()[pp.meshPoints()[pointi]]
730  << " points through the surrounding patch faces" << endl;
731  return false;
732  }
733  }
734  else
735  {
736  //? Displacement small but in wrong direction. Would probably be ok.
737  }
738  }
739  return true;
740 }
741 
742 
743 void Foam::snappySnapDriver::freezeExposedPoints
744 (
745  const meshRefinement& meshRefiner,
746  const word& fzName, // faceZone name
747  const word& pzName, // pointZone name
748  const indirectPrimitivePatch& outside,
749  vectorField& outsideDisp
750 )
751 {
752  const fvMesh& mesh = meshRefiner.mesh();
753  const pointZoneMesh& pointZones = mesh.pointZones();
754 
755  bitSet isFrozenPoint(mesh.nPoints());
756 
757  // Add frozen points
758  const label pointZonei = pointZones.findZoneID(pzName);
759  if (pointZonei != -1)
760  {
761  isFrozenPoint.set(pointZones[pointZonei]);
762  }
763 
764  // Add (inside) points of frozen faces
765  const faceZoneMesh& faceZones = mesh.faceZones();
766  const label faceZonei = faceZones.findZoneID(fzName);
767  if (faceZonei != -1)
768  {
770  (
771  UIndirectList<face>(mesh.faces(), faceZones[faceZonei]),
772  mesh.points()
773  );
774 
775  // Count number of faces per edge
776  const labelList nEdgeFaces(meshRefiner.countEdgeFaces(pp));
777 
778  // Freeze all internal points
779  forAll(nEdgeFaces, edgei)
780  {
781  if (nEdgeFaces[edgei] != 1)
782  {
783  const edge& e = pp.edges()[edgei];
784  isFrozenPoint.set(pp.meshPoints()[e[0]]);
785  isFrozenPoint.set(pp.meshPoints()[e[1]]);
786  }
787  }
788  }
789 
791  (
792  mesh,
793  isFrozenPoint,
794  orEqOp<unsigned int>(),
795  0u
796  );
797 
798  for (const label pointi : isFrozenPoint)
799  {
800  const auto iter = outside.meshPointMap().find(pointi);
801  if (iter.good())
802  {
803  outsideDisp[iter.val()] = Zero;
804  }
805  }
806 }
807 
808 
809 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
810 
811 Foam::snappySnapDriver::snappySnapDriver
812 (
813  meshRefinement& meshRefiner,
814  const labelList& globalToMasterPatch,
815  const labelList& globalToSlavePatch,
816  const bool dryRun
817 )
818 :
819  meshRefiner_(meshRefiner),
820  globalToMasterPatch_(globalToMasterPatch),
821  globalToSlavePatch_(globalToSlavePatch),
822  dryRun_(dryRun)
823 {}
824 
825 
826 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
827 
829 (
830  const fvMesh& mesh,
831  const snapParameters& snapParams,
833 )
834 {
835  const edgeList& edges = pp.edges();
836  const labelListList& pointEdges = pp.pointEdges();
837  const pointField& localPoints = pp.localPoints();
838 
839  scalarField maxEdgeLen(localPoints.size(), -GREAT);
840 
841  forAll(pointEdges, pointi)
842  {
843  const labelList& pEdges = pointEdges[pointi];
844 
845  forAll(pEdges, pEdgei)
846  {
847  const edge& e = edges[pEdges[pEdgei]];
848 
849  scalar len = e.mag(localPoints);
850 
851  maxEdgeLen[pointi] = max(maxEdgeLen[pointi], len);
852  }
853  }
854 
856  (
857  mesh,
858  pp.meshPoints(),
859  maxEdgeLen,
860  maxEqOp<scalar>(), // combine op
861  -GREAT // null value
862  );
864  return scalarField(snapParams.snapTol()*maxEdgeLen);
865 }
866 
867 
869 (
870  const meshRefinement& meshRefiner,
871  const snapParameters& snapParams,
872  const label nInitErrors,
873  const List<labelPair>& baffles,
874  motionSmoother& meshMover
875 )
876 {
877  addProfiling(smooth, "snappyHexMesh::snap::smoothing");
878  const fvMesh& mesh = meshRefiner.mesh();
879 
880  labelList checkFaces;
881 
882  if (snapParams.nSmoothInternal() > 0)
883  {
884  Info<< "Smoothing patch and internal points ..." << endl;
885  }
886  else
887  {
888  Info<< "Smoothing patch points ..." << endl;
889  }
890 
891  vectorField& pointDisp = meshMover.pointDisplacement().primitiveFieldRef();
892 
893  for
894  (
895  label smoothIter = 0;
896  smoothIter < snapParams.nSmoothPatch();
897  smoothIter++
898  )
899  {
900  Info<< "Smoothing iteration " << smoothIter << endl;
901  checkFaces.setSize(mesh.nFaces());
902  forAll(checkFaces, facei)
903  {
904  checkFaces[facei] = facei;
905  }
906 
907  // If enabled smooth the internal points
908  if (snapParams.nSmoothInternal() > smoothIter)
909  {
910  // Override values on internal points on refinement interfaces
911  pointDisp = smoothInternalDisplacement(meshRefiner, meshMover);
912  }
913 
914  // Smooth the patch points
915  pointField patchDisp(smoothPatchDisplacement(meshMover, baffles));
916  //pointField patchDisp
917  //(
918  // smoothLambdaMuPatchDisplacement(meshMover, baffles)
919  //);
920 
921  // Take over patch displacement as boundary condition on
922  // pointDisplacement
923  meshMover.setDisplacement(patchDisp);
924 
925  // Start off from current mesh.points()
926  meshMover.correct();
927 
928  scalar oldErrorReduction = -1;
929 
930  for (label snapIter = 0; snapIter < 2*snapParams.nSnap(); snapIter++)
931  {
932  Info<< nl << "Scaling iteration " << snapIter << endl;
933 
934  if (snapIter == snapParams.nSnap())
935  {
936  Info<< "Displacement scaling for error reduction set to 0."
937  << endl;
938  oldErrorReduction = meshMover.setErrorReduction(0.0);
939  }
940 
941  // Try to adapt mesh to obtain displacement by smoothly
942  // decreasing displacement at error locations.
943  if (meshMover.scaleMesh(checkFaces, baffles, true, nInitErrors))
944  {
945  Info<< "Successfully moved mesh" << endl;
946  break;
947  }
948  }
949 
950  if (oldErrorReduction >= 0)
951  {
952  meshMover.setErrorReduction(oldErrorReduction);
953  }
954  Info<< endl;
955  }
956 
957 
958  // The current mesh is the starting mesh to smooth from.
959  meshMover.correct();
960 
962  {
963  const_cast<Time&>(mesh.time())++;
964  Info<< "Writing patch smoothed mesh to time "
965  << meshRefiner.timeName() << '.' << endl;
966  meshRefiner.write
967  (
970  (
973  ),
974  mesh.time().path()/meshRefiner.timeName()
975  );
976  Info<< "Dumped mesh in = "
977  << mesh.time().cpuTimeIncrement() << " s\n" << nl << endl;
978  }
979 
980  Info<< "Patch points smoothed in = "
981  << mesh.time().cpuTimeIncrement() << " s\n" << nl << endl;
982 }
983 
984 
985 // Get (pp-local) indices of points that are both on zone and on patched surface
986 void Foam::snappySnapDriver::getZoneSurfacePoints
987 (
988  const fvMesh& mesh,
989  const indirectPrimitivePatch& pp,
990  const word& zoneName,
991 
992  bitSet& pointOnZone
993 )
994 {
995  label zonei = mesh.faceZones().findZoneID(zoneName);
996 
997  if (zonei == -1)
998  {
1000  << "Cannot find zone " << zoneName
1001  << exit(FatalError);
1002  }
1003 
1004  const faceZone& fZone = mesh.faceZones()[zonei];
1005 
1006 
1007  // Could use PrimitivePatch & localFaces to extract points but might just
1008  // as well do it ourselves.
1009 
1010  forAll(fZone, i)
1011  {
1012  const face& f = mesh.faces()[fZone[i]];
1013 
1014  forAll(f, fp)
1015  {
1016  label meshPointi = f[fp];
1017 
1018  const auto iter = pp.meshPointMap().cfind(meshPointi);
1019 
1020  if (iter.good())
1021  {
1022  const label pointi = iter.val();
1023  pointOnZone[pointi] = true;
1024  }
1025  }
1026  }
1027 }
1028 
1029 
1031 (
1032  const fvMesh& mesh,
1033  const indirectPrimitivePatch& pp
1034 )
1035 {
1036  const labelListList& pointFaces = pp.pointFaces();
1037 
1038  Field<weightedPosition> avgBoundary
1039  (
1040  pointFaces.size(),
1042  );
1043 
1044  forAll(pointFaces, pointi)
1045  {
1046  const labelList& pFaces = pointFaces[pointi];
1047 
1048  avgBoundary[pointi].first() = pFaces.size();
1049  forAll(pFaces, pfi)
1050  {
1051  label facei = pFaces[pfi];
1052  label own = mesh.faceOwner()[pp.addressing()[facei]];
1053  avgBoundary[pointi].second() += mesh.cellCentres()[own];
1054  }
1055  }
1056 
1057  // Add coupled contributions
1059 
1060  tmp<pointField> tavgBoundary(new pointField(avgBoundary.size()));
1061  weightedPosition::getPoints(avgBoundary, tavgBoundary.ref());
1062 
1063  return tavgBoundary;
1064 }
1065 
1066 
1067 //Foam::tmp<Foam::scalarField> Foam::snappySnapDriver::calcEdgeLen
1068 //(
1069 // const indirectPrimitivePatch& pp
1070 //) const
1071 //{
1072 // // Get local edge length based on refinement level
1073 // // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1074 // // (Ripped from snappyLayerDriver)
1075 //
1076 // tmp<scalarField> tedgeLen(new scalarField(pp.nPoints()));
1077 // scalarField& edgeLen = tedgeLen();
1078 // {
1079 // const fvMesh& mesh = meshRefiner_.mesh();
1080 // const scalar edge0Len = meshRefiner_.meshCutter().level0EdgeLength();
1081 // const labelList& cellLevel = meshRefiner_.meshCutter().cellLevel();
1082 //
1083 // labelList maxPointLevel(pp.nPoints(), labelMin);
1084 //
1085 // forAll(pp, i)
1086 // {
1087 // label ownLevel = cellLevel[mesh.faceOwner()[pp.addressing()[i]]];
1088 // const face& f = pp.localFaces()[i];
1089 // forAll(f, fp)
1090 // {
1091 // maxPointLevel[f[fp]] = max(maxPointLevel[f[fp]], ownLevel);
1092 // }
1093 // }
1094 //
1095 // syncTools::syncPointList
1096 // (
1097 // mesh,
1098 // pp.meshPoints(),
1099 // maxPointLevel,
1100 // maxEqOp<label>(),
1101 // labelMin // null value
1102 // );
1103 //
1104 //
1105 // forAll(maxPointLevel, pointi)
1106 // {
1107 // // Find undistorted edge size for this level.
1108 // edgeLen[pointi] = edge0Len/(1<<maxPointLevel[pointi]);
1109 // }
1110 // }
1111 // return tedgeLen;
1112 //}
1113 
1114 
1116 (
1117  const scalar planarCos,
1118  const indirectPrimitivePatch& pp,
1119  const pointField& nearestPoint,
1120  const vectorField& nearestNormal,
1121 
1122  vectorField& disp
1123 ) const
1124 {
1125  Info<< "Detecting near surfaces ..." << endl;
1126 
1127  const pointField& localPoints = pp.localPoints();
1128  const labelList& meshPoints = pp.meshPoints();
1129  const refinementSurfaces& surfaces = meshRefiner_.surfaces();
1130  const fvMesh& mesh = meshRefiner_.mesh();
1131 
1133  //const scalarField edgeLen(calcEdgeLen(pp));
1134  //
1137  //
1138  //{
1139  // const vector n = normalised(vector::one);
1140  //
1141  // pointField start(14*pp.nPoints());
1142  // pointField end(start.size());
1143  //
1144  // label rayi = 0;
1145  // forAll(localPoints, pointi)
1146  // {
1147  // const point& pt = localPoints[pointi];
1148  //
1149  // // Along coordinate axes
1150  //
1151  // {
1152  // start[rayi] = pt;
1153  // point& endPt = end[rayi++];
1154  // endPt = pt;
1155  // endPt.x() -= edgeLen[pointi];
1156  // }
1157  // {
1158  // start[rayi] = pt;
1159  // point& endPt = end[rayi++];
1160  // endPt = pt;
1161  // endPt.x() += edgeLen[pointi];
1162  // }
1163  // {
1164  // start[rayi] = pt;
1165  // point& endPt = end[rayi++];
1166  // endPt = pt;
1167  // endPt.y() -= edgeLen[pointi];
1168  // }
1169  // {
1170  // start[rayi] = pt;
1171  // point& endPt = end[rayi++];
1172  // endPt = pt;
1173  // endPt.y() += edgeLen[pointi];
1174  // }
1175  // {
1176  // start[rayi] = pt;
1177  // point& endPt = end[rayi++];
1178  // endPt = pt;
1179  // endPt.z() -= edgeLen[pointi];
1180  // }
1181  // {
1182  // start[rayi] = pt;
1183  // point& endPt = end[rayi++];
1184  // endPt = pt;
1185  // endPt.z() += edgeLen[pointi];
1186  // }
1187  //
1188  // // At 45 degrees
1189  //
1190  // const vector vec(edgeLen[pointi]*n);
1191  //
1192  // {
1193  // start[rayi] = pt;
1194  // point& endPt = end[rayi++];
1195  // endPt = pt;
1196  // endPt.x() += vec.x();
1197  // endPt.y() += vec.y();
1198  // endPt.z() += vec.z();
1199  // }
1200  // {
1201  // start[rayi] = pt;
1202  // point& endPt = end[rayi++];
1203  // endPt = pt;
1204  // endPt.x() -= vec.x();
1205  // endPt.y() += vec.y();
1206  // endPt.z() += vec.z();
1207  // }
1208  // {
1209  // start[rayi] = pt;
1210  // point& endPt = end[rayi++];
1211  // endPt = pt;
1212  // endPt.x() += vec.x();
1213  // endPt.y() -= vec.y();
1214  // endPt.z() += vec.z();
1215  // }
1216  // {
1217  // start[rayi] = pt;
1218  // point& endPt = end[rayi++];
1219  // endPt = pt;
1220  // endPt.x() -= vec.x();
1221  // endPt.y() -= vec.y();
1222  // endPt.z() += vec.z();
1223  // }
1224  // {
1225  // start[rayi] = pt;
1226  // point& endPt = end[rayi++];
1227  // endPt = pt;
1228  // endPt.x() += vec.x();
1229  // endPt.y() += vec.y();
1230  // endPt.z() -= vec.z();
1231  // }
1232  // {
1233  // start[rayi] = pt;
1234  // point& endPt = end[rayi++];
1235  // endPt = pt;
1236  // endPt.x() -= vec.x();
1237  // endPt.y() += vec.y();
1238  // endPt.z() -= vec.z();
1239  // }
1240  // {
1241  // start[rayi] = pt;
1242  // point& endPt = end[rayi++];
1243  // endPt = pt;
1244  // endPt.x() += vec.x();
1245  // endPt.y() -= vec.y();
1246  // endPt.z() -= vec.z();
1247  // }
1248  // {
1249  // start[rayi] = pt;
1250  // point& endPt = end[rayi++];
1251  // endPt = pt;
1252  // endPt.x() -= vec.x();
1253  // endPt.y() -= vec.y();
1254  // endPt.z() -= vec.z();
1255  // }
1256  // }
1257  //
1258  // labelList surface1;
1259  // List<pointIndexHit> hit1;
1260  // labelList region1;
1261  // vectorField normal1;
1262  //
1263  // labelList surface2;
1264  // List<pointIndexHit> hit2;
1265  // labelList region2;
1266  // vectorField normal2;
1267  // surfaces.findNearestIntersection
1268  // (
1269  // unzonedSurfaces, // surfacesToTest,
1270  // start,
1271  // end,
1272  //
1273  // surface1,
1274  // hit1,
1275  // region1,
1276  // normal1,
1277  //
1278  // surface2,
1279  // hit2,
1280  // region2,
1281  // normal2
1282  // );
1283  //
1284  // // All intersections
1285  // {
1286  // OBJstream str
1287  // (
1288  // mesh.time().path()
1289  // / "surfaceHits_" + meshRefiner_.timeName() + ".obj"
1290  // );
1291  //
1292  // Info<< "Dumping intersections with rays to " << str.name()
1293  // << endl;
1294  //
1295  // forAll(hit1, i)
1296  // {
1297  // if (hit1[i].hit())
1298  // {
1299  // str.writeLine(start[i], hit1[i].point());
1300  // }
1301  // if (hit2[i].hit())
1302  // {
1303  // str.writeLine(start[i], hit2[i].point());
1304  // }
1305  // }
1306  // }
1307  //
1308  // // Co-planar intersections
1309  // {
1310  // OBJstream str
1311  // (
1312  // mesh.time().path()
1313  // / "coplanarHits_" + meshRefiner_.timeName() + ".obj"
1314  // );
1315  //
1316  // Info<< "Dumping intersections with co-planar surfaces to "
1317  // << str.name() << endl;
1318  //
1319  // forAll(localPoints, pointi)
1320  // {
1321  // bool hasNormal = false;
1322  // point surfPointA;
1323  // vector surfNormalA;
1324  // point surfPointB;
1325  // vector surfNormalB;
1326  //
1327  // bool isCoplanar = false;
1328  //
1329  // label rayi = 14*pointi;
1330  // for (label i = 0; i < 14; i++)
1331  // {
1332  // if (hit1[rayi].hit())
1333  // {
1334  // const point& pt = hit1[rayi].point();
1335  // const vector& n = normal1[rayi];
1336  //
1337  // if (!hasNormal)
1338  // {
1339  // hasNormal = true;
1340  // surfPointA = pt;
1341  // surfNormalA = n;
1342  // }
1343  // else
1344  // {
1345  // if
1346  // (
1347  // meshRefiner_.isGap
1348  // (
1349  // planarCos,
1350  // surfPointA,
1351  // surfNormalA,
1352  // pt,
1353  // n
1354  // )
1355  // )
1356  // {
1357  // isCoplanar = true;
1358  // surfPointB = pt;
1359  // surfNormalB = n;
1360  // break;
1361  // }
1362  // }
1363  // }
1364  // if (hit2[rayi].hit())
1365  // {
1366  // const point& pt = hit2[rayi].point();
1367  // const vector& n = normal2[rayi];
1368  //
1369  // if (!hasNormal)
1370  // {
1371  // hasNormal = true;
1372  // surfPointA = pt;
1373  // surfNormalA = n;
1374  // }
1375  // else
1376  // {
1377  // if
1378  // (
1379  // meshRefiner_.isGap
1380  // (
1381  // planarCos,
1382  // surfPointA,
1383  // surfNormalA,
1384  // pt,
1385  // n
1386  // )
1387  // )
1388  // {
1389  // isCoplanar = true;
1390  // surfPointB = pt;
1391  // surfNormalB = n;
1392  // break;
1393  // }
1394  // }
1395  // }
1396  //
1397  // rayi++;
1398  // }
1399  //
1400  // if (isCoplanar)
1401  // {
1402  // str.writeLine(surfPointA, surfPointB);
1403  // }
1404  // }
1405  // }
1406  //}
1407 
1408 
1409  const pointField avgCc(avgCellCentres(mesh, pp));
1410 
1411  // Construct rays through localPoints to beyond cell centre
1412  pointField start(pp.nPoints());
1413  pointField end(pp.nPoints());
1414  forAll(localPoints, pointi)
1415  {
1416  const point& pt = localPoints[pointi];
1417  const vector d = 2*(avgCc[pointi]-pt);
1418  start[pointi] = pt - d;
1419  end[pointi] = pt + d;
1420  }
1421 
1422 
1423  autoPtr<OBJstream> gapStr;
1425  {
1426  gapStr.reset
1427  (
1428  new OBJstream
1429  (
1430  mesh.time().path()
1431  / "detectNearSurfaces_" + meshRefiner_.timeName() + ".obj"
1432  )
1433  );
1434  }
1435 
1436 
1437  const bitSet isPatchMasterPoint
1438  (
1440  (
1441  mesh,
1442  meshPoints
1443  )
1444  );
1445 
1446  label nOverride = 0;
1447 
1448  // 1. All points to non-interface surfaces
1449  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1450  {
1451  const labelList unzonedSurfaces =
1453  (
1454  meshRefiner_.surfaces().surfZones()
1455  );
1456 
1457  // Do intersection test
1458  labelList surface1;
1459  List<pointIndexHit> hit1;
1460  labelList region1;
1461  vectorField normal1;
1462 
1463  labelList surface2;
1464  List<pointIndexHit> hit2;
1465  labelList region2;
1466  vectorField normal2;
1467  surfaces.findNearestIntersection
1468  (
1469  unzonedSurfaces,
1470  start,
1471  end,
1472 
1473  surface1,
1474  hit1,
1475  region1,
1476  normal1,
1477 
1478  surface2,
1479  hit2,
1480  region2,
1481  normal2
1482  );
1483 
1484 
1485  forAll(localPoints, pointi)
1486  {
1487  // Current location
1488  const point& pt = localPoints[pointi];
1489 
1490  bool override = false;
1491 
1492  //if (hit1[pointi].hit())
1493  //{
1494  // if
1495  // (
1496  // meshRefiner_.isGap
1497  // (
1498  // planarCos,
1499  // nearestPoint[pointi],
1500  // nearestNormal[pointi],
1501  // hit1[pointi].point(),
1502  // normal1[pointi]
1503  // )
1504  // )
1505  // {
1506  // disp[pointi] = hit1[pointi].point()-pt;
1507  // override = true;
1508  // }
1509  //}
1510  //if (hit2[pointi].hit())
1511  //{
1512  // if
1513  // (
1514  // meshRefiner_.isGap
1515  // (
1516  // planarCos,
1517  // nearestPoint[pointi],
1518  // nearestNormal[pointi],
1519  // hit2[pointi].point(),
1520  // normal2[pointi]
1521  // )
1522  // )
1523  // {
1524  // disp[pointi] = hit2[pointi].point()-pt;
1525  // override = true;
1526  // }
1527  //}
1528 
1529  if (hit1[pointi].hit() && hit2[pointi].hit())
1530  {
1531  if
1532  (
1533  meshRefiner_.isGap
1534  (
1535  planarCos,
1536  hit1[pointi].point(),
1537  normal1[pointi],
1538  hit2[pointi].point(),
1539  normal2[pointi]
1540  )
1541  )
1542  {
1543  // TBD: check if the attraction (to nearest) would attract
1544  // good enough and not override attraction
1545 
1546  if (gapStr)
1547  {
1548  gapStr().writeLine(pt, hit2[pointi].point());
1549  }
1550 
1551  // Choose hit2 : nearest to end point (so inside the domain)
1552  disp[pointi] = hit2[pointi].point()-pt;
1553  override = true;
1554  }
1555  }
1556 
1557  if (override && isPatchMasterPoint[pointi])
1558  {
1559  nOverride++;
1560  }
1561  }
1562  }
1563 
1564 
1565  // 2. All points on zones to their respective surface
1566  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1567 
1568  {
1569  // Surfaces with zone information
1570  const PtrList<surfaceZonesInfo>& surfZones = surfaces.surfZones();
1571 
1572  const labelList zonedSurfaces = surfaceZonesInfo::getNamedSurfaces
1573  (
1574  surfZones
1575  );
1576 
1577  forAll(zonedSurfaces, i)
1578  {
1579  label zoneSurfi = zonedSurfaces[i];
1580  const labelList surfacesToTest(1, zoneSurfi);
1581 
1582  const wordList& faceZoneNames =
1583  surfZones[zoneSurfi].faceZoneNames();
1584  forAll(faceZoneNames, namei)
1585  {
1586  const word& faceZoneName = faceZoneNames[namei];
1587 
1588  // Get indices of points both on faceZone and on pp.
1589  bitSet pointOnZone(pp.nPoints());
1590  getZoneSurfacePoints
1591  (
1592  mesh,
1593  pp,
1594  faceZoneName,
1595  pointOnZone
1596  );
1597  const labelList zonePointIndices(pointOnZone.toc());
1598 
1599  // Do intersection test
1600  labelList surface1;
1601  List<pointIndexHit> hit1;
1602  labelList region1;
1603  vectorField normal1;
1604 
1605  labelList surface2;
1606  List<pointIndexHit> hit2;
1607  labelList region2;
1608  vectorField normal2;
1609  surfaces.findNearestIntersection
1610  (
1611  surfacesToTest,
1612  pointField(start, zonePointIndices),
1613  pointField(end, zonePointIndices),
1614 
1615  surface1,
1616  hit1,
1617  region1,
1618  normal1,
1619 
1620  surface2,
1621  hit2,
1622  region2,
1623  normal2
1624  );
1625 
1626 
1627  forAll(hit1, i)
1628  {
1629  label pointi = zonePointIndices[i];
1630 
1631  // Current location
1632  const point& pt = localPoints[pointi];
1633 
1634  bool override = false;
1635 
1636  //if (hit1[i].hit())
1637  //{
1638  // if
1639  // (
1640  // meshRefiner_.isGap
1641  // (
1642  // planarCos,
1643  // nearestPoint[pointi],
1644  // nearestNormal[pointi],
1645  // hit1[i].point(),
1646  // normal1[i]
1647  // )
1648  // )
1649  // {
1650  // disp[pointi] = hit1[i].point()-pt;
1651  // override = true;
1652  // }
1653  //}
1654  //if (hit2[i].hit())
1655  //{
1656  // if
1657  // (
1658  // meshRefiner_.isGap
1659  // (
1660  // planarCos,
1661  // nearestPoint[pointi],
1662  // nearestNormal[pointi],
1663  // hit2[i].point(),
1664  // normal2[i]
1665  // )
1666  // )
1667  // {
1668  // disp[pointi] = hit2[i].point()-pt;
1669  // override = true;
1670  // }
1671  //}
1672 
1673  if (hit1[i].hit() && hit2[i].hit())
1674  {
1675  if
1676  (
1677  meshRefiner_.isGap
1678  (
1679  planarCos,
1680  hit1[i].point(),
1681  normal1[i],
1682  hit2[i].point(),
1683  normal2[i]
1684  )
1685  )
1686  {
1687  if (gapStr)
1688  {
1689  gapStr().writeLine(pt, hit2[i].point());
1690  }
1691 
1692  disp[pointi] = hit2[i].point()-pt;
1693  override = true;
1694  }
1695  }
1696 
1697  if (override && isPatchMasterPoint[pointi])
1698  {
1699  nOverride++;
1700  }
1701  }
1702  }
1703  }
1704  }
1705 
1706  Info<< "Overriding nearest with intersection of close gaps at "
1707  << returnReduce(nOverride, sumOp<label>())
1708  << " out of " << returnReduce(pp.nPoints(), sumOp<label>())
1709  << " points." << endl;
1710 }
1711 
1712 
1713 void Foam::snappySnapDriver::calcNearestSurface
1714 (
1715  const refinementSurfaces& surfaces,
1716 
1717  const labelList& surfacesToTest,
1718  const labelListList& regionsToTest,
1719 
1720  const pointField& localPoints,
1721  const labelList& zonePointIndices,
1722 
1723  scalarField& minSnapDist,
1724  labelList& snapSurf,
1725  vectorField& patchDisp,
1726 
1727  // Optional: nearest point, normal
1728  pointField& nearestPoint,
1729  vectorField& nearestNormal
1730 )
1731 {
1732  // Find nearest for points both on faceZone and pp.
1733  List<pointIndexHit> hitInfo;
1734  labelList hitSurface;
1735 
1736  if (nearestNormal.size() == localPoints.size())
1737  {
1738  labelList hitRegion;
1739  vectorField hitNormal;
1740  surfaces.findNearestRegion
1741  (
1742  surfacesToTest,
1743  regionsToTest,
1744 
1745  pointField(localPoints, zonePointIndices),
1746  sqr(scalarField(minSnapDist, zonePointIndices)),
1747 
1748  hitSurface,
1749  hitInfo,
1750  hitRegion,
1751  hitNormal
1752  );
1753 
1754  forAll(hitInfo, i)
1755  {
1756  if (hitInfo[i].hit())
1757  {
1758  label pointi = zonePointIndices[i];
1759  nearestPoint[pointi] = hitInfo[i].point();
1760  nearestNormal[pointi] = hitNormal[i];
1761  }
1762  }
1763  }
1764  else
1765  {
1766  surfaces.findNearest
1767  (
1768  surfacesToTest,
1769  regionsToTest,
1770 
1771  pointField(localPoints, zonePointIndices),
1772  sqr(scalarField(minSnapDist, zonePointIndices)),
1773 
1774  hitSurface,
1775  hitInfo
1776  );
1777  }
1778 
1779  forAll(hitInfo, i)
1780  {
1781  if (hitInfo[i].hit())
1782  {
1783  label pointi = zonePointIndices[i];
1784 
1785  patchDisp[pointi] = hitInfo[i].point() - localPoints[pointi];
1786  minSnapDist[pointi] = mag(patchDisp[pointi]);
1787  snapSurf[pointi] = hitSurface[i];
1788  }
1789  }
1790 }
1791 
1792 
1793 Foam::vectorField Foam::snappySnapDriver::calcNearestSurface
1794 (
1795  const bool strictRegionSnap,
1796  const meshRefinement& meshRefiner,
1797  const labelList& globalToMasterPatch,
1798  const labelList& globalToSlavePatch,
1799  const scalarField& snapDist,
1800  const indirectPrimitivePatch& pp,
1801  pointField& nearestPoint,
1802  vectorField& nearestNormal
1803 )
1804 {
1805  Info<< "Calculating patchDisplacement as distance to nearest surface"
1806  << " point ..." << endl;
1807  if (strictRegionSnap)
1808  {
1809  Info<< " non-zone points : attract to local region on surface only"
1810  << nl
1811  << " zone points : attract to local region on surface only"
1812  << nl
1813  << endl;
1814  }
1815  else
1816  {
1817  Info<< " non-zone points :"
1818  << " attract to nearest of all non-zone surfaces"
1819  << nl
1820  << " zone points : attract to zone surface only" << nl
1821  << endl;
1822  }
1823 
1824 
1825  const pointField& localPoints = pp.localPoints();
1826  const refinementSurfaces& surfaces = meshRefiner.surfaces();
1827  const fvMesh& mesh = meshRefiner.mesh();
1828 
1829  // Displacement per patch point
1830  vectorField patchDisp(localPoints.size(), Zero);
1831 
1832  if (returnReduceOr(localPoints.size()))
1833  {
1834  // Current surface snapped to. Used to check whether points have been
1835  // snapped at all
1836  labelList snapSurf(localPoints.size(), -1);
1837 
1838  // Current best snap distance (since point might be on multiple
1839  // regions)
1840  scalarField minSnapDist(snapDist);
1841 
1842 
1843  if (strictRegionSnap)
1844  {
1845  // Attract patch points to same region only
1846 
1847  forAll(surfaces.surfaces(), surfi)
1848  {
1849  label geomi = surfaces.surfaces()[surfi];
1850  label nRegions = surfaces.geometry()[geomi].regions().size();
1851 
1852  const labelList surfacesToTest(1, surfi);
1853 
1854  for (label regioni = 0; regioni < nRegions; regioni++)
1855  {
1856  label globali = surfaces.globalRegion(surfi, regioni);
1857  label masterPatchi = globalToMasterPatch[globali];
1858 
1859  // Get indices of points both on patch and on pp
1860  labelList zonePointIndices
1861  (
1862  getFacePoints
1863  (
1864  pp,
1865  mesh.boundaryMesh()[masterPatchi]
1866  )
1867  );
1868 
1869  calcNearestSurface
1870  (
1871  surfaces,
1872 
1873  surfacesToTest,
1874  labelListList(1, labelList(1, regioni)), //regionsToTest
1875 
1876  localPoints,
1877  zonePointIndices,
1878 
1879  minSnapDist,
1880  snapSurf,
1881  patchDisp,
1882 
1883  // Optional: nearest point, normal
1884  nearestPoint,
1885  nearestNormal
1886  );
1887 
1888  if (globalToSlavePatch[globali] != masterPatchi)
1889  {
1890  label slavePatchi = globalToSlavePatch[globali];
1891 
1892  // Get indices of points both on patch and on pp
1893  labelList zonePointIndices
1894  (
1895  getFacePoints
1896  (
1897  pp,
1898  mesh.boundaryMesh()[slavePatchi]
1899  )
1900  );
1901 
1902  calcNearestSurface
1903  (
1904  surfaces,
1905 
1906  surfacesToTest,
1907  labelListList(1, labelList(1, regioni)),
1908 
1909  localPoints,
1910  zonePointIndices,
1911 
1912  minSnapDist,
1913  snapSurf,
1914  patchDisp,
1915 
1916  // Optional: nearest point, normal
1917  nearestPoint,
1918  nearestNormal
1919  );
1920  }
1921  }
1922  }
1923  }
1924  else
1925  {
1926  // Divide surfaces into zoned and unzoned
1927  const labelList unzonedSurfaces =
1929  (
1930  meshRefiner.surfaces().surfZones()
1931  );
1932 
1933 
1934  // 1. All points to non-interface surfaces
1935  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1936 
1937  List<pointIndexHit> hitInfo;
1938  labelList hitSurface;
1939 
1940  if (nearestNormal.size() == localPoints.size())
1941  {
1942  labelList hitRegion;
1943  vectorField hitNormal;
1944  surfaces.findNearestRegion
1945  (
1946  unzonedSurfaces,
1947  localPoints,
1948  sqr(snapDist),
1949  hitSurface,
1950  hitInfo,
1951  hitRegion,
1952  hitNormal
1953  );
1954 
1955  forAll(hitInfo, pointi)
1956  {
1957  if (hitInfo[pointi].hit())
1958  {
1959  nearestPoint[pointi] = hitInfo[pointi].point();
1960  nearestNormal[pointi] = hitNormal[pointi];
1961  }
1962  }
1963  }
1964  else
1965  {
1966  surfaces.findNearest
1967  (
1968  unzonedSurfaces,
1969  localPoints,
1970  sqr(snapDist), // sqr of attract distance
1971  hitSurface,
1972  hitInfo
1973  );
1974  }
1975 
1976  forAll(hitInfo, pointi)
1977  {
1978  if (hitInfo[pointi].hit())
1979  {
1980  patchDisp[pointi] =
1981  hitInfo[pointi].point()
1982  - localPoints[pointi];
1983 
1984  snapSurf[pointi] = hitSurface[pointi];
1985  }
1986  }
1987 
1988 
1989  const labelList zonedSurfaces = surfaceZonesInfo::getNamedSurfaces
1990  (
1991  meshRefiner.surfaces().surfZones()
1992  );
1993 
1994 
1995  // 2. All points on zones to their respective surface
1996  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1997  // (ignoring faceZone subdivision)
1998 
1999  // Surfaces with zone information
2000  const PtrList<surfaceZonesInfo>& surfZones = surfaces.surfZones();
2001 
2002  forAll(zonedSurfaces, i)
2003  {
2004  label surfi = zonedSurfaces[i];
2005  const labelList surfacesToTest(1, surfi);
2006  const label geomi = surfaces.surfaces()[surfi];
2007  const label nRegions =
2008  surfaces.geometry()[geomi].regions().size();
2009 
2010  const wordList& faceZoneNames =
2011  surfZones[surfi].faceZoneNames();
2012 
2013  // Get indices of points both on any faceZone and on pp.
2014  bitSet pointOnZone(pp.nPoints());
2015  forAll(faceZoneNames, locali)
2016  {
2017  getZoneSurfacePoints
2018  (
2019  mesh,
2020  pp,
2021  faceZoneNames[locali],
2022  pointOnZone
2023  );
2024  }
2025  const labelList zonePointIndices(pointOnZone.toc());
2026 
2027  calcNearestSurface
2028  (
2029  surfaces,
2030 
2031  surfacesToTest,
2032  labelListList(1, identity(nRegions)),
2033 
2034  localPoints,
2035  zonePointIndices,
2036 
2037  minSnapDist,
2038  snapSurf,
2039  patchDisp,
2040 
2041  // Optional: nearest point, normal
2042  nearestPoint,
2043  nearestNormal
2044  );
2045  }
2046  }
2047 
2048 
2049  // Check if all points are being snapped
2050  forAll(snapSurf, pointi)
2051  {
2052  if (snapSurf[pointi] == -1)
2053  {
2054  static label nWarn = 0;
2055 
2056  if (nWarn < 100)
2057  {
2059  << "For point:" << pointi
2060  << " coordinate:" << localPoints[pointi]
2061  << " did not find any surface within:"
2062  << minSnapDist[pointi] << " metre." << endl;
2063  nWarn++;
2064  if (nWarn == 100)
2065  {
2067  << "Reached warning limit " << nWarn
2068  << ". Suppressing further warnings." << endl;
2069  }
2070  }
2071  }
2072  }
2073 
2074  {
2075  const bitSet isPatchMasterPoint
2076  (
2078  (
2079  mesh,
2080  pp.meshPoints()
2081  )
2082  );
2083 
2084  scalarField magDisp(mag(patchDisp));
2085 
2086  Info<< "Wanted displacement : average:"
2087  << meshRefinement::gAverage(isPatchMasterPoint, magDisp)
2088  << " min:" << gMin(magDisp)
2089  << " max:" << gMax(magDisp) << endl;
2090  }
2091  }
2092 
2093  Info<< "Calculated surface displacement in = "
2094  << mesh.time().cpuTimeIncrement() << " s\n" << nl << endl;
2095 
2096 
2097  // Limit amount of movement. Can not happen for triSurfaceMesh but
2098  // can happen for some analytical shapes?
2099  forAll(patchDisp, patchPointi)
2100  {
2101  scalar magDisp = mag(patchDisp[patchPointi]);
2102 
2103  if (magDisp > snapDist[patchPointi])
2104  {
2105  patchDisp[patchPointi] *= snapDist[patchPointi] / magDisp;
2106 
2107  Pout<< "Limiting displacement for " << patchPointi
2108  << " from " << magDisp << " to " << snapDist[patchPointi]
2109  << endl;
2110  }
2111  }
2112 
2113  // Points on zones in one domain but only present as point on other
2114  // will not do condition 2 on all. Sync explicitly.
2116  (
2117  mesh,
2118  pp.meshPoints(),
2119  patchDisp,
2120  minMagSqrEqOp<point>(), // combine op
2121  vector(GREAT, GREAT, GREAT) // null value (note: cannot use VGREAT)
2122  );
2124  return patchDisp;
2125 }
2126 
2127 
2129 (
2130  const snapParameters& snapParams,
2131  motionSmoother& meshMover
2132 ) const
2133 {
2134  if (dryRun_)
2135  {
2136  return;
2137  }
2138 
2139  const fvMesh& mesh = meshRefiner_.mesh();
2140  const indirectPrimitivePatch& pp = meshMover.patch();
2141 
2142  Info<< "Smoothing displacement ..." << endl;
2143 
2144  // Set edge diffusivity as inverse of distance to patch
2145  scalarField edgeGamma(1.0/(edgePatchDist(meshMover.pMesh(), pp) + SMALL));
2146  //scalarField edgeGamma(mesh.nEdges(), 1.0);
2147  //scalarField edgeGamma(wallGamma(mesh, pp, 10, 1));
2148 
2149  // Get displacement field
2150  pointVectorField& disp = meshMover.displacement();
2151 
2152  for (label iter = 0; iter < snapParams.nSmoothDispl(); iter++)
2153  {
2154  if ((iter % 10) == 0)
2155  {
2156  Info<< "Iteration " << iter << endl;
2157  }
2158  pointVectorField oldDisp(disp);
2159  meshMover.smooth(oldDisp, edgeGamma, disp);
2160  }
2161  Info<< "Displacement smoothed in = "
2162  << mesh.time().cpuTimeIncrement() << " s\n" << nl << endl;
2163 
2165  {
2166  const_cast<Time&>(mesh.time())++;
2167  Info<< "Writing smoothed mesh to time " << meshRefiner_.timeName()
2168  << endl;
2169 
2170  // Moving mesh creates meshPhi. Can be cleared out by a mesh.clearOut
2171  // but this will also delete all pointMesh but not pointFields which
2172  // gives an illegal situation.
2173 
2174  meshRefiner_.write
2175  (
2178  (
2181  ),
2182  mesh.time().path()/meshRefiner_.timeName()
2183  );
2184  Info<< "Writing displacement field ..." << endl;
2185  disp.write();
2186  tmp<pointScalarField> magDisp(mag(disp));
2187  magDisp().write();
2188 
2189  Info<< "Writing actual patch displacement ..." << endl;
2190  vectorField actualPatchDisp(disp, pp.meshPoints());
2191  dumpMove
2192  (
2193  mesh.time().path()
2194  / "actualPatchDisplacement_" + meshRefiner_.timeName() + ".obj",
2195  pp.localPoints(),
2196  pp.localPoints() + actualPatchDisp
2197  );
2198  }
2199 }
2200 
2201 
2203 (
2204  const snapParameters& snapParams,
2205  const label nInitErrors,
2206  const List<labelPair>& baffles,
2207  motionSmoother& meshMover
2208 )
2209 {
2210  addProfiling(scale, "snappyHexMesh::snap::scale");
2211  const fvMesh& mesh = meshRefiner_.mesh();
2212 
2213  // Relax displacement until correct mesh
2214  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2215  labelList checkFaces(identity(mesh.nFaces()));
2216 
2217  scalar oldErrorReduction = -1;
2218 
2219  bool meshOk = false;
2220 
2221  Info<< "Moving mesh ..." << endl;
2222  for (label iter = 0; iter < 2*snapParams.nSnap(); iter++)
2223  {
2224  Info<< nl << "Iteration " << iter << endl;
2225 
2226  if (iter == snapParams.nSnap())
2227  {
2228  Info<< "Displacement scaling for error reduction set to 0." << endl;
2229  oldErrorReduction = meshMover.setErrorReduction(0.0);
2230  }
2231 
2232  meshOk = meshMover.scaleMesh(checkFaces, baffles, true, nInitErrors);
2233 
2234  if (meshOk)
2235  {
2236  Info<< "Successfully moved mesh" << endl;
2237  break;
2238  }
2240  {
2241  const_cast<Time&>(mesh.time())++;
2242  Info<< "Writing scaled mesh to time " << meshRefiner_.timeName()
2243  << endl;
2244  mesh.write();
2245 
2246  Info<< "Writing displacement field ..." << endl;
2247  meshMover.displacement().write();
2248  tmp<pointScalarField> magDisp(mag(meshMover.displacement()));
2249  magDisp().write();
2250  }
2251  }
2252 
2253  if (oldErrorReduction >= 0)
2254  {
2255  meshMover.setErrorReduction(oldErrorReduction);
2256  }
2257  Info<< "Moved mesh in = "
2258  << mesh.time().cpuTimeIncrement() << " s\n" << nl << endl;
2259 
2260  return meshOk;
2261 }
2262 
2263 
2264 // After snapping: correct patching according to nearest surface.
2265 // Code is very similar to calcNearestSurface.
2266 // - calculate face-wise snap distance as max of point-wise
2267 // - calculate face-wise nearest surface point
2268 // - repatch face according to patch for surface point.
2270 (
2271  const snapParameters& snapParams,
2272  const labelList& adaptPatchIDs,
2273  const labelList& preserveFaces
2274 )
2275 {
2276  const fvMesh& mesh = meshRefiner_.mesh();
2277  const refinementSurfaces& surfaces = meshRefiner_.surfaces();
2278 
2279  Info<< "Repatching faces according to nearest surface ..." << endl;
2280 
2281  // Get the labels of added patches.
2283  (
2285  (
2286  mesh,
2287  adaptPatchIDs
2288  )
2289  );
2290  indirectPrimitivePatch& pp = ppPtr();
2291 
2292  // Divide surfaces into zoned and unzoned
2293  labelList zonedSurfaces =
2295  labelList unzonedSurfaces =
2297 
2298 
2299  // Faces that do not move
2300  bitSet isZonedFace(mesh.nFaces());
2301  {
2302  // 1. Preserve faces in preserveFaces list
2303  forAll(preserveFaces, facei)
2304  {
2305  if (preserveFaces[facei] != -1)
2306  {
2307  isZonedFace.set(facei);
2308  }
2309  }
2310 
2311  // 2. All faces on zoned surfaces
2312  const PtrList<surfaceZonesInfo>& surfZones = surfaces.surfZones();
2313  const faceZoneMesh& fZones = mesh.faceZones();
2314 
2315  forAll(zonedSurfaces, i)
2316  {
2317  const label zoneSurfi = zonedSurfaces[i];
2318  const wordList& fZoneNames = surfZones[zoneSurfi].faceZoneNames();
2319  forAll(fZoneNames, i)
2320  {
2321  const faceZone& fZone = fZones[fZoneNames[i]];
2322  isZonedFace.set(fZone);
2323  }
2324  }
2325  }
2326 
2327 
2328  // Determine per pp face which patch it should be in
2329  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2330 
2331  // Patch that face should be in
2332  labelList closestPatch(pp.size(), -1);
2333  {
2334  // face snap distance as max of point snap distance
2335  scalarField faceSnapDist(pp.size(), -GREAT);
2336  {
2337  // Distance to attract to nearest feature on surface
2338  const scalarField snapDist
2339  (
2340  calcSnapDistance
2341  (
2342  mesh,
2343  snapParams,
2344  pp
2345  )
2346  );
2347 
2348  const faceList& localFaces = pp.localFaces();
2349 
2350  forAll(localFaces, facei)
2351  {
2352  const face& f = localFaces[facei];
2353 
2354  forAll(f, fp)
2355  {
2356  faceSnapDist[facei] = max
2357  (
2358  faceSnapDist[facei],
2359  snapDist[f[fp]]
2360  );
2361  }
2362  }
2363  }
2364 
2365  pointField localFaceCentres(mesh.faceCentres(), pp.addressing());
2366 
2367  // Get nearest surface and region
2368  labelList hitSurface;
2369  labelList hitRegion;
2370  surfaces.findNearestRegion
2371  (
2372  unzonedSurfaces,
2373  localFaceCentres,
2374  sqr(faceSnapDist), // sqr of attract distance
2375  hitSurface,
2376  hitRegion
2377  );
2378 
2379  // Get patch
2380  forAll(pp, i)
2381  {
2382  label facei = pp.addressing()[i];
2383 
2384  if (hitSurface[i] != -1 && !isZonedFace.test(facei))
2385  {
2386  closestPatch[i] = globalToMasterPatch_
2387  [
2388  surfaces.globalRegion
2389  (
2390  hitSurface[i],
2391  hitRegion[i]
2392  )
2393  ];
2394  }
2395  }
2396  }
2397 
2398 
2399  // Change those faces for which there is a different closest patch
2400  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2401 
2402  labelList ownPatch(mesh.nFaces(), -1);
2403  labelList neiPatch(mesh.nFaces(), -1);
2404 
2406 
2407  forAll(patches, patchi)
2408  {
2409  const polyPatch& pp = patches[patchi];
2410 
2411  forAll(pp, i)
2412  {
2413  ownPatch[pp.start()+i] = patchi;
2414  neiPatch[pp.start()+i] = patchi;
2415  }
2416  }
2417 
2418  label nChanged = 0;
2419  forAll(closestPatch, i)
2420  {
2421  label facei = pp.addressing()[i];
2422 
2423  if (closestPatch[i] != -1 && closestPatch[i] != ownPatch[facei])
2424  {
2425  ownPatch[facei] = closestPatch[i];
2426  neiPatch[facei] = closestPatch[i];
2427  nChanged++;
2428  }
2429  }
2430 
2431  Info<< "Repatched " << returnReduce(nChanged, sumOp<label>())
2432  << " faces in = " << mesh.time().cpuTimeIncrement() << " s\n" << nl
2433  << endl;
2434 
2435  return meshRefiner_.createBaffles(ownPatch, neiPatch);
2436 }
2437 
2438 
2439 void Foam::snappySnapDriver::detectWarpedFaces
2440 (
2441  const scalar featureCos,
2442  const indirectPrimitivePatch& pp,
2443 
2444  DynamicList<label>& splitFaces,
2445  DynamicList<labelPair>& splits
2446 ) const
2447 {
2448  const fvMesh& mesh = meshRefiner_.mesh();
2449  const faceList& localFaces = pp.localFaces();
2450  const pointField& localPoints = pp.localPoints();
2451  const labelList& bFaces = pp.addressing();
2452 
2453  splitFaces.clear();
2454  splitFaces.setCapacity(bFaces.size());
2455  splits.clear();
2456  splits.setCapacity(bFaces.size());
2457 
2458  // Determine parallel consistent normals on points
2459  const vectorField pointNormals(PatchTools::pointNormals(mesh, pp));
2460 
2461  face f0(4);
2462  face f1(4);
2463 
2464  forAll(localFaces, facei)
2465  {
2466  const face& f = localFaces[facei];
2467 
2468  if (f.size() >= 4)
2469  {
2470  // See if splitting face across diagonal would make two faces with
2471  // biggish normal angle
2472 
2473  labelPair minDiag(-1, -1);
2474  scalar minCos(GREAT);
2475 
2476  for (label startFp = 0; startFp < f.size()-2; startFp++)
2477  {
2478  label minFp = f.rcIndex(startFp);
2479 
2480  for
2481  (
2482  label endFp = f.fcIndex(f.fcIndex(startFp));
2483  endFp < f.size() && endFp != minFp;
2484  endFp++
2485  )
2486  {
2487  // Form two faces
2488  f0.setSize(endFp-startFp+1);
2489  label i0 = 0;
2490  for (label fp = startFp; fp <= endFp; fp++)
2491  {
2492  f0[i0++] = f[fp];
2493  }
2494  f1.setSize(f.size()+2-f0.size());
2495  label i1 = 0;
2496  for (label fp = endFp; fp != startFp; fp = f.fcIndex(fp))
2497  {
2498  f1[i1++] = f[fp];
2499  }
2500  f1[i1++] = f[startFp];
2501 
2502  //Info<< "Splitting face:" << f << " into f0:" << f0
2503  // << " f1:" << f1 << endl;
2504 
2505  const vector n0 = f0.areaNormal(localPoints);
2506  const scalar n0Mag = mag(n0);
2507 
2508  const vector n1 = f1.areaNormal(localPoints);
2509  const scalar n1Mag = mag(n1);
2510 
2511  if (n0Mag > ROOTVSMALL && n1Mag > ROOTVSMALL)
2512  {
2513  scalar cosAngle = (n0/n0Mag) & (n1/n1Mag);
2514  if (cosAngle < minCos)
2515  {
2516  minCos = cosAngle;
2517  minDiag = labelPair(startFp, endFp);
2518  }
2519  }
2520  }
2521  }
2522 
2523 
2524  if (minCos < featureCos)
2525  {
2526  splitFaces.append(bFaces[facei]);
2527  splits.append(minDiag);
2528  }
2529  }
2530  }
2531 }
2532 
2533 
2534 Foam::labelList Foam::snappySnapDriver::getInternalOrBaffleDuplicateFace() const
2535 {
2536  const fvMesh& mesh = meshRefiner_.mesh();
2537 
2538  labelList internalOrBaffleFaceZones;
2539  {
2540  List<surfaceZonesInfo::faceZoneType> fzTypes(2);
2541  fzTypes[0] = surfaceZonesInfo::INTERNAL;
2542  fzTypes[1] = surfaceZonesInfo::BAFFLE;
2543  internalOrBaffleFaceZones = meshRefiner_.getZones(fzTypes);
2544  }
2545 
2546  List<labelPair> baffles
2547  (
2548  meshRefiner_.subsetBaffles
2549  (
2550  mesh,
2551  internalOrBaffleFaceZones,
2553  )
2554  );
2555 
2556  labelList faceToDuplicate(mesh.nFaces(), -1);
2557  forAll(baffles, i)
2558  {
2559  const labelPair& p = baffles[i];
2560  faceToDuplicate[p[0]] = p[1];
2561  faceToDuplicate[p[1]] = p[0];
2562  }
2564  return faceToDuplicate;
2565 }
2566 
2567 
2569 (
2570  const dictionary& snapDict,
2571  const dictionary& motionDict,
2572  const meshRefinement::FaceMergeType mergeType,
2573  const scalar featureCos,
2574  const scalar planarAngle,
2575  const snapParameters& snapParams
2576 )
2577 {
2578  addProfiling(snap, "snappyHexMesh::snap");
2579  fvMesh& mesh = meshRefiner_.mesh();
2580 
2581  Info<< nl
2582  << "Morphing phase" << nl
2583  << "--------------" << nl
2584  << endl;
2585 
2586  // faceZone handling
2587  // ~~~~~~~~~~~~~~~~~
2588  //
2589  // We convert all faceZones into baffles during snapping so we can use
2590  // a standard mesh motion (except for the mesh checking which for baffles
2591  // created from internal faces should check across the baffles). The state
2592  // is stored in two variables:
2593  // baffles : pairs of boundary faces
2594  // duplicateFace : from mesh face to its baffle colleague (or -1 for
2595  // normal faces)
2596  // There are three types of faceZones according to the faceType property:
2597  //
2598  // internal
2599  // --------
2600  // - baffles: need to be checked across
2601  // - duplicateFace: from face to duplicate face. Contains
2602  // all faces on faceZone to prevents merging patch faces.
2603  //
2604  // baffle
2605  // ------
2606  // - baffles: no need to be checked across
2607  // - duplicateFace: contains all faces on faceZone to prevent
2608  // merging patch faces.
2609  //
2610  // boundary
2611  // --------
2612  // - baffles: no need to be checked across. Also points get duplicated
2613  // so will no longer be baffles
2614  // - duplicateFace: contains no faces on faceZone since both sides can
2615  // merge faces independently.
2616 
2617 
2618 
2619  // faceZones of type internal
2620  const labelList internalFaceZones
2621  (
2622  meshRefiner_.getZones
2623  (
2625  (
2626  1,
2628  )
2629  )
2630  );
2631 
2632 
2633  // Create baffles (pairs of faces that share the same points)
2634  // Baffles stored as owner and neighbour face that have been created.
2635  {
2636  List<labelPair> baffles;
2637  labelList originatingFaceZone;
2638  meshRefiner_.createZoneBaffles
2639  (
2640  identity(mesh.faceZones().size()),
2641  baffles,
2642  originatingFaceZone
2643  );
2644  }
2645 
2646  // Duplicate points on faceZones of type boundary
2647  meshRefiner_.dupNonManifoldBoundaryPoints();
2648 
2649 
2650  bool doFeatures = false;
2651  label nFeatIter = 1;
2652  if (snapParams.nFeatureSnap() > 0)
2653  {
2654  doFeatures = true;
2655 
2656  if (!dryRun_)
2657  {
2658  nFeatIter = snapParams.nFeatureSnap();
2659  }
2660 
2661  Info<< "Snapping to features in " << nFeatIter
2662  << " iterations ..." << endl;
2663  }
2664 
2665 
2666  bool meshOk = false;
2667 
2668 
2669  // Get the labels of added patches.
2670  labelList adaptPatchIDs(meshRefiner_.meshedPatches());
2671 
2672 
2673 
2674  {
2675  autoPtr<indirectPrimitivePatch> ppPtr
2676  (
2678  (
2679  mesh,
2680  adaptPatchIDs
2681  )
2682  );
2683 
2684 
2685  // Distance to attract to nearest feature on surface
2686  scalarField snapDist(calcSnapDistance(mesh, snapParams, ppPtr()));
2687 
2688 
2689  // Construct iterative mesh mover.
2690  Info<< "Constructing mesh displacer ..." << endl;
2691  Info<< "Using mesh parameters " << motionDict << nl << endl;
2692 
2693  autoPtr<motionSmoother> meshMoverPtr
2694  (
2695  new motionSmoother
2696  (
2697  mesh,
2698  ppPtr(),
2699  adaptPatchIDs,
2701  (
2703  adaptPatchIDs
2704  ),
2705  motionDict,
2706  dryRun_
2707  )
2708  );
2709 
2710 
2711  // Check initial mesh
2712  Info<< "Checking initial mesh ..." << endl;
2713  labelHashSet wrongFaces(mesh.nFaces()/100);
2714  motionSmoother::checkMesh(false, mesh, motionDict, wrongFaces, dryRun_);
2715  const label nInitErrors = returnReduce
2716  (
2717  wrongFaces.size(),
2718  sumOp<label>()
2719  );
2720 
2721  Info<< "Detected " << nInitErrors << " illegal faces"
2722  << " (concave, zero area or negative cell pyramid volume)"
2723  << endl;
2724 
2725 
2726  Info<< "Checked initial mesh in = "
2727  << mesh.time().cpuTimeIncrement() << " s\n" << nl << endl;
2728 
2729  // Extract baffles across internal faceZones (for checking mesh quality
2730  // across
2731  labelPairList internalBaffles
2732  (
2733  meshRefiner_.subsetBaffles
2734  (
2735  mesh,
2736  internalFaceZones,
2738  )
2739  );
2740 
2741 
2742 
2743  // Pre-smooth patch vertices (so before determining nearest)
2744  preSmoothPatch
2745  (
2746  meshRefiner_,
2747  snapParams,
2748  nInitErrors,
2749  internalBaffles,
2750  meshMoverPtr()
2751  );
2752 
2753  // TBD. Include re-patching?
2754 
2755 
2756  //- Only if in feature attraction mode:
2757  // Nearest feature
2758  vectorField patchAttraction;
2759  // Constraints at feature
2760  List<pointConstraint> patchConstraints;
2761 
2762 
2763  //- Any faces to split
2764  DynamicList<label> splitFaces;
2765  //- Indices in face to split across
2766  DynamicList<labelPair> splits;
2767 
2768 
2769  for (label iter = 0; iter < nFeatIter; iter++)
2770  {
2771  Info<< nl
2772  << "Morph iteration " << iter << nl
2773  << "-----------------" << endl;
2774 
2775  // Splitting iteration?
2776  bool doSplit = false;
2777  if
2778  (
2779  doFeatures
2780  && snapParams.nFaceSplitInterval() > 0
2781  && (
2782  (iter == nFeatIter-1)
2783  || (iter > 0 && (iter%snapParams.nFaceSplitInterval()) == 0)
2784  )
2785  )
2786  {
2787  doSplit = true;
2788  }
2789 
2790 
2791 
2792  indirectPrimitivePatch& pp = ppPtr();
2793  motionSmoother& meshMover = meshMoverPtr();
2794 
2795 
2796  // Calculate displacement at every patch point if we need it:
2797  // - if automatic near-surface detection
2798  // - if face splitting active
2799  pointField nearestPoint;
2800  vectorField nearestNormal;
2801 
2802  if (snapParams.detectNearSurfacesSnap() || doSplit)
2803  {
2804  nearestPoint.setSize(pp.nPoints(), vector::max);
2805  nearestNormal.setSize(pp.nPoints(), Zero);
2806  }
2807 
2808  vectorField disp = calcNearestSurface
2809  (
2810  snapParams.strictRegionSnap(), // attract points to region only
2811  meshRefiner_,
2812  globalToMasterPatch_, // for if strictRegionSnap
2813  globalToSlavePatch_, // for if strictRegionSnap
2814  snapDist,
2815  pp,
2816 
2817  nearestPoint,
2818  nearestNormal
2819  );
2820 
2821 
2822  // Override displacement at thin gaps
2823  if (snapParams.detectNearSurfacesSnap())
2824  {
2825  detectNearSurfaces
2826  (
2827  Foam::cos(degToRad(planarAngle)),// planar cos for gaps
2828  pp,
2829  nearestPoint, // surfacepoint from nearest test
2830  nearestNormal, // surfacenormal from nearest test
2831 
2832  disp
2833  );
2834  }
2835 
2836  // Override displacement with feature edge attempt
2837  if (doFeatures)
2838  {
2839  splitFaces.clear();
2840  splits.clear();
2841  disp = calcNearestSurfaceFeature
2842  (
2843  snapParams,
2844  !doSplit, // alignMeshEdges
2845  iter,
2846  featureCos,
2847  scalar(iter+1)/nFeatIter,
2848 
2849  snapDist,
2850  disp,
2851  nearestNormal,
2852  meshMover,
2853 
2854  patchAttraction,
2855  patchConstraints,
2856 
2857  splitFaces,
2858  splits
2859  );
2860  }
2861 
2862  // Check for displacement being outwards.
2863  outwardsDisplacement(pp, disp);
2864 
2865  // Freeze points on exposed points/faces
2866  freezeExposedPoints
2867  (
2868  meshRefiner_,
2869  "frozenFaces", // faceZone name
2870  "frozenPoints", // pointZone name
2871  pp,
2872  disp
2873  );
2874 
2875  // Set initial distribution of displacement field (on patches)
2876  // from patchDisp and make displacement consistent with b.c.
2877  // on displacement pointVectorField.
2878  meshMover.setDisplacement(disp);
2879 
2880 
2882  {
2883  dumpMove
2884  (
2885  mesh.time().path()
2886  / "patchDisplacement_" + name(iter) + ".obj",
2887  pp.localPoints(),
2888  pp.localPoints() + disp
2889  );
2890  }
2891 
2892  // Get smoothly varying internal displacement field.
2893  smoothDisplacement(snapParams, meshMover);
2894 
2895  // Apply internal displacement to mesh.
2896  meshOk = scaleMesh
2897  (
2898  snapParams,
2899  nInitErrors,
2900  internalBaffles,
2901  meshMover
2902  );
2903 
2904  if (!meshOk)
2905  {
2907  << "Did not successfully snap mesh."
2908  << " Continuing to snap to resolve easy" << nl
2909  << " surfaces but the"
2910  << " resulting mesh will not satisfy your quality"
2911  << " constraints" << nl << endl;
2912  }
2913 
2915  {
2916  const_cast<Time&>(mesh.time())++;
2917  Info<< "Writing scaled mesh to time "
2918  << meshRefiner_.timeName() << endl;
2919  meshRefiner_.write
2920  (
2923  (
2926  ),
2927  mesh.time().path()/meshRefiner_.timeName()
2928  );
2929  Info<< "Writing displacement field ..." << endl;
2930  meshMover.displacement().write();
2931  tmp<pointScalarField> magDisp
2932  (
2933  mag(meshMover.displacement())
2934  );
2935  magDisp().write();
2936  }
2937 
2938  // Use current mesh as base mesh
2939  meshMover.correct();
2940 
2941 
2942 
2943  // See if any faces need splitting
2944  label nTotalSplit = returnReduce(splitFaces.size(), sumOp<label>());
2945  if (nTotalSplit && doSplit)
2946  {
2947  // Filter out baffle faces from faceZones of type
2948  // internal/baffle
2949 
2950  labelList duplicateFace(getInternalOrBaffleDuplicateFace());
2951 
2952  {
2953  labelList oldSplitFaces(std::move(splitFaces));
2954  List<labelPair> oldSplits(std::move(splits));
2955  forAll(oldSplitFaces, i)
2956  {
2957  if (duplicateFace[oldSplitFaces[i]] == -1)
2958  {
2959  splitFaces.append(oldSplitFaces[i]);
2960  splits.append(oldSplits[i]);
2961  }
2962  }
2963  nTotalSplit = returnReduce
2964  (
2965  splitFaces.size(),
2966  sumOp<label>()
2967  );
2968  }
2969 
2970  // Update mesh
2971  meshRefiner_.splitFacesUndo
2972  (
2973  splitFaces,
2974  splits,
2975  motionDict,
2976 
2977  duplicateFace,
2978  internalBaffles
2979  );
2980 
2981  // Redo meshMover
2982  meshMoverPtr.clear();
2983  ppPtr.clear();
2984 
2985  // Update mesh mover
2986  ppPtr = meshRefinement::makePatch(mesh, adaptPatchIDs);
2987  meshMoverPtr.reset
2988  (
2989  new motionSmoother
2990  (
2991  mesh,
2992  ppPtr(),
2993  adaptPatchIDs,
2995  (
2997  adaptPatchIDs
2998  ),
2999  motionDict,
3000  dryRun_
3001  )
3002  );
3003 
3004  // Update snapping distance
3005  snapDist = calcSnapDistance(mesh, snapParams, ppPtr());
3006 
3007 
3009  {
3010  const_cast<Time&>(mesh.time())++;
3011  Info<< "Writing split-faces mesh to time "
3012  << meshRefiner_.timeName() << endl;
3013  meshRefiner_.write
3014  (
3017  (
3020  ),
3021  mesh.time().path()/meshRefiner_.timeName()
3022  );
3023  }
3024  }
3025 
3026 
3028  {
3029  forAll(internalBaffles, i)
3030  {
3031  const labelPair& p = internalBaffles[i];
3032  const point& fc0 = mesh.faceCentres()[p[0]];
3033  const point& fc1 = mesh.faceCentres()[p[1]];
3034 
3035  if (mag(fc0-fc1) > meshRefiner_.mergeDistance())
3036  {
3038  << "Separated baffles : f0:" << p[0]
3039  << " centre:" << fc0
3040  << " f1:" << p[1] << " centre:" << fc1
3041  << " distance:" << mag(fc0-fc1)
3042  << exit(FatalError);
3043  }
3044  }
3045  }
3046  }
3047  }
3048 
3049 
3050  // Merge any introduced baffles (from faceZones of faceType 'internal')
3051  {
3052  autoPtr<mapPolyMesh> mapPtr = meshRefiner_.mergeZoneBaffles
3053  (
3054  true, // internal zones
3055  false // baffle zones
3056  );
3057 
3058  if (mapPtr)
3059  {
3061  {
3062  const_cast<Time&>(mesh.time())++;
3063  Info<< "Writing baffle-merged mesh to time "
3064  << meshRefiner_.timeName() << endl;
3065  meshRefiner_.write
3066  (
3069  (
3072  ),
3073  meshRefiner_.timeName()
3074  );
3075  }
3076  }
3077  }
3078 
3079  // Repatch faces according to nearest. Do not repatch baffle faces.
3080  {
3081  labelList duplicateFace(getInternalOrBaffleDuplicateFace());
3082 
3083  repatchToSurface(snapParams, adaptPatchIDs, duplicateFace);
3084  }
3085 
3086  if
3087  (
3088  mergeType == meshRefinement::FaceMergeType::GEOMETRIC
3089  || mergeType == meshRefinement::FaceMergeType::IGNOREPATCH
3090  )
3091  {
3092  labelList duplicateFace(getInternalOrBaffleDuplicateFace());
3093 
3094  // Repatching might have caused faces to be on same patch and hence
3095  // mergeable so try again to merge coplanar faces. Do not merge baffle
3096  // faces to ensure they both stay the same.
3097  label nChanged = meshRefiner_.mergePatchFacesUndo
3098  (
3099  featureCos, // minCos
3100  featureCos, // concaveCos
3101  meshRefiner_.meshedPatches(),
3102  motionDict,
3103  duplicateFace, // faces not to merge
3104  mergeType
3105  );
3106 
3107  nChanged += meshRefiner_.mergeEdgesUndo(featureCos, motionDict);
3108 
3109  if (nChanged > 0 && debug & meshRefinement::MESH)
3110  {
3111  const_cast<Time&>(mesh.time())++;
3112  Info<< "Writing patchFace merged mesh to time "
3113  << meshRefiner_.timeName() << endl;
3114  meshRefiner_.write
3115  (
3118  (
3121  ),
3122  meshRefiner_.timeName()
3123  );
3124  }
3125  }
3126 
3128  {
3129  const_cast<Time&>(mesh.time())++;
3130  }
3131 }
3132 
3133 
3134 // ************************************************************************* //
Given a displacement moves the mesh by scaling the displacement back until there are no more mesh err...
faceListList boundary
const labelListList & pointEdges() const
Return point-edge addressing.
Switch detectNearSurfacesSnap() const
Override attraction to nearest with intersection location.
static bool checkMesh(const bool report, const polyMesh &mesh, const dictionary &dict, labelHashSet &wrongFaces, const bool dryRun=false)
Check mesh with mesh settings in dict. Collects incorrect faces.
label nPoints() const
Number of points supporting patch faces.
const List< face_type > & localFaces() const
Return patch faces addressing into local point list.
reference val() const
Const access to referenced object (value)
Definition: HashTable.H:1185
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
ZoneMesh< faceZone, polyMesh > faceZoneMesh
A ZoneMesh with the type faceZone.
fileName path() const
Return path = rootPath/caseName. Same as TimePaths::path()
Definition: Time.H:503
label nSmoothDispl() const
Number of mesh displacement smoothing iterations.
const labelList & surfaces() const
bool scaleMesh(const snapParameters &snapParams, const label nInitErrors, const List< labelPair > &baffles, motionSmoother &)
Do the hard work: move the mesh according to displacement,.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:68
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
const Field< point_type > & localPoints() const
Return pointField of points in patch.
static bitSet getMasterFaces(const polyMesh &mesh)
Get per face whether it is uncoupled or a master of a coupled set of faces.
Definition: syncTools.C:119
label nPoints() const noexcept
Number of mesh points.
void doSnap(const dictionary &snapDict, const dictionary &motionDict, const meshRefinement::FaceMergeType mergeType, const scalar featureCos, const scalar planarAngle, const snapParameters &snapParams)
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
Type gMin(const FieldField< Field, Type > &f)
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1122
static const weightedPosition zero
static const pointMesh & New(const polyMesh &mesh, Args &&... args)
Get existing or create a new MeshObject. Registered with typeName.
Definition: MeshObject.C:53
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
std::enable_if< std::is_same< bool, TypeT >::value, bool >::type set(const label i, bool val=true)
A bitSet::set() method for a list of bool.
Definition: List.H:489
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Unit conversion functions.
bool scaleMesh(labelList &checkFaces, const bool smoothMesh=true, const label nAllow=0)
Move mesh with given scale. Return true if mesh ok or has.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
const pointMesh & pMesh() const
Reference to pointMesh.
label nSmoothPatch() const
Number of patch smoothing iterations before finding.
void detectNearSurfaces(const scalar planarCos, const indirectPrimitivePatch &, const pointField &nearestPoint, const vectorField &nearestNormal, vectorField &disp) const
Per patch point override displacement if in gap situation.
T & first()
Access first element of the list, position [0].
Definition: UList.H:853
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
static writeType writeLevel()
Get/set write level.
GeometricField< vector, pointPatchField, pointMesh > pointVectorField
void setCapacity(const label len)
Alter the size of the underlying storage.
Definition: DynamicListI.H:303
List< labelPair > labelPairList
List of labelPair.
Definition: labelPair.H:33
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of a point.
Definition: meshTools.C:196
static void getPoints(const UList< weightedPosition > &in, List< point > &out)
Get points.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
label nFaces() const noexcept
Number of mesh faces.
PrimitivePatch< IndirectList< face >, const pointField & > indirectPrimitivePatch
A PrimitivePatch with an IndirectList for the faces, const reference for the point field...
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
List< labelList > labelListList
List of labelList.
Definition: labelList.H:38
void findNearest(const labelList &surfacesToTest, const pointField &samples, const scalarField &nearestDistSqr, labelList &surfaces, List< pointIndexHit > &) const
Find nearest point on surfaces.
label nFeatureSnap() const
const labelList & meshPoints() const
Return labelList of mesh points in patch.
void smooth(volScalarField &field, const scalar coeff)
Definition: fvcSmooth.C:37
const fvMesh & mesh() const
Reference to mesh.
static labelPairList findDuplicateFacePairs(const polyMesh &)
Helper routine to find all baffles (two boundary faces.
label nTotalPoints() const noexcept
Return total number of points in decomposed mesh. Not compensated for duplicate points! ...
Container for data on surfaces used for surface-driven refinement. Contains all the data about the le...
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1078
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
void correct()
Take over existing mesh position.
label fcIndex(const label i) const noexcept
The forward circular index. The next index in the list which returns to the first at the end of the l...
Definition: UListI.H:97
label mergePoints(const PointList &points, labelList &pointToUnique, labelList &uniquePoints, const scalar mergeTol=SMALL, const bool verbose=false)
Calculate merge mapping, preserving the original point order. All points closer/equal mergeTol are to...
void reset(T *p=nullptr) noexcept
Delete managed object and set to new given pointer.
Definition: autoPtrI.H:37
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition: HashSet.H:85
const Field< point_type > & faceNormals() const
Return face unit normals for patch.
List< face > faceList
List of faces.
Definition: faceListFwd.H:39
A list of faces which address into the list of points.
const Map< label > & meshPointMap() const
Mesh point map.
#define addProfiling(Name, Descr)
Define profiling trigger with specified name and description string.
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
void setSize(const label n)
Alias for resize()
Definition: List.H:316
scalar setErrorReduction(const scalar)
Set the errorReduction (by how much to scale the displacement.
FaceMergeType
Enumeration for what to do with co-planar patch faces on a single.
dynamicFvMesh & mesh
dimensionedScalar cos(const dimensionedScalar &ds)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
const pointField & points
void findNearestRegion(const labelList &surfacesToTest, const pointField &samples, const scalarField &nearestDistSqr, labelList &hitSurface, labelList &hitRegion) const
Find nearest point on surfaces. Return surface and region on.
label nFaceSplitInterval() const
label findZoneID(const word &zoneName) const
Find zone index by name, return -1 if not found.
Definition: ZoneMesh.C:629
const_iterator cfind(const Key &key) const
Find and return an const_iterator set at the hashed entry.
Definition: HashTableI.H:113
An edge is a list of two vertex labels. This can correspond to a directed graph edge or an edge on a ...
Definition: edge.H:59
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
Definition: labelLists.C:44
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition: polyMesh.H:608
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Simple container to keep together snap specific information.
autoPtr< mapPolyMesh > repatchToSurface(const snapParameters &snapParams, const labelList &adaptPatchIDs, const labelList &preserveFaces)
Repatch faces according to surface nearest the face centre.
label size() const noexcept
The number of entries in the list.
Definition: UPtrListI.H:106
const Field< point_type > & faceCentres() const
Return face centres for patch.
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1116
static void setDisplacement(const labelList &patchIDs, const indirectPrimitivePatch &pp, pointField &patchDisp, pointVectorField &displacement)
Set displacement field from displacement on patch points.
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
void smoothDisplacement(const snapParameters &snapParams, motionSmoother &) const
Smooth the displacement field to the internal.
const globalMeshData & globalData() const
Return parallel info (demand-driven)
Definition: polyMesh.C:1311
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
static void preSmoothPatch(const meshRefinement &meshRefiner, const snapParameters &snapParams, const label nInitErrors, const List< labelPair > &baffles, motionSmoother &)
Smooth the mesh (patch and internal) to increase visibility.
label nInternalFaces() const noexcept
Number of internal faces.
Vector< scalar > vector
Definition: vector.H:57
ZoneMesh< pointZone, polyMesh > pointZoneMesh
A ZoneMesh with the type pointZone.
const Field< point_type > & points() const noexcept
Return reference to global points.
void append(const T &val)
Copy append an element to the end of this list.
Definition: DynamicList.H:584
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1103
const PtrList< surfaceZonesInfo > & surfZones() const
const vectorField & cellCentres() const
bool visNormal(const vector &n, const vectorField &faceNormals, const labelList &faceLabels)
Check if n is in same direction as normals of all faceLabels.
Definition: meshTools.C:30
static labelList getUnnamedSurfaces(const PtrList< surfaceZonesInfo > &surfList)
Get indices of unnamed surfaces (surfaces without faceZoneName)
errorManip< error > abort(error &err)
Definition: errorManip.H:139
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie...
Definition: PtrList.H:159
label nEdges() const
Number of mesh edges.
static void syncPointList(const polyMesh &mesh, List< T > &pointValues, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronize values on all mesh points.
const labelListList & pointCells() const
virtual bool write(const bool writeOnProc=true) const
Write mesh using IO settings from time.
Definition: fvMesh.C:1102
void smooth(const GeometricField< Type, pointPatchField, pointMesh > &fld, const scalarField &edgeWeight, GeometricField< Type, pointPatchField, pointMesh > &newFld) const
Fully explicit smoothing of fields (not positions)
word timeName() const
Replacement for Time::timeName() that returns oldInstance (if overwrite_)
static scalarField calcSnapDistance(const fvMesh &mesh, const snapParameters &snapParams, const indirectPrimitivePatch &)
Calculate edge length per patch point.
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO...
const labelListList & pointFaces() const
Return point-face addressing.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
An OFstream that keeps track of vertices and provides convenience output methods for OBJ files...
Definition: OBJstream.H:55
int debug
Static debugging option.
Pair< label > labelPair
A pair of labels.
Definition: Pair.H:51
Type gMax(const FieldField< Field, Type > &f)
constexpr auto end(C &c) -> decltype(c.end())
Return iterator to the end of the container c.
Definition: stdFoam.H:201
const faceZoneMesh & faceZones() const noexcept
Return face zone mesh.
Definition: polyMesh.H:670
defineTypeNameAndDebug(combustionModel, 0)
const indirectPrimitivePatch & patch() const
Reference to patch.
static tmp< pointField > avgCellCentres(const fvMesh &mesh, const indirectPrimitivePatch &)
Helper: calculate average cell centre per point.
void clear() noexcept
Clear the addressed list, i.e. set the size to zero.
Definition: DynamicListI.H:405
Geometric merging of points. See below.
labelList f(nPoints)
static autoPtr< indirectPrimitivePatch > makePatch(const polyMesh &, const labelList &)
Create patch from set of patches.
Info<< "Finished reading KIVA file"<< endl;cellShapeList cellShapes(nPoints);labelList cellZoning(nPoints, -1);const cellModel &hex=cellModel::ref(cellModel::HEX);labelList hexLabels(8);label activeCells=0;labelList pointMap(nPoints);forAll(pointMap, i){ pointMap[i]=i;}for(label i=0;i< nPoints;i++){ if(f[i] > 0.0) { hexLabels[0]=i;hexLabels[1]=i1tab[i];hexLabels[2]=i3tab[i1tab[i]];hexLabels[3]=i3tab[i];hexLabels[4]=i8tab[i];hexLabels[5]=i1tab[i8tab[i]];hexLabels[6]=i3tab[i1tab[i8tab[i]]];hexLabels[7]=i3tab[i8tab[i]];cellShapes[activeCells].reset(hex, hexLabels);edgeList edges=cellShapes[activeCells].edges();forAll(edges, ei) { if(edges[ei].mag(points)< SMALL) { label start=pointMap[edges[ei].start()];while(start !=pointMap[start]) { start=pointMap[start];} label end=pointMap[edges[ei].end()];while(end !=pointMap[end]) { end=pointMap[end];} label minLabel=min(start, end);pointMap[start]=pointMap[end]=minLabel;} } cellZoning[activeCells]=idreg[i];activeCells++;}}cellShapes.setSize(activeCells);cellZoning.setSize(activeCells);forAll(cellShapes, celli){ cellShape &cs=cellShapes[celli];forAll(cs, i) { cs[i]=pointMap[cs[i]];} cs.collapse();}label bcIDs[11]={-1, 0, 2, 4, -1, 5, -1, 6, 7, 8, 9};const label nBCs=12;const word *kivaPatchTypes[nBCs]={ &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &symmetryPolyPatch::typeName, &wedgePolyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &symmetryPolyPatch::typeName, &oldCyclicPolyPatch::typeName};enum patchTypeNames{ PISTON, VALVE, LINER, CYLINDERHEAD, AXIS, WEDGE, INFLOW, OUTFLOW, PRESIN, PRESOUT, SYMMETRYPLANE, CYCLIC};const char *kivaPatchNames[nBCs]={ "piston", "valve", "liner", "cylinderHead", "axis", "wedge", "inflow", "outflow", "presin", "presout", "symmetryPlane", "cyclic"};List< SLList< face > > pFaces[nBCs]
Definition: readKivaGrid.H:235
double cpuTimeIncrement() const
Return CPU time (in seconds) since last call to cpuTimeIncrement()
Definition: cpuTimePosix.C:80
messageStream Warning
Warning stream (stdout output on master, null elsewhere), with additional &#39;FOAM Warning&#39; header text...
Helper class which maintains intersections of (changing) mesh with (static) surfaces.
pointVectorField & pointDisplacement()
Return reference to the point motion displacement field.
static void swapBoundaryCellList(const polyMesh &mesh, const UList< T > &cellData, List< T > &neighbourCellData)
Swap to obtain neighbour cell values for all boundary faces.
const vectorField & faceCentres() const
const pointZoneMesh & pointZones() const noexcept
Return point zone mesh.
Definition: polyMesh.H:662
List< word > wordList
List of word.
Definition: fileName.H:59
label nSnap() const
Maximum number of snapping relaxation iterations. Should stop.
PrimitivePatch< UIndirectList< face >, const pointField & > uindirectPrimitivePatch
A PrimitivePatch with UIndirectList for the faces, const reference for the point field.
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field values.
scalar snapTol() const
Relative distance for points to be attracted by surface.
vector point
Point is a vector.
Definition: point.H:37
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:59
#define WarningInFunction
Report a warning using Foam::Warning.
Switch strictRegionSnap() const
Attract point to corresponding surface region only.
label newPointi
Definition: readKivaGrid.H:496
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
static bitSet getMasterPoints(const polyMesh &mesh, const labelList &meshPoints)
Determine master point for subset of points. If coupled.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
const polyBoundaryMesh & patches
label globalRegion(const label surfI, const label regionI) const
From surface and region on surface to global region.
static labelList getNamedSurfaces(const PtrList< surfaceZonesInfo > &surfList)
Get indices of named surfaces (surfaces with faceZoneName)
messageStream Info
Information stream (stdout output on master, null elsewhere)
label nTotalFaces() const noexcept
Return total number of faces in decomposed mesh. Not compensated for duplicate faces! ...
SubField< vector > subField
Declare type of subField.
Definition: Field.H:128
writeType
Enumeration for what to write. Used as a bit-pattern.
Field< vector > vectorField
Specialisation of Field<T> for vector.
debugType
Enumeration for what to debug. Used as a bit-pattern.
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:60
static void syncPoints(const polyMesh &mesh, List< weightedPosition > &)
Synchronisation for mesh point positions.
bool write() const
Write mesh and all data.
List< label > labelList
A List of labels.
Definition: List.H:62
static tmp< pointVectorField > makeDisplacementField(const pointMesh &pMesh, const labelList &adaptPatchIDs)
Helper function to make a pointVectorField with correct.
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
label nSmoothInternal() const
Number of internal point smoothing iterations (combined with.
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:69
bool coupled
static T gAverage(const bitSet &isMasterElem, const UList< T > &values)
Helper: calculate average.
void findNearestIntersection(const labelList &surfacesToTest, const pointField &start, const pointField &end, labelList &surface1, List< pointIndexHit > &hit1, labelList &region1, labelList &surface2, List< pointIndexHit > &hit2, labelList &region2) const
Find intersection nearest to the endpoints. surface1,2 are.
label rcIndex(const label i) const noexcept
The reverse circular index. The previous index in the list which returns to the last at the beginning...
Definition: UListI.H:104
static tmp< pointField > pointNormals(const polyMesh &, const PrimitivePatch< FaceList, PointField > &, const bitSet &flipMap=bitSet::null())
Return parallel consistent point normals for patches using mesh points.
bool returnReduceOr(const bool value, const label comm=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
constexpr scalar degToRad(const scalar deg) noexcept
Conversion from degrees to radians.
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Namespace for OpenFOAM.
pointVectorField & displacement()
Reference to displacement field.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127