JuliaF64LinearAlgebraΒΆ
jla64.spad line 1 [edit on github]
Linear Algebra functions computed using Julia and its algorithms. 64 bits version.
- conditionNumber: (JuliaFloat64Matrix, JuliaFloat64) -> JuliaFloat64
conditionNumber(m, p)computes thep-condition number ofm.
- conditionNumber: JuliaFloat64Matrix -> JuliaFloat64
conditionNumber(m)computes the condition number ofm.
- condSkeel: JuliaFloat64Matrix -> JuliaFloat64
condSkeel(m)computes the Skeel condition number ofm.
- eigen!: JuliaFloat64Matrix -> Record(values: JuliaComplexF64Vector, vectors: JuliaComplexF64Matrix)
eigen!(m)computes the spectral decomposition ofmbut overwritesmto save memory space.
- eigen: JuliaFloat64Matrix -> Record(values: JuliaComplexF64Vector, vectors: JuliaComplexF64Matrix)
eigen(m)computes the spectral decomposition ofm.
- eigenSystem!: JuliaFloat64Matrix -> Record(values: JuliaComplexF64Vector, leftVectors: JuliaFloat64Matrix, rightVectors: JuliaFloat64Matrix)
eigenSystem!(m)computes the spectral decomposition ofmbut overwritesmto save memory space. If thej-th eigenvalue (values) is real, then the left eigenvectorsu(j) = column(lefVectors,j), thej-th column of lefVectors. If thej-th and (j+1)-steigenvalues form a complex conjugate pair, then the left eigenvectors areu(j) = column(lefVectors,j) + %i*column(lefVectors,j+1) andu(j+1) = column(lefVectors,j) - %i*column((lefVectors,j+1). This applieas also to righVectors.
- eigenSystem: JuliaFloat64Matrix -> Record(values: JuliaComplexF64Vector, leftVectors: JuliaFloat64Matrix, rightVectors: JuliaFloat64Matrix)
eigenSystem(m)computes the spectral decomposition ofm. If thej-th eigenvalue (values) is real, then the left eigenvectorsu(j) = column(lefVectors,j), thej-th column of lefVectors. If thej-th and (j+1)-steigenvalues form a complex conjugate pair, then the left eigenvectors areu(j) = column(lefVectors,j) + %i*column(lefVectors,j+1) andu(j+1) = column(lefVectors,j) - %i*column((lefVectors,j+1). This applieas also to righVectors.
- eigvals!: JuliaFloat64Matrix -> JuliaComplexF64Vector
eigvals!(m)returns the eigen values ofmbut overwritesmto save memory space.
- eigvals: JuliaFloat64Matrix -> JuliaComplexF64Vector
eigvals(m)returns the eigen values ofm.
- eigvecs: JuliaFloat64Matrix -> JuliaComplexF64Matrix
eigvecs(m)returns the eigen vectors ofm.
- exp: JuliaFloat64Matrix -> JuliaFloat64Matrix
exp(m)returns the matrix exponential ofm.
- jlPeakFlops: () -> JuliaFloat64
jlPeakFlops()returns the peak flop rate using matrix multiplication. You can modify the number of threads used or the BLAS/LAPACK libraries used to see if that fits your needs.
- log: JuliaFloat64Matrix -> JuliaComplexF64Matrix
log(m)tries to compute the principal matrix logarithm ofm. Otherwise, returns a non pricipal matrix logarithm ofmif possible.
- logDeterminant: JuliaFloat64Matrix -> JuliaFloat64
logDeterminant(m)computes the logarithm of the determinant ofm, possibly with more accuracy and avoding nder/overflow.
- lu!: JuliaFloat64Matrix -> Record(LU: JuliaFloat64Matrix, ipiv: JuliaInt64Vector)
lu!(m)computes the LU factorisation ofminm.
- lu: JuliaFloat64Matrix -> Record(LU: JuliaFloat64Matrix, L: JuliaFloat64Matrix, U: JuliaFloat64Matrix, ipiv: JuliaInt64Vector)
lu(m)computes the LU factorisation ofm.
- luReorder!: (JuliaFloat64Matrix, JuliaInt64Vector) -> JuliaFloat64Matrix
luOrder(mat, ipiv) returns mat in-place reordered with ipiv pivot indices.
- luReorder: (JuliaFloat64Matrix, JuliaInt64Vector) -> JuliaFloat64Matrix
luOrder(mat, ipiv) returns a copy of mat reordered with ipiv pivot indices.
- mpInverse: JuliaFloat64Matrix -> JuliaFloat64Matrix
mpInverse(m)returns the Moore-Penrose pseudo inverse ofm.
- norm: (JuliaFloat64Matrix, JuliaFloat64) -> JuliaFloat64
norm(m,p)computes thep-norm ofm.
- norm: (JuliaFloat64Vector, JuliaFloat64) -> JuliaFloat64
norm(v,p)computes thp-norm ofv.
- norm: JuliaFloat64Matrix -> JuliaFloat64
norm(m)computes the 2-norm ofm, also known as the Frobenius norm.
- norm: JuliaFloat64Vector -> JuliaFloat64
norm(v)computes the 2-norm ofv.
- normalize!: JuliaFloat64Matrix -> JuliaFloat64Matrix
normalize!(m)destructively normalizemsuch that its norm equals to 1.
- normalize!: JuliaFloat64Vector -> JuliaFloat64Vector
normalize!(v)destructively normalizevsuch that norm(v) equals to 1.
- normalize: JuliaFloat64Matrix -> JuliaFloat64Matrix
normalize(m)returns normalizedmsuch that its norm equals to 1.
- normalize: JuliaFloat64Vector -> JuliaFloat64Vector
normalize(v)returns normalizedvsuch that its norm equals to 1.
- operatorNorm: (JuliaFloat64Matrix, JuliaFloat64) -> JuliaFloat64
operatorNorm(m,p)computes the operator norm ofminduced by the vectorp-norm.
- operatorNorm: JuliaFloat64Matrix -> JuliaFloat64
operatorNorm(m)computes the operator norm ofminduced by the vector 2-norm.
- rank!: (JuliaFloat64Matrix, JuliaFloat64) -> NonNegativeInteger
rank!(m, tol)computes rank ofm. Counts singular value with magnitude greater than tol but overwritesmto save memory space.
- rank: (JuliaFloat64Matrix, JuliaFloat64) -> NonNegativeInteger
rank(m, tol)computes rank ofm. Counts singular value with magnitude greater than tol.
- solve!: (JuliaFloat64Matrix, JuliaFloat64Matrix) -> JuliaFloat64Matrix
solve!(A,B)solves the matrix equation A*X=B. OverwritesBwith matrixXand returnsX.
- solve: (JuliaFloat64Matrix, JuliaFloat64Matrix) -> JuliaFloat64Matrix
solve(A,B)solves the matrix equation A*X=B, and returnsX.
- sqrt: JuliaFloat64Matrix -> JuliaComplexF64Matrix
sqrt(m)returns the principal square root ofm.
- svd!: JuliaFloat64Matrix -> Record(U: JuliaFloat64Matrix, sv: JuliaFloat64Vector, Vt: JuliaFloat64Matrix)
svd!(m)is the same assvd(m) but overwites a to save memory space.
- svd: JuliaFloat64Matrix -> Record(U: JuliaFloat64Matrix, sv: JuliaFloat64Vector, Vt: JuliaFloat64Matrix)
svd(m)computes the singular value decompositionSVDofmsuch thatSVD.U* diagonalMatrix(sv) *SVD.Vt=m.
- svdvals!: JuliaFloat64Matrix -> JuliaFloat64Vector
svdvals!(m)returns the singular values ofmbut overwritesmto save memory space.
- svdvals: JuliaFloat64Matrix -> JuliaFloat64Vector
svdvals(m)returns the singular values ofm.
- trace: JuliaFloat64Matrix -> JuliaFloat64
trace(m)computes the trace ofm.
- tril!: JuliaFloat64Matrix -> JuliaFloat64Matrix
tril!(m)overwritesmwith its upper triangular matrix counterpart. Returnsm.
- tril: JuliaFloat64Matrix -> JuliaFloat64Matrix
tril(m)returns the lower triangular matrix ofm
- triu!: JuliaFloat64Matrix -> JuliaFloat64Matrix
triu!(m)overwritesmwith its upper triangular matrix counterpart. Returnsm.
- triu: JuliaFloat64Matrix -> JuliaFloat64Matrix
triu(m)returns the upper triangular matrix ofm.