Newer
Older
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
!pbl_dry_contrgradient.f90
!>@brief implements mass-flux contr-gradient parametrization
module pbl_dry_contrgradient
implicit none
type pblContrGradDataType
real, allocatable :: wu(:), &
theta_up(:), &
thetav_up(:), &
ua_up(:), &
va_up(:), &
qv_up(:), &
entr(:)
real :: theta_up0, thetav_up0, qv_up0
real :: ua_up0, va_up0
real :: ua_upH, va_upH
real :: theta_upH, thetav_upH, qv_upH
end type pblContrGradDataType
real, parameter :: c_entr = 0.5 !< dry entrainment parameterization
real, parameter :: b1 = 1.0; !< vertical velocity in updraft Bernulli equation constant
real, parameter :: b2 = 2.0; !< vertical velocity in updraft Bernulli equation constant
public allocate_cbl, deallocate_cbl
public get_entrainment, get_up, get_wu, apply_cntrg
contains
subroutine allocate_cbl(cbl, kmax)
integer, intent(in):: kmax
type(pblContrGradDataType), intent(out):: cbl
allocate(cbl%entr(kmax))
allocate(cbl%theta_up(kmax))
allocate(cbl%thetav_up(kmax))
allocate(cbl%qv_up(kmax))
allocate(cbl%ua_up(kmax))
allocate(cbl%va_up(kmax))
allocate(cbl%wu(kmax))
end subroutine allocate_cbl
subroutine deallocate_cbl(cbl)
deallocate(cbl%entr)
deallocate(cbl%theta_up)
deallocate(cbl%thetav_up)
deallocate(cbl%qv_up)
deallocate(cbl%ua_up)
deallocate(cbl%va_up)
deallocate(cbl%wu)
end subroutine deallocate_cbl
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
use scm_state_data, only : stateBLDataType
use pbl_turb_data, only : turbBLDataType
use phys_fluid, only: fluidParamsDataType
use pbl_grid, only : pblgridDataType
implicit none
type(pblContrGradDataType) :: cbl
type(stateBLDataType), intent(inout):: bl
type(fluidParamsDataType), intent(in) :: fluid
type(pblgridDataType), intent(in) :: grid
integer k
real dz
do k = grid%kmax, 2, -1
dz =(grid%z_cell(k-1) - grid%z_cell(k))
if (grid%z_cell(k) > bl%hpbla_diag) then
cbl%entr(k) = 0.0
else
cbl%entr(k) = c_entr * ( 1/(grid%z_cell(k) + dz) &
+ 1/(bl%hpbla_diag - grid%z_cell(k) + dz));
end if
end do
cbl%entr(1) = 0
end subroutine get_entrainment
subroutine get_up(cbl, bl, fluid, grid)
use scm_state_data, only : stateBLDataType
use pbl_turb_data, only : turbBLDataType
use phys_fluid, only: fluidParamsDataType
use pbl_grid, only : pblgridDataType
implicit none
type(pblContrGradDataType) :: cbl
type(stateBLDataType), intent(inout):: bl
type(fluidParamsDataType), intent(in) :: fluid
type(pblgridDataType), intent(in) :: grid
real phys_beta, wstr, dz
integer k, kmax
kmax = grid%kmax
phys_beta = 0.5 * fluid%g /(bl%theta_v(kmax) + bl%theta_v(kmax-1))
wstr = (phys_beta * (-bl%surf%hs/bl%rho(kmax)) &
* bl%hpbla_diag)**(0.33333)
cbl%thetav_up0 = (-bl%surf%hs/bl%rho(kmax)) / wstr
cbl%theta_up0 = (-bl%surf%hs/bl%rho(kmax)) / wstr
cbl%qv_up0 = (-bl%surf%es/bl%rho(kmax)) / wstr
cbl%ua_up0 = (-bl%surf%cm2u * bl%u(kmax)/bl%rho(kmax)) / wstr
cbl%va_up0 = (-bl%surf%cm2u * bl%v(kmax)/bl%rho(kmax)) / wstr
cbl%thetav_up(kmax) = cbl%thetav_up0
cbl%theta_up(kmax) = cbl%theta_up0
cbl%qv_up(kmax) = cbl%qv_up0
cbl%ua_up(kmax) = cbl%ua_up0
cbl%va_up(kmax) = cbl%va_up0
do k = kmax-1 ,1,-1
if (grid%z_cell(k) > bl%hpbla_diag ) then
cbl%theta_up(k) = 0.0
cbl%thetav_up(k) = 0.0
cbl%qv_up(k) = 0
else
dz = (grid%z_cell(k) - grid%z_cell(k+1))
cbl%theta_up(k) = cbl%theta_up(k + 1) &
- dz * cbl%entr(k+1) * cbl%theta_up(k+1)
cbl%thetav_up(k) = cbl%thetav_up(k + 1) &
- dz * cbl%entr(k+1) * (cbl%thetav_up(k+1))
cbl%qv_up(k) = cbl%qv_up(k + 1) &
- dz * cbl%entr(k+1) * (cbl%qv_up(k+1))
cbl%ua_up(k) = cbl%ua_up(k + 1) &
- dz * cbl%entr(k+1) * (cbl%ua_up(k+1))
cbl%ua_up(k) = cbl%va_up(k + 1) &
- dz * cbl%entr(k+1) * (cbl%va_up(k+1))
end if
end do
do k = kmax ,1,-1
cbl%thetav_up(k) = cbl%thetav_up(k) + bl%theta_v(k)
cbl%theta_up(k) = cbl%theta_up(k) + bl%theta(k)
cbl%qv_up(k) = cbl%qv_up(k) + bl%qv(k)
cbl%ua_up(k) = cbl%ua_up(k) + bl%u(k)
cbl%va_up(k) = cbl%va_up(k) + bl%v(k)
end do
end subroutine get_up
subroutine get_wu(cbl, bl, fluid, grid)
use scm_state_data, only : stateBLDataType
use pbl_turb_data, only : turbBLDataType
use phys_fluid, only: fluidParamsDataType
use pbl_grid, only : pblgridDataType
implicit none
type(pblContrGradDataType) :: cbl
type(stateBLDataType), intent(inout):: bl
type(fluidParamsDataType), intent(in) :: fluid
type(pblgridDataType), intent(in) :: grid
real wu, phys_beta, umod, ustr, rhs, a, b, dz
integer k, kmax
kmax = grid%kmax
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
umod = sqrt(bl%u(kmax)*bl%u(kmax) + bl%v(kmax)*bl%v(kmax))
ustr = sqrt(bl%surf%cm2u /bl%rho(kmax) * umod)
phys_beta = 0.5 * fluid%g /(bl%theta_v(kmax) + bl%theta_v(kmax-1))
cbl%wu(kmax) = 2.46 * ((ustr)**3 &
- 0.4 * phys_beta * bl%surf%hs &
* (grid%z_cell(k))/bl%rho(kmax))**(1.0/3.0) &
* sqrt(1.0 - grid%z_cell(k) / bl%hpbla_diag)
do k = kmax - 1, 1, -1
dz = (grid%z_cell(K) - grid%z_cell(K+1))
a = 0.0
if (abs(cbl%entr(k)) >= 1.0e-8 .and. abs(cbl%wu(k+1)) > 1.0e-6 ) then
a = b1 * cbl%entr(k) * cbl%wu(k+1) * dz
b = b2 * 9.81 * (dz/cbl%wu(k+1)) &
* (cbl%thetav_up(k) - bl%theta_v(k))/bl%theta_v(k);
rhs = cbl%wu(k+1) - a + b;
if (rhs > 0) then
cbl%wu(k) = rhs;
else
cbl%wu(k) = 0.0;
endif
else
cbl%wu(k) = 0.0;
endif
enddo
end subroutine get_wu
subroutine apply_cntrg(cbl, bl, fluid, grid, dt)
use scm_state_data, only : stateBLDataType
use pbl_turb_data, only : turbBLDataType
use phys_fluid, only: fluidParamsDataType
use pbl_grid, only : pblgridDataType
implicit none
type(stateBLDataType), intent(inout):: bl
type(fluidParamsDataType), intent(in) :: fluid
type(pblgridDataType), intent(in) :: grid
real, intent(in) :: dt
real, dimension(grid%kmax) :: rhs_up, rhs
real a0w0
integer kmax, kpbl, k
kmax = grid%kmax
kpbl = bl%kpbl
a0w0 = 0.15 * cbl%wu(kpbl)
call get_entrainment(cbl, bl, fluid, grid)
call get_up(cbl, bl, fluid, grid)
call get_wu(cbl, bl, fluid, grid)
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
!apply upwind
rhs_up(:) = 0
do k = kmax-1, 1, -1
if (abs(cbl%wu(k)) > 1.0e-8) then
rhs_up(k) = 0.1 * (cbl%wu(k)) * &
(cbl%theta_up(k+1) - cbl%theta_up(k) + bl%theta(k+1) - bl%theta(k)) &
/ (grid%z_cell(k) - grid%z_cell(k+1))
else
rhs_up(k)=0
end if
end do
bl%theta = bl%theta + dt * rhs_up
rhs_up = 0
do k = kmax-1, 1, -1
if (abs(cbl%wu(k)) > 1.0e-7) then
rhs_up(k) = 0.1 * (cbl%wu(k)) * &
(cbl%qv_up(k+1) - cbl%qv_up(k) &
+ bl%qv(k+1) - bl%qv(k)) / (grid%z_cell(K) -grid%z_cell(K+1))
else
rhs_up(k)=0
end if
end do
bl%qv(:) = bl%qv(:) + dt * rhs_up(:)
end subroutine apply_cntrg
end module pbl_dry_contrgradient