1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | foam-extend: Open Source CFD
4 \\ / O peration | Version: 3.2
5 \\ / A nd | Web: http://www.foam-extend.org
6 \\/ M anipulation | For copyright notice see file Copyright
7 -------------------------------------------------------------------------------
9 This file is part of foam-extend.
11 foam-extend 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 3 of the License, or (at your
14 option) any later version.
16 foam-extend is distributed in the hope that it will be useful, but
17 WITHOUT ANY WARRANTY; without even the implied warranty of
18 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 General Public License for more details.
21 You should have received a copy of the GNU General Public License
22 along with foam-extend. If not, see <http://www.gnu.org/licenses/>.
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.
35 \*---------------------------------------------------------------------------*/
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 /*---------------------------------------------------------------------------*\
48 Class smoothDelta Declaration
49 \*---------------------------------------------------------------------------*/
57 //- Raw geometric delta, to be used in smoothing
58 autoPtr<LESdelta> geometricDelta_;
61 scalar maxDeltaRatio_;
64 // Private Member Functions
66 //- Disallow default bitwise copy construct and assignment
67 smoothDelta(const smoothDelta&);
68 void operator=(const smoothDelta&);
70 // Calculate the delta values
73 //- Private member class used by mesh-wave to propagate the delta-ratio
78 // Private Member Functions
80 //- Update. Gets information from neighbouring face/cell and
81 // uses this to update itself (if nessecary) and return true.
92 // Static data members
95 static scalar maxDeltaRatio;
103 //- Construct from origin, yStar, distance
104 inline deltaData(const scalar delta);
117 // Needed by FaceCellWave
119 //- Check whether origin has been changed at all or
120 // still contains original (invalid) value.
121 inline bool valid() const;
123 //- Check for identical geometrical data.
124 // Used for cyclics checking.
125 inline bool sameGeometry
132 //- Convert any absolute coordinates into relative to
133 // (patch)face centre
134 inline void leaveDomain
138 const label patchFaceI,
139 const point& faceCentre
142 //- Reverse of leaveDomain
143 inline void enterDomain
147 const label patchFaceI,
148 const point& faceCentre
151 //- Apply rotation matrix to any coordinates
152 inline void transform
158 //- Influence of neighbouring face.
159 inline bool updateCell
162 const label thisCellI,
163 const label neighbourFaceI,
164 const deltaData& neighbourInfo,
168 //- Influence of neighbouring cell.
169 inline bool updateFace
172 const label thisFaceI,
173 const label neighbourCellI,
174 const deltaData& neighbourInfo,
178 //- Influence of different value on same face.
179 inline bool updateFace
182 const label thisFaceI,
183 const deltaData& neighbourInfo,
189 // Needed for List IO
190 inline bool operator==(const deltaData&) const;
192 inline bool operator!=(const deltaData&) const;
194 // IOstream Operators
196 friend Ostream& operator<<
199 const deltaData& wDist
202 return os << wDist.delta_;
205 friend Istream& operator>>(Istream& is, deltaData& wDist)
207 return is >> wDist.delta_;
214 const polyMesh& mesh,
215 const volScalarField& delta,
216 dynamicLabelList& changedFaces,
217 DynamicList<deltaData>& changedFacesInfo
223 //- Runtime type information
229 //- Construct from name, mesh and IOdictionary
239 virtual ~smoothDelta()
245 //- Read the LESdelta dictionary
246 virtual void read(const dictionary&);
249 virtual void correct();
253 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
255 } // End namespace Foam
257 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
259 #include "smoothDeltaDeltaDataI.H"
261 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
265 // ************************************************************************* //