gonum.org/v1/gonum@v0.14.0/mat/svd.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 "gonum.org/v1/gonum/blas/blas64" 9 "gonum.org/v1/gonum/lapack" 10 "gonum.org/v1/gonum/lapack/lapack64" 11 ) 12 13 const badRcond = "mat: invalid rcond value" 14 15 // SVD is a type for creating and using the Singular Value Decomposition 16 // of a matrix. 17 type SVD struct { 18 kind SVDKind 19 20 s []float64 21 u blas64.General 22 vt blas64.General 23 } 24 25 // SVDKind specifies the treatment of singular vectors during an SVD 26 // factorization. 27 type SVDKind int 28 29 const ( 30 // SVDNone specifies that no singular vectors should be computed during 31 // the decomposition. 32 SVDNone SVDKind = 0 33 34 // SVDThinU specifies the thin decomposition for U should be computed. 35 SVDThinU SVDKind = 1 << (iota - 1) 36 // SVDFullU specifies the full decomposition for U should be computed. 37 SVDFullU 38 // SVDThinV specifies the thin decomposition for V should be computed. 39 SVDThinV 40 // SVDFullV specifies the full decomposition for V should be computed. 41 SVDFullV 42 43 // SVDThin is a convenience value for computing both thin vectors. 44 SVDThin SVDKind = SVDThinU | SVDThinV 45 // SVDFull is a convenience value for computing both full vectors. 46 SVDFull SVDKind = SVDFullU | SVDFullV 47 ) 48 49 // succFact returns whether the receiver contains a successful factorization. 50 func (svd *SVD) succFact() bool { 51 return len(svd.s) != 0 52 } 53 54 // Factorize computes the singular value decomposition (SVD) of the input matrix A. 55 // The singular values of A are computed in all cases, while the singular 56 // vectors are optionally computed depending on the input kind. 57 // 58 // The full singular value decomposition (kind == SVDFull) is a factorization 59 // of an m×n matrix A of the form 60 // 61 // A = U * Σ * Vᵀ 62 // 63 // where Σ is an m×n diagonal matrix, U is an m×m orthogonal matrix, and V is an 64 // n×n orthogonal matrix. The diagonal elements of Σ are the singular values of A. 65 // The first min(m,n) columns of U and V are, respectively, the left and right 66 // singular vectors of A. 67 // 68 // Significant storage space can be saved by using the thin representation of 69 // the SVD (kind == SVDThin) instead of the full SVD, especially if 70 // m >> n or m << n. The thin SVD finds 71 // 72 // A = U~ * Σ * V~ᵀ 73 // 74 // where U~ is of size m×min(m,n), Σ is a diagonal matrix of size min(m,n)×min(m,n) 75 // and V~ is of size n×min(m,n). 76 // 77 // Factorize returns whether the decomposition succeeded. If the decomposition 78 // failed, routines that require a successful factorization will panic. 79 func (svd *SVD) Factorize(a Matrix, kind SVDKind) (ok bool) { 80 // kill previous factorization 81 svd.s = svd.s[:0] 82 svd.kind = kind 83 84 m, n := a.Dims() 85 var jobU, jobVT lapack.SVDJob 86 87 // TODO(btracey): This code should be modified to have the smaller 88 // matrix written in-place into aCopy when the lapack/native/dgesvd 89 // implementation is complete. 90 switch { 91 case kind&SVDFullU != 0: 92 jobU = lapack.SVDAll 93 svd.u = blas64.General{ 94 Rows: m, 95 Cols: m, 96 Stride: m, 97 Data: use(svd.u.Data, m*m), 98 } 99 case kind&SVDThinU != 0: 100 jobU = lapack.SVDStore 101 svd.u = blas64.General{ 102 Rows: m, 103 Cols: min(m, n), 104 Stride: min(m, n), 105 Data: use(svd.u.Data, m*min(m, n)), 106 } 107 default: 108 jobU = lapack.SVDNone 109 } 110 switch { 111 case kind&SVDFullV != 0: 112 svd.vt = blas64.General{ 113 Rows: n, 114 Cols: n, 115 Stride: n, 116 Data: use(svd.vt.Data, n*n), 117 } 118 jobVT = lapack.SVDAll 119 case kind&SVDThinV != 0: 120 svd.vt = blas64.General{ 121 Rows: min(m, n), 122 Cols: n, 123 Stride: n, 124 Data: use(svd.vt.Data, min(m, n)*n), 125 } 126 jobVT = lapack.SVDStore 127 default: 128 jobVT = lapack.SVDNone 129 } 130 131 // A is destroyed on call, so copy the matrix. 132 aCopy := DenseCopyOf(a) 133 svd.kind = kind 134 svd.s = use(svd.s, min(m, n)) 135 136 work := []float64{0} 137 lapack64.Gesvd(jobU, jobVT, aCopy.mat, svd.u, svd.vt, svd.s, work, -1) 138 work = getFloat64s(int(work[0]), false) 139 ok = lapack64.Gesvd(jobU, jobVT, aCopy.mat, svd.u, svd.vt, svd.s, work, len(work)) 140 putFloat64s(work) 141 if !ok { 142 svd.kind = 0 143 } 144 return ok 145 } 146 147 // Kind returns the SVDKind of the decomposition. If no decomposition has been 148 // computed, Kind returns -1. 149 func (svd *SVD) Kind() SVDKind { 150 if !svd.succFact() { 151 return -1 152 } 153 return svd.kind 154 } 155 156 // Rank returns the rank of A based on the count of singular values greater than 157 // rcond scaled by the largest singular value. 158 // Rank will panic if the receiver does not contain a successful factorization or 159 // rcond is negative. 160 func (svd *SVD) Rank(rcond float64) int { 161 if rcond < 0 { 162 panic(badRcond) 163 } 164 if !svd.succFact() { 165 panic(badFact) 166 } 167 s0 := svd.s[0] 168 for i, v := range svd.s { 169 if v <= rcond*s0 { 170 return i 171 } 172 } 173 return len(svd.s) 174 } 175 176 // Cond returns the 2-norm condition number for the factorized matrix. Cond will 177 // panic if the receiver does not contain a successful factorization. 178 func (svd *SVD) Cond() float64 { 179 if !svd.succFact() { 180 panic(badFact) 181 } 182 return svd.s[0] / svd.s[len(svd.s)-1] 183 } 184 185 // Values returns the singular values of the factorized matrix in descending order. 186 // 187 // If the input slice is non-nil, the values will be stored in-place into 188 // the slice. In this case, the slice must have length min(m,n), and Values will 189 // panic with ErrSliceLengthMismatch otherwise. If the input slice is nil, a new 190 // slice of the appropriate length will be allocated and returned. 191 // 192 // Values will panic if the receiver does not contain a successful factorization. 193 func (svd *SVD) Values(s []float64) []float64 { 194 if !svd.succFact() { 195 panic(badFact) 196 } 197 if s == nil { 198 s = make([]float64, len(svd.s)) 199 } 200 if len(s) != len(svd.s) { 201 panic(ErrSliceLengthMismatch) 202 } 203 copy(s, svd.s) 204 return s 205 } 206 207 // UTo extracts the matrix U from the singular value decomposition. The first 208 // min(m,n) columns are the left singular vectors and correspond to the singular 209 // values as returned from SVD.Values. 210 // 211 // If dst is empty, UTo will resize dst to be m×m if the full U was computed 212 // and size m×min(m,n) if the thin U was computed. When dst is non-empty, then 213 // UTo will panic if dst is not the appropriate size. UTo will also panic if 214 // the receiver does not contain a successful factorization, or if U was 215 // not computed during factorization. 216 func (svd *SVD) UTo(dst *Dense) { 217 if !svd.succFact() { 218 panic(badFact) 219 } 220 kind := svd.kind 221 if kind&SVDThinU == 0 && kind&SVDFullU == 0 { 222 panic("svd: u not computed during factorization") 223 } 224 r := svd.u.Rows 225 c := svd.u.Cols 226 if dst.IsEmpty() { 227 dst.ReuseAs(r, c) 228 } else { 229 r2, c2 := dst.Dims() 230 if r != r2 || c != c2 { 231 panic(ErrShape) 232 } 233 } 234 235 tmp := &Dense{ 236 mat: svd.u, 237 capRows: r, 238 capCols: c, 239 } 240 dst.Copy(tmp) 241 } 242 243 // VTo extracts the matrix V from the singular value decomposition. The first 244 // min(m,n) columns are the right singular vectors and correspond to the singular 245 // values as returned from SVD.Values. 246 // 247 // If dst is empty, VTo will resize dst to be n×n if the full V was computed 248 // and size n×min(m,n) if the thin V was computed. When dst is non-empty, then 249 // VTo will panic if dst is not the appropriate size. VTo will also panic if 250 // the receiver does not contain a successful factorization, or if V was 251 // not computed during factorization. 252 func (svd *SVD) VTo(dst *Dense) { 253 if !svd.succFact() { 254 panic(badFact) 255 } 256 kind := svd.kind 257 if kind&SVDThinV == 0 && kind&SVDFullV == 0 { 258 panic("svd: v not computed during factorization") 259 } 260 r := svd.vt.Rows 261 c := svd.vt.Cols 262 if dst.IsEmpty() { 263 dst.ReuseAs(c, r) 264 } else { 265 r2, c2 := dst.Dims() 266 if c != r2 || r != c2 { 267 panic(ErrShape) 268 } 269 } 270 271 tmp := &Dense{ 272 mat: svd.vt, 273 capRows: r, 274 capCols: c, 275 } 276 dst.Copy(tmp.T()) 277 } 278 279 // SolveTo calculates the minimum-norm solution to a linear least squares problem 280 // 281 // minimize over n-element vectors x: |b - A*x|_2 and |x|_2 282 // 283 // where b is a given m-element vector, using the SVD of m×n matrix A stored in 284 // the receiver. A may be rank-deficient, that is, the given effective rank can be 285 // 286 // rank ≤ min(m,n) 287 // 288 // The rank can be computed using SVD.Rank. 289 // 290 // Several right-hand side vectors b and solution vectors x can be handled in a 291 // single call. Vectors b are stored in the columns of the m×k matrix B and the 292 // resulting vectors x will be stored in the columns of dst. dst must be either 293 // empty or have the size equal to n×k. 294 // 295 // The decomposition must have been factorized computing both the U and V 296 // singular vectors. 297 // 298 // SolveTo returns the residuals calculated from the complete SVD. For this 299 // value to be valid the factorization must have been performed with at least 300 // SVDFullU. 301 func (svd *SVD) SolveTo(dst *Dense, b Matrix, rank int) []float64 { 302 if !svd.succFact() { 303 panic(badFact) 304 } 305 if rank < 1 || len(svd.s) < rank { 306 panic("svd: rank out of range") 307 } 308 kind := svd.kind 309 if kind&SVDThinU == 0 && kind&SVDFullU == 0 { 310 panic("svd: u not computed during factorization") 311 } 312 if kind&SVDThinV == 0 && kind&SVDFullV == 0 { 313 panic("svd: v not computed during factorization") 314 } 315 316 u := Dense{ 317 mat: svd.u, 318 capRows: svd.u.Rows, 319 capCols: svd.u.Cols, 320 } 321 vt := Dense{ 322 mat: svd.vt, 323 capRows: svd.vt.Rows, 324 capCols: svd.vt.Cols, 325 } 326 s := svd.s[:rank] 327 328 _, bc := b.Dims() 329 c := getDenseWorkspace(svd.u.Cols, bc, false) 330 defer putDenseWorkspace(c) 331 c.Mul(u.T(), b) 332 333 y := getDenseWorkspace(rank, bc, false) 334 defer putDenseWorkspace(y) 335 y.DivElem(c.slice(0, rank, 0, bc), repVector{vec: s, cols: bc}) 336 dst.Mul(vt.slice(0, rank, 0, svd.vt.Cols).T(), y) 337 338 res := make([]float64, bc) 339 if rank < svd.u.Cols { 340 c = c.slice(len(s), svd.u.Cols, 0, bc) 341 for j := range res { 342 col := c.ColView(j) 343 res[j] = Dot(col, col) 344 } 345 } 346 return res 347 } 348 349 type repVector struct { 350 vec []float64 351 cols int 352 } 353 354 func (m repVector) Dims() (r, c int) { return len(m.vec), m.cols } 355 func (m repVector) At(i, j int) float64 { 356 if i < 0 || len(m.vec) <= i || j < 0 || m.cols <= j { 357 panic(ErrIndexOutOfRange.string) // Panic with string to prevent mat.Error recovery. 358 } 359 return m.vec[i] 360 } 361 func (m repVector) T() Matrix { return Transpose{m} } 362 363 // SolveVecTo calculates the minimum-norm solution to a linear least squares problem 364 // 365 // minimize over n-element vectors x: |b - A*x|_2 and |x|_2 366 // 367 // where b is a given m-element vector, using the SVD of m×n matrix A stored in 368 // the receiver. A may be rank-deficient, that is, the given effective rank can be 369 // 370 // rank ≤ min(m,n) 371 // 372 // The rank can be computed using SVD.Rank. 373 // 374 // The resulting vector x will be stored in dst. dst must be either empty or 375 // have length equal to n. 376 // 377 // The decomposition must have been factorized computing both the U and V 378 // singular vectors. 379 // 380 // SolveVecTo returns the residuals calculated from the complete SVD. For this 381 // value to be valid the factorization must have been performed with at least 382 // SVDFullU. 383 func (svd *SVD) SolveVecTo(dst *VecDense, b Vector, rank int) float64 { 384 if !svd.succFact() { 385 panic(badFact) 386 } 387 if rank < 1 || len(svd.s) < rank { 388 panic("svd: rank out of range") 389 } 390 kind := svd.kind 391 if kind&SVDThinU == 0 && kind&SVDFullU == 0 { 392 panic("svd: u not computed during factorization") 393 } 394 if kind&SVDThinV == 0 && kind&SVDFullV == 0 { 395 panic("svd: v not computed during factorization") 396 } 397 398 u := Dense{ 399 mat: svd.u, 400 capRows: svd.u.Rows, 401 capCols: svd.u.Cols, 402 } 403 vt := Dense{ 404 mat: svd.vt, 405 capRows: svd.vt.Rows, 406 capCols: svd.vt.Cols, 407 } 408 s := svd.s[:rank] 409 410 c := getVecDenseWorkspace(svd.u.Cols, false) 411 defer putVecDenseWorkspace(c) 412 c.MulVec(u.T(), b) 413 414 y := getVecDenseWorkspace(rank, false) 415 defer putVecDenseWorkspace(y) 416 y.DivElemVec(c.sliceVec(0, rank), NewVecDense(rank, s)) 417 dst.MulVec(vt.slice(0, rank, 0, svd.vt.Cols).T(), y) 418 419 var res float64 420 if rank < c.Len() { 421 c = c.sliceVec(rank, c.Len()) 422 res = Dot(c, c) 423 } 424 return res 425 }