pyku.geo#

pyku.geo.align_georeferencing(ds, ref, tolerance=None)[source]#

Align georeferencing of one dataset onto the other.

Parameters:
  • ds (xarray.Dataset) – The dataset to be aligned.

  • ref (xarray.Dataset) – The reference dataset for alignment.

  • tolerance (float) – Absolute tolerance for projection yx coordinates.

Raises:

Exception – If the y and x projection coordinates are outside the tolerance and the georeferencing cannot be aligned.

Returns:

Dataset with georeferencing aligned to reference dataset.

Return type:

xarray.Dataset

pyku.geo.apply_georeferencing(ds, area_def)[source]#

Apply/force georeferencing to dataset. This function is usefull if data with a known georeferencing are available but the georeferencing is not correct/broken.

Parameters:
  • ds (xarray.Dataset) – The input dataset.

  • area_def (pyresample.AreaDefinition) – The area definition.

Returns:

Dataset containing projection coordinates, geographical coordinates according to the georeferincing dataset.

Return type:

xarray.Dataset

Example

Exemplary data are loaded and regridded to a default pyku projection. Then the georeferencing is broken for the example.

In [1]: %%time
   ...: import pyku
   ...: # Get data and project to HYRAS 5km resolution grid
   ...: # -------------------------------------------------
   ...: 
   ...: ds = pyku.resources.get_test_data('global_data')
   ...: ds = ds.pyku.project('HYR-LAEA-5')
   ...: 
   ...: # Break georeferencing
   ...: # --------------------
   ...: 
   ...: broken = ds.drop_vars(['crs', 'y', 'x'])
   ...: print("Broken dataset")
   ...: print(broken)
   ...: 
Broken dataset
<xarray.Dataset> Size: 84MB
Dimensions:    (time: 366, y: 220, x: 258, bnds: 2)
Coordinates:
  * time       (time) datetime64[ns] 3kB 1980-01-01T12:00:00 ... 1980-12-31T1...
    height     float64 8B 2.0
    lat        (y, x) float64 454kB 55.13 55.13 55.14 ... 45.09 45.08 45.08
    lon        (y, x) float64 454kB 1.95 2.028 2.106 2.184 ... 19.7 19.76 19.83
Dimensions without coordinates: y, x, bnds
Data variables:
    tas        (time, y, x) float32 83MB dask.array<chunksize=(366, 220, 258), meta=np.ndarray>
    time_bnds  (time, bnds) datetime64[ns] 6kB dask.array<chunksize=(366, 2), meta=np.ndarray>
Attributes: (12/46)
    Conventions:            CF-1.7 CMIP-6.2
    activity_id:            CMIP
    branch_method:          standard
    branch_time_in_child:   0.0
    branch_time_in_parent:  65744.0
    contact:                cmip6-data@ec-earth.org
    ...                     ...
    variant_label:          r1i1p1f1
    license:                CMIP6 model data produced by EC-Earth-Consortium ...
    cmor_version:           3.4.0
    tracking_id:            hdl:21.14100/712ab96c-590d-4a68-9efa-013ecb56371f
    history:                2019-06-06T10:58:51Z ; CMOR rewrote data to be co...
    CORDEX_domain:          undefined
CPU times: user 254 ms, sys: 284 ms, total: 539 ms
Wall time: 495 ms

The georeferencing is then repaired by applying the georeferencing to the broken data.

In [2]: %%time
   ...: import pyku.geo as geo
   ...: repaired = geo.apply_georeferencing(broken, 'HYR-LAEA-5')
   ...: print(repaired)
   ...: 
<xarray.Dataset> Size: 84MB
Dimensions:    (time: 366, y: 220, x: 258, bnds: 2)
Coordinates:
  * time       (time) datetime64[ns] 3kB 1980-01-01T12:00:00 ... 1980-12-31T1...
  * y          (y) float64 2kB 3.586e+06 3.582e+06 ... 2.496e+06 2.492e+06
  * x          (x) float64 2kB 3.808e+06 3.814e+06 ... 5.088e+06 5.094e+06
    height     float64 8B 2.0
    lat        (y, x) float64 454kB 55.13 55.13 55.14 ... 45.09 45.08 45.08
    lon        (y, x) float64 454kB 1.95 2.028 2.106 2.184 ... 19.7 19.76 19.83
Dimensions without coordinates: bnds
Data variables:
    tas        (time, y, x) float32 83MB dask.array<chunksize=(366, 220, 258), meta=np.ndarray>
    time_bnds  (time, bnds) datetime64[ns] 6kB dask.array<chunksize=(366, 2), meta=np.ndarray>
    crs        int32 4B 1
Attributes: (12/46)
    Conventions:            CF-1.7 CMIP-6.2
    activity_id:            CMIP
    branch_method:          standard
    branch_time_in_child:   0.0
    branch_time_in_parent:  65744.0
    contact:                cmip6-data@ec-earth.org
    ...                     ...
    variant_label:          r1i1p1f1
    license:                CMIP6 model data produced by EC-Earth-Consortium ...
    cmor_version:           3.4.0
    tracking_id:            hdl:21.14100/712ab96c-590d-4a68-9efa-013ecb56371f
    history:                2019-06-06T10:58:51Z ; CMOR rewrote data to be co...
    CORDEX_domain:          undefined
CPU times: user 181 ms, sys: 19.4 ms, total: 200 ms
Wall time: 152 ms
pyku.geo.are_longitudes_wrapped(ds)[source]#

Check if the longitudes in the dataset are within the range [-180, 180).

Parameters:

ds (xarray.Dataset) – The input dataset containing longitude values.

Returns:

True if all longitudes are within the range [-180, 180), otherwise False.

Return type:

bool

Example

In [1]: import pyku
   ...: ds = pyku.resources.get_test_data('global_data')
   ...: ds.pyku.are_longitudes_wrapped()
   ...: 
Out[1]: False
pyku.geo.are_yx_projection_coordinates_strictly_monotonic(ds)[source]#

Checks if y projection coordinates are strictly increasing or strictly decreasing

Parameters:

ds (xarray.Dataset) – The input dataset.

Returns:

True if y and x coordinates are strictly monotonic.

Return type:

bool

Example

In [1]: import pyku
   ...: ds = pyku.resources.get_test_data('model_data')
   ...: ds.pyku.are_yx_projection_coordinates_strictly_monotonic()
   ...: 
Out[1]: True
pyku.geo.get_area_at_true_scale(ds)[source]#

Get the area of a single pixel at true scale. Given the dataset resolution and area definition, the area of a pixel centered at true scale is determined.

Parameters:

ds (xarray.Dataset) – Georeferenced dataset

Returns:

Area at true scale.

Return type:

float

Example

In [1]: import xarray
   ...: import pyku.geo as geo
   ...: ds = pyku.resources.get_test_data('hyras')
   ...: geo.get_area_at_true_scale(ds)
   ...: 
Out[1]: 25000001.923909195
pyku.geo.get_area_def(ds)[source]#

Get area definition.

Parameters:

ds (xarray.Dataset) – The input dataset.

Returns:

The area definition.

Return type:

pyresample.geometry.AreaDefinition

Example

In [1]: import xarray, pyku
   ...: ds = pyku.resources.get_test_data('model_data')
   ...: print(ds.pyku.get_area_def())
   ...: 
Area ID: rotated_pole
Description: rotated_pole
Projection: {'datum': 'WGS84', 'lon_0': '18', 'no_defs': 'None', 'o_lat_p': '39.25', 'o_lon_p': '0', 'o_proj': 'longlat', 'proj': 'ob_tran', 'type': 'crs'}
Number of columns: 424
Number of rows: 412
Area extent: (-28.43, 21.89, 18.21, -23.43)
pyku.geo.get_areas_cf_definitions()[source]#

Get the CF conform area definitions included in pyku

Returns:

Dictionary of definitions of CF conform area definitions

Return type:

dict

Example

In [1]: import pyku.geo as geo
   ...: geo.get_areas_cf_definitions().get('EUR-11').get('crs_cf')
   ...: 
pyku.geo.get_areas_definitions()[source]#

Get all default areas definitions included in pyku.

Returns:

Definitions of all areas included in pyku

Return type:

dict

Example

In [1]: import pyku.geo as geo
   ...: geo.get_areas_definitions().keys()
   ...: 
Out[1]: dict_keys(['EUR-44', 'EUR-11', 'EUR-6', 'GER-0275', 'EPISODES', 'radolan_wx', 'DWD_RADOLAN_1100x900', 'DWD_RADOLAN_900x900', 'HOSTRADA-1', 'linet_eqc', 'euradcom_ps', 'euradcom_eqc', 'world_360x180', 'world_540x270', 'world_7200x3600', 'world_3600x1800', 'world_43200x21600', 'seamless_europe', 'seamless_world', 'seamless_germany', 'HYR-LCC-1', 'HYR-LCC-3', 'HYR-LCC-5', 'HYR-LCC-125', 'HYR-LCC-50', 'HYR-LAEA-1', 'HYR-LAEA-2', 'HYR-LAEA-3', 'HYR-LAEA-5', 'HYR-LAEA-12.5', 'HYR-LAEA-50', 'HYR-LAEA-100', 'HYR-GER-LAEA-1', 'HYR-GER-LAEA-2', 'HYR-GER-LAEA-3', 'HYR-GER-LAEA-5', 'HYR-GER-LAEA-12.5', 'HYR-GER-LAEA-50', 'HYR-GER-LAEA-100', 'cerra'])
pyku.geo.get_georeferencing(ds)[source]#

Return georeferencing from dataset.

The returned dataset contains projection coordinates x/y, geographical coordinates lat/lon.

Parameters:

ds (xarray.Dataset) – The input dataset.

Returns:

Dataset containing projection coordinates, geographical coordinates.

Return type:

xarray.Dataset

Notes

Accepted names for geographic longitude, latitude, projection coordinates x and y are defined in pyku under etc/metadata.yaml

pyku.geo.get_lonlats(ds, dtype='ndarray')[source]#

Return lon lats in dataset

Parameters:
  • ds (xarray.Dataset) – The input dataset.

  • dtype (str) – Output type, either ‘ndarray’ for a tuple of ndarrays, or xarray.Dataset for an xarray dataset

Returns:

(lons, lats). If no latidudes or longitudes are found in the dataset, None is returned.

Return type:

Tuple[numpy.ndarray]

Example

In [1]: import pyku
   ...: 
   ...: ds = pyku.resources.get_test_data('hyras-tas-monthly')
   ...: 
   ...: ds.pyku.get_lonlats()
   ...: 
Out[1]: 
(array([[ 6.05492766,  6.06807811,  6.08122878, ..., 14.77308078,
         14.78621582, 14.79935059],
        [ 6.05424991,  6.06740261,  6.08055552, ..., 14.77389986,
         14.78703713, 14.80017414],
        [ 6.05357185,  6.0667268 ,  6.07988197, ..., 14.7747193 ,
         14.78785881, 14.80099806],
        ...,
        [ 5.30969451,  5.32531316,  5.34093224, ..., 15.67337877,
         15.68896785, 15.70455642],
        [ 5.30865437,  5.32427647,  5.33989898, ..., 15.67463484,
         15.69022733, 15.70581932],
        [ 5.30761368,  5.32323922,  5.33886517, ..., 15.67589158,
         15.69148749, 15.70708289]]),
 array([[47.1123055 , 47.1127836 , 47.11326011, ..., 47.07897204,
         47.07839197, 47.07781031],
        [47.12129323, 47.12177142, 47.12224802, ..., 47.08795317,
         47.08737298, 47.0867912 ],
        [47.13028091, 47.1307592 , 47.13123589, ..., 47.09693424,
         47.09635394, 47.09577204],
        ...,
        [55.07115969, 55.07174272, 55.07232382, ..., 55.03051822,
         55.02981111, 55.02910206],
        [55.08012406, 55.08070724, 55.08128848, ..., 55.03947233,
         55.03876505, 55.03805582],
        [55.08908843, 55.08967176, 55.09025315, ..., 55.04842644,
         55.04771897, 55.04700957]]))
pyku.geo.get_nx(ds)[source]#

Get number of pixels in the x/lon direction.

Parameters:

ds (xarray.Dataset) – The input dataset.

Returns:

Number of pixels in the x/lon direction.

Return type:

int

Example

In [1]: import pyku
   ...: ds = pyku.resources.get_test_data('hyras-tas-monthly')
   ...: ds.pyku.get_nx()
   ...: 
Out[1]: 665
pyku.geo.get_ny(ds)[source]#

Get number of pixels in the y/lat direction

Parameters:

ds (xarray.Dataset) – The input dataset.

Returns:

Number of pixels y/lat direction.

Return type:

int

Raises:

Exception – if the dimension of projection coordinate y is more than 1.

Example

In [1]: pyku
   ...: 
   ...: ds = pyku.resources.get_test_data('hyras-tas-monthly')
   ...: 
   ...: ds.pyku.get_ny()
   ...: 
Out[1]: 890
pyku.geo.get_yx(ds, dtype='ndarray')[source]#

Return y x projection coordinates in dataset.

Parameters:
  • ds (xarray.Dataset) – The input dataset.

  • dtype (str) – Output type, either ndarray for a tuple of ndarrays, or xr.Dataset.

Returns:

(y, x). If no y or x projection

coordinates are found in the dataset, None is returned.

Return type:

tuple[numpy.ndarray

pyku.geo.get_yx_area_extent(ds)[source]#

Get area exent in projection coordinates

Parameters:

ds (xarray.Dataset) – The input dataset.

Returns:

(lower_left_x, lower_left_y, upper_right_x, upper_right_y)

Return type:

tuple

pyku.geo.is_georeferencing_sorted(ds)[source]#

Check if the georeferencing coordinates in the dataset are sorted in a standard geospatial order: x-coordinates (longitude/easting) increasing from left to right and y-coordinates (latitude/northing) increasing from bottom to top.

Parameters:

ds (xarray.Dataset) – The input dataset containing georeferencing coordinates.

Returns:

True if the georeferencing coordinates are sorted as described, False otherwise.

Return type:

bool

Example

In [1]: import pyku
   ...: ds = pyku.resources.get_test_data('model_data')
   ...: ds.pyku.is_georeferencing_sorted()
   ...: 
Out[1]: False
pyku.geo.list_standard_areas()[source]#

List pyku standard area definitions

Returns:

List of pyku standard areas

Return type:

List[str]

pyku.geo.load_area_def(projection_name, area_file=None)[source]#

Load area definition from projection file and projection name.

Parameters:
  • projection_name (str) – The projection name.

  • area_file (str) – File containing projection definition

Returns:

Area definition

Return type:

pyresample.geometry.AreaDefinition

Example

In [1]: import pyku.geo as geo
   ...: geo.load_area_def('HYR-LAEA-5')
   ...: 
Out[1]: 
Area ID: HYR-LAEA-5
Description: HYR-LAEA-5
Projection: {'ellps': 'GRS80', 'lat_0': '52', 'lon_0': '10', 'no_defs': 'None', 'proj': 'laea', 'towgs84': '[0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0]', 'type': 'crs', 'units': 'm', 'x_0': '4321000', 'y_0': '3210000'}
Number of columns: 258
Number of rows: 220
Area extent: (3806000.0, 2489000.0, 5096000.0, 3589000.0)
pyku.geo.project(ds, area_def=None, roi=1000000.0, method='nearest_neighbor', use_dask=True, area_file=None, keep_mask=False, power_parameter=None)[source]#

Project dataset to target projection.

Parameters:
  • ds (xarray.Dataset) – The input dataset.

  • area_def (pyresample.geometry.AreaDefinition, str) – Target projection can be given as a pyresample.AreaDefinition, as a predefined projection identified by a string (e.g. ‘EUR-44’), or as an xarray.Dataset for swath resampling (Experimental).

  • roi (int) – Radius of influence for nearest-neighbor calculations (default 1000 km); it should be at least five times the data resolution to avoid data gaps, though lower values can be used to increase calculation speed.

  • method (str) –

    Specifies the resampling method to use. Options include:

    • nearest_neighbor (pyresample)

    • bilinear (pyresample)

    • idw (Inverse Distance Weighting, pyresample)

    • conservative (ESMF)

    The default method is nearest_neighbor. When using the idw method, the power_parameter must be specified. The resampling can be performed using either the pyresample or the ESMF library. The specific method and library can be explicitly chosen by passing one of the following options:

    • pyresample_nearest_neighbor

    • pyresample_bilinear

    • pyresample_idw

    • esmf_bilinear

    • esmf_conservative

    • esmf_conservative_normed

    • esmf_patch

    • esmf_nearest_s2d

    • esmf_nearest_d2s

    • pyresample_bilinear_swath_to_grid

    • pyresample_bilinear_swath_to_grid_legacy

    • pyresample_bilinear_swath_to_swath_legacy

    • pyresample_nearest_neighbor_legacy

    Depending on the choice, the resampling will be executed with the corresponding library and method.

  • keep_mask (bool) – If the data is projected to a region larger than the original, the nearest_neighbor algorithm alwas returns the value of the nearest data point and does not return NaNs! Hence in that case the radius of influence shall be reduced to set values to NaN beyond the roi, or the mask is calculated from the input and applied to the output.

  • use_dask (bool) – Whether to use dask for block resampling. Defaults to True.

  • area_file (str) – Path to area file where predefined projections which are defined by a string are defined.

  • power_parameter (float) – Power parameter for Inverse Distance Weighting (IDW). Defaults to 2.5. For other interpolation methods, this parameter is ignored.

Returns:

Dataset in target projection.

Return type:

xarray.Dataset

Examples

In [1]: import pyku
   ...: 
   ...: # Select first timestep and wrap longitudes of global data
   ...: ds = (
   ...:     pyku.resources.get_test_data('global_data')
   ...:     .isel(time=0)
   ...: )
   ...: 
   ...: # Project to HYRAS grid
   ...: projected = ds.pyku.project(area_def='EUR-44')
   ...: 
   ...: # Plot
   ...: projected.ana.one_map(var='tas', crs='EUR-44')
   ...: 
_images/geo_project_standard.png
In [2]: import pyku
   ...: 
   ...: # Select first timestep and wrap longitudes of global data
   ...: ds = (
   ...:     pyku.resources.get_test_data('global_data')
   ...:     .isel(time=0)
   ...: )
   ...: 
   ...: # Project to HYRAS grid
   ...: projected = ds.pyku.project(
   ...:    area_def='HYR-LAEA-5',
   ...:    method='bilinear'
   ...: )
   ...: 
   ...: # Plot
   ...: projected.ana.one_map(var='tas', crs='HYR-LAEA-5')
   ...: 
_images/geo_project_bilinear.png
In [3]: import pyku
   ...: 
   ...: # Select first timestep and wrap longitudes of global data
   ...: ds = (
   ...:    pyku.resources.get_test_data('global_data')
   ...:    .isel(time=0)
   ...: )
   ...: 
   ...: # Project to HYRAS grid
   ...: projected = ds.pyku.project(
   ...:    area_def='DWD_RADOLAN_900x900',
   ...:    method='idw',
   ...:    roi=500000,
   ...:    power_parameter=1
   ...: )
   ...: 
   ...: # Plot
   ...: projected.ana.one_map(var='tas', crs='DWD_RADOLAN_900x900')
   ...: 
_images/geo_project_idw.png
pyku.geo.select_area_extent(ds, lower_left_lat=None, lower_left_lon=None, upper_right_lat=None, upper_right_lon=None)[source]#

Select an area extent in geographic coordinates, defined by the latitude and longitude of the lower-left and upper-right corners.

Since only indexed coordinates can be selected in xarray, and these are typically the y and x projection coordinates rather than geographic latitudes and longitudes, this function allows selection based on geographic coordinates, which are generally not indexed.

Parameters:
  • ds (xarray.Dataset) – The input dataset.

  • lower_left_lat (float) – Latitude of the lower left corner.

  • lower_left_lon (float) – Longitude of the lower left corner.

  • upper_right_lat (float) – Latitude of the upper right corner.

  • upper_right_lon (float) – Longitude of the upper right corner.

Returns:

Dataset within the selected area extent.

Return type:

xarray.Dataset

Example

In [1]: import pyku
   ...: ds = pyku.resources.get_test_data('hyras')
   ...: print(ds.pyku.get_area_def())
   ...: 
   ...: ds = ds.pyku.select_area_extent(
   ...:     lower_left_lat=48.0,
   ...:     lower_left_lon=5.0,
   ...:     upper_right_lat=51.0,
   ...:     upper_right_lon=15.0,
   ...: )
   ...: 
   ...: ds.isel(time=0).ana.one_map(var='tas')
   ...: 
Area ID: crs
Description: crs
Projection: {'ellps': 'GRS80', 'lat_0': '52', 'lon_0': '10', 'no_defs': 'None', 'proj': 'laea', 'type': 'crs', 'units': 'm', 'x_0': '4321000', 'y_0': '3210000'}
Number of columns: 133
Number of rows: 178
Area extent: (4021000.0, 2674000.0, 4686000.0, 3564000.0)
pyku.geo.select_neighborhood(ds, lat, lon, roi=10000, neighbours=1000, crop=False)[source]#

Select all the given number of neighbours within the radius of influence (roi) around the selected point given in geographic coordinates.

Parameters:
  • ds (xarray.Dataset) – The input dataset.

  • lat (float) – Geographic latitude.

  • lon (float) – Geographic longitude.

  • roi (float) – Radius of influence, given in meters.

  • neighbours (int) – Number of neighbours to be selected within radius of influence. Only the neighbours within the radius of influence are returned. If the number of neigbours within the radius of influence is smaller than the value given, only a subset is returned.

  • crop (bool) – Crop dataset with selected neighborhood.

Returns:

Dataset within the selected neighborhood.

Return type:

xarray.Dataset

Example

In [1]: import pyku
   ...: ds = pyku.resources.get_test_data('hyras')
   ...: selection = ds.pyku.select_neighborhood(
   ...:     lat=51,
   ...:     lon=10,
   ...:     roi=50000,
   ...:     neighbours=1000,
   ...:     crop=False
   ...: )
   ...: selection.ana.mean_map(var='tas')
   ...: 
_images/geo_select_neighborhood.png
pyku.geo.set_latlon_bounds(ds)[source]#

Determine and set geographic coordinate bounds.

Parameters:

ds (xarray.Dataset) – The input dataset.

Returns:

The dataset with geographic coordinate bounds.

Return type:

xarray.Dataset

Example

In [1]: import pyku
   ...: ds = pyku.resources.get_test_data('hyras-tas-monthly')
   ...: ds.pyku.set_latlon_bounds()
   ...: 
Out[1]: 
<xarray.Dataset> Size: 162MB
Dimensions:             (time: 24, y: 890, x: 665, bnds: 2, bounds: 4)
Coordinates:
  * time                (time) datetime64[ns] 192B 1961-01-01 ... 1962-12-01
  * y                   (y) float64 7kB 2.674e+06 2.676e+06 ... 3.564e+06
  * x                   (x) float64 5kB 4.022e+06 4.022e+06 ... 4.686e+06
    lat                 (y, x) float64 5MB dask.array<chunksize=(890, 665), meta=np.ndarray>
    lon                 (y, x) float64 5MB dask.array<chunksize=(890, 665), meta=np.ndarray>
    lon_bounds          (y, x, bounds) float64 19MB 6.049 6.062 ... 15.72 15.7
    lat_bounds          (y, x, bounds) float64 19MB 47.11 47.11 ... 55.05 55.05
Dimensions without coordinates: bnds, bounds
Data variables:
    crs                 (time) int32 96B dask.array<chunksize=(24,), meta=np.ndarray>
    number_of_stations  (time) float64 192B dask.array<chunksize=(24,), meta=np.ndarray>
    tas                 (time, y, x) float64 114MB dask.array<chunksize=(24, 890, 665), meta=np.ndarray>
    time_bnds           (time, bnds) datetime64[ns] 384B dask.array<chunksize=(24, 2), meta=np.ndarray>
    x_bnds              (time, x, bnds) float64 255kB dask.array<chunksize=(24, 665, 2), meta=np.ndarray>
    y_bnds              (time, y, bnds) float64 342kB dask.array<chunksize=(24, 890, 2), meta=np.ndarray>
Attributes: (12/22)
    source:                 surface observations
    institution:            Deutscher Wetterdienst (DWD)
    Conventions:            CF-1.11
    title:                  gridded_temperature_dataset_(HYRAS TAS)
    realization:            v6-1
    project_id:             HYRAS
    ...                     ...
    ConventionsURL:         http://cfconventions.org/Data/cf-conventions/cf-c...
    license:                The HYRAS data, produced by DWD, is licensed unde...
    comment:                Please be aware that the parameters are stored as...
    filename:               tas_hyras_1_1961_v6-1_de_monmean.nc
    creation_date:          2025-07-04T06:55:30Z
    unique_dataset_id:      DWD_HYRAS_DE_tas_v6-1_1961_mon_159ccae4-6f61-462a...
pyku.geo.set_spatial_weights(ds, how='area')[source]#

Set area weights for a dataset.

Parameters:
  • ds (xarray.Dataset) – Dataset which should contain latitude (‘lat’) as a coordinate.

  • how (str) – Method to determine weights. Default is set to ‘area’. The gridcell areas will be calculated according to the lat and lon bounds. Optionally one can select weighting by cosinus of latitude how=’cos’.

Returns:

Dataset with an added ‘weights’ variable containing the calculated area weights.

Return type:

xarray.Dataset

Example

In [1]: import xarray, pyku
   ...: import pyku.analyse as analyse
   ...: ds = pyku.resources.get_test_data('global_data')
   ...: ds = ds.pyku.project('world_360x180')
   ...: ds = ds.pyku.set_spatial_weights(how='area')
   ...: analyse.one_map(ds.weights)
   ...: 
_images/geo_spatial_weights1.png
In [2]: import xarray, pyku
   ...: import pyku.analyse as analyse
   ...: ds = pyku.resources.get_test_data('global_data')
   ...: ds = ds.pyku.project('EUR-44')
   ...: ds = ds.pyku.set_spatial_weights(how='area')
   ...: analyse.one_map(ds.weights, crs='EUR-44')
   ...: 
_images/geo_spatial_weights2.png
pyku.geo.sort_georeferencing(ds)[source]#

Sort georeferencing. After this operation, the data are arranged from left to right and top to bottom.

While Pyku can handle any arranging of the data, following this recommended structure ensures a more standardized data layout, reducing the likelihood of encountering edge cases.

Parameters:

ds (xarray.Dataset) – The input data.

Returns:

Data with georeferencing sorted from left to right and top to bottom.

Return type:

xarray.Dataset

Example

In [1]: import matplotlib.pyplot as plt
   ...: import pyku
   ...: 
   ...: plt.close('all')
   ...: 
   ...: # Get coordinate reference system of original data
   ...: # ------------------------------------------------
   ...: 
   ...: ccrs_original = (
   ...:     pyku.resources.get_test_data('model_data')
   ...:     .pyku.get_area_def()
   ...:     .to_cartopy_crs()
   ...: )
   ...: 
   ...: # Get coordinate reference system of sorted data
   ...: # ----------------------------------------------
   ...: 
   ...: ccrs_sorted = (
   ...:     pyku.resources.get_test_data('model_data')
   ...:     .pyku.sort_georeferencing()
   ...:     .pyku.get_area_def()
   ...:     .to_cartopy_crs()
   ...: )
   ...: 
   ...: # Plot side-by-side
   ...: # -----------------
   ...: 
   ...: plt.figure(figsize=(10, 5))
   ...: 
   ...: ax1 = plt.subplot(2, 2, 1, projection=ccrs_original)
   ...: ax1.set_title('original')
   ...: ax1.set_global()
   ...: ax1.coastlines('auto')
   ...: ax1.gridlines()
   ...: 
   ...: ax2 = plt.subplot(2, 1, 1, projection=ccrs_sorted)
   ...: ax2.set_title('sorted')
   ...: ax2.set_global()
   ...: ax2.coastlines('auto')
   ...: ax2.gridlines()
   ...: 
Out[1]: <cartopy.mpl.gridliner.Gridliner at 0x7f53b8e60cb0>
_images/geo_sort_georeferencing.png
pyku.geo.wrap_longitudes(ds)[source]#

Wrap longitudes to [-180, +180[, and sort by increasing longitudes.

Parameters:

ds (xarray.Dataset) – The input data

Returns:

Data with longitudes wrapped to [-180,180[, and sorted by increasing longitude.

Return type:

xarray.Dataset

Example

In [1]: import pyku
   ...: ds = pyku.resources.get_test_data('global_data')
   ...: ds.pyku.wrap_longitudes()
   ...: 
Out[1]: 
<xarray.Dataset> Size: 192MB
Dimensions:    (lat: 256, bnds: 2, lon: 512, time: 366)
Coordinates:
  * lat        (lat) float64 2kB -89.46 -88.77 -88.07 ... 88.07 88.77 89.46
  * lon        (lon) float64 4kB -180.0 -179.3 -178.6 ... 177.9 178.6 179.3
  * time       (time) datetime64[ns] 3kB 1980-01-01T12:00:00 ... 1980-12-31T1...
    height     float64 8B ...
Dimensions without coordinates: bnds
Data variables:
    lat_bnds   (lat, bnds) float64 4kB ...
    lon_bnds   (lon, bnds) float64 8kB ...
    tas        (time, lat, lon) float32 192MB ...
    time_bnds  (time, bnds) datetime64[ns] 6kB ...
Attributes: (12/45)
    Conventions:            CF-1.7 CMIP-6.2
    activity_id:            CMIP
    branch_method:          standard
    branch_time_in_child:   0.0
    branch_time_in_parent:  65744.0
    contact:                cmip6-data@ec-earth.org
    ...                     ...
    variable_id:            tas
    variant_label:          r1i1p1f1
    license:                CMIP6 model data produced by EC-Earth-Consortium ...
    cmor_version:           3.4.0
    tracking_id:            hdl:21.14100/712ab96c-590d-4a68-9efa-013ecb56371f
    history:                2019-06-06T10:58:51Z ; CMOR rewrote data to be co...