pyretis.simulation package

This package defines different simulations for use with PyRETIS.

The different simulations are defined as objects which inherit from the base Simulation object defined in simulation.py. The simulation object defines as simulation as a series of tasks to be executed, typically at each step of the simulation. These tasks may produce results which can be outputted to the user in some way.

Package structure

Modules

md_simulation.py (pyretis.simulation.md_simulation)
Defines simulation classes for molecular dynamics simulations.
mc_simulation.py (pyretis.simulation.mc_simulation)
Define simulation classes for Monte Carlo simulations.
path_simulation.py (pyretis.simulation.path_simulation)
Defines simulation classes for path simulations.
simulation.py (pyretis.simulation.simulation)
Defines the Simulation class which is the base class for simulations.
simulation_task.py (pyretis.simulation.simulation_task)
Defines classes for the handling of simulation tasks.

Important classes defined in this package

Simulation (Simulation)
The base class for simulations.
SimulationTask (SimulationTask)
A class for creating tasks for simulations.
SimulationSingleTIS (SimulationSingleTIS)
A class for running a TIS simulation for a single ensemble.
SimulationRETIS (SimulationRETIS)
A class for running a RETIS simulation for a set of ensembles.

pyretis.simulation.mc_simulation module

Definition of simulation objects for Monte Carlo simulations.

This module defines some classes and functions for performing Monte Carlo simulations.

Important classes defined here

UmbrellaWindowSimulation (UmbrellaWindowSimulation)
Defines a simulation for performing umbrella window simulations. Several umbrella window simulations can be joined to perform a umbrella simulation.
class pyretis.simulation.mc_simulation.UmbrellaWindowSimulation(system, umbrella, overlap, maxdx, rgen=None, mincycle=0, startcycle=0)[source]

Bases: pyretis.simulation.simulation.Simulation

This class defines an Umbrella simulation.

The Umbrella simulation is a special case of the simulation class with settings to simplify the execution of the umbrella simulation.

system

The system to act on.

Type:object like System
umbrella

The umbrella window.

Type:list = [float, float]
overlap

The positions that must be crossed before the simulation is done.

Type:float
startcycle

The current simulation cycle.

Type:int
mincycle

The MINIMUM number of cycles to perform.

Type:int
rgen

Object to use for random number generation.

Type:object like RandomGenerator
maxdx

Maximum allowed displacement in the Monte Carlo step.

Type:float
__init__(system, umbrella, overlap, maxdx, rgen=None, mincycle=0, startcycle=0)[source]

Initialise the umbrella simulation simulation.

Parameters:
  • system (object like System) – The system to act on.
  • umbrella (list, [float, float]) – The umbrella window to consider.
  • overlap (float) – The position we have to cross before the simulation is done.
  • maxdx (float) – Defines the maximum movement allowed in the Monte Carlo moves.
  • rgen (object like RandomGenerator) – Object to use for random number generation.
  • mincycle (int, optional) – The MINIMUM number of cycles to perform. Note that in the base Simulation class this is the MAXIMUM number of cycles to perform. The meaning is redefined in this class by overriding self.simulation_finished.
  • startcycle (int, optional) – The current simulation cycle, i.e. where we start.
__str__()[source]

Return some info about the simulation as a string.

is_finished()[source]

Check if the simulation is done.

In the umbrella simulation, the simulation is finished when we cycle is larger than maxcycle and all particles have crossed self.overlap.

Returns:out (boolean) – True if the simulation is finished, False otherwise.
simulation_type = 'umbrella-window'

pyretis.simulation.md_simulation module

Definitions of simulation objects for molecular dynamics simulations.

This module contains definitions of classes for performing molecular dynamics simulations.

Important classes defined here

SimulationNVE (SimulationNVE)
Definition of a simple NVE simulation. The engine used for this simulation must have dynamics equal to NVE.
SimulationMD (SimulationMD)
Definition of a simulation for running somply MD.
SimulationMDFlux (SimulationMDFlux)
Definition of a simulation for determining the initial flux. This is used for calculating rates in TIS simulations.
class pyretis.simulation.md_simulation.SimulationMD(system, engine, order_function=None, steps=0, startcycle=0)[source]

Bases: pyretis.simulation.simulation.Simulation

A generic MD simulation.

This class is used to define a simple MD simulation.

system

This is the system the simulation will act on.

Type:object like System
engine

The engine to use for integrating the equations of motion.

Type:object like EngineBase
order_function

A class that can be used to calculate an order parameter, if needed.

Type:object like OrderParameter
__init__(system, engine, order_function=None, steps=0, startcycle=0)[source]

Initialise the MD simulation.

Here we just add variables and do not do any other setup.

Parameters:
  • system (object like System) – This is the system we are investigating.
  • engine (object like EngineBase) – This is the integrator that is used to propagate the system in time.
  • order_function (object like OrderParameter, optional) – A class that can be used to calculate an order parameter, if needed.
  • steps (int, optional) – The number of simulation steps to perform.
  • startcycle (int, optional) – The cycle we start the simulation on.
__str__()[source]

Return a string with info about the simulation.

run()[source]

Run the MD simulation.

Yields:results (dict) – The results from a single step in the simulation.
simulation_output = [{'type': 'energy', 'name': 'md-energy-file'}, {'type': 'thermo-file', 'name': 'md-thermo-file'}, {'type': 'traj-xyz', 'name': 'md-traj-file'}, {'type': 'thermo-screen', 'name': 'md-thermo-screen'}, {'type': 'order', 'name': 'md-order-file'}]
simulation_type = 'md'
class pyretis.simulation.md_simulation.SimulationNVE(system, engine, order_function=None, steps=0, startcycle=0)[source]

Bases: pyretis.simulation.md_simulation.SimulationMD

A MD NVE simulation class.

This class is used to define a NVE simulation. Compared with the SimulationMD we here require that the engine supports NVE dynamics.

__init__(system, engine, order_function=None, steps=0, startcycle=0)[source]

Initialise the NVE simulation object.

Here we will set up the tasks that are to be performed in the simulation, such as the integration and thermodynamics calculation(s).

Parameters:
  • system (object like System) – This is the system we are investigating.
  • engine (object like EngineBase) – This is the integrator that is used to propagate the system in time.
  • order_function (object like OrderParameter, optional) – A class that can be used to calculate an order parameter, if needed.
  • steps (int, optional) – The number of simulation steps to perform.
  • startcycle (int, optional) – The cycle we start the simulation on.
__str__()[source]

Return a string with info about the simulation.

simulation_output = [{'type': 'energy', 'name': 'nve-energy-file'}, {'type': 'thermo-file', 'name': 'nve-thermo-file'}, {'type': 'traj-xyz', 'name': 'nve-traj-file'}, {'type': 'thermo-screen', 'name': 'nve-thermo-screen'}, {'type': 'order', 'name': 'nve-order-file'}]
simulation_type = 'md-nve'
step()[source]

Run a single simulation step.

class pyretis.simulation.md_simulation.SimulationMDFlux(system, order_function, engine, interfaces, steps=0, startcycle=0)[source]

Bases: pyretis.simulation.md_simulation.SimulationMD

A simulation for obtaining the initial flux for TIS.

This class is used to define a MD simulation where the goal is to calculate crossings in order to obtain the initial flux for a TIS calculation.

interfaces

These floats define the interfaces used in the crossing calculation.

Type:list of floats
leftside_prev

These are used to store the previous positions with respect to the interfaces.

Type:list of booleans or None
__init__(system, order_function, engine, interfaces, steps=0, startcycle=0)[source]

Initialise the MD-Flux simulation object.

Parameters:
  • system (object like System) – This is the system we are investigating
  • order_function (object like OrderParameter) – The class used for calculating the order parameters.
  • engine (object like EngineBase) – This is the integrator that is used to propagate the system in time.
  • interfaces (list of floats) – These define the interfaces for which we will check the crossing(s).
  • steps (int, optional) – The number of steps to perform.
  • startcycle (int, optional) – The cycle we start the simulation on.
__str__()[source]

Return a string with info about the simulation.

load_restart_info(info)[source]

Load the restart information.

restart_info()[source]

Return restart info.

Here we report the cycle number and the random number generator status.

run()[source]

Run the MD simulation.

Yields:results (dict) – The results from a single step in the simulation.
simulation_output = [{'type': 'energy', 'name': 'flux-energy-file'}, {'type': 'traj-xyz', 'name': 'flux-traj-file'}, {'type': 'thermo-screen', 'name': 'flux-thermo-screen'}, {'type': 'order', 'name': 'flux-order-file'}, {'type': 'cross', 'name': 'flux-cross-file'}]
simulation_type = 'md-flux'

pyretis.simulation.path_simulation module

Definitions of simulation objects for path sampling simulations.

This module defines simulations for performing path sampling simulations.

Important classes defined here

PathSimulation (PathSimulation)
The base class for path simulations.
SimulationSingleTIS (SimulationSingleTIS)
Definition of a TIS simulation for a single path ensemble.
SimulationRETIS (SimulationRETIS)
Definition of a RETIS simulation.
class pyretis.simulation.path_simulation.PathSimulation(system, order_function, engine, path_ensembles, settings, rgen=None, steps=0, startcycle=0)[source]

Bases: pyretis.simulation.simulation.Simulation

A base class for TIS/RETIS simulations.

engine

This is the integrator that is used to propagate the system in time.

Type:object like EngineBase
path_ensembles

This is used for storing results for the different path ensembles.

Type:list of objects like PathEnsemble
rgen

This is a random generator used for the generation of paths.

Type:object like RandomGenerator
system

This is the system the simulation will act on.

Type:object like System
settings

A dictionary with TIS and RETIS settings. We expect that we can find settings['tis'] and possibly settings['retis']. For settings['tis'] we further expect to find the keys:

  • aimless: Determines if we should do aimless shooting (True) or not (False).
  • sigma_v: Scale used for non-aimless shooting.
  • seed: A integer seed for the random generator used for the path ensemble we are simulating here.

Note that the pyretis.core.tis.make_tis_step_ensemble() method will make use of additional keys from settings['tis']. Please see this method for further details. For the settings['retis'] we expect to find the following keys:

  • swapfreq: The frequency for swapping moves.
  • relative_shoots: If we should do relative shooting for the ensembles.
  • nullmoves: Should we perform null moves.
  • swapsimul: Should we just swap a single pair or several pairs.
Type:dict
required_settings

This is just a list of the settings that the simulation requires. Here it is used as a check to see that we have all we need to set up the simulation.

Type:tuple of strings
__init__(system, order_function, engine, path_ensembles, settings, rgen=None, steps=0, startcycle=0)[source]

Initialise the path simulation object.

Parameters:
  • system (object like System) – This is the system we are investigating.
  • order_function (object like OrderParameter) – The object used for calculating the order parameter.
  • engine (object like EngineBase) – This is the integrator that is used to propagate the system in time.
  • path_ensembles (list of objects like PathEnsemble) – This is used for storing results for the different path ensembles.
  • rgen (object like RandomGenerator) – This object is the random generator to use in the simulation.
  • settings (dict of dicts) – A dictionary with TIS and RETIS settings.
  • steps (int, optional) – The number of simulation steps to perform.
  • startcycle (int, optional) – The cycle we start the simulation on.
create_output_tasks(settings, progress=False)[source]

Create output tasks for the simulation.

This method will generate output tasks based on the tasks listed in simulation_output.

Parameters:
  • settings (dict) – These are the simulation settings.
  • progress (boolean) – For some simulations, the user may select to display a progress bar, we then need to disable the screen output.
initiate(settings)[source]

Initialise the path simulation.

Parameters:settings (dictionary) – The simulation settings.
name = 'Generic path simulation'
required_settings = ('tis',)
restart_info()[source]

Return restart info.

The restart info for the path simulation includes the state of the random number generator(s).

simulation_output = [{'type': 'pathensemble', 'name': 'path-ensemble', 'result': ('pathensemble-{}',)}, {'type': 'path-order', 'name': 'path-ensemble-order', 'result': ('path-{}', 'status-{}')}, {'type': 'path-traj-{}', 'name': 'path-ensemble-traj', 'result': ('path-{}', 'status-{}', 'pathensemble-{}')}, {'type': 'path-energy', 'name': 'path-ensemble-energy', 'result': ('path-{}', 'status-{}')}]
simulation_type = 'generic-path'
write_restart(now=False)[source]

Create a restart file.

Parameters:now (boolean, optional) – If True, the output file will be written irrespective of the step number.
class pyretis.simulation.path_simulation.SimulationSingleTIS(system, order_function, engine, path_ensemble, settings, rgen=None, steps=0, startcycle=0)[source]

Bases: pyretis.simulation.path_simulation.PathSimulation

A single TIS simulation.

This class is used to define a TIS simulation where the goal is to calculate crossing probabilities for a single path ensemble.

path_ensemble

This is used for storing results for the simulation. Note that we also have the path_ensembles attribute defined by the parent class. The attribute path_ensemble used here, is meant to reflect that this class is intended for simulating a single ensemble only.

Type:object like PathEnsemble
__init__(system, order_function, engine, path_ensemble, settings, rgen=None, steps=0, startcycle=0)[source]

Initialise the TIS simulation object.

Parameters:
  • system (object like System) – This is the system we are investigating.
  • order_function (object like OrderParameter) – The object used for calculating the order parameter.
  • engine (object like EngineBase) – This is the integrator that is used to propagate the system in time.
  • 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.
  • rgen (object like RandomGenerator) – This is the random generator to use in the simulation.
  • settings (dict) – This dict contains settings for the simulation.
  • steps (int, optional) – The number of simulation steps to perform.
  • startcycle (int, optional) – The cycle we start the simulation on.
__str__()[source]

Just a small function to return some info about the simulation.

name = 'Single TIS simulation'
required_settings = ('tis',)
simulation_output = [{'type': 'pathensemble', 'name': 'path-ensemble', 'result': ('pathensemble-{}',)}, {'type': 'path-order', 'name': 'path-ensemble-order', 'result': ('path-{}', 'status-{}')}, {'type': 'path-traj-{}', 'name': 'path-ensemble-traj', 'result': ('path-{}', 'status-{}', 'pathensemble-{}')}, {'type': 'path-energy', 'name': 'path-ensemble-energy', 'result': ('path-{}', 'status-{}')}, {'type': 'pathensemble-screen', 'name': 'path-ensemble-screen', 'result': ('pathensemble-{}',)}]
simulation_type = 'tis'
step()[source]

Perform a TIS simulation step.

Returns:out (dict) – This list contains the results of the TIS step.
class pyretis.simulation.path_simulation.SimulationRETIS(system, order_function, engine, path_ensembles, settings, rgen=None, steps=0, startcycle=0)[source]

Bases: pyretis.simulation.path_simulation.PathSimulation

A RETIS simulation.

This class is used to define a RETIS simulation where the goal is to calculate crossing probabilities for several path ensembles.

The attributes are documented in the parent class, please see: PathSimulation.

__init__(system, order_function, engine, path_ensembles, settings, rgen=None, steps=0, startcycle=0)[source]

Initialise the RETIS simulation object.

Parameters:
  • system (object like System) – This is the system we are investigating.
  • order_function (object like OrderParameter) – The object used for calculating the order parameter.
  • engine (object like EngineBase) – This is the integrator that is used to propagate the system in time.
  • path_ensembles (list of objects like PathEnsemble) – This is used for storing results for the different path ensembles.
  • rgen (object like RandomGenerator) – This object is the random generator to use in the simulation.
  • settings (dict) – A dictionary with settings for TIS and RETIS.
  • steps (int, optional) – The number of simulation steps to perform.
  • startcycle (int, optional) – The cycle we start the simulation on.
__str__()[source]

Just a small function to return some info about the simulation.

name = 'RETIS simulation'
required_settings = ('tis', 'retis')
simulation_output = [{'type': 'pathensemble', 'name': 'path-ensemble', 'result': ('pathensemble-{}',)}, {'type': 'path-order', 'name': 'path-ensemble-order', 'result': ('path-{}', 'status-{}')}, {'type': 'path-traj-{}', 'name': 'path-ensemble-traj', 'result': ('path-{}', 'status-{}', 'pathensemble-{}')}, {'type': 'path-energy', 'name': 'path-ensemble-energy', 'result': ('path-{}', 'status-{}')}, {'type': 'pathensemble-retis-screen', 'name': 'path-ensemble-retis-screen', 'result': ('pathensemble-{}',)}]
simulation_type = 'retis'
step()[source]

Perform a RETIS simulation step.

Returns:out (dict) – This list contains the results of the defined tasks.

pyretis.simulation.simulation module

Definitions of generic simulation objects.

This module defines the generic simulation object. This is the base class for all other simulations.

Important classes defined here

Simulation (Simulation)
A class defining a generic simulation.
class pyretis.simulation.simulation.Simulation(steps=0, startcycle=0)[source]

Bases: object

This class defines a generic simulation.

cycle

This dictionary stores information about the number of cycles. The keywords are:

  • end: Represents the cycle number where the simulation should end.
  • step: The current cycle number.
  • start: The cycle number we started at.
  • stepno: The number of cycles we have performed to arrive at cycle number given by cycle[‘step’]. Note that cycle[‘stepno’] might be different from cycle[‘step’] since cycle[‘start’] might be != 0.
Type:dict of integers
exe_dir

The path we are running the simulation from.

Type:string
restart_freq

The frequency for creating restart files.

Type:integer
first_step

True if the first step has not been executed yet.

Type:boolean
system

This is the system the simulation will act on.

Type:object like System
simulation_output

This list defines the output tasks associated with the simulation.

Type:list of dicts
simulation_type

An identifier for the simulation.

Type:string
tasks

This is the list of simulation tasks to execute.

Type:list of objects like SimulationTask
__init__(steps=0, startcycle=0)[source]

Initialise the simulation object.

Parameters:
  • steps (int, optional) – The number of simulation steps to perform.
  • startcycle (int, optional) – The cycle we start the simulation on.
__str__()[source]

Just a small function to return some info about the simulation.

add_task(task, position=None)[source]

Add a new simulation task.

A task can still be added manually by simply appending to py:attr:.tasks. This function will, however, do some checks so that the task added can be executed.

Parameters:
  • task (dict) – A dict defining the task. A task is represented by an object of type SimulationTask with some additional settings on how to store the output and when to execute the task. The keywords in the dict defining the task are:
    • func: Callable, this is a function to execute in the task.
    • args: List, with arguments for the function.
    • kwargs: Dict, with the keyword arguments for the function.
    • when: Dict, which defines when the task should be executed.
    • first: Boolean, determines if the task should be executed on the initial step, i.e. before the full simulation starts.
    • result: String, used to label the result.
  • position (int, optional) – Can be used to give the tasks a specific position in the task list.
create_output_tasks(settings, progress=False)[source]

Create output tasks for the simulation.

This method will generate output tasks based on the tasks listed in simulation_output.

Parameters:
  • settings (dict) – These are the simulation settings.
  • progress (boolean) – For some simulations, the user may select to display a progress bar, we then need to disable the screen output.
execute_tasks()[source]

Execute all the tasks in sequential order.

Returns:results (dict) – The results from the different tasks (if any).
extend_cycles(steps)[source]

Extend a simulation with the given number of steps.

Parameters:steps (int) – The number of steps to extend the simulation with.
Returns:out (None) – Returns None but modifies self.cycle.
is_finished()[source]

Determine if the simulation is finished.

In this object, the simulation is done if the current step number is larger than the end cycle. Note that the number of steps performed is dependent on the value of self.cycle[‘start’].

Returns:out (boolean) – True if the simulation is finished, False otherwise.
load_restart_info(info)[source]

Load restart information.

Note, we do not change the end property here as we probably are extending a simulation.

Parameters:info (dict) – The dictionary with the restart information.
restart_info()[source]

Return information which can be used to restart the simulation.

run()[source]

Run a simulation.

The intended usage is for simulations where all tasks have been defined in self.tasks.

Note

This function will just run the tasks via executing step() In general, this is probably too generic for the simulation you want, if you are creating a custom simulation. Please consider customizing the run() (or the step()) method of your simulation class.

Yields:out (dict) – This dictionary contains the results from the simulation.
set_up_output(settings, progress=False)[source]

Set up output from the simulation.

This includes the predefined output tasks, but also output related to the restart file(s).

Parameters:
  • settings (dict) – These are the simulation settings.
  • progress (boolean) – For some simulations, the user may select to display a progress bar, we then need to disable the screen output.
simulation_output = []
simulation_type = 'generic'
soft_exit()[source]

Force simulation to stop at the current step.

step()[source]

Execute a simulation step.

Here, the tasks in tasks will be executed sequentially.

Returns:out (dict) – This dictionary contains the results of the defined tasks.

Note

This function will have ‘side effects’ and update/change the state of other attached variables such as the system or other variables that are not explicitly shown. This is intended and the behavior is defined by the tasks in tasks.

write_restart(now=False)[source]

Create a restart file.

Parameters:now (boolean, optional) – If True, the output file will be written irrespective of the step number.

pyretis.simulation.simulation_task module

Definition of a class for simulation tasks.

Important classes defined here

SimulationTask (SimulationTask)
A class representing a simulation task.
SimulationTaskList (SimulationTaskList)
A class for representing a list of simulation tasks. This class defines functionality for adding tasks from a dictionary description.
class pyretis.simulation.simulation_task.SimulationTask(function, args=None, kwargs=None, when=None, result=None, first=False)[source]

Bases: pyretis.inout.simulationio.Task

Representation of simulation tasks.

This class defines a task object. A task is executed at specific points, at regular intervals etc. in a simulation. A task will typically provide a result, but it does not need to. It can simply just alter the state of the passed argument(s).

function

The function to execute.

Type:function
when

Determines when the task should be executed.

Type:dict
args

List of arguments to the function.

Type:list
kwargs

The keyword arguments to the function.

Type:dict
first

True if this task should be executed before the first step of the simulation.

Type:boolean
result

This is a label for the result created by the task.

Type:string
__call__(step)[source]

Execute the task.

Parameters:

step (dict of ints) – The keys are:

  • ‘step’: the current cycle number.
  • ‘start’: the cycle number at the start.
  • ‘stepno’: the number of cycles we have performed so far.
Returns:

out (unknown type) – The result of self.execute(step).

__init__(function, args=None, kwargs=None, when=None, result=None, first=False)[source]

Initialise the task.

Parameters:
  • function (callable) – The function to execute.
  • args (list, optional) – List of arguments to the function.
  • kwargs (dict, optional) – The keyword arguments to the function.
  • when (dict, optional) – Determines if the task should be executed.
  • result (string, optional) – This is a label for the result created by the task.
  • first (boolean, optional) – True if this task should be executed before the first step of the simulation.
__str__()[source]

Output info about the task.

execute(step)[source]

Execute the task.

Parameters:

step (dict of ints) – The keys are:

  • ‘step’: the current cycle number.
  • ‘start’: the cycle number at the start.
  • ‘stepno’: the number of cycles we have performed so far.
Returns:

out (unknown type) – The result of running self.function.

property result

Return the result label.

run_first()[source]

Return True if task should be executed before first step.

task_dict()[source]

Return a dict representing the task.

pyretis.simulation.simulation_task._check_args(function, given_args=None, given_kwargs=None)[source]

Check consistency for function and the given (keyword) arguments.

Here we assume that the arguments are given in a list and that the keyword arguments are given as a dictionary. The function inspect.getargspec is used to check the input function.

Parameters:
  • function (callable) – The function we will inspect.
  • given_args (list, optional) – A list of the arguments to pass to the function. ‘self’ will not be considered here since it passed implicitly.
  • given_kwargs (dict, optional) – A dictionary with keyword arguments.
Returns:

out (boolean) – False if there is some inconsistencies, i.e. when the calling of the given function will probably fail. True otherwise.