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