Skip to content

FBC3 KeyValuePair annotations + metaid #429

@pascalaldo

Description

@pascalaldo

Hi,

I recently picked up the implementation of SBML L3 FBC v3 in cobrapy (see opencobra/cobrapy#1440). All of the specification implemented so far was nicely available in libsbml, so thank you for that! The only issue I ran into is with KeyValuePair objects:
The issue is with fbc-21501 and fbc-21502 of the spec, which state that KeyValuePair can have sboTerm and metaid attributes, and an annotation subobject. This is indeed implemented in the libsbml code by having KeyValuePair inherit from SBase. However, the code to convert this into XML is not present. I think the issue lies within the functions at FbcSBasePlugin::writeKeyValuePairsAnnotation and KeyValuePair::toXML(). The latter indeed seems to lack any call to export SBase attributes etc., whilst the first lacks the ability to write the attributes metaid and sboTerm, as per fbc-21509 (we're not using fbc-21509 in cobrapy, so I do not really care about that, just noticed it).

The following code illustrates the issue:

import libsbml
from logging import getLogger
logger = getLogger(__name__)

def check(value: int, message: str) -> bool:
    """Check the libsbml return value and prints message if something happened.

    If 'value' is None, prints an error message constructed using
      'message' and then exits with status code 1. If 'value' is an integer,
      it assumes it is a libSBML return status code. If the code value is
      LIBSBML_OPERATION_SUCCESS, returns without further action; if it is not,
      prints an error message constructed using 'message' along with text from
      libSBML explaining the meaning of the code, and exits with status code 1.
    """
    valid = True
    if value is None:
        logger.error(f"Error: LibSBML returned a null value trying to <{message}>.")
        valid = False
    elif isinstance(value, int):
        if value != libsbml.LIBSBML_OPERATION_SUCCESS:
            logger.error(f"Error encountered trying to '{message}'.")
            logger.error(
                f"LibSBML returned error code {str(value)}: "
                f"{libsbml.OperationReturnValue_toString(value).strip()}"
            )
            valid = False

    return valid


TEST_SBML = """<?xml version="1.0" encoding="UTF-8"?>
<sbml xmlns="http://www.sbml.org/sbml/level3/version2/core" xmlns:fbc="http://www.sbml.org/sbml/level3/version1/fbc/version3" level="3" version="2" fbc:required="false">
<model>
    <annotation>
        <listOfKeyValuePairs xmlns="http://sbml.org/fbc/keyvaluepair">
            <keyValuePair metaid="meta_kvp_1" id="kvp_1" key="cobra_test" value="true">
                <annotation>
                    <rdf:RDF xmlns:rdf="http://www.w3.org/1999/02/22-rdf-syntax-ns#" xmlns:bqmodel="http://biomodels.net/model-qualifiers/" xmlns:bqbiol="http://biomodels.net/biology-qualifiers/">
                        <rdf:Description rdf:about="#meta_kvp_1">
                            <bqbiol:is>
                                <rdf:Bag>
                                    <rdf:li rdf:resource="http://identifiers.org/bigg.model/e_coli_core" />
                                </rdf:Bag>
                            </bqbiol:is>
                        </rdf:Description>
                    </rdf:RDF>
                </annotation>
            </keyValuePair>
        </listOfKeyValuePairs>
    </annotation>
</model>
</sbml>"""


def test_read():
    print("Reading nested metadata example.")
    doc: libsbml.SBMLDocument = libsbml.readSBMLFromString(TEST_SBML)

    model: "libsbml.Model" = doc.getModel()
    model_fbc: "libsbml.FbcModelPlugin" = model.getPlugin("fbc")

    for kvp in model_fbc.getListOfKeyValuePairs():
        cvterms = kvp.getCVTerms()
        for cv in cvterms:
            print(f"KeyValuePair->CVTerm->Resource[0]: {cv.getResourceURI(0)}")


def test_write():
    print("# Writing nested metadata example.")
    # create document with new fbc version 3 namespace
    sbmlns: libsbml.SBMLNamespaces = libsbml.SBMLNamespaces(3, 2)
    sbmlns.addPkgNamespace("fbc", 3)
    doc: libsbml.SBMLDocument = libsbml.SBMLDocument(sbmlns)
    doc_fbc: libsbml.FbcSBMLDocumentPlugin = doc.getPlugin("fbc")
    doc_fbc.setRequired(False)

    model: libsbml.Model = doc.createModel()
    model_fbc: libsbml.FbcModelPlugin = model.getPlugin("fbc")
    model.setMetaId("meta_model_e_coli_core")

    kvp: libsbml.KeyValuePair = model_fbc.createKeyValuePair()
    kvp.setKey("cobra_test")
    kvp.setValue("true")
    kvp.setId("kvp_1")
    check(kvp.setMetaId("meta_kvp_1"), "Setting metaid")

    cv: "libsbml.CVTerm" = libsbml.CVTerm()
    check(cv.setQualifierType(libsbml.BIOLOGICAL_QUALIFIER), "Set qualifier type")
    check(cv.setBiologicalQualifierType(libsbml.BQB_IS), "Set qualifier")
    check(
        cv.addResource("http://identifiers.org/bigg.model/e_coli_core"), "Add resource"
    )

    # check(model.addCVTerm(cv, newBag=True), "Adding CVTerm to model")
    check(kvp.addCVTerm(cv, newBag=True), "Adding cvterm")
    cvterms_check = kvp.getCVTerms()
    for cv_check in cvterms_check:
        print(f"KeyValuePair->CVTerm->Resource[0]: {cv_check.getResourceURI(0)}")

    sbml_str: str = libsbml.writeSBMLToString(doc)
    print("#####")
    print(sbml_str)
    print("#####")


if __name__ == "__main__":
    test_read()
    test_write()

The output of this is:

Reading nested metadata example.
KeyValuePair->CVTerm->Resource[0]: http://identifiers.org/bigg.model/e_coli_core
# Writing nested metadata example.
KeyValuePair->CVTerm->Resource[0]: http://identifiers.org/bigg.model/e_coli_core
#####
<?xml version="1.0" encoding="UTF-8"?>
<sbml xmlns="http://www.sbml.org/sbml/level3/version2/core" xmlns:fbc="http://www.sbml.org/sbml/level3/version1/fbc/version3" level="3" version="2" fbc:required="false">
  <model metaid="meta_model_e_coli_core">
    <annotation>
      <listOfKeyValuePairs xmlns="http://sbml.org/fbc/keyvaluepair">
        <keyValuePair id="kvp_1" key="cobra_test" value="true"/>
      </listOfKeyValuePairs>
    </annotation>
  </model>
</sbml>

#####

This shows that setting annotations (cvterms) to a KeyValuePair works fine and reading it from an SBML file also works. It's just writing them to xml that is missing.
Hope you can help with this, let me know if you need more info.

Metadata

Metadata

Assignees

No one assigned

    Labels

    No labels
    No labels

    Type

    No type

    Projects

    No projects

    Milestone

    No milestone

    Relationships

    None yet

    Development

    No branches or pull requests

    Issue actions