Skip to content

Commit

Permalink
a fix to cfunction
Browse files Browse the repository at this point in the history
Some function-valued argument slots will now be specialized as `Any`,
and we need to allow those to match Refs in cfunction.
This was triggered by the Gtk.jl tests.

undo some changes that aren't actually needed anymore

restore some commented-out tests
  • Loading branch information
JeffBezanson committed Jan 28, 2016
1 parent 0b5e87a commit b43c5e0
Show file tree
Hide file tree
Showing 15 changed files with 82 additions and 51 deletions.
16 changes: 8 additions & 8 deletions base/LineEdit.jl
Original file line number Diff line number Diff line change
Expand Up @@ -317,7 +317,7 @@ end
# Edit functionality
is_non_word_char(c) = c in " \t\n\"\\'`@\$><=:;|&{}()[].,+-*/?%^~"

function reset_key_repeats(f, s::MIState)
function reset_key_repeats(f::Function, s::MIState)
key_repeats_sav = s.key_repeats
try
s.key_repeats = 0
Expand Down Expand Up @@ -715,7 +715,7 @@ function normalize_keys(keymap::Dict)
return ret
end

function add_nested_key!(keymap::Dict, key, value::ANY; override = false)
function add_nested_key!(keymap::Dict, key, value; override = false)
i = start(key)
while !done(key, i)
c, i = next(key, i)
Expand Down Expand Up @@ -747,7 +747,7 @@ immutable KeyAlias
KeyAlias(seq) = new(normalize_key(seq))
end

match_input(k, s, term, cs, keymap) = (update_key_repeats(s, cs); return keymap_fcn(k, s, ByteString(cs)))
match_input(k::Function, s, term, cs, keymap) = (update_key_repeats(s, cs); return keymap_fcn(k, s, ByteString(cs)))
match_input(k::Void, s, term, cs, keymap) = (s,p) -> return :ok
match_input(k::KeyAlias, s, term, cs, keymap) = match_input(keymap, s, IOBuffer(k.seq), Char[], keymap)
function match_input(k::Dict, s, term=terminal(s), cs=Char[], keymap = k)
Expand All @@ -762,7 +762,7 @@ function match_input(k::Dict, s, term=terminal(s), cs=Char[], keymap = k)
end

keymap_fcn(f::Void, s, c) = (s, p) -> return :ok
function keymap_fcn(f::ANY, s, c)
function keymap_fcn(f::Function, s, c)
return (s, p) -> begin
r = f(s, p, c)
if isa(r, Symbol)
Expand Down Expand Up @@ -858,7 +858,7 @@ function add_specialisations(dict, subdict, level)
end
end

postprocess!(others::ANY) = nothing
postprocess!(others) = nothing
function postprocess!(dict::Dict)
# needs to be done first for every branch
if haskey(dict, '\0')
Expand All @@ -885,7 +885,7 @@ end
# source is the keymap specified by the user (with normalized keys)
function keymap_merge(target,source)
ret = copy(target)
direct_keys = filter((k,v::ANY) -> isa(v, Union{Function, KeyAlias, Void}), source)
direct_keys = filter((k,v) -> isa(v, Union{Function, KeyAlias, Void}), source)
# first direct entries
for key in keys(direct_keys)
add_nested_key!(ret, key, source[key]; override = true)
Expand Down Expand Up @@ -1103,7 +1103,7 @@ function reset_state(s::PrefixSearchState)
reset_state(s.histprompt.hp)
end

function transition(f, s::PrefixSearchState, mode)
function transition(f::Function, s::PrefixSearchState, mode)
if isdefined(s,:mi)
transition(f,s.mi,mode)
end
Expand Down Expand Up @@ -1533,7 +1533,7 @@ activate(m::ModalInterface, s::MIState, termbuf, term::TextTerminal) =
activate(s.current_mode, s, termbuf, term)

commit_changes(t::UnixTerminal, termbuf) = write(t, takebuf_array(termbuf.out_stream))
function transition(f, s::MIState, mode)
function transition(f::Function, s::MIState, mode)
if mode == :abort
s.aborted = true
return
Expand Down
2 changes: 1 addition & 1 deletion base/REPLCompletions.jl
Original file line number Diff line number Diff line change
Expand Up @@ -10,7 +10,7 @@ function completes_global(x, name)
return startswith(x, name) && !('#' in x)
end

function filtered_mod_names(ffunc, mod::Module, name::AbstractString, all::Bool=false, imported::Bool=false)
function filtered_mod_names(ffunc::Function, mod::Module, name::AbstractString, all::Bool=false, imported::Bool=false)
ssyms = names(mod, all, imported)
filter!(ffunc, ssyms)
syms = UTF8String[string(s) for s in ssyms]
Expand Down
2 changes: 1 addition & 1 deletion base/array.jl
Original file line number Diff line number Diff line change
Expand Up @@ -750,7 +750,7 @@ function findprev(testf::Function, A, start)
end
findlast(testf::Function, A) = findprev(testf, A, length(A))

function find(testf, A::AbstractArray)
function find(testf::Function, A::AbstractArray)
# use a dynamic-length array to store the indexes, then copy to a non-padded
# array for the return
tmpI = Array(Int, 0)
Expand Down
2 changes: 1 addition & 1 deletion base/client.jl
Original file line number Diff line number Diff line change
Expand Up @@ -126,7 +126,7 @@ end
syntax_deprecation_warnings(warn::Bool) =
ccall(:jl_parse_depwarn, Cint, (Cint,), warn) == 1

function syntax_deprecation_warnings(f, warn::Bool)
function syntax_deprecation_warnings(f::Function, warn::Bool)
prev = syntax_deprecation_warnings(warn)
try
f()
Expand Down
2 changes: 1 addition & 1 deletion base/dates/adjusters.jl
Original file line number Diff line number Diff line change
Expand Up @@ -130,7 +130,7 @@ lastdayofquarter(dt::DateTime) = DateTime(lastdayofquarter(Date(dt)))

# Temporal Adjusters
immutable DateFunction
f
f::Function
# validate boolean, single-arg inner constructor
function DateFunction(f::ANY,negate::Bool,dt::TimeType)
try
Expand Down
2 changes: 1 addition & 1 deletion base/deprecated.jl
Original file line number Diff line number Diff line change
Expand Up @@ -978,7 +978,7 @@ export fieldoffsets
# 14766
@deprecate write(io::IO, p::Ptr, nb::Integer) unsafe_write(io, p, nb)

@deprecate isgeneric(f) true
@deprecate isgeneric(f) isa(f,Function)

# need to do this manually since the front end deprecates method defs of `call`
const call = @eval function(f,args...)
Expand Down
4 changes: 2 additions & 2 deletions base/file.jl
Original file line number Diff line number Diff line change
Expand Up @@ -50,7 +50,7 @@ cd() = cd(homedir())
systemerror(:close, ccall(:close,Int32,(Int32,),fd) != 0)
end
end
@windows_only function cd(f, dir::AbstractString)
@windows_only function cd(f::Function, dir::AbstractString)
old = pwd()
try
cd(dir)
Expand All @@ -59,7 +59,7 @@ end
cd(old)
end
end
cd(f) = cd(f, homedir())
cd(f::Function) = cd(f, homedir())

function mkdir(path::AbstractString, mode::Unsigned=0o777)
@unix_only ret = ccall(:mkdir, Int32, (Cstring,UInt32), path, mode)
Expand Down
2 changes: 1 addition & 1 deletion base/tuple.jl
Original file line number Diff line number Diff line change
Expand Up @@ -41,7 +41,7 @@ eltype{T,_}(::Type{NTuple{_,T}}) = T

## mapping ##

ntuple(f, n::Integer) =
ntuple(f::Function, n::Integer) =
n<=0 ? () :
n==1 ? (f(1),) :
n==2 ? (f(1),f(2),) :
Expand Down
2 changes: 1 addition & 1 deletion base/util.jl
Original file line number Diff line number Diff line change
Expand Up @@ -311,7 +311,7 @@ end

## printing with color ##

function with_output_color(f, color::Symbol, io::IO, args...)
function with_output_color(f::Function, color::Symbol, io::IO, args...)
buf = IOBuffer()
have_color && print(buf, get(text_colors, color, color_normal))
try f(buf, args...)
Expand Down
25 changes: 16 additions & 9 deletions src/codegen.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -1170,9 +1170,16 @@ static Function *jl_cfunction_object(jl_function_t *f, jl_value_t *rt, jl_tuplet

jl_lambda_info_t *li = jl_get_specialization1((jl_tupletype_t*)sigt, NULL);
if (li != NULL) {
if (!jl_types_equal((jl_value_t*)li->specTypes, sigt)) {
jl_errorf("cfunction: type signature of %s does not match specification",
jl_symbol_name(li->name));
for(i=1; i < nargs+1; i++) {
jl_value_t *speci = jl_nth_slot_type(li->specTypes, i);
jl_value_t *sigi = jl_nth_slot_type((jl_tupletype_t*)sigt, i);
if ((isref & (2<<(i-1))) && speci == (jl_value_t*)jl_any_type) {
// specialized for Any => can accept any Ref
}
else if (!jl_types_equal(speci, sigi)) {
jl_errorf("cfunction: type signature of %s does not match specification",
jl_symbol_name(li->name));
}
}
jl_value_t *astrt = li->rettype;
if (rt != NULL) {
Expand Down Expand Up @@ -2864,8 +2871,8 @@ static Value *emit_jlcall(Value *theFptr, Value *theF, jl_value_t **args,
return emit_jlcall(theFptr, theF, argStart, nargs, ctx);
}

static jl_cgval_t emit_call_function_object(jl_lambda_info_t *li, jl_cgval_t theF, Value *theFptr,
jl_value_t **args, size_t nargs, jl_codectx_t *ctx)
static jl_cgval_t emit_call_function_object(jl_lambda_info_t *li, const jl_cgval_t &theF, Value *theFptr,
jl_value_t **args, size_t nargs, jl_value_t *callexpr, jl_codectx_t *ctx)
{
if (li->functionObjects.specFunctionObject != NULL) {
// emit specialized call site
Expand Down Expand Up @@ -2941,7 +2948,8 @@ static jl_cgval_t emit_call_function_object(jl_lambda_info_t *li, jl_cgval_t the
call->setAttributes(cf->getAttributes());
return sret ? mark_julia_slot(result, jlretty) : mark_julia_type(call, retboxed, jlretty);
}
return mark_julia_type(emit_jlcall(theFptr, boxed(theF,ctx), &args[1], nargs, ctx), true, jl_any_type); // (typ will be patched up by caller)
return mark_julia_type(emit_jlcall(theFptr, boxed(theF,ctx), &args[1], nargs, ctx), true,
expr_type(callexpr, ctx));
}

static jl_cgval_t emit_call(jl_value_t **args, size_t arglen, jl_codectx_t *ctx, jl_value_t *expr)
Expand Down Expand Up @@ -3007,13 +3015,12 @@ static jl_cgval_t emit_call(jl_value_t **args, size_t arglen, jl_codectx_t *ctx,
!jl_is_leaf_type(f)) {
jl_add_linfo_root(ctx->linfo, f);
}
fval = mark_julia_type(literal_pointer_val((jl_value_t*)f), true, jl_typeof(f));
fval = mark_julia_const((jl_value_t*)f);
}
else {
fval = emit_expr(args[0], ctx);
}
// TODO jb/functions - mark with expr_type(expr,ctx) ??
result = emit_call_function_object(li, fval, theFptr, args, nargs, ctx);
result = emit_call_function_object(li, fval, theFptr, args, nargs, expr, ctx);
ctx->gc.argDepth = argStart;
JL_GC_POP();
return result;
Expand Down
3 changes: 2 additions & 1 deletion src/gf.c
Original file line number Diff line number Diff line change
Expand Up @@ -978,7 +978,8 @@ JL_DLLEXPORT jl_lambda_info_t *jl_instantiate_staged(jl_methlist_t *m, jl_tuplet
}
ex = oldast;
}
func = (jl_lambda_info_t*)jl_toplevel_eval_in_warn(m->func->module, (jl_value_t*)ex, 1); // need to eval macros in the right module, but not give a warning for the `eval` call unless that results in a call to `eval`
// need to eval macros in the right module, but not give a warning for the `eval` call unless that results in a call to `eval`
func = (jl_lambda_info_t*)jl_toplevel_eval_in_warn(m->func->module, (jl_value_t*)ex, 1);
func->name = m->func->name;
JL_GC_POP();
return func;
Expand Down
31 changes: 23 additions & 8 deletions src/julia-syntax.scm
Original file line number Diff line number Diff line change
Expand Up @@ -424,7 +424,17 @@
;; call with keyword args pre-sorted - original method code goes here
,(method-def-expr-
mangled sparams
`((|::| ,mangled (call (top typeof) ,mangled)) ,@vars ,@restkw ,@pargl ,@vararg)
`((|::| ,mangled (call (top typeof) ,mangled)) ,@vars ,@restkw
;; strip type off function self argument if not needed for a static param.
;; then it is ok for cl-convert to move this definition above the original def.
,(if (decl? (car pargl))
(if (any (lambda (sp)
(expr-contains-eq (sparam-name sp) (caddr (car pargl))))
positional-sparams)
(car pargl)
(decl-var (car pargl)))
(car pargl))
,@(cdr pargl) ,@vararg)
`(block
,@(if (null? lno) '()
;; TODO jb/functions get a better `name` for functions specified by type
Expand Down Expand Up @@ -3004,18 +3014,19 @@ So far only the second case can actually occur.
(cl-convert (cadddr lam2) 'anon lam2 (table) #f interp)))
,(last e))))
(else
(let* ((newlam (to-goto-form
(renumber-jlgensym
(convert-lambda lam2 '|#anon| #t))))
(vi (lam:vinfo newlam))
(let* ((exprs (lift-toplevel (convert-lambda lam2 '|#anon| #t)))
(top-stmts (cdr exprs))
(newlam (to-goto-form (renumber-jlgensym (car exprs))))
(vi (lam:vinfo newlam))
;; insert `list` expression heads to make the lambda vinfo
;; lists quotable
(newlam `(lambda (list ,@(cadr newlam))
(list (list ,@(map (lambda (l) (cons 'list l))
(car vi)))
(list ,@(cadr vi)) ,(caddr vi) (list ,@(cadddr vi)))
,@(cdddr newlam))))
`(block
`(toplevel-butlast
,@top-stmts
,@sp-inits
(method ,name ,(cl-convert sig fname lam namemap toplevel interp)
,(julia-expand-macros `(quote ,newlam))
Expand Down Expand Up @@ -3055,7 +3066,9 @@ So far only the second case can actually occur.
(let* ((var-exprs (map (lambda (v)
(let ((cv (assq v (cadr (lam:vinfo lam)))))
(if cv
`(call (top getfield) ,fname (inert ,v))
(if interp
`($ (call (top QuoteNode) ,v))
`(call (top getfield) ,fname (inert ,v)))
v)))
cvs))
(P (filter identity (map (lambda (v ve)
Expand Down Expand Up @@ -3095,7 +3108,9 @@ So far only the second case can actually occur.
(block
,@(map-cl-convert (cdr (lam:body e)) 'anon
(lambda-optimize-vars! e)
(table) #t interp))))
(table)
(null? (cadr e)) ;; only toplevel thunks have 0 args
interp))))
(else (cons (car e)
(map-cl-convert (cdr e) fname lam namemap toplevel interp))))))))

Expand Down
5 changes: 4 additions & 1 deletion src/toplevel.c
Original file line number Diff line number Diff line change
Expand Up @@ -746,8 +746,11 @@ JL_DLLEXPORT void jl_method_def(jl_svec_t *argdata, jl_lambda_info_t *f, jl_valu
assert(jl_is_svec(tvars));
assert(jl_nparams(argtypes)>0);

if (jl_is_tuple_type(argtypes) && jl_nparams(argtypes) > 0 && !jl_is_type(jl_tparam0(argtypes)))
jl_error("function type in method definition is not a type");
jl_value_t *ftype = jl_first_argument_datatype((jl_value_t*)argtypes);
if (!(jl_is_type_type(ftype) ||
if (ftype == NULL ||
!(jl_is_type_type(ftype) ||
(jl_is_datatype(ftype) &&
(!((jl_datatype_t*)ftype)->abstract || jl_is_leaf_type(ftype)) &&
((jl_datatype_t*)ftype)->name->mt != NULL)))
Expand Down
29 changes: 14 additions & 15 deletions test/core.jl
Original file line number Diff line number Diff line change
Expand Up @@ -354,10 +354,10 @@ sptest1{T,S}(x::T, y::S) = 43
sptest2{T}(x::T) = T
@test is(sptest2(:a),Symbol)

#sptest3{T}(x::T) = y->T
#let m = sptest3(:a)
# @test is(m(0),Symbol)
#end
sptest3{T}(x::T) = y->T
let m = sptest3(:a)
@test is(m(0),Symbol)
end

sptest4{T}(x::T, y::T) = 42
sptest4{T}(x::T, y) = 44
Expand Down Expand Up @@ -812,8 +812,7 @@ let
local my_func, a, c
my_func{T}(P::Vector{T}, Q::Vector{T}) = 0
my_func{T}(x::T, P::Vector{T}) = 1
# todo: this gives an ambiguity warning
#my_func{T}(P::Vector{T}, x::T) = 2
my_func{T}(P::Vector{T}, x::T) = 2
a = Int[3]
c = Vector[a]

Expand All @@ -835,13 +834,13 @@ let
boor(x) = 0
boor(x::Union) = 1
@test boor(StridedArray) == 0
#@test boor(StridedArray.body) == 1
@test boor(StridedArray.body) == 1

# issue #1202
foor(x::Union) = 1
@test_throws MethodError foor(StridedArray)
@test foor(StridedArray.body) == 1
#@test_throws MethodError foor(StridedArray)
@test_throws MethodError foor(StridedArray)
end

# issue #1153
Expand Down Expand Up @@ -894,9 +893,9 @@ let
@test a2 == [101,102,909]
end

#@test unsafe_pointer_to_objref(ccall(:jl_call1, Ptr{Void}, (Any,Any),
# x -> x+1, 314158)) == 314159
#@test unsafe_pointer_to_objref(pointer_from_objref(e+pi)) == e+pi
@test unsafe_pointer_to_objref(ccall(:jl_call1, Ptr{Void}, (Any,Any),
x -> x+1, 314158)) == 314159
@test unsafe_pointer_to_objref(pointer_from_objref(e+pi)) == e+pi

let
local a, aa
Expand Down Expand Up @@ -1067,10 +1066,10 @@ let
end

# issue #2169
#let
# i2169{T}(a::Array{T}) = typemin(T)
# @test invoke(i2169, Tuple{Array} ,Int8[1]) === Int8(-128)
#end
let
i2169{T}(a::Array{T}) = typemin(T)
@test invoke(i2169, Tuple{Array} ,Int8[1]) === Int8(-128)
end

# issue #2365
type B2365{T}
Expand Down
6 changes: 6 additions & 0 deletions test/staged.jl
Original file line number Diff line number Diff line change
Expand Up @@ -151,3 +151,9 @@ module TestGeneratedThrow
cfunction(foo,Void,())
end
end

# @generated functions including inner functions
@generated function _g_f_with_inner(x)
:(y->y)
end
@test (_g_f_with_inner(1))(8) == 8

0 comments on commit b43c5e0

Please sign in to comment.