Skip to contents

This class defines the SEIaImIsRD model, an extension of the SEIRD model. The model shows how populations of susceptible, exposed, infectious and recovered individuals evolve over time, in which the infectious individuals are subgrouped into compartments according to different severity of symptoms, i.e. asymptomatic, mild and severe.

Usage

# S4 method for SEIaImIsRD
initial_conditions(object)

# S4 method for SEIaImIsRD
transmission_parameters(object)

# S4 method for SEIaImIsRD
initial_conditions(object) <- value

# S4 method for SEIaImIsRD
transmission_parameters(object) <- value

# S4 method for SEIaImIsRD
run(object, times, solve_method = "lsoda")

# S4 method for SEIaImIsRD
R0(model)

# S4 method for SEIaImIsRD
ode_structure_diagram(model)

Arguments

object

An object of class SEIaImIsRD

value

A numeric named list of values for transmission parameters: beta - a named list of the effective contact rate from each infected group (i.e. rate at which an infected individual exposes susceptible), each element with value in [0,1] kappa - rate of progression from exposed to infectious (the reciprocal is the incubation period), p_symptom - a named list of the probability of exposed individuals moving into each of the different infected groups, Here, only probabilities for the mild (p_symptom.mild) and severe (p_symptom.severe) groups need be specified and the asymptomatic probability is the remainder (1 - p_symptom.mild - p_symptom.severe). Thus we require p_symptom.mild + p_symptom.severe <= 1. gamma - a list of the rate of removal of each infected group (i.e. recovery rate of an infected individual), each element with value in [0,1] mu - a list of the rate of disease-caused mortality of each infected group, each element with value in [0,1]

times

A list of time points of the simulation period

solve_method

A string indicating which ode integrator to use. Default is set to 'lsoda'

model

an SEIaImIsRD model

Value

An object of class SEIaImIsRD with initial population and parameters An object of class SEIaImIsRD with initial population and parameters A list of two dataframes: one with the time steps, time series of S, E, I_asymptomatic, I_mild, I_severe and R population fractions, and one with the time steps, time series of incidences and deaths population fraction An R0 value An ODE-compartmental structure diagram object of class html

Details

Notes:

  1. Total initial population size is normalised to 1.

  2. The current model does not include natural death or birth.

Methods (by generic)

  • initial_conditions: Retrieves initial conditions of SEIaImIsRD model.

  • transmission_parameters: Retrieves transmission parameters of SEIaImIsRD model.

  • initial_conditions<-: Setter method for initial population sizes (in fraction) of the SEIaImIsRD model.

  • transmission_parameters<-: Setter method for transmission parameters of the SEIaImIsRD model.

  • run: Solves the ode system.

    $$\frac{dS(t)}{dt} = -S(t)(\beta.asymptomatic I_{asymptomatic}(t) + \beta.mild I_{mild}(t) + \beta.severe I_{severe}(t))$$ $$\frac{dE(t)}{dt} = S(t)(\beta.asymptomatic I_{asymptomatic}(t) + \beta.mild I_{mild}(t) + \beta.severe I_{severe}(t)) - \kappa E(t)$$ $$\frac{dI_{asymptomatic}(t)}{dt} = p_symptom.asymptomatic \kappa E(t) - (\gamma.asymptomatic + \mu.asymptomatic) I_{asymptomatic}(t)$$ $$\frac{dI_{mild}(t)}{dt} = p_symptom.mild\kappa E(t) - (\gamma.mild + \mu.mild) I_{mild}(t)$$ $$\frac{dI_{severe}(t)}{dt} = p_symptom.severe\kappa E(t) - (\gamma.severe + \mu.severe) I_{severe}(t)$$ $$\frac{dR(t)}{dt} = \gamma.asymptomatic I_{asymptomatic}(t) + \gamma.mild I_{mild}(t) + \gamma.severe I_{severe}(t)$$ $$\frac{dD(t)}{dt} = \mu.asymptomatic I_{asymptomatic}(t) + \mu.mild I_{mild}(t) + \mu.severe I_{severe}(t)$$

  • R0: Calculate the basic reproduction number (R_0) of the system using the next generation matrix approach.

  • ode_structure_diagram: Prints a compartmental diagram for the SEIaImIsRD model

Slots

initial_condition_names

list of names of initial conditions (characters). Default is list("S0", "E0", "I_asymptomatic0", "I_mild0", "I_severe0", "R0", "D0").

transmission_parameter_names

list of names of transmission parameters (characters). Default is list("beta", "kappa", "p_symptom", "gamma", "mu").

initial_conditions

list of values for initial conditions (double).

transmission_parameters

list of values for transmission parameters (double).

output_names

list of compartments name which are used by the model and incidence.

See also

https://www.ncbi.nlm.nih.gov/pmc/articles/PMC6002118/pdf/main.pdf mathematical details of the next generation matrix approach.

The R0 parameter is given by: $$R_0 = \rho(FV^{-1})$$, where $$F=\frac{\partial F_i{(x_0)}}{\partial x_j}$$ and $$V=\frac{\partial V_i{(x_0)}}{\partial x_j}$$ and rho represents the spectral radius of a matrix (i.e. the largest absolute of eigenvalue). The F_i are the new infections, while the V_i transfers of infections from one compartment to another. x_0 is the disease-free equilibrium state.