Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

new vertical coord decode and reordering #315

Merged
merged 22 commits into from
Mar 29, 2022
Merged
Show file tree
Hide file tree
Changes from 15 commits
Commits
Show all changes
22 commits
Select commit Hold shift + click to select a range
7ff8814
new vertical coord decode and reordering
kthyng Mar 25, 2022
3b5197d
[pre-commit.ci] auto fixes from pre-commit.com hooks
pre-commit-ci[bot] Mar 26, 2022
57299a9
Update doc/whats-new.rst
kthyng Mar 27, 2022
02f1250
Update cf_xarray/accessor.py
kthyng Mar 27, 2022
8fafedc
preferentially use outnames over prefix
kthyng Mar 27, 2022
d5da68e
[pre-commit.ci] auto fixes from pre-commit.com hooks
pre-commit-ci[bot] Mar 27, 2022
429349d
added in transpose statement that does not work yet
kthyng Mar 27, 2022
1720fe3
Merge branch 'add_ocean_sigma_coordinate' of github.com:kthyng/cf-xar…
kthyng Mar 27, 2022
9c461c5
[pre-commit.ci] auto fixes from pre-commit.com hooks
pre-commit-ci[bot] Mar 27, 2022
4ede9a5
removed reordering stuff that did not work
kthyng Mar 28, 2022
df216bb
merged
kthyng Mar 28, 2022
8f69bc1
provide a default name now
kthyng Mar 28, 2022
5bf1582
I think the tests are back how they were
kthyng Mar 28, 2022
6af2ebd
trying to follow dcherian suggestions
kthyng Mar 28, 2022
41c8acf
small precommit change
kthyng Mar 28, 2022
c89d1f0
Update cf_xarray/tests/test_accessor.py
kthyng Mar 29, 2022
6043932
small fix from review
kthyng Mar 29, 2022
11f8029
[pre-commit.ci] auto fixes from pre-commit.com hooks
pre-commit-ci[bot] Mar 29, 2022
02a27a0
[pre-commit.ci] auto fixes from pre-commit.com hooks
pre-commit-ci[bot] Mar 29, 2022
b100452
Update cf_xarray/accessor.py
dcherian Mar 29, 2022
64380b9
Update cf_xarray/tests/test_accessor.py
dcherian Mar 29, 2022
6fbe702
fixed tests I think
kthyng Mar 29, 2022
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
41 changes: 35 additions & 6 deletions cf_xarray/accessor.py
Original file line number Diff line number Diff line change
Expand Up @@ -2190,12 +2190,15 @@ def bounds_to_vertices(
)
return obj

def decode_vertical_coords(self, prefix="z"):
def decode_vertical_coords(self, *, outnames=None, prefix=None):
"""
Decode parameterized vertical coordinates in place.

Parameters
----------
outnames : dict, optional
Keys of outnames are the input sigma/s coordinate variable name and
the values are the name to use for the associated vertical coordinate.
prefix : str, optional
Prefix for newly created z variables.
E.g. ``s_rho`` becomes ``z_rho``
Expand All @@ -2209,7 +2212,8 @@ def decode_vertical_coords(self, prefix="z"):
Will only decode when the ``formula_terms`` and ``standard_name`` attributes
are set on the parameter (e.g ``s_rho`` )

Currently only supports ``ocean_s_coordinate_g1`` and ``ocean_s_coordinate_g2``.
Currently only supports ``ocean_s_coordinate_g1``, ``ocean_s_coordinate_g2``,
and ``ocean_sigma_coordinate``.

.. warning::
Very lightly tested. Please double check the results.
Expand All @@ -2223,12 +2227,29 @@ def decode_vertical_coords(self, prefix="z"):
requirements = {
"ocean_s_coordinate_g1": {"depth_c", "depth", "s", "C", "eta"},
"ocean_s_coordinate_g2": {"depth_c", "depth", "s", "C", "eta"},
"ocean_sigma_coordinate": {"sigma", "eta", "depth"},
}

allterms = self.formula_terms
for dim in allterms:
suffix = dim.split("_")
zname = f"{prefix}_" + "_".join(suffix[1:])
if prefix is None:
# set outnames here
# need a default name option
if outnames is None:
zname = f"z_{dim}"
else:
kthyng marked this conversation as resolved.
Show resolved Hide resolved
try:
zname = outnames[dim]
except KeyError:
print("Your `outnames` need to include a key of `dim`.")

else:
warnings.warn(
dcherian marked this conversation as resolved.
Show resolved Hide resolved
"`prefix` is being deprecated; use `outnames` instead.",
DeprecationWarning,
)
suffix = dim.split("_")
zname = f"{prefix}_" + "_".join(suffix[1:])

if "standard_name" not in ds[dim].attrs:
continue
Expand All @@ -2254,23 +2275,31 @@ def decode_vertical_coords(self, prefix="z"):
terms["depth_c"] * terms["s"]
+ (terms["depth"] - terms["depth_c"]) * terms["C"]
)

# z(n,k,j,i) = S(k,j,i) + eta(n,j,i) * (1 + S(k,j,i) / depth(j,i))
ds.coords[zname] = S + terms["eta"] * (1 + S / terms["depth"])
ztemp = S + terms["eta"] * (1 + S / terms["depth"])

elif stdname == "ocean_s_coordinate_g2":
# make sure all necessary terms are present in terms
# (depth_c * s(k) + depth(j,i) * C(k)) / (depth_c + depth(j,i))
S = (terms["depth_c"] * terms["s"] + terms["depth"] * terms["C"]) / (
terms["depth_c"] + terms["depth"]
)

# z(n,k,j,i) = eta(n,j,i) + (eta(n,j,i) + depth(j,i)) * S(k,j,i)
ds.coords[zname] = terms["eta"] + (terms["eta"] + terms["depth"]) * S
ztemp = terms["eta"] + (terms["eta"] + terms["depth"]) * S

elif stdname == "ocean_sigma_coordinate":
dcherian marked this conversation as resolved.
Show resolved Hide resolved
# z(n,k,j,i) = eta(n,j,i) + sigma(k)*(depth(j,i)+eta(n,j,i))
ztemp = terms["eta"] + terms["sigma"] * (terms["depth"] + terms["eta"])

else:
raise NotImplementedError(
f"Coordinate function for {stdname!r} not implemented yet. Contributions welcome!"
)

ds.coords[zname] = ztemp


@xr.register_dataarray_accessor("cf")
class CFDataArrayAccessor(CFAccessor):
Expand Down
23 changes: 23 additions & 0 deletions cf_xarray/datasets.py
Original file line number Diff line number Diff line change
Expand Up @@ -16,6 +16,29 @@
ds_no_attrs[variable].attrs = {}


# POM dataset
pomds = xr.Dataset()
pomds["sigma"] = (
# fmt: off
"sigma",
[-0.983333, -0.95 , -0.916667, -0.883333, -0.85 , -0.816667,
-0.783333, -0.75 , -0.716667, -0.683333, -0.65 , -0.616667,
-0.583333, -0.55 , -0.516667, -0.483333, -0.45 , -0.416667,
-0.383333, -0.35 , -0.316667, -0.283333, -0.25 , -0.216667,
-0.183333, -0.15 , -0.116667, -0.083333, -0.05 , -0.016667],
# fmt: on
{
"units": "sigma_level",
"long_name": "Sigma Stretched Vertical Coordinate at Nodes",
"positive": "down",
"standard_name": "ocean_sigma_coordinate",
"formula_terms": "sigma: sigma eta: zeta depth: depth",
}
)
pomds["depth"] = 175.0
pomds["zeta"] = ("ocean_time", [-0.155356, -0.127435])


popds = xr.Dataset()
popds.coords["TLONG"] = (
("nlat", "nlon"),
Expand Down
17 changes: 15 additions & 2 deletions cf_xarray/tests/test_accessor.py
Original file line number Diff line number Diff line change
Expand Up @@ -24,6 +24,7 @@
forecast,
mollwds,
multiple,
pomds,
popds,
romsds,
vert,
Expand Down Expand Up @@ -1029,15 +1030,16 @@ def test_param_vcoord_ocean_s_coord():
romsds.hc + romsds.h
)
expected = romsds.zeta + (romsds.zeta + romsds.h) * Zo_rho
romsds.cf.decode_vertical_coords()
romsds.cf.decode_vertical_coords(outnames={"s_rho": "z_rho"})
assert_allclose(
romsds.z_rho.reset_coords(drop=True), expected.reset_coords(drop=True)
)

romsds.s_rho.attrs["standard_name"] = "ocean_s_coordinate_g1"
Zo_rho = romsds.hc * (romsds.s_rho - romsds.Cs_r) + romsds.Cs_r * romsds.h

expected = Zo_rho + romsds.zeta * (1 + Zo_rho / romsds.h)
romsds.cf.decode_vertical_coords()
romsds.cf.decode_vertical_coords(outnames={"s_rho": "z_rho"})
assert_allclose(
romsds.z_rho.reset_coords(drop=True), expected.reset_coords(drop=True)
)
Expand All @@ -1056,6 +1058,17 @@ def test_param_vcoord_ocean_s_coord():
copy.cf.decode_vertical_coords()


def test_param_vcoord_ocean_sigma_coordinate():
expected = pomds.zeta + pomds.sigma * (pomds.depth + pomds.zeta)
pomds.cf.decode_vertical_coords(outnames={"sigma": "z"})
assert_allclose(pomds.z.reset_coords(drop=True), expected.reset_coords(drop=True))

copy = pomds.copy(deep=True)
del copy["zeta"]
with pytest.raises(KeyError):
dcherian marked this conversation as resolved.
Show resolved Hide resolved
copy.cf.decode_vertical_coords()
kthyng marked this conversation as resolved.
Show resolved Hide resolved


def test_formula_terms():
srhoterms = {
"s": "s_rho",
Expand Down
4 changes: 2 additions & 2 deletions doc/parametricz.md
Original file line number Diff line number Diff line change
Expand Up @@ -21,7 +21,7 @@ xr.set_options(display_expand_data=False)

# Parametric Vertical Coordinates

`cf_xarray` supports decoding [parametric vertical coordinates](http://cfconventions.org/Data/cf-conventions/cf-conventions-1.8/cf-conventions.html#parametric-vertical-coordinate) encoded in the `formula_terms` attribute using {py:meth}`Dataset.cf.decode_vertical_coords`. Right now, only the two ocean s-coordinates are supported, but support for the [rest](http://cfconventions.org/Data/cf-conventions/cf-conventions-1.8/cf-conventions.html#parametric-v-coord) should be easy to add (Pull Requests are very welcome!).
`cf_xarray` supports decoding [parametric vertical coordinates](http://cfconventions.org/Data/cf-conventions/cf-conventions-1.8/cf-conventions.html#parametric-vertical-coordinate) encoded in the `formula_terms` attribute using {py:meth}`Dataset.cf.decode_vertical_coords`. Right now, only the two ocean s-coordinates and `ocean_sigma_coordinate` are supported, but support for the [rest](http://cfconventions.org/Data/cf-conventions/cf-conventions-1.8/cf-conventions.html#parametric-v-coord) should be easy to add (Pull Requests are very welcome!).

## Decoding parametric coordinates
```{code-cell}
Expand All @@ -33,7 +33,7 @@ romsds
Now we decode the vertical coordinates **in-place**. Note the new `z_rho` variable. `cf_xarray` sees that `s_rho` has a `formula_terms` attribute, looks up the right formula using `s_rho.attrs["standard_name"]` and computes a new vertical coordinate variable.

```{code-cell}
romsds.cf.decode_vertical_coords() # adds new z_rho variable
romsds.cf.decode_vertical_coords(outnames={'s_rho': 'z_rho'}) # adds new z_rho variable
romsds.z_rho
```

Expand Down
6 changes: 6 additions & 0 deletions doc/whats-new.rst
Original file line number Diff line number Diff line change
Expand Up @@ -2,6 +2,12 @@

What's New
----------

v 0.7.1 (unreleased)
====================
- added another type of vertical coordinate to decode: ``ocean_sigma_coordinate``. By `Kristen Thyng`_.


v0.7.0 (January 24, 2022)
=========================
- Many improvements to autoguessing for plotting. By `Deepak Cherian`_
Expand Down