pylj.md

Functions related to the molecular dynamics components of pylj.

pylj.md.calculate_msd(particles, initial_particles, box_length)[source]

Determines the mean squared displacement of the particles in the system. :param particles: Information about the particles. :type particles: util.particle_dt, array_like :param initial_particles: Information about the initial state of the particles. :type initial_particles: util.particle_dt, array_like :param box_length: Size of the cell vector. :type box_length: float

Returns:Mean squared deviation for the particles at the given timestep.
Return type:float
pylj.md.calculate_temperature(particles, mass)[source]

Determine the instantaneous temperature of the system. :param particles: Information about the particles. :type particles: util.particle_dt, array_like

Returns:Calculated instantaneous simulation temperature.
Return type:float
pylj.md.compute_energy(particles, box_length, cut_off, constants, forcefield)[source]

Calculates the total energy of the simulation. :param particles: Information about the particles. :type particles: util.particle_dt, array_like :param box_length: Length of a single dimension of the simulation square, in Angstrom. :type box_length: float :param cut_off: The distance greater than which the energies between particles is

taken as zero.
Parameters:
  • constants (float, array_like (optional)) – The constants associated with the particular forcefield used, e.g. for the function forcefields.lennard_jones, theses are [A, B]
  • forcefield (function (optional)) – The particular forcefield to be used to find the energy and forces.
Returns:

  • util.particle_dt, array_like – Information about particles, with updated accelerations and forces.
  • float, array_like – Current distances between pairs of particles in the simulation.
  • float, array_like – Current energies between pairs of particles in the simulation.

pylj.md.compute_force(particles, box_length, cut_off, constants, forcefield, mass)[source]

Calculates the forces and therefore the accelerations on each of the particles in the simulation. :param particles: Information about the particles. :type particles: util.particle_dt, array_like :param box_length: Length of a single dimension of the simulation square, in Angstrom. :type box_length: float :param cut_off: The distance greater than which the forces between particles is taken

as zero.
Parameters:
  • constants (float, array_like (optional)) – The constants associated with the particular forcefield used, e.g. for the function forcefields.lennard_jones, theses are [A, B]
  • forcefield (function (optional)) – The particular forcefield to be used to find the energy and forces.
  • mass (float (optional)) – The mass of the particle being simulated (units of atomic mass units).
Returns:

  • util.particle_dt, array_like – Information about particles, with updated accelerations and forces.
  • float, array_like – Current distances between pairs of particles in the simulation.
  • float, array_like – Current forces between pairs of particles in the simulation.
  • float, array_like – Current energies between pairs of particles in the simulation.

pylj.md.heat_bath(particles, temperature_sample, bath_temperature)[source]

Rescales the velocities of the particles in the system to control the temperature of the simulation. Thereby allowing for an NVT ensemble. The velocities are rescaled according the following relationship, .. math:

v_{\text{new}} = v_{\text{old}} \times
\sqrt{\frac{T_{\text{desired}}}{\bar{T}}}
Parameters:
  • particles (util.particle_dt, array_like) – Information about the particles.
  • temperature_sample (float, array_like) – The temperature at each timestep in the simulation.
  • bath_temp (float) – The desired temperature of the simulation.
Returns:

Information about the particles with new, rescaled velocities.

Return type:

util.particle_dt, array_like

pylj.md.initialise(number_of_particles, temperature, box_length, init_conf, timestep_length=1e-14, mass=39.948, constants=[1.363e-134, 9.273e-78], forcefield=CPUDispatcher(<function lennard_jones>))[source]

Initialise the particle positions (this can be either as a square or random arrangement), velocities (based on the temperature defined, and calculate the initial forces/accelerations.

Parameters:
  • number_of_particles (int) – Number of particles to simulate.
  • temperature (float) – Initial temperature of the particles, in Kelvin.
  • box_length (float) – Length of a single dimension of the simulation square, in Angstrom.
  • init_conf (string, optional) – The way that the particles are initially positioned. Should be one of: - ‘square’ - ‘random’
  • timestep_length (float (optional)) – Length for each Velocity-Verlet integration step, in seconds.
  • mass (float (optional)) – The mass of the particles being simulated.
  • constants (float, array_like (optional)) – The values of the constants for the forcefield used.
  • forcefield (function (optional)) – The particular forcefield to be used to find the energy and forces.
Returns:

System information.

Return type:

System

pylj.md.initialize(number_particles, temperature, box_length, init_conf, timestep_length=1e-14)[source]

Maps to the md.initialise function to account for US english spelling.

pylj.md.sample(particles, box_length, initial_particles, system)[source]

Sample parameters of interest in the simulation. :param particles: Information about the particles. :type particles: util.particle_dt, array_like :param box_length: Length of a single dimension of the simulation square, in Angstrom. :type box_length: float :param initial_particles: Information about the initial particle conformation. :type initial_particles: util.particle_dt, array-like :param system: Details about the whole system :type system: System

Returns:Details about the whole system, with the new temperature, pressure, msd, and force appended to the appropriate arrays.
Return type:System
pylj.md.update_positions(positions, old_positions, velocities, accelerations, timestep_length, box_length)[source]

Update the particle positions using the Velocity-Verlet integrator. :param positions: Where N is the number of particles, and the first row are the x

positions and the second row the y positions.
Parameters:
  • old_positions ((2, N) array_like) – Where N is the number of particles, and the first row are the previous x positions and the second row are the y positions.
  • velocities ((2, N) array_like) – Where N is the number of particles, and the first row are the x velocities and the second row the y velocities.
  • accelerations ((2, N) array_like) – Where N is the number of particles, and the first row are the x accelerations and the second row the y accelerations.
  • timestep_length (float) – Length for each Velocity-Verlet integration step, in seconds.
  • box_length (float) – Length of a single dimension of the simulation square, in Angstrom.
Returns:

Updated positions.

Return type:

(2, N) array_like

pylj.md.update_velocities(velocities, accelerations_old, accelerations_new, timestep_length)[source]

Update the particle velocities using the Velocity-Verlet algoritm. :param velocities: Where N is the number of particles, and the first row are the x

velocities and the second row the y velocities.
Parameters:
  • accelerations ((2, N) array_like) – Where N is the number of particles, and the first row are the x accelerations and the second row the y accelerations.
  • timestep_length (float) – Length for each Velocity-Verlet integration step, in seconds.
Returns:

Updated velocities.

Return type:

(2, N) array_like

pylj.md.velocity_verlet[source]

Uses the Velocity-Verlet integrator to move forward in time. The Updates the particles positions and velocities in terms of the Velocity Verlet algorithm. Also calculates the instanteous temperature, pressure, and force and appends these to the appropriate system array. :param particles: Information about the particles. :type particles: util.particle_dt, array_like :param timestep_length: Length for each Velocity-Verlet integration step, in seconds. :type timestep_length: float :param box_length: Length of a single dimension of the simulation square, in Angstrom. :type box_length: float

Returns:Information about the particles, with new positions and velocities.
Return type:util.particle_dt, array_like