#![deny(rustdoc::broken_intra_doc_links, rustdoc::private_intra_doc_links)]
#[cfg(any(
feature = "intel-mkl-static", //decprecated
feature = "intel-mkl-dynamic", //decprecated
feature = "intel-mkl-dynamic-lp64-iomp",
feature = "intel-mkl-dynamic-lp64-seq",
feature = "intel-mkl-static-ilp64-iomp",
feature = "intel-mkl-static-lp64-iomp",
feature = "intel-mkl-dynamic-ilp64-iomp",
feature = "intel-mkl-static-ilp64-seq",
feature = "intel-mkl-static-lp64-seq"
))]
extern crate intel_mkl_src as _src;
#[cfg(any(feature = "openblas-system", feature = "openblas-static"))]
extern crate openblas_src as _src;
#[cfg(any(feature = "netlib-system", feature = "netlib-static"))]
extern crate netlib_src as _src;
pub mod alloc;
pub mod cholesky;
pub mod eig;
pub mod eig_generalized;
pub mod eigh;
pub mod eigh_generalized;
pub mod error;
pub mod flags;
pub mod layout;
pub mod least_squares;
pub mod opnorm;
pub mod qr;
pub mod rcond;
pub mod solve;
pub mod solveh;
pub mod svd;
pub mod svddc;
pub mod triangular;
pub mod tridiagonal;
pub use crate::eig_generalized::GeneralizedEigenvalue;
pub use self::flags::*;
pub use self::least_squares::LeastSquaresOwned;
pub use self::svd::{SvdOwned, SvdRef};
pub use self::tridiagonal::{LUFactorizedTridiagonal, Tridiagonal};
use self::{alloc::*, error::*, layout::*};
use cauchy::*;
use std::mem::MaybeUninit;
pub type Pivot = Vec<i32>;
#[cfg_attr(doc, katexit::katexit)]
pub trait Lapack: Scalar {
fn eig(
calc_v: bool,
l: MatrixLayout,
a: &mut [Self],
) -> Result<(Vec<Self::Complex>, Vec<Self::Complex>)>;
fn eig_generalized(
calc_v: bool,
l: MatrixLayout,
a: &mut [Self],
b: &mut [Self],
thresh_opt: Option<Self::Real>,
) -> Result<(
Vec<GeneralizedEigenvalue<Self::Complex>>,
Vec<Self::Complex>,
)>;
fn eigh(
calc_eigenvec: bool,
layout: MatrixLayout,
uplo: UPLO,
a: &mut [Self],
) -> Result<Vec<Self::Real>>;
fn eigh_generalized(
calc_eigenvec: bool,
layout: MatrixLayout,
uplo: UPLO,
a: &mut [Self],
b: &mut [Self],
) -> Result<Vec<Self::Real>>;
fn householder(l: MatrixLayout, a: &mut [Self]) -> Result<Vec<Self>>;
fn q(l: MatrixLayout, a: &mut [Self], tau: &[Self]) -> Result<()>;
fn qr(l: MatrixLayout, a: &mut [Self]) -> Result<Vec<Self>>;
fn svd(l: MatrixLayout, calc_u: bool, calc_vt: bool, a: &mut [Self]) -> Result<SvdOwned<Self>>;
fn svddc(layout: MatrixLayout, jobz: JobSvd, a: &mut [Self]) -> Result<SvdOwned<Self>>;
fn least_squares(
a_layout: MatrixLayout,
a: &mut [Self],
b: &mut [Self],
) -> Result<LeastSquaresOwned<Self>>;
fn least_squares_nrhs(
a_layout: MatrixLayout,
a: &mut [Self],
b_layout: MatrixLayout,
b: &mut [Self],
) -> Result<LeastSquaresOwned<Self>>;
fn lu(l: MatrixLayout, a: &mut [Self]) -> Result<Pivot>;
fn inv(l: MatrixLayout, a: &mut [Self], p: &Pivot) -> Result<()>;
fn solve(l: MatrixLayout, t: Transpose, a: &[Self], p: &Pivot, b: &mut [Self]) -> Result<()>;
fn bk(l: MatrixLayout, uplo: UPLO, a: &mut [Self]) -> Result<Pivot>;
fn invh(l: MatrixLayout, uplo: UPLO, a: &mut [Self], ipiv: &Pivot) -> Result<()>;
fn solveh(l: MatrixLayout, uplo: UPLO, a: &[Self], ipiv: &Pivot, b: &mut [Self]) -> Result<()>;
fn cholesky(l: MatrixLayout, uplo: UPLO, a: &mut [Self]) -> Result<()>;
fn inv_cholesky(l: MatrixLayout, uplo: UPLO, a: &mut [Self]) -> Result<()>;
fn solve_cholesky(l: MatrixLayout, uplo: UPLO, a: &[Self], b: &mut [Self]) -> Result<()>;
fn rcond(l: MatrixLayout, a: &[Self], anorm: Self::Real) -> Result<Self::Real>;
fn opnorm(t: NormType, l: MatrixLayout, a: &[Self]) -> Self::Real;
fn solve_triangular(
al: MatrixLayout,
bl: MatrixLayout,
uplo: UPLO,
d: Diag,
a: &[Self],
b: &mut [Self],
) -> Result<()>;
fn lu_tridiagonal(a: Tridiagonal<Self>) -> Result<LUFactorizedTridiagonal<Self>>;
fn rcond_tridiagonal(lu: &LUFactorizedTridiagonal<Self>) -> Result<Self::Real>;
fn solve_tridiagonal(
lu: &LUFactorizedTridiagonal<Self>,
bl: MatrixLayout,
t: Transpose,
b: &mut [Self],
) -> Result<()>;
}
macro_rules! impl_lapack {
($s:ty) => {
impl Lapack for $s {
fn eig(
calc_v: bool,
l: MatrixLayout,
a: &mut [Self],
) -> Result<(Vec<Self::Complex>, Vec<Self::Complex>)> {
use eig::*;
let work = EigWork::<$s>::new(calc_v, l)?;
let EigOwned { eigs, vr, vl } = work.eval(a)?;
Ok((eigs, vr.or(vl).unwrap_or_default()))
}
fn eig_generalized(
calc_v: bool,
l: MatrixLayout,
a: &mut [Self],
b: &mut [Self],
thresh_opt: Option<Self::Real>,
) -> Result<(
Vec<GeneralizedEigenvalue<Self::Complex>>,
Vec<Self::Complex>,
)> {
use eig_generalized::*;
let work = EigGeneralizedWork::<$s>::new(calc_v, l)?;
let eig_generalized_owned = work.eval(a, b)?;
let eigs = eig_generalized_owned.calc_eigs(thresh_opt);
let vr = eig_generalized_owned.vr;
let vl = eig_generalized_owned.vl;
Ok((eigs, vr.or(vl).unwrap_or_default()))
}
fn eigh(
calc_eigenvec: bool,
layout: MatrixLayout,
uplo: UPLO,
a: &mut [Self],
) -> Result<Vec<Self::Real>> {
use eigh::*;
let work = EighWork::<$s>::new(calc_eigenvec, layout)?;
work.eval(uplo, a)
}
fn eigh_generalized(
calc_eigenvec: bool,
layout: MatrixLayout,
uplo: UPLO,
a: &mut [Self],
b: &mut [Self],
) -> Result<Vec<Self::Real>> {
use eigh_generalized::*;
let work = EighGeneralizedWork::<$s>::new(calc_eigenvec, layout)?;
work.eval(uplo, a, b)
}
fn householder(l: MatrixLayout, a: &mut [Self]) -> Result<Vec<Self>> {
use qr::*;
let work = HouseholderWork::<$s>::new(l)?;
work.eval(a)
}
fn q(l: MatrixLayout, a: &mut [Self], tau: &[Self]) -> Result<()> {
use qr::*;
let mut work = QWork::<$s>::new(l)?;
work.calc(a, tau)?;
Ok(())
}
fn qr(l: MatrixLayout, a: &mut [Self]) -> Result<Vec<Self>> {
let tau = Self::householder(l, a)?;
let r = Vec::from(&*a);
Self::q(l, a, &tau)?;
Ok(r)
}
fn svd(
l: MatrixLayout,
calc_u: bool,
calc_vt: bool,
a: &mut [Self],
) -> Result<SvdOwned<Self>> {
use svd::*;
let work = SvdWork::<$s>::new(l, calc_u, calc_vt)?;
work.eval(a)
}
fn svddc(layout: MatrixLayout, jobz: JobSvd, a: &mut [Self]) -> Result<SvdOwned<Self>> {
use svddc::*;
let work = SvdDcWork::<$s>::new(layout, jobz)?;
work.eval(a)
}
fn least_squares(
l: MatrixLayout,
a: &mut [Self],
b: &mut [Self],
) -> Result<LeastSquaresOwned<Self>> {
let b_layout = l.resized(b.len() as i32, 1);
Self::least_squares_nrhs(l, a, b_layout, b)
}
fn least_squares_nrhs(
a_layout: MatrixLayout,
a: &mut [Self],
b_layout: MatrixLayout,
b: &mut [Self],
) -> Result<LeastSquaresOwned<Self>> {
use least_squares::*;
let work = LeastSquaresWork::<$s>::new(a_layout, b_layout)?;
work.eval(a, b)
}
fn lu(l: MatrixLayout, a: &mut [Self]) -> Result<Pivot> {
use solve::*;
LuImpl::lu(l, a)
}
fn inv(l: MatrixLayout, a: &mut [Self], p: &Pivot) -> Result<()> {
use solve::*;
let mut work = InvWork::<$s>::new(l)?;
work.calc(a, p)?;
Ok(())
}
fn solve(
l: MatrixLayout,
t: Transpose,
a: &[Self],
p: &Pivot,
b: &mut [Self],
) -> Result<()> {
use solve::*;
SolveImpl::solve(l, t, a, p, b)
}
fn bk(l: MatrixLayout, uplo: UPLO, a: &mut [Self]) -> Result<Pivot> {
use solveh::*;
let work = BkWork::<$s>::new(l)?;
work.eval(uplo, a)
}
fn invh(l: MatrixLayout, uplo: UPLO, a: &mut [Self], ipiv: &Pivot) -> Result<()> {
use solveh::*;
let mut work = InvhWork::<$s>::new(l)?;
work.calc(uplo, a, ipiv)
}
fn solveh(
l: MatrixLayout,
uplo: UPLO,
a: &[Self],
ipiv: &Pivot,
b: &mut [Self],
) -> Result<()> {
use solveh::*;
SolvehImpl::solveh(l, uplo, a, ipiv, b)
}
fn cholesky(l: MatrixLayout, uplo: UPLO, a: &mut [Self]) -> Result<()> {
use cholesky::*;
CholeskyImpl::cholesky(l, uplo, a)
}
fn inv_cholesky(l: MatrixLayout, uplo: UPLO, a: &mut [Self]) -> Result<()> {
use cholesky::*;
InvCholeskyImpl::inv_cholesky(l, uplo, a)
}
fn solve_cholesky(
l: MatrixLayout,
uplo: UPLO,
a: &[Self],
b: &mut [Self],
) -> Result<()> {
use cholesky::*;
SolveCholeskyImpl::solve_cholesky(l, uplo, a, b)
}
fn rcond(l: MatrixLayout, a: &[Self], anorm: Self::Real) -> Result<Self::Real> {
use rcond::*;
let mut work = RcondWork::<$s>::new(l);
work.calc(a, anorm)
}
fn opnorm(t: NormType, l: MatrixLayout, a: &[Self]) -> Self::Real {
use opnorm::*;
let mut work = OperatorNormWork::<$s>::new(t, l);
work.calc(a)
}
fn solve_triangular(
al: MatrixLayout,
bl: MatrixLayout,
uplo: UPLO,
d: Diag,
a: &[Self],
b: &mut [Self],
) -> Result<()> {
use triangular::*;
SolveTriangularImpl::solve_triangular(al, bl, uplo, d, a, b)
}
fn lu_tridiagonal(a: Tridiagonal<Self>) -> Result<LUFactorizedTridiagonal<Self>> {
use tridiagonal::*;
let work = LuTridiagonalWork::<$s>::new(a.l);
work.eval(a)
}
fn rcond_tridiagonal(lu: &LUFactorizedTridiagonal<Self>) -> Result<Self::Real> {
use tridiagonal::*;
let mut work = RcondTridiagonalWork::<$s>::new(lu.a.l);
work.calc(lu)
}
fn solve_tridiagonal(
lu: &LUFactorizedTridiagonal<Self>,
bl: MatrixLayout,
t: Transpose,
b: &mut [Self],
) -> Result<()> {
use tridiagonal::*;
SolveTridiagonalImpl::solve_tridiagonal(lu, bl, t, b)
}
}
};
}
impl_lapack!(c64);
impl_lapack!(c32);
impl_lapack!(f64);
impl_lapack!(f32);