Public Documentation
Documentation for EpiLatentModels.jl's public interface.
See the Internals section of the manual for internal package docs covering all submodules.
Contents
Index
EpiAware.EpiLatentModelsEpiAware.EpiLatentModels.AREpiAware.EpiLatentModels.BroadcastLatentModelEpiAware.EpiLatentModels.CombineLatentModelsEpiAware.EpiLatentModels.ConcatLatentModelsEpiAware.EpiLatentModels.DiffLatentModelEpiAware.EpiLatentModels.FixedInterceptEpiAware.EpiLatentModels.HierarchicalNormalEpiAware.EpiLatentModels.IIDEpiAware.EpiLatentModels.InterceptEpiAware.EpiLatentModels.NullEpiAware.EpiLatentModels.PrefixLatentModelEpiAware.EpiLatentModels.RandomWalkEpiAware.EpiLatentModels.RecordExpectedLatentEpiAware.EpiLatentModels.RepeatBlockEpiAware.EpiLatentModels.RepeatEachEpiAware.EpiLatentModels.SEIRParamsEpiAware.EpiLatentModels.SIRParamsEpiAware.EpiLatentModels.TransformLatentModelEpiAware.EpiAwareBase.broadcast_ruleEpiAware.EpiAwareBase.broadcast_ruleEpiAware.EpiLatentModels.arimaEpiAware.EpiLatentModels.armaEpiAware.EpiLatentModels.broadcast_dayofweekEpiAware.EpiLatentModels.broadcast_weeklyEpiAware.EpiLatentModels.equal_dimensions
Public API
EpiAware.EpiLatentModels — ModuleModule for defining latent models.
EpiAware.EpiLatentModels.AR — Typestruct AR{D<:Distributions.Sampleable, I<:Distributions.Sampleable, P<:Int64, E<:AbstractTuringLatentModel} <: AbstractTuringLatentModelThe autoregressive (AR) model struct.
Constructors
AR(damp_prior::Sampleable, init_prior::Sampleable; p::Int = 1, ϵ_t::AbstractTuringLatentModel = HierarchicalNormal()): Constructs an AR model with the specified prior distributions for damping coefficients and initial conditions. The order of the AR model can also be specified.AR(; damp_priors::Vector{D} = [truncated(Normal(0.0, 0.05), 0, 1)], init_priors::Vector{I} = [Normal()], ϵ_t::AbstractTuringLatentModel = HierarchicalNormal()) where {D <: Sampleable, I <: Sampleable}: Constructs an AR model with the specified prior distributions for damping coefficients and initial conditions. The order of the AR model is determined by the length of thedamp_priorsvector.
Examples
using Distributions, Turing, EpiAware
ar = AR()
ar
nothing
# outputmdl = generate_latent(ar, 10)
mdl()
nothing
# outputrand(mdl)
nothing
# outputFields
damp_prior::Distributions.Sampleable: Prior distribution for the damping coefficients.init_prior::Distributions.Sampleable: Prior distribution for the initial conditionsp::Int64: Order of the AR model.ϵ_t::AbstractTuringLatentModel: Prior distribution for the error term.
EpiAware.EpiLatentModels.BroadcastLatentModel — Typestruct BroadcastLatentModel{M<:AbstractTuringLatentModel, P<:Integer, B<:AbstractBroadcastRule} <: AbstractTuringLatentModelThe BroadcastLatentModel struct represents a latent model that supports broadcasting of latent periods.
Constructors
BroadcastLatentModel(;model::M; period::Int, broadcast_rule::B): Constructs aBroadcastLatentModelwith the givenmodel,period, andbroadcast_rule.BroadcastLatentModel(model::M, period::Int, broadcast_rule::B): An alternative constructor that allows themodel,period, andbroadcast_ruleto be specified without keyword arguments.
Examples
using EpiAware, Turing
each_model = BroadcastLatentModel(RandomWalk(), 7, RepeatEach())
gen_each_model = generate_latent(each_model, 10)
rand(gen_each_model)
block_model = BroadcastLatentModel(RandomWalk(), 3, RepeatBlock())
gen_block_model = generate_latent(block_model, 10)
rand(gen_block_model)Fields
model::AbstractTuringLatentModel: The underlying latent model.period::Integer: The period of the broadcast.broadcast_rule::AbstractBroadcastRule: The broadcast rule to be applied.
EpiAware.EpiLatentModels.CombineLatentModels — Typestruct CombineLatentModels{M<:(AbstractVector{<:AbstractTuringLatentModel}), P<:(AbstractVector{<:String})} <: AbstractTuringLatentModelThe CombineLatentModels struct.
This struct is used to combine multiple latent models into a single latent model. If a prefix is supplied wraps each model with PrefixLatentModel.
Constructors
CombineLatentModels(models::M, prefixes::P) where {M <: AbstractVector{<:AbstractTuringLatentModel}, P <: AbstractVector{<:String}}: Constructs aCombineLatentModelsinstance with specified models and prefixes, ensuring that there are at least two models and the number of models and prefixes are equal.CombineLatentModels(models::M) where {M <: AbstractVector{<:AbstractTuringLatentModel}}: Constructs aCombineLatentModelsinstance with specified models, automatically generating prefixes for each model. The
automatic prefixes are of the form Combine.1, Combine.2, etc.
Examples
using EpiAware, Distributions
combined_model = CombineLatentModels([Intercept(Normal(2, 0.2)), AR()])
latent_model = generate_latent(combined_model, 10)
latent_model()Fields
models::AbstractVector{<:AbstractTuringLatentModel}: A vector of latent modelsprefixes::AbstractVector{<:String}: A vector of prefixes for the latent models
EpiAware.EpiLatentModels.ConcatLatentModels — Typestruct ConcatLatentModels{M<:(AbstractVector{<:AbstractTuringLatentModel}), N<:Int64, F<:Function, P<:(AbstractVector{<:String})} <: AbstractTuringLatentModelThe ConcatLatentModels struct.
This struct is used to concatenate multiple latent models into a single latent model.
Constructors
ConcatLatentModels(models::M, no_models::I, dimension_adaptor::F, prefixes::P) where {M <: AbstractVector{<:AbstractTuringLatentModel}, I <: Int, F <: Function, P <: AbstractVector{String}}: Constructs aConcatLatentModelsinstance with specified models, number of models, dimension adaptor, and prefixes.ConcatLatentModels(models::M, dimension_adaptor::F; prefixes::P = "Concat." * string.(1:length(models))) where {M <: AbstractVector{<:AbstractTuringLatentModel}, F <: Function}: Constructs aConcatLatentModelsinstance with specified models and dimension adaptor. The number of models is automatically determined as are the prefixes (of the formConcat.1,Concat.2, etc.) by default.ConcatLatentModels(models::M; dimension_adaptor::Function, prefixes::P) where {M <: AbstractVector{<:AbstractTuringLatentModel}, P <: AbstractVector{String}}: Constructs aConcatLatentModelsinstance with specified models, dimension adaptor, prefixes, and automatically determines the number of models.The default dimension adaptor isequal_dimensions. The default prefixes are of the formConcat.1,Concat.2, etc.ConcatLatentModels(; models::M, dimension_adaptor::Function, prefixes::P) where {M <: AbstractVector{<:AbstractTuringLatentModel}, P <: AbstractVector{String}}: Constructs aConcatLatentModelsinstance with specified models, dimension adaptor, prefixes, and automatically determines the number of models. The default dimension adaptor isequal_dimensions. The default prefixes are of the formConcat.1,Concat.2, etc.
Examples
using EpiAware, Distributions
combined_model = ConcatLatentModels([Intercept(Normal(2, 0.2)), AR()])
latent_model = generate_latent(combined_model, 10)
latent_model()Fields
models::AbstractVector{<:AbstractTuringLatentModel}: A vector of latent modelsno_models::Int64: The number of models in the collectiondimension_adaptor::Function: The dimension function for the latent variables. By default this divides the number of latent variables by the number of models and returns a vector of dimensions rounding up the first element and rounding down the rest.prefixes::AbstractVector{<:String}: A vector of prefixes for the latent models
EpiAware.EpiLatentModels.DiffLatentModel — Typestruct DiffLatentModel{M<:AbstractTuringLatentModel, P<:Distributions.Distribution} <: AbstractTuringLatentModelModel the latent process as a d-fold differenced version of another process.
Mathematical specification
Let $\Delta$ be the differencing operator. If $\tilde{Z}_t$ is a realisation of undifferenced latent model supplied to DiffLatentModel, then the differenced process is given by,
\[\Delta^{(d)} Z_t = \tilde{Z}_t, \quad t = d+1, \ldots.\]
We can recover $Z_t$ by applying the inverse differencing operator $\Delta^{-1}$, which corresponds to the cumulative sum operator cumsum in Julia, d-times. The d initial terms $Z_1, \ldots, Z_d$ are inferred.
Constructors
DiffLatentModel(latent_model, init_prior_distribution::Distribution; d::Int)Constructs aDiffLatentModelford-fold differencing withlatent_modelas the undifferenced latent process. All initial terms have common priorinit_prior_distribution.DiffLatentModel(;model, init_priors::Vector{D} where {D <: Distribution})Constructs aDiffLatentModelford-fold differencing withlatent_modelas the undifferenced latent process. Thedinitial terms have priors given by the vectorinit_priors, thereforelength(init_priors)setsd.
Example usage with generate_latent
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_tFields
model::AbstractTuringLatentModel: Underlying latent model for the differenced processinit_prior::Distributions.Distribution: The prior distribution for the initial latent variables.d::Int64: Number of times differenced.
EpiAware.EpiLatentModels.FixedIntercept — Typestruct FixedIntercept{F<:Real} <: AbstractTuringInterceptA variant of the Intercept struct that represents a fixed intercept value for a latent model.
Constructors
FixedIntercept(intercept): Constructs aFixedInterceptinstance with the specified intercept value.FixedIntercept(; intercept): Constructs aFixedInterceptinstance with the specified intercept value using named arguments.
Examples
using EpiAware
fi = FixedIntercept(2.0)
fi_model = generate_latent(fi, 10)
fi_model()
nothing
# outputFields
intercept::Real
EpiAware.EpiLatentModels.HierarchicalNormal — Typestruct HierarchicalNormal{R<:Real, D<:Distributions.Sampleable, M<:Bool} <: AbstractTuringLatentModelThe HierarchicalNormal struct represents a non-centered hierarchical normal distribution.
Constructors
HierarchicalNormal(mean, std_prior): Constructs aHierarchicalNormalinstance with the specified mean and standard deviation prior.HierarchicalNormal(; mean = 0.0, std_prior = truncated(Normal(0,0.1), 0, Inf)): Constructs aHierarchicalNormalinstance with the specified mean and standard deviation prior using named arguments and with default values.HierarchicalNormal(std_prior): Constructs aHierarchicalNormalinstance with the specified standard deviation prior.HierarchicalNormal(mean, std_prior): Constructs aHierarchicalNormalinstance with the specified mean and standard deviation prior.
Examples
using Distributions, Turing, EpiAware
hn = HierarchicalNormal()
mdl = generate_latent(hn, 10)
mdl()
rand(mdl)
nothing
# outputFields
mean::Real: Mean of the normal distribution.std_prior::Distributions.Sampleable: Prior distribution for the standard deviation.add_mean::Bool: Flag to indicate if mean should be added (false when mean = 0)
EpiAware.EpiLatentModels.IID — Typestruct IID{D<:Distributions.Sampleable} <: AbstractTuringLatentModelModel latent process $\epsilon_t$ as independent and identically distributed random variables.
Mathematical specification
The IID process $\epsilon_t$ is specified as a sequence of independent and identically distributed random variables,
\[\epsilon_t \sim \text{Prior}, \quad t = 1, 2, \ldots\]
where Prior is the specified distribution.
Constructors
IID(prior::Distribution = Normal(0, 1)): Create an IID model with the specified prior distribution.
Examples
using EpiAware, Distributions
model = IID(Normal(0, 1))
# output
IID{Normal{Float64}}(Distributions.Normal{Float64}(μ=0.0, σ=1.0))
idd = generate_latent(model, 10)
idd()
nothing
# outputFields
ϵ_t::Distributions.Sampleable
EpiAware.EpiLatentModels.Intercept — Typestruct Intercept{D<:Distributions.Sampleable} <: AbstractTuringInterceptThe Intercept struct is used to model the intercept of a latent process. It broadcasts a single intercept value to a length n latent process.
Constructors
- Intercept(intercept_prior)
- Intercept(; intercept_prior)
Examples
using Distributions, Turing, EpiAware
int = Intercept(Normal(0, 1))
int
# output
Intercept{Normal{Float64}}(Distributions.Normal{Float64}(μ=0.0, σ=1.0))mdl = generate_latent(int, 10)
mdl()
nothingrand(mdl)
nothing
# outputFields
intercept_prior::Distributions.Sampleable: Prior distribution for the intercept.
EpiAware.EpiLatentModels.Null — Typestruct Null <: AbstractTuringLatentModelA null model struct. This struct is used to indicate that no latent variables are to be generated.
Fields
EpiAware.EpiLatentModels.PrefixLatentModel — Typestruct PrefixLatentModel{M<:AbstractTuringLatentModel, P<:String} <: AbstractTuringLatentModelGenerate a latent model with a prefix. A lightweight wrapper around `EpiAwareUtils.prefix_submodel`.
# Constructors
- `PrefixLatentModel(model::M, prefix::P)`: Create a `PrefixLatentModel` with the latent model `model` and the prefix `prefix`.
- `PrefixLatentModel(; model::M, prefix::P)`: Create a `PrefixLatentModel` with the latent model `model` and the prefix `prefix`.
# Examples
```julia
using EpiAware
latent_model = PrefixLatentModel(model = HierarchicalNormal(), prefix = "Test")
mdl = generate_latent(latent_model, 10)
rand(mdl)
```Fields
model::AbstractTuringLatentModel: The latent modelprefix::String: The prefix for the latent model
EpiAware.EpiLatentModels.RandomWalk — Typestruct RandomWalk{D<:Distributions.Sampleable, E<:AbstractTuringLatentModel} <: AbstractTuringLatentModelModel latent process $Z_t$ as a random walk.
Mathematical specification
The random walk $Z_t$ is specified as a parameteric transformation of the white noise sequence $(\epsilon_t)_{t\geq 1}$,
\[Z_t = Z_0 + \sigma \sum_{i = 1}^t \epsilon_t\]
Constructing a random walk requires specifying:
- An
init_prioras a prior for $Z_0$. Default isNormal(). - A
std_priorfor $\sigma$. The default is HalfNormal with a mean of 0.25. - An
ϵ_tprior for the white noise sequence. The default isIID(Normal()).
Constructors
RandomWalk(init_prior::Sampleable, ϵ_t::AbstractTuringLatentModel): Constructs a random walk model with the specified prior distributions for the initial condition and white noise sequence.RandomWalk(; init_prior::Sampleable = Normal(), ϵ_t::AbstractTuringLatentModel = HierarchicalNormal()): Constructs a random walk model with the specified prior distributions for the initial condition and white noise sequence.
Example usage
using Distributions, Turing, EpiAware
rw = RandomWalk()
rw
nothing
# outputmdl = generate_latent(rw, 10)
mdl()
nothing
# outputrand(mdl)
nothingFields
init_prior::Distributions.Sampleableϵ_t::AbstractTuringLatentModel
EpiAware.EpiLatentModels.RecordExpectedLatent — Typestruct RecordExpectedLatent{M<:AbstractTuringLatentModel} <: AbstractTuringLatentModelRecord a variable (using the Turing := syntax) in a latent model.
# Fields
- `model::AbstractTuringLatentModel`: The latent model to dispatch to.
# Constructors
- `RecordExpectedLatent(model::AbstractTuringLatentModel)`: Record the expected latent vector from the model as `exp_latent`.
# Examples
```julia
using EpiAware, Turing
mdl = RecordExpectedLatent(FixedIntercept(0.1))
gen_latent = generate_latent(mdl, 1)
sample(gen_latent, Prior(), 10)
```Fields
model::AbstractTuringLatentModel
EpiAware.EpiLatentModels.RepeatBlock — Typestruct RepeatBlock <: AbstractBroadcastRuleRepeatBlock is a struct that represents a broadcasting rule. It is a subtype of AbstractBroadcastRule.
It repeats the latent process in blocks of size period. An example of this rule is to repeat the latent process in blocks of size 7 to model a weekly process (though for this we also provide the broadcast_weekly helper function).
Examples
using EpiAware
rule = RepeatBlock()
latent = [1, 2, 3, 4, 5]
n = 10
period = 2
broadcast_rule(rule, latent, n, period)Fields
EpiAware.EpiLatentModels.RepeatEach — Typestruct RepeatEach <: AbstractBroadcastRuleRepeatEach is a struct that represents a broadcasting rule. It is a subtype of AbstractBroadcastRule.
It repeats the latent process at each period. An example of this rule is to repeat the latent process at each day of the week (though for this we also provide the dayofweek helper function).
Examples
using EpiAware
rule = RepeatEach()
latent = [1, 2]
n = 10
period = 2
broadcast_rule(rule, latent, n, period)Fields
EpiAware.EpiLatentModels.SEIRParams — Typestruct SEIRParams{P<:SciMLBase.ODEProblem, D<:Distributions.Sampleable, E<:Distributions.Sampleable, F<:Distributions.Sampleable, G<:Distributions.Sampleable} <: AbstractTuringLatentModelA structure representing the SEIR (Susceptible-Exposed-Infectious-Recovered) model and priors for the infectiousness and recovery rate parameters.
Constructors
SEIRParams(; tspan, infectiousness::Distribution, incubation_rate::Distribution, recovery_rate::Distribution, initial_prop_infected::Distribution):
Construct a SEIRParams object with the specified time span for ODE solving, infectiousness prior, incubation rate prior and recovery rate prior.
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 γ.
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)
)
nothingFields
prob::SciMLBase.ODEProblem: The ODE problem instance for the SIR model.infectiousness::Distributions.Sampleable: Prior distribution for the infectiousness parameter.incubation_rate::Distributions.Sampleable: Prior distribution for the incubation rate parameter.recovery_rate::Distributions.Sampleable: Prior distribution for the recovery rate parameter.initial_prop_infected::Distributions.Sampleable: Prior distribution for initial proportion of the population that is infected.
EpiAware.EpiLatentModels.SIRParams — Typestruct SIRParams{P<:SciMLBase.ODEProblem, D<:Distributions.Sampleable, E<:Distributions.Sampleable, F<:Distributions.Sampleable} <: AbstractTuringLatentModelA structure representing the SIR (Susceptible-Infectious-Recovered) model and priors for the infectiousness and recovery rate parameters.
Constructors
SIRParams(; tspan, infectiousness::Distribution, recovery_rate::Distribution, initial_prop_infected::Distribution):
Construct an SIRParams object with the specified time span for ODE solving, infectiousness prior, and recovery rate prior.
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)
)
nothingFields
prob::SciMLBase.ODEProblem: The ODE problem instance for the SIR model.infectiousness::Distributions.Sampleable: Prior distribution for the infectiousness parameter.recovery_rate::Distributions.Sampleable: Prior distribution for the recovery rate parameter.initial_prop_infected::Distributions.Sampleable: Prior distribution for initial proportion of the population that is infected.
EpiAware.EpiLatentModels.TransformLatentModel — Typestruct TransformLatentModel{M<:AbstractTuringLatentModel, F<:Function} <: AbstractTuringLatentModelThe TransformLatentModel struct represents a latent model that applies a transformation function to the latent variables generated by another latent model.
Constructors
TransformLatentModel(model, transform): Constructs aTransformLatentModelinstance with the specified latent model and transformation function.TransformLatentModel(; model, transform): Constructs aTransformLatentModelinstance with the specified latent model and transformation function using named arguments.
Example
using EpiAware, Distributions
trans = TransformLatentModel(Intercept(Normal(2, 0.2)), x -> x .|> exp)
trans_model = generate_latent(trans, 5)
trans_model()Fields
model::AbstractTuringLatentModel: The latent model to transform.transform::Function: The transformation function.
EpiAware.EpiAwareBase.broadcast_rule — Methodbroadcast_rule(_::RepeatBlock, latent, n, period) -> Any
broadcast_rule is a function that applies the RepeatBlock rule to the latent process latent to generate n samples.
Arguments
rule::RepeatBlock: The broadcasting rule.latent::Vector: The latent process.n: The number of samples to generate.period: The period of the broadcast.
Returns
latent: The generated broadcasted latent periods.
EpiAware.EpiAwareBase.broadcast_rule — Methodbroadcast_rule(_::RepeatEach, latent, n, period) -> Any
broadcast_rule is a function that applies the RepeatEach rule to the latent process latent to generate n samples.
Arguments
rule::RepeatEach: The broadcasting rule.latent::Vector: The latent process.n: The number of samples to generate.period: The period of the broadcast.
Returns
latent: The generated broadcasted latent periods.
EpiAware.EpiLatentModels.arima — Methodarima(; ar_init, diff_init, damp, θ, ϵ_t) -> DiffLatentModel
Define an ARIMA model by wrapping define_arma and applying differencing via DiffLatentModel.
Arguments
ar_init: Prior distribution for AR initial conditions. A vector of distributions.diff_init: Prior distribution for differencing initial conditions. A vector of distributions.θ: Prior distribution for MA coefficients. A vector of distributions.damp: Prior distribution for AR damping coefficients. A vector of distributions.ϵ_t: Distribution of the error term. Default isHierarchicalNormal().
Returns
An ARIMA model consisting of AR and MA components with differencing applied.
Example
using EpiAware, Distributions
ARIMA = arima()
arima_model = generate_latent(ARIMA, 10)
arima_model()
nothingEpiAware.EpiLatentModels.arma — Methodarma(
;
init,
damp,
θ,
ϵ_t
) -> Union{AR{D, I, Int64, MA{Distributions.Product{Distributions.Continuous, Distributions.Truncated{Distributions.Normal{Float64}, Distributions.Continuous, Float64, Float64, Float64}, FillArrays.Fill{Distributions.Truncated{Distributions.Normal{Float64}, Distributions.Continuous, Float64, Float64, Float64}, 1, Tuple{Base.OneTo{Int64}}}}, Int64, HierarchicalNormal{Float64, Distributions.Truncated{Distributions.Normal{Float64}, Distributions.Continuous, Float64, Float64, Float64}, Bool}}} where {D<:Distributions.Sampleable, I<:Distributions.Sampleable}, AR{D, I, Int64, MA{Distributions.Product{Distributions.Continuous, Distributions.Truncated{Distributions.Normal{Float64}, Distributions.Continuous, Float64, Float64, Float64}, Vector{Distributions.Truncated{Distributions.Normal{Float64}, Distributions.Continuous, Float64, Float64, Float64}}}, Int64, HierarchicalNormal{Float64, Distributions.Truncated{Distributions.Normal{Float64}, Distributions.Continuous, Float64, Float64, Float64}, Bool}}} where {D<:Distributions.Sampleable, I<:Distributions.Sampleable}}
Define an ARMA model using AR and MA components.
Arguments
init: Prior distribution for AR initial conditions. A vector of distributions.θ: Prior distribution for MA coefficients. A vector of distributions.damp: Prior distribution for AR damping coefficients. A vector of distributions.ϵ_t: Distribution of the error term. Default isHierarchicalNormal().
Returns
An AR model with an MA model as its error term, effectively creating an ARMA model.
Example
using EpiAware, Distributions
ARMA = arma(;
θ = [truncated(Normal(0.0, 0.02), -1, 1)],
damp = [truncated(Normal(0.0, 0.02), 0, 1)]
)
arma_model = generate_latent(ARMA, 10)
arma_model()
nothing
# outputEpiAware.EpiLatentModels.broadcast_dayofweek — Methodbroadcast_dayofweek(
model::AbstractTuringLatentModel;
link
) -> BroadcastLatentModel{TransformLatentModel{M, EpiAware.EpiLatentModels.var"#44#46"}, Int64, RepeatEach} where M<:AbstractTuringLatentModel
Constructs a BroadcastLatentModel appropriate for modelling the day of the week for a given AbstractTuringLatentModel.
Arguments
model::AbstractTuringLatentModel: The latent model to be repeated.link::Function: The link function to transform the latent model before broadcasting
to periodic weekly. Default is x -> 7 * softmax(x) which implements constraint of the sum week effects to be 7.
Returns
BroadcastLatentModel: The broadcast latent model.
EpiAware.EpiLatentModels.broadcast_weekly — Methodbroadcast_weekly(
model::AbstractTuringLatentModel
) -> BroadcastLatentModel{<:AbstractTuringLatentModel, Int64, RepeatBlock}
Constructs a BroadcastLatentModel appropriate for modelling piecewise constant weekly processes for a given AbstractTuringLatentModel.
Arguments
model::AbstractTuringLatentModel: The latent model to be repeated.
Returns
BroadcastLatentModel: The broadcast latent model.
EpiAware.EpiLatentModels.equal_dimensions — Methodequal_dimensions(n::Int64, m::Int64) -> Vector{Int64}
Return a vector of dimensions that are equal or as close as possible, given the total number of elements n and the number of dimensions m. The default dimension adaptor for ConcatLatentModels.
Arguments
n::Int: The total number of elements.m::Int: The number of dimensions.
Returns
dims::AbstractVector{Int}: A vector of dimensions, where the first element is the ceiling ofn / mand the remaining elements are the floor ofn / m.