forked from ESCOMP/atmospheric_physics
-
Notifications
You must be signed in to change notification settings - Fork 0
/
kessler.F90
306 lines (264 loc) · 12.7 KB
/
kessler.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
module kessler
use ccpp_kinds, only: kind_phys
implicit none
private
save
public :: kessler_run ! Main routine
public :: kessler_init ! init routine
public :: kessler_timestep_init ! init timestep routine
! Private module data (constants set at initialization)
real(kind_phys) :: rd ! gas constant for dry air, J/(kgK)
real(kind_phys) :: cp ! heat capacity at constant pressure, J/(kgK)
real(kind_phys) :: lv ! latent heat of vaporization, J/kg
real(kind_phys) :: psl ! reference pressure at sea level, mb
real(kind_phys) :: rhoqr ! density of liquid water, kg/m^3
CONTAINS
!> \section arg_table_kessler_init Argument Table
!! \htmlinclude kessler_init.html
subroutine kessler_init(rd_in, cp_in, lv_in, psl_in, rhoqr_in, errmsg, errflg)
! Set physical constants to be consistent with calling model
real(kind_phys), intent(in) :: rd_in ! gas constant for dry air, J/(kgK)
real(kind_phys), intent(in) :: cp_in ! heat capacity at constant pres., J/(kgK)
real(kind_phys), intent(in) :: lv_in ! latent heat of vaporization, J/kg
real(kind_phys), intent(in) :: psl_in ! reference pressure at sea level, mb
real(kind_phys), intent(in) :: rhoqr_in ! density of liquid water, kg/m^3
character(len=512), intent(out) :: errmsg
integer, intent(out) :: errflg
errmsg = ''
errflg = 0
rd = rd_in
cp = cp_in
lv = lv_in
psl = psl_in/100._kind_phys
rhoqr = rhoqr_in
end subroutine kessler_init
!> \section arg_table_kessler_timestep_init Argument Table
!! \htmlinclude kessler_timestep_init.html
subroutine kessler_timestep_init(ncol, nz, pdel, pdeldry, qv, qc, qr, errmsg, errflg)
use state_converters, only : wet_to_dry_run
! Dummy arguments
integer, intent(in) :: ncol
integer, intent(in) :: nz
real(kind_phys), intent(in) :: pdel(:,:)
real(kind_phys), intent(in) :: pdeldry(:,:)
real(kind_phys), intent(inout) :: qv(:,:)
real(kind_phys), intent(inout) :: qc(:,:)
real(kind_phys), intent(inout) :: qr(:,:)
character(len=*), intent(out) :: errmsg
integer, intent(out) :: errflg
! Local variables
integer :: k
errflg = 0
errmsg = ''
!!XXgoldyXX ==> cacraigucar: Why is this here instead of having
! a wet_to_dry_timestep_init?
call wet_to_dry_run(ncol, nz, pdel, pdeldry, qv, qc, qr, errmsg, errflg)
end subroutine kessler_timestep_init
!-----------------------------------------------------------------------
!
! Version: 2.0
!
! Date: January 22nd, 2015
!
! Change log:
! v2 - Added sub-cycling of rain sedimentation so as not to violate
! CFL condition.
!
! The KESSLER subroutine implements the Kessler (1969) microphysics
! parameterization as described by Soong and Ogura (1973) and Klemp
! and Wilhelmson (1978, KW). KESSLER is called at the end of each
! time step and makes the final adjustments to the potential
! temperature and moisture variables due to microphysical processes
! occurring during that time step. KESSLER is called once for each
! vertical column of grid cells. Increments are computed and added
! into the respective variables. The Kessler scheme contains three
! moisture categories: water vapor, cloud water (liquid water that
! moves with the flow), and rain water (liquid water that falls
! relative to the surrounding air). There are no ice categories.
! Variables in the column are ordered from the surface to the top.
!
! SUBROUTINE KESSLER(theta, qv, qc, qr, rho, pk, dt, z, nz, rainnc)
!
! Input variables:
! temp - temperature (K)
! qv - water vapor mixing ratio (gm/gm)
! qc - cloud water mixing ratio (gm/gm)
! qr - rain water mixing ratio (gm/gm)
! rho - dry air density (not mean state as in KW) (kg/m^3)
! pk - Exner function (not mean state as in KW) (p/p0)**(R/cp)
! dt - time step (s)
! z - heights of thermodynamic levels in the grid column (m)
! nz - number of thermodynamic levels in the column
! precl - Precipitation rate (m_water/s)
!
! Output variables:
! Increments are added into t, qv, qc, qr, and rainnc which are
! returned to the routine from which KESSLER was called. To obtain
! the total precip qt, after calling the KESSLER routine, compute:
!
! qt = sum over surface grid cells of (rainnc * cell area) (kg)
! [here, the conversion to kg uses (10^3 kg/m^3)*(10^-3 m/mm) = 1]
!
!
! Authors: Paul Ullrich
! University of California, Davis
! Email: [email protected]
!
! Based on a code by Joseph Klemp
! (National Center for Atmospheric Research)
!
! Reference:
!
! Klemp, J. B., W. C. Skamarock, W. C., and S.-H. Park, 2015:
! Idealized Global Nonhydrostatic Atmospheric Test Cases on a Reduced
! Radius Sphere. Journal of Advances in Modeling Earth Systems.
! doi:10.1002/2015MS000435
!
!=======================================================================
!> \section arg_table_kessler_run Argument Table
!! \htmlinclude kessler_run.html
subroutine kessler_run(ncol, nz, dt, rho, z, pk, theta, &
qv, qc, qr, precl, errmsg, errflg)
!------------------------------------------------
! Input / output parameters
!------------------------------------------------
integer, intent(in) :: ncol ! Number of columns
integer, intent(in) :: nz ! Number of vertical levels
real(kind_phys), intent(in) :: dt ! Time step (s)
real(kind_phys), intent(in) :: rho(:,:) ! Dry air density (kg/m^3)
real(kind_phys), intent(in) :: z(:,:) ! Heights of thermo. levels (m)
real(kind_phys), intent(in) :: pk(:,:) ! Exner function (p/p0)**(R/cp)
real(kind_phys), intent(inout) :: theta(:,:) ! temperature (K)
real(kind_phys), intent(inout) :: qv(:,:) ! Water vapor mixing ratio (gm/gm)
real(kind_phys), intent(inout) :: qc(:,:) ! Cloud water mixing ratio (gm/gm)
real(kind_phys), intent(inout) :: qr(:,:) ! Rain water mixing ratio (gm/gm)
real(kind_phys), intent(out) :: precl(:) ! Precipi tation rate (m_water / s)
character(len=*), intent(out) :: errmsg
integer, intent(out) :: errflg
!------------------------------------------------
! Local variables
!------------------------------------------------
real(kind_phys) :: r(nz), rhalf(nz), velqr(nz), sed(nz), pc(nz)
real(kind_phys) :: f5, f2x, xk, ern, qrprod, prod, qvs, dt_max, dt0
integer :: col, klev, rainsplit, nt
! Initialize output variables
precl = 0._kind_phys
errmsg = ''
errflg = 0
! Check inputs
if (dt <= 0._kind_phys) then
write(errmsg,*) 'KESSLER called with nonpositive dt'
errflg = 1
return
end if
!------------------------------------------------
! Begin calculation
!------------------------------------------------
f2x = 17.27_kind_phys
f5 = 237.3_kind_phys * f2x * lv / cp
xk = .2875_kind_phys ! kappa (r/cp)
! Loop through columns
do col = 1, ncol
do klev = 1, nz
r(klev) = 0.001_kind_phys * rho(col, klev)
rhalf(klev) = sqrt(rho(col, 1) / rho(col, klev))
pc(klev) = 3.8_kind_phys / (pk(col, klev)**(1._kind_phys/xk)*psl)
!
! if qr is (round-off) negative then the computation of
! velqr triggers floating point exception error when running
! in debugging mode with NAG
!
qr(col,klev) = MAX(qr(col,klev),0.0_kind_phys)
!
! Liquid water terminal velocity (m/s) following KW eq. 2.15
velqr(klev) = 36.34_kind_phys * rhalf(klev) * &
(qr(col, klev) * r(klev))**0.1364_kind_phys
end do
! Maximum time step size in accordance with CFL condition
dt_max = dt
do klev = 1, nz - 1
! NB: Original test for velqr /= 0 numerically unstable
if (abs(velqr(klev)) > 1.0E-12_kind_phys) then
dt_max = min(dt_max, 0.8_kind_phys*(z(col, klev+1) - &
z(col, klev)) / velqr(klev))
end if
end do
! Number of subcycles
rainsplit = ceiling(dt / dt_max)
if (rainsplit < 1) then
write(errmsg, *) 'KESSLER: bad rainsplit ',dt,dt_max,rainsplit
errflg = 1
return
end if
dt0 = dt / real(rainsplit, kind_phys)
! Subcycle through rain process
nt = 1
do while (nt <= rainsplit)
! Precipitation rate (m/s)
precl(col) = precl(col) + rho(col, 1) * qr(col, 1) * &
velqr(1) / rhoqr
! Sedimentation term using upstream differencing
do klev = 1, nz-1
sed(klev) = dt0 * &
((r(klev+1) * qr(col, klev+1) * velqr(klev+1)) - &
(r(klev) * qr(col, klev) * velqr(klev))) / &
(r(klev) * (z(col, klev+1) - z(col, klev)))
end do
sed(nz) = -dt0 * qr(col, nz) * velqr(nz) / &
(0.5_kind_phys * (z(col, nz)-z(col, nz-1)))
! Adjustment terms
do klev = 1, nz
! Autoconversion and accretion rates following KW eq. 2.13a,b
qrprod = qc(col, klev) - (qc(col, klev) - dt0 * &
max(.001_kind_phys * (qc(col, klev)-.001_kind_phys), &
0._kind_phys)) / &
(1._kind_phys + dt0 * 2.2_kind_phys * &
qr(col, klev)**.875_kind_phys)
qc(col, klev) = max(qc(col, klev) - qrprod, 0._kind_phys)
qr(col, klev) = max(qr(col, klev) + qrprod + sed(klev), &
0._kind_phys)
! Saturation vapor mixing ratio (gm/gm) following KW eq. 2.11
qvs = pc(klev) * exp(f2x*(pk(col, klev)*theta(col, klev) - 273._kind_phys) / (pk(col, klev)*theta(col, klev) &
- 36._kind_phys))
prod = (qv(col, klev) - qvs) / (1._kind_phys + qvs*f5 / (pk(col, klev)*theta(col, klev) - 36._kind_phys)**2)
! Evaporation rate following KW eq. 2.14a,b
ern = min(dt0 * (((1.6_kind_phys + 124.9_kind_phys*(r(klev)*qr(col, klev))**.2046_kind_phys) * &
(r(klev) * qr(col, klev))**.525_kind_phys) / &
(2550000._kind_phys * pc(klev) / (3.8_kind_phys*qvs) + 540000._kind_phys)) * &
(dim(qvs,qv(col, klev)) / (r(klev)*qvs)), &
max(-prod-qc(col, klev),0._kind_phys),qr(col, klev))
! Saturation adjustment following KW eq. 3.10
theta(col, klev)= theta(col, klev) + (lv / (cp * pk(col, klev)) * (max(prod,-qc(col, klev)) - ern))
qv(col, klev) = max(qv(col, klev) - max(prod, -qc(col, klev)) + ern, 0._kind_phys)
qc(col, klev) = qc(col, klev) + max(prod, -qc(col, klev))
qr(col, klev) = qr(col, klev) - ern
end do
! Recalculate liquid water terminal velocity
if (nt /= rainsplit) then
do klev = 1, nz
velqr(klev) = 36.34_kind_phys * rhalf(klev) * (qr(col, klev)*r(klev))**0.1364_kind_phys
end do
!
! recompute rainsplit since velqr has changed
!
do klev = 1, nz - 1
if (abs(velqr(klev)) > 1.0E-12_kind_phys) then
dt_max = min(dt_max, 0.8_kind_phys*(z(col, klev+1) - z(col, klev)) / velqr(klev))
end if
end do
! Number of subcycles
rainsplit = ceiling(dt / dt_max)
if (rainsplit < 1) then
write(errmsg, *) 'KESSLER: bad rainsplit ',dt,dt_max,rainsplit
errflg = 1
return
end if
dt0 = dt / real(rainsplit, kind_phys)
end if
nt=nt+1
end do
precl(col) = precl(col) / real(rainsplit, kind_phys)
end do ! column loop
end subroutine kessler_run
!=======================================================================
end module kessler