5. Running the simulation

To run the integrator call the method espressomd.integrate.Integrator.run():

system.integrator.run(number_of_steps, recalc_forces=False, reuse_forces=False)

where number_of_steps is the number of time steps the integrator should perform. The two main integration schemes of ESPResSo are the Velocity Verlet algorithm and an adaption of the Velocity Verlet algorithm to simulate an NPT ensemble. A steepest descent implementation is also available.

5.1. Velocity Verlet algorithm

espressomd.integrate.Integrator.set_vv()

The equations of motion for the trajectory of point-like particles read

\[\begin{split}\dot v_i(t) = F_i(\{x_j\},v_i,t)/m_i \\ \dot x_i(t) = v_i(t),\end{split}\]

where \(x_i\), \(v_i\), \(m_i\) are position, velocity and mass of particle \(i\) and \(F_i(\{x_j\},v_i,t)\) the forces acting on it. These forces comprise all interactions with other particles and external fields as well as non-deterministic contributions described in Thermostats.

For numerical integration, this equation is discretized to the following steps ([Rap04] eqs. 3.5.8 - 3.5.10):

  1. Calculate the velocity at the half step

\[v(t+dt/2) = v(t) + \frac{F(x(t),v(t-dt/2),t)}{m} dt/2\]
  1. Calculate the new position

\[x(t+dt) = x(t) + v(t+dt/2) dt\]
  1. Calculate the force based on the new position

\[F = F(x(t+dt), v(t+dt/2), t+dt)\]
  1. Calculate the new velocity

\[v(t+dt) = v(t+dt/2) + \frac{F(x(t+dt),t+dt)}{m} dt/2\]

Note that this implementation of the Velocity Verlet algorithm reuses forces in step 1. That is, they are computed once in step 3, but used twice, in step 4 and in step 1 of the next iteration. In the first time step after setting up, there are no forces present yet. Therefore, ESPResSo has to compute them before the first time step. That has two consequences: first, random forces are redrawn, resulting in a narrower distribution of the random forces, which we compensate by stretching. Second, coupling forces of e.g. the lattice Boltzmann fluid cannot be computed and are therefore lacking in the first half time step. In order to minimize these effects, ESPResSo has a quite conservative heuristics to decide whether a change makes it necessary to recompute forces before the first time step. Therefore, calling 100 times espressomd.integrate.Integrator.run() with steps=1 does the same as with steps=100, apart from some small calling overhead.

However, for checkpointing, there is no way for ESPResSo to tell that the forces that you read back in actually match the parameters that are set. Therefore, ESPResSo would recompute the forces before the first time step, which makes it essentially impossible to checkpoint LB simulations, where it is vital to keep the coupling forces. To work around this, there is an additional parameter reuse_forces, which tells integrate to not recalculate the forces for the first time step, but use that the values still stored with the particles. Use this only if you are absolutely sure that the forces stored match your current setup!

The opposite problem occurs when timing interactions: In this case, one would like to recompute the forces, despite the fact that they are already correctly calculated. To this aim, the option recalc_forces can be used to enforce force recalculation.

5.2. Isotropic NPT thermostat

set_isotropic_npt()

As the NPT thermostat alters the way the equations of motion are integrated, it is discussed here and only a brief summary is given in Thermostats.

To activate the NPT integrator, use set_isotropic_npt() with parameters:

  • ext_pressure: (float) The external pressure as float variable.

  • piston: (float) The mass of the applied piston as float variable.

  • direction: [int, int, int] Flags to enable/disable box dimensions to be subject to fluctuations. By default, all directions are enabled.

Additionally, a NPT thermostat has to be set by set_npt() with parameters:

  • kT: (float) Thermal energy of the heat bath

  • gamma0: (float) Friction coefficient of the bath

  • gammav: (float) Artificial friction coefficient for the volume fluctuations.

A code snippet would look like:

import espressomd

system = espressomd.System()
system.thermostat.set_npt(kT=1.0, gamma0=1.0, gammav=1.0)
system.integrator.set_isotropic_npt(ext_pressure=1.0, piston=1.0)

The physical meaning of these parameters is described below:

The relaxation towards a desired pressure \(P\) (parameter ext_pressure) is enabled by treating the box volume \(V\) as a degree of freedom with corresponding momentum \(\Pi = Q\dot{V}\), where \(Q\) (parameter piston) is an artificial piston mass. Which box dimensions are affected to change the volume can be controlled by a list of boolean flags for parameter direction. An additional energy \(H_V = 1/(2Q)\Pi + PV\) associated with the volume is postulated. This results in a “force” on the box such that

\[\dot{\Pi} = \mathcal{P} - P\]

where

\[\mathcal{P} = \frac{1}{Vd} \sum_{i,j} f_{ij}x_{ij} + \frac{1}{Vd} \sum_i m_i v_i^2\]

Here \(\mathcal{P}\) is the instantaneous pressure, \(d\) the dimension of the system (number of flags set by direction), \(f_{ij}\) the short range interaction force between particles \(i\) and \(j\) and \(x_{ij}= x_j - x_i\).

In addition to this deterministic force, a friction \(-\frac{\gamma^V}{Q}\Pi(t)\) and noise \(\sqrt{k_B T \gamma^V} \eta(t)\) are added for the box volume dynamics and the particle dynamics. . This introduces three new parameters: The friction coefficient for the box \(\gamma^V\) (parameter gammav), the friction coefficient of the particles \(\gamma^0\) (parameter gamma0) and the thermal energy \(k_BT\) (parameter kT). For a discussion of these terms and their discretisation, see Langevin thermostat, which uses the same approach, but only for particles. As a result of box geometry changes, the particle positions and velocities have to be rescaled during integration.

The discretisation consists of the following steps (see [KDunweg99] for a full derivation of the algorithm):

  1. Calculate the particle velocities at the half step

\[v'(t+dt/2) = v(t) + \frac{F(x(t),v(t-dt/2),t)}{m} dt/2\]
  1. Calculate the instantaneous pressure and “volume momentum”

\[\mathcal{P} = \mathcal{P}(x(t),V(t),f(x(t)), v'(t+dt/2))\]
\[\Pi(t+dt/2) = \Pi(t) + (\mathcal{P}-P) dt/2 -\frac{\gamma^V}{Q}\Pi(t) dt/2 + \sqrt{k_B T \gamma^V dt} \overline{\eta}\]
  1. Calculate box volume and scaling parameter \(L\) at half step and full step, scale the simulation box accordingly

\[V(t+dt/2) = V(t) + \frac{\Pi(t+dt/2)}{Q} dt/2\]
\[L(t+dt/2) = V(t+dt/2)^{1/d}\]
\[V(t+dt) = V(t+dt/2) + \frac{\Pi(t+dt/2)}{Q} dt/2\]
\[L(t+dt) = V(t+dt)^{1/d}\]
  1. Update particle positions and scale velocities

\[x(t+dt) = \frac{L(t+dt)}{L(t)} \left[ x(t) + \frac{L^2(t)}{L^2(t+dt/2)} v(t+dt/2) dt \right]\]
\[v(t+dt/2) = \frac{L(t)}{L(t+dt)} v'(t+dt/2)\]
  1. Calculate forces, instantaneous pressure and “volume momentum”

\[F = F(x(t+dt),v(t+dt/2),t)\]
\[\mathcal{P} = \mathcal{P}(x(t+dt),V(t+dt),f(x(t+dt)), v(t+dt/2))\]
\[\Pi(t+dt) = \Pi(t+dt/2) + (\mathcal{P}-P) dt/2 -\frac{\gamma^V}{Q}\Pi(t+dt/2) dt/2 + \sqrt{k_B T \gamma^V dt} \overline{\eta}\]
  1. Update the velocities

\[v(t+dt) = v(t+dt/2) + \frac{F(t+dt)}{m} dt/2\]

Notes:

  • The NPT algorithm is only tested for all 3 directions enabled for scaling. Usage of direction is considered an experimental feature.

  • In step 4, only those coordinates are scaled for which direction is set.

  • The for the instantaneous pressure the same limitations of applicability hold as described in Pressure.

  • The particle forces \(F\) include interactions as well as a friction and noise term analogous to the terms in the Langevin thermostat.

  • The particle forces are only calculated in step 5 and then reused in step 1 of the next iteration. See Velocity Verlet algorithm for the implications of that.

5.3. Rotational degrees of freedom and particle anisotropy

When the feature ROTATION is compiled in, particles not only have a position, but also an orientation that changes with an angular velocity. A torque on a particle leads to a change in angular velocity depending on the particles rotational inertia. The property espressomd.particle_data.ParticleHandle.rinertia has to be specified as the three eigenvalues of the particles rotational inertia tensor.

The rotational degrees of freedom are also integrated using a velocity Verlet scheme. The implementation is based on a quaternion representation of the particle orientation and described in [Ome98].

When the Langevin thermostat is enabled, the rotational degrees of freedom are also thermalized.

Whether or not rotational degrees of freedom are propagated, is controlled on a per-particle and per-axis level, where the axes are the Cartesian axes of the particle in its body-fixed frame. It is important to note that starting from version 4.0 and unlike in earlier versions of ESPResSo, the particles’ rotation is disabled by default. In this way, just compiling in the ROTATION feature no longer changes the physics of the system.

The rotation of a particle is controlled via the espressomd.particle_data.ParticleHandle.rotation property. E.g., the following code adds a particle with rotation enabled on the x axis:

import espressomd
s = espressomd.System()
s.part.add(pos=(0, 0, 0), rotation=(1, 0, 0))

Notes:

  • The orientation of a particle is stored as a quaternion in the espressomd.particle_data.ParticleHandle.quat property. For a value of (1,0,0,0), the body and space frames coincide.

  • The space-frame direction of the particle’s z-axis in its body frame is accessible through the espressomd.particle_data.ParticleHandle.director property.

  • Any other vector can be converted from body to space fixed frame using the espressomd.particle_data.ParticleHandle.convert_vector_body_to_space method.

  • When DIPOLES are compiled in, the particles dipole moment is always co-aligned with the z-axis in the body-fixed frame.

  • Changing the particles dipole moment and director will re-orient the particle such that its z-axis in space frame is aligned parallel to the given vector. No guarantees are made for the other two axes after setting the director or the dipole moment.

The following particle properties are related to rotation:

5.4. Steepest descent

espressomd.integrate.Integrator.set_steepest_descent()

This feature is used to propagate each particle by a small distance parallel to the force acting on it. When only conservative forces for which a potential exists are in use, this is equivalent to a steepest descent energy minimization. A common application is removing overlap between randomly placed particles.

Please note that the behavior is undefined if a thermostat is activated. It runs a simple steepest descent algorithm:

\[\vec{r}_{i+1} = \vec{r}_i + \min(\gamma \vec{F}_i, \vec{r}_{\text{max_displacement}}),\]

while the maximal force is bigger than f_max or for at most max_steps times. The energy is relaxed by gamma, while the change per coordinate per step is limited to max_displacement. The combination of gamma and max_displacement can be used to get a poor man’s adaptive update. Rotational degrees of freedom are treated similarly: each particle is rotated around an axis parallel to the torque acting on the particle. Please be aware of the fact that this needs not to converge to a local minimum in periodic boundary conditions. Translational and rotational coordinates that are fixed using the fix and rotation attribute of particles are not altered.

Usage example:

system.integrator.set_steepest_descent(
    f_max=0, gamma=0.1, max_displacement=0.1)
system.integrator.run(20)
system.integrator.set_vv()  # to switch back to velocity verlet