Initial commit for version 2.0.x patch release
[OpenFOAM-2.0.x.git] / src / turbulenceModels / LES / LESdeltas / smoothDelta / smoothDelta.H
blob1dc93f2d68a4b3ea662b6d836df41faacbb87292
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2004-2010 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
13     the Free Software Foundation, either version 3 of the License, or
14     (at your 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, see <http://www.gnu.org/licenses/>.
24 Class
25     Foam::smoothDelta
27 Description
28     Smoothed delta which takes a given simple geometric delta and applies
29     smoothing to it such that the ratio of deltas between two cells is no
30     larger than a specified amount, typically 1.15.
32 SourceFiles
33     smoothDelta.C
35 \*---------------------------------------------------------------------------*/
37 #ifndef smoothDelta_H
38 #define smoothDelta_H
40 #include "LESdelta.H"
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44 namespace Foam
47 /*---------------------------------------------------------------------------*\
48                            Class smoothDelta Declaration
49 \*---------------------------------------------------------------------------*/
51 class smoothDelta
53     public LESdelta
55 public:
57     //- Public member class used by mesh-wave to propagate the delta-ratio
58     class deltaData
59     {
60         scalar delta_;
62         // Private Member Functions
64             //- Update. Gets information from neighbouring face/cell and
65             //  uses this to update itself (if nessecary) and return true.
66             template<class TrackingData>
67             inline bool update
68             (
69                 const deltaData& w2,
70                 const scalar scale,
71                 const scalar tol,
72                 TrackingData& td
73             );
76     public:
78         // Constructors
80             //- Construct null
81             inline deltaData();
83             //- Construct from origin, yStar, distance
84             inline deltaData(const scalar delta);
87         // Member Functions
89             // Access
91             scalar delta() const
92             {
93                 return delta_;
94             }
97             // Needed by FaceCellWave
99                 //- Check whether origin has been changed at all or
100                 //  still contains original (invalid) value.
101                 template<class TrackingData>
102                 inline bool valid(TrackingData& td) const;
104                 //- Check for identical geometrical data.
105                 //  Used for cyclics checking.
106                 template<class TrackingData>
107                 inline bool sameGeometry
108                 (
109                     const polyMesh&,
110                     const deltaData&,
111                     const scalar,
112                     TrackingData& td
113                 ) const;
115                 //- Convert any absolute coordinates into relative to
116                 //  (patch)face centre
117                 template<class TrackingData>
118                 inline void leaveDomain
119                 (
120                     const polyMesh&,
121                     const polyPatch&,
122                     const label patchFaceI,
123                     const point& faceCentre,
124                     TrackingData& td
125                 );
127                 //- Reverse of leaveDomain
128                 template<class TrackingData>
129                 inline void enterDomain
130                 (
131                     const polyMesh&,
132                     const polyPatch&,
133                     const label patchFaceI,
134                     const point& faceCentre,
135                     TrackingData& td
136                 );
138                 //- Apply rotation matrix to any coordinates
139                 template<class TrackingData>
140                 inline void transform
141                 (
142                     const polyMesh&,
143                     const tensor&,
144                     TrackingData& td
145                 );
147                 //- Influence of neighbouring face.
148                 template<class TrackingData>
149                 inline bool updateCell
150                 (
151                     const polyMesh&,
152                     const label thisCellI,
153                     const label neighbourFaceI,
154                     const deltaData& neighbourInfo,
155                     const scalar tol,
156                     TrackingData& td
157                 );
159                 //- Influence of neighbouring cell.
160                 template<class TrackingData>
161                 inline bool updateFace
162                 (
163                     const polyMesh&,
164                     const label thisFaceI,
165                     const label neighbourCellI,
166                     const deltaData& neighbourInfo,
167                     const scalar tol,
168                     TrackingData& td
169                 );
171                 //- Influence of different value on same face.
172                 template<class TrackingData>
173                 inline bool updateFace
174                 (
175                     const polyMesh&,
176                     const label thisFaceI,
177                     const deltaData& neighbourInfo,
178                     const scalar tol,
179                     TrackingData& td
180                 );
182                 //- Same (like operator==)
183                 template<class TrackingData>
184                 inline bool equal(const deltaData&, TrackingData& td) const;
186             // Member Operators
188                 // Needed for List IO
189                 inline bool operator==(const deltaData&) const;
191                 inline bool operator!=(const deltaData&) const;
193             // IOstream Operators
195                 friend Ostream& operator<<
196                 (
197                     Ostream& os,
198                     const deltaData& wDist
199                 )
200                 {
201                     return os << wDist.delta_;
202                 }
204                 friend Istream& operator>>(Istream& is, deltaData& wDist)
205                 {
206                     return is >> wDist.delta_;
207                 }
208     };
211 private:
213     // Private data
215         autoPtr<LESdelta> geometricDelta_;
216         scalar maxDeltaRatio_;
219     // Private Member Functions
221         //- Disallow default bitwise copy construct and assignment
222         smoothDelta(const smoothDelta&);
223         void operator=(const smoothDelta&);
225         // Calculate the delta values
226         void calcDelta();
229         void setChangedFaces
230         (
231             const polyMesh& mesh,
232             const volScalarField& delta,
233             DynamicList<label>& changedFaces,
234             DynamicList<deltaData>& changedFacesInfo
235         );
238 public:
240     //- Runtime type information
241     TypeName("smooth");
244     // Constructors
246         //- Construct from name, mesh and IOdictionary
247         smoothDelta
248         (
249             const word& name,
250             const fvMesh& mesh,
251             const dictionary&
252         );
255     //- Destructor
256     virtual ~smoothDelta()
257     {}
260     // Member Functions
262         //- Read the LESdelta dictionary
263         virtual void read(const dictionary&);
265         // Correct values
266         virtual void correct();
270 //- Data associated with deltaData type are contiguous
271 template<>
272 inline bool contiguous<smoothDelta::deltaData>()
274     return true;
278 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
280 } // End namespace Foam
282 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
284 #include "smoothDeltaDeltaDataI.H"
286 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
288 #endif
290 // ************************************************************************* //