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