Sampling Schemes
Functions
FreeBird.SamplingSchemes
— ModuleSamplingSchemes
Module for defining the sampling schemes.
FreeBird.SamplingSchemes.LatticeNestedSamplingParameters
— Typemutable struct LatticeNestedSamplingParameters <: SamplingParameters
The LatticeNestedSamplingParameters
struct represents the parameters used in the lattice nested sampling scheme.
Fields
mc_steps::Int64
: The number of total Monte Carlo moves to perform.energy_perturbation::Float64
: The energy perturbation used in the sampling process.fail_count::Int64
: The number of failed MC moves in a row.allowed_fail_count::Int64
: The maximum number of failed MC moves allowed before resetting the step size.random_seed::Int64
: The seed for the random number generator.
FreeBird.SamplingSchemes.MCMixedMoves
— Typestruct MCMixedMoves <: MCRoutine
A type for generating a new walker by performing random walks and swapping atoms. Currently, it is intended to use this routine for multi-component systems. The actual number of random walks and swaps to perform is determined by the weights of the fields walks_freq
and swaps_freq
. For example, if walks_freq=4
and swaps_freq=1
, then the probability of performing a random walk is 4/5, and the probability of performing a swap is 1/5.
Fields
walks_freq::Int
: The frequency of random walks to perform.swaps_freq::Int
: The frequency of atom swaps to perform.
FreeBird.SamplingSchemes.MCNewSample
— Typestruct MCNewSample <: MCRoutine
A type for generating a new walker from a random configuration. Currently, it is intended to use this routine for lattice gas systems.
FreeBird.SamplingSchemes.MCRandomWalkClone
— Typestruct MCRandomWalkClone <: MCRoutine
A type for generating a new walker by cloning an existing walker and performing a random walk for decorrelation.
FreeBird.SamplingSchemes.MCRandomWalkMaxE
— Typestruct MCRandomWalkMaxE <: MCRoutine
A type for generating a new walker by performing a random walk for decorrelation on the highest-energy walker.
FreeBird.SamplingSchemes.MCRejectionSampling
— Typestruct MCRejectionSampling <: MCRoutine
A type for generating a new walker by performing rejection sampling. Currently, it is intended to use this routine for lattice gas systems.
FreeBird.SamplingSchemes.MCRoutine
— Typeabstract type MCRoutine
An abstract type representing a Monte Carlo routine.
Currently, the following concrete types are supported:
MCRandomWalkMaxE
: A type for generating a new walker by performing a random walk for decorrelation on the
highest-energy walker.
MCRandomWalkClone
: A type for generating a new walker by cloning an existing walker and performing a random walk
for decorrelation.
MCNewSample
: A type for generating a new walker from a random configuration. Currently, it is intended to use
this routine for lattice gas systems.
MCMixedMoves
: A type for generating a new walker by performing random walks and swapping atoms. Currently, it is
intended to use this routine for multi-component systems. The actual number of random walks and swaps to perform is determined by the weights of the fields walks_freq
and swaps_freq
. See MCMixedMoves
.
MCRejectionSampling
: A type for generating a new walker by performing rejection sampling. Currently, it is intended
to use this routine for lattice gas systems.
FreeBird.SamplingSchemes.MetropolisMCParameters
— TypeMetropolisMCParameters <: SamplingParameters
Parameters for the Metropolis Monte Carlo algorithm.
Fields
temperature::Float64
: The temperature of the system.equilibrium_steps::Int64
: The number of steps to equilibrate the system.sampling_steps::Int64
: The number of steps to sample the system.step_size::Float64
: The step size for the random walk (for atomistic systems).step_size_lo::Float64
: The lower bound of the step size.step_size_up::Float64
: The upper bound of the step size.accept_range::Tuple{Float64, Float64}
: The range of acceptance rates for adjusting the step size.
e.g. (0.25, 0.75) means that the step size will decrease if the acceptance rate is below 0.25 and increase if it is above 0.75.
random_seed::Int64
: The seed for the random number generator.
FreeBird.SamplingSchemes.NestedSamplingParameters
— Typemutable struct NestedSamplingParameters <: SamplingParameters
The NestedSamplingParameters
struct represents the parameters used in the nested sampling scheme.
Fields
mc_steps::Int64
: The number of total Monte Carlo moves to perform.initial_step_size::Float64
: The initial step size, which is the fallback step size if MC routine fails to accept a move.step_size::Float64
: The on-the-fly step size used in the sampling process.step_size_lo::Float64
: The lower bound of the step size.step_size_up::Float64
: The upper bound of the step size.accept_range::Tuple{Float64, Float64}
: The range of acceptance rates for adjusting the step size.
e.g. (0.25, 0.75) means that the step size will decrease if the acceptance rate is below 0.25 and increase if it is above 0.75.
fail_count::Int64
: The number of failed MC moves in a row.allowed_fail_count::Int64
: The maximum number of failed MC moves allowed before resetting the step size.random_seed::Int64
: The seed for the random number generator.
FreeBird.SamplingSchemes.SamplingParameters
— Typestruct NestedSamplingParameters
The NestedSamplingParameters
struct represents the parameters for various sampling algorithm.
FreeBird.SamplingSchemes.WangLandauParameters
— TypeWangLandauParameters
A structure to hold the parameters for the Wang-Landau sampling scheme.
Fields
num_steps::Int64
: The number of Monte Carlo steps.flatness_criterion::Float64
: The criterion for flatness of the histogram.f_initial::Float64
: The initial modification factor.f_min::Float64
: The minimum modification factor.energy_bins::Vector{Float64}
: The pre-supplied energy bins.max_iter::Int64
: The maximum number of iterations in each flatness check.step_size::Float64
: The step size for the random walk (for atomistic systems).random_seed::Int64
: The seed for the random number generator.
FreeBird.SamplingSchemes.WangLandauParameters
— MethodWangLandauParameters(;
num_steps::Int64=100,
flatness_criterion::Float64=0.8,
f_initial::Float64=Float64(MathConstants.e),
f_min::Float64=exp(1e-8),
energy_min::Float64=0.0,
energy_max::Float64=1.0,
num_energy_bins::Int64=100,
max_iter::Int64=1000,
step_size::Float64=0.01,
random_seed::Int64=1234
)
Create a WangLandauParameters
object with the specified parameters.
Arguments
num_steps::Int64
: The number of Monte Carlo steps.flatness_criterion::Float64
: The criterion for flatness of the histogram.f_initial::Float64
: The initial modification factor.f_min::Float64
: The minimum modification factor.energy_min::Float64
: The minimum energy.energy_max::Float64
: The maximum energy.num_energy_bins::Int64
: The number of energy bins.max_iter::Int64
: The maximum number of iterations in each flatness check.random_seed::Int64
: The seed for the random number generator.
Returns
WangLandauParameters
: The parameters for the Wang-Landau sampling scheme.
FreeBird.SamplingSchemes.adjust_step_size
— Methodadjust_step_size(params::SamplingParameters, rate::Float64)
Adjusts the step size of the sampling algorithm based on the acceptance rate. The step size is increased by 10% if the acceptance rate is above the upper limit of the range, and decreased by 10% if the acceptance rate is below the lower limit of the range.
Arguments
params::SamplingParameters
: The parameters of the sampling algorithm.rate::Float64
: The acceptance rate of the algorithm.range::Tuple{Float64, Float64}
: The range of acceptance rates for adjusting the step size. Default is (0.25, 0.75).
Returns
params::SamplingParameters
: The updated parameters with adjusted step size.
FreeBird.SamplingSchemes.exact_enumeration
— Methodexact_enumeration(lattice::SLattice{G}, cutoff_radii::Tuple{Float64, Float64}, h::ClassicalHamiltonian) where G
Enumerate all possible configurations of a lattice system and compute the energy of each configuration.
Arguments
lattice::SLattice{G}
: The (starting) lattice system to enumerate. All possible configurations will be generated from this lattice system.h::ClassicalHamiltonian
: The Hamiltonian containing the on-site and nearest-neighbor interaction energies.
Returns
DataFrame
: A DataFrame containing the energy and configuration of each configuration.LatticeGasWalkers
: A collection of lattice walkers for each configuration.
FreeBird.SamplingSchemes.monte_carlo_sampling
— Methodmonte_carlo_sampling(
lattice::AbstractLattice,
h::ClassicalHamiltonian,
mc_params::MetropolisMCParameters
)
Perform the Metropolis Monte Carlo sampling algorithm for a range of temperatures.
Note: The Boltzmann constant is set to 8.617333262e-5 eV K$^{-1}$. Thus, the units of the temperature should be in Kelvin, and the units of the energy should be in eV (defined in the Hamiltonian).
Arguments
lattice::AbstractLattice
: The initial lattice configuration.h::ClassicalHamiltonian
: The Hamiltonian containing the on-site and nearest-neighbor interaction energies.mc_params::MetropolisMCParameters
: The parameters for the Metropolis Monte Carlo algorithm.
Returns
energies::Vector{Float64}
: The energies of the system at each temperature.configs::Vector{typeof(lattice)}
: The configurations of the system at each temperature.cvs::Vector{Float64}
: The heat capacities of the system at each temperature.acceptance_rates::Vector{Float64}
: The acceptance rates of the system at each temperature.
FreeBird.SamplingSchemes.monte_carlo_sampling
— Methodmonte_carlo_sampling(
at::AtomWalker,
lj::LennardJonesParametersSets,
mc_params::MetropolisMCParameters
)
Perform the Metropolis Monte Carlo sampling algorithm for a range of temperatures.
Note: The Boltzmann constant is set to 8.617333262e-5 eV K$^{-1}$. Thus, the units of the temperature should be in Kelvin, and the units of the energy should be in eV.
Arguments
at::AtomWalker
: The initial atom walker configuration.lj::LennardJonesParametersSets
: The Lennard-Jones parameters.mc_params::MetropolisMCParameters
: The parameters for the Metropolis Monte Carlo algorithm.
Returns
energies::Vector{Float64}
: The energies of the system at each temperature.configs::Vector{typeof(at)}
: The configurations of the system at each temperature.cvs::Vector{Float64}
: The heat capacities of the system at each temperature.acceptance_rates::Vector{Float64}
: The acceptance rates of the system at each temperature.
FreeBird.SamplingSchemes.nested_sampling
— Methodnested_sampling(liveset::AtomWalkers, ns_params::NestedSamplingParameters, n_steps::Int64, mc_routine::MCRoutine; args...)
Perform a nested sampling loop for a given number of steps.
Arguments
liveset::AtomWalkers
: The initial set of walkers.ns_params::NestedSamplingParameters
: The parameters for nested sampling.n_steps::Int64
: The number of steps to perform.mc_routine::MCRoutine
: The Monte Carlo routine to use.
Returns
df
: A DataFrame containing the iteration number and maximum energy for each step.liveset
: The updated set of walkers.ns_params
: The updated nested sampling parameters.
FreeBird.SamplingSchemes.nested_sampling
— Methodnested_sampling(liveset::LatticeGasWalkers, ns_params::LatticeNestedSamplingParameters, n_steps::Int64, mc_routine::MCRoutine, save_strategy::DataSavingStrategy)
Perform a nested sampling loop on a lattice gas system for a given number of steps.
Arguments
liveset::LatticeGasWalkers
: The initial set of walkers.ns_params::LatticeNestedSamplingParameters
: The parameters for nested sampling.n_steps::Int64
: The number of steps to perform.mc_routine::MCRoutine
: The Monte Carlo routine to use.save_strategy::DataSavingStrategy
: The strategy for saving data.
Returns
df
: A DataFrame containing the iteration number and maximum energy for each step.liveset
: The updated set of walkers.ns_params
: The updated nested sampling parameters.
FreeBird.SamplingSchemes.nested_sampling_step!
— Methodnested_sampling_step!(liveset::AtomWalkers, ns_params::NestedSamplingParameters, mc_routine::MCMixedMoves)
Perform a single step of the nested sampling algorithm using the Monte Carlo mixed moves routine.
Arguments
liveset::AtomWalkers
: The set of atom walkers.ns_params::NestedSamplingParameters
: The parameters for nested sampling.mc_routine::MCMixedMoves
: The Monte Carlo mixed moves routine.
Returns
iter
: The iteration number after the step.emax
: The highest energy recorded during the step.liveset
: The updated set of atom walkers.ns_params
: The updated nested sampling parameters.
FreeBird.SamplingSchemes.nested_sampling_step!
— Methodnested_sampling_step!(liveset::AtomWalkers, ns_params::NestedSamplingParameters, mc_routine::MCRoutine)
Perform a single step of the nested sampling algorithm using the Monte Carlo random walk routine.
Arguments
liveset::AtomWalkers
: The set of atom walkers.ns_params::NestedSamplingParameters
: The parameters for nested sampling.mc_routine::MCRoutine
: The Monte Carlo routine for generating new samples. SeeMCRoutine
.
Returns
iter
: The iteration number after the step.emax
: The highest energy recorded during the step.liveset
: The updated set of atom walkers.ns_params
: The updated nested sampling parameters.
FreeBird.SamplingSchemes.nested_sampling_step!
— Methodnested_sampling_step!(liveset::LatticeGasWalkers, ns_params::LatticeNestedSamplingParameters, mc_routine::MCNewSample)
Perform a single step of the nested sampling algorithm.
This function takes a liveset
of lattice gas walkers, ns_params
containing the parameters for nested sampling, and mc_routine
representing the Monte Carlo routine for generating new samples. It performs a single step of the nested sampling algorithm by updating the liveset with a new walker.
Arguments
liveset::LatticeGasWalkers
: The liveset of lattice gas walkers.ns_params::LatticeNestedSamplingParameters
: The parameters for nested sampling.mc_routine::MCNewSample
: The Monte Carlo routine for generating new samples.
Returns
iter
: The iteration number of the liveset after the step.emax
: The maximum energy of the liveset after the step.liveset::LatticeGasWalkers
: The updated liveset after the step.ns_params::LatticeNestedSamplingParameters
: The updated nested sampling parameters after the step.
FreeBird.SamplingSchemes.nested_sampling_step!
— Methodnested_sampling_step!(liveset::LatticeGasWalkers, ns_params::LatticeNestedSamplingParameters, mc_routine::MCRoutine)
Perform a single step of the nested sampling algorithm.
This function takes a liveset
of lattice gas walkers, ns_params
containing the parameters for nested sampling, and mc_routine
representing the Monte Carlo routine for generating new samples. It performs a single step of the nested sampling algorithm by updating the liveset with a new walker.
Arguments
liveset::LatticeGasWalkers
: The liveset of lattice gas walkers.ns_params::LatticeNestedSamplingParameters
: The parameters for nested sampling.mc_routine::MCRoutine
: The Monte Carlo routine for generating new samples.
Returns
iter
: The iteration number of the liveset after the step.emax
: The maximum energy of the liveset after the step.
FreeBird.SamplingSchemes.nvt_monte_carlo
— Methodnvt_monte_carlo(
lattice::AbstractLattice,
h::ClassicalHamiltonian,
temperature::Float64,
num_steps::Int64,
random_seed::Int64
)
Perform the NVT Monte Carlo algorithm to sample the lattice configurations.
Note: The Boltzmann constant is set to 8.617333262e-5 eV K$^{-1}$. Thus, the units of the temperature should be in Kelvin, and the units of the energy should be in eV (defined in the Hamiltonian).
Arguments
lattice::AbstractLattice
: The initial lattice configuration.h::ClassicalHamiltonian
: The Hamiltonian containing the on-site and nearest-neighbor interaction energies.temperature::Float64
: The temperature of the system.num_steps::Int64
: The number of Monte Carlo steps.random_seed::Int64
: The seed for the random number generator.
Returns
energies::Vector{Float64}
: The energies of the system at each step.configurations::Vector{typeof(lattice)}
: The configurations of the system at each step.accepted_steps::Int64
: The number of accepted steps.
FreeBird.SamplingSchemes.sort_by_energy!
— Methodsort_by_energy!(liveset::LJAtomWalkers)
Sorts the walkers in the liveset by their energy in descending order.
Arguments
liveset::LJAtomWalkers
: The liveset of walkers to be sorted.
Returns
liveset::LJAtomWalkers
: The sorted liveset.
FreeBird.SamplingSchemes.update_iter!
— Methodupdate_iter!(liveset::AtomWalkers)
Update the iteration count for each walker in the liveset.
Arguments
liveset::AtomWalkers
: The set of walkers to update.
FreeBird.SamplingSchemes.wang_landau
— Methodwang_landau(
lattice::AbstractLattice,
h::ClassicalHamiltonian,
wl_params::WangLandauParameters
)
wang_landau(
walker::AtomWalker,
lj::LennardJonesParametersSets,
wl_params::WangLandauParameters
)
Perform the Wang-Landau sampling scheme for a lattice or an atomistic system.
Arguments
lattice::AbstractLattice
/walker::AtomWalker
: The initial lattice/atomistic configuration.h::ClassicalHamiltonian
/lj::LennardJonesParametersSets
: The Hamiltonian parameters for the lattice/atomistic system.wl_params::WangLandauParameters
: The parameters for the Wang-Landau sampling scheme.
Returns
df::DataFrame
/energies::Vector{Float64}
: The energies of the system at each step.configs::Vector{AbstractLattice}
/configs::Vector{AtomWalker}
: The configurations of the system at each step.wl_params::WangLandauParameters
: The parameters for the Wang-Landau sampling scheme.S::Vector{Float64}
: The entropy of the system.H::Vector{Int64}
: The histogram of the system.