cholesky.go 19 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674
  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^T*U, we get for the condition number κ that
  52. // κ(A) := |A| |A^-1| = |U^T*U| |A^-1| <= |U^T| |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. // SetFromU sets the Cholesky decomposition from the given triangular matrix.
  141. // SetFromU panics if t is not upper triangular. Note that t is copied into,
  142. // not stored inside, the receiver.
  143. func (c *Cholesky) SetFromU(t *TriDense) {
  144. n, kind := t.Triangle()
  145. if kind != Upper {
  146. panic("cholesky: matrix must be upper triangular")
  147. }
  148. if c.chol == nil {
  149. c.chol = NewTriDense(n, Upper, nil)
  150. } else {
  151. c.chol = NewTriDense(n, Upper, use(c.chol.mat.Data, n*n))
  152. }
  153. c.chol.Copy(t)
  154. c.updateCond(-1)
  155. }
  156. // Clone makes a copy of the input Cholesky into the receiver, overwriting the
  157. // previous value of the receiver. Clone does not place any restrictions on receiver
  158. // shape. Clone panics if the input Cholesky is not the result of a valid decomposition.
  159. func (c *Cholesky) Clone(chol *Cholesky) {
  160. if !chol.valid() {
  161. panic(badCholesky)
  162. }
  163. n := chol.Symmetric()
  164. if c.chol == nil {
  165. c.chol = NewTriDense(n, Upper, nil)
  166. } else {
  167. c.chol = NewTriDense(n, Upper, use(c.chol.mat.Data, n*n))
  168. }
  169. c.chol.Copy(chol.chol)
  170. c.cond = chol.cond
  171. }
  172. // Det returns the determinant of the matrix that has been factorized.
  173. func (c *Cholesky) Det() float64 {
  174. if !c.valid() {
  175. panic(badCholesky)
  176. }
  177. return math.Exp(c.LogDet())
  178. }
  179. // LogDet returns the log of the determinant of the matrix that has been factorized.
  180. func (c *Cholesky) LogDet() float64 {
  181. if !c.valid() {
  182. panic(badCholesky)
  183. }
  184. var det float64
  185. for i := 0; i < c.chol.mat.N; i++ {
  186. det += 2 * math.Log(c.chol.mat.Data[i*c.chol.mat.Stride+i])
  187. }
  188. return det
  189. }
  190. // SolveTo finds the matrix X that solves A * X = B where A is represented
  191. // by the Cholesky decomposition. The result is stored in-place into dst.
  192. func (c *Cholesky) SolveTo(dst *Dense, b Matrix) error {
  193. if !c.valid() {
  194. panic(badCholesky)
  195. }
  196. n := c.chol.mat.N
  197. bm, bn := b.Dims()
  198. if n != bm {
  199. panic(ErrShape)
  200. }
  201. dst.reuseAs(bm, bn)
  202. if b != dst {
  203. dst.Copy(b)
  204. }
  205. lapack64.Potrs(c.chol.mat, dst.mat)
  206. if c.cond > ConditionTolerance {
  207. return Condition(c.cond)
  208. }
  209. return nil
  210. }
  211. // SolveCholTo finds the matrix X that solves A * X = B where A and B are represented
  212. // by their Cholesky decompositions a and b. The result is stored in-place into
  213. // dst.
  214. func (a *Cholesky) SolveCholTo(dst *Dense, b *Cholesky) error {
  215. if !a.valid() || !b.valid() {
  216. panic(badCholesky)
  217. }
  218. bn := b.chol.mat.N
  219. if a.chol.mat.N != bn {
  220. panic(ErrShape)
  221. }
  222. dst.reuseAsZeroed(bn, bn)
  223. dst.Copy(b.chol.T())
  224. blas64.Trsm(blas.Left, blas.Trans, 1, a.chol.mat, dst.mat)
  225. blas64.Trsm(blas.Left, blas.NoTrans, 1, a.chol.mat, dst.mat)
  226. blas64.Trmm(blas.Right, blas.NoTrans, 1, b.chol.mat, dst.mat)
  227. if a.cond > ConditionTolerance {
  228. return Condition(a.cond)
  229. }
  230. return nil
  231. }
  232. // SolveVecTo finds the vector X that solves A * x = b where A is represented
  233. // by the Cholesky decomposition. The result is stored in-place into
  234. // dst.
  235. func (c *Cholesky) SolveVecTo(dst *VecDense, b Vector) error {
  236. if !c.valid() {
  237. panic(badCholesky)
  238. }
  239. n := c.chol.mat.N
  240. if br, bc := b.Dims(); br != n || bc != 1 {
  241. panic(ErrShape)
  242. }
  243. switch rv := b.(type) {
  244. default:
  245. dst.reuseAs(n)
  246. return c.SolveTo(dst.asDense(), b)
  247. case RawVectorer:
  248. bmat := rv.RawVector()
  249. if dst != b {
  250. dst.checkOverlap(bmat)
  251. }
  252. dst.reuseAs(n)
  253. if dst != b {
  254. dst.CopyVec(b)
  255. }
  256. lapack64.Potrs(c.chol.mat, dst.asGeneral())
  257. if c.cond > ConditionTolerance {
  258. return Condition(c.cond)
  259. }
  260. return nil
  261. }
  262. }
  263. // RawU returns the Triangular matrix used to store the Cholesky decomposition of
  264. // the original matrix A. The returned matrix should not be modified. If it is
  265. // modified, the decomposition is invalid and should not be used.
  266. func (c *Cholesky) RawU() Triangular {
  267. return c.chol
  268. }
  269. // UTo extracts the n×n upper triangular matrix U from a Cholesky
  270. // decomposition into dst and returns the result. If dst is nil a new
  271. // TriDense is allocated.
  272. // A = U^T * U.
  273. func (c *Cholesky) UTo(dst *TriDense) *TriDense {
  274. if !c.valid() {
  275. panic(badCholesky)
  276. }
  277. n := c.chol.mat.N
  278. if dst == nil {
  279. dst = NewTriDense(n, Upper, make([]float64, n*n))
  280. } else {
  281. dst.reuseAs(n, Upper)
  282. }
  283. dst.Copy(c.chol)
  284. return dst
  285. }
  286. // LTo extracts the n×n lower triangular matrix L from a Cholesky
  287. // decomposition into dst and returns the result. If dst is nil a new
  288. // TriDense is allocated.
  289. // A = L * L^T.
  290. func (c *Cholesky) LTo(dst *TriDense) *TriDense {
  291. if !c.valid() {
  292. panic(badCholesky)
  293. }
  294. n := c.chol.mat.N
  295. if dst == nil {
  296. dst = NewTriDense(n, Lower, make([]float64, n*n))
  297. } else {
  298. dst.reuseAs(n, Lower)
  299. }
  300. dst.Copy(c.chol.TTri())
  301. return dst
  302. }
  303. // ToSym reconstructs the original positive definite matrix given its
  304. // Cholesky decomposition into dst and returns the result. If dst is nil
  305. // a new SymDense is allocated.
  306. func (c *Cholesky) ToSym(dst *SymDense) *SymDense {
  307. if !c.valid() {
  308. panic(badCholesky)
  309. }
  310. n := c.chol.mat.N
  311. if dst == nil {
  312. dst = NewSymDense(n, nil)
  313. } else {
  314. dst.reuseAs(n)
  315. }
  316. // Create a TriDense representing the Cholesky factor U with dst's
  317. // backing slice.
  318. // Operations on u are reflected in s.
  319. u := &TriDense{
  320. mat: blas64.Triangular{
  321. Uplo: blas.Upper,
  322. Diag: blas.NonUnit,
  323. N: n,
  324. Data: dst.mat.Data,
  325. Stride: dst.mat.Stride,
  326. },
  327. cap: n,
  328. }
  329. u.Copy(c.chol)
  330. // Compute the product U^T*U using the algorithm from LAPACK/TESTING/LIN/dpot01.f
  331. a := u.mat.Data
  332. lda := u.mat.Stride
  333. bi := blas64.Implementation()
  334. for k := n - 1; k >= 0; k-- {
  335. a[k*lda+k] = bi.Ddot(k+1, a[k:], lda, a[k:], lda)
  336. if k > 0 {
  337. bi.Dtrmv(blas.Upper, blas.Trans, blas.NonUnit, k, a, lda, a[k:], lda)
  338. }
  339. }
  340. return dst
  341. }
  342. // InverseTo computes the inverse of the matrix represented by its Cholesky
  343. // factorization and stores the result into s. If the factorized
  344. // matrix is ill-conditioned, a Condition error will be returned.
  345. // Note that matrix inversion is numerically unstable, and should generally be
  346. // avoided where possible, for example by using the Solve routines.
  347. func (c *Cholesky) InverseTo(s *SymDense) error {
  348. if !c.valid() {
  349. panic(badCholesky)
  350. }
  351. s.reuseAs(c.chol.mat.N)
  352. // Create a TriDense representing the Cholesky factor U with the backing
  353. // slice from s.
  354. // Operations on u are reflected in s.
  355. u := &TriDense{
  356. mat: blas64.Triangular{
  357. Uplo: blas.Upper,
  358. Diag: blas.NonUnit,
  359. N: s.mat.N,
  360. Data: s.mat.Data,
  361. Stride: s.mat.Stride,
  362. },
  363. cap: s.mat.N,
  364. }
  365. u.Copy(c.chol)
  366. _, ok := lapack64.Potri(u.mat)
  367. if !ok {
  368. return Condition(math.Inf(1))
  369. }
  370. if c.cond > ConditionTolerance {
  371. return Condition(c.cond)
  372. }
  373. return nil
  374. }
  375. // Scale multiplies the original matrix A by a positive constant using
  376. // its Cholesky decomposition, storing the result in-place into the receiver.
  377. // That is, if the original Cholesky factorization is
  378. // U^T * U = A
  379. // the updated factorization is
  380. // U'^T * U' = f A = A'
  381. // Scale panics if the constant is non-positive, or if the receiver is non-zero
  382. // and is of a different size from the input.
  383. func (c *Cholesky) Scale(f float64, orig *Cholesky) {
  384. if !orig.valid() {
  385. panic(badCholesky)
  386. }
  387. if f <= 0 {
  388. panic("cholesky: scaling by a non-positive constant")
  389. }
  390. n := orig.Symmetric()
  391. if c.chol == nil {
  392. c.chol = NewTriDense(n, Upper, nil)
  393. } else if c.chol.mat.N != n {
  394. panic(ErrShape)
  395. }
  396. c.chol.ScaleTri(math.Sqrt(f), orig.chol)
  397. c.cond = orig.cond // Scaling by a positive constant does not change the condition number.
  398. }
  399. // ExtendVecSym computes the Cholesky decomposition of the original matrix A,
  400. // whose Cholesky decomposition is in a, extended by a the n×1 vector v according to
  401. // [A w]
  402. // [w' k]
  403. // where k = v[n-1] and w = v[:n-1]. The result is stored into the receiver.
  404. // In order for the updated matrix to be positive definite, it must be the case
  405. // that k > w' A^-1 w. If this condition does not hold then ExtendVecSym will
  406. // return false and the receiver will not be updated.
  407. //
  408. // ExtendVecSym will panic if v.Len() != a.Symmetric()+1 or if a does not contain
  409. // a valid decomposition.
  410. func (c *Cholesky) ExtendVecSym(a *Cholesky, v Vector) (ok bool) {
  411. n := a.Symmetric()
  412. if v.Len() != n+1 {
  413. panic(badSliceLength)
  414. }
  415. if !a.valid() {
  416. panic(badCholesky)
  417. }
  418. // The algorithm is commented here, but see also
  419. // https://math.stackexchange.com/questions/955874/cholesky-factor-when-adding-a-row-and-column-to-already-factorized-matrix
  420. // We have A and want to compute the Cholesky of
  421. // [A w]
  422. // [w' k]
  423. // We want
  424. // [U c]
  425. // [0 d]
  426. // to be the updated Cholesky, and so it must be that
  427. // [A w] = [U' 0] [U c]
  428. // [w' k] [c' d] [0 d]
  429. // Thus, we need
  430. // 1) A = U'U (true by the original decomposition being valid),
  431. // 2) U' * c = w => c = U'^-1 w
  432. // 3) c'*c + d'*d = k => d = sqrt(k-c'*c)
  433. // First, compute c = U'^-1 a
  434. // TODO(btracey): Replace this with CopyVec when issue 167 is fixed.
  435. w := NewVecDense(n, nil)
  436. for i := 0; i < n; i++ {
  437. w.SetVec(i, v.At(i, 0))
  438. }
  439. k := v.At(n, 0)
  440. var t VecDense
  441. t.SolveVec(a.chol.T(), w)
  442. dot := Dot(&t, &t)
  443. if dot >= k {
  444. return false
  445. }
  446. d := math.Sqrt(k - dot)
  447. newU := NewTriDense(n+1, Upper, nil)
  448. newU.Copy(a.chol)
  449. for i := 0; i < n; i++ {
  450. newU.SetTri(i, n, t.At(i, 0))
  451. }
  452. newU.SetTri(n, n, d)
  453. c.chol = newU
  454. c.updateCond(-1)
  455. return true
  456. }
  457. // SymRankOne performs a rank-1 update of the original matrix A and refactorizes
  458. // its Cholesky factorization, storing the result into the receiver. That is, if
  459. // in the original Cholesky factorization
  460. // U^T * U = A,
  461. // in the updated factorization
  462. // U'^T * U' = A + alpha * x * x^T = A'.
  463. //
  464. // Note that when alpha is negative, the updating problem may be ill-conditioned
  465. // and the results may be inaccurate, or the updated matrix A' may not be
  466. // positive definite and not have a Cholesky factorization. SymRankOne returns
  467. // whether the updated matrix A' is positive definite.
  468. //
  469. // SymRankOne updates a Cholesky factorization in O(n²) time. The Cholesky
  470. // factorization computation from scratch is O(n³).
  471. func (c *Cholesky) SymRankOne(orig *Cholesky, alpha float64, x Vector) (ok bool) {
  472. if !orig.valid() {
  473. panic(badCholesky)
  474. }
  475. n := orig.Symmetric()
  476. if r, c := x.Dims(); r != n || c != 1 {
  477. panic(ErrShape)
  478. }
  479. if orig != c {
  480. if c.chol == nil {
  481. c.chol = NewTriDense(n, Upper, nil)
  482. } else if c.chol.mat.N != n {
  483. panic(ErrShape)
  484. }
  485. c.chol.Copy(orig.chol)
  486. }
  487. if alpha == 0 {
  488. return true
  489. }
  490. // Algorithms for updating and downdating the Cholesky factorization are
  491. // described, for example, in
  492. // - J. J. Dongarra, J. R. Bunch, C. B. Moler, G. W. Stewart: LINPACK
  493. // Users' Guide. SIAM (1979), pages 10.10--10.14
  494. // or
  495. // - P. E. Gill, G. H. Golub, W. Murray, and M. A. Saunders: Methods for
  496. // modifying matrix factorizations. Mathematics of Computation 28(126)
  497. // (1974), Method C3 on page 521
  498. //
  499. // The implementation is based on LINPACK code
  500. // http://www.netlib.org/linpack/dchud.f
  501. // http://www.netlib.org/linpack/dchdd.f
  502. // and
  503. // https://icl.cs.utk.edu/lapack-forum/viewtopic.php?f=2&t=2646
  504. //
  505. // According to http://icl.cs.utk.edu/lapack-forum/archives/lapack/msg00301.html
  506. // LINPACK is released under BSD license.
  507. //
  508. // See also:
  509. // - M. A. Saunders: Large-scale Linear Programming Using the Cholesky
  510. // Factorization. Technical Report Stanford University (1972)
  511. // http://i.stanford.edu/pub/cstr/reports/cs/tr/72/252/CS-TR-72-252.pdf
  512. // - Matthias Seeger: Low rank updates for the Cholesky decomposition.
  513. // EPFL Technical Report 161468 (2004)
  514. // http://infoscience.epfl.ch/record/161468
  515. work := getFloats(n, false)
  516. defer putFloats(work)
  517. var xmat blas64.Vector
  518. if rv, ok := x.(RawVectorer); ok {
  519. xmat = rv.RawVector()
  520. } else {
  521. var tmp *VecDense
  522. tmp.CopyVec(x)
  523. xmat = tmp.RawVector()
  524. }
  525. blas64.Copy(xmat, blas64.Vector{N: n, Data: work, Inc: 1})
  526. if alpha > 0 {
  527. // Compute rank-1 update.
  528. if alpha != 1 {
  529. blas64.Scal(math.Sqrt(alpha), blas64.Vector{N: n, Data: work, Inc: 1})
  530. }
  531. umat := c.chol.mat
  532. stride := umat.Stride
  533. for i := 0; i < n; i++ {
  534. // Compute parameters of the Givens matrix that zeroes
  535. // the i-th element of x.
  536. c, s, r, _ := blas64.Rotg(umat.Data[i*stride+i], work[i])
  537. if r < 0 {
  538. // Multiply by -1 to have positive diagonal
  539. // elemnts.
  540. r *= -1
  541. c *= -1
  542. s *= -1
  543. }
  544. umat.Data[i*stride+i] = r
  545. if i < n-1 {
  546. // Multiply the extended factorization matrix by
  547. // the Givens matrix from the left. Only
  548. // the i-th row and x are modified.
  549. blas64.Rot(
  550. blas64.Vector{N: n - i - 1, Data: umat.Data[i*stride+i+1 : i*stride+n], Inc: 1},
  551. blas64.Vector{N: n - i - 1, Data: work[i+1 : n], Inc: 1},
  552. c, s)
  553. }
  554. }
  555. c.updateCond(-1)
  556. return true
  557. }
  558. // Compute rank-1 downdate.
  559. alpha = math.Sqrt(-alpha)
  560. if alpha != 1 {
  561. blas64.Scal(alpha, blas64.Vector{N: n, Data: work, Inc: 1})
  562. }
  563. // Solve U^T * p = x storing the result into work.
  564. ok = lapack64.Trtrs(blas.Trans, c.chol.RawTriangular(), blas64.General{
  565. Rows: n,
  566. Cols: 1,
  567. Stride: 1,
  568. Data: work,
  569. })
  570. if !ok {
  571. // The original matrix is singular. Should not happen, because
  572. // the factorization is valid.
  573. panic(badCholesky)
  574. }
  575. norm := blas64.Nrm2(blas64.Vector{N: n, Data: work, Inc: 1})
  576. if norm >= 1 {
  577. // The updated matrix is not positive definite.
  578. return false
  579. }
  580. norm = math.Sqrt((1 + norm) * (1 - norm))
  581. cos := getFloats(n, false)
  582. defer putFloats(cos)
  583. sin := getFloats(n, false)
  584. defer putFloats(sin)
  585. for i := n - 1; i >= 0; i-- {
  586. // Compute parameters of Givens matrices that zero elements of p
  587. // backwards.
  588. cos[i], sin[i], norm, _ = blas64.Rotg(norm, work[i])
  589. if norm < 0 {
  590. norm *= -1
  591. cos[i] *= -1
  592. sin[i] *= -1
  593. }
  594. }
  595. umat := c.chol.mat
  596. stride := umat.Stride
  597. for i := n - 1; i >= 0; i-- {
  598. work[i] = 0
  599. // Apply Givens matrices to U.
  600. // TODO(vladimir-ch): Use workspace to avoid modifying the
  601. // receiver in case an invalid factorization is created.
  602. blas64.Rot(
  603. blas64.Vector{N: n - i, Data: work[i:n], Inc: 1},
  604. blas64.Vector{N: n - i, Data: umat.Data[i*stride+i : i*stride+n], Inc: 1},
  605. cos[i], sin[i])
  606. if umat.Data[i*stride+i] == 0 {
  607. // The matrix is singular (may rarely happen due to
  608. // floating-point effects?).
  609. ok = false
  610. } else if umat.Data[i*stride+i] < 0 {
  611. // Diagonal elements should be positive. If it happens
  612. // that on the i-th row the diagonal is negative,
  613. // multiply U from the left by an identity matrix that
  614. // has -1 on the i-th row.
  615. blas64.Scal(-1, blas64.Vector{N: n - i, Data: umat.Data[i*stride+i : i*stride+n], Inc: 1})
  616. }
  617. }
  618. if ok {
  619. c.updateCond(-1)
  620. } else {
  621. c.Reset()
  622. }
  623. return ok
  624. }
  625. func (c *Cholesky) valid() bool {
  626. return c.chol != nil && !c.chol.IsZero()
  627. }