forked from lijiext/lammps
197 lines
8.6 KiB
Plaintext
197 lines
8.6 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
|
|
|
|
fix qeq/point command :h3
|
|
fix qeq/shielded command :h3
|
|
fix qeq/slater command :h3
|
|
fix qeq/dynamic command :h3
|
|
|
|
[Syntax:]
|
|
|
|
fix ID group-ID style Nevery cutoff tolerance maxiter qfile :pre
|
|
|
|
ID, group-ID are documented in "fix"_fix.html command
|
|
style = {qeq/point} or {qeq/shielded} or {qeq/slater} or {qeq/dynamic}
|
|
Nevery = perform charge equilibration every this many steps
|
|
cutoff = global cutoff for charge-charge interactions (distance unit)
|
|
tolerance = precision to which charges will be equilibrated
|
|
maxiter = maximum iterations to perform charge equilibration
|
|
qfile = a filename with QEq parameters :ul
|
|
|
|
[Examples:]
|
|
|
|
fix 1 all qeq/point 1 10 1.0e-6 200 param.qeq1
|
|
fix 1 qeq qeq/shielded 1 8 1.0e-6 100 param.qeq2
|
|
fix 1 all qeq/slater 5 10 1.0e-6 100 params
|
|
fix 1 qeq qeq/dynamic 1 12 1.0e-3 100 my_qeq :pre
|
|
|
|
[Description:]
|
|
|
|
Perform the charge equilibration (QEq) method as described in "(Rappe
|
|
and Goddard)"_#Rappe and formulated in "(Nakano)"_#Nakano (also known
|
|
as the matrix inversion method) and in "(Rick and Stuart)"_#Rick (also
|
|
known as the extended Lagrangian method) based on the
|
|
electronegativity equilization principle.
|
|
|
|
These fixes can be used with any "pair style"_pair_style.html in
|
|
LAMMPS, so long as per-atom charges are defined. The most typical
|
|
use-case is in conjunction with a "pair style"_pair_style.html that
|
|
performs charge equilibration periodically (e.g. every timestep), such
|
|
as the ReaxFF or Streitz-Mintmire potential (the latter is not yet
|
|
implemented in LAMMPS). But these fixes can also be used with
|
|
potentials that normally assume per-atom charges are fixed, e.g. a
|
|
"Buckingham"_pair_buck.html or "LJ/Coulombic"_pair_lj.html potential.
|
|
|
|
Because the charge equilibration calculation is effectively
|
|
independent of the pair style, these fixes can also be used to perform
|
|
a one-time assignment of charges to atoms. For example, you could
|
|
define the QEq fix, perform a zero-timestep run via the "run"_run.html
|
|
command without any pair style defined which would set per-atom
|
|
charges (based on the current atom configuration), then remove the fix
|
|
via the "unfix"_unfix.html command before performing further dynamics.
|
|
|
|
IMPORTANT NOTE: Computing and using charge values different from
|
|
published values defined for a fixed-charge potential like Buckingham
|
|
or CHARMM or AMBER, can have a strong effect on energies and forces,
|
|
and produces a different model than the published versions.
|
|
|
|
IMPORTANT NOTE: The "fix qeq/comb"_fix_qeq_comb.html command must
|
|
still be used to perform charge equliibration with the "COMB
|
|
potential"_pair_comb.html. The "fix qeq/reax"_fix_qeq_reax.html
|
|
command can be used to perform charge equilibration with the "ReaxFF
|
|
force field"_pair_reax_c.html, although fix qeq/shielded yields the
|
|
same results as fix qeq/reax if {Nevery}, {cutoff}, and {tolerance} are
|
|
the same. Eventually the fix qeq/reax command will be deprecated.
|
|
|
|
The QEq method minimizes the electrostatic energy of the system (or
|
|
equalizes the derivative of energy with respect to charge of all the
|
|
atoms) by adjusting the partial charge on individual atoms based on
|
|
interactions with their neighbors within {cutoff}. It reqires a few
|
|
parameters, in {metal} units, for each atom type which provided in a
|
|
file specified by {qfile}. The file has the following format
|
|
|
|
1 chi eta gamma zeta qcore
|
|
2 chi eta gamma zeta qcore
|
|
...
|
|
Ntype chi eta gamma zeta qcore :pre
|
|
|
|
There is one line per atom type with the following parameters.
|
|
Only a subset of the parameters is used by each QEq style as descibed
|
|
below, thus the others can be set to 0.0 if desired.
|
|
|
|
{chi} = electronegativity in energy units
|
|
{eta} = self-Coulomb potential in energy units
|
|
{gamma} = shielded Coulomb constant defined by "ReaxFF force field"_#vanDuin in distance units
|
|
{zeta} = Slater type orbital exponent defined by the "Streitz-Mintmire"_#Streitz potential in reverse distance units
|
|
{qcore} = charge of the nucleus defined by the "Streitz-Mintmire potential"_#Streitz potential in charge units :ul
|
|
|
|
The {qeq/point} style describes partial charges on atoms as point
|
|
charges. Interaction between a pair of charged particles is 1/r,
|
|
which is the simplest description of the interaction between charges.
|
|
Only the {chi} and {eta} parameters from the {qfile} file are used.
|
|
Note that Coulomb catastrophe can occur if repulsion between the pair
|
|
of charged particles is too weak. This style solves partial charges
|
|
on atoms via the matrix inversion method. A tolerance of 1.0e-6 is
|
|
usually a good number.
|
|
|
|
The {qeq/shielded} style describes partial charges on atoms also as
|
|
point charges, but uses a shielded Coulomb potential to describe the
|
|
interaction between a pair of charged particles. Interaction through
|
|
the shielded Coulomb is given by equation (13) of the "ReaxFF force
|
|
field"_#vanDuin paper. The shielding accounts for charge overlap
|
|
between charged particles at small separation. This style is the same
|
|
as "fix qeq/reax"_fix_qeq_reax.html, and can be used with "pair_style
|
|
reax/c"_pair_reax_c.html. Only the {chi}, {eta}, and {gamma}
|
|
parameters from the {qfile} file are used. This style solves partial
|
|
charges on atoms via the matrix inversion method. A tolerance of
|
|
1.0e-6 is usually a good number.
|
|
|
|
The {qeq/slater} style describes partial charges on atoms as spherical
|
|
charge densities centered around atoms via the Slater 1{s} orbital, so
|
|
that the interaction between a pair of charged particles is the
|
|
product of two Slater 1{s} orbitals. The expression for the Slater
|
|
1{s} orbital is given under equation (6) of the
|
|
"Streitz-Mintmire"_#Streitz paper. Only the {chi}, {eta}, {zeta}, and
|
|
{qcore} parameters from the {qfile} file are used. This style solves
|
|
partial charges on atoms via the matrix inversion method. A tolerance
|
|
of 1.0e-6 is usually a good number.
|
|
|
|
The {qeq/dynamic} style describes partial charges on atoms as point
|
|
charges that interact through 1/r, but the extended Lagrangian method
|
|
is used to solve partial charges on atoms. Only the {chi} and {eta}
|
|
parameters from the {qfile} file are used. Note that Coulomb
|
|
catastrophe can occur if repulsion between the pair of charged
|
|
particles is too weak. A tolerance of 1.0e-3 is usually a good
|
|
number.
|
|
|
|
Note that {qeq/point}, {qeq/shielded}, and {qeq/slater} describe
|
|
different charge models, whereas the matrix inversion method and the
|
|
extended Lagrangian method ({qeq/dynamic}) are different solvers.
|
|
|
|
Note that the {qeq/point} and the {qeq/dynamic} styles both describe
|
|
charges as point charges that interact through 1/r relationship, but
|
|
solve partial charges on atoms using different solvers. Styles
|
|
{qeq/point} and the {qeq/dynamic} should yield comparable results if
|
|
the QEq parameters and {Nevery}, {cutoff}, and {tolerance} are the
|
|
same. Style {qeq/point} is typically faster, but {qeq/dynamic} scales
|
|
better on larger sizes.
|
|
|
|
IMPORTANT NOTE: To avoid the evaluation of the derivative of charge
|
|
with respect to position, which is typically ill-defined, the system
|
|
should have a zero net charge.
|
|
|
|
IMPORTANT NOTE: Developing QEq parameters (chi, eta, gamma, zeta, and
|
|
qcore) is an "art". Charges on atoms are not guaranteed to
|
|
equilibrate with arbitrary choices of these parameters. We do not
|
|
develop these QEq paramters. See the examples/qeq directory for some
|
|
examples.
|
|
|
|
[Restart, fix_modify, output, run start/stop, minimize info:]
|
|
|
|
No information about these fixes is written to "binary restart
|
|
files"_restart.html. No global scalar or vector or per-atom
|
|
quantities are stored by these fixes for access by various "output
|
|
commands"_Section_howto.html#howto_15. No parameter of these fixes
|
|
can be used with the {start/stop} keywords of the "run"_run.html
|
|
command.
|
|
|
|
Thexe fixes are invoked during "energy minimization"_minimize.html.
|
|
|
|
[Restrictions:]
|
|
|
|
These fixes are part of the QEQ package. They are only enabled if
|
|
LAMMPS was built with that package. See the "Making
|
|
LAMMPS"_Section_start.html#start_3 section for more info.
|
|
|
|
[Related commands:]
|
|
|
|
"fix qeq/reax"_fix_qeq_reax.html, "fix qeq/comb"_fix_qeq_comb.html
|
|
|
|
[Default:] none
|
|
|
|
:line
|
|
|
|
:link(Rappe)
|
|
[(Rappe and Goddard)] A. K. Rappe and W. A. Goddard III, J Physical
|
|
Chemistry, 95, 3358-3363 (1991).
|
|
|
|
:link(Nakano)
|
|
[(Nakano)] A. Nakano, Computer Physics Communications, 104, 59-69 (1997).
|
|
|
|
:link(Rick)
|
|
[(Rick and Stuart)] S. W. Rick, S. J. Stuart, B. J. Berne, J Chemical Physics
|
|
101, 16141 (1994).
|
|
|
|
:link(Streitz)
|
|
[(Streitz-Mintmire)] F. H. Streitz, J. W. Mintmire, Physical Review B, 50,
|
|
16, 11996 (1994)
|
|
|
|
:link(vanDuin)
|
|
[(ReaxFF)] A. C. T. van Duin, S. Dasgupta, F. Lorant, W. A. Goddard III, J
|
|
Physical Chemistry, 105, 9396-9049 (2001)
|