-
Notifications
You must be signed in to change notification settings - Fork 34
/
decimateGrid.py
337 lines (270 loc) · 13.2 KB
/
decimateGrid.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
326
327
328
329
330
331
332
333
334
from netCDF4 import Dataset
from netCDF4 import num2date
import numpy as np
import time
import os
import shutil
import calculateGRDMetrics
def createGrid(grdROMS,infilename,outfilename,decimate):
shutil.copy2(infilename, outfilename)
f1 = Dataset(outfilename, mode='w', format='NETCDF4')
f1.description="This is a grid file for ROMS - KINO project"
f1.history = 'Created (decimated) from IMR 800M grid' + time.ctime(time.time())
f1.source = 'Trond Kristiansen ([email protected])'
f1.type='NetCDF4 classic created using MODEL2ROMS - https://github.com/trondkr/model2roms'
if not int(grdROMS.xi_rho) % 2 == 0:
deltaXI=1
else:
deltaXI=0
if not int(grdROMS.eta_rho) % 2 == 0:
deltaETA=1
startindex=0;endindex=-2
else:
deltaETA=0
startindex=0;endindex=-1
print('Old dimensions were (XI_RHO, ETA_RHO) : %ix%i'%(int(grdROMS.xi_rho)-deltaXI,int(grdROMS.eta_rho)-deltaETA))
print('Old dimensions were (XI_U, ETA_U) : %ix%i'%(int(grdROMS.xi_u)-deltaXI,int(grdROMS.eta_u)-deltaETA))
print('Old dimensions were (XI_V, ETA_V) : %ix%i'%(int(grdROMS.xi_v)-deltaXI,int(grdROMS.eta_v)-deltaETA))
print('Old dimensions were (XI_PSI, ETA_PSI) : %ix%i'%(int(grdROMS.xi_psi)-deltaXI,int(grdROMS.eta_psi)-deltaETA))
print('New dimensions will be (XI_RHO, ETA_RHO) : %ix%i'%(((int(grdROMS.xi_rho)-deltaXI)/decimate),((int(grdROMS.eta_rho)-deltaETA)/decimate)))
print('New dimensions will be (XI_U, ETA_U) : %ix%i'%(((int(grdROMS.xi_u)-deltaXI)/decimate),((int(grdROMS.eta_u)-deltaETA)/decimate)))
print('New dimensions will be (XI_V, ETA_V) : %ix%i'%(((int(grdROMS.xi_v)-deltaXI)/decimate),((int(grdROMS.eta_v)-deltaETA)/decimate)))
print('New dimensions will be (XI_PSI, ETA_PSI) : %ix%i'%(((int(grdROMS.xi_psi)-deltaXI)/decimate),((int(grdROMS.eta_psi)-deltaETA)/decimate)))
# Define dimensions
xi_rho=(int(grdROMS.xi_rho)-deltaXI)/decimate
xi_vert=((int(grdROMS.xi_rho)-deltaXI)/decimate) + 1
eta_rho=(int(grdROMS.eta_rho)-deltaETA)/decimate
eta_vert=((int(grdROMS.eta_rho)-deltaETA)/decimate) + 1
xi_u=(int(grdROMS.xi_u)-deltaXI)/decimate
eta_u=(int(grdROMS.eta_u)-deltaETA)/decimate
xi_v=(int(grdROMS.xi_v)-deltaXI)/decimate
eta_v=(int(grdROMS.eta_v)-deltaETA)/decimate
xi_psi=(int(grdROMS.xi_psi)-deltaXI)/decimate
eta_psi=(int(grdROMS.eta_psi)-deltaETA)/decimate
s_rho=int(len(grdROMS.s_rho))
s_w=int(len(grdROMS.s_w))
f1.createDimension('xi_rho', xi_rho)
f1.createDimension('eta_rho', eta_rho)
f1.createDimension('xi_u', xi_u)
f1.createDimension('eta_u', eta_u)
f1.createDimension('xi_v', xi_v)
f1.createDimension('eta_v', eta_v)
f1.createDimension('xi_psi', xi_psi)
f1.createDimension('eta_psi', eta_psi)
f1.createDimension('xi_vert', xi_vert)
f1.createDimension('eta_vert', eta_vert)
f1.createDimension('s_rho', s_rho)
f1.createDimension('s_w', s_w)
f1.createDimension('bath', None)
vnc = f1.createVariable('lon_rho', 'd', ('eta_rho','xi_rho',),zlib=False)
vnc.long_name = 'Longitude at RHO points'
vnc.units = 'degree_east'
vnc[:,:] = grdROMS.lon_rho[startindex:endindex:decimate,startindex:endindex:decimate]
vnc = f1.createVariable('lat_rho', 'd', ('eta_rho','xi_rho',),zlib=False)
vnc.long_name = 'Latitude at RHO points'
vnc.units = 'degree_north'
vnc[:,:] = grdROMS.lat_rho[startindex:endindex:decimate,startindex:endindex:decimate]
vnc = f1.createVariable('lon_u', 'd', ('eta_u','xi_u',),zlib=False)
vnc.long_name = 'Longitude at U points'
vnc.units = 'degree_east'
vnc[:,:] = grdROMS.lon_u[startindex:endindex:decimate,startindex:endindex:decimate]
vnc = f1.createVariable('lat_u', 'd', ('eta_u','xi_u',),zlib=False)
vnc.long_name = 'Latitude at U points'
vnc.units = 'degree_north'
vnc[:,:] = grdROMS.lat_u[startindex:endindex:decimate,startindex:endindex:decimate]
vnc = f1.createVariable('lon_v', 'd', ('eta_v','xi_v',),zlib=False)
vnc.long_name = 'Longitude at V points'
vnc.units = 'degree_east'
vnc[:,:] = grdROMS.lon_v[startindex:endindex:decimate,startindex:endindex:decimate]
vnc = f1.createVariable('lat_v', 'd', ('eta_v','xi_v',),zlib=False)
vnc.long_name = 'Latitude at V points'
vnc.units = 'degree_north'
vnc[:,:] = grdROMS.lat_v[startindex:endindex:decimate,startindex:endindex:decimate]
vnc = f1.createVariable('lat_psi', 'd', ('eta_psi','xi_psi',),zlib=False)
vnc.long_name = 'Latitude at PSI points'
vnc.units = 'degree_north'
vnc[:,:] = grdROMS.lat_psi[startindex:endindex:decimate,startindex:endindex:decimate]
vnc = f1.createVariable('lon_psi', 'd', ('eta_psi','xi_psi',),zlib=False)
vnc.long_name = 'Longitude at PSI points'
vnc.units = 'degree_east'
vnc[:,:] = grdROMS.lon_psi[startindex:endindex:decimate,startindex:endindex:decimate]
vnc = f1.createVariable('lon_vert', 'd', ('eta_vert','xi_vert',),zlib=False)
vnc.long_name = 'Longitude at vertices points'
vnc.units = 'degree_east'
vnc[:,:] = grdROMS.lon_vert[::decimate,::decimate]
vnc = f1.createVariable('lat_vert', 'd', ('eta_vert','xi_vert',),zlib=False)
vnc.long_name = 'Latitude at vertices points'
vnc.units = 'degree_north'
vnc[:,:] = grdROMS.lat_vert[::decimate,::decimate]
vnc = f1.createVariable('x_rho', 'd', ('eta_rho','xi_rho',),zlib=False)
vnc.long_name = 'X location of RHO points'
vnc.units = 'meter'
vnc[:,:] = grdROMS.x_rho[startindex:endindex:decimate,startindex:endindex:decimate]
vnc = f1.createVariable('y_rho', 'd', ('eta_rho','xi_rho',),zlib=False)
vnc.long_name = 'Y location of RHO points'
vnc.units = 'meter'
vnc[:,:] = grdROMS.y_rho[startindex:endindex:decimate,startindex:endindex:decimate]
vnc = f1.createVariable('x_u', 'd', ('eta_u','xi_u',),zlib=False)
vnc.long_name = 'X location of U points'
vnc.units = 'meter'
vnc[:,:] = grdROMS.x_u[startindex:endindex:decimate,startindex:endindex:decimate]
vnc = f1.createVariable('y_u', 'd', ('eta_u','xi_u',),zlib=False)
vnc.long_name = 'Y location of U points'
vnc.units = 'meter'
vnc[:,:] = grdROMS.y_u[startindex:endindex:decimate,startindex:endindex:decimate]
vnc = f1.createVariable('x_v', 'd', ('eta_v','xi_v',),zlib=False)
vnc.long_name = 'X location of V points'
vnc.units = 'meter'
vnc[:,:] = grdROMS.x_v[startindex:endindex:decimate,startindex:endindex:decimate]
vnc = f1.createVariable('y_v', 'd', ('eta_v','xi_v',),zlib=False)
vnc.long_name = 'Y location of V points'
vnc.units = 'meter'
vnc[:,:] = grdROMS.y_v[startindex:endindex:decimate,startindex:endindex:decimate]
vnc = f1.createVariable('x_psi', 'd', ('eta_psi','xi_psi',),zlib=False)
vnc.long_name = 'X location of PSI points'
vnc.units = 'meter'
vnc[:,:] = grdROMS.x_psi[startindex:endindex:decimate,startindex:endindex:decimate]
vnc = f1.createVariable('y_psi', 'd', ('eta_psi','xi_psi',),zlib=False)
vnc.long_name = 'Y location of PSI points'
vnc.units = 'meter'
vnc[:,:] = grdROMS.y_psi[startindex:endindex:decimate,startindex:endindex:decimate]
vnc = f1.createVariable('x_vert', 'd', ('eta_vert','xi_vert',),zlib=False)
vnc.long_name = 'X location of VERT points'
vnc.units = 'meter'
vnc[:,:] = grdROMS.x_vert[::decimate,::decimate]
vnc = f1.createVariable('y_vert', 'd', ('eta_vert','xi_vert',),zlib=False)
vnc.long_name = 'Y location of VERT points'
vnc.units = 'meter'
vnc[:,:] = grdROMS.y_vert[::decimate,::decimate]
vnc = f1.createVariable('h', 'd', ('eta_rho','xi_rho',),zlib=False)
vnc.long_name = 'Final bathymetry at RHO points'
vnc.units = 'meter'
vnc.field = "bath, scalar"
vnc[:,:] = grdROMS.h[startindex:endindex:decimate,startindex:endindex:decimate]
vnc = f1.createVariable('hraw', 'd', ('bath','eta_rho','xi_rho',),zlib=False)
vnc.long_name = 'Raw bathymetry at RHO-points'
vnc.units = 'meter'
vnc.field = "bath, scalar"
vnc[:,:,:] = grdROMS.hraw[:,startindex:endindex:decimate,startindex:endindex:decimate]
vnc = f1.createVariable('f', 'd', ('eta_rho','xi_rho',),zlib=False)
vnc.long_name = 'Coriolis parameter at RHO points'
vnc.units = 'second-1'
vnc.field = "Coriolis, scalar"
vnc[:,:] = grdROMS.f[startindex:endindex:decimate,startindex:endindex:decimate]
# Calculate grid metrics
dndx,dmde,pm,pn = calculateGRDMetrics.calculateGridMetrics(grdROMS,True,decimate,startindex,endindex)
vnc = f1.createVariable('pm', 'd', ('eta_rho','xi_rho',),zlib=False)
vnc.long_name = 'curvilinear coordinate metric in XI'
vnc.units = 'meter-1'
vnc.field = "pm, scalar"
vnc[:,:] = pm
print("Average DX in meters: ",1/np.average(pm,axis=None))
vnc = f1.createVariable('pn', 'd', ('eta_rho','xi_rho',),zlib=False)
vnc.long_name = 'curvilinear coordinate metric in ETA'
vnc.units = 'meter-1'
vnc.field = "pn, scalar"
vnc[:,:] = pn
print("Average DY in meters: ",1/np.average(pn,axis=None))
vnc = f1.createVariable('dmde', 'd', ('eta_rho','xi_rho',),zlib=False)
vnc.long_name = 'XI derivative of inverse metric factor pn'
vnc.units = 'meter'
vnc.field = "dmde, scalar"
vnc[:,:] = dmde
vnc = f1.createVariable('dndx', 'd', ('eta_rho','xi_rho',),zlib=False)
vnc.long_name = 'ETA derivative of inverse metric factor pm'
vnc.units = 'meter'
vnc.field = "dndx, scalar"
vnc[:,:] = dndx
vnc=f1.createVariable('angle','d',('eta_rho','xi_rho'),zlib=False)
vnc.long_name = "angle between xi axis and east"
vnc.units = "radian"
vnc[:,:]=grdROMS.angle[startindex:endindex:decimate,startindex:endindex:decimate]
vnc=f1.createVariable('mask_rho','d',('eta_rho', 'xi_rho'),zlib=False)
vnc.long_name = "mask on RHO-points"
vnc.option_0 = "land"
vnc.option_1 = "water"
vnc.FillValue = 1.0
values=grdROMS.mask_rho[startindex:endindex:decimate,startindex:endindex:decimate]
infile="/Users/trondkr/Projects/KINO/GRID/mask_change.txt"
fmask=open(infile,"r")
lines=fmask.readlines()
import string
for line in lines:
l=string.split(line," ")
i=int(float(l[0].strip()))
j=int(float(l[1].strip()))
m=float(l[2].strip())
print("Changing %s %s from %s to %s"%(j,i,values[j,i],m))
values[j,i]=m
vnc[:,:]=values[:,:]
vnc=f1.createVariable('mask_u','d',('eta_u', 'xi_u'),zlib=False)
vnc.long_name = "mask on U-points"
vnc.option_0 = "land"
vnc.option_1 = "water"
vnc.FillValue = 1.0
vnc[:,:]=grdROMS.mask_u[startindex:endindex:decimate,startindex:endindex:decimate]
vnc=f1.createVariable('mask_v','d',('eta_v', 'xi_v'),zlib=False)
vnc.long_name = "mask on V-points"
vnc.option_0 = "land"
vnc.option_1 = "water"
vnc.FillValue = 1.0
vnc[:,:]=grdROMS.mask_v[startindex:endindex:decimate,startindex:endindex:decimate]
vnc=f1.createVariable('mask_psi','d',('eta_psi', 'xi_psi'),zlib=False)
vnc.long_name = "mask on PSI-points"
vnc.option_0 = "land"
vnc.option_1 = "water"
vnc.FillValue = 1.0
vnc[:,:]=grdROMS.mask_psi[startindex:endindex:decimate,startindex:endindex:decimate]
vnc = f1.createVariable('s_rho', 'd', ('s_rho',), zlib=False)
vnc.long_name = "S-coordinate at RHO-points"
vnc.valid_min = -1.
vnc.valid_max = 0.
vnc.field = "s_rho, scalar"
vnc[:] = grdROMS.s_rho[:]
vnc = f1.createVariable('s_w', 'd', ('s_w',), zlib=False)
vnc.long_name = "S-coordinate at W-points"
vnc.valid_min = -1.
vnc.valid_max = 0.
vnc.field = "s_w, scalar"
vnc[:] = grdROMS.s_w[:]
vnc = f1.createVariable('Cs_r', 'd', ('s_rho',), zlib=False)
vnc.long_name = "S-coordinate stretching curves at RHO-points"
vnc.valid_min = -1.
vnc.valid_max = 0.
vnc.field = "s_rho, scalar"
vnc[:] = grdROMS.Cs_rho[:]
vnc = f1.createVariable('Cs_w', 'd', ('s_w',), zlib=False)
vnc.long_name = "S-coordinate stretching curves at W-points"
vnc.valid_min = -1.
vnc.valid_max = 0.
vnc.field = "s_w, scalar"
vnc[:] = grdROMS.Cs_w[:]
vnc = f1.createVariable('hc', 'd')
vnc.long_name = "S-coordinate parameter, critical depth";
vnc.units = "meter"
vnc[:] = grdROMS.hc
vnc = f1.createVariable('xl', 'd')
vnc.long_name = "domain length in the XI-direction";
vnc.units = "meter"
vnc[:] = grdROMS.xl
vnc = f1.createVariable('el', 'd')
vnc.long_name = "domain length in the ETA-direction";
vnc.units = "meter"
vnc[:] = grdROMS.el
vnc = f1.createVariable('Tcline', 'd')
vnc.long_name = "S-coordinate surface/bottom layer width";
vnc.units = "meter"
vnc[:] = grdROMS.tcline
vnc = f1.createVariable('theta_s', 'd')
vnc.long_name = "S-coordinate surface control parameter";
vnc[:] = grdROMS.theta_s
vnc = f1.createVariable('theta_b', 'd')
vnc.long_name = "S-coordinate bottom control parameter";
vnc[:] = grdROMS.theta_b
vnc = f1.createVariable('spherical', 'c')
vnc.long_name = "Grid type logical switch";
vnc.flag_values = "T, F" ;
vnc.flag_meanings = "spherical Cartesian" ;
vnc[:] = grdROMS.spherical
f1.close()
print("Creating new decimated grid file: %s"%(outfilename))