Skip to content

Commit 92977df

Browse files
committed
Added member and time dimension
1 parent 2e0aca2 commit 92977df

File tree

1 file changed

+69
-10
lines changed

1 file changed

+69
-10
lines changed

pysteps/xarray_helpers.py

Lines changed: 69 additions & 10 deletions
Original file line numberDiff line numberDiff line change
@@ -83,6 +83,7 @@ def convert_input_to_xarray_dataset(
8383
precip: np.ndarray,
8484
quality: np.ndarray | None,
8585
metadata: dict[str, str | float | None],
86+
startdate: datetime | None = None,
8687
) -> xr.Dataset:
8788
"""
8889
Read a precip, quality, metadata tuple as returned by the importers
@@ -99,6 +100,8 @@ def convert_input_to_xarray_dataset(
99100
metadata: dict
100101
Metadata dictionary containing the attributes described in the
101102
documentation of :py:mod:`pysteps.io.importers`.
103+
startdate: datetime, None
104+
Datetime object containing the start date and time for the nowcast
102105
103106
Returns
104107
-------
@@ -107,7 +110,31 @@ def convert_input_to_xarray_dataset(
107110
108111
"""
109112
var_name, attrs = cf_parameters_from_unit(metadata["unit"])
110-
h, w = precip.shape
113+
114+
dims = None
115+
timesteps = None
116+
ens_number = None
117+
118+
if precip.ndim == 4:
119+
ens_number, timesteps, h, w = precip.shape
120+
dims = ["ens_number", "time", "y", "x"]
121+
122+
if startdate is None:
123+
raise Exception("startdate missing")
124+
125+
elif precip.ndim == 3:
126+
timesteps, h, w = precip.shape
127+
dims = ["time", "y", "x"]
128+
129+
if startdate is None:
130+
raise Exception("startdate missing")
131+
132+
elif precip.ndim == 2:
133+
h, w = precip.shape
134+
dims = ["y", "x"]
135+
else:
136+
raise Exception(f"Precip field shape: {precip.shape} not supported")
137+
111138
x_r = np.linspace(metadata["x1"], metadata["x2"], w + 1)[:-1]
112139
x_r += 0.5 * (x_r[1] - x_r[0])
113140
y_r = np.linspace(metadata["y1"], metadata["y2"], h + 1)[:-1]
@@ -142,25 +169,37 @@ def convert_input_to_xarray_dataset(
142169

143170
data_vars = {
144171
var_name: (
145-
["y", "x"],
172+
dims,
146173
precip,
147174
{
148175
"units": attrs["units"],
149176
"standard_name": attrs["standard_name"],
150177
"long_name": attrs["long_name"],
151178
"grid_mapping": "projection",
152-
"transform": metadata["transform"],
153-
"accutime": metadata["accutime"],
154-
"threshold": metadata["threshold"],
155-
"zerovalue": metadata["zerovalue"],
156-
"zr_a": metadata["zr_a"],
157-
"zr_b": metadata["zr_b"],
179+
"xpixelsize": xpixelsize,
180+
"ypixelsize": ypixelsize,
158181
},
159182
)
160183
}
184+
185+
metadata_keys = [
186+
"transform",
187+
"accutime",
188+
"threshold",
189+
"zerovalue",
190+
"zr_a",
191+
"zr_b",
192+
"cartesian_unit",
193+
"yorigin",
194+
]
195+
196+
for metadata_field in metadata_keys:
197+
if metadata_field in metadata:
198+
data_vars[var_name][2][metadata_field] = metadata[metadata_field]
199+
161200
if quality is not None:
162201
data_vars["quality"] = (
163-
["y", "x"],
202+
dims,
164203
quality,
165204
{
166205
"units": "1",
@@ -210,6 +249,26 @@ def convert_input_to_xarray_dataset(
210249
},
211250
),
212251
}
252+
253+
if ens_number is not None:
254+
coords["ens_number"] = (
255+
["ens_number"],
256+
list(range(1, ens_number + 1, 1)),
257+
{
258+
"long_name": "ensemble member",
259+
"standard_name": "realization",
260+
"units": "",
261+
},
262+
)
263+
264+
if timesteps is not None:
265+
startdate_str = datetime.strftime(startdate, "%Y-%m-%d %H:%M:%S")
266+
267+
coords["time"] = (
268+
["time"],
269+
list(range(1, timesteps + 1, 1)),
270+
{"long_name": "forecast time", "units": "seconds since %s" % startdate_str},
271+
)
213272
if grid_mapping_var_name is not None:
214273
coords[grid_mapping_name] = (
215274
[],
@@ -223,7 +282,7 @@ def convert_input_to_xarray_dataset(
223282
"precip_var": var_name,
224283
}
225284
dataset = xr.Dataset(data_vars=data_vars, coords=coords, attrs=attrs)
226-
return dataset.sortby(["y", "x"])
285+
return dataset.sortby(dims)
227286

228287

229288
def convert_output_to_xarray_dataset(

0 commit comments

Comments
 (0)