PAPRECA hybrid off-lattice kMC/MD simulator
2.0.0 (17 September 2024)
|
This page contains information regarding PAPRECA commands. For details related to LAMMPS commands please visit this documentation page.
Important note:: PAPRECA will always ignore text to the right of a "#" character. You can use "#" characters to add comments to your input files.
LAMMPS fix command to be used in the LAMMPS input file. This command has to be utilized after defining the LAMMPS simulation box in your LAMMPS input file. The papreca fix for LAMMPS initializes and updates two neighbors lists (a half and a full) and facilitates predefined event detection and execution. Certain kMC events (e.g., diffusion, deposition) require interference (i.e., collision) checks. Therefore, a full neighbor list has to be used with such events as it guarantees that the collision can be detected when scanning through the neighbor list of either atom. On the other hand, the discovery of other kMC events (e.g., bond formation) is significantly more efficient when half lists are used (since those lists each pair of atoms once). Effectively, this is why two neighbor lists are required for PAPRECA. Please refer to the LAMMPS documentation for more information regarding neighbor lists.
Note 1: If you plan to use special_bonds in your simulation refrain from setting ANY of the special_bonds to zero. Setting a special_bond to zero eliminates the (1-2,1-3, or 1-4) neighbors from the neighbor lists. Please use a double number beyond the accuracy limits of a C++ double instead of zero (e.g., use "special_bonds lj 1e-100 1.0 1.0 coul 1e-100 1.0 1.0" in your input file instead of "special_bonds lj 0.0 1.0 1.0 coul 0.0 1.0 1.0" to include the 1-2 neighbors). Once again, this does not affect the computational efficiency of the MD stage but includes additional neighbor pairs in the neighbor list.
Note 2: Consider using neigh_modify command with the "every 1", "delay 0", and "check yes" options in your LAMMPS input file (where applicable). This forces LAMMPS to update the neighbors lists on every timestep and increases the accuracy of your PAPRECA run (in the expense of computational efficiency).
Sets the total number of kMC stages for the PAPRECA run.
Note: This is a mandatory command. The PAPRECA simulation will not start unless the total number of kMC stages is set.
Sets the total number of kMC stages per MD stage for the PAPRECA run.
Note: This is a mandatory command. The PAPRECA simulation will not start unless the total number of kMC stages per MD stage is set.
Defines the random seed to initialize the random number generator. PAPRECA uses the following random number generator: RANMAR (the algorithm of Marsaglia, Zaman, Tsang)[1].
Note: This is a mandatory command. The PAPRECA simulation will not start unless the random seed is set.
[1] James, F. "A review of pseudorandom number generators." Computer Physics Communications, vol. 60, 1990
Defines the options related to the initialization of "sigma" (as in the Lennard-Jones potential, see here ) values for the PAPRECA simulation. Note that, for events requiring collision tests (e.g., deposition events) the sigma value defines the maximum acceptable distance between two species (i.e., if d ≤ sigma PAPRECA assumes that an interference exists and that atom collides) [1] , [2].
For style manual, the sigma values have to be manually declared in the PAPRECA input file (using sequential init_sigma command) for all types.
For style LAMMPS, the sigma values are retrieved from LAMMPS. Note that, it is necessary that your pair_style of choice uses sigma values (e.g., lj/cut/ and that all sigma values are properly set in the LAMMPS input file.
The mix keyword performs geometric or arithmetic mixing on already defined sigma values. The mixing is performed after the last init_sigma command line in the PAPRECA input (or right before the start of the PAPRECA run if the init_sigma command has not been used). Hence, it can be very useful if you wish to solely define the diagonal sigma terms. See this LAMMPS documentation page for more information regarding the mixing styles.
Note: The sigma values for all pairs of atom types (diagonal and cross-terms) have to be defined before the start of a PAPRECA run. The PAPRECA run will start normally even if the sigmas of all pairs of atom types have not been defined. However, the PAPRECA run will abort as soon as an interference check is performed for pairs of atom types of unknown sigma values.
[1] Ntioudis, S., et al. "PAPRECA: A parallel hybrid off-lattice kinetic Monte Carlo/molecular dynamics simulator", Journal of Open Source Software, 9(98), 6714 (2024). https://doi.org/10.21105/joss.06714
[2] Ntioudis, S., et al. "A hybrid off-lattice kinetic Monte Carlo/molecular dynamics method for amorphous thin film growth", Computational Materials Science, 229, 112421 (2023). https://doi.org/10.1016/j.commatsci.2023.112421
Defines a "sigma" value (as in the Lennard-Jones potential, see here ). Note that, for events requiring collision tests (e.g., deposition events) the sigma value defines the maximum acceptable distance between two species (i.e., if d &le sigma PAPRECA assumes that an interference exists and that atom collides) [1], [2].
Note 1: The sigma values for all pairs of atom types (diagonal and cross-terms) have to be defined before the start of a PAPRECA run. The PAPRECA run will start normally even if the sigmas of all pairs of atom types have not been defined. However, the PAPRECA run will abort as soon as an interference check is performed for pairs of atom types of unknown sigma values.
Note 2: A sigmas_options command has to be present in your PAPRECA input file before you use this command.
Note 3: This command can be used as many times as necessary in order to define all the sigmas for the PAPRECA simulation. When used multiple times, previously defined sigma values will be overwritten.
Note 4: PAPRECA will not check if the defined sigma corresponds to atom types that exist in the simulation (i.e., defined in the LAMMPS input file).
[1] Ntioudis, S., et al. "PAPRECA: A parallel hybrid off-lattice kinetic Monte Carlo/molecular dynamics simulator", Journal of Open Source Software, 9(98), 6714 (2024). https://doi.org/10.21105/joss.06714
[2] Ntioudis, S., et al. "A hybrid off-lattice kinetic Monte Carlo/molecular dynamics method for amorphous thin film growth", Computational Materials Science, 229, 112421 (2023). https://doi.org/10.1016/j.commatsci.2023.112421
Defines the mobile atom types for the MD stage of the PAPRECA simulation
Note 1: This is a mandatory command. The PAPRECA simulation will not start unless the fluid atom types have been defined.
Note 2: PAPRECA will not check if the defined mobile atom types are consistent with the mobile types of the MD (LAMMPS) stage.
Defines the frozen atom types for the MD stage PAPRECA simulation.
Note 1: PAPRECA will not check if the defined frozen atom types are consistent with the mobile types of the MD (LAMMPS) stage.
Note 2: This command is not mandatory. However, you must use this command if certain atom types are defined as frozen in your LAMMPS input file (e.g., if you used "fix freeze_atoms frozen setforce 0.0 0.0 0.0" ). See the LAMMPS documentation page for more information.
PAPRECA will start with an empty frozen_atomtypes list if the frozen-atomtypes command is not used.
Sets a desired ending time for the PAPRECA run, i.e., the PAPRECA run will stop when the simulation time is equal to or greater than time_end. This command does not replace kMC_steps command. It simply places a time limit on top of kMC_steps command.
time_end = std::numeric_limits< double >::max( ). This is the maximum limit of a double number in C++. Therefore, you will never reach that number.
Sets up a height calculation for the PAPRECA simulation. If this command is used, the height across the z-coordinate is calculated before every kMC stage. This command was designed to dynamically capture the height of thin films.
For style mass_bins the height calculation is performed as follows: firstly, the whole simulation box is divided into x-y segments. The x-y segments span from on side to another in the x- and y- directions and have a thickness (in the z-direction) of bin_width. Secondly, the total mass of each x-y segment is calculated. Finally, a running sum of bin masses is calculated (starting from the lowermost x-y segment) and the film height is defined as the height of the first mass bin for which the running mass sum is greater that or equal to M × cutoff (where M is the total mass of the system at the current PAPRECA step).
No height calculations are performed if the user does not include this command in the PAPRECA input file.
This command does NOT set up a predefined desorption kMC event template. This command was created to assist film growth studies and is designed to delete all atoms whose z-coordinate is equal to or greater than height (as set in the command).
For styles "gather_all" and "gather_local" results must be identical. However, the choice of style might affect the efficiency of the PAPRECA run due to different implementations. See the relevant C++ function documentation for more information about these 2 different approaches: deleteDesorbedAtoms(). Also, for styles "gather_all" and "gather_local" atoms bonded to deleted atoms are also deleted. For example, if the deletion height is set to 30 and an atom is above 30 (LAMMPS length units) and bonded to another atom whose z-coordinate is 29, then, both atoms will be deleted. Bonded atoms are deleted to prevent "bond atoms missing from proc %d" errors (see LAMMPS documentation page) within the MD stage of PAPRECA run. If the "max" keyword is used, then the maximum number of atoms that can be deleted at once becomes "N".
For style "LAMMPS_region" a wrapper (see deleteAtomsInBoxRegion()) around the region and delete_atoms commands are used to delete atoms above the deletion height. Note that, unlike the "gather_all" and "gather_local" styles, the "LAMMPS_region" style will simply delete all bonded interactions (i.e., bond, angles, dihedrals, and impropers) associated with the deleted atoms, but not the bonded atoms to deleted atoms. Consider This as it may lead to instabilities to do sudden system energy change.
Note: Prior to using this command you must set up a height calculation for your PAPRECA simulation (i.e., a height_calculation command must be present in your PAPRECA above the current command line).
No atom deletions (desorptions) are performed if the user does not include this command in the PAPRECA input file. Also, for styles "gather_all" and "gather_local", if the desorption command is used without the "max" keyword then any number of atoms can be deleted at once.
Defines a minimization command to be executed within the MD stage of the PAPRECA run and before simulating the MD trajectory. This command might be helpful to relax the system to the closest Potential Energy Surface (PES) valley and avoid instabilities [1], [2] within the MD trajectory (e.g., "bond atoms missing from proc %d" errors (see LAMMPS documentation page). Please see the relevant LAMMPS documentation page for the minimize command for more information.
When the "yes" keyword is utilized, the user has to provide a valid LAMMPS minimization command. Note that, PAPRECA will not check if the command is valid before the start of the simulation but will probably abort during runtime. Moreover, note that any valid LAMMPS command can be passed to PAPRECA with this command.
Note: A LAMMPS minimization style has to be defined in the LAMMPS input file (see min_style for more information regarding minimization styles).
No minimization prior to the MD trajectory is performed if this command is not included in the PAPRECA input file.
[1] Ntioudis, S., et al. "PAPRECA: A parallel hybrid off-lattice kinetic Monte Carlo/molecular dynamics simulator", Journal of Open Source Software, 9(98), 6714 (2024). https://doi.org/10.21105/joss.06714
[2] Ntioudis, S., et al. "A hybrid off-lattice kinetic Monte Carlo/molecular dynamics method for amorphous thin film growth", Computational Materials Science, 229, 112421 (2023). https://doi.org/10.1016/j.commatsci.2023.112421
Defines a minimization command to be executed within the MD stage of the PAPRECA run and right after simulating the MD trajectory. This command might be helpful to relax the system to the closest Potential Energy Surface (PES) valley and avoid instabilities [1], [2] within the MD trajectory (e.g., "bond atoms missing from proc %d" errors (see LAMMPS documentation page). Please see the relevant LAMMPS documentation page for the minimize command for more information.
When the "yes" keyword is utilized, the user has to provide a valid LAMMPS minimization command. Note that, PAPRECA will not check if the command is valid before the start of the simulation but will probably abort during runtime. Moreover, note that any valid LAMMPS command can be passed to PAPRECA with this command.
Note: A LAMMPS minimization style has to be defined in the LAMMPS input file (see min_style for more information regarding minimization styles).
No minimization after the MD trajectory is performed if this command is not included in the PAPRECA input file.
[1] Ntioudis, S., et al. "PAPRECA: A parallel hybrid off-lattice kinetic Monte Carlo/molecular dynamics simulator", Journal of Open Source Software, 9(98), 6714 (2024). https://doi.org/10.21105/joss.06714
[2] Ntioudis, S., et al. "A hybrid off-lattice kinetic Monte Carlo/molecular dynamics method for amorphous thin film growth", Computational Materials Science, 229, 112421 (2023). https://doi.org/10.1016/j.commatsci.2023.112421
Sets the duration of the MD stage of the PAPRECA run. This represents how many MD timesteps will be simulated by LAMMPS during each MD stage. The trajectory duration can be 0 (means that PAPRECA will perform a pure off-lattice kMC run).
trajectory_duration = 0.
Create a predefined bond-break template (see PAPRECA::PredefinedReaction and PAPRECA::BondBreak) for the kMC stage of the PAPRECA run. Note that, for the kMC event to function properly, you have to explicitly declare the bond_style (typically harmonic) along with the relevant bond_coeff in your LAMMPS input file.
On every kMC stage, PAPRECA will attempt to remove bonds with a given probability (based on the chosen rate). Bonds are removed from the system by calling the PAPRECA::deleteBond() LAMMPS wrapper function and the executeBondBreak() function of the papreca.cpp driver code.
You can provide the bond-breaking rate manually or input the activation energy, attempt frequency, and temperature of that kMC event to obtain the corresponding rate from the Arrhenius equation (see rates_calc.h rates_calc.cpp, and PAPRECA::getRateFromArrhenius() ).
When the catalyzed keyword is used, the bond-breaking event can be selected and executed (within the kMC stage of the PAPRECA run) only if at least one atom of the specified atom type is present in the full neighbor list of the parent atom (i.e., the atom on which the bond-breaking was discovered). Note that, PAPRECA will not check if the provided catalyzing atom type is valid (i.e., exists in your simulation and has been defined in the LAMMPS input file). Also, note that the bond-breaking probabilities will be influenced by any settings related to the building and updating of the neighbor lists. For example, if your pair_style cutoff is too small, then fewer "catalyzing" types will be in the neighborhood of the parent atom.
When the limit keyword is used, bond-breaking events are only valid if the current distance between bonded atoms obeys the following inequality: (1-length_perc) * length_equil <= (1+length_perc) * length_equil. Careful, to avoid bond-missing errors it is suggested that the length_equil variable is set to the equilibrium bond length as defined in the LAMMPS input file.
Note 1: In the current version each pair of atom types (e.g., type 1 and type 2) and their corresponding bond (e.g., bond type 5 for atom types 1 and 2) are allowed to be associated with only one predefined bond-breaking template.
Note 2: PAPRECA uses a dummy group with id=1 to perform bond deletions (see this LAMMPS wrapper function for more information PAPRECA::deleteBond()), which means that you should NOT use this bond id to define some another bonded interaction. Please see the "kmc.lmp" file located in ./Examples/Phosphate Film Growth from TCP on Fe110/ for an example demonstrating how you can define multiple bond types.
Create a predefined bond-formation template (see PAPRECA::PredefinedBondForm and PAPRECA::BondForm) for the kMC stage of the PAPRECA run. Note that, for the kMC event to function properly, you have to explicitly declare the bond_style (typically harmonic) along with the relevant bond_coeff in your LAMMPS input file. PAPRECA does not permit the formation of two bonds between exactly the same atoms. If you wish to create a "double" bond consider defining an additional bond type whose coefficients are representative of a "double" bond.
On every kMC stage, PAPRECA will attempt to create the bond (as defined in the template) with a given probability (based on the chosen rate). Bonds are created by calling PAPRECA::formBond() LAMMPS wrapper function and the executeBondForm() function of the papreca.cpp driver code.
You can provide the bond-formation rate manually or input the activation energy, attempt frequency, and temperature of that kMC event to obtain the corresponding rate from the Arrhenius equation (see rates_calc.h rates_calc.cpp, and PAPRECA::getRateFromArrhenius() ).
For more information regarding the catalyzed option please see create_BondBreak command.
Note: In the current version each pair of atom types (e.g., type 1 and type 2) and their corresponding bond (e.g., bond type 5 for atom types 1 and 2) are allowed to be associated with only one predefined bond-forming template.
Sets a maximum number of bonds for a specific atom type. This command can be useful to prevent "overbonding" in PAPRECA studies involving bond-formation events.
The maximum number of bonds of a species is unlimited if the user does not include this command in the PAPRECA input file.
Sets a maximum number of bonds of a specific bond type for a specific atom type. This command can be useful to prevent "overbonding" in PAPRECA studies involving bond-formation events.
The maximum number of bonds of a specific bond type for a specific atom species is unlimited if the user does not include this command in the PAPRECA input file.
Create a predefined deposition template (see PAPRECA::PredefinedDeposition and PAPRECA::Deposition) for the kMC stage of the PAPRECA run. This command will only work properly if you have previously defined a molecule template (i.e., you have prepared a separate molecule file and you have used the molecule command in your LAMMPS input file. See this part of the LAMMPS documentation for more information.
On every kMC stage, PAPRECA will attempt to insert a molecule in the system with a given probability (based on the chosen rate). Molecules are inserted in the system by calling the PAPRECA::insertMolecule() LAMMPS wrapper function and the executeDeposition() function of the papreca.cpp driver code. Note that, the deposition candidate coordinates coincide with the geometric center of the molecule (as defined in the LAMMPS input file).
For sticking_coeff = constant the user must select a constant sticking coefficient value. The sticking coefficient value will remain unchanged (as set) throughout the simulation. Conversely, if sticking_coeff = variable, then the sticking coefficient is dynamically calculated by dividing the number of available (collision-free) deposition sites by the total number of sites (occupied and collision-free) [1] , [2]. Note that, if you have defined multiple deposition events with the same adsorbate_name (e.g., mmmTCP), then PAPRECA will assume that the adsorption sites of all such events are identical and will assign/calculate an identical sticking coefficient.
You can provide the deposition rate manually or input the activation energy, attempt frequency, and temperature of that kMC event to obtain the corresponding rate from the Arrhenius equation (see rates_calc.h rates_calc.cpp, and PAPRECA::getRateFromArrhenius() ). For predefined deposition events the rate can also be calculated from the kinetic theory of gases (Hertz-Knudsen equation). Please see this function PAPRECA::getDepoRateFromHertzKnudsen() and here for a brief theoretical background.
Note: In the current version each parent atom (e.g., type 1) is allowed to be associated with only one predefined deposition template. However, different atom types can be associated with the same adsorbate_name (in separate deposition templates). See the phosphates example located in (./Examples/Phosphate Film Growth from TCP on Fe110).
If the sticking_coeff keyword is not used, then the predefined deposition event will be initialized with variable sticking coefficients.
[1] Ntioudis, S., et al. "PAPRECA: A parallel hybrid off-lattice kinetic Monte Carlo/molecular dynamics simulator", Journal of Open Source Software, 9(98), 6714 (2024). https://doi.org/10.21105/joss.06714
[2] Ntioudis, S., et al. "A hybrid off-lattice kinetic Monte Carlo/molecular dynamics method for amorphous thin film growth", Computational Materials Science, 229, 112421 (2023). https://doi.org/10.1016/j.commatsci.2023.112421
When random_depovecs are not activated (i.e., random_depovecs = no) the center of mass (COM) of any deposited molecule is located directly above the parent atom (i.e., the atom on which searches for kMC events are performed) and at a distance of depo_offset (see create_Deposition command).
When random_depovecs are activated (i.e., random_depovecs = yes) the COM of any deposited molecule is located on the surface of the upper hemisphere (i.e., along the +z-direction) of a sphere centered on the parent atom. The exact position of the COM is determined randomly (by drawing a random number, see getDepoPointCandidateCoords() for more information).
Note: In the current version deposition candidates are always placed above the parent atom. This means that you cannot use this version to place a molecule underneath the parent atom. Hence, you cannot get a film growing towards -z.
random_depovecs = no
Set limits for discovering deposition events (see create_Deposition command) in the kMC stage of the PAPRECA run. Note that, those limitations are applied to all deposition events defined in the simulation. Specifically, $PAPRECA scans for deposition events between film_height - height_scan and film_height + height scan and rejects deposition candidates above film_height + height_reject.
Note: Prior to using this command you must set up a height calculation for your PAPRECA simulation (i.e., a height_calculation command must be present in your PAPRECA above the current command line).
Deposition scans are not limited to a distance above/below the current film height if the user does not include this command in the PAPRECA input file.
Create a predefined diffusion hop template (see PAPRECA::PredefinedDiffusionHop and PAPRECA::Diffusion) for the kMC stage of the PAPRECA run. A "displacive" diffusion template (i.e., displacive = "yes" ), means that the parent atom moves from its initial position to the position of the vacant site. Conversely, a "non-displacive" diffusion template (i.e., displacive = "no" ), means that the parent atom does not move and that a new atom is inserted on the vacant site. Note that, in the current version of PAPRECA, the inserted atom has zero charge and zero velocity.
On every kMC stage, PAPRECA will attempt to move ("displacive" templates) or insert ("non-displacive" templates) an atom in the system with a given probability (based on the chosen rate). Atoms are diffused in the system by calling the PAPRECA::diffuseAtom() LAMMPS wrapper function and the executeDiffusion() function of the papreca.cpp driver code.
If the custom template "Fe_4PO4neib" is used (see example above for syntax), then diffusion events require at least four PO4 structures in their neighborhood to be valid. A PO4 structure consists of Phosphorus atom which is bonded (with an explicit bonded interaction) to four Oxygen atoms. The "Fe_4PO4neib" custom template performs searches on the neighbor list of the parent atom. Also, note that the diffusion probabilities will be influenced by any settings related to the building and updating of the neighbor lists. For example, if your pair_style cutoff is too small, then fewer PO4 structures will be in the neighborhood of the parent atom. The "Fe_4PO4neib" custom template was created to cover the needs of a very specific application related to the formation and growth of thin film from tricresyl phosphate (TCP) molecules on an iron Fe110 surface [1] , [2].
You can provide the diffusion rate manually or input the activation energy, attempt frequency, and temperature of that kMC event to obtain the corresponding rate from the Arrhenius equation (see rates_calc.h rates_calc.cpp, and PAPRECA::getRateFromArrhenius() ).
Note: In the current version each parent atom (e.g., type 1) is allowed to be associated with only one predefined diffusion template.
[1] Ntioudis, S., et al. "PAPRECA: A parallel hybrid off-lattice kinetic Monte Carlo/molecular dynamics simulator", Journal of Open Source Software, 9(98), 6714 (2024). https://doi.org/10.21105/joss.06714
[2] Ntioudis, S., et al. "A hybrid off-lattice kinetic Monte Carlo/molecular dynamics method for amorphous thin film growth", Computational Materials Science, 229, 112421 (2023). https://doi.org/10.1016/j.commatsci.2023.112421
When the "2D" option is used diffusion sites (i.e., locations towards which the atoms migrate) can only be located above the parent atom (i.e., the atom on which searches for kMC events are performed). This means that atoms will solely diffuse towards the +z-direction. When the "3D" option is used diffusion sites can be located above as well as below the parent atom.
When random_diffvecs are not activated (i.e., random_diffvecs = no) diffusion sites will be directly above/below (depending on whether the "2D" or "3D" options were used) and at a distance of diff_dist (see create_DiffusionHop command).
When random_diffvecs are activated (i.e., random_diffvecs = yes ) diffusion sites will be located on the surface of a sphere centered on the parent atom. The exact position of the COM is selected randomly (by drawing a random number, see getDiffPointCandidateCoords() for more information).
random_diffvecs = no
Create a monoatomic desorption template (see PAPRECA::PredefinedMonoatomicDesorption and PAPRECA::MonoatomicDesorption) for the kMC stage of the PAPRECA run. Note that, in the current PAPRECA version monoatomic desorption can only occur if the candidate atom is lone (i.e., has no explicit bonds with any other system atom).
On every kMC stage, PAPRECA will attempt to delete a single atom from the system with a given probability (based on the chosen rate). Atoms are desorbed in the system by calling the PAPRECA::deletAtoms() LAMMPS wrapper function and the executeMonoatomicDesorption() function of the papreca.cpp driver code.
You can provide the monoatomic desorption rate manually or input the activation energy, attempt frequency, and temperature of that kMC event to obtain the corresponding rate from the Arrhenius equation (see rates_calc.h rates_calc.cpp, and PAPRECA::getRateFromArrhenius() ).
Note: In the current version each parent atom (e.g., type 1) is allowed to be associated with only one monoatomic desorption template.
Generates a file named "heightVtime.log" in the current directory of the PAPRECA run. The HeightVtime file comprises 2 columns. The first column lists the time (in seconds) of a PAPRECA step. The second column lists the height of the system (along the z-coordinate) in units consistent with the LAMMPS length units. PAPRECA will append to the HeightVtime file every N steps.
If this command is not used in your PAPRECA input file, no "heightVtime.log" file will be generated in your run PAPRECA run directory.
Generates a file named "surface_coverage.log" in the current directory of the PAPRECA run. The SurfaceCoverage file comprises 2 columns. The first column lists the time (in seconds) of a PAPRECA step. The second column lists the surface coverage, i.e., the number of occupied sites divided by the total number of sites in the system. Note that, the surface coverage is linked to the sticking coefficient of the deposition event (see create_Deposition command).
PAPRECA will append to the SurfaceCoverage file every N steps.
If this command is not used in your PAPRECA input file, no "surface_coverage.log" file will be generated in your run PAPRECA run directory.
Generates a file named distributions.log every N steps, in the current directory of the PAPRECA run. The first column of the distributions lists the height (in length units as in LAMMPS). The remaining columns list the number of atoms of each type (in each height-bin).
For more information regarding the x-y bins and the domain segmentation please refer to height_calculation command.
If this command is not used in your PAPRECA input file, no "distributions.log" files will be generated in your run PAPRECA run directory.
Generates a file named "execTimes.log" in the current directory of the PAPRECA run. The execTimes file lists the minimum, average, and maximum wall times for the kMC, and MD stages. Note that the minimum, average, and maximum times will be equal if you choose to run the simulation on a single processor.
PAPRECA will append to the execTimes.log file every N steps.
Note: The total walltime (at the bottom of the execTimes.log file) is the sum of all the average times. Note that, if you choose a print frequency (N) different than 1, the reported total wall time will be smaller than the actual walltime. Nevertheless, the total walltime will be printed (by LAMMPS) in the console at the end of the run.
If this command is not used in your PAPRECA input file, no "execTimes.log" file will be generated in your run PAPRECA run directory.
Dump a LAMMPS restart file every N PAPRECA steps, in the current directory. The LAMMPS restart file can be used to restart a PAPRECA run. To restart a PAPRECA run, create a new LAMMPS input where you read a dumped LAMMPS restart file (i.e., read_restart and start your run with that LAMMPS input file and the already used PAPRECA input file.
See the PAPRECA::dumpRestart() function for a brief explanation regarding the reason why the dumping of restart files has to be controlled by PAPRECA (and not directly from the LAMMPS input file).
If this command is not used in your PAPRECA input file, no restart files will be generated in your run PAPRECA run directory.