md.pair.aniso¶
Overview
Anistropic LJ force. |
|
Base class anisotropic pair force. |
|
Screened dipole-dipole pair forces. |
|
Gay-Berne anisotropic pair force. |
|
Patchy pair potentials. |
|
Modulate |
|
Modulate |
|
Modulate |
|
Modulate |
|
Modulate |
|
Modulate |
|
Modulate |
Details
Anisotropic pair forces.
Anisotropic pair force classes apply a force, torque, and virial on every particle in the simulation state commensurate with the potential energy:
AnisotropicPair
applies cuttoffs, exclusions, and assigns per particle
energies and virials in the same manner as hoomd.md.pair.Pair
AnisotropicPair
does not support the 'xplor'
shifting mode or the r_on
parameter.
- class hoomd.md.pair.aniso.ALJ(nlist, default_r_cut=None)¶
Bases:
AnisotropicPair
Anistropic LJ force.
- Parameters:
nlist (hoomd.md.nlist.NeighborList) – Neighbor list
default_r_cut (float) – Default cutoff radius \([length]\).
ALJ
computes the Lennard-Jones force between anisotropic particles as described in Ramasubramani, V. et al. 2020, using the formula:\[U(r, r_c) = U_0(r) + U_c(r_c)\]The first term is the central interaction \(U_0\), the standard center-center interaction between two Lennard-Jones particles with center-center distance \(r\). The second term is the contact interaction \(U_c\), computed from the smallest distance between the surfaces of the two shapes \(r_c\). The central and contact interactions are defined as follows:
\[ \begin{align}\begin{aligned}&U_0(r) = 4 \varepsilon \left[ \left( \frac{\sigma}{r} \right)^{12} - \left( \frac{\sigma}{r} \right)^{6} \right]\\&U_c(r_c) = 4 \varepsilon_c(\varepsilon) \left[ \left( \frac{\sigma_c}{r_c} \right)^{12} - \left( \frac{\sigma_c}{r_c} \right)^{6} \right]\end{aligned}\end{align} \]where \(\varepsilon\) (
epsilon
) affects strength of both the central and contact interactions, \(\varepsilon_c\) is an energy coefficient set proportional to \(\varepsilon\) to preserve the shape, \(\sigma\) is the interaction distance of the central term computed as the average of \(\sigma_i\) (sigma_i
) and \(\sigma_j\) (sigma_j
). Lastly,ALJ
uses the contact ratios \(\beta_i\) (contact_ratio_i
) and \(\beta_j\) (contact_ratio_j
) to compute the contact sigma \(\sigma_c\) as follows:\[ \begin{align}\begin{aligned}\sigma_c &= \frac{1}{2} \left[\sigma_{ci} + \sigma_{cj} \right]\\\sigma_{ci} &= \beta_i \cdot \sigma_i\\\sigma_{cj} &= \beta_j \cdot \sigma_j\end{aligned}\end{align} \]The total potential energy is therefore the sum of two interactions, a central Lennard-Jones potential and a radially-shifted Lennard-Jones potential where the shift is anisotropic and depends on the extent of the shape in each direction.
Each term has an independent cutoff at which the energy is set to zero. The behavior of these cutoffs is dependent on whether a user requires LJ or Weeks-Chandler-Anderson (WCA)-like (repulsive-only) behavior. This behavior is controlled using the
alpha
parameter, which can take on the following values:¶ center-center repulsive only
center-center full-range
contact-contact repulsive only
alpha = 0
alpha = 1
contact-contact full-range
alpha = 2
alpha = 3
For polytopes, computing interactions using a single contact point leads to significant instabilities in the torques because the contact point can jump from one end of a face to another in an arbitrarily small time interval. To ameliorate this, the ALJ potential performs a local averaging over all the features associated with the closest simplices on two polytopes. This averaging can be turned off by setting the
average_simplices
key for the type pair toFalse
.Specifying only
rounding_radii
creates an ellipsoid, while specifying onlyvertices
creates a convex polytope (setvertices
andfaces
to empty lists to create the ellipsoid).Important
The repulsive part of the contact interaction \(U_c(r_c)\) prevents two
ALJ
particles from approaching closely, effectively rounding the shape by a radius \(\sigma_c\). For this reason, the shape written bytype_shapes
includes the rounding due torounding_radii
and that due to \(\sigma_c\).Choosing
r_cut
:Set
r_cut
for each pair of particle types so thatALJ
can compute interactions for all possible relative placements and orientations of the particles. The farthest apart two particles can be while still interacting depends on the value ofalpha
.In the following list, the first argument to the \(\max\) function is for the center-center interaction. The second argument is for the contact-contact interaction, where \(R_i\) is the radius of the shape’s minimal origin-centered bounding sphere of the particle with type \(i\).
Let \(\lambda_{min} = 2^{1/6}\) be the position of the potential energy minimum of the Lennard-Jones potential and \(\lambda_{cut}^{attractive}\) be a larger value, such as 2.5 (typically used in isotropic LJ systems).
For alpha=0:
\[\begin{split}r_{\mathrm{cut},ij} = \max \bigg( & \frac{\lambda_{min}}{2} (\sigma_i + \sigma_j), \\ & R_i + R_j + R_{\mathrm{rounding},i} + R_{\mathrm{rounding},j} + \frac{\lambda_{min}}{2} (\beta_i \cdot \sigma_i + \beta_j \cdot \sigma_j) \bigg)\end{split}\]For alpha=1:
\[\begin{split}r_{\mathrm{cut},ij} = \max \bigg( & \frac{\lambda_{cut}^{attractive}}{2} (\sigma_i + \sigma_j), \\ & R_i + R_j + R_{\mathrm{rounding},i} + R_{\mathrm{rounding},j}+ \frac{\lambda_{min}}{2} (\beta_i \cdot \sigma_i + \beta_j \cdot \sigma_j) \bigg)\end{split}\]For alpha=2:
\[\begin{split}r_{\mathrm{cut},ij} = \max \bigg( & \frac{\lambda_{min}}{2} (\sigma_i + \sigma_j)), \\ & R_i + R_j + R_{\mathrm{rounding},i} + R_{\mathrm{rounding},j} + \frac{\lambda_{cut}^{attractive}}{2} (\beta_i \cdot \sigma_i + \beta_j \cdot \sigma_j) \bigg)\end{split}\]For alpha=3:
\[\begin{split}r_{\mathrm{cut},ij} = \max \bigg( & \frac{\lambda_{cut}^{attractive}}{2} (\sigma_i + \sigma_j), \\ & R_i + R_j + R_{\mathrm{rounding},i} + R_{\mathrm{rounding},j} + \frac{\lambda_{cut}^{attractive}}{2} (\beta_i \cdot \sigma_i + \beta_j \cdot \sigma_j) \bigg)\end{split}\]
Warning
Changing dimension in a simulation will invalidate this force and will lead to error or unrealistic behavior.
- params¶
The ALJ potential parameters. The dictionary has the following keys:
epsilon
(float
, required) - base energy scale \(\varepsilon\) \([energy]\).sigma_i
(float
, required) - the insphere diameter of the first particle type, \(\sigma_i\) \([length]\).sigma_j
(float
, required) - the insphere diameter of the second particle type, \(\sigma_j\) \([length]\).alpha
(int
, required) - Integer 0-3 indicating whether or not to include the attractive component of the interaction (see above for details).contact_ratio_i
(float
, optional) - \(\beta_i\), the ratio of the contact sphere diameter of the first type withsigma_i
. Defaults to 0.15.contact_ratio_j
(float
, optional) - \(\beta_j\), the ratio of the contact sphere diameter of the second type withsigma_j
. Defaults to 0.15.average_simplices
(bool
, optional) - Whether to average over simplices. Defaults toTrue
. See class documentation for more information.
Type:
hoomd.data.typeparam.TypeParameter
[tuple
[particle_types
,particle_types
],dict
]
Note
While the evaluation of the potential is symmetric with respect to the potential parameter labels
i
andj
, the parameters which physically represent a specific particle type must appear in all sets of pair parameters which include that particle type.- shape¶
The shape of a given type. The dictionary has the following keys per type:
vertices
(list
[tuple
[float
,float
,float
]], required) - The vertices of a convex polytope in 2 or 3 dimensions. The third dimension in 2D is ignored.rounding_radii
(tuple
[float
,float
,float
] orfloat
) - The semimajor axes of a rounding ellipsoid \(R_{\mathrm{rounding},i}\). If a single value is specified, the rounding ellipsoid is a sphere. Defaults to (0.0, 0.0, 0.0).faces
(list
[list
[int
]], required) - The faces of the polyhedron specified as a list of list of integers. The indices corresponding to the vertices must be ordered counterclockwise with respect to the face normal vector pointing outward from the origin.
Type:
hoomd.data.typeparam.TypeParameter
[particle_types
,dict
]
Example:
nl = hoomd.md.nlist.Cell(buffer=0.4) alj = hoomd.md.pair.aniso.ALJ(nl, r_cut=2.5) cube_verts = [(-0.5, -0.5, -0.5), (-0.5, -0.5, 0.5), (-0.5, 0.5, -0.5), (-0.5, 0.5, 0.5), (0.5, -0.5, -0.5), (0.5, -0.5, 0.5), (0.5, 0.5, -0.5), (0.5, 0.5, 0.5)]; cube_faces = [[0, 2, 6], [6, 4, 0], [5, 0, 4], [5,1,0], [5,4,6], [5,6,7], [3,2,0], [3,0,1], [3,6,2], [3,7,6], [3,1,5], [3,5,7]] alj.params[("A", "A")] = dict(epsilon=2.0, sigma_i=1.0, sigma_j=1.0, alpha=1, ) alj.shape["A"] = dict(vertices=cube_verts, faces=cube_faces)
The following example shows how to easily get the faces, with vertex indices properly ordered, for a shape with known vertices by using the coxeter package:
Example:
import coxeter nl = hoomd.md.nlist.Cell(buffer=0.4) alj = hoomd.md.pair.aniso.ALJ(nl, r_cut=2.5) cube_verts = [[-0.5, -0.5, -0.5], [-0.5, -0.5, 0.5], [-0.5, 0.5, -0.5], [-0.5, 0.5, 0.5], [0.5, -0.5, -0.5], [0.5, -0.5, 0.5], [0.5, 0.5, -0.5], [0.5, 0.5, 0.5]] cube = coxeter.shapes.ConvexPolyhedron(cube_verts) alj.params[("A", "A")] = dict(epsilon=2.0, sigma_i=1.0, sigma_j=1.0, alpha=1, ) alj.shape["A"] = dict(vertices=cube.vertices, faces=cube.faces)
- property type_shapes¶
The shape specification for use with GSD files for visualization.
This is not meant to be used for access to shape information in Python. See the attribute
shape
for programatic assess. Use this property to log shape for visualization and storage through the GSD file type.(
Loggable
: category=”object”)
- class hoomd.md.pair.aniso.AnisotropicPair(nlist, default_r_cut=None, mode='none')¶
Bases:
Pair
Base class anisotropic pair force.
AnisotropicPair
is the base class for all anisotropic pair forces.Warning
This class should not be instantiated by users. The class can be used for
isinstance
orissubclass
checks.- Parameters:
nlist (hoomd.md.nlist.NeighborList) – The neighbor list.
default_r_cut (
float
, optional) – The default cutoff for the potential, defaults toNone
which means no cutoff \([\mathrm{length}]\).mode (
str
, optional) – the energy shifting mode, defaults to “none”.
- class hoomd.md.pair.aniso.Dipole(nlist, default_r_cut=None)¶
Bases:
AnisotropicPair
Screened dipole-dipole pair forces.
- Parameters:
nlist (hoomd.md.nlist.NeighborList) – Neighbor list
default_r_cut (float) – Default cutoff radius \([\mathrm{length}]\).
Dipole
computes the (screened) interaction between pairs of particles with dipoles and electrostatic charges:\[ \begin{align}\begin{aligned}U &= U_{dd} + U_{de} + U_{ee}\\U_{dd} &= A e^{-\kappa r} \left(\frac{\vec{\mu_i}\cdot\vec{\mu_j}}{r^3} - 3\frac{(\vec{\mu_i}\cdot \vec{r_{ji}}) (\vec{\mu_j}\cdot \vec{r_{ji}})} {r^5} \right)\\U_{de} &= A e^{-\kappa r} \left(\frac{(\vec{\mu_j}\cdot \vec{r_{ji}})q_i}{r^3} - \frac{(\vec{\mu_i}\cdot \vec{r_{ji}})q_j}{r^3} \right)\\U_{ee} &= A e^{-\kappa r} \frac{q_i q_j}{r}\end{aligned}\end{align} \]Note
All units are documented electronic dipole moments. However,
Dipole
can also be used to represent magnetic dipoles.Example:
nl = nlist.Cell() dipole = md.pair.ansio.Dipole(nl, default_r_cut=3.0) dipole.params[('A', 'B')] = dict(A=1.0, kappa=4.0) dipole.mu['A'] = (4.0, 1.0, 0.0)
- params¶
The dipole potential parameters. The dictionary has the following keys:
A
(float
, required) - \(A\) - electrostatic energy scale (default: 1.0) \([\mathrm{energy} \cdot \mathrm{length} \cdot \mathrm{charge}^{-2}]\)kappa
(float
, required) - \(\kappa\) - inverse screening length \([\mathrm{length}^{-1}]\)
Type:
TypeParameter
[tuple
[particle_type
,particle_type
],dict
]
- class hoomd.md.pair.aniso.GayBerne(nlist, default_r_cut=None, mode='none')¶
Bases:
AnisotropicPair
Gay-Berne anisotropic pair force.
- Parameters:
nlist (hoomd.md.nlist.NeighborList) – Neighbor list
default_r_cut (float) – Default cutoff radius \([\mathrm{length}]\).
mode (str) – energy shifting/smoothing mode.
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 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)
- 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
]
- class hoomd.md.pair.aniso.Patchy(nlist, default_r_cut=None, mode='none')¶
Bases:
AnisotropicPair
Patchy pair potentials.
Patchy
combines an isotropicPair
potential with an orientation dependent modulation function \(f\). UsePatchy
with an attractive potential to create localized sticky patches on the surface of a particle. Use it with a repulsive potential to create localized bumps.Patchy
computes both forces and torques on particles.Note
Patchy
provides no interaction when there are no patches or particles are oriented such that \(f = 0\). UsePatchy
along with a repulsive isotropicPair
potential to prevent particles from passing through each other.The specific form of the patchy pair potential between particles \(i\) and \(j\) is:
\[U(r_{ij}, \mathbf{q}_i, \mathbf{q}_j) = \sum_{m=1}^{N_{\mathrm{patches},i}} \sum_{n=1}^{N_{\mathrm{patches},j}} f(\theta_{m,i}, \alpha, \omega) f(\theta_{n,j}, \alpha, \omega) U_{\mathrm{pair}}(r_{ij})\]where \(U_{\mathrm{pair}}(r_{ij})\) is the isotropic pair potential and \(f\) is an orientation-dependent factor of the patchy spherical cap half-angle \(\alpha\) and patch steepness \(\omega\):
\[\begin{split}\begin{align} f(\theta, \alpha, \omega) &= \frac{\big(1+e^{-\omega (\cos{\theta} - \cos{\alpha}) }\big)^{-1} - f_{min}}{f_{max} - f_{min}}\\ f_{max} &= \big( 1 + e^{-\omega (1 - \cos{\alpha}) } \big)^{-1} \\ f_{min} &= \big( 1 + e^{-\omega (-1 - \cos{\alpha}) } \big)^{-1} \\ \end{align}\end{split}\]directors
sets the locations of the patches in the local reference frame of the particle.Patchy
rotates the local director \(\vec{d}\) by the particle’s orientation and computes the \(cos(\theta)\) in \(f\) as the cosine of the angle between the director and \(\vec{r}_{ij}\): \(cos(\theta_i) = \mathbf{q} \hat{d} \mathbf{q}^* \cdot \hat{r}_{ij}\) and \(cos(\theta_j) = \mathbf{q} \hat{d} \mathbf{q}^* \cdot -\hat{r}_{ij}\).\(\alpha\) and \(\omega\) control the shape of the patch:
See also
hoomd.hpmc.pair.AngularStep
provides a similar functional form for HPMC, except that \(f\) is a step function.Warning
This class should not be instantiated by users. The class can be used for
isinstance
orissubclass
checks.- params¶
The Patchy potential parameters unique to each pair of particle types. The dictionary has the following keys:
envelope_params
(dict
, required)pair_params
(dict
, required) - passed to isotropic potential (see subclasses).
Example:
envelope_params = {'alpha': math.pi/4, 'omega': 30} patchy.params[('A', 'A')] = dict(pair_params=pair_params, envelope_params=envelope_params)
Type:
TypeParameter
[tuple
[particle_type
,particle_type
],dict
]
- directors¶
List of vectors pointing to patch centers, by particle type (normalized when set). When a particle type does not have patches, set an empty list.
Type:
TypeParameter
[particle_type
,list
[tuple
[float
,float
,float
]]Examples:
patchy.directors['A'] = [(1,0,0), (1,1,1)]
patchy.directors['A'] = []
- class hoomd.md.pair.aniso.PatchyExpandedGaussian(nlist, default_r_cut=None, mode='none')¶
Bases:
Patchy
Modulate
hoomd.md.pair.ExpandedGaussian
with angular patches.- Parameters:
nlist (hoomd.md.nlist.NeighborList) – Neighbor list
default_r_cut (float) – Default cutoff radius \([\mathrm{length}]\).
mode (str) – energy shifting/smoothing mode.
Example:
gauss_params=dict(epsilon=1, sigma=1, delta=0.5) envelope_params=dict(alpha=math.pi/2, omega=40) patchy_expanded_gaussian = hoomd.md.pair.aniso.PatchyExpandedGaussian( nlist=neighbor_list, default_r_cut=3.0) patchy_expanded_gaussian.params[('A', 'A')] = dict( pair_params=gauss_params, envelope_params=envelope_params) patchy_expanded_gaussian.directors['A'] = [(1,0,0), (1,1,1)] simulation.operations.integrator.forces = [patchy_expanded_gaussian]
- params¶
The Patchy potential parameters unique to each pair of particle types. The dictionary has the following keys:
envelope_params
(dict
, required)pair_params
(dict
, required) - passed tomd.pair.ExpandedGaussian.params
.
Type:
TypeParameter
[tuple
[particle_type
,particle_type
],dict
]
Inherited from:
Patchy
- directors¶
List of vectors pointing to patch centers, by particle type.
- class hoomd.md.pair.aniso.PatchyExpandedLJ(nlist, default_r_cut=None, mode='none')¶
Bases:
Patchy
Modulate
hoomd.md.pair.ExpandedLJ
with angular patches.- Parameters:
nlist (hoomd.md.nlist.NeighborList) – Neighbor list
default_r_cut (float) – Default cutoff radius \([\mathrm{length}]\).
mode (str) – energy shifting/smoothing mode.
lj_params=dict(epsilon=1, sigma=1) envelope_params=dict(alpha=math.pi/2, omega=20) patchylj = hoomd.md.pair.aniso.PatchyLJ(nlist=neighbor_list, default_r_cut=3.0) patchylj.params[('A', 'A')] = dict(pair_params=lj_params, envelope_params=envelope_params) patchylj.directors['A'] = [(1,0,0)] simulation.operation.integrator.forces = [patchylj]
- params¶
The Patchy potential parameters unique to each pair of particle types. The dictionary has the following keys:
envelope_params
(dict
, required)pair_params
(dict
, required) - passed tomd.pair.ExpandedLJ.params
.
Type:
TypeParameter
[tuple
[particle_type
,particle_type
],dict
]
Inherited from:
Patchy
- directors¶
List of vectors pointing to patch centers, by particle type.
- class hoomd.md.pair.aniso.PatchyExpandedMie(nlist, default_r_cut=None, mode='none')¶
Bases:
Patchy
Modulate
hoomd.md.pair.ExpandedMie
with angular patches.- Parameters:
nlist (hoomd.md.nlist.NeighborList) – Neighbor list
default_r_cut (float) – Default cutoff radius \([\mathrm{length}]\).
mode (str) – energy shifting/smoothing mode.
Example:
expanded_mie_params = dict(epsilon=1, sigma=1, n=15, m=10, delta=1) envelope_params = dict(alpha=math.pi/3, omega=20) patchy_expanded_mie = hoomd.md.pair.aniso.PatchyExpandedMie( nlist=neighbor_list, default_r_cut=3.0) patchy_expanded_mie.params[('A', 'A')] = dict( pair_params=expanded_mie_params envelope_params=envelope_params) patchy_expanded_mie.directors['A'] = [(1,0,0)] simulation.operations.integrator.forces = [patchy_expanded_mie]
- params¶
The Patchy potential parameters unique to each pair of particle types. The dictionary has the following keys:
envelope_params
(dict
, required)pair_params
(dict
, required) - passed tomd.pair.ExpandedMie.params
.
Type:
TypeParameter
[tuple
[particle_type
,particle_type
],dict
]
Inherited from:
Patchy
- directors¶
List of vectors pointing to patch centers, by particle type.
- class hoomd.md.pair.aniso.PatchyGaussian(nlist, default_r_cut=None, mode='none')¶
Bases:
Patchy
Modulate
hoomd.md.pair.Gaussian
with angular patches.- Parameters:
nlist (hoomd.md.nlist.NeighborList) – Neighbor list
default_r_cut (float) – Default cutoff radius \([\mathrm{length}]\).
mode (str) – energy shifting/smoothing mode.
Example:
gauss_params=dict(epsilon=1, sigma=1) envelope_params=dict(alpha=math.pi/4, omega=30) patchy_gaussian = hoomd.md.pair.aniso.PatchyGaussian(nlist=neighbor_list, default_r_cut=3.0) patchy_gaussian.params[('A', 'A')] = dict(pair_params=gauss_params, envelope_params=envelope_params) patchy_gaussian.directors['A'] = [(1,0,0)] simulation.operations.integrator.forces = [patchy_gaussian]
- params¶
The Patchy potential parameters unique to each pair of particle types. The dictionary has the following keys:
envelope_params
(dict
, required)pair_params
(dict
, required) - passed tomd.pair.ExpandedMie.params
.
Type:
TypeParameter
[tuple
[particle_type
,particle_type
],dict
]
Inherited from:
Patchy
- directors¶
List of vectors pointing to patch centers, by particle type.
- class hoomd.md.pair.aniso.PatchyLJ(nlist, default_r_cut=None, mode='none')¶
Bases:
Patchy
Modulate
hoomd.md.pair.LJ
with angular patches.- Parameters:
nlist (hoomd.md.nlist.NeighborList) – Neighbor list
default_r_cut (float) – Default cutoff radius \([\mathrm{length}]\).
mode (str) – energy shifting/smoothing mode.
Example:
lj_params = dict(epsilon=1, sigma=1) envelope_params=dict(alpha=math.pi/2, omega=20) patchylj = hoomd.md.pair.aniso.PatchyLJ(nlist=neighbor_list, default_r_cut=3.0) patchylj.params[('A', 'A')] = dict(pair_params=lj_params, envelope_params=envelope_params) patchylj.directors['A'] = [(1,0,0)] simulation.operations.integrator.forces = [patchylj]
- params¶
The Patchy potential parameters unique to each pair of particle types. The dictionary has the following keys:
envelope_params
(dict
, required)pair_params
(dict
, required) - passed tomd.pair.LJ.params
.
Type:
TypeParameter
[tuple
[particle_type
,particle_type
],dict
]
Inherited from:
Patchy
- directors¶
List of vectors pointing to patch centers, by particle type.
- class hoomd.md.pair.aniso.PatchyMie(nlist, default_r_cut=None, mode='none')¶
Bases:
Patchy
Modulate
hoomd.md.pair.Mie
with angular patches.- Parameters:
nlist (hoomd.md.nlist.NeighborList) – Neighbor list
default_r_cut (float) – Default cutoff radius \([\mathrm{length}]\).
mode (str) – energy shifting/smoothing mode.
Example:
mie_params = dict(epsilon=1, sigma=1, n=15, m=10) envelope_params = dict(alpha=math.pi/3, omega=20) patchy_mie = hoomd.md.pair.aniso.PatchyMie(nlist=neighbor_list, default_r_cut=3.0) patchy_mie.params[('A', 'A')] = dict(pair_params=mie_params envelope_params = envelope_params) patchy_mie.directors['A'] = [(1,0,0)] simulation.operations.integrator.forces = [patchy_mie]
- params¶
The Patchy potential parameters unique to each pair of particle types. The dictionary has the following keys:
envelope_params
(dict
, required)pair_params
(dict
, required) - passed tomd.pair.Mie.params
.
Type:
TypeParameter
[tuple
[particle_type
,particle_type
],dict
]
Inherited from:
Patchy
- directors¶
List of vectors pointing to patch centers, by particle type.
- class hoomd.md.pair.aniso.PatchyYukawa(nlist, default_r_cut=None, mode='none')¶
Bases:
Patchy
Modulate
hoomd.md.pair.Yukawa
with angular patches.- Parameters:
nlist (hoomd.md.nlist.NeighborList) – Neighbor list
default_r_cut (float) – Default cutoff radius \([\mathrm{length}]\).
mode (str) – energy shifting/smoothing mode.
Example:
yukawa_params = dict(epsilon=1, kappa=10) envelope_params = dict(alpha=math.pi/4, omega=25) patchy_yukawa = hoomd.md.pair.aniso.PatchyYukawa(nlist=neighbor_list, default_r_cut=5.0) patchy_yukawa.params[('A', 'A')] = dict(pair_params=yukawa_params envelope_params=envelope_params) patchy_yukawa.directors['A'] = [(1,0,0)] simulation.operations.integrator.forces = [patchy_yukawa]
- params¶
The Patchy potential parameters unique to each pair of particle types. The dictionary has the following keys:
envelope_params
(dict
, required)pair_params
(dict
, required) - passed tomd.pair.Yukawa.params
.
Type:
TypeParameter
[tuple
[particle_type
,particle_type
],dict
]
Inherited from:
Patchy
- directors¶
List of vectors pointing to patch centers, by particle type.