1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243
|
<HTML>
<CENTER><A HREF = "http://lammps.sandia.gov">LAMMPS WWW Site</A> - <A HREF = "Manual.html">LAMMPS Documentation</A> - <A HREF = "Section_commands.html#comm">LAMMPS Commands</A>
</CENTER>
<HR>
<H3>pair_style resquared command
</H3>
<H3>pair_style resquared/gpu command
</H3>
<H3>pair_style resquared/omp command
</H3>
<P><B>Syntax:</B>
</P>
<PRE>pair_style resquared cutoff
</PRE>
<UL><LI>cutoff = global cutoff for interactions (distance units)
</UL>
<P><B>Examples:</B>
</P>
<PRE>pair_style resquared 10.0
pair_coeff * * 1.0 1.0 1.7 3.4 3.4 1.0 1.0 1.0
</PRE>
<P><B>Description:</B>
</P>
<P>Style <I>resquared</I> computes the RE-squared anisotropic interaction
<A HREF = "#Everaers">(Everaers)</A>, <A HREF = "#Babadi">(Babadi)</A> between pairs of
ellipsoidal and/or spherical Lennard-Jones particles. For ellipsoidal
interactions, the potential considers the ellipsoid as being comprised
of small spheres of size sigma. LJ particles are a single sphere of
size sigma. The distinction is made to allow the pair style to make
efficient calculations of ellipsoid/solvent interactions.
</P>
<P>Details for the equations used are given in the references below and
in <A HREF = "PDF/pair_resquared_extra.pdf">this supplementary document</A>.
</P>
<P>Use of this pair style requires the NVE, NVT, or NPT fixes with the
<I>asphere</I> extension (e.g. <A HREF = "fix_nve_asphere.html">fix nve/asphere</A>) in
order to integrate particle rotation. Additionally, <A HREF = "atom_style.html">atom_style
ellipsoid</A> should be used since it defines the
rotational state and the size and shape of each ellipsoidal particle.
</P>
<P>The following coefficients must be defined for each pair of atoms
types via the <A HREF = "pair_coeff.html">pair_coeff</A> command as in the examples
above, or in the data file or restart files read by the
<A HREF = "read_data.html">read_data</A> or <A HREF = "read_restart.html">read_restart</A>
commands:
</P>
<UL><LI>A12 = Energy Prefactor/Hamaker constant (energy units)
<LI>sigma = atomic interaction diameter (distance units)
<LI>epsilon_i_a = relative well depth of type I for side-to-side interactions
<LI>epsilon_i_b = relative well depth of type I for face-to-face interactions
<LI>epsilon_i_c = relative well depth of type I for end-to-end interactions
<LI>epsilon_j_a = relative well depth of type J for side-to-side interactions
<LI>epsilon_j_b = relative well depth of type J for face-to-face interactions
<LI>epsilon_j_c = relative well depth of type J for end-to-end interactions
<LI>cutoff (distance units)
</UL>
<P>The parameters used depend on the type of the interacting particles,
i.e. ellipsoids or LJ spheres. The type of a particle is determined
by the diameters specified for its 3 shape paramters. If all 3 shape
parameters = 0.0, then the particle is treated as an LJ sphere. The
epsilon_i_* or epsilon_j_* parameters are ignored for LJ spheres. If
the 3 shape paraemters are > 0.0, then the particle is treated as an
ellipsoid (even if the 3 parameters are equal to each other).
</P>
<P>A12 specifies the energy prefactor which depends on the types of the
two interacting particles.
</P>
<P>For ellipsoid/ellipsoid interactions, the interaction is computed by
the formulas in the supplementary docuement referenced above. A12 is
the Hamaker constant as described in <A HREF = "#Everaers">(Everaers)</A>. In LJ
units:
</P>
<CENTER><IMG SRC = "Eqs/pair_resquared.jpg">
</CENTER>
<P>where rho gives the number density of the spherical particles
composing the ellipsoids and epsilon_LJ determines the interaction
strength of the spherical particles.
</P>
<P>For ellipsoid/LJ sphere interactions, the interaction is also computed
by the formulas in the supplementary docuement referenced above. A12
has a modifed form (see <A HREF = "PDF/pair_resquared_extra.pdf">here</A> for
details):
</P>
<CENTER><IMG SRC = "Eqs/pair_resquared2.jpg">
</CENTER>
<P>For ellipsoid/LJ sphere interactions, a correction to the distance-
of-closest approach equation has been implemented to reduce the error
from two particles of disparate sizes; see <A HREF = "PDF/pair_resquared_extra.pdf">this supplementary
document</A>.
</P>
<P>For LJ sphere/LJ sphere interactions, the interaction is computed
using the standard Lennard-Jones formula, which is much cheaper to
compute than the ellipsoidal formulas. A12 is used as epsilon in the
standard LJ formula:
</P>
<CENTER><IMG SRC = "Eqs/pair_resquared3.jpg">
</CENTER>
<P>and the specified <I>sigma</I> is used as the sigma in the standard LJ
formula.
</P>
<P>When one of both of the interacting particles are ellipsoids, then
<I>sigma</I> specifies the diameter of the continuous distribution of
constituent particles within each ellipsoid used to model the
RE-squared potential. Note that this is a different meaning for
<I>sigma</I> than the <A HREF = "pair_gayberne.html">pair_style gayberne</A> potential
uses.
</P>
<P>The epsilon_i and epsilon_j coefficients are defined for atom types,
not for pairs of atom types. Thus, in a series of pair_coeff
commands, they only need to be specified once for each atom type.
</P>
<P>Specifically, if any of epsilon_i_a, epsilon_i_b, epsilon_i_c are
non-zero, the three values are assigned to atom type I. If all the
epsilon_i values are zero, they are ignored. If any of epsilon_j_a,
epsilon_j_b, epsilon_j_c are non-zero, the three values are assigned
to atom type J. If all three epsilon_i values are zero, they are
ignored. Thus the typical way to define the epsilon_i and epsilon_j
coefficients is to list their values in "pair_coeff I J" commands when
I = J, but set them to 0.0 when I != J. If you do list them when I !=
J, you should insure they are consistent with their values in other
pair_coeff commands.
</P>
<P>Note that if this potential is being used as a sub-style of
<A HREF = "pair_hybrid.html">pair_style hybrid</A>, and there is no "pair_coeff I I"
setting made for RE-squared for a particular type I (because I-I
interactions are computed by another hybrid pair potential), then you
still need to insure the epsilon a,b,c coefficients are assigned to
that type in a "pair_coeff I J" command.
</P>
<P>For large uniform molecules it has been shown that the epsilon_*_*
energy parameters are approximately representable in terms of local
contact curvatures <A HREF = "#Everaers">(Everaers)</A>:
</P>
<CENTER><IMG SRC = "Eqs/pair_resquared4.jpg">
</CENTER>
<P>where a, b, and c give the particle diameters.
</P>
<P>The last coefficient is optional. If not specified, the global cutoff
specified in the pair_style command is used.
</P>
<HR>
<P>Styles with a <I>cuda</I>, <I>gpu</I>, <I>omp</I>, or <I>opt</I> suffix are functionally
the same as the corresponding style without the suffix. They have
been optimized to run faster, depending on your available hardware, as
discussed in <A HREF = "Section_accelerate.html">Section_accelerate</A> of the
manual. The accelerated styles take the same arguments and should
produce the same results, except for round-off and precision issues.
</P>
<P>These accelerated styles are part of the USER-CUDA, GPU, USER-OMP and OPT
packages, respectively. They are only enabled if LAMMPS was built with
those packages. See the <A HREF = "Section_start.html#start_3">Making LAMMPS</A>
section for more info.
</P>
<P>You can specify the accelerated styles explicitly in your input script
by including their suffix, or you can use the <A HREF = "Section_start.html#start_7">-suffix command-line
switch</A> when you invoke LAMMPS, or you can
use the <A HREF = "suffix.html">suffix</A> command in your input script.
</P>
<P>See <A HREF = "Section_accelerate.html">Section_accelerate</A> of the manual for
more instructions on how to use the accelerated styles effectively.
</P>
<HR>
<P><B>Mixing, shift, table, tail correction, restart, rRESPA info</B>:
</P>
<P>For atom type pairs I,J and I != J, the epsilon and sigma coefficients
and cutoff distance can be mixed, but only for sphere pairs. The
default mix value is <I>geometric</I>. See the "pair_modify" command for
details. Other type pairs cannot be mixed, due to the different
meanings of the energy prefactors used to calculate the interactions
and the implicit dependence of the ellipsoid-sphere interaction on the
equation for the Hamaker constant presented here. Mixing of sigma and
epsilon followed by calculation of the energy prefactors using the
equations above is recommended.
</P>
<P>This pair styles supports the <A HREF = "pair_modify.html">pair_modify</A> shift
option for the energy of the Lennard-Jones portion of the pair
interaction, but only for sphere-sphere interactions. There is no
shifting performed for ellipsoidal interactions due to the anisotropic
dependence of the interaction.
</P>
<P>The <A HREF = "pair_modify.html">pair_modify</A> table option is not relevant
for this pair style.
</P>
<P>This pair style does not support the <A HREF = "pair_modify.html">pair_modify</A>
tail option for adding long-range tail corrections to energy and
pressure.
</P>
<P>This pair style writes its information to <A HREF = "restart.html">binary restart
files</A>, so pair_style and pair_coeff commands do not need
to be specified in an input script that reads a restart file.
</P>
<P>This pair style can only be used via the <I>pair</I> keyword of the
<A HREF = "run_style.html">run_style respa</A> command. It does not support the
<I>inner</I>, <I>middle</I>, <I>outer</I> keywords of the <A HREF = "run_style.html">run_style
command</A>.
</P>
<HR>
<P><B>Restrictions:</B>
</P>
<P>This style is part of the ASPHERE package. It is only enabled if
LAMMPS was built with that package. See the <A HREF = "Section_start.html#start_3">Making
LAMMPS</A> section for more info.
</P>
<P>This pair style requires that atoms be ellipsoids as defined by the
<A HREF = "atom_style.html">atom_style ellipsoid</A> command.
</P>
<P>Particles acted on by the potential can be finite-size aspherical or
spherical particles, or point particles. Spherical particles have all
3 of their shape parameters equal to each other. Point particles have
all 3 of their shape parameters equal to 0.0.
</P>
<P>The distance-of-closest-approach approximation used by LAMMPS becomes
less accurate when high-aspect ratio ellipsoids are used.
</P>
<P><B>Related commands:</B>
</P>
<P><A HREF = "pair_coeff.html">pair_coeff</A>, <A HREF = "fix_nve_asphere.html">fix nve/asphere</A>,
<A HREF = "compute_temp_asphere.html">compute temp/asphere</A>, <A HREF = "pair_gayberne.html">pair_style
gayberne</A>
</P>
<P><B>Default:</B> none
</P>
<HR>
<A NAME = "Everaers"></A>
<P><B>(Everaers)</B> Everaers and Ejtehadi, Phys Rev E, 67, 041710 (2003).
</P>
<A NAME = "Babadi"></A>
<P><B>(Berardi)</B> Babadi, Ejtehadi, Everaers, J Comp Phys, 219, 770-779 (2006).
</P>
</HTML>
|