lammps/doc/fix_ave_chunk.html

455 lines
24 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>fix ave/chunk command
</H3>
<P><B>Syntax:</B>
</P>
<PRE>fix ID group-ID ave/chunk Nevery Nrepeat Nfreq chunkID value1 value2 ... keyword args ...
</PRE>
<UL><LI>ID, group-ID are documented in <A HREF = "fix.html">fix</A> command
<LI>ave/chunk = style name of this fix command
<LI>Nevery = use input values every this many timesteps
<LI>Nrepeat = # of times to use input values for calculating averages
<LI>Nfreq = calculate averages every this many timesteps
<LI>chunkID = ID of <A HREF = "compute_chunk_atom.html">compute chunk/atom</A> command
<LI>one or more input values can be listed
<LI>value = vx, vy, vz, fx, fy, fz, density/mass, density/number, temp, c_ID, c_ID[I], f_ID, f_ID[I], v_name
<PRE> vx,vy,vz,fx,fy,fz = atom attribute (velocity, force component)
density/number, density/mass = number or mass density
temp = temperature
c_ID = per-atom vector calculated by a compute with ID
c_ID[I] = Ith column of per-atom array calculated by a compute with ID
f_ID = per-atom vector calculated by a fix with ID
f_ID[I] = Ith column of per-atom array calculated by a fix with ID
v_name = per-atom vector calculated by an atom-style variable with name
</PRE>
<LI>zero or more keyword/arg pairs may be appended
<LI>keyword = <I>norm</I> or <I>ave</I> or <I>bias</I> or <I>adof</I> or <I>cdof</I> or <I>file</I> or <I>overwrite</I> or <I>title1</I> or <I>title2</I> or <I>title3</I>
<PRE> <I>norm</I> arg = <I>all</I> or <I>sample</I> or <I>none</I> = how output on <I>Nfreq</I> steps is normalized
all = output is sum of atoms across all <I>Nrepeat</I> samples, divided by atom count
sample = output is sum of <I>Nrepeat</I> sample averages, divided by <I>Nrepeat</I>
none = output is sum of <I>Nrepeat</I> sums, divided by <I>Nrepeat</I>
<I>ave</I> args = <I>one</I> or <I>running</I> or <I>window M</I>
one = output new average value every Nfreq steps
running = output cumulative average of all previous Nfreq steps
window M = output average of M most recent Nfreq steps
<I>bias</I> arg = bias-ID
bias-ID = ID of a temperature compute that removes a velocity bias for temperature calculation
<I>adof</I> value = dof_per_atom
dof_per_atom = define this many degrees-of-freedom per atom for temperature calculation
<I>cdof</I> value = dof_per_chunk
dof_per_chunk = define this many degrees-of-freedom per chunk for temperature calculation
<I>file</I> arg = filename
filename = file to write results to
<I>overwrite</I> arg = none = overwrite output file with only latest output
<I>title1</I> arg = string
string = text to print as 1st line of output file
<I>title2</I> arg = string
string = text to print as 2nd line of output file
<I>title3</I> arg = string
string = text to print as 3rd line of output file
</PRE>
</UL>
<P><B>Examples:</B>
</P>
<PRE>fix 1 all ave/chunk 10000 1 10000 binchunk c_myCentro title1 "My output values"
fix 1 flow ave/chunk 100 10 1000 molchunk vx vz norm sample file vel.profile
fix 1 flow ave/chunk 100 5 1000 binchunk density/mass ave running
fix 1 flow ave/chunk 100 5 1000 binchunk density/mass ave running
</PRE>
<P><B>IMPORTANT NOTE:</B>
</P>
<P>If you are trying to replace an older fix ave/spatial command with the
newer, more flexible fix ave/chunk and <A HREF = "compute_chunk_atom.html">compute
chunk/atom</A> commands, you simply need to split
the fix ave/spatial arguments across the two new commands. For
example, this command:
</P>
<PRE>fix 1 flow ave/spatial 100 10 1000 y 0.0 1.0 vx vz norm sample file vel.profile
</PRE>
<P>could be replaced by:
</P>
<PRE>compute cc1 flow chunk/atom bin/1d y 0.0 1.0
fix 1 flow ave/chunk 100 10 1000 cc1 vx vz norm sample file vel.profile
</PRE>
<P><B>Description:</B>
</P>
<P>Use one or more per-atom vectors as inputs every few timesteps, sum
the values over the atoms in each chunk at each timestep, then average
the per-chunk values over longer timescales. The resulting chunk
averages can be used by other <A HREF = "Section_howto.html#howto_15">output
commands</A> such as <A HREF = "thermo_style.html">thermo_style
custom</A>, and can also be written to a file.
</P>
<P>In LAMMPS, chunks are collections of atoms defined by a <A HREF = "compute_chunk_atom.html">compute
chunk/atom</A> command, which assigns each atom
to a single chunk (or no chunk). The ID for this command is specified
as chunkID. For example, a single chunk could be the atoms in a
molecule or atoms in a spatial bin. See the <A HREF = "compute_chunk_atom.html">compute
chunk/atom</A> doc page and "<A HREF = "Section_howto.html#howto_23">Section_howto
23</A> for details of how chunks can be
defined and examples of how they can be used to measure properties of
a system.
</P>
<P>Note that only atoms in the specified group contribute to the summing
and averaging calculations. The <A HREF = "compute_chunk_atom.html">compute
chunk/atom</A> command defines its own group as
well as an optional region. Atoms will have a chunk ID = 0, meaning
they belong to no chunk, if they are not in that group or region.
Thus you can specify the "all" group for this command if you simply
want to use the chunk definitions provided by chunkID.
</P>
<P>Each specified per-atom value can be an atom attribute (position,
velocity, force component), a mass or number density, or the result of
a <A HREF = "compute.html">compute</A> or <A HREF = "fix.html">fix</A> or the evaluation of an
atom-style <A HREF = "variable.html">variable</A>. In the latter cases, the
compute, fix, or variable must produce a per-atom quantity, not a
global quantity. Note that the <A HREF = "compute_property_atom.html">compute
property/atom</A> command provides access to
any attribute defined and stored by atoms. If you wish to
time-average global quantities from a compute, fix, or variable, then
see the <A HREF = "fix_ave_time.html">fix ave/time</A> command.
</P>
<P><A HREF = "compute.html">Computes</A> that produce per-atom quantities are those
which have the word <I>atom</I> in their style name. See the doc pages for
individual <A HREF = "fix.html">fixes</A> to determine which ones produce per-atom
quantities. <A HREF = "variable.html">Variables</A> of style <I>atom</I> are the only
ones that can be used with this fix since all other styles of variable
produce global quantities.
</P>
<P>The per-atom values of each input vector are summed and averaged
independently of the per-atom values in other input vectors.
</P>
<P>IMPORTANT NOTE: This fix works by creating an array of size <I>Nchunk</I>
by Nvalues on each processor. <I>Nchunk</I> is the number of chunks which
is defined by the <A HREF = "doc/compute_chunk_atom.html">compute chunk/atom</A>
command. Nvalues is the number of input values specified. Each
processor loops over its atoms, tallying its values to the appropriate
chunk. Then the entire array is summed across all processors. This
means that using a large number of chunks will incur an overhead in
memory and computational cost (summing across processors), so be
careful to define a reasonable number of chunks.
</P>
<HR>
<P>The <I>Nevery</I>, <I>Nrepeat</I>, and <I>Nfreq</I> arguments specify on what
timesteps the input values will be accessed and contribute to the
average. The final averaged quantities are generated on timesteps
that are a multiples of <I>Nfreq</I>. The average is over <I>Nrepeat</I>
quantities, computed in the preceding portion of the simulation every
<I>Nevery</I> timesteps. <I>Nfreq</I> must be a multiple of <I>Nevery</I> and
<I>Nevery</I> must be non-zero even if <I>Nrepeat</I> is 1. Also, the timesteps
contributing to the average value cannot overlap, i.e. Nfreq >
(Nrepeat-1)*Nevery is required.
</P>
<P>For example, if Nevery=2, Nrepeat=6, and Nfreq=100, then values on
timesteps 90,92,94,96,98,100 will be used to compute the final average
on timestep 100. Similarly for timesteps 190,192,194,196,198,200 on
timestep 200, etc. If Nrepeat=1 and Nfreq = 100, then no time
averaging is done; values are simply generated on timesteps
100,200,etc.
</P>
<P>Each input value can also be averaged over the atoms in each chunk.
The way the averaging is done across the <I>Nrepeat</I> timesteps to
produce output on the <I>Nfreq</I> timesteps, and across multiple <I>Nfreq</I>
outputs, is determined by the <I>norm</I> and <I>ave</I> keyword settings, as
discussed below.
</P>
<P>IMPORTANT NOTE: To perform per-chunk averaging within a <I>Nfreq</I> time
window, the number of chunks <I>Nchunk</I> defined by the <A HREF = "compute_chunk_atom.html">compute
chunk/atom</A> command must remain constant. If
the <I>ave</I> keyword is set to <I>running</I> or <I>window</I> then <I>Nchunk</I> must
remain constant for the duration of the simulation. This fix forces
the chunk/atom compute specified by chunkID to hold <I>Nchunk</I> constant
for the appropriate time windows, by not allowing it to re-calcualte
<I>Nchunk</I>, which can also affect how it assigns chunk IDs to atoms.
More details are given on the <A HREF = "compute_chunk_atom.html">compute
chunk/atom</A> doc page.
</P>
<HR>
<P>The atom attribute values (vx,vy,vz,fx,fy,fz) are self-explanatory.
As noted above, any other atom attributes can be used as input values
to this fix by using the <A HREF = "compute_property_atom.html">compute
property/atom</A> command and then specifying
an input value from that compute.
</P>
<P>The <I>density/number</I> value means the number density is computed for
each chunk, i.e. number/volume. The <I>density/mass</I> value means the
mass density is computed for each chunk, i.e. total-mass/volume. The
output values are in units of 1/volume or density (mass/volume). See
the <A HREF = "units.html">units</A> command doc page for the definition of density
for each choice of units, e.g. gram/cm^3. If the chunks defined by
the <A HREF = "compute_chunk_atom.html">compute chunk/atom</A> command are spatial
bins, the volume is the bin volume. Otherwise it is the volume of the
entire simulation box.
</P>
<P>The <I>temp</I> value means the temperature is computed for each chunk, by
the formula KE = DOF/2 k T, where KE = total kinetic energy of the
chunk of atoms (sum of 1/2 m v^2), DOF = the total number of degrees
of freedom for all atoms in the chunk, k = Boltzmann constant, and T =
temperature.
</P>
<P>The DOF is calculated as N*adof + cdof, where N = number of atoms in
the chunk, adof = degrees of freedom per atom, and cdof = degrees of
freedom per chunk. By default adof = 2 or 3 = dimensionality of
system, as set via the <A HREF = "dimension.html">dimension</A> command, and cdof =
0.0. This gives the usual formula for temperature.
</P>
<P>Note that currently this temperature only includes translational
degrees of freedom for each atom. No rotational degrees of freedom
are included for finite-size particles. Also no degrees of freedom
are subtracted for any velocity bias or constraints that are applied,
such as <A HREF = "compute_temp_partial.html">compute temp/partial</A>, or <A HREF = "fix_shake.html">fix
shake</A> or <A HREF = "fix_rigid.html">fix rigid</A>. This is because
those degrees of freedom (e.g. a constrained bond) could apply to sets
of atoms that are both included and excluded from a specific chunk,
and hence the concept is somewhat ill-defined. In some cases, you can
use the <I>adof</I> and <I>cdof</I> keywords to adjust the calculated degress of
freedom appropriately, as explained below.
</P>
<P>Also note that a bias can be subtracted from atom velocities before
they are used in the above formula for KE, by using the <I>bias</I>
keyword. This allows, for example, a thermal temperature to be
computed after removal of a flow velocity profile.
</P>
<P>Note that the per-chunk temperature calculated by this fix and the
<A HREF = "compute_temp_chunk.html">compute temp/chunk</A> command can be different.
The compute calculates the temperature for each chunk for a single
snapshot. This fix can do that but can also time average those values
over many snapshots, or it can compute a temperature as if the atoms
in the chunk on different timesteps were collected together as one set
of atoms to calculate their temperature. The compute allows the
center-of-mass velocity of each chunk to be subtracted before
calculating the temperature; this fix does not.
</P>
<P>If a value begins with "c_", a compute ID must follow which has been
previously defined in the input script. If no bracketed integer is
appended, the per-atom vector calculated by the compute is used. If a
bracketed integer is appended, the Ith column of the per-atom array
calculated by the compute is used. Users can also write code for
their own compute styles and <A HREF = "Section_modify.html">add them to LAMMPS</A>.
</P>
<P>If a value begins with "f_", a fix ID must follow which has been
previously defined in the input script. If no bracketed integer is
appended, the per-atom vector calculated by the fix is used. If a
bracketed integer is appended, the Ith column of the per-atom array
calculated by the fix is used. Note that some fixes only produce
their values on certain timesteps, which must be compatible with
<I>Nevery</I>, else an error results. Users can also write code for their
own fix styles and <A HREF = "Section_modify.html">add them to LAMMPS</A>.
</P>
<P>If a value begins with "v_", a variable name must follow which has
been previously defined in the input script. Variables of style
<I>atom</I> can reference thermodynamic keywords and various per-atom
attributes, or invoke other computes, fixes, or variables when they
are evaluated, so this is a very general means of generating per-atom
quantities to average within chunks.
</P>
<HR>
<P>Additional optional keywords also affect the operation of this fix
and its outputs.
</P>
<P>The <I>norm</I> keyword affects how averaging is done for the per-chunk
values that are output every <I>Nfreq</I> timesteps.
</P>
<P>It the <I>norm</I> setting is <I>all</I>, which is the default, a chunk value is
summed over all atoms in all <I>Nrepeat</I> samples, as is the count of
atoms in the chunk. The averaged output value for the chunk on the
<I>Nfreq</I> timesteps is Total-sum / Total-count. In other words it is an
average over atoms across the entire <I>Nfreq</I> timescale.
</P>
<P>If the <I>norm</I> setting is <I>sample</I>, the chunk value is summed over atoms
for each sample, as is the count, and an "average sample value" is
computed for each sample, i.e. Sample-sum / Sample-count. The outuput
value for the chunk on the <I>Nfreq</I> timesteps is the average of the
<I>Nrepeat</I> "average sample values", i.e. the sum of <I>Nrepeat</I> "average
sample values" divided by <I>Nrepeat</I>. In other words it is an average
of an average.
</P>
<P>If the <I>norm</I> setting is <I>none</I>, a similar computation as for the
<I>sample</I> seting is done, except the individual "average sample values"
are "summed sample values". A summed sample value is simply the chunk
value summed over atoms in the sample, without dividing by the number
of atoms in the sample. The outuput value for the chunk on the
<I>Nfreq</I> timesteps is the average of the <I>Nrepeat</I> "summed sample
values", i.e. the sum of <I>Nrepeat</I> "summed sample values" divided by
<I>Nrepeat</I>.
</P>
<P>The <I>ave</I> keyword determines how the per-chunk values produced every
<I>Nfreq</I> steps are averaged with values produced on previous steps that
were multiples of <I>Nfreq</I>, before they are accessed by another output
command or written to a file.
</P>
<P>If the <I>ave</I> setting is <I>one</I>, which is the default, then the chunk
values produced on timesteps that are multiples of <I>Nfreq</I> are
independent of each other; they are output as-is without further
averaging.
</P>
<P>If the <I>ave</I> setting is <I>running</I>, then the chunk values produced on
timesteps that are multiples of <I>Nfreq</I> are summed and averaged in a
cumulative sense before being output. Each output chunk value is thus
the average of the chunk value produced on that timestep with all
preceding values for the same chunk. This running average begins when
the fix is defined; it can only be restarted by deleting the fix via
the <A HREF = "unfix.html">unfix</A> command, or re-defining the fix by
re-specifying it.
</P>
<P>If the <I>ave</I> setting is <I>window</I>, then the chunk values produced on
timesteps that are multiples of <I>Nfreq</I> are summed and averaged within
a moving "window" of time, so that the last M values for the same
chunk are used to produce the output. E.g. if M = 3 and Nfreq = 1000,
then the output on step 10000 will be the average of the individual
chunk values on steps 8000,9000,10000. Outputs on early steps will
average over less than M values if they are not available.
</P>
<P>The <I>bias</I> keyword specifies the ID of a temperature compute that
removes a "bias" velocity from each atom, specified as <I>bias-ID</I>. It
is only used when the <I>temp</I> value is calculated, to compute the
thermal temperature of each chunk after the translational kinetic
energy components have been altered in a prescribed way, e.g. to
remove a flow velocity profile. See the doc pages for individual
computes that calculate a temperature to see which ones implement a
bias.
</P>
<P>The <I>adof</I> and <I>cdof</I> keywords define the values used in the degree of
freedom (DOF) formula described above for for temperature calculation
for each chunk. They are only used when the <I>temp</I> value is
calculated. They can be used to calculate a more appropriate
temperature for some kinds of chunks. Here are 3 examples.
</P>
<P>If spatially binned chunks contain some number of water molecules and
<A HREF = "fix_shake.html">fix shake</A> is used to make each molecule rigid, then
you could calculate a temperature with 6 degrees of freedom (DOF) (3
translational, 3 rotational) per molecule by setting <I>adof</I> to 2.0.
If <A HREF = "compute_temp_partial.html">compute temp/partial</A> is used with the
<I>bias</I> keyword to only allow the x component of velocity to contribute
to the temperature, then <I>adof</I> = 1.0 would be appropriate. If each
chunk consists of a large molecule, with some number of its bonds
constrained by <A HREF = "fix_shake.html">fix shake</A> or the entire molecule by
<A HREF = "fix_rigid.html">fix rigid/small</A>, then <I>cdof</I> could be set to the
remaining degrees of freedom for the entire molecule (entire chunk in
this case).
</P>
<P>The <I>file</I> keyword allows a filename to be specified. Every <I>Nfreq</I>
timesteps, a section of chunk info will be written to a text file in
the following format. A line with the timestep and number of chunks
is written. Then one line per chunk is written, containing the chunk
ID (1-Nchunk), an optional original ID value, optional coordinate
values for chunks that represent spatial bins, the number of atoms in
the chunk, and one or more calculated values. More explanation of the
optional values is given below. The number of values in each line
corresponds to the number of values specified in the fix ave/chunk
command. The number of atoms and the value(s) are summed or average
quantities, as explained above.
</P>
<P>The <I>overwrite</I> keyword will continuously overwrite the output file
with the latest output, so that it only contains one timestep worth of
output. This option can only be used with the <I>ave running</I> setting.
</P>
<P>The <I>title1</I> and <I>title2</I> and <I>title3</I> keywords allow specification of
the strings that will be printed as the first 3 lines of the output
file, assuming the <I>file</I> keyword was used. LAMMPS uses default
values for each of these, so they do not need to be specified.
</P>
<P>By default, these header lines are as follows:
</P>
<PRE># Chunk-averaged data for fix ID and group name
# Timestep Number-of-chunks
# Chunk (OrigID) (Coord1) (Coord2) (Coord3) Ncount value1 value2 ...
</PRE>
<P>In the first line, ID and name are replaced with the fix-ID and group
name. The second line describes the two values that are printed at
the first of each section of output. In the third line the values are
replaced with the appropriate value names, e.g. fx or c_myCompute<B>2</B>.
</P>
<P>The words in parenthesis only appear with corresponding columns if the
chunk style specified for the <A HREF = "compute_chunk_atom.html">compute
chunk/atom</A> command supports them. The OrigID
column is only used if the <I>compress</I> keyword was set to <I>yes</I> for the
<A HREF = "compute_chunk_atom.html">compute chunk/atom</A> command. This means that
the original chunk IDs (e.g. molecule IDs) will have been compressed
to remove chunk IDs with no atoms assigned to them. Thus a compresed
chunk ID of 3 may correspond to an original chunk ID or molecule ID of
415. The OrigID column will list 415 for the 3rd chunk.
</P>
<P>The CoordN columns only appear if a <I>binning</I> style was used in the
<A HREF = "compute_chunk_atom.html">compute chunk/atom</A> command. For <I>bin/1d</I>,
<I>bin/2d</I>, and <I>bin/3d</I> styles the column values are the center point
of the bin in the corresponding dimension. Just Coord1 is used for
<I>bin/1d</I>, Coord2 is added for <I>bin/2d</I>, Coord3 is added for <I>bin/3d</I>.
</P>
<P>Note that if the value of the <I>units</I> keyword used in the <A HREF = "compute_chunk_atom.html">compute
chunk/atom command</A> is <I>box</I> or <I>lattice</I>, the
coordinate values will be in distance <A HREF = "units.html">units</A>. If the
value of the <I>units</I> keyword is <I>reduced</I>, the coordinate values will
be in unitless reduced units (0-1).
</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>. None of the <A HREF = "fix_modify.html">fix_modify</A> options
are relevant to this fix.
</P>
<P>This fix computes a global array of values which can be accessed by
various <A HREF = "Section_howto.html#howto_15">output commands</A>. The values can
only be accessed on timesteps that are multiples of <I>Nfreq</I> since that
is when averaging is performed. The global array has # of rows =
the number of chunks <I>Nchunk</I> as calculated by the specified <A HREF = "compute_chunk_atom.html">compute
chunk/atom</A> command. The # of columns =
M+1+Nvalues, where M = 1 to 4, depending on whether the optional
columns for OrigID and CoordN are used, as explained above.
Following the optional columns, the next column contains the count of
atoms in the chunk, and the remaining columns are the Nvalue
quantities. When the array is accessed with a row I that exceeds the
current number of chunks, than a 0.0 is returned by the fix instead of
an error, since the number of chunks can vary as a simulation runs
depending on how that value is computed by the compute chunk/atom
command.
</P>
<P>The array values calculated by this fix are treated as "intensive",
since they are typically already normalized by the count of atoms in
each chunk.
</P>
<P>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> none
</P>
<P><B>Related commands:</B>
</P>
<P><A HREF = "compute.html">compute</A>, <A HREF = "fix_ave_atom.html">fix ave/atom</A>, <A HREF = "fix_ave_histo.html">fix
ave/histo</A>, <A HREF = "fix_ave_time.html">fix ave/time</A>,
<A HREF = "variable.html">variable</A>, <A HREF = "fix_ave_correlate.html">fix ave/correlate</A>
</P>
<P><B>Default:</B>
</P>
<P>The option defaults are norm = all, ave = one, bias = none, no file output, and
title 1,2,3 = strings as described above.
</P>
</HTML>