pyresample API

pyresample.geometry

Classes for geometry operations

class pyresample.geometry.AreaDefinition(area_id, name, proj_id, proj_dict, x_size, y_size, area_extent, nprocs=1, lons=None, lats=None, dtype=<type 'numpy.float64'>)

Holds definition of an area.

Parameters:
  • area_id (str) – ID of area
  • name (str) – Name of area
  • proj_id (str) – ID of projection
  • proj_dict (dict) – Dictionary with Proj.4 parameters
  • x_size (int) – x dimension in number of pixels
  • y_size (int) – y dimension in number of pixels
  • area_extent (list) – Area extent as a list (LL_x, LL_y, UR_x, UR_y)
  • nprocs (int, optional) – Number of processor cores to be used
  • lons (numpy array, optional) – Grid lons
  • lats (numpy array, optional) – Grid lats
area_id

str – ID of area

name

str – Name of area

proj_id

str – ID of projection

proj_dict

dict – Dictionary with Proj.4 parameters

x_size

int – x dimension in number of pixels

y_size

int – y dimension in number of pixels

shape

tuple – Corresponding array shape as (rows, cols)

size

int – Number of points in grid

area_extent

tuple – Area extent as a tuple (LL_x, LL_y, UR_x, UR_y)

area_extent_ll

tuple – Area extent in lons lats as a tuple (LL_lon, LL_lat, UR_lon, UR_lat)

pixel_size_x

float – Pixel width in projection units

pixel_size_y

float – Pixel height in projection units

pixel_upper_left

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

pixel_offset_x

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

pixel_offset_y

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

proj4_string

str – Projection defined as Proj.4 string

lons

object – Grid lons

lats

object – Grid lats

cartesian_coords

object – Grid cartesian coordinates

projection_x_coords

object – Grid projection x coordinate

projection_y_coords

object – Grid projection y coordinate

colrow2lonlat(cols, rows)

Return longitudes and latitudes 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).

get_lonlat(row, col)

Retrieves lon and lat values of single point in area grid

Parameters:
  • row (int) –
  • col (int) –
Returns:

(lon, lat)

Return type:

tuple of floats

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

Returns lon and lat arrays of area.

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 result the result. Requires data_slice to be None
Returns:

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

Return type:

tuple of numpy arrays

get_proj_coords(data_slice=None, cache=False, dtype=None)

Get projection coordinates of grid

Parameters:
  • data_slice (slice object, optional) – Calculate only coordinates for specified slice
  • cache (bool, optional) – Store result the result. Requires data_slice to be None
Returns:

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

Return type:

tuple of numpy arrays

get_xy_from_lonlat(lon, lat)

Retrieve 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.

Input:

lon : point or sequence (list or array) of longitudes lat : point or sequence (list or array) of latitudes

Returns:

(x, y) : tuple of integer points/arrays

get_xy_from_proj_coords(xm_, ym_)

Retrieve closest x and y coordinates (column, row indices) for a location specified with projection coordinates (xm_,ym_) in meters. A ValueError is raised, if the return point is outside the area domain. If xm_,ym_ is a tuple of sequences of projection coordinates, a tuple of masked arrays are returned.

Input:

xm_ : point or sequence (list or array) of x-coordinates in m (map projection) ym_ : point or sequence (list or array) of y-coordinates in m (map projection)

Returns:

(x, y) : tuple of integer points/arrays

lonlat2colrow(lons, lats)

Return image columns and rows for the given longitudes and latitudes. Both scalars and arrays are supported. Same as get_xy_from_lonlat, renamed for convenience.

outer_boundary_corners

Returns the lon,lat of the outer edges of the corner points

proj4_string

Returns projection definition as Proj.4 string

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

Base class for geometry definitions

corners

Returns the corners of the current area.

get_area()

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

get_area_extent_for_subset(row_LR, col_LR, row_UL, col_UL)

Retrieves area_extent for a subdomain rows are counted from upper left to lower left 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
: list
Area extent as a list (LL_x, LL_y, UR_x, UR_y) of the subset
Author:

Ulrich Hamann

get_boundary_lonlats()

Returns 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_lonlat(row, col)

Retrieve lon and lat of single pixel

Parameters:
  • row (int) –
  • col (int) –
Returns:

(lon, lat)

Return type:

tuple of floats

get_lonlats(data_slice=None, **kwargs)

Base method for lon lat retrieval with slicing

intersection(other)

Returns 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
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)

Tests 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
class pyresample.geometry.Boundary(side1, side2, side3, side4)

Container for geometry boundary. Labelling starts in upper left corner and proceeds clockwise

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

Base class for geometry definitions defined by lons and lats only

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

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

tuple – Grid shape as (rows, cols)

size

int – Number of elements in grid

lons

object – Grid lons

lats

object – Grid lats

cartesian_coords

object – Grid cartesian coordinates

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

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.
shape

tuple – Swath shape

size

int – Number of elements in swath

ndims

int – Swath dimensions

lons

object – Swath lons

lats

object – Swath lats

cartesian_coords

object – Swath cartesian coordinates

pyresample.image

Handles resampling of images with assigned geometry definitions

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

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

numpy array – Image data

geo_def

object – Geometry definition

fill_value

int or None – Resample result fill value

nprocs

int – Number of processor cores to be used for geometry operations

get_array_from_linesample(row_indices, col_indices)

Samples 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)

Base method for resampling from preprocessed data.

resample(target_geo_def)

Base method for resampling

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

Holds image with geometry definition. Allows nearest neighbour resampling 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

numpy array – Image data

geo_def

object – Geometry definition

radius_of_influence

float – Cut off distance in meters

epsilon

float – Allowed uncertainty in meters

fill_value

int or None – Resample result fill value

reduce_data

bool – Perform coarse data reduction before resampling

nprocs

int – Number of processor cores to be used

segments

int or None – Number of segments to use when resampling

resample(target_geo_def)

Resamples 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)

Holds image with area definition. ‘ 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

numpy array – Image data

geo_def

object – Area definition as AreaDefinition object

fill_value

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

nprocs

int – Number of processor cores to be used

segments

int or None – Number of segments to use when resampling

resample(target_area_def)

Resamples 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.grid

Resample image from one projection to another using nearest neighbour method in cartesian projection coordinate systems

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

Samples 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)

Samples from image based on lon lat arrays using nearest neighbour method in cartesian projection coordinate systems.

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)

Returns 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)

Resamples 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.kd_tree

Handles reprojection of geolocated data. Several types of resampling are supported

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)

Returns 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
  • 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 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 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 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.kd_tree.which_kdtree()

Returns the name of the kdtree used for resampling

pyresample.utils

Utility functions for pyresample

exception pyresample.utils.AreaNotFound

Exception raised when specified are is no found in file

pyresample.utils.fwhm2sigma(fwhm)

Calculate sigma for gauss function from FWHM (3 dB level)

Parameters:fwhm (float) – FWHM of gauss function (3 dB level of beam footprint)
Returns:sigma – sigma for use in resampling gauss function
Return type:float
pyresample.utils.generate_nearest_neighbour_linesample_arrays(source_area_def, target_area_def, radius_of_influence, nprocs=1)

Generate linesample arrays for nearest neighbour grid resampling

Parameters:
  • source_area_def (object) – Source area definition as AreaDefinition object
  • target_area_def (object) – Target area definition as AreaDefinition object
  • radius_of_influence (float) – Cut off distance in meters
  • nprocs (int, optional) – Number of processor cores to be used
Returns:

(row_indices, col_indices)

Return type:

tuple of numpy arrays

pyresample.utils.generate_quick_linesample_arrays(source_area_def, target_area_def, nprocs=1)

Generate linesample arrays for quick grid resampling

Parameters:
  • source_area_def (object) – Source area definition as AreaDefinition object
  • target_area_def (object) – Target area definition as AreaDefinition object
  • nprocs (int, optional) – Number of processor cores to be used
Returns:

(row_indices, col_indices)

Return type:

tuple of numpy arrays

pyresample.utils.get_area_def(area_id, area_name, proj_id, proj4_args, x_size, y_size, area_extent)

Construct AreaDefinition object from arguments

Parameters:
  • area_id (str) – ID of area
  • proj_id (str) – ID of projection
  • area_name (str) – Description of area
  • proj4_args (list or str) – Proj4 arguments as list of arguments or string
  • x_size (int) – Number of pixel in x dimension
  • y_size (int) – Number of pixel in y dimension
  • area_extent (list) – Area extent as a list of ints (LL_x, LL_y, UR_x, UR_y)
Returns:

area_def – AreaDefinition object

Return type:

object

pyresample.utils.load_area(area_file_name, *regions)

Load area(s) from area file

Parameters:
  • area_file_name (str) – Path to area definition file
  • 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:

object or list

Raises:

AreaNotFound: – If a specified area name is not found

pyresample.utils.parse_area_file(area_file_name, *regions)

Parse area information from area file

Parameters:
  • area_file_name (str) – Path to area definition file
  • 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.utils.wrap_longitudes(lons)

Wrap longitudes to the [-180:+180[ validity range (preserves dtype)

Parameters:lons (numpy array) – Longitudes in degrees
Returns:lons – Longitudes wrapped into [-180:+180[ validity range
Return type:numpy array

pyresample.data_reduce

Reduce data sets based on geographical information

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

Calculates 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)

Calculates 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)

Makes 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)

Makes 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)

Makes 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.plot

plot.area_def2basemap(area_def, **kwargs)

Get Basemap object from AreaDefinition

Parameters:
  • area_def (object) – geometry.AreaDefinition object
  • **kwargs (Keyword arguments) – Additional initialization arguments for Basemap
Returns:

bmap

Return type:

Basemap object

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
plot.save_quicklook(filename, area_def, data, vmin=None, vmax=None, label='Variable (units)', num_meridians=45, num_parallels=10, coast_res='c', backend='AGG')

Display 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
  • backend (str, optional) – matplotlib backend to use’
plot.show_quicklook(area_def, data, vmin=None, vmax=None, label='Variable (units)', num_meridians=45, num_parallels=10, coast_res='c')

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.ewa