Flag Variables

cf_xarray has some support for flag variables, including flag masks.

Flag Values

import cf_xarray
import xarray as xr

da = xr.DataArray(
    [1, 2, 1, 1, 2, 2, 3, 3, 3, 3],
    dims=("time",),
    attrs={
        "flag_values": [1, 2, 3],
        "flag_meanings": "atlantic_ocean pacific_ocean indian_ocean",
        "standard_name": "region",
    },
  name="region",
)
da.cf
Flag Variable:
       Flag Meanings:   atlantic_ocean:     1 
                         pacific_ocean:     2 
                          indian_ocean:     3 

Coordinates:
             CF Axes:   X, Y, Z, T: n/a

      CF Coordinates:   longitude, latitude, vertical, time: n/a

       Cell Measures:   area, volume: n/a

      Standard Names:   n/a

              Bounds:   n/a

       Grid Mappings:   n/a

Now you can perform meaningful boolean comparisons that take advantage of the flag_meanings attribute:

# compare to da == 1
da.cf == "atlantic_ocean"
<xarray.DataArray 'region' (time: 10)> Size: 10B
array([ True, False,  True,  True, False, False, False, False, False,
       False])
Dimensions without coordinates: time

Similarly with membership tests using DataArray.cf.isin()

# compare to da.isin([2, 3])
da.cf.isin(["indian_ocean", "pacific_ocean"])
<xarray.DataArray 'region' (time: 10)> Size: 10B
array([False,  True, False, False,  True,  True,  True,  True,  True,
        True])
Dimensions without coordinates: time

You can also check whether a DataArray has the appropriate attributes to be recognized as a flag variable using DataArray.cf.is_flag_variable()

da.cf.is_flag_variable
True

Flag Masks

Warning

Interpreting flag masks is very lightly tested. Please double-check the results and open an issue or pull request to suggest improvements.

Load an example dataset:

from cf_xarray.datasets import flag_indep

flag_indep
<xarray.DataArray 'flag_var' (time: 8)> Size: 8B
array([0, 1, 2, 3, 4, 5, 6, 7], dtype=uint8)
Dimensions without coordinates: time
Attributes:
    flag_masks:     [1, 2, 4]
    flag_meanings:  flag_1 flag_2 flag_4
    standard_name:  flag_independent
flag_indep.cf
Flag Variable:
       Flag Meanings:   flag_1:        1  / Bit: .......1
                        flag_2:        2  / Bit: ......1.
                        flag_4:        4  / Bit: .....1..

Coordinates:
             CF Axes:   X, Y, Z, T: n/a

      CF Coordinates:   longitude, latitude, vertical, time: n/a

       Cell Measures:   area, volume: n/a

      Standard Names:   n/a

              Bounds:   n/a

       Grid Mappings:   n/a
flag_indep.cf == "flag_1"
<xarray.DataArray 'flag_var' (time: 8)> Size: 8B
array([False,  True, False,  True, False,  True, False,  True])
Dimensions without coordinates: time

And isin

flag_indep.cf.isin(["flag_1", "flag_3"])
<xarray.DataArray 'flag_var' (time: 8)> Size: 8B
array([False,  True, False,  True, False,  True, False,  True])
Dimensions without coordinates: time

Combined masks and values

Warning

Interpreting a mix of flag masks and flag values is very lightly tested. Please double-check the results and open an issue or pull request to suggest improvements.

Load an example dataset:

from cf_xarray.datasets import flag_mix

flag_mix
<xarray.DataArray 'flag_var' (time: 8)> Size: 8B
array([ 4,  8, 13,  5, 10, 14,  7,  3], dtype=uint8)
Dimensions without coordinates: time
Attributes:
    flag_values:    [1, 2, 4, 8, 12]
    flag_masks:     [1, 2, 12, 12, 12]
    flag_meanings:  flag_1 flag_2 flag_3 flag_4 flag_5
    standard_name:  flag_mix
flag_mix.cf
Flag Variable:
       Flag Meanings:   flag_1:        1  / Bit: .......1
                        flag_2:        2  / Bit: ......1.
                        flag_3:       12 & 4  / Bit: ....01..
                        flag_4:       12 & 8  / Bit: ....10..
                        flag_5:       12 & 12  / Bit: ....11..

Coordinates:
             CF Axes:   X, Y, Z, T: n/a

      CF Coordinates:   longitude, latitude, vertical, time: n/a

       Cell Measures:   area, volume: n/a

      Standard Names:   n/a

              Bounds:   n/a

       Grid Mappings:   n/a
flag_mix.cf == 'flag_4'
<xarray.DataArray 'flag_var' (time: 8)> Size: 8B
array([False,  True, False, False,  True, False, False, False])
Dimensions without coordinates: time