-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathsol_eddy_visco.py
354 lines (283 loc) · 10.7 KB
/
sol_eddy_visco.py
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
344
345
346
347
348
349
350
351
352
353
354
"""
In turbulence_models.py several models, including SA, k-epsilon, k-omega SST,
v2-f, and cess are implemented for RANS simulation
"""
import numpy as np
from subfunc.subfuncs import sol_eqn
def kepXu(mesh, u, k, e, vv, r, mu, ReTau):
n = mesh.nPoints
d = np.minimum(mesh.y, mesh.y[-1] - mesh.y) #wall distance, y_n
# yplus = d * retau
ystar = d * np.sqrt(k)/mu
# Model constants
# cmu = 0.09
sigk = 1.0
sige = 1.3
Ce1 = 1.44
Ce2 = 1.92
## vv
vv =k *((7.19e-3)*ystar-(4.33e-5)*np.power(ystar,2.0)+(8.8e-8)*np.power(ystar,3.0))
ystar_vv=np.sqrt(vv)*d/mu
# Model functions
lmu=0.544*d/(1+5.025e-4*pow(ystar_vv,1.65))
leps=8.8*d/(1+10/ystar_vv +5.15e-2*ystar_vv)
mut =np.sqrt(vv)*lmu
e=np.sqrt(vv)*k/leps
# Turbulent production: Pk = mut*dudy^2
Pk = mut * np.power(mesh.ddy @ u, 2)
# ---------------------------------------------------------------------
# e-equation
mueff = mu + mut / sige
fs = fd = f2 = np.ones(n)
# diffusion matrix: mueff*d2()/dy2 + dmueff/dy d()/dy
A = np.einsum('i,ij->ij', mueff * fd, mesh.d2dy2) \
+ np.einsum('i,ij->ij', (mesh.ddy @ mueff) * fd, mesh.ddy)
# Left-hand-side, implicitly treated source term
np.fill_diagonal(A, A.diagonal() - Ce2 * f2 * r * e / k / fs)
# Right-hand-side
b = -e[1:-1] / k[1:-1] * Ce1 * Pk[1:-1]
# Wall boundary conditions
e[0] = mu[0] / r[0] * k[1] / np.power(d[1], 2)
e[-1] = mu[-1] / r[-1] * k[-2] / np.power(d[-2], 2)
# Solve eps equation
e = sol_eqn(e * fs, A, b, 0.8) / fs
e[1:-1] = np.maximum(e[1:-1], 1.e-12)
# ---------------------------------------------------------------------
# k-equation
# effective viscosity
mueff = mu + mut / sigk
fs = fd = np.ones(n)
# diffusion matrix: mueff*d2()/dy2 + dmueff/dy d()/dy
A = np.einsum('i,ij->ij', mueff * fd, mesh.d2dy2) \
+ np.einsum('i,ij->ij', (mesh.ddy @ mueff) * fd, mesh.ddy)
# implicitly treated source term
np.fill_diagonal(A, A.diagonal() - r * e / k / fs)
# Right-hand-side
b = -Pk[1:-1]
# Wall boundary conditions
k[0] = k[-1] = 0.0
# Solve TKE
k = sol_eqn(k * fs, A, b, 0.7) / fs
k[1:-1] = np.maximum(k[1:-1], 1.e-12)
return mut, k, e, Pk
def kep(mesh, u, k, e, r, mu):
n = mesh.nPoints
d = np.minimum(mesh.y, mesh.y[-1] - mesh.y)
# yplus = d * retau
# Model constants
cmu = 0.09
sigk = 1.0
sige = 1.3
Ce1 = 1.44
Ce2 = 1.92
# Model functions
mut = cmu * r / e * np.power(k, 2)
mut[1:-1] = np.minimum(np.maximum(mut[1:-1], 1.0e-10), 100.0)
# Turbulent production: Pk = mut*dudy^2
Pk = mut * np.power(mesh.ddy @ u, 2)
# ---------------------------------------------------------------------
# e-equation
# effective viscosity
mueff = mu + mut / sige
fs = fd = np.ones(n)
# diffusion matrix: mueff*d2()/dy2 + dmueff/dy d()/dy
A = np.einsum('i,ij->ij', mueff * fd, mesh.d2dy2) \
+ np.einsum('i,ij->ij', (mesh.ddy @ mueff) * fd, mesh.ddy)
# Left-hand-side, implicitly treated source term
np.fill_diagonal(A, A.diagonal() - Ce2 * r * e / k / fs)
# Right-hand-side
b = -e[1:-1] / k[1:-1] * Ce1 * Pk[1:-1]
# Wall boundary conditions
e[0] = mu[0] / r[0] * k[1] / np.power(d[1], 2)
e[-1] = mu[-1] / r[-1] * k[-2] / np.power(d[-2], 2)
# Solve eps equation
e = sol_eqn(e * fs, A, b, 0.8) / fs
e[1:-1] = np.maximum(e[1:-1], 1.e-12)
# ---------------------------------------------------------------------
# k-equation
# effective viscosity
mueff = mu + mut / sigk
fs = fd = np.ones(n)
# diffusion matrix: mueff*d2()/dy2 + dmueff/dy d()/dy
A = np.einsum('i,ij->ij', mueff * fd, mesh.d2dy2) \
+ np.einsum('i,ij->ij', (mesh.ddy @ mueff) * fd, mesh.ddy)
# implicitly treated source term
np.fill_diagonal(A, A.diagonal() - r * e / k / fs)
# Right-hand-side
b = -Pk[1:-1]
# Wall boundary conditions
k[0] = k[-1] = 0.0
# Solve TKE
k = sol_eqn(k * fs, A, b, 0.7) / fs
k[1:-1] = np.maximum(k[1:-1], 1.e-12)
return mut, k, e, Pk
def sst(mesh, u, k, om, r, mu ):
n = mesh.nPoints
# model constants
sigma_k1 = 0.85
sigma_k2 = 1.0
sigma_om1 = 0.5
sigma_om2 = 0.856
beta_1 = 0.075
beta_2 = 0.0828
betaStar = 0.09
a1 = 0.31
# alfa_1 = beta_1 / betaStar - sigma_om1 * 0.41 ** 2.0 / betaStar ** 0.5
# alfa_2 = beta_2 / betaStar - sigma_om2 * 0.41 ** 2.0 / betaStar ** 0.5
alfa_1 = 5.0/9.0 #beta_1 / betaStar - sigma_om1 * 0.41 ** 2.0 / betaStar ** 0.5
alfa_2 = 0.44 #beta_2 / betaStar - sigma_om2 * 0.41 ** 2.0 / betaStar ** 0.5
# Relaxation factors
underrelaxK = 0.6
underrelaxOm = 0.4
# required gradients
dkdy = mesh.ddy @ k
domdy = mesh.ddy @ om
wallDist = np.minimum(mesh.y, mesh.y[-1] - mesh.y)
wallDist = np.maximum(wallDist, 1.0e-8)
# VortRate = StrainRate in fully developed channel
strMag = np.absolute(mesh.ddy @ u)
# Blending functions
CDkom = 2.0 * sigma_om2 * r / om * dkdy * domdy
gamma1 = 500.0 * mu / (r * om * wallDist * wallDist)
gamma2 = 4.0 * sigma_om2 * r * k / (wallDist * wallDist * np.maximum(CDkom, 1.0e-20))
gamma3 = np.sqrt(k) / (betaStar * om * wallDist)
gamma = np.minimum(np.maximum(gamma1, gamma3), gamma2)
bF1 = np.tanh(np.power(gamma, 4.0))
gamma = np.maximum(2.0 * gamma3, gamma1)
bF2 = np.tanh(np.power(gamma, 2.0))
# more model constants
alfa = alfa_1 * bF1 + (1 - bF1) * alfa_2
beta = beta_1 * bF1 + (1 - bF1) * beta_2
sigma_k = sigma_k1 * bF1 + (1 - bF1) * sigma_k2
sigma_om = sigma_om1 * bF1 + (1 - bF1) * sigma_om2
# Eddy viscosity
zeta = np.minimum(1.0 / om, a1 / (strMag * bF2))
mut = r * k * zeta
mut = np.minimum(np.maximum(mut, 0.0), 100.0)
# ---------------------------------------------------------------------
# om-equation
# effective viscosity
mueff = mu + sigma_om * mut
fs = np.ones(n)
# diffusion matrix: mueff*d2()/dy2 + dmueff/dy d()/dy
A = np.einsum('i,ij->ij', mueff, mesh.d2dy2) \
+ np.einsum('i,ij->ij', mesh.ddy @ mueff, mesh.ddy)
# implicitly treated source term
np.fill_diagonal(A, A.diagonal() - beta * r * om / fs)
# Right-hand-side
b = -alfa[1:-1] * r[1:-1] * strMag[1:-1] * strMag[1:-1] - (1 - bF1[1:-1]) * CDkom[1:-1]
# Wall boundary conditions
om[0] = 60.0 * mu[0] / beta_1 / r[0] / wallDist[2] / wallDist[1]
om[-1] = 60.0 * mu[-1] / beta_1 / r[-1] / wallDist[-2] / wallDist[-2]
# Solve
om = sol_eqn(om * fs, A, b, underrelaxOm) / fs
om[1:-1] = np.maximum(om[1:-1], 1.e-12)
# ---------------------------------------------------------------------
# k-equation
# effective viscosity
mueff = mu + sigma_k * mut
fs = fd = np.ones(n)
# diffusion matrix: mueff*d2()/dy2 + dmueff/dy d()/dy
A = np.einsum('i,ij->ij', mueff * fd, mesh.d2dy2) \
+ np.einsum('i,ij->ij', (mesh.ddy @ mueff) * fd, mesh.ddy)
# implicitly treated source term
np.fill_diagonal(A, A.diagonal() - betaStar * r * om / fs)
# Right-hand-side
Pk = np.minimum(mut * strMag * strMag, 20 * betaStar * k * r * om)
b = -Pk[1:-1]
# Wall boundary conditions
k[0] = k[-1] = 0.0
# Solve
k = sol_eqn(k * fs, A, b, underrelaxK) / fs
k[1:-1] = np.maximum(k[1:-1], 1.e-12)
e = betaStar*k*om
return mut, k, e, Pk, om
def v2f(mesh, u, k, e, v2, r, mu):
n = mesh.nPoints
f = np.zeros(n)
# Model constants
cmu = 0.22
sigk = 1.0
sige = 1.3
Ce2 = 1.9
Ct = 6
Cl = 0.23
Ceta = 70
C1 = 1.4
# C1 = 2.2
C2 = 0.3
# Relaxation factors
underrelaxK = 0.6
underrelaxE = 0.6
underrelaxV2 = 0.6
# Time and length scales, eddy viscosity and turbulent production
Tt = np.maximum(k / e, Ct * np.power(mu / (r * e), 0.5))
Lt = Cl * np.maximum(np.power(k, 1.5) / e, Ceta * np.power(np.power(mu / r, 3) / e, 0.25))
mut = np.maximum(cmu * r * v2 * Tt, 0.0)
Pk = mut * np.power(mesh.ddy @ u, 2.0)
# ---------------------------------------------------------------------
# f-equation
# implicitly treated source term
A = np.einsum('i,ij->ij', Lt * Lt, mesh.d2dy2)
np.fill_diagonal(A, A.diagonal() - 1.0)
# Right-hand-side
vok = v2[1:-1] / k[1:-1]
rhsf = ((C1 - 6) * vok - 2 / 3 * (C1 - 1)) / Tt[1:-1] - C2 * Pk[1:-1] / (r[1:-1] * k[1:-1])
# Solve
f = sol_eqn(f, A, rhsf, 1)
f[1:-1] = np.maximum(f[1:-1], 1.e-12)
# ---------------------------------------------------------------------
# v2-equation:
# effective viscosity and pre-factors for compressibility implementation
mueff = mu + mut
fs = fd = np.ones(n)
# diffusion matrix: mueff*d2()/dy2 + dmueff/dy d()/dy
A = np.einsum('i,ij->ij', mueff * fd, mesh.d2dy2) \
+ np.einsum('i,ij->ij', (mesh.ddy @ mueff) * fd, mesh.ddy)
# implicitly treated source term
np.fill_diagonal(A, A.diagonal() - 6.0 * r * e / k / fs)
# Right-hand-side
b = -r[1:-1] * k[1:-1] * f[1:-1]
# Wall boundary conditions
v2[0] = v2[-1] = 0.0
# Solve
v2 = sol_eqn(v2 * fs, A, b, underrelaxV2) / fs
v2[1:-1] = np.maximum(v2[1:-1], 1.e-12)
# ---------------------------------------------------------------------
# e-equation
# effective viscosity
mueff = mu + mut / sige
fs = np.ones(n)
fd = np.ones(n)
# diffusion matrix: mueff*d2()/dy2 + dmueff/dy d()/dy
A = np.einsum('i,ij->ij', mueff * fd, mesh.d2dy2) \
+ np.einsum('i,ij->ij', (mesh.ddy @ mueff) * fd, mesh.ddy)
# implicitly treated source term
np.fill_diagonal(A, A.diagonal() - Ce2 / Tt * r / fs)
# Right-hand-side
Ce1 = 1.4 * (1 + 0.045 * np.sqrt(k[1:-1] / v2[1:-1]))
b = -1 / Tt[1:-1] * Ce1 * Pk[1:-1]
# Wall boundary conditions
e[0] = mu[0] * k[1] / r[0] / np.power(mesh.y[1] - mesh.y[0], 2)
e[-1] = mu[-1] * k[-2] / r[-1] / np.power(mesh.y[-1] - mesh.y[-2], 2)
# Solve
e = sol_eqn(e * fs, A, b, underrelaxE) / fs
e[1:-1] = np.maximum(e[1:-1], 1.e-12)
# ---------------------------------------------------------------------
# k-equation
# effective viscosity
mueff = mu + mut / sigk
fs = fd = np.ones(n)
# diffusion matrix: mueff*d2()/dy2 + dmueff/dy d()/dy
A = np.einsum('i,ij->ij', mueff * fd, mesh.d2dy2) \
+ np.einsum('i,ij->ij', (mesh.ddy @ mueff) * fd, mesh.ddy)
# implicitly treated source term
np.fill_diagonal(A, A.diagonal() - r * e / k / fs)
# Right-hand-side
b = -Pk[1:-1]
# Wall boundary conditions
k[0] = k[-1] = 0.0
# Solve
k = sol_eqn(k * fs, A, b, underrelaxK) / fs
k[1:-1] = np.maximum(k[1:-1], 1.e-12)
return mut, k, e, Pk, v2