Loading...
Searching...
No Matches
snappySnapDriver.H
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2015 OpenFOAM Foundation
9 Copyright (C) 2020 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27Class
28 Foam::snappySnapDriver
29
30Description
31 All to do with snapping to surface
32
33SourceFiles
34 snappySnapDriver.C
35 snappySnapDriverFeature.C
36
37\*---------------------------------------------------------------------------*/
38
39#ifndef snappySnapDriver_H
40#define snappySnapDriver_H
41
42#include "meshRefinement.H"
43#include "DynamicField.H"
44#include "pointConstraint.H"
45
46// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47
48namespace Foam
49{
50
51// Forward declaration of classes
52class motionSmoother;
54class snapParameters;
55class pointConstraint;
56class layerParameters;
60/*---------------------------------------------------------------------------*\
61 Class snappySnapDriver Declaration
62\*---------------------------------------------------------------------------*/
63
64class snappySnapDriver
65{
66 // Private data
67
68 //- Mesh+surface
69 meshRefinement& meshRefiner_;
70
71 //- From global surface region to master side patch
72 const labelList globalToMasterPatch_;
73
74 //- From global surface region to slave side patch
75 const labelList globalToSlavePatch_;
76
77 //- Are we operating in test mode?
78 const bool dryRun_;
79
80
81 // Private Member Functions
82
83
84 // Snapping
85
86 //- Top-level: snap onto surface & features and add buffer layer
87 void doSnapBufferLayers
88 (
89 const dictionary& snapDict,
90 const dictionary& motionDict,
91 const meshRefinement::FaceMergeType mergeType,
92 const scalar featureCos,
93 const scalar planarAngle,
94 const snapParameters& snapParams
95 );
96
97 //- Calculates (geometric) shared points
98 // Requires bitSet to be sized and initialised
99 static label getCollocatedPoints
100 (
101 const scalar tol,
102 const pointField&,
103 bitSet&
104 );
105
106 //- Calculate displacement (over all mesh points) to move points
107 // to average of connected cell centres
108 static tmp<pointField> smoothInternalDisplacement
109 (
110 const meshRefinement& meshRefiner,
111 const motionSmoother&
112 );
113
114 //- Calculate displacement per patch point to smooth out patch.
115 // Quite complicated in determining which points to move where.
116 static tmp<pointField> smoothPatchDisplacement
117 (
118 const motionSmoother&,
119 const List<labelPair>&
120 );
121
122 static tmp<pointField> avg
123 (
124 const polyMesh& mesh,
125 const bitSet& isMasterPoint,
127 const pointField&
128 );
129
130 //- Calculate displacement per patch point. Wip.
131 static tmp<pointField> smoothLambdaMuPatchDisplacement
132 (
133 const polyMesh& mesh,
135 const List<labelPair>& baffles
136 );
137
138 tmp<scalarField> wantedThickness
139 (
141 const scalar cellSizeFraction
142 ) const;
143
144 //- Create/update pointMesh
145 const pointMesh& makePointMesh
146 (
149 const word& allEdgePatchName,
150 const word& allPointPatchName
151 ) const;
152
153 autoPtr<displacementMotionSolver> makeMotionSolver
154 (
155 const pointMesh& pMesh,
156 const dictionary& snapDict,
157 const labelList& adaptPatchIDs
158 ) const;
159
160 void setDisplacement
161 (
163 const pointField& patchDisp,
164 const labelList& adaptPatchIDs,
165 const pointField& points0,
167 );
168
169 autoPtr<mapPolyMesh> addBufferLayers
170 (
172
173 //const scalar cellSizeFraction,
174 const pointField& thickness,
175 // Layer mesh modifier
176 addPatchCellLayer& addLayer
177 );
178
179 //- Check that face zones are synced
180 void checkCoupledFaceZones() const;
181
182 //- Per edge distance to patch
183 static tmp<scalarField> edgePatchDist
184 (
185 const pointMesh&,
187 );
188
189 //- Write displacement as .obj file.
190 static void dumpMove
191 (
192 const fileName&,
193 const pointField&,
194 const pointField&
195 );
196
197 //- Check displacement is outwards pointing
198 static bool outwardsDisplacement
199 (
201 const vectorField&
202 );
203
204 //- Freeze points on pointZone or (inside of) faceZone
205 static void freezeExposedPoints
206 (
207 const meshRefinement& meshRefiner,
208 const word& fzName, // faceZone name
209 const word& pzName, // pointZone name
210 const indirectPrimitivePatch& outside,
211 vectorField& patchDisp
212 );
213
214 //- Detect warpage
215 void detectWarpedFaces
216 (
217 const scalar featureCos,
219
220 DynamicList<label>& splitFaces,
222 ) const;
223
224 //- Get per face -1 or label of opposite face if on internal/baffle
225 // faceZone
226 labelList getInternalOrBaffleDuplicateFace() const;
227
228 //- Get points both on patch and facezone.
229 static void getZoneSurfacePoints
230 (
231 const fvMesh& mesh,
233 const word& zoneName,
234
235 bitSet& pointOnZone
236 );
237
238 //- Get points both on patch and facezone.
239 template<class FaceList>
240 static labelList getFacePoints
241 (
243 const FaceList& faces
244 );
245
246 //- Per patch point calculate point on nearest surface.
247 // Return displacement of patch points.
248 static void calcNearestSurface
249 (
250 const refinementSurfaces& surfaces,
251
252 const labelList& surfacesToTest,
253 const labelListList& regionsToTest,
254
255 const pointField& localPoints,
256 const labelList& zonePointIndices,
257
258 scalarField& minSnapDist,
259 labelList& snapSurf,
260 vectorField& patchDisp,
261
262 // Optional: nearest point, normal
263 pointField& nearestPoint,
264 vectorField& nearestNormal
265 );
266
267
268 // Feature line snapping
269
270 //- Is point on two feature edges that make a largish angle?
271 bool isFeaturePoint
272 (
273 const scalar featureCos,
275 const bitSet& isFeatureEdge,
276 const label pointi
277 ) const;
278
279 void smoothAndConstrain
280 (
281 const bitSet& isMasterEdge,
283 const labelList& meshEdges,
284 const List<pointConstraint>& constraints,
285 vectorField& disp
286 ) const;
287 //void smoothAndConstrain2
288 //(
289 // const bool applyConstraints,
290 // const indirectPrimitivePatch& pp,
291 // const List<pointConstraint>& constraints,
292 // vectorField& disp
293 //) const;
294 void calcNearest
295 (
296 const label iter,
298 vectorField& pointDisp,
299 vectorField& pointSurfaceNormal,
300 vectorField& pointRotation
301 ) const;
302 void calcNearestFace
303 (
304 const label iter,
306 const pointField& ppLocalPoints,
307 const scalarField& faceSnapDist,
308 vectorField& faceDisp,
309 vectorField& faceSurfaceNormal,
310 labelList& faceSurfaceRegion
311 //vectorField& faceRotation
312 ) const;
313 void calcNearestFacePointProperties
314 (
315 const label iter,
317 const pointField& ppLocalPoints,
318
319 const vectorField& faceDisp,
320 const vectorField& faceSurfaceNormal,
321 const labelList& faceSurfaceRegion,
322
323 List<List<point>>& pointFaceSurfNormals,
324 List<List<point>>& pointFaceDisp,
325 List<List<point>>& pointFaceCentres,
326 List<labelList>& pointFacePatchID
327 ) const;
328 void correctAttraction
329 (
330 const UList<point>& surfacePoints,
331 const UList<label>& surfaceCounts,
332 const point& edgePt,
333 const vector& edgeNormal, // normalised normal
334 const point& pt,
335 vector& edgeOffset // offset from pt to point on edge
336 ) const;
337
338
339 //- For any reverse (so from feature back to mesh) attraction:
340 // add attraction if diagonal points on face attracted
341 void stringFeatureEdges
342 (
343 const label iter,
344 const scalar featureCos,
345
347 const pointField& ppLocalPoints,
348 const scalarField& snapDist,
349
350 const vectorField& rawPatchAttraction,
351 const List<pointConstraint>& rawPatchConstraints,
352
353 vectorField& patchAttraction,
354 List<pointConstraint>& patchConstraints
355 ) const;
356
357 //- Remove constraints of points next to multi-patch points
358 // to give a bit more freedom of the mesh to conform to the
359 // multi-patch points. Bit dodgy for simple cases.
360 void releasePointsNextToMultiPatch
361 (
362 const label iter,
363 const scalar featureCos,
364
366 const pointField& ppLocalPoints,
367 const scalarField& snapDist,
368
369 const List<List<point>>& pointFaceCentres,
370 const labelListList& pointFacePatchID,
371
372 const vectorField& rawPatchAttraction,
373 const List<pointConstraint>& rawPatchConstraints,
374
375 vectorField& patchAttraction,
376 List<pointConstraint>& patchConstraints
377 ) const;
378
379 //- Detect any diagonal attraction. Returns indices in face
380 // or (-1, -1) if none
381 labelPair findDiagonalAttraction
382 (
384 const vectorField& patchAttraction,
385 const List<pointConstraint>& patchConstraints,
386 const label facei
387 ) const;
388
389 scalar pyrVol
390 (
392 const vectorField& featureAttraction,
393 const face& localF,
394 const point& cc
395 ) const;
396 void facePoints
397 (
399 const vectorField& featureAttraction,
400 const vectorField& nearestAttraction,
401 const face& f,
403 ) const;
404 scalar pyrVol
405 (
407 const vectorField& featureAttraction,
408 const vectorField& nearestAttraction,
409 const face& localF,
410 const point& cc
411 ) const;
412 Tuple2<point, vector> centreAndNormal
413 (
415 const vectorField& featureAttraction,
416 const vectorField& nearestAttraction,
417 const face& localF
418 ) const;
419 bool isSplitAlignedWithFeature
420 (
421 const scalar featureCos,
422 const point& newPt0,
423 const pointConstraint& pc0,
424 const point& newPt1,
425 const pointConstraint& pc1
426 ) const;
427 bool isConcave
428 (
429 const point& c0,
430 const vector& area0,
431 const point& c1,
432 const vector& area1,
433 const scalar concaveCos
434 ) const;
435 labelPair findDiagonalAttraction
436 (
437 const scalar featureCos,
438 const scalar concaveCos,
439 const scalar minAreaFraction,
441 const pointField& ppLocalPoints,
442 const vectorField& patchAttraction,
443 const List<pointConstraint>& patchConstraints,
444 const vectorField& nearestAttraction,
445 const vectorField& nearestNormal,
446 const label faceI,
447
449 DynamicField<point>& points1
450 ) const;
451
452 //- Do all logic on whether to add face cut to diagonal
453 // attraction
454 void splitDiagonals
455 (
456 const scalar featureCos,
457 const scalar concaveCos,
458 const scalar minAreaFraction,
459
461 const pointField& ppLocalPoints,
462 const vectorField& nearestAttraction,
463 const vectorField& nearestNormal,
464 const List<labelList>& pointFacePatchID,
465
466 vectorField& patchAttraction,
467 List<pointConstraint>& patchConstraints,
468 DynamicList<label>& splitFaces,
470 DynamicList<labelPair>& splitPatches
471 ) const;
472
473 //- Avoid attraction across face diagonal since would
474 // cause face squeeze
475 void avoidDiagonalAttraction
476 (
477 const label iter,
478 const scalar featureCos,
480 const pointField& ppLocalPoints,
481 vectorField& patchAttraction,
482 List<pointConstraint>& patchConstraints
483 ) const;
484
485 //- Write some stats about constraints
486 void writeStats
487 (
489 const bitSet& isMasterPoint,
490 const List<pointConstraint>& patchConstraints
491 ) const;
492
493 //- Return hit if on multiple points
494 pointIndexHit findMultiPatchPoint
495 (
496 const point& pt,
497 const labelList& patchIDs,
498 const List<point>& faceCentres
499 ) const;
500
501 //- Return hit if faces-on-the-same-normalplane are on multiple
502 // patches
503 pointIndexHit findMultiPatchPoint
504 (
505 const point& pt,
506 const labelList& pfPatchID,
507 const UList<vector>& surfaceNormals,
508 const labelList& faceToNormalBin
509 ) const;
510
511 //- Return index of similar normal
512 label findNormal
513 (
514 const scalar featureCos,
515 const vector& faceSurfaceNormal,
516 const UList<vector>& surfaceNormals
517 ) const;
518
519 //- Determine attraction and constraints for single point
520 // using sampled surrounding of the point
521 void featureAttractionUsingReconstruction
522 (
523 const label iter,
524 const scalar featureCos,
525
527 const pointField& ppLocalPoints,
528 const scalarField& snapDist,
529 const vectorField& nearestDisp,
530 const label pointi,
531
532 const List<List<point>>& pointFaceSurfNormals,
533 const List<List<point>>& pointFaceDisp,
534 const List<List<point>>& pointFaceCentres,
535 const labelListList& pointFacePatchID,
536
537 DynamicList<point>& surfacePoints,
538 DynamicList<vector>& surfaceNormals,
539 labelList& faceToNormalBin,
540
541 vector& patchAttraction,
542 pointConstraint& patchConstraint
543 ) const;
544
545 //- Determine attraction and constraints for all points
546 // using sampled surrounding of the point
547 void featureAttractionUsingReconstruction
548 (
549 const label iter,
550 const scalar featureCos,
552 const pointField& ppLocalPoints,
553 const scalarField& snapDist,
554 const vectorField& nearestDisp,
555
556 const List<List<point>>& pointFaceSurfNormals,
557 const List<List<point>>& pointFaceDisp,
558 const List<List<point>>& pointFaceCentres,
559 const labelListList& pointFacePatchID,
560
561 vectorField& patchAttraction,
562 List<pointConstraint>& patchConstraints
563 ) const;
564
565 //- Determine geometric features and attraction to equivalent
566 // surface features
567 void determineFeatures
568 (
569 const label iter,
570 const scalar featureCos,
571 const bool multiRegionFeatureSnap,
572 const bool strictRegionFeatureSnap, // special feat-point
573
575 const pointField& ppLocalPoints,
576 const scalarField& snapDist,
577 const vectorField& nearestDisp,
578
579 const List<List<point>>& pointFaceSurfNormals,
580 const List<List<point>>& pointFaceDisp,
581 const List<List<point>>& pointFaceCentres,
582 const labelListList& pointFacePatchID,
583
584 List<labelList>& pointAttractor,
586 // Feature-edge to pp point
587 List<List<DynamicList<point>>>& edgeAttractors,
588 List<List<DynamicList<pointConstraint>>>& edgeConstraints,
589 vectorField& patchAttraction,
590 List<pointConstraint>& patchConstraints
591 ) const;
592
593 //- Determine features originating from bafles and
594 // and add attraction to equivalent surface features
595 void determineBaffleFeatures
596 (
597 const label iter,
598 const bool baffleFeaturePoints,
599 const scalar featureCos,
600
602 const pointField& ppLocalPoints,
603 const scalarField& snapDist,
604
605 // Feature-point to pp point
606 List<labelList>& pointAttractor,
608 // Feature-edge to pp point
609 List<List<DynamicList<point>>>& edgeAttractors,
610 List<List<DynamicList<pointConstraint>>>& edgeConstraints,
611 // pp point to nearest feature
612 vectorField& patchAttraction,
613 List<pointConstraint>& patchConstraints
614 ) const;
615 void reverseAttractMeshPoints
616 (
617 const label iter,
618
620 const pointField& ppLocalPoints,
621 const scalarField& snapDist,
622
623 // Feature-point to pp point
624 const List<labelList>& pointAttractor,
626 // Feature-edge to pp point
627 const List<List<DynamicList<point>>>& edgeAttractors,
629
630 const vectorField& rawPatchAttraction,
631 const List<pointConstraint>& rawPatchConstraints,
632
633 // pp point to nearest feature
634 vectorField& patchAttraction,
635 List<pointConstraint>& patchConstraints
636 ) const;
637
638 //- Find point on nearest feature edge (within searchDist).
639 // Return point and feature
640 // and store feature-edge to mesh-point and vice versa
641 Tuple2<label, pointIndexHit> findNearFeatureEdge
642 (
643 const bool isRegionEdge,
644
646 const pointField& ppLocalPoints,
647 const scalarField& snapDist,
648 const label pointi,
649 const point& estimatedPt,
650
655 ) const;
656
657 //- Find nearest feature point (within searchDist).
658 // Return feature point
659 // and store feature-point to mesh-point and vice versa.
660 // If another mesh point already referring to this feature
661 // point and further away, reset that one to a near feature
662 // edge (using findNearFeatureEdge above)
663 Tuple2<label, pointIndexHit> findNearFeaturePoint
664 (
665 const bool isRegionEdge,
666
668 const pointField& ppLocalPoints,
669 const scalarField& snapDist,
670 const label pointi,
671 const point& estimatedPt,
672
673 // Feature-point to pp point
674 List<labelList>& pointAttractor,
676 // Feature-edge to pp point
677 List<List<DynamicList<point>>>& edgeAttractors,
678 List<List<DynamicList<pointConstraint>>>& edgeConstraints,
679 // pp point to nearest feature
680 vectorField& patchAttraction,
681 List<pointConstraint>& patchConstraints
682 ) const;
683
684 void featureAttractionUsingFeatureEdges
685 (
686 const label iter,
687 const bool multiRegionFeatureSnap,
688 const bool strictRegionFeatureSnap,
689
690 const bool detectBaffles,
691 const bool baffleFeaturePoints,
692 const bool releasePoints,
693 const bool stringFeatures,
694 const bool avoidDiagonal,
695
696 const scalar featureCos,
697
699 const pointField& ppLocalPoints,
700 const scalarField& snapDist,
701 const vectorField& nearestDisp,
702 const vectorField& nearestNormal,
703
704 const List<List<point>>& pointFaceSurfNormals,
705 const List<List<point>>& pointFaceDisp,
706 const List<List<point>>& pointFaceCentres,
707 const labelListList& pointFacePatchID,
708
709 vectorField& patchAttraction,
710 List<pointConstraint>& patchConstraints
711 ) const;
712
713 void preventFaceSqueeze
714 (
715 const label iter,
716 const scalar featureCos,
718 const pointField& ppLocalPoints,
719 const scalarField& snapDist,
720 const vectorField& nearestAttraction,
721
722 vectorField& patchAttraction,
723 List<pointConstraint>& patchConstraints
724 ) const;
725
726 //- Top level feature attraction routine. Gets given
727 // displacement to nearest surface in nearestDisp
728 // and calculates new displacement taking into account
729 // features
730 vectorField calcNearestSurfaceFeature
731 (
732 const snapParameters& snapParams,
733 const bool alignMeshEdges,
734 const bool strictRegionFeatureSnap, //points on >=3 patches
735 const label iter,
736 const scalar featureCos,
737 const scalar featureAttract,
738 const scalarField& snapDist,
739 const vectorField& nearestDisp,
740 const vectorField& nearestNormal,
742 const pointField& ppLocalPoints,
743 vectorField& patchAttraction,
744 List<pointConstraint>& patchConstraints,
745
746 DynamicList<label>& splitFaces,
748 DynamicList<labelPair>& splitPatches
749 ) const;
750
751
752 //- No copy construct
753 snappySnapDriver(const snappySnapDriver&) = delete;
754
755 //- No copy assignment
756 void operator=(const snappySnapDriver&) = delete;
757
758
759public:
760
761 //- Runtime type information
762 ClassName("snappySnapDriver");
763
764
765 // Constructors
766
767 //- Construct from components
768 snappySnapDriver
769 (
770 meshRefinement& meshRefiner,
771 const labelList& globalToMasterPatch,
772 const labelList& globalToSlavePatch,
773 const bool dryRun = false
774 );
775
776
777 // Member Functions
778
779 // Snapping
780
781 //- Merge baffles.
783
784 //- Calculate edge length per patch point.
786 (
787 const fvMesh& mesh,
788 const snapParameters& snapParams,
790 );
791
792 //- Smooth the mesh (patch and internal) to increase visibility
793 // of surface points (on castellated mesh) w.r.t. surface.
794 static void preSmoothPatch
795 (
796 const meshRefinement& meshRefiner,
797 const snapParameters& snapParams,
798 const label nInitErrors,
799 const List<labelPair>& baffles,
801 );
802
803 //- Helper: calculate average cell centre per point
805 (
806 const fvMesh& mesh,
808 );
809
810 //- Per patch point override displacement if in gap situation
812 (
813 const scalar planarCos,
815 const pointField& ppLocalPoints,
816 const pointField& nearestPoint,
817 const vectorField& nearestNormal,
818 vectorField& disp
819 ) const;
820
821 //- Per patch point calculate point on nearest surface. Set as
822 // boundary conditions of motionSmoother displacement field. Return
823 // displacement of patch points.
824 static vectorField calcNearestSurface
825 (
826 const bool strictRegionSnap,
827 const meshRefinement& meshRefiner,
828 const labelList& globalToMasterPatch,
829 const labelList& globalToSlavePatch,
831 const pointField& ppLocalPoints,
832 const scalarField& snapDist,
833 pointField& nearestPoint,
834 vectorField& nearestNormal
835 );
836
840 //static vectorField calcNearestLocalSurface
841 //(
842 // const meshRefinement& meshRefiner,
843 // const scalarField& snapDist,
844 // const indirectPrimitivePatch&
845 //);
846
847 //- Smooth the displacement field to the internal.
849 (
850 const snapParameters& snapParams,
852 ) const;
853
854 //- Do the hard work: move the mesh according to displacement,
855 // locally relax the displacement. Return true if ended up with
856 // correct mesh, false if not.
857 bool scaleMesh
858 (
859 const snapParameters& snapParams,
860 const label nInitErrors,
861 const List<labelPair>& baffles,
863 );
864
865 //- Repatch faces according to surface nearest the face centre
866 // - calculate face-wise snap distance as max of point-wise
867 // - calculate face-wise nearest surface point
868 // - repatch face according to patch for surface point.
870 (
871 const snapParameters& snapParams,
872 const labelList& adaptPatchIDs,
873 const labelList& preserveFaces
874 );
875
876 //- Snap onto surface & features
877 void doSnap
878 (
879 const dictionary& snapDict,
880 const dictionary& motionDict,
881 const meshRefinement::FaceMergeType mergeType,
882 const scalar featureCos,
883 const scalar planarAngle,
884 const snapParameters& snapParams
885 );
886};
887
888
889// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
890
891} // End namespace Foam
892
893// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
894
895#ifdef NoRepository
897#endif
898
899// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
900
901#endif
902
903// ************************************************************************* //
Info<< nl;Info<< "Write faMesh in vtk format:"<< nl;{ vtk::uindirectPatchWriter writer(aMesh.patch(), fileName(aMesh.time().globalPath()/vtkBaseFileName));writer.writeGeometry();globalIndex procAddr(aMesh.nFaces());labelList cellIDs;if(UPstream::master()) { cellIDs.resize(procAddr.totalSize());for(const labelRange &range :procAddr.ranges()) { auto slice=cellIDs.slice(range);slice=identity(range);} } writer.beginCellData(4);writer.writeProcIDs();writer.write("cellID", cellIDs);writer.write("area", aMesh.S().field());writer.write("normal", aMesh.faceAreaNormals());writer.beginPointData(1);writer.write("normal", aMesh.pointAreaNormals());Info<< " "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edges")));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
labelList patchIDs
Dynamically sized Field. Similar to DynamicList, but inheriting from a Field instead of a List.
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
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 2-tuple for storing two objects of dissimilar types. The container is similar in purpose to std::pa...
Definition Tuple2.H:51
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition UList.H:89
Adds layers of cells to outside (or inside) of polyMesh. Can optionally create stand-alone extruded m...
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition bitSet.H:61
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Virtual base class for displacement motion solver.
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
Simple container to keep together layer specific information.
Helper class which maintains intersections of (changing) mesh with (static) surfaces.
FaceMergeType
Enumeration for what to do with co-planar patch faces on a single.
Given a displacement moves the mesh by scaling the displacement back until there are no more mesh err...
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
Simple container to keep together refinement specific information.
Container for data on surfaces used for surface-driven refinement. Contains all the data about the le...
Simple container to keep together snap specific information.
void detectNearSurfaces(const scalar planarCos, const indirectPrimitivePatch &pp, const pointField &ppLocalPoints, const pointField &nearestPoint, const vectorField &nearestNormal, vectorField &disp) const
Per patch point override displacement if in gap situation.
autoPtr< mapPolyMesh > mergeZoneBaffles(const List< labelPair > &)
Merge baffles.
bool scaleMesh(const snapParameters &snapParams, const label nInitErrors, const List< labelPair > &baffles, motionSmoother &)
Do the hard work: move the mesh according to displacement,.
ClassName("snappySnapDriver")
Runtime type information.
static void preSmoothPatch(const meshRefinement &meshRefiner, const snapParameters &snapParams, const label nInitErrors, const List< labelPair > &baffles, motionSmoother &)
Smooth the mesh (patch and internal) to increase visibility.
static tmp< pointField > avgCellCentres(const fvMesh &mesh, const indirectPrimitivePatch &)
Helper: calculate average cell centre per point.
void doSnap(const dictionary &snapDict, const dictionary &motionDict, const meshRefinement::FaceMergeType mergeType, const scalar featureCos, const scalar planarAngle, const snapParameters &snapParams)
Snap onto surface & features.
void smoothDisplacement(const snapParameters &snapParams, motionSmoother &) const
Smooth the displacement field to the internal.
autoPtr< mapPolyMesh > repatchToSurface(const snapParameters &snapParams, const labelList &adaptPatchIDs, const labelList &preserveFaces)
Repatch faces according to surface nearest the face centre.
static scalarField calcSnapDistance(const fvMesh &mesh, const snapParameters &snapParams, const indirectPrimitivePatch &)
Calculate edge length per patch point.
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
#define ClassName(TypeNameString)
Add typeName information from argument TypeNameString to a class.
Definition className.H:74
dynamicFvMesh & mesh
const pointField & points
Namespace for OpenFOAM.
Pair< label > labelPair
A pair of labels.
Definition Pair.H:54
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
List< label > labelList
A List of labels.
Definition List.H:62
List< pointConstraint > pointConstraintList
List of pointConstraint.
GeometricField< vector, pointPatchField, pointMesh > pointVectorField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
PrimitivePatch< IndirectList< face >, const pointField & > indirectPrimitivePatch
A PrimitivePatch with an IndirectList for the faces, const reference for the point field.
Field< vector > vectorField
Specialisation of Field<T> for vector.
vector point
Point is a vector.
Definition point.H:37
vectorField pointField
pointField is a vectorField.
Vector< scalar > vector
Definition vector.H:57
PointIndexHit< point > pointIndexHit
A PointIndexHit with a 3D point.
labelList f(nPoints)
pointField points0(pointIOField(IOobject("points", mesh.time().constant(), polyMesh::meshSubDir, mesh, IOobject::MUST_READ, IOobject::NO_WRITE, IOobject::NO_REGISTER)))