Skip to content

Commit

Permalink
Implement sortslices, deprecate sortrows/sortcols (#28332)
Browse files Browse the repository at this point in the history
As discussed on triage, `sortslices` is the higher dimensional extension
of `sortrows`/`sortcols`. The dimensions being specified are the dimensions
(and for higher dimensions the order of the dimensions) to slice along. See
the help text for an example of the higher dimensional behavior. Deprecate
sortrows/sortcols in favor of sortslices.
  • Loading branch information
Keno committed Jul 30, 2018
1 parent 0fb00c9 commit 33e8b62
Show file tree
Hide file tree
Showing 9 changed files with 185 additions and 83 deletions.
2 changes: 2 additions & 0 deletions NEWS.md
Original file line number Diff line number Diff line change
Expand Up @@ -1350,6 +1350,8 @@ Deprecated or removed

* `realmin`/`realmax` are deprecated in favor of `floatmin`/`floatmax` ([#28302]).

* `sortrows`/`sortcols` have been deprecated in favor of the more general `sortslices`.

Command-line option changes
---------------------------

Expand Down
3 changes: 3 additions & 0 deletions base/deprecated.jl
Original file line number Diff line number Diff line change
Expand Up @@ -1777,6 +1777,9 @@ end
@deprecate realmin floatmin
@deprecate realmax floatmax

@deprecate sortrows(A::AbstractMatrix; kws...) sortslices(A, dims=1, kws...)
@deprecate sortcols(A::AbstractMatrix; kws...) sortslices(A, dims=2, kws...)

# END 0.7 deprecations

# BEGIN 1.0 deprecations
Expand Down
3 changes: 1 addition & 2 deletions base/exports.jl
Original file line number Diff line number Diff line change
Expand Up @@ -419,10 +419,9 @@ export
selectdim,
sort!,
sort,
sortcols,
sortperm,
sortperm!,
sortrows,
sortslices,
dropdims,
step,
stride,
Expand Down
156 changes: 156 additions & 0 deletions base/multidimensional.jl
Original file line number Diff line number Diff line change
Expand Up @@ -1496,3 +1496,159 @@ end
function Base.showarg(io::IO, r::Iterators.Pairs{<:CartesianIndex, <:Any, <:Any, T}, toplevel) where T<:AbstractVector
print(io, "pairs(IndexCartesian(), ::$T)")
end

## sortslices

"""
sortslices(A; dims, alg::Algorithm=DEFAULT_UNSTABLE, lt=isless, by=identity, rev::Bool=false, order::Ordering=Forward)
Sort slices of an array `A`. The required keyword argument `dims` must
be either an integer or a tuple of integers. It specifies the
dimension(s) over which the slices are sorted.
E.g., if `A` is a matrix, `dims=1` will sort rows, `dims=2` will sort columns.
Note that the default comparison function on one dimensional slices sorts
lexicographically.
For the remaining keyword arguments, see the documentation of [`sort!`](@ref).
# Examples
```jldoctest
julia> sortslices([7 3 5; -1 6 4; 9 -2 8], dims=1) # Sort rows
3×3 Array{Int64,2}:
-1 6 4
7 3 5
9 -2 8
julia> sortslices([7 3 5; -1 6 4; 9 -2 8], dims=1, lt=(x,y)->isless(x[2],y[2]))
3×3 Array{Int64,2}:
9 -2 8
7 3 5
-1 6 4
julia> sortslices([7 3 5; -1 6 4; 9 -2 8], dims=1, rev=true)
3×3 Array{Int64,2}:
9 -2 8
7 3 5
-1 6 4
julia> sortslices([7 3 5; 6 -1 -4; 9 -2 8], dims=2) # Sort columns
3×3 Array{Int64,2}:
3 5 7
-1 -4 6
-2 8 9
julia> sortslices([7 3 5; 6 -1 -4; 9 -2 8], dims=2, alg=InsertionSort, lt=(x,y)->isless(x[2],y[2]))
3×3 Array{Int64,2}:
5 3 7
-4 -1 6
8 -2 9
julia> sortslices([7 3 5; 6 -1 -4; 9 -2 8], dims=2, rev=true)
3×3 Array{Int64,2}:
7 5 3
6 -4 -1
9 8 -2
```
# Higher dimensions
`sortslices` extends naturally to higher dimensions. E.g., if `A` is a
a 2x2x2 array, `sortslices(A, dims=3)` will sort slices within the 3rd dimension,
passing the 2x2 slices `A[:, :, 1]` and `A[:, :, 2]` to the comparison function.
Note that while there is no default order on higher-dimensional slices, you may
use the `by` or `lt` keyword argument to specify such an order.
If `dims` is a tuple, the order of the dimensions in `dims` is
relevant and specifies the linear order of the slices. E.g., if `A` is three
dimensional and `dims` is `(1, 2)`, the orderings of the first two dimensions
are re-arranged such such that the slices (of the remaining third dimension) are sorted.
If `dims` is `(2, 1)` instead, the same slices will be taken,
but the result order will be row-major instead.
# Higher dimensional examples
```
julia> A = permutedims(reshape([4 3; 2 1; 'A' 'B'; 'C' 'D'], (2, 2, 2)), (1, 3, 2))
2×2×2 Array{Any,3}:
[:, :, 1] =
4 3
2 1
[:, :, 2] =
'A' 'B'
'C' 'D'
julia> sortslices(A, dims=(1,2))
2×2×2 Array{Any,3}:
[:, :, 1] =
1 3
2 4
[:, :, 2] =
'D' 'B'
'C' 'A'
julia> sortslices(A, dims=(2,1))
2×2×2 Array{Any,3}:
[:, :, 1] =
1 2
3 4
[:, :, 2] =
'D' 'C'
'B' 'A'
julia> sortslices(reshape([5; 4; 3; 2; 1], (1,1,5)), dims=3, by=x->x[1,1])
1×1×5 Array{Int64,3}:
[:, :, 1] =
1
[:, :, 2] =
2
[:, :, 3] =
3
[:, :, 4] =
4
[:, :, 5] =
5
```
"""
function sortslices(A::AbstractArray; dims::Union{Integer, Tuple{Vararg{Integer}}}, kws...)
_sortslices(A, Val{dims}(); kws...)
end

# Works around inference's lack of ability to recognize partial constness
struct DimSelector{dims, T}
A::T
end
DimSelector{dims}(x::T) where {dims, T} = DimSelector{dims, T}(x)
(ds::DimSelector{dims, T})(i) where {dims, T} = i in dims ? axes(ds.A, i) : (:,)

_negdims(n, dims) = filter(i->!(i in dims), 1:n)

function compute_itspace(A, ::Val{dims}) where {dims}
negdims = _negdims(ndims(A), dims)
axs = Iterators.product(ntuple(DimSelector{dims}(A), ndims(A))...)
vec(permutedims(collect(axs), (dims..., negdims...)))
end

function _sortslices(A::AbstractArray, d::Val{dims}; kws...) where dims
itspace = compute_itspace(A, d)
vecs = map(its->view(A, its...), itspace)
p = sortperm(vecs; kws...)
if ndims(A) == 2 && isa(dims, Integer) && isa(A, Array)
# At the moment, the performance of the generic version is subpar
# (about 5x slower). Hardcode a fast-path until we're able to
# optimize this.
return dims == 1 ? A[p, :] : A[:, p]
else
B = similar(A)
for (x, its) in zip(p, itspace)
B[its...] = vecs[x]
end
B
end
end
71 changes: 0 additions & 71 deletions base/sort.jl
Original file line number Diff line number Diff line change
Expand Up @@ -37,8 +37,6 @@ export # also exported by Base
partialsort!,
partialsortperm,
partialsortperm!,
sortrows,
sortcols,
# algorithms:
InsertionSort,
QuickSort,
Expand Down Expand Up @@ -933,75 +931,6 @@ end
Av
end


"""
sortrows(A; alg::Algorithm=DEFAULT_UNSTABLE, lt=isless, by=identity, rev::Bool=false, order::Ordering=Forward)
Sort the rows of matrix `A` lexicographically.
See [`sort!`](@ref) for a description of possible
keyword arguments.
# Examples
```jldoctest
julia> sortrows([7 3 5; -1 6 4; 9 -2 8])
3×3 Array{Int64,2}:
-1 6 4
7 3 5
9 -2 8
julia> sortrows([7 3 5; -1 6 4; 9 -2 8], lt=(x,y)->isless(x[2],y[2]))
3×3 Array{Int64,2}:
9 -2 8
7 3 5
-1 6 4
julia> sortrows([7 3 5; -1 6 4; 9 -2 8], rev=true)
3×3 Array{Int64,2}:
9 -2 8
7 3 5
-1 6 4
```
"""
function sortrows(A::AbstractMatrix; kws...)
rows = [view(A, i, :) for i in axes(A,1)]
p = sortperm(rows; kws...)
A[p,:]
end

"""
sortcols(A; alg::Algorithm=DEFAULT_UNSTABLE, lt=isless, by=identity, rev::Bool=false, order::Ordering=Forward)
Sort the columns of matrix `A` lexicographically.
See [`sort!`](@ref) for a description of possible
keyword arguments.
# Examples
```jldoctest
julia> sortcols([7 3 5; 6 -1 -4; 9 -2 8])
3×3 Array{Int64,2}:
3 5 7
-1 -4 6
-2 8 9
julia> sortcols([7 3 5; 6 -1 -4; 9 -2 8], alg=InsertionSort, lt=(x,y)->isless(x[2],y[2]))
3×3 Array{Int64,2}:
5 3 7
-4 -1 6
8 -2 9
julia> sortcols([7 3 5; 6 -1 -4; 9 -2 8], rev=true)
3×3 Array{Int64,2}:
7 5 3
6 -4 -1
9 8 -2
```
"""
function sortcols(A::AbstractMatrix; kws...)
cols = [view(A, :, i) for i in axes(A,2)]
p = sortperm(cols; kws...)
A[:,p]
end

## fast clever sorting for floats ##

module Float
Expand Down
3 changes: 1 addition & 2 deletions doc/src/base/sort.md
Original file line number Diff line number Diff line change
Expand Up @@ -111,8 +111,7 @@ Base.sort!
Base.sort
Base.sortperm
Base.Sort.sortperm!
Base.Sort.sortrows
Base.Sort.sortcols
Base.Sort.sortslices
```

## Order-Related Functions
Expand Down
24 changes: 19 additions & 5 deletions test/arrayops.jl
Original file line number Diff line number Diff line change
Expand Up @@ -659,7 +659,7 @@ let A, B, C, D
# 10 repeats of each row
B = A[shuffle!(repeat(1:10, 10)), :]
C = unique(B, dims=1)
@test sortrows(C) == sortrows(A)
@test sortslices(C, dims=1) == sortslices(A, dims=1)
@test unique(B, dims=2) == B
@test unique(B', dims=2)' == C

Expand Down Expand Up @@ -1173,11 +1173,11 @@ end
@testset "sort on arrays" begin
local a = rand(3,3)

asr = sortrows(a)
asr = sortslices(a, dims=1)
@test isless(asr[1,:],asr[2,:])
@test isless(asr[2,:],asr[3,:])

asc = sortcols(a)
asc = sortslices(a, dims=2)
@test isless(asc[:,1],asc[:,2])
@test isless(asc[:,2],asc[:,3])

Expand All @@ -1187,11 +1187,11 @@ end
@test m == zeros(3, 4)
@test o == fill(1, 3, 4)

asr = sortrows(a, rev=true)
asr = sortslices(a, dims=1, rev=true)
@test isless(asr[2,:],asr[1,:])
@test isless(asr[3,:],asr[2,:])

asc = sortcols(a, rev=true)
asc = sortslices(a, dims=2, rev=true)
@test isless(asc[:,2],asc[:,1])
@test isless(asc[:,3],asc[:,2])

Expand Down Expand Up @@ -1223,6 +1223,20 @@ end
@test all(bs[:,:,1] .<= bs[:,:,2])
end

@testset "higher dimensional sortslices" begin
A = permutedims(reshape([4 3; 2 1; 'A' 'B'; 'C' 'D'], (2, 2, 2)), (1, 3, 2))
@test sortslices(A, dims=(1, 2)) ==
permutedims(reshape([1 3; 2 4; 'D' 'B'; 'C' 'A'], (2, 2, 2)), (1, 3, 2))
@test sortslices(A, dims=(2, 1)) ==
permutedims(reshape([1 2; 3 4; 'D' 'C'; 'B' 'A'], (2, 2, 2)), (1, 3, 2))
B = reshape(1:8, (2,2,2))
@test sortslices(B, dims=(3,1))[:, :, 1] == [
1 3;
5 7
]
@test sortslices(B, dims=(1,3)) == B
end

@testset "fill" begin
@test fill!(Float64[1.0], -0.0)[1] === -0.0
A = fill(1.,3,3)
Expand Down
4 changes: 2 additions & 2 deletions test/offsetarray.jl
Original file line number Diff line number Diff line change
Expand Up @@ -434,8 +434,8 @@ amin, amax = extrema(parent(A))
@test unique(A, dims=2) == OffsetArray(parent(A), first(axes(A, 1)) - 1, 0)
v = OffsetArray(rand(8), (-2,))
@test sort(v) == OffsetArray(sort(parent(v)), v.offsets)
@test sortrows(A) == OffsetArray(sortrows(parent(A)), A.offsets)
@test sortcols(A) == OffsetArray(sortcols(parent(A)), A.offsets)
@test sortslices(A, dims=1) == OffsetArray(sortslices(parent(A), dims=1), A.offsets)
@test sortslices(A, dims=2) == OffsetArray(sortslices(parent(A), dims=2), A.offsets)
@test sort(A, dims=1) == OffsetArray(sort(parent(A), dims=1), A.offsets)
@test sort(A, dims=2) == OffsetArray(sort(parent(A), dims=2), A.offsets)

Expand Down
2 changes: 1 addition & 1 deletion test/testhelpers/OffsetArrays.jl
Original file line number Diff line number Diff line change
Expand Up @@ -45,7 +45,7 @@ Base.eachindex(::IndexLinear, A::OffsetVector) = axes(A, 1)
_indices(::Tuple{}, ::Tuple{}) = ()
Base.axes1(A::OffsetArray{T,0}) where {T} = Base.Slice(1:1) # we only need to specialize this one

const OffsetAxis = Union{Integer, UnitRange, Base.Slice{<:UnitRange}}
const OffsetAxis = Union{Integer, UnitRange, Base.Slice{<:UnitRange}, Base.OneTo}
function Base.similar(A::OffsetArray, T::Type, dims::Dims)
B = similar(parent(A), T, dims)
end
Expand Down

1 comment on commit 33e8b62

@nanosoldier
Copy link
Collaborator

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Executing the daily benchmark build, I will reply here when finished:

@nanosoldier runbenchmarks(ALL, isdaily = true)

Please sign in to comment.