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