Loading...
Searching...
No Matches
hexMeshSmootherMotionSolver.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) 2021,2024-2025 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
30#include "syncTools.H"
31#include "pointConstraints.H"
32#include "unitConversion.H"
33#include "OBJstream.H"
34#include "PatchTools.H"
35//#include "geometricElementTransformPointSmoother.H"
36#include "pointList.H"
37#include "vectorList.H"
38#include "meshPointPatch.H"
39#include "pointSmoother.H"
40#include "fvMesh.H"
44// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
45
46namespace Foam
47{
48 defineTypeNameAndDebug(hexMeshSmootherMotionSolver, 0);
49
51 (
53 hexMeshSmootherMotionSolver,
55 );
56
58 (
60 hexMeshSmootherMotionSolver,
61 displacement
62 );
63}
64
65
66// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
67
68Foam::labelList Foam::hexMeshSmootherMotionSolver::nonConstraintPatches
69(
70 const polyMesh& mesh
71)
72{
73 // Get list of all non-constraint patches. These are the ones where
74 // laplacian smoothing is applied.
75
76 const auto& pbm = mesh.boundaryMesh();
77
78 DynamicList<label> patchIDs(pbm.size());
79 for (const auto& pp : pbm)
80 {
81 if (!polyPatch::constraintType(pp.type()))
82 {
83 patchIDs.append(pp.index());
84 }
85 }
86 return patchIDs;
87}
88
89
91Foam::hexMeshSmootherMotionSolver::makePatch
92(
93 const polyMesh& mesh,
94 const labelList& patchIDs,
95 const labelList& zoneIDs,
96 const pointField& points0
97)
98{
99 // Mark all faces
100 bitSet isPatchFace(mesh.nFaces());
101
102 // Mark all boundary faces (or just patchIDs?)
103 for (const label patchi : patchIDs)
104 {
105 const polyPatch& pp = mesh.boundaryMesh()[patchi];
106 isPatchFace.set(pp.range());
107 }
108
109 const auto& fzs = mesh.faceZones();
110 for (const label zonei : zoneIDs)
111 {
112 isPatchFace.set(fzs[zonei]);
113 }
114
115 syncTools::syncFaceList(mesh, isPatchFace, orEqOp<unsigned int>());
116
117 const labelList patchFaces(isPatchFace.sortedToc());
118
119 return autoPtr<indirectPrimitivePatch>::New
120 (
121 IndirectList<face>(mesh.faces(), patchFaces),
122 points0
123 );
124}
125
126
127void Foam::hexMeshSmootherMotionSolver::checkMesh
128(
129 const pointField& currentPoints,
130 const vectorField& fCtrs,
131 const vectorField& fAreas,
132 const vectorField& cellCtrs,
133 const scalarField& cellVols,
134 labelHashSet& markedFaces,
135 bitSet& markedPoints
136) const
137{
138 // Replacement for motionSmootherAlgo::checkMesh. Adds to markedFaces
139 // any faces that are insufficient quality
140
141
142 markedFaces.clear();
143 markedPoints = false;
144
145
146 /*
147
148 tmp<scalarField> tminCellQ
149 (
150 pointSmoothers::geometricElementTransformPointSmoother::cellQuality
151 (
152 mesh(),
153 currentPoints
154 )
155 );
156 const scalarField& minCellQ = tminCellQ();
157
158 markedPoints.setSize(mesh().nPoints());
159
160 labelHashSet set;
161 DynamicList<label> storage;
162
163 for (label facei = 0; facei < mesh().nFaces(); facei++)
164 {
165 const label own = mesh().faceOwner()[facei];
166 if (minCellQ[own] < VSMALL)
167 {
168 markedFaces.insert(facei);
169 markedPoints.set(mesh().cellPoints(own, set, storage));
170 }
171 else if
172 (
173 mesh().isInternalFace(facei)
174 && minCellQ[mesh().faceNeighbour()[facei]] < VSMALL
175 )
176 {
177 markedFaces.insert(facei);
178 markedPoints.set
179 (
180 mesh().cellPoints(mesh().faceNeighbour()[facei], set, storage)
181 );
182 }
183 }
184 */
185
186 // Get measure of face quality
187 tmp<scalarField> tfaceQ
188 (
189 pointSmoother_->faceQuality
190 (
191 currentPoints,
192 fCtrs,
193 fAreas,
194 cellCtrs,
196 )
197 );
198 const auto& faceQ = tfaceQ();
199
200 markedPoints.setSize(mesh().nPoints());
201 forAll(faceQ, facei)
202 {
203 if (faceQ[facei] < VSMALL)
204 {
205 markedFaces.insert(facei);
206 markedPoints.set(mesh().faces()[facei]);
207 }
208 }
209
210
211 syncTools::syncPointList
212 (
213 mesh(),
214 markedPoints,
215 orEqOp<unsigned int>(),
216 0U
217 );
218
219 // Par sync. TBD.
220 {
221 bitSet isMarkedFace(mesh().nFaces());
222 isMarkedFace.set(markedFaces.toc());
223 syncTools::syncFaceList
224 (
225 mesh(),
226 isMarkedFace,
227 orEqOp<unsigned int>()
228 );
229 markedFaces.insert(isMarkedFace.toc());
230 }
231}
232
233
234//void Foam::hexMeshSmootherMotionSolver::constrainDisplacement
235//(
236// pointField& points
237//) const
238//{
239// // Make sure the points obey the boundary conditions
240// // on pointDisplacement
241//
242// // Update pointDisplacement for suppled points
243// pointDisplacement_.primitiveFieldRef() = points-points0();
244// const pointConstraints& pcs =
245// pointConstraints::New(pointDisplacement_.mesh());
246// pcs.constrainDisplacement(pointDisplacement_, false);
248// points = points0()+pointDisplacement();
249//}
250
251
252bool Foam::hexMeshSmootherMotionSolver::relax
253(
254 const scalarList& relaxationFactors,
255 const bitSet& pointsToRelax,
256 const pointField& initialPoints,
257 const pointField& wantedPoints,
258 pointField& relaxedPoints,
259 labelList& relaxationLevel
260) const
261{
262 // Find relaxation level that makes mesh quality acceptable. Gets given
263 // initial set of mesh points
264
265 relaxedPoints = wantedPoints;
266
267 {
268 vectorField fCtrs(mesh().nFaces());
269 vectorField fAreas(mesh().nFaces());
270 vectorField cellCtrs(mesh().nCells());
271 scalarField cellVols(mesh().nCells());
272
273 // Calculate mesh quantities with new locations
274 const auto& geom =
275 reinterpret_cast<const fvMesh&>(mesh()).geometry();
276 geom.updateGeom
277 (
278 relaxedPoints,
279 mesh().points(), // old points (for avoiding recalculation)
280 fCtrs,
281 fAreas,
282 cellCtrs,
284 );
285
286 // Check the modified face quality. Marks faces with insufficient
287 // quality.
288 labelHashSet markedFaces;
289 bitSet markedPoints;
290 checkMesh
291 (
292 relaxedPoints,
293 fCtrs,
294 fAreas,
295 cellCtrs,
296 cellVols,
297 markedFaces,
298 markedPoints
299 );
300 if (debug)
301 {
302 Pout<< "** hexMeshSmootherMotionSolver::relax : errorfaces:"
303 << markedFaces.size()
304 << " errorpoints:" << markedPoints.count() << endl;
305 }
306 }
307
308
309 // Create a list of relaxation levels
310 // -1 indicates a point which is not to be moved
311 // 0 is the starting value for a moving point
312 relaxationLevel.setSize(mesh().nPoints());
313 relaxationLevel = -1;
314 for (const label pointi : pointsToRelax)
315 {
316 relaxationLevel[pointi] = 0;
317 }
318
319 syncTools::syncPointList
320 (
321 mesh(),
322 relaxationLevel,
323 maxEqOp<label>(),
324 label(-1)
325 );
326
327 vectorField fCtrs(mesh().nFaces());
328 vectorField fAreas(mesh().nFaces());
329 vectorField cellCtrs(mesh().nCells());
330 scalarField cellVols(mesh().nCells());
331
332 // Loop whilst relaxation levels are being incremented
333 bool complete(false);
334 while (!complete)
335 {
336 //Info<< " Moving "
337 // << countZeroOrPos(relaxationFactors.size(), relaxationLevel)
338 // << " points" << endl;
339
340 // Calculate current points (relaxationLevel >= 0)
341 forAll(relaxationLevel, pointi)
342 {
343 if (relaxationLevel[pointi] >= 0)
344 {
345 const scalar x
346 (
347 relaxationFactors[relaxationLevel[pointi]]
348 );
349
350 relaxedPoints[pointi] =
351 (1 - x)*initialPoints[pointi]
352 + x*wantedPoints[pointi];
353 }
354 }
355
356 // Make sure the relaxed points still obey the boundary conditions
357 // on pointDisplacement. Note: could do this afterwards but better
358 // as soon as possible so we pick it up in the checkMesh
359 //constrainDisplacement(relaxedPoints);
360
361
362 // Calculate mesh quantities with new locations
363 const auto& geom =
364 reinterpret_cast<const fvMesh&>(mesh()).geometry();
365 geom.updateGeom
366 (
367 relaxedPoints,
368 mesh().points(), // old points
369 fCtrs,
370 fAreas,
371 cellCtrs,
373 );
374
375 // Check the modified face quality. Marks faces with insufficient
376 // quality.
377 labelHashSet markedFaces;
378 bitSet markedPoints;
379 checkMesh
380 (
381 relaxedPoints,
382 fCtrs,
383 fAreas,
384 cellCtrs,
385 cellVols,
386 markedFaces,
387 markedPoints
388 );
389 //Pout<< " checkMesh : errorfaces:" << markedFaces.size()
390 // << " errorpoints:" << markedPoints.count() << endl;
391
392 complete = true;
393 for (const label pointi : markedPoints)
394 {
395 if (relaxationLevel[pointi] < relaxationFactors.size() - 1)
396 {
397 ++relaxationLevel[pointi];
398 complete = false;
399 }
400 }
401
402 //Info<< " After adjustment:"
403 // << countZeroOrPos(relaxationFactors.size(), relaxationLevel)
404 // << " points relaxed" << endl;
405
406
407 // Synchronise convergence
408 UPstream::reduceAnd(complete);
409
410 // Synchronise relaxation levels
411 syncTools::syncPointList
412 (
413 mesh(),
414 relaxationLevel,
415 maxEqOp<label>(),
416 label(0)
417 );
418 }
419
420 // Check for convergence
421 const label count(countPos(relaxationLevel));
422 const bool converged(count == 0);
423
424 //if (converged)
425 //{
426 // Info<< "... Converged" << endl << endl;
427 //}
428 //else
429 //{
430 // Info<< "... Not converged" << endl << endl;
431 //}
432
433 return converged;
434}
435Foam::label Foam::hexMeshSmootherMotionSolver::countPos
436(
437 const labelList& elems
438) const
439{
440 label n = 0;
441 for (const label elem : elems)
442 {
443 if (elem > 0)
444 {
445 n++;
446 }
447 }
448 return returnReduce(n, sumOp<label>());
449}
450
451
452Foam::labelList Foam::hexMeshSmootherMotionSolver::countZeroOrPos
453(
454 const label size,
455 const labelList& elems
456) const
457{
458 labelList n(size, Zero);
459 for (const label elem : elems)
460 {
461 if (elem >= 0)
462 {
463 n[elem]++;
464 }
465 }
466
467 Pstream::listReduce(n, sumOp<label>());
468 return n;
469}
470
471
472void Foam::hexMeshSmootherMotionSolver::select
473(
474 const labelUList& lst,
475 const label val,
476 bitSet& isVal
477) const
478{
479 isVal.set(lst.size());
480 isVal = false;
481 forAll(lst, i)
482 {
483 isVal[i] = (lst[i] == val);
484 }
485}
486
487
488void Foam::hexMeshSmootherMotionSolver::laplaceSmooth
489(
490 const label type,
491 const pointField& initialPoints,
492 pointField& newPoints
493) const
494{
495 if (initialPoints.size() != mesh().nPoints())
496 {
497 FatalErrorInFunction << "mesh().nPoints:" << mesh().nPoints()
498 << " initial:" << initialPoints.size() << exit(FatalError);
499 }
500
501 newPoints.setSize(initialPoints.size());
502 newPoints = Zero;
503 labelList n(initialPoints.size(), 0);
504
505 DynamicList<label> storage;
506 forAll(pointTypes_, pointi)
507 {
508 if (pointTypes_[pointi] == INTERIOR)
509 {
510 const labelList& pPoints = mesh().pointPoints(pointi, storage);
511 for (const label otherPointi : pPoints)
512 {
513 if (isMasterPoint_[otherPointi])
514 {
515 newPoints[pointi] += initialPoints[otherPointi];
516 n[pointi]++;
517 }
518 }
519 //Pout<< "Moving internal point " << initialPoints[pointi]
520 // << " to average " << newPoints[pointi]/n[pointi]
521 // << " of " << n[pointi] << " points" << endl;
522 }
523 }
524
525 // Combine
526 syncTools::syncPointList
527 (
528 mesh(),
529 n,
530 plusEqOp<label>(),
531 label(0)
532 );
533 syncTools::syncPointList
534 (
535 mesh(),
536 newPoints,
537 plusEqOp<vector>(),
538 vector::zero
539 );
540 forAll(newPoints, pointi)
541 {
542 if (n[pointi] == 0)
543 {
544 // This can happen if not interior point
545 newPoints[pointi] = initialPoints[pointi];
546 //Pout<< "Not Moving boundary point " << newPoints[pointi] << endl;
547 }
548 else
549 {
550 newPoints[pointi] /= n[pointi];
551 //Pout<< "Moving internal point " << initialPoints[pointi]
552 // << " to " << newPoints[pointi] << endl;
553 }
554 }
555}
556void Foam::hexMeshSmootherMotionSolver::featLaplaceSmooth
557(
558 const indirectPrimitivePatch& pp,
559 const pointField& initialPoints,
560 pointField& newPoints
561) const
562{
563 if (initialPoints.size() != pp.nPoints())
564 {
565 FatalErrorInFunction << "pp.nPoints:" << pp.nPoints()
566 << " initial:" << initialPoints.size() << exit(FatalError);
567 }
568
569 newPoints.setSize(pp.nPoints());
570 newPoints = Zero;
571 labelList n(pp.nPoints(), 0);
572
573 const edgeList& edges = pp.edges();
574 const labelListList& pointEdges = pp.pointEdges();
575 const labelList& meshPoints = pp.meshPoints();
576
577 forAll(pointEdges, pointi)
578 {
579 const label myConstraint = pointTypes_[meshPoints[pointi]];
580 if (myConstraint != INTERIOR) // pp points should never be interior
581 {
582 const labelList& pEdges = pointEdges[pointi];
583 //Pout<< "For boundary point:" << initialPoints[pointi]
584 // << endl;
585
586 for (const label edgei : pEdges)
587 {
588 const label otherPointi = edges[edgei].otherVertex(pointi);
589 const label otherMeshPointi = meshPoints[otherPointi];
590 const label otherConstraint = pointTypes_[otherMeshPointi];
591
592 if
593 (
594 (otherConstraint != INTERIOR) // Should not happen
595 && (myConstraint <= otherConstraint)
596 && isMasterPoint_[otherMeshPointi]
597 )
598 {
599 //Pout<< " summing boundary point:"
600 // << initialPoints[otherPointi] << endl;
601
602 newPoints[pointi] += initialPoints[otherPointi];
603 n[pointi]++;
604 }
605 }
606 }
607 }
608
609 // Combine
610 syncTools::syncPointList
611 (
612 mesh(),
613 meshPoints,
614 n,
615 plusEqOp<label>(),
616 label(0)
617 );
618 syncTools::syncPointList
619 (
620 mesh(),
621 meshPoints,
622 newPoints,
623 plusEqOp<vector>(),
624 vector::zero
625 );
626
627 forAll(newPoints, pointi)
628 {
629 if (n[pointi] == 0)
630 {
631 // This can happen if surface point surrounded by feature points
632 // only.
633 newPoints[pointi] = initialPoints[pointi];
634 //Pout<< "Not Moving boundary point " << newPoints[pointi] << endl;
635 }
636 else
637 {
638 newPoints[pointi] /= n[pointi];
639 //Pout<< "Moving surface point " << initialPoints[pointi]
640 // << " to average " << newPoints[pointi]
641 // << " of " << n[pointi] << " points" << endl;
642 }
643 }
644}
645void Foam::hexMeshSmootherMotionSolver::snapBoundaryPoints
646(
647 const scalar scale,
648 const pointField& initialPoints,
649 pointField& newPoints
650) const
651{
652 if (initialPoints.size() != pointDisplacement_.mesh().size())
653 {
655 << "mesh.nPoints():" << pointDisplacement_.mesh().size()
656 << " initial:" << initialPoints.size() << exit(FatalError);
657 }
658 const indirectPrimitivePatch& bnd0 = bnd0Ptr_();
659 const labelList& mp = bnd0.meshPoints();
660
661 // Save old point location
662 const vectorField bndPoints(initialPoints, mp);
663
664 // Update pointDisplacement_ to be consistent with mesh points being set to
665 // initialPoints. This makes sure that the snapping is done using the
666 // initialPoints as starting point
667 pointDisplacement_.primitiveFieldRef() = initialPoints-points0();
668 // 'snap' using boundary conditions
669 pointDisplacement_.correctBoundaryConditions();
670
671 // Calculate new position
672 newPoints = points0() + pointDisplacement().internalField();
673
674 if (scale < 1.0)
675 {
676 // Underrelax
677 vectorField d(newPoints, mp);
678 d -= bndPoints;
679 d *= scale;
680 d += bndPoints;
681 UIndirectList<point>(newPoints, mp) = d;
682 }
683}
684
685//void Foam::hexMeshSmootherMotionSolver::writeOBJ
686//(
687// const fileName& name,
688// const pointField& p0,
689// const pointField& p1
690//) const
691//{
692// OBJstream os(mesh().time().path()/name);
693// forAll(p0, pointi)
694// {
695// os.write(linePointRef(p0[pointi], p1[pointi]));
696// }
697// Pout<< "Dumped to " << os.name() << endl;
698//}
699//void Foam::hexMeshSmootherMotionSolver::writeOBJ
700//(
701// const fileName& name,
702// const UIndirectList<face>& pp,
703// const pointField& points
704//) const
705//{
706// const faceList fcs(pp);
707//
708// OBJstream os(mesh().time().path()/name);
709// os.write(fcs, points, false);
710// Pout<< "Dumped faces to " << os.name() << endl;
711//}
712
713
714void Foam::hexMeshSmootherMotionSolver::emptyCorrectPoints
715(
716 pointVectorField& pointDisplacement
717) const
718{
719 // Assume empty point patches are already in correct location
720 // so knock out any off-plane displacement.
721 auto& fld = pointDisplacement.primitiveFieldRef();
722 for (const auto& ppf : pointDisplacement.boundaryField())
723 {
724 if (isA<emptyPointPatchVectorField>(ppf))
725 {
726 const auto& mp = ppf.patch().meshPoints();
727 forAll(mp, i)
728 {
729 pointConstraint pc;
730 ppf.patch().applyConstraint(i, pc);
731 fld[mp[i]] = pc.constrainDisplacement(fld[mp[i]]);
732 }
733 }
734 }
735
736 pointField wantedPoints(points0() + fld);
737 twoDCorrectPoints(wantedPoints);
738 fld = wantedPoints-points0();
739}
740
741
742// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
743
744Foam::hexMeshSmootherMotionSolver::
745hexMeshSmootherMotionSolver
746(
747 const polyMesh& mesh,
748 const IOdictionary& dict
749)
750:
751 displacementMotionSolver(mesh, dict, typeName),
752 pointSmoother_(pointSmoother::New(mesh, coeffDict())),
753 nPointSmootherIter_
754 (
755 coeffDict().get<label>("nPointSmootherIter")
756 ),
757 relaxationFactors_(coeffDict().get<scalarList>("relaxationFactors")),
758 relaxationLevel_(mesh.nPoints(), 0),
759 relaxedPoints_(mesh.points()),
760 //surfacesDict_(coeffDict().subDict("geometry")),
761 //featureAngle_(coeffDict().get<scalar>("featureAngle")),
762 //snapPatches_
763 //(
764 // mesh.boundaryMesh().patchSet
765 // (
766 // coeffDict().get<wordRes>("patches")
767 // ).sortedToc()
768 //),
769 //snapZones_
770 //(
771 // mesh.faceZones().indices
772 // (
773 // coeffDict().get<wordRes>("faceZones")
774 // )
775 //),
776 snapScale_(Function1<scalar>::New("snapScale", coeffDict())),
777 isMasterPoint_(syncTools::getMasterPoints(mesh)),
778 // Create big primitivePatch for all outside and any features on it
779 //bnd0Ptr_(makePatch(mesh, snapPatches_, snapZones_, points0()))
780 bnd0Ptr_
781 (
782 makePatch
783 (
784 mesh,
785 nonConstraintPatches(mesh),
786 labelList::null(),
787 points0()
788 )
789 )
790{
791// findSurfaces();
792//
793// // Do multi-patch constraints
794// // ~~~~~~~~~~~~~~~~~~~~~~~~~~
795//
796// const scalar featureEdgeCos(Foam::cos(featureAngle_));
797// const scalar featurePointCos(featureEdgeCos);
798//
799// const indirectPrimitivePatch& bnd0 = bnd0Ptr_();
800//
801// calcConstraints
802// (
803// featureEdgeCos,
804// featurePointCos,
805// bnd0,
806// bnd0EdgeConstraints_,
807// bnd0PointConstraints_
808// );
809
811
812 pointTypes_.setSize(pMesh.size());
813 pointTypes_ = INTERIOR;
814
815 for (const auto& pp : pMesh.boundary())
816 {
818 {
819 const auto& mp = pp.meshPoints();
820 UIndirectList<label>(pointTypes_, mp) = pointType::SURFACE;
821 }
822 }
823
824 // Override with any explicit constraint boundaries
825 for (const auto& pp : pMesh.boundary())
826 {
827 const auto* meshPointPtr = isA<meshPointPatch>(pp);
828 if (meshPointPtr)
829 {
830 const auto& constraints = meshPointPtr->constraints();
831 const auto& mp = meshPointPtr->meshPoints();
832
833 forAll(mp, i)
834 {
835 pointTypes_[mp[i]] = pointType(constraints[i].first());
836 }
837 }
838 }
839
840 // Make sure coupled points agree. Max constraint wins.
842 (
843 mesh,
844 pointTypes_,
846 label(0)
847 );
848
849 bitSet isVal;
850 select(pointTypes_, POINT, isVal);
851 const label nFeatPoint = returnReduce(isVal.count(), sumOp<label>());
852 select(pointTypes_, EDGE, isVal);
853 const label nFeatEdge = returnReduce(isVal.count(), sumOp<label>());
854 select(pointTypes_, SURFACE, isVal);
855 const label nSurface = returnReduce(isVal.count(), sumOp<label>());
856 select(pointTypes_, INTERIOR, isVal);
857 const label nInternal = returnReduce(isVal.count(), sumOp<label>());
858 Info<< "Attraction:" << nl
859 << " feature point:" << nFeatPoint << nl
860 << " feature edge :" << nFeatEdge << nl
861 << " surface :" << nSurface << nl
862 << " none :" << nInternal
863 << endl;
864}
865
866
867Foam::hexMeshSmootherMotionSolver::
868hexMeshSmootherMotionSolver
869(
870 const polyMesh& mesh,
871 const IOdictionary& dict,
872 const pointVectorField& pointDisplacement,
873 const pointIOField& points0
874)
875:
876 displacementMotionSolver(mesh, dict, pointDisplacement, points0, typeName),
877 pointSmoother_(pointSmoother::New(mesh, coeffDict())),
878 //pointSmoother_
879 //(
880 // pointSmoother::New
881 // (
882 // coeffDict(),
883 // displacementMotionSolver::pointDisplacement()
884 // )
885 //),
886 nPointSmootherIter_
887 (
888 coeffDict().get<label>("nPointSmootherIter")
889 ),
890 relaxationFactors_(coeffDict().get<scalarList>("relaxationFactors")),
891 relaxationLevel_(mesh.nPoints(), 0),
892 relaxedPoints_(mesh.points()),
893 //surfacesDict_(coeffDict().subDict("geometry")),
894 //featureAngle_(coeffDict().get<scalar>("featureAngle")),
895 //snapPatches_
896 //(
897 // mesh.boundaryMesh().patchSet
898 // (
899 // coeffDict().get<wordReList>("patches")
900 // ).sortedToc()
901 //),
902 //snapZones_
903 //(
904 // mesh.faceZones().indices
905 // (
906 // coeffDict().get<wordRes>("faceZones")
907 // )
908 //),
909 snapScale_(Function1<scalar>::New("snapScale", coeffDict())),
910 isMasterPoint_(syncTools::getMasterPoints(mesh)),
911 // Create big primitivePatch for all outside and any features on it
912 //bnd0Ptr_(makePatch(mesh, snapPatches_, snapZones_, points0))
913 bnd0Ptr_
914 (
915 makePatch
916 (
917 mesh,
918 nonConstraintPatches(mesh),
919 labelList::null(),
920 points0
921 )
922 )
923{
924// findSurfaces();
925//
926// const scalar featureEdgeCos(Foam::cos(featureAngle_));
927// const scalar featurePointCos(featureEdgeCos);
928//
929// const indirectPrimitivePatch& bnd0 = bnd0Ptr_();
930//
931// calcConstraints
932// (
933// featureEdgeCos,
934// featurePointCos,
935// bnd0,
936// bnd0EdgeConstraints_,
937// bnd0PointConstraints_
938// );
939
941
942 pointTypes_.setSize(mesh.nPoints());
943 pointTypes_ = INTERIOR;
944
945 for (const auto& pp : pMesh.boundary())
946 {
947 if (!isA<meshPointPatch>(pp) && !pp.coupled())
948 {
949 const auto& mp = pp.meshPoints();
950 UIndirectList<label>(pointTypes_, mp) = pointType::SURFACE;
951 }
952 }
953
954 // Override with any explicit constraint boundaries
955 for (const auto& pp : pMesh.boundary())
956 {
957 const auto* meshPointPtr = isA<meshPointPatch>(pp);
958 if (meshPointPtr)
959 {
960 const auto& constraints = meshPointPtr->constraints();
961 const auto& mp = meshPointPtr->meshPoints();
962
963 forAll(mp, i)
964 {
965 pointTypes_[mp[i]] = pointType(constraints[i].first());
966 }
967 }
968 }
969
970 // Make sure coupled points agree. Max constraint wins.
972 (
973 mesh,
974 pointTypes_,
976 label(0)
977 );
978
979 bitSet isVal;
980 select(pointTypes_, POINT, isVal);
981 const label nFeatPoint = returnReduce(isVal.count(), sumOp<label>());
982 select(pointTypes_, EDGE, isVal);
983 const label nFeatEdge = returnReduce(isVal.count(), sumOp<label>());
984 select(pointTypes_, SURFACE, isVal);
985 const label nSurface = returnReduce(isVal.count(), sumOp<label>());
986 select(pointTypes_, INTERIOR, isVal);
987 const label nInternal = returnReduce(isVal.count(), sumOp<label>());
988 Info<< "Attraction:" << nl
989 << " feature point:" << nFeatPoint << nl
990 << " feature edge :" << nFeatEdge << nl
991 << " surface :" << nSurface << nl
992 << " none :" << nInternal
993 << endl;
994}
995
996
997// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
998
999Foam::hexMeshSmootherMotionSolver::
1000~hexMeshSmootherMotionSolver()
1001{}
1002
1003
1004// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
1005
1007Foam::hexMeshSmootherMotionSolver::curPoints() const
1008{
1009 //Note: twoDCorrect already done by ::solve
1010 return relaxedPoints_;
1011}
1012
1013
1014void Foam::hexMeshSmootherMotionSolver::solve()
1015{
1016 // Update any internal storage for current state
1017 movePoints(mesh().points());
1018
1019 // No updating of bc since we don't want to do snapping yet - is done
1020 // later on
1021 //pointDisplacement().boundaryFieldRef().updateCoeffs();
1022
1023 const indirectPrimitivePatch& bnd0 = bnd0Ptr_();
1024 const labelList& mp = bnd0.meshPoints();
1025
1026 // Points on boundary
1027 bitSet isBndPoint(mesh().nPoints(), false);
1028 isBndPoint.set(mp);
1029
1030 // Points not on boundary
1031 bitSet isInternalPoint;
1032 select(pointTypes_, INTERIOR, isInternalPoint);
1033
1034
1035 const pointField& initialPoints = mesh().points();
1036
1037 // Wanted locations - starts off from current mesh points. Note : could
1038 // use relaxedPoints_ but this should be equal to current mesh points unless
1039 // we have multiple motion solvers ...
1040 pointField movedPoints(initialPoints);
1041
1042
1043 // -1 indicates a point which is not to be moved
1044 // 0 is the starting value for a moving point
1045 const label nRelaxed = countPos(relaxationLevel_);
1046 //Pout<< "Starting relaxed:" << nRelaxed << endl;
1047
1048 if (nRelaxed > 0)
1049 {
1050 // MeshSmoother::snapSmoothing()
1051 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1052 // - set initialPoint from relaxedPoint
1053 // - smooth internal points
1054 // - and relax point motion (= adapt relaxationLevel_)
1055 // - snap boundary points
1056 // - and relax point motion (= adapt relaxationLevel_)
1057 // - mark all points that are not snapped (i.e. snapped but relaxed)
1058 // - smooth all boundary points (using features)
1059 // - and relax point motion (= adapt relaxationLevel_)
1060
1061 // Laplace smoothing of interior (=unconstrained) points
1062 laplaceSmooth(INTERIOR, initialPoints, movedPoints);
1063
1064 // Apply constraints
1065 //constrainDisplacement(movedPoints);
1066
1067 // Determine relaxation level to move points
1068 relax
1069 (
1070 relaxationFactors_,
1071 isInternalPoint,
1072 initialPoints, // starting location
1073 movedPoints, // wanted location
1074 relaxedPoints_, // inbetween location without errors
1075 relaxationLevel_
1076 );
1077 //Pout<< "After laplaceSmooth:"
1078 // << countZeroOrPos(relaxationFactors_.size(), relaxationLevel_)
1079 // << endl;
1080 //writeOBJ
1081 //(
1082 // "laplaceSmooth_relax_" + mesh().time().timeName() + ".obj",
1083 // bnd0,
1084 // relaxedPoints_
1085 //);
1086
1087
1088 // Snap boundary points
1089 const scalar scale = snapScale_->value(mesh().time().timeIndex());
1090 if (scale > 0)
1091 {
1092 // snap to surface (=apply boundary conditions)
1093 snapBoundaryPoints(scale, initialPoints, movedPoints);
1094
1095 // Apply constraints
1096 //constrainDisplacement(movedPoints);
1097
1098 relax
1099 (
1100 relaxationFactors_,
1101 isBndPoint,
1102 initialPoints,
1103 movedPoints,
1104 relaxedPoints_,
1105 relaxationLevel_
1106 );
1107 //Pout<< "After snapping:"
1108 // << countZeroOrPos(relaxationFactors_.size(), relaxationLevel_)
1109 // << endl;
1110 //writeOBJ("snap_relax.obj", initialPoints, relaxedPoints_);
1111
1112 // Now relaxedPoints_ with relaxationLevel_ 0 are perfectly snapped
1113 // (only applicable for bnd0 points)
1114 }
1115
1116 // Laplace smoothing of (now snapped&relaxed) boundary points. Use
1117 // average of surrounding boundary points of same type only
1118 pointField bndMovedPoints;
1119 featLaplaceSmooth
1120 (
1121 bnd0,
1122 pointField(relaxedPoints_, mp),
1123 bndMovedPoints
1124 );
1125 UIndirectList<point>(movedPoints, mp) = bndMovedPoints;
1126
1127 // Apply constraints
1128 //constrainDisplacement(movedPoints);
1129 //writeOBJ("featLaplaceSmooth.obj", initialPoints, movedPoints);
1130
1131 relax
1132 (
1133 relaxationFactors_,
1134 isBndPoint,
1135 initialPoints,
1136 movedPoints,
1137 relaxedPoints_,
1138 relaxationLevel_
1139 );
1140 }
1141 else
1142 {
1143 // MeshSmoother::GETMeSmoothing()
1144 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1145 // - set initialPoint from relaxedPoint
1146 // - smooth all boundary points (using features)
1147 // - and relax point motion (= adapt relaxationLevel_)
1148 // - snap boundary points
1149 // - and relax point motion (= adapt relaxationLevel_)
1150 // - do all points GETMe
1151 // - and relax point motion (= adapt relaxationLevel_)
1152
1153
1154 // Laplace smoothing of boundary points. Use average of
1155 // surrounding boundary points of same type only. TBD: make consistent
1156 // wrt coupled points not on boundary patch.
1157 pointField bndMovedPoints;
1158 featLaplaceSmooth
1159 (
1160 bnd0,
1161 pointField(relaxedPoints_, mp),
1162 bndMovedPoints
1163 );
1164 UIndirectList<point>(movedPoints, mp) = bndMovedPoints;
1165 //writeOBJ
1166 //(
1167 // "featLaplaceSmooth_unconstrain_"
1168 // + mesh().time().timeName()
1169 // + ".obj",
1170 // pointField(initialPoints, mp),
1171 // bndMovedPoints
1172 //);
1173
1174 // Apply constraints
1175 //constrainDisplacement(movedPoints);
1176
1177 relax
1178 (
1179 relaxationFactors_,
1180 isBndPoint,
1181 initialPoints,
1182 movedPoints,
1183 relaxedPoints_,
1184 relaxationLevel_
1185 );
1186 //writeOBJ
1187 //(
1188 // "featLaplaceSmooth_relax"
1189 // + mesh().time().timeName()
1190 // + ".obj",
1191 // initialPoints,
1192 // relaxedPoints_
1193 //);
1194
1195 // Snap boundary points
1196 const scalar scale = snapScale_->value(mesh().time().timeIndex());
1197 if (scale > 0)
1198 {
1199 // snap to surface (=apply boundary conditions)
1200 snapBoundaryPoints(scale, relaxedPoints_, movedPoints);
1201
1202 // Apply constraints
1203 //constrainDisplacement(movedPoints);
1204 //writeOBJ("snap.obj", initialPoints, movedPoints);
1205
1206 relax
1207 (
1208 relaxationFactors_,
1209 isBndPoint,
1210 initialPoints,
1211 movedPoints,
1212 relaxedPoints_,
1213 relaxationLevel_
1214 );
1215 //writeOBJ("snap_relax.obj", initialPoints, relaxedPoints_);
1216 }
1217
1218 vectorField fCtrs(mesh().nFaces());
1219 vectorField fAreas(mesh().nFaces());
1220 vectorField cellCtrs(mesh().nCells());
1221 scalarField cellVols(mesh().nCells());
1222
1223 movedPoints = relaxedPoints_;
1224
1225 for(label i = 0; i < nPointSmootherIter_; i ++)
1226 {
1227 // Starting from current points do smoothing
1228 // Calculate mesh quantities with new locations
1229
1230 const auto& geom =
1231 reinterpret_cast<const fvMesh&>(mesh()).geometry();
1232 geom.updateGeom
1233 (
1234 movedPoints,
1235 mesh().points(), // old points
1236 fCtrs,
1237 fAreas,
1238 cellCtrs,
1239 cellVols
1240 );
1241
1242 //- Smooth point positions (returned as pointDisplacement w.r.t.
1243 // points0)
1244 pointSmoother_->update
1245 (
1246 identity(mesh().nFaces()),
1247 points0(),
1248 movedPoints,
1249 fCtrs,
1250 fAreas,
1251 cellCtrs,
1252 cellVols,
1253 pointDisplacement_
1254 );
1255 // Keep points on empty patches
1256 emptyCorrectPoints(pointDisplacement());
1257
1258 // Update moving points
1259 movedPoints = points0() + pointDisplacement().internalField();
1260
1262 //snapBoundaryPoints(scale, movedPoints, movedPoints);
1263 }
1264
1265 // snap to surface (=apply boundary conditions)
1266 //snapBoundaryPoints(scale, movedPoints, movedPoints);
1267
1268 //writeOBJ
1269 //(
1270 // "GETMeSmoothing_snapped_"
1271 // + mesh().time().timeName()
1272 // + ".obj",
1273 // pointField(initialPoints, mp),
1274 // pointField(movedPoints, mp)
1275 //);
1276
1277 relax
1278 (
1279 relaxationFactors_,
1280 bitSet(mesh().nPoints(), true),
1281 initialPoints,
1282 movedPoints,
1283 relaxedPoints_,
1284 relaxationLevel_
1285 );
1286 //writeOBJ("GETMeSmoothing_relax.obj", initialPoints, relaxedPoints_);
1287 }
1288}
1289
1290
1291// ************************************************************************* //
label n
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Info<< nl;Info<< "Write faMesh in vtk format:"<< nl;{ vtk::uindirectPatchWriter writer(aMesh.patch(), fileName(aMesh.time().globalPath()/vtkBaseFileName));writer.writeGeometry();globalIndex procAddr(aMesh.nFaces());labelList cellIDs;if(UPstream::master()) { cellIDs.resize(procAddr.totalSize());for(const labelRange &range :procAddr.ranges()) { auto slice=cellIDs.slice(range);slice=identity(range);} } writer.beginCellData(4);writer.writeProcIDs();writer.write("cellID", cellIDs);writer.write("area", aMesh.S().field());writer.write("normal", aMesh.faceAreaNormals());writer.beginPointData(1);writer.write("normal", aMesh.pointAreaNormals());Info<< " "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edges")));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
labelList patchIDs
const polyBoundaryMesh & pbm
@ READ_IF_PRESENT
Reading is optional [identical to LAZY_READ].
static FOAM_NO_DANGLING_REFERENCE const pointMesh & New(const polyMesh &mesh, Args &&... args)
const labelList & meshPoints() const
Return labelList of mesh points in patch.
A List with indirect addressing. Like IndirectList but does not store addressing.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition bitSet.H:61
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Virtual base class for displacement motion solver.
Virtual base class for mesh motion solver.
Mesh representing a set of points created from polyMesh.
Definition pointMesh.H:49
label nPoints() const noexcept
Number of mesh points.
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
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
U
Definition pEqn.H:72
UEqn relax()
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
const pointField & points
label nPoints
const labelIOList & zoneIDs
Definition correctPhi.H:59
return returnReduce(nRefine-oldNRefine, sumOp< label >())
List< bool > select(const label n, const labelUList &locations)
Construct a selection list of bools (all false) with the given pre-size, subsequently add specified l...
Definition BitOps.C:139
unsigned int count(const UList< bool > &bools, const bool val=true)
Count number of 'true' entries.
Definition BitOps.H:73
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
const dimensionedScalar mp
Proton mass.
Namespace for OpenFOAM.
List< edge > edgeList
List of edge.
Definition edgeList.H:32
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
List< label > labelList
A List of labels.
Definition List.H:62
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition HashSet.H:85
messageStream Info
Information stream (stdout output on master, null elsewhere).
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.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
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...
Field< vector > vectorField
Specialisation of Field<T> for vector.
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
vectorField pointField
pointField is a vectorField.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
List< scalar > scalarList
List of scalar.
Definition scalarList.H:32
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
label timeIndex
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
const scalarField & cellVols
Unit conversion functions.
pointField points0(pointIOField(IOobject("points", mesh.time().constant(), polyMesh::meshSubDir, mesh, IOobject::MUST_READ, IOobject::NO_WRITE, IOobject::NO_REGISTER)))