ldl/lib.rs
1extern crate num_traits;
2
3#[cfg(test)]
4mod tests;
5
6use num_traits::{Bounded, Float, FromPrimitive, NumAssignOps, PrimInt, ToPrimitive};
7
8#[derive(Debug, PartialEq, Clone)]
9pub enum Marker {
10 Unused,
11 Used,
12}
13
14/// Compute the elimination tree for a quasidefinite matrix
15/// in compressed sparse column form, where the input matrix is
16/// assumed to contain data for the upper triangular part of `A` only,
17/// and there are no duplicate indices.
18///
19/// Returns an elimination tree for the factorization `A = LDL^T` and a
20/// count of the nonzeros in each column of `L` that are strictly below the
21/// diagonal.
22///
23/// Does not allocate. It is assumed that the arrays `work`, `l_nz`, and
24/// `etree` will be allocated with a number of elements equal to `n`.
25///
26/// The data in (`n`,`a_p`,`a_i`) are from a square matrix `A` in CSC format, and
27/// should include the upper triangular part of `A` only.
28///
29/// This function is only intended for factorisation of QD matrices specified
30/// by their upper triangular part. An error is returned if any column has
31/// data below the diagonal or is completely empty.
32///
33/// For matrices with a non-empty column but a zero on the corresponding diagonal,
34/// this function will *not* return an error, as it may still be possible to factor
35/// such a matrix in LDL form. No promises are made in this case though...
36///
37/// # Arguments
38///
39/// * `n` - number of columns in CSC matrix `A` (assumed square)
40/// * `a_p` - column pointers (size `n`+1) for columns of `A`
41/// * `a_i` - row indices of `A`. Has `a_p[n]` elements
42/// * `work` - work vector (size `n`) (no meaning on return)
43/// * `l_nz` - count of nonzeros in each column of `L` (size `n`) below diagonal
44/// * `etree` - elimination tree (size `n`)
45///
46/// # Returns
47///
48/// Sum of `Lnz` (i.e. total nonzeros in `L` below diagonal).
49/// Returns 1 if the input is not triu or has an empty column.
50/// Returns 2 if the return value overflows.
51pub fn etree<I: PrimInt + NumAssignOps + Bounded + ToPrimitive + FromPrimitive>(
52 n: I,
53 a_p: &[I],
54 a_i: &[I],
55 work: &mut [I],
56 l_nz: &mut [I],
57 etree: &mut [Option<I>],
58) -> Result<I, I> {
59 for i in 0..n.to_usize().unwrap() {
60 // Zero out Lnz and work. Set all etree values to unknown.
61 work[i] = I::zero();
62 l_nz[i] = I::zero();
63 etree[i] = None;
64
65 // Abort if A doesn't have at least one entry
66 // one entry in every column.
67 if a_p[i] == a_p[i + 1] {
68 return Err(I::one());
69 }
70 }
71
72 for j in 0..n.to_usize().unwrap() {
73 work[j] = I::from_usize(j).unwrap();
74 for p in a_p[j].to_usize().unwrap()..a_p[j + 1].to_usize().unwrap() {
75 let mut i = a_i[p].to_usize().unwrap();
76 if i > j {
77 // Abort if entries on lower triangle.
78 return Err(I::one());
79 }
80 while work[i].to_usize().unwrap() != j {
81 if etree[i].is_none() {
82 etree[i] = Some(I::from_usize(j).unwrap());
83 }
84 l_nz[i] += I::one(); // nonzeros in this column
85 work[i] = I::from_usize(j).unwrap();
86 i = etree[i].unwrap().to_usize().unwrap();
87 }
88 }
89 }
90
91 // Compute the total nonzeros in L. This much
92 // space is required to store Li and Lx. Return
93 // error code -2 if the nonzero count will overflow
94 // its integer type.
95 let mut sum_l_nz = I::zero();
96 for i in 0..n.to_usize().unwrap() {
97 if sum_l_nz > I::max_value() - l_nz[i] {
98 return Err(I::from(2).unwrap());
99 } else {
100 sum_l_nz += l_nz[i];
101 }
102 }
103
104 Ok(sum_l_nz)
105}
106
107/// Compute an LDL decomposition for a quasidefinite matrix
108/// in compressed sparse column form, where the input matrix is
109/// assumed to contain data for the upper triangular part of `A` only,
110/// and there are no duplicate indices.
111///
112/// Returns factors `L`, `D` and `Dinv = 1./D`.
113///
114/// Does not allocate. It is assumed that `L` will be a compressed
115/// sparse column matrix with data (`n`,`l_p`,`l_i`,`l_x`) with sufficient space
116/// allocated, with a number of nonzeros equal to the count given
117/// as a return value by `etree`.
118///
119/// # Arguments
120///
121/// * `n` - number of columns in `L` and `A` (both square)
122/// * `a_p` - column pointers (size `n`+1) for columns of `A` (not modified)
123/// * `a_i` - row indices of `A`. Has `a_p[n]` elements (not modified)
124/// * `a_x` - data of `A`. Has `a_p[n]` elements (not modified)
125/// * `l_p` - column pointers (size `n`+1) for columns of `L`
126/// * `l_i` - row indices of `L`. Has `l_p[n]` elements
127/// * `l_x` - data of `L`. Has `l_p[n]` elements
128/// * `d` - vectorized factor `D`. Length is `n`
129/// * `d_inv` - reciprocal of `D`. Length is `n`
130/// * `l_nz` - count of nonzeros in each column of `L` below diagonal,
131/// as given by `etree()` (not modified)
132/// * `etree` - elimination tree as as given by `etree()` (not modified)
133/// * `bwork` - working array of bools. Length is `n`
134/// * `iwork` - working array of integers. Length is 3*`n`
135/// * `fwork` - working array of floats. Length is `n`
136///
137/// # Returns
138///
139/// Returns a count of the number of positive elements in `D`.
140/// Returns 1 and exits immediately if any element of `D` evaluates
141/// exactly to zero (matrix is not quasidefinite or otherwise LDL factorisable)
142pub fn factor<F: Float + NumAssignOps, I: PrimInt + NumAssignOps + FromPrimitive + ToPrimitive>(
143 n: I,
144 a_p: &[I],
145 a_i: &[I],
146 a_x: &[F],
147 l_p: &mut [I],
148 l_i: &mut [I],
149 l_x: &mut [F],
150 d: &mut [F],
151 d_inv: &mut [F],
152 l_nz: &[I],
153 etree: &[Option<I>],
154 bwork: &mut [Marker],
155 iwork: &mut [I],
156 fwork: &mut [F],
157) -> Result<I, I> {
158 let un = n.to_usize().unwrap();
159
160 let mut nnz_y: usize;
161 let mut bidx: usize;
162 let mut cidx: usize;
163 let mut next_idx: Option<usize>;
164 let mut nnz_e: usize;
165 let mut tmp_idx: usize;
166
167 let mut positive_values_in_d = I::zero();
168
169 // Partition working memory into pieces.
170 let y_markers = bwork;
171 let (y_idx, iwork) = iwork.split_at_mut(un);
172 let (elim_buffer, iwork) = iwork.split_at_mut(un);
173 let (l_next_space_in_col, _) = iwork.split_at_mut(un);
174 let y_vals = fwork;
175
176 l_p[0] = I::zero(); // first column starts at index zero
177
178 for i in 0..un {
179 // Compute L column indices.
180 l_p[i + 1] = l_p[i] + l_nz[i]; // cumsum, total at the end
181
182 // Set all Yidx to be 'unused' initially
183 // in each column of L, the next available space
184 // to start is just the first space in the column
185 y_markers[i] = Marker::Unused;
186 y_vals[i] = F::zero();
187 d[i] = F::zero();
188 l_next_space_in_col[i] = l_p[i];
189 }
190
191 // First element of the diagonal D.
192 d[0] = a_x[0];
193 if d[0] == F::zero() {
194 return Err(I::one());
195 }
196 if d[0] > F::zero() {
197 positive_values_in_d += I::one();
198 }
199 d_inv[0] = F::one() / d[0];
200
201 // Start from 1 here. The upper LH corner is trivially 0
202 // in L b/c we are only computing the subdiagonal elements.
203 for k in 1..un {
204 // NB : For each k, we compute a solution to
205 // y = L(0:(k-1),0:k-1))\b, where b is the kth
206 // column of A that sits above the diagonal.
207 // The solution y is then the kth row of L,
208 // with an implied '1' at the diagonal entry.
209
210 // Number of nonzeros in this row of L.
211 nnz_y = 0; // Number of elements in this row.
212
213 // This loop determines where nonzeros
214 // will go in the kth row of L, but doesn't
215 // compute the actual values.
216 tmp_idx = a_p[k + 1].to_usize().unwrap();
217
218 for i in a_p[k].to_usize().unwrap()..tmp_idx {
219 bidx = a_i[i].to_usize().unwrap(); // We are working on this element of b.
220
221 // Initialize D[k] as the element of this column
222 // corresponding to the diagonal place. Don't use
223 // this element as part of the elimination step
224 // that computes the k^th row of L.
225 if bidx == k {
226 d[k] = a_x[i];
227 continue;
228 }
229
230 y_vals[bidx] = a_x[i]; // initialise y(bidx) = b(bidx)
231
232 // Use the forward elimination tree to figure
233 // out which elements must be eliminated after
234 // this element of b.
235 next_idx = Some(bidx);
236
237 if y_markers[next_idx.unwrap()] == Marker::Unused {
238 // This y term not already visited.
239
240 y_markers[next_idx.unwrap()] = Marker::Used; // I touched this one.
241 elim_buffer[0] = I::from(next_idx.unwrap()).unwrap(); // It goes at the start of the current list.
242 nnz_e = 1; // Length of unvisited elimination path from here.
243
244 next_idx = etree[bidx].map(|next| next.to_usize().unwrap());
245
246 while next_idx.is_some() && next_idx.unwrap() < k {
247 if y_markers[next_idx.unwrap()] == Marker::Used {
248 break;
249 }
250
251 y_markers[next_idx.unwrap()] = Marker::Used; // I touched this one.
252 elim_buffer[nnz_e] = I::from_usize(next_idx.unwrap()).unwrap(); // It goes in the current list.
253 nnz_e += 1; // The list is one longer than before.
254
255 // One step further along tree.
256 next_idx = etree[next_idx.unwrap()].map(|next| next.to_usize().unwrap());
257 }
258
259 // Now I put the buffered elimination list into
260 // my current ordering in reverse order.
261 while nnz_e != 0 {
262 // yIdx[nnzY++] = elim_buffer[--nnzE];
263 nnz_e -= 1;
264 y_idx[nnz_y] = elim_buffer[nnz_e];
265 nnz_y += 1;
266 }
267 }
268 }
269
270 // This for loop places nonzeros values in the k^th row.
271 // let mut i: isize = nnz_y as isize - 1;
272 // while i >= 0 {
273 if nnz_y > 0 {
274 for i in (0..=(nnz_y - 1)).rev() {
275 // for(i = (nnzY-1); i >=0; i--){
276
277 // Which column are we working on?
278 cidx = y_idx[i as usize].to_usize().unwrap();
279
280 // Loop along the elements in this
281 // column of L and subtract to solve to y.
282 tmp_idx = l_next_space_in_col[cidx].to_usize().unwrap();
283 let y_vals_cidx = y_vals[cidx];
284 for j in l_p[cidx].to_usize().unwrap()..tmp_idx {
285 y_vals[l_i[j].to_usize().unwrap()] -= l_x[j] * y_vals_cidx;
286 }
287
288 // Now I have the cidx^th element of y = L\b.
289 // so compute the corresponding element of
290 // this row of L and put it into the right place.
291 l_i[tmp_idx] = I::from_usize(k).unwrap();
292 l_x[tmp_idx] = y_vals_cidx * d_inv[cidx];
293
294 // D[k] -= yVals[cidx]*yVals[cidx]*Dinv[cidx];
295 d[k] -= y_vals_cidx * l_x[tmp_idx];
296 l_next_space_in_col[cidx] += I::one();
297
298 // Reset the yvalues and indices back to zero and UNUSED
299 // once I'm done with them.
300 y_vals[cidx] = F::zero();
301 y_markers[cidx] = Marker::Unused;
302
303 // i -= 1;
304 }
305 }
306
307 // Maintain a count of the positive entries
308 // in D. If we hit a zero, we can't factor
309 // this matrix, so abort
310 if d[k] == F::zero() {
311 return Err(I::one());
312 }
313 if d[k] > F::zero() {
314 positive_values_in_d += I::one();
315 }
316
317 // Compute the inverse of the diagonal.
318 d_inv[k] = F::one() / d[k];
319 }
320
321 Ok(positive_values_in_d)
322}
323
324/// Solves `LDL'x = b`
325///
326/// It is assumed that `L` will be a compressed
327/// sparse column matrix with data (`n`,`l_p`,`l_i`,`l_x`).
328///
329/// # Arguments
330///
331/// * `n` - number of columns in `L`
332/// * `l_p` - column pointers (size `n`+1) for columns of `L`
333/// * `l_i` - row indices of `L`. Has `l_p[n]` elements
334/// * `l_x` - data of `L`. Has `l_p[n]` elements
335/// * `d_inv` - reciprocal of `D`. Length is `n`
336/// * `x` - initialized to `b`. Equal to `x` on return
337pub fn solve<F: Float + NumAssignOps, I: PrimInt>(
338 n: I,
339 l_p: &[I],
340 l_i: &[I],
341 l_x: &[F],
342 d_inv: &[F],
343 x: &mut [F],
344) {
345 lsolve(n, l_p, l_i, l_x, x);
346 for i in 0..n.to_usize().unwrap() {
347 x[i] *= d_inv[i];
348 }
349 ltsolve(n, l_p, l_i, l_x, x);
350}
351
352/// Solves `(L+I)x = b`
353///
354/// It is assumed that `L` will be a compressed
355/// sparse column matrix with data (`n`,`l_p`,`l_i`,`l_x`).
356///
357/// # Arguments
358///
359/// * `n` - number of columns in `L`
360/// * `l_p` - column pointers (size `n`+1) for columns of `L`
361/// * `l_i` - row indices of `L`. Has `l_p[n]` elements
362/// * `l_x` - data of `L`. Has `l_p[n]` elements
363/// * `x` - initialized to `b`. Equal to `x` on return
364pub fn lsolve<F: Float + NumAssignOps, I: PrimInt>(
365 n: I,
366 l_p: &[I],
367 l_i: &[I],
368 l_x: &[F],
369 x: &mut [F],
370) {
371 for i in 0..n.to_usize().unwrap() {
372 let val = x[i];
373 for j in l_p[i].to_usize().unwrap()..l_p[i + 1].to_usize().unwrap() {
374 x[l_i[j].to_usize().unwrap()] -= l_x[j] * val;
375 }
376 }
377}
378
379/// Solves `(L+I)'x = b`
380///
381/// It is assumed that `L` will be a compressed
382/// sparse column matrix with data (`n`,`l_p`,`l_i`,`l_x`).
383///
384/// # Arguments
385///
386/// * `n` - number of columns in `L`
387/// * `l_p` - column pointers (size n+1) for columns of `L`
388/// * `l_i` - row indices of `L`. Has `l_p[n]` elements
389/// * `l_x` - data of `L`. Has `l_p[n]` elements
390/// * `x` - initialized to `b`. Equal to `x` on return
391pub fn ltsolve<F: Float + NumAssignOps, I: PrimInt>(
392 n: I,
393 l_p: &[I],
394 l_i: &[I],
395 l_x: &[F],
396 x: &mut [F],
397) {
398 for i in (0..=n.to_usize().unwrap() - 1).rev() {
399 //for(i = n-1; i>=0; i--){
400 let mut val = x[i];
401 for j in l_p[i].to_usize().unwrap()..l_p[i + 1].to_usize().unwrap() {
402 val -= l_x[j] * x[l_i[j].to_usize().unwrap()];
403 }
404 x[i] = val;
405 }
406}
407
408pub fn factor_solve<
409 F: Float + NumAssignOps,
410 I: PrimInt + NumAssignOps + FromPrimitive + ToPrimitive + Clone,
411>(
412 a_n: I,
413 a_p: &[I],
414 a_i: &[I],
415 a_x: &[F],
416 b: &mut [F],
417) -> Result<(), I> {
418 let un = a_n.to_usize().unwrap();
419
420 let l_n = a_n;
421
422 // Pre-factorisation memory allocations //
423
424 // These can happen *before* the etree is calculated
425 // since the sizes are not sparsity pattern specific.
426
427 // For the elimination tree.
428 let mut etree: Vec<Option<I>> = vec![None; un];
429 let mut l_nz: Vec<I> = vec![I::zero(); un];
430
431 // For the L factors. Li and Lx are sparsity dependent
432 // so must be done after the etree is constructed.
433 let mut l_p: Vec<I> = vec![I::zero(); un + 1];
434 let mut d: Vec<F> = vec![F::zero(); un];
435 let mut d_inv: Vec<F> = vec![F::zero(); un];
436
437 // Working memory. Note that both the etree and factor
438 // calls requires a working vector of int, with
439 // the factor function requiring 3*An elements and the
440 // etree only An elements. Just allocate the larger
441 // amount here and use it in both places.
442 let mut iwork: Vec<I> = vec![I::zero(); 3 * un];
443 let mut bwork: Vec<Marker> = vec![Marker::Unused; un];
444 let mut fwork: Vec<F> = vec![F::zero(); un];
445
446 // Elimination tree calculation //
447
448 let sum_l_nz = crate::etree(a_n, a_p, a_i, &mut iwork, &mut l_nz, &mut etree)?;
449
450 // LDL factorisation //
451
452 let mut l_i: Vec<I> = vec![I::zero(); sum_l_nz.to_usize().unwrap()];
453 let mut l_x: Vec<F> = vec![F::zero(); sum_l_nz.to_usize().unwrap()];
454
455 factor(
456 a_n, a_p, a_i, a_x, &mut l_p, &mut l_i, &mut l_x, &mut d, &mut d_inv, &l_nz, &etree,
457 &mut bwork, &mut iwork, &mut fwork,
458 )?;
459
460 // Solve //
461
462 solve(l_n, &l_p, &l_i, &l_x, &d_inv, b);
463
464 Ok(())
465}