GayBerne

class hoomd.md.pair.aniso.GayBerne(nlist, default_r_cut=None, mode='none')

Bases: AnisotropicPair

Gay-Berne anisotropic pair force.

Parameters:

GayBerne computes the Gay-Berne anisotropic pair force on every particle in the simulation state. This version of the Gay-Berne force supports identical pairs of uniaxial ellipsoids, with orientation-independent energy-well depth. The potential comes from the following paper Allen et. al. 2006.

\[\begin{split}U(\vec r, \vec e_i, \vec e_j) = \begin{cases} 4 \varepsilon \left[ \zeta^{-12} - \zeta^{-6} \right] & \zeta < \zeta_{\mathrm{cut}} \\ 0 & \zeta \ge \zeta_{\mathrm{cut}} \\ \end{cases}\end{split}\]

where

\[ \begin{align}\begin{aligned}\zeta &= \left(\frac{r-\sigma+\sigma_{\mathrm{min}}} {\sigma_{\mathrm{min}}}\right),\\\sigma^{-2} &= \frac{1}{2} \hat{\vec{r}} \cdot \vec{H^{-1}} \cdot \hat{\vec{r}},\\\vec{H} &= 2 \ell_\perp^2 \vec{1} + (\ell_\parallel^2 - \ell_\perp^2) (\vec{e_i} \otimes \vec{e_i} + \vec{e_j} \otimes \vec{e_j}),\end{aligned}\end{align} \]

and \(\sigma_{\mathrm{min}} = 2 \min(\ell_\perp, \ell_\parallel)\). The parallel direction is aligned with z axis in the particle’s reference frame.

The cut-off parameter \(r_{\mathrm{cut}}\) is defined for two particles oriented parallel along the long axis, i.e. \(\zeta_{\mathrm{cut}} = \left(\frac{r-\sigma_{\mathrm{max}} + \sigma_{\mathrm{min}}}{\sigma_{\mathrm{min}}}\right)\) where \(\sigma_{\mathrm{max}} = 2 \max(\ell_\perp, \ell_\parallel)\) .

The quantities \(\ell_\parallel\) and \(\ell_\perp\) denote the semi-axis lengths parallel and perpendicular to particle orientation.

Example:

nl = nlist.Cell()
gay_berne = md.pair.aniso.GayBerne(nlist=nl, default_r_cut=2.5)
gay_berne.params[('A', 'A')] = dict(epsilon=1.0, lperp=0.45, lpar=0.5)
gay_berne.r_cut[('A', 'B')] = 2 ** (1.0 / 6.0)

Members inherited from AutotunedObject:

property kernel_parameters

Kernel parameters. Read more...

property is_tuning_complete

Check if kernel parameter tuning is complete. Read more...

tune_kernel_parameters()

Start tuning kernel parameters. Read more...


Members inherited from Force:

additional_energy

Additional energy term. Read more...

additional_virial

Additional virial tensor term \(W_\mathrm{additional}\). Read more...

cpu_local_force_arrays

Local force arrays on the CPU. Read more...

energies

Energy contribution \(U_i\) from each particle. Read more...

energy

The potential energy \(U\) of the system from this force. Read more...

forces

The force \(\vec{F}_i\) applied to each particle. Read more...

gpu_local_force_arrays

Local force arrays on the GPU. Read more...

torques

The torque \(\vec{\tau}_i\) applied to each particle. Read more...

virials

Virial tensor contribution \(W_i\) from each particle. Read more...


Members inherited from Pair:

nlist

Neighbor list used to compute the pair force. Read more...

mode

Energy smoothing/cutoff mode. Read more...

r_cut

Cuttoff radius beyond which the energy and force are 0. Read more...

r_on

Radius at which the XPLOR smoothing function starts. Read more...

compute_energy()

Compute the energy between two sets of particles. Read more...


Members defined in GayBerne:

params

The Gay-Berne potential parameters. The dictionary has the following keys:

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

  • lperp (float, required) - \(\ell_\perp\) \([\mathrm{length}]\)

  • lpar (float, required) - \(\ell_\parallel\) \([\mathrm{length}]\)

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

property type_shapes

Get all the types of shapes in the current simulation.

Example

>>> gay_berne.type_shapes
[{'type': 'Ellipsoid', 'a': 1.0, 'b': 1.0, 'c': 1.5}]
Returns:

A list of dictionaries, one for each particle type in the system.

(Loggable: category=”object”)