github.com/jingcheng-WU/gonum@v0.9.1-0.20210323123734-f1a2a11a8f7b/mat/cholesky.go (about) 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 5 package mat 6 7 import ( 8 "math" 9 10 "github.com/jingcheng-WU/gonum/blas" 11 "github.com/jingcheng-WU/gonum/blas/blas64" 12 "github.com/jingcheng-WU/gonum/lapack/lapack64" 13 ) 14 15 const ( 16 badTriangle = "mat: invalid triangle" 17 badCholesky = "mat: invalid Cholesky factorization" 18 ) 19 20 var ( 21 _ Matrix = (*Cholesky)(nil) 22 _ Symmetric = (*Cholesky)(nil) 23 24 _ Matrix = (*BandCholesky)(nil) 25 _ Symmetric = (*BandCholesky)(nil) 26 _ Banded = (*BandCholesky)(nil) 27 _ SymBanded = (*BandCholesky)(nil) 28 ) 29 30 // Cholesky is a symmetric positive definite matrix represented by its 31 // Cholesky decomposition. 32 // 33 // The decomposition can be constructed using the Factorize method. The 34 // factorization itself can be extracted using the UTo or LTo methods, and the 35 // original symmetric matrix can be recovered with ToSym. 36 // 37 // Note that this matrix representation is useful for certain operations, in 38 // particular finding solutions to linear equations. It is very inefficient 39 // at other operations, in particular At is slow. 40 // 41 // Cholesky methods may only be called on a value that has been successfully 42 // initialized by a call to Factorize that has returned true. Calls to methods 43 // of an unsuccessful Cholesky factorization will panic. 44 type Cholesky struct { 45 // The chol pointer must never be retained as a pointer outside the Cholesky 46 // struct, either by returning chol outside the struct or by setting it to 47 // a pointer coming from outside. The same prohibition applies to the data 48 // slice within chol. 49 chol *TriDense 50 cond float64 51 } 52 53 // updateCond updates the condition number of the Cholesky decomposition. If 54 // norm > 0, then that norm is used as the norm of the original matrix A, otherwise 55 // the norm is estimated from the decomposition. 56 func (c *Cholesky) updateCond(norm float64) { 57 n := c.chol.mat.N 58 work := getFloats(3*n, false) 59 defer putFloats(work) 60 if norm < 0 { 61 // This is an approximation. By the definition of a norm, 62 // |AB| <= |A| |B|. 63 // Since A = Uᵀ*U, we get for the condition number κ that 64 // κ(A) := |A| |A^-1| = |Uᵀ*U| |A^-1| <= |Uᵀ| |U| |A^-1|, 65 // so this will overestimate the condition number somewhat. 66 // The norm of the original factorized matrix cannot be stored 67 // because of update possibilities. 68 unorm := lapack64.Lantr(CondNorm, c.chol.mat, work) 69 lnorm := lapack64.Lantr(CondNormTrans, c.chol.mat, work) 70 norm = unorm * lnorm 71 } 72 sym := c.chol.asSymBlas() 73 iwork := getInts(n, false) 74 v := lapack64.Pocon(sym, norm, work, iwork) 75 putInts(iwork) 76 c.cond = 1 / v 77 } 78 79 // Dims returns the dimensions of the matrix. 80 func (ch *Cholesky) Dims() (r, c int) { 81 if !ch.valid() { 82 panic(badCholesky) 83 } 84 r, c = ch.chol.Dims() 85 return r, c 86 } 87 88 // At returns the element at row i, column j. 89 func (c *Cholesky) At(i, j int) float64 { 90 if !c.valid() { 91 panic(badCholesky) 92 } 93 n := c.Symmetric() 94 if uint(i) >= uint(n) { 95 panic(ErrRowAccess) 96 } 97 if uint(j) >= uint(n) { 98 panic(ErrColAccess) 99 } 100 101 var val float64 102 for k := 0; k <= min(i, j); k++ { 103 val += c.chol.at(k, i) * c.chol.at(k, j) 104 } 105 return val 106 } 107 108 // T returns the receiver, the transpose of a symmetric matrix. 109 func (c *Cholesky) T() Matrix { 110 return c 111 } 112 113 // Symmetric implements the Symmetric interface and returns the number of rows 114 // in the matrix (this is also the number of columns). 115 func (c *Cholesky) Symmetric() int { 116 r, _ := c.chol.Dims() 117 return r 118 } 119 120 // Cond returns the condition number of the factorized matrix. 121 func (c *Cholesky) Cond() float64 { 122 if !c.valid() { 123 panic(badCholesky) 124 } 125 return c.cond 126 } 127 128 // Factorize calculates the Cholesky decomposition of the matrix A and returns 129 // whether the matrix is positive definite. If Factorize returns false, the 130 // factorization must not be used. 131 func (c *Cholesky) Factorize(a Symmetric) (ok bool) { 132 n := a.Symmetric() 133 if c.chol == nil { 134 c.chol = NewTriDense(n, Upper, nil) 135 } else { 136 c.chol.Reset() 137 c.chol.reuseAsNonZeroed(n, Upper) 138 } 139 copySymIntoTriangle(c.chol, a) 140 141 sym := c.chol.asSymBlas() 142 work := getFloats(c.chol.mat.N, false) 143 norm := lapack64.Lansy(CondNorm, sym, work) 144 putFloats(work) 145 _, ok = lapack64.Potrf(sym) 146 if ok { 147 c.updateCond(norm) 148 } else { 149 c.Reset() 150 } 151 return ok 152 } 153 154 // Reset resets the factorization so that it can be reused as the receiver of a 155 // dimensionally restricted operation. 156 func (c *Cholesky) Reset() { 157 if c.chol != nil { 158 c.chol.Reset() 159 } 160 c.cond = math.Inf(1) 161 } 162 163 // IsEmpty returns whether the receiver is empty. Empty matrices can be the 164 // receiver for size-restricted operations. The receiver can be emptied using 165 // Reset. 166 func (c *Cholesky) IsEmpty() bool { 167 return c.chol == nil || c.chol.IsEmpty() 168 } 169 170 // SetFromU sets the Cholesky decomposition from the given triangular matrix. 171 // SetFromU panics if t is not upper triangular. If the receiver is empty it 172 // is resized to be n×n, the size of t. If dst is non-empty, SetFromU panics 173 // if c is not of size n×n. Note that t is copied into, not stored inside, the 174 // receiver. 175 func (c *Cholesky) SetFromU(t Triangular) { 176 n, kind := t.Triangle() 177 if kind != Upper { 178 panic("cholesky: matrix must be upper triangular") 179 } 180 if c.chol == nil { 181 c.chol = NewTriDense(n, Upper, nil) 182 } else { 183 c.chol.reuseAsNonZeroed(n, Upper) 184 } 185 c.chol.Copy(t) 186 c.updateCond(-1) 187 } 188 189 // Clone makes a copy of the input Cholesky into the receiver, overwriting the 190 // previous value of the receiver. Clone does not place any restrictions on receiver 191 // shape. Clone panics if the input Cholesky is not the result of a valid decomposition. 192 func (c *Cholesky) Clone(chol *Cholesky) { 193 if !chol.valid() { 194 panic(badCholesky) 195 } 196 n := chol.Symmetric() 197 if c.chol == nil { 198 c.chol = NewTriDense(n, Upper, nil) 199 } else { 200 c.chol = NewTriDense(n, Upper, use(c.chol.mat.Data, n*n)) 201 } 202 c.chol.Copy(chol.chol) 203 c.cond = chol.cond 204 } 205 206 // Det returns the determinant of the matrix that has been factorized. 207 func (c *Cholesky) Det() float64 { 208 if !c.valid() { 209 panic(badCholesky) 210 } 211 return math.Exp(c.LogDet()) 212 } 213 214 // LogDet returns the log of the determinant of the matrix that has been factorized. 215 func (c *Cholesky) LogDet() float64 { 216 if !c.valid() { 217 panic(badCholesky) 218 } 219 var det float64 220 for i := 0; i < c.chol.mat.N; i++ { 221 det += 2 * math.Log(c.chol.mat.Data[i*c.chol.mat.Stride+i]) 222 } 223 return det 224 } 225 226 // SolveTo finds the matrix X that solves A * X = B where A is represented 227 // by the Cholesky decomposition. The result is stored in-place into dst. 228 // If the Cholesky decomposition is singular or near-singular a Condition error 229 // is returned. See the documentation for Condition for more information. 230 func (c *Cholesky) SolveTo(dst *Dense, b Matrix) error { 231 if !c.valid() { 232 panic(badCholesky) 233 } 234 n := c.chol.mat.N 235 bm, bn := b.Dims() 236 if n != bm { 237 panic(ErrShape) 238 } 239 240 dst.reuseAsNonZeroed(bm, bn) 241 if b != dst { 242 dst.Copy(b) 243 } 244 lapack64.Potrs(c.chol.mat, dst.mat) 245 if c.cond > ConditionTolerance { 246 return Condition(c.cond) 247 } 248 return nil 249 } 250 251 // SolveCholTo finds the matrix X that solves A * X = B where A and B are represented 252 // by their Cholesky decompositions a and b. The result is stored in-place into 253 // dst. 254 // If the Cholesky decomposition is singular or near-singular a Condition error 255 // is returned. See the documentation for Condition for more information. 256 func (a *Cholesky) SolveCholTo(dst *Dense, b *Cholesky) error { 257 if !a.valid() || !b.valid() { 258 panic(badCholesky) 259 } 260 bn := b.chol.mat.N 261 if a.chol.mat.N != bn { 262 panic(ErrShape) 263 } 264 265 dst.reuseAsZeroed(bn, bn) 266 dst.Copy(b.chol.T()) 267 blas64.Trsm(blas.Left, blas.Trans, 1, a.chol.mat, dst.mat) 268 blas64.Trsm(blas.Left, blas.NoTrans, 1, a.chol.mat, dst.mat) 269 blas64.Trmm(blas.Right, blas.NoTrans, 1, b.chol.mat, dst.mat) 270 if a.cond > ConditionTolerance { 271 return Condition(a.cond) 272 } 273 return nil 274 } 275 276 // SolveVecTo finds the vector x that solves A * x = b where A is represented 277 // by the Cholesky decomposition. The result is stored in-place into 278 // dst. 279 // If the Cholesky decomposition is singular or near-singular a Condition error 280 // is returned. See the documentation for Condition for more information. 281 func (c *Cholesky) SolveVecTo(dst *VecDense, b Vector) error { 282 if !c.valid() { 283 panic(badCholesky) 284 } 285 n := c.chol.mat.N 286 if br, bc := b.Dims(); br != n || bc != 1 { 287 panic(ErrShape) 288 } 289 switch rv := b.(type) { 290 default: 291 dst.reuseAsNonZeroed(n) 292 return c.SolveTo(dst.asDense(), b) 293 case RawVectorer: 294 bmat := rv.RawVector() 295 if dst != b { 296 dst.checkOverlap(bmat) 297 } 298 dst.reuseAsNonZeroed(n) 299 if dst != b { 300 dst.CopyVec(b) 301 } 302 lapack64.Potrs(c.chol.mat, dst.asGeneral()) 303 if c.cond > ConditionTolerance { 304 return Condition(c.cond) 305 } 306 return nil 307 } 308 } 309 310 // RawU returns the Triangular matrix used to store the Cholesky decomposition of 311 // the original matrix A. The returned matrix should not be modified. If it is 312 // modified, the decomposition is invalid and should not be used. 313 func (c *Cholesky) RawU() Triangular { 314 return c.chol 315 } 316 317 // UTo stores into dst the n×n upper triangular matrix U from a Cholesky 318 // decomposition 319 // A = Uᵀ * U. 320 // If dst is empty, it is resized to be an n×n upper triangular matrix. When dst 321 // is non-empty, UTo panics if dst is not n×n or not Upper. UTo will also panic 322 // if the receiver does not contain a successful factorization. 323 func (c *Cholesky) UTo(dst *TriDense) { 324 if !c.valid() { 325 panic(badCholesky) 326 } 327 n := c.chol.mat.N 328 if dst.IsEmpty() { 329 dst.ReuseAsTri(n, Upper) 330 } else { 331 n2, kind := dst.Triangle() 332 if n != n2 { 333 panic(ErrShape) 334 } 335 if kind != Upper { 336 panic(ErrTriangle) 337 } 338 } 339 dst.Copy(c.chol) 340 } 341 342 // LTo stores into dst the n×n lower triangular matrix L from a Cholesky 343 // decomposition 344 // A = L * Lᵀ. 345 // If dst is empty, it is resized to be an n×n lower triangular matrix. When dst 346 // is non-empty, LTo panics if dst is not n×n or not Lower. LTo will also panic 347 // if the receiver does not contain a successful factorization. 348 func (c *Cholesky) LTo(dst *TriDense) { 349 if !c.valid() { 350 panic(badCholesky) 351 } 352 n := c.chol.mat.N 353 if dst.IsEmpty() { 354 dst.ReuseAsTri(n, Lower) 355 } else { 356 n2, kind := dst.Triangle() 357 if n != n2 { 358 panic(ErrShape) 359 } 360 if kind != Lower { 361 panic(ErrTriangle) 362 } 363 } 364 dst.Copy(c.chol.TTri()) 365 } 366 367 // ToSym reconstructs the original positive definite matrix from its 368 // Cholesky decomposition, storing the result into dst. If dst is 369 // empty it is resized to be n×n. If dst is non-empty, ToSym panics 370 // if dst is not of size n×n. ToSym will also panic if the receiver 371 // does not contain a successful factorization. 372 func (c *Cholesky) ToSym(dst *SymDense) { 373 if !c.valid() { 374 panic(badCholesky) 375 } 376 n := c.chol.mat.N 377 if dst.IsEmpty() { 378 dst.ReuseAsSym(n) 379 } else { 380 n2 := dst.Symmetric() 381 if n != n2 { 382 panic(ErrShape) 383 } 384 } 385 // Create a TriDense representing the Cholesky factor U with dst's 386 // backing slice. 387 // Operations on u are reflected in s. 388 u := &TriDense{ 389 mat: blas64.Triangular{ 390 Uplo: blas.Upper, 391 Diag: blas.NonUnit, 392 N: n, 393 Data: dst.mat.Data, 394 Stride: dst.mat.Stride, 395 }, 396 cap: n, 397 } 398 u.Copy(c.chol) 399 // Compute the product Uᵀ*U using the algorithm from LAPACK/TESTING/LIN/dpot01.f 400 a := u.mat.Data 401 lda := u.mat.Stride 402 bi := blas64.Implementation() 403 for k := n - 1; k >= 0; k-- { 404 a[k*lda+k] = bi.Ddot(k+1, a[k:], lda, a[k:], lda) 405 if k > 0 { 406 bi.Dtrmv(blas.Upper, blas.Trans, blas.NonUnit, k, a, lda, a[k:], lda) 407 } 408 } 409 } 410 411 // InverseTo computes the inverse of the matrix represented by its Cholesky 412 // factorization and stores the result into s. If the factorized 413 // matrix is ill-conditioned, a Condition error will be returned. 414 // Note that matrix inversion is numerically unstable, and should generally be 415 // avoided where possible, for example by using the Solve routines. 416 func (c *Cholesky) InverseTo(dst *SymDense) error { 417 if !c.valid() { 418 panic(badCholesky) 419 } 420 dst.reuseAsNonZeroed(c.chol.mat.N) 421 // Create a TriDense representing the Cholesky factor U with the backing 422 // slice from dst. 423 // Operations on u are reflected in dst. 424 u := &TriDense{ 425 mat: blas64.Triangular{ 426 Uplo: blas.Upper, 427 Diag: blas.NonUnit, 428 N: dst.mat.N, 429 Data: dst.mat.Data, 430 Stride: dst.mat.Stride, 431 }, 432 cap: dst.mat.N, 433 } 434 u.Copy(c.chol) 435 436 _, ok := lapack64.Potri(u.mat) 437 if !ok { 438 return Condition(math.Inf(1)) 439 } 440 if c.cond > ConditionTolerance { 441 return Condition(c.cond) 442 } 443 return nil 444 } 445 446 // Scale multiplies the original matrix A by a positive constant using 447 // its Cholesky decomposition, storing the result in-place into the receiver. 448 // That is, if the original Cholesky factorization is 449 // Uᵀ * U = A 450 // the updated factorization is 451 // U'ᵀ * U' = f A = A' 452 // Scale panics if the constant is non-positive, or if the receiver is non-empty 453 // and is of a different size from the input. 454 func (c *Cholesky) Scale(f float64, orig *Cholesky) { 455 if !orig.valid() { 456 panic(badCholesky) 457 } 458 if f <= 0 { 459 panic("cholesky: scaling by a non-positive constant") 460 } 461 n := orig.Symmetric() 462 if c.chol == nil { 463 c.chol = NewTriDense(n, Upper, nil) 464 } else if c.chol.mat.N != n { 465 panic(ErrShape) 466 } 467 c.chol.ScaleTri(math.Sqrt(f), orig.chol) 468 c.cond = orig.cond // Scaling by a positive constant does not change the condition number. 469 } 470 471 // ExtendVecSym computes the Cholesky decomposition of the original matrix A, 472 // whose Cholesky decomposition is in a, extended by a the n×1 vector v according to 473 // [A w] 474 // [w' k] 475 // where k = v[n-1] and w = v[:n-1]. The result is stored into the receiver. 476 // In order for the updated matrix to be positive definite, it must be the case 477 // that k > w' A^-1 w. If this condition does not hold then ExtendVecSym will 478 // return false and the receiver will not be updated. 479 // 480 // ExtendVecSym will panic if v.Len() != a.Symmetric()+1 or if a does not contain 481 // a valid decomposition. 482 func (c *Cholesky) ExtendVecSym(a *Cholesky, v Vector) (ok bool) { 483 n := a.Symmetric() 484 485 if v.Len() != n+1 { 486 panic(badSliceLength) 487 } 488 if !a.valid() { 489 panic(badCholesky) 490 } 491 492 // The algorithm is commented here, but see also 493 // https://math.stackexchange.com/questions/955874/cholesky-factor-when-adding-a-row-and-column-to-already-factorized-matrix 494 // We have A and want to compute the Cholesky of 495 // [A w] 496 // [w' k] 497 // We want 498 // [U c] 499 // [0 d] 500 // to be the updated Cholesky, and so it must be that 501 // [A w] = [U' 0] [U c] 502 // [w' k] [c' d] [0 d] 503 // Thus, we need 504 // 1) A = U'U (true by the original decomposition being valid), 505 // 2) U' * c = w => c = U'^-1 w 506 // 3) c'*c + d'*d = k => d = sqrt(k-c'*c) 507 508 // First, compute c = U'^-1 a 509 w := NewVecDense(n, nil) 510 w.CopyVec(v) 511 k := v.At(n, 0) 512 513 var t VecDense 514 _ = t.SolveVec(a.chol.T(), w) 515 516 dot := Dot(&t, &t) 517 if dot >= k { 518 return false 519 } 520 d := math.Sqrt(k - dot) 521 522 newU := NewTriDense(n+1, Upper, nil) 523 newU.Copy(a.chol) 524 for i := 0; i < n; i++ { 525 newU.SetTri(i, n, t.At(i, 0)) 526 } 527 newU.SetTri(n, n, d) 528 c.chol = newU 529 c.updateCond(-1) 530 return true 531 } 532 533 // SymRankOne performs a rank-1 update of the original matrix A and refactorizes 534 // its Cholesky factorization, storing the result into the receiver. That is, if 535 // in the original Cholesky factorization 536 // Uᵀ * U = A, 537 // in the updated factorization 538 // U'ᵀ * U' = A + alpha * x * xᵀ = A'. 539 // 540 // Note that when alpha is negative, the updating problem may be ill-conditioned 541 // and the results may be inaccurate, or the updated matrix A' may not be 542 // positive definite and not have a Cholesky factorization. SymRankOne returns 543 // whether the updated matrix A' is positive definite. If the update fails 544 // the receiver is left unchanged. 545 // 546 // SymRankOne updates a Cholesky factorization in O(n²) time. The Cholesky 547 // factorization computation from scratch is O(n³). 548 func (c *Cholesky) SymRankOne(orig *Cholesky, alpha float64, x Vector) (ok bool) { 549 if !orig.valid() { 550 panic(badCholesky) 551 } 552 n := orig.Symmetric() 553 if r, c := x.Dims(); r != n || c != 1 { 554 panic(ErrShape) 555 } 556 if orig != c { 557 if c.chol == nil { 558 c.chol = NewTriDense(n, Upper, nil) 559 } else if c.chol.mat.N != n { 560 panic(ErrShape) 561 } 562 c.chol.Copy(orig.chol) 563 } 564 565 if alpha == 0 { 566 return true 567 } 568 569 // Algorithms for updating and downdating the Cholesky factorization are 570 // described, for example, in 571 // - J. J. Dongarra, J. R. Bunch, C. B. Moler, G. W. Stewart: LINPACK 572 // Users' Guide. SIAM (1979), pages 10.10--10.14 573 // or 574 // - P. E. Gill, G. H. Golub, W. Murray, and M. A. Saunders: Methods for 575 // modifying matrix factorizations. Mathematics of Computation 28(126) 576 // (1974), Method C3 on page 521 577 // 578 // The implementation is based on LINPACK code 579 // http://www.netlib.org/linpack/dchud.f 580 // http://www.netlib.org/linpack/dchdd.f 581 // and 582 // https://icl.cs.utk.edu/lapack-forum/viewtopic.php?f=2&t=2646 583 // 584 // According to http://icl.cs.utk.edu/lapack-forum/archives/lapack/msg00301.html 585 // LINPACK is released under BSD license. 586 // 587 // See also: 588 // - M. A. Saunders: Large-scale Linear Programming Using the Cholesky 589 // Factorization. Technical Report Stanford University (1972) 590 // http://i.stanford.edu/pub/cstr/reports/cs/tr/72/252/CS-TR-72-252.pdf 591 // - Matthias Seeger: Low rank updates for the Cholesky decomposition. 592 // EPFL Technical Report 161468 (2004) 593 // http://infoscience.epfl.ch/record/161468 594 595 work := getFloats(n, false) 596 defer putFloats(work) 597 var xmat blas64.Vector 598 if rv, ok := x.(RawVectorer); ok { 599 xmat = rv.RawVector() 600 } else { 601 var tmp *VecDense 602 tmp.CopyVec(x) 603 xmat = tmp.RawVector() 604 } 605 blas64.Copy(xmat, blas64.Vector{N: n, Data: work, Inc: 1}) 606 607 if alpha > 0 { 608 // Compute rank-1 update. 609 if alpha != 1 { 610 blas64.Scal(math.Sqrt(alpha), blas64.Vector{N: n, Data: work, Inc: 1}) 611 } 612 umat := c.chol.mat 613 stride := umat.Stride 614 for i := 0; i < n; i++ { 615 // Compute parameters of the Givens matrix that zeroes 616 // the i-th element of x. 617 c, s, r, _ := blas64.Rotg(umat.Data[i*stride+i], work[i]) 618 if r < 0 { 619 // Multiply by -1 to have positive diagonal 620 // elemnts. 621 r *= -1 622 c *= -1 623 s *= -1 624 } 625 umat.Data[i*stride+i] = r 626 if i < n-1 { 627 // Multiply the extended factorization matrix by 628 // the Givens matrix from the left. Only 629 // the i-th row and x are modified. 630 blas64.Rot( 631 blas64.Vector{N: n - i - 1, Data: umat.Data[i*stride+i+1 : i*stride+n], Inc: 1}, 632 blas64.Vector{N: n - i - 1, Data: work[i+1 : n], Inc: 1}, 633 c, s) 634 } 635 } 636 c.updateCond(-1) 637 return true 638 } 639 640 // Compute rank-1 downdate. 641 alpha = math.Sqrt(-alpha) 642 if alpha != 1 { 643 blas64.Scal(alpha, blas64.Vector{N: n, Data: work, Inc: 1}) 644 } 645 // Solve Uᵀ * p = x storing the result into work. 646 ok = lapack64.Trtrs(blas.Trans, c.chol.RawTriangular(), blas64.General{ 647 Rows: n, 648 Cols: 1, 649 Stride: 1, 650 Data: work, 651 }) 652 if !ok { 653 // The original matrix is singular. Should not happen, because 654 // the factorization is valid. 655 panic(badCholesky) 656 } 657 norm := blas64.Nrm2(blas64.Vector{N: n, Data: work, Inc: 1}) 658 if norm >= 1 { 659 // The updated matrix is not positive definite. 660 return false 661 } 662 norm = math.Sqrt((1 + norm) * (1 - norm)) 663 cos := getFloats(n, false) 664 defer putFloats(cos) 665 sin := getFloats(n, false) 666 defer putFloats(sin) 667 for i := n - 1; i >= 0; i-- { 668 // Compute parameters of Givens matrices that zero elements of p 669 // backwards. 670 cos[i], sin[i], norm, _ = blas64.Rotg(norm, work[i]) 671 if norm < 0 { 672 norm *= -1 673 cos[i] *= -1 674 sin[i] *= -1 675 } 676 } 677 workMat := getWorkspaceTri(c.chol.mat.N, c.chol.triKind(), false) 678 defer putWorkspaceTri(workMat) 679 workMat.Copy(c.chol) 680 umat := workMat.mat 681 stride := workMat.mat.Stride 682 for i := n - 1; i >= 0; i-- { 683 work[i] = 0 684 // Apply Givens matrices to U. 685 blas64.Rot( 686 blas64.Vector{N: n - i, Data: work[i:n], Inc: 1}, 687 blas64.Vector{N: n - i, Data: umat.Data[i*stride+i : i*stride+n], Inc: 1}, 688 cos[i], sin[i]) 689 if umat.Data[i*stride+i] == 0 { 690 // The matrix is singular (may rarely happen due to 691 // floating-point effects?). 692 ok = false 693 } else if umat.Data[i*stride+i] < 0 { 694 // Diagonal elements should be positive. If it happens 695 // that on the i-th row the diagonal is negative, 696 // multiply U from the left by an identity matrix that 697 // has -1 on the i-th row. 698 blas64.Scal(-1, blas64.Vector{N: n - i, Data: umat.Data[i*stride+i : i*stride+n], Inc: 1}) 699 } 700 } 701 if ok { 702 c.chol.Copy(workMat) 703 c.updateCond(-1) 704 } 705 return ok 706 } 707 708 func (c *Cholesky) valid() bool { 709 return c.chol != nil && !c.chol.IsEmpty() 710 } 711 712 // BandCholesky is a symmetric positive-definite band matrix represented by its 713 // Cholesky decomposition. 714 // 715 // Note that this matrix representation is useful for certain operations, in 716 // particular finding solutions to linear equations. It is very inefficient at 717 // other operations, in particular At is slow. 718 // 719 // BandCholesky methods may only be called on a value that has been successfully 720 // initialized by a call to Factorize that has returned true. Calls to methods 721 // of an unsuccessful Cholesky factorization will panic. 722 type BandCholesky struct { 723 // The chol pointer must never be retained as a pointer outside the Cholesky 724 // struct, either by returning chol outside the struct or by setting it to 725 // a pointer coming from outside. The same prohibition applies to the data 726 // slice within chol. 727 chol *TriBandDense 728 cond float64 729 } 730 731 // Factorize calculates the Cholesky decomposition of the matrix A and returns 732 // whether the matrix is positive definite. If Factorize returns false, the 733 // factorization must not be used. 734 func (ch *BandCholesky) Factorize(a SymBanded) (ok bool) { 735 n, k := a.SymBand() 736 if ch.chol == nil { 737 ch.chol = NewTriBandDense(n, k, Upper, nil) 738 } else { 739 ch.chol.Reset() 740 ch.chol.ReuseAsTriBand(n, k, Upper) 741 } 742 copySymBandIntoTriBand(ch.chol, a) 743 cSym := blas64.SymmetricBand{ 744 Uplo: blas.Upper, 745 N: n, 746 K: k, 747 Data: ch.chol.RawTriBand().Data, 748 Stride: ch.chol.RawTriBand().Stride, 749 } 750 _, ok = lapack64.Pbtrf(cSym) 751 if !ok { 752 ch.Reset() 753 return false 754 } 755 work := getFloats(3*n, false) 756 iwork := getInts(n, false) 757 aNorm := lapack64.Lansb(CondNorm, cSym, work) 758 ch.cond = 1 / lapack64.Pbcon(cSym, aNorm, work, iwork) 759 putInts(iwork) 760 putFloats(work) 761 return true 762 } 763 764 // SolveTo finds the matrix X that solves A * X = B where A is represented by 765 // the Cholesky decomposition. The result is stored in-place into dst. 766 // If the Cholesky decomposition is singular or near-singular a Condition error 767 // is returned. See the documentation for Condition for more information. 768 func (ch *BandCholesky) SolveTo(dst *Dense, b Matrix) error { 769 if !ch.valid() { 770 panic(badCholesky) 771 } 772 br, bc := b.Dims() 773 if br != ch.chol.mat.N { 774 panic(ErrShape) 775 } 776 dst.reuseAsNonZeroed(br, bc) 777 if b != dst { 778 dst.Copy(b) 779 } 780 lapack64.Pbtrs(ch.chol.mat, dst.mat) 781 if ch.cond > ConditionTolerance { 782 return Condition(ch.cond) 783 } 784 return nil 785 } 786 787 // SolveVecTo finds the vector x that solves A * x = b where A is represented by 788 // the Cholesky decomposition. The result is stored in-place into dst. 789 // If the Cholesky decomposition is singular or near-singular a Condition error 790 // is returned. See the documentation for Condition for more information. 791 func (ch *BandCholesky) SolveVecTo(dst *VecDense, b Vector) error { 792 if !ch.valid() { 793 panic(badCholesky) 794 } 795 n := ch.chol.mat.N 796 if br, bc := b.Dims(); br != n || bc != 1 { 797 panic(ErrShape) 798 } 799 if b, ok := b.(RawVectorer); ok && dst != b { 800 dst.checkOverlap(b.RawVector()) 801 } 802 dst.reuseAsNonZeroed(n) 803 if dst != b { 804 dst.CopyVec(b) 805 } 806 lapack64.Pbtrs(ch.chol.mat, dst.asGeneral()) 807 if ch.cond > ConditionTolerance { 808 return Condition(ch.cond) 809 } 810 return nil 811 } 812 813 // Cond returns the condition number of the factorized matrix. 814 func (ch *BandCholesky) Cond() float64 { 815 if !ch.valid() { 816 panic(badCholesky) 817 } 818 return ch.cond 819 } 820 821 // Reset resets the factorization so that it can be reused as the receiver of 822 // a dimensionally restricted operation. 823 func (ch *BandCholesky) Reset() { 824 if ch.chol != nil { 825 ch.chol.Reset() 826 } 827 ch.cond = math.Inf(1) 828 } 829 830 // Dims returns the dimensions of the matrix. 831 func (ch *BandCholesky) Dims() (r, c int) { 832 if !ch.valid() { 833 panic(badCholesky) 834 } 835 r, c = ch.chol.Dims() 836 return r, c 837 } 838 839 // At returns the element at row i, column j. 840 func (ch *BandCholesky) At(i, j int) float64 { 841 if !ch.valid() { 842 panic(badCholesky) 843 } 844 n, k, _ := ch.chol.TriBand() 845 if uint(i) >= uint(n) { 846 panic(ErrRowAccess) 847 } 848 if uint(j) >= uint(n) { 849 panic(ErrColAccess) 850 } 851 852 if i > j { 853 i, j = j, i 854 } 855 if j-i > k { 856 return 0 857 } 858 var aij float64 859 for k := max(0, j-k); k <= i; k++ { 860 aij += ch.chol.at(k, i) * ch.chol.at(k, j) 861 } 862 return aij 863 } 864 865 // T returns the receiver, the transpose of a symmetric matrix. 866 func (ch *BandCholesky) T() Matrix { 867 return ch 868 } 869 870 // TBand returns the receiver, the transpose of a symmetric band matrix. 871 func (ch *BandCholesky) TBand() Banded { 872 return ch 873 } 874 875 // Symmetric implements the Symmetric interface and returns the number of rows 876 // in the matrix (this is also the number of columns). 877 func (ch *BandCholesky) Symmetric() int { 878 n, _ := ch.chol.Triangle() 879 return n 880 } 881 882 // Bandwidth returns the lower and upper bandwidth values for the matrix. 883 // The total bandwidth of the matrix is kl+ku+1. 884 func (ch *BandCholesky) Bandwidth() (kl, ku int) { 885 _, k, _ := ch.chol.TriBand() 886 return k, k 887 } 888 889 // SymBand returns the number of rows/columns in the matrix, and the size of the 890 // bandwidth. The total bandwidth of the matrix is 2*k+1. 891 func (ch *BandCholesky) SymBand() (n, k int) { 892 n, k, _ = ch.chol.TriBand() 893 return n, k 894 } 895 896 // IsEmpty returns whether the receiver is empty. Empty matrices can be the 897 // receiver for dimensionally restricted operations. The receiver can be emptied 898 // using Reset. 899 func (ch *BandCholesky) IsEmpty() bool { 900 return ch == nil || ch.chol.IsEmpty() 901 } 902 903 // Det returns the determinant of the matrix that has been factorized. 904 func (ch *BandCholesky) Det() float64 { 905 if !ch.valid() { 906 panic(badCholesky) 907 } 908 return math.Exp(ch.LogDet()) 909 } 910 911 // LogDet returns the log of the determinant of the matrix that has been factorized. 912 func (ch *BandCholesky) LogDet() float64 { 913 if !ch.valid() { 914 panic(badCholesky) 915 } 916 var det float64 917 for i := 0; i < ch.chol.mat.N; i++ { 918 det += 2 * math.Log(ch.chol.mat.Data[i*ch.chol.mat.Stride]) 919 } 920 return det 921 } 922 923 func (ch *BandCholesky) valid() bool { 924 return ch.chol != nil && !ch.chol.IsEmpty() 925 }