1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
11 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13 * Copyright (c) 2001-2004, The GROMACS development team,
14 * check out http://www.gromacs.org for more information.
16 * This program is free software; you can redistribute it and/or
17 * modify it under the terms of the GNU General Public License
18 * as published by the Free Software Foundation; either version 2
19 * of the License, or (at your option) any later version.
21 * If you want to redistribute modifications, please consider that
22 * scientific software is very special. Version control is crucial -
23 * bugs must be traceable. We will be happy to consider code for
24 * inclusion in the official distribution, but derived work must not
25 * be called official GROMACS. Details are found in the README & COPYING
26 * files - if they are missing, get the official version at www.gromacs.org.
28 * To help us fund GROMACS development, we humbly ask that you cite
29 * the papers on the package - you can find them in the top README file.
31 * For more info, check our website at http://www.gromacs.org
34 * Green Red Orange Magenta Azure Cyan Skyblue
62 static void low_print_data(FILE *fp, real time, rvec x[], int n, atom_id *index,
63 gmx_bool bDim[], const char *sffmt)
67 fprintf(fp, " %g", time);
68 for (i = 0; i < n; i++)
78 for (d = 0; d < DIM; d++)
82 fprintf(fp, sffmt, x[ii][d]);
87 fprintf(fp, sffmt, norm(x[ii]));
93 static void average_data(rvec x[], rvec xav[], real *mass,
94 int ngrps, int isize[], atom_id **index)
99 double sum[DIM], mtot;
101 for (g = 0; g < ngrps; g++)
106 for (i = 0; i < isize[g]; i++)
112 svmul(m, x[ind], tmp);
113 for (d = 0; d < DIM; d++)
121 for (d = 0; d < DIM; d++)
129 for (d = 0; d < DIM; d++)
131 xav[g][d] = sum[d]/mtot;
136 /* mass=NULL, so these are forces: sum only (do not average) */
137 for (d = 0; d < DIM; d++)
145 static void print_data(FILE *fp, real time, rvec x[], real *mass, gmx_bool bCom,
146 int ngrps, int isize[], atom_id **index, gmx_bool bDim[],
149 static rvec *xav = NULL;
157 average_data(x, xav, mass, ngrps, isize, index);
158 low_print_data(fp, time, xav, ngrps, NULL, bDim, sffmt);
162 low_print_data(fp, time, x, isize[0], index[0], bDim, sffmt);
166 static void write_trx_x(t_trxstatus *status, t_trxframe *fr, real *mass, gmx_bool bCom,
167 int ngrps, int isize[], atom_id **index)
169 static rvec *xav = NULL;
170 static t_atoms *atoms = NULL;
183 snew(atoms->atom, ngrps);
185 /* Note that atom and residue names will be the ones
186 * of the first atom in each group.
188 for (i = 0; i < ngrps; i++)
190 atoms->atom[i] = fr->atoms->atom[index[i][0]];
191 atoms->atomname[i] = fr->atoms->atomname[index[i][0]];
194 average_data(fr->x, xav, mass, ngrps, isize, index);
196 fr_av.natoms = ngrps;
199 write_trxframe(status, &fr_av, NULL);
203 write_trxframe_indexed(status, fr, isize[0], index[0], NULL);
207 static void make_legend(FILE *fp, int ngrps, int isize, atom_id index[],
208 char **name, gmx_bool bCom, gmx_bool bMol, gmx_bool bDim[],
209 const output_env_t oenv)
212 const char *dimtxt[] = { " X", " Y", " Z", "" };
226 for (i = 0; i < n; i++)
228 for (d = 0; d <= DIM; d++)
232 snew(leg[j], STRLEN);
235 sprintf(leg[j], "mol %d%s", index[i]+1, dimtxt[d]);
239 sprintf(leg[j], "%s%s", name[i], dimtxt[d]);
243 sprintf(leg[j], "atom %d%s", index[i]+1, dimtxt[d]);
249 xvgr_legend(fp, j, (const char**)leg, oenv);
251 for (i = 0; i < j; i++)
258 static real ekrot(rvec x[], rvec v[], real mass[], int isize, atom_id index[])
260 static real **TCM = NULL, **L;
261 double tm, m0, lxx, lxy, lxz, lyy, lyz, lzz, ekrot;
270 for (i = 0; i < DIM; i++)
275 for (i = 0; i < DIM; i++)
285 for (i = 0; i < isize; i++)
290 cprod(x[j], v[j], a0);
291 for (m = 0; (m < DIM); m++)
293 xcm[m] += m0*x[j][m]; /* c.o.m. position */
294 vcm[m] += m0*v[j][m]; /* c.o.m. velocity */
295 acm[m] += m0*a0[m]; /* rotational velocity around c.o.m. */
298 dcprod(xcm, vcm, b0);
299 for (m = 0; (m < DIM); m++)
306 lxx = lxy = lxz = lyy = lyz = lzz = 0.0;
307 for (i = 0; i < isize; i++)
311 for (m = 0; m < DIM; m++)
313 dx[m] = x[j][m] - xcm[m];
315 lxx += dx[XX]*dx[XX]*m0;
316 lxy += dx[XX]*dx[YY]*m0;
317 lxz += dx[XX]*dx[ZZ]*m0;
318 lyy += dx[YY]*dx[YY]*m0;
319 lyz += dx[YY]*dx[ZZ]*m0;
320 lzz += dx[ZZ]*dx[ZZ]*m0;
323 L[XX][XX] = lyy + lzz;
327 L[YY][YY] = lxx + lzz;
331 L[ZZ][ZZ] = lxx + lyy;
333 m_inv_gen(L, DIM, TCM);
335 /* Compute omega (hoeksnelheid) */
338 for (m = 0; m < DIM; m++)
340 for (n = 0; n < DIM; n++)
342 ocm[m] += TCM[m][n]*acm[n];
344 ekrot += 0.5*ocm[m]*acm[m];
350 static real ektrans(rvec v[], real mass[], int isize, atom_id index[])
358 for (i = 0; i < isize; i++)
361 for (d = 0; d < DIM; d++)
363 mvcom[d] += mass[j]*v[j][d];
368 return dnorm2(mvcom)/(mtot*2);
371 static real temp(rvec v[], real mass[], int isize, atom_id index[])
377 for (i = 0; i < isize; i++)
380 ekin2 += mass[j]*norm2(v[j]);
383 return ekin2/(3*isize*BOLTZ);
386 static void remove_jump(matrix box, int natoms, rvec xp[], rvec x[])
391 for (d = 0; d < DIM; d++)
393 hbox[d] = 0.5*box[d][d];
395 for (i = 0; i < natoms; i++)
397 for (m = DIM-1; m >= 0; m--)
399 while (x[i][m]-xp[i][m] <= -hbox[m])
401 for (d = 0; d <= m; d++)
403 x[i][d] += box[m][d];
406 while (x[i][m]-xp[i][m] > hbox[m])
408 for (d = 0; d <= m; d++)
410 x[i][d] -= box[m][d];
417 static void write_pdb_bfac(const char *fname, const char *xname,
418 const char *title, t_atoms *atoms, int ePBC, matrix box,
419 int isize, atom_id *index, int nfr_x, rvec *x,
420 int nfr_v, rvec *sum,
421 gmx_bool bDim[], real scale_factor,
422 const output_env_t oenv)
425 real max, len2, scale;
430 if ((nfr_x == 0) || (nfr_v == 0))
432 fprintf(stderr, "No frames found for %s, will not write %s\n",
437 fprintf(stderr, "Used %d frames for %s\n", nfr_x, "coordinates");
438 fprintf(stderr, "Used %d frames for %s\n", nfr_v, title);
443 for (i = 0; i < DIM; i++)
457 for (i = 0; i < isize; i++)
459 svmul(scale, sum[index[i]], sum[index[i]]);
462 fp = xvgropen(xname, title, "Atom", "", oenv);
463 for (i = 0; i < isize; i++)
465 fprintf(fp, "%-5d %10.3f %10.3f %10.3f\n", 1+i,
466 sum[index[i]][XX], sum[index[i]][YY], sum[index[i]][ZZ]);
471 for (i = 0; i < isize; i++)
474 for (m = 0; m < DIM; m++)
476 if (bDim[m] || bDim[DIM])
478 len2 += sqr(sum[index[i]][m]);
487 if (scale_factor != 0)
489 scale = scale_factor;
499 scale = 10.0/sqrt(max);
503 printf("Maximum %s is %g on atom %d %s, res. %s %d\n",
504 title, sqrt(max), maxi+1, *(atoms->atomname[maxi]),
505 *(atoms->resinfo[atoms->atom[maxi].resind].name),
506 atoms->resinfo[atoms->atom[maxi].resind].nr);
508 if (atoms->pdbinfo == NULL)
510 snew(atoms->pdbinfo, atoms->nr);
514 for (i = 0; i < isize; i++)
517 for (m = 0; m < DIM; m++)
519 if (bDim[m] || bDim[DIM])
521 len2 += sqr(sum[index[i]][m]);
524 atoms->pdbinfo[index[i]].bfac = sqrt(len2)*scale;
529 for (i = 0; i < isize; i++)
531 atoms->pdbinfo[index[i]].bfac = sum[index[i]][onedim]*scale;
534 write_sto_conf_indexed(fname, title, atoms, x, NULL, ePBC, box, isize, index);
538 static void update_histo(int gnx, atom_id index[], rvec v[],
539 int *nhisto, int **histo, real binwidth)
547 for (i = 0; (i < gnx); i++)
549 vn = norm(v[index[i]]);
550 vnmax = max(vn, vnmax);
553 *nhisto = 1+(vnmax/binwidth);
554 snew(*histo, *nhisto);
556 for (i = 0; (i < gnx); i++)
558 vn = norm(v[index[i]]);
563 fprintf(stderr, "Extending histogram from %d to %d\n", *nhisto, nnn);
566 for (m = *nhisto; (m < nnn); m++)
576 static void print_histo(const char *fn, int nhisto, int histo[], real binwidth,
577 const output_env_t oenv)
582 fp = xvgropen(fn, "Velocity distribution", "V (nm/ps)", "arbitrary units",
584 for (i = 0; (i < nhisto); i++)
586 fprintf(fp, "%10.3e %10d\n", i*binwidth, histo[i]);
591 int gmx_traj(int argc, char *argv[])
593 const char *desc[] = {
594 "[TT]g_traj[tt] plots coordinates, velocities, forces and/or the box.",
595 "With [TT]-com[tt] the coordinates, velocities and forces are",
596 "calculated for the center of mass of each group.",
597 "When [TT]-mol[tt] is set, the numbers in the index file are",
598 "interpreted as molecule numbers and the same procedure as with",
599 "[TT]-com[tt] is used for each molecule.[PAR]",
600 "Option [TT]-ot[tt] plots the temperature of each group,",
601 "provided velocities are present in the trajectory file.",
602 "No corrections are made for constrained degrees of freedom!",
603 "This implies [TT]-com[tt].[PAR]",
604 "Options [TT]-ekt[tt] and [TT]-ekr[tt] plot the translational and",
605 "rotational kinetic energy of each group,",
606 "provided velocities are present in the trajectory file.",
607 "This implies [TT]-com[tt].[PAR]",
608 "Options [TT]-cv[tt] and [TT]-cf[tt] write the average velocities",
609 "and average forces as temperature factors to a [TT].pdb[tt] file with",
610 "the average coordinates or the coordinates at [TT]-ctime[tt].",
611 "The temperature factors are scaled such that the maximum is 10.",
612 "The scaling can be changed with the option [TT]-scale[tt].",
613 "To get the velocities or forces of one",
614 "frame set both [TT]-b[tt] and [TT]-e[tt] to the time of",
615 "desired frame. When averaging over frames you might need to use",
616 "the [TT]-nojump[tt] option to obtain the correct average coordinates.",
617 "If you select either of these option the average force and velocity",
618 "for each atom are written to an [TT].xvg[tt] file as well",
619 "(specified with [TT]-av[tt] or [TT]-af[tt]).[PAR]",
620 "Option [TT]-vd[tt] computes a velocity distribution, i.e. the",
621 "norm of the vector is plotted. In addition in the same graph",
622 "the kinetic energy distribution is given."
624 static gmx_bool bMol = FALSE, bCom = FALSE, bPBC = TRUE, bNoJump = FALSE;
625 static gmx_bool bX = TRUE, bY = TRUE, bZ = TRUE, bNorm = FALSE, bFP = FALSE;
626 static int ngroups = 1;
627 static real ctime = -1, scale = 0, binwidth = 1;
629 { "-com", FALSE, etBOOL, {&bCom},
630 "Plot data for the com of each group" },
631 { "-pbc", FALSE, etBOOL, {&bPBC},
632 "Make molecules whole for COM" },
633 { "-mol", FALSE, etBOOL, {&bMol},
634 "Index contains molecule numbers iso atom numbers" },
635 { "-nojump", FALSE, etBOOL, {&bNoJump},
636 "Remove jumps of atoms across the box" },
637 { "-x", FALSE, etBOOL, {&bX},
638 "Plot X-component" },
639 { "-y", FALSE, etBOOL, {&bY},
640 "Plot Y-component" },
641 { "-z", FALSE, etBOOL, {&bZ},
642 "Plot Z-component" },
643 { "-ng", FALSE, etINT, {&ngroups},
644 "Number of groups to consider" },
645 { "-len", FALSE, etBOOL, {&bNorm},
646 "Plot vector length" },
647 { "-fp", FALSE, etBOOL, {&bFP},
648 "Full precision output" },
649 { "-bin", FALSE, etREAL, {&binwidth},
650 "Binwidth for velocity histogram (nm/ps)" },
651 { "-ctime", FALSE, etREAL, {&ctime},
652 "Use frame at this time for x in [TT]-cv[tt] and [TT]-cf[tt] instead of the average x" },
653 { "-scale", FALSE, etREAL, {&scale},
654 "Scale factor for [TT].pdb[tt] output, 0 is autoscale" }
656 FILE *outx = NULL, *outv = NULL, *outf = NULL, *outb = NULL, *outt = NULL;
657 FILE *outekt = NULL, *outekr = NULL;
663 t_trxframe fr, frout;
664 int flags, nvhisto = 0, *vhisto = NULL;
665 rvec *xtop, *xp = NULL;
666 rvec *sumx = NULL, *sumv = NULL, *sumf = NULL;
669 t_trxstatus *status_out = NULL;
670 gmx_rmpbc_t gpbc = NULL;
672 int nr_xfr, nr_vfr, nr_ffr;
675 atom_id **index0, **index;
678 gmx_bool bTop, bOX, bOXT, bOV, bOF, bOB, bOT, bEKT, bEKR, bCV, bCF;
679 gmx_bool bDim[4], bDum[4], bVD;
680 char *sffmt, sffmt6[1024];
681 const char *box_leg[6] = { "XX", "YY", "ZZ", "YX", "ZX", "ZY" };
685 { efTRX, "-f", NULL, ffREAD },
686 { efTPS, NULL, NULL, ffREAD },
687 { efNDX, NULL, NULL, ffOPTRD },
688 { efXVG, "-ox", "coord.xvg", ffOPTWR },
689 { efTRX, "-oxt", "coord.xtc", ffOPTWR },
690 { efXVG, "-ov", "veloc.xvg", ffOPTWR },
691 { efXVG, "-of", "force.xvg", ffOPTWR },
692 { efXVG, "-ob", "box.xvg", ffOPTWR },
693 { efXVG, "-ot", "temp.xvg", ffOPTWR },
694 { efXVG, "-ekt", "ektrans.xvg", ffOPTWR },
695 { efXVG, "-ekr", "ekrot.xvg", ffOPTWR },
696 { efXVG, "-vd", "veldist.xvg", ffOPTWR },
697 { efPDB, "-cv", "veloc.pdb", ffOPTWR },
698 { efPDB, "-cf", "force.pdb", ffOPTWR },
699 { efXVG, "-av", "all_veloc.xvg", ffOPTWR },
700 { efXVG, "-af", "all_force.xvg", ffOPTWR }
702 #define NFILE asize(fnm)
704 parse_common_args(&argc, argv,
705 PCA_CAN_TIME | PCA_TIME_UNIT | PCA_CAN_VIEW | PCA_BE_NICE,
706 NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL, &oenv);
710 fprintf(stderr, "Interpreting indexfile entries as molecules.\n"
711 "Using center of mass.\n");
714 bOX = opt2bSet("-ox", NFILE, fnm);
715 bOXT = opt2bSet("-oxt", NFILE, fnm);
716 bOV = opt2bSet("-ov", NFILE, fnm);
717 bOF = opt2bSet("-of", NFILE, fnm);
718 bOB = opt2bSet("-ob", NFILE, fnm);
719 bOT = opt2bSet("-ot", NFILE, fnm);
720 bEKT = opt2bSet("-ekt", NFILE, fnm);
721 bEKR = opt2bSet("-ekr", NFILE, fnm);
722 bCV = opt2bSet("-cv", NFILE, fnm) || opt2bSet("-av", NFILE, fnm);
723 bCF = opt2bSet("-cf", NFILE, fnm) || opt2bSet("-af", NFILE, fnm);
724 bVD = opt2bSet("-vd", NFILE, fnm) || opt2parg_bSet("-bin", asize(pa), pa);
725 if (bMol || bOT || bEKT || bEKR)
737 sffmt = "\t" gmx_real_fullprecision_pfmt;
743 sprintf(sffmt6, "%s%s%s%s%s%s", sffmt, sffmt, sffmt, sffmt, sffmt, sffmt);
745 bTop = read_tps_conf(ftp2fn(efTPS, NFILE, fnm), title, &top, &ePBC,
747 bCom && (bOX || bOXT || bOV || bOT || bEKT || bEKR));
749 if ((bMol || bCV || bCF) && !bTop)
751 gmx_fatal(FARGS, "Need a run input file for option -mol, -cv or -cf");
756 indexfn = ftp2fn(efNDX, NFILE, fnm);
760 indexfn = ftp2fn_null(efNDX, NFILE, fnm);
763 if (!(bCom && !bMol))
767 snew(grpname, ngroups);
768 snew(isize0, ngroups);
769 snew(index0, ngroups);
770 get_index(&(top.atoms), indexfn, ngroups, isize0, index0, grpname);
777 snew(isize, ngroups);
778 snew(index, ngroups);
779 for (i = 0; i < ngroups; i++)
781 if (index0[0][i] < 0 || index0[0][i] >= mols->nr)
783 gmx_fatal(FARGS, "Molecule index (%d) is out of range (%d-%d)",
784 index0[0][i]+1, 1, mols->nr);
786 isize[i] = atndx[index0[0][i]+1] - atndx[index0[0][i]];
787 snew(index[i], isize[i]);
788 for (j = 0; j < isize[i]; j++)
790 index[i][j] = atndx[index0[0][i]] + j;
801 snew(mass, top.atoms.nr);
802 for (i = 0; i < top.atoms.nr; i++)
804 mass[i] = top.atoms.atom[i].m;
815 flags = flags | TRX_READ_X;
816 outx = xvgropen(opt2fn("-ox", NFILE, fnm),
817 bCom ? "Center of mass" : "Coordinate",
818 output_env_get_xvgr_tlabel(oenv), "Coordinate (nm)", oenv);
819 make_legend(outx, ngroups, isize0[0], index0[0], grpname, bCom, bMol, bDim, oenv);
823 flags = flags | TRX_READ_X;
824 status_out = open_trx(opt2fn("-oxt", NFILE, fnm), "w");
828 flags = flags | TRX_READ_V;
829 outv = xvgropen(opt2fn("-ov", NFILE, fnm),
830 bCom ? "Center of mass velocity" : "Velocity",
831 output_env_get_xvgr_tlabel(oenv), "Velocity (nm/ps)", oenv);
832 make_legend(outv, ngroups, isize0[0], index0[0], grpname, bCom, bMol, bDim, oenv);
836 flags = flags | TRX_READ_F;
837 outf = xvgropen(opt2fn("-of", NFILE, fnm), "Force",
838 output_env_get_xvgr_tlabel(oenv), "Force (kJ mol\\S-1\\N nm\\S-1\\N)",
840 make_legend(outf, ngroups, isize0[0], index0[0], grpname, bCom, bMol, bDim, oenv);
844 outb = xvgropen(opt2fn("-ob", NFILE, fnm), "Box vector elements",
845 output_env_get_xvgr_tlabel(oenv), "(nm)", oenv);
847 xvgr_legend(outb, 6, box_leg, oenv);
855 flags = flags | TRX_READ_V;
856 outt = xvgropen(opt2fn("-ot", NFILE, fnm), "Temperature",
857 output_env_get_xvgr_tlabel(oenv), "(K)", oenv);
858 make_legend(outt, ngroups, isize[0], index[0], grpname, bCom, bMol, bDum, oenv);
866 flags = flags | TRX_READ_V;
867 outekt = xvgropen(opt2fn("-ekt", NFILE, fnm), "Center of mass translation",
868 output_env_get_xvgr_tlabel(oenv), "Energy (kJ mol\\S-1\\N)", oenv);
869 make_legend(outekt, ngroups, isize[0], index[0], grpname, bCom, bMol, bDum, oenv);
877 flags = flags | TRX_READ_X | TRX_READ_V;
878 outekr = xvgropen(opt2fn("-ekr", NFILE, fnm), "Center of mass rotation",
879 output_env_get_xvgr_tlabel(oenv), "Energy (kJ mol\\S-1\\N)", oenv);
880 make_legend(outekr, ngroups, isize[0], index[0], grpname, bCom, bMol, bDum, oenv);
884 flags = flags | TRX_READ_V;
888 flags = flags | TRX_READ_X | TRX_READ_V;
892 flags = flags | TRX_READ_X | TRX_READ_F;
894 if ((flags == 0) && !bOB)
896 fprintf(stderr, "Please select one or more output file options\n");
900 read_first_frame(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &fr, flags);
904 snew(sumx, fr.natoms);
908 snew(sumv, fr.natoms);
912 snew(sumf, fr.natoms);
920 gpbc = gmx_rmpbc_init(&top.idef, ePBC, fr.natoms);
925 time = output_env_conv_time(oenv, fr.time);
927 if (fr.bX && bNoJump && fr.bBox)
931 remove_jump(fr.box, fr.natoms, xp, fr.x);
937 for (i = 0; i < fr.natoms; i++)
939 copy_rvec(fr.x[i], xp[i]);
943 if (fr.bX && bCom && bPBC)
945 gmx_rmpbc_trxfr(gpbc, &fr);
950 update_histo(isize[0], index[0], fr.v, &nvhisto, &vhisto, binwidth);
955 print_data(outx, time, fr.x, mass, bCom, ngroups, isize, index, bDim, sffmt);
962 frout.atoms = &top.atoms;
965 write_trx_x(status_out, &frout, mass, bCom, ngroups, isize, index);
969 print_data(outv, time, fr.v, mass, bCom, ngroups, isize, index, bDim, sffmt);
973 print_data(outf, time, fr.f, NULL, bCom, ngroups, isize, index, bDim, sffmt);
977 fprintf(outb, "\t%g", fr.time);
978 fprintf(outb, sffmt6,
979 fr.box[XX][XX], fr.box[YY][YY], fr.box[ZZ][ZZ],
980 fr.box[YY][XX], fr.box[ZZ][XX], fr.box[ZZ][YY]);
985 fprintf(outt, " %g", time);
986 for (i = 0; i < ngroups; i++)
988 fprintf(outt, sffmt, temp(fr.v, mass, isize[i], index[i]));
994 fprintf(outekt, " %g", time);
995 for (i = 0; i < ngroups; i++)
997 fprintf(outekt, sffmt, ektrans(fr.v, mass, isize[i], index[i]));
999 fprintf(outekt, "\n");
1001 if (bEKR && fr.bX && fr.bV)
1003 fprintf(outekr, " %g", time);
1004 for (i = 0; i < ngroups; i++)
1006 fprintf(outekr, sffmt, ekrot(fr.x, fr.v, mass, isize[i], index[i]));
1008 fprintf(outekr, "\n");
1010 if ((bCV || bCF) && fr.bX &&
1011 (ctime < 0 || (fr.time >= ctime*0.999999 &&
1012 fr.time <= ctime*1.000001)))
1014 for (i = 0; i < fr.natoms; i++)
1016 rvec_inc(sumx[i], fr.x[i]);
1022 for (i = 0; i < fr.natoms; i++)
1024 rvec_inc(sumv[i], fr.v[i]);
1030 for (i = 0; i < fr.natoms; i++)
1032 rvec_inc(sumf[i], fr.f[i]);
1038 while (read_next_frame(oenv, status, &fr));
1042 gmx_rmpbc_done(gpbc);
1045 /* clean up a bit */
1054 close_trx(status_out);
1083 print_histo(opt2fn("-vd", NFILE, fnm), nvhisto, vhisto, binwidth, oenv);
1090 if (ePBC != epbcNONE && !bNoJump)
1092 fprintf(stderr, "\nWARNING: More than one frame was used for option -cv or -cf\n"
1093 "If atoms jump across the box you should use the -nojump or -ctime option\n\n");
1095 for (i = 0; i < isize[0]; i++)
1097 svmul(1.0/nr_xfr, sumx[index[0][i]], sumx[index[0][i]]);
1100 else if (nr_xfr == 0)
1102 fprintf(stderr, "\nWARNING: No coordinate frames found for option -cv or -cf\n\n");
1107 write_pdb_bfac(opt2fn("-cv", NFILE, fnm),
1108 opt2fn("-av", NFILE, fnm), "average velocity", &(top.atoms),
1109 ePBC, topbox, isize[0], index[0], nr_xfr, sumx,
1110 nr_vfr, sumv, bDim, scale, oenv);
1114 write_pdb_bfac(opt2fn("-cf", NFILE, fnm),
1115 opt2fn("-af", NFILE, fnm), "average force", &(top.atoms),
1116 ePBC, topbox, isize[0], index[0], nr_xfr, sumx,
1117 nr_ffr, sumf, bDim, scale, oenv);
1121 view_all(oenv, NFILE, fnm);