Reference¶
Command Line Interface (CLI)¶
- Entry point module for the command-line
interface. The kmcos executable should be on the program path, import this modules main function and run it.
To call kmcos command as you would from the shell, use
kmcos.cli.main('...')
Every command can be shortened as long as it is non-ambiguous, e.g.
kmcos ex <xml-file>
instead of
kmcos export <xml-file>
etc.
List of commands¶
kmcos benchmark- Run 1 mio. kMC steps on model in current directory and report runtime.
kmcos buildBuild kmc_model.so from *f90 files in the current directory.
- Additional Parameters ::
- -d/–debug
- Turn on assertion statements in F90 code
- -n/–no-compiler-optimization
- Do not send optimizing flags to compiler.
kmcos edit <xml-file> (deprecated)- Open the kmcos xml-file in a GUI to edit the model.
kmcos export <xml-file> [<export-path>]Take a kmcos xml-file and export all generated source code to the export-path. There try to build the kmc_model.so.
Additional Parameters
-s/--source-only Export source only and don't build binary -b/--backend (local_smart|lat_int) Choose backend. Default is "local_smart". lat_int is EXPERIMENTAL and not made for production, yet. -d/--debug Turn on assertion statements in F90 code. (Only active in compile step) --acf Build the modules base_acf.f90 and proclist_acf.f90. Default is false. This both modules contain functions to calculate ACF (autocorrelation function) and MSD (mean squared displacement). -n/--no-compiler-optimization Do not send optimizing flags to compiler.
kmcos help <command>- Print usage information for the given command.
kmcos help all- Display documentation for all commands.
kmcos import <xml-file>- Take a kmcos xml-file and open an ipython shell with the project_tree imported as kmc_model.
kmcos rebuildExport code and rebuild binary module from XML information included in kmc_settings.py in current directory.
- Additional Parameters ::
- -d/–debug
- Turn on assertion statements in F90 code
kmcos run- Open an interactive shell and create a KMC_Model in it
- run == shell
kmcos settings-export <xml-file> [<export-path>]- Take a kmcos xml-file and export kmc_settings.py to the export-path.
kmcos shell- Open an interactive shell and create a KMC_Model in it
- run == shell
kmcos version- Print version number and exit.
kmcos viewTake a kmc_model.so and kmc_settings.py in the same directory and start to simulate the model visually.
- Additional Parameters ::
- -v/–steps-per-frame <number>
- Number of steps per frame
kmcos xml- Print xml representation of model to stdout
Data Types¶
kmcos.types¶
Holds all the data models used in kmcos.
-
class
kmcos.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 overall structure is the following as is also displayed in the editor GUI.
Project:
- Meta - Parameters - Lattice(s) - Species - Processes
-
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.
- color (str) – Color of species in editor GUI (#ffffff hex-type specification).
- representation (str) – ase.atoms.Atoms constructor describing species geometry.
- tags (str) – Tags of species (space separated string).
-
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_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
-
validate_model()¶ Run various consistency and completeness test of the model to make sure we have a minimally complete model.
-
-
class
kmcos.types.Meta(*args, **kwargs)¶ Class holding the meta-information about the kMC project
-
class
kmcos.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’
-
class
kmcos.types.LayerList(**kwargs)¶ A list of layers
Parameters: - cell (np.array (3x3)) – Size of unit-cell.
- default_layer (str.) – name of default 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.
-
class
kmcos.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: [])
-
class
kmcos.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.
-
class
kmcos.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.
- color (str) – Color of species in editor GUI (#ffffff hex-type specification).
- representation (str) – ase.atoms.Atoms constructor describing species geometry.
- tags (str) – Tags of species (space separated string).
-
class
kmcos.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.
-
class
kmcos.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: - coord (Coord) – Relative Coord (generated by
LayerList.generate_coord()orLattice.generate_coord_set()). - species (str) – Name of species.
- coord (Coord) – Relative Coord (generated by
-
class
kmcos.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.
kmcos.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.
-
kmcos.io.export_source(project_tree, export_dir=None, code_generator=None, options=None, accelerated=False)¶ Export a kmcos 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.
export_source is the central feature of the kmcos approach. In order to generate different backend solvers, additional candidates of this methods could be implemented.
-
kmcos.io.export_xml(project_tree, filename=None)¶ Writes a project to an XML file.
-
class
kmcos.io.ProcListWriter(data, dir)¶ Write the different parts of Fortran 90 code needed to run a kMC model.
-
write_proclist(smart=True, code_generator='local_smart', accelerated=False)¶ Write the proclist.f90 module, i.e. the rules which make up the kMC process list.
-
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.
-
Runtime frontend¶
kmcos.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 kmcos.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
kmcos.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 kmcos.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)
-
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 – Steps to run model before sampling (.ie. to reach steady-state). (Default: 1e8) :type init_steps: int :param sample_steps: Number of steps to sample over (Default: 1e8) :type sample_steps: int :param cores: Number of parallel processes to launch. :type cores: int :param samples: Number of samples. Use more samples if precise coverages are needed (Default: 1). :type samples: int
-
-
class
kmcos.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.
-
class
kmcos.run.PressureParameter(*args, **kwargs)¶ Create a grid of p in [p_min, p_max] such that ln({p}) is a regular grid.
-
class
kmcos.run.TemperatureParameter(*args, **kwargs)¶ Create a grid of p in [T_min, T_max] such that ({T})**(-1) is a regular grid.
-
class
kmcos.run.LinearParameter(*args, **kwargs)¶ Create a regular grid between min and max.
-
class
kmcos.run.LogParameter(*args, **kwargs)¶ Create a log grid between 10^min and 10^max (like np.logspace)
kmcos.run.KMC_Model (typical usage model.___)¶
-
class
kmcos.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, 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
model._put([0,0,0,model.lattice.lattice_bridge], model.proclist.co]) # puts a CO molecule at the `bridge` site of the lower left unit cell model._put([1,0,0,model.lattice.lattice_bridge], model.proclist.co ]) # puts a CO molecule at the `bridge` site one to the right # ... many more 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, check the site_name in kmc_settings.py.
Ex: site_names = [‘simple_cubic_hollow’] set site name as ‘sg.model.lattice.simple_cubic_hollow’
Ex: site_names = [‘ruo2_bridge’, ‘ruo2_cus’, ‘ruo2_Burrowed’] set site name as ‘sg.model.lattice.ruo2_bridge’ or ‘sg.model.lattice.ruo2_cus’ or ‘sg.model.lattice.ruo2_Burrowed’
- To see all the available species names, check the species_tags in 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.
-
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_steps(n=10000, progress=False)¶ Propagate the model n steps.
Parameters: n (int) – Number of steps to run (Default: 10000)
-
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)¶ 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
-
get_atoms(geometry=True, tag=None, reset_time_overrun=False)¶ Return an ASE 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 ASE object of current configuration (Default: True).
-
get_backend()¶ Return name of backend that model was compiled with.
Return type: str
-
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_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_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.
-
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
-
plot_configuration(filename='', directory='./exported_configurations', resolution=150, scale=20, representation='spatial', plot_settings={})¶ 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.
-
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.
Examples
model.put([0,0,0,model.lattice.site], model.proclist.co ]) # puts a CO molecule at the `bridge` site # of the lower left unit cell
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, check the site_name in kmc_settings.py.
Ex: site_names = [‘simple_cubic_hollow’] set site name as ‘sg.model.lattice.simple_cubic_hollow’
Ex: site_names = [‘ruo2_bridge’, ‘ruo2_cus’, ‘ruo2_Burrowed’] set site name as ‘sg.model.lattice.ruo2_bridge’ or ‘sg.model.lattice.ruo2_cus’ or ‘sg.model.lattice.ruo2_Burrowed’
- To see all the available species names, check the species_tags in kmc_settings.py.
- Ex: species_tags = {
- “CO”:, “O”:, “empty”:, }
-
run()¶ Runs the model indefinitely. To control the simulations, model must have been initialized with proper Queues.
-
start()¶ Start child process
-
view(scaleA=None)¶ Start current model in live view mode.
-
xml()¶ Returns the XML representation that this model was created from.
Return type: str
-
-
class
kmcos.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 (kmcos 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.
-
-
class
kmcos.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
kmcos.view¶
kmcos.cli¶
Entry point module for the command-line interface. The kmcos executable should be on the program path, import this modules main function and run it.
To call kmcos command as you would from the shell, use
kmcos.cli.main('...')
Every command can be shortened as long as it is non-ambiguous, e.g.
kmcos ex <xml-file>
instead of
kmcos export <xml-file>
etc.
You may also use syntax kmcos.export(”…”) for any cli command.
-
kmcos.cli.main(args=None)¶ The CLI main entry point function.
The optional argument args, can be used to directly supply command line argument like
$ kmcos <args>
otherwise args will be taken from STDIN.
kmcos.utils¶
Several utility functions that do not seem to fit somewhere else.
-
kmcos.utils.build(options)¶ Build binary with f2py binding from complete set of source file in the current directory.
-
kmcos.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.
-
kmcos.utils.get_ase_constructor(atoms)¶ Return the ASE constructor string for atoms.
-
kmcos.utils.split_sequence(seq, size)¶ Take a list and a number n and return list divided into n sublists of roughly equal size.
-
kmcos.utils.write_py(fileobj, images, **kwargs)¶ Write a ASE atoms construction string for images into fileobj.
kmcos kMC project DTD¶
The central storage and exchange format is XML. XML was chosen over JSON, pickle or alike because it still seems as 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 kmcos file is validated against the DTD below.
<!ELEMENT kmc (meta?,species_list?,parameter_list?, lattice, process_list?,output_list?)>
<!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)*>
<!ATTLIST process
name CDATA #REQUIRED
rate_constant CDATA #REQUIRED
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 output_list (output)*>
<!ELEMENT output EMPTY>
<!ATTLIST output
item CDATA #REQUIRED
>
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 kmcos.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 kmcos.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¶
kmcos/base¶
The base kMC module, which implements the kMC method on a
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. Letbe the rate constants
the number of available sites, and
the accumulated rates, then
is calculated according to
.
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.
procpositive integer number that represents the process to be added.sitepositive 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_nameidentifier of this simulation, used as name of punch filevolumethe total number of sitesnr_of_procthe 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.
acondition that is supposed to hold truermessage 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
procinteger representing the requested process.siteinteger representing the requested site.canwriteable 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.
procpositive integer that states the processsitepositive 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_procRandom real number fromthat selects the next process
ran_siteRandom real number fromthat selects the next site
procReturn integersiteReturn integer
base/get_accum_rate¶
Return accumulated rate at a given process.
proc_nrinteger representing the requested process.return_accum_ratewriteable real, where the requested accumulated rate will be stored.
base/get_avail_site¶
Return field from the avail_sites database
proc_nrinteger representing the requested process.fieldinteger for the site at questionswitch1 or 2 for site or storage location
base/get_integ_rate¶
Return integrated rate at a given process.
proc_nrinteger representing the requested process.return_integ_ratewriteable real, where the requested integrated rate will be stored.
base/get_kmc_step¶
Return the current kmc_step
kmc_stepWriteable integer
base/get_kmc_time¶
Returns current kmc_time as rdouble real as defined in kind_values.f90.
return_kmc_timewriteable real, where the kmc_time will be stored.
base/get_kmc_time_step¶
Returns current kmc_time_step (the time increment).
return_kmc_stepwriteable real, where the kmc_time_step will be stored.
base/get_kmc_volume¶
Return the total number of sites.
volumeWriteable integer.
base/get_nrofsites¶
Return how many sites are available for a certain process. Usually used for debugging
procinteger representing the requested processreturn_nrofsiteswriteable integer, where nr of sites gets stored
base/get_procstat¶
Return process counter for process proc as integer.
procinteger representing the requested process.return_procstatwriteable integer, where the process counter will be stored.
base/get_rate¶
Return rate of given process.
proc_nrinteger representing the requested process.return_ratewriteable real, where the requested rate will be stored.
base/get_species¶
Return the species that occupies site.
siteinteger representing the site
base/get_system_name¶
Return the systems name, that was specified with base/allocate_system
system_nameWriteable string of type character(len=200).
base/get_walltime¶
Return the current walltime.
return_walltimewriteable real where the walltime will be stored.
base/increment_procstat¶
Increment the process counter for process proc by one.
procinteger 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. Letthe integrated rates,
be the rate constants,
the number of available sites during kMC-time interval i,
the corresponding timesteps then
at the time
is calculated according to
.
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!
arrreal array of type rsingle (kind_values.f90) in monotonically (not strictly) increasing ordervaluereal 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_namestring of 200 characters which will make the reload_system look for a file called ./<system_name>.reloadreloadedlogical 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.
siteinteger representing the siteold_speciesinteger representing the species to be removednew_speciesinteger 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.
siteinteger representing the requested site.speciesinteger 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 valueIn 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.
newreadable 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_nThe process number as defined in the corresponding proclist_ module.ratethe rate in
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_nameReadable 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_timeRandom real number
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.
kmcos/lattice¶
Implements the mappings between the real space lattice and the 1-D lattice, which kmcos/base operates on. Furthermore replicates all geometry specific functions of kmcos/base in terms of lattice coordinates. Using this module each site can be addressed with 4-tuple(i, j, k, n)wherei, j, kdefine the unit cell andnthe 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
siteinteger array of size (4) a lattice coordinate
lattice/calculate_nr2lattice¶
Maps a continuous set of of integers
to a 4-tuple representing a lattice coordinate
nrinteger 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 python
import kmc_settings import kmcos.run
kmc_settings.simulation_size = 9, 9, 4
- with kmcos.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.
kmcos/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
nkMC 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/get_kmc_step¶
Determines next step without executing it.
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/init¶
Allocates the system and initializes all sites in the given layer.
input_system_sizenumber of unit cell per axis.system_nameidentifier for reload file.layerinitial layer.no_banner[optional] if True no copyright is issued.
proclist/initialize_state¶
Initialize all sites and book-keeping array for the given layer.
layerinteger representing layer
proclist/run_proc_nr¶
Runs process
procon sitenr_site.
procinteger representing the process numbernr_siteinteger representing the site
libkmc/kind_values¶
This module offers kind_values for commonly used intrinsic types in a platform independent way.
lat_int¶
kmcos/base¶
The base kMC module, which implements the kMC method on a
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. Letbe the rate constants
the number of available sites, and
the accumulated rates, then
is calculated according to
.
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.
procpositive integer number that represents the process to be added.sitepositive 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_nameidentifier of this simulation, used as name of punch filevolumethe total number of sitesnr_of_procthe 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.
acondition that is supposed to hold truermessage 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
procinteger representing the requested process.siteinteger representing the requested site.canwriteable 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.
procpositive integer that states the processsitepositive 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_procRandom real number fromthat selects the next process
ran_siteRandom real number fromthat selects the next site
procReturn integersiteReturn integer
base/get_accum_rate¶
Return accumulated rate at a given process.
proc_nrinteger representing the requested process.return_accum_ratewriteable real, where the requested accumulated rate will be stored.
base/get_avail_site¶
Return field from the avail_sites database
proc_nrinteger representing the requested process.fieldinteger for the site at questionswitch1 or 2 for site or storage location
base/get_integ_rate¶
Return integrated rate at a given process.
proc_nrinteger representing the requested process.return_integ_ratewriteable real, where the requested integrated rate will be stored.
base/get_kmc_step¶
Return the current kmc_step
kmc_stepWriteable integer
base/get_kmc_time¶
Returns current kmc_time as rdouble real as defined in kind_values.f90.
return_kmc_timewriteable real, where the kmc_time will be stored.
base/get_kmc_time_step¶
Returns current kmc_time_step (the time increment).
return_kmc_stepwriteable integer, where the kmc_time_step will be stored.
base/get_kmc_volume¶
Return the total number of sites.
volumeWriteable integer.
base/get_nrofsites¶
Return how many sites are available for a certain process. Usually used for debugging
procinteger representing the requested processreturn_nrofsiteswriteable integer, where nr of sites gets stored
base/get_procstat¶
Return process counter for process proc as integer.
procinteger representing the requested process.return_procstatwriteable integer, where the process counter will be stored.
base/get_rate¶
Return rate of given process.
proc_nrinteger representing the requested process.return_ratewriteable real, where the requested rate will be stored.
base/get_species¶
Return the species that occupies site.
siteinteger representing the site
base/get_system_name¶
Return the systems name, that was specified with base/allocate_system
system_nameWriteable string of type character(len=200).
base/get_walltime¶
Return the current walltime.
return_walltimewriteable real where the walltime will be stored.
base/increment_procstat¶
Increment the process counter for process proc by one.
procinteger 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. Letthe integrated rates,
be the rate constants,
the number of available sites during kMC-time interval i,
the corresponding timesteps then
at the time
is calculated according to
.
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!
arrreal array of type rsingle (kind_values.f90) in monotonically (not strictly) increasing ordervaluereal 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_namestring of 200 characters which will make the reload_system look for a file called ./<system_name>.reloadreloadedlogical 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.
siteinteger representing the siteold_speciesinteger representing the species to be removednew_speciesinteger 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.
siteinteger representing the requested site.speciesinteger 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 valueIn 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.
newreadable 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_nThe process number as defined in the corresponding proclist_ module.ratethe rate in
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_nameReadable 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_timeRandom real number
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.
kmcos/lattice¶
Implements the mappings between the real space lattice and the 1-D lattice, which kmcos/base operates on. Furthermore replicates all geometry specific functions of kmcos/base in terms of lattice coordinates. Using this module each site can be addressed with 4-tuple(i, j, k, n)wherei, j, kdefine the unit cell andnthe 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
siteinteger array of size (4) a lattice coordinate
lattice/calculate_nr2lattice¶
Maps a continuous set of of integers
to a 4-tuple representing a lattice coordinate
nrinteger 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 python
import kmc_settings import kmcos.run
kmc_settings.simulation_size = 9, 9, 4
- with kmcos.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.
kmcos/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
nkMC 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/get_kmc_step¶
Determines next step without executing it.
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/init¶
Allocates the system and initializes all sites in the given layer.
input_system_sizenumber of unit cell per axis.system_nameidentifier for reload file.layerinitial layer.no_banner[optional] if True no copyright is issued.
proclist/initialize_state¶
Initialize all sites and book-keeping array for the given layer.
layerinteger representing layer
libkmc/kind_values¶
This module offers kind_values for commonly used intrinsic types in a platform independent way.
otf¶
kmcos/base¶
The base kMC module, which implements the kMC method on a
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.
procpositive integer number that represents the process to be added.sitepositive 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_nameidentifier of this simulation, used as name of punch filevolumethe total number of sitesnr_of_procthe 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.
acondition that is supposed to hold truermessage 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
procinteger representing the requested process.siteinteger representing the requested site.canwriteable 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.
procpositive integer that states the processsitepositive 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_procRandom real number fromthat selects the next process
ran_siteRandom real number fromthat selects the next site
procReturn integersiteReturn integer
base/get_accum_rate¶
Return accumulated rate at a given process.
proc_nrinteger representing the requested process.return_accum_ratewriteable real, where the requested accumulated rate will be stored.
base/get_avail_site¶
Return field from the avail_sites database
proc_nrinteger representing the requested process.fieldinteger for the site at questionswitch1 or 2 for site or storage location
base/get_integ_rate¶
Return integrated rate at a given process.
proc_nrinteger representing the requested process.return_integ_ratewriteable real, where the requested integrated rate will be stored.
base/get_kmc_step¶
Return the current kmc_step
kmc_stepWriteable integer
base/get_kmc_time¶
Returns current kmc_time as rdouble real as defined in kind_values.f90.
return_kmc_timewriteable real, where the kmc_time will be stored.
base/get_kmc_time_step¶
Returns current kmc_time_step (the time increment).
return_kmc_stepwriteable integer, where the kmc_time_step will be stored.
base/get_kmc_volume¶
Return the total number of sites.
volumeWriteable integer.
base/get_nrofsites¶
Return how many sites are available for a certain process. Usually used for debugging
procinteger representing the requested processreturn_nrofsiteswriteable integer, where nr of sites gets stored
base/get_procstat¶
Return process counter for process proc as integer.
procinteger representing the requested process.return_procstatwriteable integer, where the process counter will be stored.
base/get_rate¶
Return rate of given process.
proc_nrinteger representing the requested process.return_ratewriteable real, where the requested rate will be stored.
base/get_species¶
Return the species that occupies site.
siteinteger representing the site
base/get_system_name¶
Return the systems name, that was specified with base/allocate_system
system_nameWriteable string of type character(len=200).
base/get_walltime¶
Return the current walltime.
return_walltimewriteable real where the walltime will be stored.
base/increment_procstat¶
Increment the process counter for process proc by one.
procinteger 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. Letthe integrated rates,
be the rate constants,
the number of available sites during kMC-time interval i,
the corresponding timesteps then
at the time
is calculated according to .
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!
arrreal array of type rsingle (kind_values.f90) in monotonically (not strictly) increasing ordervaluereal 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/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_namestring of 200 characters which will make the reload_system look for a file called ./<system_name>.reloadreloadedlogical 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.
siteinteger representing the siteold_speciesinteger representing the species to be removednew_speciesinteger 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.
siteinteger representing the requested site.speciesinteger 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 valueIn 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.
newreadable 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_nThe process number as defined in the corresponding proclist_ module.ratethe rate in
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_nameReadable 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_timeRandom real number
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
!
procpositive integer number that represents the process whose rate is changed.sitepositive integer number that represents the site for the processratepositive real number that represents the updated rate
base/volume¶
Total number of sites.
base/walltime¶
Total CPU time spent on this simulation.
kmcos/lattice¶
Implements the mappings between the real space lattice and the 1-D lattice, which kmcos/base operates on. Furthermore replicates all geometry specific functions of kmcos/base in terms of lattice coordinates. Using this module each site can be addressed with 4-tuple(i, j, k, n)wherei, j, kdefine the unit cell andnthe 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
siteinteger array of size (4) a lattice coordinate
lattice/calculate_nr2lattice¶
Maps a continuous set of of integers
to a 4-tuple representing a lattice coordinate
nrinteger 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 python
import kmc_settings import kmcos.run
kmc_settings.simulation_size = 9, 9, 4
- with kmcos.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.
kmcos/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
nkMC 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/get_kmc_step¶
Determines next step without executing it.
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/init¶
Allocates the system and initializes all sites in the given layer.
input_system_sizenumber of unit cell per axis.system_nameidentifier for reload file.layerinitial layer.no_banner[optional] if True no copyright is issued.
proclist/initialize_state¶
Initialize all sites and book-keeping array for the given layer.
layerinteger representing layer
proclist/run_proc_nr¶
Runs process
procon sitenr_site.
procinteger representing the process numbernr_siteinteger representing the site
libkmc/kind_values¶
This module offers kind_values for commonly used intrinsic types in a platform independent way.
lattice. Virtually any lattice kMC model can be build on top of this.
The methods offered are:
be the rate constants
the number of available sites, and
the accumulated rates, then
is calculated according to
.
that selects the next process

the number of available sites
during kMC-time interval i,
the corresponding
timesteps then
at the time
is calculated according to
.
![\in [0,1]](../_images/math/0268e0a93e47bc316c3d24060f782557873cd7bd.png)
![\in [1,volume]](../_images/math/ea7b78ca3fd084d7c3ad123e8f96329c4706f363.png)