SBML_ODESOLVER_API cvodeSettings_t *CvodeSettings_clone(cvodeSettings_t *set)
{
  int i;
  cvodeSettings_t *clone;
  ASSIGN_NEW_MEMORY(clone, struct cvodeSettings, NULL);
    
  /* Setting SBML ODE Solver integration parameters */
  CvodeSettings_setErrors(clone, set->Error, set->RError, set->Mxstep);
  CvodeSettings_setSwitches(clone, set->UseJacobian, set->Indefinitely,
			    set->HaltOnEvent, set->SteadyState,
			    set->StoreResults,
			    set->Sensitivity, set->SensMethod);

  CvodeSettings_setMethod(clone, set->CvodeMethod, set->MaxOrder);
  CvodeSettings_setIterMethod(clone, set->IterMethod);
  
  /* Unless indefinite integration is chosen, generate a TimePoints array  */
  if  ( !clone->Indefinitely ) {    
    ASSIGN_NEW_MEMORY_BLOCK(clone->TimePoints,clone->PrintStep+1,double,NULL);
    /* copy TimePoint array */
    for ( i=0; i<=clone->PrintStep; i++ ) {
      clone->TimePoints[i] = set->TimePoints[i];
    }
  }
SBML_ODESOLVER_API cvodeSettings_t *CvodeSettings_createWith(double Time, int PrintStep, double Error, double RError, int Mxstep, int UseJacobian, int Indefinitely, int HaltOnEvent, int SteadyState, int StoreResults, int Sensitivity, int SensMethod)
{

  cvodeSettings_t *set;
  ASSIGN_NEW_MEMORY(set, struct cvodeSettings, NULL);

  /* 1. Setting SBML ODE Solver integration parameters */
  CvodeSettings_setErrors(set, Error, RError, Mxstep);
  /* set non-linear solver defaults (BDF, Newton, max.order 5*/
  set->CvodeMethod = 0;
  set->IterMethod = 0;
  set->MaxOrder = 5;
  CvodeSettings_setSwitches(set, UseJacobian, Indefinitely,
			    HaltOnEvent, SteadyState, StoreResults,
			    Sensitivity, SensMethod);

  /* 2. Setting Requested Time Series */
  /* Unless indefinite integration, generate a TimePoints array  */
  if  ( !Indefinitely ) {    
    /* ... generate default TimePoint array */
    CvodeSettings_setTime(set, Time, PrintStep);
  }
  return set;
}
int
main (int argc, char *argv[]){
  int i, j;
  char *model, *parameter, *reaction;
  double start, end, steps, value;
  double time ;
  double printstep;

  /* libSBML types */
  SBMLDocument_t *d;
  SBMLReader_t *sr;

  /* SOSlib types */
  cvodeSettings_t *set;
  varySettings_t *vs;
  SBMLResultsMatrix_t *resM;
  SBMLResults_t *results;
 
  /* parsing command-line arguments */
  if (argc < 8 ) {
    fprintf(stderr,
	    "usage %s sbml-model-file simulation-time time-steps"
	    " start-value end-value step-number parameter-id"
	    " [optional reaction-id]\n",
            argv[0]);
    exit(EXIT_FAILURE);
  }
  model = argv[1];
  time = atof(argv[2]);
  printstep = atoi(argv[3]);
  
  parameter = argv[7];
  start = atof(argv[4]);
  end = atof(argv[5]);
  steps = atoi(argv[6]);
  
  if ( argc > 8 ) {
      reaction = argv[8];
  }
  else{
    reaction = NULL;
  }
  
  printf("### Varying parameter %s (reaction %s) from %f to %f in %f steps\n",
	 parameter, reaction, start, end, steps);
  
  /* parsing the SBML model with libSBML */
  sr = SBMLReader_create();
  d = SBMLReader_readSBML(sr, model);
  SBMLReader_free(sr);  

  /* Setting SBML ODE Solver integration parameters with default values */
  set = CvodeSettings_create();
  /* resetting the values we need */
  CvodeSettings_setTime(set, time, printstep);
  CvodeSettings_setErrors(set, 1e-18, 1e-10, 10000);
  CvodeSettings_setSwitches(set, 1, 0, 1, 1, 1, 0, 0); 
  CvodeSettings_setSteadyState(set, 1); 

  /* Setting SBML Ode Solver batch integration parameters */
  vs = VarySettings_allocate(1, steps+1);
  VarySettings_addParameter(vs, parameter, reaction, start, end);
  VarySettings_dump(vs);

  /* calling the SBML ODE Solver Batch function,
     and retrieving SBMLResults */
  resM = SBML_odeSolverBatch(d, set, vs);

  if ( resM == NULL ) {
    printf("### Parameter variation not succesful!\n");
    return(0);
  }
    /* we don't need these anymore */
  CvodeSettings_free(set);  
  SBMLDocument_free(d);
  VarySettings_free(vs);

  results = resM->results[0][0];

  for ( i=0; i<resM->i; i++ ) {
    for ( j=0; j<resM->j; j++ ) {
      results = SBMLResultsMatrix_getResults(resM, i, j);
      printf("### RESULTS Parameter %d, Step %d \n", i+1, j+1);
      /* printing results only for species*/
      SBMLResults_dumpSpecies(results);
    }
  }

  /* SolverError_dumpAndClearErrors(); */
  SBMLResultsMatrix_free(resM);
  return (EXIT_SUCCESS);  
}
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;
}
int doit(void)
{
    int i ;
    cvodeSettings_t *settings ;
    variableIndex_t *s1, *s2;
    integratorInstance_t *integratorInstanceA;
    integratorInstance_t *integratorInstanceB;

    odeModel_t *model = ODEModel_createFromFile("events-2-events-1-assignment-l2.xml");
    RETURN_ON_ERRORS_WITH(1);

    assert(ODEModel_hasVariable(model, "S1"));
    assert(ODEModel_hasVariable(model, "S1"));
    assert(!ODEModel_hasVariable(model, "foobar"));

    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 1 will be ignored and Time =
       0.1 will be used as timestep for infinite integration */
    CvodeSettings_setTime(settings, .01, 1);

    /* Setting Cvode Parameters: absolute and relative tolerances and
       maximal internal step */
    CvodeSettings_setErrors(settings, 1e-18, 1e-14, 500);

    
    /* Setting Integration Switches: see documentation or
       example simpleIntegratorInstance.c for details on
       the passed values */
    CvodeSettings_setSwitches(settings, 1, 1, 0, 0, 0, 0);

    integratorInstanceA = IntegratorInstance_create(model, settings);
    RETURN_ON_ERRORS_WITH(1);

    integratorInstanceB = IntegratorInstance_create(model, settings);
    RETURN_ON_ERRORS_WITH(1);

    DumpState(integratorInstanceA, integratorInstanceB, s1, s2);

    for (i=0; i != 500; i++)
    {
        IntegratorInstance_integrateOneStep(integratorInstanceA);
        RETURN_ON_ERRORS_WITH(1);

        IntegratorInstance_integrateOneStep(integratorInstanceB);
        RETURN_ON_ERRORS_WITH(1);

        DumpState(integratorInstanceA, integratorInstanceB, s1, s2);

    }

    IntegratorInstance_free(integratorInstanceA);
    IntegratorInstance_free(integratorInstanceB);
    VariableIndex_free(s1);
    VariableIndex_free(s2);
    ODEModel_free(model);
    CvodeSettings_free(settings);

    return 0;
}
int doit(void)
{
    int i ;
    double value;
    cvodeSettings_t *settings ;
    variableIndex_t *s1, *s2;
    integratorInstance_t *integratorInstanceA;
    integratorInstance_t *integratorInstanceB;

    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");
    RETURN_ON_ERRORS_WITH(1);

    /* 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 1 will be ignored and Time =
       0.1 will be used as timestep for infinite integration */
    CvodeSettings_setTime(settings, .1, 1);

    /* Setting Cvode Parameters: absolute and relative tolerances and
       maximal internal step */
    CvodeSettings_setErrors(settings, 1e-18, 1e-14, 500);

    
    /* Setting Integration Switches: see documentation or
       example simpleIntegratorInstance.c for details on
       the passed values */
    CvodeSettings_setSwitches(settings, 1, 1, 0, 0, 0, 0, 0);

    /* use compiled model */
    CvodeSettings_setCompileFunctions(settings, 0);

    /* Generate two independent integrator instances from the same
       odeModel and cvodeSettings */
    integratorInstanceA = IntegratorInstance_create(model, settings);
    RETURN_ON_ERRORS_WITH(1);

    integratorInstanceB = IntegratorInstance_create(model, settings);
    RETURN_ON_ERRORS_WITH(1);

    DumpState(integratorInstanceA, integratorInstanceB, s1, s2);

    for (i=0; i != 30; i++)
    {

        /* run integrations A and B */
        IntegratorInstance_integrateOneStep(integratorInstanceA);
        RETURN_ON_ERRORS_WITH(1);

        IntegratorInstance_integrateOneStep(integratorInstanceB);
        RETURN_ON_ERRORS_WITH(1);

	/* While variable s1 from integration A is between 1e-15 and
	   1e-16, set variables s1 and s2 in integration B to some
	   value.  This function also takes care of creating and
	   freeing CVODE solver structures when ODE variables are
	   changed!
	*/
	value = IntegratorInstance_getVariableValue(integratorInstanceA, s1);
        if ( value < 1.e-15 && value >  1.e-16 )
        {
            IntegratorInstance_setVariableValue(integratorInstanceB,
						s1,1.5e-15);
            IntegratorInstance_setVariableValue(integratorInstanceB,
						s2,1.5e-15);
        }

        DumpState(integratorInstanceA, integratorInstanceB, s1, s2);

    }

    IntegratorInstance_free(integratorInstanceA);
    IntegratorInstance_free(integratorInstanceB);
    VariableIndex_free(s1);
    VariableIndex_free(s2);
    ODEModel_free(model);
    CvodeSettings_free(settings);

    return 0;
}