gonum.org/v1/gonum@v0.14.0/lapack/lapack.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 lapack 6 7 import "gonum.org/v1/gonum/blas" 8 9 // Complex128 defines the public complex128 LAPACK API supported by gonum/lapack. 10 type Complex128 interface{} 11 12 // Float64 defines the public float64 LAPACK API supported by gonum/lapack. 13 type Float64 interface { 14 Dgecon(norm MatrixNorm, n int, a []float64, lda int, anorm float64, work []float64, iwork []int) float64 15 Dgeev(jobvl LeftEVJob, jobvr RightEVJob, n int, a []float64, lda int, wr, wi []float64, vl []float64, ldvl int, vr []float64, ldvr int, work []float64, lwork int) (first int) 16 Dgels(trans blas.Transpose, m, n, nrhs int, a []float64, lda int, b []float64, ldb int, work []float64, lwork int) bool 17 Dgelqf(m, n int, a []float64, lda int, tau, work []float64, lwork int) 18 Dgeqrf(m, n int, a []float64, lda int, tau, work []float64, lwork int) 19 Dgesvd(jobU, jobVT SVDJob, m, n int, a []float64, lda int, s, u []float64, ldu int, vt []float64, ldvt int, work []float64, lwork int) (ok bool) 20 Dgetrf(m, n int, a []float64, lda int, ipiv []int) (ok bool) 21 Dgetri(n int, a []float64, lda int, ipiv []int, work []float64, lwork int) (ok bool) 22 Dgetrs(trans blas.Transpose, n, nrhs int, a []float64, lda int, ipiv []int, b []float64, ldb int) 23 Dggsvd3(jobU, jobV, jobQ GSVDJob, m, n, p int, a []float64, lda int, b []float64, ldb int, alpha, beta, u []float64, ldu int, v []float64, ldv int, q []float64, ldq int, work []float64, lwork int, iwork []int) (k, l int, ok bool) 24 Dlantr(norm MatrixNorm, uplo blas.Uplo, diag blas.Diag, m, n int, a []float64, lda int, work []float64) float64 25 Dlange(norm MatrixNorm, m, n int, a []float64, lda int, work []float64) float64 26 Dlansy(norm MatrixNorm, uplo blas.Uplo, n int, a []float64, lda int, work []float64) float64 27 Dlapmr(forward bool, m, n int, x []float64, ldx int, k []int) 28 Dlapmt(forward bool, m, n int, x []float64, ldx int, k []int) 29 Dormqr(side blas.Side, trans blas.Transpose, m, n, k int, a []float64, lda int, tau, c []float64, ldc int, work []float64, lwork int) 30 Dormlq(side blas.Side, trans blas.Transpose, m, n, k int, a []float64, lda int, tau, c []float64, ldc int, work []float64, lwork int) 31 Dpbcon(uplo blas.Uplo, n, kd int, ab []float64, ldab int, anorm float64, work []float64, iwork []int) float64 32 Dpbtrf(uplo blas.Uplo, n, kd int, ab []float64, ldab int) (ok bool) 33 Dpbtrs(uplo blas.Uplo, n, kd, nrhs int, ab []float64, ldab int, b []float64, ldb int) 34 Dpocon(uplo blas.Uplo, n int, a []float64, lda int, anorm float64, work []float64, iwork []int) float64 35 Dpotrf(ul blas.Uplo, n int, a []float64, lda int) (ok bool) 36 Dpotri(ul blas.Uplo, n int, a []float64, lda int) (ok bool) 37 Dpotrs(ul blas.Uplo, n, nrhs int, a []float64, lda int, b []float64, ldb int) 38 Dpstrf(uplo blas.Uplo, n int, a []float64, lda int, piv []int, tol float64, work []float64) (rank int, ok bool) 39 Dsyev(jobz EVJob, uplo blas.Uplo, n int, a []float64, lda int, w, work []float64, lwork int) (ok bool) 40 Dtbtrs(uplo blas.Uplo, trans blas.Transpose, diag blas.Diag, n, kd, nrhs int, a []float64, lda int, b []float64, ldb int) (ok bool) 41 Dtrcon(norm MatrixNorm, uplo blas.Uplo, diag blas.Diag, n int, a []float64, lda int, work []float64, iwork []int) float64 42 Dtrtri(uplo blas.Uplo, diag blas.Diag, n int, a []float64, lda int) (ok bool) 43 Dtrtrs(uplo blas.Uplo, trans blas.Transpose, diag blas.Diag, n, nrhs int, a []float64, lda int, b []float64, ldb int) (ok bool) 44 } 45 46 // Direct specifies the direction of the multiplication for the Householder matrix. 47 type Direct byte 48 49 const ( 50 Forward Direct = 'F' // Reflectors are right-multiplied, H_0 * H_1 * ... * H_{k-1}. 51 Backward Direct = 'B' // Reflectors are left-multiplied, H_{k-1} * ... * H_1 * H_0. 52 ) 53 54 // Sort is the sorting order. 55 type Sort byte 56 57 const ( 58 SortIncreasing Sort = 'I' 59 SortDecreasing Sort = 'D' 60 ) 61 62 // StoreV indicates the storage direction of elementary reflectors. 63 type StoreV byte 64 65 const ( 66 ColumnWise StoreV = 'C' // Reflector stored in a column of the matrix. 67 RowWise StoreV = 'R' // Reflector stored in a row of the matrix. 68 ) 69 70 // MatrixNorm represents the kind of matrix norm to compute. 71 type MatrixNorm byte 72 73 const ( 74 MaxAbs MatrixNorm = 'M' // max(abs(A(i,j))) 75 MaxColumnSum MatrixNorm = 'O' // Maximum absolute column sum (one norm) 76 MaxRowSum MatrixNorm = 'I' // Maximum absolute row sum (infinity norm) 77 Frobenius MatrixNorm = 'F' // Frobenius norm (sqrt of sum of squares) 78 ) 79 80 // MatrixType represents the kind of matrix represented in the data. 81 type MatrixType byte 82 83 const ( 84 General MatrixType = 'G' // A general dense matrix. 85 UpperTri MatrixType = 'U' // An upper triangular matrix. 86 LowerTri MatrixType = 'L' // A lower triangular matrix. 87 ) 88 89 // Pivot specifies the pivot type for plane rotations. 90 type Pivot byte 91 92 const ( 93 Variable Pivot = 'V' 94 Top Pivot = 'T' 95 Bottom Pivot = 'B' 96 ) 97 98 // ApplyOrtho specifies which orthogonal matrix is applied in Dormbr. 99 type ApplyOrtho byte 100 101 const ( 102 ApplyP ApplyOrtho = 'P' // Apply P or Pᵀ. 103 ApplyQ ApplyOrtho = 'Q' // Apply Q or Qᵀ. 104 ) 105 106 // GenOrtho specifies which orthogonal matrix is generated in Dorgbr. 107 type GenOrtho byte 108 109 const ( 110 GeneratePT GenOrtho = 'P' // Generate Pᵀ. 111 GenerateQ GenOrtho = 'Q' // Generate Q. 112 ) 113 114 // SVDJob specifies the singular vector computation type for SVD. 115 type SVDJob byte 116 117 const ( 118 SVDAll SVDJob = 'A' // Compute all columns of the orthogonal matrix U or V. 119 SVDStore SVDJob = 'S' // Compute the singular vectors and store them in the orthogonal matrix U or V. 120 SVDOverwrite SVDJob = 'O' // Compute the singular vectors and overwrite them on the input matrix A. 121 SVDNone SVDJob = 'N' // Do not compute singular vectors. 122 ) 123 124 // GSVDJob specifies the singular vector computation type for Generalized SVD. 125 type GSVDJob byte 126 127 const ( 128 GSVDU GSVDJob = 'U' // Compute orthogonal matrix U. 129 GSVDV GSVDJob = 'V' // Compute orthogonal matrix V. 130 GSVDQ GSVDJob = 'Q' // Compute orthogonal matrix Q. 131 GSVDUnit GSVDJob = 'I' // Use unit-initialized matrix. 132 GSVDNone GSVDJob = 'N' // Do not compute orthogonal matrix. 133 ) 134 135 // EVComp specifies how eigenvectors are computed in Dsteqr. 136 type EVComp byte 137 138 const ( 139 EVOrig EVComp = 'V' // Compute eigenvectors of the original symmetric matrix. 140 EVTridiag EVComp = 'I' // Compute eigenvectors of the tridiagonal matrix. 141 EVCompNone EVComp = 'N' // Do not compute eigenvectors. 142 ) 143 144 // EVJob specifies whether eigenvectors are computed in Dsyev. 145 type EVJob byte 146 147 const ( 148 EVCompute EVJob = 'V' // Compute eigenvectors. 149 EVNone EVJob = 'N' // Do not compute eigenvectors. 150 ) 151 152 // LeftEVJob specifies whether left eigenvectors are computed in Dgeev. 153 type LeftEVJob byte 154 155 const ( 156 LeftEVCompute LeftEVJob = 'V' // Compute left eigenvectors. 157 LeftEVNone LeftEVJob = 'N' // Do not compute left eigenvectors. 158 ) 159 160 // RightEVJob specifies whether right eigenvectors are computed in Dgeev. 161 type RightEVJob byte 162 163 const ( 164 RightEVCompute RightEVJob = 'V' // Compute right eigenvectors. 165 RightEVNone RightEVJob = 'N' // Do not compute right eigenvectors. 166 ) 167 168 // BalanceJob specifies matrix balancing operation. 169 type BalanceJob byte 170 171 const ( 172 Permute BalanceJob = 'P' 173 Scale BalanceJob = 'S' 174 PermuteScale BalanceJob = 'B' 175 BalanceNone BalanceJob = 'N' 176 ) 177 178 // SchurJob specifies whether the Schur form is computed in Dhseqr. 179 type SchurJob byte 180 181 const ( 182 EigenvaluesOnly SchurJob = 'E' 183 EigenvaluesAndSchur SchurJob = 'S' 184 ) 185 186 // SchurComp specifies whether and how the Schur vectors are computed in Dhseqr. 187 type SchurComp byte 188 189 const ( 190 SchurOrig SchurComp = 'V' // Compute Schur vectors of the original matrix. 191 SchurHess SchurComp = 'I' // Compute Schur vectors of the upper Hessenberg matrix. 192 SchurNone SchurComp = 'N' // Do not compute Schur vectors. 193 ) 194 195 // UpdateSchurComp specifies whether the matrix of Schur vectors is updated in Dtrexc. 196 type UpdateSchurComp byte 197 198 const ( 199 UpdateSchur UpdateSchurComp = 'V' // Update the matrix of Schur vectors. 200 UpdateSchurNone UpdateSchurComp = 'N' // Do not update the matrix of Schur vectors. 201 ) 202 203 // EVSide specifies what eigenvectors are computed in Dtrevc3. 204 type EVSide byte 205 206 const ( 207 EVRight EVSide = 'R' // Compute only right eigenvectors. 208 EVLeft EVSide = 'L' // Compute only left eigenvectors. 209 EVBoth EVSide = 'B' // Compute both right and left eigenvectors. 210 ) 211 212 // EVHowMany specifies which eigenvectors are computed in Dtrevc3 and how. 213 type EVHowMany byte 214 215 const ( 216 EVAll EVHowMany = 'A' // Compute all right and/or left eigenvectors. 217 EVAllMulQ EVHowMany = 'B' // Compute all right and/or left eigenvectors multiplied by an input matrix. 218 EVSelected EVHowMany = 'S' // Compute selected right and/or left eigenvectors. 219 ) 220 221 // MaximizeNormX specifies the heuristic method for computing a contribution to 222 // the reciprocal Dif-estimate in Dlatdf. 223 type MaximizeNormXJob byte 224 225 const ( 226 LocalLookAhead MaximizeNormXJob = 0 // Solve Z*x=h-f where h is a vector of ±1. 227 NormalizedNullVector MaximizeNormXJob = 2 // Compute an approximate null-vector e of Z, normalize e and solve Z*x=±e-f. 228 )