faer_ext/
lib.rs

1#![cfg_attr(docsrs, feature(doc_cfg))]
2
3/// Conversions from external library matrix views into `faer` types.
4pub trait IntoFaer {
5    type Faer;
6    fn into_faer(self) -> Self::Faer;
7}
8
9#[cfg(feature = "nalgebra")]
10#[cfg_attr(docsrs, doc(cfg(feature = "nalgebra")))]
11/// Conversions from external library matrix views into `nalgebra` types.
12pub trait IntoNalgebra {
13    type Nalgebra;
14    fn into_nalgebra(self) -> Self::Nalgebra;
15}
16
17#[cfg(feature = "ndarray")]
18#[cfg_attr(docsrs, doc(cfg(feature = "ndarray")))]
19/// Conversions from external library matrix views into `ndarray` types.
20pub trait IntoNdarray {
21    type Ndarray;
22    fn into_ndarray(self) -> Self::Ndarray;
23}
24
25#[cfg(feature = "nalgebra")]
26#[cfg_attr(docsrs, doc(cfg(feature = "nalgebra")))]
27const _: () = {
28    use faer::prelude::*;
29    use nalgebra::{Dim, Dyn, MatrixView, MatrixViewMut, ViewStorage, ViewStorageMut};
30
31    impl<'a, T, R: Dim, C: Dim, RStride: Dim, CStride: Dim> IntoFaer
32        for MatrixView<'a, T, R, C, RStride, CStride>
33    {
34        type Faer = MatRef<'a, T>;
35
36        #[track_caller]
37        fn into_faer(self) -> Self::Faer {
38            let nrows = self.nrows();
39            let ncols = self.ncols();
40            let strides = self.strides();
41            let ptr = self.as_ptr();
42            unsafe {
43                faer::MatRef::from_raw_parts(
44                    ptr,
45                    nrows,
46                    ncols,
47                    strides.0.try_into().unwrap(),
48                    strides.1.try_into().unwrap(),
49                )
50            }
51        }
52    }
53
54    impl<'a, T, R: Dim, C: Dim, RStride: Dim, CStride: Dim> IntoFaer
55        for MatrixViewMut<'a, T, R, C, RStride, CStride>
56    {
57        type Faer = MatMut<'a, T>;
58
59        #[track_caller]
60        fn into_faer(self) -> Self::Faer {
61            let nrows = self.nrows();
62            let ncols = self.ncols();
63            let strides = self.strides();
64            let ptr = { self }.as_mut_ptr();
65            unsafe {
66                faer::MatMut::from_raw_parts_mut(
67                    ptr,
68                    nrows,
69                    ncols,
70                    strides.0.try_into().unwrap(),
71                    strides.1.try_into().unwrap(),
72                )
73            }
74        }
75    }
76
77    impl<'a, T> IntoNalgebra for MatRef<'a, T> {
78        type Nalgebra = MatrixView<'a, T, Dyn, Dyn, Dyn, Dyn>;
79
80        #[track_caller]
81        fn into_nalgebra(self) -> Self::Nalgebra {
82            let nrows = self.nrows();
83            let ncols = self.ncols();
84            let row_stride = self.row_stride();
85            let col_stride = self.col_stride();
86            let ptr = self.as_ptr();
87            unsafe {
88                MatrixView::<'_, T, Dyn, Dyn, Dyn, Dyn>::from_data(ViewStorage::<
89                    '_,
90                    T,
91                    Dyn,
92                    Dyn,
93                    Dyn,
94                    Dyn,
95                >::from_raw_parts(
96                    ptr,
97                    (Dyn(nrows), Dyn(ncols)),
98                    (
99                        Dyn(row_stride.try_into().unwrap()),
100                        Dyn(col_stride.try_into().unwrap()),
101                    ),
102                ))
103            }
104        }
105    }
106
107    impl<'a, T> IntoNalgebra for MatMut<'a, T> {
108        type Nalgebra = MatrixViewMut<'a, T, Dyn, Dyn, Dyn, Dyn>;
109
110        #[track_caller]
111        fn into_nalgebra(self) -> Self::Nalgebra {
112            let nrows = self.nrows();
113            let ncols = self.ncols();
114            let row_stride = self.row_stride();
115            let col_stride = self.col_stride();
116            let ptr = self.as_ptr_mut();
117            unsafe {
118                MatrixViewMut::<'_, T, Dyn, Dyn, Dyn, Dyn>::from_data(ViewStorageMut::<
119                    '_,
120                    T,
121                    Dyn,
122                    Dyn,
123                    Dyn,
124                    Dyn,
125                >::from_raw_parts(
126                    ptr,
127                    (Dyn(nrows), Dyn(ncols)),
128                    (
129                        Dyn(row_stride.try_into().unwrap()),
130                        Dyn(col_stride.try_into().unwrap()),
131                    ),
132                ))
133            }
134        }
135    }
136};
137
138#[cfg(feature = "ndarray")]
139#[cfg_attr(docsrs, doc(cfg(feature = "ndarray")))]
140const _: () = {
141    use faer::prelude::*;
142    use ndarray::{ArrayView, ArrayViewMut, Ix2, ShapeBuilder};
143
144    impl<'a, T> IntoFaer for ArrayView<'a, T, Ix2> {
145        type Faer = MatRef<'a, T>;
146
147        #[track_caller]
148        fn into_faer(self) -> Self::Faer {
149            let nrows = self.nrows();
150            let ncols = self.ncols();
151            let strides: [isize; 2] = self.strides().try_into().unwrap();
152            let ptr = self.as_ptr();
153            unsafe { faer::MatRef::from_raw_parts(ptr, nrows, ncols, strides[0], strides[1]) }
154        }
155    }
156
157    impl<'a, T> IntoFaer for ArrayViewMut<'a, T, Ix2> {
158        type Faer = MatMut<'a, T>;
159
160        #[track_caller]
161        fn into_faer(self) -> Self::Faer {
162            let nrows = self.nrows();
163            let ncols = self.ncols();
164            let strides: [isize; 2] = self.strides().try_into().unwrap();
165            let ptr = { self }.as_mut_ptr();
166            unsafe { faer::MatMut::from_raw_parts_mut(ptr, nrows, ncols, strides[0], strides[1]) }
167        }
168    }
169
170    impl<'a, T> IntoNdarray for MatRef<'a, T> {
171        type Ndarray = ArrayView<'a, T, Ix2>;
172
173        #[track_caller]
174        fn into_ndarray(self) -> Self::Ndarray {
175            let nrows = self.nrows();
176            let ncols = self.ncols();
177            let row_stride: usize = self.row_stride().try_into().unwrap();
178            let col_stride: usize = self.col_stride().try_into().unwrap();
179            let ptr = self.as_ptr();
180            unsafe {
181                ArrayView::<'_, T, Ix2>::from_shape_ptr(
182                    (nrows, ncols).strides((row_stride, col_stride)),
183                    ptr,
184                )
185            }
186        }
187    }
188
189    impl<'a, T> IntoNdarray for MatMut<'a, T> {
190        type Ndarray = ArrayViewMut<'a, T, Ix2>;
191
192        #[track_caller]
193        fn into_ndarray(self) -> Self::Ndarray {
194            let nrows = self.nrows();
195            let ncols = self.ncols();
196            let row_stride: usize = self.row_stride().try_into().unwrap();
197            let col_stride: usize = self.col_stride().try_into().unwrap();
198            let ptr = self.as_ptr_mut();
199            unsafe {
200                ArrayViewMut::<'_, T, Ix2>::from_shape_ptr(
201                    (nrows, ncols).strides((row_stride, col_stride)),
202                    ptr,
203                )
204            }
205        }
206    }
207};
208
209#[cfg(all(feature = "nalgebra", feature = "ndarray"))]
210#[cfg_attr(docsrs, doc(cfg(all(feature = "nalgebra", feature = "ndarray"))))]
211const _: () = {
212    use nalgebra::{Dim, Dyn, MatrixView, MatrixViewMut, ViewStorage, ViewStorageMut};
213    use ndarray::{ArrayView, ArrayViewMut, Ix2, ShapeBuilder};
214
215    impl<'a, T> IntoNalgebra for ArrayView<'a, T, Ix2> {
216        type Nalgebra = MatrixView<'a, T, Dyn, Dyn, Dyn, Dyn>;
217
218        #[track_caller]
219        fn into_nalgebra(self) -> Self::Nalgebra {
220            let nrows = self.nrows();
221            let ncols = self.ncols();
222            let [row_stride, col_stride]: [isize; 2] = self.strides().try_into().unwrap();
223            let ptr = self.as_ptr();
224
225            unsafe {
226                MatrixView::<'_, T, Dyn, Dyn, Dyn, Dyn>::from_data(ViewStorage::<
227                    '_,
228                    T,
229                    Dyn,
230                    Dyn,
231                    Dyn,
232                    Dyn,
233                >::from_raw_parts(
234                    ptr,
235                    (Dyn(nrows), Dyn(ncols)),
236                    (
237                        Dyn(row_stride.try_into().unwrap()),
238                        Dyn(col_stride.try_into().unwrap()),
239                    ),
240                ))
241            }
242        }
243    }
244    impl<'a, T> IntoNalgebra for ArrayViewMut<'a, T, Ix2> {
245        type Nalgebra = MatrixViewMut<'a, T, Dyn, Dyn, Dyn, Dyn>;
246
247        #[track_caller]
248        fn into_nalgebra(self) -> Self::Nalgebra {
249            let nrows = self.nrows();
250            let ncols = self.ncols();
251            let [row_stride, col_stride]: [isize; 2] = self.strides().try_into().unwrap();
252            let ptr = { self }.as_mut_ptr();
253
254            unsafe {
255                MatrixViewMut::<'_, T, Dyn, Dyn, Dyn, Dyn>::from_data(ViewStorageMut::<
256                    '_,
257                    T,
258                    Dyn,
259                    Dyn,
260                    Dyn,
261                    Dyn,
262                >::from_raw_parts(
263                    ptr,
264                    (Dyn(nrows), Dyn(ncols)),
265                    (
266                        Dyn(row_stride.try_into().unwrap()),
267                        Dyn(col_stride.try_into().unwrap()),
268                    ),
269                ))
270            }
271        }
272    }
273
274    impl<'a, T, R: Dim, C: Dim, RStride: Dim, CStride: Dim> IntoNdarray
275        for MatrixView<'a, T, R, C, RStride, CStride>
276    {
277        type Ndarray = ArrayView<'a, T, Ix2>;
278
279        #[track_caller]
280        fn into_ndarray(self) -> Self::Ndarray {
281            let nrows = self.nrows();
282            let ncols = self.ncols();
283            let (row_stride, col_stride) = self.strides();
284            let ptr = self.as_ptr();
285
286            unsafe {
287                ArrayView::<'_, T, Ix2>::from_shape_ptr(
288                    (nrows, ncols).strides((row_stride, col_stride)),
289                    ptr,
290                )
291            }
292        }
293    }
294    impl<'a, T, R: Dim, C: Dim, RStride: Dim, CStride: Dim> IntoNdarray
295        for MatrixViewMut<'a, T, R, C, RStride, CStride>
296    {
297        type Ndarray = ArrayViewMut<'a, T, Ix2>;
298
299        #[track_caller]
300        fn into_ndarray(self) -> Self::Ndarray {
301            let nrows = self.nrows();
302            let ncols = self.ncols();
303            let (row_stride, col_stride) = self.strides();
304            let ptr = { self }.as_mut_ptr();
305
306            unsafe {
307                ArrayViewMut::<'_, T, Ix2>::from_shape_ptr(
308                    (nrows, ncols).strides((row_stride, col_stride)),
309                    ptr,
310                )
311            }
312        }
313    }
314};
315
316#[cfg(feature = "polars")]
317#[cfg_attr(docsrs, doc(cfg(feature = "polars")))]
318pub mod polars {
319    use faer::Mat;
320    use polars::prelude::*;
321
322    pub trait Frame {
323        fn is_valid(self) -> PolarsResult<LazyFrame>;
324    }
325
326    impl Frame for LazyFrame {
327        fn is_valid(self) -> PolarsResult<LazyFrame> {
328            let test_dtypes: bool = self
329                .clone()
330                .limit(0)
331                .collect()?
332                .dtypes()
333                .into_iter()
334                .map(|e| {
335                    matches!(
336                        e,
337                        DataType::UInt8
338                            | DataType::UInt16
339                            | DataType::UInt32
340                            | DataType::UInt64
341                            | DataType::Int8
342                            | DataType::Int16
343                            | DataType::Int32
344                            | DataType::Int64
345                            | DataType::Float32
346                            | DataType::Float64
347                    )
348                })
349                .all(|e| e);
350            let test_no_nulls: bool = self
351                .clone()
352                .null_count()
353                .cast_all(DataType::UInt64, true)
354                .with_column(
355                    fold_exprs(
356                        lit(0).cast(DataType::UInt64),
357                        |acc, x| Ok(Some(acc + x)),
358                        [col("*")],
359                    )
360                    .alias("sum"),
361                )
362                .select(&[col("sum")])
363                .collect()?
364                .column("sum")?
365                .u64()?
366                .into_iter()
367                .map(|e| e.eq(&Some(0u64)))
368                .collect::<Vec<_>>()[0];
369            match (test_dtypes, test_no_nulls) {
370                (true, true) => Ok(self),
371                (false, true) => Err(PolarsError::InvalidOperation(
372                    "frame contains non-numerical data".into(),
373                )),
374                (true, false) => Err(PolarsError::InvalidOperation(
375                    "frame contains null entries".into(),
376                )),
377                (false, false) => Err(PolarsError::InvalidOperation(
378                    "frame contains non-numerical data and null entries".into(),
379                )),
380            }
381        }
382    }
383
384    macro_rules! polars_impl {
385        ($ty: ident, $dtype: ident, $fn_name: ident) => {
386            /// Converts a `polars` lazyframe into a [`Mat`].
387            ///
388            /// Note that this function expects that the frame passed "looks like"
389            /// a numerical array and all values will be cast to either f32 or f64
390            /// prior to building [`Mat`].
391            ///
392            /// Passing a frame with either non-numerical column data or null
393            /// entries will result in a error. Users are expected to reolve
394            /// these issues in `polars` prior calling this function.
395            #[cfg(feature = "polars")]
396            #[cfg_attr(docsrs, doc(cfg(feature = "polars")))]
397            pub fn $fn_name(
398                frame: impl Frame,
399            ) -> PolarsResult<Mat<$ty>> {
400                use core::{iter::zip, mem::MaybeUninit};
401
402                fn implementation(
403                    lf: LazyFrame,
404                ) -> PolarsResult<Mat<$ty>> {
405                    let df = lf
406                        .select(&[col("*").cast(DataType::$dtype)])
407                        .collect()?;
408
409                    let nrows = df.height();
410                    let ncols = df.get_column_names().len();
411
412                    let mut out = Mat::<$ty>::with_capacity(df.height(), df.get_column_names().len());
413                    let stride = out.col_stride() as usize;
414                    let out_ptr = out.as_ptr_mut();
415
416                    df.get_column_names().iter()
417                        .enumerate()
418                        .try_for_each(|(j, col)| -> PolarsResult<()> {
419                            let mut row_start = 0usize;
420
421                            // SAFETY: this is safe since we allocated enough space for `ncols` columns and
422                            // `nrows` rows
423                            let out_col = unsafe {
424                                core::slice::from_raw_parts_mut(
425                                    out_ptr.wrapping_add(stride * j) as *mut MaybeUninit<$ty>,
426                                    nrows,
427                                )
428                            };
429
430                            df.column(col)?.$ty()?.downcast_iter().try_for_each(
431                                |chunk| -> PolarsResult<()> {
432                                    let len = chunk.len();
433                                    if len == 0 {
434                                        return Ok(());
435                                    }
436
437                                    match row_start.checked_add(len) {
438                                        Some(next_row_start) => {
439                                            if next_row_start <= nrows {
440                                                let mut out_slice = &mut out_col[row_start..next_row_start];
441                                                let mut values = chunk.values_iter().as_slice();
442                                                let validity = chunk.validity();
443
444                                                assert_eq!(values.len(), len);
445
446                                                match validity {
447                                                    Some(bitmap) => {
448                                                        let (mut bytes, offset, bitmap_len) = bitmap.as_slice();
449                                                        assert_eq!(bitmap_len, len);
450                                                        const BITS_PER_BYTE: usize = 8;
451
452                                                        if offset > 0 {
453                                                            let first_byte_len = Ord::min(len, 8 - offset);
454
455                                                            let (out_prefix, out_suffix) = out_slice.split_at_mut(first_byte_len);
456                                                            let (values_prefix, values_suffix) = values.split_at(first_byte_len);
457
458                                                            for (out_elem, value_elem) in zip(
459                                                                out_prefix,
460                                                                values_prefix,
461                                                            ) {
462                                                                *out_elem = MaybeUninit::new(*value_elem)
463                                                            }
464
465                                                            bytes = &bytes[1..];
466                                                            values = values_suffix;
467                                                            out_slice = out_suffix;
468                                                        }
469
470                                                        if bytes.len() > 0 {
471                                                            for (out_slice8, values8) in zip(
472                                                                out_slice.chunks_exact_mut(BITS_PER_BYTE),
473                                                                values.chunks_exact(BITS_PER_BYTE),
474                                                            ) {
475                                                                for (out_elem, value_elem) in zip(out_slice8, values8) {
476                                                                    *out_elem = MaybeUninit::new(*value_elem);
477                                                                }
478                                                            }
479
480                                                            for (out_elem, value_elem) in zip(
481                                                                out_slice.chunks_exact_mut(BITS_PER_BYTE).into_remainder(),
482                                                                values.chunks_exact(BITS_PER_BYTE).remainder(),
483                                                            ) {
484                                                                *out_elem = MaybeUninit::new(*value_elem);
485                                                            }
486                                                        }
487                                                    }
488                                                    None => {
489                                                        // SAFETY: T and MaybeUninit<T> have the same layout
490                                                        // NOTE: This state should not be reachable
491                                                        let values = unsafe {
492                                                            core::slice::from_raw_parts(
493                                                                values.as_ptr() as *const MaybeUninit<$ty>,
494                                                                values.len(),
495                                                            )
496                                                        };
497                                                        out_slice.copy_from_slice(values);
498                                                    }
499                                                }
500
501                                                row_start = next_row_start;
502                                                Ok(())
503                                            } else {
504                                                Err(PolarsError::ShapeMismatch(
505                                                    format!("too many values in column {col}").into(),
506                                                ))
507                                            }
508                                        }
509                                        None => Err(PolarsError::ShapeMismatch(
510                                            format!("too many values in column {col}").into(),
511                                        )),
512                                    }
513                                },
514                            )?;
515
516                            if row_start < nrows {
517                                Err(PolarsError::ShapeMismatch(
518                                    format!("not enough values in column {col} (column has {row_start} values, while dataframe has {nrows} rows)").into(),
519                                ))
520                            } else {
521                                Ok(())
522                            }
523                        })?;
524
525                    // SAFETY: we initialized every `ncols` columns, and each one was initialized with `nrows`
526                    // elements
527                    unsafe { out.set_dims(nrows, ncols) };
528
529                    Ok(out)
530                }
531
532                implementation(frame.is_valid()?)
533            }
534        };
535    }
536
537    polars_impl!(f32, Float32, polars_to_faer_f32);
538    polars_impl!(f64, Float64, polars_to_faer_f64);
539}
540
541#[cfg(feature = "numpy")]
542#[cfg_attr(docsrs, doc(cfg(feature = "numpy")))]
543const _: () = {
544    use faer::prelude::*;
545    use numpy::{PyReadonlyArray1, PyReadonlyArray2, PyArrayMethods};
546    use numpy::Element;
547
548    impl<'a, T: Element + 'a> IntoFaer for PyReadonlyArray2<'a, T> {
549        type Faer = MatRef<'a, T>;
550
551        #[track_caller]
552        fn into_faer(self) -> Self::Faer {
553            let raw_arr = self.as_raw_array();
554            let nrows = raw_arr.nrows();
555            let ncols = raw_arr.ncols();
556            let strides: [isize; 2] = raw_arr.strides().try_into().unwrap();
557            let ptr = raw_arr.as_ptr();
558            unsafe { MatRef::from_raw_parts(ptr, nrows, ncols, strides[0], strides[1]) }
559        }
560    }
561
562    impl<'a, T: Element + 'a> IntoFaer for PyReadonlyArray1<'a, T> {
563        type Faer = ColRef<'a, T>;
564
565        #[track_caller]
566        fn into_faer(self) -> Self::Faer {
567            let raw_arr = self.as_raw_array();
568            let nrows = raw_arr.len();
569            let strides: [isize; 1] = raw_arr.strides().try_into().unwrap();
570            let ptr = raw_arr.as_ptr();
571            unsafe { ColRef::from_raw_parts(ptr, nrows, strides[0]) }
572        }
573    }
574
575};
576
577#[cfg(test)]
578mod tests {
579    #![allow(unused_imports)]
580    #![allow(non_snake_case)]
581
582    use super::*;
583    use faer::mat;
584    use faer::prelude::*;
585
586    #[cfg(feature = "ndarray")]
587    #[test]
588    fn test_ext_ndarray() {
589        let mut I_faer = Mat::<f32>::identity(8, 7);
590        let mut I_ndarray = ndarray::Array2::<f32>::zeros([8, 7]);
591        I_ndarray.diag_mut().fill(1.0);
592
593        assert_eq!(I_ndarray.view().into_faer(), I_faer);
594        assert_eq!(I_faer.as_ref().into_ndarray(), I_ndarray);
595
596        assert_eq!(I_ndarray.view_mut().into_faer(), I_faer);
597        assert_eq!(I_faer.as_mut().into_ndarray(), I_ndarray);
598    }
599
600    #[cfg(feature = "nalgebra")]
601    #[test]
602    fn test_ext_nalgebra() {
603        let mut I_faer = Mat::<f32>::identity(8, 7);
604        let mut I_nalgebra = nalgebra::DMatrix::<f32>::identity(8, 7);
605
606        assert_eq!(I_nalgebra.view_range(.., ..).into_faer(), I_faer);
607        assert_eq!(I_faer.as_ref().into_nalgebra(), I_nalgebra);
608
609        assert_eq!(I_nalgebra.view_range_mut(.., ..).into_faer(), I_faer);
610        assert_eq!(I_faer.as_mut().into_nalgebra(), I_nalgebra);
611    }
612
613    #[cfg(feature = "polars")]
614    #[test]
615    fn test_polars_pos() {
616        use crate::polars::{polars_to_faer_f32, polars_to_faer_f64};
617        #[rustfmt::skip]
618        use ::polars::prelude::*;
619
620        let s0: Series = Series::new("a", [1, 2, 3]);
621        let s1: Series = Series::new("b", [10, 11, 12]);
622
623        let lf = DataFrame::new(vec![s0, s1]).unwrap().lazy();
624
625        let arr_32 = polars_to_faer_f32(lf.clone()).unwrap();
626        let arr_64 = polars_to_faer_f64(lf).unwrap();
627
628        let expected_32 = mat![[1f32, 10f32], [2f32, 11f32], [3f32, 12f32]];
629        let expected_64 = mat![[1f64, 10f64], [2f64, 11f64], [3f64, 12f64]];
630
631        assert_eq!(arr_32, expected_32);
632        assert_eq!(arr_64, expected_64);
633    }
634
635    #[cfg(feature = "polars")]
636    #[test]
637    #[should_panic(expected = "frame contains null entries")]
638    fn test_polars_neg_32_null() {
639        use crate::polars::polars_to_faer_f32;
640        #[rustfmt::skip]
641        use ::polars::prelude::*;
642
643        let s0: Series = Series::new("a", [1, 2, 3]);
644        let s1: Series = Series::new("b", [Some(10), Some(11), None]);
645
646        let lf = DataFrame::new(vec![s0, s1]).unwrap().lazy();
647
648        polars_to_faer_f32(lf).unwrap();
649    }
650
651    #[cfg(feature = "polars")]
652    #[test]
653    #[should_panic(expected = "frame contains non-numerical data")]
654    fn test_polars_neg_32_strl() {
655        use crate::polars::polars_to_faer_f32;
656        #[rustfmt::skip]
657        use ::polars::prelude::*;
658
659        let s0: Series = Series::new("a", [1, 2, 3]);
660        let s1: Series = Series::new("b", ["fish", "dog", "crocodile"]);
661
662        let lf = DataFrame::new(vec![s0, s1]).unwrap().lazy();
663
664        polars_to_faer_f32(lf).unwrap();
665    }
666
667    #[cfg(feature = "polars")]
668    #[test]
669    #[should_panic(expected = "frame contains non-numerical data and null entries")]
670    fn test_polars_neg_32_combo() {
671        use crate::polars::polars_to_faer_f32;
672        #[rustfmt::skip]
673        use ::polars::prelude::*;
674
675        let s0: Series = Series::new("a", [1, 2, 3]);
676        let s1: Series = Series::new("b", [Some(10), Some(11), None]);
677        let s2: Series = Series::new("c", [Some("fish"), Some("dog"), None]);
678
679        let lf = DataFrame::new(vec![s0, s1, s2]).unwrap().lazy();
680
681        polars_to_faer_f32(lf).unwrap();
682    }
683
684    #[cfg(feature = "polars")]
685    #[test]
686    #[should_panic(expected = "frame contains null entries")]
687    fn test_polars_neg_64_null() {
688        use crate::polars::polars_to_faer_f64;
689        #[rustfmt::skip]
690        use ::polars::prelude::*;
691
692        let s0: Series = Series::new("a", [1, 2, 3]);
693        let s1: Series = Series::new("b", [Some(10), Some(11), None]);
694
695        let lf = DataFrame::new(vec![s0, s1]).unwrap().lazy();
696
697        polars_to_faer_f64(lf).unwrap();
698    }
699
700    #[cfg(feature = "polars")]
701    #[test]
702    #[should_panic(expected = "frame contains non-numerical data")]
703    fn test_polars_neg_64_strl() {
704        use crate::polars::polars_to_faer_f64;
705        #[rustfmt::skip]
706        use ::polars::prelude::*;
707
708        let s0: Series = Series::new("a", [1, 2, 3]);
709        let s1: Series = Series::new("b", ["fish", "dog", "crocodile"]);
710
711        let lf = DataFrame::new(vec![s0, s1]).unwrap().lazy();
712
713        polars_to_faer_f64(lf).unwrap();
714    }
715
716    #[cfg(feature = "polars")]
717    #[test]
718    #[should_panic(expected = "frame contains non-numerical data and null entries")]
719    fn test_polars_neg_64_combo() {
720        use crate::polars::polars_to_faer_f64;
721        #[rustfmt::skip]
722        use ::polars::prelude::*;
723
724        let s0: Series = Series::new("a", [1, 2, 3]);
725        let s1: Series = Series::new("b", [Some(10), Some(11), None]);
726        let s2: Series = Series::new("c", [Some("fish"), Some("dog"), None]);
727
728        let lf = DataFrame::new(vec![s0, s1, s2]).unwrap().lazy();
729
730        polars_to_faer_f64(lf).unwrap();
731    }
732
733    #[cfg(feature = "numpy")]
734    #[test]
735    fn test_ext_numpy() {
736        use pyo3::{Bound, Python};
737        use numpy::{PyArray1, PyReadonlyArray1, PyReadonlyArray2, PyArrayMethods};
738        use numpy::ndarray::array;
739        use numpy::pyarray;
740        use pyo3::prelude::*;
741        pyo3::prepare_freethreaded_python();
742
743        Python::with_gil(|py| {
744            let arr_1 = PyArray1::from_vec(py, vec![1., 0., 1., 0.]);
745            let py_array1: PyReadonlyArray1<f64> = arr_1.readonly();
746            let expected_f64: Mat<f64> = mat![[1., 0., 1., 0.]];
747            assert_eq!(py_array1.into_faer(), expected_f64.transpose().col(0));
748
749            let arr_2 = pyarray![py, [1., 0.], [0., 1.]];
750            let py_array2: PyReadonlyArray2<f64> = arr_2.readonly();
751            let expected_f64 = Mat::<f64>::identity(2, 2);
752            assert_eq!(py_array2.into_faer(), expected_f64.as_ref());
753
754            let arr_1 = PyArray1::from_vec(py, vec![1., 0., 1., 0.]);
755            let py_array1: PyReadonlyArray1<f32> = arr_1.readonly();
756            let expected_f32: Mat<f32> = mat![[1., 0., 1., 0.]];
757            assert_eq!(py_array1.into_faer(), expected_f32.transpose().col(0));
758
759            let arr_2 = pyarray![py, [1., 0.], [0., 1.]];
760            let py_array2: PyReadonlyArray2<f32> = arr_2.readonly();
761            let expected_f32 = Mat::<f32>::identity(2, 2);
762            assert_eq!(py_array2.into_faer(), expected_f32.as_ref());
763        });
764    }
765}