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, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
43 #include "gromacs/commandline/pargs.h"
50 #include "gromacs/fileio/futil.h"
54 #include "gromacs/fileio/tpxio.h"
55 #include "gromacs/fileio/trxio.h"
59 #include "gromacs/fileio/confio.h"
63 static void low_print_data(FILE *fp, real time, rvec x[], int n, atom_id *index,
64 gmx_bool bDim[], const char *sffmt)
68 fprintf(fp, " %g", time);
69 for (i = 0; i < n; i++)
79 for (d = 0; d < DIM; d++)
83 fprintf(fp, sffmt, x[ii][d]);
88 fprintf(fp, sffmt, norm(x[ii]));
94 static void average_data(rvec x[], rvec xav[], real *mass,
95 int ngrps, int isize[], atom_id **index)
100 double sum[DIM], mtot;
102 for (g = 0; g < ngrps; g++)
107 for (i = 0; i < isize[g]; i++)
113 svmul(m, x[ind], tmp);
114 for (d = 0; d < DIM; d++)
122 for (d = 0; d < DIM; d++)
130 for (d = 0; d < DIM; d++)
132 xav[g][d] = sum[d]/mtot;
137 /* mass=NULL, so these are forces: sum only (do not average) */
138 for (d = 0; d < DIM; d++)
146 static void print_data(FILE *fp, real time, rvec x[], real *mass, gmx_bool bCom,
147 int ngrps, int isize[], atom_id **index, gmx_bool bDim[],
150 static rvec *xav = NULL;
158 average_data(x, xav, mass, ngrps, isize, index);
159 low_print_data(fp, time, xav, ngrps, NULL, bDim, sffmt);
163 low_print_data(fp, time, x, isize[0], index[0], bDim, sffmt);
167 static void write_trx_x(t_trxstatus *status, t_trxframe *fr, real *mass, gmx_bool bCom,
168 int ngrps, int isize[], atom_id **index)
170 static rvec *xav = NULL;
171 static t_atoms *atoms = NULL;
184 snew(atoms->atom, ngrps);
186 /* Note that atom and residue names will be the ones
187 * of the first atom in each group.
189 for (i = 0; i < ngrps; i++)
191 atoms->atom[i] = fr->atoms->atom[index[i][0]];
192 atoms->atomname[i] = fr->atoms->atomname[index[i][0]];
195 average_data(fr->x, xav, mass, ngrps, isize, index);
197 fr_av.natoms = ngrps;
200 write_trxframe(status, &fr_av, NULL);
204 write_trxframe_indexed(status, fr, isize[0], index[0], NULL);
208 static void make_legend(FILE *fp, int ngrps, int isize, atom_id index[],
209 char **name, gmx_bool bCom, gmx_bool bMol, gmx_bool bDim[],
210 const output_env_t oenv)
213 const char *dimtxt[] = { " X", " Y", " Z", "" };
227 for (i = 0; i < n; i++)
229 for (d = 0; d <= DIM; d++)
233 snew(leg[j], STRLEN);
236 sprintf(leg[j], "mol %d%s", index[i]+1, dimtxt[d]);
240 sprintf(leg[j], "%s%s", name[i], dimtxt[d]);
244 sprintf(leg[j], "atom %d%s", index[i]+1, dimtxt[d]);
250 xvgr_legend(fp, j, (const char**)leg, oenv);
252 for (i = 0; i < j; i++)
259 static real ekrot(rvec x[], rvec v[], real mass[], int isize, atom_id index[])
261 static real **TCM = NULL, **L;
262 double tm, m0, lxx, lxy, lxz, lyy, lyz, lzz, ekrot;
271 for (i = 0; i < DIM; i++)
276 for (i = 0; i < DIM; i++)
286 for (i = 0; i < isize; i++)
291 cprod(x[j], v[j], a0);
292 for (m = 0; (m < DIM); m++)
294 xcm[m] += m0*x[j][m]; /* c.o.m. position */
295 vcm[m] += m0*v[j][m]; /* c.o.m. velocity */
296 acm[m] += m0*a0[m]; /* rotational velocity around c.o.m. */
299 dcprod(xcm, vcm, b0);
300 for (m = 0; (m < DIM); m++)
307 lxx = lxy = lxz = lyy = lyz = lzz = 0.0;
308 for (i = 0; i < isize; i++)
312 for (m = 0; m < DIM; m++)
314 dx[m] = x[j][m] - xcm[m];
316 lxx += dx[XX]*dx[XX]*m0;
317 lxy += dx[XX]*dx[YY]*m0;
318 lxz += dx[XX]*dx[ZZ]*m0;
319 lyy += dx[YY]*dx[YY]*m0;
320 lyz += dx[YY]*dx[ZZ]*m0;
321 lzz += dx[ZZ]*dx[ZZ]*m0;
324 L[XX][XX] = lyy + lzz;
328 L[YY][YY] = lxx + lzz;
332 L[ZZ][ZZ] = lxx + lyy;
334 m_inv_gen(L, DIM, TCM);
336 /* Compute omega (hoeksnelheid) */
339 for (m = 0; m < DIM; m++)
341 for (n = 0; n < DIM; n++)
343 ocm[m] += TCM[m][n]*acm[n];
345 ekrot += 0.5*ocm[m]*acm[m];
351 static real ektrans(rvec v[], real mass[], int isize, atom_id index[])
359 for (i = 0; i < isize; i++)
362 for (d = 0; d < DIM; d++)
364 mvcom[d] += mass[j]*v[j][d];
369 return dnorm2(mvcom)/(mtot*2);
372 static real temp(rvec v[], real mass[], int isize, atom_id index[])
378 for (i = 0; i < isize; i++)
381 ekin2 += mass[j]*norm2(v[j]);
384 return ekin2/(3*isize*BOLTZ);
387 static void remove_jump(matrix box, int natoms, rvec xp[], rvec x[])
392 for (d = 0; d < DIM; d++)
394 hbox[d] = 0.5*box[d][d];
396 for (i = 0; i < natoms; i++)
398 for (m = DIM-1; m >= 0; m--)
400 while (x[i][m]-xp[i][m] <= -hbox[m])
402 for (d = 0; d <= m; d++)
404 x[i][d] += box[m][d];
407 while (x[i][m]-xp[i][m] > hbox[m])
409 for (d = 0; d <= m; d++)
411 x[i][d] -= box[m][d];
418 static void write_pdb_bfac(const char *fname, const char *xname,
419 const char *title, t_atoms *atoms, int ePBC, matrix box,
420 int isize, atom_id *index, int nfr_x, rvec *x,
421 int nfr_v, rvec *sum,
422 gmx_bool bDim[], real scale_factor,
423 const output_env_t oenv)
426 real max, len2, scale;
431 if ((nfr_x == 0) || (nfr_v == 0))
433 fprintf(stderr, "No frames found for %s, will not write %s\n",
438 fprintf(stderr, "Used %d frames for %s\n", nfr_x, "coordinates");
439 fprintf(stderr, "Used %d frames for %s\n", nfr_v, title);
444 for (i = 0; i < DIM; i++)
458 for (i = 0; i < isize; i++)
460 svmul(scale, sum[index[i]], sum[index[i]]);
463 fp = xvgropen(xname, title, "Atom", "", oenv);
464 for (i = 0; i < isize; i++)
466 fprintf(fp, "%-5d %10.3f %10.3f %10.3f\n", 1+i,
467 sum[index[i]][XX], sum[index[i]][YY], sum[index[i]][ZZ]);
472 for (i = 0; i < isize; i++)
475 for (m = 0; m < DIM; m++)
477 if (bDim[m] || bDim[DIM])
479 len2 += sqr(sum[index[i]][m]);
488 if (scale_factor != 0)
490 scale = scale_factor;
500 scale = 10.0/sqrt(max);
504 printf("Maximum %s is %g on atom %d %s, res. %s %d\n",
505 title, sqrt(max), maxi+1, *(atoms->atomname[maxi]),
506 *(atoms->resinfo[atoms->atom[maxi].resind].name),
507 atoms->resinfo[atoms->atom[maxi].resind].nr);
509 if (atoms->pdbinfo == NULL)
511 snew(atoms->pdbinfo, atoms->nr);
515 for (i = 0; i < isize; i++)
518 for (m = 0; m < DIM; m++)
520 if (bDim[m] || bDim[DIM])
522 len2 += sqr(sum[index[i]][m]);
525 atoms->pdbinfo[index[i]].bfac = sqrt(len2)*scale;
530 for (i = 0; i < isize; i++)
532 atoms->pdbinfo[index[i]].bfac = sum[index[i]][onedim]*scale;
535 write_sto_conf_indexed(fname, title, atoms, x, NULL, ePBC, box, isize, index);
539 static void update_histo(int gnx, atom_id index[], rvec v[],
540 int *nhisto, int **histo, real binwidth)
548 for (i = 0; (i < gnx); i++)
550 vn = norm(v[index[i]]);
551 vnmax = max(vn, vnmax);
554 *nhisto = 1+(vnmax/binwidth);
555 snew(*histo, *nhisto);
557 for (i = 0; (i < gnx); i++)
559 vn = norm(v[index[i]]);
564 fprintf(stderr, "Extending histogram from %d to %d\n", *nhisto, nnn);
567 for (m = *nhisto; (m < nnn); m++)
577 static void print_histo(const char *fn, int nhisto, int histo[], real binwidth,
578 const output_env_t oenv)
583 fp = xvgropen(fn, "Velocity distribution", "V (nm/ps)", "arbitrary units",
585 for (i = 0; (i < nhisto); i++)
587 fprintf(fp, "%10.3e %10d\n", i*binwidth, histo[i]);
592 int gmx_traj(int argc, char *argv[])
594 const char *desc[] = {
595 "[THISMODULE] plots coordinates, velocities, forces and/or the box.",
596 "With [TT]-com[tt] the coordinates, velocities and forces are",
597 "calculated for the center of mass of each group.",
598 "When [TT]-mol[tt] is set, the numbers in the index file are",
599 "interpreted as molecule numbers and the same procedure as with",
600 "[TT]-com[tt] is used for each molecule.[PAR]",
601 "Option [TT]-ot[tt] plots the temperature of each group,",
602 "provided velocities are present in the trajectory file.",
603 "No corrections are made for constrained degrees of freedom!",
604 "This implies [TT]-com[tt].[PAR]",
605 "Options [TT]-ekt[tt] and [TT]-ekr[tt] plot the translational and",
606 "rotational kinetic energy of each group,",
607 "provided velocities are present in the trajectory file.",
608 "This implies [TT]-com[tt].[PAR]",
609 "Options [TT]-cv[tt] and [TT]-cf[tt] write the average velocities",
610 "and average forces as temperature factors to a [TT].pdb[tt] file with",
611 "the average coordinates or the coordinates at [TT]-ctime[tt].",
612 "The temperature factors are scaled such that the maximum is 10.",
613 "The scaling can be changed with the option [TT]-scale[tt].",
614 "To get the velocities or forces of one",
615 "frame set both [TT]-b[tt] and [TT]-e[tt] to the time of",
616 "desired frame. When averaging over frames you might need to use",
617 "the [TT]-nojump[tt] option to obtain the correct average coordinates.",
618 "If you select either of these option the average force and velocity",
619 "for each atom are written to an [TT].xvg[tt] file as well",
620 "(specified with [TT]-av[tt] or [TT]-af[tt]).[PAR]",
621 "Option [TT]-vd[tt] computes a velocity distribution, i.e. the",
622 "norm of the vector is plotted. In addition in the same graph",
623 "the kinetic energy distribution is given."
625 static gmx_bool bMol = FALSE, bCom = FALSE, bPBC = TRUE, bNoJump = FALSE;
626 static gmx_bool bX = TRUE, bY = TRUE, bZ = TRUE, bNorm = FALSE, bFP = FALSE;
627 static int ngroups = 1;
628 static real ctime = -1, scale = 0, binwidth = 1;
630 { "-com", FALSE, etBOOL, {&bCom},
631 "Plot data for the com of each group" },
632 { "-pbc", FALSE, etBOOL, {&bPBC},
633 "Make molecules whole for COM" },
634 { "-mol", FALSE, etBOOL, {&bMol},
635 "Index contains molecule numbers iso atom numbers" },
636 { "-nojump", FALSE, etBOOL, {&bNoJump},
637 "Remove jumps of atoms across the box" },
638 { "-x", FALSE, etBOOL, {&bX},
639 "Plot X-component" },
640 { "-y", FALSE, etBOOL, {&bY},
641 "Plot Y-component" },
642 { "-z", FALSE, etBOOL, {&bZ},
643 "Plot Z-component" },
644 { "-ng", FALSE, etINT, {&ngroups},
645 "Number of groups to consider" },
646 { "-len", FALSE, etBOOL, {&bNorm},
647 "Plot vector length" },
648 { "-fp", FALSE, etBOOL, {&bFP},
649 "Full precision output" },
650 { "-bin", FALSE, etREAL, {&binwidth},
651 "Binwidth for velocity histogram (nm/ps)" },
652 { "-ctime", FALSE, etREAL, {&ctime},
653 "Use frame at this time for x in [TT]-cv[tt] and [TT]-cf[tt] instead of the average x" },
654 { "-scale", FALSE, etREAL, {&scale},
655 "Scale factor for [TT].pdb[tt] output, 0 is autoscale" }
657 FILE *outx = NULL, *outv = NULL, *outf = NULL, *outb = NULL, *outt = NULL;
658 FILE *outekt = NULL, *outekr = NULL;
664 t_trxframe fr, frout;
665 int flags, nvhisto = 0, *vhisto = NULL;
666 rvec *xtop, *xp = NULL;
667 rvec *sumx = NULL, *sumv = NULL, *sumf = NULL;
670 t_trxstatus *status_out = NULL;
671 gmx_rmpbc_t gpbc = NULL;
673 int nr_xfr, nr_vfr, nr_ffr;
676 atom_id **index0, **index;
679 gmx_bool bTop, bOX, bOXT, bOV, bOF, bOB, bOT, bEKT, bEKR, bCV, bCF;
680 gmx_bool bDim[4], bDum[4], bVD;
681 char *sffmt, sffmt6[1024];
682 const char *box_leg[6] = { "XX", "YY", "ZZ", "YX", "ZX", "ZY" };
686 { efTRX, "-f", NULL, ffREAD },
687 { efTPS, NULL, NULL, ffREAD },
688 { efNDX, NULL, NULL, ffOPTRD },
689 { efXVG, "-ox", "coord.xvg", ffOPTWR },
690 { efTRX, "-oxt", "coord.xtc", ffOPTWR },
691 { efXVG, "-ov", "veloc.xvg", ffOPTWR },
692 { efXVG, "-of", "force.xvg", ffOPTWR },
693 { efXVG, "-ob", "box.xvg", ffOPTWR },
694 { efXVG, "-ot", "temp.xvg", ffOPTWR },
695 { efXVG, "-ekt", "ektrans.xvg", ffOPTWR },
696 { efXVG, "-ekr", "ekrot.xvg", ffOPTWR },
697 { efXVG, "-vd", "veldist.xvg", ffOPTWR },
698 { efPDB, "-cv", "veloc.pdb", ffOPTWR },
699 { efPDB, "-cf", "force.pdb", ffOPTWR },
700 { efXVG, "-av", "all_veloc.xvg", ffOPTWR },
701 { efXVG, "-af", "all_force.xvg", ffOPTWR }
703 #define NFILE asize(fnm)
705 if (!parse_common_args(&argc, argv,
706 PCA_CAN_TIME | PCA_TIME_UNIT | PCA_CAN_VIEW | PCA_BE_NICE,
707 NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL, &oenv))
714 fprintf(stderr, "Interpreting indexfile entries as molecules.\n"
715 "Using center of mass.\n");
718 bOX = opt2bSet("-ox", NFILE, fnm);
719 bOXT = opt2bSet("-oxt", NFILE, fnm);
720 bOV = opt2bSet("-ov", NFILE, fnm);
721 bOF = opt2bSet("-of", NFILE, fnm);
722 bOB = opt2bSet("-ob", NFILE, fnm);
723 bOT = opt2bSet("-ot", NFILE, fnm);
724 bEKT = opt2bSet("-ekt", NFILE, fnm);
725 bEKR = opt2bSet("-ekr", NFILE, fnm);
726 bCV = opt2bSet("-cv", NFILE, fnm) || opt2bSet("-av", NFILE, fnm);
727 bCF = opt2bSet("-cf", NFILE, fnm) || opt2bSet("-af", NFILE, fnm);
728 bVD = opt2bSet("-vd", NFILE, fnm) || opt2parg_bSet("-bin", asize(pa), pa);
729 if (bMol || bOT || bEKT || bEKR)
741 sffmt = "\t" gmx_real_fullprecision_pfmt;
747 sprintf(sffmt6, "%s%s%s%s%s%s", sffmt, sffmt, sffmt, sffmt, sffmt, sffmt);
749 bTop = read_tps_conf(ftp2fn(efTPS, NFILE, fnm), title, &top, &ePBC,
751 bCom && (bOX || bOXT || bOV || bOT || bEKT || bEKR));
753 if ((bMol || bCV || bCF) && !bTop)
755 gmx_fatal(FARGS, "Need a run input file for option -mol, -cv or -cf");
760 indexfn = ftp2fn(efNDX, NFILE, fnm);
764 indexfn = ftp2fn_null(efNDX, NFILE, fnm);
767 if (!(bCom && !bMol))
771 snew(grpname, ngroups);
772 snew(isize0, ngroups);
773 snew(index0, ngroups);
774 get_index(&(top.atoms), indexfn, ngroups, isize0, index0, grpname);
781 snew(isize, ngroups);
782 snew(index, ngroups);
783 for (i = 0; i < ngroups; i++)
785 if (index0[0][i] < 0 || index0[0][i] >= mols->nr)
787 gmx_fatal(FARGS, "Molecule index (%d) is out of range (%d-%d)",
788 index0[0][i]+1, 1, mols->nr);
790 isize[i] = atndx[index0[0][i]+1] - atndx[index0[0][i]];
791 snew(index[i], isize[i]);
792 for (j = 0; j < isize[i]; j++)
794 index[i][j] = atndx[index0[0][i]] + j;
805 snew(mass, top.atoms.nr);
806 for (i = 0; i < top.atoms.nr; i++)
808 mass[i] = top.atoms.atom[i].m;
819 flags = flags | TRX_READ_X;
820 outx = xvgropen(opt2fn("-ox", NFILE, fnm),
821 bCom ? "Center of mass" : "Coordinate",
822 output_env_get_xvgr_tlabel(oenv), "Coordinate (nm)", oenv);
823 make_legend(outx, ngroups, isize0[0], index0[0], grpname, bCom, bMol, bDim, oenv);
827 flags = flags | TRX_READ_X;
828 status_out = open_trx(opt2fn("-oxt", NFILE, fnm), "w");
832 flags = flags | TRX_READ_V;
833 outv = xvgropen(opt2fn("-ov", NFILE, fnm),
834 bCom ? "Center of mass velocity" : "Velocity",
835 output_env_get_xvgr_tlabel(oenv), "Velocity (nm/ps)", oenv);
836 make_legend(outv, ngroups, isize0[0], index0[0], grpname, bCom, bMol, bDim, oenv);
840 flags = flags | TRX_READ_F;
841 outf = xvgropen(opt2fn("-of", NFILE, fnm), "Force",
842 output_env_get_xvgr_tlabel(oenv), "Force (kJ mol\\S-1\\N nm\\S-1\\N)",
844 make_legend(outf, ngroups, isize0[0], index0[0], grpname, bCom, bMol, bDim, oenv);
848 outb = xvgropen(opt2fn("-ob", NFILE, fnm), "Box vector elements",
849 output_env_get_xvgr_tlabel(oenv), "(nm)", oenv);
851 xvgr_legend(outb, 6, box_leg, oenv);
859 flags = flags | TRX_READ_V;
860 outt = xvgropen(opt2fn("-ot", NFILE, fnm), "Temperature",
861 output_env_get_xvgr_tlabel(oenv), "(K)", oenv);
862 make_legend(outt, ngroups, isize[0], index[0], grpname, bCom, bMol, bDum, oenv);
870 flags = flags | TRX_READ_V;
871 outekt = xvgropen(opt2fn("-ekt", NFILE, fnm), "Center of mass translation",
872 output_env_get_xvgr_tlabel(oenv), "Energy (kJ mol\\S-1\\N)", oenv);
873 make_legend(outekt, ngroups, isize[0], index[0], grpname, bCom, bMol, bDum, oenv);
881 flags = flags | TRX_READ_X | TRX_READ_V;
882 outekr = xvgropen(opt2fn("-ekr", NFILE, fnm), "Center of mass rotation",
883 output_env_get_xvgr_tlabel(oenv), "Energy (kJ mol\\S-1\\N)", oenv);
884 make_legend(outekr, ngroups, isize[0], index[0], grpname, bCom, bMol, bDum, oenv);
888 flags = flags | TRX_READ_V;
892 flags = flags | TRX_READ_X | TRX_READ_V;
896 flags = flags | TRX_READ_X | TRX_READ_F;
898 if ((flags == 0) && !bOB)
900 fprintf(stderr, "Please select one or more output file options\n");
904 read_first_frame(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &fr, flags);
908 snew(sumx, fr.natoms);
912 snew(sumv, fr.natoms);
916 snew(sumf, fr.natoms);
924 gpbc = gmx_rmpbc_init(&top.idef, ePBC, fr.natoms);
929 time = output_env_conv_time(oenv, fr.time);
931 if (fr.bX && bNoJump && fr.bBox)
935 remove_jump(fr.box, fr.natoms, xp, fr.x);
941 for (i = 0; i < fr.natoms; i++)
943 copy_rvec(fr.x[i], xp[i]);
947 if (fr.bX && bCom && bPBC)
949 gmx_rmpbc_trxfr(gpbc, &fr);
954 update_histo(isize[0], index[0], fr.v, &nvhisto, &vhisto, binwidth);
959 print_data(outx, time, fr.x, mass, bCom, ngroups, isize, index, bDim, sffmt);
966 frout.atoms = &top.atoms;
969 write_trx_x(status_out, &frout, mass, bCom, ngroups, isize, index);
973 print_data(outv, time, fr.v, mass, bCom, ngroups, isize, index, bDim, sffmt);
977 print_data(outf, time, fr.f, NULL, bCom, ngroups, isize, index, bDim, sffmt);
981 fprintf(outb, "\t%g", fr.time);
982 fprintf(outb, sffmt6,
983 fr.box[XX][XX], fr.box[YY][YY], fr.box[ZZ][ZZ],
984 fr.box[YY][XX], fr.box[ZZ][XX], fr.box[ZZ][YY]);
989 fprintf(outt, " %g", time);
990 for (i = 0; i < ngroups; i++)
992 fprintf(outt, sffmt, temp(fr.v, mass, isize[i], index[i]));
998 fprintf(outekt, " %g", time);
999 for (i = 0; i < ngroups; i++)
1001 fprintf(outekt, sffmt, ektrans(fr.v, mass, isize[i], index[i]));
1003 fprintf(outekt, "\n");
1005 if (bEKR && fr.bX && fr.bV)
1007 fprintf(outekr, " %g", time);
1008 for (i = 0; i < ngroups; i++)
1010 fprintf(outekr, sffmt, ekrot(fr.x, fr.v, mass, isize[i], index[i]));
1012 fprintf(outekr, "\n");
1014 if ((bCV || bCF) && fr.bX &&
1015 (ctime < 0 || (fr.time >= ctime*0.999999 &&
1016 fr.time <= ctime*1.000001)))
1018 for (i = 0; i < fr.natoms; i++)
1020 rvec_inc(sumx[i], fr.x[i]);
1026 for (i = 0; i < fr.natoms; i++)
1028 rvec_inc(sumv[i], fr.v[i]);
1034 for (i = 0; i < fr.natoms; i++)
1036 rvec_inc(sumf[i], fr.f[i]);
1042 while (read_next_frame(oenv, status, &fr));
1046 gmx_rmpbc_done(gpbc);
1049 /* clean up a bit */
1058 close_trx(status_out);
1078 gmx_ffclose(outekt);
1082 gmx_ffclose(outekr);
1087 print_histo(opt2fn("-vd", NFILE, fnm), nvhisto, vhisto, binwidth, oenv);
1094 if (ePBC != epbcNONE && !bNoJump)
1096 fprintf(stderr, "\nWARNING: More than one frame was used for option -cv or -cf\n"
1097 "If atoms jump across the box you should use the -nojump or -ctime option\n\n");
1099 for (i = 0; i < isize[0]; i++)
1101 svmul(1.0/nr_xfr, sumx[index[0][i]], sumx[index[0][i]]);
1104 else if (nr_xfr == 0)
1106 fprintf(stderr, "\nWARNING: No coordinate frames found for option -cv or -cf\n\n");
1111 write_pdb_bfac(opt2fn("-cv", NFILE, fnm),
1112 opt2fn("-av", NFILE, fnm), "average velocity", &(top.atoms),
1113 ePBC, topbox, isize[0], index[0], nr_xfr, sumx,
1114 nr_vfr, sumv, bDim, scale, oenv);
1118 write_pdb_bfac(opt2fn("-cf", NFILE, fnm),
1119 opt2fn("-af", NFILE, fnm), "average force", &(top.atoms),
1120 ePBC, topbox, isize[0], index[0], nr_xfr, sumx,
1121 nr_ffr, sumf, bDim, scale, oenv);
1125 view_all(oenv, NFILE, fnm);