Skip to content

Commit 69e407b

Browse files
committed
Based on discussions in #2212.
Minimal invasive changes to the Factorization objects for 0.1. All existing functionality is retained, with only a few renamings for clarity. Replaces the *d routines with *fact routines: lud ->lufact chold -> cholfact cholpd-> cholpfact qrd ->qrfact qrpd -> qrpfact Replaces * and Ac_mul_B on QRDense objects: qmulQR performs Q*A qTmulQR performs Q'*A Updates to exports, documentation, tests, and deprecations. This is a self-contained commit that can be reverted if necessary close #2062
1 parent ce4c54d commit 69e407b

File tree

5 files changed

+100
-85
lines changed

5 files changed

+100
-85
lines changed

base/deprecated.jl

+5
Original file line numberDiff line numberDiff line change
@@ -86,6 +86,11 @@ end
8686
@deprecate cov_pearson cov
8787
@deprecate areduce reducedim
8888
@deprecate tmpnam tempname
89+
@deprecate lud lufact
90+
@deprecate chold cholfact
91+
@deprecate cholpd cholpfact
92+
@deprecate qrd qrfact
93+
@deprecate qrpd qrpfact
8994

9095
export grow!
9196
function grow!(a, d)

base/exports.jl

+12-10
Original file line numberDiff line numberDiff line change
@@ -517,10 +517,10 @@ export
517517

518518
# linear algebra
519519
chol,
520-
chold,
521-
chold!,
522-
cholpd,
523-
cholpd!,
520+
cholfact,
521+
cholfact!,
522+
cholpfact,
523+
cholpfact!,
524524
cond,
525525
cross,
526526
ctranspose,
@@ -549,18 +549,20 @@ export
549549
ldltd,
550550
linreg,
551551
lu,
552-
lud,
553-
lud!,
552+
lufact,
553+
lufact!,
554554
norm,
555555
normfro,
556556
null,
557557
pinv,
558558
qr,
559-
qrd!,
560-
qrd,
559+
qrfact!,
560+
qrfact,
561561
qrp,
562-
qrpd!,
563-
qrpd,
562+
qrpfact!,
563+
qrpfact,
564+
qmulQR,
565+
qTmulQR,
564566
randsym,
565567
rank,
566568
rref,

base/linalg_dense.jl

+45-45
Original file line numberDiff line numberDiff line change
@@ -425,14 +425,14 @@ function inv(C::CholeskyDense)
425425
end
426426

427427
## Should these functions check that the matrix is Hermitian?
428-
chold!{T<:BlasFloat}(A::Matrix{T}, UL::BlasChar) = CholeskyDense{T}(A, UL)
429-
chold!{T<:BlasFloat}(A::Matrix{T}) = chold!(A, 'U')
430-
chold{T<:BlasFloat}(A::Matrix{T}, UL::BlasChar) = chold!(copy(A), UL)
431-
chold{T<:Number}(A::Matrix{T}, UL::BlasChar) = chold(float64(A), UL)
432-
chold{T<:Number}(A::Matrix{T}) = chold(A, 'U')
428+
cholfact!{T<:BlasFloat}(A::Matrix{T}, UL::BlasChar) = CholeskyDense{T}(A, UL)
429+
cholfact!{T<:BlasFloat}(A::Matrix{T}) = cholfact!(A, 'U')
430+
cholfact{T<:BlasFloat}(A::Matrix{T}, UL::BlasChar) = cholfact!(copy(A), UL)
431+
cholfact{T<:Number}(A::Matrix{T}, UL::BlasChar) = cholfact(float64(A), UL)
432+
cholfact{T<:Number}(A::Matrix{T}) = cholfact(A, 'U')
433433

434434
## Matlab (and R) compatible
435-
chol(A::Matrix, UL::BlasChar) = factors(chold(A, UL))
435+
chol(A::Matrix, UL::BlasChar) = factors(cholfact(A, UL))
436436
chol(A::Matrix) = chol(A, 'U')
437437
chol(x::Number) = imag(x) == 0 && real(x) > 0 ? sqrt(x) : error("Argument not positive-definite")
438438

@@ -489,18 +489,18 @@ function inv(C::CholeskyPivotedDense)
489489
end
490490

491491
## Should these functions check that the matrix is Hermitian?
492-
cholpd!{T<:BlasFloat}(A::Matrix{T}, UL::BlasChar, tol::Real) = CholeskyPivotedDense{T}(A, UL, tol)
493-
cholpd!{T<:BlasFloat}(A::Matrix{T}, UL::BlasChar) = cholpd!(A, UL, -1.)
494-
cholpd!{T<:BlasFloat}(A::Matrix{T}, tol::Real) = cholpd!(A, 'U', tol)
495-
cholpd!{T<:BlasFloat}(A::Matrix{T}) = cholpd!(A, 'U', -1.)
496-
cholpd{T<:Number}(A::Matrix{T}, UL::BlasChar, tol::Real) = cholpd(float64(A), UL, tol)
497-
cholpd{T<:Number}(A::Matrix{T}, UL::BlasChar) = cholpd(float64(A), UL, -1.)
498-
cholpd{T<:Number}(A::Matrix{T}, tol::Real) = cholpd(float64(A), true, tol)
499-
cholpd{T<:Number}(A::Matrix{T}) = cholpd(float64(A), 'U', -1.)
500-
cholpd{T<:BlasFloat}(A::Matrix{T}, UL::BlasChar, tol::Real) = cholpd!(copy(A), UL, tol)
501-
cholpd{T<:BlasFloat}(A::Matrix{T}, UL::BlasChar) = cholpd!(copy(A), UL, -1.)
502-
cholpd{T<:BlasFloat}(A::Matrix{T}, tol::Real) = cholpd!(copy(A), 'U', tol)
503-
cholpd{T<:BlasFloat}(A::Matrix{T}) = cholpd!(copy(A), 'U', -1.)
492+
cholpfact!{T<:BlasFloat}(A::Matrix{T}, UL::BlasChar, tol::Real) = CholeskyPivotedDense{T}(A, UL, tol)
493+
cholpfact!{T<:BlasFloat}(A::Matrix{T}, UL::BlasChar) = cholpfact!(A, UL, -1.)
494+
cholpfact!{T<:BlasFloat}(A::Matrix{T}, tol::Real) = cholpfact!(A, 'U', tol)
495+
cholpfact!{T<:BlasFloat}(A::Matrix{T}) = cholpfact!(A, 'U', -1.)
496+
cholpfact{T<:Number}(A::Matrix{T}, UL::BlasChar, tol::Real) = cholpfact(float64(A), UL, tol)
497+
cholpfact{T<:Number}(A::Matrix{T}, UL::BlasChar) = cholpfact(float64(A), UL, -1.)
498+
cholpfact{T<:Number}(A::Matrix{T}, tol::Real) = cholpfact(float64(A), true, tol)
499+
cholpfact{T<:Number}(A::Matrix{T}) = cholpfact(float64(A), 'U', -1.)
500+
cholpfact{T<:BlasFloat}(A::Matrix{T}, UL::BlasChar, tol::Real) = cholpfact!(copy(A), UL, tol)
501+
cholpfact{T<:BlasFloat}(A::Matrix{T}, UL::BlasChar) = cholpfact!(copy(A), UL, -1.)
502+
cholpfact{T<:BlasFloat}(A::Matrix{T}, tol::Real) = cholpfact!(copy(A), 'U', tol)
503+
cholpfact{T<:BlasFloat}(A::Matrix{T}) = cholpfact!(copy(A), 'U', -1.)
504504

505505
type LUDense{T} <: Factorization{T}
506506
lu::Matrix{T}
@@ -530,16 +530,16 @@ function factors{T}(lu::LUDense{T})
530530
L, U, P
531531
end
532532

533-
function lud!{T<:BlasFloat}(A::Matrix{T})
533+
function lufact!{T<:BlasFloat}(A::Matrix{T})
534534
lu, ipiv, info = LAPACK.getrf!(A)
535535
LUDense{T}(lu, ipiv, info)
536536
end
537537

538-
lud{T<:BlasFloat}(A::Matrix{T}) = lud!(copy(A))
539-
lud{T<:Number}(A::Matrix{T}) = lud(float64(A))
538+
lufact{T<:BlasFloat}(A::Matrix{T}) = lufact!(copy(A))
539+
lufact{T<:Number}(A::Matrix{T}) = lufact(float64(A))
540540

541541
## Matlab-compatible
542-
lu{T<:Number}(A::Matrix{T}) = factors(lud(A))
542+
lu{T<:Number}(A::Matrix{T}) = factors(lufact(A))
543543
lu(x::Number) = (one(x), x, [1])
544544

545545
function det{T}(lu::LUDense{T})
@@ -571,31 +571,31 @@ end
571571
size(A::QRDense) = size(A.hh)
572572
size(A::QRDense,n) = size(A.hh,n)
573573

574-
qrd!{T<:BlasFloat}(A::StridedMatrix{T}) = QRDense{T}(LAPACK.geqrf!(A)...)
575-
qrd{T<:BlasFloat}(A::StridedMatrix{T}) = qrd!(copy(A))
576-
qrd{T<:Real}(A::StridedMatrix{T}) = qrd(float64(A))
574+
qrfact!{T<:BlasFloat}(A::StridedMatrix{T}) = QRDense{T}(LAPACK.geqrf!(A)...)
575+
qrfact{T<:BlasFloat}(A::StridedMatrix{T}) = qrfact!(copy(A))
576+
qrfact{T<:Real}(A::StridedMatrix{T}) = qrfact(float64(A))
577577

578-
function factors{T<:BlasFloat}(qrd::QRDense{T})
579-
aa = copy(qrd.hh)
578+
function factors{T<:BlasFloat}(qrfact::QRDense{T})
579+
aa = copy(qrfact.hh)
580580
R = triu(aa[1:min(size(aa)),:]) # must be *before* call to orgqr!
581-
LAPACK.orgqr!(aa, qrd.tau, size(aa,2)), R
581+
LAPACK.orgqr!(aa, qrfact.tau, size(aa,2)), R
582582
end
583583

584-
qr{T<:Number}(x::StridedMatrix{T}) = factors(qrd(x))
584+
qr{T<:Number}(x::StridedMatrix{T}) = factors(qrfact(x))
585585
qr(x::Number) = (one(x), x)
586586

587587
## Multiplication by Q from the QR decomposition
588-
(*){T<:BlasFloat}(A::QRDense{T}, B::StridedVecOrMat{T}) =
588+
qmulQR{T<:BlasFloat}(A::QRDense{T}, B::StridedVecOrMat{T}) =
589589
LAPACK.ormqr!('L', 'N', A.hh, size(A.hh,2), A.tau, copy(B))
590590

591591
## Multiplication by Q' from the QR decomposition
592-
Ac_mul_B{T<:BlasFloat}(A::QRDense{T}, B::StridedVecOrMat{T}) =
592+
qTmulQR{T<:BlasFloat}(A::QRDense{T}, B::StridedVecOrMat{T}) =
593593
LAPACK.ormqr!('L', iscomplex(A.tau)?'C':'T', A.hh, size(A.hh,2), A.tau, copy(B))
594594

595595
## Least squares solution. Should be more careful about cases with m < n
596596
function (\){T<:BlasFloat}(A::QRDense{T}, B::StridedVecOrMat{T})
597597
n = length(A.tau)
598-
ans, info = LAPACK.trtrs!('U','N','N',A.hh[1:n,:],(A'*B)[1:n,:])
598+
ans, info = LAPACK.trtrs!('U','N','N',A.hh[1:n,:],(qTmulQR(A,B))[1:n,:])
599599
if info > 0; throw(LAPACK.SingularException(info)); end
600600
return ans
601601
end
@@ -615,28 +615,28 @@ end
615615
size(x::QRPivotedDense) = size(x.hh)
616616
size(x::QRPivotedDense,d) = size(x.hh,d)
617617
## Multiplication by Q from the QR decomposition
618-
(*){T<:BlasFloat}(A::QRPivotedDense{T}, B::StridedVecOrMat{T}) =
618+
qmulQR{T<:BlasFloat}(A::QRPivotedDense{T}, B::StridedVecOrMat{T}) =
619619
LAPACK.ormqr!('L', 'N', A.hh, size(A,2), A.tau, copy(B))
620620
## Multiplication by Q' from the QR decomposition
621-
Ac_mul_B{T<:BlasFloat}(A::QRPivotedDense{T}, B::StridedVecOrMat{T}) =
621+
qTmulQR{T<:BlasFloat}(A::QRPivotedDense{T}, B::StridedVecOrMat{T}) =
622622
LAPACK.ormqr!('L', iscomplex(A.tau)?'C':'T', A.hh, size(A,2), A.tau, copy(B))
623623

624-
qrpd!{T<:BlasFloat}(A::StridedMatrix{T}) = QRPivotedDense{T}(LAPACK.geqp3!(A)...)
625-
qrpd{T<:BlasFloat}(A::StridedMatrix{T}) = qrpd!(copy(A))
626-
qrpd{T<:Real}(x::StridedMatrix{T}) = qrpd(float64(x))
624+
qrpfact!{T<:BlasFloat}(A::StridedMatrix{T}) = QRPivotedDense{T}(LAPACK.geqp3!(A)...)
625+
qrpfact{T<:BlasFloat}(A::StridedMatrix{T}) = qrpfact!(copy(A))
626+
qrpfact{T<:Real}(x::StridedMatrix{T}) = qrpfact(float64(x))
627627

628628
function factors{T<:BlasFloat}(x::QRPivotedDense{T})
629629
aa = copy(x.hh)
630630
R = triu(aa[1:min(size(aa)),:])
631631
LAPACK.orgqr!(aa, x.tau, size(aa,2)), R, x.jpvt
632632
end
633633

634-
qrp{T<:BlasFloat}(x::StridedMatrix{T}) = factors(qrpd(x))
634+
qrp{T<:BlasFloat}(x::StridedMatrix{T}) = factors(qrpfact(x))
635635
qrp{T<:Real}(x::StridedMatrix{T}) = qrp(float64(x))
636636

637637
function (\){T<:BlasFloat}(A::QRPivotedDense{T}, B::StridedVecOrMat{T})
638638
n = length(A.tau)
639-
x, info = LAPACK.trtrs!('U','N','N',A.hh[1:n,:],(A'*B)[1:n,:])
639+
x, info = LAPACK.trtrs!('U','N','N',A.hh[1:n,:],(qTmulQR(A,B))[1:n,:])
640640
if info > 0; throw(LAPACK.SingularException(info)); end
641641
isa(B, Vector) ? x[invperm(A.jpvt)] : x[:,invperm(A.jpvt)]
642642
end
@@ -674,7 +674,7 @@ function det(A::Matrix)
674674
m, n = size(A)
675675
if m != n; error("det only defined for square matrices"); end
676676
if istriu(A) | istril(A); return prod(diag(A)); end
677-
return det(lud(A))
677+
return det(lufact(A))
678678
end
679679
det(x::Number) = x
680680

@@ -921,7 +921,7 @@ function cond(A::StridedMatrix, p)
921921
elseif p == 1 || p == Inf
922922
m, n = size(A)
923923
if m != n; error("Use 2-norm for non-square matrices"); end
924-
cnd = 1 / LAPACK.gecon!(p == 1 ? '1' : 'I', lud(A).lu, norm(A, p))
924+
cnd = 1 / LAPACK.gecon!(p == 1 ? '1' : 'I', lufact(A).lu, norm(A, p))
925925
else
926926
error("Norm type must be 1, 2 or Inf")
927927
end
@@ -1217,16 +1217,16 @@ end
12171217

12181218
#show(io, lu::LUTridiagonal) = print(io, "LU decomposition of ", summary(lu.lu))
12191219

1220-
lud!{T}(A::Tridiagonal{T}) = LUTridiagonal{T}(LAPACK.gttrf!(A.dl,A.d,A.du)...)
1221-
lud{T}(A::Tridiagonal{T}) = LUTridiagonal{T}(LAPACK.gttrf!(copy(A.dl),copy(A.d),copy(A.du))...)
1222-
lu(A::Tridiagonal) = factors(lud(A))
1220+
lufact!{T}(A::Tridiagonal{T}) = LUTridiagonal{T}(LAPACK.gttrf!(A.dl,A.d,A.du)...)
1221+
lufact{T}(A::Tridiagonal{T}) = LUTridiagonal{T}(LAPACK.gttrf!(copy(A.dl),copy(A.d),copy(A.du))...)
1222+
lu(A::Tridiagonal) = factors(lufact(A))
12231223

12241224
function det{T}(lu::LUTridiagonal{T})
12251225
n = length(lu.d)
12261226
prod(lu.d) * (bool(sum(lu.ipiv .!= 1:n) % 2) ? -one(T) : one(T))
12271227
end
12281228

1229-
det(A::Tridiagonal) = det(lud(A))
1229+
det(A::Tridiagonal) = det(lufact(A))
12301230

12311231
(\){T<:BlasFloat}(lu::LUTridiagonal{T}, B::StridedVecOrMat{T}) =
12321232
LAPACK.gttrs!('N', lu.dl, lu.d, lu.du, lu.du2, lu.ipiv, copy(B))

doc/stdlib/base.rst

+29-21
Original file line numberDiff line numberDiff line change
@@ -1860,57 +1860,65 @@ Linear algebra functions in Julia are largely implemented by calling functions f
18601860

18611861
Compute the LU factorization of ``A``, such that ``A[P,:] = L*U``.
18621862

1863-
.. function:: lud(A) -> LUDense
1863+
.. function:: lufact(A) -> LUDense
18641864

1865-
Compute the LU factorization of ``A`` and return a ``LUDense`` object. ``factors(lud(A))`` returns the triangular matrices containing the factorization. The following functions are available for ``LUDense`` objects: ``size``, ``factors``, ``\``, ``inv``, ``det``.
1865+
Compute the LU factorization of ``A`` and return a ``LUDense`` object. ``factors(lufact(A))`` returns the triangular matrices containing the factorization. The following functions are available for ``LUDense`` objects: ``size``, ``factors``, ``\``, ``inv``, ``det``.
18661866

1867-
.. function:: lud!(A) -> LUDense
1867+
.. function:: lufact!(A) -> LUDense
18681868

1869-
``lud!`` is the same as ``lud`` but saves space by overwriting the input A, instead of creating a copy.
1869+
``lufact!`` is the same as ``lufact`` but saves space by overwriting the input A, instead of creating a copy.
18701870

18711871
.. function:: chol(A, [LU]) -> F
18721872

18731873
Compute Cholesky factorization of a symmetric positive-definite matrix ``A`` and return the matrix ``F``. If ``LU`` is ``L`` (Lower), ``A = L*L'``. If ``LU`` is ``U`` (Upper), ``A = R'*R``.
18741874

1875-
.. function:: chold(A, [LU]) -> CholeskyDense
1875+
.. function:: cholfact(A, [LU]) -> CholeskyDense
18761876

1877-
Compute the Cholesky factorization of a symmetric positive-definite matrix ``A`` and return a ``CholeskyDense`` object. ``LU`` may be 'L' for using the lower part or 'U' for the upper part. The default is to use 'U'. ``factors(chold(A))`` returns the triangular matrix containing the factorization. The following functions are available for ``CholeskyDense`` objects: ``size``, ``factors``, ``\``, ``inv``, ``det``. A ``LAPACK.PosDefException`` error is thrown in case the matrix is not positive definite.
1877+
Compute the Cholesky factorization of a symmetric positive-definite matrix ``A`` and return a ``CholeskyDense`` object. ``LU`` may be 'L' for using the lower part or 'U' for the upper part. The default is to use 'U'. ``factors(cholfact(A))`` returns the triangular matrix containing the factorization. The following functions are available for ``CholeskyDense`` objects: ``size``, ``factors``, ``\``, ``inv``, ``det``. A ``LAPACK.PosDefException`` error is thrown in case the matrix is not positive definite.
18781878

1879-
.. function: chold!(A, [LU]) -> CholeskyDense
1879+
.. function: cholfact!(A, [LU]) -> CholeskyDense
18801880
1881-
``chold!`` is the same as ``chold`` but saves space by overwriting the input A, instead of creating a copy.
1881+
``cholfact!`` is the same as ``cholfact`` but saves space by overwriting the input A, instead of creating a copy.
18821882
1883-
.. function:: cholpd(A, [LU]) -> CholeskyPivotedDense
1883+
.. function:: cholpfact(A, [LU]) -> CholeskyPivotedDense
18841884

1885-
Compute the pivoted Cholesky factorization of a symmetric positive semi-definite matrix ``A`` and return a ``CholeskyDensePivoted`` object. ``LU`` may be 'L' for using the lower part or 'U' for the upper part. The default is to use 'U'. ``factors(cholpd(A))`` returns the triangular matrix containing the factorization. The following functions are available for ``CholeskyDensePivoted`` objects: ``size``, ``factors``, ``\``, ``inv``, ``det``. A ``LAPACK.RankDeficientException`` error is thrown in case the matrix is rank deficient.
1885+
Compute the pivoted Cholesky factorization of a symmetric positive semi-definite matrix ``A`` and return a ``CholeskyDensePivoted`` object. ``LU`` may be 'L' for using the lower part or 'U' for the upper part. The default is to use 'U'. ``factors(cholpfact(A))`` returns the triangular matrix containing the factorization. The following functions are available for ``CholeskyDensePivoted`` objects: ``size``, ``factors``, ``\``, ``inv``, ``det``. A ``LAPACK.RankDeficientException`` error is thrown in case the matrix is rank deficient.
18861886

1887-
.. function:: cholpd!(A, [LU]) -> CholeskyPivotedDense
1887+
.. function:: cholpfact!(A, [LU]) -> CholeskyPivotedDense
18881888

1889-
``cholpd!`` is the same as ``cholpd`` but saves space by overwriting the input A, instead of creating a copy.
1889+
``cholpfact!`` is the same as ``cholpfact`` but saves space by overwriting the input A, instead of creating a copy.
18901890

18911891
.. function:: qr(A) -> Q, R
18921892

18931893
Compute the QR factorization of ``A`` such that ``A = Q*R``. Also see ``qrd``.
18941894

1895-
.. function:: qrd(A)
1895+
.. function:: qrfact(A)
18961896

1897-
Compute the QR factorization of ``A`` and return a ``QRDense`` object. ``factors(qrd(A))`` returns ``Q`` and ``R``. The following functions are available for ``QRDense`` objects: ``size``, ``factors``, ``*``, ``Ac_mul_B``, ``\``.
1897+
Compute the QR factorization of ``A`` and return a ``QRDense`` object. ``factors(qrfact(A))`` returns ``Q`` and ``R``. The following functions are available for ``QRDense`` objects: ``size``, ``factors``, ``qmulQR``, ``qTmulQR``, ``\``.
18981898

1899-
.. function:: qrd!(A)
1899+
.. function:: qrfact!(A)
19001900

1901-
``qrd!`` is the same as ``qrd`` but saves space by overwriting the input A, instead of creating a copy.
1901+
``qrfact!`` is the same as ``qrfact`` but saves space by overwriting the input A, instead of creating a copy.
19021902

19031903
.. function:: qrp(A) -> Q, R, P
19041904

1905-
Compute the QR factorization of ``A`` with pivoting, such that ``A*I[:,P] = Q*R``, where ``I`` is the identity matrix. Also see ``qrpd``.
1905+
Compute the QR factorization of ``A`` with pivoting, such that ``A*I[:,P] = Q*R``, where ``I`` is the identity matrix. Also see ``qrpfact``.
19061906

1907-
.. function:: qrpd(A) -> QRPivotedDense
1907+
.. function:: qrpfact(A) -> QRPivotedDense
19081908

1909-
Compute the QR factorization of ``A`` with pivoting and return a ``QRDensePivoted`` object. ``factors(qrpd(A))`` returns ``Q`` and ``R``. The following functions are available for ``QRDensePivoted`` objects: ``size``, ``factors``, ``*``, ``Ac_mul_B``, ``\``.
1909+
Compute the QR factorization of ``A`` with pivoting and return a ``QRDensePivoted`` object. ``factors(qrpfact(A))`` returns ``Q`` and ``R``. The following functions are available for ``QRDensePivoted`` objects: ``size``, ``factors``, ``qmulQR``, ``qTmulQR``, ``\``.
19101910

1911-
.. function:: qrpd!(A) -> QRPivotedDense
1911+
.. function:: qrpfact!(A) -> QRPivotedDense
19121912

1913-
``qrpd!`` is the same as ``qrpd`` but saves space by overwriting the input A, instead of creating a copy.
1913+
``qrpfact!`` is the same as ``qrpfact`` but saves space by overwriting the input A, instead of creating a copy.
1914+
1915+
.. function:: qmulQR(QR, A)
1916+
1917+
Perform Q*A efficiently, where Q is a an orthogonal matrix defined as the product of k elementary reflectors from the QR decomposition.
1918+
1919+
.. function:: qTmulQR(QR, A)
1920+
1921+
Perform Q'*A efficiently, where Q is a an orthogonal matrix defined as the product of k elementary reflectors from the QR decomposition.
19141922
19151923
.. function:: eig(A) -> D, V
19161924

test/linalg.jl

+9-9
Original file line numberDiff line numberDiff line change
@@ -8,18 +8,18 @@ for elty in (Float32, Float64, Complex64, Complex128)
88
apd = a'*a # symmetric positive-definite
99
b = convert(Vector{elty}, b)
1010

11-
capd = chold(apd) # upper Cholesky factor
11+
capd = cholfact(apd) # upper Cholesky factor
1212
r = factors(capd)
1313
@assert_approx_eq r'*r apd
1414
@assert_approx_eq b apd * (capd\b)
1515
@assert_approx_eq apd * inv(capd) eye(elty, n)
1616
@assert_approx_eq a*(capd\(a'*b)) b # least squares soln for square a
1717
@assert_approx_eq det(capd) det(apd)
1818

19-
l = factors(chold(apd, 'L')) # lower Cholesky factor
19+
l = factors(cholfact(apd, 'L')) # lower Cholesky factor
2020
@assert_approx_eq l*l' apd
2121

22-
cpapd = cholpd(apd) # pivoted Choleksy decomposition
22+
cpapd = cholpfact(apd) # pivoted Choleksy decomposition
2323
@test rank(cpapd) == n
2424
@test all(diff(diag(real(cpapd.LR))).<=0.) # diagonal should be non-increasing
2525
@assert_approx_eq b apd * (cpapd\b)
@@ -32,7 +32,7 @@ for elty in (Float32, Float64, Complex64, Complex128)
3232
@assert_approx_eq inv(bc2) * apd eye(elty, n)
3333
@assert_approx_eq apd * (bc2\b) b
3434

35-
lua = lud(a) # LU decomposition
35+
lua = lufact(a) # LU decomposition
3636
l,u,p = lu(a)
3737
L,U,P = factors(lua)
3838
@test l == L && u == U && p == P
@@ -41,18 +41,18 @@ for elty in (Float32, Float64, Complex64, Complex128)
4141
@assert_approx_eq a * inv(lua) eye(elty, n)
4242
@assert_approx_eq a*(lua\b) b
4343

44-
qra = qrd(a) # QR decomposition
44+
qra = qrfact(a) # QR decomposition
4545
q,r = factors(qra)
4646
@assert_approx_eq q'*q eye(elty, n)
4747
@assert_approx_eq q*q' eye(elty, n)
4848
Q,R = qr(a)
4949
@test q == Q && r == R
5050
@assert_approx_eq q*r a
51-
@assert_approx_eq qra*b Q*b
52-
@assert_approx_eq qra'*b Q'*b
51+
@assert_approx_eq qmulQR(qra,b) Q*b
52+
@assert_approx_eq qTmulQR(qra,b) Q'*b
5353
@assert_approx_eq a*(qra\b) b
5454

55-
qrpa = qrpd(a) # pivoted QR decomposition
55+
qrpa = qrpfact(a) # pivoted QR decomposition
5656
q,r,p = factors(qrpa)
5757
@assert_approx_eq q'*q eye(elty, n)
5858
@assert_approx_eq q*q' eye(elty, n)
@@ -294,7 +294,7 @@ for elty in (Float32, Float64, Complex64, Complex128)
294294
@assert_approx_eq solve(T,v) invFv
295295
B = convert(Matrix{elty}, B)
296296
@assert_approx_eq solve(T, B) F\B
297-
Tlu = lud(T)
297+
Tlu = lufact(T)
298298
x = Tlu\v
299299
@assert_approx_eq x invFv
300300
@assert_approx_eq det(T) det(F)

0 commit comments

Comments
 (0)