202eda2b1084569ee8361f9808b3926044cdff73
[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, 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 = nullptr;
158
159     natoms = read_first_x(oenv, &status, trxfn, &t, &x, box);
160
161     check_index(nullptr, n, index, nullptr, 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 (nullptr != top)
175     {
176         gpbc = gmx_rmpbc_init(&top->idef, ePBC, natoms);
177     }
178
179     bFirst = TRUE;
180     do
181     {
182         if (nullptr != 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 (nullptr != 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).c_str(),
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 == nullptr)
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 = nullptr, **maxdres = nullptr;
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 = nullptr;
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) : nullptr;
362     atm    = afile ? gmx_ffopen(afile, "w") : nullptr;
363     trxout = xfile ? open_trx(xfile, "w") : nullptr;
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 && output_env_get_print_xvgr_codes(oenv) )
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
429     if (nres)
430     {
431         snew(mindres, ng-1);
432         snew(maxdres, ng-1);
433         for (i = 1; i < ng; i++)
434         {
435             snew(mindres[i-1], nres);
436             snew(maxdres[i-1], nres);
437             for (j = 0; j < nres; j++)
438             {
439                 mindres[i-1][j] = 1e6;
440             }
441             /* maxdres[*][*] is already 0 */
442         }
443     }
444     bFirst = TRUE;
445     do
446     {
447         if (bSplit && !bFirst && std::abs(t/output_env_get_time_factor(oenv)) < 1e-5)
448         {
449             fprintf(dist, "%s\n", output_env_get_print_xvgr_codes(oenv) ? "&" : "");
450             if (num)
451             {
452                 fprintf(num, "%s\n", output_env_get_print_xvgr_codes(oenv) ? "&" : "");
453             }
454             if (atm)
455             {
456                 fprintf(atm, "%s\n", output_env_get_print_xvgr_codes(oenv) ? "&" : "");
457             }
458         }
459         fprintf(dist, "%12e", output_env_conv_time(oenv, t));
460         if (num)
461         {
462             fprintf(num, "%12e", output_env_conv_time(oenv, t));
463         }
464
465         if (bMat)
466         {
467             if (ng == 1)
468             {
469                 calc_dist(rcut, bPBC, ePBC, box, x0, gnx[0], gnx[0], index[0], index[0], bGroup,
470                           &dmin, &dmax, &nmin, &nmax, &min1, &min2, &max1, &max2);
471                 fprintf(dist, "  %12e", bMin ? dmin : dmax);
472                 if (num)
473                 {
474                     fprintf(num, "  %8d", bMin ? nmin : nmax);
475                 }
476             }
477             else
478             {
479                 for (i = 0; (i < ng-1); i++)
480                 {
481                     for (k = i+1; (k < ng); k++)
482                     {
483                         calc_dist(rcut, bPBC, ePBC, box, x0, gnx[i], gnx[k], index[i], index[k],
484                                   bGroup, &dmin, &dmax, &nmin, &nmax, &min1, &min2, &max1, &max2);
485                         fprintf(dist, "  %12e", bMin ? dmin : dmax);
486                         if (num)
487                         {
488                             fprintf(num, "  %8d", bMin ? nmin : nmax);
489                         }
490                     }
491                 }
492             }
493         }
494         else
495         {
496             for (i = 1; (i < ng); i++)
497             {
498                 calc_dist(rcut, bPBC, ePBC, box, x0, gnx[0], gnx[i], index[0], index[i], bGroup,
499                           &dmin, &dmax, &nmin, &nmax, &min1, &min2, &max1, &max2);
500                 fprintf(dist, "  %12e", bMin ? dmin : dmax);
501                 if (num)
502                 {
503                     fprintf(num, "  %8d", bMin ? nmin : nmax);
504                 }
505                 if (nres)
506                 {
507                     for (j = 0; j < nres; j++)
508                     {
509                         calc_dist(rcut, bPBC, ePBC, box, x0, residue[j+1]-residue[j], gnx[i],
510                                   &(index[0][residue[j]]), index[i], bGroup,
511                                   &dmin, &dmax, &nmin, &nmax, &min1r, &min2r, &max1r, &max2r);
512                         mindres[i-1][j] = std::min(mindres[i-1][j], dmin);
513                         maxdres[i-1][j] = std::max(maxdres[i-1][j], dmax);
514                     }
515                 }
516             }
517         }
518         fprintf(dist, "\n");
519         if (num)
520         {
521             fprintf(num, "\n");
522         }
523         if ( (bMin ? min1 : max1) != -1)
524         {
525             if (atm)
526             {
527                 fprintf(atm, "%12e  %12d  %12d\n",
528                         output_env_conv_time(oenv, t), 1+(bMin ? min1 : max1),
529                         1+(bMin ? min2 : max2));
530             }
531         }
532
533         if (trxout)
534         {
535             oindex[0] = bMin ? min1 : max1;
536             oindex[1] = bMin ? min2 : max2;
537             write_trx(trxout, 2, oindex, atoms, i, t, box, x0, nullptr, nullptr);
538         }
539         bFirst = FALSE;
540         /*dmin should be minimum distance for residue and group*/
541         if (bEachResEachTime)
542         {
543             fprintf(respertime, "%12e", t);
544             for (i = 1; i < ng; i++)
545             {
546                 for (j = 0; j < nres; j++)
547                 {
548                     fprintf(respertime, " %7g", bMin ? mindres[i-1][j] : maxdres[i-1][j]);
549                     /*reset distances for next time point*/
550                     mindres[i-1][j] = 1e6;
551                     maxdres[i-1][j] = 0;
552                 }
553             }
554             fprintf(respertime, "\n");
555         }
556     }
557     while (read_next_x(oenv, status, &t, x0, box));
558
559     close_trx(status);
560     xvgrclose(dist);
561     if (num)
562     {
563         xvgrclose(num);
564     }
565     if (atm)
566     {
567         gmx_ffclose(atm);
568     }
569     if (trxout)
570     {
571         close_trx(trxout);
572     }
573     if (respertime)
574     {
575         xvgrclose(respertime);
576     }
577
578     if (nres && !bEachResEachTime)
579     {
580         FILE *res;
581
582         sprintf(buf, "%simum Distance", bMin ? "Min" : "Max");
583         res = xvgropen(rfile, buf, "Residue (#)", "Distance (nm)", oenv);
584         xvgr_legend(res, ng-1, (const char**)leg, oenv);
585         for (j = 0; j < nres; j++)
586         {
587             fprintf(res, "%4d", j+1);
588             for (i = 1; i < ng; i++)
589             {
590                 fprintf(res, " %7g", bMin ? mindres[i-1][j] : maxdres[i-1][j]);
591             }
592             fprintf(res, "\n");
593         }
594         xvgrclose(res);
595     }
596
597     sfree(x0);
598 }
599
600 int find_residues(const t_atoms *atoms, int n, int index[], int **resindex)
601 {
602     int  i;
603     int  nres      = 0, resnr, presnr = 0;
604     bool presFound = false;
605     int *residx;
606
607     /* build index of first atom numbers for each residue */
608     snew(residx, atoms->nres+1);
609     for (i = 0; i < n; i++)
610     {
611         resnr = atoms->atom[index[i]].resind;
612         if (!presFound || resnr != presnr)
613         {
614             residx[nres] = i;
615             nres++;
616             presnr    = resnr;
617             presFound = true;
618         }
619     }
620     if (debug)
621     {
622         printf("Found %d residues out of %d (%d/%d atoms)\n",
623                nres, atoms->nres, atoms->nr, n);
624     }
625     srenew(residx, nres+1);
626     /* mark end of last residue */
627     residx[nres] = n;
628     *resindex    = residx;
629     return nres;
630 }
631
632 void dump_res(FILE *out, int nres, int *resindex, int index[])
633 {
634     int i, j;
635
636     for (i = 0; i < nres-1; i++)
637     {
638         fprintf(out, "Res %d (%d):", i, resindex[i+1]-resindex[i]);
639         for (j = resindex[i]; j < resindex[i+1]; j++)
640         {
641             fprintf(out, " %d(%d)", j, index[j]);
642         }
643         fprintf(out, "\n");
644     }
645 }
646
647 int gmx_mindist(int argc, char *argv[])
648 {
649     const char       *desc[] = {
650         "[THISMODULE] computes the distance between one group and a number of",
651         "other groups. Both the minimum distance",
652         "(between any pair of atoms from the respective groups)",
653         "and the number of contacts within a given",
654         "distance are written to two separate output files.",
655         "With the [TT]-group[tt] option a contact of an atom in another group",
656         "with multiple atoms in the first group is counted as one contact",
657         "instead of as multiple contacts.",
658         "With [TT]-or[tt], minimum distances to each residue in the first",
659         "group are determined and plotted as a function of residue number.[PAR]",
660         "With option [TT]-pi[tt] the minimum distance of a group to its",
661         "periodic image is plotted. This is useful for checking if a protein",
662         "has seen its periodic image during a simulation. Only one shift in",
663         "each direction is considered, giving a total of 26 shifts.",
664         "It also plots the maximum distance within the group and the lengths",
665         "of the three box vectors.[PAR]",
666         "Also [gmx-distance] and [gmx-pairdist] calculate distances."
667     };
668
669     static gmx_bool   bMat             = FALSE, bPI = FALSE, bSplit = FALSE, bMax = FALSE, bPBC = TRUE;
670     static gmx_bool   bGroup           = FALSE;
671     static real       rcutoff          = 0.6;
672     static int        ng               = 1;
673     static gmx_bool   bEachResEachTime = FALSE, bPrintResName = FALSE;
674     t_pargs           pa[]             = {
675         { "-matrix", FALSE, etBOOL, {&bMat},
676           "Calculate half a matrix of group-group distances" },
677         { "-max",    FALSE, etBOOL, {&bMax},
678           "Calculate *maximum* distance instead of minimum" },
679         { "-d",      FALSE, etREAL, {&rcutoff},
680           "Distance for contacts" },
681         { "-group",      FALSE, etBOOL, {&bGroup},
682           "Count contacts with multiple atoms in the first group as one" },
683         { "-pi",     FALSE, etBOOL, {&bPI},
684           "Calculate minimum distance with periodic images" },
685         { "-split",  FALSE, etBOOL, {&bSplit},
686           "Split graph where time is zero" },
687         { "-ng",       FALSE, etINT, {&ng},
688           "Number of secondary groups to compute distance to a central group" },
689         { "-pbc",    FALSE, etBOOL, {&bPBC},
690           "Take periodic boundary conditions into account" },
691         { "-respertime",  FALSE, etBOOL, {&bEachResEachTime},
692           "When writing per-residue distances, write distance for each time point" },
693         { "-printresname",  FALSE, etBOOL, {&bPrintResName},
694           "Write residue names" }
695     };
696     gmx_output_env_t *oenv;
697     t_topology       *top  = nullptr;
698     int               ePBC = -1;
699     rvec             *x;
700     matrix            box;
701     gmx_bool          bTop = FALSE;
702
703     int               i, nres = 0;
704     const char       *trxfnm, *tpsfnm, *ndxfnm, *distfnm, *numfnm, *atmfnm, *oxfnm, *resfnm;
705     char            **grpname;
706     int              *gnx;
707     int             **index, *residues = nullptr;
708     t_filenm          fnm[] = {
709         { efTRX, "-f",  nullptr,      ffREAD },
710         { efTPS,  nullptr, nullptr,      ffOPTRD },
711         { efNDX,  nullptr, nullptr,      ffOPTRD },
712         { efXVG, "-od", "mindist",  ffWRITE },
713         { efXVG, "-on", "numcont",  ffOPTWR },
714         { efOUT, "-o", "atm-pair", ffOPTWR },
715         { efTRO, "-ox", "mindist",  ffOPTWR },
716         { efXVG, "-or", "mindistres", ffOPTWR }
717     };
718 #define NFILE asize(fnm)
719
720     if (!parse_common_args(&argc, argv,
721                            PCA_CAN_VIEW | PCA_CAN_TIME | PCA_TIME_UNIT,
722                            NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, nullptr, &oenv))
723     {
724         return 0;
725     }
726
727     trxfnm  = ftp2fn(efTRX, NFILE, fnm);
728     ndxfnm  = ftp2fn_null(efNDX, NFILE, fnm);
729     distfnm = opt2fn("-od", NFILE, fnm);
730     numfnm  = opt2fn_null("-on", NFILE, fnm);
731     atmfnm  = ftp2fn_null(efOUT, NFILE, fnm);
732     oxfnm   = opt2fn_null("-ox", NFILE, fnm);
733     resfnm  = opt2fn_null("-or", NFILE, fnm);
734     if (bPI || resfnm != nullptr)
735     {
736         /* We need a tps file */
737         tpsfnm = ftp2fn(efTPS, NFILE, fnm);
738     }
739     else
740     {
741         tpsfnm = ftp2fn_null(efTPS, NFILE, fnm);
742     }
743
744     if (!tpsfnm && !ndxfnm)
745     {
746         gmx_fatal(FARGS, "You have to specify either the index file or a tpr file");
747     }
748
749     if (bPI)
750     {
751         ng = 1;
752         fprintf(stderr, "Choose a group for distance calculation\n");
753     }
754     else if (!bMat)
755     {
756         ng++;
757     }
758
759     snew(gnx, ng);
760     snew(index, ng);
761     snew(grpname, ng);
762
763     if (tpsfnm || resfnm || !ndxfnm)
764     {
765         snew(top, 1);
766         bTop = read_tps_conf(tpsfnm, top, &ePBC, &x, nullptr, box, FALSE);
767         if (bPI && !bTop)
768         {
769             printf("\nWARNING: Without a run input file a trajectory with broken molecules will not give the correct periodic image distance\n\n");
770         }
771     }
772     get_index(top ? &(top->atoms) : nullptr, ndxfnm, ng, gnx, index, grpname);
773
774     if (bMat && (ng == 1))
775     {
776         ng = gnx[0];
777         printf("Special case: making distance matrix between all atoms in group %s\n",
778                grpname[0]);
779         srenew(gnx, ng);
780         srenew(index, ng);
781         srenew(grpname, ng);
782         for (i = 1; (i < ng); i++)
783         {
784             gnx[i]      = 1;
785             grpname[i]  = grpname[0];
786             snew(index[i], 1);
787             index[i][0] = index[0][i];
788         }
789         gnx[0] = 1;
790     }
791
792     if (resfnm)
793     {
794         GMX_RELEASE_ASSERT(top != nullptr, "top pointer cannot be NULL when finding residues");
795         nres = find_residues(&(top->atoms), gnx[0], index[0], &residues);
796
797         if (debug)
798         {
799             dump_res(debug, nres, residues, index[0]);
800         }
801     }
802
803     if (bPI)
804     {
805         periodic_mindist_plot(trxfnm, distfnm, top, ePBC, gnx[0], index[0], bSplit, oenv);
806     }
807     else
808     {
809         dist_plot(trxfnm, atmfnm, distfnm, numfnm, resfnm, oxfnm,
810                   rcutoff, bMat, top ? &(top->atoms) : nullptr,
811                   ng, index, gnx, grpname, bSplit, !bMax, nres, residues, bPBC, ePBC,
812                   bGroup, bEachResEachTime, bPrintResName, oenv);
813     }
814
815     do_view(oenv, distfnm, "-nxy");
816     if (!bPI)
817     {
818         do_view(oenv, numfnm, "-nxy");
819     }
820
821     return 0;
822 }