END_TEST

START_TEST ( test_EventAssignment )
{
  EventAssignment* ea = new EventAssignment(2, 4);
  
  fail_unless (!(ea->hasRequiredAttributes()));

  ea->setVariable("ea");

  fail_unless (ea->hasRequiredAttributes());

  delete ea;
}
END_TEST

START_TEST ( test_EventAssignment )
{
  EventAssignment* ea = new EventAssignment(2, 4);
  
  fail_unless (!(ea->hasRequiredElements()));

  ea->setMath(SBML_parseFormula("fd"));

  fail_unless (ea->hasRequiredElements());

  delete ea;
}
END_TEST


START_TEST ( test_EventAssignment_parent_create )
{
    Event *e = new Event(2, 4);

    EventAssignment *ea = e->createEventAssignment();

    ListOf *lo = e->getListOfEventAssignments();

    fail_unless(lo == e->getEventAssignment(0)->getParentSBMLObject());
    fail_unless(lo == ea->getParentSBMLObject());
    fail_unless(e == lo->getParentSBMLObject());
}
END_TEST


START_TEST ( test_Event_parent_NULL )
{
    SBMLDocument *d = new SBMLDocument(2, 4);
    Model *m = d->createModel();
    Event *c = m->createEvent();
    EventAssignment *ea = c->createEventAssignment();
    Trigger *t = new Trigger(2, 4);
    t->setMath(new ASTNode());
    Delay *dy = new Delay(2, 4);
    dy->setMath(new ASTNode());
    c->setTrigger(t);
    c->setDelay(dy);

    fail_unless(c->getAncestorOfType(SBML_MODEL) == m);
    fail_unless(c->getTrigger()->getParentSBMLObject() == c);
    fail_unless (c->getDelay()->getSBMLDocument() == d);
    fail_unless(ea->getAncestorOfType(SBML_EVENT) == c);

    Event *c1 = c->clone();
    delete d;

    fail_unless(c1->getAncestorOfType(SBML_MODEL) == NULL);
    fail_unless(c1->getParentSBMLObject() == NULL);
    fail_unless (c1->getSBMLDocument() == NULL);

    fail_unless(c1->getEventAssignment(0)->getAncestorOfType(SBML_MODEL) == NULL);
    fail_unless(c1->getEventAssignment(0)->getAncestorOfType(SBML_EVENT) == c1);
    fail_unless(c1->getEventAssignment(0)->getParentSBMLObject() != NULL);
    fail_unless(c1->getEventAssignment(0)->getSBMLDocument() == NULL);

    fail_unless(c1->getTrigger()->getAncestorOfType(SBML_MODEL) == NULL);
    fail_unless(c1->getTrigger()->getAncestorOfType(SBML_EVENT) == c1);
    fail_unless(c1->getTrigger()->getParentSBMLObject() != NULL);
    fail_unless(c1->getTrigger()->getSBMLDocument() == NULL);

    fail_unless(c1->getDelay()->getAncestorOfType(SBML_MODEL) == NULL);
    fail_unless(c1->getDelay()->getAncestorOfType(SBML_EVENT) == c1);
    fail_unless(c1->getDelay()->getParentSBMLObject() != NULL);
    fail_unless(c1->getDelay()->getSBMLDocument() == NULL);

    delete c1;
}
END_TEST


START_TEST ( test_EventAssignment_parent_add )
{
    Event *e = new Event(2, 4);
    EventAssignment *ea = new EventAssignment(2, 4);
    ea->setVariable("c");
    ea->setMath(SBML_parseFormula("K+L"));

    e->addEventAssignment(ea);

    delete ea;

    ListOf *lo = e->getListOfEventAssignments();

    fail_unless(lo == e->getEventAssignment(0)->getParentSBMLObject());
    fail_unless(e == lo->getParentSBMLObject());
}
示例#6
0
int Submodel::convertTimeAndExtentWith(const ASTNode* tcf, const ASTNode* xcf, const ASTNode* klmod)
{
  if (tcf==NULL && xcf==NULL) return LIBSBML_OPERATION_SUCCESS;
  Model* model = getInstantiation();
  if (model==NULL) {
    //getInstantiation sets its own error messages.
    return LIBSBML_OPERATION_FAILED;
  }
  ASTNode tcftimes(AST_TIMES);
  ASTNode tcfdiv(AST_DIVIDE);
  if (tcf != NULL) {
    tcftimes.addChild(tcf->deepCopy());
    tcfdiv.addChild(tcf->deepCopy());
  }
  ASTNode rxndivide(AST_DIVIDE);
  if (klmod != NULL) {
    ASTNode rxnref(AST_NAME);
    rxndivide.addChild(rxnref.deepCopy());
    rxndivide.addChild(klmod->deepCopy());
  }
  List* allElements = model->getAllElements();
  for (ListIterator iter = allElements->begin(); iter != allElements->end(); ++iter)
  {
    SBase* element = static_cast<SBase*>(*iter);
    assert(element != NULL);
    ASTNode* ast1 = NULL;
    ASTNode* ast2 = NULL;
    Constraint* constraint = NULL;
    Delay* delay = NULL;
    EventAssignment* ea = NULL;
    InitialAssignment* ia = NULL;
    KineticLaw* kl = NULL;
    Priority* priority = NULL;
    RateRule* rrule = NULL;
    Rule* rule = NULL;
    Submodel* submodel = NULL;
    Trigger* trigger = NULL;
    string cf = "";
    //Reaction math will be converted below, in the bits with the kinetic law.  But because of that, we need to handle references *to* the reaction:  even if it has no kinetic law, the units have changed, and this needs to be reflected by the flattening routine.
    if (rxndivide.getNumChildren() != 0 && element->getTypeCode()==SBML_REACTION && element->isSetId()) {
      rxndivide.getChild(0)->setName(element->getId().c_str());
      for (ListIterator iter = allElements->begin(); iter != allElements->end(); ++iter)
      {
        SBase* subelement = static_cast<SBase*>(*iter);
        subelement->replaceSIDWithFunction(element->getId(), &rxndivide);
      }
    }

    //Submodels need their timeConversionFactor and extentConversionFactor attributes converted.  We're moving top-down, so all we need to do here is fix the conversion factor attributes themselves, pointing them to new parameters if need be.
    if ((tcf !=NULL || xcf != NULL) && element->getTypeCode()==SBML_COMP_SUBMODEL) {
      submodel = static_cast<Submodel*>(element);
      if (tcf != NULL) {
        if (submodel->isSetTimeConversionFactor()) {
          createNewConversionFactor(cf, tcf, submodel->getTimeConversionFactor(), model);
          submodel->setTimeConversionFactor(cf);
        }
        else {
          submodel->setTimeConversionFactor(tcf->getName());
        }
      }
      if (xcf != NULL) {
        if (submodel->isSetExtentConversionFactor()) {
          createNewConversionFactor(cf, xcf, submodel->getExtentConversionFactor(), model);
          submodel->setExtentConversionFactor(cf);
        }
        else {
          submodel->setExtentConversionFactor(xcf->getName());
        }
      }
    }
    if (tcf==NULL) {
      if (klmod !=NULL && element->getTypeCode()==SBML_KINETIC_LAW) {
        kl = static_cast<KineticLaw*>(element);
        if (kl->isSetMath()) {
          ast1 = new ASTNode(AST_TIMES);
          ast1->addChild(klmod->deepCopy());
          ast1->addChild(kl->getMath()->deepCopy());
          kl->setMath(ast1);
          delete ast1;
        }
      }
    }
    else {
      // All math 'time' and 'delay' csymbols must still be converted.
      // Also, several constructs are modified directly.
      switch(element->getTypeCode()) {
        //This would be a WHOLE LOT SIMPLER if there was a 'hasMath' class in libsbml.  But even so, it would have to
        // handle the kinetic laws, rate rules, and delays separately.
      case SBML_KINETIC_LAW:
        //Kinetic laws are multiplied by 'klmod'.
        kl = static_cast<KineticLaw*>(element);
        ast1 = kl->getMath()->deepCopy();
        convertCSymbols(ast1, &tcfdiv, &tcftimes);
        if (klmod !=NULL) {
          kl = static_cast<KineticLaw*>(element);
          if (kl->isSetMath()) {
            ast2 = new ASTNode(AST_TIMES);
            ast2->addChild(klmod->deepCopy());
            ast2->addChild(ast1);
            kl->setMath(ast2);
            delete ast2;
          }
        }
        else {
          kl->setMath(ast1);
          delete ast1;
        }
        break;
      case SBML_DELAY:
        //Delays are multiplied by the time conversion factor.
        delay = static_cast<Delay*>(element);
        if (delay->isSetMath()) {
          ast1 = delay->getMath()->deepCopy();
          convertCSymbols(ast1, &tcfdiv, &tcftimes);
          tcftimes.addChild(ast1);
          delay->setMath(&tcftimes);
          tcftimes.removeChild(1);
          delete ast1;
        }
        break;
      case SBML_RATE_RULE:
        //Rate rules are divided by the time conversion factor.
        rrule = static_cast<RateRule*>(element);
        if (rrule->isSetMath()) {
          ast1 = rrule->getMath()->deepCopy();
          tcfdiv.insertChild(0, ast1);
          rrule->setMath(&tcfdiv);
          tcfdiv.removeChild(0);
          delete ast1;
        }
        //Fall through to:
      case SBML_ASSIGNMENT_RULE:
      case SBML_ALGEBRAIC_RULE:
        //Rules in general need csymbols converted.
        rule = static_cast<Rule*>(element);
        if (rule->isSetMath()) {
          ast1 = rule->getMath()->deepCopy();
          convertCSymbols(ast1, &tcfdiv, &tcftimes);
          rule->setMath(ast1);
          delete ast1;
        }
        break;
      case SBML_EVENT_ASSIGNMENT:
        //Event assignments need csymbols converted.
        ea = static_cast<EventAssignment*>(element);
        if (ea->isSetMath()) {
          ast1 = ea->getMath()->deepCopy();
          convertCSymbols(ast1, &tcfdiv, &tcftimes);
          ea->setMath(ast1);
          delete ast1;
        }
        break;
      case SBML_INITIAL_ASSIGNMENT:
        //Initial assignments need csymbols converted.
        ia = static_cast<InitialAssignment*>(element);
        if (ia->isSetMath()) {
          ast1 = ia->getMath()->deepCopy();
          convertCSymbols(ast1, &tcfdiv, &tcftimes);
          ia->setMath(ast1);
          delete ast1;
        }
        break;
      case SBML_CONSTRAINT:
        //Constraints need csymbols converted.
        constraint = static_cast<Constraint*>(element);
        if (constraint->isSetMath()) {
          ast1 = constraint->getMath()->deepCopy();
          convertCSymbols(ast1, &tcfdiv, &tcftimes);
          constraint->setMath(ast1);
          delete ast1;
        }
        break;
      case SBML_PRIORITY:
        //Priorities need csymbols converted.
        priority = static_cast<Priority*>(element);
        if (priority->isSetMath()) {
          ast1 = priority->getMath()->deepCopy();
          convertCSymbols(ast1, &tcfdiv, &tcftimes);
          priority->setMath(ast1);
          delete ast1;
        }
        break;
      case SBML_TRIGGER:
        //Triggers need csymbols converted.
        trigger = static_cast<Trigger*>(element);
        if (trigger->isSetMath()) {
          ast1 = trigger->getMath()->deepCopy();
          convertCSymbols(ast1, &tcfdiv, &tcftimes);
          trigger->setMath(ast1);
          delete ast1;
        }
        break;
      default:
        //Do nothing!  If we wanted to call a plugin routine, this would be the place.  The only other alternative is to #ifdef some code in here that deals with the math-containing package objects explicitly.  Which might be the best option, all told.
        break;
      }
    }
  }

  delete allElements;

  return LIBSBML_OPERATION_SUCCESS;
}
示例#7
0
int
main (int argc, char* argv[])
{
  if (argc != 2)
  {
    cout << endl << "Usage: printNotes filename" << endl << endl;
    return 1;
  }

  unsigned int i,j;
  const char* filename   = argv[1];
  SBMLDocument* document;
  SBMLReader reader;

  document = reader.readSBML(filename);

  unsigned int errors = document->getNumErrors();

  cout << endl;
  cout << "filename: " << filename << endl;
  cout << endl;

  if(errors > 0)
  {
    document->printErrors(cerr);
    delete document;

    return errors;
  }

  /* Model */

  Model* m = document->getModel();
  printNotes(m);

  for(i=0; i < m->getNumReactions(); i++)
  {
    Reaction* re = m->getReaction(i);
    printNotes(re);

    /* SpeciesReference (Reacatant) */

    for(j=0; j < re->getNumReactants(); j++)
    {
      SpeciesReference* rt = re->getReactant(j);
      if (rt->isSetNotes()) cout << "   ";
      printNotes(rt, (rt->isSetSpecies() ? rt->getSpecies() : std::string("")) );
    }

    /* SpeciesReference (Product) */

    for(j=0; j < re->getNumProducts(); j++)
    {
      SpeciesReference* rt = re->getProduct(j);
      if (rt->isSetNotes()) cout << "   ";
      printNotes(rt, (rt->isSetSpecies() ? rt->getSpecies() : std::string("")) );
    }

    /* ModifierSpeciesReference (Modifier) */

    for(j=0; j < re->getNumModifiers(); j++)
    {
      ModifierSpeciesReference* md = re->getModifier(j);
      if (md->isSetNotes()) cout << "   ";
      printNotes(md, (md->isSetSpecies() ? md->getSpecies() : std::string("")) );
    }

    /* Kineticlaw */

    if(re->isSetKineticLaw())
    {
      KineticLaw* kl = re->getKineticLaw();
      if (kl->isSetNotes()) cout << "   ";
      printNotes(kl);

      /* Parameter */

      for(j=0; j < kl->getNumParameters(); j++)
      {
        Parameter* pa = kl->getParameter(j);
        if (pa->isSetNotes()) cout << "      ";
        printNotes(pa);
      }
    }

  }

  /* Species */

  for(i=0; i < m->getNumSpecies(); i++)
  {
    Species* sp = m->getSpecies(i);
    printNotes(sp);
  }

  /* Compartment */

  for(i=0; i < m->getNumCompartments(); i++)
  {
    Compartment* sp = m->getCompartment(i);
    printNotes(sp);
  }

  /* FunctionDefinition */

  for(i=0; i < m->getNumFunctionDefinitions(); i++)
  {
    FunctionDefinition* sp = m->getFunctionDefinition(i);
    printNotes(sp);
  }

  /* UnitDefinition */

  for(i=0; i < m->getNumUnitDefinitions(); i++)
  {
    UnitDefinition* sp = m->getUnitDefinition(i);
    printNotes(sp);
  }

  /* Parameter */

  for(i=0; i < m->getNumParameters(); i++)
  {
    Parameter* sp = m->getParameter(i);
    printNotes(sp);
  }

  /* Rule */

  for(i=0; i < m->getNumRules(); i++)
  {
    Rule* sp = m->getRule(i);
    printNotes(sp);
  }

  /* InitialAssignment */

  for(i=0; i < m->getNumInitialAssignments(); i++)
  {
    InitialAssignment* sp = m->getInitialAssignment(i);
    printNotes(sp);
  }

  /* Event */

  for(i=0; i < m->getNumEvents(); i++)
  {
    Event* sp = m->getEvent(i);
    printNotes(sp);

    /* Trigger */

    if(sp->isSetTrigger())
    {
      const Trigger* tg = sp->getTrigger();
      if (tg->isSetNotes()) cout << "   ";
      printNotes(const_cast<Trigger*>(tg));
    }

    /* Delay */

    if(sp->isSetDelay())
    {
      const Delay* dl = sp->getDelay();
      if (dl->isSetNotes()) cout << "   ";
      printNotes(const_cast<Delay*>(dl));
    }

    /* EventAssignment */

    for(j=0; j < sp->getNumEventAssignments(); j++)
    {
      EventAssignment* ea = sp->getEventAssignment(j);
      if (ea->isSetNotes()) cout << "   ";
      printNotes(ea);
    }
  }

  /* SpeciesType */

  for(i=0; i < m->getNumSpeciesTypes(); i++)
  {
    SpeciesType* sp = m->getSpeciesType(i);
    printNotes(sp);
  }

  /* Constraint */

  for(i=0; i < m->getNumConstraints(); i++)
  {
    Constraint* sp = m->getConstraint(i);
    printNotes(sp);
  }

  delete document;
  return errors;
}
示例#8
0
void Module::CreateSBMLModel()
{
  Model* sbmlmod = m_sbml.createModel();
  sbmlmod->setId(m_modulename);
  sbmlmod->setName(m_modulename);
  sbmlmod->setNotes("<body xmlns=\"http://www.w3.org/1999/xhtml\"><p> Originally created by libAntimony " VERSION_STRING " (using libSBML " LIBSBML_DOTTED_VERSION ") </p></body>");
  char cc = g_registry.GetCC();
  //User-defined functions
  for (size_t uf=0; uf<g_registry.GetNumUserFunctions(); uf++) {
    const UserFunction* userfunction = g_registry.GetNthUserFunction(uf);
    assert(userfunction != NULL);
    FunctionDefinition* fd = sbmlmod->createFunctionDefinition();
    fd->setId(userfunction->GetModuleName());
    ASTNode* math = parseStringToASTNode(userfunction->ToSBMLString());
    fd->setMath(math);
    delete math;
  }
  //Compartments
  Compartment* defaultCompartment = sbmlmod->createCompartment();
  defaultCompartment->setId(DEFAULTCOMP);
  defaultCompartment->setConstant(true);
  defaultCompartment->setSize(1);
  defaultCompartment->setSBOTerm(410); //The 'implicit compartment'
  size_t numcomps = GetNumVariablesOfType(allCompartments);
  for (size_t comp=0; comp<numcomps; comp++) {
    const Variable* compartment = GetNthVariableOfType(allCompartments, comp);
    Compartment* sbmlcomp = sbmlmod->createCompartment();
    sbmlcomp->setId(compartment->GetNameDelimitedBy(cc));
    if (compartment->GetDisplayName() != "") {
      sbmlcomp->setName(compartment->GetDisplayName());
    }
    sbmlcomp->setConstant(compartment->GetIsConst());
    formula_type ftype = compartment->GetFormulaType();
    assert (ftype == formulaINITIAL || ftype==formulaASSIGNMENT || ftype==formulaRATE);
    if (ftype != formulaINITIAL) {
      sbmlcomp->setConstant(false);
    }
    const Formula* formula = compartment->GetFormula();
    if (formula->IsDouble()) {
      sbmlcomp->setSize(atof(formula->ToSBMLString().c_str()));
    }
    SetAssignmentFor(sbmlmod, compartment);
  }

  //Species
  size_t numspecies = GetNumVariablesOfType(allSpecies);
  for (size_t spec=0; spec < numspecies; spec++) {
    const Variable* species = GetNthVariableOfType(allSpecies, spec);
    Species* sbmlspecies = sbmlmod->createSpecies();
    sbmlspecies->setId(species->GetNameDelimitedBy(cc));
    if (species->GetDisplayName() != "") {
      sbmlspecies->setName(species->GetDisplayName());
    }
    sbmlspecies->setConstant(false); //There's no need to try to distinguish between const and var for species.
    if (species->GetIsConst()) {
      sbmlspecies->setBoundaryCondition(true);
    }
    else {
      sbmlspecies->setBoundaryCondition(false);
    }
    const Variable* compartment = species->GetCompartment();
    if (compartment == NULL) {
      sbmlspecies->setCompartment(defaultCompartment->getId());
    }
    else {
      sbmlspecies->setCompartment(compartment->GetNameDelimitedBy(cc));
    }
    const Formula* formula = species->GetFormula();
    if (formula->IsDouble()) {
      sbmlspecies->setInitialConcentration(atof(formula->ToSBMLString().c_str()));
    }
    else if (formula->IsAmountIn(species->GetCompartment())) {
      sbmlspecies->setInitialAmount(formula->ToAmount());
    }
    SetAssignmentFor(sbmlmod, species);
  }

  //Formulas
  size_t numforms = GetNumVariablesOfType(allFormulas);
  for (size_t form=0; form < numforms; form++) {
    const Variable* formvar = GetNthVariableOfType(allFormulas, form);
    const Formula*  formula = formvar->GetFormula();
    Parameter* param = sbmlmod->createParameter();
    param->setId(formvar->GetNameDelimitedBy(cc));
    if (formvar->GetDisplayName() != "") {
      param->setName(formvar->GetDisplayName());
    }
    param->setConstant(formvar->GetIsConst());
    if (formula->IsDouble()) {
      param->setValue(atof(formula->ToSBMLString().c_str()));
    }
    SetAssignmentFor(sbmlmod, formvar);
    formula_type ftype = formvar->GetFormulaType();
    assert (ftype == formulaINITIAL || ftype==formulaASSIGNMENT || ftype==formulaRATE);
    if (ftype != formulaINITIAL) {
      param->setConstant(false);
    }
  }

  //Reactions
  size_t numrxns = GetNumVariablesOfType(allReactions);
  for (size_t rxn=0; rxn < numrxns; rxn++) {
    const Variable* rxnvar = GetNthVariableOfType(allReactions, rxn);
    const AntimonyReaction* reaction = rxnvar->GetReaction();
    if (reaction->IsEmpty()) {
      continue; //Reactions that involve no species are illegal in SBML.
    }
    Reaction* sbmlrxn = sbmlmod->createReaction();
    sbmlrxn->setId(rxnvar->GetNameDelimitedBy(cc));
    if (rxnvar->GetDisplayName() != "") {
      sbmlrxn->setName(rxnvar->GetDisplayName());
    }
    if (reaction->GetType() == rdBecomes) {
      sbmlrxn->setReversible(true);
    }
    else {
      assert(reaction->GetType() == rdBecomesIrreversibly);
      sbmlrxn->setReversible(false);
    }
    const Formula* formula = reaction->GetFormula();
    string formstring = formula->ToSBMLString(rxnvar->GetStrandVars());
    if (!formula->IsEmpty()) {
      KineticLaw* kl = sbmlmod->createKineticLaw();
      ASTNode* math = parseStringToASTNode(formstring);
      kl->setMath(math);
      delete math;
    }
    const ReactantList* left = reaction->GetLeft();
    for (size_t lnum=0; lnum<left->Size(); lnum++) {
      const Variable* nthleft = left->GetNthReactant(lnum);
      double nthstoich = left->GetStoichiometryFor(lnum);
      SpeciesReference* sr = sbmlmod->createReactant();
      sr->setSpecies(nthleft->GetNameDelimitedBy(cc));
      sr->setStoichiometry(nthstoich);
    }
    const ReactantList* right = reaction->GetRight();
    for (size_t rnum=0; rnum<right->Size(); rnum++) {
      const Variable* nthright = right->GetNthReactant(rnum);
      double nthstoich = right->GetStoichiometryFor(rnum);
      SpeciesReference* sr = sbmlmod->createProduct();
      sr->setSpecies(nthright->GetNameDelimitedBy(cc));
      sr->setStoichiometry(nthstoich);
    }
    //Find 'modifiers' and add them.
    vector<const Variable*> subvars = formula->GetVariablesFrom(formstring, m_modulename);
    for (size_t v=0; v<subvars.size(); v++) {
      if (subvars[v] != NULL && subvars[v]->GetType() == varSpeciesUndef) {
        if (left->GetStoichiometryFor(subvars[v]) == 0 &&
            right->GetStoichiometryFor(subvars[v]) == 0) {
          ModifierSpeciesReference* msr = sbmlmod->createModifier();
          msr->setSpecies(subvars[v]->GetNameDelimitedBy(cc));
        }
      }
    }
  }

  //Events
  size_t numevents = GetNumVariablesOfType(allEvents);
  for (size_t ev=0; ev < numevents; ev++) {
    const Variable* eventvar = GetNthVariableOfType(allEvents, ev);
    const AntimonyEvent* event = eventvar->GetEvent();
    Event* sbmlevent = sbmlmod->createEvent();
    sbmlevent->setId(eventvar->GetNameDelimitedBy(cc));
    if (eventvar->GetDisplayName() != "") {
      sbmlevent->setName(eventvar->GetDisplayName());
    }
    Trigger* trig = sbmlevent->createTrigger();
    ASTNode* ASTtrig = parseStringToASTNode(event->GetTrigger()->ToSBMLString());
    trig->setMath(ASTtrig);
    delete ASTtrig;
    const Formula* delay = event->GetDelay();
    if (!delay->IsEmpty()) {
      ASTtrig = parseStringToASTNode(delay->ToSBMLString());
      Delay* sbmldelay = sbmlevent->createDelay();
      sbmldelay->setMath(ASTtrig);
      delete ASTtrig;
    }
      
    long numasnts = static_cast<long>(event->GetNumAssignments());
    for (long asnt=numasnts-1; asnt>=0; asnt--) {
      //events are stored in reverse order.  Don't ask...
      EventAssignment* sbmlasnt = sbmlmod->createEventAssignment();
      sbmlasnt->setVariable(event->GetNthAssignmentVariableName(asnt, cc));
      ASTNode* ASTasnt = parseStringToASTNode(event->GetNthAssignmentFormulaString(asnt, '_', true));
      sbmlasnt->setMath(ASTasnt);
      delete ASTasnt;
    }
  }

  //Interactions
  size_t numinteractions = GetNumVariablesOfType(allInteractions);
  for (size_t irxn=0; irxn<numinteractions; irxn++) {
    const Variable* arxnvar = GetNthVariableOfType(allInteractions, irxn);
    const AntimonyReaction* arxn = arxnvar->GetReaction();
    Reaction* rxn = sbmlmod->getReaction(arxn->GetRight()->GetNthReactant(0)->GetNameDelimitedBy(cc));
    if (rxn != NULL) {
      for (size_t interactor=0; interactor<arxn->GetLeft()->Size(); interactor++) {
        ModifierSpeciesReference* msr = rxn->createModifier();
        msr->setSpecies(arxn->GetLeft()->GetNthReactant(interactor)->GetNameDelimitedBy(cc));
        msr->setName(arxnvar->GetNameDelimitedBy(cc));
      }
    }
  }

  //Unknown variables (turn into parameters)
  size_t numunknown = GetNumVariablesOfType(allUnknown);
  for (size_t form=0; form < numunknown; form++) {
    const Variable* formvar = GetNthVariableOfType(allUnknown, form);
    Parameter* param = sbmlmod->createParameter();
    param->setId(formvar->GetNameDelimitedBy(cc));
    if (formvar->GetDisplayName() != "") {
      param->setName(formvar->GetDisplayName());
    }
    switch(formvar->GetConstType()) {
    case constVAR:
      param->setConstant(true);
      break;
    case constCONST:
      param->setConstant(false);
      break;
    case constDEFAULT:
      break;
    }
  }
}
示例#9
0
void
CompIdBase::checkId (const EventAssignment& x)
{
  if (x.isSetVariable()) doCheckId(x.getVariable(), x);
}
示例#10
0
LIBSBML_CPP_NAMESPACE_USE

int
main (int argc, char *argv[])
{
    if (argc != 2)
    {
        cout << endl << "Usage: printUnits filename" << endl << endl;
        return 1;
    }

    const char* filename   = argv[1];
    SBMLDocument* document = readSBML(filename);

    if (document->getNumErrors() > 0)
    {
        cerr << "Encountered the following SBML errors:" << endl;
        document->printErrors(cerr);
        return 1;
    }

    Model* model = document->getModel();

    if (model == 0)
    {
        cout << "No model present." << endl;
        return 1;
    }

    unsigned int i,j;
    for (i = 0; i < model->getNumSpecies(); i++)
    {
        Species* s = model->getSpecies(i);
        cout << "Species " << i << ": "
             << UnitDefinition::printUnits(s->getDerivedUnitDefinition()) << endl;
    }

    for (i = 0; i < model->getNumCompartments(); i++)
    {
        Compartment *c = model->getCompartment(i);
        cout << "Compartment " << i << ": "
             << UnitDefinition::printUnits(c->getDerivedUnitDefinition())
             << endl;
    }

    for (i = 0; i < model->getNumParameters(); i++)
    {
        Parameter *p = model->getParameter(i);
        cout << "Parameter " << i << ": "
             << UnitDefinition::printUnits(p->getDerivedUnitDefinition())
             << endl;
    }


    for (i = 0; i < model->getNumInitialAssignments(); i++)
    {
        InitialAssignment *ia = model->getInitialAssignment(i);
        cout << "InitialAssignment " << i << ": "
             << UnitDefinition::printUnits(ia->getDerivedUnitDefinition()) << endl;
        cout << "        undeclared units: ";
        cout << (ia->containsUndeclaredUnits() ? "yes\n" : "no\n");
    }

    for (i = 0; i < model->getNumEvents(); i++)
    {
        Event *e = model->getEvent(i);
        cout << "Event " << i << ": " << endl;

        if (e->isSetDelay())
        {
            cout << "Delay: "
                 << UnitDefinition::printUnits(e->getDelay()->getDerivedUnitDefinition()) << endl;
            cout << "        undeclared units: ";
            cout << (e->getDelay()->containsUndeclaredUnits() ? "yes\n" : "no\n");
        }

        for (j = 0; j < e->getNumEventAssignments(); j++)
        {
            EventAssignment *ea = e->getEventAssignment(j);
            cout << "EventAssignment " << j << ": "
                 << UnitDefinition::printUnits(ea->getDerivedUnitDefinition()) << endl;
            cout << "        undeclared units: ";
            cout << (ea->containsUndeclaredUnits() ? "yes\n" : "no\n");
        }
    }

    for (i = 0; i < model->getNumReactions(); i++)
    {
        Reaction *r = model->getReaction(i);

        cout << "Reaction " << i << ": " << endl;

        if (r->isSetKineticLaw())
        {
            cout << "Kinetic Law: "
                 << UnitDefinition::printUnits(r->getKineticLaw()->getDerivedUnitDefinition()) << endl;
            cout << "        undeclared units: ";
            cout << (r->getKineticLaw()->containsUndeclaredUnits() ? "yes\n" : "no\n");
        }

        for (j = 0; j < r->getNumReactants(); j++)
        {
            SpeciesReference *sr = r->getReactant(j);

            if (sr->isSetStoichiometryMath())
            {
                cout << "Reactant stoichiometryMath" << j << ": "
                     << UnitDefinition::printUnits(sr->getStoichiometryMath()->getDerivedUnitDefinition()) << endl;
                cout << "        undeclared units: ";
                cout << (sr->getStoichiometryMath()->containsUndeclaredUnits() ? "yes\n" : "no\n");
            }
        }

        for (j = 0; j < r->getNumProducts(); j++)
        {
            SpeciesReference *sr = r->getProduct(j);

            if (sr->isSetStoichiometryMath())
            {
                cout << "Product stoichiometryMath" << j << ": "
                     << UnitDefinition::printUnits(sr->getStoichiometryMath()->getDerivedUnitDefinition()) << endl;
                cout << "        undeclared units: ";
                cout << (sr->getStoichiometryMath()->containsUndeclaredUnits() ? "yes\n" : "no\n");
            }
        }
    }

    for (i = 0; i < model->getNumRules(); i++)
    {
        Rule *r = model->getRule(i);
        cout << "Rule " << i << ": "
             << UnitDefinition::printUnits(r->getDerivedUnitDefinition()) << endl;
        cout << "        undeclared units: ";
        cout << (r->containsUndeclaredUnits() ? "yes\n" : "no\n");
    }

    delete document;
    return 0;
}