Skip to content

Commit e94e03d

Browse files
committed
Document for Cholesky decomposition
1 parent 2625168 commit e94e03d

File tree

2 files changed

+56
-8
lines changed

2 files changed

+56
-8
lines changed

lax/src/cholesky.rs

Lines changed: 54 additions & 7 deletions
Original file line numberDiff line numberDiff line change
@@ -1,21 +1,68 @@
1-
//! Cholesky decomposition
2-
31
use super::*;
42
use crate::{error::*, layout::*};
53
use cauchy::*;
64

5+
#[cfg_attr(doc, katexit::katexit)]
6+
/// Solve symmetric/hermite positive-definite linear equations using Cholesky decomposition
7+
///
8+
/// For a given positive definite matrix $A$,
9+
/// Cholesky decomposition is described as $A = U^T U$ or $A = LL^T$ where
10+
///
11+
/// - $L$ is lower matrix
12+
/// - $U$ is upper matrix
13+
///
14+
/// This is designed as two step computation according to LAPACK API
15+
///
16+
/// 1. Factorize input matrix $A$ into $L$ or $U$
17+
/// 2. Solve linear equation $Ax = b$ or compute inverse matrix $A^{-1}$
18+
/// using $U$ or $L$.
719
pub trait Cholesky_: Sized {
8-
/// Cholesky: wrapper of `*potrf`
20+
/// Compute Cholesky decomposition $A = U^T U$ or $A = L L^T$ according to [UPLO]
21+
///
22+
/// LAPACK correspondance
23+
/// ----------------------
24+
///
25+
/// | f32 | f64 | c32 | c64 |
26+
/// |:---------|:---------|:---------|:---------|
27+
/// | [spotrf] | [dpotrf] | [cpotrf] | [zpotrf] |
28+
///
29+
/// [spotrf]: https://netlib.org/lapack/explore-html/d8/db2/group__real_p_ocomputational_gaaf31db7ab15b4f4ba527a3d31a15a58e.html
30+
/// [dpotrf]: https://netlib.org/lapack/explore-html/d1/d7a/group__double_p_ocomputational_ga2f55f604a6003d03b5cd4a0adcfb74d6.html
31+
/// [cpotrf]: https://netlib.org/lapack/explore-html/d8/d6c/group__variants_p_ocomputational_ga4e85f48dbd837ccbbf76aa077f33de19.html
32+
/// [zpotrf]: https://netlib.org/lapack/explore-html/d3/d8d/group__complex16_p_ocomputational_ga93e22b682170873efb50df5a79c5e4eb.html
933
///
10-
/// **Warning: Only the portion of `a` corresponding to `UPLO` is written.**
1134
fn cholesky(l: MatrixLayout, uplo: UPLO, a: &mut [Self]) -> Result<()>;
1235

13-
/// Wrapper of `*potri`
36+
/// Compute inverse matrix $A^{-1}$ using $U$ or $L$
37+
///
38+
/// LAPACK correspondance
39+
/// ----------------------
40+
///
41+
/// | f32 | f64 | c32 | c64 |
42+
/// |:---------|:---------|:---------|:---------|
43+
/// | [spotri] | [dpotri] | [cpotri] | [zpotri] |
44+
///
45+
/// [spotri]: https://netlib.org/lapack/explore-html/d8/db2/group__real_p_ocomputational_ga4c381894bb34b1583fcc0dceafc5bea1.html
46+
/// [dpotri]: https://netlib.org/lapack/explore-html/d1/d7a/group__double_p_ocomputational_ga9dfc04beae56a3b1c1f75eebc838c14c.html
47+
/// [cpotri]: https://netlib.org/lapack/explore-html/d6/df6/group__complex_p_ocomputational_ga52b8da4d314abefaee93dd5c1ed7739e.html
48+
/// [zpotri]: https://netlib.org/lapack/explore-html/d3/d8d/group__complex16_p_ocomputational_gaf37e3b8bbacd3332e83ffb3f1018bcf1.html
1449
///
15-
/// **Warning: Only the portion of `a` corresponding to `UPLO` is written.**
1650
fn inv_cholesky(l: MatrixLayout, uplo: UPLO, a: &mut [Self]) -> Result<()>;
1751

18-
/// Wrapper of `*potrs`
52+
/// Solve linear equation $Ax = b$ using $U$ or $L$
53+
///
54+
/// LAPACK correspondance
55+
/// ----------------------
56+
///
57+
/// | f32 | f64 | c32 | c64 |
58+
/// |:---------|:---------|:---------|:---------|
59+
/// | [spotrs] | [dpotrs] | [cpotrs] | [zpotrs] |
60+
///
61+
/// [spotrs]: https://netlib.org/lapack/explore-html/d8/db2/group__real_p_ocomputational_gaf5cc1531aa5ffe706533fbca343d55dd.html
62+
/// [dpotrs]: https://netlib.org/lapack/explore-html/d1/d7a/group__double_p_ocomputational_ga167aa0166c4ce726385f65e4ab05e7c1.html
63+
/// [cpotrs]: https://netlib.org/lapack/explore-html/d6/df6/group__complex_p_ocomputational_gad9052b4b70569dfd6e8943971c9b38b2.html
64+
/// [zpotrs]: https://netlib.org/lapack/explore-html/d3/d8d/group__complex16_p_ocomputational_gaa2116ea574b01efda584dff0b74c9fcd.html
65+
///
1966
fn solve_cholesky(l: MatrixLayout, uplo: UPLO, a: &[Self], b: &mut [Self]) -> Result<()>;
2067
}
2168

lax/src/lib.rs

Lines changed: 2 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -49,7 +49,8 @@
4949
//! According to the property input metrix, several types of triangular factorization are used:
5050
//!
5151
//! - [Solve_] trait provides methods for LU-decomposition for general matrix.
52-
//! - [Solveh_] triat provides methods for Bunch-Kaufman diagonal pivoting method for symmetric/hermite indefinite matrix
52+
//! - [Solveh_] triat provides methods for Bunch-Kaufman diagonal pivoting method for symmetric/hermite indefinite matrix.
53+
//! - [Cholesky_] triat provides methods for Cholesky decomposition for symmetric/hermite positive dinite matrix.
5354
//!
5455
//! Eigenvalue Problem
5556
//! -------------------

0 commit comments

Comments
 (0)