pyretis.core

This package defines the core PyRETIS tools.

The core tools are intended to define classes which are used in simulations.

Package structure

Modules

__init__.py
Import core functions from the modules.
box.py (pyretis.core.box)
Definition of the simulation box class.
common.py (pyretis.core.common)
Some common core methods, for instance for initiating classes.
montecarlo.py (pyretis.core.montecarlo)
This module defines methods for performing Monte Carlo moves.
particlefunctions.py (pyretis.core.particlefunctions)
Functions that operate on (a selection of) particles, for instance calculation of the kinetic temperature, pressure, momentum etc.
particles.py (pyretis.core.particles)
Definition of the particle class which is used to represent a collection of particles.
pathensemble.py (pyretis.core.pathensemble)
Definition of a class for a collection of paths (i.e. a path ensemble).
path.py (pyretis.core.path)
This module defines functions and classes for paths.
properties.py (pyretis.core.properties)
This module defines a class for a generic property.
random_gen.py (pyretis.core.random_gen)
This module defines a class for generating random numbers.
retis.py (pyretis.core.retis)
Module defining methods for performing replica exchange transition interface sampling.
system.py (pyretis.core.system)
This module defines the system class which connects different parts (for instance box, forcefield and particles) into a single structure.
tis.py (pyretis.core.tis)
This module contains methods used in the transition interface sampling algorithm.
units.py (pyretis.core.units)
This module defines conversion between units.

Important classes defined in this package

BoxBase (BoxBase)
A base class for a simulation box. The box will also handle the periodic boundaries.
System (System)
A class which defines the system we are working with. This class contain a lot of information and is used to group the information into a structure which the simulations will make use of. Typically the system will contain a reference to a box, a list of particles and also a force field.
Particles (Particles)
A class defining a list of particles. This will contain the positions, velocities and forces for the particles.
Path (Path)
A class representing a path. The path contains snapshots with some additional information (energies and order parameters).
PathEnsemble (PathEnsemble)
A class representing a collection of paths. The path ensemble will not store the full trajectories, but only a simplified representation.
PathEnsembleExt (PathEnsembleExt)
A class representing external path ensembles. This handles additional bookkeeping related to external paths.
RandomGenerator (RandomGenerator)
A class for generating random numbers.

pyretis.core.box module

Definition of a class for a simulation box.

The simulation box handles the periodic boundaries if needed. It is typically referenced via the System class, i.e. as System.box.

Important classes defined here

BoxBase (BoxBase)
Class for a simulation box.
RectangularBox (RectangularBox)
Class representing a rectangular simulation box.
TriclinicBox (TriclinicBox)
Class representing a triclinic simulation box.

Examples

>>> from pyretis.core.box import create_box
>>> box = create_box(cell=[10, 10, 10], periodic=[True, False, True])
pyretis.core.box._cos(angle)[source]

Return cosine of an angle.

We also check if the angle is close to 90.0 and if so, we return a zero.

Parameters:angle (float) – The angle in degrees.
Returns:out (float) – The cosine of the angle.
pyretis.core.box._get_low_high_length(low, high, length, periodic)[source]

Determine box cell parameters from input.

This method will consider the following cases:

  1. We are given low, high and length.
  2. We are given high and length, and determine the low values.
  3. We are given low and length, and determine the high values.
  4. We are given length, assume low to be zero and determine high.
  5. We are given low and high, and determine the length.
  6. We are given just high, assume low to be zero and determine the length.
  7. We are just given low, and assume high and length to be infinite.
  8. We are given none of the values, and assume all to infinite.
Parameters:
  • low (numpy.array or None) – The lower bounds for the box.
  • high (numpy.array or None) – The upper bounds for the box.
  • length (numpy.array or None) – The lengths of the box.
  • periodic (list of boolean or None) – We will assume a periodic box for dimensions where this list is True.
Returns:

  • out[0] (numpy.array) – The updated lower bounds for the box.
  • out[1] (numpy.array) – The updated upper bounds for the box.
  • out[2] (numpy.array) – The updated lengths of the box.
  • out[3] (list of boolean) – The updated periodic settings for the box.

pyretis.core.box.angles_from_box_matrix(box_matrix)[source]

Obtain angles and lengths from a given box matrix.

Parameters:box_matrix (np.array) – The box matrix (2D).
Returns:
  • length (numpy.array) – 1D array, the box-lengths on form [a, b, c].
  • alpha (float) – The alpha angle, in degrees.
  • beta (float) – The beta angle, in degrees.
  • gamma (float) – The gamma angle, in degrees.
pyretis.core.box.box_matrix_to_list(matrix, full=False)[source]

Return a list representation of the box matrix.

This method ensures correct ordering of the elements for PyRETIS: xx, yy, zz, xy, xz, yx, yz, zx, zy.

Parameters:
  • matrix (numpy.array) – A matrix (2D) representing the box.
  • full (boolean, optional) – Return a full set of parameters (9) if set to True. If False, and we need 3 or fewer parameters (i.e. the other 6 are zero) we will only return the 3 non-zero ones.
Returns:

out (list) – A list with the box-parametres.

pyretis.core.box.box_vector_angles(length, alpha, beta, gamma)[source]

Obtain the box matrix from given lengths and angles.

Parameters:
  • length (numpy.array) – 1D array, the box-lengths on form [a, b, c].
  • alpha (float) – The alpha angle, in degrees.
  • beta (float) – The beta angle, in degrees.
  • gamma (float) – The gamma angle, in degrees.
Returns:

out (numpy.array, 2D) – The (upper triangular) box matrix.

pyretis.core.box.check_consistency(low, high, length)[source]

Check that given box bounds are consistent.

Parameters:
  • low (numpy.array) – The lower bounds for the box.
  • high (numpy.array) – The upper bounds for the box.
  • length (numpy.array) – The lengths of the box.
pyretis.core.box.create_box(low=None, high=None, cell=None, periodic=None)[source]

Set up and create a box.

Parameters:
  • low (numpy.array, optional) – 1D array containing the lower bounds of the cell.
  • high (numpy.array, optional) – 1D array containing the higher bounds of the cell.
  • cell (numpy.array, optional) – 1D array, a flattened version of the simulation box matrix. This array is expected to contain 1, 2, 3, 6 or 9 items. These are the xx, yy, zz, xy, xz, yx, yz, zx, zy elements, respectively.
  • periodic (list of boolean, optional) – If periodic[i] then we should apply periodic boundaries to dimension i.
Returns:

out (object like BoxBase) – The object representing the simulation box.

pyretis.core.common module

Definition of some common methods that might be useful.

Important methods defined here

bit_fat_comparer (:py:func`.big_fat_comparer`)
Method to compare two nested list/dictionaries.
compare_objects (:py:func`.compare_objects`)
Method to compare two PyRETIS objects.
counter (:py:func`.counter`)
Function to count the number of iterations.
crossing_counter (:py:func`.crossing_counter`)
Function to count the crossing of a path on an interface.
crossing_finder (:py:func`.crossing_finder`)
Function to get the shooting points of the crossing of a path on an interface.
import_from (import_from())
A method to dynamically import method/classes etc. from user specified modules.
inspect_function (inspect_function())
A method to obtain information about arguments, keyword arguments for functions.
initiate_instance (initiate_instance())
Method to initiate a class with optional arguments.
generic_factory (generic_factory())
Create instances of classes based on settings.
null_move (:py:func`.compare_objects`)
Method to do not move.
priority_checker (:py:func`.priority_checker`)
Method to check ensemble to prioritize.
relative_shoots_select (:py:func`.compare_objects`)
Method to select the shooting ensemble.
segments_counter (:py:func`.segments_counter`)
Function that counts the number of segments between two interfaces.
trim_path_between_interfaces (:py:func`.trim_path_between_interfaces`)
Function to trim a path between interfaces.
select_and_trim_a_segment (:py:func`.select_and_trim_a_segment`)
Function to trim a path between interfaces plus the two external points.
compute_weight (compute_weight())
A method to compute the statistical weight of a path generated by a Stone Skipping and Wire Fencing move.
soft_partial_exit (:py:func`.soft_partial_exit`)
Function that check the presence of the EXIT file, and kindly stops the iterator.
pyretis.core.common._arg_kind(arg)[source]

Determine kind for a given argument.

This method will help inspect_function() to determine the correct kind for arguments.

Parameters:arg (object like inspect.Parameter) – The argument we will determine the type of.
Returns:out (string) – A string we use for determine the kind.
pyretis.core.common._pick_out_arg_kwargs(klass, settings)[source]

Pick out arguments for a class from settings.

Parameters:
  • klass (class) – The class to initiate.
  • settings (dict) – Positional and keyword arguments to pass to klass.__init__().
Returns:

  • out[0] (list) – A list of the positional arguments.
  • out[1] (dict) – The keyword arguments.

pyretis.core.common.big_fat_comparer(any1, any2, hard=False)[source]

Check if two dictionary are the same, regardless their complexity.

Parameters:
  • any1 (anything)
  • any2 (anything)
  • hard (boolean, optional) – Raise ValueError if any1 and any2 are different
Returns:

out (boolean) – True if any1 = any2, false otherwise

pyretis.core.common.compute_weight(path, interfaces, move)[source]

Compute the High Acceptance path weight after a MC move.

This function computes the weights that will be used in the computation of the P cross. This trick allows the use of the High Acceptance version of Stone Skipping or Wire Fencing, allowing the acceptance of B to A paths. The drawback is that swapping moves needs to account also for this different weights. The weight 1 will be returned for a path not generated by SS or WF.

Parameters:
  • path (object like PathBase) – This is the input path which will be checked.
  • interfaces (list/tuple of floats) – These are the interface positions of the form [left, middle, right].
  • move (string, optional) – The MC move to compute the weights for.
Returns:

out[0] (float) – The weight of the path.

pyretis.core.common.counter()[source]

Return how many times this function is called.

pyretis.core.common.crossing_counter(path, interface)[source]

Count the crossing to an interfaces.

Method to count the crosses of a path over an interface.

Parameters:
  • path (object like PathBase) – Input path which will be trimmed.
  • interface (float) – The position of the interface.
Returns:

cnt (integer) – Number of crossing of the given interface.

pyretis.core.common.crossing_finder(path, interface, last_frame=False)[source]

Find the crossing to an interfaces.

Method to select the crosses of a path over an interface.

Parameters:
  • path (object like PathBase) – Input path which will be trimmed.
  • interface (float) – Interface position.
  • last_frame (boolean, optional) – Determines if the last crossing will be selected or not.
Returns:

ph1, ph2 (snapshots) – Snapshots to define the randomly picked crossing, one right before and one right after the interface.

pyretis.core.common.generic_factory(settings, object_map, name='generic')[source]

Create instances of classes based on settings.

This method is intended as a semi-generic factory for creating instances of different objects based on simulation input settings. The input settings define what classes should be created and the object_map defines a mapping between settings and the class.

Parameters:
  • settings (dict) – This defines how we set up and select the order parameter.
  • object_map (dict) – Definitions on how to initiate the different classes.
  • name (string, optional) – Short name for the object type. Only used for error messages.
Returns:

out (instance of a class) – The created object, in case we were successful. Otherwise we return none.

pyretis.core.common.import_from(module_path, function_name)[source]

Import a method/class from a module.

This method will dynamically import a specified method/object from a module and return it. If the module can not be imported or if we can’t find the method/class in the module we will raise exceptions.

Parameters:
  • module_path (string) – The path/filename to load from.
  • function_name (string) – The name of the method/class to load.
Returns:

out (object) – The thing we managed to import.

pyretis.core.common.initiate_instance(klass, settings)[source]

Initialise a class with optional arguments.

Parameters:
  • klass (class) – The class to initiate.
  • settings (dict) – Positional and keyword arguments to pass to klass.__init__().
Returns:

out (instance of klass) – Here, we just return the initiated instance of the given class.

pyretis.core.common.inspect_function(function)[source]

Return arguments/kwargs of a given function.

This method is intended for use where we are checking that we can call a certain function. This method will return arguments and keyword arguments a function expects. This method may be fragile - we assume here that we are not really interested in args and kwargs and we do not look for more information about these here.

Parameters:function (callable) – The function to inspect.
Returns:out (dict) – A dict with the arguments, the following keys are defined:
  • args : list of the positional arguments
  • kwargs : list of keyword arguments
  • varargs : list of arguments
  • keywords : list of keyword arguments
pyretis.core.common.null_move(ensemble, cycle)[source]

Perform a null move for a path ensemble.

The null move simply consist of accepting the last accepted path again.

Parameters:
  • ensemble (dict, it contains:) –

    • path_ensembleobject like PathEnsemble

      This is the path ensemble to update with the null move

  • cycle (integer) – The current cycle number

Returns:

  • out[0] (boolean) – Should the path be accepted or not? Here, it’s always True since the null move is always accepted.
  • out[1] (object like PathBase) – The unchanged path.
  • out[2] (string) – The status will here be ‘ACC’ since we just accept the last accepted path again in this move.

pyretis.core.common.priority_checker(ensembles, settings)[source]

Determine the shooting ensemble during a RETIS simulation.

Here we check whether to do priority shooting or not. If True, we either shoot from the ensemble with the fewest paths or ensemble [0^-] if all ensembles have the same no. of paths.

Parameters:
  • ensembles (list of dictionaries of objects) – Lit of dict of ensembles we are using in a path method.
  • settings (dict) – This dict contains the settings for the RETIS method.
Returns:

out[0] (list) – Returns a list of boolean dictating whether certain ensembles are to be skipped or not.

pyretis.core.common.relative_shoots_select(ensembles, rgen, relative)[source]

Randomly select the ensemble for ‘relative’ shooting moves.

Here we select the ensemble to do the shooting in based on relative probabilities. We draw a random number between 0 and 1 which is used to select the ensemble.

Parameters:
  • ensembles (list of objects like PathEnsemble) – This is a list of the ensembles we are using in the RETIS method.
  • rgen (object like RandomGenerator) – This is a random generator. Here we assume that we can call rgen.rand() to draw random uniform numbers.
  • relative (list of floats) – These are the relative probabilities for the ensembles. We assume here that these numbers are normalised.
Returns:

  • out[0] (integer) – The index of the path ensemble to shoot in.
  • out[1] (object like PathEnsemble) – The selected path ensemble for shooting.

pyretis.core.common.segments_counter(path, interface_l, interface_r, reverse=False)[source]

Count the directional segment between interfaces.

Method to count the number of the directional segments of the path, along the orderp, that connect FROM interface_l TO interface_r.

Parameters:
  • path (object like PathBase) – This is the input path which segments will be counted.
  • interface_r (float) – This is the position of the RIGHT interface.
  • interface_l (float) – This is the position of the LEFT interface.
  • reverse (boolean, optional) – Check on a reversed path.
Returns:

n_segments (integer) – Segment counter

pyretis.core.common.select_and_trim_a_segment(path, interface_l, interface_r, segment_to_pick=None)[source]

Cut a directional segment from interface_l to interface_r.

It keeps what is within the range [interface_l interface_r) AND the snapshots just after/before the interface.

Parameters:
  • path (object like PathBase) – This is the input path which will be trimmed.
  • interface_r (float) – This is the position of the RIGHT interface.
  • interface_l (float) – This is the position of the LEFT interface.
  • segment_to_pick (integer (n.b. it starts from 0)) – This is the segment to be selected, None = random
Returns:

segment (a path segment composed only the snapshots for which) – orderp is between interface_r and interface_l and the ones right after/before the interfaces.

pyretis.core.common.soft_partial_exit(exe_path='')[source]

Check the presence of the EXIT file.

Parameters:exe_path (string, optional) – Path for the EXIT file.
Returns:out (boolean) – True if EXIT is present. False if EXIT in not present.
pyretis.core.common.trim_path_between_interfaces(path, interface_l, interface_r)[source]

Cut a path between the two interfaces.

The method cut a path and keeps only what is within the range (interface_l interface_r). -Be careful, it can provide multiple discontinuous segments- =Be carefull2 consider if you need to make this check left inclusive (as the ensemble should be left inclusive)

Parameters:
  • path (object like PathBase) – This is the input path which will be trimmed.
  • interface_r (float) – This is the position of the RIGHT interface.
  • interface_l (float) – This is the position of the LEFT interface.
Returns:

new_path (object like PathBase) – This is the output trimmed path.

pyretis.core.montecarlo module

Module for Monte Carlo Algorithms and other “random” functions.

In this module, Monte Carlo Algorithms are defined. Note that some derived or “random” functions are also defined in the TIS module.

Important methods defined here

metropolis_accept_reject (metropolis_accept_reject())
Accept/reject an energy change according to the metropolis rule.
max_displace_step (max_displace_step())
Monte Carlo routine for displacing particles. It will select and displace one particle randomly.
pyretis.core.montecarlo.max_displace_step(rgen, system, maxdx=0.1, idx=None)[source]

Monte Carlo routine for displacing particles.

It selects and displaces one particle randomly. If the move is accepted, the new positions and energy are returned. Otherwise, the move is rejected and the old positions and potential energy is returned. The function accept_reject is used to accept/reject the move.

Parameters:
  • rgen (object like RandomGenerator) – The random number generator.
  • system (object like System) – The system object to operate on
  • maxdx (float, optional) – The maximum displacement (default is 0.1).
  • idx (int, optional) – Index of the particle to displace. If idx is not given, the particle is chosen randomly.
Returns:

out (boolean) – The outcome of applying the function accept_reject to the system and trial position.

pyretis.core.montecarlo.metropolis_accept_reject(rgen, system, deltae)[source]

Accept/reject a energy change according to the metropolis rule.

FIXME: Check if metropolis really is a good name here.

Parameters:
  • rgen (object like RandomGenerator) – The random number generator.
  • system (object like System) – The system object we are investigating. This is used to access the beta factor.
  • deltae (float) – The change in energy.
Returns:

out (boolean) – True if the move is accepted, False otherwise.

Notes

An overflow is possible when using numpy.exp() here. This can, for instance, happen in an umbrella simulation where the bias potential is infinite or very large. Right now, this is just ignored.

pyretis.core.particlefunctions module

This file contains functions that act on (a selection of) particles.

These functions are intended for calculating particle properties as the kinetic temperature, pressure etc.

Important methods defined here

atomic_kinetic_energy_tensor (atomic_kinetic_energy_tensor())
Return the kinetic energy tensor for each atom in a selection of particles.
calculate_kinetic_energy (calculate_kinetic_energy())
Return the kinetic energy of a collection of particles.
calculate_kinetic_energy_tensor (calculate_kinetic_energy_tensor())
Return the kinetic energy tensor for a selection of particles.
calculate_kinetic_temperature (calculate_kinetic_temperature())
Return the kinetic temperature of a collection of particles.
calculate_linear_momentum (calculate_linear_momentum())
Calculates the linear momentum of a collection of particles.
calculate_pressure_from_temp (calculate_pressure_from_temp())
Return the scalar pressure using the temperature and the virial.
calculate_pressure_tensor (calculate_pressure_tensor())
Return the pressure tensor, obtained from the virial and the kinetic energy tensor.
calculate_scalar_pressure (calculate_scalar_pressure())
Return the scalar pressure (from the trace of the pressure tensor).
calculate_thermo (calculate_thermo())
Calculate and return several “thermodynamic” properties as the potential, kinetic and total energies per particle, the temperature, the pressure and the momentum.
calculate_thermo_path (calculate_thermo_path())
Calculate and return some thermodynamic properties. This method is similar to the calculate_thermo, however, it is simpler and calculates fewer quantities.
kinetic_energy (kinetic_energy())
Return the kinetic energy for velocities and masses given as numpy arrays.
kinetic_temperature (kinetic_temperature())
Return the temperature for velocities and masses given as numpy arrays.
reset_momentum (reset_momentum())
Set linear momentum (for a selection of particles) to zero.
pyretis.core.particlefunctions._get_vel_mass(particles, selection=None)[source]

Return velocity and mass for a selection.

This is just for convenience since we are using this selection a lot.

Parameters:
  • particles (object like Particles) – This object represents the particles.
  • selection (list of integers, optional) – A list with indices of particles to use in the calculation.
Returns:

  • out[0] (numpy.array) – The velocities corresponding to the selection.
  • out[1] (numpy.array, optional) – The masses corresponding to the selection.

pyretis.core.particlefunctions.atomic_kinetic_energy_tensor(particles, selection=None)[source]

Return kinetic energy tensors for a particle selection.

Parameters:
  • particles (object like Particles) – This object represents the particles.
  • selection (list of integers, optional) – A list with indices of particles to use in the calculation.
Returns:

kin (numpy.array) – A numpy array with dimensionality equal to (len(selection), dim, dim) where dim is the number of dimensions used in the velocities. kin[i] contains the kinetic energy tensor formed by the outer product of mol[selection][i] and vel[selection][i]. The sum of the tensor should equal the output from calculate_kinetic_energy_tensor.

pyretis.core.particlefunctions.calculate_kinetic_energy(particles, selection=None, kin_tensor=None)[source]

Return the kinetic energy of a collection of particles.

Parameters:
  • particles (object like Particles) – This object represents the particles.
  • selection (list of integers, optional) – A list with indices of particles to use in the calculation.
  • kin_tensor (numpy.array, optional) – If kinetic_tensor is not given, the kinetic energy tensor will be calculated.
Returns:

  • out[0] (float) – The scalar kinetic energy.
  • out[1] (numpy.array) – The kinetic energy tensor.

pyretis.core.particlefunctions.calculate_kinetic_energy_tensor(particles, selection=None)[source]

Return the kinetic energy tensor for a selection of particles.

The tensor is formed as the outer product of the velocities.

Parameters:
  • particles (object like Particles) – This object represents the particles.
  • selection (list of integers, optional) – A list with indices of particles to use in the calculation.
Returns:

out (numpy.array) – The kinetic energy tensor. Dimensionality is equal to (dim, dim) where dim is the number of dimensions used in the velocities. The trace gives the kinetic energy.

pyretis.core.particlefunctions.calculate_kinetic_temperature(particles, boltzmann, dof=None, selection=None, kin_tensor=None)[source]

Return the kinetic temperature of a collection of particles.

Parameters:
  • particles (object like Particles) – This object represents the particles.
  • boltzmann (float) – This is the Boltzmann factor/constant in correct units.
  • dof (list of floats, optional) – The degrees of freedom to subtract. Its shape should be equal to the number of dimensions.
  • selection (list of integers, optional) – A list with indices of particles to use in the calculation.
  • kin_tensor (numpy.array optional) – The kinetic energy tensor. If the kinetic energy tensor is not given, it will be recalculated here.
Returns:

  • out[0] (numpy.array) – Array with the same size as the kinetic energy. It contains the temperature in each spatial dimension.
  • out[1] (float) – The temperature averaged over all dimensions.
  • out[2] (numpy.array) – The kinetic energy tensor.

pyretis.core.particlefunctions.calculate_linear_momentum(particles, selection=None)[source]

Calculate the linear momentum for a collection of particles.

Parameters:
  • particles (object like Particles) – This object represents the particles.
  • selection (list of integers, optional) – A list with indices of particles to use in the calculation.
Returns:

out (numpy.array) – The array contains the linear momentum for each dimension.

pyretis.core.particlefunctions.calculate_pressure_from_temp(particles, dim, boltzmann, volume, dof=None)[source]

Evaluate the scalar pressure.

The scalar pressure is here calculated using the temperature and the degrees of freedom.

Parameters:
  • particles (object like Particles) – This object represents the particles.
  • dim (int) – This is the dimensionality of the system. Typically provided by system.get_dim().
  • boltzmann (float) – This is the Boltzmann factor/constant in correct units. Typically it can be supplied by system.get_boltzmann().
  • volume (float) – This is the volume ‘occupied’ by the particles. It can typically be obtained by a box.calculate_volume()
  • dof (list of floats, optional) – The degrees of freedom to subtract. Its shape should be equal to the number of dimensions.
Returns:

  • out[0] (float) – Pressure times volume.
  • out[1] (float) – The pressure.

Note

This function may possibly be removed - it does not appear to be very useful right now.

pyretis.core.particlefunctions.calculate_pressure_tensor(particles, volume, kin_tensor=None)[source]

Calculate the pressure tensor.

The pressure tensor is obtained from the virial the kinetic energy tensor.

Parameters:
  • particles (object like Particles) – This object represents the particles.
  • volume (float) – This is the volume ‘occupied’ by the particles. It can typically be obtained by a box.calculate_volume().
  • kin_tensor (numpy.array, optional) – The kinetic energy tensor. If kin_tensor is not given, it will be calculated here.
Returns:

out (numpy.array) – The symmetric pressure tensor, dimensions (dim, dim), where dim = the number of dimensions considered in the simulation.

pyretis.core.particlefunctions.calculate_scalar_pressure(particles, volume, dim, press_tensor=None, kin_tensor=None)[source]

Evaluate the scalar pressure using the pressure tensor.

Parameters:
  • particles (object like Particles) – This object represents the particles.
  • volume (float) – This is the volume ‘occupied’ by the particles. It can typically be obtained by a box.calculate_volume().
  • dim (int) – This is the dimensionality of the system. Typically provided by system.get_dim()
  • press_tensor (numpy.array, optional) – If press_tensor is not given, the pressure tensor will be calculated here.
  • kin_tensor (numpy.array, optional) – If kin_tensor is not given, the kinetic energy tensor will be calculated here.
Returns:

out (float) – The scalar pressure averaged over the diagonal components of the pressure tensor.

pyretis.core.particlefunctions.calculate_thermo(system, dof=None, dim=None, volume=None, vpot=None)[source]

Calculate and return several thermodynamic properties.

The calculated properties are the potential, kinetic and total energies per particle, the temperature, the pressure and the momentum.

Parameters:
  • system (object like System) – This object is used to access the particles and the box.
  • dof (list of floats, optional) – The degrees of freedom.
  • dim (float, optional) – The dimensionality of, typically provided with a system.get_dim().
  • volume (float, optional) – This is the volume ‘occupied’ by the particles. It can typically be obtained by a system.box.calculate_volume().
  • vpot (float, optional) – The potential energy of the particles.
Returns:

out (dict) – This dict contains the float that is calculated in this routine.

pyretis.core.particlefunctions.calculate_thermo_path(system)[source]

Calculate and return several thermodynamic properties.

The calculated properties are the potential, kinetic and total energies for the system and the current temperature. The name here calculate_thermo_path just indicates that this function is useful in connection with path sampling simulations, i.e. it just calculates a few energy terms.

Parameters:system (object like System) – This object is used to access the particles and the box.
Returns:out (dict) – This dict contains the float that is calculated in this routine.
pyretis.core.particlefunctions.kinetic_energy(vel, mass)[source]

Obtain the kinetic energy for given velocities and masses.

Parameters:
  • vel (numpy.array) – The velocities
  • mass (numpy.array) – The masses. This is assumed to be a column vector.
Returns:

  • out[0] (float) – The kinetic energy
  • out[1] (numpy.array) – The kinetic energy tensor.

pyretis.core.particlefunctions.kinetic_temperature(vel, mass, boltzmann, dof=None, kin_tensor=None)[source]

Return the kinetic temperature given velocities and masses.

This method does not work on a particle object, but rather with numpy arrays. That is, it is intended for use when we can’t rely on the particle object.

Parameters:
  • vel (numpy.array) – The velocities
  • mass (numpy.array) – The masses. This is assumed to be a column vector.
  • boltzmann (float) – This is the Boltzmann factor/constant in correct units.
  • dof (list of floats, optional) – The degrees of freedom to subtract. Its shape should be equal to the number of dimensions.
  • kin_tensor (numpy.array optional) – The kinetic energy tensor. If the kinetic energy tensor is not given, it will be recalculated here.
Returns:

  • out[0] (numpy.array) – Array with the same size as the kinetic energy. It contains the temperature in each spatial dimension.
  • out[1] (float) – The temperature averaged over all dimensions.
  • out[2] (numpy.array) – The kinetic energy tensor.

pyretis.core.particlefunctions.reset_momentum(particles, selection=None, dim=None)[source]

Set the linear momentum of a selection of particles to zero.

Parameters:
  • particles (object like Particles) – This object represents the particles.
  • selection (list of integers, optional) – A list with indices of particles to use in the calculation.
  • dim (list or None, optional) – If dim is None, the momentum will be reset for ALL dimensions. Otherwise, it will only be applied to the dimensions where dim is True.
Returns:

out (None) – Returns None and modifies velocities of the selected particles.

pyretis.core.particles module

This module defines a class, representing a collection of particles.

The class for particles is, in reality, a simplistic particle list which stores positions, velocities, masses etc. and is used for representing the particles in the simulations.

Important classes defined here

Particles (Particles)
Class for a list of particles.
ParticlesExt (ParticlesExt)
Class for an external particle list.
class pyretis.core.particles.Particles(dim=1)[source]

Bases: object

Base class for a collection of particles.

This is a simple particle list. It stores the positions, velocities, forces, masses (and inverse masses) and type information for a set of particles. This class also defines a method for iterating over pairs of particles, which could be useful for implementing neighbour lists. In this particular class, this method will just define an all-pairs list.

npart

Number of particles.

Type:integer
pos

Positions of the particles.

Type:numpy.array
vel

Velocities of the particles.

Type:numpy.array
force

Forces on the particles.

Type:numpy.array
virial

The current virial for the particles.

Type:numpy.array
mass

Masses of the particles.

Type:numpy.array
imass

Inverse masses, 1.0 / self.mass.

Type:numpy.array
name

A name for the particle. This may be used as short text describing the particle.

Type:list of strings
ptype

A type for the particle. Particles with identical ptype are of the same kind.

Type:numpy.array of integers
dim

This variable is the dimensionality of the particle list. This should be derived from the box. For some functions, it is convenient to be able to access the dimensionality directly from the particle list. It is therefore set as an attribute here.

Type:int
vpot

The potential energy of the particles.

Type:float
ekin

The kinetic energy of the particles.

Type:float
__eq__(other)[source]

Compare two particle states.

__init__(dim=1)[source]

Initialise the Particle list.

Here we just create an empty particle list.

Parameters:dim (integer, optional) – The number of dimensions we are considering for positions, velocities and forces.
__iter__()[source]

Iterate over the particles.

This function will yield the properties of the different particles.

Yields:out (dict) – The information in self.pos, self.vel, … etc.
__str__()[source]

Print out basic info about the particle list.

_copy_attribute(attr, copy_function)[source]

Copy an attribute.

Parameters:
  • attr (string) – The attribute to copy.
  • copy_function (callable) – The method to use for copying the attribute.
Returns:

out (object) – A copy of the selected attribute.

add_particle(pos, vel, force, mass=1.0, name='?', ptype=0)[source]

Add a particle to the system.

Parameters:
  • pos (numpy.array) – Positions of new particle.
  • vel (numpy.array) – Velocities of new particle.
  • force (numpy.array) – Forces on the new particle.
  • mass (float, optional) – The mass of the particle.
  • name (string, optional) – The name of the particle.
  • ptype (integer, optional) – The particle type.
copy()[source]

Return a copy of the particle state.

Returns:out (object like Particles) – A copy of the current state of the particles.
empty_list()[source]

Reset the particle list.

This will delete all particles in the list and set other variables to None.

Note

This is almost self.__init__ repeated. The reason for this is simply that we want to define all attributes in self.__init__ and not get any ‘surprise attributes’ defined elsewhere. Also, note that the dimensionality (self.dim) is not changed in this method.

get_force()[source]

Return (a copy of) the forces.

get_pos()[source]

Return (a copy of) positions.

get_selection(properties, selection=None)[source]

Return selected properties for a selection of particles.

Parameters:
  • properties (list of strings) – The strings represent the properties to return.
  • selection (list with indices to return, optional) – If a selection is not given, data for all particles are returned.
Returns:

out (list) – A list with the properties in the order they were asked for in the properties argument.

get_vel()[source]

Return (a copy of) the velocities.

load_restart_info(info)[source]

Load restart information.

Parameters:info (dict) – Dictionary with the settings to load.
pairs()[source]

Iterate over all pairs of particles.

Yields:
  • out[0] (integer) – The index for the first particle in the pair.
  • out[1] (integer) – The index for the second particle in the pair.
  • out[2] (integer) – The particle type of the first particle.
  • out[3] (integer) – The particle type of the second particle.
particle_type = 'internal'
restart_info()[source]

Generate information for saving a restart file.

reverse_velocities()[source]

Reverse the velocities in the system.

set_force(force)[source]

Set the forces for the particles.

Parameters:force (numpy.array) – The forces to set.
set_pos(pos)[source]

Set positions for the particles.

Parameters:pos (numpy.array) – The positions to set.
set_vel(vel)[source]

Set velocities for the particles.

Parameters:vel (numpy.array) – The velocities to set.
class pyretis.core.particles.ParticlesExt(dim=1)[source]

Bases: Particles

A particle list, when positions and velocities are stored in files.

config

The file name and index in this file for the configuration the particle list is representing.

Type:tuple
vel_rev

If this is True, the velocities in the file represeting the configuration will have to be reversed before they are used.

Type:boolean
__init__(dim=1)[source]

Create an empty ParticleExt list.

Parameters:dim (integer, optional) – The number of dimensions we are considering for positions, velocities and forces.
__str__()[source]

Print out basic info about the particle list.

add_particle(pos, vel, force, mass=1.0, name='?', ptype=0)[source]

Add a particle to the system.

Parameters:
  • pos (tuple) – Positions of new particle.
  • vel (boolean) – Velocities of new particle.
  • force (tuple) – Forces on the new particle.
  • mass (float, optional) – The mass of the particle.
  • name (string, optional) – The name of the particle.
  • ptype (integer, optional) – The particle type.
empty_list()[source]

Just empty the list.

get_pos()[source]

Just return the positions of the particles.

get_vel()[source]

Return info about the velocities.

particle_type = 'external'
reverse_velocities()[source]

Reverse the velocities in the system.

set_pos(pos)[source]

Set positions for the particles.

This will copy the input positions, for this class, the input positions are assumed to be a file name with a corresponding integer which determines the index for the positions in the file for cases where the file contains several snapshots.

Parameters:pos (tuple of (string, int)) – The positions to set, this represents the file name and the index for the frame in this file.
set_vel(rev_vel)[source]

Set velocities for the particles.

Here we store information which tells if the velocities should be reversed or not.

Parameters:rev_vel (boolean) – The velocities to set. If True, the velocities should be reversed before used.

pyretis.core.path module

Classes and functions for paths.

The classes and functions defined in this module are useful for representing paths.

Important classes defined here

PathBase (PathBase)
A base class for paths.
Path (Path)
Class for a generic path that stores all possible information.

Important methods defined here

paste_paths (paste_paths())
Function for joining two paths, one is in a backward time direction and the other is in the forward time direction.
class pyretis.core.path.Path(rgen=None, maxlen=None, time_origin=0)[source]

Bases: PathBase

A path where the full trajectory is stored in memory.

This class represents a path. A path consists of a series of consecutive snapshots (the trajectory) with the corresponding order parameter. Here we store all information for all phase points on the path.

empty_path(**kwargs)[source]

Return an empty path of same class as the current one.

Returns:out (object like PathBase) – A new empty path.
get_shooting_point(criteria='rnd', interfaces=None)[source]

Return a shooting point from the path.

This will simply draw a shooting point from the path at random. All points can be selected with equal probability with the exception of the end points which are not considered.

Parameters:
  • criteria (string, optional) – The criteria to select the shooting point: ‘rnd’: random, except the first and last point, standard sh. ‘exp’: selection towards low density region.
  • list/tuple of floats, optional – These are the interface positions of the form [left, middle, right].
Returns:

  • out[0] (object like System) – The phase point we selected.
  • out[1] (int) – The shooting point index.

load_restart_info(info)[source]

Set up the path using restart information.

restart_info()[source]

Return a dictionary with restart information.

class pyretis.core.path.PathBase(rgen=None, maxlen=None, time_origin=0)[source]

Bases: object

Base class for representation of paths.

This class represents a path. A path consists of a series of consecutive snapshots (the trajectory) with the corresponding order parameter.

generated

This contains information on how the path was generated. generated[0] : string, as defined in the variable _GENERATED generated[1:] : additional information: For generated[0] == 'sh' the additional information is the index of the shooting point on the old path, the new path and the corresponding order parameter.

Type:tuple
maxlen

This is the maximum path length. Some algorithms require this to be set. Others don’t, which is indicated by setting maxlen equal to None.

Type:int
ordermax

This is the (current) maximum order parameter for the path. ordermax[0] is the value, ordermax[1] is the index in self.path.

Type:tuple
ordermin

This is the (current) minimum order parameter for the path. ordermin[0] is the value, ordermin[1] is the index in self.path.

Type:tuple
phasepoints

The phase points the path is made up of.

Type:list of objects like System
rgen

This is the random generator that will be used for the paths that required random numbers.

Type:object like RandomGenerator
time_origin

This is the location of the phase point path[0] relative to its parent. This might be useful for plotting.

Type:int
status

The status of the path. The possibilities are defined in the variable _STATUS.

Type:str or None
weight

The statistical weight of the path.

Type:real
__eq__(other)[source]

Check if two paths are equal.

__iadd__(other)[source]

Add path data to a path from another path, i.e. self += other.

This will simply append the phase points from other.

Parameters:other (object of type Path) – The object to add path data from.
Returns:self (object of type Path) – The updated path object.
__init__(rgen=None, maxlen=None, time_origin=0)[source]

Initialise the PathBase object.

Parameters:
  • rgen (object like RandomGenerator, optional) – This is the random generator that will be used.
  • maxlen (int, optional) – This is the max-length of the path. The default value, None, is just a path of arbitrary length.
  • time_origin (int, optional) – This can be used to store the shooting point of a parent trajectory.
__ne__(other)[source]

Check if two paths are not equal.

__str__()[source]

Return a simple string representation of the Path.

append(phasepoint)[source]

Append a new phase point to the path.

Parameters:out (object like System) – The system information we add to the path.
check_interfaces(interfaces)[source]

Check current status of the path.

Get the current status of the path with respect to the interfaces. This is intended to determine if we have crossed certain interfaces or not.

Parameters:interfaces (list of floats) – This list is assumed to contain the three interface values left, middle and right.
Returns:
  • out[0] (str, ‘L’ or ‘R’ or None) – Start condition: did the trajectory start at the left (‘L’) or right (‘R’) interface.
  • out[1] (str, ‘L’ or ‘R’ or None) – Ending condition: did the trajectory end at the left (‘L’) or right (‘R’) interface or None of them.
  • out[2] str, ‘M’ or ‘’* – ‘M’ if the middle interface is crossed, ‘*’ otherwise.
  • out[3] (list of boolean) – These values are given by ordermin < interfaces[i] <= ordermax.
copy()[source]

Return a copy of the path.

delete(idx)[source]

Remove a phase point from the path.

Parameters:idx (integer) – The index of the frame to remove.
abstract empty_path(**kwargs)[source]

Return an empty path of same class as the current one.

This function is intended to spawn sibling paths that share some properties and also some characteristics of the current path. The idea here is that a path of a certain class should only be able to create paths of the same class.

Returns:out (object like PathBase) – A new empty path.
get_end_point(left, right=None)[source]

Return the end point of the path as a string.

The end point is either to the left of the left interface or to the right of the right interface, or somewhere in between.

Parameters:
  • left (float) – The left interface.
  • right (float, optional) – The right interface, equal to left if not specified.
Returns:

out (string) – A string representing where the end point is (‘L’ - left, ‘R’ - right or None).

get_move()[source]

Return the move used to generate the path.

get_path_data(status, interfaces)[source]

Return information about the path.

This information can be stored in a object like PathEnsemble.

Parameters:
  • status (string) – This represents the current status of the path.
  • interfaces (list) – These are just the interfaces we are currently considering.
abstract get_shooting_point()[source]

Return a shooting point from the path.

Returns:
  • phasepoint (object like System) – A phase point which will be the state to shoot from.
  • idx (int) – The index of the shooting point.
get_start_point(left, right=None)[source]

Return the start point of the path as a string.

The start point is either to the left of the left interface or to the right of the right interface.

Parameters:
  • left (float) – The left interface.
  • right (float, optional) – The right interface, equal to left if not specified.
Returns:

out (string) – A string representing where the start point is (‘L’ - left, ‘R’ - right or None).

property length

Compute the length of the path.

abstract load_restart_info(info)[source]

Read a dictionary with restart information.

property ordermax

Compute the maximum order parameter of the path.

property ordermin

Compute the minimum order parameter of the path.

abstract restart_info()[source]

Return a dictionary with restart information.

reverse(order_function=False, rev_v=True)[source]

Reverse a path and return the reverse path as a new path.

This will reverse a path and return the reversed path as a new object like PathBase object.

Returns:
  • new_path (object like PathBase) – The time reversed path.
  • order_function (object like OrderParameter, optional) – The method to use to re-calculate the order parameter, if it is velocity dependent.
  • rev_v (boolean, optional) – If True, also the velocities are reversed, if False, the velocities for each frame are not altered.
static reverse_velocities(system)[source]

Reverse the velocities in the phase points.

set_move(move)[source]

Update the path move.

The path move is a short string that represents how the path was generated. It should preferably match one of the moves defined in _GENERATED.

Parameters:move (string) – A short description of the move.
sorting(key, reverse=False)[source]

Re-order the phase points according to the given key.

Parameters:
  • key (string) – The attribute we will sort according to.
  • reverse (boolean, optional) – If this is False, the sorting is from big to small.
Yields:

out (object like System) – The ordered phase points from the path.

success(target_interface)[source]

Check if the path is successful.

The check is based on the maximum order parameter and the value of target_interface. It is successful if the maximum order parameter is greater than target_interface.

Parameters:target_interface (float) – The value for which the path is successful, i.e. the “target_interface” interface.
update_energies(ekin, vpot)[source]

Update the energies for the phase points.

This method is useful in cases where the energies are read from external engines and returned as a list of floats.

Parameters:
  • ekin (list of floats) – The kinetic energies to set.
  • vpot (list of floats) – The potential energies to set.
pyretis.core.path.check_crossing(cycle, orderp, interfaces, leftside_prev)[source]

Check if we have crossed an interface during the last step.

This function is useful for checking if an interface was crossed from the previous step till the current one. This is for instance used in the MD simulations for the initial flux. If will use a variable to store the previous positions with respect to the interfaces and check if interfaces were crossed here.

Parameters:
  • cycle (int) – This is the current simulation cycle number.
  • orderp (float) – The current order parameter.
  • interfaces (list of floats) – These are the interfaces to check.
  • leftside_prev (list of booleans or None) – These are used to store the previous positions with respect to the interfaces.
Returns:

  • leftside_curr (list of booleans) – These are the updated positions with respect to the interfaces.
  • cross (list of tuples) – If a certain interface is crossed, a tuple will be added to this list. The tuple is of form (cycle number, interface number, direction) where the direction is ‘-’ for a crossing in the negative direction and ‘+’ for a crossing in the positive direction.

pyretis.core.path.paste_paths(path_back, path_forw, overlap=True, maxlen=None)[source]

Merge a backward with a forward path into a new path.

The resulting path is equal to the two paths stacked, in correct time. Note that the ordering is important here so that: paste_paths(path1, path2) != paste_paths(path2, path1).

There are two things we need to take care of here:

  • path_back must be iterated in reverse (it is assumed to be a backward trajectory).
  • we may have to remove one point in path2 (if the paths overlap).
Parameters:
  • path_back (object like PathBase) – This is the backward trajectory.
  • path_forw (object like PathBase) – This is the forward trajectory.
  • overlap (boolean, optional) – If True, path_back and path_forw have a common starting-point, that is, the first point in path_forw is identical to the first point in path_back. In time-space, this means that the first point in path_forw is identical to the last point in path_back (the backward and forward path started at the same location in space).
  • maxlen (float, optional) – This is the maximum length for the new path. If it’s not given, it will just be set to the largest of the maxlen of the two given paths.

Note

Some information about the path will not be set here. This must be set elsewhere. This includes how the path was generated (path.generated) and the status of the path (path.status).

pyretis.core.pathensemble module

Classes and functions for path ensembles.

The classes and functions defined in this module are useful for representing path ensembles.

Important classes defined here

PathEnsemble (PathEnsemble)
Class for defining path ensembles.
PathEnsembleExt (PathEnsembleExt)
Class for defining path ensembles when we are working with paths stored on disk and not in memory only.
class pyretis.core.pathensemble.PathEnsemble(ensemble_number, interfaces, rgen=None, maxpath=10000, exe_dir=None)[source]

Bases: object

Representation of a path ensemble.

This class represents a collection of paths in a path ensemble. In general, paths may be “long and complicated” so here, we really just store a simplified abstraction of the path, which is obtained by the Path.get_path_data() function for a given Path object. The returned dictionary is stored in the list PathEnsemble.paths. The only full path we store is the last accepted path. This is convenient for the RETIS method where paths may be swapped between path ensembles.

ensemble_number

This integer is used to represent the path ensemble, for RETIS simulations it’s useful to identify the path ensemble. The path ensembles are numbered sequentially 0, 1, 2, etc. This corresponds to [0^-], [0^+], [1^+], etc.

Type:integer
ensemble_name

A string which can be used for printing the ensemble name. This is of form [0^-], [0^+], [1^+], etc.

Type:string
ensemble_name_simple

A string with a simpler representation of the ensemble name, can be used for creating output files etc.

Type:string
interfaces

Interfaces, specified with the values for the order parameters: [left, middle, right].

Type:list of floats
paths

This list contains the stored information for the paths. Here we only store the data returned by calling the get_path_data() function of the Path object.

Type:list
nstats

This dict store some statistics for the path ensemble. The keys are:

  • npath : The number of paths stored.
  • nshoot : The number of accepted paths generated by shooting.
  • ACC, BWI, … : Number of paths with given status (from _STATUS).
Type:dict of ints
maxpath

The maximum number of paths to store.

Type:int
last_path

This is the last accepted path.

Type:object like PathBase
__eq__(other)[source]

Check if two path_ensemble are equal.

__init__(ensemble_number, interfaces, rgen=None, maxpath=10000, exe_dir=None)[source]

Initialise the PathEnsemble object.

Parameters:
  • ensemble_number (integer) – An integer used to identify the ensemble.
  • interfaces (list of floats) – These are the interfaces specified with the values for the order parameters: [left, middle, right].
  • rgen (object like RandomGenerator, optional) – The random generator that will be used for the paths that required random numbers.
  • maxpath (integer, optional) – The maximum number of paths to store information for in memory. Note, that this will not influence the analysis as long as you are using the output files when running the analysis.
  • exe_dir (string, optional) – The base folder where the simulation was executed from. This is used to set up output directories for the path ensemble.
__ne__(other)[source]

Check if two paths are not equal.

__str__()[source]

Return a string with some info about the path ensemble.

add_path_data(path, status, cycle=0)[source]

Append data from the given path to self.path_data.

This will add the data from a given` path` to the list path data for this ensemble. If will also update self.last_path if the given path is accepted.

Parameters:
  • path (object like PathBase) – This is the object to store data from.
  • status (string) – This is the status of the path. Note that the path object also has a status property. However, this one might not be set, for instance when the path is just None. We therefore use status here as a parameter.
  • cycle (int, optional) – The current cycle number.
clear_generate()[source]

Remove all the files in an ensemble/generate/ directory.

This is toggled on by adding ‘remove_generate = True’ to the retis.rst input file, under the [simulation] section. TODO: Combining this with the retis_swap_zero move causes problems

during a ‘save-cycle’. As the files of 000/generate are removed, it will try to copy 000/generate/second_last.gro to 001/traj/…, but is already deleted. I did a quick-fix by not deleting GRO or G96 files. This only works for GROMACS, and a more general fix will be needed for other external engines. Update: I cannot recreate the problem anymore. I’m thinking now that it had to do with the fact that rejected swapping moves did not return the swapped paths of the previous accepted cycle. Therefore, the pointers were not updated correctly. As this has been fixed during unittesting, I think this is okay now.
copy_path_to_generate(path, _prefix=None)[source]

Copy a path for temporary storing.

Parameters:path (object like PathBase) – The path to copy.
Returns:out (object like PathBase) – The copy of the path.
directories()[source]

Yield the directories PyRETIS should make.

get_acceptance_rate()[source]

Return acceptance rate for the path ensemble.

The acceptance rate is obtained as the fraction of accepted paths to the total number of paths in the path ensemble. This will only consider the paths that are currently stored in self.paths.

Returns:out (float) – The acceptance rate.
get_accepted()[source]

Yield accepted paths from the path ensemble.

This function will return an iterator useful for iterating over accepted paths only. In the path ensemble we store both accepted and rejected paths. This function will loop over all paths stored and yield the accepted paths the correct number of times.

get_paths()[source]

Yield the different paths stored in the path ensemble.

It is included here in order to have a simple compatibility between the PathEnsemble object and the py:class:.PathEnsembleFile object. This is useful for the analysis.

Yields:out (dict) – This is the dictionary representing the path data.
load_restart_info(info, cycle=0)[source]

Load restart information.

Parameters:
  • info (dict) – A dictionary with the restart information.
  • cycle (integer, optional) – The current simulation cycle.
move_path_to_generate(_path, _prefix=None)[source]

Move a path for temporary storing.

reset_data()[source]

Erase the stored data in the path ensemble.

It can be used in combination with flushing the data to a file in order to periodically write and empty the amount of data stored in memory.

Notes

We do not reset self.last_path as this might be used in the RETIS function.

restart_info()[source]

Return a dictionary with restart information.

store_path(path)[source]

Store a new accepted path in the path ensemble.

Parameters:path (object like PathBase) – The path we are going to store.
Returns:None, but we update self.last_path.
update_directories(path)[source]

Update directory names.

This method will not create new directories, but it will update the directory names.

Parameters:path (string) – The base path to set.
class pyretis.core.pathensemble.PathEnsembleExt(ensemble_number, interfaces, rgen=None, maxpath=10000, exe_dir=None)[source]

Bases: PathEnsemble

Representation of a path ensemble.

This class is similar to PathEnsemble but it is made to work with external paths. That is, some extra file handling is done when accepting a path.

static _copy_path(path, target_dir, prefix=None)[source]

Copy a path to a given target directory.

Parameters:
  • path (object like PathBase) – This is the path object we are going to copy.
  • target_dir (string) – The location where we are copying the path to.
  • prefix (string, optional) – To give a prefix to the name of copied files.
Returns:

out (object like py:class:.PathBase) – A copy of the input path.

static _move_path(path, target_dir, prefix=None)[source]

Move a path to a given target directory.

Parameters:
  • path (object like PathBase) – This is the path object we are going to move.
  • target_dir (string) – The location where we are moving the path to.
  • prefix (string, optional) – To give a prefix to the name of moved files.
copy_path_to_generate(path, pref=None)[source]

Copy a path for temporary storing.

list_superfluous()[source]

List files in accepted directory that we do not need.

load_restart_info(info, cycle=0)[source]

Load restart for external path.

move_path_to_generate(path, prefix=None)[source]

Move a path for temporary storing.

store_path(path)[source]

Store a path by explicitly moving it.

Parameters:path (object like PathBase) – This is the path object we are going to store.
pyretis.core.pathensemble._generate_file_names(path, target_dir, prefix=None)[source]

Generate new file names for moving copying paths.

Parameters:
  • path (object like PathBase) – This is the path object we are going to store.
  • target_dir (string) – The location where we are moving the path to.
  • prefix (string, optional) – The prefix can be used to prefix the name of the files.
Returns:

  • out[0] (list) – A list with new file names.
  • out[1] (dict) – A dict which defines the unique “source -> destination” for copy/move operations.

pyretis.core.pathensemble.generate_ensemble_name(ensemble_number, zero_pad=3)[source]

Generate a simple name for an ensemble.

The simple name will have a format like 01, 001, 0001 etc. and it is used to name the path ensemble and the output directory.

Parameters:
  • ensemble_number (int) – The number representing the ensemble.
  • zero_pad (int, optional) – The number of zeros to use for padding the name.
Returns:

out (string) – The ensemble name.

pyretis.core.pathensemble.get_path_ensemble_class(ensemble_type)[source]

Return the path ensemble class consistent with the given engine.

Parameters:ensemble_type (string) – The type of ensemble we are requesting.

pyretis.core.properties module

This file contains a class for a generic property.

class pyretis.core.properties.Property(desc='')[source]

Bases: object

A generic numerical value with standard deviation and average.

A generic object to store values obtained during a simulation. It will maintain the mean and variance as values are added using Property.add(val)

desc

Description of the property.

Type:string
nval

Number values stored.

Type:integer
mean

The current mean.

Type:float
delta2

Helper variable used for calculating the variance.

Type:float
variance

The current variance

Type:float
val

Store all values.

Type:list
Parameters:desc (string, optional) – Used to set the attribute desc.

Examples

>>> from pyretis.core.properties import Property
>>> ener = Property(desc='Energy of the particle(s)')
>>> ener.add(42.0)
>>> ener.add(12.220)
>>> ener.add(99.22)
>>> ener.mean
__init__(desc='')[source]

Initialise the property.

Parameters:desc (string, optional) – Description of the object.
__str__()[source]

Return string description of the property.

add(val)[source]

Add a value to the property & update the mean and variance.

Parameters:val (float or another type (list, numpy.array)) – The value to add.
Returns:out (None) – Returns None but updates the mean and variance.
dump_to_file(filename)[source]

Dump the contents in self.val to a file.

Parameters:filename (string) – Name/path of file to write.

Note

Consider if this should be moved/deleted and just replaced with a function from a more general input-output module

update_mean_and_variance()[source]

Calculate the mean and variance on the fly.

Source: http://en.wikipedia.org/wiki/Algorithms_for_calculating_variance

Returns:out (None) – Returns None but updates the mean and variance.

Note

Consider if this should be moved/deleted and just replaced with a function from the analysis package.

pyretis.core.random_gen module

This module handles random number generation.

It derives most of the random number procedures from RandomState in numpy.random and defines a class which used RandomState to generate pseudo-random numbers.

Important classes defined here

RandomGeneratorBase (RandomGeneratorBase)
An abstract base class for random number generators.
RandomGenerator (RandomGenerator)
A class representing a random number generator.
ReservoirSampler (ReservoirSampler)
A class for reservoir sampling.
MockRandomGenerator (MockRandomGenerator)
A non-random number generator which is useful for testing.
class pyretis.core.random_gen.MockRandomGenerator(seed=0, norm_shift=False)[source]

Bases: RandomGeneratorBase

A mock random generator, useful only for testing.

This class represents a random generator that can be used for testing algorithms. It will simply return numbers from a small list of pseudo-random numbers. This class is only useful for testing algorithms on different systems. It should NEVER be used for actual production runs!

__init__(seed=0, norm_shift=False)[source]

Initialise the mock random number generator.

Here, we set up some predefined random number which we will use as a pool for the generation.

Parameters:
  • seed (int, optional) – An integer used for seeding the generator if needed.
  • norm_shift (boolean) – If True is will ensure that the fake ‘normal distribution’ is shifted to get the right mean.
choice(a, size=None, replace=True, p=None)[source]

Choose random samples.

Parameters:
  • a (iterable) – The group to choose from.
  • size (int, optional) – The number of samples to choose.
  • replace (bool, optional) – If to replace samples between draws, default True.
  • p (iterable, optional) – The probabilities for every option in a, default is an uniform dirstribution.
Returns:

choice (array-like) – The picked choices.

get_state()[source]

Return current state.

multivariate_normal(mean, cov, cho=None, size=1)[source]

Draw numbers from a multi-variate distribution.

This is an attempt on speeding up the call of RandomState.multivariate_normal if we need to call it over and over again. Such repeated calling will do an SVD repeatedly, which is wasteful. In this function, this transform can be supplied and it is only estimated if it’s not explicitly given.

Parameters:
  • mean (numpy array (1D, 2)) – Mean of the N-dimensional array
  • cov (numpy array (2D, (2, 2))) – Covariance matrix of the distribution.
  • cho (numpy.array (2D, (2, 2)), optional) – Cholesky factorization of the covariance. If not given, it will be calculated here.
  • size (int, optional) – The number of samples to do.
Returns:

out (float or numpy.array of floats size) – The random numbers that are drawn here.

See also

numpy.random.multivariate_normal

normal(loc=0.0, scale=1.0, size=None)[source]

Return values from a normal distribution.

Parameters:
  • loc (float, optional) – The mean of the distribution
  • scale (float, optional) – The standard deviation of the distribution
  • size (int, tuple of ints, optional) – Output shape, i.e. how many values to generate. Default is None which is just a single value.
Returns:

out (float, numpy.array of floats) – The random numbers generated.

Note

This is part of the Mock-random number generator. Hence, it won’t provide a true normal distribution, though the mean is set to the value of loc.

rand(shape=1)[source]

Draw random numbers in [0, 1).

Parameters:shape (int) – The number of numbers to draw.
Returns:out (float) – Pseudo-random number in [0, 1).
random_integers(low, high)[source]

Return random integers in [low, high].

Parameters:
  • low (int) – This is the lower limit
  • high (int) – This is the upper limit
Returns:

out (int) – This is a pseudo-random integer in [low, high].

set_state(state)[source]

Set current state.

class pyretis.core.random_gen.RandomGenerator(seed=0)[source]

Bases: RandomGeneratorBase

A random number generator from numpy.

This class that defines a random number generator. It will use numpy.random.RandomState for the actual generation and we refer to the numpy documentation [1].

seed

A seed for the pseudo-random generator.

Type:int
rgen

This is a container for the Mersenne Twister pseudo-random number generator as implemented in numpy [13].

Type:object like numpy.random.RandomState

References

[13]The NumPy documentation on RandomState, http://docs.scipy.org/doc/numpy/reference/generated/numpy.random.RandomState.html
__init__(seed=0)[source]

Initialise the random number generator.

If a seed is given, the random number generator will be seeded.

Parameters:seed (int, optional) – An integer used for seeding the generator if needed.
choice(a, size=None, replace=True, p=None)[source]

Chooses random samples.

Parameters:
  • a (iterable, int) – The group to choose from, an int is converted to range(a).
  • size (int, optional) – The number of samples to choose.
  • replace (bool, optional) – If to replace samples between draws, default True.
  • p (iterable, optional) – The probabilities for every option in a, default is an uniform dirstribution.
Returns:

choice (array-like) – The picked choices.

get_state()[source]

Return current state.

multivariate_normal(mean, cov, cho=None, size=1)[source]

Draw numbers from a multi-variate distribution.

This is an attempt on speeding up the call of RandomState.multivariate_normal if we need to call it over and over again. Such repeated calling will do an SVD repeatedly, which is wasteful. In this function, this transform can be supplied and it is only estimated if it’s not explicitly given.

Parameters:
  • mean (numpy array (1D, 2)) – Mean of the N-dimensional array
  • cov (numpy array (2D, (2, 2))) – Covariance matrix of the distribution.
  • cho (numpy.array (2D, (2, 2)), optional) – Cholesky factorization of the covariance. If not given, it will be calculated here.
  • size (int, optional) – The number of samples to do.
Returns:

out (float or numpy.array of floats size) – The random numbers that are drawn here.

See also

numpy.random.multivariate_normal

normal(loc=0.0, scale=1.0, size=None)[source]

Return values from a normal distribution.

Parameters:
  • loc (float, optional) – The mean of the distribution
  • scale (float, optional) – The standard deviation of the distribution
  • size (int, tuple of ints, optional) – Output shape, i.e. how many values to generate. Default is None which is just a single value.
Returns:

out (float, numpy.array of floats) – The random numbers generated.

rand(shape=1)[source]

Draw random numbers in [0, 1).

Parameters:shape (int, optional) – The number of numbers to draw
Returns:out (float) – Pseudo-random number in [0, 1)

Note

Here, we will just draw a list of numbers and not for an arbitrary shape.

random_integers(low, high)[source]

Draw random integers in [low, high].

Parameters:
  • low (int) – This is the lower limit
  • high (int) – This is the upper limit
Returns:

out (int) – The pseudo-random integers in [low, high].

Note

np.random.randint(low, high) is defined as drawing from low (inclusive) to high (exclusive).

set_state(state)[source]

Set state for random generator.

class pyretis.core.random_gen.RandomGeneratorBase(seed=0)[source]

Bases: object

A base class for random number generators.

This is a base class for random number generators. It does not actually implement a generator.

seed

A seed for the generator

Type:int
__init__(seed=0)[source]

Initialise the random number generator.

Parameters:seed (int, optional) – An integer used for seeding the generator if needed.
draw_maxwellian_velocities(system, sigma_v=None)[source]

Draw numbers from a Gaussian distribution.

Parameters:
  • system (object like System) – This is used to determine the temperature parameter(s) and the shape (number of particles and dimensionality)
  • sigma_v (numpy.array, optional) – The standard deviation in velocity, one for each particle. If it’s not given it will be estimated.
generate_maxwellian_velocities(particles, boltzmann, temperature, dof, selection=None, momentum=True)[source]

Generate velocities from a Maxwell distribution.

The velocities are drawn to match a given temperature and this function can be applied to a subset of the particles.

The generation is done in three steps:

  1. We generate velocities from a standard normal distribution.
  2. We scale the velocity of particle i with 1.0/sqrt(mass_i) and reset the momentum.
  3. We scale the velocities to the set temperature.
Parameters:
  • particles (object like Particles) – These are the particles to set the velocity of.
  • boltzmann (float) – The Boltzmann factor in correct units.
  • temperature (float) – The desired temperature. Typically, system.temperature[‘set’] will be used here.
  • dof (list of floats, optional) – The degrees of freedom to subtract. Its shape should be equal to the number of dimensions.
  • selection (list of ints, optional) – A list with indices of the particles to consider. Can be used to only apply it to a selection of particles
  • momentum (boolean, optional) – If true, we will reset the momentum.
Returns:

out (None) – Returns None but modifies velocities of the selected particles.

abstract get_state()[source]

Return info about the current state.

abstract multivariate_normal(mean, cov, cho=None, size=1)[source]

Draw numbers from a multi-variate distribution.

Parameters:
  • mean (numpy array (1D, 2)) – Mean of the N-dimensional array
  • cov (numpy array (2D, (2, 2))) – Covariance matrix of the distribution.
  • cho (numpy.array (2D, (2, 2)), optional) – Cholesky factorization of the covariance. If not given, it will be calculated here.
  • size (int, optional) – The number of samples to do.
Returns:

out (float or numpy.array of floats size) – The random numbers that are drawn here.

abstract normal(loc=0.0, scale=1.0, size=None)[source]

Return values from a normal distribution.

Parameters:
  • loc (float, optional) – The mean of the distribution
  • scale (float, optional) – The standard deviation of the distribution
  • size (int, tuple of ints, optional) – Output shape, i.e. how many values to generate. Default is None which is just a single value.
Returns:

out (float, numpy.array of floats) – The random numbers generated.

abstract rand(shape=1)[source]

Draw random numbers in [0, 1).

Parameters:shape (int, optional) – The number of numbers to draw
Returns:out (float) – Pseudo-random number in [0, 1)
abstract random_integers(low, high)[source]

Draw random integers in [low, high].

Parameters:
  • low (int) – This is the lower limit
  • high (int) – This is the upper limit
Returns:

out (int) – The pseudo-random integers in [low, high].

abstract set_state(state)[source]

Set state for random generator.

class pyretis.core.random_gen.ReservoirSampler(seed=0, length=10, rgen=None)[source]

Bases: object

A class representing a reservoir sampler.

The reservoir sampler will maintain a list of k items drawn randomly from a set of N > k items. The list is created and maintained so that we only need to store k`items This is useful when `N is very large or when storing all N items require a lot of memory. The algorithm is described by Knuth [14] but here we do a variation, so that each item may be picked several times.

rgen

This is a container for the Mersenne Twister pseudo-random number generator as implemented in numpy, see the documentation of RandomGenerator.

Type:object like numpy.random.RandomState
items

The number of items seen so far, i.e. the current N.

Type:integer
reservoir

The items we have stored.

Type:list
length

The maximum number of items to store in the reservoir (i.e. k).

Type:integer
ret_idx

This is the index of the item to return if we are requesting items from the reservoir.

Type:integer

References

[14]The Art of Computer Programming.
__init__(seed=0, length=10, rgen=None)[source]

Initialise the reservoir.

Parameters:
  • seed (int, optional) – An integer used for seeding the generator.
  • length (int, optional) – The maximum number of items to store.
  • rgen (object like RandomGenerator, optional) – In case we want to re-use a random generator object. If this is specified, the parameter seed is ignored.
append(new_item)[source]

Try to add an item to the reservoir.

Parameters:new_item (any type) – This is the item we try to add to the reservoir.
get_item()[source]

Return the next item from the reservoir.

Returns:out (any type) – Returns an item from the reservoir.
pyretis.core.random_gen.create_random_generator(settings=None)[source]

Create a random generator from given settings.

Parameters:settings (dict, optional) – This is the dict used for creating the random generator. Currently, we will actually just look for a seed value.
Returns:out (object like RandomGenerator) – The random generator created.

pyretis.core.retis module

This module contains functions for RETIS.

This module defines functions that are needed to perform Replica Exchange Transition Interface Sampling (RETIS). The RETIS algorithm was first described by van Erp [RETIS].

Methods defined here

make_retis_step (make_retis_step())
Function to select and execute the RETIS move.
retis_tis_moves (retis_tis_moves())
Function to execute the TIS steps in the RETIS algorithm.
retis_moves (retis_moves())
Function to perform RETIS swapping moves - it selects what scheme to use, i.e. [0^-] <-> [0^+], [1^+] <-> [2^+], ... or [0^+] <-> [1^+], [2^+] <-> [3^+], ....
retis_swap_wrapper (retis_swap_wrapper())
The function that actually swaps two path ensembles. This function decides which swap will be done: it can be a typical RETIS swap (retis_swap), or a swap between two PPTIS ensembles with limited memory (prretis_swap).
retis_swap (retis_swap())
The function that actually swaps two path ensembles. A swap between two typical RETIS ensembles, not PPTIS ensembles.
repptis_swap (repptis_swap())
The function that actually swaps two path ensembles. A swap between two PPTIS ensembles with truncated paths, so this is not the typical RETIS swap.
retis_swap_zero (retis_swap_zero())
The function that performs the swapping for the [0^-] <-> [0^+] swap.
high_acc_swap (high_acc_wap())
The function coputes if a path generated via SS can be accepted for swapping in accordance to super detail balance.

References

[RETIS]T. S. van Erp, Phys. Rev. Lett. 98, 26830 (2007), http://dx.doi.org/10.1103/PhysRevLett.98.268301
pyretis.core.retis.make_retis_step(ensembles, rgen, settings, cycle)[source]

Determine and execute the appropriate RETIS move.

Here we will determine what kind of RETIS moves we should do. We have two options:

  1. Do the RETIS swapping moves. This is done by calling retis_moves().
  2. Do TIS moves, either for all ensembles or for just one, based on values of relative shoot frequencies. This is done by calling retis_tis_moves().

This function will just determine and execute the appropriate move (1 or 2) based on the given swapping frequencies in the settings and drawing a random number from the random number generator rgen.

Parameters:
  • ensembles (list of dicts of objects) – This is a list of the ensembles we are using in the RETIS method. It contains:
    • path_ensemble: object like PathEnsemble This is used for storing results for the simulation. It is also used for defining the interfaces for this simulation.
    • system: object like System System is used here since we need access to the temperature and to the particle list.
    • order_function: object like OrderParameter The class used for calculating the order parameter(s).
    • engine: object like EngineBase The engine to use for propagating a path.
  • rgen (object like RandomGenerator) – This is a random generator. Here we assume that we can call rgen.rand() to draw random uniform numbers.
  • settings (dict) – Contains the settings for the simulation.
  • cycle (integer) – The current cycle number.
Returns:

out (generator) – This generator yields the results after performing the RETIS moves.

pyretis.core.retis.null_move(ensemble, cycle)[source]

Perform a null move for a path ensemble.

The null move simply consist of accepting the last accepted path again.

Parameters:
  • ensemble (dict, it contains:) –

    • path_ensembleobject like PathEnsemble

      This is the path ensemble to update with the null move

  • cycle (integer) – The current cycle number

Returns:

  • out[0] (boolean) – Should the path be accepted or not? Here, it’s always True since the null move is always accepted.
  • out[1] (object like PathBase) – The unchanged path.
  • out[2] (string) – The status will here be ‘ACC’ since we just accept the last accepted path again in this move.

pyretis.core.retis.retis_moves(ensembles, rgen, settings, cycle)[source]

Perform RETIS moves on the given ensembles.

This function will perform RETIS moves on the given ensembles. First we have two strategies based on settings[‘retis’][‘swapsimul’]:

  1. If settings[‘retis’][‘swapsimul’] is True we will perform several swaps, either [0^-] <-> [0^+], [1^+] <-> [2^+], ... or [0^+] <-> [1^+], [2^+] <-> [3^+], .... Which one of these two swap options we use is determined randomly and they have equal probability.
  2. If settings[‘retis’][‘swapsimul’] is False we will just perform one swap for randomly chosen ensembles, i.e. we pick a random ensemble and try to swap with the ensemble to the right. Here we may also perform null moves if the settings[‘retis’][‘nullmove’] specifies so.
Parameters:
  • ensembles (list of dicts of objects) – Each contains:
    • path_ensemble: object like PathEnsemble This is used for storing results for the simulation. It is also used for defining the interfaces for this simulation.
    • system: object like System System is used here since we need access to the temperature and to the particle list
    • order_function: object like OrderParameter The class used for calculating the order parameter(s).
    • engine: object like EngineBase The engine to use for propagating a path.
  • rgen (object like RandomGenerator) – This is a random generator. Here we assume that we can call rgen.rand() to draw random uniform numbers.
  • settings (dict) – This dict contains the settings for the RETIS method.
  • cycle (integer) – The current cycle number.
Yields:

out (dict) – This dictionary contains the result of the RETIS moves.

pyretis.core.retis.retis_swap(ensembles, idx, settings, cycle)[source]

Perform a RETIS swapping move for two ensembles.

These ensembles are not PPTIS ensembles.

The RETIS swapping move will attempt to swap accepted paths between two ensembles in the hope that path from [i^+] is an acceptable path for [(i+1)^+] as well. We have two cases:

  1. If we try to swap between [0^-] and [0^+] we need to integrate the equations of motion.
  2. Otherwise, we can just swap and accept if the path from [i^+] is an acceptable path for [(i+1)^+]. The path from [(i+1)^+] is always acceptable for [i^+] (by construction).
Parameters:
  • ensembles (list of dictionaries of objects) – This is a list of the ensembles we are using in the RETIS method. Each one contains:
    • path_ensemble: object like PathEnsemble This is used for storing results for the simulation. It is also used for defining the interfaces for this simulation.
    • system: object like System System is used here since we need access to the temperature and to the particle list
    • order_function: object like OrderParameter The class used for calculating the order parameter(s).
    • engine: object like EngineBase The engine to use for propagating a path.
  • idx (integer) – Definition of what path ensembles to swap. We will swap ensembles[idx] with ensembles[idx+1]. If idx == 0 we have case 1) defined above.
  • settings (dict) – This dict contains the settings for the RETIS method.
  • cycle (integer) – Current cycle number.
Returns:

  • out[0] (boolean) – Should the path be accepted or not?
  • out[1] (list of object like PathBase) – The trial paths.
  • out[2] (string) – The status for the trial paths.

Note

Note that path.generated is NOT updated here. This is because we are just swapping references and not the paths. In case the swap is rejected updating this would invalidate the last accepted path.

pyretis.core.retis.retis_swap_zero(ensembles, settings, cycle)[source]

Perform the RETIS swapping for [0^-] <-> [0^+] swaps.

The RETIS swapping move for ensembles [0^-] and [0^+] requires some extra integration. Here we are generating new paths for [0^-] and [0^+] in the following way:

  1. For [0^-] we take the initial point in [0^+] and integrate backward in time. This is merged with the second point in [0^+] to give the final path. The initial point in [0^+] starts to the left of the interface and the second point is on the right side - i.e. the path will cross the interface at the end points. If we let the last point in [0^+] be called A_0 and the second last point B, and we let A_1, A_2, ... be the points on the backward trajectory generated from A_0 then the final path will be made up of the points [..., A_2, A_1, A_0, B]. Here, B will be on the right side of the interface and the first point of the path will also be on the right side.
  2. For [0^+] we take the last point of [0^-] and use that as an initial point to generate a new trajectory for [0^+] by integration forward in time. We also include the second last point of the [0^-] trajectory which is on the left side of the interface. We let the second last point be B (this is on the left side of the interface), the last point A_0 and the points generated from A_0 we denote by A_1, A_2, .... Then the resulting path will be [B, A_0, A_1, A_2, ...]. Here, B will be on the left side of the interface and the last point of the path will also be on the left side of the interface.
Parameters:
  • ensembles (list of dictionaries of objects) – This is a list of the ensembles we are using in the RETIS method. It contains:
    • path_ensemble: object like PathEnsemble This is used for storing results for the simulation. It is also used for defining the interfaces for this simulation.
    • system: object like System System is used here since we need access to the temperature and to the particle list.
    • order_function: object like OrderParameter The class used for calculating the order parameter(s).
    • engine: object like EngineBase The engine to use for propagating a path.
  • settings (dict) – This dict contains the settings for the RETIS method.
  • cycle (integer) – The current cycle number.
Returns:

out (string) – The result of the swapping move.

pyretis.core.system module

Module defining the system class.

The system class is used to group together many important objects in PyRETIS, for instance, the particles, force field etc.

Important classes defined here

System (System)
A class representing a system. A system object defines the system the simulation acts and contains information about particles, force fields etc.
class pyretis.core.system.System(units='lj', box=None, temperature=None)[source]

Bases: object

This class defines a generic system for simulations.

box

Defines the simulation box.

Type:object like Box
temperature

This dictionary contains information on the temperature. The following information is stored:

  • set: The set temperature, T, (if any).
  • beta: The derived property 1.0/(k_B*T).
  • dof: Information about the degrees of freedom for the system.
Type:dict
order

The order parameter(s) for the current state of the system (if they have been calculated).

Type:tuple
particles

Defines the particle list which represents the particles and the properties of the particles (positions, velocities, forces etc.).

Type:object like Particles
post_setup

This list contains extra functions that should be called when preparing to run a simulation. This is typically functions that should only be called after the system is fully set up. The tuples should correspond to (‘function’, args) where such that system.function(*args) can be called.

Type:list of tuples
forcefield

Defines the force field to use and implements the actual force and potential calculation.

Type:object like ForceField
units

Units to use for the system/simulation. Should match the defined units in pyretis.core.units.

Type:string
__eq__(other)[source]

Compare two system objects.

__init__(units='lj', box=None, temperature=None)[source]

Initialise the system.

Parameters:
  • units (string, optional) – The system of units to use in the simulation box.
  • box (object like Box, optional) – This variable represents the simulation box. It is used to determine the number of dimensions.
  • temperature (float, optional) – The (desired) temperature of the system, if applicable.

Note

self.temperature is defined as a dictionary. This is just because it’s convenient to include information about the degrees of freedom of the system here.

__ne__(other)[source]

Check if two systems are not equal.

__str__()[source]

Just print some basic info about the system.

_adjust_dof_according_to_box()[source]

Adjust the dof according to the box connected to the system.

For each ‘True’ in the periodic settings of the box, we subtract one degree of freedom for that dimension.

add_particle(pos, vel=None, force=None, mass=1.0, name='?', ptype=0)[source]

Add a particle to the system.

Parameters:
  • pos (numpy.array,) – Position of the particle.
  • vel (numpy.array, optional) – The velocity of the particle. If not given numpy.zeros will be used.
  • force (numpy.array, optional) – Force on the particle. If not given np.zeros will be used.
  • mass (float, optional) – Mass of the particle, the default is 1.0.
  • name (string, optional) – Name of the particle, the default is ‘?’.
  • ptype (integer, optional) – Particle type, the default is 0.
Returns:

out (None) – Does not return anything, but updates particles.

adjust_dof(dof)[source]

Adjust the degrees of freedom to neglect in the system.

Parameters:dof (numpy.array) – The degrees of freedom to neglect, in addition to the ones we already have neglected.
calculate_beta(temperature=None)[source]

Return the so-called beta factor for the system.

Beta is defined as \beta = 1/(k_\text{B} \times T) where k_\text{B} is the Boltzmann constant and the temperature T is either specified in the parameters or assumed equal to the set temperature of the system.

Parameters:temperature (float, optional) – The temperature of the system. If the temperature is not given, self.temperature will be used.
Returns:out (float) – The calculated beta factor, or None if no temperature data is available.
calculate_temperature()[source]

Calculate the temperature of the system.

It is included here for convenience since the degrees of freedom are easily accessible and it’s a very common calculation to perform, even though it might be cleaner to include it as a particle function.

Returns:out (float) – The temperature of the system.
copy()[source]

Return a copy of the system.

This copy is useful for storing snapshots obtained during a simulation.

Returns:out (object like System) – A copy of the system.
evaluate_force()[source]

Evaluate forces on the particles.

Returns:
  • out[1] (numpy.array) – Forces on the particles.
  • out[2] (numpy.array) – The virial.

Note

This function will not update the forces, just calculate them. Use self.force to update the forces.

evaluate_potential()[source]

Evaluate the potential energy.

Returns:out (float) – The potential energy.

Note

This function will not update the potential, but it will just return its value for the (possibly given) configuration. The function self.potential can be used to update the potential for the particles in the system.

evaluate_potential_and_force()[source]

Evaluate the potential and/or the force.

Returns:
  • out[1] (float) – The potential energy.
  • out[2] (numpy.array) – Forces on the particles.
  • out[3] (numpy.array) – The virial.

Note

This function will not update the forces/potential energy for the particles. To update these, call self.potential_and_force.

extra_setup()[source]

Perform extra set-up for the system.

The extra set-up will typically be tasks that can only be performed after the system is fully set-up, for instance after the force field is properly defined.

force()[source]

Update the forces and the virial.

The update is done by calling self._evaluate_potential_force.

Returns:
  • out[1] (numpy.array) – Forces on the particles. Note that self.particles.force will also be updated.
  • out[2] (numpy.array) – The virial. Note that self.particles.virial will be updated.
generate_velocities(rgen=None, seed=0, momentum=True, temperature=None, distribution='maxwell')[source]

Set velocities for the particles according to a given temperature.

The temperature can be specified, or it can be taken from self.temperature[‘set’].

Parameters:
  • rgen (string, optional) – This string can be used to select a particular random generator. Typically this is only useful for testing.
  • seed (int, optional) – The seed for the random generator.
  • momentum (boolean, optional) – Determines if the momentum should be reset.
  • temperature (float, optional) – The desired temperature to set.
  • distribution (str, optional) – Selects a distribution for generating the velocities.
Returns:

out (None) – Does not return anything, but updates system.particles.vel.

get_boltzmann()[source]

Return the Boltzmann constant in correct units for the system.

Returns:out (float) – The Boltzmann constant.
get_dim()[source]

Return the dimensionality of the system.

The value is obtained from the box. In other words, it is the box object that defines the dimensionality of the system.

Returns:out (integer) – The number of dimensions of the system.
load_restart_info(info)[source]

Load restart information.

Parameters:info (dict) – The dictionary with the restart information, should be similar to the dict produced by restart_info().
potential()[source]

Update the potential energy.

Returns:out (float) – The potential energy.
potential_and_force()[source]

Update the potential energy and forces.

The potential in self.particles.vpot and the forces in self.particles.force are here updated by calling forcefield.evaluate_potential_force().

Returns:
  • out[1] (float) – The potential energy, note self.particles.vpot is also updated.
  • out[2] (numpy.array) – Forces on the particles. Note that self.particles.force will also be updated.
  • out[3] (numpy.array) – The virial. Note that self.particles.virial will also be updated.
rescale_velocities(energy, external=False)[source]

Re-scale the kinetic energy to a given total energy.

Parameters:
  • energy (float) – The desired energy.
  • energy (boolean, optional) – If True, self.particles.vpot will be used as the potential energy.
Returns:

None, but updates the velocities of the particles.

restart_info()[source]

Return a dictionary with restart information.

update_box(length)[source]

Update the system box, create if needed.

Parameters:length (numpy.array, list or iterable.) – The box vectors represented as a list.

pyretis.core.tis module

This file contains functions used in TIS.

This module defines the functions needed to perform TIS simulations (algorithms implemented as described by van Erp et al. [TIS] ) with advanced shooting moves, such as Stone Skipping and Web Throwing (algorithms implemented as described by Riccardi et al. [SS+WT] ).

Methods defined here

make_tis_step (make_tis_step())
A method that will perform a single TIS step.
make_tis_step_ensemble (make_tis_step_ensemble())
A method to perform a TIS step for a path ensemble. It will handle adding of the path to a path ensemble object.
shoot (shoot())
A method that will perform a shooting move.
select_shoot (select_shoot())
A method that randomly selects the shooting point.
wire_fencing (wire_fencing())
A method that will do a wire fencing shooting move.
stone_skipping (stone_skipping())
A method that will do a stone skipping shooting move.
web_throwing (web_throwing())
A method that will do a web throwing shooting move.
extender (extender())
A method that will extend a path to target interfaces.
time_reversal (time_reversal())
A method for performing the time reversal move.
ss_wt_wf_acceptance (ss_wt_wf_acceptance())
A method to check the acceptance of a newly generated path according to super detailed balance rules for Stone Skipping, Wire Fencing (basic and High Acceptance version) and Web Throwing.

References

[SS+WT]E. Riccardi, O. Dahlen, T. S. van Erp, J. Phys. Chem. letters, 8, 18, 4456, (2017), https://doi.org/10.1021/acs.jpclett.7b01617
[TIS]T. S. van Erp, D. Moroni, P. G. Bolhuis, J. Chem. Phys. 118, 7762 (2003), https://dx.doi.org/10.1063%2F1.1562614
pyretis.core.tis.extender(source_seg, ensemble, tis_set, start_cond=('R', 'L'))[source]

Extend a path to the given interfaces.

This function will perform the web throwing move from an initial path.

Parameters:
  • source_seg (object like PathBase) – This is the input path which will be prolonged.
  • ensemble (dictionary of objects) – It contains:
    • order_function: object like OrderParameter The class used for obtaining the order parameter(s).
    • engine: object like EngineBase The engine to use for propagating a path.
    • system: object like System System is used here since we need access to the temperature and to the particle list.
    • interfaces: list of floats These are the interface positions on form [left, middle, right].
  • tis_set (dict) – This contains the settings for TIS. Keys used here:
    • aimless: boolean, is the shooting aimless or not?
    • allowmaxlength: boolean, should paths be allowed to reach maximum length?
    • maxlength: integer, maximum allowed length of paths.
    • high_accept: boolean, the option for High Acceptance WF and SS.
  • start_cond (string, optional) – The starting condition for the current ensemble, ‘L’eft or ‘R’ight. (‘L’, ‘R’), the default option, implies no directional difference.
Returns:

  • out[0] (boolean) – True if the path can be accepted.
  • out[1] (object like PathBase) – Returns the generated path.
  • out[2] (string) – Status of the path, this is one of the strings defined in path._STATUS.

pyretis.core.tis.make_tis(ensembles, rgen, settings, cycle)[source]

Perform a TIS step in a given path ensemble.

This function will execute the TIS steps. When used in the RETIS method, three different options can be given:

  1. If relative_shoots is given in the input settings as a list of probability for each ensemble to be pick, then an ensemble will be ranomly picked and TIS performed on. For all the other ensembles we again have two options based on the given settings[‘nullmoves’]:

    1. Do a ‘null move’ in all other ensembles.
    2. Do nothing for all other ensembles.

    Performing the null move in an ensemble will simply just accept the previously accepted path in that ensemble again.

  2. If relative_shoots is not given in the input settings, TIS moves will be executed for all path ensembles.

Parameters:
  • ensembles (list of dictionaries of objects) – This is a list of the ensembles we are using in the TIS method with their properties. They contain:
    • path_ensemble: object like PathEnsemble This is used for storing results for the simulation. It is also used for defining the interfaces for this simulation.
    • system: object like System System is used here since we need access to the temperature and to the particle list
    • order_function: object like OrderParameter The class used for calculating the order parameter(s).
    • engine: object like EngineBase The engine to use for propagating a path.
  • rgen (object like RandomGenerator) – This is the random generator that will be used.
  • settings (dict) – This dictionary contains the TIS settings.
  • cycle (int) – The current cycle number.
Yields:

out (dict) – This dictionary contains the results of the TIS move(s).

pyretis.core.tis.make_tis_step(ensemble, tis_settings, start_cond)[source]

Perform a TIS step and generate a new path.

The new path will be generated from the input path, either by performing a time-reversal move or by a shooting move. This is determined pseudo-randomly by drawing a random number from a uniform distribution using the given random generator.

Parameters:
  • ensemble (dictionary of objects) – It contains:
    • path_ensemble: object like PathEnsemble This is the path ensemble to perform the TIS step for.
    • path: object like PathBase This is the input path which will be used for generating a new path.
    • system: object like System System is used here since we need access to the temperature and to the particle list.
    • order_function: object like OrderParameter The class used for obtaining the order parameter(s).
    • engine: object like EngineBase The engine to use for propagating a path.
    • rgen: object like RandomGenerator This is the random generator that will be used.
    • interfaces: list of floats These are the interface positions on form [left, middle, right].
  • tis_settings (dict) – This dictionary contains the settings for the TIS method. Here we explicitly use:
    • freq: float, the frequency of how often we should do time reversal moves.
    • shooting_move: string, the label of the shooting move to perform.
  • start_cond (string) – The starting condition for the path. This is determined by the ensemble we are generating for - it is ‘R’ight or ‘L’eft.
Returns:

  • out[0] (boolean) – True if the new path can be accepted.
  • out[1] (object like PathBase) – The generated path.
  • out[2] (string) – The status of the path.

pyretis.core.tis.make_tis_step_ensemble(ensemble, settings, cycle)[source]

Perform a TIS step for a given path ensemble.

This function will run make_tis_step for the given path_ensemble. It will handle adding of the path. This function is intended for convenience when working with path ensembles. If we are using the path ensemble [0^-] then the start condition should be ‘R’ for right.

Parameters:
  • ensemble (dictionary of objects) – It contains:
    • path_ensemble: object like PathEnsemble This is used for storing results for the simulation. It is also used for defining the interfaces for this simulation.
    • system: object like System System is used here since we need access to the temperature and to the particle list.
    • order_function: object like OrderParameter The class used for obtaining the order parameter(s).
    • engine: object like EngineBase The engine to use for propagating a path.
    • rgen: object like RandomGenerator This is the random generator that will be used.
  • settings (dict) – This dictionary contains the ensemble settings.
  • cycle (int) – The current cycle number.
Returns:

  • out[0] (boolean) – True if the new path can be accepted, False otherwise.
  • out[1] (object like PathBase) – The generated path.
  • out[2] (string) – The status of the path.

pyretis.core.tis.paste_paths(path_back, path_forw, overlap=True, maxlen=None)[source]

Merge a backward with a forward path into a new path.

The resulting path is equal to the two paths stacked, in correct time. Note that the ordering is important here so that: paste_paths(path1, path2) != paste_paths(path2, path1).

There are two things we need to take care of here:

  • path_back must be iterated in reverse (it is assumed to be a backward trajectory).
  • we may have to remove one point in path2 (if the paths overlap).
Parameters:
  • path_back (object like PathBase) – This is the backward trajectory.
  • path_forw (object like PathBase) – This is the forward trajectory.
  • overlap (boolean, optional) – If True, path_back and path_forw have a common starting-point, that is, the first point in path_forw is identical to the first point in path_back. In time-space, this means that the first point in path_forw is identical to the last point in path_back (the backward and forward path started at the same location in space).
  • maxlen (float, optional) – This is the maximum length for the new path. If it’s not given, it will just be set to the largest of the maxlen of the two given paths.

Note

Some information about the path will not be set here. This must be set elsewhere. This includes how the path was generated (path.generated) and the status of the path (path.status).

pyretis.core.tis.select_shoot(ensemble, tis_settings, start_cond)[source]

Select the shooting move to generate a new path.

The new path will be generated from the input path, either by performing a normal shooting or web-throwing. This is determined pseudo-randomly by drawing a random number from a uniform distribution using the given random generator.

Parameters:
  • ensemble (dictionary of objects) – It contains:
    • path_ensemble: object like PathEnsemble This is the path ensemble to perform the TIS step for.
    • path: object like PathBase This is the input path which will be used for generating a new path.
    • system: object like System System is used here since we need access to the temperature and to the particle list.
    • order_function: object like OrderParameter The class used for obtaining the order parameter(s).
    • engine: object like EngineBase The engine to use for propagating a path.
    • rgen: object like RandomGenerator This is the random generator that will be used.
    • interfaces: list of floats These are the interface positions on form [left, middle, right].
  • tis_settings (dict) – This dictionary contains the settings for the TIS method. Here we explicitly use:
    • freq: float, the frequency of how often we should do time reversal moves.
    • shooting_move: string, the label of the shooting move to perform.
  • start_cond (string) – The starting condition for the path. This is determined by the ensemble we are generating for - it is ‘R’ight or ‘L’eft.
Returns:

  • out[0] (boolean) – True if the new path can be accepted.
  • out[1] (object like PathBase) – The generated path.
  • out[2] (string) – The status of the path.

pyretis.core.tis.shoot(ensemble, tis_settings, start_cond, shooting_point=None)[source]

Perform a shooting-move.

This function will perform the shooting move from a randomly selected time-slice.

Parameters:
  • ensemble (dict) – It contains:
    • path_ensemble: object like PathEnsemble This is the path ensemble to perform the TIS step for.
    • system: object like System System is used here since we need access to the temperature and to the particle list.
    • order_function: object like OrderParameter The class used for obtaining the order parameter(s).
    • engine: object like EngineBase The engine to use for propagating a path.
    • rgen: object like RandomGenerator This is the random generator that will be used.
    • interfaces: list of floats These are the interface positions on form [left, middle, right].
  • tis_settings (dict) – This contains the settings for TIS:
    • aimless: boolean, is the shooting aimless or not?
    • allowmaxlength: boolean, should paths be allowed to reach maximum length?
    • maxlength: integer, maximum allowed length of paths.
  • start_cond (string) – The starting condition for the current ensemble, ‘L’eft or ‘R’ight.
  • shooting_point (object like System, optional) – If given, it is the shooting point from which the path is generated.
Returns:

  • out[0] (boolean) – True if the path can be accepted.
  • out[1] (object like PathBase) – Returns the generated path.
  • out[2] (string) – Status of the path, this is one of the strings defined in path._STATUS.

pyretis.core.tis.ss_wt_wf_acceptance(trial_path, ensemble, tis_settings, start_cond='L')[source]

Weights, possibly reverses and accept/rejects generated SS/WT/WFpaths.

Parameters:
  • trial_path (object like PathBase) – This is the new path that will obtain weights, and might be reversed and accepted.
  • ensemble (dict) – It contains:
    • path_ensemble: object like PathEnsemble This is the path ensemble to perform the TIS step for.
    • order_function: object like OrderParameter The class used for obtaining the order parameter(s).
    • rgen: object like RandomGenerator This is the random generator that will be used.
    • interfaces : list/tuple of floats These are the interface positions of the form [left, middle, right].
  • tis_settings (dict) – This contains the settings for TIS. KEys used here;
    • high_accept : boolean, the option for High Acceptance SS/WF.
    • shooting_move : string, the label of the shooting move to perform.
  • start_cond (string, optional) – The starting condition for the current ensemble, ‘L’eft or ‘R’ight.
Returns:

  • out[0] (boolean) – True if the path can be accepted.
  • out[1] (object like PathBase) – Returns the weighed and possibly reversed path.

pyretis.core.tis.stone_skipping(ensemble, tis_settings, start_cond)[source]

Perform a stone_skipping move.

This function will perform the famous stone skipping move from an initial path.

Parameters:
  • ensemble (dict) – It contains:
    • path_ensemble: object like PathEnsemble This is the path ensemble to perform the TIS step for.
    • path: object like PathBase This is the input path which will be used for generating a new path.
    • system: object like System System is used here since we need access to the temperature and to the particle list.
    • order_function: object like OrderParameter The class used for obtaining the order parameter(s).
    • engine: object like EngineBase The engine to use for propagating a path.
    • rgen: object like RandomGenerator This is the random generator that will be used.
    • interfaces: list of floats These are the interface positions on form [left, middle, right].
  • tis_settings (dict) – This contains the settings for TIS. Keys used here:
    • aimless: boolean, is the shooting aimless or not?
    • allowmaxlength: boolean, should paths be allowed to reach maximum length?
    • maxlength: integer, maximum allowed length of paths.
    • high_accept: boolean, the option for High Acceptance SS.
  • start_cond (string) – The starting condition for the current ensemble, ‘L’eft or ‘R’ight.
Returns:

  • out[0] (boolean) – True if the path can be accepted.
  • out[1] (object like PathBase) – Returns the generated path.
  • out[2] (string) – Status of the path, this is one of the strings defined in path._STATUS.

pyretis.core.tis.time_reversal(path, order_function, interfaces, start_condition)[source]

Perform a time-reversal move.

Parameters:
  • path (object like PathBase) – This is the input path which will be used for generating a new path.
  • order_function (object like OrderParameter) – The class used for obtaining the order parameter(s).
  • interfaces (list/tuple of floats) – These are the interface positions of the form [left, middle, right].
  • start_condition (string) – The starting condition, ‘L’eft or ‘R’ight.
Returns:

  • out[0] (boolean) – True if the path can be accepted.
  • out[1] (object like PathBase) – The time reversed path.
  • out[2] (string) – Status of the path, this is one of the strings defined in .path._STATUS.

pyretis.core.tis.web_throwing(ensemble, tis_set, start_cond='L')[source]

Perform a web_throwing move.

This function performs the great web throwing move from an initial path.

Parameters:
  • ensemble (dict) – It contains:
    • path_ensemble: object like PathEnsemble This is the path ensemble to perform the TIS step for.
    • system: object like System System is used here since we need access to the temperature and to the particle list.
    • order_function: object like OrderParameter The class used for obtaining the order parameter(s).
    • engine: object like EngineBase The engine to use for propagating a path.
    • rgen: object like RandomGenerator This is the random generator that will be used.
    • interfaces: list of floats These are the interface positions on form [left, middle, right].
  • tis_set (dict) – This contains the settings for TIS. Keys used here:
    • aimless: boolean, is the shooting aimless or not?
    • allowmaxlength: boolean, should paths be allowed to reach maximum length?
    • maxlength: integer, maximum allowed length of paths.
  • start_cond (string, optional) – The starting condition for the current ensemble, ‘L’eft or ‘R’ight.
Returns:

  • out[0] (boolean) – True if the path can be accepted.
  • out[1] (object like PathBase) – Returns the generated path.
  • out[2] (string) – Status of the path, this is one of the strings defined in path._STATUS.

pyretis.core.tis.wire_fencing(ensemble, tis_settings, start_cond)[source]

Perform a wire_fencing move.

This function will perform the non famous wire fencing move from an initial path.

Parameters:
  • ensemble (dict) – It contains:
    • path_ensemble: object like PathEnsemble This is the path ensemble to perform the TIS step for.
    • system: object like System System is used here since we need access to the temperature and to the particle list.
    • order_function: object like OrderParameter The class used for obtaining the order parameter(s).
    • engine: object like EngineBase The engine to use for propagating a path.
    • rgen: object like RandomGenerator This is the random generator that will be used.
    • interfaces: list of floats These are the interface positions on form [left, middle, right].
  • tis_settings (dict) – This contains the settings for TIS. Keys used here:
    • aimless: boolean, is the shooting aimless or not?
    • allowmaxlength: boolean, should paths be allowed to reach maximum length?
    • maxlength: integer, maximum allowed length of paths.
    • high_accept: boolean, the option for High Acceptance WF.
  • start_cond (string) – The starting condition for the current ensemble, ‘L’eft or ‘R’ight.
Returns:

  • out[0] (boolean) – True if the path can be accepted.
  • out[1] (object like PathBase) – Returns the generated path.
  • out[2] (string) – Status of the path, this is one of the strings defined in path._STATUS.

pyretis.core.units module

This module defines natural constants and unit conversions.

This module defines some natural constants and conversions between units which can be used by the PyRETIS program. The natural constants are mainly used for conversions but it is also used to define the Boltzmann constant for internal use in PyRETIS. The unit conversions are mainly useful for input/output.

All numerical values are from the National Institute of Standards and Technology and can be accessed through a web interface http://physics.nist.gov/constants or in plain text. [7]

Internally, all computations are carried out in units which are defined by a length scale, an energy scale and a mass scale. The time scale is defined by these choices.

Charges are typically given (in the input) in units of the electron charge. The internal unit for charge is not yet implemented, but a choice is here to include the factor \frac{1}{\sqrt{4\pi\varepsilon_0}}. An internal calculation of q_1 q_2 will then include coulombs constant in the correct units.

The different sets of unit systems are described below in the section on unit systems.

Natural constants

The keys for CONSTANTS defines the natural constant and its units, for instance CONSTANTS['kB']['J/K'] is the Boltzmann constants in units of Joule per Kelvin. The currently defined natural constants are:

  • kB : The Boltzmann constant. [1]
  • NA : The Avogadro constant. [2]
  • e : The elementary charge. [3]
  • c0 : The velocity of light in vacuum. [4]
  • mu0: Vacuum permeability. [5]
  • e0: Vacuum permittivity (or permittivity of free space or electric constant). [6]

Unit conversions

For defining the different unit conversions a simple set of base conversions are defined. These represent some common units that are convenient for input and output. For each dimension [12] we define some units and the conversion between these. The base units are:

Charge:
  • e: Electron charge.
  • C: Coulomb.
Energy:
  • kcal: Kilocalorie.
  • kcal/mol: Kilocalorie per mol.
  • J: Joule.
  • J/mol: Joule per mol.
  • kJ/mol: Kilojoule per mol.
  • eV: Electronvolt.
  • hartree: Hartree (atomic unit of energy).
Force:
  • N: Newton.
  • pN: Piconewton.
  • dyn: Dyne.
Length:
  • A: Ångström.
  • nm: Nanometre.
  • bohr: Bohr radius.
  • m: Meter.
Mass:
  • g/mol: Grams per mol, numerically equal to the atomic mass unit.
  • g: Gram.
  • kg: Kilogram.
Pressure:
  • Pa: Pascal.
  • bar: Bar.
  • atm: Atmosphere.
Temperature:
  • K: Kelvin.
Time:
  • s: Second.
  • ps: Picosecond.
  • fs: Femtosecond
  • ns: Nanosecond.
  • us: Microsecond.
  • ms: Millisecond.
Velocity:
  • m/s: Meter per second.
  • nm/ps: Nanometer per picosecond.
  • A/fs: Ångström per femtosecond.
  • A/ps: Ångström per picosecond.

Unit conversions and internal systems of units

The following system of units are defined for PyRETIS:

  • lj: A Lennard-Jones type of units.
  • real: A system of units similar to the LAMMPS unit real. [8]
  • metal: A system of units similar to the LAMMPS unit metal. [8]
  • au: Atomic units. [9]
  • electron: A system of units similar to the LAMMPS unit electron. [8]
  • si: A system of units similar to the LAMMPS unit si. [8]
  • gromacs: A system of units similar to the units used by GROMACS. [10]
  • reduced: A reduced system of units.
  • cp2k: A system of units consistent with CP2K.

The defining units for the Lennard-Jones units (lj) are typically based on the Lennard-Jones parameters for one of the components, e.g. \varepsilon, \sigma and the atomic mass of argon (119.8 kB, 3.405 Å, 39.948 g/mol). [11] The defining units for the other systems are given in the table below:

Table 50 Defining units for energy systems
System name Energy unit Length unit Mass unit
real 1 kcal/mol 1 Å 1 g/mol
metal 1 eV 1 Å 1 g/mol
au 1 hartree 1 bohr 9.10938291e-31 kg
electron 1 hartree 1 bohr 1 g/mol
si 1 J 1 m 1 kg
gromacs 1 kJ/mol 1 nm 1 g/mol
cp2k 1 hartree 1 Å 1 g/mol

These units are also used for the input and define the time unit. Further, all system of units expect an input temperature in Kelvin (K) and all systems, with the exception of si, expect a charge in units of electron charges. The si system uses here Coulomb as its unit of charge. The time units for the different energy systems are given in the table below.

Table 51 Time units for different systems
System name Time unit
real 48.8882129084 fs
metal 10.1805056505 fs
au 0.0241888423521 fs
electron 1.03274987345 fs
si 1.0 s
gromacs 1.0 ps
cp2k 0.0457102889683 fs

The interpretation here is that if you are for instance using the system real and would like to have a time step equal to 0.5 fs, then the input time step to PyRETIS should be 0.5 fs / 48.8882129084 fs. NOTE: When using external engines, PyRETIS will not do any assumptions on the input time/length etc and simply assume that the input values are given in correct units for the external engine. In this case, the only time PyRETIS will make use of the unit system is when the Boltzmann constant is used together with energies. That is, the Boltzmann constant must be in units consistent with the energy output from the external engine.

References and footnotes

[1](1, 2) https://en.wikipedia.org/wiki/Boltzmann_constant
[2]https://en.wikipedia.org/wiki/Avogadro_constant
[3]https://en.wikipedia.org/wiki/Elementary_charge
[4]https://en.wikipedia.org/wiki/Speed_of_light
[5]https://en.wikipedia.org/wiki/Vacuum_permeability
[6]https://en.wikipedia.org/wiki/Vacuum_permittivity
[7]National Institute of Standards and Technology, http://physics.nist.gov/cuu/Constants/Table/allascii.txt
[8](1, 2, 3, 4) The LAMMPS manual, http://lammps.sandia.gov/doc/units.html
[9]https://en.wikipedia.org/wiki/Atomic_units
[10]The GROMACS manual, tables 2.1 and 2.2 on page. 8, http://manual.gromacs.org/documentation/5.1.1/manual-5.1.1.pdf
[11]Rowley et al., J. Comput. Phys., vol. 17, pp. 401-414, 1975, doi: http://dx.doi.org/10.1016/0021-9991
[12]Note that ‘dimension’ here is, strictly speaking, not a true dimension, for instance, we define conversions for the dimension velocity which in reality is composed of the dimensions length and time.

Examples

>>> from pyretis.core.units import CONVERT  
>>> print(CONVERT['length'])
{('A', 'nm'): 0.1, ('A', 'bohr'): 1.8897261254578281, ('A', 'm'): 1e-10}
>>> from pyretis.core.units import create_conversion_factors
>>> create_conversion_factors('lj', length=(3.405, 'A'), energy=(119.8, 'kB'),
...                           mass=(39.948, 'g/mol'), charge='e')
>>> print(CONVERT['length']['bohr', 'nm'])
0.052917721067
>>> print(CONVERT['length']['lj', 'nm'])
0.3405
>>> print(CONVERT['length']['bohr', 'lj'])
0.155411809301
>>> create_conversion_factors('cgs', length=(0.01, 'm'),
...                           energy=(1.0e-7, 'J'),
...                           mass=(1.0, 'g'), charge='e')
>>> print(round(CONVERT['force']['cgs', 'dyn'], 2))
1.0
>>> print(round(CONVERT['time']['cgs', 's'], 2))
1.0
pyretis.core.units._add_conversion_and_inverse(conv_dict, value, unit1, unit2)[source]

Add a specific conversion and it’s inverse.

This function is mainly here to ensure that we don’t forget to add the inverse conversions.

Parameters:
  • conv_dict (dict) – This is where we store the conversion.
  • value (float) – The conversion factor to add
  • unit1 (string) – The unit we are converting from.
  • unit2 (string) – The unit we are converting to.
Returns:

None, but updates the given conv_dict.

pyretis.core.units._check_input_unit(unit, dim, input_unit)[source]

Check input units for create_conversion_factors().

Parameters:
  • unit (string) – Name for the unit system we are dealing with.
  • dim (string) – The dimension we are looking at, typically ‘length’, ‘mass’ or ‘energy’.
  • input_unit (tuple) – This is the input unit on the form (value, string) where the value is the numerical value and the string the unit, e.g. (1.0, nm).
Returns:

out (tuple) – The input_unit if it passes the tests, otherwise an exception will be raised. If the input_unit is None the default values from UNIT_SYSTEMS will be returned if they have been defined.

Raises:

ValueError – If the unit in input_unit is unknown or malformed.

pyretis.core.units._examples()[source]

Just some examples of usage.

pyretis.core.units._generate_conversion_for_dim(conv_dict, dim, unit)[source]

Generate conversion factors for the specified dimension.

It will generate all conversions from a specified unit for a given dimension considering all other units defined in UNITS.

Parameters:
  • conv_dict (dict) – A dictionary with conversions which we wish to update.
  • dim (string) – The dimension to consider
  • unit (string) – The unit we create conversions for.
Returns:

None, but updates the given conv_dict

pyretis.core.units.bfs_convert(conversions, unit_from, unit_to)[source]

Generate unit conversion between the provided units.

The unit conversion can be obtained given that a “path” between these units exist. This path is obtained by a Breadth-first search.

Parameters:
  • conversions (dictionary) – The unit conversion as convert[quantity].
  • unit_from (string) – Starting unit.
  • unit_to (string) – Target unit.
Returns:

  • out[0] (tuple) – A tuple containing the two units: (unit_from, unit_to).
  • out[1] (float) – The conversion factor.
  • out[2] (list of tuples) – The ‘path’ between the units, i.e. the traversal from unit_from to unit_to. out[2][i] gives the (unit_from, unit_to) tuple for step i in the conversion.

pyretis.core.units.convert_bases(dimension)[source]

Create all conversions between base units.

This function will generate all conversions between base units defined in a UNITS[dimension] dictionary. It assumes that one of the bases have been used to defined conversions to all other bases.

Parameters:dimension (string) – The dimension to convert for.
pyretis.core.units.create_conversion_factors(unit, length=None, energy=None, mass=None, charge=None)[source]

Set up conversion factors for a system of units.

Parameters:
  • unit (string) – A name for the system of units
  • length (tuple, optional) – This is the length unit given as (float, string) where the float is the numerical value and the string the unit, e.g. (1.0, nm).
  • energy (tuple, optional) – This is the energy unit given as (float, string) where the float is the numerical value and the string the unit, e.g. (1.0, eV).
  • mass (tuple, optional) – This is the mass unit given as (float, string) where the float is the numerical value and the string the unit, e.g. (1.0, g/mol).
  • charge (string, optional) – This is the unit of charge given as a string, e.g. ‘e’ or ‘C’.
Returns:

  • None but will update CONVERT so that the conversion factors are
  • available.

pyretis.core.units.generate_conversion_factors(unit, distance, energy, mass, charge='e')[source]

Create conversions for a system of units from fundamental units.

This will create a system of units from the three fundamental units distance, energy and mass.

Parameters:
  • unit (string) – This is a label for the unit
  • distance (tuple) – This is the distance unit. The form is assumed to be (value, unit) where the unit is one of the known distance units, ‘nm’, ‘A’, ‘m’.
  • energy (tuple) – This is the energy unit. The form is assumed to be (value, unit) where unit is one of the known energy units, ‘J’, ‘kcal’, ‘kcal/mol’, ‘kb’.
  • mass (tuple) – This is the mass unit. The form is assumed to be (value, unit) where unit is one of the known mass units, ‘g/mol’, ‘kg’, ‘g’.
  • charge (string, optional) – This selects the base charge. It can be ‘C’ or ‘e’ for Coulomb or the electron charge. This will determine how we treat Coulomb’s constant.
pyretis.core.units.generate_inverse(conversions)[source]

Generate all simple inverse conversions.

A simple inverse conversion is something we can obtain by doing a 1 / unit type of conversion.

Parameters:conversions (dictionary) – The unit conversions as convert[quantity]. Note that this variable will be updated by this function.
Returns:out (None) – Will not return anything, but will update the given parameter conversions.
pyretis.core.units.print_table(unit, system=False)[source]

Print out tables with conversion factors.

This is a table in rst format which is useful for displaying conversions.

Parameters:
  • unit (string) – The unit to print out conversions for.
  • system (boolean, optional) – Determines if information for system conversions should be printed.
Returns:

out (None) – Does not return anything, but will print to the screen.

pyretis.core.units.read_conversions(filename='units.txt', select_units=None)[source]

Load conversion factors from a file.

This will load unit conversions from a file.

Parameters:
  • filename (string, optional) – The file to load units from.
  • select_units (string, optional) – If select_units is different from None, it can be used to pick out only conversions for a specific system of units, e.g. ‘real’ or ‘gromacs’, … etc.
Returns:

out (dict) – A dictionary with the conversions.

pyretis.core.units.units_from_settings(settings)[source]

Set up units from given input settings.

Parameters:settings (dict) – A dict defining the units.
Returns:msg (string) – Just a string with some information about the units created. This can be used for printing out some info to the user.
pyretis.core.units.write_conversions(filename='units.txt')[source]

Print out the information currently stored in CONVERT.

This function is intended for creating a big list of conversion factors that can be included in this script for defining unit conversions.

Parameters:filename (string, optional) – The file to write units to.
Returns:out (None) – Will not return anything, but writes the given file.
pyretis.core.units.CONSTANTS = Dictionary with natural constants

A dict containing the natural constants. Natural constants can be accessed with CONSTANTS[‘kB’][‘eV/K’] etc.

pyretis.core.units.CONVERT = Dictionary with conversion factors

A dictionary with conversion factors. It is used to convert between different units, e.g. CONVERT[‘length’][‘bohr’, ‘nm] will convert from the length unit bohr to the length unit nm.

pyretis.core.units.DIMENSIONS = Dictionary with known dimensions

A dictionary with the known dimensions. Note that not all of these are true dimensions, for instance, we are using velocity as a dimension here. This is just because it is convenient to use this to get conversion factors for velocities.

pyretis.core.units.UNITS = Dictionary with known base units

A dictionary of sets. Each set defines the known base unit for a dimension, i.e. UNITS[‘length’] is the set with known base units for the length: UNITS[‘length’] = {‘A’, ‘nm’, ‘bohr’, ‘m’}

pyretis.core.units.UNIT_SYSTEMS = Dictionary with definitions of systems of units

A dictionary containing basic information about the different unit systems. E.g. UNIT_SYSTEMS[‘lj’][‘length’] contains the length unit for the ‘lj’ unit system.