TimeStampedForceBiasMonteCarlo

TimeStampedForceBiasMonteCarlo(configuration, constraints=None, trajectory_filename=None, steps=None, log_interval=None, method=None, hook_functions=None, write_velocities=None, write_forces=None, write_stresses=None, trajectory_interval=None, trajectory_object_id=None, pre_step_hook=None, post_step_hook=None, measurement_hook=None)

Function for performing a time-stamped force-bias Monte Carlo simulation.

Parameters:
  • configuration (MoleculeConfiguration | BulkConfiguration | DeviceConfiguration | SurfaceConfiguration | MDTrajectory) – The initial configuration or a trajectory. If a trajectory is given the last configuration will be used as the initial configuration for this MC run and the trajectory will be appended to. If a trajectory is given, then trajectory_interval and trajectory_object_id should be None.

  • constraints (list of int | list of constraint objects.) – The list of atomic indices, denoting fixed atoms, or constraint objects. FixAtomConstraints, FixCenterOfMass, RigidBody and FixStrain constraint objects are supported
    Default: [].

  • trajectory_filename (str | None) – The filename of the file to be used for storing the trajectory, or None if no trajectory should be written.
    Default: None.

  • steps (int) – The number of steps to take in the simulation.
    Default: 50.

  • log_interval (int) – The resolution used in saving steps to a trajectory file ad writing the log file, where a value of 1 results in all steps being saved, and e.g. a value of 2 resulting in every other step being discarded.
    Default: 1.

  • method (ForceBiasMontCarlo | ForceBiasMonteCarloNPTBerendsen) – The Monte Carlo method used for the simulation.
    Default: ForceBiasMonteCarlo.

  • hook_functions (HookFunctions) – An optional HookFunctions object, which contains user-defined functions or a list of functions which will be called just before (pre-step) or just after (post-step) the forces evaluation.
    Default: None.

  • write_velocities (bool) – Write the velocities to the trajectory file every log_interval steps. Since, the time-stamped force-bias Monte Carlo algorithm does not use velocities explicitly, zero velocities will be written.
    Default: False.

  • write_forces (bool) – Write the forces to the trajectory file every log_interval steps.
    Default: True.

  • write_stresses (bool) – Write the stress to the trajectory file every log_interval steps.
    Default: False.

  • trajectory_interval (int) – The resolution used in saving steps to a trajectory file. A value of 1 results in all steps being saved; a value of 2 results in every second step being saved; etc. If configuration is a trajectory (i.e. this is a restart MC simulation) this parameter will use the same value as was used in the previous trajectory.
    Default: The same value as log_interval.

  • trajectory_object_id (str | None) – The object id of the trajectory written to trajectory_filename. If a value of None is given, then an object id will be chosen automatically.
    Default: None.

  • pre_step_hook (function | list of functions | None) – An optional user-defined function or a list of functions which will be called just before the forces evaluation. The signature of the function requires the arguments (step, time, configuration, forces, stress). The return status is ignored. Unhandled exceptions will terminate the simulation. If a list is given the functions will be called in the given order.
    Default: None.

  • post_step_hook (function | list of functions | None) – An optional user-defined function or a list of functions which will be called just after the forces evaluation. The signature of the function requires the arguments (step, time, configuration, forces, stress). The return status is ignored. Unhandled exceptions will terminate the simulation. If a list is given the functions will be called in the given order.
    Default: None.

  • measurement_hook (function | list of functions | None) – An optional user-defined function or a list of functions which will be called at the end of each MD step after all constraints have been applied. The signature of the function requires the arguments (step, time, configuration). Optional arguments include, forces, local_forces (for distributed MD), stress, trajectory (the MD trajectory), temperature, pressure, potential_energy, and/or kinetic_energy. The return value should be a dictionary that maps string keys to values. The values may be numbers, numpy arrays, or PhysicalQuantities. These values are stored on the MDTrajectory and may be accessed using the measurement method. Unhandled exceptions will terminate the Molecular Dynamics evaluation. If a list is given the functions will be called in the given order.
    Default: None.

Returns:

The Monte Carlo trajectory.

Return type:

MDTrajectory

Usage Examples

Perform a time-stamped force-bias Monte Carlo run of 500 steps of a bulk silicon crystal, using the Stillinger-Weber potential:

# -------------------------------------------------------------
# Bulk configuration
# -------------------------------------------------------------

# Set up lattice
vector_a = [5.4306, 0.0, 0.0]*Angstrom
vector_b = [0.0, 5.4306, 0.0]*Angstrom
vector_c = [0.0, 0.0, 5.4306]*Angstrom
lattice = UnitCell(vector_a, vector_b, vector_c)

# Define elements
elements = [Silicon, Silicon, Silicon, Silicon, Silicon, Silicon, Silicon,
            Silicon]

# Define coordinates
fractional_coordinates = [[0.0 ,  0.0 ,  0.0 ],
                          [0.25,  0.25,  0.25],
                          [0.5 ,  0.5 ,  0.0 ],
                          [0.75,  0.75,  0.25],
                          [0.5 ,  0.0 ,  0.5 ],
                          [0.75,  0.25,  0.75],
                          [0.0 ,  0.5 ,  0.5 ],
                          [0.25,  0.75,  0.75]]

# Set up configuration
bulk_configuration = BulkConfiguration(
    bravais_lattice=lattice,
    elements=elements,
    fractional_coordinates=fractional_coordinates
    )

# -------------------------------------------------------------
# Calculator
# -------------------------------------------------------------

potentialSet = StillingerWeber_Si_1985()
calculator = TremoloXCalculator(parameters=potentialSet)

bulk_configuration.setCalculator(calculator)
bulk_configuration.update()

# -------------------------------------------------------------
# Time-Stamped Force-Bias Monte Carlo
# -------------------------------------------------------------

method = ForceBiasMonteCarlo(
    reservoir_temperature=500.0*Kelvin,
    max_atom_displacement=0.3*Ang,
)

mc_trajectory = TimeStampedForceBiasMonteCarlo(
    bulk_configuration,
    constraints=[],
    trajectory_filename='tfmc_trajectory.hdf5',
    steps=500,
    log_interval=50,
    method=method
)

bulk_configuration = mc_trajectory.lastImage()


tfmc_example.py

Notes

This function performs a time-stamped force-bias MonteCarlo (TFMC) simulation on the given configuration (see [1] for a summary on the TFMC method). The TFMC algorithm essentially consists of a Monte Carlo sampling of the energy landscape at a given temperature.

For details about the Monte Carlo algorithm and how the simulation parameters are specified see ForceBiasMonteCarlo.

Each Monte Carlo step in the TFMC algorithm is a possible realization of how the atoms might move in the potential energy landscape, starting from a given configuration. Therefore the simulation bears some similarity to Molecular Dynamics, and it is even possible to estimate a time step by which subsequent images are separated on average [2]. The time step in TFMC is typically a factor of 2-50 larger than in Molecular Dynamics simulations, however, one must keep in mind that it is not very rigorously defined and primarily suitable for long-timescale-kinetics.

The TimeStampedForceBiasMonteCarlo function returns an MDTrajectory object. After the simulation has completed, you can extract various properties of the stored snapshots, such as coordinates, or forces, by using the MDTrajectory class methods.

Similar to MolecularDynamics, you can perform custom operations on the configuration during the simulation by using the pre_step_hook or post_step_hook functionality.

The pre-step-hook one is invoked immediately before the force and stress calculations, which allows you to modify atomic positions or the cell vectors.

The post-step-hook is invoked after the force and stress calculation, which means that you can modify the forces and stress, e.g. by adding a bias potential.

Suitable functions must have the signature (step, time, configuration, forces, stress) where step denotes the integration step number, time is the current simulation time, and configuration is the current configuration, forces is the PhysicalQuantity array holding the current forces, and stress is the stress tensor. Note, that you need to modify these two arrays in-place, as the return status of the function is ignored. You can use the class methods of the configuration object (see e.g. BulkConfiguration) to perform operations on the atoms.

Moreover, it is possible to give a list of hook functions, which will be called sequentially at each TFMC step.

Only FixAtomConstraints are supported in TFMC simulations.