Skip to content
Open
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
84 changes: 68 additions & 16 deletions base/reducedim.jl
Original file line number Diff line number Diff line change
Expand Up @@ -240,7 +240,6 @@ Extract first entry of slices of array A into existing array R.
copyfirst!(R::AbstractArray, A::AbstractArray) = mapfirst!(identity, R, A)

function mapfirst!(f::F, R::AbstractArray, A::AbstractArray{<:Any,N}) where {N, F}
lsiz = check_reducedims(R, A)
t = _firstreducedslice(axes(R), axes(A))
map!(f, R, view(A, t...))
end
Expand Down Expand Up @@ -1014,21 +1013,25 @@ for (fname, op) in [(:sum, :add_sum), (:prod, :mul_prod),
end

##### findmin & findmax #####
# `iterate` based fallback for compatibility.
# The initial values of Rval are not used if the corresponding indices in Rind are 0.
#
function findminmax!(f, op, Rval, Rind, A::AbstractArray{T,N}) where {T,N}
function _findminmax!(f, op, Rval, Rind, A::AbstractArray{T,N}, keys, init = true) where {T,N}
(isempty(Rval) || isempty(A)) && return Rval, Rind
lsiz = check_reducedims(Rval, A)
for i = 1:N
axes(Rval, i) == axes(Rind, i) || throw(DimensionMismatch("Find-reduction: outputs must have the same indices"))
end
zi = zero(eltype(keys))
zi in keys && throw(ArgumentError(LazyString("`keys` containing ", zi, " is not supported!")))
if init
fill!(Rind, zi)
fill!(Rval, f(first(A)))
end
# If we're reducing along dimension 1, for efficiency we can make use of a temporary.
# Otherwise, keep the result in Rval/Rind so that we traverse A in storage order.
indsAt, indsRt = safe_tail(axes(A)), safe_tail(axes(Rval))
keep, Idefault = Broadcast.shapeindexer(indsRt)
ks = keys(A)
y = iterate(ks)
zi = zero(eltype(ks))
y = iterate(keys)
if reducedim1(Rval, A)
i1 = first(axes1(Rval))
for IA in CartesianIndices(indsAt)
Expand All @@ -1042,7 +1045,7 @@ function findminmax!(f, op, Rval, Rind, A::AbstractArray{T,N}) where {T,N}
tmpRv = tmpAv
tmpRi = k
end
y = iterate(ks, kss)
y = iterate(keys, kss)
end
@inbounds Rval[i1,IR] = tmpRv
@inbounds Rind[i1,IR] = tmpRi
Expand All @@ -1059,7 +1062,56 @@ function findminmax!(f, op, Rval, Rind, A::AbstractArray{T,N}) where {T,N}
@inbounds Rval[i,IR] = tmpAv
@inbounds Rind[i,IR] = k
end
y = iterate(ks, kss)
y = iterate(keys, kss)
end
end
end
Rval, Rind
end

# Optimized fallback for `keys` with the same axes, e.g. `LinearIndices`, `CartesianIndices`.
# We initialize `Rval`/`Rind` via `map/copyfirst!` to support non-1 based `A`.
function _findminmax!(f, op, Rval, Rind, A::AbstractArray{T,N}, keys::AbstractArray{<:Any,N}, init = true) where {T,N}
axes(keys) == axes(A) || return @invoke _findminmax!(f, op, Rval, Rind, A::AbstractArray, keys::Any, init)
(isempty(Rval) || isempty(A)) && return Rval, Rind
lsiz = check_reducedims(Rval, A)
for i = 1:N
axes(Rval, i) == axes(Rind, i) || throw(DimensionMismatch("Find-reduction: outputs must have the same indices"))
end
if init
copyfirst!(Rind, keys)
mapfirst!(f, Rval, A)
end
# If we're reducing along dimension 1, for efficiency we can make use of a temporary.
# Otherwise, keep the result in Rval/Rind so that we traverse A in storage order.
indsAt, indsRt = safe_tail(axes(A)), safe_tail(axes(Rval))
keep, Idefault = Broadcast.shapeindexer(indsRt)
if reducedim1(Rval, A)
i1 = first(axes1(Rval))
@inbounds for IA in CartesianIndices(indsAt)
IR = Broadcast.newindex(IA, keep, Idefault)
tmpRv = Rval[i1,IR]
tmpRi = Rind[i1,IR]
for i in axes(A,1)
tmpAv = f(A[i,IA])
if op(tmpRv, tmpAv)
tmpRv = tmpAv
tmpRi = keys[i,IA]
end
end
Rval[i1,IR] = tmpRv
Rind[i1,IR] = tmpRi
end
else
@inbounds for IA in CartesianIndices(indsAt)
IR = Broadcast.newindex(IA, keep, Idefault)
for i in axes(A, 1)
tmpAv = f(A[i,IA])
tmpRv = Rval[i,IR]
if op(tmpRv, tmpAv)
Rval[i,IR] = tmpAv
Rind[i,IR] = keys[i,IA]
end
end
end
end
Expand All @@ -1077,7 +1129,7 @@ $(_DOCS_ALIASING_WARNING)
"""
function findmin!(rval::AbstractArray, rind::AbstractArray, A::AbstractArray;
init::Bool=true)
findminmax!(identity, isgreater, init && !isempty(A) ? fill!(rval, first(A)) : rval, fill!(rind,zero(eltype(keys(A)))), A)
_findminmax!(identity, isgreater, rval, rind, A, keys(A), init)
end

"""
Expand Down Expand Up @@ -1133,9 +1185,9 @@ function _findmin(f, A, region::D) where {D}
end
similar(A, promote_op(f, eltype(A)), ri), zeros(eltype(keys(A)), ri)
else
fA = f(first(A))
findminmax!(f, isgreater, fill!(similar(A, _findminmax_inittype(f, A), ri), fA),
zeros(eltype(keys(A)), ri), A)
rval = similar(A, _findminmax_inittype(f, A), ri)
rind = similar(A, eltype(keys(A)), ri)
_findminmax!(f, isgreater, rval, rind, A, keys(A))
end
end

Expand All @@ -1150,7 +1202,7 @@ $(_DOCS_ALIASING_WARNING)
"""
function findmax!(rval::AbstractArray, rind::AbstractArray, A::AbstractArray;
init::Bool=true)
findminmax!(identity, isless, init && !isempty(A) ? fill!(rval, first(A)) : rval, fill!(rind,zero(eltype(keys(A)))), A)
_findminmax!(identity, isless, rval, rind, A, keys(A), init)
end

"""
Expand Down Expand Up @@ -1206,9 +1258,9 @@ function _findmax(f, A, region::D) where {D}
end
similar(A, promote_op(f, eltype(A)), ri), zeros(eltype(keys(A)), ri)
else
fA = f(first(A))
findminmax!(f, isless, fill!(similar(A, _findminmax_inittype(f, A), ri), fA),
zeros(eltype(keys(A)), ri), A)
rval = similar(A, _findminmax_inittype(f, A), ri)
rind = similar(A, eltype(keys(A)), ri)
_findminmax!(f, isless, rval, rind, A, keys(A))
end
end

Expand Down
8 changes: 8 additions & 0 deletions test/offsetarray.jl
Original file line number Diff line number Diff line change
Expand Up @@ -896,6 +896,14 @@ end
@test CartesianIndices(A) == CartesianIndices(B)
end

# issue #38660
@testset "`findmin/max` for OffsetArray" begin
ov = OffsetVector([-1, 1], 0:1)
@test @inferred(findmin(ov; dims = 1)) .|> first == (-1, 0)
ov = OffsetVector([-1, 1], -1:0)
@test @inferred(findmax(ov; dims = 1)) .|> first == (1, 0)
end

@testset "overflowing show" begin
A = OffsetArray(repeat([1], 1), typemax(Int)-1)
b = IOBuffer(maxsize=10)
Expand Down
65 changes: 40 additions & 25 deletions test/reducedim.jl
Original file line number Diff line number Diff line change
Expand Up @@ -294,33 +294,48 @@ end
end
end

struct WithIteratorsKeys{T,N} <: AbstractArray{T,N}
data::Array{T,N}
end
Base.size(a::WithIteratorsKeys) = size(a.data)
Base.getindex(a::WithIteratorsKeys, inds...) = a.data[inds...]
struct IteratorsKeys{T,A<:AbstractArray{T}}
iter::A
IteratorsKeys(iter) = new{eltype(iter),typeof(iter)}(iter)
end
Base.iterate(a::IteratorsKeys, state...) = Base.iterate(a.iter, state...)
Base.keys(a::WithIteratorsKeys) = IteratorsKeys(keys(a.data))
Base.eltype(::IteratorsKeys{T}) where {T} = T

# findmin/findmax function arguments: output type inference
@testset "findmin/findmax output type inference" begin
A = ["1" "22"; "333" "4444"]
for (tup, rval, rind) in [((1,), [1 2], [CartesianIndex(1, 1) CartesianIndex(1, 2)]),
((2,), reshape([1, 3], 2, 1), reshape([CartesianIndex(1, 1), CartesianIndex(2, 1)], 2, 1)),
((1,2), fill(1,1,1), fill(CartesianIndex(1,1),1,1))]
rval′, rind′ = findmin(length, A, dims=tup)
@test (rval, rind) == (rval′, rind′)
@test typeof(rval′) == Matrix{Int}
end
for (tup, rval, rind) in [((1,), [3 4], [CartesianIndex(2, 1) CartesianIndex(2, 2)]),
((2,), reshape([2, 4], 2, 1), reshape([CartesianIndex(1, 2), CartesianIndex(2, 2)], 2, 1)),
((1,2), fill(4,1,1), fill(CartesianIndex(2,2),1,1))]
rval′, rind′ = findmax(length, A, dims=tup)
@test (rval, rind) == (rval′, rind′)
@test typeof(rval) == Matrix{Int}
end
B = [1.5 1.0; 5.5 6.0]
for (tup, rval, rind) in [((1,), [3//2 1//1], [CartesianIndex(1, 1) CartesianIndex(1, 2)]),
((2,), reshape([1//1, 11//2], 2, 1), reshape([CartesianIndex(1, 2), CartesianIndex(2, 1)], 2, 1)),
((1,2), fill(1//1,1,1), fill(CartesianIndex(1,2),1,1))]
rval′, rind′ = findmin(Rational, B, dims=tup)
@test (rval, rind) == (rval′, rind′)
@test typeof(rval) == Matrix{Rational{Int}}
rval′, rind′ = findmin(Rational ∘ abs ∘ complex, B, dims=tup)
@test (rval, rind) == (rval′, rind′)
@test typeof(rval) == Matrix{Rational{Int}}
for wrapper in (identity, WithIteratorsKeys)
A = wrapper(["1" "22"; "333" "4444"])
for (tup, rval, rind) in [((1,), [1 2], [CartesianIndex(1, 1) CartesianIndex(1, 2)]),
((2,), reshape([1, 3], 2, 1), reshape([CartesianIndex(1, 1), CartesianIndex(2, 1)], 2, 1)),
((1,2), fill(1,1,1), fill(CartesianIndex(1,1),1,1))]
rval′, rind′ = @inferred findmin(length, A, dims=tup)
@test (rval, rind) == (rval′, rind′)
@test typeof(rval′) == Matrix{Int}
end
for (tup, rval, rind) in [((1,), [3 4], [CartesianIndex(2, 1) CartesianIndex(2, 2)]),
((2,), reshape([2, 4], 2, 1), reshape([CartesianIndex(1, 2), CartesianIndex(2, 2)], 2, 1)),
((1,2), fill(4,1,1), fill(CartesianIndex(2,2),1,1))]
rval′, rind′ = @inferred findmax(length, A, dims=tup)
@test (rval, rind) == (rval′, rind′)
@test typeof(rval) == Matrix{Int}
end
B = wrapper([1.5 1.0; 5.5 6.0])
for (tup, rval, rind) in [((1,), [3//2 1//1], [CartesianIndex(1, 1) CartesianIndex(1, 2)]),
((2,), reshape([1//1, 11//2], 2, 1), reshape([CartesianIndex(1, 2), CartesianIndex(2, 1)], 2, 1)),
((1,2), fill(1//1,1,1), fill(CartesianIndex(1,2),1,1))]
rval′, rind′ = @inferred findmin(Rational, B, dims=tup)
@test (rval, rind) == (rval′, rind′)
@test typeof(rval) == Matrix{Rational{Int}}
rval′, rind′ = @inferred findmin(Rational ∘ abs ∘ complex, B, dims=tup)
@test (rval, rind) == (rval′, rind′)
@test typeof(rval) == Matrix{Rational{Int}}
end
end
end

Expand Down