Merge branch release-4-6 into master
[alexxy/gromacs.git] / src / gromacs / gmxana / gmx_mindist.c
1 /*
2  *
3  *                This source code is part of
4  *
5  *                 G   R   O   M   A   C   S
6  *
7  *          GROningen MAchine for Chemical Simulations
8  *
9  *                        VERSION 3.2.0
10  * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12  * Copyright (c) 2001-2004, The GROMACS development team,
13  * check out http://www.gromacs.org for more information.
14
15  * This program is free software; you can redistribute it and/or
16  * modify it under the terms of the GNU General Public License
17  * as published by the Free Software Foundation; either version 2
18  * of the License, or (at your option) any later version.
19  *
20  * If you want to redistribute modifications, please consider that
21  * scientific software is very special. Version control is crucial -
22  * bugs must be traceable. We will be happy to consider code for
23  * inclusion in the official distribution, but derived work must not
24  * be called official GROMACS. Details are found in the README & COPYING
25  * files - if they are missing, get the official version at www.gromacs.org.
26  *
27  * To help us fund GROMACS development, we humbly ask that you cite
28  * the papers on the package - you can find them in the top README file.
29  *
30  * For more info, check our website at http://www.gromacs.org
31  *
32  * And Hey:
33  * Green Red Orange Magenta Azure Cyan Skyblue
34  */
35 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
38
39 #include <math.h>
40 #include <stdlib.h>
41
42 #include "sysstuff.h"
43 #include <string.h>
44 #include "typedefs.h"
45 #include "smalloc.h"
46 #include "macros.h"
47 #include "vec.h"
48 #include "xvgr.h"
49 #include "pbc.h"
50 #include "futil.h"
51 #include "statutil.h"
52 #include "index.h"
53 #include "tpxio.h"
54 #include "rmpbc.h"
55 #include "xtcio.h"
56 #include "gmx_ana.h"
57
58
59 static void periodic_dist(matrix box, rvec x[], int n, atom_id index[],
60                           real *rmin, real *rmax, int *min_ind)
61 {
62 #define NSHIFT 26
63     int  sx, sy, sz, i, j, s;
64     real sqr_box, r2min, r2max, r2;
65     rvec shift[NSHIFT], d0, d;
66
67     sqr_box = sqr(min(norm(box[XX]), min(norm(box[YY]), norm(box[ZZ]))));
68
69     s = 0;
70     for (sz = -1; sz <= 1; sz++)
71     {
72         for (sy = -1; sy <= 1; sy++)
73         {
74             for (sx = -1; sx <= 1; sx++)
75             {
76                 if (sx != 0 || sy != 0 || sz != 0)
77                 {
78                     for (i = 0; i < DIM; i++)
79                     {
80                         shift[s][i] = sx*box[XX][i]+sy*box[YY][i]+sz*box[ZZ][i];
81                     }
82                     s++;
83                 }
84             }
85         }
86     }
87
88     r2min = sqr_box;
89     r2max = 0;
90
91     for (i = 0; i < n; i++)
92     {
93         for (j = i+1; j < n; j++)
94         {
95             rvec_sub(x[index[i]], x[index[j]], d0);
96             r2 = norm2(d0);
97             if (r2 > r2max)
98             {
99                 r2max = r2;
100             }
101             for (s = 0; s < NSHIFT; s++)
102             {
103                 rvec_add(d0, shift[s], d);
104                 r2 = norm2(d);
105                 if (r2 < r2min)
106                 {
107                     r2min      = r2;
108                     min_ind[0] = i;
109                     min_ind[1] = j;
110                 }
111             }
112         }
113     }
114
115     *rmin = sqrt(r2min);
116     *rmax = sqrt(r2max);
117 }
118
119 static void periodic_mindist_plot(const char *trxfn, const char *outfn,
120                                   t_topology *top, int ePBC,
121                                   int n, atom_id index[], gmx_bool bSplit,
122                                   const output_env_t oenv)
123 {
124     FILE        *out;
125     const char  *leg[5] = { "min per.", "max int.", "box1", "box2", "box3" };
126     t_trxstatus *status;
127     real         t;
128     rvec        *x;
129     matrix       box;
130     int          natoms, ind_min[2] = {0, 0}, ind_mini = 0, ind_minj = 0;
131     real         r, rmin, rmax, rmint, tmint;
132     gmx_bool     bFirst;
133     gmx_rmpbc_t  gpbc = NULL;
134
135     natoms = read_first_x(oenv, &status, trxfn, &t, &x, box);
136
137     check_index(NULL, n, index, NULL, natoms);
138
139     out = xvgropen(outfn, "Minimum distance to periodic image",
140                    output_env_get_time_label(oenv), "Distance (nm)", oenv);
141     if (output_env_get_print_xvgr_codes(oenv))
142     {
143         fprintf(out, "@ subtitle \"and maximum internal distance\"\n");
144     }
145     xvgr_legend(out, 5, leg, oenv);
146
147     rmint = box[XX][XX];
148     tmint = 0;
149
150     if (NULL != top)
151     {
152         gpbc = gmx_rmpbc_init(&top->idef, ePBC, natoms);
153     }
154
155     bFirst = TRUE;
156     do
157     {
158         if (NULL != top)
159         {
160             gmx_rmpbc(gpbc, natoms, box, x);
161         }
162
163         periodic_dist(box, x, n, index, &rmin, &rmax, ind_min);
164         if (rmin < rmint)
165         {
166             rmint    = rmin;
167             tmint    = t;
168             ind_mini = ind_min[0];
169             ind_minj = ind_min[1];
170         }
171         if (bSplit && !bFirst && abs(t/output_env_get_time_factor(oenv)) < 1e-5)
172         {
173             fprintf(out, "&\n");
174         }
175         fprintf(out, "\t%g\t%6.3f %6.3f %6.3f %6.3f %6.3f\n",
176                 output_env_conv_time(oenv, t), rmin, rmax, norm(box[0]), norm(box[1]), norm(box[2]));
177         bFirst = FALSE;
178     }
179     while (read_next_x(oenv, status, &t, x, box));
180
181     if (NULL != top)
182     {
183         gmx_rmpbc_done(gpbc);
184     }
185
186     ffclose(out);
187
188     fprintf(stdout,
189             "\nThe shortest periodic distance is %g (nm) at time %g (%s),\n"
190             "between atoms %d and %d\n",
191             rmint, output_env_conv_time(oenv, tmint), output_env_get_time_unit(oenv),
192             index[ind_mini]+1, index[ind_minj]+1);
193 }
194
195 static void calc_dist(real rcut, gmx_bool bPBC, int ePBC, matrix box, rvec x[],
196                       int nx1, int nx2, atom_id index1[], atom_id index2[],
197                       gmx_bool bGroup,
198                       real *rmin, real *rmax, int *nmin, int *nmax,
199                       int *ixmin, int *jxmin, int *ixmax, int *jxmax)
200 {
201     int      i, j, i0 = 0, j1;
202     int      ix, jx;
203     atom_id *index3;
204     rvec     dx;
205     real     r2, rmin2, rmax2, rcut2;
206     t_pbc    pbc;
207     int      nmin_j, nmax_j;
208
209     *ixmin = -1;
210     *jxmin = -1;
211     *ixmax = -1;
212     *jxmax = -1;
213     *nmin  = 0;
214     *nmax  = 0;
215
216     rcut2 = sqr(rcut);
217
218     /* Must init pbc every step because of pressure coupling */
219     if (bPBC)
220     {
221         set_pbc(&pbc, ePBC, box);
222     }
223     if (index2)
224     {
225         i0     = 0;
226         j1     = nx2;
227         index3 = index2;
228     }
229     else
230     {
231         j1     = nx1;
232         index3 = index1;
233     }
234
235     rmin2 = 1e12;
236     rmax2 = -1e12;
237
238     for (j = 0; (j < j1); j++)
239     {
240         jx = index3[j];
241         if (index2 == NULL)
242         {
243             i0 = j + 1;
244         }
245         nmin_j = 0;
246         nmax_j = 0;
247         for (i = i0; (i < nx1); i++)
248         {
249             ix = index1[i];
250             if (ix != jx)
251             {
252                 if (bPBC)
253                 {
254                     pbc_dx(&pbc, x[ix], x[jx], dx);
255                 }
256                 else
257                 {
258                     rvec_sub(x[ix], x[jx], dx);
259                 }
260                 r2 = iprod(dx, dx);
261                 if (r2 < rmin2)
262                 {
263                     rmin2  = r2;
264                     *ixmin = ix;
265                     *jxmin = jx;
266                 }
267                 if (r2 > rmax2)
268                 {
269                     rmax2  = r2;
270                     *ixmax = ix;
271                     *jxmax = jx;
272                 }
273                 if (r2 <= rcut2)
274                 {
275                     nmin_j++;
276                 }
277                 else if (r2 > rcut2)
278                 {
279                     nmax_j++;
280                 }
281             }
282         }
283         if (bGroup)
284         {
285             if (nmin_j > 0)
286             {
287                 (*nmin)++;
288             }
289             if (nmax_j > 0)
290             {
291                 (*nmax)++;
292             }
293         }
294         else
295         {
296             *nmin += nmin_j;
297             *nmax += nmax_j;
298         }
299     }
300     *rmin = sqrt(rmin2);
301     *rmax = sqrt(rmax2);
302 }
303
304 void dist_plot(const char *fn, const char *afile, const char *dfile,
305                const char *nfile, const char *rfile, const char *xfile,
306                real rcut, gmx_bool bMat, t_atoms *atoms,
307                int ng, atom_id *index[], int gnx[], char *grpn[], gmx_bool bSplit,
308                gmx_bool bMin, int nres, atom_id *residue, gmx_bool bPBC, int ePBC,
309                gmx_bool bGroup, gmx_bool bEachResEachTime, gmx_bool bPrintResName,
310                const output_env_t oenv)
311 {
312     FILE            *atm, *dist, *num;
313     t_trxstatus     *trxout;
314     char             buf[256];
315     char           **leg;
316     real             t, dmin, dmax, **mindres = NULL, **maxdres = NULL;
317     int              nmin, nmax;
318     t_trxstatus     *status;
319     int              i = -1, j, k, natoms;
320     int              min1, min2, max1, max2, min1r, min2r, max1r, max2r;
321     atom_id          oindex[2];
322     rvec            *x0;
323     matrix           box;
324     t_trxframe       frout;
325     gmx_bool         bFirst;
326     FILE            *respertime = NULL;
327
328     if ((natoms = read_first_x(oenv, &status, fn, &t, &x0, box)) == 0)
329     {
330         gmx_fatal(FARGS, "Could not read coordinates from statusfile\n");
331     }
332
333     sprintf(buf, "%simum Distance", bMin ? "Min" : "Max");
334     dist = xvgropen(dfile, buf, output_env_get_time_label(oenv), "Distance (nm)", oenv);
335     sprintf(buf, "Number of Contacts %s %g nm", bMin ? "<" : ">", rcut);
336     num    = nfile ? xvgropen(nfile, buf, output_env_get_time_label(oenv), "Number", oenv) : NULL;
337     atm    = afile ? ffopen(afile, "w") : NULL;
338     trxout = xfile ? open_trx(xfile, "w") : NULL;
339
340     if (bMat)
341     {
342         if (ng == 1)
343         {
344             snew(leg, 1);
345             sprintf(buf, "Internal in %s", grpn[0]);
346             leg[0] = strdup(buf);
347             xvgr_legend(dist, 0, (const char**)leg, oenv);
348             if (num)
349             {
350                 xvgr_legend(num, 0, (const char**)leg, oenv);
351             }
352         }
353         else
354         {
355             snew(leg, (ng*(ng-1))/2);
356             for (i = j = 0; (i < ng-1); i++)
357             {
358                 for (k = i+1; (k < ng); k++, j++)
359                 {
360                     sprintf(buf, "%s-%s", grpn[i], grpn[k]);
361                     leg[j] = strdup(buf);
362                 }
363             }
364             xvgr_legend(dist, j, (const char**)leg, oenv);
365             if (num)
366             {
367                 xvgr_legend(num, j, (const char**)leg, oenv);
368             }
369         }
370     }
371     else
372     {
373         snew(leg, ng-1);
374         for (i = 0; (i < ng-1); i++)
375         {
376             sprintf(buf, "%s-%s", grpn[0], grpn[i+1]);
377             leg[i] = strdup(buf);
378         }
379         xvgr_legend(dist, ng-1, (const char**)leg, oenv);
380         if (num)
381         {
382             xvgr_legend(num, ng-1, (const char**)leg, oenv);
383         }
384     }
385
386     if (bEachResEachTime)
387     {
388         sprintf(buf, "%simum Distance", bMin ? "Min" : "Max");
389         respertime = xvgropen(rfile, buf, output_env_get_time_label(oenv), "Distance (nm)", oenv);
390         xvgr_legend(respertime, ng-1, (const char**)leg, oenv);
391         if (bPrintResName)
392         {
393             fprintf(respertime, "# ");
394         }
395         for (j = 0; j < nres; j++)
396         {
397             fprintf(respertime, "%s%d ", *(atoms->resinfo[atoms->atom[index[0][residue[j]]].resind].name), atoms->atom[index[0][residue[j]]].resind);
398         }
399         fprintf(respertime, "\n");
400
401     }
402
403     j = 0;
404     if (nres)
405     {
406         snew(mindres, ng-1);
407         snew(maxdres, ng-1);
408         for (i = 1; i < ng; i++)
409         {
410             snew(mindres[i-1], nres);
411             snew(maxdres[i-1], nres);
412             for (j = 0; j < nres; j++)
413             {
414                 mindres[i-1][j] = 1e6;
415             }
416             /* maxdres[*][*] is already 0 */
417         }
418     }
419     bFirst = TRUE;
420     do
421     {
422         if (bSplit && !bFirst && abs(t/output_env_get_time_factor(oenv)) < 1e-5)
423         {
424             fprintf(dist, "&\n");
425             if (num)
426             {
427                 fprintf(num, "&\n");
428             }
429             if (atm)
430             {
431                 fprintf(atm, "&\n");
432             }
433         }
434         fprintf(dist, "%12e", output_env_conv_time(oenv, t));
435         if (num)
436         {
437             fprintf(num, "%12e", output_env_conv_time(oenv, t));
438         }
439
440         if (bMat)
441         {
442             if (ng == 1)
443             {
444                 calc_dist(rcut, bPBC, ePBC, box, x0, gnx[0], gnx[0], index[0], index[0], bGroup,
445                           &dmin, &dmax, &nmin, &nmax, &min1, &min2, &max1, &max2);
446                 fprintf(dist, "  %12e", bMin ? dmin : dmax);
447                 if (num)
448                 {
449                     fprintf(num, "  %8d", bMin ? nmin : nmax);
450                 }
451             }
452             else
453             {
454                 for (i = 0; (i < ng-1); i++)
455                 {
456                     for (k = i+1; (k < ng); k++)
457                     {
458                         calc_dist(rcut, bPBC, ePBC, box, x0, gnx[i], gnx[k], index[i], index[k],
459                                   bGroup, &dmin, &dmax, &nmin, &nmax, &min1, &min2, &max1, &max2);
460                         fprintf(dist, "  %12e", bMin ? dmin : dmax);
461                         if (num)
462                         {
463                             fprintf(num, "  %8d", bMin ? nmin : nmax);
464                         }
465                     }
466                 }
467             }
468         }
469         else
470         {
471             for (i = 1; (i < ng); i++)
472             {
473                 calc_dist(rcut, bPBC, ePBC, box, x0, gnx[0], gnx[i], index[0], index[i], bGroup,
474                           &dmin, &dmax, &nmin, &nmax, &min1, &min2, &max1, &max2);
475                 fprintf(dist, "  %12e", bMin ? dmin : dmax);
476                 if (num)
477                 {
478                     fprintf(num, "  %8d", bMin ? nmin : nmax);
479                 }
480                 if (nres)
481                 {
482                     for (j = 0; j < nres; j++)
483                     {
484                         calc_dist(rcut, bPBC, ePBC, box, x0, residue[j+1]-residue[j], gnx[i],
485                                   &(index[0][residue[j]]), index[i], bGroup,
486                                   &dmin, &dmax, &nmin, &nmax, &min1r, &min2r, &max1r, &max2r);
487                         mindres[i-1][j] = min(mindres[i-1][j], dmin);
488                         maxdres[i-1][j] = max(maxdres[i-1][j], dmax);
489                     }
490                 }
491             }
492         }
493         fprintf(dist, "\n");
494         if (num)
495         {
496             fprintf(num, "\n");
497         }
498         if ( (bMin ? min1 : max1) != -1)
499         {
500             if (atm)
501             {
502                 fprintf(atm, "%12e  %12d  %12d\n",
503                         output_env_conv_time(oenv, t), 1+(bMin ? min1 : max1),
504                         1+(bMin ? min2 : max2));
505             }
506         }
507
508         if (trxout)
509         {
510             oindex[0] = bMin ? min1 : max1;
511             oindex[1] = bMin ? min2 : max2;
512             write_trx(trxout, 2, oindex, atoms, i, t, box, x0, NULL, NULL);
513         }
514         bFirst = FALSE;
515         /*dmin should be minimum distance for residue and group*/
516         if (bEachResEachTime)
517         {
518             fprintf(respertime, "%12e", t);
519             for (i = 1; i < ng; i++)
520             {
521                 for (j = 0; j < nres; j++)
522                 {
523                     fprintf(respertime, " %7g", bMin ? mindres[i-1][j] : maxdres[i-1][j]);
524                     /*reset distances for next time point*/
525                     mindres[i-1][j] = 1e6;
526                     maxdres[i-1][j] = 0;
527                 }
528             }
529             fprintf(respertime, "\n");
530         }
531     }
532     while (read_next_x(oenv, status, &t, x0, box));
533
534     close_trj(status);
535     ffclose(dist);
536     if (num)
537     {
538         ffclose(num);
539     }
540     if (atm)
541     {
542         ffclose(atm);
543     }
544     if (trxout)
545     {
546         close_trx(trxout);
547     }
548
549     if (nres && !bEachResEachTime)
550     {
551         FILE *res;
552
553         sprintf(buf, "%simum Distance", bMin ? "Min" : "Max");
554         res = xvgropen(rfile, buf, "Residue (#)", "Distance (nm)", oenv);
555         xvgr_legend(res, ng-1, (const char**)leg, oenv);
556         for (j = 0; j < nres; j++)
557         {
558             fprintf(res, "%4d", j+1);
559             for (i = 1; i < ng; i++)
560             {
561                 fprintf(res, " %7g", bMin ? mindres[i-1][j] : maxdres[i-1][j]);
562             }
563             fprintf(res, "\n");
564         }
565     }
566
567     sfree(x0);
568 }
569
570 int find_residues(t_atoms *atoms, int n, atom_id index[], atom_id **resindex)
571 {
572     int  i;
573     int  nres = 0, resnr, presnr;
574     int *residx;
575
576     /* build index of first atom numbers for each residue */
577     presnr = NOTSET;
578     snew(residx, atoms->nres+1);
579     for (i = 0; i < n; i++)
580     {
581         resnr = atoms->atom[index[i]].resind;
582         if (resnr != presnr)
583         {
584             residx[nres] = i;
585             nres++;
586             presnr = resnr;
587         }
588     }
589     if (debug)
590     {
591         printf("Found %d residues out of %d (%d/%d atoms)\n",
592                nres, atoms->nres, atoms->nr, n);
593     }
594     srenew(residx, nres+1);
595     /* mark end of last residue */
596     residx[nres] = n;
597     *resindex    = residx;
598     return nres;
599 }
600
601 void dump_res(FILE *out, int nres, atom_id *resindex, atom_id index[])
602 {
603     int i, j;
604
605     for (i = 0; i < nres-1; i++)
606     {
607         fprintf(out, "Res %d (%d):", i, resindex[i+1]-resindex[i]);
608         for (j = resindex[i]; j < resindex[i+1]; j++)
609         {
610             fprintf(out, " %d(%d)", j, index[j]);
611         }
612         fprintf(out, "\n");
613     }
614 }
615
616 int gmx_mindist(int argc, char *argv[])
617 {
618     const char     *desc[] = {
619         "[TT]g_mindist[tt] computes the distance between one group and a number of",
620         "other groups. Both the minimum distance",
621         "(between any pair of atoms from the respective groups)",
622         "and the number of contacts within a given",
623         "distance are written to two separate output files.",
624         "With the [TT]-group[tt] option a contact of an atom in another group",
625         "with multiple atoms in the first group is counted as one contact",
626         "instead of as multiple contacts.",
627         "With [TT]-or[tt], minimum distances to each residue in the first",
628         "group are determined and plotted as a function of residue number.[PAR]",
629         "With option [TT]-pi[tt] the minimum distance of a group to its",
630         "periodic image is plotted. This is useful for checking if a protein",
631         "has seen its periodic image during a simulation. Only one shift in",
632         "each direction is considered, giving a total of 26 shifts.",
633         "It also plots the maximum distance within the group and the lengths",
634         "of the three box vectors.[PAR]",
635         "Other programs that calculate distances are [TT]g_dist[tt]",
636         "and [TT]g_bond[tt]."
637     };
638     const char     *bugs[] = {
639         "The [TT]-pi[tt] option is very slow."
640     };
641
642     static gmx_bool bMat             = FALSE, bPI = FALSE, bSplit = FALSE, bMax = FALSE, bPBC = TRUE;
643     static gmx_bool bGroup           = FALSE;
644     static real     rcutoff          = 0.6;
645     static int      ng               = 1;
646     static gmx_bool bEachResEachTime = FALSE, bPrintResName = FALSE;
647     t_pargs         pa[]             = {
648         { "-matrix", FALSE, etBOOL, {&bMat},
649           "Calculate half a matrix of group-group distances" },
650         { "-max",    FALSE, etBOOL, {&bMax},
651           "Calculate *maximum* distance instead of minimum" },
652         { "-d",      FALSE, etREAL, {&rcutoff},
653           "Distance for contacts" },
654         { "-group",      FALSE, etBOOL, {&bGroup},
655           "Count contacts with multiple atoms in the first group as one" },
656         { "-pi",     FALSE, etBOOL, {&bPI},
657           "Calculate minimum distance with periodic images" },
658         { "-split",  FALSE, etBOOL, {&bSplit},
659           "Split graph where time is zero" },
660         { "-ng",       FALSE, etINT, {&ng},
661           "Number of secondary groups to compute distance to a central group" },
662         { "-pbc",    FALSE, etBOOL, {&bPBC},
663           "Take periodic boundary conditions into account" },
664         { "-respertime",  FALSE, etBOOL, {&bEachResEachTime},
665           "When writing per-residue distances, write distance for each time point" },
666         { "-printresname",  FALSE, etBOOL, {&bPrintResName},
667           "Write residue names" }
668     };
669     output_env_t    oenv;
670     t_topology     *top  = NULL;
671     int             ePBC = -1;
672     char            title[256];
673     real            t;
674     rvec           *x;
675     matrix          box;
676     gmx_bool        bTop = FALSE;
677
678     FILE           *atm;
679     int             i, j, nres = 0;
680     const char     *trxfnm, *tpsfnm, *ndxfnm, *distfnm, *numfnm, *atmfnm, *oxfnm, *resfnm;
681     char          **grpname;
682     int            *gnx;
683     atom_id       **index, *residues = NULL;
684     t_filenm        fnm[] = {
685         { efTRX, "-f",  NULL,      ffREAD },
686         { efTPS,  NULL, NULL,      ffOPTRD },
687         { efNDX,  NULL, NULL,      ffOPTRD },
688         { efXVG, "-od", "mindist",  ffWRITE },
689         { efXVG, "-on", "numcont",  ffOPTWR },
690         { efOUT, "-o", "atm-pair", ffOPTWR },
691         { efTRO, "-ox", "mindist",  ffOPTWR },
692         { efXVG, "-or", "mindistres", ffOPTWR }
693     };
694 #define NFILE asize(fnm)
695
696     parse_common_args(&argc, argv,
697                       PCA_CAN_VIEW | PCA_CAN_TIME | PCA_TIME_UNIT | PCA_BE_NICE,
698                       NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL, &oenv);
699
700     trxfnm  = ftp2fn(efTRX, NFILE, fnm);
701     ndxfnm  = ftp2fn_null(efNDX, NFILE, fnm);
702     distfnm = opt2fn("-od", NFILE, fnm);
703     numfnm  = opt2fn_null("-on", NFILE, fnm);
704     atmfnm  = ftp2fn_null(efOUT, NFILE, fnm);
705     oxfnm   = opt2fn_null("-ox", NFILE, fnm);
706     resfnm  = opt2fn_null("-or", NFILE, fnm);
707     if (bPI || resfnm != NULL)
708     {
709         /* We need a tps file */
710         tpsfnm = ftp2fn(efTPS, NFILE, fnm);
711     }
712     else
713     {
714         tpsfnm = ftp2fn_null(efTPS, NFILE, fnm);
715     }
716
717     if (!tpsfnm && !ndxfnm)
718     {
719         gmx_fatal(FARGS, "You have to specify either the index file or a tpr file");
720     }
721
722     if (bPI)
723     {
724         ng = 1;
725         fprintf(stderr, "Choose a group for distance calculation\n");
726     }
727     else if (!bMat)
728     {
729         ng++;
730     }
731
732     snew(gnx, ng);
733     snew(index, ng);
734     snew(grpname, ng);
735
736     if (tpsfnm || resfnm || !ndxfnm)
737     {
738         snew(top, 1);
739         bTop = read_tps_conf(tpsfnm, title, top, &ePBC, &x, NULL, box, FALSE);
740         if (bPI && !bTop)
741         {
742             printf("\nWARNING: Without a run input file a trajectory with broken molecules will not give the correct periodic image distance\n\n");
743         }
744     }
745     get_index(top ? &(top->atoms) : NULL, ndxfnm, ng, gnx, index, grpname);
746
747     if (bMat && (ng == 1))
748     {
749         ng = gnx[0];
750         printf("Special case: making distance matrix between all atoms in group %s\n",
751                grpname[0]);
752         srenew(gnx, ng);
753         srenew(index, ng);
754         srenew(grpname, ng);
755         for (i = 1; (i < ng); i++)
756         {
757             gnx[i]      = 1;
758             grpname[i]  = grpname[0];
759             snew(index[i], 1);
760             index[i][0] = index[0][i];
761         }
762         gnx[0] = 1;
763     }
764
765     if (resfnm)
766     {
767         nres = find_residues(top ? &(top->atoms) : NULL,
768                              gnx[0], index[0], &residues);
769         if (debug)
770         {
771             dump_res(debug, nres, residues, index[0]);
772         }
773     }
774
775     if (bPI)
776     {
777         periodic_mindist_plot(trxfnm, distfnm, top, ePBC, gnx[0], index[0], bSplit, oenv);
778     }
779     else
780     {
781         dist_plot(trxfnm, atmfnm, distfnm, numfnm, resfnm, oxfnm,
782                   rcutoff, bMat, top ? &(top->atoms) : NULL,
783                   ng, index, gnx, grpname, bSplit, !bMax, nres, residues, bPBC, ePBC,
784                   bGroup, bEachResEachTime, bPrintResName, oenv);
785     }
786
787     do_view(oenv, distfnm, "-nxy");
788     if (!bPI)
789     {
790         do_view(oenv, numfnm, "-nxy");
791     }
792
793     return 0;
794 }