API Reference

Cutout

class atlite.Cutout(path, **cutoutparams)

Cutout base class.

This class builds the starting point for most atlite functionalities.

availabilitymatrix(shapes, excluder, nprocesses=None, disable_progressbar=False)

Compute the eligible share within cutout cells in the overlap with shapes.

For parallel calculation (nprocesses not None) the excluder must not be initialized and all raster references must be strings. Otherwise processes are colliding when reading from one common rasterio.DatasetReader.

Parameters
  • cutout (atlite.Cutout) – Cutout which the availability matrix is aligned to.

  • shapes (geopandas.Series/geopandas.DataFrame) – Geometries for which the availabilities are calculated.

  • excluder (atlite.gis.ExclusionContainer) – Container of all meta data or objects which to exclude, i.e. rasters and geometries.

  • nprocesses (int, optional) – Number of processes to use for calculating the matrix. The paralle- lization can heavily boost the calculation speed. The default is None.

  • disable_progressbar (bool, optional) – Disable the progressbar if nprocesses is not None. Then the map function instead of the imap function is used for the multiprocessing pool. This speeds up the calculation.

Returns

availabilities – DataArray of shape (|shapes|, |y|, |x|) containing all the eligible share of cutout cell (x,y) in the overlap with shape i.

Return type

xr.DataArray

Notes

The rasterio (or GDAL) average downsampling returns different results dependent on how the target raster (the cutout raster) is spanned. Either it is spanned from the top left going downwards, e.g. Affine(0.25, 0, 0, 0, -0.25, 50), or starting in the lower left corner and going up, e.g. Affine(0.25, 0, 0, 0, 0.25, 50). Here we stick to the top down version which is why we use cutout.transform_r and flipping the y-axis in the end.

property bounds

Total bounds of the area covered by the cutout (x, y, X, Y).

convert_and_aggregate(convert_func, matrix=None, index=None, layout=None, shapes=None, shapes_crs=4326, per_unit=False, return_capacity=False, capacity_factor=False, show_progress=True, dask_kwargs={}, **convert_kwds)

Convert and aggregate a weather-based renewable generation time-series.

NOTE: Not meant to be used by the user him or herself. Rather it is a gateway function that is called by all the individual time-series generation functions like pv and wind. Thus, all its parameters are also available from these.

Parameters
  • matrix (N x S - xr.DataArray or sp.sparse.csr_matrix or None) – If given, it is used to aggregate the grid cells to buses. N is the number of buses, S the number of spatial coordinates, in the order of cutout.grid.

  • index (pd.Index) – Index of Buses.

  • layout (X x Y - xr.DataArray) – The capacity to be build in each of the grid_cells.

  • shapes (list or pd.Series of shapely.geometry.Polygon) – If given, matrix is constructed as indicatormatrix of the polygons, its index determines the bus index on the time-series.

  • shapes_crs (pyproj.CRS or compatible) – If different to the map crs of the cutout, the shapes are transformed to match cutout.crs (defaults to EPSG:4326).

  • per_unit (boolean) – Returns the time-series in per-unit units, instead of in MW (defaults to False).

  • return_capacity (boolean) – Additionally returns the installed capacity at each bus corresponding to layout (defaults to False).

  • capacity_factor (boolean) – If True, the static capacity factor of the chosen resource for each grid cell is computed.

  • show_progress (boolean, default True) – Whether to show a progress bar.

  • dask_kwargs (dict, default {}) – Dict with keyword arguments passed to dask.compute.

  • convert_func (Function) – Callback like convert_wind, convert_pv

Returns

  • resource (xr.DataArray) – Time-series of renewable generation aggregated to buses, if matrix or equivalents are provided else the total sum of generated energy.

  • units (xr.DataArray (optional)) – The installed units per bus in MW corresponding to layout (only if return_capacity is True).

property extent

Total extent of the area covered by the cutout (x, X, y, Y).

heat_demand(threshold=15.0, a=1.0, constant=0.0, hour_shift=0.0, **params)

Convert outside temperature into daily heat demand using the degree-day approximation.

Since “daily average temperature” means different things in different time zones and since xarray coordinates do not handle time zones gracefully like pd.DateTimeIndex, you can provide an hour_shift to redefine when the day starts.

E.g. for Moscow in winter, hour_shift = 4, for New York in winter, hour_shift = -5

This time shift applies across the entire spatial scope of ds for all times. More fine-grained control will be built in a some point, i.e. space- and time-dependent time zones.

WARNING: Because the original data is provided every month, at the month boundaries there is untidiness if you use a time shift. The resulting xarray will have duplicates in the index for the parts of the day in each month at the boundary. You will have to re-average these based on the number of hours in each month for the duplicated day.

Parameters
  • threshold (float) – Outside temperature in degrees Celsius above which there is no heat demand.

  • a (float) – Linear factor relating heat demand to outside temperature.

  • constant (float) – Constant part of heat demand that does not depend on outside temperature (e.g. due to water heating).

  • hour_shift (float) – Time shift relative to UTC for taking daily average

Note

You can also specify all of the general conversion arguments documented in the convert_and_aggregate function.

hydro(plants, hydrobasins, flowspeed=1, weight_with_height=False, show_progress=True, **kwargs)

Compute inflow time-series for plants by aggregating over catchment basins from hydrobasins

Parameters
  • plants (pd.DataFrame) – Run-of-river plants or dams with lon, lat columns.

  • hydrobasins (str|gpd.GeoDataFrame) – Filename or GeoDataFrame of one level of the HydroBASINS dataset.

  • flowspeed (float) – Average speed of water flows to estimate the water travel time from basin to plant (default: 1 m/s).

  • weight_with_height (bool) – Whether surface runoff should be weighted by potential height (probably better for coarser resolution).

  • show_progress (bool) – Whether to display progressbars.

References

[1] Liu, Hailiang, et al. “A validated high-resolution hydro power time-series model for energy systems analysis.” arXiv preprint arXiv:1901.08476 (2019).

[2] Lehner, B., Grill G. (2013): Global river hydrography and network routing: baseline data and new approaches to study the world’s large river systems. Hydrological Processes, 27(15): 2171–2186. Data is available at www.hydrosheds.org.

indicatormatrix(shapes, shapes_crs=4326)

Compute the indicatormatrix.

The indicatormatrix I[i,j] is a sparse representation of the ratio of the area in orig[j] lying in dest[i], where orig and dest are collections of polygons, i.e.

A value of I[i,j] = 1 indicates that the shape orig[j] is fully contained in shape dest[j].

Note that the polygons must be in the same crs.

Parameters

shapes (Collection of shapely polygons) –

Returns

I – Indicatormatrix

Return type

sp.sparse.lil_matrix

layout_from_capacity_list(data, col='Capacity')

Get a capacity layout aligned to the cutout based on a capacity list.

Parameters
  • data (pandas.DataFrame) – Capacity list with columns ‘x’, ‘y’ and col. Each capacity entry is added to the grid cell intersecting with the coordinate (x,y).

  • col (str, optional) – Name of the column with capacity values. The default is ‘Capacity’.

Returns

Capacity layout with dimensions ‘x’ and ‘y’ indicating the total capacity placed within one grid cell.

Return type

xr.DataArray

Example

>>> import atlite
>>> import powerplantmatching as pm
>>> data = pm.data.OPSD_VRE_country('DE')
>>> data = (data.query('Fueltype == "Solar"')
            .rename(columns={'lon':'x', 'lat':'y'}))
>>> cutout = atlite.Cutout('Germany', x = slice(-5, 15), y = slice(40, 55),
               time='2013-06-01', module='era5')
>>> cutout.prepare(features=['influx', 'temperature'])
>>> layout = cutout.layout_from_capacity_list(data)
>>> pv = cutout.pv('CdTe', 'latitude_optimal', layout=layout)
>>> pv.plot()
merge(other, path=None, **kwargs)

Merge two cutouts into a single cutout.

Parameters
  • other (atlite.Cutout) – Other cutout to merge.

  • path (str | path-like) – File where to store the merged cutout. Defaults to a temporary file.

  • **kwargs – Keyword arguments passed to xarray.merge().

Returns

merged – Merged cutout.

Return type

Cutout

prepare(features=None, tmpdir=None, overwrite=False)

Prepare all or a selection of features in a cutout.

This function loads the feature data of a cutout, e.g. influx or runoff. When not specifying the feature argument, all available features will be loaded. The function compares the variables which are already included in the cutout with the available variables of the modules specified by the cutout. It detects missing variables and stores them into the netcdf file of the cutout.

Parameters
  • cutout (atlite.Cutout) –

  • features (str/list, optional) – Feature(s) to be prepared. The default slice(None) results in all available features.

  • tmpdir (str/Path, optional) – Directory in which temporary files (for example retrieved ERA5 netcdf files) are stored. If set, the directory will not be deleted and the intermediate files can be examined.

  • overwrite (bool, optional) – Whether to overwrite variables which are already included in the cutout. The default is False.

Returns

cutout – Cutout with prepared data. The variables are stored in cutout.data.

Return type

atlite.Cutout

pv(panel, orientation, clearsky_model=None, **params)

Convert downward-shortwave, upward-shortwave radiation flux and ambient temperature into a pv generation time-series.

Parameters
  • panel (str or dict) – Panel config dictionary with the parameters for the electrical model in [3]. Alternatively, name of yaml file stored in atlite.config.solarpanel_dir.

  • orientation (str, dict or callback) – Panel orientation can be chosen from either ‘latitude_optimal’, a constant orientation {‘slope’: 0.0, ‘azimuth’: 0.0} or a callback function with the same signature as the callbacks generated by the ‘atlite.pv.orientation.make_*’ functions.

  • clearsky_model (str or None) – Either the ‘simple’ or the ‘enhanced’ Reindl clearsky model. The default choice of None will choose dependending on data availability, since the ‘enhanced’ model also incorporates ambient air temperature and relative humidity.

Returns

pv – Time-series or capacity factors based on additional general conversion arguments.

Return type

xr.DataArray

Note

You can also specify all of the general conversion arguments documented in the convert_and_aggregate function.

References

[1] Soteris A. Kalogirou. Solar Energy Engineering: Processes and Systems, pages 49–117,469–516. Academic Press, 2009. ISBN 0123745012.

[2] D.T. Reindl, W.A. Beckman, and J.A. Duffie. Diffuse fraction correla- tions. Solar Energy, 45(1):1 – 7, 1990.

[3] Hans Georg Beyer, Gerd Heilscher and Stefan Bofinger. A Robust Model for the MPP Performance of Different Types of PV-Modules Applied for the Performance Check of Grid Connected Systems, Freiburg, June 2004. Eurosun (ISES Europe Solar Congress).

sel(path=None, bounds=None, buffer=0, **kwargs)

Select parts of the cutout.

Parameters
  • path (str | path-like) – File where to store the sub-cutout. Defaults to a temporary file.

  • bounds (GeoSeries.bounds | DataFrame, optional) – The outer bounds of the cutout or as a DataFrame containing (min.long, min.lat, max.long, max.lat).

  • buffer (float, optional) – Buffer around the bounds. The default is 0.

  • **kwargs – Passed to xr.Dataset.sel for data selection.

Returns

selected – Selected cutout.

Return type

Cutout

solar_thermal(orientation={'azimuth': 180.0, 'slope': 45.0}, trigon_model='simple', clearsky_model='simple', c0=0.8, c1=3.0, t_store=80.0, **params)

Convert downward short-wave radiation flux and outside temperature into time series for solar thermal collectors.

Mathematical model and defaults for c0, c1 based on model in [1].

Parameters
  • cutout (cutout) –

  • orientation (dict or str or function) – Panel orientation with slope and azimuth (units of degrees), or ‘latitude_optimal’.

  • trigon_model (str) – Type of trigonometry model

  • clearsky_model (str or None) – Type of clearsky model for diffuse irradiation. Either ‘simple’ or ‘enhanced’.

  • c0 (float) – Parameters for model in [1] (defaults to 0.8 and 3., respectively)

  • c1 (float) – Parameters for model in [1] (defaults to 0.8 and 3., respectively)

  • t_store (float) – Store temperature in degree Celsius

Note

You can also specify all of the general conversion arguments documented in the convert_and_aggregate function.

References

[1] Henning and Palzer, Renewable and Sustainable Energy Reviews 30 (2014) 1003-1018

to_file(fn=None)

Save cutout to a netcdf file.

Parameters

fn (str | path-like) – File name where to store the cutout, defaults to cutout.path.

property transform

Get the affine transform of the cutout.

property transform_r

Get the affine transform of the cutout with reverse y-order.

uniform_layout()

Get a uniform capacity layout for all grid cells.

wind(turbine, smooth=False, **params)

Generate wind generation time-series

Extrapolates 10m wind speed with monthly surface roughness to hub height and evaluates the power curve.

Parameters
  • turbine (str or dict) – A turbineconfig dictionary with the keys ‘hub_height’ for the hub height and ‘V’, ‘POW’ defining the power curve. Alternatively a str refering to a local or remote turbine configuration as accepted by atlite.resource.get_windturbineconfig().

  • smooth (bool or dict) – If True smooth power curve with a gaussian kernel as determined for the Danish wind fleet to Delta_v = 1.27 and sigma = 2.29. A dict allows to tune these values.

Note

You can also specify all of the general conversion arguments documented in the convert_and_aggregate function.

References

[1] Andresen G B, Søndergaard A A and Greiner M 2015 Energy 93, Part 1 1074 – 1088. doi:10.1016/j.energy.2015.09.071

Data

Management of data retrieval and structure.

atlite.data.available_features(module=None)

Inspect the available features of all or a selection of modules.

Parameters

module (str/list, optional) – Module name(s) which to inspect. The default None will result in all modules

Returns

A Series of all variables. The MultiIndex indicated which module provides the variable and with which feature name the variable can be obtained.

Return type

pd.Series

atlite.data.cutout_prepare(cutout, features=None, tmpdir=None, overwrite=False)

Prepare all or a selection of features in a cutout.

This function loads the feature data of a cutout, e.g. influx or runoff. When not specifying the feature argument, all available features will be loaded. The function compares the variables which are already included in the cutout with the available variables of the modules specified by the cutout. It detects missing variables and stores them into the netcdf file of the cutout.

Parameters
  • cutout (atlite.Cutout) –

  • features (str/list, optional) – Feature(s) to be prepared. The default slice(None) results in all available features.

  • tmpdir (str/Path, optional) – Directory in which temporary files (for example retrieved ERA5 netcdf files) are stored. If set, the directory will not be deleted and the intermediate files can be examined.

  • overwrite (bool, optional) – Whether to overwrite variables which are already included in the cutout. The default is False.

Returns

cutout – Cutout with prepared data. The variables are stored in cutout.data.

Return type

atlite.Cutout

atlite.data.get_features(cutout, module, features, tmpdir=None)

Load the feature data for a given module.

This get the data for a set of features from a module. All modules in atlite.datasets are allowed.

atlite.data.maybe_remove_tmpdir(func)

Use this wrapper to make tempfile deletion compatible with windows machines.

atlite.data.non_bool_dict(d)

Convert bool to int for netCDF4 storing

Convert

All functions for converting weather data into energy system model data.

atlite.convert.convert_and_aggregate(cutout, convert_func, matrix=None, index=None, layout=None, shapes=None, shapes_crs=4326, per_unit=False, return_capacity=False, capacity_factor=False, show_progress=True, dask_kwargs={}, **convert_kwds)

Convert and aggregate a weather-based renewable generation time-series.

NOTE: Not meant to be used by the user him or herself. Rather it is a gateway function that is called by all the individual time-series generation functions like pv and wind. Thus, all its parameters are also available from these.

Parameters
  • matrix (N x S - xr.DataArray or sp.sparse.csr_matrix or None) – If given, it is used to aggregate the grid cells to buses. N is the number of buses, S the number of spatial coordinates, in the order of cutout.grid.

  • index (pd.Index) – Index of Buses.

  • layout (X x Y - xr.DataArray) – The capacity to be build in each of the grid_cells.

  • shapes (list or pd.Series of shapely.geometry.Polygon) – If given, matrix is constructed as indicatormatrix of the polygons, its index determines the bus index on the time-series.

  • shapes_crs (pyproj.CRS or compatible) – If different to the map crs of the cutout, the shapes are transformed to match cutout.crs (defaults to EPSG:4326).

  • per_unit (boolean) – Returns the time-series in per-unit units, instead of in MW (defaults to False).

  • return_capacity (boolean) – Additionally returns the installed capacity at each bus corresponding to layout (defaults to False).

  • capacity_factor (boolean) – If True, the static capacity factor of the chosen resource for each grid cell is computed.

  • show_progress (boolean, default True) – Whether to show a progress bar.

  • dask_kwargs (dict, default {}) – Dict with keyword arguments passed to dask.compute.

  • convert_func (Function) – Callback like convert_wind, convert_pv

Returns

  • resource (xr.DataArray) – Time-series of renewable generation aggregated to buses, if matrix or equivalents are provided else the total sum of generated energy.

  • units (xr.DataArray (optional)) – The installed units per bus in MW corresponding to layout (only if return_capacity is True).

atlite.convert.convert_soil_temperature(ds)

Return soil temperature (useful for e.g. heat pump T-dependent coefficient of performance).

atlite.convert.convert_temperature(ds)

Return outside temperature (useful for e.g. heat pump T-dependent coefficient of performance).

atlite.convert.convert_wind(ds, turbine)

Convert wind speeds for turbine to wind energy generation.

atlite.convert.heat_demand(cutout, threshold=15.0, a=1.0, constant=0.0, hour_shift=0.0, **params)

Convert outside temperature into daily heat demand using the degree-day approximation.

Since “daily average temperature” means different things in different time zones and since xarray coordinates do not handle time zones gracefully like pd.DateTimeIndex, you can provide an hour_shift to redefine when the day starts.

E.g. for Moscow in winter, hour_shift = 4, for New York in winter, hour_shift = -5

This time shift applies across the entire spatial scope of ds for all times. More fine-grained control will be built in a some point, i.e. space- and time-dependent time zones.

WARNING: Because the original data is provided every month, at the month boundaries there is untidiness if you use a time shift. The resulting xarray will have duplicates in the index for the parts of the day in each month at the boundary. You will have to re-average these based on the number of hours in each month for the duplicated day.

Parameters
  • threshold (float) – Outside temperature in degrees Celsius above which there is no heat demand.

  • a (float) – Linear factor relating heat demand to outside temperature.

  • constant (float) – Constant part of heat demand that does not depend on outside temperature (e.g. due to water heating).

  • hour_shift (float) – Time shift relative to UTC for taking daily average

Note

You can also specify all of the general conversion arguments documented in the convert_and_aggregate function.

atlite.convert.hydro(cutout, plants, hydrobasins, flowspeed=1, weight_with_height=False, show_progress=True, **kwargs)

Compute inflow time-series for plants by aggregating over catchment basins from hydrobasins

Parameters
  • plants (pd.DataFrame) – Run-of-river plants or dams with lon, lat columns.

  • hydrobasins (str|gpd.GeoDataFrame) – Filename or GeoDataFrame of one level of the HydroBASINS dataset.

  • flowspeed (float) – Average speed of water flows to estimate the water travel time from basin to plant (default: 1 m/s).

  • weight_with_height (bool) – Whether surface runoff should be weighted by potential height (probably better for coarser resolution).

  • show_progress (bool) – Whether to display progressbars.

References

[1] Liu, Hailiang, et al. “A validated high-resolution hydro power time-series model for energy systems analysis.” arXiv preprint arXiv:1901.08476 (2019).

[2] Lehner, B., Grill G. (2013): Global river hydrography and network routing: baseline data and new approaches to study the world’s large river systems. Hydrological Processes, 27(15): 2171–2186. Data is available at www.hydrosheds.org.

atlite.convert.maybe_progressbar(ds, show_progress, **kwargs)

Load a xr.dataset with dask arrays either with or without progressbar.

atlite.convert.pv(cutout, panel, orientation, clearsky_model=None, **params)

Convert downward-shortwave, upward-shortwave radiation flux and ambient temperature into a pv generation time-series.

Parameters
  • panel (str or dict) – Panel config dictionary with the parameters for the electrical model in [3]. Alternatively, name of yaml file stored in atlite.config.solarpanel_dir.

  • orientation (str, dict or callback) – Panel orientation can be chosen from either ‘latitude_optimal’, a constant orientation {‘slope’: 0.0, ‘azimuth’: 0.0} or a callback function with the same signature as the callbacks generated by the ‘atlite.pv.orientation.make_*’ functions.

  • clearsky_model (str or None) – Either the ‘simple’ or the ‘enhanced’ Reindl clearsky model. The default choice of None will choose dependending on data availability, since the ‘enhanced’ model also incorporates ambient air temperature and relative humidity.

Returns

pv – Time-series or capacity factors based on additional general conversion arguments.

Return type

xr.DataArray

Note

You can also specify all of the general conversion arguments documented in the convert_and_aggregate function.

References

[1] Soteris A. Kalogirou. Solar Energy Engineering: Processes and Systems, pages 49–117,469–516. Academic Press, 2009. ISBN 0123745012.

[2] D.T. Reindl, W.A. Beckman, and J.A. Duffie. Diffuse fraction correla- tions. Solar Energy, 45(1):1 – 7, 1990.

[3] Hans Georg Beyer, Gerd Heilscher and Stefan Bofinger. A Robust Model for the MPP Performance of Different Types of PV-Modules Applied for the Performance Check of Grid Connected Systems, Freiburg, June 2004. Eurosun (ISES Europe Solar Congress).

atlite.convert.solar_thermal(cutout, orientation={'azimuth': 180.0, 'slope': 45.0}, trigon_model='simple', clearsky_model='simple', c0=0.8, c1=3.0, t_store=80.0, **params)

Convert downward short-wave radiation flux and outside temperature into time series for solar thermal collectors.

Mathematical model and defaults for c0, c1 based on model in [1].

Parameters
  • cutout (cutout) –

  • orientation (dict or str or function) – Panel orientation with slope and azimuth (units of degrees), or ‘latitude_optimal’.

  • trigon_model (str) – Type of trigonometry model

  • clearsky_model (str or None) – Type of clearsky model for diffuse irradiation. Either ‘simple’ or ‘enhanced’.

  • c0 (float) – Parameters for model in [1] (defaults to 0.8 and 3., respectively)

  • c1 (float) – Parameters for model in [1] (defaults to 0.8 and 3., respectively)

  • t_store (float) – Store temperature in degree Celsius

Note

You can also specify all of the general conversion arguments documented in the convert_and_aggregate function.

References

[1] Henning and Palzer, Renewable and Sustainable Energy Reviews 30 (2014) 1003-1018

atlite.convert.wind(cutout, turbine, smooth=False, **params)

Generate wind generation time-series

Extrapolates 10m wind speed with monthly surface roughness to hub height and evaluates the power curve.

Parameters
  • turbine (str or dict) – A turbineconfig dictionary with the keys ‘hub_height’ for the hub height and ‘V’, ‘POW’ defining the power curve. Alternatively a str refering to a local or remote turbine configuration as accepted by atlite.resource.get_windturbineconfig().

  • smooth (bool or dict) – If True smooth power curve with a gaussian kernel as determined for the Danish wind fleet to Delta_v = 1.27 and sigma = 2.29. A dict allows to tune these values.

Note

You can also specify all of the general conversion arguments documented in the convert_and_aggregate function.

References

[1] Andresen G B, Søndergaard A A and Greiner M 2015 Energy 93, Part 1 1074 – 1088. doi:10.1016/j.energy.2015.09.071

Resource

Module for providing access to external ressources, like windturbine or pv panel configurations.

atlite.resource.get_oedb_windturbineconfig(search=None, **search_params)

Download a windturbine configuration from the OEDB database.

Download the configuration of a windturbine model from the OEDB database into the local ‘windturbine_dir’. The OEDB database can be viewed here: https://openenergy-platform.org/dataedit/view/supply/wind_turbine_library (2019-07-22) Only one turbine configuration is downloaded at a time, if the search parameters yield an ambigious result, no data is downloaded.

Parameters
  • search (int|str) – Smart search parameter, if int use as model id, if str look in name or turbine_type

  • **search_params (dict) – Recognized arguments are ‘id’, ‘name’, ‘turbine_type’ and ‘manufacturer’

Returns

turbineconfig – The turbine configuration in the format from ‘atlite.ressource.get_turbineconf(name)’.

Return type

dict

Example

>>> get_oedb_windturbineconfig(10)
{'V': ..., 'POW': ..., ...}
>>> get_oedb_windturbineconfig(name="E-53/800", manufacturer="Enercon")
{'V': ..., 'POW': ..., ...}
atlite.resource.get_solarpanelconfig(panel)

Load the ‘panel’.yaml file from local disk and provide a solar panel dict.

atlite.resource.get_windturbineconfig(turbine)

Load the wind ‘turbine’ configuration.

The configuration can either be one from local storage, then ‘turbine’ is considered part of the file base name ‘<turbine>.yaml’ in config.windturbine_dir. Alternatively the configuration can be downloaded from the Open Energy Database (OEDB), in which case ‘turbine’ is a dictionary used for selecting a turbine from the database.

Parameters

turbine (str) – Name of the local turbine file. Alternatively a dict for selecting a turbine from the Open Energy Database, in this case the key ‘source’ should be contained. For all other key arguments to retrieve the matching turbine, see atlite.resource.download_windturbineconfig() for details.

atlite.resource.windturbine_smooth(turbine, params=None)

Smooth the powercurve in turbine with a gaussian kernel

Parameters
  • turbine (dict) – Turbine config with at least V and POW

  • params (dict) – Allows adjusting fleet availability eta, mean Delta_v and stdev sigma. Defaults to values from Andresen’s paper: 0.95, 1.27 and 2.29, respectively.

Returns

turbine – Turbine config with a smoothed power curve

Return type

dict

References

G. B. Andresen, A. A. Søndergaard, M. Greiner, Validation of Danish wind time series from a new global renewable energy atlas for energy system analysis, Energy 93, Part 1 (2015) 1074–1088.

Wind

Functions for use in conjunction with wind data generation.

atlite.wind.extrapolate_wind_speed(ds, to_height, from_height=None)

Extrapolate the wind speed from a given height above ground to another.

If ds already contains a key refering to wind speeds at the desired to_height, no conversion is done and the wind speeds are directly returned.

Extrapolation of the wind speed follows the logarithmic law as desribed in [1].

Parameters
  • ds (xarray.Dataset) – Dataset containing the wind speed time-series at ‘from_height’ with key ‘wnd{height:d}m’ and the surface orography with key ‘roughness’ at the geographic locations of the wind speeds.

  • from_height (int) – (Optional) Height (m) from which the wind speeds are interpolated to ‘to_height’. If not provided, the closest height to ‘to_height’ is selected.

  • to_height (int|float) – Height (m) to which the wind speeds are extrapolated to.

Returns

da – DataArray containing the extrapolated wind speeds. Name of the DataArray is ‘wnd{to_height:d}’.

Return type

xarray.DataArray

References

[1] Equation (2) in Andresen, G. et al (2015): ‘Validation of Danish wind time series from a new global renewable energy atlas for energy system analysis’.

[2] https://en.wikipedia.org/w/index.php?title=Roughness_length&oldid=862127433, Retrieved 2019-02-15.

GIS

Functions for Geographic Information System.

class atlite.gis.ExclusionContainer(crs=3035, res=100)

Container for exclusion objects and meta data.

add_geometry(geometry, buffer=0, invert=False)

Register a collection of geometries to the ExclusionContainer.

Parameters
  • geometry (str/path/geopandas.GeoDataFrame) – Path to geometries or geometries which to exclude.

  • buffer (float, optional) – Buffer around the excluded areas in units of ExclusionContainer.crs. The default is 0.

  • invert (bool, optional) – Whether to exclude (False) or include (True) the specified areas of the geometries. The default is False.

add_raster(raster, codes=None, buffer=0, invert=False, nodata=255, allow_no_overlap=False, crs=None)

Register a raster to the ExclusionContainer.

Parameters
  • raster (str/rasterio.DatasetReader) – Raster or path to raster which to exclude.

  • codes (int/list/function, optional) – Codes in the raster which to exclude. Can be a callable function which takes the mask (np.array) as argument and performs a elementwise condition (must not change the shape). The default is 1.

  • buffer (int, optional) – Buffer around the excluded areas in units of ExclusionContainer.crs. Use this to create a buffer around the excluded/included area. The default is 0.

  • invert (bool, optional) – Whether to exclude (False) or include (True) the specified areas of the raster. The default is False.

  • allow_no_overlap – Allow that a raster and a shape (for which the raster will be used as a mask) do not overlap. In this case an array with only nodata is returned.

  • crs (rasterio.CRS/EPSG) – CRS of the raster. Specify this if the raster has invalid crs.

property all_closed

Check whether all files in the raster container are closed.

property all_open

Check whether all files in the raster container are open.

open_files()

Open rasters and load geometries.

atlite.gis.compute_availabilitymatrix(cutout, shapes, excluder, nprocesses=None, disable_progressbar=False)

Compute the eligible share within cutout cells in the overlap with shapes.

For parallel calculation (nprocesses not None) the excluder must not be initialized and all raster references must be strings. Otherwise processes are colliding when reading from one common rasterio.DatasetReader.

Parameters
  • cutout (atlite.Cutout) – Cutout which the availability matrix is aligned to.

  • shapes (geopandas.Series/geopandas.DataFrame) – Geometries for which the availabilities are calculated.

  • excluder (atlite.gis.ExclusionContainer) – Container of all meta data or objects which to exclude, i.e. rasters and geometries.

  • nprocesses (int, optional) – Number of processes to use for calculating the matrix. The paralle- lization can heavily boost the calculation speed. The default is None.

  • disable_progressbar (bool, optional) – Disable the progressbar if nprocesses is not None. Then the map function instead of the imap function is used for the multiprocessing pool. This speeds up the calculation.

Returns

availabilities – DataArray of shape (|shapes|, |y|, |x|) containing all the eligible share of cutout cell (x,y) in the overlap with shape i.

Return type

xr.DataArray

Notes

The rasterio (or GDAL) average downsampling returns different results dependent on how the target raster (the cutout raster) is spanned. Either it is spanned from the top left going downwards, e.g. Affine(0.25, 0, 0, 0, -0.25, 50), or starting in the lower left corner and going up, e.g. Affine(0.25, 0, 0, 0, 0.25, 50). Here we stick to the top down version which is why we use cutout.transform_r and flipping the y-axis in the end.

atlite.gis.compute_indicatormatrix(orig, dest, orig_crs=4326, dest_crs=4326)

Compute the indicatormatrix.

The indicatormatrix I[i,j] is a sparse representation of the ratio of the area in orig[j] lying in dest[i], where orig and dest are collections of polygons, i.e.

A value of I[i,j] = 1 indicates that the shape orig[j] is fully contained in shape dest[j].

Note that the polygons must be in the same crs.

Parameters
  • orig (Collection of shapely polygons) –

  • dest (Collection of shapely polygons) –

Returns

I – Indicatormatrix

Return type

sp.sparse.lil_matrix

atlite.gis.get_coords(x, y, time, dx=0.25, dy=0.25, dt='h', **kwargs)

Create an cutout coordinate system on the basis of slices and step sizes.

Parameters
  • x (slice) – Numerical slices with lower and upper bound of the x dimension.

  • y (slice) – Numerical slices with lower and upper bound of the y dimension.

  • time (slice) – Slice with strings with lower and upper bound of the time dimension.

  • dx (float, optional) – Step size of the x coordinate. The default is 0.25.

  • dy (float, optional) – Step size of the y coordinate. The default is 0.25.

  • dt (str, optional) – Frequency of the time coordinate. The default is ‘h’. Valid are all pandas offset aliases.

Returns

ds – Dataset with x, y and time variables, representing the whole coordinate system.

Return type

xarray.Dataset

atlite.gis.maybe_swap_spatial_dims(ds, namex='x', namey='y')

Swap order of spatial dimensions according to atlite concention.

atlite.gis.pad_extent(src, src_transform, dst_transform, src_crs, dst_crs, **kwargs)

Pad the extent of src by an equivalent of one cell of the target raster.

This ensures that the array is large enough to not be treated as nodata in all cells of the destination raster. If src.ndim > 2, the function expects the last two dimensions to be y,x. Additional keyword arguments are used in np.pad().

atlite.gis.padded_transform_and_shape(bounds, res)

Get the (transform, shape) tuple of a raster with resolution res and bounds bounds.

atlite.gis.projected_mask(raster, geom, transform=None, shape=None, crs=None, allow_no_overlap=False, **kwargs)

Load a mask and optionally project it to target resolution and shape.

atlite.gis.regrid(ds, dimx, dimy, **kwargs)

Interpolate Dataset or DataArray ds to a new grid, using rasterio’s reproject facility.

See also: https://mapbox.github.io/rasterio/topics/resampling.html

Parameters
  • ds (xr.Dataset|xr.DataArray) – N-dim data on a spatial grid

  • dimx (pd.Index) – New x-coordinates in destination crs. dimx.name MUST refer to x-coord of ds.

  • dimy (pd.Index) – New y-coordinates in destination crs. dimy.name MUST refer to y-coord of ds.

  • **kwargs – Arguments passed to rio.wrap.reproject; of note: - resampling is one of gis.Resampling.{average,cubic,bilinear,nearest} - src_crs, dst_crs define the different crs (default: EPSG 4326, ie latlong)

atlite.gis.reproject(shapes, p1, p2)

Project a collection of shapes from one crs to another.

atlite.gis.reproject_shapes(shapes, crs1, crs2)

Project a collection of shapes from one crs to another.

atlite.gis.shape_availability(geometry, excluder)

Compute the eligible area in one or more geometries.

Parameters
  • geometry (geopandas.Series) – Geometry of which the eligible area is computed. If the series contains more than one geometry, the eligble area of the combined geometries is computed.

  • excluder (atlite.gis.ExclusionContainer) – Container of all meta data or objects which to exclude, i.e. rasters and geometries.

Returns

  • masked (np.array) – Mask whith eligible raster cells indicated by 1 and excluded cells by 0.

  • transform (rasterion.Affine) – Affine transform of the mask.

atlite.gis.shape_availability_reprojected(geometry, excluder, dst_transform, dst_crs, dst_shape)

Compute and reproject the eligible area of one or more geometries.

The function executes shape_availability and reprojects the calculated mask onto a new raster defined by (dst_transform, dst_crs, dst_shape). Before reprojecting, the function pads the mask such all non-nodata data points are projected in full cells of the target raster. The ensures that all data within the mask are projected correclty (GDAL inherent ‘problem’).

geometrygeopandas.Series

Geometry in which the eligible area is computed. If the series contains more than one geometry, the eligble area of the combined geometries is computed.

excluderatlite.gis.ExclusionContainer

Container of all meta data or objects which to exclude, i.e. rasters and geometries.

dst_transformrasterio.Affine

Transform of the target raster.

dst_crsrasterio.CRS/proj.CRS

CRS of the target raster.

dst_shapetuple

Shape of the target raster.

maskednp.array

Average share of available area per grid cell. 0 indicates excluded, 1 is fully included.

transformrasterio.Affine

Affine transform of the mask.

atlite.gis.spdiag(v)

Create a sparse diagonal matrix from a 1-dimensional array.

Utils

General utility functions for internal use.

class atlite.utils.CachedAttribute(method, name=None, doc=None)

Computes attribute value and caches it in the instance. From the Python Cookbook (Denis Otkidach) This decorator allows you to create a property which can be computed once and accessed many times. Sort of like memoization.

class atlite.utils.arrowdict

A subclass of dict, which allows you to get items in the dict using the attribute syntax!

atlite.utils.migrate_from_cutout_directory(old_cutout_dir, path)

Convert an old style cutout directory to new style netcdf file