ReadSave
ParallelTemperingMonteCarlo.ReadSave — Module
CheckpointModule designed to save relevant parameters and configurations throughout the simulation to allow restarting.
ParallelTemperingMonteCarlo.ReadSave.build_states — Method
build_states(mc_params::MCParams, ensemble::AbstractEnsemble, temp::TempGrid, potential::AbstractPotential)For use initialising states and outputs when NOT restarting, but beginning from files. Builds empty Output struct named results and a vector of mc_states using either: one configuration stored in config.data OR a series of configurations stored in config.i. NB if config.i doesn't exist the default will be config.1. In this way states can be initialised with different starting configurations.
ParallelTemperingMonteCarlo.ReadSave.checkpoint — Method
checkpoint(index::Int, mc_states::MCStateVector, results::Output, ensemble::AbstractEnsemble, rdfsave::Bool)
checkpoint(index::Int, mc_states::MCStateVector, results::Output, ensemble::NPT, rdfsave::Bool)Function to save relevant information about the current state of the system at step index. Saves the configurations in each mc_state save_configs as well as the histograms stored in results. Optionally stores the volume histograms if using the NPT ensemble and the radial distribution functions if desired.
ParallelTemperingMonteCarlo.ReadSave.checkpoint_config — Method
checkpoint_config(savefile::SaveFile, state::MCState{T, N, BC, Ptype, Etype}) where {T, N, BC <: SphericalBC, Ptype, Etype}
checkpoint_config(savefile::SaveFile, state::MCState{T, N, BC, Ptype, Etype}) where {T, N, BC <: CubicBC, Ptype, Etype}
checkpoint_config(savefile::SaveFile, state::MCState{T, N, BC, Ptype, Etype}) where {T, N, BC <: RhombicBC, Ptype, Etype}Function writes a single config in the standard xyz format. N atoms, the comment line contains the boundary condition information (implemented for Spherical BC and both types of Periodic BC) as well as max_displ information determining the stepsize used at the current step of the monte carlo simulation. The comment row is followed by 1 as a placeholder for the atom type to be implemented in future and the positions x,y,z in order.
ParallelTemperingMonteCarlo.ReadSave.read_checkpoint_config — Method
read_checkpoint_config(xyzdata)Function designed to take a single xyz-style checkpoint file and return the configuration and max displacement data associated with a saved mc_state. This is used to reconstruct an MC simulation from checkpoints.
ParallelTemperingMonteCarlo.ReadSave.read_config — Method
read_config(xyzdata)Designed to read in one xyz-style file with one configuration and return this for starting a simulation from files without restarting.
ParallelTemperingMonteCarlo.ReadSave.read_init — Method
read_init()Function to reinitialise the fixed parameters of the MC simulation as saved by the save_init function.
ParallelTemperingMonteCarlo.ReadSave.read_params — Method
read_params(paramsvec)
read_params(paramsvec,restart)Function to turn a delimited paramsvec into an MCParams and TempGrid struct. Second method includes a bool restart to determine whether or not to set eq_cycles=0 or 0.2*cycles if not restarting
ParallelTemperingMonteCarlo.ReadSave.readensemble — Method
readensemble(ensemblevec)Function to convert delimited file contents ensemblevec and convert them into an ensemble.
ParallelTemperingMonteCarlo.ReadSave.readpotential — Method
readpotential(potinfovec)Function to convert delimited file contents potinfovec into a potential. Implemented for:
ParallelTemperingMonteCarlo.ReadSave.rebuild_states — Method
rebuild_states(n_traj::Int, ensemble::AbstractEnsemble, temps::TempGrid, potential::AbstractPotential)Function to rebuild the MCStates vector and results struct from checkpoint information. The ensemble temps and potential along with n_traj are reconstructed elsewhere, but required to accurately recreate the states.
ParallelTemperingMonteCarlo.ReadSave.save_configs — Method
save_configs(mc_states::MCStateVector)Function to save the configuration of each state in a vector of mc_states. Writes each configuration according to checkpoint_config into a file config.i where i indicates the order of the states.
ParallelTemperingMonteCarlo.ReadSave.save_histparams — Method
save_histparams(results)Initialises and populates a data file containing the information necessary to interpret histogram data.
ParallelTemperingMonteCarlo.ReadSave.save_init — Method
save_init(potential,ensemble,params,temp)Function to write all static parameters into a single parameters file. If a params file does not exist, it is created in ./checkpoint as params.data. This contains the mc_params ensemble and potential data.
ParallelTemperingMonteCarlo.ReadSave.setresults — Method
setresults(histparams,histdata,histv_data,r2data)Function to re-initialise the results struct on restarting a simulation.
ParallelTemperingMonteCarlo.ReadSave.writeensemble — Method
writeensemble(savefile::SaveFile, ensemble::NVT)
writeensemble(savefile::SaveFile, ensemble::NPT)
writeensemble(savefile::SaveFile, ensemble::NNVT)Function to write the ensemble data into a savefile including the move types. First method is for the NVT ensemble which does not include volume changes, second method is NPT ensemble and does inclue volume moves.
ParallelTemperingMonteCarlo.ReadSave.writeparams — Method
writeparams(savefile,params,temp)Function to write the mc_params and temp data into a savefile. These are static parameters that define how the simulation is to proceed such as the number of cycles, trajectories and the temperatures to be covered.
ParallelTemperingMonteCarlo.ReadSave.writepotential — Method
writepotential(savefile::SaveFile, potential::Ptype) where Ptype <: AbstractDimerPotential
writepotential(savefile::SaveFile, potential::Ptype) where Ptype <: AbstractDimerPotentialB
writepotential(savefile::SaveFile, potential::Ptype) where Ptype <: EmbeddedAtomPotential
writepotential(savefile::SaveFile, potential::Ptype) where Ptype <: AbstractMachineLearningPotentialFunction to write potential surface information into savefile. implemented methods are the Embedded Atom Model, Extended Lennard-Jones and ELJ in Magnetic Field. This does not work for machine learning potentials.