int doIt(void)
{
    int i ;
    cvodeSettings_t *settings, *set2;
    variableIndex_t *s1, *s2;
    integratorInstance_t *integratorInstance;

    odeModel_t *model =
      ODEModel_createFromFile("basic-model1-forward-l2.xml");
    
    RETURN_ON_ERRORS_WITH(1)

    s1 = ODEModel_getVariableIndex(model, "S1");
    s2 = ODEModel_getVariableIndex(model, "S2");

    /* Creating settings with default values */
    settings = CvodeSettings_create();
    
    /* Setting end time to .1, number of time steps to 1 and NULL
       instead of an optional predefined time series (double *); due
       to Indefinitely == 1, Printstep 5 will be ignored and Time =
       0.1 will be used as timestep for infinite integration */
    CvodeSettings_setTime(settings, .1, 5);

    /* Setting Cvode Parameters: absolute tolerance, relative
       tolerance and maximal internal step, respectively */
    CvodeSettings_setErrors(settings, 1e-20, 1e-14, 500);
    
    /* Setting Integration Switches */
    CvodeSettings_setMethod(settings, 1, 12);
    CvodeSettings_setJacobian(settings, 1);
    CvodeSettings_setIndefinitely(settings, 1);
    CvodeSettings_setHaltOnEvent(settings, 0);
    CvodeSettings_setSteadyState(settings, 0);
    CvodeSettings_setStoreResults(settings, 0);
    CvodeSettings_setSensitivity(settings, 0);

    /* first integration run */
    printf("\nFIRST INTEGRATION RUN WITH:\n");
    CvodeSettings_dump(settings);
    integratorInstance = IntegratorInstance_create(model, settings);
    RETURN_ON_ERRORS_WITH(1);

    /* print variable (ODE, assignments) and constant names */
    printf("#time  ");
    IntegratorInstance_dumpNames(integratorInstance);
    /* print initial conditions and parameters */
    IntegratorInstance_dumpData(integratorInstance);
    
    for (i=0; i != 12; i++)
    {
        /* setting the next time step, uncomment the following
	   function call to see its effect */
     /* IntegratorInstance_setNextTimeStep(integratorInstance, 0.2*(i+1)); */

        IntegratorInstance_integrateOneStep(integratorInstance);
        RETURN_ON_ERRORS_WITH(1);
	/* print current data */
	IntegratorInstance_dumpData(integratorInstance);
    }

    /* now, let's try again, with different settings */

    set2 = CvodeSettings_create();
    /* as Indefinitely will be set to 0, a finite integration
       to time 0.24 in 6 steps will be run */
    CvodeSettings_setTime(set2, 1.2, 6);
    CvodeSettings_setErrors(set2, 1e-16, 1e-14, 500);
    /* switches can be set all together, same order as above */
    CvodeSettings_setSwitches(set2, 1, 0, 0, 0, 0, 0, 0);
    
    printf("\nNOW, LET'S TRY AGAIN WITH DIFFERENT SETTINGS:\n");
    CvodeSettings_dump(set2);
    
    IntegratorInstance_set(integratorInstance, set2);

    IntegratorInstance_dumpData(integratorInstance);  

    while( !IntegratorInstance_timeCourseCompleted(integratorInstance) ) {

        IntegratorInstance_integrateOneStep(integratorInstance);
        RETURN_ON_ERRORS_WITH(1);

 	IntegratorInstance_dumpData(integratorInstance);
    }
    
    printf("\n\nFINISHED SUCCESSFULLY!\n");
    printf("Please, note the different values e.g. at time 1.2.\n");
    printf("The values for the first run a more exact, due to the much\n");
    printf("lower error tolerances. The error tolerances have to be\n");
    printf("adapted to the ranges of each model!!\n\n");

    IntegratorInstance_free(integratorInstance);
    ODEModel_free(model);
    VariableIndex_free(s1);
    VariableIndex_free(s2);
    CvodeSettings_free(settings);
    CvodeSettings_free(set2);
    
    return 0;
}
Beispiel #2
0
int
main (int argc, char *argv[]){
  int i;
  
  odeModel_t *om;
  cvodeSettings_t *set;
  integratorInstance_t *ii;
  variableIndex_t *p;
  int flag;
   
  /* Setting SBML ODE Solver integration parameters  */
  set = CvodeSettings_create();
  CvodeSettings_setTime(set, 1000, 10);
  CvodeSettings_setErrors(set, 1e-15, 1e-8, 1e4);
  CvodeSettings_setMethod(set, 0, 5);
  /*   CvodeSettings_setStoreResults(set, 0); */
  CvodeSettings_setJacobian(set, 1); /* for testing only */
  CvodeSettings_setCompileFunctions(set, 0); /* for testing only */

 
  /* creating the odeModel */
  om = ODEModel_createFromFile("MAPK.xml");
  ii = IntegratorInstance_create(om, set);
  

  /* ACTIVATE SENSITIVITY ANALYSIS */
  CvodeSettings_setSensitivity(set, 1);
  /* 0: simultaneous 1: staggered, 2: staggered1
     see CVODES user guide for details */
  CvodeSettings_setSensMethod(set, 0);

  /* ACTIVATE ADJOINT ANALYSIS */
  CvodeSettings_setDoAdj(set);
  CvodeSettings_setAdjTime(set, 1000, 100);
  CvodeSettings_setAdjErrors(set, 1e-15, 1e-8);
  CvodeSettings_setnSaveSteps(set, 1000);

  printf("Try 3 integrations with selected parameters/ICs!\n");
  char *sensIDTest[4];  
  sensIDTest[0] = "MAPK";
  sensIDTest[2] = "MAPK_P";
  sensIDTest[1] = "K1";
  sensIDTest[3] = "Ki";
  CvodeSettings_setSensParams(set, sensIDTest, 4);

  fprintf(stdout, "\n\nReading in linear objective function from: 'MAPK.linobjfun'\n");
  fprintf(stdout, "Demonstration of forward/adjoint sensitivity (near) equivalence. \n\n");

  /* Initially, only linear objective is present */
  flag = IntegratorInstance_setLinearObjectiveFunction(ii, "MAPK.linobjfun");
  if (flag!=1)
    return(EXIT_FAILURE);  

  /* reset integrator to new settings */
  IntegratorInstance_reset(ii);
  
  /* get a parameter for which we will check sensitivities */
  p = ODEModel_getVariableIndex(om, "K1");  
   
  i = 0;
  while ( i < 4 ) {

   /*  if ( i == 2) break; */
    /* Set nonlinear objective function after 2 loops  */
    if ( i == 2)
    {
      fprintf(stdout, "\nReading in nonlinear objective now: 'MAPK.objfun'\n\n");
      flag = IntegratorInstance_setObjectiveFunction(ii, "MAPK.objfun");
      if (flag!=1)
	return(EXIT_FAILURE);
    }

    IntegratorInstance_reset(ii);
    
    while( !IntegratorInstance_timeCourseCompleted(ii) )
     if ( !IntegratorInstance_integrateOneStep(ii) )
       break;
    printf("\nIntegration time was %g\n\n",
	 IntegratorInstance_getIntegrationTime(ii));
     

    /*  IntegratorInstance_dumpData(ii); */
    printf("Param default: %s\n", ODEModel_getVariableName(om, p));
    IntegratorInstance_dumpPSensitivities(ii, p);

    flag = IntegratorInstance_CVODEQuad(ii);
    if (flag!=1)
	return(EXIT_FAILURE);

    if ( i < 2)
      fprintf(stdout, "\n### Printing Forward Sensitivities\n");
    else
      fprintf(stdout, "\n### Printing Objective Function (since nonlinear objective is present)\n");  

     
    flag = IntegratorInstance_printQuad(ii, stdout);


    if (flag!=1)
	return(EXIT_FAILURE);

     /* Now go into adjoint phase */   
    IntegratorInstance_resetAdjPhase(ii);
    /* Adjoint phase */
    while( !IntegratorInstance_timeCourseCompleted(ii) )
      if ( !IntegratorInstance_integrateOneStep(ii) )
	break;
    printf("\nIntegration time was %g\n\n",
	   IntegratorInstance_getIntegrationTime(ii));
    
    /* Print out adjoint soln */
    IntegratorInstance_dumpAdjData(ii);

    /* adjoint quadrature */
    flag = IntegratorInstance_CVODEQuad(ii);
     if (flag!=1) 
      return(EXIT_FAILURE);

    fprintf(stdout, "\n### Printing Adjoint Sensitivities: int_0^T <df/dp, psi> dt\n");
    flag = IntegratorInstance_printQuad(ii, stdout);
    if (flag!=1)
	return(EXIT_FAILURE); 


    fprintf(stdout, "\n############# DONE RUN NUMBER %d  #############\n", i); 
    i++;
  }
  
  /*   VariableIndex_free(y); */
  VariableIndex_free(p);
  /* now we have the results and can free the inputs */
  IntegratorInstance_free(ii);
  CvodeSettings_free(set);
  ODEModel_free(om);

  return (EXIT_SUCCESS);  
}
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);  
}
Beispiel #4
0
int
main (int argc, char *argv[])
{
  int neq, nass, nconst;
  char *x[3];
  double x0[3];
  char *ode[2];
  ASTNode_t *ast[2];
  
  double time, rtol, atol;
  double printstep;

  /* SBML model containing events */
  Model_t *events;
  
  /* SOSlib types */
  odeModel_t *om;
  cvodeSettings_t *set;
  integratorInstance_t *ii;
  variableIndex_t *vi;

  /* parsing input ODEs and settings */
  /* time and error settings */
  time = 0.5;
  rtol = atol = 1.0E-9;
  printstep = 10;
  /* variables */
  neq = 2;
  x[0] = "C_cy"; 
  x0[0] = 0.0;
  ode[0] = "((150.0 * (3.8 - (p * D_cy) - (p * C_cy)) * (1.0 - (p * C_cy))) - (9.0 * C_cy))";
  x[1] = "D_cy";
  x0[1] = 9.0;
  ode[1] = "(3.8 - (3.0 * D_cy) - (p * D_cy) - (p * C_cy))";
  /* parameters */
  nconst = 1;
  x[2] = "p";
  x0[2] = 0.2;
  /* assignments */
  nass = 0;
  /* SBML model containing events */
  events = NULL;
  
  /* creating ODE ASTs from strings */
  ast[0] = SBML_parseFormula(ode[0]);
  ast[1] = SBML_parseFormula(ode[1]);

  /* end of input parsing */
  
  /* Setting SBML ODE Solver integration parameters  */
  set = CvodeSettings_createWithTime(time, printstep);
  CvodeSettings_setErrors(set, atol, rtol, 1e4);
  CvodeSettings_setStoreResults(set, 0);

  /* activating default sensitivity */
  CvodeSettings_setSensitivity(set, 1);
  
  /* creating odeModel */
  /* `events' could be an SBML model containing events */
  om = ODEModel_createFromODEs(ast, neq, nass, nconst, x, x0, events);

  /* creating integrator */
  ii = IntegratorInstance_create(om,set);

  /* integrate */
  printf("integrating:\n");
  while( !IntegratorInstance_timeCourseCompleted(ii) )
  {
    if ( !IntegratorInstance_integrateOneStep(ii) )
    {
      SolverError_dump();
      break;
    }
    else
      IntegratorInstance_dumpData(ii);    
  }

  /* Specific data interfaces: */
  /* names are stored in the odeModel and values 
     are independently stored in the integrator: you can
     get different values from different integrator instances
     built from the same model */
  vi = ODEModel_getVariableIndex(om, "C_cy");
  printf("\nVariable %s has final value of %g at time %g\n\n",
	 ODEModel_getVariableName(om, vi),
	 IntegratorInstance_getVariableValue(ii, vi),
	 IntegratorInstance_getTime(ii));
  VariableIndex_free(vi);

  vi = ODEModel_getVariableIndex(om, "p");
  printf("Sensitivies to %s:\n",  ODEModel_getVariableName(om, vi));
  IntegratorInstance_dumpPSensitivities(ii, vi);
  printf("\n\n");
  VariableIndex_free(vi);

  /* finished, free stuff */
  /* free ASTs */
  ASTNode_free(ast[0]);
  ASTNode_free(ast[1]);
  /* free solver structures */
  IntegratorInstance_free(ii);
  ODEModel_free(om);  
  CvodeSettings_free(set);
  
  return (EXIT_SUCCESS);  
}
Beispiel #5
0
int
main (int argc, char *argv[])
{
    int i, j, k;
    int neq, nsens;
    odeModel_t *om;
    odeSense_t *os;
    cvodeSettings_t *set;
    integratorInstance_t *ii;

    variableIndex_t *p1, *p2, *p3, *vi;
    double *weights;

    /* Setting SBML ODE Solver integration parameters  */
    set = CvodeSettings_create();
    CvodeSettings_setTime(set, 30, 10);
    CvodeSettings_setErrors(set, 1e-15, 1e-10, 1e9);
    CvodeSettings_setMethod(set, 0, 5);
    /*   CvodeSettings_setStoreResults(set, 0); */
    CvodeSettings_setJacobian(set, 1);         /* for testing only */
    CvodeSettings_setCompileFunctions(set, 0); /* for testing only */
    /* CvodeSettings_dump(set); */
    CvodeSettings_setFIM(set);                 /* ACTIVATE FIM */

    /* creating the odeModel */
    om = ODEModel_createFromFile("MAPK.xml");
    ii = IntegratorInstance_create(om, set);

    printf("\nFirst try a few integrations without sensitivity\n");
    IntegratorInstance_dumpNames(ii);
    printf("\n");
    for ( i=0; i<2; i++ )
    {
        printf("Run #%d:\n", i);
        IntegratorInstance_integrate(ii);
        /*     IntegratorInstance_dumpData(ii); */
        IntegratorInstance_printResults(ii, stderr);
        IntegratorInstance_reset(ii);
        printf("finished\n\n");
    }

    /* ACTIVATE SENSITIVITY ANALYSIS */

    CvodeSettings_setSensitivity(set, 1);
    /* 0: simultaneous 1: staggered, 2: staggered1
       see CVODES user guide for details */
    CvodeSettings_setSensMethod(set, 0);

    /* reset integrator to new settings */
    IntegratorInstance_reset(ii);

    printf("Now Activate Sensitivity\n\n");

    os = IntegratorInstance_getSensitivityModel(ii);
    nsens = ODESense_getNsens(os);

    printf("nsens  = %i\n\n", nsens);
    for ( i=0; i<nsens; i++ ) {
        vi = ODESense_getSensParamIndexByNum(os, i);
        printf("%s\n", ODEModel_getVariableName(om, vi) );
        VariableIndex_free(vi);
    }

    printf("\n");
    printf("sensitivities calculated for all constants\n");

    p1 = ODESense_getSensParamIndexByNum(os, 1);
    p2 = ODESense_getSensParamIndexByNum(os, 2);
    p3 = ODESense_getSensParamIndexByNum(os, 3);
    printf("sensitivities printed for constants %s, %s, %s\n\n",
           ODEModel_getVariableName(om, p1),
           ODEModel_getVariableName(om, p2),
           ODEModel_getVariableName(om, p3));

    /* create non-default weights for computation of FIM */
    /* weights should be extracted from objective function! */

    neq = ODEModel_getNeq(om);

    ASSIGN_NEW_MEMORY_BLOCK(weights, neq, double, 0);
    for ( i=0; i<neq; i++ )
        weights[i] = 1.;

    /* set weights (to non-default values) */
    IntegratorInstance_setFIMweights(ii, weights, neq);

    /* *** *** *** *** *** *** discrete data *** *** *** *** *** *** *** */
    CvodeSettings_setDiscreteObservation(set);
    printf("DISCRETE DATA\n\n");

    i = 0;
    while ( i < 2 )
    {
        printf("Run #%d:\n", i);
        while( !IntegratorInstance_timeCourseCompleted(ii) )
        {
            IntegratorInstance_dumpPSensitivities(ii, p1);
            IntegratorInstance_dumpPSensitivities(ii, p2);
            IntegratorInstance_dumpPSensitivities(ii, p3);
            if ( !IntegratorInstance_integrateOneStep(ii) )
                break;
        }
        IntegratorInstance_dumpPSensitivities(ii, p1);
        IntegratorInstance_dumpPSensitivities(ii, p2);
        IntegratorInstance_dumpPSensitivities(ii, p3);

        fprintf(stderr, "FIM =\n");
        for ( j=0; j<nsens; j++ )
        {
            for ( k=0; k<nsens; k++ )
                fprintf(stderr, "%g ", IntegratorInstance_getFIM(ii,j,k));
            fprintf(stderr, "\n");
        }
        fprintf(stderr, "\n");

        IntegratorInstance_reset(ii);

        i++;
    }

    /* *** *** *** *** *** *** continuous data *** *** *** *** *** *** *** */
    CvodeSettings_unsetDiscreteObservation(set);
    printf("CONTINUOUS DATA\n\n");

    i = 0;
    while ( i < 2 )
    {
        printf("Run #%d:\n", i);
        while( !IntegratorInstance_timeCourseCompleted(ii) )
        {
            IntegratorInstance_dumpPSensitivities(ii, p1);
            IntegratorInstance_dumpPSensitivities(ii, p2);
            IntegratorInstance_dumpPSensitivities(ii, p3);
            if ( !IntegratorInstance_integrateOneStep(ii) )
                break;
        }
        IntegratorInstance_dumpPSensitivities(ii, p1);
        IntegratorInstance_dumpPSensitivities(ii, p2);
        IntegratorInstance_dumpPSensitivities(ii, p3);

        /* calculate FIM */
        IntegratorInstance_CVODEQuad(ii);

        fprintf(stderr, "FIM =\n");
        for ( j=0; j<nsens; j++ )
        {
            for ( k=0; k<nsens; k++ )
                fprintf(stderr, "%g ", IntegratorInstance_getFIM(ii,j,k));
            fprintf(stderr, "\n");
        }
        fprintf(stderr, "\n");

        IntegratorInstance_reset(ii);

        i++;
    }

    fprintf(stderr, "finished\n");

    /*   VariableIndex_free(y); */
    VariableIndex_free(p1);
    VariableIndex_free(p2);
    VariableIndex_free(p3);
    free(weights);

    /* now we have the results and can free the inputs */
    IntegratorInstance_free(ii);
    CvodeSettings_free(set);
    ODEModel_free(om);

    return (EXIT_SUCCESS);
}