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, 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.
47 #include "gromacs/utility/smalloc.h"
52 #include "gromacs/fileio/futil.h"
53 #include "gromacs/commandline/pargs.h"
55 #include "gromacs/fileio/tpxio.h"
56 #include "gromacs/fileio/trxio.h"
61 static void periodic_dist(matrix box, rvec x[], int n, atom_id index[],
62 real *rmin, real *rmax, int *min_ind)
65 int sx, sy, sz, i, j, s;
66 real sqr_box, r2min, r2max, r2;
67 rvec shift[NSHIFT], d0, d;
69 sqr_box = sqr(min(norm(box[XX]), min(norm(box[YY]), norm(box[ZZ]))));
72 for (sz = -1; sz <= 1; sz++)
74 for (sy = -1; sy <= 1; sy++)
76 for (sx = -1; sx <= 1; sx++)
78 if (sx != 0 || sy != 0 || sz != 0)
80 for (i = 0; i < DIM; i++)
82 shift[s][i] = sx*box[XX][i]+sy*box[YY][i]+sz*box[ZZ][i];
93 for (i = 0; i < n; i++)
95 for (j = i+1; j < n; j++)
97 rvec_sub(x[index[i]], x[index[j]], d0);
103 for (s = 0; s < NSHIFT; s++)
105 rvec_add(d0, shift[s], d);
121 static void periodic_mindist_plot(const char *trxfn, const char *outfn,
122 t_topology *top, int ePBC,
123 int n, atom_id index[], gmx_bool bSplit,
124 const output_env_t oenv)
127 const char *leg[5] = { "min per.", "max int.", "box1", "box2", "box3" };
132 int natoms, ind_min[2] = {0, 0}, ind_mini = 0, ind_minj = 0;
133 real r, rmin, rmax, rmint, tmint;
135 gmx_rmpbc_t gpbc = NULL;
137 natoms = read_first_x(oenv, &status, trxfn, &t, &x, box);
139 check_index(NULL, n, index, NULL, natoms);
141 out = xvgropen(outfn, "Minimum distance to periodic image",
142 output_env_get_time_label(oenv), "Distance (nm)", oenv);
143 if (output_env_get_print_xvgr_codes(oenv))
145 fprintf(out, "@ subtitle \"and maximum internal distance\"\n");
147 xvgr_legend(out, 5, leg, oenv);
154 gpbc = gmx_rmpbc_init(&top->idef, ePBC, natoms);
162 gmx_rmpbc(gpbc, natoms, box, x);
165 periodic_dist(box, x, n, index, &rmin, &rmax, ind_min);
170 ind_mini = ind_min[0];
171 ind_minj = ind_min[1];
173 if (bSplit && !bFirst && abs(t/output_env_get_time_factor(oenv)) < 1e-5)
177 fprintf(out, "\t%g\t%6.3f %6.3f %6.3f %6.3f %6.3f\n",
178 output_env_conv_time(oenv, t), rmin, rmax, norm(box[0]), norm(box[1]), norm(box[2]));
181 while (read_next_x(oenv, status, &t, x, box));
185 gmx_rmpbc_done(gpbc);
191 "\nThe shortest periodic distance is %g (nm) at time %g (%s),\n"
192 "between atoms %d and %d\n",
193 rmint, output_env_conv_time(oenv, tmint), output_env_get_time_unit(oenv),
194 index[ind_mini]+1, index[ind_minj]+1);
197 static void calc_dist(real rcut, gmx_bool bPBC, int ePBC, matrix box, rvec x[],
198 int nx1, int nx2, atom_id index1[], atom_id index2[],
200 real *rmin, real *rmax, int *nmin, int *nmax,
201 int *ixmin, int *jxmin, int *ixmax, int *jxmax)
203 int i, j, i0 = 0, j1;
207 real r2, rmin2, rmax2, rcut2;
220 /* Must init pbc every step because of pressure coupling */
223 set_pbc(&pbc, ePBC, box);
240 for (j = 0; (j < j1); j++)
249 for (i = i0; (i < nx1); i++)
256 pbc_dx(&pbc, x[ix], x[jx], dx);
260 rvec_sub(x[ix], x[jx], dx);
306 void dist_plot(const char *fn, const char *afile, const char *dfile,
307 const char *nfile, const char *rfile, const char *xfile,
308 real rcut, gmx_bool bMat, t_atoms *atoms,
309 int ng, atom_id *index[], int gnx[], char *grpn[], gmx_bool bSplit,
310 gmx_bool bMin, int nres, atom_id *residue, gmx_bool bPBC, int ePBC,
311 gmx_bool bGroup, gmx_bool bEachResEachTime, gmx_bool bPrintResName,
312 const output_env_t oenv)
314 FILE *atm, *dist, *num;
318 real t, dmin, dmax, **mindres = NULL, **maxdres = NULL;
321 int i = -1, j, k, natoms;
322 int min1, min2, max1, max2, min1r, min2r, max1r, max2r;
328 FILE *respertime = NULL;
330 if ((natoms = read_first_x(oenv, &status, fn, &t, &x0, box)) == 0)
332 gmx_fatal(FARGS, "Could not read coordinates from statusfile\n");
335 sprintf(buf, "%simum Distance", bMin ? "Min" : "Max");
336 dist = xvgropen(dfile, buf, output_env_get_time_label(oenv), "Distance (nm)", oenv);
337 sprintf(buf, "Number of Contacts %s %g nm", bMin ? "<" : ">", rcut);
338 num = nfile ? xvgropen(nfile, buf, output_env_get_time_label(oenv), "Number", oenv) : NULL;
339 atm = afile ? gmx_ffopen(afile, "w") : NULL;
340 trxout = xfile ? open_trx(xfile, "w") : NULL;
347 sprintf(buf, "Internal in %s", grpn[0]);
348 leg[0] = strdup(buf);
349 xvgr_legend(dist, 0, (const char**)leg, oenv);
352 xvgr_legend(num, 0, (const char**)leg, oenv);
357 snew(leg, (ng*(ng-1))/2);
358 for (i = j = 0; (i < ng-1); i++)
360 for (k = i+1; (k < ng); k++, j++)
362 sprintf(buf, "%s-%s", grpn[i], grpn[k]);
363 leg[j] = strdup(buf);
366 xvgr_legend(dist, j, (const char**)leg, oenv);
369 xvgr_legend(num, j, (const char**)leg, oenv);
376 for (i = 0; (i < ng-1); i++)
378 sprintf(buf, "%s-%s", grpn[0], grpn[i+1]);
379 leg[i] = strdup(buf);
381 xvgr_legend(dist, ng-1, (const char**)leg, oenv);
384 xvgr_legend(num, ng-1, (const char**)leg, oenv);
388 if (bEachResEachTime)
390 sprintf(buf, "%simum Distance", bMin ? "Min" : "Max");
391 respertime = xvgropen(rfile, buf, output_env_get_time_label(oenv), "Distance (nm)", oenv);
392 xvgr_legend(respertime, ng-1, (const char**)leg, oenv);
395 fprintf(respertime, "# ");
397 for (j = 0; j < nres; j++)
399 fprintf(respertime, "%s%d ", *(atoms->resinfo[atoms->atom[index[0][residue[j]]].resind].name), atoms->atom[index[0][residue[j]]].resind);
401 fprintf(respertime, "\n");
410 for (i = 1; i < ng; i++)
412 snew(mindres[i-1], nres);
413 snew(maxdres[i-1], nres);
414 for (j = 0; j < nres; j++)
416 mindres[i-1][j] = 1e6;
418 /* maxdres[*][*] is already 0 */
424 if (bSplit && !bFirst && abs(t/output_env_get_time_factor(oenv)) < 1e-5)
426 fprintf(dist, "&\n");
436 fprintf(dist, "%12e", output_env_conv_time(oenv, t));
439 fprintf(num, "%12e", output_env_conv_time(oenv, t));
446 calc_dist(rcut, bPBC, ePBC, box, x0, gnx[0], gnx[0], index[0], index[0], bGroup,
447 &dmin, &dmax, &nmin, &nmax, &min1, &min2, &max1, &max2);
448 fprintf(dist, " %12e", bMin ? dmin : dmax);
451 fprintf(num, " %8d", bMin ? nmin : nmax);
456 for (i = 0; (i < ng-1); i++)
458 for (k = i+1; (k < ng); k++)
460 calc_dist(rcut, bPBC, ePBC, box, x0, gnx[i], gnx[k], index[i], index[k],
461 bGroup, &dmin, &dmax, &nmin, &nmax, &min1, &min2, &max1, &max2);
462 fprintf(dist, " %12e", bMin ? dmin : dmax);
465 fprintf(num, " %8d", bMin ? nmin : nmax);
473 for (i = 1; (i < ng); i++)
475 calc_dist(rcut, bPBC, ePBC, box, x0, gnx[0], gnx[i], index[0], index[i], bGroup,
476 &dmin, &dmax, &nmin, &nmax, &min1, &min2, &max1, &max2);
477 fprintf(dist, " %12e", bMin ? dmin : dmax);
480 fprintf(num, " %8d", bMin ? nmin : nmax);
484 for (j = 0; j < nres; j++)
486 calc_dist(rcut, bPBC, ePBC, box, x0, residue[j+1]-residue[j], gnx[i],
487 &(index[0][residue[j]]), index[i], bGroup,
488 &dmin, &dmax, &nmin, &nmax, &min1r, &min2r, &max1r, &max2r);
489 mindres[i-1][j] = min(mindres[i-1][j], dmin);
490 maxdres[i-1][j] = max(maxdres[i-1][j], dmax);
500 if ( (bMin ? min1 : max1) != -1)
504 fprintf(atm, "%12e %12d %12d\n",
505 output_env_conv_time(oenv, t), 1+(bMin ? min1 : max1),
506 1+(bMin ? min2 : max2));
512 oindex[0] = bMin ? min1 : max1;
513 oindex[1] = bMin ? min2 : max2;
514 write_trx(trxout, 2, oindex, atoms, i, t, box, x0, NULL, NULL);
517 /*dmin should be minimum distance for residue and group*/
518 if (bEachResEachTime)
520 fprintf(respertime, "%12e", t);
521 for (i = 1; i < ng; i++)
523 for (j = 0; j < nres; j++)
525 fprintf(respertime, " %7g", bMin ? mindres[i-1][j] : maxdres[i-1][j]);
526 /*reset distances for next time point*/
527 mindres[i-1][j] = 1e6;
531 fprintf(respertime, "\n");
534 while (read_next_x(oenv, status, &t, x0, box));
551 if (nres && !bEachResEachTime)
555 sprintf(buf, "%simum Distance", bMin ? "Min" : "Max");
556 res = xvgropen(rfile, buf, "Residue (#)", "Distance (nm)", oenv);
557 xvgr_legend(res, ng-1, (const char**)leg, oenv);
558 for (j = 0; j < nres; j++)
560 fprintf(res, "%4d", j+1);
561 for (i = 1; i < ng; i++)
563 fprintf(res, " %7g", bMin ? mindres[i-1][j] : maxdres[i-1][j]);
572 int find_residues(t_atoms *atoms, int n, atom_id index[], atom_id **resindex)
575 int nres = 0, resnr, presnr;
578 /* build index of first atom numbers for each residue */
580 snew(residx, atoms->nres+1);
581 for (i = 0; i < n; i++)
583 resnr = atoms->atom[index[i]].resind;
593 printf("Found %d residues out of %d (%d/%d atoms)\n",
594 nres, atoms->nres, atoms->nr, n);
596 srenew(residx, nres+1);
597 /* mark end of last residue */
603 void dump_res(FILE *out, int nres, atom_id *resindex, atom_id index[])
607 for (i = 0; i < nres-1; i++)
609 fprintf(out, "Res %d (%d):", i, resindex[i+1]-resindex[i]);
610 for (j = resindex[i]; j < resindex[i+1]; j++)
612 fprintf(out, " %d(%d)", j, index[j]);
618 int gmx_mindist(int argc, char *argv[])
620 const char *desc[] = {
621 "[THISMODULE] computes the distance between one group and a number of",
622 "other groups. Both the minimum distance",
623 "(between any pair of atoms from the respective groups)",
624 "and the number of contacts within a given",
625 "distance are written to two separate output files.",
626 "With the [TT]-group[tt] option a contact of an atom in another group",
627 "with multiple atoms in the first group is counted as one contact",
628 "instead of as multiple contacts.",
629 "With [TT]-or[tt], minimum distances to each residue in the first",
630 "group are determined and plotted as a function of residue number.[PAR]",
631 "With option [TT]-pi[tt] the minimum distance of a group to its",
632 "periodic image is plotted. This is useful for checking if a protein",
633 "has seen its periodic image during a simulation. Only one shift in",
634 "each direction is considered, giving a total of 26 shifts.",
635 "It also plots the maximum distance within the group and the lengths",
636 "of the three box vectors.[PAR]",
637 "Also [gmx-distance] calculates distances."
639 const char *bugs[] = {
640 "The [TT]-pi[tt] option is very slow."
643 static gmx_bool bMat = FALSE, bPI = FALSE, bSplit = FALSE, bMax = FALSE, bPBC = TRUE;
644 static gmx_bool bGroup = FALSE;
645 static real rcutoff = 0.6;
647 static gmx_bool bEachResEachTime = FALSE, bPrintResName = FALSE;
649 { "-matrix", FALSE, etBOOL, {&bMat},
650 "Calculate half a matrix of group-group distances" },
651 { "-max", FALSE, etBOOL, {&bMax},
652 "Calculate *maximum* distance instead of minimum" },
653 { "-d", FALSE, etREAL, {&rcutoff},
654 "Distance for contacts" },
655 { "-group", FALSE, etBOOL, {&bGroup},
656 "Count contacts with multiple atoms in the first group as one" },
657 { "-pi", FALSE, etBOOL, {&bPI},
658 "Calculate minimum distance with periodic images" },
659 { "-split", FALSE, etBOOL, {&bSplit},
660 "Split graph where time is zero" },
661 { "-ng", FALSE, etINT, {&ng},
662 "Number of secondary groups to compute distance to a central group" },
663 { "-pbc", FALSE, etBOOL, {&bPBC},
664 "Take periodic boundary conditions into account" },
665 { "-respertime", FALSE, etBOOL, {&bEachResEachTime},
666 "When writing per-residue distances, write distance for each time point" },
667 { "-printresname", FALSE, etBOOL, {&bPrintResName},
668 "Write residue names" }
671 t_topology *top = NULL;
677 gmx_bool bTop = FALSE;
681 const char *trxfnm, *tpsfnm, *ndxfnm, *distfnm, *numfnm, *atmfnm, *oxfnm, *resfnm;
684 atom_id **index, *residues = NULL;
686 { efTRX, "-f", NULL, ffREAD },
687 { efTPS, NULL, NULL, ffOPTRD },
688 { efNDX, NULL, NULL, ffOPTRD },
689 { efXVG, "-od", "mindist", ffWRITE },
690 { efXVG, "-on", "numcont", ffOPTWR },
691 { efOUT, "-o", "atm-pair", ffOPTWR },
692 { efTRO, "-ox", "mindist", ffOPTWR },
693 { efXVG, "-or", "mindistres", ffOPTWR }
695 #define NFILE asize(fnm)
697 if (!parse_common_args(&argc, argv,
698 PCA_CAN_VIEW | PCA_CAN_TIME | PCA_TIME_UNIT | PCA_BE_NICE,
699 NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL, &oenv))
704 trxfnm = ftp2fn(efTRX, NFILE, fnm);
705 ndxfnm = ftp2fn_null(efNDX, NFILE, fnm);
706 distfnm = opt2fn("-od", NFILE, fnm);
707 numfnm = opt2fn_null("-on", NFILE, fnm);
708 atmfnm = ftp2fn_null(efOUT, NFILE, fnm);
709 oxfnm = opt2fn_null("-ox", NFILE, fnm);
710 resfnm = opt2fn_null("-or", NFILE, fnm);
711 if (bPI || resfnm != NULL)
713 /* We need a tps file */
714 tpsfnm = ftp2fn(efTPS, NFILE, fnm);
718 tpsfnm = ftp2fn_null(efTPS, NFILE, fnm);
721 if (!tpsfnm && !ndxfnm)
723 gmx_fatal(FARGS, "You have to specify either the index file or a tpr file");
729 fprintf(stderr, "Choose a group for distance calculation\n");
740 if (tpsfnm || resfnm || !ndxfnm)
743 bTop = read_tps_conf(tpsfnm, title, top, &ePBC, &x, NULL, box, FALSE);
746 printf("\nWARNING: Without a run input file a trajectory with broken molecules will not give the correct periodic image distance\n\n");
749 get_index(top ? &(top->atoms) : NULL, ndxfnm, ng, gnx, index, grpname);
751 if (bMat && (ng == 1))
754 printf("Special case: making distance matrix between all atoms in group %s\n",
759 for (i = 1; (i < ng); i++)
762 grpname[i] = grpname[0];
764 index[i][0] = index[0][i];
771 nres = find_residues(top ? &(top->atoms) : NULL,
772 gnx[0], index[0], &residues);
775 dump_res(debug, nres, residues, index[0]);
781 periodic_mindist_plot(trxfnm, distfnm, top, ePBC, gnx[0], index[0], bSplit, oenv);
785 dist_plot(trxfnm, atmfnm, distfnm, numfnm, resfnm, oxfnm,
786 rcutoff, bMat, top ? &(top->atoms) : NULL,
787 ng, index, gnx, grpname, bSplit, !bMax, nres, residues, bPBC, ePBC,
788 bGroup, bEachResEachTime, bPrintResName, oenv);
791 do_view(oenv, distfnm, "-nxy");
794 do_view(oenv, numfnm, "-nxy");