Skip to content

Commit

Permalink
got model five working
Browse files Browse the repository at this point in the history
there are some untested things in here
  • Loading branch information
skeating committed Sep 7, 2024
1 parent a5e453d commit c19d5b0
Show file tree
Hide file tree
Showing 12 changed files with 1,530 additions and 21 deletions.
77 changes: 76 additions & 1 deletion src/sbml/conversion/SBMLRateRuleConverter.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -196,6 +196,74 @@ SBMLRateRuleConverter::matchesProperties(const ConversionProperties &props) cons
}


void
SBMLRateRuleConverter::catchAnomalies()
{
// if a rate rule is actually assigned to a variable that is itself
// the subject of an assignment rule we want to replace that math
// before we do any analysis
AssignmentRule* assignment;
Model* model = mDocument->getModel();
for (unsigned int i = 0; i < mDocument->getModel()->getNumRules(); ++i)
{
Rule* rule = model->getRule(i);
if (rule->getTypeCode() != SBML_RATE_RULE)
{
continue;
}
ASTNode* math = const_cast<ASTNode*>(rule->getMath());
if (math->getType() == AST_NAME)
{
assignment = model->getAssignmentRuleByVariable(math->getName());
if (assignment == NULL)
{
continue;
}

const ASTNode* assignment_math = assignment->getMath();
math->setType(assignment_math->getType());
for (unsigned int n = 0; n < assignment_math->getNumChildren(); ++n)
{
math->addChild(assignment_math->getChild(n)->deepCopy());
}
}

}


//if (this->getType() == AST_NAME)
//{
// // here the node is a variable
// std::string name = this->getName();
// const SBase* parent = this->getParentSBMLObject();
// const Model* model = static_cast<const Model*>(parent->getAncestorOfType(SBML_MODEL));
// const AssignmentRule* rule = model->getAssignmentRuleByVariable(name);

// // and the variable is assigned by rule
// // we want to replace this node by the rule
// if (rule != NULL)
// {
// const ASTNode* math = rule->getMath();
// // if this is a name we just want to replace
// ASTNodeType_t type = math->getType();
// if (type == AST_NAME)
//
// {
// this->setName(math->getName());
// }
// else
// {
// this->setType(math->getType());
// for (unsigned int i = 0; i < math->getNumChildren(); ++i)
// {
// this->addChild(math->getChild(i)->deepCopy());
// }
// }

// }
//}
}

int
SBMLRateRuleConverter::convert()
{
Expand All @@ -205,6 +273,7 @@ SBMLRateRuleConverter::convert()
{
return returnValue;
}
catchAnomalies();

// Fages algo 3.6 Steps 1-2
populateODEinfo();
Expand Down Expand Up @@ -376,7 +445,12 @@ SBMLRateRuleConverter::determineCoefficient(ASTNode* ode, unsigned int termN, do
// first child should be a number
// take it out of the term
// it will be used as a coefficient
if (term1->getType() == AST_TIMES && term1->getNumChildren() > 0)
if (term1->getType() == AST_REAL)
{
coeff = term1->getValue();
found = true;
}
else if (term1->getType() == AST_TIMES && term1->getNumChildren() > 0)
{
if (term1->getChild(0)->isNumber())
{
Expand Down Expand Up @@ -647,6 +721,7 @@ SBMLRateRuleConverter::addToTerms(ASTNode* node)
delete term;
return;
}
term->setValue(1.0);
}

if (mTerms.size() == 0)
Expand Down
2 changes: 2 additions & 0 deletions src/sbml/conversion/SBMLRateRuleConverter.h
Original file line number Diff line number Diff line change
Expand Up @@ -244,6 +244,8 @@ class LIBSBML_EXTERN SBMLRateRuleConverter : public SBMLConverter
void createReactions();
void removeRules();

void catchAnomalies();


// member variables populated during analysis
pairODEs mODEs;
Expand Down
70 changes: 51 additions & 19 deletions src/sbml/conversion/test/TestSBMLRateRuleConverter.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -905,31 +905,63 @@ START_TEST(test_model_valid_51)
}
END_TEST



START_TEST(test_model_valid_52)
{
ConversionProperties props;
props.addOption("inferReactions", true);

SBMLConverter* converter = new SBMLRateRuleConverter();
converter->setProperties(&props);

std::string filename(TestDataDirectory);
filename += "valid_52_rr.xml";
std::string filename1(TestDataDirectory);
filename1 += "valid_52_bio.xml";

SBMLDocument* d = readSBMLFromFile(filename.c_str());

converter->setDocument(d);
fail_unless(converter->convert() == LIBSBML_OPERATION_SUCCESS);

SBMLDocument* d1 = readSBMLFromFile(filename1.c_str());
std::string out = writeSBMLToStdString(d);
std::string expected = writeSBMLToStdString(d1);

fail_unless(equals(expected.c_str(), out.c_str()));

delete converter;
delete d;
delete d1;
}
END_TEST

Suite *
create_suite_TestSBMLRateRuleConverter (void)
{
Suite *suite = suite_create("SBMLRateRuleConverter");
TCase *tcase = tcase_create("SBMLRateRuleConverter");

/*tcase_add_test(tcase, test_conversion_raterule_converter_invalid);
tcase_add_test(tcase, test_conversion_raterule_converter);
tcase_add_test(tcase, test_conversion_raterule_converter_non_standard_stoichiometry);
tcase_add_test(tcase, test_conversion_raterule_converter_hidden_variable);
tcase_add_test(tcase, test_crash_converter);
tcase_add_test(tcase, test_model);
tcase_add_test(tcase, test_model1);
tcase_add_test(tcase, test_model2);
tcase_add_test(tcase, test_model3);
tcase_add_test(tcase, test_model4);
tcase_add_test(tcase, test_model5);
tcase_add_test(tcase, test_model6);*/
//tcase_add_test(tcase, test_model_valid_01); // this one works
//tcase_add_test(tcase, test_model_valid_02); // this one works
//tcase_add_test(tcase, test_model_valid_03); // this one replaces the kinetic law'b' with the assignment rule for b
//tcase_add_test(tcase, test_model_valid_04); // this one works
// tcase_add_test(tcase, test_model_valid_05); // this one puts the wrong version number
tcase_add_test(tcase, test_model_valid_51); // this one works although puts the kinetic law in a different order
// tcase_add_test(tcase, test_conversion_raterule_converter_invalid);
// tcase_add_test(tcase, test_conversion_raterule_converter);
// tcase_add_test(tcase, test_conversion_raterule_converter_non_standard_stoichiometry);
// tcase_add_test(tcase, test_conversion_raterule_converter_hidden_variable);
// tcase_add_test(tcase, test_crash_converter);
// tcase_add_test(tcase, test_model);
// tcase_add_test(tcase, test_model1);
// tcase_add_test(tcase, test_model2);
// tcase_add_test(tcase, test_model3);
// tcase_add_test(tcase, test_model4);
// tcase_add_test(tcase, test_model5);
// tcase_add_test(tcase, test_model6);
// tcase_add_test(tcase, test_model_valid_01); // this one works
// tcase_add_test(tcase, test_model_valid_02); // this one works
// tcase_add_test(tcase, test_model_valid_03); // this one replaces the kinetic law'b' with the assignment rule for b
// tcase_add_test(tcase, test_model_valid_04); // this one works
tcase_add_test(tcase, test_model_valid_05); // this one works
//tcase_add_test(tcase, test_model_valid_51); // this one works although puts the kinetic law in a different order
// tcase_add_test(tcase, test_model_valid_52); // this one works although puts the kinetic law in a different order



Expand Down
2 changes: 1 addition & 1 deletion src/sbml/conversion/test/test-data/valid_05_rr.xml
Original file line number Diff line number Diff line change
@@ -1,4 +1,4 @@
<?xml version="1.0" encoding="UTF-8"?>!
<?xml version="1.0" encoding="UTF-8"?>
<sbml xmlns="http://www.sbml.org/sbml/level3/version1/core" level="3" version="1">
<model>
<listOfCompartments>
Expand Down
183 changes: 183 additions & 0 deletions src/sbml/conversion/test/test-data/valid_52_bio.xml
Original file line number Diff line number Diff line change
@@ -0,0 +1,183 @@
<?xml version="1.0" encoding="UTF-8"?>
<sbml xmlns="http://www.sbml.org/sbml/level3/version1/core" level="3" version="1">
<model>
<listOfCompartments>
<compartment id="compartmentOne" spatialDimensions="3" size="1" constant="true"/>
</listOfCompartments>
<listOfSpecies>
<species id="SimData_0" name="SimData_0" compartment="compartmentOne" initialConcentration="1" hasOnlySubstanceUnits="false" boundaryCondition="false" constant="false"/>
<species id="SimData_1" name="SimData_1" compartment="compartmentOne" initialConcentration="1" hasOnlySubstanceUnits="false" boundaryCondition="false" constant="false"/>
<species id="SimData_2" name="SimData_2" compartment="compartmentOne" initialConcentration="1" hasOnlySubstanceUnits="false" boundaryCondition="false" constant="false"/>
</listOfSpecies>
<listOfParameters>
<parameter id="I" name="I" value="1" constant="true"/>
<parameter id="k_1" name="k_1" value="4.2416" constant="true"/>
<parameter id="k_2" name="k_2" value="5.9816" constant="true"/>
<parameter id="k_3" name="k_3" value="0.1009" constant="true"/>
<parameter id="k_4" name="k_4" value="1.1549" constant="true"/>
<parameter id="k_5" name="k_5" value="1.3618" constant="true"/>
<parameter id="k_6" name="k_6" value="1.4219" constant="true"/>
<parameter id="k_7" name="k_7" value="0.0051" constant="true"/>
<parameter id="k_8" name="k_8" value="0.0972" constant="true"/>
<parameter id="k_9" name="k_9" value="0.0012" constant="true"/>
<parameter id="k_10" name="k_10" value="56.8583" constant="true"/>
<parameter id="k_11" name="k_11" value="0.0111" constant="true"/>
<parameter id="k_12" name="k_12" value="0.0014" constant="true"/>
</listOfParameters>
<listOfReactions>
<reaction id="J1" reversible="false" fast="false">
<listOfReactants>
<speciesReference species="SimData_0" stoichiometry="1" constant="true"/>
<speciesReference species="SimData_2" stoichiometry="1" constant="true"/>
</listOfReactants>
<listOfProducts>
<speciesReference species="SimData_0" stoichiometry="1" constant="true"/>
<speciesReference species="SimData_2" stoichiometry="2" constant="true"/>
</listOfProducts>
<kineticLaw>
<math xmlns="http://www.w3.org/1998/Math/MathML">
<apply>
<times/>
<ci> SimData_0 </ci>
<ci> SimData_2 </ci>
<apply>
<divide/>
<ci> k_5 </ci>
<apply>
<plus/>
<ci> SimData_2 </ci>
<ci> k_11 </ci>
</apply>
</apply>
</apply>
</math>
</kineticLaw>
</reaction>
<reaction id="J2" reversible="false" fast="false">
<listOfReactants>
<speciesReference species="SimData_1" stoichiometry="1" constant="true"/>
<speciesReference species="SimData_2" stoichiometry="1" constant="true"/>
</listOfReactants>
<listOfProducts>
<speciesReference species="SimData_1" stoichiometry="1" constant="true"/>
</listOfProducts>
<kineticLaw>
<math xmlns="http://www.w3.org/1998/Math/MathML">
<apply>
<times/>
<ci> SimData_1 </ci>
<ci> SimData_2 </ci>
<apply>
<divide/>
<ci> k_6 </ci>
<apply>
<plus/>
<ci> SimData_2 </ci>
<ci> k_12 </ci>
</apply>
</apply>
</apply>
</math>
</kineticLaw>
</reaction>
<reaction id="J3" reversible="false" fast="false">
<listOfReactants>
<speciesReference species="SimData_0" stoichiometry="1" constant="true"/>
<speciesReference species="SimData_1" stoichiometry="1" constant="true"/>
</listOfReactants>
<listOfProducts>
<speciesReference species="SimData_0" stoichiometry="1" constant="true"/>
<speciesReference species="SimData_1" stoichiometry="2" constant="true"/>
</listOfProducts>
<kineticLaw>
<math xmlns="http://www.w3.org/1998/Math/MathML">
<apply>
<times/>
<ci> SimData_0 </ci>
<ci> SimData_1 </ci>
<apply>
<divide/>
<ci> k_3 </ci>
<apply>
<plus/>
<ci> SimData_1 </ci>
<ci> k_9 </ci>
</apply>
</apply>
</apply>
</math>
</kineticLaw>
</reaction>
<reaction id="J4" reversible="false" fast="false">
<listOfReactants>
<speciesReference species="SimData_1" stoichiometry="1" constant="true"/>
</listOfReactants>
<kineticLaw>
<math xmlns="http://www.w3.org/1998/Math/MathML">
<apply>
<times/>
<ci> SimData_1 </ci>
<apply>
<divide/>
<ci> k_4 </ci>
<apply>
<plus/>
<ci> SimData_1 </ci>
<ci> k_10 </ci>
</apply>
</apply>
</apply>
</math>
</kineticLaw>
</reaction>
<reaction id="J5" reversible="false" fast="false">
<listOfReactants>
<speciesReference species="SimData_0" stoichiometry="1" constant="true"/>
</listOfReactants>
<listOfProducts>
<speciesReference species="SimData_0" stoichiometry="2" constant="true"/>
</listOfProducts>
<kineticLaw>
<math xmlns="http://www.w3.org/1998/Math/MathML">
<apply>
<times/>
<ci> k_1 </ci>
<ci> SimData_0 </ci>
<apply>
<divide/>
<ci> I </ci>
<apply>
<plus/>
<ci> SimData_0 </ci>
<ci> k_7 </ci>
</apply>
</apply>
</apply>
</math>
</kineticLaw>
</reaction>
<reaction id="J6" reversible="false" fast="false">
<listOfReactants>
<speciesReference species="SimData_0" stoichiometry="1" constant="true"/>
</listOfReactants>
<kineticLaw>
<math xmlns="http://www.w3.org/1998/Math/MathML">
<apply>
<times/>
<ci> SimData_0 </ci>
<apply>
<divide/>
<ci> k_2 </ci>
<apply>
<plus/>
<ci> SimData_0 </ci>
<ci> k_8 </ci>
</apply>
</apply>
</apply>
</math>
</kineticLaw>
</reaction>
</listOfReactions>
</model>
</sbml>
Loading

0 comments on commit c19d5b0

Please sign in to comment.