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