model.h Source File
|
CPP API
|
ode_secir/model.h
Go to the documentation of this file.
63 class Model : public FlowModel<FP, InfectionState, Populations<FP, AgeGroup, InfectionState>, Parameters<FP>, Flows>
65 using Base = FlowModel<FP, InfectionState, mio::Populations<FP, AgeGroup, InfectionState>, Parameters<FP>, Flows>;
82 void get_flows(Eigen::Ref<const Eigen::VectorX<FP>> pop, Eigen::Ref<const Eigen::VectorX<FP>> y, FP t,
104 size_t INSCi = this->populations.get_flat_index({i, InfectionState::InfectedNoSymptomsConfirmed});
106 size_t ISyCi = this->populations.get_flat_index({i, InfectionState::InfectedSymptomsConfirmed});
119 //symptomatic are less well quarantined when testing and tracing is overwhelmed so they infect more people
135 pop[Sj] + pop[Ej] + pop[INSj] + pop[ISyj] + pop[ISevj] + pop[ICrj] + pop[Rj]; // without died people
143 flows[this->template get_flat_flow_index<InfectionState::Susceptible, InfectionState::Exposed>({i})] +=
149 icu_occupancy, 0.90 * params.template get<ICUCapacity<FP>>(), params.template get<ICUCapacity<FP>>(),
152 FP deathsPerSevereAdjusted = params.template get<CriticalPerSevere<FP>>()[i] - criticalPerSevereAdjusted;
155 flows[this->template get_flat_flow_index<InfectionState::Exposed, InfectionState::InfectedNoSymptoms>(
163 flows[this->template get_flat_flow_index<InfectionState::InfectedNoSymptoms, InfectionState::Recovered>(
178 flows[this->template get_flat_flow_index<InfectionState::InfectedSymptoms, InfectionState::InfectedSevere>(
181 flows[this->template get_flat_flow_index<InfectionState::InfectedSymptoms, InfectionState::Recovered>(
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})] =
201 flows[this->template get_flat_flow_index<InfectionState::InfectedSevere, InfectionState::Dead>({i})] =
205 flows[this->template get_flat_flow_index<InfectionState::InfectedCritical, InfectionState::Dead>({i})] =
206 params.template get<DeathsPerCritical<FP>>()[i] / params.template get<TimeInfectedCritical<FP>>()[i] *
208 flows[this->template get_flat_flow_index<InfectionState::InfectedCritical, InfectionState::Recovered>(
258 FP get_infections_relative(const Simulation<FP, Base>& model, FP t, const Eigen::Ref<const Eigen::VectorX<FP>>& y);
295 FP delay_npi_implementation; // delay which is needed to implement a NPI that is criterion-dependent
297 const FP dt = dyn_npis.get_thresholds().size() > 0 ? FP(dyn_npis.get_interval().get()) : FP(tmax);
345 std::pair<FP, SimulationTime<FP>> m_dynamic_npi = {-std::numeric_limits<FP>::max(), mio::SimulationTime<FP>(0)};
364 return mio::simulate<FP, Model<FP>, Simulation<FP>>(t0, tmax, dt, model, std::move(integrator_core));
389 FP get_infections_relative(const Simulation<FP, Base>& sim, FP /* t*/, const Eigen::Ref<const Eigen::VectorX<FP>>& y)
415 return mio::failure(mio::StatusCode::OutOfRange, "t_idx is not a valid index for the TimeSeries");
420 //The infected compartments in the SECIR Model are: Exposed, Carrier, Infected, Hospitalized, ICU in respective agegroups
425 //F encodes new Infections and V encodes transition times in the next-generation matrix calculation of R_t
448 ((sim.get_model().parameters.template get<StartDay<FP>>() + sim.get_result().get_time(t_idx)) / 182.5 +
450 ContactMatrixGroup<FP> const& contact_matrix = sim.get_model().parameters.template get<ContactPatterns<FP>>();
479 (params.template get<TestAndTraceCapacity<FP>>()) * params.template get<TestAndTraceCapacityMaxRisk<FP>>(),
550 (cont_freq_eff(i, j) * params.template get<RelativeTransmissionNoSymptoms<FP>>()[(mio::AgeGroup)j] *
556 F(i, j + 2 * num_groups) = sim.get_result().get_value(t_idx)[sim.get_model().populations.get_flat_index(
566 icu_occupancy, 0.90 * params.template get<ICUCapacity<FP>>(), params.template get<ICUCapacity<FP>>(),
571 V(i + num_groups, i + num_groups) = 1 / params.template get<TimeInfectedNoSymptoms<FP>>()[(mio::AgeGroup)i];
592 Eigen::MatrixX<ScalarType> NextGenMatrix(num_infected_compartments * num_groups, 5 * num_groups);
641 auto time_late = std::distance(times.begin(), std::lower_bound(times.begin(), times.end(), t_value));
663 auto get_mobility_factors(const Simulation<FP, Base>& sim, FP /*t*/, const Eigen::Ref<const Eigen::VectorX<FP>>& y)
667 auto&& p_asymp = params.template get<RecoveredPerInfectedNoSymptoms<FP>>().array().template cast<FP>();
668 auto&& p_inf = params.template get<RiskOfInfectionFromSymptomatic<FP>>().array().template cast<FP>();
669 auto&& p_inf_max = params.template get<MaxRiskOfInfectionFromSymptomatic<FP>>().array().template cast<FP>();
676 ((1 - p_asymp) / params.template get<TimeInfectedNoSymptoms<FP>>().array().template cast<FP>() * y_INS.array())
679 auto test_and_trace_capacity_max_risk = FP(params.template get<TestAndTraceCapacityMaxRisk<FP>>());
680 auto riskFromInfectedSymptomatic = smoother_cosine<FP>(test_and_trace_required, test_and_trace_capacity,
686 slice(factors, {Eigen::Index(InfectionState::InfectedSymptoms), Eigen::Index(size_t(params.get_num_groups())),
693 auto test_commuters(Simulation<FP, Base>& sim, Eigen::Ref<Eigen::VectorX<FP>> mobile_population, FP time)
703 auto INSCi = model.populations.get_flat_index({i, InfectionState::InfectedNoSymptomsConfirmed});
707 //put detected commuters in their own compartment so they don't contribute to infections in their home node
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
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 ¶ms)
Definition: ode_secir/model.h:71
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:231
void serialize(IOContext &io) const
serialize this.
Definition: ode_secir/model.h:219
Parameters of an age-resolved SECIR/SECIHURD model.
Definition: ode_secir/parameters.h:399
specialization of compartment model simulation for secir models.
Definition: ode_secir/model.h:267
Eigen::Ref< Eigen::VectorX< FP > > advance(FP tmax)
advance simulation to tmax.
Definition: ode_secir/model.h:288
Simulation(mio::osecir::Model< FP > const &model, FP t0=0., FP dt=0.1)
construct a simulation.
Definition: ode_secir/model.h:275
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
@ InfectedNoSymptomsConfirmed
@ InfectedSevere
@ Susceptible
@ InfectedCritical
@ Recovered
@ InfectedSymptomsConfirmed
@ Exposed
@ InfectedSymptoms
@ InfectedNoSymptoms
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:611
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:628
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
@ OutOfRange
@ UnknownError
boost::outcome_v2::in_place_type_t< T > Tag
Type that is used for overload resolution.
Definition: io.h:407
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:294
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
The AgeGroup struct is used as a dynamically sized tag for all age dependent categories.
Definition: age_group.h:32
mio::Populations< FP, AgeGroup, InfectionState > Populations
Definition: compartmental_model.h:62
Parameters< FP > ParameterSet
Definition: compartmental_model.h:63
A Flow defines a possible transition between two Compartments in a FlowModel.
Definition: flow.h:33
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 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:372
capacity to test and trace contacts of infected for quarantine per day.
Definition: ode_secir/parameters.h:356
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
Generated by