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