Skip to content

Commit 5b2e3e7

Browse files
sam0410simonbyrne
authored andcommittedDec 11, 2018
Adding rtol and atol for pinv and nullspace (#29998)
* Add code for rtol and atol * add tests * resolve comment * fix typo * fix tests * add news.md item * Not deprecated yet. * Tweak docs slightly * typo from diff [skip ci]
1 parent 217d330 commit 5b2e3e7

File tree

4 files changed

+46
-23
lines changed

4 files changed

+46
-23
lines changed
 

‎NEWS.md

+1
Original file line numberDiff line numberDiff line change
@@ -110,6 +110,7 @@ Standard library changes
110110
* `mul!`, `rmul!` and `lmul!` methods for `UniformScaling` ([#29506]).
111111
* `Symmetric` and `Hermitian` matrices now preserve the wrapper when scaled with a number ([#29469]).
112112
* Exponentiation operator `^` now supports raising an `Irrational` to an `AbstractMatrix` power ([#29782]).
113+
* Added keyword arguments `rtol`, `atol` to `pinv`, `nullspace` and `rank` ([#29998], [#29926]).
113114

114115
#### Random
115116
* `randperm` and `randcycle` now use the type of their argument to determine the element type of

‎stdlib/LinearAlgebra/src/dense.jl

+34-23
Original file line numberDiff line numberDiff line change
@@ -1237,19 +1237,21 @@ factorize(A::Transpose) = transpose(factorize(parent(A)))
12371237
## Moore-Penrose pseudoinverse
12381238

12391239
"""
1240-
pinv(M[, rtol::Real])
1240+
pinv(M; atol::Real=0, rtol::Real=atol>0 ? 0 : n*ϵ)
1241+
pinv(M, rtol::Real) = pinv(M; rtol=rtol) # to be deprecated in Julia 2.0
12411242
12421243
Computes the Moore-Penrose pseudoinverse.
12431244
12441245
For matrices `M` with floating point elements, it is convenient to compute
12451246
the pseudoinverse by inverting only singular values greater than
1246-
`rtol * maximum(svdvals(M))`.
1247+
`max(atol, rtol*σ₁)` where `σ₁` is the largest singular value of `M`.
12471248
1248-
The optimal choice of `rtol` varies both with the value of `M` and the intended application
1249-
of the pseudoinverse. The default value of `rtol` is
1250-
`eps(real(float(one(eltype(M)))))*minimum(size(M))`, which is essentially machine epsilon
1251-
for the real part of a matrix element multiplied by the larger matrix dimension. For
1252-
inverting dense ill-conditioned matrices in a least-squares sense,
1249+
The optimal choice of absolute (`atol`) and relative tolerance (`rtol`) varies
1250+
both with the value of `M` and the intended application of the pseudoinverse.
1251+
The default relative tolerance is `n*ϵ`, where `n` is the size of the smallest
1252+
dimension of `M`, and `ϵ` is the [`eps`](@ref) of the element type of `M`.
1253+
1254+
For inverting dense ill-conditioned matrices in a least-squares sense,
12531255
`rtol = sqrt(eps(real(float(one(eltype(M))))))` is recommended.
12541256
12551257
For more information, see [^issue8859], [^B96], [^S84], [^KY88].
@@ -1280,7 +1282,7 @@ julia> M * N
12801282
12811283
[^KY88]: Konstantinos Konstantinides and Kung Yao, "Statistical analysis of effective singular values in matrix rank determination", IEEE Transactions on Acoustics, Speech and Signal Processing, 36(5), 1988, 757-763. [doi:10.1109/29.1585](https://doi.org/10.1109/29.1585)
12821284
"""
1283-
function pinv(A::AbstractMatrix{T}, rtol::Real) where T
1285+
function pinv(A::AbstractMatrix{T}; atol::Real = 0.0, rtol::Real = (eps(real(float(one(T))))*min(size(A)...))*iszero(atol)) where T
12841286
m, n = size(A)
12851287
Tout = typeof(zero(T)/sqrt(one(T) + one(T)))
12861288
if m == 0 || n == 0
@@ -1289,9 +1291,10 @@ function pinv(A::AbstractMatrix{T}, rtol::Real) where T
12891291
if istril(A)
12901292
if istriu(A)
12911293
maxabsA = maximum(abs.(diag(A)))
1294+
tol = max(rtol*maxabsA, atol)
12921295
B = zeros(Tout, n, m)
12931296
for i = 1:min(m, n)
1294-
if abs(A[i,i]) > rtol*maxabsA
1297+
if abs(A[i,i]) > tol
12951298
Aii = inv(A[i,i])
12961299
if isfinite(Aii)
12971300
B[i,i] = Aii
@@ -1302,17 +1305,14 @@ function pinv(A::AbstractMatrix{T}, rtol::Real) where T
13021305
end
13031306
end
13041307
SVD = svd(A, full = false)
1308+
tol = max(rtol*maximum(SVD.S), atol)
13051309
Stype = eltype(SVD.S)
13061310
Sinv = zeros(Stype, length(SVD.S))
1307-
index = SVD.S .> rtol*maximum(SVD.S)
1311+
index = SVD.S .> tol
13081312
Sinv[index] = one(Stype) ./ SVD.S[index]
13091313
Sinv[findall(.!isfinite.(Sinv))] .= zero(Stype)
13101314
return SVD.Vt' * (Diagonal(Sinv) * SVD.U')
13111315
end
1312-
function pinv(A::AbstractMatrix{T}) where T
1313-
rtol = eps(real(float(one(T))))*min(size(A)...)
1314-
return pinv(A, rtol)
1315-
end
13161316
function pinv(x::Number)
13171317
xi = inv(x)
13181318
return ifelse(isfinite(xi), xi, zero(xi))
@@ -1321,13 +1321,16 @@ end
13211321
## Basis for null space
13221322

13231323
"""
1324-
nullspace(M[, rtol::Real])
1324+
nullspace(M; atol::Real=0, rtol::Rea=atol>0 ? 0 : n*ϵ)
1325+
nullspace(M, rtol::Real) = nullspace(M; rtol=rtol) # to be deprecated in Julia 2.0
13251326
13261327
Computes a basis for the nullspace of `M` by including the singular
1327-
vectors of A whose singular have magnitude are greater than `rtol*σ₁`,
1328-
where `σ₁` is `A`'s largest singular values. By default, the value of
1329-
`rtol` is the smallest dimension of `A` multiplied by the [`eps`](@ref)
1330-
of the [`eltype`](@ref) of `A`.
1328+
vectors of A whose singular have magnitude are greater than `max(atol, rtol*σ₁)`,
1329+
where `σ₁` is `M`'s largest singularvalue.
1330+
1331+
By default, the relative tolerance `rtol` is `n*ϵ`, where `n`
1332+
is the size of the smallest dimension of `M`, and `ϵ` is the [`eps`](@ref) of
1333+
the element type of `M`.
13311334
13321335
# Examples
13331336
```jldoctest
@@ -1343,21 +1346,29 @@ julia> nullspace(M)
13431346
0.0
13441347
1.0
13451348
1346-
julia> nullspace(M, 2)
1349+
julia> nullspace(M, rtol=3)
13471350
3×3 Array{Float64,2}:
13481351
0.0 1.0 0.0
13491352
1.0 0.0 0.0
13501353
0.0 0.0 1.0
1354+
1355+
julia> nullspace(M, atol=0.95)
1356+
3×1 Array{Float64,2}:
1357+
0.0
1358+
0.0
1359+
1.0
13511360
```
13521361
"""
1353-
function nullspace(A::AbstractMatrix, rtol::Real = min(size(A)...)*eps(real(float(one(eltype(A))))))
1362+
function nullspace(A::AbstractMatrix; atol::Real = 0.0, rtol::Real = (min(size(A)...)*eps(real(float(one(eltype(A))))))*iszero(atol))
13541363
m, n = size(A)
13551364
(m == 0 || n == 0) && return Matrix{eltype(A)}(I, n, n)
13561365
SVD = svd(A, full=true)
1357-
indstart = sum(s -> s .> SVD.S[1]*rtol, SVD.S) + 1
1366+
tol = max(atol, SVD.S[1]*rtol)
1367+
indstart = sum(s -> s .> tol, SVD.S) + 1
13581368
return copy(SVD.Vt[indstart:end,:]')
13591369
end
1360-
nullspace(a::AbstractVector, rtol::Real = min(size(a)...)*eps(real(float(one(eltype(a)))))) = nullspace(reshape(a, length(a), 1), rtol)
1370+
1371+
nullspace(A::AbstractVector; atol::Real = 0.0, rtol::Real = (min(size(A)...)*eps(real(float(one(eltype(A))))))*iszero(atol)) = nullspace(reshape(A, length(A), 1), rtol= rtol, atol= atol)
13611372

13621373
"""
13631374
cond(M, p::Real=2)

‎stdlib/LinearAlgebra/src/deprecated.jl

+3
Original file line numberDiff line numberDiff line change
@@ -2,3 +2,6 @@
22

33
# To be deprecated in 2.0
44
rank(A::AbstractMatrix, tol::Real) = rank(A,rtol=tol)
5+
nullspace(A::AbstractVector, tol::Real) = nullspace(reshape(A, length(A), 1), rtol= tol)
6+
nullspace(A::AbstractMatrix, tol::Real) = nullspace(A, rtol=tol)
7+
pinv(A::AbstractMatrix{T}, tol::Real) where T = pinv(A, rtol=tol)

‎stdlib/LinearAlgebra/test/dense.jl

+8
Original file line numberDiff line numberDiff line change
@@ -72,6 +72,8 @@ bimg = randn(n,2)/2
7272
@test norm(a[:,1:n1]'a15null,Inf) zero(eltya) atol=300ε
7373
@test norm(a15null'a[:,1:n1],Inf) zero(eltya) atol=400ε
7474
@test size(nullspace(b), 2) == 0
75+
@test size(nullspace(b, rtol=0.001), 2) == 0
76+
@test size(nullspace(b, atol=100*εb), 2) == 0
7577
@test size(nullspace(b, 100*εb), 2) == 0
7678
@test nullspace(zeros(eltya,n)) == Matrix(I, 1, 1)
7779
@test nullspace(zeros(eltya,n), 0.1) == Matrix(I, 1, 1)
@@ -82,6 +84,12 @@ bimg = randn(n,2)/2
8284
end
8385
end # for eltyb
8486

87+
@testset "Test pinv (rtol, atol)" begin
88+
M = [1 0 0; 0 1 0; 0 0 0]
89+
@test pinv(M,atol=1)== zeros(3,3)
90+
@test pinv(M,rtol=0.5)== M
91+
end
92+
8593
for (a, a2) in ((copy(ainit), copy(ainit2)), (view(ainit, 1:n, 1:n), view(ainit2, 1:n, 1:n)))
8694
@testset "Test pinv" begin
8795
pinva15 = pinv(a[:,1:n1])

0 commit comments

Comments
 (0)
Please sign in to comment.