Linear Algebra¶
Linear algebra functions in Julia are largely implemented by calling functions from LAPACK. Sparse factorizations call functions from SuiteSparse.
-
*(A, B) Matrix multiplication
-
\(A, B) Matrix division using a polyalgorithm. For input matrices
AandB, the resultXis such thatA*X == BwhenAis square. The solver that is used depends upon the structure ofA. A direct solver is used for upper- or lower triangularA. For HermitianA(equivalent to symmetricAfor non-complexA) theBunchKaufmanfactorization is used. Otherwise an LU factorization is used. For rectangularAthe result is the minimum-norm least squares solution computed by reducingAto bidiagonal form and solving the bidiagonal least squares problem. For sparse, squareAthe LU factorization (from UMFPACK) is used.
-
dot(x, y)¶ Compute the dot product. For complex vectors, the first vector is conjugated.
-
cross(x, y)¶ Compute the cross product of two 3-vectors.
-
rref(A)¶ Compute the reduced row echelon form of the matrix A.
-
factorize(A)¶ Compute a convenient factorization (including LU, Cholesky, Bunch-Kaufman, Triangular) of A, based upon the type of the input matrix. The return value can then be reused for efficient solving of multiple systems. For example:
A=factorize(A); x=A\\b; y=A\\C.
-
factorize!(A)¶ factorize!is the same asfactorize(), but saves space by overwriting the inputA, instead of creating a copy.
-
lu(A) → L, U, p¶ Compute the LU factorization of
A, such thatA[p,:] = L*U.
-
lufact(A[, pivot=true]) → F¶ Compute the LU factorization of
A. The return type ofFdepends on the type ofA. In most cases, ifAis a subtypeSof AbstractMatrix with an element typeT`supporting+,-,*and/the return type isLU{T,S{T}}. If pivoting is chosen (default) the element type should also supportabsand<. WhenAis sparse and have element of typeFloat32,Float64,Complex{Float32}, orComplex{Float64}the return type isUmfpackLU. Some examples are shown in the table below.Type of input AType of output FRelationship between FandAMatrix()LUF[:L]*F[:U] == A[F[:p], :]Tridiagonal()LU{T,Tridiagonal{T}}N/A SparseMatrixCSC()UmfpackLUF[:L]*F[:U] == Rs .* A[F[:p], F[:q]]The individual components of the factorization
Fcan be accessed by indexing:Component Description LULU{T,Tridiagonal{T}}UmfpackLUF[:L]L(lower triangular) part ofLU✓ ✓ F[:U]U(upper triangular) part ofLU✓ ✓ F[:p](right) permutation Vector✓ ✓ F[:P](right) permutation Matrix✓ F[:q]left permutation Vector✓ F[:Rs]Vectorof scaling factors✓ F[:(:)](L,U,p,q,Rs)components✓ Supported function LULU{T,Tridiagonal{T}}UmfpackLU/✓ \✓ ✓ ✓ cond✓ ✓ det✓ ✓ ✓ size✓ ✓
-
lufact!(A) → LU¶ lufact!is the same aslufact(), but saves space by overwriting the input A, instead of creating a copy. For sparseAthenzvalfield is not overwritten but the index fields,colptrandrowvalare decremented in place, converting from 1-based indices to 0-based indices.
-
chol(A[, LU]) → F¶ Compute the Cholesky factorization of a symmetric positive definite matrix
Aand return the matrixF. IfLUis:L(Lower),A = L*L'. IfLUis:U(Upper),A = R'*R.
-
cholfact(A, [LU,][pivot=false,][tol=-1.0]) → Cholesky¶ Compute the Cholesky factorization of a dense symmetric positive (semi)definite matrix
Aand return either aCholeskyifpivot=falseorCholeskyPivotedifpivot=true.LUmay be:Lfor using the lower part or:Ufor the upper part. The default is to use:U. The triangular matrix can be obtained from the factorizationFwith:F[:L]andF[:U]. The following functions are available forCholeskyobjects:size,\,inv,det. ForCholeskyPivotedthere is also defined arank. Ifpivot=falseaPosDefExceptionexception is thrown in case the matrix is not positive definite. The argumenttoldetermines the tolerance for determining the rank. For negative values, the tolerance is the machine precision.
-
cholfact(A[, ll]) → CholmodFactor Compute the sparse Cholesky factorization of a sparse matrix
A. IfAis Hermitian its Cholesky factor is determined. IfAis not Hermitian the Cholesky factor ofA*A'is determined. A fill-reducing permutation is used. Methods forsize,solve,\,findn_nzs,diag,detandlogdet. One of the solve methods includes an integer argument that can be used to solve systems involving parts of the factorization only. The optional boolean argument,lldetermines whether the factorization returned is of theA[p,p] = L*L'form, whereLis lower triangular orA[p,p] = L*Diagonal(D)*L'form whereLis unit lower triangular andDis a non-negative vector. The default is LDL.
-
cholfact!(A, [LU,][pivot=false,][tol=-1.0]) → Cholesky¶ cholfact!is the same ascholfact(), but saves space by overwriting the inputA, instead of creating a copy.
-
ldltfact(A) → LDLtFactorization¶ Compute a factorization of a positive definite matrix
Asuch thatA=L*Diagonal(d)*L'whereLis a unit lower triangular matrix anddis a vector with non-negative elements.
-
qr(A, [pivot=false,][thin=true]) → Q, R, [p]¶ Compute the (pivoted) QR factorization of
Asuch that eitherA = Q*RorA[:,p] = Q*R. Also seeqrfact. The default is to compute a thin factorization. Note thatRis not extended with zeros when the fullQis requested.
-
qrfact(A[, pivot=false]) → F¶ Computes the QR factorization of
A. The return type ofFdepends on the element type ofAand whether pivoting is specified (withpivot=true).Return type eltype(A)pivotRelationship between FandAQRnot BlasFloateither A==F[:Q]*F[:R]QRCompactWYBlasFloatfalseA==F[:Q]*F[:R]QRPivotedBlasFloattrueA[:,F[:p]]==F[:Q]*F[:R]BlasFloatrefers to any of:Float32,Float64,Complex64orComplex128.The individual components of the factorization
Fcan be accessed by indexing:Component Description QRQRCompactWYQRPivotedF[:Q]Q(orthogonal/unitary) part ofQR✓ ( QRPackedQ)✓ ( QRCompactWYQ)✓ ( QRPackedQ)F[:R]R(upper right triangular) part ofQR✓ ✓ ✓ F[:p]pivot Vector✓ F[:P](pivot) permutation Matrix✓ The following functions are available for the
QRobjects:size,\. WhenAis rectangular,\will return a least squares solution and if the solution is not unique, the one with smallest norm is returned.Multiplication with respect to either thin or full
Qis allowed, i.e. bothF[:Q]*F[:R]andF[:Q]*Aare supported. AQmatrix can be converted into a regular matrix withfull()which has a named argumentthin.注解
qrfactreturns multiple types because LAPACK uses several representations that minimize the memory storage requirements of products of Householder elementary reflectors, so that theQandRmatrices can be stored compactly rather as two separate dense matrices.The data contained in
QRorQRPivotedcan be used to construct theQRPackedQtype, which is a compact representation of the rotation matrix:\[Q = \prod_{i=1}^{\min(m,n)} (I - \tau_i v_i v_i^T)\]where \(\tau_i\) is the scale factor and \(v_i\) is the projection vector associated with the \(i^{th}\) Householder elementary reflector.
The data contained in
QRCompactWYcan be used to construct theQRCompactWYQtype, which is a compact representation of the rotation matrix\[Q = I + Y T Y^T\]where
Yis \(m \times r\) lower trapezoidal andTis \(r \times r\) upper triangular. The compact WY representation [Schreiber1989] is not to be confused with the older, WY representation [Bischof1987]. (The LAPACK documentation usesVin lieu ofY.)[Bischof1987] C Bischof and C Van Loan, The WY representation for products of Householder matrices, SIAM J Sci Stat Comput 8 (1987), s2-s13. doi:10.1137/0908009 [Schreiber1989] R Schreiber and C Van Loan, A storage-efficient WY representation for products of Householder transformations, SIAM J Sci Stat Comput 10 (1989), 53-57. doi:10.1137/0910005
-
qrfact!(A[, pivot=false])¶ qrfact!is the same asqrfact(), but saves space by overwriting the inputA, instead of creating a copy.
-
bkfact(A) → BunchKaufman¶ Compute the Bunch-Kaufman [Bunch1977] factorization of a real symmetric or complex Hermitian matrix
Aand return aBunchKaufmanobject. The following functions are available forBunchKaufmanobjects:size,\,inv,issym,ishermitian.
| [Bunch1977] | J R Bunch and L Kaufman, Some stable methods for calculating inertia and solving symmetric linear systems, Mathematics of Computation 31:137 (1977), 163-179. url. |
-
bkfact!(A) → BunchKaufman¶ bkfact!is the same asbkfact(), but saves space by overwriting the inputA, instead of creating a copy.
-
sqrtm(A)¶ Compute the matrix square root of
A. IfB = sqrtm(A), thenB*B == Awithin roundoff error.sqrtmuses a polyalgorithm, computing the matrix square root using Schur factorizations (schurfact()) unless it detects the matrix to be Hermitian or real symmetric, in which case it computes the matrix square root from an eigendecomposition (eigfact()). In the latter situation for positive definite matrices, the matrix square root hasRealelements, otherwise it hasComplexelements.
-
eig(A,[irange,][vl,][vu,][permute=true,][scale=true]) → D, V¶ Compute eigenvalues and eigenvectors of
A. Seeeigfact()for details on thebalancekeyword argument.julia> eig([1.0 0.0 0.0; 0.0 3.0 0.0; 0.0 0.0 18.0]) ([1.0,3.0,18.0], 3x3 Array{Float64,2}: 1.0 0.0 0.0 0.0 1.0 0.0 0.0 0.0 1.0)
eigis a wrapper aroundeigfact(), extracting all parts of the factorization to a tuple; where possible, usingeigfact()is recommended.
-
eig(A, B) → D, V Computes generalized eigenvalues and vectors of
Awith respect toB.eigis a wrapper aroundeigfact(), extracting all parts of the factorization to a tuple; where possible, usingeigfact()is recommended.
-
eigvals(A,[irange,][vl,][vu])¶ Returns the eigenvalues of
A. IfAisSymmetric(),Hermitian()orSymTridiagonal(), it is possible to calculate only a subset of the eigenvalues by specifying either aUnitRange()irangecovering indices of the sorted eigenvalues, or a pairvlandvufor the lower and upper boundaries of the eigenvalues.For general non-symmetric matrices it is possible to specify how the matrix is balanced before the eigenvector calculation. The option
permute=truepermutes the matrix to become closer to upper triangular, andscale=truescales the matrix by its diagonal elements to make rows and columns more equal in norm. The default istruefor both options.
-
eigmax(A)¶ Returns the largest eigenvalue of
A.
-
eigmin(A)¶ Returns the smallest eigenvalue of
A.
-
eigvecs(A, [eigvals,][permute=true,][scale=true])¶ Returns the eigenvectors of
A. Thepermuteandscalekeywords are the same as foreigfact().For
SymTridiagonal()matrices, if the optional vector of eigenvalueseigvalsis specified, returns the specific corresponding eigenvectors.
-
eigfact(A,[il,][iu,][vl,][vu,][permute=true,][scale=true])¶ Compute the eigenvalue decomposition of
Aand return anEigenobject. IfFis the factorization object, the eigenvalues can be accessed withF[:values]and the eigenvectors withF[:vectors]. The following functions are available forEigenobjects:inv,det.If
AisSymmetric,HermitianorSymTridiagonal, it is possible to calculate only a subset of the eigenvalues by specifying either a UnitRange`irangecovering indices of the sorted eigenvalues or a pairvlandvufor the lower and upper boundaries of the eigenvalues.For general non-symmetric matrices it is possible to specify how the matrix is balanced before the eigenvector calculation. The option
permute=truepermutes the matrix to become closer to upper triangular, andscale=truescales the matrix by its diagonal elements to make rows and columns more equal in norm. The default istruefor both options.
-
eigfact(A, B) Compute the generalized eigenvalue decomposition of
AandBand return anGeneralizedEigenobject. IfFis the factorization object, the eigenvalues can be accessed withF[:values]and the eigenvectors withF[:vectors].
-
eigfact!(A[, B])¶ eigfact!is the same aseigfact(), but saves space by overwriting the input A (and B), instead of creating a copy.
-
hessfact(A)¶ Compute the Hessenberg decomposition of
Aand return aHessenbergobject. IfFis the factorization object, the unitary matrix can be accessed withF[:Q]and the Hessenberg matrix withF[:H]. WhenQis extracted, the resulting type is theHessenbergQobject, and may be converted to a regular matrix withfull().
-
hessfact!(A)¶ hessfact!is the same ashessfact(), but saves space by overwriting the input A, instead of creating a copy.
-
schurfact(A) → Schur¶ Computes the Schur factorization of the matrix
A. The (quasi) triangular Schur factor can be obtained from theSchurobjectFwith eitherF[:Schur]orF[:T]and the unitary/orthogonal Schur vectors can be obtained withF[:vectors]orF[:Z]such thatA=F[:vectors]*F[:Schur]*F[:vectors]'. The eigenvalues ofAcan be obtained withF[:values].
-
schurfact!(A)¶ Computer the Schur factorization of
A, overwritingAin the process. Seeschurfact()
-
schur(A) → Schur[:T], Schur[:Z], Schur[:values]¶ See
schurfact()
-
schurfact(A, B) → GeneralizedSchur Computes the Generalized Schur (or QZ) factorization of the matrices
AandB. The (quasi) triangular Schur factors can be obtained from theSchurobjectFwithF[:S]andF[:T], the left unitary/orthogonal Schur vectors can be obtained withF[:left]orF[:Q]and the right unitary/orthogonal Schur vectors can be obtained withF[:right]orF[:Z]such thatA=F[:left]*F[:S]*F[:right]'andB=F[:left]*F[:T]*F[:right]'. The generalized eigenvalues ofAandBcan be obtained withF[:alpha]./F[:beta].
-
schur(A, B) → GeneralizedSchur[:S], GeneralizedSchur[:T], GeneralizedSchur[:Q], GeneralizedSchur[:Z] See
schurfact()
-
svdfact(A[, thin=true]) → SVD¶ Compute the Singular Value Decomposition (SVD) of
Aand return anSVDobject.U,S,VandVtcan be obtained from the factorizationFwithF[:U],F[:S],F[:V]andF[:Vt], such thatA = U*diagm(S)*Vt. Ifthinistrue, an economy mode decomposition is returned. The algorithm producesVtand henceVtis more efficient to extract thanV. The default is to produce a thin decomposition.
-
svdfact!(A[, thin=true]) → SVD¶ svdfact!is the same assvdfact(), but saves space by overwriting the input A, instead of creating a copy. Ifthinistrue, an economy mode decomposition is returned. The default is to produce a thin decomposition.
-
svd(A[, thin=true]) → U, S, V¶ Wrapper around
svdfactextracting all parts the factorization to a tuple. Direct use ofsvdfactis therefore generally more efficient. Computes the SVD of A, returningU, vectorS, andVsuch thatA == U*diagm(S)*V'. Ifthinistrue, an economy mode decomposition is returned. The default is to produce a thin decomposition.
-
svdvals(A)¶ Returns the singular values of
A.
-
svdvals!(A)¶ Returns the singular values of
A, while saving space by overwriting the input.
-
svdfact(A, B) → GeneralizedSVD Compute the generalized SVD of
AandB, returning aGeneralizedSVDFactorization objectF, such thatA = F[:U]*F[:D1]*F[:R0]*F[:Q]'andB = F[:V]*F[:D2]*F[:R0]*F[:Q]'.
-
svd(A, B) → U, V, Q, D1, D2, R0 Wrapper around
svdfactextracting all parts the factorization to a tuple. Direct use ofsvdfactis therefore generally more efficient. The function returns the generalized SVD ofAandB, returningU,V,Q,D1,D2, andR0such thatA = U*D1*R0*Q'andB = V*D2*R0*Q'.
-
svdvals(A, B) Return only the singular values from the generalized singular value decomposition of
AandB.
-
triu(M)¶ Upper triangle of a matrix.
-
triu!(M)¶ Upper triangle of a matrix, overwriting
Min the process.
-
tril(M)¶ Lower triangle of a matrix.
-
tril!(M)¶ Lower triangle of a matrix, overwriting
Min the process.
-
diagind(M[, k])¶ A
Rangegiving the indices of thek-th diagonal of the matrixM.
-
diag(M[, k])¶ The
k-th diagonal of a matrix, as a vector. Usediagmto construct a diagonal matrix.
-
diagm(v[, k])¶ Construct a diagonal matrix and place
von thek-th diagonal.
-
scale(A, b)¶
-
scale(b, A) Scale an array
Aby a scalarb, returning a new array.If
Ais a matrix andbis a vector, thenscale(A,b)scales each columniofAbyb[i](similar toA*diagm(b)), whilescale(b,A)scales each rowiofAbyb[i](similar todiagm(b)*A), returning a new array.Note: for large
A,scalecan be much faster thanA .* borb .* A, due to the use of BLAS.
-
scale!(A, b)¶
-
scale!(b, A) Scale an array
Aby a scalarb, similar toscale()but overwritingAin-place.If
Ais a matrix andbis a vector, thenscale!(A,b)scales each columniofAbyb[i](similar toA*diagm(b)), whilescale!(b,A)scales each rowiofAbyb[i](similar todiagm(b)*A), again operating in-place onA.
-
Tridiagonal(dl, d, du)¶ Construct a tridiagonal matrix from the lower diagonal, diagonal, and upper diagonal, respectively. The result is of type
Tridiagonaland provides efficient specialized linear solvers, but may be converted into a regular matrix withfull().
-
Bidiagonal(dv, ev, isupper)¶ Constructs an upper (
isupper=true) or lower (isupper=false) bidiagonal matrix using the given diagonal (dv) and off-diagonal (ev) vectors. The result is of typeBidiagonaland provides efficient specialized linear solvers, but may be converted into a regular matrix withfull().
-
SymTridiagonal(d, du)¶ Construct a real symmetric tridiagonal matrix from the diagonal and upper diagonal, respectively. The result is of type
SymTridiagonaland provides efficient specialized eigensolvers, but may be converted into a regular matrix withfull().
-
Woodbury(A, U, C, V)¶ Construct a matrix in a form suitable for applying the Woodbury matrix identity.
-
rank(M)¶ Compute the rank of a matrix.
-
norm(A[, p])¶ Compute the
p-norm of a vector or the operator norm of a matrixA, defaulting to thep=2-norm.For vectors,
pcan assume any numeric value (even though not all values produce a mathematically valid vector norm). In particular,norm(A, Inf)returns the largest value inabs(A), whereasnorm(A, -Inf)returns the smallest.For matrices, valid values of
pare1,2, orInf. (Note that for sparse matrices,p=2is currently not implemented.) Usevecnorm()to compute the Frobenius norm.
-
vecnorm(A[, p])¶ For any iterable container
A(including arrays of any dimension) of numbers, compute thep-norm (defaulting top=2) as ifAwere a vector of the corresponding length.For example, if
Ais a matrix andp=2, then this is equivalent to the Frobenius norm.
-
cond(M[, p])¶ Condition number of the matrix
M, computed using the operatorp-norm. Valid values forpare1,2(default), orInf.
-
condskeel(M[, x, p])¶ - \[\begin{split}\kappa_S(M, p) & = \left\Vert \left\vert M \right\vert \left\vert M^{-1} \right\vert \right\Vert_p \\ \kappa_S(M, x, p) & = \left\Vert \left\vert M \right\vert \left\vert M^{-1} \right\vert \left\vert x \right\vert \right\Vert_p\end{split}\]
Skeel condition number \(\kappa_S\) of the matrix
M, optionally with respect to the vectorx, as computed using the operatorp-norm.pisInfby default, if not provided. Valid values forpare1,2, orInf.This quantity is also known in the literature as the Bauer condition number, relative condition number, or componentwise relative condition number.
-
trace(M)¶ Matrix trace
-
det(M)¶ Matrix determinant
-
logdet(M)¶ Log of matrix determinant. Equivalent to
log(det(M)), but may provide increased accuracy and/or speed.
-
inv(M)¶ Matrix inverse
-
pinv(M)¶ Moore-Penrose pseudoinverse
-
null(M)¶ Basis for nullspace of
M.
-
repmat(A, n, m)¶ Construct a matrix by repeating the given matrix
ntimes in dimension 1 andmtimes in dimension 2.
-
repeat(A, inner = Int[], outer = Int[])¶ Construct an array by repeating the entries of
A. The i-th element ofinnerspecifies the number of times that the individual entries of the i-th dimension ofAshould be repeated. The i-th element ofouterspecifies the number of times that a slice along the i-th dimension ofAshould be repeated.
-
kron(A, B)¶ Kronecker tensor product of two vectors or two matrices.
-
blkdiag(A...)¶ Concatenate matrices block-diagonally. Currently only implemented for sparse matrices.
-
linreg(x, y) → [a; b]¶ Linear Regression. Returns
aandbsuch thata+b*xis the closest line to the given points(x,y). In other words, this function determines parameters[a, b]that minimize the squared error betweenyanda+b*x.Example:
using PyPlot; x = float([1:12]) y = [5.5; 6.3; 7.6; 8.8; 10.9; 11.79; 13.48; 15.02; 17.77; 20.81; 22.0; 22.99] a, b = linreg(x,y) # Linear regression plot(x, y, "o") # Plot (x,y) points plot(x, [a+b*i for i in x]) # Plot the line determined by the linear regression
-
linreg(x, y, w) Weighted least-squares linear regression.
-
expm(A)¶ Matrix exponential.
-
lyap(A, C)¶ Computes the solution
Xto the continuous Lyapunov equationAX + XA' + C = 0, where no eigenvalue ofAhas a zero real part and no two eigenvalues are negative complex conjugates of each other.
-
sylvester(A, B, C)¶ Computes the solution
Xto the Sylvester equationAX + XB + C = 0, whereA,BandChave compatible dimensions andAand-Bhave no eigenvalues with equal real part.
-
issym(A) → Bool¶ Test whether a matrix is symmetric.
-
isposdef(A) → Bool¶ Test whether a matrix is positive definite.
-
isposdef!(A) → Bool¶ Test whether a matrix is positive definite, overwriting
Ain the processes.
-
istril(A) → Bool¶ Test whether a matrix is lower triangular.
-
istriu(A) → Bool¶ Test whether a matrix is upper triangular.
-
ishermitian(A) → Bool¶ Test whether a matrix is Hermitian.
-
transpose(A)¶ The transposition operator (
.').
-
ctranspose(A)¶ The conjugate transposition operator (
').
-
eigs(A, [B, ]; nev=6, which="LM", tol=0.0, maxiter=1000, sigma=nothing, ritzvec=true, v0=zeros((0, ))) -> (d, [v, ]nconv, niter, nmult, resid)¶ eigscomputes eigenvaluesdofAusing Lanczos or Arnoldi iterations for real symmetric or general nonsymmetric matrices respectively. IfBis provided, the generalized eigen-problem is solved. The following keyword arguments are supported:nev: Number of eigenvaluesncv: Number of Krylov vectors used in the computation; should satisfynev+1 <= ncv <= nfor real symmetric problems andnev+2 <= ncv <= nfor other problems; default isncv = max(20,2*nev+1).which: type of eigenvalues to compute. See the note below.whichtype of eigenvalues :LMeigenvalues of largest magnitude (default) :SMeigenvalues of smallest magnitude :LReigenvalues of largest real part :SReigenvalues of smallest real part :LIeigenvalues of largest imaginary part (nonsymmetric or complex Aonly):SIeigenvalues of smallest imaginary part (nonsymmetric or complex Aonly):BEcompute half of the eigenvalues from each end of the spectrum, biased in favor of the high end. (real symmetric Aonly)tol: tolerance (\(tol \le 0.0\) defaults toDLAMCH('EPS'))maxiter: Maximum number of iterations (default = 300)sigma: Specifies the level shift used in inverse iteration. Ifnothing(default), defaults to ordinary (forward) iterations. Otherwise, find eigenvalues close tosigmausing shift and invert iterations.ritzvec: Returns the Ritz vectorsv(eigenvectors) iftruev0: starting vector from which to start the iterations
eigsreturns thenevrequested eigenvalues ind, the corresponding Ritz vectorsv(only ifritzvec=true), the number of converged eigenvaluesnconv, the number of iterationsniterand the number of matrix vector multiplicationsnmult, as well as the final residual vectorresid.注解
The
sigmaandwhichkeywords interact: the description of eigenvalues searched for bywhichdo _not_ necessarily refer to the eigenvalues ofA, but rather the linear operator constructed by the specification of the iteration mode implied bysigma.sigmaiteration mode whichrefers to eigenvalues ofnothingordinary (forward) \(A\) real or complex inverse with level shift sigma\((A - \sigma I )^{-1}\)
-
peakflops(n; parallel=false)¶ peakflopscomputes the peak flop rate of the computer by using BLAS double precisiongemm!(). By default, if no arguments are specified, it multiplies a matrix of sizen x n, wheren = 2000. If the underlying BLAS is using multiple threads, higher flop rates are realized. The number of BLAS threads can be set withblas_set_num_threads(n).If the keyword argument
parallelis set totrue,peakflopsis run in parallel on all the worker processors. The flop rate of the entire parallel computer is returned. When running in parallel, only 1 BLAS thread is used. The argumentnstill refers to the size of the problem that is solved on each processor.
BLAS Functions¶
This module provides wrappers for some of the BLAS functions for
linear algebra. Those BLAS functions that overwrite one of the input
arrays have names ending in '!'.
Usually a function has 4 methods defined, one each for Float64,
Float32, Complex128 and Complex64 arrays.
-
dot(n, X, incx, Y, incy)¶ Dot product of two vectors consisting of
nelements of arrayXwith strideincxandnelements of arrayYwith strideincy.
-
dotu(n, X, incx, Y, incy)¶ Dot function for two complex vectors.
-
dotc(n, X, incx, U, incy)¶ Dot function for two complex vectors conjugating the first vector.
-
blascopy!(n, X, incx, Y, incy)¶ Copy
nelements of arrayXwith strideincxto arrayYwith strideincy. ReturnsY.
-
nrm2(n, X, incx)¶ 2-norm of a vector consisting of
nelements of arrayXwith strideincx.
-
asum(n, X, incx)¶ sum of the absolute values of the first
nelements of arrayXwith strideincx.
-
axpy!(n, a, X, incx, Y, incy)¶ Overwrite
Ywitha*X + Y. ReturnsY.
-
scal!(n, a, X, incx)¶ Overwrite
Xwitha*X. ReturnsX.
-
scal(n, a, X, incx)¶ Returns
a*X.
-
syrk!(uplo, trans, alpha, A, beta, C)¶ Rank-k update of the symmetric matrix
Casalpha*A*A.' + beta*Coralpha*A.'*A + beta*Caccording to whethertransis ‘N’ or ‘T’. Whenuplois ‘U’ the upper triangle ofCis updated (‘L’ for lower triangle). ReturnsC.
-
syrk(uplo, trans, alpha, A)¶ Returns either the upper triangle or the lower triangle, according to
uplo(‘U’ or ‘L’), ofalpha*A*A.'oralpha*A.'*A, according totrans(‘N’ or ‘T’).
-
herk!(uplo, trans, alpha, A, beta, C)¶ Methods for complex arrays only. Rank-k update of the Hermitian matrix
Casalpha*A*A' + beta*Coralpha*A'*A + beta*Caccording to whethertransis ‘N’ or ‘T’. Whenuplois ‘U’ the upper triangle ofCis updated (‘L’ for lower triangle). ReturnsC.
-
herk(uplo, trans, alpha, A)¶ Methods for complex arrays only. Returns either the upper triangle or the lower triangle, according to
uplo(‘U’ or ‘L’), ofalpha*A*A'oralpha*A'*A, according totrans(‘N’ or ‘T’).
-
gbmv!(trans, m, kl, ku, alpha, A, x, beta, y)¶ Update vector
yasalpha*A*x + beta*yoralpha*A'*x + beta*yaccording totrans(‘N’ or ‘T’). The matrixAis a general band matrix of dimensionmbysize(A,2)withklsub-diagonals andkusuper-diagonals. Returns the updatedy.
-
gbmv(trans, m, kl, ku, alpha, A, x, beta, y)¶ Returns
alpha*A*xoralpha*A'*xaccording totrans(‘N’ or ‘T’). The matrixAis a general band matrix of dimensionmbysize(A,2)withklsub-diagonals andkusuper-diagonals.
-
sbmv!(uplo, k, alpha, A, x, beta, y)¶ Update vector
yasalpha*A*x + beta*ywhereAis a a symmetric band matrix of ordersize(A,2)withksuper-diagonals stored in the argumentA. The storage layout forAis described the reference BLAS module, level-2 BLAS at http://www.netlib.org/lapack/explore-html/.Returns the updated
y.
-
sbmv(uplo, k, alpha, A, x)¶ Returns
alpha*A*xwhereAis a symmetric band matrix of ordersize(A,2)withksuper-diagonals stored in the argumentA.
-
sbmv(uplo, k, A, x) Returns
A*xwhereAis a symmetric band matrix of ordersize(A,2)withksuper-diagonals stored in the argumentA.
-
gemm!(tA, tB, alpha, A, B, beta, C)¶ Update
Casalpha*A*B + beta*Cor the other three variants according totA(transposeA) andtB. Returns the updatedC.
-
gemm(tA, tB, alpha, A, B)¶ Returns
alpha*A*Bor the other three variants according totA(transposeA) andtB.
-
gemm(tA, tB, A, B) Returns
A*Bor the other three variants according totA(transposeA) andtB.
-
gemv!(tA, alpha, A, x, beta, y)¶ Update the vector
yasalpha*A*x + beta*xoralpha*A'x + beta*xaccording totA(transposeA). Returns the updatedy.
-
gemv(tA, alpha, A, x)¶ Returns
alpha*A*xoralpha*A'xaccording totA(transposeA).
-
gemv(tA, A, x) Returns
A*xorA'xaccording totA(transposeA).
-
symm!(side, ul, alpha, A, B, beta, C)¶ Update
Casalpha*A*B + beta*Coralpha*B*A + beta*Caccording toside.Ais assumed to be symmetric. Only theultriangle ofAis used. Returns the updatedC.
-
symm(side, ul, alpha, A, B)¶ Returns
alpha*A*Boralpha*B*Aaccording toside.Ais assumed to be symmetric. Only theultriangle ofAis used.
-
symm(side, ul, A, B) Returns
A*BorB*Aaccording toside.Ais assumed to be symmetric. Only theultriangle ofAis used.
-
symm(tA, tB, alpha, A, B) Returns
alpha*A*Bor the other three variants according totA(transposeA) andtB.
-
symv!(ul, alpha, A, x, beta, y)¶ Update the vector
yasalpha*A*y + beta*y.Ais assumed to be symmetric. Only theultriangle ofAis used. Returns the updatedy.
-
symv(ul, alpha, A, x)¶ Returns
alpha*A*x.Ais assumed to be symmetric. Only theultriangle ofAis used.
-
symv(ul, A, x) Returns
A*x.Ais assumed to be symmetric. Only theultriangle ofAis used.
-
trmm!(side, ul, tA, dA, alpha, A, B)¶ Update
Basalpha*A*Bor one of the other three variants determined byside(A on left or right) andtA(transpose A). Only theultriangle ofAis used.dAindicates ifAis unit-triangular (the diagonal is assumed to be all ones). Returns the updatedB.
-
trmm(side, ul, tA, dA, alpha, A, B)¶ Returns
alpha*A*Bor one of the other three variants determined byside(A on left or right) andtA(transpose A). Only theultriangle ofAis used.dAindicates ifAis unit-triangular (the diagonal is assumed to be all ones).
-
trsm!(side, ul, tA, dA, alpha, A, B)¶ Overwrite
Bwith the solution toA*X = alpha*Bor one of the other three variants determined byside(A on left or right ofX) andtA(transpose A). Only theultriangle ofAis used.dAindicates ifAis unit-triangular (the diagonal is assumed to be all ones). Returns the updatedB.
-
trsm(side, ul, tA, dA, alpha, A, B)¶ Returns the solution to
A*X = alpha*Bor one of the other three variants determined byside(A on left or right ofX) andtA(transpose A). Only theultriangle ofAis used.dAindicates ifAis unit-triangular (the diagonal is assumed to be all ones).
-
trmv!(side, ul, tA, dA, alpha, A, b)¶ Update
basalpha*A*bor one of the other three variants determined byside(A on left or right) andtA(transpose A). Only theultriangle ofAis used.dAindicates ifAis unit-triangular (the diagonal is assumed to be all ones). Returns the updatedb.
-
trmv(side, ul, tA, dA, alpha, A, b)¶ Returns
alpha*A*bor one of the other three variants determined byside(A on left or right) andtA(transpose A). Only theultriangle ofAis used.dAindicates ifAis unit-triangular (the diagonal is assumed to be all ones).
-
trsv!(side, ul, tA, dA, alpha, A, b)¶ Overwrite
bwith the solution toA*X = alpha*bor one of the other three variants determined byside(A on left or right ofX) andtA(transpose A). Only theultriangle ofAis used.dAindicates ifAis unit-triangular (the diagonal is assumed to be all ones). Returns the updatedb.
-
trsv(side, ul, tA, dA, alpha, A, b)¶ Returns the solution to
A*X = alpha*bor one of the other three variants determined byside(A on left or right ofX) andtA(transpose A). Only theultriangle ofAis used.dAindicates ifAis unit-triangular (the diagonal is assumed to be all ones).
-
blas_set_num_threads(n)¶ Set the number of threads the BLAS library should use.