lammps/doc/doc2/fix_pour.html

243 lines
12 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>fix pour command
</H3>
<P><B>Syntax:</B>
</P>
<PRE>fix ID group-ID pour N type seed keyword values ...
</PRE>
<UL><LI>ID, group-ID are documented in <A HREF = "fix.html">fix</A> command
<LI>pour = style name of this fix command
<LI>N = # of atoms to insert
<LI>type = atom type to assign to inserted atoms (offset for molecule insertion)
<LI>seed = random # seed (positive integer)
<LI>one or more keyword/value pairs may be appended to args
<LI>keyword = <I>region</I> or <I>diam</I> or <I>dens</I> or <I>vol</I> or <I>rate</I> or <I>vel</I> or <I>mol</I> or <I>rigid</I> or <I>shake</I>
<PRE> <I>region</I> value = region-ID
region-ID = ID of region to use as insertion volume
<I>diam</I> values = dstyle args
dstyle = <I>one</I> or <I>range</I> or <I>poly</I>
<I>one</I> args = D
D = single diameter for inserted particles (distance units)
<I>range</I> args = Dlo Dhi
Dlo,Dhi = range of diameters for inserted particles (distance units)
<I>poly</I> args = Npoly D1 P1 D2 P2 ...
Npoly = # of (D,P) pairs
D1,D2,... = diameter for subset of inserted particles (distance units)
P1,P2,... = percentage of inserted particles with this diameter (0-1)
<I>vol</I> values = fraction Nattempt
fraction = desired volume fraction for filling insertion volume
Nattempt = max # of insertion attempts per atom
<I>rate</I> value = V
V = z velocity (3d) or y velocity (2d) at which
insertion volume moves (velocity units)
<I>vel</I> values (3d) = vxlo vxhi vylo vyhi vz
<I>vel</I> values (2d) = vxlo vxhi vy
vxlo,vxhi = range of x velocities for inserted particles (velocity units)
vylo,vyhi = range of y velocities for inserted particles (velocity units)
vz = z velocity (3d) assigned to inserted particles (velocity units)
vy = y velocity (2d) assigned to inserted particles (velocity units)
<I>mol</I> value = template-ID
template-ID = ID of molecule template specified in a separate <A HREF = "molecule.html">molecule</A> command
<I>molfrac</I> values = f1 f2 ... fN
f1 to fN = relative probability of creating each of N molecules in template-ID
<I>rigid</I> value = fix-ID
fix-ID = ID of <A HREF = "fix_rigid.html">fix rigid/small</A> command
<I>shake</I> value = fix-ID
fix-ID = ID of <A HREF = "fix_shake.html">fix shake</A> command
</PRE>
</UL>
<P><B>Examples:</B>
</P>
<PRE>fix 3 all pour 1000 2 29494 region myblock
fix 2 all pour 10000 1 19985583 region disk vol 0.33 100 rate 1.0 diam range 0.9 1.1
fix 2 all pour 10000 1 19985583 region disk diam poly 2 0.7 0.4 1.5 0.6
fix ins all pour 500 1 4767548 vol 0.8 10 region slab mol object rigid myRigid
</PRE>
<P><B>Description:</B>
</P>
<P>Insert finite-size particles or molecules into the simulation box
every few timesteps within a specified region until N particles or
molecules have been inserted. This is typically used to model the
pouring of granular particles into a container under the influence of
gravity. For the remainder of this doc page, a single inserted atom
or molecule is referred to as a "particle".
</P>
<P>If inserted particles are individual atoms, they are assigned the
specified atom type. If they are molecules, the type of each atom in
the inserted molecule is specified in the file read by the
<A HREF = "molecule.html">molecule</A> command, and those values are added to the
specified atom type. E.g. if the file specifies atom types 1,2,3, and
those are the atom types you want for inserted molecules, then specify
<I>type</I> = 0. If you specify <I>type</I> = 2, the in the inserted molecule
will have atom types 3,4,5.
</P>
<P>All atoms in the inserted particle are assigned to two groups: the
default group "all" and the group specified in the fix pour command
(which can also be "all").
</P>
<P>This command must use the <I>region</I> keyword to define an insertion
volume. The specified region must have been previously defined with a
<A HREF = "region.html">region</A> command. It must be of type <I>block</I> or a z-axis
<I>cylinder</I> and must be defined with side = <I>in</I>. The cylinder style
of region can only be used with 3d simulations.
</P>
<P>Individual atoms are inserted, unless the <I>mol</I> keyword is used. It
specifies a <I>template-ID</I> previously defined using the
<A HREF = "molecule.html">molecule</A> command, which reads a file that defines the
molecule. The coordinates, atom types, center-of-mass, moments of
inertia, etc, as well as any bond/angle/etc and special neighbor
information for the molecule can be specified in the molecule file.
See the <A HREF = "molecule.html">molecule</A> command for details. The only
settings required to be in this file are the coordinates and types of
atoms in the molecule.
</P>
<P>If the molecule template contains more than one molecule, the relative
probability of depositing each molecule can be specified by the
<I>molfrac</I> keyword. N relative probablities, each from 0.0 to 1.0, are
specified, where N is the number of molecules in the template. Each
time a molecule is inserted, a random number is used to sample from
the list of relative probabilities. The N values must sum to 1.0.
</P>
<P>If you wish to insert molecules via the <I>mol</I> keyword, that will be
treated as rigid bodies, use the <I>rigid</I> keyword, specifying as its
value the ID of a separate <A HREF = "fix_rigid_small.html">fix rigid/small</A>
command which also appears in your input script.
</P>
<P>If you wish to insert molecules via the <I>mol</I> keyword, that will have
their bonds or angles constrained via SHAKE, use the <I>shake</I> keyword,
specifying as its value the ID of a separate <A HREF = "fix_shake.html">fix
shake</A> command which also appears in your input script.
</P>
<P>Each timestep particles are inserted, they are placed randomly inside
the insertion volume so as to mimic a stream of poured particles. If
they are molecules they are also oriented randomly. Each atom in the
particle is tested for overlaps with existing particles, including
effects due to periodic boundary conditions if applicable. If an
overlap is detected, another random insertion attempt is made; see the
<I>vol</I> keyword discussion below. The larger the volume of the
insertion region, the more particles that can be inserted at any one
timestep. Particles are inserted again after enough time has elapsed
that the previously inserted particles fall out of the insertion
volume under the influence of gravity. Insertions continue every so
many timesteps until the desired # of particles has been inserted.
</P>
<P>All other keywords are optional with defaults as shown below.
</P>
<P>The <I>diam</I> option is only used when inserting atoms and specifes the
diameters of inserted particles. There are 3 styles: <I>one</I>, <I>range</I>,
or <I>poly</I>. For <I>one</I>, all particles will have diameter <I>D</I>. For
<I>range</I>, the diameter of each particle will be chosen randomly and
uniformly between the specified <I>Dlo</I> and <I>Dhi</I> bounds. For <I>poly</I>, a
series of <I>Npoly</I> diameters is specified. For each diameter a
percentage value from 0.0 to 1.0 is also specified. The <I>Npoly</I>
percentages must sum to 1.0. For the example shown above with "diam 2
0.7 0.4 1.5 0.6", all inserted particles will have a diameter of 0.7
or 1.5. 40% of the particles will be small; 60% will be large.
</P>
<P>Note that for molecule insertion, the diameters of individual atoms in
the molecule can be specified in the file read by the
<A HREF = "molecule.html">molecule</A> command. If not specified, the diameter of
each atom in the molecule has a default diameter of 1.0.
</P>
<P>The <I>dens</I> and <I>vel</I> options enable inserted particles to have a range
of densities or xy velocities. The specific values for a particular
inserted particle will be chosen randomly and uniformly between the
specified bounds. The <I>vz</I> or <I>vy</I> value for option <I>vel</I> assigns a
z-velocity (3d) or y-velocity (2d) to each inserted particle.
</P>
<P>The <I>vol</I> option specifies what volume fraction of the insertion
volume will be filled with particles. For particles with a size
specified by the <I>diam range</I> keyword, they are assumed to all be of
maximum diamter <I>Dhi</I> for purposes of computing their contribution to
the volume fraction.
</P>
<P>The higher the volume fraction value, the more particles are inserted
each timestep. Since inserted particles cannot overlap, the maximum
volume fraction should be no higher than about 0.6. Each timestep
particles are inserted, LAMMPS will make up to a total of M tries to
insert the new particles without overlaps, where M = # of inserted
particles * Nattempt. If LAMMPS is unsuccessful at completing all
insertions, it prints a warning.
</P>
<P>The <I>rate</I> option moves the insertion volume in the z direction (3d)
or y direction (2d). This enables pouring particles from a
successively higher height over time.
</P>
<P>IMPORTANT NOTE: If you are monitoring the temperature of a system
where the particle count is changing due to adding particles, you
typically should use the <A HREF = "compute_modify.html">compute_modify dynamic
yes</A> command for the temperature compute you are
using.
</P>
<HR>
<P><B>Restart, fix_modify, output, run start/stop, minimize info:</B>
</P>
<P>No information about this fix is written to <A HREF = "restart.html">binary restart
files</A>. This means you must be careful when restarting a
pouring simulation, when the restart file was written in the middle of
the pouring operation. Specifically, you should use a new fix pour
command in the input script for the restarted simulation that
continues the operation. You will need to adjust the arguments of the
original fix pour command to do this.
</P>
<P>Also note that because the state of the random number generator is not
saved in restart files, you cannot do "exact" restarts with this fix,
where the simulation continues on the same as if no restart had taken
place. However, in a statistical sense, a restarted simulation should
produce the same behavior if you adjust the fix pour parameters
appropriately.
</P>
<P>None of the <A HREF = "fix_modify.html">fix_modify</A> options are relevant to this
fix. No global or per-atom quantities are stored by this fix for
access by various <A HREF = "Section_howto.html#howto_15">output commands</A>. No
parameter of this fix can be used with the <I>start/stop</I> keywords of
the <A HREF = "run.html">run</A> command. This fix is not invoked during <A HREF = "minimize.html">energy
minimization</A>.
</P>
<P><B>Restrictions:</B>
</P>
<P>This fix is part of the GRANULAR 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>For 3d simulations, a gravity fix in the -z direction must be defined
for use in conjunction with this fix. For 2d simulations, gravity
must be defined in the -y direction.
</P>
<P>The specified insertion region cannot be a "dynamic" region, as
defined by the <A HREF = "region.html">region</A> command.
</P>
<P><B>Related commands:</B>
</P>
<P><A HREF = "fix_deposit.html">fix_deposit</A>, <A HREF = "fix_gravity.html">fix_gravity</A>,
<A HREF = "region.html">region</A>
</P>
<P><B>Default:</B>
</P>
<P>Insertions are performed for individual particles, i.e. no <I>mol</I>
setting is defined. If the <I>mol</I> keyword is used, the default for
<I>molfrac</I> is an equal probabilities for all molecules in the template.
Additional option defaults are diam = one 1.0, dens = 1.0 1.0, vol =
0.25 50, rate = 0.0, vel = 0.0 0.0 0.0 0.0 0.0.
</P>
</HTML>