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