epygram.geometries.GaussGeometry — Gauss Geometry class

Contains the classes for 3D geometries of fields.


class epygram.geometries.GaussGeometry.GaussGeometry(name, grid, dimensions, vcoordinate, position_on_horizontal_grid='__unknown__', geoid=None)[source]

Bases: epygram.geometries.AbstractGeometry.Geometry

Handles the geometry for a Global Gauss grid 3-Dimensions Field.

Parameters
  • name – Name of geometrical type of representation of points on the Globe. Name must be among [‘rotated_reduced_gauss’, ‘reduced_gauss’, ‘regular_gauss’

  • grid – Handles description of the horizontal grid.

  • dimensions – Handles grid dimensions.

  • vcoordinate – Handles vertical geometry parameters.

  • position_on_horizontal_grid

    Position of points w/r to the horizontal. among: [‘upper-right’, ‘upper-left’,

    ’lower-left’, ‘lower-right’, ‘center-left’, ‘center-right’, ‘lower-center’, ‘upper-center’, ‘center’, ‘__unknown__’]

  • geoid – To specify geoid shape.

default_cartopy_CRS()

Note

Requires plugin: with_cartopy (config.activate_plugins)

Create a cartopy.crs appropriate to the Geometry.

distance_to_nearest_neighbour_ij(i, j)[source]

Returns the distance to the nearest point of (i,j) point

Parameters
  • i – X index of point in the 2D matrix of gridpoints

  • j – Y index of point in the 2D matrix of gridpoints

distance_to_nearest_neighbour_ll(lon, lat)[source]

Returns the local resolution at the nearest point of lon/lat. It’s the distance between this point and its closest neighbour.

Parameters
  • lon – longitude of the point in degrees

  • lat – latitude of the point in degrees

fill_maskedvalues(data, fill_value=None)[source]

Returns a copy of data with ‘real’ masked values (i.e. not those linked to reduced Gauss) filled with fill_value. data must be already 4D for simplicity reasons.

get_datashape(force_dimZ=None, dimT=None, d4=False, **_)[source]

Returns the data shape according to the geometry.

Parameters
  • force_dimZ – if supplied, force the Z dimension instead of that of the vertical geometry

  • dimT – if supplied, is the time dimension to be added to the data shape

  • d4

    • if True, shape is 4D (need to specify dimT)

    • if False, shape is 3D if dimZ > 1 else 2D

get_lonlat_grid(position=None, d4=False, nb_validities=0, force_longitudes=None, **_)[source]

Returns a tuple of two tables containing one the longitude of each point, the other the latitude, with 2D shape.

Shape of 2D data in Gauss grids:

  • grid[0, 0:Nj] is first (Northern) band of latitude, masked after Nj = number of longitudes for latitude j

  • grid[-1, 0:Nj] is last (Southern) band of latitude (idem).

Parameters
  • position – position of lonlat grid with respect to the model cell. Defaults to self.position_on_horizontal_grid.

  • d4

    • if True, returned values are shaped in a 4 dimensions array

    • if False, shape of returned values is determined with respect to geometry.

      d4=True requires nb_validities > 0

  • nb_validities – number of validities represented in data values

  • force_longitudes – if ‘positive’, the longitudes will be forced positive if ‘]-180,180]’, the longitudes will be in the ]-180, 180] interval

property gridpoints_number

Returns the number of gridpoints of the grid.

horizontally_flattened(data)[source]

Returns a copy of data with horizontal dimensions flattened and compressed (cf. numpy.ma.masked_array.compressed). data must be 4D for simplicity reasons.

ij2ll(i, j, position=None)[source]

Return the (lon, lat) coordinates of point (i,j), in degrees.

Parameters
  • i – X index of point in the 2D matrix of gridpoints

  • j – Y index of point in the 2D matrix of gridpoints

  • position – lat lon position to return with respect to the model cell. Defaults to self.position_on_horizontal_grid.

property isglobal
Returns

True if geometry is global

ll2ij(lon, lat, position=None)[source]

Return the (i, j) coordinates in the 2D matrix of gridpoints, of the gridpoint nearest to (lon, lat).

Parameters
  • lon – longitude of point in degrees

  • lat – latitude of point in degrees

  • position – lat lon position to return with respect to the model cell. Defaults to self.position_on_horizontal_grid.

map_factor(lon, lat)[source]

Returns the map factor at the given longitude/latitude(s)

Parameters
  • lon – longitude in degrees

  • lat – latitude in degrees

map_factor_field(position=None)[source]

Returns a new field whose data is the map factor over the field.

Parameters

position – grid position with respect to the model cell. Defaults to self.position_on_horizontal_grid.

meridian_resolution_j(j)[source]

Returns the average meridian resolution at longitude circle number j.

nearest_points(lon, lat, request, position=None, external_distance=None, squeeze=True)[source]

Returns a list of the (i, j) position of the points needed to perform an interpolation.

Parameters
  • lon – longitude of point in degrees.

  • lat – latitude of point in degrees.

  • request

    criteria for selecting the points, among: * {‘n’:’1’} - the nearest point * {‘n’:’2*2’} - the 2*2 square points around the position * {‘n’:’4*4’} - the 4*4 square points around the position * {‘n’:’N*N’} - the N*N square points around the position: N must be even * {‘radius’:xxxx, ‘shape’:’square’} - the points which are xxxx metres

    around the position in X or Y direction

    • {‘radius’:xxxx, ‘shape’:’circle’} - the points within xxxx metres around the position. (default shape == circle)

  • position – position in the model cell of the lat lon position. Defaults to self.position_on_horizontal_grid.

  • external_distance – can be a dict containing the target point value and an external field on the same grid as self, to which the distance is computed within the 4 horizontally nearest points; e.g. {‘target_value’:4810, ‘external_field’:a_3DField_with_same_geometry}. If so, the nearest point is selected with distance = |target_value - external_field.data|

  • squeeze – True to suppress useless dimensions

Return type

general output form is [list, list, …, list] with as many list items as the length of lon/lat. Each list item is of the form [tuple, tuple, …, tuple] with as many tuples as the request implies. A tuple represents one of the nearest points associated with one value taken from lon/lat. Each tuple as the form (i, j).

Dimensions with a length of one are removed except if squeeze is False. If squeeze is True and if request implies only one nearest point, the list item of the general output form is replaced by the tuple item; if length of lon/lat is one, the output is directly the list item of the general output form. Hence, if length of lon/lat is one and the request implies only one point, the output is a tuple.

In case of a simple square request, output is actually an array. Otherwise, the output is as described (it cannot be an array because the number of nearest points can vary with the entry point).

point_is_inside_domain_ij(i, j, margin=- 0.1)[source]

Returns True if the point(s) of lon/lat coordinates is(are) inside the field.

Parameters
  • i – X index of point in the 2D matrix of gridpoints

  • j – Y index of point in the 2D matrix of gridpoints

  • margin – DEPRECATED

point_is_inside_domain_ll(lon, *_, **__)[source]

Returns True if the point(s) of lon/lat coordinates is(are) inside the field. This is always the case in Gauss grids, no real meaning.

Parameters
  • lon – longitude of the point in degrees

  • lat – latitude of the point in degrees

  • margin – considers the point inside if at least ‘margin’ points far from the border. The -0.1 default is a safety for precision errors.

  • position – position of the grid with respect to the model cell. Defaults to self.position_on_horizontal_grid.

reproject_wind_on_lonlat(u, v, lon, lat, map_factor_correction=True, reverse=False)[source]

Reprojects a wind vector (u, v) on rotated/stretched sphere onto real sphere, i.e. with components on true zonal/meridian axes.

Parameters
  • u – the u == zonal-on-the-grid component of wind

  • v – the v == meridian-on-the-grid component of wind

  • lon – longitudes of points in degrees

  • lat – latitudes of points in degrees

  • map_factor_correction – applies a correction of magnitude due to map factor.

  • reverse – if True, apply the reverse reprojection.

lon/lat are coordinates on real sphere.

reshape_data(data, first_dimension=None, d4=False)[source]

Returns a 2D data (horizontal dimensions) reshaped from 1D, according to geometry.

Parameters
  • data – the 1D data (or 3D with a T and Z dimensions, or 2D with either a T/Z dimension, to be specified), of dimension concording with geometry. In case data is 3D, T must be first dimension and Z the second.

  • first_dimension – in case data is 2D, specify what is the first dimension of data among (‘T’, ‘Z’)

  • d4

    • if True, returned values are shaped in a 4 dimensions array

    • if False, shape of returned values is determined with respect to geometry

resolution_field(direction='meridian')[source]

Returns a field whose values are the local resolution in m.

Parameters

direction – among (‘zonal’, ‘meridian’), direction in which the resolution is computed.

resolution_field_from_stretching(*args, **kwargs)

Returns a field which values are the local resolution computed as the nominal resolution stretched locally by the map factor.

resolution_j(j)[source]

Returns the average meridian resolution at longitude circle number j.

resolution_ll(lon, lat)[source]

Returns the average meridian resolution (worst directional resolution) at point position.

Parameters
  • lon – longitude of the point in degrees

  • lat – latitude of the point in degrees

zonal_resolution_j(j)[source]

Returns the average zonal resolution at longitude circle number j.