cholesky.go 19 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698
  1. // Copyright ©2013 The Gonum Authors. All rights reserved.
  2. // Use of this source code is governed by a BSD-style
  3. // license that can be found in the LICENSE file.
  4. package mat
  5. import (
  6. "math"
  7. "gonum.org/v1/gonum/blas"
  8. "gonum.org/v1/gonum/blas/blas64"
  9. "gonum.org/v1/gonum/lapack/lapack64"
  10. )
  11. const (
  12. badTriangle = "mat: invalid triangle"
  13. badCholesky = "mat: invalid Cholesky factorization"
  14. )
  15. var (
  16. _ Matrix = (*Cholesky)(nil)
  17. _ Symmetric = (*Cholesky)(nil)
  18. )
  19. // Cholesky is a symmetric positive definite matrix represented by its
  20. // Cholesky decomposition.
  21. //
  22. // The decomposition can be constructed using the Factorize method. The
  23. // factorization itself can be extracted using the UTo or LTo methods, and the
  24. // original symmetric matrix can be recovered with ToSym.
  25. //
  26. // Note that this matrix representation is useful for certain operations, in
  27. // particular finding solutions to linear equations. It is very inefficient
  28. // at other operations, in particular At is slow.
  29. //
  30. // Cholesky methods may only be called on a value that has been successfully
  31. // initialized by a call to Factorize that has returned true. Calls to methods
  32. // of an unsuccessful Cholesky factorization will panic.
  33. type Cholesky struct {
  34. // The chol pointer must never be retained as a pointer outside the Cholesky
  35. // struct, either by returning chol outside the struct or by setting it to
  36. // a pointer coming from outside. The same prohibition applies to the data
  37. // slice within chol.
  38. chol *TriDense
  39. cond float64
  40. }
  41. // updateCond updates the condition number of the Cholesky decomposition. If
  42. // norm > 0, then that norm is used as the norm of the original matrix A, otherwise
  43. // the norm is estimated from the decomposition.
  44. func (c *Cholesky) updateCond(norm float64) {
  45. n := c.chol.mat.N
  46. work := getFloats(3*n, false)
  47. defer putFloats(work)
  48. if norm < 0 {
  49. // This is an approximation. By the definition of a norm,
  50. // |AB| <= |A| |B|.
  51. // Since A = Uᵀ*U, we get for the condition number κ that
  52. // κ(A) := |A| |A^-1| = |Uᵀ*U| |A^-1| <= |Uᵀ| |U| |A^-1|,
  53. // so this will overestimate the condition number somewhat.
  54. // The norm of the original factorized matrix cannot be stored
  55. // because of update possibilities.
  56. unorm := lapack64.Lantr(CondNorm, c.chol.mat, work)
  57. lnorm := lapack64.Lantr(CondNormTrans, c.chol.mat, work)
  58. norm = unorm * lnorm
  59. }
  60. sym := c.chol.asSymBlas()
  61. iwork := getInts(n, false)
  62. v := lapack64.Pocon(sym, norm, work, iwork)
  63. putInts(iwork)
  64. c.cond = 1 / v
  65. }
  66. // Dims returns the dimensions of the matrix.
  67. func (ch *Cholesky) Dims() (r, c int) {
  68. if !ch.valid() {
  69. panic(badCholesky)
  70. }
  71. r, c = ch.chol.Dims()
  72. return r, c
  73. }
  74. // At returns the element at row i, column j.
  75. func (c *Cholesky) At(i, j int) float64 {
  76. if !c.valid() {
  77. panic(badCholesky)
  78. }
  79. n := c.Symmetric()
  80. if uint(i) >= uint(n) {
  81. panic(ErrRowAccess)
  82. }
  83. if uint(j) >= uint(n) {
  84. panic(ErrColAccess)
  85. }
  86. var val float64
  87. for k := 0; k <= min(i, j); k++ {
  88. val += c.chol.at(k, i) * c.chol.at(k, j)
  89. }
  90. return val
  91. }
  92. // T returns the the receiver, the transpose of a symmetric matrix.
  93. func (c *Cholesky) T() Matrix {
  94. return c
  95. }
  96. // Symmetric implements the Symmetric interface and returns the number of rows
  97. // in the matrix (this is also the number of columns).
  98. func (c *Cholesky) Symmetric() int {
  99. r, _ := c.chol.Dims()
  100. return r
  101. }
  102. // Cond returns the condition number of the factorized matrix.
  103. func (c *Cholesky) Cond() float64 {
  104. if !c.valid() {
  105. panic(badCholesky)
  106. }
  107. return c.cond
  108. }
  109. // Factorize calculates the Cholesky decomposition of the matrix A and returns
  110. // whether the matrix is positive definite. If Factorize returns false, the
  111. // factorization must not be used.
  112. func (c *Cholesky) Factorize(a Symmetric) (ok bool) {
  113. n := a.Symmetric()
  114. if c.chol == nil {
  115. c.chol = NewTriDense(n, Upper, nil)
  116. } else {
  117. c.chol = NewTriDense(n, Upper, use(c.chol.mat.Data, n*n))
  118. }
  119. copySymIntoTriangle(c.chol, a)
  120. sym := c.chol.asSymBlas()
  121. work := getFloats(c.chol.mat.N, false)
  122. norm := lapack64.Lansy(CondNorm, sym, work)
  123. putFloats(work)
  124. _, ok = lapack64.Potrf(sym)
  125. if ok {
  126. c.updateCond(norm)
  127. } else {
  128. c.Reset()
  129. }
  130. return ok
  131. }
  132. // Reset resets the factorization so that it can be reused as the receiver of a
  133. // dimensionally restricted operation.
  134. func (c *Cholesky) Reset() {
  135. if c.chol != nil {
  136. c.chol.Reset()
  137. }
  138. c.cond = math.Inf(1)
  139. }
  140. // IsEmpty returns whether the receiver is empty. Empty matrices can be the
  141. // receiver for size-restricted operations. The receiver can be emptied using
  142. // Reset.
  143. func (c *Cholesky) IsEmpty() bool {
  144. return c.chol == nil || c.chol.IsEmpty()
  145. }
  146. // SetFromU sets the Cholesky decomposition from the given triangular matrix.
  147. // SetFromU panics if t is not upper triangular. If the receiver is empty it
  148. // is resized to be n×n, the size of t. If dst is non-empty, SetFromU panics
  149. // if c is not of size n×n. Note that t is copied into, not stored inside, the
  150. // receiver.
  151. func (c *Cholesky) SetFromU(t Triangular) {
  152. n, kind := t.Triangle()
  153. if kind != Upper {
  154. panic("cholesky: matrix must be upper triangular")
  155. }
  156. if c.chol == nil {
  157. c.chol = NewTriDense(n, Upper, nil)
  158. } else {
  159. c.chol.reuseAsNonZeroed(n, Upper)
  160. }
  161. c.chol.Copy(t)
  162. c.updateCond(-1)
  163. }
  164. // Clone makes a copy of the input Cholesky into the receiver, overwriting the
  165. // previous value of the receiver. Clone does not place any restrictions on receiver
  166. // shape. Clone panics if the input Cholesky is not the result of a valid decomposition.
  167. func (c *Cholesky) Clone(chol *Cholesky) {
  168. if !chol.valid() {
  169. panic(badCholesky)
  170. }
  171. n := chol.Symmetric()
  172. if c.chol == nil {
  173. c.chol = NewTriDense(n, Upper, nil)
  174. } else {
  175. c.chol = NewTriDense(n, Upper, use(c.chol.mat.Data, n*n))
  176. }
  177. c.chol.Copy(chol.chol)
  178. c.cond = chol.cond
  179. }
  180. // Det returns the determinant of the matrix that has been factorized.
  181. func (c *Cholesky) Det() float64 {
  182. if !c.valid() {
  183. panic(badCholesky)
  184. }
  185. return math.Exp(c.LogDet())
  186. }
  187. // LogDet returns the log of the determinant of the matrix that has been factorized.
  188. func (c *Cholesky) LogDet() float64 {
  189. if !c.valid() {
  190. panic(badCholesky)
  191. }
  192. var det float64
  193. for i := 0; i < c.chol.mat.N; i++ {
  194. det += 2 * math.Log(c.chol.mat.Data[i*c.chol.mat.Stride+i])
  195. }
  196. return det
  197. }
  198. // SolveTo finds the matrix X that solves A * X = B where A is represented
  199. // by the Cholesky decomposition. The result is stored in-place into dst.
  200. func (c *Cholesky) SolveTo(dst *Dense, b Matrix) error {
  201. if !c.valid() {
  202. panic(badCholesky)
  203. }
  204. n := c.chol.mat.N
  205. bm, bn := b.Dims()
  206. if n != bm {
  207. panic(ErrShape)
  208. }
  209. dst.reuseAsNonZeroed(bm, bn)
  210. if b != dst {
  211. dst.Copy(b)
  212. }
  213. lapack64.Potrs(c.chol.mat, dst.mat)
  214. if c.cond > ConditionTolerance {
  215. return Condition(c.cond)
  216. }
  217. return nil
  218. }
  219. // SolveCholTo finds the matrix X that solves A * X = B where A and B are represented
  220. // by their Cholesky decompositions a and b. The result is stored in-place into
  221. // dst.
  222. func (a *Cholesky) SolveCholTo(dst *Dense, b *Cholesky) error {
  223. if !a.valid() || !b.valid() {
  224. panic(badCholesky)
  225. }
  226. bn := b.chol.mat.N
  227. if a.chol.mat.N != bn {
  228. panic(ErrShape)
  229. }
  230. dst.reuseAsZeroed(bn, bn)
  231. dst.Copy(b.chol.T())
  232. blas64.Trsm(blas.Left, blas.Trans, 1, a.chol.mat, dst.mat)
  233. blas64.Trsm(blas.Left, blas.NoTrans, 1, a.chol.mat, dst.mat)
  234. blas64.Trmm(blas.Right, blas.NoTrans, 1, b.chol.mat, dst.mat)
  235. if a.cond > ConditionTolerance {
  236. return Condition(a.cond)
  237. }
  238. return nil
  239. }
  240. // SolveVecTo finds the vector X that solves A * x = b where A is represented
  241. // by the Cholesky decomposition. The result is stored in-place into
  242. // dst.
  243. func (c *Cholesky) SolveVecTo(dst *VecDense, b Vector) error {
  244. if !c.valid() {
  245. panic(badCholesky)
  246. }
  247. n := c.chol.mat.N
  248. if br, bc := b.Dims(); br != n || bc != 1 {
  249. panic(ErrShape)
  250. }
  251. switch rv := b.(type) {
  252. default:
  253. dst.reuseAsNonZeroed(n)
  254. return c.SolveTo(dst.asDense(), b)
  255. case RawVectorer:
  256. bmat := rv.RawVector()
  257. if dst != b {
  258. dst.checkOverlap(bmat)
  259. }
  260. dst.reuseAsNonZeroed(n)
  261. if dst != b {
  262. dst.CopyVec(b)
  263. }
  264. lapack64.Potrs(c.chol.mat, dst.asGeneral())
  265. if c.cond > ConditionTolerance {
  266. return Condition(c.cond)
  267. }
  268. return nil
  269. }
  270. }
  271. // RawU returns the Triangular matrix used to store the Cholesky decomposition of
  272. // the original matrix A. The returned matrix should not be modified. If it is
  273. // modified, the decomposition is invalid and should not be used.
  274. func (c *Cholesky) RawU() Triangular {
  275. return c.chol
  276. }
  277. // UTo stores into dst the n×n upper triangular matrix U from a Cholesky
  278. // decomposition
  279. // A = Uᵀ * U.
  280. // If dst is empty, it is resized to be an n×n upper triangular matrix. When dst
  281. // is non-empty, UTo panics if dst is not n×n or not Upper. UTo will also panic
  282. // if the receiver does not contain a successful factorization.
  283. func (c *Cholesky) UTo(dst *TriDense) {
  284. if !c.valid() {
  285. panic(badCholesky)
  286. }
  287. n := c.chol.mat.N
  288. if dst.IsEmpty() {
  289. dst.ReuseAsTri(n, Upper)
  290. } else {
  291. n2, kind := dst.Triangle()
  292. if n != n2 {
  293. panic(ErrShape)
  294. }
  295. if kind != Upper {
  296. panic(ErrTriangle)
  297. }
  298. }
  299. dst.Copy(c.chol)
  300. }
  301. // LTo stores into dst the n×n lower triangular matrix L from a Cholesky
  302. // decomposition
  303. // A = L * Lᵀ.
  304. // If dst is empty, it is resized to be an n×n lower triangular matrix. When dst
  305. // is non-empty, LTo panics if dst is not n×n or not Lower. LTo will also panic
  306. // if the receiver does not contain a successful factorization.
  307. func (c *Cholesky) LTo(dst *TriDense) {
  308. if !c.valid() {
  309. panic(badCholesky)
  310. }
  311. n := c.chol.mat.N
  312. if dst.IsEmpty() {
  313. dst.ReuseAsTri(n, Lower)
  314. } else {
  315. n2, kind := dst.Triangle()
  316. if n != n2 {
  317. panic(ErrShape)
  318. }
  319. if kind != Lower {
  320. panic(ErrTriangle)
  321. }
  322. }
  323. dst.Copy(c.chol.TTri())
  324. }
  325. // ToSym reconstructs the original positive definite matrix from its
  326. // Cholesky decomposition, storing the result into dst. If dst is
  327. // empty it is resized to be n×n. If dst is non-empty, ToSym panics
  328. // if dst is not of size n×n. ToSym will also panic if the receiver
  329. // does not contain a successful factorization.
  330. func (c *Cholesky) ToSym(dst *SymDense) {
  331. if !c.valid() {
  332. panic(badCholesky)
  333. }
  334. n := c.chol.mat.N
  335. if dst.IsEmpty() {
  336. dst.ReuseAsSym(n)
  337. } else {
  338. n2 := dst.Symmetric()
  339. if n != n2 {
  340. panic(ErrShape)
  341. }
  342. }
  343. // Create a TriDense representing the Cholesky factor U with dst's
  344. // backing slice.
  345. // Operations on u are reflected in s.
  346. u := &TriDense{
  347. mat: blas64.Triangular{
  348. Uplo: blas.Upper,
  349. Diag: blas.NonUnit,
  350. N: n,
  351. Data: dst.mat.Data,
  352. Stride: dst.mat.Stride,
  353. },
  354. cap: n,
  355. }
  356. u.Copy(c.chol)
  357. // Compute the product Uᵀ*U using the algorithm from LAPACK/TESTING/LIN/dpot01.f
  358. a := u.mat.Data
  359. lda := u.mat.Stride
  360. bi := blas64.Implementation()
  361. for k := n - 1; k >= 0; k-- {
  362. a[k*lda+k] = bi.Ddot(k+1, a[k:], lda, a[k:], lda)
  363. if k > 0 {
  364. bi.Dtrmv(blas.Upper, blas.Trans, blas.NonUnit, k, a, lda, a[k:], lda)
  365. }
  366. }
  367. }
  368. // InverseTo computes the inverse of the matrix represented by its Cholesky
  369. // factorization and stores the result into s. If the factorized
  370. // matrix is ill-conditioned, a Condition error will be returned.
  371. // Note that matrix inversion is numerically unstable, and should generally be
  372. // avoided where possible, for example by using the Solve routines.
  373. func (c *Cholesky) InverseTo(s *SymDense) error {
  374. if !c.valid() {
  375. panic(badCholesky)
  376. }
  377. s.reuseAsNonZeroed(c.chol.mat.N)
  378. // Create a TriDense representing the Cholesky factor U with the backing
  379. // slice from s.
  380. // Operations on u are reflected in s.
  381. u := &TriDense{
  382. mat: blas64.Triangular{
  383. Uplo: blas.Upper,
  384. Diag: blas.NonUnit,
  385. N: s.mat.N,
  386. Data: s.mat.Data,
  387. Stride: s.mat.Stride,
  388. },
  389. cap: s.mat.N,
  390. }
  391. u.Copy(c.chol)
  392. _, ok := lapack64.Potri(u.mat)
  393. if !ok {
  394. return Condition(math.Inf(1))
  395. }
  396. if c.cond > ConditionTolerance {
  397. return Condition(c.cond)
  398. }
  399. return nil
  400. }
  401. // Scale multiplies the original matrix A by a positive constant using
  402. // its Cholesky decomposition, storing the result in-place into the receiver.
  403. // That is, if the original Cholesky factorization is
  404. // Uᵀ * U = A
  405. // the updated factorization is
  406. // U'ᵀ * U' = f A = A'
  407. // Scale panics if the constant is non-positive, or if the receiver is non-empty
  408. // and is of a different size from the input.
  409. func (c *Cholesky) Scale(f float64, orig *Cholesky) {
  410. if !orig.valid() {
  411. panic(badCholesky)
  412. }
  413. if f <= 0 {
  414. panic("cholesky: scaling by a non-positive constant")
  415. }
  416. n := orig.Symmetric()
  417. if c.chol == nil {
  418. c.chol = NewTriDense(n, Upper, nil)
  419. } else if c.chol.mat.N != n {
  420. panic(ErrShape)
  421. }
  422. c.chol.ScaleTri(math.Sqrt(f), orig.chol)
  423. c.cond = orig.cond // Scaling by a positive constant does not change the condition number.
  424. }
  425. // ExtendVecSym computes the Cholesky decomposition of the original matrix A,
  426. // whose Cholesky decomposition is in a, extended by a the n×1 vector v according to
  427. // [A w]
  428. // [w' k]
  429. // where k = v[n-1] and w = v[:n-1]. The result is stored into the receiver.
  430. // In order for the updated matrix to be positive definite, it must be the case
  431. // that k > w' A^-1 w. If this condition does not hold then ExtendVecSym will
  432. // return false and the receiver will not be updated.
  433. //
  434. // ExtendVecSym will panic if v.Len() != a.Symmetric()+1 or if a does not contain
  435. // a valid decomposition.
  436. func (c *Cholesky) ExtendVecSym(a *Cholesky, v Vector) (ok bool) {
  437. n := a.Symmetric()
  438. if v.Len() != n+1 {
  439. panic(badSliceLength)
  440. }
  441. if !a.valid() {
  442. panic(badCholesky)
  443. }
  444. // The algorithm is commented here, but see also
  445. // https://math.stackexchange.com/questions/955874/cholesky-factor-when-adding-a-row-and-column-to-already-factorized-matrix
  446. // We have A and want to compute the Cholesky of
  447. // [A w]
  448. // [w' k]
  449. // We want
  450. // [U c]
  451. // [0 d]
  452. // to be the updated Cholesky, and so it must be that
  453. // [A w] = [U' 0] [U c]
  454. // [w' k] [c' d] [0 d]
  455. // Thus, we need
  456. // 1) A = U'U (true by the original decomposition being valid),
  457. // 2) U' * c = w => c = U'^-1 w
  458. // 3) c'*c + d'*d = k => d = sqrt(k-c'*c)
  459. // First, compute c = U'^-1 a
  460. w := NewVecDense(n, nil)
  461. w.CopyVec(v)
  462. k := v.At(n, 0)
  463. var t VecDense
  464. _ = t.SolveVec(a.chol.T(), w)
  465. dot := Dot(&t, &t)
  466. if dot >= k {
  467. return false
  468. }
  469. d := math.Sqrt(k - dot)
  470. newU := NewTriDense(n+1, Upper, nil)
  471. newU.Copy(a.chol)
  472. for i := 0; i < n; i++ {
  473. newU.SetTri(i, n, t.At(i, 0))
  474. }
  475. newU.SetTri(n, n, d)
  476. c.chol = newU
  477. c.updateCond(-1)
  478. return true
  479. }
  480. // SymRankOne performs a rank-1 update of the original matrix A and refactorizes
  481. // its Cholesky factorization, storing the result into the receiver. That is, if
  482. // in the original Cholesky factorization
  483. // Uᵀ * U = A,
  484. // in the updated factorization
  485. // U'ᵀ * U' = A + alpha * x * xᵀ = A'.
  486. //
  487. // Note that when alpha is negative, the updating problem may be ill-conditioned
  488. // and the results may be inaccurate, or the updated matrix A' may not be
  489. // positive definite and not have a Cholesky factorization. SymRankOne returns
  490. // whether the updated matrix A' is positive definite. If the update fails
  491. // the receiver is left unchanged.
  492. //
  493. // SymRankOne updates a Cholesky factorization in O(n²) time. The Cholesky
  494. // factorization computation from scratch is O(n³).
  495. func (c *Cholesky) SymRankOne(orig *Cholesky, alpha float64, x Vector) (ok bool) {
  496. if !orig.valid() {
  497. panic(badCholesky)
  498. }
  499. n := orig.Symmetric()
  500. if r, c := x.Dims(); r != n || c != 1 {
  501. panic(ErrShape)
  502. }
  503. if orig != c {
  504. if c.chol == nil {
  505. c.chol = NewTriDense(n, Upper, nil)
  506. } else if c.chol.mat.N != n {
  507. panic(ErrShape)
  508. }
  509. c.chol.Copy(orig.chol)
  510. }
  511. if alpha == 0 {
  512. return true
  513. }
  514. // Algorithms for updating and downdating the Cholesky factorization are
  515. // described, for example, in
  516. // - J. J. Dongarra, J. R. Bunch, C. B. Moler, G. W. Stewart: LINPACK
  517. // Users' Guide. SIAM (1979), pages 10.10--10.14
  518. // or
  519. // - P. E. Gill, G. H. Golub, W. Murray, and M. A. Saunders: Methods for
  520. // modifying matrix factorizations. Mathematics of Computation 28(126)
  521. // (1974), Method C3 on page 521
  522. //
  523. // The implementation is based on LINPACK code
  524. // http://www.netlib.org/linpack/dchud.f
  525. // http://www.netlib.org/linpack/dchdd.f
  526. // and
  527. // https://icl.cs.utk.edu/lapack-forum/viewtopic.php?f=2&t=2646
  528. //
  529. // According to http://icl.cs.utk.edu/lapack-forum/archives/lapack/msg00301.html
  530. // LINPACK is released under BSD license.
  531. //
  532. // See also:
  533. // - M. A. Saunders: Large-scale Linear Programming Using the Cholesky
  534. // Factorization. Technical Report Stanford University (1972)
  535. // http://i.stanford.edu/pub/cstr/reports/cs/tr/72/252/CS-TR-72-252.pdf
  536. // - Matthias Seeger: Low rank updates for the Cholesky decomposition.
  537. // EPFL Technical Report 161468 (2004)
  538. // http://infoscience.epfl.ch/record/161468
  539. work := getFloats(n, false)
  540. defer putFloats(work)
  541. var xmat blas64.Vector
  542. if rv, ok := x.(RawVectorer); ok {
  543. xmat = rv.RawVector()
  544. } else {
  545. var tmp *VecDense
  546. tmp.CopyVec(x)
  547. xmat = tmp.RawVector()
  548. }
  549. blas64.Copy(xmat, blas64.Vector{N: n, Data: work, Inc: 1})
  550. if alpha > 0 {
  551. // Compute rank-1 update.
  552. if alpha != 1 {
  553. blas64.Scal(math.Sqrt(alpha), blas64.Vector{N: n, Data: work, Inc: 1})
  554. }
  555. umat := c.chol.mat
  556. stride := umat.Stride
  557. for i := 0; i < n; i++ {
  558. // Compute parameters of the Givens matrix that zeroes
  559. // the i-th element of x.
  560. c, s, r, _ := blas64.Rotg(umat.Data[i*stride+i], work[i])
  561. if r < 0 {
  562. // Multiply by -1 to have positive diagonal
  563. // elemnts.
  564. r *= -1
  565. c *= -1
  566. s *= -1
  567. }
  568. umat.Data[i*stride+i] = r
  569. if i < n-1 {
  570. // Multiply the extended factorization matrix by
  571. // the Givens matrix from the left. Only
  572. // the i-th row and x are modified.
  573. blas64.Rot(
  574. blas64.Vector{N: n - i - 1, Data: umat.Data[i*stride+i+1 : i*stride+n], Inc: 1},
  575. blas64.Vector{N: n - i - 1, Data: work[i+1 : n], Inc: 1},
  576. c, s)
  577. }
  578. }
  579. c.updateCond(-1)
  580. return true
  581. }
  582. // Compute rank-1 downdate.
  583. alpha = math.Sqrt(-alpha)
  584. if alpha != 1 {
  585. blas64.Scal(alpha, blas64.Vector{N: n, Data: work, Inc: 1})
  586. }
  587. // Solve Uᵀ * p = x storing the result into work.
  588. ok = lapack64.Trtrs(blas.Trans, c.chol.RawTriangular(), blas64.General{
  589. Rows: n,
  590. Cols: 1,
  591. Stride: 1,
  592. Data: work,
  593. })
  594. if !ok {
  595. // The original matrix is singular. Should not happen, because
  596. // the factorization is valid.
  597. panic(badCholesky)
  598. }
  599. norm := blas64.Nrm2(blas64.Vector{N: n, Data: work, Inc: 1})
  600. if norm >= 1 {
  601. // The updated matrix is not positive definite.
  602. return false
  603. }
  604. norm = math.Sqrt((1 + norm) * (1 - norm))
  605. cos := getFloats(n, false)
  606. defer putFloats(cos)
  607. sin := getFloats(n, false)
  608. defer putFloats(sin)
  609. for i := n - 1; i >= 0; i-- {
  610. // Compute parameters of Givens matrices that zero elements of p
  611. // backwards.
  612. cos[i], sin[i], norm, _ = blas64.Rotg(norm, work[i])
  613. if norm < 0 {
  614. norm *= -1
  615. cos[i] *= -1
  616. sin[i] *= -1
  617. }
  618. }
  619. workMat := getWorkspaceTri(c.chol.mat.N, c.chol.triKind(), false)
  620. defer putWorkspaceTri(workMat)
  621. workMat.Copy(c.chol)
  622. umat := workMat.mat
  623. stride := workMat.mat.Stride
  624. for i := n - 1; i >= 0; i-- {
  625. work[i] = 0
  626. // Apply Givens matrices to U.
  627. blas64.Rot(
  628. blas64.Vector{N: n - i, Data: work[i:n], Inc: 1},
  629. blas64.Vector{N: n - i, Data: umat.Data[i*stride+i : i*stride+n], Inc: 1},
  630. cos[i], sin[i])
  631. if umat.Data[i*stride+i] == 0 {
  632. // The matrix is singular (may rarely happen due to
  633. // floating-point effects?).
  634. ok = false
  635. } else if umat.Data[i*stride+i] < 0 {
  636. // Diagonal elements should be positive. If it happens
  637. // that on the i-th row the diagonal is negative,
  638. // multiply U from the left by an identity matrix that
  639. // has -1 on the i-th row.
  640. blas64.Scal(-1, blas64.Vector{N: n - i, Data: umat.Data[i*stride+i : i*stride+n], Inc: 1})
  641. }
  642. }
  643. if ok {
  644. c.chol.Copy(workMat)
  645. c.updateCond(-1)
  646. }
  647. return ok
  648. }
  649. func (c *Cholesky) valid() bool {
  650. return c.chol != nil && !c.chol.IsEmpty()
  651. }