rrRoadRunner

class RoadRunner

The main RoadRunner class.

The RoadRunner class is responsible for loading and simulating SBML models.

MemoryManagment: Any pointer returned by a get… method is owned by the RoadRunner object and does NOT have to be deleted.

Public Functions

RoadRunner(unsigned int level = 3, unsigned int version = 2)

create an empty RoadRunner object with given SBML level and version.

RoadRunner(const std::string &uriOrSBML, const Dictionary *options = 0)

load an sbml document from anywhere.

same arguments as load.

If options is not null, then the RoadRunner::computeAndAssignConservationLaws flag is set to whatever value is specified in the options struct.

Parameters:
  • uriOrSBML – a URI, local path or sbml document contents.

  • options – an options struct, if null, default values are used.

RoadRunner(const std::string &compiler, const std::string &tempDir, const std::string &supportCodeDir)

All three of the RoadRunner options default to the empty std::string, in this case, the default values are used.

Parameters:
  • compiler – If LLVM build is enabled, the compiler defaults to LLVM.

  • tempDir – If the old external C compiler is used, this is the where the C files are written to.

  • supportCodeDir – If the old external C compiler is used, this is the location where roadrunner C include files are.

RoadRunner(const RoadRunner &rr)

RoadRunner copy constructor.

void operator=(const RoadRunner &rr)

RoadRunner copy constructor.

virtual ~RoadRunner()

free any memory this class allocated

int getInstanceID()

When there are multiple instances of RoadRunner, this is the instance id.

int getInstanceCount()

Number of currently running RoadRunner instances.

std::string getInfo()

information about the current state of this object.

double getCurrentTime()

return the current time as a double

class Compiler *getCompiler()

The Compiler that the ModelGenerator is using to compile / interpret sbml code.

void setCompiler(const std::string &compiler)

Set the name of the externa compiler to use.

Some ModelGenerators may have no use for this value.

An exception is raised if the std::string is invalid.

Integrator *getIntegrator()

get a pointer to the integrator which is currently being used to time evolve the system.

SteadyStateSolver *getSteadyStateSolver()

get a pointer to the current steady state solver

SensitivitySolver *getSensitivitySolver()

get a pointer to the current sensitivities solver

Integrator *getIntegratorByName(const std::string &name)

Get an integrator by name without switching the current integrator.

SteadyStateSolver *getSteadyStateSolverByName(const std::string &name)

Get an SteadyStateSolver by name without switching the current solver.

SensitivitySolver *getSensitivitySolverByName(const std::string &name)

Get an SensitivitySolver by name without switching the current solver.

Integrator *makeIntegrator(const std::string &name)

Make an integrator for the given name.

SteadyStateSolver *makeSteadyStateSolver(const std::string &name)

Make an SteadyStateSolver* for the given name.

SensitivitySolver *makeSensitivitySolver(const std::string &name)

Make SensitivitySolver for the given name.

std::vector<std::string> getExistingIntegratorNames()

Return a list of the names of all existing integrators.

std::vector<std::string> getExistingSteadyStateSolverNames()

Return a list of the names of all existing integrators.

std::vector<std::string> getExistingSensitivitySolverNames()

Return a list of the names of all existing sensitivity solvers.

void setIntegrator(const std::string &name)

set the current Integrator to

Parameters:

name;

void setSteadyStateSolver(const std::string &name)

set the current SteadyStateSolver to

use getRegisteredSteadyStateSolverNames to get a list of available names

Parameters:

name;

void setSensitivitySolver(const std::string &name)

set the current SensitivitySolver to

use getRegisteredSensitivitySolverNames to get a list of available sensitivity solver names

Parameters:

name;

bool integratorExists(const std::string &name)

returns true if integrator

Use getRegisteredIntegratorNames to get a list of available integrator solver names.

Parameters:

name – exists.

bool steadyStateSolverExists(const std::string &name)

returns true if SteadyStateSolver

Parameters:

name – exists.

bool sensitivitySolverExists(const std::string &name)

returns true if SensitivitySolver

Parameters:

name – exists.

bool isModelLoaded()
std::string getModelName()

returns the model name if a model is loaded, empty std::string otherwise.

void setModelName(const std::string &name)

sets the model name if a model is loaded.

std::string getModelId()

returns the model id if a model is loaded, empty std::string otherwise.

void setModelId(const std::string &id)

sets the model id if a model is loaded.

bool clearModel()

Clears the currently loaded model and all associated memory.

Deletes jitted code and libStruct data

Returns:

True if memory was freed, false if no model was loaded in the first place

double oneStep(double currentTime, double stepSize, bool reset = true)

Carry out a single integration step using a stepsize as indicated in the method call.

Parameters:
  • currentTime – current time

  • stepSize – the step size

  • reset – whether to reset the model (defaults to true).

double internalOneStep(double currentTime, double stepSize, bool reset = true)

Carry out a single internal solver integration step.

Parameters:
  • currentTime – current time

  • stepSize – the step size

  • reset – whether to reset the model (defaults to true).

const ls::DoubleMatrix *simulate(const SimulateOptions *options = 0)

Simulate the current SBML model.

If options is null, then the current simulation settings (start time, end time, n steps) are used. If options is not null, then the current simulation settings are set to the values specified by options and they are used.

The options Dictionary

may contain a large number of options, for a complete list of all available options for each integrator type,

For example, to perform a simulation from time 0 to 10 with 1000 steps, using a stiff integtator, you would:

RoadRunner rr = RoadRunner("someFile.xml");
SimulateOptions opt = rr.getSimulateOptions();
opt.start = 0;
opt.duration = 10;
opt.steps = 1000;
const DoubleMatrix *result = rr.simulate(&opt);

See also

IntegratorFactory,

See also

SimulateOptions.

Similarly, options specific to a particular integrator, such as the ‘seed’ option with the Gillespie integrator, this is set via the ‘setIntegrator’ “integrator” key, i.e.

RoadRunner rr = RoadRunner("someFile.xml");
rr.setIntegrator("gillespie");
SimulateOptions opt;
opt.start = 0;
opt.duration = 10;
opt.steps = 1000;
opt.setItem("stiff", true);
opt.setItem("seed", 12345);
const DoubleMatrix *result = rr.simulate(&opt);

Here, the “integrator” specifies the integrator to use. The “stiff” key is only used by the deterministic solvers, and it is safely ignored by the stochastic solvers. Also, the “seed” sets the random seed that the integrator uses. For more information about all of the available options for each integrator,

If one wants to not store the result matrix in memory and instead write it to a file during simulation, one can set the output_file option. When output file is nonempty,

simulate() will write to its path once every Config::K_ROWS_PER_WRITE rows are generated and clear the matrix. So an empty result matrix is returned, and the last simulation results are not stored.

See also

IntegratorFactory::getIntegratorOptions”.

Throws:

an – std::exception if any options are invalid.

Returns:

a borrowed reference to a DoubleMatrix object if successful. The matrix will be empty if output_file is specified and nonempty.

const ls::DoubleMatrix *simulate(double start, double stop, int points)

simulate the model using currently set integrator

Parameters:
  • start – starting time to simulate

  • stop – what time point does the simulation end?

  • points – how many points to output (one greater than the number of steps to take).

const ls::DoubleMatrix *simulate(const std::vector<double> &times)

simulate the model using currently set integrator

Parameters:

times – a vector of all the time outputs desired.

Matrix3D<double, double> timeSeriesSensitivities(double start, double stop, int num, std::vector<std::string> params = std::vector<std::string>(), std::vector<std::string> species = std::vector<std::string>(), int k = 0)

simulate a timeseries with sensitivities from start to step with num data points.

Matrix3D indexed by time. Each element of the 3D matrix is a Matrix<double> with rows and columns parameters and model variables respectively. The parameter k determines the kth order derivative of the sensitivity information that will be returned

Parameters:
  • start – starting time for time series simulation

  • stop – last time point for time series simulation

  • num – number of data points to simulate. Determines Z of Matrix3D.

  • params – vector of parameters that you want sensitivity for. When empty (default), compute sensitivities for all parameters vs all variables.

  • species – vector of species to include in the results Default is empty, in which case all species will be included. All species are selected during solving and slicing only occurs at the end.

  • k – (default 0) return the kth other derivative of the sensitivity data.

void saveState(std::string filename, char opt = 'b')

similar to saveStateS but save data to file

Parameters:
  • filename – the filename to save the data to.

  • opt – either ‘b’ (default) or ‘r’. The latter is used for debugging.

std::stringstream *saveStateS(char opt = 'b')

save state as binary to a stringstream so it can be loaded again later.

Parameters:

opt – either ‘b’ (default) or ‘r’. The latter is used for debugging.

Returns:

new reference to a stringstream - the caller is responsible for reclaiming memory Importantly, if the stream is read into a new RoadRunner instance via RoadRunner::loadStateS then the user no longer needs to handle the reference manually because RoadRunner::loadStateS does it after loading the state.

void loadState(const std::string &filename)

Loads a roadrunner instance saved by saveState with the ‘b’ option.

void loadStateS(std::stringstream *state)

load state from a

The stingstream pointer should be heap allocated and generated from RoadRunner::saveStateS. It is an error to use a stringstream that was generated in any other way. The memory associated with the stringstream is automatically cleaned up after loading the state.

Parameters:

state – stringstream that was produced by RoadRunner::saveStateS.

const ls::DoubleMatrix *getSimulationData() const

RoadRunner keeps a copy of the simulation data around until the next call to simulate.

This matrix can be obtained here.

void setSimulateOptions(const SimulateOptions &settings)
SimulateOptions &getSimulateOptions()

get a reference to the SimulateOptions that were set either by setSimulateOptions or simulate.

RoadRunnerOptions &getOptions()

Get a reference to the options that determine how this class should behave.

These are general options. For options specific for loading or simulationg,

See also

getSimulateOptions.

void setOptions(const RoadRunnerOptions&)
std::string getSBML(int level = 0, int version = 0)

get the originally loaded sbml document as a std::string.

This may optionally up or downconvert the document to a different version, if the level and version arguments are non-zero.

If both arguments are zero, then the document is left alone and the

std::string getCurrentSBML(int level = 0, int version = 0)

Returns the SBML with the current model parameters.

This is different than getSBML which returns the orginal SBML.

This may optionally up or downconvert the document to a different version, if the level and version arguments are non-zero.

If both arguments are zero, then the document is left alone and the original version is returned.

void reset()

Picks up default options from config.

The current default is SelectionRecord::TIME | SelectionRecord::RATE | SelectionRecord::FLOATING.

void reset(int options)

reset the model accordign to a bitfield specified by the SelectionRecord::SelectionType values.

Note, this would make more sense as an unsigned, however SWIG has issues mapping an enum to an unsigned, but seems to std::map enums just fine to an int.

For example, to reset the floating species, time and rate rule values:

rr.reset(SelectionRecord::TIME | SelectionRecord::RATE | SelectionRecord::FLOATING);

Parameters:

options – a bitmask made from the SelectionRecord::SelectionTypes values.

void resetSelectionLists()

resets time course and steady state selection lists to defaults

void changeInitialConditions(const std::vector<double> &ic)

set the floating species initial concentrations.

equivalent to ExecutableModel::reset, then ExecutableModel::setFloatingSpeciesConcentrations

ExecutableModel *getModel()

get a pointer to the ExecutableModel owned by the RoadRunner object.

void load(const std::string &uriOrSBML, const Dictionary *options = 0)

load an sbml document from anywhere.

If options is not null, then the RoadRunner::computeAndAssignConservationLaws flag is set to whatever value is specified in the options struct.

WARNING, will no longer return a value in next version.

Parameters:
  • uriOrSBML – a URI, local path or sbml document contents.

  • options – an options struct, if null, default values are used.

rr::SelectionRecord createSelection(const std::string &str)

create a selection record.

This record can be used to select values.

std::vector<rr::SelectionRecord> &getSelections()

Returns the currently selected columns that will be returned by calls to simulate() or simulateEx(,,).

double getValue(const std::string &sel)

Creates a new selection based on the selection std::string, and returns the value it queries.

double getValue(const SelectionRecord &record)
void setSelections(const std::vector<std::string> &selections)
void setSelections(const std::vector<rr::SelectionRecord> &selections)
std::vector<double> getSelectedValues()

returns the values selected with SimulateOptions for the current model time / timestep”)

void getIds(int types, std::list<std::string> &ids)

populates a given list with all the ids that this class can accept.

std::vector<std::string> getIndependentFloatingSpeciesIds()

Gets the ids for all independent floating species.

Independent means (in this case) that the species is unique up to a conserved cycle, i.e. each conserved cycle counts for only one species. See Sauro, Systems Biology: Introduction to Pathway Modeling, 1st ed. pp. 60.

Author

JKM

std::vector<std::string> getDependentFloatingSpeciesIds()

Gets the ids for all dependent floating species.

See getIndependentFloatingSpeciesIds for an explanation of independent vs. dependent.

Author

JKM

std::vector<std::string> getFloatingSpeciesConcentrationIds()

Gets the ids for all floating species concentrations.

Author

KC

std::vector<std::string> getFloatingSpeciesInitialConcentrationIds()

Gets the ids for all initial floating species concentrations.

See getIndependentFloatingSpeciesIds for an explanation of independent vs. dependent.

int getSupportedIdTypes()

returns a bit field of the ids that this class supports.

void setValue(const std::string &id, double value)

sets the value coresponding to the given selection std::string

raises an exception in the selection std::string is invalid.

ls::DoubleMatrix getFloatingSpeciesAmountsNamedArray()

Returns the floating species amounts as a named array.

Author

KC

ls::DoubleMatrix getFloatingSpeciesConcentrationsNamedArray()

Returns the floating species concentrations as a named array.

Author

KC

ls::DoubleMatrix getBoundarySpeciesAmountsNamedArray()

Returns the boundary species amounts as a named array.

Author

KC

ls::DoubleMatrix getBoundarySpeciesConcentrationsNamedArray()

Returns the boundary species concentrations as a named array.

Author

KC

std::vector<double> getIndependentFloatingSpeciesAmountsV()

Returns the independent floating species amounts as an array.

Author

Claire

std::vector<double> getDependentFloatingSpeciesAmountsV()

Returns the dependent floating species amounts as an array.

Author

Claire

std::vector<double> getIndependentFloatingSpeciesConcentrationsV()

Returns the independent floating species concentrations as an array.

Author

Claire

std::vector<double> getDependentFloatingSpeciesConcentrationsV()

Returns the dependent floating species concentrations as an array.

Author

Claire

ls::DoubleMatrix getIndependentFloatingSpeciesAmountsNamedArray()

Returns the independent floating species amounts as a named array.

Author

Claire

ls::DoubleMatrix getDependentFloatingSpeciesAmountsNamedArray()

Returns the dependent floating species amounts as a named array.

Author

Claire

ls::DoubleMatrix getIndependentFloatingSpeciesConcentrationsNamedArray()

Returns the independent floating concentrations amounts as a named array.

Author

Claire

ls::DoubleMatrix getDependentFloatingSpeciesConcentrationsNamedArray()

Returns the dependent floating concentrations amounts as a named array.

Author

Claire

std::vector<double> getRatesOfChange()

Returns the rate of change of the floating species as an array.

Author

KC

ls::DoubleMatrix getRatesOfChangeNamedArray()

Returns the rate of change of the floating species as a named array.

Author

KC

std::vector<double> getIndependentRatesOfChange()

Returns the rate of change of the independent floating species as an array.

Author

KC

ls::DoubleMatrix getIndependentRatesOfChangeNamedArray()

Returns the rate of change of the independent floating species as a named array.

Author

KC

std::vector<double> getDependentRatesOfChange()

Returns the rate of change of the dependent floating species as an array.

Author

KC

ls::DoubleMatrix getDependentRatesOfChangeNamedArray()

Returns the rate of change of the dependent floating species as a named array.

Author

KC

ls::DoubleMatrix getFullJacobian()

compute the full Jacobian at the current operating point

ls::DoubleMatrix getFullReorderedJacobian()
ls::DoubleMatrix getReducedJacobian(double h = -1.0)

Compute the reduced Jacobian at the current operating point.

Parameters:

h – The step sized used for central difference method. If negative, the default value from the config file is used.

std::vector<std::complex<double>> getFullEigenValues()

Returns the eigenvalues of the full jacobian.

This returns a vector of Complex numbers. In the python bindings, the values are complex only if there are non-zero imaginary parts of the values.

std::vector<std::complex<double>> getReducedEigenValues()

Returns the eigenvalues of the reduced jacobian.

This returns a vector of Complex numbers. In the python bindings, the values are complex only if there are non-zero imaginary parts of the values.

ls::DoubleMatrix getFullEigenValuesNamedArray()

Returns the eigenvalues of the full jacobian as a named array.

This returns an N x 2 matrix where the first column is the real values and the second is the imaginary part. The rows are labeled with the corresponding species ids, and the columns are labeled ‘real’ and ‘imaginary’.

ls::DoubleMatrix getReducedEigenValuesNamedArray()

Returns the eigenvalues of the reduced jacobian as a named array.

This returns an N x 2 matrix where the first column is the real values and the second is the imaginary part. The rows are labeled with the corresponding species ids, and the columns are labeled ‘real’ and ‘imaginary’.

ls::DoubleMatrix getLinkMatrix()
ls::DoubleMatrix getNrMatrix()

get the reduced stochiometry matrix.

If conservation conversion is enabled, this is the matrix that coresponds to the independent species.

A synonym for getReducedStoichiometryMatrix().

ls::DoubleMatrix getKMatrix()

Returns the K matrix (right nullspace of Nr) The K matrix has the structure, [I K0]’.

ls::DoubleMatrix getReducedStoichiometryMatrix()

get the reduced stochiometry matrix.

If conservation conversion is enabled, this is the matrix that coresponds to the independent species.

A synonym for getNrMatrix().

ls::DoubleMatrix getFullStoichiometryMatrix()

Get the stoichiometry matrix that coresponds to the full model, even it it was converted via conservation conversion.

ls::DoubleMatrix getExtendedStoichiometryMatrix()
ls::DoubleMatrix getL0Matrix()
ls::DoubleMatrix getConservationMatrix()
ls::DoubleMatrix getUnscaledConcentrationControlCoefficientMatrix()
ls::DoubleMatrix getScaledConcentrationControlCoefficientMatrix()
ls::DoubleMatrix getUnscaledFluxControlCoefficientMatrix()
ls::DoubleMatrix getScaledFluxControlCoefficientMatrix()
std::vector<std::string> getEigenValueIds()

returns the list of floating species, but with a “eigen(…)” std::string wrapped around them.

double getUnscaledParameterElasticity(const std::string &reactionName, const std::string &parameterName)

Returns the unscaled elasticity for a named reaction with respect to a named parameter.

ls::DoubleMatrix getFrequencyResponse(double startFrequency, int numberOfDecades, int numberOfPoints, const std::string &parameterName, const std::string &variableName, bool useDB, bool useHz)
void setConservedMoietyAnalysis(bool value)

This method turns on / off the computation and adherence to conservation laws.

bool getConservedMoietyAnalysis()

is conservation analysis enabled.

This is set

double getDiffStepSize() const

Set the differential step size used in MCA routines like getCC.

Author

JKM

void setDiffStepSize(double val)

Set the differential step size used in MCA routines like getCC.

Author

JKM

double getSteadyStateThreshold() const

Get the steady state threshold used in getCC.

In the MCA routines, RoadRunner will keep trying to converge to a steady state until this threshold is met

Author

JKM

void setSteadyStateThreshold(double val)

Set the steady state threshold used in getCC.

In the MCA routines, RoadRunner will keep trying to converge to a steady state until this threshold is met

Author

JKM

double getFluxThreshold() const

Get the flux threshold used in getCC.

In the MCA routines, RoadRunner will keep trying to converge to a steady state until this threshold is met

void setFluxThreshold(double val)

Set the flux threshold used in getCC.

In getCC, RoadRunner will return 0.0 for flux control coefficients if the

double getuCC(const std::string &variableName, const std::string &parameterName)

Get unscaled control coefficient with respect to a global parameter.

variableName must be either a reaction or floating species.

parameterName must be eithe a global parameter, boundary species, or conserved sum.

double getCC(const std::string &variableName, const std::string &parameterName)

Get scaled control coefficient with respect to a global parameter.

The variableName must be either a reaction id, or a floating species id.

The parameterName must be either a global parameter, boundary species, or conserved sum.

double getuEE(const std::string &reactionName, const std::string &parameterName)

Get unscaled elasticity coefficient with respect to a global parameter or species.

double getuEE(const std::string &reactionName, const std::string &parameterName, bool computeSteadystate)

Get unscaled elasticity coefficient with respect to a global parameter or species.

Optionally the model is brought to steady state after the computation.

double getEE(const std::string &reactionName, const std::string &parameterName)

Get scaled elasticity coefficient with respect to a global parameter or species.

double getEE(const std::string &reactionName, const std::string &parameterName, bool computeSteadyState)

Get scaled elasticity coefficient with respect to a global parameter or species.

Optionally the model is brought to steady state after the computation.

ls::DoubleMatrix getUnscaledElasticityMatrix()

Compute the unscaled species elasticity matrix at the current operating point.

ls::DoubleMatrix getScaledElasticityMatrix()

Compute the unscaled elasticity matrix at the current operating point.

double getScaledFloatingSpeciesElasticity(const std::string &reactionName, const std::string &speciesName)

Compute the scaled elasticity for a given reaction and given species.

double getUnscaledSpeciesElasticity(int reactionId, int speciesIndex)

Get a single species elasticity value IMPORTANT: Assumes that the reaction rates have been precomputed at the operating point !!

void addSpeciesConcentration(const std::string &sid, const std::string &compartment, double initConcentration, bool hasOnlySubstanceUnits = false, bool boundaryCondition = false, const std::string &substanceUnits = "", bool forceRegenerate = true)

Add a species to the current model and set its concentration.

Parameters:
  • sid – the ID of the species to be added

  • compartment – the compartment of the species to be added

  • initConcentration – the initial concentration of the species to be added

  • hasOnlySubstanceUnits – the boolean value of the hasOnlySubstanceUnits attribute, indicating whether the species is to be used as an amount when encountered in formulas elsewhere in the model.

  • boundaryCondition – the boolean value of the boundaryCondition attribute, indicating whether the species is a boundary species.

  • substanceUnits – the substance unit of the species to be added

  • forceRegenerate – a boolean value to indicate if the model is regenerated after this function call default value is true to regenerate model after each call of editing function to save time for editing for multiple times, one could set this flag to true only in the last call of editing

void addSpeciesAmount(const std::string &sid, const std::string &compartment, double initAmount = 0, bool hasOnlySubstanceUnits = false, bool boundaryCondition = false, const std::string &substanceUnits = "", bool forceRegenerate = true)

Add a species to the current model.

Parameters:
  • sid – the ID of the species to be added

  • compartment – the compartment of the species to be added

  • initAmount – the initial amount of the species to be added

  • hasOnlySubstanceUnits – the boolean value of the hasOnlySubstanceUnits attribute, indicating whether the species is to be used as an amount when encountered in formulas elsewhere in the model.

  • boundaryCondition – the boolean value of the boundaryCondition attribute, indicating whether the species is a boundary species.

  • substanceUnits – the substance unit of the species to be added

  • forceRegenerate – a boolean value to indicate if the model is regenerated after this function call default value is true to regenerate model after each call of editing function to save time for editing for multiple times, one could set this flag to true only in the last call of editing

void removeSpecies(const std::string &sid, bool forceRegenerate = true)
void setBoundary(const std::string &sid, bool boundaryCondition, bool forceRegenerate = true)

Set the boundary condition of an existing species.

By default, the boundary condition for a species is false.

Parameters:
  • sid – the ID of the species

  • boundaryCondition – boundary condition to be set

  • forceRegenerate – a boolean value to indicate if the model is regenerated after this function call default value is true to regenerate model after each call of editing function to save time for editing for multiple times, one could set this flag to true only in the last call of editing

void setHasOnlySubstanceUnits(const std::string &sid, bool hasOnlySubstanceUnits, bool forceRegenerate = true)

Set the hasOnlySubstanceUnits attribute for an existing species.

Parameters:
  • sid – the ID of a species

  • hasOnlySubstanceUnits – the value of hasOnlySubstanceUnits attribute to be set

  • forceRegenerate – a boolean value to indicate if the model is regenerated after this function call default value is true to regenerate model after each call of editing function to save time for editing for multiple times, one could set this flag to true only in the last call of editing

bool getHasOnlySubstanceUnits(const std::string &sid)

Get the hasOnlySubstanceUnits attribute for an existing species.

Parameters:

sid – the ID of a species

void setInitAmount(const std::string &sid, double initAmount, bool forceRegenerate = true)

Set initial amount for an existing species.

Initial amount/concentration set before will be unset.

Parameters:
  • sid – the ID of a species

  • initAmount – the initial amount to be set

  • forceRegenerate – a boolean value to indicate if the model is regenerated after this function call default value is true to regenerate model after each call of editing function to save time for editing for multiple times, one could set this flag to true only in the last call of editing

void setInitConcentration(const std::string &sid, double initConcentration, bool forceRegenerate = true)

Set initial concentration for an existing species.

Initial amount/concentration set before will be unset.

Parameters:
  • sid – the ID of a species

  • initConcentration – the initial concentration to be set

  • forceRegenerate – a boolean value to indicate if the model is regenerated after this function call default value is true to regenerate model after each call of editing function to save time for editing for multiple times, one could set this flag to true only in the last call of editing

void setConstant(const std::string &sid, bool constant, bool forceRegenerate = true)

Set the constant attribute for an existing species/ parameter/ compartment By default, the constant attribute is false.

Parameters:
  • sid – the ID of a species/ parameter/ compartment

  • constant – the constant attribute to be set

  • forceRegenerate – a boolean value to indicate if the model is regenerated after this function call default value is true to regenerate model after each call of editing function to save time for editing for multiple times, one could set this flag to true only in the last call of editing

void addReaction(const std::string &sbmlRep, bool forceRegenerate = true)
void addReaction(const std::string &rid, std::vector<std::string> reactants, std::vector<std::string> products, const std::string &kineticLaw, bool forceRegenerate = true)
void removeReaction(const std::string &rid, bool deleteUnusedParameters = false, bool forceRegenerate = true)

Remove a reaction from the current model.

Parameters:
  • rid – the ID of the reaction to be removed

  • deleteUnusedParameters; – a boolean value to indicate whether to remove unused parameters from the model as well.

  • forceRegenerate – a boolean value to indicate if the model is regenerated after this function call default value is true to regenerate model after each call of editing function to save time for editing for multiple times, one could set this flag to true only in the last call of editing

void setReversible(const std::string &rid, bool reversible, bool forceRegenerate = true)
void setKineticLaw(const std::string &rid, const std::string &kineticLaw, bool forceRegenerate = true)
std::string getKineticLaw(const std::string &rid)

Get the kinetic law of an existing reaction in the current model.

Parameters:

rid – the ID of reaction to be modified

void addParameter(const std::string &pid, double value, bool forceRegenerate = true)

Add a parameter to the current model.

Parameters:
  • pid – the ID of the parameter to be added

  • value – the value of the parameter to be added

  • forceRegenerate – a boolean value to indicate if the model is regenerated after this function call default value is true to regenerate model after each call of editing function to save time for editing for multiple times, one could set this flag to true only in the last call of editing

void removeParameter(const std::string &pid, bool forceRegenerate = true)

Remove a parameter from the current model.

Parameters:
  • pid – the ID of the parameter to be removed

  • forceRegenerate – a boolean value to indicate if the model is regenerated after this function call default value is true to regenerate model after each call of editing function to save time for editing for multiple times, one could set this flag to true only in the last call of editing

void addCompartment(const std::string &cid, double initVolume, bool forceRegenerate = true)

Add a compartment to the current model.

Parameters:
  • cid – the ID of the compartment to be added

  • initVolume – the initial volume of the compartment to be added

  • forceRegenerate – a boolean value to indicate if the model is regenerated after this function call default value is true to regenerate model after each call of editing function to save time for editing for multiple times, one could set this flag to true only in the last call of editing

void removeCompartment(const std::string &cid, bool forceRegenerate = true)

Remove a compartment from the current model.

Parameters:
  • cid – the ID of the compartment to be removed

  • forceRegenerate – a boolean value to indicate if the model is regenerated after this function call default value is true to regenerate model after each call of editing function to save time for editing for multiple times, one could set this flag to true only in the last call of editing

void addAssignmentRule(const std::string &vid, const std::string &formula, bool forceRegenerate = true)
void addRateRule(const std::string &vid, const std::string &formula, bool forceRegenerate = true)
void removeRules(const std::string &vid, bool useInitialValue = false, bool forceRegenerate = true)

Remove rules related to given variable from the current model.

Parameters:
  • vid – ID of variable that rules assign formula to

  • useInitialValue – a boolean value to indicate whether using initial value as current value when the assignment rule of a variable is removed otherwise, it will keep the value of assignment rule as current value

  • forceRegenerate – a boolean value to indicate if the model is regenerated after this function call default value is true to regenerate model after each call of editing function to save time for editing for multiple times, one could set this flag to true only in the last call of editing

void addInitialAssignment(const std::string &vid, const std::string &formula, bool forceRegenerate = true)
void removeInitialAssignment(const std::string &vid, bool forceRegenerate = true, bool errIfNotExist = true)

Remove initial assignment for a symbol from the current model.

Parameters:
  • vid – ID of the symbol

  • forceRegenerate – a boolean value to indicate if the model is regenerated after this function call default value is true to regenerate model after each call of editing function to save time for editing for multiple times, one could set this flag to true only in the last call of editing

void addEvent(const std::string &eid, bool useValuesFromTriggerTime, const std::string &trigger, bool forceRegenerate = true)
void addTrigger(const std::string &eid, const std::string &trigger, bool forceRegenerate = true)
void setPersistent(const std::string &eid, bool persistent, bool forceRegenerate = true)
void setTriggerInitialValue(const std::string &eid, bool initValue, bool forceRegenerate = true)
void addPriority(const std::string &eid, const std::string &priority, bool forceRegenerate = true)
void addDelay(const std::string &eid, const std::string &delay, bool forceRegenerate = true)
void addEventAssignment(const std::string &eid, const std::string &vid, const std::string &formula, bool forceRegenerate = true)
void removeEventAssignments(const std::string &eid, const std::string &vid, bool forceRegenerate = true)

Remove event assignments for given variable from an existing event.

Parameters:
  • eid – the ID of the event

  • vid – the ID of the vairable of the event assignments

  • forceRegenerate – a boolean value to indicate if the model is regenerated after this function call default value is true to regenerate model after each call of editing function to save time for editing for multiple times, one could set this flag to true only in the last call of editing

void removeEvent(const std::string &eid, bool forceRegenerate = true)

Remove an event from the current model.

Parameters:
  • eid – the ID of the event to be removed

  • forceRegenerate – a boolean value to indicate if the model is regenerated after this function call default value is true to regenerate model after each call of editing function to save time for editing for multiple times, one could set this flag to true only in the last call of editing

void validateCurrentSBML()

Validate the current SBML.

double mcaSteadyState()
double steadyState(Dictionary *dict = 0)

Compute the steady state of the model, returns the sum of squares of the solution.

The steady state solver and whatever options it needs may be specified via the given dictionary. For a list of all available steady state solvers,

Parameters:

dict – a pointer to a dictionary which has the steady state options. May be NULL, in this case the existing options are used.

ls::DoubleMatrix steadyStateNamedArray(const Dictionary *dict = 0)

Like steadyState but returns a named array of the steady state values.

std::vector<rr::SelectionRecord> &getSteadyStateSelections()

returns the current set of steady state selections.

void setSteadyStateSelections(const std::vector<std::string> &steadyStateSelections)

parses the given list of strings and generates selections records which will be used for the steady state selections.

void setSteadyStateSelections(const std::vector<rr::SelectionRecord> &steadyStateSelections)

makes a copy of an existing list of selection records.

These will be saved and used for selection values in getSteadyStateValues().

std::vector<double> getSteadyStateValues()

Performs a steady state calculation (evolves the system to a steady state), then calculates and returns the set of values specifed by the steady state selections.

std::vector<std::string> getSteadyStateSelectionStrings() const

Returns a std::vector of the steady state selection strings.

ls::DoubleMatrix getSteadyStateValuesNamedArray()

Like getSteadyStateValues but returns a named array.

void regenerateModel(bool forceRegenerate = true, bool reset = false)

regenerate a model, this recompiles the internal runtime model.

This would be used after a series of addSpecies, addReaction calls and the model is then ready to be compiled into executable form.

void setBoundarySpeciesByIndex(const int &index, const double &value)

Alias for this function on the child model object.

int getNumberOfIndependentSpecies()

Alias for this function on the child model object.

int getNumberOfDependentSpecies()

Alias for this function on the child model object.

std::vector<std::string> getGlobalParameterIds()

Alias for this function on the child model object.

use ExecutableModel::getGlobalParameterIds

std::vector<std::string> getBoundarySpeciesIds()

Alias for this function on the child model object.

std::vector<std::string> getAssignmentRuleIds()

Get the Ids of the elements assigned by assignment rules.

std::vector<std::string> getRateRuleIds()

Get the Ids of the elements assigned by rate rules.

std::vector<std::string> getInitialAssignmentIds()

Get the Ids of the elements with an initial assignment.

std::vector<std::string> getBoundarySpeciesConcentrationIds()

Gets the ids for all boundary species concentrations.

double getBoundarySpeciesByIndex(const int &index)

Alias for this function on the child model object.

double getGlobalParameterByIndex(const int &index)

Alias for this function on the child model object.

use ExecutableModel::getGlobalParameterValues

void setGlobalParameterByName(const std::string &param, double value)

set the parameter with id

Author

ciaran welsh

Parameters:
  • param – to

  • value

double getGlobalParameterByName(const std::string &param)

get the value of global parameter with id

Author

ciaran welsh

Parameters:

param

std::vector<double> getConservedMoietyValues()

Returns the sum of each conserved cycle.

Author

MTK, JKM

std::vector<std::string> getConservedMoietyIds()
void setSeed(long int seed, bool resetModel = true)

Set the value of Config::RANDOM_SEED.

int64_t getSeed(const std::string &integratorName = "")

Returns the value of Config::RANDOM_SEED.

void resetSeed()

Reset seed’s value.

Set the value of Config::RANDOM_SEED to -1

int getNumberOfReactions()

Alias for this function on the child model object.

double getReactionRate(const int &index)

Alias for this function on the child model object.

double getRateOfChange(const int &index)

Alias for this function on the child model object.

std::vector<std::string> getRateOfChangeIds()

Alias for this function on the child model object.

int getNumberOfCompartments()

Alias for this function on the child model object.

void setCompartmentByIndex(const int &index, const double &value)

Alias for this function on the child model object.

double getCompartmentByIndex(const int &index)

Alias for this function on the child model object.

std::vector<std::string> getCompartmentIds()

Alias for this function on the child model object.

int getNumberOfBoundarySpecies()

Alias for this function on the child model object.

void setBoundarySpeciesConcentrations(const std::vector<double> &values)

Alias for this function on the child model object.

void setBoundarySpeciesAmounts(const std::vector<double> &values)

Alias for this function on the child model object.

int getNumberOfFloatingSpecies()

Alias for this function on the child model object.

double getFloatingSpeciesByIndex(int index)

Alias for this function on the child model object.

void setFloatingSpeciesByIndex(int index, double value)

Alias for this function on the child model object.

std::vector<double> getFloatingSpeciesConcentrationsV()

Alias for this function on the child model object.

std::vector<double> getFloatingSpeciesAmountsV()

Alias for this function on the child model object.

std::vector<double> getBoundarySpeciesConcentrationsV()

Alias for this function on the child model object.

std::vector<double> getBoundarySpeciesAmountsV()

Alias for this function on the child model object.

std::vector<double> getFloatingSpeciesInitialConcentrations()

Alias for this function on the child model object.

void setFloatingSpeciesConcentrations(const std::vector<double> &values)

Alias for this function on the child model object.

void setFloatingSpeciesInitialConcentrationByIndex(const int &index, const double &value)

Alias for this function on the child model object.

void setFloatingSpeciesInitialConcentrations(const std::vector<double> &values)

Alias for this function on the child model object.

std::vector<std::string> getFloatingSpeciesIds()

Alias for this function on the child model object.

std::vector<std::string> getFloatingSpeciesInitialConditionIds()

Alias for this function on the child model object.

size_t getNumberOfGlobalParameters()

Alias for this function on the child model object.

void setGlobalParameterByIndex(const int index, const double value)

Alias for this function on the child model object.

std::vector<double> getGlobalParameterValues()

Alias for this function on the child model object.

use ExecutableModel::getGlobalParameterValues

void evalModel()

Updates the model based on all recent changes.

std::vector<double> getReactionRates()

Alias for this function on the child model object.

use ExecutableModel::getReactionRates

std::vector<std::string> getReactionIds()

Alias for this function on the child model object.

returns a list of reaction ids obtained from ExecutableModel::getReactionId

Public Static Functions

static std::string getParamPromotedSBML(const std::string &sArg)

Given an sbml document or path to an sbml document, this method moves all the local parameters to global parameters.

static std::vector<std::string> getRegisteredIntegratorNames()

Return a std::vector of the names of all registered integrators.

static std::vector<std::string> getRegisteredSteadyStateSolverNames()

Return a std::vector of the names of all registered SteadyStateSolver names.

static std::vector<std::string> getRegisteredSensitivitySolverNames()

Return a std::vector of the names of all registered sensitivity solvers.

static void registerSolvers()

Ensures all integrators and steady state solvers are registered.

static std::string getExtendedVersionInfo()

getVersion plus info about dependent libs versions.

Private Types

enum VariableType

Values:

enumerator vtSpecies
enumerator vtFlux
enum JacobianMode

Values:

enumerator JACOBIAN_FULL
enumerator JACOBIAN_REDUCED

Private Functions

void fixDependentSpeciesValues(int except, double *ref)
void initLLVM()

calls llvm initialization routines in a thread safe way

size_t createDefaultSteadyStateSelectionList()
size_t createDefaultTimeCourseSelectionList()

PRIVATE.

void getSelectedValues(ls::DoubleMatrix &results, int nRow, double currentTime)

copies the current selection values into the n’th row of the given matrix

void getSelectedValues(std::vector<double> &results, double currentTime)

copies the current selection values into the given std::vector.

bool populateResult()
double getNthSelectedOutput(size_t index, double currentTime)
bool isParameterUsed(const std::string &sid)
void getAllVariables(const libsbml::ASTNode *node, std::set<std::string> &ids)
int getTimeRowIndex()

Get the row index of the time variable in the output array (returns -1 if time is not selected)

double getVariableValue(const VariableType variableType, const int variableIndex)
ls::LibStructural *getLibStruct()

the LibStruct is normally null, only created on demand here.

bool createDefaultSelectionLists()

If the specified integrator does not exist, create it, and point the integrator pointer to it.

size_t createTimeCourseSelectionList()

creates a selection list from the amounts / conc / variables ivars of the SimulationOptions struct.

std::vector<SelectionRecord> getSelectionList()
void applySimulateOptions()

The simulateOptions may be changed at any time.

This updates the integrators just before they are used with the potentially changed options.

std::vector<std::complex<double>> getEigenValues(JacobianMode mode)
ls::DoubleMatrix getEigenValuesNamedArray(JacobianMode mode)
void checkID(const std::string &functionName, const std::string &sid)
void parseSpecies(const std::string &species, double *stoichiometry, char **sid)
void removeVariable(const std::string &sid)
bool hasVariable(const libsbml::ASTNode *node, const std::string &sid)
void getSpeciesIdsFromAST(const libsbml::ASTNode *node, std::vector<std::string> &species)
void getSpeciesIdsFromAST(const libsbml::ASTNode *node, std::vector<std::string> &species, std::vector<std::string> &speciesNames)
void saveSelectionVector(std::ostream&, std::vector<SelectionRecord>&)
void loadSelectionVector(std::istream&, std::vector<SelectionRecord>&)

Private Members

class RoadRunnerImpl *impl

private implementation class, can only access if inside the implementation file.

const int fileMagicNumber = 0xAD6F52
const int dataVersionNumber

Private Static Attributes

static bool llvmInitialized = false

PRIVATE.

True once the llvm initialization routines have been run. This is necessary because they are not threadsafe.

static bool solversRegistered = false