From: Rainer M. <ra...@tb...> - 2006-11-15 15:56:14
|
Dear Dr. Fernando, I just realized that you might want to change local parameters given in kinetic laws: No, this is currently not possible, as local parameters are replaced with their values during construction of the odeModel. SOSlib only keeps parameters of the global parameter list as odeModel parameters. You would need to `globalize' local parameters of your model yourself, via libSBML. In the SOSlib file src/odeSolver.c you will find an example how to do that, the function: static int globalizeParameter(Model_t *m, char *id, char *rid) at line 309 of the 1.6.0 release, where id is the parameter ID and rid is the ID of the reaction, where this local parameter is defined and used. We are currently not handling local parameters explicitly because it is much easier to do this via libSBML before passing a model to SOSlib, but we don't want to change an applications' input model. But I will think about offering our globalizeParameter function for the SOSlib API. You could do this by hand, if you want: Just move the function declaration to src/sbmlsolver/odeSolver.h and remove the `static'. Rainer Rainer Machne wrote: > Dear Dr. Fernando, > > yes you can change any parameter or variable (in SBML: global > parameters, species, compartments) between iterative integration runs > and even during an integration (e.g. to handle external events). > > First, you need to get a pointer to structure variableIndex_t for your > parameter from the odeModel structure: > > odeModel = ODEModel_createFromFile("model.xml"); > vi = ODEModel_getVariableIndex(odeModel, "paramX"); > /* where paramX is the SBML ID of the parameter you want to set */ > > Then you can use this variableIndex to get or set the value of this > parameter or variable via the integratorInstance: > > IntegratorInstance_setVariableValue(ii, vi, 55); > > There is one caveat with this function: > > If you run multiple integrations wiht one integratorInstance, you need > to do an > IntegratorInstance_reset(ii); > before the next integration run, and this will reset your initial > conditions and parameters to the original values in the input SBML file! > And thus, currently the setVariableValue function MUST be called AFTER > the reset function! (We hope to provide additional reset functions with > the next release, which keep current or manually reset values.) > > See > http://www.tbi.univie.ac.at/~raim/odeSolver/doc/api/group__variableIndex.html > for functions to get the variableIndex from the odeModel, and > http://www.tbi.univie.ac.at/~raim/odeSolver/doc/api/group__integrator.html#ga25 > for setting variable/parameter values in the integrator. > > See > examples/ChangingIntegratorInstance.c > for an example program which changes variables during an integration. > > Please, keep asking if anything is unclear. > > Thanks, > Rainer > > > > SourceForge.net wrote: > >>Support Requests item #1597040, was opened at 2006-11-15 07:00 >>Message generated for change (Tracker Item Submitted) made by Item Submitter >>You can respond by visiting: >>https://sourceforge.net/tracker/?func=detail&atid=744812&aid=1597040&group_id=139893 >> >>Please note that this message will contain a full copy of the comment thread, >>including the initial issue submission, for this request, >>not just the latest update. >>Category: None >>Group: None >>Status: Open >>Priority: 5 >>Private: No >>Submitted By: Nobody/Anonymous (nobody) >>Assigned to: Nobody/Anonymous (nobody) >>Summary: Altering the parameters of an ODE object. >> >>Initial Comment: >>Is there any way of altering the rate constants once an SBML file has been read in? This is quite important for using evolutionary algorithms, because a significant amount of time is wasted reading the new SBML file into the ODEsolver. >> >>Yours, >>Dr. C. Fernando >>Birmingham University >>UK. >> >> >>---------------------------------------------------------------------- >> >>You can respond by visiting: >>https://sourceforge.net/tracker/?func=detail&atid=744812&aid=1597040&group_id=139893 >> >>------------------------------------------------------------------------- >>Take Surveys. Earn Cash. Influence the Future of IT >>Join SourceForge.net's Techsay panel and you'll get the chance to share your >>opinions on IT & business topics through brief surveys - and earn cash >>http://www.techsay.com/default.php?page=join.php&p=sourceforge&CID=DEVDEV >>_______________________________________________ >>sbmlsolver-devel mailing list >>sbm...@li... >>https://lists.sourceforge.net/lists/listinfo/sbmlsolver-devel > > > > ------------------------------------------------------------------------- > Take Surveys. Earn Cash. Influence the Future of IT > Join SourceForge.net's Techsay panel and you'll get the chance to share your > opinions on IT & business topics through brief surveys - and earn cash > http://www.techsay.com/default.php?page=join.php&p=sourceforge&CID=DEVDEV > _______________________________________________ > sbmlsolver-discuss mailing list > sbm...@li... > https://lists.sourceforge.net/lists/listinfo/sbmlsolver-discuss > > ------------------------------------------------------------------------- > Take Surveys. Earn Cash. Influence the Future of IT > Join SourceForge.net's Techsay panel and you'll get the chance to share your > opinions on IT & business topics through brief surveys - and earn cash > http://www.techsay.com/default.php?page=join.php&p=sourceforge&CID=DEVDEV > _______________________________________________ > sbmlsolver-devel mailing list > sbm...@li... > https://lists.sourceforge.net/lists/listinfo/sbmlsolver-devel |