libRoadRunner C++ API  1.0.0
 All Classes Functions Variables Enumerations Enumerator Pages
rrRoadRunner.h
1 #ifndef rrRoadRunnerH
2 #define rrRoadRunnerH
3 
4 #include "rr-libstruct/lsMatrix.h"
5 #include "rrVariableType.h"
6 #include "rrParameterType.h"
7 #include "rrSelectionRecord.h"
8 #include "rrRoadRunnerData.h"
9 #include "rrConstants.h"
10 #include "rrRoadRunnerOptions.h"
11 #include "Configurable.h"
12 
13 #include <string>
14 #include <vector>
15 #include <list>
16 
17 namespace ls
18 {
19 class LibStructural;
20 }
21 
22 namespace rr
23 {
24 
25 class ModelGenerator;
26 class SBMLModelSimulation;
27 class ExecutableModel;
28 class Integrator;
29 
38 class RR_DECLSPEC RoadRunner : public Configurable
39 {
40 
41 public:
42 
53  RoadRunner(const std::string& compiler = "", const std::string& tempDir = "",
54  const std::string& supportCodeDir = "");
55 
59  virtual ~RoadRunner();
60 
64  int getInstanceID();
65 
69  int getInstanceCount();
70 
75  static std::string getParamPromotedSBML(const std::string& sArg);
76 
80  std::string getInfo();
81 
85  class Compiler* getCompiler();
86 
91  bool setCompiler(const std::string& compiler);
92 
97  Integrator* getIntegrator();
98 
99  bool isModelLoaded();
100 
104  std::string getModelName();
105 
106  bool unLoadModel();
107 
116  bool setTempFileFolder(const std::string& folder);
117 
125  std::string getTempFolder();
126 
132  double oneStep(double currentTime, double stepSize, bool reset = true);
133 
145  const RoadRunnerData *simulate(const SimulateOptions* options = 0);
146 
152  RoadRunnerData *getSimulationResult();
153 
154  void setSimulateOptions(const SimulateOptions& settings);
155 
160  SimulateOptions& getSimulateOptions();
161 
165  std::string getSBML();
166 
171  void reset();
172 
181  void changeInitialConditions(const std::vector<double>& ic);
182 
186  ModelGenerator* getModelGenerator();
187 
191  ExecutableModel* getModel();
192 
202  bool load(const std::string& uriOrSBML,
203  const LoadSBMLOptions* options = 0);
204 
209  std::vector<double> getReactionRates();
210 
215  std::vector<double> getRatesOfChange();
216 
223  std::vector<std::string> getReactionIds();
224 
231  virtual _xmlNode *createConfigNode();
232 
238  virtual void loadConfig(const _xmlDoc* doc);
239 
247  std::string getConfigurationXML();
248 
254  void setConfigurationXML(const std::string& xml);
255 
256 
261  void correctMaxStep();
262 
263 /************************ Selection Ids Species Section ***********************/
264 #if (1) /**********************************************************************/
265 /******************************************************************************/
266 
270  rr::SelectionRecord createSelection(const std::string& str);
271 
276  std::vector<rr::SelectionRecord>& getSelections();
277 
282  double getValue(const std::string& sel);
283 
284  double getValue(const SelectionRecord& record);
285 
286 
287  void setSelections(const std::vector<std::string>& selections);
288 
289  void setSelections(const std::vector<rr::SelectionRecord>& selections);
290 
291 
295  std::vector<double> getSelectedValues();
296 
300  void getIds(int types, std::list<std::string> &ids);
301 
305  int getSupportedIdTypes();
306 
307 
311  bool setValue(const std::string& id, double value);
312 
313 /************************ End Selection Ids Species Section *******************/
314 #endif /***********************************************************************/
315 /******************************************************************************/
316 
320  ls::DoubleMatrix getFullJacobian();
321 
322  ls::DoubleMatrix getFullReorderedJacobian();
323 
327  ls::DoubleMatrix getReducedJacobian();
328 
332  ls::DoubleMatrix getEigenvalues();
333 
334  std::vector<Complex> getEigenvaluesCpx();
335 
336  ls::DoubleMatrix* getLinkMatrix();
337  ls::DoubleMatrix* getNrMatrix();
338  ls::DoubleMatrix* getL0Matrix();
339 
340 
341  ls::DoubleMatrix getConservationMatrix();
342  ls::DoubleMatrix getUnscaledConcentrationControlCoefficientMatrix();
343  ls::DoubleMatrix getScaledConcentrationControlCoefficientMatrix();
344  ls::DoubleMatrix getUnscaledFluxControlCoefficientMatrix();
345  ls::DoubleMatrix getScaledFluxControlCoefficientMatrix();
346  int getNumberOfDependentSpecies();
347  int getNumberOfIndependentSpecies();
348 
349 
354  std::vector<std::string> getEigenvalueIds();
355 
360  double getUnscaledParameterElasticity(const string& reactionName,
361  const string& parameterName);
362 
363 
364  Matrix<double> getFrequencyResponse(double startFrequency,
365  int numberOfDecades, int numberOfPoints,
366  const string& parameterName, const string& variableName,
367  bool useDB, bool useHz);
368 
372  void setConservedMoietyAnalysis(bool value);
373 
377  bool getConservedMoietyAnalysis();
378 
382  std::string getCurrentSBML();
383 
388  int getNumberOfReactions();
389 
394  double getReactionRate(const int& index);
395 
404  double getRateOfChange(const int& index);
405 
409  std::vector<std::string> getRateOfChangeIds();
410 
411 
412  std::vector<std::string> getConservedMoietyIds();
413 
414  std::vector<double> getConservedMoietyValues();
415 
416  int getNumberOfCompartments();
417 
423  void setCompartmentByIndex(const int& index, const double& value);
424 
430  double getCompartmentByIndex(const int& index);
431 
436  std::vector<std::string> getCompartmentIds();
437 
443  int getNumberOfBoundarySpecies();
444 
449  void setBoundarySpeciesByIndex(const int& index, const double& value);
450 
455  double getBoundarySpeciesByIndex(const int& index);
456 
461  std::vector<double> getBoundarySpeciesConcentrations();
462 
467  void setBoundarySpeciesConcentrations(const std::vector<double>& values);
468 
473  std::vector<std::string> getBoundarySpeciesIds();
474 
479  int getNumberOfFloatingSpecies();
480 
486  double getFloatingSpeciesByIndex(int index);
487 
492  void setFloatingSpeciesByIndex(int index, double value);
493 
498  std::vector<double> getFloatingSpeciesConcentrations();
499 
504  std::vector<double> getFloatingSpeciesInitialConcentrations();
505 
510  void setFloatingSpeciesConcentrations(const std::vector<double>& values);
511 
516  void setFloatingSpeciesInitialConcentrationByIndex(const int& index,
517  const double& value);
518 
523  void setFloatingSpeciesInitialConcentrations(const std::vector<double>& values);
524 
529  std::vector<std::string> getFloatingSpeciesIds();
530 
535  std::vector<std::string> getFloatingSpeciesInitialConditionIds();
536 
541  int getNumberOfGlobalParameters();
542 
547  void setGlobalParameterByIndex(const int index, const double value);
548 
553  double getGlobalParameterByIndex(const int& index);
554 
559  std::vector<double> getGlobalParameterValues();
560 
565  std::vector<std::string> getGlobalParameterIds();
566 
573  void evalModel();
574 
578  static std::string getExtendedVersionInfo();
579 
580 
589  double getuCC(const std::string& variableName, const std::string& parameterName);
590 
599  double getCC(const std::string& variableName, const std::string& parameterName);
600 
604  double getuEE(const std::string& reactionName, const std::string& parameterName);
605 
610  double getuEE(const std::string& reactionName, const std::string& parameterName,
611  bool computeSteadystate);
612 
616  double getEE(const std::string& reactionName, const std::string& parameterName);
617 
622  double getEE(const std::string& reactionName, const std::string& parameterName,
623  bool computeSteadyState);
624 
628  ls::DoubleMatrix getUnscaledElasticityMatrix();
629 
633  ls::DoubleMatrix getScaledElasticityMatrix();
634 
638  double getScaledFloatingSpeciesElasticity(const std::string& reactionName,
639  const std::string& speciesName);
640 
646  double getUnscaledSpeciesElasticity(int reactionId, int speciesIndex);
647 
648 
649  /******************************* Steady State Section *************************/
650  #if (1) /**********************************************************************/
651  /******************************************************************************/
652 
657  double steadyState();
658 
662  std::vector<rr::SelectionRecord>& getSteadyStateSelections();
663 
668  void setSteadyStateSelections(const std::vector<std::string>&
669  steadyStateSelections);
670 
675  void setSteadyStateSelections(const std::vector<rr::SelectionRecord>&
676  steadyStateSelections);
677 
683  std::vector<double> getSteadyStateValues();
684 
685  /******************************* End Steady State Section *********************/
686  #endif /***********************************************************************/
687  /******************************************************************************/
688 
689 private:
690  static int mInstanceCount;
691  int mInstanceID;
692  bool mUseKinsol;
693  const double mDiffStepSize;
694 
695  const double mSteadyStateThreshold;
696  ls::DoubleMatrix mRawRoadRunnerData;
697  RoadRunnerData mRoadRunnerData;
698 
699  std::string mCurrentSBMLFileName;
700 
705  class CvodeInterface *mCVode;
706  std::vector<SelectionRecord> mSelectionList;
707 
711  ModelGenerator *mModelGenerator;
712 
716  bool mComputeAndAssignConservationLaws;
717 
718  std::vector<SelectionRecord> mSteadyStateSelection;
719 
720  ExecutableModel* mModel;
721 
722  std::string mCurrentSBML;
723 
727  LibStructural* mLS;
728 
729  SimulateOptions mSettings;
730 
731 
732  int createDefaultSteadyStateSelectionList();
733  int createDefaultTimeCourseSelectionList();
734 
735  void addNthOutputToResult(ls::DoubleMatrix& results, int nRow,
736  double dCurrentTime);
737  bool populateResult();
738 
739 
740  double getNthSelectedOutput(const int& index, const double& dCurrentTime);
741 
742  double getVariableValue(const VariableType::VariableType variableType,
743  const int variableIndex);
744 
745 
746 
747  std::string createModelName(const std::string& mCurrentSBMLFileName);
748 
752  LibStructural* getLibStruct();
753 
754  bool initializeModel();
755 
756  bool createDefaultSelectionLists();
757 
762  int createTimeCourseSelectionList();
763 
769  void setParameterValue(const ParameterType::ParameterType parameterType,
770  const int parameterIndex, const double value);
771 
772  double getParameterValue(const ParameterType::ParameterType parameterType,
773  const int parameterIndex);
774 
778  void changeParameter(ParameterType::ParameterType parameterType,
779  int reactionIndex, int parameterIndex, double originalValue,
780  double increment);
781 
782 
783  std::vector<SelectionRecord> getSelectionList();
784 
790  std::string configurationXML;
791 
792 
793  friend class aFinalizer;
794 };
795 
796 }
797 
798 #endif
799 
Definition: rrModelGenerator.h:18
Definition: rrCompiler.h:25
Definition: rrRoadRunnerOptions.h:30
Definition: rrExecutableModel.h:26
Definition: rrRoadRunner.h:38
Definition: rrSelectionRecord.h:15
Definition: Configurable.h:65
Definition: rrRoadRunnerOptions.h:154