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