forked from lijiext/lammps
345 lines
14 KiB
HTML
345 lines
14 KiB
HTML
<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 lj/cut command
|
|
</H3>
|
|
<H3>pair_style lj/cut/cuda command
|
|
</H3>
|
|
<H3>pair_style lj/cut/experimental/cuda command
|
|
</H3>
|
|
<H3>pair_style lj/cut/gpu command
|
|
</H3>
|
|
<H3>pair_style lj/cut/opt command
|
|
</H3>
|
|
<H3>pair_style lj/cut/omp command
|
|
</H3>
|
|
<H3>pair_style lj/cut/coul/cut command
|
|
</H3>
|
|
<H3>pair_style lj/cut/coul/cut/cuda command
|
|
</H3>
|
|
<H3>pair_style lj/cut/coul/cut/gpu command
|
|
</H3>
|
|
<H3>pair_style lj/cut/coul/cut/omp command
|
|
</H3>
|
|
<H3>pair_style lj/cut/coul/debye command
|
|
</H3>
|
|
<H3>pair_style lj/cut/coul/debye/cuda command
|
|
</H3>
|
|
<H3>pair_style lj/cut/coul/debye/gpu command
|
|
</H3>
|
|
<H3>pair_style lj/cut/coul/debye/omp command
|
|
</H3>
|
|
<H3>pair_style lj/cut/coul/dsf command
|
|
</H3>
|
|
<H3>pair_style lj/cut/coul/dsf/gpu command
|
|
</H3>
|
|
<H3>pair_style lj/cut/coul/long command
|
|
</H3>
|
|
<H3>pair_style lj/cut/coul/long/cuda command
|
|
</H3>
|
|
<H3>pair_style lj/cut/coul/long/gpu command
|
|
</H3>
|
|
<H3>pair_style lj/cut/coul/long/opt command
|
|
</H3>
|
|
<H3>pair_style lj/cut/coul/long/omp command
|
|
</H3>
|
|
<H3>pair_style lj/cut/coul/msm command
|
|
</H3>
|
|
<H3>pair_style lj/cut/coul/msm/omp command
|
|
</H3>
|
|
<H3>pair_style lj/cut/tip4p/cut command
|
|
</H3>
|
|
<H3>pair_style lj/cut/tip4p/cut/omp command
|
|
</H3>
|
|
<H3>pair_style lj/cut/tip4p/long command
|
|
</H3>
|
|
<H3>pair_style lj/cut/tip4p/long/omp command
|
|
</H3>
|
|
<H3>pair_style lj/cut/tip4p/long/opt command
|
|
</H3>
|
|
<P><B>Syntax:</B>
|
|
</P>
|
|
<PRE>pair_style style args
|
|
</PRE>
|
|
<UL><LI>style = <I>lj/cut</I> or <I>lj/cut/coul/cut</I> or <I>lj/cut/coul/debye</I> or <I>lj/cut/coul/dsf</I> or <I>lj/cut/coul/long</I> or <I>lj/cut/coul/msm</I> or <I>lj/cut/tip4p/long</I>
|
|
<LI>args = list of arguments for a particular style
|
|
</UL>
|
|
<PRE> <I>lj/cut</I> args = cutoff
|
|
cutoff = global cutoff for Lennard Jones interactions (distance units)
|
|
<I>lj/cut/coul/cut</I> args = cutoff (cutoff2)
|
|
cutoff = global cutoff for LJ (and Coulombic if only 1 arg) (distance units)
|
|
cutoff2 = global cutoff for Coulombic (optional) (distance units)
|
|
<I>lj/cut/coul/debye</I> args = kappa cutoff (cutoff2)
|
|
kappa = inverse of the Debye length (inverse distance units)
|
|
cutoff = global cutoff for LJ (and Coulombic if only 1 arg) (distance units)
|
|
cutoff2 = global cutoff for Coulombic (optional) (distance units)
|
|
<I>lj/cut/coul/dsf</I> args = alpha cutoff (cutoff2)
|
|
alpha = damping parameter (inverse distance units)
|
|
cutoff = global cutoff for LJ (and Coulombic if only 1 arg) (distance units)
|
|
cutoff2 = global cutoff for Coulombic (distance units)
|
|
<I>lj/cut/coul/long</I> args = cutoff (cutoff2)
|
|
cutoff = global cutoff for LJ (and Coulombic if only 1 arg) (distance units)
|
|
cutoff2 = global cutoff for Coulombic (optional) (distance units)
|
|
<I>lj/cut/coul/msm</I> args = cutoff (cutoff2)
|
|
cutoff = global cutoff for LJ (and Coulombic if only 1 arg) (distance units)
|
|
cutoff2 = global cutoff for Coulombic (optional) (distance units)
|
|
<I>lj/cut/tip4p/cut</I> args = otype htype btype atype qdist cutoff (cutoff2)
|
|
otype,htype = atom types for TIP4P O and H
|
|
btype,atype = bond and angle types for TIP4P waters
|
|
qdist = distance from O atom to massless charge (distance units)
|
|
cutoff = global cutoff for LJ (and Coulombic if only 1 arg) (distance units)
|
|
cutoff2 = global cutoff for Coulombic (optional) (distance units)
|
|
<I>lj/cut/tip4p/long</I> args = otype htype btype atype qdist cutoff (cutoff2)
|
|
otype,htype = atom types for TIP4P O and H
|
|
btype,atype = bond and angle types for TIP4P waters
|
|
qdist = distance from O atom to massless charge (distance units)
|
|
cutoff = global cutoff for LJ (and Coulombic if only 1 arg) (distance units)
|
|
cutoff2 = global cutoff for Coulombic (optional) (distance units)
|
|
</PRE>
|
|
<P><B>Examples:</B>
|
|
</P>
|
|
<PRE>pair_style lj/cut 2.5
|
|
pair_coeff * * 1 1
|
|
pair_coeff 1 1 1 1.1 2.8
|
|
</PRE>
|
|
<PRE>pair_style lj/cut/coul/cut 10.0
|
|
pair_style lj/cut/coul/cut 10.0 8.0
|
|
pair_coeff * * 100.0 3.0
|
|
pair_coeff 1 1 100.0 3.5 9.0
|
|
pair_coeff 1 1 100.0 3.5 9.0 9.0
|
|
</PRE>
|
|
<PRE>pair_style lj/cut/coul/debye 1.5 3.0
|
|
pair_style lj/cut/coul/debye 1.5 2.5 5.0
|
|
pair_coeff * * 1.0 1.0
|
|
pair_coeff 1 1 1.0 1.5 2.5
|
|
pair_coeff 1 1 1.0 1.5 2.5 5.0
|
|
</PRE>
|
|
<PRE>pair_style lj/cut/coul/dsf 0.05 2.5 10.0
|
|
pair_coeff * * 1.0 1.0
|
|
pair_coeff 1 1 1.0 1.0 2.5
|
|
</PRE>
|
|
<PRE>pair_style lj/cut/coul/long 10.0
|
|
pair_style lj/cut/coul/long 10.0 8.0
|
|
pair_coeff * * 100.0 3.0
|
|
pair_coeff 1 1 100.0 3.5 9.0
|
|
</PRE>
|
|
<PRE>pair_style lj/cut/coul/msm 10.0
|
|
pair_style lj/cut/coul/msm 10.0 8.0
|
|
pair_coeff * * 100.0 3.0
|
|
pair_coeff 1 1 100.0 3.5 9.0
|
|
</PRE>
|
|
<PRE>pair_style lj/cut/tip4p/cut 1 2 7 8 0.15 12.0
|
|
pair_style lj/cut/tip4p/cut 1 2 7 8 0.15 12.0 10.0
|
|
pair_coeff * * 100.0 3.0
|
|
pair_coeff 1 1 100.0 3.5 9.0
|
|
</PRE>
|
|
<PRE>pair_style lj/cut/tip4p/long 1 2 7 8 0.15 12.0
|
|
pair_style lj/cut/tip4p/long 1 2 7 8 0.15 12.0 10.0
|
|
pair_coeff * * 100.0 3.0
|
|
pair_coeff 1 1 100.0 3.5 9.0
|
|
</PRE>
|
|
<P><B>Description:</B>
|
|
</P>
|
|
<P>The <I>lj/cut</I> styles compute the standard 12/6 Lennard-Jones potential,
|
|
given by
|
|
</P>
|
|
<CENTER><IMG SRC = "Eqs/pair_lj.jpg">
|
|
</CENTER>
|
|
<P>Rc is the cutoff.
|
|
</P>
|
|
<P>Style <I>lj/cut/coul/cut</I> adds a Coulombic pairwise interaction given by
|
|
</P>
|
|
<CENTER><IMG SRC = "Eqs/pair_coulomb.jpg">
|
|
</CENTER>
|
|
<P>where C is an energy-conversion constant, Qi and Qj are the charges on
|
|
the 2 atoms, and epsilon is the dielectric constant which can be set
|
|
by the <A HREF = "dielectric.html">dielectric</A> command. If one cutoff is
|
|
specified in the pair_style command, it is used for both the LJ and
|
|
Coulombic terms. If two cutoffs are specified, they are used as
|
|
cutoffs for the LJ and Coulombic terms respectively.
|
|
</P>
|
|
<P>Style <I>lj/cut/coul/debye</I> adds an additional exp() damping factor
|
|
to the Coulombic term, given by
|
|
</P>
|
|
<CENTER><IMG SRC = "Eqs/pair_debye.jpg">
|
|
</CENTER>
|
|
<P>where kappa is the inverse of the Debye length. This potential is
|
|
another way to mimic the screening effect of a polar solvent.
|
|
</P>
|
|
<P>Style <I>lj/cut/coul/dsf</I> computes the Coulombic term via the damped
|
|
shifted force model described in <A HREF = "#Fennell">Fennell</A>, given by:
|
|
</P>
|
|
<CENTER><IMG SRC = "Eqs/pair_coul_dsf.jpg">
|
|
</CENTER>
|
|
<P>where <I>alpha</I> is the damping parameter and erfc() is the complementary
|
|
error-function. This potential is essentially a short-range,
|
|
spherically-truncated, charge-neutralized, shifted, pairwise <I>1/r</I>
|
|
summation. The potential is based on Wolf summation, proposed as an
|
|
alternative to Ewald summation for condensed phase systems where
|
|
charge screening causes electrostatic interactions to become
|
|
effectively short-ranged. In order for the electrostatic sum to be
|
|
absolutely convergent, charge neutralization within the cutoff radius
|
|
is enforced by shifting the potential through placement of image
|
|
charges on the cutoff sphere. Convergence can often be improved by
|
|
setting <I>alpha</I> to a small non-zero value.
|
|
</P>
|
|
<P>Styles <I>lj/cut/coul/long</I> and <I>lj/cut/coul/msm</I> compute the same
|
|
Coulombic interactions as style <I>lj/cut/coul/cut</I> except that an
|
|
additional damping factor is applied to the Coulombic term so it can
|
|
be used in conjunction with the <A HREF = "kspace_style.html">kspace_style</A>
|
|
command and its <I>ewald</I> or <I>pppm</I> option. The Coulombic cutoff
|
|
specified for this style means that pairwise interactions within this
|
|
distance are computed directly; interactions outside that distance are
|
|
computed in reciprocal space.
|
|
</P>
|
|
<P>Styles <I>lj/cut/tip4p/cut</I> and <I>lj/cut/tip4p/long</I> implement the TIP4P
|
|
water model of <A HREF = "#Jorgensen">(Jorgensen)</A>, which introduces a massless
|
|
site located a short distance away from the oxygen atom along the
|
|
bisector of the HOH angle. The atomic types of the oxygen and
|
|
hydrogen atoms, the bond and angle types for OH and HOH interactions,
|
|
and the distance to the massless charge site are specified as
|
|
pair_style arguments. Style <I>lj/cut/tip4p/cut</I> uses a cutoff for
|
|
Coulomb interactions; style <I>lj/cut/tip4p/long</I> is for use with a
|
|
long-range Coulombic solver (Ewald or PPPM).
|
|
</P>
|
|
<P>IMPORTANT NOTE: For each TIP4P water molecule in your system, the atom
|
|
IDs for the O and 2 H atoms must be consecutive, with the O atom
|
|
first. This is to enable LAMMPS to "find" the 2 H atoms associated
|
|
with each O atom. For example, if the atom ID of an O atom in a TIP4P
|
|
water molecule is 500, then its 2 H atoms must have IDs 501 and 502.
|
|
</P>
|
|
<P>See the <A HREF = "Section_howto.html#howto_8">howto section</A> for more
|
|
information on how to use the TIP4P pair styles and lists of
|
|
parameters to set. Note that the neighobr list cutoff for Coulomb
|
|
interactions is effectively extended by a distance 2*qdist when using
|
|
the TIP4P pair style, to account for the offset distance of the
|
|
fictitious charges on O atoms in water molecules. Thus it is
|
|
typically best in an efficiency sense to use a LJ cutoff >= Coulomb
|
|
cutoff + 2*qdist, to shrink the size of the neighbor list. This leads
|
|
to slightly larger cost for the long-range calculation, so you can
|
|
test the trade-off for your model.
|
|
</P>
|
|
<P>For all of the <I>lj/cut</I> pair styles, 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, or by mixing as
|
|
described below:
|
|
</P>
|
|
<UL><LI>epsilon (energy units)
|
|
<LI>sigma (distance units)
|
|
<LI>cutoff1 (distance units)
|
|
<LI>cutoff2 (distance units)
|
|
</UL>
|
|
<P>Note that sigma is defined in the LJ formula as the zero-crossing
|
|
distance for the potential, not as the energy minimum at 2^(1/6)
|
|
sigma.
|
|
</P>
|
|
<P>The latter 2 coefficients are optional. If not specified, the global
|
|
LJ and Coulombic cutoffs specified in the pair_style command are used.
|
|
If only one cutoff is specified, it is used as the cutoff for both LJ
|
|
and Coulombic interactions for this type pair. If both coefficients
|
|
are specified, they are used as the LJ and Coulombic cutoffs for this
|
|
type pair. You cannot specify 2 cutoffs for style <I>lj/cut</I>, since it
|
|
has no Coulombic terms.
|
|
</P>
|
|
<P>For <I>lj/cut/coul/long</I> and <I>lj/cut/coul/msm</I> and <I>lj/cut/tip4p/cut</I>
|
|
and <I>lj/cut/tip4p/long</I> only the LJ cutoff can be specified since a
|
|
Coulombic cutoff cannot be specified for an individual I,J type pair.
|
|
All type pairs use the same global Coulombic cutoff specified in the
|
|
pair_style command.
|
|
</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 for all of the lj/cut pair styles can be mixed.
|
|
The default mix value is <I>geometric</I>. See the "pair_modify" command
|
|
for details.
|
|
</P>
|
|
<P>All of the <I>lj/cut</I> pair styles support the
|
|
<A HREF = "pair_modify.html">pair_modify</A> shift option for the energy of the
|
|
Lennard-Jones portion of the pair interaction.
|
|
</P>
|
|
<P>The <I>lj/cut/coul/long</I> and <I>lj/cut/tip4p/long</I> pair styles support the
|
|
<A HREF = "pair_modify.html">pair_modify</A> table option since they can tabulate
|
|
the short-range portion of the long-range Coulombic interaction.
|
|
</P>
|
|
<P>All of the <I>lj/cut</I> pair styles support the
|
|
<A HREF = "pair_modify.html">pair_modify</A> tail option for adding a long-range
|
|
tail correction to the energy and pressure for the Lennard-Jones
|
|
portion of the pair interaction.
|
|
</P>
|
|
<P>All of the <I>lj/cut</I> pair styles write their 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>The <I>lj/cut</I> and <I>lj/cut/coul/long</I> pair styles support the use of the
|
|
<I>inner</I>, <I>middle</I>, and <I>outer</I> keywords of the <A HREF = "run_style.html">run_style
|
|
respa</A> command, meaning the pairwise forces can be
|
|
partitioned by distance at different levels of the rRESPA hierarchy.
|
|
The other styles only support the <I>pair</I> keyword of run_style respa.
|
|
See the <A HREF = "run_style.html">run_style</A> command for details.
|
|
</P>
|
|
<HR>
|
|
|
|
<P><B>Restrictions:</B>
|
|
</P>
|
|
<P>The <I>lj/cut/coul/long</I> and <I>lj/cut/tip4p/long</I> styles are part of the
|
|
KSPACE package. The <I>lj/cut/tip4p/cut</I> style is part of the MOLECULE
|
|
package. These styles 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. Note that the KSPACE and MOLECULE packages are
|
|
installed by default.
|
|
</P>
|
|
<P><B>Related commands:</B>
|
|
</P>
|
|
<P><A HREF = "pair_coeff.html">pair_coeff</A>
|
|
</P>
|
|
<P><B>Default:</B> none
|
|
</P>
|
|
<HR>
|
|
|
|
<A NAME = "Jorgensen"></A>
|
|
|
|
<P><B>(Jorgensen)</B> Jorgensen, Chandrasekhar, Madura, Impey, Klein, J Chem
|
|
Phys, 79, 926 (1983).
|
|
</P>
|
|
<A NAME = "Fennell"></A>
|
|
|
|
<P><B>(Fennell)</B> C. J. Fennell, J. D. Gezelter, J Chem Phys, 124,
|
|
234104 (2006).
|
|
</P>
|
|
</HTML>
|