Region masks#
POP includes a default region mask as a component of the grid information. This is often not super
relevant for analyses. pop_tools
provides several alternative region masks; these are demostrated here.
Import packages#
%matplotlib inline
import matplotlib.colors as colors
import matplotlib.pyplot as plt
import numpy as np
import xarray as xr
import pop_tools
Load POP grid as xarray.Dataset
#
grid_name = 'POP_gx1v7'
ds = pop_tools.get_grid(grid_name)
ds
Downloading file 'inputdata/ocn/pop/gx1v7/grid/horiz_grid_20010402.ieeer8' from 'https://svn-ccsm-inputdata.cgd.ucar.edu/trunk/inputdata/ocn/pop/gx1v7/grid/horiz_grid_20010402.ieeer8' to '/home/docs/.pop_tools'.
/home/docs/checkouts/readthedocs.org/user_builds/pop-tools/conda/latest/lib/python3.9/site-packages/urllib3/connectionpool.py:1095: InsecureRequestWarning: Unverified HTTPS request is being made to host 'svn-ccsm-inputdata.cgd.ucar.edu'. Adding certificate verification is strongly advised. See: https://urllib3.readthedocs.io/en/latest/advanced-usage.html#tls-warnings
warnings.warn(
0%| | 0.00/6.88M [00:00<?, ?B/s]
1%|β | 66.6k/6.88M [00:00<00:12, 563kB/s]
9%|ββββ | 607k/6.88M [00:00<00:02, 2.89MB/s]
55%|βββββββββββββββββββββ | 3.77M/6.88M [00:00<00:00, 14.8MB/s]
100%|βββββββββββββββββββββββββββββββββββββ| 6.87M/6.88M [00:00<00:00, 20.7MB/s]
0%| | 0.00/6.88M [00:00<?, ?B/s]
100%|βββββββββββββββββββββββββββββββββββββ| 6.88M/6.88M [00:00<00:00, 2.96GB/s]
Downloading file 'inputdata/ocn/pop/gx1v7/grid/topography_20161215.ieeei4' from 'https://svn-ccsm-inputdata.cgd.ucar.edu/trunk/inputdata/ocn/pop/gx1v7/grid/topography_20161215.ieeei4' to '/home/docs/.pop_tools'.
/home/docs/checkouts/readthedocs.org/user_builds/pop-tools/conda/latest/lib/python3.9/site-packages/urllib3/connectionpool.py:1095: InsecureRequestWarning: Unverified HTTPS request is being made to host 'svn-ccsm-inputdata.cgd.ucar.edu'. Adding certificate verification is strongly advised. See: https://urllib3.readthedocs.io/en/latest/advanced-usage.html#tls-warnings
warnings.warn(
0%| | 0.00/492k [00:00<?, ?B/s]
17%|βββββββ | 82.9k/492k [00:00<00:00, 700kB/s]
0%| | 0.00/492k [00:00<?, ?B/s]
100%|ββββββββββββββββββββββββββββββββββββββββ| 492k/492k [00:00<00:00, 617MB/s]
Downloading file 'inputdata/ocn/pop/gx1v7/grid/region_mask_20151008.ieeei4' from 'https://svn-ccsm-inputdata.cgd.ucar.edu/trunk/inputdata/ocn/pop/gx1v7/grid/region_mask_20151008.ieeei4' to '/home/docs/.pop_tools'.
/home/docs/checkouts/readthedocs.org/user_builds/pop-tools/conda/latest/lib/python3.9/site-packages/urllib3/connectionpool.py:1095: InsecureRequestWarning: Unverified HTTPS request is being made to host 'svn-ccsm-inputdata.cgd.ucar.edu'. Adding certificate verification is strongly advised. See: https://urllib3.readthedocs.io/en/latest/advanced-usage.html#tls-warnings
warnings.warn(
0%| | 0.00/492k [00:00<?, ?B/s]
20%|ββββββββ | 99.3k/492k [00:00<00:00, 831kB/s]
0%| | 0.00/492k [00:00<?, ?B/s]
100%|ββββββββββββββββββββββββββββββββββββββββ| 492k/492k [00:00<00:00, 610MB/s]
<xarray.Dataset> Dimensions: (nlat: 384, nlon: 320, z_t: 60, z_w: 60, z_w_bot: 60, nreg: 13) Coordinates: * z_t (z_t) float64 500.0 1.5e+03 2.5e+03 ... 5.125e+05 5.375e+05 * z_w (z_w) float64 0.0 1e+03 2e+03 3e+03 ... 4.75e+05 5e+05 5.25e+05 * z_w_bot (z_w_bot) float64 1e+03 2e+03 3e+03 ... 5e+05 5.25e+05 5.5e+05 * nreg (nreg) int64 0 1 2 3 4 5 6 7 8 9 10 11 12 Dimensions without coordinates: nlat, nlon Data variables: (12/15) TLAT (nlat, nlon) float64 -79.22 -79.22 -79.22 ... 72.2 72.19 72.19 TLONG (nlat, nlon) float64 320.6 321.7 322.8 ... 318.9 319.4 319.8 ULAT (nlat, nlon) float64 -78.95 -78.95 -78.95 ... 72.42 72.41 72.41 ULONG (nlat, nlon) float64 321.1 322.3 323.4 ... 319.2 319.6 320.0 DXT (nlat, nlon) float64 1.894e+06 1.893e+06 ... 1.473e+06 DYT (nlat, nlon) float64 5.94e+06 5.94e+06 ... 5.046e+06 5.046e+06 ... ... UAREA (nlat, nlon) float64 1.423e+13 1.423e+13 ... 7.639e+12 KMT (nlat, nlon) int32 0 0 0 0 0 0 0 0 0 0 ... 0 0 0 0 0 0 0 0 0 0 REGION_MASK (nlat, nlon) int32 0 0 0 0 0 0 0 0 0 0 ... 0 0 0 0 0 0 0 0 0 0 dz (z_t) float64 1e+03 1e+03 1e+03 ... 2.499e+04 2.5e+04 2.5e+04 region_name (nreg) <U21 'Black Sea' 'Baltic Sea' ... 'Hudson Bay' region_val (nreg) int64 -13 -12 -5 1 2 3 4 6 7 8 9 10 11 Attributes: lateral_dims: [384, 320] vertical_dims: 60 vert_grid_file: gx1v7_vert_grid horiz_grid_fname: inputdata/ocn/pop/gx1v7/grid/horiz_grid_20010402.ieeer8 topography_fname: inputdata/ocn/pop/gx1v7/grid/topography_20161215.ieeei4 region_mask_fname: inputdata/ocn/pop/gx1v7/grid/region_mask_20151008.ieeei4 type: dipole title: POP_gx1v7 grid
Plot default REGION_MASK
#
The default REGION_MASK
is a 2-D array with unique integer values for each region. Negative integers denote βmarginal seas,β which are not directly connected to the ocean.
regions = np.array(np.unique(ds.REGION_MASK))
regions
array([-13, -12, -5, 0, 1, 2, 3, 4, 6, 7, 8, 9, 10,
11], dtype=int32)
ds.REGION_MASK.plot.contourf(levels=regions, cmap='tab20');

More useful region masks#
Itβs often more useful to define a region mask as a 3-D array of zeros and ones, where the first dimension is region
; this permits overlapping regions and is convenient for computation because the mask can be applied by multiplication, which yields a region
dimension via broadcasting.
pop_tools
supports converting the default REGION_MASK
to this type of mask thru the region_mask_3d
function.
mask3d = pop_tools.region_mask_3d(grid_name, mask_name='default')
mask3d
<xarray.DataArray (region: 13, nlat: 384, nlon: 320)> array([[[0., 0., 0., ..., 0., 0., 0.], [0., 0., 0., ..., 0., 0., 0.], [0., 0., 0., ..., 0., 0., 0.], ..., [0., 0., 0., ..., 0., 0., 0.], [0., 0., 0., ..., 0., 0., 0.], [0., 0., 0., ..., 0., 0., 0.]], [[0., 0., 0., ..., 0., 0., 0.], [0., 0., 0., ..., 0., 0., 0.], [0., 0., 0., ..., 0., 0., 0.], ..., [0., 0., 0., ..., 0., 0., 0.], [0., 0., 0., ..., 0., 0., 0.], [0., 0., 0., ..., 0., 0., 0.]], [[0., 0., 0., ..., 0., 0., 0.], [0., 0., 0., ..., 0., 0., 0.], [0., 0., 0., ..., 0., 0., 0.], ..., ... ..., [0., 0., 0., ..., 0., 0., 0.], [0., 0., 0., ..., 0., 0., 0.], [0., 0., 0., ..., 0., 0., 0.]], [[0., 0., 0., ..., 0., 0., 0.], [0., 0., 0., ..., 0., 0., 0.], [0., 0., 0., ..., 0., 0., 0.], ..., [0., 0., 0., ..., 0., 0., 0.], [0., 0., 0., ..., 0., 0., 0.], [0., 0., 0., ..., 0., 0., 0.]], [[0., 0., 0., ..., 0., 0., 0.], [0., 0., 0., ..., 0., 0., 0.], [0., 0., 0., ..., 0., 0., 0.], ..., [0., 0., 0., ..., 0., 0., 0.], [0., 0., 0., ..., 0., 0., 0.], [0., 0., 0., ..., 0., 0., 0.]]]) Coordinates: * region (region) <U21 'Black Sea' 'Baltic Sea' ... 'Hudson Bay' Dimensions without coordinates: nlat, nlon Attributes: mask_name: default
This mask includes a coordinate variable with the region names.
mask3d.region
<xarray.DataArray 'region' (region: 13)> array(['Black Sea', 'Baltic Sea', 'Red Sea', 'Southern Ocean', 'Pacific Ocean', 'Indian Ocean', 'Persian Gulf', 'Atlantic Ocean', 'Mediterranean Sea', 'Lab. Sea & Baffin Bay', 'GIN Seas', 'Arctic Ocean', 'Hudson Bay'], dtype='<U21') Coordinates: * region (region) <U21 'Black Sea' 'Baltic Sea' ... 'Hudson Bay'
A particular region can be selected by name.
mask3d.sel(region='Southern Ocean').plot();

To visualize all the regions, we can define a help plotting function,
def visualize_mask(mask3d):
nregion = len(mask3d.region)
# mask out land
mask3d = mask3d.where(ds.KMT > 0)
# visualize the regions
ncol = int(np.sqrt(nregion))
nrow = int(nregion / ncol) + min(1, nregion % ncol)
fig, ax = plt.subplots(nrow, ncol, figsize=(4 * ncol, 3 * nrow), constrained_layout=True)
for i, region in enumerate(mask3d.region.values):
plt.axes(ax.ravel()[i])
mask3d.sel(region=region).plot()
# delete the unused axes
for i in range(nregion, ncol * nrow):
fig.delaxes(ax.ravel()[i])
fig.suptitle(f'Mask name = {mask3d.mask_name}', fontsize=16)
and apply it to the default mask created above.
visualize_mask(mask3d)

Alternative region masks#
Other useful region masks are pre-defined in the package. list_region_masks
returns a list of pre-defined masks.
region_masks = pop_tools.list_region_masks(grid_name)
region_masks
['lat-range-basin', 'Pacific-Indian-Atlantic']
We can visualize all of these using the helper function above.
for region_mask in region_masks:
mask3d = pop_tools.region_mask_3d(grid_name, mask_name=region_mask)
visualize_mask(mask3d)


To illustrated how regions cover the global domain, including with overlap, we can sum over the region
dimension.
mask3d = pop_tools.region_mask_3d(grid_name, mask_name='lat-range-basin')
mask3d.sum('region').plot();

User defined region masks#
Finally, it is also possible to make a region mask on the fly by building a dictionary containing the defining logic. region_mask_3d
accepts a region_defs
argument. This is a dictionary of the following form.
region_defs = {region1_name: list_of_criteria_dicts_1,
region2_name: list_of_criteria_dicts_2,...}
The list_of_criteria_dicts
are lists of dictionaries; each must include the keys βmatchβ or βboundsβ. For instance:
list_of_criteria_dicts_1 = [{'match': {'REGION_MASK': [1, 2, 3, 6]},
'bounds': {'TLAT': [-90., -30.]}}]
will return a mask where the default REGION_MASK
matches the specified values and TLAT
falls between the specified bounds. Multiple entries in the list_of_criteria_dicts
are applied with an βorβ condition.
Hereβs an example region mask generated for the North Atlantic Subpolar and Subtropical Gyres.
region_defs = {
'NAtl-STG': [
{'match': {'REGION_MASK': [6]}, 'bounds': {'TLAT': [32.0, 42.0], 'TLONG': [310.0, 350.0]}}
],
'NAtl-SPG': [
{'match': {'REGION_MASK': [6]}, 'bounds': {'TLAT': [50.0, 60.0], 'TLONG': [310.0, 350.0]}}
],
}
mask3d = pop_tools.region_mask_3d(grid_name, region_defs=region_defs, mask_name='N. Atlantic Gyres')
visualize_mask(mask3d)

%load_ext watermark
%watermark -d -iv -m -g -h
Compiler : GCC 11.3.0
OS : Linux
Release : 5.15.0-1004-aws
Machine : x86_64
Processor : x86_64
CPU cores : 2
Architecture: 64bit
Hostname: build-21213814-project-451810-pop-tools
Git hash: d3c80c0576ae4838c0e04a0157734eb0c977e613
pop_tools : 2023.3.0.post2+dirty
xarray : 2023.6.0
matplotlib: 3.7.1
numpy : 1.24.4