Loading...
Searching...
No Matches
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,2025 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{
52/*---------------------------------------------------------------------------*\
53 Class smoothDelta Declaration
54\*---------------------------------------------------------------------------*/
55
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
98
99 scalar delta() const
100 {
101 return delta_;
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 //- Interpolate between two values (lerp). Returns true if
193 //- causes changes. Not sure if needs to be specialised between
194 //- face and cell and what index is needed...
195 template<class TrackingData>
196 inline bool interpolate
197 (
198 const polyMesh&,
199 const point& pt,
200 const label i0,
201 const deltaData& f0,
202 const label i1,
203 const deltaData& f1,
204 const scalar weight,
205 const scalar tol,
206 TrackingData& td
207 );
208
209
210 // Member Operators
211
212 //- Test for equality
213 inline bool operator==(const deltaData&) const;
214
215 //- Test for inequality
216 inline bool operator!=(const deltaData&) const;
217
218
219 // IOstream Operators
220
221 friend Ostream& operator<<(Ostream& os, const deltaData& rhs)
222 {
223 return os << rhs.delta_;
224 }
225
226 friend Istream& operator>>(Istream& is, deltaData& rhs)
227 {
228 return is >> rhs.delta_;
229 }
230 };
231
232
233private:
234
235 // Private Data
236
237 autoPtr<LESdelta> geometricDelta_;
238
239 scalar maxDeltaRatio_;
240
241
242 // Private Member Functions
243
244 //- No copy construct
245 smoothDelta(const smoothDelta&) = delete;
246
247 //- No copy assignment
248 void operator=(const smoothDelta&) = delete;
249
250 // Calculate the delta values
251 void calcDelta();
252
253 //- Fill changedFaces (with face labels) and changedFacesInfo
254 // (with delta).
255 // This is the initial set of faces from which to start the waves.
256 // Since there might be lots of places with delta jumps we can follow
257 // various strategies for this initial 'seed'.
258 // - start from single cell/face and let FaceCellWave pick up all
259 // others from there. might be quite a few waves before everything
260 // settles.
261 // - start from all faces. Lots of initial transfers.
262 // We do something in between:
263 // - start from all faces where there is a jump. Since we cannot easily
264 // determine this across coupled patches (cyclic, processor)
265 // introduce all faces of these and let FaceCellWave sort it out.
266 void setChangedFaces
267 (
268 const polyMesh& mesh,
269 const volScalarField& delta,
270 DynamicList<label>& changedFaces,
271 DynamicList<deltaData>& changedFacesInfo
272 );
273
274public:
275
276 //- Declare type-name, virtual type (with debug switch)
277 TypeName("smooth");
278
279
280 // Constructors
281
282 //- Construct from name, turbulenceModel and dictionary
283 smoothDelta
284 (
285 const word& name,
287 const dictionary&
288 );
289
290
291 //- Destructor
292 virtual ~smoothDelta() = default;
293
294
295 // Member Functions
296
297 //- Read the LESdelta dictionary
298 virtual void read(const dictionary& dict);
299
300 // Correct values
301 virtual void correct();
302};
303
304
305// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
306
307} // End namespace LESModels
308
309
310// * * * * * * * * * * * * * * * * * Traits * * * * * * * * * * * * * * * * //
311
312//- Contiguous data for deltaData
313template<>
314struct is_contiguous<LESModels::smoothDelta::deltaData> : std::true_type {};
315
316//- Interpolation - used in e.g. cyclicAMI interpolation
318(
321 const scalar t
322)
323{
324 return LESModels::smoothDelta::deltaData(lerp(a.delta(), b.delta(), t));
325}
326
327
328// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
329
330} // End namespace Foam
331
332// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
333
335
336// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
337
338#endif
339
340// ************************************************************************* //
scalar delta
#define w2
Definition blockCreate.C:28
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition Istream.H:60
Public class used by mesh-wave to propagate the delta-ratio.
Definition smoothDelta.H:61
bool operator==(const deltaData &) const
Test for equality.
bool interpolate(const polyMesh &, const point &pt, const label i0, const deltaData &f0, const label i1, const deltaData &f1, const scalar weight, const scalar tol, TrackingData &td)
Interpolate between two values (lerp). Returns true if causes changes. Not sure if needs to be specia...
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)
friend Istream & operator>>(Istream &is, deltaData &rhs)
void leaveDomain(const polyMesh &, const polyPatch &, const label patchFacei, const point &faceCentre, TrackingData &td)
Convert any absolute coordinates into relative to (patch)face centre.
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.
LESdelta(const LESdelta &)=delete
No copy construct.
const turbulenceModel & turbulence() const
Return turbulenceModel reference.
Definition LESdelta.H:146
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
Smoothed delta which takes a given simple geometric delta and applies smoothing to it such that the r...
Abstract base class for turbulence models (RAS, LES and laminar).
A class for handling words, derived from Foam::string.
Definition word.H:66
mesh update()
dynamicFvMesh & mesh
OBJstream os(runTime.globalPath()/outputName)
wallPoints::trackData td(isBlockedFace, regionToBlockSize)
Namespace for LES SGS models.
Namespace for OpenFOAM.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Tensor< scalar > tensor
Definition symmTensor.H:57
vector point
Point is a vector.
Definition point.H:37
void rhs(fvMatrix< typename Expr::value_type > &m, const Expr &expression)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
dimensioned< Type > lerp(const dimensioned< Type > &a, const dimensioned< Type > &b, const scalar t)
dictionary dict
volScalarField & b
A template class to specify that a data type can be considered as being contiguous in memory.
Definition contiguous.H:70
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition typeInfo.H:68