diff --git a/doc/src/fix_qeq.rst b/doc/src/fix_qeq.rst index c655076ce8..7ea79a89ce 100644 --- a/doc/src/fix_qeq.rst +++ b/doc/src/fix_qeq.rst @@ -34,13 +34,14 @@ Syntax * maxiter = maximum iterations to perform charge equilibration * qfile = a filename with QEq parameters or *coul/streitz* or *reax/c* * zero or more keyword/value pairs may be appended -* keyword = *alpha* or *qdamp* or *qstep* +* keyword = *alpha* or *qdamp* or *qstep* or *warn* .. parsed-literal:: *alpha* value = Slater type orbital exponent (qeq/slater only) *qdamp* value = damping factor for damped dynamics charge solver (qeq/dynamic and qeq/fire only) *qstep* value = time step size for damped dynamics charge solver (qeq/dynamic and qeq/fire only) + *warn* value = do (=yes) or do not (=no) print a warning when the maximum number of iterations is reached Examples """""""" @@ -101,8 +102,8 @@ The QEq method minimizes the electrostatic energy of the system (or equalizes the derivative of energy with respect to charge of all the atoms) by adjusting the partial charge on individual atoms based on interactions with their neighbors within *cutoff*\ . It requires a few -parameters, in *metal* units, for each atom type which provided in a -file specified by *qfile*\ . The file has the following format +parameters in the appropriate units for each atom type which are read +from a file specified by *qfile*\ . The file has the following format .. parsed-literal:: @@ -112,7 +113,7 @@ file specified by *qfile*\ . The file has the following format Ntype chi eta gamma zeta qcore There have to be parameters given for every atom type. Wildcard entries -are possible using the same syntax as elsewhere in LAMMPS +are possible using the same type range syntax as for "coeff" commands (i.e., n\*m, n\*, \*m, \*). Later entries will overwrite previous ones. Empty lines or any text following the pound sign (#) are ignored. Each line starts with the atom type followed by five parameters. @@ -126,6 +127,14 @@ entries per line are required. * *zeta* = Slater type orbital exponent defined by the :ref:`Streitz-Mintmire ` potential in reverse distance units * *qcore* = charge of the nucleus defined by the :ref:`Streitz-Mintmire potential ` potential in charge units +The fix qeq styles will print a warning if the charges are not +equilibrated within *tolerance* by *maxiter* steps, unless the +*warn* keyword is used with "no" as argument. This latter option +may be useful for testing and benchmarking purposes, as it allows +to use a fixed number of QEq iterations when *tolerance* is set +to a small enough value to alway reach the *maxiter* limit. Turning +off warnings will avoid the excessive output in that case. + The *qeq/point* style describes partial charges on atoms as point charges. Interaction between a pair of charged particles is 1/r, which is the simplest description of the interaction between charges. @@ -229,7 +238,7 @@ Related commands Default """"""" -none +warn yes ---------- diff --git a/examples/qeq/buck.inc b/examples/qeq/buck.inc new file mode 100644 index 0000000000..af3fe3dfb6 --- /dev/null +++ b/examples/qeq/buck.inc @@ -0,0 +1,22 @@ + +kspace_style pppm 1e-6 + +neighbor 1.0 bin +neigh_modify delay 0 every 1 check yes + +group type1 type 1 +compute charge1 type1 property/atom q +compute q1 type1 reduce ave c_charge1 +group type2 type 2 +compute charge2 type2 property/atom q +compute q2 type2 reduce ave c_charge2 +variable qtot equal count(type1)*c_q1+count(type2)*c_q2 +variable nqeq equal f_2 + +thermo_style custom step pe c_q1 c_q2 v_qtot v_nqeq +thermo 10 +thermo_modify format line "%4d %12.9g %12.8f %12.8f %16.12f %6.0f" + +timestep 0.0001 + +fix 1 all nve diff --git a/examples/qeq/in.qeq.buck b/examples/qeq/in.qeq.buck index 2258a9e2ae..4b74b44186 100644 --- a/examples/qeq/in.qeq.buck +++ b/examples/qeq/in.qeq.buck @@ -1,41 +1,68 @@ # This example demonstrates the use of various fix qeq variants with -# that defines and uses charges, in this case pair_style buck/coul/long +# a pair style using charges, in this case pair_style buck/coul/long units metal atom_style charge read_data data.aC -replicate 2 2 2 +#replicate 2 2 2 pair_style buck/coul/long 12.0 pair_coeff 2 2 1388.77 .3623188 175.0 pair_coeff 1 2 18003 .2052124 133.5381 pair_coeff 1 1 0 .1 0 -kspace_style ewald 1e-6 -neighbor 1.0 bin -neigh_modify delay 0 every 1 check yes +fix 2 all qeq/shielded 1 10 1.0e-20 10 param.qeq2 -group type1 type 1 -compute charge1 type1 property/atom q -compute q1 type1 reduce ave c_charge1 -group type2 type 2 -compute charge2 type2 property/atom q -compute q2 type2 reduce ave c_charge2 -variable qtot equal count(type1)*c_q1+count(type2)*c_q2 - -thermo_style custom step pe c_q1 c_q2 v_qtot -thermo 10 - -timestep 0.0001 +include buck.inc velocity all create 300.0 1281937 -fix 1 all nve +run 0 post no -#fix 2 all qeq/point 1 10 1.0e-6 100 param.qeq2 -#fix 2 all qeq/shielded 1 10 1.0e-6 100 param.qeq2 -#fix 2 all qeq/slater 1 10 1.0e-6 100 param.qeq2 -#fix 2 all qeq/dynamic 1 10 1.0e-4 100 param.qeq2 -fix 2 all qeq/fire 1 10 1.0e-4 100 param.qeq2 +write_restart qeq.restart -run 100 +clear + +print "Using fix qeq/point" +read_restart qeq.restart +fix 2 all qeq/point 1 10 1.0e-6 100 param.qeq2 +include buck.inc + +run 100 + +clear + +print "Using fix qeq/shielded" +read_restart qeq.restart +fix 2 all qeq/shielded 1 10 1.0e-6 100 param.qeq2 +include buck.inc + +run 100 + + +clear + +print "Using fix qeq/slater" +read_restart qeq.restart +fix 2 all qeq/shielded 1 10 1.0e-6 100 param.qeq2 +include buck.inc + +run 100 + +clear + +print "Using fix qeq/dynamic" +read_restart qeq.restart +fix 2 all qeq/dynamic 1 10 1.0e-3 100 param.qeq2 +include buck.inc + +run 100 + +clear + +print "Using fix qeq/fire" +read_restart qeq.restart +fix 2 all qeq/fire 1 10 1.0e-3 100 param.qeq2 +include buck.inc + +run 100 diff --git a/examples/qeq/log.20Apr21.qeq.buck.g++.1 b/examples/qeq/log.20Apr21.qeq.buck.g++.1 new file mode 100644 index 0000000000..d5b41867b1 --- /dev/null +++ b/examples/qeq/log.20Apr21.qeq.buck.g++.1 @@ -0,0 +1,650 @@ +LAMMPS (8 Apr 2021) + using 1 OpenMP thread(s) per MPI task +# This example demonstrates the use of various fix qeq variants with +# a pair style using charges, in this case pair_style buck/coul/long + +units metal +atom_style charge + +read_data data.aC +Reading data file ... + orthogonal box = (0.0000000 0.0000000 0.0000000) to (25.158320 25.158320 28.020256) + 1 by 1 by 1 MPI processor grid + reading atoms ... + 1200 atoms + read_data CPU = 0.009 seconds +#replicate 2 2 2 + +pair_style buck/coul/long 12.0 +pair_coeff 2 2 1388.77 .3623188 175.0 +pair_coeff 1 2 18003 .2052124 133.5381 +pair_coeff 1 1 0 .1 0 + +fix 2 all qeq/shielded 1 10 1.0e-20 10 param.qeq2 + +include buck.inc + +kspace_style pppm 1e-6 + +neighbor 1.0 bin +neigh_modify delay 0 every 1 check yes + +group type1 type 1 +400 atoms in group type1 +compute charge1 type1 property/atom q +compute q1 type1 reduce ave c_charge1 +group type2 type 2 +800 atoms in group type2 +compute charge2 type2 property/atom q +compute q2 type2 reduce ave c_charge2 +variable qtot equal count(type1)*c_q1+count(type2)*c_q2 +variable nqeq equal f_2 + +thermo_style custom step pe c_q1 c_q2 v_qtot v_nqeq +thermo 10 +thermo_modify format line "%4d %12.9g %12.8f %12.8f %16.12f %6.0f" + +timestep 0.0001 + +fix 1 all nve + +velocity all create 300.0 1281937 +run 0 post no +PPPM initialization ... + using 12-bit tables for long-range coulomb (src/kspace.cpp:339) + G vector (1/distance) = 0.30705229 + grid = 48 48 54 + stencil order = 5 + estimated absolute RMS force accuracy = 1.8909403e-05 + estimated relative force accuracy = 1.3131854e-06 + using double precision FFTW3 + 3d grid and FFT values/proc = 184525 124416 +Neighbor list info ... + update every 1 steps, delay 0 steps, check yes + max neighbors/atom: 2000, page size: 100000 + master list distance cutoff = 13 + ghost atom cutoff = 13 + binsize = 6.5, bins = 4 4 5 + 2 neighbor lists, perpetual/occasional/extra = 2 0 0 + (1) pair buck/coul/long, perpetual, half/full from (2) + attributes: half, newton on + pair build: halffull/newton + stencil: none + bin: none + (2) fix qeq/shielded, perpetual + attributes: full, newton on + pair build: full/bin/atomonly + stencil: full/bin/3d + bin: standard +WARNING: Fix qeq CG convergence failed (4.357840257025601e-19) after 10 iterations at step 0 (src/QEQ/fix_qeq.cpp:410) +WARNING: Fix qeq CG convergence failed (5.274094378414531e-18) after 10 iterations at step 0 (src/QEQ/fix_qeq.cpp:410) +Per MPI rank memory allocation (min/avg/max) = 38.75 | 38.75 | 38.75 Mbytes +Step PotEng c_q1 c_q2 v_qtot v_nqeq + 0 -2879.00327 0.76536977 -0.38268489 0.000000000000 10 +Loop time of 1.66893e-06 on 1 procs for 0 steps with 1200 atoms + + +write_restart qeq.restart +System init for write_restart ... +PPPM initialization ... + using 12-bit tables for long-range coulomb (src/kspace.cpp:339) + G vector (1/distance) = 0.27644401 + grid = 27 27 30 + stencil order = 5 + estimated absolute RMS force accuracy = 1.4502702e-05 + estimated relative force accuracy = 1.0071569e-06 + using double precision FFTW3 + 3d grid and FFT values/proc = 42772 21870 + +clear + using 1 OpenMP thread(s) per MPI task + +print "Using fix qeq/point" +Using fix qeq/point +read_restart qeq.restart +Reading restart file ... + restart file = 8 Apr 2021, LAMMPS = 8 Apr 2021 + restoring atom style charge from restart + orthogonal box = (0.0000000 0.0000000 0.0000000) to (25.158320 25.158320 28.020256) + 1 by 1 by 1 MPI processor grid + restoring pair style buck/coul/long from restart + 1200 atoms + read_restart CPU = 0.001 seconds +fix 2 all qeq/point 1 10 1.0e-6 100 param.qeq2 +include buck.inc + +kspace_style pppm 1e-6 + +neighbor 1.0 bin +neigh_modify delay 0 every 1 check yes + +group type1 type 1 +400 atoms in group type1 +compute charge1 type1 property/atom q +compute q1 type1 reduce ave c_charge1 +group type2 type 2 +800 atoms in group type2 +compute charge2 type2 property/atom q +compute q2 type2 reduce ave c_charge2 +variable qtot equal count(type1)*c_q1+count(type2)*c_q2 +variable nqeq equal f_2 + +thermo_style custom step pe c_q1 c_q2 v_qtot v_nqeq +thermo 10 +thermo_modify format line "%4d %12.9g %12.8f %12.8f %16.12f %6.0f" + +timestep 0.0001 + +fix 1 all nve + +run 100 +PPPM initialization ... + using 12-bit tables for long-range coulomb (src/kspace.cpp:339) + G vector (1/distance) = 0.27644401 + grid = 27 27 30 + stencil order = 5 + estimated absolute RMS force accuracy = 1.4502702e-05 + estimated relative force accuracy = 1.0071569e-06 + using double precision FFTW3 + 3d grid and FFT values/proc = 42772 21870 +Neighbor list info ... + update every 1 steps, delay 0 steps, check yes + max neighbors/atom: 2000, page size: 100000 + master list distance cutoff = 13 + ghost atom cutoff = 13 + binsize = 6.5, bins = 4 4 5 + 2 neighbor lists, perpetual/occasional/extra = 2 0 0 + (1) pair buck/coul/long, perpetual, half/full from (2) + attributes: half, newton on + pair build: halffull/newton + stencil: none + bin: none + (2) fix qeq/point, perpetual + attributes: full, newton on + pair build: full/bin/atomonly + stencil: full/bin/3d + bin: standard +Per MPI rank memory allocation (min/avg/max) = 24.69 | 24.69 | 24.69 Mbytes +Step PotEng c_q1 c_q2 v_qtot v_nqeq + 0 -3432.17988 0.85228288 -0.42614144 -0.000000000000 3 + 10 -3452.03328 0.85475605 -0.42737803 -0.000000000000 8 + 20 -3497.57515 0.85994936 -0.42997468 0.000000000000 8 + 30 -3568.22095 0.86767937 -0.43383969 0.000000000001 8 + 40 -3633.24956 0.87335551 -0.43667775 0.000000000000 8 + 50 -3700.10219 0.87805056 -0.43902528 0.000000000000 8 + 60 -3784.36769 0.88402303 -0.44201151 -0.000000000000 8 + 70 -3877.51378 0.89008950 -0.44504475 0.000000000000 8 + 80 -3965.29722 0.89431515 -0.44715757 -0.000000000001 8 + 90 -4048.36764 0.89698588 -0.44849294 0.000000000000 8 + 100 -4118.65809 0.89719102 -0.44859551 0.000000000000 8 +Loop time of 11.5935 on 1 procs for 100 steps with 1200 atoms + +Performance: 0.075 ns/day, 322.041 hours/ns, 8.626 timesteps/s +99.9% 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.8257 | 2.8257 | 2.8257 | 0.0 | 24.37 +Kspace | 1.2136 | 1.2136 | 1.2136 | 0.0 | 10.47 +Neigh | 0 | 0 | 0 | 0.0 | 0.00 +Comm | 0.015541 | 0.015541 | 0.015541 | 0.0 | 0.13 +Output | 0.0014489 | 0.0014489 | 0.0014489 | 0.0 | 0.01 +Modify | 7.5351 | 7.5351 | 7.5351 | 0.0 | 64.99 +Other | | 0.00206 | | | 0.02 + +Nlocal: 1200.00 ave 1200 max 1200 min +Histogram: 1 0 0 0 0 0 0 0 0 0 +Nghost: 8100.00 ave 8100 max 8100 min +Histogram: 1 0 0 0 0 0 0 0 0 0 +Neighs: 367600.0 ave 367600 max 367600 min +Histogram: 1 0 0 0 0 0 0 0 0 0 +FullNghs: 735200.0 ave 735200 max 735200 min +Histogram: 1 0 0 0 0 0 0 0 0 0 + +Total # of neighbors = 735200 +Ave neighs/atom = 612.66667 +Neighbor list builds = 0 +Dangerous builds = 0 + +clear + using 1 OpenMP thread(s) per MPI task + +print "Using fix qeq/shielded" +Using fix qeq/shielded +read_restart qeq.restart +Reading restart file ... + restart file = 8 Apr 2021, LAMMPS = 8 Apr 2021 + restoring atom style charge from restart + orthogonal box = (0.0000000 0.0000000 0.0000000) to (25.158320 25.158320 28.020256) + 1 by 1 by 1 MPI processor grid + restoring pair style buck/coul/long from restart + 1200 atoms + read_restart CPU = 0.001 seconds +fix 2 all qeq/shielded 1 10 1.0e-6 100 param.qeq2 +include buck.inc + +kspace_style pppm 1e-6 + +neighbor 1.0 bin +neigh_modify delay 0 every 1 check yes + +group type1 type 1 +400 atoms in group type1 +compute charge1 type1 property/atom q +compute q1 type1 reduce ave c_charge1 +group type2 type 2 +800 atoms in group type2 +compute charge2 type2 property/atom q +compute q2 type2 reduce ave c_charge2 +variable qtot equal count(type1)*c_q1+count(type2)*c_q2 +variable nqeq equal f_2 + +thermo_style custom step pe c_q1 c_q2 v_qtot v_nqeq +thermo 10 +thermo_modify format line "%4d %12.9g %12.8f %12.8f %16.12f %6.0f" + +timestep 0.0001 + +fix 1 all nve + +run 100 +PPPM initialization ... + using 12-bit tables for long-range coulomb (src/kspace.cpp:339) + G vector (1/distance) = 0.27644401 + grid = 27 27 30 + stencil order = 5 + estimated absolute RMS force accuracy = 1.4502702e-05 + estimated relative force accuracy = 1.0071569e-06 + using double precision FFTW3 + 3d grid and FFT values/proc = 42772 21870 +Neighbor list info ... + update every 1 steps, delay 0 steps, check yes + max neighbors/atom: 2000, page size: 100000 + master list distance cutoff = 13 + ghost atom cutoff = 13 + binsize = 6.5, bins = 4 4 5 + 2 neighbor lists, perpetual/occasional/extra = 2 0 0 + (1) pair buck/coul/long, perpetual, half/full from (2) + attributes: half, newton on + pair build: halffull/newton + stencil: none + bin: none + (2) fix qeq/shielded, perpetual + attributes: full, newton on + pair build: full/bin/atomonly + stencil: full/bin/3d + bin: standard +Per MPI rank memory allocation (min/avg/max) = 24.69 | 24.69 | 24.69 Mbytes +Step PotEng c_q1 c_q2 v_qtot v_nqeq + 0 -2879.00309 0.76536977 -0.38268489 -0.000000000000 3 + 10 -2882.50998 0.76536972 -0.38268486 0.000000000000 2 + 20 -2893.89472 0.76536950 -0.38268475 0.000000000000 2 + 30 -2913.6181 0.76536875 -0.38268438 0.000000000001 1 + 40 -2942.24129 0.76536939 -0.38268470 -0.000000000001 1 + 50 -2980.18817 0.76536780 -0.38268390 0.000000000000 2 + 60 -3027.60957 0.76536804 -0.38268402 0.000000000000 2 + 70 -3084.12552 0.76536573 -0.38268287 0.000000000000 2 + 80 -3148.8697 0.76536550 -0.38268275 0.000000000001 1 + 90 -3220.43086 0.76536380 -0.38268190 0.000000000000 2 + 100 -3297.0618 0.76536251 -0.38268126 0.000000000000 2 +Loop time of 7.93936 on 1 procs for 100 steps with 1200 atoms + +Performance: 0.109 ns/day, 220.538 hours/ns, 12.595 timesteps/s +99.9% 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.8061 | 2.8061 | 2.8061 | 0.0 | 35.34 +Kspace | 1.2176 | 1.2176 | 1.2176 | 0.0 | 15.34 +Neigh | 0 | 0 | 0 | 0.0 | 0.00 +Comm | 0.015528 | 0.015528 | 0.015528 | 0.0 | 0.20 +Output | 0.0014365 | 0.0014365 | 0.0014365 | 0.0 | 0.02 +Modify | 3.8966 | 3.8966 | 3.8966 | 0.0 | 49.08 +Other | | 0.002076 | | | 0.03 + +Nlocal: 1200.00 ave 1200 max 1200 min +Histogram: 1 0 0 0 0 0 0 0 0 0 +Nghost: 8100.00 ave 8100 max 8100 min +Histogram: 1 0 0 0 0 0 0 0 0 0 +Neighs: 367600.0 ave 367600 max 367600 min +Histogram: 1 0 0 0 0 0 0 0 0 0 +FullNghs: 735200.0 ave 735200 max 735200 min +Histogram: 1 0 0 0 0 0 0 0 0 0 + +Total # of neighbors = 735200 +Ave neighs/atom = 612.66667 +Neighbor list builds = 0 +Dangerous builds = 0 + + +clear + using 1 OpenMP thread(s) per MPI task + +print "Using fix qeq/slater" +Using fix qeq/slater +read_restart qeq.restart +Reading restart file ... + restart file = 8 Apr 2021, LAMMPS = 8 Apr 2021 + restoring atom style charge from restart + orthogonal box = (0.0000000 0.0000000 0.0000000) to (25.158320 25.158320 28.020256) + 1 by 1 by 1 MPI processor grid + restoring pair style buck/coul/long from restart + 1200 atoms + read_restart CPU = 0.001 seconds +fix 2 all qeq/shielded 1 10 1.0e-6 100 param.qeq2 +include buck.inc + +kspace_style pppm 1e-6 + +neighbor 1.0 bin +neigh_modify delay 0 every 1 check yes + +group type1 type 1 +400 atoms in group type1 +compute charge1 type1 property/atom q +compute q1 type1 reduce ave c_charge1 +group type2 type 2 +800 atoms in group type2 +compute charge2 type2 property/atom q +compute q2 type2 reduce ave c_charge2 +variable qtot equal count(type1)*c_q1+count(type2)*c_q2 +variable nqeq equal f_2 + +thermo_style custom step pe c_q1 c_q2 v_qtot v_nqeq +thermo 10 +thermo_modify format line "%4d %12.9g %12.8f %12.8f %16.12f %6.0f" + +timestep 0.0001 + +fix 1 all nve + +run 100 +PPPM initialization ... + using 12-bit tables for long-range coulomb (src/kspace.cpp:339) + G vector (1/distance) = 0.27644401 + grid = 27 27 30 + stencil order = 5 + estimated absolute RMS force accuracy = 1.4502702e-05 + estimated relative force accuracy = 1.0071569e-06 + using double precision FFTW3 + 3d grid and FFT values/proc = 42772 21870 +Neighbor list info ... + update every 1 steps, delay 0 steps, check yes + max neighbors/atom: 2000, page size: 100000 + master list distance cutoff = 13 + ghost atom cutoff = 13 + binsize = 6.5, bins = 4 4 5 + 2 neighbor lists, perpetual/occasional/extra = 2 0 0 + (1) pair buck/coul/long, perpetual, half/full from (2) + attributes: half, newton on + pair build: halffull/newton + stencil: none + bin: none + (2) fix qeq/shielded, perpetual + attributes: full, newton on + pair build: full/bin/atomonly + stencil: full/bin/3d + bin: standard +Per MPI rank memory allocation (min/avg/max) = 24.69 | 24.69 | 24.69 Mbytes +Step PotEng c_q1 c_q2 v_qtot v_nqeq + 0 -2879.00309 0.76536977 -0.38268489 -0.000000000000 3 + 10 -2882.50998 0.76536972 -0.38268486 0.000000000000 2 + 20 -2893.89472 0.76536950 -0.38268475 0.000000000000 2 + 30 -2913.6181 0.76536875 -0.38268438 0.000000000001 1 + 40 -2942.24129 0.76536939 -0.38268470 -0.000000000001 1 + 50 -2980.18817 0.76536780 -0.38268390 0.000000000000 2 + 60 -3027.60957 0.76536804 -0.38268402 0.000000000000 2 + 70 -3084.12552 0.76536573 -0.38268287 0.000000000000 2 + 80 -3148.8697 0.76536550 -0.38268275 0.000000000001 1 + 90 -3220.43086 0.76536380 -0.38268190 0.000000000000 2 + 100 -3297.0618 0.76536251 -0.38268126 0.000000000000 2 +Loop time of 7.9652 on 1 procs for 100 steps with 1200 atoms + +Performance: 0.108 ns/day, 221.256 hours/ns, 12.555 timesteps/s +99.9% 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.809 | 2.809 | 2.809 | 0.0 | 35.27 +Kspace | 1.2214 | 1.2214 | 1.2214 | 0.0 | 15.33 +Neigh | 0 | 0 | 0 | 0.0 | 0.00 +Comm | 0.015635 | 0.015635 | 0.015635 | 0.0 | 0.20 +Output | 0.0014393 | 0.0014393 | 0.0014393 | 0.0 | 0.02 +Modify | 3.9157 | 3.9157 | 3.9157 | 0.0 | 49.16 +Other | | 0.002091 | | | 0.03 + +Nlocal: 1200.00 ave 1200 max 1200 min +Histogram: 1 0 0 0 0 0 0 0 0 0 +Nghost: 8100.00 ave 8100 max 8100 min +Histogram: 1 0 0 0 0 0 0 0 0 0 +Neighs: 367600.0 ave 367600 max 367600 min +Histogram: 1 0 0 0 0 0 0 0 0 0 +FullNghs: 735200.0 ave 735200 max 735200 min +Histogram: 1 0 0 0 0 0 0 0 0 0 + +Total # of neighbors = 735200 +Ave neighs/atom = 612.66667 +Neighbor list builds = 0 +Dangerous builds = 0 + +clear + using 1 OpenMP thread(s) per MPI task + +print "Using fix qeq/dynamic" +Using fix qeq/dynamic +read_restart qeq.restart +Reading restart file ... + restart file = 8 Apr 2021, LAMMPS = 8 Apr 2021 + restoring atom style charge from restart + orthogonal box = (0.0000000 0.0000000 0.0000000) to (25.158320 25.158320 28.020256) + 1 by 1 by 1 MPI processor grid + restoring pair style buck/coul/long from restart + 1200 atoms + read_restart CPU = 0.001 seconds +fix 2 all qeq/dynamic 1 10 1.0e-3 100 param.qeq2 +include buck.inc + +kspace_style pppm 1e-6 + +neighbor 1.0 bin +neigh_modify delay 0 every 1 check yes + +group type1 type 1 +400 atoms in group type1 +compute charge1 type1 property/atom q +compute q1 type1 reduce ave c_charge1 +group type2 type 2 +800 atoms in group type2 +compute charge2 type2 property/atom q +compute q2 type2 reduce ave c_charge2 +variable qtot equal count(type1)*c_q1+count(type2)*c_q2 +variable nqeq equal f_2 + +thermo_style custom step pe c_q1 c_q2 v_qtot v_nqeq +thermo 10 +thermo_modify format line "%4d %12.9g %12.8f %12.8f %16.12f %6.0f" + +timestep 0.0001 + +fix 1 all nve + +run 100 +PPPM initialization ... + using 12-bit tables for long-range coulomb (src/kspace.cpp:339) + G vector (1/distance) = 0.27644401 + grid = 27 27 30 + stencil order = 5 + estimated absolute RMS force accuracy = 1.4502702e-05 + estimated relative force accuracy = 1.0071569e-06 + using double precision FFTW3 + 3d grid and FFT values/proc = 42772 21870 +Neighbor list info ... + update every 1 steps, delay 0 steps, check yes + max neighbors/atom: 2000, page size: 100000 + master list distance cutoff = 13 + ghost atom cutoff = 13 + binsize = 6.5, bins = 4 4 5 + 2 neighbor lists, perpetual/occasional/extra = 2 0 0 + (1) pair buck/coul/long, perpetual + attributes: half, newton on + pair build: half/bin/atomonly/newton + stencil: half/bin/3d/newton + bin: standard + (2) fix qeq/dynamic, perpetual, copy from (1) + attributes: half, newton on + pair build: copy + stencil: none + bin: none +Per MPI rank memory allocation (min/avg/max) = 17.87 | 17.87 | 17.87 Mbytes +Step PotEng c_q1 c_q2 v_qtot v_nqeq + 0 -3432.38094 0.85231286 -0.42615643 0.000000000001 43 + 10 -3452.05217 0.85475894 -0.42737947 -0.000000000001 17 + 20 -3497.8643 0.85999180 -0.42999590 -0.000000000007 22 + 30 -3568.53169 0.86772479 -0.43386239 -0.000000000006 22 + 40 -3633.43753 0.87338291 -0.43669146 -0.000000000006 22 + 50 -3700.27953 0.87807632 -0.43903816 -0.000000000005 22 + 60 -3784.4004 0.88402822 -0.44201411 0.000000000002 17 + 70 -3877.73706 0.89012201 -0.44506100 0.000000000002 22 + 80 -3965.36111 0.89432486 -0.44716243 0.000000000008 17 + 90 -4048.57901 0.89701688 -0.44850844 -0.000000000004 22 + 100 -4118.62736 0.89718691 -0.44859346 -0.000000000026 17 +Loop time of 18.5333 on 1 procs for 100 steps with 1200 atoms + +Performance: 0.047 ns/day, 514.815 hours/ns, 5.396 timesteps/s +99.9% 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.8268 | 2.8268 | 2.8268 | 0.0 | 15.25 +Kspace | 1.2138 | 1.2138 | 1.2138 | 0.0 | 6.55 +Neigh | 0 | 0 | 0 | 0.0 | 0.00 +Comm | 0.015407 | 0.015407 | 0.015407 | 0.0 | 0.08 +Output | 0.0014303 | 0.0014303 | 0.0014303 | 0.0 | 0.01 +Modify | 14.474 | 14.474 | 14.474 | 0.0 | 78.10 +Other | | 0.001973 | | | 0.01 + +Nlocal: 1200.00 ave 1200 max 1200 min +Histogram: 1 0 0 0 0 0 0 0 0 0 +Nghost: 8100.00 ave 8100 max 8100 min +Histogram: 1 0 0 0 0 0 0 0 0 0 +Neighs: 367600.0 ave 367600 max 367600 min +Histogram: 1 0 0 0 0 0 0 0 0 0 + +Total # of neighbors = 367600 +Ave neighs/atom = 306.33333 +Neighbor list builds = 0 +Dangerous builds = 0 + +clear + using 1 OpenMP thread(s) per MPI task + +print "Using fix qeq/fire" +Using fix qeq/fire +read_restart qeq.restart +Reading restart file ... + restart file = 8 Apr 2021, LAMMPS = 8 Apr 2021 + restoring atom style charge from restart + orthogonal box = (0.0000000 0.0000000 0.0000000) to (25.158320 25.158320 28.020256) + 1 by 1 by 1 MPI processor grid + restoring pair style buck/coul/long from restart + 1200 atoms + read_restart CPU = 0.001 seconds +fix 2 all qeq/fire 1 10 1.0e-3 100 param.qeq2 +include buck.inc + +kspace_style pppm 1e-6 + +neighbor 1.0 bin +neigh_modify delay 0 every 1 check yes + +group type1 type 1 +400 atoms in group type1 +compute charge1 type1 property/atom q +compute q1 type1 reduce ave c_charge1 +group type2 type 2 +800 atoms in group type2 +compute charge2 type2 property/atom q +compute q2 type2 reduce ave c_charge2 +variable qtot equal count(type1)*c_q1+count(type2)*c_q2 +variable nqeq equal f_2 + +thermo_style custom step pe c_q1 c_q2 v_qtot v_nqeq +thermo 10 +thermo_modify format line "%4d %12.9g %12.8f %12.8f %16.12f %6.0f" + +timestep 0.0001 + +fix 1 all nve + +run 100 +PPPM initialization ... + using 12-bit tables for long-range coulomb (src/kspace.cpp:339) + G vector (1/distance) = 0.27644401 + grid = 27 27 30 + stencil order = 5 + estimated absolute RMS force accuracy = 1.4502702e-05 + estimated relative force accuracy = 1.0071569e-06 + using double precision FFTW3 + 3d grid and FFT values/proc = 42772 21870 +Neighbor list info ... + update every 1 steps, delay 0 steps, check yes + max neighbors/atom: 2000, page size: 100000 + master list distance cutoff = 13 + ghost atom cutoff = 13 + binsize = 6.5, bins = 4 4 5 + 2 neighbor lists, perpetual/occasional/extra = 2 0 0 + (1) pair buck/coul/long, perpetual + attributes: half, newton on + pair build: half/bin/atomonly/newton + stencil: half/bin/3d/newton + bin: standard + (2) fix qeq/fire, perpetual, copy from (1) + attributes: half, newton on + pair build: copy + stencil: none + bin: none +Per MPI rank memory allocation (min/avg/max) = 17.87 | 17.87 | 17.87 Mbytes +Step PotEng c_q1 c_q2 v_qtot v_nqeq + 0 -3432.06113 0.85226679 -0.42613339 0.000000000004 37 + 10 -3452.0494 0.85475813 -0.42737906 0.000000000001 10 + 20 -3497.83503 0.85998739 -0.42999370 0.000000000003 13 + 30 -3568.47507 0.86771599 -0.43385799 0.000000000004 13 + 40 -3633.35368 0.87337029 -0.43668514 0.000000000004 13 + 50 -3700.15601 0.87805847 -0.43902924 0.000000000005 13 + 60 -3784.32042 0.88401635 -0.44200818 0.000000000000 11 + 70 -3877.59818 0.89010162 -0.44505081 0.000000000000 13 + 80 -3965.28426 0.89431356 -0.44715678 0.000000000000 11 + 90 -4048.3338 0.89698069 -0.44849034 0.000000000001 13 + 100 -4118.63638 0.89718818 -0.44859409 0.000000000003 12 +Loop time of 13.0492 on 1 procs for 100 steps with 1200 atoms + +Performance: 0.066 ns/day, 362.479 hours/ns, 7.663 timesteps/s +99.9% 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.7996 | 2.7996 | 2.7996 | 0.0 | 21.45 +Kspace | 1.2141 | 1.2141 | 1.2141 | 0.0 | 9.30 +Neigh | 0 | 0 | 0 | 0.0 | 0.00 +Comm | 0.015527 | 0.015527 | 0.015527 | 0.0 | 0.12 +Output | 0.0014405 | 0.0014405 | 0.0014405 | 0.0 | 0.01 +Modify | 9.0166 | 9.0166 | 9.0166 | 0.0 | 69.10 +Other | | 0.001981 | | | 0.02 + +Nlocal: 1200.00 ave 1200 max 1200 min +Histogram: 1 0 0 0 0 0 0 0 0 0 +Nghost: 8100.00 ave 8100 max 8100 min +Histogram: 1 0 0 0 0 0 0 0 0 0 +Neighs: 367600.0 ave 367600 max 367600 min +Histogram: 1 0 0 0 0 0 0 0 0 0 + +Total # of neighbors = 367600 +Ave neighs/atom = 306.33333 +Neighbor list builds = 0 +Dangerous builds = 0 +Total wall time: 0:01:00 diff --git a/examples/qeq/log.20Apr21.qeq.buck.g++.4 b/examples/qeq/log.20Apr21.qeq.buck.g++.4 new file mode 100644 index 0000000000..cafa94bae9 --- /dev/null +++ b/examples/qeq/log.20Apr21.qeq.buck.g++.4 @@ -0,0 +1,650 @@ +LAMMPS (8 Apr 2021) + using 1 OpenMP thread(s) per MPI task +# This example demonstrates the use of various fix qeq variants with +# a pair style using charges, in this case pair_style buck/coul/long + +units metal +atom_style charge + +read_data data.aC +Reading data file ... + orthogonal box = (0.0000000 0.0000000 0.0000000) to (25.158320 25.158320 28.020256) + 1 by 2 by 2 MPI processor grid + reading atoms ... + 1200 atoms + read_data CPU = 0.009 seconds +#replicate 2 2 2 + +pair_style buck/coul/long 12.0 +pair_coeff 2 2 1388.77 .3623188 175.0 +pair_coeff 1 2 18003 .2052124 133.5381 +pair_coeff 1 1 0 .1 0 + +fix 2 all qeq/shielded 1 10 1.0e-20 10 param.qeq2 + +include buck.inc + +kspace_style pppm 1e-6 + +neighbor 1.0 bin +neigh_modify delay 0 every 1 check yes + +group type1 type 1 +400 atoms in group type1 +compute charge1 type1 property/atom q +compute q1 type1 reduce ave c_charge1 +group type2 type 2 +800 atoms in group type2 +compute charge2 type2 property/atom q +compute q2 type2 reduce ave c_charge2 +variable qtot equal count(type1)*c_q1+count(type2)*c_q2 +variable nqeq equal f_2 + +thermo_style custom step pe c_q1 c_q2 v_qtot v_nqeq +thermo 10 +thermo_modify format line "%4d %12.9g %12.8f %12.8f %16.12f %6.0f" + +timestep 0.0001 + +fix 1 all nve + +velocity all create 300.0 1281937 +run 0 post no +PPPM initialization ... + using 12-bit tables for long-range coulomb (src/kspace.cpp:339) + G vector (1/distance) = 0.30705229 + grid = 48 48 54 + stencil order = 5 + estimated absolute RMS force accuracy = 1.8909403e-05 + estimated relative force accuracy = 1.3131854e-06 + using double precision FFTW3 + 3d grid and FFT values/proc = 57970 32256 +Neighbor list info ... + update every 1 steps, delay 0 steps, check yes + max neighbors/atom: 2000, page size: 100000 + master list distance cutoff = 13 + ghost atom cutoff = 13 + binsize = 6.5, bins = 4 4 5 + 2 neighbor lists, perpetual/occasional/extra = 2 0 0 + (1) pair buck/coul/long, perpetual, half/full from (2) + attributes: half, newton on + pair build: halffull/newton + stencil: none + bin: none + (2) fix qeq/shielded, perpetual + attributes: full, newton on + pair build: full/bin/atomonly + stencil: full/bin/3d + bin: standard +WARNING: Fix qeq CG convergence failed (4.299911728887494e-19) after 10 iterations at step 0 (src/QEQ/fix_qeq.cpp:410) +WARNING: Fix qeq CG convergence failed (5.273380778822746e-18) after 10 iterations at step 0 (src/QEQ/fix_qeq.cpp:410) +Per MPI rank memory allocation (min/avg/max) = 14.97 | 15.02 | 15.08 Mbytes +Step PotEng c_q1 c_q2 v_qtot v_nqeq + 0 -2879.00327 0.76536977 -0.38268489 0.000000000000 10 +Loop time of 3.33786e-06 on 4 procs for 0 steps with 1200 atoms + + +write_restart qeq.restart +System init for write_restart ... +PPPM initialization ... + using 12-bit tables for long-range coulomb (src/kspace.cpp:339) + G vector (1/distance) = 0.27644401 + grid = 27 27 30 + stencil order = 5 + estimated absolute RMS force accuracy = 1.4502702e-05 + estimated relative force accuracy = 1.0071569e-06 + using double precision FFTW3 + 3d grid and FFT values/proc = 14960 5832 + +clear + using 1 OpenMP thread(s) per MPI task + +print "Using fix qeq/point" +Using fix qeq/point +read_restart qeq.restart +Reading restart file ... + restart file = 8 Apr 2021, LAMMPS = 8 Apr 2021 + restoring atom style charge from restart + orthogonal box = (0.0000000 0.0000000 0.0000000) to (25.158320 25.158320 28.020256) + 1 by 2 by 2 MPI processor grid + restoring pair style buck/coul/long from restart + 1200 atoms + read_restart CPU = 0.001 seconds +fix 2 all qeq/point 1 10 1.0e-6 100 param.qeq2 +include buck.inc + +kspace_style pppm 1e-6 + +neighbor 1.0 bin +neigh_modify delay 0 every 1 check yes + +group type1 type 1 +400 atoms in group type1 +compute charge1 type1 property/atom q +compute q1 type1 reduce ave c_charge1 +group type2 type 2 +800 atoms in group type2 +compute charge2 type2 property/atom q +compute q2 type2 reduce ave c_charge2 +variable qtot equal count(type1)*c_q1+count(type2)*c_q2 +variable nqeq equal f_2 + +thermo_style custom step pe c_q1 c_q2 v_qtot v_nqeq +thermo 10 +thermo_modify format line "%4d %12.9g %12.8f %12.8f %16.12f %6.0f" + +timestep 0.0001 + +fix 1 all nve + +run 100 +PPPM initialization ... + using 12-bit tables for long-range coulomb (src/kspace.cpp:339) + G vector (1/distance) = 0.27644401 + grid = 27 27 30 + stencil order = 5 + estimated absolute RMS force accuracy = 1.4502702e-05 + estimated relative force accuracy = 1.0071569e-06 + using double precision FFTW3 + 3d grid and FFT values/proc = 14960 5832 +Neighbor list info ... + update every 1 steps, delay 0 steps, check yes + max neighbors/atom: 2000, page size: 100000 + master list distance cutoff = 13 + ghost atom cutoff = 13 + binsize = 6.5, bins = 4 4 5 + 2 neighbor lists, perpetual/occasional/extra = 2 0 0 + (1) pair buck/coul/long, perpetual, half/full from (2) + attributes: half, newton on + pair build: halffull/newton + stencil: none + bin: none + (2) fix qeq/point, perpetual + attributes: full, newton on + pair build: full/bin/atomonly + stencil: full/bin/3d + bin: standard +Per MPI rank memory allocation (min/avg/max) = 11.10 | 11.14 | 11.16 Mbytes +Step PotEng c_q1 c_q2 v_qtot v_nqeq + 0 -3432.17988 0.85228288 -0.42614144 -0.000000000000 3 + 10 -3452.03328 0.85475605 -0.42737803 -0.000000000000 8 + 20 -3497.57515 0.85994936 -0.42997468 0.000000000000 8 + 30 -3568.22095 0.86767937 -0.43383969 0.000000000000 8 + 40 -3633.24956 0.87335551 -0.43667775 -0.000000000000 8 + 50 -3700.10219 0.87805056 -0.43902528 0.000000000000 8 + 60 -3784.36769 0.88402303 -0.44201151 0.000000000000 8 + 70 -3877.51378 0.89008950 -0.44504475 0.000000000000 8 + 80 -3965.29722 0.89431515 -0.44715757 0.000000000000 8 + 90 -4048.36764 0.89698588 -0.44849294 -0.000000000000 8 + 100 -4118.65809 0.89719102 -0.44859551 0.000000000000 8 +Loop time of 3.30911 on 4 procs for 100 steps with 1200 atoms + +Performance: 0.261 ns/day, 91.920 hours/ns, 30.220 timesteps/s +99.0% 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.67613 | 0.68904 | 0.71562 | 1.9 | 20.82 +Kspace | 0.36056 | 0.3881 | 0.39892 | 2.6 | 11.73 +Neigh | 0 | 0 | 0 | 0.0 | 0.00 +Comm | 0.013339 | 0.017982 | 0.019974 | 2.0 | 0.54 +Output | 0.0006721 | 0.00099713 | 0.0019572 | 0.0 | 0.03 +Modify | 2.2109 | 2.211 | 2.211 | 0.0 | 66.81 +Other | | 0.002041 | | | 0.06 + +Nlocal: 300.000 ave 300 max 300 min +Histogram: 4 0 0 0 0 0 0 0 0 0 +Nghost: 4875.00 ave 4880 max 4870 min +Histogram: 2 0 0 0 0 0 0 0 0 2 +Neighs: 91900.0 ave 91900 max 91900 min +Histogram: 4 0 0 0 0 0 0 0 0 0 +FullNghs: 183800.0 ave 183800 max 183800 min +Histogram: 4 0 0 0 0 0 0 0 0 0 + +Total # of neighbors = 735200 +Ave neighs/atom = 612.66667 +Neighbor list builds = 0 +Dangerous builds = 0 + +clear + using 1 OpenMP thread(s) per MPI task + +print "Using fix qeq/shielded" +Using fix qeq/shielded +read_restart qeq.restart +Reading restart file ... + restart file = 8 Apr 2021, LAMMPS = 8 Apr 2021 + restoring atom style charge from restart + orthogonal box = (0.0000000 0.0000000 0.0000000) to (25.158320 25.158320 28.020256) + 1 by 2 by 2 MPI processor grid + restoring pair style buck/coul/long from restart + 1200 atoms + read_restart CPU = 0.006 seconds +fix 2 all qeq/shielded 1 10 1.0e-6 100 param.qeq2 +include buck.inc + +kspace_style pppm 1e-6 + +neighbor 1.0 bin +neigh_modify delay 0 every 1 check yes + +group type1 type 1 +400 atoms in group type1 +compute charge1 type1 property/atom q +compute q1 type1 reduce ave c_charge1 +group type2 type 2 +800 atoms in group type2 +compute charge2 type2 property/atom q +compute q2 type2 reduce ave c_charge2 +variable qtot equal count(type1)*c_q1+count(type2)*c_q2 +variable nqeq equal f_2 + +thermo_style custom step pe c_q1 c_q2 v_qtot v_nqeq +thermo 10 +thermo_modify format line "%4d %12.9g %12.8f %12.8f %16.12f %6.0f" + +timestep 0.0001 + +fix 1 all nve + +run 100 +PPPM initialization ... + using 12-bit tables for long-range coulomb (src/kspace.cpp:339) + G vector (1/distance) = 0.27644401 + grid = 27 27 30 + stencil order = 5 + estimated absolute RMS force accuracy = 1.4502702e-05 + estimated relative force accuracy = 1.0071569e-06 + using double precision FFTW3 + 3d grid and FFT values/proc = 14960 5832 +Neighbor list info ... + update every 1 steps, delay 0 steps, check yes + max neighbors/atom: 2000, page size: 100000 + master list distance cutoff = 13 + ghost atom cutoff = 13 + binsize = 6.5, bins = 4 4 5 + 2 neighbor lists, perpetual/occasional/extra = 2 0 0 + (1) pair buck/coul/long, perpetual, half/full from (2) + attributes: half, newton on + pair build: halffull/newton + stencil: none + bin: none + (2) fix qeq/shielded, perpetual + attributes: full, newton on + pair build: full/bin/atomonly + stencil: full/bin/3d + bin: standard +Per MPI rank memory allocation (min/avg/max) = 11.10 | 11.14 | 11.16 Mbytes +Step PotEng c_q1 c_q2 v_qtot v_nqeq + 0 -2879.00309 0.76536977 -0.38268489 0.000000000000 3 + 10 -2882.50998 0.76536972 -0.38268486 0.000000000000 2 + 20 -2893.89472 0.76536950 -0.38268475 -0.000000000000 2 + 30 -2913.6181 0.76536875 -0.38268438 -0.000000000000 1 + 40 -2942.24129 0.76536939 -0.38268470 0.000000000000 1 + 50 -2980.18817 0.76536780 -0.38268390 0.000000000000 2 + 60 -3027.60957 0.76536804 -0.38268402 -0.000000000000 2 + 70 -3084.12552 0.76536573 -0.38268287 0.000000000000 2 + 80 -3148.8697 0.76536550 -0.38268275 0.000000000000 1 + 90 -3220.43086 0.76536380 -0.38268190 0.000000000000 2 + 100 -3297.0618 0.76536251 -0.38268126 0.000000000000 2 +Loop time of 2.25559 on 4 procs for 100 steps with 1200 atoms + +Performance: 0.383 ns/day, 62.655 hours/ns, 44.334 timesteps/s +97.9% 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.67442 | 0.69181 | 0.70907 | 2.0 | 30.67 +Kspace | 0.39381 | 0.41151 | 0.43023 | 2.6 | 18.24 +Neigh | 0 | 0 | 0 | 0.0 | 0.00 +Comm | 0.012851 | 0.01426 | 0.015146 | 0.7 | 0.63 +Output | 0.00066686 | 0.00098681 | 0.0019395 | 0.0 | 0.04 +Modify | 1.1349 | 1.135 | 1.135 | 0.0 | 50.32 +Other | | 0.002035 | | | 0.09 + +Nlocal: 300.000 ave 300 max 300 min +Histogram: 4 0 0 0 0 0 0 0 0 0 +Nghost: 4875.00 ave 4880 max 4870 min +Histogram: 2 0 0 0 0 0 0 0 0 2 +Neighs: 91900.0 ave 91900 max 91900 min +Histogram: 4 0 0 0 0 0 0 0 0 0 +FullNghs: 183800.0 ave 183800 max 183800 min +Histogram: 4 0 0 0 0 0 0 0 0 0 + +Total # of neighbors = 735200 +Ave neighs/atom = 612.66667 +Neighbor list builds = 0 +Dangerous builds = 0 + + +clear + using 1 OpenMP thread(s) per MPI task + +print "Using fix qeq/slater" +Using fix qeq/slater +read_restart qeq.restart +Reading restart file ... + restart file = 8 Apr 2021, LAMMPS = 8 Apr 2021 + restoring atom style charge from restart + orthogonal box = (0.0000000 0.0000000 0.0000000) to (25.158320 25.158320 28.020256) + 1 by 2 by 2 MPI processor grid + restoring pair style buck/coul/long from restart + 1200 atoms + read_restart CPU = 0.012 seconds +fix 2 all qeq/shielded 1 10 1.0e-6 100 param.qeq2 +include buck.inc + +kspace_style pppm 1e-6 + +neighbor 1.0 bin +neigh_modify delay 0 every 1 check yes + +group type1 type 1 +400 atoms in group type1 +compute charge1 type1 property/atom q +compute q1 type1 reduce ave c_charge1 +group type2 type 2 +800 atoms in group type2 +compute charge2 type2 property/atom q +compute q2 type2 reduce ave c_charge2 +variable qtot equal count(type1)*c_q1+count(type2)*c_q2 +variable nqeq equal f_2 + +thermo_style custom step pe c_q1 c_q2 v_qtot v_nqeq +thermo 10 +thermo_modify format line "%4d %12.9g %12.8f %12.8f %16.12f %6.0f" + +timestep 0.0001 + +fix 1 all nve + +run 100 +PPPM initialization ... + using 12-bit tables for long-range coulomb (src/kspace.cpp:339) + G vector (1/distance) = 0.27644401 + grid = 27 27 30 + stencil order = 5 + estimated absolute RMS force accuracy = 1.4502702e-05 + estimated relative force accuracy = 1.0071569e-06 + using double precision FFTW3 + 3d grid and FFT values/proc = 14960 5832 +Neighbor list info ... + update every 1 steps, delay 0 steps, check yes + max neighbors/atom: 2000, page size: 100000 + master list distance cutoff = 13 + ghost atom cutoff = 13 + binsize = 6.5, bins = 4 4 5 + 2 neighbor lists, perpetual/occasional/extra = 2 0 0 + (1) pair buck/coul/long, perpetual, half/full from (2) + attributes: half, newton on + pair build: halffull/newton + stencil: none + bin: none + (2) fix qeq/shielded, perpetual + attributes: full, newton on + pair build: full/bin/atomonly + stencil: full/bin/3d + bin: standard +Per MPI rank memory allocation (min/avg/max) = 11.10 | 11.14 | 11.16 Mbytes +Step PotEng c_q1 c_q2 v_qtot v_nqeq + 0 -2879.00309 0.76536977 -0.38268489 0.000000000000 3 + 10 -2882.50998 0.76536972 -0.38268486 0.000000000000 2 + 20 -2893.89472 0.76536950 -0.38268475 -0.000000000000 2 + 30 -2913.6181 0.76536875 -0.38268438 -0.000000000000 1 + 40 -2942.24129 0.76536939 -0.38268470 0.000000000000 1 + 50 -2980.18817 0.76536780 -0.38268390 0.000000000000 2 + 60 -3027.60957 0.76536804 -0.38268402 -0.000000000000 2 + 70 -3084.12552 0.76536573 -0.38268287 0.000000000000 2 + 80 -3148.8697 0.76536550 -0.38268275 0.000000000000 1 + 90 -3220.43086 0.76536380 -0.38268190 0.000000000000 2 + 100 -3297.0618 0.76536251 -0.38268126 0.000000000000 2 +Loop time of 2.39249 on 4 procs for 100 steps with 1200 atoms + +Performance: 0.361 ns/day, 66.458 hours/ns, 41.797 timesteps/s +96.0% 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.6751 | 0.70301 | 0.71919 | 2.1 | 29.38 +Kspace | 0.45569 | 0.47315 | 0.49885 | 2.6 | 19.78 +Neigh | 0 | 0 | 0 | 0.0 | 0.00 +Comm | 0.012967 | 0.018681 | 0.020909 | 2.4 | 0.78 +Output | 0.00066733 | 0.00099397 | 0.0019579 | 0.0 | 0.04 +Modify | 1.1945 | 1.1946 | 1.1947 | 0.0 | 49.93 +Other | | 0.002046 | | | 0.09 + +Nlocal: 300.000 ave 300 max 300 min +Histogram: 4 0 0 0 0 0 0 0 0 0 +Nghost: 4875.00 ave 4880 max 4870 min +Histogram: 2 0 0 0 0 0 0 0 0 2 +Neighs: 91900.0 ave 91900 max 91900 min +Histogram: 4 0 0 0 0 0 0 0 0 0 +FullNghs: 183800.0 ave 183800 max 183800 min +Histogram: 4 0 0 0 0 0 0 0 0 0 + +Total # of neighbors = 735200 +Ave neighs/atom = 612.66667 +Neighbor list builds = 0 +Dangerous builds = 0 + +clear + using 1 OpenMP thread(s) per MPI task + +print "Using fix qeq/dynamic" +Using fix qeq/dynamic +read_restart qeq.restart +Reading restart file ... + restart file = 8 Apr 2021, LAMMPS = 8 Apr 2021 + restoring atom style charge from restart + orthogonal box = (0.0000000 0.0000000 0.0000000) to (25.158320 25.158320 28.020256) + 1 by 2 by 2 MPI processor grid + restoring pair style buck/coul/long from restart + 1200 atoms + read_restart CPU = 0.002 seconds +fix 2 all qeq/dynamic 1 10 1.0e-3 100 param.qeq2 +include buck.inc + +kspace_style pppm 1e-6 + +neighbor 1.0 bin +neigh_modify delay 0 every 1 check yes + +group type1 type 1 +400 atoms in group type1 +compute charge1 type1 property/atom q +compute q1 type1 reduce ave c_charge1 +group type2 type 2 +800 atoms in group type2 +compute charge2 type2 property/atom q +compute q2 type2 reduce ave c_charge2 +variable qtot equal count(type1)*c_q1+count(type2)*c_q2 +variable nqeq equal f_2 + +thermo_style custom step pe c_q1 c_q2 v_qtot v_nqeq +thermo 10 +thermo_modify format line "%4d %12.9g %12.8f %12.8f %16.12f %6.0f" + +timestep 0.0001 + +fix 1 all nve + +run 100 +PPPM initialization ... + using 12-bit tables for long-range coulomb (src/kspace.cpp:339) + G vector (1/distance) = 0.27644401 + grid = 27 27 30 + stencil order = 5 + estimated absolute RMS force accuracy = 1.4502702e-05 + estimated relative force accuracy = 1.0071569e-06 + using double precision FFTW3 + 3d grid and FFT values/proc = 14960 5832 +Neighbor list info ... + update every 1 steps, delay 0 steps, check yes + max neighbors/atom: 2000, page size: 100000 + master list distance cutoff = 13 + ghost atom cutoff = 13 + binsize = 6.5, bins = 4 4 5 + 2 neighbor lists, perpetual/occasional/extra = 2 0 0 + (1) pair buck/coul/long, perpetual + attributes: half, newton on + pair build: half/bin/atomonly/newton + stencil: half/bin/3d/newton + bin: standard + (2) fix qeq/dynamic, perpetual, copy from (1) + attributes: half, newton on + pair build: copy + stencil: none + bin: none +Per MPI rank memory allocation (min/avg/max) = 9.195 | 9.246 | 9.278 Mbytes +Step PotEng c_q1 c_q2 v_qtot v_nqeq + 0 -3432.38094 0.85231286 -0.42615643 -0.000000000001 43 + 10 -3452.05217 0.85475894 -0.42737947 -0.000000000003 17 + 20 -3497.8643 0.85999180 -0.42999590 0.000000000000 22 + 30 -3568.53169 0.86772479 -0.43386239 -0.000000000000 22 + 40 -3633.43753 0.87338291 -0.43669146 0.000000000006 22 + 50 -3700.27953 0.87807632 -0.43903816 0.000000000003 22 + 60 -3784.4004 0.88402822 -0.44201411 0.000000000009 17 + 70 -3877.73706 0.89012201 -0.44506100 0.000000000010 22 + 80 -3965.36111 0.89432486 -0.44716243 0.000000000011 17 + 90 -4048.57901 0.89701688 -0.44850844 0.000000000012 22 + 100 -4118.62736 0.89718691 -0.44859346 0.000000000013 17 +Loop time of 5.27704 on 4 procs for 100 steps with 1200 atoms + +Performance: 0.164 ns/day, 146.584 hours/ns, 18.950 timesteps/s +98.5% 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.68437 | 0.69096 | 0.69826 | 0.7 | 13.09 +Kspace | 0.38484 | 0.38941 | 0.39524 | 0.7 | 7.38 +Neigh | 0 | 0 | 0 | 0.0 | 0.00 +Comm | 0.012609 | 0.01529 | 0.016842 | 1.3 | 0.29 +Output | 0.00067735 | 0.0010006 | 0.0019588 | 1.7 | 0.02 +Modify | 4.1783 | 4.1783 | 4.1784 | 0.0 | 79.18 +Other | | 0.002027 | | | 0.04 + +Nlocal: 300.000 ave 300 max 300 min +Histogram: 4 0 0 0 0 0 0 0 0 0 +Nghost: 4875.00 ave 4880 max 4870 min +Histogram: 2 0 0 0 0 0 0 0 0 2 +Neighs: 91900.0 ave 93081 max 90719 min +Histogram: 2 0 0 0 0 0 0 0 0 2 + +Total # of neighbors = 367600 +Ave neighs/atom = 306.33333 +Neighbor list builds = 0 +Dangerous builds = 0 + +clear + using 1 OpenMP thread(s) per MPI task + +print "Using fix qeq/fire" +Using fix qeq/fire +read_restart qeq.restart +Reading restart file ... + restart file = 8 Apr 2021, LAMMPS = 8 Apr 2021 + restoring atom style charge from restart + orthogonal box = (0.0000000 0.0000000 0.0000000) to (25.158320 25.158320 28.020256) + 1 by 2 by 2 MPI processor grid + restoring pair style buck/coul/long from restart + 1200 atoms + read_restart CPU = 0.001 seconds +fix 2 all qeq/fire 1 10 1.0e-3 100 param.qeq2 +include buck.inc + +kspace_style pppm 1e-6 + +neighbor 1.0 bin +neigh_modify delay 0 every 1 check yes + +group type1 type 1 +400 atoms in group type1 +compute charge1 type1 property/atom q +compute q1 type1 reduce ave c_charge1 +group type2 type 2 +800 atoms in group type2 +compute charge2 type2 property/atom q +compute q2 type2 reduce ave c_charge2 +variable qtot equal count(type1)*c_q1+count(type2)*c_q2 +variable nqeq equal f_2 + +thermo_style custom step pe c_q1 c_q2 v_qtot v_nqeq +thermo 10 +thermo_modify format line "%4d %12.9g %12.8f %12.8f %16.12f %6.0f" + +timestep 0.0001 + +fix 1 all nve + +run 100 +PPPM initialization ... + using 12-bit tables for long-range coulomb (src/kspace.cpp:339) + G vector (1/distance) = 0.27644401 + grid = 27 27 30 + stencil order = 5 + estimated absolute RMS force accuracy = 1.4502702e-05 + estimated relative force accuracy = 1.0071569e-06 + using double precision FFTW3 + 3d grid and FFT values/proc = 14960 5832 +Neighbor list info ... + update every 1 steps, delay 0 steps, check yes + max neighbors/atom: 2000, page size: 100000 + master list distance cutoff = 13 + ghost atom cutoff = 13 + binsize = 6.5, bins = 4 4 5 + 2 neighbor lists, perpetual/occasional/extra = 2 0 0 + (1) pair buck/coul/long, perpetual + attributes: half, newton on + pair build: half/bin/atomonly/newton + stencil: half/bin/3d/newton + bin: standard + (2) fix qeq/fire, perpetual, copy from (1) + attributes: half, newton on + pair build: copy + stencil: none + bin: none +Per MPI rank memory allocation (min/avg/max) = 9.195 | 9.246 | 9.278 Mbytes +Step PotEng c_q1 c_q2 v_qtot v_nqeq + 0 -3432.05316 0.85226679 -0.42613339 0.000000000001 37 + 10 -3452.04937 0.85475813 -0.42737906 0.000000000001 10 + 20 -3497.83659 0.85998739 -0.42999370 0.000000000002 13 + 30 -3568.47793 0.86771599 -0.43385799 0.000000000002 13 + 40 -3633.35326 0.87337029 -0.43668514 0.000000000002 13 + 50 -3700.16079 0.87805847 -0.43902924 0.000000000000 13 + 60 -3784.31906 0.88401635 -0.44200818 -0.000000000001 11 + 70 -3877.60163 0.89010162 -0.44505081 -0.000000000000 13 + 80 -3965.28179 0.89431356 -0.44715678 0.000000000001 11 + 90 -4048.33861 0.89698069 -0.44849034 0.000000000001 13 + 100 -4118.63861 0.89718818 -0.44859409 0.000000000002 12 +Loop time of 3.88026 on 4 procs for 100 steps with 1200 atoms + +Performance: 0.223 ns/day, 107.785 hours/ns, 25.771 timesteps/s +98.0% 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.68424 | 0.69912 | 0.73572 | 2.5 | 18.02 +Kspace | 0.38093 | 0.41715 | 0.43168 | 3.2 | 10.75 +Neigh | 0 | 0 | 0 | 0.0 | 0.00 +Comm | 0.012711 | 0.013318 | 0.014003 | 0.4 | 0.34 +Output | 0.00066566 | 0.00098735 | 0.0019317 | 0.0 | 0.03 +Modify | 2.7477 | 2.7477 | 2.7477 | 0.0 | 70.81 +Other | | 0.002004 | | | 0.05 + +Nlocal: 300.000 ave 300 max 300 min +Histogram: 4 0 0 0 0 0 0 0 0 0 +Nghost: 4875.00 ave 4880 max 4870 min +Histogram: 2 0 0 0 0 0 0 0 0 2 +Neighs: 91900.0 ave 93081 max 90719 min +Histogram: 2 0 0 0 0 0 0 0 0 2 + +Total # of neighbors = 367600 +Ave neighs/atom = 306.33333 +Neighbor list builds = 0 +Dangerous builds = 0 +Total wall time: 0:00:17 diff --git a/examples/qeq/log.27Nov18.qeq.buck.g++.1 b/examples/qeq/log.27Nov18.qeq.buck.g++.1 deleted file mode 100644 index 4d5225ccc3..0000000000 --- a/examples/qeq/log.27Nov18.qeq.buck.g++.1 +++ /dev/null @@ -1,118 +0,0 @@ -LAMMPS (27 Nov 2018) - using 1 OpenMP thread(s) per MPI task -# This example demonstrates the use of various fix qeq variants with -# that defines and uses charges, in this case pair_style buck/coul/long - -units metal -atom_style charge - -read_data data.aC - orthogonal box = (0 0 0) to (25.1583 25.1583 28.0203) - 1 by 1 by 1 MPI processor grid - reading atoms ... - 1200 atoms -replicate 2 2 2 - orthogonal box = (0 0 0) to (50.3166 50.3166 56.0405) - 1 by 1 by 1 MPI processor grid - 9600 atoms - Time spent = 0.00114894 secs - -pair_style buck/coul/long 12.0 -pair_coeff 2 2 1388.77 .3623188 175.0 -pair_coeff 1 2 18003 .2052124 133.5381 -pair_coeff 1 1 0 .1 0 -kspace_style ewald 1e-6 - -neighbor 1.0 bin -neigh_modify delay 0 every 1 check yes - -group type1 type 1 -3200 atoms in group type1 -compute charge1 type1 property/atom q -compute q1 type1 reduce ave c_charge1 -group type2 type 2 -6400 atoms in group type2 -compute charge2 type2 property/atom q -compute q2 type2 reduce ave c_charge2 -variable qtot equal count(type1)*c_q1+count(type2)*c_q2 - -thermo_style custom step pe c_q1 c_q2 v_qtot spcpu -thermo 10 - -timestep 0.0001 - -velocity all create 300.0 1281937 -fix 1 all nve - -#fix 2 all qeq/point 1 10 1.0e-6 100 param.qeq2 -#fix 2 all qeq/shielded 1 10 1.0e-6 100 param.qeq2 -#fix 2 all qeq/slater 1 10 1.0e-6 100 param.qeq2 -#fix 2 all qeq/dynamic 1 10 1.0e-4 100 param.qeq2 -fix 2 all qeq/fire 1 10 1.0e-4 100 param.qeq2 - -run 100 -Ewald initialization ... - using 12-bit tables for long-range coulomb (src/kspace.cpp:321) - G vector (1/distance) = 0.305064 - estimated absolute RMS force accuracy = 2.07629e-05 - estimated relative force accuracy = 1.44191e-06 - KSpace vectors: actual max1d max3d = 13556 20 34460 - kxmax kymax kzmax = 18 18 20 -Neighbor list info ... - update every 1 steps, delay 0 steps, check yes - max neighbors/atom: 2000, page size: 100000 - master list distance cutoff = 13 - ghost atom cutoff = 13 - binsize = 6.5, bins = 8 8 9 - 2 neighbor lists, perpetual/occasional/extra = 2 0 0 - (1) pair buck/coul/long, perpetual - attributes: half, newton on - pair build: half/bin/atomonly/newton - stencil: half/bin/3d/newton - bin: standard - (2) fix qeq/fire, perpetual, copy from (1) - attributes: half, newton on - pair build: copy - stencil: none - bin: none -Per MPI rank memory allocation (min/avg/max) = 134 | 134 | 134 Mbytes -Step PotEng c_q1 c_q2 v_qtot S/CPU - 0 -27457.219 0.85227886 -0.42613943 -2.1827873e-10 0 - 10 -27626.057 0.85486228 -0.42743114 -2.0372681e-10 0.64313877 - 20 -27975.085 0.85968531 -0.42984266 -1.036824e-10 0.55119179 - 30 -28552.628 0.86755661 -0.4337783 1.3051249e-10 0.53160643 - 40 -29133.643 0.87426387 -0.43713193 1.1368684e-10 0.53075341 - 50 -29697.011 0.8794039 -0.43970195 1.200533e-10 0.52358127 - 60 -30342.001 0.88478594 -0.44239297 6.002665e-11 0.5366762 - 70 -31081.138 0.8906973 -0.44534865 -4.7293724e-11 0.55904546 - 80 -31792.732 0.89506635 -0.44753317 -4.3200998e-11 0.59606079 - 90 -32424.749 0.89714841 -0.44857421 -1.1596057e-10 0.58047419 - 100 -32998.353 0.89755721 -0.44877861 -1.0231815e-10 0.59444001 -Loop time of 177.79 on 1 procs for 100 steps with 9600 atoms - -Performance: 0.005 ns/day, 4938.612 hours/ns, 0.562 timesteps/s -99.8% CPU use with 1 MPI tasks x 1 OpenMP threads - -MPI task timing breakdown: -Section | min time | avg time | max time |%varavg| %total ---------------------------------------------------------------- -Pair | 11.518 | 11.518 | 11.518 | 0.0 | 6.48 -Kspace | 107.37 | 107.37 | 107.37 | 0.0 | 60.39 -Neigh | 0 | 0 | 0 | 0.0 | 0.00 -Comm | 0.019721 | 0.019721 | 0.019721 | 0.0 | 0.01 -Output | 0.002218 | 0.002218 | 0.002218 | 0.0 | 0.00 -Modify | 58.869 | 58.869 | 58.869 | 0.0 | 33.11 -Other | | 0.007197 | | | 0.00 - -Nlocal: 9600 ave 9600 max 9600 min -Histogram: 1 0 0 0 0 0 0 0 0 0 -Nghost: 22125 ave 22125 max 22125 min -Histogram: 1 0 0 0 0 0 0 0 0 0 -Neighs: 2.9408e+06 ave 2.9408e+06 max 2.9408e+06 min -Histogram: 1 0 0 0 0 0 0 0 0 0 - -Total # of neighbors = 2940800 -Ave neighs/atom = 306.333 -Neighbor list builds = 0 -Dangerous builds = 0 -Total wall time: 0:03:01 diff --git a/examples/qeq/log.27Nov18.qeq.buck.g++.4 b/examples/qeq/log.27Nov18.qeq.buck.g++.4 deleted file mode 100644 index 947c3caeaf..0000000000 --- a/examples/qeq/log.27Nov18.qeq.buck.g++.4 +++ /dev/null @@ -1,118 +0,0 @@ -LAMMPS (27 Nov 2018) - using 1 OpenMP thread(s) per MPI task -# This example demonstrates the use of various fix qeq variants with -# that defines and uses charges, in this case pair_style buck/coul/long - -units metal -atom_style charge - -read_data data.aC - orthogonal box = (0 0 0) to (25.1583 25.1583 28.0203) - 1 by 2 by 2 MPI processor grid - reading atoms ... - 1200 atoms -replicate 2 2 2 - orthogonal box = (0 0 0) to (50.3166 50.3166 56.0405) - 1 by 2 by 2 MPI processor grid - 9600 atoms - Time spent = 0.000675201 secs - -pair_style buck/coul/long 12.0 -pair_coeff 2 2 1388.77 .3623188 175.0 -pair_coeff 1 2 18003 .2052124 133.5381 -pair_coeff 1 1 0 .1 0 -kspace_style ewald 1e-6 - -neighbor 1.0 bin -neigh_modify delay 0 every 1 check yes - -group type1 type 1 -3200 atoms in group type1 -compute charge1 type1 property/atom q -compute q1 type1 reduce ave c_charge1 -group type2 type 2 -6400 atoms in group type2 -compute charge2 type2 property/atom q -compute q2 type2 reduce ave c_charge2 -variable qtot equal count(type1)*c_q1+count(type2)*c_q2 - -thermo_style custom step pe c_q1 c_q2 v_qtot spcpu -thermo 10 - -timestep 0.0001 - -velocity all create 300.0 1281937 -fix 1 all nve - -#fix 2 all qeq/point 1 10 1.0e-6 100 param.qeq2 -#fix 2 all qeq/shielded 1 10 1.0e-6 100 param.qeq2 -#fix 2 all qeq/slater 1 10 1.0e-6 100 param.qeq2 -#fix 2 all qeq/dynamic 1 10 1.0e-4 100 param.qeq2 -fix 2 all qeq/fire 1 10 1.0e-4 100 param.qeq2 - -run 100 -Ewald initialization ... - using 12-bit tables for long-range coulomb (src/kspace.cpp:321) - G vector (1/distance) = 0.305064 - estimated absolute RMS force accuracy = 2.07629e-05 - estimated relative force accuracy = 1.44191e-06 - KSpace vectors: actual max1d max3d = 13556 20 34460 - kxmax kymax kzmax = 18 18 20 -Neighbor list info ... - update every 1 steps, delay 0 steps, check yes - max neighbors/atom: 2000, page size: 100000 - master list distance cutoff = 13 - ghost atom cutoff = 13 - binsize = 6.5, bins = 8 8 9 - 2 neighbor lists, perpetual/occasional/extra = 2 0 0 - (1) pair buck/coul/long, perpetual - attributes: half, newton on - pair build: half/bin/atomonly/newton - stencil: half/bin/3d/newton - bin: standard - (2) fix qeq/fire, perpetual, copy from (1) - attributes: half, newton on - pair build: copy - stencil: none - bin: none -Per MPI rank memory allocation (min/avg/max) = 53.06 | 53.13 | 53.21 Mbytes -Step PotEng c_q1 c_q2 v_qtot S/CPU - 0 -27457.215 0.85227886 -0.42613943 2.1373125e-11 0 - 10 -27626.057 0.85486228 -0.42743114 3.0468073e-11 2.4245312 - 20 -27975.085 0.85968531 -0.42984266 1.0095391e-10 2.0185316 - 30 -28552.627 0.86755661 -0.4337783 1.3096724e-10 1.9605335 - 40 -29133.643 0.87426387 -0.43713193 1.5279511e-10 1.9624139 - 50 -29697.01 0.8794039 -0.43970195 1.6461854e-10 1.8113263 - 60 -30342 0.88478594 -0.44239297 1.7826096e-10 1.9537722 - 70 -31081.139 0.89069733 -0.44534866 1.4733814e-10 2.058406 - 80 -31792.732 0.89506635 -0.44753317 1.3824319e-10 2.2160813 - 90 -32424.752 0.89714841 -0.44857421 1.2914825e-10 2.0952145 - 100 -32998.353 0.89755721 -0.44877861 1.4824764e-10 2.1292486 -Loop time of 48.7541 on 4 procs for 100 steps with 9600 atoms - -Performance: 0.018 ns/day, 1354.281 hours/ns, 2.051 timesteps/s -97.4% CPU use with 4 MPI tasks x 1 OpenMP threads - -MPI task timing breakdown: -Section | min time | avg time | max time |%varavg| %total ---------------------------------------------------------------- -Pair | 2.9747 | 3.0315 | 3.0758 | 2.1 | 6.22 -Kspace | 27.873 | 28.264 | 28.63 | 5.3 | 57.97 -Neigh | 0 | 0 | 0 | 0.0 | 0.00 -Comm | 0.53835 | 0.8523 | 1.2286 | 28.2 | 1.75 -Output | 0.0012984 | 0.001591 | 0.0024178 | 1.2 | 0.00 -Modify | 16.58 | 16.59 | 16.601 | 0.3 | 34.03 -Other | | 0.01409 | | | 0.03 - -Nlocal: 2400 ave 2400 max 2400 min -Histogram: 4 0 0 0 0 0 0 0 0 0 -Nghost: 11550 ave 11550 max 11550 min -Histogram: 4 0 0 0 0 0 0 0 0 0 -Neighs: 735200 ave 740758 max 729642 min -Histogram: 2 0 0 0 0 0 0 0 0 2 - -Total # of neighbors = 2940800 -Ave neighs/atom = 306.333 -Neighbor list builds = 0 -Dangerous builds = 0 -Total wall time: 0:00:49 diff --git a/src/QEQ/fix_qeq.cpp b/src/QEQ/fix_qeq.cpp index 4504865303..56491228fd 100644 --- a/src/QEQ/fix_qeq.cpp +++ b/src/QEQ/fix_qeq.cpp @@ -54,10 +54,15 @@ FixQEq::FixQEq(LAMMPS *lmp, int narg, char **arg) : { if (narg < 8) error->all(FLERR,"Illegal fix qeq command"); + scalar_flag = 1; + extscalar = 0; + nevery = utils::inumeric(FLERR,arg[3],false,lmp); cutoff = utils::numeric(FLERR,arg[4],false,lmp); tolerance = utils::numeric(FLERR,arg[5],false,lmp); maxiter = utils::inumeric(FLERR,arg[6],false,lmp); + maxwarn = 1; + matvecs = 0; // check for sane arguments if ((nevery <= 0) || (cutoff <= 0.0) || (tolerance <= 0.0) || (maxiter <= 0)) @@ -123,7 +128,6 @@ FixQEq::FixQEq(LAMMPS *lmp, int narg, char **arg) : } else { read_file(arg[7]); } - } /* ---------------------------------------------------------------------- */ @@ -277,6 +281,13 @@ void FixQEq::reallocate_matrix() /* ---------------------------------------------------------------------- */ +double FixQEq::compute_scalar() +{ + return matvecs; +} + +/* ---------------------------------------------------------------------- */ + void FixQEq::init_list(int /*id*/, NeighList *ptr) { list = ptr; @@ -395,7 +406,7 @@ int FixQEq::CG( double *b, double *x ) vector_sum( d, 1., p, beta, d, inum ); } - if ((comm->me == 0) && (loop >= maxiter)) + if ((comm->me == 0) && maxwarn && (loop >= maxiter)) error->warning(FLERR,fmt::format("Fix qeq CG convergence failed ({}) " "after {} iterations at step {}", sqrt(sig_new)/b_norm,loop, diff --git a/src/QEQ/fix_qeq.h b/src/QEQ/fix_qeq.h index 9a3087840b..a11a41841e 100644 --- a/src/QEQ/fix_qeq.h +++ b/src/QEQ/fix_qeq.h @@ -35,6 +35,8 @@ class FixQEq : public Fix { void pre_force_respa(int, int, int); void min_pre_force(int); + virtual double compute_scalar(); + // derived child classes must provide these functions virtual void init() = 0; @@ -64,7 +66,8 @@ class FixQEq : public Fix { double swa, swb; // lower/upper Taper cutoff radius double Tap[8]; // Taper function double tolerance; // tolerance for the norm of the rel residual in CG - int maxiter; // maximum number of QEq iterations + int maxiter; // maximum number of QEq iterations + int maxwarn; // print warning when max iterations was reached double cutoff, cutoff_sq; // neighbor cutoff double *chi,*eta,*gamma,*zeta,*zcore; // qeq parameters diff --git a/src/QEQ/fix_qeq_dynamic.cpp b/src/QEQ/fix_qeq_dynamic.cpp index 397393b786..5af7a4c9a0 100644 --- a/src/QEQ/fix_qeq_dynamic.cpp +++ b/src/QEQ/fix_qeq_dynamic.cpp @@ -53,6 +53,12 @@ FixQEqDynamic::FixQEqDynamic(LAMMPS *lmp, int narg, char **arg) : if (iarg+2 > narg) error->all(FLERR,"Illegal fix qeq/dynamic command"); qstep = atof(arg[iarg+1]); iarg += 2; + } else if (strcmp(arg[iarg],"warn") == 0) { + if (iarg+2 > narg) error->all(FLERR,"Illegal fix qeq/dynamic command"); + if (strcmp(arg[iarg+1],"no") == 0) maxwarn = 0; + else if (strcmp(arg[iarg+1],"yes") == 0) maxwarn = 1; + else error->all(FLERR,"Illegal fix qeq/dynamic command"); + iarg += 2; } else error->all(FLERR,"Illegal fix qeq/dynamic command"); } } @@ -145,7 +151,7 @@ void FixQEqDynamic::pre_force(int /*vflag*/) MPI_Allreduce(&enegmax,&enegmaxall,1,MPI_DOUBLE,MPI_MAX,world); enegmax = enegmaxall; - if (enegchk <= tolerance && enegmax <= 100.0*tolerance) break; + if ((enegchk <= tolerance) && (enegmax <= 100.0*tolerance)) break; for (ii = 0; ii < inum; ii++) { i = ilist[ii]; @@ -153,8 +159,9 @@ void FixQEqDynamic::pre_force(int /*vflag*/) q1[i] += qf[i]*dtq2 - qdamp*q1[i]; } } + matvecs = iloop; - if ((comm->me == 0) && (iloop >= maxiter)) + if ((comm->me == 0) && maxwarn && (iloop >= maxiter)) error->warning(FLERR,fmt::format("Charges did not converge at step " "{}: {}",update->ntimestep,enegchk)); diff --git a/src/QEQ/fix_qeq_fire.cpp b/src/QEQ/fix_qeq_fire.cpp index 0bdf65dc18..0e89ee3d17 100644 --- a/src/QEQ/fix_qeq_fire.cpp +++ b/src/QEQ/fix_qeq_fire.cpp @@ -63,6 +63,12 @@ FixQEqFire::FixQEqFire(LAMMPS *lmp, int narg, char **arg) : if (iarg+2 > narg) error->all(FLERR,"Illegal fix qeq/fire command"); qstep = atof(arg[iarg+1]); iarg += 2; + } else if (strcmp(arg[iarg],"warn") == 0) { + if (iarg+2 > narg) error->all(FLERR,"Illegal fix qeq/fire command"); + if (strcmp(arg[iarg+1],"no") == 0) maxwarn = 0; + else if (strcmp(arg[iarg+1],"yes") == 0) maxwarn = 1; + else error->all(FLERR,"Illegal fix qeq/fire command"); + iarg += 2; } else error->all(FLERR,"Illegal fix qeq/fire command"); } } @@ -213,8 +219,9 @@ void FixQEqFire::pre_force(int /*vflag*/) if (enegchk < tolerance) break; } + matvecs = iloop; - if ((comm->me == 0) && (iloop >= maxiter)) + if ((comm->me == 0) && maxwarn && (iloop >= maxiter)) error->warning(FLERR,fmt::format("Charges did not converge at step " "{}: {}",update->ntimestep,enegchk)); diff --git a/src/QEQ/fix_qeq_point.cpp b/src/QEQ/fix_qeq_point.cpp index ac31f906e0..8d60585773 100644 --- a/src/QEQ/fix_qeq_point.cpp +++ b/src/QEQ/fix_qeq_point.cpp @@ -38,7 +38,15 @@ using namespace LAMMPS_NS; /* ---------------------------------------------------------------------- */ FixQEqPoint::FixQEqPoint(LAMMPS *lmp, int narg, char **arg) : - FixQEq(lmp, narg, arg) {} + FixQEq(lmp, narg, arg) { + if (narg == 10) { + if (strcmp(arg[8],"warn") == 0) { + if (strcmp(arg[9],"no") == 0) maxwarn = 0; + else if (strcmp(arg[9],"yes") == 0) maxwarn = 1; + else error->all(FLERR,"Illegal fix qeq/point command"); + } else error->all(FLERR,"Illegal fix qeq/point command"); + } else if (narg > 8) error->all(FLERR,"Illegal fix qeq/point command"); +} /* ---------------------------------------------------------------------- */ @@ -80,6 +88,7 @@ void FixQEqPoint::pre_force(int /*vflag*/) init_matvec(); matvecs = CG(b_s, s); // CG on s - parallel matvecs += CG(b_t, t); // CG on t - parallel + matvecs /= 2; calculate_Q(); if (force->kspace) force->kspace->qsum_qsq(); diff --git a/src/QEQ/fix_qeq_shielded.cpp b/src/QEQ/fix_qeq_shielded.cpp index ad6202abd8..351700057b 100644 --- a/src/QEQ/fix_qeq_shielded.cpp +++ b/src/QEQ/fix_qeq_shielded.cpp @@ -40,6 +40,13 @@ using namespace LAMMPS_NS; FixQEqShielded::FixQEqShielded(LAMMPS *lmp, int narg, char **arg) : FixQEq(lmp, narg, arg) { + if (narg == 10) { + if (strcmp(arg[8],"warn") == 0) { + if (strcmp(arg[9],"no") == 0) maxwarn = 0; + else if (strcmp(arg[9],"yes") == 0) maxwarn = 1; + else error->all(FLERR,"Illegal fix qeq/shielded command"); + } else error->all(FLERR,"Illegal fix qeq/shielded command"); + } else if (narg > 8) error->all(FLERR,"Illegal fix qeq/shielded command"); if (reax_flag) extract_reax(); } @@ -143,6 +150,7 @@ void FixQEqShielded::pre_force(int /*vflag*/) init_matvec(); matvecs = CG(b_s, s); // CG on s - parallel matvecs += CG(b_t, t); // CG on t - parallel + matvecs /= 2; calculate_Q(); if (force->kspace) force->kspace->qsum_qsq(); diff --git a/src/QEQ/fix_qeq_slater.cpp b/src/QEQ/fix_qeq_slater.cpp index 326d71c93b..da0d2090cf 100644 --- a/src/QEQ/fix_qeq_slater.cpp +++ b/src/QEQ/fix_qeq_slater.cpp @@ -52,6 +52,12 @@ FixQEqSlater::FixQEqSlater(LAMMPS *lmp, int narg, char **arg) : if (iarg+2 > narg) error->all(FLERR,"Illegal fix qeq/slater command"); alpha = atof(arg[iarg+1]); iarg += 2; + } else if (strcmp(arg[iarg],"warn") == 0) { + if (iarg+2 > narg) error->all(FLERR,"Illegal fix qeq/slater command"); + if (strcmp(arg[iarg+1],"no") == 0) maxwarn = 0; + else if (strcmp(arg[iarg+1],"yes") == 0) maxwarn = 1; + else error->all(FLERR,"Illegal fix qeq/slater command"); + iarg += 2; } else error->all(FLERR,"Illegal fix qeq/slater command"); } @@ -120,6 +126,7 @@ void FixQEqSlater::pre_force(int /*vflag*/) init_matvec(); matvecs = CG(b_s, s); // CG on s - parallel matvecs += CG(b_t, t); // CG on t - parallel + matvecs /= 2; calculate_Q(); if (force->kspace) force->kspace->qsum_qsq();