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.
42 #include "gromacs/utility/smalloc.h"
45 #include "gromacs/commandline/pargs.h"
46 #include "gromacs/fileio/tpxio.h"
47 #include "gromacs/fileio/trxio.h"
48 #include "gromacs/fileio/strdb.h"
49 #include "gromacs/math/vec.h"
51 #include "gromacs/topology/index.h"
52 #include "gromacs/pbcutil/pbc.h"
53 #include "gromacs/fileio/xvgr.h"
55 #include "gromacs/utility/futil.h"
56 #include "gromacs/fileio/matio.h"
60 static void calc_dist(int nind, atom_id index[], rvec x[], int ePBC, matrix box,
69 set_pbc(&pbc, ePBC, box);
70 for (i = 0; (i < nind-1); i++)
73 for (j = i+1; (j < nind); j++)
75 pbc_dx(&pbc, xi, x[index[j]], dx);
77 d[i][j] = sqrt(temp2);
82 static void calc_dist_tot(int nind, atom_id index[], rvec x[],
84 real **d, real **dtot, real **dtot2,
85 gmx_bool bNMR, real **dtot1_3, real **dtot1_6)
89 real temp, temp2, temp1_3;
93 set_pbc(&pbc, ePBC, box);
94 for (i = 0; (i < nind-1); i++)
97 for (j = i+1; (j < nind); j++)
99 pbc_dx(&pbc, xi, x[index[j]], dx);
104 dtot2[i][j] += temp2;
107 temp1_3 = 1.0/(temp*temp2);
108 dtot1_3[i][j] += temp1_3;
109 dtot1_6[i][j] += temp1_3*temp1_3;
115 static void calc_nmr(int nind, int nframes, real **dtot1_3, real **dtot1_6,
116 real *max1_3, real *max1_6)
119 real temp1_3, temp1_6;
121 for (i = 0; (i < nind-1); i++)
123 for (j = i+1; (j < nind); j++)
125 temp1_3 = pow(dtot1_3[i][j]/nframes, -1.0/3.0);
126 temp1_6 = pow(dtot1_6[i][j]/nframes, -1.0/6.0);
127 if (temp1_3 > *max1_3)
131 if (temp1_6 > *max1_6)
135 dtot1_3[i][j] = temp1_3;
136 dtot1_6[i][j] = temp1_6;
137 dtot1_3[j][i] = temp1_3;
138 dtot1_6[j][i] = temp1_6;
143 static char Hnum[] = "123";
168 static int read_equiv(const char *eq_fn, t_equiv ***equivptr)
171 char line[STRLEN], resname[10], atomname[10], *lp;
172 int neq, na, n, resnr;
175 fp = gmx_ffopen(eq_fn, "r");
178 while (get_a_line(fp, line, STRLEN))
181 /* this is not efficient, but I'm lazy */
182 srenew(equiv, neq+1);
185 if (sscanf(lp, "%s %n", atomname, &n) == 1)
189 equiv[neq][0].nname = strdup(atomname);
190 while (sscanf(lp, "%d %s %s %n", &resnr, resname, atomname, &n) == 3)
192 /* this is not efficient, but I'm lazy (again) */
193 srenew(equiv[neq], na+1);
194 equiv[neq][na].rnr = resnr-1;
195 equiv[neq][na].rname = strdup(resname);
196 equiv[neq][na].aname = strdup(atomname);
199 equiv[neq][na].nname = NULL;
205 /* make empty element as flag for end of array */
206 srenew(equiv[neq], na+1);
207 equiv[neq][na].rnr = NOTSET;
208 equiv[neq][na].rname = NULL;
209 equiv[neq][na].aname = NULL;
221 static void dump_equiv(FILE *out, int neq, t_equiv **equiv)
225 fprintf(out, "Dumping equivalent list\n");
226 for (i = 0; i < neq; i++)
228 fprintf(out, "%s", equiv[i][0].nname);
229 for (j = 0; equiv[i][j].rnr != NOTSET; j++)
231 fprintf(out, " %d %s %s",
232 equiv[i][j].rnr, equiv[i][j].rname, equiv[i][j].aname);
238 static gmx_bool is_equiv(int neq, t_equiv **equiv, char **nname,
239 int rnr1, char *rname1, char *aname1,
240 int rnr2, char *rname2, char *aname2)
246 /* we can terminate each loop when bFound is true! */
247 for (i = 0; i < neq && !bFound; i++)
249 /* find first atom */
250 for (j = 0; equiv[i][j].rnr != NOTSET && !bFound; j++)
252 bFound = ( equiv[i][j].rnr == rnr1 &&
253 strcmp(equiv[i][j].rname, rname1) == 0 &&
254 strcmp(equiv[i][j].aname, aname1) == 0 );
258 /* find second atom */
260 for (j = 0; equiv[i][j].rnr != NOTSET && !bFound; j++)
262 bFound = ( equiv[i][j].rnr == rnr2 &&
263 strcmp(equiv[i][j].rname, rname2) == 0 &&
264 strcmp(equiv[i][j].aname, aname2) == 0 );
270 *nname = strdup(equiv[i-1][0].nname);
276 static int analyze_noe_equivalent(const char *eq_fn,
277 t_atoms *atoms, int isize, atom_id *index,
279 atom_id *noe_index, t_noe_gr *noe_gr)
282 int i, j, anmil, anmjl, rnri, rnrj, gi, groupnr, neq;
283 char *anmi, *anmj, **nnm;
284 gmx_bool bMatch, bEquiv;
292 neq = read_equiv(eq_fn, &equiv);
295 dump_equiv(debug, neq, equiv);
305 for (i = 0; i < isize; i++)
307 if (equiv && i < isize-1)
309 /* check explicit list of equivalent atoms */
313 rnri = atoms->atom[index[i]].resind;
314 rnrj = atoms->atom[index[j]].resind;
316 is_equiv(neq, equiv, &nnm[i],
317 rnri, *atoms->resinfo[rnri].name, *atoms->atomname[index[i]],
318 rnrj, *atoms->resinfo[rnrj].name, *atoms->atomname[index[j]]);
319 if (nnm[i] && bEquiv)
321 nnm[j] = strdup(nnm[i]);
325 /* set index for matching atom */
326 noe_index[i] = groupnr;
327 /* skip matching atom */
331 while (bEquiv && i < isize-1);
339 /* look for triplets of consecutive atoms with name XX?,
340 X are any number of letters or digits and ? goes from 1 to 3
341 This is supposed to cover all CH3 groups and the like */
342 anmi = *atoms->atomname[index[i]];
343 anmil = strlen(anmi);
344 bMatch = i <= isize-3 && anmi[anmil-1] == '1';
347 for (j = 1; j < 3; j++)
349 anmj = *atoms->atomname[index[i+j]];
350 anmjl = strlen(anmj);
351 bMatch = bMatch && ( anmil == anmjl && anmj[anmjl-1] == Hnum[j] &&
352 strncmp(anmi, anmj, anmil-1) == 0 );
355 /* set index for this atom */
356 noe_index[i] = groupnr;
359 /* set index for next two matching atoms */
360 for (j = 1; j < 3; j++)
362 noe_index[i+j] = groupnr;
364 /* skip matching atoms */
373 /* make index without looking for equivalent atoms */
374 for (i = 0; i < isize; i++)
380 noe_index[isize] = groupnr;
385 for (i = 0; i < isize; i++)
387 rnri = atoms->atom[index[i]].resind;
388 fprintf(debug, "%s %s %d -> %s\n", *atoms->atomname[index[i]],
389 *atoms->resinfo[rnri].name, rnri, nnm[i] ? nnm[i] : "");
393 for (i = 0; i < isize; i++)
396 if (!noe_gr[gi].aname)
399 noe_gr[gi].anr = index[i];
402 noe_gr[gi].aname = strdup(nnm[i]);
406 noe_gr[gi].aname = strdup(*atoms->atomname[index[i]]);
407 if (noe_index[i] == noe_index[i+1])
409 noe_gr[gi].aname[strlen(noe_gr[gi].aname)-1] = '*';
412 noe_gr[gi].rnr = atoms->atom[index[i]].resind;
413 noe_gr[gi].rname = strdup(*atoms->resinfo[noe_gr[gi].rnr].name);
414 /* dump group definitions */
417 fprintf(debug, "%d %d %d %d %s %s %d\n", i, gi,
418 noe_gr[gi].ianr, noe_gr[gi].anr, noe_gr[gi].aname,
419 noe_gr[gi].rname, noe_gr[gi].rnr);
423 for (i = 0; i < isize; i++)
432 /* #define NSCALE 3 */
433 /* static char *noe_scale[NSCALE+1] = { */
434 /* "strong", "medium", "weak", "none" */
438 static char *noe2scale(real r3, real r6, real rmax)
441 static char buf[NSCALE+1];
443 /* r goes from 0 to rmax
444 NSCALE*r/rmax goes from 0 to NSCALE
445 NSCALE - NSCALE*r/rmax goes from NSCALE to 0 */
446 s3 = NSCALE - min(NSCALE, (int)(NSCALE*r3/rmax));
447 s6 = NSCALE - min(NSCALE, (int)(NSCALE*r6/rmax));
449 for (i = 0; i < s3; i++)
462 static void calc_noe(int isize, atom_id *noe_index,
463 real **dtot1_3, real **dtot1_6, int gnr, t_noe **noe)
467 /* make half matrix of noe-group distances from atom distances */
468 for (i = 0; i < isize; i++)
471 for (j = i; j < isize; j++)
475 noe[gi][gj].i_3 += pow(dtot1_3[i][j], -3);
476 noe[gi][gj].i_6 += pow(dtot1_6[i][j], -6);
481 for (i = 0; i < gnr; i++)
483 for (j = i+1; j < gnr; j++)
485 noe[i][j].r_3 = pow(noe[i][j].i_3/noe[i][j].nr, -1.0/3.0);
486 noe[i][j].r_6 = pow(noe[i][j].i_6/noe[i][j].nr, -1.0/6.0);
487 noe[j][i] = noe[i][j];
492 static void write_noe(FILE *fp, int gnr, t_noe **noe, t_noe_gr *noe_gr, real rmax)
495 real r3, r6, min3, min6;
496 char buf[10], b3[10], b6[10];
501 ";%4s %3s %4s %4s%3s %4s %4s %4s %4s%3s %5s %5s %8s %2s %2s %s\n",
502 "ianr", "anr", "anm", "rnm", "rnr", "ianr", "anr", "anm", "rnm", "rnr",
503 "1/r^3", "1/r^6", "intnsty", "Dr", "Da", "scale");
504 for (i = 0; i < gnr; i++)
507 for (j = i+1; j < gnr; j++)
512 min3 = min(r3, min3);
513 min6 = min(r6, min6);
514 if (r3 < rmax || r6 < rmax)
516 if (grj.rnr == gri.rnr)
518 sprintf(buf, "%2d", grj.anr-gri.anr);
526 sprintf(b3, "%-5.3f", r3);
534 sprintf(b6, "%-5.3f", r6);
541 "%4d %4d %4s %4s%3d %4d %4d %4s %4s%3d %5s %5s %8d %2d %2s %s\n",
542 gri.ianr+1, gri.anr+1, gri.aname, gri.rname, gri.rnr+1,
543 grj.ianr+1, grj.anr+1, grj.aname, grj.rname, grj.rnr+1,
544 b3, b6, (int)(noe[i][j].i_6+0.5), grj.rnr-gri.rnr, buf,
545 noe2scale(r3, r6, rmax));
549 #define MINI ((i == 3) ? min3 : min6)
550 for (i = 3; i <= 6; i += 3)
554 fprintf(stdout, "NOTE: no 1/r^%d averaged distances found below %g, "
555 "smallest was %g\n", i, rmax, MINI );
559 fprintf(stdout, "Smallest 1/r^%d averaged distance was %g\n", i, MINI );
565 static void calc_rms(int nind, int nframes,
566 real **dtot, real **dtot2,
567 real **rmsmat, real *rmsmax,
568 real **rmscmat, real *rmscmax,
569 real **meanmat, real *meanmax)
572 real mean, mean2, rms, rmsc;
573 /* N.B. dtot and dtot2 contain the total distance and the total squared
574 * distance respectively, BUT they return RMS and the scaled RMS resp.
581 for (i = 0; (i < nind-1); i++)
583 for (j = i+1; (j < nind); j++)
585 mean = dtot[i][j]/nframes;
586 mean2 = dtot2[i][j]/nframes;
587 rms = sqrt(max(0, mean2-mean*mean));
601 meanmat[i][j] = meanmat[j][i] = mean;
602 rmsmat[i][j] = rmsmat[j][i] = rms;
603 rmscmat[i][j] = rmscmat[j][i] = rmsc;
608 real rms_diff(int natom, real **d, real **d_r)
614 for (i = 0; (i < natom-1); i++)
616 for (j = i+1; (j < natom); j++)
618 r = d[i][j]-d_r[i][j];
622 r2 /= (natom*(natom-1))/2;
627 int gmx_rmsdist(int argc, char *argv[])
629 const char *desc[] = {
630 "[THISMODULE] computes the root mean square deviation of atom distances,",
631 "which has the advantage that no fit is needed like in standard RMS",
632 "deviation as computed by [gmx-rms].",
633 "The reference structure is taken from the structure file.",
634 "The RMSD at time t is calculated as the RMS",
635 "of the differences in distance between atom-pairs in the reference",
636 "structure and the structure at time t.[PAR]",
637 "[THISMODULE] can also produce matrices of the rms distances, rms distances",
638 "scaled with the mean distance and the mean distances and matrices with",
639 "NMR averaged distances (1/r^3 and 1/r^6 averaging). Finally, lists",
640 "of atom pairs with 1/r^3 and 1/r^6 averaged distance below the",
641 "maximum distance ([TT]-max[tt], which will default to 0.6 in this case)",
642 "can be generated, by default averaging over equivalent hydrogens",
643 "(all triplets of hydrogens named *[123]). Additionally a list of",
644 "equivalent atoms can be supplied ([TT]-equiv[tt]), each line containing",
645 "a set of equivalent atoms specified as residue number and name and",
646 "atom name; e.g.:[PAR]",
647 "[TT]HB* 3 SER HB1 3 SER HB2[tt][PAR]",
648 "Residue and atom names must exactly match those in the structure",
649 "file, including case. Specifying non-sequential atoms is undefined."
653 int natom, i, j, teller, gi, gj;
665 atom_id *index, *noe_index;
667 real **d_r, **d, **dtot, **dtot2, **mean, **rms, **rmsc, *resnr;
668 real **dtot1_3 = NULL, **dtot1_6 = NULL;
669 real rmsnow, meanmax, rmsmax, rmscmax;
671 t_noe_gr *noe_gr = NULL;
675 gmx_bool bRMS, bScale, bMean, bNOE, bNMR3, bNMR6, bNMR;
677 static int nlevels = 40;
678 static real scalemax = -1.0;
679 static gmx_bool bSumH = TRUE;
680 static gmx_bool bPBC = TRUE;
684 { "-nlevels", FALSE, etINT, {&nlevels},
685 "Discretize RMS in this number of levels" },
686 { "-max", FALSE, etREAL, {&scalemax},
687 "Maximum level in matrices" },
688 { "-sumh", FALSE, etBOOL, {&bSumH},
689 "Average distance over equivalent hydrogens" },
690 { "-pbc", FALSE, etBOOL, {&bPBC},
691 "Use periodic boundary conditions when computing distances" }
694 { efTRX, "-f", NULL, ffREAD },
695 { efTPS, NULL, NULL, ffREAD },
696 { efNDX, NULL, NULL, ffOPTRD },
697 { efDAT, "-equiv", "equiv", ffOPTRD },
698 { efXVG, NULL, "distrmsd", ffWRITE },
699 { efXPM, "-rms", "rmsdist", ffOPTWR },
700 { efXPM, "-scl", "rmsscale", ffOPTWR },
701 { efXPM, "-mean", "rmsmean", ffOPTWR },
702 { efXPM, "-nmr3", "nmr3", ffOPTWR },
703 { efXPM, "-nmr6", "nmr6", ffOPTWR },
704 { efDAT, "-noe", "noe", ffOPTWR },
706 #define NFILE asize(fnm)
708 if (!parse_common_args(&argc, argv, PCA_CAN_VIEW | PCA_CAN_TIME | PCA_BE_NICE,
709 NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL, &oenv))
714 bRMS = opt2bSet("-rms", NFILE, fnm);
715 bScale = opt2bSet("-scl", NFILE, fnm);
716 bMean = opt2bSet("-mean", NFILE, fnm);
717 bNOE = opt2bSet("-noe", NFILE, fnm);
718 bNMR3 = opt2bSet("-nmr3", NFILE, fnm);
719 bNMR6 = opt2bSet("-nmr6", NFILE, fnm);
720 bNMR = bNMR3 || bNMR6 || bNOE;
726 if (bNOE && scalemax < 0)
729 fprintf(stderr, "WARNING: using -noe without -max "
730 "makes no sense, setting -max to %g\n\n", scalemax);
733 /* get topology and index */
734 read_tps_conf(ftp2fn(efTPS, NFILE, fnm), buf, &top, &ePBC, &x, NULL, box, FALSE);
740 atoms = &(top.atoms);
742 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &isize, &index, &grpname);
744 /* initialize arrays */
750 snew(dtot1_3, isize);
751 snew(dtot1_6, isize);
758 for (i = 0; (i < isize); i++)
761 snew(dtot[i], isize);
762 snew(dtot2[i], isize);
765 snew(dtot1_3[i], isize);
766 snew(dtot1_6[i], isize);
768 snew(mean[i], isize);
770 snew(rmsc[i], isize);
776 calc_dist(isize, index, x, ePBC, box, d_r);
779 /*open output files*/
780 fp = xvgropen(ftp2fn(efXVG, NFILE, fnm), "RMS Deviation", "Time (ps)", "RMSD (nm)",
782 if (output_env_get_print_xvgr_codes(oenv))
784 fprintf(fp, "@ subtitle \"of distances between %s atoms\"\n", grpname);
788 natom = read_first_x(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &t, &x, box);
793 calc_dist_tot(isize, index, x, ePBC, box, d, dtot, dtot2, bNMR, dtot1_3, dtot1_6);
795 rmsnow = rms_diff(isize, d, d_r);
796 fprintf(fp, "%g %g\n", t, rmsnow);
799 while (read_next_x(oenv, status, &t, x, box));
800 fprintf(stderr, "\n");
806 calc_rms(isize, teller, dtot, dtot2, mean, &meanmax, rms, &rmsmax, rmsc, &rmscmax);
807 fprintf(stderr, "rmsmax = %g, rmscmax = %g\n", rmsmax, rmscmax);
811 calc_nmr(isize, teller, dtot1_3, dtot1_6, &max1_3, &max1_6);
825 /* make list of noe atom groups */
826 snew(noe_index, isize+1);
828 gnr = analyze_noe_equivalent(opt2fn_null("-equiv", NFILE, fnm),
829 atoms, isize, index, bSumH, noe_index, noe_gr);
830 fprintf(stdout, "Found %d non-equivalent atom-groups in %d atoms\n",
832 /* make half matrix of of noe-group distances from atom distances */
834 for (i = 0; i < gnr; i++)
838 calc_noe(isize, noe_index, dtot1_3, dtot1_6, gnr, noe);
841 rlo.r = 1.0, rlo.g = 1.0, rlo.b = 1.0;
842 rhi.r = 0.0, rhi.g = 0.0, rhi.b = 0.0;
846 write_xpm(opt2FILE("-rms", NFILE, fnm, "w"), 0,
847 "RMS of distance", "RMS (nm)", "Atom Index", "Atom Index",
848 isize, isize, resnr, resnr, rms, 0.0, rmsmax, rlo, rhi, &nlevels);
853 write_xpm(opt2FILE("-scl", NFILE, fnm, "w"), 0,
854 "Relative RMS", "RMS", "Atom Index", "Atom Index",
855 isize, isize, resnr, resnr, rmsc, 0.0, rmscmax, rlo, rhi, &nlevels);
860 write_xpm(opt2FILE("-mean", NFILE, fnm, "w"), 0,
861 "Mean Distance", "Distance (nm)", "Atom Index", "Atom Index",
862 isize, isize, resnr, resnr, mean, 0.0, meanmax, rlo, rhi, &nlevels);
867 write_xpm(opt2FILE("-nmr3", NFILE, fnm, "w"), 0, "1/r^3 averaged distances",
868 "Distance (nm)", "Atom Index", "Atom Index",
869 isize, isize, resnr, resnr, dtot1_3, 0.0, max1_3, rlo, rhi, &nlevels);
873 write_xpm(opt2FILE("-nmr6", NFILE, fnm, "w"), 0, "1/r^6 averaged distances",
874 "Distance (nm)", "Atom Index", "Atom Index",
875 isize, isize, resnr, resnr, dtot1_6, 0.0, max1_6, rlo, rhi, &nlevels);
880 write_noe(opt2FILE("-noe", NFILE, fnm, "w"), gnr, noe, noe_gr, scalemax);
883 do_view(oenv, ftp2fn(efXVG, NFILE, fnm), NULL);