Microsimulation API
cervical.cpp
Go to the documentation of this file.
1 #include <microsimulation.h>
2 
3 #include <numeric>
4 
5 namespace {
6 
7  using namespace std;
8  using namespace Rcpp;
9  using namespace ssim;
10 
11  // declarations
12 
13  enum hpv_t {LR_HPV,HPV_16,HPV_18,Other_HR_HPV};
14 
15  enum state_t {Normal, HPV, CIN1, CIN23, LocalCancer, RegionalCancer, DistantCancer, Death};
16 
17  enum event_t {toHPV, toCIN1, toNormal, toCIN23, toNoCIN, toLocalCancer, toRegionalCancer,
18  toDistantCancer, toUtility, toUtilityChange, toOtherDeath, toCancerDeath};
19 
20  // forward declarations
21  class Person;
22  class HPV_infection;
23 
24  // namespace FullState {
25  // typedef std::tuple<short,short,short,bool,double> Type;
26  // enum Fields {state, ext_grade, dx, psa_ge_3, cohort};
27  // // string names[5] = {"state","ext_grade","dx","psa_ge_3","cohort"};
28  // }
29  // EventReport<FullState::Type,short,double> report;
30  // typedef pair<string,double> CostKey;
31  // CostReport<CostKey> costs;
32 
33  typedef std::tuple<hpv_t,state_t,state_t> H_key;
34  typedef map<H_key,NumericInterpolate> H_t;
35  typedef vector<HPV_infection *> infections_t;
36 
37  class SimInput {
38  public:
39  bool debug;
40 
41  H_t H;
42  Rng * rngNh, * rngOther, * rngScreen, * rngTreatment;
43  Rpexp rmu0;
44 
45  NumericVector parameter;
46  LogicalVector bparameter;
47 
48  SimInput() {
49  debug = false;
50  rngNh = new Rng();
51  rngOther = new Rng();
52  rngScreen = new Rng();
53  rngTreatment = new Rng();
54  }
55 
56  virtual ~SimInput() {
57  delete rngNh;
58  delete rngOther;
59  delete rngScreen;
60  delete rngTreatment;
61  }
62  };
63 
64  class SimOutput {
65  public:
66  SimpleReport<double> lifeHistories;
67  SimpleReport<double> outParameters;
68  };
69 
70  // typedef Table<std::tuple<double,double,int>,double> TablePrtx; // Age, DxY, G
71  // typedef Table<std::tuple<int,double,double,int>,double> TablePradt;
72  // typedef Table<pair<double,double>,double> TableBiopsyCompliance;
73  // typedef Table<pair<double,double>,double> TableDDD; // as per TableBiopsyCompliance
74  // typedef map<int,NumericInterpolate> H_dist_t;
75  // typedef map<pair<double,int>,NumericInterpolate> H_local_t;
76  // TablePrtx prtxCM, prtxRP;
77  // TablePradt pradt;
78  // TableBiopsyCompliance tableOpportunisticBiopsyCompliance, tableFormalBiopsyCompliance;
79  // TableDDD rescreen_shape, rescreen_scale, rescreen_cure;
80  // NumericInterpolate interp_prob_grade7;
81  // H_dist_t H_dist;
82  // H_local_t H_local;
83  // set<double,greater<double> > H_local_age_set;
84 
85 
86  class cMessageByHPV : public cMessage {
87  public:
88  cMessageByHPV(event_t event, hpv_t hpv) : cMessage(event), hpv(hpv) { }
89  hpv_t hpv;
90  };
91 
92  inline bool cMessageHPVPred(const ssim::Event* e, const hpv_t hpv) {
93  const cMessageByHPV * msg = dynamic_cast<const cMessageByHPV *>(e);
94  return (msg != 0 && msg->hpv == hpv);
95  }
96 
97  inline bool cMessageAnyHPVPred(const ssim::Event* e) {
98  return (dynamic_cast<const cMessageByHPV *>(e) != 0);
99  }
100 
101  class cMessageUtilityChange : public cMessage {
102  public:
103  cMessageUtilityChange(double change) : cMessage(toUtilityChange), change(change) { }
104  double change;
105  };
106 
107  class cMessageUtility : public cMessage {
108  public:
109  cMessageUtility(double utility) : cMessage(toUtility), utility(utility) { }
110  double utility;
111  };
112 
113  // utilities
114  template<class T>
115  T bounds(T x, T a, T b) {
116  return (x<a)?a:((x>b)?b:x);
117  }
118  template<class T>
119  T max(T left, T right) {
120  return (left<right) ? right : left;
121  }
122 
123  class HPV_infection : public cProcess
124  {
125  public:
126  hpv_t hpv;
127  state_t state;
128  bool ever_infected, immunity;
129  double p_immunity;
130  Person * person;
131  SimInput * in;
132  SimOutput * out;
133  HPV_infection(hpv_t hpv = LR_HPV, double p_immunity = 1.0) :
134  hpv(hpv), p_immunity(p_immunity) {}
135  double event_time(state_t from, state_t to);
136  void init();
137  virtual void handleMessage(const cMessage* msg);
138  void removeEvents() {
139  cancel([hpv](const cMessage * msg) { return this->hpv == hpv; });
140  }
141  };
142 
143  class Person : public cProcess
144  {
145  public:
146  SimInput * in;
147  SimOutput * out;
148  int id;
149  double cohort, utility;
150  bool dx;
151  state_t state;
152  infections_t infections;
153  Person(SimInput * in = NULL, SimOutput * out = NULL, const int id = 0, const double cohort = 1950) :
154  in(in), out(out), id(id), cohort(cohort), utility(1.0) { }
155  void init();
156  virtual void handleMessage(const cMessage* msg);
157  void link(HPV_infection * infection) {
158  infection->out = out;
159  infection->in = in;
160  infections.push_back(infection);
161  infection->person=this;
162  }
163  state_t max_infection_state();
164  void removeHPVevents() {
165  // for (infections_t::iterator it = infections.begin(); it!=infections.end(); ++it)
166  // (*it)->removeEvents();
167  Sim::ignore_event(cMessageAnyHPVPred);
168  }
169  // void add_costs(string item);
170  // void scheduleUtilityChange(double at, std::string category, bool transient = true,
171  // double sign = -1.0);
172  };
173 
174  // state_t Person::max_infection_state() {
175  // int max_value = (int) Normal;
176  // for (vector<HPV_infection *>::iterator it = infections.begin(); it!=infections.end(); ++it)
177  // max_value = max<int>((int) (*it)->state, max_value);
178  // return (state_t) max_value;
179  // }
180 
181  double HPV_infection::event_time(state_t from, state_t to) {
182  return in->H[H_key(hpv,from,to)].invert(-log(R::runif(0.0,1.0)),now());
183  }
184 
185  void HPV_infection::init() {
186  // initialise state
187  state = Normal;
188  ever_infected = false;
189  immunity = false;
190  // determine event times
191  in->rngNh->set();
192  scheduleAt(event_time(Normal,HPV), new cMessageByHPV(toHPV,hpv));
193  }
194 
195  void HPV_infection::handleMessage(const cMessage* msg) {
196 
197  in->rngNh->set();
198 
199  if (person->id < in->parameter["nLifeHistories"]) { // only record up to the first n individuals
200  out->lifeHistories.record("id",person->id);
201  out->lifeHistories.record("time",now());
202  out->lifeHistories.record("hpv",hpv);
203  out->lifeHistories.record("event",msg->kind);
204  }
205 
206  switch(msg->kind) {
207  case toNormal:
208  state = Normal;
209  removeEvents();
210  if (!immunity)
211  scheduleAt(event_time(Normal,HPV),
212  new cMessageByHPV(toHPV,hpv));
213  break;
214  case toHPV:
215  state = HPV;
216  ever_infected = true;
217  removeEvents();
218  scheduleAt(event_time(HPV,Normal),
219  new cMessageByHPV(toNormal,hpv));
220  scheduleAt(event_time(HPV,CIN1),
221  new cMessageByHPV(toCIN1,hpv));
222  scheduleAt(event_time(HPV,CIN23),
223  new cMessageByHPV(toCIN23,hpv));
224  break;
225  case toCIN1:
226  state = CIN1;
227  removeEvents();
228  scheduleAt(event_time(CIN1,Normal),
229  new cMessageByHPV(toNormal,hpv));
230  if (!immunity)
231  scheduleAt(event_time(CIN1,HPV),
232  new cMessageByHPV(toHPV,hpv));
233  scheduleAt(event_time(CIN1,CIN23),
234  new cMessageByHPV(toCIN23,hpv));
235  break;
236  case toCIN23:
237  state = CIN23;
238  removeEvents();
239  scheduleAt(event_time(CIN1,Normal),
240  new cMessageByHPV(toNormal,hpv));
241  if (!immunity)
242  scheduleAt(event_time(CIN1,HPV),
243  new cMessageByHPV(toHPV,hpv));
244  scheduleAt(event_time(CIN1,CIN23),
245  new cMessageByHPV(toCIN23,hpv));
246  break;
247  case toLocalCancer:
248  state = LocalCancer;
249  person->removeHPVevents();
250  person->scheduleAt(now(), new cMessage(toLocalCancer));
251  break;
252  };
253  }
254 
258  // void Person::add_costs(string item) {
259  // costs.add(CostKey(item,cohort),now(),cost_parameters[item]);
260  // }
261 
266  // void Person::scheduleUtilityChange(double at, std::string category, bool transient, double sign) {
267  // scheduleAt(at, new cMessageUtilityChange(sign*utility_estimates[category]));
268  // if (transient) {
269  // scheduleAt(at + utility_duration[category],
270  // new cMessageUtilityChange(-sign*utility_estimates[category]));
271  // }
272  // }
273 
277  void Person::init() {
278 
279  // declarations
280 
281  // initialise state variables
282  dx = false;
283 
284  // determine event times
285  in->rngNh->set();
286 
287  // schedule natural history events
288  double aoc = in->rmu0.rand(R::runif(0.0,1.0));
289  scheduleAt(aoc, toOtherDeath);
290 
291  // schedule screening events
292  // rngScreen->set();
293 
294  in->rngNh->set();
295 
296  // // utilities
297  // // (i) set initial baseline utility
298  // utility = 0.98;
299  // // (ii) schedule changes in the baseline by age
300  // scheduleAt(20.0, new cMessageUtility(0.97));
301  // scheduleAt(40.0, new cMessageUtility(0.96));
302  // scheduleAt(60.0, new cMessageUtility(0.95));
303  // scheduleAt(80.0, new cMessageUtility(0.91));
304  // // Mark, which reference do/should we use for this?
305 
306  // record some parameters using SimpleReport - too many for a tuple
307  if (id < in->parameter["nLifeHistories"]) {
308  out->outParameters.record("id",double(id));
309  out->outParameters.record("aoc",aoc);
310  out->outParameters.record("cohort",cohort);
311  }
312  }
313 
317 void Person::handleMessage(const cMessage* msg) {
318 
319  // by default, use the natural history RNG
320  in->rngNh->set();
321 
322  // declarations
323 
324  // reports
325  // if (bparameter["full_report"])
326  // report.add(FullState::Type(state, ext_grade, dx, psa>=3.0, cohort), msg->kind, previousEventTime, age, utility);
327  if (id < in->parameter["nLifeHistories"]) { // only record up to the first n individuals
328  out->lifeHistories.record("id",id);
329  out->lifeHistories.record("time",now());
330  out->lifeHistories.record("hpv",-1.0);
331  out->lifeHistories.record("event",msg->kind);
332  }
333 
334  // handle messages by kind
335 
336  switch(msg->kind) {
337 
338  case toCancerDeath:
339  // add_costs("CancerDeath"); // cost for death, should this be zero???
340  if (id < in->parameter["nLifeHistories"]) {
341  out->outParameters.record("dx",dx);
342  out->outParameters.record("age_d",now());
343  out->outParameters.record("cancer_death",1.0);
344  }
345  Sim::stop_simulation();
346  break;
347 
348  case toOtherDeath:
349  if (id < in->parameter["nLifeHistories"]) {
350  out->outParameters.record("dx",dx);
351  out->outParameters.record("age_d",now());
352  out->outParameters.record("cancer_death",0.0);
353  }
354  Sim::stop_simulation();
355  break;
356 
357  case toLocalCancer:
358  state = LocalCancer;
359  // progression?
360  break;
361 
362  case toRegionalCancer:
363  state = RegionalCancer;
364  // progression?
365  break;
366 
367  case toDistantCancer:
368  state = DistantCancer;
369  // progression?
370  break;
371 
372 
373  // case toUtility:
374  // {
375  // const cMessageUtility * msgUtility;
376  // if ((msgUtility = dynamic_cast<const cMessageUtility *>(msg)) != 0) {
377  // utility = msgUtility->utility;
378  // } else {
379  // REprintf("Could not cast to cMessageUtility.");
380  // }
381  // } break;
382 
383  // case toUtilityChange:
384  // {
385  // const cMessageUtilityChange * msgUtilityChange;
386  // if ((msgUtilityChange = dynamic_cast<const cMessageUtilityChange *>(msg)) != 0) {
387  // utility += msgUtilityChange->change;
388  // } else {
389  // REprintf("Could not cast to cMessageUtilityChange.");
390  // }
391  // } break;
392 
393  default:
394  REprintf("No valid kind of event: %i\n",msg->kind);
395  break;
396 
397  } // switch
398 
399 } // handleMessage()
400 
401 
402 RcppExport SEXP callCervical(SEXP parmsIn) {
403 
404  // declarations
405  Person person;
406  HPV_infection lr_hpv, hpv_16, hpv_18, other_hr_hpv;
407  SimInput in;
408  SimOutput out;
409 
410  in.rngNh->set();
411 
412  // read in the parameters
413  List parms(parmsIn);
414  List tables = parms["tables"];
415  in.parameter = parms["parameter"];
416  in.bparameter = parms["bparameter"]; // scalar bools
417  List otherParameters = parms["otherParameters"];
418  in.debug = as<bool>(in.bparameter["debug"]);
419  NumericVector mu0 = as<NumericVector>(otherParameters["mu0"]);
420  // cost_parameters = as<NumericVector>(otherParameters["cost_parameters"]);
421  // utility_estimates = as<NumericVector>(otherParameters["utility_estimates"]);
422  // utility_duration = as<NumericVector>(otherParameters["utility_duration"]);
423 
424  {
425  in.H.clear();
426  DataFrame df_H = as<DataFrame>(tables["H"]); // age,hpv,from_state,to_state,survival
427  // extract the columns
428  IntegerVector hpv = df_H["hpv"],
429  from_state = df_H["from_state"],
430  to_state = df_H["to_state"];
431  NumericVector age = df_H["age"],
432  survival = df_H["survival"];
433  typedef pair<double,double> dpair;
434  // push to the map values and set of ages
435  for (int i=0; i<age.size(); ++i) {
436  in.H[H_key((hpv_t) hpv[i], (state_t) from_state[i], (state_t) to_state[i])].push_back(dpair(age[i],-log(survival[i])));
437  }
438  // prepare the map values for lookup
439  for (H_t::iterator it = in.H.begin(); it != in.H.end(); it++)
440  it->second.prepare();
441  // to use: event_time(from_state,to_state)
442  }
443 
444  NumericVector cohort = as<NumericVector>(parms["cohort"]); // at present, this is the only chuck-specific data
445 
446  // set up the parameters
447  double ages0[106];
448  std::iota(ages0, ages0+106, 0.0);
449  in.rmu0 = Rpexp(&mu0[0], ages0, 106);
450  vector<double> ages(101);
451  std::iota(ages.begin(), ages.end(), 0.0);
452  ages.push_back(1.0e+6);
453 
454  // re-set the output objects
455  // report.clear();
456  // costs.clear();
457  out.outParameters.clear();
458  out.lifeHistories.clear();
459 
460  // report.discountRate = parameter["discountRate.effectiveness"];
461  // report.setPartition(ages);
462  // costs.discountRate = parameter["discountRate.costs"];
463  // costs.setPartition(ages);
464 
465  // main loop
466  for (int i = 0; i < in.parameter["n"]; ++i) {
467  // define infections
468  int id = i+as<int>(in.parameter["firstID"]);
469  lr_hpv = HPV_infection(LR_HPV);
470  hpv_16 = HPV_infection(HPV_16);
471  hpv_18 = HPV_infection(HPV_18);
472  other_hr_hpv = HPV_infection(Other_HR_HPV);
473  // define person
474  person = Person(&in,&out,id,cohort[i]);
475  // link infections to the person
476  person.link(&lr_hpv);
477  person.link(&hpv_16);
478  person.link(&hpv_18);
479  person.link(&other_hr_hpv);
480  // create processes
481  Sim::create_process(&person);
482  Sim::create_process(&lr_hpv);
483  Sim::create_process(&hpv_16);
484  Sim::create_process(&hpv_18);
485  Sim::create_process(&other_hr_hpv);
486  // run!!
487  Sim::run_simulation();
488  // tidy up
489  Sim::clear();
490  in.rngNh->nextSubstream();
491  in.rngOther->nextSubstream();
492  in.rngScreen->nextSubstream();
493  in.rngTreatment->nextSubstream();
494  if (i % 10000 == 0) Rcpp::checkUserInterrupt(); /* be polite -- did the user hit ctrl-C? */
495  }
496 
497  // output
498  // TODO: clean up these objects in C++ (cf. R)
499  return List::create(
500  // _("costs") = costs.wrap(), // CostReport
501  // _("summary") = report.wrap(), // EventReport
502  _("lifeHistories") = out.lifeHistories.wrap(), // SimpleReport<double>
503  _("parameters") = out.outParameters.wrap() // SimpleReport<double>
504  );
505 }
506 
507 } // anonymous namespace
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
it
GNU Free Documentation License March Inc Temple MA USA Everyone is permitted to copy and distribute verbatim copies of this license but changing it is not allowed PREAMBLE The purpose of this License is to make a or other written document free in the sense of with or without modifying it
Definition: fdl.txt:15
ssim::cProcess
cProcess class for OMNET++ API compatibility.
Definition: microsimulation.h:314
illnessDeath::toOtherDeath
@ toOtherDeath
Definition: illness-death.cpp:11
ssim::now
Time now()
now() function for compatibility with C++SIM
Definition: microsimulation.cc:9
illnessDeath::state_t
state_t
Definition: illness-death.cpp:9
ssim::Event
basic event in the simulation.
Definition: ssim.h:106
ssim::SimpleReport
SimpleReport class for collecting data for homogeneous fields of type T with string names.
Definition: microsimulation.h:1099
Death
@ Death
Definition: person-r-20121231.cc:41
Person::init
void init()
Definition: person-r-20121231.cc:98
ssim::Rpexp
Rpexp is a random number generator class for piecewise constant hazards. Given time lower bounds t an...
Definition: microsimulation.h:491
illnessDeath::toCancerDeath
@ toCancerDeath
Definition: illness-death.cpp:11
ssim
name space for the Siena simulator.
Definition: microsimulation.cc:3
illnessDeath::event_t
event_t
Definition: illness-death.cpp:11
Rcpp
Definition: microsimulation.h:101
ssim::Rng
Definition: microsimulation.h:548
microsimulation.h
ssim::cMessage::kind
short kind
Definition: microsimulation.h:222