2006-09-22 00:22:34 +08:00
|
|
|
<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>atom_modify command
|
|
|
|
</H3>
|
|
|
|
<P><B>Syntax:</B>
|
|
|
|
</P>
|
2010-01-14 04:56:56 +08:00
|
|
|
<PRE>atom_modify keyword values ...
|
2006-09-22 00:22:34 +08:00
|
|
|
</PRE>
|
|
|
|
<UL><LI>one or more keyword/value pairs may be appended
|
|
|
|
|
2010-01-14 04:56:56 +08:00
|
|
|
<LI>keyword = <I>map</I> or <I>first</I> or <I>sort</I>
|
2006-09-22 00:22:34 +08:00
|
|
|
|
2008-02-16 03:09:18 +08:00
|
|
|
<PRE> <I>map</I> value = <I>array</I> or <I>hash</I>
|
2010-01-14 04:56:56 +08:00
|
|
|
<I>first</I> value = group-ID = group whose atoms will appear first in internal atom lists
|
|
|
|
<I>sort</I> values = Nfreq binsize
|
|
|
|
Nfreq = sort atoms spatially every this many time steps
|
|
|
|
binsize = bin size for spatial sorting (distance units)
|
2006-09-22 00:22:34 +08:00
|
|
|
</PRE>
|
|
|
|
|
|
|
|
</UL>
|
|
|
|
<P><B>Examples:</B>
|
|
|
|
</P>
|
2008-02-16 03:09:18 +08:00
|
|
|
<PRE>atom_modify map hash
|
2010-01-14 04:56:56 +08:00
|
|
|
atom_modify map array sort 10000 2.0
|
2008-02-16 03:09:18 +08:00
|
|
|
atom_modify first colloid
|
2006-09-22 00:22:34 +08:00
|
|
|
</PRE>
|
|
|
|
<P><B>Description:</B>
|
|
|
|
</P>
|
|
|
|
<P>Modify properties of the atom style selected within LAMMPS.
|
|
|
|
</P>
|
|
|
|
<P>The <I>map</I> keyword determines how atom ID lookup is done for molecular
|
|
|
|
problems. Lookups are performed by bond (angle, etc) routines in
|
|
|
|
LAMMPS to find the local atom index associated with a global atom ID.
|
|
|
|
When the <I>array</I> value is used, each processor stores a lookup table
|
|
|
|
of length N, where N is the total # of atoms in the system. This is
|
|
|
|
the fastest method for most simulations, but a processor can run out
|
|
|
|
of memory to store the table for very large simulations. The <I>hash</I>
|
|
|
|
value uses a hash table to perform the lookups. This method can be
|
|
|
|
slightly slower than the <I>array</I> method, but its memory cost is
|
|
|
|
proportional to N/P on each processor, where P is the total number of
|
|
|
|
processors running the simulation.
|
|
|
|
</P>
|
2008-02-16 03:09:18 +08:00
|
|
|
<P>The <I>first</I> keyword allows a <A HREF = "group.html">group</A> to be specified whose
|
|
|
|
atoms will be maintained as the first atoms in each processor's list
|
|
|
|
of owned atoms. This in only useful when the specified group is a
|
|
|
|
small fraction of all the atoms, and there are other operations LAMMPS
|
|
|
|
is performing that will be sped-up significantly by being able to loop
|
|
|
|
over the smaller set of atoms. Otherwise the reordering required by
|
|
|
|
this option will be a net slow-down. The <A HREF = "neigh_modify.html">neigh_modify
|
|
|
|
include</A> and <A HREF = "communicate.html">communicate group</A>
|
|
|
|
commands are two examples of commands that require this setting to
|
|
|
|
work efficiently. Several <A HREF = "fix.html">fixes</A>, most notably time
|
|
|
|
integration fixes like <A HREF = "fix_nve.html">fix nve</A>, also take advantage of
|
|
|
|
this setting if the group they operate on is the group specified by
|
|
|
|
this command.
|
|
|
|
</P>
|
2010-01-14 04:56:56 +08:00
|
|
|
<P>It is OK to use the <I>first</I> keyword with a group that has not yet been
|
|
|
|
defined, e.g. to use the atom_modify command at the beginning of your
|
|
|
|
input script. LAMMPS will check that the group exists before the
|
|
|
|
first simulation is performed.
|
|
|
|
</P>
|
|
|
|
<P>The <I>sort</I> keyword turns on a spatial sorting or reordering of atoms
|
|
|
|
within each processor's sub-domain every <I>Nfreq</I> timesteps. This can
|
|
|
|
improve cache performance and thus speed=up a LAMMPS simulation, as
|
|
|
|
discussed in a paper by <A HREF = "#Meloni">(Meloni)</A>. In tests we have done,
|
|
|
|
the amount of speed-up can range from zero to 3-4x. It is typically
|
|
|
|
more effective at speeding up simulations of liquids as opposed to
|
|
|
|
solids.
|
|
|
|
</P>
|
|
|
|
<P>Reordering is peformed every <I>Nfreq</I> timesteps during a dynamics run
|
|
|
|
or iterations during a minimization. More precisely, reordering
|
|
|
|
occurs at the first reneighboring that occurs after the target
|
|
|
|
timestep. The reordering is performed locally by each processor,
|
|
|
|
using bins of the specified <I>binsize</I>. If <I>binsize</I> is set to 0.0,
|
|
|
|
then a binsize equal to half the <A HREF = "neighbor.html">neighbor</A> cutoff
|
|
|
|
distance (force cutoff plus skin distance) is used, which is a
|
|
|
|
reasonable value. After the atoms have been binned, they are
|
|
|
|
reordered so that atoms in the same bin are adjacent to each other in
|
|
|
|
the processor's 1d list of atoms.
|
|
|
|
</P>
|
|
|
|
<P>The goal of this procedure is for atoms be near each other in the
|
|
|
|
processor's 1d list of atoms that are also near to each other
|
|
|
|
spatially. This can improve cache performance when pairwise
|
|
|
|
intereractions and neighbor lists are computed. Note that if bins are
|
|
|
|
too small, there will be few atoms/bin. Likewise if bins are too
|
|
|
|
large, there will be many atoms/bin. In both cases, the goal of cache
|
|
|
|
locality can be undermined.
|
|
|
|
</P>
|
|
|
|
<P>IMPORTANT NOTE: Running a simulation with sorting on versus off should
|
|
|
|
not change the simulation results in a statistical sense. However,
|
|
|
|
reordering will induce round-off differences, which will lead to
|
|
|
|
diverging trajectories when comparing two simluations. Various
|
|
|
|
commands, particularly those which use random numbers (e.g. <A HREF = "velocity.html">velocity
|
|
|
|
create</A>, and <A HREF = "fix_langevin.html">fix langevin</A>), may
|
|
|
|
generate different results (but statistically identical) depending on
|
|
|
|
the order in which they process atoms. The order of atoms in a
|
|
|
|
<A HREF = "dump.html">dump</A> file will also change if sorting is enabled.
|
2008-02-16 03:09:18 +08:00
|
|
|
</P>
|
2006-09-22 00:22:34 +08:00
|
|
|
<P><B>Restrictions:</B>
|
|
|
|
</P>
|
2010-01-14 04:56:56 +08:00
|
|
|
<P>The map keyword can only be used before the simulation box is defined
|
|
|
|
by a <A HREF = "read_data.html">read_data</A> or <A HREF = "create_box.html">create_box</A>
|
|
|
|
command.
|
|
|
|
</P>
|
2010-01-14 04:57:58 +08:00
|
|
|
<P>The <I>first</I> and <I>sort</I> options cannot be used together. Since sorting
|
|
|
|
is on by default, it will be turned off if the <I>first</I> keyword is
|
|
|
|
used.
|
2006-09-22 00:22:34 +08:00
|
|
|
</P>
|
|
|
|
<P><B>Related commands:</B> none
|
|
|
|
</P>
|
|
|
|
<P><B>Default:</B>
|
|
|
|
</P>
|
2008-02-16 03:09:18 +08:00
|
|
|
<P>By default, atomic (non-molecular) problems do not allocate maps. For
|
2009-03-19 01:53:44 +08:00
|
|
|
molecular problems, the option default is map = array. By default, a
|
2010-01-14 04:56:56 +08:00
|
|
|
"first" group is not defined. By default, sorting is enabled with a
|
|
|
|
frequency of 1000 and a binsize of 0.0, which means the neighbor
|
|
|
|
cutoff will be used to set the bin size.
|
|
|
|
</P>
|
|
|
|
<HR>
|
|
|
|
|
|
|
|
<A NAME = "Meloni"></A>
|
|
|
|
|
|
|
|
<P><B>(Meloni)</B> Meloni and Rasati, J Chem Phys, 126, 121102 (2007).
|
2006-09-22 00:22:34 +08:00
|
|
|
</P>
|
|
|
|
</HTML>
|