Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

support SBML writing #639

Merged
merged 3 commits into from
Jul 26, 2022
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension


Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
2 changes: 1 addition & 1 deletion .github/workflows/pr-format.yml
Original file line number Diff line number Diff line change
Expand Up @@ -57,6 +57,6 @@ jobs:
":x: Auto-formatting triggered by [this comment](${{ github.event.comment.html_url }}) failed, perhaps someone pushed to the PR in the meantime?"
fi
else
then gh pr comment ${{ github.event.issue.number }} --body \
gh pr comment ${{ github.event.issue.number }} --body \
":sunny: Auto-formatting triggered by [this comment](${{ github.event.comment.html_url }}) succeeded, but the code was already formatted correctly."
fi
2 changes: 1 addition & 1 deletion Project.toml
Original file line number Diff line number Diff line change
Expand Up @@ -28,7 +28,7 @@ MAT = "0.10"
MacroTools = "0.5.6"
OSQP = "0.6"
OrderedCollections = "1.4"
SBML = "0.10.0"
SBML = "~1.1"
StableRNGs = "1.0"
Tulip = "0.7.0, 0.8.0, 0.9.2"
julia = "1.5"
Expand Down
54 changes: 38 additions & 16 deletions src/base/types/SBMLModel.jl
Original file line number Diff line number Diff line change
Expand Up @@ -126,9 +126,10 @@ function _sbml_export_annotation(annotation)::Maybe{String}
if isnothing(annotation) || isempty(annotation)
nothing
elseif length(annotation) != 1 || first(annotation).first != ""
@_io_log @warn "Possible data loss: annotation cannot be exported to SBML" annotation
@_io_log @warn "Possible data loss: multiple annotations converted to text for SBML" annotation
join(["$k: $v" for (k, v) in annotation], "\n")
else
@_io_log @warn "Possible data loss: trying to represent annotation in SBML is unlikely to work " annotation
first(annotation).second
end
end
Expand All @@ -142,11 +143,12 @@ Return the stoichiometry of reaction with ID `rid`.
"""
function reaction_stoichiometry(m::SBMLModel, rid::String)::Dict{String,Float64}
s = Dict{String,Float64}()
for (mid, x) in m.sbml.reactions[rid].reactants
s[mid] = get(s, mid, 0.0) - x
default1(x) = isnothing(x) ? 1 : x
for sr in m.sbml.reactions[rid].reactants
s[sr.species] = get(s, sr.species, 0.0) - default1(sr.stoichiometry)
end
for (mid, x) in m.sbml.reactions[rid].products
s[mid] = get(s, mid, 0.0) + x
for sr in m.sbml.reactions[rid].products
s[sr.species] = get(s, sr.species, 0.0) + default1(sr.stoichiometry)
end
return s
end
Expand Down Expand Up @@ -186,37 +188,52 @@ function Base.convert(::Type{SBMLModel}, mm::MetabolicModel)
rxns = reactions(mm)
stoi = stoichiometry(mm)
(lbs, ubs) = bounds(mm)
comps = _default.("", metabolite_compartment.(Ref(mm), mets))
comps = _default.("compartment", metabolite_compartment.(Ref(mm), mets))
compss = Set(comps)

return SBMLModel(
SBML.Model(
compartments = Dict(comp => SBML.Compartment() for comp in compss),
compartments = Dict(
comp => SBML.Compartment(constant = true) for comp in compss
),
species = Dict(
mid => SBML.Species(
name = metabolite_name(mm, mid),
compartment = _default("", comps[mi]),
compartment = _default("compartment", comps[mi]),
formula = metabolite_formula(mm, mid),
charge = metabolite_charge(mm, mid),
constant = false,
boundary_condition = false,
only_substance_units = false,
notes = _sbml_export_notes(metabolite_notes(mm, mid)),
annotation = _sbml_export_annotation(metabolite_annotations(mm, mid)),
) for (mi, mid) in enumerate(mets)
),
reactions = Dict(
rid => SBML.Reaction(
name = reaction_name(mm, rid),
reactants = Dict(
mets[i] => -stoi[i, ri] for
reactants = [
SBML.SpeciesReference(
species = mets[i],
stoichiometry = -stoi[i, ri],
constant = true,
) for
i in SparseArrays.nonzeroinds(stoi[:, ri]) if stoi[i, ri] <= 0
),
products = Dict(
mets[i] => stoi[i, ri] for
],
products = [
SBML.SpeciesReference(
species = mets[i],
stoichiometry = stoi[i, ri],
constant = true,
) for
i in SparseArrays.nonzeroinds(stoi[:, ri]) if stoi[i, ri] > 0
),
],
kinetic_parameters = Dict(
"LOWER_BOUND" => SBML.Parameter(value = lbs[ri]),
"UPPER_BOUND" => SBML.Parameter(value = ubs[ri]),
),
lower_bound = "LOWER_BOUND",
upper_bound = "UPPER_BOUND",
gene_product_association = _maybemap(
x -> _unparse_grr(SBML.GeneProductAssociation, x),
reaction_gene_association(mm, rid),
Expand All @@ -228,13 +245,18 @@ function Base.convert(::Type{SBMLModel}, mm::MetabolicModel)
),
gene_products = Dict(
gid => SBML.GeneProduct(
label = gid,
name = gene_name(mm, gid),
notes = _sbml_export_notes(gene_notes(mm, gid)),
annotation = _sbml_export_annotation(gene_annotations(mm, gid)),
) for gid in genes(mm)
),
objective = Dict(
rid => oc for (rid, oc) in zip(rxns, objective(mm)) if oc != 0
active_objective = "objective",
objectives = Dict(
"objective" => SBML.Objective(
"maximize",
Dict(rid => oc for (rid, oc) in zip(rxns, objective(mm)) if oc != 0),
),
),
),
)
Expand Down
4 changes: 2 additions & 2 deletions src/base/utils/Reaction.jl
Original file line number Diff line number Diff line change
Expand Up @@ -12,7 +12,7 @@ See also: [`reaction_mass_balanced`](@ref)
function check_duplicate_reaction(
crxn::Reaction,
rxns::OrderedDict{String,Reaction};
only_metabolites=true
only_metabolites = true,
)
for (k, rxn) in rxns
if rxn.id != crxn.id # skip if same ID
Expand Down Expand Up @@ -120,7 +120,7 @@ julia> stoichiometry_string(req; format_id = x -> x[1:end-2])
"coa + pyr = for + accoa"
```
"""
function stoichiometry_string(req; format_id=x -> x)
function stoichiometry_string(req; format_id = x -> x)
count_prefix(n) = abs(n) == 1 ? "" : string(abs(n), " ")
substrates =
join((string(count_prefix(n), format_id(met)) for (met, n) in req if n < 0), " + ")
Expand Down
2 changes: 1 addition & 1 deletion src/base/utils/gene_associations.jl
Original file line number Diff line number Diff line change
Expand Up @@ -30,7 +30,7 @@ Convert a GeneAssociation to the corresponding `SBML.jl` structure.
function _unparse_grr(
::Type{SBML.GeneProductAssociation},
x::GeneAssociation,
)::SBML.GeneAssociation
)::SBML.GeneProductAssociation
SBML.GPAOr([SBML.GPAAnd([SBML.GPARef(j) for j in i]) for i in x])
end

Expand Down
2 changes: 2 additions & 0 deletions src/io/io.jl
Original file line number Diff line number Diff line change
Expand Up @@ -52,6 +52,8 @@ Currently, these model types are supported:
function save_model(model::MetabolicModel, file_name::String)
if endswith(file_name, ".json")
return save_json_model(model, file_name)
elseif endswith(file_name, ".xml")
return save_sbml_model(model, file_name)
elseif endswith(file_name, ".mat")
return save_mat_model(model, file_name)
else
Expand Down
16 changes: 16 additions & 0 deletions src/io/sbml.jl
Original file line number Diff line number Diff line change
Expand Up @@ -7,3 +7,19 @@ Load and return a SBML XML model in `file_name`.
function load_sbml_model(file_name::String)::SBMLModel
return SBMLModel(SBML.readSBML(file_name))
end

"""
write_sbml_model(model::MetabolicModel, file_name::String)

Write a given SBML model to `file_name`.
"""
function save_sbml_model(model::MetabolicModel, file_name::String)
m =
typeof(model) == SBMLModel ? model :
begin
@_io_log @warn "Automatically converting $(typeof(model)) to SBMLModel for saving, information may be lost."
convert(SBMLModel, model)
end

SBML.writeSBML(m.sbml, file_name)
end
21 changes: 18 additions & 3 deletions test/base/types/SBMLModel.jl
Original file line number Diff line number Diff line change
Expand Up @@ -7,10 +7,25 @@
@test Set(reactions(sbmlm)) == Set(reactions(sbmlm2))
@test Set(reactions(sbmlm)) == Set(reactions(sm))
@test Set(metabolites(sbmlm)) == Set(metabolites(sbmlm2))
sp(x) = x.species
@test all([
sbmlm.sbml.reactions[i].reactants == sbmlm2.sbml.reactions[i].reactants &&
sbmlm.sbml.reactions[i].products == sbmlm2.sbml.reactions[i].products for
i in reactions(sbmlm2)
issetequal(
sp.(sbmlm.sbml.reactions[i].reactants),
sp.(sbmlm2.sbml.reactions[i].reactants),
) && issetequal(
sp.(sbmlm.sbml.reactions[i].products),
sp.(sbmlm2.sbml.reactions[i].products),
) for i in reactions(sbmlm2)
])
st(x) = isnothing(x.stoichiometry) ? 1.0 : x.stoichiometry
@test all([
issetequal(
st.(sbmlm.sbml.reactions[i].reactants),
st.(sbmlm2.sbml.reactions[i].reactants),
) && issetequal(
st.(sbmlm.sbml.reactions[i].products),
st.(sbmlm2.sbml.reactions[i].products),
) for i in reactions(sbmlm2)
])
end

Expand Down
8 changes: 8 additions & 0 deletions test/io/json.jl
Original file line number Diff line number Diff line change
Expand Up @@ -14,3 +14,11 @@
@test sum(jlbs) == sum(slbs)
@test sum(jubs) == sum(subs)
end

@testset "Save JSON model" begin
model = load_model(CoreModel, model_paths["e_coli_core.json"])
testpath = tmpfile("modeltest.json")
save_model(model, testpath)
wrote = convert(CoreModel, load_json_model(testpath))
@test isequal(model, wrote)
end
2 changes: 1 addition & 1 deletion test/io/mat.jl
Original file line number Diff line number Diff line change
Expand Up @@ -6,8 +6,8 @@
end

@testset "Save MAT model" begin
testpath = tmpfile("iJR904-clone.mat")
loaded = load_model(CoreModel, model_paths["iJR904.mat"])
testpath = tmpfile("iJR904-clone.mat")
save_model(loaded, testpath)
wrote = load_model(CoreModel, testpath)
@test wrote isa CoreModel
Expand Down
8 changes: 8 additions & 0 deletions test/io/sbml.jl
Original file line number Diff line number Diff line change
Expand Up @@ -16,3 +16,11 @@
cm = convert(CoreModelCoupled, sbmlm)
@test n_coupling_constraints(cm) == 0
end

@testset "Save SBML model" begin
model = load_model(CoreModel, model_paths["e_coli_core.xml"])
testpath = tmpfile("modeltest.xml")
save_model(convert(SBMLModel, model), testpath)
wrote = convert(CoreModel, load_sbml_model(testpath))
@test isequal(model, wrote)
end