autoSnapDriver.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 | Copyright (C) 2011-2015 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Class
25  Foam::autoSnapDriver
26 
27 Description
28  All to do with snapping to surface
29 
30 SourceFiles
31  autoSnapDriver.C
32  autoSnapDriverFeature.C
33 
34 \*---------------------------------------------------------------------------*/
35 
36 #ifndef autoSnapDriver_H
37 #define autoSnapDriver_H
38 
39 #include "meshRefinement.H"
40 #include "DynamicField.H"
41 
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 
44 namespace Foam
45 {
46 
47 // Forward declaration of classes
48 class motionSmoother;
49 class refinementParameters;
50 class snapParameters;
51 class pointConstraint;
52 
53 /*---------------------------------------------------------------------------*\
54  Class autoSnapDriver Declaration
55 \*---------------------------------------------------------------------------*/
56 
57 class autoSnapDriver
58 {
59  // Private data
60 
61  //- Mesh+surface
63 
64  //- From global surface region to master side patch
66 
67  //- From global surface region to slave side patch
69 
70 
71  // Private Member Functions
72 
73 
74  // Snapping
75 
76  //- Calculates (geometric) shared points
78  (
79  const scalar tol,
80  const pointField&,
82  );
83 
84  //- Calculate displacement (over all mesh points) to move points
85  // to average of connected cell centres
87  (
88  const meshRefinement& meshRefiner,
89  const motionSmoother&
90  );
91 
92  //- Calculate displacement per patch point to smooth out patch.
93  // Quite complicated in determining which points to move where.
95  (
96  const motionSmoother&,
97  const List<labelPair>&
98  );
99 
100  static tmp<pointField> avg
101  (
102  const indirectPrimitivePatch&,
103  const pointField&
104  );
105 
106  //- Calculate displacement per patch point. Wip.
108  (
109  const motionSmoother& meshMover,
110  const List<labelPair>& baffles
111  );
112 
113 
114  //- Check that face zones are synced
115  void checkCoupledFaceZones() const;
116 
117  //- Per edge distance to patch
119  (
120  const pointMesh&,
122  );
123 
124  //- Write displacement as .obj file.
125  static void dumpMove
126  (
127  const fileName&,
128  const pointField&,
129  const pointField&
130  );
131 
132  //- Check displacement is outwards pointing
133  static bool outwardsDisplacement
134  (
135  const indirectPrimitivePatch&,
136  const vectorField&
137  );
138 
139  //- Detect warpage
140  void detectWarpedFaces
141  (
142  const scalar featureCos,
143  const indirectPrimitivePatch& pp,
144 
145  DynamicList<label>& splitFaces,
146  DynamicList<labelPair>& splits
147  ) const;
148 
149  //- Get per face -1 or label of opposite face if on internal/baffle
150  // faceZone
152 
153  //- Get points both on patch and facezone.
155  (
156  const fvMesh& mesh,
157  const indirectPrimitivePatch&,
158  const word& zoneName
159  );
160 
161  //- Get points both on patch and facezone.
162  template<class FaceList>
163  static labelList getFacePoints
164  (
165  const indirectPrimitivePatch& pp,
166  const FaceList& faces
167  );
168 
169  //- Per patch point calculate point on nearest surface.
170  // Return displacement of patch points.
171  static void calcNearestSurface
172  (
173  const refinementSurfaces& surfaces,
174 
175  const labelList& surfacesToTest,
176  const labelListList& regionsToTest,
177 
178  const pointField& localPoints,
179  const labelList& zonePointIndices,
180 
181  scalarField& minSnapDist,
182  labelList& snapSurf,
183  vectorField& patchDisp,
184 
185  // Optional: nearest point, normal
186  pointField& nearestPoint,
187  vectorField& nearestNormal
188  );
189 
190 
191  // Feature line snapping
192 
193  //- Is point on two feature edges that make a largish angle?
194  bool isFeaturePoint
195  (
196  const scalar featureCos,
197  const indirectPrimitivePatch& pp,
198  const PackedBoolList& isFeatureEdge,
199  const label pointI
200  ) const;
201 
202  void smoothAndConstrain
203  (
204  const PackedBoolList& isMasterEdge,
205  const indirectPrimitivePatch& pp,
206  const labelList& meshEdges,
207  const List<pointConstraint>& constraints,
208  vectorField& disp
209  ) const;
210  //void smoothAndConstrain2
211  //(
212  // const bool applyConstraints,
213  // const indirectPrimitivePatch& pp,
214  // const List<pointConstraint>& constraints,
215  // vectorField& disp
216  //) const;
217  void calcNearest
218  (
219  const label iter,
220  const indirectPrimitivePatch& pp,
221  vectorField& pointDisp,
222  vectorField& pointSurfaceNormal,
223  vectorField& pointRotation
224  ) const;
225  void calcNearestFace
226  (
227  const label iter,
228  const indirectPrimitivePatch& pp,
229  const scalarField& faceSnapDist,
230  vectorField& faceDisp,
231  vectorField& faceSurfaceNormal,
232  labelList& faceSurfaceRegion
233  //vectorField& faceRotation
234  ) const;
236  (
237  const label iter,
238  const indirectPrimitivePatch& pp,
239 
240  const vectorField& faceDisp,
241  const vectorField& faceSurfaceNormal,
242  const labelList& faceSurfaceRegion,
243 
244  List<List<point> >& pointFaceSurfNormals,
245  List<List<point> >& pointFaceDisp,
246  List<List<point> >& pointFaceCentres,
247  List<labelList>& pointFacePatchID
248  ) const;
249  void correctAttraction
250  (
251  const DynamicList<point>& surfacePoints,
252  const DynamicList<label>& surfaceCounts,
253  const point& edgePt,
254  const vector& edgeNormal, // normalised normal
255  const point& pt,
256  vector& edgeOffset // offset from pt to point on edge
257  ) const;
258 
259 
260  //- For any reverse (so from feature back to mesh) attraction:
261  // add attraction if diagonal points on face attracted
262  void stringFeatureEdges
263  (
264  const label iter,
265  const scalar featureCos,
266 
267  const indirectPrimitivePatch& pp,
268  const scalarField& snapDist,
269 
270  const vectorField& rawPatchAttraction,
271  const List<pointConstraint>& rawPatchConstraints,
272 
273  vectorField& patchAttraction,
274  List<pointConstraint>& patchConstraints
275  ) const;
276 
277  //- Remove constraints of points next to multi-patch points
278  // to give a bit more freedom of the mesh to conform to the
279  // multi-patch points. Bit dodgy for simple cases.
281  (
282  const label iter,
283  const scalar featureCos,
284 
285  const indirectPrimitivePatch& pp,
286  const scalarField& snapDist,
287 
288  const List<List<point> >& pointFaceCentres,
289  const labelListList& pointFacePatchID,
290 
291  const vectorField& rawPatchAttraction,
292  const List<pointConstraint>& rawPatchConstraints,
293 
294  vectorField& patchAttraction,
295  List<pointConstraint>& patchConstraints
296  ) const;
297 
298  //- Detect any diagonal attraction. Returns indices in face
299  // or (-1, -1) if none
301  (
302  const indirectPrimitivePatch& pp,
303  const vectorField& patchAttraction,
304  const List<pointConstraint>& patchConstraints,
305  const label faceI
306  ) const;
307 
308  scalar pyrVol
309  (
310  const indirectPrimitivePatch& pp,
311  const vectorField& featureAttraction,
312  const face& localF,
313  const point& cc
314  ) const;
315  void facePoints
316  (
317  const indirectPrimitivePatch& pp,
318  const vectorField& featureAttraction,
319  const vectorField& nearestAttraction,
320  const face& f,
322  ) const;
323  scalar pyrVol
324  (
325  const indirectPrimitivePatch& pp,
326  const vectorField& featureAttraction,
327  const vectorField& nearestAttraction,
328  const face& localF,
329  const point& cc
330  ) const;
332  (
333  const indirectPrimitivePatch& pp,
334  const vectorField& featureAttraction,
335  const vectorField& nearestAttraction,
336  const face& localF
337  ) const;
339  (
340  const scalar featureCos,
341  const point& newPt0,
342  const pointConstraint& pc0,
343  const point& newPt1,
344  const pointConstraint& pc1
345  ) const;
346  bool isConcave
347  (
348  const point& c0,
349  const vector& area0,
350  const point& c1,
351  const vector& area1,
352  const scalar concaveCos
353  ) const;
355  (
356  const scalar featureCos,
357  const scalar concaveCos,
358  const scalar minAreaFraction,
359  const indirectPrimitivePatch& pp,
360  const vectorField& patchAttraction,
361  const List<pointConstraint>& patchConstraints,
362  const vectorField& nearestAttraction,
363  const vectorField& nearestNormal,
364  const label faceI,
365 
366  DynamicField<point>& points0,
367  DynamicField<point>& points1
368  ) const;
369 
370  //- Do all logic on whether to add face cut to diagonal
371  // attraction
372  void splitDiagonals
373  (
374  const scalar featureCos,
375  const scalar concaveCos,
376  const scalar minAreaFraction,
377 
378  const indirectPrimitivePatch& pp,
379  const vectorField& nearestAttraction,
380  const vectorField& nearestNormal,
381 
382  vectorField& patchAttraction,
383  List<pointConstraint>& patchConstraints,
384  DynamicList<label>& splitFaces,
385  DynamicList<labelPair>& splits
386  ) const;
387 
388  //- Avoid attraction across face diagonal since would
389  // cause face squeeze
391  (
392  const label iter,
393  const scalar featureCos,
394  const indirectPrimitivePatch& pp,
395  vectorField& patchAttraction,
396  List<pointConstraint>& patchConstraints
397  ) const;
398 
399  //- Write some stats about constraints
400  void writeStats
401  (
402  const indirectPrimitivePatch& pp,
403  const PackedBoolList& isMasterPoint,
404  const List<pointConstraint>& patchConstraints
405  ) const;
406 
407  //- Return hit if on multiple points
409  (
410  const point& pt,
411  const labelList& patchIDs,
412  const List<point>& faceCentres
413  ) const;
414 
415  //- Return hit if faces-on-the-same-normalplane are on multiple
416  // patches
418  (
419  const point& pt,
420  const labelList& pfPatchID,
421  const DynamicList<vector>& surfaceNormals,
422  const labelList& faceToNormalBin
423  ) const;
424 
425  //- Return index of similar normal
427  (
428  const scalar featureCos,
429  const vector& faceSurfaceNormal,
430  const DynamicList<vector>& surfaceNormals
431  ) const;
432 
433  //- Determine attraction and constraints for single point
434  // using sampled surrounding of the point
436  (
437  const label iter,
438  const scalar featureCos,
439 
440  const indirectPrimitivePatch& pp,
441  const scalarField& snapDist,
442  const vectorField& nearestDisp,
443  const label pointI,
444 
445  const List<List<point> >& pointFaceSurfNormals,
446  const List<List<point> >& pointFaceDisp,
447  const List<List<point> >& pointFaceCentres,
448  const labelListList& pointFacePatchID,
449 
450  DynamicList<point>& surfacePoints,
451  DynamicList<vector>& surfaceNormals,
452  labelList& faceToNormalBin,
453 
454  vector& patchAttraction,
455  pointConstraint& patchConstraint
456  ) const;
457 
458  //- Determine attraction and constraints for all points
459  // using sampled surrounding of the point
461  (
462  const label iter,
463  const scalar featureCos,
464  const indirectPrimitivePatch& pp,
465  const scalarField& snapDist,
466  const vectorField& nearestDisp,
467 
468  const List<List<point> >& pointFaceSurfNormals,
469  const List<List<point> >& pointFaceDisp,
470  const List<List<point> >& pointFaceCentres,
471  const labelListList& pointFacePatchID,
472 
473  vectorField& patchAttraction,
474  List<pointConstraint>& patchConstraints
475  ) const;
476 
477  //- Determine geometric features and attraction to equivalent
478  // surface features
479  void determineFeatures
480  (
481  const label iter,
482  const scalar featureCos,
483  const bool multiRegionFeatureSnap,
484 
485  const indirectPrimitivePatch&,
486  const scalarField& snapDist,
487  const vectorField& nearestDisp,
488 
489  const List<List<point> >& pointFaceSurfNormals,
490  const List<List<point> >& pointFaceDisp,
491  const List<List<point> >& pointFaceCentres,
492  const labelListList& pointFacePatchID,
493 
494  List<labelList>& pointAttractor,
496  // Feature-edge to pp point
497  List<List<DynamicList<point> > >& edgeAttractors,
498  List<List<DynamicList<pointConstraint> > >& edgeConstraints,
499  vectorField& patchAttraction,
500  List<pointConstraint>& patchConstraints
501  ) const;
502 
503  //- Determine features originating from bafles and
504  // and add attraction to equivalent surface features
506  (
507  const label iter,
508  const bool baffleFeaturePoints,
509  const scalar featureCos,
510 
511  const indirectPrimitivePatch& pp,
512  const scalarField& snapDist,
513 
514  // Feature-point to pp point
515  List<labelList>& pointAttractor,
517  // Feature-edge to pp point
518  List<List<DynamicList<point> > >& edgeAttractors,
519  List<List<DynamicList<pointConstraint> > >& edgeConstraints,
520  // pp point to nearest feature
521  vectorField& patchAttraction,
522  List<pointConstraint>& patchConstraints
523  ) const;
525  (
526  const label iter,
527 
528  const indirectPrimitivePatch& pp,
529  const scalarField& snapDist,
530 
531  // Feature-point to pp point
532  const List<labelList>& pointAttractor,
534  // Feature-edge to pp point
535  const List<List<DynamicList<point> > >& edgeAttractors,
537 
538  const vectorField& rawPatchAttraction,
539  const List<pointConstraint>& rawPatchConstraints,
540 
541  // pp point to nearest feature
542  vectorField& patchAttraction,
543  List<pointConstraint>& patchConstraints
544  ) const;
545 
546  //- Find point on nearest feature edge (within searchDist).
547  // Return point and feature
548  // and store feature-edge to mesh-point and vice versa
550  (
551  const bool isRegionEdge,
552 
553  const indirectPrimitivePatch& pp,
554  const scalarField& snapDist,
555  const label pointI,
556  const point& estimatedPt,
557 
560  vectorField&,
562  ) const;
563 
564  //- Find nearest feature point (within searchDist).
565  // Return feature point
566  // and store feature-point to mesh-point and vice versa.
567  // If another mesh point already referring to this feature
568  // point and further away, reset that one to a near feature
569  // edge (using findNearFeatureEdge above)
571  (
572  const bool isRegionEdge,
573 
574  const indirectPrimitivePatch& pp,
575  const scalarField& snapDist,
576  const label pointI,
577  const point& estimatedPt,
578 
579  // Feature-point to pp point
580  List<labelList>& pointAttractor,
582  // Feature-edge to pp point
583  List<List<DynamicList<point> > >& edgeAttractors,
584  List<List<DynamicList<pointConstraint> > >& edgeConstraints,
585  // pp point to nearest feature
586  vectorField& patchAttraction,
587  List<pointConstraint>& patchConstraints
588  ) const;
589 
591  (
592  const label iter,
593  const bool multiRegionFeatureSnap,
594 
595  const bool detectBaffles,
596  const bool baffleFeaturePoints,
597  const bool releasePoints,
598  const bool stringFeatures,
599  const bool avoidDiagonal,
600 
601  const scalar featureCos,
602 
603  const indirectPrimitivePatch& pp,
604  const scalarField& snapDist,
605  const vectorField& nearestDisp,
606  const vectorField& nearestNormal,
607 
608  const List<List<point> >& pointFaceSurfNormals,
609  const List<List<point> >& pointFaceDisp,
610  const List<List<point> >& pointFaceCentres,
611  const labelListList& pointFacePatchID,
612 
613  vectorField& patchAttraction,
614  List<pointConstraint>& patchConstraints
615  ) const;
616 
617  void preventFaceSqueeze
618  (
619  const label iter,
620  const scalar featureCos,
621  const indirectPrimitivePatch& pp,
622  const scalarField& snapDist,
623  const vectorField& nearestAttraction,
624 
625  vectorField& patchAttraction,
626  List<pointConstraint>& patchConstraints
627  ) const;
628 
629  //- Top level feature attraction routine. Gets given
630  // displacement to nearest surface in nearestDisp
631  // and calculates new displacement taking into account
632  // features
634  (
635  const snapParameters& snapParams,
636  const bool alignMeshEdges,
637  const label iter,
638  const scalar featureCos,
639  const scalar featureAttract,
640  const scalarField& snapDist,
641  const vectorField& nearestDisp,
642  const vectorField& nearestNormal,
643  motionSmoother& meshMover,
644  vectorField& patchAttraction,
645  List<pointConstraint>& patchConstraints,
646 
647  DynamicList<label>& splitFaces,
648  DynamicList<labelPair>& splits
649  ) const;
650 
651 
652  //- Disallow default bitwise copy construct
654 
655  //- Disallow default bitwise assignment
656  void operator=(const autoSnapDriver&);
657 
658 
659 public:
660 
661  //- Runtime type information
662  ClassName("autoSnapDriver");
663 
664 
665  // Constructors
666 
667  //- Construct from components
669  (
670  meshRefinement& meshRefiner,
671  const labelList& globalToMasterPatch,
672  const labelList& globalToSlavePatch
673  );
674 
675 
676  // Member Functions
677 
678  // Snapping
679 
680  //- Merge baffles.
682 
683  //- Calculate edge length per patch point.
685  (
686  const fvMesh& mesh,
687  const snapParameters& snapParams,
689  );
690 
691  //- Smooth the mesh (patch and internal) to increase visibility
692  // of surface points (on castellated mesh) w.r.t. surface.
693  static void preSmoothPatch
694  (
695  const meshRefinement& meshRefiner,
696  const snapParameters& snapParams,
697  const label nInitErrors,
698  const List<labelPair>& baffles,
700  );
701 
702  //- Helper: calculate average cell centre per point
704  (
705  const fvMesh& mesh,
707  );
708 
709  //- Per patch point override displacement if in gap situation
710  void detectNearSurfaces
711  (
712  const scalar planarCos,
713  const indirectPrimitivePatch&,
714  const pointField& nearestPoint,
715  const vectorField& nearestNormal,
716  vectorField& disp
717  ) const;
718 
719  //- Per patch point calculate point on nearest surface. Set as
720  // boundary conditions of motionSmoother displacement field. Return
721  // displacement of patch points.
723  (
724  const bool strictRegionSnap,
725  const meshRefinement& meshRefiner,
726  const labelList& globalToMasterPatch,
727  const labelList& globalToSlavePatch,
728  const scalarField& snapDist,
729  const indirectPrimitivePatch&,
730  pointField& nearestPoint,
731  vectorField& nearestNormal
732  );
733 
737  //static vectorField calcNearestLocalSurface
738  //(
739  // const meshRefinement& meshRefiner,
740  // const scalarField& snapDist,
741  // const indirectPrimitivePatch&
742  //);
743 
744  //- Smooth the displacement field to the internal.
745  void smoothDisplacement
746  (
747  const snapParameters& snapParams,
749  ) const;
750 
751  //- Do the hard work: move the mesh according to displacement,
752  // locally relax the displacement. Return true if ended up with
753  // correct mesh, false if not.
754  bool scaleMesh
755  (
756  const snapParameters& snapParams,
757  const label nInitErrors,
758  const List<labelPair>& baffles,
760  );
761 
762  //- Repatch faces according to surface nearest the face centre
764  (
765  const snapParameters& snapParams,
766  const labelList& adaptPatchIDs,
767  const labelList& preserveFaces
768  );
769 
770  void doSnap
771  (
772  const dictionary& snapDict,
773  const dictionary& motionDict,
774  const bool mergePatchFaces,
775  const scalar featureCos,
776  const scalar planarAngle,
777  const snapParameters& snapParams
778  );
779 
780 };
781 
782 
783 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
784 
785 } // End namespace Foam
786 
787 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
788 
789 #ifdef NoRepository
790 # include "autoSnapDriverTemplates.C"
791 #endif
792 
793 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
794 
795 #endif
796 
797 // ************************************************************************* //
Foam::autoSnapDriver::preSmoothPatch
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.
Definition: autoSnapDriver.C:827
Foam::autoSnapDriver::globalToSlavePatch_
const labelList globalToSlavePatch_
From global surface region to slave side patch.
Definition: autoSnapDriver.H:67
Foam::autoSnapDriver::getZoneSurfacePoints
static labelList getZoneSurfacePoints(const fvMesh &mesh, const indirectPrimitivePatch &, const word &zoneName)
Get points both on patch and facezone.
Definition: autoSnapDriver.C:943
Foam::PackedBoolList
A bit-packed bool list.
Definition: PackedBoolList.H:63
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::autoSnapDriver::isConcave
bool isConcave(const point &c0, const vector &area0, const point &c1, const vector &area1, const scalar concaveCos) const
Definition: autoSnapDriverFeature.C:1616
Foam::fileName
A class for handling file names.
Definition: fileName.H:69
Foam::autoSnapDriver::smoothLambdaMuPatchDisplacement
static pointField smoothLambdaMuPatchDisplacement(const motionSmoother &meshMover, const List< labelPair > &baffles)
Calculate displacement per patch point. Wip.
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::autoSnapDriver::determineFeatures
void determineFeatures(const label iter, const scalar featureCos, const bool multiRegionFeatureSnap, const indirectPrimitivePatch &, const scalarField &snapDist, const vectorField &nearestDisp, const List< List< point > > &pointFaceSurfNormals, const List< List< point > > &pointFaceDisp, const List< List< point > > &pointFaceCentres, const labelListList &pointFacePatchID, List< labelList > &pointAttractor, List< List< pointConstraint > > &pointConstraints, List< List< DynamicList< point > > > &edgeAttractors, List< List< DynamicList< pointConstraint > > > &edgeConstraints, vectorField &patchAttraction, List< pointConstraint > &patchConstraints) const
Determine geometric features and attraction to equivalent.
Definition: autoSnapDriverFeature.C:2242
Foam::DynamicList< label >
Foam::autoSnapDriver::meshRefiner_
meshRefinement & meshRefiner_
Mesh+surface.
Definition: autoSnapDriver.H:61
Foam::autoSnapDriver::scaleMesh
bool scaleMesh(const snapParameters &snapParams, const label nInitErrors, const List< labelPair > &baffles, motionSmoother &)
Do the hard work: move the mesh according to displacement,.
Definition: autoSnapDriver.C:2164
Foam::autoSnapDriver::autoSnapDriver
autoSnapDriver(const autoSnapDriver &)
Disallow default bitwise copy construct.
Foam::pointConstraint
Accumulates point constraints through successive applications of the applyConstraint function.
Definition: pointConstraint.H:56
Foam::autoSnapDriver::smoothAndConstrain
void smoothAndConstrain(const PackedBoolList &isMasterEdge, const indirectPrimitivePatch &pp, const labelList &meshEdges, const List< pointConstraint > &constraints, vectorField &disp) const
Definition: autoSnapDriverFeature.C:130
Foam::autoSnapDriver::calcNearestFacePointProperties
void calcNearestFacePointProperties(const label iter, const indirectPrimitivePatch &pp, const vectorField &faceDisp, const vectorField &faceSurfaceNormal, const labelList &faceSurfaceRegion, List< List< point > > &pointFaceSurfNormals, List< List< point > > &pointFaceDisp, List< List< point > > &pointFaceCentres, List< labelList > &pointFacePatchID) const
Definition: autoSnapDriverFeature.C:447
Foam::autoSnapDriver::getInternalOrBaffleDuplicateFace
labelList getInternalOrBaffleDuplicateFace() const
Get per face -1 or label of opposite face if on internal/baffle.
Definition: autoSnapDriver.C:2493
Foam::autoSnapDriver::findNormal
label findNormal(const scalar featureCos, const vector &faceSurfaceNormal, const DynamicList< vector > &surfaceNormals) const
Return index of similar normal.
Definition: autoSnapDriverFeature.C:753
Foam::autoSnapDriver::isSplitAlignedWithFeature
bool isSplitAlignedWithFeature(const scalar featureCos, const point &newPt0, const pointConstraint &pc0, const point &newPt1, const pointConstraint &pc1) const
Definition: autoSnapDriverFeature.C:1576
Foam::autoSnapDriver::mergeZoneBaffles
autoPtr< mapPolyMesh > mergeZoneBaffles(const List< labelPair > &)
Merge baffles.
Foam::autoSnapDriver::avoidDiagonalAttraction
void avoidDiagonalAttraction(const label iter, const scalar featureCos, const indirectPrimitivePatch &pp, vectorField &patchAttraction, List< pointConstraint > &patchConstraints) const
Avoid attraction across face diagonal since would.
Definition: autoSnapDriverFeature.C:1970
Foam::help::mergePatchFaces
faceList mergePatchFaces(const List< DynList< label > > &pfcs, const pointField &polyPoints)
Definition: helperFunctionsGeometryQueriesI.H:232
Foam::autoSnapDriver::dumpMove
static void dumpMove(const fileName &, const pointField &, const pointField &)
Write displacement as .obj file.
Definition: autoSnapDriver.C:700
Foam::autoSnapDriver::findDiagonalAttraction
labelPair findDiagonalAttraction(const indirectPrimitivePatch &pp, const vectorField &patchAttraction, const List< pointConstraint > &patchConstraints, const label faceI) const
Detect any diagonal attraction. Returns indices in face.
Definition: autoSnapDriverFeature.C:1533
Foam::autoSnapDriver::globalToMasterPatch_
const labelList globalToMasterPatch_
From global surface region to master side patch.
Definition: autoSnapDriver.H:64
Foam::autoSnapDriver::checkCoupledFaceZones
void checkCoupledFaceZones() const
Check that face zones are synced.
Foam::DynamicField
Dynamically sized Field.
Definition: DynamicField.H:49
Foam::autoSnapDriver::findMultiPatchPoint
pointIndexHit findMultiPatchPoint(const point &pt, const labelList &patchIDs, const List< point > &faceCentres) const
Return hit if on multiple points.
Definition: autoSnapDriverFeature.C:729
Foam::autoSnapDriver::centreAndNormal
Tuple2< point, vector > centreAndNormal(const indirectPrimitivePatch &pp, const vectorField &featureAttraction, const vectorField &nearestAttraction, const face &localF) const
Foam::autoSnapDriver::writeStats
void writeStats(const indirectPrimitivePatch &pp, const PackedBoolList &isMasterPoint, const List< pointConstraint > &patchConstraints) const
Write some stats about constraints.
Definition: autoSnapDriverFeature.C:866
Foam::snapParameters
Simple container to keep together snap specific information.
Definition: snapParameters.H:50
Foam::autoSnapDriver::splitDiagonals
void splitDiagonals(const scalar featureCos, const scalar concaveCos, const scalar minAreaFraction, const indirectPrimitivePatch &pp, const vectorField &nearestAttraction, const vectorField &nearestNormal, vectorField &patchAttraction, List< pointConstraint > &patchConstraints, DynamicList< label > &splitFaces, DynamicList< labelPair > &splits) const
Do all logic on whether to add face cut to diagonal.
Definition: autoSnapDriverFeature.C:1885
Foam::PointIndexHit
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:53
Foam::constant::physicoChemical::c1
const dimensionedScalar c1
First radiation constant: default SI units: [W/m2].
Foam::autoSnapDriver::calcSnapDistance
static scalarField calcSnapDistance(const fvMesh &mesh, const snapParameters &snapParams, const indirectPrimitivePatch &)
Calculate edge length per patch point.
Definition: autoSnapDriver.C:787
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::autoSnapDriver::reverseAttractMeshPoints
void reverseAttractMeshPoints(const label iter, const indirectPrimitivePatch &pp, const scalarField &snapDist, const List< labelList > &pointAttractor, const List< List< pointConstraint > > &pointConstraints, const List< List< DynamicList< point > > > &edgeAttractors, const List< List< DynamicList< pointConstraint > > > &, const vectorField &rawPatchAttraction, const List< pointConstraint > &rawPatchConstraints, vectorField &patchAttraction, List< pointConstraint > &patchConstraints) const
Definition: autoSnapDriverFeature.C:3112
Foam::motionSmoother
Given a displacement moves the mesh by scaling the displacement back until there are no more mesh err...
Definition: motionSmoother.H:89
Foam::autoSnapDriver::findNearFeaturePoint
Tuple2< label, pointIndexHit > findNearFeaturePoint(const bool isRegionEdge, const indirectPrimitivePatch &pp, const scalarField &snapDist, const label pointI, const point &estimatedPt, List< labelList > &pointAttractor, List< List< pointConstraint > > &pointConstraints, List< List< DynamicList< point > > > &edgeAttractors, List< List< DynamicList< pointConstraint > > > &edgeConstraints, vectorField &patchAttraction, List< pointConstraint > &patchConstraints) const
Find nearest feature point (within searchDist).
Definition: autoSnapDriverFeature.C:2136
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
meshRefinement.H
Foam::autoSnapDriver::outwardsDisplacement
static bool outwardsDisplacement(const indirectPrimitivePatch &, const vectorField &)
Check displacement is outwards pointing.
Definition: autoSnapDriver.C:729
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::pointMesh
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:48
Foam::autoSnapDriver::calcNearest
void calcNearest(const label iter, const indirectPrimitivePatch &pp, vectorField &pointDisp, vectorField &pointSurfaceNormal, vectorField &pointRotation) const
Foam::autoSnapDriver::smoothDisplacement
void smoothDisplacement(const snapParameters &snapParams, motionSmoother &) const
Smooth the displacement field to the internal.
Definition: autoSnapDriver.C:2095
autoSnapDriverTemplates.C
Foam::autoSnapDriver::avgCellCentres
static tmp< pointField > avgCellCentres(const fvMesh &mesh, const indirectPrimitivePatch &)
Helper: calculate average cell centre per point.
Definition: autoSnapDriver.C:990
Foam::pointConstraints
Application of (multi-)patch point contraints.
Definition: pointConstraints.H:62
Foam::autoSnapDriver::featureAttractionUsingFeatureEdges
void featureAttractionUsingFeatureEdges(const label iter, const bool multiRegionFeatureSnap, const bool detectBaffles, const bool baffleFeaturePoints, const bool releasePoints, const bool stringFeatures, const bool avoidDiagonal, const scalar featureCos, const indirectPrimitivePatch &pp, const scalarField &snapDist, const vectorField &nearestDisp, const vectorField &nearestNormal, const List< List< point > > &pointFaceSurfNormals, const List< List< point > > &pointFaceDisp, const List< List< point > > &pointFaceCentres, const labelListList &pointFacePatchID, vectorField &patchAttraction, List< pointConstraint > &patchConstraints) const
Definition: autoSnapDriverFeature.C:3373
Foam::autoSnapDriver::determineBaffleFeatures
void determineBaffleFeatures(const label iter, const bool baffleFeaturePoints, const scalar featureCos, const indirectPrimitivePatch &pp, const scalarField &snapDist, List< labelList > &pointAttractor, List< List< pointConstraint > > &pointConstraints, List< List< DynamicList< point > > > &edgeAttractors, List< List< DynamicList< pointConstraint > > > &edgeConstraints, vectorField &patchAttraction, List< pointConstraint > &patchConstraints) const
Determine features originating from bafles and.
Definition: autoSnapDriverFeature.C:2784
Foam::autoPtr
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:117
Foam::autoSnapDriver::facePoints
void facePoints(const indirectPrimitivePatch &pp, const vectorField &featureAttraction, const vectorField &nearestAttraction, const face &f, DynamicField< point > &points) const
Foam::autoSnapDriver::findNearFeatureEdge
Tuple2< label, pointIndexHit > findNearFeatureEdge(const bool isRegionEdge, const indirectPrimitivePatch &pp, const scalarField &snapDist, const label pointI, const point &estimatedPt, List< List< DynamicList< point > > > &, List< List< DynamicList< pointConstraint > > > &, vectorField &, List< pointConstraint > &) const
Find point on nearest feature edge (within searchDist).
Definition: autoSnapDriverFeature.C:2069
DynamicField.H
Foam::autoSnapDriver
All to do with snapping to surface.
Definition: autoSnapDriver.H:56
Foam::autoSnapDriver::detectWarpedFaces
void detectWarpedFaces(const scalar featureCos, const indirectPrimitivePatch &pp, DynamicList< label > &splitFaces, DynamicList< labelPair > &splits) const
Detect warpage.
Definition: autoSnapDriver.C:2400
Foam::autoSnapDriver::smoothInternalDisplacement
static tmp< pointField > smoothInternalDisplacement(const meshRefinement &meshRefiner, const motionSmoother &)
Calculate displacement (over all mesh points) to move points.
Definition: autoSnapDriver.C:124
Foam::Pair
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: contiguous.H:49
Foam::autoSnapDriver::correctAttraction
void correctAttraction(const DynamicList< point > &surfacePoints, const DynamicList< label > &surfaceCounts, const point &edgePt, const vector &edgeNormal, const point &pt, vector &edgeOffset) const
Definition: autoSnapDriverFeature.C:691
Foam::autoSnapDriver::smoothPatchDisplacement
static tmp< pointField > smoothPatchDisplacement(const motionSmoother &, const List< labelPair > &)
Calculate displacement per patch point to smooth out patch.
Definition: autoSnapDriver.C:304
Foam::autoSnapDriver::featureAttractionUsingReconstruction
void featureAttractionUsingReconstruction(const label iter, const scalar featureCos, const indirectPrimitivePatch &pp, const scalarField &snapDist, const vectorField &nearestDisp, const label pointI, const List< List< point > > &pointFaceSurfNormals, const List< List< point > > &pointFaceDisp, const List< List< point > > &pointFaceCentres, const labelListList &pointFacePatchID, DynamicList< point > &surfacePoints, DynamicList< vector > &surfaceNormals, labelList &faceToNormalBin, vector &patchAttraction, pointConstraint &patchConstraint) const
Determine attraction and constraints for single point.
Definition: autoSnapDriverFeature.C:914
f
labelList f(nPoints)
Foam::autoSnapDriver::ClassName
ClassName("autoSnapDriver")
Runtime type information.
Foam::Vector< scalar >
Foam::autoSnapDriver::calcNearestSurface
static void calcNearestSurface(const refinementSurfaces &surfaces, const labelList &surfacesToTest, const labelListList &regionsToTest, const pointField &localPoints, const labelList &zonePointIndices, scalarField &minSnapDist, labelList &snapSurf, vectorField &patchDisp, pointField &nearestPoint, vectorField &nearestNormal)
Per patch point calculate point on nearest surface.
Definition: autoSnapDriver.C:1692
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::autoSnapDriver::calcNearestFace
void calcNearestFace(const label iter, const indirectPrimitivePatch &pp, const scalarField &faceSnapDist, vectorField &faceDisp, vectorField &faceSurfaceNormal, labelList &faceSurfaceRegion) const
Definition: autoSnapDriverFeature.C:221
Foam::meshRefinement
Helper class which maintains intersections of (changing) mesh with (static) surfaces.
Definition: meshRefinement.H:82
Foam::autoSnapDriver::pyrVol
scalar pyrVol(const indirectPrimitivePatch &pp, const vectorField &featureAttraction, const face &localF, const point &cc) const
Foam::autoSnapDriver::stringFeatureEdges
void stringFeatureEdges(const label iter, const scalar featureCos, const indirectPrimitivePatch &pp, const scalarField &snapDist, const vectorField &rawPatchAttraction, const List< pointConstraint > &rawPatchConstraints, vectorField &patchAttraction, List< pointConstraint > &patchConstraints) const
For any reverse (so from feature back to mesh) attraction:
Definition: autoSnapDriverFeature.C:1213
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::autoSnapDriver::calcNearestSurfaceFeature
vectorField calcNearestSurfaceFeature(const snapParameters &snapParams, const bool alignMeshEdges, const label iter, const scalar featureCos, const scalar featureAttract, const scalarField &snapDist, const vectorField &nearestDisp, const vectorField &nearestNormal, motionSmoother &meshMover, vectorField &patchAttraction, List< pointConstraint > &patchConstraints, DynamicList< label > &splitFaces, DynamicList< labelPair > &splits) const
Top level feature attraction routine. Gets given.
Definition: autoSnapDriverFeature.C:3807
Foam::autoSnapDriver::getCollocatedPoints
static label getCollocatedPoints(const scalar tol, const pointField &, PackedBoolList &)
Calculates (geometric) shared points.
Definition: autoSnapDriver.C:64
Foam::autoSnapDriver::repatchToSurface
autoPtr< mapPolyMesh > repatchToSurface(const snapParameters &snapParams, const labelList &adaptPatchIDs, const labelList &preserveFaces)
Repatch faces according to surface nearest the face centre.
Definition: autoSnapDriver.C:2230
Foam::autoSnapDriver::edgePatchDist
static tmp< scalarField > edgePatchDist(const pointMesh &, const indirectPrimitivePatch &)
Per edge distance to patch.
Definition: autoSnapDriver.C:654
Foam::autoSnapDriver::avg
static tmp< pointField > avg(const indirectPrimitivePatch &, const pointField &)
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
Foam::autoSnapDriver::isFeaturePoint
bool isFeaturePoint(const scalar featureCos, const indirectPrimitivePatch &pp, const PackedBoolList &isFeatureEdge, const label pointI) const
Is point on two feature edges that make a largish angle?
Definition: autoSnapDriverFeature.C:69
Foam::autoSnapDriver::operator=
void operator=(const autoSnapDriver &)
Disallow default bitwise assignment.
Foam::autoSnapDriver::releasePointsNextToMultiPatch
void releasePointsNextToMultiPatch(const label iter, const scalar featureCos, const indirectPrimitivePatch &pp, const scalarField &snapDist, const List< List< point > > &pointFaceCentres, const labelListList &pointFacePatchID, const vectorField &rawPatchAttraction, const List< pointConstraint > &rawPatchConstraints, vectorField &patchAttraction, List< pointConstraint > &patchConstraints) const
Remove constraints of points next to multi-patch points.
Definition: autoSnapDriverFeature.C:1402
Foam::Tuple2
A 2-tuple for storing two objects of different types.
Definition: Tuple2.H:47
Foam::autoSnapDriver::doSnap
void doSnap(const dictionary &snapDict, const dictionary &motionDict, const bool mergePatchFaces, const scalar featureCos, const scalar planarAngle, const snapParameters &snapParams)
Definition: autoSnapDriver.C:2528
Foam::autoSnapDriver::preventFaceSqueeze
void preventFaceSqueeze(const label iter, const scalar featureCos, const indirectPrimitivePatch &pp, const scalarField &snapDist, const vectorField &nearestAttraction, vectorField &patchAttraction, List< pointConstraint > &patchConstraints) const
Definition: autoSnapDriverFeature.C:3671
Foam::autoSnapDriver::detectNearSurfaces
void detectNearSurfaces(const scalar planarCos, const indirectPrimitivePatch &, const pointField &nearestPoint, const vectorField &nearestNormal, vectorField &disp) const
Per patch point override displacement if in gap situation.
Definition: autoSnapDriver.C:1094
Foam::refinementSurfaces
Container for data on surfaces used for surface-driven refinement. Contains all the data about the le...
Definition: refinementSurfaces.H:60
Foam::PrimitivePatch
A list of faces which address into the list of points.
Definition: PrimitivePatchTemplate.H:88
Foam::autoSnapDriver::getFacePoints
static labelList getFacePoints(const indirectPrimitivePatch &pp, const FaceList &faces)
Get points both on patch and facezone.