Loading...
Searching...
No Matches
snappySnapDriverBufferLayers.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) 2024 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
26Description
27 All to do with snapping to the surface with buffer layer
28
29\*----------------------------------------------------------------------------*/
30
31#include "snappySnapDriver.H"
32#include "polyTopoChange.H"
33#include "syncTools.H"
34#include "fvMesh.H"
35#include "Time.H"
36#include "OBJstream.H"
37#include "mapPolyMesh.H"
38#include "snapParameters.H"
39#include "unitConversion.H"
40#include "PatchTools.H"
41#include "profiling.H"
42#include "addPatchCellLayer.H"
43#include "snappyLayerDriver.H"
44#include "weightedPosition.H"
45
46#include "localPointRegion.H"
47#include "pointConstraints.H"
49#include "meshPointPatch.H"
50#include "processorPointPatch.H"
51#include "dummyTransform.H"
52#include "faceSet.H"
53#include "motionSmoother.H"
54#include "tetDecomposer.H"
55#include "tetMatcher.H"
56
57// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
58
59Foam::tmp<Foam::pointField> Foam::snappySnapDriver::avg
60(
61 const polyMesh& mesh,
62 const bitSet& isMasterPoint,
63 const indirectPrimitivePatch& pp,
64 const pointField& localPoints
65)
66{
67 const labelListList& pointEdges = pp.pointEdges();
68 const labelList& meshPoints = pp.meshPoints();
69 const edgeList& edges = pp.edges();
70
71 Field<weightedPosition> wps
72 (
73 pointEdges.size(),
75 );
76
77 // Calculate sum of all contributions (so not positions)
78 forAll(pointEdges, verti)
79 {
80 weightedPosition& wp = wps[verti];
81 for (const label edgei : pointEdges[verti])
82 {
83 const label otherVerti = edges[edgei].otherVertex(verti);
84
85 if (isMasterPoint[meshPoints[otherVerti]])
86 {
87 wp.first() += 1.0;
88 wp.second() += localPoints[otherVerti];
89 }
90 }
91 }
92
93 weightedPosition::syncPoints(mesh, meshPoints, wps);
94
95 tmp<pointField> tavg(new pointField(wps.size()));
96 pointField& avg = tavg.ref();
97
98 forAll(wps, verti)
99 {
100 const weightedPosition& wp = wps[verti];
101
102 if (mag(wp.first()) < VSMALL)
103 {
104 // Set to zero?
105 avg[verti] = Zero;
106 }
107 else
108 {
109 avg[verti] = wp.second()/wp.first();
110 }
111 }
112
113 return tavg;
114}
115
116
117Foam::tmp<Foam::pointField>
118Foam::snappySnapDriver::smoothLambdaMuPatchDisplacement
119(
120 const polyMesh& mesh,
122 const List<labelPair>& baffles
123)
124{
125 const bitSet isMasterPoint(syncTools::getMasterPoints(mesh));
126
127 pointField newLocalPoints(pp.localPoints());
128
129 const label iters = 90;
130 const scalar lambda = 0.33;
131 const scalar mu = 0.34;
132
133 for (label iter = 0; iter < iters; iter++)
134 {
135 // Lambda
136 newLocalPoints =
137 (1 - lambda)*newLocalPoints
138 + lambda*avg(mesh, isMasterPoint, pp, newLocalPoints);
139
140 // Mu
141 newLocalPoints =
142 (1 + mu)*newLocalPoints
143 - mu*avg(mesh, isMasterPoint, pp, newLocalPoints);
144 }
145 return newLocalPoints-pp.localPoints();
146}
147
148
149Foam::tmp<Foam::scalarField>
150Foam::snappySnapDriver::wantedThickness
151(
153 const scalar cellSizeFraction
154) const
155{
156 fvMesh& mesh = meshRefiner_.mesh();
157
158 const labelList& cellLevel = meshRefiner_.meshCutter().cellLevel();
159 const labelList& owner = mesh.faceOwner();
160
161 // Undistorted edge length
162 const scalar edge0Len =
163 meshRefiner_.meshCutter().level0EdgeLength();
164
165
167 scalarField& thickness = tthickness.ref();
168
169 labelList maxPointLevel(pp.nPoints(), labelMin);
170
171 forAll(pp, i)
172 {
173 label ownLevel = cellLevel[owner[pp.addressing()[i]]];
174
175 const face& f = pp.localFaces()[i];
176
177 forAll(f, fp)
178 {
179 maxPointLevel[f[fp]] = max(maxPointLevel[f[fp]], ownLevel);
180 }
181 }
182
184 (
185 mesh,
186 pp.meshPoints(),
187 maxPointLevel,
189 labelMin // null value
190 );
191
192 forAll(thickness, pointi)
193 {
194 const scalar edgeLen = edge0Len/(1<<maxPointLevel[pointi]);
195 thickness[pointi] = cellSizeFraction*edgeLen;
196 }
197 return tthickness;
198}
199
200
201//const Foam::pointMesh& Foam::snappySnapDriver::makePointMesh
202//(
203// const indirectPrimitivePatch& pp,
204// const pointConstraintList& pointConstraints,
205// const word& allEdgePatchName,
206// const word& allPointPatchName
207//) const
208//{
209// fvMesh& mesh = meshRefiner_.mesh();
210//
211// if (pointConstraints.size() != pp.nPoints())
212// {
213// FatalErrorInFunction<< "pointConstraints:" << pointConstraints.size()
214// << " pp:" << pp.nPoints() << exit(FatalError);
215// }
216//
217//
218// // Expand pointConstraints to all meshPoints
219// pointConstraintList meshPointConstraints(mesh.nPoints());
220// UIndirectList<pointConstraint>(meshPointConstraints, pp.meshPoints()) =
221// pointConstraints;
222//
223// const auto combineConstraints = [&]
224// (
225// pointConstraint& x,
226// const pointConstraint& y
227// )
228// {
229// x.combine(y);
230// };
231//
232//
233// syncTools::syncPointList
234// (
235// mesh,
236// meshPointConstraints,
237// combineConstraints,
238// pointConstraint(),
239// dummyTransform()
240// );
241//
242//
243// // Sort point constraints:
244// // - 0 : does not happen?
245// // - 1 : attract to surface
246// // - 2 : attract to feature edge
247// // - 3 : attract to feature point
248//
249// DynamicList<label> featEdgeMeshPoints(pointConstraints.size());
250// DynamicList<label> featPointMeshPoints(pointConstraints.size());
251// forAll(meshPointConstraints, pointi)
252// {
253// if (meshPointConstraints[pointi].first() == 2)
254// {
255// featEdgeMeshPoints.append(pointi);
256// }
257// else if (meshPointConstraints[pointi].first() == 3)
258// {
259// featPointMeshPoints.append(pointi);
260// }
261// }
262//
263// // Lookup / construct (from polyPatches) the pointMesh
264// auto& pMesh = pointMesh::New(meshRefiner_.mesh());
265//
266// pointBoundaryMesh& pointBm =
267// const_cast<pointBoundaryMesh&>(pMesh.boundary());
268//
269// // Check if already has constraint patches
270// label edgePatchi = pointBm.findPatchID(allEdgePatchName);
271// if (edgePatchi != -1)
272// {
273// // Delete patch. TBD: clear patchGroup
274// pointBm.set(edgePatchi, nullptr);
275// }
276// label pointPatchi = pointBm.findPatchID(allPointPatchName);
277// if (pointPatchi != -1)
278// {
279// // Delete patch. TBD: clear patchGroup
280// pointBm.set(pointPatchi, nullptr);
281// }
282//
283// // Add additional point patches in order:
284// // - polyPatch based
285// // - featEdge-based constraints
286// // - featPoint-based constraints (so can override edge-based constraints)
287// // - processor boundaries (or should this be all constraint patches, e.g.
288// // symmetry plane. Note: hopefully this is already handled in the
289// // feature extraction ...)
290//
291// if (returnReduce(featEdgeMeshPoints.size(), sumOp<label>()))
292// {
293// if (edgePatchi != -1)
294// {
295// // Override patch
296// const_cast<pointBoundaryMesh&>(pointBm).set
297// (
298// edgePatchi,
299// new meshPointPatch
300// (
301// allEdgePatchName,
302// featEdgeMeshPoints,
303// List<pointConstraint>
304// (
305// meshPointConstraints,
306// featEdgeMeshPoints
307// ),
308// edgePatchi,
309// pointBm,
310// meshPointPatch::typeName
311// )
312// );
313// }
314// else
315// {
316// // Append
317// const_cast<pointBoundaryMesh&>(pointBm).push_back
318// (
319// new meshPointPatch
320// (
321// allEdgePatchName,
322// featEdgeMeshPoints,
323// List<pointConstraint>
324// (
325// meshPointConstraints,
326// featEdgeMeshPoints
327// ),
328// pointBm.size(),
329// pointBm,
330// meshPointPatch::typeName
331// )
332// );
333// }
334// }
335// if (returnReduce(featPointMeshPoints.size(), sumOp<label>()))
336// {
337// if (pointPatchi != -1)
338// {
339// // Override patch
340// const_cast<pointBoundaryMesh&>(pointBm).set
341// (
342// pointPatchi,
343// new meshPointPatch
344// (
345// allPointPatchName,
346// featPointMeshPoints,
347// List<pointConstraint>
348// (
349// meshPointConstraints,
350// featPointMeshPoints
351// ),
352// pointPatchi,
353// pointBm,
354// meshPointPatch::typeName
355// )
356// );
357// }
358// else
359// {
360// // Append
361// const_cast<pointBoundaryMesh&>(pointBm).push_back
362// (
363// new meshPointPatch
364// (
365// allPointPatchName,
366// featPointMeshPoints,
367// List<pointConstraint>
368// (
369// meshPointConstraints,
370// featPointMeshPoints
371// ),
372// pointBm.size(),
373// pointBm,
374// meshPointPatch::typeName
375// )
376// );
377// }
378// }
379//
380// // Shuffle into order
381// labelList oldToNew(pointBm.size());
382// label newPatchi = 0;
383// forAll(pointBm, patchi)
384// {
385// if (!isA<processorPointPatch>(pointBm[patchi]))
386// {
387// oldToNew[patchi] = newPatchi++;
388// }
389// }
390// forAll(pointBm, patchi)
391// {
392// if (isA<processorPointPatch>(pointBm[patchi]))
393// {
394// oldToNew[patchi] = newPatchi++;
395// }
396// }
397// pointBm.reorder(oldToNew, true);
398//
399// return pMesh;
400//}
401
402
403Foam::autoPtr<Foam::displacementMotionSolver>
404Foam::snappySnapDriver::makeMotionSolver
405(
406 const pointMesh& pMesh,
407 const dictionary& snapDict,
408 const labelList& adaptPatchIDs
409// const pointConstraintList& pointConstraints,
410) const
411{
412 fvMesh& mesh = meshRefiner_.mesh();
413
414 tmp<pointVectorField> tallDisp
415 (
417 (
418 pMesh,
419 adaptPatchIDs
420 )
421 );
422
423 // Make sure the pointDisplacement is not registered (since
424 // displacementMotionSolver itself holds it)
425 tallDisp.ref().checkOut();
426
428 (
430 (
431 snapDict.get<word>("solver"),
432 mesh,
434 (
436 (
437 "motionSolverDict",
438 pMesh.thisDb().time().constant(),
439 pMesh.thisDb(),
442 false
443 ),
444 snapDict
445 ),
446 tallDisp(),
448 (
450 (
451 "points0",
452 pMesh.thisDb().time().constant(),
454 pMesh.thisDb(),
457 false
458 ),
459 mesh.points()
460 )
461 )
462 );
463 return motionPtr;
464}
465
466
467void Foam::snappySnapDriver::setDisplacement
468(
470 const pointField& patchDisp, // displacement w.r.t. current mesh
471 const labelList& adaptPatchIDs,
472 const pointField& points0,
473 pointVectorField& fld // displacement w.r.t. points0
474)
475{
476 const pointMesh& pMesh = fld.mesh();
477 const pointField& points = pMesh().points();
478 const labelList& meshPoints = pp.meshPoints();
479
480 if
481 (
482 (points0.size() != points.size())
483 || (points0.size() != fld.size())
484 || (points0.size() != pMesh.size())
485 || (meshPoints.size() != patchDisp.size())
486 )
487 {
489 << "Sizing :"
490 << " points0.size():" << points0.size()
491 << " points.size():" << points.size()
492 << " fld.size():" << fld.size()
493 << " patchDisp.size():" << patchDisp.size()
494 << " meshPoints.size():" << meshPoints.size()
495 << " mesh.nPoints():" << pMesh.size()
496 << exit(FatalError);
497 }
498
499
500 // Problem is that the patchDisp might not be consistent (parallel etc)
501 // across shared points so expand to mesh points
502
503 pointField meshDisp(pMesh.size(), Zero);
504
505 forAll(meshPoints, patchPointi)
506 {
507 const label meshPointi = meshPoints[patchPointi];
508
509 meshDisp[meshPointi] =
510 patchDisp[patchPointi]+points[meshPointi]-points0[meshPointi];
511 }
512
513 // Assign to bc
514 pointVectorField::Boundary& bfld = fld.boundaryFieldRef();
515 for (const label patchi : adaptPatchIDs)
516 {
517 bfld[patchi] == pointField(meshDisp, bfld[patchi].patch().meshPoints());
518 }
519
520 // Apply multi-patch constraints. Problem: most patches originate from
521 // meshing so are type 'wall'. The pointConstraints are only on any points
522 // remaining from the starting mesh.
523 const pointConstraints& pcs = pointConstraints::New(pMesh);
524 pcs.constrainDisplacement(fld, true);
525}
526
527
528Foam::autoPtr<Foam::mapPolyMesh> Foam::snappySnapDriver::addBufferLayers
529(
531 const pointField& thickness,
532 // Layer mesh modifier
533 addPatchCellLayer& addLayer
534)
535{
536 fvMesh& mesh = meshRefiner_.mesh();
537
538 // Introduce single layer of cells. Straight from snappyLayerDriver
539
540 // Global face indices engine
541 const globalIndex globalFaces(mesh.nFaces());
542
543 // Determine extrudePatch.edgeFaces in global numbering (so across
544 // coupled patches). This is used only to string up edges
545 // between coupled
546 // faces (all edges between same (global)face indices get extruded).
547 labelListList edgeGlobalFaces
548 (
550 (
551 mesh,
552 globalFaces,
553 pp
554 )
555 );
556
557
558
559 // Use global edge - face connectivity to
560 // - disable any non-manifold extrusions
561 // - boundary edges can be extruded - is handled by addPatchCellLayer
562 // - what about setting numLayers to 0 to disable layers? Should that
563 // affect buffer layer addition? Guess not since buffer layer is
564 // to aid snapping ...
565 labelList nFaceLayers(pp.size(), 1);
566 labelList nPointLayers(pp.nPoints(), 1);
567 forAll(edgeGlobalFaces, edgei)
568 {
569 if (edgeGlobalFaces[edgei].size() > 2)
570 {
571 const edge& e = pp.edges()[edgei];
572 const edge meshE(pp.meshPoints()[e[0]], pp.meshPoints()[e[1]]);
573
574 nPointLayers[e[0]] = 0;
575 nPointLayers[e[1]] = 0;
576 }
577 //else if (edgeGlobalFaces[edgei].size() == 1)
578 //{
579 // const edge& e = pp.edges()[edgei];
580 // const edge meshE(pp.meshPoints()[e[0]], pp.meshPoints()[e[1]]);
581 //
582 // nPointLayers[e[0]] = 0;
583 // nPointLayers[e[1]] = 0;
584 //}
585 else if (edgeGlobalFaces[edgei].size() == 0)
586 {
587 const edge& e = pp.edges()[edgei];
588 const edge meshE(pp.meshPoints()[e[0]], pp.meshPoints()[e[1]]);
589 FatalErrorInFunction << "Edge:" << meshE
590 << " At:" << meshE.line(mesh.points())
591 << " has no faces!" << exit(FatalError);
592 }
593 }
594
596 (
597 mesh,
598 pp.meshPoints(),
599 nPointLayers,
601 labelMax // null value
602 );
603
604 forAll(pp.localFaces(), facei)
605 {
606 const face& f = pp.localFaces()[facei];
607 const UIndirectList<label> pLayers(nPointLayers, f);
608 if (!pLayers.found(label(1)))
609 {
610 nFaceLayers[facei] = 0;
611 }
612 }
613
614
615
616 // Determine patches for extruded boundary edges. Adds any
617 // additional processor patches (since extruding coupled edge can cause
618 // additional connectivity)
619
620 labelList edgePatchID;
621 labelList edgeZoneID;
622 boolList edgeFlip;
623 labelList inflateFaceID;
625 (
626 meshRefiner_,
627 globalFaces,
628 edgeGlobalFaces,
629 pp,
630
631 edgePatchID,
632 edgeZoneID,
633 edgeFlip,
634 inflateFaceID
635 );
636
637
638 // Mesh topo change engine. Insert current mesh.
639 polyTopoChange meshMod(mesh);
640
641
642 // Add topo regardless of whether extrudeStatus is extruderemove.
643 // Not add layer if patchDisp is zero.
644
645 addLayer.setRefinement
646 (
647 globalFaces,
648 edgeGlobalFaces,
649
650 scalarField(pp.nPoints(), 1), // expansion ratio
651 pp,
652 bitSet(pp.size()), // no flip
653
654 edgePatchID, // boundary patch for extruded boundary edges
655 edgeZoneID, // zone for extruded edges
656 edgeFlip,
657 inflateFaceID,
658
659
660 labelList(0), // exposed patchIDs, not used for adding layers
661 nFaceLayers,
662 nPointLayers,
663 thickness, //patchAttraction,
664
665 meshMod
666 );
667
668 // Apply the stored topo changes to the current mesh.
669 autoPtr<mapPolyMesh> mapPtr = meshMod.changeMesh(mesh, false);
670 mapPolyMesh& map = *mapPtr;
671
672 // Update fields
673 mesh.updateMesh(map);
674
675 // Move mesh (since morphing does not do this)
676 if (map.hasMotionPoints())
677 {
678 mesh.movePoints(map.preMotionPoints());
679 }
680 else
681 {
682 // Hack to remove meshPhi - mapped incorrectly. TBD.
683 mesh.clearOut();
684 }
685
686 // Reset the instance for if in overwrite mode
687 mesh.setInstance(meshRefiner_.timeName());
688
689 // Update numbering on layer
690 addLayer.updateMesh
691 (
692 map,
693 identity(pp.size()),
695 );
696
697 // Update intersections
698 const labelListList addedCells(addLayer.addedCells());
699 bitSet isChangedFace(mesh.nFaces());
700 for (const labelList& faceToCells : addedCells)
701 {
702 for (const label celli : faceToCells)
703 {
704 isChangedFace.set(mesh.cells()[celli]);
705 }
706 }
707 meshRefiner_.updateMesh(map, isChangedFace.toc());
708
710 {
711 const_cast<Time&>(mesh.time())++;
712 Info<< "Writing mesh-with-layer to time "
713 << meshRefiner_.timeName() << endl;
714 meshRefiner_.write
715 (
718 (
721 ),
722 meshRefiner_.timeName()
723 );
724 }
725 return mapPtr;
726}
727
728
729void Foam::snappySnapDriver::doSnapBufferLayers
730(
731 const dictionary& snapDict,
732 const dictionary& motionDict,
733 const meshRefinement::FaceMergeType mergeType,
734 const scalar featureCos,
735 const scalar planarAngle,
736 const snapParameters& snapParams
737)
738{
739 addProfiling(snap, "snappyHexMesh::snap");
740 fvMesh& mesh = meshRefiner_.mesh();
741
742 Info<< nl
743 << "Morphing phase" << nl
744 << "--------------" << nl
745 << endl;
746
747
748 // Name of pointPatch for any points attracted to feature edges
749 //const word allEdgePatchName("boundaryEdges");
750 // Name of pointPatch for any points attracted to feature points
751 //const word allPointPatchName("boundaryPoints");
752
753
754 // faceZone handling
755 // ~~~~~~~~~~~~~~~~~
756 //
757 // We convert all faceZones into baffles during snapping so we can use
758 // a standard mesh motion (except for the mesh checking which for baffles
759 // created from internal faces should check across the baffles). The state
760 // is stored in two variables:
761 // baffles : pairs of boundary faces
762 // duplicateFace : from mesh face to its baffle colleague (or -1 for
763 // normal faces)
764 // There are three types of faceZones according to the faceType property:
765 //
766 // internal
767 // --------
768 // - baffles: need to be checked across
769 // - duplicateFace: from face to duplicate face. Contains
770 // all faces on faceZone to prevents merging patch faces.
771 //
772 // baffle
773 // ------
774 // - baffles: no need to be checked across
775 // - duplicateFace: contains all faces on faceZone to prevent
776 // merging patch faces.
777 //
778 // boundary
779 // --------
780 // - baffles: no need to be checked across. Also points get duplicated
781 // so will no longer be baffles
782 // - duplicateFace: contains no faces on faceZone since both sides can
783 // merge faces independently.
784
785
786
787 // faceZones of type internal
788 const labelList internalFaceZones
789 (
790 meshRefiner_.getZones
791 (
793 (
794 1,
796 )
797 )
798 );
799
800
801 // Create baffles (pairs of faces that share the same points)
802 // Baffles stored as owner and neighbour face that have been created.
803 {
804 List<labelPair> baffles;
805 labelList originatingFaceZone;
806 meshRefiner_.createZoneBaffles
807 (
809 baffles,
810 originatingFaceZone
811 );
812 }
813
814
815 // Duplicate points on faceZones of type boundary
816 meshRefiner_.dupNonManifoldBoundaryPoints();
817
818 // Extract baffles across internal faceZones (for checking mesh quality
819 // across)
820 labelPairList internalBaffles
821 (
822 meshRefiner_.subsetBaffles
823 (
824 mesh,
825 internalFaceZones,
827 )
828 );
829
830
831 bool doFeatures = true; //false;
832 label nFeatIter = 1;
833 if (snapParams.nFeatureSnap() > 0)
834 {
835 doFeatures = true;
836
837 if (!dryRun_)
838 {
839 nFeatIter = snapParams.nFeatureSnap();
840 }
841
842 Info<< "Snapping to features in " << nFeatIter
843 << " iterations ..." << endl;
844 }
845
846 // Get the labels of added patches.
847 const labelList adaptPatchIDs(meshRefiner_.meshedPatches());
848
850 (
852 (
853 mesh,
854 adaptPatchIDs
855 )
856 );
857
858
859 if (debug)
860 {
861 const fileName dir(mesh.time().path());
862 mkDir(dir);
863 OBJstream str
864 (
865 dir
866 / "pp_initial_" + meshRefiner_.timeName() + ".obj"
867 );
868 str.write
869 (
870 ppPtr().localFaces(),
871 ppPtr().localPoints(),
872 false
873 );
874 }
875
876
877 // Maximum distance to attract to nearest feature on surface
878 scalarField snapDist(calcSnapDistance(mesh, snapParams, ppPtr()));
879
880
881 // Smooth patch points. Equivalent of preSmoothPatch but using mesh
882 // motion solver.
883 vectorField patchDisp
884 (
885 smoothLambdaMuPatchDisplacement
886 (
887 mesh,
888 ppPtr(),
889 internalBaffles
890 )
891 );
892 pointField ppLocalPoints(ppPtr().localPoints()+patchDisp);
893
894 const bool smoothInternal = true;
895 if (smoothInternal)
896 {
897 // Create pointMesh with the correct patches
898 //const pointMesh& pMesh = makePointMesh
899 //(
900 // ppPtr(),
901 // patchConstraints,
902 // allEdgePatchName,
903 // allPointPatchName
904 //);
905 const pointMesh& pMesh = pointMesh::New(mesh);
906
908 (
909 makeMotionSolver
910 (
911 pMesh,
912 snapDict,
913 adaptPatchIDs
914 //patchConstraints
915 )
916 );
917
918 // Insert as bc to motionSmoother. Note that this takes displacement
919 // relative to points0
920 setDisplacement
921 (
922 ppPtr(),
923 patchDisp,
924 adaptPatchIDs,
925 motionPtr().points0(),
926 motionPtr().pointDisplacement()
927 );
928 // Solve internal displacement
929 tmp<pointField> tnewPoints(motionPtr->newPoints());
930
931 // Move points
932 mesh.movePoints(tnewPoints);
933
934 // Update pp for new mesh points. Ok as long as we also update geometry
935 // and wanted displacement (usually zero if mesh motion has succeeded)
936 ppLocalPoints = pointField(mesh.points(), ppPtr().meshPoints());
937 patchDisp -= (ppLocalPoints-ppPtr().localPoints());
938
940 {
941 const_cast<Time&>(mesh.time())++;
942 Info<< "Writing smoothed mesh to time "
943 << meshRefiner_.timeName() << endl;
944 meshRefiner_.write
945 (
948 (
951 ),
952 meshRefiner_.timeName()
953 );
954 }
955 }
956
957
958 // Construct iterative mesh mover.
959 Info<< "Constructing mesh displacer ..." << endl;
960 Info<< "Using mesh parameters " << motionDict << nl << endl;
961
963 (
964 makeMotionSolver
965 (
967 snapDict,
968 adaptPatchIDs
969 //patchConstraints
970 )
971 );
972 autoPtr<motionSmoother> meshMoverPtr
973 (
975 (
976 mesh,
977 ppPtr(),
978 adaptPatchIDs,
979 motionPtr().pointDisplacement(),
980 motionDict,
981 dryRun_
982 )
983 );
984
985
986 // Check initial mesh
987 Info<< "Checking initial mesh ..." << endl;
988 labelHashSet wrongFaces(mesh.nFaces()/100);
989 motionSmoother::checkMesh(false, mesh, motionDict, wrongFaces, dryRun_);
990 const label nInitErrors = returnReduce
991 (
992 wrongFaces.size(),
994 );
995
996 Info<< "Detected " << nInitErrors << " illegal faces"
997 << " (concave, zero area or negative cell pyramid volume)"
998 << endl;
999
1000
1001 Info<< "Checked initial mesh in = "
1002 << mesh.time().cpuTimeIncrement() << " s\n" << nl << endl;
1003
1004
1005
1006 //- Only if in feature attraction mode:
1007 // Point on nearest feature
1008 vectorField patchFeaturePoint(ppPtr().nPoints(), Zero);
1009 // Constraints at feature
1010 List<pointConstraint> patchConstraints(ppPtr().nPoints());
1011
1012 for (label iter = 0; iter < nFeatIter; iter++)
1013 {
1014 Info<< nl
1015 << "Morph iteration " << iter << nl
1016 << "-----------------" << endl;
1017
1018
1019 // Calculate displacement at every patch point if we need it:
1020 // - if automatic near-surface detection
1021 // - if face splitting active
1022 pointField nearestPoint(ppPtr().nPoints(), vector::max);
1023 vectorField nearestNormal(ppPtr().nPoints(), Zero);
1024
1025
1026 const bool strictRegionSnap
1027 (
1028 iter < nFeatIter/2
1029 ? snapParams.strictRegionSnap()
1030 : Switch(true)
1031 );
1032
1033 vectorField disp = calcNearestSurface
1034 (
1035 strictRegionSnap, // attract points to region only
1036 meshRefiner_,
1037 globalToMasterPatch_, // for if strictRegionSnap
1038 globalToSlavePatch_, // for if strictRegionSnap
1039 ppPtr(),
1040 ppLocalPoints,
1041 snapDist, // max snap distance
1042
1043 nearestPoint,
1044 nearestNormal
1045 );
1046
1047 // Override displacement at thin gaps
1048 if (snapParams.detectNearSurfacesSnap())
1049 {
1050 detectNearSurfaces
1051 (
1052 Foam::cos(degToRad(planarAngle)),// planar cos for gaps
1053 ppPtr(),
1054 ppLocalPoints,
1055 nearestPoint, // surfacepoint from nearest test
1056 nearestNormal, // surfacenormal from nearest test
1057
1058 disp
1059 );
1060 }
1061
1062 // Override displacement with feature edge attempt
1063 if (doFeatures)
1064 {
1065 //- Any faces to split
1066 DynamicList<label> splitFaces;
1067 //- Indices in face to split across
1069 //- Patches for split face
1070 DynamicList<labelPair> splitPatches;
1071
1072 // Offset to project to nearest feature. Use in combination with
1073 // patchConstraints.
1074 vectorField patchAttraction;
1075
1076 disp = calcNearestSurfaceFeature
1077 (
1078 snapParams,
1079 false, // no alignMeshEdges
1080 true, // check >=3 patch points
1081
1082 iter,
1083 featureCos,
1084 scalar(iter+1)/nFeatIter,
1085
1086 snapDist,
1087 disp, // nearest surface
1088 nearestNormal,
1089 ppPtr(),
1090 ppLocalPoints,
1091
1092 patchAttraction, // offset wrt ppLocalPoints to nearest
1093 // feature edge/point
1094 patchConstraints, // feature type + constraint
1095
1096 splitFaces,
1097 splits,
1098 splitPatches
1099 );
1100
1101
1102 // Freeze points on exposed points/faces
1103 freezeExposedPoints
1104 (
1105 meshRefiner_,
1106 "frozenFaces", // faceZone name
1107 "frozenPoints", // pointZone name
1108 ppPtr(),
1109 disp
1110 );
1111
1112 patchFeaturePoint = ppLocalPoints+patchAttraction;
1113
1114 if (debug)
1115 {
1116 OBJstream str
1117 (
1118 mesh.time().path()
1119 / "calcNearestSurfaceFeature"
1120 + meshRefiner_.timeName()
1121 + ".obj"
1122 );
1123 forAll(ppLocalPoints, pointi)
1124 {
1125 const point& pt = ppLocalPoints[pointi];
1126 str.write(linePointRef(pt, pt+disp[pointi]));
1127 }
1128 }
1129
1130
1131
1132 // Split any faces:
1133 // - does not move/add any points
1134 // - tries to move new faces to correct patch
1135
1136 if (returnReduce(splitFaces.size(), sumOp<label>()))
1137 {
1138 polyTopoChange meshMod(mesh);
1139
1140 // Insert the mesh changes
1141 meshRefiner_.doSplitFaces
1142 (
1143 splitFaces,
1144 splits,
1145 splitPatches,
1146 meshMod
1147 );
1148
1149 // Save old meshPoints before changing mesh
1150 const Map<label> oldMeshPointMap(ppPtr->meshPointMap());
1151 Pout<< "old pp points:" << ppPtr->nPoints()
1152 << " oldMeshPointMap:" << oldMeshPointMap.size()
1153 << endl;
1154
1155 // Remove any unnecessary fields
1156 meshMoverPtr.clear();
1157 motionPtr.clear();
1158 ppPtr.clear();
1159 mesh.clearOut();
1160 mesh.moving(false);
1161
1162 // Change the mesh (no inflation)
1163 autoPtr<mapPolyMesh> mapPtr = meshMod.changeMesh(mesh, false);
1164 mapPolyMesh& map = *mapPtr;
1165
1166 // Update fields
1167 mesh.updateMesh(map);
1168
1169 // Move mesh (since morphing might not do this)
1170 if (map.hasMotionPoints())
1171 {
1172 mesh.movePoints(map.preMotionPoints());
1173 }
1174 else
1175 {
1176 mesh.clearOut();
1177 }
1178
1179 // Reset the instance for if in overwrite mode
1180 mesh.setInstance(meshRefiner_.timeName());
1181 meshRefiner_.setInstance(mesh.facesInstance());
1182
1183 // Update intersections on split faces
1184 {
1185 DynamicList<label> changedFaces(splitFaces.size());
1186 Map<label> splitFacesMap(splitFaces.size());
1187 forAll(splitFaces, i)
1188 {
1189 splitFacesMap.insert(splitFaces[i], i);
1190 }
1191 forAll(map.faceMap(), facei)
1192 {
1193 if (splitFacesMap.find(map.faceMap()[facei]))
1194 {
1195 changedFaces.append(facei);
1196 }
1197 }
1198 // Update intersections on changed faces
1199 meshRefiner_.updateMesh
1200 (
1201 map,
1202 meshRefiner_.growFaceCellFace(changedFaces)
1203 );
1204 }
1205
1207 {
1208 const_cast<Time&>(mesh.time())++;
1209 Info<< "Writing split-faces mesh to time "
1210 << meshRefiner_.timeName() << endl;
1211 meshRefiner_.write
1212 (
1215 (
1218 ),
1219 mesh.time().path()/meshRefiner_.timeName()
1220 );
1221 }
1222
1223
1224
1225 // Update local mesh data
1226 // ~~~~~~~~~~~~~~~~~~~~~~
1227
1228 Info<< "Updating for face-splitting" << endl;
1229
1230 // baffles
1231 forAll(internalBaffles, i)
1232 {
1233 labelPair& baffle = internalBaffles[i];
1234 baffle.first() = map.reverseFaceMap()[baffle.first()];
1235 baffle.second() = map.reverseFaceMap()[baffle.second()];
1236 }
1237
1238 // re-do patch (since faces might have been split)
1239 ppPtr = meshRefinement::makePatch(mesh, adaptPatchIDs);
1240 motionPtr = makeMotionSolver
1241 (
1243 snapDict,
1244 adaptPatchIDs
1245 //patchConstraints
1246 );
1247 meshMoverPtr.reset
1248 (
1249 new motionSmoother
1250 (
1251 mesh,
1252 ppPtr(),
1253 adaptPatchIDs,
1254 motionPtr().pointDisplacement(),
1255 motionDict,
1256 dryRun_
1257 )
1258 );
1259
1260 const auto& mp = ppPtr->meshPoints();
1261 // pointMap (new-to-old) for pp points. Note: no points changed
1262 // but local point ordering might have changed since faces
1263 // split.
1264 labelList ppMap(mp.size());
1265 forAll(mp, i)
1266 {
1267 ppMap[i] = oldMeshPointMap[mp[i]];
1268 }
1269 // patchDisp
1270 meshRefinement::updateList(ppMap, vector::zero, patchDisp);
1271 // snapDist
1272 meshRefinement::updateList(ppMap, scalar(0), snapDist);
1273 // patchFeaturePoint
1275 (
1276 ppMap,
1278 patchFeaturePoint
1279 );
1280 // patchConstraints
1282 (
1283 ppMap,
1285 patchConstraints
1286 );
1287// // nearestPoint
1288// meshRefinement::updateList(ppMap, vector::zero, nearestPoint);
1289// // nearestNormal
1290// meshRefinement::updateList(ppMap, vector::zero, nearestNormal);
1291 // disp
1293 // ppLocalPoints
1294 meshRefinement::updateList(ppMap, vector::zero, ppLocalPoints);
1295
1296 Info<< "DONE Updating for face-splitting" << endl;
1297 }
1298 }
1299
1300
1301 // Attract/slide the mesh a bit
1302 {
1307 //const pointMesh& pMesh = makePointMesh
1308 //(
1309 // ppPtr(),
1310 // patchConstraints,
1311 // allEdgePatchName,
1312 // allPointPatchName
1313 //);
1314 //const pointMesh& pMesh = pointMesh::New(mesh);
1315 //
1316 //autoPtr<displacementMotionSolver> motionPtr
1317 //(
1318 // makeMotionSolver
1319 // (
1320 // pMesh,
1321 // snapDict,
1322 // adaptPatchIDs
1323 // //patchConstraints
1324 // )
1325 //);
1326
1327 // Insert as bc to motionSmoother. Note that this takes displacement
1328 // relative to points0
1329 setDisplacement
1330 (
1331 ppPtr(),
1332 disp,
1333 adaptPatchIDs,
1334 motionPtr().points0(),
1335 motionPtr().pointDisplacement()
1336 );
1337
1338 // Solve internal displacement
1339 tmp<pointField> tnewPoints(motionPtr->newPoints());
1340
1341 // Move points
1342 if (false)
1343 {
1344 // 1. Directly move points
1345 mesh.movePoints(tnewPoints);
1346 // Optional? Only geometry used is ppLocalPoints which we keep
1347 // and update 'by hand'.
1348 //ppPtr().movePoints(tnewPoints);
1349 }
1350 else
1351 {
1352 // 2. Use motionSmoother
1353 // Set initial distribution of displacement field (on patches)
1354 // from patchDisp and make displacement consistent with b.c.
1355 // on displacement pointVectorField.
1356 const vectorField newDisp(tnewPoints()-meshMoverPtr().oldPoints());
1357
1358 meshMoverPtr().displacement().vectorField::operator=(newDisp);
1359 meshMoverPtr().setDisplacement(disp);
1360
1361 // Apply internal displacement to mesh.
1362 const bool meshOk = scaleMesh
1363 (
1364 snapParams,
1365 nInitErrors,
1366 internalBaffles,
1367 meshMoverPtr()
1368 );
1369
1370 if (!meshOk)
1371 {
1373 << "Did not successfully snap mesh."
1374 << " Continuing to snap to resolve easy" << nl
1375 << " surfaces but the"
1376 << " resulting mesh will not satisfy your quality"
1377 << " constraints" << nl << endl;
1378 }
1379
1380 // Use current mesh as base mesh
1381 meshMoverPtr().correct();
1382 }
1383
1384 // Update pp for new mesh points. Ok as long as we also update
1385 // geometry and wanted displacement (usually zero if mesh motion
1386 // has succeeded)
1387 ppLocalPoints = pointField(mesh.points(), ppPtr().meshPoints());
1388 patchDisp -= (ppLocalPoints-ppPtr().localPoints());
1389
1391 {
1392 const_cast<Time&>(mesh.time())++;
1393 Info<< "Writing partially moved mesh to time "
1394 << meshRefiner_.timeName() << endl;
1395 meshRefiner_.write
1396 (
1399 (
1402 ),
1403 meshRefiner_.timeName()
1404 );
1405 }
1406 }
1407
1408 // Split problematic cells
1409 {
1410 Info<< nl << "Checking moved mesh ..." << endl;
1411 faceSet wrongFaces(mesh, "wrongFaces", mesh.nFaces()/1000);
1413 (
1414 false,
1415 mesh,
1416 motionDict,
1417 identity(mesh.nFaces()),
1418 internalBaffles,
1419 wrongFaces,
1420 false // dryRun_
1421 );
1422 const label nWrong(returnReduce(wrongFaces.size(), sumOp<label>()));
1423 Info<< "Detected " << nWrong
1424 << " illegal faces"
1425 << " (concave, zero area or negative cell pyramid volume)"
1426 << endl;
1427
1428 //if (nWrong)
1429 if (false)
1430 {
1431 bitSet decomposeCell(mesh.nCells());
1432 for (const label facei : wrongFaces)
1433 {
1434 const label own = mesh.faceOwner()[facei];
1435 if (!tetMatcher::test(mesh, own))
1436 {
1437 decomposeCell.set(own);
1438 }
1439 if (mesh.isInternalFace(facei))
1440 {
1441 const label nei = mesh.faceNeighbour()[facei];
1442 if (!tetMatcher::test(mesh, nei))
1443 {
1444 decomposeCell.set(nei);
1445 }
1446 }
1447 }
1448 Pout<< "spliyyinG :" << decomposeCell.count() << " cells"
1449 << endl;
1450
1451 tetDecomposer tetDecomp(mesh);
1452 polyTopoChange meshMod(mesh);
1453 tetDecomp.setRefinement
1454 (
1456 decomposeCell,
1457 meshMod
1458 );
1459
1460 // Save old meshPoints before changing mesh
1461 const Map<label> oldMeshPointMap(ppPtr->meshPointMap());
1462
1463 Pout<< "old pp points:" << ppPtr->nPoints()
1464 << " oldMeshPointMap:" << oldMeshPointMap.size()
1465 << endl;
1466
1467 // Remove any unnecessary fields
1468 meshMoverPtr.clear();
1469 motionPtr.clear();
1470 ppPtr.clear();
1471 mesh.clearOut();
1472 mesh.moving(false);
1473
1474 // Change the mesh (no inflation)
1475 autoPtr<mapPolyMesh> mapPtr = meshMod.changeMesh(mesh, false);
1476 mapPolyMesh& map = *mapPtr;
1477
1478 // Update fields
1479 mesh.updateMesh(map);
1480
1481 // Move mesh (since morphing does not do this)
1482 if (map.hasMotionPoints())
1483 {
1484 mesh.movePoints(map.preMotionPoints());
1485 }
1486 else
1487 {
1488 mesh.clearOut();
1489 }
1490
1491 // Reset the instance for if in overwrite mode
1492 mesh.setInstance(meshRefiner_.timeName());
1493 meshRefiner_.setInstance(mesh.facesInstance());
1494
1495 //- Update numbering on tet-decomposition engine
1496 tetDecomp.updateMesh(map);
1497
1498 bitSet isChangedFace(mesh.nFaces());
1499 forAll(map.cellMap(), celli)
1500 {
1501 if (decomposeCell[map.cellMap()[celli]])
1502 {
1503 isChangedFace.set(mesh.cells()[celli]);
1504 }
1505 }
1507 (
1508 mesh,
1509 isChangedFace,
1511 );
1512
1513 Pout<< "isChangedFace :" << decomposeCell.count() << " faces"
1514 << endl;
1515
1516
1517 // Update intersection info
1518 meshRefiner_.updateMesh(map, isChangedFace.toc());
1519
1520
1522 {
1523 const_cast<Time&>(mesh.time())++;
1524 Info<< "Writing tet-decomp mesh to time "
1525 << meshRefiner_.timeName() << endl;
1526 meshRefiner_.write
1527 (
1530 (
1533 ),
1534 mesh.time().path()/meshRefiner_.timeName()
1535 );
1536 }
1537
1538
1539 // Update local mesh data
1540 // ~~~~~~~~~~~~~~~~~~~~~~
1541
1542 Info<< "Updating for tet-decomp" << endl;
1543
1544 // baffles
1545 forAll(internalBaffles, i)
1546 {
1547 labelPair& baffle = internalBaffles[i];
1548 baffle.first() = map.reverseFaceMap()[baffle.first()];
1549 baffle.second() = map.reverseFaceMap()[baffle.second()];
1550 }
1551
1552 // re-do patch (since faces might have been split)
1553 ppPtr = meshRefinement::makePatch(mesh, adaptPatchIDs);
1554 Pout<< "new pp points:" << ppPtr->nPoints()
1555 << " new meshPointMap:" << ppPtr->meshPointMap().size()
1556 << endl;
1557 const auto& mp = ppPtr->meshPoints();
1558 // pointMap (new-to-old) for pp points. Might have new
1559 // face-centre points - these get mapped from any point on
1560 // the originating face.
1561 labelList ppMap(mp.size(), -1);
1562 forAll(mp, i)
1563 {
1564 const label oldMeshPointi = map.pointMap()[mp[i]];
1565 const auto mpFnd = oldMeshPointMap.find(oldMeshPointi);
1566 if (mpFnd)
1567 {
1568 ppMap[i] = mpFnd();
1569 }
1570 }
1571
1572 // patchDisp
1573 meshRefinement::updateList(ppMap, vector::zero, patchDisp);
1574 // snapDist
1575 meshRefinement::updateList(ppMap, scalar(0), snapDist);
1576 // patchFeaturePoint
1578 (
1579 ppMap,
1581 patchFeaturePoint
1582 );
1583 // patchConstraints
1585 (
1586 ppMap,
1588 patchConstraints
1589 );
1590// // nearestPoint
1591// meshRefinement::updateList(ppMap, vector::zero, nearestPoint);
1592// // nearestNormal
1593// meshRefinement::updateList(ppMap, vector::zero, nearestNormal);
1594// // disp
1595// meshRefinement::updateList(ppMap, vector::zero, disp);
1596
1597 // Maximum distance to attract to nearest feature on surface
1598 snapDist = calcSnapDistance(mesh, snapParams, ppPtr());
1599 ppLocalPoints = pointField(mesh.points(), ppPtr().meshPoints());
1600
1601 Info<< "DONE Updating for tet-decomp" << endl;
1602 }
1603 }
1604 }
1605
1606//XXXXXX
1607/*
1608 {
1609 // Introduce single layer of cells. Straight from snappyLayerDriver
1610
1611 // Global face indices engine
1612 const globalIndex globalFaces(mesh.nFaces());
1613
1614 // Determine extrudePatch.edgeFaces in global numbering (so across
1615 // coupled patches). This is used only to string up edges
1616 // between coupled
1617 // faces (all edges between same (global)face indices get extruded).
1618 labelListList edgeGlobalFaces
1619 (
1620 addPatchCellLayer::globalEdgeFaces
1621 (
1622 mesh,
1623 globalFaces,
1624 ppPtr()
1625 )
1626 );
1627
1628 // Determine patches for extruded boundary edges. Calculates if any
1629 // additional processor patches need to be constructed.
1630
1631 labelList edgePatchID;
1632 labelList edgeZoneID;
1633 boolList edgeFlip;
1634 labelList inflateFaceID;
1635 snappyLayerDriver::determineSidePatches
1636 (
1637 meshRefiner_,
1638 globalFaces,
1639 edgeGlobalFaces,
1640 ppPtr(),
1641
1642 edgePatchID,
1643 edgeZoneID,
1644 edgeFlip,
1645 inflateFaceID
1646 );
1647
1648
1649 // Mesh topo change engine. Insert current mesh.
1650 polyTopoChange meshMod(mesh);
1651
1652 // Layer mesh modifier
1653 addPatchCellLayer addLayer(mesh);
1654
1655 // Extrude very thin layer of cells
1656 pointField extrusion(PatchTools::pointNormals(mesh, ppPtr()));
1657 const tmp<scalarField> thickness
1658 (
1659 wantedThickness
1660 (
1661 ppPtr(),
1662 1e-3 // cellSizeFraction
1663 )
1664 );
1665 extrusion *= thickness;
1666
1667 // Add topo regardless of whether extrudeStatus is extruderemove.
1668 // Not add layer if patchDisp is zero.
1669
1670 const label ppNFaces = ppPtr().size();
1671 const label ppNPoints = ppPtr().nPoints();
1672
1673 addLayer.setRefinement
1674 (
1675 globalFaces,
1676 edgeGlobalFaces,
1677
1678 scalarField(ppNPoints, 1), // expansion ratio
1679 ppPtr(),
1680 bitSet(ppPtr().size()), // no flip
1681
1682 edgePatchID, // boundary patch for extruded boundary edges
1683 edgeZoneID, // zone for extruded edges
1684 edgeFlip,
1685 inflateFaceID,
1686
1687
1688 labelList(0), // exposed patchIDs, not used for adding layers
1689 labelList(ppNFaces, 1),
1690 labelList(ppNPoints, 1),
1691 extrusion, //patchAttraction,
1692
1693 meshMod
1694 );
1695
1696 // Save old meshPoints before changing mesh
1697 const Map<label> oldMeshPointMap(ppPtr->meshPointMap());
1698 Pout<< "old pp points:" << ppPtr->nPoints()
1699 << " oldMeshPointMap:" << oldMeshPointMap.size()
1700 << endl;
1701
1702 // Remove any unnecessary fields
1703 meshMoverPtr.clear();
1704 motionPtr.clear();
1705 ppPtr.clear();
1706 mesh.clearOut();
1707 mesh.moving(false);
1708
1709 // Apply the stored topo changes to the current mesh.
1710 autoPtr<mapPolyMesh> mapPtr = meshMod.changeMesh(mesh, false);
1711 mapPolyMesh& map = *mapPtr;
1712
1713 // Update fields
1714 mesh.updateMesh(map);
1715
1716 // Move mesh (since morphing does not do this)
1717 if (map.hasMotionPoints())
1718 {
1719 mesh.movePoints(map.preMotionPoints());
1720 }
1721 else
1722 {
1723 // Hack to remove meshPhi - mapped incorrectly. TBD.
1724 mesh.clearOut();
1725 }
1726
1727 // Reset the instance for if in overwrite mode
1728 mesh.setInstance(meshRefiner_.timeName());
1729
1730 // Re-do the patch
1731 ppPtr = meshRefinement::makePatch(mesh, adaptPatchIDs);
1732 Pout<< "new pp points:" << ppPtr->nPoints()
1733 << " new meshPointMap:" << ppPtr->meshPointMap().size()
1734 << endl;
1735 // Map the old-patch-point data to the new patch points
1736 // pointMap (new-to-old) for new pp points
1737 labelList ppMap(ppPtr->nPoints(), -1);
1738 {
1739 const labelListList& oldAddedPoints = addLayer.addedPoints();
1740 forAll(oldAddedPoints, oldPatchPointi)
1741 {
1742 const label oldPointi = oldAddedPoints[oldPatchPointi].last();
1743 const label newPointi = map.reversePointMap()[oldPointi];
1744 const label newPatchPointi = ppPtr().meshPointMap()[newPointi];
1745
1746 ppMap[newPatchPointi] = oldPatchPointi;
1747 }
1748 }
1749
1750 // Update attraction
1751 meshRefinement::updateList(ppMap, vector::zero, patchDisp);
1752 // snapDist
1753 meshRefinement::updateList(ppMap, scalar(0), snapDist);
1754 // patchFeaturePoint
1755 meshRefinement::updateList
1756 (
1757 ppMap,
1758 vector::zero,
1759 patchFeaturePoint
1760 );
1761 // patchConstraints
1762 meshRefinement::updateList
1763 (
1764 ppMap,
1765 pointConstraint(),
1766 patchConstraints
1767 );
1768 // ppLocalPoints
1769 //meshRefinement::updateList(ppMap, vector::zero, ppLocalPoints);
1770 ppLocalPoints = pointField(mesh.points(), ppPtr().meshPoints());
1771
1772
1773 // Update numbering on layer
1774 addLayer.updateMesh
1775 (
1776 map,
1777 identity(ppNFaces),
1778 identity(ppNPoints)
1779 );
1780
1781 // Update intersections
1782 const labelListList addedCells(addLayer.addedCells());
1783 bitSet isChangedFace(mesh.nFaces());
1784 for (const labelList& faceToCells : addedCells)
1785 {
1786 for (const label celli : faceToCells)
1787 {
1788 isChangedFace.set(mesh.cells()[celli]);
1789 }
1790 }
1791 meshRefiner_.updateMesh(map, isChangedFace.toc());
1792
1793
1794 if (debug & meshRefinement::MESH)
1795 {
1796 const_cast<Time&>(mesh.time())++;
1797 Info<< "Writing mesh-with-layer to time "
1798 << meshRefiner_.timeName() << endl;
1799 meshRefiner_.write
1800 (
1801 meshRefinement::debugType(debug),
1802 meshRefinement::writeType
1803 (
1804 meshRefinement::writeLevel()
1805 | meshRefinement::WRITEMESH
1806 ),
1807 meshRefiner_.timeName()
1808 );
1809 }
1810 if (debug)
1811 {
1812 OBJstream str
1813 (
1814 mesh.time().path()
1815 / "new_projection"
1816 + meshRefiner_.timeName()
1817 + ".obj"
1818 );
1819 forAll(ppLocalPoints, pointi)
1820 {
1821 const point& pt = ppLocalPoints[pointi];
1822 str.write(linePointRef(pt, patchFeaturePoint[pointi]));
1823 }
1824 Pout<< "** writing mapped attraction to " << str.name() << endl;
1825 }
1826 }
1827*/
1828//XXXX
1829 {
1830 // Layer mesh modifier
1831 addPatchCellLayer addLayer(mesh);
1832
1833 // Save old mesh points (to construct the new-to-old local patch points)
1834 const Map<label> oldMeshPointMap(ppPtr->meshPointMap());
1835
1836 meshMoverPtr.clear();
1837 motionPtr.clear();
1838
1839 const pointField thickness
1840 (
1841 wantedThickness(ppPtr(), 1e-3) //cellSizeFraction
1842 * PatchTools::pointNormals(mesh, ppPtr())
1843 );
1844 autoPtr<mapPolyMesh> mapPtr = addBufferLayers
1845 (
1846 ppPtr(),
1847 thickness,
1848 addLayer
1849 );
1850
1851 // Re-do the patch
1852 ppPtr = meshRefinement::makePatch(mesh, adaptPatchIDs);
1853
1854 // Map the old-patch-point data to the new patch points
1855 // pointMap (new-to-old) for new pp points
1856 labelList ppMap(ppPtr->nPoints(), -1);
1857 {
1858 const labelListList& addedPoints = addLayer.addedPoints();
1859 forAll(addedPoints, oldPatchPointi)
1860 {
1861 const label newPointi = addedPoints[oldPatchPointi].last();
1862 const label newPatchPointi = ppPtr().meshPointMap()[newPointi];
1863 ppMap[newPatchPointi] = oldPatchPointi;
1864 }
1865 }
1866
1867 // Update attraction
1868 meshRefinement::updateList(ppMap, vector::zero, patchDisp);
1869 // snapDist
1870 meshRefinement::updateList(ppMap, scalar(0), snapDist);
1871 // patchFeaturePoint
1873 (
1874 ppMap,
1876 patchFeaturePoint
1877 );
1878 // patchConstraints
1880 (
1881 ppMap,
1883 patchConstraints
1884 );
1885 // ppLocalPoints
1886 ppLocalPoints = pointField(mesh.points(), ppPtr().meshPoints());
1887 }
1888//XXXXXX
1889 const bool snapToGeometry = true;
1890 if (snapToGeometry)
1891 {
1892 // Create pointMesh with the correct patches
1893 // - points on single polyPatches stay as is
1894 // - points on two polyPatches go to allEdgePatchName
1895 // - points on >two polyPatches go to allPointPatchName
1896 //const pointMesh& pMesh = makePointMesh
1897 //(
1898 // ppPtr(),
1899 // patchConstraints,
1900 // allEdgePatchName,
1901 // allPointPatchName
1902 //);
1903 const pointMesh& pMesh = pointMesh::New(mesh);
1904
1906 (
1907 makeMotionSolver
1908 (
1909 pMesh,
1910 snapDict,
1911 adaptPatchIDs
1912 //patchConstraints
1913 )
1914 );
1915
1916 // Insert as bc to motionSmoother. Note that this takes displacement
1917 // relative to points0
1918 setDisplacement
1919 (
1920 ppPtr(),
1921 patchFeaturePoint-ppLocalPoints,
1922 adaptPatchIDs,
1923 motionPtr().points0(),
1924 motionPtr().pointDisplacement()
1925 );
1926
1927
1928 if (debug)
1929 {
1930 OBJstream str
1931 (
1932 mesh.time().path()
1933 / "buffer_layer_new_projection"
1934 + meshRefiner_.timeName()
1935 + ".obj"
1936 );
1937 forAll(ppLocalPoints, pointi)
1938 {
1939 const point& pt = ppLocalPoints[pointi];
1940 str.write(linePointRef(pt, patchFeaturePoint[pointi]));
1941 }
1942 Pout<< "** writing mapped attraction to " << str.name() << endl;
1943 }
1944
1945
1946 // Solve internal displacement
1947 tmp<pointField> tnewPoints(motionPtr->newPoints());
1948
1949 // Move points
1950 mesh.movePoints(tnewPoints);
1951
1952 // Update pp for new mesh points. Ok as long as we also update geometry
1953 // and wanted displacement (usually zero if mesh motion has succeeded)
1954 ppLocalPoints = pointField(mesh.points(), ppPtr().meshPoints());
1955 patchDisp -= (ppLocalPoints-ppPtr().localPoints());
1956
1958 {
1959 const_cast<Time&>(mesh.time())++;
1960 Info<< "Writing smoothed LAYER mesh to time "
1961 << meshRefiner_.timeName() << endl;
1962 meshRefiner_.write
1963 (
1966 (
1969 ),
1970 meshRefiner_.timeName()
1971 );
1972 }
1973 }
1974
1975
1976 // Merge any introduced baffles (from faceZones of faceType 'internal')
1977 {
1978 autoPtr<mapPolyMesh> mapPtr = meshRefiner_.mergeZoneBaffles
1979 (
1980 true, // internal zones
1981 false // baffle zones
1982 );
1983
1984 if (mapPtr.valid())
1985 {
1987 {
1988 const_cast<Time&>(mesh.time())++;
1989 Info<< "Writing baffle-merged mesh to time "
1990 << meshRefiner_.timeName() << endl;
1991 meshRefiner_.write
1992 (
1995 (
1998 ),
1999 meshRefiner_.timeName()
2000 );
2001 }
2002 }
2003 }
2004
2005 // Repatch faces according to nearest. Do not repatch baffle faces.
2006 {
2007 labelList duplicateFace(getInternalOrBaffleDuplicateFace());
2008
2009 repatchToSurface(snapParams, adaptPatchIDs, duplicateFace);
2010 }
2011
2013 {
2014 const_cast<Time&>(mesh.time())++;
2015 }
2016}
2017
2018
2019// ************************************************************************* //
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())
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
GeometricBoundaryField< vector, pointPatchField, pointMesh > Boundary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
@ 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 Addr & addressing() const noexcept
The addressing used for the list.
label size() const noexcept
The number of elements in 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
A HashTable to objects of type <T> with a label key.
Definition Map.H:54
static FOAM_NO_DANGLING_REFERENCE const pointConstraints & New(const pointMesh &mesh, Args &&... args)
const Mesh & mesh() const noexcept
Reference to the mesh.
Definition MeshObject.H:257
An OFstream that keeps track of vertices and provides convenience output methods for OBJ files.
Definition OBJstream.H:58
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.
label nPoints() const
Number of points supporting patch faces.
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
const labelList & meshPoints() const
Return labelList of mesh points in patch.
const Field< point_type > & localPoints() const
Return pointField of points in patch.
const List< face_type > & localFaces() const
Return patch faces addressing into local point list.
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition Switch.H:81
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
A List with indirect addressing. Like IndirectList but does not store addressing.
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 zero
static const Form max
Adds layers of cells to outside (or inside) of polyMesh. Can optionally create stand-alone extruded m...
static labelListList globalEdgeFaces(const polyMesh &, const globalIndex &globalFaces, const indirectPrimitivePatch &pp, const bitSet &orientation)
Per patch edge the pp faces (in global indices) using it.
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
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
static autoPtr< displacementMotionSolver > New(const word &solverTypeName, const polyMesh &, const IOdictionary &, const pointVectorField &pointDisplacement, const pointIOField &points0)
Select constructed from polyMesh, dictionary and components.
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
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 void updateMesh(const mapPolyMesh &mpm)
Update mesh corresponding to the given map.
Definition fvMesh.C:960
virtual void movePoints(const pointField &)
Move points, returns volumes swept by faces in motion.
Definition fvMesh.C:884
void clearOut(const bool isMeshUpdate=false)
Clear all geometry and addressing.
Definition fvMesh.C:227
Calculates a non-overlapping list of offsets based on an input size (eg, number of cells) from differ...
Definition globalIndex.H:77
static labelPairList findDuplicateFacePairs(const polyMesh &)
Helper routine to find all baffles (two boundary faces.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
writeType
Enumeration for what to write. Used as a bit-pattern.
debugType
Enumeration for what to debug. Used as a bit-pattern.
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.
static void updateList(const labelList &newToOld, const T &nullValue, List< T > &elems)
Helper: reorder list according to map.
static writeType writeLevel()
Get/set write level.
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.
Given a displacement moves the mesh by scaling the displacement back until there are no more mesh err...
static const weightedPosition zero
Accumulates point constraints through successive applications of the applyConstraint function.
Application of (multi-)patch point constraints.
Mesh representing a set of points created from polyMesh.
Definition pointMesh.H:49
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
void setInstance(const fileName &instance, const IOobjectOption::writeOption wOpt=IOobject::AUTO_WRITE)
Set the instance for mesh files.
Definition polyMeshIO.C:29
const faceZoneMesh & faceZones() const noexcept
Return face zone mesh.
Definition polyMesh.H:671
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition polyMesh.C:844
bool moving() const noexcept
Is mesh moving.
Definition polyMesh.H:732
virtual const labelList & faceOwner() const
Return face owner.
Definition polyMesh.C:1101
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition polyMesh.C:1107
virtual const pointField & points() const
Return raw points.
Definition polyMesh.C:1063
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh").
Definition polyMesh.H:411
Direct mesh changes based on v1.3 polyTopoChange syntax.
bool isInternalFace(const label faceIndex) const noexcept
Return true if given face label is internal to the mesh.
label nCells() const noexcept
Number of mesh cells.
label nFaces() const noexcept
Number of mesh faces.
const cellList & cells() const
Simple container to keep together snap specific information.
static void determineSidePatches(meshRefinement &meshRefiner, const globalIndex &globalFaces, const labelListList &edgeGlobalFaces, const indirectPrimitivePatch &pp, labelList &edgePatchID, labelList &edgeZoneID, boolList &edgeFlip, labelList &inflateFaceID)
Helper: see what zones and patches edges should be extruded into.
static bitSet getMasterPoints(const polyMesh &mesh)
Get per point whether it is uncoupled or a master of a coupled set of points.
Definition syncTools.C:61
static void syncFaceList(const polyMesh &mesh, UList< T > &faceValues, const CombineOp &cop, const bool parRun=UPstream::parRun())
Synchronize values on all mesh faces.
Definition syncTools.H:465
static void syncPointList(const polyMesh &mesh, List< T > &pointValues, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronize values on all mesh points.
Decomposes polyMesh into tets (or pyramids).
static bool test(const UList< face > &faces)
Test if given list of faces satisfies criteria for TET. (4 tri).
Definition tetMatcher.C:82
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
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
dynamicFvMesh & mesh
Dummy transform to be used with syncTools.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
const pointField & points
label nPoints
#define WarningInFunction
Report a warning using Foam::Warning.
const dimensionedScalar mp
Proton mass.
Namespace for handling debugging switches.
Definition debug.C:45
const std::string patch
OpenFOAM patch number as a std::string.
List< edge > edgeList
List of edge.
Definition edgeList.H:32
Pair< label > labelPair
A pair of labels.
Definition Pair.H:54
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
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
bool mkDir(const fileName &pathName, mode_t mode=0777)
Make a directory and return an error if it could not be created.
Definition POSIX.C:616
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition HashSet.H:85
constexpr label labelMin
Definition label.H:54
messageStream Info
Information stream (stdout output on master, null elsewhere).
GeometricField< vector, pointPatchField, pointMesh > pointVectorField
vectorIOField pointIOField
pointIOField is a vectorIOField.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
constexpr label labelMax
Definition label.H:55
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
constexpr scalar degToRad() noexcept
Multiplication factor for degrees to radians conversion.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
line< point, const point & > linePointRef
A line using referred points.
Definition line.H:66
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.
vector point
Point is a vector.
Definition point.H:37
List< bool > boolList
A List of bools.
Definition List.H:60
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...
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
dimensionedScalar cos(const dimensionedScalar &ds)
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
labelList f(nPoints)
loopControl iters(runTime, aMesh.solutionDict(), "solution")
volScalarField & e
dimensionedScalar lambda("lambda", dimTime/sqr(dimLength), laminarTransport)
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Unit conversion functions.
pointField points0(pointIOField(IOobject("points", mesh.time().constant(), polyMesh::meshSubDir, mesh, IOobject::MUST_READ, IOobject::NO_WRITE, IOobject::NO_REGISTER)))