Skip to content

Commit bbf5584

Browse files
stevengjandreasnoack
authored andcommitted
isapprox: test max(atol,rtol*...) rather than atol+rtol*... (#22742)
* isapprox now tests max(atol,rtol*...) rather than atol+rtol*... * move news to breaking section * in isapprox, make default rtol=0 if atol>0 is specified * consistent isapprox for UniformScaling * whoops, remove old code * deprecate 2-arg rtoldefault
1 parent 55ab3e8 commit bbf5584

File tree

6 files changed

+32
-16
lines changed

6 files changed

+32
-16
lines changed

NEWS.md

+4
Original file line numberDiff line numberDiff line change
@@ -120,6 +120,10 @@ This section lists changes that do not have deprecation warnings.
120120
`Bidiagonal{T,V<:AbstractVector{T}}` and `SymTridiagonal{T,V<:AbstractVector{T}}`
121121
respectively ([#22718], [#22925], [#23035]).
122122

123+
* `isapprox(x,y)` now tests `norm(x-y) <= max(atol, rtol*max(norm(x), norm(y)))`
124+
rather than `norm(x-y) <= atol + ...`, and `rtol` defaults to zero
125+
if an `atol > 0` is specified ([#22742]).
126+
123127
* Spaces are no longer allowed between `@` and the name of a macro in a macro call ([#22868]).
124128

125129
* Juxtaposition of a non-literal with a macro call (`x@macro`) is no longer valid syntax ([#22868]).

base/deprecated.jl

+2
Original file line numberDiff line numberDiff line change
@@ -1661,6 +1661,8 @@ end
16611661
@deprecate num2hex(x::Union{Float16,Float32,Float64}) hex(reintepret(Unsigned, x), sizeof(x)*2)
16621662
@deprecate num2hex(n::Integer) hex(n, sizeof(n)*2)
16631663

1664+
# PR #22742: change in isapprox semantics
1665+
@deprecate rtoldefault(x,y) rtoldefault(x,y,0) false
16641666

16651667
# END 0.7 deprecations
16661668

base/floatfuncs.jl

+12-8
Original file line numberDiff line numberDiff line change
@@ -189,15 +189,16 @@ end
189189

190190
# isapprox: approximate equality of numbers
191191
"""
192-
isapprox(x, y; rtol::Real=sqrt(eps), atol::Real=0, nans::Bool=false, norm::Function)
192+
isapprox(x, y; rtol::Real=atol>0 ? √eps : 0, atol::Real=0, nans::Bool=false, norm::Function)
193193
194-
Inexact equality comparison: `true` if `norm(x-y) <= atol + rtol*max(norm(x), norm(y))`. The
194+
Inexact equality comparison: `true` if `norm(x-y) <= max(atol, rtol*max(norm(x), norm(y)))`. The
195195
default `atol` is zero and the default `rtol` depends on the types of `x` and `y`. The keyword
196196
argument `nans` determines whether or not NaN values are considered equal (defaults to false).
197197
198-
For real or complex floating-point values, `rtol` defaults to
199-
`sqrt(eps(typeof(real(x-y))))`. This corresponds to requiring equality of about half of the
200-
significand digits. For other types, `rtol` defaults to zero.
198+
For real or complex floating-point values, if an `atol > 0` is not specified, `rtol` defaults to
199+
the square root of [`eps`](@ref) of the type of `x` or `y`, whichever is bigger (least precise).
200+
This corresponds to requiring equality of about half of the significand digits. Otherwise,
201+
e.g. for integer arguments or if an `atol > 0` is supplied, `rtol` defaults to zero.
201202
202203
`x` and `y` may also be arrays of numbers, in which case `norm` defaults to `vecnorm` but
203204
may be changed by passing a `norm::Function` keyword argument. (For numbers, `norm` is the
@@ -220,8 +221,8 @@ julia> isapprox([10.0^9, 1.0], [10.0^9, 2.0])
220221
true
221222
```
222223
"""
223-
function isapprox(x::Number, y::Number; rtol::Real=rtoldefault(x,y), atol::Real=0, nans::Bool=false)
224-
x == y || (isfinite(x) && isfinite(y) && abs(x-y) <= atol + rtol*max(abs(x), abs(y))) || (nans && isnan(x) && isnan(y))
224+
function isapprox(x::Number, y::Number; atol::Real=0, rtol::Real=rtoldefault(x,y,atol), nans::Bool=false)
225+
x == y || (isfinite(x) && isfinite(y) && abs(x-y) <= max(atol, rtol*max(abs(x), abs(y)))) || (nans && isnan(x) && isnan(y))
225226
end
226227

227228
const = isapprox
@@ -230,7 +231,10 @@ const ≈ = isapprox
230231
# default tolerance arguments
231232
rtoldefault(::Type{T}) where {T<:AbstractFloat} = sqrt(eps(T))
232233
rtoldefault(::Type{<:Real}) = 0
233-
rtoldefault(x::Union{T,Type{T}}, y::Union{S,Type{S}}) where {T<:Number,S<:Number} = max(rtoldefault(real(T)),rtoldefault(real(S)))
234+
function rtoldefault(x::Union{T,Type{T}}, y::Union{S,Type{S}}, atol::Real) where {T<:Number,S<:Number}
235+
rtol = max(rtoldefault(real(T)),rtoldefault(real(S)))
236+
return atol > 0 ? zero(rtol) : rtol
237+
end
234238

235239
# fused multiply-add
236240
fma_libm(x::Float32, y::Float32, z::Float32) =

base/linalg/generic.jl

+4-3
Original file line numberDiff line numberDiff line change
@@ -1296,11 +1296,12 @@ promote_leaf_eltypes(x::Union{AbstractArray,Tuple}) = mapreduce(promote_leaf_elt
12961296
# Supports nested arrays; e.g., for `a = [[1,2, [3,4]], 5.0, [6im, [7.0, 8.0]]]`
12971297
# `a ≈ a` is `true`.
12981298
function isapprox(x::AbstractArray, y::AbstractArray;
1299-
rtol::Real=Base.rtoldefault(promote_leaf_eltypes(x),promote_leaf_eltypes(y)),
1300-
atol::Real=0, nans::Bool=false, norm::Function=vecnorm)
1299+
atol::Real=0,
1300+
rtol::Real=Base.rtoldefault(promote_leaf_eltypes(x),promote_leaf_eltypes(y),atol),
1301+
nans::Bool=false, norm::Function=vecnorm)
13011302
d = norm(x - y)
13021303
if isfinite(d)
1303-
return d <= atol + rtol*max(norm(x), norm(y))
1304+
return d <= max(atol, rtol*max(norm(x), norm(y)))
13041305
else
13051306
# Fall back to a component-wise approximate comparison
13061307
return all(ab -> isapprox(ab[1], ab[2]; rtol=rtol, atol=atol, nans=nans), zip(x, y))

base/linalg/uniformscaling.jl

+5-5
Original file line numberDiff line numberDiff line change
@@ -197,15 +197,16 @@ broadcast(::typeof(/), J::UniformScaling,x::Number) = UniformScaling(J.λ/x)
197197
==(J1::UniformScaling,J2::UniformScaling) = (J1.λ == J2.λ)
198198

199199
function isapprox(J1::UniformScaling{T}, J2::UniformScaling{S};
200-
rtol::Real=Base.rtoldefault(T,S), atol::Real=0, nans::Bool=false) where {T<:Number,S<:Number}
200+
atol::Real=0, rtol::Real=Base.rtoldefault(T,S,atol), nans::Bool=false) where {T<:Number,S<:Number}
201201
isapprox(J1.λ, J2.λ, rtol=rtol, atol=atol, nans=nans)
202202
end
203203
function isapprox(J::UniformScaling,A::AbstractMatrix;
204-
rtol::Real=rtoldefault(promote_leaf_eltypes(A),eltype(J)),
205-
atol::Real=0, nans::Bool=false, norm::Function=vecnorm)
204+
atol::Real=0,
205+
rtol::Real=rtoldefault(promote_leaf_eltypes(A),eltype(J),atol),
206+
nans::Bool=false, norm::Function=vecnorm)
206207
n = checksquare(A)
207208
Jnorm = norm === vecnorm ? abs(J.λ)*sqrt(n) : (norm === Base.norm ? abs(J.λ) : norm(diagm(fill(J.λ, n))))
208-
return norm(A - J) <= atol + rtol*max(norm(A), Jnorm)
209+
return norm(A - J) <= max(atol, rtol*max(norm(A), Jnorm))
209210
end
210211
isapprox(A::AbstractMatrix,J::UniformScaling;kwargs...) = isapprox(J,A;kwargs...)
211212

@@ -335,4 +336,3 @@ UniformScaling{Float64}
335336
```
336337
"""
337338
chol(J::UniformScaling, args...) = ((C, info) = _chol!(J, nothing); @assertposdef C info)
338-

test/math.jl

+5
Original file line numberDiff line numberDiff line change
@@ -669,6 +669,11 @@ end
669669
@test exp10(Float16(1.0)) === Float16(exp10(1.0))
670670
end
671671

672+
# #22742: updated isapprox semantics
673+
@test !isapprox(1.0, 1.0+1e-12, atol=1e-14)
674+
@test isapprox(1.0, 1.0+0.5*sqrt(eps(1.0)))
675+
@test !isapprox(1.0, 1.0+1.5*sqrt(eps(1.0)), atol=sqrt(eps(1.0)))
676+
672677
# test AbstractFloat fallback pr22716
673678
struct Float22716{T<:AbstractFloat} <: AbstractFloat
674679
x::T

0 commit comments

Comments
 (0)