Add replacements for pbc enumerations
[gromacs.git] / src / gromacs / gmxana / gmx_mindist.cpp
blob52cd5148f0fbe05c0ce5c177348c2599b7302df5
1 /*
2 * This file is part of the GROMACS molecular simulation package.
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,2019, 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.
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.
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.
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.
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.
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.
37 #include "gmxpre.h"
39 #include <cmath>
40 #include <cstdlib>
41 #include <cstring>
43 #include <algorithm>
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"
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)
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;
75 sqr_box = std::min(norm2(box[XX]), norm2(box[YY]));
76 if (ePBC == epbcXYZ)
78 sqr_box = std::min(sqr_box, norm2(box[ZZ]));
79 nsz = 1;
81 else if (ePBC == epbcXY)
83 nsz = 0;
85 else
87 gmx_fatal(FARGS, "pbc = %s is not supported by g_mindist",
88 epbc_names[ePBC]);
91 nshift = 0;
92 for (sz = -nsz; sz <= nsz; sz++)
94 for (sy = -1; sy <= 1; sy++)
96 for (sx = -1; sx <= 1; sx++)
98 if (sx != 0 || sy != 0 || sz != 0)
100 for (i = 0; i < DIM; i++)
102 shift[nshift][i] =
103 sx*box[XX][i] + sy*box[YY][i] + sz*box[ZZ][i];
105 nshift++;
111 r2min = sqr_box;
112 r2max = 0;
114 for (i = 0; i < n; i++)
116 for (j = i+1; j < n; j++)
118 rvec_sub(x[index[i]], x[index[j]], d0);
119 r2 = norm2(d0);
120 if (r2 > r2max)
122 r2max = r2;
124 for (s = 0; s < nshift; s++)
126 rvec_add(d0, shift[s], d);
127 r2 = norm2(d);
128 if (r2 < r2min)
130 r2min = r2;
131 min_ind[0] = i;
132 min_ind[1] = j;
138 *rmin = std::sqrt(r2min);
139 *rmax = std::sqrt(r2max);
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)
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;
158 natoms = read_first_x(oenv, &status, trxfn, &t, &x, box);
160 check_index(nullptr, n, index, nullptr, natoms);
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))
166 fprintf(out, "@ subtitle \"and maximum internal distance\"\n");
168 xvgr_legend(out, 5, leg, oenv);
170 rmint = box[XX][XX];
171 tmint = 0;
173 if (nullptr != top)
175 gpbc = gmx_rmpbc_init(&top->idef, ePBC, natoms);
178 bFirst = TRUE;
181 if (nullptr != top)
183 gmx_rmpbc(gpbc, natoms, box, x);
186 periodic_dist(ePBC, box, x, n, index, &rmin, &rmax, ind_min);
187 if (rmin < rmint)
189 rmint = rmin;
190 tmint = t;
191 ind_mini = ind_min[0];
192 ind_minj = ind_min[1];
194 if (bSplit && !bFirst && std::abs(t/output_env_get_time_factor(oenv)) < 1e-5)
196 fprintf(out, "%s\n", output_env_get_print_xvgr_codes(oenv) ? "&" : "");
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;
202 while (read_next_x(oenv, status, &t, x, box));
204 if (nullptr != top)
206 gmx_rmpbc_done(gpbc);
209 xvgrclose(out);
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);
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)
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;
232 *ixmin = -1;
233 *jxmin = -1;
234 *ixmax = -1;
235 *jxmax = -1;
236 *nmin = 0;
237 *nmax = 0;
239 rcut2 = gmx::square(rcut);
241 /* Must init pbc every step because of pressure coupling */
242 if (bPBC)
244 set_pbc(&pbc, ePBC, box);
246 if (index2)
248 i0 = 0;
249 j1 = nx2;
250 index3 = index2;
252 else
254 j1 = nx1;
255 index3 = index1;
258 rmin2 = 1e12;
259 rmax2 = -1e12;
261 for (j = 0; (j < j1); j++)
263 jx = index3[j];
264 if (index2 == nullptr)
266 i0 = j + 1;
268 nmin_j = 0;
269 nmax_j = 0;
270 for (i = i0; (i < nx1); i++)
272 ix = index1[i];
273 if (ix != jx)
275 if (bPBC)
277 pbc_dx(&pbc, x[ix], x[jx], dx);
279 else
281 rvec_sub(x[ix], x[jx], dx);
283 r2 = iprod(dx, dx);
284 if (r2 < rmin2)
286 rmin2 = r2;
287 *ixmin = ix;
288 *jxmin = jx;
290 if (r2 > rmax2)
292 rmax2 = r2;
293 *ixmax = ix;
294 *jxmax = jx;
296 if (r2 <= rcut2)
298 nmin_j++;
300 else
302 nmax_j++;
306 if (bGroup)
308 if (nmin_j > 0)
310 (*nmin)++;
312 if (nmax_j > 0)
314 (*nmax)++;
317 else
319 *nmin += nmin_j;
320 *nmax += nmax_j;
323 *rmin = std::sqrt(rmin2);
324 *rmax = std::sqrt(rmax2);
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)
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;
352 if (read_first_x(oenv, &status, fn, &t, &x0, box) == 0)
354 gmx_fatal(FARGS, "Could not read coordinates from statusfile\n");
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;
364 if (bMat)
366 if (ng == 1)
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)
374 xvgr_legend(num, 0, leg, oenv);
377 else
379 snew(leg, (ng*(ng-1))/2);
380 for (i = j = 0; (i < ng-1); i++)
382 for (k = i+1; (k < ng); k++, j++)
384 sprintf(buf, "%s-%s", grpn[i], grpn[k]);
385 leg[j] = gmx_strdup(buf);
388 xvgr_legend(dist, j, leg, oenv);
389 if (num)
391 xvgr_legend(num, j, leg, oenv);
395 else
397 snew(leg, ng-1);
398 for (i = 0; (i < ng-1); i++)
400 sprintf(buf, "%s-%s", grpn[0], grpn[i+1]);
401 leg[i] = gmx_strdup(buf);
403 xvgr_legend(dist, ng-1, leg, oenv);
404 if (num)
406 xvgr_legend(num, ng-1, leg, oenv);
410 if (bEachResEachTime)
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) )
417 fprintf(respertime, "# ");
419 for (j = 0; j < nres; j++)
421 fprintf(respertime, "%s%d ", *(atoms->resinfo[atoms->atom[index[0][residue[j]]].resind].name), atoms->atom[index[0][residue[j]]].resind);
423 fprintf(respertime, "\n");
428 if (nres)
430 snew(mindres, ng-1);
431 snew(maxdres, ng-1);
432 for (i = 1; i < ng; i++)
434 snew(mindres[i-1], nres);
435 snew(maxdres[i-1], nres);
436 for (j = 0; j < nres; j++)
438 mindres[i-1][j] = 1e6;
440 /* maxdres[*][*] is already 0 */
443 bFirst = TRUE;
446 if (bSplit && !bFirst && std::abs(t/output_env_get_time_factor(oenv)) < 1e-5)
448 fprintf(dist, "%s\n", output_env_get_print_xvgr_codes(oenv) ? "&" : "");
449 if (num)
451 fprintf(num, "%s\n", output_env_get_print_xvgr_codes(oenv) ? "&" : "");
453 if (atm)
455 fprintf(atm, "%s\n", output_env_get_print_xvgr_codes(oenv) ? "&" : "");
458 fprintf(dist, "%12e", output_env_conv_time(oenv, t));
459 if (num)
461 fprintf(num, "%12e", output_env_conv_time(oenv, t));
464 if (bMat)
466 if (ng == 1)
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)
473 fprintf(num, " %8d", bMin ? nmin : nmax);
476 else
478 for (i = 0; (i < ng-1); i++)
480 for (k = i+1; (k < ng); k++)
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)
487 fprintf(num, " %8d", bMin ? nmin : nmax);
493 else
495 for (i = 1; (i < ng); i++)
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)
502 fprintf(num, " %8d", bMin ? nmin : nmax);
504 if (nres)
506 for (j = 0; j < nres; j++)
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);
517 fprintf(dist, "\n");
518 if (num)
520 fprintf(num, "\n");
522 if ( (bMin ? min1 : max1) != -1)
524 if (atm)
526 fprintf(atm, "%12e %12d %12d\n",
527 output_env_conv_time(oenv, t), 1+(bMin ? min1 : max1),
528 1+(bMin ? min2 : max2));
532 if (trxout)
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);
538 bFirst = FALSE;
539 /*dmin should be minimum distance for residue and group*/
540 if (bEachResEachTime)
542 fprintf(respertime, "%12e", t);
543 for (i = 1; i < ng; i++)
545 for (j = 0; j < nres; j++)
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;
553 fprintf(respertime, "\n");
556 while (read_next_x(oenv, status, &t, x0, box));
558 close_trx(status);
559 xvgrclose(dist);
560 if (num)
562 xvgrclose(num);
564 if (atm)
566 gmx_ffclose(atm);
568 if (trxout)
570 close_trx(trxout);
572 if (respertime)
574 xvgrclose(respertime);
577 if (nres && !bEachResEachTime)
579 FILE *res;
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++)
586 fprintf(res, "%4d", j+1);
587 for (i = 1; i < ng; i++)
589 fprintf(res, " %7g", bMin ? mindres[i-1][j] : maxdres[i-1][j]);
591 fprintf(res, "\n");
593 xvgrclose(res);
596 if (x0)
598 sfree(x0);
601 int freeLeg = bMat ? (ng == 1 ? 1 : (ng*(ng-1))/2) : ng - 1;
602 for (int i = 0; i < freeLeg; i++)
604 sfree(leg[i]);
606 sfree(leg);
609 static int find_residues(const t_atoms *atoms, int n, const int index[], int **resindex)
611 int i;
612 int nres = 0, resnr, presnr = 0;
613 bool presFound = false;
614 int *residx;
616 /* build index of first atom numbers for each residue */
617 snew(residx, atoms->nres+1);
618 for (i = 0; i < n; i++)
620 resnr = atoms->atom[index[i]].resind;
621 if (!presFound || resnr != presnr)
623 residx[nres] = i;
624 nres++;
625 presnr = resnr;
626 presFound = true;
629 if (debug)
631 printf("Found %d residues out of %d (%d/%d atoms)\n",
632 nres, atoms->nres, atoms->nr, n);
634 srenew(residx, nres+1);
635 /* mark end of last residue */
636 residx[nres] = n;
637 *resindex = residx;
638 return nres;
641 static void dump_res(FILE *out, int nres, int *resindex, int index[])
643 int i, j;
645 for (i = 0; i < nres-1; i++)
647 fprintf(out, "Res %d (%d):", i, resindex[i+1]-resindex[i]);
648 for (j = resindex[i]; j < resindex[i+1]; j++)
650 fprintf(out, " %d(%d)", j, index[j]);
652 fprintf(out, "\n");
656 int gmx_mindist(int argc, char *argv[])
658 const char *desc[] = {
659 "[THISMODULE] computes the distance between one group and a number of",
660 "other groups. Both the minimum distance",
661 "(between any pair of atoms from the respective groups)",
662 "and the number of contacts within a given",
663 "distance are written to two separate output files.",
664 "With the [TT]-group[tt] option a contact of an atom in another group",
665 "with multiple atoms in the first group is counted as one contact",
666 "instead of as multiple contacts.",
667 "With [TT]-or[tt], minimum distances to each residue in the first",
668 "group are determined and plotted as a function of residue number.[PAR]",
669 "With option [TT]-pi[tt] the minimum distance of a group to its",
670 "periodic image is plotted. This is useful for checking if a protein",
671 "has seen its periodic image during a simulation. Only one shift in",
672 "each direction is considered, giving a total of 26 shifts. Note",
673 "that periodicity information is required from the file supplied with",
674 "with [TT]-s[tt], either as a .tpr file or a .pdb file with CRYST1 fields.",
675 "It also plots the maximum distance within the group and the lengths",
676 "of the three box vectors.[PAR]",
677 "Also [gmx-distance] and [gmx-pairdist] calculate distances."
680 gmx_bool bMat = FALSE, bPI = FALSE, bSplit = FALSE, bMax = FALSE, bPBC = TRUE;
681 gmx_bool bGroup = FALSE;
682 real rcutoff = 0.6;
683 int ng = 1;
684 gmx_bool bEachResEachTime = FALSE, bPrintResName = FALSE;
685 t_pargs pa[] = {
686 { "-matrix", FALSE, etBOOL, {&bMat},
687 "Calculate half a matrix of group-group distances" },
688 { "-max", FALSE, etBOOL, {&bMax},
689 "Calculate *maximum* distance instead of minimum" },
690 { "-d", FALSE, etREAL, {&rcutoff},
691 "Distance for contacts" },
692 { "-group", FALSE, etBOOL, {&bGroup},
693 "Count contacts with multiple atoms in the first group as one" },
694 { "-pi", FALSE, etBOOL, {&bPI},
695 "Calculate minimum distance with periodic images" },
696 { "-split", FALSE, etBOOL, {&bSplit},
697 "Split graph where time is zero" },
698 { "-ng", FALSE, etINT, {&ng},
699 "Number of secondary groups to compute distance to a central group" },
700 { "-pbc", FALSE, etBOOL, {&bPBC},
701 "Take periodic boundary conditions into account" },
702 { "-respertime", FALSE, etBOOL, {&bEachResEachTime},
703 "When writing per-residue distances, write distance for each time point" },
704 { "-printresname", FALSE, etBOOL, {&bPrintResName},
705 "Write residue names" }
707 gmx_output_env_t *oenv;
708 t_topology *top = nullptr;
709 int ePBC = -1;
710 rvec *x = nullptr;
711 matrix box;
712 gmx_bool bTop = FALSE;
714 int i, nres = 0;
715 const char *trxfnm, *tpsfnm, *ndxfnm, *distfnm, *numfnm, *atmfnm, *oxfnm, *resfnm;
716 char **grpname;
717 int *gnx;
718 int **index, *residues = nullptr;
719 t_filenm fnm[] = {
720 { efTRX, "-f", nullptr, ffREAD },
721 { efTPS, nullptr, nullptr, ffOPTRD },
722 { efNDX, nullptr, nullptr, ffOPTRD },
723 { efXVG, "-od", "mindist", ffWRITE },
724 { efXVG, "-on", "numcont", ffOPTWR },
725 { efOUT, "-o", "atm-pair", ffOPTWR },
726 { efTRO, "-ox", "mindist", ffOPTWR },
727 { efXVG, "-or", "mindistres", ffOPTWR }
729 #define NFILE asize(fnm)
731 if (!parse_common_args(&argc, argv,
732 PCA_CAN_VIEW | PCA_CAN_TIME | PCA_TIME_UNIT,
733 NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, nullptr, &oenv))
735 return 0;
738 trxfnm = ftp2fn(efTRX, NFILE, fnm);
739 ndxfnm = ftp2fn_null(efNDX, NFILE, fnm);
740 distfnm = opt2fn("-od", NFILE, fnm);
741 numfnm = opt2fn_null("-on", NFILE, fnm);
742 atmfnm = ftp2fn_null(efOUT, NFILE, fnm);
743 oxfnm = opt2fn_null("-ox", NFILE, fnm);
744 resfnm = opt2fn_null("-or", NFILE, fnm);
745 if (bPI || resfnm != nullptr)
747 /* We need a tps file */
748 tpsfnm = ftp2fn(efTPS, NFILE, fnm);
750 else
752 tpsfnm = ftp2fn_null(efTPS, NFILE, fnm);
755 if (!tpsfnm && !ndxfnm)
757 gmx_fatal(FARGS, "You have to specify either the index file or a tpr file");
760 if (bPI)
762 ng = 1;
763 fprintf(stderr, "Choose a group for distance calculation\n");
765 else if (!bMat)
767 ng++;
770 snew(gnx, ng);
771 snew(index, ng);
772 snew(grpname, ng);
774 if (tpsfnm || resfnm || !ndxfnm)
776 snew(top, 1);
777 bTop = read_tps_conf(tpsfnm, top, &ePBC, &x, nullptr, box, FALSE);
778 if (bPI && !bTop)
780 printf("\nWARNING: Without a run input file a trajectory with broken molecules will not give the correct periodic image distance\n\n");
783 get_index(top ? &(top->atoms) : nullptr, ndxfnm, ng, gnx, index, grpname);
785 if (bMat && (ng == 1))
787 ng = gnx[0];
788 printf("Special case: making distance matrix between all atoms in group %s\n",
789 grpname[0]);
790 srenew(gnx, ng);
791 srenew(index, ng);
792 srenew(grpname, ng);
793 for (i = 1; (i < ng); i++)
795 gnx[i] = 1;
796 grpname[i] = grpname[0];
797 snew(index[i], 1);
798 index[i][0] = index[0][i];
800 gnx[0] = 1;
803 if (resfnm)
805 GMX_RELEASE_ASSERT(top != nullptr, "top pointer cannot be NULL when finding residues");
806 nres = find_residues(&(top->atoms), gnx[0], index[0], &residues);
808 if (debug)
810 dump_res(debug, nres, residues, index[0]);
813 else if (bEachResEachTime || bPrintResName)
815 gmx_fatal(FARGS, "Option -or needs to be set to print residues");
818 if (bPI)
820 periodic_mindist_plot(trxfnm, distfnm, top, ePBC, gnx[0], index[0], bSplit, oenv);
822 else
824 dist_plot(trxfnm, atmfnm, distfnm, numfnm, resfnm, oxfnm,
825 rcutoff, bMat, top ? &(top->atoms) : nullptr,
826 ng, index, gnx, grpname, bSplit, !bMax, nres, residues, bPBC, ePBC,
827 bGroup, bEachResEachTime, bPrintResName, oenv);
830 do_view(oenv, distfnm, "-nxy");
831 if (!bPI)
833 do_view(oenv, numfnm, "-nxy");
836 output_env_done(oenv);
837 done_top(top);
838 for (int i = 0; i < ng; i++)
840 sfree(index[i]);
842 sfree(index);
843 sfree(gnx);
844 sfree(x);
845 sfree(grpname);
846 sfree(top);
848 return 0;