Loading...
Searching...
No Matches
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-2025 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#include "addPatchCellLayer.H"
55#include "snappyLayerDriver.H"
56#include "IOmanip.H"
58// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
59
60namespace Foam
61{
62
64
65} // End namespace Foam
66
67
68// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
69
70// Calculate geometrically collocated points, Requires bitSet to be
71// sized and initialised!
72Foam::label Foam::snappySnapDriver::getCollocatedPoints
73(
74 const scalar tol,
75 const pointField& points,
76 bitSet& isCollocatedPoint
77)
78{
79 labelList pointMap;
80 label nUnique = Foam::mergePoints
81 (
82 points, // points
83 tol, // mergeTol
84 false, // verbose
85 pointMap
86 );
87 bool hasMerged = (nUnique < points.size());
88
89 if (!returnReduceOr(hasMerged))
90 {
91 return 0;
92 }
93
94 // Determine which merged points are referenced more than once
95 label nCollocated = 0;
96
97 // Per old point the newPoint. Or -1 (not set yet) or -2 (already seen
98 // twice)
99 labelList firstOldPoint(nUnique, -1);
100 forAll(pointMap, oldPointi)
101 {
102 label newPointi = pointMap[oldPointi];
103
104 if (firstOldPoint[newPointi] == -1)
105 {
106 // First use of oldPointi. Store.
107 firstOldPoint[newPointi] = oldPointi;
108 }
109 else if (firstOldPoint[newPointi] == -2)
110 {
111 // Third or more reference of oldPointi -> non-manifold
112 isCollocatedPoint.set(oldPointi);
113 nCollocated++;
114 }
115 else
116 {
117 // Second reference of oldPointi -> non-manifold
118 isCollocatedPoint.set(firstOldPoint[newPointi]);
119 nCollocated++;
120
121 isCollocatedPoint.set(oldPointi);
122 nCollocated++;
123
124 // Mark with special value to save checking next time round
125 firstOldPoint[newPointi] = -2;
126 }
127 }
128 return returnReduce(nCollocated, sumOp<label>());
129}
130
131
132Foam::tmp<Foam::pointField> Foam::snappySnapDriver::smoothInternalDisplacement
133(
134 const meshRefinement& meshRefiner,
135 const motionSmoother& meshMover
136)
137{
138 const indirectPrimitivePatch& pp = meshMover.patch();
139 const polyMesh& mesh = meshMover.mesh();
140
141 // Get neighbour refinement
142 const hexRef8& cutter = meshRefiner.meshCutter();
143 const labelList& cellLevel = cutter.cellLevel();
144
145
146 // Get the faces on the boundary
147 bitSet isFront(mesh.nFaces(), pp.addressing());
148
149 // Walk out from the surface a bit. Poor man's FaceCellWave.
150 // Commented out for now - not sure if needed and if so how much
151 //for (label iter = 0; iter < 2; iter++)
152 //{
153 // bitSet newIsFront(mesh.nFaces());
154 //
155 // forAll(isFront, facei)
156 // {
157 // if (isFront.test(facei))
158 // {
159 // label own = mesh.faceOwner()[facei];
160 // const cell& ownFaces = mesh.cells()[own];
161 // newIsFront.set(ownFaces);
162 //
163 // if (mesh.isInternalFace(facei))
164 // {
165 // label nei = mesh.faceNeighbour()[facei];
166 // const cell& neiFaces = mesh.cells()[nei];
167 // newIsFront.set(neiFaces);
168 // }
169 // }
170 // }
171 //
172 // syncTools::syncFaceList
173 // (
174 // mesh,
175 // newIsFront,
176 // orEqOp<unsigned int>()
177 // );
178 //
179 // isFront = newIsFront;
180 //}
181
182 // Mark all points on faces
183 // - not on the boundary
184 // - inbetween differing refinement levels
185 bitSet isMovingPoint(mesh.nPoints());
186
187 label nInterface = 0;
188
189 for (label facei = 0; facei < mesh.nInternalFaces(); facei++)
190 {
191 label ownLevel = cellLevel[mesh.faceOwner()[facei]];
192 label neiLevel = cellLevel[mesh.faceNeighbour()[facei]];
193
194 if (!isFront.test(facei) && ownLevel != neiLevel)
195 {
196 const face& f = mesh.faces()[facei];
197 isMovingPoint.set(f);
198
199 ++nInterface;
200 }
201 }
202
203 labelList neiCellLevel;
204 syncTools::swapBoundaryCellList(mesh, cellLevel, neiCellLevel);
205
206 for (label facei = mesh.nInternalFaces(); facei < mesh.nFaces(); facei++)
207 {
208 label ownLevel = cellLevel[mesh.faceOwner()[facei]];
209 label neiLevel = neiCellLevel[facei-mesh.nInternalFaces()];
210
211 if (!isFront.test(facei) && ownLevel != neiLevel)
212 {
213 const face& f = mesh.faces()[facei];
214 isMovingPoint.set(f);
215
216 ++nInterface;
217 }
218 }
219
220 if (debug)
221 {
222 Info<< "Found " << returnReduce(nInterface, sumOp<label>())
223 << " faces out of " << mesh.globalData().nTotalFaces()
224 << " inbetween refinement regions." << endl;
225 }
226
227 // Make sure that points that are coupled to a moving point are marked
228 // as well
230
231 // Unmark any point on the boundary. If we're doing zero iterations of
232 // face-cell wave we might have coupled points not being unmarked.
233 isMovingPoint.unset(pp.meshPoints());
234
235 // Make sure that points that are coupled to meshPoints but not on a patch
236 // are unmarked as well
238
239
240 // Calculate average of connected cells
241 Field<weightedPosition> sumLocation
242 (
243 mesh.nPoints(),
245 );
246
247 forAll(isMovingPoint, pointi)
248 {
249 if (isMovingPoint.test(pointi))
250 {
251 const labelList& pCells = mesh.pointCells(pointi);
252
253 sumLocation[pointi].first() = pCells.size();
254 for (const label celli : pCells)
255 {
256 sumLocation[pointi].second() += mesh.cellCentres()[celli];
257 }
258 }
259 }
260
261 // Add coupled contributions
263
264 auto tdisplacement = tmp<pointField>::New(mesh.nPoints(), Zero);
265 auto& displacement = tdisplacement.ref();
266
267 label nAdapted = 0;
268
269 forAll(displacement, pointi)
270 {
271 const weightedPosition& wp = sumLocation[pointi];
272 if (mag(wp.first()) > VSMALL)
273 {
274 displacement[pointi] =
275 wp.second()/wp.first()
276 - mesh.points()[pointi];
277 nAdapted++;
278 }
279 }
280
281 Info<< "Smoothing " << returnReduce(nAdapted, sumOp<label>())
282 << " points inbetween refinement regions."
283 << endl;
284
285 return tdisplacement;
286}
287
288
289// Calculate displacement as average of patch points.
290Foam::tmp<Foam::pointField> Foam::snappySnapDriver::smoothPatchDisplacement
291(
292 const motionSmoother& meshMover,
293 const List<labelPair>& baffles
294)
295{
296 const indirectPrimitivePatch& pp = meshMover.patch();
297
298 // Calculate geometrically non-manifold points on the patch to be moved.
299 bitSet nonManifoldPoint(pp.nPoints());
300 label nNonManifoldPoints = getCollocatedPoints
301 (
302 SMALL,
303 pp.localPoints(),
304 nonManifoldPoint
305 );
306 Info<< "Found " << nNonManifoldPoints << " non-manifold point(s)."
307 << endl;
308
309
310 // Average points
311 // ~~~~~~~~~~~~~~
312
313 // We determine three points:
314 // - average of (centres of) connected patch faces
315 // - average of (centres of) connected internal mesh faces
316 // - as fallback: centre of any connected cell
317 // so we can do something moderately sensible for non/manifold points.
318
319 // Note: the averages are calculated properly parallel. This is
320 // necessary to get the points shared by processors correct.
321
322
323 const labelListList& pointFaces = pp.pointFaces();
324 const labelList& meshPoints = pp.meshPoints();
325 const pointField& points = pp.points();
326 const polyMesh& mesh = meshMover.mesh();
327
328 // Get labels of faces to count (master of coupled faces and baffle pairs)
330
331 {
332 forAll(baffles, i)
333 {
334 label f0 = baffles[i].first();
335 label f1 = baffles[i].second();
336
337 if (isMasterFace.test(f0))
338 {
339 // Make f1 a slave
340 isMasterFace.unset(f1);
341 }
342 else if (isMasterFace.test(f1))
343 {
344 isMasterFace.unset(f0);
345 }
346 else
347 {
349 << "Both sides of baffle consisting of faces " << f0
350 << " and " << f1 << " are already slave faces."
351 << abort(FatalError);
352 }
353 }
354 }
355
356
357 // Get average position of boundary face centres
358 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
359
360 Field<weightedPosition> avgBoundary
361 (
362 pointFaces.size(),
364 );
365 {
366 forAll(pointFaces, patchPointi)
367 {
368 const labelList& pFaces = pointFaces[patchPointi];
369
370 forAll(pFaces, pfi)
371 {
372 label facei = pFaces[pfi];
373
374 if (isMasterFace.test(pp.addressing()[facei]))
375 {
376 avgBoundary[patchPointi].first() += 1.0;
377 avgBoundary[patchPointi].second() +=
378 pp[facei].centre(points);
379 }
380 }
381 }
382
383 // Add coupled contributions
385
386 // Normalise
387 forAll(avgBoundary, i)
388 {
389 // Note: what if there is no master boundary face?
390 if (mag(avgBoundary[i].first()) > VSMALL)
391 {
392 avgBoundary[i].second() /= avgBoundary[i].first();
393 }
394 }
395 }
396
397
398 // Get average position of internal face centres
399 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
400
401 Field<weightedPosition> avgInternal;
402 {
404 (
405 mesh.nPoints(),
407 );
408
409 // Note: no use of pointFaces
410 const faceList& faces = mesh.faces();
411
412 for (label facei = 0; facei < mesh.nInternalFaces(); facei++)
413 {
414 const face& f = faces[facei];
415 const point& fc = mesh.faceCentres()[facei];
416
417 forAll(f, fp)
418 {
419 weightedPosition& wp = globalSum[f[fp]];
420 wp.first() += 1.0;
421 wp.second() += fc;
422 }
423 }
424
425 // Count coupled faces as internal ones (but only once)
427
428 forAll(patches, patchi)
429 {
430 if
431 (
432 patches[patchi].coupled()
433 && refCast<const coupledPolyPatch>(patches[patchi]).owner()
434 )
435 {
436 const coupledPolyPatch& pp =
438
439 const vectorField::subField faceCentres = pp.faceCentres();
440
441 forAll(pp, i)
442 {
443 const face& f = pp[i];
444 const point& fc = faceCentres[i];
445
446 forAll(f, fp)
447 {
448 weightedPosition& wp = globalSum[f[fp]];
449 wp.first() += 1.0;
450 wp.second() += fc;
451 }
452 }
453 }
454 }
455
456 // Add coupled contributions
458
459 avgInternal.setSize(meshPoints.size());
460
461 forAll(avgInternal, patchPointi)
462 {
463 label meshPointi = meshPoints[patchPointi];
464 const weightedPosition& wp = globalSum[meshPointi];
465
466 avgInternal[patchPointi].first() = wp.first();
467 if (mag(wp.first()) < VSMALL)
468 {
469 // Set to zero?
470 avgInternal[patchPointi].second() = wp.second();
471 }
472 else
473 {
474 avgInternal[patchPointi].second() = wp.second()/wp.first();
475 }
476 }
477 }
478
479
480 // Precalculate any cell using mesh point (replacement of pointCells()[])
481 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
482
483 labelList anyCell(mesh.nPoints(), -1);
484 forAll(mesh.faceOwner(), facei)
485 {
486 label own = mesh.faceOwner()[facei];
487 const face& f = mesh.faces()[facei];
488
489 forAll(f, fp)
490 {
491 anyCell[f[fp]] = own;
492 }
493 }
494
495
496 // Displacement to calculate.
497 auto tpatchDisp = tmp<pointField>::New(meshPoints.size(), Zero);
498 auto& patchDisp = tpatchDisp.ref();
499
500 forAll(pointFaces, i)
501 {
502 label meshPointi = meshPoints[i];
503 const point& currentPos = pp.points()[meshPointi];
504
505 // Now we have the two average points and their counts:
506 // avgBoundary and avgInternal
507 // Do some blending between the two.
508 // Note: the following section has some reasoning behind it but the
509 // blending factors can be experimented with.
510
511 const weightedPosition& internal = avgInternal[i];
512 const weightedPosition& boundary = avgBoundary[i];
513
514 point newPos;
515
516 if (!nonManifoldPoint.test(i))
517 {
518 // Points that are manifold. Weight the internal and boundary
519 // by their number of faces and blend with
520 scalar internalBlend = 0.1;
521 scalar blend = 0.1;
522
523 point avgPos =
524 (
525 internalBlend*internal.first()*internal.second()
526 +(1-internalBlend)*boundary.first()*boundary.second()
527 )
528 / (
529 internalBlend*internal.first()
530 +(1-internalBlend)*boundary.first()
531 );
532
533 newPos = (1-blend)*avgPos + blend*currentPos;
534 }
535 else if (internal.first() == 0)
536 {
537 // Non-manifold without internal faces. Use any connected cell
538 // as internal point instead. Use precalculated any cell to avoid
539 // e.g. pointCells()[meshPointi][0]
540
541 const point& cc = mesh.cellCentres()[anyCell[meshPointi]];
542
543 scalar cellCBlend = 0.8;
544 scalar blend = 0.1;
545
546 point avgPos = (1-cellCBlend)*boundary.second() + cellCBlend*cc;
547
548 newPos = (1-blend)*avgPos + blend*currentPos;
549 }
550 else
551 {
552 // Non-manifold point with internal faces connected to them
553 scalar internalBlend = 0.9;
554 scalar blend = 0.1;
555
556 point avgPos =
557 internalBlend*internal.second()
558 + (1-internalBlend)*boundary.second();
559
560 newPos = (1-blend)*avgPos + blend*currentPos;
561 }
562
563 patchDisp[i] = newPos - currentPos;
564 }
565
566 return tpatchDisp;
567}
568//XXXXXXX
569//Foam::tmp<Foam::pointField> Foam::snappySnapDriver::avg
570//(
571// const indirectPrimitivePatch& pp,
572// const pointField& localPoints
573//)
574//{
575// const labelListList& pointEdges = pp.pointEdges();
576// const edgeList& edges = pp.edges();
577//
578// auto tavg = tmp<pointField>::New(pointEdges.size(), Zero);
579// auto& avg = tavg.ref();
580//
581// forAll(pointEdges, verti)
582// {
583// vector& avgPos = avg[verti];
584//
585// const labelList& pEdges = pointEdges[verti];
586//
587// forAll(pEdges, myEdgei)
588// {
589// const edge& e = edges[pEdges[myEdgei]];
590//
591// label otherVerti = e.otherVertex(verti);
592//
593// avgPos += localPoints[otherVerti];
594// }
595//
596// avgPos /= pEdges.size();
597// }
598// return tavg;
599//}
600//Foam::tmp<Foam::pointField>
601//Foam::snappySnapDriver::smoothLambdaMuPatchDisplacement
602//(
603// const motionSmoother& meshMover,
604// const List<labelPair>& baffles
605//)
606//{
607// const indirectPrimitivePatch& pp = meshMover.patch();
608// pointField newLocalPoints(pp.localPoints());
609//
610// const label iters = 90;
611// const scalar lambda = 0.33;
612// const scalar mu = 0.34;
613//
614// for (label iter = 0; iter < iters; iter++)
615// {
616// // Lambda
617// newLocalPoints =
618// (1 - lambda)*newLocalPoints
619// + lambda*avg(pp, newLocalPoints);
620//
621// // Mu
622// newLocalPoints =
623// (1 + mu)*newLocalPoints
624// - mu*avg(pp, newLocalPoints);
625// }
626// return newLocalPoints-pp.localPoints();
627//}
628//XXXXXXX
629
630
631Foam::tmp<Foam::scalarField> Foam::snappySnapDriver::edgePatchDist
632(
633 const pointMesh& pMesh,
635)
636{
637 const polyMesh& mesh = pMesh();
638
639 // Set initial changed points to all the patch points
640 List<pointEdgePoint> wallInfo(pp.nPoints());
641
642 forAll(pp.localPoints(), ppi)
643 {
644 wallInfo[ppi] = pointEdgePoint(pp.localPoints()[ppi], 0.0);
645 }
646
647 // Current info on points
648 List<pointEdgePoint> allPointInfo(mesh.nPoints());
649
650 // Current info on edges
651 List<pointEdgePoint> allEdgeInfo(mesh.nEdges());
652
654 (
655 mesh,
656 pp.meshPoints(),
657 wallInfo,
658
659 allPointInfo,
660 allEdgeInfo,
661 mesh.globalData().nTotalPoints() // max iterations
662 );
663
664 // Copy edge values into scalarField
665 auto tedgeDist = tmp<scalarField>::New(mesh.nEdges());
666 auto& edgeDist = tedgeDist.ref();
667
668 forAll(allEdgeInfo, edgei)
669 {
670 edgeDist[edgei] = Foam::sqrt(allEdgeInfo[edgei].distSqr());
671 }
672
673 return tedgeDist;
674}
675
676
677void Foam::snappySnapDriver::dumpMove
678(
679 const fileName& fName,
680 const pointField& meshPts,
681 const pointField& surfPts
682)
683{
684 // Dump direction of growth into file
685 Info<< "Dumping move direction to " << fName << endl;
686
687 OFstream nearestStream(fName);
688
689 label verti = 0;
690
691 forAll(meshPts, pti)
692 {
693 meshTools::writeOBJ(nearestStream, meshPts[pti]);
694 verti++;
695
696 meshTools::writeOBJ(nearestStream, surfPts[pti]);
697 verti++;
698
699 nearestStream<< "l " << verti-1 << ' ' << verti << nl;
700 }
701}
702
703
704// Check whether all displacement vectors point outwards of patch. Return true
705// if so.
706bool Foam::snappySnapDriver::outwardsDisplacement
707(
709 const vectorField& patchDisp
710)
711{
712 const vectorField& faceNormals = pp.faceNormals();
713 const labelListList& pointFaces = pp.pointFaces();
714
715 forAll(pointFaces, pointi)
716 {
717 const labelList& pFaces = pointFaces[pointi];
718
719 vector disp(patchDisp[pointi]);
720
721 scalar magDisp = mag(disp);
722
723 if (magDisp > SMALL)
724 {
725 disp /= magDisp;
726
727 bool outwards = meshTools::visNormal(disp, faceNormals, pFaces);
728
729 if (!outwards)
730 {
731 Warning<< "Displacement " << patchDisp[pointi]
732 << " at mesh point " << pp.meshPoints()[pointi]
733 << " coord " << pp.points()[pp.meshPoints()[pointi]]
734 << " points through the surrounding patch faces" << endl;
735 return false;
736 }
737 }
738 else
739 {
740 //? Displacement small but in wrong direction. Would probably be ok.
741 }
742 }
743 return true;
744}
745
746
747void Foam::snappySnapDriver::freezeExposedPoints
748(
749 const meshRefinement& meshRefiner,
750 const word& fzName, // faceZone name
751 const word& pzName, // pointZone name
752 const indirectPrimitivePatch& outside,
753 vectorField& outsideDisp
754)
755{
756 const fvMesh& mesh = meshRefiner.mesh();
757 const pointZoneMesh& pointZones = mesh.pointZones();
758
759 bitSet isFrozenPoint(mesh.nPoints());
760
761 // Add frozen points
762 const label pointZonei = pointZones.findZoneID(pzName);
763 if (pointZonei != -1)
764 {
765 isFrozenPoint.set(pointZones[pointZonei]);
766 }
767
768 // Add (inside) points of frozen faces
769 const faceZoneMesh& faceZones = mesh.faceZones();
770 const label faceZonei = faceZones.findZoneID(fzName);
771 if (faceZonei != -1)
772 {
774 (
775 UIndirectList<face>(mesh.faces(), faceZones[faceZonei]),
776 mesh.points()
777 );
778
779 // Count number of faces per edge
780 const labelList nEdgeFaces(meshRefiner.countEdgeFaces(pp));
781
782 // Freeze all internal points
783 forAll(nEdgeFaces, edgei)
784 {
785 if (nEdgeFaces[edgei] != 1)
786 {
787 const edge& e = pp.edges()[edgei];
788 isFrozenPoint.set(pp.meshPoints()[e[0]]);
789 isFrozenPoint.set(pp.meshPoints()[e[1]]);
790 }
791 }
792 }
793
795 (
796 mesh,
797 isFrozenPoint,
799 0u
800 );
801
802 for (const label pointi : isFrozenPoint)
803 {
804 const auto iter = outside.meshPointMap().find(pointi);
805 if (iter.good())
806 {
807 outsideDisp[iter.val()] = Zero;
808 }
810}
811
812
813// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
814
815Foam::snappySnapDriver::snappySnapDriver
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,
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
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,
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.good())
1025 {
1026 const label pointi = iter.val();
1027 pointOnZone[pointi] = true;
1029 }
1030 }
1031}
1032
1033
1034Foam::tmp<Foam::pointField> Foam::snappySnapDriver::avgCellCentres
1035(
1036 const fvMesh& mesh,
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
1063
1064 auto tavgBoundary = tmp<pointField>::New(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// auto tedgeLen = tmp<scalarField>::New(pp.nPoints());
1081// auto& edgeLen = tedgeLen.ref();
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,
1123 const pointField& localPoints,
1124 const pointField& nearestPoint,
1125 const vectorField& nearestNormal,
1126
1127 vectorField& disp
1128) const
1129{
1130 Info<< "Detecting near surfaces ..." << endl;
1131
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.writeLine(start[i], hit1[i].point());
1304 // }
1305 // if (hit2[i].hit())
1306 // {
1307 // str.writeLine(start[i], hit2[i].point());
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].point();
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].point();
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.writeLine(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].point(),
1506 // normal1[pointi]
1507 // )
1508 // )
1509 // {
1510 // disp[pointi] = hit1[pointi].point()-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].point(),
1524 // normal2[pointi]
1525 // )
1526 // )
1527 // {
1528 // disp[pointi] = hit2[pointi].point()-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].point(),
1541 normal1[pointi],
1542 hit2[pointi].point(),
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 gapStr().writeLine(pt, hit2[pointi].point());
1553 }
1554
1555 // Choose hit2 : nearest to end point (so inside the domain)
1556 disp[pointi] = hit2[pointi].point()-pt;
1557 override = true;
1558 }
1559 }
1560
1561 if (override && isPatchMasterPoint[pointi])
1562 {
1563 nOverride++;
1564 }
1565 }
1566 }
1567
1568
1569 // 2. All points on zones to their respective surface
1570 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1571
1572 {
1573 // Surfaces with zone information
1574 const PtrList<surfaceZonesInfo>& surfZones = surfaces.surfZones();
1575
1576 const labelList zonedSurfaces = surfaceZonesInfo::getNamedSurfaces
1577 (
1578 surfZones
1579 );
1580
1581 forAll(zonedSurfaces, i)
1582 {
1583 label zoneSurfi = zonedSurfaces[i];
1584 const labelList surfacesToTest(1, zoneSurfi);
1585
1586 const wordList& faceZoneNames =
1587 surfZones[zoneSurfi].faceZoneNames();
1588 forAll(faceZoneNames, namei)
1589 {
1590 const word& faceZoneName = faceZoneNames[namei];
1591
1592 // Get indices of points both on faceZone and on pp.
1593 bitSet pointOnZone(pp.nPoints());
1594 getZoneSurfacePoints
1595 (
1596 mesh,
1597 pp,
1598 faceZoneName,
1599 pointOnZone
1600 );
1601 const labelList zonePointIndices(pointOnZone.toc());
1602
1603 // Do intersection test
1604 labelList surface1;
1606 labelList region1;
1607 vectorField normal1;
1608
1609 labelList surface2;
1611 labelList region2;
1612 vectorField normal2;
1614 (
1615 surfacesToTest,
1616 pointField(start, zonePointIndices),
1617 pointField(end, zonePointIndices),
1618
1619 surface1,
1620 hit1,
1621 region1,
1622 normal1,
1623
1624 surface2,
1625 hit2,
1626 region2,
1627 normal2
1628 );
1629
1630
1631 forAll(hit1, i)
1632 {
1633 label pointi = zonePointIndices[i];
1634
1635 // Current location
1636 const point& pt = localPoints[pointi];
1637
1638 bool override = false;
1639
1640 //if (hit1[i].hit())
1641 //{
1642 // if
1643 // (
1644 // meshRefiner_.isGap
1645 // (
1646 // planarCos,
1647 // nearestPoint[pointi],
1648 // nearestNormal[pointi],
1649 // hit1[i].point(),
1650 // normal1[i]
1651 // )
1652 // )
1653 // {
1654 // disp[pointi] = hit1[i].point()-pt;
1655 // override = true;
1656 // }
1657 //}
1658 //if (hit2[i].hit())
1659 //{
1660 // if
1661 // (
1662 // meshRefiner_.isGap
1663 // (
1664 // planarCos,
1665 // nearestPoint[pointi],
1666 // nearestNormal[pointi],
1667 // hit2[i].point(),
1668 // normal2[i]
1669 // )
1670 // )
1671 // {
1672 // disp[pointi] = hit2[i].point()-pt;
1673 // override = true;
1674 // }
1675 //}
1676
1677 if (hit1[i].hit() && hit2[i].hit())
1678 {
1679 if
1680 (
1681 meshRefiner_.isGap
1682 (
1683 planarCos,
1684 hit1[i].point(),
1685 normal1[i],
1686 hit2[i].point(),
1687 normal2[i]
1688 )
1689 )
1690 {
1691 if (gapStr)
1692 {
1693 gapStr().writeLine(pt, hit2[i].point());
1694 }
1695
1696 disp[pointi] = hit2[i].point()-pt;
1697 override = true;
1698 }
1699 }
1700
1701 if (override && isPatchMasterPoint[pointi])
1702 {
1703 nOverride++;
1704 }
1705 }
1706 }
1707 }
1708 }
1709
1710 Info<< "Overriding nearest with intersection of close gaps at "
1711 << returnReduce(nOverride, sumOp<label>())
1712 << " out of " << returnReduce(pp.nPoints(), sumOp<label>())
1713 << " points." << endl;
1714}
1715
1716
1717void Foam::snappySnapDriver::calcNearestSurface
1718(
1719 const refinementSurfaces& surfaces,
1720
1721 const labelList& surfacesToTest,
1722 const labelListList& regionsToTest,
1723
1724 const pointField& localPoints,
1725 const labelList& zonePointIndices,
1726
1727 scalarField& minSnapDist,
1728 labelList& snapSurf,
1729 vectorField& patchDisp,
1730
1731 // Optional: nearest point, normal
1732 pointField& nearestPoint,
1733 vectorField& nearestNormal
1734)
1735{
1736 // Find nearest for points both on faceZone and pp.
1737 List<pointIndexHit> hitInfo;
1738 labelList hitSurface;
1739
1740 if (nearestNormal.size() == localPoints.size())
1741 {
1742 labelList hitRegion;
1743 vectorField hitNormal;
1744 surfaces.findNearestRegion
1745 (
1746 surfacesToTest,
1747 regionsToTest,
1748
1749 pointField(localPoints, zonePointIndices),
1750 sqr(scalarField(minSnapDist, zonePointIndices)),
1751
1752 hitSurface,
1753 hitInfo,
1754 hitRegion,
1755 hitNormal
1756 );
1757
1758 forAll(hitInfo, i)
1759 {
1760 if (hitInfo[i].hit())
1761 {
1762 label pointi = zonePointIndices[i];
1763 nearestPoint[pointi] = hitInfo[i].point();
1764 nearestNormal[pointi] = hitNormal[i];
1765 }
1766 }
1767 }
1768 else
1769 {
1770 surfaces.findNearest
1771 (
1772 surfacesToTest,
1773 regionsToTest,
1774
1775 pointField(localPoints, zonePointIndices),
1776 sqr(scalarField(minSnapDist, zonePointIndices)),
1777
1778 hitSurface,
1779 hitInfo
1780 );
1781 }
1782
1783 forAll(hitInfo, i)
1784 {
1785 if (hitInfo[i].hit())
1786 {
1787 label pointi = zonePointIndices[i];
1788
1789 patchDisp[pointi] = hitInfo[i].point() - localPoints[pointi];
1790 minSnapDist[pointi] = mag(patchDisp[pointi]);
1791 snapSurf[pointi] = hitSurface[i];
1792 }
1793 }
1794}
1795
1796
1797Foam::vectorField Foam::snappySnapDriver::calcNearestSurface
1798(
1799 const bool strictRegionSnap,
1800 const meshRefinement& meshRefiner,
1801 const labelList& globalToMasterPatch,
1802 const labelList& globalToSlavePatch,
1804 const pointField& localPoints,
1805 const scalarField& snapDist,
1806 pointField& nearestPoint,
1807 vectorField& nearestNormal
1808)
1809{
1810 Info<< "Calculating patchDisplacement as distance to nearest surface"
1811 << " point ..." << endl;
1812 if (strictRegionSnap)
1813 {
1814 Info<< " non-zone points : attract to local region on surface only"
1815 << nl
1816 << " zone points : attract to local region on surface only"
1817 << nl
1818 << endl;
1819 }
1820 else
1821 {
1822 Info<< " non-zone points :"
1823 << " attract to nearest of all non-zone surfaces"
1824 << nl
1825 << " zone points : attract to zone surface only" << nl
1826 << endl;
1827 }
1828
1829
1830 const refinementSurfaces& surfaces = meshRefiner.surfaces();
1831 const fvMesh& mesh = meshRefiner.mesh();
1832
1833 // Displacement per patch point
1834 vectorField patchDisp(localPoints.size(), Zero);
1835
1836 if (returnReduceOr(localPoints.size()))
1837 {
1838 // Current surface snapped to. Used to check whether points have been
1839 // snapped at all
1840 labelList snapSurf(localPoints.size(), -1);
1841
1842 // Current best snap distance (since point might be on multiple
1843 // regions)
1844 scalarField minSnapDist(snapDist);
1845
1846
1847 if (strictRegionSnap)
1848 {
1849 // Attract patch points to same region only
1850
1851 forAll(surfaces.surfaces(), surfi)
1852 {
1853 label geomi = surfaces.surfaces()[surfi];
1854 label nRegions = surfaces.geometry()[geomi].regions().size();
1855
1856 const labelList surfacesToTest(1, surfi);
1857
1858 for (label regioni = 0; regioni < nRegions; regioni++)
1859 {
1860 label globali = surfaces.globalRegion(surfi, regioni);
1861 label masterPatchi = globalToMasterPatch[globali];
1862
1863 // Get indices of points both on patch and on pp
1864 labelList zonePointIndices
1865 (
1866 getFacePoints
1867 (
1868 pp,
1869 mesh.boundaryMesh()[masterPatchi]
1870 )
1871 );
1872
1873 calcNearestSurface
1874 (
1875 surfaces,
1876
1877 surfacesToTest,
1878 labelListList(1, labelList(1, regioni)), //regionsToTest
1879
1880 localPoints,
1881 zonePointIndices,
1882
1883 minSnapDist,
1884 snapSurf,
1885 patchDisp,
1886
1887 // Optional: nearest point, normal
1888 nearestPoint,
1889 nearestNormal
1890 );
1891
1892 if (globalToSlavePatch[globali] != masterPatchi)
1893 {
1894 label slavePatchi = globalToSlavePatch[globali];
1895
1896 // Get indices of points both on patch and on pp
1897 labelList zonePointIndices
1898 (
1899 getFacePoints
1900 (
1901 pp,
1902 mesh.boundaryMesh()[slavePatchi]
1903 )
1904 );
1905
1906 calcNearestSurface
1907 (
1908 surfaces,
1909
1910 surfacesToTest,
1911 labelListList(1, labelList(1, regioni)),
1912
1913 localPoints,
1914 zonePointIndices,
1915
1916 minSnapDist,
1917 snapSurf,
1918 patchDisp,
1919
1920 // Optional: nearest point, normal
1921 nearestPoint,
1922 nearestNormal
1923 );
1924 }
1925 }
1926 }
1927 }
1928 else
1929 {
1930 // Divide surfaces into zoned and unzoned
1931 const labelList unzonedSurfaces =
1933 (
1934 meshRefiner.surfaces().surfZones()
1935 );
1936
1937
1938 // 1. All points to non-interface surfaces
1939 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1940
1941 List<pointIndexHit> hitInfo;
1942 labelList hitSurface;
1943
1944 if (nearestNormal.size() == localPoints.size())
1945 {
1946 labelList hitRegion;
1947 vectorField hitNormal;
1948 surfaces.findNearestRegion
1949 (
1950 unzonedSurfaces,
1951 localPoints,
1952 sqr(snapDist),
1953 hitSurface,
1954 hitInfo,
1955 hitRegion,
1956 hitNormal
1957 );
1958
1959 forAll(hitInfo, pointi)
1960 {
1961 if (hitInfo[pointi].hit())
1962 {
1963 nearestPoint[pointi] = hitInfo[pointi].point();
1964 nearestNormal[pointi] = hitNormal[pointi];
1965 }
1966 }
1967 }
1968 else
1969 {
1970 surfaces.findNearest
1971 (
1972 unzonedSurfaces,
1973 localPoints,
1974 sqr(snapDist), // sqr of attract distance
1975 hitSurface,
1976 hitInfo
1977 );
1978 }
1979
1980 forAll(hitInfo, pointi)
1981 {
1982 if (hitInfo[pointi].hit())
1983 {
1984 patchDisp[pointi] =
1985 hitInfo[pointi].point()
1986 - localPoints[pointi];
1987
1988 snapSurf[pointi] = hitSurface[pointi];
1989 }
1990 }
1991
1992
1993 const labelList zonedSurfaces = surfaceZonesInfo::getNamedSurfaces
1994 (
1995 meshRefiner.surfaces().surfZones()
1996 );
1997
1998
1999 // 2. All points on zones to their respective surface
2000 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2001 // (ignoring faceZone subdivision)
2002
2003 // Surfaces with zone information
2004 const PtrList<surfaceZonesInfo>& surfZones = surfaces.surfZones();
2005
2006 forAll(zonedSurfaces, i)
2007 {
2008 label surfi = zonedSurfaces[i];
2009 const labelList surfacesToTest(1, surfi);
2010 const label geomi = surfaces.surfaces()[surfi];
2011 const label nRegions =
2012 surfaces.geometry()[geomi].regions().size();
2013
2014 const wordList& faceZoneNames =
2015 surfZones[surfi].faceZoneNames();
2016
2017 // Get indices of points both on any faceZone and on pp.
2018 bitSet pointOnZone(pp.nPoints());
2019 forAll(faceZoneNames, locali)
2020 {
2021 getZoneSurfacePoints
2022 (
2023 mesh,
2024 pp,
2025 faceZoneNames[locali],
2026 pointOnZone
2027 );
2028 }
2029 const labelList zonePointIndices(pointOnZone.toc());
2030
2031 calcNearestSurface
2032 (
2033 surfaces,
2034
2035 surfacesToTest,
2036 labelListList(1, identity(nRegions)),
2037
2038 localPoints,
2039 zonePointIndices,
2040
2041 minSnapDist,
2042 snapSurf,
2043 patchDisp,
2044
2045 // Optional: nearest point, normal
2046 nearestPoint,
2047 nearestNormal
2048 );
2049 }
2050 }
2051
2052
2053 // Check if all points are being snapped
2054 forAll(snapSurf, pointi)
2055 {
2056 if (snapSurf[pointi] == -1)
2057 {
2058 static label nWarn = 0;
2059
2060 if (nWarn < 100)
2061 {
2063 << "For point:" << pointi
2064 << " coordinate:" << localPoints[pointi]
2065 << " did not find any surface within:"
2066 << minSnapDist[pointi] << " metre." << endl;
2067 nWarn++;
2068 if (nWarn == 100)
2069 {
2071 << "Reached warning limit " << nWarn
2072 << ". Suppressing further warnings." << endl;
2073 }
2074 }
2075 }
2076 }
2077
2078 {
2079 const bitSet isPatchMasterPoint
2080 (
2082 (
2083 mesh,
2084 pp.meshPoints()
2085 )
2086 );
2087
2088 scalarField magDisp(mag(patchDisp));
2089
2090 auto limits = gMinMax(magDisp);
2091
2092 Info<< "Wanted displacement : average:"
2093 << meshRefinement::gAverage(isPatchMasterPoint, magDisp)
2094 << " min:" << limits.min()
2095 << " max:" << limits.max() << 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)
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
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 }
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
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 (
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
2411 const polyBoundaryMesh& patches = mesh.boundaryMesh();
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,
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 {
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];
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 if (meshRefiner_.meshType() == meshRefinement::CASTELLATEDBUFFERLAYER2)
2585 {
2586 // Buffer-layer replacement for this routine
2587
2588 doSnapBufferLayers
2589 (
2590 snapDict,
2591 motionDict,
2592 mergeType,
2593 featureCos,
2594 planarAngle,
2595 snapParams
2596 );
2597 return;
2598 }
2599
2600 addProfiling(snap, "snappyHexMesh::snap");
2601 fvMesh& mesh = meshRefiner_.mesh();
2602
2603 Info<< nl
2604 << "Morphing phase" << nl
2605 << "--------------" << nl
2606 << endl;
2607
2608 // faceZone handling
2609 // ~~~~~~~~~~~~~~~~~
2610 //
2611 // We convert all faceZones into baffles during snapping so we can use
2612 // a standard mesh motion (except for the mesh checking which for baffles
2613 // created from internal faces should check across the baffles). The state
2614 // is stored in two variables:
2615 // baffles : pairs of boundary faces
2616 // duplicateFace : from mesh face to its baffle colleague (or -1 for
2617 // normal faces)
2618 // There are three types of faceZones according to the faceType property:
2619 //
2620 // internal
2621 // --------
2622 // - baffles: need to be checked across
2623 // - duplicateFace: from face to duplicate face. Contains
2624 // all faces on faceZone to prevents merging patch faces.
2625 //
2626 // baffle
2627 // ------
2628 // - baffles: no need to be checked across
2629 // - duplicateFace: contains all faces on faceZone to prevent
2630 // merging patch faces.
2631 //
2632 // boundary
2633 // --------
2634 // - baffles: no need to be checked across. Also points get duplicated
2635 // so will no longer be baffles
2636 // - duplicateFace: contains no faces on faceZone since both sides can
2637 // merge faces independently.
2638
2639
2640
2641 // Get labels of patches where optional buffer layers are added
2642 DynamicList<label> bufPatchIDs;
2643 if (meshRefiner_.meshType() == meshRefinement::CASTELLATEDBUFFERLAYER)
2644 {
2645 bufPatchIDs.setCapacity(globalToMasterPatch_.size());
2646
2647 const auto& addLayers =
2648 meshRefiner_.surfaces().addBufferLayers();
2649
2650 // Normal patches
2651 forAll(globalToMasterPatch_, globalRegioni)
2652 {
2653 if (addLayers[globalRegioni])
2654 {
2655 bufPatchIDs.push_uniq(globalToMasterPatch_[globalRegioni]);
2656 bufPatchIDs.push_uniq(globalToSlavePatch_[globalRegioni]);
2657 }
2658 }
2659
2660 // Temporary patches from faceZones
2661 for (const auto& fz : mesh.faceZones())
2662 {
2663 label mpI, spI;
2665 if (meshRefiner_.getFaceZoneInfo(fz.name(), mpI, spI, type))
2666 {
2667 bufPatchIDs.push_uniq(mpI);
2668 bufPatchIDs.push_uniq(spI);
2669 }
2670 }
2671 }
2672
2673
2674
2675 // faceZones of type internal
2676 const labelList internalFaceZones
2677 (
2678 meshRefiner_.getZones
2679 (
2681 (
2682 1,
2684 )
2685 )
2686 );
2687
2688 // Create baffles (pairs of faces that share the same points)
2689 // Baffles stored as owner and neighbour face that have been created.
2690 List<labelPair> baffles;
2691 {
2692 labelList originatingFaceZone;
2693 meshRefiner_.createZoneBaffles
2694 (
2696 baffles,
2697 originatingFaceZone
2698 );
2699 }
2700
2701 // Duplicate points on faceZones of type boundary. Renumber baffles
2702 // (probably not necessary - faceIDs should not change)
2703 {
2704 autoPtr<mapPolyMesh> map = meshRefiner_.dupNonManifoldBoundaryPoints();
2705 if (map)
2706 {
2707 const labelList& reverseFaceMap = map->reverseFaceMap();
2708 forAll(baffles, i)
2709 {
2710 label f0 = reverseFaceMap[baffles[i].first()];
2711 label f1 = reverseFaceMap[baffles[i].second()];
2712 baffles[i] = labelPair(f0, f1);
2713 }
2714 }
2715 }
2716
2717
2718 bool doFeatures = false;
2719 label nFeatIter = 1;
2720 if (snapParams.nFeatureSnap() > 0)
2721 {
2722 doFeatures = true;
2723
2724 if (!dryRun_)
2725 {
2726 nFeatIter = snapParams.nFeatureSnap();
2727 }
2728
2729 Info<< "Snapping to features in " << nFeatIter
2730 << " iterations ..." << endl;
2731 }
2732
2733
2734 bool meshOk = false;
2735
2736
2737 // Get the labels of added patches.
2738 const labelList adaptPatchIDs(meshRefiner_.meshedPatches());
2739
2740
2741
2742 {
2744 (
2746 (
2747 mesh,
2748 adaptPatchIDs
2749 )
2750 );
2751
2752 // Distance to attract to nearest feature on surface
2753 scalarField snapDist(calcSnapDistance(mesh, snapParams, ppPtr()));
2754
2755
2756 // Construct iterative mesh mover.
2757 Info<< "Constructing mesh displacer ..." << endl;
2758 Info<< "Using mesh parameters " << motionDict << nl << endl;
2759
2760 autoPtr<motionSmoother> meshMoverPtr
2761 (
2762 new motionSmoother
2763 (
2764 mesh,
2765 ppPtr(),
2766 adaptPatchIDs,
2768 (
2770 adaptPatchIDs
2771 ),
2772 motionDict,
2773 dryRun_
2774 )
2775 );
2776
2777
2778 // Check initial mesh
2779 Info<< "Checking initial mesh ..." << endl;
2780 labelHashSet wrongFaces(mesh.nFaces()/100);
2781 motionSmoother::checkMesh(false, mesh, motionDict, wrongFaces, dryRun_);
2782 const label nInitErrors = returnReduce
2783 (
2784 wrongFaces.size(),
2785 sumOp<label>()
2786 );
2787
2788 Info<< "Detected " << nInitErrors << " illegal faces"
2789 << " (concave, zero area or negative cell pyramid volume)"
2790 << endl;
2791
2792
2793 Info<< "Checked initial mesh in = "
2794 << mesh.time().cpuTimeIncrement() << " s\n" << nl << endl;
2795
2796 // Extract baffles across internal faceZones (for checking mesh quality
2797 // across)
2798 labelPairList internalBaffles
2799 (
2800 meshRefiner_.subsetBaffles
2801 (
2802 mesh,
2803 internalFaceZones,
2805 )
2806 );
2807
2808
2809
2810 // Pre-smooth patch vertices (so before determining nearest)
2811 preSmoothPatch
2812 (
2813 meshRefiner_,
2814 snapParams,
2815 nInitErrors,
2816 internalBaffles,
2817 meshMoverPtr()
2818 );
2819
2820 // Reset moving flag in case we do any topo changes
2821 mesh.moving(false);
2822
2823
2824 // Optionally add buffer layers
2825 if (meshRefiner_.meshType() == meshRefinement::CASTELLATEDBUFFERLAYER)
2826 {
2829 //autoPtr<mapPolyMesh> mapPtr = meshRefiner_.mergeZoneBaffles
2830 //(
2831 // true, // internal zones
2832 // false // baffle zones
2833 //);
2834 //
2835 //if (mapPtr)
2836 //{
2837 // if (debug & meshRefinement::MESH)
2838 // {
2839 // const_cast<Time&>(mesh.time())++;
2840 // Info<< "Writing baffle-merged mesh to time "
2841 // << meshRefiner_.timeName() << endl;
2842 // meshRefiner_.write
2843 // (
2844 // meshRefinement::debugType(debug),
2845 // meshRefinement::writeType
2846 // (
2847 // meshRefinement::writeLevel()
2848 // | meshRefinement::WRITEMESH
2849 // ),
2850 // meshRefiner_.timeName()
2851 // );
2852 // }
2853 //}
2854
2856 //Info<< "Adding buffer layers ..." << endl;
2857 //
2858 //{
2859 // // Remove references to pp
2860 // meshMoverPtr.clear();
2861 //
2862 // const labelList meshFaces
2863 // (
2864 // mesh.faceZones().selection
2865 // (
2866 // internalFaceZones
2867 // ).sortedToc()
2868 // );
2869 // ppPtr.reset
2870 // (
2871 // new indirectPrimitivePatch
2872 // (
2873 // IndirectList<face>
2874 // (
2875 // mesh.faces(),
2876 // meshFaces
2877 // ),
2878 // mesh.points()
2879 // )
2880 // );
2881 // const pointField thickness
2882 // (
2883 // wantedThickness(ppPtr(), 1e-1) //cellSizeFraction
2884 // * PatchTools::pointNormals(mesh, ppPtr())
2885 // );
2886 //
2887 // // Layer mesh modifier
2888 // // - use intrusion, not extrusion
2889 // addPatchCellLayer addLayer(mesh, true, false);
2890 //
2891 // // Do mesh changes : introduce point, faces, cells
2892 // autoPtr<mapPolyMesh> mapPtr = addBufferLayers
2893 // (
2894 // ppPtr(),
2895 // -thickness, //1e-3, //cellSizeFraction,
2896 // addLayer
2897 // );
2898 //
2899 // // Update numbering on baffles, pointToMaster. Note: uses
2900 // // geometric tolerance - could be avoided since now we have
2901 // // addPatchCellLayer still intact. However would like to avoid
2902 // // use of addPatchCellLayer altogether.
2903 // if (mapPtr)
2904 // {
2905 // // Invalidate extrusion (face numbering might have changed)
2906 // ppPtr.clear();
2907 //
2908 // labelList dummyPointToMaster;
2909 // snappyLayerDriver::mapFaceZonePoints
2910 // (
2911 // meshRefiner_,
2912 // mapPtr(),
2913 // internalBaffles,
2914 // dummyPointToMaster
2915 // );
2916 //
2917 // if (debug & meshRefinement::MESH)
2918 // {
2919 // const_cast<Time&>(mesh.time())++;
2920 // Info<< "Writing INTERNAL ZONE buffer layer mesh"
2921 // << " to time " << meshRefiner_.timeName() << endl;
2922 // meshRefiner_.write
2923 // (
2924 // meshRefinement::debugType(debug),
2925 // meshRefinement::writeType
2926 // (
2927 // meshRefinement::writeLevel()
2928 // | meshRefinement::WRITEMESH
2929 // ),
2930 // meshRefiner_.timeName()
2931 // );
2932 // }
2933 // }
2934 //}
2935
2937 //{
2938 // labelList originatingFaceZone;
2939 // meshRefiner_.createZoneBaffles
2940 // (
2941 // internalFaceZones,
2942 // baffles,
2943 // originatingFaceZone
2944 // );
2945 //}
2946
2947
2948 Info<< "Adding buffer layers ..." << endl;
2949
2950 // Remove references to pp
2951 meshMoverPtr.clear();
2952 // Invalidate extrusion (face numbering might have changed)
2953 ppPtr.clear();
2954
2955 // Note: all the way up to addBufferLayers can probably be replaced
2956 // once addPatchCellLayer can add to both sides of faceZone ...
2957
2958 // Duplicate points on faceZones that layers are added to
2959 labelList pointToMaster;
2960 {
2961 labelList numLayers(mesh.boundaryMesh().size(), 0);
2962 UIndirectList<label>(numLayers, bufPatchIDs) = 1;
2963
2964 autoPtr<mapPolyMesh> mapPtr =
2966 (
2967 meshRefiner_,
2968 bufPatchIDs, // patch indices
2969 numLayers, // num layers per patch
2970 baffles,
2971 pointToMaster
2972 );
2973 if (mapPtr)
2974 {
2975 // Update numbering on any baffles
2977 (
2978 mesh,
2979 mapPtr().faceMap(),
2980 internalBaffles
2981 );
2982
2983
2986 //if (debug & meshRefinement::MESH)
2987 //{
2988 // pointField newPoints(mesh.nPoints());
2989 // forAll(newPoints, pointi)
2990 // {
2991 // const auto& pCells = mesh.pointCells()[pointi];
2992 //
2993 // point avg(Zero);
2994 // for (const label celli : pCells)
2995 // {
2996 // avg += mesh.cellCentres()[celli];
2997 // }
2998 // avg /= pCells.size();
2999 //
3000 // newPoints[pointi] =
3001 // 0.5*mesh.points()[pointi]
3002 // + 0.5*avg;
3003 // }
3004 // mesh.movePoints(newPoints);
3005 //
3006 // const_cast<Time&>(mesh.time())++;
3007 // Info<< "Writing DEBUG-shrunk layer mesh to time "
3008 // << meshRefiner_.timeName() << endl;
3009 // meshRefiner_.write
3010 // (
3011 // meshRefinement::debugType(debug),
3012 // meshRefinement::writeType
3013 // (
3014 // meshRefinement::writeLevel()
3015 // | meshRefinement::WRITEMESH
3016 // ),
3017 // meshRefiner_.timeName()
3018 // );
3019 //}
3020 }
3021 }
3022
3023 //- Not needed: shrinking of mesh since now using intrusion ...
3024 // // Move mesh back with thickness. Two purposes:
3025 // // - avoid mapFaceZonePoints below merging points extraneously
3026 // // (does not use addPatchCellLayer structure; uses geometric
3027 // // tolerance)
3028 // // - see what is happening
3029 // {
3030 // pointField newPoints(mesh.points());
3031 // const auto& mp = ppPtr().meshPoints();
3032 // forAll(mp, i)
3033 // {
3034 // newPoints[mp[i]] -= thickness[i];
3035 // }
3036 // mesh.movePoints(newPoints);
3037 // ppPtr().movePoints(mesh.points());
3038 //
3039 // if (debug & meshRefinement::MESH)
3040 // {
3041 // const_cast<Time&>(mesh.time())++;
3042 // Info<< "Writing shrunk buffer layer mesh to time "
3043 // << meshRefiner_.timeName() << endl;
3044 // meshRefiner_.write
3045 // (
3046 // meshRefinement::debugType(debug),
3047 // meshRefinement::writeType
3048 // (
3049 // meshRefinement::writeLevel()
3050 // | meshRefinement::WRITEMESH
3051 // ),
3052 // meshRefiner_.timeName()
3053 // );
3054 // }
3055 // }
3056
3057
3058 {
3059 // Layer mesh modifier
3060 // - use intrusion, not extrusion
3061 addPatchCellLayer addLayer(mesh, true, false);
3062
3063 // Redo pp
3065 (
3066 meshRefinement::makePatch(mesh, bufPatchIDs)
3067 );
3068
3069 pointField thickness
3070 (
3071 wantedThickness(bufPatchPtr(), 1e-1) //cellSizeFraction
3072 * PatchTools::pointNormals(mesh, bufPatchPtr())
3073 );
3074
3075 // Make sure to adhere to constraints
3076 {
3077
3078 const pointMesh& pMesh = pointMesh::New(mesh);
3079 const labelList& mp = bufPatchPtr().meshPoints();
3080
3082 (
3084 (
3085 pMesh,
3086 adaptPatchIDs
3087 )
3088 );
3089 // Set internal field
3090 UIndirectList<point>(tdisp.ref(), mp) = thickness;
3091 // Take over onto boundary field. Needed since constraint
3092 // patch might be before the fixedValue patches so its
3093 // value gets overwritten
3094 for (auto& ppf : tdisp.ref().boundaryFieldRef())
3095 {
3096 ppf == ppf.patchInternalField();
3097 }
3098
3099 // Adhere to multi-point constraints
3100 const pointConstraints& pcs = pointConstraints::New(pMesh);
3101 pcs.constrainDisplacement(tdisp.ref(), false);
3102
3103 thickness = UIndirectList<point>(tdisp(), mp);
3104 }
3105
3106
3107
3108 // Print a bit
3109 {
3110 // See snappyLayerDriver::calculateLayerThickness.
3111 const auto& pbm = mesh.boundaryMesh();
3112 label maxLen = 0;
3113 for (const label patchi : bufPatchIDs)
3114 {
3115 maxLen = max(maxLen, label(pbm[patchi].name().size()));
3116 }
3117
3118 const int oldPrecision = Info.stream().precision();
3119
3120 Info<< nl
3121 << setf(ios_base::left) << setw(maxLen) << "patch"
3122 << setw(0) << " faces layers thickness[m]" << nl
3123 << setf(ios_base::left) << setw(maxLen) << "-----"
3124 << setw(0) << " ----- ------ ------------" << endl;
3125
3126 for (const label patchi : bufPatchIDs)
3127 {
3128 Info<< setf(ios_base::left) << setw(maxLen)
3129 << pbm[patchi].name() << setprecision(3)
3130 << " " << setw(8)
3131 << returnReduce(pbm[patchi].size(), sumOp<scalar>())
3132 << " " << setw(6) << 1
3133 << " " << setw(8) << gAverage(mag(thickness))
3134 << endl;
3135 }
3136 Info<< setprecision(oldPrecision) << endl;
3137 }
3138
3139
3140 // Do mesh changes : introduce point, faces, cells
3141 autoPtr<mapPolyMesh> mapPtr = addBufferLayers
3142 (
3143 bufPatchPtr(),
3144 -thickness, //1e-3, //cellSizeFraction,
3145 addLayer
3146 );
3147
3148 // Update numbering on baffles, pointToMaster. Note: uses
3149 // geometric tolerance - could be avoided since now we have
3150 // addPatchCellLayer still intact. However would like to avoid
3151 // use of addPatchCellLayer altogether.
3152 if (mapPtr)
3153 {
3154 // Invalidate extrusion (face numbering might have changed)
3155 ppPtr.clear();
3156
3158 (
3159 meshRefiner_,
3160 mapPtr(),
3161 internalBaffles,
3162 pointToMaster
3163 );
3164 }
3165 }
3166
3167
3168
3169 {
3170 // Merge duplicated points (this creates illegal cells -
3171 // hopefully they will be smoothed out)
3172 autoPtr<mapPolyMesh> mapPtr =
3173 meshRefiner_.mergePoints(pointToMaster);
3174 if (mapPtr)
3175 {
3176 // Invalidate extrusion (point numbering might have changed)
3177 ppPtr.clear();
3178
3179 // Update numbering on any baffles
3181 (
3182 mesh,
3183 mapPtr().faceMap(),
3184 internalBaffles
3185 );
3186 // Extract baffles across internal faceZones
3187 internalBaffles = meshRefinement::subsetBaffles
3188 (
3189 mesh,
3190 internalFaceZones,
3191 internalBaffles
3192 );
3193
3195 {
3196 const_cast<Time&>(mesh.time())++;
3197 Info<< "Writing merged points buffer layer mesh"
3198 << " to time " << meshRefiner_.timeName() << endl;
3199 meshRefiner_.write
3200 (
3203 (
3206 ),
3207 meshRefiner_.timeName()
3208 );
3209 }
3210 }
3211 }
3212
3213
3214 Info<< "Inflating buffer layers ..." << endl;
3215
3216 const pointMesh& pMesh = pointMesh::New(mesh);
3217
3218 {
3220 (
3221 makeMotionSolver
3222 (
3223 pMesh,
3224 snapDict,
3225 bufPatchIDs
3226 //patchConstraints
3227 )
3228 );
3229
3230 // Solve internal displacement
3231 tmp<pointField> tnewPoints(motionPtr->newPoints());
3232
3233 // Move points
3234 mesh.movePoints(tnewPoints);
3235
3236 // Reset moving flag to avoid problems with topo changes
3237 mesh.moving(false);
3238
3240 {
3241 const_cast<Time&>(mesh.time())++;
3242 Info<< "Writing smoothed buffer layer mesh to time "
3243 << meshRefiner_.timeName() << endl;
3244 meshRefiner_.write
3245 (
3248 (
3251 ),
3252 meshRefiner_.timeName()
3253 );
3254 }
3255 }
3256
3257 // Update mesh mover
3258 ppPtr = meshRefinement::makePatch(mesh, adaptPatchIDs);
3259
3260 // Update distance to attract to nearest feature on surface
3261 snapDist = calcSnapDistance(mesh, snapParams, ppPtr());
3262
3263 meshMoverPtr.reset
3264 (
3265 new motionSmoother
3266 (
3267 mesh,
3268 ppPtr(),
3269 adaptPatchIDs,
3271 (
3273 adaptPatchIDs
3274 ),
3275 motionDict,
3276 dryRun_
3277 )
3278 );
3279 }
3280
3281
3282 //- Only if in feature attraction mode:
3283 // Nearest feature
3284 vectorField patchAttraction;
3285 // Constraints at feature
3286 List<pointConstraint> patchConstraints;
3287
3288
3289 //- Any faces to split
3290 DynamicList<label> splitFaces;
3291 //- Indices in face to split across
3293 //- Patch for both sides of the face
3294 DynamicList<labelPair> splitPatches;
3295
3296
3297 for (label iter = 0; iter < nFeatIter; iter++)
3298 {
3299 Info<< nl
3300 << "Morph iteration " << iter << nl
3301 << "-----------------" << endl;
3302
3303 // Splitting iteration?
3304 bool doSplit = false;
3305 if
3306 (
3307 doFeatures
3308 && snapParams.nFaceSplitInterval() > 0
3309 && (
3310 (iter == nFeatIter-1)
3311 || (iter > 0 && (iter%snapParams.nFaceSplitInterval()) == 0)
3312 )
3313 )
3314 {
3315 doSplit = true;
3316 }
3317
3318
3319
3320 indirectPrimitivePatch& pp = ppPtr();
3321 motionSmoother& meshMover = meshMoverPtr();
3322
3323
3324 // Calculate displacement at every patch point if we need it:
3325 // - if automatic near-surface detection
3326 // - if face splitting active
3327 pointField nearestPoint;
3328 vectorField nearestNormal;
3329
3330 if (snapParams.detectNearSurfacesSnap() || doSplit)
3331 {
3332 nearestPoint.setSize(pp.nPoints(), vector::max);
3333 nearestNormal.setSize(pp.nPoints(), Zero);
3334 }
3335
3336 vectorField disp = calcNearestSurface
3337 (
3338 snapParams.strictRegionSnap(), // attract points to region only
3339 meshRefiner_,
3340 globalToMasterPatch_, // for if strictRegionSnap
3341 globalToSlavePatch_, // for if strictRegionSnap
3342 pp,
3343 pp.localPoints(),
3344 snapDist,
3345
3346 nearestPoint,
3347 nearestNormal
3348 );
3349
3350
3351 // Override displacement at thin gaps
3352 if (snapParams.detectNearSurfacesSnap())
3353 {
3354 detectNearSurfaces
3355 (
3356 Foam::cos(degToRad(planarAngle)),// planar cos for gaps
3357 pp,
3358 pp.localPoints(),
3359 nearestPoint, // surfacepoint from nearest test
3360 nearestNormal, // surfacenormal from nearest test
3361
3362 disp
3363 );
3364 }
3365
3366 // Override displacement with feature edge attempt
3367 if (doFeatures)
3368 {
3369 splitFaces.clear();
3370 splits.clear();
3371 splitPatches.clear();
3372 disp = calcNearestSurfaceFeature
3373 (
3374 snapParams,
3375 !doSplit, // alignMeshEdges
3376 false, // no special handling for >=3 patch points
3377 iter,
3378 featureCos,
3379 scalar(iter+1)/nFeatIter,
3380
3381 snapDist,
3382 disp,
3383 nearestNormal,
3384 pp,
3385 pp.localPoints(),
3386
3387 patchAttraction,
3388 patchConstraints,
3389
3390 splitFaces,
3391 splits,
3392 splitPatches
3393 );
3394 }
3395
3396 // Check for displacement being outwards.
3397 outwardsDisplacement(pp, disp);
3398
3399 // Freeze points on exposed points/faces
3400 freezeExposedPoints
3401 (
3402 meshRefiner_,
3403 "frozenFaces", // faceZone name
3404 "frozenPoints", // pointZone name
3405 pp,
3406 disp
3407 );
3408
3409 // Set initial distribution of displacement field (on patches)
3410 // from patchDisp and make displacement consistent with b.c.
3411 // on displacement pointVectorField.
3412 meshMover.setDisplacement(disp);
3413
3414
3416 {
3417 dumpMove
3418 (
3419 mesh.time().path()
3420 / "patchDisplacement_" + name(iter) + ".obj",
3421 pp.localPoints(),
3422 pp.localPoints() + disp
3423 );
3424 }
3425
3426 // Get smoothly varying internal displacement field.
3427 smoothDisplacement(snapParams, meshMover);
3428
3429 // Apply internal displacement to mesh.
3430 meshOk = scaleMesh
3431 (
3432 snapParams,
3433 nInitErrors,
3434 internalBaffles,
3435 meshMover
3436 );
3437
3438 if (!meshOk)
3439 {
3441 << "Did not successfully snap mesh."
3442 << " Continuing to snap to resolve easy" << nl
3443 << " surfaces but the"
3444 << " resulting mesh will not satisfy your quality"
3445 << " constraints" << nl << endl;
3446 }
3447
3449 {
3450 const_cast<Time&>(mesh.time())++;
3451 Info<< "Writing scaled mesh to time "
3452 << meshRefiner_.timeName() << endl;
3453 meshRefiner_.write
3454 (
3457 (
3460 ),
3461 mesh.time().path()/meshRefiner_.timeName()
3462 );
3463 Info<< "Writing displacement field ..." << endl;
3464 meshMover.displacement().write();
3465 tmp<pointScalarField> magDisp
3466 (
3467 mag(meshMover.displacement())
3468 );
3469 magDisp().write();
3470 }
3471
3472 // Use current mesh as base mesh
3473 meshMover.correct();
3474
3475
3476
3477 // See if any faces need splitting
3478 label nTotalSplit = returnReduce(splitFaces.size(), sumOp<label>());
3479 if (nTotalSplit && doSplit)
3480 {
3481 // Filter out baffle faces from faceZones of type
3482 // internal/baffle
3483
3484 labelList duplicateFace(getInternalOrBaffleDuplicateFace());
3485
3486 {
3487 labelList oldSplitFaces(std::move(splitFaces));
3488 List<labelPair> oldSplits(std::move(splits));
3489 List<labelPair> oldSplitPatches(std::move(splitPatches));
3490 forAll(oldSplitFaces, i)
3491 {
3492 if (duplicateFace[oldSplitFaces[i]] == -1)
3493 {
3494 splitFaces.append(oldSplitFaces[i]);
3495 splits.append(oldSplits[i]);
3496 splitPatches.append(oldSplitPatches[i]);
3497 }
3498 }
3499 nTotalSplit = returnReduce
3500 (
3501 splitFaces.size(),
3502 sumOp<label>()
3503 );
3504 }
3505
3506 // Update mesh
3507
3508 // Reset moving flag to avoid meshPhi problems with topo changes
3509 mesh.moving(false);
3510
3511 meshRefiner_.splitFacesUndo
3512 (
3513 splitFaces,
3514 splits,
3515 splitPatches,
3516 motionDict,
3517
3518 duplicateFace,
3519 internalBaffles
3520 );
3521
3522 // Redo meshMover
3523 meshMoverPtr.clear();
3524 ppPtr.clear();
3525
3526 // Update mesh mover
3527 ppPtr = meshRefinement::makePatch(mesh, adaptPatchIDs);
3528 meshMoverPtr.reset
3529 (
3530 new motionSmoother
3531 (
3532 mesh,
3533 ppPtr(),
3534 adaptPatchIDs,
3536 (
3538 adaptPatchIDs
3539 ),
3540 motionDict,
3541 dryRun_
3542 )
3543 );
3544
3545 // Update snapping distance
3546 snapDist = calcSnapDistance(mesh, snapParams, ppPtr());
3547
3548
3550 {
3551 const_cast<Time&>(mesh.time())++;
3552 Info<< "Writing split-faces mesh to time "
3553 << meshRefiner_.timeName() << endl;
3554 meshRefiner_.write
3555 (
3558 (
3561 ),
3562 mesh.time().path()/meshRefiner_.timeName()
3563 );
3564 }
3565 }
3566
3567
3569 {
3570 forAll(internalBaffles, i)
3571 {
3572 const labelPair& p = internalBaffles[i];
3573 const point& fc0 = mesh.faceCentres()[p[0]];
3574 const point& fc1 = mesh.faceCentres()[p[1]];
3575
3576 if (mag(fc0-fc1) > meshRefiner_.mergeDistance())
3577 {
3579 << "Separated baffles : f0:" << p[0]
3580 << " centre:" << fc0
3581 << " f1:" << p[1] << " centre:" << fc1
3582 << " distance:" << mag(fc0-fc1)
3583 << exit(FatalError);
3584 }
3585 }
3586 }
3587 }
3588 }
3589
3590
3591 // Reset moving flag to avoid any meshPhi mapping problems
3592 mesh.moving(false);
3593
3594
3595 // Merge any introduced baffles (from faceZones of faceType 'internal')
3596 {
3597 autoPtr<mapPolyMesh> mapPtr = meshRefiner_.mergeZoneBaffles
3598 (
3599 true, // internal zones
3600 false // baffle zones
3601 );
3602
3603 if (mapPtr)
3604 {
3606 {
3607 const_cast<Time&>(mesh.time())++;
3608 Info<< "Writing baffle-merged mesh to time "
3609 << meshRefiner_.timeName() << endl;
3610 meshRefiner_.write
3611 (
3614 (
3617 ),
3618 meshRefiner_.timeName()
3619 );
3620 }
3621 }
3622 }
3623
3624 // Repatch faces according to nearest. Do not repatch baffle faces.
3625 {
3626 labelList duplicateFace(getInternalOrBaffleDuplicateFace());
3627
3628 repatchToSurface(snapParams, adaptPatchIDs, duplicateFace);
3629 }
3630
3631 if
3632 (
3635 )
3636 {
3637 labelList duplicateFace(getInternalOrBaffleDuplicateFace());
3638
3639 // Repatching might have caused faces to be on same patch and hence
3640 // mergeable so try again to merge coplanar faces. Do not merge baffle
3641 // faces to ensure they both stay the same.
3642 label nChanged = meshRefiner_.mergePatchFacesUndo
3643 (
3644 featureCos, // minCos
3645 featureCos, // concaveCos
3646 meshRefiner_.meshedPatches(),
3647 motionDict,
3648 duplicateFace, // faces not to merge
3649 mergeType
3650 );
3651
3652 nChanged += meshRefiner_.mergeEdgesUndo(featureCos, motionDict);
3653
3654 if (nChanged > 0 && debug & meshRefinement::MESH)
3655 {
3656 const_cast<Time&>(mesh.time())++;
3657 Info<< "Writing patchFace merged mesh to time "
3658 << meshRefiner_.timeName() << endl;
3659 meshRefiner_.write
3660 (
3663 (
3666 ),
3667 meshRefiner_.timeName()
3668 );
3669 }
3670 }
3671
3673 {
3674 const_cast<Time&>(mesh.time())++;
3675 }
3676}
3677
3678
3679// ************************************************************************* //
Istream and Ostream manipulators taking arguments.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
const polyBoundaryMesh & pbm
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
void clear() noexcept
Clear the addressed list, i.e. set the size to zero.
void append(const T &val)
Copy append an element to the end of this list.
label push_uniq(const T &val)
Append an element if not already in the list.
void setCapacity(const label len)
Alter the size of the underlying storage.
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
SubField< vector > subField
Definition Field.H:183
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field values.
const_iterator cfind(const Key &key) const
Find and return an const_iterator set at the hashed entry.
Definition HashTableI.H:113
label size() const noexcept
The number of elements in table.
Definition HashTable.H:358
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
const Addr & addressing() const noexcept
The addressing used for the list.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
bool set(const label i, bool val=true)
A bitSet::set() method for a list of bool.
Definition List.H:469
void setSize(label n)
Alias for resize().
Definition List.H:536
static FOAM_NO_DANGLING_REFERENCE const pointMesh & New(const polyMesh &mesh, Args &&... args)
An OFstream that keeps track of vertices and provides convenience output methods for OBJ files.
Definition OBJstream.H:58
Output to file stream as an OSstream, normally using std::ofstream for the actual output.
Definition OFstream.H:75
static tmp< pointField > pointNormals(const polyMesh &, const PrimitivePatch< FaceList, PointField > &, const bitSet &flipMap=bitSet::null())
Return parallel consistent point normals for patches using mesh points.
Wave propagation of information through grid. Every iteration information goes through one layer of e...
label nPoints() const
Number of points supporting patch faces.
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 Field< point_type > & faceNormals() const
Return face unit normals for patch.
const Field< point_type > & points() const noexcept
Return reference to global points.
const Field< point_type > & faceCentres() const
Return face centres for patch.
const labelListList & pointFaces() const
Return point-face addressing.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
fileName path() const
The path for the case = rootPath/caseName.
Definition TimePathsI.H:102
const T1 & first() const noexcept
Access the first element.
Definition Tuple2.H:132
A List with indirect addressing. Like IndirectList but does not store addressing.
T & first()
Access first element of the list, position [0].
Definition UList.H:957
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
label size() const noexcept
The number of entries in the list.
Definition UPtrListI.H:106
static const Form max
label findZoneID(const word &zoneName) const
Find zone index by name, return -1 if not found.
Definition ZoneMesh.C:757
Adds layers of cells to outside (or inside) of polyMesh. Can optionally create stand-alone extruded m...
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
void clear() noexcept
Same as reset(nullptr).
Definition autoPtr.H:255
void reset(T *p=nullptr) noexcept
Delete managed object and set to new given pointer.
Definition autoPtrI.H:37
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition bitSet.H:61
void set(const bitSet &bitset)
Set specified bits from another bitset.
Definition bitSetI.H:502
bool test(const label pos) const
Test for True value at specified position, never auto-vivify entries.
Definition bitSet.H:334
labelList toc() const
The indices of the on bits as a sorted labelList.
Definition bitSet.C:476
The coupledPolyPatch is an abstract base class for patches that couple regions of the computational d...
double cpuTimeIncrement() const
Return CPU time [seconds] since last call to cpuTimeIncrement(), resetCpuTimeIncrement().
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
An edge is a list of two vertex labels. This can correspond to a directed graph edge or an edge on a ...
Definition edge.H:62
A subset of mesh faces organised as a primitive patch.
Definition faceZone.H:63
A face is a list of labels corresponding to mesh vertices.
Definition face.H:71
A class for handling file names.
Definition fileName.H:75
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
const Time & time() const
Return the top-level database.
Definition fvMesh.H:360
virtual bool write(const bool writeOnProc=true) const
Write mesh using IO settings from time.
Definition fvMesh.C:1068
virtual void movePoints(const pointField &)
Move points, returns volumes swept by faces in motion.
Definition fvMesh.C:884
label nTotalFaces() const noexcept
Total global number of mesh faces. Not compensated for duplicate faces!
label nTotalPoints() const noexcept
Total global number of mesh points. Not compensated for duplicate points!
Refinement of (split) hexes using polyTopoChange.
Definition hexRef8.H:63
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.
static void mapBaffles(const polyMesh &mesh, const labelList &faceMap, List< labelPair > &baffles)
Map baffles after layer addition. Gets new-to-old face map.
static List< labelPair > subsetBaffles(const polyMesh &mesh, const labelList &zoneIDs, const List< labelPair > &baffles)
Subset baffles according to zones.
writeType
Enumeration for what to write. Used as a bit-pattern.
const refinementSurfaces & surfaces() const
Reference to surface search engines.
word timeName() const
Replacement for Time::timeName() that returns oldInstance (if overwrite_).
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.
static bool checkMesh(const bool report, const polyMesh &mesh, const dictionary &dict, labelHashSet &wrongFaces, const bool dryRun=false)
Check mesh with mesh settings in dict. Collects incorrect faces.
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...
A traits class, which is primarily used for primitives and vector-space.
Definition pTraits.H:64
Application of (multi-)patch point constraints.
void constrainDisplacement(pointVectorField &displacement, const bool overrideValue=false) const
Apply boundary conditions (single-patch constraints),.
Holds information regarding nearest wall point. Used in PointEdgeWave. (so not standard FaceCellWave)...
Mesh representing a set of points created from polyMesh.
Definition pointMesh.H:49
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh,...
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition polyMesh.H:609
const faceZoneMesh & faceZones() const noexcept
Return face zone mesh.
Definition polyMesh.H:671
virtual const faceList & faces() const
Return raw faces.
Definition polyMesh.C:1088
bool moving() const noexcept
Is mesh moving.
Definition polyMesh.H:732
virtual const labelList & faceOwner() const
Return face owner.
Definition polyMesh.C:1101
const globalMeshData & globalData() const
Return parallel info (demand-driven).
Definition polyMesh.C:1296
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition polyMesh.C:1107
const pointZoneMesh & pointZones() const noexcept
Return point zone mesh.
Definition polyMesh.H:663
virtual const pointField & points() const
Return raw points.
Definition polyMesh.C:1063
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
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 writeOnProc=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.
static autoPtr< mapPolyMesh > dupFaceZonePoints(meshRefinement &meshRefiner, const labelList &patchIDs, const labelList &numLayers, List< labelPair > baffles, labelList &pointToMaster)
Duplicate points on faceZones with layers. Re-used when adding buffer layers. Can be made private aga...
static void mapFaceZonePoints(meshRefinement &meshRefiner, const mapPolyMesh &map, labelPairList &baffles, labelList &pointToMaster)
Map numbering after adding cell layers.
All to do with snapping to surface.
void detectNearSurfaces(const scalar planarCos, const indirectPrimitivePatch &pp, const pointField &ppLocalPoints, 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)
Snap onto surface & features.
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).
faceZoneType
What to do with faceZone faces.
static labelList getNamedSurfaces(const PtrList< surfaceZonesInfo > &surfList)
Get indices of named surfaces (surfaces with faceZoneName).
static void swapBoundaryCellList(const polyMesh &mesh, const UList< T > &cellData, List< T > &neighbourCellData, const bool parRun=UPstream::parRun())
Extract and swap to obtain neighbour cell values for all boundary faces.
static bitSet getMasterFaces(const polyMesh &mesh)
Get per face whether it is uncoupled or a master of a coupled set of faces.
Definition syncTools.C:119
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:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition tmpI.H:235
Wrapper for position + weight to be used in e.g. averaging.
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:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
volScalarField & p
auto limits
Definition setRDeltaT.H:186
const polyBoundaryMesh & patches
faceListList boundary
bool coupled
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
auto & name
const pointField & points
Geometric merging of points. See below.
#define WarningInFunction
Report a warning using Foam::Warning.
const dimensionedScalar mp
Proton mass.
Namespace for handling debugging switches.
Definition debug.C:45
const wordList internal
Standard volume internal 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:30
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of a point.
Definition meshTools.C:196
Namespace for OpenFOAM.
List< edge > edgeList
List of edge.
Definition edgeList.H:32
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
Type gAverage(const FieldField< Field, Type > &f, const label comm)
The global arithmetic average of a FieldField.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
Pair< label > labelPair
A pair of labels.
Definition Pair.H:54
List< word > wordList
List of word.
Definition fileName.H:60
bool returnReduceOr(const bool value, const int communicator=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
ZoneMesh< pointZone, polyMesh > pointZoneMesh
A ZoneMesh with pointZone content on a polyMesh.
List< labelPair > labelPairList
List of labelPair.
Definition labelPair.H:33
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
List< label > labelList
A List of labels.
Definition List.H:62
dimensionedSymmTensor sqr(const dimensionedVector &dv)
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition HashSet.H:85
Smanip< std::ios_base::fmtflags > setf(std::ios_base::fmtflags flags)
Definition IOmanip.H:169
messageStream Info
Information stream (stdout output on master, null elsewhere).
ZoneMesh< faceZone, polyMesh > faceZoneMesh
A ZoneMesh with faceZone content on a polyMesh.
List< face > faceList
List of faces.
Definition faceListFwd.H:41
GeometricField< vector, pointPatchField, pointMesh > pointVectorField
Omanip< int > setprecision(const int i)
Definition IOmanip.H:205
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
T returnReduce(const T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
Omanip< int > setw(const int i)
Definition IOmanip.H:199
label mergePoints(const PointList &points, labelList &pointToUnique, labelList &uniquePoints, const scalar mergeTol=SMALL, const bool verbose=false)
Calculate merge mapping, preserving the original point order. All points closer/equal mergeTol are to...
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
constexpr scalar degToRad() noexcept
Multiplication factor for degrees to radians conversion.
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
PrimitivePatch< IndirectList< face >, const pointField & > indirectPrimitivePatch
A PrimitivePatch with an IndirectList for the faces, const reference for the point field.
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
errorManip< error > abort(error &err)
Definition errorManip.H:139
Field< vector > vectorField
Specialisation of Field<T> for vector.
vector point
Point is a vector.
Definition point.H:37
MinMax< Type > gMinMax(const FieldField< Field, Type > &f)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
PrimitivePatch< UIndirectList< face >, const pointField & > uindirectPrimitivePatch
A PrimitivePatch with UIndirectList for the faces, const reference for the point field.
vectorField pointField
pointField is a vectorField.
Vector< scalar > vector
Definition vector.H:57
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
dimensionedScalar cos(const dimensionedScalar &ds)
messageStream Warning
Warning stream (stdout output on master, null elsewhere), with additional 'FOAM Warning' header text.
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
#define addProfiling(Name,...)
Define profiling trigger with specified name and description string. The description is generated by ...
label newPointi
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];} pointMap[start]=pointMap[end]=Foam::min(start, end);} } 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};constexpr 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< DynamicList< face > > pFaces[nBCs]
labelList f(nPoints)
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Unit conversion functions.