smoothDelta.H
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2016 OpenFOAM Foundation
9 Copyright (C) 2016-2020 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27Class
28 Foam::smoothDelta
29
30Description
31 Smoothed delta which takes a given simple geometric delta and applies
32 smoothing to it such that the ratio of deltas between two cells is no
33 larger than a specified amount, typically 1.15.
34
35SourceFiles
36 smoothDelta.C
37
38\*---------------------------------------------------------------------------*/
39
40#ifndef smoothDelta_H
41#define smoothDelta_H
42
43#include "LESdelta.H"
44
45// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46
47namespace Foam
48{
49namespace LESModels
50{
51
52/*---------------------------------------------------------------------------*\
53 Class smoothDelta Declaration
54\*---------------------------------------------------------------------------*/
56class smoothDelta
57:
58 public LESdelta
59{
60public:
61
62 //- Public class used by mesh-wave to propagate the delta-ratio
63 class deltaData
64 {
65 // Private Data
66
67 scalar delta_;
68
69
70 // Private Member Functions
71
72 //- Update gets information from neighbouring face/cell and
73 //- uses this to update itself (if necessary) and return true.
74 template<class TrackingData>
75 inline bool update
76 (
77 const deltaData& w2,
78 const scalar scale,
79 const scalar tol,
80 TrackingData& td
81 );
82
83
84 public:
85
86 // Constructors
87
88 //- Default construct
89 inline deltaData();
90
91 //- Construct from delta value
92 inline deltaData(const scalar delta);
93
94
95 // Member Functions
96
97 // Access
99 scalar delta() const
100 {
101 return delta_;
102 }
103
104
105 // Needed by FaceCellWave
106
107 //- Changed or contains original (invalid) value
108 template<class TrackingData>
109 inline bool valid(TrackingData& td) const;
110
111 //- Check for identical geometrical data (eg, cyclics checking)
112 template<class TrackingData>
113 inline bool sameGeometry
114 (
115 const polyMesh&,
116 const deltaData&,
117 const scalar,
118 TrackingData& td
119 ) const;
120
121 //- Convert any absolute coordinates into relative to
122 //- (patch)face centre
123 template<class TrackingData>
124 inline void leaveDomain
125 (
126 const polyMesh&,
127 const polyPatch&,
128 const label patchFacei,
129 const point& faceCentre,
130 TrackingData& td
131 );
132
133 //- Reverse of leaveDomain
134 template<class TrackingData>
135 inline void enterDomain
136 (
137 const polyMesh&,
138 const polyPatch&,
139 const label patchFacei,
140 const point& faceCentre,
141 TrackingData& td
142 );
143
144 //- Apply rotation matrix to any coordinates
145 template<class TrackingData>
146 inline void transform
147 (
148 const polyMesh&,
149 const tensor&,
150 TrackingData& td
151 );
152
153 //- Influence of neighbouring face.
154 template<class TrackingData>
155 inline bool updateCell
156 (
157 const polyMesh&,
158 const label thisCelli,
159 const label neighbourFacei,
160 const deltaData& neighbourInfo,
161 const scalar tol,
162 TrackingData& td
163 );
164
165 //- Influence of neighbouring cell.
166 template<class TrackingData>
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 TrackingData& td
175 );
176
177 //- Influence of different value on same face.
178 template<class TrackingData>
179 inline bool updateFace
180 (
181 const polyMesh&,
182 const label thisFacei,
183 const deltaData& neighbourInfo,
184 const scalar tol,
185 TrackingData& td
186 );
187
188 //- Test for equality, with TrackingData
189 template<class TrackingData>
190 inline bool equal(const deltaData&, TrackingData& td) const;
191
192
193 // Member Operators
194
195 //- Test for equality
196 inline bool operator==(const deltaData&) const;
197
198 //- Test for inequality
199 inline bool operator!=(const deltaData&) const;
200
201
202 // IOstream Operators
204 friend Ostream& operator<<(Ostream& os, const deltaData& rhs)
205 {
206 return os << rhs.delta_;
207 }
209 friend Istream& operator>>(Istream& is, deltaData& rhs)
210 {
211 return is >> rhs.delta_;
212 }
213 };
214
215
216private:
217
218 // Private Data
219
220 autoPtr<LESdelta> geometricDelta_;
221
222 scalar maxDeltaRatio_;
223
224
225 // Private Member Functions
226
227 //- No copy construct
228 smoothDelta(const smoothDelta&) = delete;
229
230 //- No copy assignment
231 void operator=(const smoothDelta&) = delete;
232
233 // Calculate the delta values
234 void calcDelta();
235
236 //- Fill changedFaces (with face labels) and changedFacesInfo
237 // (with delta).
238 // This is the initial set of faces from which to start the waves.
239 // Since there might be lots of places with delta jumps we can follow
240 // various strategies for this initial 'seed'.
241 // - start from single cell/face and let FaceCellWave pick up all
242 // others from there. might be quite a few waves before everything
243 // settles.
244 // - start from all faces. Lots of initial transfers.
245 // We do something in between:
246 // - start from all faces where there is a jump. Since we cannot easily
247 // determine this across coupled patches (cyclic, processor)
248 // introduce all faces of these and let FaceCellWave sort it out.
249 void setChangedFaces
250 (
251 const polyMesh& mesh,
252 const volScalarField& delta,
253 DynamicList<label>& changedFaces,
254 DynamicList<deltaData>& changedFacesInfo
255 );
256
257public:
258
259 //- Declare type-name, virtual type (with debug switch)
260 TypeName("smooth");
261
262
263 // Constructors
264
265 //- Construct from name, turbulenceModel and dictionary
267 (
268 const word& name,
270 const dictionary&
271 );
272
273
274 //- Destructor
275 virtual ~smoothDelta() = default;
276
277
278 // Member Functions
279
280 //- Read the LESdelta dictionary
281 virtual void read(const dictionary& dict);
282
283 // Correct values
284 virtual void correct();
285};
286
287
288// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
289
290} // End namespace LESModels
291
292
293// * * * * * * * * * * * * * * * * * Traits * * * * * * * * * * * * * * * * //
294
295//- Contiguous data for deltaData
296template<>
297struct is_contiguous<LESModels::smoothDelta::deltaData> : std::true_type {};
298
299
300// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
301
302} // End namespace Foam
303
304// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
305
307
308// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
309
310#endif
311
312// ************************************************************************* //
scalar delta
#define w2
Definition: blockCreate.C:35
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:72
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:64
Public class used by mesh-wave to propagate the delta-ratio.
Definition: smoothDelta.H:63
bool operator==(const deltaData &) const
Test for equality.
void transform(const polyMesh &, const tensor &, TrackingData &td)
Apply rotation matrix to any coordinates.
void enterDomain(const polyMesh &, const polyPatch &, const label patchFacei, const point &faceCentre, TrackingData &td)
Reverse of leaveDomain.
bool equal(const deltaData &, TrackingData &td) const
Test for equality, with TrackingData.
bool sameGeometry(const polyMesh &, const deltaData &, const scalar, TrackingData &td) const
Check for identical geometrical data (eg, cyclics checking)
bool operator!=(const deltaData &) const
Test for inequality.
bool updateCell(const polyMesh &, const label thisCelli, const label neighbourFacei, const deltaData &neighbourInfo, const scalar tol, TrackingData &td)
Influence of neighbouring face.
bool valid(TrackingData &td) const
Changed or contains original (invalid) value.
friend Ostream & operator<<(Ostream &os, const deltaData &rhs)
Definition: smoothDelta.H:203
friend Istream & operator>>(Istream &is, deltaData &rhs)
Definition: smoothDelta.H:208
void leaveDomain(const polyMesh &, const polyPatch &, const label patchFacei, const point &faceCentre, TrackingData &td)
bool updateFace(const polyMesh &, const label thisFacei, const label neighbourCelli, const deltaData &neighbourInfo, const scalar tol, TrackingData &td)
Influence of neighbouring cell.
virtual ~smoothDelta()=default
Destructor.
TypeName("smooth")
Declare type-name, virtual type (with debug switch)
virtual void read(const dictionary &dict)
Read the LESdelta dictionary.
Definition: smoothDelta.C:176
Abstract base class for LES deltas.
Definition: LESdelta.H:54
const turbulenceModel & turbulence() const
Return turbulenceModel reference.
Definition: LESdelta.H:134
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:75
Abstract base class for turbulence models (RAS, LES and laminar).
A class for handling words, derived from Foam::string.
Definition: word.H:68
mesh update()
dynamicFvMesh & mesh
OBJstream os(runTime.globalPath()/outputName)
Namespace for OpenFOAM.
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
dictionary dict
A template class to specify that a data type can be considered as being contiguous in memory.
Definition: contiguous.H:78
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73