Main Page | Modules | Data Structures | Directories | File List | Data Fields | Globals

integratorInstance.h

Go to the documentation of this file.
00001 /*
00002   Last changed Time-stamp: <2005-12-15 20:37:17 raim>
00003   $Id: integratorInstance.h,v 1.21 2005/12/15 19:54:06 raimc Exp $ 
00004 */
00005 /* 
00006  *
00007  * This library is free software; you can redistribute it and/or modify it
00008  * under the terms of the GNU Lesser General Public License as published
00009  * by the Free Software Foundation; either version 2.1 of the License, or
00010  * any later version.
00011  *
00012  * This library is distributed in the hope that it will be useful, but
00013  * WITHOUT ANY WARRANTY, WITHOUT EVEN THE IMPLIED WARRANTY OF
00014  * MERCHANTABILITY OR FITNESS FOR A PARTICULAR PURPOSE. The software and
00015  * documentation provided hereunder is on an "as is" basis, and the
00016  * authors have no obligations to provide maintenance, support,
00017  * updates, enhancements or modifications.  In no event shall the
00018  * authors be liable to any party for direct, indirect, special,
00019  * incidental or consequential damages, including lost profits, arising
00020  * out of the use of this software and its documentation, even if the
00021  * authors have been advised of the possibility of such damage.  See
00022  * the GNU Lesser General Public License for more details.
00023  *
00024  * You should have received a copy of the GNU Lesser General Public License
00025  * along with this library; if not, write to the Free Software Foundation,
00026  * Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA.
00027  *
00028  * The original code contained here was initially developed by:
00029  *
00030  *     Andrew Finney
00031  *
00032  * Contributor(s):
00033  *     Rainer Machne
00034  */
00035 
00036 #ifndef _INTEGRATORINSTANCE_H_
00037 #define _INTEGRATORINSTANCE_H_
00038 
00039 /* Header Files for CVODE */
00040 #include "nvector_serial.h"
00041 
00042 #include "sbmlsolver/exportdefs.h"
00043 #include "sbmlsolver/integratorSettings.h"
00044 #include "sbmlsolver/cvodedatatype.h"
00045 #include "sbmlsolver/odeModel.h"
00046 #include "sbmlsolver/cvodedata.h"
00047 
00048 #ifdef __cplusplus
00049 extern "C" {
00050 #endif
00051 
00052   typedef struct cvodeSolver cvodeSolver_t;
00053   typedef struct integratorInstance integratorInstance_t ;
00054 
00056   struct cvodeSolver
00057   {
00058     
00059     double t, tout, t0;
00060     int iout, nout;  
00061     realtype reltol, atol1;
00062     N_Vector y, abstol, senstol;
00063     void *cvode_mem; 
00064     int nsens;
00065     N_Vector *yS;    
00066     N_Vector dy;     
00067   };
00068 
00069 
00071   struct integratorInstance
00072   {
00075     odeModel_t *om;
00078     cvodeSettings_t *opt;
00081     cvodeData_t *data;
00083     cvodeSolver_t *solver;
00085     cvodeResults_t *results; 
00086   };
00087   
00088   /* common to all solvers */
00089   SBML_ODESOLVER_API integratorInstance_t *IntegratorInstance_create(odeModel_t *, cvodeSettings_t *);
00090   SBML_ODESOLVER_API int IntegratorInstance_set(integratorInstance_t *, cvodeSettings_t *);
00091   SBML_ODESOLVER_API int IntegratorInstance_reset(integratorInstance_t *);
00092   SBML_ODESOLVER_API cvodeSettings_t *IntegratorInstance_getSettings(integratorInstance_t *);
00093   SBML_ODESOLVER_API void IntegratorInstance_copyVariableState(integratorInstance_t *target, integratorInstance_t *source);
00094   SBML_ODESOLVER_API double IntegratorInstance_getTime(integratorInstance_t *);
00095   SBML_ODESOLVER_API double IntegratorInstance_getVariableValue(integratorInstance_t *, variableIndex_t *);
00096   SBML_ODESOLVER_API double IntegratorInstance_getSensitivity(integratorInstance_t *, variableIndex_t *y,  variableIndex_t *p);
00097   SBML_ODESOLVER_API int IntegratorInstance_setNextTimeStep(integratorInstance_t *, double);
00098   SBML_ODESOLVER_API void IntegratorInstance_dumpNames(integratorInstance_t *);
00099   SBML_ODESOLVER_API void IntegratorInstance_dumpData(integratorInstance_t *);
00100   SBML_ODESOLVER_API void IntegratorInstance_dumpYSensitivities(integratorInstance_t *, variableIndex_t *);
00101   SBML_ODESOLVER_API void IntegratorInstance_dumpPSensitivities(integratorInstance_t *, variableIndex_t *);
00102   SBML_ODESOLVER_API cvodeData_t *IntegratorInstance_getData(integratorInstance_t *);
00103   SBML_ODESOLVER_API int IntegratorInstance_integrate(integratorInstance_t *);
00104   SBML_ODESOLVER_API int IntegratorInstance_checkTrigger(integratorInstance_t *);
00105   SBML_ODESOLVER_API int IntegratorInstance_checkSteadyState(integratorInstance_t *);
00106   SBML_ODESOLVER_API int IntegratorInstance_timeCourseCompleted(integratorInstance_t *);
00107   SBML_ODESOLVER_API cvodeResults_t *IntegratorInstance_createResults(integratorInstance_t *);
00108   SBML_ODESOLVER_API int IntegratorInstance_updateModel(integratorInstance_t*);
00109   SBML_ODESOLVER_API int IntegratorInstance_simpleOneStep(integratorInstance_t *);
00110   
00111   /* these functions contain solver specific switches and need to be adapted
00112      for any new solver, and so does the local
00113      integratorInstance_initialiyeSolverStructures */    
00114   SBML_ODESOLVER_API void IntegratorInstance_setVariableValue(integratorInstance_t *, variableIndex_t *, double);
00115   SBML_ODESOLVER_API int IntegratorInstance_integrateOneStep(integratorInstance_t *);
00116   SBML_ODESOLVER_API void IntegratorInstance_dumpSolver(integratorInstance_t *);
00117   SBML_ODESOLVER_API void IntegratorInstance_free(integratorInstance_t *);
00118   SBML_ODESOLVER_API int IntegratorInstance_handleError(integratorInstance_t *);
00119   SBML_ODESOLVER_API void IntegratorInstance_printStatistics(integratorInstance_t *, FILE *f);
00120   
00121   
00122 #ifdef __cplusplus
00123 }
00124 #endif
00125 
00126 /* default function for data update, event and steady state handling,
00127    result storage and loop variables; to be used by solver
00128    specific ...OneStep functions */
00129 int IntegratorInstance_updateData(integratorInstance_t *);
00130 
00131 #endif

Generated on Wed Dec 21 18:10:07 2005 for SBML ODE Solver Library API by  doxygen 1.4.4