faer/
lib.rs

1//! `faer` is a general-purpose linear algebra library for rust, with a focus on high performance
2//! for algebraic operations on medium/large matrices, as well as matrix decompositions
3//!
4//! most of the high-level functionality in this library is provided through associated functions in
5//! its vocabulary types: [`Mat`]/[`MatRef`]/[`MatMut`]
6//!
7//! `faer` is recommended for applications that handle medium to large dense matrices, and its
8//! design is not well suited for applications that operate mostly on low dimensional vectors and
9//! matrices such as computer graphics or game development. for such applications, `nalgebra` and
10//! `cgmath` may be better suited
11//!
12//! # basic usage
13//!
14//! [`Mat`] is a resizable matrix type with dynamic capacity, which can be created using
15//! [`Mat::new`] to produce an empty $0\times 0$ matrix, [`Mat::zeros`] to create a rectangular
16//! matrix filled with zeros, [`Mat::identity`] to create an identity matrix, or [`Mat::from_fn`]
17//! for the more general case
18//!
19//! Given a `&Mat<T>` (resp. `&mut Mat<T>`), a [`MatRef<'_, T>`](MatRef) (resp. [`MatMut<'_,
20//! T>`](MatMut)) can be created by calling [`Mat::as_ref`] (resp. [`Mat::as_mut`]), which allow
21//! for more flexibility than `Mat` in that they allow slicing ([`MatRef::get`]) and splitting
22//! ([`MatRef::split_at`])
23//!
24//! `MatRef` and `MatMut` are lightweight view objects. the former can be copied freely while the
25//! latter has move and reborrow semantics, as described in its documentation
26//!
27//! most of the matrix operations can be used through the corresponding math operators: `+` for
28//! matrix addition, `-` for subtraction, `*` for either scalar or matrix multiplication depending
29//! on the types of the operands.
30//!
31//! ## example
32//! ```
33//! use faer::{Mat, Scale, mat};
34//!
35//! let a = mat![
36//! 	[1.0, 5.0, 9.0], //
37//! 	[2.0, 6.0, 10.0],
38//! 	[3.0, 7.0, 11.0],
39//! 	[4.0, 8.0, 12.0f64],
40//! ];
41//!
42//! let b = Mat::from_fn(4, 3, |i, j| (i + j) as f64);
43//!
44//! let add = &a + &b;
45//! let sub = &a - &b;
46//! let scale = Scale(3.0) * &a;
47//! let mul = &a * b.transpose();
48//!
49//! let a00 = a[(0, 0)];
50//! ```
51//!
52//! # matrix decompositions
53//! `faer` provides a variety of matrix factorizations, each with its own advantages and drawbacks:
54//!
55//! ## $LL^\top$ decomposition
56//! [`Mat::llt`] decomposes a self-adjoint positive definite matrix $A$ such that
57//! $$A = LL^H,$$
58//! where $L$ is a lower triangular matrix. this decomposition is highly efficient and has good
59//! stability properties
60//!
61//! [an implementation for sparse matrices is also available](sparse::linalg::solvers::Llt)
62//!
63//! ## $LBL^\top$ decomposition
64//! [`Mat::lblt`] decomposes a self-adjoint (possibly indefinite) matrix $A$ such that
65//! $$P A P^\top = LBL^H,$$
66//! where $P$ is a permutation matrix, $L$ is a lower triangular matrix, and $B$ is a block
67//! diagonal matrix, with $1 \times 1$ or $2 \times 2$ diagonal blocks.
68//! this decomposition is efficient and has good stability properties
69//!
70//! ## $LU$ decomposition with partial pivoting
71//! [`Mat::partial_piv_lu`] decomposes a square invertible matrix $A$ into a lower triangular
72//! matrix $L$, a unit upper triangular matrix $U$, and a permutation matrix $P$, such that
73//! $$PA = LU$$
74//! it is used by default for computing the determinant, and is generally the recommended method
75//! for solving a square linear system or computing the inverse of a matrix (although we generally
76//! recommend using a [`faer::linalg::solvers::Solve`](crate::linalg::solvers::Solve) instead of
77//! computing the inverse explicitly)
78//!
79//! [an implementation for sparse matrices is also available](sparse::linalg::solvers::Lu)
80//!
81//! ## $LU$ decomposition with full pivoting
82//! [`Mat::full_piv_lu`] decomposes a generic rectangular matrix $A$ into a lower triangular
83//! matrix $L$, a unit upper triangular matrix $U$, and permutation matrices $P$ and $Q$, such that
84//! $$PAQ^\top = LU$$
85//! it can be more stable than the LU decomposition with partial pivoting, in exchange for being
86//! more computationally expensive
87//!
88//! ## $QR$ decomposition
89//! [`Mat::qr`] decomposes a matrix $A$ into the product $$A = QR,$$
90//! where $Q$ is a unitary matrix, and $R$ is an upper trapezoidal matrix. it is often used for
91//! solving least squares problems
92//!
93//! [an implementation for sparse matrices is also available](sparse::linalg::solvers::Qr)
94//!
95//! ## $QR$ decomposition with column pivoting
96//! ([`Mat::col_piv_qr`]) decomposes a matrix $A$ into the product $$AP^\top = QR,$$
97//! where $P$ is a permutation matrix, $Q$ is a unitary matrix, and $R$ is an upper trapezoidal
98//! matrix
99//!
100//! it is slower than the version with no pivoting, in exchange for being more numerically stable
101//! for rank-deficient matrices
102//!
103//! ## singular value decomposition
104//! the SVD of a matrix $A$ of shape $(m, n)$ is a decomposition into three components $U$, $S$,
105//! and $V$, such that:
106//!
107//! - $U$ has shape $(m, m)$ and is a unitary matrix,
108//! - $V$ has shape $(n, n)$ and is a unitary matrix,
109//! - $S$ has shape $(m, n)$ and is zero everywhere except the main diagonal, with nonnegative
110//! diagonal values in nonincreasing order,
111//! - and finally:
112//!
113//! $$A = U S V^H$$
114//!
115//! the SVD is provided in two forms: either the full matrices $U$ and $V$ are computed, using
116//! [`Mat::svd`], or only their first $\min(m, n)$ columns are computed, using
117//! [`Mat::thin_svd`]
118//!
119//! if only the singular values (elements of $S$) are desired, they can be obtained in
120//! nonincreasing order using [`Mat::singular_values`]
121//!
122//! ## eigendecomposition
123//! **note**: the order of the eigenvalues is currently unspecified and may be changed in a future
124//! release
125//!
126//! the eigenvalue decomposition of a square matrix $A$ of shape $(n, n)$ is a decomposition into
127//! two components $U$, $S$:
128//!
129//! - $U$ has shape $(n, n)$ and is invertible,
130//! - $S$ has shape $(n, n)$ and is a diagonal matrix,
131//! - and finally:
132//!
133//! $$A = U S U^{-1}$$
134//!
135//! if $A$ is self-adjoint, then $U$ can be made unitary ($U^{-1} = U^H$), and $S$ is real valued.
136//! additionally, the eigenvalues are sorted in nondecreasing order
137//!
138//! Depending on the domain of the input matrix and whether it is self-adjoint, multiple methods
139//! are provided to compute the eigendecomposition:
140//! * [`Mat::self_adjoint_eigen`] can be used with either real or complex matrices,
141//! producing an eigendecomposition of the same type,
142//! * [`Mat::eigen`] can be used with real or complex matrices, but always produces complex values.
143//!
144//! if only the eigenvalues (elements of $S$) are desired, they can be obtained using
145//! [`Mat::self_adjoint_eigenvalues`] (nondecreasing order), [`Mat::eigenvalues`]
146//! with the same conditions described above.
147//!
148//! # crate features
149//!
150//! - `std`: enabled by default. links with the standard library to enable additional features such
151//!   as cpu feature detection at runtime
152//! - `rayon`: enabled by default. enables the `rayon` parallel backend and enables global
153//!   parallelism by default
154//! - `serde`: Enables serialization and deserialization of [`Mat`]
155//! - `npy`: enables conversions to/from numpy's matrix file format
156//! - `perf-warn`: produces performance warnings when matrix operations are called with suboptimal
157//! data layout
158//! - `nightly`: requires the nightly compiler. enables experimental simd features such as avx512
159
160#![cfg_attr(not(feature = "std"), no_std)]
161#![allow(non_snake_case)]
162#![warn(missing_docs)]
163#![warn(rustdoc::broken_intra_doc_links)]
164
165extern crate alloc;
166#[cfg(feature = "std")]
167extern crate std;
168
169/// see: [`generativity::make_guard`]
170#[macro_export]
171macro_rules! make_guard {
172    ($($name:ident),* $(,)?) => {$(
173        #[allow(unused_unsafe)]
174        let $name = unsafe { extern crate generativity; ::generativity::Id::new() };
175        #[allow(unused, unused_unsafe)]
176        let lifetime_brand = unsafe { extern crate generativity; ::generativity::LifetimeBrand::new(&$name) };
177        #[allow(unused_unsafe)]
178        let $name = unsafe { extern crate generativity; ::generativity::Guard::new($name) };
179    )*};
180}
181
182macro_rules! repeat_n {
183	($e: expr, $n: expr) => {
184		iter::repeat_n($e, $n)
185	};
186}
187
188macro_rules! try_const {
189	($e: expr) => {
190		::pulp::try_const! { $e }
191	};
192}
193
194use core::num::NonZeroUsize;
195use core::sync::atomic::AtomicUsize;
196use equator::{assert, debug_assert};
197use faer_traits::*;
198
199macro_rules! auto {
200	($ty: ty) => {
201		$crate::Auto::<$ty>::auto()
202	};
203}
204
205macro_rules! dispatch {
206	($imp: expr, $ty: ident, $T: ty $(,)?) => {
207		if try_const! { <$T>::IS_NATIVE_C32 } {
208			unsafe { transmute(<ComplexImpl<f32> as ComplexField>::Arch::default().dispatch(transmute::<_, $ty<ComplexImpl<f32>>>($imp))) }
209		} else if try_const! { <$T>::IS_NATIVE_C64 } {
210			unsafe { transmute(<ComplexImpl<f64> as ComplexField>::Arch::default().dispatch(transmute::<_, $ty<ComplexImpl<f64>>>($imp))) }
211		} else {
212			<$T>::Arch::default().dispatch($imp)
213		}
214	};
215}
216
217macro_rules! stack_mat {
218	($name: ident, $m: expr, $n: expr, $A: expr, $N: expr, $T: ty $(,)?) => {
219		let mut __tmp = {
220			#[repr(align(64))]
221			struct __Col<T, const A: usize>([T; A]);
222			struct __Mat<T, const A: usize, const N: usize>([__Col<T, A>; N]);
223
224			core::mem::MaybeUninit::<__Mat<$T, $A, $N>>::uninit()
225		};
226		let __stack = MemStack::new_any(core::slice::from_mut(&mut __tmp));
227		let mut $name = $crate::linalg::temp_mat_zeroed::<$T, _, _>($m, $n, __stack).0;
228		let mut $name = $name.as_mat_mut();
229	};
230
231	($name: ident, $m: expr, $n: expr,  $T: ty $(,)?) => {
232		stack_mat!($name, $m, $n, $m, $n, $T)
233	};
234}
235
236#[macro_export]
237#[doc(hidden)]
238macro_rules! __dbg {
239    () => {
240        std::eprintln!("[{}:{}:{}]", std::file!(), std::line!(), std::column!())
241    };
242    ($val:expr $(,)?) => {
243        match $val {
244            tmp => {
245                std::eprintln!("[{}:{}:{}] {} = {:16.12?}",
246                    std::file!(), std::line!(), std::column!(), std::stringify!($val), &tmp);
247                tmp
248            }
249        }
250    };
251    ($($val:expr),+ $(,)?) => {
252        ($($crate::__dbg!($val)),+,)
253    };
254}
255
256#[cfg(feature = "perf-warn")]
257#[macro_export]
258#[doc(hidden)]
259macro_rules! __perf_warn {
260	($name: ident) => {{
261		#[inline(always)]
262		#[allow(non_snake_case)]
263		fn $name() -> &'static ::core::sync::atomic::AtomicBool {
264			static $name: ::core::sync::atomic::AtomicBool = ::core::sync::atomic::AtomicBool::new(false);
265			&$name
266		}
267		::core::matches!(
268			$name().compare_exchange(
269				false,
270				true,
271				::core::sync::atomic::Ordering::Relaxed,
272				::core::sync::atomic::Ordering::Relaxed,
273			),
274			Ok(_)
275		)
276	}};
277}
278
279#[doc(hidden)]
280#[macro_export]
281macro_rules! with_dim {
282	($name: ident, $value: expr $(,)?) => {
283		let __val__ = $value;
284		$crate::make_guard!($name);
285		let $name = $crate::utils::bound::Dim::new(__val__, $name);
286	};
287
288	({$(let $name: ident = $value: expr;)*}) => {$(
289		let __val__ = $value;
290		$crate::make_guard!($name);
291		let $name = $crate::utils::bound::Dim::new(__val__, $name);
292	)*};
293}
294
295/// zips together matrix of the same size, so that coefficient-wise operations can be performed on
296/// their elements.
297///
298/// # note
299/// the order in which the matrix elements are traversed is unspecified.
300///
301/// # example
302/// ```
303/// use faer::{Mat, mat, unzip, zip};
304///
305/// let nrows = 2;
306/// let ncols = 3;
307///
308/// let a = mat![[1.0, 3.0, 5.0], [2.0, 4.0, 6.0]];
309/// let b = mat![[7.0, 9.0, 11.0], [8.0, 10.0, 12.0]];
310/// let mut sum = Mat::<f64>::zeros(nrows, ncols);
311///
312/// zip!(&mut sum, &a, &b).for_each(|unzip!(sum, a, b)| {
313/// 	*sum = a + b;
314/// });
315///
316/// for i in 0..nrows {
317/// 	for j in 0..ncols {
318/// 		assert_eq!(sum[(i, j)], a[(i, j)] + b[(i, j)]);
319/// 	}
320/// }
321/// ```
322#[macro_export]
323macro_rules! zip {
324    ($head: expr $(,)?) => {
325        $crate::linalg::zip::LastEq($crate::linalg::zip::IntoView::into_view($head), ::core::marker::PhantomData)
326    };
327
328    ($head: expr, $($tail: expr),* $(,)?) => {
329        $crate::linalg::zip::ZipEq::new($crate::linalg::zip::IntoView::into_view($head), $crate::zip!($($tail,)*))
330    };
331}
332
333/// used to undo the zipping by the [`zip!`] macro.
334///
335/// # example
336/// ```
337/// use faer::{Mat, mat, unzip, zip};
338///
339/// let nrows = 2;
340/// let ncols = 3;
341///
342/// let a = mat![[1.0, 3.0, 5.0], [2.0, 4.0, 6.0]];
343/// let b = mat![[7.0, 9.0, 11.0], [8.0, 10.0, 12.0]];
344/// let mut sum = Mat::<f64>::zeros(nrows, ncols);
345///
346/// zip!(&mut sum, &a, &b).for_each(|unzip!(sum, a, b)| {
347/// 	*sum = a + b;
348/// });
349///
350/// for i in 0..nrows {
351/// 	for j in 0..ncols {
352/// 		assert_eq!(sum[(i, j)], a[(i, j)] + b[(i, j)]);
353/// 	}
354/// }
355/// ```
356#[macro_export]
357macro_rules! unzip {
358    ($head: pat $(,)?) => {
359        $crate::linalg::zip::Last($head)
360    };
361
362    ($head: pat, $($tail: pat),* $(,)?) => {
363        $crate::linalg::zip::Zip($head, $crate::unzip!($($tail,)*))
364    };
365}
366
367#[macro_export]
368#[doc(hidden)]
369macro_rules! __transpose_impl {
370    ([$([$($col:expr),*])*] $($v:expr;)* ) => {
371        [$([$($col,)*],)* [$($v,)*]]
372    };
373    ([$([$($col:expr),*])*] $($v0:expr, $($v:expr),* ;)*) => {
374        $crate::__transpose_impl!([$([$($col),*])* [$($v0),*]] $($($v),* ;)*)
375    };
376}
377
378/// creates a [`Mat`] containing the arguments.
379///
380/// ```
381/// use faer::mat;
382///
383/// let matrix = mat![
384/// 	[1.0, 5.0, 9.0], //
385/// 	[2.0, 6.0, 10.0],
386/// 	[3.0, 7.0, 11.0],
387/// 	[4.0, 8.0, 12.0f64],
388/// ];
389///
390/// assert_eq!(matrix[(0, 0)], 1.0);
391/// assert_eq!(matrix[(1, 0)], 2.0);
392/// assert_eq!(matrix[(2, 0)], 3.0);
393/// assert_eq!(matrix[(3, 0)], 4.0);
394///
395/// assert_eq!(matrix[(0, 1)], 5.0);
396/// assert_eq!(matrix[(1, 1)], 6.0);
397/// assert_eq!(matrix[(2, 1)], 7.0);
398/// assert_eq!(matrix[(3, 1)], 8.0);
399///
400/// assert_eq!(matrix[(0, 2)], 9.0);
401/// assert_eq!(matrix[(1, 2)], 10.0);
402/// assert_eq!(matrix[(2, 2)], 11.0);
403/// assert_eq!(matrix[(3, 2)], 12.0);
404/// ```
405#[macro_export]
406macro_rules! mat {
407    () => {
408        {
409            compile_error!("number of columns in the matrix is ambiguous");
410        }
411    };
412
413    ($([$($v:expr),* $(,)?] ),* $(,)?) => {
414        {
415            let __data = ::core::mem::ManuallyDrop::new($crate::__transpose_impl!([] $($($v),* ;)*));
416            let __data = &*__data;
417            let __ncols = __data.len();
418            let __nrows = (*__data.get(0).unwrap()).len();
419
420            #[allow(unused_unsafe)]
421            unsafe {
422                $crate::mat::Mat::from_fn(__nrows, __ncols, |i, j| ::core::ptr::from_ref(&__data[j][i]).read())
423            }
424        }
425    };
426}
427
428/// creates a [`col::Col`] containing the arguments
429///
430/// ```
431/// use faer::col;
432///
433/// let col_vec = col![3.0, 5.0, 7.0, 9.0];
434///
435/// assert_eq!(col_vec[0], 3.0);
436/// assert_eq!(col_vec[1], 5.0);
437/// assert_eq!(col_vec[2], 7.0);
438/// assert_eq!(col_vec[3], 9.0);
439/// ```
440#[macro_export]
441macro_rules! col {
442    ($($v: expr),* $(,)?) => {
443        {
444            let __data = ::core::mem::ManuallyDrop::new([$($v,)*]);
445            let __data = &*__data;
446            let __len = __data.len();
447
448            #[allow(unused_unsafe)]
449            unsafe {
450                $crate::col::Col::from_fn(__len, |i| ::core::ptr::from_ref(&__data[i]).read())
451            }
452        }
453    };
454}
455
456/// creates a [`row::Row`] containing the arguments
457///
458/// ```
459/// use faer::row;
460///
461/// let row_vec = row![3.0, 5.0, 7.0, 9.0];
462///
463/// assert_eq!(row_vec[0], 3.0);
464/// assert_eq!(row_vec[1], 5.0);
465/// assert_eq!(row_vec[2], 7.0);
466/// assert_eq!(row_vec[3], 9.0);
467/// ```
468#[macro_export]
469macro_rules! row {
470    ($($v: expr),* $(,)?) => {
471        {
472            let __data = ::core::mem::ManuallyDrop::new([$($v,)*]);
473            let __data = &*__data;
474            let __len = __data.len();
475
476            #[allow(unused_unsafe)]
477            unsafe {
478                $crate::row::Row::from_fn(__len, |i| ::core::ptr::from_ref(&__data[i]).read())
479            }
480        }
481    };
482}
483
484/// convenience function to concatenate a nested list of matrices into a single
485/// big ['Mat']. concatonation pattern follows the numpy.block convention that
486/// each sub-list must have an equal number of columns (net) but the boundaries
487/// do not need to align. in other words, this sort of thing:
488/// ```notcode
489///   AAAbb
490///   AAAbb
491///   cDDDD
492/// ```
493/// is perfectly acceptable
494#[doc(hidden)]
495#[track_caller]
496pub fn concat_impl<T: ComplexField>(blocks: &[&[(mat::MatRef<'_, T>, Conj)]]) -> mat::Mat<T> {
497	#[inline(always)]
498	fn count_total_columns<T: ComplexField>(block_row: &[(mat::MatRef<'_, T>, Conj)]) -> usize {
499		let mut out: usize = 0;
500		for (elem, _) in block_row.iter() {
501			out += elem.ncols();
502		}
503		out
504	}
505
506	#[inline(always)]
507	#[track_caller]
508	fn count_rows<T: ComplexField>(block_row: &[(mat::MatRef<'_, T>, Conj)]) -> usize {
509		let mut out: usize = 0;
510		for (i, (e, _)) in block_row.iter().enumerate() {
511			if i == 0 {
512				out = e.nrows();
513			} else {
514				assert!(e.nrows() == out);
515			}
516		}
517		out
518	}
519
520	// get size of result while doing checks
521	let mut n: usize = 0;
522	let mut m: usize = 0;
523	for row in blocks.iter() {
524		n += count_rows(row);
525	}
526	for (i, row) in blocks.iter().enumerate() {
527		let cols = count_total_columns(row);
528		if i == 0 {
529			m = cols;
530		} else {
531			assert!(cols == m);
532		}
533	}
534
535	let mut mat = mat::Mat::<T>::zeros(n, m);
536	let mut ni: usize = 0;
537	let mut mj: usize;
538	for row in blocks.iter() {
539		mj = 0;
540
541		for (elem, conj) in row.iter() {
542			let mut dst = mat.as_mut().submatrix_mut(ni, mj, elem.nrows(), elem.ncols());
543			if *conj == Conj::No {
544				dst.copy_from(elem);
545			} else {
546				dst.copy_from(elem.conjugate());
547			}
548			mj += elem.ncols();
549		}
550		ni += row[0].0.nrows();
551	}
552
553	mat
554}
555
556/// concatenates the matrices in each row horizontally,
557/// then concatenates the results vertically
558///
559/// `concat![[a0, a1, a2], [b1, b2]]` results in the matrix
560/// ```notcode
561/// [a0 | a1 | a2][b0 | b1]
562/// ```
563#[macro_export]
564macro_rules! concat {
565    () => {
566        {
567            compile_error!("number of columns in the matrix is ambiguous");
568        }
569    };
570
571    ($([$($v:expr),* $(,)?] ),* $(,)?) => {
572        {
573            $crate::concat_impl(&[$(&[$(($v).as_ref().__canonicalize(),)*],)*])
574        }
575    };
576}
577
578/// helper utilities
579pub mod utils;
580
581/// diagonal matrix
582pub mod diag;
583/// rectangular matrix
584pub mod mat;
585/// permutation matrix
586pub mod perm;
587
588/// column vector
589pub mod col;
590/// row vector
591pub mod row;
592
593pub mod linalg;
594#[path = "./operator/mod.rs"]
595pub mod matrix_free;
596pub mod sparse;
597
598/// de-serialization from common matrix file formats
599#[cfg(feature = "std")]
600pub mod io;
601
602#[cfg(feature = "serde")]
603mod serde;
604
605/// native unsigned integer type
606pub trait Index: traits::IndexCore + traits::Index + seal::Seal {}
607impl<T: faer_traits::Index<Signed: seal::Seal> + seal::Seal> Index for T {}
608
609mod seal {
610	pub trait Seal {}
611	impl<T: faer_traits::Seal> Seal for T {}
612	impl Seal for crate::utils::bound::Dim<'_> {}
613	impl<I: crate::Index> Seal for crate::utils::bound::Idx<'_, I> {}
614	impl<I: crate::Index> Seal for crate::utils::bound::IdxInc<'_, I> {}
615	impl<I: crate::Index> Seal for crate::utils::bound::MaybeIdx<'_, I> {}
616	impl<I: crate::Index> Seal for crate::utils::bound::IdxIncOne<I> {}
617	impl<I: crate::Index> Seal for crate::utils::bound::MaybeIdxOne<I> {}
618	impl Seal for crate::utils::bound::One {}
619	impl Seal for crate::utils::bound::Zero {}
620	impl Seal for crate::ContiguousFwd {}
621	impl Seal for crate::ContiguousBwd {}
622}
623
624/// sealed trait for types that can be created from "unbound" values, as long as their
625/// struct preconditions are upheld
626pub trait Unbind<I = usize>: Send + Sync + Copy + core::fmt::Debug + seal::Seal {
627	/// creates new value
628	/// # safety
629	/// safety invariants must be upheld
630	unsafe fn new_unbound(idx: I) -> Self;
631
632	/// returns the unbound value, unconstrained by safety invariants
633	fn unbound(self) -> I;
634}
635
636/// type that can be used to index into a range
637pub type Idx<Dim, I = usize> = <Dim as ShapeIdx>::Idx<I>;
638/// type that can be used to partition a range
639pub type IdxInc<Dim, I = usize> = <Dim as ShapeIdx>::IdxInc<I>;
640/// either an index or a negative value
641pub type MaybeIdx<Dim, I = usize> = <Dim as ShapeIdx>::MaybeIdx<I>;
642
643/// base trait for [`Shape`]
644pub trait ShapeIdx {
645	/// type that can be used to index into a range
646	type Idx<I: Index>: Unbind<I> + Ord + Eq;
647	/// type that can be used to partition a range
648	type IdxInc<I: Index>: Unbind<I> + Ord + Eq + From<Idx<Self, I>>;
649	/// either an index or a negative value
650	type MaybeIdx<I: Index>: Unbind<I::Signed> + Ord + Eq;
651}
652
653/// matrix dimension
654pub trait Shape: Unbind + Ord + ShapeIdx<Idx<usize>: Ord + Eq + PartialOrd<Self>, IdxInc<usize>: Ord + Eq + PartialOrd<Self>> {
655	/// whether the types involved have any safety invariants
656	const IS_BOUND: bool = true;
657
658	/// bind the current value using a invariant lifetime guard
659	#[inline]
660	fn bind<'n>(self, guard: generativity::Guard<'n>) -> utils::bound::Dim<'n> {
661		utils::bound::Dim::new(self.unbound(), guard)
662	}
663
664	/// cast a slice of bound values to unbound values
665	#[inline]
666	fn cast_idx_slice<I: Index>(slice: &[Idx<Self, I>]) -> &[I] {
667		unsafe { core::slice::from_raw_parts(slice.as_ptr() as _, slice.len()) }
668	}
669
670	/// cast a slice of bound values to unbound values
671	#[inline]
672	fn cast_idx_inc_slice<I: Index>(slice: &[IdxInc<Self, I>]) -> &[I] {
673		unsafe { core::slice::from_raw_parts(slice.as_ptr() as _, slice.len()) }
674	}
675
676	/// returns the index `0`, which is always valid
677	#[inline(always)]
678	fn start() -> IdxInc<Self> {
679		unsafe { IdxInc::<Self>::new_unbound(0) }
680	}
681
682	/// returns the incremented value, as an inclusive index
683	#[inline(always)]
684	fn next(idx: Idx<Self>) -> IdxInc<Self> {
685		unsafe { IdxInc::<Self>::new_unbound(idx.unbound() + 1) }
686	}
687
688	/// returns the last value, equal to the dimension
689	#[inline(always)]
690	fn end(self) -> IdxInc<Self> {
691		unsafe { IdxInc::<Self>::new_unbound(self.unbound()) }
692	}
693
694	/// checks if the index is valid, returning `Some(_)` in that case
695	#[inline(always)]
696	fn idx(self, idx: usize) -> Option<Idx<Self>> {
697		if idx < self.unbound() {
698			Some(unsafe { Idx::<Self>::new_unbound(idx) })
699		} else {
700			None
701		}
702	}
703
704	/// checks if the index is valid, returning `Some(_)` in that case
705	#[inline(always)]
706	fn idx_inc(self, idx: usize) -> Option<IdxInc<Self>> {
707		if idx <= self.unbound() {
708			Some(unsafe { IdxInc::<Self>::new_unbound(idx) })
709		} else {
710			None
711		}
712	}
713
714	/// checks if the index is valid, and panics otherwise
715	#[inline(always)]
716	fn checked_idx(self, idx: usize) -> Idx<Self> {
717		equator::assert!(idx < self.unbound());
718		unsafe { Idx::<Self>::new_unbound(idx) }
719	}
720
721	/// checks if the index is valid, and panics otherwise
722	#[inline(always)]
723	fn checked_idx_inc(self, idx: usize) -> IdxInc<Self> {
724		equator::assert!(idx <= self.unbound());
725		unsafe { IdxInc::<Self>::new_unbound(idx) }
726	}
727
728	/// assumes the index is valid
729	/// # safety
730	/// the index must be valid
731	#[inline(always)]
732	unsafe fn unchecked_idx(self, idx: usize) -> Idx<Self> {
733		equator::debug_assert!(idx < self.unbound());
734		unsafe { Idx::<Self>::new_unbound(idx) }
735	}
736
737	/// assumes the index is valid
738	/// # safety
739	/// the index must be valid
740	#[inline(always)]
741	unsafe fn unchecked_idx_inc(self, idx: usize) -> IdxInc<Self> {
742		equator::debug_assert!(idx <= self.unbound());
743		unsafe { IdxInc::<Self>::new_unbound(idx) }
744	}
745
746	/// returns an iterator over the indices between `from` and `to`
747	#[inline(always)]
748	fn indices(from: IdxInc<Self>, to: IdxInc<Self>) -> impl Clone + ExactSizeIterator + DoubleEndedIterator<Item = Idx<Self>> {
749		(from.unbound()..to.unbound()).map(
750			#[inline(always)]
751			|i| unsafe { Idx::<Self>::new_unbound(i) },
752		)
753	}
754}
755
756impl<T: Send + Sync + Copy + core::fmt::Debug + faer_traits::Seal> Unbind<T> for T {
757	#[inline(always)]
758	unsafe fn new_unbound(idx: T) -> Self {
759		idx
760	}
761
762	#[inline(always)]
763	fn unbound(self) -> T {
764		self
765	}
766}
767
768impl ShapeIdx for usize {
769	type Idx<I: Index> = I;
770	type IdxInc<I: Index> = I;
771	type MaybeIdx<I: Index> = I::Signed;
772}
773impl Shape for usize {
774	const IS_BOUND: bool = false;
775}
776
777/// stride distance between two consecutive elements along a given dimension
778pub trait Stride: seal::Seal + core::fmt::Debug + Copy + Send + Sync + 'static {
779	/// the reversed stride type
780	type Rev: Stride<Rev = Self>;
781	/// returns the reversed stride
782	fn rev(self) -> Self::Rev;
783
784	/// returns the stride in elements
785	fn element_stride(self) -> isize;
786}
787
788impl Stride for isize {
789	type Rev = Self;
790
791	#[inline(always)]
792	fn rev(self) -> Self::Rev {
793		-self
794	}
795
796	#[inline(always)]
797	fn element_stride(self) -> isize {
798		self
799	}
800}
801
802/// contiguous stride equal to `+1`
803#[derive(Copy, Clone, Debug)]
804pub struct ContiguousFwd;
805/// contiguous stride equal to `-1`
806#[derive(Copy, Clone, Debug)]
807pub struct ContiguousBwd;
808
809impl Stride for ContiguousFwd {
810	type Rev = ContiguousBwd;
811
812	#[inline(always)]
813	fn rev(self) -> Self::Rev {
814		ContiguousBwd
815	}
816
817	#[inline(always)]
818	fn element_stride(self) -> isize {
819		1
820	}
821}
822
823impl Stride for ContiguousBwd {
824	type Rev = ContiguousFwd;
825
826	#[inline(always)]
827	fn rev(self) -> Self::Rev {
828		ContiguousFwd
829	}
830
831	#[inline(always)]
832	fn element_stride(self) -> isize {
833		-1
834	}
835}
836
837/// memory allocation error
838#[derive(Copy, Clone, Debug, PartialEq, Eq)]
839pub enum TryReserveError {
840	///rRequired allocation does not fit within `isize` bytes
841	CapacityOverflow,
842	/// allocator could not provide an allocation with the requested layout
843	AllocError {
844		/// requested layout
845		layout: core::alloc::Layout,
846	},
847}
848
849/// determines whether the input should be implicitly conjugated or not
850#[derive(Copy, Clone, Debug, PartialEq, Eq, PartialOrd, Ord)]
851pub enum Conj {
852	/// no implicit conjugation
853	No,
854	/// implicit conjugation
855	Yes,
856}
857
858#[derive(Copy, Clone, Debug, PartialEq, Eq, PartialOrd, Ord)]
859pub(crate) enum DiagStatus {
860	Unit,
861	Generic,
862}
863
864/// determines whether to replace or add to the result of a matmul operatio
865#[derive(Copy, Clone, Debug, PartialEq, Eq, PartialOrd, Ord)]
866pub enum Accum {
867	/// overwrites the output buffer
868	Replace,
869	/// adds the result to the output buffer
870	Add,
871}
872
873/// determines which side of a self-adjoint matrix should be accessed
874#[derive(Copy, Clone, Debug, PartialEq, Eq, PartialOrd, Ord)]
875pub enum Side {
876	/// lower triangular half
877	Lower,
878	/// upper triangular half
879	Upper,
880}
881
882impl Conj {
883	/// returns `self == Conj::Yes`
884	#[inline]
885	pub const fn is_conj(self) -> bool {
886		matches!(self, Conj::Yes)
887	}
888
889	/// returns the composition of `self` and `other`
890	#[inline]
891	pub const fn compose(self, other: Self) -> Self {
892		match (self, other) {
893			(Conj::No, Conj::No) => Conj::No,
894			(Conj::Yes, Conj::Yes) => Conj::No,
895			(Conj::No, Conj::Yes) => Conj::Yes,
896			(Conj::Yes, Conj::No) => Conj::Yes,
897		}
898	}
899
900	/// returns `Conj::No` if `T` is the canonical representation, otherwise `Conj::Yes`
901	#[inline]
902	pub const fn get<T: Conjugate>() -> Self {
903		if T::IS_CANONICAL { Self::No } else { Self::Yes }
904	}
905
906	#[inline]
907	pub(crate) fn apply<T: Conjugate>(value: &T) -> T::Canonical {
908		let value = unsafe { &*(value as *const T as *const T::Canonical) };
909
910		if try_const! { matches!(Self::get::<T>(), Conj::Yes) } {
911			T::Canonical::conj_impl(value)
912		} else {
913			T::Canonical::copy_impl(value)
914		}
915	}
916
917	#[inline]
918	pub(crate) fn apply_rt<T: ComplexField>(self, value: &T) -> T {
919		if self.is_conj() { T::conj_impl(value) } else { T::copy_impl(value) }
920	}
921}
922
923/// determines the parallelization configuration
924#[derive(Copy, Clone, Debug, PartialEq, Eq)]
925pub enum Par {
926	/// sequential, non portable across different platforms
927	Seq,
928	/// parallelized using the global rayon threadpool, non portable across different platforms
929	#[cfg(feature = "rayon")]
930	Rayon(NonZeroUsize),
931}
932
933impl Par {
934	/// returns `Par::Rayon(nthreads)` if `nthreads` is non-zero, or
935	/// `Par::Rayon(rayon::current_num_threads())` otherwise
936	#[inline]
937	#[cfg(feature = "rayon")]
938	pub fn rayon(nthreads: usize) -> Self {
939		if nthreads == 0 {
940			Self::Rayon(NonZeroUsize::new(rayon::current_num_threads()).unwrap())
941		} else {
942			Self::Rayon(NonZeroUsize::new(nthreads).unwrap())
943		}
944	}
945
946	/// the number of threads that should ideally execute an operation with the given parallelism
947	#[inline]
948	pub fn degree(&self) -> usize {
949		utils::thread::parallelism_degree(*self)
950	}
951}
952
953#[allow(non_camel_case_types)]
954/// `Complex<f32>`
955pub type c32 = traits::c32;
956#[allow(non_camel_case_types)]
957/// `Complex<f64>`
958pub type c64 = traits::c64;
959#[allow(non_camel_case_types)]
960/// `Complex<f64>`
961pub type cx128 = traits::cx128;
962#[allow(non_camel_case_types)]
963/// `Complex<f64>`
964pub type fx128 = traits::fx128;
965
966pub use col::{Col, ColMut, ColRef};
967pub use mat::{Mat, MatMut, MatRef};
968pub use row::{Row, RowMut, RowRef};
969
970#[allow(unused_imports, dead_code)]
971mod internal_prelude {
972	pub trait DivCeil: Sized {
973		fn msrv_div_ceil(self, rhs: Self) -> Self;
974		fn msrv_next_multiple_of(self, rhs: Self) -> Self;
975		fn msrv_checked_next_multiple_of(self, rhs: Self) -> Option<Self>;
976	}
977
978	impl DivCeil for usize {
979		#[inline]
980		fn msrv_div_ceil(self, rhs: Self) -> Self {
981			let d = self / rhs;
982			let r = self % rhs;
983			if r > 0 { d + 1 } else { d }
984		}
985
986		#[inline]
987		fn msrv_next_multiple_of(self, rhs: Self) -> Self {
988			match self % rhs {
989				0 => self,
990				r => self + (rhs - r),
991			}
992		}
993
994		#[inline]
995		fn msrv_checked_next_multiple_of(self, rhs: Self) -> Option<Self> {
996			{
997				match self.checked_rem(rhs)? {
998					0 => Some(self),
999					r => self.checked_add(rhs - r),
1000				}
1001			}
1002		}
1003	}
1004
1005	#[cfg(test)]
1006	pub(crate) use std::dbg;
1007	#[cfg(test)]
1008	pub(crate) use {alloc::boxed::Box, alloc::vec, alloc::vec::Vec};
1009
1010	pub use faer_traits::{ComplexImpl, ComplexImplConj, Symbolic};
1011
1012	pub(crate) use crate::col::{Col, ColMut, ColRef};
1013	pub(crate) use crate::diag::{Diag, DiagMut, DiagRef};
1014	pub(crate) use crate::hacks::transmute;
1015	pub(crate) use crate::linalg::{self, temp_mat_scratch, temp_mat_uninit, temp_mat_zeroed};
1016	pub(crate) use crate::mat::{AsMat, AsMatMut, AsMatRef, Mat, MatMut, MatRef};
1017	pub(crate) use crate::perm::{Perm, PermRef};
1018	pub(crate) use crate::prelude::*;
1019	pub(crate) use crate::row::{AsRowMut, AsRowRef, Row, RowMut, RowRef};
1020	pub(crate) use crate::utils::bound::{Array, Dim, Idx, IdxInc, MaybeIdx};
1021	pub(crate) use crate::utils::simd::SimdCtx;
1022	pub(crate) use crate::{Auto, NonExhaustive, Side, Spec};
1023
1024	pub use num_complex::Complex;
1025
1026	pub use faer_macros::math;
1027	pub use faer_traits::math_utils::*;
1028	pub use faer_traits::{ComplexField, Conjugate, Index, IndexCore, Real, RealField, SignedIndex, SimdArch};
1029
1030	#[inline]
1031	pub fn simd_align(i: usize) -> usize {
1032		i.wrapping_neg()
1033	}
1034
1035	pub(crate) use crate::{Accum, Conj, ContiguousBwd, ContiguousFwd, DiagStatus, Par, Shape, Stride, Unbind, unzip, zip};
1036
1037	pub use {unzip as uz, zip as z};
1038
1039	pub use crate::make_guard;
1040	pub use dyn_stack::{MemStack, StackReq};
1041	pub use equator::{assert, assert as Assert, debug_assert, debug_assert as DebugAssert};
1042	pub use reborrow::*;
1043}
1044
1045#[allow(unused_imports)]
1046pub(crate) mod internal_prelude_sp {
1047	pub(crate) use crate::internal_prelude::*;
1048	pub(crate) use crate::sparse::{
1049		FaerError, NONE, Pair, SparseColMat, SparseColMatMut, SparseColMatRef, SparseRowMat, SparseRowMatMut, SparseRowMatRef, SymbolicSparseColMat,
1050		SymbolicSparseColMatRef, SymbolicSparseRowMat, SymbolicSparseRowMatRef, Triplet, csc_numeric, csc_symbolic, csr_numeric, csr_symbolic,
1051		linalg as linalg_sp, try_collect, try_zeroed, windows2,
1052	};
1053	pub(crate) use core::cell::Cell;
1054	pub(crate) use core::iter;
1055	pub(crate) use dyn_stack::MemBuffer;
1056}
1057
1058/// useful imports for general usage of the library
1059pub mod prelude {
1060	pub use reborrow::{IntoConst, Reborrow, ReborrowMut};
1061
1062	pub use super::{Par, Scale, c32, c64, col, mat, row, unzip, zip};
1063	pub use col::{Col, ColMut, ColRef};
1064	pub use mat::{Mat, MatMut, MatRef};
1065	pub use row::{Row, RowMut, RowRef};
1066
1067	#[cfg(feature = "linalg")]
1068	pub use super::linalg::solvers::{DenseSolve, Solve, SolveLstsq};
1069
1070	#[cfg(feature = "sparse")]
1071	pub use super::prelude_sp::*;
1072
1073	/// see [`Default`]
1074	#[inline]
1075	pub fn default<T: Default>() -> T {
1076		Default::default()
1077	}
1078}
1079
1080#[cfg(feature = "sparse")]
1081mod prelude_sp {
1082	use crate::sparse;
1083
1084	pub use sparse::{SparseColMat, SparseColMatMut, SparseColMatRef, SparseRowMat, SparseRowMatMut, SparseRowMatRef};
1085}
1086
1087/// scaling factor for multiplying matrices.
1088#[derive(Copy, Clone, Debug)]
1089#[repr(transparent)]
1090pub struct Scale<T>(pub T);
1091impl<T> Scale<T> {
1092	/// create a reference to a scaling factor from a reference to a value.
1093	#[inline(always)]
1094	pub fn from_ref(value: &T) -> &Self {
1095		unsafe { &*(value as *const T as *const Self) }
1096	}
1097
1098	/// create a mutable reference to a scaling factor from a mutable reference to a value.
1099	#[inline(always)]
1100	pub fn from_mut(value: &mut T) -> &mut Self {
1101		unsafe { &mut *(value as *mut T as *mut Self) }
1102	}
1103}
1104
1105/// 0: disabled
1106/// 1: `Seq`
1107/// n >= 2: `Rayon(n - 2)`
1108///
1109/// default: `Rayon(0)`
1110static GLOBAL_PARALLELISM: AtomicUsize = {
1111	#[cfg(all(not(miri), feature = "rayon"))]
1112	{
1113		AtomicUsize::new(2)
1114	}
1115	#[cfg(not(all(not(miri), feature = "rayon")))]
1116	{
1117		AtomicUsize::new(1)
1118	}
1119};
1120
1121/// causes functions that access global parallelism settings to panic.
1122pub fn disable_global_parallelism() {
1123	GLOBAL_PARALLELISM.store(0, core::sync::atomic::Ordering::Relaxed);
1124}
1125
1126/// sets the global parallelism settings.
1127pub fn set_global_parallelism(par: Par) {
1128	let value = match par {
1129		Par::Seq => 1,
1130		#[cfg(feature = "rayon")]
1131		Par::Rayon(n) => n.get().saturating_add(2),
1132	};
1133	GLOBAL_PARALLELISM.store(value, core::sync::atomic::Ordering::Relaxed);
1134}
1135
1136/// gets the global parallelism settings.
1137///
1138/// # panics
1139/// panics if global parallelism is disabled.
1140#[track_caller]
1141pub fn get_global_parallelism() -> Par {
1142	let value = GLOBAL_PARALLELISM.load(core::sync::atomic::Ordering::Relaxed);
1143	match value {
1144		0 => panic!("Global parallelism is disabled."),
1145		1 => Par::Seq,
1146		#[cfg(feature = "rayon")]
1147		n => Par::rayon(n - 2),
1148		#[cfg(not(feature = "rayon"))]
1149		_ => unreachable!(),
1150	}
1151}
1152
1153#[doc(hidden)]
1154pub mod hacks;
1155
1156/// statistics and randomness functionality
1157pub mod stats;
1158
1159mod non_exhaustive {
1160	#[doc(hidden)]
1161	#[derive(Debug, Copy, Clone, PartialEq, Eq)]
1162	#[repr(transparent)]
1163	pub struct NonExhaustive(pub(crate) ());
1164}
1165pub(crate) use non_exhaustive::NonExhaustive;
1166
1167/// like `Default`, but with an extra type parameter so that algorithm hyperparameters can be tuned
1168/// per scalar type.
1169pub trait Auto<T> {
1170	/// returns the default value for the type `T`
1171	fn auto() -> Self;
1172}
1173
1174/// implements [`Default`] based on `Config`'s [`Auto`] implementation for the type `T`.
1175pub struct Spec<Config, T> {
1176	/// wrapped config value
1177	pub config: Config,
1178	__marker: core::marker::PhantomData<fn() -> T>,
1179}
1180
1181impl<Config, T> core::ops::Deref for Spec<Config, T> {
1182	type Target = Config;
1183
1184	#[inline]
1185	fn deref(&self) -> &Self::Target {
1186		&self.config
1187	}
1188}
1189
1190impl<Config, T> core::ops::DerefMut for Spec<Config, T> {
1191	#[inline]
1192	fn deref_mut(&mut self) -> &mut Self::Target {
1193		&mut self.config
1194	}
1195}
1196
1197impl<Config: Copy, T> Copy for Spec<Config, T> {}
1198impl<Config: Clone, T> Clone for Spec<Config, T> {
1199	#[inline]
1200	fn clone(&self) -> Self {
1201		Self::new(self.config.clone())
1202	}
1203}
1204impl<Config: core::fmt::Debug, T> core::fmt::Debug for Spec<Config, T> {
1205	#[inline]
1206	fn fmt(&self, f: &mut core::fmt::Formatter<'_>) -> core::fmt::Result {
1207		self.config.fmt(f)
1208	}
1209}
1210
1211impl<Config, T> Spec<Config, T> {
1212	/// wraps the given config value
1213	#[inline]
1214	pub fn new(config: Config) -> Self {
1215		Spec {
1216			config,
1217			__marker: core::marker::PhantomData,
1218		}
1219	}
1220}
1221
1222impl<T, Config> From<Config> for Spec<Config, T> {
1223	#[inline]
1224	fn from(config: Config) -> Self {
1225		Spec {
1226			config,
1227			__marker: core::marker::PhantomData,
1228		}
1229	}
1230}
1231
1232impl<T, Config: Auto<T>> Default for Spec<Config, T> {
1233	#[inline]
1234	fn default() -> Self {
1235		Spec {
1236			config: Auto::<T>::auto(),
1237			__marker: core::marker::PhantomData,
1238		}
1239	}
1240}
1241
1242mod into_range {
1243	use super::*;
1244
1245	pub trait IntoRange<I> {
1246		type Len<N: Shape>: Shape;
1247
1248		fn into_range(self, min: I, max: I) -> core::ops::Range<I>;
1249	}
1250
1251	impl<I> IntoRange<I> for core::ops::Range<I> {
1252		type Len<N: Shape> = usize;
1253
1254		#[inline]
1255		fn into_range(self, _: I, _: I) -> core::ops::Range<I> {
1256			self
1257		}
1258	}
1259	impl<I> IntoRange<I> for core::ops::RangeFrom<I> {
1260		type Len<N: Shape> = usize;
1261
1262		#[inline]
1263		fn into_range(self, _: I, max: I) -> core::ops::Range<I> {
1264			self.start..max
1265		}
1266	}
1267	impl<I> IntoRange<I> for core::ops::RangeTo<I> {
1268		type Len<N: Shape> = usize;
1269
1270		#[inline]
1271		fn into_range(self, min: I, _: I) -> core::ops::Range<I> {
1272			min..self.end
1273		}
1274	}
1275	impl<I> IntoRange<I> for core::ops::RangeFull {
1276		type Len<N: Shape> = N;
1277
1278		#[inline]
1279		fn into_range(self, min: I, max: I) -> core::ops::Range<I> {
1280			min..max
1281		}
1282	}
1283}
1284
1285mod sort;
1286
1287pub extern crate dyn_stack;
1288pub extern crate faer_traits as traits;
1289pub extern crate reborrow;