Skip to content

Commit ca0f2a0

Browse files
committed
Add .cf.formula_terms
Closes #211
1 parent 64e8695 commit ca0f2a0

File tree

4 files changed

+53
-14
lines changed

4 files changed

+53
-14
lines changed

cf_xarray/accessor.py

Lines changed: 33 additions & 14 deletions
Original file line numberDiff line numberDiff line change
@@ -1,6 +1,7 @@
11
import functools
22
import inspect
33
import itertools
4+
import re
45
import warnings
56
from collections import ChainMap
67
from typing import (
@@ -1480,6 +1481,16 @@ def __getitem__(self, key: Union[str, List[str]]) -> Union[DataArray, Dataset]:
14801481
"""
14811482
return _getitem(self, key)
14821483

1484+
@property
1485+
def formula_terms(self) -> Dict[str, Dict[str, str]]:
1486+
"""
1487+
Property that returns a dictionary
1488+
{parametric_coord_name: {standard_term_name: variable_name}}
1489+
"""
1490+
return {
1491+
dim: self._obj[dim].cf.formula_terms for dim in _get_dims(self._obj, "Z")
1492+
}
1493+
14831494
def get_bounds(self, key: str) -> DataArray:
14841495
"""
14851496
Get bounds variable corresponding to key.
@@ -1668,36 +1679,29 @@ def decode_vertical_coords(self, prefix="z"):
16681679
.. warning::
16691680
Very lightly tested. Please double check the results.
16701681
"""
1671-
import re
1672-
16731682
ds = self._obj
1674-
dims = _get_dims(ds, "Z")
16751683

16761684
requirements = {
16771685
"ocean_s_coordinate_g1": {"depth_c", "depth", "s", "C", "eta"},
16781686
"ocean_s_coordinate_g2": {"depth_c", "depth", "s", "C", "eta"},
16791687
}
16801688

1681-
for dim in dims:
1689+
allterms = self.formula_terms
1690+
for dim in allterms:
16821691
suffix = dim.split("_")
16831692
zname = f"{prefix}_" + "_".join(suffix[1:])
16841693

1685-
if (
1686-
"formula_terms" not in ds[dim].attrs
1687-
or "standard_name" not in ds[dim].attrs
1688-
):
1694+
if "standard_name" not in ds[dim].attrs:
16891695
continue
1690-
1691-
formula_terms = ds[dim].attrs["formula_terms"]
16921696
stdname = ds[dim].attrs["standard_name"]
16931697

16941698
# map "standard" formula term names to actual variable names
16951699
terms = {}
1696-
for mapping in re.sub(": ", ":", formula_terms).split(" "):
1697-
key, value = mapping.split(":")
1700+
for key, value in allterms[dim].items():
16981701
if value not in ds:
16991702
raise KeyError(
1700-
f"Variable {value!r} is required to decode coordinate for {dim} but it is absent in the Dataset."
1703+
f"Variable {value!r} is required to decode coordinate for {dim!r}"
1704+
" but it is absent in the Dataset."
17011705
)
17021706
terms[key] = ds[value]
17031707

@@ -1725,12 +1729,27 @@ def decode_vertical_coords(self, prefix="z"):
17251729

17261730
else:
17271731
raise NotImplementedError(
1728-
f"Coordinate function for {stdname} not implemented yet. Contributions welcome!"
1732+
f"Coordinate function for {stdname!r} not implemented yet. Contributions welcome!"
17291733
)
17301734

17311735

17321736
@xr.register_dataarray_accessor("cf")
17331737
class CFDataArrayAccessor(CFAccessor):
1738+
@property
1739+
def formula_terms(self) -> Dict[str, str]:
1740+
"""
1741+
Property that returns a dictionary
1742+
{parametric_coord_name: {standard_term_name: variable_name}}
1743+
"""
1744+
da = self._obj
1745+
dims = _single(_get_dims)(da, "Z")[0]
1746+
terms = {}
1747+
formula_terms = da[dims].attrs.get("formula_terms", "")
1748+
for mapping in re.sub(r"\s*:\s*", ":", formula_terms).split():
1749+
key, value = mapping.split(":")
1750+
terms[key] = value
1751+
return terms
1752+
17341753
def __getitem__(self, key: Union[str, List[str]]) -> DataArray:
17351754
"""
17361755
Index into a DataArray making use of CF attributes.

cf_xarray/tests/test_accessor.py

Lines changed: 16 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -915,6 +915,22 @@ def test_param_vcoord_ocean_s_coord():
915915
copy.cf.decode_vertical_coords()
916916

917917

918+
def test_formula_terms():
919+
srhoterms = {
920+
"s": "s_rho",
921+
"C": "Cs_r",
922+
"eta": "zeta",
923+
"depth": "h",
924+
"depth_c": "hc",
925+
}
926+
assert romsds.cf.formula_terms == {"s_rho": srhoterms}
927+
assert romsds["temp"].cf.formula_terms == srhoterms
928+
assert romsds["s_rho"].cf.formula_terms == srhoterms
929+
with pytest.raises(KeyError):
930+
# x,y,t variable
931+
romsds["zeta"].cf.formula_terms
932+
933+
918934
def test_standard_name_mapper():
919935
da = xr.DataArray(
920936
np.arange(6),

doc/api.rst

Lines changed: 2 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -31,6 +31,7 @@ Attributes
3131
DataArray.cf.cell_measures
3232
DataArray.cf.coordinates
3333
DataArray.cf.standard_names
34+
DataArray.cf.formula_terms
3435
DataArray.cf.plot
3536

3637

@@ -65,6 +66,7 @@ Attributes
6566
Dataset.cf.axes
6667
Dataset.cf.cell_measures
6768
Dataset.cf.coordinates
69+
Dataset.cf.formula_terms
6870
Dataset.cf.standard_names
6971

7072
.. _dsmeth:

doc/whats-new.rst

Lines changed: 2 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -6,6 +6,8 @@ What's New
66
v0.5.2 (unreleased)
77
===================
88

9+
- Added :py:attr:`DataArray.cf.formula_terms` and :py:attr:`Dataset.cf.formula_terms`.
10+
By `Deepak Cherian`_.
911
- :py:meth:`DataArray.cf.differentiate` and :py:meth:`Dataset.cf.differentiate` can optionally correct
1012
sign of the derivative by interpreting the ``"positive"`` attribute. By `Deepak Cherian`_.
1113

0 commit comments

Comments
 (0)