<PRE>fix ID group-ID <I>shake/rattle</I> tol iter N constraint values ... keyword value ...
</PRE>
<UL><LI>ID, group-ID are documented in <A HREF = "fix.html">fix</A> command
<LI>shake or rattle = style name of this fix command
<LI>tol = accuracy tolerance of SHAKE solution
<LI>iter = max # of iterations in each SHAKE solution
<LI>N = print SHAKE statistics every this many timesteps (0 = never)
<LI>one or more constraint/value pairs are appended
<LI>constraint = <I>b</I> or <I>a</I> or <I>t</I> or <I>m</I>
<PRE> <I>b</I> values = one or more bond types
<I>a</I> values = one or more angle types
<I>t</I> values = one or more atom types
<I>m</I> value = one or more mass values
</PRE>
<LI>zero or more keyword/value pairs may be appended
<LI>keyword = <I>mol</I>
<PRE> <I>mol</I> value = template-ID
template-ID = ID of molecule template specified in a separate <A HREF = "molecule.html">molecule</A> command
</PRE>
</UL>
<P><B>Examples:</B>
</P>
<PRE>fix 1 sub shake 0.0001 20 10 b 4 19 a 3 5 2
fix 1 sub shake 0.0001 20 10 t 5 6 m 1.0 a 31
fix 1 sub shake 0.0001 20 10 t 5 6 m 1.0 a 31 mol myMol
</PRE>
<PRE>fix 1 sub rattle 0.0001 20 10 t 5 6 m 1.0 a 31
fix 1 sub rattle 0.0001 20 10 t 5 6 m 1.0 a 31 mol myMol
</PRE>
<P><B>Description:</B>
</P>
<P>Apply bond and angle constraints to specified bonds and angles in the
simulation. This typically enables a longer timestep.
</P>
<P><B>SHAKE vs. RATTLE:</B>
</P>
<P>The SHAKE algorithm was invented for schemes such as Verlet, where only the coordinates
are integrated and the velocities are approximated as finite differences to the trajectories (<A HREF = "#Ryckaert">Ryckaert et al. (1977)</A>).
If the velocities are integrated explicitly, for example with velocity Verlet,
a second set of constraining forces is required in order to eliminate velocity components along the bonds (<A HREF = "#Andersen">Andersen (1983)</A>).
</P>
<P>In order to formulate the individual constraints for SHAKE and RATTLE, let us focus on a single molecule whose
bonds are to be constrained. Let image(Eqs/fix_rattle_ri.jpg) and image(Eqs/fix_rattle_vi.jpg) be the position and velocity of atom <I>i</I> at time <I>n</I>, for <I>i</I>=1,...,<I>N</I>, where <I>N</I> is the number of sites of our reference molecule. The distance vector between sites <I>i</I> and <I>j</I> is given by image(Eqs/fix_rattle_rij.jpg). The constraints can then be formulated as
<P>The SHAKE algorithm satisfies the first condition, i.e. the sites at time <I>n+1</I> will have the desired separations image(Eqs/fix_rattle_dij.jpg) immediately after the coordinates are integrated.
If we also enforce the second condition, the velocity components along the bonds will vanish.
RATTLE satisfies both conditions and fix rattle is implemented in a way that it uses fix shake for dealing with the coordinate constraints. Therefore all the information about SHAKE below is also relevant for RATTLE.
</P>
<P><B>SHAKE:</B>
</P>
<P>Each timestep the specified bonds and angles are reset to their
equilibrium lengths and angular values via the well-known SHAKE
algorithm. This is done by applying an additional constraint force so
that the new positions preserve the desired atom separations. The
equations for the additional force are solved via an iterative method
that typically converges to an accurate solution in a few iterations.
The desired tolerance (e.g. 1.0e-4 = 1 part in 10000) and maximum # of
iterations are specified as arguments. Setting the N argument will
print statistics to the screen and log file about regarding the
lengths of bonds and angles that are being constrained. Small delta
values mean SHAKE is doing a good job.
</P>
<P>In LAMMPS, only small clusters of atoms can be constrained. This is
so the constraint calculation for a cluster can be performed by a
single processor, to enable good parallel performance. A cluster is
defined as a central atom connected to others in the cluster by
constrained bonds. LAMMPS allows for the following kinds of clusters
to be constrained: one central atom bonded to 1 or 2 or 3 atoms, or
one central atom bonded to 2 others and the angle between the 3 atoms
also constrained. This means water molecules or CH2 or CH3 groups may
be constrained, but not all the C-C backbone bonds of a long polymer
chain.
</P>
<P>The <I>b</I> constraint lists bond types that will be constrained. The <I>t</I>
constraint lists atom types. All bonds connected to an atom of the
specified type will be constrained. The <I>m</I> constraint lists atom
masses. All bonds connected to atoms of the specified masses will be
constrained (within a fudge factor of MASSDELTA specified in
fix_shake.cpp). The <I>a</I> constraint lists angle types. If both bonds
in the angle are constrained then the angle will also be constrained
if its type is in the list.
</P>
<P>For all constraints, a particular bond is only constrained if both
atoms in the bond are in the group specified with the SHAKE fix.
</P>
<P>The degrees-of-freedom removed by SHAKE bonds and angles are accounted
for in temperature and pressure computations. Similarly, the SHAKE
contribution to the pressure of the system (virial) is also accounted
for.
</P>
<P>IMPORTANT NOTE: This command works by using the current forces on
atoms to caculate an additional constraint force which when added will
leave the atoms in positions that satisfy the SHAKE constraints
(e.g. bond length) after the next time integration step. If you