Struct numpy::array::PyArray

source ·
pub struct PyArray<T, D>(/* private fields */);
Expand description

A safe, statically-typed wrapper for NumPy’s ndarray class.

§Memory location

These methods transfers ownership of the Rust allocation into a suitable Python object and uses the memory as the internal buffer backing the NumPy array.

Please note that some destructive methods like resize will fail when used with this kind of array as NumPy cannot reallocate the internal buffer.

These methods allocate memory in Python’s private heap via NumPy’s API.

In both cases, PyArray is managed by Python so it can neither be moved from nor deallocated manually.

§References

Like new, all constructor methods of PyArray return a shared reference &PyArray instead of an owned value. This design follows PyO3’s ownership concept, i.e. the return value is GIL-bound owning reference into Python’s heap.

§Element type and dimensionality

PyArray has two type parametes T and D. T represents the type of its elements, e.g. f32 or [PyObject]. D represents its dimensionality, e.g Ix2 or IxDyn.

Element types are Rust types which implement the Element trait. Dimensions are represented by the ndarray::Dimension trait.

Typically, Ix1, Ix2, ... are used for fixed dimensionality arrays, and IxDyn is used for dynamic dimensionality arrays. Type aliases for combining PyArray with these types are provided, e.g. PyArray1 or PyArrayDyn.

To specify concrete dimension like 3×4×5, types which implement the ndarray::IntoDimension trait are used. Typically, this means arrays like [3, 4, 5] or tuples like (3, 4, 5).

§Example

use numpy::{PyArray, PyArrayMethods};
use ndarray::{array, Array};
use pyo3::Python;

Python::with_gil(|py| {
    let pyarray = PyArray::arange_bound(py, 0., 4., 1.).reshape([2, 2]).unwrap();
    let array = array![[3., 4.], [5., 6.]];

    assert_eq!(
        array.dot(&pyarray.readonly().as_array()),
        array![[8., 15.], [12., 23.]]
    );
});

Implementations§

source§

impl<T, D> PyArray<T, D>

source

pub fn as_untyped(&self) -> &PyUntypedArray

Access an untyped representation of this array.

source

pub fn to_owned(&self) -> Py<Self>

👎Deprecated since 0.21.0: use Bound::unbind() instead

Turn &PyArray<T,D> into Py<PyArray<T,D>>, i.e. a pointer into Python’s heap which is independent of the GIL lifetime.

This method can be used to avoid lifetime annotations of function arguments or return values.

§Example
use numpy::{PyArray1, PyArrayMethods};
use pyo3::{Py, Python};

let array: Py<PyArray1<f64>> = Python::with_gil(|py| {
    PyArray1::zeros_bound(py, 5, false).unbind()
});

Python::with_gil(|py| {
    assert_eq!(array.bind(py).readonly().as_slice().unwrap(), [0.0; 5]);
});
source

pub unsafe fn from_owned_ptr<'py>( py: Python<'py>, ptr: *mut PyObject ) -> &'py Self

Constructs a reference to a PyArray from a raw pointer to a Python object.

§Safety

This is a wrapper around [pyo3::FromPyPointer::from_owned_ptr_or_opt] and inherits its safety contract.

source

pub unsafe fn from_borrowed_ptr<'py>( py: Python<'py>, ptr: *mut PyObject ) -> &'py Self

Constructs a reference to a PyArray from a raw point to a Python object.

§Safety

This is a wrapper around [pyo3::FromPyPointer::from_borrowed_ptr_or_opt] and inherits its safety contract.

source

pub fn data(&self) -> *mut T

Returns a pointer to the first element of the array.

source§

impl<T: Element, D: Dimension> PyArray<T, D>

source

pub fn dims(&self) -> D

Same as shape, but returns D instead of &[usize].

source

pub unsafe fn new<'py, ID>(py: Python<'py>, dims: ID, is_fortran: bool) -> &Self
where ID: IntoDimension<Dim = D>,

👎Deprecated since 0.21.0: will be replaced by PyArray::new_bound in the future

Deprecated form of PyArray<T, D>::new_bound

§Safety

Same as PyArray<T, D>::new_bound

source

pub unsafe fn new_bound<'py, ID>( py: Python<'py>, dims: ID, is_fortran: bool ) -> Bound<'py, Self>
where ID: IntoDimension<Dim = D>,

Creates a new uninitialized NumPy array.

If is_fortran is true, then it has Fortran/column-major order, otherwise it has C/row-major order.

§Safety

The returned array will always be safe to be dropped as the elements must either be trivially copyable (as indicated by <T as Element>::IS_COPY) or be pointers into Python’s heap, which NumPy will automatically zero-initialize.

However, the elements themselves will not be valid and should be initialized manually using raw pointers obtained via uget_raw. Before that, all methods which produce references to the elements invoke undefined behaviour. In particular, zero-initialized pointers are not valid instances of PyObject.

§Example
use numpy::prelude::*;
use numpy::PyArray3;
use pyo3::Python;

Python::with_gil(|py| {
    let arr = unsafe {
        let arr = PyArray3::<i32>::new_bound(py, [4, 5, 6], false);

        for i in 0..4 {
            for j in 0..5 {
                for k in 0..6 {
                    arr.uget_raw([i, j, k]).write((i * j * k) as i32);
                }
            }
        }

        arr
    };

    assert_eq!(arr.shape(), &[4, 5, 6]);
});
source

pub unsafe fn borrow_from_array<'py, S>( array: &ArrayBase<S, D>, container: &'py PyAny ) -> &'py Self
where S: Data<Elem = T>,

👎Deprecated since 0.21.0: will be replaced by PyArray::borrow_from_array_bound in the future
source

pub unsafe fn borrow_from_array_bound<'py, S>( array: &ArrayBase<S, D>, container: Bound<'py, PyAny> ) -> Bound<'py, Self>
where S: Data<Elem = T>,

Creates a NumPy array backed by array and ties its ownership to the Python object container.

§Safety

container is set as a base object of the returned array which must not be dropped until container is dropped. Furthermore, array must not be reallocated from the time this method is called and until container is dropped.

§Example
#[pyclass]
struct Owner {
    array: Array1<f64>,
}

#[pymethods]
impl Owner {
    #[getter]
    fn array<'py>(this: Bound<'py, Self>) -> Bound<'py, PyArray1<f64>> {
        let array = &this.borrow().array;

        // SAFETY: The memory backing `array` will stay valid as long as this object is alive
        // as we do not modify `array` in any way which would cause it to be reallocated.
        unsafe { PyArray1::borrow_from_array_bound(array, this.into_any()) }
    }
}
source

pub fn zeros<'py, ID>(py: Python<'py>, dims: ID, is_fortran: bool) -> &Self
where ID: IntoDimension<Dim = D>,

👎Deprecated since 0.21.0: will be replaced by PyArray::zeros_bound in the future

Deprecated form of PyArray<T, D>::zeros_bound

source

pub fn zeros_bound<ID>( py: Python<'_>, dims: ID, is_fortran: bool ) -> Bound<'_, Self>
where ID: IntoDimension<Dim = D>,

Construct a new NumPy array filled with zeros.

If is_fortran is true, then it has Fortran/column-major order, otherwise it has C/row-major order.

For arrays of Python objects, this will fill the array with valid pointers to zero-valued Python integer objects.

See also numpy.zeros and PyArray_Zeros.

§Example
use numpy::{PyArray2, PyArrayMethods};
use pyo3::Python;

Python::with_gil(|py| {
    let pyarray = PyArray2::<usize>::zeros_bound(py, [2, 2], true);

    assert_eq!(pyarray.readonly().as_slice().unwrap(), [0; 4]);
});
source

pub unsafe fn as_slice(&self) -> Result<&[T], NotContiguousError>

Returns an immutable view of the internal data as a slice.

§Safety

Calling this method is undefined behaviour if the underlying array is aliased mutably by other instances of PyArray or concurrently modified by Python or other native code.

Please consider the safe alternative PyReadonlyArray::as_slice.

source

pub unsafe fn as_slice_mut(&self) -> Result<&mut [T], NotContiguousError>

Returns a mutable view of the internal data as a slice.

§Safety

Calling this method is undefined behaviour if the underlying array is aliased immutably or mutably by other instances of PyArray or concurrently modified by Python or other native code.

Please consider the safe alternative PyReadwriteArray::as_slice_mut.

source

pub fn from_owned_array<'py>(py: Python<'py>, arr: Array<T, D>) -> &'py Self

👎Deprecated since 0.21.0: will be replaced by PyArray::from_owned_array_bound in the future
source

pub fn from_owned_array_bound( py: Python<'_>, arr: Array<T, D> ) -> Bound<'_, Self>

Constructs a NumPy from an ndarray::Array

This method uses the internal Vec of the ndarray::Array as the base object of the NumPy array.

§Example
use numpy::{PyArray, PyArrayMethods};
use ndarray::array;
use pyo3::Python;

Python::with_gil(|py| {
    let pyarray = PyArray::from_owned_array_bound(py, array![[1, 2], [3, 4]]);

    assert_eq!(pyarray.readonly().as_array(), array![[1, 2], [3, 4]]);
});
source

pub unsafe fn get(&self, index: impl NpyIndex<Dim = D>) -> Option<&T>

Get a reference of the specified element if the given index is valid.

§Safety

Calling this method is undefined behaviour if the underlying array is aliased mutably by other instances of PyArray or concurrently modified by Python or other native code.

Consider using safe alternatives like PyReadonlyArray::get.

§Example
use numpy::{PyArray, PyArrayMethods};
use pyo3::Python;

Python::with_gil(|py| {
    let pyarray = PyArray::arange_bound(py, 0, 16, 1).reshape([2, 2, 4]).unwrap();

    assert_eq!(unsafe { *pyarray.get([1, 0, 3]).unwrap() }, 11);
});
source

pub unsafe fn get_mut(&self, index: impl NpyIndex<Dim = D>) -> Option<&mut T>

Same as get, but returns Option<&mut T>.

§Safety

Calling this method is undefined behaviour if the underlying array is aliased immutably or mutably by other instances of PyArray or concurrently modified by Python or other native code.

Consider using safe alternatives like PyReadwriteArray::get_mut.

§Example
use numpy::{PyArray, PyArrayMethods};
use pyo3::Python;

Python::with_gil(|py| {
    let pyarray = PyArray::arange_bound(py, 0, 16, 1).reshape([2, 2, 4]).unwrap();

    unsafe {
        *pyarray.get_mut([1, 0, 3]).unwrap() = 42;
    }

    assert_eq!(unsafe { *pyarray.get([1, 0, 3]).unwrap() }, 42);
});
source

pub unsafe fn uget<Idx>(&self, index: Idx) -> &T
where Idx: NpyIndex<Dim = D>,

Get an immutable reference of the specified element, without checking the given index.

See NpyIndex for what types can be used as the index.

§Safety

Passing an invalid index is undefined behavior. The element must also have been initialized and all other references to it is must also be shared.

See PyReadonlyArray::get for a safe alternative.

§Example
use numpy::{PyArray, PyArrayMethods};
use pyo3::Python;

Python::with_gil(|py| {
    let pyarray = PyArray::arange_bound(py, 0, 16, 1).reshape([2, 2, 4]).unwrap();

    assert_eq!(unsafe { *pyarray.uget([1, 0, 3]) }, 11);
});
source

pub unsafe fn uget_mut<Idx>(&self, index: Idx) -> &mut T
where Idx: NpyIndex<Dim = D>,

Same as uget, but returns &mut T.

§Safety

Passing an invalid index is undefined behavior. The element must also have been initialized and other references to it must not exist.

See PyReadwriteArray::get_mut for a safe alternative.

source

pub unsafe fn uget_raw<Idx>(&self, index: Idx) -> *mut T
where Idx: NpyIndex<Dim = D>,

Same as uget, but returns *mut T.

§Safety

Passing an invalid index is undefined behavior.

source

pub fn get_owned<Idx>(&self, index: Idx) -> Option<T>
where Idx: NpyIndex<Dim = D>,

Get a copy of the specified element in the array.

See NpyIndex for what types can be used as the index.

§Example
use numpy::{PyArray, PyArrayMethods};
use pyo3::Python;

Python::with_gil(|py| {
    let pyarray = PyArray::arange_bound(py, 0, 16, 1).reshape([2, 2, 4]).unwrap();

    assert_eq!(pyarray.get_owned([1, 0, 3]), Some(11));
});
source

pub fn to_dyn(&self) -> &PyArray<T, IxDyn>

Turn an array with fixed dimensionality into one with dynamic dimensionality.

source

pub fn to_vec(&self) -> Result<Vec<T>, NotContiguousError>

Returns a copy of the internal data of the array as a Vec.

Fails if the internal array is not contiguous. See also as_slice.

§Example
use numpy::PyArray2;
use pyo3::Python;

Python::with_gil(|py| {
    let pyarray= py
        .eval("__import__('numpy').array([[0, 1], [2, 3]], dtype='int64')", None, None)
        .unwrap()
        .downcast::<PyArray2<i64>>()
        .unwrap();

    assert_eq!(pyarray.to_vec().unwrap(), vec![0, 1, 2, 3]);
});
source

pub fn from_array<'py, S>(py: Python<'py>, arr: &ArrayBase<S, D>) -> &'py Self
where S: Data<Elem = T>,

👎Deprecated since 0.21.0: will be replaced by PyArray::from_array_bound in the future
source

pub fn from_array_bound<'py, S>( py: Python<'py>, arr: &ArrayBase<S, D> ) -> Bound<'py, Self>
where S: Data<Elem = T>,

Construct a NumPy array from a ndarray::ArrayBase.

This method allocates memory in Python’s heap via the NumPy API, and then copies all elements of the array there.

§Example
use numpy::{PyArray, PyArrayMethods};
use ndarray::array;
use pyo3::Python;

Python::with_gil(|py| {
    let pyarray = PyArray::from_array_bound(py, &array![[1, 2], [3, 4]]);

    assert_eq!(pyarray.readonly().as_array(), array![[1, 2], [3, 4]]);
});
source

pub fn try_readonly(&self) -> Result<PyReadonlyArray<'_, T, D>, BorrowError>

Get an immutable borrow of the NumPy array

source

pub fn readonly(&self) -> PyReadonlyArray<'_, T, D>

Get an immutable borrow of the NumPy array

§Panics

Panics if the allocation backing the array is currently mutably borrowed.

For a non-panicking variant, use try_readonly.

source

pub fn try_readwrite(&self) -> Result<PyReadwriteArray<'_, T, D>, BorrowError>

Get a mutable borrow of the NumPy array

source

pub fn readwrite(&self) -> PyReadwriteArray<'_, T, D>

Get a mutable borrow of the NumPy array

§Panics

Panics if the allocation backing the array is currently borrowed or if the array is flagged as not writeable.

For a non-panicking variant, use try_readwrite.

source

pub unsafe fn as_array(&self) -> ArrayView<'_, T, D>

Returns an ArrayView of the internal array.

See also PyReadonlyArray::as_array.

§Safety

Calling this method invalidates all exclusive references to the internal data, e.g. &mut [T] or ArrayViewMut.

source

pub unsafe fn as_array_mut(&self) -> ArrayViewMut<'_, T, D>

Returns an ArrayViewMut of the internal array.

See also PyReadwriteArray::as_array_mut.

§Safety

Calling this method invalidates all other references to the internal data, e.g. ArrayView or ArrayViewMut.

source

pub fn as_raw_array(&self) -> RawArrayView<T, D>

Returns the internal array as RawArrayView enabling element access via raw pointers

source

pub fn as_raw_array_mut(&self) -> RawArrayViewMut<T, D>

Returns the internal array as RawArrayViewMut enabling element access via raw pointers

source

pub fn to_owned_array(&self) -> Array<T, D>

Get a copy of the array as an ndarray::Array.

§Example
use numpy::{PyArray, PyArrayMethods};
use ndarray::array;
use pyo3::Python;

Python::with_gil(|py| {
    let pyarray = PyArray::arange_bound(py, 0, 4, 1).reshape([2, 2]).unwrap();

    assert_eq!(
        pyarray.to_owned_array(),
        array![[0, 1], [2, 3]]
    )
});
source§

impl<N, D> PyArray<N, D>
where N: Scalar + Element, D: Dimension,

source

pub unsafe fn try_as_matrix<R, C, RStride, CStride>( &self ) -> Option<MatrixView<'_, N, R, C, RStride, CStride>>
where R: Dim, C: Dim, RStride: Dim, CStride: Dim,

Try to convert this array into a nalgebra::MatrixView using the given shape and strides.

See PyReadonlyArray::try_as_matrix for a discussion of the memory layout requirements.

§Safety

Calling this method invalidates all exclusive references to the internal data, e.g. ArrayViewMut or MatrixSliceMut.

source

pub unsafe fn try_as_matrix_mut<R, C, RStride, CStride>( &self ) -> Option<MatrixViewMut<'_, N, R, C, RStride, CStride>>
where R: Dim, C: Dim, RStride: Dim, CStride: Dim,

Try to convert this array into a nalgebra::MatrixViewMut using the given shape and strides.

See PyReadonlyArray::try_as_matrix for a discussion of the memory layout requirements.

§Safety

Calling this method invalidates all other references to the internal data, e.g. ArrayView, MatrixSlice, ArrayViewMut or MatrixSliceMut.

source§

impl<D: Dimension> PyArray<PyObject, D>

source

pub fn from_owned_object_array<'py, T>( py: Python<'py>, arr: Array<Py<T>, D> ) -> &'py Self

👎Deprecated since 0.21.0: will be replaced by PyArray::from_owned_object_array_bound in the future
source

pub fn from_owned_object_array_bound<T>( py: Python<'_>, arr: Array<Py<T>, D> ) -> Bound<'_, Self>

Construct a NumPy array containing objects stored in a ndarray::Array

This method uses the internal Vec of the ndarray::Array as the base object of the NumPy array.

§Example
use ndarray::array;
use pyo3::{pyclass, Py, Python};
use numpy::{PyArray, PyArrayMethods};

#[pyclass]
struct CustomElement {
    foo: i32,
    bar: f64,
}

Python::with_gil(|py| {
    let array = array![
        Py::new(py, CustomElement {
            foo: 1,
            bar: 2.0,
        }).unwrap(),
        Py::new(py, CustomElement {
            foo: 3,
            bar: 4.0,
        }).unwrap(),
    ];

    let pyarray = PyArray::from_owned_object_array_bound(py, array);

    assert!(pyarray.readonly().as_array().get(0).unwrap().as_ref(py).is_instance_of::<CustomElement>());
});
source§

impl<T: Copy + Element> PyArray<T, Ix0>

source

pub fn item(&self) -> T

Get the single element of a zero-dimensional array.

See inner for an example.

source§

impl<T: Element> PyArray<T, Ix1>

source

pub fn from_slice<'py>(py: Python<'py>, slice: &[T]) -> &'py Self

👎Deprecated since 0.21.0: will be replaced by PyArray::from_slice_bound in the future
source

pub fn from_slice_bound<'py>(py: Python<'py>, slice: &[T]) -> Bound<'py, Self>

Construct a one-dimensional array from a slice.

§Example
use numpy::{PyArray, PyArrayMethods};
use pyo3::Python;

Python::with_gil(|py| {
    let slice = &[1, 2, 3, 4, 5];
    let pyarray = PyArray::from_slice_bound(py, slice);
    assert_eq!(pyarray.readonly().as_slice().unwrap(), &[1, 2, 3, 4, 5]);
});
source

pub fn from_vec<'py>(py: Python<'py>, vec: Vec<T>) -> &'py Self

👎Deprecated since 0.21.0: will be replaced by PyArray::from_vec_bound in the future
source

pub fn from_vec_bound<'py>(py: Python<'py>, vec: Vec<T>) -> Bound<'py, Self>

Construct a one-dimensional array from a Vec<T>.

§Example
use numpy::{PyArray, PyArrayMethods};
use pyo3::Python;

Python::with_gil(|py| {
    let vec = vec![1, 2, 3, 4, 5];
    let pyarray = PyArray::from_vec_bound(py, vec);
    assert_eq!(pyarray.readonly().as_slice().unwrap(), &[1, 2, 3, 4, 5]);
});
source

pub fn from_iter<'py, I>(py: Python<'py>, iter: I) -> &'py Self
where I: IntoIterator<Item = T>,

👎Deprecated since 0.21.0: will be replaced by PyArray::from_iter_bound in the future
source

pub fn from_iter_bound<I>(py: Python<'_>, iter: I) -> Bound<'_, Self>
where I: IntoIterator<Item = T>,

Construct a one-dimensional array from an Iterator.

If no reliable size_hint is available, this method can allocate memory multiple times, which can hurt performance.

§Example
use numpy::{PyArray, PyArrayMethods};
use pyo3::Python;

Python::with_gil(|py| {
    let pyarray = PyArray::from_iter_bound(py, "abcde".chars().map(u32::from));
    assert_eq!(pyarray.readonly().as_slice().unwrap(), &[97, 98, 99, 100, 101]);
});
source§

impl<T: Element> PyArray<T, Ix2>

source

pub fn from_vec2<'py>( py: Python<'py>, v: &[Vec<T>] ) -> Result<&'py Self, FromVecError>

👎Deprecated since 0.21.0: will be replaced by PyArray::from_vec2_bound in the future
source

pub fn from_vec2_bound<'py>( py: Python<'py>, v: &[Vec<T>] ) -> Result<Bound<'py, Self>, FromVecError>

Construct a two-dimension array from a Vec<Vec<T>>.

This function checks all dimensions of the inner vectors and returns an error if they are not all equal.

§Example
use numpy::{PyArray, PyArrayMethods};
use pyo3::Python;
use ndarray::array;

Python::with_gil(|py| {
    let vec2 = vec![vec![11, 12], vec![21, 22]];
    let pyarray = PyArray::from_vec2_bound(py, &vec2).unwrap();
    assert_eq!(pyarray.readonly().as_array(), array![[11, 12], [21, 22]]);

    let ragged_vec2 = vec![vec![11, 12], vec![21]];
    assert!(PyArray::from_vec2_bound(py, &ragged_vec2).is_err());
});
source§

impl<T: Element> PyArray<T, Ix3>

source

pub fn from_vec3<'py>( py: Python<'py>, v: &[Vec<Vec<T>>] ) -> Result<&'py Self, FromVecError>

👎Deprecated since 0.21.0: will be replaced by PyArray::from_vec3_bound in the future
source

pub fn from_vec3_bound<'py>( py: Python<'py>, v: &[Vec<Vec<T>>] ) -> Result<Bound<'py, Self>, FromVecError>

Construct a three-dimensional array from a Vec<Vec<Vec<T>>>.

This function checks all dimensions of the inner vectors and returns an error if they are not all equal.

§Example
use numpy::{PyArray, PyArrayMethods};
use pyo3::Python;
use ndarray::array;

Python::with_gil(|py| {
    let vec3 = vec![
        vec![vec![111, 112], vec![121, 122]],
        vec![vec![211, 212], vec![221, 222]],
    ];
    let pyarray = PyArray::from_vec3_bound(py, &vec3).unwrap();
    assert_eq!(
        pyarray.readonly().as_array(),
        array![[[111, 112], [121, 122]], [[211, 212], [221, 222]]]
    );

    let ragged_vec3 = vec![
        vec![vec![111, 112], vec![121, 122]],
        vec![vec![211], vec![221, 222]],
    ];
    assert!(PyArray::from_vec3_bound(py, &ragged_vec3).is_err());
});
source§

impl<T: Element, D> PyArray<T, D>

source

pub fn copy_to<U: Element>(&self, other: &PyArray<U, D>) -> PyResult<()>

Copies self into other, performing a data type conversion if necessary.

See also PyArray_CopyInto.

§Example
use numpy::{PyArray, PyArrayMethods};
use pyo3::Python;

Python::with_gil(|py| {
    let pyarray_f = PyArray::arange_bound(py, 2.0, 5.0, 1.0);
    let pyarray_i = unsafe { PyArray::<i64, _>::new_bound(py, [3], false) };

    assert!(pyarray_f.copy_to(&pyarray_i).is_ok());

    assert_eq!(pyarray_i.readonly().as_slice().unwrap(), &[2, 3, 4]);
});
source

pub fn cast<'py, U: Element>( &'py self, is_fortran: bool ) -> PyResult<&'py PyArray<U, D>>

Cast the PyArray<T> to PyArray<U>, by allocating a new array.

See also PyArray_CastToType.

§Example
use numpy::{PyArray, PyArrayMethods};
use pyo3::Python;

Python::with_gil(|py| {
    let pyarray_f = PyArray::arange_bound(py, 2.0, 5.0, 1.0);

    let pyarray_i = pyarray_f.cast::<i32>(false).unwrap();

    assert_eq!(pyarray_i.readonly().as_slice().unwrap(), &[2, 3, 4]);
});
source

pub fn reshape_with_order<'py, ID: IntoDimension>( &'py self, dims: ID, order: NPY_ORDER ) -> PyResult<&'py PyArray<T, ID::Dim>>

Construct a new array which has same values as self, but has different dimensions specified by dims and a possibly different memory order specified by order.

See also numpy.reshape and PyArray_Newshape.

§Example
use numpy::prelude::*;
use numpy::{npyffi::NPY_ORDER, PyArray};
use pyo3::Python;
use ndarray::array;

Python::with_gil(|py| {
    let array =
        PyArray::from_iter_bound(py, 0..9).reshape_with_order([3, 3], NPY_ORDER::NPY_FORTRANORDER).unwrap();

    assert_eq!(array.readonly().as_array(), array![[0, 3, 6], [1, 4, 7], [2, 5, 8]]);
    assert!(array.is_fortran_contiguous());

    assert!(array.reshape([5]).is_err());
});
source

pub fn reshape<'py, ID: IntoDimension>( &'py self, dims: ID ) -> PyResult<&'py PyArray<T, ID::Dim>>

Special case of reshape_with_order which keeps the memory order the same.

source

pub unsafe fn resize<ID: IntoDimension>(&self, dims: ID) -> PyResult<()>

Extends or truncates the dimensions of an array.

This method works only on contiguous arrays. Missing elements will be initialized as if calling zeros.

See also ndarray.resize and PyArray_Resize.

§Safety

There should be no outstanding references (shared or exclusive) into the array as this method might re-allocate it and thereby invalidate all pointers into it.

§Example
use numpy::prelude::*;
use numpy::PyArray;
use pyo3::Python;

Python::with_gil(|py| {
    let pyarray = PyArray::<f64, _>::zeros_bound(py, (10, 10), false);
    assert_eq!(pyarray.shape(), [10, 10]);

    unsafe {
        pyarray.resize((100, 100)).unwrap();
    }
    assert_eq!(pyarray.shape(), [100, 100]);
});
source§

impl<T: Element + AsPrimitive<f64>> PyArray<T, Ix1>

source

pub fn arange<'py>(py: Python<'py>, start: T, stop: T, step: T) -> &Self

👎Deprecated since 0.21.0: will be replaced by PyArray::arange_bound in the future

Deprecated form of PyArray<T, Ix1>::arange_bound

source

pub fn arange_bound<'py>( py: Python<'py>, start: T, stop: T, step: T ) -> Bound<'py, Self>

Return evenly spaced values within a given interval.

See numpy.arange for the Python API and PyArray_Arange for the C API.

§Example
use numpy::{PyArray, PyArrayMethods};
use pyo3::Python;

Python::with_gil(|py| {
    let pyarray = PyArray::arange_bound(py, 2.0, 4.0, 0.5);
    assert_eq!(pyarray.readonly().as_slice().unwrap(), &[2.0, 2.5, 3.0, 3.5]);

    let pyarray = PyArray::arange_bound(py, -2, 4, 3);
    assert_eq!(pyarray.readonly().as_slice().unwrap(), &[-2, 1]);
});

Methods from Deref<Target = PyUntypedArray>§

source

pub fn as_array_ptr(&self) -> *mut PyArrayObject

Returns a raw pointer to the underlying PyArrayObject.

source

pub fn dtype(&self) -> &PyArrayDescr

Returns the dtype of the array.

See also ndarray.dtype and PyArray_DTYPE.

§Example
use numpy::prelude::*;
use numpy::{dtype_bound, PyArray};
use pyo3::Python;

Python::with_gil(|py| {
   let array = PyArray::from_vec_bound(py, vec![1_i32, 2, 3]);

   assert!(array.dtype().is_equiv_to(&dtype_bound::<i32>(py)));
});
source

pub fn is_contiguous(&self) -> bool

Returns true if the internal data of the array is contiguous, indepedently of whether C-style/row-major or Fortran-style/column-major.

§Example
use numpy::{PyArray1, PyUntypedArrayMethods};
use pyo3::{types::{IntoPyDict, PyAnyMethods}, Python};

Python::with_gil(|py| {
    let array = PyArray1::arange_bound(py, 0, 10, 1);
    assert!(array.is_contiguous());

    let view = py
        .eval_bound("array[::2]", None, Some(&[("array", array)].into_py_dict_bound(py)))
        .unwrap()
        .downcast_into::<PyArray1<i32>>()
        .unwrap();
    assert!(!view.is_contiguous());
});
source

pub fn is_fortran_contiguous(&self) -> bool

Returns true if the internal data of the array is Fortran-style/column-major contiguous.

source

pub fn is_c_contiguous(&self) -> bool

Returns true if the internal data of the array is C-style/row-major contiguous.

source

pub fn ndim(&self) -> usize

Returns the number of dimensions of the array.

See also ndarray.ndim and PyArray_NDIM.

§Example
use numpy::{PyArray3, PyUntypedArrayMethods};
use pyo3::Python;

Python::with_gil(|py| {
    let arr = PyArray3::<f64>::zeros_bound(py, [4, 5, 6], false);

    assert_eq!(arr.ndim(), 3);
});
source

pub fn strides(&self) -> &[isize]

Returns a slice indicating how many bytes to advance when iterating along each axis.

See also ndarray.strides and PyArray_STRIDES.

§Example
use numpy::{PyArray3, PyUntypedArrayMethods};
use pyo3::Python;

Python::with_gil(|py| {
    let arr = PyArray3::<f64>::zeros_bound(py, [4, 5, 6], false);

    assert_eq!(arr.strides(), &[240, 48, 8]);
});
source

pub fn shape(&self) -> &[usize]

Returns a slice which contains dimmensions of the array.

See also [ndarray.shape][ndaray-shape] and PyArray_DIMS.

§Example
use numpy::{PyArray3, PyUntypedArrayMethods};
use pyo3::Python;

Python::with_gil(|py| {
    let arr = PyArray3::<f64>::zeros_bound(py, [4, 5, 6], false);

    assert_eq!(arr.shape(), &[4, 5, 6]);
});
source

pub fn len(&self) -> usize

Calculates the total number of elements in the array.

source

pub fn is_empty(&self) -> bool

Returns true if the there are no elements in the array.

Methods from Deref<Target = PyAny>§

pub fn is<T>(&self, other: &T) -> bool
where T: AsPyPointer,

Returns whether self and other point to the same object. To compare the equality of two objects (the == operator), use eq.

This is equivalent to the Python expression self is other.

pub fn hasattr<N>(&self, attr_name: N) -> Result<bool, PyErr>
where N: IntoPy<Py<PyString>>,

Determines whether this object has the given attribute.

This is equivalent to the Python expression hasattr(self, attr_name).

To avoid repeated temporary allocations of Python strings, the [intern!] macro can be used to intern attr_name.

§Example: intern!ing the attribute name
#[pyfunction]
fn has_version(sys: &Bound<'_, PyModule>) -> PyResult<bool> {
    sys.hasattr(intern!(sys.py(), "version"))
}

pub fn getattr<N>(&self, attr_name: N) -> Result<&PyAny, PyErr>
where N: IntoPy<Py<PyString>>,

Retrieves an attribute value.

This is equivalent to the Python expression self.attr_name.

To avoid repeated temporary allocations of Python strings, the [intern!] macro can be used to intern attr_name.

§Example: intern!ing the attribute name
#[pyfunction]
fn version<'py>(sys: &Bound<'py, PyModule>) -> PyResult<Bound<'py, PyAny>> {
    sys.getattr(intern!(sys.py(), "version"))
}

pub fn setattr<N, V>(&self, attr_name: N, value: V) -> Result<(), PyErr>
where N: IntoPy<Py<PyString>>, V: ToPyObject,

Sets an attribute value.

This is equivalent to the Python expression self.attr_name = value.

To avoid repeated temporary allocations of Python strings, the [intern!] macro can be used to intern name.

§Example: intern!ing the attribute name
#[pyfunction]
fn set_answer(ob: &Bound<'_, PyAny>) -> PyResult<()> {
    ob.setattr(intern!(ob.py(), "answer"), 42)
}

pub fn delattr<N>(&self, attr_name: N) -> Result<(), PyErr>
where N: IntoPy<Py<PyString>>,

Deletes an attribute.

This is equivalent to the Python statement del self.attr_name.

To avoid repeated temporary allocations of Python strings, the [intern!] macro can be used to intern attr_name.

pub fn compare<O>(&self, other: O) -> Result<Ordering, PyErr>
where O: ToPyObject,

Returns an Ordering between self and other.

This is equivalent to the following Python code:

if self == other:
    return Equal
elif a < b:
    return Less
elif a > b:
    return Greater
else:
    raise TypeError("PyAny::compare(): All comparisons returned false")
§Examples
use pyo3::prelude::*;
use pyo3::types::PyFloat;
use std::cmp::Ordering;

Python::with_gil(|py| -> PyResult<()> {
    let a = PyFloat::new_bound(py, 0_f64);
    let b = PyFloat::new_bound(py, 42_f64);
    assert_eq!(a.compare(b)?, Ordering::Less);
    Ok(())
})?;

It will return PyErr for values that cannot be compared:

use pyo3::prelude::*;
use pyo3::types::{PyFloat, PyString};

Python::with_gil(|py| -> PyResult<()> {
    let a = PyFloat::new_bound(py, 0_f64);
    let b = PyString::new_bound(py, "zero");
    assert!(a.compare(b).is_err());
    Ok(())
})?;

pub fn rich_compare<O>( &self, other: O, compare_op: CompareOp ) -> Result<&PyAny, PyErr>
where O: ToPyObject,

Tests whether two Python objects obey a given [CompareOp].

lt, le, eq, ne, gt and ge are the specialized versions of this function.

Depending on the value of compare_op, this is equivalent to one of the following Python expressions:

compare_opPython expression
[CompareOp::Eq]self == other
[CompareOp::Ne]self != other
[CompareOp::Lt]self < other
[CompareOp::Le]self <= other
[CompareOp::Gt]self > other
[CompareOp::Ge]self >= other
§Examples
use pyo3::class::basic::CompareOp;
use pyo3::prelude::*;
use pyo3::types::PyInt;

Python::with_gil(|py| -> PyResult<()> {
    let a: Bound<'_, PyInt> = 0_u8.into_py(py).into_bound(py).downcast_into()?;
    let b: Bound<'_, PyInt> = 42_u8.into_py(py).into_bound(py).downcast_into()?;
    assert!(a.rich_compare(b, CompareOp::Le)?.is_truthy()?);
    Ok(())
})?;

pub fn lt<O>(&self, other: O) -> Result<bool, PyErr>
where O: ToPyObject,

Tests whether this object is less than another.

This is equivalent to the Python expression self < other.

pub fn le<O>(&self, other: O) -> Result<bool, PyErr>
where O: ToPyObject,

Tests whether this object is less than or equal to another.

This is equivalent to the Python expression self <= other.

pub fn eq<O>(&self, other: O) -> Result<bool, PyErr>
where O: ToPyObject,

Tests whether this object is equal to another.

This is equivalent to the Python expression self == other.

pub fn ne<O>(&self, other: O) -> Result<bool, PyErr>
where O: ToPyObject,

Tests whether this object is not equal to another.

This is equivalent to the Python expression self != other.

pub fn gt<O>(&self, other: O) -> Result<bool, PyErr>
where O: ToPyObject,

Tests whether this object is greater than another.

This is equivalent to the Python expression self > other.

pub fn ge<O>(&self, other: O) -> Result<bool, PyErr>
where O: ToPyObject,

Tests whether this object is greater than or equal to another.

This is equivalent to the Python expression self >= other.

pub fn is_callable(&self) -> bool

Determines whether this object appears callable.

This is equivalent to Python’s callable() function.

§Examples
use pyo3::prelude::*;

Python::with_gil(|py| -> PyResult<()> {
    let builtins = PyModule::import_bound(py, "builtins")?;
    let print = builtins.getattr("print")?;
    assert!(print.is_callable());
    Ok(())
})?;

This is equivalent to the Python statement assert callable(print).

Note that unless an API needs to distinguish between callable and non-callable objects, there is no point in checking for callability. Instead, it is better to just do the call and handle potential exceptions.

pub fn call( &self, args: impl IntoPy<Py<PyTuple>>, kwargs: Option<&PyDict> ) -> Result<&PyAny, PyErr>

Calls the object.

This is equivalent to the Python expression self(*args, **kwargs).

§Examples
use pyo3::prelude::*;
use pyo3::types::PyDict;

const CODE: &str = r#"
def function(*args, **kwargs):
    assert args == ("hello",)
    assert kwargs == {"cruel": "world"}
    return "called with args and kwargs"
"#;

Python::with_gil(|py| {
    let module = PyModule::from_code_bound(py, CODE, "", "")?;
    let fun = module.getattr("function")?;
    let args = ("hello",);
    let kwargs = PyDict::new_bound(py);
    kwargs.set_item("cruel", "world")?;
    let result = fun.call(args, Some(&kwargs))?;
    assert_eq!(result.extract::<String>()?, "called with args and kwargs");
    Ok(())
})

pub fn call0(&self) -> Result<&PyAny, PyErr>

Calls the object without arguments.

This is equivalent to the Python expression self().

§Examples
use pyo3::prelude::*;

Python::with_gil(|py| -> PyResult<()> {
    let module = PyModule::import_bound(py, "builtins")?;
    let help = module.getattr("help")?;
    help.call0()?;
    Ok(())
})?;

This is equivalent to the Python expression help().

pub fn call1(&self, args: impl IntoPy<Py<PyTuple>>) -> Result<&PyAny, PyErr>

Calls the object with only positional arguments.

This is equivalent to the Python expression self(*args).

§Examples
use pyo3::prelude::*;

const CODE: &str = r#"
def function(*args, **kwargs):
    assert args == ("hello",)
    assert kwargs == {}
    return "called with args"
"#;

Python::with_gil(|py| {
    let module = PyModule::from_code_bound(py, CODE, "", "")?;
    let fun = module.getattr("function")?;
    let args = ("hello",);
    let result = fun.call1(args)?;
    assert_eq!(result.extract::<String>()?, "called with args");
    Ok(())
})

pub fn call_method<N, A>( &self, name: N, args: A, kwargs: Option<&PyDict> ) -> Result<&PyAny, PyErr>
where N: IntoPy<Py<PyString>>, A: IntoPy<Py<PyTuple>>,

Calls a method on the object.

This is equivalent to the Python expression self.name(*args, **kwargs).

To avoid repeated temporary allocations of Python strings, the [intern!] macro can be used to intern name.

§Examples
use pyo3::prelude::*;
use pyo3::types::PyDict;

const CODE: &str = r#"
class A:
    def method(self, *args, **kwargs):
        assert args == ("hello",)
        assert kwargs == {"cruel": "world"}
        return "called with args and kwargs"
a = A()
"#;

Python::with_gil(|py| {
    let module = PyModule::from_code_bound(py, CODE, "", "")?;
    let instance = module.getattr("a")?;
    let args = ("hello",);
    let kwargs = PyDict::new_bound(py);
    kwargs.set_item("cruel", "world")?;
    let result = instance.call_method("method", args, Some(&kwargs))?;
    assert_eq!(result.extract::<String>()?, "called with args and kwargs");
    Ok(())
})

pub fn call_method0<N>(&self, name: N) -> Result<&PyAny, PyErr>
where N: IntoPy<Py<PyString>>,

Calls a method on the object without arguments.

This is equivalent to the Python expression self.name().

To avoid repeated temporary allocations of Python strings, the [intern!] macro can be used to intern name.

§Examples
use pyo3::prelude::*;

const CODE: &str = r#"
class A:
    def method(self, *args, **kwargs):
        assert args == ()
        assert kwargs == {}
        return "called with no arguments"
a = A()
"#;

Python::with_gil(|py| {
    let module = PyModule::from_code_bound(py, CODE, "", "")?;
    let instance = module.getattr("a")?;
    let result = instance.call_method0("method")?;
    assert_eq!(result.extract::<String>()?, "called with no arguments");
    Ok(())
})

pub fn call_method1<N, A>(&self, name: N, args: A) -> Result<&PyAny, PyErr>
where N: IntoPy<Py<PyString>>, A: IntoPy<Py<PyTuple>>,

Calls a method on the object with only positional arguments.

This is equivalent to the Python expression self.name(*args).

To avoid repeated temporary allocations of Python strings, the [intern!] macro can be used to intern name.

§Examples
use pyo3::prelude::*;

const CODE: &str = r#"
class A:
    def method(self, *args, **kwargs):
        assert args == ("hello",)
        assert kwargs == {}
        return "called with args"
a = A()
"#;

Python::with_gil(|py| {
    let module = PyModule::from_code_bound(py, CODE, "", "")?;
    let instance = module.getattr("a")?;
    let args = ("hello",);
    let result = instance.call_method1("method", args)?;
    assert_eq!(result.extract::<String>()?, "called with args");
    Ok(())
})

pub fn is_true(&self) -> Result<bool, PyErr>

👎Deprecated since 0.21.0: use .is_truthy() instead

Returns whether the object is considered to be true.

This is equivalent to the Python expression bool(self).

pub fn is_truthy(&self) -> Result<bool, PyErr>

Returns whether the object is considered to be true.

This applies truth value testing equivalent to the Python expression bool(self).

pub fn is_none(&self) -> bool

Returns whether the object is considered to be None.

This is equivalent to the Python expression self is None.

pub fn is_ellipsis(&self) -> bool

👎Deprecated since 0.20.0: use .is(py.Ellipsis()) instead

Returns whether the object is Ellipsis, e.g. ....

This is equivalent to the Python expression self is ....

pub fn is_empty(&self) -> Result<bool, PyErr>

Returns true if the sequence or mapping has a length of 0.

This is equivalent to the Python expression len(self) == 0.

pub fn get_item<K>(&self, key: K) -> Result<&PyAny, PyErr>
where K: ToPyObject,

Gets an item from the collection.

This is equivalent to the Python expression self[key].

pub fn set_item<K, V>(&self, key: K, value: V) -> Result<(), PyErr>
where K: ToPyObject, V: ToPyObject,

Sets a collection item value.

This is equivalent to the Python expression self[key] = value.

pub fn del_item<K>(&self, key: K) -> Result<(), PyErr>
where K: ToPyObject,

Deletes an item from the collection.

This is equivalent to the Python expression del self[key].

pub fn iter(&self) -> Result<&PyIterator, PyErr>

Takes an object and returns an iterator for it.

This is typically a new iterator but if the argument is an iterator, this returns itself.

pub fn get_type(&self) -> &PyType

Returns the Python type object for this object’s type.

pub fn get_type_ptr(&self) -> *mut PyTypeObject

Returns the Python type pointer for this object.

pub fn downcast<T>(&self) -> Result<&T, PyDowncastError<'_>>
where T: PyTypeCheck<AsRefTarget = T>,

Downcast this PyAny to a concrete Python type or pyclass.

Note that you can often avoid downcasting yourself by just specifying the desired type in function or method signatures. However, manual downcasting is sometimes necessary.

For extracting a Rust-only type, see PyAny::extract.

§Example: Downcasting to a specific Python object
use pyo3::prelude::*;
use pyo3::types::{PyDict, PyList};

Python::with_gil(|py| {
    let dict = PyDict::new_bound(py);
    assert!(dict.is_instance_of::<PyAny>());
    let any = dict.as_any();

    assert!(any.downcast::<PyDict>().is_ok());
    assert!(any.downcast::<PyList>().is_err());
});
§Example: Getting a reference to a pyclass

This is useful if you want to mutate a PyObject that might actually be a pyclass.

use pyo3::prelude::*;

#[pyclass]
struct Class {
    i: i32,
}

Python::with_gil(|py| {
    let class = Py::new(py, Class { i: 0 }).unwrap().into_bound(py).into_any();

    let class_bound: &Bound<'_, Class> = class.downcast()?;

    class_bound.borrow_mut().i += 1;

    // Alternatively you can get a `PyRefMut` directly
    let class_ref: PyRefMut<'_, Class> = class.extract()?;
    assert_eq!(class_ref.i, 1);
    Ok(())
})

pub fn downcast_exact<T>(&self) -> Result<&T, PyDowncastError<'_>>
where T: PyTypeInfo<AsRefTarget = T>,

Downcast this PyAny to a concrete Python type or pyclass (but not a subclass of it).

It is almost always better to use [PyAny::downcast] because it accounts for Python subtyping. Use this method only when you do not want to allow subtypes.

The advantage of this method over [PyAny::downcast] is that it is faster. The implementation of downcast_exact uses the equivalent of the Python expression type(self) is T, whereas downcast uses isinstance(self, T).

For extracting a Rust-only type, see PyAny::extract.

§Example: Downcasting to a specific Python object but not a subtype
use pyo3::prelude::*;
use pyo3::types::{PyBool, PyLong};

Python::with_gil(|py| {
    let b = PyBool::new_bound(py, true);
    assert!(b.is_instance_of::<PyBool>());
    let any: &Bound<'_, PyAny> = b.as_any();

    // `bool` is a subtype of `int`, so `downcast` will accept a `bool` as an `int`
    // but `downcast_exact` will not.
    assert!(any.downcast::<PyLong>().is_ok());
    assert!(any.downcast_exact::<PyLong>().is_err());

    assert!(any.downcast_exact::<PyBool>().is_ok());
});

pub unsafe fn downcast_unchecked<T>(&self) -> &T
where T: HasPyGilRef<AsRefTarget = T>,

Converts this PyAny to a concrete Python type without checking validity.

§Safety

Callers must ensure that the type is valid or risk type confusion.

pub fn extract<'py, D>(&'py self) -> Result<D, PyErr>
where D: FromPyObjectBound<'py, 'py>,

Extracts some type from the Python object.

This is a wrapper function around FromPyObject::extract().

pub fn get_refcnt(&self) -> isize

Returns the reference count for the Python object.

pub fn repr(&self) -> Result<&PyString, PyErr>

Computes the “repr” representation of self.

This is equivalent to the Python expression repr(self).

pub fn str(&self) -> Result<&PyString, PyErr>

Computes the “str” representation of self.

This is equivalent to the Python expression str(self).

pub fn hash(&self) -> Result<isize, PyErr>

Retrieves the hash code of self.

This is equivalent to the Python expression hash(self).

pub fn len(&self) -> Result<usize, PyErr>

Returns the length of the sequence or mapping.

This is equivalent to the Python expression len(self).

pub fn dir(&self) -> &PyList

Returns the list of attributes of this object.

This is equivalent to the Python expression dir(self).

pub fn is_instance(&self, ty: &PyAny) -> Result<bool, PyErr>

Checks whether this object is an instance of type ty.

This is equivalent to the Python expression isinstance(self, ty).

pub fn is_exact_instance(&self, ty: &PyAny) -> bool

Checks whether this object is an instance of exactly type ty (not a subclass).

This is equivalent to the Python expression type(self) is ty.

pub fn is_instance_of<T>(&self) -> bool
where T: PyTypeInfo,

Checks whether this object is an instance of type T.

This is equivalent to the Python expression isinstance(self, T), if the type T is known at compile time.

pub fn is_exact_instance_of<T>(&self) -> bool
where T: PyTypeInfo,

Checks whether this object is an instance of exactly type T.

This is equivalent to the Python expression type(self) is T, if the type T is known at compile time.

pub fn contains<V>(&self, value: V) -> Result<bool, PyErr>
where V: ToPyObject,

Determines if self contains value.

This is equivalent to the Python expression value in self.

pub fn py(&self) -> Python<'_>

Returns a GIL marker constrained to the lifetime of this type.

pub fn as_ptr(&self) -> *mut PyObject

Returns the raw FFI pointer represented by self.

§Safety

Callers are responsible for ensuring that the pointer does not outlive self.

The reference is borrowed; callers should not decrease the reference count when they are finished with the pointer.

pub fn into_ptr(&self) -> *mut PyObject

Returns an owned raw FFI pointer represented by self.

§Safety

The reference is owned; when finished the caller should either transfer ownership of the pointer or decrease the reference count (e.g. with pyo3::ffi::Py_DecRef).

pub fn py_super(&self) -> Result<&PySuper, PyErr>

Return a proxy object that delegates method calls to a parent or sibling class of type.

This is equivalent to the Python expression super()

Trait Implementations§

source§

impl<T, D> AsPyPointer for PyArray<T, D>

source§

fn as_ptr(&self) -> *mut PyObject

Returns the underlying FFI pointer as a borrowed pointer.
source§

impl<T, D> AsRef<PyAny> for PyArray<T, D>

source§

fn as_ref(&self) -> &PyAny

Converts this type into a shared reference of the (usually inferred) input type.
source§

impl<T, D> Debug for PyArray<T, D>

source§

fn fmt(&self, f: &mut Formatter<'_>) -> Result<(), Error>

Formats the value using the given formatter. Read more
source§

impl<T, D> Deref for PyArray<T, D>

§

type Target = PyUntypedArray

The resulting type after dereferencing.
source§

fn deref(&self) -> &Self::Target

Dereferences the value.
source§

impl<T, D> Display for PyArray<T, D>

source§

fn fmt(&self, f: &mut Formatter<'_>) -> Result<(), Error>

Formats the value using the given formatter. Read more
source§

impl<'a, T, D> From<&'a PyArray<T, D>> for &'a PyAny

source§

fn from(ob: &'a PyArray<T, D>) -> Self

Converts to this type from the input type.
source§

impl<T, D> From<&PyArray<T, D>> for Py<PyArray<T, D>>

source§

fn from(other: &PyArray<T, D>) -> Self

Converts to this type from the input type.
source§

impl<'py, T: Element, D: Dimension> FromPyObject<'py> for &'py PyArray<T, D>

source§

fn extract_bound(ob: &Bound<'py, PyAny>) -> PyResult<Self>

Extracts Self from the bound smart pointer obj. Read more
§

fn extract(ob: &'py PyAny) -> Result<Self, PyErr>

Extracts Self from the source GIL Ref obj. Read more
source§

impl<T, D> IntoPy<Py<PyAny>> for PyArray<T, D>

source§

fn into_py<'py>(self, py: Python<'py>) -> PyObject

Performs the conversion.
source§

impl<T, D> IntoPy<Py<PyArray<T, D>>> for &PyArray<T, D>

source§

fn into_py<'py>(self, py: Python<'py>) -> Py<PyArray<T, D>>

Performs the conversion.
source§

impl<T, D> PyNativeType for PyArray<T, D>

§

type AsRefSource = PyArray<T, D>

The form of this which is stored inside a Py<T> smart pointer.
§

fn as_borrowed(&self) -> Borrowed<'_, '_, Self::AsRefSource>

Cast &self to a Borrowed smart pointer. Read more
§

fn py(&self) -> Python<'_>

Returns a GIL marker constrained to the lifetime of this type.
§

unsafe fn unchecked_downcast(obj: &PyAny) -> &Self

Cast &PyAny to &Self without no type checking. Read more
source§

impl<T: Element, D: Dimension> PyTypeInfo for PyArray<T, D>

source§

const NAME: &'static str = "PyArray<T, D>"

Class name.
source§

const MODULE: Option<&'static str> = _

Module name, if any.
source§

fn type_object_raw<'py>(py: Python<'py>) -> *mut PyTypeObject

Returns the PyTypeObject instance for this type.
source§

fn is_type_of_bound(ob: &Bound<'_, PyAny>) -> bool

Checks if object is an instance of this type or a subclass of this type.
§

fn type_object(py: Python<'_>) -> &PyType

Returns the safe abstraction over the type object.
§

fn type_object_bound(py: Python<'_>) -> Bound<'_, PyType>

Returns the safe abstraction over the type object.
§

fn is_type_of(object: &PyAny) -> bool

Checks if object is an instance of this type or a subclass of this type.
§

fn is_exact_type_of(object: &PyAny) -> bool

Checks if object is an instance of this type.
§

fn is_exact_type_of_bound(object: &Bound<'_, PyAny>) -> bool

Checks if object is an instance of this type.
source§

impl<T, D> ToPyObject for PyArray<T, D>

source§

fn to_object(&self, py: Python<'_>) -> PyObject

Converts self into a Python object.
source§

impl<T, D> DerefToPyAny for PyArray<T, D>

Auto Trait Implementations§

§

impl<T, D> !RefUnwindSafe for PyArray<T, D>

§

impl<T, D> !Send for PyArray<T, D>

§

impl<T, D> !Sync for PyArray<T, D>

§

impl<T, D> Unpin for PyArray<T, D>
where D: Unpin, T: Unpin,

§

impl<T, D> UnwindSafe for PyArray<T, D>
where D: UnwindSafe, T: UnwindSafe,

Blanket Implementations§

source§

impl<T> Any for T
where T: 'static + ?Sized,

source§

fn type_id(&self) -> TypeId

Gets the TypeId of self. Read more
source§

impl<T> Borrow<T> for T
where T: ?Sized,

source§

fn borrow(&self) -> &T

Immutably borrows from an owned value. Read more
source§

impl<T> BorrowMut<T> for T
where T: ?Sized,

source§

fn borrow_mut(&mut self) -> &mut T

Mutably borrows from an owned value. Read more
source§

impl<T> From<T> for T

source§

fn from(t: T) -> T

Returns the argument unchanged.

§

impl<'p, T> FromPyPointer<'p> for T
where T: 'p + PyNativeType,

§

unsafe fn from_owned_ptr_or_opt( py: Python<'p>, ptr: *mut PyObject ) -> Option<&'p T>

👎Deprecated since 0.21.0
Convert from an arbitrary PyObject. Read more
§

unsafe fn from_borrowed_ptr_or_opt( _py: Python<'p>, ptr: *mut PyObject ) -> Option<&'p T>

👎Deprecated since 0.21.0
Convert from an arbitrary borrowed PyObject. Read more
§

unsafe fn from_owned_ptr_or_panic( py: Python<'p>, ptr: *mut PyObject ) -> &'p Self

👎Deprecated since 0.21.0
Convert from an arbitrary PyObject or panic. Read more
§

unsafe fn from_owned_ptr(py: Python<'p>, ptr: *mut PyObject) -> &'p Self

👎Deprecated since 0.21.0
Convert from an arbitrary PyObject or panic. Read more
§

unsafe fn from_owned_ptr_or_err( py: Python<'p>, ptr: *mut PyObject ) -> Result<&'p Self, PyErr>

👎Deprecated since 0.21.0
Convert from an arbitrary PyObject. Read more
§

unsafe fn from_borrowed_ptr_or_panic( py: Python<'p>, ptr: *mut PyObject ) -> &'p Self

👎Deprecated since 0.21.0
Convert from an arbitrary borrowed PyObject. Read more
§

unsafe fn from_borrowed_ptr(py: Python<'p>, ptr: *mut PyObject) -> &'p Self

👎Deprecated since 0.21.0
Convert from an arbitrary borrowed PyObject. Read more
§

unsafe fn from_borrowed_ptr_or_err( py: Python<'p>, ptr: *mut PyObject ) -> Result<&'p Self, PyErr>

👎Deprecated since 0.21.0
Convert from an arbitrary borrowed PyObject. Read more
§

impl<T> HasPyGilRef for T
where T: PyNativeType,

§

type AsRefTarget = T

Utility type to make Py::as_ref work.
source§

impl<T, U> Into<U> for T
where U: From<T>,

source§

fn into(self) -> U

Calls U::from(self).

That is, this conversion is whatever the implementation of From<T> for U chooses to do.

§

impl<'v, T> PyTryFrom<'v> for T
where T: PyTypeInfo<AsRefTarget = T> + PyNativeType,

§

fn try_from<V>(value: V) -> Result<&'v T, PyDowncastError<'v>>
where V: Into<&'v PyAny>,

👎Deprecated since 0.21.0: use value.downcast::<T>() instead of T::try_from(value)
Cast from a concrete Python object type to PyObject.
§

fn try_from_exact<V>(value: V) -> Result<&'v T, PyDowncastError<'v>>
where V: Into<&'v PyAny>,

👎Deprecated since 0.21.0: use value.downcast_exact::<T>() instead of T::try_from_exact(value)
Cast from a concrete Python object type to PyObject. With exact type check.
§

unsafe fn try_from_unchecked<V>(value: V) -> &'v T
where V: Into<&'v PyAny>,

👎Deprecated since 0.21.0: use value.downcast_unchecked::<T>() instead of T::try_from_unchecked(value)
Cast a PyAny to a specific type of PyObject. The caller must have already verified the reference is for this type. Read more
§

impl<T> PyTypeCheck for T
where T: PyTypeInfo,

§

const NAME: &'static str = <T as PyTypeInfo>::NAME

Name of self. This is used in error messages, for example.
§

fn type_check(object: &Bound<'_, PyAny>) -> bool

Checks if object is an instance of Self, which may include a subtype. Read more
source§

impl<T> Same for T

§

type Output = T

Should always be Self
§

impl<SS, SP> SupersetOf<SS> for SP
where SS: SubsetOf<SP>,

§

fn to_subset(&self) -> Option<SS>

The inverse inclusion map: attempts to construct self from the equivalent element of its superset. Read more
§

fn is_in_subset(&self) -> bool

Checks if self is actually part of its subset T (and can be converted to it).
§

fn to_subset_unchecked(&self) -> SS

Use with care! Same as self.to_subset but without any property checks. Always succeeds.
§

fn from_subset(element: &SS) -> SP

The inclusion map: converts self to the equivalent element of its superset.
source§

impl<T> ToString for T
where T: Display + ?Sized,

source§

default fn to_string(&self) -> String

Converts the given value to a String. Read more
source§

impl<T, U> TryFrom<U> for T
where U: Into<T>,

§

type Error = Infallible

The type returned in the event of a conversion error.
source§

fn try_from(value: U) -> Result<T, <T as TryFrom<U>>::Error>

Performs the conversion.
source§

impl<T, U> TryInto<U> for T
where U: TryFrom<T>,

§

type Error = <U as TryFrom<T>>::Error

The type returned in the event of a conversion error.
source§

fn try_into(self) -> Result<U, <U as TryFrom<T>>::Error>

Performs the conversion.