Beispiel #1
0
/**
   returns a pointer to the function named 'symbol' in the given 'code'
*/
void *CompiledCode_getFunction(compiled_code_t *code, const char *symbol)
{
  void *result = NULL;

#ifdef WIN32
  
  result = GetProcAddress(code->dllHandle, symbol);

  if (result)
    return result ;

  SolverError_storeLastWin32Error("");
  result = NULL;

#else /* default case gcc */

  char* returnvalue = NULL;

  /* Clear any existing error */
  returnvalue = dlerror();
  result = dlsym(code->dllHandle, symbol);

  returnvalue = dlerror();
  if ( returnvalue != NULL )
    SolverError_error(FATAL_ERROR_TYPE, SOLVER_ERROR_DL_SYMBOL_UNDEFINED,
		      "dlsym(): couldn't get symbol %s from shared library %s",
		      symbol, code->dllFileName);    

#endif /* end WIN32 */
  
  return (result);
}
SBML_ODESOLVER_API int IntegratorInstance_checkTrigger(integratorInstance_t *engine)
{  
    int i, j, fired;
    ASTNode_t *trigger, *assignment;
    Event_t *e;
    EventAssignment_t *ea;
    variableIndex_t *vi;

    cvodeSettings_t *opt = engine->opt;
    cvodeData_t *data = engine->data;
    odeModel_t *om = engine->om;

    fired = 0;

    for ( i=0; i<Model_getNumEvents(om->simple); i++ ) {
      e = Model_getEvent(om->simple, i);
      trigger = (ASTNode_t *) Event_getTrigger(e);
      if ( data->trigger[i] == 0 && evaluateAST(trigger, data) ) {

	if (opt->HaltOnEvent)
	  SolverError_error(ERROR_ERROR_TYPE,
			    SOLVER_ERROR_EVENT_TRIGGER_FIRED,
			    "Event Trigger %d (%s) fired at time %g. "
			    "Aborting simulation.",
			    i, SBML_formulaToString(trigger),
			    data->currenttime);
    /* removed AMF 08/11/05
	else 
	  SolverError_error(WARNING_ERROR_TYPE,
			    SOLVER_ERROR_EVENT_TRIGGER_FIRED,
			    "Event Trigger %d (%s) fired at time %g. ",
			    i, SBML_formulaToString(trigger),
			    data->currenttime); */
	fired++;
	data->trigger[i] = 1;      
	for ( j=0; j<Event_getNumEventAssignments(e); j++ ) {
	  ea = Event_getEventAssignment(e, j);
	  assignment = (ASTNode_t *) EventAssignment_getMath(ea);
	  vi = ODEModel_getVariableIndex(om,
					 EventAssignment_getVariable(ea));
	  IntegratorInstance_setVariableValue(engine, vi,
					      evaluateAST(assignment, data));
	  VariableIndex_free(vi);
	}
      }
      else {
	data->trigger[i] = 0;
      }
    }

    return fired;

}
Beispiel #3
0
/*
  force to load an SBML file
*/
static SBMLDocument_t *loadFile()
{
    char *filename;
    SBMLDocument_t *d;

    while(1){
        printf("Please enter a filename: ");
        filename = get_line(stdin);
        filename = util_trim(filename);
        if ( (strlen(filename)) == 0 ) {
            printf("No filename found.\n\n");
        }
        else {
            if ( (d = parseModelWithArguments(filename)) == 0 ) {
                if ( Opt.Validate ) 
                    SolverError_error(
                        WARNING_ERROR_TYPE,
                        SOLVER_ERROR_MAKE_SURE_SCHEMA_IS_ON_PATH,
                        "Please make sure that path >%s< contains"
                        "the correct SBML schema for validation."
                        "Or try running without validation.", Opt.SchemaPath);

                SolverError_error(
                    ERROR_ERROR_TYPE,
                    SOLVER_ERROR_CANNOT_PARSE_MODEL,
                    "Can't parse Model >%s<", filename);
                SolverError_dumpAndClearErrors();
            }
            else {
                printf("SBML file %s successfully loaded.\n", filename);
                return d;
            }
        }
        free(filename);
    }
    return NULL;
}
SBML_ODESOLVER_API void IntegratorInstance_copyVariableState(integratorInstance_t *target, integratorInstance_t *source)
{
    int i;
    cvodeData_t *targetData = target->data;
    cvodeData_t *sourceData = source->data;
    odeModel_t *model = target->om;

    if (model == source->om)
    {
        for ( i=0; i<sourceData->nvalues; i++ )
            targetData->value[i] = sourceData->value[i];
    }
    else
        SolverError_error(
            ERROR_ERROR_TYPE,
            SOLVER_ERROR_ATTEMPTING_TO_COPY_VARIABLE_STATE_BETWEEN_INSTANCES_OF_DIFFERENT_MODELS,
            "Attempting to copy variable state between instances of "
	    "different models");
}
/* NOTE: provisional steady state finding! */
SBML_ODESOLVER_API int IntegratorInstance_checkSteadyState(integratorInstance_t *engine)
{
  int i;
  double dy_mean, dy_var, dy_std;
  cvodeData_t *data = engine->data;
  odeModel_t *om = engine->om;
  
  /* calculate the mean and standard deviation of rates of change and
     store in cvodeData_t * */
  dy_mean = 0.0;
  dy_var = 0.0;
  dy_std = 0.0;
  
  for ( i=0; i<om->neq; i++ ) {
    dy_mean += fabs(evaluateAST(om->ode[i],data));
  }
  dy_mean = dy_mean / om->neq;
  for ( i=0; i<om->neq; i++ ) {
    dy_var += MySQR(evaluateAST(om->ode[i],data) - dy_mean);
  }
  dy_var = dy_var / (om->neq -1);
  dy_std = MySQRT(dy_var);

  /* stop integrator if mean + std of rates of change are lower than
     1e-11 */
  if ( (dy_mean + dy_std) < 1e-11 ) {
    data->steadystate = 1;
    SolverError_error(WARNING_ERROR_TYPE,
		      SOLVER_MESSAGE_STEADYSTATE_FOUND,
		      "Steady state found. "
		      "Simulation aborted at %g seconds. "
		      "Mean of rates: %g, std %g",
		      data->currenttime, dy_mean, dy_std);
    return(1) ;
  }
  else {
    data->steadystate = 0;
    return(0);
  }  
}
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;
}
Beispiel #7
0
compiled_code_t *Compiler_compile_win32_tcc(const char *sourceCode)
{
  compiled_code_t *code = NULL;
  char tempDir[MAX_PATH+1];
  TCHAR tccFileName[MAX_PATH]; 
  int i;
  int result;
  char *cFileName;
  char *dllFileName;
  char *outFileName;
  FILE *cFile;
  char command[4*MAX_PATH];
  char *dllFileNameDot ;
  HMODULE dllHandle, solverHandle ;
#ifdef _DEBUG
  char *solverFileName = "SBML_odeSolverD.dll";
#else
  char *solverFileName = "SBML_odeSolver.dll";
#endif

  /*printf("Source code:\n%s\n", sourceCode);*/

  /* avoid creating files in current directory if environment
     variables not set */
  if (!GetTempPath(MAX_PATH+1, tempDir))
  {
    SolverError_storeLastWin32Error("Trying to find out location of system temp directory");
    return NULL;
  }
  
  solverHandle = GetModuleHandle(solverFileName);
    
  if (!solverHandle)
  {
    SolverError_storeLastWin32Error("Trying to get handle of solver dll");
    return NULL;
  }
  
  /* compute tcc path from the path to this dll */
  if( !GetModuleFileName( solverHandle, tccFileName, MAX_PATH ) )
  {
    SolverError_storeLastWin32Error("Trying find location of the soslib dll");
    return NULL ;
  }

  for (i = strlen(tccFileName); i != -1 && tccFileName[i] != '\\'; i--);

  tccFileName[i + 1] = '\0';
  strcat(tccFileName, "tcc\\tcc.exe");

  cFileName = tempnam(tempDir, "temp_soslib_c_file");
  dllFileName = tempnam(tempDir, "temp_soslib_dll");
  outFileName = tempnam(tempDir, "temp_soslib_compilation_output");
  cFile = fopen(cFileName, "w");

  if (!cFile)
  {
    SolverError_storeLastWin32Error("Unable to open C source file for write");
    return NULL;
  }

  fprintf(cFile, "%s", sourceCode);
  fclose(cFile);

  sprintf(command, "%s -o %s -shared %s > %s",
	  tccFileName, dllFileName, cFileName, outFileName);

  /*printf("Command:\n%s\n", command);*/

  result = system(command);

  if (result == -1)
  {
    SolverError_storeLastWin32Error("Whilst running compile command");
    remove(cFileName);
    free(cFileName);
    return NULL ;
  }
  else if (result != 0)
  {
    SolverError_error(ERROR_ERROR_TYPE, SOLVER_ERROR_COMPILATION_FAILED,
		      "Compile command failed - returned %d", result);
    remove(cFileName);
    free(cFileName);
    return NULL ;
  }
  
  remove(cFileName);
  free(cFileName);
  remove(outFileName);
  free(outFileName);

  ASSIGN_NEW_MEMORY_BLOCK(dllFileNameDot, (strlen(dllFileName) +2), char, NULL);
  
  strcpy(dllFileNameDot, dllFileName);
  strcat(dllFileNameDot, ".");

  dllHandle = LoadLibrary(dllFileNameDot);
  free(dllFileNameDot);

  if (!dllHandle)
  {
    SolverError_storeLastWin32Error("While loading compiled dll");
    return NULL;
  }
  
  ASSIGN_NEW_MEMORY(code, compiled_code_t, NULL);

  code->dllHandle = dllHandle ;
  code->dllFileName = dllFileName;

  return (code);
}
Beispiel #8
0
/**
   Returns a pointer to code that is compiled from the given source code
*/
compiled_code_t *Compiler_compile_with_gcc(const char *sourceCode)
{
  compiled_code_t *code = NULL;
  char gccFileName[MAX_PATH+1] = "g++";
  int result;
  char *tmpFileName = NULL;
  char *cFileName   = NULL;
  char *dllFileName = NULL;
  char *oFileName   = NULL;
  FILE *cFile;
  char command[4*MAX_PATH];
  void *dllHandle;
  
  /* generate a unique temprorary filename template */
  ASSIGN_NEW_MEMORY_BLOCK(tmpFileName, (MAX_PATH+1), char, NULL);
  tmpFileName = tmpnam(tmpFileName);
  
#ifdef _DEBUG
  Warn(NULL,"Temporary File Name is %s\n", tmpFileName);
#endif
  
  /* generate needed file names from the template*/
  ASSIGN_NEW_MEMORY_BLOCK(cFileName, (strlen(tmpFileName)+3), char, NULL);
  strcpy(cFileName, tmpFileName);
  strcat(cFileName, ".c");
  ASSIGN_NEW_MEMORY_BLOCK(oFileName, (strlen(tmpFileName)+3), char, NULL);
  strcpy(oFileName, tmpFileName);
  strcat(oFileName, ".o");  
  ASSIGN_NEW_MEMORY_BLOCK(dllFileName,
			  (strlen(tmpFileName)+strlen(SHAREDLIBEXT)+1),
			  char, NULL);
  strcpy(dllFileName, tmpFileName);
  strcat(dllFileName, SHAREDLIBEXT);

  /* open file and dump source code to it */
  cFile = fopen(cFileName, "w");
  
  if (!cFile)
  {
    SolverError_error(WARNING_ERROR_TYPE, SOLVER_ERROR_OPEN_FILE,
		      "Could not open file %s - %s!",
		      cFileName, strerror(errno));  
    free(cFileName);
    free(oFileName);
    free(dllFileName);
    return NULL;
  }

  fprintf(cFile, "%s", sourceCode);
  fclose(cFile);

  /* construct command for compiling */
#if defined (__APPLE__) && defined (__MACH__)
  sprintf(command,
	  "%s -I%s -I%s %s -I../src -pipe -O -dynamiclib -fPIC -o %s %s -L../src -L%s -L%s -L%s -lODES -lsbml -lm",
 	  gccFileName,
	  SOSLIB_CPPFLAGS, /* changed order: SOSLIB first */
	  SUNDIALS_CPPFLAGS,
	  SBML_CPPFLAGS,
	  dllFileName,
	  cFileName,
	  SUNDIALS_LDFLAGS,
	  SBML_LDFLAGS,
	  SOSLIB_LDFLAGS);
#else
  sprintf(command,
	  "%s -I%s -I%s %s -I../src -pipe -O -shared -fPIC -o %s %s -L../src -L%s -L%s -L%s -lODES -lsbml -lm",
 	  gccFileName,
	  SOSLIB_CPPFLAGS, /* changed order: SOSLIB first */
	  SUNDIALS_CPPFLAGS,
	  SBML_CPPFLAGS,
	  dllFileName,
	  cFileName,
	  SUNDIALS_LDFLAGS,
	  SBML_LDFLAGS,
	  SOSLIB_LDFLAGS);
#endif
 
#ifdef _DEBUG
  Warn(NULL, "Command: %s\n", command);
#endif
  
  /* compile source to shared library */
  result = system(command);

  /* clean up compilation intermediates */
  free(tmpFileName);
  remove(cFileName);
  free(cFileName);
  remove(oFileName);
  free(oFileName);

  /* handle possible errors */
  if (result != 0)
  {
    if (result == -1)
      SolverError_error(WARNING_ERROR_TYPE, SOLVER_ERROR_GCC_FORK_FAILED,
			"forking gcc compiler subprocess failed!");
    else 
      SolverError_error(WARNING_ERROR_TYPE, SOLVER_ERROR_COMPILATION_FAILED,
			"compiling failed with errno %d - %s!",
			result, strerror(result));    
    return (NULL);
  }
  


  /* load shared library */
  dllHandle = dlopen(dllFileName, RTLD_LAZY);
  if (dllHandle == NULL)
  {
    SolverError_error(WARNING_ERROR_TYPE, SOLVER_ERROR_DL_LOAD_FAILED,
		      "loading shared library %s failed %d - %s!",
		      dllFileName, errno, strerror(errno));
    SolverError_dumpAndClearErrors();
    return (NULL);
  }

  ASSIGN_NEW_MEMORY(code, compiled_code_t, NULL);
  code->dllHandle   = dllHandle;
  code->dllFileName = dllFileName;

  return (code);
}
Beispiel #9
0
SBML_ODESOLVER_API double evaluateAST(ASTNode_t *n, cvodeData_t *data)
{
  int i, j, childnum;
  int found, datafound;
  int true;
  time_series_t *ts=data->model->time_series;
  double findtol=1e-5;
 
  ASTNodeType_t type;
  /* ASTNode_t **child; */

  /* value* are for intermediate values. */
  double value1, value2, value3, result;

  if ( n == NULL )
  {
    SolverError_error(FATAL_ERROR_TYPE,
		      SOLVER_ERROR_AST_UNKNOWN_NODE_TYPE,
		      "evaluateAST: empty Abstract Syntax Tree (AST).");
    return (0);
  }
  if ( ASTNode_isUnknown(n) )
  {
    SolverError_error(FATAL_ERROR_TYPE,
		      SOLVER_ERROR_AST_UNKNOWN_NODE_TYPE,
		      "evaluateAST: unknown ASTNode type");
  }
  result = 0;

  childnum = ASTNode_getNumChildren(n);
  type = ASTNode_getType(n);
  switch(type)
  {
  case AST_INTEGER:
    result = (double) ASTNode_getInteger(n);      
    break;
  case AST_REAL:
    result = ASTNode_getReal(n);
    break;
  case AST_REAL_E:
    result = ASTNode_getReal(n);
    break;      
  case AST_RATIONAL:
    result = ASTNode_getReal(n) ;
    break;
	
  case AST_NAME:
    /** VARIABLES:

    find the value of the variable in the data->value
    array. SOSlib's extension to libSBML's AST allows to add the
    index of the variable in this array to AST_NAME
    (ASTIndexName). If the ASTNode is not indexed, its array
    index is searched via the data->model->names array, which
    corresponds to the data->value array. For nodes with name
    `Time' or `time' the data->currenttime is returned.  If no
    value is found a fatal error is produced. */
    found = 0;
    if ( ASTNode_isSetIndex(n) )
    {
      if ( ASTNode_isSetData(n) )
      {

        /* if continuous data is observed, obtain interpolated result */  
        if ( (data->model->discrete_observation_data != 1) || (data->model->compute_vector_v != 1) )
	{
	  result = call(ASTNode_getIndex(n),
			data->currenttime, ts);

	}
	else  /* if discrete data is observed, simply obtain value from time_series */
	{
          datafound = 0;
          i = data->TimeSeriesIndex;
        
	    if ( fabs(data->currenttime - ts->time[i]) < findtol )
	    {    
	      result = ts->data[ASTNode_getIndex(n)][i];
              datafound++;
	    }
	
          if ( datafound != 1)
	  {	            
	    SolverError_error(FATAL_ERROR_TYPE,
			SOLVER_ERROR_AST_EVALUATION_FAILED_DISCRETE_DATA,
			"use of discrete time series data failed; none or several time points matching current time");
	    result = 0;
           /*  break;  */
	  }
          else found = 1;

	}
      }
      else
      {
	/* majority case: just return the
	   according value from data->values
	   from the index stored by SOSlib
	   ASTIndexNameNode sub-class of libSBML's ASTNode */
	result = data->value[ASTNode_getIndex(n)];
  
      }
      
       found++;
    }

    if ( found == 0 )
    {
      for ( j=0; j<data->nvalues; j++ )
      {
	if ( (strcmp(ASTNode_getName(n),data->model->names[j]) == 0) )
	{
	  
	  result = data->value[j];
	  found++;
	}
      }
    }

    if ( found == 0 )
    {
      SolverError_error(FATAL_ERROR_TYPE,
			SOLVER_ERROR_AST_EVALUATION_FAILED_MISSING_VALUE,
			"No value found for AST_NAME %s . Defaults to Zero "
			"to avoid program crash", ASTNode_getName(n));
      result = 0;
    }
    break;

  case AST_FUNCTION_DELAY:
    SolverError_error(FATAL_ERROR_TYPE,
		      SOLVER_ERROR_AST_EVALUATION_FAILED_DELAY,
		      "Solving ODEs with Delay is not implemented. "
		      "Defaults to 0 to avoid program crash");
    result = 0.0;
    break;
  case AST_NAME_TIME:
    result = (double) data->currenttime;
    break;
	
  case AST_CONSTANT_E:
    /** exp(1) is used to adjust exponentiale to machine precision */
    result = exp(1.);
    break;
  case AST_CONSTANT_FALSE:
    result = 0.0;
    break;
  case AST_CONSTANT_PI:
    /** pi = 4 * atan 1  is used to adjust Pi to machine precision */
    result = 4.*atan(1.);
    break;
  case AST_CONSTANT_TRUE:
    result = 1.0;
    break;

  case AST_PLUS:
    result = 0.0;
    for ( i=0; i<childnum; i++) 
      result += evaluateAST(child(n,i),data);   
    break;      
  case AST_MINUS:
    if ( childnum<2 )
      result = - (evaluateAST(child(n,0),data));
    else
      result = evaluateAST(child(n,0),data) - evaluateAST(child(n,1),data);
    break;
  case AST_TIMES:
    result = 1.0;
    for ( i=0; i<childnum; i++)
    {
      result *= evaluateAST(child(n,i),data);
      if (result == 0.0) break; /* small optimization by skipping the rest of children */
    }
    break;
  case AST_DIVIDE:
    result = evaluateAST(child(n,0),data) / evaluateAST(child(n,1),data);
    break;
  case AST_POWER:
    result = pow(evaluateAST(child(n,0),data),evaluateAST(child(n,1),data));
    break;
  case AST_LAMBDA:
    SolverError_error(FATAL_ERROR_TYPE,
		      SOLVER_ERROR_AST_EVALUATION_FAILED_LAMBDA,
		      "Lambda can only be used in SBML function definitions."
		      " Defaults to 0 to avoid program crash");
    result = 0.0;
    break;
    /** FUNCTIONS: */
  case AST_FUNCTION:
    /**  Evaluate external functions, if it was set with
	 setUserDefinedFunction */      
    if ( UsrDefFunc == NULL )
    {
      SolverError_error(FATAL_ERROR_TYPE,
			SOLVER_ERROR_AST_EVALUATION_FAILED_FUNCTION,
			"The function %s() has not been defined "
			"in the SBML input model or as an externally "
			"supplied function. Defaults to 0 to avoid "
			"program crash",
			ASTNode_getName(n));
      result = 0.0;
    }
    else
    {
      double *func_vals = NULL;
      ASSIGN_NEW_MEMORY_BLOCK(func_vals, childnum+1, double, 0);
      for ( i=0; i<childnum; i++ ) 
	func_vals[i] = evaluateAST(child(n,i), data);      
      result = UsrDefFunc((char *)ASTNode_getName(n), childnum, func_vals);
      free(func_vals);
    }
    break;
  case AST_FUNCTION_ABS:
    result = fabs(evaluateAST(child(n,0),data));
    break;
  case AST_FUNCTION_ARCCOS:
    result = acos(evaluateAST(child(n,0),data)) ;
    break;
  case AST_FUNCTION_ARCCOSH:
    result = aCosh(evaluateAST(child(n,0),data));
    break;
  case AST_FUNCTION_ARCCOT:
    /** arccot x =  arctan (1 / x) */
    result = atan(1./ evaluateAST(child(n,0),data));
    break;
  case AST_FUNCTION_ARCCOTH:
    /** arccoth x = 1/2 * ln((x+1)/(x-1)) */
    value1 = evaluateAST(child(n,0),data);
    result = log((value1 + 1) / (value1 - 1))/2;
    break;
  case AST_FUNCTION_ARCCSC:
    /** arccsc(x) = Arctan(1 / sqrt((x - 1)(x + 1))) */
    value1 = evaluateAST(child(n,0),data);
    result = atan(1/sqrt((value1-1)*(value1+1)));
    break;
  case AST_FUNCTION_ARCCSCH:
    /** arccsch(x) = ln((1/x) + sqrt((1/(x*x)) + 1)) */
    value1 = evaluateAST(child(n,0),data);
    result = log(1/value1 + sqrt(1/(value1*value1)+1));
    break;
  case AST_FUNCTION_ARCSEC:
    /** arcsec(x) = arccos(1/x) */
    result = acos(1/evaluateAST(child(n,0),data));
    break;
  case AST_FUNCTION_ARCSECH:
    /* arcsech(x) = arccosh(1/x) */
    result = aCosh( 1. /  evaluateAST(child(n,0),data));
    break;
  case AST_FUNCTION_ARCSIN:
    result = asin(evaluateAST(child(n,0),data));
    break;
  case AST_FUNCTION_ARCSINH:
    result = aSinh(evaluateAST(child(n,0),data));
    break;
  case AST_FUNCTION_ARCTAN:
    result = atan(evaluateAST(child(n,0),data));
    break;
  case AST_FUNCTION_ARCTANH:
    result = aTanh(evaluateAST(child(n,0),data));
    break;
  case AST_FUNCTION_CEILING:
    result = ceil(evaluateAST(child(n,0),data));
    break;
  case AST_FUNCTION_COS:
    result = cos(evaluateAST(child(n,0),data));
    break;
  case AST_FUNCTION_COSH:
    result = cosh(evaluateAST(child(n,0),data));
    break;
  case AST_FUNCTION_COT:
    /** cot x = 1 / tan x */
    result = (1./tan(evaluateAST(child(n,0),data)));
    break;
  case AST_FUNCTION_COTH:
    /** coth x = cosh x / sinh x */
    result = 1/tanh(evaluateAST(child(n,0),data));
    break;  
  case AST_FUNCTION_CSC:
    /** csc x = 1 / sin x */
    result = (1./sin(evaluateAST(child(n,0),data)));
    break;
  case AST_FUNCTION_CSCH:
    /** csch x = 1 / sinh x  */
    result = (1./sinh(evaluateAST(child(n,0),data)));
    break;
  case AST_FUNCTION_EXP:
    result = exp(evaluateAST(child(n,0),data));
    break;
  case AST_FUNCTION_FACTORIAL:
    {
      int j;
      value1 = evaluateAST(child(n,0),data);
      j = floor(value1);
      if ( value1 != j )
	SolverError_error(FATAL_ERROR_TYPE,
			  SOLVER_ERROR_AST_EVALUATION_FAILED_FLOAT_FACTORIAL,
			  "The factorial is only implemented."
			  "for integer values. The floor value of the "
			  "passed float is used for calculation!");

      for(result=1;j>1;--j)
	result *= j;
    }
    break;
  case AST_FUNCTION_FLOOR:
    result = floor(evaluateAST(child(n,0),data));
    break;
  case AST_FUNCTION_LN:
    result = log(evaluateAST(child(n,0),data));
    break;
  case AST_FUNCTION_LOG:
    /** log(x,y) = log10(y)/log10(x) (where x is the base)  */
    result = log10(evaluateAST(child(n,1),data)) /
      log10(evaluateAST(child(n,0),data));
    break;
  case AST_FUNCTION_PIECEWISE:
    /** Piecewise: */
    found = 0;
    /** Go through n pieces with 2 AST children for each piece, */
    for ( i=0; i<(childnum-1); i=i+2 )
    {
      if ( evaluateAST(child(n, i+1), data) )
      {
	found++;
	result = evaluateAST(child(n, i), data);
      }
    }
    /** odd number of AST children: if no piece was true, otherwise remains */
    /* i should be equal to childnum for even number piecewise AST
       and equal to childnum-1 for odd numbered piecewise ASTs */
    if ( i == childnum-1 && found == 0 )
    {
      found++;
      result = evaluateAST(child(n, i), data);
    }
    if ( found == 0 )
      SolverError_error(FATAL_ERROR_TYPE,
			SOLVER_ERROR_AST_EVALUATION_FAILED_PIECEWISE,
			"Piecewise function failed; no true piece");
    if ( found > 1 )
      SolverError_error(FATAL_ERROR_TYPE,
			SOLVER_ERROR_AST_EVALUATION_FAILED_PIECEWISE,
			"Piecewise function failed; several true pieces");
    break;
  case AST_FUNCTION_POWER:
    result = pow(evaluateAST(child(n,0),data),evaluateAST(child(n,1),data));
    break;
  case AST_FUNCTION_ROOT:
    /*!!! ALSO do this in compiled code */
    value1 = evaluateAST(child(n,1),data);
    value2 = evaluateAST(child(n,0),data);
    value3 = floor(value2);
    /* for odd root degrees, negative numbers are OK */
    if ( value2 == value3 ) /* check whether degree is integer */
    {
      if ( (value1 < 0) && ((int)value2 % 2 != 0) )
	result = - pow(fabs(value1), 1./value2);
      else
	result = pow(value1, 1./value2);	
    }
    else
      result = pow(value1, 1./value2);
    break;
  case AST_FUNCTION_SEC:
    /** sec x = 1 / cos x */
    result = 1./cos(evaluateAST(child(n,0),data));
    break;
  case AST_FUNCTION_SECH:
    /** sech x = 1 / cosh x */
    result = 1./cosh(evaluateAST(child(n,0),data));
    break;
  case AST_FUNCTION_SIN:
    result = sin(evaluateAST(child(n,0),data));
    break;
  case AST_FUNCTION_SINH:
    result = sinh(evaluateAST(child(n,0),data));
    break;
  case AST_FUNCTION_TAN:
    result = tan(evaluateAST(child(n,0),data));
    break;
  case AST_FUNCTION_TANH:
    result = tanh(evaluateAST(child(n,0),data));
    break;

  case AST_LOGICAL_AND:
    /** AND: all children must be true */
    true = 0;
    for ( i=0; i<childnum; i++ ) true += evaluateAST(child(n,i),data);
    if ( true == childnum ) result = 1.0;
    else result = 0.0;
    break;
  case AST_LOGICAL_NOT:
    result = (double) (!(evaluateAST(child(n,0),data)));
    break;
  case AST_LOGICAL_OR:
    /** OR: at least one child must be true */
    true = 0;
    for ( i=0; i<childnum; i++ ) true += evaluateAST(child(n,i),data);    
    if ( true > 0 ) result = 1.0;
    else result = 0.0;
    break;
  case AST_LOGICAL_XOR:
    /* n-ary: true if an odd number of children is true */
    true = 0;
    for ( i=0; i<childnum; i++ ) true += evaluateAST(child(n,i),data);    
    if ( true % 2 != 0 ) result = 1.0;
    else result = 0.0;
    break;
    /* !!! check n-ary definitions for relational operators !!! */
  case AST_RELATIONAL_EQ:
    /** n-ary: all children must be equal */
    result = 1.0;
    if (childnum <= 1) break;
    value1 = evaluateAST(child(n,0),data);
    for ( i=1; i<childnum; i++ )
      if ( value1 != evaluateAST(child(n,i),data) )
      {
        result = 0.0;
        break;
      }
    break;
  case AST_RELATIONAL_GEQ:
    /** n-ary: each child must be greater than or equal to the following */
    result = 1.0;
    if (childnum <= 1) break;
    value2 = evaluateAST(child(n,0),data);
    for ( i=1; i<childnum; i++ )
    {
      value1 = value2;
      value2 = evaluateAST(child(n,i),data);
      if (value1 < value2)
      {
        result = 0.0;
        break;
      }
    }
    break;
  case AST_RELATIONAL_GT:
    /** n-ary: each child must be greater than the following */
    result = 1.0;
    if (childnum <= 1) break;
    value2 = evaluateAST(child(n,0),data);
    for ( i=1; i<childnum; i++ )
    {
      value1 = value2;
      value2 = evaluateAST(child(n,i),data);
      if (value1 <= value2)
      {
        result = 0.0;
        break;
      }
    }
    break;
  case AST_RELATIONAL_LEQ:
    /** n-ary: each child must be lower than or equal to the following */
    result = 1.0;
    if (childnum <= 1) break;
    value2 = evaluateAST(child(n,0),data);
    for ( i=1; i<childnum; i++ )
    {
      value1 = value2;
      value2 = evaluateAST(child(n,i),data);
      if (value1 > value2)
      {
        result = 0.0;
        break;
      }
    }
    break;
  case AST_RELATIONAL_LT :
    /* n-ary: each child must be lower than the following */
    result = 1.0;
    if (childnum <= 1) break;
    value2 = evaluateAST(child(n,0),data);
    for ( i=1; i<childnum; i++ )
    {
      value1 = value2;
      value2 = evaluateAST(child(n,i),data);
      if (value1 >= value2)
      {
        result = 0.0;
        break;
      }
    }
    break;
  case AST_RELATIONAL_NEQ:
    /* binary "not equal" */
    result = (evaluateAST(child(n, 0), data) != evaluateAST(child(n, 1), data));
    break;
  default:
    SolverError_error(FATAL_ERROR_TYPE,
                      SOLVER_ERROR_AST_UNKNOWN_NODE_TYPE,
                      "evaluateAST: unknown ASTNode type: %d", type);
    result = 0;
    break;
  }
  
  return result;
}
SBML_ODESOLVER_API int drawJacoby(cvodeData_t *data, char *file, char *format) {

#if !USE_GRAPHVIZ

  SolverError_error(
		    WARNING_ERROR_TYPE,
		    SOLVER_ERROR_NO_GRAPHVIZ,
		    "odeSolver has been compiled without GRAPHIZ functionality. ",
		    "Graphs are printed to stdout in the graphviz' .dot format.");

  drawJacobyTxt(data, file);

#else

  int i, j;
  GVC_t *gvc;
  Agraph_t *g;
  Agnode_t *r;
  Agnode_t *s;  
  Agedge_t *e;
  Agsym_t *a;
  char name[WORDSIZE];
  char label[WORDSIZE];  
  char *output[3];
  char *command = "dot";
  char *formatopt;
  char *outfile;
  

  /* setting name of outfile */
  ASSIGN_NEW_MEMORY_BLOCK(outfile, strlen(file)+ strlen(format)+7, char, 0);
  sprintf(outfile, "-o%s_jm.%s", file, format);
  
  /* setting output format */
  ASSIGN_NEW_MEMORY_BLOCK(formatopt, strlen(format)+3, char, 0);
  sprintf(formatopt, "-T%s", format); 

  /* construct command-line */
  output[0] = command;
  output[1] = formatopt;
  output[2] = outfile;
  output[3] = NULL;
    
  /* set up renderer context */
  gvc = (GVC_t *) gvContext();
#if GRAPHVIZ_MAJOR_VERSION == 2 && GRAPHVIZ_MINOR_VERSION < 4
  dotneato_initialize(gvc, 3, output);
#elif GRAPHVIZ_MAJOR_VERSION == 2 && GRAPHVIZ_MINOR_VERSION == 4
  parse_args(gvc, 3, output);
#elif GRAPHVIZ_MAJOR_VERSION == 2 && GRAPHVIZ_MINOR_VERSION >= 6 || GRAPHVIZ_MAJOR_VERSION >= 3
  gvParseArgs(gvc, 3, output);
#endif

  g = agopen("G", AGDIGRAPH);

  /* avoid overlapping nodes, for graph embedding by neato */
  a = agraphattr(g, "overlap", "");
  agxset(g, a->index, "scale");

  /* set graph label */
  if ( Model_isSetName(data->model->m) )
    sprintf(label, "%s at time %g",  Model_getName(data->model->m),
	    data->currenttime);
  else if ( Model_isSetId(data->model->m) )
    sprintf(label, "%s at time %g",  Model_getId(data->model->m),
	    data->currenttime);
  else
    sprintf(label, "label=\"at time %g\";\n", data->currenttime);

  a = agraphattr(g, "label", "");
  agxset(g, a->index, label);
  
  /*
    Set edges from species A to species B if the
    corresponding entry in the jacobian ((d[B]/dt)/d[A])
    is not '0'. Set edge color 'red' and arrowhead 'tee'
    if negative.
  */

  for ( i=0; i<data->model->neq; i++ ) {
    for ( j=0; j<data->model->neq; j++ ) {
      if ( evaluateAST(data->model->jacob[i][j], data) != 0 ) {
	
	sprintf(name, "%s", data->model->names[j]);
	r = agnode(g,name);
	agset(r, "label", data->model->names[j]);

	sprintf(label, "%s.htm", data->model->names[j]);
	a = agnodeattr(g, "URL", "");
	agxset(r, a->index, label);
	
	sprintf(name,"%s", data->model->names[i]);
	s = agnode(g,name);
	agset(s, "label", data->model->names[i]);

	sprintf(label, "%s.htm", data->model->names[i]);	
	a = agnodeattr(g, "URL", "");
	agxset(s, a->index, label);
	
	e = agedge(g,r,s);

	a = agedgeattr(g, "label", "");
	sprintf(name, "%g",  evaluateAST(data->model->jacob[i][j], data)); 
	agxset (e, a->index, name);


	
	if ( evaluateAST(data->model->jacob[i][j], data) < 0 ) {
	  a = agedgeattr(g, "arrowhead", "");
	  agxset(e, a->index, "tee");
	  a = agedgeattr(g, "color", "");
	  agxset(e, a->index, "red"); 	    
	}	
      }
    }
  }
  
  /* Compute a layout */
#if GRAPHVIZ_MAJOR_VERSION == 2 && GRAPHVIZ_MINOR_VERSION <= 2
  gvBindContext(gvc, g);
  dot_layout(g);
#elif GRAPHVIZ_MAJOR_VERSION == 2 && GRAPHVIZ_MINOR_VERSION == 4
  gvlayout_layout(gvc, g);
#elif GRAPHVIZ_MAJOR_VERSION == 2 && GRAPHVIZ_MINOR_VERSION >= 6 || GRAPHVIZ_MAJOR_VERSION >= 3
  gvLayoutJobs(gvc, g);
#endif
  
  /* Write the graph according to -T and -o options */
#if GRAPHVIZ_MAJOR_VERSION == 2 && GRAPHVIZ_MINOR_VERSION <= 2
  dotneato_write(gvc);
#elif GRAPHVIZ_MAJOR_VERSION == 2 && GRAPHVIZ_MINOR_VERSION == 4
  emit_jobs(gvc, g);
#elif GRAPHVIZ_MAJOR_VERSION == 2 && GRAPHVIZ_MINOR_VERSION >= 6 || GRAPHVIZ_MAJOR_VERSION >= 3
  gvRenderJobs(gvc, g);
#endif
  
  /* Clean out layout data */
#if GRAPHVIZ_MAJOR_VERSION == 2 && GRAPHVIZ_MINOR_VERSION <= 2
  dot_cleanup(g);
#elif GRAPHVIZ_MAJOR_VERSION == 2 && GRAPHVIZ_MINOR_VERSION == 4
  gvlayout_cleanup(gvc, g);
#elif GRAPHVIZ_MAJOR_VERSION == 2 && GRAPHVIZ_MINOR_VERSION >= 6 || GRAPHVIZ_MAJOR_VERSION >= 3
  gvFreeLayout(gvc, g);
#endif
  
  /* Free graph structures */
#if GRAPHVIZ_MAJOR_VERSION == 2 && GRAPHVIZ_MINOR_VERSION <= 2
  dot_cleanup(g);
#endif
  agclose(g);

  /* Clean up output file and errors */
#if GRAPHVIZ_MAJOR_VERSION == 2 && GRAPHVIZ_MINOR_VERSION <= 2
  gvFREEcontext(gvc);
  dotneato_eof(gvc);
#elif GRAPHVIZ_MAJOR_VERSION == 2 && GRAPHVIZ_MINOR_VERSION == 4
  dotneato_terminate(gvc);
#elif (GRAPHVIZ_MAJOR_VERSION == 2 && GRAPHVIZ_MINOR_VERSION >= 6) || GRAPHVIZ_MAJOR_VERSION >= 3
  gvFreeContext(gvc);
#endif

  xfree(formatopt);
  xfree(outfile);

#endif

  return 1;
}
SBML_ODESOLVER_API int drawModel(Model_t *m, char* file, char *format) {
  
#if !USE_GRAPHVIZ

  SolverError_error(
		    WARNING_ERROR_TYPE,
		    SOLVER_ERROR_NO_GRAPHVIZ,
		    "odeSolver has been compiled without GRAPHIZ functionality. ",
		    "Graphs are printed to stdout in the graphviz' .dot format.");
  drawModelTxt(m, file);
  
#else

  GVC_t *gvc;
  Agraph_t *g;
  Agnode_t *r;
  Agnode_t *s;  
  Agedge_t *e;
  Agsym_t *a;
  Species_t *sp;
  Reaction_t *re;
  const ASTNode_t *math;  
  SpeciesReference_t *sref;
  ModifierSpeciesReference_t *mref;
  char *output[4];
  char *command = "dot";
  char *formatopt;
  char *outfile;
  int i,j;
  int reversible;
  char name[WORDSIZE];
  char label[WORDSIZE];

  /* setting name of outfile */
  ASSIGN_NEW_MEMORY_BLOCK(outfile, strlen(file)+ strlen(format)+7, char, 0);
  sprintf(outfile, "-o%s_rn.%s", file, format);

  /* setting output format */
  ASSIGN_NEW_MEMORY_BLOCK(formatopt, strlen(format)+3, char, 0);
  sprintf(formatopt, "-T%s", format);

  /* construct command-line */
  output[0] = command;
  output[1] = formatopt;
  output[2] = outfile;
  output[3] = NULL;

  /* set up renderer context */
  gvc = (GVC_t *) gvContext();
#if GRAPHVIZ_MAJOR_VERSION == 2 && GRAPHVIZ_MINOR_VERSION < 4
  dotneato_initialize(gvc, 3, output);
#elif GRAPHVIZ_MAJOR_VERSION == 2 && GRAPHVIZ_MINOR_VERSION == 4
  parse_args(gvc, 3, output);
#elif GRAPHVIZ_MAJOR_VERSION == 2 && GRAPHVIZ_MINOR_VERSION >= 6 || GRAPHVIZ_MAJOR_VERSION >= 3
  gvParseArgs(gvc, 3, output);  
#endif  

  g = agopen("G", AGDIGRAPH);
  
  /* avoid overlapping nodes, for graph embedding by neato */ 
  a = agraphattr(g, "overlap", "");
  agxset(g, a->index, "scale");

  for ( i=0; i<Model_getNumReactions(m); i++ ) {

    re = Model_getReaction(m,i);
    reversible = Reaction_getReversible(re);
    sprintf(name, "%s", Reaction_getId(re));
    r = agnode(g,name);
    a = agnodeattr(g, "shape", "ellipse");    
    agxset(r, a->index, "box");
    
    sprintf(label, "%s", Reaction_isSetName(re) ?
	    Reaction_getName(re) : Reaction_getId(re));
    agset(r, "label", label);
    
    sprintf(label, "%s.htm", Reaction_getId(re));
    a = agnodeattr(g, "URL", "");
    agxset(r, a->index, label);
    
    for ( j=0; j<Reaction_getNumModifiers(re); j++ ) {

      mref = Reaction_getModifier(re,j);
      sp = Model_getSpeciesById(m, ModifierSpeciesReference_getSpecies(mref));
      
      sprintf(name,"%s", Species_getId(sp));
      s = agnode(g,name);
      sprintf(label, "%s", Species_isSetName(sp) ? 
	      Species_getName(sp) : Species_getId(sp));
      agset(s, "label", label);

      if ( Species_getBoundaryCondition(sp) ) {
	a = agnodeattr(g, "color", "");
	agxset(s, a->index, "blue");
      }
      if ( Species_getConstant(sp) ) {
	a = agnodeattr(g, "color", "");
	agxset(s, a->index, "green4");
      }

      sprintf(label, "%s.htm", Species_getId(sp));
      a = agnodeattr(g, "URL", "");
      agxset(s, a->index, label);
	
      e = agedge(g,s,r);
      a = agedgeattr(g, "style", "");
      agxset(e, a->index, "dashed");
      a = agedgeattr(g, "arrowhead", "");
      agxset(e, a->index, "odot");
    }

    for ( j=0; j<Reaction_getNumReactants(re); j++ ) {

      sref = Reaction_getReactant(re,j);
      sp = Model_getSpeciesById(m, SpeciesReference_getSpecies(sref));
      
      sprintf(name,"%s", Species_getId(sp));
      s = agnode(g, name);
      sprintf(label, "%s", Species_isSetName(sp) ? 
	      Species_getName(sp) : Species_getId(sp));
      agset(s, "label", label);

      if ( Species_getBoundaryCondition(sp) ) {
	a = agnodeattr(g, "color", "");
	agxset(s, a->index, "blue");
      }
      if ( Species_getConstant(sp) ) {
	a = agnodeattr(g, "color", "");
	agxset(s, a->index, "green4");
      }

      sprintf(label, "%s.htm", Species_getId(sp));
      a = agnodeattr(g, "URL", "");
      agxset(s, a->index, label);
      
      e = agedge(g,s,r);
      a = agedgeattr(g, "label", "");
      
      if ( (SpeciesReference_isSetStoichiometryMath(sref)) ) {
	math = SpeciesReference_getStoichiometryMath(sref);
	if ( (strcmp(SBML_formulaToString(math),"1") !=
	      0) ) {
	  agxset (e, a->index, SBML_formulaToString(math));
	}
      }
      else {
	if ( SpeciesReference_getStoichiometry(sref) != 1 ) {
	  sprintf(name, "%g", SpeciesReference_getStoichiometry(sref));
	  agxset (e, a->index, name);
	}
      }
      if ( reversible == 1 ) {
	a = agedgeattr(g, "arrowtail", "");
	agxset(e, a->index, "onormal");
      }      
    }
    
    for ( j=0; j<Reaction_getNumProducts(re); j++ ) {
      sref = Reaction_getProduct(re,j);
      sp = Model_getSpeciesById(m, SpeciesReference_getSpecies(sref));
      sprintf(name,"%s", Species_getId(sp));
      s = agnode(g,name);
      sprintf(label, "%s", Species_isSetName(sp) ? 
	      Species_getName(sp) : Species_getId(sp));
      agset(s, "label", label);

      if ( Species_getBoundaryCondition(sp) ) {
	a = agnodeattr(g, "color", "");
	agxset(s, a->index, "blue");
      }
      if ( Species_getConstant(sp) ) {
	a = agnodeattr(g, "color", "");
	agxset(s, a->index, "green4");
      }

      sprintf(label, "%s.htm", Species_getId(sp));
      a = agnodeattr(g, "URL", "");
      agxset(s, a->index, label);
            
      e = agedge(g,r,s);
      a = agedgeattr(g, "label", "");
      if ( SpeciesReference_isSetStoichiometryMath(sref) ) {
	math = SpeciesReference_getStoichiometryMath(sref);
	if ( (strcmp(SBML_formulaToString(math),"1") !=
	      0) ) {
	  agxset (e, a->index, SBML_formulaToString(math));
	}
      }
      else {
	if ( SpeciesReference_getStoichiometry(sref) != 1 ) {
	  sprintf(name, "%g",SpeciesReference_getStoichiometry(sref));
	  agxset (e, a->index,name);
	}
      }
      if ( reversible == 1 ) {
	a = agedgeattr(g, "arrowtail", "");
	agxset(e, a->index, "onormal");
      }      
    }   
  }

  /* Compute a layout */
#if GRAPHVIZ_MAJOR_VERSION == 2 && GRAPHVIZ_MINOR_VERSION <= 2
  gvBindContext(gvc, g);
  dot_layout(g);
#elif GRAPHVIZ_MAJOR_VERSION == 2 && GRAPHVIZ_MINOR_VERSION == 4
  gvlayout_layout(gvc, g);
#elif GRAPHVIZ_MAJOR_VERSION == 2 && GRAPHVIZ_MINOR_VERSION >= 6 || GRAPHVIZ_MAJOR_VERSION >= 3
  gvLayoutJobs(gvc, g);
#endif

  /* Write the graph according to -T and -o options */
#if GRAPHVIZ_MAJOR_VERSION == 2 && GRAPHVIZ_MINOR_VERSION <= 2
  dotneato_write(gvc);
#elif GRAPHVIZ_MAJOR_VERSION == 2 && GRAPHVIZ_MINOR_VERSION == 4
  emit_jobs(gvc, g);
#elif GRAPHVIZ_MAJOR_VERSION == 2 && GRAPHVIZ_MINOR_VERSION >= 6 || GRAPHVIZ_MAJOR_VERSION >= 3
  gvRenderJobs(gvc, g);
#endif
  
  /* Clean out layout data */
#if GRAPHVIZ_MAJOR_VERSION == 2 && GRAPHVIZ_MINOR_VERSION <= 2
  dot_cleanup(g);
#elif GRAPHVIZ_MAJOR_VERSION == 2 && GRAPHVIZ_MINOR_VERSION == 4
  gvlayout_cleanup(gvc, g);
#elif GRAPHVIZ_MAJOR_VERSION == 2 && GRAPHVIZ_MINOR_VERSION >= 6 || GRAPHVIZ_MAJOR_VERSION >= 3
  gvFreeLayout(gvc, g);
#endif
  
  /* Free graph structures */
#if GRAPHVIZ_MAJOR_VERSION == 2 && GRAPHVIZ_MINOR_VERSION <= 2
  dot_cleanup(g);
#else
  agclose(g);
#endif

  /* Clean up output file and errors */
#if GRAPHVIZ_MAJOR_VERSION == 2 && GRAPHVIZ_MINOR_VERSION <= 2
  gvFREEcontext(gvc);
  dotneato_eof(gvc);
#elif GRAPHVIZ_MAJOR_VERSION == 2 && GRAPHVIZ_MINOR_VERSION == 4
  dotneato_terminate(gvc);
#elif GRAPHVIZ_MAJOR_VERSION == 2 && GRAPHVIZ_MINOR_VERSION >= 6 || GRAPHVIZ_MAJOR_VERSION >= 3
  gvFreeContext(gvc); 
#endif  

  xfree(formatopt);  
  xfree(outfile);
  
#endif
  
  return 1;

}
Beispiel #12
0
SBML_ODESOLVER_API int IntegratorInstance_idaOneStep(integratorInstance_t *engine)
{
    int i, flag;
    realtype *ydata = NULL;
    
    cvodeSolver_t *solver = engine->solver;
    cvodeData_t *data = engine->data;
/*     cvodeSettings_t *opt = engine->opt; */
/*     cvodeResults_t *results = engine->results; */
    odeModel_t *om = engine->om;
    
    /* !!!! calling CVODE !!!! */
    flag = -1; /* IDASolver(solver->cvode_mem, solver->tout, &(solver->t),
		  solver->y, solver->dy, IDA_NORMAL); */

    if ( flag != IDA_SUCCESS )
      {
        static const char *message[] =
	  {
	    /*  0 IDA_SUCCESS */
	    "Success",
	    /**/
	    /*  1 IDA_ROOT_RETURN */
	    /*   "CVode succeeded, and found one or more roots" */
	    /*  2 IDA_TSTOP_RETURN */
	    /*   "CVode succeeded and returned at tstop" */
	    /**/
	    /* -1 IDA_MEM_NULL -1 (old CVODE_NO_MEM) */
	    "The cvode_mem argument was NULL",
	    /* -2 IDA_ILL_INPUT */
	    "One of the inputs to CVode is illegal. This "
	    "includes the situation when a component of the "
	    "error weight vectors becomes < 0 during "
	    "internal time-stepping. The ILL_INPUT flag "
	    "will also be returned if the linear solver "
	    "routine CV--- (called by the user after "
	    "calling CVodeMalloc) failed to set one of the "
	    "linear solver-related fields in cvode_mem or "
	    "if the linear solver's init routine failed. In "
	    "any case, the user should see the printed "
	    "error message for more details.",
	    /* -3 IDA_NO_MALLOC */
	    "cvode_mem was not allocated",
	    /* -4 IDA_TOO_MUCH_WORK */
	    "The solver took %g internal steps but could not "
	    "compute variable values for time %g",
	    /* -5 IDA_TOO_MUCH_ACC */
	    "The solver could not satisfy the accuracy " 
	    "requested for some internal step.",
	    /* -6 IDA_ERR_FAILURE */
	    "Error test failures occurred too many times "
	    "during one internal time step or "
	    "occurred with |h| = hmin.",
	    /* -7 IDA_CONV_FAILURE */
	    "Convergence test failures occurred too many "
	    "times during one internal time step or occurred "
	    "with |h| = hmin.",
	    /* -8 IDA_LINIT_FAIL */
	    "CVode -- Initial Setup: "
	    "The linear solver's init routine failed.",
	    /* -9 IDA_LSETUP_FAIL */
	    "The linear solver's setup routine failed in an "
	    "unrecoverable manner.",
	    /* -10 IDA_LSOLVE_FAIL */
	    "The linear solver's solve routine failed in an "
	    "unrecoverable manner.",
	    /* -11 IDA_MEM_FAIL */
	    "A memory allocation failed. "
	    "(including an attempt to increase maxord)",
	    /* -12 IDA_RTFUNC_NULL */
	    "nrtfn > 0 but g = NULL.",
	    /* -13 IDA_NO_SLDET */
	    "CVodeGetNumStabLimOrderReds -- Illegal attempt "
	    "to call without enabling SLDET.",
	    /* -14 IDA_BAD_K */
	    "CVodeGetDky -- Illegal value for k.",
	    /* -15 IDA_BAD_T */
	    "CVodeGetDky -- Illegal value for t.",
	    /* -16 IDA_BAD_DKY */
	    "CVodeGetDky -- dky = NULL illegal.",
	    /* -17 IDA_PDATA_NULL */
	    "???",
	  };
	    
	SolverError_error(
			  ERROR_ERROR_TYPE,
			  flag,
			  (abs(flag) < (int)NUMBER_OF_ELEMENTS(message)) ? message[abs(flag)] : "???",
			  solver->tout);
	SolverError_error(
			  WARNING_ERROR_TYPE,
			  SOLVER_ERROR_INTEGRATION_NOT_SUCCESSFUL,
			  "Integration not successful. Results are not complete.");

	return 0 ; /* Error - stop integration*/
      }
    
    ydata = NV_DATA_S(solver->y);

    
    /* update cvodeData time dependent variables */    
    for ( i=0; i<om->neq; i++ )
      data->value[i] = ydata[i];

    /* update rest of data with internal default function */
    return IntegratorInstance_updateData(engine);

}
Beispiel #13
0
/* creates CVODES structures and fills cvodeSolver 
   return 1 => success
   return 0 => failure
*/
int
IntegratorInstance_createIDASolverStructures(integratorInstance_t *engine)
{
  int i, flag, neq, nalg;
  realtype *ydata, *abstoldata, *dydata;
  
  odeModel_t *om = engine->om;
  cvodeData_t *data = engine->data;
  cvodeSolver_t *solver = engine->solver;
  cvodeSettings_t *opt = engine->opt;
  
  neq = engine->om->neq;   /* number of ODEs */
  nalg = engine->om->nalg; /* number of algebraic constraints */
  
  /* construct jacobian, if wanted and not yet existing */
  if ( opt->UseJacobian && om->jacob == NULL ) 
    /* reset UseJacobian option, depending on success */
    engine->UseJacobian = ODEModel_constructJacobian(om);
  else if ( !opt->UseJacobian )
  {
    /* free jacobian from former runs (not necessary, frees also
       unsuccessful jacobians from former runs ) */
    ODEModel_freeJacobian(om);
    SolverError_error(WARNING_ERROR_TYPE,
		      SOLVER_ERROR_MODEL_NOT_SIMPLIFIED,
		      "Jacobian matrix construction skipped.");
    engine->UseJacobian = om->jacobian;
  }
  /* construct algebraic `Jacobian' (or do that in constructJacobian */
  
  /* CVODESolverStructures from former runs must be freed */
  if ( engine->run > 1 )
    IntegratorInstance_freeIDASolverStructures(engine);
  
  
    /*
     * Allocate y, abstol vectors
     */
  solver->y = N_VNew_Serial(neq + nalg);
  CVODE_HANDLE_ERROR((void *)solver->y, "N_VNew_Serial for vector y", 0);
  
  solver->dy = N_VNew_Serial(neq + nalg);
  CVODE_HANDLE_ERROR((void *)solver->dy, "N_VNew_Serial for vector dy", 0);
    
  solver->abstol = N_VNew_Serial(neq + nalg);
  CVODE_HANDLE_ERROR((void *)solver->abstol,
		     "N_VNew_Serial for vector abstol", 0);
  
  /*
   * Initialize y, abstol vectors
   */
  ydata      = NV_DATA_S(solver->y);
  abstoldata = NV_DATA_S(solver->abstol);
  dydata     = NV_DATA_S(solver->dy);
  
  for ( i=0; i<neq; i++ )
  {
    /* Set initial value vector components of y and y' */
    ydata[i] = data->value[i];
    /* Set absolute tolerance vector components,
       currently the same absolute error is used for all y */ 
    abstoldata[i] = opt->Error;
    dydata[i] = evaluateAST(om->ode[i], data);
  }
  /* set initial value vector components for algebraic rule variables  */
    
  /* scalar relative tolerance: the same for all y */
  solver->reltol = opt->RError;

  /*
   * Call IDACreate to create the solver memory:
   *
   */
  solver->cvode_mem = IDACreate();
  CVODE_HANDLE_ERROR((void *)(solver->cvode_mem), "IDACreate", 0);

  /*
   * Call IDAInit to initialize the integrator memory:
   *
   * cvode_mem  pointer to the CVode memory block returned by CVodeCreate
   * fRes         user's right hand side function
   * t0         initial value of time
   * y          the dependent variable vector
   * dy         the ODE value vector
   */
  flag = IDAInit(solver->cvode_mem, fRes, solver->t0, solver->y,
                 solver->dy);
  CVODE_HANDLE_ERROR(&flag, "IDAInit", 1);
  /*
   * specify scalar relative and vector absolute tolerances
   * reltol     the scalar relative tolerance
   * abstol     pointer to the absolute tolerance vector
   */
  flag = IDASVtolerances(solver->cvode_mem, solver->reltol, solver->abstol);
  CVODE_HANDLE_ERROR(&flag, "IDASVtolerances", 1);

  /* 
   * Link the main integrator with data for right-hand side function
   */ 
  flag = IDASetUserData(solver->cvode_mem, engine->data);
  CVODE_HANDLE_ERROR(&flag, "IDASetUserData", 1);
    
  /*
   * Link the main integrator with the IDADENSE linear solver
   */
  flag = IDADense(solver->cvode_mem, neq);
  CVODE_HANDLE_ERROR(&flag, "IDADense", 1);


  /*
   * Set the routine used by the IDADense linear solver
   * to approximate the Jacobian matrix to ...
   */
  if ( opt->UseJacobian == 1 ) {
    /* ... user-supplied routine JacRes : put JacRes instead of NULL
       when working */
    flag = IDADlsSetDenseJacFn(solver->cvode_mem, JacRes);
    CVODE_HANDLE_ERROR(&flag, "IDADlsSetDenseJacFn", 1);
  }
     
  return 1; /* OK */
}
Beispiel #14
0
void storeSBMLError(errorType_t type, const ParseMessage_t *pm)
{
    SolverError_error(type, ParseMessage_getId(pm),
		      (char*) ParseMessage_getMessage(pm)); 
}