Internal Documentation
Documentation for EpiLatentModels.jl's internal interface.
Contents
Index
EpiAware.EpiLatentModels._seir_jac_prototypeEpiAware.EpiLatentModels._sir_jac_prototypeEpiAware.EpiLatentModels.ARStepEpiAware.EpiLatentModels.ARStepEpiAware.EpiLatentModels.MAStepEpiAware.EpiLatentModels.MAStepEpiAware.EpiLatentModels.RWStepEpiAware.EpiLatentModels.RWStepEpiAware.EpiAwareBase.broadcast_nEpiAware.EpiAwareBase.broadcast_nEpiAware.EpiAwareBase.generate_latentEpiAware.EpiAwareBase.generate_latentEpiAware.EpiAwareBase.generate_latentEpiAware.EpiAwareBase.generate_latentEpiAware.EpiAwareBase.generate_latentEpiAware.EpiAwareBase.generate_latentEpiAware.EpiAwareBase.generate_latentEpiAware.EpiAwareBase.generate_latentEpiAware.EpiAwareBase.generate_latentEpiAware.EpiAwareBase.generate_latentEpiAware.EpiAwareBase.generate_latentEpiAware.EpiAwareBase.generate_latentEpiAware.EpiAwareBase.generate_latentEpiAware.EpiAwareBase.generate_latentEpiAware.EpiAwareBase.generate_latentEpiAware.EpiAwareUtils.get_stateEpiAware.EpiLatentModels._expand_distEpiAware.EpiLatentModels._seir_functionEpiAware.EpiLatentModels._seir_jacEpiAware.EpiLatentModels._seir_vfEpiAware.EpiLatentModels._sir_functionEpiAware.EpiLatentModels._sir_jacEpiAware.EpiLatentModels._sir_vf
Internal API
EpiAware.EpiLatentModels._seir_jac_prototype — ConstantSparse Jacobian matrix prototype for the basic SEIR model written in density/per-capita form.
EpiAware.EpiLatentModels._sir_jac_prototype — ConstantSparse Jacobian matrix prototype for the basic SIR model written in density/per-capita form.
EpiAware.EpiLatentModels.ARStep — Typestruct ARStep{D<:(AbstractVector{<:Real})} <: AbstractAccumulationStepThe autoregressive (AR) step function struct
Fields
damp_AR::AbstractVector{<:Real}
EpiAware.EpiLatentModels.ARStep — MethodThe autoregressive (AR) step function for use with accumulate_scan.
EpiAware.EpiLatentModels.MAStep — Typestruct MAStep{C<:(AbstractVector{<:Real})} <: AbstractAccumulationStepThe moving average (MA) step function struct
Fields
θ::AbstractVector{<:Real}
EpiAware.EpiLatentModels.MAStep — MethodThe moving average (MA) step function for use with accumulate_scan.
EpiAware.EpiLatentModels.RWStep — Typestruct RWStep <: AbstractAccumulationStepThe random walk (RW) step function struct
Fields
EpiAware.EpiLatentModels.RWStep — MethodThe random walk (RW) step function for use with accumulate_scan.
EpiAware.EpiAwareBase.broadcast_n — Methodbroadcast_n(_::RepeatBlock, n, period) -> Any
A function that returns the length of the latent periods to generate using the RepeatBlock rule which is equal n divided by the period and rounded up to the nearest integer.
Arguments
rule::RepeatBlock: The broadcasting rule.n: The number of samples to generate.period: The period of the broadcast.
EpiAware.EpiAwareBase.broadcast_n — Methodbroadcast_n(_::RepeatEach, n, period) -> Any
A function that returns the length of the latent periods to generate using the RepeatEach rule which is equal to the period.
Arguments
rule::RepeatEach: The broadcasting rule.n: The number of samples to generate.period: The period of the broadcast.
Returns
m: The length of the latent periods to generate.
EpiAware.EpiAwareBase.generate_latent — Methodgenerate_latent(latent_model::AR, n) -> Any
Generate a latent AR series.
Arguments
latent_model::AR: The AR model.n::Int: The length of the AR series.
Returns
ar::Vector{Float64}: The generated AR series.
Notes
- The length of
damp_priorandinit_priormust be the same. nmust be longer than the order of the autoregressive process.
EpiAware.EpiAwareBase.generate_latent — Methodgenerate_latent(model::BroadcastLatentModel, n) -> Any
Generates latent periods using the specified model and n number of samples.
Arguments
model::BroadcastLatentModel: The broadcast latent model.n::Any: The number of samples to generate.
Returns
broadcasted_latent: The generated broadcasted latent periods.
EpiAware.EpiAwareBase.generate_latent — Methodgenerate_latent(
latent_models::CombineLatentModels,
n
) -> Any
Generate latent variables using a combination of multiple latent models.
Arguments
latent_models::CombineLatentModels: An instance of theCombineLatentModelstype representing the collection of latent models.n: The number of latent variables to generate.
Returns
- The combined latent variables generated from all the models.
Example
EpiAware.EpiAwareBase.generate_latent — Methodgenerate_latent(latent_models::ConcatLatentModels, n) -> Any
Generate latent variables by concatenating multiple latent models.
Arguments
latent_models::ConcatLatentModels: An instance of theConcatLatentModelstype representing the collection of latent models.n: The number of latent variables to generate.
Returns
concatenated_latents: The combined latent variables generated from all the models.latent_aux: A tuple containing the auxiliary latent variables generated from each individual model.
EpiAware.EpiAwareBase.generate_latent — Methodgenerate_latent(latent_model::DiffLatentModel, n) -> Any
Generate a Turing model for n-step latent process $Z_t$ using a differenced latent model defined by latent_model.
Arguments
latent_model::DiffLatentModel: The differential latent model.n: The length of the latent variables.
Turing model specifications
Sampled random variables
latent_init: The initial latent process variables.- Other random variables defined by
model<:AbstractTuringLatentModelfield of the undifferenced model.
Generated quantities
- A tuple containing the generated latent process as its first argument and a
NamedTupleof sampled auxiliary variables as second argument.
Example usage with DiffLatentModel model constructor
generate_latent can be used to construct a Turing model for the differenced latent process. In this example, the underlying undifferenced process is a RandomWalk model.
First, we construct a RandomWalk struct with an initial value prior and a step size standard deviation prior.
using Distributions, EpiAware
rw = RandomWalk(Normal(0.0, 1.0), truncated(Normal(0.0, 0.05), 0.0, Inf))Then, we can use DiffLatentModel to construct a DiffLatentModel for d-fold differenced process with rw as the undifferenced latent process.
We have two constructor options for DiffLatentModel. The first option is to supply a common prior distribution for the initial terms and specify d as follows:
diff_model = DiffLatentModel(rw, Normal(); d = 2)Or we can supply a vector of priors for the initial terms and d is inferred as follows:
diff_model2 = DiffLatentModel(;undiffmodel = rw, init_priors = [Normal(), Normal()])Then, we can use generate_latent to construct a Turing model for the differenced latent process generating a length n process,
# Construct a Turing model
n = 100
difference_mdl = generate_latent(diff_model, n)Now we can use the Turing PPL API to sample underlying parameters and generate the unobserved latent process.
#Sample random parameters from prior
θ = rand(difference_mdl)
#Get a sampled latent process as a generated quantity from the model
(Z_t, _) = generated_quantities(difference_mdl, θ)
Z_tEpiAware.EpiAwareBase.generate_latent — Methodgenerate_latent(latent_model::FixedIntercept, n) -> Any
Generate a latent intercept series with a fixed intercept value.
Arguments
latent_model::FixedIntercept: The fixed intercept latent model.n: The number of latent variables to generate.
Returns
latent_vars: An array of lengthnfilled with the fixed intercept value.
EpiAware.EpiAwareBase.generate_latent — Methodgenerate_latent(obs_model::HierarchicalNormal, n) -> Any
function EpiAwareBase.generate_latent(obs_model::HierarchicalNormal, n)Generate latent variables from the hierarchical normal distribution.
Arguments
obs_model::HierarchicalNormal: The hierarchical normal distribution model.n: Number of latent variables to generate.
Returns
η_t: Generated latent variables.
EpiAware.EpiAwareBase.generate_latent — Methodgenerate_latent(model::IID, n) -> Any
Generate latent variables from the IID (Independent and Identically Distributed) model.
Arguments
model::IID: The IID model.n: Number of latent variables to generate.
Returns
ϵ_t: Generated latent variables.
EpiAware.EpiAwareBase.generate_latent — Methodgenerate_latent(latent_model::Intercept, n) -> Any
Generate a latent intercept series.
Arguments
latent_model::Intercept: The intercept model.n::Int: The length of the intercept series.
Returns
intercept::Vector{Float64}: The generated intercept series.
EpiAware.EpiAwareBase.generate_latent — Methodgenerate_latent(latent_model::MA, n) -> Any
Generate a latent MA series.
Arguments
latent_model::MA: The MA model.n::Int: The length of the MA series.
Returns
ma::Vector{Float64}: The generated MA series.
Notes
nmust be longer than the order of the moving average process.
EpiAware.EpiAwareBase.generate_latent — Methodgenerate_latent(latent_model::Null, n) -> Any
Generates nothing as latent variables for the given latent_model of type Null.
Example
using EpiAware
null = Null()
null_mdl = generate_latent(null, 10)
isnothing(null_mdl())
# output
true
EpiAware.EpiAwareBase.generate_latent — Methodgenerate_latent(latent_model::RandomWalk, n) -> Any
Generate a latent RW series using accumulate_scan.
Arguments
latent_model::RandomWalk: The RandomWalk model.n::Int: The length of the RW series.
Returns
rw::Vector{Float64}: The generated RW series.
Notes
nmust be greater than 0.
EpiAware.EpiAwareBase.generate_latent — Methodgenerate_latent(params::SEIRParams, n) -> Any
Generates the initial parameters and initial conditions for the basic SEIR model.
SEIR model
\[\begin{aligned} \frac{dS}{dt} &= -\beta SI \\ \frac{dE}{dt} &= \beta SI - \alpha E \\ \frac{dI}{dt} &= \alpha E - \gamma I \\ \frac{dR}{dt} &= \gamma I \end{aligned}\]
Where S is the proportion of the population that is susceptible, E is the proportion of the population that is exposed, I is the proportion of the population that is infected and R is the proportion of the population that is recovered. The parameters are the infectiousness β, the incubation rate α and the recovery rate γ.
Initial conditions
For this version of the SEIR model we sample the initial proportion of the population that is infected (exposed or infectious). The proportion of the infected group that is exposed is α / (α + γ) and the proportion of the infected group that is infectious is γ / (α + γ). The reason for this is that these are the equilibrium proportions in a constant incidence environment.
Example
using EpiAware, OrdinaryDiffEq, Distributions
# Create an instance of SIRParams
seirparams = SEIRParams(
tspan = (0.0, 30.0),
infectiousness = LogNormal(log(0.3), 0.05),
incubation_rate = LogNormal(log(0.1), 0.05),
recovery_rate = LogNormal(log(0.1), 0.05),
initial_prop_infected = Beta(1, 99)
)
seirparam_mdl = generate_latent(seirparams, nothing)
# Sample the parameters of SEIR model
sampled_params = rand(seirparam_mdl)
nothingReturns
A tuple (u0, p) where:
u0: A vector representing the initial state of the system[S₀, E₀, I₀, R₀]whereS₀
is the initial proportion of susceptible individuals, E₀ is the initial proportion of exposed individuals,I₀ is the initial proportion of infected individuals, and R₀ is the initial proportion of recovered individuals.
p: A vector containing the parameters[β, α, γ]whereβis the infectiousness rate,
α is the incubation rate, and γ is the recovery rate.
EpiAware.EpiAwareBase.generate_latent — Methodgenerate_latent(params::SIRParams, n) -> Any
Generates the initial parameters and initial conditions for the basic SIR model.
SIR model
\[\begin{aligned} \frac{dS}{dt} &= -\beta SI \\ \frac{dI}{dt} &= \beta SI - \gamma I \\ \frac{dR}{dt} &= \gamma I \end{aligned}\]
Where S is the proportion of the population that is susceptible, I is the proportion of the population that is infected and R is the proportion of the population that is recovered. The parameters are the infectiousness β and the recovery rate γ.
Example
using EpiAware, OrdinaryDiffEq, Distributions
# Create an instance of SIRParams
sirparams = SIRParams(
tspan = (0.0, 30.0),
infectiousness = LogNormal(log(0.3), 0.05),
recovery_rate = LogNormal(log(0.1), 0.05),
initial_prop_infected = Beta(1, 99)
)
sirparam_mdl = generate_latent(sirparams, nothing)
#sample the parameters of SIR model
sampled_params = rand(sirparam_mdl)
nothingReturns
- A tuple
(u0, p)where:u0: A vector representing the initial state of the system[S₀, I₀, R₀]whereS₀is the initial proportion of susceptible individuals,I₀is the initial proportion of infected individuals, andR₀is the initial proportion of recovered individuals.p: A vector containing the parameters[β, γ]whereβis the infectiousness rate andγis the recovery rate.
EpiAware.EpiAwareBase.generate_latent — Methodgenerate_latent(model::TransformLatentModel, n) -> Any
generate_latent(model::TransformLatentModel, n)Generate latent variables using the specified TransformLatentModel.
Arguments
model::TransformLatentModel: TheTransformLatentModelto generate latent variables from.n: The number of latent variables to generate.
Returns
- The transformed latent variables.
EpiAware.EpiAwareUtils.get_state — Methodget_state(
acc_step::EpiAware.EpiLatentModels.MAStep,
initial_state,
state
) -> Any
The MA step function method for get_state.
EpiAware.EpiLatentModels._expand_dist — Method_expand_dist(
dist::Vector{D} where D<:Distributions.Distribution
) -> Any
Internal function to expand a vector of distributions into a product distribution.
Implementation note
If all distributions in the input vector dist are the same, it returns a filled distribution using filldist. Otherwise, it returns an array distribution using arraydist.
EpiAware.EpiLatentModels._seir_function — FunctionInternal function for the ODE function of the basic SIR model written in density/per-capita form. The function passes vector field and Jacobian functions to the ODE solver.
EpiAware.EpiLatentModels._seir_jac — Method_seir_jac(J, u, p, t)
Internal function for the Jacobian of the basic SEIR model written in density/per-capita form. The function is used to define the ODE problem for the SEIR model. The Jacobian is used to speed up the solution of the ODE problem when using a stiff solver.
EpiAware.EpiLatentModels._seir_vf — Method_seir_vf(du, u, p, t)
Internal function for the vector field of a basic SEIR model written in density/per-capita form. The function is used to define the ODE problem for the SIR model.
EpiAware.EpiLatentModels._sir_function — FunctionInternal function for the ODE function of the basic SIR model written in density/per-capita form. The function passes vector field and Jacobian functions to the ODE solver.
EpiAware.EpiLatentModels._sir_jac — Method_sir_jac(J, u, p, t)
Internal function for the Jacobian of the basic SIR model written in density/per-capita form. The function is used to define the ODE problem for the SIR model. The Jacobian is used to speed up the solution of the ODE problem when using a stiff solver.
EpiAware.EpiLatentModels._sir_vf — Method_sir_vf(du, u, p, t)
Internal function for the vector field of a basic SIR model written in density/per-capita form. The function is used to define the ODE problem for the SIR model.