API reference

Here is a reference of SuSMoST programming interface generated from source code docstrings.

SuSMoST main module

susmost.cell.universal_cell_size(cell_size)

Converts cell_size into 3x3 matrix from different simpler forms: n -> [n,n,1] [n,m] -> [n,m,1] [n,m,k] -> np.diag([n,m,k]) [[a,b],[c,d]] -> [[a,b,0],[c,d,0],[0,0,1]]

class susmost.LatticeTask(site_state_types, unit_cell, states, edges_dict, IM, INF_E, precission, max_dimensions)[source]

Class for description of a lattice model

Constructor of LatticeTask

Parameters:
site_state_types:

list of types of possible states of lattice sites as SiteStateType objects. Usually correspond to considered adsorption complexes, including an empty site.

unit_cell:

unit cell of the lattice as Cell object

states:

list of possible states of lattice sites as CellState objects. Usually correspond to various orientations of adsorption complexes from site_state_types

property edges_array

Get edges considered in the model as a two-dimensional array of size (edges_count, 3)

property edges_count

Get number of different edges in the lattice, that are considered in the model

get_ads_energy(ac_name)[source]

Returns current energy of adsorption of the specified adsorption complex

Parameters:
ac_name:

name of adsorption complex

Returns:

Сurrent energy of adsorption of the specified adsorption complex

get_property(property_name, ac_name=None)[source]

Returns current value of the property of the specified adsorption complex

Parameters:
property_name:

name of property

ac_name:

name of adsorption complex

Returns:

Current value of the property of the specified adsorption complex

property interactions_array

Get paiwise interactions of the model as a three-dimensional array of size (states_count, edges_count, states_count).

Returns:
result [i,j,k]:

interaction energy between states i-th and k-th separated by j-th edge

set_ads_energy(ac_name, energy)[source]

Sets adsorption energy for adsorption complex by it’s name

Parameters:
ac_name:

name of adsorption complex

energy:

new adsorption energy of the adsorption complex

Returns:

None

set_property(property_name, values, default_value=0)[source]

Sets a property for all adsoprtion complexes

Parameters:
property_name:

name of the property

values:

dictionary with AC names as keys and new property values as values

default_value:

value to be set for ACs not included in :values:, default value - :0:.

Returns:

None

property states_count

Get number of possible cell states

property zero_coverage_state_index

Get index of the state with zero coverage, that is empty state

class susmost.SiteStateType(name, configuration=None, ads_energy=None, xyz_mark=None, **kwargs)[source]

Descriptor of a lattice site state.

Parameters:
name:

name of the state

configuration:

geometry of the state. Internal coordinates <https://en.wikipedia.org/wiki/Z-matrix_(chemistry)> of lattice sites spanned by the state. See functions make_dimer_zmatrix(), make_trimer_zmatrix()`` etc

ads_energy:

energy of the state, can be consdidered as chemical potential of the surface or adsorption energy.

xyz_mark:

list of chemical element symbols, one element for each site of the configuration. Used for visualization purposes only. For example for trajectory generation by mc.run() function.

**kwargs:

other properties of the state. All states of the model must have the same set of properties. New properties must be registered with function register_property(). Property ``coverage``is registered by default.

susmost.joined_cells_lattice_task(lat_task, supercell, cell_self_pbc=[False, False, False], E_limit=None)

Transforms LatticeTask task object by joining lattice unit cells in a supercell. Similar to real-space renormalization approach.

Parameters:
lat_task:

input LatticeTask task object

supercell:

size of a new lattice unit cell in input unit cells. See susmost.cell.universal_cell_size() for possible formats

cell_self_pbc:

periodic boundary conditions along each of three axis in an output LatticeTask

E_limit:

limit of internal energy of the state, by default lat_task.INF_E/2

Returns:

Renormalized LatticeTask object

susmost.load_lattice_task(samples_dir, INF_E=1000000.0, precission=6, max_dimensions=[-1, -1, 0])

Creates LatticeTask object from data in samples_dir, populated with make_ac_samples() function.

Parameters:
samples_dir:

directory populated with make_ac_samples() function and with energies file. energies file must consit of row of the following format: {sample file name} = {energy}, where {sample file name} is a name of every sample_*.xyz file in the directory, {energy} - is an overall energy of the sample, computed with DFT, empirical potentials, fitted to experiments, etc.

INF_E:

energy value that will be considered as infinite in a lattice model

precission:

precission of edges comparison on lattice graph generation

max_dimensions:

list of 3 integers, maximal size of the considered system in unit cells, -1 - mean infinite size. Meant to be used for simulaion of confined surfaces, i.e. steps.

Returns:

LatticeTask object.

susmost.make_ac_samples(samples_dir, ac_fns, r_cut_off, max_dimensions=[-1, -1, 0], sym_prec=0.1, borders=None)

Create samples of adsorption complexes, pairs of adsorption complexes and interaction file.

Parameters:
samples_dir:

Directory to be created and populated with output files (see bellow)

ac_fns:

List of names of xyz-files with sketches of adsorption complexes and a clean adsorbate (surface) slab. The first element of the list should be a clean surface slab. All files must be in Extended XYZ format with the following additional parameters specified:

  1. Lattice - cell vectors in a standart format of nine numbers.

  2. ac_origin - 3-dimensional coordinates of the unit cell origin

  3. size - 2-dimensional integer size of the slab

  4. sym_prec - recommended precission for symmetries search by spglib

Files with adsorption complexes should provide information for distinguishing adsorbate (surface) and adsorbent atoms. There are several ways to do it:

  • With key-value pairs in a comments section of XYZ file:

    1. surface_zlevel - value of z-coordinate that separates adsorbate atoms (above) and adsorbent atoms (bellow). Important: surface_zlevel is meassured with respect to ac_origin. It means that if ac_origin="0 0 3" and surface_zlevel=5, then all atoms with z-coordinate less than 3 + 5 = 8 are considered as surface atoms.

    2. adsorbate_indices - list of integers that enumerate indices of adsorbate atoms.

    3. adsorbate_elements - list of chemical symbols of adsorbate atoms. Can be used only if sets of adsorbate and adsorbent elements do not intersect.

  • With additional atom-wise parameter adsorption - a symbolic parameter, that should be specified for each atom. All atoms with adsorption parameter equal to s are considered as surface atoms. Important: structure of atom-wise parameters in Extended XYZ format is specified as key-value pair Parameters in a a comments section of XYZ file.

r_cut_off:

cutoff distance between surface unit cells - used for computation of sample slab sizes and maximal distance of pairwise interactions.

max_dimensions:

list of 3 integers, maximal size of the considered system in unit cells, -1 - mean infinite size. Meant to be used for simulaion of confined surfaces, i.e. steps.

sym_prec:

default symmetry detection precission. Used if no sym_prec parameter found in the adsorption complex XYZ file. The higher value means more tollerant symmetry detection.

borders:

(‘border’, x, y) or (‘fixed’, x1, y1, x2, y2) Size of samples in unit cells. (‘border’, x, y) gives variable size of samples with fixed distances from adsorption complexes to slab edges. (‘fixed’, x1, y1, x2, y2) gives fixed size samples that can be non-collinear to unit cell - (x1,y1) and (x2,y2) are lattice vectors of samples in unit cell coordinates.

Returns:

None

Output files:
  1. sample_empty_unit_cell.xyz - xyz-file with surface unit cell

  2. sample_ac_{AC-name}_{index}.xyz - xyz-files with adsorption complexes on undisturbed surface slabs of different sizes, here:

    • {AC-name} - name of adsorption complex taken from ac_name property of input xyz-file,

    • {index} - indexes that enumerate different sizes of slabs

  3. sample_{index}.xyz - xyz-files with pairs of adsorption complexes.

  4. interactions - matrices of pairwise interactions - to be used by load_lattice_task() function

  5. marginal_samples - list of pairwise samples that should have zero energy of lateral interactions, because some of their rotations fall out of the cutoff radius

  6. adsorbate_state_{index}.xyz - xyz-files for visualization of states of lattice in a lattice model

susmost.make_tensor(lt, beta, scheme='TRG', dtype=<class 'float'>)
Parameters:
lt:

LatticeTask object

beta:

Inverse temperature

scheme:

‘TRG’ (default) - memory-efficient for solve_TRG; ‘TM’ - memory-efficient for solve_TM; ‘generic’ - memory-expensive but universal (for both solve_TRG and solve_TM), supports ad hoc manybody interactions, sometimes more CPU-efficient than ‘TRG’; ‘raw’ - for debug purposes.

Returns:

Object that represents IRF tensor or corresponding tensor network.

susmost.nearest_int(cc1, cc2, int_en=0.0, int_r=0.5, h=0.5, INF_E=1000000.0, dop_int=0.0, empty_name='Empty')[source]

Calculates the energy of interaction between nearest neighbors. If name of the molecule is ‘Empty’, then there is no any interactions with it.

Parameters:
cc1,cc2:

Coordinates of first and second molecules.

int_en:

Simple interaction energy between elements of the molecule. Default is 0.0.

int_r:

The radius of the interaction. If the distance R between the molecules is h <= R <= int_r, then there is interaction. Default is 0.5.

h:

Radius of a solid sphere. If the distance between the elements of the molecules is smaller or equal, then an infinitely strong repulsion arises. Default is 0.5.

INF_E:

The value of an infinitely strong repulsion. Default is 1e6.

dop_int:

Complex interaction energy between elements of the molecule. Has the form dop_int = {(‘mol_1_name’,’mol_2_name’):energy1,…} Default is 0.0.

empty_name:

Name of the element with properties of empty element. Default is “Empty”.

Returns:

Energy of interaction between molecules.

susmost.normal_lattice_task(site_state_types, unit_cell, unit_cell_atoms, interaction_func, r_cut_off, INF_E, precission=6, max_dimensions=[-1, -1, 0])

Creates LatticeTask object for adstract lattice model.

Parameters:
site_state_types:

list of SiteStateType objects

unit_cell:

list of unit cell basis vectors, specified as 3x3 array

unit_cell_atoms:

list of 3-d coordinates of atoms in a unit cell

interaction_func:

function that returns pair interaction for a pair of site states specified in carthesian coordinates

r_cut_off:

cutoff distance between surface unit cells - used for computation of sample slab sizes and maximal distance of pairwise interactions

INF_E:

energy value that will be considered as infinite in a lattice model

precission:

precission of edges comparison on lattice graph generation

max_dimensions:

list of 3 integers, maximal size of the considered system in unit cells, -1 - mean infinite size. Meant to be used for simulaion of confined surfaces, i.e. steps.

Returns:

LatticeTask object

susmost.phash(x)[source]

Persistent hash

Parameters:
x:

any object with defined repr(x)

Returns:

hashlib.md5(repr(x).encode('utf-8')).hexdigest()

susmost.register_property(name, cmp_func='==', merge_func=None, join_func=None)[source]

Register new state property

Parameters:
name:

name of the property

susmost.solve_TRG(W, max_n=64, max_value=1e-08, atol_lnZ=1e-09, contractions_count=30, cuda=False, use_rndsvd=False, MIZK_nsqr_threshold='auto-memory', MIZK_r='auto-memory', rndsvd_oversampling=2.0, rndsvd_iterations=4, memory_limit='auto', verbose=False)
Parameters:
W:

tensor representation

max_n:

number of singular values to keep on reduced singular value decompositions (SVD’s)

max_value:

maximal singular value that can be thrown out on reduced SVD’s

atol_lnZ:

convergence tollerance of contraction iterations

contractions_count:

maximum number of TRG contraction iterations

cuda:

use CUDA, should be used only with use_rndsvd=True

use_rndsvd:

use randomized SVD with Levin-Nave, otherwise :func:scipy.linalg.svd is used.

MIZK_nsqr_threshold:

n^2 threshold for usage of four 3-leg tensors decomposition (4x3) instead of a single 4-leg tensor. See 10.1103/PhysRevE.97.033310 by Morita-Igarashi-Zhao-Kawashima (MIZK) for details. If W.shape[0]*W.shape[1] is greater than MIZK_nsqr_threshold then the MIZK decopmosition is used, else a single 4-leg tensor as by Levin-Nave will be created. Should not affect final result if use_rndsvd=False. By default calculated automatically from memory_limit.

MIZK_r:

MIZK block size for contraction of 4x3 tensor decomposition. Smaller values slow down computation, but reduce memory usage and vice versa. By default calculated automatically from memory_limit.

rndsvd_oversampling:

oversampling for randomized SVD

rndsvd_iterations:

number of power iterations for randomized SVD

memory_limit:

maximum memory to be used. By default is psutil.virtual_memory().available - 1 GB. For CUDA should be specified.

verbose:

print intermediate ln(Z) and trace values, timings, used algorithms, tensor dimension sizes etc.

Returns:

Natural logarithm of partition function per lattice cell ln(Z).

SuSMoST Monte Carlo module

class susmost.mc.Lattice(lattice_task=None, lattice_size=None, metropolis=None)[source]

python -c ‘import numpy as np; from susmost import mc, load_lattice_task; lt = load_lattice_task(“tmp_00”); lt.set_ads_energy(“atop-triang”, -1.); m = mc.make_metropolis(lt, 10, [10.], 1.); mc.run(m, 10); print(m.cells, lt.states[0].state_type.__dict__); l = mc.Lattice(metropolis=m); print(l.calc_fp()); print (l.calc_lateral_energy_per_site()), print (m.curE.internal); print (lt.site_state_types); print (“PPP”,l.calc_property_per_site(“ads_energy”)) ‘

susmost.mc.make_metropolis(lat_task, supercell, T, k_B, param_names=None, precoverage=None, seed=None, pbc_flags=[True, True, False])[source]

Create Metropolis object for specified LatticeTask and ReGraph objects

Parameters:
lat_task:

LatticeTask object

supercell:

size of the lattice. See susmost.cell.universal_cell_size() for possible formats

T:

list of temperatures in Kalvin units for replica exchange (aka parallel tempering), number of temperatures must be equal to number of MPI processes

k_B:

Boltzmann constant, used as energy scale parameter

param_names:

list of names of states properties that will be sampled in Metropolis iterations as simulaion parameters; default value - [‘coverage’]

precoverage:

index of the state that will be used for initialization of the lattice, by default the state with zero coverage parameter is used

seed:

random number generator for Metropolis simulaion

pbc_flags:

list or tuple of three boolean flags, that enables/disables Periodic Boundary Conditions along each axis; default value - [True, True, False]

Returns:

Metropolis object

class susmost.mc.param_limit((object)arg1, (object)param_index, (object)lower_bound, (object)upper_bound, (object)enabled)

Descriptor of a parameter limit

C++ signature :

void __init__(_object*,int,double,double,bool)

check((param_limit)arg1, (object)arg2) bool :
C++ signature :

bool check(param_limit_t {lvalue},double)

susmost.mc.run(m, log_periods_cnt, log_period_steps=None, log_callback=None, params_period_steps=None, relaxation_steps=None, bcast_full_params_log=True, traj_fns='auto', recalc_curs=True)[source]

Perform Metropolis iterations

Parameters:
m:

Metropolis object

log_periods_cnt:

number of log periods to perform

log_period_steps:

number of Metropolis iterations in a single log period. By default m.cells_count value is used so log period is equivalent to conventional Monte Carlo step notion

log_callback:

function to be called after each log period, log_callback function must take m as a parameter and can return any object as result. List of all objects returned from log_callback will be returned as a result of overall run() call. By default is None, so nothing is logged.

params_period_steps:

period for sampling of energy and state properties (see notes on make_metropolis() function). Samples are stored in full_params_log attribute of m object. Meassured in Metropolis iterations. Default value (log_periods_cnt * log_period_steps)//1000 + 1 so up to 1000 entries will be added to full_params_log.

relaxation_steps:

number of Metropolis iterations to be skipped before logging and parameters sampling. Default value log_periods_cnt * log_period_steps

bcast_full_params_log:

set to False if parameters samples should not be broadcasted to all MPI processes. By default True, so every processes receives complete sample of parameters.

traj_fns:

list of XYZ file names to be used to save snapshots of the lattice each log_period_steps iteration. traj_fns[i] is used to save samples from i-th replica (see notes on make_metropolis() function). By default “traj.{temperature}.xyz” is used for each temperature

recalc_curs:

Call recalc_curs() method before running Metropolis simulaions. recalc_curs() recomputes total energy and additive parameter values. It is computationally expensive, but necessary if cells property was modified outside of Metropolis iterations, for example due to coverage initialization.

Number of steps to execute: relaxation_steps + log_periods_cnt * log_period_steps. Number of full_params_log entries: (log_periods_cnt * log_period_steps) / params_period_steps. Number of samples in trajectory files: log_periods_cnt.

Results:

List of objects returned by log_callback function calls

After the function call m.full_params_log contains three-dimensional array with samples of simulation parameters. Parameter values are averaged over lattice sites.

m.full_params_log[i,j,k] is the value of the k-th parameter on j-th time period from i-th replica. Parameter index k correspond to the k-th item of param_names argument of make_metropolis() function. In addition to state parameters, several system-wide parameters are included into m.full_params_log. In particular:

  • m.full_params_log[i,j,-1] index of the temperature of this sample

  • m.full_params_log[i,j,-2] total energy of the system per one lattice cell

  • m.full_params_log[i,j,-3] lateral interactions energy (total energy minus adsoprtion energy) per lattice cell

  • m.full_params_log[i,j,-4] acceptance rate of Metropolis iterations

  • m.full_params_log[i,j,-5] KMC r0 parameter, i.e. maximal step rate, should be constant through correct simulaion

  • m.full_params_log[i,j,-6] free adsorption energy per lattice cell

  • m.full_params_log[i,j,-7] ground state (minimal) adsorption energy per lattice cell

  • m.full_params_log[i,j,-8] entropy of adsoprtion complexes per lattice cell

  • m.full_params_log[i,j,-9] mean internal energy (enthalpy) of adsoprtion complexes per lattice cell

  • m.full_params_log[i,j,:-9] samples of states parameters - see param_names argument of make_metropolis() function.