copasi API  0.1
Classes | Defines | Typedefs | Functions
copasi/copasi_api.cpp File Reference
#include <iostream>
#include <vector>
#include <string>
#include <set>
#include <limits>
#include <QString>
#include <QStringList>
#include <QRegExp>
#include <QHash>
#include <QList>
#include <QPair>
#include <QFile>
#include "copasi_api.h"
#include "copasi/copasi.h"
#include "copasi/report/CCopasiRootContainer.h"
#include "copasi/CopasiDataModel/CCopasiDataModel.h"
#include "copasi/model/CModel.h"
#include "copasi/model/CCompartment.h"
#include "copasi/model/CMetab.h"
#include "copasi/model/CReaction.h"
#include "copasi/model/CChemEq.h"
#include "copasi/model/CModelValue.h"
#include "copasi/function/CFunctionDB.h"
#include "copasi/function/CFunction.h"
#include "copasi/function/CEvaluationTree.h"
#include "copasi/report/CReport.h"
#include "copasi/report/CReportDefinition.h"
#include "copasi/report/CReportDefinitionVector.h"
#include "copasi/trajectory/CTrajectoryTask.h"
#include "copasi/trajectory/CTrajectoryMethod.h"
#include "copasi/trajectory/CTrajectoryProblem.h"
#include "copasi/scan/CScanTask.h"
#include "copasi/scan/CScanMethod.h"
#include "copasi/scan/CScanProblem.h"
#include "copasi/trajectory/CTimeSeries.h"
#include "copasi/steadystate/CSteadyStateTask.h"
#include "copasi/steadystate/CMCATask.h"
#include "copasi/steadystate/CMCAMethod.h"
#include "copasi/elementaryFluxModes/CFluxMode.h"
#include "copasi/elementaryFluxModes/CEFMTask.h"
#include "copasi/elementaryFluxModes/CEFMProblem.h"
#include "copasi/commandline/COptions.h"
#include "copasi/report/CCopasiContainer.h"
#include "copasi/parameterFitting/CFitTask.h"
#include "copasi/parameterFitting/CFitMethod.h"
#include "copasi/parameterFitting/CFitProblem.h"
#include "copasi/parameterFitting/CFitItem.h"
#include "copasi/parameterFitting/CExperimentSet.h"
#include "copasi/parameterFitting/CExperiment.h"
#include "copasi/parameterFitting/CExperimentObjectMap.h"
#include "copasi/report/CKeyFactory.h"
#include "GASStateGA.h"
#include "GA1DArrayGenome.h"
#include "libstructural.h"
#include "matrix.h"
#include "muParserDef.h"
#include "muParser.h"
#include "muParserInt.h"
#include "mtrand.h"
#include "src/antimony_api.h"

Classes

struct  CopasiPtr
struct  GAData

Defines

#define COPASI_MAIN   1
#define LIB_EXPORTS   1

Typedefs

typedef QHash< QString, CopasiPtrCQHash
typedef GA1DArrayGenome< float > RealGenome

Functions

void copasi_init ()
void copasi_end ()
 initialize copasi -- MUST BE CALLED before calling any other functions
int cSetAssignmentRuleHelper (copasi_model, CMetab *, const char *)
int copasi_cleanup_assignments (copasi_model model, bool doWhile=false)
void cRemoveModel (copasi_model model)
 remove a model
void clearCopasiModel (copasi_model model)
copasi_model cCreateModel (const char *name)
 create a model
void cCreateSpecies (copasi_compartment compartment, const char *name, double iv)
 add a species to the model
copasi_compartment cCreateCompartment (copasi_model model, const char *name, double volume)
 create compartment
int cSetValue (copasi_model 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.
void cSetVolume (copasi_model model, const char *name, double vol)
 set a volume of compartment
void cSetConcentration (copasi_model model, const char *name, double conc)
 set a species as boundary or floating (will remove any assignment rules)
int cSetGlobalParameter (copasi_model model, const char *name, double value)
 set the value of an existing global parameter or create a new global parameter
void cSetBoundarySpecies (copasi_model model, const char *name, int isBoundary)
 set a species as boundary or floating (will remove any assignment rules)
int cSetAssignmentRule (copasi_model model, const char *name, const char *formula)
 set the assignment rule for a species (automatically assumes boundary species)
int cCreateVariable (copasi_model model, const char *name, const char *formula)
 create a new variable that is not a constant by a formula
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
copasi_reaction cCreateReaction (copasi_model model, const char *name)
 add a species or set an existing species as fixed
void cAddReactant (copasi_reaction reaction, const char *species, double stoichiometry)
 add a reactant to a reaction
void cAddProduct (copasi_reaction reaction, const char *species, double stoichiometry)
 add a product to a reaction
int cSetReactionRate (copasi_reaction reaction, const char *formula)
 set reaction rate equation
void cCompileModel (copasi_model model, int subs)
 clear all contents of a model
tc_matrix simulate (copasi_model model, double startTime, double endTime, int numSteps, CCopasiMethod::SubType method)
tc_matrix cSimulateDeterministic (copasi_model model, double startTime, double endTime, int numSteps)
 simulate using LSODA numerical integrator
tc_matrix cSimulateTauLeap (copasi_model model, double startTime, double endTime, int numSteps)
 simulate using Tau Leap stochastic algorithm
tc_matrix cSimulateStochastic (copasi_model model, double startTime, double endTime, int numSteps)
 simulate using exact stochastic algorithm
tc_matrix cSimulateHybrid (copasi_model model, double startTime, double endTime, int numSteps)
 simulate using Hybrid algorithm/deterministic algorithm
void cWriteSBMLFile (copasi_model model, const char *filename)
 save a model as an SBML file
copasi_model cReadAntimonyFile (const char *filename)
 create a model from an Antimony or SBML file
copasi_model cReadSBMLFile (const char *filename)
 create a model from an SBML file
copasi_model cReadSBMLString (const char *sbml)
 create a model from an SBML string
tc_matrix cGetJacobian (copasi_model model)
 get the Jacobian at the current state
tc_matrix cGetSteadyState2 (copasi_model model, int maxiter)
 bring the system to steady state using normal simulation
tc_matrix cGetSteadyState (copasi_model model)
 bring the system to steady state
tc_matrix cGetEigenvalues (copasi_model model)
 get the eigenvalues of the Jacobian at the current state
tc_matrix cGetUnscaledElasticities (copasi_model model)
 unscaled elasticities
tc_matrix cGetUnscaledConcentrationControlCoeffs (copasi_model model)
 unscaled elasticities
tc_matrix cGetUnscaledFluxControlCoeffs (copasi_model model)
 unscaled flux control coefficients
tc_matrix cGetScaledElasticities (copasi_model model)
 scaled elasticities
tc_matrix cGetScaledConcentrationConcentrationCoeffs (copasi_model model)
 scaled concentration control coefficients
tc_matrix cGetScaledFluxControlCoeffs (copasi_model model)
 scaled flux control coefficients
tc_matrix cGetReducedStoichiometryMatrix (copasi_model model)
 reduced stoichiometry matrix
tc_matrix cGetFullStoichiometryMatrix (copasi_model model)
 full stoichiometry matrix
tc_matrix cGetElementaryFluxModes (copasi_model model)
 elementary flux modes
void InitializeGenome (GAGenome &x)
tc_matrix cOptimize (copasi_model model, const char *objective, tc_matrix params)
 fit the model parameters to time-series data
void cSetOptimizerSize (int n)
void cSetOptimizerIterations (int n)
void cSetOptimizerMutationRate (double q)
void cSetOptimizerCrossoverRate (double q)
tc_matrix convertFromDoubleMatrix (DoubleMatrix &matrix, std::vector< std::string > &rowNames, std::vector< std::string > &colNames)
void convertToDoubleMatrix (tc_matrix m, DoubleMatrix &matrix, std::vector< std::string > &rowNames, std::vector< std::string > &colNames)
tc_matrix cGetGammaMatrix (copasi_model model)
 get Gamma matrix (i.e. conservation laws)
tc_matrix cGetKMatrix (copasi_model model)
 get K matrix (right nullspace)
tc_matrix cGetLinkMatrix (copasi_model model)
 get L matrix (left nullspace)
tc_matrix cGetK0Matrix (copasi_model model)
 get K0 matrix
tc_matrix cGetL0Matrix (copasi_model model)
 get L0 matrix

Define Documentation

#define COPASI_MAIN   1
#define LIB_EXPORTS   1

Typedef Documentation

typedef QHash< QString, CopasiPtr > CQHash
typedef GA1DArrayGenome<float> RealGenome

Function Documentation

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

add a product to a reaction

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

add a reactant to a reaction

Parameters:
copasi_reactionreaction
char* reactant
doublestoichiometry
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
copasi_compartment cCreateCompartment ( copasi_model  model,
const char *  name,
double  volume 
)

create compartment

Parameters:
char*compartment name
doublevolume
Returns:
copasi_compartment a new compartment
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
copasi_model cCreateModel ( const char *  name)

create a model

Parameters:
char*model name
Returns:
copasi_model a new copasi model
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
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)
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
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
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)
tc_matrix cGetFullStoichiometryMatrix ( copasi_model  model)

full stoichiometry matrix

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

get Gamma matrix (i.e. conservation laws)

Parameters:
copasi_modelmodel
Returns:
tc_matrix
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
tc_matrix cGetK0Matrix ( copasi_model  model)

get K0 matrix

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

get K matrix (right nullspace)

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

get L0 matrix

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

get L matrix (left nullspace)

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

reduced stoichiometry matrix

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

scaled concentration control coefficients

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

scaled elasticities

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

scaled flux control coefficients

add a compartment to the model

Parameters:
copasi_modelmodel
Returns:
tc_matrix
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
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
tc_matrix cGetUnscaledConcentrationControlCoeffs ( copasi_model  model)

unscaled elasticities

unscaled concentration control coefficients

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

unscaled elasticities

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

unscaled flux control coefficients

Parameters:
copasi_modelmodel
Returns:
tc_matrix
void clearCopasiModel ( copasi_model  model)
tc_matrix convertFromDoubleMatrix ( DoubleMatrix &  matrix,
std::vector< std::string > &  rowNames,
std::vector< std::string > &  colNames 
)
void convertToDoubleMatrix ( tc_matrix  m,
DoubleMatrix &  matrix,
std::vector< std::string > &  rowNames,
std::vector< std::string > &  colNames 
)
int copasi_cleanup_assignments ( copasi_model  model,
bool  doWhile = false 
)
void copasi_end ( )

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

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

void copasi_init ( )
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)
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
copasi_model cReadSBMLFile ( const char *  filename)

create a model from an SBML file

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

create a model from an SBML string

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

remove a model

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
int cSetAssignmentRuleHelper ( copasi_model  model,
CMetab pSpecies,
const char *  formula 
)
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)
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
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
void cSetOptimizerCrossoverRate ( double  q)
void cSetOptimizerIterations ( int  n)
void cSetOptimizerMutationRate ( double  q)
void cSetOptimizerSize ( int  n)
int cSetReactionRate ( copasi_reaction  reaction,
const char *  formula 
)

set reaction rate equation

Parameters:
copasi_reactionreaction
char*custom formula
Returns:
int success=1 failure=0
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
void cSetVolume ( copasi_model  ,
const char *  compartment,
double  volume 
)

set a volume of compartment

Parameters:
copasi_modelmodel
char* compartment name
doublevolume
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
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
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
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
void cWriteSBMLFile ( copasi_model  model,
const char *  filename 
)

save a model as an SBML file

Parameters:
copasi_modelcopasi model
char*file name
void InitializeGenome ( GAGenome &  x)
tc_matrix simulate ( copasi_model  model,
double  startTime,
double  endTime,
int  numSteps,
CCopasiMethod::SubType  method 
)
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines