AbstractWalkers

Functions

FreeBird.AbstractWalkers.AtomWalkerType
mutable struct AtomWalker

The AtomWalker struct represents a walker composed of atoms/molecules.

Fields

  • configuration::FastSystem: The configuration of the walker.
  • energy::typeof(0.0u"eV"): The energy of the walker.
  • iter::Int64: The current iteration number of the walker.
  • list_num_par::Vector{Int64}: The list of the number of particles for each component.
  • frozen::Vector{Bool}: A boolean vector indicating whether each component is frozen or not.
  • energy_frozen_part::typeof(0.0u"eV"): The energy of the frozen particles in the walker, serves as a constant energy offset to the interacting part of the system.

Constructor

  • AtomWalker(configuration::FastSystem; energy=0.0u"eV", iter=0, list_num_par=zeros(Int,C), frozen=zeros(Bool,C), energy_frozen_part=0.0u"eV"): Constructs a new AtomWalker object with the given configuration and optional parameters.
source
FreeBird.AbstractWalkers.AtomWalkerMethod
AtomWalker(configuration::FastSystem; freeze_species::Vector{Symbol}=Symbol[], merge_same_species=true)

Constructs an AtomWalker object with the given configuration.

Arguments

  • configuration::FastSystem: The configuration of the walker.
  • freeze_species::Vector{Symbol}: A vector of species to freeze.
  • merge_same_species::Bool: A boolean indicating whether to merge the same species into one component.

Returns

  • AtomWalker{C}: The constructed AtomWalker object.

Example

julia> at = FreeBirdIO.generate_multi_type_random_starting_config(10.0,[2,1,3,4,5,6];particle_types=[:H,:O,:H,:Fe,:Au,:Cl])
FastSystem(Au₅Cl₆Fe₄H₅O, periodic = FFF):
    bounding_box      : [ 5.94392        0        0;
                                0  5.94392        0;
                                0        0  5.94392]u"Å"

        .--------------.  
       /|Fel           |  
      / H   H   Cl     |  
     /  Hu   O         |  
    *   |       Au   Fe|  
    |   |FeCl        Fe|  
    |   |        Au    |  
    |   .---------Au---.  
    |  /           H  /   
    | Au Cl          /    
    |/              /     
    *--------------*      

julia> AtomWalker(at;freeze_species=[:H],merge_same_species=false)
AtomWalker{6}(FastSystem(Au₅Cl₆Fe₄H₅O, periodic = FFF, bounding_box = [[5.943921952763129, 0.0, 0.0], [0.0, 5.943921952763129, 0.0], [0.0, 0.0, 5.943921952763129]]u"Å"), 0.0 eV, 0, [2, 3, 1, 6, 4, 5], Bool[1, 1, 0, 0, 0, 0], 0.0 eV)

julia> AtomWalker(at;freeze_species=[:H],merge_same_species=true)
AtomWalker{5}(FastSystem(Au₅Cl₆Fe₄H₅O, periodic = FFF, bounding_box = [[5.943921952763129, 0.0, 0.0], [0.0, 5.943921952763129, 0.0], [0.0, 0.0, 5.943921952763129]]u"Å"), 0.0 eV, 0, [5, 1, 6, 4, 5], Bool[1, 0, 0, 0, 0], 0.0 eV)
source
FreeBird.AbstractWalkers.LatticeGeometryType
abstract type LatticeGeometry

The LatticeGeometry abstract type represents the geometry of a lattice. It has the following subtypes:

  • SquareLattice: A square lattice.
  • TriangularLattice: A triangular lattice.
  • GenericLattice: A generic lattice. Currently used for non-square and non-triangular lattices.
source
FreeBird.AbstractWalkers.LatticeWalkerType
mutable struct LatticeWalker

The LatticeWalker struct represents a walker on a 3D lattice.

Fields

  • configuration::AbstractLattice: The configuration of the walker.
  • energy::Float64: The energy of the walker.
  • iter::Int64: The current iteration number of the walker.

Constructor

LatticeWalker(configuration::AbstractLattice; energy=0.0, iter=0)

Create a new LatticeWalker with the given configuration and optional energy and iteration number.

source
FreeBird.AbstractWalkers.MLatticeType
mutable struct MLattice{C,G}

A mutable struct representing a lattice with the following fields:

  • lattice_vectors::Matrix{Float64}: The lattice vectors defining the unit cell.
  • positions::Matrix{Float64}: The positions of the lattice points.
  • basis::Vector{Tuple{Float64, Float64, Float64}}: The basis vectors within the unit cell.
  • supercell_dimensions::Tuple{Int64, Int64, Int64}: The dimensions of the supercell.
  • periodicity::Tuple{Bool, Bool, Bool}: The periodicity in each dimension.
  • components::Vector{Vector{Bool}}: The components of the lattice.
  • neighbors::Vector{Vector{Vector{Int}}}: The neighbors of each lattice point.
  • adsorptions::Vector{Bool}: The adsorption sites on the lattice.

Inner Constructor

MLattice{C,G}(
    lattice_vectors::Matrix{Float64},
    basis::Vector{Tuple{Float64, Float64, Float64}},
    supercell_dimensions::Tuple{Int64, Int64, Int64},
    periodicity::Tuple{Bool, Bool, Bool},
    components::Vector{Vector{Bool}},
    adsorptions::Vector{Bool},
    cutoff_radii::Vector{Float64},
) where {C,G}

Creates an MLattice instance with the specified parameters. The constructor performs the following steps:

  1. Validates that the number of components matches the expected value C.
  2. Computes the positions of the lattice points using lattice_positions.
  3. Computes the supercell lattice vectors.
  4. Computes the neighbors of each lattice point using compute_neighbors.

Throws an ArgumentError if the number of components does not match C.

Outer Constructors

MLattice{C,SquareLattice}(; lattice_constant::Float64=1.0,
                           basis::Vector{Tuple{Float64,Float64,Float64}}=[(0.0, 0.0, 0.0)],
                           supercell_dimensions::Tuple{Int64,Int64,Int64}=(4, 4, 1),
                           periodicity::Tuple{Bool,Bool,Bool}=(true, true, false),
                           cutoff_radii::Vector{Float64}=[1.1, 1.5],
                           components::Union{Vector{Vector{Int64}},Vector{Vector{Bool}},Symbol}=:equal,
                           adsorptions::Union{Vector{Int},Symbol}=:full)

MLattice{C,TriangularLattice}(; lattice_constant::Float64=1.0,
                              basis::Vector{Tuple{Float64,Float64,Float64}}=[(0.0, 0.0, 0.0),(1/2, sqrt(3)/2, 0.0)],
                              supercell_dimensions::Tuple{Int64,Int64,Int64}=(4, 2, 1),
                              periodicity::Tuple{Bool,Bool,Bool}=(true, true, false),
                              cutoff_radii::Vector{Float64}=[1.1, 1.5],
                              components::Union{Vector{Vector{Int64}},Vector{Vector{Bool}},Symbol}=:equal,
                              adsorptions::Union{Vector{Int},Symbol}=:full)

Constructs a square/triangular lattice with the specified parameters. The components and adsorptions arguments can be a vector of integers specifying the indices of the occupied sites, or a symbol. If components is :equal, the lattice is divided into C equal components when possible, or nearest to equal components otherwise. If adsorptions is :full, all sites are classified as adsorption sites.

Returns

  • MLattice{C,G}: A square/triangular lattice object with C components.
source
FreeBird.AbstractWalkers.check_num_componentsMethod
check_num_components(C::Int, list_num_par::Vector{Int}, frozen::Vector{Bool})

Check that the number of components matches the length of the list of number of particles and frozen particles.

Arguments

  • C::Int: The number of components.
  • list_num_par::Vector{Int}: The number of particles in each component.
  • frozen::Vector{Bool}: A vector indicating whether each component is frozen.
source
FreeBird.AbstractWalkers.lattice_positionsMethod

latticepositions(latticevectors::Matrix{Float64}, basis::Vector{Tuple{Float64, Float64, Float64}}, supercell_dimensions::Tuple{Int64, Int64, Int64})

Compute the positions of atoms in a 3D lattice.

Arguments

  • lattice_vectors::Matrix{Float64}: The lattice vectors of the system.
  • basis::Vector{Tuple{Float64, Float64, Float64}}: The basis of the system.
  • supercell_dimensions::Tuple{Int64, Int64, Int64}: The dimensions of the supercell.

Returns

  • positions::Matrix{Float64}: The positions of the atoms in the supercell.
source
FreeBird.AbstractWalkers.mlattice_setupMethod
mlattice_setup(C::Int, 
                 basis::Vector{Tuple{Float64, Float64, Float64}},
                 supercell_dimensions::Tuple{Int64, Int64, Int64},
                 components::Union{Vector{Vector{Int64}},Vector{Vector{Bool}},Symbol},
                 adsorptions::Union{Vector{Int}, Vector{Bool}, Symbol})

Setup the components and adsorptions for a lattice.

Arguments

  • C::Int: The number of components.
  • basis::Vector{Tuple{Float64, Float64, Float64}}: The basis of the lattice.
  • supercell_dimensions::Tuple{Int64, Int64, Int64}: The dimensions of the supercell.
  • components::Union{Vector{Vector{Int64}},Vector{Vector{Bool}},Symbol}: The components of the lattice.
  • adsorptions::Union{Vector{Int}, Vector{Bool}, Symbol}: The adsorption sites on the lattice.

Returns

  • lattice_comp::Vector{Vector{Bool}}: The components of the lattice.
  • lattice_adsorptions::Vector{Bool}: The adsorption sites on the lattice.
source
FreeBird.AbstractWalkers.sort_components_by_atomic_numberMethod
sort_components_by_atomic_number(at::AbstractSystem; merge_same_species=true)

Sorts the components of an AbstractSystem object at by their atomic number.

Arguments

  • at::AbstractSystem: The input AbstractSystem object.

Keyword Arguments

  • merge_same_species::Bool=true: Whether to merge components with the same species.

Returns

  • list_num_par::Vector{Int64}: A vector containing the number of each component species.
  • new_list::FastSystem: A new FastSystem object with the sorted components.

The function first extracts the atomic numbers of the components in at. If merge_same_species is true, it sorts the unique species and counts the number of each species. If merge_same_species is false, it creates a list of species and their counts. It then sorts the species and counts by atomic number. Finally, it constructs a new FastSystem object with the sorted components and returns the list of species counts and the new FastSystem object.

Examples

julia> at = FreeBirdIO.generate_multi_type_random_starting_config(10.0,[2,1,3,4,5,6];particle_types=[:H,:O,:H,:Fe,:Au,:Cl])
FastSystem(Au₅Cl₆Fe₄H₅O, periodic = FFF):
    bounding_box      : [ 5.94392        0        0;
                                0  5.94392        0;
                                0        0  5.94392]u"Å"

        .--------------.  
       /|     Cl    H  |  
      / |      Fe      |  
     /  |  Au     FeH  |  
    *   |   FeH  ACl   |  
    |   |    Cl     Au |  
    |   |            O |  
    |   .--Fe----------.  
    |  /H  Cl         /   
    | /          Au  /    
    |/Cl          Cl/     
    *--------------*      


julia> AbstractWalkers.sort_components_by_atomic_number(at; merge_same_species=false)
([2, 3, 1, 6, 4, 5], FastSystem(Au₅Cl₆Fe₄H₅O, periodic = FFF, bounding_box = [[5.943921952763129, 0.0, 0.0], [0.0, 5.943921952763129, 0.0], [0.0, 0.0, 5.943921952763129]]u"Å"))

julia> AbstractWalkers.sort_components_by_atomic_number(at)
([5, 1, 6, 4, 5], FastSystem(Au₅Cl₆Fe₄H₅O, periodic = FFF, bounding_box = [[5.943921952763129, 0.0, 0.0], [0.0, 5.943921952763129, 0.0], [0.0, 0.0, 5.943921952763129]]u"Å"))
source
FreeBird.AbstractWalkers.split_componentsMethod
split_components(at::AbstractSystem, list_num_par::Vector{Int})

Split the system into components based on the number of particles in each component.

Arguments

  • at::AbstractSystem: The system to split.
  • list_num_par::Vector{Int}: The number of particles in each component.

Returns

  • components: An array of FastSystem objects representing the components of the system.
source
FreeBird.AbstractWalkers.split_components_by_chemical_speciesMethod
split_components_by_chemical_species(at::AbstractSystem)

Split an AbstractSystem into multiple components based on the chemical species.

Arguments

  • at::AbstractSystem: The input AbstractSystem to be split.

Returns

An array of FastSystem objects, each representing a component of the input system.

Example

julia> at = FreeBirdIO.generate_multi_type_random_starting_config(10.0,[2,1,3,4,5,6];particle_types=[:H,:O,:H,:Fe,:Au,:Cl])
FastSystem(Au₅Cl₆Fe₄H₅O, periodic = FFF):
    bounding_box      : [ 5.94392        0        0;
                                0  5.94392        0;
                                0        0  5.94392]u"Å"

        .--------------.  
       /Au      Cl     |  
      / |HAu Fe        |  
     /  |     Cl Cl Cl |  
    *   |Cle           |  
    |   | Cl      H    |  
    |   |      OAuH    |  
    |FeFe-----------H--.  
    |  /          Au  /   
    | /      Au      /    
    |/              /     
    *--------------*      


julia> AbstractWalkers.split_components_by_chemical_species(at)
5-element Vector{FastSystem}:
 FastSystem(H₅, periodic = FFF, bounding_box = [[5.943921952763129, 0.0, 0.0], [0.0, 5.943921952763129, 0.0], [0.0, 0.0, 5.943921952763129]]u"Å")
 FastSystem(O, periodic = FFF, bounding_box = [[5.943921952763129, 0.0, 0.0], [0.0, 5.943921952763129, 0.0], [0.0, 0.0, 5.943921952763129]]u"Å")
 FastSystem(Cl₆, periodic = FFF, bounding_box = [[5.943921952763129, 0.0, 0.0], [0.0, 5.943921952763129, 0.0], [0.0, 0.0, 5.943921952763129]]u"Å")
 FastSystem(Fe₄, periodic = FFF, bounding_box = [[5.943921952763129, 0.0, 0.0], [0.0, 5.943921952763129, 0.0], [0.0, 0.0, 5.943921952763129]]u"Å")
 FastSystem(Au₅, periodic = FFF, bounding_box = [[5.943921952763129, 0.0, 0.0], [0.0, 5.943921952763129, 0.0], [0.0, 0.0, 5.943921952763129]]u"Å")
source
FreeBird.AbstractWalkers.split_into_subarraysMethod
split_into_subarrays(arr::AbstractVector, N::Int)

Split an array into N subarrays of approximately equal size.

Arguments

  • arr::AbstractVector: The array to split.
  • N::Int: The number of subarrays to create.

Returns

  • subarrays::Vector{Vector{eltype(arr)}}: A vector of subarrays.
source
FreeBird.AbstractWalkers.update_walker!Method
update_walker!(walker::AtomWalker, key::Symbol, value)

Update the properties of an AtomWalker object.

A convenient function that updates the value of a specific property of an AtomWalker object.

Arguments

  • walker::AtomWalker: The AtomWalker object to be updated.
  • key::Symbol: The key of the property to be updated.
  • value: The new value of the property.

Returns

  • walker::AtomWalker: The updated AtomWalker object.

Example

update_walker!(walker, :energy, 10.0u"eV")
update_walker!(walker, :iter, 1)
source