Skip to content

Commit

Permalink
fix #996, bounds check all indexes
Browse files Browse the repository at this point in the history
primitive arrayref and arrayset now support nd indexing
changing argument order of arrayset to match assign
consistently throw BoundsError for out of bounds
  • Loading branch information
JeffBezanson committed Oct 9, 2012
1 parent 6111987 commit 87e24f4
Show file tree
Hide file tree
Showing 20 changed files with 196 additions and 137 deletions.
79 changes: 29 additions & 50 deletions base/array.jl
Original file line number Diff line number Diff line change
Expand Up @@ -45,9 +45,7 @@ function copy_to_unsafe{T}(dest::Array{T}, dsto, src::Array{T}, so, N)
pointer(dest, dsto), pointer(src, so), N*sizeof(T))
else
for i=0:N-1
# NOTE: this works around the performance problem caused by all
# the doubled definitions of assign()
arrayset(dest, i+dsto, src[i+so])
arrayset(dest, src[i+so], i+dsto)
end
end
return dest
Expand Down Expand Up @@ -303,32 +301,21 @@ end

## Indexing: ref ##

ref(a::Array, i::Int) = arrayref(a,i)
ref(a::Array, i::Integer) = arrayref(a,int(i))
ref{T}(a::Array{T,0}) = arrayref(a,1)
ref{T}(a::Array{T,1}, i::Int) = arrayref(a,i)
ref{T}(a::Array{T,1}, i::Integer) = arrayref(a,int(i))
ref(a::Array{Any,1}, i::Int) = arrayref(a,i)
ref(a::Array{Any,1}, i::Integer) = arrayref(a,int(i))

# @Jeff: begin fixme issue #996
ref(A::Array, i0::Integer, i1::Integer) = A[i0 + size(A,1)*(i1-1)]
ref(A::Array, i0::Integer) = arrayref(A,int(i0))
ref(A::Array, i0::Integer, i1::Integer) = arrayref(A,int(i0),int(i1))
ref(A::Array, i0::Integer, i1::Integer, i2::Integer) =
A[i0 + size(A,1)*((i1-1) + size(A,2)*(i2-1))]
arrayref(A,int(i0),int(i1),int(i2))
ref(A::Array, i0::Integer, i1::Integer, i2::Integer, i3::Integer) =
A[i0 + size(A,1)*((i1-1) + size(A,2)*((i2-1) + size(A,3)*(i3-1)))]
arrayref(A,int(i0),int(i1),int(i2),int(i3))
ref(A::Array, i0::Integer, i1::Integer, i2::Integer, i3::Integer, i4::Integer) =
arrayref(A,int(i0),int(i1),int(i2),int(i3),int(i4))
ref(A::Array, i0::Integer, i1::Integer, i2::Integer, i3::Integer, i4::Integer, i5::Integer) =
arrayref(A,int(i0),int(i1),int(i2),int(i3),int(i4),int(i5))

function ref(A::Array, I::Integer...)
ndims = length(I)
index = I[1]
stride = 1
for k=2:ndims
stride = stride * size(A, k-1)
index += (I[k]-1) * stride
end
return A[index]
end
# end fixme issue #996
ref(A::Array, i0::Integer, i1::Integer, i2::Integer, i3::Integer, i4::Integer, i5::Integer, I::Int...) =
arrayref(A,int(i0),int(i1),int(i2),int(i3),int(i4),int(i5),I...)

# Fast copy using copy_to for Range1
function ref(A::Array, I::Range1{Int})
Expand Down Expand Up @@ -447,31 +434,23 @@ ref{T<:Integer}(A::Matrix, I::AbstractVector{T}, J::AbstractVector{Bool}) = A[I,
ref{T<:Integer}(A::Matrix, I::AbstractVector{Bool}, J::AbstractVector{T}) = A[find(I),J]

## Indexing: assign ##
# Jeff: begin fixme #996
assign(A::Array{Any}, x::ANY, i::Integer) = arrayset(A,int(i),x)
assign{T}(A::Array{T}, x, i::Integer) = arrayset(A,int(i),convert(T, x))
assign{T}(A::Array{T,0}, x) = arrayset(A,1,convert(T, x))

assign(A::Array, x, i0::Integer, i1::Integer) =
assign(A, x, i0 + size(A,1)*(i1-1))

assign(A::Array, x, i0::Integer, i1::Integer, i2::Integer) =
assign(A, x, i0 + size(A,1)*((i1-1) + size(A,2)*(i2-1)))

assign(A::Array, x, i0::Integer, i1::Integer, i2::Integer, i3::Integer) =
assign(A, x, i0 + size(A,1)*((i1-1) + size(A,2)*((i2-1) + size(A,3)*(i3-1))))

function assign(A::Array, x, I0::Integer, I::Integer...)
index = I0
stride = 1
for k=1:length(I)
stride = stride * size(A, k)
index += (I[k]-1) * stride
end
A[index] = x
return A
end
# end fixme 996
assign{T}(A::Array{T,0}, x) = arrayset(A, convert(T,x), 1)

assign(A::Array{Any}, x::ANY, i::Integer) = arrayset(A, x, int(i))

assign{T}(A::Array{T}, x, i0::Integer) = arrayset(A, convert(T,x), int(i0))
assign{T}(A::Array{T}, x, i0::Integer, i1::Integer) =
arrayset(A, convert(T,x), int(i0), int(i1))
assign{T}(A::Array{T}, x, i0::Integer, i1::Integer, i2::Integer) =
arrayset(A, convert(T,x), int(i0), int(i1), int(i2))
assign{T}(A::Array{T}, x, i0::Integer, i1::Integer, i2::Integer, i3::Integer) =
arrayset(A, convert(T,x), int(i0), int(i1), int(i2), int(i3))
assign{T}(A::Array{T}, x, i0::Integer, i1::Integer, i2::Integer, i3::Integer, i4::Integer) =
arrayset(A, convert(T,x), int(i0), int(i1), int(i2), int(i3), int(i4))
assign{T}(A::Array{T}, x, i0::Integer, i1::Integer, i2::Integer, i3::Integer, i4::Integer, i5::Integer) =
arrayset(A, convert(T,x), int(i0), int(i1), int(i2), int(i3), int(i4), int(i5))
assign{T}(A::Array{T}, x, i0::Integer, i1::Integer, i2::Integer, i3::Integer, i4::Integer, i5::Integer, I::Int...) =
arrayset(A, convert(T,x), int(i0), int(i1), int(i2), int(i3), int(i4), int(i5), I...)

function assign{T<:Integer}(A::Array, x, I::AbstractVector{T})
for i in I
Expand Down Expand Up @@ -713,7 +692,7 @@ end

function push(a::Array{Any,1}, item::ANY)
ccall(:jl_array_grow_end, Void, (Any, Uint), a, 1)
arrayset(a, length(a), item)
arrayset(a, item, length(a))
return a
end

Expand Down
2 changes: 1 addition & 1 deletion base/base.jl
Original file line number Diff line number Diff line change
Expand Up @@ -124,7 +124,7 @@ function append_any(xs...)
i = 1
for x in xs
for y in x
arrayset(out, i, y)
arrayset(out, y, i)
i += 1
end
end
Expand Down
4 changes: 2 additions & 2 deletions base/cell.jl
Original file line number Diff line number Diff line change
Expand Up @@ -2,15 +2,15 @@ function cell_1d(xs::ANY...)
n = length(xs)
a = Array(Any,n)
for i=1:n
arrayset(a,i,xs[i])
arrayset(a,xs[i],i)
end
a
end

function cell_2d(nr, nc, xs::ANY...)
a = Array(Any,nr,nc)
for i=1:(nr*nc)
arrayset(a,i,xs[i])
arrayset(a,xs[i],i)
end
a
end
Expand Down
4 changes: 1 addition & 3 deletions base/deepcopy.jl
Original file line number Diff line number Diff line change
Expand Up @@ -61,9 +61,7 @@ function _deepcopy_array_t(x, T, stackdict::ObjectIdDict)
while true
try
for i=i0:length(x)
# NOTE: this works around the performance problem caused by all
# the doubled definitions of assign()
arrayset(dest, i, deepcopy_internal(x[i], stackdict))
arrayset(dest, deepcopy_internal(x[i], stackdict), i)
end
break
catch err
Expand Down
2 changes: 1 addition & 1 deletion base/env.jl
Original file line number Diff line number Diff line change
Expand Up @@ -65,7 +65,7 @@ done(::EnvHash, i) = (ccall(:jl_environ, Any, (Int32,), i) == nothing)
function next(::EnvHash, i)
env = ccall(:jl_environ, Any, (Int32,), i)
if env == nothing
error("index out of range")
error(BoundsError)
end
env::ByteString
m = match(r"^(.*?)=(.*)$"s, env)
Expand Down
20 changes: 17 additions & 3 deletions base/inference.jl
Original file line number Diff line number Diff line change
Expand Up @@ -128,9 +128,9 @@ t_func[method_exists] = (2, 2, cmp_tfunc)
t_func[applicable] = (1, Inf, (f, args...)->Bool)
t_func[tuplelen] = (1, 1, x->Int)
t_func[arraylen] = (1, 1, x->Int)
t_func[arrayref] = (2, 2, (a,i)->(isa(a,CompositeKind) && subtype(a,Array) ?
a.parameters[1] : Any))
t_func[arrayset] = (3, 3, (a,i,v)->a)
#t_func[arrayref] = (2,Inf,(a,i...)->(isa(a,CompositeKind) && subtype(a,Array) ?
# a.parameters[1] : Any))
#t_func[arrayset] = (3, Inf, (a,v,i...)->a)
arraysize_tfunc(a, d) = Int
function arraysize_tfunc(a)
if isa(a,CompositeKind) && subtype(a,Array)
Expand Down Expand Up @@ -348,6 +348,20 @@ function builtin_tfunction(f::ANY, args::ANY, argtypes::ANY)
if is(f,tuple)
return limit_tuple_depth(argtypes)
end
if is(f,arrayset)
if length(argtypes) < 3
return None
end
return argtypes[1]
end
if is(f,arrayref)
if length(argtypes) < 2
return None
end
a = argtypes[1]
return (isa(a,CompositeKind) && subtype(a,Array) ?
a.parameters[1] : Any)
end
tf = get(t_func::ObjectIdDict, f, false)
if is(tf,false)
# struct constructor
Expand Down
2 changes: 1 addition & 1 deletion base/pcre.jl
Original file line number Diff line number Diff line change
Expand Up @@ -99,7 +99,7 @@ study(re::Array{Uint8}) = study(re, int32(0))
function exec(regex::Array{Uint8}, extra::Ptr{Void},
str::ByteString, offset::Integer, options::Integer, cap::Bool)
if offset < 0 || length(str) < offset
error("index out of range")
error(BoundsError)
end
ncap = info(regex, extra, INFO_CAPTURECOUNT, Int32)
ovec = Array(Int32, 3(ncap+1))
Expand Down
2 changes: 1 addition & 1 deletion base/regex.jl
Original file line number Diff line number Diff line change
Expand Up @@ -106,7 +106,7 @@ match(r::Regex, s::String) = match(r, s, start(s))
function search(str::ByteString, re::Regex, idx::Integer)
len = length(str)
if idx >= len+2
return idx == len+2 ? (0,0) : error("index out of range")
return idx == len+2 ? (0,0) : error(BoundsError)
end
opts = re.options & PCRE.EXECUTE_MASK
m, n = PCRE.exec(re.regex, re.extra, str, idx-1, opts, true)
Expand Down
16 changes: 8 additions & 8 deletions base/string.jl
Original file line number Diff line number Diff line change
Expand Up @@ -160,7 +160,7 @@ end
typealias Chars Union(Char,AbstractVector{Char},Set{Char})

function strchr(s::String, c::Chars, i::Integer)
if i < 1 error("index out of range") end
if i < 1 error(BoundsError) end
i = nextind(s,i-1)
while !done(s,i)
d, j = next(s,i)
Expand All @@ -179,7 +179,7 @@ function search(s::String, c::Chars, i::Integer)
if isempty(c)
return 1 <= i <= length(s)+1 ? (i,i) :
i == length(s)+2 ? (0,0) :
error("index out of range")
error(BoundsError)
end
i=strchr(s,c,i)
(i, nextind(s,i))
Expand All @@ -190,7 +190,7 @@ function search(s::String, t::String, i::Integer)
if isempty(t)
return 1 <= i <= length(s)+1 ? (i,i) :
i == length(s)+2 ? (0,0) :
error("index out of range")
error(BoundsError)
end
t1, j2 = next(t,start(t))
while true
Expand Down Expand Up @@ -342,7 +342,7 @@ SubString(s::String, i::Integer) = SubString(s, i, length(s))

function next(s::SubString, i::Int)
if i < 1 || i > s.length
error("string index out of bounds")
error(BoundsError)
end
c, i = next(s.string, i+s.offset)
c, i-s.offset
Expand All @@ -356,7 +356,7 @@ length(s::SubString) = s.length

function ref(s::String, r::Range1{Int})
if first(r) < 1 || length(s) < last(r)
error("index out of range")
error(BoundsError)
end
SubString(s, first(r), last(r))
end
Expand All @@ -373,7 +373,7 @@ strlen(s::RepString) = strlen(s.string)*s.repeat

function next(s::RepString, i::Int)
if i < 1 || i > length(s)
error("string index out of bounds")
error(BoundsError)
end
j = mod1(i,length(s.string))
c, k = next(s.string, j)
Expand Down Expand Up @@ -1223,9 +1223,9 @@ end
# find the index of the first occurrence of a value in a byte array

function memchr(a::Array{Uint8,1}, b::Integer, i::Integer)
if i < 1 error("index out of range") end
if i < 1 error(BoundsError) end
n = length(a)
if i > n return i == n+1 ? 0 : error("index out of range") end
if i > n return i == n+1 ? 0 : error(BoundsError) end
p = pointer(a)
q = ccall(:memchr, Ptr{Uint8}, (Ptr{Uint8}, Int32, Uint), p+i-1, b, n-i+1)
q == C_NULL ? 0 : int(q-p+1)
Expand Down
4 changes: 2 additions & 2 deletions extras/sparse.jl
Original file line number Diff line number Diff line change
Expand Up @@ -624,7 +624,7 @@ ref(A::SparseMatrixCSC, i::Integer) = ref(A, ind2sub(size(A),i))
ref(A::SparseMatrixCSC, I::(Integer,Integer)) = ref(A, I[1], I[2])

function ref{T}(A::SparseMatrixCSC{T}, i0::Integer, i1::Integer)
if !(1 <= i0 <= A.m && 1 <= i1 <= A.n); error("ref: index out of bounds"); end
if !(1 <= i0 <= A.m && 1 <= i1 <= A.n); error(BoundsError); end
first = A.colptr[i1]
last = A.colptr[i1+1]-1
while first <= last
Expand Down Expand Up @@ -677,7 +677,7 @@ assign(A::SparseMatrixCSC, v, i::Integer) = assign(A, v, ind2sub(size(A),i)...)
function assign{T,T_int}(A::SparseMatrixCSC{T,T_int}, v, i0::Integer, i1::Integer)
i0 = convert(T_int, i0)
i1 = convert(T_int, i1)
if !(1 <= i0 <= A.m && 1 <= i1 <= A.n); error("assign: index out of bounds"); end
if !(1 <= i0 <= A.m && 1 <= i1 <= A.n); error(BoundsError); end
v = convert(T, v)
if v == 0 #either do nothing or delete entry if it exists
first = A.colptr[i1]
Expand Down
1 change: 1 addition & 0 deletions src/alloc.c
Original file line number Diff line number Diff line change
Expand Up @@ -58,6 +58,7 @@ jl_value_t *jl_overflow_exception;
jl_value_t *jl_inexact_exception;
jl_value_t *jl_undefref_exception;
jl_value_t *jl_interrupt_exception;
jl_value_t *jl_bounds_exception;
jl_value_t *jl_memory_exception;

jl_sym_t *call_sym; jl_sym_t *dots_sym;
Expand Down
56 changes: 39 additions & 17 deletions src/array.c
Original file line number Diff line number Diff line change
Expand Up @@ -370,20 +370,46 @@ jl_value_t *jl_arrayref(jl_array_t *a, size_t i)
return elt;
}

static size_t array_nd_index(jl_array_t *a, jl_value_t **args, size_t nidxs,
char *fname)
{
size_t i=0;
if (nidxs == 1) {
if (!jl_is_long(args[0]))
jl_type_error(fname, (jl_value_t*)jl_long_type, args[0]);
i = jl_unbox_long(args[0])-1;
}
else {
size_t k, stride=1;
size_t nd = jl_array_ndims(a);
for(k=0; k < nidxs; k++) {
if (!jl_is_long(args[k]))
jl_type_error(fname, (jl_value_t*)jl_long_type, args[k]);
size_t ii = jl_unbox_long(args[k])-1;
size_t d = k>=nd ? 1 : jl_array_dim(a, k);
if (k < nidxs-1 && ii >= d) {
jl_raise(jl_bounds_exception);
}
i += ii * stride;
stride = stride * d;
}
}
if (i >= a->length) {
jl_raise(jl_bounds_exception);
}
return i;
}

JL_CALLABLE(jl_f_arrayref)
{
JL_NARGS(arrayref, 2, 2);
JL_NARGSV(arrayref, 2);
JL_TYPECHK(arrayref, array, args[0]);
JL_TYPECHK(arrayref, long, args[1]);
jl_array_t *a = (jl_array_t*)args[0];
size_t i = jl_unbox_long(args[1])-1;
if (i >= a->length) {
jl_errorf("ref array[%d]: index out of range", i+1);
}
size_t i = array_nd_index(a, &args[1], nargs-1, "arrayref");
return jl_arrayref(a, i);
}

void jl_arrayset(jl_array_t *a, size_t i, jl_value_t *rhs)
void jl_arrayset(jl_array_t *a, jl_value_t *rhs, size_t i)
{
jl_value_t *el_type = jl_tparam0(jl_typeof(a));
if (el_type != (jl_value_t*)jl_any_type) {
Expand All @@ -400,15 +426,11 @@ void jl_arrayset(jl_array_t *a, size_t i, jl_value_t *rhs)

JL_CALLABLE(jl_f_arrayset)
{
JL_NARGS(arrayset, 3, 3);
JL_NARGSV(arrayset, 3);
JL_TYPECHK(arrayset, array, args[0]);
JL_TYPECHK(arrayset, long, args[1]);
jl_array_t *b = (jl_array_t*)args[0];
size_t i = jl_unbox_long(args[1])-1;
if (i >= b->length) {
jl_errorf("assign array[%d]: index out of range", i+1);
}
jl_arrayset(b, i, args[2]);
jl_array_t *a = (jl_array_t*)args[0];
size_t i = array_nd_index(a, &args[2], nargs-2, "arrayset");
jl_arrayset(a, args[1], i);
return args[0];
}

Expand Down Expand Up @@ -454,7 +476,7 @@ void jl_array_grow_end(jl_array_t *a, size_t inc)
void jl_array_del_end(jl_array_t *a, size_t dec)
{
if (dec > a->length)
jl_error("array_del_end: index out of range");
jl_raise(jl_bounds_exception);
char *ptail = (char*)a->data + (a->length-dec)*a->elsize;
if (a->ptrarray)
memset(ptail, 0, dec*a->elsize);
Expand Down Expand Up @@ -507,7 +529,7 @@ void jl_array_del_beg(jl_array_t *a, size_t dec)
if (dec == 0)
return;
if (dec > a->length)
jl_error("array_del_beg: index out of range");
jl_raise(jl_bounds_exception);
size_t es = a->elsize;
size_t nb = dec*es;
memset(a->data, 0, nb);
Expand Down
Loading

0 comments on commit 87e24f4

Please sign in to comment.