azplugins.pair

Overview

Colloid

Colloid pair potential.

DPDGeneralWeight

Dissipative Particle Dynamics with generalized weight function.

Hertz

Hertz potential.

PerturbedLennardJones

Perturbed Lennard-Jones potential.

TwoPatchMorse

Two-patch Morse potential.

Details

Pair potentials.

class hoomd.azplugins.pair.Colloid(nlist, default_r_cut=None, default_r_on=0, mode='none')

Bases: Pair

Colloid pair potential.

Parameters:

Colloid is the effective Lennard-Jones potential obtained by integrating the Lennard-Jones potential over zero, one, or two spheres. A discussion of the application of these potentials to colloidal suspensions can be found in Grest et al.

The pair potential has three different coupling styles between particles:

  • When the radii of both particles are zero, the potential is the usual Lennard-Jones coupling:

    \[U(r) = \frac{A}{36} \left[\left(\frac{\sigma}{r}\right)^{12} - \left(\frac{\sigma}{r}\right)^6 \right]\]

    The Hamaker constant is related to the Lennard-Jones parameters as \(A = 144 \varepsilon\) (to give the usual prefactor \(4\varepsilon\)).

  • When the radius of one particle is zero and the other radius \(a\) is nonzero, the potential is the interaction between a pointlike particle and a sphere comprised of Lennard-Jones particles:

    \[U(r) = \frac{2 a^3 \sigma^3 A}{9(a^2-r^2)^3} \left[ 1 - \frac{(5 a^6 + 45 a^4 r^2 + 63 a^2 r^4 + 15 r^6) \sigma^6} {15 (a-r)^6 (a+r)^6} \right]\]

    The Hamaker constant is related to the Lennard-Jones parameters as \(A = 24 \pi \varepsilon \sigma^3 \rho\), where \(\rho\) is the number density of particles in the sphere.

  • When the radii of both particles are non-zero, the potential is the interaction between two spheres comprised of Lennard-Jones particles, given by Eqs. (16) and (17) of Everaers and Ejtehadi. The Hamaker constant is related to the Lennard-Jones parameters as \(A = 4 \pi^2 \varepsilon \sigma^6 \rho_1 \rho_2\), where \(\rho_1\) and \(\rho_2\) are the number densities of particles in each sphere.

Example:

# Explicit Solvent Model from https://doi.org/10.1063/1.5043401
nl = nlist.Cell()
colloid = pair.Colloid(default_r_cut=3.0, nlist=nl)
# standard Lennard-Jones for solvent-solvent
colloid.params[('S', 'S')] = dict(A=144.0, a_1=0, a_2=0 sigma=1.0)
# solvent-colloid
colloid.params[('S', 'C')] = dict(A=144.0, a_1=0, a_2=5.0 sigma=1.0)
colloid.r_cut[('S', 'C')] = 9.0
# colloid-colloid
colloid.params[('C','C')] = dict(A=40.0, a_1=5.0, a_2=5.0 sigma=1.0)
colloid.r_cut[('C', 'C')] = 10.581
params

The Colloid potential parameters. The dictionary has the following keys:

  • A (float, required) - Hamaker constant \(A\) \([\mathrm{energy}]\)

  • a_1 (float, required) - Radius of first particle \(a_1\) \([\mathrm{length}]\)

  • a_2 (float, required) - Radius of second particle \(a_2\) \([\mathrm{length}]\)

  • sigma (float, required) - Lennard-Jones particle size \(\sigma\) \([\mathrm{length}]\)

Type: TypeParameter [tuple [particle_type, particle_type], dict]

mode

Energy shifting/smoothing mode: "none", "shift", or "xplor".

Type: str

class hoomd.azplugins.pair.DPDGeneralWeight(nlist, kT, default_r_cut=None)

Bases: Pair

Dissipative Particle Dynamics with generalized weight function.

Parameters:

DPDGeneralWeight is a complete set of DPD forces (conservative force, dissipative force, and random force):

\[\mathbf{F} = \mathbf{F}_{\rm C} + \mathbf{F}_{\rm D} + \mathbf{F}_{\rm R}\]

The conservative force \(\mathbf{F}_{\rm C}\) is the standard form:

\[\mathbf{F}_{\rm C} = A \left(1 - \frac{r_{ij}}{r_{\rm cut}} \right) \mathbf{\hat{r}}_{ij}\]

where A is the interaction parameter and \(r_{\rm cut}\) is the cutoff radius. Here, \(r_{ij} = |\mathbf{r}_{ij}|\) and \(\mathbf{\hat{r}}_{ij} = \mathbf{r}_{ij}/r_{ij}\), where \(\mathbf{r}_{ij} = \mathbf{r}_i - \mathbf{r}_j\) is the vector from particle j to i and \(\mathbf{r}_i\) is the position of particle i. See Groot and Warren for more details.

The dissipative force \(\mathbf{F}_{\rm D}\) is:

\[\mathbf{F}_{\rm D} = -\gamma w(r_{ij}) (\mathbf{v}_{ij} \cdot \mathbf{\hat r}_{ij}) \mathbf{\hat r}_{ij}\]

and the random force \(\mathbf{F}_{\rm R}\) is:

\[\mathbf{F}_{\rm R} = \sqrt{\gamma w(r_{ij})} \theta_{ij} \mathbf{\hat r}_{ij}\]

where \(\mathbf{v}_{ij} = \mathbf{v}_i - \mathbf{v}_j\) with \(\mathbf{v}_i\) being the velocity of particle i, \(w\) is the dissipative weight function, and \(\theta_{ij}\) is a random noise term with zero mean and variance \(\langle \theta_{ij}(t) \theta_{ij}(t') \rangle = 2 k_{\rm B} T \delta(t-t')\).

Unlike hoomd.md.pair.DPD, DPDGeneralWeight uses the general weight function proposed by Fan et al.:

\[w(r) = \left(1 - \frac{r_{ij}}{r_{\rm cut}}\right)^s\]

Setting \(s = 2\) gives the standard DPD weight function.

Example:

nl = nlist.Cell()
dpd = pair.DPDGeneralWeight(nlist=nl, kT=1.0, default_r_cut=1.0)
dpd.params[('A', 'A')] = dict(A=25.0, gamma=4.5, s=0.5)
dpd.params[('A', 'B')] = dict(A=50.0, gamma=4.5, s=0.5)
dpd.params[('B', 'B')] = dict(A=25.0, gamma=4.5, s=0.5)
params

The DPDGeneralWeight potential parameters. The dictionary has the following keys:

  • A (float, required) - Repulsive parameter \(A\) \([\mathrm{force}]\)

  • gamma (float, required) - Friction parameter \(\gamma\) \([\mathrm{mass} \cdot \mathrm{time}^{-1}]\)

  • s (float, required) - Weight function exponent \(s\)

Type: TypeParameter [tuple [particle_type, particle_type], dict]

mode

Energy shifting/smoothing mode: "none".

Type: str

class hoomd.azplugins.pair.Hertz(nlist, default_r_cut=None, default_r_on=0, mode='none')

Bases: Pair

Hertz potential.

Parameters:

The Hertz potential is:

\[U(r) = \varepsilon \left( 1-\frac{ r }{ r_{\rm{cut}} } \right)^{5/2} , \quad r < r_{\rm{cut}}\]

Example:

nl = hoomd.md.nlist.cell()
hertz = azplugins.pair.Hertz(r_cut=3.0, nlist=nl)
hertz.params[('A', 'A')] = dict(epsilon=1.0)
hertz.r_cut[('A', 'B')] = 3.0
params

The Hertz potential parameters. The dictonary has the following key:

  • epsilon (float, required) - energy parameter \(\varepsilon\) \([\mathrm{energy}]\)

Type: TypeParameter [tuple [particle_type, particle_type], dict]

mode

Energy shifting/smoothing mode: "none", "shift", or "xplor".

Type: str

class hoomd.azplugins.pair.PerturbedLennardJones(nlist, default_r_cut=None, default_r_on=0, mode='none')

Bases: Pair

Perturbed Lennard-Jones potential.

Parameters:
  • nlist (hoomd.md.nlist.NeighborList) – Neighbor list.

  • default_r_cut (float) – Default cutoff radius \([\mathrm{length}]\).

  • default_r_on (float) – Default turn-on radius \([\mathrm{length}]\).

  • mode (str) – Energy shifting/smoothing mode.

PerturbedLennardJones is a Lennard-Jones perturbation potential. The potential has a purely repulsive core, and the parameter attraction_scale_factor (lambda) sets the strength of the attractive tail:

\begin{eqnarray*} U(r) &= U_{\mathrm{LJ}}(r) + (1-\lambda)\varepsilon, & r \le 2^{1/6}\sigma \\ &= \lambda U_{\mathrm{LJ}}(r), & r > 2^{1/6}\sigma \end{eqnarray*}

where \(U_{\rm LJ}\) is the standard Lennard-Jones potential (see hoomd.md.pair.LJ). When \(\lambda = 0\), \(U\) is the standard Weeks-Chandler-Anderson repulsive potential, while when \(\lambda = 1\), \(U\) is \(U_{\rm LJ}\).

Example:

nl = nlist.Cell()
perturbed_lj = azplugins.pair.PerturbedLennardJones(
    default_r_cut=3.0, nlist=nl)
perturbed_lj.params[('A', 'A')] = dict(
    epsilon=1.0, sigma=1.0, attraction_scale_factor=0.5)
perturbed_lj.r_cut[('A', 'B')] = 3.0
params

The Perturbed LJ potential parameters. The dictionary has the following keys:

  • epsilon (float, required) - energy parameter \(\varepsilon\) \([\mathrm{energy}]\)

  • sigma (float, required) - particle size \(\sigma\) \([\mathrm{length}]\)

  • attraction_scale_factor (float, required) - scale factor for attraction \(\lambda\), between 0 and 1

Type: TypeParameter [tuple [particle_type, particle_type], dict]

mode

Energy shifting/smoothing mode: "none", "shift", or "xplor".

Type: str

class hoomd.azplugins.pair.TwoPatchMorse(nlist, default_r_cut=None, mode='none')

Bases: AnisotropicPair

Two-patch Morse potential.

Parameters:

TwoPatchMorse is a Morse potential that is modulated by an orientation-dependent function.

\[U(\vec{r}_{ij}, \hat{n}_i, \hat{n}_j) = U_{\rm M}(|\vec{r}_{ij}|) \Omega(\hat{r}_{ij} \cdot \hat{n}_i) \Omega(\hat{r}_{ij} \cdot \hat{n}_j)\]

where \(U_{\rm M}\) is the potential that depends on distance

\[U_{\rm M}(r) = M_{\rm d} \left( \left[ 1 - \exp\left( -\frac{r-r_{\rm eq}}{M_r}\right) \right]^2 - 1 \right)\]

and \(\Omega(\gamma)\) depends on the orientations

\[\Omega(\gamma) = \frac{1}{1+\exp[-\omega (\gamma^2 - \alpha)]}\]

The potential can be smoothed to zero force (making it purely attractive) when \(r < r_{\rm eq}\) by making \(U_{\rm M}(r < r_{\rm eq}) = -M_{\rm d}\) with the option repulsion = False.

Here, \(\vec{r}_{ij}\) is the displacement vector between particles \(i\) and \(j\), \(|\vec{r}_{ij}|\) is the magnitude of that displacement, and \(\hat{n}_i\) is the normalized orientation vector of particle \(i\). The parameters \(M_{\rm d}\), \(M_{\rm r}\), and \(r_{\rm eq}\) control the depth, width, and position of the potential well. The parameters \(\alpha\) and \(\omega\) control the width and steepness of the orientation dependence.

Example:

nl = hoomd.md.nlist.Cell()
m2p = azplugins.pair.TwoPatchMorse(nlist=nl)
m2p.params[('A', 'A')] = dict(M_d=1.8347, M_r=0.0302, r_eq=1.0043,
    omega=20, alpha=0.50, repulsion=True)
m2p.r_cut[('A', 'A')] = 3.0
params

The two-patch Morse potential parameters. The dictionary has the following keys:

  • M_d (float, required) - \(M_{\rm d}\), controls the depth of the potential well \([\mathrm{energy}]\)

  • M_r (float, required) - \(M_r\) controls the width of the potential well \([\mathrm{length}]\)

  • r_eq (float, required) - \(r_{\rm eq}\) controls the position of the potential well \([\mathrm{length}]\)

  • omega (float, required) - \(\omega\) controls the steepness of the orientation dependence

  • alpha (float, required) - \(\alpha\) controls the width of the orientation depndence

  • repulsion (bool, required) - If True, include the repulsive part of \(U_{\rm M}\) for \(r < r_{\rm eq}\). Otherwise, set \(U_{\rm r} = -M_{\rm d}\) for \(r < r_{\rm eq}\).

Type: TypeParameter [tuple [particle_type, particle_type], dict]