Loading and saving models

using COBREXA

Getting the models reliably from the repositories

For convenience, COBREXA provides a specific function download_model to download models from repositories that also automatically uses the cached downloaded version of the model if it's already downloaded, and verifies the checksum to improve reproducibility. It will print out a warning in case the model checksum does not match the expectation:

download_model(
    "http://bigg.ucsd.edu/static/models/e_coli_core.json",
    "e_coli_core.json",
    "7bedec10576cfe935b19218dc881f3fb14f890a1871448fc19a9b4ee15b448d8",
)

download_model(
    "http://bigg.ucsd.edu/static/models/e_coli_core.xml",
    "e_coli_core.xml",
    "b4db506aeed0e434c1f5f1fdd35feda0dfe5d82badcfda0e9d1342335ab31116",
)
"e_coli_core.xml"
How do I get the model hash?

We do not need to fill in the hash values immediately – instead, it's possible to simply run the function once, and then copy the reported hash value from the warning message into the script.

Loading models

To load genome-scale metabolic models, COBREXA uses the AbstractFBCModels framework to import various kinds of models including SBML, JSON and the legacy Matlab-formatted "COBRA toolbox" models.

All models can be loaded automatically using load_model; but one must import the model-type specific packages to load the functionality. (This step is required to keep the "base" COBREXA as efficient and fast-loading as possible.)

import JSONFBCModels, SBMLFBCModels

model1 = load_model("e_coli_core.json")

model2 = load_model("e_coli_core.xml")
SBMLFBCModels.SBMLFBCModel(#= 95 reactions, 72 metabolites =#)

We can now explore the contents of the models using the AbstractFBCModels' interface:

import AbstractFBCModels as A

A.reactions(model1)

A.reactions(model2)
95-element Vector{String}:
 "R_ACALD"
 "R_ACALDt"
 "R_ACKr"
 "R_ACONTa"
 "R_ACONTb"
 "R_ACt2r"
 "R_ADK1"
 "R_AKGDH"
 "R_AKGt2r"
 "R_ALCD2x"
 ⋮
 "R_SUCCt2_2"
 "R_SUCCt3"
 "R_SUCDi"
 "R_SUCOAS"
 "R_TALA"
 "R_THD2"
 "R_TKT1"
 "R_TKT2"
 "R_TPI"

Additional extractable information can be found in the documentation of the abstract models package.

Converting model types

Normally, load_model is forced to guess the model type from the filename suffix. We can specify the model type ourselves (this also allows the users to work with non-standard file suffixes, and saves some overhead and uncertainty in the guessing process):

import JSONFBCModels: JSONFBCModel

model = load_model(JSONFBCModel, "e_coli_core.json")
JSONFBCModels.JSONFBCModel(#= 95 reactions, 72 metabolites =#)

Sometimes it is useful to convert the model data to another type, such as the SBML to a JSON model structure:

model_converted_to_json = load_model("e_coli_core.xml", JSONFBCModel)
JSONFBCModels.JSONFBCModel(#= 95 reactions, 72 metabolites =#)

Or to the "Canonical Julia model" from AbstractFBCModels:

model_in_julia_structures =
    load_model(JSONFBCModel, "e_coli_core.json", A.CanonicalModel.Model)
AbstractFBCModels.CanonicalModel.Model(
  reactions = Dict{String, AbstractFBCModels.CanonicalModel.Reaction}("ACALD" =…
  metabolites = Dict{String, AbstractFBCModels.CanonicalModel.Metabolite}("glu_…
  genes = Dict{String, AbstractFBCModels.CanonicalModel.Gene}("b4301" => Abstra…
  couplings = Dict{String, AbstractFBCModels.CanonicalModel.Coupling}(),
)
Compatibility with COBREXA v1.x

CanonicalModel is a newer, cleaned-up version of the StandardModel type used in COBREXA version 1. For all code that relied on StandardModel, the canonical one should work just as well.

The above command specifies all model types explicitly, leaving least room for guessing-based errors.

If required, it is also possible to convert all model types to each other simply by using Julia's convert:

model_in_json_structure = convert(JSONFBCModel, model_in_julia_structures)
JSONFBCModels.JSONFBCModel(#= 95 reactions, 72 metabolites =#)

Saving models

The models can be saved to file storage by using save_model:

save_model(model_converted_to_json, "e_coli_core_from_sbml.json")

Expectably, the file will contain the JSON with the model description:

println(open("e_coli_core_from_sbml.json") do f
    read(f, 100)
end |> String, "...")
{"metabolites":[{"compartment":"c","name":"3-Phospho-D-glyceroyl phosphate","formula":"C3H4O10P2","i...

Note that without the conversion, it may happen that you save the model in an unexpected format!

save_model(model_in_julia_structures, "e_coli_saved_wrongly.json")
println(open("e_coli_saved_wrongly.json") do f
    read(f, 100)
end |> String, "...")
7JL5Model|)&i!@=OZAbstractFBCModelsCanonicalModelD5DictND!Reaction...

The above code has saved the CanonicalModel in the way specified by the CanonicalModel structure – which is, in this case, a binary dump of the Julia objects, instead of the expected JSON. To prevent this, you can either specify the output type yourself:

save_model(model_in_julia_structures, "e_coli_saved_right.json", JSONFBCModel)

...or use save_converted_model to guess the model type automatically from the extension:

save_converted_model(model_in_julia_structures, "e_coli_saved_automatically_right.json")
println(open("e_coli_saved_automatically_right.json") do f
    read(f, 100)
end |> String, "...")
{"metabolites":[{"compartment":"c","name":"3-Phospho-D-glyceroyl phosphate","formula":"C3H4O10P2","i...

As with load_model, there is some overhead and uncertainty associated with save_converted_model guessing the model type from extension. For that reason, it is adviseable to rely on the guessing functionality only in interactive use in REPL, and avoid it in automated scriptage altogether.


This page was generated using Literate.jl.