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