2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015,2016,2017 by the GROMACS development team.
7 * Copyright (c) 2018,2019,2020,2021, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
47 #include "gromacs/commandline/pargs.h"
48 #include "gromacs/commandline/viewit.h"
49 #include "gromacs/fileio/confio.h"
50 #include "gromacs/fileio/trxio.h"
51 #include "gromacs/fileio/xvgr.h"
52 #include "gromacs/gmxana/gmx_ana.h"
53 #include "gromacs/linearalgebra/nrjac.h"
54 #include "gromacs/math/functions.h"
55 #include "gromacs/math/units.h"
56 #include "gromacs/math/vec.h"
57 #include "gromacs/mdtypes/md_enums.h"
58 #include "gromacs/pbcutil/pbc.h"
59 #include "gromacs/pbcutil/rmpbc.h"
60 #include "gromacs/topology/index.h"
61 #include "gromacs/topology/topology.h"
62 #include "gromacs/trajectory/trajectoryframe.h"
63 #include "gromacs/utility/arraysize.h"
64 #include "gromacs/utility/cstringutil.h"
65 #include "gromacs/utility/fatalerror.h"
66 #include "gromacs/utility/futil.h"
67 #include "gromacs/utility/smalloc.h"
68 #include "gromacs/utility/stringutil.h"
71 low_print_data(FILE* fp, real time, rvec x[], int n, const int* index, const gmx_bool bDim[], const char* sffmt)
75 fprintf(fp, " %g", time);
76 for (i = 0; i < n; i++)
86 for (d = 0; d < DIM; d++)
90 fprintf(fp, sffmt, x[ii][d]);
95 fprintf(fp, sffmt, norm(x[ii]));
101 static void average_data(rvec x[], rvec xav[], const real* mass, int ngrps, const int isize[], int** index)
106 double sum[DIM], mtot;
108 for (g = 0; g < ngrps; g++)
113 for (i = 0; i < isize[g]; i++)
119 svmul(m, x[ind], tmp);
120 for (d = 0; d < DIM; d++)
128 for (d = 0; d < DIM; d++)
136 for (d = 0; d < DIM; d++)
138 xav[g][d] = sum[d] / mtot;
143 /* mass=NULL, so these are forces: sum only (do not average) */
144 for (d = 0; d < DIM; d++)
152 static void print_data(FILE* fp,
163 static rvec* xav = nullptr;
171 average_data(x, xav, mass, ngrps, isize, index);
172 low_print_data(fp, time, xav, ngrps, nullptr, bDim, sffmt);
176 low_print_data(fp, time, x, isize[0], index[0], bDim, sffmt);
180 static void write_trx_x(t_trxstatus* status,
181 const t_trxframe* fr,
188 static rvec* xav = nullptr;
189 static t_atoms* atoms = nullptr;
200 snew(atoms->atom, ngrps);
202 /* Note that atom and residue names will be the ones
203 * of the first atom in each group.
205 for (i = 0; i < ngrps; i++)
207 atoms->atom[i] = fr->atoms->atom[index[i][0]];
208 atoms->atomname[i] = fr->atoms->atomname[index[i][0]];
211 average_data(fr->x, xav, mass, ngrps, isize, index);
213 fr_av.natoms = ngrps;
216 write_trxframe(status, &fr_av, nullptr);
220 write_trxframe_indexed(status, fr, isize[0], index[0], nullptr);
224 static void make_legend(FILE* fp,
231 const gmx_bool bDim[],
232 const gmx_output_env_t* oenv)
235 const char* dimtxt[] = { " X", " Y", " Z", "" };
249 for (i = 0; i < n; i++)
251 for (d = 0; d <= DIM; d++)
255 snew(leg[j], STRLEN);
258 sprintf(leg[j], "mol %d%s", index[i] + 1, dimtxt[d]);
262 sprintf(leg[j], "%s%s", name[i], dimtxt[d]);
266 sprintf(leg[j], "atom %d%s", index[i] + 1, dimtxt[d]);
272 xvgr_legend(fp, j, leg, oenv);
274 for (i = 0; i < j; i++)
281 static real ekrot(rvec x[], rvec v[], const real mass[], int isize, const int index[])
283 real TCM[5][5], L[5][5];
284 double tm, m0, lxx, lxy, lxz, lyy, lyz, lzz, ekrot;
294 for (i = 0; i < isize; i++)
299 cprod(x[j], v[j], a0);
300 for (m = 0; (m < DIM); m++)
302 xcm[m] += m0 * x[j][m]; /* c.o.m. position */
303 vcm[m] += m0 * v[j][m]; /* c.o.m. velocity */
304 acm[m] += m0 * a0[m]; /* rotational velocity around c.o.m. */
307 dcprod(xcm, vcm, b0);
308 for (m = 0; (m < DIM); m++)
312 acm[m] -= b0[m] / tm;
315 lxx = lxy = lxz = lyy = lyz = lzz = 0.0;
316 for (i = 0; i < isize; i++)
320 for (m = 0; m < DIM; m++)
322 dx[m] = x[j][m] - xcm[m];
324 lxx += dx[XX] * dx[XX] * m0;
325 lxy += dx[XX] * dx[YY] * m0;
326 lxz += dx[XX] * dx[ZZ] * m0;
327 lyy += dx[YY] * dx[YY] * m0;
328 lyz += dx[YY] * dx[ZZ] * m0;
329 lzz += dx[ZZ] * dx[ZZ] * m0;
332 L[XX][XX] = lyy + lzz;
336 L[YY][YY] = lxx + lzz;
340 L[ZZ][ZZ] = lxx + lyy;
342 m_inv_gen(&L[0][0], DIM, &TCM[0][0]);
344 /* Compute omega (hoeksnelheid) */
347 for (m = 0; m < DIM; m++)
349 for (n = 0; n < DIM; n++)
351 ocm[m] += TCM[m][n] * acm[n];
353 ekrot += 0.5 * ocm[m] * acm[m];
359 static real ektrans(rvec v[], const real mass[], int isize, const int index[])
367 for (i = 0; i < isize; i++)
370 for (d = 0; d < DIM; d++)
372 mvcom[d] += mass[j] * v[j][d];
377 return dnorm2(mvcom) / (mtot * 2);
380 static real temp(rvec v[], const real mass[], int isize, const int index[])
386 for (i = 0; i < isize; i++)
389 ekin2 += mass[j] * norm2(v[j]);
392 return ekin2 / (3 * isize * gmx::c_boltz);
395 static void remove_jump(matrix box, int natoms, rvec xp[], rvec x[])
400 for (d = 0; d < DIM; d++)
402 hbox[d] = 0.5 * box[d][d];
404 for (i = 0; i < natoms; i++)
406 for (m = DIM - 1; m >= 0; m--)
408 while (x[i][m] - xp[i][m] <= -hbox[m])
410 for (d = 0; d <= m; d++)
412 x[i][d] += box[m][d];
415 while (x[i][m] - xp[i][m] > hbox[m])
417 for (d = 0; d <= m; d++)
419 x[i][d] -= box[m][d];
426 static void write_pdb_bfac(const char* fname,
438 const gmx_bool bDim[],
440 const gmx_output_env_t* oenv)
443 real max, len2, scale;
447 if ((nfr_x == 0) || (nfr_v == 0))
449 fprintf(stderr, "No frames found for %s, will not write %s\n", title, fname);
453 fprintf(stderr, "Used %d frames for %s\n", nfr_x, "coordinates");
454 fprintf(stderr, "Used %d frames for %s\n", nfr_v, title);
459 for (i = 0; i < DIM; i++)
473 for (i = 0; i < isize; i++)
475 svmul(scale, sum[index[i]], sum[index[i]]);
478 fp = xvgropen(xname, title, "Atom", "Spatial component", oenv);
479 for (i = 0; i < isize; i++)
482 "%-5d %10.3f %10.3f %10.3f\n",
491 for (i = 0; i < isize; i++)
494 for (m = 0; m < DIM; m++)
496 if (bDim[m] || bDim[DIM])
498 len2 += gmx::square(sum[index[i]][m]);
507 if (scale_factor != 0)
509 scale = scale_factor;
519 scale = 10.0 / std::sqrt(max);
523 printf("Maximum %s is %g on atom %d %s, res. %s %d\n",
527 *(atoms->atomname[maxi]),
528 *(atoms->resinfo[atoms->atom[maxi].resind].name),
529 atoms->resinfo[atoms->atom[maxi].resind].nr);
531 if (atoms->pdbinfo == nullptr)
533 snew(atoms->pdbinfo, atoms->nr);
535 atoms->havePdbInfo = TRUE;
539 for (i = 0; i < isize; i++)
542 for (m = 0; m < DIM; m++)
544 if (bDim[m] || bDim[DIM])
546 len2 += gmx::square(sum[index[i]][m]);
549 atoms->pdbinfo[index[i]].bfac = std::sqrt(len2) * scale;
554 for (i = 0; i < isize; i++)
556 atoms->pdbinfo[index[i]].bfac = sum[index[i]][onedim] * scale;
559 write_sto_conf_indexed(fname, title, atoms, x, nullptr, pbcType, box, isize, index);
563 static void update_histo(int gnx, const int index[], rvec v[], int* nhisto, int** histo, real binwidth)
568 if (*histo == nullptr)
571 for (i = 0; (i < gnx); i++)
573 vn = norm(v[index[i]]);
574 vnmax = std::max(vn, vnmax);
577 *nhisto = static_cast<int>(1 + (vnmax / binwidth));
578 snew(*histo, *nhisto);
580 for (i = 0; (i < gnx); i++)
582 vn = norm(v[index[i]]);
583 in = static_cast<int>(vn / binwidth);
587 fprintf(stderr, "Extending histogram from %d to %d\n", *nhisto, nnn);
590 for (m = *nhisto; (m < nnn); m++)
600 static void print_histo(const char* fn, int nhisto, int histo[], real binwidth, const gmx_output_env_t* oenv)
605 fp = xvgropen(fn, "Velocity distribution", "V (nm/ps)", "arbitrary units", oenv);
606 for (i = 0; (i < nhisto); i++)
608 fprintf(fp, "%10.3e %10d\n", i * binwidth, histo[i]);
613 int gmx_traj(int argc, char* argv[])
615 const char* desc[] = {
616 "[THISMODULE] plots coordinates, velocities, forces and/or the box.",
617 "With [TT]-com[tt] the coordinates, velocities and forces are",
618 "calculated for the center of mass of each group.",
619 "When [TT]-mol[tt] is set, the numbers in the index file are",
620 "interpreted as molecule numbers and the same procedure as with",
621 "[TT]-com[tt] is used for each molecule.[PAR]",
622 "Option [TT]-ot[tt] plots the temperature of each group,",
623 "provided velocities are present in the trajectory file.",
624 "No corrections are made for constrained degrees of freedom!",
625 "This implies [TT]-com[tt].[PAR]",
626 "Options [TT]-ekt[tt] and [TT]-ekr[tt] plot the translational and",
627 "rotational kinetic energy of each group,",
628 "provided velocities are present in the trajectory file.",
629 "This implies [TT]-com[tt].[PAR]",
630 "Options [TT]-cv[tt] and [TT]-cf[tt] write the average velocities",
631 "and average forces as temperature factors to a [REF].pdb[ref] file with",
632 "the average coordinates or the coordinates at [TT]-ctime[tt].",
633 "The temperature factors are scaled such that the maximum is 10.",
634 "The scaling can be changed with the option [TT]-scale[tt].",
635 "To get the velocities or forces of one",
636 "frame set both [TT]-b[tt] and [TT]-e[tt] to the time of",
637 "desired frame. When averaging over frames you might need to use",
638 "the [TT]-nojump[tt] option to obtain the correct average coordinates.",
639 "If you select either of these option the average force and velocity",
640 "for each atom are written to an [REF].xvg[ref] file as well",
641 "(specified with [TT]-av[tt] or [TT]-af[tt]).[PAR]",
642 "Option [TT]-vd[tt] computes a velocity distribution, i.e. the",
643 "norm of the vector is plotted. In addition in the same graph",
644 "the kinetic energy distribution is given.",
646 "See [gmx-trajectory] for plotting similar data for selections."
648 static gmx_bool bMol = FALSE, bCom = FALSE, bPBC = TRUE, bNoJump = FALSE;
649 static gmx_bool bX = TRUE, bY = TRUE, bZ = TRUE, bNorm = FALSE, bFP = FALSE;
650 static int ngroups = 1;
651 static real ctime = -1, scale = 0, binwidth = 1;
653 { "-com", FALSE, etBOOL, { &bCom }, "Plot data for the com of each group" },
654 { "-pbc", FALSE, etBOOL, { &bPBC }, "Make molecules whole for COM" },
659 "Index contains molecule numbers instead of atom numbers" },
660 { "-nojump", FALSE, etBOOL, { &bNoJump }, "Remove jumps of atoms across the box" },
661 { "-x", FALSE, etBOOL, { &bX }, "Plot X-component" },
662 { "-y", FALSE, etBOOL, { &bY }, "Plot Y-component" },
663 { "-z", FALSE, etBOOL, { &bZ }, "Plot Z-component" },
664 { "-ng", FALSE, etINT, { &ngroups }, "Number of groups to consider" },
665 { "-len", FALSE, etBOOL, { &bNorm }, "Plot vector length" },
666 { "-fp", FALSE, etBOOL, { &bFP }, "Full precision output" },
667 { "-bin", FALSE, etREAL, { &binwidth }, "Binwidth for velocity histogram (nm/ps)" },
672 "Use frame at this time for x in [TT]-cv[tt] and [TT]-cf[tt] instead of the average x" },
677 "Scale factor for [REF].pdb[ref] output, 0 is autoscale" }
679 FILE * outx = nullptr, *outv = nullptr, *outf = nullptr, *outb = nullptr, *outt = nullptr;
680 FILE * outekt = nullptr, *outekr = nullptr;
686 int flags, nvhisto = 0, *vhisto = nullptr;
687 rvec * xtop, *xp = nullptr;
688 rvec * sumx = nullptr, *sumv = nullptr, *sumf = nullptr;
691 t_trxstatus* status_out = nullptr;
692 gmx_rmpbc_t gpbc = nullptr;
694 int nr_xfr, nr_vfr, nr_ffr;
696 int * isize0, *isize;
697 int ** index0, **index;
700 gmx_bool bTop, bOX, bOXT, bOV, bOF, bOB, bOT, bEKT, bEKR, bCV, bCF;
701 gmx_bool bDim[4], bDum[4], bVD;
703 const char* box_leg[6] = { "XX", "YY", "ZZ", "YX", "ZX", "ZY" };
704 gmx_output_env_t* oenv;
707 { efTRX, "-f", nullptr, ffREAD }, { efTPS, nullptr, nullptr, ffREAD },
708 { efNDX, nullptr, nullptr, ffOPTRD }, { efXVG, "-ox", "coord", ffOPTWR },
709 { efTRX, "-oxt", "coord", ffOPTWR }, { efXVG, "-ov", "veloc", ffOPTWR },
710 { efXVG, "-of", "force", ffOPTWR }, { efXVG, "-ob", "box", ffOPTWR },
711 { efXVG, "-ot", "temp", ffOPTWR }, { efXVG, "-ekt", "ektrans", ffOPTWR },
712 { efXVG, "-ekr", "ekrot", ffOPTWR }, { efXVG, "-vd", "veldist", ffOPTWR },
713 { efPDB, "-cv", "veloc", ffOPTWR }, { efPDB, "-cf", "force", ffOPTWR },
714 { efXVG, "-av", "all_veloc", ffOPTWR }, { efXVG, "-af", "all_force", ffOPTWR }
716 #define NFILE asize(fnm)
718 if (!parse_common_args(&argc,
720 PCA_CAN_TIME | PCA_TIME_UNIT | PCA_CAN_VIEW,
737 "Interpreting indexfile entries as molecules.\n"
738 "Using center of mass.\n");
741 bOX = opt2bSet("-ox", NFILE, fnm);
742 bOXT = opt2bSet("-oxt", NFILE, fnm);
743 bOV = opt2bSet("-ov", NFILE, fnm);
744 bOF = opt2bSet("-of", NFILE, fnm);
745 bOB = opt2bSet("-ob", NFILE, fnm);
746 bOT = opt2bSet("-ot", NFILE, fnm);
747 bEKT = opt2bSet("-ekt", NFILE, fnm);
748 bEKR = opt2bSet("-ekr", NFILE, fnm);
749 bCV = opt2bSet("-cv", NFILE, fnm) || opt2bSet("-av", NFILE, fnm);
750 bCF = opt2bSet("-cf", NFILE, fnm) || opt2bSet("-af", NFILE, fnm);
751 bVD = opt2bSet("-vd", NFILE, fnm) || opt2parg_bSet("-bin", asize(pa), pa);
752 if (bMol || bOT || bEKT || bEKR)
764 sprintf(sffmt, "\t%s", gmx_real_fullprecision_pfmt);
768 sprintf(sffmt, "\t%%g");
770 std::string sffmt6 = gmx::formatString("%s%s%s%s%s%s", sffmt, sffmt, sffmt, sffmt, sffmt, sffmt);
772 bTop = read_tps_conf(ftp2fn(efTPS, NFILE, fnm),
778 bCom && (bOX || bOXT || bOV || bOT || bEKT || bEKR));
780 if ((bMol || bCV || bCF) && !bTop)
782 gmx_fatal(FARGS, "Need a run input file for option -mol, -cv or -cf");
787 indexfn = ftp2fn(efNDX, NFILE, fnm);
791 indexfn = ftp2fn_null(efNDX, NFILE, fnm);
794 if (!(bCom && !bMol))
798 snew(grpname, ngroups);
799 snew(isize0, ngroups);
800 snew(index0, ngroups);
801 get_index(&(top.atoms), indexfn, ngroups, isize0, index0, grpname);
808 snew(isize, ngroups);
809 snew(index, ngroups);
810 for (i = 0; i < ngroups; i++)
812 if (index0[0][i] < 0 || index0[0][i] >= mols->nr)
814 gmx_fatal(FARGS, "Molecule index (%d) is out of range (%d-%d)", index0[0][i] + 1, 1, mols->nr);
816 isize[i] = atndx[index0[0][i] + 1] - atndx[index0[0][i]];
817 snew(index[i], isize[i]);
818 for (j = 0; j < isize[i]; j++)
820 index[i][j] = atndx[index0[0][i]] + j;
831 snew(mass, top.atoms.nr);
832 for (i = 0; i < top.atoms.nr; i++)
834 mass[i] = top.atoms.atom[i].m;
843 std::string label(output_env_get_xvgr_tlabel(oenv));
846 flags = flags | TRX_READ_X;
847 outx = xvgropen(opt2fn("-ox", NFILE, fnm),
848 bCom ? "Center of mass" : "Coordinate",
852 make_legend(outx, ngroups, isize0[0], index0[0], grpname, bCom, bMol, bDim, oenv);
856 flags = flags | TRX_READ_X;
857 status_out = open_trx(opt2fn("-oxt", NFILE, fnm), "w");
861 flags = flags | TRX_READ_V;
862 outv = xvgropen(opt2fn("-ov", NFILE, fnm),
863 bCom ? "Center of mass velocity" : "Velocity",
867 make_legend(outv, ngroups, isize0[0], index0[0], grpname, bCom, bMol, bDim, oenv);
871 flags = flags | TRX_READ_F;
873 opt2fn("-of", NFILE, fnm), "Force", label, "Force (kJ mol\\S-1\\N nm\\S-1\\N)", oenv);
874 make_legend(outf, ngroups, isize0[0], index0[0], grpname, bCom, bMol, bDim, oenv);
878 outb = xvgropen(opt2fn("-ob", NFILE, fnm), "Box vector elements", label, "(nm)", oenv);
880 xvgr_legend(outb, 6, box_leg, oenv);
888 flags = flags | TRX_READ_V;
889 outt = xvgropen(opt2fn("-ot", NFILE, fnm), "Temperature", label, "(K)", oenv);
890 make_legend(outt, ngroups, isize[0], index[0], grpname, bCom, bMol, bDum, oenv);
898 flags = flags | TRX_READ_V;
899 outekt = xvgropen(opt2fn("-ekt", NFILE, fnm),
900 "Center of mass translation",
902 "Energy (kJ mol\\S-1\\N)",
904 make_legend(outekt, ngroups, isize[0], index[0], grpname, bCom, bMol, bDum, oenv);
912 flags = flags | TRX_READ_X | TRX_READ_V;
913 outekr = xvgropen(opt2fn("-ekr", NFILE, fnm),
914 "Center of mass rotation",
916 "Energy (kJ mol\\S-1\\N)",
918 make_legend(outekr, ngroups, isize[0], index[0], grpname, bCom, bMol, bDum, oenv);
922 flags = flags | TRX_READ_V;
926 flags = flags | TRX_READ_X | TRX_READ_V;
930 flags = flags | TRX_READ_X | TRX_READ_F;
932 if ((flags == 0) && !bOB)
934 fprintf(stderr, "Please select one or more output file options\n");
938 read_first_frame(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &fr, flags);
941 if ((bOV || bOF) && fn2ftp(ftp2fn(efTRX, NFILE, fnm)) == efXTC)
944 "Cannot extract velocities or forces since your input XTC file does not contain "
950 snew(sumx, fr.natoms);
954 snew(sumv, fr.natoms);
958 snew(sumf, fr.natoms);
966 gpbc = gmx_rmpbc_init(&top.idef, pbcType, fr.natoms);
971 time = output_env_conv_time(oenv, fr.time);
973 if (fr.bX && bNoJump && fr.bBox)
977 remove_jump(fr.box, fr.natoms, xp, fr.x);
983 for (i = 0; i < fr.natoms; i++)
985 copy_rvec(fr.x[i], xp[i]);
989 if (fr.bX && bCom && bPBC)
991 gmx_rmpbc_trxfr(gpbc, &fr);
996 update_histo(isize[0], index[0], fr.v, &nvhisto, &vhisto, binwidth);
1001 print_data(outx, time, fr.x, mass, bCom, ngroups, isize, index, bDim, sffmt);
1005 t_trxframe frout = fr;
1008 frout.atoms = &top.atoms;
1009 frout.bAtoms = TRUE;
1013 write_trx_x(status_out, &frout, mass, bCom, ngroups, isize, index);
1017 print_data(outv, time, fr.v, mass, bCom, ngroups, isize, index, bDim, sffmt);
1021 print_data(outf, time, fr.f, nullptr, bCom, ngroups, isize, index, bDim, sffmt);
1025 fprintf(outb, "\t%g", fr.time);
1034 fprintf(outb, "\n");
1038 fprintf(outt, " %g", time);
1039 for (i = 0; i < ngroups; i++)
1041 fprintf(outt, sffmt, temp(fr.v, mass, isize[i], index[i]));
1043 fprintf(outt, "\n");
1047 fprintf(outekt, " %g", time);
1048 for (i = 0; i < ngroups; i++)
1050 fprintf(outekt, sffmt, ektrans(fr.v, mass, isize[i], index[i]));
1052 fprintf(outekt, "\n");
1054 if (bEKR && fr.bX && fr.bV)
1056 fprintf(outekr, " %g", time);
1057 for (i = 0; i < ngroups; i++)
1059 fprintf(outekr, sffmt, ekrot(fr.x, fr.v, mass, isize[i], index[i]));
1061 fprintf(outekr, "\n");
1063 if ((bCV || bCF) && fr.bX
1064 && (ctime < 0 || (fr.time >= ctime * 0.999999 && fr.time <= ctime * 1.000001)))
1066 for (i = 0; i < fr.natoms; i++)
1068 rvec_inc(sumx[i], fr.x[i]);
1074 for (i = 0; i < fr.natoms; i++)
1076 rvec_inc(sumv[i], fr.v[i]);
1082 for (i = 0; i < fr.natoms; i++)
1084 rvec_inc(sumf[i], fr.f[i]);
1089 } while (read_next_frame(oenv, status, &fr));
1091 if (gpbc != nullptr)
1093 gmx_rmpbc_done(gpbc);
1096 /* clean up a bit */
1105 close_trx(status_out);
1134 print_histo(opt2fn("-vd", NFILE, fnm), nvhisto, vhisto, binwidth, oenv);
1141 if (pbcType != PbcType::No && !bNoJump)
1144 "\nWARNING: More than one frame was used for option -cv or -cf\n"
1145 "If atoms jump across the box you should use the -nojump or -ctime "
1148 for (i = 0; i < isize[0]; i++)
1150 svmul(1.0 / nr_xfr, sumx[index[0][i]], sumx[index[0][i]]);
1153 else if (nr_xfr == 0)
1155 fprintf(stderr, "\nWARNING: No coordinate frames found for option -cv or -cf\n\n");
1160 write_pdb_bfac(opt2fn("-cv", NFILE, fnm),
1161 opt2fn("-av", NFILE, fnm),
1178 write_pdb_bfac(opt2fn("-cf", NFILE, fnm),
1179 opt2fn("-af", NFILE, fnm),
1196 view_all(oenv, NFILE, fnm);
1199 // Free index and isize only if they are distinct from index0 and isize0
1202 for (int i = 0; i < ngroups; i++)
1209 for (int i = 0; i < ngroups; i++)
1218 output_env_done(oenv);