API

This is the documentation of the Arianna module's functions, types and structures.

API Reference

Index

Types and Functions

Arianna.AriannaModule
module Arianna

Arianna is a flexible and extensible framework for Monte Carlo simulations. Instead of acting as a black-box simulator, it provides a modular structure where users define their own system and Monte Carlo "moves".

source
Arianna.ActionType
abstract type Action

Abstract type representing Monte Carlo actions/moves that can be performed on a system.

Concrete subtypes must implement:

  • sample_action!(action, policy, parameters, system, rng): Sample a new action from the policy
  • perform_action!(system, action): Apply the action to modify the system state
  • invert_action!(action, system): Invert/reverse the action
  • log_proposal_density(action, policy, parameters, system): Log probability density of proposing this action

Optional methods:

  • revert_action!(system, action): Optimized version of perform_action! that can reuse cached values
source
Arianna.AriannaSystemType
abstract type AriannaSystem

Abstract type representing a system that can be simulated using methods defined in the Arianna module.

source
Arianna.MetropolisType
Metropolis{P,R<:AbstractRNG,C<:Function} <: AriannaAlgorithm

A struct representing a Metropolis Monte Carlo algorithm.

Fields

  • pools::Vector{P}: Vector of independent pools (one for each system)
  • sweepstep::Int: Number of Monte Carlo steps per sweep
  • seed::Int: Random number seed
  • rngs::Vector{R}: Vector of random number generators
  • parallel::Bool: Flag to parallelise over different threads
  • collecter::C: Transducer to collect results from parallelised loops

Type Parameters

  • P: Type of the pool
  • R: Type of the random number generator
  • C: Type of the transducer
source
Arianna.MetropolisMethod
Metropolis(chains; pool=missing, sweepstep=1, seed=1, R=Xoshiro, parallel=false, extras...)

Create a new Metropolis instance.

Arguments

  • chains: Vector of systems to run the algorithm on
  • pool: Pool of moves to perform sweeps over
  • sweepstep: Number of Monte Carlo steps per sweep
  • seed: Random number seed
  • R: Type of the random number generator
  • parallel: Flag to parallelise over different threads
  • extras: Additional keyword arguments

Returns

  • algorithm: Metropolis instance
source
Arianna.MoveType
Move{A<:Action,P<:Policy,V<:AbstractArray,T<:AbstractFloat}

A struct representing a Monte Carlo move with an associated action, policy, and parameters.

Fields

  • action::A: The action to be performed in the move
  • policy::P: The policy used to propose actions
  • parameters::V: Parameters of the policy
  • weight::T: Weight/probability of selecting this move in a sweep
  • total_calls::Int: Total number of times this move has been attempted
  • accepted_calls::Int: Number of times this move has been accepted

Type Parameters

  • A: Type of the action
  • P: Type of the policy
  • V: Type of the parameter array
  • T: Type of the weight (floating point)
source
Arianna.MoveMethod
Move(action, policy, parameters, weight)

Create a new Move instance.

Arguments

  • action: The action to be performed in the move
  • policy: The policy used to propose actions
  • parameters: Parameters of the policy
  • weight: Weight/probability of selecting this move in a sweep
source
Arianna.PolicyType
abstract type Policy

Abstract type representing proposal policies for Monte Carlo actions.

A Policy defines how actions are sampled and their proposal probabilities are calculated. Concrete subtypes work together with specific Action types to implement the proposal mechanism.

source
Arianna.PrintTimeStepsType
PrintTimeSteps <: AriannaAlgorithm

Algorithm to display a progress bar and current timestep during simulation.

source
Arianna.SimulationType
mutable struct Simulation{S, A, VS}

A structure representing a Monte Carlo simulation.

Fields

  • chains::Vector{S}: Vector of independent Arianna systems.
  • algorithms::A: List of algorithms.
  • steps::Int: Number of MC sweeps.
  • t::Int: Current time step.
  • schedulers::VS: List of schedulers (one for each algorithm).
  • counters::Vector{Int}: Counters for the schedulers (one for each algorithm).
  • path::String: Simulation path.
  • verbose::Bool: Flag for verbose output.
source
Arianna.SimulationMethod
Simulation(chains, algorithm_list, steps; path="data", verbose=false)

Create a new Simulation instance from a list of algorithm constructors.

Arguments

  • chains: Vector of independent Arianna systems.
  • algorithm_list: List of algorithm constructors.
  • steps: Number of MC sweeps.
  • path="data": Simulation path.
  • verbose=false: Flag for verbose output.
source
Arianna.StoreBackupsType
StoreBackups <: AriannaAlgorithm

Algorithm to create backup files of system states during simulation.

Fields

  • dirs::Vector{String}: Directories for storing backup files
  • fmt::Format: Format type for output files
  • store_first::Bool: Whether to store backups at initialization
  • store_last::Bool: Whether to store backups at finalization
source
Arianna.StoreCallbacksType
StoreCallbacks{V} <: AriannaAlgorithm

Algorithm to store callback values during simulation.

Fields

  • callbacks::V: Vector of callback functions to evaluate
  • paths::Vector{String}: Paths to output files for each callback
  • files::Vector{IOStream}: File handles for writing callback values
  • store_first::Bool: Whether to store callback values at initialization
  • store_last::Bool: Whether to store callback values at finalization

Constructor

StoreCallbacks(callbacks::V, path; store_first::Bool=true, store_last::Bool=false)

Create a new StoreCallbacks instance.

Arguments

  • callbacks::V: Vector of callback functions
  • path: Base path for output files
  • store_first=true: Store values at initialization
  • store_last=false: Store values at finalization
source
Arianna.StoreLastFramesType
StoreLastFrames <: AriannaAlgorithm

Algorithm to store the final state of each system at the end of simulation.

Fields

  • paths::Vector{String}: Paths to output files
  • fmt::Format: Format type for output files
source
Arianna.StoreParametersType
StoreParameters{V<:AbstractArray} <: AriannaAlgorithm

A struct representing a parameter store.

Fields

  • paths::Vector{String}: Vector of paths to the parameter files
  • files::Vector{IOStream}: Vector of open file streams to the parameter files
  • parameters_list::V: List of parameters to store
  • store_first::Bool: Flag to store the parameters at the first step
  • store_last::Bool: Flag to store the parameters at the last step

Type Parameters

  • V: Type of the parameter array
source
Arianna.StoreParametersMethod
StoreParameters(chains; dependencies=missing, path=missing, ids=missing, store_first=true, store_last=false, extras...)

Create a new StoreParameters instance.

Arguments

  • chains: Vector of chains to store the parameters of
  • dependencies: Dependencies of the parameter store
  • path: Path to the parameter files
  • ids: IDs of the parameters to store
  • store_first: Flag to store the parameters at the first step
  • store_last: Flag to store the parameters at the last step
  • extras: Additional keyword arguments

Returns

  • algorithm: StoreParameters instance
source
Arianna.StoreTrajectoriesType
StoreTrajectories{F<:Format} <: AriannaAlgorithm

Algorithm to store system trajectories during simulation.

Fields

  • paths::Vector{String}: Paths to output files
  • files::Vector{IOStream}: File handles for writing trajectories
  • fmt::F: Format type for output files
  • store_first::Bool: Whether to store trajectories at initialization
  • store_last::Bool: Whether to store trajectories at finalization
source
Arianna.build_scheduleMethod
build_schedule(steps::Int, burn::Int, base::AbstractFloat)

Create a vector of timestep from burn to steps log-spaced with base base.

source
Arianna.build_scheduleMethod
build_schedule(steps::Int, burn::Int, Δt::Int)

Create a vector of timestep from burn to steps at intervals Δt.

source
Arianna.build_scheduleMethod
build_schedule(steps::Int, burn::Int, block::Vector{Int})

Create a vector of timestep from burn to steps with repeated blocks specified by block.

source
Arianna.callback_acceptanceMethod
callback_acceptance(simulation)

Calculate the mean acceptance rate of the Metropolis algorithm.

Arguments

  • simulation: Simulation to calculate the acceptance rate of
source
Arianna.delta_log_target_densityMethod
delta_log_target_density(x₁, x₂, system::AriannaSystem)

Calculate the change in log target density between two states.

Arguments

  • x₁: Initial state
  • x₂: Final state
  • system: System object
source
Arianna.finaliseMethod
finalise(::AriannaAlgorithm, ::Simulation)

Finalise the algorithm for the given simulation.

source
Arianna.finaliseMethod
store_lastframe(io, system, t, fmt::Format)

Store the final state of the system at time t to the given IO stream in the specified format.

source
Arianna.initialiseMethod
initialise(::AriannaAlgorithm, ::Simulation)

Initialise the algorithm for the given simulation.

source
Arianna.initialiseMethod
store_backup(io, system, t, fmt::Format)

Store a backup of the system state at time t to the given IO stream in the specified format.

source
Arianna.invert_action!Method
invert_action!(action::Action, system::AriannaSystem)

Invert/reverse an action.

Arguments

  • action: Action to invert
  • system: System the action was applied to
source
Arianna.log_proposal_densityMethod
log_proposal_density(action, policy, parameters, system::AriannaSystem)

Calculate the log probability density of proposing the given action.

Arguments

  • action: Proposed action
  • policy: Policy used for proposal
  • parameters: Parameters of the policy
  • system: System the action is applied to
source
Arianna.make_step!Method
make_step!(::Simulation, ::AriannaAlgorithm)

Perform a single step of the algorithm in the simulation.

source
Arianna.make_step!Method
make_step!(simulation::Simulation, algorithm::Metropolis)

Perform a single step of the Metropolis algorithm.

Arguments

  • simulation: Simulation to perform the step on
  • algorithm: Metropolis instance
source
Arianna.mc_step!Method
mc_step!(system::AriannaSystem, action::Action, policy::Policy, parameters::AbstractArray{T}, rng) where {T<:AbstractFloat}

Perform a single Monte Carlo step.

Arguments

  • system: System to modify
  • action: Action to perform
  • policy: Policy used for proposal
  • parameters: Parameters of the policy
  • rng: Random number generator
source
Arianna.mc_sweep!Method
mc_sweep!(system::AriannaSystem, pool, rng; mc_steps=1)

Perform a Monte Carlo sweep over a pool of moves.

Arguments

  • system: System to modify
  • pool: Pool of moves to perform sweeps over
  • rng: Random number generator
  • mc_steps: Number of Monte Carlo steps per sweep
source
Arianna.perform_action!Method
perform_action!(system::AriannaSystem, action::Action)

Apply the action to modify the system state.

Arguments

  • system: System to modify
  • action: Action to perform

Returns

A tuple of (x₁, x₂) containing the old and new states

source
Arianna.raise_errorMethod
raise_error(s)

Helper function to raise errors for unimplemented required methods.

Arguments

  • s: String describing the missing method implementation
source
Arianna.revert_action!Method
revert_action!(system, action::Action)

Optimized version of perform_action! that can reuse cached values.

Arguments

  • system: System to modify
  • action: Action to perform
source
Arianna.run!Method
run!(simulation::Simulation)

Run the Monte Carlo simulation.

Arguments

  • simulation::Simulation: The simulation instance to run.
source
Arianna.sample_action!Method
sample_action!(action::Action, policy::Policy, parameters, system::AriannaSystem, rng)

Sample a new action from the policy.

Arguments

  • action: Action to be sampled
  • policy: Policy to sample from
  • parameters: Parameters of the policy
  • system: System the action will be applied to
  • rng: Random number generator
source
Arianna.store_trajectoryMethod
store_trajectory(io, system::AriannaSystem, t, fmt::Format)

Store the system trajectory at time t to the given IO stream in the specified format.

source
Arianna.write_algorithmMethod
write_algorithm(io, algorithm::AriannaAlgorithm, scheduler)

Write a summary of the algorithm on the given IO stream.

source
Arianna.write_algorithmMethod
write_algorithm(io, algorithm::Metropolis, scheduler)

Write the algorithm to a string.

Arguments

  • io: IO stream to write to
  • algorithm: Algorithm to write
  • scheduler: Scheduler to write
source
Arianna.write_parametersMethod
write_parameters(::Policy, parameters)

Write the parameters of a policy to a string.

Arguments

  • policy: Policy to write the parameters of
  • parameters: Parameters of the policy
source