1 clim_recal.utils.xarray
clim_recal.utils.xarray
1.1 Classes
Name | Description |
---|---|
XarrayTimeSeriesCalcManager | Manage cacluations over time for .nc files. |
1.1.1 XarrayTimeSeriesCalcManager
clim_recal.utils.xarray.XarrayTimeSeriesCalcManager(self, path=climate_data_mount_path() / 'Raw/UKCP2.2/', save_folder=Path('../docs/assets/cpm-raw-medians'), sub_path=Path('latest'), variables=VariableOptions.cpm_values(), runs=RunOptions.preferred_and_first(), method_name='median', time_dim_name='time', regex=CPM_REGEX, source_folders=list())
Manage cacluations over time for .nc
files.
1.1.1.1 Attributes
Name | Type | Description |
---|---|---|
path | PathLike | Path to aggreate raw files from, following standard cpm hierarchy. |
sub_path | PathLike | None | A subpath to parse, like ‘latest’ for UKCPM2.2 . |
save_folder | PathLike | Where to save resulting summary nc files. |
variables | Sequence[str | VariableOptions] | Which variables to include |
runs | Sequence[str | RunOptions] | Which RunOptions to include. |
method_name | str | Which method to use to aggreate. Must be a standard xarray Dataset method. |
time_dim_name | str | Name of the temporal dim on the joined .nc files. |
regex | str | Check .nc paths to match and then aggregate. |
source_folders | list | List of folders to iterate over, filled via path if None. |
1.1.1.2 Examples
>>> tasmax_hads_1980_raw = getfixture('tasmax_hads_1980_raw')
>>> if not tasmax_hads_1980_raw:
... pytest.skip(mount_or_cache_doctest_skip_message)>>> tmp_save: Path = getfixture('tmp_path') / 'xarray-time-series-summary-manager'
>>> xr_var_managers = XarrayTimeSeriesCalcManager(save_folder=tmp_save)
>>> save_paths: tuple[Path, ...] = xr_var_managers.save_joined_xr_time_series(stop=2, ts_stop=2)
'05' 'tasmax' (1/2)...
Aggregating '06' 'tasmax' (2/2)...
Aggregating >>> pprint(save_paths)
'.../median-tasmax-05.nc'),
(...Path('.../median-tasmax-06.nc'))
...Path(>>> pprint(sorted(tmp_save.iterdir()))
'.../median-tasmax-05.nc'),
[...Path('.../median-tasmax-06.nc')] ...Path(
1.2 Functions
Name | Description |
---|---|
annual_group_xr_time_series | Return and plot a Dataset of time series temporally overlayed. |
apply_geo_func | Apply a Callable to netcdf_source file and export via to_netcdf . |
cftime360_to_date | Convert a Datetime360Day into a date . |
cftime_range_gen | Convert a banded time index a banded standard (Gregorian). |
check_xarray_path_and_var_name | Check and return a T_Dataset instances and included variable name. |
convert_xr_calendar | Convert cpm 360 day time series to a standard 365/366 day time series. |
cpm_check_converted | Check if cpm_xr_time_series is likely already reprojected. |
cpm_reproject_with_standard_calendar | Convert raw cpm_xr_time_series to an 365/366 days and 27700 coords. |
cpm_xarray_to_standard_calendar | Convert a CPM nc file of 360 day calendars to standard calendar. |
crop_xarray | Crop xr_time_series with crop_path shapefile . |
ensure_xr_dataset | Return xr_time_series as a xarray.Dataset instance. |
file_name_to_start_end_dates | Return dates of file name with date_format -date_format structure. |
gdal_warp_wrapper | Execute the gdalwrap function within python . |
generate_360_to_standard | Return array_to_expand 360 days expanded to 365 or 366 days. |
get_cpm_for_coord_alignment | Check if cpm_for_coord_alignment is a Dataset , process if a Path . |
hads_resample_and_reproject | Resample HADs xarray time series to 2.2km. |
interpolate_xr_ts_nans | Interpolate nan values in a Dataset time series. |
join_xr_time_series_var | Join a set of xr_time_series files chronologically. |
plot_xarray | Plot da with **kwargs to path . |
region_crop_file_name | Generate a file name for a regional crop. |
xr_reproject_crs | Reproject source_xr to target_xr coordinate structure. |
1.2.1 annual_group_xr_time_series
clim_recal.utils.xarray.annual_group_xr_time_series(joined_xr_time_series, variable_name, groupby_method='time.dayofyear', method_name='median', time_dim_name='time', regex=CPM_REGEX, start=0, stop=None, step=1, plot_path='annual-aggregated.png', time_stamp=True, **kwargs)
Return and plot a Dataset
of time series temporally overlayed.
1.2.1.1 Parameters
Name | Type | Description | Default |
---|---|---|---|
joined_xr_time_series |
T_Dataset | PathLike | Provide existing Dataset to aggregate and plot (otherwise check path ). |
required |
variable_name |
str | A variable name to specify for data expected. If none that will be extracted and checked from the files directly. | required |
groupby_method |
str | xarray method to aggregate time of xr_time_series . |
'time.dayofyear' |
method_name |
str | xarray method to calculate for plot. |
'median' |
time_dim_name |
str | Name of time dimension in passed files. | 'time' |
regex |
str | A str to filter files within path |
CPM_REGEX |
start |
int | What point to start indexing path results from. |
0 |
stop |
int | None | What point to stop indexing path restuls from. |
None |
step |
int | How many paths to jump between when iterating between stop and start . |
1 |
plot_path |
PathLike | Path to save plot to. |
'annual-aggregated.png' |
time_stamp |
bool | Whether to include a time stamp in the plot_path name. |
True |
**kwargs |
Additional parameters to pass to plot_xarray . |
{} |
1.2.1.2 Examples
>>> tasmax_cpm_1980_raw = getfixture('tasmax_cpm_1980_raw')
>>> if not tasmax_cpm_1980_raw:
... pytest.skip(mount_or_cache_doctest_skip_message)>>> tasmax_cpm_1980_raw_path = getfixture('tasmax_cpm_1980_raw_path').parents[1]
>>> results: T_Dataset = annual_group_xr_time_series(
'tasmax', stop=3)
... tasmax_cpm_1980_raw_path, >>> results
<xarray.Dataset> ...
360)
Dimensions: (dayofyear:
Coordinates:* dayofyear (dayofyear) int64 ... 1 2 3 4 5 6 7 ... 355 356 357 358 359 360
Data variables:9.2 8.95 8.408 8.747 ... 6.387 8.15 9.132 tasmax (dayofyear) float64 ...
1.2.2 apply_geo_func
clim_recal.utils.xarray.apply_geo_func(source_path, func, export_folder, new_path_name_func=None, to_netcdf=True, to_raster=False, export_path_as_output_path_kwarg=False, return_results=False, **kwargs)
Apply a Callable
to netcdf_source
file and export via to_netcdf
.
1.2.2.1 Parameters
Name | Type | Description | Default |
---|---|---|---|
source_path |
PathLike | netcdf file to apply func to. |
required |
func |
ReprojectFuncType | Callable to modify netcdf . |
required |
export_folder |
PathLike | Where to save results. | required |
new_path_name_func |
Callable[[Path], Path] | None | Callabe to generate new path to save to. |
None |
to_netcdf |
bool | Whether to call to_netcdf() method on results Dataset . |
True |
to_raster |
bool | Whether to call rio.to_raster() on results Dataset . |
False |
export_path_as_output_path_kwarg |
bool | Whether to add output_path = export_path to kwargs passed to func . Meant for cases calling gdal_warp_wrapper . |
False |
return_results |
bool | Whether to return results, which would be a Dataset or GDALDataset (the latter if gdal_warp_wrapper is used). |
False |
**kwargs |
Other parameters passed to func call. |
{} |
1.2.2.2 Returns
Type | Description |
---|---|
Either a Path to generated file or converted xarray object. |
1.2.3 cftime360_to_date
clim_recal.utils.xarray.cftime360_to_date(cf_360)
Convert a Datetime360Day
into a date
.
1.2.3.1 Examples
>>> cftime360_to_date(Datetime360Day(1980, 1, 1))
1980, 1, 1) datetime.date(
1.2.4 cftime_range_gen
clim_recal.utils.xarray.cftime_range_gen(time_data_array, **kwargs)
Convert a banded time index a banded standard (Gregorian).
1.2.5 check_xarray_path_and_var_name
clim_recal.utils.xarray.check_xarray_path_and_var_name(xr_time_series, variable_name, ignore_warnings=True)
Check and return a T_Dataset
instances and included variable name.
1.2.6 convert_xr_calendar
clim_recal.utils.xarray.convert_xr_calendar(xr_time_series, align_on=DEFAULT_CALENDAR_ALIGN, calendar=CFCalendarSTANDARD, use_cftime=False, missing_value=np.nan, interpolate_na=False, ensure_output_type_is_dataset=False, interpolate_method=DEFAULT_INTERPOLATION_METHOD, keep_crs=True, keep_attrs=True, limit=1, engine=NETCDF4_XARRAY_ENGINE, check_cftime_cols=None, cftime_range_gen_kwargs=None)
Convert cpm 360 day time series to a standard 365/366 day time series.
1.2.6.1 Notes
Short time examples (like 2 skipped out of 8 days) raises: ValueError("date_range_like was unable to generate a range as the source frequency was not inferable."
)
1.2.6.2 Parameters
Name | Type | Description | Default |
---|---|---|---|
xr_time_series |
T_DataArray | T_Dataset | PathLike | A DataArray or Dataset to convert to calendar time. |
required |
align_on |
ConvertCalendarAlignOptions | Whether and how to align calendar types. |
DEFAULT_CALENDAR_ALIGN |
calendar |
CFCalendar | Type of calendar to convert xr_time_series to. |
CFCalendarSTANDARD |
use_cftime |
bool | Whether to enforce cftime vs datetime64 time format. |
False |
missing_value |
Any | None | Missing value to populate missing date interpolations with. | np.nan |
keep_crs |
bool | Reapply initial Coordinate Reference System (CRS) after time projection. | True |
interpolate_na |
bool | Whether to apply temporal interpolation for missing values. | False |
interpolate_method |
InterpOptions | Which InterpOptions method to apply if interpolate_na is True . |
DEFAULT_INTERPOLATION_METHOD |
keep_attrs |
bool | Whether to keep all attributes on after interpolate_na |
True |
limit |
int | Limit of number of continuous missing day values allowed in interpolate_na . |
1 |
engine |
XArrayEngineType | Which XArrayEngineType to use in parsing files and operations. |
NETCDF4_XARRAY_ENGINE |
extrapolate_fill_value |
If True , then pass fill_value=extrapolate . See: * https://docs.xarray.dev/en/stable/generated/xarray.Dataset.interpolate_na.html * https://docs.scipy.org/doc/scipy/reference/generated/scipy.interpolate.interp1d.html#scipy.interpolate.interp1d |
required | |
check_cftime_cols |
tuple[str] | None | Columns to check cftime format on |
None |
cftime_range_gen_kwargs |
dict[str, Any] | None | Any kwargs to pass to cftime_range_gen |
None |
1.2.6.3 Raises
Type | Description |
---|---|
ValueError | Likely from xarray calling date_range_like . |
1.2.6.4 Returns
Type | Description |
---|---|
T_DataArrayOrSet | Converted xr_time_series to specified calendar with optional interpolation. |
1.2.6.5 Notes
Certain values may fail to interpolate in cases of 360 -> 365/366 (Gregorian) calendar. Examples include projecting CPM data, which is able to fill in measurement values (e.g. tasmax
) but the year
and month_number
variables have nan
values
1.2.6.6 Examples
2 Note a new doctest needs to be written to deal
3 with default year
vs date
parameters
>>> xr_360_to_365_datetime64: T_Dataset = convert_xr_calendar(
="date")
... xarray_spatial_4_years_360_day, align_on>>> xr_360_to_365_datetime64.sel(
=slice("1981-01-30", "1981-02-01"),
... time="Glasgow").day_360
... space<xarray.DataArray 'day_360' (time: 3)>...
Coordinates:* time (time) datetime64[ns] ...1981-01-30 1981-01-31 1981-02-01
<U10 ...'Glasgow'
space >>> xr_360_to_365_datetime64_interp: T_Dataset = convert_xr_calendar(
=True)
... xarray_spatial_4_years_360_day, interpolate_na>>> xr_360_to_365_datetime64_interp.sel(
=slice("1981-01-30", "1981-02-01"),
... time="Glasgow").day_360
... space<xarray.DataArray 'day_360' (time: 3)>...
0.23789282, 0.5356328 , 0.311945 ])
array([
Coordinates:* time (time) datetime64[ns] ...1981-01-30 1981-01-31 1981-02-01
<U10 ...'Glasgow'
space >>> convert_xr_calendar(xarray_spatial_6_days_2_skipped)
Traceback (most recent call last):
...ValueError: `date_range_like` was unable to generate a range as the source frequency was not inferable.
3.0.1 cpm_check_converted
clim_recal.utils.xarray.cpm_check_converted(cpm_xr_time_series)
Check if cpm_xr_time_series
is likely already reprojected.
3.0.1.1 Parameters
Name | Type | Description | Default |
---|---|---|---|
cpm_xr_time_series |
T_Dataset | PathLike | Dataset instance or Path to check. |
required |
3.0.1.2 Returns
Type | Description |
---|---|
True if all of the methics are True , else False |
3.0.2 cpm_reproject_with_standard_calendar
clim_recal.utils.xarray.cpm_reproject_with_standard_calendar(cpm_xr_time_series, variable_name=None, close_temp_paths=True, force=False)
Convert raw cpm_xr_time_series
to an 365/366 days and 27700 coords.
3.0.2.1 Notes
Currently makes UTM coordinate structure
3.0.2.2 Parameters
Name | Type | Description | Default |
---|---|---|---|
cpm_xr_time_series |
T_Dataset | PathLike | Dataset (or path to load as Dataset ) expected to be in raw UKCPM format, with 360 day years and a rotated coordinate system. |
required |
variable_name |
str | None | Name of variable used, usually a measure of climate change like tasmax and tasmin . |
None |
3.0.2.3 Returns
Type | Description |
---|---|
Final xarray Dataset after spatial and temporal changes. |
3.0.2.4 Examples
>>> tasmax_cpm_1980_raw = getfixture('tasmax_cpm_1980_raw')
>>> if not tasmax_cpm_1980_raw:
... pytest.skip(mount_or_cache_doctest_skip_message)>>> tasmax_cpm_1980_365_day: T_Dataset = cpm_reproject_with_standard_calendar(
=tasmax_cpm_1980_raw,
... cpm_xr_time_series="tasmax")
... variable_name.100%...
Warp: ...nc ...100%...
Translate: ...tif ..>>> tasmax_cpm_1980_365_day
<xarray.Dataset>
365, x: 493, y: 607)
Dimensions: (time:
Coordinates:* time (time) datetime64[ns]...
* x (x) float64...
* y (y) float64...
|S1...
transverse_mercator
spatial_ref int64...
Data variables:
tasmax (time, y, x) float32...12/18)
Attributes: (
...>>> tasmax_cpm_1980_365_day.dims
'time': 365, 'x': 493, 'y': 607}) Frozen...({
3.0.3 cpm_xarray_to_standard_calendar
clim_recal.utils.xarray.cpm_xarray_to_standard_calendar(cpm_xr_time_series, include_bnds_index=False)
Convert a CPM nc
file of 360 day calendars to standard calendar.
3.0.3.1 Parameters
Name | Type | Description | Default |
---|---|---|---|
cpm_xr_time_series |
T_Dataset | PathLike | A raw xarray of the form provided by CPM. |
required |
include_bnds_index |
bool | Whether to fix bnds indexing in returned Dataset . |
False |
3.0.3.2 Returns
Type | Description |
---|---|
Dataset calendar converted to standard (Gregorian). |
3.0.4 crop_xarray
clim_recal.utils.xarray.crop_xarray(xr_time_series, crop_box, **kwargs)
Crop xr_time_series
with crop_path
shapefile
.
3.0.4.1 Parameters
Name | Type | Description | Default |
---|---|---|---|
xr_time_series |
T_Dataset | PathLike | Dataset or path to netcdf file to load and crop. |
required |
crop_box |
BoundingBoxCoords | Instance of BoundingBoxCoords with coords. |
required |
kwargs |
Any additional parameters to pass to clip |
{} |
3.0.4.2 Returns
Type | Description |
---|---|
T_Dataset | Spatially cropped xr_time_series Dataset with final_crs spatial coords. |
3.0.4.3 Examples
>>> from clim_recal.utils.data import GlasgowCoordsEPSG27700
>>> from numpy.testing import assert_allclose
>>> tasmax_cpm_1980_raw = getfixture('tasmax_cpm_1980_raw')
>>> if not tasmax_cpm_1980_raw:
... pytest.skip(mount_or_cache_doctest_skip_message)>>> tasmax_cpm_1980_365_day: T_Dataset = cpm_reproject_with_standard_calendar(
=tasmax_cpm_1980_raw,
... cpm_xr_time_series="tasmax")
... variable_name>>> cropped = crop_xarray(
... tasmax_cpm_1980_365_day,=GlasgowCoordsEPSG27700)
... crop_box>>> assert_allclose(cropped.rio.bounds(),
... GlasgowCoordsEPSG27700.as_rioxarray_tuple(),=.01)
... rtol>>> tasmax_cpm_1980_365_day.sizes
'x': 529, 'y': 653, 'time': 365})
Frozen({>>> cropped.sizes
'x': 10, 'y': 8, 'time': 365}) Frozen({
3.0.5 ensure_xr_dataset
clim_recal.utils.xarray.ensure_xr_dataset(xr_time_series, default_name='to_convert')
Return xr_time_series
as a xarray.Dataset
instance.
3.0.5.1 Parameters
Name | Type | Description | Default |
---|---|---|---|
xr_time_series |
T_DataArrayOrSet | Instance to check and if necessary to convert to Dataset . |
required |
default_name |
Name to give returned Dataset if xr_time_series.name is empty. |
'to_convert' |
3.0.5.2 Returns
Type | Description |
---|---|
T_Dataset | Converted (or original) Dataset . |
3.0.5.3 Examples
>>> ensure_xr_dataset(xarray_spatial_4_days)
<xarray.Dataset>...
5, space: 3)
Dimensions: (time:
Coordinates:* time (time) datetime64[ns] ...1980-11-30 1980-12-01 ... 1980-12-04
* space (space) <U10 ...'Glasgow' 'Manchester' 'London'
Data variables:.0.5488 0.7152 ... 0.9256 0.07104 xa_template (time, space) float64 ..
3.0.6 file_name_to_start_end_dates
clim_recal.utils.xarray.file_name_to_start_end_dates(path, date_format=CLI_DATE_FORMAT_STR)
Return dates of file name with date_format
-date_format
structure.
3.0.6.1 Parameters
Name | Type | Description | Default |
---|---|---|---|
path |
PathLike | Path to file | required |
date_format |
str | Format of date for strptime |
CLI_DATE_FORMAT_STR |
3.0.6.2 Examples
The examples below are meant to demonstrate usage, and the significance of when the last date is included or not by default.
>>> from .core import date_range_generator
>>> tif_365_path: Path = (Path('some') /
'folder' /
... 'pr_rcp85_land-cpm_uk_2.2km_06_day_20761201-20771130.tif')
... >>> start_date, end_date = file_name_to_start_end_dates(tif_365_path)
>>> start_date
2076, 12, 1, 0, 0)
datetime.datetime(>>> end_date
2077, 11, 30, 0, 0)
datetime.datetime(>>> dates: tuple[date, ...] = tuple(
=start_date,
... date_range_generator(start_date=end_date,
... end_date=True))
... inclusive>>> dates[:3]
2076, 12, 1, 0, 0),
(datetime.datetime(2076, 12, 2, 0, 0),
datetime.datetime(2076, 12, 3, 0, 0))
datetime.datetime(>>> len(dates)
365
>>> tif_366_path: Path = (Path('some') /
'folder' /
... 'pr_rcp85_land-cpm_uk_2.2km_06_day_20791201-20801130.tif')
... >>> from pandas import date_range
>>> dates = date_range(*file_name_to_start_end_dates(tif_366_path))
>>> len(dates)
366
3.0.7 gdal_warp_wrapper
clim_recal.utils.xarray.gdal_warp_wrapper(input_path, output_path, output_crs=BRITISH_NATIONAL_GRID_EPSG, output_x_resolution=None, output_y_resolution=None, copy_metadata=True, return_path=True, format=GDALNetCDFFormatStr, multithread=True, warp_dict_options=DEFAULT_WARP_DICT_OPTIONS, use_tqdm_progress_bar=True, tqdm_file_name_chars=TQDM_FILE_NAME_PRINT_CHARS_INDEX, resampling_method=None, supress_warnings=True, **kwargs)
Execute the gdalwrap
function within python
.
This is following a script in the bash/
folder that uses this programme:
f=$1 # The first argument is the file to reproject
fn=${f/Raw/Reprojected_infill} # Replace Raw with Reprojected_infill in the filename
folder=`dirname $fn` # Get the folder name
mkdir -p $folder # Create the folder if it doesn't exist
gdalwarp -t_srs 'EPSG:27700' -tr 2200 2200 -r near -overwrite $f "${fn%.nc}.tif" # Reproject the file
3.0.7.1 Parameters
Name | Type | Description | Default |
---|---|---|---|
input_path |
PathLike | Path with CPRUK files to resample. srcDSOrSrcDSTab in Warp . |
required |
output_path |
PathLike | Path to save resampled input_path file(s) to destNameOrDestDS in Warp . |
required |
output_crs |
str | Coordinate system to convert input_path file(s) to. dstSRS in WarpOptions . |
BRITISH_NATIONAL_GRID_EPSG |
format |
GDALFormatsType | str | None | Format to convert input_path to in output_path . |
GDALNetCDFFormatStr |
output_x_resolution |
int | None | Resolution of x cordinates to convert input_path file(s) to. xRes in WarpOptions . |
None |
output_y_resolution |
int | None | Resolution of y cordinates to convert input_path file(s) to. yRes in WarpOptions . |
None |
copy_metadata |
bool | Whether to copy metadata when possible. | True |
return_path |
bool | Return the resulting path if True , else the new GDALDataset . |
True |
format |
GDALFormatsType | str | None | Format to write new file to. | GDALNetCDFFormatStr |
multithread |
bool | Whether to use multithread to speed up calculations. |
True |
kwargs |
Any additional parameters to pass to WarpOption . |
{} |
3.0.8 generate_360_to_standard
clim_recal.utils.xarray.generate_360_to_standard(array_to_expand)
Return array_to_expand
360 days expanded to 365 or 366 days.
This may be dropped if cpm_reproject_with_standard_calendar
is successful.
3.0.9 get_cpm_for_coord_alignment
clim_recal.utils.xarray.get_cpm_for_coord_alignment(cpm_for_coord_alignment, skip_reproject=False, cpm_regex=CPM_REGEX)
Check if cpm_for_coord_alignment
is a Dataset
, process if a Path
.
3.0.9.1 Parameters
Name | Type | Description | Default |
---|---|---|---|
cpm_for_coord_alignment |
PathLike | T_Dataset | None | Either a Path or a file or folder with a cpm file to align to or a xarray.Dataset . If a folder, the first file matching cpm_regex will be used. It will then be processed via cpm_reproject_with_standard_calendar for comparability and use alongside cpm files. |
required |
skip_reproject |
bool | Whether to skip calling cpm_reproject_with_standard_calendar . |
False |
cpm_regex |
str | A regular expression to filter suitable files if cpm_for_coord_alignment is a folder Path . |
CPM_REGEX |
3.0.9.2 Returns
Type | Description |
---|---|
An xarray.Dataset coordinate structure to align HADs coordinates. |
3.0.10 hads_resample_and_reproject
clim_recal.utils.xarray.hads_resample_and_reproject(hads_xr_time_series, variable_name, cpm_to_match, cpm_to_match_func=cpm_reproject_with_standard_calendar, x_dim_name=HADS_RAW_X_COLUMN_NAME, y_dim_name=HADS_RAW_Y_COLUMN_NAME)
Resample HADs
xarray
time series to 2.2km.
3.0.11 interpolate_xr_ts_nans
clim_recal.utils.xarray.interpolate_xr_ts_nans(xr_ts, original_xr_ts=None, check_cftime_cols=None, interpolate_method=DEFAULT_INTERPOLATION_METHOD, keep_crs=True, keep_attrs=True, limit=1, cftime_range_gen_kwargs=None, **kwargs)
Interpolate nan
values in a Dataset
time series.
3.0.11.1 Notes
For details and details of keep_attrs
, limit
and **kwargs
parameters see: https://docs.xarray.dev/en/stable/generated/xarray.DataArray.interpolate_na.html
3.0.11.2 Parameters
Name | Type | Description | Default |
---|---|---|---|
xr_ts |
T_Dataset | Dataset to interpolate via interpolate_na . Requires a time coordinate. |
required |
original_xr_ts |
T_Dataset | None | A Dataset to compare the conversion process with. If not provided, set to the original xr_ts as a reference. |
None |
check_cftime_cols |
tuple[str] | None | tuple of column names in a cftime format to check. |
None |
interpolate_method |
InterpOptions | Which of the xarray interpolation methods to use. |
DEFAULT_INTERPOLATION_METHOD |
keep_crs |
bool | Whether to ensure the original crs is kept via rio.write_crs . |
True |
keep_attrs |
bool | Passed to keep_attrs in interpolate_na . See Notes. |
True |
limit |
int | How many nan are allowed either side of data point to interpolate. See Notes. |
1 |
cftime_range_gen_kwargs |
dict[str, Any] | None | Any cftime_range_gen arguments to use with check_cftime_cols calls. |
None |
3.0.11.3 Returns
Type | Description |
---|---|
Dataset where xr_ts nan values are iterpolated with respect to the time coordinate. |
3.0.12 join_xr_time_series_var
clim_recal.utils.xarray.join_xr_time_series_var(path, variable_name=None, method_name='median', time_dim_name='time', regex=CPM_REGEX, start=0, stop=None, step=1)
Join a set of xr_time_series files chronologically.
3.0.12.1 Parameters
Name | Type | Description | Default |
---|---|---|---|
path |
PathLike | Path to collect files to process from, filtered via regex . |
required |
variable_name |
str | None | A variable name to specify for data expected. If none that will be extracted and checked from the files directly. | None |
method_name |
str | What method to use to summarise each time point results. | 'median' |
time_dim_name |
str | Name of time dimension in passed files. | 'time' |
regex |
str | A str to filter files within path |
CPM_REGEX |
start |
int | What point to start indexing path results from. |
0 |
stop |
int | None | What point to stop indexing path restuls from. |
None |
step |
int | How many paths to jump between when iterating between stop and start . |
1 |
3.0.12.2 Examples
>>> tasmax_cpm_1980_raw = getfixture('tasmax_cpm_1980_raw')
>>> if not tasmax_cpm_1980_raw:
... pytest.skip(mount_or_cache_doctest_skip_message)>>> tasmax_cpm_1980_raw_path = getfixture('tasmax_cpm_1980_raw_path').parents[1]
>>> results: T_Dataset = join_xr_time_series_var(
... tasmax_cpm_1980_raw_path,'tasmax', stop=3)
... >>> results
<xarray.Dataset> ...
1080)
Dimensions: (time:
Coordinates:* time (time) object ... 1980-12-01 12:00:00 ... 1983-11-30 12:00:00
Data variables:8.221 6.716 6.499 7.194 ... 8.456 8.153 5.501 tasmax (time) float64 ...
3.0.13 plot_xarray
clim_recal.utils.xarray.plot_xarray(da, path=None, time_stamp=False, return_path=True, **kwargs)
Plot da
with **kwargs
to path
.
3.0.13.1 Parameters
Name | Type | Description | Default |
---|---|---|---|
da |
T_DataArrayOrSet | xarray objec to plot. |
required |
path |
PathLike | None | File to write plot to. | None |
time_stamp |
bool | Whather to add a datetime str of time of writing in file name. |
False |
kwargs |
Additional parameters to pass to plot . |
{} |
3.0.13.2 Examples
>>> example_path: Path = (
'tmp_path') / 'test-path/example.png')
... getfixture(>>> image_path: Path = plot_xarray(
... xarray_spatial_4_days, example_path)>>> example_path == image_path
True
>>> example_time_stamped: Path = (
/ 'example-stamped.png')
... example_path.parent >>> timed_image_path: Path = plot_xarray(
... xarray_spatial_4_days, example_time_stamped,=True)
... time_stamp>>> example_time_stamped != timed_image_path
True
>>> print(timed_image_path)
/.../test-path/example-stamped_...-...-..._...png
3.0.14 region_crop_file_name
clim_recal.utils.xarray.region_crop_file_name(crop_region, file_name)
Generate a file name for a regional crop.
3.0.14.1 Parameters
Name | Type | Description | Default |
---|---|---|---|
crop_region |
str | RegionOptions | None | Region name to include in cropped file name. | required |
file_name |
PathLike | File name to add crop_region name to. |
required |
3.0.14.2 Examples
>>> region_crop_file_name(
'Glasgow',
... 'tasmax.nc')
... 'crop_Glasgow_tasmax.nc'
>>> region_crop_file_name(
'Glasgow',
... 'tasmax_hadukgrid_uk_2_2km_day_19800601-19800630.nc')
... 'crop_Glasgow_tasmax_hads_19800601-19800630.nc'
>>> region_crop_file_name(
'Glasgow',
... 'tasmax_rcp85_land-cpm_uk_2.2km_05_day_std_year_19861201-19871130.nc')
... 'crop_Glasgow_tasmax_cpm_05_19861201-19871130.nc'
3.0.15 xr_reproject_crs
clim_recal.utils.xarray.xr_reproject_crs(xr_time_series, x_dim_name=CPM_RAW_X_COLUMN_NAME, y_dim_name=CPM_RAW_Y_COLUMN_NAME, time_dim_name=TIME_COLUMN_NAME, variable_name=None, final_crs=BRITISH_NATIONAL_GRID_EPSG, match_xr_time_series=None, match_xr_time_series_load_func=None, match_xr_time_series_load_kwargs=None, resampling_method=DEFAULT_RESAMPLING_METHOD, nodata=np.nan, **kwargs)
Reproject source_xr
to target_xr
coordinate structure.
3.0.15.1 Parameters
Name | Type | Description | Default |
---|---|---|---|
xr_time_series |
T_Dataset | PathLike | Dataset or PathLike to load and reproject. |
required |
x_dim_name |
str | str name of x spatial dimension in xr_time_series . Default matches CPM UK projections. |
CPM_RAW_X_COLUMN_NAME |
y_dim_name |
str | str name of y spatial dimension in xr_time_series . Default matches CPM UK projections. |
CPM_RAW_Y_COLUMN_NAME |
time_dim_name |
str | str name of time dimension in xr_time_series . |
TIME_COLUMN_NAME |
variable_name |
str | None | Name of datset to apply projection to within xr_time_series . Inferred if None assuming only one data_var attribute. |
None |
final_crs |
str | Coordinate system str to project xr_time_series to. |
BRITISH_NATIONAL_GRID_EPSG |
resampling_method |
Resampling | rasterio resampling method to apply. |
DEFAULT_RESAMPLING_METHOD |
3.0.15.2 Examples
>>> tasmax_hads_1980_raw = getfixture('tasmax_hads_1980_raw')
>>> if not tasmax_hads_1980_raw:
... pytest.skip(mount_or_cache_doctest_skip_message)>>> tasmax_hads_1980_raw.dims
'time': 31,
FrozenMappingWarningOnValuesAccess({'projection_y_coordinate': 1450,
'projection_x_coordinate': 900,
'bnds': 2})
>>> tasmax_hads_2_2km: T_Dataset = xr_reproject_crs(
... tasmax_hads_1980_raw,="tasmax",
... variable_name=HADS_RAW_X_COLUMN_NAME,
... x_dim_name=HADS_RAW_Y_COLUMN_NAME,
... y_dim_name=(CPM_RESOLUTION_METERS,
... resolution
... CPM_RESOLUTION_METERS),)>>> tasmax_hads_2_2km.dims
'x': 410,
FrozenMappingWarningOnValuesAccess({'y': 660,
'time': 31})