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 void getZoneSurfacePoints
161  (
162  const fvMesh& mesh,
163  const indirectPrimitivePatch&,
164  const word& zoneName,
165 
166  bitSet& pointOnZone
167  );
168 
169  //- Get points both on patch and facezone.
170  template<class FaceList>
171  static labelList getFacePoints
172  (
173  const indirectPrimitivePatch& pp,
174  const FaceList& faces
175  );
176 
177  //- Per patch point calculate point on nearest surface.
178  // Return displacement of patch points.
179  static void calcNearestSurface
180  (
181  const refinementSurfaces& surfaces,
182 
183  const labelList& surfacesToTest,
184  const labelListList& regionsToTest,
185 
186  const pointField& localPoints,
187  const labelList& zonePointIndices,
188 
189  scalarField& minSnapDist,
190  labelList& snapSurf,
191  vectorField& patchDisp,
192 
193  // Optional: nearest point, normal
194  pointField& nearestPoint,
195  vectorField& nearestNormal
196  );
197 
198 
199  // Feature line snapping
200 
201  //- Is point on two feature edges that make a largish angle?
202  bool isFeaturePoint
203  (
204  const scalar featureCos,
205  const indirectPrimitivePatch& pp,
206  const bitSet& isFeatureEdge,
207  const label pointi
208  ) const;
209 
210  void smoothAndConstrain
211  (
212  const bitSet& isMasterEdge,
213  const indirectPrimitivePatch& pp,
214  const labelList& meshEdges,
215  const List<pointConstraint>& constraints,
216  vectorField& disp
217  ) const;
218  //void smoothAndConstrain2
219  //(
220  // const bool applyConstraints,
221  // const indirectPrimitivePatch& pp,
222  // const List<pointConstraint>& constraints,
223  // vectorField& disp
224  //) const;
225  void calcNearest
226  (
227  const label iter,
228  const indirectPrimitivePatch& pp,
229  vectorField& pointDisp,
230  vectorField& pointSurfaceNormal,
231  vectorField& pointRotation
232  ) const;
233  void calcNearestFace
234  (
235  const label iter,
236  const indirectPrimitivePatch& pp,
237  const scalarField& faceSnapDist,
238  vectorField& faceDisp,
239  vectorField& faceSurfaceNormal,
240  labelList& faceSurfaceRegion
241  //vectorField& faceRotation
242  ) const;
243  void calcNearestFacePointProperties
244  (
245  const label iter,
246  const indirectPrimitivePatch& pp,
247 
248  const vectorField& faceDisp,
249  const vectorField& faceSurfaceNormal,
250  const labelList& faceSurfaceRegion,
251 
252  List<List<point>>& pointFaceSurfNormals,
253  List<List<point>>& pointFaceDisp,
254  List<List<point>>& pointFaceCentres,
255  List<labelList>& pointFacePatchID
256  ) const;
257  void correctAttraction
258  (
259  const DynamicList<point>& surfacePoints,
260  const DynamicList<label>& surfaceCounts,
261  const point& edgePt,
262  const vector& edgeNormal, // normalised normal
263  const point& pt,
264  vector& edgeOffset // offset from pt to point on edge
265  ) const;
266 
267 
268  //- For any reverse (so from feature back to mesh) attraction:
269  // add attraction if diagonal points on face attracted
270  void stringFeatureEdges
271  (
272  const label iter,
273  const scalar featureCos,
274 
275  const indirectPrimitivePatch& pp,
276  const scalarField& snapDist,
277 
278  const vectorField& rawPatchAttraction,
279  const List<pointConstraint>& rawPatchConstraints,
280 
281  vectorField& patchAttraction,
282  List<pointConstraint>& patchConstraints
283  ) const;
284 
285  //- Remove constraints of points next to multi-patch points
286  // to give a bit more freedom of the mesh to conform to the
287  // multi-patch points. Bit dodgy for simple cases.
288  void releasePointsNextToMultiPatch
289  (
290  const label iter,
291  const scalar featureCos,
292 
293  const indirectPrimitivePatch& pp,
294  const scalarField& snapDist,
295 
296  const List<List<point>>& pointFaceCentres,
297  const labelListList& pointFacePatchID,
298 
299  const vectorField& rawPatchAttraction,
300  const List<pointConstraint>& rawPatchConstraints,
301 
302  vectorField& patchAttraction,
303  List<pointConstraint>& patchConstraints
304  ) const;
305 
306  //- Detect any diagonal attraction. Returns indices in face
307  // or (-1, -1) if none
308  labelPair findDiagonalAttraction
309  (
310  const indirectPrimitivePatch& pp,
311  const vectorField& patchAttraction,
312  const List<pointConstraint>& patchConstraints,
313  const label facei
314  ) const;
315 
316  scalar pyrVol
317  (
318  const indirectPrimitivePatch& pp,
319  const vectorField& featureAttraction,
320  const face& localF,
321  const point& cc
322  ) const;
323  void facePoints
324  (
325  const indirectPrimitivePatch& pp,
326  const vectorField& featureAttraction,
327  const vectorField& nearestAttraction,
328  const face& f,
330  ) const;
331  scalar pyrVol
332  (
333  const indirectPrimitivePatch& pp,
334  const vectorField& featureAttraction,
335  const vectorField& nearestAttraction,
336  const face& localF,
337  const point& cc
338  ) const;
339  Tuple2<point, vector> centreAndNormal
340  (
341  const indirectPrimitivePatch& pp,
342  const vectorField& featureAttraction,
343  const vectorField& nearestAttraction,
344  const face& localF
345  ) const;
346  bool isSplitAlignedWithFeature
347  (
348  const scalar featureCos,
349  const point& newPt0,
350  const pointConstraint& pc0,
351  const point& newPt1,
352  const pointConstraint& pc1
353  ) const;
354  bool isConcave
355  (
356  const point& c0,
357  const vector& area0,
358  const point& c1,
359  const vector& area1,
360  const scalar concaveCos
361  ) const;
362  labelPair findDiagonalAttraction
363  (
364  const scalar featureCos,
365  const scalar concaveCos,
366  const scalar minAreaFraction,
367  const indirectPrimitivePatch& pp,
368  const vectorField& patchAttraction,
369  const List<pointConstraint>& patchConstraints,
370  const vectorField& nearestAttraction,
371  const vectorField& nearestNormal,
372  const label faceI,
373 
375  DynamicField<point>& points1
376  ) const;
377 
378  //- Do all logic on whether to add face cut to diagonal
379  // attraction
380  void splitDiagonals
381  (
382  const scalar featureCos,
383  const scalar concaveCos,
384  const scalar minAreaFraction,
385 
386  const indirectPrimitivePatch& pp,
387  const vectorField& nearestAttraction,
388  const vectorField& nearestNormal,
389 
390  vectorField& patchAttraction,
391  List<pointConstraint>& patchConstraints,
392  DynamicList<label>& splitFaces,
393  DynamicList<labelPair>& splits
394  ) const;
395 
396  //- Avoid attraction across face diagonal since would
397  // cause face squeeze
398  void avoidDiagonalAttraction
399  (
400  const label iter,
401  const scalar featureCos,
402  const indirectPrimitivePatch& pp,
403  vectorField& patchAttraction,
404  List<pointConstraint>& patchConstraints
405  ) const;
406 
407  //- Write some stats about constraints
408  void writeStats
409  (
410  const indirectPrimitivePatch& pp,
411  const bitSet& isMasterPoint,
412  const List<pointConstraint>& patchConstraints
413  ) const;
414 
415  //- Return hit if on multiple points
416  pointIndexHit findMultiPatchPoint
417  (
418  const point& pt,
419  const labelList& patchIDs,
420  const List<point>& faceCentres
421  ) const;
422 
423  //- Return hit if faces-on-the-same-normalplane are on multiple
424  // patches
425  pointIndexHit findMultiPatchPoint
426  (
427  const point& pt,
428  const labelList& pfPatchID,
429  const DynamicList<vector>& surfaceNormals,
430  const labelList& faceToNormalBin
431  ) const;
432 
433  //- Return index of similar normal
434  label findNormal
435  (
436  const scalar featureCos,
437  const vector& faceSurfaceNormal,
438  const DynamicList<vector>& surfaceNormals
439  ) const;
440 
441  //- Determine attraction and constraints for single point
442  // using sampled surrounding of the point
443  void featureAttractionUsingReconstruction
444  (
445  const label iter,
446  const scalar featureCos,
447 
448  const indirectPrimitivePatch& pp,
449  const scalarField& snapDist,
450  const vectorField& nearestDisp,
451  const label pointi,
452 
453  const List<List<point>>& pointFaceSurfNormals,
454  const List<List<point>>& pointFaceDisp,
455  const List<List<point>>& pointFaceCentres,
456  const labelListList& pointFacePatchID,
457 
458  DynamicList<point>& surfacePoints,
459  DynamicList<vector>& surfaceNormals,
460  labelList& faceToNormalBin,
461 
462  vector& patchAttraction,
463  pointConstraint& patchConstraint
464  ) const;
465 
466  //- Determine attraction and constraints for all points
467  // using sampled surrounding of the point
468  void featureAttractionUsingReconstruction
469  (
470  const label iter,
471  const scalar featureCos,
472  const indirectPrimitivePatch& pp,
473  const scalarField& snapDist,
474  const vectorField& nearestDisp,
475 
476  const List<List<point>>& pointFaceSurfNormals,
477  const List<List<point>>& pointFaceDisp,
478  const List<List<point>>& pointFaceCentres,
479  const labelListList& pointFacePatchID,
480 
481  vectorField& patchAttraction,
482  List<pointConstraint>& patchConstraints
483  ) const;
484 
485  //- Determine geometric features and attraction to equivalent
486  // surface features
487  void determineFeatures
488  (
489  const label iter,
490  const scalar featureCos,
491  const bool multiRegionFeatureSnap,
492 
493  const indirectPrimitivePatch&,
494  const scalarField& snapDist,
495  const vectorField& nearestDisp,
496 
497  const List<List<point>>& pointFaceSurfNormals,
498  const List<List<point>>& pointFaceDisp,
499  const List<List<point>>& pointFaceCentres,
500  const labelListList& pointFacePatchID,
501 
502  List<labelList>& pointAttractor,
504  // Feature-edge to pp point
505  List<List<DynamicList<point>>>& edgeAttractors,
506  List<List<DynamicList<pointConstraint>>>& edgeConstraints,
507  vectorField& patchAttraction,
508  List<pointConstraint>& patchConstraints
509  ) const;
510 
511  //- Determine features originating from bafles and
512  // and add attraction to equivalent surface features
513  void determineBaffleFeatures
514  (
515  const label iter,
516  const bool baffleFeaturePoints,
517  const scalar featureCos,
518 
519  const indirectPrimitivePatch& pp,
520  const scalarField& snapDist,
521 
522  // Feature-point to pp point
523  List<labelList>& pointAttractor,
525  // Feature-edge to pp point
526  List<List<DynamicList<point>>>& edgeAttractors,
527  List<List<DynamicList<pointConstraint>>>& edgeConstraints,
528  // pp point to nearest feature
529  vectorField& patchAttraction,
530  List<pointConstraint>& patchConstraints
531  ) const;
532  void reverseAttractMeshPoints
533  (
534  const label iter,
535 
536  const indirectPrimitivePatch& pp,
537  const scalarField& snapDist,
538 
539  // Feature-point to pp point
540  const List<labelList>& pointAttractor,
542  // Feature-edge to pp point
543  const List<List<DynamicList<point>>>& edgeAttractors,
545 
546  const vectorField& rawPatchAttraction,
547  const List<pointConstraint>& rawPatchConstraints,
548 
549  // pp point to nearest feature
550  vectorField& patchAttraction,
551  List<pointConstraint>& patchConstraints
552  ) const;
553 
554  //- Find point on nearest feature edge (within searchDist).
555  // Return point and feature
556  // and store feature-edge to mesh-point and vice versa
557  Tuple2<label, pointIndexHit> findNearFeatureEdge
558  (
559  const bool isRegionEdge,
560 
561  const indirectPrimitivePatch& pp,
562  const scalarField& snapDist,
563  const label pointi,
564  const point& estimatedPt,
565 
568  vectorField&,
570  ) const;
571 
572  //- Find nearest feature point (within searchDist).
573  // Return feature point
574  // and store feature-point to mesh-point and vice versa.
575  // If another mesh point already referring to this feature
576  // point and further away, reset that one to a near feature
577  // edge (using findNearFeatureEdge above)
578  Tuple2<label, pointIndexHit> findNearFeaturePoint
579  (
580  const bool isRegionEdge,
581 
582  const indirectPrimitivePatch& pp,
583  const scalarField& snapDist,
584  const label pointi,
585  const point& estimatedPt,
586 
587  // Feature-point to pp point
588  List<labelList>& pointAttractor,
590  // Feature-edge to pp point
591  List<List<DynamicList<point>>>& edgeAttractors,
592  List<List<DynamicList<pointConstraint>>>& edgeConstraints,
593  // pp point to nearest feature
594  vectorField& patchAttraction,
595  List<pointConstraint>& patchConstraints
596  ) const;
597 
598  void featureAttractionUsingFeatureEdges
599  (
600  const label iter,
601  const bool multiRegionFeatureSnap,
602 
603  const bool detectBaffles,
604  const bool baffleFeaturePoints,
605  const bool releasePoints,
606  const bool stringFeatures,
607  const bool avoidDiagonal,
608 
609  const scalar featureCos,
610 
611  const indirectPrimitivePatch& pp,
612  const scalarField& snapDist,
613  const vectorField& nearestDisp,
614  const vectorField& nearestNormal,
615 
616  const List<List<point>>& pointFaceSurfNormals,
617  const List<List<point>>& pointFaceDisp,
618  const List<List<point>>& pointFaceCentres,
619  const labelListList& pointFacePatchID,
620 
621  vectorField& patchAttraction,
622  List<pointConstraint>& patchConstraints
623  ) const;
624 
625  void preventFaceSqueeze
626  (
627  const label iter,
628  const scalar featureCos,
629  const indirectPrimitivePatch& pp,
630  const scalarField& snapDist,
631  const vectorField& nearestAttraction,
632 
633  vectorField& patchAttraction,
634  List<pointConstraint>& patchConstraints
635  ) const;
636 
637  //- Top level feature attraction routine. Gets given
638  // displacement to nearest surface in nearestDisp
639  // and calculates new displacement taking into account
640  // features
641  vectorField calcNearestSurfaceFeature
642  (
643  const snapParameters& snapParams,
644  const bool alignMeshEdges,
645  const label iter,
646  const scalar featureCos,
647  const scalar featureAttract,
648  const scalarField& snapDist,
649  const vectorField& nearestDisp,
650  const vectorField& nearestNormal,
651  motionSmoother& meshMover,
652  vectorField& patchAttraction,
653  List<pointConstraint>& patchConstraints,
654 
655  DynamicList<label>& splitFaces,
656  DynamicList<labelPair>& splits
657  ) const;
658 
659 
660  //- No copy construct
661  snappySnapDriver(const snappySnapDriver&) = delete;
662 
663  //- No copy assignment
664  void operator=(const snappySnapDriver&) = delete;
665 
666 
667 public:
668 
669  //- Runtime type information
670  ClassName("snappySnapDriver");
671 
672 
673  // Constructors
674 
675  //- Construct from components
677  (
678  meshRefinement& meshRefiner,
679  const labelList& globalToMasterPatch,
680  const labelList& globalToSlavePatch,
681  const bool dryRun = false
682  );
683 
684 
685  // Member Functions
686 
687  // Snapping
688 
689  //- Merge baffles.
691 
692  //- Calculate edge length per patch point.
694  (
695  const fvMesh& mesh,
696  const snapParameters& snapParams,
698  );
699 
700  //- Smooth the mesh (patch and internal) to increase visibility
701  // of surface points (on castellated mesh) w.r.t. surface.
702  static void preSmoothPatch
703  (
704  const meshRefinement& meshRefiner,
705  const snapParameters& snapParams,
706  const label nInitErrors,
707  const List<labelPair>& baffles,
709  );
710 
711  //- Helper: calculate average cell centre per point
713  (
714  const fvMesh& mesh,
716  );
717 
718  //- Per patch point override displacement if in gap situation
719  void detectNearSurfaces
720  (
721  const scalar planarCos,
722  const indirectPrimitivePatch&,
723  const pointField& nearestPoint,
724  const vectorField& nearestNormal,
725  vectorField& disp
726  ) const;
727 
728  //- Per patch point calculate point on nearest surface. Set as
729  // boundary conditions of motionSmoother displacement field. Return
730  // displacement of patch points.
731  static vectorField calcNearestSurface
732  (
733  const bool strictRegionSnap,
734  const meshRefinement& meshRefiner,
735  const labelList& globalToMasterPatch,
736  const labelList& globalToSlavePatch,
737  const scalarField& snapDist,
738  const indirectPrimitivePatch&,
739  pointField& nearestPoint,
740  vectorField& nearestNormal
741  );
742 
746  //static vectorField calcNearestLocalSurface
747  //(
748  // const meshRefinement& meshRefiner,
749  // const scalarField& snapDist,
750  // const indirectPrimitivePatch&
751  //);
752 
753  //- Smooth the displacement field to the internal.
754  void smoothDisplacement
755  (
756  const snapParameters& snapParams,
758  ) const;
759 
760  //- Do the hard work: move the mesh according to displacement,
761  // locally relax the displacement. Return true if ended up with
762  // correct mesh, false if not.
763  bool scaleMesh
764  (
765  const snapParameters& snapParams,
766  const label nInitErrors,
767  const List<labelPair>& baffles,
769  );
770 
771  //- Repatch faces according to surface nearest the face centre
772  // - calculate face-wise snap distance as max of point-wise
773  // - calculate face-wise nearest surface point
774  // - repatch face according to patch for surface point.
776  (
777  const snapParameters& snapParams,
778  const labelList& adaptPatchIDs,
779  const labelList& preserveFaces
780  );
781 
782  void doSnap
783  (
784  const dictionary& snapDict,
785  const dictionary& motionDict,
786  const meshRefinement::FaceMergeType mergeType,
787  const scalar featureCos,
788  const scalar planarAngle,
789  const snapParameters& snapParams
790  );
791 };
792 
793 
794 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
795 
796 } // End namespace Foam
797 
798 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
799 
800 #ifdef NoRepository
801 # include "snappySnapDriverTemplates.C"
802 #endif
803 
804 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
805 
806 #endif
807 
808 // ************************************************************************* //
Foam::snappySnapDriver::smoothDisplacement
void smoothDisplacement(const snapParameters &snapParams, motionSmoother &) const
Smooth the displacement field to the internal.
Definition: snappySnapDriver.C:2066
Foam::snappySnapDriver::calcSnapDistance
static scalarField calcSnapDistance(const fvMesh &mesh, const snapParameters &snapParams, const indirectPrimitivePatch &)
Calculate edge length per patch point.
Definition: snappySnapDriver.C:764
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::fileName
A class for handling file names.
Definition: fileName.H:73
Foam::bitSet
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:63
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
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:133
Foam::DynamicField
Dynamically sized Field.
Definition: DynamicField.H:49
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:52
Foam::constant::physicoChemical::c1
const dimensionedScalar c1
First radiation constant: default SI units: [W/m2].
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:1051
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:123
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
meshRefinement.H
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::pointMesh
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:51
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:2207
Foam::pointConstraints
Application of (multi-)patch point constraints.
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:2140
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< label >
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:2506
Foam::Vector< scalar >
Foam::List< label >
Foam::meshRefinement
Helper class which maintains intersections of (changing) mesh with (static) surfaces.
Definition: meshRefinement.H:85
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:72
Foam::Tuple2
A 2-tuple for storing two objects of dissimilar types. The container is similar in purpose to std::pa...
Definition: stringOps.H:60
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:804
Foam::snappySnapDriver::avgCellCentres
static tmp< pointField > avgCellCentres(const fvMesh &mesh, const indirectPrimitivePatch &)
Helper: calculate average cell centre per point.
Definition: snappySnapDriver.C:966
Foam::PrimitivePatch
A list of faces which address into the list of points.
Definition: PrimitivePatch.H:79