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  {
280  }
281 
289  Eigen::Ref<Eigen::VectorX<FP>> advance(FP tmax)
290  {
291  using std::min;
292  auto& dyn_npis = this->get_model().parameters.template get<DynamicNPIsInfectedSymptoms<FP>>();
293 
294  // If no dynamic NPI thresholds are set, directly use the base advance
295  if (dyn_npis.get_thresholds().size() == 0) {
296  return BaseT::advance(tmax);
297  }
298 
299  auto& contact_patterns = this->get_model().parameters.template get<ContactPatterns<FP>>();
300 
301  FP delay_npi_implementation;
302  FP t = BaseT::get_result().get_last_time();
303  while (t < tmax) {
304  if (t > 0) {
305  delay_npi_implementation = FP(dyn_npis.get_implementation_delay());
306  }
307  else { // DynamicNPIs for t=0 are treated as 'from-start NPIs'. I.e., do not enforce delay.
308  delay_npi_implementation = 0;
309  }
310 
311  if (dyn_npis.get_thresholds().size() > 0) {
312  FP direc_begin = FP(dyn_npis.get_directive_begin());
313  FP direc_end = FP(dyn_npis.get_directive_end());
314  if (floating_point_greater_equal(t, direc_begin, 1e-10) && t < direc_end) {
315  auto inf_rel = get_infections_relative<FP>(*this, t, this->get_result().get_last_value()) *
316  dyn_npis.get_base_value();
317  auto exceeded_threshold = dyn_npis.get_max_exceeded_threshold(inf_rel);
318  const bool npi_expired = floating_point_greater_equal(t, FP(m_dynamic_npi.second), 1e-10);
319  if (exceeded_threshold != dyn_npis.get_thresholds().end() &&
320  (exceeded_threshold->first > m_dynamic_npi.first || npi_expired)) {
321  // old npi was weaker or is expired
322 
323  // Keep-alive: if the NPI expired but the threshold is still exceeded at the
324  // same level, renew immediately without delay to avoid a gap.
325  // Apply implementation delay only if stronger NPI needed.
326  bool is_expiry_renewal =
327  npi_expired && !(exceeded_threshold->first > m_dynamic_npi.first);
328  FP effective_delay = is_expiry_renewal ? FP(0) : delay_npi_implementation;
329  if (t + effective_delay < direc_end) {
330  auto t_start = SimulationTime<FP>(t + effective_delay);
331  // set the end to the minimum of start+duration and the end of the directive
332  auto t_end = SimulationTime<FP>(min<FP>(direc_end, FP(t_start + dyn_npis.get_duration())));
333  m_dynamic_npi = std::make_pair(exceeded_threshold->first, t_end);
334  // For new NPIs (t_start > 0): shift t_start_damping by +1 so the smooth
335  // transition window [t_start, t_start+1] lies in the future, consistent with
336  // predefined dampings. For t_start = 0 (global t0): the window [-1, 0]
337  // is in the past and no shift is needed.
338  // For keep-alive renewals (is_expiry_renewal): do not shift t_start_damping.
339  // Since t_start == t_end_damping_old, the new start damping has the same
340  // (time, level, type) as the previous entry.
341  // Therefore, the contact matrix stays constant at the NPI level with no dip.
342  auto t_start_damping = (!is_expiry_renewal && FP(t_start) > FP(0))
343  ? SimulationTime<FP>(FP(t_start) + FP(1))
344  : t_start;
345  auto t_end_damping = SimulationTime<FP>(FP(t_end) + FP(1));
346  implement_dynamic_npis<FP>(contact_patterns.get_cont_freq_mat(), exceeded_threshold->second,
347  t_start_damping, t_end_damping, [](auto& g) {
348  return make_contact_damping_matrix(g);
349  });
350  }
351  }
352  }
353  }
354 
355  FP dt_eff = min<FP>(1.0, tmax - t);
356  // Clamp step size at the NPI end time to avoid stepping over the lifting phase.
357  // This ensures the keep-alive check activates exactly at t_end, preventing a brief dip.
358  FP npi_end = FP(m_dynamic_npi.second);
359  if (t < npi_end) {
360  dt_eff = min<FP>(dt_eff, npi_end - t);
361  }
362  BaseT::advance(t + dt_eff);
363  t = t + dt_eff;
364  }
365 
366  return this->get_result().get_last_value();
367  }
368 
369 private:
370  std::pair<FP, SimulationTime<FP>> m_dynamic_npi = {-std::numeric_limits<FP>::max(), mio::SimulationTime<FP>(0)};
371 };
372 
385 template <typename FP>
386 inline auto simulate(FP t0, FP tmax, FP dt, const Model<FP>& model,
387  std::unique_ptr<OdeIntegratorCore<FP>>&& integrator_core = nullptr)
388 {
389  return mio::simulate<FP, Model<FP>, Simulation<FP>>(t0, tmax, dt, model, std::move(integrator_core));
390 }
391 
404 template <typename FP>
405 inline auto simulate_flows(FP t0, FP tmax, FP dt, const Model<FP>& model,
406  std::unique_ptr<OdeIntegratorCore<FP>>&& integrator_core = nullptr)
407 {
408  return mio::simulate_flows<FP, Model<FP>, Simulation<FP, mio::FlowSimulation<FP, Model<FP>>>>(
409  t0, tmax, dt, model, std::move(integrator_core));
410 }
411 
412 //see declaration above.
413 template <typename FP, class Base>
414 FP get_infections_relative(const Simulation<FP, Base>& sim, FP /* t*/, const Eigen::Ref<const Eigen::VectorX<FP>>& y)
415 {
416  FP sum_inf = 0;
417  for (auto i = AgeGroup(0); i < sim.get_model().parameters.get_num_groups(); ++i) {
418  sum_inf += sim.get_model().populations.get_from(y, {i, InfectionState::InfectedSymptoms});
419  }
420  FP inf_rel = sum_inf / sim.get_model().populations.get_total();
421 
422  return inf_rel;
423 }
424 
433 template <typename FP, class Base>
435 {
436  using std::abs;
437  using std::sin;
438 
439  if (!(t_idx < static_cast<size_t>(sim.get_result().get_num_time_points()))) {
440  return mio::failure(mio::StatusCode::OutOfRange, "t_idx is not a valid index for the TimeSeries");
441  }
442 
443  auto const& params = sim.get_model().parameters;
444  const size_t num_groups = (size_t)sim.get_model().parameters.get_num_groups();
445  //The infected compartments in the SECIR Model are: Exposed, Carrier, Infected, Hospitalized, ICU in respective agegroups
446  const size_t num_infected_compartments = 5;
447  const size_t total_infected_compartments = num_infected_compartments * num_groups;
448  const FP pi = std::numbers::pi_v<ScalarType>;
449 
450  //F encodes new Infections and V encodes transition times in the next-generation matrix calculation of R_t
451  Eigen::MatrixX<FP> F(total_infected_compartments, total_infected_compartments);
452  Eigen::MatrixX<FP> V(total_infected_compartments, total_infected_compartments);
453  F = Eigen::MatrixX<FP>::Zero(total_infected_compartments,
454  total_infected_compartments); //Initialize matrices F and V with zeroes
455  V = Eigen::MatrixX<FP>::Zero(total_infected_compartments, total_infected_compartments);
456 
457  FP test_and_trace_required = 0.0;
458  FP icu_occupancy = 0.0;
459  for (auto i = AgeGroup(0); i < (mio::AgeGroup)num_groups; ++i) {
460  test_and_trace_required +=
461  (1 - params.template get<RecoveredPerInfectedNoSymptoms<FP>>()[i]) /
462  params.template get<TimeInfectedNoSymptoms<FP>>()[i] *
463  sim.get_result().get_value(
464  t_idx)[sim.get_model().populations.get_flat_index({i, InfectionState::InfectedNoSymptoms})];
465  icu_occupancy += sim.get_result().get_value(
466  t_idx)[sim.get_model().populations.get_flat_index({i, InfectionState::InfectedCritical})];
467  }
468 
469  FP season_val =
470  (1 +
471  params.template get<Seasonality<FP>>() *
472  sin(pi *
473  ((sim.get_model().parameters.template get<StartDay<FP>>() + sim.get_result().get_time(t_idx)) / 182.5 +
474  0.5)));
475  ContactMatrixGroup<FP> const& contact_matrix = sim.get_model().parameters.template get<ContactPatterns<FP>>();
476 
477  Eigen::MatrixX<FP> cont_freq_eff(num_groups, num_groups);
478  Eigen::MatrixX<FP> riskFromInfectedSymptomatic_derivatives(num_groups, num_groups);
479  Eigen::VectorX<FP> divN(num_groups);
480  Eigen::VectorX<FP> riskFromInfectedSymptomatic(num_groups);
481 
482  for (mio::AgeGroup k = 0; k < (mio::AgeGroup)num_groups; k++) {
483  FP temp = sim.get_result().get_value(
484  t_idx)[sim.get_model().populations.get_flat_index({k, InfectionState::Susceptible})] +
485  sim.get_result().get_value(
486  t_idx)[sim.get_model().populations.get_flat_index({k, InfectionState::Exposed})] +
487  sim.get_result().get_value(
488  t_idx)[sim.get_model().populations.get_flat_index({k, InfectionState::InfectedNoSymptoms})] +
489  sim.get_result().get_value(
490  t_idx)[sim.get_model().populations.get_flat_index({k, InfectionState::InfectedSymptoms})] +
491  sim.get_result().get_value(
492  t_idx)[sim.get_model().populations.get_flat_index({k, InfectionState::InfectedSevere})] +
493  sim.get_result().get_value(
494  t_idx)[sim.get_model().populations.get_flat_index({k, InfectionState::InfectedCritical})] +
495  sim.get_result().get_value(
496  t_idx)[sim.get_model().populations.get_flat_index({k, InfectionState::Recovered})];
497  if (temp < Limits<FP>::zero_tolerance()) {
498  temp = 1.0;
499  }
500  divN[(size_t)k] = 1 / temp;
501 
502  riskFromInfectedSymptomatic[(size_t)k] = smoother_cosine<FP>(
503  test_and_trace_required, params.template get<TestAndTraceCapacity<FP>>(),
504  (params.template get<TestAndTraceCapacity<FP>>()) * params.template get<TestAndTraceCapacityMaxRisk<FP>>(),
505  params.template get<RiskOfInfectionFromSymptomatic<FP>>()[k],
506  params.template get<MaxRiskOfInfectionFromSymptomatic<FP>>()[k]);
507 
508  for (mio::AgeGroup l = 0; l < (mio::AgeGroup)num_groups; l++) {
509  if (test_and_trace_required < params.template get<TestAndTraceCapacity<FP>>() ||
510  test_and_trace_required > params.template get<TestAndTraceCapacityMaxRisk<FP>>() *
511  params.template get<TestAndTraceCapacity<FP>>()) {
512  riskFromInfectedSymptomatic_derivatives((size_t)k, (size_t)l) = 0;
513  }
514  else {
515  riskFromInfectedSymptomatic_derivatives((size_t)k, (size_t)l) =
516  -0.5 * pi *
517  (params.template get<MaxRiskOfInfectionFromSymptomatic<FP>>()[k] -
518  params.template get<RiskOfInfectionFromSymptomatic<FP>>()[k]) /
519  (4 * params.template get<TestAndTraceCapacity<FP>>()) *
520  (1 - params.template get<RecoveredPerInfectedNoSymptoms<FP>>()[l]) /
521  params.template get<TimeInfectedNoSymptoms<FP>>()[l] *
522  sin(pi / (4 * params.template get<TestAndTraceCapacity<FP>>()) *
523  (test_and_trace_required - params.template get<TestAndTraceCapacity<FP>>()));
524  }
525  }
526 
527  for (Eigen::Index l = 0; l < (Eigen::Index)num_groups; l++) {
528  cont_freq_eff(l, (size_t)k) =
529  season_val * contact_matrix.get_matrix_at(SimulationTime<FP>(static_cast<FP>(t_idx)))(
530  static_cast<Eigen::Index>((size_t)l), static_cast<Eigen::Index>((size_t)k));
531  }
532  }
533 
534  //Check criterion if matrix V will be invertible by checking if subblock J is invertible
535  Eigen::MatrixX<FP> J(num_groups, num_groups);
536  J = Eigen::MatrixX<FP>::Zero(num_groups, num_groups);
537  for (size_t i = 0; i < num_groups; i++) {
538  J(i, i) = 1 / (params.template get<TimeInfectedCritical<FP>>()[(mio::AgeGroup)i]);
539 
540  if (!(icu_occupancy < 0.9 * params.template get<ICUCapacity<FP>>() ||
541  icu_occupancy > params.template get<ICUCapacity<FP>>())) {
542  for (size_t j = 0; j < num_groups; j++) {
543  J(i, j) -= sim.get_result().get_value(t_idx)[sim.get_model().populations.get_flat_index(
544  {(mio::AgeGroup)i, InfectionState::InfectedSevere})] /
545  params.template get<TimeInfectedSevere<FP>>()[(mio::AgeGroup)i] * 5 *
546  params.template get<CriticalPerSevere<FP>>()[(mio::AgeGroup)i] * pi /
547  (params.template get<ICUCapacity<FP>>()) *
548  sin(pi / (0.1 * params.template get<ICUCapacity<FP>>()) *
549  (icu_occupancy - 0.9 * params.template get<ICUCapacity<FP>>()));
550  }
551  }
552  }
553 
554  //Check, if J is invertible
555  if (J.determinant() == 0) {
556  return mio::failure(mio::StatusCode::UnknownError, "Matrix V is not invertible");
557  }
558 
559  //Initialize the matrix F
560  for (size_t i = 0; i < num_groups; i++) {
561  for (size_t j = 0; j < num_groups; j++) {
562 
563  FP temp = 0;
564  for (Eigen::Index k = 0; k < (Eigen::Index)num_groups; k++) {
565  temp += cont_freq_eff(i, k) *
566  sim.get_result().get_value(t_idx)[sim.get_model().populations.get_flat_index(
567  {(mio::AgeGroup)k, InfectionState::InfectedSymptoms})] *
568  riskFromInfectedSymptomatic_derivatives(k, j) * divN[k];
569  }
570 
571  F(i, j + num_groups) =
572  sim.get_result().get_value(t_idx)[sim.get_model().populations.get_flat_index(
573  {(mio::AgeGroup)i, InfectionState::Susceptible})] *
575  (cont_freq_eff(i, j) * params.template get<RelativeTransmissionNoSymptoms<FP>>()[(mio::AgeGroup)j] *
576  divN[(size_t)j] +
577  temp);
578  }
579 
580  for (size_t j = 0; j < num_groups; j++) {
581  F(i, j + 2 * num_groups) = sim.get_result().get_value(t_idx)[sim.get_model().populations.get_flat_index(
582  {(mio::AgeGroup)i, InfectionState::Susceptible})] *
584  cont_freq_eff(i, j) * riskFromInfectedSymptomatic[(size_t)j] * divN[(size_t)j];
585  }
586  }
587 
588  //Initialize the matrix V
589  for (Eigen::Index i = 0; i < (Eigen::Index)num_groups; i++) {
590  FP criticalPerSevereAdjusted = smoother_cosine<FP>(
591  icu_occupancy, 0.90 * params.template get<ICUCapacity<FP>>(), params.template get<ICUCapacity<FP>>(),
592  params.template get<CriticalPerSevere<FP>>()[(mio::AgeGroup)i], 0);
593 
594  V(i, i) = 1 / params.template get<TimeExposed<FP>>()[(mio::AgeGroup)i];
595  V(i + num_groups, i) = -1 / params.template get<TimeExposed<FP>>()[(mio::AgeGroup)i];
596  V(i + num_groups, i + num_groups) = 1 / params.template get<TimeInfectedNoSymptoms<FP>>()[(mio::AgeGroup)i];
597  V(i + 2 * num_groups, i + num_groups) =
598  -(1 - params.template get<RecoveredPerInfectedNoSymptoms<FP>>()[(mio::AgeGroup)i]) /
599  params.template get<TimeInfectedNoSymptoms<FP>>()[(mio::AgeGroup)i];
600  V(i + 2 * num_groups, i + 2 * num_groups) =
601  (1 / params.template get<TimeInfectedSymptoms<FP>>()[(mio::AgeGroup)i]);
602  V(i + 3 * num_groups, i + 2 * num_groups) =
603  -params.template get<SeverePerInfectedSymptoms<FP>>()[(mio::AgeGroup)i] /
604  params.template get<TimeInfectedSymptoms<FP>>()[(mio::AgeGroup)i];
605  V(i + 3 * num_groups, i + 3 * num_groups) =
606  1 / (params.template get<TimeInfectedSevere<FP>>()[(mio::AgeGroup)i]);
607  V(i + 4 * num_groups, i + 3 * num_groups) =
608  -criticalPerSevereAdjusted / (params.template get<TimeInfectedSevere<FP>>()[(mio::AgeGroup)i]);
609 
610  for (size_t j = 0; j < num_groups; j++) {
611  V(i + 4 * num_groups, j + 4 * num_groups) = J(i, j);
612  }
613  }
614 
615  V = V.inverse();
616 
617  Eigen::MatrixX<ScalarType> NextGenMatrix(num_infected_compartments * num_groups, 5 * num_groups);
618  NextGenMatrix.noalias() = F * V;
619 
620  // Compute the largest eigenvalue in absolute value
621  Eigen::ComplexEigenSolver<Eigen::MatrixX<ScalarType>> ces;
622  ces.compute(NextGenMatrix);
623  FP rho = ces.eigenvalues().cwiseAbs().maxCoeff();
624 
625  return mio::success(rho);
626 }
627 
635 template <typename FP, class Base>
636 Eigen::VectorX<FP> get_reproduction_numbers(const Simulation<FP, Base>& sim)
637 {
638  Eigen::VectorX<FP> temp(sim.get_result().get_num_time_points());
639  for (int i = 0; i < sim.get_result().get_num_time_points(); i++) {
640  temp[i] = get_reproduction_number<ScalarType>((size_t)i, sim).value();
641  }
642  return temp;
643 }
644 
652 template <typename FP, class Base>
654 {
655  if (t_value < sim.get_result().get_time(0) || t_value > sim.get_result().get_last_time()) {
657  "Cannot interpolate reproduction number outside computed horizon of the TimeSeries");
658  }
659 
660  if (t_value == sim.get_result().get_time(0)) {
661  return mio::success(get_reproduction_number((size_t)0, sim).value());
662  }
663 
664  const auto& times = sim.get_result().get_times();
665 
666  auto time_late = std::distance(times.begin(), std::lower_bound(times.begin(), times.end(), t_value));
667 
668  FP y1 = get_reproduction_number(static_cast<size_t>(time_late - 1), sim).value();
669  FP y2 = get_reproduction_number(static_cast<size_t>(time_late), sim).value();
670 
671  auto result = linear_interpolation(t_value, sim.get_result().get_time(time_late - 1),
672  sim.get_result().get_time(time_late), y1, y2);
673  return mio::success(static_cast<FP>(result));
674 }
675 
687 template <typename FP, class Base = mio::Simulation<FP, Model<FP>>>
688 auto get_mobility_factors(const Simulation<FP, Base>& sim, FP /*t*/, const Eigen::Ref<const Eigen::VectorX<FP>>& y)
689 {
690  auto& params = sim.get_model().parameters;
691  //parameters as arrays
692  auto&& p_asymp = params.template get<RecoveredPerInfectedNoSymptoms<FP>>().array().template cast<FP>();
693  auto&& p_inf = params.template get<RiskOfInfectionFromSymptomatic<FP>>().array().template cast<FP>();
694  auto&& p_inf_max = params.template get<MaxRiskOfInfectionFromSymptomatic<FP>>().array().template cast<FP>();
695  //slice of InfectedNoSymptoms
696  auto y_INS = slice(y, {Eigen::Index(InfectionState::InfectedNoSymptoms),
697  Eigen::Index(size_t(params.get_num_groups())), Eigen::Index(InfectionState::Count)});
698 
699  //compute isolation, same as infection risk from main model
700  auto test_and_trace_required =
701  ((1 - p_asymp) / params.template get<TimeInfectedNoSymptoms<FP>>().array().template cast<FP>() * y_INS.array())
702  .sum();
703  auto test_and_trace_capacity = FP(params.template get<TestAndTraceCapacity<FP>>());
704  auto test_and_trace_capacity_max_risk = FP(params.template get<TestAndTraceCapacityMaxRisk<FP>>());
705  auto riskFromInfectedSymptomatic = smoother_cosine<FP>(test_and_trace_required, test_and_trace_capacity,
706  test_and_trace_capacity * test_and_trace_capacity_max_risk,
707  p_inf.matrix(), p_inf_max.matrix());
708 
709  //set factor for infected
710  auto factors = Eigen::VectorX<FP>::Ones(y.rows()).eval();
711  slice(factors, {Eigen::Index(InfectionState::InfectedSymptoms), Eigen::Index(size_t(params.get_num_groups())),
712  Eigen::Index(InfectionState::Count)})
713  .array() = riskFromInfectedSymptomatic;
714  return factors;
715 }
716 
717 template <typename FP, class Base = mio::Simulation<FP, Model<FP>>>
718 auto test_commuters(Simulation<FP, Base>& sim, Eigen::Ref<Eigen::VectorX<FP>> mobile_population, FP time)
719 {
720  auto& model = sim.get_model();
721  FP nondetection = 1.0;
722  if (time >= model.parameters.get_start_commuter_detection() &&
723  time < model.parameters.get_end_commuter_detection()) {
724  nondetection = model.parameters.get_commuter_nondetection();
725  }
726  for (auto i = AgeGroup(0); i < model.parameters.get_num_groups(); ++i) {
727  auto INSi = model.populations.get_flat_index({i, InfectionState::InfectedNoSymptoms});
728  auto INSCi = model.populations.get_flat_index({i, InfectionState::InfectedNoSymptomsConfirmed});
729  auto ISyi = model.populations.get_flat_index({i, InfectionState::InfectedSymptoms});
730  auto ISyCi = model.populations.get_flat_index({i, InfectionState::InfectedSymptomsConfirmed});
731 
732  //put detected commuters in their own compartment so they don't contribute to infections in their home node
733  sim.get_result().get_last_value()[INSi] -= mobile_population[INSi] * (1 - nondetection);
734  sim.get_result().get_last_value()[INSCi] += mobile_population[INSi] * (1 - nondetection);
735  sim.get_result().get_last_value()[ISyi] -= mobile_population[ISyi] * (1 - nondetection);
736  sim.get_result().get_last_value()[ISyCi] += mobile_population[ISyi] * (1 - nondetection);
737 
738  //reduce the number of commuters
739  mobile_population[ISyi] *= nondetection;
740  mobile_population[INSi] *= nondetection;
741  }
742 }
743 
744 } // namespace osecir
745 } // namespace mio
746 
747 #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
A class for simulating a CompartmentalModel.
Definition: memilio/compartments/simulation.h:43
const Model & get_model() const
Get a reference to the model owned and used by the simulation.
Definition: simulation_base.h:131
TimeSeries< FP > & get_result()
Returns the simulation result describing the model population in each time step.
Definition: simulation_base.h:116
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:401
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:289
Simulation(mio::osecir::Model< FP > const &model, FP t0=0., FP dt=0.1)
construct a simulation.
Definition: ode_secir/model.h:277
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:636
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:653
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 percentage of infected people of the total population in the node.
Definition: metapopulation_mobility_instant.h:483
auto failure(const IOStatus &s)
Create an object that is implicitly convertible to an error IOResult<T>.
Definition: io.h:381
boost::outcome_v2::in_place_type_t< T > Tag
Type that is used for overload resolution.
Definition: io.h:408
auto i
Definition: io.h:810
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:482
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:360
bool floating_point_greater_equal(FP v1, FP v2, FP abs_tol=0, FP rel_tol=std::numeric_limits< FP >::min())
compare two floating point values with tolerances.
Definition: floating_point.h:136
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:508
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:354
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:532
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:374
capacity to test and trace contacts of infected for quarantine per day.
Definition: ode_secir/parameters.h:358
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