initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / turbulenceModels / LES / LESdeltas / smoothDelta / smoothDelta.H
blob4b48d472f11ddecdf78df4de2264f5acb87c8cd6
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2009 OpenCFD Ltd.
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         autoPtr<LESdelta> geometricDelta_;
59         scalar maxDeltaRatio_;
62     // Private Member Functions
64         //- Disallow default bitwise copy construct and assignment
65         smoothDelta(const smoothDelta&);
66         void operator=(const smoothDelta&);
68         // Calculate the delta values
69         void calcDelta();
71         //- Private member class used by mesh-wave to propagate the delta-ratio
72         class deltaData
73         {
74             scalar delta_;
76             // Private Member Functions
78                 //- Update. Gets information from neighbouring face/cell and
79                 //  uses this to update itself (if nessecary) and return true.
80                 inline bool update
81                 (
82                     const deltaData& w2,
83                     const scalar scale,
84                     const scalar tol
85                 );
88         public:
90             // Static data members
92                 //- delta fraction
93                 static scalar maxDeltaRatio;
96                 // Constructors
98                 //- Construct null
99                 inline deltaData();
101                 //- Construct from origin, yStar, distance
102                 inline deltaData(const scalar delta);
105             // Member Functions
107                 // Access
109                 scalar delta() const
110                 {
111                     return delta_;
112                 }
115                 // Needed by FaceCellWave
117                     //- Check whether origin has been changed at all or
118                     //  still contains original (invalid) value.
119                     inline bool valid() const;
121                     //- Check for identical geometrical data.
122                     //  Used for cyclics checking.
123                     inline bool sameGeometry
124                     (
125                         const polyMesh&,
126                         const deltaData&,
127                         const scalar
128                     ) const;
130                     //- Convert any absolute coordinates into relative to
131                     //  (patch)face centre
132                     inline void leaveDomain
133                     (
134                         const polyMesh&,
135                         const polyPatch&,
136                         const label patchFaceI,
137                         const point& faceCentre
138                     );
140                     //- Reverse of leaveDomain
141                     inline void enterDomain
142                     (
143                         const polyMesh&,
144                         const polyPatch&,
145                         const label patchFaceI,
146                         const point& faceCentre
147                     );
149                     //- Apply rotation matrix to any coordinates
150                     inline void transform
151                     (
152                         const polyMesh&,
153                         const tensor&
154                     );
156                     //- Influence of neighbouring face.
157                     inline bool updateCell
158                     (
159                         const polyMesh&,
160                         const label thisCellI,
161                         const label neighbourFaceI,
162                         const deltaData& neighbourInfo,
163                         const scalar tol
164                     );
166                     //- Influence of neighbouring cell.
167                     inline bool updateFace
168                     (
169                         const polyMesh&,
170                         const label thisFaceI,
171                         const label neighbourCellI,
172                         const deltaData& neighbourInfo,
173                         const scalar tol
174                     );
176                     //- Influence of different value on same face.
177                     inline bool updateFace
178                     (
179                         const polyMesh&,
180                         const label thisFaceI,
181                         const deltaData& neighbourInfo,
182                         const scalar tol
183                     );
185                 // Member Operators
187                     // Needed for List IO
188                     inline bool operator==(const deltaData&) const;
190                     inline bool operator!=(const deltaData&) const;
192                 // IOstream Operators
194                     friend Ostream& operator<<
195                     (
196                         Ostream& os,
197                         const deltaData& wDist
198                     )
199                     {
200                         return os << wDist.delta_;
201                     }
203                     friend Istream& operator>>(Istream& is, deltaData& wDist)
204                     {
205                         return is >> wDist.delta_;
206                     }
207         };
210         void setChangedFaces
211         (
212             const polyMesh& mesh,
213             const volScalarField& delta,
214             DynamicList<label>& changedFaces,
215             DynamicList<deltaData>& changedFacesInfo
216         );
219 public:
221     //- Runtime type information
222     TypeName("smooth");
225     // Constructors
227         //- Construct from name, mesh and IOdictionary
228         smoothDelta
229         (
230             const word& name,
231             const fvMesh& mesh,
232             const dictionary&
233         );
236     //- Destructor
237     virtual ~smoothDelta()
238     {}
241     // Member Functions
243         //- Read the LESdelta dictionary
244         virtual void read(const dictionary&);
246         // Correct values
247         virtual void correct();
251 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
253 } // End namespace Foam
255 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
257 #include "smoothDeltaDeltaDataI.H"
259 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
261 #endif
263 // ************************************************************************* //