lammps/doc/pair_gayberne.txt

246 lines
10 KiB
Plaintext
Executable File

"LAMMPS WWW Site"_lws - "LAMMPS Documentation"_ld - "LAMMPS Commands"_lc :c
:link(lws,http://lammps.sandia.gov)
:link(ld,Manual.html)
:link(lc,Section_commands.html#comm)
:line
pair_style gayberne command :h3
pair_style gayberne/gpu command :h3
[Syntax:]
pair_style gayberne gamma upsilon mu cutoff :pre
pair_style gayberne/gpu gpuflag gpunum gamma upsilon mu cutoff :pre
style = {gayberne} or {gayberne/gpu}
gpuflag = {one/node} or {one/gpu} or {multi/gpu}, only used with gayberne/gpu
gpunum = gpu ID or quantity, only used with gayberne/gpu
gamma = shift for potential minimum (typically 1)
upsilon = exponent for eta orientation-dependent energy function
mu = exponent for chi orientation-dependent energy function
cutoff = global cutoff for interactions (distance units) :ul
[Examples:]
pair_style gayberne 1.0 1.0 1.0 10.0
pair_style gayberne/gpu 1.0 1.0 1.0 10.0
pair_coeff * * 1.0 1.7 1.7 3.4 3.4 1.0 1.0 1.0 :pre
[Description:]
The {gayberne} styles compute a Gay-Berne anisotropic LJ interaction
"(Berardi)"_#Berardi between pairs of ellipsoidal particles or an
ellipsoidal and spherical particle via the formulas
:c,image(Eqs/pair_gayberne.jpg)
where A1 and A2 are the transformation matrices from the simulation
box frame to the body frame and r12 is the center to center vector
between the particles. Ur controls the shifted distance dependent
interaction based on the distance of closest approach of the two
particles (h12) and the user-specified shift parameter gamma. When
both particles are spherical, the formula reduces to the usual
Lennard-Jones interaction (see details below for when Gay-Berne treats
a particle as "spherical").
Style {gayberne/gpu} is a GPU-enabled version of style {gayberne} that
should give identical answers. Depending on system size and the GPU
processor you have on your system, it may be 100x faster (for the
pairwise portion of the run time). See more details below.
For large uniform molecules it has been shown that the energy
parameters are approximately representable in terms of local contact
curvatures "(Everaers)"_#Everaers:
:c,image(Eqs/pair_gayberne2.jpg)
The variable names utilized as potential parameters are for the most
part taken from "(Everaers)"_#Everaers in order to be consistent with
its RE-squared potential fix. Details on the upsilon and mu
parameters are given "here"_gbdoc.
More details of the Gay-Berne formulation are given in the references
listed below and in "this supplementary
document"_PDF/pair_gayberne_extra.pdf.
Use of this pair style requires the NVE, NVT, or NPT fixes with the
{asphere} extension (e.g. "fix nve/asphere"_fix_nve_asphere.html) in
order to integrate particle rotation. Additionally, "atom_style
ellipsoid"_atom_style.html should be used since it defines the
rotational state of the ellipsoidal particles. The size and shape of
the ellipsoidal particles are defined by the "shape"_shape.html
command.
The following coefficients must be defined for each pair of atoms
types via the "pair_coeff"_pair_coeff.html command as in the examples
above, or in the data file or restart files read by the
"read_data"_read_data.html or "read_restart"_read_restart.html
commands, or by mixing as described below:
epsilon = well depth (energy units)
sigma = minimum effective particle radii (distance units)
epsilon_i_a = relative well depth of type I for side-to-side interactions
epsilon_i_b = relative well depth of type I for face-to-face interactions
epsilon_i_c = relative well depth of type I for end-to-end interactions
epsilon_j_a = relative well depth of type J for side-to-side interactions
epsilon_j_b = relative well depth of type J for face-to-face interactions
epsilon_j_c = relative well depth of type J for end-to-end interactions
cutoff (distance units) :ul
The last coefficient is optional. If not specified, the global
cutoff specified in the pair_style command is used.
It is typical for the Gay-Berne potential to define {sigma} as the
minimum of the 3 "shape" diameters for a I,I interaction, though this
is not required. Note that this is a different meaning for {sigma}
than the "pair_style resquared"_pair_resquared.html potential uses.
The epsilon_i and epsilon_j coefficients are actually 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.
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.
Note that if this potential is being used as a sub-style of
"pair_style hybrid"_pair_hybrid.html, and there is no "pair_coeff I I"
setting made for Gay-Berne 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.
IMPORTANT NOTE: If the epsilon a,b,c for an atom type are all 1.0, and
if the shape of the particle is spherical (see the "shape"_shape.html
command), meaning the 3 diameters are all the same, then the particle
is treated as "spherical" by the Gay-Berne potential. This is
significant because if two "spherical" particles interact, then the
simple Lennard-Jones formula is used to compute their interaction
energy/force using epsilon and sigma, which is much cheaper to compute
than the full Gay-Berne formula. Thus you should insure epsilon a,b,c
are set to 1.0 for spherical particle types and use epsilon and sigma
to specify its interaction with other spherical particles.
:line
The {gayberne/gpu} style is identical to the {gayberne} style, except
that each processor off-loads its pairwise calculations to a GPU chip.
Depending on the hardware available on your system this can provide a
significant speed-up, espcially for the relatively expensive
computations inherent in Gay-Berne interactions. See the "Running on
GPUs"_Section_start.html#2_8 section of the manual for more details
about hardware and software requirements for using GPUs.
The {gpuflag} and {gpunum} settings in the pair_style command refer to
how the GPUs on your system are configured.
Set {gpuflag} to {one/node} if you have multiple gpus on 1 CPU node of
your system. In this case, {gpunum} should be the starting GPU ID.
Set {gpuflag} to {one/node} if you have multiple gpus on 1 node
of your system. In this case, {gpunum} should be the starting GPU ID.
// set multi_gpu_mode to one/gpu to select the same gpu id on every node
// -- param is id of gpu
// set multi_gpu_mode to multi/gpu to get ma
// -- param is number of gpus per node
Additional requirements in your input script to run with style
{gayberne/gpu} are as follows:
The "newton pair"_newton.html setting must be {off}.
You should use the "neigh_modify one"_neigh_modify.html command and
set its value to something close (but slightly larger) than the number
of pairwise neighbors/atom you expect to have in your model. This is
a function of the pairwise cutoff. Note that the default for this
setting is 2000, which is much larger than most models need. Unlike
neighbor lists in LAMMPS itself, the GPU version of this pair style
uses that setting to allocate memory on the GPU for neighbor
information. If the setting is too large, it will limit the number of
atoms that can be stored on the GPU.
:line
[Mixing, shift, table, tail correction, restart, rRESPA info]:
For atom type pairs I,J and I != J, the epsilon and sigma coefficients
and cutoff distance for this pair style can be mixed. The default mix
value is {geometric}. See the "pair_modify" command for details.
This pair styles supports the "pair_modify"_pair_modify.html 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.
The "pair_modify"_pair_modify.html table option is not relevant
for this pair style.
This pair style does not support the "pair_modify"_pair_modify.html
tail option for adding long-range tail corrections to energy and
pressure.
This pair style writes its information to "binary restart
files"_restart.html, so pair_style and pair_coeff commands do not need
to be specified in an input script that reads a restart file.
This pair style can only be used via the {pair} keyword of the
"run_style respa"_run_style.html command. It does not support the
{inner}, {middle}, {outer} keywords.
:line
[Restrictions:]
The {gayberne} style is part of the "asphere" package. The
{gayberne/gpu} style is part of the "gpu" package. They are only
enabled if LAMMPS was built with the those packages. See the "Making
LAMMPS"_Section_start.html#2_3 section for more info.
This pair style requires that atoms store torque and a quaternion to
represent their orientation, as defined by the
"atom_style"_atom_style.html. It also require they store a per-type
"shape"_shape.html. The particles cannot store a per-particle
diameter.
Particles acted on by the potential can be extended aspherical or
spherical particles, or point particles.
The Gay-Berne potential does not become isotropic as r increases
"(Everaers)"_#Everaers. The distance-of-closest-approach
approximation used by LAMMPS becomes less accurate when high-aspect
ratio ellipsoids are used.
[Related commands:]
"pair_coeff"_pair_coeff.html, "fix nve/asphere"_fix_nve_asphere.html,
"compute temp/asphere"_compute_temp_asphere.html, "pair_style
resquared"_pair_resquared.html
[Default:] none
:line
:link(Everaers)
[(Everaers)] Everaers and Ejtehadi, Phys Rev E, 67, 041710 (2003).
:link(Berardi)
[(Berardi)] Berardi, Fava, Zannoni, Chem Phys Lett, 297, 8-14 (1998).
Berardi, Muccioli, Zannoni, J Chem Phys, 128, 024905 (2008).
:link(Perram)
[(Perram)] Perram and Rasmussen, Phys Rev E, 54, 6565-6572 (1996).
:link(Allen)
[(Allen)] Allen and Germano, Mol Phys 104, 3225-3235 (2006).