model.h Source File
|
CPP API
|
ode_secirvvs/model.h
Go to the documentation of this file.
46 Flow<InfectionState::InfectedNoSymptomsNaiveConfirmed, InfectionState::InfectedSymptomsNaiveConfirmed>,
47 Flow<InfectionState::InfectedNoSymptomsNaiveConfirmed, InfectionState::SusceptibleImprovedImmunity>,
51 Flow<InfectionState::InfectedSymptomsNaiveConfirmed, InfectionState::SusceptibleImprovedImmunity>,
60 Flow<InfectionState::InfectedNoSymptomsPartialImmunity, InfectionState::InfectedSymptomsPartialImmunity>,
61 Flow<InfectionState::InfectedNoSymptomsPartialImmunity, InfectionState::SusceptibleImprovedImmunity>,
62 Flow<InfectionState::InfectedNoSymptomsPartialImmunityConfirmed, InfectionState::InfectedSymptomsPartialImmunityConfirmed>,
63 Flow<InfectionState::InfectedNoSymptomsPartialImmunityConfirmed, InfectionState::SusceptibleImprovedImmunity>,
64 Flow<InfectionState::InfectedSymptomsPartialImmunity, InfectionState::InfectedSeverePartialImmunity>,
65 Flow<InfectionState::InfectedSymptomsPartialImmunity, InfectionState::SusceptibleImprovedImmunity>,
66 Flow<InfectionState::InfectedSymptomsPartialImmunityConfirmed, InfectionState::InfectedSeverePartialImmunity>,
67 Flow<InfectionState::InfectedSymptomsPartialImmunityConfirmed, InfectionState::SusceptibleImprovedImmunity>,
68 Flow<InfectionState::InfectedSeverePartialImmunity, InfectionState::InfectedCriticalPartialImmunity>,
69 Flow<InfectionState::InfectedSeverePartialImmunity, InfectionState::SusceptibleImprovedImmunity>,
72 Flow<InfectionState::InfectedCriticalPartialImmunity, InfectionState::SusceptibleImprovedImmunity>,
75 Flow<InfectionState::ExposedImprovedImmunity, InfectionState::InfectedNoSymptomsImprovedImmunity>,
76 Flow<InfectionState::InfectedNoSymptomsImprovedImmunity, InfectionState::InfectedSymptomsImprovedImmunity>,
77 Flow<InfectionState::InfectedNoSymptomsImprovedImmunity, InfectionState::SusceptibleImprovedImmunity>,
78 Flow<InfectionState::InfectedNoSymptomsImprovedImmunityConfirmed, InfectionState::InfectedSymptomsImprovedImmunityConfirmed>,
79 Flow<InfectionState::InfectedNoSymptomsImprovedImmunityConfirmed, InfectionState::SusceptibleImprovedImmunity>,
80 Flow<InfectionState::InfectedSymptomsImprovedImmunity, InfectionState::InfectedSevereImprovedImmunity>,
81 Flow<InfectionState::InfectedSymptomsImprovedImmunity, InfectionState::SusceptibleImprovedImmunity>,
82 Flow<InfectionState::InfectedSymptomsImprovedImmunityConfirmed, InfectionState::InfectedSevereImprovedImmunity>,
83 Flow<InfectionState::InfectedSymptomsImprovedImmunityConfirmed, InfectionState::SusceptibleImprovedImmunity>,
84 Flow<InfectionState::InfectedSevereImprovedImmunity, InfectionState::InfectedCriticalImprovedImmunity>,
85 Flow<InfectionState::InfectedSevereImprovedImmunity, InfectionState::SusceptibleImprovedImmunity>,
88 Flow<InfectionState::InfectedCriticalImprovedImmunity, InfectionState::SusceptibleImprovedImmunity>>;
93 : public FlowModel<FP, InfectionState, mio::Populations<FP, AgeGroup, InfectionState>, Parameters<FP>, Flows>
95 using Base = FlowModel<FP, InfectionState, mio::Populations<FP, AgeGroup, InfectionState>, Parameters<FP>, Flows>;
107 : Model(Populations({AgeGroup(num_agegroups), InfectionState::Count}), ParameterSet(AgeGroup(num_agegroups)))
111 void get_flows(Eigen::Ref<const Eigen::VectorX<FP>> pop, Eigen::Ref<const Eigen::VectorX<FP>> y, FP t,
129 this->populations.get_from(pop, {i, InfectionState::InfectedNoSymptomsPartialImmunityConfirmed}) +
130 this->populations.get_from(pop, {i, InfectionState::InfectedNoSymptomsImprovedImmunityConfirmed}));
145 size_t INSNCi = this->populations.get_flat_index({i, InfectionState::InfectedNoSymptomsNaiveConfirmed});
146 size_t ISyNCi = this->populations.get_flat_index({i, InfectionState::InfectedSymptomsNaiveConfirmed});
148 size_t SPIi = this->populations.get_flat_index({i, InfectionState::SusceptiblePartialImmunity});
150 size_t INSPIi = this->populations.get_flat_index({i, InfectionState::InfectedNoSymptomsPartialImmunity});
151 size_t ISyPIi = this->populations.get_flat_index({i, InfectionState::InfectedSymptomsPartialImmunity});
152 size_t ISevPIi = this->populations.get_flat_index({i, InfectionState::InfectedSeverePartialImmunity});
153 size_t ICrPIi = this->populations.get_flat_index({i, InfectionState::InfectedCriticalPartialImmunity});
156 this->populations.get_flat_index({i, InfectionState::InfectedNoSymptomsPartialImmunityConfirmed});
158 this->populations.get_flat_index({i, InfectionState::InfectedSymptomsPartialImmunityConfirmed});
161 size_t INSIIi = this->populations.get_flat_index({i, InfectionState::InfectedNoSymptomsImprovedImmunity});
162 size_t ISyIIi = this->populations.get_flat_index({i, InfectionState::InfectedSymptomsImprovedImmunity});
163 size_t ISevIIi = this->populations.get_flat_index({i, InfectionState::InfectedSevereImprovedImmunity});
164 size_t ICrIIi = this->populations.get_flat_index({i, InfectionState::InfectedCriticalImprovedImmunity});
167 this->populations.get_flat_index({i, InfectionState::InfectedNoSymptomsImprovedImmunityConfirmed});
169 this->populations.get_flat_index({i, InfectionState::InfectedSymptomsImprovedImmunityConfirmed});
171 size_t SIIi = this->populations.get_flat_index({i, InfectionState::SusceptibleImprovedImmunity});
185 //symptomatic are less well quarantined when testing and tracing is overwhelmed so they infect more people
206 size_t SIIj = this->populations.get_flat_index({j, InfectionState::SusceptibleImprovedImmunity});
208 size_t INSNCj = this->populations.get_flat_index({j, InfectionState::InfectedNoSymptomsNaiveConfirmed});
209 size_t ISyNCj = this->populations.get_flat_index({j, InfectionState::InfectedSymptomsNaiveConfirmed});
211 size_t SPIj = this->populations.get_flat_index({j, InfectionState::SusceptiblePartialImmunity});
215 size_t ISyPIj = this->populations.get_flat_index({j, InfectionState::InfectedSymptomsPartialImmunity});
216 size_t ISevPIj = this->populations.get_flat_index({j, InfectionState::InfectedSeverePartialImmunity});
217 size_t ICrPIj = this->populations.get_flat_index({j, InfectionState::InfectedCriticalPartialImmunity});
220 this->populations.get_flat_index({j, InfectionState::InfectedNoSymptomsPartialImmunityConfirmed});
222 this->populations.get_flat_index({j, InfectionState::InfectedSymptomsPartialImmunityConfirmed});
227 size_t ISyIIj = this->populations.get_flat_index({j, InfectionState::InfectedSymptomsImprovedImmunity});
228 size_t ISevIIj = this->populations.get_flat_index({j, InfectionState::InfectedSevereImprovedImmunity});
229 size_t ICrIIj = this->populations.get_flat_index({j, InfectionState::InfectedCriticalImprovedImmunity});
232 this->populations.get_flat_index({j, InfectionState::InfectedNoSymptomsImprovedImmunityConfirmed});
234 this->populations.get_flat_index({j, InfectionState::InfectedSymptomsImprovedImmunityConfirmed});
244 FP Nj = pop[SNj] + pop[ENj] + pop[INSNj] + pop[ISyNj] + pop[ISevNj] + pop[ICrNj] + pop[INSNCj] +
246 pop[INSPICj] + pop[ISyPICj] + pop[SIIj] + pop[EIIj] + pop[INSIIj] + pop[ISyIIj] + pop[ISevIIj] +
271 // is different for different vaccination status. This is not the case here and in addition, ICUCapacity
274 icu_occupancy, 0.90 * params.template get<ICUCapacity<FP>>(), params.template get<ICUCapacity<FP>>(),
277 FP deathsPerSevereAdjusted = params.template get<CriticalPerSevere<FP>>()[i] - criticalPerSevereAdjusted;
334 flows[this->template get_flat_flow_index<InfectionState::InfectedSevereNaive, InfectionState::DeadNaive>(
338 flows[this->template get_flat_flow_index<InfectionState::InfectedCriticalNaive, InfectionState::DeadNaive>(
364 flows[this->template get_flat_flow_index<InfectionState::InfectedNoSymptomsPartialImmunityConfirmed,
369 flows[this->template get_flat_flow_index<InfectionState::InfectedNoSymptomsPartialImmunityConfirmed,
388 flows[this->template get_flat_flow_index<InfectionState::InfectedSymptomsPartialImmunityConfirmed,
394 flows[this->template get_flat_flow_index<InfectionState::InfectedSymptomsPartialImmunityConfirmed,
403 reducInfectedSevereCriticalDeadPartialImmunity / reducInfectedSevereCriticalDeadPartialImmunity *
408 (1 - (reducInfectedSevereCriticalDeadPartialImmunity / reducInfectedSevereCriticalDeadPartialImmunity) *
414 (reducInfectedSevereCriticalDeadPartialImmunity / reducInfectedSevereCriticalDeadPartialImmunity) *
420 (reducInfectedSevereCriticalDeadPartialImmunity / reducInfectedSevereCriticalDeadPartialImmunity) *
421 params.template get<DeathsPerCritical<FP>>()[i] / params.template get<TimeInfectedCritical<FP>>()[i] *
426 (1 - (reducInfectedSevereCriticalDeadPartialImmunity / reducInfectedSevereCriticalDeadPartialImmunity) *
447 flows[this->template get_flat_flow_index<InfectionState::InfectedNoSymptomsImprovedImmunityConfirmed,
452 flows[this->template get_flat_flow_index<InfectionState::InfectedNoSymptomsImprovedImmunityConfirmed,
467 (1 - (reducInfectedSevereCriticalDeadImprovedImmunity / reducInfectedSymptomsImprovedImmunity) *
471 flows[this->template get_flat_flow_index<InfectionState::InfectedSymptomsImprovedImmunityConfirmed,
477 flows[this->template get_flat_flow_index<InfectionState::InfectedSymptomsImprovedImmunityConfirmed,
479 (1 - (reducInfectedSevereCriticalDeadImprovedImmunity / reducInfectedSymptomsImprovedImmunity) *
486 reducInfectedSevereCriticalDeadImprovedImmunity / reducInfectedSevereCriticalDeadImprovedImmunity *
492 (reducInfectedSevereCriticalDeadImprovedImmunity / reducInfectedSevereCriticalDeadImprovedImmunity) *
498 (reducInfectedSevereCriticalDeadImprovedImmunity / reducInfectedSevereCriticalDeadImprovedImmunity) *
504 (reducInfectedSevereCriticalDeadImprovedImmunity / reducInfectedSevereCriticalDeadImprovedImmunity) *
505 params.template get<DeathsPerCritical<FP>>()[i] / params.template get<TimeInfectedCritical<FP>>()[i] *
511 (reducInfectedSevereCriticalDeadImprovedImmunity / reducInfectedSevereCriticalDeadImprovedImmunity) *
560 FP get_infections_relative(const Simulation<FP, Base>& model, FP t, const Eigen::Ref<const Eigen::VectorX<FP>>& y);
595 void apply_variant(const FP t, const CustomIndexArray<UncertainValue<FP>, AgeGroup> base_infectiousness)
601 auto start_day_new_variant = this->get_model().parameters.template get<StartDayNewVariant<FP>>();
611 this->get_model().parameters.template get<TransmissionProbabilityOnContact<FP>>()[i] = new_transmission;
647 log_warning("too many first vaccinated at time {}: setting first_vacc from {} to {}", t, first_vacc,
657 log_warning("too many fully vaccinated at time {}: setting full_vacc from {} to {}", t, full_vacc,
684 // in the apply_variant function, we adjust the TransmissionProbabilityOnContact parameter. We need to store
685 // the base value to use it in the apply_variant function and also to reset the parameter after the simulation.
686 auto base_infectiousness = this->get_model().parameters.template get<TransmissionProbabilityOnContact<FP>>();
747 // reset TransmissionProbabilityOnContact. This is important for the graph simulation where the advance
749 this->get_model().parameters.template get<TransmissionProbabilityOnContact<FP>>() = base_infectiousness;
756 std::pair<FP, SimulationTime<FP>> m_dynamic_npi = {-std::numeric_limits<FP>::max(), SimulationTime<FP>(0)};
775 return mio::simulate<FP, Model<FP>, Simulation<FP>>(t0, tmax, dt, model, std::move(integrator_core));
800 FP get_infections_relative(const Simulation<FP, Base>& sim, FP /*t*/, const Eigen::Ref<const Eigen::VectorX<FP>>& y)
805 sum_inf += sim.get_model().populations.get_from(y, {i, InfectionState::InfectedSymptomsNaiveConfirmed});
806 sum_inf += sim.get_model().populations.get_from(y, {i, InfectionState::InfectedSymptomsPartialImmunity});
807 sum_inf += sim.get_model().populations.get_from(y, {i, InfectionState::InfectedSymptomsImprovedImmunity});
809 sim.get_model().populations.get_from(y, {i, InfectionState::InfectedSymptomsPartialImmunityConfirmed});
811 sim.get_model().populations.get_from(y, {i, InfectionState::InfectedSymptomsImprovedImmunityConfirmed});
829 auto get_mobility_factors(const Simulation<Base>& sim, FP /*t*/, const Eigen::Ref<const Eigen::VectorX<FP>>& y)
834 auto& p_asymp = params.template get<RecoveredPerInfectedNoSymptoms<FP>>().array().template cast<FP>();
835 auto& p_inf = params.template get<RiskOfInfectionFromSymptomatic<FP>>().array().template cast<FP>();
836 auto& p_inf_max = params.template get<MaxRiskOfInfectionFromSymptomatic<FP>>().array().template cast<FP>();
847 ((1 - p_asymp) / params.template get<TimeInfectedNoSymptoms<FP>>().array().template cast<FP>() * y_INS.array())
850 smoother_cosine<FP>(test_and_trace_required, FP(params.template get<TestAndTraceCapacity<FP>>()),
857 slice(factors, {Eigen::Index(InfectionState::InfectedSymptomsNaive), Eigen::Index(size_t(params.get_num_groups())),
870 auto test_commuters(Simulation<FP, Base>& sim, Eigen::Ref<Eigen::VectorX<FP>> mobile_population, FP time)
880 auto ISyNCi = model.populations.get_flat_index({i, InfectionState::InfectedSymptomsNaiveConfirmed});
882 auto INSNCi = model.populations.get_flat_index({i, InfectionState::InfectedNoSymptomsNaiveConfirmed});
884 auto ISPIi = model.populations.get_flat_index({i, InfectionState::InfectedSymptomsPartialImmunity});
885 auto ISPICi = model.populations.get_flat_index({i, InfectionState::InfectedSymptomsPartialImmunityConfirmed});
886 auto INSPIi = model.populations.get_flat_index({i, InfectionState::InfectedNoSymptomsPartialImmunity});
888 model.populations.get_flat_index({i, InfectionState::InfectedNoSymptomsPartialImmunityConfirmed});
890 auto ISyIIi = model.populations.get_flat_index({i, InfectionState::InfectedSymptomsImprovedImmunity});
891 auto ISyIICi = model.populations.get_flat_index({i, InfectionState::InfectedSymptomsImprovedImmunityConfirmed});
892 auto INSIIi = model.populations.get_flat_index({i, InfectionState::InfectedNoSymptomsImprovedImmunity});
894 model.populations.get_flat_index({i, InfectionState::InfectedNoSymptomsImprovedImmunityConfirmed});
896 //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 class template for an array with custom indices.
Definition: custom_index_array.h:136
A FlowModel is a CompartmentalModel defined by the flows between compartments.
Definition: flow_model.h:74
Interface class defining the integration step used in a SystemIntegrator.
Definition: integrator.h:48
decltype(auto) get_from(Arr &&y, Index const &cats) const
get_from returns the value of a flat container at the flat index corresponding to set of enum values.
Definition: populations.h:123
Represents the simulation time as an integer index.
Definition: simulation_day.h:32
Parameters of the simulation that are the same everywhere within the Model.
Definition: abm/parameters.h:764
Run the Simulation in discrete steps, evolve the Model and report results.
Definition: models/abm/simulation.h:37
Definition: ode_secirvvs/model.h:94
void serialize(IOContext &io) const
serialize this.
Definition: ode_secirvvs/model.h:522
Model(const Populations &pop, const ParameterSet ¶ms)
Definition: ode_secirvvs/model.h:101
static IOResult< Model > deserialize(IOContext &io)
deserialize an object of this class.
Definition: ode_secirvvs/model.h:534
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_secirvvs/model.h:111
Parameters of an age-resolved SECIR/SECIHURD model with paths for partial and improved immunity throu...
Definition: ode_secirvvs/parameters.h:645
specialization of compartment model simulation for the SECIRVVS model.
Definition: ode_secirvvs/model.h:569
FP m_t_last_npi_check
Definition: ode_secirvvs/model.h:755
void apply_vaccination(FP t)
Definition: ode_secirvvs/model.h:616
Eigen::Ref< Eigen::VectorX< FP > > advance(FP tmax)
advance simulation to tmax.
Definition: ode_secirvvs/model.h:674
Simulation(mio::osecirvvs::Model< FP > const &model, FP t0=0., FP dt=0.1)
construct a simulation.
Definition: ode_secirvvs/model.h:577
void apply_variant(const FP t, const CustomIndexArray< UncertainValue< FP >, AgeGroup > base_infectiousness)
Applies the effect of a new variant of a disease to the transmission probability of the model.
Definition: ode_secirvvs/model.h:595
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 double floor(const ad::internal::active_type< AD_TAPE_REAL, DATA_HANDLER_1 > &x)
Definition: ad.hpp:2451
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::binary_intermediate_aa< AD_TAPE_REAL, ad::internal::active_type< AD_TAPE_REAL, DATA_HANDLER_1 >, ad::internal::active_type< AD_TAPE_REAL, DATA_HANDLER_1 >, ad::operations::ad_pow_aa< AD_TAPE_REAL > > pow(const ad::internal::active_type< AD_TAPE_REAL, DATA_HANDLER_1 > &x1, const ad::internal::active_type< AD_TAPE_REAL, DATA_HANDLER_1 > &x2)
Definition: ad.hpp:1610
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
auto get_mobility_factors(const Simulation< Base > &sim, FP, const Eigen::Ref< const Eigen::VectorX< FP >> &y)
Get mobility factors.
Definition: ode_secirvvs/model.h:829
@ InfectedNoSymptomsNaive
@ InfectedNoSymptomsImprovedImmunity
@ InfectedNoSymptomsNaiveConfirmed
@ InfectedNoSymptomsPartialImmunity
auto test_commuters(Simulation< FP, Base > &sim, Eigen::Ref< Eigen::VectorX< FP >> mobile_population, FP time)
Definition: ode_secirvvs/model.h:870
FP get_infections_relative(const Simulation< FP, Base > &model, FP t, const Eigen::Ref< const Eigen::VectorX< FP >> &y)
get percentage of infections per total population.
Definition: ode_secirvvs/model.h:800
auto simulate(FP t0, FP tmax, FP dt, const Model< FP > &model, std::unique_ptr< OdeIntegratorCore< FP >> &&integrator_core=nullptr)
Specialization of simulate for SECIRVVS models using Simulation.
Definition: ode_secirvvs/model.h:772
auto simulate_flows(FP t0, FP tmax, FP dt, const Model< FP > &model, std::unique_ptr< OdeIntegratorCore< FP >> &&integrator_core=nullptr)
Specialization of simulate for SECIRVVS models using the FlowSimulation.
Definition: ode_secirvvs/model.h:791
A collection of classes to simplify handling of matrix shapes in meta programming.
Definition: models/abm/analyze_result.h:30
boost::outcome_v2::in_place_type_t< T > Tag
Type that is used for overload resolution.
Definition: io.h:407
void log_warning(spdlog::string_view_t fmt, const Args &... args)
Definition: logging.h:112
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 slice(V &&v, Seq< Eigen::Index > elems)
take a regular slice of a row or column vector.
Definition: eigen_util.h:107
void implement_dynamic_npis(DampingExprGroup &damping_expr_group, const std::vector< DampingSampling< FP >> &npis, SimulationTime< FP > begin, SimulationTime< FP > end, MakeMatrix &&make_matrix)
implement dynamic NPIs for a time span.
Definition: dynamic_npis.h:279
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
The AgeGroup struct is used as a dynamically sized tag for all age dependent categories.
Definition: age_group.h:32
Populations populations
Definition: compartmental_model.h:156
mio::Populations< FP, AgeGroup, InfectionState > Populations
Definition: compartmental_model.h:62
Parameters< FP > ParameterSet
Definition: compartmental_model.h:63
ParameterSet parameters
Definition: compartmental_model.h:157
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
Multiplier for the test and trace capacity to determine when it is considered overloaded by symptomat...
Definition: ode_secirvvs/parameters.h:144
capacity to test and trace contacts of infected for quarantine per day.
Definition: ode_secirvvs/parameters.h:112
the (mean) time in day unit for asymptomatic cases that are infected but have not yet developed sympt...
Definition: ode_secirvvs/parameters.h:225
Generated by