-
Notifications
You must be signed in to change notification settings - Fork 2
/
Copy path__init__.py
219 lines (195 loc) · 7.27 KB
/
__init__.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
import xarray as xr
import numpy as np
from lib.log import log
import os
from datetime import timedelta, datetime
from cli.applications.croco.post_process_v1.depth_functions import z_levels
from cli.applications.croco.post_process_v1.functions import (
hour_rounder,
u2rho_4d,
v2rho_4d,
)
import subprocess
# All dates in the CROCO output are represented
# in seconds from 1 Jan 2000 (i.e. the reference date)
REFERENCE_DATE = datetime(2000, 1, 1, 0, 0, 0)
# Model variables use the dimensions time (time from reference date),
# eta_rho (lat) and xi_rho (lon). We are changing eta_rho and xi_rho
# from grid points to real lat and lon data.
def post_process_v1(args):
id = args.id
grid = os.path.abspath(args.grid)
input = os.path.abspath(args.input)
output = os.path.abspath(args.output)
log("Running CROCO output post-processing (v1)")
log("CONFIG::id", id)
log("CONFIG::input", input)
log("CONFIG::grid", grid)
log("CONFIG::output", output)
# Ensure the directory for the specified output exists
os.makedirs(os.path.dirname(output), exist_ok=True)
data = xr.open_dataset(input)
data_grid = xr.open_dataset(grid)
# Dimensions that need to be transformed
time = data.time.values # Time steps
lon_rho = data.lon_rho.values # Longitude (4326)
lat_rho = data.lat_rho.values # Latitude (4326)
# Convert time to human readable
time_steps = []
for t in time:
date_now = REFERENCE_DATE + timedelta(seconds=np.float64(t))
date_round = hour_rounder(date_now)
time_steps.append(date_round)
# Get the run_date from the CROCO output NetCDF file
run_date = time_steps[119].strftime("%Y%m%d")
log("CONFIG::run_date", run_date)
# Variables used in the visualisations
temperature = data.temp.values
salt = data.salt.values
ssh = data.zeta.values # Sea-surface height
u = data.u.values # East-West velocity
v = data.v.values # North-South velocity
h = data.h.values # Bathymetry elevation model
# Variables used to calculate depth levels
# CROCO uses these params to determine how to deform the grid
s_rho = data.s_rho.values # Vertical levels
theta_s = data.theta_s
theta_b = data.theta_b
# Variables used to calculate depth levels from grid (bathymetry)
h = data_grid.h.values
# Convert u and v current components to the rho grid (otherwise these values are on cell edges and not the center of the cells)
# use the function u2rho_4d and v2rho_4d
log("Normalizing u/v model output variables")
u_rho = u2rho_4d(u)
v_rho = v2rho_4d(v)
log("Masking 0 values (land) from variables")
temperature[np.where(temperature == 0)] = np.nan
salt[np.where(salt == 0)] = np.nan
u[np.where(u == 0)] = np.nan
v[np.where(v == 0)] = np.nan
# Variables hard coded set during model configuration
# Relative to each model
# Critical depth. Convergence of surface layers with bottom layers due to grid squeezing
hc = data.hc.values
N = np.shape(data.s_rho)[0]
type_coordinate = "rho"
vtransform = (
2 if data.VertCoordType == "NEW" else 1 if data.VertCoordType == "OLD" else -1
)
if vtransform == -1:
raise Exception("Unexpected value for vtransform (" + vtransform + ")")
# m_rho refers to the depth level in meters
log("Converting depth levels to depth in meters")
m_rho = np.zeros(np.shape(temperature))
for x in np.arange(np.size(temperature, 0)):
depth_temp = z_levels(
h, ssh[x, :, :], theta_s, theta_b, hc, N, type_coordinate, vtransform
)
m_rho[x, ::] = depth_temp
# Create new xarray dataset with selected variables
log("Generating dataset")
data_out = xr.Dataset(
attrs={
"description": "CROCO output from algoa Bay model transformed lon/lat/depth/time",
"model_name": id,
"run_date": run_date,
},
data_vars={
"temperature": xr.Variable(
["time", "depth", "lat", "lon"],
temperature,
{
"long_name": data.temp.long_name,
"units": data.temp.units,
"description": "Temperature at grid points",
},
),
"salt": xr.Variable(
["time", "depth", "lat", "lon"],
salt,
{
"long_name": data.salt.long_name,
"units": data.salt.units,
"description": "Salinity at grid points",
},
),
"u": xr.Variable(
["time", "depth", "lat", "lon"],
u_rho,
{
"long_name": data.u.long_name,
"units": data.u.units,
"description": "@MATT TODO at grid points",
},
),
"v": xr.Variable(
["time", "depth", "lat", "lon"],
v_rho,
{
"long_name": data.v.long_name,
"units": data.v.units,
"description": "@MATT TODO at grid points",
},
),
"m_rho": xr.Variable(
["time", "depth", "lat", "lon"],
m_rho,
{
"description": "Depth level of grid points in meters, centred in grid cells"
},
),
"h": xr.Variable(
["lat", "lon"],
h,
{
"long_name": data.h.long_name,
"units": data.h.units,
"description": "Bathymetry elevation at grid XY points",
},
),
},
coords={
"lon_rho": xr.Variable(
["lat", "lon"],
lon_rho,
{
"long_name": data.lon_rho.long_name,
"units": data.lon_rho.units,
"description": "Longitude coordinate values of curvilinear grid cells",
},
),
"lat_rho": xr.Variable(
["lat", "lon"],
lat_rho,
{
"long_name": data.lat_rho.long_name,
"units": data.lat_rho.units,
"description": "Latitude coordinate values of curvilinear grid cells",
},
),
"depth": xr.Variable(
["depth"], s_rho, {"description": "Depth levels (grid levels)"}
),
"time": xr.Variable(
["time"],
time_steps,
{"description": "Time steps in hours"},
),
},
)
encoding = {
"temperature": {"dtype": "float32"},
"salt": {"dtype": "float32"},
"u": {"dtype": "float32"},
"v": {"dtype": "float32"},
"m_rho": {"dtype": "float32"},
"lon_rho": {"dtype": "float32"},
"lat_rho": {"dtype": "float32"},
"depth": {"dtype": "u2", "_FillValue": 65535},
"time": {"dtype": "i4"},
"h": {"dtype": "float32"},
}
log("Writing NetCDF file")
data_out.to_netcdf(output, encoding=encoding, mode="w")
subprocess.call(["chmod", "-R", "775", output])
log("Done!")