|
copasi API
0.1
|
#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, CopasiPtr > | CQHash |
| 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 COPASI_MAIN 1 |
| #define LIB_EXPORTS 1 |
| typedef GA1DArrayGenome<float> RealGenome |
| void cAddProduct | ( | copasi_reaction | reaction, |
| const char * | species, | ||
| double | stoichiometry | ||
| ) |
add a product to a reaction
| copasi_reaction | reaction |
| char | * product |
| double | stoichiometry |
| void cAddReactant | ( | copasi_reaction | reaction, |
| const char * | species, | ||
| double | stoichiometry | ||
| ) |
add a reactant to a reaction
| copasi_reaction | reaction |
| char | * reactant |
| double | stoichiometry |
| 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.
| copasi_model | model |
| int | substitute nested assignments |
| copasi_compartment cCreateCompartment | ( | copasi_model | model, |
| const char * | name, | ||
| double | volume | ||
| ) |
create compartment
| char* | compartment name |
| double | volume |
| 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_model | model |
| char | * event name |
| char | * trigger |
| char | * response: name of variable or species |
| char* | response: assignment formula |
| copasi_model cCreateModel | ( | const char * | name | ) |
| copasi_reaction cCreateReaction | ( | copasi_model | model, |
| const char * | name | ||
| ) |
add a species or set an existing species as fixed
| copasi_model | model |
| char* | species name |
| void cCreateSpecies | ( | copasi_compartment | compartment, |
| const char * | name, | ||
| double | initialValue | ||
| ) |
add a species to the model
| copasi_compartment | model |
| char* | species name |
| double | initial 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
| copasi_model | model |
| char* | name of new variable |
| char* | formula |
| tc_matrix cGetEigenvalues | ( | copasi_model | model | ) |
get the eigenvalues of the Jacobian at the current state
| copasi_model | model |
| tc_matrix cGetElementaryFluxModes | ( | copasi_model | model | ) |
elementary flux modes
| copasi_model | model |
| tc_matrix cGetFullStoichiometryMatrix | ( | copasi_model | model | ) |
| tc_matrix cGetGammaMatrix | ( | copasi_model | model | ) |
| tc_matrix cGetJacobian | ( | copasi_model | model | ) |
get the Jacobian at the current state
| copasi_model | model |
| tc_matrix cGetK0Matrix | ( | copasi_model | model | ) |
| tc_matrix cGetKMatrix | ( | copasi_model | model | ) |
| tc_matrix cGetL0Matrix | ( | copasi_model | model | ) |
| tc_matrix cGetLinkMatrix | ( | copasi_model | model | ) |
| tc_matrix cGetReducedStoichiometryMatrix | ( | copasi_model | model | ) |
| tc_matrix cGetScaledConcentrationConcentrationCoeffs | ( | copasi_model | model | ) |
| tc_matrix cGetScaledElasticities | ( | copasi_model | model | ) |
| TCAPIEXPORT tc_matrix cGetScaledFluxControlCoeffs | ( | copasi_model | model | ) |
scaled flux control coefficients
add a compartment to the model
| copasi_model | model |
| tc_matrix cGetSteadyState | ( | copasi_model | model | ) |
bring the system to steady state
| copasi_model | model |
| tc_matrix cGetSteadyState2 | ( | copasi_model | model, |
| int | iter | ||
| ) |
bring the system to steady state using normal simulation
| copasi_model | model |
| int | max iterations (each iteration doubles the time duration) |
| tc_matrix cGetUnscaledConcentrationControlCoeffs | ( | copasi_model | model | ) |
unscaled elasticities
unscaled concentration control coefficients
| copasi_model | model |
| tc_matrix cGetUnscaledElasticities | ( | copasi_model | model | ) |
| tc_matrix cGetUnscaledFluxControlCoeffs | ( | copasi_model | model | ) |
| 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
| copasi_model | model |
| char | * filename (tab separated) |
| tc_matrix | parameters 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
| copasi_model | model |
| char | * objective function or filename |
| tc_matrix | parameter initial values and min and max values (3 columns) |
| copasi_model cReadAntimonyFile | ( | const char * | filename | ) |
create a model from an Antimony or SBML file
| char* | file name |
| copasi_model cReadSBMLFile | ( | const char * | filename | ) |
| copasi_model cReadSBMLString | ( | const char * | sbml | ) |
create a model from an SBML string
| char* | SBML string |
| 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)
| copasi_model | model |
| char | * species name |
| char* | formula, use 0 to remove assignment rule |
| 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)
| copasi_model | model |
| char | * name |
| int | boundary = 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)
| copasi_model | model |
| char | * species name |
| double | concentration 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
| copasi_model | model |
| char* | parameter name |
| double | 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
| copasi_reaction | reaction |
| char* | custom formula |
| 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.
| copasi_model | model |
| char | * name |
| double | value |
| void cSetVolume | ( | copasi_model | , |
| const char * | compartment, | ||
| double | volume | ||
| ) |
set a volume of compartment
| copasi_model | model |
| char | * compartment name |
| double | volume |
| tc_matrix cSimulateDeterministic | ( | copasi_model | model, |
| double | startTime, | ||
| double | endTime, | ||
| int | numSteps | ||
| ) |
simulate using LSODA numerical integrator
| copasi_model | model |
| double | start time |
| double | end time |
| int | number of steps in the output |
| tc_matrix cSimulateHybrid | ( | copasi_model | model, |
| double | startTime, | ||
| double | endTime, | ||
| int | numSteps | ||
| ) |
simulate using Hybrid algorithm/deterministic algorithm
| copasi_model | model |
| double | start time |
| double | end time |
| int | number of steps in the output |
| tc_matrix cSimulateStochastic | ( | copasi_model | model, |
| double | startTime, | ||
| double | endTime, | ||
| int | numSteps | ||
| ) |
simulate using exact stochastic algorithm
| copasi_model | model |
| double | start time |
| double | end time |
| int | number of steps in the output |
| tc_matrix cSimulateTauLeap | ( | copasi_model | model, |
| double | startTime, | ||
| double | endTime, | ||
| int | numSteps | ||
| ) |
simulate using Tau Leap stochastic algorithm
| copasi_model | model |
| double | start time |
| double | end time |
| int | number of steps in the output |
| void cWriteSBMLFile | ( | copasi_model | model, |
| const char * | filename | ||
| ) |
save a model as an SBML file
| copasi_model | copasi model |
| char* | file name |
| void InitializeGenome | ( | GAGenome & | x | ) |
| tc_matrix simulate | ( | copasi_model | model, |
| double | startTime, | ||
| double | endTime, | ||
| int | numSteps, | ||
| CCopasiMethod::SubType | method | ||
| ) |
1.7.5.1