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 | Copyright (C) 2011-2016 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  static label getCollocatedPoints
76  (
77  const scalar tol,
78  const pointField&,
80  );
81 
82  //- Calculate displacement per patch point to smooth out patch.
83  // Quite complicated in determining which points to move where.
84  static pointField smoothPatchDisplacement
85  (
86  const motionSmoother&,
87  const List<labelPair>&
88  );
89 
90  //- Check that face zones are synced
91  void checkCoupledFaceZones() const;
92 
93  //- Per edge distance to patch
94  static tmp<scalarField> edgePatchDist
95  (
96  const pointMesh&,
98  );
99 
100  //- Write displacement as .obj file.
101  static void dumpMove
102  (
103  const fileName&,
104  const pointField&,
105  const pointField&
106  );
107 
108  //- Check displacement is outwards pointing
109  static bool outwardsDisplacement
110  (
111  const indirectPrimitivePatch&,
112  const vectorField&
113  );
114 
115  //- Detect warpage
116  void detectWarpedFaces
117  (
118  const scalar featureCos,
119  const indirectPrimitivePatch& pp,
120 
121  DynamicList<label>& splitFaces,
122  DynamicList<labelPair>& splits
123  ) const;
124 
125  // Feature line snapping
126 
127  //- Is point on two feature edges that make a largish angle?
128  bool isFeaturePoint
129  (
130  const scalar featureCos,
131  const indirectPrimitivePatch& pp,
132  const PackedBoolList& isFeatureEdge,
133  const label pointi
134  ) const;
135 
136  void smoothAndConstrain
137  (
138  const PackedBoolList& isMasterEdge,
139  const indirectPrimitivePatch& pp,
140  const labelList& meshEdges,
141  const List<pointConstraint>& constraints,
142  vectorField& disp
143  ) const;
144 
145  void calcNearest
146  (
147  const label iter,
148  const indirectPrimitivePatch& pp,
149  vectorField& pointDisp,
150  vectorField& pointSurfaceNormal,
151  vectorField& pointRotation
152  ) const;
153 
154  void calcNearestFace
155  (
156  const label iter,
157  const indirectPrimitivePatch& pp,
158  const scalarField& faceSnapDist,
159  vectorField& faceDisp,
160  vectorField& faceSurfaceNormal,
161  labelList& faceSurfaceRegion,
162  vectorField& faceRotation
163  ) const;
164 
165  //- Collect (possibly remote) per point data of all surrounding
166  // faces
167  // - faceSurfaceNormal
168  // - faceDisp
169  // - faceCentres&faceNormal
170  void calcNearestFacePointProperties
171  (
172  const label iter,
173  const indirectPrimitivePatch& pp,
174 
175  const vectorField& faceDisp,
176  const vectorField& faceSurfaceNormal,
177  const labelList& faceSurfaceRegion,
178 
179  List<List<point>>& pointFaceSurfNormals,
180  List<List<point>>& pointFaceDisp,
181  List<List<point>>& pointFaceCentres,
182  List<labelList>& pointFacePatchID
183  ) const;
184 
185  //- Gets passed in offset to nearest point on feature
186  // edge. Calculates if the point has a different number of
187  // faces on either side of the feature and if so attracts the
188  // point to that non-dominant plane.
189  void correctAttraction
190  (
191  const DynamicList<point>& surfacePoints,
192  const DynamicList<label>& surfaceCounts,
193  const point& edgePt,
194  const vector& edgeNormal, // normalised normal
195  const point& pt,
196  vector& edgeOffset // offset from pt to point on edge
197  ) const;
198 
199 
200  //- For any reverse (so from feature back to mesh) attraction:
201  // add attraction if diagonal points on face attracted
202  void stringFeatureEdges
203  (
204  const label iter,
205  const scalar featureCos,
206 
207  const indirectPrimitivePatch& pp,
208  const scalarField& snapDist,
209 
210  const vectorField& rawPatchAttraction,
211  const List<pointConstraint>& rawPatchConstraints,
212 
213  vectorField& patchAttraction,
214  List<pointConstraint>& patchConstraints
215  ) const;
216 
217  //- Remove constraints of points next to multi-patch points
218  // to give a bit more freedom of the mesh to conform to the
219  // multi-patch points. Bit dodgy for simple cases.
220  void releasePointsNextToMultiPatch
221  (
222  const label iter,
223  const scalar featureCos,
224 
225  const indirectPrimitivePatch& pp,
226  const scalarField& snapDist,
227 
228  const List<List<point>>& pointFaceCentres,
229  const labelListList& pointFacePatchID,
230 
231  const vectorField& rawPatchAttraction,
232  const List<pointConstraint>& rawPatchConstraints,
233 
234  vectorField& patchAttraction,
235  List<pointConstraint>& patchConstraints
236  ) const;
237 
238  //- Detect any diagonal attraction. Returns indices in face
239  // or (-1, -1) if none
240  labelPair findDiagonalAttraction
241  (
242  const indirectPrimitivePatch& pp,
243  const vectorField& patchAttraction,
244  const List<pointConstraint>& patchConstraints,
245  const label facei
246  ) const;
247 
248  //- Avoid attraction across face diagonal since would
249  // cause face squeeze
250  void avoidDiagonalAttraction
251  (
252  const label iter,
253  const scalar featureCos,
254  const indirectPrimitivePatch& pp,
255  vectorField& patchAttraction,
256  List<pointConstraint>& patchConstraints
257  ) const;
258 
259  //- Return hit if on multiple points
260  pointIndexHit findMultiPatchPoint
261  (
262  const point& pt,
263  const labelList& patchIDs,
264  const List<point>& faceCentres
265  ) const;
266 
267  //- Return hit if faces-on-the-same-normalplane are on multiple
268  // patches
269  // - false, index=-1 : single patch
270  // - true , index=0 : multiple patches but on different
271  // normals planes (so geometric feature
272  // edge is also a region edge)
273  // - true , index=1 : multiple patches on same normals plane
274  // i.e. flat region edge
275  pointIndexHit findMultiPatchPoint
276  (
277  const point& pt,
278  const labelList& pfPatchID,
279  const DynamicList<vector>& surfaceNormals,
280  const labelList& faceToNormalBin
281  ) const;
282 
283  //- Return index of similar normal
284  label findNormal
285  (
286  const scalar featureCos,
287  const vector& faceSurfaceNormal,
288  const DynamicList<vector>& surfaceNormals
289  ) const;
290 
291  //- Determine attraction and constraints for single point
292  // using sampled surrounding of the point
293  void featureAttractionUsingReconstruction
294  (
295  const label iter,
296  const scalar featureCos,
297 
298  const indirectPrimitivePatch& pp,
299  const scalarField& snapDist,
300  const vectorField& nearestDisp,
301  const label pointi,
302 
303  const List<List<point>>& pointFaceSurfNormals,
304  const List<List<point>>& pointFaceDisp,
305  const List<List<point>>& pointFaceCentres,
306  const labelListList& pointFacePatchID,
307 
308  DynamicList<point>& surfacePoints,
309  DynamicList<vector>& surfaceNormals,
310  labelList& faceToNormalBin,
311 
312  vector& patchAttraction,
313  pointConstraint& patchConstraint
314  ) const;
315 
316  //- Determine attraction and constraints for all points
317  // using sampled surrounding of the point
318  void featureAttractionUsingReconstruction
319  (
320  const label iter,
321  const bool avoidSnapProblems,
322  const scalar featureCos,
323  const indirectPrimitivePatch& pp,
324  const scalarField& snapDist,
325  const vectorField& nearestDisp,
326 
327  const List<List<point>>& pointFaceSurfNormals,
328  const List<List<point>>& pointFaceDisp,
329  const List<List<point>>& pointFaceCentres,
330  const labelListList& pointFacePatchID,
331 
332  vectorField& patchAttraction,
333  List<pointConstraint>& patchConstraints
334  ) const;
335 
336  //- Determine geometric features and attraction to equivalent
337  // surface features
338  void determineFeatures
339  (
340  const label iter,
341  const scalar featureCos,
342  const bool multiRegionFeatureSnap,
343 
344  const indirectPrimitivePatch&,
345  const scalarField& snapDist,
346  const vectorField& nearestDisp,
347 
348  const List<List<point>>& pointFaceSurfNormals,
349  const List<List<point>>& pointFaceDisp,
350  const List<List<point>>& pointFaceCentres,
351  const labelListList& pointFacePatchID,
352 
353  List<labelList>& pointAttractor,
355  // Feature-edge to pp point
356  List<List<DynamicList<point>>>& edgeAttractors,
357  List<List<DynamicList<pointConstraint>>>& edgeConstraints,
358  vectorField& patchAttraction,
359  List<pointConstraint>& patchConstraints
360  ) const;
361 
362  //- Determine features originating from bafles and
363  // and add attraction to equivalent surface features
364  void determineBaffleFeatures
365  (
366  const label iter,
367  const scalar featureCos,
368 
369  const indirectPrimitivePatch& pp,
370  const scalarField& snapDist,
371 
372  // Feature-point to pp point
373  List<labelList>& pointAttractor,
375  // Feature-edge to pp point
376  List<List<DynamicList<point>>>& edgeAttractors,
377  List<List<DynamicList<pointConstraint>>>& edgeConstraints,
378  // pp point to nearest feature
379  vectorField& patchAttraction,
380  List<pointConstraint>& patchConstraints
381  ) const;
382 
383  void reverseAttractMeshPoints
384  (
385  const label iter,
386 
387  const indirectPrimitivePatch& pp,
388  const scalarField& snapDist,
389 
390  // Feature-point to pp point
391  const List<labelList>& pointAttractor,
393  // Feature-edge to pp point
394  const List<List<DynamicList<point>>>& edgeAttractors,
396 
397  const vectorField& rawPatchAttraction,
398  const List<pointConstraint>& rawPatchConstraints,
399 
400  // pp point to nearest feature
401  vectorField& patchAttraction,
402  List<pointConstraint>& patchConstraints
403  ) const;
404 
405  //- Find point on nearest feature edge (within searchDist).
406  // Return point and feature
407  // and store feature-edge to mesh-point and vice versa
408  Tuple2<label, pointIndexHit> findNearFeatureEdge
409  (
410  const bool isRegionEdge,
411 
412  const indirectPrimitivePatch& pp,
413  const scalarField& snapDist,
414  const label pointi,
415  const point& estimatedPt,
416 
419  vectorField&,
421  ) const;
422 
423  //- Find nearest feature point (within searchDist).
424  // Return feature point
425  // and store feature-point to mesh-point and vice versa.
426  // If another mesh point already referring to this feature
427  // point and further away, reset that one to a near feature
428  // edge (using findNearFeatureEdge above)
429  Tuple2<label, pointIndexHit> findNearFeaturePoint
430  (
431  const bool isRegionEdge,
432 
433  const indirectPrimitivePatch& pp,
434  const scalarField& snapDist,
435  const label pointi,
436  const point& estimatedPt,
437 
438  // Feature-point to pp point
439  List<labelList>& pointAttractor,
441  // Feature-edge to pp point
442  List<List<DynamicList<point>>>& edgeAttractors,
443  List<List<DynamicList<pointConstraint>>>& edgeConstraints,
444  // pp point to nearest feature
445  vectorField& patchAttraction,
446  List<pointConstraint>& patchConstraints
447  ) const;
448 
449  void featureAttractionUsingFeatureEdges
450  (
451  const label iter,
452  const bool avoidSnapProblems,
453  const scalar featureCos,
454  const bool multiRegionFeatureSnap,
455  const indirectPrimitivePatch& pp,
456  const scalarField& snapDist,
457  const vectorField& nearestDisp,
458 
459  const List<List<point>>& pointFaceSurfNormals,
460  const List<List<point>>& pointFaceDisp,
461  const List<List<point>>& pointFaceCentres,
462  const labelListList& pointFacePatchID,
463 
464  vectorField& patchAttraction,
465  List<pointConstraint>& patchConstraints
466  ) const;
467 
468  void preventFaceSqueeze
469  (
470  const label iter,
471  const scalar featureCos,
472  const indirectPrimitivePatch& pp,
473  const scalarField& snapDist,
474 
475  vectorField& patchAttraction,
476  List<pointConstraint>& patchConstraints
477  ) const;
478 
479  //- Top level feature attraction routine. Gets given
480  // displacement to nearest surface in nearestDisp
481  // and calculates new displacement taking into account
482  // features
483  vectorField calcNearestSurfaceFeature
484  (
485  const snapParameters& snapParams,
486  const bool avoidSnapProblems,
487  const label iter,
488  const scalar featureCos,
489  const scalar featureAttract,
490  const scalarField& snapDist,
491  const vectorField& nearestDisp,
492  motionSmoother& meshMover,
493  vectorField& patchAttraction,
494  List<pointConstraint>& patchConstraints
495  ) const;
496 
497 
498  //- Disallow default bitwise copy construct
500 
501  //- Disallow default bitwise assignment
502  void operator=(const snappySnapDriver&);
503 
504 
505 public:
506 
507  //- Runtime type information
508  ClassName("snappySnapDriver");
509 
510 
511  // Constructors
512 
513  //- Construct from components
515  (
516  meshRefinement& meshRefiner,
517  const labelList& globalToMasterPatch,
518  const labelList& globalToSlavePatch
519  );
520 
521 
522  // Member Functions
523 
524  // Snapping
525 
526  //- Merge baffles.
528 
529  //- Calculate edge length per patch point.
531  (
532  const fvMesh& mesh,
533  const snapParameters& snapParams,
535  );
536 
537  //- Smooth the mesh (patch and internal) to increase visibility
538  // of surface points (on castellated mesh) w.r.t. surface.
539  static void preSmoothPatch
540  (
541  const meshRefinement& meshRefiner,
542  const snapParameters& snapParams,
543  const label nInitErrors,
544  const List<labelPair>& baffles,
546  );
547 
548  //- Get points both on patch and facezone.
550  (
551  const fvMesh& mesh,
552  const indirectPrimitivePatch&,
553  const word& zoneName
554  );
555 
556  //- Helper: calculate average cell centre per point
558  (
559  const fvMesh& mesh,
561  );
562 
563  //- Per patch point override displacement if in gap situation
564  void detectNearSurfaces
565  (
566  const scalar planarCos,
567  const indirectPrimitivePatch&,
568  const pointField& nearestPoint,
569  const vectorField& nearestNormal,
570  vectorField& disp
571  ) const;
572 
573  //- Per patch point calculate point on nearest surface. Set as
574  // boundary conditions of motionSmoother displacement field. Return
575  // displacement of patch points.
577  (
578  const meshRefinement& meshRefiner,
579  const scalarField& snapDist,
580  const indirectPrimitivePatch&,
581  pointField& nearestPoint,
582  vectorField& nearestNormal
583  );
584 
585  //- Smooth the displacement field to the internal.
586  void smoothDisplacement
587  (
588  const snapParameters& snapParams,
590  ) const;
591 
592  //- Do the hard work: move the mesh according to displacement,
593  // locally relax the displacement. Return true if ended up with
594  // correct mesh, false if not.
595  bool scaleMesh
596  (
597  const snapParameters& snapParams,
598  const label nInitErrors,
599  const List<labelPair>& baffles,
601  );
602 
603  //- Repatch faces according to surface nearest the face centre
605  (
606  const snapParameters& snapParams,
607  const labelList& adaptPatchIDs,
608  const labelList& preserveFaces
609  );
610 
611  void doSnap
612  (
613  const dictionary& snapDict,
614  const dictionary& motionDict,
615  const scalar featureCos,
616  const scalar planarAngle,
617  const snapParameters& snapParams
618  );
619 };
620 
621 
622 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
623 
624 } // End namespace Foam
625 
626 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
627 
628 #endif
629 
630 // ************************************************************************* //
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: Tuple2.H:47
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)
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 contraints.
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.
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.
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.
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:53
A class for managing temporary objects.
Definition: PtrList.H:54
void smoothDisplacement(const snapParameters &snapParams, motionSmoother &) const
Smooth the displacement field to the internal.
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.