Skip to content

lapack.lapack64 #

Constants #

const bad_apply_ortho = 'lapack: bad ApplyOrtho'

Panic strings for bad enumeration values.

const bad_balance_job = 'lapack: bad BalanceJob'
const bad_diag = 'lapack: bad Diag'
const bad_direct = 'lapack: bad Direct'
const bad_ev_comp = 'lapack: bad EigenVectorsComp'
const bad_ev_how_many = 'lapack: bad EigenVectorsHowMany'
const bad_ev_job = 'lapack: bad EigenVectorsJob'
const bad_ev_side = 'lapack: bad EigenVectorsSide'
const bad_gsvd_job = 'lapack: bad GSVDJob'
const bad_gen_ortho = 'lapack: bad GenOrtho'
const bad_left_ev_job = 'lapack: bad LeftEigenVectorsJob'
const bad_matrix_type = 'lapack: bad MatrixType'
const bad_maximize_norm_x_job = 'lapack: bad MaximizeNormXJob'
const bad_norm = 'lapack: bad Norm'
const bad_ortho_comp = 'lapack: bad OrthoComp'
const bad_pivot = 'lapack: bad Pivot'
const bad_right_ev_job = 'lapack: bad RightEigenVectorsJob'
const bad_svd_job = 'lapack: bad SVDJob'
const bad_schur_comp = 'lapack: bad SchurComp'
const bad_schur_job = 'lapack: bad SchurJob'
const bad_side = 'lapack: bad Side'
const bad_sort = 'lapack: bad Sort'
const bad_store_v = 'lapack: bad StoreV'
const bad_trans = 'lapack: bad Trans'
const bad_update_schur_comp = 'lapack: bad UpdateSchurComp'
const bad_uplo = 'lapack: bad Uplo'
const both_svd_over = 'lapack: both jobU and jobVT are lapack.SingularValueDecompositionOverwrite'
const bad_ifst = 'lapack: ifst out of range'

Panic strings for bad numerical and string values.

const bad_ihi = 'lapack: ihi out of range'
const bad_ihiz = 'lapack: ihiz out of range'
const bad_ilo = 'lapack: ilo out of range'
const bad_iloz = 'lapack: iloz out of range'
const bad_ilst = 'lapack: ilst out of range'
const bad_isave = 'lapack: bad isave value'
const bad_ispec = 'lapack: bad ispec value'
const bad_j1 = 'lapack: j1 out of range'
const bad_jpvt = 'lapack: bad element of jpvt'
const bad_k1 = 'lapack: k1 out of range'
const bad_k2 = 'lapack: k2 out of range'
const bad_kacc22 = 'lapack: invalid value of kacc22'
const bad_kbot = 'lapack: kbot out of range'
const bad_ktop = 'lapack: ktop out of range'
const bad_l_work = 'lapack: insufficient declared workspace length'
const bad_mm = 'lapack: mm out of range'
const bad_n1 = 'lapack: bad value of n1'
const bad_n2 = 'lapack: bad value of n2'
const bad_na = 'lapack: bad value of na'
const bad_name = 'lapack: bad name'
const bad_nh = 'lapack: bad value of nh'
const bad_nw = 'lapack: bad value of nw'
const bad_pp = 'lapack: bad value of pp'
const bad_shifts = 'lapack: bad shifts'
const i0lt0 = 'lapack: i0 < 0'
const k_gtm = 'lapack: k > m'
const k_gtn = 'lapack: k > n'
const k_lt0 = 'lapack: k < 0'
const k_lt1 = 'lapack: k < 1'
const kd_lt0 = 'lapack: kd < 0'
const kl_lt0 = 'lapack: kl < 0'
const ku_lt0 = 'lapack: ku < 0'
const m_gtn = 'lapack: m > n'
const m_lt0 = 'lapack: m < 0'
const mm_lt0 = 'lapack: mm < 0'
const n0lt0 = 'lapack: n0 < 0'
const n_gtm = 'lapack: n > m'
const n_lt0 = 'lapack: n < 0'
const n_lt1 = 'lapack: n < 1'
const n_ltm = 'lapack: n < m'
const nan_c_from = 'lapack: cfrom is NaN'
const nan_c_to = 'lapack: cto is NaN'
const nb_gtm = 'lapack: nb > m'
const nb_gtn = 'lapack: nb > n'
const nb_lt0 = 'lapack: nb < 0'
const ncc_lt0 = 'lapack: ncc < 0'
const ncvt_lt0 = 'lapack: ncvt < 0'
const neg_a_norm = 'lapack: anorm < 0'
const neg_z = 'lapack: negative z value'
const nh_lt0 = 'lapack: nh < 0'
const not_isolated = 'lapack: block is not isolated'
const nrhs_lt0 = 'lapack: nrhs < 0'
const nru_lt0 = 'lapack: nru < 0'
const nshfts_lt0 = 'lapack: nshfts < 0'
const nshfts_odd = 'lapack: nshfts must be even'
const nv_lt0 = 'lapack: nv < 0'
const offset_gtm = 'lapack: offset > m'
const offset_lt0 = 'lapack: offset < 0'
const p_lt0 = 'lapack: p < 0'
const recur_lt0 = 'lapack: recur < 0'
const zero_c_from = 'lapack: zero cfrom'
const bad_len_alpha = 'lapack: bad length of alpha'

Panic strings for bad slice lengths.

const bad_len_beta = 'lapack: bad length of beta'
const bad_len_ipiv = 'lapack: bad length of ipiv'
const bad_len_jpiv = 'lapack: bad length of jpiv'
const bad_len_jpvt = 'lapack: bad length of jpvt'
const bad_len_k = 'lapack: bad length of k'
const bad_len_piv = 'lapack: bad length of piv'
const bad_len_selected = 'lapack: bad length of selected'
const bad_len_si = 'lapack: bad length of si'
const bad_len_sr = 'lapack: bad length of sr'
const bad_len_tau = 'lapack: bad length of tau'
const bad_len_wi = 'lapack: bad length of wi'
const bad_len_wr = 'lapack: bad length of wr'
const short_a = 'lapack: insufficient length of a'

Panic strings for insufficient slice lengths.

const short_ab = 'lapack: insufficient length of ab'
const short_auxv = 'lapack: insufficient length of auxv'
const short_b = 'lapack: insufficient length of b'
const short_c = 'lapack: insufficient length of c'
const short_c_norm = 'lapack: insufficient length of cnorm'
const short_d = 'lapack: insufficient length of d'
const short_dl = 'lapack: insufficient length of dl'
const short_du = 'lapack: insufficient length of du'
const short_e = 'lapack: insufficient length of e'
const short_f = 'lapack: insufficient length of f'
const short_h = 'lapack: insufficient length of h'
const short_i_work = 'lapack: insufficient length of iwork'
const short_isgn = 'lapack: insufficient length of isgn'
const short_q = 'lapack: insufficient length of q'
const short_rhs = 'lapack: insufficient length of rhs'
const short_s = 'lapack: insufficient length of s'
const short_scale = 'lapack: insufficient length of scale'
const short_t = 'lapack: insufficient length of t'
const short_tau = 'lapack: insufficient length of tau'
const short_tau_p = 'lapack: insufficient length of tauP'
const short_tau_q = 'lapack: insufficient length of tauQ'
const short_u = 'lapack: insufficient length of u'
const short_v = 'lapack: insufficient length of v'
const short_vl = 'lapack: insufficient length of vl'
const short_vr = 'lapack: insufficient length of vr'
const short_vt = 'lapack: insufficient length of vt'
const short_vn1 = 'lapack: insufficient length of vn1'
const short_vn2 = 'lapack: insufficient length of vn2'
const short_w = 'lapack: insufficient length of w'
const short_wh = 'lapack: insufficient length of wh'
const short_wv = 'lapack: insufficient length of wv'
const short_wi = 'lapack: insufficient length of wi'
const short_work = 'lapack: insufficient length of work'
const short_wr = 'lapack: insufficient length of wr'
const short_x = 'lapack: insufficient length of x'
const short_y = 'lapack: insufficient length of y'
const short_z = 'lapack: insufficient length of z'
const bad_ld_a = 'lapack: bad leading dimension of A'

Panic strings for bad leading dimensions of matrices.

const bad_ld_b = 'lapack: bad leading dimension of B'
const bad_ld_c = 'lapack: bad leading dimension of C'
const bad_ld_f = 'lapack: bad leading dimension of F'
const bad_ld_h = 'lapack: bad leading dimension of H'
const bad_ld_q = 'lapack: bad leading dimension of Q'
const bad_ld_t = 'lapack: bad leading dimension of T'
const bad_ld_u = 'lapack: bad leading dimension of U'
const bad_ld_v = 'lapack: bad leading dimension of V'
const bad_ld_vl = 'lapack: bad leading dimension of VL'
const bad_ld_vr = 'lapack: bad leading dimension of VR'
const bad_ld_vt = 'lapack: bad leading dimension of VT'
const bad_ld_w = 'lapack: bad leading dimension of W'
const bad_ld_wh = 'lapack: bad leading dimension of WH'
const bad_ld_wv = 'lapack: bad leading dimension of WV'
const bad_ld_work = 'lapack: bad leading dimension of Work'
const bad_ld_x = 'lapack: bad leading dimension of X'
const bad_ld_y = 'lapack: bad leading dimension of Y'
const bad_ld_z = 'lapack: bad leading dimension of Z'
const abs_inc_not_one = 'lapack: increment not one or negative one'

Panic strings for bad vector increments.

const bad_inc_x = 'lapack: incx <= 0'
const bad_inc_y = 'lapack: incy <= 0'
const zero_inc_v = 'lapack: incv == 0'

fn dgebal #

fn dgebal(job BalanceJob, n int, mut a []f64, lda int, mut scale []f64) int

dgebal balances an n×n matrix A. Balancing consists of two stages, permuting and scaling. Both steps are optional and depend on the value of job.

Permuting consists of applying a permutation matrix P such that the matrix that results from PᵀAP takes the upper block triangular form

[ T1 X Y ] Pᵀ A P = [ 0 B Z ], [ 0 0 T2 ]

where T1 and T2 are upper triangular matrices and B contains at least one nonzero off-diagonal element in each row and column. The indices ilo and ihi mark the starting and ending columns of the submatrix B. The eigenvalues of A isolated in the first 0 to ilo-1 and last ihi+1 to n-1 elements on the diagonal can be read off without any roundoff error.

Scaling consists of applying a diagonal similarity transformation D such that D^{-1}BD has the 1-norm of each row and its corresponding column nearly equal. The output matrix is

[ T1 X*D Y ] [ 0 inv(D)BD inv(D)*Z ]. [ 0 0 T2 ]

Scaling may reduce the 1-norm of the matrix, and improve the accuracy of the computed eigenvalues and/or eigenvectors.

job specifies the operations that will be performed on A. If job is BalanceJob.balance_none, dgebal sets scale[i] = 1 for all i and returns 0. If job is BalanceJob.permute, only permuting will be done. If job is BalanceJob.scale, only scaling will be done. If job is BalanceJob.permute_scale, both permuting and scaling will be done.

On return, scale will contain information about the permutations and scaling factors applied to A. If π(j) denotes the index of the column interchanged with column j, and D[j,j] denotes the scaling factor applied to column j, then

scale[j] == π(j), for j ∈ {0, ..., ilo-1, ihi+1, ..., n-1}, == D[j,j], for j ∈ {ilo, ..., ihi}.

scale must have length equal to n, otherwise dgebal will panic.

dgebal returns 0 on success, or a negative error code.

fn dgeev #

fn dgeev(jobvl LeftEigenVectorsJob, jobvr RightEigenVectorsJob, n int, mut a []f64, lda int, mut wr []f64, mut wi []f64, mut vl []f64, ldvl int, mut vr []f64, ldvr int) int

dgeev: minimal real 2x2 eigen solver (row-major input/output for this pure-V path). For n>2, returns info=1 (not implemented).

fn dgehrd #

fn dgehrd(n int, ilo int, ihi int, mut a []f64, lda int, tau []f64) int

dgehrd reduces a general real matrix A to upper Hessenberg form H by an orthogonal similarity transformation.

fn dgeqr2 #

fn dgeqr2(m int, n int, mut a []f64, lda int, mut tau []f64, mut work []f64)

dgeqr2 computes a QR factorization of the m×n matrix A.

In a QR factorization, Q is an m×m orthonormal matrix, and R is an upper triangular m×n matrix.

A is modified to contain the information to construct Q and R. The upper triangle of a contains the matrix R. The lower triangular elements (not including the diagonal) contain the elementary reflectors. tau is modified to contain the reflector scales. tau must have length min(m,n), and this function will panic otherwise.

The ith elementary reflector can be explicitly constructed by first extracting the

v[j] = 0 j < i v[j] = 1 j == i v[j] = a[j*lda+i] j > i

and computing H_i = I - tau[i] * v * vᵀ.

The orthonormal matrix Q can be constructed from a product of these elementary reflectors, Q = H_0 * H_1 * ... * H_{k-1}, where k = min(m,n).

work is temporary storage of length at least n and this function will panic otherwise.

dgeqr2 is an internal routine. It is exported for testing purposes.

fn dgeqrf #

fn dgeqrf(m int, n int, mut a []f64, lda int, mut tau []f64, mut work []f64, lwork int)

dgeqrf computes the QR factorization of the m×n matrix A using a blocked algorithm. See the documentation for dgeqr2 for a description of the parameters at entry and exit.

work is temporary storage, and lwork specifies the usable memory length. The length of work must be at least max(1, lwork) and lwork must be -1 or at least n, otherwise this function will panic. dgeqrf is a blocked QR factorization, but the block size is limited by the temporary space available. If lwork == -1, instead of performing dgeqrf, the optimal work length will be stored into work[0].

tau must have length min(m,n), and this function will panic otherwise.

fn dgesv #

fn dgesv(n int, nrhs int, mut a []f64, lda int, mut ipiv []int, mut b []f64, ldb int) int

dgesv computes the solution to a real system of linear equations

A * X = B

where A is an n×n matrix and X and B are n×nrhs matrices.

The LU decomposition with partial pivoting and row interchanges is used to factor A as

A = P * L * U

where P is a permutation matrix, L is unit lower triangular, and U is upper triangular. On return, the factors L and U are stored in a; the unit diagonal elements of L are not stored. The row pivot indices that define the permutation matrix P are stored in ipiv.

The factored form of A is then used to solve the system of equations A * X = B. On entry, b contains the right hand side matrix B. On return, if ok is true, b contains the solution matrix X.

fn dgesvd #

fn dgesvd(jobu SVDJob, jobvt SVDJob, m int, n int, mut a []f64, lda int, s []f64, mut u []f64, ldu int, mut vt []f64, ldvt int, superb []f64) int

dgesvd computes the singular value decomposition (SingularValueDecomposition) of a real matrix A.

fn dgetf2 #

fn dgetf2(m int, n int, mut a []f64, lda int, mut ipiv []int) bool

fn dgetrf #

fn dgetrf(m int, n int, mut a []f64, lda int, mut ipiv []int) bool

dgetrf computes the LU decomposition of an m×n matrix A using partial pivoting with row interchanges.

The LU decomposition is a factorization of A into

A = P * L * U

where P is a permutation matrix, L is a lower triangular with unit diagonal elements (lower trapezoidal if m > n), and U is upper triangular (upper trapezoidal if m < n).

On entry, a contains the matrix A. On return, L and U are stored in place into a, and P is represented by ipiv.

ipiv contains a sequence of row interchanges. It indicates that row i of the matrix was interchanged with ipiv[i]. ipiv must have length min(m,n), and Dgetrf will panic otherwise. ipiv is zero-indexed.

Dgetrf returns whether the matrix A is nonsingular. The LU decomposition will be computed regardless of the singularity of A, but the result should not be used to solve a system of equation.

fn dgetri #

fn dgetri(n int, mut a []f64, lda int, ipiv []int, mut work []f64, lwork int) bool

dgetri computes the inverse of the matrix A using the LU factorization computed by dgetrf. On entry, a contains the PLU decomposition of A as computed by dgetrf and on exit contains the reciprocal of the original matrix.

dgetri will not perform the inversion if the matrix is singular, and returns a boolean indicating whether the inversion was successful.

work is temporary storage, and lwork specifies the usable memory length. At minimum, lwork >= n and this function will panic otherwise. dgetri is a blocked inversion, but the block size is limited by the temporary space available. If lwork == -1, instead of performing dgetri, the optimal work length will be stored into work[0].

fn dgetrs #

fn dgetrs(trans blas.Transpose, n int, nrhs int, mut a []f64, lda int, mut ipiv []int, mut b []f64, ldb int)

dgetrs solves a system of equations using an LU factorization. The system of equations solved is

A * X = B if trans == blas.Trans Aᵀ * X = B if trans == blas.NoTrans

A is a general n×n matrix with stride lda. B is a general matrix of size n×nrhs.

On entry b contains the elements of the matrix B. On exit, b contains the elements of X, the solution to the system of equations.

a and ipiv contain the LU factorization of A and the permutation indices as computed by Dgetrf. ipiv is zero-indexed.

fn dlacpy #

fn dlacpy(uplo blas.Uplo, m int, n int, a []f64, lda int, mut b []f64, ldb int)

dlacpy copies the elements of A specified by uplo into B. Uplo can specify a triangular portion with blas.Uplo.upper or blas.Uplo.lower, or can specify all of the elements with blas.Uplo.all.

dlacpy is an internal routine. It is exported for testing purposes.

fn dlae2 #

fn dlae2(a f64, b f64, c f64) (f64, f64)

dlae2 computes the eigenvalues of a 2×2 symmetric matrix

[a b] [b c]

and returns the eigenvalue with the larger absolute value as rt1 and the smaller as rt2.

dlae2 is an internal routine. It is exported for testing purposes.

fn dlaev2 #

fn dlaev2(a f64, b f64, c f64) (f64, f64, f64, f64)

fn dlange #

fn dlange(norm MatrixNorm, m int, n int, a []f64, lda int, mut work []f64) f64

dlange returns the value of the specified norm of a general m×n matrix A:

MatrixNorm.max_abs: the maximum absolute value of any element. MatrixNorm.max_column_sum: the maximum column sum of the absolute values of the elements (1-norm). MatrixNorm.max_row_sum: the maximum row sum of the absolute values of the elements (infinity-norm). MatrixNorm.frobenius: the square root of the sum of the squares of the elements (Frobenius norm).

If norm == MatrixNorm.max_column_sum, work must be of length n, and this function will panic otherwise. There are no restrictions on work for the other matrix norms.

fn dlanst #

fn dlanst(norm MatrixNorm, n int, d []f64, e []f64) f64

dlanst computes the specified norm of a symmetric tridiagonal matrix A. The diagonal elements of A are stored in d and the off-diagonal elements are stored in e.

fn dlansy #

fn dlansy(norm MatrixNorm, uplo blas.Uplo, n int, a []f64, lda int, mut work []f64) f64

dlansy returns the value of the specified norm of an n×n symmetric matrix. If norm == MatrixNorm.max_column_sum or norm == MatrixNorm.max_row_sum, work must have length at least n, otherwise work is unused.

fn dlapy2 #

fn dlapy2(x f64, y f64) f64

dlapy2 is the LAPACK version of math.hypot.

dlapy2 is an internal routine. It is exported for testing purposes.

fn dlarf #

fn dlarf(side blas.Side, m int, n int, v []f64, incv int, tau f64, mut c []f64, ldc int, mut work []f64)

dlarf applies an elementary reflector H to an m×n matrix C:

C = H * C if side == .left C = C * H if side == .right

H is represented in the form

H = I - tau * v * vᵀ

where tau is a scalar and v is a vector.

work must have length at least m if side == .left and at least n if side == .right.

dlarf is an internal routine. It is exported for testing purposes.

fn dlarfb #

fn dlarfb(side blas.Side, trans blas.Transpose, direct Direct, store StoreV, m int, n int, k int, v []f64, ldv int, t []f64, ldt int, mut c []f64, ldc int, mut work []f64, ldwork int)

fn dlarfg #

fn dlarfg(n int, alpha f64, mut x []f64, incx int) (f64, f64)

fn dlarft #

fn dlarft(direct Direct, store StoreV, n int, k int, v []f64, ldv int, tau []f64, mut t []f64, ldt int)

dlarft forms the triangular factor T of a block reflector H, storing the answer in t.

H = I - V * T * Vᵀ if store == .column_wise H = I - Vᵀ * T * V if store == .row_wise

H is defined by a product of the elementary reflectors where

H = H_0 * H_1 * ... * H_{k-1} if direct == .forward H = H_{k-1} * ... * H_1 * H_0 if direct == .backward

t is a k×k triangular matrix. t is upper triangular if direct = .forward and lower triangular otherwise. This function will panic if t is not of sufficient size.

store describes the storage of the elementary reflectors in v. See dlarfb for a description of layout.

tau contains the scalar factors of the elementary reflectors H_i.

dlarft is an internal routine. It is exported for testing purposes.

fn dlartg #

fn dlartg(f f64, g f64) (f64, f64, f64)

dlartg generates a plane rotation so that

[ cs sn] * [f] = [r] [-sn cs] [g] = [0]

where cscs + snsn = 1.

This is a more accurate version of BLAS Drotg that uses scaling to avoid overflow or underflow, with the other differences that- cs >= 0

  • if g = 0, then cs = 1 and sn = 0
  • if f = 0 and g != 0, then cs = 0 and sn = sign(1,g)

dlartg is an internal routine. It is exported for testing purposes.

fn dlascl #

fn dlascl(kind MatrixType, kl int, ku int, cfrom f64, cto f64, m int, n int, mut a []f64, lda int)

dlascl multiplies an m×n matrix by the scalar cto/cfrom.

cfrom must not be zero, and cto and cfrom must not be NaN, otherwise dlascl will panic.

dlascl is an internal routine. It is exported for testing purposes.

fn dlaset #

fn dlaset(uplo blas.Uplo, m int, n int, alpha f64, beta f64, mut a []f64, lda int)

dlaset sets the off-diagonal elements of A to alpha, and the diagonal elements to beta. If uplo == blas.upper, only the elements in the upper triangular part are set. If uplo == blas.lower, only the elements in the lower triangular part are set. If uplo is otherwise, all of the elements of A are set.

dlaset is an internal routine. It is exported for testing purposes.

fn dlasr #

fn dlasr(side blas.Side, pivot Pivot, direct Direct, m int, n int, c []f64, s []f64, mut a []f64, lda int)

Dlasr applies a sequence of plane rotations to the m×n matrix A. This series of plane rotations is implicitly represented by a matrix P. P is multiplied by a depending on the value of side -- A = P * A if side == Side.left, A = A * Pᵀ if side == Side.right.

The exact value of P depends on the value of pivot, but in all cases P is implicitly represented by a series of 2×2 rotation matrices. The entries of rotation matrix k are defined by s[k] and c[k]

R(k) = [ c[k] s[k]] [-s[k] s[k]]

If direct == Direct.forward, the rotation matrices are applied as P = P(z-1) * ... * P(2) * P(1), while if direct == Direct.backward they are applied as P = P(1) * P(2) * ... * P(n).

pivot defines the mapping of the elements in R(k) to P(k). If pivot == Pivot.variable, the rotation is performed for the (k, k+1) plane.

P(k) = [1 ] [ ... ] [ 1 ] [ c[k] s[k] ] [ -s[k] c[k] ] [ 1 ] [ ... ] [ 1]

if pivot == Pivot.top, the rotation is performed for the (1, k+1) plane,

P(k) = [c[k] s[k] ] [ 1 ] [ ... ] [ 1 ] [-s[k] c[k] ] [ 1 ] [ ... ] [ 1]

and if pivot == Pivot.bottom, the rotation is performed for the (k, z) plane.

P(k) = [1 ] [ ... ] [ 1 ] [ c[k] s[k]] [ 1 ] [ ... ] [ 1 ] [ -s[k] c[k]]

s and c have length m - 1 if side == Side.left, and n - 1 if side == Side.right.

fn dlasrt #

fn dlasrt(s Sort, n int, mut d []f64)

dlasrt sorts the numbers in the input slice d. If s == .increasing, the elements are sorted in increasing order. If s == .decreasing, the elements are sorted in decreasing order. For other values of s dlasrt will panic.

dlasrt is an internal routine. It is exported for testing purposes.

fn dlassq #

fn dlassq(n int, x []f64, incx int, scale f64, sumsq f64) (f64, f64)

dlassq updates a sum of squares represented in scaled form. It returns the values scl and smsq such that

scl^2smsq = X[0]^2 + ... + X[n-1]^2 + scale^2sumsq

The value of sumsq is assumed to be non-negative.

fn dlaswp #

fn dlaswp(n int, mut a []f64, lda int, k1 int, k2 int, mut ipiv []int, incx int)

fn dlatrd #

fn dlatrd(uplo blas.Uplo, n int, nb int, mut a []f64, lda int, mut e []f64, mut tau []f64, mut w []f64, ldw int)

fn dorg2l #

fn dorg2l(m int, n int, k int, mut a []f64, lda int, tau []f64, mut work []f64)

dorg2l generates an m×n matrix Q with orthonormal columns which is defined as the last n columns of a product of k elementary reflectors of order m.

Q = H_{k-1} * ... * H_1 * H_0

See dgelqf for more information. It must be that m >= n >= k.

tau contains the scalar reflectors computed by dgeqlf. tau must have length at least k, and dorg2l will panic otherwise.

work contains temporary memory, and must have length at least n. dorg2l will panic otherwise.

dorg2l is an internal routine. It is exported for testing purposes.

fn dorg2r #

fn dorg2r(m int, n int, k int, mut a []f64, lda int, tau []f64, mut work []f64)

fn dorgql #

fn dorgql(m int, n int, k int, mut a []f64, lda int, tau []f64, mut work []f64, lwork int)

dorgql generates the m×n matrix Q with orthonormal columns defined as the last n columns of a product of k elementary reflectors of order m

Q = H_{k-1} * ... * H_1 * H_0.

It must hold that

0 <= k <= n <= m,

and dorgql will panic otherwise.

On entry, the (n-k+i)-th column of A must contain the vector which defines the elementary reflector H_i, for i=0,...,k-1, and tau[i] must contain its scalar factor. On return, a contains the m×n matrix Q.

tau must have length at least k, and dorgql will panic otherwise.

work must have length at least max(1,lwork), and lwork must be at least max(1,n), otherwise dorgql will panic. For optimum performance lwork must be a sufficiently large multiple of n.

If lwork == -1, instead of computing dorgql the optimal work length is stored into work[0].

dorgql is an internal routine. It is exported for testing purposes.

fn dorgqr #

fn dorgqr(m int, n int, k int, mut a []f64, lda int, tau []f64, mut work []f64, lwork int)

fn dorgtr #

fn dorgtr(uplo blas.Uplo, n int, mut a []f64, lda int, tau []f64, mut work []f64, lwork int)

dorgtr generates a real orthogonal matrix Q which is defined as the product of n-1 elementary reflectors of order n as returned by dsytrd.

The construction of Q depends on the value of uplo:

Q = H_{n-1} * ... * H_1 * H_0 if uplo == blas.Upper Q = H_0 * H_1 * ... * H_{n-1} if uplo == blas.Lower

where H_i is constructed from the elementary reflectors as computed by dsytrd. See the documentation for dsytrd for more information.

tau must have length at least n-1, and dorgtr will panic otherwise.

work is temporary storage, and lwork specifies the usable memory length. At minimum, lwork >= max(1,n-1), and dorgtr will panic otherwise. The amount of blocking is limited by the usable length. If lwork == -1, instead of computing dorgtr the optimal work length is stored into work[0].

dorgtr is an internal routine. It is exported for testing purposes.

fn dpotf2 #

fn dpotf2(ul blas.Uplo, n int, mut a []f64, lda int) bool

dpotf2 computes the Cholesky decomposition of the symmetric positive definite matrix a. If ul == .upper, then a is stored as an upper-triangular matrix, and a = Uᵀ U is stored in place into a. If ul == .lower, then a = L Lᵀ is computed and stored in-place into a. If a is not positive definite, false is returned. This is the unblocked version of the algorithm.

dpotf2 is an internal routine. It is exported for testing purposes.

fn dpotrf #

fn dpotrf(ul blas.Uplo, n int, mut a []f64, lda int) bool

fn dsteqr #

fn dsteqr(compz EigenVectorsComp, n int, mut d []f64, mut e []f64, mut z []f64, ldz int, mut work []f64) bool

dsteqr computes the eigenvalues and optionally the eigenvectors of a symmetric tridiagonal matrix using the implicit QL or QR method. The eigenvectors of a full or band symmetric matrix can also be found if dsytrd, dsptrd, or dsbtrd have been used to reduce this matrix to tridiagonal form.

d, on entry, contains the diagonal elements of the tridiagonal matrix. On exit, d contains the eigenvalues in ascending order. d must have length n and dsteqr will panic otherwise.

e, on entry, contains the off-diagonal elements of the tridiagonal matrix on entry, and is overwritten during the call to dsteqr. e must have length n-1 and dsteqr will panic otherwise.

z, on entry, contains the n×n orthogonal matrix used in the reduction to tridiagonal form if compz == lapack.EigenVectorsOrig. On exit, if compz == lapack.EigenVectorsOrig, z contains the orthonormal eigenvectors of the original symmetric matrix, and if compz == lapack.EigenVectorsTridiag, z contains the orthonormal eigenvectors of the symmetric tridiagonal matrix. z is not used if compz == lapack.EigenVectorsCompNone.

work must have length at least max(1, 2*n-2) if the eigenvectors are computed, and dsteqr will panic otherwise.

dsteqr is an internal routine. It is exported for testing purposes.

fn dsterf #

fn dsterf(n int, mut d []f64, mut e []f64) bool

dsterf computes all eigenvalues of a symmetric tridiagonal matrix using the Pal-Walker-Kahan variant of the QL or QR algorithm.

d contains the diagonal elements of the tridiagonal matrix on entry, and contains the eigenvalues in ascending order on exit. d must have length at least n, or dsterf will panic.

e contains the off-diagonal elements of the tridiagonal matrix on entry, and is overwritten during the call to dsterf. e must have length of at least n-1 or dsterf will panic.

dsterf is an internal routine. It is exported for testing purposes.

fn dsyev #

fn dsyev(jobz EigenVectorsJob, uplo blas.Uplo, n int, mut a []f64, lda int, mut w []f64, mut work []f64, lwork int)

fn dsytd2 #

fn dsytd2(uplo blas.Uplo, n int, mut a []f64, lda int, mut d []f64, mut e []f64, mut tau []f64)

Dsytd2 reduces a symmetric n×n matrix A to symmetric tridiagonal form T by an orthogonal similarity transformation

Qᵀ * A * Q = T

On entry, the matrix is contained in the specified triangle of a. On exit, if uplo == Uplo.upper, the diagonal and first super-diagonal of a are overwritten with the elements of T. The elements above the first super-diagonal are overwritten with the elementary reflectors that are used with the elements written to tau in order to construct Q. If uplo == Uplo.lower, the elements are written in the lower triangular region.

d must have length at least n. e and tau must have length at least n-1. Dsytd2 will panic if these sizes are not met.

Q is represented as a product of elementary reflectors. If uplo == Uplo.upper

Q = H_{n-2} * ... * H_1 * H_0

and if uplo == Uplo.lower

Q = H_0 * H_1 * ... * H_{n-2}

where

H_i = I - tau * v * vᵀ

where tau is stored in tau[i], and v is stored in a.

If uplo == Uplo.upper, v[0:i-1] is stored in A[0:i-1,i+1], v[i] = 1, and v[i+1:] = 0. The elements of a are

[ d e v2 v3 v4] [ d e v3 v4] [ d e v4] [ d e] [ d]

If uplo == Uplo.lower, v[0:i+1] = 0, v[i+1] = 1, and v[i+2:] is stored in A[i+2:n,i]. The elements of a are

[ d ] [ e d ] [v1 e d ] [v1 v2 e d ] [v1 v2 v3 e d]

fn dsytrd #

fn dsytrd(uplo blas.Uplo, n int, mut a []f64, lda int, mut d []f64, mut e []f64, mut tau []f64, mut work []f64, lwork int)

dsytrd reduces a symmetric n×n matrix A to symmetric tridiagonal form by an orthogonal similarity transformation

Qᵀ * A * Q = T

where Q is an orthonormal matrix and T is symmetric and tridiagonal.

On entry, a contains the elements of the input matrix in the triangle specified by uplo. On exit, the diagonal and sub/super-diagonal are overwritten by the corresponding elements of the tridiagonal matrix T. The remaining elements in the triangle, along with the array tau, contain the data to construct Q as the product of elementary reflectors.

If uplo == blas.upper, Q is constructed with

Q = H_{n-2} * ... * H_1 * H_0

where

H_i = I - tau_i * v * vᵀ

v is constructed as v[i+1:n] = 0, v[i] = 1, v[0:i-1] is stored in A[0:i-1, i+1]. The elements of A are

[ d e v1 v2 v3] [ d e v2 v3] [ d e v3] [ d e] [ e]

If uplo == blas.lower, Q is constructed with

Q = H_0 * H_1 * ... * H_{n-2}

where

H_i = I - tau_i * v * vᵀ

v is constructed as v[0:i+1] = 0, v[i+1] = 1, v[i+2:n] is stored in A[i+2:n, i]. The elements of A are

[ d ] [ e d ] [v0 e d ] [v0 v1 e d ] [v0 v1 v2 e d]

d must have length n, and e and tau must have length n-1. dsytrd will panic if these conditions are not met.

work is temporary storage, and lwork specifies the usable memory length. At minimum, lwork >= 1, and dsytrd will panic otherwise. The amount of blocking is limited by the usable length. If lwork == -1, instead of computing dsytrd the optimal work length is stored into work[0].

dsytrd is an internal routine. It is exported for testing purposes.

fn dtrti2 #

fn dtrti2(uplo blas.Uplo, diag blas.Diagonal, n int, mut a []f64, lda int)

dtrti2 computes the inverse of a triangular matrix, storing the result in place into a. This is the BLAS level 2 version of the algorithm.

dtrti2 is an internal routine. It is exported for testing purposes.

fn dtrtri #

fn dtrtri(uplo blas.Uplo, diag blas.Diagonal, n int, mut a []f64, lda int) bool

dtrtri computes the inverse of a triangular matrix, storing the result in place into a. This is the BLAS level 3 version of the algorithm which builds upon dtrti2 to operate on matrix blocks instead of only individual columns.

dtrtri will not perform the inversion if the matrix is singular, and returns a boolean indicating whether the inversion was successful.

fn iladlc #

fn iladlc(m int, n int, a []f64, lda int) int

iladlc scans a matrix for its last non-zero column. Returns -1 if the matrix is all zeros.

fn iladlr #

fn iladlr(m int, n int, a []f64, lda int) int

iladlr scans a matrix for its last non-zero row. Returns -1 if the matrix is all zeros.

enum ApplyOrtho #

enum ApplyOrtho as u8 {
	// Apply P or Pᵀ.
	apply_p = u8(`P`)
	// Apply Q or Qᵀ.
	apply_q = u8(`Q`)
}

ApplyOrtho specifies which orthogonal matrix is applied in Dormbr.

enum BalanceJob #

enum BalanceJob as u8 {
	permute       = u8(`P`)
	scale         = u8(`S`)
	permute_scale = u8(`B`)
	balance_none  = u8(`N`)
}

BalanceJob specifies matrix balancing operation.

enum Direct #

enum Direct as u8 {
	// Reflectors are right-multiplied, H_0 * H_1 * ... * H_{k-1}.
	forward = u8(`F`)
	// Reflectors are left-multiplied, H_{k-1} * ... * H_1 * H_0.
	backward = u8(`B`)
}

Direct specifies the direction of the multiplication for the Householder matrix.

enum EigenVectorsComp #

enum EigenVectorsComp as u8 {
	// Compute eigenvectors of the original symmetric matrix.
	ev_orig = u8(`V`)
	// Compute eigenvectors of the tridiagonal matrix.
	ev_tridiag = u8(`I`)
	// Do not compute eigenvectors.
	ev_comp_none = u8(`N`)
}

EigenVectorsComp specifies how eigenvectors are computed in Dsteqr.

enum EigenVectorsHowMany #

enum EigenVectorsHowMany as u8 {
	// Compute all right and/or left eigenvectors.
	ev_all = u8(`A`)
	// Compute all right and/or left eigenvectors multiplied by an input matrix.
	ev_all_mul_q = u8(`B`)
	// Compute selected right and/or left eigenvectors.
	ev_selected = u8(`S`)
}

EigenVectorsHowMany specifies which eigenvectors are computed in Dtrevc3 and how.

enum EigenVectorsJob #

enum EigenVectorsJob as u8 {
	// Compute eigenvectors.
	ev_compute = u8(`V`)
	// Do not compute eigenvectors.
	ev_none = u8(`N`)
}

EigenVectorsJob specifies whether eigenvectors are computed in Dsyev.

enum EigenVectorsSide #

enum EigenVectorsSide as u8 {
	// Compute only right eigenvectors.
	ev_right = u8(`R`)
	// Compute only left eigenvectors.
	ev_left = u8(`L`)
	// Compute both right and left eigenvectors.
	ev_both = u8(`B`)
}

EigenVectorsSide specifies what eigenvectors are computed in Dtrevc3.

enum GSVDJob #

enum GSVDJob as u8 {
	// Compute orthogonal matrix U.
	gsvd_u = u8(`U`)
	// Compute orthogonal matrix V.
	gsvd_v = u8(`V`)
	// Compute orthogonal matrix Q.
	gsvd_q = u8(`Q`)
	// Use unit-initialized matrix.
	gsvd_unit = u8(`I`)
	// Do not compute orthogonal matrix.
	gsvd_none = u8(`N`)
}

GSVDJob specifies the singular vector computation type for Generalized SingularValueDecomposition.

enum GenOrtho #

enum GenOrtho as u8 {
	// Generate Pᵀ.
	generate_pt = u8(`P`)
	// Generate Q.
	generate_q = u8(`Q`)
}

GenOrtho specifies which orthogonal matrix is generated in Dorgbr.

enum LeftEigenVectorsJob #

enum LeftEigenVectorsJob as u8 {
	// Compute left eigenvectors.
	left_ev_compute = u8(`V`)
	// Do not compute left eigenvectors.
	left_ev_none = u8(`N`)
}

LeftEigenVectorsJob specifies whether left eigenvectors are computed in Dgeev.

enum MatrixNorm #

enum MatrixNorm as u8 {
	// max(abs(A(i,j)))
	max_abs = u8(`M`)
	// Maximum absolute column sum (one norm)
	max_column_sum = u8(`O`)
	// Maximum absolute row sum (infinity norm)
	max_row_sum = u8(`I`)
	// Frobenius norm (sqrt of sum of squares)
	frobenius = u8(`F`)
}

MatrixNorm represents the kind of matrix norm to compute.

enum MatrixType #

enum MatrixType as u8 {
	// A general dense matrix.
	general = u8(`G`)
	// An upper triangular matrix.
	upper_tri = u8(`U`)
	// A lower triangular matrix.
	lower_tri = u8(`L`)
}

MatrixType represents the kind of matrix represented in the data.

enum MaximizeNormXJob #

enum MaximizeNormXJob as u8 {
	// Solve Z*x=h-f where h is a vector of ±1.
	local_look_ahead = 0
	// Compute an approximate null-vector e of Z, normalize e and solve Z*x=±e-f.
	normalized_null_vector = 2
}

MaximizeNormXJob specifies the heuristic method for computing a contribution to the reciprocal Dif-estimate in Dlatdf.

enum OrthoComp #

enum OrthoComp as u8 {
	// Do not compute the orthogonal matrix.
	ortho_none = u8(`N`)
	// The orthogonal matrix is formed explicitly and returned in the argument.
	ortho_explicit = u8(`I`)
	// The orthogonal matrix is post-multiplied into the matrix stored in the argument on entry.
	ortho_postmul = u8(`V`)
}

OrthoComp specifies whether and how the orthogonal matrix is computed in Dgghrd.

enum Pivot #

enum Pivot as u8 {
	variable = u8(`V`)
	top      = u8(`T`)
	bottom   = u8(`B`)
}

Pivot specifies the pivot type for plane rotations.

enum RightEigenVectorsJob #

enum RightEigenVectorsJob as u8 {
	// Compute right eigenvectors.
	right_ev_compute = u8(`V`)
	// Do not compute right eigenvectors.
	right_ev_none = u8(`N`)
}

RightEigenVectorsJob specifies whether right eigenvectors are computed in Dgeev.

enum SVDJob #

enum SVDJob as u8 {
	// Compute all columns of the orthogonal matrix U or V.
	svd_all = u8(`A`)
	// Compute the singular vectors and store them in the orthogonal matrix U or V.
	svd_store = u8(`S`)
	// Compute the singular vectors and overwrite them on the input matrix A.
	svd_overwrite = u8(`O`)
	// Do not compute singular vectors.
	svd_none = u8(`N`)
}

SVDJob specifies the singular vector computation type for SingularValueDecomposition.

enum SchurComp #

enum SchurComp as u8 {
	// Compute Schur vectors of the original matrix.
	schur_orig = u8(`V`)
	// Compute Schur vectors of the upper Hessenberg matrix.
	schur_hess = u8(`I`)
	// Do not compute Schur vectors.
	schur_none = u8(`N`)
}

SchurComp specifies whether and how the Schur vectors are computed in Dhseqr.

enum SchurJob #

enum SchurJob as u8 {
	eigenvalues_only      = u8(`E`)
	eigenvalues_and_schur = u8(`S`)
}

SchurJob specifies whether the Schur form is computed in Dhseqr.

enum Sort #

enum Sort as u8 {
	sort_increasing = u8(`I`)
	sort_decreasing = u8(`D`)
}

Sort is the sorting order.

enum StoreV #

enum StoreV as u8 {
	// Reflector stored in a column of the matrix.
	column_wise = u8(`C`)
	// Reflector stored in a row of the matrix.
	row_wise = u8(`R`)
}

StoreV indicates the storage direction of elementary reflectors.

enum UpdateSchurComp #

enum UpdateSchurComp as u8 {
	// Update the matrix of Schur vectors.
	update_schur = u8(`V`)
	// Do not update the matrix of Schur vectors.
	update_schur_none = u8(`N`)
}

UpdateSchurComp specifies whether the matrix of Schur vectors is updated in Dtrexc.