2
<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>
12
<H3>fix ave/chunk command
16
<PRE>fix ID group-ID ave/chunk Nevery Nrepeat Nfreq chunkID value1 value2 ... keyword args ...
18
<UL><LI>ID, group-ID are documented in <A HREF = "fix.html">fix</A> command
20
<LI>ave/chunk = style name of this fix command
22
<LI>Nevery = use input values every this many timesteps
24
<LI>Nrepeat = # of times to use input values for calculating averages
26
<LI>Nfreq = calculate averages every this many timesteps
28
<LI>chunkID = ID of <A HREF = "compute_chunk_atom.html">compute chunk/atom</A> command
30
<LI>one or more input values can be listed
32
<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
34
<PRE> vx,vy,vz,fx,fy,fz = atom attribute (velocity, force component)
35
density/number, density/mass = number or mass density
37
c_ID = per-atom vector calculated by a compute with ID
38
c_ID[I] = Ith column of per-atom array calculated by a compute with ID
39
f_ID = per-atom vector calculated by a fix with ID
40
f_ID[I] = Ith column of per-atom array calculated by a fix with ID
41
v_name = per-atom vector calculated by an atom-style variable with name
43
<LI>zero or more keyword/arg pairs may be appended
45
<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>
47
<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
48
all = output is sum of atoms across all <I>Nrepeat</I> samples, divided by atom count
49
sample = output is sum of <I>Nrepeat</I> sample averages, divided by <I>Nrepeat</I>
50
none = output is sum of <I>Nrepeat</I> sums, divided by <I>Nrepeat</I>
51
<I>ave</I> args = <I>one</I> or <I>running</I> or <I>window M</I>
52
one = output new average value every Nfreq steps
53
running = output cumulative average of all previous Nfreq steps
54
window M = output average of M most recent Nfreq steps
55
<I>bias</I> arg = bias-ID
56
bias-ID = ID of a temperature compute that removes a velocity bias for temperature calculation
57
<I>adof</I> value = dof_per_atom
58
dof_per_atom = define this many degrees-of-freedom per atom for temperature calculation
59
<I>cdof</I> value = dof_per_chunk
60
dof_per_chunk = define this many degrees-of-freedom per chunk for temperature calculation
61
<I>file</I> arg = filename
62
filename = file to write results to
63
<I>overwrite</I> arg = none = overwrite output file with only latest output
64
<I>title1</I> arg = string
65
string = text to print as 1st line of output file
66
<I>title2</I> arg = string
67
string = text to print as 2nd line of output file
68
<I>title3</I> arg = string
69
string = text to print as 3rd line of output file
75
<PRE>fix 1 all ave/chunk 10000 1 10000 binchunk c_myCentro title1 "My output values"
76
fix 1 flow ave/chunk 100 10 1000 molchunk vx vz norm sample file vel.profile
77
fix 1 flow ave/chunk 100 5 1000 binchunk density/mass ave running
78
fix 1 flow ave/chunk 100 5 1000 binchunk density/mass ave running
80
<P><B>IMPORTANT NOTE:</B>
82
<P>If you are trying to replace an older fix ave/spatial command with the
83
newer, more flexible fix ave/chunk and <A HREF = "compute_chunk_atom.html">compute
84
chunk/atom</A> commands, you simply need to split
85
the fix ave/spatial arguments across the two new commands. For
86
example, this command:
88
<PRE>fix 1 flow ave/spatial 100 10 1000 y 0.0 1.0 vx vz norm sample file vel.profile
90
<P>could be replaced by:
92
<PRE>compute cc1 flow chunk/atom bin/1d y 0.0 1.0
93
fix 1 flow ave/chunk 100 10 1000 cc1 vx vz norm sample file vel.profile
95
<P><B>Description:</B>
97
<P>Use one or more per-atom vectors as inputs every few timesteps, sum
98
the values over the atoms in each chunk at each timestep, then average
99
the per-chunk values over longer timescales. The resulting chunk
100
averages can be used by other <A HREF = "Section_howto.html#howto_15">output
101
commands</A> such as <A HREF = "thermo_style.html">thermo_style
102
custom</A>, and can also be written to a file.
104
<P>In LAMMPS, chunks are collections of atoms defined by a <A HREF = "compute_chunk_atom.html">compute
105
chunk/atom</A> command, which assigns each atom
106
to a single chunk (or no chunk). The ID for this command is specified
107
as chunkID. For example, a single chunk could be the atoms in a
108
molecule or atoms in a spatial bin. See the <A HREF = "compute_chunk_atom.html">compute
109
chunk/atom</A> doc page and "<A HREF = "Section_howto.html#howto_23">Section_howto
110
23</A> for details of how chunks can be
111
defined and examples of how they can be used to measure properties of
114
<P>Note that only atoms in the specified group contribute to the summing
115
and averaging calculations. The <A HREF = "compute_chunk_atom.html">compute
116
chunk/atom</A> command defines its own group as
117
well as an optional region. Atoms will have a chunk ID = 0, meaning
118
they belong to no chunk, if they are not in that group or region.
119
Thus you can specify the "all" group for this command if you simply
120
want to use the chunk definitions provided by chunkID.
122
<P>Each specified per-atom value can be an atom attribute (position,
123
velocity, force component), a mass or number density, or the result of
124
a <A HREF = "compute.html">compute</A> or <A HREF = "fix.html">fix</A> or the evaluation of an
125
atom-style <A HREF = "variable.html">variable</A>. In the latter cases, the
126
compute, fix, or variable must produce a per-atom quantity, not a
127
global quantity. Note that the <A HREF = "compute_property_atom.html">compute
128
property/atom</A> command provides access to
129
any attribute defined and stored by atoms. If you wish to
130
time-average global quantities from a compute, fix, or variable, then
131
see the <A HREF = "fix_ave_time.html">fix ave/time</A> command.
133
<P><A HREF = "compute.html">Computes</A> that produce per-atom quantities are those
134
which have the word <I>atom</I> in their style name. See the doc pages for
135
individual <A HREF = "fix.html">fixes</A> to determine which ones produce per-atom
136
quantities. <A HREF = "variable.html">Variables</A> of style <I>atom</I> are the only
137
ones that can be used with this fix since all other styles of variable
138
produce global quantities.
140
<P>The per-atom values of each input vector are summed and averaged
141
independently of the per-atom values in other input vectors.
143
<P>IMPORTANT NOTE: This fix works by creating an array of size <I>Nchunk</I>
144
by Nvalues on each processor. <I>Nchunk</I> is the number of chunks which
145
is defined by the <A HREF = "doc/compute_chunk_atom.html">compute chunk/atom</A>
146
command. Nvalues is the number of input values specified. Each
147
processor loops over its atoms, tallying its values to the appropriate
148
chunk. Then the entire array is summed across all processors. This
149
means that using a large number of chunks will incur an overhead in
150
memory and computational cost (summing across processors), so be
151
careful to define a reasonable number of chunks.
155
<P>The <I>Nevery</I>, <I>Nrepeat</I>, and <I>Nfreq</I> arguments specify on what
156
timesteps the input values will be accessed and contribute to the
157
average. The final averaged quantities are generated on timesteps
158
that are a multiples of <I>Nfreq</I>. The average is over <I>Nrepeat</I>
159
quantities, computed in the preceding portion of the simulation every
160
<I>Nevery</I> timesteps. <I>Nfreq</I> must be a multiple of <I>Nevery</I> and
161
<I>Nevery</I> must be non-zero even if <I>Nrepeat</I> is 1. Also, the timesteps
162
contributing to the average value cannot overlap, i.e. Nfreq >
163
(Nrepeat-1)*Nevery is required.
165
<P>For example, if Nevery=2, Nrepeat=6, and Nfreq=100, then values on
166
timesteps 90,92,94,96,98,100 will be used to compute the final average
167
on timestep 100. Similarly for timesteps 190,192,194,196,198,200 on
168
timestep 200, etc. If Nrepeat=1 and Nfreq = 100, then no time
169
averaging is done; values are simply generated on timesteps
172
<P>Each input value can also be averaged over the atoms in each chunk.
173
The way the averaging is done across the <I>Nrepeat</I> timesteps to
174
produce output on the <I>Nfreq</I> timesteps, and across multiple <I>Nfreq</I>
175
outputs, is determined by the <I>norm</I> and <I>ave</I> keyword settings, as
178
<P>IMPORTANT NOTE: To perform per-chunk averaging within a <I>Nfreq</I> time
179
window, the number of chunks <I>Nchunk</I> defined by the <A HREF = "compute_chunk_atom.html">compute
180
chunk/atom</A> command must remain constant. If
181
the <I>ave</I> keyword is set to <I>running</I> or <I>window</I> then <I>Nchunk</I> must
182
remain constant for the duration of the simulation. This fix forces
183
the chunk/atom compute specified by chunkID to hold <I>Nchunk</I> constant
184
for the appropriate time windows, by not allowing it to re-calcualte
185
<I>Nchunk</I>, which can also affect how it assigns chunk IDs to atoms.
186
More details are given on the <A HREF = "compute_chunk_atom.html">compute
187
chunk/atom</A> doc page.
191
<P>The atom attribute values (vx,vy,vz,fx,fy,fz) are self-explanatory.
192
As noted above, any other atom attributes can be used as input values
193
to this fix by using the <A HREF = "compute_property_atom.html">compute
194
property/atom</A> command and then specifying
195
an input value from that compute.
197
<P>The <I>density/number</I> value means the number density is computed for
198
each chunk, i.e. number/volume. The <I>density/mass</I> value means the
199
mass density is computed for each chunk, i.e. total-mass/volume. The
200
output values are in units of 1/volume or density (mass/volume). See
201
the <A HREF = "units.html">units</A> command doc page for the definition of density
202
for each choice of units, e.g. gram/cm^3. If the chunks defined by
203
the <A HREF = "compute_chunk_atom.html">compute chunk/atom</A> command are spatial
204
bins, the volume is the bin volume. Otherwise it is the volume of the
205
entire simulation box.
207
<P>The <I>temp</I> value means the temperature is computed for each chunk, by
208
the formula KE = DOF/2 k T, where KE = total kinetic energy of the
209
chunk of atoms (sum of 1/2 m v^2), DOF = the total number of degrees
210
of freedom for all atoms in the chunk, k = Boltzmann constant, and T =
213
<P>The DOF is calculated as N*adof + cdof, where N = number of atoms in
214
the chunk, adof = degrees of freedom per atom, and cdof = degrees of
215
freedom per chunk. By default adof = 2 or 3 = dimensionality of
216
system, as set via the <A HREF = "dimension.html">dimension</A> command, and cdof =
217
0.0. This gives the usual formula for temperature.
219
<P>Note that currently this temperature only includes translational
220
degrees of freedom for each atom. No rotational degrees of freedom
221
are included for finite-size particles. Also no degrees of freedom
222
are subtracted for any velocity bias or constraints that are applied,
223
such as <A HREF = "compute_temp_partial.html">compute temp/partial</A>, or <A HREF = "fix_shake.html">fix
224
shake</A> or <A HREF = "fix_rigid.html">fix rigid</A>. This is because
225
those degrees of freedom (e.g. a constrained bond) could apply to sets
226
of atoms that are both included and excluded from a specific chunk,
227
and hence the concept is somewhat ill-defined. In some cases, you can
228
use the <I>adof</I> and <I>cdof</I> keywords to adjust the calculated degress of
229
freedom appropriately, as explained below.
231
<P>Also note that a bias can be subtracted from atom velocities before
232
they are used in the above formula for KE, by using the <I>bias</I>
233
keyword. This allows, for example, a thermal temperature to be
234
computed after removal of a flow velocity profile.
236
<P>Note that the per-chunk temperature calculated by this fix and the
237
<A HREF = "compute_temp_chunk.html">compute temp/chunk</A> command can be different.
238
The compute calculates the temperature for each chunk for a single
239
snapshot. This fix can do that but can also time average those values
240
over many snapshots, or it can compute a temperature as if the atoms
241
in the chunk on different timesteps were collected together as one set
242
of atoms to calculate their temperature. The compute allows the
243
center-of-mass velocity of each chunk to be subtracted before
244
calculating the temperature; this fix does not.
246
<P>If a value begins with "c_", a compute ID must follow which has been
247
previously defined in the input script. If no bracketed integer is
248
appended, the per-atom vector calculated by the compute is used. If a
249
bracketed integer is appended, the Ith column of the per-atom array
250
calculated by the compute is used. Users can also write code for
251
their own compute styles and <A HREF = "Section_modify.html">add them to LAMMPS</A>.
253
<P>If a value begins with "f_", a fix ID must follow which has been
254
previously defined in the input script. If no bracketed integer is
255
appended, the per-atom vector calculated by the fix is used. If a
256
bracketed integer is appended, the Ith column of the per-atom array
257
calculated by the fix is used. Note that some fixes only produce
258
their values on certain timesteps, which must be compatible with
259
<I>Nevery</I>, else an error results. Users can also write code for their
260
own fix styles and <A HREF = "Section_modify.html">add them to LAMMPS</A>.
262
<P>If a value begins with "v_", a variable name must follow which has
263
been previously defined in the input script. Variables of style
264
<I>atom</I> can reference thermodynamic keywords and various per-atom
265
attributes, or invoke other computes, fixes, or variables when they
266
are evaluated, so this is a very general means of generating per-atom
267
quantities to average within chunks.
271
<P>Additional optional keywords also affect the operation of this fix
274
<P>The <I>norm</I> keyword affects how averaging is done for the per-chunk
275
values that are output every <I>Nfreq</I> timesteps.
277
<P>It the <I>norm</I> setting is <I>all</I>, which is the default, a chunk value is
278
summed over all atoms in all <I>Nrepeat</I> samples, as is the count of
279
atoms in the chunk. The averaged output value for the chunk on the
280
<I>Nfreq</I> timesteps is Total-sum / Total-count. In other words it is an
281
average over atoms across the entire <I>Nfreq</I> timescale.
283
<P>If the <I>norm</I> setting is <I>sample</I>, the chunk value is summed over atoms
284
for each sample, as is the count, and an "average sample value" is
285
computed for each sample, i.e. Sample-sum / Sample-count. The outuput
286
value for the chunk on the <I>Nfreq</I> timesteps is the average of the
287
<I>Nrepeat</I> "average sample values", i.e. the sum of <I>Nrepeat</I> "average
288
sample values" divided by <I>Nrepeat</I>. In other words it is an average
291
<P>If the <I>norm</I> setting is <I>none</I>, a similar computation as for the
292
<I>sample</I> seting is done, except the individual "average sample values"
293
are "summed sample values". A summed sample value is simply the chunk
294
value summed over atoms in the sample, without dividing by the number
295
of atoms in the sample. The outuput value for the chunk on the
296
<I>Nfreq</I> timesteps is the average of the <I>Nrepeat</I> "summed sample
297
values", i.e. the sum of <I>Nrepeat</I> "summed sample values" divided by
300
<P>The <I>ave</I> keyword determines how the per-chunk values produced every
301
<I>Nfreq</I> steps are averaged with values produced on previous steps that
302
were multiples of <I>Nfreq</I>, before they are accessed by another output
303
command or written to a file.
305
<P>If the <I>ave</I> setting is <I>one</I>, which is the default, then the chunk
306
values produced on timesteps that are multiples of <I>Nfreq</I> are
307
independent of each other; they are output as-is without further
310
<P>If the <I>ave</I> setting is <I>running</I>, then the chunk values produced on
311
timesteps that are multiples of <I>Nfreq</I> are summed and averaged in a
312
cumulative sense before being output. Each output chunk value is thus
313
the average of the chunk value produced on that timestep with all
314
preceding values for the same chunk. This running average begins when
315
the fix is defined; it can only be restarted by deleting the fix via
316
the <A HREF = "unfix.html">unfix</A> command, or re-defining the fix by
319
<P>If the <I>ave</I> setting is <I>window</I>, then the chunk values produced on
320
timesteps that are multiples of <I>Nfreq</I> are summed and averaged within
321
a moving "window" of time, so that the last M values for the same
322
chunk are used to produce the output. E.g. if M = 3 and Nfreq = 1000,
323
then the output on step 10000 will be the average of the individual
324
chunk values on steps 8000,9000,10000. Outputs on early steps will
325
average over less than M values if they are not available.
327
<P>The <I>bias</I> keyword specifies the ID of a temperature compute that
328
removes a "bias" velocity from each atom, specified as <I>bias-ID</I>. It
329
is only used when the <I>temp</I> value is calculated, to compute the
330
thermal temperature of each chunk after the translational kinetic
331
energy components have been altered in a prescribed way, e.g. to
332
remove a flow velocity profile. See the doc pages for individual
333
computes that calculate a temperature to see which ones implement a
336
<P>The <I>adof</I> and <I>cdof</I> keywords define the values used in the degree of
337
freedom (DOF) formula described above for for temperature calculation
338
for each chunk. They are only used when the <I>temp</I> value is
339
calculated. They can be used to calculate a more appropriate
340
temperature for some kinds of chunks. Here are 3 examples.
342
<P>If spatially binned chunks contain some number of water molecules and
343
<A HREF = "fix_shake.html">fix shake</A> is used to make each molecule rigid, then
344
you could calculate a temperature with 6 degrees of freedom (DOF) (3
345
translational, 3 rotational) per molecule by setting <I>adof</I> to 2.0.
346
If <A HREF = "compute_temp_partial.html">compute temp/partial</A> is used with the
347
<I>bias</I> keyword to only allow the x component of velocity to contribute
348
to the temperature, then <I>adof</I> = 1.0 would be appropriate. If each
349
chunk consists of a large molecule, with some number of its bonds
350
constrained by <A HREF = "fix_shake.html">fix shake</A> or the entire molecule by
351
<A HREF = "fix_rigid.html">fix rigid/small</A>, then <I>cdof</I> could be set to the
352
remaining degrees of freedom for the entire molecule (entire chunk in
355
<P>The <I>file</I> keyword allows a filename to be specified. Every <I>Nfreq</I>
356
timesteps, a section of chunk info will be written to a text file in
357
the following format. A line with the timestep and number of chunks
358
is written. Then one line per chunk is written, containing the chunk
359
ID (1-Nchunk), an optional original ID value, optional coordinate
360
values for chunks that represent spatial bins, the number of atoms in
361
the chunk, and one or more calculated values. More explanation of the
362
optional values is given below. The number of values in each line
363
corresponds to the number of values specified in the fix ave/chunk
364
command. The number of atoms and the value(s) are summed or average
365
quantities, as explained above.
367
<P>The <I>overwrite</I> keyword will continuously overwrite the output file
368
with the latest output, so that it only contains one timestep worth of
369
output. This option can only be used with the <I>ave running</I> setting.
371
<P>The <I>title1</I> and <I>title2</I> and <I>title3</I> keywords allow specification of
372
the strings that will be printed as the first 3 lines of the output
373
file, assuming the <I>file</I> keyword was used. LAMMPS uses default
374
values for each of these, so they do not need to be specified.
376
<P>By default, these header lines are as follows:
378
<PRE># Chunk-averaged data for fix ID and group name
379
# Timestep Number-of-chunks
380
# Chunk (OrigID) (Coord1) (Coord2) (Coord3) Ncount value1 value2 ...
382
<P>In the first line, ID and name are replaced with the fix-ID and group
383
name. The second line describes the two values that are printed at
384
the first of each section of output. In the third line the values are
385
replaced with the appropriate value names, e.g. fx or c_myCompute<B>2</B>.
387
<P>The words in parenthesis only appear with corresponding columns if the
388
chunk style specified for the <A HREF = "compute_chunk_atom.html">compute
389
chunk/atom</A> command supports them. The OrigID
390
column is only used if the <I>compress</I> keyword was set to <I>yes</I> for the
391
<A HREF = "compute_chunk_atom.html">compute chunk/atom</A> command. This means that
392
the original chunk IDs (e.g. molecule IDs) will have been compressed
393
to remove chunk IDs with no atoms assigned to them. Thus a compresed
394
chunk ID of 3 may correspond to an original chunk ID or molecule ID of
395
415. The OrigID column will list 415 for the 3rd chunk.
397
<P>The CoordN columns only appear if a <I>binning</I> style was used in the
398
<A HREF = "compute_chunk_atom.html">compute chunk/atom</A> command. For <I>bin/1d</I>,
399
<I>bin/2d</I>, and <I>bin/3d</I> styles the column values are the center point
400
of the bin in the corresponding dimension. Just Coord1 is used for
401
<I>bin/1d</I>, Coord2 is added for <I>bin/2d</I>, Coord3 is added for <I>bin/3d</I>.
403
<P>Note that if the value of the <I>units</I> keyword used in the <A HREF = "compute_chunk_atom.html">compute
404
chunk/atom command</A> is <I>box</I> or <I>lattice</I>, the
405
coordinate values will be in distance <A HREF = "units.html">units</A>. If the
406
value of the <I>units</I> keyword is <I>reduced</I>, the coordinate values will
407
be in unitless reduced units (0-1).
411
<P><B>Restart, fix_modify, output, run start/stop, minimize info:</B>
413
<P>No information about this fix is written to <A HREF = "restart.html">binary restart
414
files</A>. None of the <A HREF = "fix_modify.html">fix_modify</A> options
415
are relevant to this fix.
417
<P>This fix computes a global array of values which can be accessed by
418
various <A HREF = "Section_howto.html#howto_15">output commands</A>. The values can
419
only be accessed on timesteps that are multiples of <I>Nfreq</I> since that
420
is when averaging is performed. The global array has # of rows =
421
the number of chunks <I>Nchunk</I> as calculated by the specified <A HREF = "compute_chunk_atom.html">compute
422
chunk/atom</A> command. The # of columns =
423
M+1+Nvalues, where M = 1 to 4, depending on whether the optional
424
columns for OrigID and CoordN are used, as explained above.
425
Following the optional columns, the next column contains the count of
426
atoms in the chunk, and the remaining columns are the Nvalue
427
quantities. When the array is accessed with a row I that exceeds the
428
current number of chunks, than a 0.0 is returned by the fix instead of
429
an error, since the number of chunks can vary as a simulation runs
430
depending on how that value is computed by the compute chunk/atom
433
<P>The array values calculated by this fix are treated as "intensive",
434
since they are typically already normalized by the count of atoms in
437
<P>No parameter of this fix can be used with the <I>start/stop</I> keywords of
438
the <A HREF = "run.html">run</A> command. This fix is not invoked during <A HREF = "minimize.html">energy
441
<P><B>Restrictions:</B> none
443
<P><B>Related commands:</B>
445
<P><A HREF = "compute.html">compute</A>, <A HREF = "fix_ave_atom.html">fix ave/atom</A>, <A HREF = "fix_ave_histo.html">fix
446
ave/histo</A>, <A HREF = "fix_ave_time.html">fix ave/time</A>,
447
<A HREF = "variable.html">variable</A>, <A HREF = "fix_ave_correlate.html">fix ave/correlate</A>
451
<P>The option defaults are norm = all, ave = one, bias = none, no file output, and
452
title 1,2,3 = strings as described above.