Remove sysstuff.h
[alexxy/gromacs.git] / src / gromacs / gmxana / gmx_traj.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 #ifdef HAVE_CONFIG_H
38 #include <config.h>
39 #endif
40
41 #include <math.h>
42 #include <stdlib.h>
43 #include <string.h>
44
45 #include "gromacs/commandline/pargs.h"
46 #include "typedefs.h"
47 #include "gromacs/utility/smalloc.h"
48 #include "macros.h"
49 #include "vec.h"
50 #include "pbc.h"
51 #include "gromacs/utility/futil.h"
52 #include "index.h"
53 #include "mshift.h"
54 #include "xvgr.h"
55 #include "gromacs/fileio/tpxio.h"
56 #include "gromacs/fileio/trxio.h"
57 #include "rmpbc.h"
58 #include "physics.h"
59 #include "gromacs/fileio/confio.h"
60 #include "gmx_ana.h"
61
62 #include "gromacs/linearalgebra/nrjac.h"
63 #include "gromacs/utility/fatalerror.h"
64
65 static void low_print_data(FILE *fp, real time, rvec x[], int n, atom_id *index,
66                            gmx_bool bDim[], const char *sffmt)
67 {
68     int i, ii, d;
69
70     fprintf(fp, " %g", time);
71     for (i = 0; i < n; i++)
72     {
73         if (index != NULL)
74         {
75             ii = index[i];
76         }
77         else
78         {
79             ii = i;
80         }
81         for (d = 0; d < DIM; d++)
82         {
83             if (bDim[d])
84             {
85                 fprintf(fp, sffmt, x[ii][d]);
86             }
87         }
88         if (bDim[DIM])
89         {
90             fprintf(fp, sffmt, norm(x[ii]));
91         }
92     }
93     fprintf(fp, "\n");
94 }
95
96 static void average_data(rvec x[], rvec xav[], real *mass,
97                          int ngrps, int isize[], atom_id **index)
98 {
99     int    g, i, ind, d;
100     real   m;
101     rvec   tmp;
102     double sum[DIM], mtot;
103
104     for (g = 0; g < ngrps; g++)
105     {
106         clear_dvec(sum);
107         clear_rvec(xav[g]);
108         mtot = 0;
109         for (i = 0; i < isize[g]; i++)
110         {
111             ind = index[g][i];
112             if (mass != NULL)
113             {
114                 m = mass[ind];
115                 svmul(m, x[ind], tmp);
116                 for (d = 0; d < DIM; d++)
117                 {
118                     sum[d] += tmp[d];
119                 }
120                 mtot += m;
121             }
122             else
123             {
124                 for (d = 0; d < DIM; d++)
125                 {
126                     sum[d] += x[ind][d];
127                 }
128             }
129         }
130         if (mass != NULL)
131         {
132             for (d = 0; d < DIM; d++)
133             {
134                 xav[g][d] = sum[d]/mtot;
135             }
136         }
137         else
138         {
139             /* mass=NULL, so these are forces: sum only (do not average) */
140             for (d = 0; d < DIM; d++)
141             {
142                 xav[g][d] = sum[d];
143             }
144         }
145     }
146 }
147
148 static void print_data(FILE *fp, real time, rvec x[], real *mass, gmx_bool bCom,
149                        int ngrps, int isize[], atom_id **index, gmx_bool bDim[],
150                        const char *sffmt)
151 {
152     static rvec *xav = NULL;
153
154     if (bCom)
155     {
156         if (xav == NULL)
157         {
158             snew(xav, ngrps);
159         }
160         average_data(x, xav, mass, ngrps, isize, index);
161         low_print_data(fp, time, xav, ngrps, NULL, bDim, sffmt);
162     }
163     else
164     {
165         low_print_data(fp, time, x, isize[0], index[0], bDim, sffmt);
166     }
167 }
168
169 static void write_trx_x(t_trxstatus *status, t_trxframe *fr, real *mass, gmx_bool bCom,
170                         int ngrps, int isize[], atom_id **index)
171 {
172     static rvec    *xav   = NULL;
173     static t_atoms *atoms = NULL;
174     t_trxframe      fr_av;
175     int             i;
176
177     fr->bV = FALSE;
178     fr->bF = FALSE;
179     if (bCom)
180     {
181         if (xav == NULL)
182         {
183             snew(xav, ngrps);
184             snew(atoms, 1);
185             *atoms = *fr->atoms;
186             snew(atoms->atom, ngrps);
187             atoms->nr = ngrps;
188             /* Note that atom and residue names will be the ones
189              * of the first atom in each group.
190              */
191             for (i = 0; i < ngrps; i++)
192             {
193                 atoms->atom[i]     = fr->atoms->atom[index[i][0]];
194                 atoms->atomname[i] = fr->atoms->atomname[index[i][0]];
195             }
196         }
197         average_data(fr->x, xav, mass, ngrps, isize, index);
198         fr_av        = *fr;
199         fr_av.natoms = ngrps;
200         fr_av.atoms  = atoms;
201         fr_av.x      = xav;
202         write_trxframe(status, &fr_av, NULL);
203     }
204     else
205     {
206         write_trxframe_indexed(status, fr, isize[0], index[0], NULL);
207     }
208 }
209
210 static void make_legend(FILE *fp, int ngrps, int isize, atom_id index[],
211                         char **name, gmx_bool bCom, gmx_bool bMol, gmx_bool bDim[],
212                         const output_env_t oenv)
213 {
214     char      **leg;
215     const char *dimtxt[] = { " X", " Y", " Z", "" };
216     int         n, i, j, d;
217
218     if (bCom)
219     {
220         n = ngrps;
221     }
222     else
223     {
224         n = isize;
225     }
226
227     snew(leg, 4*n);
228     j = 0;
229     for (i = 0; i < n; i++)
230     {
231         for (d = 0; d <= DIM; d++)
232         {
233             if (bDim[d])
234             {
235                 snew(leg[j], STRLEN);
236                 if (bMol)
237                 {
238                     sprintf(leg[j], "mol %d%s", index[i]+1, dimtxt[d]);
239                 }
240                 else if (bCom)
241                 {
242                     sprintf(leg[j], "%s%s", name[i], dimtxt[d]);
243                 }
244                 else
245                 {
246                     sprintf(leg[j], "atom %d%s", index[i]+1, dimtxt[d]);
247                 }
248                 j++;
249             }
250         }
251     }
252     xvgr_legend(fp, j, (const char**)leg, oenv);
253
254     for (i = 0; i < j; i++)
255     {
256         sfree(leg[i]);
257     }
258     sfree(leg);
259 }
260
261 static real ekrot(rvec x[], rvec v[], real mass[], int isize, atom_id index[])
262 {
263     static real **TCM = NULL, **L;
264     double        tm, m0, lxx, lxy, lxz, lyy, lyz, lzz, ekrot;
265     rvec          a0, ocm;
266     dvec          dx, b0;
267     dvec          xcm, vcm, acm;
268     int           i, j, m, n;
269
270     if (TCM == NULL)
271     {
272         snew(TCM, DIM);
273         for (i = 0; i < DIM; i++)
274         {
275             snew(TCM[i], DIM);
276         }
277         snew(L, DIM);
278         for (i = 0; i < DIM; i++)
279         {
280             snew(L[i], DIM);
281         }
282     }
283
284     clear_dvec(xcm);
285     clear_dvec(vcm);
286     clear_dvec(acm);
287     tm = 0.0;
288     for (i = 0; i < isize; i++)
289     {
290         j   = index[i];
291         m0  = mass[j];
292         tm += m0;
293         cprod(x[j], v[j], a0);
294         for (m = 0; (m < DIM); m++)
295         {
296             xcm[m] += m0*x[j][m]; /* c.o.m. position */
297             vcm[m] += m0*v[j][m]; /* c.o.m. velocity */
298             acm[m] += m0*a0[m];   /* rotational velocity around c.o.m. */
299         }
300     }
301     dcprod(xcm, vcm, b0);
302     for (m = 0; (m < DIM); m++)
303     {
304         xcm[m] /= tm;
305         vcm[m] /= tm;
306         acm[m] -= b0[m]/tm;
307     }
308
309     lxx = lxy = lxz = lyy = lyz = lzz = 0.0;
310     for (i = 0; i < isize; i++)
311     {
312         j  = index[i];
313         m0 = mass[j];
314         for (m = 0; m < DIM; m++)
315         {
316             dx[m] = x[j][m] - xcm[m];
317         }
318         lxx += dx[XX]*dx[XX]*m0;
319         lxy += dx[XX]*dx[YY]*m0;
320         lxz += dx[XX]*dx[ZZ]*m0;
321         lyy += dx[YY]*dx[YY]*m0;
322         lyz += dx[YY]*dx[ZZ]*m0;
323         lzz += dx[ZZ]*dx[ZZ]*m0;
324     }
325
326     L[XX][XX] =  lyy + lzz;
327     L[YY][XX] = -lxy;
328     L[ZZ][XX] = -lxz;
329     L[XX][YY] = -lxy;
330     L[YY][YY] =  lxx + lzz;
331     L[ZZ][YY] = -lyz;
332     L[XX][ZZ] = -lxz;
333     L[YY][ZZ] = -lyz;
334     L[ZZ][ZZ] =  lxx + lyy;
335
336     m_inv_gen(L, DIM, TCM);
337
338     /* Compute omega (hoeksnelheid) */
339     clear_rvec(ocm);
340     ekrot = 0;
341     for (m = 0; m < DIM; m++)
342     {
343         for (n = 0; n < DIM; n++)
344         {
345             ocm[m] += TCM[m][n]*acm[n];
346         }
347         ekrot += 0.5*ocm[m]*acm[m];
348     }
349
350     return ekrot;
351 }
352
353 static real ektrans(rvec v[], real mass[], int isize, atom_id index[])
354 {
355     dvec   mvcom;
356     double mtot;
357     int    i, j, d;
358
359     clear_dvec(mvcom);
360     mtot = 0;
361     for (i = 0; i < isize; i++)
362     {
363         j = index[i];
364         for (d = 0; d < DIM; d++)
365         {
366             mvcom[d] += mass[j]*v[j][d];
367         }
368         mtot += mass[j];
369     }
370
371     return dnorm2(mvcom)/(mtot*2);
372 }
373
374 static real temp(rvec v[], real mass[], int isize, atom_id index[])
375 {
376     double ekin2;
377     int    i, j;
378
379     ekin2 = 0;
380     for (i = 0; i < isize; i++)
381     {
382         j      = index[i];
383         ekin2 += mass[j]*norm2(v[j]);
384     }
385
386     return ekin2/(3*isize*BOLTZ);
387 }
388
389 static void remove_jump(matrix box, int natoms, rvec xp[], rvec x[])
390 {
391     rvec hbox;
392     int  d, i, m;
393
394     for (d = 0; d < DIM; d++)
395     {
396         hbox[d] = 0.5*box[d][d];
397     }
398     for (i = 0; i < natoms; i++)
399     {
400         for (m = DIM-1; m >= 0; m--)
401         {
402             while (x[i][m]-xp[i][m] <= -hbox[m])
403             {
404                 for (d = 0; d <= m; d++)
405                 {
406                     x[i][d] += box[m][d];
407                 }
408             }
409             while (x[i][m]-xp[i][m] > hbox[m])
410             {
411                 for (d = 0; d <= m; d++)
412                 {
413                     x[i][d] -= box[m][d];
414                 }
415             }
416         }
417     }
418 }
419
420 static void write_pdb_bfac(const char *fname, const char *xname,
421                            const char *title, t_atoms *atoms, int ePBC, matrix box,
422                            int isize, atom_id *index, int nfr_x, rvec *x,
423                            int nfr_v, rvec *sum,
424                            gmx_bool bDim[], real scale_factor,
425                            const output_env_t oenv)
426 {
427     FILE       *fp;
428     real        max, len2, scale;
429     atom_id     maxi;
430     int         i, m, onedim;
431     gmx_bool    bOne;
432
433     if ((nfr_x == 0) || (nfr_v == 0))
434     {
435         fprintf(stderr, "No frames found for %s, will not write %s\n",
436                 title, fname);
437     }
438     else
439     {
440         fprintf(stderr, "Used %d frames for %s\n", nfr_x, "coordinates");
441         fprintf(stderr, "Used %d frames for %s\n", nfr_v, title);
442         onedim = -1;
443         if (!bDim[DIM])
444         {
445             m = 0;
446             for (i = 0; i < DIM; i++)
447             {
448                 if (bDim[i])
449                 {
450                     onedim = i;
451                     m++;
452                 }
453             }
454             if (m != 1)
455             {
456                 onedim = -1;
457             }
458         }
459         scale = 1.0/nfr_v;
460         for (i = 0; i < isize; i++)
461         {
462             svmul(scale, sum[index[i]], sum[index[i]]);
463         }
464
465         fp = xvgropen(xname, title, "Atom", "", oenv);
466         for (i = 0; i < isize; i++)
467         {
468             fprintf(fp, "%-5d  %10.3f  %10.3f  %10.3f\n", 1+i,
469                     sum[index[i]][XX], sum[index[i]][YY], sum[index[i]][ZZ]);
470         }
471         gmx_ffclose(fp);
472         max  = 0;
473         maxi = 0;
474         for (i = 0; i < isize; i++)
475         {
476             len2 = 0;
477             for (m = 0; m < DIM; m++)
478             {
479                 if (bDim[m] || bDim[DIM])
480                 {
481                     len2 += sqr(sum[index[i]][m]);
482                 }
483             }
484             if (len2 > max)
485             {
486                 max  = len2;
487                 maxi = index[i];
488             }
489         }
490         if (scale_factor != 0)
491         {
492             scale = scale_factor;
493         }
494         else
495         {
496             if (max == 0)
497             {
498                 scale = 1;
499             }
500             else
501             {
502                 scale = 10.0/sqrt(max);
503             }
504         }
505
506         printf("Maximum %s is %g on atom %d %s, res. %s %d\n",
507                title, sqrt(max), maxi+1, *(atoms->atomname[maxi]),
508                *(atoms->resinfo[atoms->atom[maxi].resind].name),
509                atoms->resinfo[atoms->atom[maxi].resind].nr);
510
511         if (atoms->pdbinfo == NULL)
512         {
513             snew(atoms->pdbinfo, atoms->nr);
514         }
515         if (onedim == -1)
516         {
517             for (i = 0; i < isize; i++)
518             {
519                 len2 = 0;
520                 for (m = 0; m < DIM; m++)
521                 {
522                     if (bDim[m] || bDim[DIM])
523                     {
524                         len2 += sqr(sum[index[i]][m]);
525                     }
526                 }
527                 atoms->pdbinfo[index[i]].bfac = sqrt(len2)*scale;
528             }
529         }
530         else
531         {
532             for (i = 0; i < isize; i++)
533             {
534                 atoms->pdbinfo[index[i]].bfac = sum[index[i]][onedim]*scale;
535             }
536         }
537         write_sto_conf_indexed(fname, title, atoms, x, NULL, ePBC, box, isize, index);
538     }
539 }
540
541 static void update_histo(int gnx, atom_id index[], rvec v[],
542                          int *nhisto, int **histo, real binwidth)
543 {
544     int  i, m, in, nnn;
545     real vn, vnmax;
546
547     if (*histo == NULL)
548     {
549         vnmax = 0;
550         for (i = 0; (i < gnx); i++)
551         {
552             vn    = norm(v[index[i]]);
553             vnmax = max(vn, vnmax);
554         }
555         vnmax  *= 2;
556         *nhisto = 1+(vnmax/binwidth);
557         snew(*histo, *nhisto);
558     }
559     for (i = 0; (i < gnx); i++)
560     {
561         vn = norm(v[index[i]]);
562         in = vn/binwidth;
563         if (in >= *nhisto)
564         {
565             nnn = in+100;
566             fprintf(stderr, "Extending histogram from %d to %d\n", *nhisto, nnn);
567
568             srenew(*histo, nnn);
569             for (m = *nhisto; (m < nnn); m++)
570             {
571                 (*histo)[m] = 0;
572             }
573             *nhisto = nnn;
574         }
575         (*histo)[in]++;
576     }
577 }
578
579 static void print_histo(const char *fn, int nhisto, int histo[], real binwidth,
580                         const output_env_t oenv)
581 {
582     FILE *fp;
583     int   i;
584
585     fp = xvgropen(fn, "Velocity distribution", "V (nm/ps)", "arbitrary units",
586                   oenv);
587     for (i = 0; (i < nhisto); i++)
588     {
589         fprintf(fp, "%10.3e  %10d\n", i*binwidth, histo[i]);
590     }
591     gmx_ffclose(fp);
592 }
593
594 int gmx_traj(int argc, char *argv[])
595 {
596     const char     *desc[] = {
597         "[THISMODULE] plots coordinates, velocities, forces and/or the box.",
598         "With [TT]-com[tt] the coordinates, velocities and forces are",
599         "calculated for the center of mass of each group.",
600         "When [TT]-mol[tt] is set, the numbers in the index file are",
601         "interpreted as molecule numbers and the same procedure as with",
602         "[TT]-com[tt] is used for each molecule.[PAR]",
603         "Option [TT]-ot[tt] plots the temperature of each group,",
604         "provided velocities are present in the trajectory file.",
605         "No corrections are made for constrained degrees of freedom!",
606         "This implies [TT]-com[tt].[PAR]",
607         "Options [TT]-ekt[tt] and [TT]-ekr[tt] plot the translational and",
608         "rotational kinetic energy of each group,",
609         "provided velocities are present in the trajectory file.",
610         "This implies [TT]-com[tt].[PAR]",
611         "Options [TT]-cv[tt] and [TT]-cf[tt] write the average velocities",
612         "and average forces as temperature factors to a [TT].pdb[tt] file with",
613         "the average coordinates or the coordinates at [TT]-ctime[tt].",
614         "The temperature factors are scaled such that the maximum is 10.",
615         "The scaling can be changed with the option [TT]-scale[tt].",
616         "To get the velocities or forces of one",
617         "frame set both [TT]-b[tt] and [TT]-e[tt] to the time of",
618         "desired frame. When averaging over frames you might need to use",
619         "the [TT]-nojump[tt] option to obtain the correct average coordinates.",
620         "If you select either of these option the average force and velocity",
621         "for each atom are written to an [TT].xvg[tt] file as well",
622         "(specified with [TT]-av[tt] or [TT]-af[tt]).[PAR]",
623         "Option [TT]-vd[tt] computes a velocity distribution, i.e. the",
624         "norm of the vector is plotted. In addition in the same graph",
625         "the kinetic energy distribution is given."
626     };
627     static gmx_bool bMol    = FALSE, bCom = FALSE, bPBC = TRUE, bNoJump = FALSE;
628     static gmx_bool bX      = TRUE, bY = TRUE, bZ = TRUE, bNorm = FALSE, bFP = FALSE;
629     static int      ngroups = 1;
630     static real     ctime   = -1, scale = 0, binwidth = 1;
631     t_pargs         pa[]    = {
632         { "-com", FALSE, etBOOL, {&bCom},
633           "Plot data for the com of each group" },
634         { "-pbc", FALSE, etBOOL, {&bPBC},
635           "Make molecules whole for COM" },
636         { "-mol", FALSE, etBOOL, {&bMol},
637           "Index contains molecule numbers iso atom numbers" },
638         { "-nojump", FALSE, etBOOL, {&bNoJump},
639           "Remove jumps of atoms across the box" },
640         { "-x", FALSE, etBOOL, {&bX},
641           "Plot X-component" },
642         { "-y", FALSE, etBOOL, {&bY},
643           "Plot Y-component" },
644         { "-z", FALSE, etBOOL, {&bZ},
645           "Plot Z-component" },
646         { "-ng",       FALSE, etINT, {&ngroups},
647           "Number of groups to consider" },
648         { "-len", FALSE, etBOOL, {&bNorm},
649           "Plot vector length" },
650         { "-fp", FALSE, etBOOL, {&bFP},
651           "Full precision output" },
652         { "-bin", FALSE, etREAL, {&binwidth},
653           "Binwidth for velocity histogram (nm/ps)" },
654         { "-ctime", FALSE, etREAL, {&ctime},
655           "Use frame at this time for x in [TT]-cv[tt] and [TT]-cf[tt] instead of the average x" },
656         { "-scale", FALSE, etREAL, {&scale},
657           "Scale factor for [TT].pdb[tt] output, 0 is autoscale" }
658     };
659     FILE           *outx   = NULL, *outv = NULL, *outf = NULL, *outb = NULL, *outt = NULL;
660     FILE           *outekt = NULL, *outekr = NULL;
661     t_topology      top;
662     int             ePBC;
663     real           *mass, time;
664     char            title[STRLEN];
665     const char     *indexfn;
666     t_trxframe      fr, frout;
667     int             flags, nvhisto = 0, *vhisto = NULL;
668     rvec           *xtop, *xp = NULL;
669     rvec           *sumx = NULL, *sumv = NULL, *sumf = NULL;
670     matrix          topbox;
671     t_trxstatus    *status;
672     t_trxstatus    *status_out = NULL;
673     gmx_rmpbc_t     gpbc       = NULL;
674     int             i, j, n;
675     int             nr_xfr, nr_vfr, nr_ffr;
676     char          **grpname;
677     int            *isize0, *isize;
678     atom_id       **index0, **index;
679     atom_id        *atndx;
680     t_block        *mols;
681     gmx_bool        bTop, bOX, bOXT, bOV, bOF, bOB, bOT, bEKT, bEKR, bCV, bCF;
682     gmx_bool        bDim[4], bDum[4], bVD;
683     char           *sffmt, sffmt6[1024];
684     const char     *box_leg[6] = { "XX", "YY", "ZZ", "YX", "ZX", "ZY" };
685     output_env_t    oenv;
686
687     t_filenm        fnm[] = {
688         { efTRX, "-f", NULL, ffREAD },
689         { efTPS, NULL, NULL, ffREAD },
690         { efNDX, NULL, NULL, ffOPTRD },
691         { efXVG, "-ox", "coord.xvg", ffOPTWR },
692         { efTRX, "-oxt", "coord.xtc", ffOPTWR },
693         { efXVG, "-ov", "veloc.xvg", ffOPTWR },
694         { efXVG, "-of", "force.xvg", ffOPTWR },
695         { efXVG, "-ob", "box.xvg",   ffOPTWR },
696         { efXVG, "-ot", "temp.xvg",  ffOPTWR },
697         { efXVG, "-ekt", "ektrans.xvg", ffOPTWR },
698         { efXVG, "-ekr", "ekrot.xvg", ffOPTWR },
699         { efXVG, "-vd", "veldist.xvg", ffOPTWR },
700         { efPDB, "-cv", "veloc.pdb", ffOPTWR },
701         { efPDB, "-cf", "force.pdb", ffOPTWR },
702         { efXVG, "-av", "all_veloc.xvg", ffOPTWR },
703         { efXVG, "-af", "all_force.xvg", ffOPTWR }
704     };
705 #define NFILE asize(fnm)
706
707     if (!parse_common_args(&argc, argv,
708                            PCA_CAN_TIME | PCA_TIME_UNIT | PCA_CAN_VIEW | PCA_BE_NICE,
709                            NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL, &oenv))
710     {
711         return 0;
712     }
713
714     if (bMol)
715     {
716         fprintf(stderr, "Interpreting indexfile entries as molecules.\n"
717                 "Using center of mass.\n");
718     }
719
720     bOX  = opt2bSet("-ox", NFILE, fnm);
721     bOXT = opt2bSet("-oxt", NFILE, fnm);
722     bOV  = opt2bSet("-ov", NFILE, fnm);
723     bOF  = opt2bSet("-of", NFILE, fnm);
724     bOB  = opt2bSet("-ob", NFILE, fnm);
725     bOT  = opt2bSet("-ot", NFILE, fnm);
726     bEKT = opt2bSet("-ekt", NFILE, fnm);
727     bEKR = opt2bSet("-ekr", NFILE, fnm);
728     bCV  = opt2bSet("-cv", NFILE, fnm) || opt2bSet("-av", NFILE, fnm);
729     bCF  = opt2bSet("-cf", NFILE, fnm) || opt2bSet("-af", NFILE, fnm);
730     bVD  = opt2bSet("-vd", NFILE, fnm) || opt2parg_bSet("-bin", asize(pa), pa);
731     if (bMol || bOT || bEKT || bEKR)
732     {
733         bCom = TRUE;
734     }
735
736     bDim[XX]  = bX;
737     bDim[YY]  = bY;
738     bDim[ZZ]  = bZ;
739     bDim[DIM] = bNorm;
740
741     if (bFP)
742     {
743         sffmt = "\t" gmx_real_fullprecision_pfmt;
744     }
745     else
746     {
747         sffmt = "\t%g";
748     }
749     sprintf(sffmt6, "%s%s%s%s%s%s", sffmt, sffmt, sffmt, sffmt, sffmt, sffmt);
750
751     bTop = read_tps_conf(ftp2fn(efTPS, NFILE, fnm), title, &top, &ePBC,
752                          &xtop, NULL, topbox,
753                          bCom && (bOX || bOXT || bOV || bOT || bEKT || bEKR));
754     sfree(xtop);
755     if ((bMol || bCV || bCF) && !bTop)
756     {
757         gmx_fatal(FARGS, "Need a run input file for option -mol, -cv or -cf");
758     }
759
760     if (bMol)
761     {
762         indexfn = ftp2fn(efNDX, NFILE, fnm);
763     }
764     else
765     {
766         indexfn = ftp2fn_null(efNDX, NFILE, fnm);
767     }
768
769     if (!(bCom && !bMol))
770     {
771         ngroups = 1;
772     }
773     snew(grpname, ngroups);
774     snew(isize0, ngroups);
775     snew(index0, ngroups);
776     get_index(&(top.atoms), indexfn, ngroups, isize0, index0, grpname);
777
778     if (bMol)
779     {
780         mols    = &(top.mols);
781         atndx   = mols->index;
782         ngroups = isize0[0];
783         snew(isize, ngroups);
784         snew(index, ngroups);
785         for (i = 0; i < ngroups; i++)
786         {
787             if (index0[0][i] < 0 || index0[0][i] >= mols->nr)
788             {
789                 gmx_fatal(FARGS, "Molecule index (%d) is out of range (%d-%d)",
790                           index0[0][i]+1, 1, mols->nr);
791             }
792             isize[i] = atndx[index0[0][i]+1] - atndx[index0[0][i]];
793             snew(index[i], isize[i]);
794             for (j = 0; j < isize[i]; j++)
795             {
796                 index[i][j] = atndx[index0[0][i]] + j;
797             }
798         }
799     }
800     else
801     {
802         isize = isize0;
803         index = index0;
804     }
805     if (bCom)
806     {
807         snew(mass, top.atoms.nr);
808         for (i = 0; i < top.atoms.nr; i++)
809         {
810             mass[i] = top.atoms.atom[i].m;
811         }
812     }
813     else
814     {
815         mass = NULL;
816     }
817
818     flags = 0;
819     if (bOX)
820     {
821         flags = flags | TRX_READ_X;
822         outx  = xvgropen(opt2fn("-ox", NFILE, fnm),
823                          bCom ? "Center of mass" : "Coordinate",
824                          output_env_get_xvgr_tlabel(oenv), "Coordinate (nm)", oenv);
825         make_legend(outx, ngroups, isize0[0], index0[0], grpname, bCom, bMol, bDim, oenv);
826     }
827     if (bOXT)
828     {
829         flags      = flags | TRX_READ_X;
830         status_out = open_trx(opt2fn("-oxt", NFILE, fnm), "w");
831     }
832     if (bOV)
833     {
834         flags = flags | TRX_READ_V;
835         outv  = xvgropen(opt2fn("-ov", NFILE, fnm),
836                          bCom ? "Center of mass velocity" : "Velocity",
837                          output_env_get_xvgr_tlabel(oenv), "Velocity (nm/ps)", oenv);
838         make_legend(outv, ngroups, isize0[0], index0[0], grpname, bCom, bMol, bDim, oenv);
839     }
840     if (bOF)
841     {
842         flags = flags | TRX_READ_F;
843         outf  = xvgropen(opt2fn("-of", NFILE, fnm), "Force",
844                          output_env_get_xvgr_tlabel(oenv), "Force (kJ mol\\S-1\\N nm\\S-1\\N)",
845                          oenv);
846         make_legend(outf, ngroups, isize0[0], index0[0], grpname, bCom, bMol, bDim, oenv);
847     }
848     if (bOB)
849     {
850         outb = xvgropen(opt2fn("-ob", NFILE, fnm), "Box vector elements",
851                         output_env_get_xvgr_tlabel(oenv), "(nm)", oenv);
852
853         xvgr_legend(outb, 6, box_leg, oenv);
854     }
855     if (bOT)
856     {
857         bDum[XX]  = FALSE;
858         bDum[YY]  = FALSE;
859         bDum[ZZ]  = FALSE;
860         bDum[DIM] = TRUE;
861         flags     = flags | TRX_READ_V;
862         outt      = xvgropen(opt2fn("-ot", NFILE, fnm), "Temperature",
863                              output_env_get_xvgr_tlabel(oenv), "(K)", oenv);
864         make_legend(outt, ngroups, isize[0], index[0], grpname, bCom, bMol, bDum, oenv);
865     }
866     if (bEKT)
867     {
868         bDum[XX]  = FALSE;
869         bDum[YY]  = FALSE;
870         bDum[ZZ]  = FALSE;
871         bDum[DIM] = TRUE;
872         flags     = flags | TRX_READ_V;
873         outekt    = xvgropen(opt2fn("-ekt", NFILE, fnm), "Center of mass translation",
874                              output_env_get_xvgr_tlabel(oenv), "Energy (kJ mol\\S-1\\N)", oenv);
875         make_legend(outekt, ngroups, isize[0], index[0], grpname, bCom, bMol, bDum, oenv);
876     }
877     if (bEKR)
878     {
879         bDum[XX]  = FALSE;
880         bDum[YY]  = FALSE;
881         bDum[ZZ]  = FALSE;
882         bDum[DIM] = TRUE;
883         flags     = flags | TRX_READ_X | TRX_READ_V;
884         outekr    = xvgropen(opt2fn("-ekr", NFILE, fnm), "Center of mass rotation",
885                              output_env_get_xvgr_tlabel(oenv), "Energy (kJ mol\\S-1\\N)", oenv);
886         make_legend(outekr, ngroups, isize[0], index[0], grpname, bCom, bMol, bDum, oenv);
887     }
888     if (bVD)
889     {
890         flags = flags | TRX_READ_V;
891     }
892     if (bCV)
893     {
894         flags = flags | TRX_READ_X | TRX_READ_V;
895     }
896     if (bCF)
897     {
898         flags = flags | TRX_READ_X | TRX_READ_F;
899     }
900     if ((flags == 0) && !bOB)
901     {
902         fprintf(stderr, "Please select one or more output file options\n");
903         exit(0);
904     }
905
906     read_first_frame(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &fr, flags);
907
908     if (bCV || bCF)
909     {
910         snew(sumx, fr.natoms);
911     }
912     if (bCV)
913     {
914         snew(sumv, fr.natoms);
915     }
916     if (bCF)
917     {
918         snew(sumf, fr.natoms);
919     }
920     nr_xfr = 0;
921     nr_vfr = 0;
922     nr_ffr = 0;
923
924     if (bCom && bPBC)
925     {
926         gpbc = gmx_rmpbc_init(&top.idef, ePBC, fr.natoms);
927     }
928
929     do
930     {
931         time = output_env_conv_time(oenv, fr.time);
932
933         if (fr.bX && bNoJump && fr.bBox)
934         {
935             if (xp)
936             {
937                 remove_jump(fr.box, fr.natoms, xp, fr.x);
938             }
939             else
940             {
941                 snew(xp, fr.natoms);
942             }
943             for (i = 0; i < fr.natoms; i++)
944             {
945                 copy_rvec(fr.x[i], xp[i]);
946             }
947         }
948
949         if (fr.bX && bCom && bPBC)
950         {
951             gmx_rmpbc_trxfr(gpbc, &fr);
952         }
953
954         if (bVD && fr.bV)
955         {
956             update_histo(isize[0], index[0], fr.v, &nvhisto, &vhisto, binwidth);
957         }
958
959         if (bOX && fr.bX)
960         {
961             print_data(outx, time, fr.x, mass, bCom, ngroups, isize, index, bDim, sffmt);
962         }
963         if (bOXT && fr.bX)
964         {
965             frout = fr;
966             if (!frout.bAtoms)
967             {
968                 frout.atoms  = &top.atoms;
969                 frout.bAtoms = TRUE;
970             }
971             write_trx_x(status_out, &frout, mass, bCom, ngroups, isize, index);
972         }
973         if (bOV && fr.bV)
974         {
975             print_data(outv, time, fr.v, mass, bCom, ngroups, isize, index, bDim, sffmt);
976         }
977         if (bOF && fr.bF)
978         {
979             print_data(outf, time, fr.f, NULL, bCom, ngroups, isize, index, bDim, sffmt);
980         }
981         if (bOB && fr.bBox)
982         {
983             fprintf(outb, "\t%g", fr.time);
984             fprintf(outb, sffmt6,
985                     fr.box[XX][XX], fr.box[YY][YY], fr.box[ZZ][ZZ],
986                     fr.box[YY][XX], fr.box[ZZ][XX], fr.box[ZZ][YY]);
987             fprintf(outb, "\n");
988         }
989         if (bOT && fr.bV)
990         {
991             fprintf(outt, " %g", time);
992             for (i = 0; i < ngroups; i++)
993             {
994                 fprintf(outt, sffmt, temp(fr.v, mass, isize[i], index[i]));
995             }
996             fprintf(outt, "\n");
997         }
998         if (bEKT && fr.bV)
999         {
1000             fprintf(outekt, " %g", time);
1001             for (i = 0; i < ngroups; i++)
1002             {
1003                 fprintf(outekt, sffmt, ektrans(fr.v, mass, isize[i], index[i]));
1004             }
1005             fprintf(outekt, "\n");
1006         }
1007         if (bEKR && fr.bX && fr.bV)
1008         {
1009             fprintf(outekr, " %g", time);
1010             for (i = 0; i < ngroups; i++)
1011             {
1012                 fprintf(outekr, sffmt, ekrot(fr.x, fr.v, mass, isize[i], index[i]));
1013             }
1014             fprintf(outekr, "\n");
1015         }
1016         if ((bCV || bCF) && fr.bX &&
1017             (ctime < 0 || (fr.time >= ctime*0.999999 &&
1018                            fr.time <= ctime*1.000001)))
1019         {
1020             for (i = 0; i < fr.natoms; i++)
1021             {
1022                 rvec_inc(sumx[i], fr.x[i]);
1023             }
1024             nr_xfr++;
1025         }
1026         if (bCV && fr.bV)
1027         {
1028             for (i = 0; i < fr.natoms; i++)
1029             {
1030                 rvec_inc(sumv[i], fr.v[i]);
1031             }
1032             nr_vfr++;
1033         }
1034         if (bCF && fr.bF)
1035         {
1036             for (i = 0; i < fr.natoms; i++)
1037             {
1038                 rvec_inc(sumf[i], fr.f[i]);
1039             }
1040             nr_ffr++;
1041         }
1042
1043     }
1044     while (read_next_frame(oenv, status, &fr));
1045
1046     if (gpbc != NULL)
1047     {
1048         gmx_rmpbc_done(gpbc);
1049     }
1050
1051     /* clean up a bit */
1052     close_trj(status);
1053
1054     if (bOX)
1055     {
1056         gmx_ffclose(outx);
1057     }
1058     if (bOXT)
1059     {
1060         close_trx(status_out);
1061     }
1062     if (bOV)
1063     {
1064         gmx_ffclose(outv);
1065     }
1066     if (bOF)
1067     {
1068         gmx_ffclose(outf);
1069     }
1070     if (bOB)
1071     {
1072         gmx_ffclose(outb);
1073     }
1074     if (bOT)
1075     {
1076         gmx_ffclose(outt);
1077     }
1078     if (bEKT)
1079     {
1080         gmx_ffclose(outekt);
1081     }
1082     if (bEKR)
1083     {
1084         gmx_ffclose(outekr);
1085     }
1086
1087     if (bVD)
1088     {
1089         print_histo(opt2fn("-vd", NFILE, fnm), nvhisto, vhisto, binwidth, oenv);
1090     }
1091
1092     if (bCV || bCF)
1093     {
1094         if (nr_xfr > 1)
1095         {
1096             if (ePBC != epbcNONE && !bNoJump)
1097             {
1098                 fprintf(stderr, "\nWARNING: More than one frame was used for option -cv or -cf\n"
1099                         "If atoms jump across the box you should use the -nojump or -ctime option\n\n");
1100             }
1101             for (i = 0; i < isize[0]; i++)
1102             {
1103                 svmul(1.0/nr_xfr, sumx[index[0][i]], sumx[index[0][i]]);
1104             }
1105         }
1106         else if (nr_xfr == 0)
1107         {
1108             fprintf(stderr, "\nWARNING: No coordinate frames found for option -cv or -cf\n\n");
1109         }
1110     }
1111     if (bCV)
1112     {
1113         write_pdb_bfac(opt2fn("-cv", NFILE, fnm),
1114                        opt2fn("-av", NFILE, fnm), "average velocity", &(top.atoms),
1115                        ePBC, topbox, isize[0], index[0], nr_xfr, sumx,
1116                        nr_vfr, sumv, bDim, scale, oenv);
1117     }
1118     if (bCF)
1119     {
1120         write_pdb_bfac(opt2fn("-cf", NFILE, fnm),
1121                        opt2fn("-af", NFILE, fnm), "average force", &(top.atoms),
1122                        ePBC, topbox, isize[0], index[0], nr_xfr, sumx,
1123                        nr_ffr, sumf, bDim, scale, oenv);
1124     }
1125
1126     /* view it */
1127     view_all(oenv, NFILE, fnm);
1128
1129     return 0;
1130 }