MCRun
ParallelTemperingMonteCarlo.MCRun.acc_test! — Method
acc_test!(mc_state::MCState,ensemble::Etype,movetype::String) where Etype <: AbstractEnsembleacc_test! function now significantly contracted as a method of calculating the metropolis condition, comparing it to a random variable and if the condition is met using the swap_config! function to exchange the current mc_state with the internally defined new variables. ensemble and movetype dictate the exact calculation of the metropolis condition, and the internal potential_variables within the mcstates dictate how [`swapconfig!`](@ref) operates.
ParallelTemperingMonteCarlo.MCRun.check_e_bounds — Method
check_e_bounds(energy::Number, ebounds::VorS)Function to determine if an energy value is greater than or less than the min/max, used in equilibration cycle.
ParallelTemperingMonteCarlo.MCRun.equilibration — Method
equilibration(mc_states::MCStateVector, move_strat::MoveStrategy{N, E}, mc_params::MCParams, pot::Ptype, ensemble::Etype, n_steps::Int, results::Output, restart::Bool) where {N, E}While initialisation sets mc_states, params etc. we require something to thermalise our simulation and set the histograms. This function is mostly a wrapper for the equilibration_cycle! function that optionally removes the thermalisation from restart.
N.B. Restart is currently non-functional, do not try use it
ParallelTemperingMonteCarlo.MCRun.equilibration_cycle! — Method
equilibration_cycle!(mc_states::MCStateVector, move_strat::MoveStrategy{N, E}, mc_params::MCParams, pot::Ptype, ensemble::Etype, n_steps::Int, results::Output) where {N, E}Function to thermalise a set of mc_states ensuring that the number of equilibration cycles defined in mc_params are completed without updating the results before initialising the results struct according to the maximum and minimum energy determined throughout the equilibration cycle.
ParallelTemperingMonteCarlo.MCRun.get_energy! — Method
get_energy!(mc_state::MCState{T,N,BC,P,E},pot::PType,movetype::String) where PType <: AbstractPotential where {T,N,BC,P<:PotentialVariables,E<:NVTVariables}
get_energy!(mc_state::MCState{T,N,BC,P,E},pot::PType,movetype::String) where PType <: AbstractPotential where {T,N,BC,P<:PotentialVariables,E<:NPTVariables}
get_energy!(mc_state::MCState{T,N,BC,P,E},pot::PType,movetype::String) where PType <: AbstractPotential where {T,N,BC,P<:AbstractPotentialVariables,E<:NNVTVariables}Curry function designed to separate energy calculations into their respective ensembles and move types. Currently implemented for:
- NVT ensemble without r_cut
- NPT ensemble with r_cut
- NNVT ensemble for multiple-species atomsParallelTemperingMonteCarlo.MCRun.mc_cycle! — Method
mc_cycle!(mc_states::MCStateVector, move_strat::MoveStrategy{N, E}, mc_params::MCParams, pot::Ptype, ensemble::Etype, n_steps::Int, index::Int) where {N, E}
mc_cycle!(mc_states::MCStateVector, move_strat::MoveStrategy{N, E}, mc_params::MCParams, pot::Ptype, ensemble::Etype, n_steps::Int, results::Output, idx::Int, rdfsave::Bool) where {N, E}Basic function utilised by the simulation. For each of the n_steps run a single mc_step! on the mc_states according to pot, move_strat and ensemble, then complete the parallel_tempering_exchange! and update_step_size!.
Second method includes the sampling_step! which updates the results struct. The first method is used by the equilibration_cycle! and therefore does not update the results struct.
ParallelTemperingMonteCarlo.MCRun.mc_move! — Method
mc_move!(mc_state::MCState,move_strat::MoveStrategy{N,E},pot::Ptype,ensemble::Etype) where Ptype <: AbstractPotential where Etype <: AbstractEnsemble where {N,E}Basic move for one mc_state according to a move_strat dictating the types of moves allowed within the ensemble when moving across a pot defining the PES.
- Calculates an index for the move
- Generates either a volume or atom move depending on
movestrat[index] - Calculates energy based on the pot and new move
- Tests acc and swaps if relevant
ParallelTemperingMonteCarlo.MCRun.mc_step! — Method
mc_step!(mc_states::MCStateVector, move_strat::MoveStrategy{N, E}, pot::Ptype, ensemble::Etype, n_steps::Int) where {N, E}Distributes each state in mc_state to the mc_move! function in accordance with a move_strat, ensemble and pot.
ParallelTemperingMonteCarlo.MCRun.ptmc_run! — Method
(ptmc_run!(mc_params::MCParams, temp::TempGrid, start_config::Config, potential::Ptype, ensemble::Etype; rdfsave = false, restart = false, save = false, workingdirectory = pwd()) where Ptype <: AbstractPotential) where Etype <: AbstractEnsemble
ptmc_run!(restart::Bool; rdfsave = false, save = 1000, eq_cycles = 0.2)Main call for the ptmc program. Given mc_params dictating the number of cycles etc. the temps containing the temperature and beta values we aim to simulate, an initial start_config and the potential and ensemble we run a complete simulation, explicitly outputting the mc_states and results structs.
- Second method:
The second method relies on a series of checkpoint files -see Checkpoint module ReadSave- to autoinitialise an MC cycle. Still accepts restart as an argument to indicate whether this is a clean start with configs or a restart from a checkpoint at a given index.
- kwargs currently implemented are:
rdfsave::Bool: tells the simulation whether or not to generate and save radial distribution functions (a resource intensive step) – set to falserestart::Bool: tells the simulation whether or not we are beginning from a partially complete simulation - set false for method one.acc::Vector: sets the min and max acceptance rates used to adjust stepsize for the simulation - set [0.4 0.6] for a target of 40-60% acceptancesave::BoolorInt: tells the simulation whether to write checkpoints - set false for no save or integer expressing save frequency
ParallelTemperingMonteCarlo.MCRun.reset_counters — Method
reset_counters(state::MCState)After equilibration this resets the count stats to zero
ParallelTemperingMonteCarlo.MCRun.swap_atom_config! — Method
swap_atom_config!(mc_state::MCState{T, N, BC, P, E}, i_atom::Int, trial_pos::PositionVector) where {T, N, BC, P <: AbstractPotentialVariables, E <: AbstractEnsembleVariables}ParallelTemperingMonteCarlo.MCRun.swap_config! — Method
swap_config!(mc_state::MCState{T, N, BC, P, E}, movetype::String) where {T, N, BC, P <: AbstractPotentialVariables, E <: AbstractEnsembleVariables}Function for replacing the MC state and potential variables values with the updated values when metropolis condition is met. Implemented for the following move_type:
atommovevolumemoveswapmovefor atom swaps
All methods also call the swap_vars! function which distributes the appropriate mc_states.potential_variables values into the current mc_state struct.
ParallelTemperingMonteCarlo.MCRun.swap_config_v! — Method
swap_config_v!(mc_state::MCState, bc::CubicBC, trial_config::Config, new_dist2_mat::Matrix{N}, en_vec_new::VorS, new_en_tot::Number) where N <: Number
swap_config_v!(mc_state::MCState, bc::RhombicBC, trial_config::Config, new_dist2_mat::Matrix{N}, en_vec_new::VorS, new_en_tot::Number) where N <: NumberSwaps mc_states and ensemble variables in case of accepted volume move for NPT ensemble. Implemented for CubicBC and RhombicBC
ParallelTemperingMonteCarlo.MCRun.swap_move_config! — Method
swap_move_config!(mc_state::MCState,indices::VorS)Function designed to exchange relevant variables when swapping an atom. Accepts the mc_state and the atom indices and exchanges atom indices[1] with atom indices[2]
ParallelTemperingMonteCarlo.MCRun.swap_vars! — Method
swap_vars!(i_atom::Int, potential_variables::V) where V <: DimerPotentialVariables
swap_vars!(i_atom::Int, potential_variables::ELJPotentialBVariables)
swap_vars!(i_atom::Int, potential_variables::EmbeddedAtomVariables)
swap_vars!(i_atom::Int, potential_variables::NNPVariables)
swap_vars!(i_atom::Int, potential_variables::NNPVariables2a)Called by swap_atom_config! function; takes the appropriate potential_variables that are specific to the potential energy surface under consideration and replaces the current values with the new values such as:
- Under magnetic fields, the new tan matrix replaces the old
- In the EAM, we replace the rho and phi vectors with the new updated versions
- Using an NNP we require the new G matrix and F matrix to replace the old versions.
Implemented for potential variables: