ENH: patchCloud: return pTraits<Type>::max for unfound points
[OpenFOAM-1.7.x.git] / src / turbulenceModels / LES / LESdeltas / smoothDelta / smoothDelta.H
blob627c23801ec3ec401796dee308b119521b0eda4b
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-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     // Private data
57         autoPtr<LESdelta> geometricDelta_;
58         scalar maxDeltaRatio_;
61     // Private Member Functions
63         //- Disallow default bitwise copy construct and assignment
64         smoothDelta(const smoothDelta&);
65         void operator=(const smoothDelta&);
67         // Calculate the delta values
68         void calcDelta();
70         //- Private member class used by mesh-wave to propagate the delta-ratio
71         class deltaData
72         {
73             scalar delta_;
75             // Private Member Functions
77                 //- Update. Gets information from neighbouring face/cell and
78                 //  uses this to update itself (if nessecary) and return true.
79                 inline bool update
80                 (
81                     const deltaData& w2,
82                     const scalar scale,
83                     const scalar tol
84                 );
87         public:
89             // Static data members
91                 //- delta fraction
92                 static scalar maxDeltaRatio;
95                 // Constructors
97                 //- Construct null
98                 inline deltaData();
100                 //- Construct from origin, yStar, distance
101                 inline deltaData(const scalar delta);
104             // Member Functions
106                 // Access
108                 scalar delta() const
109                 {
110                     return delta_;
111                 }
114                 // Needed by FaceCellWave
116                     //- Check whether origin has been changed at all or
117                     //  still contains original (invalid) value.
118                     inline bool valid() const;
120                     //- Check for identical geometrical data.
121                     //  Used for cyclics checking.
122                     inline bool sameGeometry
123                     (
124                         const polyMesh&,
125                         const deltaData&,
126                         const scalar
127                     ) const;
129                     //- Convert any absolute coordinates into relative to
130                     //  (patch)face centre
131                     inline void leaveDomain
132                     (
133                         const polyMesh&,
134                         const polyPatch&,
135                         const label patchFaceI,
136                         const point& faceCentre
137                     );
139                     //- Reverse of leaveDomain
140                     inline void enterDomain
141                     (
142                         const polyMesh&,
143                         const polyPatch&,
144                         const label patchFaceI,
145                         const point& faceCentre
146                     );
148                     //- Apply rotation matrix to any coordinates
149                     inline void transform
150                     (
151                         const polyMesh&,
152                         const tensor&
153                     );
155                     //- Influence of neighbouring face.
156                     inline bool updateCell
157                     (
158                         const polyMesh&,
159                         const label thisCellI,
160                         const label neighbourFaceI,
161                         const deltaData& neighbourInfo,
162                         const scalar tol
163                     );
165                     //- Influence of neighbouring cell.
166                     inline bool updateFace
167                     (
168                         const polyMesh&,
169                         const label thisFaceI,
170                         const label neighbourCellI,
171                         const deltaData& neighbourInfo,
172                         const scalar tol
173                     );
175                     //- Influence of different value on same face.
176                     inline bool updateFace
177                     (
178                         const polyMesh&,
179                         const label thisFaceI,
180                         const deltaData& neighbourInfo,
181                         const scalar tol
182                     );
184                 // Member Operators
186                     // Needed for List IO
187                     inline bool operator==(const deltaData&) const;
189                     inline bool operator!=(const deltaData&) const;
191                 // IOstream Operators
193                     friend Ostream& operator<<
194                     (
195                         Ostream& os,
196                         const deltaData& wDist
197                     )
198                     {
199                         return os << wDist.delta_;
200                     }
202                     friend Istream& operator>>(Istream& is, deltaData& wDist)
203                     {
204                         return is >> wDist.delta_;
205                     }
206         };
209         void setChangedFaces
210         (
211             const polyMesh& mesh,
212             const volScalarField& delta,
213             DynamicList<label>& changedFaces,
214             DynamicList<deltaData>& changedFacesInfo
215         );
218 public:
220     //- Runtime type information
221     TypeName("smooth");
224     // Constructors
226         //- Construct from name, mesh and IOdictionary
227         smoothDelta
228         (
229             const word& name,
230             const fvMesh& mesh,
231             const dictionary&
232         );
235     //- Destructor
236     virtual ~smoothDelta()
237     {}
240     // Member Functions
242         //- Read the LESdelta dictionary
243         virtual void read(const dictionary&);
245         // Correct values
246         virtual void correct();
250 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
252 } // End namespace Foam
254 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
256 #include "smoothDeltaDeltaDataI.H"
258 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
260 #endif
262 // ************************************************************************* //