Loading...
Searching...
No Matches
motionSmootherAlgo.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-2017 OpenFOAM Foundation
9 Copyright (C) 2018-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
27\*---------------------------------------------------------------------------*/
28
29#include "motionSmootherAlgo.H"
30#include "twoDPointCorrector.H"
31#include "faceSet.H"
32#include "pointSet.H"
34#include "pointConstraints.H"
35#include "syncTools.H"
36#include "meshTools.H"
37#include "OFstream.H"
38
39// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40
41namespace Foam
42{
44}
45
46
47// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
48
49void Foam::motionSmootherAlgo::testSyncPositions
50(
51 const pointField& fld,
52 const scalar maxMag
53) const
54{
55 pointField syncedFld(fld);
56
58 (
59 mesh_,
60 syncedFld,
61 minEqOp<point>(), // combine op
62 point(GREAT,GREAT,GREAT) // null
63 );
64
65 forAll(syncedFld, i)
66 {
67 if (mag(syncedFld[i] - fld[i]) > maxMag)
68 {
70 << "On point " << i << " point:" << fld[i]
71 << " synchronised point:" << syncedFld[i]
72 << abort(FatalError);
73 }
74 }
75}
76
77
78void Foam::motionSmootherAlgo::checkFld(const pointScalarField& fld)
79{
80 forAll(fld, pointi)
81 {
82 const scalar val = fld[pointi];
83
84 if ((val > -GREAT) && (val < GREAT))
85 {}
86 else
87 {
89 << "Problem : point:" << pointi << " value:" << val
90 << abort(FatalError);
91 }
92 }
93}
94
95
96Foam::labelHashSet Foam::motionSmootherAlgo::getPoints
97(
99) const
100{
101 labelHashSet usedPoints(mesh_.nPoints()/100);
102
103 for (const label faceId : faceLabels)
104 {
105 usedPoints.insert(mesh_.faces()[faceId]);
106 }
107
108 return usedPoints;
109}
110
111
112Foam::tmp<Foam::scalarField> Foam::motionSmootherAlgo::calcEdgeWeights
113(
114 const pointField& points
115) const
116{
117 const edgeList& edges = mesh_.edges();
118
119 auto twght = tmp<scalarField>::New(edges.size());
120 auto& wght = twght.ref();
121
122 forAll(edges, edgei)
123 {
124 wght[edgei] = 1.0/(edges[edgei].mag(points)+SMALL);
125 }
126 return twght;
127}
128
129
130// Smooth on selected points (usually patch points)
131void Foam::motionSmootherAlgo::minSmooth
132(
133 const scalarField& edgeWeights,
134 const bitSet& isAffectedPoint,
135 const labelList& meshPoints,
136 const pointScalarField& fld,
137 pointScalarField& newFld
138) const
139{
140 tmp<pointScalarField> tavgFld = avg
141 (
142 fld,
143 edgeWeights //scalarField(mesh_.nEdges(), 1.0) // uniform weighting
144 );
145 const pointScalarField& avgFld = tavgFld();
146
147 for (const label pointi : meshPoints)
148 {
149 if (isAffectedPoint.test(pointi))
150 {
151 newFld[pointi] = min
152 (
153 fld[pointi],
154 0.5*fld[pointi] + 0.5*avgFld[pointi]
155 );
156 }
157 }
158
159 // Single and multi-patch constraints
160 pointConstraints::New(pMesh()).constrain(newFld, false);
161}
162
163
164// Smooth on all internal points
165void Foam::motionSmootherAlgo::minSmooth
166(
167 const scalarField& edgeWeights,
168 const bitSet& isAffectedPoint,
169 const pointScalarField& fld,
170 pointScalarField& newFld
171) const
172{
173 tmp<pointScalarField> tavgFld = avg
174 (
175 fld,
176 edgeWeights //scalarField(mesh_.nEdges(), 1.0) // uniform weighting
177 );
178 const pointScalarField& avgFld = tavgFld();
179
180 forAll(fld, pointi)
181 {
182 if (isAffectedPoint.test(pointi) && isInternalPoint_.test(pointi))
183 {
184 newFld[pointi] = min
185 (
186 fld[pointi],
187 0.5*fld[pointi] + 0.5*avgFld[pointi]
188 );
189 }
190 }
191
192 // Single and multi-patch constraints
193 pointConstraints::New(pMesh()).constrain(newFld, false);
194
195}
196
197
198// Scale on all internal points
199void Foam::motionSmootherAlgo::scaleField
200(
202 const scalar scale,
204) const
205{
206 for (const label pointi : pointLabels)
207 {
208 if (isInternalPoint_.test(pointi))
209 {
210 fld[pointi] *= scale;
211 }
212 }
213
214 // Single and multi-patch constraints
215 pointConstraints::New(pMesh()).constrain(fld, false);
216}
217
218
219// Scale on selected points (usually patch points)
220void Foam::motionSmootherAlgo::scaleField
221(
222 const labelList& meshPoints,
224 const scalar scale,
226) const
227{
228 for (const label pointi : meshPoints)
229 {
230 if (pointLabels.found(pointi))
231 {
232 fld[pointi] *= scale;
233 }
234 }
235}
236
237
238// Lower on internal points
239void Foam::motionSmootherAlgo::subtractField
240(
242 const scalar f,
244) const
245{
246 for (const label pointi : pointLabels)
247 {
248 if (isInternalPoint_.test(pointi))
249 {
250 fld[pointi] = max(0.0, fld[pointi]-f);
251 }
252 }
253
254 // Single and multi-patch constraints
256}
257
258
259// Scale on selected points (usually patch points)
260void Foam::motionSmootherAlgo::subtractField
261(
262 const labelList& meshPoints,
264 const scalar f,
266) const
267{
268 for (const label pointi : meshPoints)
269 {
270 if (pointLabels.found(pointi))
271 {
272 fld[pointi] = max(0.0, fld[pointi]-f);
273 }
274 }
275}
276
277
278void Foam::motionSmootherAlgo::getAffectedFacesAndPoints
279(
280 const label nPointIter,
281 const faceSet& wrongFaces,
282
283 labelList& affectedFaces,
284 bitSet& isAffectedPoint
285) const
286{
287 isAffectedPoint.reset();
288 isAffectedPoint.resize(mesh_.nPoints());
289
290 faceSet nbrFaces(mesh_, "checkFaces", wrongFaces);
291
292 // Find possible points influenced by nPointIter iterations of
293 // scaling and smoothing by doing pointCellpoint walk.
294 // Also update faces-to-be-checked to extend one layer beyond the points
295 // that will get updated.
296
297 for (label i = 0; i < nPointIter; i++)
298 {
299 pointSet nbrPoints(mesh_, "grownPoints", getPoints(nbrFaces));
300
301 for (const label pointi : nbrPoints)
302 {
303 for (const label celli : mesh_.pointCells(pointi))
304 {
305 nbrFaces.set(mesh_.cells()[celli]); // set multiple
306 }
307 }
308 nbrFaces.sync(mesh_);
309
310 if (i == nPointIter - 2)
311 {
312 for (const label facei : nbrFaces)
313 {
314 isAffectedPoint.set(mesh_.faces()[facei]); // set multiple
315 }
316 }
317 }
319 affectedFaces = nbrFaces.toc();
320}
321
322
323// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
324
325Foam::motionSmootherAlgo::motionSmootherAlgo
326(
327 polyMesh& mesh,
328 pointMesh& pMesh,
330 pointVectorField& displacement,
331 pointScalarField& scale,
332 pointField& oldPoints,
333 const labelList& adaptPatchIDs,
334 const dictionary& paramDict,
335 const bool dryRun
336)
337:
338 mesh_(mesh),
339 pMesh_(pMesh),
340 pp_(pp),
341 displacement_(displacement),
342 scale_(scale),
343 oldPoints_(oldPoints),
344 adaptPatchIDs_(adaptPatchIDs),
345 paramDict_(paramDict),
346 dryRun_(dryRun),
347 isInternalPoint_(mesh_.nPoints(), true)
349 updateMesh();
350}
351
352
353// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
358
359// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
362{
363 return mesh_;
364}
365
368{
369 return pMesh_;
370}
371
374{
375 return pp_;
376}
377
380{
381 return adaptPatchIDs_;
382}
383
386{
387 return paramDict_;
388}
389
390
392{
393 oldPoints_ = mesh_.points();
394
395 scale_ = 1.0;
397 // No need to update twoDmotion corrector since only holds edge labels
398 // which will remain the same as before. So unless the mesh was distorted
399 // severely outside of motionSmootherAlgo there will be no need.
400}
401
402
404(
405 const labelList& patchIDs,
406 pointVectorField& displacement
407)
408{
409 pointVectorField::Boundary& displacementBf =
410 displacement.boundaryFieldRef();
411
412 // Make consistent with non-adapted bc's by evaluating those now and
413 // resetting the displacement from the values.
414 // Note that we're just doing a correctBoundaryConditions with
415 // fixedValue bc's first.
416 {
417 const labelHashSet adaptPatchSet(patchIDs);
418 const label startOfRequests = UPstream::nRequests();
419 for (auto& pfld : displacementBf)
420 {
421 if (!adaptPatchSet.found(pfld.patch().index()))
422 {
423 pfld.initEvaluate(UPstream::commsTypes::nonBlocking);
424 }
425 }
426 // Wait for outstanding requests (non-blocking)
427 UPstream::waitRequests(startOfRequests);
428
429 for (auto& pfld : displacementBf)
430 {
431 const label patchi = pfld.patch().index();
432 if (!adaptPatchSet.found(patchi))
433 {
435 }
436 }
437 }
438
439 // Multi-patch constraints. Takes priority over any bc.
440 pointConstraints::New(displacement.mesh()).constrainCorners(displacement);
441
442 // Combine any coupled points. This makes sure that the ppMeshPoint
443 // value gets taken over to any coupled points. Note use of min to
444 // make sure zeroFixedValue 'wins'.
446 (
447 displacement.mesh()(),
448 displacement,
449 minMagSqrEqOp<vector>(), // combine op
450 vector::uniform(GREAT) // null value
451 );
452
453 // Adapt the fixedValue bc's (i.e. copy internal point data to
454 // boundaryField for all affected patches) to take the changes caused
455 // by multi-corner constraints into account.
456 for (const label patchi : patchIDs)
457 {
458 displacementBf[patchi] == displacementBf[patchi].patchInternalField();
459 }
460}
461
464{
465 setDisplacementPatchFields(adaptPatchIDs_, displacement_);
466}
467
468
470(
471 const labelList& patchIDs,
473 pointField& patchDisp,
474 pointVectorField& displacementFld
475)
476{
477 const polyMesh& mesh = displacementFld.mesh()();
478
479 // See comment in .H file about shared points.
480 // We want to disallow effect of loose coupled points - we only
481 // want to see effect of proper fixedValue boundary conditions
482
483 const labelList& cppMeshPoints =
484 mesh.globalData().coupledPatch().meshPoints();
485
486 const labelList& ppMeshPoints = pp.meshPoints();
487
488 pointField& displacement = displacementFld.internalFieldRef();
489
490 // - zero out patch points ('master') or points coupled to them ('slave')
491 // - assign master values
492 // - sync so the master values get transferred to the slaves
493 // - take over these values onto the fixedValue bcs. This removes
494 // any inconsistencies (since all bcs get copied from a single, internal
495 // value)
496
497 // Knock out displacement on points which are on pp or coupled
498 // to them since we want 'proper' values from displacement to take
499 // precedence.
500 bitSet isPatchPoint(mesh.nPoints(), ppMeshPoints);
501 const bitSet oldPatchPoint(isPatchPoint);
502 {
503 // Expand patch points to include points coupled to them
505 (
506 mesh,
507 isPatchPoint,
509 0u
510 );
511
512 for (const label pointi : cppMeshPoints)
513 {
514 if (isPatchPoint.test(pointi))
515 {
516 displacement[pointi] = vector::uniform(GREAT);
517 }
518 }
519 }
520
521 // Set internal point data from displacement on combined patch points.
522 forAll(ppMeshPoints, patchPointi)
523 {
524 displacement[ppMeshPoints[patchPointi]] = patchDisp[patchPointi];
525 }
526
527
528 // Adapt the fixedValue bc's (i.e. copy internal point data to
529 // boundaryField for all affected patches)
530 setDisplacementPatchFields(patchIDs, displacementFld);
531
532 pointVectorField::Boundary& displacementBf =
533 displacementFld.boundaryFieldRef();
534 for (const label patchi : patchIDs)
535 {
536 displacementBf[patchi] == displacementBf[patchi].patchInternalField();
537 }
538
539
540 if (debug)
541 {
542 OFstream str(mesh.db().path()/"changedPoints.obj");
543 label nVerts = 0;
544 forAll(ppMeshPoints, patchPointi)
545 {
546 const vector& newDisp = displacement[ppMeshPoints[patchPointi]];
547
548 if (mag(newDisp-patchDisp[patchPointi]) > SMALL)
549 {
550 const point& pt = mesh.points()[ppMeshPoints[patchPointi]];
551
552 meshTools::writeOBJ(str, pt);
553 nVerts++;
554 // Pout<< "Point:" << pt
555 // << " oldDisp:" << patchDisp[patchPointi]
556 // << " newDisp:" << newDisp << endl;
557 }
558 }
559 Pout<< "Written " << nVerts << " points that are changed to file "
560 << str.name() << endl;
561 }
562
563 // Now reset input displacement
564 forAll(ppMeshPoints, patchPointi)
565 {
566 patchDisp[patchPointi] = displacement[ppMeshPoints[patchPointi]];
567 }
568}
569
570
572{
573 setDisplacement(adaptPatchIDs_, pp_, patchDisp, displacement_);
574}
575
576
577// correctBoundaryConditions with fixedValue bc's first.
579(
580 pointVectorField& displacement
581) const
582{
583 const labelHashSet adaptPatchSet(adaptPatchIDs_);
584
585 //const lduSchedule& patchSchedule = mesh_.globalData().patchSchedule();
586
587 auto& displacementBf = displacement.boundaryFieldRef();
588
589 // 1. evaluate on adaptPatches
590 {
591 const label startOfRequests = UPstream::nRequests();
592 for (const label patchi : adaptPatchIDs_)
593 {
594 displacementBf[patchi].initEvaluate
595 (
597 );
598 }
599 // Wait for outstanding requests (non-blocking)
600 UPstream::waitRequests(startOfRequests);
601
602 for (const label patchi : adaptPatchIDs_)
603 {
604 displacementBf[patchi].evaluate(UPstream::commsTypes::nonBlocking);
605 }
606 }
607
608 // 2. evaluate on non-AdaptPatches
609 {
610 const label startOfRequests = UPstream::nRequests();
611 for (auto& pfld : displacementBf)
612 {
613 const label patchi = pfld.patch().index();
614 if (!adaptPatchSet.found(patchi))
615 {
616 pfld.initEvaluate(UPstream::commsTypes::nonBlocking);
617 }
618 }
619 // Wait for outstanding requests (non-blocking)
620 UPstream::waitRequests(startOfRequests);
621
622 for (auto& pfld : displacementBf)
623 {
624 const label patchi = pfld.patch().index();
625 if (!adaptPatchSet.found(patchi))
626 {
628 }
629 }
630 }
631
632 // Multi-patch constraints
633 pointConstraints::New(displacement.mesh()).constrainCorners(displacement);
634
635 // Correct for problems introduced by corner constraints
637 (
638 mesh_,
639 displacement,
640 maxMagSqrEqOp<vector>(), // combine op
641 vector::zero // null value
642 );
643}
644
645
647{
648 // Correct for 2-D motion
649 const twoDPointCorrector& twoDCorrector = twoDPointCorrector::New(mesh_);
650
651 if (twoDCorrector.required())
652 {
653 Info<< "Correcting 2-D mesh motion";
654
655 if (mesh_.globalData().parallel())
656 {
658 << "2D mesh-motion probably not correct in parallel" << endl;
659 }
660
661 // We do not want to move 3D planes so project all points onto those
662 const pointField& oldPoints = mesh_.points();
663 const edgeList& edges = mesh_.edges();
664
665 const labelList& neIndices = twoDCorrector.normalEdgeIndices();
666 const vector& pn = twoDCorrector.planeNormal();
667
668 for (const label edgei : neIndices)
669 forAll(neIndices, i)
670 {
671 const edge& e = edges[edgei];
672
673 point& pStart = newPoints[e.start()];
674
675 pStart += pn*(pn & (oldPoints[e.start()] - pStart));
676
677 point& pEnd = newPoints[e.end()];
678
679 pEnd += pn*(pn & (oldPoints[e.end()] - pEnd));
680 }
681
682 // Correct tangentially
683 twoDCorrector.correctPoints(newPoints);
684 Info<< " ...done" << endl;
685 }
686
687 if (debug)
688 {
689 Pout<< "motionSmootherAlgo::modifyMotionPoints :"
690 << " testing sync of newPoints."
691 << endl;
692 testSyncPositions(newPoints, 1e-6*mesh_.bounds().mag());
693 }
694}
695
696
698{
699 // Make sure to clear out tetPtIs since used in checks (sometimes, should
700 // really check)
701 mesh_.clearTetBasePtIs();
702 pp_.movePoints(mesh_.points());
703}
704
705
707(
708 const scalar errorReduction
709)
710{
711 scalar old = paramDict_.get<scalar>("errorReduction");
712
713 paramDict_.remove("errorReduction");
714 paramDict_.add("errorReduction", errorReduction);
715
716 return old;
717}
718
719
721(
722 labelList& checkFaces,
723 const bool smoothMesh,
724 const label nAllowableErrors
725)
726{
727 List<labelPair> emptyBaffles;
728 return scaleMesh
729 (
730 checkFaces,
731 emptyBaffles,
732 smoothMesh,
733 nAllowableErrors
734 );
735}
736
737
739(
740 labelList& checkFaces,
741 const List<labelPair>& baffles,
742 const bool smoothMesh,
743 const label nAllowableErrors
744)
745{
746 return scaleMesh
747 (
748 checkFaces,
749 baffles,
750 paramDict_,
751 paramDict_,
752 smoothMesh,
753 nAllowableErrors
754 );
755}
756
757
759{
760 // Set newPoints as old + scale*displacement
761
762 // Create overall displacement with same b.c.s as displacement_
763 wordList actualPatchTypes;
764 {
765 const pointBoundaryMesh& pbm = displacement_.mesh().boundary();
766 actualPatchTypes.setSize(pbm.size());
767 forAll(pbm, patchi)
768 {
769 actualPatchTypes[patchi] = pbm[patchi].type();
770 }
771 }
772
773 wordList actualPatchFieldTypes;
774 {
775 const pointVectorField::Boundary& pfld =
776 displacement_.boundaryField();
777 actualPatchFieldTypes.setSize(pfld.size());
778 forAll(pfld, patchi)
779 {
780 if (isA<fixedValuePointPatchField<vector>>(pfld[patchi]))
781 {
782 // Get rid of funny
783 actualPatchFieldTypes[patchi] =
785 }
786 else
787 {
788 actualPatchFieldTypes[patchi] = pfld[patchi].type();
789 }
790 }
791 }
792
793 pointVectorField totalDisplacement
794 (
796 (
797 "totalDisplacement",
798 mesh_.time().timeName(),
799 mesh_,
803 ),
804 scale_*displacement_,
805 actualPatchFieldTypes,
806 actualPatchTypes
807 );
808 correctBoundaryConditions(totalDisplacement);
809
810 if (debug)
811 {
812 Pout<< "scaleMesh : testing sync of totalDisplacement" << endl;
813 testSyncField
814 (
815 totalDisplacement,
817 vector::zero, // null value
818 1e-6*mesh_.bounds().mag()
819 );
820 }
821
822 tmp<pointField> tnewPoints(oldPoints_ + totalDisplacement.primitiveField());
823
824 // Correct for 2-D motion
825 modifyMotionPoints(tnewPoints.ref());
826
827 return tnewPoints;
828}
829
830
832(
833 labelList& checkFaces,
834 const List<labelPair>& baffles,
835 const dictionary& paramDict,
836 const dictionary& meshQualityDict,
837 const bool smoothMesh,
838 const label nAllowableErrors
839)
840{
841 if (!smoothMesh && adaptPatchIDs_.empty())
842 {
844 << "You specified both no movement on the internal mesh points"
845 << " (smoothMesh = false)" << nl
846 << "and no movement on the patch (adaptPatchIDs is empty)" << nl
847 << "Hence nothing to adapt."
848 << exit(FatalError);
849 }
850
851 if (debug)
852 {
853 // Had a problem with patches moved non-synced. Check transformations.
854 const polyBoundaryMesh& patches = mesh_.boundaryMesh();
855
856 Pout<< "Entering scaleMesh : coupled patches:" << endl;
857 forAll(patches, patchi)
858 {
859 if (patches[patchi].coupled())
860 {
861 const coupledPolyPatch& pp =
863
864 Pout<< '\t' << patchi << '\t' << pp.name()
865 << " parallel:" << pp.parallel()
866 << " separated:" << pp.separated()
867 << " forwardT:" << pp.forwardT().size()
868 << endl;
869 }
870 }
871 }
872
873 const scalar errorReduction = get<scalar>
874 (
875 paramDict, "errorReduction", dryRun_, keyType::REGEX_RECURSIVE
876 );
877 const label nSmoothScale = get<label>
878 (
879 paramDict, "nSmoothScale", dryRun_, keyType::REGEX_RECURSIVE
880 );
881
882 // Note: displacement_ should already be synced already from setDisplacement
883 // but just to make sure.
885 (
886 mesh_,
887 displacement_,
889 vector::zero // null value
890 );
891
892 {
893 auto limits = gMinMax(scale_.primitiveField());
894 Info<< "Moving mesh using displacement scaling :"
895 << " min:" << limits.min()
896 << " max:" << limits.max()
897 << endl;
898 }
899
900 // Get points using current displacement and scale. Optionally 2D corrected.
901 pointField newPoints(curPoints());
902
903 // Move. No need to do 2D correction since curPoints already done that.
904 mesh_.movePoints(newPoints);
905 movePoints();
906
907 // Check. Returns parallel number of incorrect faces.
908 faceSet wrongFaces(mesh_, "wrongFaces", mesh_.nFaces()/100+100);
909 checkMesh
910 (
911 false,
912 mesh_,
913 meshQualityDict,
914 checkFaces,
915 baffles,
916 wrongFaces,
917 dryRun_
918 );
919
920 if (returnReduce(wrongFaces.size(), sumOp<label>()) <= nAllowableErrors)
921 {
922 return true;
923 }
924 else
925 {
926 // Sync across coupled faces by extending the set.
927 wrongFaces.sync(mesh_);
928
929 // Special case:
930 // if errorReduction is set to zero, extend wrongFaces
931 // to face-Cell-faces to ensure quick return to previously valid mesh
932
933 if (mag(errorReduction) < SMALL)
934 {
935 labelHashSet newWrongFaces(wrongFaces);
936 for (const label facei : wrongFaces)
937 {
938 const label own = mesh_.faceOwner()[facei];
939 const cell& ownFaces = mesh_.cells()[own];
940
941 newWrongFaces.insert(ownFaces);
942
943 if (facei < mesh_.nInternalFaces())
944 {
945 const label nei = mesh_.faceNeighbour()[facei];
946 const cell& neiFaces = mesh_.cells()[nei];
947
948 newWrongFaces.insert(neiFaces);
949 }
950 }
951 wrongFaces.transfer(newWrongFaces);
952 wrongFaces.sync(mesh_);
953 }
954
955
956 // Find out points used by wrong faces and scale displacement.
957 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
958
959 pointSet usedPoints(mesh_, "usedPoints", getPoints(wrongFaces));
960 usedPoints.sync(mesh_);
961
962
963 // Grow a few layers to determine
964 // - points to be smoothed
965 // - faces to be checked in next iteration
966 bitSet isAffectedPoint(mesh_.nPoints());
967 getAffectedFacesAndPoints
968 (
969 nSmoothScale, // smoothing iterations
970 wrongFaces, // error faces
971 checkFaces,
972 isAffectedPoint
973 );
974
975 if (debug)
976 {
977 Pout<< "Faces in error:" << wrongFaces.size()
978 << " with points:" << usedPoints.size()
979 << endl;
980 }
981
982 if (adaptPatchIDs_.size())
983 {
984 // Scale conflicting patch points
985 scaleField(pp_.meshPoints(), usedPoints, errorReduction, scale_);
986 //subtractField(pp_.meshPoints(), usedPoints, 0.1, scale_);
987 }
988 if (smoothMesh)
989 {
990 // Scale conflicting internal points
991 scaleField(usedPoints, errorReduction, scale_);
992 //subtractField(usedPoints, 0.1, scale_);
993 }
994
995 scalarField eWeights(calcEdgeWeights(oldPoints_));
996
997 for (label i = 0; i < nSmoothScale; ++i)
998 {
999 if (adaptPatchIDs_.size())
1000 {
1001 // Smooth patch values
1002 pointScalarField oldScale(scale_);
1003 minSmooth
1004 (
1005 eWeights,
1006 isAffectedPoint,
1007 pp_.meshPoints(),
1008 oldScale,
1009 scale_
1010 );
1011 checkFld(scale_);
1012 }
1013 if (smoothMesh)
1014 {
1015 // Smooth internal values
1016 pointScalarField oldScale(scale_);
1017 minSmooth(eWeights, isAffectedPoint, oldScale, scale_);
1018 checkFld(scale_);
1019 }
1020 }
1021
1023 (
1024 mesh_,
1025 scale_,
1027 -GREAT // null value
1028 );
1029
1030
1031 if (debug)
1032 {
1033 auto limits = gMinMax(scale_);
1034
1035 Pout<< "scale_ after smoothing :"
1036 << " min:" << limits.min()
1037 << " max:" << limits.max()
1038 << endl;
1040
1041 return false;
1042 }
1043}
1044
1045
1047{
1048 const pointBoundaryMesh& patches = pMesh_.boundary();
1049
1050 // Check whether displacement has fixed value b.c. on adaptPatchID
1051 for (const label patchi : adaptPatchIDs_)
1052 {
1053 if
1054 (
1056 (
1057 displacement_.boundaryField()[patchi]
1058 )
1059 )
1060 {
1062 << "Patch " << patches[patchi].name()
1063 << " has wrong boundary condition "
1064 << displacement_.boundaryField()[patchi].type()
1065 << " on field " << displacement_.name() << nl
1066 << "Only type allowed is "
1067 << fixedValuePointPatchVectorField::typeName
1068 << exit(FatalError);
1069 }
1070 }
1071
1072
1073 // Determine internal points. Note that for twoD there are no internal
1074 // points so we use the points of adaptPatchIDs instead
1075
1076 isInternalPoint_.unset(pp_.meshPoints()); // unset multiple
1077
1078 // Calculate master edge addressing
1079 isMasterEdge_ = syncTools::getMasterEdges(mesh_);
1080}
1081
1082
1083// ************************************************************************* //
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))
labelList faceLabels(nFaceLabels)
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
labelList patchIDs
const polyBoundaryMesh & pbm
const Mesh & mesh() const noexcept
Return const reference to mesh.
static const char *const typeName
Typename for Field.
Definition Field.H:93
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
GeometricBoundaryField< vector, pointPatchField, pointMesh > Boundary
const Internal::FieldType & primitiveField() const noexcept
Return a const-reference to the internal field values.
Internal & internalFieldRef(const bool updateAccessTime=true)
Return a reference to the dimensioned internal field.
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition HashSet.H:229
bool found(const Key &key) const
Same as contains().
Definition HashTable.H:1370
label size() const noexcept
The number of elements in table.
Definition HashTable.H:358
@ NO_REGISTER
Do not request registration (bool: false).
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
const objectRegistry & db() const noexcept
Return the local objectRegistry.
Definition IOobject.C:450
fileName path() const
The complete path for the object (with instance, local,...).
Definition IOobject.C:500
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
void setSize(label n)
Alias for resize().
Definition List.H:536
static FOAM_NO_DANGLING_REFERENCE const pointConstraints & New(const pointMesh &mesh, Args &&... args)
Output to file stream as an OSstream, normally using std::ofstream for the actual output.
Definition OFstream.H:75
virtual const fileName & name() const override
Read/write access to the name of the stream.
Definition OSstream.H:134
bool found(const T &val, label pos=0) const
Same as contains().
Definition UList.H:983
static label nRequests() noexcept
Number of outstanding requests (on the internal list of requests).
@ nonBlocking
"nonBlocking" (immediate) : (MPI_Isend, MPI_Irecv)
Definition UPstream.H:84
static void waitRequests()
Wait for all requests to finish.
Definition UPstream.H:2497
label size() const noexcept
The number of entries in the list.
Definition UPtrListI.H:106
static const Form zero
static Form uniform(const Cmpt &s)
Return a VectorSpace with all elements = s.
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition bitSet.H:61
bool test(const label pos) const
Test for True value at specified position, never auto-vivify entries.
Definition bitSet.H:334
A cell is defined as a list of faces with extra functionality.
Definition cell.H:56
The coupledPolyPatch is an abstract base class for patches that couple regions of the computational d...
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 list of face labels.
Definition faceSet.H:50
virtual void sync(const polyMesh &mesh)
Sync faceSet across coupled patches.
Definition faceSet.C:152
A FixedValue boundary condition for pointField.
@ REGEX_RECURSIVE
Definition keyType.H:88
Given a displacement moves the mesh by scaling the displacement back until there are no more mesh err...
void setDisplacementPatchFields()
Set patch fields on displacement to be consistent with.
tmp< pointField > curPoints() const
Get the current points (oldPoints+scale*displacement).
scalar setErrorReduction(const scalar)
Set the errorReduction (by how much to scale the displacement.
const labelList & adaptPatchIDs() const
Patch labels that are being adapted.
static void setDisplacement(const labelList &patchIDs, const indirectPrimitivePatch &pp, pointField &patchDisp, pointVectorField &displacement)
Set displacement field from displacement on patch points.
static Type get(const dictionary &dict, const word &keyword, const bool noExit, enum keyType::option matchOpt, const Type &defaultValue=Zero)
Wrapper around dictionary::get which does not exit.
void correct()
Take over existing mesh position.
const polyMesh & mesh() const
Reference to mesh.
void movePoints()
Update for new mesh geometry.
const indirectPrimitivePatch & patch() const
Reference to patch.
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.
const dictionary & paramDict() const
void correctBoundaryConditions(pointVectorField &) const
Special correctBoundaryConditions which evaluates fixedValue.
static void setDisplacementPatchFields(const labelList &patchIDs, pointVectorField &pointDisplacement)
Set patch fields on patchIDs to be consistent with.
bool scaleMesh(labelList &checkFaces, const bool smoothMesh=true, const label nAllow=0)
Move mesh with given scale. Return true if mesh ok or has.
void updateMesh()
Update for new mesh topology.
void modifyMotionPoints(pointField &newPoints) const
Apply optional point constraint (2d correction).
A pointBoundaryMesh is a pointPatch list with registered IO, a reference to the associated pointMesh,...
void constrain(GeometricField< Type, pointPatchField, pointMesh > &pf, const bool overrideValue=false) const
Apply boundary conditions (single-patch constraints) and.
void constrainCorners(GeometricField< Type, pointPatchField, pointMesh > &pf) const
Apply patch-patch constraints only.
Mesh representing a set of points created from polyMesh.
Definition pointMesh.H:49
A set of point labels.
Definition pointSet.H:50
virtual void sync(const polyMesh &mesh)
Sync set across coupled patches. Adds coupled points to set.
Definition pointSet.C:153
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
virtual const pointField & points() const
Return raw points.
Definition polyMesh.C:1063
static bitSet getMasterEdges(const polyMesh &mesh)
Get per edge whether it is uncoupled or a master of a coupled set of edges.
Definition syncTools.C:90
static void syncPointPositions(const polyMesh &mesh, List< point > &positions, const CombineOp &cop, const point &nullValue)
Synchronize locations on all mesh points.
Definition syncTools.H:240
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
Class applies a two-dimensional correction to mesh motion point field.
bool required() const
Is 2D correction required, i.e. is the mesh a wedge or slab.
void correctPoints(pointField &p) const
Correct motion points.
const labelList & normalEdgeIndices() const
Return indices of normal edges.
const vector & planeNormal() const
Return plane normal.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
auto limits
Definition setRDeltaT.H:186
const polyBoundaryMesh & patches
bool coupled
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
const pointField & points
label nPoints
label faceId(-1)
#define WarningInFunction
Report a warning using Foam::Warning.
Namespace for handling debugging switches.
Definition debug.C:45
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
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
List< word > wordList
List of word.
Definition fileName.H:60
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
GeometricField< scalar, pointPatchField, pointMesh > pointScalarField
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).
GeometricField< vector, pointPatchField, pointMesh > pointVectorField
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
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.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
errorManip< error > abort(error &err)
Definition errorManip.H:139
vector point
Point is a vector.
Definition point.H:37
MinMax< Type > gMinMax(const FieldField< Field, Type > &f)
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
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
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
labelList f(nPoints)
labelList pointLabels(nPoints, -1)
cellMask correctBoundaryConditions()
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299