pub mod data;
use data::{Nmrc, Numeric, Sprs, Symb};
#[derive(Copy, Clone, Debug)]
pub enum Error {
NotPositiveDefinite,
NoPivot,
}
impl std::fmt::Display for Error {
fn fmt(&self, f: &mut std::fmt::Formatter<'_>) -> std::fmt::Result {
match self {
Self::NoPivot => write!(f, "Could not find a pivot"),
Self::NotPositiveDefinite => write!(f, "Could not complete Cholesky factorization. Please provide a positive definite matrix"),
}
}
}
impl std::error::Error for Error {}
pub fn add<T: Numeric<T>>(a: &Sprs<T>, b: &Sprs<T>, alpha: T, beta: T) -> Sprs<T> {
let mut nz = 0;
let m = a.m;
let n = b.n;
let anz = a.p[a.n] as usize;
let bnz = b.p[n] as usize;
let mut w = vec![0; m];
let mut x = vec![T::zero(); m];
let mut c = Sprs::zeros(m, n, anz + bnz);
for j in 0..n {
c.p[j] = nz as isize; nz = scatter(a, j, alpha, &mut w[..], &mut x[..], j + 1, &mut c, nz); nz = scatter(b, j, beta, &mut w[..], &mut x[..], j + 1, &mut c, nz);
for p in c.p[j] as usize..nz {
c.x[p] = x[c.i[p]];
}
}
c.p[n] = nz as isize;
c.quick_trim();
c
}
#[must_use]
pub fn chol<T: Numeric<T>>(a: &Sprs<T>, s: &mut Symb) -> Result<Nmrc<T>, Error> {
let mut top;
let mut d;
let mut lki;
let mut i;
let n = a.n;
let mut n_mat = Nmrc::new();
let mut w = vec![0; 3 * n];
let ws = n; let wc = 2 * n; let mut x = vec![T::zero(); n];
let c = match s.pinv {
Some(_) => symperm(a, &s.pinv),
None => a.clone(),
};
n_mat.l = Sprs::zeros(n, n, s.cp[n] as usize);
for k in 0..n {
w[wc + k] = s.cp[k]; n_mat.l.p[k] = w[wc + k];
x[k] = T::zero(); w[k] = k as isize; top = ereach(&c, k, &s.parent[..], ws, &mut w[..], &mut x[..], n); d = x[k]; x[k] = T::zero();
while top < n {
i = w[ws + top] as usize; lki = x[i] / n_mat.l.x[n_mat.l.p[i] as usize]; x[i] = T::zero(); for p in (n_mat.l.p[i] + 1)..w[wc + i] as isize {
x[n_mat.l.i[p as usize]] -= n_mat.l.x[p as usize] * lki;
}
d -= lki * lki; let p = w[wc + i] as usize;
w[wc + i] += 1;
n_mat.l.i[p] = k; n_mat.l.x[p] = lki;
top += 1;
}
if d <= T::zero() {
return Err(Error::NotPositiveDefinite);
}
let p = w[wc + k];
w[wc + k] += 1;
n_mat.l.i[p as usize] = k; n_mat.l.x[p as usize] = T::powf(d, 0.5);
}
n_mat.l.p[n] = s.cp[n];
Ok(n_mat)
}
#[must_use]
pub fn cholsol<T: Numeric<T>>(a: &Sprs<T>, b: &mut [T], order: i8) -> Result<(), Error> {
let n = a.n;
let mut s = schol(a, order); let n_mat = chol(a, &mut s)?; let mut x = vec![T::zero(); n];
ipvec(n, &s.pinv, b, &mut x[..]); lsolve(&n_mat.l, &mut x); ltsolve(&n_mat.l, &mut x); pvec(n, &s.pinv, &x[..], &mut b[..]);
Ok(())
}
pub fn gaxpy<T: Numeric<T>>(a_mat: &Sprs<T>, x: &[T], y: &[T]) -> Vec<T> {
let mut r = y.to_vec();
for (j, &x_j) in x.iter().enumerate().take(a_mat.n) {
for p in a_mat.p[j]..a_mat.p[j + 1] {
let p_u = p as usize;
r[a_mat.i[p_u]] += a_mat.x[p_u] * x_j;
}
}
r
}
pub fn lsolve<T: Numeric<T>>(l: &Sprs<T>, x: &mut [T]) {
for j in 0..l.n {
x[j] /= l.x[l.p[j] as usize];
for p in (l.p[j] + 1) as usize..l.p[j + 1] as usize {
x[l.i[p]] -= l.x[p] * x[j];
}
}
}
pub fn ltsolve<T: Numeric<T>>(l: &Sprs<T>, x: &mut [T]) {
for j in (0..l.n).rev() {
for p in (l.p[j] + 1) as usize..l.p[j + 1] as usize {
x[j] -= l.x[p] * x[l.i[p]];
}
x[j] /= l.x[l.p[j] as usize];
}
}
#[must_use]
pub fn lu<T: Numeric<T>>(a: &Sprs<T>, s: &mut Symb, tol: T) -> Result<Nmrc<T>, Error> {
let n = a.n;
let mut col;
let mut top;
let mut ipiv;
let mut a_f;
let mut t;
let mut pivot;
let mut x = vec![T::zero(); n];
let mut xi = vec![0; 2 * n];
let mut n_mat = Nmrc {
l: Sprs::zeros(n, n, s.lnz), u: Sprs::zeros(n, n, s.unz),
pinv: Some(vec![0; n]),
b: Vec::new(),
};
x[0..n].fill(T::zero()); n_mat.pinv.as_mut().unwrap()[0..n].fill(-1); n_mat.l.p[0..n + 1].fill(0);
s.lnz = 0;
s.unz = 0;
for k in 0..n {
n_mat.l.p[k] = s.lnz as isize; n_mat.u.p[k] = s.unz as isize;
if s.lnz + n > n_mat.l.nzmax {
let nsz = 2 * n_mat.l.nzmax + n;
n_mat.l.nzmax = nsz;
n_mat.l.i.resize(nsz, 0);
n_mat.l.x.resize(nsz, T::zero());
}
if s.unz + n > n_mat.u.nzmax {
let nsz = 2 * n_mat.u.nzmax + n;
n_mat.u.nzmax = nsz;
n_mat.u.i.resize(nsz, 0);
n_mat.u.x.resize(nsz, T::zero());
}
col = s.q.as_ref().map_or(k, |q| q[k] as usize);
top = splsolve(&mut n_mat.l, a, col, &mut xi[..], &mut x[..], &n_mat.pinv);
ipiv = -1;
a_f = -T::one();
for &i in xi[top..n].iter() {
let i = i as usize;
if n_mat.pinv.as_ref().unwrap()[i] < 0 {
t = T::abs(x[i]);
if t > a_f {
a_f = t; ipiv = i as isize;
}
} else {
n_mat.u.i[s.unz] = n_mat.pinv.as_ref().unwrap()[i] as usize;
n_mat.u.x[s.unz] = x[i];
s.unz += 1;
}
}
if ipiv == -1 || a_f <= T::zero() {
return Err(Error::NoPivot);
}
if n_mat.pinv.as_ref().unwrap()[col] < 0 && T::abs(x[col]) >= a_f * tol {
ipiv = col as isize;
}
pivot = x[ipiv as usize]; n_mat.u.i[s.unz] = k; n_mat.u.x[s.unz] = pivot;
s.unz += 1;
n_mat.pinv.as_mut().unwrap()[ipiv as usize] = k as isize; n_mat.l.i[s.lnz] = ipiv as usize; n_mat.l.x[s.lnz] = T::one();
s.lnz += 1;
for &i in xi[top..n].iter() {
let i = i as usize;
if n_mat.pinv.as_ref().unwrap()[i] < 0 {
n_mat.l.i[s.lnz] = i; n_mat.l.x[s.lnz] = x[i] / pivot; s.lnz += 1
}
x[i] = T::zero(); }
}
n_mat.l.p[n] = s.lnz as isize;
n_mat.u.p[n] = s.unz as isize;
for p in 0..s.lnz {
n_mat.l.i[p] = n_mat.pinv.as_ref().unwrap()[n_mat.l.i[p]] as usize;
}
n_mat.l.quick_trim();
n_mat.u.quick_trim();
Ok(n_mat)
}
#[must_use]
pub fn lusol<T: Numeric<T>>(a: &Sprs<T>, b: &mut [T], order: i8, tol: T) -> Result<(), Error> {
let mut x = vec![T::zero(); a.n];
let mut s;
s = sqr(a, order, false); let n = lu(a, &mut s, tol)?;
ipvec(a.n, &n.pinv, b, &mut x[..]); lsolve(&n.l, &mut x); usolve(&n.u, &mut x[..]); ipvec(a.n, &s.q, &x[..], &mut b[..]); Ok(())
}
pub fn multiply<T: Numeric<T>>(a: &Sprs<T>, b: &Sprs<T>) -> Sprs<T> {
let mut nz = 0;
let mut w = vec![0; a.m];
let mut x = vec![T::zero(); a.m];
let mut c = Sprs::zeros(a.m, b.n, 2 * (a.p[a.n] + b.p[b.n]) as usize + a.m);
for j in 0..b.n {
if nz + a.m > c.nzmax {
let nsz = 2 * c.nzmax + a.m;
c.nzmax = nsz;
c.i.resize(nsz, 0);
c.x.resize(nsz, T::zero());
}
c.p[j] = nz as isize; for p in b.p[j]..b.p[j + 1] {
nz = scatter(
a,
b.i[p as usize],
b.x[p as usize],
&mut w[..],
&mut x[..],
j + 1,
&mut c,
nz,
);
}
for p in c.p[j] as usize..nz {
c.x[p] = x[c.i[p]];
}
}
c.p[b.n] = nz as isize;
c.quick_trim();
c
}
pub fn norm<T: Numeric<T>>(a: &Sprs<T>) -> T {
let mut norm_r = T::zero();
for j in 0..a.n {
let mut s = T::zero();
for p in a.p[j] as usize..a.p[j + 1] as usize {
s += a.x[p].abs();
}
norm_r = T::max(norm_r, s);
}
norm_r
}
pub fn qr<T: Numeric<T>>(a: &Sprs<T>, s: &Symb) -> Nmrc<T> {
let mut p1;
let mut top;
let mut col;
let mut i;
let m = a.m;
let n = a.n;
let mut vnz = s.lnz;
let mut rnz = s.unz;
let mut v = Sprs::zeros(s.m2, n, vnz); let mut r = Sprs::zeros(s.m2, n, rnz);
let leftmost_p = m + n; let mut w = vec![0; s.m2 + n];
let ws = s.m2; let mut x = vec![T::zero(); s.m2];
let mut n_mat = Nmrc::new();
let mut beta = vec![T::zero(); n];
w[0..s.m2].fill(-1); rnz = 0;
vnz = 0;
for k in 0..n {
r.p[k] = rnz as isize; v.p[k] = vnz as isize; p1 = vnz;
w[k] = k as isize; v.i[vnz] = k;
vnz += 1;
top = n;
col = s.q.as_ref().map_or(k, |q| q[k] as usize);
for p in a.p[col]..a.p[col + 1] {
i = s.pinv.as_ref().unwrap()[leftmost_p + a.i[p as usize]]; let mut len = 0;
while w[i as usize] != k as isize {
w[ws + len] = i;
len += 1;
w[i as usize] = k as isize;
i = s.parent[i as usize];
}
for j in 1..(len + 1) {
top -= 1;
w[ws + top] = w[ws + len - j]; }
i = s.pinv.as_ref().unwrap()[a.i[p as usize]]; x[i as usize] = a.x[p as usize]; if i > k as isize && w[i as usize] < k as isize {
v.i[vnz] = i as usize; vnz += 1;
w[i as usize] = k as isize;
}
}
for p in top..n {
i = w[ws + p]; happly(&v, i as usize, beta[i as usize], &mut x[..]); r.i[rnz] = i as usize; r.x[rnz] = x[i as usize];
rnz += 1;
x[i as usize] = T::zero();
if s.parent[i as usize] == k as isize {
vnz = scatter_no_x(i as usize, &mut w[..], k, &mut v, vnz);
}
}
for p in p1..vnz {
v.x[p] = x[v.i[p]];
x[v.i[p]] = T::zero();
}
r.i[rnz] = k; r.x[rnz] = house(&mut v.x[..], Some(p1), &mut beta[..], Some(k), vnz - p1); rnz += 1;
}
r.p[n] = rnz as isize; v.p[n] = vnz as isize;
n_mat.l = v;
n_mat.u = r;
n_mat.b = beta;
n_mat
}
pub fn qrsol<T: Numeric<T>>(a: &Sprs<T>, b: &mut [T], order: i8) {
let n = a.n;
let m = a.m;
if m >= n {
let s = sqr(a, order, true); let n_mat = qr(a, &s); let mut x = vec![T::zero(); s.m2];
ipvec(m, &s.pinv, b, &mut x[..]); for k in 0..n {
happly(&n_mat.l, k, n_mat.b[k], &mut x[..]);
}
usolve(&n_mat.u, &mut x[..]); ipvec(n, &s.q, &x[..], &mut b[..]); } else {
let at = transpose(a); let s = sqr(&at, order, true); let n_mat = qr(&at, &s); let mut x = vec![T::zero(); s.m2];
pvec(m, &s.q, b, &mut x[..]); utsolve(&n_mat.u, &mut x[..]); for k in (0..m).rev() {
happly(&n_mat.l, k, n_mat.b[k], &mut x[..]);
}
pvec(n, &s.pinv, &x[..], &mut b[..]); }
}
pub fn schol<T: Numeric<T>>(a: &Sprs<T>, order: i8) -> Symb {
let n = a.n;
let mut s = Symb::new(); let p = amd(a, order); s.pinv = pinvert(&p, n); drop(p);
let c_mat = symperm(a, &s.pinv); s.parent = etree(&c_mat, false); let post = post::<T>(n, &s.parent[..]); let mut c = counts(&c_mat, &s.parent[..], &post[..], false); drop(post);
drop(c_mat);
s.cp = vec![0; n + 1]; s.unz = cumsum(&mut s.cp[..], &mut c[..], n);
s.lnz = s.unz;
drop(c);
s
}
pub fn scpmat<T: Numeric<T>>(alpha: T, a: &Sprs<T>) -> Sprs<T> {
let mut c = Sprs::<T>::new();
c.m = a.m;
c.n = a.n;
c.nzmax = a.nzmax;
c.p = a.p.clone();
c.i = a.i.clone();
c.x = a.x.iter().map(|x| *x + alpha).collect();
c
}
pub fn scxmat<T: Numeric<T>>(alpha: T, a: &Sprs<T>) -> Sprs<T> {
let mut c = Sprs::<T>::new();
c.m = a.m;
c.n = a.n;
c.nzmax = a.nzmax;
c.p = a.p.clone();
c.i = a.i.clone();
c.x = a.x.iter().map(|x| *x * alpha).collect();
c
}
pub fn sprs_print<T: Numeric<T>>(a: &Sprs<T>, brief: bool) {
let m = a.m;
let n = a.n;
let nzmax = a.nzmax;
println!(
"{}-by-{}, nzmax: {} nnz: {}, 1-norm: {}",
m,
n,
nzmax,
a.p[n],
norm(a)
);
for j in 0..n {
println!(
" col {} : locations {} to {}",
j,
a.p[j],
a.p[j + 1] - 1
);
for p in a.p[j]..a.p[j + 1] {
println!(" {} : {}", a.i[p as usize], a.x[p as usize]);
if brief && p > 20 {
println!(" ...");
return;
}
}
}
}
pub fn sqr<T: Numeric<T>>(a: &Sprs<T>, order: i8, qr: bool) -> Symb {
let mut s = Symb::new();
let pst;
s.q = amd(a, order); if qr {
let c = if order >= 0 {
permute(a, &None, &s.q)
} else {
a.clone()
};
s.parent = etree(&c, true); pst = post::<T>(a.n, &s.parent[..]);
s.cp = counts(&c, &s.parent[..], &pst[..], true); s.pinv = vcount(&c, &s.parent[..], &mut s.m2, &mut s.lnz);
s.unz = 0;
for k in 0..a.n {
s.unz += s.cp[k] as usize;
}
} else {
s.unz = 4 * a.p[a.n] as usize + a.n; s.lnz = s.unz; }
s
}
pub fn transpose<T: Numeric<T>>(a: &Sprs<T>) -> Sprs<T> {
let mut q;
let mut w = vec![0; a.m];
let mut c = Sprs::zeros(a.n, a.m, a.p[a.n] as usize);
for p in 0..a.p[a.n] as usize {
w[a.i[p]] += 1; }
cumsum(&mut c.p[..], &mut w[..], a.m); for j in 0..a.n {
for p in a.p[j] as usize..a.p[j + 1] as usize {
q = w[a.i[p]] as usize;
c.i[q] = j; c.x[q] = a.x[p];
w[a.i[p]] += 1;
}
}
c
}
pub fn usolve<T: Numeric<T>>(u: &Sprs<T>, x: &mut [T]) {
for j in (0..u.n).rev() {
x[j] /= u.x[(u.p[j + 1] - 1) as usize];
for p in u.p[j]..u.p[j + 1] - 1 {
x[u.i[p as usize]] -= u.x[p as usize] * x[j];
}
}
}
pub fn utsolve<T: Numeric<T>>(u: &Sprs<T>, x: &mut [T]) {
for j in 0..u.n {
for p in u.p[j] as usize..(u.p[j + 1] - 1) as usize {
x[j] -= u.x[p] * x[u.i[p]];
}
x[j] /= u.x[(u.p[j + 1] - 1) as usize];
}
}
fn amd<T: Numeric<T>>(a: &Sprs<T>, order: i8) -> Option<Vec<isize>> {
let mut dense;
let mut c;
let mut nel = 0;
let mut mindeg = 0;
let mut elenk;
let mut nvk;
let mut p;
let mut p2;
let mut pk1;
let mut pk2;
let mut e;
let mut pj;
let mut ln;
let mut nvi;
let mut i;
let mut mark_v;
let mut lemax = 0;
let mut eln;
let mut wnvi;
let mut p1;
let mut pn;
let mut h;
let mut d;
let mut dext;
let mut p3;
let mut p4;
let mut j;
let mut nvj;
let mut jlast;
if order < 0 {
return None;
}
let mut at = transpose(a); let m = a.m;
let n = a.n;
dense = std::cmp::max(16, (10. * f32::sqrt(n as f32)) as isize); dense = std::cmp::min((n - 2) as isize, dense);
if order == 0 && n == m {
c = add(a, &at, T::zero(), T::zero()); } else if order == 1 {
p2 = 0; for j in 0..m {
p = at.p[j]; at.p[j] = p2; if at.p[j + 1] - p > dense {
continue; }
while p < at.p[j + 1] {
at.i[p2 as usize] = at.i[p as usize];
p2 += 1;
p += 1;
}
}
at.p[m] = p2; let a2 = transpose(&at); c = multiply(&at, &a2); } else {
c = multiply(&at, a); }
drop(at);
let mut p_v = vec![0; n + 1]; let mut ww = vec![0; 8 * (n + 1)]; let len = 0; let nv = n + 1; let next = 2 * (n + 1); let head = 3 * (n + 1); let elen = 4 * (n + 1); let degree = 5 * (n + 1); let w = 6 * (n + 1); let hhead = 7 * (n + 1); let last = 0;
fkeep(&mut c, &diag); let mut cnz = c.p[n];
let nsz = cnz as usize + cnz as usize / 5 + 2 * n;
c.nzmax = nsz;
c.i.resize(nsz, 0);
c.x.resize(nsz, T::zero());
for k in 0..n {
ww[len + k] = c.p[k + 1] - c.p[k];
}
ww[len + n] = 0;
for i in 0..=n {
ww[head + i] = -1; p_v[last + i] = -1;
ww[next + i] = -1;
ww[hhead + i] = -1; ww[nv + i] = 1; ww[w + i] = 1; ww[elen + i] = 0; ww[degree + i] = ww[len + i]; }
mark_v = wclear(0, 0, &mut ww[..], w, n); ww[elen + n] = -2; c.p[n] = -1; ww[w + n] = 0;
for i in 0..n {
let d = ww[degree + i];
if d == 0 {
ww[elen + i] = -2; nel += 1;
c.p[i] = -1; ww[w + i] = 0;
} else if d > dense {
ww[nv + i] = 0; ww[elen + i] = -1; nel += 1;
c.p[i] = flip(n as isize);
ww[nv + n] += 1;
} else {
if ww[(head as isize + d) as usize] != -1 {
let wt = ww[(head as isize + d) as usize];
p_v[(last as isize + wt) as usize] = i as isize;
}
ww[next + i] = ww[(head as isize + d) as usize]; ww[(head as isize + d) as usize] = i as isize;
}
}
while nel < n {
let mut k;
loop {
k = ww[head + mindeg];
if !(mindeg < n && k == -1) {
break;
}
mindeg += 1;
}
if ww[(next as isize + k) as usize] != -1 {
let wt = ww[(next as isize + k) as usize];
p_v[(last as isize + wt) as usize] = -1;
}
ww[head + mindeg] = ww[(next as isize + k) as usize]; elenk = ww[(elen as isize + k) as usize]; nvk = ww[(nv as isize + k) as usize]; nel += nvk as usize;
if elenk > 0 && (cnz + mindeg as isize) as usize >= c.nzmax {
for j in 0..n {
p = c.p[j];
if p >= 0 {
c.p[j] = c.i[p as usize] as isize; c.i[p as usize] = flip(j as isize) as usize; }
}
let mut q = 0;
p = 0;
while p < cnz {
let j = flip(c.i[p as usize] as isize);
p += 1;
if j >= 0 {
c.i[q] = c.p[j as usize] as usize; c.p[j as usize] = q as isize; q += 1;
for _ in 0..ww[(len as isize + j) as usize] - 1 {
c.i[q] = c.i[p as usize];
q += 1;
p += 1;
}
}
}
cnz = q as isize; }
let mut dk = 0;
ww[(nv as isize + k) as usize] = -nvk; p = c.p[k as usize];
if elenk == 0 {
pk1 = p;
} else {
pk1 = cnz;
}
pk2 = pk1;
for k1 in 1..=(elenk + 1) as usize {
if k1 > elenk as usize {
e = k; pj = p; ln = ww[(len as isize + k) as usize] - elenk; } else {
e = c.i[p as usize] as isize; p += 1;
pj = c.p[e as usize];
ln = ww[(len as isize + e) as usize]; }
for _ in 1..=ln {
i = c.i[pj as usize] as isize;
pj += 1;
nvi = ww[(nv as isize + i) as usize];
if nvi <= 0 {
continue; }
dk += nvi; ww[(nv as isize + i) as usize] = -nvi; c.i[pk2 as usize] = i as usize; pk2 += 1;
if ww[(next as isize + i) as usize] != -1 {
let wt = ww[(next as isize + i) as usize];
p_v[(last as isize + wt) as usize] = p_v[last + i as usize];
}
if p_v[(last as isize + i) as usize] != -1 {
let wt = p_v[(last as isize + i) as usize];
ww[(next as isize + wt) as usize] = ww[(next as isize + i) as usize];
} else {
let wt = ww[degree + i as usize];
ww[(head as isize + wt) as usize] = ww[next + i as usize];
}
}
if e != k {
c.p[e as usize] = flip(k); ww[(w as isize + e) as usize] = 0; }
}
if elenk != 0 {
cnz = pk2; }
ww[(degree as isize + k) as usize] = dk; c.p[k as usize] = pk1; ww[(len as isize + k) as usize] = pk2 - pk1;
ww[(elen as isize + k) as usize] = -2;
mark_v = wclear(mark_v, lemax, &mut ww[..], w, n); for pk in pk1..pk2 {
i = c.i[pk as usize] as isize;
eln = ww[(elen as isize + i) as usize];
if eln <= 0 {
continue; }
nvi = -ww[(nv as isize + i) as usize]; wnvi = mark_v - nvi;
for p in c.p[i as usize] as usize..=(c.p[i as usize] + eln - 1) as usize {
e = c.i[p] as isize;
if ww[(w as isize + e) as usize] >= mark_v {
ww[(w as isize + e) as usize] -= nvi; } else if ww[(w as isize + e) as usize] != 0 {
ww[(w as isize + e) as usize] = ww[(degree as isize + e) as usize] + wnvi;
}
}
}
for pk in pk1..pk2 {
i = c.i[pk as usize] as isize; p1 = c.p[i as usize];
p2 = p1 + ww[(elen as isize + i) as usize] - 1;
pn = p1;
h = 0;
d = 0;
for p in p1..=p2 {
e = c.i[p as usize] as isize;
if ww[(w as isize + e) as usize] != 0 {
dext = ww[(w as isize + e) as usize] - mark_v; if dext > 0 {
d += dext; c.i[pn as usize] = e as usize; pn += 1;
h += e as usize; } else {
c.p[e as usize] = flip(k); ww[(w as isize + e) as usize] = 0; }
}
}
ww[(elen as isize + i) as usize] = pn - p1 + 1; p3 = pn;
p4 = p1 + ww[(len as isize + i) as usize];
for p in p2 + 1..p4 {
j = c.i[p as usize] as isize;
nvj = ww[(nv as isize + j) as usize];
if nvj <= 0 {
continue; }
d += nvj; c.i[pn as usize] = j as usize; pn += 1;
h += j as usize; }
if d == 0 {
c.p[i as usize] = flip(k); nvi = -ww[(nv as isize + i) as usize];
dk -= nvi; nvk += nvi; nel += nvi as usize;
ww[(nv as isize + i) as usize] = 0;
ww[(elen as isize + i) as usize] = -1; } else {
ww[(degree as isize + i) as usize] =
std::cmp::min(ww[(degree as isize + i) as usize], d); c.i[pn as usize] = c.i[p3 as usize]; c.i[p3 as usize] = c.i[p1 as usize]; c.i[p1 as usize] = k as usize; ww[(len as isize + i) as usize] = pn - p1 + 1; h %= n; ww[(next as isize + i) as usize] = ww[hhead + h]; ww[hhead + h] = i;
p_v[(last as isize + i) as usize] = h as isize; }
} ww[(degree as isize + k) as usize] = dk; lemax = std::cmp::max(lemax, dk);
mark_v = wclear(mark_v + lemax, lemax, &mut ww[..], w, n);
for pk in pk1..pk2 {
i = c.i[pk as usize] as isize;
if ww[(nv as isize + i) as usize] >= 0 {
continue; }
h = p_v[(last as isize + i) as usize] as usize; i = ww[hhead + h];
ww[hhead + h] = -1;
while i != -1 && ww[(next as isize + i) as usize] != -1 {
ln = ww[(len as isize + i) as usize];
eln = ww[(elen as isize + i) as usize];
for p in c.p[i as usize] + 1..=c.p[i as usize] + ln - 1 {
ww[w + c.i[p as usize]] = mark_v;
}
jlast = i;
let mut ok;
j = ww[(next as isize + i) as usize];
while j != -1 {
ok = ww[(len as isize + j) as usize] == ln
&& ww[(elen as isize + j) as usize] == eln;
p = c.p[j as usize] + 1;
while ok && p < c.p[j as usize] + ln {
if ww[w + c.i[p as usize]] != mark_v {
ok = false;
}
p += 1;
}
if ok {
c.p[j as usize] = flip(i); ww[(nv as isize + i) as usize] += ww[(nv as isize + j) as usize];
ww[(nv as isize + j) as usize] = 0;
ww[(elen as isize + j) as usize] = -1; j = ww[(next as isize + j) as usize]; ww[(next as isize + jlast) as usize] = j;
} else {
jlast = j; j = ww[(next as isize + j) as usize];
}
}
i = ww[(next as isize + i) as usize];
mark_v += 1;
}
}
p = pk1;
for pk in pk1..pk2 {
i = c.i[pk as usize] as isize;
nvi = -ww[(nv as isize + i) as usize];
if nvi <= 0 {
continue; }
ww[(nv as isize + i) as usize] = nvi; d = ww[(degree as isize + i) as usize] + dk - nvi; d = std::cmp::min(d, n as isize - nel as isize - nvi);
if ww[(head as isize + d) as usize] != -1 {
let wt = ww[(head as isize + d) as usize];
p_v[(last as isize + wt) as usize] = i;
}
ww[(next as isize + i) as usize] = ww[(head as isize + d) as usize]; p_v[(last as isize + i) as usize] = -1;
ww[(head as isize + d) as usize] = i;
mindeg = std::cmp::min(mindeg, d as usize); ww[(degree as isize + i) as usize] = d;
c.i[p as usize] = i as usize; p += 1;
}
ww[(nv as isize + k) as usize] = nvk; ww[(len as isize + k) as usize] = p - pk1;
if ww[(len as isize + k) as usize] == 0 {
c.p[k as usize] = -1; ww[(w as isize + k) as usize] = 0; }
if elenk != 0 {
cnz = p; }
}
for i in 0..n {
c.p[i] = flip(c.p[i]); }
for j in 0..=n {
ww[head + j] = -1;
}
for j in (0..=n).rev() {
if ww[nv + j] > 0 {
continue; }
ww[next + j] = ww[(head as isize + c.p[j]) as usize]; ww[(head as isize + c.p[j]) as usize] = j as isize;
}
for e in (0..=n).rev() {
if ww[nv + e] <= 0 {
continue; }
if c.p[e] != -1 {
ww[next + e] = ww[(head as isize + c.p[e]) as usize]; ww[(head as isize + c.p[e]) as usize] = e as isize;
}
}
let mut k = 0;
for i in 0..=n {
if c.p[i] == -1 {
k = tdfs::<T>(i as isize, k, &mut ww[..], head, next, &mut p_v[..], w);
}
}
Some(p_v)
}
fn cedge(
j: isize,
i: isize,
w: &mut [isize],
first: usize,
maxfirst: usize,
delta_colcount: &mut [isize],
prevleaf: usize,
ancestor: usize,
) {
let mut q;
let mut s;
let mut sparent;
if i <= j || w[(first as isize + j) as usize] <= w[(maxfirst as isize + i) as usize] {
return;
}
w[(maxfirst as isize + i) as usize] = w[(first as isize + j) as usize]; let jprev = w[(prevleaf as isize + i) as usize]; delta_colcount[j as usize] += 1; if jprev != -1 {
q = jprev;
while q != w[(ancestor as isize + q) as usize] {
q = w[(ancestor as isize + q) as usize];
}
s = jprev;
while s != q {
sparent = w[(ancestor as isize + s) as usize]; w[(ancestor as isize + s) as usize] = q;
s = sparent;
}
delta_colcount[q as usize] -= 1; }
w[(prevleaf as isize + i) as usize] = j; }
fn counts<T: Numeric<T>>(a: &Sprs<T>, parent: &[isize], post: &[isize], ata: bool) -> Vec<isize> {
let m = a.m;
let n = a.n;
let s = if ata { 4 * n + n + m + 1 } else { 4 * n };
let mut w = vec![0; s]; let first = 3 * n; let ancestor = 0; let maxfirst = n; let prevleaf = 2 * n; let head = 4 * n; let next = 5 * n + 1; let mut delta_colcount = vec![0; n]; let mut j;
let mut k;
let at = transpose(a);
w[0..s].fill(-1); for (k, &p) in post.iter().enumerate() {
j = p;
delta_colcount[j as usize] = match w[(first as isize + j) as usize] {
-1 => 1,
_ => 0,
};
while j != -1 && w[(first as isize + j) as usize] == -1 {
w[(first as isize + j) as usize] = k as isize;
j = parent[j as usize];
}
}
if ata {
for k in 0..n {
w[post[k] as usize] = k as isize; }
for i in 0..m {
k = n; for p in at.p[i]..at.p[i + 1] {
k = std::cmp::min(k, w[at.i[p as usize]] as usize);
}
w[next + i] = w[head + k]; w[head + k] = i as isize;
}
}
for i in 0..n {
w[ancestor + i] = i as isize; }
for k in 0..n {
j = post[k]; if parent[j as usize] != -1 {
delta_colcount[parent[j as usize] as usize] -= 1; }
if ata {
let mut ii = w[head + k];
while ii != -1 {
for p in at.p[ii as usize]..at.p[ii as usize + 1] {
cedge(
j,
at.i[p as usize] as isize,
&mut w[..],
first,
maxfirst,
&mut delta_colcount[..],
prevleaf,
ancestor,
);
}
ii = w[(next as isize + ii) as usize];
}
} else {
for p in at.p[j as usize]..at.p[j as usize + 1] {
cedge(
j,
at.i[p as usize] as isize,
&mut w[..],
first,
maxfirst,
&mut delta_colcount[..],
prevleaf,
ancestor,
);
}
}
if parent[j as usize] != -1 {
w[(ancestor as isize + j) as usize] = parent[j as usize];
}
}
for j in 0..n {
if parent[j] != -1 {
delta_colcount[parent[j] as usize] += delta_colcount[j];
}
}
delta_colcount
}
fn cumsum(p: &mut [isize], c: &mut [isize], n: usize) -> usize {
let mut nz = 0;
for (p_i, c_i) in p.iter_mut().zip(c.iter_mut()).take(n) {
*p_i = nz;
nz += *c_i;
*c_i = *p_i;
}
p[n] = nz;
nz as usize
}
fn dfs<T: Numeric<T>>(
j: usize,
l: &mut Sprs<T>,
top: usize,
xi: &mut [isize],
pstack_i: &usize,
pinv: &Option<Vec<isize>>,
) -> usize {
let mut i;
let mut j = j;
let mut jnew;
let mut head = 0;
let mut done;
let mut p2;
let mut top = top;
xi[0] = j as isize; while head >= 0 {
j = xi[head as usize] as usize; if pinv.is_some() {
jnew = pinv.as_ref().unwrap()[j];
} else {
jnew = j as isize;
}
if !marked(&l.p[..], j) {
mark(&mut l.p[..], j); if jnew < 0 {
xi[pstack_i + head as usize] = 0;
} else {
xi[pstack_i + head as usize] = unflip(l.p[jnew as usize]);
}
}
done = true; if jnew < 0 {
p2 = 0;
} else {
p2 = unflip(l.p[(jnew + 1) as usize]);
}
for p in xi[pstack_i + head as usize]..p2 {
i = l.i[p as usize]; if marked(&l.p[..], i) {
continue; }
xi[pstack_i + head as usize] = p; head += 1;
xi[head as usize] = i as isize; done = false; break; }
if done {
head -= 1; top -= 1;
xi[top] = j as isize; }
}
top
}
fn diag<T: Numeric<T>>(i: isize, j: isize, _: T) -> bool {
i != j
}
fn ereach<T: Numeric<T>>(
a: &Sprs<T>,
k: usize,
parent: &[isize],
s: usize,
w: &mut [isize],
x: &mut [T],
top: usize,
) -> usize {
let mut top = top;
let mut i;
let mut len;
for p in a.p[k]..a.p[k + 1] {
i = a.i[p as usize] as isize; if i > k as isize {
continue; }
x[i as usize] = a.x[p as usize]; len = 0;
while w[i as usize] != k as isize {
w[s + len] = i; len += 1;
w[i as usize] = k as isize;
i = parent[i as usize];
}
for j in 1..(len + 1) {
top -= 1;
w[s + top] = w[s + len - j]; }
}
top }
fn etree<T: Numeric<T>>(a: &Sprs<T>, ata: bool) -> Vec<isize> {
let mut parent = vec![0; a.n];
let mut i;
let mut inext;
let mut w = if ata {
vec![0; a.n + a.m]
} else {
vec![0; a.n]
};
let ancestor = 0;
let prev = ancestor + a.n;
if ata {
w[prev..prev + a.m].fill(-1);
}
for k in 0..a.n {
parent[k] = -1; w[ancestor + k] = -1; for p in a.p[k] as usize..a.p[k + 1] as usize {
if ata {
i = w[prev + a.i[p]];
} else {
i = a.i[p] as isize;
}
while i != -1 && i < k as isize {
inext = w[(ancestor as isize + i) as usize]; w[(ancestor as isize + i) as usize] = k as isize; if inext == -1 {
parent[i as usize] = k as isize; }
i = inext
}
if ata {
w[prev + a.i[p]] = k as isize;
}
}
}
parent
}
fn fkeep<T: Numeric<T>>(a: &mut Sprs<T>, f: &dyn Fn(isize, isize, T) -> bool) -> isize {
let mut p;
let mut nz = 0;
let n = a.n;
for j in 0..n {
p = a.p[j]; a.p[j] = nz; while p < a.p[j + 1] {
if f(a.i[p as usize] as isize, j as isize, a.x[p as usize]) {
a.x[nz as usize] = a.x[p as usize]; a.i[nz as usize] = a.i[p as usize];
nz += 1;
}
p += 1;
}
}
a.p[n] = nz;
nz
}
fn happly<T: Numeric<T>>(v: &Sprs<T>, i: usize, beta: T, x: &mut [T]) {
let mut tau = T::zero();
for p in v.p[i]..v.p[i + 1] {
tau += v.x[p as usize] * x[v.i[p as usize]];
}
tau *= beta; for p in v.p[i]..v.p[i + 1] {
x[v.i[p as usize]] -= v.x[p as usize] * tau;
}
}
fn house<T: Numeric<T>>(
x: &mut [T],
xp: Option<usize>,
beta: &mut [T],
betap: Option<usize>,
n: usize,
) -> T {
let s;
let xp = xp.unwrap_or(0);
let betap = betap.unwrap_or(0);
let sigma: T = (1..n).map(|i| x[i + xp] * x[i + xp]).sum();
if sigma != T::zero() {
s = (x[xp] * x[xp] + sigma).sqrt(); x[xp] = if x[xp] <= T::zero() {
x[xp] - s
} else {
-sigma / (x[xp] + s)
};
beta[betap] = T::one() / (-s * x[xp]);
} else {
s = x[xp].abs(); beta[betap] = if x[xp] <= T::zero() {
T::one() + T::one()
} else {
T::zero()
};
x[xp] = T::one();
}
s
}
fn ipvec<T: Numeric<T>>(n: usize, p: &Option<Vec<isize>>, b: &[T], x: &mut [T]) {
for k in 0..n {
if p.is_some() {
x[p.as_ref().unwrap()[k] as usize] = b[k];
} else {
x[k] = b[k];
}
}
}
fn permute<T: Numeric<T>>(
a: &Sprs<T>,
pinv: &Option<Vec<isize>>,
q: &Option<Vec<isize>>,
) -> Sprs<T> {
let mut j;
let mut nz = 0;
let mut c = Sprs::zeros(a.m, a.n, a.p[a.n] as usize);
for k in 0..a.n {
c.p[k] = nz as isize; if q.is_some() {
j = q.as_ref().unwrap()[k] as usize;
} else {
j = k;
}
for p in a.p[j] as usize..a.p[j + 1] as usize {
c.x[nz] = a.x[p]; if pinv.is_some() {
c.i[nz] = pinv.as_ref().unwrap()[a.i[p]] as usize;
} else {
c.i[nz] = a.i[p];
}
nz += 1;
}
}
c.p[a.n] = nz as isize;
c
}
fn pinvert(p: &Option<Vec<isize>>, n: usize) -> Option<Vec<isize>> {
if p.is_none() {
return None;
}
let mut pinv = vec![0; n]; for k in 0..n {
pinv[p.as_ref().unwrap()[k] as usize] = k as isize; }
Some(pinv)
}
fn post<T: Numeric<T>>(n: usize, parent: &[isize]) -> Vec<isize> {
let mut k = 0;
let mut post = vec![0; n]; let mut w = vec![0; 3 * n]; let head = 0; let next = head + n; let stack = head + 2 * n;
for j in 0..n {
w[head + j] = -1; }
for j in (0..n).rev() {
if parent[j] == -1 {
continue; }
w[next + j] = w[(head as isize + parent[j]) as usize]; w[(head as isize + parent[j]) as usize] = j as isize;
}
for (j, par) in parent.iter().enumerate().take(n) {
if *par != -1 {
continue; }
k = tdfs::<T>(j as isize, k, &mut w[..], head, next, &mut post[..], stack);
}
post
}
fn pvec<T: Numeric<T>>(n: usize, p: &Option<Vec<isize>>, b: &[T], x: &mut [T]) {
for (k, x_k) in x.iter_mut().enumerate().take(n) {
*x_k = match p {
Some(p) => b[p[k] as usize],
None => b[k],
};
}
}
fn reach<T: Numeric<T>>(
l: &mut Sprs<T>,
b: &Sprs<T>,
k: usize,
xi: &mut [isize],
pinv: &Option<Vec<isize>>,
) -> usize {
let mut top = l.n;
for p in b.p[k] as usize..b.p[k + 1] as usize {
if !marked(&l.p[..], b.i[p]) {
let n = l.n;
top = dfs(b.i[p], l, top, &mut xi[..], &n, pinv);
}
}
for i in xi.iter().take(l.n).skip(top) {
mark(&mut l.p[..], *i as usize); }
top
}
fn scatter<T: Numeric<T>>(
a: &Sprs<T>,
j: usize,
beta: T,
w: &mut [isize],
x: &mut [T],
mark: usize,
c: &mut Sprs<T>,
nz: usize,
) -> usize {
let mut i;
let mut nzo = nz;
for p in a.p[j] as usize..a.p[j + 1] as usize {
i = a.i[p]; if w[i] < mark as isize {
w[i] = mark as isize; c.i[nzo] = i; nzo += 1;
x[i] = beta * a.x[p]; } else {
x[i] += beta * a.x[p]; }
}
nzo
}
fn scatter_no_x<T: Numeric<T>>(
j: usize,
w: &mut [isize],
mark: usize,
c: &mut Sprs<T>,
nz: usize,
) -> usize {
let mut i;
let mut nzo = nz;
for p in c.p[j] as usize..c.p[j + 1] as usize {
i = c.i[p]; if w[i] < mark as isize {
w[i] = mark as isize; c.i[nzo] = i; nzo += 1;
}
}
nzo
}
fn splsolve<T: Numeric<T>>(
l: &mut Sprs<T>,
b: &Sprs<T>,
k: usize,
xi: &mut [isize],
x: &mut [T],
pinv: &Option<Vec<isize>>,
) -> usize {
let mut jnew;
let top = reach(l, b, k, &mut xi[..], pinv);
for p in top..l.n {
x[xi[p] as usize] = T::zero(); }
for p in b.p[k] as usize..b.p[k + 1] as usize {
x[b.i[p]] = b.x[p]; }
for j in xi.iter().take(l.n).skip(top) {
let j_u = *j as usize; jnew = match pinv {
Some(pinv) => pinv[j_u], None => *j, };
if jnew < 0 {
continue; }
for p in (l.p[jnew as usize] + 1) as usize..l.p[jnew as usize + 1] as usize {
x[l.i[p]] -= l.x[p] * x[j_u]; }
}
top }
fn symperm<T: Numeric<T>>(a: &Sprs<T>, pinv: &Option<Vec<isize>>) -> Sprs<T> {
let n = a.n;
let mut i;
let mut i2;
let mut j2;
let mut q;
let mut c = Sprs::zeros(n, n, a.p[n] as usize);
let mut w = vec![0; n];
for j in 0..n {
j2 = pinv.as_ref().map_or(j, |pinv| pinv[j] as usize);
for p in a.p[j]..a.p[j + 1] {
i = a.i[p as usize];
if i > j {
continue; }
i2 = pinv.as_ref().map_or(i, |pinv| pinv[i] as usize); w[std::cmp::max(i2, j2)] += 1; }
}
cumsum(&mut c.p[..], &mut w[..], n); for j in 0..n {
j2 = pinv.as_ref().map_or(j, |pinv| pinv[j] as usize); for p in a.p[j]..a.p[j + 1] {
i = a.i[p as usize];
if i > j {
continue; }
i2 = pinv.as_ref().map_or(i, |pinv| pinv[i] as usize); q = w[std::cmp::max(i2, j2)] as usize;
w[std::cmp::max(i2, j2)] += 1;
c.i[q] = std::cmp::min(i2, j2);
c.x[q] = a.x[p as usize];
}
}
c
}
fn tdfs<T: Numeric<T>>(
j: isize,
k: isize,
ww: &mut [isize],
head: usize,
next: usize,
post: &mut [isize],
stack: usize,
) -> isize {
let mut i;
let mut p;
let mut top = 0;
let mut k = k;
ww[stack] = j; while top >= 0 {
p = ww[(stack as isize + top) as usize]; i = ww[(head as isize + p) as usize]; match i {
-1 => {
top -= 1; post[k as usize] = p; k += 1;
}
_ => {
ww[(head as isize + p) as usize] = ww[(next as isize + i) as usize]; top += 1;
ww[(stack as isize + top) as usize] = i; }
}
}
k
}
fn vcount<T: Numeric<T>>(
a: &Sprs<T>,
parent: &[isize],
m2: &mut usize,
vnz: &mut usize,
) -> Option<Vec<isize>> {
let n = a.n;
let m = a.m;
let mut pinv: Vec<isize> = vec![0; 2 * m + n];
let leftmost = m + n; let mut w = vec![0; m + 3 * n];
let next = 0;
let head = m;
let tail = m + n;
let nque = m + 2 * n;
w[head..head + n].fill(-1); w[tail..tail + n].fill(-1);
w[nque..nque + n].fill(0);
pinv[leftmost..leftmost + m].fill(-1);
for k in (0..n).rev() {
for p in a.p[k]..a.p[k + 1] {
pinv[leftmost + a.i[p as usize]] = k as isize; }
}
let mut k;
for i in (0..m).rev() {
pinv[i] = -1; k = pinv[leftmost + i];
if k == -1 {
continue; }
if w[(nque as isize + k) as usize] == 0 {
w[(tail as isize + k) as usize] = i as isize; }
w[(nque as isize + k) as usize] += 1;
w[next + i] = w[(head as isize + k) as usize]; w[(head as isize + k) as usize] = i as isize;
}
*vnz = 0;
*m2 = m;
let mut i;
for k in 0..n {
i = w[head + k]; *vnz += 1; if i < 0 {
i = *m2 as isize; *m2 += 1;
}
pinv[i as usize] = k as isize; w[nque + k] -= 1;
if w[nque + k] <= 0 {
continue; }
*vnz += w[nque + k] as usize; let pa = parent[k]; if pa != -1 {
if w[(nque as isize + pa) as usize] == 0 {
w[(tail as isize + pa) as usize] = w[tail + k];
}
let tw = w[tail + k];
w[(next as isize + tw) as usize] = w[(head as isize + pa) as usize];
w[(head as isize + pa) as usize] = w[(next as isize + i) as usize];
w[(nque as isize + pa) as usize] += w[nque + k];
}
}
let mut k = n;
for p in pinv.iter_mut().take(m) {
if *p < 0 {
*p = k as isize;
k += 1;
}
}
Some(pinv)
}
fn wclear(mark_v: isize, lemax: isize, ww: &mut [isize], w: usize, n: usize) -> isize {
let mut mark = mark_v;
if mark < 2 || (mark + lemax < 0) {
for k in 0..n {
if ww[w + k] != 0 {
ww[w + k] = 1;
}
}
mark = 2;
}
mark }
#[inline]
fn flip(i: isize) -> isize {
-(i) - 2
}
#[inline]
fn unflip(i: isize) -> isize {
if i < 0 {
flip(i)
} else {
i
}
}
#[inline]
fn marked(ap: &[isize], j: usize) -> bool {
ap[j] < 0
}
#[inline]
fn mark(ap: &mut [isize], j: usize) {
ap[j] = flip(ap[j])
}