forked from o2msolutions/DFSEMInlet
-
Notifications
You must be signed in to change notification settings - Fork 0
/
turbulentDFSEMInletFvPatchVectorField.H
343 lines (253 loc) · 9.93 KB
/
turbulentDFSEMInletFvPatchVectorField.H
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
// Backport of turbulentDFSEMInlet BC to OpenFOAM 5.x
// Copyright (C) 2015 OpenFOAM Foundation
// Copyright (C) 2016-2017 OpenCFD Ltd.
// Copyright (C) 2019 Alexey Matveichev
// License: GPLv3
// Description
// Velocity boundary condition including synthesised eddies for use with LES
// and DES turbulent flows.
//
// Reference:
// \verbatim
// Poletto, R., Craft, T., and Revell, A.,
// "A New Divergence Free Synthetic Eddy Method for the Reproduction
// of Inlet Flow Conditions for LES",
// Flow Turbulence Combust (2013) 91:519-539
// \endverbatim
//
// Reynolds stress, velocity and turbulence length scale values can either
// be sepcified directly, or mapped. If mapping, the values should be
// entered in the same form as the timeVaryingMappedFixedValue condition,
// except that no interpolation in time is supported. These should be
// located in the directory:
//
// \$FOAM_CASE/constant/boundaryData/\<patchName\>/points
// \$FOAM_CASE/constant/boundaryData/\<patchName\>/0/\{R|U|L\}
//
//
// Usage
// \table
// Property | Description | Required | Default value
// value | Restart value | yes |
// delta | Local limiting length scale | yes |
// R | Reynolds stress field | no |
// U | Velocity field | no |
// L | Turbulence length scale field | no |
// d | Eddy density(fill fraction) | no | 1
// kappa | Von Karman constant | no | 0.41
// mapMethod | Method to map reference values | no | planarInterpolation
// perturb | Point perturbation for interpolation | no | 1e-5
// writeEddies | Flag to write eddies as OBJ file | no | no
// \endtable
//
// Note
// - The \c delta value typically represents a channel half-height
// - For R, U and L specification: if the entry is not user input, it is
// assumed that the data will be mapped
#ifndef turbulentDFSEMInletFvPatchVectorField_H
#define turbulentDFSEMInletFvPatchVectorField_H
#include "fixedValueFvPatchFields.H"
#include "Random.H"
#include "eddy.H"
#include "pointIndexHit.H"
#include "instantList.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
class pointToPointPlanarInterpolation;
/*---------------------------------------------------------------------------*\
Class turbulentDFSEMInletFvPatchVectorField Declaration
\*---------------------------------------------------------------------------*/
class turbulentDFSEMInletFvPatchVectorField
:
public fixedValueFvPatchVectorField
{
// Private data
//- Maximum number of attempts when seeding eddies
static label seedIterMax_;
//- Typical length scale, e.g. half channel height
const scalar delta_;
//- Ratio of sum of eddy volumes to eddy box volume; default = 1
const scalar d_;
//- Von Karman constant
const scalar kappa_;
//- Global numbering for faces
mutable autoPtr<globalIndex> globalFacesPtr_;
// Table reading for patch inlet flow properties
//- Fraction of perturbation (fraction of bounding box) to add
scalar perturb_;
//- Interpolation scheme to use
word mapMethod_;
//- 2D interpolation (for 'planarInterpolation' mapMethod)
mutable autoPtr<pointToPointPlanarInterpolation> mapperPtr_;
//- Flag to identify to interpolate the R field
bool interpolateR_;
//- Reynolds stress tensor
symmTensorField R_;
//- Flag to identify to interpolate the L field
bool interpolateL_;
//- Length scale
scalarField L_;
//- Flag to identify to interpolate the U field
bool interpolateU_;
//- Inlet velocity
vectorField U_;
//- Mean inlet velocity
vector UMean_;
// Patch information
//- Patch area - total across all processors
scalar patchArea_;
//- Decomposed patch faces as a list of triangles
faceList triFace_;
//- Addressing from per triangle to patch face
labelList triToFace_;
//- Cumulative triangle area per triangle face
scalarList triCumulativeMagSf_;
//- Cumulative area fractions per processor
scalarList sumTriMagSf_;
//- List of eddies
List<eddy> eddies_;
//- Minimum number of cells required to resolve an eddy
label nCellPerEddy_;
//- Patch normal into the domain
vector patchNormal_;
//- Eddy box volume
scalar v0_;
//- Random number generator
Random rndGen_;
//- Length scale per patch face
scalarField sigmax_;
//- Maximum length scale (across all processors)
scalar maxSigmaX_;
//- Global number of eddies
label nEddy_;
//- Current time index (used for updating)
label curTimeIndex_;
//- Patch bounds (local processor)
boundBox patchBounds_;
//- Single processor contains all eddies (flag)
bool singleProc_;
//- Flag to write the eddies to file
bool writeEddies_;
// Private Member Functions
//- Initialise info for patch point search
void initialisePatch();
//- Initialise the eddy box
void initialiseEddyBox();
//- Set a new eddy position
pointIndexHit setNewPosition(const bool global);
//- Initialise eddies
void initialiseEddies();
//- Convect the eddies
void convectEddies(const scalar deltaT);
//- Calculate the velocity fluctuation at a point
vector uDashEddy(const List<eddy>& eddies, const point& globalX) const;
//- Helper function to interpolate values from the boundary data or
// read from dictionary
template<class Type>
tmp<Field<Type>> interpolateOrRead
(
const word& fieldName,
const dictionary& dict,
bool& interpolateField
) const;
//- Helper function to interpolate values from the boundary data
template<class Type>
tmp<Field<Type>> interpolateBoundaryData
(
const word& fieldName
) const;
//- Write Lumley coefficients to file
void writeLumleyCoeffs() const;
//- Write eddy info in OBJ format
void writeEddyOBJ() const;
//- Return a reference to the patch mapper object
const pointToPointPlanarInterpolation& patchMapper() const;
//- Return eddies from remote processors that interact with local
// processor
void calcOverlappingProcEddies
(
List<List<eddy>>& overlappingEddies
) const;
public:
//- Runtime type information
TypeName("turbulentDFSEMInlet");
// Constructors
//- Construct from patch and internal field
turbulentDFSEMInletFvPatchVectorField
(
const fvPatch&,
const DimensionedField<vector, volMesh>&
);
//- Construct from patch, internal field and dictionary
turbulentDFSEMInletFvPatchVectorField
(
const fvPatch&,
const DimensionedField<vector, volMesh>&,
const dictionary&
);
//- Construct by mapping given turbulentDFSEMInletFvPatchVectorField
// onto a new patch
turbulentDFSEMInletFvPatchVectorField
(
const turbulentDFSEMInletFvPatchVectorField&,
const fvPatch&,
const DimensionedField<vector, volMesh>&,
const fvPatchFieldMapper&
);
//- Construct as copy
turbulentDFSEMInletFvPatchVectorField
(
const turbulentDFSEMInletFvPatchVectorField&
);
//- Construct and return a clone
virtual tmp<fvPatchVectorField> clone() const
{
return tmp<fvPatchVectorField>
(
new turbulentDFSEMInletFvPatchVectorField(*this)
);
}
//- Construct as copy setting internal field reference
turbulentDFSEMInletFvPatchVectorField
(
const turbulentDFSEMInletFvPatchVectorField&,
const DimensionedField<vector, volMesh>&
);
//- Construct and return a clone setting internal field reference
virtual tmp<fvPatchVectorField> clone
(
const DimensionedField<vector, volMesh>& iF
) const
{
return tmp<fvPatchVectorField>
(
new turbulentDFSEMInletFvPatchVectorField(*this, iF)
);
}
virtual ~turbulentDFSEMInletFvPatchVectorField();
// Member functions
// Mapping functions
//- Map (and resize as needed) from self given a mapping object
virtual void autoMap(const fvPatchFieldMapper& m);
//- Reverse map the given fvPatchField onto this fvPatchField
virtual void rmap
(
const fvPatchVectorField& ptf,
const labelList& addr
);
// Evaluation functions
//- Update the coefficients associated with the patch field
virtual void updateCoeffs();
//- Write
virtual void write(Ostream&) const;
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#ifdef NoRepository
#include "turbulentDFSEMInletFvPatchVectorFieldTemplates.C"
#endif
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //