lammps/doc/kspace_style.html

190 lines
7.8 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>kspace_style command
</H3>
<P><B>Syntax:</B>
</P>
<PRE>kspace_style style value
</PRE>
<UL><LI>style = <I>none</I> or <I>ewald</I> or <I>pppm</I> or <I>pppm/cg</I> or <I>pppm/tip4p</I> or <I>ewald/n</I> or <I>pppm/gpu</I>
<PRE> <I>none</I> value = none
<I>ewald</I> value = precision
precision = desired accuracy
<I>pppm</I> value = precision
precision = desired accuracy
<I>pppm/cg</I> value = precision (smallq)
precision = desired accuracy
smallq = cutoff for charges to be considered (optional) (charge units)
<I>pppm/tip4p</I> value = precision
precision = desired accuracy
<I>ewald/n</I> value = precision
precision = desired accuracy
<I>pppm/gpu</I> value = precision
precision = desired accuracy
</PRE>
</UL>
<P><B>Examples:</B>
</P>
<PRE>kspace_style pppm 1.0e-4
kspace_style pppm/cg 1.0e-5 1.0e-6
kspace_style none
</PRE>
<P><B>Description:</B>
</P>
<P>Define a K-space solver for LAMMPS to use each timestep to compute
long-range Coulombic interactions or long-range 1/r^N interactions.
When such a solver is used in conjunction with an appropriate pair
style, the cutoff for Coulombic or other 1/r^N interactions is
effectively infinite; each charge in the system interacts with charges
in an infinite array of periodic images of the simulation domain.
</P>
<P>The <I>ewald</I> style performs a standard Ewald summation as described in
any solid-state physics text.
</P>
<P>The <I>pppm</I> style invokes a particle-particle particle-mesh solver
<A HREF = "#Hockney">(Hockney)</A> which maps atom charge to a 3d mesh, uses 3d FFTs
to solve Poisson's equation on the mesh, then interpolates electric
fields on the mesh points back to the atoms. It is closely related to
the particle-mesh Ewald technique (PME) <A HREF = "#Darden">(Darden)</A> used in
AMBER and CHARMM. The cost of traditional Ewald summation scales as
N^(3/2) where N is the number of atoms in the system. The PPPM solver
scales as Nlog(N) due to the FFTs, so it is almost always a faster
choice <A HREF = "#Pollock">(Pollock)</A>.
</P>
<P>The <I>pppm/cg</I> style is identical to the <I>pppm</I> style except that it
has an optimization for systems where most particles are uncharged.
The optional <I>smallq</I> argument defines the cutoff for the absolute
charge value which determines whether a particle is considered charged
or not. Its default value is 1.0e-5.
</P>
<P>The <I>pppm/tip4p</I> style is identical to the <I>pppm</I> style except that it
adds a charge at the massless 4th site in each TIP4P water molecule.
It should be used with <A HREF = "pair_style.html">pair styles</A> with a
<I>long/tip4p</I> in their style name.
</P>
<P>The <I>ewald/n</I> style augments <I>ewald</I> by adding long-range dispersion
sum capabilities for 1/r^N potentials and is useful for simulation of
interfaces <A HREF = "#Veld">(Veld)</A>. It also performs standard coulombic Ewald
summations, but in a more efficient manner than the <I>ewald</I> style.
The 1/r^N capability means that Lennard-Jones or Buckingham potentials
can be used with <I>ewald/n</I> without a cutoff, i.e. they become full
long-range potentials.
</P>
<P>Currently, only the <I>ewald/n</I> style can be used with non-orthogonal
(triclinic symmetry) simulation boxes.
</P>
<P>Note that the PPPM styles can be used with single-precision FFTs by
using the compiler switch -DFFT_SINGLE for the FFT_INC setting in your
lo-level Makefile. This setting also changes some of the PPPM
operations (e.g. mapping charge to mesh and interpolating electric
fields to particles) to be performed in single precision. This option
can speed-up long-range calulations, particularly in parallel or on
GPUs. The use of the -DFFT_SINGLE flag is discussed in <A HREF = "Section_start.html#start_2_4">this
section</A> of the manual.
</P>
<HR>
<P>When a kspace style is used, a pair style that includes the
short-range correction to the pairwise Coulombic or other 1/r^N forces
must also be selected. For Coulombic interactions, these styles are
ones that have a <I>coul/long</I> in their style name. For 1/r^6
dispersion forces in a Lennard-Jones or Buckingham potential, see the
<A HREF = "pair_lj_coul.html">pair_style lj/coul</A> or <A HREF = "pair_buck_coul.html">pair_style
buck/coul</A> commands.
</P>
<P>A precision value of 1.0e-4 means one part in 10000. This setting is
used in conjunction with the pairwise cutoff to determine the number
of K-space vectors for style <I>ewald</I> or the FFT grid size for style
<I>pppm</I>.
</P>
<P>See the <A HREF = "kspace_modify.html">kspace_modify</A> command for additional
options of the K-space solvers that can be set.
</P>
<HR>
<P>Styles with a <I>cuda</I>, <I>gpu</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">this section</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>More specifically, the <I>pppm/gpu</I> style performs charge assignment and
force interpolation calculations on the GPU. These processes are
performed either in single or double precision, depending on whether
the -DFFT_SINGLE setting was specified in your lo-level Makefile, as
discussed above. The FFTs themselves are still calculated on the CPU.
If <I>pppm/gpu</I> is used with a GPU-enabled pair style, part of the PPPM
calculation can be performed concurrently on the GPU while other
calculations for non-bonded and bonded force calculation are performed
on the CPU.
</P>
<P>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 <A HREF = "Section_start.html#start_3">Making LAMMPS</A>
section for more info.
</P>
<P>See <A HREF = "Section_accelerate.html">this section</A> of the manual for more
instructions on how to use the accelerated styles effectively.
</P>
<P><B>Restrictions:</B>
</P>
<P>A simulation must be 3d and periodic in all dimensions to use an Ewald
or PPPM solver. The only exception is if the slab option is set with
<A HREF = "kspace_modify.html">kspace_modify</A>, in which case the xy dimensions
must be periodic and the z dimension must be non-periodic.
</P>
<P>Kspace styles are part of the KSPACE package. They are 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>The <I>ewald/n</I> style is part of the USER-EWALDN 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>When using a long-range pairwise TIP4P potential, you must use kspace
style <I>pppm/tip4p</I> and vice versa.
</P>
<P><B>Related commands:</B>
</P>
<P><A HREF = "kspace_modify.html">kspace_modify</A>, <A HREF = "pair_lj.html">pair_style
lj/cut/coul/long</A>, <A HREF = "pair_charmm.html">pair_style
lj/charmm/coul/long</A>, <A HREF = "pair_lj_coul.html">pair_style
lj/coul</A>, <A HREF = "pair_buck.html">pair_style buck/coul/long</A>
</P>
<P><B>Default:</B>
</P>
<PRE>kspace_style none
</PRE>
<HR>
<A NAME = "Darden"></A>
<P><B>(Darden)</B> Darden, York, Pedersen, J Chem Phys, 98, 10089 (1993).
</P>
<A NAME = "Hockney"></A>
<P><B>(Hockney)</B> Hockney and Eastwood, Computer Simulation Using Particles,
Adam Hilger, NY (1989).
</P>
<A NAME = "Pollock"></A>
<P><B>(Pollock)</B> Pollock and Glosli, Comp Phys Comm, 95, 93 (1996).
</P>
<A NAME = "Veld"></A>
<P><B>(Veld)</B> In 't Veld, Ismail, Grest, J Chem Phys, in press (2007).
</P>
</HTML>