copasi API  0.1
Classes | Public Types | Public Member Functions | Static Public Attributes
CModel Class Reference

#include <CModel.h>

Inheritance diagram for CModel:
CModelEntity CCopasiContainer CAnnotation CCopasiObject

List of all members.

Classes

class  CLinkMatrixView

Public Types

enum  VolumeUnit {
  dimensionlessVolume = 0, m3, l, ml,
  microl, nl, pl, fl
}
enum  AreaUnit {
  dimensionlessArea = 0, m2, dm2, cm2,
  mm2, microm2, nm2, pm2,
  fm2
}
enum  LengthUnit {
  dimensionlessLength = 0, m, dm, cm,
  mm, microm, nm, pm,
  fm
}
enum  TimeUnit {
  dimensionlessTime = 0, d, h, min,
  s, ms, micros, ns,
  ps, fs, OldMinute
}
enum  QuantityUnit {
  dimensionlessQuantity = 0, Mol, mMol, microMol,
  nMol, pMol, fMol, number,
  OldXML
}
enum  ModelType { deterministic = 0, stochastic }

Public Member Functions

 CModel (CCopasiContainer *pParent)
 ~CModel ()
virtual std::string getChildObjectUnits (const CCopasiObject *pObject) const
void cleanup ()
bool convert2NonReversible ()
C_INT32 load (CReadConfig &configBuffer)
void initializeMetabolites ()
void setCompileFlag (bool flag=true)
bool compileIfNecessary (CProcessReport *pProcessReport)
bool forceCompile (CProcessReport *pProcessReport)
void compileDefaultMetabInitialValueDependencies ()
void buildStoi ()
void buildLinkZero ()
void buildRedStoi ()
void buildMoieties ()
const CCopasiVector< CMetab > & getMetabolites () const
CCopasiVector< CMetab > & getMetabolites ()
const CCopasiVector< CMetab > & getMetabolitesX () const
CCopasiVector< CMetab > & getMetabolitesX ()
unsigned C_INT32 getNumMetabs () const
unsigned C_INT32 getNumVariableMetabs () const
unsigned C_INT32 getNumODEMetabs () const
unsigned C_INT32 getNumAssignmentMetabs () const
unsigned C_INT32 getNumIndependentReactionMetabs () const
unsigned C_INT32 getNumDependentReactionMetabs () const
const CCopasiVectorN
< CModelValue > & 
getModelValues () const
CCopasiVectorN< CModelValue > & getModelValues ()
unsigned C_INT32 getNumModelValues () const
CCopasiVectorNS< CReaction > & getReactions ()
const CCopasiVectorNS
< CReaction > & 
getReactions () const
unsigned C_INT32 getTotSteps () const
const CVector< C_FLOAT64 > & getParticleFlux () const
CCopasiVectorN< CEvent > & getEvents ()
const CCopasiVectorN< CEvent > & getEvents () const
const std::string & getKey () const
bool setTitle (const std::string &title)
void setInitialTime (const C_FLOAT64 &time)
const C_FLOAT64 & getInitialTime () const
void setTime (const C_FLOAT64 &time)
const C_FLOAT64 & getTime () const
CCopasiVectorNS< CCompartment > & getCompartments ()
const CCopasiVectorNS
< CCompartment > & 
getCompartments () const
const CMatrix< C_FLOAT64 > & getStoi () const
const CMatrix< C_FLOAT64 > & getRedStoi () const
const CMatrix< C_FLOAT64 > & getStoiReordered () const
const CCopasiVector< CMoiety > & getMoieties () const
unsigned C_INT32 findMetabByName (const std::string &Target) const
unsigned C_INT32 findMoiety (const std::string &Target) const
const CLinkMatrixViewgetL () const
const CMatrix< C_FLOAT64 > & getL0 () const
void applyInitialValues ()
bool updateInitialValues ()
const CStategetInitialState () const
const CStategetState () const
void setInitialState (const CState &state)
void setState (const CState &state)
void updateSimulatedValues (const bool &updateMoieties)
void updateNonSimulatedValues (void)
void calculateDerivatives (C_FLOAT64 *derivatives)
void calculateDerivativesX (C_FLOAT64 *derivativesX)
void calculateElasticityMatrix (const C_FLOAT64 &factor, const C_FLOAT64 &resolution)
void calculateJacobian (CMatrix< C_FLOAT64 > &jacobian, const C_FLOAT64 &derivationFactor, const C_FLOAT64 &resolution)
void calculateJacobianX (CMatrix< C_FLOAT64 > &jacobianX, const C_FLOAT64 &derivationFactor, const C_FLOAT64 &resolution)
C_FLOAT64 calculateDivergence () const
bool setVolumeUnit (const std::string &name)
bool setVolumeUnit (const CModel::VolumeUnit &unit)
std::string getVolumeUnitName () const
CModel::VolumeUnit getVolumeUnitEnum () const
bool setAreaUnit (const std::string &name)
bool setAreaUnit (const CModel::AreaUnit &unit)
std::string getAreaUnitName () const
CModel::AreaUnit getAreaUnitEnum () const
bool setLengthUnit (const std::string &name)
bool setLengthUnit (const CModel::LengthUnit &unit)
std::string getLengthUnitName () const
CModel::LengthUnit getLengthUnitEnum () const
bool setTimeUnit (const std::string &name)
bool setTimeUnit (const CModel::TimeUnit &unit)
std::string getTimeUnitName () const
CModel::TimeUnit getTimeUnitEnum () const
bool setQuantityUnit (const std::string &name)
bool setQuantityUnit (const CModel::QuantityUnit &unit)
std::string getQuantityUnitName () const
std::string getQuantityUnitOldXMLName () const
CModel::QuantityUnit getQuantityUnitEnum () const
void setModelType (const ModelType &modelType)
const ModelTypegetModelType () const
void setAvogadro (const C_FLOAT64 &avogadro)
const C_FLOAT64 & getAvogadro () const
const C_FLOAT64 & getQuantity2NumberFactor () const
const C_FLOAT64 & getNumber2QuantityFactor () const
CMetabcreateMetabolite (const std::string &name, const std::string &compartment, const C_FLOAT64 &iconc=1.0, const CModelEntity::Status &status=CModelEntity::REACTIONS)
bool removeMetabolite (const std::string &key, const bool &recursive=true)
bool removeMetabolite (const unsigned C_INT32 index, const bool &recursive=true)
bool removeMetabolite (const CMetab *pMetabolite, const bool &recursive=true)
bool appendDependentModelObjects (const std::set< const CCopasiObject * > &candidates, std::set< const CCopasiObject * > &dependentReactions, std::set< const CCopasiObject * > &dependentMetabolites, std::set< const CCopasiObject * > &dependentCompartments, std::set< const CCopasiObject * > &dependentModelValues, std::set< const CCopasiObject * > &dependentEvents) const
bool appendDependentReactions (std::set< const CCopasiObject * > candidates, std::set< const CCopasiObject * > &dependents) const
bool appendDependentMetabolites (std::set< const CCopasiObject * > candidates, std::set< const CCopasiObject * > &dependents) const
bool appendDependentCompartments (std::set< const CCopasiObject * > candidates, std::set< const CCopasiObject * > &dependents) const
bool appendDependentModelValues (std::set< const CCopasiObject * > candidates, std::set< const CCopasiObject * > &dependents) const
bool appendDependentEvents (std::set< const CCopasiObject * > candidates, std::set< const CCopasiObject * > &dependents) const
void removeDependentModelObjects (const std::set< const CCopasiObject * > &deletedObjects)
CCompartmentcreateCompartment (const std::string &name, const C_FLOAT64 &volume=1.0)
bool removeCompartment (const unsigned C_INT32 index, const bool &recursive=true)
bool removeCompartment (const std::string &key, const bool &recursive=true)
bool removeCompartment (const CCompartment *pCompartment, const bool &recursive=true)
CReactioncreateReaction (const std::string &name)
bool removeReaction (const CReaction *pReaction, const bool &recursive=true)
bool removeReaction (const std::string &key, const bool &recursive=true)
bool removeReaction (const unsigned C_INT32 index, const bool &recursive=true)
bool removeLocalReactionParameter (const std::string &key, const bool &recursive=true)
CEventcreateEvent (const std::string &name)
bool removeEvent (const unsigned C_INT32 index, const bool &recursive=true)
bool removeEvent (const std::string &key, const bool &recursive=true)
bool removeEvent (const CEvent *pEvent, const bool &recursive=true)
void synchronizeEventOrder (const CEvent *pEvent, const unsigned C_INT32 newOrder)
CModelValuecreateModelValue (const std::string &name, const C_FLOAT64 &value=0.0)
bool removeModelValue (const CModelValue *pModelValue, const bool &recursive=true)
bool removeModelValue (const std::string &key, const bool &recursive=true)
bool removeModelValue (const unsigned C_INT32 index, const bool &recursive=true)
const CVector< unsigned C_INT32 > & getMetabolitePermutation () const
const CStateTemplategetStateTemplate () const
CStateTemplategetStateTemplate ()
const std::set< const
CCopasiObject * > & 
getUptoDateObjects () const
const std::vector< Refresh * > & getListOfInitialRefreshes () const
const std::vector< Refresh * > & getListOfSimulatedRefreshes () const
const std::vector< Refresh * > & getListOfConstantRefreshes () const
const std::vector< Refresh * > & getListOfNonSimulatedRefreshes () const
bool hasReversibleReaction () const
std::string suitableForStochasticSimulation () const
const bool & isAutonomous () const
std::vector< Refresh * > buildInitialRefreshSequence (std::set< const CCopasiObject * > &changedObjects)
CVector< C_FLOAT64 > initializeAtolVector (const C_FLOAT64 &baseTolerance, const bool &reducedModel) const
std::string printParameterOverview ()
std::string getTimeUnitsDisplayString () const
std::string getFrequencyUnitsDisplayString () const
std::string getVolumeUnitsDisplayString () const
std::string getVolumeRateUnitsDisplayString () const
std::string getConcentrationUnitsDisplayString () const
std::string getConcentrationRateUnitsDisplayString () const
std::string getQuantityRateUnitsDisplayString () const
void evaluateRoots (CVectorCore< C_FLOAT64 > &rootValues, const bool &ignoreDiscrete)
bool processQueue (const C_FLOAT64 &time, const bool &equality, CProcessQueue::resolveSimultaneousAssignments pResolveSimultaneousAssignments)
void processRoots (const C_FLOAT64 &time, const bool &equality, const CVector< C_INT > &roots)
const C_FLOAT64 & getProcessQueueExecutionTime () const
size_t getNumRoots () const
void calculateRootDerivatives (CVector< C_FLOAT64 > &rootDerivatives)
const CVector
< CMathTrigger::CRootFinder * > & 
getRootFinders () const

Static Public Attributes

static const char * VolumeUnitNames []
static const char * AreaUnitNames []
static const char * LengthUnitNames []
static const char * TimeUnitNames []
static const char * QuantityUnitOldXMLNames []
static const char * QuantityUnitNames []
static const char * ModelTypeNames []

Member Enumeration Documentation

Enum of valid area units

Enumerator:
dimensionlessArea 
m2 
dm2 
cm2 
mm2 
microm2 
nm2 
pm2 
fm2 

Enum of valid length units

Enumerator:
dimensionlessLength 
m 
dm 
cm 
mm 
microm 
nm 
pm 
fm 

Enum of valid model types.

Enumerator:
deterministic 
stochastic 

Enum of valid quantity units

Enumerator:
dimensionlessQuantity 
Mol 
mMol 
microMol 
nMol 
pMol 
fMol 
number 
OldXML 

Enum of valid time units

Enumerator:
dimensionlessTime 
d 
h 
min 
s 
ms 
micros 
ns 
ps 
fs 
OldMinute 

Enum of valid volume units

Enumerator:
dimensionlessVolume 
m3 
l 
ml 
microl 
nl 
pl 
fl 

Constructor & Destructor Documentation

CModel::CModel ( CCopasiContainer pParent)

constructor

CModel::~CModel ( )

Destructor


Member Function Documentation

bool CModel::appendDependentCompartments ( std::set< const CCopasiObject * >  candidates,
std::set< const CCopasiObject * > &  dependents 
) const

Appends pointers to compartments which are dependent on the candidates to the list.

Parameters:
std::set<const CCopasiObject * > candidates
std::set<const CCopasiObject * > & dependents
Returns:
bool objectsAppended
bool CModel::appendDependentEvents ( std::set< const CCopasiObject * >  candidates,
std::set< const CCopasiObject * > &  dependents 
) const

Appends a pointers to events which are dependent on the candidates to the list.

Parameters:
std::set<const CCopasiObject * > candidates
std::set<const CCopasiObject * > & dependents
Returns:
bool objectsAppended
bool CModel::appendDependentMetabolites ( std::set< const CCopasiObject * >  candidates,
std::set< const CCopasiObject * > &  dependents 
) const

Appends pointers to metabolites which are dependent on the candidates to the list.

Parameters:
std::set<const CCopasiObject * > candidates
std::set<const CCopasiObject * > & dependents
Returns:
bool objectsAppended
bool CModel::appendDependentModelObjects ( const std::set< const CCopasiObject * > &  candidates,
std::set< const CCopasiObject * > &  dependentReactions,
std::set< const CCopasiObject * > &  dependentMetabolites,
std::set< const CCopasiObject * > &  dependentCompartments,
std::set< const CCopasiObject * > &  dependentModelValues,
std::set< const CCopasiObject * > &  dependentEvents 
) const

Appends pointers to all model objects, which are dependent on the candidates to appropriate lists.

Parameters:
conststd::set< const CCopasiObject * > & candidates
std::set<const CCopasiObject * > & dependentReactions
std::set<const CCopasiObject * > & dependentMetabolites
std::set<const CCopasiObject * > & dependentCompartments
std::set<const CCopasiObject * > & dependentModelValues
std::set<const CCopasiObject * > & dependentEvents
Returns:
bool objectsAppended
bool CModel::appendDependentModelValues ( std::set< const CCopasiObject * >  candidates,
std::set< const CCopasiObject * > &  dependents 
) const

Appends a pointers to model values which are dependent on the candidates to the list.

Parameters:
std::set<const CCopasiObject * > candidates
std::set<const CCopasiObject * > & dependents
Returns:
bool objectsAppended
bool CModel::appendDependentReactions ( std::set< const CCopasiObject * >  candidates,
std::set< const CCopasiObject * > &  dependents 
) const

Appends pointers to reactions which are dependent on the candidates to the list.

Parameters:
std::set<const CCopasiObject * > candidates
std::set<const CCopasiObject * > & dependents
Returns:
bool objectsAppended
void CModel::applyInitialValues ( )

initialize all values of the model with their initial values

std::vector< Refresh * > CModel::buildInitialRefreshSequence ( std::set< const CCopasiObject * > &  changedObjects)

Build the update sequence used to calculate all initial values depending on the changed objects. For metabolites the initial particle number is updated by default unless itself is in the list of changed objects. In that case the initial concentration is updated.

Parameters:
std::set<const CCopasiObject * > & changedObjects
Returns:
std::vector< Refresh * > initialRefreshSequence
void CModel::buildLinkZero ( )

Build the core of the link matrix L

Parameters:
constCMatrix< C_FLOAT64 > & LU
void CModel::buildMoieties ( )

Build the moieties based on the LU decomposition

void CModel::buildRedStoi ( )

Build the core of the link matrix L

Parameters:
constCMatrix< C_FLOAT64 > & LU Build the Reduced Stoichiometry Matrix from the LU decomposition
void CModel::buildStoi ( )

Build the Stoichiometry Matrix from the chemical equations of the steps

void CModel::calculateDerivatives ( C_FLOAT64 *  derivatives)

Calculate the changes of all model quantities determined by ODEs for the model in the current state. The parameter derivatives must at least provide space for state->getVariableNumberSize() double &param C_FLOAT64 * derivatives (output)

void CModel::calculateDerivativesX ( C_FLOAT64 *  derivativesX)

Calculate the changes of all model quantities determined by ODEs for the reduced model in the current state. The parameter derivatives must at least provide space for state->getDependentNumberSize() double &param C_FLOAT64 * derivatives (output)

C_FLOAT64 CModel::calculateDivergence ( ) const

Calculates the divergence for the current state. calculateElasticityMatrix() needs to be called before. It makes only sense to use this method if the Jacobian is not also calculated. In this case it would be more efficient to use the trace of the Jacobian

void CModel::calculateElasticityMatrix ( const C_FLOAT64 &  factor,
const C_FLOAT64 &  resolution 
)

Calculates the elasticity matrix of the model for the current state and stores it in the provided matrix.

Parameters:
CMatrix<C_FLOAT64 > & elasticityMatrix
constC_FLOAT64 & factor,
constC_FLOAT64 & resolution
void CModel::calculateJacobian ( CMatrix< C_FLOAT64 > &  jacobian,
const C_FLOAT64 &  derivationFactor,
const C_FLOAT64 &  resolution 
)

Calculates the jacobian of the full model for the current state and stores it in the provided matrix. calculateElasticityMatrix() needs to be called before.

Parameters:
CMatrix<C_FLOAT64 > & jacobian
constC_FLOAT64 & derivationFactor,
constC_FLOAT64 & resolution
void CModel::calculateJacobianX ( CMatrix< C_FLOAT64 > &  jacobianX,
const C_FLOAT64 &  derivationFactor,
const C_FLOAT64 &  resolution 
)

Calculates the Jacobian of the reduced model for the current state and stores it in the provided matrix. calculateElasticityMatrix() needs to be called before.

Parameters:
constC_FLOAT64 & derivationFactor,
constC_FLOAT64 & resolution
CMatrix<C_FLOAT64 > & jacobianX
void CModel::calculateRootDerivatives ( CVector< C_FLOAT64 > &  rootDerivatives)

Calculate the time derivative of all roots

Parameters:
CVector<C_FLOAT64 > & rootDerivatives
void CModel::cleanup ( )

Cleanup

void CModel::compileDefaultMetabInitialValueDependencies ( )

Compile the default initial value dependencies, which is that the initial concentration is updated.

bool CModel::compileIfNecessary ( CProcessReport pProcessReport)

Compile the model if necessary

Parameters:
CProcessReport*pProcessReport
Returns:
bool success
bool CModel::convert2NonReversible ( )

Converts the set of reactions to a set of reactions where all reactions are irreversible.

CCompartment * CModel::createCompartment ( const std::string &  name,
const C_FLOAT64 &  volume = 1.0 
)

Add a compartment to the model

Parameters:
conststd::string &name
constC_FLOAT64 & volume (default 1.0)
Returns:
bool success (false if failed)
CEvent * CModel::createEvent ( const std::string &  name)

Add a new event to the model

Parameters:
conststd::string &name
Returns:
bool success (false if failed)
CMetab * CModel::createMetabolite ( const std::string &  name,
const std::string &  compartment,
const C_FLOAT64 &  iconc = 1.0,
const CModelEntity::Status status = CModelEntity::REACTIONS 
)

Add a metabolite to the model

Parameters:
conststd::string & name
conststd::string & compartment
constC_FLOAT64 & iconc (default 1.0)
constCMetab::Status & status (default CMetab::METAB_VARIABL)
Returns:
bool success (false if failed)
See also:
CMetab for more information
CModelValue * CModel::createModelValue ( const std::string &  name,
const C_FLOAT64 &  value = 0.0 
)

Add a non concentration value to the model

Parameters:
conststd::string &name
constC_FLOAT64 & value (default 0.0)
CReaction * CModel::createReaction ( const std::string &  name)

Add a new reaction to the model

Parameters:
conststd::string &name
Returns:
bool success (false if failed)
void CModel::evaluateRoots ( CVectorCore< C_FLOAT64 > &  rootValues,
const bool &  ignoreDiscrete 
)

Evaluate all root values for the current state of the model. If ignoreDiscrete is true discrete roots evaluate to 1.0.

Parameters:
CVectorCore<C_FLOAT64 > & rootValues
constbool & ignoreDiscrete
unsigned C_INT32 CModel::findMetabByName ( const std::string &  Target) const

Returns the index of the metab

Parameters:
conststd::string & Target
Returns:
index

Returns the index of the metab

unsigned C_INT32 CModel::findMoiety ( const std::string &  Target) const

Returns the index of the moiety

Parameters:
conststd::string & Target
Returns:
index

Returns the index of the Moiety

bool CModel::forceCompile ( CProcessReport pProcessReport)

Force a compile the model.

Parameters:
CProcessReport*pProcessReport
Returns:
bool success
CModel::AreaUnit CModel::getAreaUnitEnum ( ) const

Get the unit for areas

Returns:
CModel::AreaUnit areaUnit
std::string CModel::getAreaUnitName ( ) const

Get the unit for areas

Returns:
std::string areaUnit
const C_FLOAT64 & CModel::getAvogadro ( ) const

Retrieve the Avogadro number.

std::string CModel::getChildObjectUnits ( const CCopasiObject pObject) const [virtual]

Retrieve the units of the child object.

Returns:
std::string units

Reimplemented from CCopasiContainer.

CCopasiVectorNS< CCompartment > & CModel::getCompartments ( )

Return the compartments of this model

Returns:
CCopasiVectorNS < CCompartment > *
const CCopasiVectorNS< CCompartment > & CModel::getCompartments ( ) const

Return the compartments of this model

Returns:
const CCopasiVectorNS < CCompartment > *
std::string CModel::getConcentrationRateUnitsDisplayString ( ) const

Retrieve the concentration rate units

Returns:
std::string concentrationRateUnits
std::string CModel::getConcentrationUnitsDisplayString ( ) const

Retrieve the concentration units

Returns:
std::string concentrationUnits
CCopasiVectorN< CEvent > & CModel::getEvents ( )

Return the vector of events

const CCopasiVectorN< CEvent > & CModel::getEvents ( ) const

Return the vector of events

std::string CModel::getFrequencyUnitsDisplayString ( ) const

Retrieve the frequency units

Returns:
std::string frequencyUnits
const CState & CModel::getInitialState ( ) const

Get the current state of the model, i.e., all current model quantities.

Returns:
const CState & initialState
const C_FLOAT64 & CModel::getInitialTime ( ) const

Retrieve the initial time

Returns:
const C_FLOAT64 & time
const std::string & CModel::getKey ( ) const [virtual]

Return the key of this model

Returns:
string key

Reimplemented from CModelEntity.

const CModel::CLinkMatrixView & CModel::getL ( ) const

Get the LU decomposition matrix of this model

Returns:
const TNT::Matrix < C_FLOAT64 > & LU Get the link matrix L of the relation: Stoi = L * RedStoi
const CLinkMatrixView L
const CMatrix< C_FLOAT64 > & CModel::getL0 ( ) const

Get the relevant portion of the link matrix L.

Returns:
const CMatrix< C_FLOAT64 > & L0
CModel::LengthUnit CModel::getLengthUnitEnum ( ) const

Get the unit for lengths

Returns:
CModel::LengthUnit lengthUnit
std::string CModel::getLengthUnitName ( ) const

Get the unit for lengths

Returns:
std::string lengthUnit
const std::vector< Refresh * > & CModel::getListOfConstantRefreshes ( ) const

Retrieve the sequence of refresh calls to be executed for updating the constant values

Returns:
const std::vector< Refresh * > & constantRefreshSequence
const std::vector< Refresh * > & CModel::getListOfInitialRefreshes ( ) const

Retrieve the sequence of refresh calls to be executed for updating the initial values

Returns:
const std::vector< Refresh * > & initialRefreshSequence
const std::vector< Refresh * > & CModel::getListOfNonSimulatedRefreshes ( ) const

Retrieve the sequence of refresh calls to be executed for updating the non simulated values

Returns:
const std::vector< Refresh * > & nonsimulatedRefreshSequence
const std::vector< Refresh * > & CModel::getListOfSimulatedRefreshes ( ) const

Retrieve the sequence of refresh calls to be executed for updating the simulated values

Returns:
const std::vector< Refresh * > & simulatedRefreshSequence
const CVector< unsigned C_INT32 > & CModel::getMetabolitePermutation ( ) const

Retrieve the metabolite permutation vector

Returns:
CVector< unsigned C_INT32 > & permutation
const CCopasiVector< CMetab > & CModel::getMetabolites ( ) const

Return the metabolites of this model

Returns:
CCopasiVectorN< CMetab > & metabolites
CCopasiVector< CMetab > & CModel::getMetabolites ( )
const CCopasiVector< CMetab > & CModel::getMetabolitesX ( ) const

Retrieves the vector of metabolites at it is used in the reduced model.

Returns:
const CCopasiVectorN< CMetab > &metabolites
CCopasiVector< CMetab > & CModel::getMetabolitesX ( )
const CModel::ModelType & CModel::getModelType ( ) const

Retrieve the type of the model.

Returns:
const ModelType & modelType
const CCopasiVectorN< CModelValue > & CModel::getModelValues ( ) const

Return the non concentration values of this model

Returns:
CCopasiVectorN< CModelValue > & values
CCopasiVectorN< CModelValue > & CModel::getModelValues ( )
const CCopasiVector< CMoiety > & CModel::getMoieties ( ) const

Return the mMoieties of this model

Returns:
CCopasiVectorN < CMoiety > &
unsigned C_INT32 CModel::getNumAssignmentMetabs ( ) const

Get the number of metabolites determined by assigments

Returns:
unsigned C_INT32 dimension
const C_FLOAT64 & CModel::getNumber2QuantityFactor ( ) const

Get the conversion factor number -> quantity

unsigned C_INT32 CModel::getNumDependentReactionMetabs ( ) const

Get the number of dependent metabolites determined by reactions

Returns:
unsigned C_INT32 dimension
unsigned C_INT32 CModel::getNumIndependentReactionMetabs ( ) const

Get the number of independent metabolites determined by reactions

Returns:
unsigned C_INT32 dimension
unsigned C_INT32 CModel::getNumMetabs ( ) const

Get the number of total metabolites

Returns:
C_INT32 totMetab
unsigned C_INT32 CModel::getNumModelValues ( ) const

Get the number of non concentration values

Returns:
C_INT32
unsigned C_INT32 CModel::getNumODEMetabs ( ) const

Get the number of metabolites determined by ODEs

Returns:
unsigned C_INT32 dimension
size_t CModel::getNumRoots ( ) const

Retrieve the number of roots used in checking for discontinuities.

Returns:
size_t numRoots
unsigned C_INT32 CModel::getNumVariableMetabs ( ) const

Get the number of variable metabolites

Returns:
unsigned C_INT32 totMetab
const CVector< C_FLOAT64 > & CModel::getParticleFlux ( ) const

Retrieve the vector of particle fluxes from the model

Returns:
const CVector< C_FLOAT64 > & particleFlux
const C_FLOAT64 & CModel::getProcessQueueExecutionTime ( ) const

Retrieve the next execution time scheduled in the process queue

Returns:
const C_FLOAT64 & processQueueExecutionTime
const C_FLOAT64 & CModel::getQuantity2NumberFactor ( ) const

Get the conversion factor quantity -> number

std::string CModel::getQuantityRateUnitsDisplayString ( ) const

Retrieve the quantity rate units

Returns:
std::string quantityRateUnits
CModel::QuantityUnit CModel::getQuantityUnitEnum ( ) const

Get the unit for quantities

Returns:
CModel::QuantityUnit quantityUnit
std::string CModel::getQuantityUnitName ( ) const

Get the unit for quantities

Returns:
std::string quantityUnit
std::string CModel::getQuantityUnitOldXMLName ( ) const
CCopasiVectorNS< CReaction > & CModel::getReactions ( )

Return the vector of reactions

Returns:
CCopasiVectorNS <CReaction> & reactions
const CCopasiVectorNS< CReaction > & CModel::getReactions ( ) const

Return the vector of reactions

Returns:
const CCopasiVectorS <CReaction> & reactions
const CMatrix< C_FLOAT64 > & CModel::getRedStoi ( ) const

Get the Reduced Stoichiometry Matrix of this Model

const CVector< CMathTrigger::CRootFinder * > & CModel::getRootFinders ( ) const

Retrieve a vector of root finders

Returns:
const CVector< CMathTrigger::CRootFinder * > & rootFinders
const CState & CModel::getState ( ) const

Get the current state of the model, i.e., all current model quantities.

Returns:
const CState & currentState
const CStateTemplate & CModel::getStateTemplate ( ) const

Retrieve the state template

Returns:
const CModel::CStateTemplate & stateTemplate
CStateTemplate & CModel::getStateTemplate ( )

Retrieve the state template

Returns:
const CModel::CStateTemplate & stateTemplate
const CMatrix< C_FLOAT64 > & CModel::getStoi ( ) const

Get the Stoichiometry Matrix of this Model

const CMatrix< C_FLOAT64 > & CModel::getStoiReordered ( ) const

Get the reordered stoichiometry matrix of this model

const C_FLOAT64 & CModel::getTime ( ) const

Retrieve the actual model time

Returns:
const C_FLOAT64 & time
CModel::TimeUnit CModel::getTimeUnitEnum ( ) const

Get the unit for time

Returns:
CModel::TimeUnit timeUnit
std::string CModel::getTimeUnitName ( ) const

Get the unit for time

Returns:
std::string timeUnit
std::string CModel::getTimeUnitsDisplayString ( ) const

Retrieve the time units

Returns:
std::string timeUnits
unsigned C_INT32 CModel::getTotSteps ( ) const

Get the total steps

Returns:
unsigned C_INT32 total steps;
const std::set< const CCopasiObject * > & CModel::getUptoDateObjects ( ) const

Retrieve the list of objects which are up to date after a call to apply assignment.

Returns:
const std::set< const CCopasiObject * > & uptoDateObjects
std::string CModel::getVolumeRateUnitsDisplayString ( ) const

Retrieve the area units

Returns:
std::string volumeUnits Retrieve the volume units
std::string volumeUnits Retrieve the volume rate units
std::string volumeRateUnits
CModel::VolumeUnit CModel::getVolumeUnitEnum ( ) const

Get the unit for volumes

Returns:
CModel::VolumeUnit volumeUnit
std::string CModel::getVolumeUnitName ( ) const

Get the unit for volumes

Returns:
std::string volumeUnit
std::string CModel::getVolumeUnitsDisplayString ( ) const

Retrieve the volume units

Returns:
std::string volumeUnits
bool CModel::hasReversibleReaction ( ) const

Check whether the model contains reversible reactions

Returns:
bool hasReversibleReaction
CVector< C_FLOAT64 > CModel::initializeAtolVector ( const C_FLOAT64 &  baseTolerance,
const bool &  reducedModel 
) const

Initialize a vector of individual absolute tolerances

Parameters:
constC_FLOAT64 & baseTolerance
constbool & reducedModel
Returns:
CVector< C_FLOAT64 > absoluteTolerances
void CModel::initializeMetabolites ( )

This function must be called to initialize the vector of Metabolites after finishing adding metabolites to compartments.

const bool & CModel::isAutonomous ( ) const

Check whether the model is autonomous

Returns:
const bool &isAutonomous
C_INT32 CModel::load ( CReadConfig configBuffer)

Loads an object with data coming from a CReadConfig object. (CReadConfig object reads an input stream)

Parameters:
pconfigbufferreference to a CReadConfig object.
Returns:
Fail
std::string CModel::printParameterOverview ( )

generates a string that contains a text description of all model parameters (initial values and reaction parameters)

bool CModel::processQueue ( const C_FLOAT64 &  time,
const bool &  equality,
CProcessQueue::resolveSimultaneousAssignments  pResolveSimultaneousAssignments 
)

Process events scheduled at the given which a are checked for equality or not

Parameters:
constC_FLOAT64 & time
constbool & equality
CProcessQueue::resolveSimultaneousAssignmentspResolveSimultaneousAssignments
Returns:
bool stateChanged
void CModel::processRoots ( const C_FLOAT64 &  time,
const bool &  equality,
const CVector< C_INT > &  roots 
)

Check whether the roots which have value 1 lead to firing of events and schedule them if needed.

Parameters:
constC_FLOAT64 & time
constbool & equality
constCVector< C_INT > & roots
bool CModel::removeCompartment ( const unsigned C_INT32  index,
const bool &  recursive = true 
)

Remove a Compartment from the model

bool CModel::removeCompartment ( const std::string &  key,
const bool &  recursive = true 
)

Remove a Compartment from the model

bool CModel::removeCompartment ( const CCompartment pCompartment,
const bool &  recursive = true 
)

Remove a Compartment from the model

void CModel::removeDependentModelObjects ( const std::set< const CCopasiObject * > &  deletedObjects)

Remove all model objects which depend on the deleted objects

Parameters:
conststd::set<const CCopasiObject*> & deletedObjects
bool CModel::removeEvent ( const unsigned C_INT32  index,
const bool &  recursive = true 
)

Remove an event from the model

Parameters:
constunsigned C_INT32 index
constbool & recursive (default: true)
Returns:
bool success
bool CModel::removeEvent ( const std::string &  key,
const bool &  recursive = true 
)

Remove an event from the model

Parameters:
conststd::string & key
constbool & recursive (default: true)
Returns:
bool success
bool CModel::removeEvent ( const CEvent pEvent,
const bool &  recursive = true 
)

Remove an event from the model

Parameters:
constCEvent * pEvent
constbool & recursive (default: true)
Returns:
bool success
bool CModel::removeLocalReactionParameter ( const std::string &  key,
const bool &  recursive = true 
)

Remove a local reaction parameter from the model

bool CModel::removeMetabolite ( const std::string &  key,
const bool &  recursive = true 
)
bool CModel::removeMetabolite ( const unsigned C_INT32  index,
const bool &  recursive = true 
)
bool CModel::removeMetabolite ( const CMetab pMetabolite,
const bool &  recursive = true 
)
bool CModel::removeModelValue ( const CModelValue pModelValue,
const bool &  recursive = true 
)
bool CModel::removeModelValue ( const std::string &  key,
const bool &  recursive = true 
)
bool CModel::removeModelValue ( const unsigned C_INT32  index,
const bool &  recursive = true 
)
bool CModel::removeReaction ( const CReaction pReaction,
const bool &  recursive = true 
)

Remove a reaction from the model using its pointer

bool CModel::removeReaction ( const std::string &  key,
const bool &  recursive = true 
)

Remove a reaction from the model using its key

bool CModel::removeReaction ( const unsigned C_INT32  index,
const bool &  recursive = true 
)

Remove a reaction from the model using its index

bool CModel::setAreaUnit ( const std::string &  name)

Set the unit for areas.

Parameters:
conststd::string & name
Returns:
bool success
bool CModel::setAreaUnit ( const CModel::AreaUnit unit)

Set the unit for areas.

Parameters:
constCModel::AreaUnit & unit
Returns:
bool success
void CModel::setAvogadro ( const C_FLOAT64 &  avogadro)

Set the Avogadro number used for the model.

Parameters:
constC_FLOAT64 & avogadro
void CModel::setCompileFlag ( bool  flag = true)

This must be called whenever something is changed in the model that would make it necessary to recalculate the matrix decomposition

void CModel::setInitialState ( const CState state)

Set all initial model quantities to the one given by the state and updates moieties and metabolite concentrations.

Parameters:
constCState & state
void CModel::setInitialTime ( const C_FLOAT64 &  time)

Set the start time for modeling

Parameters:
constC_FLOAT64 & time
bool CModel::setLengthUnit ( const std::string &  name)

Set the unit for lengths.

Parameters:
conststd::string & name
Returns:
bool success
bool CModel::setLengthUnit ( const CModel::LengthUnit unit)

Set the unit for lengths.

Parameters:
constCModel::LengthUnit & unit
Returns:
bool success
void CModel::setModelType ( const ModelType modelType)

Set the type of the model

Parameters:
constModelType & modelType
bool CModel::setQuantityUnit ( const std::string &  name)

Set the unit for quantities. If COPASI recognizes the unit the conversion factors are set accordingly and true is returned.

Parameters:
conststd::string & name
Returns:
bool success
bool CModel::setQuantityUnit ( const CModel::QuantityUnit unit)

Set the unit for quantities. If COPASI recognizes the unit the conversion factors are set accordingly and true is returned.

Parameters:
constCModel::QuantityUnit & unit
Returns:
bool success
void CModel::setState ( const CState state)

Set all independent current model quantities to the one given by the state.

Parameters:
constCState & state
void CModel::setTime ( const C_FLOAT64 &  time)

Set the actual model time

Parameters:
constC_FLOAT64 & time
bool CModel::setTimeUnit ( const std::string &  name)

Set the unit for time. If COPASI recognizes the unit the conversion factors are set accordingly and true is returned.

Parameters:
conststd::string & name
Returns:
bool success
bool CModel::setTimeUnit ( const CModel::TimeUnit unit)

Set the unit for time. If COPASI recognizes the unit the conversion factors are set accordingly and true is returned.

Parameters:
constconst CModel::TimeUnit & unit
Returns:
bool success
bool CModel::setTitle ( const std::string &  title)

Return a pointer to the current time Set the title of this model

Parameters:
const string &title title for this model
bool CModel::setVolumeUnit ( const std::string &  name)

Set the unit for volumes. If COPASI recognizes the unit the conversion factors are set accordingly and true is returned.

Parameters:
conststd::string & name
Returns:
bool success
bool CModel::setVolumeUnit ( const CModel::VolumeUnit unit)

Set the unit for volumes. If COPASI recognizes the unit the conversion factors are set accordingly and true is returned.

Parameters:
constCModel::VolumeUnit & unit
Returns:
bool success
std::string CModel::suitableForStochasticSimulation ( ) const

check if the model is suitable for stochastic simulation

void CModel::synchronizeEventOrder ( const CEvent pEvent,
const unsigned C_INT32  newOrder 
)

Synchronize the order of other events effected by the change of the given event

Parameters:
constCEvent * pEvent
constunsigned C_INT32 newOrder
bool CModel::updateInitialValues ( )

Update all initial values.

Returns:
bool success
void CModel::updateNonSimulatedValues ( void  )

Calling this method after updateSimulatedValues assure that all model values even those not needed for simulation are consistent with the current state

void CModel::updateSimulatedValues ( const bool &  updateMoieties)

This method calculates all values needed for simulation based on the current current state. If updateMoities is true the particle numbers of dependent metabolites of type REACTION are calculated otherwise they are assumed to be synchronized.

Parameters:
constbool & updateMoieties

Member Data Documentation

const char * CModel::AreaUnitNames [static]
Initial value:
  {"dimensionless", "m\xc2\xb2", "dm\xc2\xb2", "cm\xc2\xb2", "mm\xc2\xb2", "\xc2\xb5m\xc2\xb2", "nm\xc2\xb2", "pm\xc2\xb2", "fm\xc2\xb2", NULL}

String representation of valid area units

const char * CModel::LengthUnitNames [static]
Initial value:
  {"dimensionless", "m", "dm", "cm", "mm", "\xc2\xb5m", "nm", "pm", "fm", NULL}

String representation of valid length units

const char * CModel::ModelTypeNames [static]
Initial value:
  {"deterministic", "stochastic", NULL}

String representation of the valid model types.

const char * CModel::QuantityUnitNames [static]
Initial value:
  {"dimensionless", "mol", "mmol", "\xc2\xb5mol", "nmol", "pmol", "fmol", "#", NULL}

String representation of valid quantity units

const char * CModel::QuantityUnitOldXMLNames [static]
Initial value:
  {"dimensionless", "Mol", "mMol", "\xc2\xb5Mol", "nMol", "pMol", "fMol", "#", NULL}

String representation of valid quantity units as used in old (up to Build 18) COPASI files

const char * CModel::TimeUnitNames [static]
Initial value:
  {"dimensionless", "d", "h", "min", "s", "ms", "\xc2\xb5s", "ns", "ps", "fs", NULL}

String representation of valid time units

const char * CModel::VolumeUnitNames [static]
Initial value:
  {"dimensionless", "m\xc2\xb3", "l", "ml", "\xc2\xb5l", "nl", "pl", "fl", NULL}

String representation of valid volume units


The documentation for this class was generated from the following files:
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines