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, 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.
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)
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]));
78 sqr_box = std::min(sqr_box, norm2(box[ZZ]));
81 else if (ePBC == epbcXY)
87 gmx_fatal(FARGS, "pbc = %s is not supported by g_mindist",
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++)
103 sx*box[XX][i] + sy*box[YY][i] + sz*box[ZZ][i];
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);
124 for (s = 0; s < nshift; s++)
126 rvec_add(d0, shift[s], d);
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)
148 const char *leg[5] = { "min per.", "max int.", "box1", "box2", "box3" };
153 int natoms, ind_min[2] = {0, 0}, ind_mini = 0, ind_minj = 0;
154 real rmin, rmax, rmint, tmint;
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);
175 gpbc = gmx_rmpbc_init(&top->idef, ePBC, natoms);
183 gmx_rmpbc(gpbc, natoms, box, x);
186 periodic_dist(ePBC, box, x, n, index, &rmin, &rmax, ind_min);
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]));
202 while (read_next_x(oenv, status, &t, x, box));
206 gmx_rmpbc_done(gpbc);
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[],
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;
228 real r2, rmin2, rmax2, rcut2;
239 rcut2 = gmx::square(rcut);
241 /* Must init pbc every step because of pressure coupling */
244 set_pbc(&pbc, ePBC, box);
261 for (j = 0; (j < j1); j++)
264 if (index2 == nullptr)
270 for (i = i0; (i < nx1); i++)
277 pbc_dx(&pbc, x[ix], x[jx], dx);
281 rvec_sub(x[ix], x[jx], dx);
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;
339 real t, dmin, dmax, **mindres = nullptr, **maxdres = nullptr;
343 int min2, max2, min1r, min2r, max1r, max2r;
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;
369 sprintf(buf, "Internal in %s", grpn[0]);
370 leg[0] = gmx_strdup(buf);
371 xvgr_legend(dist, 0, leg, oenv);
374 xvgr_legend(num, 0, leg, oenv);
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);
391 xvgr_legend(num, j, leg, oenv);
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);
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");
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 */
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) ? "&" : "");
451 fprintf(num, "%s\n", output_env_get_print_xvgr_codes(oenv) ? "&" : "");
455 fprintf(atm, "%s\n", output_env_get_print_xvgr_codes(oenv) ? "&" : "");
458 fprintf(dist, "%12e", output_env_conv_time(oenv, t));
461 fprintf(num, "%12e", output_env_conv_time(oenv, t));
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);
473 fprintf(num, " %8d", bMin ? nmin : nmax);
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);
487 fprintf(num, " %8d", bMin ? nmin : nmax);
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);
502 fprintf(num, " %8d", bMin ? nmin : nmax);
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);
522 if ( (bMin ? min1 : max1) != -1)
526 fprintf(atm, "%12e %12d %12d\n",
527 output_env_conv_time(oenv, t), 1+(bMin ? min1 : max1),
528 1+(bMin ? min2 : max2));
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);
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;
553 fprintf(respertime, "\n");
556 while (read_next_x(oenv, status, &t, x0, box));
574 xvgrclose(respertime);
577 if (nres && !bEachResEachTime)
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]);
599 static int find_residues(const t_atoms *atoms, int n, const int index[], int **resindex)
602 int nres = 0, resnr, presnr = 0;
603 bool presFound = false;
606 /* build index of first atom numbers for each residue */
607 snew(residx, atoms->nres+1);
608 for (i = 0; i < n; i++)
610 resnr = atoms->atom[index[i]].resind;
611 if (!presFound || resnr != presnr)
621 printf("Found %d residues out of %d (%d/%d atoms)\n",
622 nres, atoms->nres, atoms->nr, n);
624 srenew(residx, nres+1);
625 /* mark end of last residue */
631 static void dump_res(FILE *out, int nres, int *resindex, int index[])
635 for (i = 0; i < nres-1; i++)
637 fprintf(out, "Res %d (%d):", i, resindex[i+1]-resindex[i]);
638 for (j = resindex[i]; j < resindex[i+1]; j++)
640 fprintf(out, " %d(%d)", j, index[j]);
646 int gmx_mindist(int argc, char *argv[])
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."
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;
674 static gmx_bool bEachResEachTime = FALSE, bPrintResName = FALSE;
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" }
697 gmx_output_env_t *oenv;
698 t_topology *top = nullptr;
702 gmx_bool bTop = FALSE;
705 const char *trxfnm, *tpsfnm, *ndxfnm, *distfnm, *numfnm, *atmfnm, *oxfnm, *resfnm;
708 int **index, *residues = nullptr;
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 }
719 #define NFILE asize(fnm)
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))
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)
737 /* We need a tps file */
738 tpsfnm = ftp2fn(efTPS, NFILE, fnm);
742 tpsfnm = ftp2fn_null(efTPS, NFILE, fnm);
745 if (!tpsfnm && !ndxfnm)
747 gmx_fatal(FARGS, "You have to specify either the index file or a tpr file");
753 fprintf(stderr, "Choose a group for distance calculation\n");
764 if (tpsfnm || resfnm || !ndxfnm)
767 bTop = read_tps_conf(tpsfnm, top, &ePBC, &x, nullptr, box, FALSE);
770 printf("\nWARNING: Without a run input file a trajectory with broken molecules will not give the correct periodic image distance\n\n");
773 get_index(top ? &(top->atoms) : nullptr, ndxfnm, ng, gnx, index, grpname);
775 if (bMat && (ng == 1))
778 printf("Special case: making distance matrix between all atoms in group %s\n",
783 for (i = 1; (i < ng); i++)
786 grpname[i] = grpname[0];
788 index[i][0] = index[0][i];
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);
800 dump_res(debug, nres, residues, index[0]);
803 else if (bEachResEachTime || bPrintResName)
805 gmx_fatal(FARGS, "Option -or needs to be set to print residues");
810 periodic_mindist_plot(trxfnm, distfnm, top, ePBC, gnx[0], index[0], bSplit, oenv);
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);
820 do_view(oenv, distfnm, "-nxy");
823 do_view(oenv, numfnm, "-nxy");