gcpy.grid

Module containing variables and functions that define and manipulate GEOS-Chem horizontal and vertical grids

Functions

calc_delta_lon(lon_edge)

Compute grid cell longitude widths from an edge vector.

calc_rectilinear_grid_area(lon_edge, lat_edge)

Compute grid cell areas (in m2) for a rectilinear grid.

calc_rectilinear_lat_edge(lat_stride, ...)

Compute latitude edge vector for a rectilinear grid.

calc_rectilinear_lon_edge(lon_stride, ...)

Compute longitude edge vector for a rectilinear grid.

call_make_grid(res, gridtype[, in_extent, ...])

Create a mask with NaN values removed from an input array

cartesian_to_latlon(x, y, z[, ret_xyz])

Convert a cartesian coordinate to latitude/longitude coordinates.

cartesian_to_spherical(x, y, z)

Convert cartesian coordinates to spherical in the form (theta, phi[, r]) with the origin remaining at the center of the original spherical coordinate system.

convert_lev_to_pres(dataset, pmid, pedge[, ...])

Convert lev dimension to pressure in a GEOS-Chem dataset.

csgrid_gmao(res[, offset])

Return cubed-sphere coordinates with GMAO face orientation.

get_grid_extents(data[, edges])

Get min and max lat and lon from an input GEOS-Chem xarray dataset or grid dict.

get_ilev_coord([n_lev, AP_edge, BP_edge, ...])

Returns the eta values (defined as (A/P0) + B) at vertical level edges.

get_ind_of_pres(dataset, pres)

Get index of pressure level that contains the requested pressure value.

get_input_res(data)

Returns resolution of dataset passed to compare_single_level or compare_zonal_means.

get_lev_coord([n_lev, AP_edge, BP_edge, ...])

Returns the eta values (defined as (A/P0) + B) at vertical level midpoints.

get_nearest_model_data(gc_data, lon_value, ...)

Returns GEOS-Chem model data at the grid box closest to a given (lat, lon) location, dispatching to the appropriate function depending on the grid type.

get_nearest_model_data_cs(gc_data, ...[, ...])

Returns GEOS-Chem model data (on a cubed-sphere grid) at the grid box closest to a given (lat, lon) location.

get_nearest_model_data_ll(gc_data, ...[, ...])

Returns GEOS-Chem model data (on a lat/lon grid) at the grid box closest to a given (lat, lon) location.

get_pressure_indices(pedge, pres_range)

Get indices where edge pressure values are within a given pressure range.

get_troposphere_mask(ds)

Returns a mask array for picking out the tropospheric grid boxes.

get_vert_grid(dataset[, AP, BP, p_sfc])

Determine vertical grid of input dataset.

latlon_to_cartesian(lon, lat)

Convert latitude/longitude coordinates along the unit sphere to cartesian coordinates defined by a vector pointing from the sphere's center to its surface.

make_grid_cs(csres)

Creates a cubed-sphere grid description.

make_grid_ll(llres[, in_extent, out_extent])

Creates a lat/lon grid description.

make_grid_sg(csres, stretch_factor, ...)

Creates a stretched-grid grid description.

pad_pressure_edges(pedge_ind, max_ind, pmid_len)

Add outer indices to edge pressure index list.

rotate_sphere_3d(theta, phi, r, rot_ang[, ...])

Rotate a spherical coordinate in the form (theta, phi[, r]) about the indicated axis, rot_axis.

spherical_to_cartesian(theta, phi[, r])

Convert spherical coordinates in the form (theta, phi[, r]) to cartesian, with the origin at the center of the original spherical coordinate system.

Classes

CSGrid(c[, offset])

Generator for cubed-sphere grid geometries.

VertGrid([AP, BP, p_sfc])

Class that defines a vertical grid given the Ap and Bp grid parameters and surface pressure.

gcpy.grid.get_troposphere_mask(ds)[source]

Returns a mask array for picking out the tropospheric grid boxes.

Parameters:

ds (xarray.Dataset) – Dataset containing certain met field variables (i.e. Met_TropLev, Met_BXHEIGHT).

Returns:

tropmask – Tropospheric mask. False denotes grid boxes that are in the troposphere and True in the stratosphere (as per Python masking logic).

Return type:

numpy.ndarray

gcpy.grid.get_input_res(data)[source]

Returns resolution of dataset passed to compare_single_level or compare_zonal_means.

Parameters:

data (xarray.Dataset) – Input GEOS-Chem dataset.

Returns:

  • res (str or int) – Lat/lon res of the form ‘latresxlonres’ or cubed-sphere resolution.

  • gridtype (str) – ‘ll’ for lat/lon or ‘cs’ for cubed-sphere.

gcpy.grid.call_make_grid(res, gridtype, in_extent=None, out_extent=None, sg_params=None)[source]

Create a mask with NaN values removed from an input array

Parameters:
  • res (str or int) – Resolution of grid (format ‘latxlon’ or csres).

  • gridtype (str) – ‘ll’ for lat/lon or ‘cs’ for cubed-sphere.

  • in_extent (list of float, optional) – Describes minimum and maximum latitude and longitude of input data in the format [minlon, maxlon, minlat, maxlat]. Default value: GLOBAL_LL_EXTENT, i.e. [-180, 180, -90, 90]

  • out_extent (list of float, optional) – Desired minimum and maximum latitude and longitude of output grid in the format [minlon, maxlon, minlat, maxlat]. Default value: GLOBAL_LL_EXTENT, i.e. [-180, 180, -90, 90]

  • sg_params (list, optional) – Desired stretched-grid parameters in the format [stretch_factor, target_longitude, target_latitude]. Will trigger stretched-grid creation if not default values. Default value: NO_STRETCH_SG_PARAMS (no stretching, [1, 170, -90])

Returns:

  • grid (dict) – The created grid.

  • grid_list (list of dict or None) – List of grids if gridtype is ‘cs’, otherwise None.

gcpy.grid.get_grid_extents(data, edges=True)[source]

Get min and max lat and lon from an input GEOS-Chem xarray dataset or grid dict.

Parameters:
  • data (xarray.Dataset or dict) – A GEOS-Chem dataset or a grid dict.

  • edges (bool, optional) – Whether grid extents should use cell edges instead of centers. Default value: True

Returns:

  • minlon (float) – Minimum longitude of data grid.

  • maxlon (float) – Maximum longitude of data grid.

  • minlat (float) – Minimum latitude of data grid.

  • maxlat (float) – Maximum latitude of data grid.

gcpy.grid.get_vert_grid(dataset, AP=None, BP=None, p_sfc=1013.25)[source]

Determine vertical grid of input dataset.

Parameters:
  • dataset (xr.Dataset) – A GEOS-Chem output dataset.

  • AP (list-like, optional) – Hybrid grid parameter A (hPa).

  • BP (list-like, optional) – Hybrid grid parameter B (unitless).

  • p_sfc (float, optional) – Surface pressure in hPa. Default value: 1013.25

Returns:

  • pedge (numpy.ndarray) – Edge pressure values for vertical grid.

  • p_mid (numpy.ndarray) – Midpoint pressure values for vertical grid.

  • nlev (int) – Number of levels in vertical grid.

gcpy.grid.get_ilev_coord(n_lev=72, AP_edge=None, BP_edge=None, top_down=False, gchp_indices=False)[source]

Returns the eta values (defined as (A/P0) + B) at vertical level edges. These are used to define the “ilev” netCDF coordinate variable.

Parameters:
  • n_lev (int, optional) – Number of levels in the grid. Default value: 72

  • AP_edge (list-like, optional) – Hybrid grid parameter A (hPa), with values placed on level edges. If not specified, values from the _GEOS_72L_AP array in this module will be used.

  • BP_edge (list-like, optional) – Hybrid grid parameter B (unitless), with values placed on level edges. If not specified, values from the _GEOS_72L_BP array in this module will be used.

  • top_down (bool, optional) – Set this to True if the eta coordinate will be arranged from top-of-atm downward (True) or from the surface upward (False). Default value: False

  • gchp_indices (bool, optional) – Set this to True to return an array of indices (as is used in GCHP files). Default value: False

Returns:

ilev – List of eta values at vertical grid edges.

Return type:

numpy.ndarray

gcpy.grid.get_lev_coord(n_lev=72, AP_edge=None, BP_edge=None, top_down=False, gchp_indices=False)[source]

Returns the eta values (defined as (A/P0) + B) at vertical level midpoints. These are used to define the “lev” netCDF coordinate variable.

Parameters:
  • n_lev (int, optional) – Number of levels in the grid. Default value: 72

  • AP_edge (list-like, optional) – Hybrid grid parameter A (hPa), with values placed on level edges. If not specified, values from the _GEOS_72L_AP array in this module will be used.

  • BP_edge (list-like, optional) – Hybrid grid parameter B (unitless), with values placed on level edges. If not specified, values from the _GEOS_72L_BP array in this module will be used.

  • top_down (bool, optional) – Set this to true if the eta coordinate will be arranged from top-of-atm downward (True) or from the surface upward (False). Default value: False

  • gchp_indices (bool, optional) – Set this to True to return an array of indices (as is used in GCHP files). Default value: False

Returns:

lev – List of eta values at vertical grid midpoints.

Return type:

numpy.ndarray

gcpy.grid.get_pressure_indices(pedge, pres_range)[source]

Get indices where edge pressure values are within a given pressure range.

Parameters:
  • pedge (numpy.ndarray) – Edge pressure values for the vertical grid.

  • pres_range (list of float) – Contains minimum and maximum pressure.

Returns:

indices – Indices where edge pressure values are within the given pressure range.

Return type:

numpy.ndarray

gcpy.grid.pad_pressure_edges(pedge_ind, max_ind, pmid_len)[source]

Add outer indices to edge pressure index list.

Parameters:
  • pedge_ind (list) – List of edge pressure indices.

  • max_ind (int) – Maximum index.

  • pmid_len (int) – Length of pmid which should not be exceeded by indices.

Returns:

pedge_ind – List of edge pressure indices, possibly with new minimum and maximum indices.

Return type:

list

gcpy.grid.get_ind_of_pres(dataset, pres)[source]

Get index of pressure level that contains the requested pressure value.

Parameters:
Returns:

index – Index of level in dataset that corresponds to requested pressure.

Return type:

int

gcpy.grid.convert_lev_to_pres(dataset, pmid, pedge, lev_type='pmid')[source]

Convert lev dimension to pressure in a GEOS-Chem dataset.

Parameters:
  • dataset (xarray.Dataset) – GEOS-Chem dataset.

  • pmid (numpy.ndarray) – Midpoint pressure values.

  • pedge (numpy.ndarray) – Edge pressure values.

  • lev_type (str, optional) – Denote whether lev is ‘pedge’ or ‘pmid’ if grid is not 72/73 or 47/48 levels. Default value: ‘pmid’

Returns:

dataset – Input dataset with “lev” dimension values replaced with pressure values.

Return type:

xarray.Dataset

class gcpy.grid.VertGrid(AP=None, BP=None, p_sfc=1013.25)[source]

Bases: object

Class that defines a vertical grid given the Ap and Bp grid parameters and surface pressure.

Parameters:
  • AP (list-like) – Hybrid-grid A parameter (hPa).

  • BP (list-like) – Hybrid-grid B parameter (unitless).

  • p_sfc (float, optional) – Surface pressure in hPa. Default value: 1013.25

p_edge()[source]

Compute pressure at grid cell edges.

Returns:

p_edge – Edge pressure values (hPa) computed as AP + BP * p_sfc.

Return type:

numpy.ndarray

p_mid()[source]

Compute pressure at grid cell midpoints.

Returns:

p_mid – Midpoint pressure values (hPa) computed as the average of adjacent edge pressure values.

Return type:

numpy.ndarray

gcpy.grid.make_grid_ll(llres, in_extent=None, out_extent=None)[source]

Creates a lat/lon grid description.

Parameters:
  • llres (str) – Lat/lon resolution in ‘latxlon’ format (e.g. ‘4x5’).

  • in_extent (list of float, optional) – Describes minimum and maximum latitude and longitude of initial grid in the format [minlon, maxlon, minlat, maxlat]. Default value: [-180, 180, -90, 90]

  • out_extent (list of float, optional) – Describes minimum and maximum latitude and longitude of target grid in the format [minlon, maxlon, minlat, maxlat]. Needed to trim extent of input data. Default value: [] (assumes value of in_extent)

Returns:

llgrid

Dict grid description of format {‘lat’lat midpoints,

’lon’ : lon midpoints, ‘lat_b’ : lat edges, ‘lon_b’ : lon edges}.

Return type:

dict

gcpy.grid.make_grid_cs(csres)[source]

Creates a cubed-sphere grid description.

Parameters:

csres (int) – Cubed-sphere resolution of target grid.

Returns:

  • csgrid (dict) – Dict of format {‘lat’: lat midpoints, ‘lon’: lon midpoints, ‘lat_b’: lat edges, ‘lon_b’: lon edges} where each value has an extra face dimension of length 6.

  • csgrid_list (list of dict) – List of dicts separated by face index.

gcpy.grid.make_grid_sg(csres, stretch_factor, target_lon, target_lat)[source]

Creates a stretched-grid grid description.

Parameters:
  • csres (int) – Cubed-sphere resolution of target grid.

  • stretch_factor (float) – Stretch factor of target grid.

  • target_lon (float) – Target stretching longitude of target grid.

  • target_lat (float) – Target stretching latitude of target grid.

Returns:

  • csgrid (dict) – Dict of format {‘lat’: lat midpoints, ‘lon’: lon midpoints, ‘lat_b’: lat edges, ‘lon_b’: lon edges} where each value has an extra face dimension of length 6.

  • csgrid_list (list of dict) – List of dicts separated by face index.

gcpy.grid.calc_rectilinear_lon_edge(lon_stride, center_at_180)[source]

Compute longitude edge vector for a rectilinear grid.

Parameters:
  • lon_stride (float) – Stride length in degrees. For example, for a standard GEOS-Chem Classic 4x5 grid, lon_stride would be 5.

  • center_at_180 (bool) – Whether or not the grid should have a cell center at 180 degrees (i.e. on the date line). If true, the first grid cell is centered on the date line; if false, the first grid edge is on the date line.

Returns:

lon_edge – Longitudes of cell edges in degrees East.

Return type:

numpy.ndarray

Notes

All values are forced to be between [-180,180]. For a grid with N cells in each band, N+1 edges will be returned, with the first and last value being duplicates.

gcpy.grid.calc_rectilinear_lat_edge(lat_stride, half_polar_grid)[source]

Compute latitude edge vector for a rectilinear grid.

Parameters:
  • lat_stride (float) – Stride length in degrees. For example, for a standard GEOS-Chem Classic 4x5 grid, lat_stride would be 4.

  • half_polar_grid (bool) – Whether or not the grid should be “half-polar” (i.e. bands at poles are half the size). In either case the grid will start and end at -/+ 90, but when half_polar_grid is True, the first and last bands will have a width of 1/2 the normal lat_stride.

Returns:

lat_edge – Latitudes of cell edges in degrees North.

Return type:

numpy.ndarray

Notes

All values are forced to be between [-90,90]. For a grid with N cells in each band, N+1 edges will be returned, with the first and last value being duplicates.

gcpy.grid.calc_rectilinear_grid_area(lon_edge, lat_edge)[source]

Compute grid cell areas (in m2) for a rectilinear grid.

Parameters:
  • lon_edge (numpy.ndarray) – Grid box longitude edges (in degrees east).

  • lat_edge (numpy.ndarray) – Grid box latitude edges (in degrees north).

Returns:

grid_area – Array of grid box areas in m2.

Return type:

numpy.ndarray

gcpy.grid.calc_delta_lon(lon_edge)[source]

Compute grid cell longitude widths from an edge vector.

Parameters:

lon_edge (numpy.ndarray) – Vector of longitude edges, in degrees East.

Returns:

lon_delta – Width of each cell in degrees East.

Return type:

numpy.ndarray

Notes

Accounts for looping over the domain.

gcpy.grid.csgrid_gmao(res, offset=-10)[source]

Return cubed-sphere coordinates with GMAO face orientation.

Parameters:
  • res (int) – Cubed-sphere resolution.

  • offset (float, optional) – Degrees to offset grid in the longitudinal direction. Default value: -10

Returns:

grid – Dictionary with keys ‘lon’, ‘lat’, ‘lon_b’, ‘lat_b’ containing the cubed-sphere grid coordinates.

Return type:

dict

Notes

This function was originally written by Jiawei Zhuange and included in package cubedsphere: https://github.com/JiaweiZhuang/cubedsphere

class gcpy.grid.CSGrid(c, offset=None)[source]

Bases: object

Generator for cubed-sphere grid geometries.

CSGrid computes the latitutde and longitudes of cell centers and edges on a cubed-sphere grid, providing a way to retrieve these geometries on-the-fly if your model output data does not include them.

lon_center, lat_center

Lat/lon coordinates for each cell center along the cubed-sphere mesh.

Type:

numpy.ndarray

lon_edge, lat_edge

Lat/lon coordinates for the midpoint of the edges separating each element on the cubed-sphere mesh.

Type:

numpy.ndarray

xyz_center, xyz_edge

As above, except coordinates are projected into a 3D cartesian space with common origin to the original lat/lon coordinate system, assuming a unit sphere.

Type:

numpy.ndarray

Notes

This class was originally written by Jiawei Zhuange and included in package cubedsphere: https://github.com/JiaweiZhuang/cubedsphere

Parameters:
  • c (int) –

    Number of edges along each cubed-sphere edge.

    C

    Lat/Lon Resolution

    24

    4 deg x 5 deg

    48,45

    2 deg x 2.5 deg

    96,90

    1 deg x 1.25 deg

    192,180

    0.5 deg x 0.625 deg

    384,360

    0.25 deg x 0.3125 deg

    720

    0.125 deg x 0.15625 deg

  • offset (float, optional) – Degrees to offset the first faces’ edge in the latitudinal direction. If not passed, then the western edge of the first face will align with the prime meridian.

Notes

This function was originally written by Jiawei Zhuange and included in package cubedsphere: https://github.com/JiaweiZhuang/cubedsphere

__init__(c, offset=None)[source]
Parameters:
  • c (int) –

    Number of edges along each cubed-sphere edge.

    C

    Lat/Lon Resolution

    24

    4 deg x 5 deg

    48,45

    2 deg x 2.5 deg

    96,90

    1 deg x 1.25 deg

    192,180

    0.5 deg x 0.625 deg

    384,360

    0.25 deg x 0.3125 deg

    720

    0.125 deg x 0.15625 deg

  • offset (float, optional) – Degrees to offset the first faces’ edge in the latitudinal direction. If not passed, then the western edge of the first face will align with the prime meridian.

Notes

This function was originally written by Jiawei Zhuange and included in package cubedsphere: https://github.com/JiaweiZhuang/cubedsphere

gcpy.grid.latlon_to_cartesian(lon, lat)[source]

Convert latitude/longitude coordinates along the unit sphere to cartesian coordinates defined by a vector pointing from the sphere’s center to its surface.

Parameters:
  • lon (float) – Longitude coordinate in radians.

  • lat (float) – Latitude coordinate in radians.

Returns:

x, y, z – Cartesian coordinates on the unit sphere.

Return type:

float

Notes

This function was originally written by Jiawei Zhuange and included in package cubedsphere: https://github.com/JiaweiZhuang/cubedsphere

gcpy.grid.cartesian_to_latlon(x, y, z, ret_xyz=False)[source]

Convert a cartesian coordinate to latitude/longitude coordinates. Optionally return the original cartesian coordinate as a tuple.

Parameters:
  • x (float) – Cartesian x-coordinate.

  • y (float) – Cartesian y-coordinate.

  • z (float) – Cartesian z-coordinate.

  • ret_xyz (bool, optional) – If True, also return the normalized cartesian coordinate vector. Default value: False

Returns:

  • lon (float) – Longitude in radians.

  • lat (float) – Latitude in radians.

  • xyz (numpy.ndarray, optional) – Normalized cartesian coordinate vector. Only returned if ret_xyz is True.

Notes

This function was originally written by Jiawei Zhuange and included in package cubedsphere: https://github.com/JiaweiZhuang/cubedsphere

gcpy.grid.spherical_to_cartesian(theta, phi, r=1)[source]

Convert spherical coordinates in the form (theta, phi[, r]) to cartesian, with the origin at the center of the original spherical coordinate system.

Parameters:
  • theta (float) – Azimuthal angle in radians.

  • phi (float) – Polar angle in radians.

  • r (float, optional) – Radius. Default value: 1

Returns:

x, y, z – Cartesian coordinates.

Return type:

float

Notes

This function was originally written by Jiawei Zhuange and included in package cubedsphere: https://github.com/JiaweiZhuang/cubedsphere

gcpy.grid.cartesian_to_spherical(x, y, z)[source]

Convert cartesian coordinates to spherical in the form (theta, phi[, r]) with the origin remaining at the center of the original spherical coordinate system.

Parameters:
  • x (float) – Cartesian x-coordinate.

  • y (float) – Cartesian y-coordinate.

  • z (float) – Cartesian z-coordinate.

Returns:

  • theta (float) – Azimuthal angle in radians.

  • phi (float) – Polar angle in radians.

  • r (float) – Radius.

Notes

This function was originally written by Jiawei Zhuange and included in package cubedsphere: https://github.com/JiaweiZhuang/cubedsphere

gcpy.grid.rotate_sphere_3d(theta, phi, r, rot_ang, rot_axis='x')[source]

Rotate a spherical coordinate in the form (theta, phi[, r]) about the indicated axis, rot_axis.

This method accomplishes the rotation by projecting to a cartesian coordinate system and performing a solid body rotation around the requested axis.

Parameters:
  • theta (float) – Azimuthal angle in radians.

  • phi (float) – Polar angle in radians.

  • r (float) – Radius.

  • rot_ang (float) – Rotation angle in radians.

  • rot_axis (str, optional) – Axis about which to rotate (‘x’, ‘y’, or ‘z’). Default value: ‘x’

Returns:

  • theta_new (float) – New azimuthal angle in radians after rotation.

  • phi_new (float) – New polar angle in radians after rotation.

  • r_new (float) – New radius after rotation.

Notes

This function was originally written by Jiawei Zhuange and included in package cubedsphere: https://github.com/JiaweiZhuang/cubedsphere

gcpy.grid.get_nearest_model_data_cs(gc_data, gc_cs_grid, lon_value, lat_value, varlist=None)[source]

Returns GEOS-Chem model data (on a cubed-sphere grid) at the grid box closest to a given (lat, lon) location.

Parameters:
  • gc_data (xarray.DataArray or xarray.Dataset) – GEOS-Chem model data for a single variable.

  • gc_cs_grid (xarray.Dataset) – Coordinate arrays defining the cubed-sphere grid.

  • lon_value (float) – Longitude at the location of interest.

  • lat_value (float) – Latitude at the location of interest.

  • varlist (list of str, optional) – List of data variables to include in the output.

Returns:

dataframe – Model data closest to the observation site.

Return type:

pandas.DataFrame

gcpy.grid.get_nearest_model_data_ll(gc_data, lon_value, lat_value, varlist=None)[source]

Returns GEOS-Chem model data (on a lat/lon grid) at the grid box closest to a given (lat, lon) location.

Parameters:
  • gc_data (xarray.DataArray or xarray.Dataset) – GEOS-Chem model data.

  • lon_value (float) – Longitude at the location of interest.

  • lat_value (float) – Latitude at the location of interest.

  • varlist (list of str, optional) – List of data variables to include in the output.

Returns:

dataframe – Model data closest to the observation site.

Return type:

pandas.DataFrame

gcpy.grid.get_nearest_model_data(gc_data, lon_value, lat_value, gc_cs_grid=None, varlist=None)[source]

Returns GEOS-Chem model data at the grid box closest to a given (lat, lon) location, dispatching to the appropriate function depending on the grid type.

Parameters:
  • gc_data (xarray.DataArray or xarray.Dataset) – GEOS-Chem data for a single variable.

  • lon_value (float) – Longitude at the location of interest.

  • lat_value (float) – Latitude at the location of interest.

  • gc_cs_grid (xarray.Dataset, optional) – Coordinate arrays defining the cubed-sphere grid. This can be obtained as the output of function gcpy.util.extract_data().

  • varlist (list of str, optional) – List of data variables to include in the output.

Returns:

dataframe – Model data closest to the observation site.

Return type:

pandas.DataFrame