-
Notifications
You must be signed in to change notification settings - Fork 11
/
Copy pathEDMF_Updrafts.pyx
425 lines (367 loc) · 19.2 KB
/
EDMF_Updrafts.pyx
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
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
#!python
#cython: boundscheck=False
#cython: wraparound=False
#cython: initializedcheck=True
#cython: cdivision=False
import numpy as np
include "parameters.pxi"
from thermodynamic_functions cimport *
from microphysics_functions cimport *
import cython
cimport Grid
cimport ReferenceState
from Variables cimport GridMeanVariables
from NetCDFIO cimport NetCDFIO_Stats
from EDMF_Environment cimport EnvironmentVariables
from libc.math cimport fmax, fmin
import pylab as plt
cdef class UpdraftVariable:
def __init__(self, nu, nz, loc, kind, name, units):
self.values = np.zeros((nu,nz),dtype=np.double, order='c')
self.old = np.zeros((nu,nz),dtype=np.double, order='c') # needed for prognostic updrafts
self.new = np.zeros((nu,nz),dtype=np.double, order='c') # needed for prognostic updrafts
self.tendencies = np.zeros((nu,nz),dtype=np.double, order='c')
self.flux = np.zeros((nu,nz),dtype=np.double, order='c')
self.bulkvalues = np.zeros((nz,), dtype=np.double, order = 'c')
if loc != 'half' and loc != 'full':
print('Invalid location setting for variable! Must be half or full')
self.loc = loc
if kind != 'scalar' and kind != 'velocity':
print ('Invalid kind setting for variable! Must be scalar or velocity')
self.kind = kind
self.name = name
self.units = units
cpdef set_bcs(self,Grid.Grid Gr):
cdef:
Py_ssize_t i,k
Py_ssize_t start_low = Gr.gw - 1
Py_ssize_t start_high = Gr.nzg - Gr.gw - 1
n_updrafts = np.shape(self.values)[0]
if self.name == 'w':
for i in xrange(n_updrafts):
self.values[i,start_high] = 0.0
self.values[i,start_low] = 0.0
for k in xrange(1,Gr.gw):
self.values[i,start_high+ k] = -self.values[i,start_high - k ]
self.values[i,start_low- k] = -self.values[i,start_low + k ]
else:
for k in xrange(Gr.gw):
for i in xrange(n_updrafts):
self.values[i,start_high + k +1] = self.values[i,start_high - k]
self.values[i,start_low - k] = self.values[i,start_low + 1 + k]
return
cdef class UpdraftVariables:
def __init__(self, nu, namelist, paramlist, Grid.Grid Gr):
self.Gr = Gr
self.n_updrafts = nu
cdef:
Py_ssize_t nzg = Gr.nzg
Py_ssize_t i, k
self.W = UpdraftVariable(nu, nzg, 'full', 'velocity', 'w','m/s' )
self.Area = UpdraftVariable(nu, nzg, 'full', 'scalar', 'area_fraction','[-]' )
self.QT = UpdraftVariable(nu, nzg, 'half', 'scalar', 'qt','kg/kg' )
self.QL = UpdraftVariable(nu, nzg, 'half', 'scalar', 'ql','kg/kg' )
self.QR = UpdraftVariable(nu, nzg, 'half', 'scalar', 'qr','kg/kg' )
if namelist['thermodynamics']['thermal_variable'] == 'entropy':
self.H = UpdraftVariable(nu, nzg, 'half', 'scalar', 's','J/kg/K' )
elif namelist['thermodynamics']['thermal_variable'] == 'thetal':
self.H = UpdraftVariable(nu, nzg, 'half', 'scalar', 'thetal','K' )
self.THL = UpdraftVariable(nu, nzg, 'half', 'scalar', 'thetal', 'K')
self.T = UpdraftVariable(nu, nzg, 'half', 'scalar', 'temperature','K' )
self.B = UpdraftVariable(nu, nzg, 'half', 'scalar', 'buoyancy','m^2/s^3' )
if namelist['turbulence']['scheme'] == 'EDMF_PrognosticTKE':
try:
use_steady_updrafts = namelist['turbulence']['EDMF_PrognosticTKE']['use_steady_updrafts']
except:
use_steady_updrafts = False
if use_steady_updrafts:
self.prognostic = False
else:
self.prognostic = True
self.updraft_fraction = paramlist['turbulence']['EDMF_PrognosticTKE']['surface_area']
else:
self.prognostic = False
self.updraft_fraction = paramlist['turbulence']['EDMF_BulkSteady']['surface_area']
self.cloud_base = np.zeros((nu,), dtype=np.double, order='c')
self.cloud_top = np.zeros((nu,), dtype=np.double, order='c')
self.cloud_cover = np.zeros((nu,), dtype=np.double, order='c')
return
cpdef initialize(self, GridMeanVariables GMV):
cdef:
Py_ssize_t i,k
Py_ssize_t gw = self.Gr.gw
double dz = self.Gr.dz
with nogil:
for i in xrange(self.n_updrafts):
for k in xrange(self.Gr.nzg):
self.W.values[i,k] = 0.0
# Simple treatment for now, revise when multiple updraft closures
# become more well defined
if self.prognostic:
self.Area.values[i,k] = 0.0 #self.updraft_fraction/self.n_updrafts
else:
self.Area.values[i,k] = self.updraft_fraction/self.n_updrafts
self.QT.values[i,k] = GMV.QT.values[k]
self.QL.values[i,k] = GMV.QL.values[k]
self.QR.values[i,k] = GMV.QR.values[k]
self.H.values[i,k] = GMV.H.values[k]
self.T.values[i,k] = GMV.T.values[k]
self.B.values[i,k] = 0.0
self.Area.values[i,gw] = self.updraft_fraction/self.n_updrafts
self.QT.set_bcs(self.Gr)
self.QR.set_bcs(self.Gr)
self.H.set_bcs(self.Gr)
return
cpdef initialize_io(self, NetCDFIO_Stats Stats):
Stats.add_profile('updraft_area')
Stats.add_profile('updraft_w')
Stats.add_profile('updraft_qt')
Stats.add_profile('updraft_ql')
Stats.add_profile('updraft_qr')
if self.H.name == 'thetal':
Stats.add_profile('updraft_thetal')
else:
# Stats.add_profile('updraft_thetal')
Stats.add_profile('updraft_s')
Stats.add_profile('updraft_temperature')
Stats.add_profile('updraft_buoyancy')
Stats.add_ts('updraft_cloud_cover')
Stats.add_ts('updraft_cloud_base')
Stats.add_ts('updraft_cloud_top')
return
cpdef set_means(self, GridMeanVariables GMV):
cdef:
Py_ssize_t i, k
self.Area.bulkvalues = np.sum(self.Area.values,axis=0)
self.W.bulkvalues[:] = 0.0
self.QT.bulkvalues[:] = 0.0
self.QL.bulkvalues[:] = 0.0
self.QR.bulkvalues[:] = 0.0
self.H.bulkvalues[:] = 0.0
self.T.bulkvalues[:] = 0.0
self.B.bulkvalues[:] = 0.0
with nogil:
for k in xrange(self.Gr.gw, self.Gr.nzg-self.Gr.gw):
if self.Area.bulkvalues[k] > 1.0e-20:
for i in xrange(self.n_updrafts):
self.QT.bulkvalues[k] += self.Area.values[i,k] * self.QT.values[i,k]/self.Area.bulkvalues[k]
self.QL.bulkvalues[k] += self.Area.values[i,k] * self.QL.values[i,k]/self.Area.bulkvalues[k]
self.QR.bulkvalues[k] += self.Area.values[i,k] * self.QR.values[i,k]/self.Area.bulkvalues[k]
self.H.bulkvalues[k] += self.Area.values[i,k] * self.H.values[i,k]/self.Area.bulkvalues[k]
self.T.bulkvalues[k] += self.Area.values[i,k] * self.T.values[i,k]/self.Area.bulkvalues[k]
self.B.bulkvalues[k] += self.Area.values[i,k] * self.B.values[i,k]/self.Area.bulkvalues[k]
self.W.bulkvalues[k] += ((self.Area.values[i,k] + self.Area.values[i,k+1]) * self.W.values[i,k]
/(self.Area.bulkvalues[k] + self.Area.bulkvalues[k+1]))
else:
self.QT.bulkvalues[k] = GMV.QT.values[k]
self.QR.bulkvalues[k] = GMV.QR.values[k]
self.QL.bulkvalues[k] = 0.0
self.H.bulkvalues[k] = GMV.H.values[k]
self.T.bulkvalues[k] = GMV.T.values[k]
self.B.bulkvalues[k] = 0.0
self.W.bulkvalues[k] = 0.0
return
# quick utility to set "new" arrays with values in the "values" arrays
cpdef set_new_with_values(self):
with nogil:
for i in xrange(self.n_updrafts):
for k in xrange(self.Gr.nzg):
self.W.new[i,k] = self.W.values[i,k]
self.Area.new[i,k] = self.Area.values[i,k]
self.QT.new[i,k] = self.QT.values[i,k]
self.QL.new[i,k] = self.QL.values[i,k]
self.QR.new[i,k] = self.QR.values[i,k]
self.H.new[i,k] = self.H.values[i,k]
self.THL.new[i,k] = self.THL.values[i,k]
self.T.new[i,k] = self.T.values[i,k]
self.B.new[i,k] = self.B.values[i,k]
return
# quick utility to set "new" arrays with values in the "values" arrays
cpdef set_old_with_values(self):
with nogil:
for i in xrange(self.n_updrafts):
for k in xrange(self.Gr.nzg):
self.W.old[i,k] = self.W.values[i,k]
self.Area.old[i,k] = self.Area.values[i,k]
self.QT.old[i,k] = self.QT.values[i,k]
self.QL.old[i,k] = self.QL.values[i,k]
self.QR.old[i,k] = self.QR.values[i,k]
self.H.old[i,k] = self.H.values[i,k]
self.THL.old[i,k] = self.THL.values[i,k]
self.T.old[i,k] = self.T.values[i,k]
self.B.old[i,k] = self.B.values[i,k]
return
# quick utility to set "tmp" arrays with values in the "new" arrays
cpdef set_values_with_new(self):
with nogil:
for i in xrange(self.n_updrafts):
for k in xrange(self.Gr.nzg):
self.W.values[i,k] = self.W.new[i,k]
self.Area.values[i,k] = self.Area.new[i,k]
self.QT.values[i,k] = self.QT.new[i,k]
self.QL.values[i,k] = self.QL.new[i,k]
self.QR.values[i,k] = self.QR.new[i,k]
self.H.values[i,k] = self.H.new[i,k]
self.THL.values[i,k] = self.THL.new[i,k]
self.T.values[i,k] = self.T.new[i,k]
self.B.values[i,k] = self.B.new[i,k]
return
cpdef io(self, NetCDFIO_Stats Stats):
Stats.write_profile('updraft_area', self.Area.bulkvalues[self.Gr.gw:self.Gr.nzg-self.Gr.gw])
Stats.write_profile('updraft_w', self.W.bulkvalues[self.Gr.gw:self.Gr.nzg-self.Gr.gw])
Stats.write_profile('updraft_qt', self.QT.bulkvalues[self.Gr.gw:self.Gr.nzg-self.Gr.gw])
Stats.write_profile('updraft_ql', self.QL.bulkvalues[self.Gr.gw:self.Gr.nzg-self.Gr.gw])
Stats.write_profile('updraft_qr', self.QR.bulkvalues[self.Gr.gw:self.Gr.nzg-self.Gr.gw])
if self.H.name == 'thetal':
Stats.write_profile('updraft_thetal', self.H.bulkvalues[self.Gr.gw:self.Gr.nzg-self.Gr.gw])
else:
Stats.write_profile('updraft_s', self.H.bulkvalues[self.Gr.gw:self.Gr.nzg-self.Gr.gw])
#Stats.write_profile('updraft_thetal', self.THL.bulkvalues[self.Gr.gw:self.Gr.nzg-self.Gr.gw])
Stats.write_profile('updraft_temperature', self.T.bulkvalues[self.Gr.gw:self.Gr.nzg-self.Gr.gw])
Stats.write_profile('updraft_buoyancy', self.B.bulkvalues[self.Gr.gw:self.Gr.nzg-self.Gr.gw])
self.get_cloud_base_top_cover()
# Note definition of cloud cover : each updraft is associated with a cloud cover equal to the maximum
# area fraction of the updraft where ql > 0. Each updraft is assumed to have maximum overlap with respect to
# itself (i.e. no consideration of tilting due to shear) while the updraft classes are assumed to have no overlap
# at all. Thus total updraft cover is the sum of each updraft's cover
Stats.write_ts('updraft_cloud_cover', np.sum(self.cloud_cover))
Stats.write_ts('updraft_cloud_base', np.amin(self.cloud_base))
Stats.write_ts('updraft_cloud_top', np.amax(self.cloud_top))
return
cpdef get_cloud_base_top_cover(self):
cdef Py_ssize_t i, k
for i in xrange(self.n_updrafts):
# Todo check the setting of ghost point z_half
self.cloud_base[i] = self.Gr.z_half[self.Gr.nzg-self.Gr.gw-1]
self.cloud_top[i] = 0.0
self.cloud_cover[i] = 0.0
for k in xrange(self.Gr.gw,self.Gr.nzg-self.Gr.gw):
if self.QL.values[i,k] > 1e-8 and self.Area.values[i,k] > 1e-3:
self.cloud_base[i] = fmin(self.cloud_base[i], self.Gr.z_half[k])
self.cloud_top[i] = fmax(self.cloud_top[i], self.Gr.z_half[k])
self.cloud_cover[i] = fmax(self.cloud_cover[i], self.Area.values[i,k])
return
cdef class UpdraftThermodynamics:
def __init__(self, n_updraft, Grid.Grid Gr, ReferenceState.ReferenceState Ref, UpdraftVariables UpdVar):
self.Gr = Gr
self.Ref = Ref
self.n_updraft = n_updraft
if UpdVar.H.name == 's':
self.t_to_prog_fp = t_to_entropy_c
self.prog_to_t_fp = eos_first_guess_entropy
elif UpdVar.H.name == 'thetal':
self.t_to_prog_fp = t_to_thetali_c
self.prog_to_t_fp = eos_first_guess_thetal
return
cpdef satadjust(self, UpdraftVariables UpdVar):
#Update T, QL
cdef:
Py_ssize_t k, i
eos_struct sa
with nogil:
for i in xrange(self.n_updraft):
for k in xrange(self.Gr.nzg):
sa = eos(self.t_to_prog_fp,self.prog_to_t_fp, self.Ref.p0_half[k],
UpdVar.QT.values[i,k], UpdVar.H.values[i,k])
UpdVar.QL.values[i,k] = sa.ql
UpdVar.T.values[i,k] = sa.T
return
cpdef buoyancy(self, UpdraftVariables UpdVar, EnvironmentVariables EnvVar,GridMeanVariables GMV, bint extrap):
cdef:
Py_ssize_t k, i
double alpha, qv, qt, t, h
Py_ssize_t gw = self.Gr.gw
UpdVar.Area.bulkvalues = np.sum(UpdVar.Area.values,axis=0)
if not extrap:
with nogil:
for i in xrange(self.n_updraft):
for k in xrange(self.Gr.nzg):
qv = UpdVar.QT.values[i,k] - UpdVar.QL.values[i,k]
alpha = alpha_c(self.Ref.p0_half[k], UpdVar.T.values[i,k], UpdVar.QT.values[i,k], qv)
UpdVar.B.values[i,k] = buoyancy_c(self.Ref.alpha0_half[k], alpha) #- GMV.B.values[k]
else:
with nogil:
for i in xrange(self.n_updraft):
for k in xrange(self.Gr.gw, self.Gr.nzg-self.Gr.gw):
if UpdVar.Area.values[i,k] > 1e-3:
qt = UpdVar.QT.values[i,k]
qv = UpdVar.QT.values[i,k] - UpdVar.QL.values[i,k]
h = UpdVar.H.values[i,k]
t = UpdVar.T.values[i,k]
alpha = alpha_c(self.Ref.p0_half[k], t, qt, qv)
UpdVar.B.values[i,k] = buoyancy_c(self.Ref.alpha0_half[k], alpha)
else:
sa = eos(self.t_to_prog_fp, self.prog_to_t_fp, self.Ref.p0_half[k],
qt, h)
qt -= sa.ql
qv = qt
t = sa.T
alpha = alpha_c(self.Ref.p0_half[k], t, qt, qv)
UpdVar.B.values[i,k] = buoyancy_c(self.Ref.alpha0_half[k], alpha)
with nogil:
for k in xrange(self.Gr.gw, self.Gr.nzg-self.Gr.gw):
GMV.B.values[k] = (1.0 - UpdVar.Area.bulkvalues[k]) * EnvVar.B.values[k]
for i in xrange(self.n_updraft):
GMV.B.values[k] += UpdVar.Area.values[i,k] * UpdVar.B.values[i,k]
for i in xrange(self.n_updraft):
UpdVar.B.values[i,k] -= GMV.B.values[k]
EnvVar.B.values[k] -= GMV.B.values[k]
return
#Implements a simple "microphysics" that clips excess humidity above a user-specified level
cdef class UpdraftMicrophysics:
def __init__(self, paramlist, n_updraft, Grid.Grid Gr, ReferenceState.ReferenceState Ref):
self.Gr = Gr
self.Ref = Ref
self.n_updraft = n_updraft
self.max_supersaturation = paramlist['turbulence']['updraft_microphysics']['max_supersaturation']
self.prec_source_h = np.zeros((n_updraft, Gr.nzg), dtype=np.double, order='c')
self.prec_source_qt = np.zeros((n_updraft, Gr.nzg), dtype=np.double, order='c')
self.prec_source_h_tot = np.zeros((Gr.nzg,), dtype=np.double, order='c')
self.prec_source_qt_tot = np.zeros((Gr.nzg,), dtype=np.double, order='c')
return
cpdef compute_sources(self, UpdraftVariables UpdVar):
"""
Compute precipitation source terms for QT, QR and H
"""
cdef:
Py_ssize_t k, i
double tmp_qr
with nogil:
for i in xrange(self.n_updraft):
for k in xrange(self.Gr.nzg):
tmp_qr = acnv_instant(UpdVar.QL.values[i,k], UpdVar.QT.values[i,k], self.max_supersaturation,\
UpdVar.T.values[i,k], self.Ref.p0_half[k])
self.prec_source_qt[i,k] = -tmp_qr
self.prec_source_h[i,k] = rain_source_to_thetal(self.Ref.p0_half[k], UpdVar.T.values[i,k],\
UpdVar.QT.values[i,k], UpdVar.QL.values[i,k], 0.0, tmp_qr)
#TODO assumes no ice
self.prec_source_h_tot = np.sum(np.multiply(self.prec_source_h, UpdVar.Area.values), axis=0)
self.prec_source_qt_tot = np.sum(np.multiply(self.prec_source_qt, UpdVar.Area.values), axis=0)
return
cpdef update_updraftvars(self, UpdraftVariables UpdVar):
"""
Apply precipitation source terms to QL, QR and H
"""
cdef:
Py_ssize_t k, i
with nogil:
for i in xrange(self.n_updraft):
for k in xrange(self.Gr.nzg):
UpdVar.QT.values[i,k] += self.prec_source_qt[i,k]
UpdVar.QL.values[i,k] += self.prec_source_qt[i,k]
UpdVar.QR.values[i,k] -= self.prec_source_qt[i,k]
UpdVar.H.values[i,k] += self.prec_source_h[i,k]
return
cdef void compute_update_combined_local_thetal(self, double p0, double T, double *qt, double *ql, double *qr, double *h,
Py_ssize_t i, Py_ssize_t k) nogil :
# Language note: array indexing must be used to dereference pointers in Cython. * notation (C-style dereferencing)
# is reserved for packing tuples
tmp_qr = acnv_instant(ql[0], qt[0], self.max_supersaturation, T, p0)
self.prec_source_qt[i,k] = -tmp_qr
self.prec_source_h[i,k] = rain_source_to_thetal(p0, T, qt[0], ql[0], 0.0, tmp_qr)
#TODO - assumes no ice
qt[0] += self.prec_source_qt[i,k]
ql[0] += self.prec_source_qt[i,k]
qr[0] -= self.prec_source_qt[i,k]
h[0] += self.prec_source_h[i,k]
return