easyclimate.interp¶
Submodules¶
Functions¶
|
Regridding regular or lat-lon grid data. |
|
Computes Barnes interpolation for observation values |
|
Rust-accelerated Barnes interpolation wrapper (lon/lat in degrees). |
|
Computes Barnes interpolation on the unit sphere \(S^2\) (spherical metric) |
|
Rust-accelerated Barnes interpolation on the unit sphere \(S^2\) (spherical metric). |
|
Interpolating variables from the model levels (e.g., sigma vertical coordinate) to the standard pressure levels by 1-D function. |
|
Interpolating variables from the pressure levels (e.g., sigma vertical coordinate) to the altitude levels by 1-D function. |
Interpolating variables from the pressure levels (e.g., sigma vertical coordinate) to the altitude levels using easyclimate-rust. |
|
|
Interpolate values from a regular grid/mesh to specific point locations. |
|
Interpolate gridded data to specific point locations with time series preservation. |
|
Interpolate atmospheric data from Community Atmosphere Model (CAM) hybrid sigma-pressure coordinates to pressure levels. |
|
Interpolate atmospheric data from Community Atmosphere Model (CAM) hybrid sigma-pressure |
|
Interpolates data at multidimensional pressure levels to constant pressure coordinates and uses an ECMWF formulation to extrapolate values below ground. |
Package Contents¶
- easyclimate.interp.interp_mesh2mesh(data_input: xarray.DataArray | xarray.Dataset, target_grid: xarray.DataArray | xarray.Dataset, lon_dim: str = 'lon', lat_dim: str = 'lat', method: Literal['linear', 'nearest', 'cubic', 'conservative'] = 'linear')¶
Regridding regular or lat-lon grid data.
Parameters¶
- data_input
xarray.DataArrayorxarray.Dataset The spatio-temporal data to be calculated.
- target_grid:
xarray.DataArrayorxarray.Dataset Target grid to be regridding.
xarray.DataArrayversion sampletarget_grid = xr.DataArray( dims=('lat', 'lon'), coords={'lat': np.arange(-89, 89, 3) + 1 / 1.0, 'lon': np.arange(-180, 180, 3) + 1 / 1.0} )
xarray.Datasetversion sampletarget_grid = xr.Dataset() target_grid['lat'] = np.arange(-89, 89, 3) + 1 / 1.0 target_grid['lon'] = np.arange(-180, 180, 3) + 1 / 1.0
- lon_dim:
str, default: lon. Longitude coordinate dimension name. By default extracting is applied over the lon dimension.
- lat_dim:
str, default: lat. Latitude coordinate dimension name. By default extracting is applied over the lat dimension.
- method:
str, default: linear. The methods of regridding.
linear: linear, bilinear, or higher dimensional linear interpolation.
nearest: nearest-neighbor regridding.
cubic: cubic spline regridding.
conservative: conservative regridding.
Reference¶
- data_input
- easyclimate.interp.interp_spatial_barnes(data: pandas.DataFrame, var_name: str, *, lon_dim: str = 'lon', lat_dim: str = 'lat', grid_res_deg: float = 0.25, auto_domain: bool = True, buffer_deg: float = 1.0, point: tuple[float, float] | None = None, grid_x: float | None = None, grid_y: float | None = None, sigma_deg: float | None = None, influence_radius_deg: float | None = None, influence_k_sigma: float = 4.0, num_iter: int = 2, method: str = 'optimized_convolution', max_dist_sigma: float | None = None, min_weight: float = 1e-06, mask_radius_deg: float | None = None, missing_sentinels: tuple[float, Ellipsis] = (9999.9, 9999.0, -9999.0, -9999.9))¶
Computes Barnes interpolation for observation values
var_namesampled at irregular locations indata(lon/lat in degrees), using Gaussian weights with width parametersigma_degand returning a regular lon/lat grid as anxarray.DataArray.Barnes interpolation is widely used in meteorology and geosciences to remodel irregular point observations into a smooth gridded field. It can be written as
\[\begin{split}f(\\boldsymbol{x})=\\frac{\\sum_{k=1}^N f_k\\cdot w_k(\\boldsymbol{x})}{\\sum_{k=1}^N w_k(\\boldsymbol{x})}\end{split}\]with Gaussian weights
\[\begin{split}w_k(\\boldsymbol{x})=\\text{e}^{-\\frac{1}{2\\sigma^2}\\left|x-\\boldsymbol{x}_k\\right|^2}\end{split}\]Naive computation of Barnes interpolation leads to an algorithmic complexity of \(O(N \\times W \\times H)\), where \(N\) is the number of sample points and \(W \\times H\) the size of the underlying grid.
For sufficiently large \(n\) (in general in the range from 3 to 6) a good approximation of Barnes interpolation with a reduced complexity \(O(N + W \\times H)\) can be obtained by the convolutional expression
\[\begin{split}f(\\boldsymbol{x})\\approx \\frac{ (\\sum_{k=1}^{N}f_k\\cdot\\delta_{\\boldsymbol{x}_k}) * ( r_n^{*n[x]}(x)\\cdot r_n^{*n[y]}(y) ) }{ ( \\sum_{k=1}^{N} \\delta_{\\boldsymbol{x}_k} ) * ( r_{n}^{*n[x]}(x)\\cdot r_{n}^{*n[y]}(y) ) }\end{split}\]where \(\\delta\) is the Dirac impulse function and \(r(.)\) an elementary rectangular function of a specific length that depends on \(\\sigma\) and \(n\).
This wrapper is degree-based:
grid_res_deg/sigma_deg/influence_radius_deg/mask_radius_degare all in degrees.Internally converts to the fast-barnes “grid units” used by the core implementation:
\[\begin{split}\\sigma_{\\text{grid}} = \\sigma_{\deg} / \\Delta_{\\deg}, \\qquad \\text{max_dist_sigma} = R_{\\deg} / \\sigma_{\\deg}.\end{split}\]If
sigma_degis not provided, it is estimated from station density using a simple spacing heuristic.Optional masking keeps only grid points that have at least one station within
mask_radius_degusing a radius-search KD-tree.
Parameters¶
- data
pandas.DataFrame Input table containing at least
lon_dim,lat_dimandvar_namecolumns. Any rows with NaN/Inf in these columns are dropped after cleaning.There should be a similar structure as follows
lon
lat
qff
-3.73
56.33
995.1
2.64
47.05
1012.5
…
…
…
Note
Data points should contain longitude (lon), latitude (lat) and data variables (the above data variable name is qff).
- var_name
str Name of the variable column to interpolate. This should match the one in the parameter data.
- lon_dim, lat_dim
str, optional Column names for longitude/latitude (degrees). Defaults are
"lon"and"lat".- grid_res_deg
float, optional Output grid spacing in degrees. Must be > 0. Default is 0.25.
- auto_domain
bool, optional If True (default), the interpolation domain is set to the data bounding box expanded by
buffer_degon each side.- buffer_deg
float, optional Padding (degrees) added around the data bounding box when
auto_domain=True.- pointtuple(float, float), optional
Lower-left corner (lon0, lat0) of the output grid (degrees). Required when
auto_domain=False.- grid_x, grid_y
float, optional Domain size in degrees in x (lon) / y (lat). Required when
auto_domain=False.- sigma_deg
float, optional Gaussian width in degrees. If None, an automatic estimate based on station density is used.
- influence_radius_deg
float, optional Radius of influence in degrees. If None, uses
influence_k_sigma * sigma_deg.- influence_k_sigma
float, optional Multiplier used when
influence_radius_degis None. Default is 4.0.- num_iter
int, optional Number of self-convolutions used by convolution-based methods. Must be >= 1. The number of performed self-convolutions of the underlying rect-kernel. Applies only if method is ‘optimized_convolution’ or ‘convolution’. The default is 2. Applies only to Convol interpolations: one of 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 20, 50.
- method
str, optional Passed to
fastbarnes_new.interpolation.barnes(). Common values include"optimized_convolution","convolution","radius","naive". The possible implementations that can be chosen are ‘naive’ for the straightforward implementation (algorithm A from paper), ‘radius’ to consider only sample points within a specific radius of influence, both with an algorithmic complexity of \(O(N \\times W \\times H)\). The choice ‘convolution’ implements algorithm B specified in the paper and ‘optimized_convolution’ is its optimization by appending tail values to the rectangular kernel. The latter two algorithms reduce the complexity down to \(O(N + W \\times H)\).- max_dist_sigma
float, optional Maximum distance (in units of
sigma) for which interpolation is computed. If None, usesinfluence_radius_deg / sigma_deg.- min_weight
float, optional Minimum Gaussian weight threshold used by radius-based methods in the backend.
- mask_radius_deg
float, optional If provided (or if None defaults to
influence_radius_deg), grid points farther than this radius (degrees) from all stations are set to NaN.- missing_sentinelstuple(float, …), optional
Values treated as missing in
var_nameand replaced by NaN before cleaning.
Returns¶
xarray.DataArrayInterpolated field on a regular lat/lon grid. Dimensions are
(lat_dim, lon_dim). The returned DataArray includes useful metadata inattrs(grid spacing, sigma, influence radius, domain definition, etc.).
See also
Zürcher, B. K.: Fast approximate Barnes interpolation: illustrated by Python-Numba implementation fast-barnes-py v1.0, Geosci. Model Dev., 16, 1697–1711, https://doi.org/10.5194/gmd-16-1697-2023, 2023.
- easyclimate.interp.interp_spatial_barnes_rs(data: pandas.DataFrame, var_name: str, *, lon_dim: str = 'lon', lat_dim: str = 'lat', grid_res_deg: float = 0.25, auto_domain: bool = True, buffer_deg: float = 1.0, point: tuple[float, float] | None = None, grid_x: float | None = None, grid_y: float | None = None, sigma_deg: float | None = None, influence_radius_deg: float | None = None, influence_k_sigma: float = 4.0, num_iter: int = 2, method: str = 'optimized_convolution', max_dist_sigma: float | None = None, min_weight: float = 1e-06, mask_radius_deg: float | None = None, missing_sentinels: tuple[float, Ellipsis] = (9999.9, 9999.0, -9999.0, -9999.9))¶
Rust-accelerated Barnes interpolation wrapper (lon/lat in degrees).
This function is API-compatible with
interp_spatial_barnes()but calls the Rust backend (easyclimate_rust._easyclimate_rust.barnes()) for the core interpolation step, and (optionally) uses the Rust radius mask (easyclimate_rust._easyclimate_rust.radius_mask_2d()) to avoid showing extrapolation far from stations.Mathematically, the method targets the same Barnes analysis:
\[f(\boldsymbol{x})=\frac{\sum_{k=1}^N f_k\cdot w_k(\boldsymbol{x})}{\sum_{k=1}^N w_k(\boldsymbol{x})}\]with Gaussian weights
\[w_k(\boldsymbol{x})=\text{e}^{-\frac{1}{2\sigma^2}\left|x-\boldsymbol{x}_k\right|^2}\]Note
All user-facing spatial parameters are in degrees (same conversion rules as
interp_spatial_barnes()).sigma_degauto-estimation uses a station-density spacing heuristic.If
mask_radius_degis enabled, masking is performed in Rust on the regular grid.
Parameters¶
Identical to
interp_spatial_barnes().Returns¶
xarray.DataArrayInterpolated field on a regular lat/lon grid. Dimensions are
(lat_dim, lon_dim). The returned DataArray includes useful metadata inattrs(grid spacing, sigma, influence radius, domain definition, etc.).
See also
interp_spatial_barnes()(pure Python/NumPy backend)Zürcher, B. K.: Fast approximate Barnes interpolation: illustrated by Python-Numba implementation fast-barnes-py v1.0, Geosci. Model Dev., 16, 1697–1711, https://doi.org/10.5194/gmd-16-1697-2023, 2023.
- easyclimate.interp.interp_spatial_barnesS2(data: pandas.DataFrame, var_name: str, *, lon_dim: str = 'lon', lat_dim: str = 'lat', grid_res_deg: float = 0.25, auto_domain: bool = True, buffer_deg: float = 1.0, point: tuple[float, float] | None = None, grid_x: float | None = None, grid_y: float | None = None, sigma_deg: float | None = None, influence_radius_deg: float | None = None, influence_k_sigma: float = 4.0, num_iter: int = 4, method: Literal['optimized_convolution_S2', 'naive_S2'] = 'optimized_convolution_S2', max_dist_sigma: float | None = None, resample: bool = True, mask_radius_deg: float | None = None, missing_sentinels: tuple[float, Ellipsis] = (9999.9, 9999.0, -9999.0, -9999.9), lambert_proj=None, lambert_grid=None, auto_proj: bool = True) xarray.DataArray¶
Computes Barnes interpolation on the unit sphere \(S^2\) (spherical metric) for observations given in lon/lat degrees, returning a regular lon/lat grid as an
xarray.DataArray.Compared to planar Barnes interpolation, the \(S^2\) variant uses spherical geometry internally (implemented by
fastbarnes_new.interpolationS2.barnes_S2()). For performance, the optimized method relies on a Lambert projection and a projected grid; this wrapper forwards Lambert options and performs early validation so errors are raised in Python (more readable than deep backend errors).Parameters¶
- data
pandas.DataFrame Input table containing at least
lon_dim,lat_dimandvar_namecolumns. Any rows with NaN/Inf in these columns are dropped after cleaning.There should be a similar structure as follows
lon
lat
qff
-3.73
56.33
995.1
2.64
47.05
1012.5
…
…
…
Note
Data points should contain longitude (lon), latitude (lat) and data variables (the above data variable name is qff).
- var_name
str Name of the variable column to interpolate. This should match the one in the parameter data.
- lon_dim, lat_dim
str, optional Column names for longitude/latitude (degrees). Defaults are
"lon"and"lat".- grid_res_deg
float, optional Output grid spacing in degrees. Must be > 0. Default is 0.25.
- auto_domain
bool, optional If True (default), the interpolation domain is set to the data bounding box expanded by
buffer_degon each side.- buffer_deg
float, optional Padding (degrees) added around the data bounding box when
auto_domain=True.- pointtuple(float, float), optional
Lower-left corner (lon0, lat0) of the output grid (degrees). Required when
auto_domain=False.- grid_x, grid_y
float, optional Domain size in degrees in x (lon) / y (lat). Required when
auto_domain=False.- sigma_deg
float, optional Gaussian width in degrees. If None, an automatic estimate based on station density is used.
- influence_radius_deg
float, optional Radius of influence in degrees. If None, uses
influence_k_sigma * sigma_deg.- influence_k_sigma
float, optional Multiplier used when
influence_radius_degis None. Default is 4.0.- num_iter
int, optional Number of self-convolutions used by convolution-based methods. Must be >= 1. The number of performed self-convolutions of the underlying rect-kernel. Applies only if method is ‘optimized_convolution’ or ‘convolution’. The default is 2. Applies only to Convol interpolations: one of 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 20, 50.
- method{‘optimized_convolution_S2’, ‘naive_S2’}, default: ‘optimized_convolution_S2’.
Designates the Barnes interpolation method to be used. The possible implementations that can be chosen are ‘naive_S2’ for the straightforward implementation (algorithm A from the paper) with an algorithmic complexity of \(O(N \times W \times H)\). The choice ‘optimized_convolution_S2’ implements the optimized algorithm B specified in the paper by appending tail values to the rectangular kernel. The latter algorithm has a reduced complexity of \(O(N + W \times H)\). The default is ‘optimized_convolution_S2’.
- max_dist_sigma
float, optional Maximum distance (in units of
sigma) for which interpolation is computed. If None, usesinfluence_radius_deg / sigma_deg.- min_weight
float, optional Minimum Gaussian weight threshold used by radius-based methods in the backend.
- mask_radius_deg
float, optional If provided (or if None defaults to
influence_radius_deg), grid points farther than this radius (degrees) from all stations are set to NaN.- missing_sentinelstuple(float, …), optional
Values treated as missing in
var_nameand replaced by NaN before cleaning.- resample
bool, default True Passed through to the S2 backend (controls internal resampling behavior).
- lambert_proj, lambert_grid, auto_proj
Lambert projection configuration forwarded to
barnes_S2(). Ifauto_proj=Trueand nolambert_projis provided, the backend will infer a projection; domains crossing the equator are not supported by the optimized S2 convolution path (split into hemispheres or pass an explicit projection/grid).
Returns¶
xarray.DataArrayInterpolated spherical Barnes field on a regular lat/lon grid, with metadata in
attrs.
See also
interp_spatial_barnes()(planar metric)interp_spatial_barnesS2_rs()(Rust backend on \(S^2\))Zürcher, B. K.: Fast approximate Barnes interpolation: illustrated by Python-Numba implementation fast-barnes-py v1.0, Geosci. Model Dev., 16, 1697–1711, https://doi.org/10.5194/gmd-16-1697-2023, 2023.
- data
- easyclimate.interp.interp_spatial_barnesS2_rs(data: pandas.DataFrame, var_name: str, *, lon_dim: str = 'lon', lat_dim: str = 'lat', grid_res_deg: float = 0.25, auto_domain: bool = True, buffer_deg: float = 1.0, point: tuple[float, float] | None = None, grid_x: float | None = None, grid_y: float | None = None, sigma_deg: float | None = None, influence_radius_deg: float | None = None, influence_k_sigma: float = 4.0, num_iter: int = 4, method: Literal['optimized_convolution_S2', 'naive_S2'] = 'optimized_convolution_S2', max_dist_sigma: float | None = None, resample: bool = True, mask_radius_deg: float | None = None, missing_sentinels: tuple[float, Ellipsis] = (9999.9, 9999.0, -9999.0, -9999.9), lambert_proj=None, lambert_grid=None, auto_proj: bool = True) xarray.DataArray¶
Rust-accelerated Barnes interpolation on the unit sphere \(S^2\) (spherical metric).
This function mirrors
interp_spatial_barnesS2()but uses the Rust backend (easyclimate_rust._easyclimate_rust.barnes_s2()) for the S2 interpolation. It keeps the same degree-based, user-friendly API and returns the samexarray.DataArraylayout.Parameters¶
Identical to
interp_spatial_barnesS2().Returns¶
xarray.DataArrayInterpolated spherical Barnes field on a regular lat/lon grid, with metadata in
attrs.
See also
interp_spatial_barnes()(planar metric)interp_spatial_barnesS2_rs()(Rust backend on \(S^2\))Zürcher, B. K.: Fast approximate Barnes interpolation: illustrated by Python-Numba implementation fast-barnes-py v1.0, Geosci. Model Dev., 16, 1697–1711, https://doi.org/10.5194/gmd-16-1697-2023, 2023.
- easyclimate.interp.interp1d_vertical_model2pressure(pressure_data: xarray.DataArray, variable_data: xarray.DataArray, vertical_input_dim: str, vertical_output_dim: str, vertical_output_level: list[int | float], kind: str = 'linear', bounds_error=None, fill_value=np.nan, assume_sorted: bool = False) xarray.DataArray¶
Interpolating variables from the model levels (e.g., sigma vertical coordinate) to the standard pressure levels by 1-D function.
- pressure_data:
xarray.DataArray. The pressure on each model level.
- variable_data:
xarray.DataArray. variable to be interpolated on model levels.
- vertical_input_dim:
str. The name of the model levels dimension.
- vertical_output_dim:
str. The name of the standard pressure levels dimension, often assigned the value ‘plev’ (Pa) or ‘lev’ (hPa).
- vertical_output_level:
list[int | float]. Customized standard pressure levels to be interpolated, e.g., [5000, 10000, 15000, 20000, 25000, 30000, 40000, 50000, 60000, 70000, 85000, 92500, 100000].
Note
The unit of vertical_output_level shuold be same as pressure_data, e.g., the unit of pressure_data is Pa, and the unit of vertical_output_level shuold be Pa.
- kind:
strorint, default: ‘linear’. Specifies the kind of interpolation as a string or as an integer specifying the order of the spline interpolator to use. The string has to be one of ‘linear’, ‘nearest’, ‘nearest-up’, ‘zero’, ‘slinear’, ‘quadratic’, ‘cubic’, ‘previous’, or ‘next’. ‘zero’, ‘slinear’, ‘quadratic’ and ‘cubic’ refer to a spline interpolation of zeroth, first, second or third order; ‘previous’ and ‘next’ simply return the previous or next value of the point; ‘nearest-up’ and ‘nearest’ differ when interpolating half-integers (e.g. 0.5, 1.5) in that ‘nearest-up’ rounds up and ‘nearest’ rounds down.
- bounds_error:
bool, default: None. If True, a ValueError is raised any time interpolation is attempted on a value outside of the range of pressure_data (where extrapolation is necessary). If False, out of bounds values are assigned fill_value. By default, an error is raised unless fill_value=”extrapolate”.
- fill_value: array-like or (array-like, array_like) or “extrapolate”, default: np.nan.
If a ndarray (or float), this value will be used to fill in for requested points outside of the data range. If not provided, then the default is NaN. The array-like must broadcast properly to the dimensions of the non-interpolation axes.
If a two-element tuple, then the first element is used as a fill value for x_new < x[0] and the second element is used for x_new > x[-1]. Anything that is not a 2-element tuple (e.g., list or ndarray, regardless of shape) is taken to be a single array-like argument meant to be used for both bounds as below, above = fill_value, fill_value. Using a two-element tuple or ndarray requires bounds_error=False.
If “extrapolate”, then points outside the data range will be extrapolated.
- assume_sortedbool:
bool, default: False. If False, values of x can be in any order and they are sorted first. If True, pressure_data has to be an array of monotonically increasing values.
See also
- pressure_data:
- easyclimate.interp.interp1d_vertical_pressure2altitude(z_data: xarray.DataArray, variable_data: xarray.DataArray, target_heights: numpy.array, vertical_input_dim: str, vertical_output_dim: str, kind: str = 'cubic', bounds_error=None, fill_value='extrapolate', assume_sorted: bool = False) xarray.DataArray¶
Interpolating variables from the pressure levels (e.g., sigma vertical coordinate) to the altitude levels by 1-D function.
- z_data:
xarray.DataArray. The vertical atmospheric geopotential height on each pressure level.
- variable_data:
xarray.DataArray. variable to be interpolated on model levels.
Note
The shape of z_data z_data and variable_data shuold be the same.
- target_heights:
numpy.array. Altitude interpolation sampling range. e.g., np.linspace(0, 15000, 100).
Note
The unit of target_heights shuold be same as z_data, e.g., the unit of target_heights is meter, and the unit of z_data shuold be meter.
- vertical_input_dim:
str. The name of the standard pressure levels dimension, often assigned the value ‘plev’ (Pa) or ‘lev’ (hPa).
Note
The vertical_input_dim z_data and variable_data shuold be the same.
- vertical_output_dim:
str. The name of the altitude dimension, often assigned the value ‘altitude’.
- kind:
strorint, default: ‘cubic’. Specifies the kind of interpolation as a string or as an integer specifying the order of the spline interpolator to use. The string has to be one of ‘linear’, ‘nearest’, ‘nearest-up’, ‘zero’, ‘slinear’, ‘quadratic’, ‘cubic’, ‘previous’, or ‘next’. ‘zero’, ‘slinear’, ‘quadratic’ and ‘cubic’ refer to a spline interpolation of zeroth, first, second or third order; ‘previous’ and ‘next’ simply return the previous or next value of the point; ‘nearest-up’ and ‘nearest’ differ when interpolating half-integers (e.g. 0.5, 1.5) in that ‘nearest-up’ rounds up and ‘nearest’ rounds down.
- bounds_error:
bool, default: None. If True, a ValueError is raised any time interpolation is attempted on a value outside of the range of pressure_data (where extrapolation is necessary). If False, out of bounds values are assigned fill_value. By default, an error is raised unless fill_value=”extrapolate”.
- fill_value: array-like or (array-like, array_like) or “extrapolate”, default: extrapolate.
If a ndarray (or float), this value will be used to fill in for requested points outside of the data range. If not provided, then the default is NaN. The array-like must broadcast properly to the dimensions of the non-interpolation axes.
If a two-element tuple, then the first element is used as a fill value for x_new < x[0] and the second element is used for x_new > x[-1]. Anything that is not a 2-element tuple (e.g., list or ndarray, regardless of shape) is taken to be a single array-like argument meant to be used for both bounds as below, above = fill_value, fill_value. Using a two-element tuple or ndarray requires bounds_error=False.
If “extrapolate”, then points outside the data range will be extrapolated.
- assume_sortedbool:
bool, default: False. If False, values of x can be in any order and they are sorted first. If True, pressure_data has to be an array of monotonically increasing values.
See also
- z_data:
- easyclimate.interp.interp1d_vertical_pressure2altitude_linear_rs(z_data: xarray.DataArray, variable_data: xarray.DataArray, target_heights: numpy.array, vertical_input_dim: str, vertical_output_dim: str) xarray.DataArray¶
Interpolating variables from the pressure levels (e.g., sigma vertical coordinate) to the altitude levels using easyclimate-rust.
- z_data:
xarray.DataArray. The vertical atmospheric geopotential height on each pressure level.
- variable_data:
xarray.DataArray. variable to be interpolated on model levels.
Note
The shape of z_data z_data and variable_data shuold be the same.
- target_heights:
numpy.array. Altitude interpolation sampling range. e.g., np.linspace(0, 15000, 100).
Note
The unit of target_heights shuold be same as z_data, e.g., the unit of target_heights is meter, and the unit of z_data shuold be meter.
- vertical_input_dim:
str. The name of the standard pressure levels dimension, often assigned the value ‘plev’ (Pa) or ‘lev’ (hPa).
Note
The vertical_input_dim z_data and variable_data shuold be the same.
- vertical_output_dim:
str. The name of the altitude dimension, often assigned the value ‘altitude’.
See also
- z_data:
- easyclimate.interp.interp_mesh2point(data_input: xarray.DataArray, df: pandas.DataFrame, lon_dim_mesh: str = 'lon', lat_dim_mesh: str = 'lat', lon_dim_df: str = 'lon', lat_dim_df: str = 'lat', method: Literal['linear', 'nearest', 'slinear', 'cubic', 'quintic', 'pchip'] = 'linear')¶
Interpolate values from a regular grid/mesh to specific point locations.
Parameters¶
- data_input
xarray.DataArray. Input grid data with latitude and longitude dimensions.
- df
pandas.DataFrame. DataFrame containing point locations with latitude and longitude columns.
- lon_dim_mesh
str, optional Name of the longitude dimension in the input grid, by default
"lon".- lat_dim_mesh
str, optional Name of the latitude dimension in the input grid, by default
"lat".- lon_dim_df
str, optional Name of the longitude column in the DataFrame, by default
"lon".- lat_dim_df
str, optional Name of the latitude column in the DataFrame, by default
"lat".- method
str, optional Interpolation method to use. Options are:
"linear": bilinear interpolation (default)"nearest": nearest neighbor"slinear": spline linear"cubic": spline cubic"quintic": spline quintic"pchip": piecewise cubic Hermite interpolating polynomial
Returns¶
pandas.DataFrame.Original DataFrame with an additional column
"interpolated_value"containing the interpolated values. Points outside the grid range will have NaN values.
Raises¶
- ValueError
If no points in the DataFrame fall within the grid’s spatial extent.
Examples¶
>>> import xarray as xr >>> import pandas as pd >>> # Create sample grid data >>> lats = np.linspace(-90, 90, 181) >>> lons = np.linspace(-180, 180, 361) >>> data = xr.DataArray(np.random.rand(181, 361), dims=['lat', 'lon'], ... coords={'lat': lats, 'lon': lons}) >>> # Create sample points >>> points = pd.DataFrame({'lat': [45.5, 30.2], 'lon': [-120.3, 150.7]}) >>> # Interpolate >>> result = interp_mesh2point(data, points)
- data_input
- easyclimate.interp.interp_mesh2point_withtime(data_input: xarray.DataArray, stations_df: pandas.DataFrame, lon_dim_df: str, lat_dim_df: str, station_dim_df: str, lon_dim_mesh: str = 'lon', lat_dim_mesh: str = 'lat', time_dim_mesh: str = 'time', method: Literal['linear', 'nearest', 'slinear', 'cubic', 'quintic', 'pchip'] = 'linear') xarray.DataArray¶
Interpolate gridded data to specific point locations with time series preservation.
Parameters¶
- data_input
xarray.DataArray. Input grid data with time, latitude and longitude dimensions.
- stations_df
pandas.DataFrame. DataFrame containing station locations with ID, latitude and longitude.
- lon_dim_df
str. Name of the longitude column in the DataFrame.
- lat_dim_df
str. Name of the latitude column in the DataFrame.
- station_dim_df
str. Name of the station ID column in the DataFrame.
- lon_dim_mesh
str, default: “lon” Name of the longitude dimension in the input grid.
- lat_dim_mesh
str, default: “lat” Name of the latitude dimension in the input grid.
- time_dim_mesh
str, default: “time” Name of the time dimension in the input grid.
- method
str, optional Interpolation method:
“linear”: bilinear interpolation (default)
“nearest”: nearest neighbor
“slinear”: spline linear
“cubic”: spline cubic
“quintic”: spline quintic
“pchip”: piecewise cubic Hermite interpolating polynomial
Returns¶
xarray.Dataset.xarray.Datasetwith dimensions(station, time)containing:Interpolated values
Station coordinates
Examples¶
>>> import xarray as xr >>> import pandas as pd >>> import numpy as np >>> # Create sample grid >>> times = pd.date_range("2020-01-01", periods=5) >>> lats = np.linspace(-45, -10, 100) >>> lons = np.linspace(110, 156, 120) >>> data = xr.DataArray( ... np.random.rand(5, 100, 120), ... dims=["time", "lat", "lon"], ... coords={"time": times, "lat": lats, "lon": lons} ... ) >>> # Create stations >>> stations = pd.DataFrame({ ... "station_id_col": [1001, 1002], ... "lat_col": [-15.5, -20.3], ... "lon_col": [125.5, 130.2] ... }) >>> # Interpolate >>> result = interp_mesh2point_withtime( ... data, stations, ... lon_dim_df="lon_col", lat_dim_df="lat_col", station_dim_df="station_id_col" ... )
- data_input
- easyclimate.interp.interp_vinth2p_dp(data_input: xarray.DataArray, surface_pressure_data: xarray.DataArray, surface_pressure_data_units: Literal['hPa', 'Pa', 'mbar'], hybrid_A_coefficients: xarray.DataArray, hybrid_B_coefficients: xarray.DataArray, vertical_output_level: list[int | float], vertical_input_dim: str, vertical_output_dim: str, vertical_output_dim_units: str, interp_method: Literal['linear', 'log', 'loglog'] = 'linear', extrapolation: bool = False, lon_dim: str = 'lon', lat_dim: str = 'lat', time_dim: str = 'time', p0_hPa=1000.0) xarray.DataArray¶
Interpolate atmospheric data from Community Atmosphere Model (CAM) hybrid sigma-pressure coordinates to pressure levels.
This function performs vertical interpolation of atmospheric data (typically temperature) from hybrid sigma-pressure coordinates to specified pressure levels using the NCL’s
vinth2palgorithm implemented in Fortran.Parameters¶
- data_input
xarray.DataArray. Input 3D temperature field on hybrid levels with dimensions, e.g., (time, lev, lat, lon).
- surface_pressure_data
xarray.DataArray. Surface pressure field with dimensions, e.g., (time, lat, lon).
- surface_pressure_data_unitsLiteral[“hPa”, “Pa”, “mbar”]
Units of the surface pressure data.
- hybrid_A_coefficients
xarray.DataArray. Hybrid A coefficients (pressure term) for the model levels.
- hybrid_B_coefficients
xarray.DataArray. Hybrid B coefficients (sigma term) for the model levels.
- vertical_output_levellist[
int] List of target pressure levels for interpolation.
- vertical_input_dim
str. Name of the vertical dimension in the input data.
- vertical_output_dim
str. Name to use for the vertical dimension in the output data.
- vertical_output_dim_units
str. Units for the output pressure levels (must be convertible to hPa).
- interp_methodLiteral[“linear”, “log”, “loglog”], optional
Interpolation method:
"linear": Linear interpolation"log": Logarithmic interpolation"loglog": Log-log interpolation
Default is
"linear".- extrapolation
bool., optional Whether to extrapolate below the lowest model level when needed. Default is
False.- lon_dim
str., optional Name of the longitude dimension. Default is
"lon".- lat_dim
str., optional Name of the latitude dimension. Default is
"lat".- time_dim:
str, default: time. The time coordinate dimension name. Default is
"time".- p0_hPa
float., optional Reference pressure in hPa for hybrid level calculation. Default is
1000.0hPa.
Returns¶
xarray.DataArrayInterpolated data on pressure levels with dimensions, e.g., (time, plev, lat, lon), where plev corresponds to vertical_output_level.
Examples¶
>>> # Interpolate temperature to pressure levels >>> interp_vinth2p_dp( ... data_input=temp_data, ... surface_pressure_data=psfc_data, ... surface_pressure_data_units="Pa", ... hybrid_A_coefficients=hyam, ... hybrid_B_coefficients=hybm, ... vertical_output_level=[1000, 850, 700, 500, 300], ... vertical_input_dim="lev", ... vertical_output_dim="plev", ... vertical_output_dim_units="hPa", ... interp_method="log" ... )
- data_input
- easyclimate.interp.interp_vinth2p_ecmwf(data_input: xarray.DataArray, surface_pressure_data: xarray.DataArray, surface_pressure_data_units: Literal['hPa', 'Pa', 'mbar'], hybrid_A_coefficients: xarray.DataArray, hybrid_B_coefficients: xarray.DataArray, vertical_output_level: list[int | float], vertical_input_dim: str, vertical_output_dim: str, vertical_output_dim_units: str, variable_flag: Literal['T', 'Z', 'other'], temperature_bottom_data: xarray.DataArray | None = None, surface_geopotential_data: xarray.DataArray | None = None, interp_method: Literal['linear', 'log', 'loglog'] = 'linear', extrapolation: bool = True, lon_dim: str = 'lon', lat_dim: str = 'lat', time_dim: str = 'time', p0_hPa: float = 1000.0) xarray.DataArray¶
Interpolate atmospheric data from Community Atmosphere Model (CAM) hybrid sigma-pressure coordinates to pressure levels using ECMWF extrapolation methods.
This function performs vertical interpolation of atmospheric data from hybrid sigma-pressure coordinates to specified pressure levels using the NCL’s
vinth2p_ecmwfalgorithm implemented in Fortran. It supports ECMWF-specific extrapolation for temperature (‘T’) and geopotential height (‘Z’) below the lowest hybrid level.Parameters¶
- data_input
xarray.DataArray Input 3D field (e.g., temperature or geopotential) on hybrid levels with dimensions, e.g., (time, lev, lat, lon).
- surface_pressure_data
xarray.DataArray Surface pressure field with dimensions, e.g., (time, lat, lon).
- surface_pressure_data_units
Literal["hPa", "Pa", "mbar"] Units of the surface pressure data.
- hybrid_A_coefficients
xarray.DataArray Hybrid A coefficients (pressure term) for the model levels.
- hybrid_B_coefficients
xarray.DataArray Hybrid B coefficients (sigma term) for the model levels.
- vertical_output_levellist[
int|float] List of target pressure levels for interpolation (in specified units).
- vertical_input_dim
str Name of the vertical dimension in the input data.
- vertical_output_dim
str Name to use for the vertical dimension in the output data.
- vertical_output_dim_units
str Units for the output pressure levels (must be convertible to hPa).
- variable_flag
Literal["T", "Z", "other"] Indicates the type of variable being interpolated: - “T”: Temperature (uses ECMWF extrapolation if enabled) - “Z”: Geopotential height (uses ECMWF extrapolation if enabled) - “other”: Any other variable (uses lowest level value for extrapolation)
- temperature_bottom_dataOptional[
xarray.DataArray], optional Temperature at the lowest model level (required for ‘Z’ extrapolation). Dimensions, e.g., (time, lat, lon). Default is None.
- surface_geopotential_dataOptional[
xarray.DataArray], optional Surface geopotential (required for ‘T’ or ‘Z’ extrapolation). Dimensions, e.g., (time, lat, lon). Default is None.
- interp_method
Literal["linear", "log", "loglog"], optional Interpolation method:
"linear": Linear interpolation"log": Logarithmic interpolation"loglog": Log-log interpolation
Default is
"linear".- extrapolation
bool, optional Whether to extrapolate below the lowest model level when needed. Default is False.
- lon_dim
str, optional Name of the longitude dimension. Default is “lon”.
- lat_dim
str, optional Name of the latitude dimension. Default is “lat”.
- time_dim:
str, default: time. The time coordinate dimension name. Default is “time”.
- p0_hPa
float, optional Reference pressure in hPa for hybrid level calculation. Default is
1000.0 hPa.
Returns¶
xarray.DataArrayInterpolated data on pressure levels with dimensions, e.g., (time, plev, lat, lon), where plev corresponds to
vertical_output_level.
Notes¶
The hybrid level pressure is calculated as: \(P = A \cdot p_0 + B \cdot \mathrm{psfc}\)
Output pressure levels are converted to hPa internally for calculations.
Missing values are converted to NaN in the output.
ECMWF extrapolation is applied only when
extrapolation=Trueand variable_flag is ‘T’ or ‘Z’.For ‘T’ or ‘Z’, tbot and phis must be provided when
extrapolation=True.
Examples¶
>>> # Interpolate temperature to pressure levels with ECMWF extrapolation >>> interp_vinth2p_ecmwf( ... data_input=temp_data, ... surface_pressure_data=psfc_data, ... surface_pressure_data_units="Pa", ... hybrid_A_coefficients=hyam, ... hybrid_B_coefficients=hybm, ... vertical_output_level=[1000, 850, 700, 500, 300], ... vertical_input_dim="lev", ... vertical_output_dim="plev", ... vertical_output_dim_units="hPa", ... variable_flag="T", ... temperature_bottom_data=tbot_data, ... surface_geopotential_data=phis_data, ... interp_method="log", ... extrapolation=True ... )
- data_input
- easyclimate.interp.interp_vintp2p_ecmwf(data_input: xarray.DataArray, pressure_data: xarray.DataArray, pressure_data_units: Literal['hPa', 'Pa', 'mbar'], surface_pressure_data: xarray.DataArray, surface_pressure_data_units: Literal['hPa', 'Pa', 'mbar'], vertical_output_level: list[int | float], vertical_input_dim: str, vertical_output_dim: str, vertical_output_dim_units: str, variable_flag: Literal['T', 'Z', 'other'], temperature_bottom_data: xarray.DataArray | None = None, surface_geopotential_data: xarray.DataArray | None = None, interp_method: Literal['linear', 'log', 'loglog'] = 'linear', extrapolation: bool = False, lon_dim: str = 'lon', lat_dim: str = 'lat', time_dim: str = 'time') xarray.DataArray¶
Interpolates data at multidimensional pressure levels to constant pressure coordinates and uses an ECMWF formulation to extrapolate values below ground.
This function performs vertical interpolation of atmospheric data from input pressure levels to specified output pressure levels using the NCL’s vintp2p_ecmwf algorithm implemented in Fortran. It supports ECMWF-specific extrapolation for temperature (‘T’) and geopotential height (‘Z’) below the lowest pressure level.
Parameters¶
- data_input
xarray.DataArray Input 3D field (e.g., temperature or geopotential) on pressure levels with dimensions, e.g., (time, lev, lat, lon).
- pressure_data
xarray.DataArray 3D pressure field corresponding to the input data levels, with dimensions, e.g., (time, lev, lat, lon).
- pressure_data_unitsLiteral[“hPa”, “Pa”, “mbar”]
Units of the pressure data.
- surface_pressure_data
xarray.DataArray Surface pressure field with dimensions, e.g., (time, lat, lon).
- surface_pressure_data_unitsLiteral[“hPa”, “Pa”, “mbar”]
Units of the surface pressure data.
- vertical_output_levellist[
int|float] List of target pressure levels for interpolation (in specified units).
- vertical_input_dim
str Name of the vertical dimension in the input data.
- vertical_output_dim
str Name to use for the vertical dimension in the output data.
- vertical_output_dim_units
str Units for the output pressure levels (must be convertible to hPa).
- variable_flag
Literal["T", "Z", "other"] Indicates the type of variable being interpolated: - “T”: Temperature (uses ECMWF extrapolation if enabled) - “Z”: Geopotential height (uses ECMWF extrapolation if enabled) - “other”: Any other variable (uses lowest level value for extrapolation)
- temperature_bottom_dataOptional[
xarray.DataArray], optional Temperature at the lowest model level (required for ‘Z’ extrapolation). Dimensions, e.g., (time, lat, lon). Default is None.
- surface_geopotential_dataOptional[
xarray.DataArray], optional Surface geopotential (required for ‘T’ or ‘Z’ extrapolation). Dimensions, e.g., (time, lat, lon). Default is None.
- interp_method
Literal["linear", "log", "loglog"], optional Interpolation method:
“linear”: Linear interpolation
“log”: Logarithmic interpolation
“loglog”: Log-log interpolation
Default is “linear”.
- extrapolation
bool, optional Whether to extrapolate below the lowest pressure level when needed. Default is False.
- lon_dim
str, optional Name of the longitude dimension. Default is “lon”.
- lat_dim
str, optional Name of the latitude dimension. Default is “lat”.
- time_dim:
str, default: time. The time coordinate dimension name. Default is “time”.
Returns¶
xarray.DataArrayInterpolated data on pressure levels with dimensions, e.g., (time, plev, lat, lon), where plev corresponds to vertical_output_level.
Notes¶
The input pressure levels are provided directly via
pressure_data.Output pressure levels and surface pressure are converted to hPa internally for calculations.
Missing values are converted to NaN in the output.
ECMWF extrapolation is applied only when
extrapolation=Trueand variable_flag is ‘T’ or ‘Z’.For ‘T’ or ‘Z’,
temperature_bottom_dataandsurface_geopotential_datamust be provided whenextrapolation=True.
Examples¶
>>> # Interpolate temperature to pressure levels with ECMWF extrapolation >>> interp_vintp2p_ecmwf( ... data=temp_data, ... pressure_data=pres_data, ... pressure_data_units="Pa", ... surface_pressure_data=psfc_data, ... surface_pressure_data_units="Pa", ... vertical_output_level=[1000, 850, 700, 500, 300], ... vertical_input_dim="lev", ... vertical_output_dim="plev", ... vertical_output_dim_units="hPa", ... variable_flag="T", ... temperature_bottom_data=tbot_data, ... surface_geopotential_data=phis_data, ... interp_method="log", ... extrapolation=True ... )
- data_input