[][src]Struct numpy::array::PyArray

pub struct PyArray<T, D>(_, _, _);

A safe, static-typed interface for NumPy ndarray.

Memory location

These methods don't allocate memory and use Box<[T]> as a internal buffer.

Please take care that you cannot use some destructive methods like resize, for this kind of array.

These methods allocate memory in Python's private heap.

In both cases, PyArray is managed by Python GC. So you can neither retrieve it nor deallocate it manually.

Reference

Like new, all constractor methods of PyArray returns &PyArray.

This design follows pyo3's ownership concept.

Data type and Dimension

PyArray has 2 type parametes T and D. T represents its data type like f32, and D represents its dimension.

All data types you can use implements Element.

Dimensions are represented by ndarray's Dimension trait.

Typically, you can use Ix1, Ix2, .. for fixed size arrays, and use IxDyn for dynamic dimensioned arrays. They're re-exported from ndarray crate.

You can also use various type aliases we provide, like PyArray1 or PyArrayDyn.

To specify concrete dimension like 3×4×5, you can use types which implements ndarray's IntoDimension trait. Typically, you can use array(e.g. [3, 4, 5]) or tuple(e.g. (3, 4, 5)) as a dimension.

Example

use pyo3::{GILGuard, Python};
use numpy::PyArray;
use ndarray::Array;
let gil = Python::acquire_gil();
let pyarray = PyArray::arange(gil.python(), 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

impl<T, D> PyArray<T, D>[src]

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

Gets a raw PyArrayObject pointer.

pub fn readonly(&self) -> PyReadonlyArray<T, D>[src]

Returns a temporally unwriteable reference of the array.

pub fn is_contiguous(&self) -> bool[src]

Returns true if the internal data of the array is C-style contiguous (default of numpy and ndarray) or Fortran-style contiguous.

Example

use pyo3::types::IntoPyDict;
let gil = pyo3::Python::acquire_gil();
let py = gil.python();
let array = numpy::PyArray::arange(py, 0, 1, 10);
assert!(array.is_contiguous());
let locals = [("np", numpy::get_array_module(py).unwrap())].into_py_dict(py);
let not_contiguous: &numpy::PyArray1<f32> = py
    .eval("np.zeros((3, 5))[::2, 4]", Some(locals), None)
    .unwrap()
    .downcast()
    .unwrap();
assert!(!not_contiguous.is_contiguous());

pub fn is_fortran_contiguous(&self) -> bool[src]

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

pub fn is_c_contiguous(&self) -> bool[src]

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

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

Get Py<PyArray> from &PyArray, which is the owned wrapper of PyObject.

You can use this method when you have to avoid lifetime annotation to your function args or return types, like used with pyo3's pymethod.

Example

use pyo3::{GILGuard, Python, Py, AsPyRef};
use numpy::PyArray1;
fn return_py_array() -> Py<PyArray1<i32>> {
   let gil = Python::acquire_gil();
   let array = PyArray1::zeros(gil.python(), [5], false);
   array.to_owned()
}
let gil = Python::acquire_gil();
let array = return_py_array();
assert_eq!(array.as_ref(gil.python()).readonly().as_slice().unwrap(), &[0, 0, 0, 0, 0]);

pub unsafe fn from_owned_ptr(py: Python, ptr: *mut PyObject) -> &Self[src]

Constructs PyArray from raw python object without incrementing reference counts.

pub unsafe fn from_borrowed_ptr(py: Python, ptr: *mut PyObject) -> &Self[src]

Constructs PyArray from raw python object and increments reference counts.

pub fn ndim(&self) -> usize[src]

Returns the number of dimensions in the array.

Same as numpy.ndarray.ndim

Example

use numpy::PyArray3;
let gil = pyo3::Python::acquire_gil();
let arr = PyArray3::<f64>::new(gil.python(), [4, 5, 6], false);
assert_eq!(arr.ndim(), 3);

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

Returns a slice which contains how many bytes you need to jump to the next row.

Same as numpy.ndarray.strides

Example

use numpy::PyArray3;
let gil = pyo3::Python::acquire_gil();
let arr = PyArray3::<f64>::new(gil.python(), [4, 5, 6], false);
assert_eq!(arr.strides(), &[240, 48, 8]);

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

Returns a slice which contains dimmensions of the array.

Same as numpy.ndarray.shape

Example

use numpy::PyArray3;
let gil = pyo3::Python::acquire_gil();
let arr = PyArray3::<f64>::new(gil.python(), [4, 5, 6], false);
assert_eq!(arr.shape(), &[4, 5, 6]);

pub fn len(&self) -> usize[src]

Calcurates the total number of elements in the array.

pub fn is_empty(&self) -> bool[src]

impl<T: Element, D: Dimension> PyArray<T, D>[src]

pub fn dims(&self) -> D[src]

Same as shape, but returns D

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

Creates a new uninitialized PyArray in python heap.

If is_fortran == true, returns Fortran-order array. Else, returns C-order array.

Example

use numpy::PyArray3;
let gil = pyo3::Python::acquire_gil();
let pyarray = PyArray3::<i32>::new(gil.python(), [4, 5, 6], false);
assert_eq!(pyarray.shape(), &[4, 5, 6]);

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

Construct a new nd-dimensional array filled with 0.

If is_fortran is true, then a fortran order array is created, otherwise a C-order array is created.

See also PyArray_Zeros

Example

use numpy::PyArray2;
let gil = pyo3::Python::acquire_gil();
let pyarray: &PyArray2<usize> = PyArray2::zeros(gil.python(), [2, 2], false);
assert_eq!(pyarray.readonly().as_array(), array![[0, 0], [0, 0]]);

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

Returns the immutable view of the internal data of PyArray as slice. Please consider the use of PyReadonlyArray::as_slice instead of this.

Safety

If the internal array is not readonly and can be mutated from Python code, holding the slice might cause undefined behavior.

pub fn as_cell_slice(&self) -> Result<&[Cell<T>], NotContiguousError>[src]

Returns the view of the internal data of PyArray as &[Cell<T>].

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

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

Safety

If another reference to the internal data exists(e.g., &[T] or ArrayView), it might cause undefined behavior.

In such case, please consider the use of as_cell_slice,

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

Construct PyArray from ndarray::Array.

This method uses internal Vec of ndarray::Array as numpy array.

Example

use numpy::PyArray;
let gil = pyo3::Python::acquire_gil();
let pyarray = PyArray::from_owned_array(gil.python(), array![[1, 2], [3, 4]]);
assert_eq!(pyarray.readonly().as_array(), array![[1, 2], [3, 4]]);

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

Get an immutable reference of a specified element, with checking the passed index is valid.

See NpyIndex for what types you can use as index.

If you pass an invalid index to this function, it returns None.

Example

use numpy::PyArray;
let gil = pyo3::Python::acquire_gil();
let arr = PyArray::arange(gil.python(), 0, 16, 1).reshape([2, 2, 4]).unwrap();
assert_eq!(*arr.get([1, 0, 3]).unwrap(), 11);
assert!(arr.get([2, 0, 3]).is_none());

For fixed dimension arrays, passing an index with invalid dimension causes compile error.

This example deliberately fails to compile
use numpy::PyArray;
let gil = pyo3::Python::acquire_gil();
let arr = PyArray::arange(gil.python(), 0, 16, 1).reshape([2, 2, 4]).unwrap();
let a = arr.get([1, 2]); // Compile Error!

However, for dinamic arrays, we cannot raise a compile error and just returns None.

use numpy::PyArray;
let gil = pyo3::Python::acquire_gil();
let arr = PyArray::arange(gil.python(), 0, 16, 1).reshape([2, 2, 4]).unwrap();
let arr = arr.to_dyn();
assert!(arr.get([1, 2].as_ref()).is_none());

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

Get an immutable reference of a specified element, without checking the passed index is valid.

See NpyIndex for what types you can use as index.

Passing an invalid index can cause undefined behavior(mostly SIGSEGV).

Example

use numpy::PyArray;
let gil = pyo3::Python::acquire_gil();
let arr = PyArray::arange(gil.python(), 0, 16, 1).reshape([2, 2, 4]).unwrap();
assert_eq!(unsafe { *arr.uget([1, 0, 3]) }, 11);

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

Same as uget, but returns &mut T.

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

Get dynamic dimensioned array from fixed dimension array.

See get for usage.

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

Returns the copy of the internal data of PyArray to Vec.

Returns ErrorKind::NotContiguous if the internal array is not contiguous. See also as_slice

Example

use numpy::PyArray2;
use pyo3::types::IntoPyDict;
let gil = pyo3::Python::acquire_gil();
let py = gil.python();
let locals = [("np", numpy::get_array_module(py).unwrap())].into_py_dict(py);
let array: &PyArray2<i64> = py
    .eval("np.array([[0, 1], [2, 3]], dtype='int64')", Some(locals), None)
    .unwrap()
    .downcast()
    .unwrap();
assert_eq!(array.to_vec().unwrap(), vec![0, 1, 2, 3]);

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

Construct PyArray from ndarray::ArrayBase.

This method allocates memory in Python's heap via numpy api, and then copies all elements of the array there.

Example

use numpy::PyArray;
let gil = pyo3::Python::acquire_gil();
let pyarray = PyArray::from_array(gil.python(), &array![[1, 2], [3, 4]]);
assert_eq!(pyarray.readonly().as_array(), array![[1, 2], [3, 4]]);

pub unsafe fn as_array(&self) -> ArrayView<T, D>[src]

Get the immutable view of the internal data of PyArray, as ndarray::ArrayView.

Safety

If the internal array is not readonly and can be mutated from Python code, holding the ArrayView might cause undefined behavior.

pub unsafe fn as_array_mut(&self) -> ArrayViewMut<T, D>[src]

Returns the internal array as ArrayViewMut. See also as_array.

Safety

If another reference to the internal data exists(e.g., &[T] or ArrayView), it might cause undefined behavior.

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

Get a copy of PyArray as ndarray::Array.

Example

use numpy::PyArray;
let gil = pyo3::Python::acquire_gil();
let py_array = PyArray::arange(gil.python(), 0, 4, 1).reshape([2, 2]).unwrap();
assert_eq!(
    py_array.to_owned_array(),
    array![[0, 1], [2, 3]]
)

impl<T: Element> PyArray<T, Ix1>[src]

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

Construct one-dimension PyArray from slice.

Example

use numpy::PyArray;
let gil = pyo3::Python::acquire_gil();
let array = [1, 2, 3, 4, 5];
let pyarray = PyArray::from_slice(gil.python(), &array);
assert_eq!(pyarray.readonly().as_slice().unwrap(), &[1, 2, 3, 4, 5]);

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

Construct one-dimension PyArray from Vec.

Example

use numpy::PyArray;
let gil = pyo3::Python::acquire_gil();
let vec = vec![1, 2, 3, 4, 5];
let pyarray = PyArray::from_vec(gil.python(), vec);
assert_eq!(pyarray.readonly().as_slice().unwrap(), &[1, 2, 3, 4, 5]);

pub fn from_exact_iter(
    py: Python,
    iter: impl ExactSizeIterator<Item = T>
) -> &Self
[src]

Construct one-dimension PyArray from a type which implements ExactSizeIterator.

Example

use numpy::PyArray;
use std::collections::BTreeSet;
let gil = pyo3::Python::acquire_gil();
let vec = vec![1, 2, 3, 4, 5];
let pyarray = PyArray::from_iter(gil.python(), vec.iter().map(|&x| x));
assert_eq!(pyarray.readonly().as_slice().unwrap(), &[1, 2, 3, 4, 5]);

pub fn from_iter(py: Python, iter: impl IntoIterator<Item = T>) -> &Self[src]

Construct one-dimension PyArray from a type which implements IntoIterator.

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

Example

use numpy::PyArray;
use std::collections::BTreeSet;
let gil = pyo3::Python::acquire_gil();
let set: BTreeSet<u32> = [4, 3, 2, 5, 1].into_iter().cloned().collect();
let pyarray = PyArray::from_iter(gil.python(), set);
assert_eq!(pyarray.readonly().as_slice().unwrap(), &[1, 2, 3, 4, 5]);

pub fn resize(&self, new_elems: usize) -> PyResult<()>[src]

Extends or trancates the length of 1 dimension PyArray.

Example

use numpy::PyArray;
let gil = pyo3::Python::acquire_gil();
let pyarray = PyArray::arange(gil.python(), 0, 10, 1);
assert_eq!(pyarray.len(), 10);
pyarray.resize(100).unwrap();
assert_eq!(pyarray.len(), 100);

impl<T: Element> PyArray<T, Ix2>[src]

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

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

This function checks all dimension of inner vec, and if there's any vec where its dimension differs from others, it returns ArrayCastError.

Example

use numpy::PyArray;
let gil = pyo3::Python::acquire_gil();
let vec2 = vec![vec![1, 2, 3]; 2];
let pyarray = PyArray::from_vec2(gil.python(), &vec2).unwrap();
assert_eq!(pyarray.readonly().as_array(), array![[1, 2, 3], [1, 2, 3]]);
assert!(PyArray::from_vec2(gil.python(), &[vec![1], vec![2, 3]]).is_err());

impl<T: Element> PyArray<T, Ix3>[src]

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

Construct a three-dimension PyArray from Vec<Vec<Vec<T>>>.

This function checks all dimension of inner vec, and if there's any vec where its dimension differs from others, it returns error.

Example

use numpy::PyArray;
let gil = pyo3::Python::acquire_gil();
let vec2 = vec![vec![vec![1, 2]; 2]; 2];
let pyarray = PyArray::from_vec3(gil.python(), &vec2).unwrap();
assert_eq!(
    pyarray.readonly().as_array(),
    array![[[1, 2], [1, 2]], [[1, 2], [1, 2]]]
);
assert!(PyArray::from_vec3(gil.python(), &[vec![vec![1], vec![]]]).is_err());

impl<T: Element, D> PyArray<T, D>[src]

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

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

Example

use numpy::PyArray;
let gil = pyo3::Python::acquire_gil();
let pyarray_f = PyArray::arange(gil.python(), 2.0, 5.0, 1.0);
let pyarray_i = PyArray::<i64, _>::new(gil.python(), [3], false);
assert!(pyarray_f.copy_to(pyarray_i).is_ok());
assert_eq!(pyarray_i.readonly().as_slice().unwrap(), &[2, 3, 4]);

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

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

Example

use numpy::PyArray;
let gil = pyo3::Python::acquire_gil();
let pyarray_f = PyArray::arange(gil.python(), 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]);

pub fn reshape<'py, ID, D2>(
    &'py self,
    dims: ID
) -> PyResult<&'py PyArray<T, D2>> where
    ID: IntoDimension<Dim = D2>,
    D2: Dimension
[src]

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

Since a returned array can contain a same pointer as self, we highly recommend to drop an old array, if this method returns Ok.

Example

use numpy::PyArray;
let gil = pyo3::Python::acquire_gil();
let array = PyArray::from_exact_iter(gil.python(), 0..9);
let array = array.reshape([3, 3]).unwrap();
assert_eq!(array.readonly().as_array(), array![[0, 1, 2], [3, 4, 5], [6, 7, 8]]);
assert!(array.reshape([5]).is_err());

pub fn reshape_with_order<'py, ID, D2>(
    &'py self,
    dims: ID,
    order: NPY_ORDER
) -> PyResult<&'py PyArray<T, D2>> where
    ID: IntoDimension<Dim = D2>,
    D2: Dimension
[src]

Same as reshape, but you can change the order of returned matrix.

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

pub fn arange(py: Python, start: T, stop: T, step: T) -> &Self[src]

Return evenly spaced values within a given interval. Same as numpy.arange.

See also PyArray_Arange.

Example

use numpy::PyArray;
let gil = pyo3::Python::acquire_gil();
let pyarray = PyArray::arange(gil.python(), 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(gil.python(), -2, 4, 3);
assert_eq!(pyarray.readonly().as_slice().unwrap(), &[-2, 1]);

Trait Implementations

impl<T, D> AsPyPointer for PyArray<T, D>[src]

fn as_ptr(&self) -> *mut PyObject[src]

Gets the underlying FFI pointer, returns a borrowed pointer.

impl<T, D> AsRef<PyAny> for PyArray<T, D>[src]

impl<'py, T, D> AsRef<PyArray<T, D>> for PyReadonlyArray<'py, T, D>[src]

impl<T, D> Debug for PyArray<T, D>[src]

impl<T, D> Deref for PyArray<T, D>[src]

type Target = PyAny

The resulting type after dereferencing.

impl<T, D> Display for PyArray<T, D>[src]

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

impl<'py, T, D> From<&'py PyArray<T, D>> for PyReadonlyArray<'py, T, D>[src]

impl<'a, T: Element, D: Dimension> FromPyObject<'a> for &'a PyArray<T, D>[src]

impl<T, D> IntoPy<PyObject> for PyArray<T, D>[src]

impl<T, D> PartialEq<PyArray<T, D>> for PyArray<T, D>[src]

impl<T, D> PyLayout<PyArray<T, D>> for PyArrayObject[src]

impl<T, D> PyNativeType for PyArray<T, D>[src]

impl<T, D> PySizedLayout<PyArray<T, D>> for PyArrayObject[src]

impl<T, D> PyTypeInfo for PyArray<T, D>[src]

type Type = ()

Type of objects to store in PyObject struct

type BaseType = PyAny

Base class

type Layout = PyArrayObject

Layout

type BaseLayout = PyObject

Layout of Basetype.

type Initializer = PyNativeTypeInitializer<Self>

Initializer for layout

type AsRefTarget = Self

Utility type to make AsPyRef work

impl<T, D> ToPyObject for PyArray<T, D>[src]

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

impl<T> Any for T where
    T: 'static + ?Sized
[src]

impl<T> Borrow<T> for T where
    T: ?Sized
[src]

impl<T> BorrowMut<T> for T where
    T: ?Sized
[src]

impl<T> From<T> for T[src]

impl<T> FromPy<T> for T[src]

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

impl<T, U> Into<U> for T where
    U: From<T>, 
[src]

impl<T, U> IntoPy<U> for T where
    U: FromPy<T>, 
[src]

impl<'v, T> PyTryFrom<'v> for T where
    T: PyTypeInfo + PyNativeType
[src]

impl<T> PyTypeObject for T where
    T: PyTypeInfo
[src]

impl<T> ToBorrowedObject for T where
    T: ToPyObject
[src]

impl<T> ToString for T where
    T: Display + ?Sized
[src]

impl<T, U> TryFrom<U> for T where
    U: Into<T>, 
[src]

type Error = Infallible

The type returned in the event of a conversion error.

impl<T, U> TryInto<U> for T where
    U: TryFrom<T>, 
[src]

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

The type returned in the event of a conversion error.