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.
45 #include "gromacs/commandline/pargs.h"
47 #include "gromacs/utility/smalloc.h"
49 #include "gromacs/math/vec.h"
51 #include "gromacs/utility/futil.h"
54 #include "gromacs/fileio/xvgr.h"
56 #include "gromacs/fileio/tpxio.h"
57 #include "gromacs/fileio/trxio.h"
60 #include "gromacs/fileio/confio.h"
63 #include "gromacs/linearalgebra/nrjac.h"
64 #include "gromacs/utility/fatalerror.h"
66 static void low_print_data(FILE *fp, real time, rvec x[], int n, atom_id *index,
67 gmx_bool bDim[], const char *sffmt)
71 fprintf(fp, " %g", time);
72 for (i = 0; i < n; i++)
82 for (d = 0; d < DIM; d++)
86 fprintf(fp, sffmt, x[ii][d]);
91 fprintf(fp, sffmt, norm(x[ii]));
97 static void average_data(rvec x[], rvec xav[], real *mass,
98 int ngrps, int isize[], atom_id **index)
103 double sum[DIM], mtot;
105 for (g = 0; g < ngrps; g++)
110 for (i = 0; i < isize[g]; i++)
116 svmul(m, x[ind], tmp);
117 for (d = 0; d < DIM; d++)
125 for (d = 0; d < DIM; d++)
133 for (d = 0; d < DIM; d++)
135 xav[g][d] = sum[d]/mtot;
140 /* mass=NULL, so these are forces: sum only (do not average) */
141 for (d = 0; d < DIM; d++)
149 static void print_data(FILE *fp, real time, rvec x[], real *mass, gmx_bool bCom,
150 int ngrps, int isize[], atom_id **index, gmx_bool bDim[],
153 static rvec *xav = NULL;
161 average_data(x, xav, mass, ngrps, isize, index);
162 low_print_data(fp, time, xav, ngrps, NULL, bDim, sffmt);
166 low_print_data(fp, time, x, isize[0], index[0], bDim, sffmt);
170 static void write_trx_x(t_trxstatus *status, t_trxframe *fr, real *mass, gmx_bool bCom,
171 int ngrps, int isize[], atom_id **index)
173 static rvec *xav = NULL;
174 static t_atoms *atoms = NULL;
187 snew(atoms->atom, ngrps);
189 /* Note that atom and residue names will be the ones
190 * of the first atom in each group.
192 for (i = 0; i < ngrps; i++)
194 atoms->atom[i] = fr->atoms->atom[index[i][0]];
195 atoms->atomname[i] = fr->atoms->atomname[index[i][0]];
198 average_data(fr->x, xav, mass, ngrps, isize, index);
200 fr_av.natoms = ngrps;
203 write_trxframe(status, &fr_av, NULL);
207 write_trxframe_indexed(status, fr, isize[0], index[0], NULL);
211 static void make_legend(FILE *fp, int ngrps, int isize, atom_id index[],
212 char **name, gmx_bool bCom, gmx_bool bMol, gmx_bool bDim[],
213 const output_env_t oenv)
216 const char *dimtxt[] = { " X", " Y", " Z", "" };
230 for (i = 0; i < n; i++)
232 for (d = 0; d <= DIM; d++)
236 snew(leg[j], STRLEN);
239 sprintf(leg[j], "mol %d%s", index[i]+1, dimtxt[d]);
243 sprintf(leg[j], "%s%s", name[i], dimtxt[d]);
247 sprintf(leg[j], "atom %d%s", index[i]+1, dimtxt[d]);
253 xvgr_legend(fp, j, (const char**)leg, oenv);
255 for (i = 0; i < j; i++)
262 static real ekrot(rvec x[], rvec v[], real mass[], int isize, atom_id index[])
264 static real **TCM = NULL, **L;
265 double tm, m0, lxx, lxy, lxz, lyy, lyz, lzz, ekrot;
274 for (i = 0; i < DIM; i++)
279 for (i = 0; i < DIM; i++)
289 for (i = 0; i < isize; i++)
294 cprod(x[j], v[j], a0);
295 for (m = 0; (m < DIM); m++)
297 xcm[m] += m0*x[j][m]; /* c.o.m. position */
298 vcm[m] += m0*v[j][m]; /* c.o.m. velocity */
299 acm[m] += m0*a0[m]; /* rotational velocity around c.o.m. */
302 dcprod(xcm, vcm, b0);
303 for (m = 0; (m < DIM); m++)
310 lxx = lxy = lxz = lyy = lyz = lzz = 0.0;
311 for (i = 0; i < isize; i++)
315 for (m = 0; m < DIM; m++)
317 dx[m] = x[j][m] - xcm[m];
319 lxx += dx[XX]*dx[XX]*m0;
320 lxy += dx[XX]*dx[YY]*m0;
321 lxz += dx[XX]*dx[ZZ]*m0;
322 lyy += dx[YY]*dx[YY]*m0;
323 lyz += dx[YY]*dx[ZZ]*m0;
324 lzz += dx[ZZ]*dx[ZZ]*m0;
327 L[XX][XX] = lyy + lzz;
331 L[YY][YY] = lxx + lzz;
335 L[ZZ][ZZ] = lxx + lyy;
337 m_inv_gen(L, DIM, TCM);
339 /* Compute omega (hoeksnelheid) */
342 for (m = 0; m < DIM; m++)
344 for (n = 0; n < DIM; n++)
346 ocm[m] += TCM[m][n]*acm[n];
348 ekrot += 0.5*ocm[m]*acm[m];
354 static real ektrans(rvec v[], real mass[], int isize, atom_id index[])
362 for (i = 0; i < isize; i++)
365 for (d = 0; d < DIM; d++)
367 mvcom[d] += mass[j]*v[j][d];
372 return dnorm2(mvcom)/(mtot*2);
375 static real temp(rvec v[], real mass[], int isize, atom_id index[])
381 for (i = 0; i < isize; i++)
384 ekin2 += mass[j]*norm2(v[j]);
387 return ekin2/(3*isize*BOLTZ);
390 static void remove_jump(matrix box, int natoms, rvec xp[], rvec x[])
395 for (d = 0; d < DIM; d++)
397 hbox[d] = 0.5*box[d][d];
399 for (i = 0; i < natoms; i++)
401 for (m = DIM-1; m >= 0; m--)
403 while (x[i][m]-xp[i][m] <= -hbox[m])
405 for (d = 0; d <= m; d++)
407 x[i][d] += box[m][d];
410 while (x[i][m]-xp[i][m] > hbox[m])
412 for (d = 0; d <= m; d++)
414 x[i][d] -= box[m][d];
421 static void write_pdb_bfac(const char *fname, const char *xname,
422 const char *title, t_atoms *atoms, int ePBC, matrix box,
423 int isize, atom_id *index, int nfr_x, rvec *x,
424 int nfr_v, rvec *sum,
425 gmx_bool bDim[], real scale_factor,
426 const output_env_t oenv)
429 real max, len2, scale;
434 if ((nfr_x == 0) || (nfr_v == 0))
436 fprintf(stderr, "No frames found for %s, will not write %s\n",
441 fprintf(stderr, "Used %d frames for %s\n", nfr_x, "coordinates");
442 fprintf(stderr, "Used %d frames for %s\n", nfr_v, title);
447 for (i = 0; i < DIM; i++)
461 for (i = 0; i < isize; i++)
463 svmul(scale, sum[index[i]], sum[index[i]]);
466 fp = xvgropen(xname, title, "Atom", "", oenv);
467 for (i = 0; i < isize; i++)
469 fprintf(fp, "%-5d %10.3f %10.3f %10.3f\n", 1+i,
470 sum[index[i]][XX], sum[index[i]][YY], sum[index[i]][ZZ]);
475 for (i = 0; i < isize; i++)
478 for (m = 0; m < DIM; m++)
480 if (bDim[m] || bDim[DIM])
482 len2 += sqr(sum[index[i]][m]);
491 if (scale_factor != 0)
493 scale = scale_factor;
503 scale = 10.0/sqrt(max);
507 printf("Maximum %s is %g on atom %d %s, res. %s %d\n",
508 title, sqrt(max), maxi+1, *(atoms->atomname[maxi]),
509 *(atoms->resinfo[atoms->atom[maxi].resind].name),
510 atoms->resinfo[atoms->atom[maxi].resind].nr);
512 if (atoms->pdbinfo == NULL)
514 snew(atoms->pdbinfo, atoms->nr);
518 for (i = 0; i < isize; i++)
521 for (m = 0; m < DIM; m++)
523 if (bDim[m] || bDim[DIM])
525 len2 += sqr(sum[index[i]][m]);
528 atoms->pdbinfo[index[i]].bfac = sqrt(len2)*scale;
533 for (i = 0; i < isize; i++)
535 atoms->pdbinfo[index[i]].bfac = sum[index[i]][onedim]*scale;
538 write_sto_conf_indexed(fname, title, atoms, x, NULL, ePBC, box, isize, index);
542 static void update_histo(int gnx, atom_id index[], rvec v[],
543 int *nhisto, int **histo, real binwidth)
551 for (i = 0; (i < gnx); i++)
553 vn = norm(v[index[i]]);
554 vnmax = max(vn, vnmax);
557 *nhisto = 1+(vnmax/binwidth);
558 snew(*histo, *nhisto);
560 for (i = 0; (i < gnx); i++)
562 vn = norm(v[index[i]]);
567 fprintf(stderr, "Extending histogram from %d to %d\n", *nhisto, nnn);
570 for (m = *nhisto; (m < nnn); m++)
580 static void print_histo(const char *fn, int nhisto, int histo[], real binwidth,
581 const output_env_t oenv)
586 fp = xvgropen(fn, "Velocity distribution", "V (nm/ps)", "arbitrary units",
588 for (i = 0; (i < nhisto); i++)
590 fprintf(fp, "%10.3e %10d\n", i*binwidth, histo[i]);
595 int gmx_traj(int argc, char *argv[])
597 const char *desc[] = {
598 "[THISMODULE] plots coordinates, velocities, forces and/or the box.",
599 "With [TT]-com[tt] the coordinates, velocities and forces are",
600 "calculated for the center of mass of each group.",
601 "When [TT]-mol[tt] is set, the numbers in the index file are",
602 "interpreted as molecule numbers and the same procedure as with",
603 "[TT]-com[tt] is used for each molecule.[PAR]",
604 "Option [TT]-ot[tt] plots the temperature of each group,",
605 "provided velocities are present in the trajectory file.",
606 "No corrections are made for constrained degrees of freedom!",
607 "This implies [TT]-com[tt].[PAR]",
608 "Options [TT]-ekt[tt] and [TT]-ekr[tt] plot the translational and",
609 "rotational kinetic energy of each group,",
610 "provided velocities are present in the trajectory file.",
611 "This implies [TT]-com[tt].[PAR]",
612 "Options [TT]-cv[tt] and [TT]-cf[tt] write the average velocities",
613 "and average forces as temperature factors to a [TT].pdb[tt] file with",
614 "the average coordinates or the coordinates at [TT]-ctime[tt].",
615 "The temperature factors are scaled such that the maximum is 10.",
616 "The scaling can be changed with the option [TT]-scale[tt].",
617 "To get the velocities or forces of one",
618 "frame set both [TT]-b[tt] and [TT]-e[tt] to the time of",
619 "desired frame. When averaging over frames you might need to use",
620 "the [TT]-nojump[tt] option to obtain the correct average coordinates.",
621 "If you select either of these option the average force and velocity",
622 "for each atom are written to an [TT].xvg[tt] file as well",
623 "(specified with [TT]-av[tt] or [TT]-af[tt]).[PAR]",
624 "Option [TT]-vd[tt] computes a velocity distribution, i.e. the",
625 "norm of the vector is plotted. In addition in the same graph",
626 "the kinetic energy distribution is given."
628 static gmx_bool bMol = FALSE, bCom = FALSE, bPBC = TRUE, bNoJump = FALSE;
629 static gmx_bool bX = TRUE, bY = TRUE, bZ = TRUE, bNorm = FALSE, bFP = FALSE;
630 static int ngroups = 1;
631 static real ctime = -1, scale = 0, binwidth = 1;
633 { "-com", FALSE, etBOOL, {&bCom},
634 "Plot data for the com of each group" },
635 { "-pbc", FALSE, etBOOL, {&bPBC},
636 "Make molecules whole for COM" },
637 { "-mol", FALSE, etBOOL, {&bMol},
638 "Index contains molecule numbers iso atom numbers" },
639 { "-nojump", FALSE, etBOOL, {&bNoJump},
640 "Remove jumps of atoms across the box" },
641 { "-x", FALSE, etBOOL, {&bX},
642 "Plot X-component" },
643 { "-y", FALSE, etBOOL, {&bY},
644 "Plot Y-component" },
645 { "-z", FALSE, etBOOL, {&bZ},
646 "Plot Z-component" },
647 { "-ng", FALSE, etINT, {&ngroups},
648 "Number of groups to consider" },
649 { "-len", FALSE, etBOOL, {&bNorm},
650 "Plot vector length" },
651 { "-fp", FALSE, etBOOL, {&bFP},
652 "Full precision output" },
653 { "-bin", FALSE, etREAL, {&binwidth},
654 "Binwidth for velocity histogram (nm/ps)" },
655 { "-ctime", FALSE, etREAL, {&ctime},
656 "Use frame at this time for x in [TT]-cv[tt] and [TT]-cf[tt] instead of the average x" },
657 { "-scale", FALSE, etREAL, {&scale},
658 "Scale factor for [TT].pdb[tt] output, 0 is autoscale" }
660 FILE *outx = NULL, *outv = NULL, *outf = NULL, *outb = NULL, *outt = NULL;
661 FILE *outekt = NULL, *outekr = NULL;
667 t_trxframe fr, frout;
668 int flags, nvhisto = 0, *vhisto = NULL;
669 rvec *xtop, *xp = NULL;
670 rvec *sumx = NULL, *sumv = NULL, *sumf = NULL;
673 t_trxstatus *status_out = NULL;
674 gmx_rmpbc_t gpbc = NULL;
676 int nr_xfr, nr_vfr, nr_ffr;
679 atom_id **index0, **index;
682 gmx_bool bTop, bOX, bOXT, bOV, bOF, bOB, bOT, bEKT, bEKR, bCV, bCF;
683 gmx_bool bDim[4], bDum[4], bVD;
684 char *sffmt, sffmt6[1024];
685 const char *box_leg[6] = { "XX", "YY", "ZZ", "YX", "ZX", "ZY" };
689 { efTRX, "-f", NULL, ffREAD },
690 { efTPS, NULL, NULL, ffREAD },
691 { efNDX, NULL, NULL, ffOPTRD },
692 { efXVG, "-ox", "coord.xvg", ffOPTWR },
693 { efTRX, "-oxt", "coord.xtc", ffOPTWR },
694 { efXVG, "-ov", "veloc.xvg", ffOPTWR },
695 { efXVG, "-of", "force.xvg", ffOPTWR },
696 { efXVG, "-ob", "box.xvg", ffOPTWR },
697 { efXVG, "-ot", "temp.xvg", ffOPTWR },
698 { efXVG, "-ekt", "ektrans.xvg", ffOPTWR },
699 { efXVG, "-ekr", "ekrot.xvg", ffOPTWR },
700 { efXVG, "-vd", "veldist.xvg", ffOPTWR },
701 { efPDB, "-cv", "veloc.pdb", ffOPTWR },
702 { efPDB, "-cf", "force.pdb", ffOPTWR },
703 { efXVG, "-av", "all_veloc.xvg", ffOPTWR },
704 { efXVG, "-af", "all_force.xvg", ffOPTWR }
706 #define NFILE asize(fnm)
708 if (!parse_common_args(&argc, argv,
709 PCA_CAN_TIME | PCA_TIME_UNIT | PCA_CAN_VIEW | PCA_BE_NICE,
710 NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL, &oenv))
717 fprintf(stderr, "Interpreting indexfile entries as molecules.\n"
718 "Using center of mass.\n");
721 bOX = opt2bSet("-ox", NFILE, fnm);
722 bOXT = opt2bSet("-oxt", NFILE, fnm);
723 bOV = opt2bSet("-ov", NFILE, fnm);
724 bOF = opt2bSet("-of", NFILE, fnm);
725 bOB = opt2bSet("-ob", NFILE, fnm);
726 bOT = opt2bSet("-ot", NFILE, fnm);
727 bEKT = opt2bSet("-ekt", NFILE, fnm);
728 bEKR = opt2bSet("-ekr", NFILE, fnm);
729 bCV = opt2bSet("-cv", NFILE, fnm) || opt2bSet("-av", NFILE, fnm);
730 bCF = opt2bSet("-cf", NFILE, fnm) || opt2bSet("-af", NFILE, fnm);
731 bVD = opt2bSet("-vd", NFILE, fnm) || opt2parg_bSet("-bin", asize(pa), pa);
732 if (bMol || bOT || bEKT || bEKR)
744 sffmt = "\t" gmx_real_fullprecision_pfmt;
750 sprintf(sffmt6, "%s%s%s%s%s%s", sffmt, sffmt, sffmt, sffmt, sffmt, sffmt);
752 bTop = read_tps_conf(ftp2fn(efTPS, NFILE, fnm), title, &top, &ePBC,
754 bCom && (bOX || bOXT || bOV || bOT || bEKT || bEKR));
756 if ((bMol || bCV || bCF) && !bTop)
758 gmx_fatal(FARGS, "Need a run input file for option -mol, -cv or -cf");
763 indexfn = ftp2fn(efNDX, NFILE, fnm);
767 indexfn = ftp2fn_null(efNDX, NFILE, fnm);
770 if (!(bCom && !bMol))
774 snew(grpname, ngroups);
775 snew(isize0, ngroups);
776 snew(index0, ngroups);
777 get_index(&(top.atoms), indexfn, ngroups, isize0, index0, grpname);
784 snew(isize, ngroups);
785 snew(index, ngroups);
786 for (i = 0; i < ngroups; i++)
788 if (index0[0][i] < 0 || index0[0][i] >= mols->nr)
790 gmx_fatal(FARGS, "Molecule index (%d) is out of range (%d-%d)",
791 index0[0][i]+1, 1, mols->nr);
793 isize[i] = atndx[index0[0][i]+1] - atndx[index0[0][i]];
794 snew(index[i], isize[i]);
795 for (j = 0; j < isize[i]; j++)
797 index[i][j] = atndx[index0[0][i]] + j;
808 snew(mass, top.atoms.nr);
809 for (i = 0; i < top.atoms.nr; i++)
811 mass[i] = top.atoms.atom[i].m;
822 flags = flags | TRX_READ_X;
823 outx = xvgropen(opt2fn("-ox", NFILE, fnm),
824 bCom ? "Center of mass" : "Coordinate",
825 output_env_get_xvgr_tlabel(oenv), "Coordinate (nm)", oenv);
826 make_legend(outx, ngroups, isize0[0], index0[0], grpname, bCom, bMol, bDim, oenv);
830 flags = flags | TRX_READ_X;
831 status_out = open_trx(opt2fn("-oxt", NFILE, fnm), "w");
835 flags = flags | TRX_READ_V;
836 outv = xvgropen(opt2fn("-ov", NFILE, fnm),
837 bCom ? "Center of mass velocity" : "Velocity",
838 output_env_get_xvgr_tlabel(oenv), "Velocity (nm/ps)", oenv);
839 make_legend(outv, ngroups, isize0[0], index0[0], grpname, bCom, bMol, bDim, oenv);
843 flags = flags | TRX_READ_F;
844 outf = xvgropen(opt2fn("-of", NFILE, fnm), "Force",
845 output_env_get_xvgr_tlabel(oenv), "Force (kJ mol\\S-1\\N nm\\S-1\\N)",
847 make_legend(outf, ngroups, isize0[0], index0[0], grpname, bCom, bMol, bDim, oenv);
851 outb = xvgropen(opt2fn("-ob", NFILE, fnm), "Box vector elements",
852 output_env_get_xvgr_tlabel(oenv), "(nm)", oenv);
854 xvgr_legend(outb, 6, box_leg, oenv);
862 flags = flags | TRX_READ_V;
863 outt = xvgropen(opt2fn("-ot", NFILE, fnm), "Temperature",
864 output_env_get_xvgr_tlabel(oenv), "(K)", oenv);
865 make_legend(outt, ngroups, isize[0], index[0], grpname, bCom, bMol, bDum, oenv);
873 flags = flags | TRX_READ_V;
874 outekt = xvgropen(opt2fn("-ekt", NFILE, fnm), "Center of mass translation",
875 output_env_get_xvgr_tlabel(oenv), "Energy (kJ mol\\S-1\\N)", oenv);
876 make_legend(outekt, ngroups, isize[0], index[0], grpname, bCom, bMol, bDum, oenv);
884 flags = flags | TRX_READ_X | TRX_READ_V;
885 outekr = xvgropen(opt2fn("-ekr", NFILE, fnm), "Center of mass rotation",
886 output_env_get_xvgr_tlabel(oenv), "Energy (kJ mol\\S-1\\N)", oenv);
887 make_legend(outekr, ngroups, isize[0], index[0], grpname, bCom, bMol, bDum, oenv);
891 flags = flags | TRX_READ_V;
895 flags = flags | TRX_READ_X | TRX_READ_V;
899 flags = flags | TRX_READ_X | TRX_READ_F;
901 if ((flags == 0) && !bOB)
903 fprintf(stderr, "Please select one or more output file options\n");
907 read_first_frame(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &fr, flags);
911 snew(sumx, fr.natoms);
915 snew(sumv, fr.natoms);
919 snew(sumf, fr.natoms);
927 gpbc = gmx_rmpbc_init(&top.idef, ePBC, fr.natoms);
932 time = output_env_conv_time(oenv, fr.time);
934 if (fr.bX && bNoJump && fr.bBox)
938 remove_jump(fr.box, fr.natoms, xp, fr.x);
944 for (i = 0; i < fr.natoms; i++)
946 copy_rvec(fr.x[i], xp[i]);
950 if (fr.bX && bCom && bPBC)
952 gmx_rmpbc_trxfr(gpbc, &fr);
957 update_histo(isize[0], index[0], fr.v, &nvhisto, &vhisto, binwidth);
962 print_data(outx, time, fr.x, mass, bCom, ngroups, isize, index, bDim, sffmt);
969 frout.atoms = &top.atoms;
972 write_trx_x(status_out, &frout, mass, bCom, ngroups, isize, index);
976 print_data(outv, time, fr.v, mass, bCom, ngroups, isize, index, bDim, sffmt);
980 print_data(outf, time, fr.f, NULL, bCom, ngroups, isize, index, bDim, sffmt);
984 fprintf(outb, "\t%g", fr.time);
985 fprintf(outb, sffmt6,
986 fr.box[XX][XX], fr.box[YY][YY], fr.box[ZZ][ZZ],
987 fr.box[YY][XX], fr.box[ZZ][XX], fr.box[ZZ][YY]);
992 fprintf(outt, " %g", time);
993 for (i = 0; i < ngroups; i++)
995 fprintf(outt, sffmt, temp(fr.v, mass, isize[i], index[i]));
1001 fprintf(outekt, " %g", time);
1002 for (i = 0; i < ngroups; i++)
1004 fprintf(outekt, sffmt, ektrans(fr.v, mass, isize[i], index[i]));
1006 fprintf(outekt, "\n");
1008 if (bEKR && fr.bX && fr.bV)
1010 fprintf(outekr, " %g", time);
1011 for (i = 0; i < ngroups; i++)
1013 fprintf(outekr, sffmt, ekrot(fr.x, fr.v, mass, isize[i], index[i]));
1015 fprintf(outekr, "\n");
1017 if ((bCV || bCF) && fr.bX &&
1018 (ctime < 0 || (fr.time >= ctime*0.999999 &&
1019 fr.time <= ctime*1.000001)))
1021 for (i = 0; i < fr.natoms; i++)
1023 rvec_inc(sumx[i], fr.x[i]);
1029 for (i = 0; i < fr.natoms; i++)
1031 rvec_inc(sumv[i], fr.v[i]);
1037 for (i = 0; i < fr.natoms; i++)
1039 rvec_inc(sumf[i], fr.f[i]);
1045 while (read_next_frame(oenv, status, &fr));
1049 gmx_rmpbc_done(gpbc);
1052 /* clean up a bit */
1061 close_trx(status_out);
1081 gmx_ffclose(outekt);
1085 gmx_ffclose(outekr);
1090 print_histo(opt2fn("-vd", NFILE, fnm), nvhisto, vhisto, binwidth, oenv);
1097 if (ePBC != epbcNONE && !bNoJump)
1099 fprintf(stderr, "\nWARNING: More than one frame was used for option -cv or -cf\n"
1100 "If atoms jump across the box you should use the -nojump or -ctime option\n\n");
1102 for (i = 0; i < isize[0]; i++)
1104 svmul(1.0/nr_xfr, sumx[index[0][i]], sumx[index[0][i]]);
1107 else if (nr_xfr == 0)
1109 fprintf(stderr, "\nWARNING: No coordinate frames found for option -cv or -cf\n\n");
1114 write_pdb_bfac(opt2fn("-cv", NFILE, fnm),
1115 opt2fn("-av", NFILE, fnm), "average velocity", &(top.atoms),
1116 ePBC, topbox, isize[0], index[0], nr_xfr, sumx,
1117 nr_vfr, sumv, bDim, scale, oenv);
1121 write_pdb_bfac(opt2fn("-cf", NFILE, fnm),
1122 opt2fn("-af", NFILE, fnm), "average force", &(top.atoms),
1123 ePBC, topbox, isize[0], index[0], nr_xfr, sumx,
1124 nr_ffr, sumf, bDim, scale, oenv);
1128 view_all(oenv, NFILE, fnm);