Running the Model From Runfiles
Running the Model–The API Way
Normally, one uses python runfiles to do simulations. However, it is convenient to initially run
commands interactively for learning purposes. The simplest thing to do is to start the model from
within a compiled model directory using kmos3 shell
(or kmos3 run
). That command will
start an IPython shell, allowing one to skip the below commands within a Python script
#!/usr/bin/env python3
from kmos3.run import KMC_Model
model = KMC_Model()
and just interact directly with model. It is often a good idea to use
%logstart some_scriptname.py
as your first command in the IPython shell to save what you have typed for later use.
When using a runfile, the starting banner can be turned off by using
model = KMC_Model(print_rates=False, banner=False)
Now that you have got a model object, you can do some KMC steps
model.do_steps(100000)
which will run 100,000 KMC steps.
Let’s say you want to change the temperature and the CO partial pressure of the model you can type
model.parameters.T = 550
model.parameters.p_COgas = 0.5
and all rate constants are instantly updated. In order get a quick overview of the current settings you can issue, e.g.,
print(model.parameters)
print(model.rate_constants)
or just
print(model)
Now an instantiated and configured model has mainly two functions: run KMC steps and report its current configuration. To analyze the current state, use
atoms = model.get_atoms()
Note
If you want to fetch data from the current state without actually visualizing the geometry, you
can speed up the get_atoms()
call using
atoms = model.get_atoms(geometry=False)
This will return an ASE atoms object of the current state of the system, including additional piggy-backed data such as the occupation, the turnover frequency (TOF), as well as the simulated time and the number of KMC steps simulated
model.get_occupation_header()
atoms.occupation
model.get_tof_header()
atoms.tof_data
atoms.kmc_time
atoms.kmc_step
If one wants to know what the next KMC step will be and at which site, without executing the step, one can use
model.get_next_kmc_step()
These quantities are often sufficient when running and simulating a catalyst surface, but of course the model could be expanded to more observables. The Fortran modules base, lattice, and proclist are attributes of the model instance, so feel free to explore the model instance, e.g., using IPython and the tabulator key (<TAB>)
model.base.<TAB>
model.lattice.<TAB>
model.proclist.<TAB>
The occupation is a 2-dimensional array which contains the occupation for each surface site divided by the number of unit cells. The first slot denotes the species and the second slot denotes the surface site, i.e.,
occupation = model.get_atoms().occupation
occupation[species, site-1]
So given there is a hydrogen species in the model, the occupation of hydrogen across all site types can be accessed like
hydrogen_occupation = occupation[model.proclist.hydrogen]
To access the coverage of one surface site, we have to remember to subtract 1, when using the the built-in constants, like so
hollow_occupation = occupation[:, model.lattice.hollow-1]
Lastly it is important to call
model.deallocate()
once the simulation is finished as this frees the memory allocated by the Fortan modules. This is particularly necessary if you want to run more than one simulation in one script.
Generate Grids of Sampled Data
For some KMC applications you may need a large number of data points across a range of external parameters (phase diagrams, microkinetic models). For this case there is the convenience class ModelRunner to work with
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)
ScanKinetics().run(init_steps=1e8, sample_steps=1e8, cores=4)
This script generates data points over the specified range(s). The TemperatureParameter
generates uniform grids over 1/T and the PressureParameter
uniform grids over log(p). The
script can be run synchronously over many cores as long as the cores can access the same file
system. You have to test whether the equilibration steps before sampling (init_steps) as well as
the batch size (sample_steps) is sufficient.
Manipulating the Model Species at Runtime
To change species on the lattice at the start of a simulation or at any other time during the
simulation, one can use the put
command. There are several syntaxes to use the put
command
model.put(site=[x,y,z,n], model.proclist.<species>)
Where n
and <species>
are the site type and species, respectively. For example
model.put([0,0,0,model.lattice.ruo2_bridge], model.proclist.co)
model.put([0,0,0,"ruo2_bridge"], "model.proclist.co")
model.put([0,0,0,2], 1)
Note
The n
is has indexing starting from 1 (there is no 0 for n
), whereas the <species>
indexing starts at 0.
If changing many sites at once, the above command is quite inefficient, since each put
call
adjusts the book-keeping database. To avoid these database updates, you can use the _put
method, like so
model._put(...)
model._put(...)
...
model._adjust_database()
Note
After using _put
, one must remember to call _adjust_database()
before executing any next
steps. Otherwise, the database of available processes will not match the species leading to an
incorrect continuation of the KMC simulation and a probable crash after some steps.
Saving and Reloading the State of the Simulation
To set the whole configuration of the lattice once, you can save and load it with the following commands
model.dump_config("YourConfigurationName")
model.load_config("YourConfigurationName")
While it is not necessary for a regular user to know, those commands use the following internal commands as part of how they function
# saving the configuration uses:
config = model._get_configuration()
# loading configuration uses:
model._set_configuration(config)
model._adjust_database()
However, simply saving and loading the configuration will not allow you to exactly reproduce the simulation where it left off. To achieve this, you must also save and reload the state of the pseudo-random number generator
# This list can be saved as a pickle or in a text file.
PRNG_state = model.proclist.get_seed().tolist()
# This command takes the PRNG_state as a list and inputs it into the model instance.
model.proclist.put_seed(PRNG_state)
By saving both the configuration and the pseudo-random number generator state, one can continue a simulation on the same trajectory (provided that the parameters are set to the same values as when the simulation stopped). The snapshots module includes methods for saving and loading the configuration, pseudo-random number generator state, and parameters. A single command to save all aspects of the simulation and reload the simulation where it was left off is planned to be added to kmos3 and the tutorials.
Todo
Add this functionality and tutorial instructions.
Running Models in Parallel
Due to the global clock in KMC there seems to be no simple and efficient way to parallelize a KMC simulation. Kmos3 cannot parallelize a single system over processors. However, one can run several kmos3 instances in parallel, which might accelerate sampling or efficiently check for steady state conditions.
However, in many applications it is still useful to run several models separately at once, for example to scan some set of parameters. This kind of problem can be considered embarrassingly parallel since it requires no communication between the different runs.
This is made very simple through the multiprocessing module, which is part of the Python standard library. Then, besides kmos3, we need to import multiprocessing in the Python script
from multiprocessing import Process
from numpy import linspace
from kmos3.run import KMC_Model
and let’s say you wanted to scan a range of temperatures while keeping all other parameters constant. You first define a function that takes a set of temperatures and runs an independent simulation for each
def run_temperatures(temperatures):
for T in temperatures:
model = KMC_Model()
model.parameters.T = T
model.do_steps(100000)
# do some evaluation
model.deallocate()
In order to split our full range of input parameters, one can use a utility function
from kmos3.utils import split_sequence
All that is left to do, is to define the input parameters, split the list and start subprocesses for each sub-list
if __name__ == '__main__':
temperatures = linspace(300, 600, 50)
nproc = 8
for temperatures in split_sequence(temperatures, nproc):
p = Process(target=run_temperatures, args=(temperatures, ))
p.start()
by executing the Python script.