cmomy.xCentralMoments

Contents

cmomy.xCentralMoments#

class cmomy.xCentralMoments(data, *, mom_ndim=1, copy=None, order=None, dtype=None, fastpath=False)[source]#

Bases: CentralMomentsABC[FloatT, DataArray]

Wrapper of xarray.DataArray based central moments data.

Parameters:
  • data (xarray.DataArray) – Central moments array.

  • mom_ndim ({1, 2}) – Value indicates if moments (mom_ndim = 1) or comoments (mom_ndim=2).

  • copy (bool, optional) – If True, copy the data. If None or False, attempt to use view. Note that False values will be converted to None for numpy versions >2.0. This will be changed to reflect the new behavior of the copy parameter to numpy.array() when the minimum numpy version >2.0.

  • order ({"C", "F", "A", "K"}, optional) – Order argument to numpy.asarray().

  • dtype (dtype) – Optional dtype for output data.

  • fastpath (bool) – For internal use.

Notes

Base data has the form

\[\begin{split}{\rm data}[..., i, j] = \begin{cases} \text{weight} & i = j = 0 \\ \langle x \rangle & i = 1, j = 0 \\ \langle (x - \langle x \rangle^i) (y - \langle y \rangle^j) \rangle & i + j > 0 \end{cases}\end{split}\]

Methods

assign_attrs(*args, **kwargs)

Assign attributes to data and return new object.

assign_coords([coords])

Assign coordinates to data and return new object.

assign_weight(weight[, copy])

Create object with updated weights

astype(dtype, *[, order, casting, subok, copy])

Underlying data cast to specified type

block(block_size, *[, dim, axis, block_dim, ...])

type block_size:

int | None

cmom()

Central moments.

copy()

Create a new object with copy of data.

drop_vars(names, *[, errors, _reorder, ...])

Interface to xarray.DataArray.drop_vars()

fill([value])

Fill data with value.

from_centralmoments(obj, *[, dims, attrs, ...])

Create and xCentralMoments object from CentralMoments.

from_raw(raw, *, mom_ndim)

Create object from raw moment data.

from_resample_vals(x, *y, mom, freq[, ...])

Create from resample observations/values.

from_vals(x, *y, mom[, weight, axis, dim, ...])

Create from observations/values.

isel([indexers, drop, missing_dims, ...])

Select subset of data by position.

mean([dim_combined, coords_combined])

Return mean/first moment(s) of data.

moments_to_comoments(*, mom[, mom_dims, ...])

Convert moments (mom_ndim=1) to comoments (mom_ndim=2).

new_like([data, copy, order, verify, dtype])

Create new object like self, with new data.

pipe(func_or_method, *args[, _reorder, ...])

Apply func_or_method to underlying data and wrap results in class:cmomy.xCentralMoments object.

push_data(data, *[, order, parallel])

Push data object to moments.

push_datas(datas, *[, axis, dim, order, ...])

Push and reduce multiple average central moments.

push_val(x, *y[, weight, order, parallel])

Push single sample to central moments.

push_vals(x, *y[, weight, axis, dim, order, ...])

Push multiple samples to central moments.

randsamp_freq(*[, axis, dim, nrep, nsamp, ...])

Interface to resample.randsamp_freq()

reduce(*[, by, axis, order, parallel, ...])

Create new object reduce along axis.

rename([new_name_or_name_dict])

Rename object.

resample_and_reduce(*, freq[, dim, axis, ...])

Bootstrap resample and reduce.

reset_index(dims_or_levels[, drop, ...])

Interface to xarray.DataArray.reset_index().

rmom()

Raw moments.

sel([indexers, method, tolerance, drop, ...])

Select subset of data.

set_index([indexes, append, _reorder, ...])

Interface to xarray.DataArray.set_index()

set_values(values)

Set values.

stack([dimensions, _reorder, _copy, _order, ...])

Stack dimensions.

std()

Standard deviation.

swap_dims([dims_dict, _reorder, _copy, ...])

Interface to xarray.DataArray.swap_dims().

to_c([copy])

Alias to to_centralmoments()

to_centralmoments([copy])

Create a CentralMoments object from xCentralMoments.

to_dataarray()

Underlying xarray.DataArray.

to_numpy()

Access to numpy array underlying class.

to_raw(*[, weight])

Raw moments accumulation array.

to_values()

Access underlying values

transpose(*dims[, transpose_coords, ...])

Transpose dimensions of data.

unstack([dim, fill_value, sparse, _reorder, ...])

Unstack dimensions.

var([dim_combined, coords_combined])

Return variance (second central moment) of data.

weight()

Weight data.

zero()

Zero out underlying data.

zeros(*, mom[, val_shape, dtype, order, ...])

Create a new base object.

zeros_like()

Create new empty object like self.

Attributes

attrs

Attributes of values.

coords

Coordinates of values.

data

Accessor to numpy array underlying data.

dims

Dimensions of values.

dtype

self.data.dtype.

indexes

Indexes of values.

mom

Number of moments.

mom_dims

Names of moment dimensions.

mom_ndim

Length of moments.

mom_shape

Shape of moments part.

name

Name of values.

ndim

self.data.ndim.

shape

self.data.shape.

sizes

Sizes of values.

val_dims

Names of value dimensions.

val_ndim

Number of value dimensions.

val_shape

Shape of values dimensions.

values

Access underlying central moments array.

Dunder Methods

__add__(b)

Add objects to new object.

__getitem__(key)

Get new object by indexing.

__iadd__(b)

Self adder.

__imul__(scale)

Inplace multiply.

__isub__(b)

Inplace subtraction.

__mul__(scale)

New object with weight scaled by scale.

__sub__(b)

Subtract objects.

set_values(values)[source]#

Set values.

to_values()[source]#

Access underlying values

to_dataarray()[source]#

Underlying xarray.DataArray.

property attrs#

Attributes of values.

property dims#

Dimensions of values.

property coords#

Coordinates of values.

property name#

Name of values.

property indexes#

Indexes of values.

property sizes#

Sizes of values.

property val_dims#

Names of value dimensions.

property mom_dims#

Names of moment dimensions.

new_like(data=None, *, copy=False, order=None, verify=False, dtype=None)[source]#

Create new object like self, with new data.

Parameters:
  • data (xarray.DataArray) – data for new object

  • copy (bool, optional) – If True, copy the data. If None or False, attempt to use view. Note that False values will be converted to None for numpy versions >2.0. This will be changed to reflect the new behavior of the copy parameter to numpy.array() when the minimum numpy version >2.0.

  • order ({"C", "F", "A", "K"}, optional) – Order argument to numpy.asarray().

  • verify (bool) – If True, make sure data is c-contiguous.

  • dtype (dtype) – Optional dtype for output data.

Returns:

xCentralMoments – New xCentralMoments object with zerod out data.

astype(dtype, *, order=None, casting=None, subok=None, copy=False)[source]#

Underlying data cast to specified type

Parameters:
  • dtype (str or dtype) – Typecode of data-type to cast the array data. Note that a value of None will upcast to np.float64. This is the same behaviour as asarray().

  • order ({"C", "F", "A", "K"}, optional) – Order argument to numpy.asarray().

  • casting ({'no', 'equiv', 'safe', 'same_kind', 'unsafe'}, optional) –

    Controls what kind of data casting may occur.

    • ’no’ means the data types should not be cast at all.

    • ’equiv’ means only byte-order changes are allowed.

    • ’safe’ means only casts which can preserve values are allowed.

    • ’same_kind’ means only safe casts or casts within a kind, like float64 to float32, are allowed.

    • ’unsafe’ (default) means any data conversions may be done.

  • subok (bool, optional) – If True, then sub-classes will be passed-through, otherwise the returned array will be forced to be a base-class array.

  • copy (bool, optional) – By default, astype always returns a newly allocated array. If this is set to False and the dtype requirement is satisfied, the input array is returned insteadof a copy.

Notes

Only numpy.float32 and numpy.float64 dtypes are supported.

moments_to_comoments(*, mom, mom_dims=None, keep_attrs=None)[source]#

Convert moments (mom_ndim=1) to comoments (mom_ndim=2).

Parameters:
  • mom (tuple of int) – Moments for comoments array. Pass a negative value for one of the moments to fill all available moments for that dimensions. For example, if original array has moments m (i.e., values.shape=(..., m + 1)), and pass in mom = (2, -1), then this will be transformed to mom = (2, m - 2).

  • mom_dims (hashable or tuple of hashable) – Name of moment dimensions. Defaults to ("mom_0",) for mom_ndim==1 and (mom_0, mom_1) for mom_ndim==2

  • keep_attrs ({"drop", "identical", "no_conflicts", "drop_conflicts", "override"} or bool, optional) –

    • ‘drop’ or False: empty attrs on returned xarray object.

    • ’identical’: all attrs must be the same on every object.

    • ’no_conflicts’: attrs from all objects are combined, any that have the same name must also have the same value.

    • ’drop_conflicts’: attrs from all objects are combined, any that have the same name but different values are dropped.

    • ’override’ or True: skip comparing and copy attrs from the first object to the result.

assign_weight(weight, copy=True)[source]#

Create object with updated weights

Parameters:
  • weight (array-like, optional) – Optional weights. Can be scalar, 1d array of length args[0].shape[axis] or array of same form as args[0].

  • copy (bool, default True) – If True (default), copy the underlying moments data before update. Otherwise, update weights in in place.

Returns:

output (object) – Same type as self

mean(dim_combined='variable', coords_combined=None)[source]#

Return mean/first moment(s) of data.

var(dim_combined='variable', coords_combined=None)[source]#

Return variance (second central moment) of data.

assign_coords(coords=None, **coords_kwargs)[source]#

Assign coordinates to data and return new object.

assign_attrs(*args, **kwargs)[source]#

Assign attributes to data and return new object.

rename(new_name_or_name_dict=None, **names)[source]#

Rename object.

stack(dimensions=None, *, _reorder=True, _copy=False, _order=None, _verify=False, **dimensions_kwargs)[source]#

Stack dimensions.

Returns:

output (xCentralMoments) – With dimensions stacked.

Examples

>>> import cmomy
>>> rng = cmomy.random.default_rng(0)
>>> vals = xr.DataArray(rng.random((10, 2, 3)), dims=["rec", "dim_0", "dim_1"])
>>> da = xCentralMoments.from_vals(vals, mom=2, dim="rec")
>>> da
<xCentralMoments(val_shape=(2, 3), mom=(2,))>
<xarray.DataArray (dim_0: 2, dim_1: 3, mom_0: 3)> Size: 144B
array([[[10.    ,  0.5205,  0.0452],
        [10.    ,  0.4438,  0.0734],
        [10.    ,  0.5038,  0.1153]],

       [[10.    ,  0.5238,  0.1272],
        [10.    ,  0.628 ,  0.0524],
        [10.    ,  0.412 ,  0.0865]]])
Dimensions without coordinates: dim_0, dim_1, mom_0
>>> da_stack = da.stack(z=["dim_0", "dim_1"])
>>> da_stack
<xCentralMoments(val_shape=(6,), mom=(2,))>
<xarray.DataArray (z: 6, mom_0: 3)> Size: 144B
array([[10.    ,  0.5205,  0.0452],
       [10.    ,  0.4438,  0.0734],
       [10.    ,  0.5038,  0.1153],
       [10.    ,  0.5238,  0.1272],
       [10.    ,  0.628 ,  0.0524],
       [10.    ,  0.412 ,  0.0865]])
Coordinates:
  * z        (z) object 48B MultiIndex
  * dim_0    (z) int64 48B 0 0 0 1 1 1
  * dim_1    (z) int64 48B 0 1 2 0 1 2
Dimensions without coordinates: mom_0

And unstack

>>> da_stack.unstack("z")
<xCentralMoments(val_shape=(2, 3), mom=(2,))>
<xarray.DataArray (dim_0: 2, dim_1: 3, mom_0: 3)> Size: 144B
array([[[10.    ,  0.5205,  0.0452],
        [10.    ,  0.4438,  0.0734],
        [10.    ,  0.5038,  0.1153]],

       [[10.    ,  0.5238,  0.1272],
        [10.    ,  0.628 ,  0.0524],
        [10.    ,  0.412 ,  0.0865]]])
Coordinates:
  * dim_0    (dim_0) int64 16B 0 1
  * dim_1    (dim_1) int64 24B 0 1 2
Dimensions without coordinates: mom_0
unstack(dim=None, fill_value=nan, *, sparse=False, _reorder=True, _copy=False, _order=None, _verify=False)[source]#

Unstack dimensions.

Returns:

output (xCentralMoments) – With dimensions unstacked

set_index(indexes=None, append=False, *, _reorder=True, _copy=False, _order=None, _verify=False, **indexes_kwargs)[source]#

Interface to xarray.DataArray.set_index()

Returns:

output (xCentralMoments) – With new index.

reset_index(dims_or_levels, drop=False, *, _reorder=True, _copy=False, _order=None, _verify=False)[source]#

Interface to xarray.DataArray.reset_index().

drop_vars(names, *, errors='raise', _reorder=True, _copy=False, _order=None, _verify=False)[source]#

Interface to xarray.DataArray.drop_vars()

swap_dims(dims_dict=None, _reorder=True, _copy=False, _order=None, _verify=False, **dims_kwargs)[source]#

Interface to xarray.DataArray.swap_dims().

sel(indexers=None, method=None, tolerance=None, drop=False, *, _reorder=False, _copy=False, _order=None, _verify=False, **indexers_kws)[source]#

Select subset of data.

Returns:

output (xCentralMoments) – With dimensions unstacked

Examples

>>> import cmomy
>>> rng = cmomy.random.default_rng(0)
>>> da = cmomy.CentralMoments.from_vals(
...     rng.random((10, 3)), axis=0, mom=2
... ).to_x(dims="x", coords=dict(x=list("abc")))
>>> da
<xCentralMoments(val_shape=(3,), mom=(2,))>
<xarray.DataArray (x: 3, mom_0: 3)> Size: 72B
array([[10.    ,  0.5248,  0.1106],
       [10.    ,  0.5688,  0.0689],
       [10.    ,  0.5094,  0.1198]])
Coordinates:
  * x        (x) <U1 12B 'a' 'b' 'c'
Dimensions without coordinates: mom_0

Select by value

>>> da.sel(x="a")
<xCentralMoments(val_shape=(), mom=(2,))>
<xarray.DataArray (mom_0: 3)> Size: 24B
array([10.    ,  0.5248,  0.1106])
Coordinates:
    x        <U1 4B 'a'
Dimensions without coordinates: mom_0
>>> da.sel(x=["a", "c"])
<xCentralMoments(val_shape=(2,), mom=(2,))>
<xarray.DataArray (x: 2, mom_0: 3)> Size: 48B
array([[10.    ,  0.5248,  0.1106],
       [10.    ,  0.5094,  0.1198]])
Coordinates:
  * x        (x) <U1 8B 'a' 'c'
Dimensions without coordinates: mom_0

Select by position

>>> da.isel(x=0)
<xCentralMoments(val_shape=(), mom=(2,))>
<xarray.DataArray (mom_0: 3)> Size: 24B
array([10.    ,  0.5248,  0.1106])
Coordinates:
    x        <U1 4B 'a'
Dimensions without coordinates: mom_0
>>> da.isel(x=[0, 1])
<xCentralMoments(val_shape=(2,), mom=(2,))>
<xarray.DataArray (x: 2, mom_0: 3)> Size: 48B
array([[10.    ,  0.5248,  0.1106],
       [10.    ,  0.5688,  0.0689]])
Coordinates:
  * x        (x) <U1 8B 'a' 'b'
Dimensions without coordinates: mom_0
isel(indexers=None, drop=False, missing_dims='raise', *, _reorder=False, _copy=False, _order=None, _verify=False, **indexers_kws)[source]#

Select subset of data by position.

Returns:

output (xCentralMoments) – With dimensions unstacked

transpose(*dims, transpose_coords=True, missing_dims='raise', _copy=False, _order=None, _verify=False)[source]#

Transpose dimensions of data.

Notes

mom_dims will always be put at the end of the array.

to_centralmoments(copy=False)[source]#

Create a CentralMoments object from xCentralMoments.

to_c(copy=False)[source]#

Alias to to_centralmoments()

classmethod from_centralmoments(obj, *, dims=None, attrs=None, coords=None, name=None, indexes=None, mom_dims=None, template=None, copy=False)[source]#

Create and xCentralMoments object from CentralMoments.

Parameters:
  • obj (CentralMoments) – Input object to be converted.

  • dims (hashable or sequence of hashable) –

    Dimension of resulting xarray.DataArray.

    • If len(dims) == self.ndim, then dims specifies all dimensions.

    • If len(dims) == self.val_ndim, dims = dims + mom_dims

    Default to ('dim_0', 'dim_1', ...)

  • attrs (mapping) – Attributes of output

  • coords (mapping) – Coordinates of output

  • name (hashable) – Name of output

  • indexes (Any) – indexes attribute. This is ignored.

  • template (DataArray) – If present, output will have attributes of template. Overrides other options.

  • copy (bool, optional) – If True, copy the data. If None or False, attempt to use view. Note that False values will be converted to None for numpy versions >2.0. This will be changed to reflect the new behavior of the copy parameter to numpy.array() when the minimum numpy version >2.0.

Returns:

output (xCentralMoments)

push_data(data, *, order=None, parallel=False)[source]#

Push data object to moments.

Parameters:
  • data (array-like) – Accumulation array of same form as self.data

  • order ({"C", "F", "A", "K"}, optional) – Order argument to numpy.asarray().

  • parallel (bool, default True) – flags to numba.njit

Returns:

output (xCentralMoments) – Same object with pushed data.

__add__(b)[source]#

Add objects to new object.

__getitem__(key)[source]#

Get new object by indexing.

Note that only objects with the same moment(s) shape are allowed.

If you want to extract data in general, use self.to_values()[….].

__iadd__(b)[source]#

Self adder.

__imul__(scale)[source]#

Inplace multiply.

__isub__(b)[source]#

Inplace subtraction.

__mul__(scale)[source]#

New object with weight scaled by scale.

__sub__(b)[source]#

Subtract objects.

cmom()[source]#

Central moments.

Strict central moments of the form

\[\text{cmom[..., n, m]} = \langle (x - \langle x \rangle)^n (y - \langle y \rangle)^m \rangle\]

where

\[\langle x \rangle = \sum_i w_i x_i / \sum_i w_i\]
Returns:

output (ndarray or DataArray)

copy()[source]#

Create a new object with copy of data.

Parameters:

**copy_kws – passed to parameter copy_kws in method new_like()

Returns:

output (object) – Same type as calling class. Object with same attributes as caller, but with new underlying data.

See also

new_like, zeros_like

property data#

Accessor to numpy array underlying data.

By convention data has the following meaning for the moments indexes

  • data[…,i=0,j=0], weight

  • data[…,i=1,j=0]], if only one moment index is one and all others zero, then this is the average value of the variable with unit index.

  • all other cases, the central moments <(x0-<x0>)**i0 * (x1 - <x1>)**i1 * …>

property dtype#

self.data.dtype.

fill(value=0)[source]#

Fill data with value.

Parameters:

value (scalar) – Value to insert into self.data

Returns:

self (object) – Same type as calling class. Same object as caller with data filled with values

property mom#

Number of moments.

property mom_ndim#

Length of moments.

if mom_ndim == 1, then single variable moments if mom_ndim == 2, then co-moments.

property mom_shape#

Shape of moments part.

property ndim#

self.data.ndim.

pipe(func_or_method, *args, _reorder=True, _copy=None, _order=None, _verify=False, **kwargs)[source]#

Apply func_or_method to underlying data and wrap results in class:cmomy.xCentralMoments object.

This is useful for calling any not implemented methods on numpy.ndarray or xarray.DataArray data.

Note that a ValueError will be raised if last dimension(s) do not have the same shape as self.mom_shape.

Parameters:
  • func_or_method (str or callable()) – If callable, then apply values = func_or_method(self.to_values(), *args, **kwargs). If string is passed, then values = getattr(self.to_values(), func_or_method)(*args, **kwargs).

  • *args (Any) – Extra positional arguments to func_or_method

  • _reorder (bool, default True) – If True, reorder the data such that mom_dims are last. Only applicable for DataArray like underlying data.

  • _copy (bool, default False) – If True, copy the resulting data. Otherwise, try to use a view.

  • _order (str, optional) – Array order to apply to output.

  • _verify (bool, default False)

  • **kwargs (Any) – Extra keyword arguments to func_or_method

Returns:

output (object) – New object after func_or_method is applies to self.to_values()

Notes

Use leading underscore for _order, _copy, etc to avoid possible name clashes with func_or_method.

push_datas(datas, *, axis=MISSING, dim=MISSING, order=None, parallel=None)[source]#

Push and reduce multiple average central moments.

Parameters:
  • datas (array-like, xarray.DataArray) – Collection of accumulation arrays to push onto self. This should have shape like (nrec,) + self.shape if axis=0, where nrec is the number of data objects to sum.

  • axis (int, optional) – Axis to reduce along. Note that negative values are relative to data.ndim - mom_ndim. It is assumed that the last dimensions are for moments. For example, if data.shape == (1,2,3) with mom_ndim=1, axis = -1 `` would be equivalent to ``axis = 1. Defaults to axis=-1.

  • dim (hashable) – Dimension to reduce along.

  • order ({"C", "F", "A", "K"}, optional) – Order argument to numpy.asarray().

  • parallel (bool, default True) – flags to numba.njit

Returns:

output (xCentralMoments) – Same object with pushed data.

rmom()[source]#

Raw moments.

\[\text{rmom[..., n, m]} = \langle x^n y^m \rangle\]

where

\[\langle x \rangle = \sum_i w_i x_i / \sum_i w_i\]
Returns:

raw_moments (ndarray or DataArray)

property shape#

self.data.shape.

std()[source]#

Standard deviation.

to_numpy()[source]#

Access to numpy array underlying class.

to_raw(*, weight=None)[source]#

Raw moments accumulation array.

\[\begin{split}\text{raw[..., n, m]} = \begin{cases} \text{weight} & n = m = 0 \\ \langle x^n y ^m \rangle & \text{otherwise} \end{cases}\end{split}\]

where

\[\langle x \rangle = \sum_i w_i x_i / \sum_i w_i\]
Returns:

raw (ndarray or DataArray)

property val_ndim#

Number of value dimensions.

property val_shape#

Shape of values dimensions.

That is shape less moments dimensions.

property values#

Access underlying central moments array.

weight()[source]#

Weight data.

zero()[source]#

Zero out underlying data.

Returns:

self (object) – Same type as calling class. Same object with data filled with zeros.

See also

fill

zeros_like()[source]#

Create new empty object like self.

Returns:

output (object) – Object with same attributes as caller, but with data set to zero.

See also

new_like

push_val(x, *y, weight=None, order=None, parallel=False)[source]#

Push single sample to central moments.

Parameters:
  • x (array) – Values to push onto self.

  • *y (array-like, optional) – Additional Values (needed if mom_ndim > 1)

  • weight (int, float, array-like, optional) – Weight of each sample. If scalar, broadcast w.shape to x0.shape.

  • order ({"C", "F", "A", "K"}, optional) – Order argument to numpy.asarray().

  • parallel (bool, default True) – flags to numba.njit

Returns:

output (xCentralMoments) – Same object with pushed data.

Notes

Array x0 should have same shape as self.val_shape.

push_vals(x, *y, weight=None, axis=MISSING, dim=MISSING, order=None, parallel=None)[source]#

Push multiple samples to central moments.

Parameters:
  • x (array) – Value to reduce.

  • *y (array-like) – Additional array (if self.mom_ndim == 2).

  • weight (int, float, array-like, optional) – Weight of each sample. If scalar, broadcast to x0.shape

  • axis (int) – Axis to reduce along.

  • dim (hashable) – Dimension to reduce along.

  • order ({"C", "F", "A", "K"}, optional) – Order argument to numpy.asarray().

  • parallel (bool, default True) – flags to numba.njit

Returns:

output (xCentralMoments) – Same object with pushed data.

randsamp_freq(*, axis=MISSING, dim=MISSING, nrep=None, nsamp=None, indices=None, freq=None, check=False, rng=None)[source]#

Interface to resample.randsamp_freq()

Parameters:
  • axis (int) – Axis that will be resampled. This is used to calculate ndat.

  • dim (str) – Dimension that will be resampled.

  • nrep (int) – Number of resample replicates.

  • nsamp (int) – Number of samples in a single resampled replicate. Defaults to size of data along sampled axis.

  • indices (array of int) – Array of shape (nrep, size). If passed, create freq from indices. See randsamp_freq().

  • freq (array of int) – Array of shape (nrep, size) where nrep is the number of replicates and size = self.shape[axis]. freq is the weight that each sample contributes to resamples values. See randsamp_freq()

  • check (bool, default False) – If True, perform sanity checks.

  • rng (Generator) – Random number generator object. Defaults to output of default_rng().

Returns:

freq (ndarray) – Frequency array

resample_and_reduce(*, freq, dim=MISSING, axis=MISSING, rep_dim='rep', parallel=None, order=None, keep_attrs=None)[source]#

Bootstrap resample and reduce.

Parameters:
  • freq (array of int) – Array of shape (nrep, size) where nrep is the number of replicates and size = self.shape[axis]. freq is the weight that each sample contributes to resamples values. See randsamp_freq()

  • dim (hashable) – Dimension to reduce along.

  • axis (int, optional) – Axis to reduce along. Note that negative values are relative to data.ndim - mom_ndim. It is assumed that the last dimensions are for moments. For example, if data.shape == (1,2,3) with mom_ndim=1, axis = -1 `` would be equivalent to ``axis = 1. Defaults to axis=-1.

  • rep_dim (hashable) – Name of new ‘replicated’ dimension:

  • parallel (bool, default True) – flags to numba.njit

  • order ({"C", "F", "A", "K"}, optional) – Order argument to numpy.asarray().

  • keep_attrs ({"drop", "identical", "no_conflicts", "drop_conflicts", "override"} or bool, optional) –

    • ‘drop’ or False: empty attrs on returned xarray object.

    • ’identical’: all attrs must be the same on every object.

    • ’no_conflicts’: attrs from all objects are combined, any that have the same name must also have the same value.

    • ’drop_conflicts’: attrs from all objects are combined, any that have the same name but different values are dropped.

    • ’override’ or True: skip comparing and copy attrs from the first object to the result.

Returns:

output (object) – Instance of calling class. Note that new object will have (...,shape[axis-1], shape[axis+1], ..., nrep, mom0, ...), where nrep = freq.shape[0].

See also

reduce

randsamp_freq

random frequency sample

freq_to_indices

convert frequency sample to index sample

indices_to_freq

convert index sample to frequency sample

resample_data

method to perform resampling

Examples

>>> import cmomy
>>> rng = cmomy.random.default_rng(0)
>>> da = cmomy.CentralMoments.from_vals(
...     rng.random((10, 3)),
...     mom=3,
...     axis=0,
... ).to_x(dims="rec")
>>> da
<xCentralMoments(val_shape=(3,), mom=(3,))>
<xarray.DataArray (rec: 3, mom_0: 4)> Size: 96B
array([[ 1.0000e+01,  5.2485e-01,  1.1057e-01, -4.6282e-03],
       [ 1.0000e+01,  5.6877e-01,  6.8876e-02, -1.2745e-02],
       [ 1.0000e+01,  5.0944e-01,  1.1978e-01, -1.4644e-02]])
Dimensions without coordinates: rec, mom_0

Note that for reproducible results, must set numba random seed as well

>>> freq = da.randsamp_freq(dim="rec", nrep=5)
>>> da_resamp = da.resample_and_reduce(
...     dim="rec",
...     freq=freq,
... )
>>> da_resamp
<xCentralMoments(val_shape=(5,), mom=(3,))>
<xarray.DataArray (rep: 5, mom_0: 4)> Size: 160B
array([[ 3.0000e+01,  5.0944e-01,  1.1978e-01, -1.4644e-02],
       [ 3.0000e+01,  5.3435e-01,  1.0038e-01, -1.2329e-02],
       [ 3.0000e+01,  5.2922e-01,  1.0360e-01, -1.6009e-02],
       [ 3.0000e+01,  5.5413e-01,  8.3204e-02, -1.1267e-02],
       [ 3.0000e+01,  5.4899e-01,  8.6627e-02, -1.5407e-02]])
Dimensions without coordinates: rep, mom_0

Alternatively, we can resample and reduce

>>> indices = cmomy.resample.freq_to_indices(freq)
>>> da.sel(rec=xr.DataArray(indices, dims=["rep", "rec"])).reduce(dim="rec")
<xCentralMoments(val_shape=(5,), mom=(3,))>
<xarray.DataArray (rep: 5, mom_0: 4)> Size: 160B
array([[ 3.0000e+01,  5.0944e-01,  1.1978e-01, -1.4644e-02],
       [ 3.0000e+01,  5.3435e-01,  1.0038e-01, -1.2329e-02],
       [ 3.0000e+01,  5.2922e-01,  1.0360e-01, -1.6009e-02],
       [ 3.0000e+01,  5.5413e-01,  8.3204e-02, -1.1267e-02],
       [ 3.0000e+01,  5.4899e-01,  8.6627e-02, -1.5407e-02]])
Dimensions without coordinates: rep, mom_0
reduce(*, by=None, axis=MISSING, order=None, parallel=None, coords_policy='first', dim=MISSING, group_dim=None, groups=None, keep_attrs=None)[source]#

Create new object reduce along axis.

Parameters:
  • by (ndarray or DataArray or str or iterable of str, optional) – If None, reduce over entire dim. Otherwise, reduce by group. If ndarray, use the unique values. If DataArray, use unique values and rename dimension by.name. If str or Iterable of str, Create grouper from these named coordinates.

  • axis (int, optional) – Axis to reduce along. Note that negative values are relative to data.ndim - mom_ndim. It is assumed that the last dimensions are for moments. For example, if data.shape == (1,2,3) with mom_ndim=1, axis = -1 `` would be equivalent to ``axis = 1. Defaults to axis=-1.

  • order ({"C", "F", "A", "K"}, optional) – Order argument to numpy.asarray().

  • parallel (bool, default True) – flags to numba.njit

  • coords_policy ({'first', 'last', 'group', None}) –

    Policy for handling coordinates along dim if by is specified for DataArray data. If no coordinates do nothing, otherwise use:

    • ’first’: select first value of coordinate for each block.

    • ’last’: select last value of coordinate for each block.

    • ’group’: Assign unique groups from group_idx to dim

    • None: drop any coordinates.

    Note that if coords_policy is one of first or last, parameter groups will be ignored.

  • dim (hashable) – Dimension to reduce along.

  • group_dim (str, optional) – Name of the output group dimension. Defaults to dim.

  • groups (sequence, optional) – Sequence of length by.max() + 1 to assign as coordinates for group_dim.

  • keep_attrs ({"drop", "identical", "no_conflicts", "drop_conflicts", "override"} or bool, optional) –

    • ‘drop’ or False: empty attrs on returned xarray object.

    • ’identical’: all attrs must be the same on every object.

    • ’no_conflicts’: attrs from all objects are combined, any that have the same name must also have the same value.

    • ’drop_conflicts’: attrs from all objects are combined, any that have the same name but different values are dropped.

    • ’override’ or True: skip comparing and copy attrs from the first object to the result.

Returns:

output (xCentralMoments) – If by is None, reduce all samples along axis. Otherwise, reduce for each unique value of by. In this case, output will have shape (..., shape[axis-1], shape[axis+1], ..., ngroup, mom0, ...) where ngroups = np.max(by) + 1 is the number of unique positive values in by.

Notes

This is a new feature, and subject to change. In particular, the current implementation is not smart about the output coordinates and dimensions, and is inconsistent with xarray.DataArray.groupby(). It is up the the user to manipulate the dimensions/coordinates. output dimensions and coords simplistically.

Examples

>>> import cmomy
>>> rng = cmomy.random.default_rng(0)
>>> vals = xr.DataArray(rng.random((10, 2, 3)), dims=["rec", "dim_0", "dim_1"])
>>> da = xCentralMoments.from_vals(vals, dim="rec", mom=2)
>>> da.reduce(dim="dim_0")
<xCentralMoments(val_shape=(3,), mom=(2,))>
<xarray.DataArray (dim_1: 3, mom_0: 3)> Size: 72B
array([[20.    ,  0.5221,  0.0862],
       [20.    ,  0.5359,  0.0714],
       [20.    ,  0.4579,  0.103 ]])
Dimensions without coordinates: dim_1, mom_0
block(block_size, *, dim=MISSING, axis=MISSING, block_dim=None, order=None, parallel=None, keep_attrs=None, coords_policy='first')[source]#
Parameters:
  • block_size (int) – Number of observations to include in a given block.

  • dim (hashable) – Dimension to reduce along.

  • axis (int) – Axis to reduce along.

  • block_dim (str, optional) – Name of blocked dimension. Defaults to dim.

  • order ({"C", "F", "A", "K"}, optional) – Order argument to numpy.asarray().

  • parallel (bool, default True) – flags to numba.njit

  • keep_attrs ({"drop", "identical", "no_conflicts", "drop_conflicts", "override"} or bool, optional) –

    • ‘drop’ or False: empty attrs on returned xarray object.

    • ’identical’: all attrs must be the same on every object.

    • ’no_conflicts’: attrs from all objects are combined, any that have the same name must also have the same value.

    • ’drop_conflicts’: attrs from all objects are combined, any that have the same name but different values are dropped.

    • ’override’ or True: skip comparing and copy attrs from the first object to the result.

  • coords_policy ({'first', 'last', 'group', None}) –

    Policy for handling coordinates along dim if by is specified for DataArray data. If no coordinates do nothing, otherwise use:

    • ’first’: select first value of coordinate for each block.

    • ’last’: select last value of coordinate for each block.

    • ’group’: Assign unique groups from group_idx to dim

    • None: drop any coordinates.

    Note that if coords_policy is one of first or last, parameter groups will be ignored.

Returns:

output (xCentralMoments) – Object with block averaging.

Examples

>>> import cmomy
>>> rng = cmomy.random.default_rng(0)
>>> x = rng.random((10, 10))
>>> da = cmomy.CentralMoments.from_vals(x, mom=2, axis=0).to_x()
>>> da
<xCentralMoments(val_shape=(10,), mom=(2,))>
<xarray.DataArray (dim_0: 10, mom_0: 3)> Size: 240B
array([[10.    ,  0.6247,  0.0583],
       [10.    ,  0.3938,  0.0933],
       [10.    ,  0.425 ,  0.1003],
       [10.    ,  0.5   ,  0.117 ],
       [10.    ,  0.5606,  0.0446],
       [10.    ,  0.5612,  0.0861],
       [10.    ,  0.531 ,  0.0731],
       [10.    ,  0.8403,  0.0233],
       [10.    ,  0.5097,  0.103 ],
       [10.    ,  0.5368,  0.085 ]])
Dimensions without coordinates: dim_0, mom_0
>>> da.block(block_size=5, dim="dim_0")
<xCentralMoments(val_shape=(2,), mom=(2,))>
<xarray.DataArray (dim_0: 2, mom_0: 3)> Size: 48B
array([[50.    ,  0.5008,  0.0899],
       [50.    ,  0.5958,  0.0893]])
Dimensions without coordinates: dim_0, mom_0

This is equivalent to

>>> cmomy.CentralMoments.from_vals(x.reshape(2, 50), mom=2, axis=1).to_x()
<xCentralMoments(val_shape=(2,), mom=(2,))>
<xarray.DataArray (dim_0: 2, mom_0: 3)> Size: 48B
array([[50.    ,  0.5268,  0.0849],
       [50.    ,  0.5697,  0.0979]])
Dimensions without coordinates: dim_0, mom_0

The coordinate policy can be useful to keep coordinates:

>>> da2 = da.assign_coords(dim_0=range(10))
>>> da2.block(5, dim="dim_0", coords_policy="first")
<xCentralMoments(val_shape=(2,), mom=(2,))>
<xarray.DataArray (dim_0: 2, mom_0: 3)> Size: 48B
array([[50.    ,  0.5008,  0.0899],
       [50.    ,  0.5958,  0.0893]])
Coordinates:
  * dim_0    (dim_0) int64 16B 0 5
Dimensions without coordinates: mom_0
>>> da2.block(5, dim="dim_0", coords_policy="last")
<xCentralMoments(val_shape=(2,), mom=(2,))>
<xarray.DataArray (dim_0: 2, mom_0: 3)> Size: 48B
array([[50.    ,  0.5008,  0.0899],
       [50.    ,  0.5958,  0.0893]])
Coordinates:
  * dim_0    (dim_0) int64 16B 4 9
Dimensions without coordinates: mom_0
>>> da2.block(5, dim="dim_0", coords_policy=None)
<xCentralMoments(val_shape=(2,), mom=(2,))>
<xarray.DataArray (dim_0: 2, mom_0: 3)> Size: 48B
array([[50.    ,  0.5008,  0.0899],
       [50.    ,  0.5958,  0.0893]])
Dimensions without coordinates: dim_0, mom_0
classmethod zeros(*, mom, val_shape=None, dtype=None, order=None, dims=None, mom_dims=None, attrs=None, coords=None, name=None, indexes=None, template=None)[source]#

Create a new base object.

Parameters:
  • mom (int or tuple of int) – Order or moments. If integer or length one tuple, then moments are for a single variable. If length 2 tuple, then comoments of two variables

  • val_shape (tuple) – Shape of values part of data. That is, the non-moment dimensions.

  • dtype (dtype) – Optional dtype for output data.

  • order ({"C", "F", "A", "K"}, optional) – Order argument to numpy.asarray().

  • dims (hashable or sequence of hashable) –

    Dimension of resulting xarray.DataArray.

    • If len(dims) == self.ndim, then dims specifies all dimensions.

    • If len(dims) == self.val_ndim, dims = dims + mom_dims

    Default to ('dim_0', 'dim_1', ...)

  • attrs (mapping) – Attributes of output

  • coords (mapping) – Coordinates of output

  • name (hashable) – Name of output

  • indexes (Any) – indexes attribute. This is ignored.

  • template (DataArray) – If present, output will have attributes of template. Overrides other options.

Returns:

output (xCentralMoments) – New instance with zero values.

See also

numpy.zeros

classmethod from_vals(x, *y, mom, weight=None, axis=MISSING, dim=MISSING, mom_dims=None, order=None, parallel=None, dtype=None, keep_attrs=None)[source]#

Create from observations/values.

Parameters:
  • x (DataArray) – Values to reduce.

  • *y (array-like or DataArray) – Additional values (needed if len(mom)==2).

  • weight (scalar or array-like or DataArray, optional) – Optional weight. If scalar or array, attempt to broadcast to x0.shape

  • axis (int) – Axis to reduce along.

  • dim (hashable) – Dimension to reduce along.

  • mom (int or tuple of int) – Order or moments. If integer or length one tuple, then moments are for a single variable. If length 2 tuple, then comoments of two variables

  • order ({"C", "F", "A", "K"}, optional) – Order argument to numpy.asarray().

  • parallel (bool, default True) – flags to numba.njit

  • dtype (dtype) – Optional dtype for output data.

  • mom_dims (hashable or tuple of hashable) – Name of moment dimensions. Defaults to ("mom_0",) for mom_ndim==1 and (mom_0, mom_1) for mom_ndim==2

  • keep_attrs ({"drop", "identical", "no_conflicts", "drop_conflicts", "override"} or bool, optional) –

    • ‘drop’ or False: empty attrs on returned xarray object.

    • ’identical’: all attrs must be the same on every object.

    • ’no_conflicts’: attrs from all objects are combined, any that have the same name must also have the same value.

    • ’drop_conflicts’: attrs from all objects are combined, any that have the same name but different values are dropped.

    • ’override’ or True: skip comparing and copy attrs from the first object to the result.

Returns:

output (xCentralMoments)

classmethod from_resample_vals(x, *y, mom, freq, weight=None, axis=MISSING, dim=MISSING, order=None, parallel=None, dtype=None, mom_dims=None, rep_dim='rep', keep_attrs=True)[source]#

Create from resample observations/values.

This effectively resamples x.

Parameters:
  • x (DataArray) – For moments, pass single array-like objects x=x0. For comoments, pass tuple of array-like objects x=(x0, x1).

  • *y (array-like or DataArray) – Additional values (needed if len(mom) > 1).

  • mom (int or tuple of int) – Order or moments. If integer or length one tuple, then moments are for a single variable. If length 2 tuple, then comoments of two variables

  • freq (array of int) – Array of shape (nrep, size) where nrep is the number of replicates and size = self.shape[axis]. freq is the weight that each sample contributes to resamples values. See randsamp_freq()

  • weight (array-like, optional) – Optional weights. Can be scalar, 1d array of length args[0].shape[axis] or array of same form as args[0].

  • axis (int) – Axis to reduce along.

  • dim (hashable) – Dimension to reduce along.

  • full_output (bool) – If True, also return freq array

  • order ({"C", "F", "A", "K"}, optional) – Order argument to numpy.asarray().

  • parallel (bool, default True) – flags to numba.njit

  • dtype (dtype) – Optional dtype for output data.

Returns:

out (xCentralMoments) – Instance of calling class

classmethod from_raw(raw, *, mom_ndim)[source]#

Create object from raw moment data.

raw[…, i, j] = <x**i y**j>. raw[…, 0, 0] = weight

Parameters:
  • raw (xarray.DataArray) – Raw moment array.

  • mom_ndim ({1, 2}) – Value indicates if moments (mom_ndim = 1) or comoments (mom_ndim=2).

Returns:

output (xCentralMoments)

Notes

Weights are taken from raw[...,0, 0]. Using raw moments can result in numerical issues, especially for higher moments. Use with care.