lammps/doc/pair_dpd.html

195 lines
8.0 KiB
HTML
Raw Normal View History

<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 dpd command
</H3>
<H3>pair_style dpd/omp command
</H3>
<H3>pair_style dpd/tstat command
</H3>
<H3>pair_style dpd/tstat/omp command
</H3>
<P><B>Syntax:</B>
</P>
<PRE>pair_style dpd T cutoff seed
pair_style dpd/tstat Tstart Tstop cutoff seed
</PRE>
<UL><LI>T = temperature (temperature units)
<LI>Tstart,Tstop = desired temperature at start/end of run (temperature units)
<LI>cutoff = global cutoff for DPD interactions (distance units)
<LI>seed = random # seed (positive integer)
</UL>
<P><B>Examples:</B>
</P>
<PRE>pair_style dpd 1.0 2.5 34387
pair_coeff * * 3.0 1.0
pair_coeff 1 1 3.0 1.0 1.0
</PRE>
<PRE>pair_style dpd/tstat 1.0 1.0 2.5 34387
pair_coeff * * 1.0
pair_coeff 1 1 1.0 1.0
</PRE>
<P><B>Description:</B>
</P>
<P>Style <I>dpd</I> computes a force field for dissipative particle dynamics
(DPD) following the exposition in <A HREF = "#Groot">(Groot)</A>.
</P>
<P>Style <I>dpd/tstat</I> invokes a DPD thermostat on pairwise interactions,
which is equivalent to the non-conservative portion of the DPD force
field. This thermostat can be used in conjunction with any <A HREF = "pair_style.html">pair
style</A>, and in leiu of per-particle thermostats like
<A HREF = "fix_langevin.html">fix langevin</A> or ensemble thermostats like Nose
Hoover as implemented by <A HREF = "fix_nh.html">fix nvt</A>. To use <I>dpd/stat</I>
with another pair style, use the <A HREF = "pair_hybrid.html">pair_style
hybrid/overlay</A> command to compute both the desired
pair interaction and the thermostat for each pair of particles.
</P>
<P>For style <I>dpd</I>, the force on atom I due to atom J is given as a sum
of 3 terms
</P>
<CENTER><IMG SRC = "Eqs/pair_dpd.jpg">
</CENTER>
<P>where Fc is a conservative force, Fd is a dissipative force, and Fr is
a random force. Rij is a unit vector in the direction Ri - Rj, Vij is
the vector difference in velocities of the two atoms = Vi - Vj, alpha
is a Gaussian random number with zero mean and unit variance, dt is
the timestep size, and w(r) is a weighting factor that varies between
0 and 1. Rc is the cutoff. Sigma is set equal to sqrt(2 Kb T gamma),
where Kb is the Boltzmann constant and T is the temperature parameter
in the pair_style command.
</P>
<P>For style <I>dpd/tstat</I>, the force on atom I due to atom J is the same
as the above equation, except that the conservative Fc term is
dropped. Also, during the run, T is set each timestep to a ramped
value from Tstart to Tstop.
</P>
<P>For style <I>dpd</I>, the pairwise energy associated with style <I>dpd</I> is
only due to the conservative force term Fc, and is shifted to be zero
at the cutoff distance Rc. The pairwise virial is calculated using
all 3 terms. For style <I>dpd/tstat</I> there is no pairwise energy, but
the last two terms of the formula make a contribution to the virial.
</P>
<P>For style <I>dpd</I>, 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:
</P>
<UL><LI>A (force units)
<LI>gamma (force/velocity units)
<LI>cutoff (distance units)
</UL>
<P>The last coefficient is optional. If not specified, the global DPD
cutoff is used. Note that sigma is set equal to sqrt(2 T gamma),
where T is the temperature set by the <A HREF = "pair_style.html">pair_style</A>
command so it does not need to be specified.
</P>
<P>For style <I>dpd/tstat</I>, the coefficiencts defined for each pair of
atoms types via the <A HREF = "pair_coeff.html">pair_coeff</A> command is the same,
except that A is not included.
</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>These pair styles do not support mixing. Thus, coefficients for all
I,J pairs must be specified explicitly.
</P>
<P>These pair styles do not support the <A HREF = "pair_modify.html">pair_modify</A>
shift option for the energy of the pair interaction. Note that as
discussed above, the energy due to the conservative Fc term is already
shifted to be 0.0 at the cutoff distance Rc.
</P>
<P>The <A HREF = "pair_modify.html">pair_modify</A> table option is not relevant
for these pair styles.
</P>
<P>These pair style do not support the <A HREF = "pair_modify.html">pair_modify</A>
tail option for adding long-range tail corrections to energy and
pressure.
</P>
<P>These pair styles writes 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. Note
that the user-specified random number seed is stored in the restart
file, so when a simulation is restarted, each processor will
re-initialize its random number generator the same way it did
initially. This means the random forces will be random, but will not
be the same as they would have been if the original simulation had
continued past the restart time.
</P>
<P>These pair styles can only be used via the <I>pair</I> keyword of the
<A HREF = "run_style.html">run_style respa</A> command. They do not support the
<I>inner</I>, <I>middle</I>, <I>outer</I> keywords.
</P>
<P>The <I>dpd/tstat</I> style can ramp its target temperature over multiple
runs, using the <I>start</I> and <I>stop</I> keywords of the <A HREF = "run.html">run</A>
command. See the <A HREF = "run.html">run</A> command for details of how to do
this.
</P>
<HR>
<P><B>Restrictions:</B>
</P>
<P>The default frequency for rebuilding neighbor lists is every 10 steps
(see the <A HREF = "neigh_modify.html">neigh_modify</A> command). This may be too
infrequent for style <I>dpd</I> simulations since particles move rapidly
and can overlap by large amounts. If this setting yields a non-zero
number of "dangerous" reneighborings (printed at the end of a
simulation), you should experiment with forcing reneighboring more
often and see if system energies/trajectories change.
</P>
<P>These pair styles requires you to use the <A HREF = "communicate.html">communicate vel
yes</A> option so that velocites are stored by ghost
atoms.
</P>
<P>These pair styles will not restart exactly when using the
<A HREF = "read_restart.html">read_restart</A> command, though they should provide
statistically similar results. This is because the forces they
compute depend on atom velocities. See the
<A HREF = "read_restart.html">read_restart</A> command for more details.
</P>
<P><B>Related commands:</B>
</P>
<P><A HREF = "pair_coeff.html">pair_coeff</A>, <A HREF = "fix_nh.html">fix nvt</A>, <A HREF = "fix_langevin.html">fix
langevin</A>
</P>
<P><B>Default:</B> none
</P>
<HR>
<A NAME = "Groot"></A>
<P><B>(Groot)</B> Groot and Warren, J Chem Phys, 107, 4423-35 (1997).
</P>
</HTML>