Merge remote-tracking branch 'origin/nr/multiSolverFix' into nextRelease
[foam-extend-3.2.git] / src / turbulenceModels / LES / LESdeltas / smoothDelta / smoothDelta.H
blob709c5605c7a3f3871d577d794a6415412c6ed5f0
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright held by original author
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
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 the
13     Free Software Foundation; either version 2 of the License, or (at your
14     option) any later version.
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.
21     You should have received a copy of the GNU General Public License
22     along with OpenFOAM; if not, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 Class
26     Foam::smoothDelta
28 Description
29     Smoothed delta which takes a given simple geometric delta and applies
30     smoothing to it such that the ratio of deltas between two cells is no
31     larger than a specified amount, typically 1.15.
33 SourceFiles
34     smoothDelta.C
36 \*---------------------------------------------------------------------------*/
38 #ifndef smoothDelta_H
39 #define smoothDelta_H
41 #include "LESdelta.H"
43 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45 namespace Foam
48 /*---------------------------------------------------------------------------*\
49                            Class smoothDelta Declaration
50 \*---------------------------------------------------------------------------*/
52 class smoothDelta
54     public LESdelta
56     // Private data
58         //- Raw geometric delta, to be used in smoothing
59         autoPtr<LESdelta> geometricDelta_;
61         //- Smoothing ratio
62         scalar maxDeltaRatio_;
65     // Private Member Functions
67         //- Disallow default bitwise copy construct and assignment
68         smoothDelta(const smoothDelta&);
69         void operator=(const smoothDelta&);
71         // Calculate the delta values
72         void calcDelta();
74         //- Private member class used by mesh-wave to propagate the delta-ratio
75         class deltaData
76         {
77             scalar delta_;
79             // Private Member Functions
81                 //- Update. Gets information from neighbouring face/cell and
82                 //  uses this to update itself (if nessecary) and return true.
83                 inline bool update
84                 (
85                     const deltaData& w2,
86                     const scalar scale,
87                     const scalar tol
88                 );
91         public:
93             // Static data members
95                 //- delta fraction
96                 static scalar maxDeltaRatio;
99                 // Constructors
101                 //- Construct null
102                 inline deltaData();
104                 //- Construct from origin, yStar, distance
105                 inline deltaData(const scalar delta);
108             // Member Functions
110                 // Access
112                 scalar delta() const
113                 {
114                     return delta_;
115                 }
118                 // Needed by FaceCellWave
120                     //- Check whether origin has been changed at all or
121                     //  still contains original (invalid) value.
122                     inline bool valid() const;
124                     //- Check for identical geometrical data.
125                     //  Used for cyclics checking.
126                     inline bool sameGeometry
127                     (
128                         const polyMesh&,
129                         const deltaData&,
130                         const scalar
131                     ) const;
133                     //- Convert any absolute coordinates into relative to
134                     //  (patch)face centre
135                     inline void leaveDomain
136                     (
137                         const polyMesh&,
138                         const polyPatch&,
139                         const label patchFaceI,
140                         const point& faceCentre
141                     );
143                     //- Reverse of leaveDomain
144                     inline void enterDomain
145                     (
146                         const polyMesh&,
147                         const polyPatch&,
148                         const label patchFaceI,
149                         const point& faceCentre
150                     );
152                     //- Apply rotation matrix to any coordinates
153                     inline void transform
154                     (
155                         const polyMesh&,
156                         const tensor&
157                     );
159                     //- Influence of neighbouring face.
160                     inline bool updateCell
161                     (
162                         const polyMesh&,
163                         const label thisCellI,
164                         const label neighbourFaceI,
165                         const deltaData& neighbourInfo,
166                         const scalar tol
167                     );
169                     //- Influence of neighbouring cell.
170                     inline bool updateFace
171                     (
172                         const polyMesh&,
173                         const label thisFaceI,
174                         const label neighbourCellI,
175                         const deltaData& neighbourInfo,
176                         const scalar tol
177                     );
179                     //- Influence of different value on same face.
180                     inline bool updateFace
181                     (
182                         const polyMesh&,
183                         const label thisFaceI,
184                         const deltaData& neighbourInfo,
185                         const scalar tol
186                     );
188                 // Member Operators
190                     // Needed for List IO
191                     inline bool operator==(const deltaData&) const;
193                     inline bool operator!=(const deltaData&) const;
195                 // IOstream Operators
197                     friend Ostream& operator<<
198                     (
199                         Ostream& os,
200                         const deltaData& wDist
201                     )
202                     {
203                         return os << wDist.delta_;
204                     }
206                     friend Istream& operator>>(Istream& is, deltaData& wDist)
207                     {
208                         return is >> wDist.delta_;
209                     }
210         };
213         void setChangedFaces
214         (
215             const polyMesh& mesh,
216             const volScalarField& delta,
217             DynamicList<label>& changedFaces,
218             DynamicList<deltaData>& changedFacesInfo
219         );
222 public:
224     //- Runtime type information
225     TypeName("smooth");
228     // Constructors
230         //- Construct from name, mesh and IOdictionary
231         smoothDelta
232         (
233             const word& name,
234             const fvMesh& mesh,
235             const dictionary&
236         );
239     //- Destructor
240     virtual ~smoothDelta()
241     {}
244     // Member Functions
246         //- Read the LESdelta dictionary
247         virtual void read(const dictionary&);
249         // Correct values
250         virtual void correct();
254 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
256 } // End namespace Foam
258 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
260 #include "smoothDeltaDeltaDataI.H"
262 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
264 #endif
266 // ************************************************************************* //