copasi API  0.1
Classes | Functions
copasi/copasi_api.h File Reference
#include "TC_structs.h"

Go to the source code of this file.

Classes

struct  copasi_model
 this struct is used to contain a pointer to an instance of a COPASI class More...
struct  copasi_reaction
 this struct is used to contain a pointer to an instance of a COPASI class More...
struct  copasi_compartment
 this struct is used to contain a pointer to an instance of a COPASI class More...

Functions

BEGIN_C_DECLS TCAPIEXPORT void copasi_end ()
 initialize copasi -- MUST BE CALLED before calling any other functions
TCAPIEXPORT copasi_model cCreateModel (const char *name)
 create a model
TCAPIEXPORT void cRemoveModel (copasi_model)
 remove a model
TCAPIEXPORT void cCompileModel (copasi_model model, int substitute_nested_assignments)
 clear all contents of a model
TCAPIEXPORT copasi_model cReadAntimonyFile (const char *filename)
 create a model from an Antimony or SBML file
TCAPIEXPORT copasi_model cReadSBMLFile (const char *filename)
 create a model from an SBML file
TCAPIEXPORT copasi_model cReadSBMLString (const char *sbml)
 create a model from an SBML string
TCAPIEXPORT void cWriteSBMLFile (copasi_model model, const char *filename)
 save a model as an SBML file
TCAPIEXPORT tc_matrix cGetScaledFluxControlCoeffs (copasi_model model)
 add a compartment to the model
TCAPIEXPORT copasi_compartment cCreateCompartment (copasi_model model, const char *name, double volume)
 create compartment
TCAPIEXPORT void cSetVolume (copasi_model, const char *compartment, double volume)
 set a volume of compartment
TCAPIEXPORT int cSetValue (copasi_model, const char *name, double value)
 set the concentration of a species, volume of a compartment, or value of a parameter The function will figure out which using the name (fast lookup using hashtables). If the name does not exist in the model, a new global parameter will be created.
TCAPIEXPORT void cCreateSpecies (copasi_compartment compartment, const char *name, double initialValue)
 add a species to the model
TCAPIEXPORT void cSetBoundarySpecies (copasi_model model, const char *species, int isBoundary)
 set a species as boundary or floating (will remove any assignment rules)
TCAPIEXPORT void cSetConcentration (copasi_model, const char *species, double value)
 set a species as boundary or floating (will remove any assignment rules)
TCAPIEXPORT int cSetAssignmentRule (copasi_model model, const char *species, const char *formula)
 set the assignment rule for a species (automatically assumes boundary species)
TCAPIEXPORT int cSetGlobalParameter (copasi_model model, const char *name, double value)
 set the value of an existing global parameter or create a new global parameter
TCAPIEXPORT int cCreateVariable (copasi_model model, const char *name, const char *formula)
 create a new variable that is not a constant by a formula
TCAPIEXPORT int cCreateEvent (copasi_model model, const char *name, const char *trigger, const char *variable, const char *formula)
 add a trigger and a response, where the response is defined by a target variable and an assignment formula
TCAPIEXPORT copasi_reaction cCreateReaction (copasi_model model, const char *name)
 add a species or set an existing species as fixed
TCAPIEXPORT void cAddReactant (copasi_reaction reaction, const char *species, double stoichiometry)
 add a reactant to a reaction
TCAPIEXPORT void cAddProduct (copasi_reaction reaction, const char *species, double stoichiometry)
 add a product to a reaction
TCAPIEXPORT int cSetReactionRate (copasi_reaction reaction, const char *formula)
 set reaction rate equation
TCAPIEXPORT tc_matrix cSimulateDeterministic (copasi_model model, double startTime, double endTime, int numSteps)
 simulate using LSODA numerical integrator
TCAPIEXPORT tc_matrix cSimulateStochastic (copasi_model model, double startTime, double endTime, int numSteps)
 simulate using exact stochastic algorithm
TCAPIEXPORT tc_matrix cSimulateHybrid (copasi_model model, double startTime, double endTime, int numSteps)
 simulate using Hybrid algorithm/deterministic algorithm
TCAPIEXPORT tc_matrix cSimulateTauLeap (copasi_model model, double startTime, double endTime, int numSteps)
 simulate using Tau Leap stochastic algorithm
TCAPIEXPORT tc_matrix cGetSteadyState (copasi_model model)
 bring the system to steady state
TCAPIEXPORT tc_matrix cGetSteadyState2 (copasi_model model, int iter)
 bring the system to steady state using normal simulation
TCAPIEXPORT tc_matrix cGetJacobian (copasi_model model)
 get the Jacobian at the current state
TCAPIEXPORT tc_matrix cGetEigenvalues (copasi_model model)
 get the eigenvalues of the Jacobian at the current state
TCAPIEXPORT tc_matrix cGetUnscaledElasticities (copasi_model model)
 unscaled elasticities
TCAPIEXPORT tc_matrix cGetUnscaledConcentrationControlCoeffs (copasi_model model)
 unscaled elasticities
TCAPIEXPORT tc_matrix cGetUnscaledFluxControlCoeffs (copasi_model model)
 unscaled flux control coefficients
TCAPIEXPORT tc_matrix cGetScaledElasticities (copasi_model model)
 scaled elasticities
TCAPIEXPORT tc_matrix cGetScaledConcentrationConcentrationCoeffs (copasi_model model)
 scaled concentration control coefficients
TCAPIEXPORT tc_matrix cGetFullStoichiometryMatrix (copasi_model model)
 full stoichiometry matrix
TCAPIEXPORT tc_matrix cGetReducedStoichiometryMatrix (copasi_model model)
 reduced stoichiometry matrix
TCAPIEXPORT tc_matrix cGetElementaryFluxModes (copasi_model model)
 elementary flux modes
TCAPIEXPORT tc_matrix cGetGammaMatrix (copasi_model model)
 get Gamma matrix (i.e. conservation laws)
TCAPIEXPORT tc_matrix cGetKMatrix (copasi_model model)
 get K matrix (right nullspace)
TCAPIEXPORT tc_matrix cGetK0Matrix (copasi_model model)
 get K0 matrix
TCAPIEXPORT tc_matrix cGetLinkMatrix (copasi_model model)
 get L matrix (left nullspace)
TCAPIEXPORT tc_matrix cGetL0Matrix (copasi_model model)
 get L0 matrix
TCAPIEXPORT tc_matrix cOptimize (copasi_model model, const char *objective, tc_matrix input)
 fit the model parameters to time-series data
TCAPIEXPORT void cSetOptimizerIterations (int)
TCAPIEXPORT void cSetOptimizerSize (int)
TCAPIEXPORT void cSetOptimizerMutationRate (double)
TCAPIEXPORT void cSetOptimizerCrossoverRate (double)

Function Documentation

TCAPIEXPORT void cAddProduct ( copasi_reaction  reaction,
const char *  species,
double  stoichiometry 
)

add a product to a reaction

Parameters:
copasi_reactionreaction
char* product
doublestoichiometry
TCAPIEXPORT void cAddReactant ( copasi_reaction  reaction,
const char *  species,
double  stoichiometry 
)

add a reactant to a reaction

Parameters:
copasi_reactionreaction
char* reactant
doublestoichiometry
TCAPIEXPORT void cCompileModel ( copasi_model  model,
int  substitute_nested_assignments 
)

clear all contents of a model

This function is only needed for calling COPASI methods not found in this library. This function compiles the COPASI model; it is called internally by the simulate and other anlysis functions.

Parameters:
copasi_modelmodel
intsubstitute nested assignments
TCAPIEXPORT copasi_compartment cCreateCompartment ( copasi_model  model,
const char *  name,
double  volume 
)

create compartment

Parameters:
char*compartment name
doublevolume
Returns:
copasi_compartment a new compartment
TCAPIEXPORT int cCreateEvent ( copasi_model  model,
const char *  name,
const char *  trigger,
const char *  variable,
const char *  formula 
)

add a trigger and a response, where the response is defined by a target variable and an assignment formula

Parameters:
copasi_modelmodel
char* event name
char* trigger
char* response: name of variable or species
char*response: assignment formula
Returns:
int 0=failed 1=success
TCAPIEXPORT copasi_model cCreateModel ( const char *  name)

create a model

Parameters:
char*model name
Returns:
copasi_model a new copasi model
TCAPIEXPORT copasi_reaction cCreateReaction ( copasi_model  model,
const char *  name 
)

add a species or set an existing species as fixed

Parameters:
copasi_modelmodel
char*species name
Returns:
copasi_reaction a new reaction
TCAPIEXPORT void cCreateSpecies ( copasi_compartment  compartment,
const char *  name,
double  initialValue 
)

add a species to the model

Parameters:
copasi_compartmentmodel
char*species name
doubleinitial value (concentration or count, depending on the model)
TCAPIEXPORT int cCreateVariable ( copasi_model  model,
const char *  name,
const char *  formula 
)

create a new variable that is not a constant by a formula

Parameters:
copasi_modelmodel
char*name of new variable
char*formula
Returns:
int 0=failed 1=success
TCAPIEXPORT tc_matrix cGetEigenvalues ( copasi_model  model)

get the eigenvalues of the Jacobian at the current state

Parameters:
copasi_modelmodel
Returns:
tc_matrix matrix with 1 row and n columns, each containing an eigenvalue
TCAPIEXPORT tc_matrix cGetElementaryFluxModes ( copasi_model  model)

elementary flux modes

Parameters:
copasi_modelmodel
Returns:
tc_matrix matrix with reactions as rows (with rownames) and flux modes as columns (no column names)
TCAPIEXPORT tc_matrix cGetFullStoichiometryMatrix ( copasi_model  model)

full stoichiometry matrix

Parameters:
copasi_modelmodel
Returns:
tc_matrix
TCAPIEXPORT tc_matrix cGetGammaMatrix ( copasi_model  model)

get Gamma matrix (i.e. conservation laws)

Parameters:
copasi_modelmodel
Returns:
tc_matrix
TCAPIEXPORT tc_matrix cGetJacobian ( copasi_model  model)

get the Jacobian at the current state

Parameters:
copasi_modelmodel
Returns:
tc_matrix matrix with n rows and n columns, where n = number of species
TCAPIEXPORT tc_matrix cGetK0Matrix ( copasi_model  model)

get K0 matrix

Parameters:
copasi_modelmodel
Returns:
tc_matrix
TCAPIEXPORT tc_matrix cGetKMatrix ( copasi_model  model)

get K matrix (right nullspace)

Parameters:
copasi_modelmodel
Returns:
tc_matrix
TCAPIEXPORT tc_matrix cGetL0Matrix ( copasi_model  model)

get L0 matrix

Parameters:
copasi_modelmodel
Returns:
tc_matrix
TCAPIEXPORT tc_matrix cGetLinkMatrix ( copasi_model  model)

get L matrix (left nullspace)

Parameters:
copasi_modelmodel
Returns:
tc_matrix
TCAPIEXPORT tc_matrix cGetReducedStoichiometryMatrix ( copasi_model  model)

reduced stoichiometry matrix

Parameters:
copasi_modelmodel
Returns:
tc_matrix
TCAPIEXPORT tc_matrix cGetScaledConcentrationConcentrationCoeffs ( copasi_model  model)

scaled concentration control coefficients

Parameters:
copasi_modelmodel
Returns:
tc_matrix
TCAPIEXPORT tc_matrix cGetScaledElasticities ( copasi_model  model)

scaled elasticities

Parameters:
copasi_modelmodel
Returns:
tc_matrix
TCAPIEXPORT tc_matrix cGetScaledFluxControlCoeffs ( copasi_model  model)

add a compartment to the model

scaled flux control coefficients

Parameters:
copasi_modelmodel/*! scaled flux control coefficients
copasi_modelmodel
Returns:
tc_matrix
Parameters:
copasi_modelmodel
Returns:
tc_matrix

add a compartment to the model

Parameters:
copasi_modelmodel
Returns:
tc_matrix
TCAPIEXPORT tc_matrix cGetSteadyState ( copasi_model  model)

bring the system to steady state

Parameters:
copasi_modelmodel
Returns:
tc_matrix matrix with 1 row and n columns, where n = number of species
TCAPIEXPORT tc_matrix cGetSteadyState2 ( copasi_model  model,
int  iter 
)

bring the system to steady state using normal simulation

Parameters:
copasi_modelmodel
intmax iterations (each iteration doubles the time duration)
Returns:
tc_matrix matrix with 1 row and n columns, where n = number of species
TCAPIEXPORT tc_matrix cGetUnscaledConcentrationControlCoeffs ( copasi_model  model)

unscaled elasticities

unscaled concentration control coefficients

Parameters:
copasi_modelmodel
Returns:
tc_matrix
TCAPIEXPORT tc_matrix cGetUnscaledElasticities ( copasi_model  model)

unscaled elasticities

Parameters:
copasi_modelmodel
Returns:
tc_matrix
TCAPIEXPORT tc_matrix cGetUnscaledFluxControlCoeffs ( copasi_model  model)

unscaled flux control coefficients

Parameters:
copasi_modelmodel
Returns:
tc_matrix
BEGIN_C_DECLS TCAPIEXPORT void copasi_end ( )

initialize copasi -- MUST BE CALLED before calling any other functions

destroy copasi -- MUST BE CALLED at the end of program

TCAPIEXPORT tc_matrix cOptimize ( copasi_model  model,
const char *  objective,
tc_matrix  input 
)

fit the model parameters to time-series data

Parameters:
copasi_modelmodel
char* filename (tab separated)
tc_matrixparameters to optimize. rownames should contain parameter names, column 1 contains parameter min-values, and column 2 contains parameter max values
char* pick method. Use of of the following: "GeneticAlgorithm", "LevenbergMarquardt", "SimulatedAnnealing", "NelderMead", "SRES", "ParticleSwarm", "SteepestDescent", "RandomSearch"

use genetic algorithms to generate a distribution of parameter values that satisfy an objective function or fit a data file

Parameters:
copasi_modelmodel
char* objective function or filename
tc_matrixparameter initial values and min and max values (3 columns)
TCAPIEXPORT copasi_model cReadAntimonyFile ( const char *  filename)

create a model from an Antimony or SBML file

Parameters:
char*file name
Returns:
copasi_model a new copasi model
TCAPIEXPORT copasi_model cReadSBMLFile ( const char *  filename)

create a model from an SBML file

Parameters:
char*file name
Returns:
copasi_model a new copasi model
TCAPIEXPORT copasi_model cReadSBMLString ( const char *  sbml)

create a model from an SBML string

Parameters:
char*SBML string
Returns:
copasi_model a new copasi model
TCAPIEXPORT void cRemoveModel ( copasi_model  )

remove a model

TCAPIEXPORT int cSetAssignmentRule ( copasi_model  model,
const char *  species,
const char *  formula 
)

set the assignment rule for a species (automatically assumes boundary species)

Parameters:
copasi_modelmodel
char* species name
char*formula, use 0 to remove assignment rule
Returns:
int 0=failed 1=success
TCAPIEXPORT void cSetBoundarySpecies ( copasi_model  model,
const char *  species,
int  isBoundary 
)

set a species as boundary or floating (will remove any assignment rules)

Parameters:
copasi_modelmodel
char* name
intboundary = 1, floating = 0 (default)
TCAPIEXPORT void cSetConcentration ( copasi_model  ,
const char *  species,
double  value 
)

set a species as boundary or floating (will remove any assignment rules)

Parameters:
copasi_modelmodel
char* species name
doubleconcentration or count
TCAPIEXPORT int cSetGlobalParameter ( copasi_model  model,
const char *  name,
double  value 
)

set the value of an existing global parameter or create a new global parameter

Parameters:
copasi_modelmodel
char*parameter name
doublevalue
Returns:
int 0=new value created 1=found existing value
TCAPIEXPORT void cSetOptimizerCrossoverRate ( double  )
TCAPIEXPORT void cSetOptimizerIterations ( int  )
TCAPIEXPORT void cSetOptimizerMutationRate ( double  )
TCAPIEXPORT void cSetOptimizerSize ( int  )
TCAPIEXPORT int cSetReactionRate ( copasi_reaction  reaction,
const char *  formula 
)

set reaction rate equation

Parameters:
copasi_reactionreaction
char*custom formula
Returns:
int success=1 failure=0
TCAPIEXPORT int cSetValue ( copasi_model  ,
const char *  name,
double  value 
)

set the concentration of a species, volume of a compartment, or value of a parameter The function will figure out which using the name (fast lookup using hashtables). If the name does not exist in the model, a new global parameter will be created.

Parameters:
copasi_modelmodel
char* name
doublevalue
Returns:
0 if new variable was created. 1 if existing variable was found
TCAPIEXPORT void cSetVolume ( copasi_model  ,
const char *  compartment,
double  volume 
)

set a volume of compartment

Parameters:
copasi_modelmodel
char* compartment name
doublevolume
TCAPIEXPORT tc_matrix cSimulateDeterministic ( copasi_model  model,
double  startTime,
double  endTime,
int  numSteps 
)

simulate using LSODA numerical integrator

Parameters:
copasi_modelmodel
doublestart time
doubleend time
intnumber of steps in the output
Returns:
tc_matrix matrix of concentration or particles
TCAPIEXPORT tc_matrix cSimulateHybrid ( copasi_model  model,
double  startTime,
double  endTime,
int  numSteps 
)

simulate using Hybrid algorithm/deterministic algorithm

Parameters:
copasi_modelmodel
doublestart time
doubleend time
intnumber of steps in the output
Returns:
tc_matrix matrix of concentration or particles
TCAPIEXPORT tc_matrix cSimulateStochastic ( copasi_model  model,
double  startTime,
double  endTime,
int  numSteps 
)

simulate using exact stochastic algorithm

Parameters:
copasi_modelmodel
doublestart time
doubleend time
intnumber of steps in the output
Returns:
tc_matrix matrix of concentration or particles
TCAPIEXPORT tc_matrix cSimulateTauLeap ( copasi_model  model,
double  startTime,
double  endTime,
int  numSteps 
)

simulate using Tau Leap stochastic algorithm

Parameters:
copasi_modelmodel
doublestart time
doubleend time
intnumber of steps in the output
Returns:
tc_matrix matrix of concentration or particles
TCAPIEXPORT void cWriteSBMLFile ( copasi_model  model,
const char *  filename 
)

save a model as an SBML file

Parameters:
copasi_modelcopasi model
char*file name
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines