3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2004, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
33 * Green Red Orange Magenta Azure Cyan Skyblue
59 static void periodic_dist(matrix box, rvec x[], int n, atom_id index[],
60 real *rmin, real *rmax, int *min_ind)
63 int sx, sy, sz, i, j, s;
64 real sqr_box, r2min, r2max, r2;
65 rvec shift[NSHIFT], d0, d;
67 sqr_box = sqr(min(norm(box[XX]), min(norm(box[YY]), norm(box[ZZ]))));
70 for (sz = -1; sz <= 1; sz++)
72 for (sy = -1; sy <= 1; sy++)
74 for (sx = -1; sx <= 1; sx++)
76 if (sx != 0 || sy != 0 || sz != 0)
78 for (i = 0; i < DIM; i++)
80 shift[s][i] = sx*box[XX][i]+sy*box[YY][i]+sz*box[ZZ][i];
91 for (i = 0; i < n; i++)
93 for (j = i+1; j < n; j++)
95 rvec_sub(x[index[i]], x[index[j]], d0);
101 for (s = 0; s < NSHIFT; s++)
103 rvec_add(d0, shift[s], d);
119 static void periodic_mindist_plot(const char *trxfn, const char *outfn,
120 t_topology *top, int ePBC,
121 int n, atom_id index[], gmx_bool bSplit,
122 const output_env_t oenv)
125 const char *leg[5] = { "min per.", "max int.", "box1", "box2", "box3" };
130 int natoms, ind_min[2] = {0, 0}, ind_mini = 0, ind_minj = 0;
131 real r, rmin, rmax, rmint, tmint;
133 gmx_rmpbc_t gpbc = NULL;
135 natoms = read_first_x(oenv, &status, trxfn, &t, &x, box);
137 check_index(NULL, n, index, NULL, natoms);
139 out = xvgropen(outfn, "Minimum distance to periodic image",
140 output_env_get_time_label(oenv), "Distance (nm)", oenv);
141 if (output_env_get_print_xvgr_codes(oenv))
143 fprintf(out, "@ subtitle \"and maximum internal distance\"\n");
145 xvgr_legend(out, 5, leg, oenv);
152 gpbc = gmx_rmpbc_init(&top->idef, ePBC, natoms);
160 gmx_rmpbc(gpbc, natoms, box, x);
163 periodic_dist(box, x, n, index, &rmin, &rmax, ind_min);
168 ind_mini = ind_min[0];
169 ind_minj = ind_min[1];
171 if (bSplit && !bFirst && abs(t/output_env_get_time_factor(oenv)) < 1e-5)
175 fprintf(out, "\t%g\t%6.3f %6.3f %6.3f %6.3f %6.3f\n",
176 output_env_conv_time(oenv, t), rmin, rmax, norm(box[0]), norm(box[1]), norm(box[2]));
179 while (read_next_x(oenv, status, &t, x, box));
183 gmx_rmpbc_done(gpbc);
189 "\nThe shortest periodic distance is %g (nm) at time %g (%s),\n"
190 "between atoms %d and %d\n",
191 rmint, output_env_conv_time(oenv, tmint), output_env_get_time_unit(oenv),
192 index[ind_mini]+1, index[ind_minj]+1);
195 static void calc_dist(real rcut, gmx_bool bPBC, int ePBC, matrix box, rvec x[],
196 int nx1, int nx2, atom_id index1[], atom_id index2[],
198 real *rmin, real *rmax, int *nmin, int *nmax,
199 int *ixmin, int *jxmin, int *ixmax, int *jxmax)
201 int i, j, i0 = 0, j1;
205 real r2, rmin2, rmax2, rcut2;
218 /* Must init pbc every step because of pressure coupling */
221 set_pbc(&pbc, ePBC, box);
238 for (j = 0; (j < j1); j++)
247 for (i = i0; (i < nx1); i++)
254 pbc_dx(&pbc, x[ix], x[jx], dx);
258 rvec_sub(x[ix], x[jx], dx);
304 void dist_plot(const char *fn, const char *afile, const char *dfile,
305 const char *nfile, const char *rfile, const char *xfile,
306 real rcut, gmx_bool bMat, t_atoms *atoms,
307 int ng, atom_id *index[], int gnx[], char *grpn[], gmx_bool bSplit,
308 gmx_bool bMin, int nres, atom_id *residue, gmx_bool bPBC, int ePBC,
309 gmx_bool bGroup, gmx_bool bEachResEachTime, gmx_bool bPrintResName,
310 const output_env_t oenv)
312 FILE *atm, *dist, *num;
316 real t, dmin, dmax, **mindres = NULL, **maxdres = NULL;
319 int i = -1, j, k, natoms;
320 int min1, min2, max1, max2, min1r, min2r, max1r, max2r;
326 FILE *respertime = NULL;
328 if ((natoms = read_first_x(oenv, &status, fn, &t, &x0, box)) == 0)
330 gmx_fatal(FARGS, "Could not read coordinates from statusfile\n");
333 sprintf(buf, "%simum Distance", bMin ? "Min" : "Max");
334 dist = xvgropen(dfile, buf, output_env_get_time_label(oenv), "Distance (nm)", oenv);
335 sprintf(buf, "Number of Contacts %s %g nm", bMin ? "<" : ">", rcut);
336 num = nfile ? xvgropen(nfile, buf, output_env_get_time_label(oenv), "Number", oenv) : NULL;
337 atm = afile ? ffopen(afile, "w") : NULL;
338 trxout = xfile ? open_trx(xfile, "w") : NULL;
345 sprintf(buf, "Internal in %s", grpn[0]);
346 leg[0] = strdup(buf);
347 xvgr_legend(dist, 0, (const char**)leg, oenv);
350 xvgr_legend(num, 0, (const char**)leg, oenv);
355 snew(leg, (ng*(ng-1))/2);
356 for (i = j = 0; (i < ng-1); i++)
358 for (k = i+1; (k < ng); k++, j++)
360 sprintf(buf, "%s-%s", grpn[i], grpn[k]);
361 leg[j] = strdup(buf);
364 xvgr_legend(dist, j, (const char**)leg, oenv);
367 xvgr_legend(num, j, (const char**)leg, oenv);
374 for (i = 0; (i < ng-1); i++)
376 sprintf(buf, "%s-%s", grpn[0], grpn[i+1]);
377 leg[i] = strdup(buf);
379 xvgr_legend(dist, ng-1, (const char**)leg, oenv);
382 xvgr_legend(num, ng-1, (const char**)leg, oenv);
386 if (bEachResEachTime)
388 sprintf(buf, "%simum Distance", bMin ? "Min" : "Max");
389 respertime = xvgropen(rfile, buf, output_env_get_time_label(oenv), "Distance (nm)", oenv);
390 xvgr_legend(respertime, ng-1, (const char**)leg, oenv);
393 fprintf(respertime, "# ");
395 for (j = 0; j < nres; j++)
397 fprintf(respertime, "%s%d ", *(atoms->resinfo[atoms->atom[index[0][residue[j]]].resind].name), atoms->atom[index[0][residue[j]]].resind);
399 fprintf(respertime, "\n");
408 for (i = 1; i < ng; i++)
410 snew(mindres[i-1], nres);
411 snew(maxdres[i-1], nres);
412 for (j = 0; j < nres; j++)
414 mindres[i-1][j] = 1e6;
416 /* maxdres[*][*] is already 0 */
422 if (bSplit && !bFirst && abs(t/output_env_get_time_factor(oenv)) < 1e-5)
424 fprintf(dist, "&\n");
434 fprintf(dist, "%12e", output_env_conv_time(oenv, t));
437 fprintf(num, "%12e", output_env_conv_time(oenv, t));
444 calc_dist(rcut, bPBC, ePBC, box, x0, gnx[0], gnx[0], index[0], index[0], bGroup,
445 &dmin, &dmax, &nmin, &nmax, &min1, &min2, &max1, &max2);
446 fprintf(dist, " %12e", bMin ? dmin : dmax);
449 fprintf(num, " %8d", bMin ? nmin : nmax);
454 for (i = 0; (i < ng-1); i++)
456 for (k = i+1; (k < ng); k++)
458 calc_dist(rcut, bPBC, ePBC, box, x0, gnx[i], gnx[k], index[i], index[k],
459 bGroup, &dmin, &dmax, &nmin, &nmax, &min1, &min2, &max1, &max2);
460 fprintf(dist, " %12e", bMin ? dmin : dmax);
463 fprintf(num, " %8d", bMin ? nmin : nmax);
471 for (i = 1; (i < ng); i++)
473 calc_dist(rcut, bPBC, ePBC, box, x0, gnx[0], gnx[i], index[0], index[i], bGroup,
474 &dmin, &dmax, &nmin, &nmax, &min1, &min2, &max1, &max2);
475 fprintf(dist, " %12e", bMin ? dmin : dmax);
478 fprintf(num, " %8d", bMin ? nmin : nmax);
482 for (j = 0; j < nres; j++)
484 calc_dist(rcut, bPBC, ePBC, box, x0, residue[j+1]-residue[j], gnx[i],
485 &(index[0][residue[j]]), index[i], bGroup,
486 &dmin, &dmax, &nmin, &nmax, &min1r, &min2r, &max1r, &max2r);
487 mindres[i-1][j] = min(mindres[i-1][j], dmin);
488 maxdres[i-1][j] = max(maxdres[i-1][j], dmax);
498 if ( (bMin ? min1 : max1) != -1)
502 fprintf(atm, "%12e %12d %12d\n",
503 output_env_conv_time(oenv, t), 1+(bMin ? min1 : max1),
504 1+(bMin ? min2 : max2));
510 oindex[0] = bMin ? min1 : max1;
511 oindex[1] = bMin ? min2 : max2;
512 write_trx(trxout, 2, oindex, atoms, i, t, box, x0, NULL, NULL);
515 /*dmin should be minimum distance for residue and group*/
516 if (bEachResEachTime)
518 fprintf(respertime, "%12e", t);
519 for (i = 1; i < ng; i++)
521 for (j = 0; j < nres; j++)
523 fprintf(respertime, " %7g", bMin ? mindres[i-1][j] : maxdres[i-1][j]);
524 /*reset distances for next time point*/
525 mindres[i-1][j] = 1e6;
529 fprintf(respertime, "\n");
532 while (read_next_x(oenv, status, &t, x0, box));
549 if (nres && !bEachResEachTime)
553 sprintf(buf, "%simum Distance", bMin ? "Min" : "Max");
554 res = xvgropen(rfile, buf, "Residue (#)", "Distance (nm)", oenv);
555 xvgr_legend(res, ng-1, (const char**)leg, oenv);
556 for (j = 0; j < nres; j++)
558 fprintf(res, "%4d", j+1);
559 for (i = 1; i < ng; i++)
561 fprintf(res, " %7g", bMin ? mindres[i-1][j] : maxdres[i-1][j]);
570 int find_residues(t_atoms *atoms, int n, atom_id index[], atom_id **resindex)
573 int nres = 0, resnr, presnr;
576 /* build index of first atom numbers for each residue */
578 snew(residx, atoms->nres+1);
579 for (i = 0; i < n; i++)
581 resnr = atoms->atom[index[i]].resind;
591 printf("Found %d residues out of %d (%d/%d atoms)\n",
592 nres, atoms->nres, atoms->nr, n);
594 srenew(residx, nres+1);
595 /* mark end of last residue */
601 void dump_res(FILE *out, int nres, atom_id *resindex, atom_id index[])
605 for (i = 0; i < nres-1; i++)
607 fprintf(out, "Res %d (%d):", i, resindex[i+1]-resindex[i]);
608 for (j = resindex[i]; j < resindex[i+1]; j++)
610 fprintf(out, " %d(%d)", j, index[j]);
616 int gmx_mindist(int argc, char *argv[])
618 const char *desc[] = {
619 "[TT]g_mindist[tt] computes the distance between one group and a number of",
620 "other groups. Both the minimum distance",
621 "(between any pair of atoms from the respective groups)",
622 "and the number of contacts within a given",
623 "distance are written to two separate output files.",
624 "With the [TT]-group[tt] option a contact of an atom in another group",
625 "with multiple atoms in the first group is counted as one contact",
626 "instead of as multiple contacts.",
627 "With [TT]-or[tt], minimum distances to each residue in the first",
628 "group are determined and plotted as a function of residue number.[PAR]",
629 "With option [TT]-pi[tt] the minimum distance of a group to its",
630 "periodic image is plotted. This is useful for checking if a protein",
631 "has seen its periodic image during a simulation. Only one shift in",
632 "each direction is considered, giving a total of 26 shifts.",
633 "It also plots the maximum distance within the group and the lengths",
634 "of the three box vectors.[PAR]",
635 "Other programs that calculate distances are [TT]g_dist[tt]",
636 "and [TT]g_bond[tt]."
638 const char *bugs[] = {
639 "The [TT]-pi[tt] option is very slow."
642 static gmx_bool bMat = FALSE, bPI = FALSE, bSplit = FALSE, bMax = FALSE, bPBC = TRUE;
643 static gmx_bool bGroup = FALSE;
644 static real rcutoff = 0.6;
646 static gmx_bool bEachResEachTime = FALSE, bPrintResName = FALSE;
648 { "-matrix", FALSE, etBOOL, {&bMat},
649 "Calculate half a matrix of group-group distances" },
650 { "-max", FALSE, etBOOL, {&bMax},
651 "Calculate *maximum* distance instead of minimum" },
652 { "-d", FALSE, etREAL, {&rcutoff},
653 "Distance for contacts" },
654 { "-group", FALSE, etBOOL, {&bGroup},
655 "Count contacts with multiple atoms in the first group as one" },
656 { "-pi", FALSE, etBOOL, {&bPI},
657 "Calculate minimum distance with periodic images" },
658 { "-split", FALSE, etBOOL, {&bSplit},
659 "Split graph where time is zero" },
660 { "-ng", FALSE, etINT, {&ng},
661 "Number of secondary groups to compute distance to a central group" },
662 { "-pbc", FALSE, etBOOL, {&bPBC},
663 "Take periodic boundary conditions into account" },
664 { "-respertime", FALSE, etBOOL, {&bEachResEachTime},
665 "When writing per-residue distances, write distance for each time point" },
666 { "-printresname", FALSE, etBOOL, {&bPrintResName},
667 "Write residue names" }
670 t_topology *top = NULL;
676 gmx_bool bTop = FALSE;
680 const char *trxfnm, *tpsfnm, *ndxfnm, *distfnm, *numfnm, *atmfnm, *oxfnm, *resfnm;
683 atom_id **index, *residues = NULL;
685 { efTRX, "-f", NULL, ffREAD },
686 { efTPS, NULL, NULL, ffOPTRD },
687 { efNDX, NULL, NULL, ffOPTRD },
688 { efXVG, "-od", "mindist", ffWRITE },
689 { efXVG, "-on", "numcont", ffOPTWR },
690 { efOUT, "-o", "atm-pair", ffOPTWR },
691 { efTRO, "-ox", "mindist", ffOPTWR },
692 { efXVG, "-or", "mindistres", ffOPTWR }
694 #define NFILE asize(fnm)
696 if (!parse_common_args(&argc, argv,
697 PCA_CAN_VIEW | PCA_CAN_TIME | PCA_TIME_UNIT | PCA_BE_NICE,
698 NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL, &oenv))
703 trxfnm = ftp2fn(efTRX, NFILE, fnm);
704 ndxfnm = ftp2fn_null(efNDX, NFILE, fnm);
705 distfnm = opt2fn("-od", NFILE, fnm);
706 numfnm = opt2fn_null("-on", NFILE, fnm);
707 atmfnm = ftp2fn_null(efOUT, NFILE, fnm);
708 oxfnm = opt2fn_null("-ox", NFILE, fnm);
709 resfnm = opt2fn_null("-or", NFILE, fnm);
710 if (bPI || resfnm != NULL)
712 /* We need a tps file */
713 tpsfnm = ftp2fn(efTPS, NFILE, fnm);
717 tpsfnm = ftp2fn_null(efTPS, NFILE, fnm);
720 if (!tpsfnm && !ndxfnm)
722 gmx_fatal(FARGS, "You have to specify either the index file or a tpr file");
728 fprintf(stderr, "Choose a group for distance calculation\n");
739 if (tpsfnm || resfnm || !ndxfnm)
742 bTop = read_tps_conf(tpsfnm, title, top, &ePBC, &x, NULL, box, FALSE);
745 printf("\nWARNING: Without a run input file a trajectory with broken molecules will not give the correct periodic image distance\n\n");
748 get_index(top ? &(top->atoms) : NULL, ndxfnm, ng, gnx, index, grpname);
750 if (bMat && (ng == 1))
753 printf("Special case: making distance matrix between all atoms in group %s\n",
758 for (i = 1; (i < ng); i++)
761 grpname[i] = grpname[0];
763 index[i][0] = index[0][i];
770 nres = find_residues(top ? &(top->atoms) : NULL,
771 gnx[0], index[0], &residues);
774 dump_res(debug, nres, residues, index[0]);
780 periodic_mindist_plot(trxfnm, distfnm, top, ePBC, gnx[0], index[0], bSplit, oenv);
784 dist_plot(trxfnm, atmfnm, distfnm, numfnm, resfnm, oxfnm,
785 rcutoff, bMat, top ? &(top->atoms) : NULL,
786 ng, index, gnx, grpname, bSplit, !bMax, nres, residues, bPBC, ePBC,
787 bGroup, bEachResEachTime, bPrintResName, oenv);
790 do_view(oenv, distfnm, "-nxy");
793 do_view(oenv, numfnm, "-nxy");