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