xhydro.indicators package

Indicators submodule.

Submodules

xhydro.indicators.generic module

Module to compute indicators using xclim’s build_indicator_module_from_yaml.

xhydro.indicators.generic.compute_indicators(ds: Dataset, indicators: str | PathLike | Sequence[Indicator] | Sequence[tuple[str, Indicator]] | ModuleType, *, periods: list[str] | list[list[str]] | None = None, restrict_years: bool = True, to_level: str | None = 'indicators', rechunk_input: bool = False) dict[source]

Calculate variables and indicators based on a YAML call to xclim.

The function cuts the output to be the same years as the inputs. Hence, if an indicator creates a timestep outside the original year range (e.g. the first DJF for QS-DEC), it will not appear in the output.

Parameters

dsxr.Dataset

Dataset to use for the indicators.

indicatorsstr | os.PathLike | Sequence[Indicator] | Sequence[tuple[str, Indicator]] | ModuleType

Path to a YAML file that instructs on how to calculate missing variables. Can also be only the « stem », if translations and custom indices are implemented. Can be the indicator module directly, or a sequence of indicators or a sequence of tuples (indicator name, indicator) as returned by iter_indicators().

periodslist of str or list of lists of str, optional

Either [start, end] or list of [start, end] of continuous periods over which to compute the indicators. This is needed when the time axis of ds contains some jumps in time. If None, the dataset will be considered continuous.

restrict_yearsbool

If True, cut the time axis to be within the same years as the input. This is mostly useful for frequencies that do not start in January, such as QS-DEC. In that instance, xclim would start on previous_year-12-01 (DJF), with a NaN. restrict_years will cut that first timestep. This should have no effect on YS and MS indicators.

to_levelstr, optional

The processing level to assign to the output. If None, the processing level of the inputs is preserved.

rechunk_inputbool

If True, the dataset is rechunked with flox.xarray.rechunk_for_blockwise() according to the resampling frequency of the indicator. Each rechunking is done once per frequency with xscen.utils.rechunk_for_resample().

Returns

dict

Dictionary (keys = timedeltas) with indicators separated by temporal resolution.

See Also

xclim.indicators, xclim.core.indicator.build_indicator_module_from_yaml

xhydro.indicators.generic.compute_volume(da: DataArray, *, out_units: str = 'm3', attrs: dict | None = None) DataArray[source]

Compute the volume of water from a streamflow variable, keeping the same frequency.

Parameters

daxr.DataArray

Streamflow variable.

out_unitsstr

Output units. Defaults to « m3 ».

attrsdict, optional

Attributes to add to the output variable. Default attributes for « long_name », « units », « cell_methods » and « description » will be added if not provided.

Returns

xr.DataArray

Volume of water.

xhydro.indicators.generic.get_yearly_op(ds, op, *, input_var: str = 'streamflow', window: int = 1, timeargs: dict | None = None, missing: str = 'skip', missing_options: dict | None = None, interpolate_na: bool = False) Dataset[source]

Compute yearly operations on a variable.

Parameters

dsxr.Dataset

Dataset containing the variable to compute the operation on.

opstr

Operation to compute. One of [« max », « min », « mean », « sum »].

input_varstr

Name of the input variable. Defaults to « streamflow ».

windowint

Size of the rolling window. A « mean » operation is performed on the rolling window before the call to xclim. This parameter cannot be used with the « sum » operation.

timeargsdict, optional

Dictionary of time arguments for the operation. Keys are the name of the period that will be added to the results (e.g. « winter », « summer », « annual »). Values are up to two dictionaries, with both being optional. The first is {“freq”: str}, where str is a frequency supported by xarray (e.g. « YS », « YS-JAN », « YS-DEC »). It needs to be a yearly frequency. Defaults to « YS-JAN ». The second is an indexer as supported by xclim.core.calendar.select_time(). Defaults to {}, which means the whole year. See xclim.core.calendar.select_time() for more information. Examples: {« winter »: {« freq »: « YS-DEC », « date_bounds »: [« 12-01 », « 02-28 »]}}, {« jan »: {« freq »: « YS », « month »: 1}}, {« annual »: {}}.

missingstr

How to handle missing values. One of « skip », « any », « at_least_n », « pct », « wmo ». See xclim.core.missing() for more information.

missing_optionsdict, optional

Dictionary of options for the missing values” method. See xclim.core.missing() for more information.

interpolate_nabool

Whether to interpolate missing values before computing the operation. Only used with the « sum » operation. Defaults to False.

Returns

xr.Dataset

Dataset containing the computed operations, with one variable per indexer. The name of the variable follows the pattern {input_var}{window}_{op}_{indexer}.

Notes

If you want to perform a frequency analysis on a frequency that is finer than annual, simply use multiple timeargs (e.g. 1 per month) to create multiple distinct variables.

xhydro.indicators.pmp module

Module to compute Probable Maximum Precipitation.

xhydro.indicators.pmp.compute_spring_and_summer_mask(snw: DataArray, thresh: str = '1 cm', window_wint_start: int = 14, window_wint_end: int = 45, spr_start: int = 60, spr_end: int = 30, freq: str = 'YS-JUL')[source]

Create a mask that defines the spring and summer seasons based on the snow water equivalent.

Parameters

snwxarray.DataArray

Snow water equivalent. Must be a length (e.g. « mm ») or a mass (e.g. « kg m-2 »).

threshQuantified

Threshold snow thickness to define the start and end of winter.

window_wint_startint

Minimum number of days with snow depth above or equal to threshold to define the start of winter.

window_wint_endint

Maximum number of days with snow depth below or equal to threshold to define the end of winter.

spr_startint

Number of days before the end of winter to define the start of spring.

spr_endint

Number of days after the end of winter to define the end of spring.

freqstr

Frequency of the time axis (annual frequency). Defaults to « YS-JUL ».

Returns

xr.Dataset

Dataset with two DataArrays (mask_spring and mask_summer), with values of 1 where the spring and summer criteria are met and 0 where they are not.

xhydro.indicators.pmp.major_precipitation_events(da: DataArray, windows: list[int], quantile: float = 0.9)[source]

Get precipitation events that exceed a given quantile for a given time step accumulation. Based on Clavet-Gaumont et al. (2017).

Parameters

daxr.DataArray

DataArray containing the precipitation values.

windowslist of int

List of the number of time steps to accumulate precipitation.

quantilefloat

Threshold that limits the events to those that exceed this quantile. Defaults to 0.9.

Returns

xr.DataArray

Masked DataArray containing the major precipitation events.

Notes

https://doi.org/10.1016/j.ejrh.2017.07.003

xhydro.indicators.pmp.precipitable_water(hus: DataArray, zg: DataArray, orog: DataArray, windows: list[int] | None = None, beta_func: bool = True, add_pre_lay: bool = False)[source]

Compute precipitable water based on Clavet-Gaumont et al. (2017) and Rousseau et al (2014).

Parameters

husxr.DataArray

Specific humidity. Must have a pressure level (plev) dimension.

zgxr.DataArray

Geopotential height. Must have a pressure level (plev) dimension.

orogxr.DataArray

Surface altitude.

windowslist of int, optional

Duration of the event in time steps. Defaults to [1]. Note that an additional time step will be added to the window size to account for antecedent conditions.

beta_funcbool, optional

If True, use the beta function proposed by Boer (1982) to get the pressure layers above the topography. If False, the surface altitude is used as the lower boundary of the atmosphere assuming that the surface altitude and the geopotential height are virtually identical at low altitudes.

add_pre_laybool, optional

If True, add the pressure layer between the surface and the lowest pressure level (e.g., at sea level). If False, only the pressure layers between the lowest and highest pressure level are considered.

Returns

xr.DataArray

Precipitable water.

Notes

1) The precipitable water of an event is defined as the maximum precipitable water found during the entire duration of the event, extending up to one time step before the start of the event. Thus, the rolling operation made using windows is a maximum, not a sum.

  1. beta_func = True and add_pre_lay = False follow Clavet-Gaumont et al. (2017) and Rousseau et al (2014).

  2. https://doi.org/10.1016/j.ejrh.2017.07.003 https://doi.org/10.1016/j.jhydrol.2014.10.053 https://doi.org/10.1175/1520-0493(1982)110<1801:DEIIC>2.0.CO;2

xhydro.indicators.pmp.precipitable_water_100y(pw: DataArray, dist: str, method: str, mf: float = 0.2, rebuild_time: bool = True)[source]

Compute the 100-year return period of precipitable water for each month. Based on Clavet-Gaumont et al. (2017).

Parameters

pwxr.DataArray

DataArray containing the precipitable water.

diststr

Probability distributions.

method{« ML » or « MLE », « MM », « PWM », « APP »}

Fitting method, either maximum likelihood (ML or MLE), method of moments (MM) or approximate method (APP). Can also be the probability weighted moments (PWM), also called L-Moments, if a compatible dist object is passed.

mffloat

Maximum majoration factor of the 100-year event compared to the maximum of the timeseries. Used as an upper limit for the frequency analysis.

rebuild_timebool

Whether or not to reconstruct a timeseries with the same time dimensions as pw.

Returns

xr.DataArray

Precipitable water for a 100-year return period.

Notes

https://doi.org/10.1016/j.ejrh.2017.07.003

xhydro.indicators.pmp.spatial_average_storm_configurations(da, radius)[source]

Compute the spatial average for different storm configurations proposed by Clavet-Gaumont et al. (2017).

Parameters

daxr.DataArray

DataArray containing the precipitation values.

radiusfloat

Maximum radius of the storm.

Returns

xr.DataSet

DataSet contaning the spatial averages for all the storm configurations. The y and x coordinates indicate the location of the storm. This location is determined by n//2, where n is the total number of cells for both the rows and columns in the configuration, and // represents floor division.

Notes

https://doi.org/10.1016/j.ejrh.2017.07.003.