lapack64.go 28 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613
  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. package lapack64
  5. import (
  6. "gonum.org/v1/gonum/blas"
  7. "gonum.org/v1/gonum/blas/blas64"
  8. "gonum.org/v1/gonum/lapack"
  9. "gonum.org/v1/gonum/lapack/gonum"
  10. )
  11. var lapack64 lapack.Float64 = gonum.Implementation{}
  12. // Use sets the LAPACK float64 implementation to be used by subsequent BLAS calls.
  13. // The default implementation is native.Implementation.
  14. func Use(l lapack.Float64) {
  15. lapack64 = l
  16. }
  17. func max(a, b int) int {
  18. if a > b {
  19. return a
  20. }
  21. return b
  22. }
  23. // Potrf computes the Cholesky factorization of a.
  24. // The factorization has the form
  25. // A = Uᵀ * U if a.Uplo == blas.Upper, or
  26. // A = L * Lᵀ if a.Uplo == blas.Lower,
  27. // where U is an upper triangular matrix and L is lower triangular.
  28. // The triangular matrix is returned in t, and the underlying data between
  29. // a and t is shared. The returned bool indicates whether a is positive
  30. // definite and the factorization could be finished.
  31. func Potrf(a blas64.Symmetric) (t blas64.Triangular, ok bool) {
  32. ok = lapack64.Dpotrf(a.Uplo, a.N, a.Data, max(1, a.Stride))
  33. t.Uplo = a.Uplo
  34. t.N = a.N
  35. t.Data = a.Data
  36. t.Stride = a.Stride
  37. t.Diag = blas.NonUnit
  38. return
  39. }
  40. // Potri computes the inverse of a real symmetric positive definite matrix A
  41. // using its Cholesky factorization.
  42. //
  43. // On entry, t contains the triangular factor U or L from the Cholesky
  44. // factorization A = Uᵀ*U or A = L*Lᵀ, as computed by Potrf.
  45. //
  46. // On return, the upper or lower triangle of the (symmetric) inverse of A is
  47. // stored in t, overwriting the input factor U or L, and also returned in a. The
  48. // underlying data between a and t is shared.
  49. //
  50. // The returned bool indicates whether the inverse was computed successfully.
  51. func Potri(t blas64.Triangular) (a blas64.Symmetric, ok bool) {
  52. ok = lapack64.Dpotri(t.Uplo, t.N, t.Data, max(1, t.Stride))
  53. a.Uplo = t.Uplo
  54. a.N = t.N
  55. a.Data = t.Data
  56. a.Stride = t.Stride
  57. return
  58. }
  59. // Potrs solves a system of n linear equations A*X = B where A is an n×n
  60. // symmetric positive definite matrix and B is an n×nrhs matrix, using the
  61. // Cholesky factorization A = Uᵀ*U or A = L*Lᵀ. t contains the corresponding
  62. // triangular factor as returned by Potrf. On entry, B contains the right-hand
  63. // side matrix B, on return it contains the solution matrix X.
  64. func Potrs(t blas64.Triangular, b blas64.General) {
  65. lapack64.Dpotrs(t.Uplo, t.N, b.Cols, t.Data, max(1, t.Stride), b.Data, max(1, b.Stride))
  66. }
  67. // Pbtrf computes the Cholesky factorization of an n×n symmetric positive
  68. // definite band matrix
  69. // A = Uᵀ * U if a.Uplo == blas.Upper
  70. // A = L * Lᵀ if a.Uplo == blas.Lower
  71. // where U and L are upper, respectively lower, triangular band matrices.
  72. //
  73. // The triangular matrix U or L is returned in t, and the underlying data
  74. // between a and t is shared. The returned bool indicates whether A is positive
  75. // definite and the factorization could be finished.
  76. func Pbtrf(a blas64.SymmetricBand) (t blas64.TriangularBand, ok bool) {
  77. ok = lapack64.Dpbtrf(a.Uplo, a.N, a.K, a.Data, max(1, a.Stride))
  78. t.Uplo = a.Uplo
  79. t.Diag = blas.NonUnit
  80. t.N = a.N
  81. t.K = a.K
  82. t.Data = a.Data
  83. t.Stride = a.Stride
  84. return t, ok
  85. }
  86. // Pbtrs solves a system of linear equations A*X = B with an n×n symmetric
  87. // positive definite band matrix A using the Cholesky factorization
  88. // A = Uᵀ * U if t.Uplo == blas.Upper
  89. // A = L * Lᵀ if t.Uplo == blas.Lower
  90. // t contains the corresponding triangular factor as returned by Pbtrf.
  91. //
  92. // On entry, b contains the right hand side matrix B. On return, it is
  93. // overwritten with the solution matrix X.
  94. func Pbtrs(t blas64.TriangularBand, b blas64.General) {
  95. lapack64.Dpbtrs(t.Uplo, t.N, t.K, b.Cols, t.Data, max(1, t.Stride), b.Data, max(1, b.Stride))
  96. }
  97. // Gecon estimates the reciprocal of the condition number of the n×n matrix A
  98. // given the LU decomposition of the matrix. The condition number computed may
  99. // be based on the 1-norm or the ∞-norm.
  100. //
  101. // a contains the result of the LU decomposition of A as computed by Getrf.
  102. //
  103. // anorm is the corresponding 1-norm or ∞-norm of the original matrix A.
  104. //
  105. // work is a temporary data slice of length at least 4*n and Gecon will panic otherwise.
  106. //
  107. // iwork is a temporary data slice of length at least n and Gecon will panic otherwise.
  108. func Gecon(norm lapack.MatrixNorm, a blas64.General, anorm float64, work []float64, iwork []int) float64 {
  109. return lapack64.Dgecon(norm, a.Cols, a.Data, max(1, a.Stride), anorm, work, iwork)
  110. }
  111. // Gels finds a minimum-norm solution based on the matrices A and B using the
  112. // QR or LQ factorization. Gels returns false if the matrix
  113. // A is singular, and true if this solution was successfully found.
  114. //
  115. // The minimization problem solved depends on the input parameters.
  116. //
  117. // 1. If m >= n and trans == blas.NoTrans, Gels finds X such that || A*X - B||_2
  118. // is minimized.
  119. // 2. If m < n and trans == blas.NoTrans, Gels finds the minimum norm solution of
  120. // A * X = B.
  121. // 3. If m >= n and trans == blas.Trans, Gels finds the minimum norm solution of
  122. // Aᵀ * X = B.
  123. // 4. If m < n and trans == blas.Trans, Gels finds X such that || A*X - B||_2
  124. // is minimized.
  125. // Note that the least-squares solutions (cases 1 and 3) perform the minimization
  126. // per column of B. This is not the same as finding the minimum-norm matrix.
  127. //
  128. // The matrix A is a general matrix of size m×n and is modified during this call.
  129. // The input matrix B is of size max(m,n)×nrhs, and serves two purposes. On entry,
  130. // the elements of b specify the input matrix B. B has size m×nrhs if
  131. // trans == blas.NoTrans, and n×nrhs if trans == blas.Trans. On exit, the
  132. // leading submatrix of b contains the solution vectors X. If trans == blas.NoTrans,
  133. // this submatrix is of size n×nrhs, and of size m×nrhs otherwise.
  134. //
  135. // Work is temporary storage, and lwork specifies the usable memory length.
  136. // At minimum, lwork >= max(m,n) + max(m,n,nrhs), and this function will panic
  137. // otherwise. A longer work will enable blocked algorithms to be called.
  138. // In the special case that lwork == -1, work[0] will be set to the optimal working
  139. // length.
  140. func Gels(trans blas.Transpose, a blas64.General, b blas64.General, work []float64, lwork int) bool {
  141. return lapack64.Dgels(trans, a.Rows, a.Cols, b.Cols, a.Data, max(1, a.Stride), b.Data, max(1, b.Stride), work, lwork)
  142. }
  143. // Geqrf computes the QR factorization of the m×n matrix A using a blocked
  144. // algorithm. A is modified to contain the information to construct Q and R.
  145. // The upper triangle of a contains the matrix R. The lower triangular elements
  146. // (not including the diagonal) contain the elementary reflectors. tau is modified
  147. // to contain the reflector scales. tau must have length at least min(m,n), and
  148. // this function will panic otherwise.
  149. //
  150. // The ith elementary reflector can be explicitly constructed by first extracting
  151. // the
  152. // v[j] = 0 j < i
  153. // v[j] = 1 j == i
  154. // v[j] = a[j*lda+i] j > i
  155. // and computing H_i = I - tau[i] * v * vᵀ.
  156. //
  157. // The orthonormal matrix Q can be constucted from a product of these elementary
  158. // reflectors, Q = H_0 * H_1 * ... * H_{k-1}, where k = min(m,n).
  159. //
  160. // Work is temporary storage, and lwork specifies the usable memory length.
  161. // At minimum, lwork >= m and this function will panic otherwise.
  162. // Geqrf is a blocked QR factorization, but the block size is limited
  163. // by the temporary space available. If lwork == -1, instead of performing Geqrf,
  164. // the optimal work length will be stored into work[0].
  165. func Geqrf(a blas64.General, tau, work []float64, lwork int) {
  166. lapack64.Dgeqrf(a.Rows, a.Cols, a.Data, max(1, a.Stride), tau, work, lwork)
  167. }
  168. // Gelqf computes the LQ factorization of the m×n matrix A using a blocked
  169. // algorithm. A is modified to contain the information to construct L and Q. The
  170. // lower triangle of a contains the matrix L. The elements above the diagonal
  171. // and the slice tau represent the matrix Q. tau is modified to contain the
  172. // reflector scales. tau must have length at least min(m,n), and this function
  173. // will panic otherwise.
  174. //
  175. // See Geqrf for a description of the elementary reflectors and orthonormal
  176. // matrix Q. Q is constructed as a product of these elementary reflectors,
  177. // Q = H_{k-1} * ... * H_1 * H_0.
  178. //
  179. // Work is temporary storage, and lwork specifies the usable memory length.
  180. // At minimum, lwork >= m and this function will panic otherwise.
  181. // Gelqf is a blocked LQ factorization, but the block size is limited
  182. // by the temporary space available. If lwork == -1, instead of performing Gelqf,
  183. // the optimal work length will be stored into work[0].
  184. func Gelqf(a blas64.General, tau, work []float64, lwork int) {
  185. lapack64.Dgelqf(a.Rows, a.Cols, a.Data, max(1, a.Stride), tau, work, lwork)
  186. }
  187. // Gesvd computes the singular value decomposition of the input matrix A.
  188. //
  189. // The singular value decomposition is
  190. // A = U * Sigma * Vᵀ
  191. // where Sigma is an m×n diagonal matrix containing the singular values of A,
  192. // U is an m×m orthogonal matrix and V is an n×n orthogonal matrix. The first
  193. // min(m,n) columns of U and V are the left and right singular vectors of A
  194. // respectively.
  195. //
  196. // jobU and jobVT are options for computing the singular vectors. The behavior
  197. // is as follows
  198. // jobU == lapack.SVDAll All m columns of U are returned in u
  199. // jobU == lapack.SVDStore The first min(m,n) columns are returned in u
  200. // jobU == lapack.SVDOverwrite The first min(m,n) columns of U are written into a
  201. // jobU == lapack.SVDNone The columns of U are not computed.
  202. // The behavior is the same for jobVT and the rows of Vᵀ. At most one of jobU
  203. // and jobVT can equal lapack.SVDOverwrite, and Gesvd will panic otherwise.
  204. //
  205. // On entry, a contains the data for the m×n matrix A. During the call to Gesvd
  206. // the data is overwritten. On exit, A contains the appropriate singular vectors
  207. // if either job is lapack.SVDOverwrite.
  208. //
  209. // s is a slice of length at least min(m,n) and on exit contains the singular
  210. // values in decreasing order.
  211. //
  212. // u contains the left singular vectors on exit, stored columnwise. If
  213. // jobU == lapack.SVDAll, u is of size m×m. If jobU == lapack.SVDStore u is
  214. // of size m×min(m,n). If jobU == lapack.SVDOverwrite or lapack.SVDNone, u is
  215. // not used.
  216. //
  217. // vt contains the left singular vectors on exit, stored rowwise. If
  218. // jobV == lapack.SVDAll, vt is of size n×m. If jobVT == lapack.SVDStore vt is
  219. // of size min(m,n)×n. If jobVT == lapack.SVDOverwrite or lapack.SVDNone, vt is
  220. // not used.
  221. //
  222. // work is a slice for storing temporary memory, and lwork is the usable size of
  223. // the slice. lwork must be at least max(5*min(m,n), 3*min(m,n)+max(m,n)).
  224. // If lwork == -1, instead of performing Gesvd, the optimal work length will be
  225. // stored into work[0]. Gesvd will panic if the working memory has insufficient
  226. // storage.
  227. //
  228. // Gesvd returns whether the decomposition successfully completed.
  229. func Gesvd(jobU, jobVT lapack.SVDJob, a, u, vt blas64.General, s, work []float64, lwork int) (ok bool) {
  230. 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)
  231. }
  232. // Getrf computes the LU decomposition of the m×n matrix A.
  233. // The LU decomposition is a factorization of A into
  234. // A = P * L * U
  235. // where P is a permutation matrix, L is a unit lower triangular matrix, and
  236. // U is a (usually) non-unit upper triangular matrix. On exit, L and U are stored
  237. // in place into a.
  238. //
  239. // ipiv is a permutation vector. It indicates that row i of the matrix was
  240. // changed with ipiv[i]. ipiv must have length at least min(m,n), and will panic
  241. // otherwise. ipiv is zero-indexed.
  242. //
  243. // Getrf is the blocked version of the algorithm.
  244. //
  245. // Getrf returns whether the matrix A is singular. The LU decomposition will
  246. // be computed regardless of the singularity of A, but division by zero
  247. // will occur if the false is returned and the result is used to solve a
  248. // system of equations.
  249. func Getrf(a blas64.General, ipiv []int) bool {
  250. return lapack64.Dgetrf(a.Rows, a.Cols, a.Data, max(1, a.Stride), ipiv)
  251. }
  252. // Getri computes the inverse of the matrix A using the LU factorization computed
  253. // by Getrf. On entry, a contains the PLU decomposition of A as computed by
  254. // Getrf and on exit contains the reciprocal of the original matrix.
  255. //
  256. // Getri will not perform the inversion if the matrix is singular, and returns
  257. // a boolean indicating whether the inversion was successful.
  258. //
  259. // Work is temporary storage, and lwork specifies the usable memory length.
  260. // At minimum, lwork >= n and this function will panic otherwise.
  261. // Getri is a blocked inversion, but the block size is limited
  262. // by the temporary space available. If lwork == -1, instead of performing Getri,
  263. // the optimal work length will be stored into work[0].
  264. func Getri(a blas64.General, ipiv []int, work []float64, lwork int) (ok bool) {
  265. return lapack64.Dgetri(a.Cols, a.Data, max(1, a.Stride), ipiv, work, lwork)
  266. }
  267. // Getrs solves a system of equations using an LU factorization.
  268. // The system of equations solved is
  269. // A * X = B if trans == blas.Trans
  270. // Aᵀ * X = B if trans == blas.NoTrans
  271. // A is a general n×n matrix with stride lda. B is a general matrix of size n×nrhs.
  272. //
  273. // On entry b contains the elements of the matrix B. On exit, b contains the
  274. // elements of X, the solution to the system of equations.
  275. //
  276. // a and ipiv contain the LU factorization of A and the permutation indices as
  277. // computed by Getrf. ipiv is zero-indexed.
  278. func Getrs(trans blas.Transpose, a blas64.General, b blas64.General, ipiv []int) {
  279. lapack64.Dgetrs(trans, a.Cols, b.Cols, a.Data, max(1, a.Stride), ipiv, b.Data, max(1, b.Stride))
  280. }
  281. // Ggsvd3 computes the generalized singular value decomposition (GSVD)
  282. // of an m×n matrix A and p×n matrix B:
  283. // Uᵀ*A*Q = D1*[ 0 R ]
  284. //
  285. // Vᵀ*B*Q = D2*[ 0 R ]
  286. // where U, V and Q are orthogonal matrices.
  287. //
  288. // Ggsvd3 returns k and l, the dimensions of the sub-blocks. k+l
  289. // is the effective numerical rank of the (m+p)×n matrix [ Aᵀ Bᵀ ]ᵀ.
  290. // R is a (k+l)×(k+l) nonsingular upper triangular matrix, D1 and
  291. // D2 are m×(k+l) and p×(k+l) diagonal matrices and of the following
  292. // structures, respectively:
  293. //
  294. // If m-k-l >= 0,
  295. //
  296. // k l
  297. // D1 = k [ I 0 ]
  298. // l [ 0 C ]
  299. // m-k-l [ 0 0 ]
  300. //
  301. // k l
  302. // D2 = l [ 0 S ]
  303. // p-l [ 0 0 ]
  304. //
  305. // n-k-l k l
  306. // [ 0 R ] = k [ 0 R11 R12 ] k
  307. // l [ 0 0 R22 ] l
  308. //
  309. // where
  310. //
  311. // C = diag( alpha_k, ... , alpha_{k+l} ),
  312. // S = diag( beta_k, ... , beta_{k+l} ),
  313. // C^2 + S^2 = I.
  314. //
  315. // R is stored in
  316. // A[0:k+l, n-k-l:n]
  317. // on exit.
  318. //
  319. // If m-k-l < 0,
  320. //
  321. // k m-k k+l-m
  322. // D1 = k [ I 0 0 ]
  323. // m-k [ 0 C 0 ]
  324. //
  325. // k m-k k+l-m
  326. // D2 = m-k [ 0 S 0 ]
  327. // k+l-m [ 0 0 I ]
  328. // p-l [ 0 0 0 ]
  329. //
  330. // n-k-l k m-k k+l-m
  331. // [ 0 R ] = k [ 0 R11 R12 R13 ]
  332. // m-k [ 0 0 R22 R23 ]
  333. // k+l-m [ 0 0 0 R33 ]
  334. //
  335. // where
  336. // C = diag( alpha_k, ... , alpha_m ),
  337. // S = diag( beta_k, ... , beta_m ),
  338. // C^2 + S^2 = I.
  339. //
  340. // R = [ R11 R12 R13 ] is stored in A[1:m, n-k-l+1:n]
  341. // [ 0 R22 R23 ]
  342. // and R33 is stored in
  343. // B[m-k:l, n+m-k-l:n] on exit.
  344. //
  345. // Ggsvd3 computes C, S, R, and optionally the orthogonal transformation
  346. // matrices U, V and Q.
  347. //
  348. // jobU, jobV and jobQ are options for computing the orthogonal matrices. The behavior
  349. // is as follows
  350. // jobU == lapack.GSVDU Compute orthogonal matrix U
  351. // jobU == lapack.GSVDNone Do not compute orthogonal matrix.
  352. // The behavior is the same for jobV and jobQ with the exception that instead of
  353. // lapack.GSVDU these accept lapack.GSVDV and lapack.GSVDQ respectively.
  354. // The matrices U, V and Q must be m×m, p×p and n×n respectively unless the
  355. // relevant job parameter is lapack.GSVDNone.
  356. //
  357. // alpha and beta must have length n or Ggsvd3 will panic. On exit, alpha and
  358. // beta contain the generalized singular value pairs of A and B
  359. // alpha[0:k] = 1,
  360. // beta[0:k] = 0,
  361. // if m-k-l >= 0,
  362. // alpha[k:k+l] = diag(C),
  363. // beta[k:k+l] = diag(S),
  364. // if m-k-l < 0,
  365. // alpha[k:m]= C, alpha[m:k+l]= 0
  366. // beta[k:m] = S, beta[m:k+l] = 1.
  367. // if k+l < n,
  368. // alpha[k+l:n] = 0 and
  369. // beta[k+l:n] = 0.
  370. //
  371. // On exit, iwork contains the permutation required to sort alpha descending.
  372. //
  373. // iwork must have length n, work must have length at least max(1, lwork), and
  374. // lwork must be -1 or greater than n, otherwise Ggsvd3 will panic. If
  375. // lwork is -1, work[0] holds the optimal lwork on return, but Ggsvd3 does
  376. // not perform the GSVD.
  377. 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) {
  378. 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)
  379. }
  380. // Lange computes the matrix norm of the general m×n matrix A. The input norm
  381. // specifies the norm computed.
  382. // lapack.MaxAbs: the maximum absolute value of an element.
  383. // lapack.MaxColumnSum: the maximum column sum of the absolute values of the entries.
  384. // lapack.MaxRowSum: the maximum row sum of the absolute values of the entries.
  385. // lapack.Frobenius: the square root of the sum of the squares of the entries.
  386. // If norm == lapack.MaxColumnSum, work must be of length n, and this function will panic otherwise.
  387. // There are no restrictions on work for the other matrix norms.
  388. func Lange(norm lapack.MatrixNorm, a blas64.General, work []float64) float64 {
  389. return lapack64.Dlange(norm, a.Rows, a.Cols, a.Data, max(1, a.Stride), work)
  390. }
  391. // Lansy computes the specified norm of an n×n symmetric matrix. If
  392. // norm == lapack.MaxColumnSum or norm == lapack.MaxRowSum, work must have length
  393. // at least n and this function will panic otherwise.
  394. // There are no restrictions on work for the other matrix norms.
  395. func Lansy(norm lapack.MatrixNorm, a blas64.Symmetric, work []float64) float64 {
  396. return lapack64.Dlansy(norm, a.Uplo, a.N, a.Data, max(1, a.Stride), work)
  397. }
  398. // Lantr computes the specified norm of an m×n trapezoidal matrix A. If
  399. // norm == lapack.MaxColumnSum work must have length at least n and this function
  400. // will panic otherwise. There are no restrictions on work for the other matrix norms.
  401. func Lantr(norm lapack.MatrixNorm, a blas64.Triangular, work []float64) float64 {
  402. return lapack64.Dlantr(norm, a.Uplo, a.Diag, a.N, a.N, a.Data, max(1, a.Stride), work)
  403. }
  404. // Lapmt rearranges the columns of the m×n matrix X as specified by the
  405. // permutation k_0, k_1, ..., k_{n-1} of the integers 0, ..., n-1.
  406. //
  407. // If forward is true a forward permutation is performed:
  408. //
  409. // X[0:m, k[j]] is moved to X[0:m, j] for j = 0, 1, ..., n-1.
  410. //
  411. // otherwise a backward permutation is performed:
  412. //
  413. // X[0:m, j] is moved to X[0:m, k[j]] for j = 0, 1, ..., n-1.
  414. //
  415. // k must have length n, otherwise Lapmt will panic. k is zero-indexed.
  416. func Lapmt(forward bool, x blas64.General, k []int) {
  417. lapack64.Dlapmt(forward, x.Rows, x.Cols, x.Data, max(1, x.Stride), k)
  418. }
  419. // Ormlq multiplies the matrix C by the othogonal matrix Q defined by
  420. // A and tau. A and tau are as returned from Gelqf.
  421. // C = Q * C if side == blas.Left and trans == blas.NoTrans
  422. // C = Qᵀ * C if side == blas.Left and trans == blas.Trans
  423. // C = C * Q if side == blas.Right and trans == blas.NoTrans
  424. // C = C * Qᵀ if side == blas.Right and trans == blas.Trans
  425. // If side == blas.Left, A is a matrix of side k×m, and if side == blas.Right
  426. // A is of size k×n. This uses a blocked algorithm.
  427. //
  428. // Work is temporary storage, and lwork specifies the usable memory length.
  429. // At minimum, lwork >= m if side == blas.Left and lwork >= n if side == blas.Right,
  430. // and this function will panic otherwise.
  431. // Ormlq uses a block algorithm, but the block size is limited
  432. // by the temporary space available. If lwork == -1, instead of performing Ormlq,
  433. // the optimal work length will be stored into work[0].
  434. //
  435. // Tau contains the Householder scales and must have length at least k, and
  436. // this function will panic otherwise.
  437. func Ormlq(side blas.Side, trans blas.Transpose, a blas64.General, tau []float64, c blas64.General, work []float64, lwork int) {
  438. 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)
  439. }
  440. // Ormqr multiplies an m×n matrix C by an orthogonal matrix Q as
  441. // C = Q * C if side == blas.Left and trans == blas.NoTrans,
  442. // C = Qᵀ * C if side == blas.Left and trans == blas.Trans,
  443. // C = C * Q if side == blas.Right and trans == blas.NoTrans,
  444. // C = C * Qᵀ if side == blas.Right and trans == blas.Trans,
  445. // where Q is defined as the product of k elementary reflectors
  446. // Q = H_0 * H_1 * ... * H_{k-1}.
  447. //
  448. // If side == blas.Left, A is an m×k matrix and 0 <= k <= m.
  449. // If side == blas.Right, A is an n×k matrix and 0 <= k <= n.
  450. // The ith column of A contains the vector which defines the elementary
  451. // reflector H_i and tau[i] contains its scalar factor. tau must have length k
  452. // and Ormqr will panic otherwise. Geqrf returns A and tau in the required
  453. // form.
  454. //
  455. // work must have length at least max(1,lwork), and lwork must be at least n if
  456. // side == blas.Left and at least m if side == blas.Right, otherwise Ormqr will
  457. // panic.
  458. //
  459. // work is temporary storage, and lwork specifies the usable memory length. At
  460. // minimum, lwork >= m if side == blas.Left and lwork >= n if side ==
  461. // blas.Right, and this function will panic otherwise. Larger values of lwork
  462. // will generally give better performance. On return, work[0] will contain the
  463. // optimal value of lwork.
  464. //
  465. // If lwork is -1, instead of performing Ormqr, the optimal workspace size will
  466. // be stored into work[0].
  467. func Ormqr(side blas.Side, trans blas.Transpose, a blas64.General, tau []float64, c blas64.General, work []float64, lwork int) {
  468. 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)
  469. }
  470. // Pocon estimates the reciprocal of the condition number of a positive-definite
  471. // matrix A given the Cholesky decmposition of A. The condition number computed
  472. // is based on the 1-norm and the ∞-norm.
  473. //
  474. // anorm is the 1-norm and the ∞-norm of the original matrix A.
  475. //
  476. // work is a temporary data slice of length at least 3*n and Pocon will panic otherwise.
  477. //
  478. // iwork is a temporary data slice of length at least n and Pocon will panic otherwise.
  479. func Pocon(a blas64.Symmetric, anorm float64, work []float64, iwork []int) float64 {
  480. return lapack64.Dpocon(a.Uplo, a.N, a.Data, max(1, a.Stride), anorm, work, iwork)
  481. }
  482. // Syev computes all eigenvalues and, optionally, the eigenvectors of a real
  483. // symmetric matrix A.
  484. //
  485. // w contains the eigenvalues in ascending order upon return. w must have length
  486. // at least n, and Syev will panic otherwise.
  487. //
  488. // On entry, a contains the elements of the symmetric matrix A in the triangular
  489. // portion specified by uplo. If jobz == lapack.EVCompute, a contains the
  490. // orthonormal eigenvectors of A on exit, otherwise jobz must be lapack.EVNone
  491. // and on exit the specified triangular region is overwritten.
  492. //
  493. // Work is temporary storage, and lwork specifies the usable memory length. At minimum,
  494. // lwork >= 3*n-1, and Syev will panic otherwise. The amount of blocking is
  495. // limited by the usable length. If lwork == -1, instead of computing Syev the
  496. // optimal work length is stored into work[0].
  497. func Syev(jobz lapack.EVJob, a blas64.Symmetric, w, work []float64, lwork int) (ok bool) {
  498. return lapack64.Dsyev(jobz, a.Uplo, a.N, a.Data, max(1, a.Stride), w, work, lwork)
  499. }
  500. // Trcon estimates the reciprocal of the condition number of a triangular matrix A.
  501. // The condition number computed may be based on the 1-norm or the ∞-norm.
  502. //
  503. // work is a temporary data slice of length at least 3*n and Trcon will panic otherwise.
  504. //
  505. // iwork is a temporary data slice of length at least n and Trcon will panic otherwise.
  506. func Trcon(norm lapack.MatrixNorm, a blas64.Triangular, work []float64, iwork []int) float64 {
  507. return lapack64.Dtrcon(norm, a.Uplo, a.Diag, a.N, a.Data, max(1, a.Stride), work, iwork)
  508. }
  509. // Trtri computes the inverse of a triangular matrix, storing the result in place
  510. // into a.
  511. //
  512. // Trtri will not perform the inversion if the matrix is singular, and returns
  513. // a boolean indicating whether the inversion was successful.
  514. func Trtri(a blas64.Triangular) (ok bool) {
  515. return lapack64.Dtrtri(a.Uplo, a.Diag, a.N, a.Data, max(1, a.Stride))
  516. }
  517. // Trtrs solves a triangular system of the form A * X = B or Aᵀ * X = B. Trtrs
  518. // returns whether the solve completed successfully. If A is singular, no solve is performed.
  519. func Trtrs(trans blas.Transpose, a blas64.Triangular, b blas64.General) (ok bool) {
  520. return lapack64.Dtrtrs(a.Uplo, trans, a.Diag, a.N, b.Cols, a.Data, max(1, a.Stride), b.Data, max(1, b.Stride))
  521. }
  522. // Geev computes the eigenvalues and, optionally, the left and/or right
  523. // eigenvectors for an n×n real nonsymmetric matrix A.
  524. //
  525. // The right eigenvector v_j of A corresponding to an eigenvalue λ_j
  526. // is defined by
  527. // A v_j = λ_j v_j,
  528. // and the left eigenvector u_j corresponding to an eigenvalue λ_j is defined by
  529. // u_jᴴ A = λ_j u_jᴴ,
  530. // where u_jᴴ is the conjugate transpose of u_j.
  531. //
  532. // On return, A will be overwritten and the left and right eigenvectors will be
  533. // stored, respectively, in the columns of the n×n matrices VL and VR in the
  534. // same order as their eigenvalues. If the j-th eigenvalue is real, then
  535. // u_j = VL[:,j],
  536. // v_j = VR[:,j],
  537. // and if it is not real, then j and j+1 form a complex conjugate pair and the
  538. // eigenvectors can be recovered as
  539. // u_j = VL[:,j] + i*VL[:,j+1],
  540. // u_{j+1} = VL[:,j] - i*VL[:,j+1],
  541. // v_j = VR[:,j] + i*VR[:,j+1],
  542. // v_{j+1} = VR[:,j] - i*VR[:,j+1],
  543. // where i is the imaginary unit. The computed eigenvectors are normalized to
  544. // have Euclidean norm equal to 1 and largest component real.
  545. //
  546. // Left eigenvectors will be computed only if jobvl == lapack.LeftEVCompute,
  547. // otherwise jobvl must be lapack.LeftEVNone.
  548. // Right eigenvectors will be computed only if jobvr == lapack.RightEVCompute,
  549. // otherwise jobvr must be lapack.RightEVNone.
  550. // For other values of jobvl and jobvr Geev will panic.
  551. //
  552. // On return, wr and wi will contain the real and imaginary parts, respectively,
  553. // of the computed eigenvalues. Complex conjugate pairs of eigenvalues appear
  554. // consecutively with the eigenvalue having the positive imaginary part first.
  555. // wr and wi must have length n, and Geev will panic otherwise.
  556. //
  557. // work must have length at least lwork and lwork must be at least max(1,4*n) if
  558. // the left or right eigenvectors are computed, and at least max(1,3*n) if no
  559. // eigenvectors are computed. For good performance, lwork must generally be
  560. // larger. On return, optimal value of lwork will be stored in work[0].
  561. //
  562. // If lwork == -1, instead of performing Geev, the function only calculates the
  563. // optimal vaule of lwork and stores it into work[0].
  564. //
  565. // On return, first will be the index of the first valid eigenvalue.
  566. // If first == 0, all eigenvalues and eigenvectors have been computed.
  567. // If first is positive, Geev failed to compute all the eigenvalues, no
  568. // eigenvectors have been computed and wr[first:] and wi[first:] contain those
  569. // eigenvalues which have converged.
  570. func Geev(jobvl lapack.LeftEVJob, jobvr lapack.RightEVJob, a blas64.General, wr, wi []float64, vl, vr blas64.General, work []float64, lwork int) (first int) {
  571. n := a.Rows
  572. if a.Cols != n {
  573. panic("lapack64: matrix not square")
  574. }
  575. if jobvl == lapack.LeftEVCompute && (vl.Rows != n || vl.Cols != n) {
  576. panic("lapack64: bad size of VL")
  577. }
  578. if jobvr == lapack.RightEVCompute && (vr.Rows != n || vr.Cols != n) {
  579. panic("lapack64: bad size of VR")
  580. }
  581. 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)
  582. }