nalgebra_sparse/factorization/
cholesky.rs

1use crate::csc::CscMatrix;
2use crate::ops::Op;
3use crate::ops::serial::spsolve_csc_lower_triangular;
4use crate::pattern::SparsityPattern;
5use core::{iter, mem};
6use nalgebra::{DMatrix, DMatrixView, DMatrixViewMut, RealField};
7use std::fmt::{Display, Formatter};
8
9/// A symbolic sparse Cholesky factorization of a CSC matrix.
10///
11/// The symbolic factorization computes the sparsity pattern of `L`, the Cholesky factor.
12#[derive(Debug, Clone, PartialEq, Eq)]
13pub struct CscSymbolicCholesky {
14    // Pattern of the original matrix that was decomposed
15    m_pattern: SparsityPattern,
16    l_pattern: SparsityPattern,
17    // u in this context is L^T, so that M = L L^T
18    u_pattern: SparsityPattern,
19}
20
21impl CscSymbolicCholesky {
22    /// Compute the symbolic factorization for a sparsity pattern belonging to a CSC matrix.
23    ///
24    /// The sparsity pattern must be symmetric. However, this is not enforced, and it is the
25    /// responsibility of the user to ensure that this property holds.
26    ///
27    /// # Panics
28    ///
29    /// Panics if the sparsity pattern is not square.
30    pub fn factor(pattern: SparsityPattern) -> Self {
31        assert_eq!(
32            pattern.major_dim(),
33            pattern.minor_dim(),
34            "Major and minor dimensions must be the same (square matrix)."
35        );
36        let (l_pattern, u_pattern) = nonzero_pattern(&pattern);
37        Self {
38            m_pattern: pattern,
39            l_pattern,
40            u_pattern,
41        }
42    }
43
44    /// The pattern of the Cholesky factor `L`.
45    #[must_use]
46    pub fn l_pattern(&self) -> &SparsityPattern {
47        &self.l_pattern
48    }
49}
50
51/// A sparse Cholesky factorization `A = L L^T` of a [`CscMatrix`].
52///
53/// The factor `L` is a sparse, lower-triangular matrix. See the article on [Wikipedia] for
54/// more information.
55///
56/// The implementation is a port of the `CsCholesky` implementation in `nalgebra`. It is similar
57/// to Tim Davis' [`CSparse`]. The current implementation performs no fill-in reduction, and can
58/// therefore be expected to produce much too dense Cholesky factors for many matrices.
59/// It is therefore not currently recommended to use this implementation for serious projects.
60///
61/// [`CSparse`]: https://epubs.siam.org/doi/book/10.1137/1.9780898718881
62/// [Wikipedia]: https://en.wikipedia.org/wiki/Cholesky_decomposition
63// TODO: We should probably implement PartialEq/Eq, but in that case we'd probably need a
64// custom implementation, due to the need to exclude the workspace arrays
65#[derive(Debug, Clone)]
66pub struct CscCholesky<T> {
67    // Pattern of the original matrix
68    m_pattern: SparsityPattern,
69    l_factor: CscMatrix<T>,
70    u_pattern: SparsityPattern,
71    work_x: Vec<T>,
72    work_c: Vec<usize>,
73}
74
75#[derive(Debug, PartialEq, Eq, Copy, Clone)]
76#[non_exhaustive]
77/// Possible errors produced by the Cholesky factorization.
78pub enum CholeskyError {
79    /// The matrix is not positive definite.
80    NotPositiveDefinite,
81}
82
83impl Display for CholeskyError {
84    fn fmt(&self, f: &mut Formatter<'_>) -> std::fmt::Result {
85        write!(f, "Matrix is not positive definite")
86    }
87}
88
89impl std::error::Error for CholeskyError {}
90
91impl<T: RealField> CscCholesky<T> {
92    /// Computes the numerical Cholesky factorization associated with the given
93    /// symbolic factorization and the provided values.
94    ///
95    /// The values correspond to the non-zero values of the CSC matrix for which the
96    /// symbolic factorization was computed.
97    ///
98    /// # Errors
99    ///
100    /// Returns an error if the numerical factorization fails. This can occur if the matrix is not
101    /// symmetric positive definite.
102    ///
103    /// # Panics
104    ///
105    /// Panics if the number of values differ from the number of non-zeros of the sparsity pattern
106    /// of the matrix that was symbolically factored.
107    pub fn factor_numerical(
108        symbolic: CscSymbolicCholesky,
109        values: &[T],
110    ) -> Result<Self, CholeskyError> {
111        assert_eq!(
112            symbolic.l_pattern.nnz(),
113            symbolic.u_pattern.nnz(),
114            "u is just the transpose of l, so should have the same nnz"
115        );
116
117        let l_nnz = symbolic.l_pattern.nnz();
118        let l_values = vec![T::zero(); l_nnz];
119        let l_factor =
120            CscMatrix::try_from_pattern_and_values(symbolic.l_pattern, l_values).unwrap();
121
122        let (nrows, ncols) = (l_factor.nrows(), l_factor.ncols());
123
124        let mut factorization = CscCholesky {
125            m_pattern: symbolic.m_pattern,
126            l_factor,
127            u_pattern: symbolic.u_pattern,
128            work_x: vec![T::zero(); nrows],
129            // Fill with MAX so that things hopefully totally fail if values are not
130            // overwritten. Might be easier to debug this way
131            work_c: vec![usize::MAX, ncols],
132        };
133
134        factorization.refactor(values)?;
135        Ok(factorization)
136    }
137
138    /// Computes the Cholesky factorization of the provided matrix.
139    ///
140    /// The matrix must be symmetric positive definite. Symmetry is not checked, and it is up
141    /// to the user to enforce this property.
142    ///
143    /// # Errors
144    ///
145    /// Returns an error if the numerical factorization fails. This can occur if the matrix is not
146    /// symmetric positive definite.
147    ///
148    /// # Panics
149    ///
150    /// Panics if the matrix is not square.
151    pub fn factor(matrix: &CscMatrix<T>) -> Result<Self, CholeskyError> {
152        let symbolic = CscSymbolicCholesky::factor(matrix.pattern().clone());
153        Self::factor_numerical(symbolic, matrix.values())
154    }
155
156    /// Re-computes the factorization for a new set of non-zero values.
157    ///
158    /// This is useful when the values of a matrix changes, but the sparsity pattern remains
159    /// constant.
160    ///
161    /// # Errors
162    ///
163    /// Returns an error if the numerical factorization fails. This can occur if the matrix is not
164    /// symmetric positive definite.
165    ///
166    /// # Panics
167    ///
168    /// Panics if the number of values does not match the number of non-zeros in the sparsity
169    /// pattern.
170    pub fn refactor(&mut self, values: &[T]) -> Result<(), CholeskyError> {
171        self.decompose_left_looking(values)
172    }
173
174    /// Returns a reference to the Cholesky factor `L`.
175    #[must_use]
176    pub fn l(&self) -> &CscMatrix<T> {
177        &self.l_factor
178    }
179
180    /// Returns the Cholesky factor `L`.
181    pub fn take_l(self) -> CscMatrix<T> {
182        self.l_factor
183    }
184
185    /// Perform a numerical left-looking cholesky decomposition of a matrix with the same structure as the
186    /// one used to initialize `self`, but with different non-zero values provided by `values`.
187    fn decompose_left_looking(&mut self, values: &[T]) -> Result<(), CholeskyError> {
188        assert!(
189            values.len() >= self.m_pattern.nnz(),
190            // TODO: Improve error message
191            "The set of values is too small."
192        );
193
194        let n = self.l_factor.nrows();
195
196        // Reset `work_c` to the column pointers of `l`.
197        self.work_c.clear();
198        self.work_c.extend_from_slice(self.l_factor.col_offsets());
199
200        unsafe {
201            for k in 0..n {
202                // Scatter the k-th column of the original matrix with the values provided.
203                let range_begin = *self.m_pattern.major_offsets().get_unchecked(k);
204                let range_end = *self.m_pattern.major_offsets().get_unchecked(k + 1);
205                let range_k = range_begin..range_end;
206
207                *self.work_x.get_unchecked_mut(k) = T::zero();
208                for p in range_k.clone() {
209                    let irow = *self.m_pattern.minor_indices().get_unchecked(p);
210
211                    if irow >= k {
212                        *self.work_x.get_unchecked_mut(irow) = values.get_unchecked(p).clone();
213                    }
214                }
215
216                for &j in self.u_pattern.lane(k) {
217                    let factor = -self
218                        .l_factor
219                        .values()
220                        .get_unchecked(*self.work_c.get_unchecked(j))
221                        .clone();
222                    *self.work_c.get_unchecked_mut(j) += 1;
223
224                    if j < k {
225                        let col_j = self.l_factor.col(j);
226                        let col_j_entries = col_j.row_indices().iter().zip(col_j.values());
227                        for (&z, val) in col_j_entries {
228                            if z >= k {
229                                *self.work_x.get_unchecked_mut(z) += val.clone() * factor.clone();
230                            }
231                        }
232                    }
233                }
234
235                let diag = self.work_x.get_unchecked(k).clone();
236
237                if diag > T::zero() {
238                    let denom = diag.sqrt();
239
240                    {
241                        let (offsets, _, values) = self.l_factor.csc_data_mut();
242                        *values.get_unchecked_mut(*offsets.get_unchecked(k)) = denom.clone();
243                    }
244
245                    let mut col_k = self.l_factor.col_mut(k);
246                    let (col_k_rows, col_k_values) = col_k.rows_and_values_mut();
247                    let col_k_entries = col_k_rows.iter().zip(col_k_values);
248                    for (&p, val) in col_k_entries {
249                        *val = self.work_x.get_unchecked(p).clone() / denom.clone();
250                        *self.work_x.get_unchecked_mut(p) = T::zero();
251                    }
252                } else {
253                    return Err(CholeskyError::NotPositiveDefinite);
254                }
255            }
256        }
257
258        Ok(())
259    }
260
261    /// Solves the system `A X = B`, where `X` and `B` are dense matrices.
262    ///
263    /// # Panics
264    ///
265    /// Panics if `B` is not square.
266    #[must_use = "Did you mean to use solve_mut()?"]
267    pub fn solve<'a>(&'a self, b: impl Into<DMatrixView<'a, T>>) -> DMatrix<T> {
268        let b = b.into();
269        let mut output = b.clone_owned();
270        self.solve_mut(&mut output);
271        output
272    }
273
274    /// Solves the system `AX = B`, where `X` and `B` are dense matrices.
275    ///
276    /// The result is stored in-place in `b`.
277    ///
278    /// # Panics
279    ///
280    /// Panics if `b` is not square.
281    pub fn solve_mut<'a>(&'a self, b: impl Into<DMatrixViewMut<'a, T>>) {
282        let expect_msg = "If the Cholesky factorization succeeded,\
283            then the triangular solve should never fail";
284        // Solve LY = B
285        let mut y = b.into();
286        spsolve_csc_lower_triangular(Op::NoOp(self.l()), &mut y).expect(expect_msg);
287
288        // Solve L^T X = Y
289        let mut x = y;
290        spsolve_csc_lower_triangular(Op::Transpose(self.l()), &mut x).expect(expect_msg);
291    }
292}
293
294fn reach(
295    pattern: &SparsityPattern,
296    j: usize,
297    max_j: usize,
298    tree: &[usize],
299    marks: &mut Vec<bool>,
300    out: &mut Vec<usize>,
301) {
302    marks.clear();
303    marks.resize(tree.len(), false);
304
305    // TODO: avoid all those allocations.
306    let mut tmp = Vec::new();
307    let mut res = Vec::new();
308
309    for &irow in pattern.lane(j) {
310        let mut curr = irow;
311        while curr != usize::MAX && curr <= max_j && !marks[curr] {
312            marks[curr] = true;
313            tmp.push(curr);
314            curr = tree[curr];
315        }
316
317        tmp.append(&mut res);
318        mem::swap(&mut tmp, &mut res);
319    }
320
321    res.sort_unstable();
322
323    out.append(&mut res);
324}
325
326fn nonzero_pattern(m: &SparsityPattern) -> (SparsityPattern, SparsityPattern) {
327    let etree = elimination_tree(m);
328    // Note: We assume CSC, therefore rows == minor and cols == major
329    let (nrows, ncols) = (m.minor_dim(), m.major_dim());
330    let mut rows = Vec::with_capacity(m.nnz());
331    let mut col_offsets = Vec::with_capacity(ncols + 1);
332    let mut marks = Vec::new();
333
334    // NOTE: the following will actually compute the non-zero pattern of
335    // the transpose of l.
336    col_offsets.push(0);
337    for i in 0..nrows {
338        reach(m, i, i, &etree, &mut marks, &mut rows);
339        col_offsets.push(rows.len());
340    }
341
342    let u_pattern =
343        SparsityPattern::try_from_offsets_and_indices(nrows, ncols, col_offsets, rows).unwrap();
344
345    // TODO: Avoid this transpose?
346    let l_pattern = u_pattern.transpose();
347
348    (l_pattern, u_pattern)
349}
350
351fn elimination_tree(pattern: &SparsityPattern) -> Vec<usize> {
352    // Note: The pattern is assumed to of a CSC matrix, so the number of rows is
353    // given by the minor dimension
354    let nrows = pattern.minor_dim();
355    let mut forest: Vec<_> = iter::repeat(usize::MAX).take(nrows).collect();
356    let mut ancestor: Vec<_> = iter::repeat(usize::MAX).take(nrows).collect();
357
358    for k in 0..nrows {
359        for &irow in pattern.lane(k) {
360            let mut i = irow;
361
362            while i < k {
363                let i_ancestor = ancestor[i];
364                ancestor[i] = k;
365
366                if i_ancestor == usize::MAX {
367                    forest[i] = k;
368                    break;
369                }
370
371                i = i_ancestor;
372            }
373        }
374    }
375
376    forest
377}