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