Collaboration diagram for ODE Model: f(x,p,t) = dx/dt:
![]() |
Modules | |
Jacobian Matrix: J = df(x)/dx | |
Constructing and Interfacing the Jacobian matrix of an ODE system. | |
Sensitivity Matrix: P = df(x)/dp | |
Constructing and Interfacing the sensitivity matrix of an ODE system. | |
Variables + Parameters | |
Getting the variableIndex structure. | |
Functions | |
SBML_ODESOLVER_API odeModel_t * | ODEModel_create (Model_t *m) |
Create internal model odeModel from a reaction network, represented as libSBML's Model_t structure. | |
SBML_ODESOLVER_API odeModel_t * | ODEModel_createFromFile (char *sbmlFileName) |
Create internal model odeModel from an SBML file, that contains level 1 or level 2 SBML. | |
SBML_ODESOLVER_API odeModel_t * | ODEModel_createFromSBML2 (SBMLDocument_t *d) |
Create internal model odeModel_t from SBMLDocument containing a level 2 SBML model. | |
SBML_ODESOLVER_API void | ODEModel_free (odeModel_t *om) |
Frees the odeModel structures. | |
SBML_ODESOLVER_API int | ODEModel_hasVariable (odeModel_t *model, const char *symbol) |
Returns 1 if a variable or parameter with the SBML id exists in the ODEModel. | |
SBML_ODESOLVER_API int | ODEModel_getNumValues (odeModel_t *om) |
Returns the total number of values in oodeModel, equivalent to ODEModel_getNeq + ODEModel_getNumAssignments + ODEModel_getNumConstants. | |
SBML_ODESOLVER_API const char * | ODEModel_getVariableName (odeModel_t *om, variableIndex_t *vi) |
Returns the name of the variable corresponding to passed variableIndex. | |
SBML_ODESOLVER_API int | ODEModel_getNeq (odeModel_t *om) |
SBML_ODESOLVER_API int | ODEModel_getNsens (odeModel_t *om) |
Returns the number parameters for which sensitivity analysis might be requested. | |
SBML_ODESOLVER_API const ASTNode_t * | ODEModel_getOde (odeModel_t *om, variableIndex_t *vi) |
Returns the ODE from the odeModel for the variable with variableIndex vi. | |
SBML_ODESOLVER_API int | ODEModel_getNumAssignments (odeModel_t *om) |
Returns the number of variable assignments in the odeModel. | |
SBML_ODESOLVER_API const ASTNode_t * | ODEModel_getAssignment (odeModel_t *om, variableIndex_t *vi) |
Returns the assignment formula from the odeModel for the variable with variableIndex vi. | |
SBML_ODESOLVER_API int | ODEModel_getNumConstants (odeModel_t *om) |
Returns the number of constant parameters of the odeModel. | |
SBML_ODESOLVER_API int | ODEModel_getNalg (odeModel_t *om) |
Returns the number variables that are defined by an algebraic rule. | |
SBML_ODESOLVER_API void | ODEModel_dumpNames (odeModel_t *om) |
Prints the names (SBML IDs) of all model variables and parameters. | |
SBML_ODESOLVER_API const Model_t * | ODEModel_getModel (odeModel_t *om) |
Returns the SBML model that has been extracted from the input SBML model's reaction network and structures. |
The internal ODE Model (structure odeModel) can be interfaced for analytical purposes. All formulae can be retrieved as libSBML Abstract Syntax Trees (AST).
|
Create internal model odeModel from a reaction network, represented as libSBML's Model_t structure. The input model must be of SBML level2! The function at first, attempts to construct a simplified SBML model, that contains all compartments, species, parameters, events and rules of the input model, and constructs new ODEs as SBML RateRules from the reaction network of the input model. The function then creates the structure odeModel which contains variable, parameter and constant names and all formulas (ODEs and assignments) as indexed AST (iAST). This structure can be used to initialize and run several integration runs, each associated with initial conditions in cvodeData_t. Alternatively I.1a - I.1c allow to construct odeModel from higher-level data (a file, an SBML document or a reaction network model, respectively). |
|
Create internal model odeModel from an SBML file, that contains level 1 or level 2 SBML. Conversion of level 1 to level 2 models is done internally. |
|
Create internal model odeModel_t from SBMLDocument containing a level 2 SBML model.
|
|
Frees the odeModel structures.
|
|
Returns 1 if a variable or parameter with the SBML id exists in the ODEModel.
|
|
Returns the total number of values in oodeModel, equivalent to ODEModel_getNeq + ODEModel_getNumAssignments + ODEModel_getNumConstants.
|
|
Returns the name of the variable corresponding to passed variableIndex. The returned string (const char *) may NOT be changed or freed by calling applications. |
|
|
|
Returns the number parameters for which sensitivity analysis might be requested.
|
|
Returns the ODE from the odeModel for the variable with variableIndex vi. The ODE is returned as an `indexed abstract syntax tree' (iAST), which is an extension to the usual libSBML AST. Every AST_NAME type node in the tree has been replaced by an ASTIndexNameNode, that allows a O(1) retrieval of values for this node from an array of all values of the odeModel. |
|
Returns the number of variable assignments in the odeModel.
|
|
Returns the assignment formula from the odeModel for the variable with variableIndex vi. The ODE is returned as an `indexed abstract syntax tree' (iAST), which is an extension to the usual libSBML AST. Every AST_NAME type node in the tree has been replaced by an ASTIndexNameNode, that allows a O(1) retrieval of value for this node from an array of all values of the odeModel. |
|
Returns the number of constant parameters of the odeModel.
|
|
Returns the number variables that are defined by an algebraic rule. As SBML Algebraic Rules and ODE models with algebraic constraints (DAE models) can currently not be handled, this function is of no use. |
|
Prints the names (SBML IDs) of all model variables and parameters.
|
|
Returns the SBML model that has been extracted from the input SBML model's reaction network and structures. The model contains only compartments, species, parameters and SBML Rules. Changes to this model will have no effect on odeModel or integratorInstance. |