Skip to content

Commit e1a3749

Browse files
timholytkelman
authored andcommitted
Use promote_op in broadcasting and matrix multiplication
This makes it easier to support these operations on AbstractArrays with generic eltypes (cherry picked from commit 82f9a21) ref #13803 remove _pure_meta for backport also qualify Base.MulFun
1 parent b8a6d77 commit e1a3749

11 files changed

+65
-45
lines changed

base/abstractarray.jl

+9
Original file line numberDiff line numberDiff line change
@@ -1273,6 +1273,15 @@ function map_promote(f, A::AbstractArray)
12731273
return promote_to!(f, 2, dest, A)
12741274
end
12751275

1276+
# These are needed because map(eltype, As) is not inferrable
1277+
promote_eltype_op(::Any) = Bottom
1278+
promote_eltype_op{T}(op, ::AbstractArray{T}) = promote_op(op, T)
1279+
promote_eltype_op{T}(op, ::T ) = promote_op(op, T)
1280+
promote_eltype_op{R,S}(op, ::AbstractArray{R}, ::AbstractArray{S}) = promote_op(op, R, S)
1281+
promote_eltype_op{R,S}(op, ::AbstractArray{R}, ::S) = promote_op(op, R, S)
1282+
promote_eltype_op{R,S}(op, ::R, ::AbstractArray{S}) = promote_op(op, R, S)
1283+
promote_eltype_op(op, A, B, C, D...) = promote_op(op, eltype(A), promote_eltype_op(op, B, C, D...))
1284+
12761285
## 1 argument
12771286
map!{F}(f::F, A::AbstractArray) = map!(f, A, A)
12781287
function map!{F}(f::F, dest::AbstractArray, A::AbstractArray)

base/bool.jl

+2
Original file line numberDiff line numberDiff line change
@@ -59,3 +59,5 @@ fld(x::Bool, y::Bool) = div(x,y)
5959
cld(x::Bool, y::Bool) = div(x,y)
6060
rem(x::Bool, y::Bool) = y ? false : throw(DivideError())
6161
mod(x::Bool, y::Bool) = rem(x,y)
62+
63+
promote_op(op, ::Type{Bool}, ::Type{Bool}) = typeof(op(true, true))

base/broadcast.jl

+11-23
Original file line numberDiff line numberDiff line change
@@ -3,9 +3,8 @@
33
module Broadcast
44

55
using ..Cartesian
6-
import Base.promote_eltype
7-
import Base.@get!
8-
import Base.num_bit_chunks, Base._msk_end, Base.unsafe_bitgetindex
6+
using Base: promote_op, promote_eltype, promote_eltype_op, @get!, _msk_end, unsafe_bitgetindex
7+
using Base: AddFun, SubFun, MulFun, LDivFun, RDivFun, PowFun
98
import Base: .+, .-, .*, ./, .\, .//, .==, .<, .!=, .<=, .%, .<<, .>>, .^
109
export broadcast, broadcast!, broadcast_function, broadcast!_function, bitbroadcast
1110
export broadcast_getindex, broadcast_setindex!
@@ -232,7 +231,7 @@ for (Bsig, Asig, gbf, gbb) in
232231
end
233232

234233

235-
broadcast(f, As...) = broadcast!(f, Array(promote_eltype(As...), broadcast_shape(As...)), As...)
234+
broadcast(f, As...) = broadcast!(f, Array(promote_eltype_op(f, As...), broadcast_shape(As...)), As...)
236235

237236
bitbroadcast(f, As...) = broadcast!(f, BitArray(broadcast_shape(As...)), As...)
238237

@@ -286,35 +285,28 @@ end
286285

287286
## elementwise operators ##
288287

289-
.*(As::AbstractArray...) = broadcast(*, As...)
290288
.%(A::AbstractArray, B::AbstractArray) = broadcast(%, A, B)
291289
.<<(A::AbstractArray, B::AbstractArray) = broadcast(<<, A, B)
292290
.>>(A::AbstractArray, B::AbstractArray) = broadcast(>>, A, B)
293291

294-
eltype_plus(As::AbstractArray...) = promote_eltype(As...)
295-
eltype_plus(As::AbstractArray{Bool}...) = typeof(true+true)
292+
eltype_plus(As::AbstractArray...) = promote_eltype_op(AddFun(), As...)
296293

297294
.+(As::AbstractArray...) = broadcast!(+, Array(eltype_plus(As...), broadcast_shape(As...)), As...)
298295

299-
type_minus(T, S) = promote_type(T, S)
300-
type_minus(::Type{Bool}, ::Type{Bool}) = typeof(true-true)
301-
302296
function .-(A::AbstractArray, B::AbstractArray)
303-
broadcast!(-, Array(type_minus(eltype(A), eltype(B)), broadcast_shape(A,B)), A, B)
297+
broadcast!(-, Array(promote_op(SubFun(), eltype(A), eltype(B)), broadcast_shape(A,B)), A, B)
304298
end
305299

306-
type_div(T,S) = promote_type(T,S)
307-
type_div{T<:Integer,S<:Integer}(::Type{T},::Type{S}) = typeof(one(T)/one(S))
308-
type_div{T,S}(::Type{Complex{T}},::Type{Complex{S}}) = Complex{type_div(T,S)}
309-
type_div{T,S}(::Type{Complex{T}},::Type{S}) = Complex{type_div(T,S)}
310-
type_div{T,S}(::Type{T},::Type{Complex{S}}) = Complex{type_div(T,S)}
300+
eltype_mul(As::AbstractArray...) = promote_eltype_op(MulFun(), As...)
301+
302+
.*(As::AbstractArray...) = broadcast!(*, Array(eltype_mul(As...), broadcast_shape(As...)), As...)
311303

312304
function ./(A::AbstractArray, B::AbstractArray)
313-
broadcast!(/, Array(type_div(eltype(A), eltype(B)), broadcast_shape(A, B)), A, B)
305+
broadcast!(/, Array(promote_op(RDivFun(), eltype(A), eltype(B)), broadcast_shape(A, B)), A, B)
314306
end
315307

316308
function .\(A::AbstractArray, B::AbstractArray)
317-
broadcast!(\, Array(type_div(eltype(A), eltype(B)), broadcast_shape(A, B)), A, B)
309+
broadcast!(\, Array(promote_op(LDivFun(), eltype(A), eltype(B)), broadcast_shape(A, B)), A, B)
318310
end
319311

320312
typealias RatIntT{T<:Integer} Union{Type{Rational{T}},Type{T}}
@@ -327,12 +319,8 @@ function .//(A::AbstractArray, B::AbstractArray)
327319
broadcast!(//, Array(type_rdiv(eltype(A), eltype(B)), broadcast_shape(A, B)), A, B)
328320
end
329321

330-
type_pow(T,S) = promote_type(T,S)
331-
type_pow{S<:Integer}(::Type{Bool},::Type{S}) = Bool
332-
type_pow{S}(T,::Type{Rational{S}}) = type_pow(T, type_div(S, S))
333-
334322
function .^(A::AbstractArray, B::AbstractArray)
335-
broadcast!(^, Array(type_pow(eltype(A), eltype(B)), broadcast_shape(A, B)), A, B)
323+
broadcast!(^, Array(promote_op(PowFun(), eltype(A), eltype(B)), broadcast_shape(A, B)), A, B)
336324
end
337325

338326
## element-wise comparison operators returning BitArray ##

base/complex.jl

+7
Original file line numberDiff line numberDiff line change
@@ -26,6 +26,13 @@ promote_rule{T<:Real,S<:Real}(::Type{Complex{T}}, ::Type{S}) =
2626
promote_rule{T<:Real,S<:Real}(::Type{Complex{T}}, ::Type{Complex{S}}) =
2727
Complex{promote_type(T,S)}
2828

29+
promote_op{T<:Real,S<:Real}(op, ::Type{Complex{T}}, ::Type{Complex{S}}) =
30+
Complex{promote_op(op,T,S)}
31+
promote_op{T<:Real,S<:Real}(op, ::Type{Complex{T}}, ::Type{S}) =
32+
Complex{promote_op(op,T,S)}
33+
promote_op{T<:Real,S<:Real}(op, ::Type{T}, ::Type{Complex{S}}) =
34+
Complex{promote_op(op,T,S)}
35+
2936
widen{T}(::Type{Complex{T}}) = Complex{widen(T)}
3037

3138
real(z::Complex) = z.re

base/functors.jl

+4
Original file line numberDiff line numberDiff line change
@@ -120,6 +120,10 @@ EqX{T}(x::T) = EqX{T}(x)
120120

121121
call(f::EqX, y) = f.x == y
122122

123+
# More promote_op rules
124+
125+
promote_op{T<:Integer}(::PowFun, ::Type{Bool}, ::Type{T}) = Bool
126+
123127
#### Bitwise operators ####
124128

125129
# BitFunctors are functions that behave in the same bit-wise manner when applied

base/int.jl

+2
Original file line numberDiff line numberDiff line change
@@ -340,6 +340,8 @@ promote_rule(::Type{UInt128}, ::Type{Int32} ) = UInt128
340340
promote_rule(::Type{UInt128}, ::Type{Int64} ) = UInt128
341341
promote_rule(::Type{UInt128}, ::Type{Int128}) = UInt128
342342

343+
promote_op{R<:Integer,S<:Integer}(op, ::Type{R}, ::Type{S}) = typeof(op(one(R), one(S)))
344+
343345
## traits ##
344346

345347
typemin(::Type{Int8 }) = Int8(-128)

base/linalg.jl

+1
Original file line numberDiff line numberDiff line change
@@ -7,6 +7,7 @@ importall ..Base.Operators
77
import Base: USE_BLAS64, size, copy, copy_transpose!, power_by_squaring,
88
print_matrix, transpose!, unsafe_getindex, unsafe_setindex!,
99
isapprox
10+
using Base: promote_op, MulFun
1011

1112
export
1213
# Modules

base/linalg/matmul.jl

+19-19
Original file line numberDiff line numberDiff line change
@@ -38,8 +38,8 @@ function scale!(C::AbstractMatrix, b::AbstractVector, A::AbstractMatrix)
3838
end
3939
C
4040
end
41-
scale(A::Matrix, b::Vector) = scale!(similar(A, promote_type(eltype(A),eltype(b))), A, b)
42-
scale(b::Vector, A::Matrix) = scale!(similar(b, promote_type(eltype(A),eltype(b)), size(A)), b, A)
41+
scale(A::Matrix, b::Vector) = scale!(similar(A, promote_op(MulFun(),eltype(A),eltype(b))), A, b)
42+
scale(b::Vector, A::Matrix) = scale!(similar(b, promote_op(MulFun(),eltype(b),eltype(A)), size(A)), b, A)
4343

4444
# Dot products
4545

@@ -78,11 +78,11 @@ At_mul_B{T<:BlasComplex}(x::StridedVector{T}, y::StridedVector{T}) = [BLAS.dotu(
7878

7979
# Matrix-vector multiplication
8080
function (*){T<:BlasFloat,S}(A::StridedMatrix{T}, x::StridedVector{S})
81-
TS = promote_type(arithtype(T),arithtype(S))
81+
TS = promote_op(MulFun(),arithtype(T),arithtype(S))
8282
A_mul_B!(similar(x, TS, size(A,1)), A, convert(AbstractVector{TS}, x))
8383
end
8484
function (*){T,S}(A::AbstractMatrix{T}, x::AbstractVector{S})
85-
TS = promote_type(arithtype(T),arithtype(S))
85+
TS = promote_op(MulFun(),arithtype(T),arithtype(S))
8686
A_mul_B!(similar(x,TS,size(A,1)),A,x)
8787
end
8888
(*)(A::AbstractVector, B::AbstractMatrix) = reshape(A,length(A),1)*B
@@ -101,22 +101,22 @@ end
101101
A_mul_B!(y::StridedVector, A::StridedVecOrMat, x::StridedVector) = generic_matvecmul!(y, 'N', A, x)
102102

103103
function At_mul_B{T<:BlasFloat,S}(A::StridedMatrix{T}, x::StridedVector{S})
104-
TS = promote_type(arithtype(T),arithtype(S))
104+
TS = promote_op(MulFun(),arithtype(T),arithtype(S))
105105
At_mul_B!(similar(x,TS,size(A,2)), A, convert(AbstractVector{TS}, x))
106106
end
107107
function At_mul_B{T,S}(A::StridedMatrix{T}, x::StridedVector{S})
108-
TS = promote_type(arithtype(T),arithtype(S))
108+
TS = promote_op(MulFun(),arithtype(T),arithtype(S))
109109
At_mul_B!(similar(x,TS,size(A,2)), A, x)
110110
end
111111
At_mul_B!{T<:BlasFloat}(y::StridedVector{T}, A::StridedVecOrMat{T}, x::StridedVector{T}) = gemv!(y, 'T', A, x)
112112
At_mul_B!(y::StridedVector, A::StridedVecOrMat, x::StridedVector) = generic_matvecmul!(y, 'T', A, x)
113113

114114
function Ac_mul_B{T<:BlasFloat,S}(A::StridedMatrix{T}, x::StridedVector{S})
115-
TS = promote_type(arithtype(T),arithtype(S))
115+
TS = promote_op(MulFun(),arithtype(T),arithtype(S))
116116
Ac_mul_B!(similar(x,TS,size(A,2)),A,convert(AbstractVector{TS},x))
117117
end
118118
function Ac_mul_B{T,S}(A::StridedMatrix{T}, x::StridedVector{S})
119-
TS = promote_type(arithtype(T),arithtype(S))
119+
TS = promote_op(MulFun(),arithtype(T),arithtype(S))
120120
Ac_mul_B!(similar(x,TS,size(A,2)), A, x)
121121
end
122122

@@ -127,7 +127,7 @@ Ac_mul_B!(y::StridedVector, A::StridedVecOrMat, x::StridedVector) = generic_matv
127127
# Matrix-matrix multiplication
128128

129129
function (*){T,S}(A::AbstractMatrix{T}, B::StridedMatrix{S})
130-
TS = promote_type(arithtype(T), arithtype(S))
130+
TS = promote_op(MulFun(), arithtype(T), arithtype(S))
131131
A_mul_B!(similar(B, TS, (size(A,1), size(B,2))), A, B)
132132
end
133133
A_mul_B!{T<:BlasFloat}(C::StridedMatrix{T}, A::StridedVecOrMat{T}, B::StridedVecOrMat{T}) = gemm_wrapper!(C, 'N', 'N', A, B)
@@ -144,14 +144,14 @@ end
144144
A_mul_B!(C::StridedMatrix, A::StridedVecOrMat, B::StridedVecOrMat) = generic_matmatmul!(C, 'N', 'N', A, B)
145145

146146
function At_mul_B{T,S}(A::AbstractMatrix{T}, B::StridedMatrix{S})
147-
TS = promote_type(arithtype(T), arithtype(S))
147+
TS = promote_op(MulFun(),arithtype(T), arithtype(S))
148148
At_mul_B!(similar(B, TS, (size(A,2), size(B,2))), A, B)
149149
end
150150
At_mul_B!{T<:BlasFloat}(C::StridedMatrix{T}, A::StridedVecOrMat{T}, B::StridedVecOrMat{T}) = is(A,B) ? syrk_wrapper!(C, 'T', A) : gemm_wrapper!(C, 'T', 'N', A, B)
151151
At_mul_B!(C::StridedMatrix, A::StridedVecOrMat, B::StridedVecOrMat) = generic_matmatmul!(C, 'T', 'N', A, B)
152152

153153
function A_mul_Bt{T,S}(A::AbstractMatrix{T}, B::StridedMatrix{S})
154-
TS = promote_type(arithtype(T), arithtype(S))
154+
TS = promote_op(MulFun(),arithtype(T), arithtype(S))
155155
A_mul_Bt!(similar(B, TS, (size(A,1), size(B,1))), A, B)
156156
end
157157
A_mul_Bt!{T<:BlasFloat}(C::StridedMatrix{T}, A::StridedVecOrMat{T}, B::StridedVecOrMat{T}) = is(A,B) ? syrk_wrapper!(C, 'N', A) : gemm_wrapper!(C, 'N', 'T', A, B)
@@ -168,7 +168,7 @@ end
168168
A_mul_Bt!(C::StridedVecOrMat, A::StridedVecOrMat, B::StridedVecOrMat) = generic_matmatmul!(C, 'N', 'T', A, B)
169169

170170
function At_mul_Bt{T,S}(A::AbstractMatrix{T}, B::StridedVecOrMat{S})
171-
TS = promote_type(arithtype(T), arithtype(S))
171+
TS = promote_op(MulFun(),arithtype(T), arithtype(S))
172172
At_mul_Bt!(similar(B, TS, (size(A,2), size(B,1))), A, B)
173173
end
174174
At_mul_Bt!{T<:BlasFloat}(C::StridedMatrix{T}, A::StridedVecOrMat{T}, B::StridedVecOrMat{T}) = gemm_wrapper!(C, 'T', 'T', A, B)
@@ -177,7 +177,7 @@ At_mul_Bt!(C::StridedMatrix, A::StridedVecOrMat, B::StridedVecOrMat) = generic_m
177177
Ac_mul_B{T<:BlasReal}(A::StridedMatrix{T}, B::StridedMatrix{T}) = At_mul_B(A, B)
178178
Ac_mul_B!{T<:BlasReal}(C::StridedMatrix{T}, A::StridedVecOrMat{T}, B::StridedVecOrMat{T}) = At_mul_B!(C, A, B)
179179
function Ac_mul_B{T,S}(A::StridedMatrix{T}, B::StridedMatrix{S})
180-
TS = promote_type(arithtype(T), arithtype(S))
180+
TS = promote_op(MulFun(),arithtype(T), arithtype(S))
181181
Ac_mul_B!(similar(B, TS, (size(A,2), size(B,2))), A, B)
182182
end
183183
Ac_mul_B!{T<:BlasComplex}(C::StridedMatrix{T}, A::StridedVecOrMat{T}, B::StridedVecOrMat{T}) = is(A,B) ? herk_wrapper!(C,'C',A) : gemm_wrapper!(C,'C', 'N', A, B)
@@ -186,16 +186,16 @@ Ac_mul_B!(C::StridedMatrix, A::StridedVecOrMat, B::StridedVecOrMat) = generic_ma
186186
A_mul_Bc{T<:BlasFloat,S<:BlasReal}(A::StridedMatrix{T}, B::StridedMatrix{S}) = A_mul_Bt(A, B)
187187
A_mul_Bc!{T<:BlasFloat,S<:BlasReal}(C::StridedMatrix{T}, A::StridedVecOrMat{T}, B::StridedVecOrMat{S}) = A_mul_Bt!(C, A, B)
188188
function A_mul_Bc{T,S}(A::StridedMatrix{T}, B::StridedMatrix{S})
189-
TS = promote_type(arithtype(T),arithtype(S))
189+
TS = promote_op(MulFun(),arithtype(T),arithtype(S))
190190
A_mul_Bc!(similar(B,TS,(size(A,1),size(B,1))),A,B)
191191
end
192192
A_mul_Bc!{T<:BlasComplex}(C::StridedMatrix{T}, A::StridedVecOrMat{T}, B::StridedVecOrMat{T}) = is(A,B) ? herk_wrapper!(C, 'N', A) : gemm_wrapper!(C, 'N', 'C', A, B)
193193
A_mul_Bc!(C::StridedMatrix, A::StridedVecOrMat, B::StridedVecOrMat) = generic_matmatmul!(C, 'N', 'C', A, B)
194194

195-
Ac_mul_Bc{T,S}(A::AbstractMatrix{T}, B::StridedMatrix{S}) = Ac_mul_Bc!(similar(B, promote_type(arithtype(T), arithtype(S)), (size(A,2), size(B,1))), A, B)
195+
Ac_mul_Bc{T,S}(A::AbstractMatrix{T}, B::StridedMatrix{S}) = Ac_mul_Bc!(similar(B, promote_op(MulFun(),arithtype(T), arithtype(S)), (size(A,2), size(B,1))), A, B)
196196
Ac_mul_Bc!{T<:BlasFloat}(C::StridedMatrix{T}, A::StridedVecOrMat{T}, B::StridedVecOrMat{T}) = gemm_wrapper!(C, 'C', 'C', A, B)
197197
Ac_mul_Bc!(C::StridedMatrix, A::StridedVecOrMat, B::StridedVecOrMat) = generic_matmatmul!(C, 'C', 'C', A, B)
198-
Ac_mul_Bt{T,S}(A::AbstractMatrix{T}, B::StridedMatrix{S}) = Ac_mul_Bt(similar(B, promote_type(arithtype(A), arithtype(B)), (size(A,2), size(B,1))), A, B)
198+
Ac_mul_Bt{T,S}(A::AbstractMatrix{T}, B::StridedMatrix{S}) = Ac_mul_Bt(similar(B, promote_op(MulFun(),arithtype(A), arithtype(B)), (size(A,2), size(B,1))), A, B)
199199
Ac_mul_Bt!(C::StridedMatrix, A::StridedVecOrMat, B::StridedVecOrMat) = generic_matmatmul!(C, 'C', 'T', A, B)
200200

201201
# Supporting functions for matrix multiplication
@@ -420,7 +420,7 @@ end
420420
function generic_matmatmul{T,S}(tA, tB, A::AbstractVecOrMat{T}, B::AbstractMatrix{S})
421421
mA, nA = lapack_size(tA, A)
422422
mB, nB = lapack_size(tB, B)
423-
C = similar(B, promote_type(arithtype(T),arithtype(S)), mA, nB)
423+
C = similar(B, promote_op(MulFun(),arithtype(T),arithtype(S)), mA, nB)
424424
generic_matmatmul!(C, tA, tB, A, B)
425425
end
426426

@@ -625,7 +625,7 @@ end
625625

626626
# multiply 2x2 matrices
627627
function matmul2x2{T,S}(tA, tB, A::AbstractMatrix{T}, B::AbstractMatrix{S})
628-
matmul2x2!(similar(B, promote_type(T,S), 2, 2), tA, tB, A, B)
628+
matmul2x2!(similar(B, promote_op(MulFun(),T,S), 2, 2), tA, tB, A, B)
629629
end
630630

631631
function matmul2x2!{T,S,R}(C::AbstractMatrix{R}, tA, tB, A::AbstractMatrix{T}, B::AbstractMatrix{S})
@@ -654,7 +654,7 @@ end
654654

655655
# Multiply 3x3 matrices
656656
function matmul3x3{T,S}(tA, tB, A::AbstractMatrix{T}, B::AbstractMatrix{S})
657-
matmul3x3!(similar(B, promote_type(T,S), 3, 3), tA, tB, A, B)
657+
matmul3x3!(similar(B, promote_op(MulFun(),T,S), 3, 3), tA, tB, A, B)
658658
end
659659

660660
function matmul3x3!{T,S,R}(C::AbstractMatrix{R}, tA, tB, A::AbstractMatrix{T}, B::AbstractMatrix{S})

base/promotion.jl

+3
Original file line numberDiff line numberDiff line change
@@ -203,7 +203,10 @@ checked_mul(x::Integer, y::Integer) = checked_mul(promote(x,y)...)
203203
# as needed. For example, if you need to provide a custom result type
204204
# for the multiplication of two types,
205205
# promote_op{R<:MyType,S<:MyType}(::MulFun, ::Type{R}, ::Type{S}) = MyType{multype(R,S)}
206+
promote_op(::Any) = Bottom
207+
promote_op(::Any, T) = T
206208
promote_op{R,S}(::Any, ::Type{R}, ::Type{S}) = promote_type(R, S)
209+
promote_op(op, T, S, U, V...) = promote_op(op, T, promote_op(op, S, U, V...))
207210

208211
## catch-alls to prevent infinite recursion when definitions are missing ##
209212

base/sparse/sparsematrix.jl

+1-1
Original file line numberDiff line numberDiff line change
@@ -1023,7 +1023,7 @@ end # macro
10231023
(.-)(A::Number, B::SparseMatrixCSC) = A .- full(B)
10241024
( -)(A::Array , B::SparseMatrixCSC) = A - full(B)
10251025

1026-
(.*)(A::AbstractArray, B::AbstractArray) = broadcast_zpreserving(*, A, B)
1026+
(.*)(A::AbstractArray, B::AbstractArray) = broadcast_zpreserving(Base.MulFun(), A, B)
10271027
(.*)(A::SparseMatrixCSC, B::Number) = SparseMatrixCSC(A.m, A.n, copy(A.colptr), copy(A.rowval), A.nzval .* B)
10281028
(.*)(A::Number, B::SparseMatrixCSC) = SparseMatrixCSC(B.m, B.n, copy(B.colptr), copy(B.rowval), A .* B.nzval)
10291029

test/arrayops.jl

+6-2
Original file line numberDiff line numberDiff line change
@@ -1233,7 +1233,7 @@ b = rand(6,7)
12331233
# return type declarations (promote_op)
12341234
module RetTypeDecl
12351235
using Base.Test
1236-
import Base: +, *, .*
1236+
import Base: +, *, .*, zero
12371237

12381238
immutable MeterUnits{T,P} <: Number
12391239
val::T
@@ -1243,9 +1243,11 @@ module RetTypeDecl
12431243
m = MeterUnits(1.0, 1) # 1.0 meter, i.e. units of length
12441244
m2 = MeterUnits(1.0, 2) # 1.0 meter^2, i.e. units of area
12451245

1246-
(+){T}(x::MeterUnits{T,1}, y::MeterUnits{T,1}) = MeterUnits{T,1}(x.val+y.val)
1246+
(+){T,pow}(x::MeterUnits{T,pow}, y::MeterUnits{T,pow}) = MeterUnits{T,pow}(x.val+y.val)
1247+
(*){T,pow}(x::Int, y::MeterUnits{T,pow}) = MeterUnits{typeof(x*one(T)),pow}(x*y.val)
12471248
(*){T}(x::MeterUnits{T,1}, y::MeterUnits{T,1}) = MeterUnits{T,2}(x.val*y.val)
12481249
(.*){T}(x::MeterUnits{T,1}, y::MeterUnits{T,1}) = MeterUnits{T,2}(x.val*y.val)
1250+
zero{T,pow}(x::MeterUnits{T,pow}) = MeterUnits{T,pow}(zero(T))
12491251

12501252
Base.promote_op{R,S}(::Base.AddFun, ::Type{MeterUnits{R,1}}, ::Type{MeterUnits{S,1}}) = MeterUnits{promote_type(R,S),1}
12511253
Base.promote_op{R,S}(::Base.MulFun, ::Type{MeterUnits{R,1}}, ::Type{MeterUnits{S,1}}) = MeterUnits{promote_type(R,S),2}
@@ -1255,6 +1257,8 @@ module RetTypeDecl
12551257
@test @inferred([m,m]+m) == [m+m,m+m]
12561258
@test @inferred(m.*[m,m]) == [m2,m2]
12571259
@test @inferred([m,m].*m) == [m2,m2]
1260+
@test @inferred([m 2m; m m]*[m,m]) == [3m2,2m2]
1261+
@test @inferred([m m].*[m,m]) == [m2 m2; m2 m2]
12581262
end
12591263

12601264
# range, range ops

0 commit comments

Comments
 (0)