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