lattice_mc package

Submodules

lattice_mc.atom module

class lattice_mc.atom.Atom(initial_site)[source]

Bases: object

Atoms are distinguishable particles, each occupying a specific lattice site.

atom_number = 0
dr_squared()[source]

\(|dr|^2\), where \(dr\) is the total displacement vector for this Atom.

Parameters:None
Returns:\(|dr|^2\).
Return type:dr_squared (float)
reset()[source]

Reinitialise the stored displacements, number of hops, and list of sites visited for this Atom.

Parameters:None
Returns:None
site

Get or set self.site for this Atom.

lattice_mc.cluster module

class lattice_mc.cluster.Cluster(sites)[source]

Bases: object

Clusters are sets of sites.

is_neighbouring(other_cluster)[source]

logical check whether the neighbour list for cluster A includes any sites in cluster B

Parameters:other_cluster (Cluster) – the other cluster we are testing for neighbour connectivity
Returns:True if the other cluster neighbours this cluster.
Return type:(Bool)
is_periodically_contiguous()[source]

logical check whether a cluster connects with itself across the simulation periodic boundary conditions.

Parameters:none
Returns
( Bool, Bool, Bool ): Contiguity along the x, y, and z coordinate axes
merge(other_cluster)[source]

Combine two clusters into a single cluster.

Parameters:other_cluster (Cluster) – The second cluster to combine.
Returns:The combination of both clusters.
Return type:(Cluster)
remove_sites_from_neighbours(remove_labels)[source]

Removes sites from the set of neighbouring sites if these have labels in remove_labels.

Parameters:Remove_labels (List) or (Str) – List of Site labels to be removed from the cluster neighbour set.
Returns:None
sites_at_edges()[source]

Finds the six sites with the maximum and minimum coordinates along x, y, and z.

Parameters:None
Returns:In the order [ +x, -x, +y, -y, +z, -z ]
Return type:(List(List))
size()[source]

Number of sites in this cluster.

Parameters:None
Returns:The number of sites in this cluster.
Return type:(Int)

lattice_mc.global_vars module

lattice_mc.init_lattice module

lattice_mc.init_lattice.cubic_lattice(a, b, c, spacing)[source]

Generate a cubic lattice.

Parameters:
  • a (Int) – Number of lattice repeat units along x.
  • b (Int) – Number of lattice repeat units along y.
  • c (Int) – Number of lattice repeat units along z.
  • spacing (Float) – Distance between lattice sites.
Returns:

The new lattice

Return type:

(Lattice)

lattice_mc.init_lattice.honeycomb_lattice(a, b, spacing, alternating_sites=False)[source]

Generate a honeycomb lattice.

Parameters:
  • a (Int) – Number of lattice repeat units along x.
  • b (Int) – Number of lattice repeat units along y.
  • spacing (Float) – Distance between lattice sites.
  • alternating_sites (Bool, optional) – Label alternating sites with ‘A’ and ‘B’. Defaults to False.
Returns:

The new lattice

Return type:

(Lattice)

Notes

The returned lattice is 3D periodic, but all sites and edges lie in the xy plane.

lattice_mc.init_lattice.lattice_from_sites_file(site_file, cell_lengths)[source]

Generate a lattice from a sites file.

Parameters:
  • site_file (Str) – Filename for the file containing the site information.
  • cell_lengths (List(Float,Float,Float)) – A list containing the [ x, y, z ] cell lengths.
Returns:

The new lattice

Return type:

(Lattice)

Notes

The site information file format is:
<number_of_sites> (Int).
Followed by blocks of data separated by double linebreaks; one block per site.
site: <site number> (Int).
centre: <x> <y> <z> (Float,Float,Float).
neighbours: <list of site numbers of neighbouring sites> (List[Int]).
label: <site group labal> (Str).
energy: <site occupation energy> (Float).
The energy is optional, and will be set to 0.0 if not included.
Line order within each block is not meaningful.
British and American English spellings for centre|center and neighbour|neighbor are accepted.
An example file can be found in the examples directory.
lattice_mc.init_lattice.square_lattice(a, b, spacing)[source]

Generate a square lattice.

Parameters:
  • a (Int) – Number of lattice repeat units along x.
  • b (Int) – Number of lattice repeat units along y.
  • spacing (Float) – Distance between lattice sites.
Returns:

The new lattice

Return type:

(Lattice)

Notes

The returned lattice is 3D periodic, but all sites and edges lie in the xy plane.

lattice_mc.jump module

class lattice_mc.jump.Jump(initial_site, final_site, nearest_neighbour_energy=False, coordination_number_energy=False, jump_lookup_table=None)[source]

Bases: object

boltzmann_factor()[source]

Boltzmann probability factor for accepting this jump, following the Metropolis algorithm.

Parameters:None
Returns:Metropolis relative probability of accepting this jump.
Return type:(Float)
coordination_number_delta_E()[source]

Coordination-number dependent energy conrtibution to the change in system energy if this jump were accepted.

Parameters:None
Returns:delta E (coordination-number)
Return type:(Float)
delta_E()[source]

The change in system energy if this jump were accepted.

Parameters:None
Returns:delta E
Return type:(Float)
dr(cell_lengths)[source]

Particle displacement vector for this jump

Parameters:cell_lengths (np.array(x,y,z)) – Cell lengths for the orthogonal simulation cell.
Returns
(np.array(x,y,z)): dr
nearest_neighbour_delta_E()[source]

Nearest-neighbour interaction contribution to the change in system energy if this jump were accepted.

Parameters:None
Returns:delta E (nearest-neighbour)
Return type:(Float)
rate()[source]

Average rate for this jump. Calculated as ( v_0 * P_jump ).

Parameters:None
Returns:The average rate for this jump.
Return type:(Float)
relative_probability
relative_probability_from_lookup_table(jump_lookup_table)[source]

Relative probability of accepting this jump from a lookup-table.

Parameters:jump_lookup_table (LookupTable) – the lookup table to be used for this jump.
Returns:relative probability of accepting this jump.
Return type:(Float)

lattice_mc.lattice module

class lattice_mc.lattice.Lattice(sites, cell_lengths)[source]

Bases: object

Lattice class

connected_site_pairs()[source]

Returns a dictionary of all connections between pair of sites (by site label). e.g. for a linear lattice A-B-C will return:

{ 'A' : [ 'B' ], 'B' : [ 'A', 'C' ], 'C' : [ 'B' ] }
Parameters:None
Returns:A dictionary of neighbouring site types in the lattice.
Return type:site_connections (Dict{Str List[Str]})
connected_sites(site_labels=None)[source]

Searches the lattice to find sets of sites that are contiguously neighbouring. Mutually exclusive sets of contiguous sites are returned as Cluster objects.

Parameters:( (site_labels) –

obj:(List(Str)|Set(Str)|Str), optional): Labels for sites to be considered in the search. This can be a list:

[ 'A', 'B' ]

a set:

( 'A', 'B' )

or a string:

'A'.
Returns:List of Cluster objects for groups of contiguous sites.
Return type:(List(Cluster))
detached_sites(site_labels=None)[source]

Returns all sites in the lattice (optionally from the set of sites with specific labels) that are not part of a percolating network. This is determined from clusters of connected sites that do not wrap round to themselves through a periodic boundary.

Parameters:site_labels (String or List(String)) – Lables of sites to be considered.
Returns:List of sites not in a periodic percolating network.
Return type:(List(Site))
enforce_periodic_boundary_conditions()[source]

Ensure that all lattice sites are within the central periodic image of the simulation cell. Sites that are outside the central simulation cell are mapped back into this cell.

Parameters:None
Returns:None
initialise_site_lookup_table()[source]

Create a lookup table allowing sites in this lattice to be queried using self.site_lookup[n] where n is the identifying site numbe.

Parameters:None
Returns:None
is_blocked()[source]

Check whether there are any possible jumps.

Parameters:None
Returns:True if there are no possible jumps. Otherwise returns False.
Return type:(Bool)
jump()[source]

Select a jump at random from all potential jumps, then update the lattice state.

Parameters:None
Returns:None
max_site_coordination_numbers()[source]

Returns a dictionary of the maximum coordination number for each site label. e.g.:

{ 'A' : 4, 'B' : 4 }
Parameters:none
Returns:
Int)): dictionary of maxmimum coordination
number for each site label.
Return type:max_coordination_numbers (Dict(Str
occupied_site_numbers()[source]

List of site id numbers for all sites that are occupied.

Parameters:None
Returns:List of site id numbers for occupied sites.
Return type:List(Int)
occupied_sites()[source]

The set of sites occupied by atoms.

Parameters:None
Returns:List of sites that are occupied.
Return type:List(Site)
populate_sites(number_of_atoms, selected_sites=None)[source]

Populate the lattice sites with a specific number of atoms.

Parameters:
  • number_of_atoms (Int) – The number of atoms to populate the lattice sites with.
  • ( (selected_sites) – obj:List, optional): List of site labels if only some sites are to be occupied. Defaults to None.
Returns:

None

potential_jumps()[source]

All nearest-neighbour jumps not blocked by volume exclusion (i.e. from occupied to neighbouring unoccupied sites).

Parameters:None
Returns:List of possible jumps.
Return type:(List(Jump))
reset()[source]

Reset all time-dependent counters for this lattice and its constituent sites

Parameters:None
Returns:None
select_sites(site_labels)[source]

Selects sites in the lattice with specified labels.

Parameters:site_labels (List(Str)|Set(Str)|Str) – Labels of sites to select. This can be a List [ ‘A’, ‘B’ ], a Set ( ‘A’, ‘B’ ), or a String ‘A’.
Returns:List of sites with labels given by site_labels.
Return type:(List(Site))
set_cn_energies(cn_energies)[source]

Set the coordination number dependent energies for this lattice.

Parameters:(Dict(Str (cn_energies) –

Dict(Int:Float))): Dictionary of dictionaries specifying the coordination number dependent energies for each site type. e.g.:

{ 'A' : { 0 : 0.0, 1 : 1.0, 2 : 2.0 }, 'B' : { 0 : 0.0, 1 : 2.0 } }
Returns:None
set_nn_energy(delta_E)[source]

Set the lattice nearest-neighbour energy.

Parameters:delta_E (Float) – The nearest-neighbour energy E_nn.
Returns:None
set_site_energies(energies)[source]

Set the energies for every site in the lattice according to the site labels.

Parameters:(Dict(Str (energies) –

Float): Dictionary of energies for each site label, e.g.:

{ 'A' : 1.0, 'B', 0.0 }
Returns:None
site_coordination_numbers()[source]

Returns a dictionary of the coordination numbers for each site label. e.g.:

{ 'A' : { 4 }, 'B' : { 2, 4 } }
Parameters:none
Returns:
Set(Int))): dictionary of coordination
numbers for each site label.
Return type:coordination_numbers (Dict(Str
site_occupation_statistics()[source]

Average site occupation for each site type

Parameters:None
Returns:Float)): Dictionary of occupation statistics, e.g.:
{ 'A' : 2.5, 'B' : 25.3 }
Return type:(Dict(Str
site_specific_coordination_numbers()[source]

Returns a dictionary of coordination numbers for each site type.

Parameters:None
Returns:List(Int))) : Dictionary of coordination numbers for each site type, e.g.:
{ 'A' : [ 2, 4 ], 'B' : [ 2 ] }
Return type:(Dict(Str
site_with_id(number)[source]

Select the site with a specific id number.

Parameters:number (Int) – The identifying number for a specific site.
Returns:The site with id number equal to number
Return type:(Site)
transmute_sites(old_site_label, new_site_label, n_sites_to_change)[source]

Selects a random subset of sites with a specific label and gives them a different label.

Parameters:
  • old_site_label (String or List(String)) – Site label(s) of the sites to be modified..
  • new_site_label (String) – Site label to be applied to the modified sites.
  • n_sites_to_change (Int) – Number of sites to modify.
Returns:

None

update(jump)[source]

Update the lattice state by accepting a specific jump

Parameters:jump (Jump) – The jump that has been accepted.
Returns:None.
update_site_occupation_times(delta_t)[source]

Increase the time occupied for all occupied sites by delta t

Parameters:delta_t (Float) – Timestep.
Returns:None
vacant_site_numbers()[source]

List of site id numbers for all sites that are vacant.

Parameters:None
Returns:List of site id numbers for vacant sites.
Return type:List(Int)
vacant_sites()[source]

The set of sites not occupied by atoms.

Parameters:None
Returns:List of sites that are vacant.
Return type:List(Site)

lattice_mc.lattice_site module

class lattice_mc.lattice_site.Site(number, coordinates, neighbours, energy, label, cn_energies=None)[source]

Bases: object

Site class

cn_occupation_energy(delta_occupation=None)[source]

The coordination-number dependent energy for this site.

Parameters:( (delta_occupation) – obj:Dict(Str:Int), optional): A dictionary of a change in (site-type specific) coordination number, e.g. { ‘A’ : 1, ‘B’ : -1 }. If this is not None, the coordination-number dependent energy is calculated including these changes in neighbour-site occupations. Defaults to None
Returns:The coordination-number dependent energy for this site.
Return type:(Float)
index = 0
nn_occupation()[source]

The number of occupied nearest-neighbour sites.

Parameters:None
Returns:The number of occupied nearest-neighbour sites.
Return type:(Int)
set_cn_occupation_energies(cn_energies)[source]

Set the coordination-number dependent energies for this site.

Parameters:(Dict(Int (cn_energies) – Float)): Dictionary of coordination number dependent site energies, e.g. { 0 : 0.0, 1 : 0.5 }.
Returns:None
site_specific_neighbours()[source]

Returns the number of neighbouring sites, classified by site type.

Parameters:None
Returns:Int)): Dictionary of neighboring sites, classified by site label, e.g. { ‘A’ : 1, ‘B’ : 1 }.
Return type:(Dict(Str
site_specific_nn_occupation()[source]

Returns the number of occupied nearest neighbour sites, classified by site type.

Parameters:None
Returns:Int)): Dictionary of nearest-neighbour occupied site numbers, classified by site label, e.g. { ‘A’ : 2, ‘B’ : 1 }.
Return type:(Dict(Str

lattice_mc.lookup_table module

class lattice_mc.lookup_table.LookupTable(lattice, hamiltonian)[source]

Bases: object

LookupTable class

generate_nearest_neighbour_lookup_table()[source]

Construct a look-up table of relative jump probabilities for a nearest-neighbour interaction Hamiltonian.

Parameters:None.
Returns:None.
relative_probability(l1, l2, c1, c2)[source]

The relative probability for a jump between two sites with specific site types and coordination numbers.

Parameters:
  • l1 (Str) – Site label for the initial site.
  • l2 (Str) – Site label for the final site.
  • c1 (Int) – Coordination number for the initial site.
  • c2 (Int) – Coordination number for the final site.
Returns:

The relative probability of this jump occurring.

Return type:

(Float)

lattice_mc.lookup_table.metropolis(delta_E)[source]

Boltzmann probability factor for an event with an energy change delta_E, following the Metropolis algorithm.

Parameters:delta_E (Float) – The change in energy.
Returns:Metropolis relative probability for this event.
Return type:(Float)

lattice_mc.options module

class lattice_mc.options.Options[source]

Bases: object

Object for storing options for setting up and running a simulation.

read_lattice_from_file(filename)[source]

Set the filename for the sites file used to define the simulation lattice.

Parameters:filename (Str) – The filename for the sites file.
Returns:None
set_cn_energies(cn_energies)[source]

Set the coordination-number dependent energies.

Parameters:(Dict(Str (cn_energies) – Dict(Int:Float))): Dict of Dicts specifying the coordination-number dependent energies for each site type. e.g. { ‘A’ : { 0 : 0.0, 1 : 1.0 }, ‘B’ : { 0 : 0.0, 1 : 2.0 } }
Returns:None
set_cn_energy_scaling(cn_energy_scaling)[source]

Set a scaling factor for the coordination-number dependent energies.

Parameters:cn_energy_scaling (Float) – Coordination-number dependent energy scaling.
Returns:None
set_lattice_cell_lengths(cell_lengths)[source]

Set the lattice cell lengths for a simulation.

Parameters:cell_lengths (np.array(x,y,z)) – Lattice cell lengths.
Returns:None
set_nn_energy_scaling(energy_scaling)[source]

Set a scaling factor for the nearest-neighbour interaction energy.

Parameters:energy_scaling (Float) – Nearest-neighbour energy scaling.
Returns:None
set_number_of_atoms(number_of_atoms)[source]

Set the number of atoms present in the simulation.

Parameters:number_of_atoms (int) – The number of atoms.
Returns:None
set_number_of_equilibration_jumps(number_of_equilibration_jumps)[source]

Set the number of equilibration jumps for a simulation.

Parameters:number_of_equilibration_jumps (Int) – The number of equiibration jumps.
Returns:None
set_number_of_jumps(number_of_jumps)[source]

Set the on-site energies for each site type.

Parameters:(Dict(Str (site_energies) – Float)): On-site energies for each site type label. e.g. { ‘A’ : 1.0, ‘B’, -1.0 }
Returns:None
set_site_energies(site_energies)[source]

Set the on-site energies for each site type.

Parameters:(Dict(Str (site_energies) – Float)): On-site energies for each site type label. e.g. { ‘A’ : 1.0, ‘B’, -1.0 }
Returns:None

lattice_mc.simulation module

class lattice_mc.simulation.Simulation[source]

Bases: object

Simulation class

average_site_occupations

Average site occupation numbers.

Parameters:None
Returns:
Float)): Average site occupation numbers for each site label,
e.g. { ‘A’ : 12.4, ‘B’ : 231.2 }
Return type:(Dict(Str
collective_correlation

Returns the collective correlation factor, f_I

Parameters:None
Returns:The collective correlation factor, f_I.
Return type:(Float)
collective_diffusion_coefficient

Returns the collective or “jump” diffusion coefficient, D_J.

Parameters:None
Returns:The collective diffusion coefficient, D_J.
Return type:(Float)
collective_diffusion_coefficient_per_atom

The collective diffusion coefficient per atom. D_J / n_atoms.

Parameters:None
Returns:The collective diffusion coefficient per atom, D_J / n_atoms.
Return type:(Float)
define_lattice_from_file(filename, cell_lengths)[source]
Set up the simulation lattice from a file containing site data.
Uses init_lattice.lattice_from_sites_file, which defines the site file spec.
Parameters:
  • filename (Str) – sites file filename.
  • cell_lengths (List(x,y,z)) – cell lengths for the simulation cell.
Returns:

None

is_initialised()[source]

Check whether the simulation has been initialised.

Parameters:None
Returns:None
old_collective_correlation

Returns the collective correlation factor, f_I

Parameters:None
Returns:The collective correlation factor, f_I.
Return type:(Float)

Notes

This function assumes that the jump distance between sites has been normalised to a=1. If the jumps distance is not equal to 1 then the value returned by this function should be divided by a^2. Even better, use self.collective_correlation

old_tracer_correlation

Deprecated tracer correlation factor for this simulation.

Parameters:None
Returns:The tracer correlation factor, f.
Return type:(Float)

Notes

This function assumes that the jump distance between sites has been normalised to a=1. If the jump distance is not equal to 1 then the value returned by this function should be divided by a^2. Even better, use self.tracer_correlation.

reset()[source]

Reset all counters for this simulation.

Parameters:None
Returns:None
run(for_time=None)[source]

Run the simulation.

Parameters:( (for_time) – obj:Float, optional): If for_time is set, then run the simulation until a set amount of time has passed. Otherwise, run the simulation for a set number of jumps. Defaults to None.
Returns:None
set_cn_energies(cn_energies)[source]

Set the coordination-number dependent energies for this simulation.

Parameters:(Dict(Str (cn_energies) – Dict(Int:Float))): Dict of Dicts specficying the coordination-number dependent energies. e.g. { ‘A’ : { 0 : 0.0, 1: 0.5 }, ‘B’ : { 0 : -0.5, 1 : -2.0 } }
Returns:None
set_nn_energy(nn_energy)[source]

Set the nearest-neighbour energy for this simulation.

Parameters:nn_energy (Float) – nearest-neigbour energy.
Returns:None
set_number_of_atoms(n, selected_sites=None)[source]

Set the number of atoms for the simulation, and populate the simulation lattice.

Parameters:
  • n (Int) – Number of atoms for this simulation.
  • ( (selected_sites) – obj:(List|Set|String), optional): Selects a subset of site types to be populated with atoms. Defaults to None.
Returns:

None

set_number_of_equilibration_jumps(n)[source]

Set the number of equilibration jumps for this simulation.

Parameters:n (Int) – number of equilibration jumps
Returns:None
set_number_of_jumps(n)[source]

Set the number of jumps for this simulation.

Parameters:n (Int) – number of jumps
Returns:None
set_site_energies(site_energies)[source]

Set the on-site energies for this simulation.

Parameters:(Dict(Str (site_energies) – Float)): Dict of energies for each site type. e.g. { ‘A’ : 0.0,. ‘B’ : 1.0 }
Returns:None
setup_lookup_table(hamiltonian='nearest-neighbour')[source]

Create a jump-probability look-up table corresponding to the appropriate Hamiltonian.

Parameters:hamiltonian (Str, optional) – String specifying the simulation Hamiltonian. valid values are ‘nearest-neighbour’ (default) and ‘coordination_number’.
Returns:None
tracer_correlation

Tracer correlation factor, f.

Parameters:None
Returns:The tracer correlation factor, f.
Return type:(Float)
tracer_diffusion_coefficient

Tracer diffusion coefficient, D*.

Parameters:None
Returns:The tracer diffusion coefficient, D*.
Return type:(Float)

lattice_mc.species module

class lattice_mc.species.Species(atoms)[source]

Bases: object

Species class.

Contains methods that operate on sets of Atom objects

collective_correlation()[source]

Collective correlation factor, f_I, for these atoms.

Parameters:None
Returns:The collective correlation factor, f_I, for these atoms.
Return type:(Float)
collective_dr_squared()[source]

Squared sum of total displacements for these atoms.

Parameters:None
Returns:The square of the summed total displacements for these atoms.
Return type:(Float)
occupations(site_label)[source]

Number of these atoms occupying a specific site type.

Parameters:site_label (Str) – Label for the site type being considered.
Returns:Number of atoms occupying sites of type site_label.
Return type:(Int)
sites_occupied()[source]

Returns a list of sites occupied by these atoms.

Parameters:None
Returns:List of sites occupied by these atoms.
Return type:(List)
sum_dr_squared()[source]

Sum of squared total displacements for these atoms.

Parameters:None
Returns:The sum of squared total displacements for these atoms.
Return type:(Float)
summed_dr2()[source]

Sum of squared individual displacements for these atoms.

Parameters:None
Returns:The sum of squared individual displacements for these atoms.
Return type:(Float)
tracer_correlation()[source]

Tracer correlation factor, f, for these atoms.

Parameters:None
Returns:The tracer correlation factor, f, for these atoms.
Return type:(Float)

lattice_mc.transitions module

class lattice_mc.transitions.Transitions(jumps)[source]

Bases: object

Transitions class

Contains methods that operate on sets of Jumps.

cumulative_probabilities()[source]

Cumulative sum of the relative probabilities for all possible jumps.

Parameters:None
Returns:Cumulative sum of relative jump probabilities.
Return type:(np.array)
random()[source]

Select a jump at random with appropriate relative probabilities.

Parameters:None
Returns:The randomly selected Jump.
Return type:(Jump)
time_to_jump()[source]

The timestep until the next jump.

Parameters:None
Returns:The timestep until the next jump.
Return type:(Float)

Module contents