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