forked from ImperialCollegeLondon/WInc3D
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathABL.f90
356 lines (294 loc) · 12.1 KB
/
ABL.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
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
!==========================================================
! Subroutines for computing the SFS from the ABL
!==========================================================
subroutine wall_sgs(ux,uy,uz,nut1,sxy1,syz1,tauwallxy,tauwallzy,wallfluxx,wallfluxy,wallfluxz)
USE MPI
USE decomp_2d
USE param
USE var, only: uxf1,uzf1, di1
implicit none
real(mytype),dimension(xsize(1),xsize(2),xsize(3)) :: ux,uy,uz, nut1
real(mytype),dimension(xsize(1),xsize(2),xsize(3)) :: wallfluxx,wallfluxy,wallfluxz
real(mytype),dimension(xsize(1),xsize(2),xsize(3)) :: wallfluxxf,wallfluxyf,wallfluxzf
real(mytype),dimension(xsize(1),xsize(2),xsize(3)) :: sxy1, syz1
real(mytype),dimension(xsize(1),xsize(3)) :: tauwallxy, tauwallzy
integer :: i,j,k,code
real(mytype) :: ut,ut1,utt,ut11, abl_vel, ABLtaux, ABLtauz, delta
real(mytype) :: ux_HAve_local, uz_HAve_local,S_HAve_local
real(mytype) :: ux_HAve, uz_HAve,S_HAve,ux12,uz12,S12
real(mytype) :: sxy_HAve_local, szy_HAve_local, nut_HAve_local, nutsxy_HAve_local, nutszy_HAve_local
real(mytype) :: sxy_HAve, szy_HAve, nut_HAve, nutsxy_HAve, nutszy_HAve
real(mytype) :: nutprimes
real(mytype) :: nuLESBar, scriptR, xi1, nuLES, TS1, TR1
real(mytype) :: CD ! drag coefficient
call filter(0.0d0)
! Filter the velocity with twice the grid scale according to Bou-zeid et al 2005
call filx(uxf1,ux,di1,fisx,fiffx,fifsx,fifwx,xsize(1),xsize(2),xsize(3),0)
call filx(uzf1,uz,di1,fisx,fiffx,fifsx,fifwx,xsize(1),xsize(2),xsize(3),0)
! Determine the shear stress using Moeng's formulation
!*****************************************************************************************
ux_HAve_local=0.
uz_HAve_local=0.
S_HAve_local=0.
if (xstart(2)==1) then
!if (nrank==0) print *, 'Max of the filtered velocity', maxval(ux), 'Max of the unfiltered velocity', maxval(ux)
do k=1,xsize(3)
do i=1,xsize(1)
ux_HAve_local=ux_HAve_local+0.5*(uxf1(i,1,k)+uxf1(i,2,k))
uz_HAve_local=uz_HAve_local+0.5*(uzf1(i,1,k)+uzf1(i,2,k))
S_HAve_local=S_HAve_local+sqrt((0.5*(uxf1(i,1,k)+uxf1(i,2,k)))**2.+(0.5*(uzf1(i,1,k)+uzf1(i,2,k)))**2.)
enddo
enddo
ux_HAve_local=ux_HAve_local/xsize(3)/xsize(1)
uz_HAve_local=uz_HAve_local/xsize(3)/xsize(1)
S_HAve_local=S_HAve_local/xsize(3)/xsize(1)
else
ux_HAve_local=0.
uz_HAve_local=0.
S_HAve_local =0.
endif
call MPI_ALLREDUCE(ux_HAve_local,ux_HAve,1,real_type,MPI_SUM,MPI_COMM_WORLD,code)
call MPI_ALLREDUCE(uz_HAve_local,uz_HAve,1,real_type,MPI_SUM,MPI_COMM_WORLD,code)
call MPI_ALLREDUCE(S_HAve_local,S_HAve,1,real_type,MPI_SUM,MPI_COMM_WORLD,code)
ux_HAve=ux_HAve/p_col
uz_HAve=uz_HAve/p_col
S_HAve= S_HAve/p_col
if (istret.ne.0) delta=(yp(2)-yp(1))/2.0
if (istret.eq.0) delta=dy/2.
! Compute the friction velocity u_shear and the drag coefficient
u_shear=k_roughness*sqrt(ux_HAve**2.+uz_HAve**2.)/log(delta/z_zero)
CD=k_roughness**2./log(delta/z_zero)**2.
if (nrank==0) then
print *, ' '
print *, ' ABL:'
print *, ' Horizontally-averaged velocity at y=1/2... ', ux_HAve,uz_Have
print *, ' Friction velocity : ', u_shear
print *, ' Drag Coefficient : ', CD
endif
!Compute the shear stresses -- only on the wall
!u_shear=ustar
wallfluxx = 0.
wallfluxy = 0.
wallfluxz = 0.
! Apply BCs locally
if (xstart(2)==1) then
do k=1,xsize(3)
do i=1,xsize(1)
ux12=0.5*(uxf1(i,1,k)+uxf1(i,2,k))
uz12=0.5*(uzf1(i,1,k)+uzf1(i,2,k))
S12=sqrt(ux12**2.+uz12**2.)
if(iwallmodel==1) then ! MOENG
tauwallxy(i,k)=-(k_roughness/log(delta/z_zero))**2.0*ux_HAve*S_HAve
tauwallzy(i,k)=-(k_roughness/log(delta/z_zero))**2.0*uz_HAve*S_HAve
else !Parlage model
tauwallxy(i,k)=-(k_roughness/log(delta/z_zero))**2.0*ux12*S12
tauwallzy(i,k)=-(k_roughness/log(delta/z_zero))**2.0*uz12*S12
endif
if(jLES.ge.2) then ! Apply second order central finite difference schemes for the near wall
wallfluxx(i,1,k) = -(-1./2.*(-2.*nut1(i,3,k)*sxy1(i,3,k))+&
2.*(-2.*nut1(i,2,k)*sxy1(i,2,k))-3./2.*tauwallxy(i,k))/(2.*delta)
wallfluxy(i,1,k) = 0.
wallfluxz(i,1,k) = -(-1./2.*(-2.*nut1(i,3,k)*syz1(i,3,k))+&
2.*(-2.*nut1(i,2,k)*syz1(i,2,k))-3./2.*tauwallzy(i,k))/(2.*delta)
else
wallfluxx(i,1,k) = -CD*S12*delta
wallfluxy(i,1,k) = 0.!
wallfluxz(i,1,k) = -CD*S12*delta
endif
enddo
enddo
endif
!*********************************************************************************************************
! Filter wallfluxes
!call filx(wallfluxxf,wallfluxx,di1,sx,vx,fiffx,fifx,ficx,fibx,fibbx,filax,&
!fiz1x,fiz2x,xsize(1),xsize(2),xsize(3),0)
!call filx(wallfluxyf,wallfluxy,di1,sx,vx,fiffx,fifx,ficx,fibx,fibbx,filax,&
!fiz1x,fiz2x,xsize(1),xsize(2),xsize(3),0)
!call filx(wallfluxzf,wallfluxz,di1,sx,vx,fiffx,fifx,ficx,fibx,fibbx,filax,&
!fiz1x,fiz2x,xsize(1),xsize(2),xsize(3),0)
!wallfluxx(:,:,:) = wallfluxxf(:,:,:)
!wallfluxy(:,:,:) = wallfluxyf(:,:,:)
!wallfluxz(:,:,:) = wallfluxzf(:,:,:)
if (nrank==0) write(*,*) 'Maximum wall shear stress for x and z', maxval(tauwallxy), maxval(tauwallzy)
if (nrank==0) write(*,*) 'Minimum wall shear stress for x and z', minval(tauwallxy), minval(tauwallzy)
if (nrank==0) write(*,*) 'Wall flux x ', maxval(wallfluxx), maxval(wallfluxz)
if (nrank==0) write(*,*) 'Wall flux z ', minval(wallfluxx), minval(wallfluxz)
return
end subroutine wall_sgs
subroutine wall_shear_stress(ux,uy,uz,di1,tauwallxy,tauwallzy)
USE MPI
USE decomp_2d
USE param
implicit none
real(mytype),dimension(xsize(1),xsize(2),xsize(3)) :: ux,uy,uz,nut1
real(mytype),dimension(xsize(1),xsize(2),xsize(3)) :: uxf,uyf,uzf
real(mytype),dimension(xsize(1),xsize(3)) :: tauwallxy, tauwallzy
real(mytype),dimension(xsize(1),xsize(2),xsize(3)) :: di1
integer :: i,j,k,code
real(mytype) :: delta
real(mytype) :: ux_HAve_local, uz_HAve_local,S_HAve_local
real(mytype) :: ux_HAve, uz_HAve,S_HAve,ux12,uz12,S12
real(mytype) :: sxy_HAve_local, szy_HAve_local, nut_HAve_local, nutsxy_HAve_local, nutszy_HAve_local
real(mytype) :: sxy_HAve, szy_HAve, nut_HAve, nutsxy_HAve, nutszy_HAve
real(mytype) :: nutprimes
real(mytype) :: nuLESBar, scriptR, xi1, nuLES, TS1, TR1
real(mytype) :: CD ! drag coefficient
call filter(0.)
call filx(uxf,ux,di1,fisx,fiffx,fifsx,fifwx,xsize(1),xsize(2),xsize(3),0)
call filx(uzf,uz,di1,fisx,fiffx,fifsx,fifwx,xsize(1),xsize(2),xsize(3),0)
! Filter the velocity with twice the grid scale according to Bou-zeid et al 2005
! Determine the shear stress using Moeng's formulation
!*****************************************************************************************
ux_HAve_local=0.
uz_HAve_local=0.
S_HAve_local=0.
if (xstart(2)==1) then
if (nrank==0) print *, 'Max of the filtered velocity', maxval(uxf), 'Max of the unfiltered velocity', maxval(ux)
do k=1,xsize(3)
do i=1,xsize(1)
ux_HAve_local=ux_HAve_local+0.5*(uxf(i,1,k)+uxf(i,2,k))
uz_HAve_local=uz_HAve_local+0.5*(uzf(i,1,k)+uzf(i,2,k))
S_HAve_local=S_HAve_local+sqrt((0.5*(uxf(i,1,k)+uxf(i,2,k)))**2.+(0.5*(uzf(i,1,k)+uzf(i,2,k)))**2.)
enddo
enddo
ux_HAve_local=ux_HAve_local/xsize(3)/xsize(1)
uz_HAve_local=uz_HAve_local/xsize(3)/xsize(1)
S_HAve_local=S_HAve_local/xsize(3)/xsize(1)
else
ux_HAve_local=0.
uz_HAve_local=0.
S_HAve_local =0.
endif
call MPI_ALLREDUCE(ux_HAve_local,ux_HAve,1,real_type,MPI_SUM,MPI_COMM_WORLD,code)
call MPI_ALLREDUCE(uz_HAve_local,uz_HAve,1,real_type,MPI_SUM,MPI_COMM_WORLD,code)
call MPI_ALLREDUCE(S_HAve_local,S_HAve,1,real_type,MPI_SUM,MPI_COMM_WORLD,code)
ux_HAve=ux_HAve/p_col
uz_HAve=uz_HAve/p_col
S_HAve= S_HAve/p_col
if (istret.ne.0) delta=(yp(2)-yp(1))/2.0
if (istret.eq.0) delta=dy/2.
! Compute the friction velocity u_shear and the drag coefficient
u_shear=k_roughness*sqrt(ux_HAve**2.+uz_HAve**2.)/log(delta/z_zero)
CD=k_roughness**2./log(delta/z_zero)**2.
if (nrank==0) then
print *, ' '
print *, ' ABL:'
print *, ' Horizontally-averaged velocity at y=1/2... ', ux_HAve,uz_Have
print *, ' Friction velocity : ', u_shear
print *, ' Drag Coefficient : ', CD
endif
!Compute the shear stresses -- only on the wall
if (xstart(2)==1) then
do k=1,xsize(3)
do i=1,xsize(1)
ux12=0.5*(uxf(i,1,k)+uxf(i,2,k))
uz12=0.5*(uzf(i,1,k)+uzf(i,2,k))
S12=sqrt(ux12**2.+uz12**2.)
if(iwallmodel==1) then ! MOENG
tauwallxy(i,k)=-u_shear**2.0*(S12*ux_HAve+S_HAve*(ux12-ux_HAve))/(S_HAve*sqrt(ux_HAve**2.+uz_HAve**2.))
tauwallzy(i,k)=-u_shear**2.0*(S12*uz_HAve+S_HAve*(uz12-uz_HAve))/(S_HAve*sqrt(ux_HAve**2.+uz_Have**2.))
else !Parlage model
tauwallxy(i,k)=-(k_roughness/log(delta/z_zero))**2.0*ux12*S12
tauwallzy(i,k)=-(k_roughness/log(delta/z_zero))**2.0*uz12*S12
endif
enddo
enddo
endif
!*********************************************************************************************************
if (nrank==0) write(*,*) 'Maximum wall shear stress for x and z', maxval(tauwallxy), maxval(tauwallzy)
if (nrank==0) write(*,*) 'Minimum wall shear stress for x and z', minval(tauwallxy), minval(tauwallzy)
return
end subroutine wall_shear_stress
subroutine tripping(tb,ta) !TRIPPING SUBROUTINE FOR ATMOSPHERIC BOUNDARY LAYERS
USE param
USE var
USE decomp_2d
USE MPI
implicit none
integer :: i,j,k
real(mytype),dimension(xsize(1),xsize(2),xsize(3)) :: tb, ta
integer :: seed0, ii, code
real(mytype) :: z_pos, randx, p_tr, b_tr, x_pos, y_pos, A_tr
!Done in X-Pencils
seed0=randomseed !Seed for random number
!A_tr=A_trip*min(1.0,0.8+real(itime)/200.0)
!xs_tr=4.0/2.853
!ys_tr=2.0/2.853
!ts_tr=4.0/2.853
!x0_tr=40.0/2.853
A_tr = 0.1*dt
if ((itime.eq.ifirst).and.(nrank.eq.0)) then
call random_seed(SIZE=ii)
call random_seed(PUT=seed0*(/ (1, i = 1, ii) /))
!DEBUG:
!call random_number(randx)
!call MPI_BCAST(randx,1,real_type,0,MPI_COMM_WORLD,code)
!write(*,*) 'RANDOM:', nrank, randx, ii
!First random generation of h_nxt
do j=1,z_modes
call random_number(randx)
h_coeff(j)=1.0*(randx-0.5)
enddo
h_coeff=h_coeff/sqrt(DBLE(z_modes))
endif
!Initialization h_nxt (always bounded by xsize(3)^2 operations)
if (itime.eq.ifirst) then
call MPI_BCAST(h_coeff,z_modes,real_type,0,MPI_COMM_WORLD,code)
nxt_itr=0
do k=1,xsize(3)
h_nxt(k)=0.0
z_pos=-zlz/2.0+(xstart(3)+(k-1)-1)*dz
do j=1,z_modes
h_nxt(k)= h_nxt(k)+h_coeff(j)*sin(2.0*pi*j*z_pos/zlz)
enddo
enddo
end if
!Time-loop
i=int(t/ts_tr)
if (i.ge.nxt_itr) then !Nxt_itr is a global variable
nxt_itr=i+1
!First random generation of h
h_i(:)=h_nxt(:)
if (nrank .eq. 0) then
do j=1,z_modes
call random_number(randx)
h_coeff(j)=1.0*(randx-0.5)
enddo
h_coeff=h_coeff/sqrt(DBLE(z_modes)) !Non-dimensionalization
end if
call MPI_BCAST(h_coeff,z_modes,real_type,0,MPI_COMM_WORLD,code)
!Initialization h_nxt (always bounded by z_steps^2 operations)
do k=1,xsize(3)
h_nxt(k)=0.0
z_pos=-zlz/2.0+(xstart(3)+(k-1)-1)*dz
do j=1,z_modes
h_nxt(k)= h_nxt(k)+h_coeff(j)*sin(2.0*pi*j*z_pos/zlz)
enddo
enddo
endif
!Time coefficient
p_tr=t/ts_tr-i
b_tr=3.0*p_tr**2-2.0*p_tr**3
!Creation of tripping velocity
do i=1,xsize(1)
x_pos=(xstart(1)+(i-1)-1)*dx
do j=1,xsize(2)
!y_pos=(xstart(2)+(j-1)-1)*dy
y_pos=yp(xstart(2)+(j-1))
do k=1,xsize(3)
!g(z)*EXP_F(X,Y)
ta(i,j,k)=((1.0-b_tr)*h_i(k)+b_tr*h_nxt(k))
!ta(i,j,k)=A_tr*exp(-((x_pos-x0_tr)/xs_tr)**2-(y_pos/ys_tr)**2)*ta(i,j,k)
ta(i,j,k)=A_tr*exp(-((x_pos-x0_tr)/xs_tr)**2-((y_pos-0.5)/ys_tr)**2)*ta(i,j,k)
tb(i,j,k)=tb(i,j,k)+ta(i,j,k)
z_pos=-zlz/2.0+(xstart(3)+(k-1)-1)*dz
! if ((((x_pos-x0_tr)**2).le.9.0e-3).and.(y_pos.le.0.0001).and.((z_pos).le.0.03))then
! open(442,file='tripping.dat',form='formatted',position='APPEND')
! write(442,*) t,ta(i,j,k)
! close(442)
! end if
enddo
enddo
enddo
return
end subroutine tripping