erlab.analysis.transform¶
Transformations.
Functions
|
Rotate an array in the plane defined by the two axes. |
|
Rotate a 2D DataArray in the plane defined by the two dimensions. |
|
Rotate a 3D DataArray in the plane defined by the two dimensions. |
|
Shifts the values of a DataArray along a single dimension. |
|
Symmetrize a DataArray along a specified dimension around a given center. |
|
Symmetrize a plane by averaging equally spaced rotations. |
- erlab.analysis.transform.rotate(darr, angle, axes=(0, 1), center=(0.0, 0.0), *, reshape=True, order=1, mode='constant', cval=np.nan, prefilter=True)[source]¶
Rotate an array in the plane defined by the two axes.
- Parameters:
darr (
DataArray) – The array to rotate.angle (
float) – The rotation angle in degrees.axes (
tupleof2 intsor2 strings, optional) – The two axes that define the plane of rotation. Default is the first two axes. If strings are provided, they must be valid dimension names in the input array.center (
tupleof2 floatsordict, optional) – The center of rotation in data coordinates. If a tuple, it is given as values along the dimensions specified inaxes. If a dict, it must have keys that correspond toaxes. Default is (0, 0).reshape (
bool, default:True) – IfTrue, the output shape is adapted so that the input array is contained completely in the output. Default isTrue.order (
int, default:1) – The order of the spline interpolation, default is 1. The order has to be in the range 0-5.mode (
str, default:"constant") – Passed toscipy.ndimage.affine_transform(). See the scipy documentation for more information.cval (
float, default:np.nan) – Passed toscipy.ndimage.affine_transform(). See the scipy documentation for more information.prefilter (
bool, default:True) – Passed toscipy.ndimage.affine_transform(). See the scipy documentation for more information.
- Returns:
darr (
xarray.DataArray) – The rotated array.- Return type:
See also
scipy.ndimage.affine_transformThe function that performs the affine transformation on the input array.
scipy.ndimage.rotateSimilar function that rotates a numpy array.
- erlab.analysis.transform.rotateinplane(data, rotate, **interp_kwargs)[source]¶
Rotate a 2D DataArray in the plane defined by the two dimensions.
Deprecated since version 2.9.0: Use
erlab.analysis.transform.rotate()instead.
- erlab.analysis.transform.rotatestackinplane(data, rotate, **interp_kwargs)[source]¶
Rotate a 3D DataArray in the plane defined by the two dimensions.
Deprecated since version 2.9.0: Use
erlab.analysis.transform.rotate()instead.
- erlab.analysis.transform.shift(darr, shift, along, *, shift_coords=False, keep_dim_order=True, assume_sorted=False, **shift_kwargs)[source]¶
Shifts the values of a DataArray along a single dimension.
The shift is applied using
scipy.ndimage.shift(), which uses spline interpolation. By default, the spline is of order 1 (linear interpolation).- Parameters:
darr (
DataArray) – The array to shift.shift (
float|DataArray) – The amount of shift to be applied along the specified dimension. Ifshiftis a DataArray, different shifts can be applied to different coordinates. The dimensions ofshiftmust be a subset of the dimensions ofdarr. For more information, see the note below. Ifshiftis afloat, the same shift is applied to all values along dimensionalong. This is equivalent to providing a 0-dimensional DataArray.along (
str) – Name of the dimension along which the shift is applied.shift_coords (
bool, default:False) – IfTrue, the coordinates of the output data will be changed so that the output contains all the values of the original data. IfFalse, the coordinates and shape of the original data will be retained, and only the data will be shifted. Defaults toFalse.keep_dim_order (
bool, default:True) – IfTrue, the output array will be transposed to match the input data. Otherwise, the axis order may change due to the application ofxarray.apply_ufunc(). Default isTrue.assume_sorted (
bool, default:False) – IfFalse, the data is sorted with respect toalongusingxarray.DataArray.sortby(). ProvidingTrueskips the sort. UseTruewhen you are already sure that the data is sorted ascending with respect toalong.**shift_kwargs – Additional keyword arguments passed onto
scipy.ndimage.shift. The default values of some parameters are different from scipy.orderis set to 1,cvalis set tonp.nan, andprefilteris set toFalse.
- Returns:
xarray.DataArray– The shifted DataArray.- Return type:
Note
All dimensions in
shiftmust be a dimension indarr.The
shiftarray values are divided by the step size along thealongdimension.NaN values in
shiftare treated as zero.
Example
>>> import xarray as xr >>> import numpy as np >>> import erlab.analysis as era >>> darr = xr.DataArray( ... np.array([[1, 2, 3], [4, 5, 6], [7, 8, 9]]).astype(float), dims=["x", "y"] ... ) >>> shift_arr = xr.DataArray([1, 0, 2], dims=["x"]) >>> shifted = era.transform.shift(darr, shift_arr, along="y") >>> print(shifted) <xarray.DataArray (x: 3, y: 3)> Size: 72B array([[nan, 1., 2.], [ 4., 5., 6.], [nan, nan, 7.]]) Dimensions without coordinates: x, y
- erlab.analysis.transform.symmetrize(darr, dim, *, center=0.0, subtract=False, mode='full', part='both', interp_kw=None)[source]¶
Symmetrize a DataArray along a specified dimension around a given center.
This function takes an input DataArray and symmetrizes its values along the specified dimension by reflecting and summing the data in regions below and above a given center.
The operation assumes that the coordinate corresponding to the dimension is evenly spaced. Internally, the function interpolates the data to a shifted coordinate grid to align with the nearest grid point, performs the reflection, and concatenates the resulting halves.
- Parameters:
darr (
DataArray) – The input xarray DataArray to be symmetrized. Its coordinate along the specified dimension must be uniformly spaced.dim (
Hashable) – The dimension along which to perform the symmetrization.center (
float, optional) – The central value about which the data is symmetrized (default is 0.0).subtract (
bool, optional) – If True, the reflected part is subtracted from the original data instead of being added, resulting in an antisymmetrized output instead of a symmetrized one. Default is False (i.e., the reflected part is added).mode (
{'valid', 'full'}, optional) – How to handle the parts of the symmetrized data that does not overlap with the original data. If ‘valid’, only the part that exists in both the original and reflected data is returned. If ‘full’, the full symmetrized data is returned. In this case, all NaN values in the part that exists in the overlapping region are replaced with 0.0.part (
{'both', 'below', 'above'}, optional) – The part of the symmetrized data to return. If ‘both’, the full symmetrized data is returned. If ‘below’, only the part below the center is returned. If ‘above’, only the part above the center is returned.interp_kw (
dict, optional) – Additional keyword arguments passed toxarray.DataArray.interp().
- Returns:
DataArray– A symmetrized DataArray where each value is the sum of its original and reflected counterpart.- Return type:
Examples
>>> import xarray as xr >>> import numpy as np >>> import erlab.analysis as era >>> # Create a sample DataArray with uniform coordinates. >>> da = xr.DataArray( ... np.array([1, 2, 3, 4, 5, 6]), dims="x", coords={"x": np.linspace(-2, 2, 6)} ... ) >>> sym_da = era.transform.symmetrize(da, dim="x", center=0.0) >>> print(sym_da) <xarray.DataArray (x: 6)> Size: 48B array([2., 4., 6., 6., 4., 2.]) Coordinates: * x (x) float64 48B -2.0 -1.2 -0.4 0.4 1.2 2.0
- erlab.analysis.transform.symmetrize_nfold(darr, fold, axes=(0, 1), center=(0.0, 0.0), *, reshape=True, order=1, mode='constant', cval=np.nan, prefilter=True)[source]¶
Symmetrize a plane by averaging equally spaced rotations.
The input is rotated in the plane defined by
axesat angles \(360° i / n\), where \(i = 0, \ldots, n - 1\), and the rotated copies are averaged on a common output grid.- Parameters:
darr (
DataArray) – The array to symmetrize.fold (
int) – The order of the rotational symmetry. Must be at least 2. For example,fold=4applies 4-fold symmetrization by averaging over the original array and arrays rotated by 90°, 180°, and 270°.axes (
tupleof2 intsor2 strings, optional) – The two axes that define the plane of rotation. Default is the first two axes. If strings are provided, they must be valid dimension names in the input array.center (
tupleof2 floatsordict, optional) – The center of rotation in data coordinates. If a tuple, it is given as values along the dimensions specified inaxes. If a dict, it must have keys that correspond toaxes. Default is (0, 0).reshape (
bool, default:True) – IfTrue, the output shape is expanded to contain the full extent of all rotated copies. IfFalse, the symmetrized result is returned on the original grid. Default isTrue.order (
int, default:1) – The order of the spline interpolation, default is 1. The order has to be in the range 0-5.mode (
str, default:"constant") – Passed toscipy.ndimage.affine_transform(). See the scipy documentation for more information.cval (
float, default:np.nan) – Passed toscipy.ndimage.affine_transform(). See the scipy documentation for more information.prefilter (
bool, default:True) – Passed toscipy.ndimage.affine_transform(). See the scipy documentation for more information.
- Returns:
darr (
xarray.DataArray) – The rotationally symmetrized array on the original or expanded grid.- Return type: