-
Notifications
You must be signed in to change notification settings - Fork 5
/
Copy pathsub_turbulence_BL.f90
337 lines (276 loc) · 9.13 KB
/
sub_turbulence_BL.f90
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
!------------------------------------------------------------------------------
subroutine turbulence_model_BL(nMesh,mBlock)
Use Global_Var
Use Flow_Var
implicit none
integer:: mBlock,nx,ny,nz,ksub,i,j,k,kflag,i1,j1,k1,i2,j2,k2,i0,j0,k0,nMesh
real(PRE_EC):: ui,vi,wi,uj,vj,wj,uk,vk,wk,ux,vx,wx,uy,vy,wy,uz,vz,wz
real(PRE_EC):: ix,iy,iz,jx,jy,jz,kx,ky,kz,x0,y0,z0
real(PRE_EC),allocatable,dimension(:,:,:):: omiga
real(PRE_EC),allocatable,dimension(:):: Amu1d,Amut1d,d1d,u1d,yy,omiga1d
integer,allocatable:: flag1(:,:,:)
Type (Block_TYPE),pointer:: B
Type (BC_MSG_TYPE),pointer:: Bc
B => Mesh(nMesh)%Block(mBlock)
nx=B%nx ; ny=B%ny; nz=B%nz
B%mu_t(:,:,:)=0.d0
! test if the block contains wall
kflag=0
do ksub=1, B%subface
if(B%bc_msg(ksub)%bc .eq. BC_WALL) kflag=1
enddo
if(Kflag .eq. 0) return ! No wall in this block
! This Block Contains Wall
allocate(omiga(0:nx,0:ny,0:nz),flag1(0:nx,0:ny,0:nz))
omiga=0.d0
flag1=0
!----- get Omiga (vorticity) omiga=sqrt(omigax**2+omigay**2+omigaz**2) at the cell's center ------
! 采用中心差分求解
!$OMP PARALLEL DO DEFAULT(PRIVATE) SHARED (nx,ny,nz,uu,v,w,omiga,B)
do k=1,nz-1
do j=1,ny-1
do i=1,nx-1
! 物理量对于计算坐标(下标)的导数
ui=uu(i+1,j,k)-uu(i-1,j,k)
vi=v(i+1,j,k)-v(i-1,j,k)
wi=w(i+1,j,k)-w(i-1,j,k)
uj=uu(i,j+1,k)-uu(i,j-1,k)
vj=v(i,j+1,k)-v(i,j-1,k)
wj=w(i,j+1,k)-w(i,j-1,k)
uk=uu(i,j,k+1)-uu(i,j,k-1)
vk=v(i,j,k+1)-v(i,j,k-1)
wk=w(i,j,k+1)-w(i,j,k-1)
ix=B%ix0(i,j,k); iy=B%iy0(i,j,k); iz=B%iz0(i,j,k)
jx=B%jx0(i,j,k); jy=B%jy0(i,j,k); jz=B%jz0(i,j,k)
kx=B%kx0(i,j,k); ky=B%ky0(i,j,k); kz=B%kz0(i,j,k)
!----对物理坐标的偏导数----------------------------------------------
ux=ui*ix+uj*jx+uk*kx
vx=vi*ix+vj*jx+vk*kx
wx=wi*ix+wj*jx+wk*kx
uy=ui*iy+uj*jy+uk*ky
vy=vi*iy+vj*jy+vk*ky
wy=wi*iy+wj*jy+wk*ky
uz=ui*iz+uj*jz+uk*kz
vz=vi*iz+vj*jz+vk*kz
wz=wi*iz+wj*jz+wk*kz
!-------------------------------------------------
omiga(i,j,k)=sqrt((wy-vz)**2+(uz-wx)**2+(vx-uy)**2)
enddo
enddo
enddo
!$OMP END PARALLEL DO
!----------------------------------------------------------------------
do ksub=1, B%subface
Bc=> B%bc_msg(ksub)
if(Bc%bc .eq. BC_WALL) then ! Wall boundary
! 沿网格线一维处理(假设网格线垂直壁面)---------------------------------------
if(Bc%face .eq. 1 .or. Bc%face .eq. 4) then ! i+ or i-
allocate(yy(nx),Amu1d(nx),Amut1d(nx),d1d(nx),u1d(nx),omiga1d(nx))
Amut1d=0.d0; Amu1d=0.d0 ! 初始化
if(Bc%face .eq. 1) then
i1=1; i2=0
else
i1=nx-1; i2=nx
endif
do k=Bc%kb,Bc%ke-1
do j=Bc%jb,Bc%je-1
x0=(B%xc(i1,j,k)+B%xc(i2,j,k))*0.5d0
y0=(B%yc(i1,j,k)+B%yc(i2,j,k))*0.5d0
z0=(B%zc(i1,j,k)+B%zc(i2,j,k))*0.5d0
do i=1,nx-1
if(Bc%face .eq. 1) then
i0=i
else
i0=nx-i
endif
yy(i0)=sqrt((B%xc(i,j,k)-x0)**2+(B%yc(i,j,k)-y0)**2+(B%zc(i,j,k)-z0)**2)
u1d(i0)=sqrt(uu(i,j,k)**2+v(i,j,k)**2+w(i,j,k)**2)
d1d(i0)=d(i,j,k)
omiga1d(i0)=omiga(i,j,k)
Amu1d(i0)=B%mu(i,j,k)
enddo
! BL模型(一维)
call BL_model_1d(nx-1,yy,Amu1d,Amut1d,d1d,u1d,omiga1d)
do i=1,nx-1
! 如果一个点处于多条线上,取粘性系数最小的值
if(Bc%face .eq. 1) then
i0=i
else
i0=nx-i
endif
if(flag1(i,j,k) .eq. 0) then
flag1(i,j,k)=1
B%mu_t(i,j,k)=Amut1d(i0)
else
B%mu_t(i,j,k)=min(B%mu_t(i,j,k),Amut1d(i0))
endif
enddo
enddo
enddo
deallocate(yy,Amu1d,Amut1d,d1d,u1d,omiga1d)
! !!! To set Amu_t=0 in the wall 设置虚网格点上的mut(-1)=-mut(1), 以保证壁面上mut=0 (mut=0.5*(mut(-1)+mut(1))
! 设置壁面第1层网格上的湍流粘性系数为0
if(Bc%face .eq. 1) then
B%mu_t(0,Bc%jb:Bc%je-1,Bc%kb:Bc%ke-1)=0.d0
B%mu_t(1,Bc%jb:Bc%je-1,Bc%kb:Bc%ke-1)=0.d0
else
B%mu_t(nx,Bc%jb:Bc%je-1,Bc%kb:Bc%ke-1)=0.d0
B%mu_t(nx-1,Bc%jb:Bc%je-1,Bc%kb:Bc%ke-1)=0.d0
endif
else if(Bc%face .eq. 2 .or. Bc%face .eq. 4) then ! face of j- or j+
allocate(yy(ny),Amu1d(ny),Amut1d(ny),d1d(ny),u1d(ny),omiga1d(ny))
Amut1d=0.d0; Amu1d=0.d0 ! 初始化
if(Bc%face .eq. 2) then
j1=1; j2=0
else
j1=ny-1; j2=ny
endif
do k=Bc%kb,Bc%ke-1
do i=Bc%ib,Bc%ie-1
x0=(B%xc(i,j1,k)+B%xc(i,j2,k))*0.5d0
y0=(B%yc(i,j1,k)+B%yc(i,j2,k))*0.5d0
z0=(B%zc(i,j1,k)+B%zc(i,j2,k))*0.5d0
do j=1,ny-1
if(Bc%face .eq. 2) then
j0=j
else
j0=ny-j
endif
yy(j0)=sqrt((B%xc(i,j,k)-x0)**2+(B%yc(i,j,k)-y0)**2+(B%zc(i,j,k)-z0)**2)
u1d(j0)=sqrt(uu(i,j,k)**2+v(i,j,k)**2+w(i,j,k)**2)
d1d(j0)=d(i,j,k)
omiga1d(j0)=omiga(i,j,k)
Amu1d(j0)=B%mu(i,j,k)
enddo
! BL模型(一维)
call BL_model_1d(ny-1,yy,Amu1d,Amut1d,d1d,u1d,omiga1d)
do j=1,ny-1
! 如果一个点处于多条线上,取粘性系数最小的值
if(Bc%face .eq. 2) then
j0=j
else
j0=ny-j
endif
if(flag1(i,j,k) .eq. 0) then
flag1(i,j,k)=1
B%mu_t(i,j,k)=Amut1d(j0)
else
B%mu_t(i,j,k)=min(B%mu_t(i,j,k),Amut1d(j0))
endif
enddo
enddo
enddo
deallocate(yy,Amu1d,Amut1d,d1d,u1d,omiga1d)
if(Bc%face .eq. 2) then
B%mu_t(Bc%ib:Bc%ie-1, 0, Bc%kb:Bc%ke-1)=0.d0
B%mu_t(Bc%ib:Bc%ie-1, 1, Bc%kb:Bc%ke-1)=0.d0
else
B%mu_t(Bc%ib:Bc%ie-1, ny, Bc%kb:Bc%ke-1)=0.d0
B%mu_t(Bc%ib:Bc%ie-1, ny-1, Bc%kb:Bc%ke-1)=0.d0
endif
else if(Bc%face .eq. 3 .or. Bc%face .eq. 6) then ! face of k- or k+
allocate(yy(nz),Amu1d(nz),Amut1d(nz),d1d(nz),u1d(nz),omiga1d(nz))
Amut1d=0.d0; Amu1d=0.d0 ! 初始化
if(Bc%face .eq. 3) then
k1=1; k2=0
else
k1=nz-1; j2=nz
endif
do j=Bc%jb,Bc%je-1
do i=Bc%ib,Bc%ie-1
x0=(B%xc(i,j,k1)+B%xc(i,j,k2))*0.5d0
y0=(B%yc(i,j,k1)+B%yc(i,j,k2))*0.5d0
z0=(B%zc(i,j,k1)+B%zc(i,j,k2))*0.5d0
do k=1,nz-1
if(Bc%face .eq. 3) then
k0=k
else
k0=nz-k
endif
yy(k0)=sqrt((B%xc(i,j,k)-x0)**2+(B%yc(i,j,k)-y0)**2+(B%zc(i,j,k)-z0)**2)
u1d(k0)=sqrt(uu(i,j,k)**2+v(i,j,k)**2+w(i,j,k)**2)
d1d(k0)=d(i,j,k)
omiga1d(k0)=omiga(i,j,k)
Amu1d(k0)=B%mu(i,j,k)
enddo
! BL模型(一维)
call BL_model_1d(nz-1,yy,Amu1d,Amut1d,d1d,u1d,omiga1d)
do k=1,nz-1
! 如果一个点处于多条线上,取粘性系数最小的值
if(Bc%face .eq. 3) then
k0=k
else
k0=nz-k
endif
if(flag1(i,j,k) .eq. 0) then
flag1(i,j,k)=1
B%mu_t(i,j,k)=Amut1d(k0)
else
B%mu_t(i,j,k)=min(B%mu_t(i,j,k),Amut1d(k0))
endif
enddo
enddo
enddo
deallocate(yy,Amu1d,Amut1d,d1d,u1d,omiga1d)
if(Bc%face .eq. 3) then
B%mu_t(Bc%ib:Bc%ie-1, Bc%jb:Bc%je-1, 0)=0.d0
B%mu_t(Bc%ib:Bc%ie-1, Bc%jb:Bc%je-1, 1)=0.d0
else
B%mu_t(Bc%ib:Bc%ie-1, Bc%jb:Bc%je-1, nz)=0.d0
B%mu_t(Bc%ib:Bc%ie-1, Bc%jb:Bc%je-1, nz-1)=0.d0
endif
endif
endif
enddo
call Amut_boundary(nMesh,mBlock)
deallocate(omiga,flag1)
end
!c------------------------------------------------------------------------
! B-L model of turbulence
! Ref: Wilox DC. Turbulence Modeling for CFD (2nd Edition), p77
subroutine BL_model_1d(ny,yy,Amu,Amu_t,d,u,omiga)
use precision_EC
implicit none
integer ny,j,Iflag
real(PRE_EC),dimension(ny) :: yy, Amu,Amu_t,d,u,omiga
real(PRE_EC),parameter:: AP=26.,Ccp=1.6,Ckleb=0.3,Cwk=0.25d0,AKT=0.4,AK=0.0168 ! Cwk= 1.d0
real(PRE_EC):: Tw,Ret,Fmax,etamax,Udif,etap,FF,Fwak,bl,Fkleb,Visti,Visto
TW=abs(Amu(1)*omiga(1))
do j=1,Ny
if(abs(Amu(j)*omiga(j)) .gt. TW) Tw=abs(Amu(j)*omiga(j))
enddo
Ret=sqrt(d(1)*TW)/Amu(1)
Fmax=0.d0 ; etamax=0.d0 ; Udif=0.d0
!----------------------------------
do j=1,Ny
if(u(j) .gt. Udif) Udif=u(j)
etap=yy(j)*Ret
FF=yy(j)*abs(omiga(j))*(1.d0-exp(-etap/AP))
if(FF.gt.Fmax) then
Fmax=FF
etamax=yy(j)
endif
! if(FF.gt.Fmax) then
! Fmax=FF
! etamax=yy(j)
! else
! goto 100 ! Find the first peak of F(y) ! 只要第1个峰值
! endif
enddo
100 continue
IFlag=0
Fwak=min(etamax*Fmax,Cwk*etamax*Udif*Udif/Fmax)
do j=1,Ny
etap=Ret*yy(j)
bl=AKT*yy(j)*(1.d0-exp(-etap/AP))
visti=d(j)*bl*bl*abs(omiga(j))
Fkleb=1.d0/(1.d0+5.5d0*(Ckleb*yy(j)/etamax)**6)
visto=AK*Ccp*d(j)*Fwak*Fkleb
if(abs(visto).lt.abs(visti)) IFlag=1
if(Iflag.eq.0) then
Amu_t(j)=visti
else
Amu_t(j)=visto
endif
enddo
end