github.com/gonum/matrix@v0.0.0-20181209220409-c518dec07be9/mat64/qr.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 mat64 6 7 import ( 8 "math" 9 10 "github.com/gonum/blas" 11 "github.com/gonum/blas/blas64" 12 "github.com/gonum/lapack/lapack64" 13 "github.com/gonum/matrix" 14 ) 15 16 // QR is a type for creating and using the QR factorization of a matrix. 17 type QR struct { 18 qr *Dense 19 tau []float64 20 cond float64 21 } 22 23 func (qr *QR) updateCond() { 24 // A = QR, where Q is orthonormal. Orthonormal multiplications do not change 25 // the condition number. Thus, ||A|| = ||Q|| ||R|| = ||R||. 26 n := qr.qr.mat.Cols 27 work := make([]float64, 3*n) 28 iwork := make([]int, n) 29 r := qr.qr.asTriDense(n, blas.NonUnit, blas.Upper) 30 v := lapack64.Trcon(matrix.CondNorm, r.mat, work, iwork) 31 qr.cond = 1 / v 32 } 33 34 // Factorize computes the QR factorization of an m×n matrix a where m >= n. The QR 35 // factorization always exists even if A is singular. 36 // 37 // The QR decomposition is a factorization of the matrix A such that A = Q * R. 38 // The matrix Q is an orthonormal m×m matrix, and R is an m×n upper triangular matrix. 39 // Q and R can be extracted from the QFromQR and RFromQR methods on Dense. 40 func (qr *QR) Factorize(a Matrix) { 41 m, n := a.Dims() 42 if m < n { 43 panic(matrix.ErrShape) 44 } 45 k := min(m, n) 46 if qr.qr == nil { 47 qr.qr = &Dense{} 48 } 49 qr.qr.Clone(a) 50 work := make([]float64, 1) 51 qr.tau = make([]float64, k) 52 lapack64.Geqrf(qr.qr.mat, qr.tau, work, -1) 53 54 work = make([]float64, int(work[0])) 55 lapack64.Geqrf(qr.qr.mat, qr.tau, work, len(work)) 56 qr.updateCond() 57 } 58 59 // TODO(btracey): Add in the "Reduced" forms for extracting the n×n orthogonal 60 // and upper triangular matrices. 61 62 // RFromQR extracts the m×n upper trapezoidal matrix from a QR decomposition. 63 func (m *Dense) RFromQR(qr *QR) { 64 r, c := qr.qr.Dims() 65 m.reuseAs(r, c) 66 67 // Disguise the QR as an upper triangular 68 t := &TriDense{ 69 mat: blas64.Triangular{ 70 N: c, 71 Stride: qr.qr.mat.Stride, 72 Data: qr.qr.mat.Data, 73 Uplo: blas.Upper, 74 Diag: blas.NonUnit, 75 }, 76 cap: qr.qr.capCols, 77 } 78 m.Copy(t) 79 80 // Zero below the triangular. 81 for i := r; i < c; i++ { 82 zero(m.mat.Data[i*m.mat.Stride : i*m.mat.Stride+c]) 83 } 84 } 85 86 // QFromQR extracts the m×m orthonormal matrix Q from a QR decomposition. 87 func (m *Dense) QFromQR(qr *QR) { 88 r, _ := qr.qr.Dims() 89 m.reuseAsZeroed(r, r) 90 91 // Set Q = I. 92 for i := 0; i < r*r; i += r + 1 { 93 m.mat.Data[i] = 1 94 } 95 96 // Construct Q from the elementary reflectors. 97 work := make([]float64, 1) 98 lapack64.Ormqr(blas.Left, blas.NoTrans, qr.qr.mat, qr.tau, m.mat, work, -1) 99 work = make([]float64, int(work[0])) 100 lapack64.Ormqr(blas.Left, blas.NoTrans, qr.qr.mat, qr.tau, m.mat, work, len(work)) 101 } 102 103 // SolveQR finds a minimum-norm solution to a system of linear equations defined 104 // by the matrices A and b, where A is an m×n matrix represented in its QR factorized 105 // form. If A is singular or near-singular a Condition error is returned. Please 106 // see the documentation for Condition for more information. 107 // 108 // The minimization problem solved depends on the input parameters. 109 // If trans == false, find X such that ||A*X - b||_2 is minimized. 110 // If trans == true, find the minimum norm solution of A^T * X = b. 111 // The solution matrix, X, is stored in place into the receiver. 112 func (m *Dense) SolveQR(qr *QR, trans bool, b Matrix) error { 113 r, c := qr.qr.Dims() 114 br, bc := b.Dims() 115 116 // The QR solve algorithm stores the result in-place into the right hand side. 117 // The storage for the answer must be large enough to hold both b and x. 118 // However, this method's receiver must be the size of x. Copy b, and then 119 // copy the result into m at the end. 120 if trans { 121 if c != br { 122 panic(matrix.ErrShape) 123 } 124 m.reuseAs(r, bc) 125 } else { 126 if r != br { 127 panic(matrix.ErrShape) 128 } 129 m.reuseAs(c, bc) 130 } 131 // Do not need to worry about overlap between m and b because x has its own 132 // independent storage. 133 x := getWorkspace(max(r, c), bc, false) 134 x.Copy(b) 135 t := qr.qr.asTriDense(qr.qr.mat.Cols, blas.NonUnit, blas.Upper).mat 136 if trans { 137 ok := lapack64.Trtrs(blas.Trans, t, x.mat) 138 if !ok { 139 return matrix.Condition(math.Inf(1)) 140 } 141 for i := c; i < r; i++ { 142 zero(x.mat.Data[i*x.mat.Stride : i*x.mat.Stride+bc]) 143 } 144 work := make([]float64, 1) 145 lapack64.Ormqr(blas.Left, blas.NoTrans, qr.qr.mat, qr.tau, x.mat, work, -1) 146 work = make([]float64, int(work[0])) 147 lapack64.Ormqr(blas.Left, blas.NoTrans, qr.qr.mat, qr.tau, x.mat, work, len(work)) 148 } else { 149 work := make([]float64, 1) 150 lapack64.Ormqr(blas.Left, blas.Trans, qr.qr.mat, qr.tau, x.mat, work, -1) 151 work = make([]float64, int(work[0])) 152 lapack64.Ormqr(blas.Left, blas.Trans, qr.qr.mat, qr.tau, x.mat, work, len(work)) 153 154 ok := lapack64.Trtrs(blas.NoTrans, t, x.mat) 155 if !ok { 156 return matrix.Condition(math.Inf(1)) 157 } 158 } 159 // M was set above to be the correct size for the result. 160 m.Copy(x) 161 putWorkspace(x) 162 if qr.cond > matrix.ConditionTolerance { 163 return matrix.Condition(qr.cond) 164 } 165 return nil 166 } 167 168 // SolveQRVec finds a minimum-norm solution to a system of linear equations. 169 // Please see Dense.SolveQR for the full documentation. 170 func (v *Vector) SolveQRVec(qr *QR, trans bool, b *Vector) error { 171 if v != b { 172 v.checkOverlap(b.mat) 173 } 174 r, c := qr.qr.Dims() 175 // The Solve implementation is non-trivial, so rather than duplicate the code, 176 // instead recast the Vectors as Dense and call the matrix code. 177 if trans { 178 v.reuseAs(r) 179 } else { 180 v.reuseAs(c) 181 } 182 return v.asDense().SolveQR(qr, trans, b.asDense()) 183 }