Microsimulation API
/home/marcle/src/R/microsimulation/src/calibperson-r.cc
Go to the documentation of this file.
1 
32 #include "microsimulation.h"
33 #include <Rcpp.h>
34 
35 namespace {
36 
37 using namespace ssim;
38 using namespace std;
39 
41 enum stage_t {DiseaseFree,Precursor,PreClinical,Clinical,Death};
42 
44 enum event_t {toPrecursor, toPreClinical, toClinical, toDeath, Count};
45 
47 string stage_names[5] = {"DiseaseFree","Precursor","PreClinical","Clinical","Death"};
48 
50  Rng * rng;
51 
53 class CalibPerson : public cProcess
54 {
55 public:
56  stage_t stage;
57  bool diseasepot;
58  double Lam1,sigm1,p2,lam2,mu3,tau3;
59  double clinTime;
60  int id;
61 
62  // static member(s)
63  static std::map<std::string, std::vector<double> > report;
64 
65  static void resetPopulation ();
66 
67  CalibPerson() {} // default constructor
68 
69  CalibPerson(double *par, int i=0) {
70  Lam1=par[0];
71  sigm1=par[1];
72  p2=par[2];
73  lam2=par[3];
74  mu3=par[4];
75  tau3=par[5];
76  id=i;
77  };
78 
79  void init();
80  virtual void handleMessage(const cMessage* msg);
81  virtual Time age() { return now(); }
82 };
83 
84 void CalibPerson::resetPopulation() {
85  report.clear();
86 }
87 
88 // initialise static member(s)
89 std::map<std::string, std::vector<double> > CalibPerson::report;
90 
94 void CalibPerson::init() {
95  if (R::runif(0,1)<p2){
96  diseasepot = true;}
97  else {
98  diseasepot = false;
99  }
100  clinTime=1000;
101  stage=DiseaseFree;
102  double lam1 = exp(R::rnorm(Lam1,sigm1));
103  scheduleAt(R::rexp(lam1), toPrecursor);
104  double x = R::runif(0,1);
105  scheduleAt((65 - 15*log(-log(x))), toDeath); //Gumbel
106  for(int i=10;i<110;i=i+10){
107  scheduleAt(i,Count);
108  }
109 }
110 
114 void CalibPerson::handleMessage(const cMessage* msg) {
115 
116  double ctime[] = {20,40,60,80};
117  int cind;
118 
119 
120  if (msg->kind == toDeath) {
121  stage=Death;
122  clinTime=std::min(clinTime,now());
123 
124  for(unsigned int i=0; i<4 ; i++){
125  if(i < report["TimeAtRisk"].size()){
126  report["TimeAtRisk"][i] += std::min(ctime[i],clinTime);
127  }
128  else {
129  report["TimeAtRisk"].push_back(std::min(ctime[i],clinTime));
130  }
131 
132  if(clinTime < ctime[i]){
133  break;
134  }
135 
136  }
137 
139  }
140 
141  else if (msg->kind == toPrecursor) {
142  stage = Precursor;
143  if (diseasepot){
144  simtime_t dwellTime = now()+ R::rexp(lam2);
145  scheduleAt(dwellTime, toPreClinical);
146  }
147  }
148 
149  else if (msg->kind == toPreClinical) {
150  stage=PreClinical;
151  simtime_t dwellTime = now()+ exp(R::rnorm(mu3,tau3*mu3));
152  scheduleAt(dwellTime, toClinical);
153  }
154 
155  else if (msg->kind == toClinical) {
156  stage=Clinical;
157  clinTime = now();
158  string stagestr = stage_names[stage];
159  }
160 
161  else if (msg->kind == Count){
162  cind = min(9,int(now()/10 - 1));
163  string stagestr = stage_names[stage];
164 
165  if(report.find(stagestr) == report.end()){ //key not found
166  report[stagestr].assign(10,0);
167  }
168  report[stagestr][cind]+=1;
169  }
170 }
171 
172 
173 extern "C" {
174 
175  RcppExport SEXP callCalibrationSimulation(SEXP parms) {
176  Rcpp::List parmsl(parms);
177  int nin = Rcpp::as<int>(parmsl["n"]);
178  std::vector<double> par = Rcpp::as<std::vector<double> >(parmsl["runpar"]);
179 
180  CalibPerson::resetPopulation();
181  rng = new Rng();
182  rng->set();
183 
184  CalibPerson::report.insert(make_pair("TimeAtRisk", std::vector<double>()));
185 
186  CalibPerson person;
187  for (int i = 0; i < nin; i++) {
188  person = CalibPerson(&par[0],0);
189  rng->nextSubstream();
190  Sim::create_process(&person);
192  Sim::clear();
193  }
194 
195  delete rng;
196 
198 
199  }
200 
201 }
202 
203 }
ssim::simtime_t
Time simtime_t
simtime_t typedef for OMNET++ API compatibility
Definition: microsimulation.h:284
callCalibrationSimulation
SEXP callCalibrationSimulation(SEXP)
ssim::Rng
Definition: microsimulation.h:393
stage_t
stage_t
enum of type of disease stage
Definition: person-r-20121231.cc:40
ssim::now
Time now()
now() function for compatibility with C++SIM
Definition: microsimulation.cc:14
ssim::cProcess
cProcess class for OMNET++ API compatibility. This provides a default for Process::process_event() th...
Definition: microsimulation.h:236
ssim::Sim::create_process
static ProcessId create_process(Process *)
creates a new process
Definition: ssim.cc:110
Death
@ Death
Definition: person-r-20121231.cc:41
Rcpp::wrap
SEXP wrap(const std::vector< std::pair< T1, T2 > > v)
microsimulation.h
ssim::Time
double Time
virtual time type
Definition: ssim.h:79
ssim::Rng::set
void set()
Definition: microsimulation.cc:31
ssim::Rng::nextSubstream
void nextSubstream()
Definition: microsimulation.h:405
ssim::cMessage::kind
short kind
Definition: microsimulation.h:199
illnessDeath::report
EventReport< short, short, double > report
Definition: illness-death.cpp:13
ssim::Sim::stop_simulation
static void stop_simulation()
stops execution of the simulation
Definition: ssim.cc:253
illnessDeath::event_t
event_t
Definition: illness-death.cpp:11
ssim::Sim::clear
static void clear()
clears out internal data structures
Definition: ssim.cc:117
ssim
name space for the Siena simulator.
Definition: microsimulation.cc:8
ssim::Sim::run_simulation
static void run_simulation()
starts execution of the simulation
Definition: ssim.cc:150
ssim::cMessage
cMessage class for OMNET++ API compatibility. This provides a heavier message class than Sim::Event,...
Definition: microsimulation.h:191