pyresample package

Subpackages

Submodules

pyresample.area_config module

Area config handling and creation utilities.

exception pyresample.area_config.AreaNotFound

Bases: KeyError

Exception raised when specified are is no found in file.

pyresample.area_config.convert_def_to_yaml(def_area_file, yaml_area_file)

Convert a legacy area def file to the yaml counter partself.

yaml_area_file will be overwritten by the operation.

pyresample.area_config.create_area_def(area_id, projection, width=None, height=None, area_extent=None, shape=None, upper_left_extent=None, center=None, resolution=None, radius=None, units=None, optimize_projection=False, **kwargs)

Create AreaDefinition from whatever information is known.

Parameters:
  • area_id (str) – ID of area

  • projection (pyproj CRS object, dict, str, int, tuple, object) – Projection parameters. This can be in any format understood by pyproj.crs.CRS.from_user_input(), such as a pyproj CRS object, proj4 dict, proj4 string, EPSG integer code, or others.

  • description (str, optional) – Description/name of area. Defaults to area_id

  • proj_id (str, optional) – ID of projection (deprecated)

  • units (str, optional) –

    Units that provided arguments should be interpreted as. This can be one of ‘deg’, ‘degrees’, ‘meters’, ‘metres’, and any parameter supported by the cs2cs -lu command. Units are determined in the following priority:

    1. units expressed with each variable through a DataArray’s attrs attribute.

    2. units passed to units

    3. units used in projection

    4. meters

  • width (str, optional) – Number of pixels in the x direction

  • height (str, optional) – Number of pixels in the y direction

  • area_extent (list, optional) – Area extent as a list (lower_left_x, lower_left_y, upper_right_x, upper_right_y)

  • shape (list, optional) – Number of pixels in the y and x direction (height, width)

  • upper_left_extent (list, optional) – Upper left corner of upper left pixel (x, y)

  • center (list, optional) – Center of projection (x, y)

  • resolution (list or float, optional) – Size of pixels: (dx, dy)

  • radius (list or float, optional) – Length from the center to the edges of the projection (dx, dy)

  • nprocs (int, optional) – Number of processor cores to be used

  • lons (numpy array, optional) – Grid lons

  • lats (numpy array, optional) – Grid lats

  • optimize_projection – Whether the projection parameters have to be optimized for a DynamicAreaDefinition.

Returns:

area – If shape and area_extent are found, an AreaDefinition object is returned. If only shape or area_extent can be found, a DynamicAreaDefinition object is returned

Return type:

pyresample.geometry.AreaDefinition or pyresample.geometry.DynamicAreaDefinition

Raises:

ValueError: – If neither shape nor area_extent could be found

Notes

  • resolution and radius can be specified with one value if dx == dy

  • If resolution and radius are provided as angles, center must be given or findable. In such a case, they represent [projection x distance from center[0] to center[0]+dx, projection y distance from center[1] to center[1]+dy]

pyresample.area_config.generate_area_def_rst_list(area_file: str) str

Create rst list of available area definitions with overview plot.

Parameters:

area_file – Path to area yaml file.

Returns:

rst list formatted string.

pyresample.area_config.get_area_def(area_id, area_name, proj_id, proj4_args, width, height, area_extent)

Construct AreaDefinition object from arguments.

Parameters:
  • area_id (str) – ID of area

  • area_name (str) – Description of area

  • proj_id (str) – ID of projection

  • proj4_args (dict, CRS, or str) – Projection information passed to pyproj’s CRS object

  • width (int) – Number of pixel in x dimension

  • height (int) – Number of pixel in y dimension

  • area_extent (list | tuple) – Area extent as a list of ints (LL_x, LL_y, UR_x, UR_y)

Returns:

area_def – AreaDefinition object

Return type:

object

pyresample.area_config.load_area(area_file_name, *regions)

Load area(s) from area file.

Parameters:
  • area_file_name (str, pathlib.Path, stream, or list thereof) – List of paths or streams. Any str or pathlib.Path will be interpreted as a path to a file. Any stream will be interpreted as containing a yaml definition file. To read directly from a string, use load_area_from_string().

  • regions (str argument list) – Regions to parse. If no regions are specified all regions in the file are returned

Returns:

area_defs – If one area name is specified a single AreaDefinition object is returned. If several area names are specified a list of AreaDefinition objects is returned

Return type:

pyresample.geometry.AreaDefinition or list

Raises:

AreaNotFound: – If a specified area name is not found

pyresample.area_config.load_area_from_string(area_strs, *regions)

Load area(s) from area strings.

Like load_area(), but load from string directly.

For the opposite (i.e. to create a YAML string from an area), use dump().

Parameters:
  • area_strs (str or List[str]) – Strings containing yaml definitions.

  • regions (str) – Regions to parse.

Returns:

area_defs – If one area name is specified a single AreaDefinition object is returned. If several area names are specified a list of AreaDefinition objects is returned

Return type:

pyresample.geometry.AreaDefinition or list

pyresample.area_config.parse_area_file(area_file_name, *regions)

Parse area information from area file.

Parameters:
  • area_file_name (str or list) – One or more paths to area definition files

  • regions (str argument list) – Regions to parse. If no regions are specified all regions in the file are returned

Returns:

area_defs – List of AreaDefinition objects

Return type:

list

Raises:

AreaNotFound: – If a specified area is not found

pyresample.boundary module

The Boundary classes.

pyresample.data_reduce module

Reduce data sets based on geographical information.

pyresample.data_reduce.get_valid_index_from_cartesian_grid(cart_grid, lons, lats, radius_of_influence)

Calculate relevant data indices using coarse data reduction of swath data by comparison with cartesian grid.

Parameters:
  • chart_grid (numpy array) – Grid of area cartesian coordinates

  • lons (numpy array) – Swath lons

  • lats (numpy array) – Swath lats

  • data (numpy array) – Swath data

  • radius_of_influence (float) – Cut off distance in meters

Returns:

valid_index – Boolean array of same size as lons and lats indicating relevant indices

Return type:

numpy array

pyresample.data_reduce.get_valid_index_from_lonlat_boundaries(boundary_lons, boundary_lats, lons, lats, radius_of_influence)

Find relevant indices from grid boundaries using the winding number theorem.

pyresample.data_reduce.get_valid_index_from_lonlat_grid(grid_lons, grid_lats, lons, lats, radius_of_influence)

Calculate relevant data indices using coarse data reduction of swath data by comparison with lon lat grid.

Parameters:
  • chart_grid (numpy array) – Grid of area cartesian coordinates

  • lons (numpy array) – Swath lons

  • lats (numpy array) – Swath lats

  • data (numpy array) – Swath data

  • radius_of_influence (float) – Cut off distance in meters

Returns:

valid_index – Boolean array of same size as lon and lat indicating relevant indices

Return type:

numpy array

pyresample.data_reduce.swath_from_cartesian_grid(cart_grid, lons, lats, data, radius_of_influence)

Make coarse data reduction of swath data by comparison with cartesian grid.

Parameters:
  • chart_grid (numpy array) – Grid of area cartesian coordinates

  • lons (numpy array) – Swath lons

  • lats (numpy array) – Swath lats

  • data (numpy array) – Swath data

  • radius_of_influence (float) – Cut off distance in meters

Returns:

(lons, lats, data) – Reduced swath data and coordinate set

Return type:

list of numpy arrays

pyresample.data_reduce.swath_from_lonlat_boundaries(boundary_lons, boundary_lats, lons, lats, data, radius_of_influence)

Make coarse data reduction of swath data by comparison with lon lat boundary.

Parameters:
  • boundary_lons (numpy array) – Grid of area lons

  • boundary_lats (numpy array) – Grid of area lats

  • lons (numpy array) – Swath lons

  • lats (numpy array) – Swath lats

  • data (numpy array) – Swath data

  • radius_of_influence (float) – Cut off distance in meters

Returns:

(lons, lats, data) – Reduced swath data and coordinate set

Return type:

list of numpy arrays

pyresample.data_reduce.swath_from_lonlat_grid(grid_lons, grid_lats, lons, lats, data, radius_of_influence)

Make coarse data reduction of swath data by comparison with lon lat grid.

Parameters:
  • grid_lons (numpy array) – Grid of area lons

  • grid_lats (numpy array) – Grid of area lats

  • lons (numpy array) – Swath lons

  • lats (numpy array) – Swath lats

  • data (numpy array) – Swath data

  • radius_of_influence (float) – Cut off distance in meters

Returns:

(lons, lats, data) – Reduced swath data and coordinate set

Return type:

list of numpy arrays

pyresample.geo_filter module

Filters based on geolocation validity.

class pyresample.geo_filter.GridFilter(area_def, filter, nprocs=1)

Bases: object

Geographic filter from a grid.

Parameters:
  • grid_ll_x (float) – Projection x coordinate of lower left corner of lower left pixel

  • grid_ll_y (float) – Projection y coordinate of lower left corner of lower left pixel

  • grid_ur_x (float) – Projection x coordinate of upper right corner of upper right pixel

  • grid_ur_y (float) – Projection y coordinate of upper right corner of upper right pixel

  • proj4_string (str) – Projection definition as a PROJ.4 string.

  • mask (numpy array) – Mask as boolean numpy array

filter(geometry_def, data)

Get coordinate definition and data where invalid lon/lats are removed.

get_valid_index(geometry_def)

Calculate valid_index array based on lons and lats.

Parameters:

geometry_def – Geometry definition (ex. SwathDefinition)

Returns:

Boolean numpy array of same shape as lons and lats

pyresample.geometry module

Classes for geometry operations.

class pyresample.geometry.AreaDefinition(area_id, description, proj_id, projection, width, height, area_extent, nprocs=1, lons=None, lats=None, dtype=<class 'numpy.float64'>)

Bases: _ProjectionDefinition

Holds definition of an area.

Parameters:
  • area_id (str) – Identifier for the area

  • description (str) – Human-readable description of the area

  • proj_id (str) – ID of projection

  • projection (dict or str or pyproj.crs.CRS) – Dictionary of PROJ parameters or string of PROJ or WKT parameters. Can also be a pyproj.crs.CRS object.

  • width (int) – x dimension in number of pixels, aka number of grid columns

  • height (int) – y dimension in number of pixels, aka number of grid rows

  • area_extent (list) – Area extent as a list (lower_left_x, lower_left_y, upper_right_x, upper_right_y)

  • nprocs (int, optional) – Number of processor cores to be used for certain calculations

area_id

Identifier for the area

Type:

str

description

Human-readable description of the area

Type:

str

proj_id

ID of projection

Type:

str

projection

Dictionary or string with Proj.4 parameters

Type:

dict or str

width

x dimension in number of pixels, aka number of grid columns

Type:

int

height

y dimension in number of pixels, aka number of grid rows

Type:

int

size

Number of points in grid

Type:

int

area_extent_ll

Area extent in lons lats as a tuple (lower_left_lon, lower_left_lat, upper_right_lon, upper_right_lat)

Type:

tuple

pixel_size_x

Pixel width in projection units

Type:

float

pixel_size_y

Pixel height in projection units

Type:

float

upper_left_extent

Coordinates (x, y) of upper left corner of upper left pixel in projection units

Type:

tuple

pixel_upper_left

Coordinates (x, y) of center of upper left pixel in projection units

Type:

tuple

pixel_offset_x

x offset between projection center and upper left corner of upper left pixel in units of pixels.

Type:

float

pixel_offset_y

y offset between projection center and upper left corner of upper left pixel in units of pixels..

Type:

float

crs

Coordinate reference system object similar to the PROJ parameters in proj_dict and proj_str. This is the preferred attribute to use when working with the pyproj library. Note, however, that this object is not thread-safe and should not be passed between threads.

Type:

pyproj.crs.CRS

crs_wkt

WellKnownText version of the CRS object. This is the preferred way of describing CRS information as a string.

Type:

str

cartesian_coords

Grid cartesian coordinates

Type:

object

aggregate(**dims)

Return an aggregated version of the area.

property area_extent

Tuple of this area’s extent (xmin, ymin, xmax, ymax).

colrow2lonlat(cols, rows)

Return lons and lats for the given image columns and rows.

Both scalars and arrays are supported. To be used with scarse data points instead of slices (see get_lonlats).

copy(**override_kwargs)

Make a copy of the current area.

This replaces the current values with anything in override_kwargs.

create_areas_def()

Generate YAML formatted representation of this area.

Deprecated. Use dump() instead.

create_areas_def_legacy()

Create area definition in legacy format.

crop_around(other_area)

Crop this area around other_area.

dump(filename=None)

Generate YAML formatted representation of this area.

For the opposite (i.e. to get an AreaDefinition from a YAML-formatted representation), see load_area_from_string().

Parameters:

filename (str or pathlib.Path or file-like object) – Yaml file location to dump the area to.

Returns:

If file is None returns yaml str

classmethod from_area_of_interest(area_id, projection, shape, center, resolution, units=None, **kwargs)

Create an AreaDefinition from center, resolution, and shape.

Parameters:
  • area_id (str) – ID of area

  • projection (dict or str) – Projection parameters as a proj4_dict or proj4_string

  • shape (list) – Number of pixels in the y and x direction (height, width)

  • center (list) – Center of projection (x, y)

  • resolution (list or float) – Size of pixels: (dx, dy). Can be specified with one value if dx == dy

  • units (str, optional) –

    Units that provided arguments should be interpreted as. This can be one of ‘deg’, ‘degrees’, ‘meters’, ‘metres’, and any parameter supported by the cs2cs -lu command. Units are determined in the following priority:

    1. units expressed with each variable through a DataArray’s attrs attribute.

    2. units passed to units

    3. units used in projection

    4. meters

  • description (str, optional) – Description/name of area. Defaults to area_id

  • proj_id (str, optional) – ID of projection

  • nprocs (int, optional) – Number of processor cores to be used

  • lons (numpy array, optional) – Grid lons

  • lats (numpy array, optional) – Grid lats

Returns:

AreaDefinition

Return type:

AreaDefinition

classmethod from_cf(cf_file, variable=None, y=None, x=None)

Create an AreaDefinition object from a netCDF/CF file.

Parameters:
  • nc_file (string or object) – path to a netCDF/CF file, or opened xarray.Dataset object

  • variable (string, optional) – name of the variable to load the AreaDefinition from If variable is None the file will be searched for valid CF area definitions

  • y (string, optional) – name of the variable to use as ‘y’ axis of the CF area definition If y is None an appropriate ‘y’ axis will be deduced from the CF file

  • x (string, optional) – name of the variable to use as ‘x’ axis of the CF area definition If x is None an appropriate ‘x’ axis will be deduced from the CF file

Returns:

AreaDefinition

Return type:

AreaDefinition

classmethod from_circle(area_id, projection, center, radius, shape=None, resolution=None, units=None, **kwargs)

Create an AreaDefinition from center, radius, and shape or from center, radius, and resolution.

Parameters:
  • area_id (str) – ID of area

  • projection (dict or str) – Projection parameters as a proj4_dict or proj4_string

  • center (list) – Center of projection (x, y)

  • radius (list or float) – Length from the center to the edges of the projection (dx, dy)

  • shape (list, optional) – Number of pixels in the y and x direction (height, width)

  • resolution (list or float, optional) – Size of pixels: (dx, dy)

  • units (str, optional) –

    Units that provided arguments should be interpreted as. This can be one of ‘deg’, ‘degrees’, ‘meters’, ‘metres’, and any parameter supported by the cs2cs -lu command. Units are determined in the following priority:

    1. units expressed with each variable through a DataArray’s attrs attribute.

    2. units passed to units

    3. units used in projection

    4. meters

  • description (str, optional) – Description/name of area. Defaults to area_id

  • proj_id (str, optional) – ID of projection

  • nprocs (int, optional) – Number of processor cores to be used

  • lons (numpy array, optional) – Grid lons

  • lats (numpy array, optional) – Grid lats

  • optimize_projection – Whether the projection parameters have to be optimized for a DynamicAreaDefinition.

Returns:

AreaDefinition or DynamicAreaDefinition – If shape or resolution are provided, an AreaDefinition object is returned. Else a DynamicAreaDefinition object is returned

Return type:

AreaDefinition or DynamicAreaDefinition

Notes

  • resolution and radius can be specified with one value if dx == dy

classmethod from_epsg(code, resolution)

Create an AreaDefinition object from an epsg code (string or int) and a resolution.

classmethod from_extent(area_id, projection, shape, area_extent, units=None, **kwargs)

Create an AreaDefinition object from area_extent and shape.

Parameters:
  • area_id (str) – ID of area

  • projection (dict or str) – Projection parameters as a proj4_dict or proj4_string

  • shape (list) – Number of pixels in the y and x direction (height, width)

  • area_extent (list) – Area extent as a list (lower_left_x, lower_left_y, upper_right_x, upper_right_y)

  • units (str, optional) –

    Units that provided arguments should be interpreted as. This can be one of ‘deg’, ‘degrees’, ‘meters’, ‘metres’, and any parameter supported by the cs2cs -lu command. Units are determined in the following priority:

    1. units expressed with each variable through a DataArray’s attrs attribute.

    2. units passed to units

    3. units used in projection

    4. meters

  • description (str, optional) – Description/name of area. Defaults to area_id

  • proj_id (str, optional) – ID of projection

  • nprocs (int, optional) – Number of processor cores to be used

  • lons (numpy array, optional) – Grid lons

  • lats (numpy array, optional) – Grid lats

Returns:

AreaDefinition

Return type:

AreaDefinition

classmethod from_ul_corner(area_id, projection, shape, upper_left_extent, resolution, units=None, **kwargs)

Create an AreaDefinition object from upper_left_extent, resolution, and shape.

Parameters:
  • area_id (str) – ID of area

  • projection (dict or str) – Projection parameters as a proj4_dict or proj4_string

  • shape (list) – Number of pixels in the y and x direction (height, width)

  • upper_left_extent (list) – Upper left corner of upper left pixel (x, y)

  • resolution (list or float) – Size of pixels in meters: (dx, dy). Can be specified with one value if dx == dy

  • units (str, optional) –

    Units that provided arguments should be interpreted as. This can be one of ‘deg’, ‘degrees’, ‘meters’, ‘metres’, and any parameter supported by the cs2cs -lu command. Units are determined in the following priority:

    1. units expressed with each variable through a DataArray’s attrs attribute.

    2. units passed to units

    3. units used in projection

    4. meters

  • description (str, optional) – Description/name of area. Defaults to area_id

  • proj_id (str, optional) – ID of projection

  • nprocs (int, optional) – Number of processor cores to be used

  • lons (numpy array, optional) – Grid lons

  • lats (numpy array, optional) – Grid lats

Returns:

AreaDefinition

Return type:

AreaDefinition

geocentric_resolution(ellps='WGS84', radius=None)

Find best estimate for overall geocentric resolution.

This method is extremely important to the results of KDTree-based resamplers like the nearest neighbor resampling. This is used to determine how far the KDTree should be queried for valid pixels before giving up (radius_of_influence). This method attempts to make a best guess at what geocentric resolution (the units used by the KDTree) represents the majority of an area.

To do this this method will:

  1. Create a vertical mid-line and a horizontal mid-line.

  2. Convert these coordinates to geocentric coordinates.

  3. Compute the distance between points along these lines.

  4. Take the histogram of each set of distances and find the bin with the most points.

  5. Take the average of the edges of that bin.

  6. Return the maximum of the vertical and horizontal bin edge averages.

get_area_slices(area_to_cover, shape_divisible_by=None)

Compute the slice to read based on an area_to_cover.

get_array_coordinates_from_lonlat(lon, lat)

Retrieve the array coordinates (float) for a given lon/lat.

If lon,lat is a tuple of sequences of longitudes and latitudes, a tuple of arrays is returned.

Parameters:
  • lon (array_like) – point or sequence of longitudes

  • lat (array_like) – point or sequence of latitudes

Returns:

the array coordinates (cols/rows)

Return type:

floats or arrays of floats

get_array_coordinates_from_projection_coordinates(xm, ym)

Find the floating-point grid cell index for a specified projection coordinate.

If xm, ym is a tuple of sequences of projection coordinates, a tuple of arrays are returned.

Parameters:
  • xm (array_like) – point or sequence of x-coordinates in meters (map projection)

  • ym (array_like) – point or sequence of y-coordinates in meters (map projection)

Returns:

the array coordinates (cols/rows)

Return type:

floats or arrays of floats

get_array_indices_from_lonlat(lon, lat)

Find the closest integer grid cell index for a given lon/lat.

If lon,lat is a point, a ValueError is raised if it is outside the area domain. If lon,lat is a tuple of sequences of longitudes and latitudes, a tuple of masked arrays are returned. The masked values are the actual row and col indexing the grid cell if the area had been big enough, or the numpy default (999999) if invalid.

Parameters:
  • lon (array_like) – point or sequence of longitudes

  • lat (array_like) – point or sequence of latitudes

Returns:

the array indices (cols/rows)

Return type:

ints or masked arrays of ints

Raises:

ValueError – if the return point is outside the area domain

get_array_indices_from_projection_coordinates(xm, ym)

Find the closest integer grid cell index for a specified projection coordinate.

If xm, ym is a point, a ValueError is raised if it is outside the area domain. If xm, ym is a tuple of sequences of projection coordinates, a tuple of masked arrays are returned.

Parameters:
  • xm (array_like) – point or sequence of x-coordinates in meters (map projection)

  • ym (array_like) – point or sequence of y-coordinates in meters (map projection)

Returns:

the array indices (cols/rows)

Return type:

ints or masked arrays of ints

Raises:

ValueError – if the return point is outside the area domain

get_edge_bbox_in_projection_coordinates(vertices_per_side: int | None = None, frequency: int | None = None)

Return the bounding box in projection coordinates.

get_lonlat(row, col)

Retrieve lon and lat values of single point in area grid.

Parameters:
Returns:

(lon, lat)

Return type:

tuple of floats

get_lonlat_from_array_coordinates(cols, rows)

Get the longitude and latitude from (floating) column and row indices.

If cols, rows is a tuple of sequences of array coordinates, a tuple of arrays is returned.

Parameters:
  • cols (array_like) – the column coordinates

  • rows (array_like) – the row coordinates

Returns:

the longitude, latitude in degrees

Return type:

floats or arrays of floats

get_lonlat_from_projection_coordinates(xm, ym)

Get the lonlat from projection coordinates.

If xm, ym is a tuple of sequences of projection coordinates, a tuple of arrays is returned.

Parameters:
  • xm (array_like) – the x projection coordinates in meters

  • ym (array_like) – the y projection coordinates in meters

Returns:

the longitude, latitude in degrees

Return type:

floats or arrays of floats

get_lonlats(nprocs=None, data_slice=None, cache=False, dtype=None, chunks=None)

Return lon and lat arrays of area.

Note that this historically this method always returns longitude/latitudes on the geodetic (unprojected) model of the Earth used by the Coordinate Reference System (CRS) for this area. However, this is not true for shifted datums. For example, a projection including a PROJ.4 parameter like +pm=180 to shift longitudes/latitudes 180 degrees, will return degrees on the +pm=0 equivalent of the geodetic CRS.

Parameters:
  • nprocs (int, optional) – Number of processor cores to be used. Defaults to the nprocs set when instantiating object

  • data_slice (slice object, optional) – Calculate only coordinates for specified slice

  • cache (bool, optional) – Store the result internally for later reuse. Requires data_slice to be None.

  • dtype (numpy.dtype, optional) – Data type of the returned arrays

  • chunks (int or tuple, optional) – Create dask arrays and use this chunk size

Returns:

(lons, lats) – Grids of area lons and and lats

Return type:

tuple of numpy arrays

get_lonlats_dask(chunks=None, dtype=None)

Get longitudes and latitudes.

get_proj_coords(data_slice=None, dtype=None, chunks=None)

Get projection coordinates of grid.

Parameters:
  • data_slice (slice object, optional) – Calculate only coordinates for specified slice

  • dtype (numpy.dtype, optional) – Data type of the returned arrays

  • chunks (int or tuple, optional) – Create dask arrays and use this chunk size

Returns:

  • (target_x, target_y) (tuple of numpy arrays) – Grids of area x- and y-coordinates in projection units

  • .. versionchanged:: 1.11.0 – Removed ‘cache’ keyword argument and add ‘chunks’ for creating dask arrays.

get_proj_coords_dask(chunks=None, dtype=None)

Get projection coordinates.

get_proj_vectors(dtype=None, chunks=None)

Calculate 1D projection coordinates for the X and Y dimension.

Parameters:
  • dtype (numpy.dtype) – Numpy data type for the returned arrays

  • chunks (int or tuple) – Return dask arrays with the chunk size specified. If this is a tuple then the first element is the Y array’s chunk size and the second is the X array’s chunk size.

Returns:

  • tuple ((X, Y) where X and Y are 1-dimensional numpy arrays)

  • The data type of the returned arrays can be controlled with the

  • dtype keyword argument. If chunks is provided then dask arrays

  • are returned instead.

get_proj_vectors_dask(chunks=None, dtype=None)

Get projection vectors.

get_projection_coordinates_from_array_coordinates(cols, rows)

Get the projection coordinate from the array coordinates.

If cols, rows is a tuple of sequences of array coordinates, a tuple of arrays is returned.

Parameters:
  • cols (array_like) – the column coordinates

  • rows (array_like) – the row coordinates

Returns:

the projection coordinates x, y in meters

Return type:

floats or arrays of floats

get_projection_coordinates_from_lonlat(lon, lat)

Get the projection coordinate from longitudes and latitudes.

If lon,lat is a tuple of sequences of longitudes and latitudes, a tuple of arrays is returned.

Parameters:
  • lon (array_like) – point or sequence of longitudes

  • lat (array_like) – point or sequence of latitudes

Returns:

the projection coordinates x, y in meters

Return type:

floats or arrays of floats

get_xy_from_lonlat(lon, lat)

Retrieve closest x and y coordinates.

Retrieve the closest x and y coordinates (column, row indices) for the specified geolocation (lon,lat) if inside area. If lon,lat is a point a ValueError is raised if the return point is outside the area domain. If lon,lat is a tuple of sequences of longitudes and latitudes, a tuple of masked arrays are returned.

Parameters:
  • lon – point or sequence (list or array) of longitudes

  • lat – point or sequence (list or array) of latitudes

Returns:

tuple of points/arrays

Return type:

(x, y)

get_xy_from_proj_coords(xm, ym)

Find closest grid cell index for a specified projection coordinate.

If xm, ym is a tuple of sequences of projection coordinates, a tuple of masked arrays are returned.

Parameters:
  • xm (list or array) – point or sequence of x-coordinates in meters (map projection)

  • ym (list or array) – point or sequence of y-coordinates in meters (map projection)

Returns:

column and row grid cell indexes as 2 scalars or arrays

Return type:

x, y

Raises:

ValueError – if the return point is outside the area domain

property is_geostationary

Whether this area is in a geostationary satellite projection or not.

lonlat2colrow(lons, lats)

Return image columns and rows for the given lons and lats.

Both scalars and arrays are supported. Same as get_xy_from_lonlat, renamed for convenience.

property name

Return area name.

property outer_boundary_corners

Return the lon,lat of the outer edges of the corner points.

property proj4_string

Return projection definition as Proj.4 string.

property proj_str

Return PROJ projection string.

This is no longer the preferred way of describing CRS information. Switch to the crs or crs_wkt properties for the most flexibility.

property projection_x_coords

Return projection X coordinates.

property projection_y_coords

Return projection Y coordinates.

property resolution

Return area resolution in X and Y direction.

to_cartopy_crs()

Convert projection to cartopy CRS object.

to_odc_geobox()

Convert AreaDefinition to ODC GeoBox.

See: https://odc-geo.readthedocs.io/en/latest/

update_hash(existing_hash: _Hash | None = None) _Hash

Update a hash, or return a new one if needed.

class pyresample.geometry.BaseDefinition(lons=None, lats=None, nprocs=1)

Bases: object

Base class for geometry definitions.

Changed in version 1.8.0: BaseDefinition no longer checks the validity of the provided longitude and latitude coordinates to improve performance. Longitude arrays are expected to be between -180 and 180 degrees, latitude -90 to 90 degrees. Use check_and_wrap() to preprocess your arrays.

boundary(vertices_per_side=None, force_clockwise=False, frequency=None)

Retrieve the AreaBoundary object.

Parameters:
  • vertices_per_side – (formerly frequency) The number of points to provide for each side. By default (None) the full width and height will be provided.

  • force_clockwise – Perform minimal checks and reordering of coordinates to ensure that the returned coordinates follow a clockwise direction. This is important for compatibility with pyresample.spherical.SphPolygon where operations depend on knowing the inside versus the outside of a polygon. These operations assume that coordinates are clockwise. Default is False.

property corners

Return the corners centroids of the current area.

get_area()

Get the area of the convex area defined by the corners of the curren area.

get_area_extent_for_subset(row_LR, col_LR, row_UL, col_UL)

Calculate extent for a subdomain of this area.

Rows are counted from upper left to lower left and columns are counted from upper left to upper right.

Parameters:
  • row_LR (int) – row of the lower right pixel

  • col_LR (int) – col of the lower right pixel

  • row_UL (int) – row of the upper left pixel

  • col_UL (int) – col of the upper left pixel

Returns:

Area extent (LL_x, LL_y, UR_x, UR_y) of the subset

Return type:

area_extent (tuple)

Author:

Ulrich Hamann

get_area_slices(area_to_cover)

Compute the slice to read based on an area_to_cover.

get_bbox_lonlats(vertices_per_side: int | None = None, force_clockwise: bool = True, frequency: int | None = None) tuple

Return the bounding box lons and lats sides.

Parameters:
  • vertices_per_side – The number of points to provide for each side. By default (None) the full width and height will be provided.

  • frequency – Deprecated, use vertices_per_side

  • force_clockwise – Perform minimal checks and reordering of coordinates to ensure that the returned coordinates follow a clockwise direction. This is important for compatibility with pyresample.spherical.SphPolygon where operations depend on knowing the inside versus the outside of a polygon. These operations assume that coordinates are clockwise. Default is True.

Returns:

Two lists of four elements each. The first list is longitude coordinates, the second latitude. Each element is a numpy array representing a specific side of the geometry. The order of the arrays is first row (index 0), last column, last row, and first column. The arrays are sliced (ordered) in a way to ensure that the coordinates follow a clockwise path. In the usual case this results in the coordinates starting in the north-west corner. In the case where the data is oriented with the first pixel (row 0, column 0) in the south-east corner, the coordinates will start in that corner. Other orientations that are detected to follow a counter-clockwise path will be reordered to provide a clockwise path in order to be compatible with other parts of pyresample (ex. pyresample.spherical.SphPolygon).

get_boundary_lonlats()

Return Boundary objects.

get_cartesian_coords(nprocs=None, data_slice=None, cache=False)

Retrieve cartesian coordinates of geometry definition.

Parameters:
  • nprocs (int, optional) – Number of processor cores to be used. Defaults to the nprocs set when instantiating object

  • data_slice (slice object, optional) – Calculate only cartesian coordnates for the defined slice

  • cache (bool, optional) – Store result the result. Requires data_slice to be None

Returns:

cartesian_coords

Return type:

numpy array

get_edge_lonlats(vertices_per_side=None, frequency=None)

Get the concatenated boundary of the current swath.

get_lonlat(row, col)

Retrieve lon and lat of single pixel.

Parameters:
Returns:

(lon, lat)

Return type:

tuple of floats

get_lonlats(data_slice=None, chunks=None, **kwargs)

Get longitude and latitude arrays representing this geometry.

Returns:

(lon, lat) – If chunks is provided then the arrays will be dask arrays with the provided chunk size. If chunks is not provided then the returned arrays are the same as the internal data types of this geometry object (numpy or dask).

Return type:

tuple of numpy arrays

get_lonlats_dask(chunks=None)

Get the lon lats as a single dask array.

intersection(other)

Return the corners of the intersection polygon of the current area with other.

Parameters:

other (object) – Instance of subclass of BaseDefinition

Returns:

(corner1, corner2, corner3, corner4)

Return type:

tuple of points

property is_geostationary

Whether this geometry is in a geostationary satellite projection or not.

overlap_rate(other)

Get how much the current area overlaps an other area.

Parameters:

other (object) – Instance of subclass of BaseDefinition

Returns:

overlap_rate

Return type:

float

overlaps(other)

Test if the current area overlaps the other area.

This is based solely on the corners of areas, assuming the boundaries to be great circles.

Parameters:

other (object) – Instance of subclass of BaseDefinition

Returns:

overlaps

Return type:

bool

update_hash(existing_hash: _Hash | None = None) _Hash

Update the hash.

class pyresample.geometry.CoordinateDefinition(lons, lats, nprocs=1)

Bases: BaseDefinition

Base class for geometry definitions defined by lons and lats only.

append(other)

Append another coordinate definition to existing one.

concatenate(other)

Concatenate coordinate definitions.

geocentric_resolution(ellps='WGS84', radius=None, nadir_factor=2)

Calculate maximum geocentric pixel resolution.

If lons is a xarray.DataArray object with a resolution attribute, this will be used instead of loading the longitude and latitude data. In this case the resolution attribute is assumed to mean the nadir resolution of a swath and will be multiplied by the nadir_factor to adjust for increases in the spatial resolution towards the limb of the swath.

Parameters:
  • ellps (str) – PROJ Ellipsoid for the Cartographic projection used as the target geocentric coordinate reference system. Default: ‘WGS84’. Ignored if radius is provided.

  • radius (float) – Spherical radius of the Earth to use instead of the definitions in ellps.

  • nadir_factor (int) – Number to multiply the nadir resolution attribute by to reflect pixel size on the limb of the swath.

Returns: Estimated maximum pixel size in meters on a geocentric

coordinate system (X, Y, Z) representing the Earth.

Raises: RuntimeError if a simple search for valid longitude/latitude

data points found no valid data points.

exception pyresample.geometry.DimensionError

Bases: ValueError

Wrap ValueError.

class pyresample.geometry.DynamicAreaDefinition(area_id=None, description=None, projection=None, width=None, height=None, area_extent=None, resolution=None, optimize_projection=False)

Bases: object

An AreaDefintion containing just a subset of the needed parameters.

The purpose of this class is to be able to adapt the area extent and shape of the area to a given set of longitudes and latitudes, such that e.g. polar satellite granules can be resampled optimally to a given projection.

Note that if the provided projection is geographic (lon/lat degrees) and the provided longitude and latitude data crosses the anti-meridian (-180/180), the resulting area will be the smallest possible in order to contain that data and avoid a large area spanning from -180 to 180 longitude. This means the resulting AreaDefinition will have a right-most X extent greater than 180 degrees. This does not apply to data crossing the north or south pole as there is no “smallest” area in this case.

area_id

The name of the area.

description

The description of the area.

projection

The dictionary or string or CRS object of projection parameters. Doesn’t have to be complete. If not complete, proj_info must be provided to freeze to “fill in” any missing parameters.

width

x dimension in number of pixels, aka number of grid columns

height

y dimension in number of pixels, aka number of grid rows

shape

Corresponding array shape as (height, width)

area_extent

The area extent of the area.

resolution

Resolution of the resulting area as (pixel_size_x, pixel_size_y) or a scalar if pixel_size_x == pixel_size_y.

optimize_projection

Whether the projection parameters have to be optimized.

compute_domain(corners: Sequence, resolution: float | tuple[float, float] | None = None, shape: tuple[int, int] | None = None, projection: CRS | dict | str | int | None = None)

Compute shape and area_extent from corners and [shape or resolution] info.

Parameters:
  • corners – 4-element sequence representing the outer corners of the region. Note that corners represents the center of pixels, while area_extent represents the edge of pixels. The four values are (xmin_corner, ymin_corner, xmax_corner, ymax_corner). If the x corners are None then the full extent (area of use) of the projection will be used. When needed, area of use is taken from the PROJ library or in the case of a geographic lon/lat projection -180/180 is used. A RuntimeError is raised if the area of use is needed (when x corners are None) and area of use can’t be determined.

  • resolution – Spatial resolution in projection units (typically meters or degrees). If not specified then shape must be provided. If a scalar then it is treated as the x and y resolution. If a tuple then x resolution is the first element, y is the second.

  • shape – Number of pixels in the area as a 2-element tuple. The first is number of rows, the second number of columns.

  • projection – PROJ.4 definition string, dictionary, integer EPSG code, or pyproj CRS object.

Note that shape is (rows, columns) and resolution is (x_size, y_size); the dimensions are flipped.

freeze(lonslats=None, resolution=None, shape=None, proj_info=None, antimeridian_mode=None)

Create an AreaDefinition from this area with help of some extra info.

Parameters:
  • lonlats (SwathDefinition or tuple) – The geographical coordinates to contain in the resulting area. A tuple should be (lons, lats).

  • resolution – the resolution of the resulting area.

  • shape – the shape of the resulting area.

  • proj_info – complementing parameters to the projection info.

  • antimeridian_mode

    How to handle lon/lat data crossing the anti-meridian of the projection. This currently only affects lon/lat geographic projections and data cases not covering the north or south pole. The possible options are:

    • ”modify_extents”: Set the X bounds to the edges of the data, but

      add 360 to the right-most bound. This has the effect of making the area coordinates continuous from the left side to the right side. However, this means that some coordinates will be outside the coordinate space of the projection. Although most PROJ and pyresample functionality can handle this there may be some edge cases.

    • ”modify_crs”: Change the prime meridian of the projection

      from 0 degrees longitude to 180 degrees longitude. This has the effect of putting the data on a continuous coordinate system. However, this means that comparing data resampled to this resulting area and an area not over the anti-meridian would be more difficult.

    • ”global_extents”: Ignore the bounds of the data and use -180/180

      degrees as the west and east bounds of the data. This will generate a large output area, but with the benefit of keeping the data on the original projection. Note that some resampling methods may produce artifacts when resampling on the edge of the area (the anti-meridian).

  • created (Shape parameters are ignored if the instance is) –

  • True. (with the optimize_projection flag set to) –

property pixel_size_x

Pixel width in projection units.

property pixel_size_y

Pixel height in projection units.

class pyresample.geometry.GridDefinition(lons, lats, nprocs=1)

Bases: CoordinateDefinition

Grid defined by lons and lats.

Parameters:
  • lons (numpy array) –

  • lats (numpy array) –

  • nprocs (int, optional) – Number of processor cores to be used for calculations.

shape

Grid shape as (rows, cols)

Type:

tuple

size

Number of elements in grid

Type:

int

lons

Grid lons

Type:

object

lats

Grid lats

Type:

object

cartesian_coords

Grid cartesian coordinates

Type:

object

exception pyresample.geometry.IncompatibleAreas

Bases: ValueError

Error when the areas to combine are not compatible.

exception pyresample.geometry.InvalidArea

Bases: ValueError

Error to be raised when an area is invalid for a given purpose.

class pyresample.geometry.StackedAreaDefinition(*definitions, **kwargs)

Bases: _ProjectionDefinition

Definition based on muliple vertically stacked AreaDefinitions.

append(definition)

Append another definition to the area.

get_lonlats(nprocs=None, data_slice=None, cache=False, dtype=None, chunks=None)

Return lon and lat arrays of the area.

get_lonlats_dask(chunks=None, dtype=None)

Return lon and lat dask arrays of the area.

property height

Return height of the area definition.

property proj4_string

Return projection definition as Proj.4 string.

property proj_str

Return projection definition as Proj.4 string.

squeeze()

Generate a single AreaDefinition if possible.

update_hash(the_hash=None)

Update the hash.

property width

Return width of the area definition.

class pyresample.geometry.SwathDefinition(lons, lats, nprocs=1, crs=None)

Bases: CoordinateDefinition

Swath defined by lons and lats.

Parameters:
  • lons (numpy array) –

  • lats (numpy array) –

  • nprocs (int, optional) – Number of processor cores to be used for calculations.

  • crs (pyproj.CRS,) – The CRS to use. longlat on WGS84 by default.

shape

Swath shape

Type:

tuple

size

Number of elements in swath

Type:

int

ndims

Swath dimensions

Type:

int

lons

Swath lons

Type:

object

lats

Swath lats

Type:

object

cartesian_coords

Swath cartesian coordinates

Type:

object

aggregate(**dims)

Aggregate the current swath definition by averaging.

For example, averaging over 2x2 windows: sd.aggregate(x=2, y=2)

compute_bb_proj_params(proj_dict)

Compute BB projection parameters.

compute_optimal_bb_area(proj_dict=None, resolution=None)

Compute the “best” bounding box area for this swath with proj_dict.

By default, the projection is Oblique Mercator (omerc in proj.4), in which case the right projection angle alpha is computed from the swath centerline. For other projections, only the appropriate center of projection and area extents are computed.

The height and width are computed so that the resolution is approximately the same across dimensions.

copy()

Copy the current swath.

pyresample.geometry.combine_area_extents_vertical(area1, area2)

Combine the area extents of areas 1 and 2.

pyresample.geometry.concatenate_area_defs(area1, area2, axis=0)

Append area2 to area1 and return the results.

pyresample.geometry.daskify_2in_2out(func)

Daskify the coordinate conversion functions.

pyresample.geometry.enclose_areas(*areas, area_id='joint-area')

Return the smallest areadefinition enclosing one or more others.

From one or more AreaDefinition objects (most usefully at least two), which shall differ only in extent, calculate the smallest AreaDefinition that encloses all. Touches only the area_extent; projection and units must be identical in all input areas and will be unchanged in the resulting area. When the input areas \(i=1..n\) have extent \((a_i, b_i, c_i, d_i)\), the resulting area will have extent \((\\min_i{a_i}, \\min_i{b_i}, \\max_i{c_i}, \\max_i{d_i})\).

Parameters:
  • *areas (AreaDefinition) – AreaDefinition objects to enclose.

  • area_id (Optional[str]) – Name of joint area, defaults to “joint-area”.

pyresample.geometry.get_array_hashable(arr)

Compute a hashable form of the array arr.

Works with numpy arrays, dask.array.Array, and xarray.DataArray.

pyresample.geometry.get_full_geostationary_bounding_box_in_proj_coords(geos_area, nb_points=50)

Get the valid boundary geos projection coordinates of the full disk.

Parameters:
  • geos_area – Geostationary area definition to get the bounding box for.

  • nb_points – Number of points on the polygon

pyresample.geometry.get_geostationary_angle_extent(geos_area)

Get the max earth (vs space) viewing angles in x and y.

pyresample.geometry.get_geostationary_bounding_box(geos_area, nb_points=50)

Get the bbox in lon/lats of the valid pixels inside geos_area.

Parameters:
  • geos_area – Geostationary area definition to get the bounding box for.

  • nb_points – Number of points on the polygon

pyresample.geometry.get_geostationary_bounding_box_in_lonlats(geos_area, nb_points=50)

Get the bbox in lon/lats of the valid pixels inside geos_area.

Parameters:
  • geos_area – Geostationary area definition to get the bounding box for.

  • nb_points – Number of points on the polygon

pyresample.geometry.get_geostationary_bounding_box_in_proj_coords(geos_area, nb_points=50)

Get the bbox in geos projection coordinates of the valid pixels inside geos_area.

Parameters:
  • geos_area – Geostationary area definition to get the bounding box for.

  • nb_points – Number of points on the polygon.

pyresample.geometry.masked_ints(func)

Return masked integer arrays when returning array indices.

pyresample.geometry.ordered_dump(data, stream=None, Dumper=<class 'yaml.dumper.Dumper'>, **kwds)

Dump the data to YAML in ordered fashion.

pyresample.geometry.preserve_scalars(func)

Preserve scalars through the coordinate conversion functions.

pyresample.grid module

Resample image from one projection to another using nearest neighbour method.

pyresample.grid.get_image_from_linesample(row_indices, col_indices, source_image, fill_value=0)

Sample from image based on index arrays.

Parameters:
  • row_indices (numpy array) – Row indices. Dimensions must match col_indices

  • col_indices (numpy array) – Col indices. Dimensions must match row_indices

  • source_image (numpy array) – Source image

  • fill_value (int or None, optional) – Set undetermined pixels to this value. If fill_value is None a masked array is returned with undetermined pixels masked

Returns:

image_data – Resampled image

Return type:

numpy array

pyresample.grid.get_image_from_lonlats(lons, lats, source_area_def, source_image_data, fill_value=0, nprocs=1)

Sample from image based on lon lat arrays using nearest neighbour method in cartesian projection coordinates.

Parameters:
  • lons (numpy array) – Lons. Dimensions must match lats

  • lats (numpy array) – Lats. Dimensions must match lons

  • source_area_def (object) – Source definition as AreaDefinition object

  • source_image_data (numpy array) – Source image data

  • fill_value (int or None, optional) – Set undetermined pixels to this value. If fill_value is None a masked array is returned with undetermined pixels masked

  • nprocs (int, optional) – Number of processor cores to be used

Returns:

image_data – Resampled image data

Return type:

numpy array

pyresample.grid.get_linesample(lons, lats, source_area_def, nprocs=1)

Return index row and col arrays for resampling.

Parameters:
  • lons (numpy array) – Lons. Dimensions must match lats

  • lats (numpy array) – Lats. Dimensions must match lons

  • source_area_def (object) – Source definition as AreaDefinition object

  • nprocs (int, optional) – Number of processor cores to be used

Returns:

(row_indices, col_indices) – Arrays for resampling area by array indexing

Return type:

tuple of numpy arrays

pyresample.grid.get_resampled_image(target_area_def, source_area_def, source_image_data, fill_value=0, nprocs=1, segments=None)

Resample image using nearest neighbour method in cartesian projection coordinate systems.

Parameters:
  • target_area_def (object) – Target definition as AreaDefinition object

  • source_area_def (object) – Source definition as AreaDefinition object

  • source_image_data (numpy array) – Source image data

  • fill_value ({int, None} optional) – Set undetermined pixels to this value. If fill_value is None a masked array is returned with undetermined pixels masked

  • nprocs (int, optional) – Number of processor cores to be used

  • segments ({int, None} optional) – Number of segments to use when resampling. If set to None an estimate will be calculated.

Returns:

image_data – Resampled image data

Return type:

numpy array

pyresample.image module

Handles resampling of images with assigned geometry definitions.

class pyresample.image.ImageContainer(image_data, geo_def, fill_value=0, nprocs=1)

Bases: object

Holds image with geometry definition. Allows indexing with linesample arrays.

Parameters:
  • image_data (numpy array) – Image data

  • geo_def (object) – Geometry definition

  • fill_value (int or None, optional) – Set undetermined pixels to this value. If fill_value is None a masked array is returned with undetermined pixels masked

  • nprocs (int, optional) – Number of processor cores to be used

image_data

Image data

Type:

numpy array

geo_def

Geometry definition

Type:

object

fill_value

Resample result fill value

Type:

int or None

nprocs

Number of processor cores to be used for geometry operations

Type:

int

get_array_from_linesample(row_indices, col_indices)

Get array sampled from image based on index arrays.

Parameters:
  • row_indices (numpy array) – Row indices. Dimensions must match col_indices

  • col_indices (numpy array) – Col indices. Dimensions must match row_indices

Returns:

image_data – Resampled image data

Return type:

numpy_array

get_array_from_neighbour_info(*args, **kwargs)

Resample from preprocessed data.

resample(target_geo_def)

Resample data to target definition.

class pyresample.image.ImageContainerBilinear(image_data, geo_def, radius_of_influence, epsilon=0, fill_value=0, reduce_data=False, nprocs=1, segments=None, neighbours=32)

Bases: ImageContainer

Holds image with geometry definition. Allows bilinear to new geometry definition.

Parameters:
  • image_data (numpy array) – Image data

  • geo_def (object) – Geometry definition

  • radius_of_influence (float) – Cut off distance in meters

  • epsilon (float, optional) – Allowed uncertainty in meters. Increasing uncertainty reduces execution time

  • fill_value (int or None, optional) – Set undetermined pixels to this value. If fill_value is None a masked array is returned with undetermined pixels masked

  • reduce_data (bool, optional) – Perform coarse data reduction before resampling in order to reduce execution time

  • nprocs (int, optional) – Number of processor cores to be used for geometry operations

  • segments (int or None) – Number of segments to use when resampling. If set to None an estimate will be calculated

image_data

Image data

Type:

numpy array

geo_def

Geometry definition

Type:

object

radius_of_influence

Cut off distance in meters

Type:

float

epsilon

Allowed uncertainty in meters

Type:

float

fill_value

Resample result fill value

Type:

int or None

reduce_data

Perform coarse data reduction before resampling

Type:

bool

nprocs

Number of processor cores to be used

Type:

int

segments

Number of segments to use when resampling

Type:

int or None

resample(target_geo_def)

Resamples image to area definition using bilinear approach.

Parameters:

target_geo_def (object) – Target geometry definition

Returns:

image_container – ImageContainerBilinear object of resampled geometry

Return type:

object

class pyresample.image.ImageContainerNearest(image_data, geo_def, radius_of_influence, epsilon=0, fill_value=0, reduce_data=True, nprocs=1, segments=None)

Bases: ImageContainer

Holds image with geometry definition. Allows nearest neighbour to new geometry definition.

Parameters:
  • image_data (numpy array) – Image data

  • geo_def (object) – Geometry definition

  • radius_of_influence (float) – Cut off distance in meters

  • epsilon (float, optional) – Allowed uncertainty in meters. Increasing uncertainty reduces execution time

  • fill_value (int or None, optional) – Set undetermined pixels to this value. If fill_value is None a masked array is returned with undetermined pixels masked

  • reduce_data (bool, optional) – Perform coarse data reduction before resampling in order to reduce execution time

  • nprocs (int, optional) – Number of processor cores to be used for geometry operations

  • segments (int or None) – Number of segments to use when resampling. If set to None an estimate will be calculated

image_data

Image data

Type:

numpy array

geo_def

Geometry definition

Type:

object

radius_of_influence

Cut off distance in meters

Type:

float

epsilon

Allowed uncertainty in meters

Type:

float

fill_value

Resample result fill value

Type:

int or None

reduce_data

Perform coarse data reduction before resampling

Type:

bool

nprocs

Number of processor cores to be used

Type:

int

segments

Number of segments to use when resampling

Type:

int or None

resample(target_geo_def)

Resample image to area definition using nearest neighbour approach.

Parameters:

target_geo_def (object) – Target geometry definition

Returns:

image_container – ImageContainerNearest object of resampled geometry

Return type:

object

class pyresample.image.ImageContainerQuick(image_data, geo_def, fill_value=0, nprocs=1, segments=None)

Bases: ImageContainer

Holds image with area definition and allows quick resampling within area.

Parameters:
  • image_data (numpy array) – Image data

  • geo_def (object) – Area definition as AreaDefinition object

  • fill_value (int or None, optional) – Set undetermined pixels to this value. If fill_value is None a masked array is returned with undetermined pixels masked

  • nprocs (int, optional) – Number of processor cores to be used for geometry operations

  • segments (int or None) – Number of segments to use when resampling. If set to None an estimate will be calculated

image_data

Image data

Type:

numpy array

geo_def

Area definition as AreaDefinition object

Type:

object

fill_value

Resample result fill value If fill_value is None a masked array is returned with undetermined pixels masked

Type:

int or None

nprocs

Number of processor cores to be used

Type:

int

segments

Number of segments to use when resampling

Type:

int or None

resample(target_area_def)

Resample image to area definition using nearest neighbour approach in projection coordinates.

Parameters:

target_area_def (object) – Target area definition as AreaDefinition object

Returns:

image_container – ImageContainerQuick object of resampled area

Return type:

object

pyresample.kd_tree module

Handle reprojection of geolocated data.

Several types of resampling are supported.

exception pyresample.kd_tree.EmptyResult

Bases: ValueError

No valid data is produced.

class pyresample.kd_tree.XArrayResamplerNN(source_geo_def, target_geo_def, radius_of_influence=None, neighbours=1, epsilon=0)

Bases: object

Resampler for Xarray DataArray objects with the nearest neighbor algorithm.

get_neighbour_info(mask=None)

Return neighbour info.

Returns:

  • (valid_input_index, valid_output_index,

  • index_array, distance_array) (tuple of numpy arrays) – Neighbour resampling info

get_sample_from_neighbour_info(data, fill_value=nan)

Get the pixels matching the target area.

This method should work for any dimensionality of the provided data array as long as the geolocation dimensions match in size and name in data.dims. Where source area definition are AreaDefinition objects the corresponding dimensions in the data should be ('y', 'x').

This method also attempts to preserve chunk sizes of dask arrays, but does require loading/sharing the fully computed source data before it can actually compute the values to write to the destination array. This can result in large memory usage for large source data arrays, but is a necessary evil until fancier indexing is supported by dask and/or pykdtree.

Parameters:
  • data (xarray.DataArray) – Source data pixels to sample

  • fill_value (float) – Output fill value when no source data is near the target pixel. When omitted, if the input data is an integer array then the maximum value for that integer type is used, but otherwise, NaN is used and can be detected in the result with res.isnull().

Returns:

The resampled array. The dtype of the array will

be the same as the input data. Pixels with no matching data from the input array will be filled (see the fill_value parameter description above).

Return type:

dask.array.Array

query_resample_kdtree(resample_kdtree, tlons, tlats, valid_oi, mask)

Query kd-tree on slice of target coordinates.

pyresample.kd_tree.get_neighbour_info(source_geo_def, target_geo_def, radius_of_influence, neighbours=8, epsilon=0, reduce_data=True, nprocs=1, segments=None)

Return neighbour info.

Parameters:
  • source_geo_def (object) – Geometry definition of source

  • target_geo_def (object) – Geometry definition of target

  • radius_of_influence (float) – Cut off distance in meters

  • neighbours (int, optional) – The number of neigbours to consider for each grid point

  • epsilon (float, optional) – Allowed uncertainty in meters. Increasing uncertainty reduces execution time

  • reduce_data (bool, optional) – Perform initial coarse reduction of source dataset in order to reduce execution time

  • nprocs (int, optional) – Number of processor cores to be used

  • segments (int or None) – Number of segments to use when resampling. If set to None an estimate will be calculated

Returns:

  • (valid_input_index, valid_output_index,

  • index_array, distance_array) (tuple of numpy arrays) – Neighbour resampling info

pyresample.kd_tree.get_sample_from_neighbour_info(resample_type, output_shape, data, valid_input_index, valid_output_index, index_array, distance_array=None, weight_funcs=None, fill_value=0, with_uncert=False)

Resamples swath based on neighbour info.

Parameters:
  • resample_type ({'nn', 'custom'}) – ‘nn’: Use nearest neighbour resampling ‘custom’: Resample based on weight_funcs

  • output_shape ((int, int)) – Shape of output as (rows, cols)

  • data (numpy array) – Source data

  • valid_input_index (numpy array) – valid_input_index from get_neighbour_info

  • valid_output_index (numpy array) – valid_output_index from get_neighbour_info

  • index_array (numpy array) – index_array from get_neighbour_info

  • distance_array (numpy array, optional) – distance_array from get_neighbour_info Not needed for ‘nn’ resample type

  • weight_funcs (list of function objects or function object, optional) – List of weight functions f(dist) to use for the weighting of each channel 1 to k. If only one channel is resampled weight_funcs is a single function object. Must be supplied when using ‘custom’ resample type

  • fill_value (int, float, numpy floating, numpy integer or None, optional) – Set undetermined pixels to this value. If fill_value is None a masked array is returned with undetermined pixels masked

Returns:

result – Source data resampled to target geometry

Return type:

numpy array

pyresample.kd_tree.resample_custom(source_geo_def, data, target_geo_def, radius_of_influence, weight_funcs, neighbours=8, epsilon=0, fill_value=0, reduce_data=True, nprocs=1, segments=None, with_uncert=False)

Resamples data using kd-tree custom radial weighting neighbour approach.

Parameters:
  • source_geo_def (object) – Geometry definition of source

  • data (numpy array) – Array of single channel data points or (source_geo_def.shape, k) array of k channels of datapoints

  • target_geo_def (object) – Geometry definition of target

  • radius_of_influence (float) – Cut off distance in meters

  • weight_funcs (list of function objects or function object) – List of weight functions f(dist) to use for the weighting of each channel 1 to k. If only one channel is resampled weight_funcs is a single function object.

  • neighbours (int, optional) – The number of neigbours to consider for each grid point

  • epsilon (float, optional) – Allowed uncertainty in meters. Increasing uncertainty reduces execution time

  • fill_value ({int, None}, optional) – Set undetermined pixels to this value. If fill_value is None a masked array is returned with undetermined pixels masked

  • reduce_data (bool, optional) – Perform initial coarse reduction of source dataset in order to reduce execution time

  • nprocs (int, optional) – Number of processor cores to be used

  • segments ({int, None}) – Number of segments to use when resampling. If set to None an estimate will be calculated

Returns:

  • data (numpy array (default)) – Source data resampled to target geometry

  • data, stddev, counts (numpy array, numpy array, numpy array (if with_uncert == True)) – Source data resampled to target geometry. Weighted standard devaition for all pixels having more than one source value Counts of number of source values used in weighting per pixel

pyresample.kd_tree.resample_gauss(source_geo_def, data, target_geo_def, radius_of_influence, sigmas, neighbours=8, epsilon=0, fill_value=0, reduce_data=True, nprocs=1, segments=None, with_uncert=False)

Resamples data using kd-tree gaussian weighting neighbour approach.

Parameters:
  • source_geo_def (object) – Geometry definition of source

  • data (numpy array) – Array of single channel data points or (source_geo_def.shape, k) array of k channels of datapoints

  • target_geo_def (object) – Geometry definition of target

  • radius_of_influence (float) – Cut off distance in meters

  • sigmas (list of floats or float) – List of sigmas to use for the gauss weighting of each channel 1 to k, w_k = exp(-dist^2/sigma_k^2). If only one channel is resampled sigmas is a single float value.

  • neighbours (int, optional) – The number of neigbours to consider for each grid point

  • epsilon (float, optional) – Allowed uncertainty in meters. Increasing uncertainty reduces execution time

  • fill_value ({int, None}, optional) – Set undetermined pixels to this value. If fill_value is None a masked array is returned with undetermined pixels masked

  • reduce_data (bool, optional) – Perform initial coarse reduction of source dataset in order to reduce execution time

  • nprocs (int, optional) – Number of processor cores to be used

  • segments (int or None) – Number of segments to use when resampling. If set to None an estimate will be calculated

  • with_uncert (bool, optional) – Calculate uncertainty estimates

Returns:

  • data (numpy array (default)) – Source data resampled to target geometry

  • data, stddev, counts (numpy array, numpy array, numpy array (if with_uncert == True)) – Source data resampled to target geometry. Weighted standard devaition for all pixels having more than one source value Counts of number of source values used in weighting per pixel

pyresample.kd_tree.resample_nearest(source_geo_def, data, target_geo_def, radius_of_influence, epsilon=0, fill_value=0, reduce_data=True, nprocs=1, segments=None)

Resamples data using kd-tree nearest neighbour approach.

Parameters:
  • source_geo_def (object) – Geometry definition of source

  • data (numpy array) – 1d array of single channel data points or (source_size, k) array of k channels of datapoints

  • target_geo_def (object) – Geometry definition of target

  • radius_of_influence (float) – Cut off distance in meters

  • epsilon (float, optional) – Allowed uncertainty in meters. Increasing uncertainty reduces execution time

  • fill_value (int, float, numpy floating, numpy integer or None, optional) – Set undetermined pixels to this value. If fill_value is None a masked array is returned with undetermined pixels masked

  • reduce_data (bool, optional) – Perform initial coarse reduction of source dataset in order to reduce execution time

  • nprocs (int, optional) – Number of processor cores to be used

  • segments (int or None) – Number of segments to use when resampling. If set to None an estimate will be calculated

Returns:

data – Source data resampled to target geometry

Return type:

numpy array

pyresample.plot module

Utility functions for quick and easy display.

pyresample.plot.area_def2basemap(area_def, **kwargs)

Get Basemap object from an AreaDefinition object.

Parameters:
  • area_def (object) – geometry.AreaDefinition object

  • **kwargs (Keyword arguments) – Additional initialization arguments for Basemap

Returns:

bmap

Return type:

Basemap object

pyresample.plot.ellps2axis(ellps_name)

Get semi-major and semi-minor axis from ellipsis definition.

Parameters:

ellps_name (str) – Standard name of ellipsis

Returns:

(a, b)

Return type:

semi-major and semi-minor axis

pyresample.plot.save_quicklook(filename, area_def, data, vmin=None, vmax=None, label='Variable (units)', num_meridians=45, num_parallels=10, coast_res='110m', cmap='RdBu_r')

Display and save default quicklook plot.

Parameters:
  • filename (str) – path to output file

  • area_def (object) – geometry.AreaDefinition object

  • data (numpy array | numpy masked array) – 2D array matching area_def. Use masked array for transparent values

  • vmin (float, optional) – Min value for luminescence scaling

  • vmax (float, optional) – Max value for luminescence scaling

  • label (str, optional) – Label for data

  • num_meridians (int, optional) – Number of meridians to plot on the globe

  • num_parallels (int, optional) – Number of parallels to plot on the globe

  • coast_res ({'c', 'l', 'i', 'h', 'f'}, optional) – Resolution of coastlines

pyresample.plot.show_quicklook(area_def, data, vmin=None, vmax=None, label='Variable (units)', num_meridians=45, num_parallels=10, coast_res='110m', cmap='RdBu_r')

Display default quicklook plot.

Parameters:
  • area_def (object) – geometry.AreaDefinition object

  • data (numpy array | numpy masked array) – 2D array matching area_def. Use masked array for transparent values

  • vmin (float, optional) – Min value for luminescence scaling

  • vmax (float, optional) – Max value for luminescence scaling

  • label (str, optional) – Label for data

  • num_meridians (int, optional) – Number of meridians to plot on the globe

  • num_parallels (int, optional) – Number of parallels to plot on the globe

  • coast_res ({'c', 'l', 'i', 'h', 'f'}, optional) – Resolution of coastlines

Returns:

bmap

Return type:

Basemap object

pyresample.resampler module

Base resampler class made for subclassing.

class pyresample.resampler.BaseResampler(source_geo_def: SwathDefinition | AreaDefinition, target_geo_def: CoordinateDefinition | AreaDefinition)

Bases: object

Base abstract resampler class.

compute(data, **kwargs)

Do the actual resampling.

This must be implemented by subclasses.

get_hash(source_geo_def=None, target_geo_def=None, **kwargs)

Get hash for the current resample with the given kwargs.

precompute(**kwargs)

Do the precomputation.

This is an optional step if the subclass wants to implement more complex features like caching or can share some calculations between multiple datasets to be processed.

resample(data, cache_dir=None, mask_area=None, **kwargs)

Resample data by calling precompute and compute methods.

Only certain resampling classes may use cache_dir and the mask provided when mask_area is True. The return value of calling the precompute method is passed as the cache_id keyword argument of the compute method, but may not be used directly for caching. It is up to the individual resampler subclasses to determine how this is used.

Parameters:
  • data (xarray.DataArray) – Data to be resampled

  • cache_dir (str) – directory to cache precomputed results (default False, optional)

  • mask_area (bool) – Mask geolocation data where data values are invalid. This should be used when data values may affect what neighbors are considered valid.

  • kwargs – Keyword arguments to pass to both the precompute and compute stages of the resampler.

Returns (xarray.DataArray): Data resampled to the target area

pyresample.resampler.crop_data_around_area(source_geo_def, src_arrays, target_geo_def)

Crop the data around the provided area.

pyresample.resampler.crop_source_area(source_geo_def, target_geo_def)

Crop a source area around the provided target area.

pyresample.resampler.resample_blocks(func, src_area, src_arrays, dst_area, dst_arrays=(), chunk_size=None, dtype=None, name=None, fill_value=None, **kwargs)

Resample dask arrays blockwise.

Resample_blocks applies a function blockwise to transform data from a source area domain to a destination area domain.

Parameters:
  • func – A callable to apply on the input data. This function is passed a block of src_arrays, dst_arrays in that order, followed by the kwargs, which include the fill_value. If the callable accepts a block_info keyword argument, block information is passed to it. Block information provides the source area, destination area, position of source and destination blocks relative to respectively src_area and dst_area.

  • src_area – a source geo definition.

  • dst_area – a destination geo definition. If the same as the source definition, a ValueError is raised.

  • src_arrays – data to use. When split into smaller bit to pass to func, they are split across the x and y dimensions, but not across the other dimensions, so all the dimensions of the smaller arrays will be using only one chunk!

  • dst_arrays – arrays to use that are already in dst_area space. If the array has more than 2 dimensions, the last two are expected to be y, x.

  • chunk_size – the chunks size(s) to use in the dst_area space. This has to be provided since it is not guaranteed that we can get this information from the other arguments. Moreover, this needs to be an iterable of k elements if the resulting array of func is to have a different number of dimensions (k) than the input array.

  • dtype – the dtype the resulting array is going to have. Has to be provided.

  • name – Name prefix of the dask tasks to be generated

  • fill_value – Desired value for any invalid values in the output array

  • kwargs – any other keyword arguments that will be passed on to func.

Principle of operations:

Resample_blocks works by iterating over chunks on the dst_area domain. For each chunk, the corresponding slice of the src_area domain is computed and the input src_arrays are cut accordingly to pass to func. To know more about how the slicing is performed, refer to the :class:Slicer class and subclasses.

Examples

To generate indices from the gradient resampler, you can apply the corresponding function with no input. Note how we provide the chunk sizes knowing that the result array with have 2 elements along a third dimension.

>>> indices_xy = resample_blocks(gradient_resampler_indices, source_geo_def, [], target_geo_def,
...                              chunk_size=(2, "auto", "auto"), dtype=float)

From these indices, to resample an array using bilinear interpolation:

>>>  resampled = resample_blocks(block_bilinear_interpolator, source_geo_def, [src_array], target_geo_def,
...                              dst_arrays=[indices_xy],
...                              chunk_size=("auto", "auto"), dtype=src_array.dtype)

pyresample.slicer module

Area and Swath Slicers.

class pyresample.slicer.AreaSlicer(area_to_crop, area_to_contain)

Bases: Slicer

A Slicer for cropping AreaDefinitions.

get_polygon_to_contain()

Get the shapely Polygon corresponding to area_to_contain in projection coordinates of area_to_crop.

get_slices_from_polygon(poly_to_contain)

Get the slices based on the polygon.

class pyresample.slicer.Slicer(area_to_crop, area_to_contain)

Bases: ABC

Abstract Slicer.

Provided an Area-to-crop and an Area-to-contain, a Slicer provides methods to find slices that enclose area-to-contain inside area-to-crop.

Example

For slicing a full-disk MSG area using a polar-stereographic area over Germany:

>>> from pyresample import slicer
>>> from satpy.resample import get_area_def
>>> msg_area = get_area_def("msg_seviri_fes_3km")
>>> germ_area = get_area_def("germ")
>>> slc = slicer.create_slicer(msg_area, germ_area)
>>> slc.get_slices()
(slice(1900, 2242, None), slice(233, 423, None))
abstract get_polygon_to_contain()

Get the shapely Polygon corresponding to area_to_contain.

get_slices()

Get the slices to crop area_to_crop enclosing area_to_contain.

abstract get_slices_from_polygon(poly)

Get the slices based on the polygon.

class pyresample.slicer.SwathSlicer(area_to_crop, area_to_contain)

Bases: Slicer

A Slicer for cropping SwathDefinitions.

get_polygon_to_contain()

Get the shapely Polygon corresponding to area_to_contain in lon/lat coordinates.

get_slices_from_polygon(poly)

Get the slices based on the polygon.

pyresample.slicer.create_slicer(area_to_crop, area_to_contain)

Create a slicer for cropping area_to_crop based on area_to_contain.

Return an AreaSlicer or a SwathSlicer based on the first area type.

pyresample.slicer.expand_slice(small_slice)

Expand slice by one.

pyresample.spherical module

Some generalized spherical functions.

base type is a numpy array of size (n, 2) (2 for lon and lats)

class pyresample.spherical.Arc(start, end)

Bases: object

An arc of the great circle between two points.

angle(other_arc)

Oriented angle between two arcs.

Returns:

Angle in radians. A straight line will be 0. A clockwise path will be a negative angle and counter-clockwise will be positive.

get_next_intersection(arcs, known_inter=None)

Get the next intersection between the current arc and arcs.

intersection(other_arc)

Return where, if the current arc and the other_arc intersect.

None is returned if there is not intersection. An arc is defined as the shortest tracks between two points.

intersections(other_arc)

Give the two intersections of the greats circles defined by the current arc and other_arc.

From http://williams.best.vwh.net/intersect.htm

intersects(other_arc)

Check if the current arc and the other_arc intersect.

An arc is defined as the shortest tracks between two points.

class pyresample.spherical.CCoordinate(cart)

Bases: object

Cartesian coordinates.

cross(point)

Get cross product with another vector.

The cross product of the same vector gives a zero vector.

dot(point)

Get dot product with another vector.

norm()

Get Euclidean norm of the vector.

normalize(inplace=False)

Normalize the vector.

If self.cart == [0,0,0], norm=0, and cart becomes [nan, nan, nan]: Note that self.cart == [0,0,0] can occurs when computing: - the cross product of the same point. - the cross product between points lying at the equator. - the cross product between points lying at the poles.

to_spherical()

Convert to SPoint/SMultiPoint object.

class pyresample.spherical.SCoordinate(lon, lat)

Bases: object

Spherical coordinates.

The lon and lat coordinates must be provided in radians.

cross2cart(point)

Compute the cross product, and convert to cartesian coordinates.

Note: - the cross product of the same point gives a zero vector. - the cross product between points lying at the equator gives a zero vector. - the cross product between points lying at the poles.

distance(point)

Get distance using Vincenty formula.

The result must be multiplied by Earth radius to obtain distance in m or km.

hdistance(point)

Get distance using Haversine formula.

The result must be multiplied by Earth radius to obtain distance in m or km.

plot(ax=None, projection_crs=None, add_coastlines=True, add_gridlines=True, add_background=True, **plot_kwargs)

Plot the point(s) using Cartopy.

Assume vertices to be in radians.

to_cart()

Convert to cartesian.

property vertices

Return point(s) radians vertices in a ndarray of shape [n,2].

property vertices_in_degrees

Return point(s) degrees vertices in a ndarray of shape [n,2].

class pyresample.spherical.SphPolygon(vertices, radius=1)

Bases: object

Spherical polygon.

Represents a polygon on a spherical geoid. Initialise with an ndarray of shape [N, 2] where the first column contains longitudes and the second column contains latitudes. The units should be in radians. The inside of the polygon is defined by the vertices being defined clockwise around it.

The optional second argument radius indicates the radius of the spherical geoid on which calculations occur.

aedges()

Get generator over the edges, in arcs of Coordinates.

area()

Find the area of a polygon.

The inside of the polygon is defined by having the vertices enumerated clockwise around it.

Uses the algorithm described in [bev1987].

[bev1987]

, Michael Bevis and Greg Cambareri, “Computing the area of a spherical polygon of arbitrary shape”, in Mathematical Geology, May 1987, Volume 19, Issue 4, pp 335-346.

Note: The article mixes up longitudes and latitudes in equation 3! Look at the fortran code appendix for the correct version.

The units are the square of the radius passed to the constructor. For example, to calculate the area in km² of a polygon near the equator of a spherical planet with a radius of 6371 km (similar to Earth):

>>> pol = SphPolygon(np.deg2rad(np.array([[0., 0.], [0., 1.], [1., 1.], [1., 0.]])),
                     radius=6371)
>>> print(pol.area())
12363.997753690213

If SphPolygon was constructed without passing any units, the result has units of square radii (i.e., the polygon containing the entire planet would have area 4π).

edges()

Get generator over the edges, in geographical coordinates.

intersection(other)

Return the intersection of this and other polygon.

inverse()

Return an inverse of the polygon.

invert()

Invert the polygon.

union(other)

Return the union of this and other polygon.

NB! If the two polygons do not overlap (they have nothing in common) None is returned.

pyresample.spherical_geometry module

Classes for spherical geometry operations.

class pyresample.spherical_geometry.Arc(start, end)

Bases: object

An arc of the great circle between two points.

angle(other_arc, snap=True)

Get oriented angle between two arcs.

Parameters:
  • other_arc (pyresample.spherical_geometry.Arc) –

  • snap (boolean) – Snap small angles to 0. Allows for detecting colinearity. Disable snapping when calculating polygon areas as it might lead to negative area values.

center_angle()

Get angle of an arc at the center of the sphere.

intersection(other_arc)

Determine the intersection point between this arc and another.

An arc is defined as the shortest tracks between two points.

intersections(other_arc)

Get the two intersections of the greats circles defined by the current arc and other_arc.

intersects(other_arc)

Determine if this arc and another arc intersect.

An arc is defined as the shortest tracks between two points.

class pyresample.spherical_geometry.Coordinate(lon=None, lat=None, x__=None, y__=None, z__=None, R__=1)

Bases: object

Point on earth in terms of lat and lon.

It expects lon,lat in degrees But self.lat and self.lon are returned in radians !

cross(point)

Get cross product with another vector.

cross2cart(point)

Compute the cross product, and convert to cartesian coordinates (assuming radius 1).

distance(point)

Get distance using Vincenty formula.

dot(point)

Get dot product with another vector.

lat = None
lon = None
norm()

Return the norm of the vector.

normalize()

Normalize the vector.

x__ = None
y__ = None
z__ = None
pyresample.spherical_geometry.get_first_intersection(b__, boundaries)

Get the first intersection on b__ with boundaries.

pyresample.spherical_geometry.get_intersections(b__, boundaries)

Get the intersections of b__ with boundaries.

Returns both the intersection coordinates and the concerned boundaries.

pyresample.spherical_geometry.get_next_intersection(p__, b__, boundaries)

Get the next intersection from the intersection of arcs p__ and b__ along segment b__ with boundaries.

pyresample.spherical_geometry.get_polygon_area(corners)

Get the area of the convex area defined by corners.

pyresample.spherical_geometry.intersection_polygon(area_corners, segment_corners)

Get the intersection polygon between two areas.

pyresample.spherical_geometry.modpi(val)

Put val between -pi and pi.

pyresample.spherical_geometry.point_inside(point, corners)

Determine if points are inside 4 corner points.

This uses great circle arcs as area boundaries.

pyresample.spherical_utils module

Functions to support the calculation of a coverage of an area by a set of spherical polygons.

It can for instance be a set of satellite overpasses to be received of a given local stations over a certain time window where we want to calculate how much of an area is covered by the onboard scanning instrument(s).

class pyresample.spherical_utils.GetNonOverlapUnions(polygons)

Bases: GetNonOverlapUnionsBaseClass

NonOverlapUnions class.

class pyresample.spherical_utils.GetNonOverlapUnionsBaseClass(geom_objects)

Bases: object

Base class to get the smallest set of union objects that does not overlap.

The objects are here Python sets of integers - but are abstracts for geometrical shapes on a sphere.

get_ids()

Get a list of identifiers identifying the gemoetry objects in each polygon union.

get_polygons()

Get a list of all non-overlapping polygon unions.

merge()

Merge all overlapping objects (sets or polygons).

pyresample.spherical_utils.check_if_two_polygons_overlap(polygon1, polygon2)

Check if two SphPolygons overlaps.

pyresample.spherical_utils.check_keys_int_or_tuple(adict)

Check if the dictionary keys are integers or tuples.

If they are not, raise a KeyError

pyresample.spherical_utils.merge_tuples(atuple)

Take a nested tuple of integers and concatenate it to a tuple of integers.

pyresample.version module

pyresample.version.get_versions()

Module contents

Pyresample package for geographic data resampling and related utilities.

pyresample.convert_def_to_yaml(def_area_file, yaml_area_file)

Convert a legacy area def file to the yaml counter partself.

yaml_area_file will be overwritten by the operation.

pyresample.create_area_def(area_id, projection, width=None, height=None, area_extent=None, shape=None, upper_left_extent=None, center=None, resolution=None, radius=None, units=None, optimize_projection=False, **kwargs)

Create AreaDefinition from whatever information is known.

Parameters:
  • area_id (str) – ID of area

  • projection (pyproj CRS object, dict, str, int, tuple, object) – Projection parameters. This can be in any format understood by pyproj.crs.CRS.from_user_input(), such as a pyproj CRS object, proj4 dict, proj4 string, EPSG integer code, or others.

  • description (str, optional) – Description/name of area. Defaults to area_id

  • proj_id (str, optional) – ID of projection (deprecated)

  • units (str, optional) –

    Units that provided arguments should be interpreted as. This can be one of ‘deg’, ‘degrees’, ‘meters’, ‘metres’, and any parameter supported by the cs2cs -lu command. Units are determined in the following priority:

    1. units expressed with each variable through a DataArray’s attrs attribute.

    2. units passed to units

    3. units used in projection

    4. meters

  • width (str, optional) – Number of pixels in the x direction

  • height (str, optional) – Number of pixels in the y direction

  • area_extent (list, optional) – Area extent as a list (lower_left_x, lower_left_y, upper_right_x, upper_right_y)

  • shape (list, optional) – Number of pixels in the y and x direction (height, width)

  • upper_left_extent (list, optional) – Upper left corner of upper left pixel (x, y)

  • center (list, optional) – Center of projection (x, y)

  • resolution (list or float, optional) – Size of pixels: (dx, dy)

  • radius (list or float, optional) – Length from the center to the edges of the projection (dx, dy)

  • nprocs (int, optional) – Number of processor cores to be used

  • lons (numpy array, optional) – Grid lons

  • lats (numpy array, optional) – Grid lats

  • optimize_projection – Whether the projection parameters have to be optimized for a DynamicAreaDefinition.

Returns:

area – If shape and area_extent are found, an AreaDefinition object is returned. If only shape or area_extent can be found, a DynamicAreaDefinition object is returned

Return type:

pyresample.geometry.AreaDefinition or pyresample.geometry.DynamicAreaDefinition

Raises:

ValueError: – If neither shape nor area_extent could be found

Notes

  • resolution and radius can be specified with one value if dx == dy

  • If resolution and radius are provided as angles, center must be given or findable. In such a case, they represent [projection x distance from center[0] to center[0]+dx, projection y distance from center[1] to center[1]+dy]

pyresample.get_area_def(area_id, area_name, proj_id, proj4_args, width, height, area_extent)

Construct AreaDefinition object from arguments.

Parameters:
  • area_id (str) – ID of area

  • area_name (str) – Description of area

  • proj_id (str) – ID of projection

  • proj4_args (dict, CRS, or str) – Projection information passed to pyproj’s CRS object

  • width (int) – Number of pixel in x dimension

  • height (int) – Number of pixel in y dimension

  • area_extent (list | tuple) – Area extent as a list of ints (LL_x, LL_y, UR_x, UR_y)

Returns:

area_def – AreaDefinition object

Return type:

object

pyresample.load_area(area_file_name, *regions)

Load area(s) from area file.

Parameters:
  • area_file_name (str, pathlib.Path, stream, or list thereof) – List of paths or streams. Any str or pathlib.Path will be interpreted as a path to a file. Any stream will be interpreted as containing a yaml definition file. To read directly from a string, use load_area_from_string().

  • regions (str argument list) – Regions to parse. If no regions are specified all regions in the file are returned

Returns:

area_defs – If one area name is specified a single AreaDefinition object is returned. If several area names are specified a list of AreaDefinition objects is returned

Return type:

pyresample.geometry.AreaDefinition or list

Raises:

AreaNotFound: – If a specified area name is not found

pyresample.parse_area_file(area_file_name, *regions)

Parse area information from area file.

Parameters:
  • area_file_name (str or list) – One or more paths to area definition files

  • regions (str argument list) – Regions to parse. If no regions are specified all regions in the file are returned

Returns:

area_defs – List of AreaDefinition objects

Return type:

list

Raises:

AreaNotFound: – If a specified area is not found