The MyAwesomePackage Module

Module Index

Detailed API

SmallCouplingDynamicCavity.EpidemicModelType
EpidemicModel

Epidemic model containing all the informations of the system.

Fields

  • Disease::SmallCouplingDynamicCavity.InfectionModel: Infection model. Currently are implemented SI, SIR, SIS and SIRS infection models.

  • G::Union{Vector{<:Graphs.AbstractGraph}, var"#s30"} where var"#s30"<:Graphs.AbstractGraph: Contact graph. It can be either an AbstractGraph (contact graph constant over time) or a Vector of AbstractGraph (time varying contact graph).

  • N::Int64: Number of nodes of the contact graph.

  • T::Int64: Number of time steps.

  • ν::Array{Float64, 3}: Infection couplings. It is a NVxNVx(T+1) Array where νᵗᵢⱼ=log(1-λᵗᵢⱼ), with λᵗᵢⱼ being the infection probability from individual i to individual j at time t.

  • obsmat::Matrix{Int8}: Observations matrix. It is a NVx(T+1) Matrix, where obsᵗᵢ is the observation of individual i at time t: it is equal to -1.0 if not observed, 0.0 if S, 1.0 if I, 2.0 if R (only for SIR and SIRS).

source
SmallCouplingDynamicCavity.EpidemicModelMethod
EpidemicModel(
    infectionmodel::TI, 
    G::TG, T::Int, 
    ν::Array{Float64, 3}, 
    obs::Matrix{Int8}) where {TI<:InfectionModel,TG<:Vector{<:AbstractGraph}}

Defines the epidemic model.

This function defines an epidemic model based on the provided infection model, contact graph, time steps, infection couplings, and observation matrix.

# Arguments
- `infectionmodel`: The infection model. Currently implemented models include SI, SIR, SIS, and SIRS infection models.
- `G`: The contact graph. It should be a T+1 vector of AbstractGraph representing the time-varying contact graph.
- `T`: The number of time-steps.
- `ν`: The infection couplings. It should be a 3-dimensional array of size NVxNVx(T+1), where νᵗᵢⱼ=log(1-λᵗᵢⱼ), with λᵗᵢⱼ being the infection probability from individual i to individual j at time t.
- `obs`: The observations obsmatrix. It should be a NVx(T+1) matrix, where obsᵗᵢ is the observation of individual i at time t: it is equal to -1.0 if not observed, 0.0 if S, 1.0 if I, 2.0 if R (only for SIR and SIRS).

# Returns
- `EpidemicModel`: An [`EpidemicModel`](@ref) object representing the defined epidemic model.
source
SmallCouplingDynamicCavity.EpidemicModelMethod
EpidemicModel(
    infectionmodel::TI, 
    G::TG, T::Int, 
    ν::Array{Float64, 3}, 
    obs::Matrix{Int8}) where {TI<:InfectionModel,TG<:AbstractGraph}

Defines the epidemic model.

This function defines an epidemic model based on the provided infection model, contact graph, time steps, infection couplings, and observation matrix.

Arguments

  • infectionmodel: The infection model. Currently implemented models include SI, SIR, SIS, and SIRS infection models.
  • G: The contact graph. It should be an AbstractGraph representing the contact graph, which is time-varying.
  • T: The number of time-steps.
  • ν: The infection couplings. It should be a 3-dimensional array of size NVxNVx(T+1), where νᵗᵢⱼ=log(1-λᵗᵢⱼ), with λᵗᵢⱼ being the infection probability from individual i to individual j at time t.
  • obs: The observations matrix. It should be a NVx(T+1) matrix, where obsᵗᵢ is the observation of individual i at time t: it is equal to -1.0 if not observed, 0.0 if S, 1.0 if I, 2.0 if R (only for SIR and SIRS).

Returns

  • EpidemicModel: An EpidemicModel object representing the defined epidemic model.
source
SmallCouplingDynamicCavity.EpidemicModelMethod
EpidemicModel(
    infectionmodel::TI, 
    G::TG, T::Int, 
    ν::Array{Float64, 3}) where {TI<:InfectionModel,TG<:Vector{<:AbstractGraph}}

Define an epidemic model.

This function defines an epidemic model based on the provided infection model, time-varying contact graph, time steps, and infection couplings. It initializes the observation matrix with zeros.

Arguments

  • infectionmodel: The infection model. Currently implemented models include SI, SIR, SIS, and SIRS infection models.
  • G: The contact graph. It should be a T+1 vector of AbstractGraph representing the time-varying contact graph.
  • T: The number of time-steps.
  • ν: The infection couplings. It should be a 3-dimensional array of size NVxNVx(T+1), where νᵗᵢⱼ=log(1-λᵗᵢⱼ), with λᵗᵢⱼ being the infection probability from individual i to individual j at time t.

Returns

  • EpidemicModel: An EpidemicModel object representing the defined epidemic model.
source
SmallCouplingDynamicCavity.EpidemicModelMethod
EpidemicModel(
    infectionmodel::TI, 
    G::TG, T::Int, 
    ν::Array{Float64, 3}) where {TI<:InfectionModel,TG<:AbstractGraph}

Define an epidemic model.

This function defines an epidemic model based on the provided infection model, contact graph, time steps, and infection couplings. It initializes the observation matrix with zeros.

Arguments

  • infectionmodel: The infection model. Currently implemented models include SI, SIR, SIS, and SIRS infection models.
  • G: The contact graph. It should be an AbstractGraph representing the contact graph, which is constant over time.
  • T: The number of time-steps.
  • ν: The infection couplings. It should be a 3-dimensional array of size NVxNVx(T+1), where νᵗᵢⱼ=log(1-λᵗᵢⱼ), with λᵗᵢⱼ being the infection probability from individual i to individual j at time t.

Returns

  • EpidemicModel: An EpidemicModel object representing the defined epidemic model.
source
SmallCouplingDynamicCavity.MarginalType
Marginal

Marginals pᵢ(xᵢ) and μᵢ.

Fields

  • i::Int64: Index of the node i.

  • m::Matrix{Float64}: (nstates)x(T+1) Matrix of marginals over time, where nstates is the number of states that the infection model has.

  • μ::Vector{Float64}: T+1 Vector of marginals μᵢ over time.

source
SmallCouplingDynamicCavity.MessageType
Message

Cavity messages mᵢⱼ and μᵢⱼ.

Fields

  • i::Int64: Index of the node i.

  • j::Int64: Index of the node j.

  • m::Vector{Float64}: T+1 Vector of messages mᵢⱼ over time.

  • μ::Vector{Float64}: T+1 Vector of messages μᵢⱼ over time.

source
SmallCouplingDynamicCavity.NodeType
Node

Type containing all the informations of a node.

Fields

  • i::Int64: Index of the node.

  • ∂::Vector{Int64}: List of neighbours. If the underlying contact graph is varying in time it is the union of all the neighbours over time.

  • ∂_idx::Dict{Int64, Int64}: Only for developers.

  • marg::SmallCouplingDynamicCavity.Marginal: Marginals of the node. It is a Marginal type.

  • cavities::Vector{SmallCouplingDynamicCavity.Message}: Cavities messages entering into the node from its neigbours. It is a vector of Message, each one corresponding to a neighbour with the same order of ∂.

  • νs::Vector{Vector{Float64}}: Infection couplings of the neighbours against the node.

  • obs::Matrix{Float64}: Observation probability matrix.

source
SmallCouplingDynamicCavity.ROC_curveMethod
ROC_curve(marg::Vector{Float64}, x::Vector{TI}) where {TI<:Integer}

Compute the Receiver Operating Characteristic (ROC) curve and the Area Under the Curve (AUC) based on inferred marginal probabilities and true configurations.

This function computes the ROC curve by sorting the marginal probabilities (marg) in decreasing order and sorting the true configuration (x) correspondingly. It then iterates over the sorted values to calculate the True Positive Rate (TPR) and False Positive Rate (FPR) at different thresholds, and computes the AUC.

Arguments

  • marg: A vector containing the marginal probabilities of each node being infected at a fixed time.
  • x: The true configuration of the nodes at the same fixed time.

Returns

  • fp_rates: Array of false positive rates corresponding to each threshold.
  • tp_rates: Array of true positive rates corresponding to each threshold.
  • auc: The Area Under the Curve (AUC) of the ROC curve.
source
SmallCouplingDynamicCavity.bethe_latticeMethod
bethe_lattice(z::Int, tmax::Int, startfrom1::Bool)

Generate a Bethe lattice (tree) with a specified degree and depth.

Arguments

  • z::Int: The degree of the Bethe lattice.
  • tmax::Int: The maximum depth of the Bethe lattice.
  • startfrom1::Bool: If true, the center of the tree is vertex 1.

Returns

  • V::Vector{Int}: A list of vertices in the Bethe lattice.
  • E::Vector{Vector{Int}}: A list of edges in the Bethe lattice.

Description

This function generates a Bethe lattice (tree) with a specified degree (z) and maximum depth (tmax). The Bethe lattice is constructed by iteratively adding nodes and edges according to the specified parameters.

If startfrom1 is true, the center of the tree is vertex 1. Otherwise, the tree is constructed starting from vertex 0.

The function returns a tuple where the first element (V) is a list of vertices, and the second element (E) is a list of edges in the Bethe lattice.

source
SmallCouplingDynamicCavity.run_SCDC!Method
run_SCDC!(
    nodes::Vector{Node{TI,TG}}, 
    model::EpidemicModel{TI,TG}, 
    γ::Float64, 
    maxiter::Vector{Int64}, 
    epsconv::Float64, 
    damp::Vector{Float64}; 
    μ_cutoff::Float64 = -Inf, 
    n_iter_nc::Int64 = 1, 
    damp_nc::Float64 = 0.0, 
    callback::Function=(x...) -> nothing, 
    rng::AbstractRNG=Xoshiro(1234)) where {TI<:InfectionModel,TG<:Union{<:AbstractGraph,Vector{<:AbstractGraph}}}

Run the SCDC algorithm for epidemic modeling. The algorithm resumes the message-passing iterations from the current state of the nodes.

This function performs SCDC inference on the specified epidemic model, using the provided evidence (likelihood) probability function, and other parameters such as the probability of being a patient zero, maximum number of iterations, convergence threshold, damping factor, etc. It iteratively updates cavity messages until convergence or until the maximum number of iterations is reached. It implements a dumping schedule for the damping factor, where the dumping factor is changed after a certain number of iterations, specified by the maxiter and damp arguments.

Arguments

  • nodes::Vector{Node{TI,TG}}: Vector of nodes in the epidemic model.
  • model::EpidemicModel{TI,TG}: The epidemic model to be used.
  • γ::Float64: A parameter for the algorithm (e.g., infection rate).
  • maxiter::Vector{Int64}: Maximum number of iterations for the algorithm.
  • epsconv::Float64: Convergence threshold for the algorithm.
  • damp::Vector{Float64}: Damping factors for the algorithm.

Keyword Arguments

  • μ_cutoff::Float64: Cutoff value for some parameter μ (default is -Inf).
  • n_iter_nc::Int64: Number of iterations for non-converging cases (default is 1).
  • damp_nc::Float64: Damping factor for non-converging cases (default is 0.0).
  • callback::Function: Callback function to be called during iterations (default does nothing).
  • rng::AbstractRNG: Random number generator (default is Xoshiro).
source
SmallCouplingDynamicCavity.run_SCDC!Method
run_SCDC!(nodes::Vector{Node{TI,TG}}, model::EpidemicModel{TI,TG}, γ::Float64, maxiter::Vector{Int64}, epsconv::Float64, damp::Vector{Float64}; μ_cutoff::Float64 = -Inf, n_iter_nc::Int64 = 1, damp_nc::Float64 = 0.0, callback::Function=(x...) -> nothing, rng::AbstractRNG=Xoshiro(1234))

Run the SCDC algorithm for epidemic modeling. The algorithm resumes the message-passing iterations from the current state of the nodes.

This function performs SCDC inference on the specified epidemic model, using the provided evidence (likelihood) probability function, and other parameters such as the probability of being a patient zero, maximum number of iterations, convergence threshold, damping factor, etc. It iteratively updates cavity messages until convergence or until the maximum number of iterations is reached. It implements a dumping schedule for the damping factor, where the dumping factor is changed after a certain number of iterations, specified by the maxiter and damp arguments.

Arguments

  • nodes::Vector{Node{TI,TG}}: Vector of nodes in the epidemic model.
  • model::EpidemicModel{TI,TG}: The epidemic model to be used.
  • γ::Float64: A parameter for the algorithm (e.g., infection rate).
  • maxiter::Vector{Int64}: Maximum number of iterations for the algorithm.
  • epsconv::Float64: Convergence threshold for the algorithm.
  • damp::Vector{Float64}: Damping factors for the algorithm.

Keyword Arguments

  • μ_cutoff::Float64: Cutoff value for some parameter μ (default is -Inf).
  • n_iter_nc::Int64: Number of iterations for non-converging cases (default is 1).
  • damp_nc::Float64: Damping factor for non-converging cases (default is 0.0).
  • callback::Function: Callback function to be called during iterations (default does nothing).
  • rng::AbstractRNG: Random number generator (default is Xoshiro).
source
SmallCouplingDynamicCavity.run_SCDCMethod
run_SCDC(
    model::EpidemicModel{TI,TG},
    obsprob::Function,
    γ::Float64,
    maxiter::Int64,
    epsconv::Float64,
    damp::Float64;
    μ_cutoff::Float64 = -Inf,
    n_iter_nc::Int64 = 1,
    damp_nc::Float64 = 0.0,
    callback::Function=(x...) -> nothing
    rng::AbstractRNG=Xoshiro(1234)) where {TI<:InfectionModel,TG<:Union{<:AbstractGraph,Vector{<:AbstractGraph}}}

Runs the Small Coupling Dynamic Cavity (SCDC) inference algorithm.

This function performs SCDC inference on the specified epidemic model, using the provided evidence (likelihood) probability function, and other parameters such as the probability of being a patient zero, maximum number of iterations, convergence threshold, damping factor, etc. It iteratively updates cavity messages until convergence or until the maximum number of iterations is reached.

Arguments

  • model: An EpidemicModel representing the epidemic model.
  • obsprob: A function representing the evidence (likelihood) probability p(O|x) of an observation O given the planted state x.
  • γ: The probability of being a patient zero.
  • maxiter: The maximum number of iterations.
  • epsconv: The convergence threshold of the algorithm.
  • damp: The damping factor of the algorithm.
  • μ_cutoff: (Optional) Lower cut-off for the values of μ.
  • n_iter_nc: (Optional) Number of iterations for non-converged messages. The messages are averaged over this number of iterations.
  • damp_nc: (Optional) Damping factor for non-converged messages.
  • callback: (Optional) A callback function to monitor the progress of the algorithm.
  • rng: (Optional) Random number generator.

Returns

  • nodes: An array of Node objects representing the updated node states after inference.
source
SmallCouplingDynamicCavity.run_SCDCMethod
run_SCDC(
    model::EpidemicModel{TI,TG},
    obsprob::Function,
    γ::Float64,
    maxiter::Vector{Int64},
    epsconv::Float64,
    damp::Vector{Float64};
    μ_cutoff::Float64 = -Inf,
    n_iter_nc::Int64 = 1,
    damp_nc::Float64 = 0.0,
    callback::Function=(x...) -> nothing,
    rng::AbstractRNG=Xoshiro(1234)) where {TI<:InfectionModel,TG<:Union{<:AbstractGraph,Vector{<:AbstractGraph}}}

Runs the Small Coupling Dynamic Cavity (SCDC) inference algorithm.

This function performs SCDC inference on the specified epidemic model, using the provided evidence (likelihood) probability function, and other parameters such as the probability of being a patient zero, maximum number of iterations, convergence threshold, damping factor, etc. It iteratively updates cavity messages until convergence or until the maximum number of iterations is reached. It implements a dumping schedule for the damping factor, where the dumping factor is changed after a certain number of iterations, specified by the maxiter and damp arguments.

Arguments

  • model: An EpidemicModel representing the epidemic model.
  • obsprob: A function representing the evidence (likelihood) probability p(O|x) of an observation O given the planted state x.
  • γ: The probability of being a patient zero.
  • maxiter: Vector of maximum number of iterations for each damping factor.
  • epsconv: The convergence threshold of the algorithm.
  • damp: Vector of damping factors used in the damping schedule.
  • μ_cutoff: (Optional) Lower cut-off for the values of μ.
  • n_iter_nc: (Optional) Number of iterations for non-converged messages. The messages are averaged over this number of iterations.
  • damp_nc: (Optional) Damping factor for non-converged messages.
  • callback: (Optional) A callback function to monitor the progress of the algorithm.
  • rng: (Optional) Random number generator.

Returns

  • nodes: An array of Node objects representing the updated node states after inference.
source
SmallCouplingDynamicCavity.SIType
struct SI <: InfectionModel
    εᵢᵗ::Array{Float64, 2} # Autoinfection probabilities
end

The SI struct represents the SI (Susceptible-Infected) infection model.

Fields

  • εᵢᵗ: An NVxT array representing the self-infection probabilities over time, where NV is the number of nodes and T is the number of time-steps. Each element εᵢᵗ[i, t] denotes the probability of node i infecting itself at time t.
source
SmallCouplingDynamicCavity.sim_epidemicsMethod
sim_epidemics(
    model::EpidemicModel{SI,TG};
    patient_zero::Union{Vector{Int},Nothing}=nothing,
    γ::Union{Float64,Nothing}=nothing) where {TG<:Union{<:AbstractGraph,Vector{<:AbstractGraph}}}

Simulates an epidemic outbreak using the SI (Susceptible-Infectious) model.

Arguments

  • model: The SI epidemic model, encapsulating information about the infection dynamics, contact graph, and other parameters.
  • patient_zero: (Optional) A vector specifying the indices of initial infected individuals. If not provided (default nothing), patient zero is selected randomly based on the probability γ.
  • γ: (Optional) The probability of being a patient zero. If patient_zero is not specified and γ is provided, patient zero is chosen randomly with probability γ. If both patient_zero and γ are not provided (default nothing), patient zero is selected randomly with equal probability for each individual.

Returns

  • A matrix representing the epidemic outbreak configuration over time. Each row corresponds to a node, and each column represents a time step. The values in the matrix indicate the state of each node at each time step: 0.0 for Susceptible (S) and 1.0 for Infected (I).
source
SmallCouplingDynamicCavity.SIRType
struct SIR <: InfectionModel
    εᵢᵗ::Array{Float64, 2} # Autoinfection probabilities
    rᵢᵗ::Array{Float64, 2} # Recovery probabilities
end

The SIR struct represents the SIR (Susceptible-Infected-Recovered) infection model.

Fields

  • εᵢᵗ: An NVxT array representing the self-infection probabilities over time, where NV is the number of nodes and T is the number of time-steps. Each element εᵢᵗ[i, t] denotes the probability of node i infecting itself at time t.
  • rᵢᵗ: An NVxT array representing the recovery probabilities over time, where NV is the number of nodes and T is the number of time-steps. Each element rᵢᵗ[i, t] denotes the probability of node i recovering from infection at time t.
source
SmallCouplingDynamicCavity.sim_epidemicsMethod
sim_epidemics(
    model::EpidemicModel{SIR,TG};
    patient_zero::Union{Vector{Int},Nothing}=nothing,
    γ::Union{Float64,Nothing}=nothing) where {TG<:Union{<:AbstractGraph,Vector{<:AbstractGraph}}}

Simulates an epidemic outbreak using the SIR (Susceptible-Infectious-Recovered) model.

Arguments

  • model: The SIR epidemic model, encapsulating information about the infection dynamics, contact graph, and other parameters.
  • patient_zero: (Optional) A vector specifying the indices of initial infected individuals. If not provided (default nothing), patient zero is selected randomly based on the probability γ.
  • γ: (Optional) The probability of being a patient zero. If patient_zero is not specified and γ is provided, patient zero is chosen randomly with probability γ. If both patient_zero and γ are not provided (default nothing), patient zero is selected randomly with equal probability for each individual.

Returns

  • A matrix representing the epidemic outbreak configuration over time. Each row corresponds to a node, and each column represents a time step. The values in the matrix indicate the state of each node at each time step: 0.0 for Susceptible (S), 1.0 for Infected (I), and 2.0 for Recovered (R).
source
SmallCouplingDynamicCavity.SISType
struct SIS <: InfectionModel
    εᵢᵗ::Array{Float64, 2} # Autoinfection probabilities
    rᵢᵗ::Array{Float64, 2} # Recovery probabilities
end

The SIS struct represents the SIS (Susceptible-Infected-Susceptible) infection model.

Fields

  • εᵢᵗ: An NVxT array representing the self-infection probabilities over time, where NV is the number of nodes and T is the number of time-steps. Each element εᵢᵗ[i, t] denotes the probability of node i infecting itself at time t.
  • rᵢᵗ: An NVxT array representing the recovery probabilities over time, where NV is the number of nodes and T is the number of time-steps. Each element rᵢᵗ[i, t] denotes the probability of node i recovering from infection at time t.
source
SmallCouplingDynamicCavity.sim_epidemicsMethod
sim_epidemics(
    model::EpidemicModel{SIS,TG};
    patient_zero::Union{Vector{Int},Nothing}=nothing,
    γ::Union{Float64,Nothing}=nothing) where {TG<:Union{<:AbstractGraph,Vector{<:AbstractGraph}}}

Simulates an epidemic outbreak using the SIS (Susceptible-Infectious-Susceptible) model.

Arguments

  • model: The SIS epidemic model, encapsulating information about the infection dynamics, contact graph, and other parameters.
  • patient_zero: (Optional) A vector specifying the indices of initial infected individuals. If not provided (default nothing), patient zero is selected randomly based on the probability γ.
  • γ: (Optional) The probability of being a patient zero. If patient_zero is not specified and γ is provided, patient zero is chosen randomly with probability γ. If both patient_zero and γ are not provided (default nothing), patient zero is selected randomly with equal probability for each individual.

Returns

  • A matrix representing the epidemic outbreak configuration over time. Each row corresponds to a node, and each column represents a time step. The values in the matrix indicate the state of each node at each time step: 0.0 for Susceptible (S) and 1.0 for Infected (I).
source
SmallCouplingDynamicCavity.SIRSType
struct SIRS <: InfectionModel
    εᵢᵗ::Array{Float64, 2} # Autoinfection probabilities
    rᵢᵗ::Array{Float64, 2} # Recovery probabilities
    σᵢᵗ::Array{Float64, 2} # Loss of immunity probabilities
end

The SIRS struct represents the SIRS (Susceptible-Infected-Recovered-Susceptible) infection model.

Fields

  • εᵢᵗ: An NVxT array representing the self-infection probabilities over time, where NV is the number of nodes and T is the number of time-steps. Each element εᵢᵗ[i, t] denotes the probability of node i infecting itself at time t.
  • rᵢᵗ: An NVxT array representing the recovery probabilities over time, where NV is the number of nodes and T is the number of time-steps. Each element rᵢᵗ[i, t] denotes the probability of node i recovering from infection at time t.
  • σᵢᵗ: An NVxT array representing the loss of immunity probabilities over time, where NV is the number of nodes and T is the number of time-steps. Each element σᵢᵗ[i, t] denotes the probability of node i losing immunity and becoming susceptible again at time t.
source
SmallCouplingDynamicCavity.sim_epidemicsMethod
sim_epidemics(
    model::EpidemicModel{SIRS,TG};
    patient_zero::Union{Vector{Int},Nothing}=nothing,
    γ::Union{Float64,Nothing}=nothing) where {TG<:Union{<:AbstractGraph,Vector{<:AbstractGraph}}}

Simulates an epidemic outbreak using the SIRS (Susceptible-Infectious-Recovered-Susceptible) model.

Arguments

  • model: The SIRS epidemic model, encapsulating information about the infection dynamics, contact graph, and other parameters.
  • patient_zero: (Optional) A vector specifying the indices of initial infected individuals. If not provided (default nothing), patient zero is selected randomly based on the probability γ.
  • γ: (Optional) The probability of being a patient zero. If patient_zero is not specified and γ is provided, patient zero is chosen randomly with probability γ. If both patient_zero and γ are not provided (default nothing), patient zero is selected randomly with equal probability for each individual.

Returns

  • A matrix representing the epidemic outbreak configuration over time. Each row corresponds to a node, and each column represents a time step. The values in the matrix indicate the state of each node at each time step: 0.0 for Susceptible (S), 1.0 for Infected (I), and 2.0 for Recovered (R).
source