lammps/doc/txt/Howto_chunk.txt

206 lines
9.1 KiB
Plaintext

"Higher level section"_Howto.html - "LAMMPS WWW Site"_lws - "LAMMPS
Documentation"_ld - "LAMMPS Commands"_lc :c
:link(lws,http://lammps.sandia.gov)
:link(ld,Manual.html)
:link(lc,Commands_all.html)
:line
Use chunks to calculate system properties :h3
In LAMMS, "chunks" are collections of atoms, as defined by the
"compute chunk/atom"_compute_chunk_atom.html command, which assigns
each atom to a chunk ID (or to no chunk at all). The number of chunks
and the assignment of chunk IDs to atoms can be static or change over
time. Examples of "chunks" are molecules or spatial bins or atoms
with similar values (e.g. coordination number or potential energy).
The per-atom chunk IDs can be used as input to two other kinds of
commands, to calculate various properties of a system:
"fix ave/chunk"_fix_ave_chunk.html
any of the "compute */chunk"_compute.html commands :ul
Here a brief overview for each of the 4 kinds of chunk-related commands
is provided. Then some examples are given of how to compute different
properties with chunk commands.
Compute chunk/atom command: :h4
This compute can assign atoms to chunks of various styles. Only atoms
in the specified group and optional specified region are assigned to a
chunk. Here are some possible chunk definitions:
atoms in same molecule | chunk ID = molecule ID |
atoms of same atom type | chunk ID = atom type |
all atoms with same atom property (charge, radius, etc) | chunk ID = output of compute property/atom |
atoms in same cluster | chunk ID = output of "compute cluster/atom"_compute_cluster_atom.html command |
atoms in same spatial bin | chunk ID = bin ID |
atoms in same rigid body | chunk ID = molecule ID used to define rigid bodies |
atoms with similar potential energy | chunk ID = output of "compute pe/atom"_compute_pe_atom.html |
atoms with same local defect structure | chunk ID = output of "compute centro/atom"_compute_centro_atom.html or "compute coord/atom"_compute_coord_atom.html command :tb(s=|,c=2)
Note that chunk IDs are integer values, so for atom properties or
computes that produce a floating point value, they will be truncated
to an integer. You could also use the compute in a variable that
scales the floating point value to spread it across multiple integers.
Spatial bins can be of various kinds, e.g. 1d bins = slabs, 2d bins =
pencils, 3d bins = boxes, spherical bins, cylindrical bins.
This compute also calculates the number of chunks {Nchunk}, which is
used by other commands to tally per-chunk data. {Nchunk} can be a
static value or change over time (e.g. the number of clusters). The
chunk ID for an individual atom can also be static (e.g. a molecule
ID), or dynamic (e.g. what spatial bin an atom is in as it moves).
Note that this compute allows the per-atom output of other
"computes"_compute.html, "fixes"_fix.html, and
"variables"_variable.html to be used to define chunk IDs for each
atom. This means you can write your own compute or fix to output a
per-atom quantity to use as chunk ID. See the "Modify"_Modify.html
doc pages for info on how to do this. You can also define a "per-atom
variable"_variable.html in the input script that uses a formula to
generate a chunk ID for each atom.
Fix ave/chunk command: :h4
This fix takes the ID of a "compute
chunk/atom"_compute_chunk_atom.html command as input. For each chunk,
it then sums one or more specified per-atom values over the atoms in
each chunk. The per-atom values can be any atom property, such as
velocity, force, charge, potential energy, kinetic energy, stress,
etc. Additional keywords are defined for per-chunk properties like
density and temperature. More generally any per-atom value generated
by other "computes"_compute.html, "fixes"_fix.html, and "per-atom
variables"_variable.html, can be summed over atoms in each chunk.
Similar to other averaging fixes, this fix allows the summed per-chunk
values to be time-averaged in various ways, and output to a file. The
fix produces a global array as output with one row of values per
chunk.
Compute */chunk commands: :h4
The following computes operate on chunks of atoms to produce per-chunk
values. Any compute whose style name ends in "/chunk" is in this
category:
"compute com/chunk"_compute_com_chunk.html
"compute gyration/chunk"_compute_gyration_chunk.html
"compute inertia/chunk"_compute_inertia_chunk.html
"compute msd/chunk"_compute_msd_chunk.html
"compute property/chunk"_compute_property_chunk.html
"compute temp/chunk"_compute_temp_chunk.html
"compute torque/chunk"_compute_vcm_chunk.html
"compute vcm/chunk"_compute_vcm_chunk.html :ul
They each take the ID of a "compute
chunk/atom"_compute_chunk_atom.html command as input. As their names
indicate, they calculate the center-of-mass, radius of gyration,
moments of inertia, mean-squared displacement, temperature, torque,
and velocity of center-of-mass for each chunk of atoms. The "compute
property/chunk"_compute_property_chunk.html command can tally the
count of atoms in each chunk and extract other per-chunk properties.
The reason these various calculations are not part of the "fix
ave/chunk command"_fix_ave_chunk.html, is that each requires a more
complicated operation than simply summing and averaging over per-atom
values in each chunk. For example, many of them require calculation
of a center of mass, which requires summing mass*position over the
atoms and then dividing by summed mass.
All of these computes produce a global vector or global array as
output, wih one or more values per chunk. The output can be used in
various ways:
As input to the "fix ave/time"_fix_ave_time.html command, which can
write the values to a file and optionally time average them. :ulb,l
As input to the "fix ave/histo"_fix_ave_histo.html command to
histogram values across chunks. E.g. a histogram of cluster sizes or
molecule diffusion rates. :l
As input to special functions of "equal-style
variables"_variable.html, like sum() and max() and ave(). E.g. to
find the largest cluster or fastest diffusing molecule or average
radius-of-gyration of a set of molecules (chunks). :l,ule
Other chunk commands: :h4
"compute chunk/spread/atom"_compute_chunk_spread_atom.html
"compute reduce/chunk"_compute_reduce_chunk.html :ul
The "compute chunk/spread/atom"_compute_chunk_spread_atom.html command
spreads per-chunk values to each atom in the chunk, producing per-atom
values as its output. This can be useful for outputting per-chunk
values to a per-atom "dump file"_dump.html. Or for using an atom's
associated chunk value in an "atom-style variable"_variable.html. Or
as input to the "fix ave/chunk"_fix_ave_chunk.html command to
spatially average per-chunk values calculated by a per-chunk compute.
The "compute reduce/chunk"_compute_reduce_chunk.html command reduces a
peratom value across the atoms in each chunk to produce a value per
chunk. When used with the "compute
chunk/spread/atom"_compute_chunk_spread_atom.html command it can
create peratom values that induce a new set of chunks with a second
"compute chunk/atom"_compute_chunk_atom.html command.
Example calculations with chunks :h4
Here are examples using chunk commands to calculate various
properties:
(1) Average velocity in each of 1000 2d spatial bins:
compute cc1 all chunk/atom bin/2d x 0.0 0.1 y lower 0.01 units reduced
fix 1 all ave/chunk 100 10 1000 cc1 vx vy file tmp.out :pre
(2) Temperature in each spatial bin, after subtracting a flow
velocity:
compute cc1 all chunk/atom bin/2d x 0.0 0.1 y lower 0.1 units reduced
compute vbias all temp/profile 1 0 0 y 10
fix 1 all ave/chunk 100 10 1000 cc1 temp bias vbias file tmp.out :pre
(3) Center of mass of each molecule:
compute cc1 all chunk/atom molecule
compute myChunk all com/chunk cc1
fix 1 all ave/time 100 1 100 c_myChunk\[*\] file tmp.out mode vector :pre
(4) Total force on each molecule and ave/max across all molecules:
compute cc1 all chunk/atom molecule
fix 1 all ave/chunk 1000 1 1000 cc1 fx fy fz file tmp.out
variable xave equal ave(f_1\[2\])
variable xmax equal max(f_1\[2\])
thermo 1000
thermo_style custom step temp v_xave v_xmax :pre
(5) Histogram of cluster sizes:
compute cluster all cluster/atom 1.0
compute cc1 all chunk/atom c_cluster compress yes
compute size all property/chunk cc1 count
fix 1 all ave/histo 100 1 100 0 20 20 c_size mode vector ave running beyond ignore file tmp.histo :pre
(6) An example for using a per-chunk value to apply per-atom forces to
compress individual polymer chains (molecules) in a mixture, is
explained on the "compute
chunk/spread/atom"_compute_chunk_spread_atom.html command doc page.
(7) An example for using one set of per-chunk values for molecule
chunks, to create a 2nd set of micelle-scale chunks (clustered
molecules, due to hydrophobicity), is explained on the "compute
chunk/reduce"_compute_reduce_chunk.html command doc page.
(8) An example for using one set of per-chunk values (dipole moment
vectors) for molecule chunks, spreading the values to each atom in
each chunk, then defining a second set of chunks as spatial bins, and
using the "fix ave/chunk"_fix_ave_chunk.html command to calculate an
average dipole moment vector for each bin. This example is explained
on the "compute chunk/spread/atom"_compute_chunk_spread_atom.html
command doc page.