Reference

Command Line Interface (CLI)

usage: kmos3 [-h] [--version]
             {benchmark,build,export,settings-export,import,rebuild,shell,run,view,xml}
             ...

Positional Arguments

parser

Possible choices: benchmark, build, export, settings-export, import, rebuild, shell, run, view, xml

Named Arguments

--version

show program’s version number and exit

Sub-commands

benchmark

Run several kMC steps on the model in the current directory and report the runtime.

kmos3 benchmark [-h] [--nsteps NSTEPS]
Named Arguments
--nsteps

number of steps to run the model (default: 1e6)

build

Build ‘kmc_model.*.so’ from the *.f90 files in the current directory.

kmos3 build [-h] [-d] [-b {local_smart,lat_int,otf}] [-n] [-f FCOMPILER]
            [-l VARIABLE_LENGTH] [-t] [--acf] [-o]
Named Arguments
-d, --debug

turn on assertion statements in the F90 code (default: False)

-b, --backend

Possible choices: local_smart, lat_int, otf

choose the backend (default: ‘local_smart’) Note: ‘lat_int’ and ‘otf’ are EXPERIMENTAL

-n, --no-compiler-optimization

do not send optimizing flags to the compiler (default: False)

-f, --fcompiler

set the fortran compiler (by default it is determined from the environment variable ‘F2PY_FCOMPILER’ or, if not set, determined by numpy)

-l, --variable-length

set the maximum length for process names within the fortran source (default: 95)

-t, --temp-acc

use the temporal acceleration scheme (default: False) Builds the modules base_acc.f90, lattice_acc.mpy, proclist_constants_acc.mpy and proclist_generic_subroutines_acc.mpy.

--acf

build the modules base_acf.f90 and proclist_acf.f90 (default: False) The two modules contain functions to calculate the autocorrelation function (ACF) and mean squared displacement (MSD).

-o, --overwrite

overwrite ‘kmc_settings.py’ and ‘kmc_model.*.so’ without confirmation (default: False)

export

Take a kmos3 xml-file and export all generated source code to the export-path. There try to build the kmc_model.*.so.

kmos3 export [-h] [-d] [-b {local_smart,lat_int,otf}] [-n] [-f FCOMPILER]
             [-l VARIABLE_LENGTH] [-t] [--acf] [-o] [--outdir OUTDIR] [-s]
             infile
Positional Arguments
infile

XML or INI model file

Named Arguments
-d, --debug

turn on assertion statements in the F90 code (default: False)

-b, --backend

Possible choices: local_smart, lat_int, otf

choose the backend (default: ‘local_smart’) Note: ‘lat_int’ and ‘otf’ are EXPERIMENTAL

-n, --no-compiler-optimization

do not send optimizing flags to the compiler (default: False)

-f, --fcompiler

set the fortran compiler (by default it is determined from the environment variable ‘F2PY_FCOMPILER’ or, if not set, determined by numpy)

-l, --variable-length

set the maximum length for process names within the fortran source (default: 95)

-t, --temp-acc

use the temporal acceleration scheme (default: False) Builds the modules base_acc.f90, lattice_acc.mpy, proclist_constants_acc.mpy and proclist_generic_subroutines_acc.mpy.

--acf

build the modules base_acf.f90 and proclist_acf.f90 (default: False) The two modules contain functions to calculate the autocorrelation function (ACF) and mean squared displacement (MSD).

-o, --overwrite

overwrite ‘kmc_settings.py’ and ‘kmc_model.*.so’ without confirmation (default: False)

--outdir

directory to build the kMC model

-s, --source-only

export source only and don’t build binary (default: False)

settings-export

Take a kmos3 xml-file and export ‘kmc_settings.py’ to the export-path.

kmos3 settings-export [-h] [--outdir OUTDIR] [-b {local_smart,lat_int,otf}]
                      infile
Positional Arguments
infile

XML or INI model file

Named Arguments
--outdir

directory to build the kMC model

-b, --backend

Possible choices: local_smart, lat_int, otf

choose the backend (default: ‘local_smart’) Note: ‘lat_int’ and ‘otf’ are EXPERIMENTAL

import

Take a kmos3 xml-file and open an ipython shell with the project_tree imported as kmc_model.

kmos3 import [-h] infile
Positional Arguments
infile

XML or INI model file

rebuild

Export code and rebuild the binary module from the ‘kmc_model.*.so’ in the current directory.

kmos3 rebuild [-h] [-d] [-n] [-f FCOMPILER]
Named Arguments
-d, --debug

turn on assertion statements in the F90 code (default: False)

-n, --no-compiler-optimization

do not send optimizing flags to the compiler (default: False)

-f, --fcompiler

set the fortran compiler (by default it is determined from the environment variable ‘F2PY_FCOMPILER’ or, if not set, determined by numpy)

shell (run)

Open an interactive shell and load the model from the curent directory as KMC_Model in it

kmos3 shell [-h]

view

Take a ‘kmc_model.*.so’ and ‘kmc_settings.py’ in the same directory and start to simulate the model visually.

kmos3 view [-h] [-v STEPS_PER_FRAME]
Named Arguments
-v, --steps-per-frame

number of steps per displayed frame (default: 50000)

xml

Print xml representation of model to stdout

kmos3 xml [-h]

Data Types

kmos3.types

Holds all the data models used in kmos3.

kmos3.types.Action

alias of ConditionAction

class kmos3.types.Bystander(**kwargs)
_shorthand()
attributes = ['coord', 'allowed_species', 'flag']
kmos3.types.Condition

alias of ConditionAction

class kmos3.types.ConditionAction(**kwargs)

Represents either a condition or an action. Since both have the same attributes we use the same class here, and just store them in different lists, depending on its role. For better readability one can also use Condition or Action which are just aliases.

Parameters:
_shorthand()
attributes = ['species', 'coord', 'implicit']
class kmos3.types.Coord(**kwargs)

Class that holds exactly one coordinate as used in the description of a process. The distinction between a Coord and a Site may seem superfluous but it is made to avoid data duplication.

Parameters:
  • name (str) – Name of coordinate.

  • offset (np.array or list) – Offset in term of unit-cells.

  • layer (str) – Name of layer.

  • tags (str) – List of tags (space separated string).

pos

pos is np.array((3, 1)) and is calculated from offset and position. Not to be set manually.

_get_genstring()
attributes = ['offset', 'name', 'layer', 'pos', 'tags']
eq_mod_offset(other)

Compares wether to coordinates are the same up to (modulo) a cell offset.

ff()

ff like ‘Fortran Form’

radd_ff()

Build term as if adding on the right, omit ‘+’ if 0 anyway (in Fortran Form :-)

rsub_ff()

Build term as if subtracting on the right, omit ‘-’ if 0 anyway (in Fortran Form :-)

site_offset_unpacked()
sort_key()
class kmos3.types.FixedObject(**kwargs)

Handy class that easily allows to define data structures that can only hold a well-defined set of fields

attributes = []
class kmos3.types.LatIntProcess(**kwargs)

A process which directly includes lateral interactions. In this model a bystander just defines a set of allowed species so, it allows for additional degrees of freedom here. Different lateral model can be accounted for through counters and placeholder in rate expression.

attributes = ['name', 'rate_constant', 'condition_list', 'action_list', 'bystanders', 'enabled', 'chemical_expression', 'tof_count']
class kmos3.types.Layer(**kwargs)

Represents one layer in a possibly multi-layer geometry.

Parameters:
  • name (str) – Name of layer.

  • sites (list) – Sites associated with this layer (Default: [])

add_site(*sites, **kwargs)

Adds a new site to a layer.

attributes = ['name', 'sites', 'active', 'color']
get_info()
get_site(site_name)
class kmos3.types.LayerList(**kwargs)

A list of layers

Parameters:
  • cell (np.array (3x3)) – Size of unit-cell.

  • default_layer (str.) – name of default layer.

attributes = ['cell', 'default_layer', 'name', 'representation', 'substrate_layer']
generate_coord(terms)

Expecting something of the form site_name.offset.layer and return a Coord object

generate_coord_set(size=[1, 1, 1], layer_name='default', site_name=None)

Generates a set of coordinates around unit cell of any desired size. By default it includes exactly all sites in the unit cell. By setting size=[2,1,1] one gets an additional set in the positive and negative x-direction.

set_representation(images)

FIXME: If there is more than one representation they should be sorted by their name!!!

class kmos3.types.Meta(*args, **kwargs)

Class holding the meta-information about the kMC project

add(attrib)
get_extra()
name = 'Meta'
setattribute(attr, value)
class kmos3.types.OutputItem(*args, **kwargs)

Not implemented yet

attributes = ['name', 'output']
class kmos3.types.OutputList

A dummy class, that will hold the values which are to be printed to logfile.

attributes = ['name']
class kmos3.types.Parameter(**kwargs)

A parameter that can be used in a rate constant expression and defined via some init file.

Parameters:
  • name (str) – The name of the parameter.

  • adjustable (bool) – Create controller in GUI.

  • min (float) – Minimum value for controller.

  • max (float) – Maximum value for controller.

  • scale (str) – Controller scale: ‘log’ or ‘lin’

attributes = ['name', 'value', 'adjustable', 'min', 'max', 'scale']
get_info()
on_adjustable__do_toggled(value)
on_name__content_changed(_)
class kmos3.types.ParameterList(**kwargs)

A list of parameters

attributes = ['name']
class kmos3.types.Process(**kwargs)

One process in a kMC process list

Parameters:
  • name (str) – Name of process.

  • rate_constant (str) – Expression for rate constant.

  • otf_rate (str) – Expression used to calculate rate on the fly using bystander’s configuration, otf backend only!.

  • condition_list (list.) – List of conditions (class Condition).

  • action_list (list.) – List of conditions (class Action).

  • bystander_list (list.) – List of bystanders (class Bystander), otf backend only!.

  • enabled (bool.) – Switch this process on or of.

  • chemical_expression (str.) – Chemical expression (i.e: A@site1 + B@site2 -> empty@site1 + AB@site2) to generate process from.

  • tof_count (dict.) – Stoichiometric factor for observable products {‘NH3’: 1, ‘H2Ogas’: 2}. Hint: avoid space in keys.

_get_max_d()
add_action(action)

Adds an action to a process

add_bystander(bystander)

Adds a bystander to a process

add_condition(condition)

Adds a conditions to a process

attributes = ['name', 'rate_constant', 'otf_rate', 'condition_list', 'action_list', 'bystander_list', 'enabled', 'chemical_expression', 'tof_count']
evaluate_rate_expression(parameters={})
executing_coord()
get_info()
class kmos3.types.ProcessFormSite(**kwargs)

This is just a little varient of the site object, with the sole difference that it has a layer attribute and is meant to be used in the process form. This separation was chosen, since the Site object as in the Project should not have a layer attribute to avoid data duplication but in the ProcessForm we need this to define processes

attributes = ['name', 'pos', 'tags', 'default_species', 'layer', 'color']
class kmos3.types.ProcessList(**kwargs)

A list of processes

attributes = ['name']
class kmos3.types.Project(model_name=None)

A Project is where (almost) everything comes together. A Project holds all other elements needed to describe one kMC Project ready to be manipulated, exported, or imported.

The Project class is primarily used in build files.

The overall structure is the following and was also displayed in the editor GUI (but that GUI is deprecated as of 2022).

Project:

- Meta
- Parameters
- Lattice(s)
- Species
- Processes
_get_etree_xml()

Produces an ElemenTree object representing the Project

_get_ini_string()

Return representation of model as can be written into a *.ini File.

_get_xml_string()

Produces an XML representation of the project data

add_layer(*layers, **kwargs)

Add a layer to the project. A Layer, or keywords that are passed to the Layer constructor are accepted.

Parameters:
  • layers (list) – List of layers.

  • cell (np.array (3x3)) – Size of unit-cell.

  • default_layer (str.) – name of default layer.

add_parameter(*parameters, **kwargs)

Add a parameter to the project. A Parameter, or keywords that are passed to the Parameter constructor are accepted.

Parameters:
  • name (str) – The name of the parameter.

  • value (float) – Default value of parameter.

  • adjustable (bool) – Create controller in GUI.

  • min (float) – Minimum value for controller.

  • max (float) – Maximum value for controller.

  • scale (str) – Controller scale: ‘log’ or ‘lin’

add_process(*processes, **kwargs)

Add a process to the project. A Process, or keywords that are passed to the Process constructor are accepted.

Parameters:
  • name (str) – Name of process.

  • rate_constant (str) – Expression for rate constant.

  • condition_list (list.) – List of conditions (class Condition).

  • action_list (list.) – List of conditions (class Action).

  • enabled (bool.) – Switch this process on or of.

  • chemical_expression (str.) – Chemical expression (i.e: A@site1 + B@site2 -> empty@site1 + AB@site2) to generate process from.

  • tof_count (dict.) – Stoichiometric factor for observable products {‘NH3’: 1, ‘H2O(gas)’: 2}. Hint: avoid space in keys.

add_site(**kwargs)

Add a site to the project. The arguments are

add_site(layer_name, site)

Parameters:
  • name (str) – Name of layer to add the site to.

  • site (Site) – Site instance to add.

add_species(*speciess, **kwargs)

Add a species to the project. A Species, or keywords that are passed to the Species constructor are accepted.

Parameters:
  • name (str) – Name of species.

  • representation (str) – ase.atoms.Atoms constructor describing species geometry.

  • tags (str) – Tags of species (space separated string).

clear_model(model_name='', backend='')
compile_model(code_generator='local_smart')
export_xml_file(filename, validate=True)
get_parameters(pattern=None)

Return list of parameters in Project.

Parameters:

pattern (str) – Pattern to fnmatch name of parameter against.

get_processes(pattern=None)

Return list of processes.

Parameters:

pattern (str) – Pattern to fnmatch name of process against.

get_speciess(pattern=None)

Return list of species in Project.

Parameters:

pattern (str) – Pattern to fnmatch name of process against.

import_file(filename)
import_ini_file(filename)
import_xml_file(filename)

Takes a filename, validates the content against kmc_project.dtd and import all fields into the current project tree

parse_and_add_process(string)

Generate and add processes using a shorthand notation like, e.g. :: process_name; species1A@coord1 + species2A@coord2 + … -> species1B@coord1 + species2A@coord2 + …; rate_constant_expression

.

Parameters:

string (str) – shorthand notation for process

parse_process(string)

Generate processes using a shorthand notation like, e.g. :: process_name; species1A@coord1 + species2A@coord2 + … -> species1B@coord1 + species2A@coord2 + …; rate_constant_expression

.

Parameters:

string (str) – shorthand notation for process

print_statistics()
save(filename='', validate=True)
save_model(filename='', validate=True)
set_meta(author=None, email=None, model_name=None, model_dimension=None, debug=None)
shorten_names(max_length=15)
validate_model()

Run various consistency and completeness test of the model to make sure we have a minimally complete model.

class kmos3.types.SingleLatIntProcess(**kwargs)

A process that corresponds to a single lateral interaction configuration. This is conceptually the same as the old condition/action model, just some conditions are now called bystanders.

attributes = ['name', 'rate_constant', 'condition_list', 'action_list', 'bystanders', 'enabled', 'chemical_expression', 'tof_count']
class kmos3.types.Site(**kwargs)

Represents one lattice site.

Parameters:
  • name (str) – Name of site.

  • pos (np.array or str) – Position within unit cell.

  • tags (str) – Tags for this site (space separated).

  • default_species (str) – Initial population for this site.

attributes = ['name', 'pos', 'tags', 'default_species', 'layer', 'color']
class kmos3.types.Species(**kwargs)

Class that represent a species such as oxygen, empty, … . Note: empty is treated just like a species.

Parameters:
  • name (str) – Name of species.

  • representation (str) – ase.atoms.Atoms constructor describing species geometry.

  • tags (str) – Tags of species (space separated string).

attributes = ['name', 'color', 'representation', 'tags']
class kmos3.types.SpeciesList(**kwargs)

A list of species

attributes = ['default_species', 'name']
kmos3.types.cmp_coords(self, other)
kmos3.types.create_kmc_model(model_name=None)
kmos3.types.parse_chemical_expression(eq, process, project_tree)

Evaluates a chemical expression ‘eq’ and adds conditions and actions accordingly. Rules are:

  • each chemical expression has the form

    conditions -> actions
    
  • each condition or action term has the form (regex)

    [$^]*SPECIES@SITE\.OFFSET\.LAYER
    
  • each SPECIES must have been defined before.

  • each SITE must have been defined before via the layer form.

  • an offset in units of units cell can be given as tuple such as ‘(0, 0)’

  • a condition or action term containing the default species, i.e. by default ‘empty’ may be omitted. However a term containing the omitted site and a species other then the default must exist on the opposite side of the expression

  • ^ and $ are special prefixes for the ‘creation’ and ‘annihilation’ of a site, respectively. In case of $ the species can be always omitted. In case of ^ the default species may be omitted. Creation and annihilation is only needed for lattice reconstructions/multi-lattice models and they only stand on the right-hand (i.e. action) side of the expression

  • white spaces may be used for readability but have no effect

Examples ::
  • oxygen@cus -> oxygen@bridge #diffusion

  • co@bridge -> co@cus.(-1,0) # diffusion

  • -> oxygen@cus + oxygen@bridge # adsorption

  • oxygen@cus + co@bridge -> # reaction

kmos3.types.parse_process(string, project_tree)
kmos3.types.prettify_xml(elem)

This function takes an XML document, which can have one or many lines and turns it into a well-breaked, nicely-indented string

kmos3.io

Features front-end import/export functions for kMC Projects. Currently import and export is supported to XML and export is supported to Fortran 90 source code.

class kmos3.io.ProcListWriter(data, dir)

Write the different parts of Fortran 90 code needed to run a kMC model.

_db_print(line, debug=False)

Write out debugging statement if requested.

_get_lat_int_groups()
_get_site_params()
_gpl_message()

Prints the GPL statement at the top of the source file

_otf_get_auxilirary_params(data)
_parse_otf_rate(expr, procname, data, indent=4)

Parses the otf_rate expression and returns the expression to be inserted into the associated get_rate subroutine. Additionally collects locally defined variables and the full set of used nr_<species>_<flag> variables in order to include them in the variable declarations in those functions

_parse_otf_rate_line(expr, procname, data, indent=4)

Parses an individual line of the otf_rate returning the processed line and a list of the nr_<species>_<flag> encountered

_write_optimal_iftree(items, indent, out)
_write_optimal_iftree_otf(items, indent, out)
write_proclist(smart=True, code_generator='local_smart', accelerated=False, options=None)

Write the proclist.f90 module, i.e. the rules which make up the kMC process list.

write_proclist_acf(smart=True, code_generator='local_smart')

Write the proclist_acf.f90 module, i.e. the routines to run the calculation of the autocorrelation function or to record the displacment..

write_proclist_acf_end(out)
write_proclist_constants(data, out, code_generator='local_smart', close_module=False, module_name='proclist', accelerated=False)
write_proclist_end(out)
write_proclist_generic_part(data, out, code_generator='local_smart', accelerated=False, options=None)
write_proclist_generic_subroutines(data, out, code_generator='local_smart', accelerated=False, options=None)
write_proclist_generic_subroutines_acf(data, out, code_generator='local_smart')
write_proclist_get_diff_sites_acf_otf(data, out)
write_proclist_get_diff_sites_acf_smart(data, out)
write_proclist_get_diff_sites_displacement_otf(data, out)
write_proclist_get_diff_sites_displacement_smart(data, out)
write_proclist_lat_int(data, out, debug=False, accelerated=False)

This a dumber version f the run_proc_nr routine. Though the source code it generates might be quite a bit smaller. On the downside, it might be a little less optimized though it is local in a very strict sense. [EXPERIMENTAL/UNFINISHED!!!]

write_proclist_lat_int_nli_caselist(data, lat_int_groups, progress_bar, out)

subroutine nli_<processname>

nli = number of lateral interaction inspect a local enviroment for a set of processes that only differ by lateral interaction and return the process number corrresponding to the present configuration. If no process is applicable an integer “0” is returned.

This version is the fastest found so far but has the problem that nr_of_species**nr_of_sites quickly runs over sys.max_int or whatever is the largest available integer for your Fortran compiler.

write_proclist_lat_int_nli_casetree(data, lat_int_groups, progress_bar)

Write out subroutines that do the following: Take a given cell and determine from a group a processes that only differ by lateral interaction which one is possible. This version writes out explicit ‘select case’-tree which is somewhat slower than the module version but can theoretically accomodate for infinitely many conditions for one elementary step.

If no process is applicable an integer “0” is returned.

write_proclist_lat_int_run_proc(data, lat_int_groups, progress_bar)

subroutine run_proc_<processname>(cell)

Performs the lattice and avail_sites updates for a given process.

write_proclist_lat_int_run_proc_nr(data, lat_int_groups, progress_bar, out)

subroutine run_proc_nr(proc, cell)

Central function at the beginning of each executed elementary step. Dispatches from the determined process number to the corresponding subroutine.

write_proclist_lat_int_touchup(lat_int_groups, out)

The touchup function

Updates the elementary steps that a cell can do given the current lattice configuration. This has to be run once for every cell to initialize the simulation book-keeping.

write_proclist_multilattice(data, out)
write_proclist_otf(data, out, separate_files=True, debug=False, accelerated=False, options=None)

Writes the proclist.f90 file for the otf backend

write_proclist_pars_otf(data, out, separate_files=False)

Writes the proclist_pars.f90 files which implements the module in charge of doing i/o from python evaluated parameters, to fortran and also handles rate constants update at fortran level

write_proclist_put_take(data, out)

HERE comes the bulk part of this code generator: the put/take/create/annihilation functions encode what all the processes we defined mean in terms updates for the geometry and the list of available processes

The updates that disable available process are pretty easy and flat so they cannot be optimized much. The updates enabling processes are more sophisticasted: most processes have more than one condition. So enabling one condition of a processes is not enough. We need to check if all the other conditions are met after this update as well. All these checks typically involve many repetitive questions, i.e. we will inquire the lattice many times about the same site. To mend this we first collect all processes that could be enabled and then use a heuristic algorithm (any theoretical computer scientist knows how to improve on this?) to construct an improved if-tree

write_proclist_run_proc_name_otf(data, out=None, separate_files=False, indent=4)

This routine implements the routines that execute an specific process. As with the local_smart backend, turning processes off is easy. For turning processes on, we reuse the same logic as in local_smart, but now working on whole processes, rather that with put/take single site routines. Aditionally, this routines must call the gr_<procname> routines, which are defined in the proclist_pars module

write_proclist_run_proc_nr_otf(data, out)
write_proclist_run_proc_nr_smart(data, out)
write_proclist_touchup(data, out)
write_proclist_touchup_otf(data, out)

The touchup function

Updates the elementary steps that a cell can do given the current lattice configuration. This has to be run once for every cell to initialize the simulation book-keeping.

write_settings(code_generator='lat_int', accelerated=False)

Write the kmc_settings.py. This contains all parameters, which can be changed on the fly and without recompilation of the Fortran 90 modules. In this function, “data” is an object analogous to what is normally in the variable “kmc_model” (within a build file), and it is a Project class object from types.py (just like kmc_model is). But this is not actually the same object, it is a fresh object made from the xml (or ini) file.

write_template(filename, target=None, options=None)
kmos3.io._casetree_dict(dictionary, indent='', out=None)

Recursively prints nested dictionaries.

kmos3.io._chop_line(outstr, line_length=100)
kmos3.io._flatten(L)
kmos3.io._most_common(L)
kmos3.io._print_dict(dictionary, indent='')

Recursively prints nested dictionaries.

kmos3.io.clear_model(model_name, backend='local_smart')

Delete all files of an existing model so that a directory is ready for exporting a new model. This function is normally not called directly from here. It is normally called by the Project class of types.py, which gets accessed by:

kmc_model = kmos3.create_kmc_model(model_name)
kmc_model.clear_model()
kmos3.io.export_source(project_tree, export_dir=None, code_generator=None, options=None, accelerated=False)

Export a kmos3 project into Fortran 90 code that can be readily compiled using f2py. The model contained in project_tree will be stored under the directory export_dir. export_dir will be created if it does not exist. The XML representation of the model will be included in the kmc_settings.py module.

The variable project_tree is a Project object (from types.py) and is analogous to the variable normally named kmc_model in a build file.

export_source is a central feature of the kmos3 approach. In order to generate different backend solvers, and allows additional candidates for kmc methods to be implemented.

kmos3.io.export_xml(project_tree, filename=None)

Writes a project to an XML file.

kmos3.io.import_xml(xml)
kmos3.io.import_xml_file(filename)

Imports and returns project from an XML file.

Runtime Frontend

kmos3.run

A front-end module to run a compiled kMC model. The actual model is imported in kmc_model.so and all parameters are stored in kmc_settings.py.

The model can be used directly like so:

from kmos3.model import KMC_Model
model = KMC_Model()

model.parameters.T = 500
model.do_steps(100000)
model.view()

which, of course can also be part of a python script.

The model can also be run in a different process using the multiprocessing module. This mode is designed for use with a GUI so that the CPU intensive kMC integration can run at full throttle without impeding the front-end. Interaction with the model happens through Queues.

class kmos3.run.KMC_Model(image_queue=None, parameter_queue=None, signal_queue=None, size=None, system_name='kmc_model', banner=True, print_rates=False, autosend=True, steps_per_frame=50000, random_seed=None, cache_file=None, buffer_parameter=None, threshold_parameter=None, sampling_steps=None, execution_steps=None, executed_rates_limit=None, save_limit=None)

API Front-end to initialize and run a kMC model using python bindings. Depending on the constructor call the model can be run either via directory calls or in a separate processes access via multiprocessing.Queues. Only one model instance can exist simultaneously per process.

_adjust_database()

Set the database of processes currently possible according to the current configuration.

_get_configuration()

Return current configuration of model.

Return type:

np.array

_put(site, new_species, reduce=False)

Works exactly like put, but without updating the database of available processes. This is faster for when one does a lot updates at once, however one must call _adjust_database afterwards.

Examples

# below puts a CO molecule at the `bridge` site of the lower left unit cell
model._put([0,0,0,model.lattice.bridge], model.proclist.co)
# below does the same:
model._put([0,0,0,"bridge"], "CO")

model._put([1,0,0,model.lattice.bridge], model.proclist.co )
# below puts a CO molecule at the `bridge` site one to the right

model._adjust_database() # Important !
Parameters:
  • site (list or np.array) – Site where to put the new species, i.e. [x, y, z, bridge]

  • new_species (str) – Name of new species.

  • reduce (bool) – Of periodic boundary conditions if site falls out site lattice (Default: False)

To see all the available site names, use model.settings.site_names, which come from kmc_settings.

Ex: site_names = [‘simple_cubic_hollow’] set site name as ‘model.lattice.simple_cubic_hollow’

Ex: site_names = [‘ruo2_bridge’] set site name as ‘model.lattice.ruo2_bridge’

To see all the available species names, one can use model.settings.species_tags which comes from kmc_settings.py.

Ex: species_tags = { “CO”:, “O”:, “empty”:, }

_set_configuration(config)

Set the current lattice configuration.

Expects a 4-dimensional array, with dimensions [X, Y, Z, N] where X, Y, Z are the lattice size and N the number of sites in each unit cell.

Parameters:

config (np.array) – Configuration to set for model. Shape of array has to match with model size.

create_configuration_plot(directory='./exported_configurations', plot_settings={}, showFigure=False, exportFigure=True, dimensionality=2)

Returns the spatial view of the kmc_model and make a graph named ‘plottedConfiguration.png,’ unless specified by ‘figure_name’ in plot_settings

‘coords’ is expected to be the results from get_species_coordinates(config, species, meshgrid = ‘cartesian’)

Ex: [[0 10 0] [0 11 0] [0 18 0] [1 6 0] [2 3 0] [2 11 0] [2 13 0]] -> This is CO positions in [x y z] [[0 0 0] [0 1 0] [0 2 0] [0 3 0] [0 4 0] [0 5 0] [0 6 0]] -> This is empty site positions in [x y z]

‘directory’ sets the directory name where the plot is saved

‘plot_settings’ is a dictionary that allows for the plot to change given the arguements EX: “y_label”: “test”, “x_label”: “test”, “legendLabel”: “Species”, “legendExport”: False, “legend”: True, “figure_name”: “Plot”, “dpi”: 220, “speciesName”: False

‘dimensionality’ is an integer (either 2 or 3 dimensional) for the number of cartesian dimensions.

deallocate()

Deallocate all arrays that are allocated by the Fortran module. This needs to be called whenever more than one simulation is started from one process.

Note that the currenty state and history of the system is lost after calling this method.

Note: explicit invocation was chosen over the __del__ method because there seems to easy portable way to control garbage collection.

do_acc_steps(n=10000, stats=True, save_exe=False, save_proc=0)

Propagate the model n steps using the temporal acceleration scheme.

Parameters:
  • n (int) – Number of steps to run (Default: 10000)

  • stats (logical) – Calculate statistics for the scaling factors

  • save_exe (logical) – Track ‘save_limit’ number of executions following the execution of the target process ‘save_proc’

  • save_proc (integer) – Process to be tracked

do_steps(n=10000, progress=False)

Propagate the model n steps.

Parameters:

n (int) – Number of steps to run (Default: 10000)

do_steps_time(t=1.0, n=10000)

Propagate the model t s, or n steps (whichever is achieved first. The n steps are intended to act as an upper limit to avoid infnite steps as well as for any other reason that the user may wish to limit the number of steps.

Parameters:
  • t (real) – Length of time (s) to run (Default: 1)

  • n (int) – Upper limit for number of steps to run (Default: 10000)

Returns the number of iterations executed.

double()

Double the size of the model in each direction and initialize larger model with current configuration in each copy.

dump_config(filename, directory='./exported_configurations')

Use numpy mechanism to store current configuration in a file.

Parameters:

filename (str) – Name of file, to write configuration to.

export_movie(filename='', directory='./exported_movies', resolution=150, scale=20, fps=1, frames=30, steps=1000000.0, representation='atomic', stitch=True)

Exports a series of atomic view snapshots of model instance to a subdirectory, creating png files in the exported_movie_images directory and then creates a .webm video file of all the images of the images into a video

‘filename’ sets the filename for the images in the image directory and the video

‘scale’ increases the size of each species in the structure (currently not working as desired)

‘resolution’ changes the dpi of the images (currently not working as desired)

‘fps’ sets how long each image will stay in the video

‘frames’ sets the total video length

‘steps’ is the amount of steps the model does between each image

export_picture(resolution=150, scale=20, filename='', **kwargs)

Gets the atoms objects of the kmc_model and returns a atomic view of the configuration and make a file named ‘atomic_view.png’ unless specified by ‘filename’ in the function’s argument

‘filename’ sets the filename for the images in the image directory and the video

‘scale’ increases the size of each species in the structure

(currently not working as desired)

‘resolution’ changes the dpi of the images (currently not working as desired)

get_atoms(geometry=True, tag=None, reset_time_overrun=False)

Return an Atoms object with additional information such as coverage and Turn-over-frequencies attached.

The additional attributes are:
  • info (extra tags assigned to species)

  • kmc_step

  • kmc_time

  • occupation

  • procstat

  • integ_rates

  • tof_data

tof_data contains previously defined TOFs in reaction per seconds per

cell sampled since the last call to get_atoms()

info can be used to better visualize similar looking molecule during

post-processing

procstat holds the number of times each process was executed since

last get_atoms() call.

Parameters:

geometry (bool) – Return object of current configuration (Default: True).

get_avail(arg)

Return available (enabled) processes or sites. If the argument is a sequence it is interpreted as a site (x, y, z, n). If it is an integer it is interpreted as a process.

Parameters:

arg (int or [int]) – type or process to query

get_backend()

Return name of backend that model was compiled with.

Return type:

str

get_buffer_parameter()
Return value:

See set function.

get_executed_rates(proc_nr)

Returns the row of executed_rates according to proc_nr

get_executed_rates_limit()
Return value:

See set function.

get_execution_steps()
Return value:

See set function.

get_global_configuration(filename_csv='', directory='./exported_configurations', export_csv=True, matrix_format='cartesian')

Gets each species and their respective coordinates and returns a 3d list that separates the coordinates of each species and EITHER returns a dictionary of the species’s name OR returns a meshgrid of all the species

‘filename_csv’ sets the name of the csv file that will be exported if ‘export_csv’ is true

‘directory’ sets the directory name where the .csv file is saved if ‘export_csv’ is true

‘export_csv’ determines whether the functions exports the species coordinates as a csv file

‘matrix_format’ has two types of options: meshgrid and cartesian. Cartesian return

as a csv with (x,y,z) format, and the meshgrid format

returns as a csv with a XX, YY format

EX: Cartesian Species Coordinates CO [0,2,0] CO [0,4,0] CO [0,5,0] empty [0,0,0] empty [0,1,0] empty [0,3,0]

EX: Meshgrid [[0, 1, 1, 0, 1, 0], [1, 1, 1, 0, 1, 0], [0, 1, 1, 0, 1, 0], [1, 1, 0, 0, 1, 0], [1, 0, 1, 1, 1, 0], [1, 0, 1, 0, 1, 0]]

Note 1: For this case, “0” is empty and “1” is CO. In general, the meshgrid can have higher numbers representing more than 2 species if htere are enough spaces in the model.

Note 2: The return value for get_global_configuration() and get_species_coordinates() have the same return values when setting matrix_format = ‘meshgrid’

get_local_configurations(configurationArray, radius=2, filename='', directory='./exported_configurations', export_files=True, unique_only=True, delimiter='|')

Takes in a meshgrid or _config object (from _get_configuration) and returns either the list of either all possible smaller meshgrids (i.e. the local configurations), or only the unique local configurations Currently, get_local_configurations is only compatible with 2D configurations.

‘meshgrid’ is a matrix with all the species

EX: Meshgrid [[0, 1, 1, 0, 1, 0], [1, 1, 1, 0, 1, 0], [0, 1, 1, 0, 1, 0], [1, 1, 0, 0, 1, 0], [1, 0, 1, 1, 1, 0], [1, 0, 1, 0, 1, 0]] Meshgrid can be obtained by called self.get_global_configuration(matrix_format=’meshgrid’)

‘radius’ is the distance from the central species to the adjacent species

EX: Radius = 1 [[0, 1, 1,], [1, 1, 1,], –> This is one of the local configurations of the meshgrid [0, 1, 1,]]

‘filename’ sets the filename for the array of local configurations as a .npy file

‘directory’ sets the directory name where the .npy file is saved

‘export_file’ is the argument that determines if the function will export the return value

as a .npy file

‘unique_only’ is the argument that determines if the function returns all possible local

configurations or the unique local configurations

Example of the function’s return value

EX: [[[0, 1, 1,], [1, 0, 1,], [0, 0, 1,]],

[[0, 1, 1,], [1, 1, 1,], [0, 1, 1,]],

[[1, 1, 1,], [1, 0, 0,], [0, 1, 1,]]]

Note that if a config object is passed in, rather than a meshgrid, then each element is a list rather than an integer. For example, the first row might be [[ [0,1,0], [1,1,0] …. ], However, the function will still work. So one can use model._get_configuration() and pass the output of that into get_local_configurations. A config object from model._get_configuration() should not be confused with a global configuration from model.get_global_configuration(), they are two representations of the global configuration but are very different in format.

#TODO: We should have an optional argument for the configurationArrayFormat which can be

“meshgrid” versus “_config” since _config is really an internal format. Then we can use the flag and the “try and except” will be a last resort in an else statement.

get_next_kmc_step()

Returns the next kmc step’s process and which site it would occur on, without taking the step. The output looks like this:

(Process model.proclist.o2_adsorption_bridge_right (13), Site (10, 19, 0, 1) [#781])

The process name and process number are shown.

For the site,the format is the unit cell position in cartesian x,y,z followed by the site type’s index (in this example, it is 1). As noted in the “_put()” function, the site type indexing starts at 1 (not at zero). One can use model.settings.site_names to see the site names, which come from kmc_settings. So a value of (10, 19, 0, 1) would mean unit cell 10,19,0 with site type model.settings.site_names[0] due to the different indexing.

get_occupation_header()

Return the names of the fields returned by self.get_atoms().occupation. Useful for the header line of an ASCII output.

get_param_header()

Return the names of field return by self.get_atoms().params. Useful for the header line of an ASCII output.

get_param_value(param)

Return the evaluated value of a parameter

get_sampling_steps()
Return value:

See set function.

get_save_limit()
Return value:

See set function.

get_saved_executions()

Returns a list with names and scaling factors of the processes executed after the target process in order of execution.

get_scaling_stats()

Returns the names of the process pairs that the scaling factors refer to, the average used value of the scaling factor, and the last set (=0 if never set) value of the scaling factor for each process pair in the temporal acceleration scheme.

get_species_coordinates(filename_csv='', directory='./exported_configurations', export_csv=True, matrix_format='cartesian')

Gets the species coordinates from config and EITHER returns a 3d array, where each sub array lists the coordinates for a single species on the surface OR returns a meshgrid of all the species

‘filename_csv’ sets the filename for the functions return value as a .csv file if ‘export_csv’ is true

‘directory’ sets the directory name where the .csv file is saved if ‘export_csv’ is true

‘matrix_format’ has two types of options: meshgrid and cartesian. Cartesian return as a csv where each row represents the coordinates for a single species, and the meshgrid format returns as a csv with a XX, YY format

EX: Cartesian

Ex: [[0 10 0] [0 11 0] [0 18 0] [1 6 0] [2 3 0] [2 11 0] [2 13 0]]

-> This is CO positions in [x y z] [[0 0 0] [0 1 0] [0 2 0] [0 3 0] [0 4 0] [0 5 0] [0 6 0]] -> This is empty site positions in [x y z]

EX: Meshgrid [[0, 1, 1, 0, 1, 0], [1, 1, 1, 0, 1, 0], [0, 1, 1, 0, 1, 0], [1, 1, 0, 0, 1, 0], [1, 0, 1, 1, 1, 0], [1, 0, 1, 0, 1, 0]]

Note 1: For this case, “0” is empty and “1” is CO. In general, the meshgrid can

have higher numbers representing more than 2 species if there are enough spaces in the model.

Note 2: The return value for get_global_configuration() and

get_species_coordinates() have the same return values when setting matrix_format = ‘meshgrid’

get_std_header()

Return commented line of field names corresponding to values returned in get_std_outdata

get_std_sampled_data(samples, sample_size, tof_method='integ', output='str', show_progress=False)

Sample an average model and return TOFs and coverages in a standardized format :

[parameters] [TOFs] [occupations] kmc_time kmc_step

Parameter tof_method allows to switch between two different methods for evaluating turn-over-frequencies. The default method procstat evaluates the procstat counter, i.e. simply the number of executed events in the simulated time interval. integ will evaluate the number of times the reaction could be evaluated in the simulated time interval based on the local configurations and the rate constant.

Credit for this latter method has to be given to Sebastian Matera for the idea and implementation.

In each case check carefully that the observable is sampled good enough!

Parameters:
  • samples – Number of batches to average coverages over.

  • sample_size (int) – Number of kMC steps in total.

  • tof_method (str) – Method of how to sample TOFs. Possible values are procrates or integ. While procrates only counts the processes actually executed, integ evaluates the configuration to estimate the actual rates. The latter can be several orders more efficient for very slow processes. Differences resulting from the two methods can be used as on estimate for the statistical error in samples.

get_threshold_parameter()
Return value:

See set function.

get_tof_header()

Return the names of the fields returned by self.get_atoms().tof_data. Useful for the header line of an ASCII output.

halve()

Halve the size of the model and initialize each site in the new model with a species randomly drawn from the sites that are reduced onto one. It is necessary that the simulation size is even.

inverse()

Print short summary of current parameters and inverse rate constants. It is advisable to include this at the beginning of every generated data file for later reconstruction.

load_config(filename, directory='./exported_configurations')

Use numpy mechanism to load configuration from a file. User must ensure that size of stored configuration is correct.

Parameters:

filename (str) – Name of file, to write configuration to.

nr2site(n)

Accepts a site index and return the site in human readable coordinates.

Parameters:

n (int) – Index of site.

Return type:

str

peek(*args, **kwargs)

Creates a static image of the model in a popup window.

pickle_export_atoms(filename='', directory='./exported_configurations')
play_ascii_movie(frames=30, steps=1, site=0, delay=0.1, species=None, hexagonal=False)

Shows a series of model snapshots in the current terminal. ‘frames’ sets the total video length ‘steps’ is the number of steps the model does between each image ‘site’ is the site of interest to animate ‘species’ is a list of species to display (default is all species)

plot_configuration(filename='', directory='./exported_configurations', resolution=150, scale=20, representation='spatial', plot_settings={}, showFigure=False, exportFigure=True)

Either calls create_configuration_plot() to create the spatial representation of the model, or calls export_picture() to create the atomic representation of the model

‘filename’ sets the filename for the plot

‘directory’ sets the directory name where the plot is saved

‘scale’ increases the size of each species in the structure (currently not working as desired)

‘resolution’ changes the dpi of the images (currently not working as desired) Note: ‘resolution’ and ‘scale’ are strictly for the atomic view

‘representation’ is an optional argument for spatial and atomic view You should specify as ‘atomic’ to see the atomic view. Leaving representation empty returns spatial view by default.

‘plot_settings’ is a dictionary that allows for the plot to change given the arguements EX: “y_label”: “test”, “x_label”: “test”, “legendLabel”: “Species”, “legendExport”: False, “legend”: True, “figure_name”: “Plot”, “dpi”: 220, “speciesName”: False

post_mortem(steps=None, propagate=False, err_code=None)

Accepts an integer and generates a post-mortem report by running that many steps and returning which process would be executed next without executing it.

Parameters:
  • steps (int) – Number of steps to run before exit occurs (Default: None).

  • propagate (bool) – Run this one more step, where error occurs (Default: False).

  • err_code (str) – Error code generated by backend if project.meta.debug > 0 at compile time.

print_accum_rate_summation(order='-rate', to_stdout=True)

Shows rate individual processes contribute to the total rate

The optional argument order can be one of: name, rate, rate_constant, nrofsites. You precede each keyword with a ‘-’, to show in decreasing order. Default: ‘-rate’. Possible values are rate, rate_constant, name, nrofsites.

print_adjustable_parameters(match=None, to_stdout=True)

Print those methods that are adjustable via the GUI.

Parameters:

pattern (str) – fname pattern to limit the parameters.

print_coverages(to_stdout=True)

Show coverages (per unit cell) for each species and site type for current configurations.

print_kmc_state(to_stdout=True)

Shows current kmc step and kmc time.

print_proc_pair_eq()

Prints the names of the forward and reverse process in each pair along with a logical for the pair that is True if the process pair is equilibrated and False if the process pair is non-equilibrated.

print_procstat(to_stdout=True)

Show how often each process has been executed.

print_scaling_factors()

Prints the names of the forward and reverse process in each pair along with the scaling factor for that pair.

print_scaling_stats()

Print the average used and last set (=0 if never set) value of the scaling factor for each process pair in the temporal acceleration scheme.

print_state_summary(order='-rate', to_stdout=True, show=False, print_parameters=False)

Show summary of current model state by showing - parameters (external, optional) - number of times each elementary process has been executed - coverage - kmc step and kmc time - fire up ASE window with current lattice configuration

procstat_normalized(match=None)

Print an overview view process names along with the number of times it has been executed divided by the current rate constant times the kmc time.

Can help to find those processes which are kinetically hindered.

Parameters:

match (str) – fname pattern to filter matching parameter name.

procstat_pprint(match=None)

Print an overview view process names along with the number of times it has been executed.

Parameters:

match (str) – fname pattern to filter matching parameter name.

put(site, new_species, reduce=False)

Puts new_species at site. The site is given by 4-entry sequence like [x, y, z, n], where the first 3 entries define the unit cell from 0 to the number of unit cells in the respective direction. And n specifies the site within the unit cell.

The database of available processes will be updated automatically. For doing many put and a single update, see the _put() function.

Examples

# below puts a CO molecule at the `bridge` site of the lower left unit cell
model.put([0,0,0,model.lattice.bridge], model.proclist.co)
# below does the same:
model.put([0,0,0,"bridge"], "CO")

model.put([1,0,0,model.lattice.bridge], model.proclist.co )
# below puts a CO molecule at the `bridge` site one to the right
Parameters:
  • site (list or np.array) – Site where to put the new species, i.e. [x, y, z, bridge]

  • new_species (str) – Name of new species.

  • reduce (bool) – Of periodic boundary conditions if site falls out site lattice (Default: False)

To see all the available site names, use model.settings.site_names, which come from kmc_settings.

Ex: site_names = [‘simple_cubic_hollow’] set site name as ‘model.lattice.simple_cubic_hollow’

Ex: site_names = [‘ruo2_bridge’] set site name as ‘model.lattice.ruo2_bridge’

To see all the available species names, one can use model.settings.species_tags which comes from kmc_settings.py.

Ex: species_tags = {

“CO”:, “O”:, “empty”:, }

rate_ratios(interactive=False)
reset()

Reset the KMC_model to the initial settings.

run()

Runs the model indefinitely. To control the simulations, model must have been initialized with proper Queues.

run_proc_nr(proc, site)
set_buffer_parameter(value=1000)
Parameters:

value (int > 1) – The value of the buffer parameter determines how many times faster equilibrated reaction steps will be than the non-equilibrated steps after scaling. It could also be called the target timescale disparity. (Default: 1000)

set_debug_level(value=0)

Set the debug level in the acceleration scheme. If a value larger than 0 is set, certain variables in base.f90 will be printed. For a value of 1 variables will be printed every time reactions are scaled or unscaled. For a value of 2 variables will be printed for every accelerated kmc step. Possible values: 0, 1, 2

set_executed_rates_limit(value=10000)

Sets the width of the executed_rates matrix. (Default: 10000) :type value: int > 1

set_execution_steps(value=200)
Parameters:

value (int > 1) – The number of steps to sample during the full simulation period before assessing the equilibrium of a process. See also threshold_parameter. This is also the number of steps of either the forward or the reverse reaction that must have been executed within the current superbasin in order for the scaling of the rate constant to be carried out. Note that changing this parameter will lead to the model being reset. (Default: 200)

set_sampling_steps(value=1000)
Parameters:

value (int > 1) – The number of steps to sample before scaling the rate constants of equilibrated processes. Note that everytime a non-equilibrated process is executed the counter of sampling steps is reset to 0. (Default: 1000)

set_save_limit(value=1000)
Parameters:

value (int > 1) – The number of executions following a target execution to save. Note that changing this parameter will lead to the model being reset. (Default: 1000)

set_threshold_parameter(value=0.2)
Parameters:

value (float > 0) – A pair of processes are flagged as equilibrated if within the last execution_steps the absolute value of the number of forward minus reverse executions divided by execution_steps is less than the threshold parameter. (Default: 0.2)

show()

Creates a static image of the model in a popup window (this is a duplicate command of ‘peek’ created for convenience).

show_ascii_picture(site, species, hexagonal=False)

Shows an ascii picture of the current configuration.

switch_surface_processes_off()

Set rate constant to zero if process has ‘diff’ or ‘react’ in the name.

switch_surface_processes_on()
view()

Start current model in live view mode.

xml()

Returns the XML representation that this model was created from.

Return type:

str

class kmos3.run.Model_Rate_Constants

Holds all rate constants currently associated with the model. To inspect the expression and current settings of it you can just call it as a function with a (glob) pattern that matches the desired processes, e.g.

model.rate_constant('*ads*')

could print all rate constants for adsorption. Given of course that ‘ads’ is part of the process name. The just get the rate constant for one specific process you can use

model.rate_constant.by_name("<process name>")

To set rate constants manually use

model.rate_constants.set("<pattern>", <rate-constant (expr.)>)
__call__(pattern=None, interactive=False, model=None)

Return rate constants.

Parameters:
  • pattern (str) – fname pattern to filter matching parameter name.

  • model (kmos3 Model) – runtime instance of kMC to extract rate constants from (optional)

by_name(proc)

Return rate constant currently set for proc

Parameters:

proc (str) – Name of process.

inverse(interactive=False)

Return inverse list of rate constants.

names(pattern=None)

Return names of processes that match pattern.

Parameters:

pattern (str) – fname pattern to filter matching parameter name.

set(pattern, rate_constant, parameters=None)

Set rate constants. Pattern can be a glob pattern, and the rate constant will be applied to all processes, where the pattern matches. The rate constant can be either a number or a rate expression.

Parameters:
  • pattern (str) – fname pattern that selects the process affected.

  • rate_constant (str or float) – Rate constant to be set.

  • parameters (list) – List of parameters to be used when evaluating expression.

class kmos3.run.Model_Parameters(print_rates=True)

Holds all user defined parameters of a model in concise form. All user defined parameters can be accessed and set as attributes, like so

model.parameters.<parameter> = X.Y
__call__(match=None, interactive=False)

Return parameters that match ‘pattern’

Parameters:

match (str) – fname pattern to filter matching parameter name.

names(pattern=None)

Return names of parameters that match pattern

Parameters:

pattern (str) – fname pattern to filter matching parameter name.

class kmos3.run.ModelRunner

Setup and initiate many runs in parallel over a regular grid of parameters. A standard type of script is given below.

To allow execution from multiple hosts connected to the same filesystem calculated points are blocked via <classname>.lock. To redo a calculation <classname>.dat and <classname>.lock should be moved out of the way

from kmos3.run import ModelRunner, PressureParameter, TemperatureParameter

class ScanKinetics(ModelRunner):
    p_O2gas = PressureParameter(1)
    T = TemperatureParameter(600)
    p_COgas = PressureParameter(min=1, max=10, steps=40)
    # ... other parameters to scan

ScanKinetics().run(init_steps=1e7, sample_steps=1e7, cores=4)
__product(*args, **kwds)

Manual implementation of itertools.product for python <= 2.5

__run_sublist(sublist, init_steps, sample_steps, samples, random_seed=None)

Run sampling run for a list of parameter-tuples.

Parameters:
  • init_steps (int) – Steps to run model before sampling (.ie. to reach steady-state).

  • sample_steps (int) – Number of steps to sample over.

  • cores (int) – Number of parallel processes to launch.

  • samples (int) – Number of samples. Use more samples if precise coverages are needed.

__split_seq(seq, size)

Split a list into n chunks of roughly equal size.

__touch(fname, times=None)

Pythonic version of Unix touch.

Parameters:
  • fname (str) – filename.

  • times (datetime timestamp) – timestamp (Default: None meaning now).

parameters = {}
plot(rcParams=None, touchup=None, filename=None, backend='Agg', suffixes=['png', 'pdf', 'eps'], variable_parameters=None, fig_width_pt=246.0, plot_tofs=None, plot_occs=None, occ_xlabel=None, occ_ylabel=None, tof_xlabel=None, tof_ylabel=None, label=None, sublabel=None, arrhenius=False)

Plot the generated data using matplotlib. By default we will try to generate publication quality output of the specified TOFs and coverages.

run(init_steps=100000000.0, sample_steps=100000000.0, cores=4, samples=1, random_seed=None)

Launch the ModelRunner instance. Creates a regular grid over all ModelParameters defined in the ModelRunner class.

Parameters:
  • init_steps (int) – Steps to run model before sampling (.ie. to reach steady-state). (Default: 1e8)

  • sample_steps (int) – Number of steps to sample over (Default: 1e8)

  • cores (int) – Number of parallel processes to launch.

  • samples (int) – Number of samples. Use more samples if precise coverages are needed (Default: 1).

runner_name = 'ModelRunner'
class kmos3.run.ModelParameter(min, max=None, steps=1, type=None, unit='')

A model parameter to be scanned. If instantiated with only one value this parameter will be fixed at this value.

Use a subclass for specific type of grid.

Parameters:
  • min (float) – Minimum value for this parameter.

  • max (float) – Maximum value for this parameter (Default: min)

  • steps (int) – Number of steps between minimum and maximum.

get_grid()
class kmos3.run.PressureParameter(*args, **kwargs)

Create a grid of p in [p_min, p_max] such that ln({p}) is a regular grid.

get_grid()
class kmos3.run.TemperatureParameter(*args, **kwargs)

Create a grid of p in [T_min, T_max] such that ({T})**(-1) is a regular grid.

get_grid()
class kmos3.run.LinearParameter(*args, **kwargs)

Create a regular grid between min and max.

get_grid()
class kmos3.run.LogParameter(*args, **kwargs)

Create a log grid between 10^min and 10^max (like np.logspace)

get_grid()

kmos3.view

Run and view a kMC model. For this to work one needs a kmc_model.(so/pyd) and a kmc_settings.py in the import path.

class kmos3.view.AutoScrollbar(*args, **kwargs)

A scroll bar that disappears/appears depending on whether the full content is visible.

Construct a scrollbar widget with the parent MASTER.

Valid resource names: activebackground, activerelief, background, bd, bg, borderwidth, command, cursor, elementborderwidth, highlightbackground, highlightcolor, highlightthickness, jump, orient, relief, repeatdelay, repeatinterval, takefocus, troughcolor, width.

set(first, last)

Set the fractional values of the slider position (upper and lower ends as value between 0 and 1).

class kmos3.view.KMC_PlotBox(master, image_queue)

The part of the viewer that displays the model’s current coverage and turnover frequencies.

on_resize(event)

Update the canvas size.

update_plots()

Update the coverage and TOF plots.

class kmos3.view.KMC_Viewer(signal_queue, parameter_queue, image_queue)

A graphical front-end to run, manipulate and view a kMC model.

exit()

Exit the viewer application cleanly killing all subprocesses before the main process.

parameter_callback(name, value)

Sent (updated) parameters to the model process.

class kmos3.view.ParamSlider(master, name, value, xmin, xmax, scale_input, parameter_callback, width=200, resolution=1000)

A horizontal slider bar allowing the user to adjust a model parameter at runtime.

Construct a scale widget with the parent MASTER.

Valid resource names: activebackground, background, bigincrement, bd, bg, borderwidth, command, cursor, digits, fg, font, foreground, from, highlightbackground, highlightcolor, highlightthickness, label, length, orient, relief, repeatdelay, repeatinterval, resolution, showvalue, sliderlength, sliderrelief, state, takefocus, tickinterval, to, troughcolor, variable, width.

linlog_scale_format()

Format a model parameter’s name for display above the slider bar.

value_changed(event)

Handle the event, that slider bar has been dragged.

kmos3.view.main(model=None, steps_per_frame=50000)

The entry point for the kmos3 viewer application. In order to run and view a model the corresponding kmc_settings.py and kmc_model.(so/pyd) must be in the current import path, e.g.

from sys import path
path.append('/path/to/model')
from kmos3.view import main
main() # launch viewer

kmos3.cli

Entry point module for the command-line interface. The kmos3 executable should be on the program path, import this modules main function and run it.

To call kmos3 commands as you would from the shell, use

kmos3.cli.main('...')

You may also use the syntax kmos3.export(”…”) for any cli command.

kmos3.cli._benchmark(args)
kmos3.cli._build(args)
kmos3.cli._export(args)
kmos3.cli._import(args)
kmos3.cli._parser()
kmos3.cli._rebuild(args)
kmos3.cli._settings_export(args)
kmos3.cli._shell(args)
kmos3.cli._view(args)
kmos3.cli._xml(args)
kmos3.cli.get_options(args=None, get_parser=False)

Get the default options and parse the the options provided on STDIN.

kmos3.cli.main(args=None)

The CLI main entry point function.

The optional argument args, can be used to directly supply command line argument like

$ kmos3 <args>

otherwise args will be taken from STDIN.

kmos3.cli.shell(banner)

Wrapper around interactive ipython shell.

kmos3.utils

Several utility functions that do not seem to fit somewhere else.

class kmos3.utils.CorrectlyNamed

Syntactic Sugar class for use with kiwi, that makes sure that the name field of the class has a name field, that always complys with the rules for variables.

on_name__validate(_, name)

Called by kiwi upon chaning a string

kmos3.utils.T_grid(T_min, T_max, n)
kmos3.utils.build(options)

Build binary with f2py binding from complete set of source file in the current directory.

kmos3.utils.col_str2tuple(hex_string)

Convenience function that turns a HTML type color into a tuple of three float between 0 and 1

kmos3.utils.col_tuple2str(tup)

Convenience function that turns a HTML type color into a tuple of three float between 0 and 1

kmos3.utils.evaluate_kind_values(infile, outfile)

Go through a given file and dynamically replace all selected_int/real_kind calls with the dynamically evaluated fortran code using only code that the function itself contains.

kmos3.utils.evaluate_template(template, escape_python=False, **kwargs)

Very simple template evaluation function using only exec and str.format()

There are two flavors of the template language, depending on whether the python parts or the template parts are escaped.

A template can use the full python syntax. Every line starts with ‘#@ ‘ is interpreted as a template line. Please use proper indentation before and note the space after ‘#@’.

The template lines are converted to TEMPLATE_LINE.format(locals()) and thefore every variable in the template line should be escape with {}.

A valid template could be

for i in range:

#@ Hello World {i}

kmos3.utils.get_ase_constructor(atoms)

Return the ASE constructor string for atoms.

kmos3.utils.jmolcolor_in_hex(i)

Return a given jmol color in hexadecimal representation.

kmos3.utils.p_grid(p_min, p_max, n)
kmos3.utils.product(*args, **kwds)

Take two lists and return iterator producing all combinations of tuples between elements of the two lists.

kmos3.utils.split_sequence(seq, size)

Take a list and a number n and return list divided into n sublists of roughly equal size.

kmos3.utils.timeit(func)

Generic timing decorator

To stop time for function call f just

from kmos3.utils import timeit
@timeit
def f():
    ...
kmos3.utils.write_py(fileobj, images, **kwargs)

Write a ASE atoms construction string for images into fileobj.

Connected Variables

The connected_variables dictionary allows to pass string-writable objects created during the model building into the runtime environment. This can be useful if a person needs access to some data structures (like lists of surrounding sites) during runtime. Dictionaries, strings, and lists can be passed. For more complex variables, one could pass the name of a pickle file. This feature is used for the surroundingSitesDict.

The basic syntax in a build_file would be as follows:

kmc_model = kmos3.create_kmc_model(model_name)
kmc_model.connected_variables['frog_list'] = [1,2,3,4]

Then, during runtime, one could do the following:

print(model.connected_variables['frog_list'])

kmos3 KMC Project DTD

The central storage and exchange format is XML. XML was chosen over JSON, pickle, or alike because it seems to be one of the most flexible and universal format with good methods to define the overall structure of the data.

One way to define an XML format is by using a document type description (DTD) and in fact at every import a kmos3 file is validated against the DTD below.

<!ELEMENT kmc (meta?,species_list?,parameter_list?, lattice, process_list?,output_list?,connected_variables?)>
  <!ATTLIST kmc
  version CDATA #REQUIRED
  >
  <!ELEMENT meta EMPTY>
  <!ATTLIST meta
    author CDATA #IMPLIED
    debug CDATA #IMPLIED
    email CDATA #IMPLIED
    model_dimension CDATA #IMPLIED
    model_name CDATA #IMPLIED
  >

  <!ELEMENT species_list (species)*>
  <!ATTLIST species_list
    default_species CDATA #IMPLIED
  >
  <!ELEMENT species EMPTY>
  <!ATTLIST species
    name CDATA #REQUIRED
    color CDATA #IMPLIED
    representation CDATA #IMPLIED
    tags CDATA #IMPLIED
  >
  <!ELEMENT parameter_list (parameter)*>
  <!ELEMENT parameter EMPTY>
  <!ATTLIST parameter
    name CDATA #REQUIRED
    value CDATA #IMPLIED
    adjustable CDATA #IMPLIED
    min CDATA #IMPLIED
    max CDATA #IMPLIED
    scale CDATA #IMPLIED
  >
  <!ELEMENT lattice (layer)*>
  <!ATTLIST lattice
    cell_size CDATA #REQUIRED
    default_layer CDATA #REQUIRED
    substrate_layer CDATA #IMPLIED
    representation CDATA #IMPLIED
  >
  <!ELEMENT layer (site)*>
  <!ATTLIST layer
    name CDATA #REQUIRED
    grid CDATA #IMPLIED
    grid_offset CDATA #IMPLIED
    color CDATA #IMPLIED
  >
  <!ELEMENT site EMPTY>
  <!ATTLIST site
    pos CDATA #REQUIRED
    type CDATA #REQUIRED
    tags CDATA #IMPLIED
    default_species CDATA #IMPLIED
  >
  <!ELEMENT process_list (process)*>
    <!ELEMENT process (condition|action|bystander)*>
      <!ATTLIST process
        name CDATA #REQUIRED
        rate_constant CDATA #REQUIRED
	 otf_rate CDATA #IMPLIED
        enabled CDATA #IMPLIED
        tof_count CDATA #IMPLIED
      >
      <!ELEMENT condition EMPTY>
      <!ATTLIST condition
        coord_name CDATA #REQUIRED
        coord_layer CDATA #REQUIRED
        coord_offset CDATA #REQUIRED
        species CDATA #REQUIRED
        implicit CDATA #IMPLIED
      >
      <!ELEMENT action EMPTY>
      <!ATTLIST action
        coord_name CDATA #REQUIRED
        coord_layer CDATA #REQUIRED
        coord_offset CDATA #REQUIRED
        species CDATA #REQUIRED
      >
      <!ELEMENT bystander EMPTY>
      <!ATTLIST bystander
        coord_name CDATA #REQUIRED
        coord_layer CDATA #REQUIRED
        coord_offset CDATA #REQUIRED
	allowed_species CDATA #REQUIRED
	flag CDATA #IMPLIED
      >

      <!ELEMENT output_list (output)*>
        <!ELEMENT output EMPTY>
        <!ATTLIST output
          item CDATA #REQUIRED
        >
  <!ELEMENT connected_variables (string)*>
  <!ATTLIST connected_variables
    connected_variables_string CDATA #IMPLIED
  >

Backends

In general, the backend includes all functions that are implemented in Fortran90, which therefore should not have to be changed by hand often. The backend is divided into three modules, which import each other in the following way

base <- lattice <- proclist

The key for this division is reusability of the code. The base module implement all aspects of the KMC code, which do not depend on the described model. Thus it “never” has to change. The latttice module basically repeats all methods of the base model in terms of lattice coordinates. Thus the lattice module only changes, when the geometry of the model changes, e.g., when you add or delete sites. The proclist module implements the process list, that is the species or states each site can have and the elementary steps. Typically that changes most often while developing a model.

The rate constants and physical parameters of the system are not implemented in the backend at all, since in the physical sense they are too high-level to justify encoding and compilation at the Fortran level and so they are typical read and parsed from a python script.

The kmos3.run.KMC_Model class implements a convenient interface for most of these functions, however all public methods (in Fortran called subroutines) and variables can also be accessed directly like so

from kmos3.run import KMC_Model
model = KMC_Model(print_rates=False, banner=False)
model.base.<TAB>
model.lattice.<TAB>
model.proclist.<TAB>

which works best in conjunction with ipython.

local_smart

kmos3/base

The base kMC module, which implements the kMC method on a d = 1 lattice. Virtually any lattice kMC model can be build on top of this. The methods offered are:

  • de/allocation of memory

  • book-keeping of the lattice configuration and all available processes

  • updating and tracking kMC time, kMC step and wall time

  • saving and reloading the current state

  • determine the process and site to be executed

base/accum_rates

Stores the accumulated rate constant multiplied with the number of sites available for that process to be used by determine_procsite. Let \mathbf{c} be the rate constants \mathbf{n} the number of available sites, and \mathbf{a} the accumulated rates, then a_{i} is calculated according to a_{i}=\sum_{j=1}^{i} c_{j} n_{j}.

base/add_proc

The main idea of this subroutine is described in del_proc. Adding one process to one capability is programmatically simpler since we can just add it to the end of the respective array in avail_sites.

  • proc positive integer number that represents the process to be added.

  • site positive integer number that represents the site to be manipulated

base/allocate_system

Allocates all book-keeping structures and stores local copies of system name and size(s):

  • systen_name identifier of this simulation, used as name of punch file

  • volume the total number of sites

  • nr_of_proc the total number of processes

base/assertion_fail

Function that shall be used by all parts of the program to print a proper message in case some assertion fails.

  • a condition that is supposed to hold true

  • r message that is printed to the poor user in case it fails

base/avail_sites

Main book-keeping array that stores for each process the sites that are available and for each site the address in this very array. The meaning of the fields are:

avail_sites(proc, field, switch)

where:

  • proc – refers to a process in the process list

  • the field within the process, but the meaning differs as explained under ‘switch’

  • switch – can be either 1 or 2 and switches between (1) the actual numbers of the sites, which are available and filled in from the left but in whatever order they come or (2) the location where the site is stored in (1).

base/can_do

Returns true if ‘site’ can do ‘proc’ right now

  • proc integer representing the requested process.

  • site integer representing the requested site.

  • can writeable boolean, where the result will be stored.

base/deallocate_system

Deallocate all allocatable arrays: avail_sites, lattice, rates, accum_rates, integ_rates, procstat.

none

base/del_proc

del_proc delete one process from the main book-keeping array avail_sites. These book-keeping operations happen in O(1) time with the help of some more book-keeping overhead. avail_sites stores for each process all sites that are available. The array for each process is filled from the left, but sites generally not ordered. With this determine_procsite can effectively pick the next site and process. On the other hand a second array (avail_sites(:,:,2) ) holds for each process and each site, the location where it is stored in avail_site(:,:,1). If a site needs to be removed this subroutine first looks up the location via avail_sites(:,:,1) and replaces it with the site that is stored as the last element for this process.

  • proc positive integer that states the process

  • site positive integer that encodes the site to be manipulated

base/determine_procsite

Expects two random numbers between 0 and 1 and determines the corresponding process and site from accum_rates and avail_sites. Technically one random number would be sufficient but to circumvent issues with wrong interval_search_real implementation or rounding errors I decided to take two random numbers:

  • ran_proc Random real number from \in[0,1] that selects the next process

  • ran_site Random real number from \in[0,1] that selects the next site

  • proc Return integer \in[1,\mathrm{nr\_of\_proc}

  • site Return integer \in [1,\mathrm{volume}

base/get_accum_rate

Return accumulated rate at a given process.

  • proc_nr integer representing the requested process.

  • return_accum_rate writeable real, where the requested accumulated rate will be stored.

base/get_avail_site

Return field from the avail_sites database

  • proc_nr integer representing the requested process.

  • field integer for the site at question

  • switch 1 or 2 for site or storage location

base/get_integ_rate

Return integrated rate at a given process.

  • proc_nr integer representing the requested process.

  • return_integ_rate writeable real, where the requested integrated rate will be stored.

base/get_kmc_step

Return the current kmc_step

  • kmc_step Writeable integer

base/get_kmc_time

Returns current kmc_time as rdouble real as defined in kind_values.f90.

  • return_kmc_time writeable real, where the kmc_time will be stored.

base/get_kmc_time_step

Returns current kmc_time_step (the time increment).

  • return_kmc_step writeable real, where the kmc_time_step will be stored.

base/get_kmc_volume

Return the total number of sites.

  • volume Writeable integer.

base/get_nrofsites

Return how many sites are available for a certain process. Usually used for debugging

  • proc integer representing the requested process

  • return_nrofsites writeable integer, where nr of sites gets stored

base/get_procstat

Return process counter for process proc as integer.

  • proc integer representing the requested process.

  • return_procstat writeable integer, where the process counter will be stored.

base/get_rate

Return rate of given process.

  • proc_nr integer representing the requested process.

  • return_rate writeable real, where the requested rate will be stored.

base/get_species

Return the species that occupies site.

  • site integer representing the site

base/get_system_name

Return the systems name, that was specified with base/allocate_system

  • system_name Writeable string of type character(len=200).

base/get_walltime

Return the current walltime.

  • return_walltime writeable real where the walltime will be stored.

base/increment_procstat

Increment the process counter for process proc by one.

  • proc integer representing the process to be increment.

base/integ_rates

Stores the time-integrated rates (non-normalized to surface area) Used to determine reaction rates, i.e. average number of reactions per unit surface and time. Let \mathbf{a} the integrated rates, \mathbf{c} be the rate constants, \mathbf{n}_i the number of available sites during kMC-time interval i, \{\Delta t_i\} the corresponding timesteps then a_{i}(t) at the time t=\sum_{i=1}\Delta t_i is calculated according to a_{i}(t)=\sum_{i=1} c_{i} n_{i}\Delta t_i.

base/interval_search_real

This is basically a standard binary search algorithm that expects an array of ascending real numbers and a scalar real and return the key of the corresponding field, with the following modification :

  • the value of the returned field is equal of larger of the given value. This is important because the given value is between 0 and the largest value in the array and otherwise the last field is never selected.

  • if two or more values in the array are identical, the function return the index of the leftmost of those field. This is important because having field with identical values means that all field except the leftmost one do not contain any sites. Refer to update_accum_rate to understand why.

  • the value of the returned field may no be zero. Therefore the index the to be equal or larger than the first non-zero field.

However: as everyone knows the binary search is trickier than it appears at first site especially real numbers. So intensive testing is suggested here!

  • arr real array of type rsingle (kind_values.f90) in monotonically (not strictly) increasing order

  • value real positive number from [0, max_arr_value]

base/kmc_step

Number of kMC steps executed.

base/kmc_time

Simulated kMC time in this run in seconds.

base/kmc_time_step

The time increment of the current kMC step.

base/lattice

Stores the actual physical lattice in a 1d array, where the value on each slot represents the species on that site.

Species constants can be conveniently defined in lattice_… and later used directly in the process list.

base/nr_of_proc

Total number of available processes.

base/nr_of_sites

Stores the number of sites available for each process.

base/procstat

Stores the total number of times each process has been executed during one simulation.

base/rates

Stores the rate constants for each process in s^-1.

base/reload_system

Restore state of simulation from *.reload file as saved by save_system(). This function also allocates the system’s memory so calling allocate_system again, will cause a runtime failure.

  • system_name string of 200 characters which will make the reload_system look for a file called ./<system_name>.reload

  • reloaded logical return variable, that is .true. reload of system could be completed successfully, and .false. otherwise.

base/replace_species

Replaces the species at a given site with new_species, given that old_species is correct, i.e. identical to the site that is already there.

  • site integer representing the site

  • old_species integer representing the species to be removed

  • new_species integer representing the species to be placed

base/reset_site

This function is a higher-level function to reset a site as if it never existed. To achieve this the species is set to null_species and all available processes are stripped from the site via del_proc.

  • site integer representing the requested site.

  • species integer representing the species that ought to be at the site, for consistency checks

base/save_system

save_system stores the entire system information in a simple ASCII filed names <system_name>.reload. All fields except avail_sites are stored in the simple scheme:

variable value

In the case of array variables, multiple values are seperated by one or more spaces, and the record is terminated with a newline. The variable avail_sites is treated slightly differently, since printed on a single line it is almost impossible to interpret from the ASCII files. Instead each process starts a new line, and the first number on the line stands for the process number and the remaining fields, hold the values.

none

base/set_kmc_step

Sets the current kmc_step

  • kmc_step Writeable integer

base/set_kmc_time

Sets current kmc_time as rdouble real as defined in kind_values.f90.

  • new readable real, that the kmc time will be set to

base/set_rate_const

Allows to set the rate constant of the process with the number proc_nr.

  • proc_n The process number as defined in the corresponding proclist_ module.

  • rate the rate in s^{-1}

base/set_system_name

Set the systems name. Useful in conjunction with base.save_system to save *.reload files under a different name than the default one.

  • system_name Readable string of type character(len=200).

base/start_time

CPU time spent in simulation at least reload.

base/system_name

Unique indentifier of this simulation to be used for restart files. This name should not contain any characters that you don’t want to have in a filename either, i.e. only [A-Za-z0-9_-].

base/update_accum_rate

Updates the vector of accum_rates.

none

base/update_clocks

Updates walltime, kmc_step and kmc_time.

  • ran_time Random real number \in [0,1]

base/update_integ_rate

Updates the vector of integ_rates.

none

base/volume

Total number of sites.

base/walltime

Total CPU time spent on this simulation.

kmos3/lattice

Implements the mappings between the real space lattice and the 1-D lattice, which kmos3/base operates on. Furthermore replicates all geometry specific functions of kmos3/base in terms of lattice coordinates. Using this module each site can be addressed with 4-tuple (i, j, k, n) where i, j, k define the unit cell and n the site within the unit cell.

lattice/allocate_system

Allocates system, fills mapping cache, and checks whether mapping is consistent

none

lattice/calculate_lattice2nr

Maps all lattice coordinates onto a continuous set of integer \in [1,volume]

  • site integer array of size (4) a lattice coordinate

lattice/calculate_nr2lattice

Maps a continuous set of of integers \in [1,volume] to a 4-tuple representing a lattice coordinate

  • nr integer representing the site index

lattice/deallocate_system

Deallocates system including mapping cache.

none

lattice/default_layer

The layer in which the model is initially in by default (only relevant for multi-lattice models).

lattice/lattice2nr

Caching array holding the mapping from index to lattice coordinate: (x, y, z, n) -> i.

lattice/model_dimension

Store the number of dimensions of this model: 1, 2, or 3

lattice/nr2lattice

Caching array holding the mapping from index to lattice coordinate: i -> (x, y, z, n).

lattice/nr_of_layers

Constant storing the number of layers (for multi-lattice models > 1)

lattice/site_positions

The positions of (adsorption) site in the unit cell in fractional coordinates.

lattice/spuck

spuck = Sites Per Unit Cell Konstant The number of sites per unit cell, i.e. for coordinate notation (x, y, n) this is the maximum value of n.

lattice/system_size

Stores the current size of the allocated system lattice (x, y, z) in an integer array. In low-dimensional system, corresponding entries will be set to 1. Note that this should be thought of as a read-only variable. Changing its value at model runtime will not the indented effect of actually changing the simulated lattice. The definitive location for custom lattice size is simulation_size in kmc_settings.py.

If the system size shall be changed programmatically, it needs to happen before the KMC_Model is instantiated and Fortran array are allocated accordingly, like to

#!/usr/bin/env python3

import kmc_settings import kmos3.run

kmc_settings.simulation_size = 9, 9, 4

with kmos3.run.KMC_Model() as model:

print(model.lattice.system_size)))`

lattice/unit_cell_size

The dimensions of the unit cell (e.g. in Angstrom) of the unit cell.

kmos3/proclist

Implements the kMC process list.

proclist/do_kmc_step

Performs exactly one kMC step. * first update clock * then configuration sampling step * last execute process

none

proclist/do_kmc_steps

Performs n kMC step. If one has to run many steps without evaluation do_kmc_steps might perform a little better.

  • first update clock

  • then configuration sampling step

  • last execute process

n : Number of steps to run

proclist/do_kmc_steps_time

Performs a variable number of KMC steps to try to match the requested simulation time as closely as possible without going over. This routine always performs at least one KMC step before terminating.

  • Determine the time step for the next process

  • If the time limit is not exceeded, update clocks, rates, execute process, etc.; otherwise, abort.

Ideally we would use state(seed_size) but that was not working, so hardcoded size.

t : Requested simulation time increment (input) n : Maximum number of steps to run (input) num_iter : the number of executed iterations (output)

proclist/get_next_kmc_step

Determines next step without executing it. However, it changes the position in the random_number sequence. The python function for model.get_next_kmc_step() should be used as it makes additional function calls to reset the random numbers. Calling model.proclist.get_next_kmc_step() is discouraged as that will call this subroutine directly and will not reset the random numbers.

none

proclist/get_occupation

Evaluate current lattice configuration and returns the normalized occupation as matrix. Different species run along the first axis and different sites run along the second.

none

proclist/get_seed
Function to retrieve the state of the random number generator to

permit reproducible restart trajectories.

  • None

proclist/init

Allocates the system and initializes all sites in the given layer.

  • input_system_size number of unit cell per axis.

  • system_name identifier for reload file.

  • layer initial layer.

  • no_banner [optional] if True no copyright is issued.

proclist/initialize_state

Initialize all sites and book-keeping array for the given layer.

  • layer integer representing layer

proclist/put_seed

Subroutine to set the state of the random number generator to permit reproducible restart trajectories.

  • state an array of integers with the state of the random number generator (input)

proclist/run_proc_nr

Runs process proc on site nr_site.

  • proc integer representing the process number

  • nr_site integer representing the site

proclist/seed_gen

Function to transform a single number into a full set of integers required for initializing the random number generator.

  • sd an integer used to seed a simple random number generator used to generate additional integers for seeding the production random number generator (input)

kmos3/kind_values

This module offers kind_values for commonly used intrinsic types in a platform independent way.

lat_int

kmos3/base

The base kMC module, which implements the kMC method on a d = 1 lattice. Virtually any lattice kMC model can be build on top of this. The methods offered are:

  • de/allocation of memory

  • book-keeping of the lattice configuration and all available processes

  • updating and tracking kMC time, kMC step and wall time

  • saving and reloading the current state

  • determine the process and site to be executed

base/accum_rates

Stores the accumulated rate constant multiplied with the number of sites available for that process to be used by determine_procsite. Let \mathbf{c} be the rate constants \mathbf{n} the number of available sites, and \mathbf{a} the accumulated rates, then a_{i} is calculated according to a_{i}=\sum_{j=1}^{i} c_{j} n_{j}.

base/add_proc

The main idea of this subroutine is described in del_proc. Adding one process to one capability is programmatically simpler since we can just add it to the end of the respective array in avail_sites.

  • proc positive integer number that represents the process to be added.

  • site positive integer number that represents the site to be manipulated

base/allocate_system

Allocates all book-keeping structures and stores local copies of system name and size(s):

  • systen_name identifier of this simulation, used as name of punch file

  • volume the total number of sites

  • nr_of_proc the total number of processes

base/assertion_fail

Function that shall be used by all parts of the program to print a proper message in case some assertion fails.

  • a condition that is supposed to hold true

  • r message that is printed to the poor user in case it fails

base/avail_sites

Main book-keeping array that stores for each process the sites that are available and for each site the address in this very array. The meaning of the fields are:

avail_sites(proc, field, switch)

where:

  • proc – refers to a process in the process list

  • the field within the process, but the meaning differs as explained under ‘switch’

  • switch – can be either 1 or 2 and switches between (1) the actual numbers of the sites, which are available and filled in from the left but in whatever order they come or (2) the location where the site is stored in (1).

base/can_do

Returns true if ‘site’ can do ‘proc’ right now

  • proc integer representing the requested process.

  • site integer representing the requested site.

  • can writeable boolean, where the result will be stored.

base/deallocate_system

Deallocate all allocatable arrays: avail_sites, lattice, rates, accum_rates, procstat.

none

base/del_proc

del_proc delete one process from the main book-keeping array avail_sites. These book-keeping operations happen in O(1) time with the help of some more book-keeping overhead. avail_sites stores for each process all sites that are available. The array for each process is filled from the left, but sites generally not ordered. With this determine_procsite can effectively pick the next site and process. On the other hand a second array (avail_sites(:,:,2) ) holds for each process and each site, the location where it is stored in avail_site(:,:,1). If a site needs to be removed this subroutine first looks up the location via avail_sites(:,:,1) and replaces it with the site that is stored as the last element for this process.

  • proc positive integer that states the process

  • site positive integer that encodes the site to be manipulated

base/determine_procsite

Expects two random numbers between 0 and 1 and determines the corresponding process and site from accum_rates and avail_sites. Technically one random number would be sufficient but to circumvent issues with wrong interval_search_real implementation or rounding errors I decided to take two random numbers:

  • ran_proc Random real number from \in[0,1] that selects the next process

  • ran_site Random real number from \in[0,1] that selects the next site

  • proc Return integer \in[1,\mathrm{nr\_of\_proc}

  • site Return integer \in [1,\mathrm{volume}

base/get_accum_rate

Return accumulated rate at a given process.

  • proc_nr integer representing the requested process.

  • return_accum_rate writeable real, where the requested accumulated rate will be stored.

base/get_avail_site

Return field from the avail_sites database

  • proc_nr integer representing the requested process.

  • field integer for the site at question

  • switch 1 or 2 for site or storage location

base/get_integ_rate

Return integrated rate at a given process.

  • proc_nr integer representing the requested process.

  • return_integ_rate writeable real, where the requested integrated rate will be stored.

base/get_kmc_step

Return the current kmc_step

  • kmc_step Writeable integer

base/get_kmc_time

Returns current kmc_time as rdouble real as defined in kind_values.f90.

  • return_kmc_time writeable real, where the kmc_time will be stored.

base/get_kmc_time_step

Returns current kmc_time_step (the time increment).

  • return_kmc_step writeable integer, where the kmc_time_step will be stored.

base/get_kmc_volume

Return the total number of sites.

  • volume Writeable integer.

base/get_nrofsites

Return how many sites are available for a certain process. Usually used for debugging

  • proc integer representing the requested process

  • return_nrofsites writeable integer, where nr of sites gets stored

base/get_procstat

Return process counter for process proc as integer.

  • proc integer representing the requested process.

  • return_procstat writeable integer, where the process counter will be stored.

base/get_rate

Return rate of given process.

  • proc_nr integer representing the requested process.

  • return_rate writeable real, where the requested rate will be stored.

base/get_species

Return the species that occupies site.

  • site integer representing the site

base/get_system_name

Return the systems name, that was specified with base/allocate_system

  • system_name Writeable string of type character(len=200).

base/get_walltime

Return the current walltime.

  • return_walltime writeable real where the walltime will be stored.

base/increment_procstat

Increment the process counter for process proc by one.

  • proc integer representing the process to be increment.

base/integ_rates

Stores the time-integrated rates (non-normalized to surface area) Used to determine reaction rates, i.e. average number of reactions per unit surface and time. Let \mathbf{a} the integrated rates, \mathbf{c} be the rate constants, \mathbf{n}_i the number of available sites during kMC-time interval i, \{\Delta t_i\} the corresponding timesteps then a_{i}(t) at the time t=\sum_{i=1}\Delta t_i is calculated according to a_{i}(t)=\sum_{i=1} c_{i} n_{i}\Delta t_i.

base/interval_search_real

This is basically a standard binary search algorithm that expects an array of ascending real numbers and a scalar real and return the key of the corresponding field, with the following modification :

  • the value of the returned field is equal of larger of the given value. This is important because the given value is between 0 and the largest value in the array and otherwise the last field is never selected.

  • if two or more values in the array are identical, the function return the index of the leftmost of those field. This is important because having field with identical values means that all field except the leftmost one do not contain any sites. Refer to update_accum_rate to understand why.

  • the value of the returned field may no be zero. Therefore the index the to be equal or larger than the first non-zero field.

However: as everyone knows the binary search is trickier than it appears at first site especially real numbers. So intensive testing is suggested here!

  • arr real array of type rsingle (kind_values.f90) in monotonically (not strictly) increasing order

  • value real positive number from [0, max_arr_value]

base/kmc_step

Number of kMC steps executed.

base/kmc_time

Simulated kMC time in this run in seconds.

base/kmc_time_step

The time increment of the current kMC step.

base/lattice

Stores the actual physical lattice in a 1d array, where the value on each slot represents the species on that site.

Species constants can be conveniently defined in lattice_… and later used directly in the process list.

base/nr_of_proc

Total number of available processes.

base/nr_of_sites

Stores the number of sites available for each process.

base/procstat

Stores the total number of times each process has been executed during one simulation.

base/rates

Stores the rate constants for each process in s^-1.

base/reload_system

Restore state of simulation from *.reload file as saved by save_system(). This function also allocates the system’s memory so calling allocate_system again, will cause a runtime failure.

  • system_name string of 200 characters which will make the reload_system look for a file called ./<system_name>.reload

  • reloaded logical return variable, that is .true. reload of system could be completed successfully, and .false. otherwise.

base/replace_species

Replaces the species at a given site with new_species, given that old_species is correct, i.e. identical to the site that is already there.

  • site integer representing the site

  • old_species integer representing the species to be removed

  • new_species integer representing the species to be placed

base/reset_site

This function is a higher-level function to reset a site as if it never existed. To achieve this the species is set to null_species and all available processes are stripped from the site via del_proc.

  • site integer representing the requested site.

  • species integer representing the species that ought to be at the site, for consistency checks

base/save_system

save_system stores the entire system information in a simple ASCII filed names <system_name>.reload. All fields except avail_sites are stored in the simple scheme:

variable value

In the case of array variables, multiple values are seperated by one or more spaces, and the record is terminated with a newline. The variable avail_sites is treated slightly differently, since printed on a single line it is almost impossible to interpret from the ASCII files. Instead each process starts a new line, and the first number on the line stands for the process number and the remaining fields, hold the values.

none

base/set_kmc_time

Sets current kmc_time as rdouble real as defined in kind_values.f90.

  • new readable real, that the kmc time will be set to

base/set_rate_const

Allows to set the rate constant of the process with the number proc_nr.

  • proc_n The process number as defined in the corresponding proclist_ module.

  • rate the rate in s^{-1}

base/set_system_name

Set the systems name. Useful in conjunction with base.save_system to save *.reload files under a different name than the default one.

  • system_name Readable string of type character(len=200).

base/start_time

CPU time spent in simulation at least reload.

base/system_name

Unique indentifier of this simulation to be used for restart files. This name should not contain any characters that you don’t want to have in a filename either, i.e. only [A-Za-z0-9_-].

base/update_accum_rate

Updates the vector of accum_rates.

none

base/update_clocks

Updates walltime, kmc_step and kmc_time.

  • ran_time Random real number \in [0,1]

base/update_integ_rate

Updates the vector of integ_rates.

none

base/volume

Total number of sites.

base/walltime

Total CPU time spent on this simulation.

kmos3/lattice

Implements the mappings between the real space lattice and the 1-D lattice, which kmos3/base operates on. Furthermore replicates all geometry specific functions of kmos3/base in terms of lattice coordinates. Using this module each site can be addressed with 4-tuple (i, j, k, n) where i, j, k define the unit cell and n the site within the unit cell.

lattice/allocate_system

Allocates system, fills mapping cache, and checks whether mapping is consistent

none

lattice/calculate_lattice2nr

Maps all lattice coordinates onto a continuous set of integer \in [1,volume]

  • site integer array of size (4) a lattice coordinate

lattice/calculate_nr2lattice

Maps a continuous set of of integers \in [1,volume] to a 4-tuple representing a lattice coordinate

  • nr integer representing the site index

lattice/deallocate_system

Deallocates system including mapping cache.

none

lattice/default_layer

The layer in which the model is initially in by default (only relevant for multi-lattice models).

lattice/lattice2nr

Caching array holding the mapping from index to lattice coordinate: (x, y, z, n) -> i.

lattice/model_dimension

Store the number of dimensions of this model: 1, 2, or 3

lattice/nr2lattice

Caching array holding the mapping from index to lattice coordinate: i -> (x, y, z, n).

lattice/nr_of_layers

Constant storing the number of layers (for multi-lattice models > 1)

lattice/site_positions

The positions of (adsorption) site in the unit cell in fractional coordinates.

lattice/spuck

spuck = Sites Per Unit Cell Konstant The number of sites per unit cell, i.e. for coordinate notation (x, y, n) this is the maximum value of n.

lattice/system_size

Stores the current size of the allocated system lattice (x, y, z) in an integer array. In low-dimensional system, corresponding entries will be set to 1. Note that this should be thought of as a read-only variable. Changing its value at model runtime will not the indented effect of actually changing the simulated lattice. The definitive location for custom lattice size is simulation_size in kmc_settings.py.

If the system size shall be changed programmatically, it needs to happen before the KMC_Model is instantiated and Fortran array are allocated accordingly, like to

#!/usr/bin/env python3

import kmc_settings import kmos3.run

kmc_settings.simulation_size = 9, 9, 4

with kmos3.run.KMC_Model() as model:

print(model.lattice.system_size)))`

lattice/unit_cell_size

The dimensions of the unit cell (e.g. in Angstrom) of the unit cell.

kmos3/proclist

Implements the kMC process list.

proclist/do_kmc_step

Performs exactly one kMC step. * first update clock * then configuration sampling step * last execute process

none

proclist/do_kmc_steps

Performs n kMC step. If one has to run many steps without evaluation do_kmc_steps might perform a little better.

  • first update clock

  • then configuration sampling step

  • last execute process

n : Number of steps to run

proclist/do_kmc_steps_time

Performs a variable number of KMC steps to try to match the requested simulation time as closely as possible without going over. This routine always performs at least one KMC step before terminating.

  • Determine the time step for the next process

  • If the time limit is not exceeded, update clocks, rates, execute process, etc.; otherwise, abort.

Ideally we would use state(seed_size) but that was not working, so hardcoded size.

t : Requested simulation time increment (input) n : Maximum number of steps to run (input) num_iter : the number of executed iterations (output)

proclist/get_next_kmc_step

Determines next step without executing it. However, it changes the position in the random_number sequence. The python function for model.get_next_kmc_step() should be used as it makes additional function calls to reset the random numbers. Calling model.proclist.get_next_kmc_step() is discouraged as that will call this subroutine directly and will not reset the random numbers.

none

proclist/get_occupation

Evaluate current lattice configuration and returns the normalized occupation as matrix. Different species run along the first axis and different sites run along the second.

none

proclist/get_seed
Function to retrieve the state of the random number generator to

permit reproducible restart trajectories.

  • None

proclist/init

Allocates the system and initializes all sites in the given layer.

  • input_system_size number of unit cell per axis.

  • system_name identifier for reload file.

  • layer initial layer.

  • no_banner [optional] if True no copyright is issued.

proclist/initialize_state

Initialize all sites and book-keeping array for the given layer.

  • layer integer representing layer

proclist/put_seed

Subroutine to set the state of the random number generator to permit reproducible restart trajectories.

  • state an array of integers with the state of the random number generator (input)

proclist/seed_gen

Function to transform a single number into a full set of integers required for initializing the random number generator.

  • sd an integer used to seed a simple random number generator used to generate additional integers for seeding the production random number generator (input)

kmos3/kind_values

This module offers kind_values for commonly used intrinsic types in a platform independent way.

otf

kmos3/base

The base kMC module, which implements the kMC method on a d = 1 lattice. Virtually any lattice kMC model can be build on top of this. The methods offered are:

  • de/allocation of memory

  • book-keeping of the lattice configuration and all available processes

  • updating and tracking kMC time, kMC step and wall time

  • saving and reloading the current state

  • determine the process and site to be executed

base/accum_rates

Stores the accumulated rate constant up to a given process number taking into account all sites in which it is possible. ###

base/accum_rates_proc

Used to store the accumulated rate associated to each process ###

base/add_proc

The main idea of this subroutine is described in del_proc. Adding one process to one capability is programmatically simpler since we can just add it to the end of the respective array in avail_sites.

  • proc positive integer number that represents the process to be added.

  • site positive integer number that represents the site to be manipulated

base/allocate_system

Allocates all book-keeping structures and stores local copies of system name and size(s):

  • systen_name identifier of this simulation, used as name of punch file

  • volume the total number of sites

  • nr_of_proc the total number of processes

base/assertion_fail

Function that shall be used by all parts of the program to print a proper message in case some assertion fails.

  • a condition that is supposed to hold true

  • r message that is printed to the poor user in case it fails

base/avail_sites

Main book-keeping array that stores for each process the sites that are available and for each site the address in this very array. The meaning of the fields are:

avail_sites(proc, field, switch)

where:

  • proc – refers to a process in the process list

  • the field within the process, but the meaning differs as explained under ‘switch’

  • switch – can be either 1 or 2 and switches between (1) the actual numbers of the sites, which are available and filled in from the left but in whatever order they come or (2) the location where the site is stored in (1).

base/can_do

Returns true if ‘site’ can do ‘proc’ right now

  • proc integer representing the requested process.

  • site integer representing the requested site.

  • can writeable boolean, where the result will be stored.

base/deallocate_system

Deallocate all allocatable arrays: avail_sites, lattice, rates, accum_rates, procstat.

none

base/debug

If > 0, prints information

base/del_proc

del_proc delete one process from the main book-keeping array avail_sites. These book-keeping operations happen in O(1) time with the help of some more book-keeping overhead. avail_sites stores for each process all sites that are available. The array for each process is filled from the left, but sites generally not ordered. With this determine_procsite can effectively pick the next site and process. On the other hand a second array (avail_sites(:,:,2) ) holds for each process and each site, the location where it is stored in avail_site(:,:,1). If a site needs to be removed this subroutine first looks up the location via avail_sites(:,:,1) and replaces it with the site that is stored as the last element for this process.

  • proc positive integer that states the process

  • site positive integer that encodes the site to be manipulated

base/determine_procsite

Expects two random numbers between 0 and 1 and determines the corresponding process and site from accum_rates and avail_sites. Technically one random number would be sufficient but to circumvent issues with wrong interval_search_real implementation or rounding errors I decided to take two random numbers:

  • ran_proc Random real number from \\in[0,1] that selects the next process

  • ran_site Random real number from \\in[0,1] that selects the next site

  • proc Return integer \\in[1,\\mathrm{nr\\_of\\_proc}

  • site Return integer \\in [1,\\mathrm{volume}

base/get_accum_rate

Return accumulated rate at a given process.

  • proc_nr integer representing the requested process.

  • return_accum_rate writeable real, where the requested accumulated rate will be stored.

base/get_avail_site

Return field from the avail_sites database

  • proc_nr integer representing the requested process.

  • field integer for the site at question

  • switch 1 or 2 for site or storage location

base/get_integ_rate

Return integrated rate at a given process.

  • proc_nr integer representing the requested process.

  • return_integ_rate writeable real, where the requested integrated rate will be stored.

base/get_kmc_step

Return the current kmc_step

  • kmc_step Writeable integer

base/get_kmc_time

Returns current kmc_time as rdouble real as defined in kind_values.f90.

  • return_kmc_time writeable real, where the kmc_time will be stored.

base/get_kmc_time_step

Returns current kmc_time_step (the time increment).

  • return_kmc_step writeable integer, where the kmc_time_step will be stored.

base/get_kmc_volume

Return the total number of sites.

  • volume Writeable integer.

base/get_nrofsites

Return how many sites are available for a certain process. Usually used for debugging

  • proc integer representing the requested process

  • return_nrofsites writeable integer, where nr of sites gets stored

base/get_procstat

Return process counter for process proc as integer.

  • proc integer representing the requested process.

  • return_procstat writeable integer, where the process counter will be stored.

base/get_rate

Return rate of given process.

  • proc_nr integer representing the requested process.

  • return_rate writeable real, where the requested rate will be stored.

base/get_rate

Return rate constant of given process.

  • proc_nr integer representing the requested process.

  • return_rate writeable real, where the requested rate will be stored.

base/get_species

Return the species that occupies site.

  • site integer representing the site

base/get_system_name

Return the systems name, that was specified with base/allocate_system

  • system_name Writeable string of type character(len=200).

base/get_walltime

Return the current walltime.

  • return_walltime writeable real where the walltime will be stored.

base/increment_procstat

Increment the process counter for process proc by one.

  • proc integer representing the process to be increment.

base/integ_rates

Stores the time-integrated rates (non-normalized to surface area) Used to determine reaction rates, i.e. average number of reactions per unit surface and time. Let \\mathbf{a} the integrated rates, \\mathbf{c} be the rate constants, \\mathbf{n}_i the number of available sites during kMC-time interval i, \\{\\Delta t_i\\} the corresponding timesteps then a_{i}(t) at the time t=\\sum_{i=1}\\Delta t_i is calculated according to a_{i}(t)=\\sum_{i=1} c_{i} n_{i}\\Delta t_i.

base/interval_search_real

This is basically a standard binary search algorithm that expects an array of ascending real numbers and a scalar real and return the key of the corresponding field, with the following modification :

  • the value of the returned field is equal or larger than given value. This is important because the given value is between 0 and the largest value in the array and otherwise the last field is never selected.

  • if two or more values in the array are identical, the function return the index of the leftmost of those field. This is important because having field with identical values means that all field except the leftmost one do not contain any sites. Refer to update_accum_rate to understand why.

  • the value of the returned field may not be zero. Therefore the index the to be equal or larger than the first non-zero field.

However: as everyone knows the binary search is trickier than it appears at first sight especially real numbers. So intensive testing is suggested here!

  • arr real array of type rsingle (kind_values.f90) in monotonically (not strictly) increasing order

  • value real positive number from [0, max_arr_value]

base/kmc_step

Number of kMC steps executed.

base/kmc_time

Simulated kMC time in this run in seconds.

base/kmc_time_step

The time increment of the current kMC step.

base/lattice

Stores the actual physical lattice in a 1d array, where the value on each slot represents the species on that site.

Species constants can be conveniently defined in lattice\_… and later used directly in the process list.

base/nr_of_proc

Total number of available processes.

base/nr_of_sites

Stores the number of sites available for each process.

base/procstat

Stores the total number of times each process has been executed during one simulation.

base/rates

Stores the rate constants for each currently possible process ordered as avail_sites(:,:,1).

base/rates

Stores the rate constants for each process in s^-1.

base/reaccumulate_rates_matrix

Performs a process wide reaccumulation of the values in the rates_matrix. To be used when some of the user parameters are updated. Expected to aleviate some of the problems arising from floating point errors

base/reaccumulate_rates_matrix_proc

Performs a reaccumulation of the values of proc in the rates_matrix. To be used when some of the user parameters are updated. Expected to aleviate some of the problems arising from floating point errors

base/reload_system

Restore state of simulation from \*.reload file as saved by save_system(). This function also allocates the system’s memory so calling allocate_system again, will cause a runtime failure.

  • system_name string of 200 characters which will make the reload_system look for a file called ./<system_name>.reload

  • reloaded logical return variable, that is .true. reload of system could be completed successfully, and .false. otherwise.

base/replace_species

Replaces the species at a given site with new_species, given that old_species is correct, i.e. identical to the site that is already there.

  • site integer representing the site

  • old_species integer representing the species to be removed

  • new_species integer representing the species to be placed

base/reset_site

This function is a higher-level function to reset a site as if it never existed. To achieve this the species is set to null_species and all available processes are stripped from the site via del_proc.

  • site integer representing the requested site.

  • species integer representing the species that ought to be at the site, for consistency checks

base/save_system

save_system stores the entire system information in a simple ASCII filed names <system_name>.reload. All fields except avail_sites are stored in the simple scheme:

variable value

In the case of array variables, multiple values are seperated by one or more spaces, and the record is terminated with a newline. The variable avail_sites is treated slightly differently, since printed on a single line it is almost impossible to interpret from the ASCII files. Instead each process starts a new line, and the first number on the line stands for the process number and the remaining fields, hold the values.

none

base/set_kmc_time

Sets current kmc_time as rdouble real as defined in kind_values.f90.

  • new readable real, that the kmc time will be set to

base/set_rate_const

Allows to set the rate constant of the process with the number proc_nr.

  • proc_n The process number as defined in the corresponding proclist\_ module.

  • rate the rate in s^{-1}

base/set_system_name

Set the systems name. Useful in conjunction with base.save_system to save \*.reload files under a different name than the default one.

  • system_name Readable string of type character(len=200).

base/start_time

CPU time spent in simulation at least reload.

base/system_name

Unique indentifier of this simulation to be used for restart files. This name should not contain any characters that you don’t want to have in a filename either, i.e. only [A-Za-z0-9\_-].

base/update_accum_rate

Updates the vector of accum_rates.

none

base/update_clocks

Updates walltime, kmc_step and kmc_time.

  • ran_time Random real number \\in [0,1]

base/update_integ_rate

Updates the vector of integ_rates.

none

base/update_rates_matrix

Updates the rates_matrix. To be used when the state of a bystander has been modified

  • proc positive integer number that represents the process whose rate is changed.

  • site positive integer number that represents the site for the process

  • rate positive real number that represents the updated rate

base/volume

Total number of sites.

base/walltime

Total CPU time spent on this simulation.

kmos3/lattice

Implements the mappings between the real space lattice and the 1-D lattice, which kmos3/base operates on. Furthermore replicates all geometry specific functions of kmos3/base in terms of lattice coordinates. Using this module each site can be addressed with 4-tuple (i, j, k, n) where i, j, k define the unit cell and n the site within the unit cell.

lattice/allocate_system

Allocates system, fills mapping cache, and checks whether mapping is consistent

none

lattice/calculate_lattice2nr

Maps all lattice coordinates onto a continuous set of integer \in [1,volume]

  • site integer array of size (4) a lattice coordinate

lattice/calculate_nr2lattice

Maps a continuous set of of integers \in [1,volume] to a 4-tuple representing a lattice coordinate

  • nr integer representing the site index

lattice/deallocate_system

Deallocates system including mapping cache.

none

lattice/default_layer

The layer in which the model is initially in by default (only relevant for multi-lattice models).

lattice/lattice2nr

Caching array holding the mapping from index to lattice coordinate: (x, y, z, n) -> i.

lattice/model_dimension

Store the number of dimensions of this model: 1, 2, or 3

lattice/nr2lattice

Caching array holding the mapping from index to lattice coordinate: i -> (x, y, z, n).

lattice/nr_of_layers

Constant storing the number of layers (for multi-lattice models > 1)

lattice/site_positions

The positions of (adsorption) site in the unit cell in fractional coordinates.

lattice/spuck

spuck = Sites Per Unit Cell Konstant The number of sites per unit cell, i.e. for coordinate notation (x, y, n) this is the maximum value of n.

lattice/system_size

Stores the current size of the allocated system lattice (x, y, z) in an integer array. In low-dimensional system, corresponding entries will be set to 1. Note that this should be thought of as a read-only variable. Changing its value at model runtime will not the indented effect of actually changing the simulated lattice. The definitive location for custom lattice size is simulation_size in kmc_settings.py.

If the system size shall be changed programmatically, it needs to happen before the KMC_Model is instantiated and Fortran array are allocated accordingly, like to

#!/usr/bin/env python3

import kmc_settings import kmos3.run

kmc_settings.simulation_size = 9, 9, 4

with kmos3.run.KMC_Model() as model:

print(model.lattice.system_size)))`

lattice/unit_cell_size

The dimensions of the unit cell (e.g. in Angstrom) of the unit cell.

kmos3/proclist

Implements the kMC process list.

proclist/do_kmc_step

Performs exactly one kMC step. * first update clock * then configuration sampling step * last execute process

none

proclist/do_kmc_steps

Performs n kMC step. If one has to run many steps without evaluation do_kmc_steps might perform a little better.

  • first update clock

  • then configuration sampling step

  • last execute process

n : Number of steps to run

proclist/do_kmc_steps_time

Performs a variable number of KMC steps to try to match the requested simulation time as closely as possible without going over. This routine always performs at least one KMC step before terminating.

  • Determine the time step for the next process

  • If the time limit is not exceeded, update clocks, rates, execute process, etc.; otherwise, abort.

Ideally we would use state(seed_size) but that was not working, so hardcoded size.

t : Requested simulation time increment (input) n : Maximum number of steps to run (input) num_iter : the number of executed iterations (output)

proclist/get_next_kmc_step

Determines next step without executing it. However, it changes the position in the random_number sequence. The python function for model.get_next_kmc_step() should be used as it makes additional function calls to reset the random numbers. Calling model.proclist.get_next_kmc_step() is discouraged as that will call this subroutine directly and will not reset the random numbers.

none

proclist/get_occupation

Evaluate current lattice configuration and returns the normalized occupation as matrix. Different species run along the first axis and different sites run along the second.

none

proclist/get_seed
Function to retrieve the state of the random number generator to

permit reproducible restart trajectories.

  • None

proclist/init

Allocates the system and initializes all sites in the given layer.

  • input_system_size number of unit cell per axis.

  • system_name identifier for reload file.

  • layer initial layer.

  • no_banner [optional] if True no copyright is issued.

proclist/initialize_state

Initialize all sites and book-keeping array for the given layer.

  • layer integer representing layer

proclist/put_seed

Subroutine to set the state of the random number generator to permit reproducible restart trajectories.

  • state an array of integers with the state of the random number generator (input)

proclist/run_proc_nr

Runs process proc on site nr_site.

  • proc integer representing the process number

  • nr_site integer representing the site

proclist/seed_gen

Function to transform a single number into a full set of integers required for initializing the random number generator.

  • sd an integer used to seed a simple random number generator used to generate additional integers for seeding the production random number generator (input)

kmos3/kind_values

This module offers kind_values for commonly used intrinsic types in a platform independent way.

Todo

Fix some functions’ docstrings.