lammps/doc/pair_gran.txt

232 lines
10 KiB
Plaintext

"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 gran/hooke command :h3
pair_style gran/cuda command :h3
pair_style gran/hooke/history command :h3
pair_style gran/hertz/history command :h3
[Syntax:]
pair_style style Kn Kt gamma_n gamma_t xmu dampflag :pre
style = {gran/hooke} or {gran/hooke/history} or {gran/hertz/history} :ulb,l
Kn = elastic constant for normal particle repulsion (force/distance units or pressure units - see discussion below) :l
Kt = elastic constant for tangential contact (force/distance units or pressure units - see discussion below) :l
gamma_n = damping coefficient for collisions in normal direction (1/time units or 1/time-distance units - see discussion below) :l
gamma_t = damping coefficient for collisions in tangential direction (1/time units or 1/time-distance units - see discussion below) :l
xmu = static yield criterion (unitless fraction between 0.0 and 1.0) :l
dampflag = 0 or 1 if tangential damping force is excluded or included :l,ule
IMPORTANT NOTE: Versions of LAMMPS before 9Jan09 had different style
names for granular force fields. This is to emphasize the fact that
the Hertzian equation has changed to model polydispersity more
accurately. A side effect of the change is that the Kn, Kt, gamma_n,
and gamma_t coefficients in the pair_style command must be specified
with different values in order to reproduce calculations made with
earlier versions of LAMMPS, even for monodisperse systems. See the
NOTE below for details.
[Examples:]
pair_style gran/hooke/history 200000.0 NULL 50.0 NULL 0.5 1
pair_style gran/hooke 200000.0 70000.0 50.0 30.0 0.5 0 :pre
[Description:]
The {gran} styles use the following formulas for the frictional force
between two granular particles, as described in
"(Brilliantov)"_#Brilliantov, "(Silbert)"_#Silbert, and
"(Zhang)"_#Zhang, when the distance r between two particles of radii
Ri and Rj is less than their contact distance d = Ri + Rj. There is
no force between the particles when r > d.
The two Hookean styles use this formula:
:c,image(Eqs/pair_gran_hooke.jpg)
The Hertzian style uses this formula:
:c,image(Eqs/pair_gran_hertz.jpg)
In both equations the first parenthesized term is the normal force
between the two particles and the second parenthesized term is the
tangential force. The normal force has 2 terms, a contact force and a
damping force. The tangential force also has 2 terms: a shear force
and a damping force. The shear force is a "history" effect that
accounts for the tangential displacement between the particles for the
duration of the time they are in contact. This term is included in
pair styles {hooke/history} and {hertz/history}, but is not included
in pair style {hooke}. The tangential damping force term is included
in all three pair styles if {dampflag} is set to 1; it is not included
if {dampflag} is set to 0.
The other quantities in the equations are as follows:
delta = d - r = overlap distance of 2 particles
Kn = elastic constant for normal contact
Kt = elastic constant for tangential contact
gamma_n = viscoelastic damping constant for normal contact
gamma_t = viscoelastic damping constant for tangential contact
m_eff = Mi Mj / (Mi + Mj) = effective mass of 2 particles of mass Mi and Mj
Delta St = tangential displacement vector between 2 spherical particles \
which is truncated to satisfy a frictional yield criterion
n_ij = unit vector along the line connecting the centers of the 2 particles
Vn = normal component of the relative velocity of the 2 particles
Vt = tangential component of the relative velocity of the 2 particles :ul
The Kn, Kt, gamma_n, and gamma_t coefficients are specified as
parameters to the pair_style command. If a NULL is used for Kt, then
a default value is used where Kt = 2/7 Kn. If a NULL is used for
gamma_t, then a default value is used where gamma_t = 1/2 gamma_n.
The interpretation and units for these 4 coefficients are different in
the Hookean versus Hertzian equations.
The Hookean model is one where the normal push-back force for two
overlapping particles is a linear function of the overlap distance.
Thus the specified Kn is in units of (force/distance). Note that this
push-back force is independent of absolute particle size (in the
monodisperse case) and of the relative sizes of the two particles (in
the polydisperse case). This model also applies to the other terms in
the force equation so that the specified gamma_n is in units of
(1/time), Kt is in units of (force/distance), and gamma_t is in units
of (1/time).
The Hertzian model is one where the normal push-back force for two
overlapping particles is proportional to the area of overlap of the
two particles, and is thus a non-linear function of overlap distance.
Thus Kn has units of force per area and is thus specified in units of
(pressure). The effects of absolute particle size (monodispersity)
and relative size (polydispersity) are captured in the radii-dependent
pre-factors. When these pre-factors are carried through to the other
terms in the force equation it means that the specified gamma_n is in
units of (1/(time*distance)), Kt is in units of (pressure), and
gamma_t is in units of (1/(time*distance)).
Note that in the Hookean case, Kn can be thought of as a linear spring
constant with units of force/distance. In the Hertzian case, Kn is
like a non-linear spring constant with units of force/area or
pressure, and as shown in the "(Zhang)"_#Zhang paper, Kn = 4G /
(3(1-nu)) where nu = the Poisson ratio, G = shear modulus = E /
(2(1+nu)), and E = Young's modulus. Similarly, Kt = 8G / (2-nu).
Thus in the Hertzian case Kn and Kt can be set to values that
corresponds to properties of the material being modeled. This is also
true in the Hookean case, except that a spring constant must be chosen
that is appropriate for the absolute size of particles in the model.
Since relative particle sizes are not accounted for, the Hookean
styles may not be a suitable model for polydisperse systems.
IMPORTANT NOTE: In versions of LAMMPS before 9Jan09, the equation for
Hertzian interactions did not include the sqrt(RiRj/Ri+Rj) term and
thus was not as accurate for polydisperse systems. For monodisperse
systems, sqrt(RiRj/Ri+Rj) is a constant factor that effectively scales
all 4 coefficients: Kn, Kt, gamma_n, gamma_t. Thus you can set the
values of these 4 coefficients appropriately in the current code to
reproduce the results of a previous Hertzian monodisperse calculation.
For example, for the common case of a monodisperse system with
particles of diameter 1, all 4 of these coefficients should now be set
2x larger than they were previously.
Xmu is also specified in the pair_style command and is the upper limit
of the tangential force through the Coulomb criterion Ft = xmu*Fn,
where Ft and Fn are the total tangential and normal force components
in the formulas above. Thus in the Hookean case, the tangential force
between 2 particles grows according to a tangential spring and
dash-pot model until Ft/Fn = xmu and is then held at Ft = Fn*xmu until
the particles lose contact. In the Hertzian case, a similar analogy
holds, though the spring is no longer linear.
For granular styles there are no additional coefficients to set for
each pair of atom types via the "pair_coeff"_pair_coeff.html command.
All settings are global and are made via the pair_style command.
However you must still use the "pair_coeff"_pair_coeff.html for all
pairs of granular atom types. For example the command
pair_coeff * * :pre
should be used if all atoms in the simulation interact via a granular
potential (i.e. one of the pair styles above is used). If a granular
potential is used as a sub-style of "pair_style
hybrid"_pair_hybrid.html, then specific atom types can be used in the
pair_coeff command to determine which atoms interact via a granular
potential.
:line
Styles with a {cuda}, {gpu}, or {opt} 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 "this section"_Section_accelerate.html of the manual.
The accelerated styles take the same arguments and should produce the
same results, except for round-off and precision issues.
These accelerated styles are part of the "user-cuda", "gpu", and "opt"
packages respectively. They are only enabled if LAMMPS was built with
those packages. See the "Making LAMMPS"_Section_start.html#2_3
section for more info.
You can specify the accelerated styles explicitly in your input script
by including their suffix, or you can use the "-suffix command-line
switch"_Section_start.html#2_6 when you invoke LAMMPS, or you can use
the "suffix"_suffix.html command in your input script.
See "this section"_Section_accelerate.html of the manual for more
instructions on how to use the accelerated styles effectively.
:line
[Mixing, shift, table, tail correction, restart, rRESPA info]:
The "pair_modify"_pair_modify.html mix, shift, table, and tail options
are not relevant for granular pair styles.
These pair styles write their information to "binary restart
files"_restart.html, so a pair_style command does not need to be
specified in an input script that reads a restart file.
These pair styles can only be used via the {pair} keyword of the
"run_style respa"_run_style.html command. They do not support the
{inner}, {middle}, {outer} keywords.
:line
[Restrictions:] none
All the granular pair styles are part of the "granular" package. It
is only enabled if LAMMPS was built with that package. See the
"Making LAMMPS"_Section_start.html#2_3 section for more info.
These pair styles require that atoms store torque and angular velocity
(omega) as defined by the "atom_style"_atom_style.html. They also
require a per-particle radius is stored. The {sphere} atom style does
all of this.
This pair style requires you to use the "communicate vel
yes"_communicate.html option so that velocites are stored by ghost
atoms.
[Related commands:]
"pair_coeff"_pair_coeff.html
[Default:] none
:line
:link(Brilliantov)
[(Brilliantov)] Brilliantov, Spahn, Hertzsch, Poschel, Phys Rev E, 53,
p 5382-5392 (1996).
:link(Silbert)
[(Silbert)] Silbert, Ertas, Grest, Halsey, Levine, Plimpton, Phys Rev
E, 64, p 051302 (2001).
:link(Zhang)
[(Zhang)] Zhang and Makse, Phys Rev E, 72, p 011301 (2005).