snappySnapDriver.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 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::snappySnapDriver
26 
27 Description
28  All to do with snapping to surface
29 
30 SourceFiles
31  snappySnapDriver.C
32  snappySnapDriverFeature.C
33 
34 \*---------------------------------------------------------------------------*/
35 
36 #ifndef snappySnapDriver_H
37 #define snappySnapDriver_H
38 
39 #include "meshRefinement.H"
40 
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 
43 namespace Foam
44 {
45 
46 // Forward declaration of classes
47 class motionSmoother;
48 class snapParameters;
49 class pointConstraint;
50 
51 /*---------------------------------------------------------------------------*\
52  Class snappySnapDriver Declaration
53 \*---------------------------------------------------------------------------*/
54 
55 class snappySnapDriver
56 {
57  // Private data
58 
59  //- Mesh+surface
60  meshRefinement& meshRefiner_;
61 
62  //- From global surface region to master side patch
63  const labelList globalToMasterPatch_;
64 
65  //- From global surface region to slave side patch
66  const labelList globalToSlavePatch_;
67 
68 
69  // Private Member Functions
70 
71 
72  // Snapping
73 
74  //- Calculates (geometric) shared points
75  // Requires PackedList to be sized and initialised
76  static label getCollocatedPoints
77  (
78  const scalar tol,
79  const pointField&,
81  );
82 
83  //- Calculate displacement per patch point to smooth out patch.
84  // Quite complicated in determining which points to move where.
85  static pointField smoothPatchDisplacement
86  (
87  const motionSmoother&,
88  const List<labelPair>&
89  );
90 
91  //- Check that face zones are synced
92  void checkCoupledFaceZones() const;
93 
94  //- Per edge distance to patch
95  static tmp<scalarField> edgePatchDist
96  (
97  const pointMesh&,
99  );
100 
101  //- Write displacement as .obj file.
102  static void dumpMove
103  (
104  const fileName&,
105  const pointField&,
106  const pointField&
107  );
108 
109  //- Check displacement is outwards pointing
110  static bool outwardsDisplacement
111  (
112  const indirectPrimitivePatch&,
113  const vectorField&
114  );
115 
116  //- Detect warpage
117  void detectWarpedFaces
118  (
119  const scalar featureCos,
120  const indirectPrimitivePatch& pp,
121 
122  DynamicList<label>& splitFaces,
123  DynamicList<labelPair>& splits
124  ) const;
125 
126  // Feature line snapping
127 
128  //- Is point on two feature edges that make a largish angle?
129  bool isFeaturePoint
130  (
131  const scalar featureCos,
132  const indirectPrimitivePatch& pp,
133  const PackedBoolList& isFeatureEdge,
134  const label pointi
135  ) const;
136 
137  void smoothAndConstrain
138  (
139  const PackedBoolList& isMasterEdge,
140  const indirectPrimitivePatch& pp,
141  const labelList& meshEdges,
142  const List<pointConstraint>& constraints,
143  vectorField& disp
144  ) const;
145 
146  void calcNearest
147  (
148  const label iter,
149  const indirectPrimitivePatch& pp,
150  vectorField& pointDisp,
151  vectorField& pointSurfaceNormal,
152  vectorField& pointRotation
153  ) const;
154 
155  void calcNearestFace
156  (
157  const label iter,
158  const indirectPrimitivePatch& pp,
159  const scalarField& faceSnapDist,
160  vectorField& faceDisp,
161  vectorField& faceSurfaceNormal,
162  labelList& faceSurfaceRegion,
163  vectorField& faceRotation
164  ) const;
165 
166  //- Collect (possibly remote) per point data of all surrounding
167  // faces
168  // - faceSurfaceNormal
169  // - faceDisp
170  // - faceCentres&faceNormal
171  void calcNearestFacePointProperties
172  (
173  const label iter,
174  const indirectPrimitivePatch& pp,
175 
176  const vectorField& faceDisp,
177  const vectorField& faceSurfaceNormal,
178  const labelList& faceSurfaceRegion,
179 
180  List<List<point>>& pointFaceSurfNormals,
181  List<List<point>>& pointFaceDisp,
182  List<List<point>>& pointFaceCentres,
183  List<labelList>& pointFacePatchID
184  ) const;
185 
186  //- Gets passed in offset to nearest point on feature
187  // edge. Calculates if the point has a different number of
188  // faces on either side of the feature and if so attracts the
189  // point to that non-dominant plane.
190  void correctAttraction
191  (
192  const DynamicList<point>& surfacePoints,
193  const DynamicList<label>& surfaceCounts,
194  const point& edgePt,
195  const vector& edgeNormal, // normalised normal
196  const point& pt,
197  vector& edgeOffset // offset from pt to point on edge
198  ) const;
199 
200 
201  //- For any reverse (so from feature back to mesh) attraction:
202  // add attraction if diagonal points on face attracted
203  void stringFeatureEdges
204  (
205  const label iter,
206  const scalar featureCos,
207 
208  const indirectPrimitivePatch& pp,
209  const scalarField& snapDist,
210 
211  const vectorField& rawPatchAttraction,
212  const List<pointConstraint>& rawPatchConstraints,
213 
214  vectorField& patchAttraction,
215  List<pointConstraint>& patchConstraints
216  ) const;
217 
218  //- Remove constraints of points next to multi-patch points
219  // to give a bit more freedom of the mesh to conform to the
220  // multi-patch points. Bit dodgy for simple cases.
221  void releasePointsNextToMultiPatch
222  (
223  const label iter,
224  const scalar featureCos,
225 
226  const indirectPrimitivePatch& pp,
227  const scalarField& snapDist,
228 
229  const List<List<point>>& pointFaceCentres,
230  const labelListList& pointFacePatchID,
231 
232  const vectorField& rawPatchAttraction,
233  const List<pointConstraint>& rawPatchConstraints,
234 
235  vectorField& patchAttraction,
236  List<pointConstraint>& patchConstraints
237  ) const;
238 
239  //- Detect any diagonal attraction. Returns indices in face
240  // or (-1, -1) if none
241  labelPair findDiagonalAttraction
242  (
243  const indirectPrimitivePatch& pp,
244  const vectorField& patchAttraction,
245  const List<pointConstraint>& patchConstraints,
246  const label facei
247  ) const;
248 
249  //- Avoid attraction across face diagonal since would
250  // cause face squeeze
251  void avoidDiagonalAttraction
252  (
253  const label iter,
254  const scalar featureCos,
255  const indirectPrimitivePatch& pp,
256  vectorField& patchAttraction,
257  List<pointConstraint>& patchConstraints
258  ) const;
259 
260  //- Return hit if on multiple points
261  pointIndexHit findMultiPatchPoint
262  (
263  const point& pt,
264  const labelList& patchIDs,
265  const List<point>& faceCentres
266  ) const;
267 
268  //- Return hit if faces-on-the-same-normalplane are on multiple
269  // patches
270  // - false, index=-1 : single patch
271  // - true , index=0 : multiple patches but on different
272  // normals planes (so geometric feature
273  // edge is also a region edge)
274  // - true , index=1 : multiple patches on same normals plane
275  // i.e. flat region edge
276  pointIndexHit findMultiPatchPoint
277  (
278  const point& pt,
279  const labelList& pfPatchID,
280  const DynamicList<vector>& surfaceNormals,
281  const labelList& faceToNormalBin
282  ) const;
283 
284  //- Return index of similar normal
285  label findNormal
286  (
287  const scalar featureCos,
288  const vector& faceSurfaceNormal,
289  const DynamicList<vector>& surfaceNormals
290  ) const;
291 
292  //- Determine attraction and constraints for single point
293  // using sampled surrounding of the point
294  void featureAttractionUsingReconstruction
295  (
296  const label iter,
297  const scalar featureCos,
298 
299  const indirectPrimitivePatch& pp,
300  const scalarField& snapDist,
301  const vectorField& nearestDisp,
302  const label pointi,
303 
304  const List<List<point>>& pointFaceSurfNormals,
305  const List<List<point>>& pointFaceDisp,
306  const List<List<point>>& pointFaceCentres,
307  const labelListList& pointFacePatchID,
308 
309  DynamicList<point>& surfacePoints,
310  DynamicList<vector>& surfaceNormals,
311  labelList& faceToNormalBin,
312 
313  vector& patchAttraction,
314  pointConstraint& patchConstraint
315  ) const;
316 
317  //- Determine attraction and constraints for all points
318  // using sampled surrounding of the point
319  void featureAttractionUsingReconstruction
320  (
321  const label iter,
322  const bool avoidSnapProblems,
323  const scalar featureCos,
324  const indirectPrimitivePatch& pp,
325  const scalarField& snapDist,
326  const vectorField& nearestDisp,
327 
328  const List<List<point>>& pointFaceSurfNormals,
329  const List<List<point>>& pointFaceDisp,
330  const List<List<point>>& pointFaceCentres,
331  const labelListList& pointFacePatchID,
332 
333  vectorField& patchAttraction,
334  List<pointConstraint>& patchConstraints
335  ) const;
336 
337  //- Determine geometric features and attraction to equivalent
338  // surface features
339  void determineFeatures
340  (
341  const label iter,
342  const scalar featureCos,
343  const bool multiRegionFeatureSnap,
344 
345  const indirectPrimitivePatch&,
346  const scalarField& snapDist,
347  const vectorField& nearestDisp,
348 
349  const List<List<point>>& pointFaceSurfNormals,
350  const List<List<point>>& pointFaceDisp,
351  const List<List<point>>& pointFaceCentres,
352  const labelListList& pointFacePatchID,
353 
354  List<labelList>& pointAttractor,
356  // Feature-edge to pp point
357  List<List<DynamicList<point>>>& edgeAttractors,
358  List<List<DynamicList<pointConstraint>>>& edgeConstraints,
359  vectorField& patchAttraction,
360  List<pointConstraint>& patchConstraints
361  ) const;
362 
363  //- Determine features originating from bafles and
364  // and add attraction to equivalent surface features
365  void determineBaffleFeatures
366  (
367  const label iter,
368  const scalar featureCos,
369 
370  const indirectPrimitivePatch& pp,
371  const scalarField& snapDist,
372 
373  // Feature-point to pp point
374  List<labelList>& pointAttractor,
376  // Feature-edge to pp point
377  List<List<DynamicList<point>>>& edgeAttractors,
378  List<List<DynamicList<pointConstraint>>>& edgeConstraints,
379  // pp point to nearest feature
380  vectorField& patchAttraction,
381  List<pointConstraint>& patchConstraints
382  ) const;
383 
384  void reverseAttractMeshPoints
385  (
386  const label iter,
387 
388  const indirectPrimitivePatch& pp,
389  const scalarField& snapDist,
390 
391  // Feature-point to pp point
392  const List<labelList>& pointAttractor,
394  // Feature-edge to pp point
395  const List<List<DynamicList<point>>>& edgeAttractors,
397 
398  const vectorField& rawPatchAttraction,
399  const List<pointConstraint>& rawPatchConstraints,
400 
401  // pp point to nearest feature
402  vectorField& patchAttraction,
403  List<pointConstraint>& patchConstraints
404  ) const;
405 
406  //- Find point on nearest feature edge (within searchDist).
407  // Return point and feature
408  // and store feature-edge to mesh-point and vice versa
409  Tuple2<label, pointIndexHit> findNearFeatureEdge
410  (
411  const bool isRegionEdge,
412 
413  const indirectPrimitivePatch& pp,
414  const scalarField& snapDist,
415  const label pointi,
416  const point& estimatedPt,
417 
420  vectorField&,
422  ) const;
423 
424  //- Find nearest feature point (within searchDist).
425  // Return feature point
426  // and store feature-point to mesh-point and vice versa.
427  // If another mesh point already referring to this feature
428  // point and further away, reset that one to a near feature
429  // edge (using findNearFeatureEdge above)
430  Tuple2<label, pointIndexHit> findNearFeaturePoint
431  (
432  const bool isRegionEdge,
433 
434  const indirectPrimitivePatch& pp,
435  const scalarField& snapDist,
436  const label pointi,
437  const point& estimatedPt,
438 
439  // Feature-point to pp point
440  List<labelList>& pointAttractor,
442  // Feature-edge to pp point
443  List<List<DynamicList<point>>>& edgeAttractors,
444  List<List<DynamicList<pointConstraint>>>& edgeConstraints,
445  // pp point to nearest feature
446  vectorField& patchAttraction,
447  List<pointConstraint>& patchConstraints
448  ) const;
449 
450  void featureAttractionUsingFeatureEdges
451  (
452  const label iter,
453  const bool avoidSnapProblems,
454  const scalar featureCos,
455  const bool multiRegionFeatureSnap,
456  const indirectPrimitivePatch& pp,
457  const scalarField& snapDist,
458  const vectorField& nearestDisp,
459 
460  const List<List<point>>& pointFaceSurfNormals,
461  const List<List<point>>& pointFaceDisp,
462  const List<List<point>>& pointFaceCentres,
463  const labelListList& pointFacePatchID,
464 
465  vectorField& patchAttraction,
466  List<pointConstraint>& patchConstraints
467  ) const;
468 
469  void preventFaceSqueeze
470  (
471  const label iter,
472  const scalar featureCos,
473  const indirectPrimitivePatch& pp,
474  const scalarField& snapDist,
475 
476  vectorField& patchAttraction,
477  List<pointConstraint>& patchConstraints
478  ) const;
479 
480  //- Top level feature attraction routine. Gets given
481  // displacement to nearest surface in nearestDisp
482  // and calculates new displacement taking into account
483  // features
484  vectorField calcNearestSurfaceFeature
485  (
486  const snapParameters& snapParams,
487  const bool avoidSnapProblems,
488  const label iter,
489  const scalar featureCos,
490  const scalar featureAttract,
491  const scalarField& snapDist,
492  const vectorField& nearestDisp,
493  motionSmoother& meshMover,
494  vectorField& patchAttraction,
495  List<pointConstraint>& patchConstraints
496  ) const;
497 
498 
499  //- Disallow default bitwise copy construct
501 
502  //- Disallow default bitwise assignment
503  void operator=(const snappySnapDriver&);
504 
505 
506 public:
507 
508  //- Runtime type information
509  ClassName("snappySnapDriver");
510 
511 
512  // Constructors
513 
514  //- Construct from components
516  (
517  meshRefinement& meshRefiner,
518  const labelList& globalToMasterPatch,
519  const labelList& globalToSlavePatch
520  );
521 
522 
523  // Member Functions
524 
525  // Snapping
526 
527  //- Merge baffles.
529 
530  //- Calculate edge length per patch point.
532  (
533  const fvMesh& mesh,
534  const snapParameters& snapParams,
536  );
537 
538  //- Smooth the mesh (patch and internal) to increase visibility
539  // of surface points (on castellated mesh) w.r.t. surface.
540  static void preSmoothPatch
541  (
542  const meshRefinement& meshRefiner,
543  const snapParameters& snapParams,
544  const label nInitErrors,
545  const List<labelPair>& baffles,
547  );
548 
549  //- Get points both on patch and facezone.
551  (
552  const fvMesh& mesh,
553  const indirectPrimitivePatch&,
554  const word& zoneName
555  );
556 
557  //- Helper: calculate average cell centre per point
559  (
560  const fvMesh& mesh,
562  );
563 
564  //- Per patch point override displacement if in gap situation
565  void detectNearSurfaces
566  (
567  const scalar planarCos,
568  const indirectPrimitivePatch&,
569  const pointField& nearestPoint,
570  const vectorField& nearestNormal,
571  vectorField& disp
572  ) const;
573 
574  //- Per patch point calculate point on nearest surface. Set as
575  // boundary conditions of motionSmoother displacement field. Return
576  // displacement of patch points.
578  (
579  const meshRefinement& meshRefiner,
580  const scalarField& snapDist,
581  const indirectPrimitivePatch&,
582  pointField& nearestPoint,
583  vectorField& nearestNormal
584  );
585 
586  //- Smooth the displacement field to the internal.
587  void smoothDisplacement
588  (
589  const snapParameters& snapParams,
591  ) const;
592 
593  //- Do the hard work: move the mesh according to displacement,
594  // locally relax the displacement. Return true if ended up with
595  // correct mesh, false if not.
596  bool scaleMesh
597  (
598  const snapParameters& snapParams,
599  const label nInitErrors,
600  const List<labelPair>& baffles,
602  );
603 
604  //- Repatch faces according to surface nearest the face centre
605  // - calculate face-wise snap distance as max of point-wise
606  // - calculate face-wise nearest surface point
607  // - repatch face according to patch for surface point.
609  (
610  const snapParameters& snapParams,
611  const labelList& adaptPatchIDs,
612  const labelList& preserveFaces
613  );
614 
615  void doSnap
616  (
617  const dictionary& snapDict,
618  const dictionary& motionDict,
619  const scalar featureCos,
620  const scalar planarAngle,
621  const snapParameters& snapParams
622  );
623 };
624 
625 
626 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
627 
628 } // End namespace Foam
629 
630 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
631 
632 #endif
633 
634 // ************************************************************************* //
Given a displacement moves the mesh by scaling the displacement back until there are no more mesh err...
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
A class for handling file names.
Definition: fileName.H:69
bool scaleMesh(const snapParameters &snapParams, const label nInitErrors, const List< labelPair > &baffles, motionSmoother &)
Do the hard work: move the mesh according to displacement,.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
A 2-tuple for storing two objects of different types.
Definition: HashTable.H:66
autoPtr< mapPolyMesh > mergeZoneBaffles(const List< labelPair > &)
Merge baffles.
void doSnap(const dictionary &snapDict, const dictionary &motionDict, const scalar featureCos, const scalar planarAngle, const snapParameters &snapParams)
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.
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:53
static labelList getZoneSurfacePoints(const fvMesh &mesh, const indirectPrimitivePatch &, const word &zoneName)
Get points both on patch and facezone.
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:48
Application of (multi-)patch point constraints.
A list of faces which address into the list of points.
ClassName("snappySnapDriver")
Runtime type information.
dynamicFvMesh & mesh
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: contiguous.H:49
A class for handling words, derived from string.
Definition: word.H:59
Simple container to keep together snap specific information.
autoPtr< mapPolyMesh > repatchToSurface(const snapParameters &snapParams, const labelList &adaptPatchIDs, const labelList &preserveFaces)
Repatch faces according to surface nearest the face centre.
void smoothDisplacement(const snapParameters &snapParams, motionSmoother &) const
Smooth the displacement field to the internal.
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.
Accumulates point constraints through successive applications of the applyConstraint function...
static scalarField calcSnapDistance(const fvMesh &mesh, const snapParameters &snapParams, const indirectPrimitivePatch &)
Calculate edge length per patch point.
static tmp< pointField > avgCellCentres(const fvMesh &mesh, const indirectPrimitivePatch &)
Helper: calculate average cell centre per point.
Helper class which maintains intersections of (changing) mesh with (static) surfaces.
All to do with snapping to surface.
A bit-packed bool list.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
A class for managing temporary objects.
Definition: PtrList.H:53
Namespace for OpenFOAM.
static vectorField calcNearestSurface(const meshRefinement &meshRefiner, const scalarField &snapDist, const indirectPrimitivePatch &, pointField &nearestPoint, vectorField &nearestNormal)
Per patch point calculate point on nearest surface. Set as.