model.h Source File

CPP API: model.h Source File
ode_secir/model.h
Go to the documentation of this file.
1 /*
2 * Copyright (C) 2020-2026 MEmilio
3 *
4 * Authors: Daniel Abele, Jan Kleinert, Martin J. Kuehn
5 *
6 * Contact: Martin J. Kuehn <Martin.Kuehn@DLR.de>
7 *
8 * Licensed under the Apache License, Version 2.0 (the "License");
9 * you may not use this file except in compliance with the License.
10 * You may obtain a copy of the License at
11 *
12 * http://www.apache.org/licenses/LICENSE-2.0
13 *
14 * Unless required by applicable law or agreed to in writing, software
15 * distributed under the License is distributed on an "AS IS" BASIS,
16 * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
17 * See the License for the specific language governing permissions and
18 * limitations under the License.
19 */
20 #ifndef ODESECIR_MODEL_H
21 #define ODESECIR_MODEL_H
22 
29 #include "ode_secir/parameters.h"
30 #include "memilio/math/smoother.h"
33 
34 #include <numbers>
35 #include <complex>
36 
37 namespace mio
38 {
39 namespace osecir
40 {
41 
42 // Create template specializations for the age resolved
43 // SECIHURD model
44 // clang-format off
60 // clang-format on
61 
62 template <typename FP>
63 class Model : public FlowModel<FP, InfectionState, Populations<FP, AgeGroup, InfectionState>, Parameters<FP>, Flows>
64 {
66 
67 public:
68  using typename Base::ParameterSet;
69  using typename Base::Populations;
70 
71  Model(const Populations& pop, const ParameterSet& params)
72  : Base(pop, params)
73  {
74  }
75 
76  Model(int num_agegroups)
78  ParameterSet(mio::AgeGroup(num_agegroups)))
79  {
80  }
81 
82  void get_flows(Eigen::Ref<const Eigen::VectorX<FP>> pop, Eigen::Ref<const Eigen::VectorX<FP>> y, FP t,
83  Eigen::Ref<Eigen::VectorX<FP>> flows) const override
84  {
85  auto const& params = this->parameters;
86  AgeGroup n_agegroups = params.get_num_groups();
87 
88  ContactMatrixGroup<FP> const& contact_matrix = params.template get<ContactPatterns<FP>>();
89 
90  FP icu_occupancy = 0.0;
91  FP test_and_trace_required = 0.0;
92  for (AgeGroup i = 0; i < n_agegroups; ++i) {
93  test_and_trace_required += (1 - params.template get<RecoveredPerInfectedNoSymptoms<FP>>()[i]) /
94  params.template get<TimeInfectedNoSymptoms<FP>>()[i] *
95  this->populations.get_from(pop, {i, InfectionState::InfectedNoSymptoms});
96  icu_occupancy += this->populations.get_from(pop, {i, InfectionState::InfectedCritical});
97  }
98 
99  for (AgeGroup i = 0; i < n_agegroups; i++) {
100 
101  size_t Si = this->populations.get_flat_index({i, InfectionState::Susceptible});
102  size_t Ei = this->populations.get_flat_index({i, InfectionState::Exposed});
103  size_t INSi = this->populations.get_flat_index({i, InfectionState::InfectedNoSymptoms});
104  size_t INSCi = this->populations.get_flat_index({i, InfectionState::InfectedNoSymptomsConfirmed});
105  size_t ISyi = this->populations.get_flat_index({i, InfectionState::InfectedSymptoms});
106  size_t ISyCi = this->populations.get_flat_index({i, InfectionState::InfectedSymptomsConfirmed});
107  size_t ISevi = this->populations.get_flat_index({i, InfectionState::InfectedSevere});
108  size_t ICri = this->populations.get_flat_index({i, InfectionState::InfectedCritical});
109 
110  for (AgeGroup j = 0; j < n_agegroups; j++) {
111  size_t Sj = this->populations.get_flat_index({j, InfectionState::Susceptible});
112  size_t Ej = this->populations.get_flat_index({j, InfectionState::Exposed});
113  size_t INSj = this->populations.get_flat_index({j, InfectionState::InfectedNoSymptoms});
114  size_t ISyj = this->populations.get_flat_index({j, InfectionState::InfectedSymptoms});
115  size_t ISevj = this->populations.get_flat_index({j, InfectionState::InfectedSevere});
116  size_t ICrj = this->populations.get_flat_index({j, InfectionState::InfectedCritical});
117  size_t Rj = this->populations.get_flat_index({j, InfectionState::Recovered});
118 
119  //symptomatic are less well quarantined when testing and tracing is overwhelmed so they infect more people
120  FP riskFromInfectedSymptomatic =
121  smoother_cosine<FP>(test_and_trace_required, params.template get<TestAndTraceCapacity<FP>>(),
122  params.template get<TestAndTraceCapacity<FP>>() *
123  params.template get<TestAndTraceCapacityMaxRisk<FP>>(),
124  params.template get<RiskOfInfectionFromSymptomatic<FP>>()[j],
125  params.template get<MaxRiskOfInfectionFromSymptomatic<FP>>()[j]);
126 
127  // effective contact rate by contact rate between groups i and j and damping j
128  FP season_val = (1 + params.template get<Seasonality<FP>>() *
129  sin(std::numbers::pi_v<ScalarType> *
130  ((params.template get<StartDay<FP>>() + t) / 182.5 + 0.5)));
131  FP cont_freq_eff =
132  season_val * contact_matrix.get_matrix_at(SimulationTime<FP>(t))(
133  static_cast<Eigen::Index>((size_t)i), static_cast<Eigen::Index>((size_t)j));
134  FP Nj =
135  pop[Sj] + pop[Ej] + pop[INSj] + pop[ISyj] + pop[ISevj] + pop[ICrj] + pop[Rj]; // without died people
136  const FP divNj = (Nj < Limits<FP>::zero_tolerance()) ? FP(0.0) : FP(1.0 / Nj);
137  FP dummy_S = y[Si] * cont_freq_eff * divNj *
138  params.template get<TransmissionProbabilityOnContact<FP>>()[i] *
139  (params.template get<RelativeTransmissionNoSymptoms<FP>>()[j] * pop[INSj] +
140  riskFromInfectedSymptomatic * pop[ISyj]);
141 
142  // Susceptible -> Exposed
143  flows[this->template get_flat_flow_index<InfectionState::Susceptible, InfectionState::Exposed>({i})] +=
144  dummy_S;
145  }
146 
147  // ICU capacity shortage is close
148  FP criticalPerSevereAdjusted = smoother_cosine<FP>(
149  icu_occupancy, 0.90 * params.template get<ICUCapacity<FP>>(), params.template get<ICUCapacity<FP>>(),
150  params.template get<CriticalPerSevere<FP>>()[i], 0.0);
151 
152  FP deathsPerSevereAdjusted = params.template get<CriticalPerSevere<FP>>()[i] - criticalPerSevereAdjusted;
153 
154  // Exposed -> InfectedNoSymptoms
155  flows[this->template get_flat_flow_index<InfectionState::Exposed, InfectionState::InfectedNoSymptoms>(
156  {i})] = (1.0 / params.template get<TimeExposed<FP>>()[i]) * y[Ei];
157 
158  // InfectedNoSymptoms -> InfectedSymptoms / Recovered
161  (1 - params.template get<RecoveredPerInfectedNoSymptoms<FP>>()[i]) *
162  (1 / params.template get<TimeInfectedNoSymptoms<FP>>()[i]) * y[INSi];
163  flows[this->template get_flat_flow_index<InfectionState::InfectedNoSymptoms, InfectionState::Recovered>(
164  {i})] = params.template get<RecoveredPerInfectedNoSymptoms<FP>>()[i] *
165  (1.0 / params.template get<TimeInfectedNoSymptoms<FP>>()[i]) * y[INSi];
166 
167  // InfectedNoSymptomsConfirmed -> InfectedSymptomsConfirmed / Recovered
170  (1.0 - params.template get<RecoveredPerInfectedNoSymptoms<FP>>()[i]) *
171  (1.0 / params.template get<TimeInfectedNoSymptoms<FP>>()[i]) * y[INSCi];
174  params.template get<RecoveredPerInfectedNoSymptoms<FP>>()[i] *
175  (1.0 / params.template get<TimeInfectedNoSymptoms<FP>>()[i]) * y[INSCi];
176 
177  // InfectedSymptoms -> InfectedSevere / Recovered
178  flows[this->template get_flat_flow_index<InfectionState::InfectedSymptoms, InfectionState::InfectedSevere>(
179  {i})] = params.template get<SeverePerInfectedSymptoms<FP>>()[i] /
180  params.template get<TimeInfectedSymptoms<FP>>()[i] * y[ISyi];
181  flows[this->template get_flat_flow_index<InfectionState::InfectedSymptoms, InfectionState::Recovered>(
182  {i})] = (1 - params.template get<SeverePerInfectedSymptoms<FP>>()[i]) /
183  params.template get<TimeInfectedSymptoms<FP>>()[i] * y[ISyi];
184 
185  // InfectedSymptomsConfirmed -> InfectedSevere / Recovered
188  params.template get<SeverePerInfectedSymptoms<FP>>()[i] /
189  params.template get<TimeInfectedSymptoms<FP>>()[i] * y[ISyCi];
192  (1.0 - params.template get<SeverePerInfectedSymptoms<FP>>()[i]) /
193  params.template get<TimeInfectedSymptoms<FP>>()[i] * y[ISyCi];
194 
195  // InfectedSevere -> InfectedCritical / Recovered / Dead
196  flows[this->template get_flat_flow_index<InfectionState::InfectedSevere, InfectionState::InfectedCritical>(
197  {i})] = criticalPerSevereAdjusted / params.template get<TimeInfectedSevere<FP>>()[i] * y[ISevi];
198  flows[this->template get_flat_flow_index<InfectionState::InfectedSevere, InfectionState::Recovered>({i})] =
199  (1.0 - params.template get<CriticalPerSevere<FP>>()[i] -
200  params.template get<DeathsPerSevere<FP>>()[i]) /
201  params.template get<TimeInfectedSevere<FP>>()[i] * y[ISevi];
202  flows[this->template get_flat_flow_index<InfectionState::InfectedSevere, InfectionState::Dead>({i})] =
203  (params.template get<DeathsPerSevere<FP>>()[i] + deathsPerSevereAdjusted) /
204  params.template get<TimeInfectedSevere<FP>>()[i] * y[ISevi];
205 
206  // InfectedCritical -> Dead / Recovered
207  flows[this->template get_flat_flow_index<InfectionState::InfectedCritical, InfectionState::Dead>({i})] =
208  params.template get<DeathsPerCritical<FP>>()[i] / params.template get<TimeInfectedCritical<FP>>()[i] *
209  y[ICri];
210  flows[this->template get_flat_flow_index<InfectionState::InfectedCritical, InfectionState::Recovered>(
211  {i})] = (1.0 - params.template get<DeathsPerCritical<FP>>()[i]) /
212  params.template get<TimeInfectedCritical<FP>>()[i] * y[ICri];
213  }
214  }
215 
220  template <class IOContext>
221  void serialize(IOContext& io) const
222  {
223  auto obj = io.create_object("Model");
224  obj.add_element("Parameters", this->parameters);
225  obj.add_element("Populations", this->populations);
226  }
227 
232  template <class IOContext>
233  static IOResult<Model> deserialize(IOContext& io)
234  {
235  auto obj = io.expect_object("Model");
236  auto par = obj.expect_element("Parameters", Tag<ParameterSet>{});
237  auto pop = obj.expect_element("Populations", Tag<Populations>{});
238  return apply(
239  io,
240  [](auto&& par_, auto&& pop_) {
241  return Model{pop_, par_};
242  },
243  par, pop);
244  }
245 };
246 
247 //forward declaration, see below.
248 template <typename FP, class BaseT = mio::Simulation<FP, Model<FP>>>
249 class Simulation;
250 
259 template <typename FP, class Base = mio::Simulation<FP, Model<FP>>>
260 FP get_infections_relative(const Simulation<FP, Base>& model, FP t, const Eigen::Ref<const Eigen::VectorX<FP>>& y);
261 
267 template <typename FP, class BaseT>
268 class Simulation : public BaseT
269 {
270 public:
277  Simulation(mio::osecir::Model<FP> const& model, FP t0 = 0., FP dt = 0.1)
278  : BaseT(model, t0, dt)
279  , m_t_last_npi_check(t0)
280  {
281  }
282 
290  Eigen::Ref<Eigen::VectorX<FP>> advance(FP tmax)
291  {
292  using std::min;
293  auto& t_end_dyn_npis = this->get_model().parameters.get_end_dynamic_npis();
294  auto& dyn_npis = this->get_model().parameters.template get<DynamicNPIsInfectedSymptoms<FP>>();
295  auto& contact_patterns = this->get_model().parameters.template get<ContactPatterns<FP>>();
296 
297  FP delay_npi_implementation; // delay which is needed to implement a NPI that is criterion-dependent
298  FP t = BaseT::get_result().get_last_time();
299  const FP dt = dyn_npis.get_thresholds().size() > 0 ? FP(dyn_npis.get_interval().get()) : FP(tmax);
300 
301  while (t < tmax) {
302  FP dt_eff = min<FP>(dt, tmax - t);
303  dt_eff = min<FP>(dt_eff, m_t_last_npi_check + dt - t);
304 
305  BaseT::advance(t + dt_eff);
306  if (t > 0) {
307  delay_npi_implementation =
308  this->get_model().parameters.template get<DynamicNPIsImplementationDelay<FP>>();
309  }
310  else { // DynamicNPIs for t=0 are 'misused' to be from-start NPIs. I.e., do not enforce delay.
311  delay_npi_implementation = 0;
312  }
313  t = t + dt_eff;
314 
315  if (dyn_npis.get_thresholds().size() > 0) {
316  if (floating_point_greater_equal<FP>(t, m_t_last_npi_check + dt)) {
317  if (t < t_end_dyn_npis) {
318  auto inf_rel = get_infections_relative<FP>(*this, t, this->get_result().get_last_value()) *
319  dyn_npis.get_base_value();
320  auto exceeded_threshold = dyn_npis.get_max_exceeded_threshold(inf_rel);
321  if (exceeded_threshold != dyn_npis.get_thresholds().end() &&
322  (exceeded_threshold->first > m_dynamic_npi.first ||
323  t > FP(m_dynamic_npi.second))) { //old npi was weaker or is expired
324 
325  auto t_start = SimulationTime<FP>(t + delay_npi_implementation);
326  auto t_end = t_start + SimulationTime<FP>(FP(dyn_npis.get_duration()));
327  m_dynamic_npi = std::make_pair(exceeded_threshold->first, t_end);
328  implement_dynamic_npis<FP>(contact_patterns.get_cont_freq_mat(), exceeded_threshold->second,
329  t_start, t_end, [](auto& g) {
330  return make_contact_damping_matrix(g);
331  });
332  }
333  }
334  m_t_last_npi_check = t;
335  }
336  }
337  else {
338  m_t_last_npi_check = t;
339  }
340  }
341 
342  return this->get_result().get_last_value();
343  }
344 
345 private:
347  std::pair<FP, SimulationTime<FP>> m_dynamic_npi = {-std::numeric_limits<FP>::max(), mio::SimulationTime<FP>(0)};
348 };
349 
362 template <typename FP>
363 inline auto simulate(FP t0, FP tmax, FP dt, const Model<FP>& model,
364  std::unique_ptr<OdeIntegratorCore<FP>>&& integrator_core = nullptr)
365 {
366  return mio::simulate<FP, Model<FP>, Simulation<FP>>(t0, tmax, dt, model, std::move(integrator_core));
367 }
368 
381 template <typename FP>
382 inline auto simulate_flows(FP t0, FP tmax, FP dt, const Model<FP>& model,
383  std::unique_ptr<OdeIntegratorCore<FP>>&& integrator_core = nullptr)
384 {
385  return mio::simulate_flows<FP, Model<FP>, Simulation<FP, mio::FlowSimulation<FP, Model<FP>>>>(
386  t0, tmax, dt, model, std::move(integrator_core));
387 }
388 
389 //see declaration above.
390 template <typename FP, class Base>
391 FP get_infections_relative(const Simulation<FP, Base>& sim, FP /* t*/, const Eigen::Ref<const Eigen::VectorX<FP>>& y)
392 {
393  FP sum_inf = 0;
394  for (auto i = AgeGroup(0); i < sim.get_model().parameters.get_num_groups(); ++i) {
395  sum_inf += sim.get_model().populations.get_from(y, {i, InfectionState::InfectedSymptoms});
396  }
397  FP inf_rel = sum_inf / sim.get_model().populations.get_total();
398 
399  return inf_rel;
400 }
401 
410 template <typename FP, class Base>
412 {
413  using std::abs;
414  using std::sin;
415 
416  if (!(t_idx < static_cast<size_t>(sim.get_result().get_num_time_points()))) {
417  return mio::failure(mio::StatusCode::OutOfRange, "t_idx is not a valid index for the TimeSeries");
418  }
419 
420  auto const& params = sim.get_model().parameters;
421  const size_t num_groups = (size_t)sim.get_model().parameters.get_num_groups();
422  //The infected compartments in the SECIR Model are: Exposed, Carrier, Infected, Hospitalized, ICU in respective agegroups
423  const size_t num_infected_compartments = 5;
424  const size_t total_infected_compartments = num_infected_compartments * num_groups;
425  const FP pi = std::numbers::pi_v<ScalarType>;
426 
427  //F encodes new Infections and V encodes transition times in the next-generation matrix calculation of R_t
428  Eigen::MatrixX<FP> F(total_infected_compartments, total_infected_compartments);
429  Eigen::MatrixX<FP> V(total_infected_compartments, total_infected_compartments);
430  F = Eigen::MatrixX<FP>::Zero(total_infected_compartments,
431  total_infected_compartments); //Initialize matrices F and V with zeroes
432  V = Eigen::MatrixX<FP>::Zero(total_infected_compartments, total_infected_compartments);
433 
434  FP test_and_trace_required = 0.0;
435  FP icu_occupancy = 0.0;
436  for (auto i = AgeGroup(0); i < (mio::AgeGroup)num_groups; ++i) {
437  test_and_trace_required +=
438  (1 - params.template get<RecoveredPerInfectedNoSymptoms<FP>>()[i]) /
439  params.template get<TimeInfectedNoSymptoms<FP>>()[i] *
440  sim.get_result().get_value(
441  t_idx)[sim.get_model().populations.get_flat_index({i, InfectionState::InfectedNoSymptoms})];
442  icu_occupancy += sim.get_result().get_value(
443  t_idx)[sim.get_model().populations.get_flat_index({i, InfectionState::InfectedCritical})];
444  }
445 
446  FP season_val =
447  (1 +
448  params.template get<Seasonality<FP>>() *
449  sin(pi *
450  ((sim.get_model().parameters.template get<StartDay<FP>>() + sim.get_result().get_time(t_idx)) / 182.5 +
451  0.5)));
452  ContactMatrixGroup<FP> const& contact_matrix = sim.get_model().parameters.template get<ContactPatterns<FP>>();
453 
454  Eigen::MatrixX<FP> cont_freq_eff(num_groups, num_groups);
455  Eigen::MatrixX<FP> riskFromInfectedSymptomatic_derivatives(num_groups, num_groups);
456  Eigen::VectorX<FP> divN(num_groups);
457  Eigen::VectorX<FP> riskFromInfectedSymptomatic(num_groups);
458 
459  for (mio::AgeGroup k = 0; k < (mio::AgeGroup)num_groups; k++) {
460  FP temp = sim.get_result().get_value(
461  t_idx)[sim.get_model().populations.get_flat_index({k, InfectionState::Susceptible})] +
462  sim.get_result().get_value(
463  t_idx)[sim.get_model().populations.get_flat_index({k, InfectionState::Exposed})] +
464  sim.get_result().get_value(
465  t_idx)[sim.get_model().populations.get_flat_index({k, InfectionState::InfectedNoSymptoms})] +
466  sim.get_result().get_value(
467  t_idx)[sim.get_model().populations.get_flat_index({k, InfectionState::InfectedSymptoms})] +
468  sim.get_result().get_value(
469  t_idx)[sim.get_model().populations.get_flat_index({k, InfectionState::InfectedSevere})] +
470  sim.get_result().get_value(
471  t_idx)[sim.get_model().populations.get_flat_index({k, InfectionState::InfectedCritical})] +
472  sim.get_result().get_value(
473  t_idx)[sim.get_model().populations.get_flat_index({k, InfectionState::Recovered})];
474  if (temp < Limits<FP>::zero_tolerance()) {
475  temp = 1.0;
476  }
477  divN[(size_t)k] = 1 / temp;
478 
479  riskFromInfectedSymptomatic[(size_t)k] = smoother_cosine<FP>(
480  test_and_trace_required, params.template get<TestAndTraceCapacity<FP>>(),
481  (params.template get<TestAndTraceCapacity<FP>>()) * params.template get<TestAndTraceCapacityMaxRisk<FP>>(),
482  params.template get<RiskOfInfectionFromSymptomatic<FP>>()[k],
483  params.template get<MaxRiskOfInfectionFromSymptomatic<FP>>()[k]);
484 
485  for (mio::AgeGroup l = 0; l < (mio::AgeGroup)num_groups; l++) {
486  if (test_and_trace_required < params.template get<TestAndTraceCapacity<FP>>() ||
487  test_and_trace_required > params.template get<TestAndTraceCapacityMaxRisk<FP>>() *
488  params.template get<TestAndTraceCapacity<FP>>()) {
489  riskFromInfectedSymptomatic_derivatives((size_t)k, (size_t)l) = 0;
490  }
491  else {
492  riskFromInfectedSymptomatic_derivatives((size_t)k, (size_t)l) =
493  -0.5 * pi *
494  (params.template get<MaxRiskOfInfectionFromSymptomatic<FP>>()[k] -
495  params.template get<RiskOfInfectionFromSymptomatic<FP>>()[k]) /
496  (4 * params.template get<TestAndTraceCapacity<FP>>()) *
497  (1 - params.template get<RecoveredPerInfectedNoSymptoms<FP>>()[l]) /
498  params.template get<TimeInfectedNoSymptoms<FP>>()[l] *
499  sin(pi / (4 * params.template get<TestAndTraceCapacity<FP>>()) *
500  (test_and_trace_required - params.template get<TestAndTraceCapacity<FP>>()));
501  }
502  }
503 
504  for (Eigen::Index l = 0; l < (Eigen::Index)num_groups; l++) {
505  cont_freq_eff(l, (size_t)k) =
506  season_val * contact_matrix.get_matrix_at(SimulationTime<FP>(static_cast<FP>(t_idx)))(
507  static_cast<Eigen::Index>((size_t)l), static_cast<Eigen::Index>((size_t)k));
508  }
509  }
510 
511  //Check criterion if matrix V will be invertible by checking if subblock J is invertible
512  Eigen::MatrixX<FP> J(num_groups, num_groups);
513  J = Eigen::MatrixX<FP>::Zero(num_groups, num_groups);
514  for (size_t i = 0; i < num_groups; i++) {
515  J(i, i) = 1 / (params.template get<TimeInfectedCritical<FP>>()[(mio::AgeGroup)i]);
516 
517  if (!(icu_occupancy < 0.9 * params.template get<ICUCapacity<FP>>() ||
518  icu_occupancy > params.template get<ICUCapacity<FP>>())) {
519  for (size_t j = 0; j < num_groups; j++) {
520  J(i, j) -= sim.get_result().get_value(t_idx)[sim.get_model().populations.get_flat_index(
521  {(mio::AgeGroup)i, InfectionState::InfectedSevere})] /
522  params.template get<TimeInfectedSevere<FP>>()[(mio::AgeGroup)i] * 5 *
523  params.template get<CriticalPerSevere<FP>>()[(mio::AgeGroup)i] * pi /
524  (params.template get<ICUCapacity<FP>>()) *
525  sin(pi / (0.1 * params.template get<ICUCapacity<FP>>()) *
526  (icu_occupancy - 0.9 * params.template get<ICUCapacity<FP>>()));
527  }
528  }
529  }
530 
531  //Check, if J is invertible
532  if (J.determinant() == 0) {
533  return mio::failure(mio::StatusCode::UnknownError, "Matrix V is not invertible");
534  }
535 
536  //Initialize the matrix F
537  for (size_t i = 0; i < num_groups; i++) {
538  for (size_t j = 0; j < num_groups; j++) {
539 
540  FP temp = 0;
541  for (Eigen::Index k = 0; k < (Eigen::Index)num_groups; k++) {
542  temp += cont_freq_eff(i, k) *
543  sim.get_result().get_value(t_idx)[sim.get_model().populations.get_flat_index(
544  {(mio::AgeGroup)k, InfectionState::InfectedSymptoms})] *
545  riskFromInfectedSymptomatic_derivatives(k, j) * divN[k];
546  }
547 
548  F(i, j + num_groups) =
549  sim.get_result().get_value(t_idx)[sim.get_model().populations.get_flat_index(
550  {(mio::AgeGroup)i, InfectionState::Susceptible})] *
552  (cont_freq_eff(i, j) * params.template get<RelativeTransmissionNoSymptoms<FP>>()[(mio::AgeGroup)j] *
553  divN[(size_t)j] +
554  temp);
555  }
556 
557  for (size_t j = 0; j < num_groups; j++) {
558  F(i, j + 2 * num_groups) = sim.get_result().get_value(t_idx)[sim.get_model().populations.get_flat_index(
559  {(mio::AgeGroup)i, InfectionState::Susceptible})] *
561  cont_freq_eff(i, j) * riskFromInfectedSymptomatic[(size_t)j] * divN[(size_t)j];
562  }
563  }
564 
565  //Initialize the matrix V
566  for (Eigen::Index i = 0; i < (Eigen::Index)num_groups; i++) {
567  FP criticalPerSevereAdjusted = smoother_cosine<FP>(
568  icu_occupancy, 0.90 * params.template get<ICUCapacity<FP>>(), params.template get<ICUCapacity<FP>>(),
569  params.template get<CriticalPerSevere<FP>>()[(mio::AgeGroup)i], 0);
570 
571  V(i, i) = 1 / params.template get<TimeExposed<FP>>()[(mio::AgeGroup)i];
572  V(i + num_groups, i) = -1 / params.template get<TimeExposed<FP>>()[(mio::AgeGroup)i];
573  V(i + num_groups, i + num_groups) = 1 / params.template get<TimeInfectedNoSymptoms<FP>>()[(mio::AgeGroup)i];
574  V(i + 2 * num_groups, i + num_groups) =
575  -(1 - params.template get<RecoveredPerInfectedNoSymptoms<FP>>()[(mio::AgeGroup)i]) /
576  params.template get<TimeInfectedNoSymptoms<FP>>()[(mio::AgeGroup)i];
577  V(i + 2 * num_groups, i + 2 * num_groups) =
578  (1 / params.template get<TimeInfectedSymptoms<FP>>()[(mio::AgeGroup)i]);
579  V(i + 3 * num_groups, i + 2 * num_groups) =
580  -params.template get<SeverePerInfectedSymptoms<FP>>()[(mio::AgeGroup)i] /
581  params.template get<TimeInfectedSymptoms<FP>>()[(mio::AgeGroup)i];
582  V(i + 3 * num_groups, i + 3 * num_groups) =
583  1 / (params.template get<TimeInfectedSevere<FP>>()[(mio::AgeGroup)i]);
584  V(i + 4 * num_groups, i + 3 * num_groups) =
585  -criticalPerSevereAdjusted / (params.template get<TimeInfectedSevere<FP>>()[(mio::AgeGroup)i]);
586 
587  for (size_t j = 0; j < num_groups; j++) {
588  V(i + 4 * num_groups, j + 4 * num_groups) = J(i, j);
589  }
590  }
591 
592  V = V.inverse();
593 
594  Eigen::MatrixX<ScalarType> NextGenMatrix(num_infected_compartments * num_groups, 5 * num_groups);
595  NextGenMatrix.noalias() = F * V;
596 
597  // Compute the largest eigenvalue in absolute value
598  Eigen::ComplexEigenSolver<Eigen::MatrixX<ScalarType>> ces;
599  ces.compute(NextGenMatrix);
600  FP rho = ces.eigenvalues().cwiseAbs().maxCoeff();
601 
602  return mio::success(rho);
603 }
604 
612 template <typename FP, class Base>
613 Eigen::VectorX<FP> get_reproduction_numbers(const Simulation<FP, Base>& sim)
614 {
615  Eigen::VectorX<FP> temp(sim.get_result().get_num_time_points());
616  for (int i = 0; i < sim.get_result().get_num_time_points(); i++) {
617  temp[i] = get_reproduction_number<ScalarType>((size_t)i, sim).value();
618  }
619  return temp;
620 }
621 
629 template <typename FP, class Base>
631 {
632  if (t_value < sim.get_result().get_time(0) || t_value > sim.get_result().get_last_time()) {
634  "Cannot interpolate reproduction number outside computed horizon of the TimeSeries");
635  }
636 
637  if (t_value == sim.get_result().get_time(0)) {
638  return mio::success(get_reproduction_number((size_t)0, sim).value());
639  }
640 
641  const auto& times = sim.get_result().get_times();
642 
643  auto time_late = std::distance(times.begin(), std::lower_bound(times.begin(), times.end(), t_value));
644 
645  FP y1 = get_reproduction_number(static_cast<size_t>(time_late - 1), sim).value();
646  FP y2 = get_reproduction_number(static_cast<size_t>(time_late), sim).value();
647 
648  auto result = linear_interpolation(t_value, sim.get_result().get_time(time_late - 1),
649  sim.get_result().get_time(time_late), y1, y2);
650  return mio::success(static_cast<FP>(result));
651 }
652 
664 template <typename FP, class Base = mio::Simulation<FP, Model<FP>>>
665 auto get_mobility_factors(const Simulation<FP, Base>& sim, FP /*t*/, const Eigen::Ref<const Eigen::VectorX<FP>>& y)
666 {
667  auto& params = sim.get_model().parameters;
668  //parameters as arrays
669  auto&& p_asymp = params.template get<RecoveredPerInfectedNoSymptoms<FP>>().array().template cast<FP>();
670  auto&& p_inf = params.template get<RiskOfInfectionFromSymptomatic<FP>>().array().template cast<FP>();
671  auto&& p_inf_max = params.template get<MaxRiskOfInfectionFromSymptomatic<FP>>().array().template cast<FP>();
672  //slice of InfectedNoSymptoms
673  auto y_INS = slice(y, {Eigen::Index(InfectionState::InfectedNoSymptoms),
674  Eigen::Index(size_t(params.get_num_groups())), Eigen::Index(InfectionState::Count)});
675 
676  //compute isolation, same as infection risk from main model
677  auto test_and_trace_required =
678  ((1 - p_asymp) / params.template get<TimeInfectedNoSymptoms<FP>>().array().template cast<FP>() * y_INS.array())
679  .sum();
680  auto test_and_trace_capacity = FP(params.template get<TestAndTraceCapacity<FP>>());
681  auto test_and_trace_capacity_max_risk = FP(params.template get<TestAndTraceCapacityMaxRisk<FP>>());
682  auto riskFromInfectedSymptomatic = smoother_cosine<FP>(test_and_trace_required, test_and_trace_capacity,
683  test_and_trace_capacity * test_and_trace_capacity_max_risk,
684  p_inf.matrix(), p_inf_max.matrix());
685 
686  //set factor for infected
687  auto factors = Eigen::VectorX<FP>::Ones(y.rows()).eval();
688  slice(factors, {Eigen::Index(InfectionState::InfectedSymptoms), Eigen::Index(size_t(params.get_num_groups())),
689  Eigen::Index(InfectionState::Count)})
690  .array() = riskFromInfectedSymptomatic;
691  return factors;
692 }
693 
694 template <typename FP, class Base = mio::Simulation<FP, Model<FP>>>
695 auto test_commuters(Simulation<FP, Base>& sim, Eigen::Ref<Eigen::VectorX<FP>> mobile_population, FP time)
696 {
697  auto& model = sim.get_model();
698  FP nondetection = 1.0;
699  if (time >= model.parameters.get_start_commuter_detection() &&
700  time < model.parameters.get_end_commuter_detection()) {
701  nondetection = model.parameters.get_commuter_nondetection();
702  }
703  for (auto i = AgeGroup(0); i < model.parameters.get_num_groups(); ++i) {
704  auto INSi = model.populations.get_flat_index({i, InfectionState::InfectedNoSymptoms});
705  auto INSCi = model.populations.get_flat_index({i, InfectionState::InfectedNoSymptomsConfirmed});
706  auto ISyi = model.populations.get_flat_index({i, InfectionState::InfectedSymptoms});
707  auto ISyCi = model.populations.get_flat_index({i, InfectionState::InfectedSymptomsConfirmed});
708 
709  //put detected commuters in their own compartment so they don't contribute to infections in their home node
710  sim.get_result().get_last_value()[INSi] -= mobile_population[INSi] * (1 - nondetection);
711  sim.get_result().get_last_value()[INSCi] += mobile_population[INSi] * (1 - nondetection);
712  sim.get_result().get_last_value()[ISyi] -= mobile_population[ISyi] * (1 - nondetection);
713  sim.get_result().get_last_value()[ISyCi] += mobile_population[ISyi] * (1 - nondetection);
714 
715  //reduce the number of commuters
716  mobile_population[ISyi] *= nondetection;
717  mobile_population[INSi] *= nondetection;
718  }
719 }
720 
721 } // namespace osecir
722 } // namespace mio
723 
724 #endif // ODESECIR_MODEL_H
represents a collection of contact frequency matrices that whose sum is the total number of contacts.
Definition: contact_matrix.h:536
A FlowModel is a CompartmentalModel defined by the flows between compartments.
Definition: flow_model.h:74
constexpr size_t get_flat_flow_index() const
A flat index into an array of flows (as computed by get_flows), if the only used category in Pop is C...
Definition: flow_model.h:201
Interface class defining the integration step used in a SystemIntegrator.
Definition: integrator.h:48
Typesafe wrapper for a floating-point simulation time value (in days).
Definition: damping.h:78
Run the Simulation in discrete steps, evolve the Model and report results.
Definition: models/abm/simulation.h:37
Model & get_model()
Get the Model that this Simulation evolves.
Definition: models/abm/simulation.h:91
Definition: ode_secir/model.h:64
Model(const Populations &pop, const ParameterSet &params)
Definition: ode_secir/model.h:71
Model(int num_agegroups)
Definition: ode_secir/model.h:76
void get_flows(Eigen::Ref< const Eigen::VectorX< FP >> pop, Eigen::Ref< const Eigen::VectorX< FP >> y, FP t, Eigen::Ref< Eigen::VectorX< FP >> flows) const override
Definition: ode_secir/model.h:82
static IOResult< Model > deserialize(IOContext &io)
deserialize an object of this class.
Definition: ode_secir/model.h:233
void serialize(IOContext &io) const
serialize this.
Definition: ode_secir/model.h:221
Parameters of an age-resolved SECIR/SECIHURD model.
Definition: ode_secir/parameters.h:417
specialization of compartment model simulation for secir models.
Definition: ode_secir/model.h:269
Eigen::Ref< Eigen::VectorX< FP > > advance(FP tmax)
advance simulation to tmax.
Definition: ode_secir/model.h:290
Simulation(mio::osecir::Model< FP > const &model, FP t0=0., FP dt=0.1)
construct a simulation.
Definition: ode_secir/model.h:277
FP m_t_last_npi_check
Definition: ode_secir/model.h:346
ad::internal::unary_intermediate< AD_TAPE_REAL, ad::internal::active_type< AD_TAPE_REAL, DATA_HANDLER_1 >, ad::operations::ad_sin< AD_TAPE_REAL > > sin(const ad::internal::active_type< AD_TAPE_REAL, DATA_HANDLER_1 > &x1)
Definition: ad.hpp:913
static min_max_return_type< ad::internal::active_type< AD_TAPE_REAL, DATA_HANDLER_1 >, ad::internal::active_type< AD_TAPE_REAL, DATA_HANDLER_1 > >::type min(const ad::internal::active_type< AD_TAPE_REAL, DATA_HANDLER_1 > &a, const ad::internal::active_type< AD_TAPE_REAL, DATA_HANDLER_1 > &b)
Definition: ad.hpp:2599
ad::internal::unary_intermediate< AD_TAPE_REAL, ad::internal::active_type< AD_TAPE_REAL, DATA_HANDLER_1 >, ad::operations::ad_fabs< AD_TAPE_REAL > > abs(const ad::internal::active_type< AD_TAPE_REAL, DATA_HANDLER_1 > &x1)
Definition: ad.hpp:1144
static min_max_return_type< ad::internal::active_type< AD_TAPE_REAL, DATA_HANDLER_1 >, ad::internal::active_type< AD_TAPE_REAL, DATA_HANDLER_1 > >::type max(const ad::internal::active_type< AD_TAPE_REAL, DATA_HANDLER_1 > &a, const ad::internal::active_type< AD_TAPE_REAL, DATA_HANDLER_1 > &b)
Definition: ad.hpp:2596
trait_value< T >::RETURN_TYPE & value(T &x)
Definition: ad.hpp:3308
Eigen::VectorX< FP > get_reproduction_numbers(const Simulation< FP, Base > &sim)
Computes the reproduction number for all time points of the Model output obtained by the Simulation.
Definition: ode_secir/model.h:613
IOResult< FP > get_reproduction_number(FP t_value, const Simulation< FP, Base > &sim)
Computes the reproduction number at a given time point of the Simulation.
Definition: ode_secir/model.h:630
A collection of classes to simplify handling of matrix shapes in meta programming.
Definition: models/abm/analyze_result.h:30
std::vector< TimeSeries< FP > > simulate_flows(FP t0, FP tmax, FP dt, Model const &model, std::unique_ptr< OdeIntegratorCore< FP >> &&integrator_core=nullptr)
Run a FlowSimulation of a FlowModel.
Definition: flow_simulation.h:111
FP get_infections_relative(const SimulationNode< FP, Sim > &node, FP t, const Eigen::Ref< const Eigen::VectorX< FP >> &y)
get the percantage of infected people of the total population in the node If dynamic NPIs are enabled...
Definition: metapopulation_mobility_instant.h:484
auto failure(const IOStatus &s)
Create an object that is implicitly convertible to an error IOResult<T>.
Definition: io.h:380
boost::outcome_v2::in_place_type_t< T > Tag
Type that is used for overload resolution.
Definition: io.h:407
auto i
Definition: io.h:809
details::ApplyResultT< F, T... > apply(IOContext &io, F f, const IOResult< T > &... rs)
Evaluate a function with zero or more unpacked IOResults as arguments.
Definition: io.h:481
auto linear_interpolation(const X &x_eval, const X &x_1, const X &x_2, const V &y1, const V &y2)
Linear interpolation between two data values.
Definition: interpolation.h:44
auto slice(V &&v, Seq< Eigen::Index > elems)
take a regular slice of a row or column vector.
Definition: eigen_util.h:107
auto success()
Create an object that is implicitly convertible to a succesful IOResult<void>.
Definition: io.h:359
TimeSeries< FP > simulate(FP t0, FP tmax, FP dt, Model const &model, std::unique_ptr< OdeIntegratorCore< FP >> &&integrator_core=nullptr)
Run a Simulation of a CompartmentalModel.
Definition: memilio/compartments/simulation.h:100
auto get_mobility_factors(const SimulationNode< FP, Sim > &node, FP t, const Eigen::Ref< const Eigen::VectorX< FP >> &y)
Get an additional mobility factor.
Definition: metapopulation_mobility_instant.h:509
constexpr std::tuple_element< I, std::tuple< Index< CategoryTags >... > >::type & get(Index< CategoryTags... > &i) noexcept
Retrieves the Index (by reference) at the Ith position of a MultiIndex.
Definition: index.h:304
boost::outcome_v2::unchecked< T, IOStatus > IOResult
Value-or-error type for operations that return a value but can fail.
Definition: io.h:353
void test_commuters(SimulationNode< FP, Sim > &node, Eigen::Ref< Eigen::VectorX< FP >> mobile_population, FP time)
Test persons when moving from their source node.
Definition: metapopulation_mobility_instant.h:533
Typesafe index representing an age group.
Definition: age_group.h:40
Populations populations
Definition: compartmental_model.h:156
mio::Populations< FP, AgeGroup, InfectionState > Populations
Definition: compartmental_model.h:62
ParameterSet parameters
Definition: compartmental_model.h:157
A Flow defines a possible transition between two Compartments in a FlowModel.
Definition: flow.h:33
Type specific limits for floating-points.
Definition: memilio/config.h:59
Collection of types. Each type is mapped to an index of type size_t.
Definition: type_list.h:32
the percentage of ICU patients per hospitalized patients in the SECIR model
Definition: ode_secir/parameters.h:276
the icu capacity in the SECIR model
Definition: ode_secir/parameters.h:80
the percentage of asymptomatic cases in the SECIR model
Definition: ode_secir/parameters.h:212
the risk of infection from symptomatic cases in the SECIR model
Definition: ode_secir/parameters.h:228
Multiplier for the test and trace capacity to determine when it is considered overloaded.
Definition: ode_secir/parameters.h:390
capacity to test and trace contacts of infected for quarantine per day.
Definition: ode_secir/parameters.h:374
the (mean) time in day unit for asymptomatic cases that are infectious but have not yet developed sym...
Definition: ode_secir/parameters.h:113
the time people are 'simply' hospitalized before returning home in the SECIR model in day unit
Definition: ode_secir/parameters.h:147
probability of getting infected from a contact
Definition: ode_secir/parameters.h:180