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