9ce1732513998f81d27c582046e1bf5767c00a1c
[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,2017,2018, 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, const 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     }
90
91     nshift = 0;
92     for (sz = -nsz; sz <= nsz; sz++)
93     {
94         for (sy = -1; sy <= 1; sy++)
95         {
96             for (sx = -1; sx <= 1; sx++)
97             {
98                 if (sx != 0 || sy != 0 || sz != 0)
99                 {
100                     for (i = 0; i < DIM; i++)
101                     {
102                         shift[nshift][i] =
103                             sx*box[XX][i] + sy*box[YY][i] + sz*box[ZZ][i];
104                     }
105                     nshift++;
106                 }
107             }
108         }
109     }
110
111     r2min = sqr_box;
112     r2max = 0;
113
114     for (i = 0; i < n; i++)
115     {
116         for (j = i+1; j < n; j++)
117         {
118             rvec_sub(x[index[i]], x[index[j]], d0);
119             r2 = norm2(d0);
120             if (r2 > r2max)
121             {
122                 r2max = r2;
123             }
124             for (s = 0; s < nshift; s++)
125             {
126                 rvec_add(d0, shift[s], d);
127                 r2 = norm2(d);
128                 if (r2 < r2min)
129                 {
130                     r2min      = r2;
131                     min_ind[0] = i;
132                     min_ind[1] = j;
133                 }
134             }
135         }
136     }
137
138     *rmin = std::sqrt(r2min);
139     *rmax = std::sqrt(r2max);
140 }
141
142 static void periodic_mindist_plot(const char *trxfn, const char *outfn,
143                                   const t_topology *top, int ePBC,
144                                   int n, int index[], gmx_bool bSplit,
145                                   const gmx_output_env_t *oenv)
146 {
147     FILE        *out;
148     const char  *leg[5] = { "min per.", "max int.", "box1", "box2", "box3" };
149     t_trxstatus *status;
150     real         t;
151     rvec        *x;
152     matrix       box;
153     int          natoms, ind_min[2] = {0, 0}, ind_mini = 0, ind_minj = 0;
154     real         rmin, rmax, rmint, tmint;
155     gmx_bool     bFirst;
156     gmx_rmpbc_t  gpbc = nullptr;
157
158     natoms = read_first_x(oenv, &status, trxfn, &t, &x, box);
159
160     check_index(nullptr, n, index, nullptr, natoms);
161
162     out = xvgropen(outfn, "Minimum distance to periodic image",
163                    output_env_get_time_label(oenv), "Distance (nm)", oenv);
164     if (output_env_get_print_xvgr_codes(oenv))
165     {
166         fprintf(out, "@ subtitle \"and maximum internal distance\"\n");
167     }
168     xvgr_legend(out, 5, leg, oenv);
169
170     rmint = box[XX][XX];
171     tmint = 0;
172
173     if (nullptr != top)
174     {
175         gpbc = gmx_rmpbc_init(&top->idef, ePBC, natoms);
176     }
177
178     bFirst = TRUE;
179     do
180     {
181         if (nullptr != top)
182         {
183             gmx_rmpbc(gpbc, natoms, box, x);
184         }
185
186         periodic_dist(ePBC, box, x, n, index, &rmin, &rmax, ind_min);
187         if (rmin < rmint)
188         {
189             rmint    = rmin;
190             tmint    = t;
191             ind_mini = ind_min[0];
192             ind_minj = ind_min[1];
193         }
194         if (bSplit && !bFirst && std::abs(t/output_env_get_time_factor(oenv)) < 1e-5)
195         {
196             fprintf(out, "%s\n", output_env_get_print_xvgr_codes(oenv) ? "&" : "");
197         }
198         fprintf(out, "\t%g\t%6.3f %6.3f %6.3f %6.3f %6.3f\n",
199                 output_env_conv_time(oenv, t), rmin, rmax, norm(box[0]), norm(box[1]), norm(box[2]));
200         bFirst = FALSE;
201     }
202     while (read_next_x(oenv, status, &t, x, box));
203
204     if (nullptr != top)
205     {
206         gmx_rmpbc_done(gpbc);
207     }
208
209     xvgrclose(out);
210
211     fprintf(stdout,
212             "\nThe shortest periodic distance is %g (nm) at time %g (%s),\n"
213             "between atoms %d and %d\n",
214             rmint, output_env_conv_time(oenv, tmint), output_env_get_time_unit(oenv).c_str(),
215             index[ind_mini]+1, index[ind_minj]+1);
216 }
217
218 static void calc_dist(real rcut, gmx_bool bPBC, int ePBC, matrix box, rvec x[],
219                       int nx1, int nx2, int index1[], int index2[],
220                       gmx_bool bGroup,
221                       real *rmin, real *rmax, int *nmin, int *nmax,
222                       int *ixmin, int *jxmin, int *ixmax, int *jxmax)
223 {
224     int      i, j, i0 = 0, j1;
225     int      ix, jx;
226     int     *index3;
227     rvec     dx;
228     real     r2, rmin2, rmax2, rcut2;
229     t_pbc    pbc;
230     int      nmin_j, nmax_j;
231
232     *ixmin = -1;
233     *jxmin = -1;
234     *ixmax = -1;
235     *jxmax = -1;
236     *nmin  = 0;
237     *nmax  = 0;
238
239     rcut2 = gmx::square(rcut);
240
241     /* Must init pbc every step because of pressure coupling */
242     if (bPBC)
243     {
244         set_pbc(&pbc, ePBC, box);
245     }
246     if (index2)
247     {
248         i0     = 0;
249         j1     = nx2;
250         index3 = index2;
251     }
252     else
253     {
254         j1     = nx1;
255         index3 = index1;
256     }
257
258     rmin2 = 1e12;
259     rmax2 = -1e12;
260
261     for (j = 0; (j < j1); j++)
262     {
263         jx = index3[j];
264         if (index2 == nullptr)
265         {
266             i0 = j + 1;
267         }
268         nmin_j = 0;
269         nmax_j = 0;
270         for (i = i0; (i < nx1); i++)
271         {
272             ix = index1[i];
273             if (ix != jx)
274             {
275                 if (bPBC)
276                 {
277                     pbc_dx(&pbc, x[ix], x[jx], dx);
278                 }
279                 else
280                 {
281                     rvec_sub(x[ix], x[jx], dx);
282                 }
283                 r2 = iprod(dx, dx);
284                 if (r2 < rmin2)
285                 {
286                     rmin2  = r2;
287                     *ixmin = ix;
288                     *jxmin = jx;
289                 }
290                 if (r2 > rmax2)
291                 {
292                     rmax2  = r2;
293                     *ixmax = ix;
294                     *jxmax = jx;
295                 }
296                 if (r2 <= rcut2)
297                 {
298                     nmin_j++;
299                 }
300                 else
301                 {
302                     nmax_j++;
303                 }
304             }
305         }
306         if (bGroup)
307         {
308             if (nmin_j > 0)
309             {
310                 (*nmin)++;
311             }
312             if (nmax_j > 0)
313             {
314                 (*nmax)++;
315             }
316         }
317         else
318         {
319             *nmin += nmin_j;
320             *nmax += nmax_j;
321         }
322     }
323     *rmin = std::sqrt(rmin2);
324     *rmax = std::sqrt(rmax2);
325 }
326
327 static void dist_plot(const char *fn, const char *afile, const char *dfile,
328                       const char *nfile, const char *rfile, const char *xfile,
329                       real rcut, gmx_bool bMat, const t_atoms *atoms,
330                       int ng, int *index[], int gnx[], char *grpn[], gmx_bool bSplit,
331                       gmx_bool bMin, int nres, int *residue, gmx_bool bPBC, int ePBC,
332                       gmx_bool bGroup, gmx_bool bEachResEachTime, gmx_bool bPrintResName,
333                       const gmx_output_env_t *oenv)
334 {
335     FILE            *atm, *dist, *num;
336     t_trxstatus     *trxout;
337     char             buf[256];
338     char           **leg;
339     real             t, dmin, dmax, **mindres = nullptr, **maxdres = nullptr;
340     int              nmin, nmax;
341     t_trxstatus     *status;
342     int              i = -1, j, k;
343     int              min2, max2, min1r, min2r, max1r, max2r;
344     int              min1 = 0;
345     int              max1 = 0;
346     int              oindex[2];
347     rvec            *x0;
348     matrix           box;
349     gmx_bool         bFirst;
350     FILE            *respertime = nullptr;
351
352     if (read_first_x(oenv, &status, fn, &t, &x0, box) == 0)
353     {
354         gmx_fatal(FARGS, "Could not read coordinates from statusfile\n");
355     }
356
357     sprintf(buf, "%simum Distance", bMin ? "Min" : "Max");
358     dist = xvgropen(dfile, buf, output_env_get_time_label(oenv), "Distance (nm)", oenv);
359     sprintf(buf, "Number of Contacts %s %g nm", bMin ? "<" : ">", rcut);
360     num    = nfile ? xvgropen(nfile, buf, output_env_get_time_label(oenv), "Number", oenv) : nullptr;
361     atm    = afile ? gmx_ffopen(afile, "w") : nullptr;
362     trxout = xfile ? open_trx(xfile, "w") : nullptr;
363
364     if (bMat)
365     {
366         if (ng == 1)
367         {
368             snew(leg, 1);
369             sprintf(buf, "Internal in %s", grpn[0]);
370             leg[0] = gmx_strdup(buf);
371             xvgr_legend(dist, 0, leg, oenv);
372             if (num)
373             {
374                 xvgr_legend(num, 0, leg, oenv);
375             }
376         }
377         else
378         {
379             snew(leg, (ng*(ng-1))/2);
380             for (i = j = 0; (i < ng-1); i++)
381             {
382                 for (k = i+1; (k < ng); k++, j++)
383                 {
384                     sprintf(buf, "%s-%s", grpn[i], grpn[k]);
385                     leg[j] = gmx_strdup(buf);
386                 }
387             }
388             xvgr_legend(dist, j, leg, oenv);
389             if (num)
390             {
391                 xvgr_legend(num, j, leg, oenv);
392             }
393         }
394     }
395     else
396     {
397         snew(leg, ng-1);
398         for (i = 0; (i < ng-1); i++)
399         {
400             sprintf(buf, "%s-%s", grpn[0], grpn[i+1]);
401             leg[i] = gmx_strdup(buf);
402         }
403         xvgr_legend(dist, ng-1, leg, oenv);
404         if (num)
405         {
406             xvgr_legend(num, ng-1, leg, oenv);
407         }
408     }
409
410     if (bEachResEachTime)
411     {
412         sprintf(buf, "%simum Distance", bMin ? "Min" : "Max");
413         respertime = xvgropen(rfile, buf, output_env_get_time_label(oenv), "Distance (nm)", oenv);
414         xvgr_legend(respertime, ng-1, leg, oenv);
415         if (bPrintResName && output_env_get_print_xvgr_codes(oenv) )
416         {
417             fprintf(respertime, "# ");
418
419             for (j = 0; j < nres; j++)
420             {
421                 fprintf(respertime, "%s%d ", *(atoms->resinfo[atoms->atom[index[0][residue[j]]].resind].name), atoms->atom[index[0][residue[j]]].resind);
422             }
423             fprintf(respertime, "\n");
424         }
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, nullptr, nullptr);
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_trx(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, 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 static int find_residues(const t_atoms *atoms, int n, const 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 static 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. Note",
663         "that periodicity information is required from the file supplied with",
664         "with [TT]-s[tt], either as a .tpr file or a .pdb file with CRYST1 fields.",
665         "It also plots the maximum distance within the group and the lengths",
666         "of the three box vectors.[PAR]",
667         "Also [gmx-distance] and [gmx-pairdist] calculate distances."
668     };
669
670     static gmx_bool   bMat             = FALSE, bPI = FALSE, bSplit = FALSE, bMax = FALSE, bPBC = TRUE;
671     static gmx_bool   bGroup           = FALSE;
672     static real       rcutoff          = 0.6;
673     static int        ng               = 1;
674     static gmx_bool   bEachResEachTime = FALSE, bPrintResName = FALSE;
675     t_pargs           pa[]             = {
676         { "-matrix", FALSE, etBOOL, {&bMat},
677           "Calculate half a matrix of group-group distances" },
678         { "-max",    FALSE, etBOOL, {&bMax},
679           "Calculate *maximum* distance instead of minimum" },
680         { "-d",      FALSE, etREAL, {&rcutoff},
681           "Distance for contacts" },
682         { "-group",      FALSE, etBOOL, {&bGroup},
683           "Count contacts with multiple atoms in the first group as one" },
684         { "-pi",     FALSE, etBOOL, {&bPI},
685           "Calculate minimum distance with periodic images" },
686         { "-split",  FALSE, etBOOL, {&bSplit},
687           "Split graph where time is zero" },
688         { "-ng",       FALSE, etINT, {&ng},
689           "Number of secondary groups to compute distance to a central group" },
690         { "-pbc",    FALSE, etBOOL, {&bPBC},
691           "Take periodic boundary conditions into account" },
692         { "-respertime",  FALSE, etBOOL, {&bEachResEachTime},
693           "When writing per-residue distances, write distance for each time point" },
694         { "-printresname",  FALSE, etBOOL, {&bPrintResName},
695           "Write residue names" }
696     };
697     gmx_output_env_t *oenv;
698     t_topology       *top  = nullptr;
699     int               ePBC = -1;
700     rvec             *x;
701     matrix            box;
702     gmx_bool          bTop = FALSE;
703
704     int               i, nres = 0;
705     const char       *trxfnm, *tpsfnm, *ndxfnm, *distfnm, *numfnm, *atmfnm, *oxfnm, *resfnm;
706     char            **grpname;
707     int              *gnx;
708     int             **index, *residues = nullptr;
709     t_filenm          fnm[] = {
710         { efTRX, "-f",  nullptr,      ffREAD },
711         { efTPS,  nullptr, nullptr,      ffOPTRD },
712         { efNDX,  nullptr, nullptr,      ffOPTRD },
713         { efXVG, "-od", "mindist",  ffWRITE },
714         { efXVG, "-on", "numcont",  ffOPTWR },
715         { efOUT, "-o", "atm-pair", ffOPTWR },
716         { efTRO, "-ox", "mindist",  ffOPTWR },
717         { efXVG, "-or", "mindistres", ffOPTWR }
718     };
719 #define NFILE asize(fnm)
720
721     if (!parse_common_args(&argc, argv,
722                            PCA_CAN_VIEW | PCA_CAN_TIME | PCA_TIME_UNIT,
723                            NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, nullptr, &oenv))
724     {
725         return 0;
726     }
727
728     trxfnm  = ftp2fn(efTRX, NFILE, fnm);
729     ndxfnm  = ftp2fn_null(efNDX, NFILE, fnm);
730     distfnm = opt2fn("-od", NFILE, fnm);
731     numfnm  = opt2fn_null("-on", NFILE, fnm);
732     atmfnm  = ftp2fn_null(efOUT, NFILE, fnm);
733     oxfnm   = opt2fn_null("-ox", NFILE, fnm);
734     resfnm  = opt2fn_null("-or", NFILE, fnm);
735     if (bPI || resfnm != nullptr)
736     {
737         /* We need a tps file */
738         tpsfnm = ftp2fn(efTPS, NFILE, fnm);
739     }
740     else
741     {
742         tpsfnm = ftp2fn_null(efTPS, NFILE, fnm);
743     }
744
745     if (!tpsfnm && !ndxfnm)
746     {
747         gmx_fatal(FARGS, "You have to specify either the index file or a tpr file");
748     }
749
750     if (bPI)
751     {
752         ng = 1;
753         fprintf(stderr, "Choose a group for distance calculation\n");
754     }
755     else if (!bMat)
756     {
757         ng++;
758     }
759
760     snew(gnx, ng);
761     snew(index, ng);
762     snew(grpname, ng);
763
764     if (tpsfnm || resfnm || !ndxfnm)
765     {
766         snew(top, 1);
767         bTop = read_tps_conf(tpsfnm, top, &ePBC, &x, nullptr, box, FALSE);
768         if (bPI && !bTop)
769         {
770             printf("\nWARNING: Without a run input file a trajectory with broken molecules will not give the correct periodic image distance\n\n");
771         }
772     }
773     get_index(top ? &(top->atoms) : nullptr, ndxfnm, ng, gnx, index, grpname);
774
775     if (bMat && (ng == 1))
776     {
777         ng = gnx[0];
778         printf("Special case: making distance matrix between all atoms in group %s\n",
779                grpname[0]);
780         srenew(gnx, ng);
781         srenew(index, ng);
782         srenew(grpname, ng);
783         for (i = 1; (i < ng); i++)
784         {
785             gnx[i]      = 1;
786             grpname[i]  = grpname[0];
787             snew(index[i], 1);
788             index[i][0] = index[0][i];
789         }
790         gnx[0] = 1;
791     }
792
793     if (resfnm)
794     {
795         GMX_RELEASE_ASSERT(top != nullptr, "top pointer cannot be NULL when finding residues");
796         nres = find_residues(&(top->atoms), gnx[0], index[0], &residues);
797
798         if (debug)
799         {
800             dump_res(debug, nres, residues, index[0]);
801         }
802     }
803     else if (bEachResEachTime || bPrintResName)
804     {
805         gmx_fatal(FARGS, "Option -or needs to be set to print residues");
806     }
807
808     if (bPI)
809     {
810         periodic_mindist_plot(trxfnm, distfnm, top, ePBC, gnx[0], index[0], bSplit, oenv);
811     }
812     else
813     {
814         dist_plot(trxfnm, atmfnm, distfnm, numfnm, resfnm, oxfnm,
815                   rcutoff, bMat, top ? &(top->atoms) : nullptr,
816                   ng, index, gnx, grpname, bSplit, !bMax, nres, residues, bPBC, ePBC,
817                   bGroup, bEachResEachTime, bPrintResName, oenv);
818     }
819
820     do_view(oenv, distfnm, "-nxy");
821     if (!bPI)
822     {
823         do_view(oenv, numfnm, "-nxy");
824     }
825
826     return 0;
827 }