int main (int argc, char *argv[]){ int i, j; char *model; double time; double printstep; /* libSBML types */ SBMLDocument_t *d; SBMLReader_t *sr; Model_t *m; /* SOSlib types */ SBMLResults_t *results; timeCourse_t *tc; cvodeSettings_t *set; /* parsing command-line arguments */ if (argc < 4 ) { fprintf(stderr, "usage %s sbml-model-file simulation-time time-steps\n", argv[0]); exit(EXIT_FAILURE); } model = argv[1]; time = atof(argv[2]); printstep = atoi(argv[3]); /* parsing the SBML model with libSBML */ sr = SBMLReader_create(); d = SBMLReader_readSBML(sr, model); SBMLReader_free(sr); /* Setting SBML ODE Solver integration parameters */ set = CvodeSettings_create(); CvodeSettings_setTime(set, time, printstep); CvodeSettings_setErrors(set, 1e-9, 1e-4, 1000); /* calling the SBML ODE Solver which returns SBMLResults */ results = SBML_odeSolver(d, set); if ( SolverError_getNum(FATAL_ERROR_TYPE) ) { printf("Integration not sucessful!\n"); SolverError_dumpAndClearErrors(); return(EXIT_FAILURE); } /* now we have the results and can free the inputs */ CvodeSettings_free(set); SBMLDocument_free(d); /* print results */ printf("### RESULTS \n"); SBMLResults_dump(results); SolverError_dumpAndClearErrors(); /* now we can also free the result structure */ SBMLResults_free(results); return (EXIT_SUCCESS); }
SBML_ODESOLVER_API int IntegratorInstance_handleError(integratorInstance_t *engine) { cvodeData_t *data; cvodeSettings_t *opt; int errorCode; if ( SolverError_getNum(ERROR_ERROR_TYPE) == 0 ) return SolverError_getLastCode(WARNING_ERROR_TYPE); errorCode = SolverError_getLastCode(ERROR_ERROR_TYPE); data = engine->data; opt = engine->opt; /* if (om->algebraic) ?? */ /* if (opt->Sensitivity) ?? */ if ( errorCode ) { /* on flag CV_CONV_FAILURE try again, but now with/without generated Jacobian matrix */ if ( errorCode == CV_CONV_FAILURE && data->run == 1 && opt->StoreResults) { SolverError_error(WARNING_ERROR_TYPE, SOLVER_MESSAGE_RERUN_WITH_OR_WO_JACOBIAN, "Rerun with %s Jacobian matrix.", opt->UseJacobian ? "CVODE's internal approximation of the" : "automatically generated"); /* integrate again */ opt->UseJacobian = !opt->UseJacobian; IntegratorInstance_reset(engine); return IntegratorInstance_integrate(engine); } } return errorCode; }
int doit(int argc, char *argv[]) { double i, j ; cvodeSettings_t *settings = CvodeSettings_create(); variableIndex_t *speciesVI, *parameterVI, *parameter2VI; integratorInstance_t *integratorInstance; char *modelStr, *parameterStr, *parameter2Str, *speciesStr; double parameter, parameterEnd, parameterStepSize, parameter2, parameter2End, parameter2StepSize, errorTolerance, relativeErrorTolerance, endtime, initCond, maxDiff, diff; int maximumIntegrationSteps; odeModel_t *model ; if (argc < 11) { fprintf( stderr, "usage %s sbml-model variable parameter1 parameter1-start parameter1-end parameter1-step parameter2 parameter2-start parameter2-end parameter2-step [endtime] [error-tolerance] [relative-error-tolerance] [maximum integration steps]\n", argv[0]); exit(0); } modelStr = argv[1]; speciesStr = argv[2]; parameterStr = argv[3]; parameter = atof(argv[4]); parameterEnd = atof(argv[5]); parameterStepSize = atof(argv[6]); parameter2Str = argv[7]; parameter2 = atof(argv[8]); parameter2End = atof(argv[9]); parameter2StepSize = atof(argv[10]); if (argc > 12) endtime = atof(argv[11]); else endtime = 10000; if (argc > 13) errorTolerance = atof(argv[12]); else errorTolerance = 1e-6; if (argc > 14) relativeErrorTolerance = atof(argv[13]); else relativeErrorTolerance = 1e-4; if (argc > 15) maximumIntegrationSteps = atoi(argv[14]); else maximumIntegrationSteps = 1e9; model = ODEModel_createFromFile(modelStr); RETURN_ON_ERRORS_WITH(1); speciesVI = ODEModel_getVariableIndex(model, speciesStr); parameterVI = ODEModel_getVariableIndex(model, parameterStr); parameter2VI = ODEModel_getVariableIndex(model, parameter2Str); RETURN_ON_ERRORS_WITH(1); /* integrate until steady state */ if ( endtime == 0 ) { /* stop integration upon a steady state */ CvodeSettings_setIndefinitely(settings, 1); CvodeSettings_setSteadyStateThreshold(settings, 1e-9); CvodeSettings_setHaltOnSteadyState(settings, 1); CvodeSettings_setTime(settings, 1, 1); } else { CvodeSettings_setHaltOnSteadyState(settings, 0); CvodeSettings_setIndefinitely(settings, 0); CvodeSettings_setTime(settings, endtime, 1000); } /* absolute tolerance in Cvode integration */ CvodeSettings_setError(settings, errorTolerance); /* relative tolerance in Cvode integration */ CvodeSettings_setRError(settings, relativeErrorTolerance); /* maximum step number for CVode integration */ CvodeSettings_setMxstep(settings, maximumIntegrationSteps); /* doesn't stop integration upon an event */ CvodeSettings_setHaltOnEvent(settings, 0); /* don't Store time course history */ CvodeSettings_setStoreResults(settings, 0); /* compile model */ CvodeSettings_setCompileFunctions(settings, 1); integratorInstance = IntegratorInstance_create(model, settings); printf("set xlabel '%s'\n", ODEModel_getVariableName(model, parameterVI)); printf("set ylabel '%s'\n", ODEModel_getVariableName(model, parameter2VI)); printf("splot '-' using 1:2:3 title '%s' with points pointtype 1 pointsize 1 palette\n", ODEModel_getVariableName(model, speciesVI) ); /* remember initial condition of observed variable */ initCond = IntegratorInstance_getVariableValue(integratorInstance, speciesVI); maxDiff = 0.0; int error = 0 ; int run = 0; for ( run=0; run<2; run++ ) { for (i = parameter; i <= parameterEnd; i += parameterStepSize) { for (j = parameter2; j <= parameter2End; j += parameter2StepSize) { /* add fourth parameter here */ IntegratorInstance_reset( integratorInstance); RETURN_ON_ERRORS_WITH(1); /* for the second run reset the initial condition of the observed variable to a multiple of the maximum observed difference of its value wrt to its initial condition */ if ( run == 1 ) IntegratorInstance_setVariableValue(integratorInstance, speciesVI, fabs(5*(initCond-maxDiff))); IntegratorInstance_setVariableValue(integratorInstance, parameterVI, i); IntegratorInstance_setVariableValue(integratorInstance, parameter2VI, j); /* printf("ic %g\t", IntegratorInstance_getVariableValue(integratorInstance, speciesVI)); */ while(!IntegratorInstance_checkSteadyState(integratorInstance) && !IntegratorInstance_timeCourseCompleted(integratorInstance) ) { IntegratorInstance_integrateOneStep(integratorInstance); if (SolverError_getNum(ERROR_ERROR_TYPE) || SolverError_getNum(FATAL_ERROR_TYPE)) { fprintf(stderr, "ERROR at parameter 1 = %g, parameter 2 = %g\n", i, j); DumpErrors(); error++; } } /* printf("end %g\n", IntegratorInstance_getVariableValue(integratorInstance, speciesVI)); */ /* check whether the final value has largest difference to initial condition */ diff = fabs(initCond - IntegratorInstance_getVariableValue(integratorInstance, speciesVI)); if ( diff > maxDiff ) maxDiff = diff; DumpState(integratorInstance, parameterVI, parameter2VI, speciesVI); } } } printf("end\n"); /* printf("end\n init %g, MAX DIFF %g, abs %g\n", initCond, maxDiff, fabs(initCond-maxDiff)); */ if ( error ) printf("\t%d errors occured\n", error); IntegratorInstance_free(integratorInstance); VariableIndex_free(parameterVI); VariableIndex_free(parameter2VI); VariableIndex_free(speciesVI); ODEModel_free(model); CvodeSettings_free(settings); return 0; }
int main (int argc, char *argv[]){ int i, j; odeModel_t *om; cvodeSettings_t *set; integratorInstance_t *ii; cvodeResults_t *results; variableIndex_t *y, *p; /* Setting SBML ODE Solver integration parameters */ set = CvodeSettings_create(); CvodeSettings_setTime(set, 300.0, 10); CvodeSettings_setErrors(set, 1e-9, 1e-4, 1e9); /* ACTIVATE SENSITIVITY ANALYSIS */ CvodeSettings_setSensitivity(set, 1); /* 0: simultaneous 1: staggered, 2: staggered1 see CVODES user guide for details */ CvodeSettings_setSensMethod(set, 0); CvodeSettings_setJacobian(set, 1); /* for testing only */ /* CvodeSettings_dump(set); */ /* creating the odeModel */ om = ODEModel_createFromFile("MAPK.xml"); /* get a parameter for which we will check sensitivities */ p = ODEModel_getVariableIndex(om, "K1"); /* calling the integrator */ ii = IntegratorInstance_create(om, set); printf("### Printing Sensitivities to %s (%g) on the fly:\n", ODEModel_getVariableName(om, p), IntegratorInstance_getVariableValue(ii, p)); printf("#time "); IntegratorInstance_dumpNames(ii); IntegratorInstance_dumpPSensitivities(ii, p); while( !IntegratorInstance_timeCourseCompleted(ii) ) { if ( !IntegratorInstance_integrateOneStep(ii) ) { break; } IntegratorInstance_dumpPSensitivities(ii, p); } VariableIndex_free(p); if ( SolverError_getNum(FATAL_ERROR_TYPE) ) { printf("Integration not sucessful!\n"); SolverError_dumpAndClearErrors(); return(EXIT_FAILURE); } y = ODEModel_getVariableIndex(om, "MAPK_P"); printf("\nLet's look at a specific ODE variable:\n"); /* print sensitivities again, but now from stored results */ printf("### RESULTS for Sensitivity Analysis for one ODE variable\n"); printf("#time Variable Sensitivity Params...\n"); printf("#time "); printf("%s ", ODEModel_getVariableName(om, y)); for ( j=0; j<ODEModel_getNsens(om); j++ ) { p = ODEModel_getSensParamIndexByNum(om, j); printf("%s ", ODEModel_getVariableName(om, p)); VariableIndex_free(p); } printf("\n"); results = IntegratorInstance_createResults(ii); for ( i=0; i<CvodeResults_getNout(results); i++ ) { printf("%g ", CvodeResults_getTime(results, i)); printf("%g ", CvodeResults_getValue(results, y, i)); for ( j=0; j<ODEModel_getNsens(om); j++ ) { p = ODEModel_getSensParamIndexByNum(om, j); printf("%g ", CvodeResults_getSensitivity(results, y, p, i)); VariableIndex_free(p); } printf("\n"); } /* drawSensitivity(ii->data, "sensitivity", "ps", 0.9); */ p = ODEModel_getVariableIndex(om, "V1"); printf("\nWhat do sensitivities mean? Let's try out!\n\n"); printf("... add 1 to %s: %g + 1 = ", ODEModel_getVariableName(om, p), IntegratorInstance_getVariableValue(ii, p)); CvodeSettings_setSensitivity(set, 0); IntegratorInstance_reset(ii); IntegratorInstance_setVariableValue(ii, p, IntegratorInstance_getVariableValue(ii,p)+1); printf("%g\n", IntegratorInstance_getVariableValue(ii, p)); printf("... and integrate again:\n\n"); CvodeResults_free(results); IntegratorInstance_integrate(ii); results = IntegratorInstance_createResults(ii); /* and print changed results */ printf("#time %s\n", ODEModel_getVariableName(om, y)); for ( i=0; i<CvodeResults_getNout(results); i++ ) { printf("%g ", CvodeResults_getTime(results, i)); printf("%g\n", CvodeResults_getValue(results, y, i)); } printf("\nSee the difference?\n"); printf("Look what happens when the sensitivity changes its sign\n"); printf("between times 180 and 210.\n\n"); VariableIndex_free(y); VariableIndex_free(p); /* now we have the results and can free the inputs */ IntegratorInstance_free(ii); CvodeSettings_free(set); CvodeResults_free(results); ODEModel_free(om); return (EXIT_SUCCESS); }