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
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);
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++)
132 xav[g][d] = sum[d]/mtot;
137 /* mass=NULL, so these are forces: sum only (do not average) */
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", "" };
229 for(d=0; d<=DIM; d++)
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);
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;
286 for(i=0; i<isize; i++)
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. */
300 for(m=0; (m<DIM); m++)
307 lxx=lxy=lxz=lyy=lyz=lzz=0.0;
308 for(i=0; i<isize; i++)
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) */
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++)
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[])
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])
404 x[i][d] += box[m][d];
407 while (x[i][m]-xp[i][m] > hbox[m])
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,
422 gmx_bool bDim[],real scale_factor,
423 const output_env_t oenv)
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);
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++)
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++)
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 "[TT]g_traj[tt] 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;
665 int flags,nvhisto=0,*vhisto=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 CopyRight(stderr,argv[0]);
706 parse_common_args(&argc,argv,
707 PCA_CAN_TIME | PCA_TIME_UNIT | PCA_CAN_VIEW | PCA_BE_NICE,
708 NFILE,fnm,asize(pa),pa,asize(desc),desc,0,NULL,&oenv);
712 fprintf(stderr,"Interpreting indexfile entries as molecules.\n"
713 "Using center of mass.\n");
716 bOX = opt2bSet("-ox",NFILE,fnm);
717 bOXT = opt2bSet("-oxt",NFILE,fnm);
718 bOV = opt2bSet("-ov",NFILE,fnm);
719 bOF = opt2bSet("-of",NFILE,fnm);
720 bOB = opt2bSet("-ob",NFILE,fnm);
721 bOT = opt2bSet("-ot",NFILE,fnm);
722 bEKT = opt2bSet("-ekt",NFILE,fnm);
723 bEKR = opt2bSet("-ekr",NFILE,fnm);
724 bCV = opt2bSet("-cv",NFILE,fnm) || opt2bSet("-av",NFILE,fnm);
725 bCF = opt2bSet("-cf",NFILE,fnm) || opt2bSet("-af",NFILE,fnm);
726 bVD = opt2bSet("-vd",NFILE,fnm) || opt2parg_bSet("-bin",asize(pa),pa);
727 if (bMol || bOT || bEKT || bEKR)
739 sffmt = "\t" gmx_real_fullprecision_pfmt;
745 sprintf(sffmt6,"%s%s%s%s%s%s",sffmt,sffmt,sffmt,sffmt,sffmt,sffmt);
747 bTop = read_tps_conf(ftp2fn(efTPS,NFILE,fnm),title,&top,&ePBC,
749 bCom && (bOX || bOXT || bOV || bOT || bEKT || bEKR));
751 if ((bMol || bCV || bCF) && !bTop)
753 gmx_fatal(FARGS,"Need a run input file for option -mol, -cv or -cf");
758 indexfn = ftp2fn(efNDX,NFILE,fnm);
762 indexfn = ftp2fn_null(efNDX,NFILE,fnm);
765 if (!(bCom && !bMol))
769 snew(grpname,ngroups);
770 snew(isize0,ngroups);
771 snew(index0,ngroups);
772 get_index(&(top.atoms),indexfn,ngroups,isize0,index0,grpname);
781 for (i=0; i<ngroups; i++)
783 if (index0[0][i] < 0 || index0[0][i] >= mols->nr)
785 gmx_fatal(FARGS,"Molecule index (%d) is out of range (%d-%d)",
786 index0[0][i]+1,1,mols->nr);
788 isize[i] = atndx[index0[0][i]+1] - atndx[index0[0][i]];
789 snew(index[i],isize[i]);
790 for(j=0; j<isize[i]; j++)
792 index[i][j] = atndx[index0[0][i]] + j;
803 snew(mass,top.atoms.nr);
804 for(i=0; i<top.atoms.nr; i++)
806 mass[i] = top.atoms.atom[i].m;
817 flags = flags | TRX_READ_X;
818 outx = xvgropen(opt2fn("-ox",NFILE,fnm),
819 bCom ? "Center of mass" : "Coordinate",
820 output_env_get_xvgr_tlabel(oenv),"Coordinate (nm)",oenv);
821 make_legend(outx,ngroups,isize0[0],index0[0],grpname,bCom,bMol,bDim,oenv);
825 flags = flags | TRX_READ_X;
826 status_out = open_trx(opt2fn("-oxt",NFILE,fnm),"w");
830 flags = flags | TRX_READ_V;
831 outv = xvgropen(opt2fn("-ov",NFILE,fnm),
832 bCom ? "Center of mass velocity" : "Velocity",
833 output_env_get_xvgr_tlabel(oenv),"Velocity (nm/ps)",oenv);
834 make_legend(outv,ngroups,isize0[0],index0[0],grpname,bCom,bMol,bDim,oenv);
837 flags = flags | TRX_READ_F;
838 outf = xvgropen(opt2fn("-of",NFILE,fnm),"Force",
839 output_env_get_xvgr_tlabel(oenv),"Force (kJ mol\\S-1\\N nm\\S-1\\N)",
841 make_legend(outf,ngroups,isize0[0],index0[0],grpname,bCom,bMol,bDim,oenv);
845 outb = xvgropen(opt2fn("-ob",NFILE,fnm),"Box vector elements",
846 output_env_get_xvgr_tlabel(oenv),"(nm)",oenv);
848 xvgr_legend(outb,6,box_leg,oenv);
856 flags = flags | TRX_READ_V;
857 outt = xvgropen(opt2fn("-ot",NFILE,fnm),"Temperature",
858 output_env_get_xvgr_tlabel(oenv),"(K)", oenv);
859 make_legend(outt,ngroups,isize[0],index[0],grpname,bCom,bMol,bDum,oenv);
867 flags = flags | TRX_READ_V;
868 outekt = xvgropen(opt2fn("-ekt",NFILE,fnm),"Center of mass translation",
869 output_env_get_xvgr_tlabel(oenv),"Energy (kJ mol\\S-1\\N)",oenv);
870 make_legend(outekt,ngroups,isize[0],index[0],grpname,bCom,bMol,bDum,oenv);
878 flags = flags | TRX_READ_X | TRX_READ_V;
879 outekr = xvgropen(opt2fn("-ekr",NFILE,fnm),"Center of mass rotation",
880 output_env_get_xvgr_tlabel(oenv),"Energy (kJ mol\\S-1\\N)",oenv);
881 make_legend(outekr,ngroups,isize[0],index[0],grpname,bCom,bMol,bDum,oenv);
885 flags = flags | TRX_READ_V;
889 flags = flags | TRX_READ_X | TRX_READ_V;
893 flags = flags | TRX_READ_X | TRX_READ_F;
895 if ((flags == 0) && !bOB)
897 fprintf(stderr,"Please select one or more output file options\n");
901 read_first_frame(oenv,&status,ftp2fn(efTRX,NFILE,fnm),&fr,flags);
905 snew(sumx,fr.natoms);
909 snew(sumv,fr.natoms);
913 snew(sumf,fr.natoms);
921 gpbc = gmx_rmpbc_init(&top.idef,ePBC,fr.natoms,fr.box);
926 time = output_env_conv_time(oenv,fr.time);
928 if (fr.bX && bNoJump && fr.bBox)
932 remove_jump(fr.box,fr.natoms,xp,fr.x);
938 for(i=0; i<fr.natoms; i++)
940 copy_rvec(fr.x[i],xp[i]);
944 if (fr.bX && bCom && bPBC)
946 gmx_rmpbc_trxfr(gpbc,&fr);
951 update_histo(isize[0],index[0],fr.v,&nvhisto,&vhisto,binwidth);
956 print_data(outx,time,fr.x,mass,bCom,ngroups,isize,index,bDim,sffmt);
963 frout.atoms = &top.atoms;
966 write_trx_x(status_out,&frout,mass,bCom,ngroups,isize,index);
970 print_data(outv,time,fr.v,mass,bCom,ngroups,isize,index,bDim,sffmt);
974 print_data(outf,time,fr.f,NULL,bCom,ngroups,isize,index,bDim,sffmt);
978 fprintf(outb,"\t%g",fr.time);
980 fr.box[XX][XX],fr.box[YY][YY],fr.box[ZZ][ZZ],
981 fr.box[YY][XX],fr.box[ZZ][XX],fr.box[ZZ][YY]);
986 fprintf(outt," %g",time);
987 for(i=0; i<ngroups; i++)
989 fprintf(outt,sffmt,temp(fr.v,mass,isize[i],index[i]));
995 fprintf(outekt," %g",time);
996 for(i=0; i<ngroups; i++)
998 fprintf(outekt,sffmt,ektrans(fr.v,mass,isize[i],index[i]));
1000 fprintf(outekt,"\n");
1002 if (bEKR && fr.bX && fr.bV)
1004 fprintf(outekr," %g",time);
1005 for(i=0; i<ngroups; i++)
1007 fprintf(outekr,sffmt,ekrot(fr.x,fr.v,mass,isize[i],index[i]));
1009 fprintf(outekr,"\n");
1011 if ((bCV || bCF) && fr.bX &&
1012 (ctime < 0 || (fr.time >= ctime*0.999999 &&
1013 fr.time <= ctime*1.000001)))
1015 for(i=0; i<fr.natoms; i++)
1017 rvec_inc(sumx[i],fr.x[i]);
1023 for(i=0; i<fr.natoms; i++)
1025 rvec_inc(sumv[i],fr.v[i]);
1031 for(i=0; i<fr.natoms; i++)
1033 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 */
1048 if (bOX) ffclose(outx);
1049 if (bOXT) close_trx(status_out);
1050 if (bOV) ffclose(outv);
1051 if (bOF) ffclose(outf);
1052 if (bOB) ffclose(outb);
1053 if (bOT) ffclose(outt);
1054 if (bEKT) ffclose(outekt);
1055 if (bEKR) ffclose(outekr);
1059 print_histo(opt2fn("-vd",NFILE,fnm),nvhisto,vhisto,binwidth,oenv);
1066 if (ePBC != epbcNONE && !bNoJump)
1068 fprintf(stderr,"\nWARNING: More than one frame was used for option -cv or -cf\n"
1069 "If atoms jump across the box you should use the -nojump or -ctime option\n\n");
1071 for(i=0; i<isize[0]; i++)
1073 svmul(1.0/nr_xfr,sumx[index[0][i]],sumx[index[0][i]]);
1076 else if (nr_xfr == 0)
1078 fprintf(stderr,"\nWARNING: No coordinate frames found for option -cv or -cf\n\n");
1083 write_pdb_bfac(opt2fn("-cv",NFILE,fnm),
1084 opt2fn("-av",NFILE,fnm),"average velocity",&(top.atoms),
1085 ePBC,topbox,isize[0],index[0],nr_xfr,sumx,
1086 nr_vfr,sumv,bDim,scale,oenv);
1090 write_pdb_bfac(opt2fn("-cf",NFILE,fnm),
1091 opt2fn("-af",NFILE,fnm),"average force",&(top.atoms),
1092 ePBC,topbox,isize[0],index[0],nr_xfr,sumx,
1093 nr_ffr,sumf,bDim,scale,oenv);
1097 view_all(oenv,NFILE, fnm);