Microsimulation API
/home/marcle/src/R/microsimulation/src/person-r.cc
Go to the documentation of this file.
1 
29 #include <microsimulation.h>
30 
31 namespace {
32 
33  using namespace std;
34  using namespace ssim;
35 
38 
42 
44 enum event_t {toDeath, toPCDeath, toLocalised, toDxLocalised,
45  toDxLocallyAdvanced,
46  toLocallyAdvanced, toMetastatic, toDxMetastatic};
47 
49 class Person : public cProcess
50 {
51 public:
52  gleason_t gleason;
53  stage_t stage;
54  bool dx;
55  int id;
56 
57  // static member(s)
58  static std::map<std::string, std::vector<double> > report;
59  static std::map<std::string, Rng *> rng;
60 
61  static void resetPopulation ();
62  //
63  Person(int i = 0) : gleason(nogleason), stage(Healthy), dx(false), id(i) {};
64  void init();
65  virtual void handleMessage(const cMessage* msg);
66  virtual Time age() { return now(); }
67 };
68 
70  report.clear();
71 }
72 
73 // initialise static member(s)
74 std::map<std::string, std::vector<double> > Person::report;
75 std::map<std::string, Rng *> Person::rng;
76 
80 double dxHR(stage_t stage) {
81  // raise error if healthy?
82  return stage==Healthy ? -1 :
83  (stage==Localised ? 1.1308 :
84  (stage==LocallyAdvanced ? 0.5900 :1.3147));
85 }
86 
90 double progressionHR(gleason_t gleason) {
91  return gleason==gleasonLt7 ? 1 :
92  (gleason==gleason7 ? 1.3874 : 1.4027 * 1.3874);
93 }
94 
98 void Person::init() {
99  rng["NH"]->set();
100  if (R::runif(0.0,1.0)<0.2241)
101  scheduleAt(R::rweibull(exp(2.3525),64.0218),toLocalised);
102  scheduleAt(R::rexp(80.0),toDeath);
103 }
104 
108 void Person::handleMessage(const cMessage* msg) {
109 
110  double pDx;
111 
112  report["id"].push_back(id);
113  report["startTime"].push_back(previousEventTime);
114  report["endTime"].push_back(now());
115  report["state"].push_back(stage);
116  report["event"].push_back(msg->kind);
117 
118  rng["NH"]->set();
119 
120  if (msg->kind == toDeath) {
121  Sim::stop_simulation();
122  }
123 
124  else if (msg->kind == toPCDeath) {
125  Sim::stop_simulation();
126  }
127 
128  else if (msg->kind == toLocalised) {
129  stage = Localised;
130  gleason = (R::runif(0.0,1.0)<0.6812) ? gleasonLt7 :
131  ((R::runif(0.0,1.0)<0.5016) ? gleason7 : gleasonGt7);
132  Time dwellTime = now()+
133  rweibullHR(exp(1.0353),19.8617,progressionHR(gleason)*
134  dxHR(stage));
135  // now separate out for different transitions
136  pDx = 1.1308/(2.1308);
137  if (R::runif(0.0,1.0)<pDx) {
138  scheduleAt(dwellTime, toDxLocalised);
139  }
140  else {
141  scheduleAt(dwellTime,toLocallyAdvanced);
142  }
143  }
144 
145  else if (msg->kind == toLocallyAdvanced) {
146  stage=LocallyAdvanced;
147  Time dwellTime = now()+
148  rweibullHR(exp(1.4404),16.3863,progressionHR(gleason)*
149  dxHR(stage));
150  // now separate out for different transitions
151  pDx = 0.5900/(1.0+0.5900);
152  if (R::runif(0.0,1.0)<pDx) {
153  scheduleAt(dwellTime, toDxLocallyAdvanced);
154  }
155  else {
156  scheduleAt(dwellTime,toMetastatic);
157  }
158  }
159 
160  else if (msg->kind == toMetastatic) {
161  stage=Metastatic;
162  Time dwellTime = now()+
163  rweibullHR(exp(1.4404),1.4242,progressionHR(gleason)*
164  dxHR(stage));
165  // now separate out for different transitions
166  pDx = 1.3147/(1.0+1.3147);
167  if (R::runif(0.0,1.0)<pDx) {
168  scheduleAt(dwellTime, toDxMetastatic);
169  }
170  else {
171  scheduleAt(dwellTime,toPCDeath); // prior to diagnosis!
172  }
173  }
174 
175  else if (msg->kind == toDxLocalised) {
176  dx=true;
177  // relative survival
178  }
179 
180  else if (msg->kind == toDxLocallyAdvanced) {
181  dx=true;
182  // relative survival
183  }
184 
185  else if (msg->kind == toDxMetastatic) {
186  dx=true;
187  // relative survival
188  };
189 
190 }
191 
192 
193 extern "C" {
194 
195  RcppExport SEXP callPersonSimulation(SEXP inseed, SEXP parms) {
196  Rcpp::List parmsl(parms);
197  Rcpp::IntegerVector inseed2(inseed);
198  int nin = Rcpp::as<int>(parmsl["n"]);
199  double seed[6];
200  for (int i=0; i<6; i++) {
201  seed[i]=inseed2[i];
202  }
203  //r_create_current_stream();
204  RngStream::SetPackageSeed(seed);
206  Person::rng["NH"] = new Rng();
207  Person::rng["S"] = new Rng();
208  Person::rng["NH"]->set();
209  Person person;
210  for (int i = 0; i < nin; i++) {
211  //Person::rng.foreach(nextSubStream);
212  Person::rng["NH"]->nextSubstream();
213  Person::rng["S"]->nextSubstream();
214  person = Person(i);
215  Sim::create_process(&person);
216  Sim::run_simulation();
217  Sim::clear();
218  }
219  // tidy up -- what needs to be deleted?
220  delete Person::rng["NH"];
221  delete Person::rng["S"];
222  // Person::rng.clear();
223  // output arguments to R
224  return Rcpp::wrap(Person::report);
225 
226  } // callPersonSimulation()
227 
228 } // extern "C"
229 
230 } // namespace person_r
231 
232 namespace {
233 
234  class VerySimple : public cProcess {
235  public:
236  void init() {
237  scheduleAt(10.0, "a message");
238  scheduleAt(11.0, "another message");
239  };
240  virtual void handleMessage(const cMessage* msg) {};
241  };
242 
243 extern "C" {
244 
245  RcppExport SEXP callSpeedTest() {
246  VerySimple simple;
247  for (int i = 0; i < 1000000; i++) {
248  simple = VerySimple();
249  Sim::create_process(&simple);
250  Sim::run_simulation();
251  Sim::clear();
252  }
253  return Rcpp::wrap(1);
254 
255  } // callSpeedTest()
256 
257 } // extern "C"
258 
259 } // anonymous namespace
260 
261 
262 
nogleason
@ nogleason
Definition: person-r-20121231.cc:37
ssim::cMessage
cMessage class for OMNET++ API compatibility. This provides a heavier message class than Sim::Event,...
Definition: microsimulation.h:211
Person::handleMessage
virtual void handleMessage(const cMessage *msg)
Definition: person-r-20121231.cc:107
gleason7
@ gleason7
Definition: person-r-20121231.cc:37
gleasonLt7
@ gleasonLt7
Definition: person-r-20121231.cc:37
ssim::cProcess
cProcess class for OMNET++ API compatibility.
Definition: microsimulation.h:314
Rcpp::wrap
SEXP wrap(const std::vector< std::pair< T1, T2 > > v)
callPersonSimulation
void callPersonSimulation(int *inseed, double *parms, int *nin, double *out, int *nout)
Definition: person-r-20121231.cc:194
stage_t
stage_t
enum of type of disease stage
Definition: person-r-20121231.cc:40
progressionHR
double progressionHR(gleason_t gleason)
Definition: person-r-20121231.cc:90
Localised
@ Localised
Definition: person-r-20121231.cc:40
Metastatic
@ Metastatic
Definition: person-r-20121231.cc:41
gleason_t
gleason_t
enum for type of Gleason score
Definition: person-r-20121231.cc:37
ssim::now
Time now()
now() function for compatibility with C++SIM
Definition: microsimulation.cc:9
Death
@ Death
Definition: person-r-20121231.cc:41
Person::init
void init()
Definition: person-r-20121231.cc:98
LocallyAdvanced
@ LocallyAdvanced
Definition: person-r-20121231.cc:40
DxMetastatic
@ DxMetastatic
Definition: person-r-20121231.cc:41
DxLocallyAdvanced
@ DxLocallyAdvanced
Definition: person-r-20121231.cc:40
gleasonGt7
@ gleasonGt7
Definition: person-r-20121231.cc:37
dxHR
double dxHR(stage_t stage)
Definition: person-r-20121231.cc:80
ssim::Time
double Time
virtual time type
Definition: ssim.h:75
ssim
name space for the Siena simulator.
Definition: microsimulation.cc:3
DxLocalised
@ DxLocalised
Definition: person-r-20121231.cc:40
Person::resetPopulation
static void resetPopulation()
Definition: person-r-20121231.cc:65
ssim::rweibullHR
double rweibullHR(double shape, double scale, double hr)
Random Weibull distribution for a given shape, scale and hazard ratio.
Definition: microsimulation.cc:5
Healthy
@ Healthy
Definition: person-r-20121231.cc:40
illnessDeath::report
EventReport< short, short, double > report
Definition: illness-death.cpp:13
illnessDeath::event_t
event_t
Definition: illness-death.cpp:11
ssim::Rng
Definition: microsimulation.h:548
microsimulation.h
ssim::cMessage::kind
short kind
Definition: microsimulation.h:222