Usage

Reading Data

We can either read mds data using

cubedsphere.open_ascii_dataset(outdir, return_grid=True, **kwargs)[source]

Wrapper that opens simulation outputs from standard mitgcm outputs.

Parameters
outdir: string

Output directory

return_grid: Boolean

Return a grid generated with xmitgcm.get_grid_from_input

**kwargs

everything else that is passed to xmitgcm

Returns
ds: xarray Dataset

Dataset of simulation output

grid: xarray Dataset

Only if return_grid is True Grid generated with xmitgcm.get_grid_from_input.

Warning

You need to use useSingleCpuIO=.TRUE. if you want to use ascii files (the default MITgcm output)

We can also read NETCDF data which has been outputed from the mnc package (depreciated).

cubedsphere.open_mnc_dataset(outdir, iternumber, fname_list=['state'])[source]

Wrapper that opens simulation outputs from mnc outputs. NOT TESTED.

Parameters
outdir: string

Output directory

iternumber: integer

iteration number of output file

fname_list: list

List of NetCDF file prefixes to read (no need to specify grid files here)

Returns
ds: xarray Dataset

Dataset of simulation output

Regridding

class cubedsphere.Regridder(ds, cs_grid, input_type='cs', d_lon=5, d_lat=4, concat_mode=False, filename='weights', method='conservative', **kwargs)[source]

Class that wraps the xESMF regridder for cs geometry.

Only two methods are possible with the cs geometry: conservative when using concat_mode=False (requires lon_b to have different shape from lon and lat_b from lat) or nearest_s2d when using concat_mode=True. Conservative regridding should be used if possible!

Notes

You can find more examples in the examples directory

Examples

>>> import cubedsphere as cs  # import cubedsphere
>>> outdir = "../run"  # specify output directory
>>> # open Dataset
>>> ds_ascii, grid = cs.open_ascii_dataset(outdir_ascii, iters='all', prefix = ["T","U","V","W"])
>>> # regrid dataset
>>> regrid = cs.Regridder(ds_ascii, grid)
>>> ds_reg = regrid()
Attributes
regridder: list or xESMF regrid object

(contains) the initialized xESMF regridder

grid: dict

output grid

__init__(ds, cs_grid, input_type='cs', d_lon=5, d_lat=4, concat_mode=False, filename='weights', method='conservative', **kwargs)[source]

Build the regridder. This step will create the output grid and weight files which will then be used to regrid the dataset.

Parameters
ds: xarray DataSet

Dataset to be regridded. Dataset must contain input grid information.

cs_grid: xarray DataSet

Dataset containing cubedsphere grid information. If input_type==”cs”: Input grid. Required! If input_type==”ll”: Output cs grid. Required!

input_type: string

Needs to be “cs” or “ll”. Will result in an error if something else is used here. If “cs”: input=cs grid -> regrid to longitude lattitude (ll) If “ll”: input=ll grid -> regrid to cs grid

d_lon: integer

Longitude step size, i.e. grid resolution. Only used if input_type==”cs”.

d_lat: integer

Latitude step size, i.e. grid resolution. Only used if input_type==”cs”.

concat_mode: boolean

use one regridding instance instead of one regridder for each face. Only used if input_type==”cs”.

filename: string

filename for weights (weights will be name filename(+ _tile{i}).nc)

method: string

Regridding method. See xe.Regridder for options.

kwargs

Optional parameters that are passed to xe.Regridder (see xe.Regridder for options).

__call__(**kwargs)[source]

Wrapper that carries out the regridding from cubedsphere to latlon.

Parameters
Returns
ds: xarray DataSet

regridded Dataset

Create XGCM grids

cubedsphere.init_grid_CS(ds=None, grid_dir=None, **kwargs)[source]

Init a xgcm grid for a cubedsphere dataset. Useful for raw datasets.

Parameters
grid_dir: string

direction where the grid can be found (optional)

ds: xarray DataSet

dataset that contains the grid (optional)

Returns
grid: xgcm.grid

xgcm grid

cubedsphere.init_grid_LL(ds=None, grid_dir=None, **kwargs)[source]

Init a xgcm grid for a latlon dataset. Useful for regridded datasets

Parameters
grid_dir: string

direction where the grid can be found (optional)

ds: xarray DataSet

dataset that contains the grid (optional)

Returns
grid: xgcm.grid

xgcm grid

Plotting

cubedsphere.plotCS(dr, ds, mask_size=None, **kwargs)[source]

A quick way to plot cubed-sphere data of resolution 6*N*N. Wrapping plotCS_quick_raw to work with xarray objects

Parameters
dr: xarray DataArray

The dimensions must be (y,x)

ds: xarray DataSet (the parent DataSet of dr)

Must contain XC, YC as coordinate variables.

mask_size: None or int

The overlap size of individual tiles. If None is chosen one might likely experience issues

**kwargs

Other keyword arguments such as cmap will be passed to plotCS_quick_raw() and subsequently to plt.pcolormesh()

Returns
ph: list

List of mappabales

cubedsphere.overplot_wind(ds_reg, U, V, stepsize=1, **kwargs)[source]

Quick and dirty function for overplotting wind of a regridded dataset

Parameters
ds_reg: xarray DataSet

regridded dataset

stepsize: integer

specify the stepsize for which wind arrows should be plotted