Commit 5f1c9ae8 authored by Miroslav Kratochvil's avatar Miroslav Kratochvil
Browse files

support parsing out function definitions

parent e924dcbe
...@@ -6,6 +6,16 @@ Helper for quickly recognizing kinds of ASTs ...@@ -6,6 +6,16 @@ Helper for quickly recognizing kinds of ASTs
""" """
ast_is(ast::VPtr, what::Symbol)::Bool = ccall(sbml(what), Cint, (VPtr,), ast) != 0 ast_is(ast::VPtr, what::Symbol)::Bool = ccall(sbml(what), Cint, (VPtr,), ast) != 0
"""
parse_math_children(ast::VPtr)::Vector{Math}
Recursively parse all children of an AST node.
"""
parse_math_children(ast::VPtr)::Vector{Math} = [
parse_math(ccall(sbml(:ASTNode_getChild), VPtr, (VPtr, Cuint), ast, i - 1)) for
i = 1:ccall(sbml(:ASTNode_getNumChildren), Cuint, (VPtr,), ast)
]
""" """
parse_math(ast::VPtr)::Math parse_math(ast::VPtr)::Math
...@@ -15,26 +25,28 @@ pointer to `ASTNode_t`. ...@@ -15,26 +25,28 @@ pointer to `ASTNode_t`.
function parse_math(ast::VPtr)::Math function parse_math(ast::VPtr)::Math
if ast_is(ast, :ASTNode_isName) || ast_is(ast, :ASTNode_isConstant) if ast_is(ast, :ASTNode_isName) || ast_is(ast, :ASTNode_isConstant)
return MathIdent(get_string(ast, :ASTNode_getName)) return MathIdent(get_string(ast, :ASTNode_getName))
elseif ast_is(ast, :ASTNode_isConstantNumber) elseif ast_is(ast, :ASTNode_isInteger)
return MathVal(ccall_sbml(:ASTNode_getValue), Cdouble, (VPtr,), ast) return MathVal(ccall(sbml(:ASTNode_getInteger), Cint, (VPtr,), ast))
elseif ast_is(ast, :ASTNode_isReal)
return MathVal(ccall(sbml(:ASTNode_getReal), Cdouble, (VPtr,), ast))
elseif ast_is(ast, :ASTNode_isFunction) elseif ast_is(ast, :ASTNode_isFunction)
return MathApply( return MathApply(get_string(ast, :ASTNode_getName), parse_math_children(ast))
get_string(ast, :ASTNode_getName),
[
parse_math(ccall(sbml(:ASTNode_getChild), VPtr, (VPtr, Cuint), ast, i - 1))
for i = 1:ccall(sbml(:ASTNode_getNumChildren), Cuint, (VPtr,), ast)
],
)
elseif ast_is(ast, :ASTNode_isOperator) elseif ast_is(ast, :ASTNode_isOperator)
return MathApply( return MathApply(
string(Char(ccall(sbml(:ASTNode_getCharacter), Cchar, (VPtr,), ast))), string(Char(ccall(sbml(:ASTNode_getCharacter), Cchar, (VPtr,), ast))),
[ parse_math_children(ast),
parse_math(ccall(sbml(:ASTNode_getChild), VPtr, (VPtr, Cuint), ast, i - 1))
for i = 1:ccall(sbml(:ASTNode_getNumChildren), Cuint, (VPtr,), ast)
],
) )
elseif ast_is(ast, :ASTNode_isLambda)
children = parse_math_children(ast)
if !isempty(children)
body = pop!(children)
return MathLambda(broadcast((x::MathIdent) -> x.id, children), body)
else
@warn "invalid function definition found"
return MathIdent("?invalid?")
end
else else
@warn "unsupported math element at $ast" @warn "unsupported math element found"
return MathIdent("?unsupported?") return MathIdent("?unsupported?")
end end
end end
...@@ -331,6 +331,23 @@ function extractModel(mdl::VPtr)::Model ...@@ -331,6 +331,23 @@ function extractModel(mdl::VPtr)::Model
end end
end end
function_definitions = Dict{String,FunctionDefinition}()
for i = 1:ccall(sbml(:Model_getNumFunctionDefinitions), Cuint, (VPtr,), mdl)
fd = ccall(sbml(:Model_getFunctionDefinition), VPtr, (VPtr, Cuint), mdl, i - 1)
def = nothing
if ccall(sbml(:FunctionDefinition_isSetMath), Cint, (VPtr,), fd) != 0
def = parse_math(ccall(sbml(:FunctionDefinition_getMath), VPtr, (VPtr,), fd))
end
function_definitions[get_string(fd, :FunctionDefinition_getId)] =
FunctionDefinition(
get_optional_string(fd, :FunctionDefinition_getName),
def,
getNotes(fd),
getAnnotation(fd),
)
end
return Model( return Model(
parameters, parameters,
units, units,
......
...@@ -82,6 +82,14 @@ struct MathApply <: Math ...@@ -82,6 +82,14 @@ struct MathApply <: Math
args::Vector{Math} args::Vector{Math}
end end
"""
Function definition (aka "lambda") in mathematical expression
"""
struct MathLambda <: Math
args::Vector{String}
body::Math
end
""" """
Reaction with stoichiometry that assigns reactants and products their relative Reaction with stoichiometry that assigns reactants and products their relative
consumption/production rates (accessible in field `stoichiometry`), lower/upper consumption/production rates (accessible in field `stoichiometry`), lower/upper
...@@ -128,6 +136,17 @@ struct GeneProduct ...@@ -128,6 +136,17 @@ struct GeneProduct
GeneProduct(na, l, no = nothing, a = nothing) = new(na, l, no, a) GeneProduct(na, l, no = nothing, a = nothing) = new(na, l, no, a)
end end
"""
Custom function definition.
"""
struct FunctionDefinition
name::Maybe{String}
body::Maybe{Math}
notes::Maybe{String}
annotation::Maybe{String}
FunctionDefinition(na, b, no = nothing, a = nothing) = new(na, b, no, a)
end
""" """
Structure that collects the model-related data. Contains `parameters`, `units`, Structure that collects the model-related data. Contains `parameters`, `units`,
`compartments`, `species` and `reactions` and `gene_products`, and additional `compartments`, `species` and `reactions` and `gene_products`, and additional
...@@ -142,7 +161,8 @@ struct Model ...@@ -142,7 +161,8 @@ struct Model
species::Dict{String,Species} species::Dict{String,Species}
reactions::Dict{String,Reaction} reactions::Dict{String,Reaction}
gene_products::Dict{String,GeneProduct} gene_products::Dict{String,GeneProduct}
function_definitions::Dict{String,FunctionDefinition}
notes::Maybe{String} notes::Maybe{String}
annotation::Maybe{String} annotation::Maybe{String}
Model(p, u, c, s, r, g, n = nothing, a = nothing) = new(p, u, c, s, r, g, n, a) Model(p, u, c, s, r, g, f, n = nothing, a = nothing) = new(p, u, c, s, r, g, f, n, a)
end end
Markdown is supported
0% or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment