numpy/
array.rs

1//! Safe interface for NumPy's [N-dimensional arrays][ndarray]
2//!
3//! [ndarray]: https://numpy.org/doc/stable/reference/arrays.ndarray.html
4
5use std::{
6    marker::PhantomData,
7    mem,
8    os::raw::{c_int, c_void},
9    ptr, slice,
10};
11
12use ndarray::{
13    Array, ArrayBase, ArrayView, ArrayViewMut, Axis, Data, Dim, Dimension, IntoDimension, Ix0, Ix1,
14    Ix2, Ix3, Ix4, Ix5, Ix6, IxDyn, RawArrayView, RawArrayViewMut, RawData, ShapeBuilder,
15    StrideShape,
16};
17use num_traits::AsPrimitive;
18use pyo3::{
19    ffi,
20    types::{DerefToPyAny, PyAnyMethods, PyModule},
21    Bound, DowncastError, Py, PyAny, PyErr, PyObject, PyResult, PyTypeInfo, Python,
22};
23
24use crate::borrow::{PyReadonlyArray, PyReadwriteArray};
25use crate::cold;
26use crate::convert::{ArrayExt, IntoPyArray, NpyIndex, ToNpyDims, ToPyArray};
27use crate::dtype::{Element, PyArrayDescrMethods};
28use crate::error::{
29    BorrowError, DimensionalityError, FromVecError, IgnoreError, NotContiguousError, TypeError,
30    DIMENSIONALITY_MISMATCH_ERR, MAX_DIMENSIONALITY_ERR,
31};
32use crate::npyffi::{self, npy_intp, NPY_ORDER, PY_ARRAY_API};
33use crate::slice_container::PySliceContainer;
34use crate::untyped_array::{PyUntypedArray, PyUntypedArrayMethods};
35
36/// A safe, statically-typed wrapper for NumPy's [`ndarray`][ndarray] class.
37///
38/// # Memory location
39///
40/// - Allocated by Rust: Constructed via [`IntoPyArray`] or
41///   [`from_vec`][Self::from_vec] or [`from_owned_array`][Self::from_owned_array].
42///
43/// These methods transfers ownership of the Rust allocation into a suitable Python object
44/// and uses the memory as the internal buffer backing the NumPy array.
45///
46/// Please note that some destructive methods like [`resize`][Self::resize] will fail
47/// when used with this kind of array as NumPy cannot reallocate the internal buffer.
48///
49/// - Allocated by NumPy: Constructed via other methods, like [`ToPyArray`] or
50///   [`from_slice`][Self::from_slice] or [`from_array`][Self::from_array].
51///
52/// These methods allocate memory in Python's private heap via NumPy's API.
53///
54/// In both cases, `PyArray` is managed by Python so it can neither be moved from
55/// nor deallocated manually.
56///
57/// # References
58///
59/// Like [`new`][Self::new], all constructor methods of `PyArray` return a shared reference `&PyArray`
60/// instead of an owned value. This design follows [PyO3's ownership concept][pyo3-memory],
61/// i.e. the return value is GIL-bound owning reference into Python's heap.
62///
63/// # Element type and dimensionality
64///
65/// `PyArray` has two type parametes `T` and `D`.
66/// `T` represents the type of its elements, e.g. [`f32`] or [`PyObject`].
67/// `D` represents its dimensionality, e.g [`Ix2`][type@Ix2] or [`IxDyn`][type@IxDyn].
68///
69/// Element types are Rust types which implement the [`Element`] trait.
70/// Dimensions are represented by the [`ndarray::Dimension`] trait.
71///
72/// Typically, `Ix1, Ix2, ...` are used for fixed dimensionality arrays,
73/// and `IxDyn` is used for dynamic dimensionality arrays. Type aliases
74/// for combining `PyArray` with these types are provided, e.g. [`PyArray1`] or [`PyArrayDyn`].
75///
76/// To specify concrete dimension like `3×4×5`, types which implement the [`ndarray::IntoDimension`]
77/// trait are used. Typically, this means arrays like `[3, 4, 5]` or tuples like `(3, 4, 5)`.
78///
79/// # Example
80///
81/// ```
82/// use numpy::{PyArray, PyArrayMethods};
83/// use ndarray::{array, Array};
84/// use pyo3::Python;
85///
86/// Python::with_gil(|py| {
87///     let pyarray = PyArray::arange(py, 0., 4., 1.).reshape([2, 2]).unwrap();
88///     let array = array![[3., 4.], [5., 6.]];
89///
90///     assert_eq!(
91///         array.dot(&pyarray.readonly().as_array()),
92///         array![[8., 15.], [12., 23.]]
93///     );
94/// });
95/// ```
96///
97/// [ndarray]: https://numpy.org/doc/stable/reference/generated/numpy.ndarray.html
98/// [pyo3-memory]: https://pyo3.rs/main/memory.html
99#[repr(transparent)]
100pub struct PyArray<T, D>(PyAny, PhantomData<T>, PhantomData<D>);
101
102/// Zero-dimensional array.
103pub type PyArray0<T> = PyArray<T, Ix0>;
104/// One-dimensional array.
105pub type PyArray1<T> = PyArray<T, Ix1>;
106/// Two-dimensional array.
107pub type PyArray2<T> = PyArray<T, Ix2>;
108/// Three-dimensional array.
109pub type PyArray3<T> = PyArray<T, Ix3>;
110/// Four-dimensional array.
111pub type PyArray4<T> = PyArray<T, Ix4>;
112/// Five-dimensional array.
113pub type PyArray5<T> = PyArray<T, Ix5>;
114/// Six-dimensional array.
115pub type PyArray6<T> = PyArray<T, Ix6>;
116/// Dynamic-dimensional array.
117pub type PyArrayDyn<T> = PyArray<T, IxDyn>;
118
119/// Returns a handle to NumPy's multiarray module.
120pub fn get_array_module<'py>(py: Python<'py>) -> PyResult<Bound<'_, PyModule>> {
121    PyModule::import(py, npyffi::array::mod_name(py)?)
122}
123
124impl<T, D> DerefToPyAny for PyArray<T, D> {}
125
126unsafe impl<T: Element, D: Dimension> PyTypeInfo for PyArray<T, D> {
127    const NAME: &'static str = "PyArray<T, D>";
128    const MODULE: Option<&'static str> = Some("numpy");
129
130    fn type_object_raw<'py>(py: Python<'py>) -> *mut ffi::PyTypeObject {
131        unsafe { npyffi::PY_ARRAY_API.get_type_object(py, npyffi::NpyTypes::PyArray_Type) }
132    }
133
134    fn is_type_of(ob: &Bound<'_, PyAny>) -> bool {
135        Self::extract::<IgnoreError>(ob).is_ok()
136    }
137}
138
139impl<T: Element, D: Dimension> PyArray<T, D> {
140    fn extract<'a, 'py, E>(ob: &'a Bound<'py, PyAny>) -> Result<&'a Bound<'py, Self>, E>
141    where
142        E: From<DowncastError<'a, 'py>> + From<DimensionalityError> + From<TypeError<'py>>,
143    {
144        // Check if the object is an array.
145        let array = unsafe {
146            if npyffi::PyArray_Check(ob.py(), ob.as_ptr()) == 0 {
147                return Err(DowncastError::new(ob, <Self as PyTypeInfo>::NAME).into());
148            }
149            ob.downcast_unchecked::<Self>()
150        };
151
152        // Check if the dimensionality matches `D`.
153        let src_ndim = array.ndim();
154        if let Some(dst_ndim) = D::NDIM {
155            if src_ndim != dst_ndim {
156                return Err(DimensionalityError::new(src_ndim, dst_ndim).into());
157            }
158        }
159
160        // Check if the element type matches `T`.
161        let src_dtype = array.dtype();
162        let dst_dtype = T::get_dtype(ob.py());
163        if !src_dtype.is_equiv_to(&dst_dtype) {
164            return Err(TypeError::new(src_dtype, dst_dtype).into());
165        }
166
167        Ok(array)
168    }
169
170    /// Creates a new uninitialized NumPy array.
171    ///
172    /// If `is_fortran` is true, then it has Fortran/column-major order,
173    /// otherwise it has C/row-major order.
174    ///
175    /// # Safety
176    ///
177    /// The returned array will always be safe to be dropped as the elements must either
178    /// be trivially copyable (as indicated by `<T as Element>::IS_COPY`) or be pointers
179    /// into Python's heap, which NumPy will automatically zero-initialize.
180    ///
181    /// However, the elements themselves will not be valid and should be initialized manually
182    /// using raw pointers obtained via [`uget_raw`][Self::uget_raw]. Before that, all methods
183    /// which produce references to the elements invoke undefined behaviour. In particular,
184    /// zero-initialized pointers are _not_ valid instances of `PyObject`.
185    ///
186    /// # Example
187    ///
188    /// ```
189    /// use numpy::prelude::*;
190    /// use numpy::PyArray3;
191    /// use pyo3::Python;
192    ///
193    /// Python::with_gil(|py| {
194    ///     let arr = unsafe {
195    ///         let arr = PyArray3::<i32>::new(py, [4, 5, 6], false);
196    ///
197    ///         for i in 0..4 {
198    ///             for j in 0..5 {
199    ///                 for k in 0..6 {
200    ///                     arr.uget_raw([i, j, k]).write((i * j * k) as i32);
201    ///                 }
202    ///             }
203    ///         }
204    ///
205    ///         arr
206    ///     };
207    ///
208    ///     assert_eq!(arr.shape(), &[4, 5, 6]);
209    /// });
210    /// ```
211    pub unsafe fn new<'py, ID>(py: Python<'py>, dims: ID, is_fortran: bool) -> Bound<'py, Self>
212    where
213        ID: IntoDimension<Dim = D>,
214    {
215        let flags = c_int::from(is_fortran);
216        Self::new_uninit(py, dims, ptr::null_mut(), flags)
217    }
218
219    /// Deprecated name for [`PyArray::new`].
220    ///
221    /// # Safety
222    /// See [`PyArray::new`].
223    #[deprecated(since = "0.23.0", note = "renamed to `PyArray::new`")]
224    #[inline]
225    pub unsafe fn new_bound<'py, ID>(
226        py: Python<'py>,
227        dims: ID,
228        is_fortran: bool,
229    ) -> Bound<'py, Self>
230    where
231        ID: IntoDimension<Dim = D>,
232    {
233        Self::new(py, dims, is_fortran)
234    }
235
236    pub(crate) unsafe fn new_uninit<'py, ID>(
237        py: Python<'py>,
238        dims: ID,
239        strides: *const npy_intp,
240        flag: c_int,
241    ) -> Bound<'py, Self>
242    where
243        ID: IntoDimension<Dim = D>,
244    {
245        let mut dims = dims.into_dimension();
246        let ptr = PY_ARRAY_API.PyArray_NewFromDescr(
247            py,
248            PY_ARRAY_API.get_type_object(py, npyffi::NpyTypes::PyArray_Type),
249            T::get_dtype(py).into_dtype_ptr(),
250            dims.ndim_cint(),
251            dims.as_dims_ptr(),
252            strides as *mut npy_intp, // strides
253            ptr::null_mut(),          // data
254            flag,                     // flag
255            ptr::null_mut(),          // obj
256        );
257
258        Bound::from_owned_ptr(py, ptr).downcast_into_unchecked()
259    }
260
261    unsafe fn new_with_data<'py, ID>(
262        py: Python<'py>,
263        dims: ID,
264        strides: *const npy_intp,
265        data_ptr: *const T,
266        container: *mut PyAny,
267    ) -> Bound<'py, Self>
268    where
269        ID: IntoDimension<Dim = D>,
270    {
271        let mut dims = dims.into_dimension();
272        let ptr = PY_ARRAY_API.PyArray_NewFromDescr(
273            py,
274            PY_ARRAY_API.get_type_object(py, npyffi::NpyTypes::PyArray_Type),
275            T::get_dtype(py).into_dtype_ptr(),
276            dims.ndim_cint(),
277            dims.as_dims_ptr(),
278            strides as *mut npy_intp,    // strides
279            data_ptr as *mut c_void,     // data
280            npyffi::NPY_ARRAY_WRITEABLE, // flag
281            ptr::null_mut(),             // obj
282        );
283
284        PY_ARRAY_API.PyArray_SetBaseObject(
285            py,
286            ptr as *mut npyffi::PyArrayObject,
287            container as *mut ffi::PyObject,
288        );
289
290        Bound::from_owned_ptr(py, ptr).downcast_into_unchecked()
291    }
292
293    pub(crate) unsafe fn from_raw_parts<'py>(
294        py: Python<'py>,
295        dims: D,
296        strides: *const npy_intp,
297        data_ptr: *const T,
298        container: PySliceContainer,
299    ) -> Bound<'py, Self> {
300        let container = Bound::new(py, container)
301            .expect("Failed to create slice container")
302            .into_ptr();
303
304        Self::new_with_data(py, dims, strides, data_ptr, container.cast())
305    }
306
307    /// Creates a NumPy array backed by `array` and ties its ownership to the Python object `container`.
308    ///
309    /// The resulting NumPy array will be writeable from Python space.  If this is undesireable, use
310    /// [PyReadwriteArray::make_nonwriteable].
311    ///
312    /// # Safety
313    ///
314    /// `container` is set as a base object of the returned array which must not be dropped until `container` is dropped.
315    /// Furthermore, `array` must not be reallocated from the time this method is called and until `container` is dropped.
316    ///
317    /// # Example
318    ///
319    /// ```rust
320    /// # use pyo3::prelude::*;
321    /// # use numpy::{ndarray::Array1, PyArray1};
322    /// #
323    /// #[pyclass]
324    /// struct Owner {
325    ///     array: Array1<f64>,
326    /// }
327    ///
328    /// #[pymethods]
329    /// impl Owner {
330    ///     #[getter]
331    ///     fn array<'py>(this: Bound<'py, Self>) -> Bound<'py, PyArray1<f64>> {
332    ///         let array = &this.borrow().array;
333    ///
334    ///         // SAFETY: The memory backing `array` will stay valid as long as this object is alive
335    ///         // as we do not modify `array` in any way which would cause it to be reallocated.
336    ///         unsafe { PyArray1::borrow_from_array(array, this.into_any()) }
337    ///     }
338    /// }
339    /// ```
340    pub unsafe fn borrow_from_array<'py, S>(
341        array: &ArrayBase<S, D>,
342        container: Bound<'py, PyAny>,
343    ) -> Bound<'py, Self>
344    where
345        S: Data<Elem = T>,
346    {
347        let (strides, dims) = (array.npy_strides(), array.raw_dim());
348        let data_ptr = array.as_ptr();
349
350        let py = container.py();
351
352        Self::new_with_data(
353            py,
354            dims,
355            strides.as_ptr(),
356            data_ptr,
357            container.into_ptr().cast(),
358        )
359    }
360
361    /// Deprecated name for [`PyArray::borrow_from_array`].
362    ///
363    /// # Safety
364    /// See [`PyArray::borrow_from_array`]
365    #[deprecated(since = "0.23.0", note = "renamed to `PyArray::borrow_from_array`")]
366    #[inline]
367    pub unsafe fn borrow_from_array_bound<'py, S>(
368        array: &ArrayBase<S, D>,
369        container: Bound<'py, PyAny>,
370    ) -> Bound<'py, Self>
371    where
372        S: Data<Elem = T>,
373    {
374        Self::borrow_from_array(array, container)
375    }
376
377    /// Construct a new NumPy array filled with zeros.
378    ///
379    /// If `is_fortran` is true, then it has Fortran/column-major order,
380    /// otherwise it has C/row-major order.
381    ///
382    /// For arrays of Python objects, this will fill the array
383    /// with valid pointers to zero-valued Python integer objects.
384    ///
385    /// See also [`numpy.zeros`][numpy-zeros] and [`PyArray_Zeros`][PyArray_Zeros].
386    ///
387    /// # Example
388    ///
389    /// ```
390    /// use numpy::{PyArray2, PyArrayMethods};
391    /// use pyo3::Python;
392    ///
393    /// Python::with_gil(|py| {
394    ///     let pyarray = PyArray2::<usize>::zeros(py, [2, 2], true);
395    ///
396    ///     assert_eq!(pyarray.readonly().as_slice().unwrap(), [0; 4]);
397    /// });
398    /// ```
399    ///
400    /// [numpy-zeros]: https://numpy.org/doc/stable/reference/generated/numpy.zeros.html
401    /// [PyArray_Zeros]: https://numpy.org/doc/stable/reference/c-api/array.html#c.PyArray_Zeros
402    pub fn zeros<ID>(py: Python<'_>, dims: ID, is_fortran: bool) -> Bound<'_, Self>
403    where
404        ID: IntoDimension<Dim = D>,
405    {
406        let mut dims = dims.into_dimension();
407        unsafe {
408            let ptr = PY_ARRAY_API.PyArray_Zeros(
409                py,
410                dims.ndim_cint(),
411                dims.as_dims_ptr(),
412                T::get_dtype(py).into_dtype_ptr(),
413                if is_fortran { -1 } else { 0 },
414            );
415            Bound::from_owned_ptr(py, ptr).downcast_into_unchecked()
416        }
417    }
418
419    /// Deprecated name for [`PyArray::zeros`].
420    #[deprecated(since = "0.23.0", note = "renamed to `PyArray::zeros`")]
421    #[inline]
422    pub fn zeros_bound<ID>(py: Python<'_>, dims: ID, is_fortran: bool) -> Bound<'_, Self>
423    where
424        ID: IntoDimension<Dim = D>,
425    {
426        Self::zeros(py, dims, is_fortran)
427    }
428
429    /// Constructs a NumPy from an [`ndarray::Array`]
430    ///
431    /// This method uses the internal [`Vec`] of the [`ndarray::Array`] as the base object of the NumPy array.
432    ///
433    /// # Example
434    ///
435    /// ```
436    /// use numpy::{PyArray, PyArrayMethods};
437    /// use ndarray::array;
438    /// use pyo3::Python;
439    ///
440    /// Python::with_gil(|py| {
441    ///     let pyarray = PyArray::from_owned_array(py, array![[1, 2], [3, 4]]);
442    ///
443    ///     assert_eq!(pyarray.readonly().as_array(), array![[1, 2], [3, 4]]);
444    /// });
445    /// ```
446    pub fn from_owned_array(py: Python<'_>, mut arr: Array<T, D>) -> Bound<'_, Self> {
447        let (strides, dims) = (arr.npy_strides(), arr.raw_dim());
448        let data_ptr = arr.as_mut_ptr();
449        unsafe {
450            Self::from_raw_parts(
451                py,
452                dims,
453                strides.as_ptr(),
454                data_ptr,
455                PySliceContainer::from(arr),
456            )
457        }
458    }
459    /// Deprecated name for [`PyArray::from_owned_array`].
460    #[deprecated(since = "0.23.0", note = "renamed to `PyArray::from_owned_array`")]
461    #[inline]
462    pub fn from_owned_array_bound(py: Python<'_>, arr: Array<T, D>) -> Bound<'_, Self> {
463        Self::from_owned_array(py, arr)
464    }
465
466    /// Construct a NumPy array from a [`ndarray::ArrayBase`].
467    ///
468    /// This method allocates memory in Python's heap via the NumPy API,
469    /// and then copies all elements of the array there.
470    ///
471    /// # Example
472    ///
473    /// ```
474    /// use numpy::{PyArray, PyArrayMethods};
475    /// use ndarray::array;
476    /// use pyo3::Python;
477    ///
478    /// Python::with_gil(|py| {
479    ///     let pyarray = PyArray::from_array(py, &array![[1, 2], [3, 4]]);
480    ///
481    ///     assert_eq!(pyarray.readonly().as_array(), array![[1, 2], [3, 4]]);
482    /// });
483    /// ```
484    pub fn from_array<'py, S>(py: Python<'py>, arr: &ArrayBase<S, D>) -> Bound<'py, Self>
485    where
486        S: Data<Elem = T>,
487    {
488        ToPyArray::to_pyarray(arr, py)
489    }
490
491    /// Deprecated name for [`PyArray::from_array`].
492    #[deprecated(since = "0.23.0", note = "renamed to `PyArray::from_array`")]
493    #[inline]
494    pub fn from_array_bound<'py, S>(py: Python<'py>, arr: &ArrayBase<S, D>) -> Bound<'py, Self>
495    where
496        S: Data<Elem = T>,
497    {
498        Self::from_array(py, arr)
499    }
500}
501
502impl<D: Dimension> PyArray<PyObject, D> {
503    /// Construct a NumPy array containing objects stored in a [`ndarray::Array`]
504    ///
505    /// This method uses the internal [`Vec`] of the [`ndarray::Array`] as the base object of the NumPy array.
506    ///
507    /// # Example
508    ///
509    /// ```
510    /// use ndarray::array;
511    /// use pyo3::{pyclass, Py, Python, types::PyAnyMethods};
512    /// use numpy::{PyArray, PyArrayMethods};
513    ///
514    /// #[pyclass]
515    /// # #[allow(dead_code)]
516    /// struct CustomElement {
517    ///     foo: i32,
518    ///     bar: f64,
519    /// }
520    ///
521    /// Python::with_gil(|py| {
522    ///     let array = array![
523    ///         Py::new(py, CustomElement {
524    ///             foo: 1,
525    ///             bar: 2.0,
526    ///         }).unwrap(),
527    ///         Py::new(py, CustomElement {
528    ///             foo: 3,
529    ///             bar: 4.0,
530    ///         }).unwrap(),
531    ///     ];
532    ///
533    ///     let pyarray = PyArray::from_owned_object_array(py, array);
534    ///
535    ///     assert!(pyarray.readonly().as_array().get(0).unwrap().bind(py).is_instance_of::<CustomElement>());
536    /// });
537    /// ```
538    pub fn from_owned_object_array<T>(py: Python<'_>, mut arr: Array<Py<T>, D>) -> Bound<'_, Self> {
539        let (strides, dims) = (arr.npy_strides(), arr.raw_dim());
540        let data_ptr = arr.as_mut_ptr() as *const PyObject;
541        unsafe {
542            Self::from_raw_parts(
543                py,
544                dims,
545                strides.as_ptr(),
546                data_ptr,
547                PySliceContainer::from(arr),
548            )
549        }
550    }
551
552    /// Deprecated name for [`PyArray::from_owned_object_array`].
553    #[deprecated(
554        since = "0.23.0",
555        note = "renamed to `PyArray::from_owned_object_array`"
556    )]
557    #[inline]
558    pub fn from_owned_object_array_bound<T>(
559        py: Python<'_>,
560        arr: Array<Py<T>, D>,
561    ) -> Bound<'_, Self> {
562        Self::from_owned_object_array(py, arr)
563    }
564}
565
566impl<T: Element> PyArray<T, Ix1> {
567    /// Construct a one-dimensional array from a [mod@slice].
568    ///
569    /// # Example
570    ///
571    /// ```
572    /// use numpy::{PyArray, PyArrayMethods};
573    /// use pyo3::Python;
574    ///
575    /// Python::with_gil(|py| {
576    ///     let slice = &[1, 2, 3, 4, 5];
577    ///     let pyarray = PyArray::from_slice(py, slice);
578    ///     assert_eq!(pyarray.readonly().as_slice().unwrap(), &[1, 2, 3, 4, 5]);
579    /// });
580    /// ```
581    pub fn from_slice<'py>(py: Python<'py>, slice: &[T]) -> Bound<'py, Self> {
582        unsafe {
583            let array = PyArray::new(py, [slice.len()], false);
584            let mut data_ptr = array.data();
585            clone_elements(py, slice, &mut data_ptr);
586            array
587        }
588    }
589
590    /// Deprecated name for [`PyArray::from_slice`].
591    #[deprecated(since = "0.23.0", note = "renamed to `PyArray::from_slice`")]
592    #[inline]
593    pub fn from_slice_bound<'py>(py: Python<'py>, slice: &[T]) -> Bound<'py, Self> {
594        Self::from_slice(py, slice)
595    }
596
597    /// Construct a one-dimensional array from a [`Vec<T>`][Vec].
598    ///
599    /// # Example
600    ///
601    /// ```
602    /// use numpy::{PyArray, PyArrayMethods};
603    /// use pyo3::Python;
604    ///
605    /// Python::with_gil(|py| {
606    ///     let vec = vec![1, 2, 3, 4, 5];
607    ///     let pyarray = PyArray::from_vec(py, vec);
608    ///     assert_eq!(pyarray.readonly().as_slice().unwrap(), &[1, 2, 3, 4, 5]);
609    /// });
610    /// ```
611    #[inline(always)]
612    pub fn from_vec<'py>(py: Python<'py>, vec: Vec<T>) -> Bound<'py, Self> {
613        vec.into_pyarray(py)
614    }
615
616    /// Deprecated name for [`PyArray::from_vec`].
617    #[deprecated(since = "0.23.0", note = "renamed to `PyArray::from_vec`")]
618    #[inline]
619    pub fn from_vec_bound<'py>(py: Python<'py>, vec: Vec<T>) -> Bound<'py, Self> {
620        Self::from_vec(py, vec)
621    }
622
623    /// Construct a one-dimensional array from an [`Iterator`].
624    ///
625    /// If no reliable [`size_hint`][Iterator::size_hint] is available,
626    /// this method can allocate memory multiple times, which can hurt performance.
627    ///
628    /// # Example
629    ///
630    /// ```
631    /// use numpy::{PyArray, PyArrayMethods};
632    /// use pyo3::Python;
633    ///
634    /// Python::with_gil(|py| {
635    ///     let pyarray = PyArray::from_iter(py, "abcde".chars().map(u32::from));
636    ///     assert_eq!(pyarray.readonly().as_slice().unwrap(), &[97, 98, 99, 100, 101]);
637    /// });
638    /// ```
639    pub fn from_iter<I>(py: Python<'_>, iter: I) -> Bound<'_, Self>
640    where
641        I: IntoIterator<Item = T>,
642    {
643        let data = iter.into_iter().collect::<Vec<_>>();
644        data.into_pyarray(py)
645    }
646
647    /// Deprecated name for [`PyArray::from_iter`].
648    #[deprecated(since = "0.23.0", note = "renamed to `PyArray::from_iter`")]
649    #[inline]
650    pub fn from_iter_bound<I>(py: Python<'_>, iter: I) -> Bound<'_, Self>
651    where
652        I: IntoIterator<Item = T>,
653    {
654        Self::from_iter(py, iter)
655    }
656}
657
658impl<T: Element> PyArray<T, Ix2> {
659    /// Construct a two-dimension array from a [`Vec<Vec<T>>`][Vec].
660    ///
661    /// This function checks all dimensions of the inner vectors and returns
662    /// an error if they are not all equal.
663    ///
664    /// # Example
665    ///
666    /// ```
667    /// use numpy::{PyArray, PyArrayMethods};
668    /// use pyo3::Python;
669    /// use ndarray::array;
670    ///
671    /// Python::with_gil(|py| {
672    ///     let vec2 = vec![vec![11, 12], vec![21, 22]];
673    ///     let pyarray = PyArray::from_vec2(py, &vec2).unwrap();
674    ///     assert_eq!(pyarray.readonly().as_array(), array![[11, 12], [21, 22]]);
675    ///
676    ///     let ragged_vec2 = vec![vec![11, 12], vec![21]];
677    ///     assert!(PyArray::from_vec2(py, &ragged_vec2).is_err());
678    /// });
679    /// ```
680    pub fn from_vec2<'py>(py: Python<'py>, v: &[Vec<T>]) -> Result<Bound<'py, Self>, FromVecError> {
681        let len2 = v.first().map_or(0, |v| v.len());
682        let dims = [v.len(), len2];
683        // SAFETY: The result of `Self::new` is always safe to drop.
684        unsafe {
685            let array = Self::new(py, dims, false);
686            let mut data_ptr = array.data();
687            for v in v {
688                if v.len() != len2 {
689                    cold();
690                    return Err(FromVecError::new(v.len(), len2));
691                }
692                clone_elements(py, v, &mut data_ptr);
693            }
694            Ok(array)
695        }
696    }
697
698    /// Deprecated name for [`PyArray::from_vec2`].
699    #[deprecated(since = "0.23.0", note = "renamed to `PyArray::from_vec2`")]
700    #[inline]
701    pub fn from_vec2_bound<'py>(
702        py: Python<'py>,
703        v: &[Vec<T>],
704    ) -> Result<Bound<'py, Self>, FromVecError> {
705        Self::from_vec2(py, v)
706    }
707}
708
709impl<T: Element> PyArray<T, Ix3> {
710    /// Construct a three-dimensional array from a [`Vec<Vec<Vec<T>>>`][Vec].
711    ///
712    /// This function checks all dimensions of the inner vectors and returns
713    /// an error if they are not all equal.
714    ///
715    /// # Example
716    ///
717    /// ```
718    /// use numpy::{PyArray, PyArrayMethods};
719    /// use pyo3::Python;
720    /// use ndarray::array;
721    ///
722    /// Python::with_gil(|py| {
723    ///     let vec3 = vec![
724    ///         vec![vec![111, 112], vec![121, 122]],
725    ///         vec![vec![211, 212], vec![221, 222]],
726    ///     ];
727    ///     let pyarray = PyArray::from_vec3(py, &vec3).unwrap();
728    ///     assert_eq!(
729    ///         pyarray.readonly().as_array(),
730    ///         array![[[111, 112], [121, 122]], [[211, 212], [221, 222]]]
731    ///     );
732    ///
733    ///     let ragged_vec3 = vec![
734    ///         vec![vec![111, 112], vec![121, 122]],
735    ///         vec![vec![211], vec![221, 222]],
736    ///     ];
737    ///     assert!(PyArray::from_vec3(py, &ragged_vec3).is_err());
738    /// });
739    /// ```
740    pub fn from_vec3<'py>(
741        py: Python<'py>,
742        v: &[Vec<Vec<T>>],
743    ) -> Result<Bound<'py, Self>, FromVecError> {
744        let len2 = v.first().map_or(0, |v| v.len());
745        let len3 = v.first().map_or(0, |v| v.first().map_or(0, |v| v.len()));
746        let dims = [v.len(), len2, len3];
747        // SAFETY: The result of `Self::new` is always safe to drop.
748        unsafe {
749            let array = Self::new(py, dims, false);
750            let mut data_ptr = array.data();
751            for v in v {
752                if v.len() != len2 {
753                    cold();
754                    return Err(FromVecError::new(v.len(), len2));
755                }
756                for v in v {
757                    if v.len() != len3 {
758                        cold();
759                        return Err(FromVecError::new(v.len(), len3));
760                    }
761                    clone_elements(py, v, &mut data_ptr);
762                }
763            }
764            Ok(array)
765        }
766    }
767
768    /// Deprecated name for [`PyArray::from_vec3`].
769    #[deprecated(since = "0.23.0", note = "renamed to `PyArray::from_vec3`")]
770    #[inline]
771    pub fn from_vec3_bound<'py>(
772        py: Python<'py>,
773        v: &[Vec<Vec<T>>],
774    ) -> Result<Bound<'py, Self>, FromVecError> {
775        Self::from_vec3(py, v)
776    }
777}
778
779impl<T: Element + AsPrimitive<f64>> PyArray<T, Ix1> {
780    /// Return evenly spaced values within a given interval.
781    ///
782    /// See [numpy.arange][numpy.arange] for the Python API and [PyArray_Arange][PyArray_Arange] for the C API.
783    ///
784    /// # Example
785    ///
786    /// ```
787    /// use numpy::{PyArray, PyArrayMethods};
788    /// use pyo3::Python;
789    ///
790    /// Python::with_gil(|py| {
791    ///     let pyarray = PyArray::arange(py, 2.0, 4.0, 0.5);
792    ///     assert_eq!(pyarray.readonly().as_slice().unwrap(), &[2.0, 2.5, 3.0, 3.5]);
793    ///
794    ///     let pyarray = PyArray::arange(py, -2, 4, 3);
795    ///     assert_eq!(pyarray.readonly().as_slice().unwrap(), &[-2, 1]);
796    /// });
797    /// ```
798    ///
799    /// [numpy.arange]: https://numpy.org/doc/stable/reference/generated/numpy.arange.html
800    /// [PyArray_Arange]: https://numpy.org/doc/stable/reference/c-api/array.html#c.PyArray_Arange
801    pub fn arange<'py>(py: Python<'py>, start: T, stop: T, step: T) -> Bound<'py, Self> {
802        unsafe {
803            let ptr = PY_ARRAY_API.PyArray_Arange(
804                py,
805                start.as_(),
806                stop.as_(),
807                step.as_(),
808                T::get_dtype(py).num(),
809            );
810            Bound::from_owned_ptr(py, ptr).downcast_into_unchecked()
811        }
812    }
813
814    /// Deprecated name for [`PyArray::arange`].
815    #[deprecated(since = "0.23.0", note = "renamed to `PyArray::arange`")]
816    #[inline]
817    pub fn arange_bound<'py>(py: Python<'py>, start: T, stop: T, step: T) -> Bound<'py, Self> {
818        Self::arange(py, start, stop, step)
819    }
820}
821
822unsafe fn clone_elements<T: Element>(py: Python<'_>, elems: &[T], data_ptr: &mut *mut T) {
823    if T::IS_COPY {
824        ptr::copy_nonoverlapping(elems.as_ptr(), *data_ptr, elems.len());
825        *data_ptr = data_ptr.add(elems.len());
826    } else {
827        for elem in elems {
828            data_ptr.write(elem.clone_ref(py));
829            *data_ptr = data_ptr.add(1);
830        }
831    }
832}
833
834/// Implementation of functionality for [`PyArray<T, D>`].
835#[doc(alias = "PyArray")]
836pub trait PyArrayMethods<'py, T, D>: PyUntypedArrayMethods<'py> {
837    /// Access an untyped representation of this array.
838    fn as_untyped(&self) -> &Bound<'py, PyUntypedArray>;
839
840    /// Returns a pointer to the first element of the array.
841    fn data(&self) -> *mut T;
842
843    /// Same as [`shape`][PyUntypedArray::shape], but returns `D` instead of `&[usize]`.
844    #[inline(always)]
845    fn dims(&self) -> D
846    where
847        D: Dimension,
848    {
849        D::from_dimension(&Dim(self.shape())).expect(DIMENSIONALITY_MISMATCH_ERR)
850    }
851
852    /// Returns an immutable view of the internal data as a slice.
853    ///
854    /// # Safety
855    ///
856    /// Calling this method is undefined behaviour if the underlying array
857    /// is aliased mutably by other instances of `PyArray`
858    /// or concurrently modified by Python or other native code.
859    ///
860    /// Please consider the safe alternative [`PyReadonlyArray::as_slice`].
861    unsafe fn as_slice(&self) -> Result<&[T], NotContiguousError>
862    where
863        T: Element,
864        D: Dimension,
865    {
866        if self.is_contiguous() {
867            Ok(slice::from_raw_parts(self.data(), self.len()))
868        } else {
869            Err(NotContiguousError)
870        }
871    }
872
873    /// Returns a mutable view of the internal data as a slice.
874    ///
875    /// # Safety
876    ///
877    /// Calling this method is undefined behaviour if the underlying array
878    /// is aliased immutably or mutably by other instances of [`PyArray`]
879    /// or concurrently modified by Python or other native code.
880    ///
881    /// Please consider the safe alternative [`PyReadwriteArray::as_slice_mut`].
882    unsafe fn as_slice_mut(&self) -> Result<&mut [T], NotContiguousError>
883    where
884        T: Element,
885        D: Dimension,
886    {
887        if self.is_contiguous() {
888            Ok(slice::from_raw_parts_mut(self.data(), self.len()))
889        } else {
890            Err(NotContiguousError)
891        }
892    }
893
894    /// Get a reference of the specified element if the given index is valid.
895    ///
896    /// # Safety
897    ///
898    /// Calling this method is undefined behaviour if the underlying array
899    /// is aliased mutably by other instances of `PyArray`
900    /// or concurrently modified by Python or other native code.
901    ///
902    /// Consider using safe alternatives like [`PyReadonlyArray::get`].
903    ///
904    /// # Example
905    ///
906    /// ```
907    /// use numpy::{PyArray, PyArrayMethods};
908    /// use pyo3::Python;
909    ///
910    /// Python::with_gil(|py| {
911    ///     let pyarray = PyArray::arange(py, 0, 16, 1).reshape([2, 2, 4]).unwrap();
912    ///
913    ///     assert_eq!(unsafe { *pyarray.get([1, 0, 3]).unwrap() }, 11);
914    /// });
915    /// ```
916    unsafe fn get(&self, index: impl NpyIndex<Dim = D>) -> Option<&T>
917    where
918        T: Element,
919        D: Dimension;
920
921    /// Same as [`get`][Self::get], but returns `Option<&mut T>`.
922    ///
923    /// # Safety
924    ///
925    /// Calling this method is undefined behaviour if the underlying array
926    /// is aliased immutably or mutably by other instances of [`PyArray`]
927    /// or concurrently modified by Python or other native code.
928    ///
929    /// Consider using safe alternatives like [`PyReadwriteArray::get_mut`].
930    ///
931    /// # Example
932    ///
933    /// ```
934    /// use numpy::{PyArray, PyArrayMethods};
935    /// use pyo3::Python;
936    ///
937    /// Python::with_gil(|py| {
938    ///     let pyarray = PyArray::arange(py, 0, 16, 1).reshape([2, 2, 4]).unwrap();
939    ///
940    ///     unsafe {
941    ///         *pyarray.get_mut([1, 0, 3]).unwrap() = 42;
942    ///     }
943    ///
944    ///     assert_eq!(unsafe { *pyarray.get([1, 0, 3]).unwrap() }, 42);
945    /// });
946    /// ```
947    unsafe fn get_mut(&self, index: impl NpyIndex<Dim = D>) -> Option<&mut T>
948    where
949        T: Element,
950        D: Dimension;
951
952    /// Get an immutable reference of the specified element,
953    /// without checking the given index.
954    ///
955    /// See [`NpyIndex`] for what types can be used as the index.
956    ///
957    /// # Safety
958    ///
959    /// Passing an invalid index is undefined behavior.
960    /// The element must also have been initialized and
961    /// all other references to it is must also be shared.
962    ///
963    /// See [`PyReadonlyArray::get`] for a safe alternative.
964    ///
965    /// # Example
966    ///
967    /// ```
968    /// use numpy::{PyArray, PyArrayMethods};
969    /// use pyo3::Python;
970    ///
971    /// Python::with_gil(|py| {
972    ///     let pyarray = PyArray::arange(py, 0, 16, 1).reshape([2, 2, 4]).unwrap();
973    ///
974    ///     assert_eq!(unsafe { *pyarray.uget([1, 0, 3]) }, 11);
975    /// });
976    /// ```
977    #[inline(always)]
978    unsafe fn uget<Idx>(&self, index: Idx) -> &T
979    where
980        T: Element,
981        D: Dimension,
982        Idx: NpyIndex<Dim = D>,
983    {
984        &*self.uget_raw(index)
985    }
986
987    /// Same as [`uget`](Self::uget), but returns `&mut T`.
988    ///
989    /// # Safety
990    ///
991    /// Passing an invalid index is undefined behavior.
992    /// The element must also have been initialized and
993    /// other references to it must not exist.
994    ///
995    /// See [`PyReadwriteArray::get_mut`] for a safe alternative.
996    #[inline(always)]
997    #[allow(clippy::mut_from_ref)]
998    unsafe fn uget_mut<Idx>(&self, index: Idx) -> &mut T
999    where
1000        T: Element,
1001        D: Dimension,
1002        Idx: NpyIndex<Dim = D>,
1003    {
1004        &mut *self.uget_raw(index)
1005    }
1006
1007    /// Same as [`uget`][Self::uget], but returns `*mut T`.
1008    ///
1009    /// # Safety
1010    ///
1011    /// Passing an invalid index is undefined behavior.
1012    #[inline(always)]
1013    unsafe fn uget_raw<Idx>(&self, index: Idx) -> *mut T
1014    where
1015        T: Element,
1016        D: Dimension,
1017        Idx: NpyIndex<Dim = D>,
1018    {
1019        let offset = index.get_unchecked::<T>(self.strides());
1020        self.data().offset(offset) as *mut _
1021    }
1022
1023    /// Get a copy of the specified element in the array.
1024    ///
1025    /// See [`NpyIndex`] for what types can be used as the index.
1026    ///
1027    /// # Example
1028    /// ```
1029    /// use numpy::{PyArray, PyArrayMethods};
1030    /// use pyo3::Python;
1031    ///
1032    /// Python::with_gil(|py| {
1033    ///     let pyarray = PyArray::arange(py, 0, 16, 1).reshape([2, 2, 4]).unwrap();
1034    ///
1035    ///     assert_eq!(pyarray.get_owned([1, 0, 3]), Some(11));
1036    /// });
1037    /// ```
1038    fn get_owned<Idx>(&self, index: Idx) -> Option<T>
1039    where
1040        T: Element,
1041        D: Dimension,
1042        Idx: NpyIndex<Dim = D>;
1043
1044    /// Turn an array with fixed dimensionality into one with dynamic dimensionality.
1045    fn to_dyn(&self) -> &Bound<'py, PyArray<T, IxDyn>>
1046    where
1047        T: Element,
1048        D: Dimension;
1049
1050    /// Returns a copy of the internal data of the array as a [`Vec`].
1051    ///
1052    /// Fails if the internal array is not contiguous. See also [`as_slice`][Self::as_slice].
1053    ///
1054    /// # Example
1055    ///
1056    /// ```
1057    /// use numpy::{PyArray2, PyArrayMethods};
1058    /// use pyo3::{Python, types::PyAnyMethods, ffi::c_str};
1059    ///
1060    /// # fn main() -> pyo3::PyResult<()> {
1061    /// Python::with_gil(|py| {
1062    ///     let pyarray= py
1063    ///         .eval(c_str!("__import__('numpy').array([[0, 1], [2, 3]], dtype='int64')"), None, None)?
1064    ///         .downcast_into::<PyArray2<i64>>()?;
1065    ///
1066    ///     assert_eq!(pyarray.to_vec()?, vec![0, 1, 2, 3]);
1067    /// #   Ok(())
1068    /// })
1069    /// # }
1070    /// ```
1071    fn to_vec(&self) -> Result<Vec<T>, NotContiguousError>
1072    where
1073        T: Element,
1074        D: Dimension;
1075
1076    /// Get an immutable borrow of the NumPy array
1077    fn try_readonly(&self) -> Result<PyReadonlyArray<'py, T, D>, BorrowError>
1078    where
1079        T: Element,
1080        D: Dimension;
1081
1082    /// Get an immutable borrow of the NumPy array
1083    ///
1084    /// # Panics
1085    ///
1086    /// Panics if the allocation backing the array is currently mutably borrowed.
1087    ///
1088    /// For a non-panicking variant, use [`try_readonly`][Self::try_readonly].
1089    fn readonly(&self) -> PyReadonlyArray<'py, T, D>
1090    where
1091        T: Element,
1092        D: Dimension,
1093    {
1094        self.try_readonly().unwrap()
1095    }
1096
1097    /// Get a mutable borrow of the NumPy array
1098    fn try_readwrite(&self) -> Result<PyReadwriteArray<'py, T, D>, BorrowError>
1099    where
1100        T: Element,
1101        D: Dimension;
1102
1103    /// Get a mutable borrow of the NumPy array
1104    ///
1105    /// # Panics
1106    ///
1107    /// Panics if the allocation backing the array is currently borrowed or
1108    /// if the array is [flagged as][flags] not writeable.
1109    ///
1110    /// For a non-panicking variant, use [`try_readwrite`][Self::try_readwrite].
1111    ///
1112    /// [flags]: https://numpy.org/doc/stable/reference/generated/numpy.ndarray.flags.html
1113    fn readwrite(&self) -> PyReadwriteArray<'py, T, D>
1114    where
1115        T: Element,
1116        D: Dimension,
1117    {
1118        self.try_readwrite().unwrap()
1119    }
1120
1121    /// Returns an [`ArrayView`] of the internal array.
1122    ///
1123    /// See also [`PyReadonlyArray::as_array`].
1124    ///
1125    /// # Safety
1126    ///
1127    /// Calling this method invalidates all exclusive references to the internal data, e.g. `&mut [T]` or `ArrayViewMut`.
1128    unsafe fn as_array(&self) -> ArrayView<'_, T, D>
1129    where
1130        T: Element,
1131        D: Dimension;
1132
1133    /// Returns an [`ArrayViewMut`] of the internal array.
1134    ///
1135    /// See also [`PyReadwriteArray::as_array_mut`].
1136    ///
1137    /// # Safety
1138    ///
1139    /// Calling this method invalidates all other references to the internal data, e.g. `ArrayView` or `ArrayViewMut`.
1140    unsafe fn as_array_mut(&self) -> ArrayViewMut<'_, T, D>
1141    where
1142        T: Element,
1143        D: Dimension;
1144
1145    /// Returns the internal array as [`RawArrayView`] enabling element access via raw pointers
1146    fn as_raw_array(&self) -> RawArrayView<T, D>
1147    where
1148        T: Element,
1149        D: Dimension;
1150
1151    /// Returns the internal array as [`RawArrayViewMut`] enabling element access via raw pointers
1152    fn as_raw_array_mut(&self) -> RawArrayViewMut<T, D>
1153    where
1154        T: Element,
1155        D: Dimension;
1156
1157    /// Get a copy of the array as an [`ndarray::Array`].
1158    ///
1159    /// # Example
1160    ///
1161    /// ```
1162    /// use numpy::{PyArray, PyArrayMethods};
1163    /// use ndarray::array;
1164    /// use pyo3::Python;
1165    ///
1166    /// Python::with_gil(|py| {
1167    ///     let pyarray = PyArray::arange(py, 0, 4, 1).reshape([2, 2]).unwrap();
1168    ///
1169    ///     assert_eq!(
1170    ///         pyarray.to_owned_array(),
1171    ///         array![[0, 1], [2, 3]]
1172    ///     )
1173    /// });
1174    /// ```
1175    fn to_owned_array(&self) -> Array<T, D>
1176    where
1177        T: Element,
1178        D: Dimension;
1179
1180    /// Copies `self` into `other`, performing a data type conversion if necessary.
1181    ///
1182    /// See also [`PyArray_CopyInto`][PyArray_CopyInto].
1183    ///
1184    /// # Example
1185    ///
1186    /// ```
1187    /// use numpy::{PyArray, PyArrayMethods};
1188    /// use pyo3::Python;
1189    ///
1190    /// Python::with_gil(|py| {
1191    ///     let pyarray_f = PyArray::arange(py, 2.0, 5.0, 1.0);
1192    ///     let pyarray_i = unsafe { PyArray::<i64, _>::new(py, [3], false) };
1193    ///
1194    ///     assert!(pyarray_f.copy_to(&pyarray_i).is_ok());
1195    ///
1196    ///     assert_eq!(pyarray_i.readonly().as_slice().unwrap(), &[2, 3, 4]);
1197    /// });
1198    /// ```
1199    ///
1200    /// [PyArray_CopyInto]: https://numpy.org/doc/stable/reference/c-api/array.html#c.PyArray_CopyInto
1201    fn copy_to<U: Element>(&self, other: &Bound<'py, PyArray<U, D>>) -> PyResult<()>
1202    where
1203        T: Element;
1204
1205    /// Cast the `PyArray<T>` to `PyArray<U>`, by allocating a new array.
1206    ///
1207    /// See also [`PyArray_CastToType`][PyArray_CastToType].
1208    ///
1209    /// # Example
1210    ///
1211    /// ```
1212    /// use numpy::{PyArray, PyArrayMethods};
1213    /// use pyo3::Python;
1214    ///
1215    /// Python::with_gil(|py| {
1216    ///     let pyarray_f = PyArray::arange(py, 2.0, 5.0, 1.0);
1217    ///
1218    ///     let pyarray_i = pyarray_f.cast::<i32>(false).unwrap();
1219    ///
1220    ///     assert_eq!(pyarray_i.readonly().as_slice().unwrap(), &[2, 3, 4]);
1221    /// });
1222    /// ```
1223    ///
1224    /// [PyArray_CastToType]: https://numpy.org/doc/stable/reference/c-api/array.html#c.PyArray_CastToType
1225    fn cast<U: Element>(&self, is_fortran: bool) -> PyResult<Bound<'py, PyArray<U, D>>>
1226    where
1227        T: Element;
1228
1229    /// A view of `self` with a different order of axes determined by `axes`.
1230    ///
1231    /// If `axes` is `None`, the order of axes is reversed which corresponds to the standard matrix transpose.
1232    ///
1233    /// See also [`numpy.transpose`][numpy-transpose] and [`PyArray_Transpose`][PyArray_Transpose].
1234    ///
1235    /// # Example
1236    ///
1237    /// ```
1238    /// use numpy::prelude::*;
1239    /// use numpy::PyArray;
1240    /// use pyo3::Python;
1241    /// use ndarray::array;
1242    ///
1243    /// Python::with_gil(|py| {
1244    ///     let array = array![[0, 1, 2], [3, 4, 5]].into_pyarray(py);
1245    ///
1246    ///     let array = array.permute(Some([1, 0])).unwrap();
1247    ///
1248    ///     assert_eq!(array.readonly().as_array(), array![[0, 3], [1, 4], [2, 5]]);
1249    /// });
1250    /// ```
1251    ///
1252    /// [numpy-transpose]: https://numpy.org/doc/stable/reference/generated/numpy.transpose.html
1253    /// [PyArray_Transpose]: https://numpy.org/doc/stable/reference/c-api/array.html#c.PyArray_Transpose
1254    fn permute<ID: IntoDimension>(&self, axes: Option<ID>) -> PyResult<Bound<'py, PyArray<T, D>>>
1255    where
1256        T: Element;
1257
1258    /// Special case of [`permute`][Self::permute] which reverses the order the axes.
1259    fn transpose(&self) -> PyResult<Bound<'py, PyArray<T, D>>>
1260    where
1261        T: Element,
1262    {
1263        self.permute::<()>(None)
1264    }
1265
1266    /// Construct a new array which has same values as `self`,
1267    /// but has different dimensions specified by `shape`
1268    /// and a possibly different memory order specified by `order`.
1269    ///
1270    /// See also [`numpy.reshape`][numpy-reshape] and [`PyArray_Newshape`][PyArray_Newshape].
1271    ///
1272    /// # Example
1273    ///
1274    /// ```
1275    /// use numpy::prelude::*;
1276    /// use numpy::{npyffi::NPY_ORDER, PyArray};
1277    /// use pyo3::Python;
1278    /// use ndarray::array;
1279    ///
1280    /// Python::with_gil(|py| {
1281    ///     let array =
1282    ///         PyArray::from_iter(py, 0..9).reshape_with_order([3, 3], NPY_ORDER::NPY_FORTRANORDER).unwrap();
1283    ///
1284    ///     assert_eq!(array.readonly().as_array(), array![[0, 3, 6], [1, 4, 7], [2, 5, 8]]);
1285    ///     assert!(array.is_fortran_contiguous());
1286    ///
1287    ///     assert!(array.reshape([5]).is_err());
1288    /// });
1289    /// ```
1290    ///
1291    /// [numpy-reshape]: https://numpy.org/doc/stable/reference/generated/numpy.reshape.html
1292    /// [PyArray_Newshape]: https://numpy.org/doc/stable/reference/c-api/array.html#c.PyArray_Newshape
1293    fn reshape_with_order<ID: IntoDimension>(
1294        &self,
1295        shape: ID,
1296        order: NPY_ORDER,
1297    ) -> PyResult<Bound<'py, PyArray<T, ID::Dim>>>
1298    where
1299        T: Element;
1300
1301    /// Special case of [`reshape_with_order`][Self::reshape_with_order] which keeps the memory order the same.
1302    #[inline(always)]
1303    fn reshape<ID: IntoDimension>(&self, shape: ID) -> PyResult<Bound<'py, PyArray<T, ID::Dim>>>
1304    where
1305        T: Element,
1306    {
1307        self.reshape_with_order(shape, NPY_ORDER::NPY_ANYORDER)
1308    }
1309
1310    /// Extends or truncates the dimensions of an array.
1311    ///
1312    /// This method works only on [contiguous][PyUntypedArrayMethods::is_contiguous] arrays.
1313    /// Missing elements will be initialized as if calling [`zeros`][PyArray::zeros].
1314    ///
1315    /// See also [`ndarray.resize`][ndarray-resize] and [`PyArray_Resize`][PyArray_Resize].
1316    ///
1317    /// # Safety
1318    ///
1319    /// There should be no outstanding references (shared or exclusive) into the array
1320    /// as this method might re-allocate it and thereby invalidate all pointers into it.
1321    ///
1322    /// # Example
1323    ///
1324    /// ```
1325    /// use numpy::prelude::*;
1326    /// use numpy::PyArray;
1327    /// use pyo3::Python;
1328    ///
1329    /// Python::with_gil(|py| {
1330    ///     let pyarray = PyArray::<f64, _>::zeros(py, (10, 10), false);
1331    ///     assert_eq!(pyarray.shape(), [10, 10]);
1332    ///
1333    ///     unsafe {
1334    ///         pyarray.resize((100, 100)).unwrap();
1335    ///     }
1336    ///     assert_eq!(pyarray.shape(), [100, 100]);
1337    /// });
1338    /// ```
1339    ///
1340    /// [ndarray-resize]: https://numpy.org/doc/stable/reference/generated/numpy.ndarray.resize.html
1341    /// [PyArray_Resize]: https://numpy.org/doc/stable/reference/c-api/array.html#c.PyArray_Resize
1342    unsafe fn resize<ID: IntoDimension>(&self, newshape: ID) -> PyResult<()>
1343    where
1344        T: Element;
1345
1346    /// Try to convert this array into a [`nalgebra::MatrixView`] using the given shape and strides.
1347    ///
1348    /// # Safety
1349    ///
1350    /// Calling this method invalidates all exclusive references to the internal data, e.g. `ArrayViewMut` or `MatrixSliceMut`.
1351    #[doc(alias = "nalgebra")]
1352    #[cfg(feature = "nalgebra")]
1353    unsafe fn try_as_matrix<R, C, RStride, CStride>(
1354        &self,
1355    ) -> Option<nalgebra::MatrixView<'_, T, R, C, RStride, CStride>>
1356    where
1357        T: nalgebra::Scalar + Element,
1358        D: Dimension,
1359        R: nalgebra::Dim,
1360        C: nalgebra::Dim,
1361        RStride: nalgebra::Dim,
1362        CStride: nalgebra::Dim;
1363
1364    /// Try to convert this array into a [`nalgebra::MatrixViewMut`] using the given shape and strides.
1365    ///
1366    /// # Safety
1367    ///
1368    /// Calling this method invalidates all other references to the internal data, e.g. `ArrayView`, `MatrixSlice`, `ArrayViewMut` or `MatrixSliceMut`.
1369    #[doc(alias = "nalgebra")]
1370    #[cfg(feature = "nalgebra")]
1371    unsafe fn try_as_matrix_mut<R, C, RStride, CStride>(
1372        &self,
1373    ) -> Option<nalgebra::MatrixViewMut<'_, T, R, C, RStride, CStride>>
1374    where
1375        T: nalgebra::Scalar + Element,
1376        D: Dimension,
1377        R: nalgebra::Dim,
1378        C: nalgebra::Dim,
1379        RStride: nalgebra::Dim,
1380        CStride: nalgebra::Dim;
1381}
1382
1383/// Implementation of functionality for [`PyArray0<T>`].
1384#[doc(alias = "PyArray", alias = "PyArray0")]
1385pub trait PyArray0Methods<'py, T>: PyArrayMethods<'py, T, Ix0> {
1386    /// Get the single element of a zero-dimensional array.
1387    ///
1388    /// See [`inner`][crate::inner] for an example.
1389    fn item(&self) -> T
1390    where
1391        T: Element + Copy,
1392    {
1393        unsafe { *self.data() }
1394    }
1395}
1396
1397#[inline(always)]
1398fn get_raw<T, D, Idx>(slf: &Bound<'_, PyArray<T, D>>, index: Idx) -> Option<*mut T>
1399where
1400    T: Element,
1401    D: Dimension,
1402    Idx: NpyIndex<Dim = D>,
1403{
1404    let offset = index.get_checked::<T>(slf.shape(), slf.strides())?;
1405    Some(unsafe { slf.data().offset(offset) })
1406}
1407
1408fn as_view<T, D, S, F>(slf: &Bound<'_, PyArray<T, D>>, from_shape_ptr: F) -> ArrayBase<S, D>
1409where
1410    T: Element,
1411    D: Dimension,
1412    S: RawData,
1413    F: FnOnce(StrideShape<D>, *mut T) -> ArrayBase<S, D>,
1414{
1415    fn inner<D: Dimension>(
1416        shape: &[usize],
1417        strides: &[isize],
1418        itemsize: usize,
1419        mut data_ptr: *mut u8,
1420    ) -> (StrideShape<D>, u32, *mut u8) {
1421        let shape = D::from_dimension(&Dim(shape)).expect(DIMENSIONALITY_MISMATCH_ERR);
1422
1423        assert!(strides.len() <= 32, "{}", MAX_DIMENSIONALITY_ERR);
1424
1425        let mut new_strides = D::zeros(strides.len());
1426        let mut inverted_axes = 0_u32;
1427
1428        for i in 0..strides.len() {
1429            // FIXME(kngwyu): Replace this hacky negative strides support with
1430            // a proper constructor, when it's implemented.
1431            // See https://github.com/rust-ndarray/ndarray/issues/842 for more.
1432            if strides[i] >= 0 {
1433                new_strides[i] = strides[i] as usize / itemsize;
1434            } else {
1435                // Move the pointer to the start position.
1436                data_ptr = unsafe { data_ptr.offset(strides[i] * (shape[i] as isize - 1)) };
1437
1438                new_strides[i] = (-strides[i]) as usize / itemsize;
1439                inverted_axes |= 1 << i;
1440            }
1441        }
1442
1443        (shape.strides(new_strides), inverted_axes, data_ptr)
1444    }
1445
1446    let (shape, mut inverted_axes, data_ptr) = inner(
1447        slf.shape(),
1448        slf.strides(),
1449        mem::size_of::<T>(),
1450        slf.data() as _,
1451    );
1452
1453    let mut array = from_shape_ptr(shape, data_ptr as _);
1454
1455    while inverted_axes != 0 {
1456        let axis = inverted_axes.trailing_zeros() as usize;
1457        inverted_axes &= !(1 << axis);
1458
1459        array.invert_axis(Axis(axis));
1460    }
1461
1462    array
1463}
1464
1465#[cfg(feature = "nalgebra")]
1466fn try_as_matrix_shape_strides<N, D, R, C, RStride, CStride>(
1467    slf: &Bound<'_, PyArray<N, D>>,
1468) -> Option<((R, C), (RStride, CStride))>
1469where
1470    N: nalgebra::Scalar + Element,
1471    D: Dimension,
1472    R: nalgebra::Dim,
1473    C: nalgebra::Dim,
1474    RStride: nalgebra::Dim,
1475    CStride: nalgebra::Dim,
1476{
1477    let ndim = slf.ndim();
1478    let shape = slf.shape();
1479    let strides = slf.strides();
1480
1481    if ndim != 1 && ndim != 2 {
1482        return None;
1483    }
1484
1485    if strides.iter().any(|strides| *strides < 0) {
1486        return None;
1487    }
1488
1489    let rows = shape[0];
1490    let cols = *shape.get(1).unwrap_or(&1);
1491
1492    if R::try_to_usize().map(|expected| rows == expected) == Some(false) {
1493        return None;
1494    }
1495
1496    if C::try_to_usize().map(|expected| cols == expected) == Some(false) {
1497        return None;
1498    }
1499
1500    let row_stride = strides[0] as usize / mem::size_of::<N>();
1501    let col_stride = strides
1502        .get(1)
1503        .map_or(rows, |stride| *stride as usize / mem::size_of::<N>());
1504
1505    if RStride::try_to_usize().map(|expected| row_stride == expected) == Some(false) {
1506        return None;
1507    }
1508
1509    if CStride::try_to_usize().map(|expected| col_stride == expected) == Some(false) {
1510        return None;
1511    }
1512
1513    let shape = (R::from_usize(rows), C::from_usize(cols));
1514
1515    let strides = (
1516        RStride::from_usize(row_stride),
1517        CStride::from_usize(col_stride),
1518    );
1519
1520    Some((shape, strides))
1521}
1522
1523impl<'py, T, D> PyArrayMethods<'py, T, D> for Bound<'py, PyArray<T, D>> {
1524    #[inline(always)]
1525    fn as_untyped(&self) -> &Bound<'py, PyUntypedArray> {
1526        unsafe { self.downcast_unchecked() }
1527    }
1528
1529    #[inline(always)]
1530    fn data(&self) -> *mut T {
1531        unsafe { (*self.as_array_ptr()).data.cast() }
1532    }
1533
1534    #[inline(always)]
1535    unsafe fn get(&self, index: impl NpyIndex<Dim = D>) -> Option<&T>
1536    where
1537        T: Element,
1538        D: Dimension,
1539    {
1540        let ptr = get_raw(self, index)?;
1541        Some(&*ptr)
1542    }
1543
1544    #[inline(always)]
1545    unsafe fn get_mut(&self, index: impl NpyIndex<Dim = D>) -> Option<&mut T>
1546    where
1547        T: Element,
1548        D: Dimension,
1549    {
1550        let ptr = get_raw(self, index)?;
1551        Some(&mut *ptr)
1552    }
1553
1554    fn get_owned<Idx>(&self, index: Idx) -> Option<T>
1555    where
1556        T: Element,
1557        D: Dimension,
1558        Idx: NpyIndex<Dim = D>,
1559    {
1560        let element = unsafe { self.get(index) };
1561        element.map(|elem| elem.clone_ref(self.py()))
1562    }
1563
1564    fn to_dyn(&self) -> &Bound<'py, PyArray<T, IxDyn>> {
1565        unsafe { self.downcast_unchecked() }
1566    }
1567
1568    fn to_vec(&self) -> Result<Vec<T>, NotContiguousError>
1569    where
1570        T: Element,
1571        D: Dimension,
1572    {
1573        let slice = unsafe { self.as_slice() };
1574        slice.map(|slc| T::vec_from_slice(self.py(), slc))
1575    }
1576
1577    fn try_readonly(&self) -> Result<PyReadonlyArray<'py, T, D>, BorrowError>
1578    where
1579        T: Element,
1580        D: Dimension,
1581    {
1582        PyReadonlyArray::try_new(self.clone())
1583    }
1584
1585    fn try_readwrite(&self) -> Result<PyReadwriteArray<'py, T, D>, BorrowError>
1586    where
1587        T: Element,
1588        D: Dimension,
1589    {
1590        PyReadwriteArray::try_new(self.clone())
1591    }
1592
1593    unsafe fn as_array(&self) -> ArrayView<'_, T, D>
1594    where
1595        T: Element,
1596        D: Dimension,
1597    {
1598        as_view(self, |shape, ptr| ArrayView::from_shape_ptr(shape, ptr))
1599    }
1600
1601    unsafe fn as_array_mut(&self) -> ArrayViewMut<'_, T, D>
1602    where
1603        T: Element,
1604        D: Dimension,
1605    {
1606        as_view(self, |shape, ptr| ArrayViewMut::from_shape_ptr(shape, ptr))
1607    }
1608
1609    fn as_raw_array(&self) -> RawArrayView<T, D>
1610    where
1611        T: Element,
1612        D: Dimension,
1613    {
1614        as_view(self, |shape, ptr| unsafe {
1615            RawArrayView::from_shape_ptr(shape, ptr)
1616        })
1617    }
1618
1619    fn as_raw_array_mut(&self) -> RawArrayViewMut<T, D>
1620    where
1621        T: Element,
1622        D: Dimension,
1623    {
1624        as_view(self, |shape, ptr| unsafe {
1625            RawArrayViewMut::from_shape_ptr(shape, ptr)
1626        })
1627    }
1628
1629    fn to_owned_array(&self) -> Array<T, D>
1630    where
1631        T: Element,
1632        D: Dimension,
1633    {
1634        let view = unsafe { self.as_array() };
1635        T::array_from_view(self.py(), view)
1636    }
1637
1638    fn copy_to<U: Element>(&self, other: &Bound<'py, PyArray<U, D>>) -> PyResult<()>
1639    where
1640        T: Element,
1641    {
1642        let self_ptr = self.as_array_ptr();
1643        let other_ptr = other.as_array_ptr();
1644        let result = unsafe { PY_ARRAY_API.PyArray_CopyInto(self.py(), other_ptr, self_ptr) };
1645        if result != -1 {
1646            Ok(())
1647        } else {
1648            Err(PyErr::fetch(self.py()))
1649        }
1650    }
1651
1652    fn cast<U: Element>(&self, is_fortran: bool) -> PyResult<Bound<'py, PyArray<U, D>>>
1653    where
1654        T: Element,
1655    {
1656        let ptr = unsafe {
1657            PY_ARRAY_API.PyArray_CastToType(
1658                self.py(),
1659                self.as_array_ptr(),
1660                U::get_dtype(self.py()).into_dtype_ptr(),
1661                if is_fortran { -1 } else { 0 },
1662            )
1663        };
1664        unsafe {
1665            Bound::from_owned_ptr_or_err(self.py(), ptr).map(|ob| ob.downcast_into_unchecked())
1666        }
1667    }
1668
1669    fn permute<ID: IntoDimension>(&self, axes: Option<ID>) -> PyResult<Bound<'py, PyArray<T, D>>> {
1670        let mut axes = axes.map(|axes| axes.into_dimension());
1671        let mut axes = axes.as_mut().map(|axes| axes.to_npy_dims());
1672        let axes = axes
1673            .as_mut()
1674            .map_or_else(ptr::null_mut, |axes| axes as *mut npyffi::PyArray_Dims);
1675
1676        let py = self.py();
1677        let ptr = unsafe { PY_ARRAY_API.PyArray_Transpose(py, self.as_array_ptr(), axes) };
1678        unsafe { Bound::from_owned_ptr_or_err(py, ptr).map(|ob| ob.downcast_into_unchecked()) }
1679    }
1680
1681    fn reshape_with_order<ID: IntoDimension>(
1682        &self,
1683        shape: ID,
1684        order: NPY_ORDER,
1685    ) -> PyResult<Bound<'py, PyArray<T, ID::Dim>>>
1686    where
1687        T: Element,
1688    {
1689        let mut shape = shape.into_dimension();
1690        let mut shape = shape.to_npy_dims();
1691
1692        let py = self.py();
1693        let ptr = unsafe {
1694            PY_ARRAY_API.PyArray_Newshape(
1695                py,
1696                self.as_array_ptr(),
1697                &mut shape as *mut npyffi::PyArray_Dims,
1698                order,
1699            )
1700        };
1701        unsafe { Bound::from_owned_ptr_or_err(py, ptr).map(|ob| ob.downcast_into_unchecked()) }
1702    }
1703
1704    unsafe fn resize<ID: IntoDimension>(&self, newshape: ID) -> PyResult<()>
1705    where
1706        T: Element,
1707    {
1708        let mut newshape = newshape.into_dimension();
1709        let mut newshape = newshape.to_npy_dims();
1710
1711        let py = self.py();
1712        let res = PY_ARRAY_API.PyArray_Resize(
1713            py,
1714            self.as_array_ptr(),
1715            &mut newshape as *mut npyffi::PyArray_Dims,
1716            1,
1717            NPY_ORDER::NPY_ANYORDER,
1718        );
1719
1720        if !res.is_null() {
1721            Ok(())
1722        } else {
1723            Err(PyErr::fetch(py))
1724        }
1725    }
1726
1727    #[cfg(feature = "nalgebra")]
1728    unsafe fn try_as_matrix<R, C, RStride, CStride>(
1729        &self,
1730    ) -> Option<nalgebra::MatrixView<'_, T, R, C, RStride, CStride>>
1731    where
1732        T: nalgebra::Scalar + Element,
1733        D: Dimension,
1734        R: nalgebra::Dim,
1735        C: nalgebra::Dim,
1736        RStride: nalgebra::Dim,
1737        CStride: nalgebra::Dim,
1738    {
1739        let (shape, strides) = try_as_matrix_shape_strides(self)?;
1740
1741        let storage = nalgebra::ViewStorage::from_raw_parts(self.data(), shape, strides);
1742
1743        Some(nalgebra::Matrix::from_data(storage))
1744    }
1745
1746    #[cfg(feature = "nalgebra")]
1747    unsafe fn try_as_matrix_mut<R, C, RStride, CStride>(
1748        &self,
1749    ) -> Option<nalgebra::MatrixViewMut<'_, T, R, C, RStride, CStride>>
1750    where
1751        T: nalgebra::Scalar + Element,
1752        D: Dimension,
1753        R: nalgebra::Dim,
1754        C: nalgebra::Dim,
1755        RStride: nalgebra::Dim,
1756        CStride: nalgebra::Dim,
1757    {
1758        let (shape, strides) = try_as_matrix_shape_strides(self)?;
1759
1760        let storage = nalgebra::ViewStorageMut::from_raw_parts(self.data(), shape, strides);
1761
1762        Some(nalgebra::Matrix::from_data(storage))
1763    }
1764}
1765
1766impl<'py, T> PyArray0Methods<'py, T> for Bound<'py, PyArray0<T>> {}
1767
1768#[cfg(test)]
1769mod tests {
1770    use super::*;
1771
1772    use ndarray::array;
1773    use pyo3::{py_run, types::PyList};
1774
1775    #[test]
1776    fn test_dyn_to_owned_array() {
1777        Python::with_gil(|py| {
1778            let array = PyArray::from_vec2(py, &[vec![1, 2], vec![3, 4]])
1779                .unwrap()
1780                .to_dyn()
1781                .to_owned_array();
1782
1783            assert_eq!(array, array![[1, 2], [3, 4]].into_dyn());
1784        });
1785    }
1786
1787    #[test]
1788    fn test_hasobject_flag() {
1789        Python::with_gil(|py| {
1790            let array: Bound<'_, PyArray<PyObject, _>> =
1791                PyArray1::from_slice(py, &[PyList::empty(py).into()]);
1792
1793            py_run!(py, array, "assert array.dtype.hasobject");
1794        });
1795    }
1796}