gonum.org/v1/gonum@v0.14.0/lapack/lapack64/lapack64.go (about) 1 // Copyright ©2015 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 lapack64 6 7 import ( 8 "gonum.org/v1/gonum/blas" 9 "gonum.org/v1/gonum/blas/blas64" 10 "gonum.org/v1/gonum/lapack" 11 "gonum.org/v1/gonum/lapack/gonum" 12 ) 13 14 var lapack64 lapack.Float64 = gonum.Implementation{} 15 16 // Use sets the LAPACK float64 implementation to be used by subsequent BLAS calls. 17 // The default implementation is native.Implementation. 18 func Use(l lapack.Float64) { 19 lapack64 = l 20 } 21 22 // Tridiagonal represents a tridiagonal matrix using its three diagonals. 23 type Tridiagonal struct { 24 N int 25 DL []float64 26 D []float64 27 DU []float64 28 } 29 30 func max(a, b int) int { 31 if a > b { 32 return a 33 } 34 return b 35 } 36 37 // Potrf computes the Cholesky factorization of a. 38 // The factorization has the form 39 // 40 // A = Uᵀ * U if a.Uplo == blas.Upper, or 41 // A = L * Lᵀ if a.Uplo == blas.Lower, 42 // 43 // where U is an upper triangular matrix and L is lower triangular. 44 // The triangular matrix is returned in t, and the underlying data between 45 // a and t is shared. The returned bool indicates whether a is positive 46 // definite and the factorization could be finished. 47 func Potrf(a blas64.Symmetric) (t blas64.Triangular, ok bool) { 48 ok = lapack64.Dpotrf(a.Uplo, a.N, a.Data, max(1, a.Stride)) 49 t.Uplo = a.Uplo 50 t.N = a.N 51 t.Data = a.Data 52 t.Stride = a.Stride 53 t.Diag = blas.NonUnit 54 return 55 } 56 57 // Potri computes the inverse of a real symmetric positive definite matrix A 58 // using its Cholesky factorization. 59 // 60 // On entry, t contains the triangular factor U or L from the Cholesky 61 // factorization A = Uᵀ*U or A = L*Lᵀ, as computed by Potrf. 62 // 63 // On return, the upper or lower triangle of the (symmetric) inverse of A is 64 // stored in t, overwriting the input factor U or L, and also returned in a. The 65 // underlying data between a and t is shared. 66 // 67 // The returned bool indicates whether the inverse was computed successfully. 68 func Potri(t blas64.Triangular) (a blas64.Symmetric, ok bool) { 69 ok = lapack64.Dpotri(t.Uplo, t.N, t.Data, max(1, t.Stride)) 70 a.Uplo = t.Uplo 71 a.N = t.N 72 a.Data = t.Data 73 a.Stride = t.Stride 74 return 75 } 76 77 // Potrs solves a system of n linear equations A*X = B where A is an n×n 78 // symmetric positive definite matrix and B is an n×nrhs matrix, using the 79 // Cholesky factorization A = Uᵀ*U or A = L*Lᵀ. t contains the corresponding 80 // triangular factor as returned by Potrf. On entry, B contains the right-hand 81 // side matrix B, on return it contains the solution matrix X. 82 func Potrs(t blas64.Triangular, b blas64.General) { 83 lapack64.Dpotrs(t.Uplo, t.N, b.Cols, t.Data, max(1, t.Stride), b.Data, max(1, b.Stride)) 84 } 85 86 // Pbcon returns an estimate of the reciprocal of the condition number (in the 87 // 1-norm) of an n×n symmetric positive definite band matrix using the Cholesky 88 // factorization 89 // 90 // A = Uᵀ*U if uplo == blas.Upper 91 // A = L*Lᵀ if uplo == blas.Lower 92 // 93 // computed by Pbtrf. The estimate is obtained for norm(inv(A)), and the 94 // reciprocal of the condition number is computed as 95 // 96 // rcond = 1 / (anorm * norm(inv(A))). 97 // 98 // The length of work must be at least 3*n and the length of iwork must be at 99 // least n. 100 func Pbcon(a blas64.SymmetricBand, anorm float64, work []float64, iwork []int) float64 { 101 return lapack64.Dpbcon(a.Uplo, a.N, a.K, a.Data, a.Stride, anorm, work, iwork) 102 } 103 104 // Pbtrf computes the Cholesky factorization of an n×n symmetric positive 105 // definite band matrix 106 // 107 // A = Uᵀ * U if a.Uplo == blas.Upper 108 // A = L * Lᵀ if a.Uplo == blas.Lower 109 // 110 // where U and L are upper, respectively lower, triangular band matrices. 111 // 112 // The triangular matrix U or L is returned in t, and the underlying data 113 // between a and t is shared. The returned bool indicates whether A is positive 114 // definite and the factorization could be finished. 115 func Pbtrf(a blas64.SymmetricBand) (t blas64.TriangularBand, ok bool) { 116 ok = lapack64.Dpbtrf(a.Uplo, a.N, a.K, a.Data, max(1, a.Stride)) 117 t.Uplo = a.Uplo 118 t.Diag = blas.NonUnit 119 t.N = a.N 120 t.K = a.K 121 t.Data = a.Data 122 t.Stride = a.Stride 123 return t, ok 124 } 125 126 // Pbtrs solves a system of linear equations A*X = B with an n×n symmetric 127 // positive definite band matrix A using the Cholesky factorization 128 // 129 // A = Uᵀ * U if t.Uplo == blas.Upper 130 // A = L * Lᵀ if t.Uplo == blas.Lower 131 // 132 // t contains the corresponding triangular factor as returned by Pbtrf. 133 // 134 // On entry, b contains the right hand side matrix B. On return, it is 135 // overwritten with the solution matrix X. 136 func Pbtrs(t blas64.TriangularBand, b blas64.General) { 137 lapack64.Dpbtrs(t.Uplo, t.N, t.K, b.Cols, t.Data, max(1, t.Stride), b.Data, max(1, b.Stride)) 138 } 139 140 // Pstrf computes the Cholesky factorization with complete pivoting of an n×n 141 // symmetric positive semidefinite matrix A. 142 // 143 // The factorization has the form 144 // 145 // Pᵀ * A * P = Uᵀ * U , if a.Uplo = blas.Upper, 146 // Pᵀ * A * P = L * Lᵀ, if a.Uplo = blas.Lower, 147 // 148 // where U is an upper triangular matrix, L is lower triangular, and P is a 149 // permutation matrix. 150 // 151 // tol is a user-defined tolerance. The algorithm terminates if the pivot is 152 // less than or equal to tol. If tol is negative, then n*eps*max(A[k,k]) will be 153 // used instead. 154 // 155 // The triangular factor U or L from the Cholesky factorization is returned in t 156 // and the underlying data between a and t is shared. P is stored on return in 157 // vector piv such that P[piv[k],k] = 1. 158 // 159 // Pstrf returns the computed rank of A and whether the factorization can be 160 // used to solve a system. Pstrf does not attempt to check that A is positive 161 // semi-definite, so if ok is false, the matrix A is either rank deficient or is 162 // not positive semidefinite. 163 // 164 // The length of piv must be n and the length of work must be at least 2*n, 165 // otherwise Pstrf will panic. 166 func Pstrf(a blas64.Symmetric, piv []int, tol float64, work []float64) (t blas64.Triangular, rank int, ok bool) { 167 rank, ok = lapack64.Dpstrf(a.Uplo, a.N, a.Data, max(1, a.Stride), piv, tol, work) 168 t.Uplo = a.Uplo 169 t.Diag = blas.NonUnit 170 t.N = a.N 171 t.Data = a.Data 172 t.Stride = a.Stride 173 return t, rank, ok 174 } 175 176 // Gecon estimates the reciprocal of the condition number of the n×n matrix A 177 // given the LU decomposition of the matrix. The condition number computed may 178 // be based on the 1-norm or the ∞-norm. 179 // 180 // a contains the result of the LU decomposition of A as computed by Getrf. 181 // 182 // anorm is the corresponding 1-norm or ∞-norm of the original matrix A. 183 // 184 // work is a temporary data slice of length at least 4*n and Gecon will panic otherwise. 185 // 186 // iwork is a temporary data slice of length at least n and Gecon will panic otherwise. 187 func Gecon(norm lapack.MatrixNorm, a blas64.General, anorm float64, work []float64, iwork []int) float64 { 188 return lapack64.Dgecon(norm, a.Cols, a.Data, max(1, a.Stride), anorm, work, iwork) 189 } 190 191 // Gels finds a minimum-norm solution based on the matrices A and B using the 192 // QR or LQ factorization. Gels returns false if the matrix 193 // A is singular, and true if this solution was successfully found. 194 // 195 // The minimization problem solved depends on the input parameters. 196 // 197 // 1. If m >= n and trans == blas.NoTrans, Gels finds X such that || A*X - B||_2 198 // is minimized. 199 // 2. If m < n and trans == blas.NoTrans, Gels finds the minimum norm solution of 200 // A * X = B. 201 // 3. If m >= n and trans == blas.Trans, Gels finds the minimum norm solution of 202 // Aᵀ * X = B. 203 // 4. If m < n and trans == blas.Trans, Gels finds X such that || A*X - B||_2 204 // is minimized. 205 // 206 // Note that the least-squares solutions (cases 1 and 3) perform the minimization 207 // per column of B. This is not the same as finding the minimum-norm matrix. 208 // 209 // The matrix A is a general matrix of size m×n and is modified during this call. 210 // The input matrix B is of size max(m,n)×nrhs, and serves two purposes. On entry, 211 // the elements of b specify the input matrix B. B has size m×nrhs if 212 // trans == blas.NoTrans, and n×nrhs if trans == blas.Trans. On exit, the 213 // leading submatrix of b contains the solution vectors X. If trans == blas.NoTrans, 214 // this submatrix is of size n×nrhs, and of size m×nrhs otherwise. 215 // 216 // Work is temporary storage, and lwork specifies the usable memory length. 217 // At minimum, lwork >= max(m,n) + max(m,n,nrhs), and this function will panic 218 // otherwise. A longer work will enable blocked algorithms to be called. 219 // In the special case that lwork == -1, work[0] will be set to the optimal working 220 // length. 221 func Gels(trans blas.Transpose, a blas64.General, b blas64.General, work []float64, lwork int) bool { 222 return lapack64.Dgels(trans, a.Rows, a.Cols, b.Cols, a.Data, max(1, a.Stride), b.Data, max(1, b.Stride), work, lwork) 223 } 224 225 // Geqrf computes the QR factorization of the m×n matrix A using a blocked 226 // algorithm. A is modified to contain the information to construct Q and R. 227 // The upper triangle of a contains the matrix R. The lower triangular elements 228 // (not including the diagonal) contain the elementary reflectors. tau is modified 229 // to contain the reflector scales. tau must have length at least min(m,n), and 230 // this function will panic otherwise. 231 // 232 // The ith elementary reflector can be explicitly constructed by first extracting 233 // the 234 // 235 // v[j] = 0 j < i 236 // v[j] = 1 j == i 237 // v[j] = a[j*lda+i] j > i 238 // 239 // and computing H_i = I - tau[i] * v * vᵀ. 240 // 241 // The orthonormal matrix Q can be constructed from a product of these elementary 242 // reflectors, Q = H_0 * H_1 * ... * H_{k-1}, where k = min(m,n). 243 // 244 // Work is temporary storage, and lwork specifies the usable memory length. 245 // At minimum, lwork >= m and this function will panic otherwise. 246 // Geqrf is a blocked QR factorization, but the block size is limited 247 // by the temporary space available. If lwork == -1, instead of performing Geqrf, 248 // the optimal work length will be stored into work[0]. 249 func Geqrf(a blas64.General, tau, work []float64, lwork int) { 250 lapack64.Dgeqrf(a.Rows, a.Cols, a.Data, max(1, a.Stride), tau, work, lwork) 251 } 252 253 // Gelqf computes the LQ factorization of the m×n matrix A using a blocked 254 // algorithm. A is modified to contain the information to construct L and Q. The 255 // lower triangle of a contains the matrix L. The elements above the diagonal 256 // and the slice tau represent the matrix Q. tau is modified to contain the 257 // reflector scales. tau must have length at least min(m,n), and this function 258 // will panic otherwise. 259 // 260 // See Geqrf for a description of the elementary reflectors and orthonormal 261 // matrix Q. Q is constructed as a product of these elementary reflectors, 262 // Q = H_{k-1} * ... * H_1 * H_0. 263 // 264 // Work is temporary storage, and lwork specifies the usable memory length. 265 // At minimum, lwork >= m and this function will panic otherwise. 266 // Gelqf is a blocked LQ factorization, but the block size is limited 267 // by the temporary space available. If lwork == -1, instead of performing Gelqf, 268 // the optimal work length will be stored into work[0]. 269 func Gelqf(a blas64.General, tau, work []float64, lwork int) { 270 lapack64.Dgelqf(a.Rows, a.Cols, a.Data, max(1, a.Stride), tau, work, lwork) 271 } 272 273 // Gesvd computes the singular value decomposition of the input matrix A. 274 // 275 // The singular value decomposition is 276 // 277 // A = U * Sigma * Vᵀ 278 // 279 // where Sigma is an m×n diagonal matrix containing the singular values of A, 280 // U is an m×m orthogonal matrix and V is an n×n orthogonal matrix. The first 281 // min(m,n) columns of U and V are the left and right singular vectors of A 282 // respectively. 283 // 284 // jobU and jobVT are options for computing the singular vectors. The behavior 285 // is as follows 286 // 287 // jobU == lapack.SVDAll All m columns of U are returned in u 288 // jobU == lapack.SVDStore The first min(m,n) columns are returned in u 289 // jobU == lapack.SVDOverwrite The first min(m,n) columns of U are written into a 290 // jobU == lapack.SVDNone The columns of U are not computed. 291 // 292 // The behavior is the same for jobVT and the rows of Vᵀ. At most one of jobU 293 // and jobVT can equal lapack.SVDOverwrite, and Gesvd will panic otherwise. 294 // 295 // On entry, a contains the data for the m×n matrix A. During the call to Gesvd 296 // the data is overwritten. On exit, A contains the appropriate singular vectors 297 // if either job is lapack.SVDOverwrite. 298 // 299 // s is a slice of length at least min(m,n) and on exit contains the singular 300 // values in decreasing order. 301 // 302 // u contains the left singular vectors on exit, stored columnwise. If 303 // jobU == lapack.SVDAll, u is of size m×m. If jobU == lapack.SVDStore u is 304 // of size m×min(m,n). If jobU == lapack.SVDOverwrite or lapack.SVDNone, u is 305 // not used. 306 // 307 // vt contains the left singular vectors on exit, stored rowwise. If 308 // jobV == lapack.SVDAll, vt is of size n×m. If jobVT == lapack.SVDStore vt is 309 // of size min(m,n)×n. If jobVT == lapack.SVDOverwrite or lapack.SVDNone, vt is 310 // not used. 311 // 312 // work is a slice for storing temporary memory, and lwork is the usable size of 313 // the slice. lwork must be at least max(5*min(m,n), 3*min(m,n)+max(m,n)). 314 // If lwork == -1, instead of performing Gesvd, the optimal work length will be 315 // stored into work[0]. Gesvd will panic if the working memory has insufficient 316 // storage. 317 // 318 // Gesvd returns whether the decomposition successfully completed. 319 func Gesvd(jobU, jobVT lapack.SVDJob, a, u, vt blas64.General, s, work []float64, lwork int) (ok bool) { 320 return lapack64.Dgesvd(jobU, jobVT, a.Rows, a.Cols, a.Data, max(1, a.Stride), s, u.Data, max(1, u.Stride), vt.Data, max(1, vt.Stride), work, lwork) 321 } 322 323 // Getrf computes the LU decomposition of the m×n matrix A. 324 // The LU decomposition is a factorization of A into 325 // 326 // A = P * L * U 327 // 328 // where P is a permutation matrix, L is a unit lower triangular matrix, and 329 // U is a (usually) non-unit upper triangular matrix. On exit, L and U are stored 330 // in place into a. 331 // 332 // ipiv is a permutation vector. It indicates that row i of the matrix was 333 // changed with ipiv[i]. ipiv must have length at least min(m,n), and will panic 334 // otherwise. ipiv is zero-indexed. 335 // 336 // Getrf is the blocked version of the algorithm. 337 // 338 // Getrf returns whether the matrix A is singular. The LU decomposition will 339 // be computed regardless of the singularity of A, but division by zero 340 // will occur if the false is returned and the result is used to solve a 341 // system of equations. 342 func Getrf(a blas64.General, ipiv []int) bool { 343 return lapack64.Dgetrf(a.Rows, a.Cols, a.Data, max(1, a.Stride), ipiv) 344 } 345 346 // Getri computes the inverse of the matrix A using the LU factorization computed 347 // by Getrf. On entry, a contains the PLU decomposition of A as computed by 348 // Getrf and on exit contains the reciprocal of the original matrix. 349 // 350 // Getri will not perform the inversion if the matrix is singular, and returns 351 // a boolean indicating whether the inversion was successful. 352 // 353 // Work is temporary storage, and lwork specifies the usable memory length. 354 // At minimum, lwork >= n and this function will panic otherwise. 355 // Getri is a blocked inversion, but the block size is limited 356 // by the temporary space available. If lwork == -1, instead of performing Getri, 357 // the optimal work length will be stored into work[0]. 358 func Getri(a blas64.General, ipiv []int, work []float64, lwork int) (ok bool) { 359 return lapack64.Dgetri(a.Cols, a.Data, max(1, a.Stride), ipiv, work, lwork) 360 } 361 362 // Getrs solves a system of equations using an LU factorization. 363 // The system of equations solved is 364 // 365 // A * X = B if trans == blas.Trans 366 // Aᵀ * X = B if trans == blas.NoTrans 367 // 368 // A is a general n×n matrix with stride lda. B is a general matrix of size n×nrhs. 369 // 370 // On entry b contains the elements of the matrix B. On exit, b contains the 371 // elements of X, the solution to the system of equations. 372 // 373 // a and ipiv contain the LU factorization of A and the permutation indices as 374 // computed by Getrf. ipiv is zero-indexed. 375 func Getrs(trans blas.Transpose, a blas64.General, b blas64.General, ipiv []int) { 376 lapack64.Dgetrs(trans, a.Cols, b.Cols, a.Data, max(1, a.Stride), ipiv, b.Data, max(1, b.Stride)) 377 } 378 379 // Ggsvd3 computes the generalized singular value decomposition (GSVD) 380 // of an m×n matrix A and p×n matrix B: 381 // 382 // Uᵀ*A*Q = D1*[ 0 R ] 383 // 384 // Vᵀ*B*Q = D2*[ 0 R ] 385 // 386 // where U, V and Q are orthogonal matrices. 387 // 388 // Ggsvd3 returns k and l, the dimensions of the sub-blocks. k+l 389 // is the effective numerical rank of the (m+p)×n matrix [ Aᵀ Bᵀ ]ᵀ. 390 // R is a (k+l)×(k+l) nonsingular upper triangular matrix, D1 and 391 // D2 are m×(k+l) and p×(k+l) diagonal matrices and of the following 392 // structures, respectively: 393 // 394 // If m-k-l >= 0, 395 // 396 // k l 397 // D1 = k [ I 0 ] 398 // l [ 0 C ] 399 // m-k-l [ 0 0 ] 400 // 401 // k l 402 // D2 = l [ 0 S ] 403 // p-l [ 0 0 ] 404 // 405 // n-k-l k l 406 // [ 0 R ] = k [ 0 R11 R12 ] k 407 // l [ 0 0 R22 ] l 408 // 409 // where 410 // 411 // C = diag( alpha_k, ... , alpha_{k+l} ), 412 // S = diag( beta_k, ... , beta_{k+l} ), 413 // C^2 + S^2 = I. 414 // 415 // R is stored in 416 // 417 // A[0:k+l, n-k-l:n] 418 // 419 // on exit. 420 // 421 // If m-k-l < 0, 422 // 423 // k m-k k+l-m 424 // D1 = k [ I 0 0 ] 425 // m-k [ 0 C 0 ] 426 // 427 // k m-k k+l-m 428 // D2 = m-k [ 0 S 0 ] 429 // k+l-m [ 0 0 I ] 430 // p-l [ 0 0 0 ] 431 // 432 // n-k-l k m-k k+l-m 433 // [ 0 R ] = k [ 0 R11 R12 R13 ] 434 // m-k [ 0 0 R22 R23 ] 435 // k+l-m [ 0 0 0 R33 ] 436 // 437 // where 438 // 439 // C = diag( alpha_k, ... , alpha_m ), 440 // S = diag( beta_k, ... , beta_m ), 441 // C^2 + S^2 = I. 442 // 443 // R = [ R11 R12 R13 ] is stored in A[1:m, n-k-l+1:n] 444 // [ 0 R22 R23 ] 445 // 446 // and R33 is stored in 447 // 448 // B[m-k:l, n+m-k-l:n] on exit. 449 // 450 // Ggsvd3 computes C, S, R, and optionally the orthogonal transformation 451 // matrices U, V and Q. 452 // 453 // jobU, jobV and jobQ are options for computing the orthogonal matrices. The behavior 454 // is as follows 455 // 456 // jobU == lapack.GSVDU Compute orthogonal matrix U 457 // jobU == lapack.GSVDNone Do not compute orthogonal matrix. 458 // 459 // The behavior is the same for jobV and jobQ with the exception that instead of 460 // lapack.GSVDU these accept lapack.GSVDV and lapack.GSVDQ respectively. 461 // The matrices U, V and Q must be m×m, p×p and n×n respectively unless the 462 // relevant job parameter is lapack.GSVDNone. 463 // 464 // alpha and beta must have length n or Ggsvd3 will panic. On exit, alpha and 465 // beta contain the generalized singular value pairs of A and B 466 // 467 // alpha[0:k] = 1, 468 // beta[0:k] = 0, 469 // 470 // if m-k-l >= 0, 471 // 472 // alpha[k:k+l] = diag(C), 473 // beta[k:k+l] = diag(S), 474 // 475 // if m-k-l < 0, 476 // 477 // alpha[k:m]= C, alpha[m:k+l]= 0 478 // beta[k:m] = S, beta[m:k+l] = 1. 479 // 480 // if k+l < n, 481 // 482 // alpha[k+l:n] = 0 and 483 // beta[k+l:n] = 0. 484 // 485 // On exit, iwork contains the permutation required to sort alpha descending. 486 // 487 // iwork must have length n, work must have length at least max(1, lwork), and 488 // lwork must be -1 or greater than n, otherwise Ggsvd3 will panic. If 489 // lwork is -1, work[0] holds the optimal lwork on return, but Ggsvd3 does 490 // not perform the GSVD. 491 func Ggsvd3(jobU, jobV, jobQ lapack.GSVDJob, a, b blas64.General, alpha, beta []float64, u, v, q blas64.General, work []float64, lwork int, iwork []int) (k, l int, ok bool) { 492 return lapack64.Dggsvd3(jobU, jobV, jobQ, a.Rows, a.Cols, b.Rows, a.Data, max(1, a.Stride), b.Data, max(1, b.Stride), alpha, beta, u.Data, max(1, u.Stride), v.Data, max(1, v.Stride), q.Data, max(1, q.Stride), work, lwork, iwork) 493 } 494 495 // Gtsv solves one of the equations 496 // 497 // A * X = B if trans == blas.NoTrans 498 // Aᵀ * X = B if trans == blas.Trans or blas.ConjTrans 499 // 500 // where A is an n×n tridiagonal matrix. It uses Gaussian elimination with 501 // partial pivoting. 502 // 503 // On entry, a contains the matrix A, on return it will be overwritten. 504 // 505 // On entry, b contains the n×nrhs right-hand side matrix B. On return, it will 506 // be overwritten. If ok is true, it will be overwritten by the solution matrix X. 507 // 508 // Gtsv returns whether the solution X has been successfully computed. 509 // 510 // Dgtsv is not part of the lapack.Float64 interface and so calls to Gtsv are 511 // always executed by the Gonum implementation. 512 func Gtsv(trans blas.Transpose, a Tridiagonal, b blas64.General) (ok bool) { 513 if trans != blas.NoTrans { 514 a.DL, a.DU = a.DU, a.DL 515 } 516 return gonum.Implementation{}.Dgtsv(a.N, b.Cols, a.DL, a.D, a.DU, b.Data, max(1, b.Stride)) 517 } 518 519 // Lagtm performs one of the matrix-matrix operations 520 // 521 // C = alpha * A * B + beta * C if trans == blas.NoTrans 522 // C = alpha * Aᵀ * B + beta * C if trans == blas.Trans or blas.ConjTrans 523 // 524 // where A is an m×m tridiagonal matrix represented by its diagonals dl, d, du, 525 // B and C are m×n dense matrices, and alpha and beta are scalars. 526 // 527 // Dlagtm is not part of the lapack.Float64 interface and so calls to Lagtm are 528 // always executed by the Gonum implementation. 529 func Lagtm(trans blas.Transpose, alpha float64, a Tridiagonal, b blas64.General, beta float64, c blas64.General) { 530 gonum.Implementation{}.Dlagtm(trans, c.Rows, c.Cols, alpha, a.DL, a.D, a.DU, b.Data, max(1, b.Stride), beta, c.Data, max(1, c.Stride)) 531 } 532 533 // Lange computes the matrix norm of the general m×n matrix A. The input norm 534 // specifies the norm computed. 535 // 536 // lapack.MaxAbs: the maximum absolute value of an element. 537 // lapack.MaxColumnSum: the maximum column sum of the absolute values of the entries. 538 // lapack.MaxRowSum: the maximum row sum of the absolute values of the entries. 539 // lapack.Frobenius: the square root of the sum of the squares of the entries. 540 // 541 // If norm == lapack.MaxColumnSum, work must be of length n, and this function will panic otherwise. 542 // There are no restrictions on work for the other matrix norms. 543 func Lange(norm lapack.MatrixNorm, a blas64.General, work []float64) float64 { 544 return lapack64.Dlange(norm, a.Rows, a.Cols, a.Data, max(1, a.Stride), work) 545 } 546 547 // Langb returns the given norm of a general m×n band matrix with kl sub-diagonals and 548 // ku super-diagonals. 549 // 550 // Dlangb is not part of the lapack.Float64 interface and so calls to Langb are always 551 // executed by the Gonum implementation. 552 func Langb(norm lapack.MatrixNorm, a blas64.Band) float64 { 553 return gonum.Implementation{}.Dlangb(norm, a.Rows, a.Cols, a.KL, a.KU, a.Data, max(1, a.Stride)) 554 } 555 556 // Langt computes the specified norm of an n×n tridiagonal matrix. 557 // 558 // Dlangt is not part of the lapack.Float64 interface and so calls to Langt are 559 // always executed by the Gonum implementation. 560 func Langt(norm lapack.MatrixNorm, a Tridiagonal) float64 { 561 return gonum.Implementation{}.Dlangt(norm, a.N, a.DL, a.D, a.DU) 562 } 563 564 // Lansb computes the specified norm of an n×n symmetric band matrix. If 565 // norm == lapack.MaxColumnSum or norm == lapack.MaxRowSum, work must have length 566 // at least n and this function will panic otherwise. 567 // There are no restrictions on work for the other matrix norms. 568 // 569 // Dlansb is not part of the lapack.Float64 interface and so calls to Lansb are always 570 // executed by the Gonum implementation. 571 func Lansb(norm lapack.MatrixNorm, a blas64.SymmetricBand, work []float64) float64 { 572 return gonum.Implementation{}.Dlansb(norm, a.Uplo, a.N, a.K, a.Data, max(1, a.Stride), work) 573 } 574 575 // Lansy computes the specified norm of an n×n symmetric matrix. If 576 // norm == lapack.MaxColumnSum or norm == lapack.MaxRowSum, work must have length 577 // at least n and this function will panic otherwise. 578 // There are no restrictions on work for the other matrix norms. 579 func Lansy(norm lapack.MatrixNorm, a blas64.Symmetric, work []float64) float64 { 580 return lapack64.Dlansy(norm, a.Uplo, a.N, a.Data, max(1, a.Stride), work) 581 } 582 583 // Lantr computes the specified norm of an m×n trapezoidal matrix A. If 584 // norm == lapack.MaxColumnSum work must have length at least n and this function 585 // will panic otherwise. There are no restrictions on work for the other matrix norms. 586 func Lantr(norm lapack.MatrixNorm, a blas64.Triangular, work []float64) float64 { 587 return lapack64.Dlantr(norm, a.Uplo, a.Diag, a.N, a.N, a.Data, max(1, a.Stride), work) 588 } 589 590 // Lantb computes the specified norm of an n×n triangular band matrix A. If 591 // norm == lapack.MaxColumnSum work must have length at least n and this function 592 // will panic otherwise. There are no restrictions on work for the other matrix 593 // norms. 594 func Lantb(norm lapack.MatrixNorm, a blas64.TriangularBand, work []float64) float64 { 595 return gonum.Implementation{}.Dlantb(norm, a.Uplo, a.Diag, a.N, a.K, a.Data, max(1, a.Stride), work) 596 } 597 598 // Lapmr rearranges the rows of the m×n matrix X as specified by the permutation 599 // k[0],k[1],...,k[m-1] of the integers 0,...,m-1. 600 // 601 // If forward is true, a forward permutation is applied: 602 // 603 // X[k[i],0:n] is moved to X[i,0:n] for i=0,1,...,m-1. 604 // 605 // If forward is false, a backward permutation is applied: 606 // 607 // X[i,0:n] is moved to X[k[i],0:n] for i=0,1,...,m-1. 608 // 609 // k must have length m, otherwise Lapmr will panic. 610 func Lapmr(forward bool, x blas64.General, k []int) { 611 lapack64.Dlapmr(forward, x.Rows, x.Cols, x.Data, max(1, x.Stride), k) 612 } 613 614 // Lapmt rearranges the columns of the m×n matrix X as specified by the 615 // permutation k_0, k_1, ..., k_{n-1} of the integers 0, ..., n-1. 616 // 617 // If forward is true a forward permutation is performed: 618 // 619 // X[0:m, k[j]] is moved to X[0:m, j] for j = 0, 1, ..., n-1. 620 // 621 // otherwise a backward permutation is performed: 622 // 623 // X[0:m, j] is moved to X[0:m, k[j]] for j = 0, 1, ..., n-1. 624 // 625 // k must have length n, otherwise Lapmt will panic. k is zero-indexed. 626 func Lapmt(forward bool, x blas64.General, k []int) { 627 lapack64.Dlapmt(forward, x.Rows, x.Cols, x.Data, max(1, x.Stride), k) 628 } 629 630 // Ormlq multiplies the matrix C by the othogonal matrix Q defined by 631 // A and tau. A and tau are as returned from Gelqf. 632 // 633 // C = Q * C if side == blas.Left and trans == blas.NoTrans 634 // C = Qᵀ * C if side == blas.Left and trans == blas.Trans 635 // C = C * Q if side == blas.Right and trans == blas.NoTrans 636 // C = C * Qᵀ if side == blas.Right and trans == blas.Trans 637 // 638 // If side == blas.Left, A is a matrix of side k×m, and if side == blas.Right 639 // A is of size k×n. This uses a blocked algorithm. 640 // 641 // Work is temporary storage, and lwork specifies the usable memory length. 642 // At minimum, lwork >= m if side == blas.Left and lwork >= n if side == blas.Right, 643 // and this function will panic otherwise. 644 // Ormlq uses a block algorithm, but the block size is limited 645 // by the temporary space available. If lwork == -1, instead of performing Ormlq, 646 // the optimal work length will be stored into work[0]. 647 // 648 // Tau contains the Householder scales and must have length at least k, and 649 // this function will panic otherwise. 650 func Ormlq(side blas.Side, trans blas.Transpose, a blas64.General, tau []float64, c blas64.General, work []float64, lwork int) { 651 lapack64.Dormlq(side, trans, c.Rows, c.Cols, a.Rows, a.Data, max(1, a.Stride), tau, c.Data, max(1, c.Stride), work, lwork) 652 } 653 654 // Ormqr multiplies an m×n matrix C by an orthogonal matrix Q as 655 // 656 // C = Q * C if side == blas.Left and trans == blas.NoTrans, 657 // C = Qᵀ * C if side == blas.Left and trans == blas.Trans, 658 // C = C * Q if side == blas.Right and trans == blas.NoTrans, 659 // C = C * Qᵀ if side == blas.Right and trans == blas.Trans, 660 // 661 // where Q is defined as the product of k elementary reflectors 662 // 663 // Q = H_0 * H_1 * ... * H_{k-1}. 664 // 665 // If side == blas.Left, A is an m×k matrix and 0 <= k <= m. 666 // If side == blas.Right, A is an n×k matrix and 0 <= k <= n. 667 // The ith column of A contains the vector which defines the elementary 668 // reflector H_i and tau[i] contains its scalar factor. tau must have length k 669 // and Ormqr will panic otherwise. Geqrf returns A and tau in the required 670 // form. 671 // 672 // work must have length at least max(1,lwork), and lwork must be at least n if 673 // side == blas.Left and at least m if side == blas.Right, otherwise Ormqr will 674 // panic. 675 // 676 // work is temporary storage, and lwork specifies the usable memory length. At 677 // minimum, lwork >= m if side == blas.Left and lwork >= n if side == 678 // blas.Right, and this function will panic otherwise. Larger values of lwork 679 // will generally give better performance. On return, work[0] will contain the 680 // optimal value of lwork. 681 // 682 // If lwork is -1, instead of performing Ormqr, the optimal workspace size will 683 // be stored into work[0]. 684 func Ormqr(side blas.Side, trans blas.Transpose, a blas64.General, tau []float64, c blas64.General, work []float64, lwork int) { 685 lapack64.Dormqr(side, trans, c.Rows, c.Cols, a.Cols, a.Data, max(1, a.Stride), tau, c.Data, max(1, c.Stride), work, lwork) 686 } 687 688 // Pocon estimates the reciprocal of the condition number of a positive-definite 689 // matrix A given the Cholesky decomposition of A. The condition number computed 690 // is based on the 1-norm and the ∞-norm. 691 // 692 // anorm is the 1-norm and the ∞-norm of the original matrix A. 693 // 694 // work is a temporary data slice of length at least 3*n and Pocon will panic otherwise. 695 // 696 // iwork is a temporary data slice of length at least n and Pocon will panic otherwise. 697 func Pocon(a blas64.Symmetric, anorm float64, work []float64, iwork []int) float64 { 698 return lapack64.Dpocon(a.Uplo, a.N, a.Data, max(1, a.Stride), anorm, work, iwork) 699 } 700 701 // Syev computes all eigenvalues and, optionally, the eigenvectors of a real 702 // symmetric matrix A. 703 // 704 // w contains the eigenvalues in ascending order upon return. w must have length 705 // at least n, and Syev will panic otherwise. 706 // 707 // On entry, a contains the elements of the symmetric matrix A in the triangular 708 // portion specified by uplo. If jobz == lapack.EVCompute, a contains the 709 // orthonormal eigenvectors of A on exit, otherwise jobz must be lapack.EVNone 710 // and on exit the specified triangular region is overwritten. 711 // 712 // Work is temporary storage, and lwork specifies the usable memory length. At minimum, 713 // lwork >= 3*n-1, and Syev will panic otherwise. The amount of blocking is 714 // limited by the usable length. If lwork == -1, instead of computing Syev the 715 // optimal work length is stored into work[0]. 716 func Syev(jobz lapack.EVJob, a blas64.Symmetric, w, work []float64, lwork int) (ok bool) { 717 return lapack64.Dsyev(jobz, a.Uplo, a.N, a.Data, max(1, a.Stride), w, work, lwork) 718 } 719 720 // Tbtrs solves a triangular system of the form 721 // 722 // A * X = B if trans == blas.NoTrans 723 // Aᵀ * X = B if trans == blas.Trans or blas.ConjTrans 724 // 725 // where A is an n×n triangular band matrix, and B is an n×nrhs matrix. 726 // 727 // Tbtrs returns whether A is non-singular. If A is singular, no solutions X 728 // are computed. 729 func Tbtrs(trans blas.Transpose, a blas64.TriangularBand, b blas64.General) (ok bool) { 730 return lapack64.Dtbtrs(a.Uplo, trans, a.Diag, a.N, a.K, b.Cols, a.Data, max(1, a.Stride), b.Data, max(1, b.Stride)) 731 } 732 733 // Trcon estimates the reciprocal of the condition number of a triangular matrix A. 734 // The condition number computed may be based on the 1-norm or the ∞-norm. 735 // 736 // work is a temporary data slice of length at least 3*n and Trcon will panic otherwise. 737 // 738 // iwork is a temporary data slice of length at least n and Trcon will panic otherwise. 739 func Trcon(norm lapack.MatrixNorm, a blas64.Triangular, work []float64, iwork []int) float64 { 740 return lapack64.Dtrcon(norm, a.Uplo, a.Diag, a.N, a.Data, max(1, a.Stride), work, iwork) 741 } 742 743 // Trtri computes the inverse of a triangular matrix, storing the result in place 744 // into a. 745 // 746 // Trtri will not perform the inversion if the matrix is singular, and returns 747 // a boolean indicating whether the inversion was successful. 748 func Trtri(a blas64.Triangular) (ok bool) { 749 return lapack64.Dtrtri(a.Uplo, a.Diag, a.N, a.Data, max(1, a.Stride)) 750 } 751 752 // Trtrs solves a triangular system of the form A * X = B or Aᵀ * X = B. Trtrs 753 // returns whether the solve completed successfully. If A is singular, no solve is performed. 754 func Trtrs(trans blas.Transpose, a blas64.Triangular, b blas64.General) (ok bool) { 755 return lapack64.Dtrtrs(a.Uplo, trans, a.Diag, a.N, b.Cols, a.Data, max(1, a.Stride), b.Data, max(1, b.Stride)) 756 } 757 758 // Geev computes the eigenvalues and, optionally, the left and/or right 759 // eigenvectors for an n×n real nonsymmetric matrix A. 760 // 761 // The right eigenvector v_j of A corresponding to an eigenvalue λ_j 762 // is defined by 763 // 764 // A v_j = λ_j v_j, 765 // 766 // and the left eigenvector u_j corresponding to an eigenvalue λ_j is defined by 767 // 768 // u_jᴴ A = λ_j u_jᴴ, 769 // 770 // where u_jᴴ is the conjugate transpose of u_j. 771 // 772 // On return, A will be overwritten and the left and right eigenvectors will be 773 // stored, respectively, in the columns of the n×n matrices VL and VR in the 774 // same order as their eigenvalues. If the j-th eigenvalue is real, then 775 // 776 // u_j = VL[:,j], 777 // v_j = VR[:,j], 778 // 779 // and if it is not real, then j and j+1 form a complex conjugate pair and the 780 // eigenvectors can be recovered as 781 // 782 // u_j = VL[:,j] + i*VL[:,j+1], 783 // u_{j+1} = VL[:,j] - i*VL[:,j+1], 784 // v_j = VR[:,j] + i*VR[:,j+1], 785 // v_{j+1} = VR[:,j] - i*VR[:,j+1], 786 // 787 // where i is the imaginary unit. The computed eigenvectors are normalized to 788 // have Euclidean norm equal to 1 and largest component real. 789 // 790 // Left eigenvectors will be computed only if jobvl == lapack.LeftEVCompute, 791 // otherwise jobvl must be lapack.LeftEVNone. 792 // Right eigenvectors will be computed only if jobvr == lapack.RightEVCompute, 793 // otherwise jobvr must be lapack.RightEVNone. 794 // For other values of jobvl and jobvr Geev will panic. 795 // 796 // On return, wr and wi will contain the real and imaginary parts, respectively, 797 // of the computed eigenvalues. Complex conjugate pairs of eigenvalues appear 798 // consecutively with the eigenvalue having the positive imaginary part first. 799 // wr and wi must have length n, and Geev will panic otherwise. 800 // 801 // work must have length at least lwork and lwork must be at least max(1,4*n) if 802 // the left or right eigenvectors are computed, and at least max(1,3*n) if no 803 // eigenvectors are computed. For good performance, lwork must generally be 804 // larger. On return, optimal value of lwork will be stored in work[0]. 805 // 806 // If lwork == -1, instead of performing Geev, the function only calculates the 807 // optimal value of lwork and stores it into work[0]. 808 // 809 // On return, first will be the index of the first valid eigenvalue. 810 // If first == 0, all eigenvalues and eigenvectors have been computed. 811 // If first is positive, Geev failed to compute all the eigenvalues, no 812 // eigenvectors have been computed and wr[first:] and wi[first:] contain those 813 // eigenvalues which have converged. 814 func Geev(jobvl lapack.LeftEVJob, jobvr lapack.RightEVJob, a blas64.General, wr, wi []float64, vl, vr blas64.General, work []float64, lwork int) (first int) { 815 n := a.Rows 816 if a.Cols != n { 817 panic("lapack64: matrix not square") 818 } 819 if jobvl == lapack.LeftEVCompute && (vl.Rows != n || vl.Cols != n) { 820 panic("lapack64: bad size of VL") 821 } 822 if jobvr == lapack.RightEVCompute && (vr.Rows != n || vr.Cols != n) { 823 panic("lapack64: bad size of VR") 824 } 825 return lapack64.Dgeev(jobvl, jobvr, n, a.Data, max(1, a.Stride), wr, wi, vl.Data, max(1, vl.Stride), vr.Data, max(1, vr.Stride), work, lwork) 826 }