-
Notifications
You must be signed in to change notification settings - Fork 20
/
gfs_pv_1.4_north_pacific_anim.py
326 lines (268 loc) · 9.56 KB
/
gfs_pv_1.4_north_pacific_anim.py
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
#
# run on python 3.7
#
# python code for some calculations related to the dynamic tropopause (DT) -
# DT pressure, DT potential temperature, 330K PV,
# and a cross-section of PV at the latitude where the tropopause is lowest -
# all based on the GFS analysis available online. As the data is accessed
# online, the program can take a while to run.
#
# the date and lat-lon range can be set below
#
# (poorly) coded by Mathew Barlow
# initial release: 14 Nov 2017
# last updated: 10 Oct 2019
#
# this code has *not* been extensively tested and has been
# awkwardly translated from other coding languages, so if you find
# any errors or have any suggestions or improvements, including for
# the plotting, please let me know at [email protected] . Thanks!
#
# Support from NSF AGS-1623912 is gratefully acknowledged
#
import numpy as np
import netCDF4
import matplotlib.pyplot as plt
import matplotlib.ticker as tick
from mpl_toolkits.mplot3d import axes3d
import cartopy.crs as ccrs
from scipy.ndimage import gaussian_filter
from cartopy.feature import NaturalEarthFeature
from cartopy.mpl.gridliner import LONGITUDE_FORMATTER, LATITUDE_FORMATTER
from cartopy.mpl.ticker import LatitudeFormatter, LongitudeFormatter
from datetime import datetime
# VALUES TO SET *************************************************
# set date, lat-lon range, and PV-value definition of tropopause
mydate='20191016'
myhour='06'
(lat1,lat2)=(20,70)
(lon1,lon2)=(120,240.1)
tpdef=2 # definition of tropopause in PVU
#****************************************************************
#constants
re=6.37e6
g=9.81
cp=1004.5
r=2*cp/7
kap=r/cp
omega=7.292e-5
pi=3.14159265
# open dataset, retreive variables, close dataset
url='https://nomads.ncep.noaa.gov/dods/gfs_0p25/gfs'+\
mydate+'/gfs_0p25_'+myhour+'z_anl'
file = netCDF4.Dataset(url)
lat_in = file.variables['lat'][:]
lon_in = file.variables['lon'][:]
lev = file.variables['lev'][:]
pres2pv_in = file.variables['pres2pv'][0,:,:]
pressfc_in = file.variables['pressfc'][0,:,:]
nlev = lev.size
nx = lon_in.size
ny = lat_in.size
u_in = np.full((nlev, ny, nx), None)
v_in = np.full((nlev, ny, nx), None)
t_in = np.full((nlev, ny, nx), None)
hgt_in = np.full((nlev, ny, nx), None)
ilev = 0
while ilev < nlev:
print(ilev)
u_in[ilev, :, :] = file.variables['ugrdprs'][0, ilev, :, :]
ilev = ilev + 1
ilev = 0
while ilev < nlev:
v_in[ilev, :, :] = file.variables['vgrdprs'][0, ilev, :, :]
ilev = ilev + 1
ilev = 0
while ilev < nlev:
t_in[ilev, :, :] = file.variables['tmpprs'][0, ilev, :, :]
ilev = ilev + 1
ilev = 0
while ilev < nlev:
hgt_in[ilev, :, :] = file.variables['hgtprs'][0, ilev, :, :]
ilev = ilev + 1
#t_in = file.variables['tmpprs'][0,:,:,:]
#u_in = file.variables['ugrdprs'][0,:,:,:]
#v_in = file.variables['vgrdprs'][0,:,:,:]
#hgt_in = file.variables['hgtprs'][0,:,:,:]
file.close()
# get array indices for lat-lon range
# specified above
iy1 = np.argmin( np.abs( lat_in - lat1 ) )
iy2 = np.argmin( np.abs( lat_in - lat2 ) )
ix1 = np.argmin( np.abs( lon_in - lon1 ) )
ix2 = np.argmin( np.abs( lon_in - lon2 ) )
# select specified lat-lon range
t=t_in[:,iy1:iy2,ix1:ix2]
lon=lon_in[ix1:ix2]
lat=lat_in[iy1:iy2]
u=u_in[:,iy1:iy2,ix1:ix2]
v=v_in[:,iy1:iy2,ix1:ix2]
hgt=hgt_in[:,iy1:iy2,ix1:ix2]
pres2pv=pres2pv_in[iy1:iy2,ix1:ix2]
pressfc=pressfc_in[iy1:iy2,ix1:ix2]
# some prep work for derivatives
xlon,ylat=np.meshgrid(lon,lat)
# define potential temperature and Coriolis parameter
theta=t*(1.E5/(lev[:,np.newaxis,np.newaxis]*100))**kap
f=2*omega*np.sin(ylat*pi/180)
lon = np.array(lon, dtype='float')
lat = np.array(lat, dtype='float')
lev = np.array(lev, dtype='float')
u = np.array(u, dtype='float')
v = np.array(v, dtype='float')
hgt = np.array(hgt, dtype='float')
pres2pv = np.array(pres2pv, dtype='float')
pressfc = np.array(pressfc, dtype='float')
theta = np.array(theta, dtype='float')
f = np.array(f, dtype='float')
# calculate derivatives
def ddp(f):
# handle unevenly-spaced levels with 2nd order
# Lagrange interpolation
# except for top and bottom, where use forward diff
lev3=lev.reshape(lev.size,1,1)*100
dpp=lev3-np.roll(lev3,-1,axis=0)
dpm=lev3-np.roll(lev3,1,axis=0)
fp=np.roll(f,-1,axis=0)
fm=np.roll(f,1,axis=0)
ddp_f=(
fm*dpp/( (dpp-dpm)*(-dpm) ) +
f*(dpp+dpm)/( dpm*dpp ) +
fp*dpm/( (dpm-dpp)*(-dpp) )
)
ddp_f[0,:,:]=(f[1,:,:]-f[0,:,:])/(lev3[1,:,:]-lev3[0,:,:])
ddp_f[-1,:,:]=(f[-1,:,:]-f[-2,:,:])/(lev3[-2,:,:]-lev3[-1,:,:])
return(ddp_f)
def ddx(f):
# use center-difference, assuming evenly spaced lon
# except for side-boundaries, where use forward diff
x=(re*np.cos(ylat*np.pi/180)*np.pi/180)*lon
x3=x.reshape(1,x.shape[0],x.shape[1])
dx3=np.roll(x3,-1,axis=2)-np.roll(x3,1,axis=2)
ddx_f=(np.roll(f,-1,axis=2)-np.roll(f,1,axis=2))/dx3
ddx_f[:,:,0]=(f[:,:,1]-f[:,:,0])/(x3[:,:,1]-x3[:,:,0])
ddx_f[:,:,-1]=(f[:,:,-2]-f[:,:,-1])/(x3[:,:,-2]-x3[:,:,-1])
return(ddx_f)
def ddy(f):
# use center-difference, assuming evenly spaced lon
# except for N/S boundaries, where use forward diff
y=(re*np.pi/180)*lat
y3=y.reshape(1,y.shape[0],1)
dy3=np.roll(y3,-1,axis=1)-np.roll(y3,1,axis=1)
ddy_f=(np.roll(f,-1,axis=1)-np.roll(f,1,axis=1))/dy3
ddy_f[:,0,:]=(f[:,1,:]-f[:,0,:])/(y3[:,1,:]-y3[:,0,:])
ddy_f[:,-1,:]=(f[:,-2,:]-f[:,-1,:])/(y3[:,-2,:]-y3[:,-1,:])
return(ddy_f)
#lev3=lev.reshape(lev.size,1,1)
#ddp_theta=np.gradient(theta,lev3*100,axis=0)
#ddx_theta=np.gradient(theta,axis=2)/dx
#ddy_theta=np.gradient(theta,axis=1)/dy
gf=1
ddp_theta=ddp(theta)
ddp_u=ddp(gaussian_filter(u,sigma=gf))
ddp_v=ddp(gaussian_filter(v,sigma=gf))
ddx_theta=ddx(theta)
ddy_theta=ddy(theta)
ddx_v=ddx(gaussian_filter(v,sigma=gf))
ddy_ucos=ddy(gaussian_filter(u,sigma=gf)*np.cos(ylat*pi/180))
# calculate contributions to PV and PV
absvort=ddx_v-(1/np.cos(ylat*pi/180))*ddy_ucos+f
pv_one=g*absvort*(-ddp_theta)
pv_two=g*(ddp_v*ddx_theta-ddp_u*ddy_theta)
pv=pv_one+pv_two
# calculate pressure of tropopause, Fortran-style (alas!)
# as well as potential temperature (theta) and height
#
# starting from 10hPa and working down, to avoid
# more complicated vertical structure higher up
#
nx=ix2-ix1+1
ny=iy2-iy1+1
nz=lev.size
nzs=np.argwhere(lev==50.0)[0,0]
tp=np.empty((ny-1,nx-1))*np.nan # initialize as undef
tp_theta=np.empty((ny-1,nx-1))*np.nan # initialize as undef
tp_hgt=np.empty((ny-1,nx-1))*np.nan # initialize as undef
for ix in range(0,nx-1):
for iy in range(0,ny-1):
for iz in range(nzs,0,-1):
if pv[iz,iy,ix]/1e-6<=tpdef:
if np.isnan(tp[iy,ix]):
tp[iy,ix]=(
(lev[iz]*(pv[iz+1,iy,ix]-tpdef*1e-6)
-lev[iz+1]*(pv[iz,iy,ix]-tpdef*1e-6))/
(pv[iz+1,iy,ix]-pv[iz,iy,ix])
)
tp_theta[iy,ix]=(
((lev[iz]-tp[iy,ix])*theta[iz+1,iy,ix]+
(tp[iy,ix]-lev[iz+1])*theta[iz,iy,ix])/
(lev[iz]-lev[iz+1])
)
tp_hgt[iy,ix]=(
((lev[iz]-tp[iy,ix])*hgt[iz+1,iy,ix]+
(tp[iy,ix]-lev[iz+1])*hgt[iz,iy,ix])/
(lev[iz]-lev[iz+1])
)
# calculate PV on the 330K isentropic surface
# (also not in a pythonic way)
nx=ix2-ix1+1
ny=iy2-iy1+1
nz=lev.size
pv330=np.empty((ny-1,nx-1))*np.nan # initialize as undef
for ix in range(0,nx-1):
for iy in range(0,ny-1):
for iz in range(nz-2,0,-1):
if theta[iz,iy,ix]>=330:
if theta[iz-1,iy,ix]<=330:
if np.isnan(pv330[iy,ix]):
pv330[iy,ix]=(
((330-theta[iz-1,iy,ix])*pv[iz,iy,ix]+
(theta[iz,iy,ix]-330)*pv[iz-1,iy,ix])/
(theta[iz,iy,ix]-theta[iz-1,iy,ix])
)
# slight smoothing of result
# (appears to work better than smoothing u,v,t first)
tp=gaussian_filter(tp,sigma=1)
tp_theta=gaussian_filter(tp_theta,sigma=1)
pv330=gaussian_filter(pv330,sigma=1)
# define spatial correlation function for testing results
def scorr(a,b):
abar=np.mean(a)
bbar=np.mean(b)
covar=sum((a-abar)*(b-bbar))
avar=sum((a-abar)**2)
bvar=sum((b-bbar)**2)
r=covar/np.sqrt(avar*bvar)
return(r)
# identify latitude of lowest tropopause
maxloc=np.argwhere(tp==np.amax(tp))
latmax=lat[maxloc[0,0]]
# now make some plots - these badly need to be improved
states = NaturalEarthFeature(category='cultural',
scale='50m', facecolor='none',
name='admin_1_states_provinces_shp')
# get date for plotting
fdate=datetime.strptime(mydate, '%Y%m%d').strftime('%d %b %Y')
plt.close(fig='all')
print('got here')
nframe=30
iframe=0
while iframe<=nframe:
plt.figure(iframe,figsize=plt.figaspect(0.5))
pressfc_smooth=gaussian_filter(pressfc,sigma=1)
ax=plt.gca(projection='3d')
surf=ax.plot_surface(xlon,ylat,tp,cmap="coolwarm",alpha=1,
rstride=1,cstride=1,
linewidth=0, antialiased=False)
ax.plot_surface(xlon,ylat,pressfc_smooth/100,color="lightgray",
rstride=1,cstride=1,
linewidth=0, antialiased=False)
ax.set_zlim(1000,100)
ax.set_xlim(lon1,lon2)
ax.set_ylim(lat1,lat2)
ax.view_init(elev=90 - iframe*90/nframe,azim=-90)
plt.title('2PVU Dynamic Tropopause over topography\n'+myhour+'Z '+fdate)
plt.colorbar(surf)
plt.savefig('goo'+ '{:04d}'.format(iframe)+'.png', bbox_inches='tight',
dpi=300)
iframe=iframe+1