AbstractLiveSets
Functions
FreeBird.AbstractLiveSets
— ModuleAbstractLiveSets
Module for defining the livesets, which are collections of walkers that are used in the sampling schemes.
FreeBird.AbstractLiveSets.LJAtomWalkers
— Typestruct LJAtomWalkers <: AtomWalkers
The LJAtomWalkers
struct represents a collection of atom walkers that interact with each other using the Lennard-Jones potential.
Fields
walkers::Vector{AtomWalker{C}}
: A vector of atom walkers, whereC
is the number of components.lj_potential::LennardJonesParametersSets
: The Lennard-Jones potential parameters. SeeLennardJonesParametersSets
.
Constructor
LJAtomWalkers(walkers::Vector{AtomWalker{C}}, lj_potential::LennardJonesParametersSets; assign_energy=true)
: Constructs a newLJAtomWalkers
object with the given walkers and Lennard-Jones potential parameters. Ifassign_energy=true
, the energy of each walker is assigned using the Lennard-Jones potential.
FreeBird.AbstractLiveSets.LatticeGasWalkers
— Typestruct LatticeGasWalkers <: LatticeWalkers
The LatticeGasWalkers
struct represents a collection of lattice walkers for a lattice gas system. It is a subtype of LatticeWalkers
.
Fields
walkers::Vector{LatticeWalker{C}}
: A vector of lattice walkers.hamiltonian::LatticeGasHamiltonian
: The lattice gas Hamiltonian associated with the walkers.
Constructors
LatticeGasWalkers(walkers::Vector{LatticeWalker{C}}, hamiltonian::LatticeGasHamiltonian; assign_energy=true, perturb_energy::Float64=0.0)
: Constructs a newLatticeGasWalkers
object with the given walkers and Hamiltonian. Ifassign_energy
istrue
, the energy of each walker is assigned using the provided Hamiltonian. The optionalperturb_energy
parameter can be used to add a small perturbation to the assigned energy.
FreeBird.AbstractLiveSets.assign_energy!
— Methodassign_energy!(walker::AtomWalker, lj::LennardJonesParametersSets)
Assigns the energy to the given walker
using the Lennard-Jones parameters lj
.
Arguments
walker::AtomWalker
: The walker object to assign the energy to.lj::LennardJonesParametersSets
: The Lennard-Jones parameters.
Returns
walker::AtomWalker
: The walker object with the assigned energy.
FreeBird.AbstractLiveSets.assign_energy!
— Methodassign_energy!(walker::LatticeWalker{C}, hamiltonian::LatticeGasHamiltonian; perturb_energy::Float64=0.0)
Assigns energy to the given walker
based on the hamiltonian
. If perturb_energy
is non-zero, a small random perturbation is added to the energy.
Arguments
walker::LatticeWalker{C}
: The walker to assign energy to.hamiltonian::LatticeGasHamiltonian
: The Hamiltonian used to calculate the energy.perturb_energy::Float64=0.0
: The amount of random perturbation to add to the energy.
Returns
walker::LatticeWalker{C}
: The walker with the assigned energy.
FreeBird.AbstractLiveSets.assign_frozen_energy!
— Methodassign_frozen_energy!(walker::AtomWalker, lj::LennardJonesParametersSets)
Assigns the frozen energy to the given walker
using the Lennard-Jones parameters lj
.
Arguments
walker::AtomWalker
: The walker object to assign the energy to.lj::LennardJonesParametersSets
: The Lennard-Jones parameters.
Returns
walker::AtomWalker
: The walker object with the assigned energy.