Merge branch 'develop' into collected-small-changes

This commit is contained in:
Axel Kohlmeyer
2025-01-09 11:17:05 -05:00
9 changed files with 1094 additions and 4 deletions

View File

@ -146,6 +146,8 @@ Lowercase directories
+-------------+------------------------------------------------------------------+
| streitz | use of Streitz/Mintmire potential with charge equilibration |
+-------------+------------------------------------------------------------------+
| stress_vcm | removing binned rigid body motion from binned stress profile |
+-------------+------------------------------------------------------------------+
| tad | temperature-accelerated dynamics of vacancy diffusion in bulk Si |
+-------------+------------------------------------------------------------------+
| threebody | regression test input for a variety of manybody potentials |

View File

@ -184,11 +184,24 @@ temp/chunk calculation to a file is to use the
The keyword/value option pairs are used in the following ways.
The *com* keyword can be used with a value of *yes* to subtract the
velocity of the center-of-mass for each chunk from the velocity of the
atoms in that chunk, before calculating either the global or per-chunk
temperature. This can be useful if the atoms are streaming or
velocity of the center-of-mass (VCM) for each chunk from the velocity of
the atoms in that chunk, before calculating either the global or per-chunk
temperature. This can be useful if the atoms are streaming or
otherwise moving collectively, and you wish to calculate only the
thermal temperature.
thermal temperature. This per-chunk VCM bias can be used in other fixes and
computes that can incorporate a temperature bias. If this compute is used
as a temperature bias in other commands then this bias is subtracted from
each atom, the command runs with the remaining thermal velocities, and
then the bias is added back in. This includes thermostatting
fixes like :doc:`fix nvt <fix_nh>`,
:doc:`fix temp/rescale <fix_temp_rescale>`,
:doc:`fix temp/berendsen <fix_temp_berendsen>`, and
:doc:`fix langevin <fix_langevin>`, and computes like
:doc:`compute stress/atom <compute_stress_atom>` and
:doc:`compute pressure <compute_pressure>`. See the input script in
examples/stress_vcm for an example of how to use the *com* keyword in
conjunction with compute stress/atom to create a stress profile of a rigid
body while removing the overall motion of the rigid body.
For the *bias* keyword, *bias-ID* refers to the ID of a temperature
compute that removes a "bias" velocity from each atom. This also

View File

@ -112,6 +112,7 @@ snap: examples for using several bundled SNAP potentials
srd: stochastic rotation dynamics (SRD) particles as solvent
steinhardt: Steinhardt-Nelson Q_l and W_l parameters usng orientorder/atom
streitz: Streitz-Mintmire potential for Al2O3
stress_vcm: removing binned rigid body motion from binned stress profile
tad: temperature-accelerated dynamics of vacancy diffusion in bulk Si
template: examples for using atom_style template and comparing to atom style molecular
tersoff: regression test input for Tersoff variants

View File

@ -0,0 +1,32 @@
README stress_vcm
=================
Contents:
- in.stress_vcm: Example script showing how to remove binned
velocities of center of mass (VCM) from stress calculations.
- stress_comparison.19Nov24.png: Plot shows the stress
calculated in bars on the y axis for each positional bin on
the x axis. Plotted are three different time steps from
stress profiles with and without the VCM removed. Plot
generated using Python.
- stress_xx.19Nov24.out: Output file generated by fix ave/time.
- log.19Nov24.stress_vcm.g++.1: LAMMPS log file with 1 proc.
- log.19Nov24.stress_vcm.g++.4: LAMMPS log file with 4 procs.
Notes:
- Running this script as-is will generate two files. A log
file with thermodynamic data and a stress_xx.out file
containing the binned stress profile with the VCM removed.
- To generate the binned stress profile without removing the
VCM then the compute stress/atom command at step three
needs the last keyword "ch_temp_vcm" to be replaced with
"NULL".
- Uncommenting the line under "Atom dump" will generate an
all atom dump file every 50 time steps containing atom ID,
type, and xyz coordinates.
- Uncommenting the lines under "Image dumps" will generate
.jpg image files every 250 timesteps.
- Uncommenting lines under "Movie dump" will generate a .avi
movie file showing timesteps every 125 timesteps.

View File

@ -0,0 +1,113 @@
# Removing Binned Velocities of Center of Mass (VCM) from Stress
# This example shows how to remove rigid body motion from
# binned stress calculations. This uses a combination of commands
# from compute chunk/atom, compute temp/chunk, compute
# stress/atom and fix ave/time. We'll show how these commands
# work in the context of a shockwave experiment on a cube of
# atoms. To shock the cube, a rectangular region of atoms is
# frozen, moved into the cube with a constant velocity along the
# x direction, and then unfrozen. As the shockwave begins
# propagating, the body of the cube also moves along the x
# direction. To better understand the stress dynamics of the
# cube we remove the velocity component belonging to the overall
# motion of each bin.
units metal
boundary p p p
atom_style atomic
lattice fcc 5.3589
processors 1 * *
# Defining regions for box and atoms.
# In this experiment an elongated simulation cell is
# defined in the x direction to allow for non-periodic
# motion of the atoms.
region box1 block -3 24 0 12 0 12 units lattice
region box2 block 0 12 0 12 0 12 units lattice
# Creating box and atoms
create_box 1 box1
create_atoms 1 region box2
mass 1 40.00
# Adding energy to the system
velocity all create 600.0 9999
pair_style lj/cut 10
pair_coeff 1 1 0.04 3.405
# Begin time integration
timestep 2e-3
fix fix_nve all nve
thermo 100
run 500
#--------------------------------------#
# Chunk, Stress, and VCM removal steps #
#--------------------------------------#
# 1. Create 20 equispaced bins sliced along the x direction.
# -"units reduced" normalizes the distance from 0.0 to 1.0
variable nbins index 20
variable fraction equal 1.0/v_nbins
variable volfrac equal 1/(vol*${fraction})
compute ch_id all chunk/atom bin/1d x lower ${fraction} units reduced
# 2. Calculate temperature bins with VCM aka COM velocities removed.
compute ch_temp_vcm all temp/chunk ch_id com yes
# 3. Compute per atom stress with VCM removed via temp-ID.
# -The velocities from specified temp-ID are used to compute stress.
# -Stress/atom units are pressure*volume! Optionally handled next step.
compute atom_stress_vcm all stress/atom ch_temp_vcm
# 4. Divide out bin volume from xx stress component.
variable stress atom -(c_atom_stress_vcm[1])/(vol*${fraction})
# 5. Sum the per atom stresses in each bin.
compute ch_stress_vcm all reduce/chunk ch_id sum v_stress
# 6. Average and output to file.
# -The average output is every 100 steps with samples collected 20 times with 5 step intervals.
fix ave_stress_vcm all ave/time 5 20 100 c_ch_stress_vcm mode vector file stress_xx.out
#--------------------------------------#
# Piston compressing along x direction
region piston block -1 1 INF INF INF INF units lattice
group piston region piston
fix fix_piston piston move linear 5 0 0 units box # strain rate ~ 8e10 1/s
thermo_style custom step temp ke pe lx ly lz pxx pyy pzz econserve
# Atom dump
# dump atom_dump all atom 50 dump.vcm
# # Image dumps
# dump 2 all image 250 image.*.jpg type type &
# axes yes 0.8 0.02 view 60 -30
# dump_modify 2 pad 1
# # Movie dump
# dump 3 all movie 125 movie.avi type type &
# axes yes 0.8 0.02 view 60 -30
# dump_modify 3 pad 1
run 500
unfix fix_piston
run 1500

View File

@ -0,0 +1,253 @@
LAMMPS (19 Nov 2024)
OMP_NUM_THREADS environment is not set. Defaulting to 1 thread. (src/comm.cpp:99)
using 1 OpenMP thread(s) per MPI task
# Removing Binned Velocities of Center of Mass (VCM) from Stress
# This example shows how to remove rigid body motion from
# binned stress calculations. This uses a combination of commands
# from compute chunk/atom, compute temp/chunk, compute
# stress/atom and fix ave/time. We'll show how these commands
# work in the context of a shockwave experiment on a cube of
# atoms. To shock the cube, a rectangular region of atoms is
# frozen, moved into the cube with a constant velocity along the
# x direction, and then unfrozen. As the shockwave begins
# propagating, the body of the cube also moves along the x
# direction. To better understand the stress dynamics of the
# cube we remove the velocity component belonging to the overall
# motion of each bin.
units metal
boundary p p p
atom_style atomic
lattice fcc 5.3589
Lattice spacing in x,y,z = 5.3589 5.3589 5.3589
processors 1 * *
# Defining regions for box and atoms.
# In this experiment an elongated simulation cell is
# defined in the x direction to allow for non-periodic
# motion of the atoms.
region box1 block -3 24 0 12 0 12 units lattice
region box2 block 0 12 0 12 0 12 units lattice
# Creating box and atoms
create_box 1 box1
Created orthogonal box = (-16.0767 0 0) to (128.6136 64.3068 64.3068)
1 by 1 by 1 MPI processor grid
create_atoms 1 region box2
Created 7200 atoms
using lattice units in orthogonal box = (-16.0767 0 0) to (128.6136 64.3068 64.3068)
create_atoms CPU = 0.002 seconds
mass 1 40.00
# Adding energy to the system
velocity all create 600.0 9999
pair_style lj/cut 10
pair_coeff 1 1 0.04 3.405
# Begin time integration
timestep 2e-3
fix fix_nve all nve
thermo 100
run 500
Generated 0 of 0 mixed pair_coeff terms from geometric mixing rule
Neighbor list info ...
update: every = 1 steps, delay = 0 steps, check = yes
max neighbors/atom: 2000, page size: 100000
master list distance cutoff = 12
ghost atom cutoff = 12
binsize = 6, bins = 25 11 11
1 neighbor lists, perpetual/occasional/extra = 1 0 0
(1) pair lj/cut, perpetual
attributes: half, newton on
pair build: half/bin/atomonly/newton
stencil: half/bin/3d
bin: standard
Per MPI rank memory allocation (min/avg/max) = 5.721 | 5.721 | 5.721 Mbytes
Step Temp E_pair E_mol TotEng Press
0 600 -2252.7567 0 -1694.4304 -974.62456
100 284.72172 -1977.4291 0 -1712.483 2453.7429
200 304.44519 -1994.7937 0 -1711.4941 1822.2699
300 304.28012 -1993.2958 0 -1710.1498 1498.3794
400 296.76492 -1985.1364 0 -1708.9836 1259.9474
500 295.00895 -1982.4224 0 -1707.9036 964.9526
Loop time of 3.01696 on 1 procs for 500 steps with 7200 atoms
Performance: 28.638 ns/day, 0.838 hours/ns, 165.730 timesteps/s, 1.193 Matom-step/s
99.4% CPU use with 1 MPI tasks x 1 OpenMP threads
MPI task timing breakdown:
Section | min time | avg time | max time |%varavg| %total
---------------------------------------------------------------
Pair | 2.8439 | 2.8439 | 2.8439 | 0.0 | 94.26
Neigh | 0.11212 | 0.11212 | 0.11212 | 0.0 | 3.72
Comm | 0.015585 | 0.015585 | 0.015585 | 0.0 | 0.52
Output | 0.003747 | 0.003747 | 0.003747 | 0.0 | 0.12
Modify | 0.026097 | 0.026097 | 0.026097 | 0.0 | 0.87
Other | | 0.01551 | | | 0.51
Nlocal: 7200 ave 7200 max 7200 min
Histogram: 1 0 0 0 0 0 0 0 0 0
Nghost: 6410 ave 6410 max 6410 min
Histogram: 1 0 0 0 0 0 0 0 0 0
Neighs: 615095 ave 615095 max 615095 min
Histogram: 1 0 0 0 0 0 0 0 0 0
Total # of neighbors = 615095
Ave neighs/atom = 85.429861
Neighbor list builds = 9
Dangerous builds = 0
#------------------------------------#
# Chunk, Stress, and VCM removal steps
#------------------------------------#
# 1. Create 20 equispaced bins sliced along the x direction.
# "units reduced" normalizes the distance from 0 to 1
variable nbins index 20
variable fraction equal 1.0/v_nbins
variable volfrac equal 1/(vol*${fraction})
variable volfrac equal 1/(vol*0.05)
compute ch_id all chunk/atom bin/1d x lower ${fraction} units reduced
compute ch_id all chunk/atom bin/1d x lower 0.05 units reduced
# 2. Calculate temperature bins with VCM aka COM velocities removed.
compute ch_temp_vcm all temp/chunk ch_id com yes
# 3. Compute per atom stress with VCM removed via temp-ID.
# The velocities from specified temp-ID are used to compute stress
# Stress/atom units are pressure*volume! Optionally handled next step.
compute atom_stress_vcm all stress/atom ch_temp_vcm
# 4. Divide out bin volume from xx stress component.
variable stress atom -(c_atom_stress_vcm[1])/(vol*${fraction})
variable stress atom -(c_atom_stress_vcm[1])/(vol*0.05)
# 5. Sum the per atom stresses in each bin.
compute ch_stress_vcm all reduce/chunk ch_id sum v_stress
# 6. Average and output to file.
# The average output is every 100 steps with samples collected 20 times with 5 step intervals
# fix ave_stress_vcm all ave/time 5 20 100 c_ch_stress_vcm mode vector file stress_xx.out
#------------------------------------#
# Piston compressing along x direction
region piston block -1 1 INF INF INF INF units lattice
group piston region piston
863 atoms in group piston
fix fix_piston piston move linear 5 0 0 units box # strain rate ~ 8e10 1/s
thermo_style custom step temp ke pe lx ly lz pxx pyy pzz econserve
# Atom dump
# dump atom_dump all atom 50 dump.vcm
# # Image dumps
# dump 2 all image 250 image.*.jpg type type # axes yes 0.8 0.02 view 60 -30
# dump_modify 2 pad 1
# # Movie dump
# dump 3 all movie 125 movie.avi type type # axes yes 0.8 0.02 view 60 -30
# dump_modify 3 pad 1
run 500
Generated 0 of 0 mixed pair_coeff terms from geometric mixing rule
WARNING: One or more atoms are time integrated more than once (src/modify.cpp:296)
Per MPI rank memory allocation (min/avg/max) = 6.975 | 6.975 | 6.975 Mbytes
Step Temp KinEng PotEng Lx Ly Lz Pxx Pyy Pzz Econserve
500 295.00895 274.51875 -1982.4224 144.6903 64.3068 64.3068 631.89976 1127.2965 1135.6616 -1707.9036
600 357.38902 332.56613 -1951.3422 144.6903 64.3068 64.3068 2236.6706 2003.2726 1943.6815 -1618.7761
700 420.30268 391.11005 -1911.8178 144.6903 64.3068 64.3068 3761.5011 3065.4699 3140.3169 -1520.7077
800 484.96279 451.27911 -1875.379 144.6903 64.3068 64.3068 5362.254 4174.4201 4166.0818 -1424.0999
900 587.78954 546.96391 -1871.217 144.6903 64.3068 64.3068 6481.4714 4875.705 4676.6083 -1324.2531
1000 684.07997 636.56636 -1868.1639 144.6903 64.3068 64.3068 7734.6158 5271.3524 5272.1276 -1231.5975
Loop time of 3.09383 on 1 procs for 500 steps with 7200 atoms
Performance: 27.927 ns/day, 0.859 hours/ns, 161.612 timesteps/s, 1.164 Matom-step/s
100.0% CPU use with 1 MPI tasks x 1 OpenMP threads
MPI task timing breakdown:
Section | min time | avg time | max time |%varavg| %total
---------------------------------------------------------------
Pair | 2.8485 | 2.8485 | 2.8485 | 0.0 | 92.07
Neigh | 0.18767 | 0.18767 | 0.18767 | 0.0 | 6.07
Comm | 0.011533 | 0.011533 | 0.011533 | 0.0 | 0.37
Output | 0.003323 | 0.003323 | 0.003323 | 0.0 | 0.11
Modify | 0.031777 | 0.031777 | 0.031777 | 0.0 | 1.03
Other | | 0.01107 | | | 0.36
Nlocal: 7200 ave 7200 max 7200 min
Histogram: 1 0 0 0 0 0 0 0 0 0
Nghost: 6409 ave 6409 max 6409 min
Histogram: 1 0 0 0 0 0 0 0 0 0
Neighs: 646408 ave 646408 max 646408 min
Histogram: 1 0 0 0 0 0 0 0 0 0
Total # of neighbors = 646408
Ave neighs/atom = 89.778889
Neighbor list builds = 15
Dangerous builds = 0
unfix fix_piston
run 1500
Generated 0 of 0 mixed pair_coeff terms from geometric mixing rule
Per MPI rank memory allocation (min/avg/max) = 6.6 | 6.6 | 6.6 Mbytes
Step Temp KinEng PotEng Lx Ly Lz Pxx Pyy Pzz Econserve
1000 684.07997 636.56636 -1868.1639 144.6903 64.3068 64.3068 7734.6158 5271.3524 5272.1276 -1231.5975
1100 710.19886 660.87113 -1894.0485 144.6903 64.3068 64.3068 8048.3485 5396.6668 5376.5956 -1233.1774
1200 717.16487 667.35331 -1901.3849 144.6903 64.3068 64.3068 8009.7984 5634.5121 5349.4113 -1234.0316
1300 710.26037 660.92837 -1894.9802 144.6903 64.3068 64.3068 8063.4125 5572.1245 5530.174 -1234.0519
1400 715.93921 666.21278 -1898.8885 144.6903 64.3068 64.3068 7752.0927 5293.5463 5322.2312 -1232.6757
1500 748.85411 696.84154 -1926.4891 144.6903 64.3068 64.3068 6030.5428 4076.8886 4012.7653 -1229.6475
1600 767.98982 714.64815 -1939.8556 144.6903 64.3068 64.3068 4200.3475 2532.5711 2530.5518 -1225.2075
1700 757.22042 704.62675 -1925.553 144.6903 64.3068 64.3068 2686.7843 1482.2796 1505.8073 -1220.9262
1800 727.30327 676.78754 -1894.6635 144.6903 64.3068 64.3068 1764.2793 781.37451 801.18668 -1217.8759
1900 688.82146 640.97853 -1856.5007 144.6903 64.3068 64.3068 1022.805 417.32394 359.74951 -1215.5221
2000 655.91228 610.35509 -1823.954 144.6903 64.3068 64.3068 551.98825 -20.148643 -56.976652 -1213.5989
2100 620.22468 577.14622 -1789.1761 144.6903 64.3068 64.3068 264.05975 -266.8323 -314.45533 -1212.0299
2200 589.13325 548.21428 -1758.9252 144.6903 64.3068 64.3068 41.369707 -533.503 -525.69401 -1210.7109
2300 563.20394 524.08593 -1733.6036 144.6903 64.3068 64.3068 -220.99189 -810.90513 -774.65084 -1209.5176
2400 540.44236 502.90528 -1711.3384 144.6903 64.3068 64.3068 -358.01508 -962.31635 -977.3253 -1208.4332
2500 523.5718 487.20648 -1694.7088 144.6903 64.3068 64.3068 -521.87444 -1152.8386 -1231.7615 -1207.5023
Loop time of 9.34327 on 1 procs for 1500 steps with 7200 atoms
Performance: 27.742 ns/day, 0.865 hours/ns, 160.543 timesteps/s, 1.156 Matom-step/s
98.5% CPU use with 1 MPI tasks x 1 OpenMP threads
MPI task timing breakdown:
Section | min time | avg time | max time |%varavg| %total
---------------------------------------------------------------
Pair | 8.4692 | 8.4692 | 8.4692 | 0.0 | 90.65
Neigh | 0.7512 | 0.7512 | 0.7512 | 0.0 | 8.04
Comm | 0.031189 | 0.031189 | 0.031189 | 0.0 | 0.33
Output | 0.010584 | 0.010584 | 0.010584 | 0.0 | 0.11
Modify | 0.053052 | 0.053052 | 0.053052 | 0.0 | 0.57
Other | | 0.02803 | | | 0.30
Nlocal: 7200 ave 7200 max 7200 min
Histogram: 1 0 0 0 0 0 0 0 0 0
Nghost: 6380 ave 6380 max 6380 min
Histogram: 1 0 0 0 0 0 0 0 0 0
Neighs: 515773 ave 515773 max 515773 min
Histogram: 1 0 0 0 0 0 0 0 0 0
Total # of neighbors = 515773
Ave neighs/atom = 71.635139
Neighbor list builds = 57
Dangerous builds = 0
Total wall time: 0:00:15

View File

@ -0,0 +1,253 @@
LAMMPS (19 Nov 2024)
OMP_NUM_THREADS environment is not set. Defaulting to 1 thread. (src/comm.cpp:99)
using 1 OpenMP thread(s) per MPI task
# Removing Binned Velocities of Center of Mass (VCM) from Stress
# This example shows how to remove rigid body motion from
# binned stress calculations. This uses a combination of commands
# from compute chunk/atom, compute temp/chunk, compute
# stress/atom and fix ave/time. We'll show how these commands
# work in the context of a shockwave experiment on a cube of
# atoms. To shock the cube, a rectangular region of atoms is
# frozen, moved into the cube with a constant velocity along the
# x direction, and then unfrozen. As the shockwave begins
# propagating, the body of the cube also moves along the x
# direction. To better understand the stress dynamics of the
# cube we remove the velocity component belonging to the overall
# motion of each bin.
units metal
boundary p p p
atom_style atomic
lattice fcc 5.3589
Lattice spacing in x,y,z = 5.3589 5.3589 5.3589
processors 1 * *
# Defining regions for box and atoms.
# In this experiment an elongated simulation cell is
# defined in the x direction to allow for non-periodic
# motion of the atoms.
region box1 block -3 24 0 12 0 12 units lattice
region box2 block 0 12 0 12 0 12 units lattice
# Creating box and atoms
create_box 1 box1
Created orthogonal box = (-16.0767 0 0) to (128.6136 64.3068 64.3068)
1 by 2 by 2 MPI processor grid
create_atoms 1 region box2
Created 7200 atoms
using lattice units in orthogonal box = (-16.0767 0 0) to (128.6136 64.3068 64.3068)
create_atoms CPU = 0.001 seconds
mass 1 40.00
# Adding energy to the system
velocity all create 600.0 9999
pair_style lj/cut 10
pair_coeff 1 1 0.04 3.405
# Begin time integration
timestep 2e-3
fix fix_nve all nve
thermo 100
run 500
Generated 0 of 0 mixed pair_coeff terms from geometric mixing rule
Neighbor list info ...
update: every = 1 steps, delay = 0 steps, check = yes
max neighbors/atom: 2000, page size: 100000
master list distance cutoff = 12
ghost atom cutoff = 12
binsize = 6, bins = 25 11 11
1 neighbor lists, perpetual/occasional/extra = 1 0 0
(1) pair lj/cut/opt, perpetual
attributes: half, newton on
pair build: half/bin/atomonly/newton
stencil: half/bin/3d
bin: standard
Per MPI rank memory allocation (min/avg/max) = 3.662 | 3.662 | 3.662 Mbytes
Step Temp E_pair E_mol TotEng Press
0 600 -2252.7567 0 -1694.4304 -974.62456
100 284.1896 -1976.961 0 -1712.5101 2462.6396
200 308.58965 -1998.6349 0 -1711.4787 1789.0033
300 300.55093 -1989.9838 0 -1710.308 1545.8576
400 297.91491 -1986.2519 0 -1709.029 1247.7121
500 294.66041 -1982.1097 0 -1707.9153 961.03073
Loop time of 0.942408 on 4 procs for 500 steps with 7200 atoms
Performance: 91.680 ns/day, 0.262 hours/ns, 530.556 timesteps/s, 3.820 Matom-step/s
82.1% CPU use with 4 MPI tasks x 1 OpenMP threads
MPI task timing breakdown:
Section | min time | avg time | max time |%varavg| %total
---------------------------------------------------------------
Pair | 0.61287 | 0.63781 | 0.65858 | 2.1 | 67.68
Neigh | 0.030246 | 0.031529 | 0.034546 | 1.0 | 3.35
Comm | 0.23074 | 0.25145 | 0.27819 | 3.7 | 26.68
Output | 0.000282 | 0.0003735 | 0.000463 | 0.0 | 0.04
Modify | 0.005566 | 0.0057635 | 0.005989 | 0.2 | 0.61
Other | | 0.01548 | | | 1.64
Nlocal: 1800 ave 1814 max 1787 min
Histogram: 1 0 1 0 0 0 0 1 0 1
Nghost: 3713.5 ave 3727 max 3699 min
Histogram: 1 0 1 0 0 0 0 1 0 1
Neighs: 153532 ave 154995 max 152312 min
Histogram: 1 0 1 0 0 1 0 0 0 1
Total # of neighbors = 614128
Ave neighs/atom = 85.295556
Neighbor list builds = 9
Dangerous builds = 0
#------------------------------------#
# Chunk, Stress, and VCM removal steps
#------------------------------------#
# 1. Create 20 equispaced bins sliced along the x direction.
# "units reduced" normalizes the distance from 0 to 1
variable nbins index 20
variable fraction equal 1.0/v_nbins
variable volfrac equal 1/(vol*${fraction})
variable volfrac equal 1/(vol*0.05)
compute ch_id all chunk/atom bin/1d x lower ${fraction} units reduced
compute ch_id all chunk/atom bin/1d x lower 0.05 units reduced
# 2. Calculate temperature bins with VCM aka COM velocities removed.
compute ch_temp_vcm all temp/chunk ch_id com yes
# 3. Compute per atom stress with VCM removed via temp-ID.
# The velocities from specified temp-ID are used to compute stress
# Stress/atom units are pressure*volume! Optionally handled next step.
compute atom_stress_vcm all stress/atom ch_temp_vcm
# 4. Divide out bin volume from xx stress component.
variable stress atom -(c_atom_stress_vcm[1])/(vol*${fraction})
variable stress atom -(c_atom_stress_vcm[1])/(vol*0.05)
# 5. Sum the per atom stresses in each bin.
compute ch_stress_vcm all reduce/chunk ch_id sum v_stress
# 6. Average and output to file.
# The average output is every 100 steps with samples collected 20 times with 5 step intervals
# fix ave_stress_vcm all ave/time 5 20 100 c_ch_stress_vcm mode vector file stress_xx.out
#------------------------------------#
# Piston compressing along x direction
region piston block -1 1 INF INF INF INF units lattice
group piston region piston
864 atoms in group piston
fix fix_piston piston move linear 5 0 0 units box # strain rate ~ 8e10 1/s
thermo_style custom step temp ke pe lx ly lz pxx pyy pzz econserve
# Atom dump
# dump atom_dump all atom 50 dump.vcm
# # Image dumps
# dump 2 all image 250 image.*.jpg type type # axes yes 0.8 0.02 view 60 -30
# dump_modify 2 pad 1
# # Movie dump
# dump 3 all movie 125 movie.avi type type # axes yes 0.8 0.02 view 60 -30
# dump_modify 3 pad 1
run 500
Generated 0 of 0 mixed pair_coeff terms from geometric mixing rule
WARNING: One or more atoms are time integrated more than once (src/modify.cpp:296)
Per MPI rank memory allocation (min/avg/max) = 4.916 | 4.916 | 4.916 Mbytes
Step Temp KinEng PotEng Lx Ly Lz Pxx Pyy Pzz Econserve
500 294.66041 274.19441 -1982.1097 144.6903 64.3068 64.3068 645.25795 1119.5337 1118.3006 -1707.9153
600 357.88641 333.02897 -1951.8158 144.6903 64.3068 64.3068 2176.0343 1929.2787 1981.8479 -1618.7869
700 418.41159 389.3503 -1912.8337 144.6903 64.3068 64.3068 3702.2875 3043.7607 3081.1607 -1523.4834
800 483.71102 450.11428 -1875.7955 144.6903 64.3068 64.3068 5254.3875 4190.9789 4158.3561 -1425.6813
900 586.0893 545.38176 -1870.9313 144.6903 64.3068 64.3068 6509.1439 4756.2216 4724.7086 -1325.5495
1000 686.32946 638.65962 -1874.811 144.6903 64.3068 64.3068 7515.1606 5193.049 5261.8688 -1236.1514
Loop time of 0.656417 on 4 procs for 500 steps with 7200 atoms
Performance: 131.624 ns/day, 0.182 hours/ns, 761.711 timesteps/s, 5.484 Matom-step/s
92.8% CPU use with 4 MPI tasks x 1 OpenMP threads
MPI task timing breakdown:
Section | min time | avg time | max time |%varavg| %total
---------------------------------------------------------------
Pair | 0.51672 | 0.52334 | 0.53259 | 0.8 | 79.73
Neigh | 0.045091 | 0.045915 | 0.047402 | 0.4 | 6.99
Comm | 0.060735 | 0.071794 | 0.079302 | 2.6 | 10.94
Output | 0.000208 | 0.000389 | 0.000926 | 0.0 | 0.06
Modify | 0.006007 | 0.0061595 | 0.00626 | 0.1 | 0.94
Other | | 0.008815 | | | 1.34
Nlocal: 1800 ave 1811 max 1785 min
Histogram: 1 0 0 1 0 0 0 0 0 2
Nghost: 3713.25 ave 3727 max 3702 min
Histogram: 2 0 0 0 0 0 0 1 0 1
Neighs: 161477 ave 162958 max 159732 min
Histogram: 1 0 0 0 1 0 0 1 0 1
Total # of neighbors = 645909
Ave neighs/atom = 89.709583
Neighbor list builds = 15
Dangerous builds = 0
unfix fix_piston
run 1500
Generated 0 of 0 mixed pair_coeff terms from geometric mixing rule
Per MPI rank memory allocation (min/avg/max) = 4.541 | 4.541 | 4.541 Mbytes
Step Temp KinEng PotEng Lx Ly Lz Pxx Pyy Pzz Econserve
1000 686.32946 638.65962 -1874.811 144.6903 64.3068 64.3068 7515.1606 5193.049 5261.8688 -1236.1514
1100 709.7333 660.43791 -1898.2844 144.6903 64.3068 64.3068 7932.8638 5334.6171 5364.5335 -1237.8465
1200 713.27253 663.73132 -1902.4588 144.6903 64.3068 64.3068 7957.2574 5500.6231 5538.0516 -1238.7275
1300 705.44796 656.45022 -1895.1575 144.6903 64.3068 64.3068 7996.7265 5584.6233 5538.2494 -1238.7072
1400 711.86463 662.42121 -1899.8416 144.6903 64.3068 64.3068 7674.2462 5292.4915 5294.5366 -1237.4204
1500 742.18946 690.63979 -1924.9562 144.6903 64.3068 64.3068 6047.915 4056.6156 4014.4446 -1234.3164
1600 762.81764 709.83522 -1939.8563 144.6903 64.3068 64.3068 4185.5873 2530.0572 2576.1943 -1230.0211
1700 754.40428 702.00621 -1927.7337 144.6903 64.3068 64.3068 2662.7604 1509.1985 1484.7252 -1225.7275
1800 721.03504 670.95468 -1893.5556 144.6903 64.3068 64.3068 1765.8783 835.89765 861.9432 -1222.6009
1900 689.64162 641.74172 -1861.8886 144.6903 64.3068 64.3068 941.58148 312.93205 409.79901 -1220.1469
2000 650.79664 605.59477 -1823.9889 144.6903 64.3068 64.3068 543.39234 28.48735 80.396505 -1218.3941
2100 616.04072 573.25286 -1790.1764 144.6903 64.3068 64.3068 308.16444 -235.20997 -248.22531 -1216.9235
2200 587.18712 546.40333 -1761.8878 144.6903 64.3068 64.3068 37.044801 -476.50396 -470.83059 -1215.4845
2300 562.84178 523.74892 -1738.2239 144.6903 64.3068 64.3068 -139.28348 -711.17273 -730.80877 -1214.475
2400 540.48362 502.94367 -1716.3529 144.6903 64.3068 64.3068 -320.98222 -951.2066 -943.93966 -1213.4093
2500 519.80431 483.70067 -1696.1896 144.6903 64.3068 64.3068 -471.61317 -1088.8457 -1131.5396 -1212.4889
Loop time of 1.97213 on 4 procs for 1500 steps with 7200 atoms
Performance: 131.431 ns/day, 0.183 hours/ns, 760.598 timesteps/s, 5.476 Matom-step/s
95.3% CPU use with 4 MPI tasks x 1 OpenMP threads
MPI task timing breakdown:
Section | min time | avg time | max time |%varavg| %total
---------------------------------------------------------------
Pair | 1.5455 | 1.5599 | 1.5723 | 0.8 | 79.10
Neigh | 0.16844 | 0.1704 | 0.17237 | 0.4 | 8.64
Comm | 0.19002 | 0.2047 | 0.22068 | 2.4 | 10.38
Output | 0.000525 | 0.0006785 | 0.001077 | 0.0 | 0.03
Modify | 0.012434 | 0.012601 | 0.012777 | 0.1 | 0.64
Other | | 0.02388 | | | 1.21
Nlocal: 1800 ave 1833 max 1776 min
Histogram: 1 0 1 0 1 0 0 0 0 1
Nghost: 3702 ave 3732 max 3674 min
Histogram: 1 0 0 1 0 0 1 0 0 1
Neighs: 129380 ave 132578 max 127003 min
Histogram: 1 0 0 2 0 0 0 0 0 1
Total # of neighbors = 517520
Ave neighs/atom = 71.877778
Neighbor list builds = 54
Dangerous builds = 0
Total wall time: 0:00:03

Binary file not shown.

After

Width:  |  Height:  |  Size: 89 KiB

View File

@ -0,0 +1,423 @@
# Time-averaged data for fix ave_stress_vcm
# TimeStep Number-of-rows
# Row c_ch_stress_vcm
600 20
1 0
2 -142.965
3 2142.79
4 12968.3
5 -336.7
6 2638.09
7 4214.83
8 3187.61
9 -488.891
10 -49.3553
11 151.373
12 -317.663
13 0
14 0
15 0
16 0
17 0
18 0
19 0
20 0
700 20
1 0
2 -14.3195
3 -1238.9
4 30664.3
5 18805.2
6 498.562
7 930.874
8 660.655
9 -266.903
10 -317.877
11 -386.989
12 -304.697
13 0
14 0
15 0
16 0
17 0
18 0
19 0
20 0
800 20
1 0
2 0
3 -1656.7
4 30424.3
5 37003.5
6 15562.5
7 -2441.9
8 -1766.09
9 272.718
10 -664.774
11 -72.6933
12 -469.765
13 0
14 0
15 0
16 0
17 0
18 0
19 0
20 0
900 20
1 0
2 0
3 -1567.21
4 24987.6
5 38068.9
6 31595
7 8864.94
8 -3423.99
9 -753.063
10 125.21
11 -50.4895
12 -172.14
13 0
14 0
15 0
16 0
17 0
18 0
19 0
20 0
1000 20
1 0
2 0
3 -893.168
4 15591.6
5 32690.6
6 30183
7 27172
8 9459.75
9 -1416.35
10 -432.731
11 444.323
12 -424.357
13 0
14 0
15 0
16 0
17 0
18 0
19 0
20 0
1100 20
1 0
2 0
3 -601.805
4 8890.79
5 23345.1
6 28529.2
7 29111.9
8 25846.2
9 7451.83
10 -1624.2
11 320.704
12 -50.9865
13 -5.50481
14 0
15 0
16 0
17 0
18 0
19 0
20 0
1200 20
1 0
2 0
3 1435.39
4 8818.29
5 7129.61
6 20281.7
7 28026.1
8 28327.7
9 26918.6
10 8277.12
11 -249.644
12 -171.806
13 -7.19065
14 0
15 0
16 0
17 0
18 0
19 0
20 0
1300 20
1 0
2 0
3 -718.118
4 3021.9
5 9010.51
6 9500.87
7 19432.8
8 27254.3
9 28638.5
10 25568.5
11 8094.66
12 -368.293
13 -2.20997
14 0
15 0
16 0
17 0
18 0
19 0
20 0
1400 20
1 0
2 0
3 -650.581
4 190.19
5 5465.38
6 7489.23
7 7575.16
8 18433.5
9 26975.3
10 28981.5
11 26987.9
12 7502.07
13 0.117312
14 0
15 0
16 0
17 0
18 0
19 0
20 0
1500 20
1 0
2 0
3 -619.311
4 561.257
5 461.5
6 4105.68
7 9272.68
8 10445.6
9 18826.1
10 25434.8
11 25653.8
12 10981.2
13 33.682
14 0
15 0
16 0
17 0
18 0
19 0
20 0
1600 20
1 0
2 0
3 -349.345
4 513.579
5 -471.384
6 1257.81
7 7122.9
8 8659.35
9 8452.08
10 16013.5
11 17091
12 5476.24
13 -136.183
14 0
15 0
16 0
17 0
18 0
19 0
20 0
1700 20
1 0
2 0
3 -273.839
4 -907.407
5 -272.136
6 594.363
7 3302.77
8 5564.07
9 8689.92
10 6446.06
11 1779.37
12 338.998
13 -171.408
14 -1.21548
15 0
16 0
17 0
18 0
19 0
20 0
1800 20
1 0
2 0
3 -164.819
4 383.877
5 -140.681
6 -10.0153
7 907.937
8 3269.05
9 5325.22
10 395.73
11 -4103.73
12 -2787.16
13 -1357.04
14 -35.2044
15 0
16 0
17 0
18 0
19 0
20 0
1900 20
1 0
2 0
3 -80.813
4 334.225
5 248.55
6 82.0566
7 207.763
8 185.714
9 -55.8635
10 -2758.51
11 -4619.33
12 -5521.92
13 -2346.36
14 -415.324
15 0
16 0
17 0
18 0
19 0
20 0
2000 20
1 0
2 0
3 -83.1832
4 264.023
5 596.087
6 40.8157
7 -267.093
8 -2288.15
9 -3387.64
10 -5566.79
11 -5640.76
12 -4925.74
13 -3096.01
14 -757.817
15 -1.13042
16 0
17 0
18 0
19 0
20 0
2100 20
1 0
2 0
3 -17.4378
4 62.1251
5 740.988
6 357.467
7 -1137.61
8 -4266.83
9 -4962.9
10 -5322.45
11 -5437.58
12 -4846.56
13 -3651.28
14 -1151.01
15 -28.3074
16 0
17 0
18 0
19 0
20 0
2200 20
1 0
2 0
3 -10.8779
4 -56.7926
5 400.261
6 -568.63
7 -2193.36
8 -3856.71
9 -6603
10 -5717.11
11 -4868.64
12 -4173.5
13 -3402.64
14 -1712.44
15 -80.6771
16 -0.123189
17 0
18 0
19 0
20 0
2300 20
1 0
2 0
3 -22.8402
4 -44.5496
5 -365.476
6 -1285.6
7 -2887.76
8 -4022.77
9 -6280.86
10 -6055.26
11 -4921.51
12 -4445.37
13 -3531.69
14 -1360.49
15 -258.99
16 0.196931
17 0
18 0
19 0
20 0
2400 20
1 0
2 0
3 -0.594396
4 -148.921
5 -1118.18
6 -2071.85
7 -3989.41
8 -4567.01
9 -4939.36
10 -5170.94
11 -4922.25
12 -4587.5
13 -3748.19
14 -1785.46
15 -460.491
16 2.54038
17 0
18 0
19 0
20 0
2500 20
1 0
2 0
3 5.64755
4 -485.854
5 -2525.68
6 -2642.35
7 -5066.15
8 -4546.03
9 -4429.45
10 -4579.15
11 -4829.56
12 -4384.77
13 -3525.99
14 -1708.9
15 -627.176
16 -23.5581
17 0
18 0
19 0
20 0