Skip to content

Commit

Permalink
This commit does not belong to any branch on this repository, and may belong to a fork outside of the repository.
Some necessary generalizations
Browse files Browse the repository at this point in the history
timholy committed Mar 21, 2016
1 parent 7325418 commit 35777df
Show file tree
Hide file tree
Showing 7 changed files with 28 additions and 24 deletions.
2 changes: 1 addition & 1 deletion base/abstractarraymath.jl
Original file line number Diff line number Diff line change
Expand Up @@ -212,7 +212,7 @@ function repmat(a::AbstractVector, m::Int)
end

# Generalized repmat
function repeat{T}(A::Array{T};
function repeat{T}(A::AbstractArray{T};
inner::Array{Int} = ones(Int, ndims(A)),
outer::Array{Int} = ones(Int, ndims(A)))
ndims_in = ndims(A)
Expand Down
6 changes: 3 additions & 3 deletions base/complex.jl
Original file line number Diff line number Diff line change
Expand Up @@ -792,7 +792,7 @@ big{T<:AbstractFloat,N}(A::AbstractArray{Complex{T},N}) = convert(AbstractArray{

promote_array_type{S<:Union{Complex, Real}, AT<:AbstractFloat}(F, ::Type{S}, ::Type{Complex{AT}}) = Complex{AT}

function complex{S<:Real,T<:Real}(A::Array{S}, B::Array{T})
function complex{S<:Real,T<:Real}(A::AbstractArray{S}, B::AbstractArray{T})
if size(A) != size(B); throw(DimensionMismatch()); end
F = similar(A, typeof(complex(zero(S),zero(T))))
for (iF, iA, iB) in zip(eachindex(F), eachindex(A), eachindex(B))
Expand All @@ -801,15 +801,15 @@ function complex{S<:Real,T<:Real}(A::Array{S}, B::Array{T})
return F
end

function complex{T<:Real}(A::Real, B::Array{T})
function complex{T<:Real}(A::Real, B::AbstractArray{T})
F = similar(B, typeof(complex(A,zero(T))))
for (iF, iB) in zip(eachindex(F), eachindex(B))
@inbounds F[iF] = complex(A, B[iB])
end
return F
end

function complex{T<:Real}(A::Array{T}, B::Real)
function complex{T<:Real}(A::AbstractArray{T}, B::Real)
F = similar(A, typeof(complex(zero(T),B)))
for (iF, iA) in zip(eachindex(F), eachindex(A))
@inbounds F[iF] = complex(A[iA], B)
Expand Down
12 changes: 6 additions & 6 deletions base/linalg/dense.jl
Original file line number Diff line number Diff line change
Expand Up @@ -150,7 +150,7 @@ function trace{T}(A::Matrix{T})
t
end

function kron{T,S}(a::Matrix{T}, b::Matrix{S})
function kron{T,S}(a::AbstractMatrix{T}, b::AbstractMatrix{S})
R = Array(promote_type(T,S), size(a,1)*size(b,1), size(a,2)*size(b,2))
m = 1
for j = 1:size(a,2), l = 1:size(b,2), i = 1:size(a,1)
Expand All @@ -163,11 +163,11 @@ function kron{T,S}(a::Matrix{T}, b::Matrix{S})
R
end

kron(a::Number, b::Union{Number, Vector, Matrix}) = a * b
kron(a::Union{Vector, Matrix}, b::Number) = a * b
kron(a::Vector, b::Vector)=vec(kron(reshape(a,length(a),1),reshape(b,length(b),1)))
kron(a::Matrix, b::Vector)=kron(a,reshape(b,length(b),1))
kron(a::Vector, b::Matrix)=kron(reshape(a,length(a),1),b)
kron(a::Number, b::Union{Number, AbstractVecOrMat}) = a * b
kron(a::AbstractVecOrMat, b::Number) = a * b
kron(a::AbstractVector, b::AbstractVector)=vec(kron(reshape(a,length(a),1),reshape(b,length(b),1)))
kron(a::AbstractMatrix, b::AbstractVector)=kron(a,reshape(b,length(b),1))
kron(a::AbstractVector, b::AbstractMatrix)=kron(reshape(a,length(a),1),b)

^(A::Matrix, p::Integer) = p < 0 ? inv(A^-p) : Base.power_by_squaring(A,p)

Expand Down
2 changes: 1 addition & 1 deletion base/pointer.jl
Original file line number Diff line number Diff line change
Expand Up @@ -26,7 +26,7 @@ cconvert(::Type{Ptr{UInt8}}, s::AbstractString) = bytestring(s)
cconvert(::Type{Ptr{Int8}}, s::AbstractString) = bytestring(s)

unsafe_convert{T}(::Type{Ptr{T}}, a::Array{T}) = ccall(:jl_array_ptr, Ptr{T}, (Any,), a)
unsafe_convert(::Type{Ptr{Void}}, a::Array) = ccall(:jl_array_ptr, Ptr{Void}, (Any,), a)
unsafe_convert{S,T}(::Type{Ptr{S}}, a::AbstractArray{T}) = convert(Ptr{S}, unsafe_convert(Ptr{T}, a))

# unsafe pointer to array conversions
pointer_to_array(p, d::Integer, own=false) = pointer_to_array(p, (d,), own)
Expand Down
15 changes: 10 additions & 5 deletions base/sparse/cholmod.jl
Original file line number Diff line number Diff line change
Expand Up @@ -957,23 +957,28 @@ end
function convert{T}(::Type{Matrix{T}}, D::Dense{T})
s = unsafe_load(D.p)
a = Array(T, s.nrow, s.ncol)
if s.d == s.nrow
unsafe_copy!(pointer(a), s.x, s.d*s.ncol)
copy!(a, D)
end
function Base.copy!(dest::AbstractArray, D::Dense)
s = unsafe_load(D.p)
if s.d == s.nrow && isa(dest, Array)
unsafe_copy!(pointer(dest), s.x, s.d*s.ncol)
else
k = 0
for j = 1:s.ncol
for i = 1:s.nrow
a[i,j] = unsafe_load(s.x, i + (j - 1)*s.d)
dest[k+=1] = unsafe_load(s.x, i + (j - 1)*s.d)
end
end
end
a
dest
end
convert{T}(::Type{Matrix}, D::Dense{T}) = convert(Matrix{T}, D)
function convert{T}(::Type{Vector{T}}, D::Dense{T})
if size(D, 2) > 1
throw(DimensionMismatch("input must be a vector but had $(size(D, 2)) columns"))
end
reshape(convert(Matrix, D), size(D, 1))
copy!(Array(T, size(D, 1)), D)
end
convert{T}(::Type{Vector}, D::Dense{T}) = convert(Vector{T}, D)

Expand Down
5 changes: 1 addition & 4 deletions base/subarray.jl
Original file line number Diff line number Diff line change
Expand Up @@ -263,10 +263,7 @@ compute_first_index(f, s, parent, dim, I::Tuple{}) = f


unsafe_convert{T,N,P<:Array,I<:Tuple{Vararg{Union{RangeIndex, NoSlice}}}}(::Type{Ptr{T}}, V::SubArray{T,N,P,I}) =
pointer(V.parent) + (first_index(V)-1)*sizeof(T)

unsafe_convert{T,N,P<:Array,I<:Tuple{Vararg{Union{RangeIndex, NoSlice}}}}(::Type{Ptr{Void}}, V::SubArray{T,N,P,I}) =
convert(Ptr{Void}, unsafe_convert(Ptr{T}, V))
unsafe_convert(Ptr{T}, V.parent) + (first_index(V)-1)*sizeof(T)

pointer(V::FastSubArray, i::Int) = pointer(V.parent, V.first_index + V.stride1*(i-1))
pointer(V::FastContiguousSubArray, i::Int) = pointer(V.parent, V.first_index + i-1)
Expand Down
10 changes: 6 additions & 4 deletions base/unicode/checkstring.jl
Original file line number Diff line number Diff line change
Expand Up @@ -56,7 +56,7 @@ Throws:
"""
function unsafe_checkstring end

function unsafe_checkstring(dat::Vector{UInt8},
function unsafe_checkstring(dat::AbstractVector{UInt8},
pos = start(dat),
endpos = endof(dat)
;
Expand Down Expand Up @@ -152,8 +152,10 @@ function unsafe_checkstring(dat::Vector{UInt8},
return totalchar, flags, num4byte, num3byte, num2byte
end

function unsafe_checkstring{T <: Union{Vector{UInt16}, Vector{UInt32}, AbstractString}}(
dat::T,
typealias AbstractString1632{Tel<:Union{UInt16,UInt32}} Union{AbstractVector{Tel}, AbstractString}

function unsafe_checkstring(
dat::AbstractString1632,
pos = start(dat),
endpos = endof(dat)
;
Expand Down Expand Up @@ -184,7 +186,7 @@ function unsafe_checkstring{T <: Union{Vector{UInt16}, Vector{UInt32}, AbstractS
ch, pos = next(dat, pos)
!is_surrogate_trail(ch) && throw(UnicodeError(UTF_ERR_NOT_TRAIL, pos, ch))
num4byte += 1
if T != Vector{UInt16}
if !(typeof(dat) <: AbstractVector{UInt16})
!accept_surrogates && throw(UnicodeError(UTF_ERR_SURROGATE, pos, ch))
flags |= UTF_SURROGATE
end
Expand Down

0 comments on commit 35777df

Please sign in to comment.