API
This is the documentation of the Arianna module's functions, types and structures.
API Reference
Index
Arianna.AriannaArianna.LinterArianna.ActionArianna.AriannaAlgorithmArianna.AriannaSystemArianna.DATArianna.FormatArianna.Linter.NodeArianna.MetropolisArianna.MetropolisArianna.MoveArianna.MoveArianna.PolicyArianna.PrintTimeStepsArianna.SimulationArianna.SimulationArianna.StoreBackupsArianna.StoreCallbacksArianna.StoreLastFramesArianna.StoreParametersArianna.StoreParametersArianna.StoreTrajectoriesArianna.TXTAbstractTrees.printnodeArianna.Linter.add_action_node!Arianna.Linter.add_policy_node!Arianna.Linter.add_system_node!Arianna.Linter.collect_defined!Arianna.Linter.collect_defined_functions!Arianna.Linter.collect_defined_structs!Arianna.Linter.collect_defined_types!Arianna.Linter.construct_treeArianna.Linter.extract_defined_functionsArianna.Linter.extract_defined_structsArianna.Linter.extract_defined_typesArianna.Linter.run_linterArianna.build_scheduleArianna.build_scheduleArianna.build_scheduleArianna.callback_acceptanceArianna.delta_log_target_densityArianna.finaliseArianna.finaliseArianna.initialiseArianna.initialiseArianna.invert_action!Arianna.log_proposal_densityArianna.make_step!Arianna.make_step!Arianna.mc_step!Arianna.mc_sweep!Arianna.perform_action!Arianna.raise_missingfunctionArianna.revert_action!Arianna.run!Arianna.sample_action!Arianna.store_trajectoryArianna.unnormalised_log_target_densityArianna.write_algorithmArianna.write_algorithmArianna.write_parameters
Types and Functions
Arianna.Arianna — Modulemodule AriannaArianna 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".
Arianna.Action — Typeabstract type ActionAbstract 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 policyperform_action!(system, action): Apply the action to modify the system stateinvert_action!(action, system): Invert/reverse the actionlog_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
Arianna.AriannaAlgorithm — Typeabstract type AriannaAlgorithmAbstract type for Simulation algorithms.
Arianna.AriannaSystem — Typeabstract type AriannaSystemAbstract type representing a system that can be simulated using methods defined in the Arianna module.
Arianna.DAT — TypeDAT <: FormatFormat type for data (.dat) file output.
Arianna.Format — TypeFormatAbstract type for output file formats.
Arianna.Metropolis — TypeMetropolis{P,R<:AbstractRNG,C<:Function} <: AriannaAlgorithmA 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 sweepseed::Int: Random number seedrngs::Vector{R}: Vector of random number generatorsparallel::Bool: Flag to parallelise over different threadscollecter::C: Transducer to collect results from parallelised loops
Type Parameters
P: Type of the poolR: Type of the random number generatorC: Type of the transducer
Arianna.Metropolis — MethodMetropolis(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 onpool: Pool of moves to perform sweeps oversweepstep: Number of Monte Carlo steps per sweepseed: Random number seedR: Type of the random number generatorparallel: Flag to parallelise over different threadsextras: Additional keyword arguments
Returns
algorithm: Metropolis instance
Arianna.Move — TypeMove{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 movepolicy::P: The policy used to propose actionsparameters::V: Parameters of the policyweight::T: Weight/probability of selecting this move in a sweeptotal_calls::Int: Total number of times this move has been attemptedaccepted_calls::Int: Number of times this move has been accepted
Type Parameters
A: Type of the actionP: Type of the policyV: Type of the parameter arrayT: Type of the weight (floating point)
Arianna.Move — MethodMove(action, policy, parameters, weight)Create a new Move instance.
Arguments
action: The action to be performed in the movepolicy: The policy used to propose actionsparameters: Parameters of the policyweight: Weight/probability of selecting this move in a sweep
Arianna.Policy — Typeabstract type PolicyAbstract 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.
Arianna.PrintTimeSteps — TypePrintTimeSteps <: AriannaAlgorithmAlgorithm to display a progress bar and current timestep during simulation.
Arianna.Simulation — Typemutable 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.
Arianna.Simulation — MethodSimulation(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.
Arianna.StoreBackups — TypeStoreBackups <: AriannaAlgorithmAlgorithm to create backup files of system states during simulation.
Fields
dirs::Vector{String}: Directories for storing backup filesfmt::Format: Format type for output filesstore_first::Bool: Whether to store backups at initializationstore_last::Bool: Whether to store backups at finalization
Arianna.StoreCallbacks — TypeStoreCallbacks{V} <: AriannaAlgorithmAlgorithm to store callback values during simulation.
Fields
callbacks::V: Vector of callback functions to evaluatepaths::Vector{String}: Paths to output files for each callbackfiles::Vector{IOStream}: File handles for writing callback valuesstore_first::Bool: Whether to store callback values at initializationstore_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 functionspath: Base path for output filesstore_first=true: Store values at initializationstore_last=false: Store values at finalization
Arianna.StoreLastFrames — TypeStoreLastFrames <: AriannaAlgorithmAlgorithm to store the final state of each system at the end of simulation.
Fields
paths::Vector{String}: Paths to output filesfmt::Format: Format type for output files
Arianna.StoreParameters — TypeStoreParameters{V<:AbstractArray} <: AriannaAlgorithmA struct representing a parameter store.
Fields
paths::Vector{String}: Vector of paths to the parameter filesfiles::Vector{IOStream}: Vector of open file streams to the parameter filesparameters_list::V: List of parameters to storestore_first::Bool: Flag to store the parameters at the first stepstore_last::Bool: Flag to store the parameters at the last step
Type Parameters
V: Type of the parameter array
Arianna.StoreParameters — MethodStoreParameters(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 ofdependencies: Dependencies of the parameter storepath: Path to the parameter filesids: IDs of the parameters to storestore_first: Flag to store the parameters at the first stepstore_last: Flag to store the parameters at the last stepextras: Additional keyword arguments
Returns
algorithm: StoreParameters instance
Arianna.StoreTrajectories — TypeStoreTrajectories{F<:Format} <: AriannaAlgorithmAlgorithm to store system trajectories during simulation.
Fields
paths::Vector{String}: Paths to output filesfiles::Vector{IOStream}: File handles for writing trajectoriesfmt::F: Format type for output filesstore_first::Bool: Whether to store trajectories at initializationstore_last::Bool: Whether to store trajectories at finalization
Arianna.TXT — TypeTXT <: FormatFormat type for text (.txt) file output.
Arianna.build_schedule — Methodbuild_schedule(steps::Int, burn::Int, base::AbstractFloat)Create a vector of timestep from burn to steps log-spaced with base base.
Arianna.build_schedule — Methodbuild_schedule(steps::Int, burn::Int, Δt::Int)Create a vector of timestep from burn to steps at intervals Δt.
Arianna.build_schedule — Methodbuild_schedule(steps::Int, burn::Int, block::Vector{Int})Create a vector of timestep from burn to steps with repeated blocks specified by block.
Arianna.callback_acceptance — Methodcallback_acceptance(simulation)Calculate the mean acceptance rate of the Metropolis algorithm.
Arguments
simulation: Simulation to calculate the acceptance rate of
Arianna.delta_log_target_density — Methoddelta_log_target_density(x₁, x₂, system::AriannaSystem)Calculate the change in log target density between two states.
Arguments
x₁: Initial statex₂: Final statesystem: System object
Arianna.finalise — Methodfinalise(::AriannaAlgorithm, ::Simulation)Finalise the algorithm for the given simulation.
Arianna.finalise — Methodstore_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.
Arianna.initialise — Methodinitialise(::AriannaAlgorithm, ::Simulation)Initialise the algorithm for the given simulation.
Arianna.initialise — Methodstore_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.
Arianna.invert_action! — Methodinvert_action!(action::Action, system::AriannaSystem)Invert/reverse an action.
Arguments
action: Action to invertsystem: System the action was applied to
Arianna.log_proposal_density — Methodlog_proposal_density(action, policy, parameters, system::AriannaSystem)Calculate the log probability density of proposing the given action.
Arguments
action: Proposed actionpolicy: Policy used for proposalparameters: Parameters of the policysystem: System the action is applied to
Arianna.make_step! — Methodmake_step!(::Simulation, ::AriannaAlgorithm)Perform a single step of the algorithm in the simulation.
Arianna.make_step! — Methodmake_step!(simulation::Simulation, algorithm::Metropolis)Perform a single step of the Metropolis algorithm.
Arguments
simulation: Simulation to perform the step onalgorithm: Metropolis instance
Arianna.mc_step! — Methodmc_step!(system::AriannaSystem, action::Action, policy::Policy, parameters::AbstractArray{T}, rng) where {T<:AbstractFloat}Perform a single Monte Carlo step.
Arguments
system: System to modifyaction: Action to performpolicy: Policy used for proposalparameters: Parameters of the policyrng: Random number generator
Arianna.mc_sweep! — Methodmc_sweep!(system::AriannaSystem, pool, rng; mc_steps=1)Perform a Monte Carlo sweep over a pool of moves.
Arguments
system: System to modifypool: Pool of moves to perform sweeps overrng: Random number generatormc_steps: Number of Monte Carlo steps per sweep
Arianna.perform_action! — Methodperform_action!(system::AriannaSystem, action::Action)Apply the action to modify the system state.
Arguments
system: System to modifyaction: Action to perform
Returns
A tuple of (x₁, x₂) containing the old and new states
Arianna.raise_missingfunction — Methodraise_missingfunction(s)Helper function to raise errors for unimplemented required methods.
Arguments
s: String describing the missing method implementation
Arianna.revert_action! — Methodrevert_action!(system, action::Action)Optimized version of perform_action! that can reuse cached values.
Arguments
system: System to modifyaction: Action to perform
Arianna.run! — Methodrun!(simulation::Simulation)Run the Monte Carlo simulation.
Arguments
simulation::Simulation: The simulation instance to run.
Arianna.sample_action! — Methodsample_action!(action::Action, policy::Policy, parameters, system::AriannaSystem, rng)Sample a new action from the policy.
Arguments
action: Action to be sampledpolicy: Policy to sample fromparameters: Parameters of the policysystem: System the action will be applied torng: Random number generator
Arianna.store_trajectory — Methodstore_trajectory(io, system::AriannaSystem, t, fmt::Format)Store the system trajectory at time t to the given IO stream in the specified format.
Arianna.unnormalised_log_target_density — Methodunnormalised_log_target_density(x, system)Calculate the unnormalized log probability density of a system state.
Arguments
x: System statesystem: System object
Arianna.write_algorithm — Methodwrite_algorithm(io, algorithm::AriannaAlgorithm, scheduler)Write a summary of the algorithm on the given IO stream.
Arianna.write_algorithm — Methodwrite_algorithm(io, algorithm::Metropolis, scheduler)Write the algorithm to a string.
Arguments
io: IO stream to write toalgorithm: Algorithm to writescheduler: Scheduler to write
Arianna.write_parameters — Methodwrite_parameters(::Policy, parameters)Write the parameters of a policy to a string.
Arguments
policy: Policy to write the parameters ofparameters: Parameters of the policy
Arianna.Linter — ModuleLinterModule for linting user's code in the Arianna framework. It checks for the presence of required structs, functions, and their definitions.
Arianna.Linter.Node — TypeNodeRepresents a node in the hierarchical Arianna system tree. Each node can represent a system, action, policy, or function, and may contain child nodes.
Fields
label::Symbol: The name or identifier of the node (e.g., function or struct name).ntype::Symbol: The type of the node (e.g.,:System,:Action,:Policy, or:Function).toptype::Symbol: The top-level category the node belongs to (typically used for system type inheritance).children::Vector{Node}: A list of childNodes representing the structure beneath this node.state::String: The linting state of the node (e.g.,"✅ Defined"or"❌ Missing").
Usage
Used to build and represent the structured tree for static analysis and linting of Arianna-based systems.
AbstractTrees.printnode — MethodAbstractTrees.printnode(io::IO, n::Node)Prints a formatted representation of a Node object to the given IO stream.
Depending on the toptype of the node, the output is styled and labeled accordingly:
- If
toptypeis:function, prints the function label and its state, padded for alignment. - If
toptypeis:AriannaSystem, prints the system label in green. - If
toptypeis:Action, prints the action label in red. - If
toptypeis:Policy, prints the policy label in blue.
Arguments
io::IO: The IO stream to print to.n::Node: The node to be printed.
Arianna.Linter.add_action_node! — Functionadd_action_node!(actions, aname, policies; atype=:Action)Adds a new action node to the actions collection.
Arguments
actions: The collection to which the new action node will be added.aname: The name of the action node.policies: A list of policy nodes to attach as children to the action node.atype: (Optional) The type of the node, defaults to:Action.
Behavior
Creates an action node with the specified name and type, attaches all required function children, and adds the provided policies as children nodes.
Returns
Modifies actions in-place by adding the new action node.
Arianna.Linter.add_policy_node! — Functionadd_policy_node!(policies, pname; ptype=:Policy)Add a new policy node to the policies collection.
Arguments
policies: The collection (e.g., dictionary or custom structure) to which the policy node will be added.pname: The name (identifier) of the policy node.ptype: The type of the policy node (default is:Policy).
Behavior
Initializes the policy node with its required function children set to "Not defined".
Returns
Modifies policies in-place by adding the new policy node.
Arianna.Linter.add_system_node! — Methodadd_system_node!(systems, sname, actions; tname=:AriannaSystem)Add a new system node to the systems collection.
Arguments
systems: The collection to which the system node will be added.sname: Symbol or string specifying the name of the system node.actions: A collection of actions to be added as children of the system node.tname: (Optional) Symbol specifying the type of the system node. Defaults to:AriannaSystem.
Behavior
Creates a system node with the specified name and type, attaches a required function child, and adds all provided actions as children nodes.
Returns
Modifies systems in-place by adding the new system node.
Arianna.Linter.collect_defined! — Methodcollect_defined!(results, node, symb::Symbol; new_types=nothing)Recursively traverses the syntax tree starting at node and collects definitions of the specified kind into results.
symb: The type of definitions to collect (:abstract,:struct, or:function).new_types: Optional set of additional struct types to recognize whensymb == :struct.
Arguments
results: A mutable collection to store the collected definitions.node: ASyntaxNodefrom CSTParser representing the current location in the syntax tree.symb::Symbol: The kind of definition to collect.new_types: (optional) Additional struct types beyond the default ones.
Returns
- Nothing. The function modifies
resultsin place.
Arianna.Linter.collect_defined_functions! — Methodcollect_defined_functions!(results, node::SyntaxNode)Helper function that collects function definitions from the given node into results.
Functions are only collected if they belong to Arianna or Arianna.PolicyGuided namespaces.
Each entry is stored as a tuple:
function_name: The name of the function as aSymbol.arity: The number of arguments.arg_types: ASet{Symbol}of the argument type names, if specified.
Arguments
results: A mutable set to store(function_name, arity, arg_types)tuples.node::SyntaxNode: The node to examine.
Returns
- Nothing. The function modifies
resultsin place.
Arianna.Linter.collect_defined_structs! — Methodcollect_defined_structs!(results, node::SyntaxNode; new_types=nothing)Helper function that collects struct definitions from the given node into results.
Structs are included only if their declared type is found in REQUIRED_STRUCTS or in the optional new_types set.
Arguments
results: A mutable set to store(struct_name, struct_type)tuples.node::SyntaxNode: The node to examine.new_types: (optional) A set of additional struct types to recognize.
Returns
- Nothing. The function modifies
resultsin place.
Arianna.Linter.collect_defined_types! — Methodcollect_defined_types!(results::Set, node)Helper function that collects abstract type definitions from the given node into results.
Only types that inherit from AriannaSystem are collected.
Arguments
results::Set: A mutable set to store the collected abstract type names.node: ASyntaxNodeto examine.
Returns
- Nothing. The function modifies
resultsin place.
Arianna.Linter.construct_tree — Methodconstruct_tree(defined_structs, defined_functions, defined_types) -> NodeBuilds a hierarchical tree representing the structure of an Arianna system from the provided sets of defined structs, functions, and abstract types.
Each system node may contain actions, which in turn may contain policies, each requiring specific functions. The function cross-references the required functions with the provided defined_functions and marks them as defined if present.
Arguments
defined_structs: A set of tuples(struct_name::Symbol, struct_type::Symbol)representing defined structs.defined_functions: A set of tuples(function_name::Symbol, arity::Int, arg_types::Set{Symbol})representing defined functions.defined_types: A set ofSymbols representing additional abstract system types.
Returns
Node: The root node of the constructed tree.
Arianna.Linter.extract_defined_functions — Methodextract_defined_functions(code::String) -> Set{Tuple{Symbol, Int, Set{Symbol}}}Parses the provided Julia source code and returns a set of tuples describing each function definition found.
function_name: The name of the function as aSymbol.arity: The number of arguments the function takes (as anInt).arg_types: ASet{Symbol}representing the types of the arguments, if specified; otherwise, the set may be empty.
Arguments
code::String: The Julia source code to analyze.
Returns
Set{Tuple{Symbol, Int, Set{Symbol}}}: A set containing tuples of function names, arities, and sets of argument type symbols.
Arianna.Linter.extract_defined_structs — Methodextract_defined_structs(code::String; new_types=nothing) -> Set{Tuple{Symbol, Symbol}}Parses the provided Julia source code and returns a set of tuples (struct_name, struct_type) for each struct definition found.
struct_name: The name of the struct as aSymbol.struct_type: The type of the struct as aSymbol(e.g.,:mutable,:struct, or any custom type innew_types).
Arguments
code::String: The Julia source code to analyze.new_types: (optional) An iterable of additional struct types to recognize beyond the defaultstructandmutable struct.
Returns
Set{Tuple{Symbol, Symbol}}: A set containing tuples of struct names and their corresponding types.
Arianna.Linter.extract_defined_types — Methodextract_defined_types(code::String) -> Set{Symbol}Analyze the provided Julia source code and return a set containing the names of all abstract types defined within it.
Arguments
code::String: A string containing Julia source code to be parsed.
Returns
Set{Symbol}: A set of symbols, each representing the name of an abstract type defined in the code.
Arianna.Linter.run_linter — Methodconstruct_tree(defined_structs, defined_functions, defined_types) -> NodeBuilds a hierarchical tree representing the structure of an Arianna system from the provided sets of defined structs, functions, and abstract types.
Each system node may contain actions, which in turn may contain policies, each requiring specific functions. The function cross-references the required functions with the provided defined_functions and marks them as defined if present.
Arguments
defined_structs: A set of tuples(struct_name::Symbol, struct_type::Symbol)representing defined structs.defined_functions: A set of tuples(function_name::Symbol, arity::Int, arg_types::Set{Symbol})representing defined functions.defined_types: A set ofSymbols representing additional abstract system types.
Returns
Node: The root node of the constructed tree.