GaussianExpansionCavityMethod

Documentation for GaussianExpansionCavityMethod.

GaussianExpansionCavityMethod.BMModelType
BMModel

A type representing the Bouchaud-Mezard model.

Fields

  • N::Integer: The number of sites.
  • K::Union{Integer, Float64}: The average number of neighbors.
  • lambdas::Vector{Float64}: The local drift terms. lambdas[i] = J |∂i|, where |∂i| is the degree of site i and J is the uniform coupling strength.
  • J::SparseMatrixCSC{Float64, Integer}: The sparse coupling matrix.
  • sigma::Float64: The multiplicative noise strength.
source
GaussianExpansionCavityMethod.BMModelEnsembleType
BMModelEnsemble

A type representing the ensemble of disordered Bouchaud-Mezard model.

Fields

  • N::Integer: The number of sites.
  • K::Union{Integer, Float64}: The average number of neighbors.
  • lambdas::Vector{Float64}: The local drift terms. lambdas[i] = J |∂i|, where |∂i| is the degree of site i and J is the uniform coupling strength.
  • gen_J::Function: The function to generate the sparse coupling matrix.
  • J_params::Vector{Float64}: The parameters for the coupling matrix generator.
  • sigma::Float64: The multiplicative noise strength.
source
GaussianExpansionCavityMethod.CavityType
Cavity

A structure representing a cavity in a graph with its associated parameters.

Fields

  • i::Int: The index of the node.
  • j::Int: The index of the neighbor.
  • mu::OffsetVector{Float64, Vector{Float64}}: The mean vector.
  • C::OffsetMatrix{Float64, Matrix{Float64}}: The autocorellation matrix.
  • R::OffsetMatrix{Float64, Matrix{Float64}}: The response matrix.

Description

The Cavity structure represents a cavity in a graph, which is a subgraph formed by removing a node and its associated edges. Each cavity has an index i, a neighbor index j, and associated cavity fields such as the mean vector mu, the autocorrelation matrix C, and the response matrix R.

Methods

  • Cavity(i::Int, j::Int, T::Int): Constructs a cavity with index i, neighbor index j, and T timesteps.
source
GaussianExpansionCavityMethod.CavityEQType
CavityEQ

A structure representing a equilibrium cavity in a graph with its associated parameters.

Fields

  • i::Int: The index of the node.
  • j::Int: The index of the neighbor.
  • C::OffsetVector{Float64, Vector{Float64}}: The autocorellation vector.

Description

The CavityEQ structure represents a equilibrium cavity in a graph, which is a subgraph formed by removing a node and its associated edges. Each cavity has an index i, a neighbor index j, and associated cavity field such as the autocorrelation vector C.

Methods

  • CavityEQ(i::Int, j::Int, T::Int): Constructs a equilibrium cavity with index i, neighbor index j, and T timesteps.
source
GaussianExpansionCavityMethod.MarginalType
Marginal

A structure representing a marginal in a graph with its associated parameters.

Fields

  • i::Int: The index of the node.
  • mu::OffsetVector{Float64, Vector{Float64}}: The mean vector.
  • C::OffsetMatrix{Float64, Matrix{Float64}}: The autocorellation matrix.
  • R::OffsetMatrix{Float64, Matrix{Float64}}: The response matrix.

Description

The Marginal structure represents a marginal in a graph, which is a subgraph formed by removing a node and its associated edges. Each marginal has an index i, and associated marginal fields such as the mean vector mu, the autocorrelation matrix C, and the response matrix R.

Methods

  • Marginal(i::Int, T::Int): Constructs a marginal with index i and T timesteps.
source
GaussianExpansionCavityMethod.MarginalEQType
MarginalEQ

A structure representing a equilibrium marginal in a graph with its associated parameters.

Fields

  • i::Int: The index of the node.
  • mu::Float64: The mean at infinity.
  • C::OffsetVector{Float64, Vector{Float64}}: The autocorellation vector.
  • R::OffsetVector{Float64, Vector{Float64}}: The response vector.

Description

The MarginalEQ structure represents a equilibrium marginal in a graph, which is a subgraph formed by removing a node and its associated edges. Each marginal has an index i, and associated marginal fields such as the mean at infinity mu, the autocorrelation vector C, and the response vector R.

Methods

  • MarginalEQ(i::Int, T::Int): Constructs an equilibrium marginal with index i and T timesteps.
source
GaussianExpansionCavityMethod.NodeType
Node

A structure representing a node in a graph with its neighbors and associated cavities and marginal.

Fields

  • i::Int: The index of the node.
  • neighs::Vector{Int}: The indices of the neighbors.
  • neighs_idxs::Dict{Int, Int}: A dictionary mapping neighbor indices to their positions in the neighs vector.
  • cavs::Vector{Cavity}: The cavities associated with the node.
  • marg::Marginal: The marginal associated with the node.
  • summu::OffsetVector{Float64, Vector{Float64}}: Internal variable.
  • sumC::OffsetMatrix{Float64, Matrix{Float64}}: Internal variable.
  • sumR::OffsetMatrix{Float64, Matrix{Float64}}: Internal variable.

Description

The Node structure represents a node in a graph, along with its neighbors and associated cavities and marginal. Each node has an index i, a vector of neighbor indices neighs, a vector of cavities cavs and a marginal marg.

Methods

  • Node(i::Int, neighs::Vector{Int}, T::Int): Constructs a node with index i, neighbors neighs, and T timesteps.
source
GaussianExpansionCavityMethod.NodeEQType
NodeEQ

A structure representing a equilibrium node in a graph with its neighbors and associated cavities and marginal.

Fields

  • i::Int: The index of the node.
  • neighs::Vector{Int}: The indices of the neighbors.
  • neighs_idxs::Dict{Int, Int}: A dictionary mapping neighbor indices to their positions in the neighs vector.
  • cavs::Vector{CavityEQ}: The equilibrium cavities associated with the node.
  • marg::MarginalEQ: The equilibrium marginal associated with the node.
  • sumC::OffsetVector{Float64, Vector{Float64}}: Internal variable.
  • sumdiffC::OffsetVector{Float64, Vector{Float64}}: Internal variable.

Description

The NodeEQ structure represents a equilibrium node in a graph, along with its neighbors and associated equilibrium cavities and marginal. Each node has an index i, a vector of neighbor indices neighs, a vector of equilibrium cavities cavs and a equilibrium marginal marg.

Methods

  • NodeEQ(i::Int, neighs::Vector{Int}, T::Int): Constructs a node with index i, neighbors neighs, and T timesteps.
source
GaussianExpansionCavityMethod.OUModelType
OUModel

A structure representing a linearly-coupled Ornstein-Uhlenbeck model with its associated parameters.

Fields

  • N::Int: The number of nodes in the model.
  • K::Union{Int, Float64}: The average number of neighbors.
  • J::SparseMatrixCSC{Float64, Int}: The coupling matrix.
  • lambdas::Vector{Float64}: The local decay rates.
  • D::Float64: The diffusion coefficient.

Description

The OUModel structure represents a linearly-coupled Ornstein-Uhlenbeck model, which is a stochastic process used to model the dynamics of a system with local decay rates and coupling between nodes. Each model has a number of nodes N, local decay rates lambdas, a coupling matrix J, and a diffusion coefficient D.

Methods

source
GaussianExpansionCavityMethod.OUModelEnsembleType
OUModelEnsemble

A type representing the ensemble of disordered Ornstein-Uhlenbeck model.

Fields

  • N::Integer: The number of sites.
  • K::Union{Integer, Float64}: The average number of neighbors.
  • gen_J::Function: The function to generate the sparse coupling matrix.
  • J_params::Vector{Float64}: The parameters for the coupling matrix generator.
  • lambdas::Vector{Float64}: The local decay rates.
  • D::Float64: The diffusion coefficient.
source
GaussianExpansionCavityMethod.Phi4ModelType
Phi4Model

A type representing the Phi^4 model.

Fields

  • N::Integer: The number of sites.
  • K::Union{Integer, Float64}: The average number of neighbors.
  • J::SparseMatrixCSC{Float64, Integer}: The sparse coupling matrix.
  • lambdas::Vector{Float64}: The on-site linear terms.
  • D::Float64: The noise strength.
  • u::Float64: The cubic perturbative constant.
source
GaussianExpansionCavityMethod.Phi4ModelEnsembleType
Phi4ModelEnsemble

A type representing the ensemble of disordered Phi^4 model.

Fields

  • N::Integer: The number of sites.
  • K::Union{Integer, Float64}: The average number of neighbors.
  • gen_J::Function: The function to generate the sparse coupling matrix.
  • J_params::Vector{Float64}: The parameters for the coupling matrix generator.
  • lambdas::Vector{Float64}: The on-site linear terms.
  • D::Float64: The noise strength.
  • u::Float64: The cubic perturbative constant.
source
GaussianExpansionCavityMethod.TwoSpinModelType
TwoSpinModel

A type representing the spherical 2-Spin model.

Fields

  • N::Integer: The number of sites.
  • K::Union{Integer, Float64}: The average number of neighbors.
  • J::SparseMatrixCSC{Float64, Integer}: The sparse coupling matrix.
  • D::Float64: The noise strength.
source
GaussianExpansionCavityMethod.TwoSpinModelEnsembleType
TwoSpinModelEnsemble

A type representing the ensemble of disordered spherical 2-Spin model.

Fields

  • N::Integer: The number of sites.
  • K::Union{Integer, Float64}: The average number of neighbors.
  • gen_J::Function: The function to generate the sparse coupling matrix.
  • J_params::Vector{Float64}: The parameters for the coupling matrix generator.
  • D::Float64: The noise strength.
source
GaussianExpansionCavityMethod.BMModelRRGMethod
BMModelRRG(N, K, J, sigma)

Construct a Bouchaud-Mezard model with a random regular graph coupling matrix with ferromagnetic interactions.

Arguments

  • N::Integer: The number of sites.
  • K::Union{Int, Float64}: The average number of neighbors.
  • J::Float64: The coupling strength.
  • sigma::Float64: The multiplicative noise strength.

Returns

  • BMModelEnsemble: The Bouchaud-Mezard model.
source
GaussianExpansionCavityMethod.OUModelRRGMethod
OUModelRRG(N, K, J, lambdas, D, x0_min, x0_max)

Construct a Ornstein-Uhlenbeck model with a random regular graph coupling matrix with ferromagnetic interactions.

Arguments

  • N::Integer: The number of sites.
  • K::Union{Int, Float64}: The average number of neighbors.
  • J::Float64: The coupling strength.
  • lambdas::Vector{Float64}: The local decay rates.
  • D::Float64: The diffusion coefficient.
  • x0_min::Float64: The minimum initial condition.
  • x0_max::Float64: The maximum initial condition.

Returns

  • OUModelEnsemble: The Ornstein-Uhlenbeck model.
source
GaussianExpansionCavityMethod.OUModelRRGMethod
OUModelRRG(N, K, J, lambda, D, x0_min, x0_max)

Construct a Ornstein-Uhlenbeck model with a random regular graph coupling matrix with ferromagnetic interactions.

Arguments

  • N::Integer: The number of sites.
  • K::Union{Int, Float64}: The average number of neighbors.
  • J::Float64: The coupling strength.
  • lambda::Float64: The uniform local decay rates.
  • D::Float64: The diffusion coefficient.
  • x0_min::Float64: The minimum initial condition.
  • x0_max::Float64: The maximum initial condition.

Returns

  • OUModelEnsemble: The Ornstein-Uhlenbeck model.
source
GaussianExpansionCavityMethod.Phi4ModelRRGMethod
Phi4ModelRRG(N, K, J, lambdas, D, u)

Construct a Phi^4 model with a random regular graph coupling matrix with ferromagnetic interactions.

Arguments

  • N::Integer: The number of sites.
  • K::Union{Int, Float64}: The average number of neighbors.
  • J::Float64: The coupling strength.
  • lambdas::Vector{Float64}: The on-site linear terms.
  • D::Float64: The noise strength.
  • u::Float64: The cubic perturbative constant.

Returns

  • Phi4ModelEnsemble: The Phi^4 model.
source
GaussianExpansionCavityMethod.Phi4ModelRRGMethod
Phi4ModelRRG(N, K, J, lambda, D, u)

Construct a Phi^4 model with a random regular graph coupling matrix with ferromagnetic interactions.

Arguments

  • N::Integer: The number of sites.
  • K::Union{Int, Float64}: The average number of neighbors.
  • J::Float64: The coupling strength.
  • lambda::Float64: The on-site linear term.
  • D::Float64: The noise strength.
  • u::Float64: The cubic perturbative constant.

Returns

  • Phi4ModelEnsemble: The Phi^4 model.
source
GaussianExpansionCavityMethod.TwoSpinModelRRGMethod
TwoSpinModelRRG(N, K, J, lambda, D)

Construct a spherical 2-Spin model with a random regular graph coupling matrix with iid symmetric couplings Jᵢⱼ = Jⱼᵢ ~ p(J).

Arguments

  • N::Integer: The number of sites.
  • K::Union{Int, Float64}: The average number of neighbors.
  • pJ::Function: The coupling distribution.
  • J_params::Vector{Float64}: The parameters for the coupling distribution.
  • D::Float64: The noise strength.

Returns

  • TwoSpinModelEnsemble: The spherical 2-Spin model.
source
GaussianExpansionCavityMethod.TwoSpinModelRRG_BimMethod
TwoSpinModelRRG_Bim(N, K, J, lambda, D)

Construct a spherical 2-Spin model with a random regular graph coupling matrix with bimodal interactions.

Arguments

  • N::Integer: The number of sites.
  • K::Union{Int, Float64}: The average number of neighbors.
  • J::Float64: The coupling strength.
  • D::Float64: The noise strength.

Returns

  • TwoSpinModelEnsemble: The spherical 2-Spin model.
source
GaussianExpansionCavityMethod.TwoSpinModelRRG_FerroMethod
TwoSpinModelRRG_Ferro(N, K, J, lambda, D)

Construct a spherical 2-Spin model with a random regular graph coupling matrix with ferromagnetic interactions.

Arguments

  • N::Integer: The number of sites.
  • K::Union{Int, Float64}: The average number of neighbors.
  • J::Float64: The coupling strength.
  • D::Float64: The noise strength.

Returns

  • TwoSpinModelEnsemble: The spherical 2-Spin model.
source
GaussianExpansionCavityMethod.compute_autocorrMethod
compute_autocorr(trajs::Matrix{Float64}, tws_idxs::Vector{Int}; time_indices=nothing, showprogress=false)

Compute the average autocorrelation of the trajectories. It computes the autocorrelation for each waiting time specified by the vector tws_idxs. For each waiting time tw, it computes the autocorrelation C(t + tw, tw) for all t in the time_indices. If time_indices is nothing, it computes the autocorrelation for all time indices after tw.

Arguments

  • trajs::Matrix{Float64}: The trajectories. Each column corresponds to a time point.
  • tws_idxs::Vector{Int}: The waiting times to compute the autocorrelation.

Keyword Arguments

  • time_indices::Union{Nothing, AbstractVector{Int}}: The time indices to compute the autocorrelation. If nothing, all time indices are used.
  • showprogress::Bool: Whether to show a progress bar. Default is false.

Returns

  • autocorrs::Vector{Matrix{Float64}}: The average autocorrelation for each waiting time. The element at index (i, l, k) is the average autocorrelation of the trajectories at discretized times l and k for the waiting time tws_idxs[i].
  • time_indices::Vector{Vector{Int}}: The time indices used to compute the autocorrelation for each waiting time.
source
GaussianExpansionCavityMethod.compute_autocorrMethod
compute_autocorr(trajs; time_indices=nothing, showprogress=false)

Compute the average autocorrelation of the trajectories.

Arguments

  • trajs::Matrix{Float64}: The trajectories. Each column corresponds to a time point.

Keyword Arguments

  • time_indices::Union{Nothing, AbstractVector{Int}}: The time indices to compute the autocorrelation. If nothing, all time indices are used.
  • showprogress::Bool: Whether to show a progress bar. Default is false.

Returns

  • autocorr::Matrix{Float64}: The average autocorrelation. The element at index (l, k) is the average autocorrelation of the trajectories at discretized times l and k.
source
GaussianExpansionCavityMethod.compute_autocorrMethod
compute_autocorr(sim::Vector{Matrix{Float64}}, tws_idxs::Vector{Int}; time_indices=nothing, showprogress=false)

Compute the average autocorrelation of the trajectories in the ensemble solution object. It computes the autocorrelation for each waiting time specified by the vector tws_idxs. For each waiting time tw, it computes the autocorrelation C(t + tw, tw) for all t in the time_indices. If time_indices is nothing, it computes the autocorrelation for all time indices after tw.

Arguments

  • sim::Vector{Matrix{Float64}}: The ensemble solution object. Each element is a matrix where each column corresponds to a time point.
  • tws_idxs::Vector{Int}: The waiting times to compute the autocorrelation.

Keyword Arguments

  • time_indices::Union{Nothing, AbstractVector{Int}}: The time indices to compute the autocorrelation. If nothing, all time indices are used.
  • showprogress::Bool: Whether to show a progress bar. Default is false.

Returns

  • autocorrs::Vector{Matrix{Float64}}: The average autocorrelation for each waiting time. The element at index (i, l, k) is the average autocorrelation of the trajectories over nodes and ensemble realizations at discretized times l and k for the waiting time tws_idxs[i].
  • time_indices::Vector{Vector{Int}}: The time indices used to compute the autocorrelation for each waiting time.
source
GaussianExpansionCavityMethod.compute_autocorrMethod
compute_autocorr(sim; time_indices=nothing, showprogress=false)

Compute the average autocorrelation of the trajectories in the ensemble solution object.

Arguments

  • sim::Vector{Matrix{Float64}}: The ensemble solution object. Each element is a matrix where each column corresponds to a time point.

Keyword Arguments

  • time_indices::Union{Nothing, AbstractVector{Int}}: The time indices to compute the autocorrelation. If nothing, all time indices are used.
  • showprogress::Bool: Whether to show a progress bar. Default is false.

Returns

  • autocorr::Matrix{Float64}: The average autocorrelation. The element at index (l, k) is the average autocorrelation of the trajectories over nodes and ensemble realizations at discretized times l and k.
  • t_idx::Vector{Int}: The time indices used to compute the autocorrelation.
source
GaussianExpansionCavityMethod.compute_autocorr_TTIMethod
compute_autocorr_TTI(trajs, teq; lag_indices=nothing, showprogress=false)

Compute the average autocorrelation of the trajectories in the stationary phase, i.e. after the transient time teq. It assumes that after the transient time, the trajectories are stationary, therefore the autocorrelation is time tranlational invariant (TTI), i.e. it only depends on the time differences C(t,t') = C(t-t').

Arguments

  • trajs::Matrix{Float64}: The trajectories. Each column corresponds to a time point.
  • teq::Int: The transient time.

Keyword Arguments

  • lag_indices::Union{Nothing, AbstractVector{Int}}: The lag indices to compute the autocorrelation. If nothing, all lag indices are used.
  • showprogress::Bool: Whether to show a progress bar. Default is false.

Returns

  • autocorr::Vector{Float64}: The average TTI autocorrelation. The element at index l is the average autocorrelation of the trajectories at discretized time difference l.
  • err_autocorr::Vector{Float64}: The error associated to the TTI autocorrelation. The element at index l is error of the average autocorrelation of the trajectories at discretized time difference l.
  • l_idx::Vector{Int}: The lag indices used to compute the autocorrelation.
source
GaussianExpansionCavityMethod.compute_autocorr_TTIMethod
compute_autocorr_TTI(sim, teq; lag_indices=nothing, showprogress=false)

Compute the average autocorrelation of the trajectories in the ensemble solution object in the stationary phase, i.e. after the transient time teq. It assumes that after the transient time, the trajectories are stationary, therefore the autocorrelation is time tranlational invariant (TTI), i.e. it only depends on the time differences C(t,t') = C(t-t').

Arguments

  • sim::Vector{Matrix{Float64}}: The ensemble solution object. Each element is a matrix where each column corresponds to a time point.
  • teq::Int: The transient time.

Keyword Arguments

  • lag_indices::Union{Nothing, AbstractVector{Int}}: The lag indices to compute the autocorrelation. If nothing, all lag indices are used.
  • showprogress::Bool: Whether to show a progress bar. Default is false.

Returns

  • autocorr::Vector{Float64}: The average TTI autocorrelation. The element at index l is the average autocorrelation of the trajectories over nodes and ensemble realizations at discretized time difference l.
  • err_autocorr::Vector{Float64}: The error associated to the TTI autocorrelation. The element at index l is error of the average autocorrelation of the trajectories over nodes and ensemble realizations at discretized time difference l.
  • l_idx::Vector{Int}: The lag indices used to compute the autocorrelation.
source
GaussianExpansionCavityMethod.compute_averagesMethod
compute_averages(nodes, model, T)

Compute the equilibrium averages over nodes of the corrlation and response functions.

Arguments

  • nodes::Vector{NodeEQ}: The nodes containing the cavity autocorrelations and the marginal fields.
  • model::OUModel: The Ornstein-Uhlenbeck model.
  • T::Int: The number of time steps.

Returns

  • C_avg::OffsetVector{Float64, Vector{Float64}}: The average correlation function over the nodes.
  • R_avg::OffsetVector{Float64, Vector{Float64}}: The average response function over the nodes.
source
GaussianExpansionCavityMethod.compute_averagesMethod
compute_averages(nodes, model, T)

Compute the averages over nodes of the corrlation and response functions.

Arguments

  • nodes::Vector{Node}: The nodes containing the cavity autocorrelations and the marginal fields.
  • model::OUModel: The Ornstein-Uhlenbeck model.
  • T::Int: The number of time steps.

Returns

  • C_avg::OffsetMatrix{Float64, Matrix{Float64}}: The average correlation function over the nodes.
  • R_avg::OffsetMatrix{Float64, Matrix{Float64}}: The average response function over the nodes.
source
GaussianExpansionCavityMethod.compute_meanMethod
compute_mean(nodes, model, T)

Compute the mean over nodes of the average degree of freedom.

Arguments

  • nodes::Vector{Node}: The nodes containing the cavity autocorrelations and the marginal fields.
  • model::OUModel: The Ornstein-Uhlenbeck model.
  • T::Int: The number of time steps.

Returns

  • mu::Float64: The mean of the average degree of freedom over the nodes.
source
GaussianExpansionCavityMethod.compute_meanstdMethod
compute_meanstd(trajs; time_indices=nothing)

Compute the mean and standard deviation of the trajectories.

arguments

  • trajs::Matrix{Float64}: The trajectories. Each column corresponds to a time point.

Keyword Arguments

  • time_indices::Union{Nothing, AbstractVector{Int}}: The time indices to compute the mean and standard deviation. If nothing, all time indices are used.

Returns

  • mean_traj::Vector{Float64}: The mean trajectory. The element at index l is the mean of the trajectories at discretized time l`.
  • std_traj::Vector{Float64}: The standard deviation trajectory. The element at index l is the standard deviation of the trajectories at discretized time l.
  • t_idx::Vector{Int}: The time indices used to compute the mean and standard deviation.
source
GaussianExpansionCavityMethod.compute_meanstdMethod
compute_meanstd(sim; time_indices=nothing)

Compute the mean and standard deviation of the trajectories in the ensemble solution object.

Arguments

  • sim::Vector{Matrix{Float64}}: The ensemble solution object. Each element is a matrix where each column corresponds to a time point.

Keyword Arguments

  • time_indices::Union{Nothing, AbstractVector{Int}}: The time indices to compute the mean and standard deviation. If nothing, all time indices are used.

Returns

  • mean_traj::Vector{Float64}: The mean trajectory. The element at index l is the mean of the trajectories over nodes and ensemble realizations at discretized time l.
  • std_traj::Vector{Float64}: The standard deviation trajectory. The element at index l is the standard deviation of the trajectories over nodes and ensemble realizations at discretized time l.
  • t_idx::Vector{Int}: The time indices used to compute the mean and standard deviation.
source
GaussianExpansionCavityMethod.compute_statsMethod
compute_stats(trajs; time_indices=nothing, showprogress=false)

Compute the mean, standard deviation and average autocorrelation of the trajectories.

Arguments

  • trajs::Matrix{Float64}: The trajectories. Each column corresponds to a time point.

Keyword Arguments

  • time_indices::Union{Nothing, AbstractVector{Int}}: The time indices to compute the mean and standard deviation. If nothing, all time indices are used.
  • showprogress::Bool: Whether to show a progress bar. Default is false.

Returns

  • mean_traj::Vector{Float64}: The mean trajectory. The element at index l is the mean of the trajectories at discretized time l.
  • std_traj::Vector{Float64}: The standard deviation trajectory. The element at index l is the standard deviation of the trajectories at discretized time l.
  • autocorr::Matrix{Float64}: The average autocorrelation. The element at index (l, k) is the average autocorrelation of the trajectories at discretized times l and k.
  • t_idx::Vector{Int}: The time indices used to compute the mean and standard deviation.
source
GaussianExpansionCavityMethod.compute_statsMethod
compute_stats(sim; time_indices=nothing, showprogress=false)

Compute the mean, standard deviation and average autocorrelation of the trajectories in the ensemble solution object.

Arguments

  • sim::Vector{Matrix{Float64}}: The ensemble solution object. Each element is a matrix where each column corresponds to a time point.

Keyword Arguments

  • time_indices::Union{Nothing, AbstractVector{Int}}: The time indices to compute the mean and standard deviation. If nothing, all time indices are used.
  • showprogress::Bool: Whether to show a progress bar. Default is false.

Returns

  • mean_traj::Vector{Float64}: The mean trajectory. The element at index l is the mean of the trajectories over nodes and ensemble realizations at discretized time l.
  • std_traj::Vector{Float64}: The standard deviation trajectory. The element at index l is the standard deviation of the trajectories over nodes and ensemble realizations at discretized time l.
  • autocorr::Matrix{Float64}: The average autocorrelation. The element at index (l, k) is the average autocorrelation of the trajectories over nodes and ensemble realizations at discretized times l and k.
  • t_idx::Vector{Int}: The time indices used to compute the mean and standard deviation.
source
GaussianExpansionCavityMethod.compute_stats_TTIMethod
compute_stats_TTI(trajs, teq; time_indices=nothing, lag_indices=nothing, showprogress=false)

Compute the mean, standard deviation and average autocorrelation of the trajectories in the stationary phase, i.e. after the transient time teq. It assumes that after the transient time, the trajectories are stationary, therefore the autocorrelation is time tranlational invariant (TTI), i.e. it only depends on the time differences C(t,t') = C(t-t').

Arguments

  • trajs::Matrix{Float64}: The trajectories. Each column corresponds to a time point.
  • teq::Int: The transient time.

Keyword Arguments

  • time_indices::Union{Nothing, AbstractVector{Int}}: The time indices to compute the mean and standard deviation. If nothing, all time indices are used.
  • lag_indices::Union{Nothing, AbstractVector{Int}}: The lag indices to compute the autocorrelation. If nothing, all lag indices are used.
  • showprogress::Bool: Whether to show a progress bar. Default is false.

Returns

  • mean_traj::Vector{Float64}: The mean TTI trajectory. The element at index l is the mean of the trajectories at discretized time l.
  • std_traj::Vector{Float64}: The standard deviation of the TTI trajectory. The element at index l is the standard deviation of the trajectories at discretized time l.
  • t_idx::Vector{Int}: The time indices used to compute the mean and standard deviation.
  • autocorr::Vector{Float64}: The average TTI autocorrelation. The element at index l is the average autocorrelation of the trajectories at discretized time difference l.
  • err_autocorr::Vector{Float64}: The error associated to the TTI autocorrelation. The element at index l is error of the average autocorrelation of the trajectories at discretized time difference l.
  • l_idx::Vector{Int}: The lag indices used to compute the autocorrelation.
source
GaussianExpansionCavityMethod.compute_stats_TTIMethod
compute_stats_TTI(sim, teq; time_indices=nothing, lag_indices=nothing, showprogress=false)

Compute the mean, standard deviation and average autocorrelation of the trajectories in the ensemble solution object in the stationary phase, i.e. after the transient time teq. It assumes that after the transient time, the trajectories are stationary, therefore the autocorrelation is time tranlational invariant (TTI), i.e. it only depends on the time differences C(t,t') = C(t-t').

Arguments

  • sim::Vector{Matrix{Float64}}: The ensemble solution object. Each element is a matrix where each column corresponds to a time point.
  • teq::Int: The transient time.

Keyword Arguments

  • time_indices::Union{Nothing, AbstractVector{Int}}: The time indices to compute the mean and standard deviation. If nothing, all time indices are used.
  • lag_indices::Union{Nothing, AbstractVector{Int}}: The lag indices to compute the autocorrelation. If nothing, all lag indices are used.
  • showprogress::Bool: Whether to show a progress bar. Default is false.

Returns

  • mean_traj::Vector{Float64}: The mean TTI trajectory. The element at index l is the mean of the trajectories over nodes and ensemble realizations at discretized time l.
  • std_traj::Vector{Float64}: The standard deviation of the TTI trajectory. The element at index l is the standard deviation of the trajectories over nodes and ensemble realizations at discretized time l.
  • t_idx::Vector{Int}: The time indices used to compute the mean and standard deviation.
  • autocorr::Vector{Float64}: The average TTI autocorrelation. The element at index l is the average autocorrelation of the trajectories over nodes and ensemble realizations at discretized time difference l.
  • err_autocorr::Vector{Float64}: The error associated to the TTI autocorrelation. The element at index l is error of the average autocorrelation of the trajectories over nodes and ensemble realizations at discretized time difference l.
  • l_idx::Vector{Int}: The lag indices used to compute the autocorrelation.
source
GaussianExpansionCavityMethod.integrate_2spin_Bim_RRGMethod
integrate_2spin_Bim_RRG(K, J, D, dt, T; backup=false, backupfile="data/RRG/backup_matrices.jld2", backupevery=1000, showprogress=false)

Integrate the disorder averaged cavity equations for a 2-spin model on a random regular graph with bimodal interactions.

Srguments

  • K::Int: The average number of neighbors.
  • J::Float64: The coupling strength.
  • D::Float64: The noise strength.
  • dt::Float64: The time step for the integration.
  • T::Int: The number of iterations to run the integration.

Keyword Arguments

  • backup::Bool: Whether to save a backup of the matrices every backupevery iterations. Default is false.
  • backupfile::String: The filename for the backup file. Default is "data/RRG/backup_matrices.jld2".
  • backupevery::Int: The number of iterations to save a backup every. Default is 1000.
  • showprogress::Bool: Whether to show a progress bar. Default is false.

Returns

  • C::OffsetArray: The disorder averaged autocorrelation C matrix.
  • R::OffsetArray: The disorder averaged response R matrix.
  • Ch::OffsetArray: The disorder averaged cavity autocorrelation Ch matrix.
  • Rh::OffsetArray: The disorder averaged cavity response Rh matrix.
  • mu::OffsetArray: The Lagrange multiplier mu array.
source
GaussianExpansionCavityMethod.run_cavityMethod
run_cavity(model, dt, T; C0=1.0, mu0=0.0, showprogress=false)

Run the cavity method for the Ornstein-Uhlenbeck model.

Arguments

  • model::OUModel: The Ornstein-Uhlenbeck model.
  • dt::Float64: The time step size.
  • T::Int: The number of time steps.

Keyword Arguments

  • mu0::Float64: The initial value of the cavity marginal field (default is 0.0).
  • C0::Float64: The initial value of the cavity autocorrelation (default is 1.0).
  • showprogress::Bool: Whether to show progress bars (default is false).

Returns

  • nodes::Vector{Node}: The updated nodes after running the cavity method, containing the cavity autocorrelations and the marginal fields.
source
GaussianExpansionCavityMethod.run_cavity_EQMethod
run_cavity_EQ(model::OUModel, dt::Float64, T::Int)

Run the equilibrium cavity method for the equilibrium Ornstein-Uhlenbeck model.

Arguments

  • model::OUModel: The Ornstein-Uhlenbeck model.
  • dt::Float64: The time step size.
  • T::Int: The number of time steps.

Keyword Arguments

  • C0::Float64: The initial value of the cavity autocorrelation (default is 1.0).
  • showprogress::Bool: Whether to show progress bars (default is false).

Returns

  • nodes::Vector{NodeEQ}: The updated nodes after running the cavity method, containing the cavity autocorrelations and the marginal fields.
source
GaussianExpansionCavityMethod.sample_2SpinMethod
sample_2Spin(model, x0, tmax, tsave; rng=Xoshiro(1234), dt=1e-4)

Sample the spherical 2-Spin model using the Euler-Maruyama solver.

Arguments

  • model::TwoSpinModel: The spherical 2-Spin model to sample.
  • x0::Vector{Float64}: The initial condition.
  • tmax::Float64: The maximum time to integrate to.
  • tsave::Vector{Float64}: The time points to save the trajectory at.

Optional arguments

  • rng::AbstractRNG: The random number generator to use (default: Xoshiro(1234)).
  • dt::Float64: The time step size (default: 1e-4).

Returns

  • tvec::Vector{Float64}: The time points.
  • traj::Matrix{Float64}: The trajectories. Each column corresponds to a time point.
  • lambda_traj::Vector{Float64}: The Lagrange multipliers at each time point.
source
GaussianExpansionCavityMethod.sample_BMMethod
sample_BM(model, x0, tmax, tsave; rng=Xoshiro(1234), diverging_threshold=1e6, dt=1e-3)

Sample the Bouchaud-Mezard model using a split-step Euler-Maruyama solver. The local part is solved exactly, while the non-local part is solved using the Euler-Maruyama method.

Arguments

  • model::BMModel: The Bouchaud-Mezard model to sample.
  • x0::Vector{Float64}: The initial condition.
  • tmax::Float64: The maximum time to integrate to.
  • tsave::Vector{Float64}: The time points to save the trajectory at.

Optional arguments

  • rng::AbstractRNG: The random number generator to use (default: Xoshiro(1234)).
  • dt::Float64: The time step size (default: 1e-3).

Returns

  • t_vals::Vector{Float64}: The time points.
  • trajectories::Matrix{Float64}: The trajectories. Each column corresponds to a time point.
source
GaussianExpansionCavityMethod.sample_OUMethod
sample_OU(model, x0, tmax, tsave; rng=Xoshiro(1234), diverging_threshold=1e6, dt=1e-3)

Sample the Ornstein-Uhlenbeck model using the Euler-Maruyama solver.

Arguments

  • model::OUModel: The Ornstein-Uhlenbeck model to sample.
  • x0::Vector{Float64}: The initial condition.
  • tmax::Float64: The maximum time to integrate to.
  • tsave::Vector{Float64}: The time points to save the trajectory at.

Optional arguments

  • rng::AbstractRNG: The random number generator to use (default: Xoshiro(1234)).
  • diverging_threshold::Float64: The threshold for detecting diverging solutions (default: 1e6).
  • dt::Float64: The time step size (default: 1e-3).

Returns

  • t_vals::Vector{Float64}: The time points.
  • trajectories::Matrix{Float64}: The trajectories. Each column corresponds to a time point.
  • sol::RODESolution: The solution object.
source
GaussianExpansionCavityMethod.sample_ensemble_2SpinMethod
sample_ensemble_2Spin(model_ensemble, x0_min, x0_max, tmax, tsave, nsample; rng=Xoshiro(1234), showprogress=false, dt=1e-4)

Sample an ensemble of spherical 2-Spin models using solvers from DifferentialEquations.jl.

Arguments

  • model_ensemble::TwoSpinModelEnsemble: The ensemble of Phi^4 models to sample.
  • x0_min::Float64: The minimum initial condition.
  • x0_max::Float64: The maximum initial condition.
  • tmax::Float64: The maximum time to integrate to.
  • tsave::Vector{Float64}: The time points to save the trajectory at.
  • nsample::Int: The number of samples to generate.

Optional arguments

  • rng::AbstractRNG: The random number generator to use (default: Xoshiro(1234)).
  • showprogress::Bool: Whether to show a progress bar (default: false).
  • dt::Float64: The time step size (default: 1e-4).

Returns

  • tvals_alls::Vector{Vector{Float64}}: The time points for each sample.
  • traj_alls::Vector{Matrix{Float64}}: The trajectories for each sample. Each column corresponds to a time point.
  • lambda_traj_alls::Vector{Vector{Float64}}: The Lagrange multipliers for each sample.
source
GaussianExpansionCavityMethod.sample_ensemble_BMMethod
sample_ensemble_BM(model_ensemble, x0_min, x0_max, tmax, tsave, nsample; rng=Xoshiro(1234), diverging_threshold=1e6, showprogress=false, dt=1e-3)

Sample an ensemble of Bouchaud-Mezard models using solvers from DifferentialEquations.jl.

Arguments

  • model_ensemble::BMModelEnsemble: The ensemble of Phi^4 models to sample.
  • x0_min::Float64: The minimum initial condition.
  • x0_max::Float64: The maximum initial condition.
  • tmax::Float64: The maximum time to integrate to.
  • tsave::Vector{Float64}: The time points to save the trajectory at.
  • nsample::Int: The number of samples to generate.

Optional arguments

  • rng::AbstractRNG: The random number generator to use (default: Xoshiro(1234)).
  • diverging_threshold::Float64: The threshold for detecting diverging solutions (default: 1e6).
  • showprogress::Bool: Whether to show a progress bar (default: false).
  • dt::Float64: The time step size (default: 1e-3).

Returns

  • tvals_alls::Vector{Vector{Float64}}: The time points for each sample.
  • traj_alls::Vector{Matrix{Float64}}: The trajectories for each sample. Each column corresponds to a time point.
source
GaussianExpansionCavityMethod.sample_ensemble_OUMethod
sample_ensemble_OU(model_ensemble, tmax, tsave, nsample; rng=Xoshiro(1234), diverging_threshold=1e6, showprogress=false, dt=1e-3)

Sample an ensemble of Ornstein-Uhlenbeck models using solvers from DifferentialEquations.jl.

Arguments

  • model_ensemble::OUModelEnsemble: The ensemble of Ornstein-Uhlenbeck models to sample.
  • tmax::Float64: The maximum time to integrate to.
  • tsave::Vector{Float64}: The time points to save the trajectory at.
  • nsample::Int: The number of samples to generate.

Optional arguments

  • rng::AbstractRNG: The random number generator to use (default: Xoshiro(1234)).
  • diverging_threshold::Float64: The threshold for detecting diverging solutions (default: 1e6).
  • showprogress::Bool: Whether to show a progress bar (default: false).
  • dt::Float64: The time step size (default: 1e-3).

Returns

  • tvals_alls::Vector{Vector{Float64}}: The time points for each sample.
  • traj_alls::Vector{Matrix{Float64}}: The trajectories for each sample. Each column corresponds to a time point.
  • sim::Vector{RODESolution}: The solution objects for each sample.
source
GaussianExpansionCavityMethod.sample_ensemble_phi4Method
sample_ensemble_phi4(model_ensemble, x0_min, x0_max, tmax, tsave, nsample; rng=Xoshiro(1234), diverging_threshold=1e6, showprogress=false, dt=1e-3)

Sample an ensemble of Phi^4 models using solvers from DifferentialEquations.jl.

Arguments

  • model_ensemble::Phi4ModelEnsemble: The ensemble of Phi^4 models to sample.
  • x0_min::Float64: The minimum initial condition.
  • x0_max::Float64: The maximum initial condition.
  • tmax::Float64: The maximum time to integrate to.
  • tsave::Vector{Float64}: The time points to save the trajectory at.
  • nsample::Int: The number of samples to generate.

Optional arguments

  • rng::AbstractRNG: The random number generator to use (default: Xoshiro(1234)).
  • diverging_threshold::Float64: The threshold for detecting diverging solutions (default: 1e6).
  • showprogress::Bool: Whether to show a progress bar (default: false).
  • dt::Float64: The time step size (default: 1e-3).

Returns

  • tvals_alls::Vector{Vector{Float64}}: The time points for each sample.
  • traj_alls::Vector{Matrix{Float64}}: The trajectories for each sample. Each column corresponds to a time point.
  • sim::Vector{RODESolution}: The solution objects for each sample.
source
GaussianExpansionCavityMethod.sample_phi4Method
sample_phi4(model, x0, tmax, tsave; rng=Xoshiro(1234), diverging_threshold=1e6, dt=1e-3)

Sample the Phi^4 model using the Euler-Maruyama solver.

Arguments

  • model::Phi4Model: The Phi^4 model to sample.
  • x0::Vector{Float64}: The initial condition.
  • tmax::Float64: The maximum time to integrate to.
  • tsave::Vector{Float64}: The time points to save the trajectory at.

Optional arguments

  • rng::AbstractRNG: The random number generator to use (default: Xoshiro(1234)).
  • diverging_threshold::Float64: The threshold for detecting diverging solutions (default: 1e6).
  • dt::Float64: The time step size (default: 1e-3).

Returns

  • t_vals::Vector{Float64}: The time points.
  • trajectories::Matrix{Float64}: The trajectories. Each column corresponds to a time point.
  • sol::RODESolution: The solution object.
source