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
45 #include "gromacs/commandline/pargs.h"
49 #include "gmx_fatal.h"
50 #include "gromacs/fileio/futil.h"
54 #include "gromacs/fileio/matio.h"
55 #include "gromacs/fileio/tpxio.h"
56 #include "gromacs/fileio/trxio.h"
62 static void norm_princ(t_atoms *atoms, int isize, atom_id *index, int natoms,
68 /* equalize principal components: */
69 /* orient principal axes, get principal components */
70 orient_princ(atoms, isize, index, natoms, x, NULL, princ);
72 /* calc our own principal components */
74 for (m = 0; m < DIM; m++)
76 for (i = 0; i < isize; i++)
78 vec[m] += sqr(x[index[i]][m]);
80 vec[m] = sqrt(vec[m] / isize);
81 /* calculate scaling constants */
82 vec[m] = 1 / (sqrt(3) * vec[m]);
85 /* scale coordinates */
86 for (i = 0; i < natoms; i++)
88 for (m = 0; m < DIM; m++)
95 int gmx_rms(int argc, char *argv[])
99 "[THISMODULE] compares two structures by computing the root mean square",
100 "deviation (RMSD), the size-independent [GRK]rho[grk] similarity parameter",
101 "([TT]rho[tt]) or the scaled [GRK]rho[grk] ([TT]rhosc[tt]), ",
102 "see Maiorov & Crippen, Proteins [BB]22[bb], 273 (1995).",
103 "This is selected by [TT]-what[tt].[PAR]"
105 "Each structure from a trajectory ([TT]-f[tt]) is compared to a",
106 "reference structure. The reference structure",
107 "is taken from the structure file ([TT]-s[tt]).[PAR]",
109 "With option [TT]-mir[tt] also a comparison with the mirror image of",
110 "the reference structure is calculated.",
111 "This is useful as a reference for 'significant' values, see",
112 "Maiorov & Crippen, Proteins [BB]22[bb], 273 (1995).[PAR]",
114 "Option [TT]-prev[tt] produces the comparison with a previous frame",
115 "the specified number of frames ago.[PAR]",
117 "Option [TT]-m[tt] produces a matrix in [TT].xpm[tt] format of",
118 "comparison values of each structure in the trajectory with respect to",
119 "each other structure. This file can be visualized with for instance",
120 "[TT]xv[tt] and can be converted to postscript with [gmx-xpm2ps].[PAR]",
122 "Option [TT]-fit[tt] controls the least-squares fitting of",
123 "the structures on top of each other: complete fit (rotation and",
124 "translation), translation only, or no fitting at all.[PAR]",
126 "Option [TT]-mw[tt] controls whether mass weighting is done or not.",
127 "If you select the option (default) and ",
128 "supply a valid [TT].tpr[tt] file masses will be taken from there, ",
129 "otherwise the masses will be deduced from the [TT]atommass.dat[tt] file in",
130 "[TT]GMXLIB[tt]. This is fine for proteins, but not",
131 "necessarily for other molecules. A default mass of 12.011 amu (carbon)",
132 "is assigned to unknown atoms. You can check whether this happend by",
133 "turning on the [TT]-debug[tt] flag and inspecting the log file.[PAR]",
135 "With [TT]-f2[tt], the 'other structures' are taken from a second",
136 "trajectory, this generates a comparison matrix of one trajectory",
137 "versus the other.[PAR]",
139 "Option [TT]-bin[tt] does a binary dump of the comparison matrix.[PAR]",
141 "Option [TT]-bm[tt] produces a matrix of average bond angle deviations",
142 "analogously to the [TT]-m[tt] option. Only bonds between atoms in the",
143 "comparison group are considered."
145 static gmx_bool bPBC = TRUE, bFitAll = TRUE, bSplit = FALSE;
146 static gmx_bool bDeltaLog = FALSE;
147 static int prev = 0, freq = 1, freq2 = 1, nlevels = 80, avl = 0;
148 static real rmsd_user_max = -1, rmsd_user_min = -1, bond_user_max = -1,
149 bond_user_min = -1, delta_maxy = 0.0;
150 /* strings and things for selecting difference method */
153 ewSel, ewRMSD, ewRho, ewRhoSc, ewNR
156 const char *what[ewNR + 1] =
157 { NULL, "rmsd", "rho", "rhosc", NULL };
158 const char *whatname[ewNR] =
159 { NULL, "RMSD", "Rho", "Rho sc" };
160 const char *whatlabel[ewNR] =
161 { NULL, "RMSD (nm)", "Rho", "Rho sc" };
162 const char *whatxvgname[ewNR] =
163 { NULL, "RMSD", "\\8r\\4", "\\8r\\4\\ssc\\N" };
164 const char *whatxvglabel[ewNR] =
165 { NULL, "RMSD (nm)", "\\8r\\4", "\\8r\\4\\ssc\\N" };
166 /* strings and things for fitting methods */
169 efSel, efFit, efReset, efNone, efNR
172 const char *fit[efNR + 1] =
173 { NULL, "rot+trans", "translation", "none", NULL };
174 const char *fitgraphlabel[efNR + 1] =
175 { NULL, "lsq fit", "translational fit", "no fit" };
177 static gmx_bool bMassWeighted = TRUE;
180 { "-what", FALSE, etENUM,
181 { what }, "Structural difference measure" },
182 { "-pbc", FALSE, etBOOL,
183 { &bPBC }, "PBC check" },
184 { "-fit", FALSE, etENUM,
185 { fit }, "Fit to reference structure" },
186 { "-prev", FALSE, etINT,
187 { &prev }, "Compare with previous frame" },
188 { "-split", FALSE, etBOOL,
189 { &bSplit }, "Split graph where time is zero" },
190 { "-fitall", FALSE, etBOOL,
191 { &bFitAll }, "HIDDENFit all pairs of structures in matrix" },
192 { "-skip", FALSE, etINT,
193 { &freq }, "Only write every nr-th frame to matrix" },
194 { "-skip2", FALSE, etINT,
195 { &freq2 }, "Only write every nr-th frame to matrix" },
196 { "-max", FALSE, etREAL,
197 { &rmsd_user_max }, "Maximum level in comparison matrix" },
198 { "-min", FALSE, etREAL,
199 { &rmsd_user_min }, "Minimum level in comparison matrix" },
200 { "-bmax", FALSE, etREAL,
201 { &bond_user_max }, "Maximum level in bond angle matrix" },
202 { "-bmin", FALSE, etREAL,
203 { &bond_user_min }, "Minimum level in bond angle matrix" },
204 { "-mw", FALSE, etBOOL,
205 { &bMassWeighted }, "Use mass weighting for superposition" },
206 { "-nlevels", FALSE, etINT,
207 { &nlevels }, "Number of levels in the matrices" },
208 { "-ng", FALSE, etINT,
209 { &nrms }, "Number of groups to compute RMS between" },
210 { "-dlog", FALSE, etBOOL,
212 "HIDDENUse a log x-axis in the delta t matrix" },
213 { "-dmax", FALSE, etREAL,
214 { &delta_maxy }, "HIDDENMaximum level in delta matrix" },
215 { "-aver", FALSE, etINT,
217 "HIDDENAverage over this distance in the RMSD matrix" }
219 int natoms_trx, natoms_trx2, natoms;
220 int i, j, k, m, teller, teller2, tel_mat, tel_mat2;
222 int maxframe = NFRAME, maxframe2 = NFRAME;
223 real t, *w_rls, *w_rms, *w_rls_m = NULL, *w_rms_m = NULL;
224 gmx_bool bNorm, bAv, bFreq2, bFile2, bMat, bBond, bDelta, bMirror, bMass;
225 gmx_bool bFit, bReset;
228 t_iatom *iatom = NULL;
231 rvec *x, *xp, *xm = NULL, **mat_x = NULL, **mat_x2, *mat_x2_j = NULL, vec1,
234 char buf[256], buf2[256];
237 real rlstot = 0, **rls, **rlsm = NULL, *time, *time2, *rlsnorm = NULL,
238 **rmsd_mat = NULL, **bond_mat = NULL, *axis, *axis2, *del_xaxis,
239 *del_yaxis, rmsd_max, rmsd_min, rmsd_avg, bond_max, bond_min, ang;
240 real **rmsdav_mat = NULL, av_tot, weight, weight_tot;
241 real **delta = NULL, delta_max, delta_scalex = 0, delta_scaley = 0,
243 int delta_xsize = 0, del_lev = 100, mx, my, abs_my;
244 gmx_bool bA1, bA2, bPrev, bTop, *bInMat = NULL;
245 int ifit, *irms, ibond = 0, *ind_bond1 = NULL, *ind_bond2 = NULL, n_ind_m =
247 atom_id *ind_fit, **ind_rms, *ind_m = NULL, *rev_ind_m = NULL, *ind_rms_m =
249 char *gn_fit, **gn_rms;
252 gmx_rmpbc_t gpbc = NULL;
256 { efTPS, NULL, NULL, ffREAD },
257 { efTRX, "-f", NULL, ffREAD },
258 { efTRX, "-f2", NULL, ffOPTRD },
259 { efNDX, NULL, NULL, ffOPTRD },
260 { efXVG, NULL, "rmsd", ffWRITE },
261 { efXVG, "-mir", "rmsdmir", ffOPTWR },
262 { efXVG, "-a", "avgrp", ffOPTWR },
263 { efXVG, "-dist", "rmsd-dist", ffOPTWR },
264 { efXPM, "-m", "rmsd", ffOPTWR },
265 { efDAT, "-bin", "rmsd", ffOPTWR },
266 { efXPM, "-bm", "bond", ffOPTWR }
268 #define NFILE asize(fnm)
270 if (!parse_common_args(&argc, argv, PCA_CAN_TIME | PCA_TIME_UNIT | PCA_CAN_VIEW
271 | PCA_BE_NICE, NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL,
276 /* parse enumerated options: */
278 if (ewhat == ewRho || ewhat == ewRhoSc)
280 please_cite(stdout, "Maiorov95");
283 bFit = efit == efFit;
284 bReset = efit == efReset;
287 bReset = TRUE; /* for fit, reset *must* be set */
294 /* mark active cmdline options */
295 bMirror = opt2bSet("-mir", NFILE, fnm); /* calc RMSD vs mirror of ref. */
296 bFile2 = opt2bSet("-f2", NFILE, fnm);
297 bMat = opt2bSet("-m", NFILE, fnm);
298 bBond = opt2bSet("-bm", NFILE, fnm);
299 bDelta = (delta_maxy > 0); /* calculate rmsd vs delta t matrix from *
300 * your RMSD matrix (hidden option */
301 bNorm = opt2bSet("-a", NFILE, fnm);
302 bFreq2 = opt2parg_bSet("-skip2", asize(pa), pa);
305 fprintf(stderr, "The number of frames to skip is <= 0. "
306 "Writing out all frames.\n\n");
313 else if (bFile2 && freq2 <= 0)
316 "The number of frames to skip in second trajectory is <= 0.\n"
317 " Writing out all frames.\n\n");
327 fprintf(stderr, "WARNING: option -skip also applies to -prev\n");
331 if (bFile2 && !bMat && !bBond)
335 "WARNING: second trajectory (-f2) useless when not calculating matrix (-m/-bm),\n"
336 " will not read from %s\n", opt2fn("-f2", NFILE,
347 "WARNING: second trajectory (-f2) useless when making delta matrix,\n"
348 " will not read from %s\n", opt2fn("-f2",
354 bTop = read_tps_conf(ftp2fn(efTPS, NFILE, fnm), buf, &top, &ePBC, &xp,
356 snew(w_rls, top.atoms.nr);
357 snew(w_rms, top.atoms.nr);
362 "WARNING: Need a run input file for bond angle matrix,\n"
363 " will not calculate bond angle matrix.\n");
369 fprintf(stderr, "Select group for %s fit\n", bFit ? "least squares"
371 get_index(&(top.atoms), ftp2fn_null(efNDX, NFILE, fnm), 1, &ifit,
381 if (bFit && ifit < 3)
383 gmx_fatal(FARGS, "Need >= 3 points to fit!\n" );
387 for (i = 0; i < ifit; i++)
391 w_rls[ind_fit[i]] = top.atoms.atom[ind_fit[i]].m;
395 w_rls[ind_fit[i]] = 1;
397 bMass = bMass || (top.atoms.atom[ind_fit[i]].m != 0);
401 fprintf(stderr, "All masses in the fit group are 0, using masses of 1\n");
402 for (i = 0; i < ifit; i++)
404 w_rls[ind_fit[i]] = 1;
418 fprintf(stderr, "Select group%s for %s calculation\n",
419 (nrms > 1) ? "s" : "", whatname[ewhat]);
420 get_index(&(top.atoms), ftp2fn_null(efNDX, NFILE, fnm),
421 nrms, irms, ind_rms, gn_rms);
425 snew(rlsnorm, irms[0]);
428 for (j = 0; j < nrms; j++)
430 snew(rls[j], maxframe);
435 for (j = 0; j < nrms; j++)
437 snew(rlsm[j], maxframe);
440 snew(time, maxframe);
441 for (j = 0; j < nrms; j++)
444 for (i = 0; i < irms[j]; i++)
448 w_rms[ind_rms[j][i]] = top.atoms.atom[ind_rms[j][i]].m;
452 w_rms[ind_rms[j][i]] = 1.0;
454 bMass = bMass || (top.atoms.atom[ind_rms[j][i]].m != 0);
458 fprintf(stderr, "All masses in group %d are 0, using masses of 1\n", j);
459 for (i = 0; i < irms[j]; i++)
461 w_rms[ind_rms[j][i]] = 1;
465 /* Prepare reference frame */
468 gpbc = gmx_rmpbc_init(&top.idef, ePBC, top.atoms.nr);
469 gmx_rmpbc(gpbc, top.atoms.nr, box, xp);
473 reset_x(ifit, ind_fit, top.atoms.nr, NULL, xp, w_rls);
477 /* generate reference structure mirror image: */
478 snew(xm, top.atoms.nr);
479 for (i = 0; i < top.atoms.nr; i++)
481 copy_rvec(xp[i], xm[i]);
482 xm[i][XX] = -xm[i][XX];
485 if (ewhat == ewRhoSc)
487 norm_princ(&top.atoms, ifit, ind_fit, top.atoms.nr, xp);
490 /* read first frame */
491 natoms_trx = read_first_x(oenv, &status, opt2fn("-f", NFILE, fnm), &t, &x, box);
492 if (natoms_trx != top.atoms.nr)
495 "\nWARNING: topology has %d atoms, whereas trajectory has %d\n",
496 top.atoms.nr, natoms_trx);
498 natoms = min(top.atoms.nr, natoms_trx);
499 if (bMat || bBond || bPrev)
505 /* With -prev we use all atoms for simplicity */
510 /* Check which atoms we need (fit/rms) */
511 snew(bInMat, natoms);
512 for (i = 0; i < ifit; i++)
514 bInMat[ind_fit[i]] = TRUE;
517 for (i = 0; i < irms[0]; i++)
519 if (!bInMat[ind_rms[0][i]])
521 bInMat[ind_rms[0][i]] = TRUE;
526 /* Make an index of needed atoms */
527 snew(ind_m, n_ind_m);
528 snew(rev_ind_m, natoms);
530 for (i = 0; i < natoms; i++)
532 if (bPrev || bInMat[i])
539 snew(w_rls_m, n_ind_m);
540 snew(ind_rms_m, irms[0]);
541 snew(w_rms_m, n_ind_m);
542 for (i = 0; i < ifit; i++)
544 w_rls_m[rev_ind_m[ind_fit[i]]] = w_rls[ind_fit[i]];
546 for (i = 0; i < irms[0]; i++)
548 ind_rms_m[i] = rev_ind_m[ind_rms[0][i]];
549 w_rms_m[ind_rms_m[i]] = w_rms[ind_rms[0][i]];
557 for (k = 0; k < F_NRE; k++)
561 iatom = top.idef.il[k].iatoms;
562 ncons += top.idef.il[k].nr/3;
565 fprintf(stderr, "Found %d bonds in topology\n", ncons);
566 snew(ind_bond1, ncons);
567 snew(ind_bond2, ncons);
569 for (k = 0; k < F_NRE; k++)
573 iatom = top.idef.il[k].iatoms;
574 ncons = top.idef.il[k].nr/3;
575 for (i = 0; i < ncons; i++)
579 for (j = 0; j < irms[0]; j++)
581 if (iatom[3*i+1] == ind_rms[0][j])
585 if (iatom[3*i+2] == ind_rms[0][j])
592 ind_bond1[ibond] = rev_ind_m[iatom[3*i+1]];
593 ind_bond2[ibond] = rev_ind_m[iatom[3*i+2]];
599 fprintf(stderr, "Using %d bonds for bond angle matrix\n", ibond);
602 gmx_fatal(FARGS, "0 bonds found");
606 /* start looping over frames: */
613 gmx_rmpbc(gpbc, natoms, box, x);
618 reset_x(ifit, ind_fit, natoms, NULL, x, w_rls);
620 if (ewhat == ewRhoSc)
622 norm_princ(&top.atoms, ifit, ind_fit, natoms, x);
627 /*do the least squares fit to original structure*/
628 do_fit(natoms, w_rls, xp, x);
631 if (teller % freq == 0)
633 /* keep frame for matrix calculation */
634 if (bMat || bBond || bPrev)
636 if (tel_mat >= NFRAME)
638 srenew(mat_x, tel_mat+1);
640 snew(mat_x[tel_mat], n_ind_m);
641 for (i = 0; i < n_ind_m; i++)
643 copy_rvec(x[ind_m[i]], mat_x[tel_mat][i]);
649 /*calculate energy of root_least_squares*/
657 for (i = 0; i < n_ind_m; i++)
659 copy_rvec(mat_x[j][i], xp[ind_m[i]]);
663 reset_x(ifit, ind_fit, natoms, NULL, xp, w_rls);
667 do_fit(natoms, w_rls, x, xp);
670 for (j = 0; (j < nrms); j++)
673 calc_similar_ind(ewhat != ewRMSD, irms[j], ind_rms[j], w_rms, x, xp);
677 for (j = 0; (j < irms[0]); j++)
680 calc_similar_ind(ewhat != ewRMSD, 1, &(ind_rms[0][j]), w_rms, x, xp);
688 /*do the least squares fit to mirror of original structure*/
689 do_fit(natoms, w_rls, xm, x);
692 for (j = 0; j < nrms; j++)
695 calc_similar_ind(ewhat != ewRMSD, irms[j], ind_rms[j], w_rms, x, xm);
698 time[teller] = output_env_conv_time(oenv, t);
701 if (teller >= maxframe)
704 srenew(time, maxframe);
705 for (j = 0; (j < nrms); j++)
707 srenew(rls[j], maxframe);
711 for (j = 0; (j < nrms); j++)
713 srenew(rlsm[j], maxframe);
718 while (read_next_x(oenv, status, &t, x, box));
723 snew(time2, maxframe2);
725 fprintf(stderr, "\nWill read second trajectory file\n");
726 snew(mat_x2, NFRAME);
728 read_first_x(oenv, &status, opt2fn("-f2", NFILE, fnm), &t, &x, box);
729 if (natoms_trx2 != natoms_trx)
732 "Second trajectory (%d atoms) does not match the first one"
733 " (%d atoms)", natoms_trx2, natoms_trx);
741 gmx_rmpbc(gpbc, natoms, box, x);
746 reset_x(ifit, ind_fit, natoms, NULL, x, w_rls);
748 if (ewhat == ewRhoSc)
750 norm_princ(&top.atoms, ifit, ind_fit, natoms, x);
755 /*do the least squares fit to original structure*/
756 do_fit(natoms, w_rls, xp, x);
759 if (teller2 % freq2 == 0)
761 /* keep frame for matrix calculation */
764 if (tel_mat2 >= NFRAME)
766 srenew(mat_x2, tel_mat2+1);
768 snew(mat_x2[tel_mat2], n_ind_m);
769 for (i = 0; i < n_ind_m; i++)
771 copy_rvec(x[ind_m[i]], mat_x2[tel_mat2][i]);
777 time2[teller2] = output_env_conv_time(oenv, t);
780 if (teller2 >= maxframe2)
783 srenew(time2, maxframe2);
786 while (read_next_x(oenv, status, &t, x, box));
796 gmx_rmpbc_done(gpbc);
800 /* calculate RMS matrix */
801 fprintf(stderr, "\n");
804 fprintf(stderr, "Building %s matrix, %dx%d elements\n",
805 whatname[ewhat], tel_mat, tel_mat2);
806 snew(rmsd_mat, tel_mat);
810 fprintf(stderr, "Building bond angle matrix, %dx%d elements\n",
812 snew(bond_mat, tel_mat);
815 snew(axis2, tel_mat2);
828 for (j = 0; j < tel_mat2; j++)
830 axis2[j] = time2[freq2*j];
836 delta_scalex = 8.0/log(2.0);
837 delta_xsize = (int)(log(tel_mat/2)*delta_scalex+0.5)+1;
841 delta_xsize = tel_mat/2;
843 delta_scaley = 1.0/delta_maxy;
844 snew(delta, delta_xsize);
845 for (j = 0; j < delta_xsize; j++)
847 snew(delta[j], del_lev+1);
851 snew(rmsdav_mat, tel_mat);
852 for (j = 0; j < tel_mat; j++)
854 snew(rmsdav_mat[j], tel_mat);
861 snew(mat_x2_j, natoms);
863 for (i = 0; i < tel_mat; i++)
865 axis[i] = time[freq*i];
866 fprintf(stderr, "\r element %5d; time %5.2f ", i, axis[i]);
869 snew(rmsd_mat[i], tel_mat2);
873 snew(bond_mat[i], tel_mat2);
875 for (j = 0; j < tel_mat2; j++)
879 for (k = 0; k < n_ind_m; k++)
881 copy_rvec(mat_x2[j][k], mat_x2_j[k]);
883 do_fit(n_ind_m, w_rls_m, mat_x[i], mat_x2_j);
887 mat_x2_j = mat_x2[j];
891 if (bFile2 || (i < j))
894 calc_similar_ind(ewhat != ewRMSD, irms[0], ind_rms_m,
895 w_rms_m, mat_x[i], mat_x2_j);
896 if (rmsd_mat[i][j] > rmsd_max)
898 rmsd_max = rmsd_mat[i][j];
900 if (rmsd_mat[i][j] < rmsd_min)
902 rmsd_min = rmsd_mat[i][j];
904 rmsd_avg += rmsd_mat[i][j];
908 rmsd_mat[i][j] = rmsd_mat[j][i];
913 if (bFile2 || (i <= j))
916 for (m = 0; m < ibond; m++)
918 rvec_sub(mat_x[i][ind_bond1[m]], mat_x[i][ind_bond2[m]], vec1);
919 rvec_sub(mat_x2_j[ind_bond1[m]], mat_x2_j[ind_bond2[m]], vec2);
920 ang += acos(cos_angle(vec1, vec2));
922 bond_mat[i][j] = ang*180.0/(M_PI*ibond);
923 if (bond_mat[i][j] > bond_max)
925 bond_max = bond_mat[i][j];
927 if (bond_mat[i][j] < bond_min)
929 bond_min = bond_mat[i][j];
934 bond_mat[i][j] = bond_mat[j][i];
941 rmsd_avg /= tel_mat*tel_mat2;
945 rmsd_avg /= tel_mat*(tel_mat - 1)/2;
947 if (bMat && (avl > 0))
952 for (j = 0; j < tel_mat-1; j++)
954 for (i = j+1; i < tel_mat; i++)
958 for (my = -avl; my <= avl; my++)
960 if ((j+my >= 0) && (j+my < tel_mat))
963 for (mx = -avl; mx <= avl; mx++)
965 if ((i+mx >= 0) && (i+mx < tel_mat))
967 weight = (real)(avl+1-max(abs(mx), abs_my));
968 av_tot += weight*rmsd_mat[i+mx][j+my];
969 weight_tot += weight;
974 rmsdav_mat[i][j] = av_tot/weight_tot;
975 rmsdav_mat[j][i] = rmsdav_mat[i][j];
976 if (rmsdav_mat[i][j] > rmsd_max)
978 rmsd_max = rmsdav_mat[i][j];
982 rmsd_mat = rmsdav_mat;
987 fprintf(stderr, "\n%s: Min %f, Max %f, Avg %f\n",
988 whatname[ewhat], rmsd_min, rmsd_max, rmsd_avg);
989 rlo.r = 1; rlo.g = 1; rlo.b = 1;
990 rhi.r = 0; rhi.g = 0; rhi.b = 0;
991 if (rmsd_user_max != -1)
993 rmsd_max = rmsd_user_max;
995 if (rmsd_user_min != -1)
997 rmsd_min = rmsd_user_min;
999 if ((rmsd_user_max != -1) || (rmsd_user_min != -1))
1001 fprintf(stderr, "Min and Max value set to resp. %f and %f\n",
1002 rmsd_min, rmsd_max);
1004 sprintf(buf, "%s %s matrix", gn_rms[0], whatname[ewhat]);
1005 write_xpm(opt2FILE("-m", NFILE, fnm, "w"), 0, buf, whatlabel[ewhat],
1006 output_env_get_time_label(oenv), output_env_get_time_label(oenv), tel_mat, tel_mat2,
1007 axis, axis2, rmsd_mat, rmsd_min, rmsd_max, rlo, rhi, &nlevels);
1008 /* Print the distribution of RMSD values */
1009 if (opt2bSet("-dist", NFILE, fnm))
1011 low_rmsd_dist(opt2fn("-dist", NFILE, fnm), rmsd_max, tel_mat, rmsd_mat, oenv);
1016 snew(delta_tot, delta_xsize);
1017 for (j = 0; j < tel_mat-1; j++)
1019 for (i = j+1; i < tel_mat; i++)
1026 mx = (int)(log(mx)*delta_scalex+0.5);
1028 my = (int)(rmsd_mat[i][j]*delta_scaley*del_lev+0.5);
1029 delta_tot[mx] += 1.0;
1030 if ((rmsd_mat[i][j] >= 0) && (rmsd_mat[i][j] <= delta_maxy))
1032 delta[mx][my] += 1.0;
1038 for (i = 0; i < delta_xsize; i++)
1040 if (delta_tot[i] > 0.0)
1042 delta_tot[i] = 1.0/delta_tot[i];
1043 for (j = 0; j <= del_lev; j++)
1045 delta[i][j] *= delta_tot[i];
1046 if (delta[i][j] > delta_max)
1048 delta_max = delta[i][j];
1053 fprintf(stderr, "Maximum in delta matrix: %f\n", delta_max);
1054 snew(del_xaxis, delta_xsize);
1055 snew(del_yaxis, del_lev+1);
1056 for (i = 0; i < delta_xsize; i++)
1058 del_xaxis[i] = axis[i]-axis[0];
1060 for (i = 0; i < del_lev+1; i++)
1062 del_yaxis[i] = delta_maxy*i/del_lev;
1064 sprintf(buf, "%s %s vs. delta t", gn_rms[0], whatname[ewhat]);
1065 fp = ffopen("delta.xpm", "w");
1066 write_xpm(fp, 0, buf, "density", output_env_get_time_label(oenv), whatlabel[ewhat],
1067 delta_xsize, del_lev+1, del_xaxis, del_yaxis,
1068 delta, 0.0, delta_max, rlo, rhi, &nlevels);
1071 if (opt2bSet("-bin", NFILE, fnm))
1073 /* NB: File must be binary if we use fwrite */
1074 fp = ftp2FILE(efDAT, NFILE, fnm, "wb");
1075 for (i = 0; i < tel_mat; i++)
1077 if (fwrite(rmsd_mat[i], sizeof(**rmsd_mat), tel_mat2, fp) != tel_mat2)
1079 gmx_fatal(FARGS, "Error writing to output file");
1087 fprintf(stderr, "\nMin. angle: %f, Max. angle: %f\n", bond_min, bond_max);
1088 if (bond_user_max != -1)
1090 bond_max = bond_user_max;
1092 if (bond_user_min != -1)
1094 bond_min = bond_user_min;
1096 if ((bond_user_max != -1) || (bond_user_min != -1))
1098 fprintf(stderr, "Bond angle Min and Max set to:\n"
1099 "Min. angle: %f, Max. angle: %f\n", bond_min, bond_max);
1101 rlo.r = 1; rlo.g = 1; rlo.b = 1;
1102 rhi.r = 0; rhi.g = 0; rhi.b = 0;
1103 sprintf(buf, "%s av. bond angle deviation", gn_rms[0]);
1104 write_xpm(opt2FILE("-bm", NFILE, fnm, "w"), 0, buf, "degrees",
1105 output_env_get_time_label(oenv), output_env_get_time_label(oenv), tel_mat, tel_mat2,
1106 axis, axis2, bond_mat, bond_min, bond_max, rlo, rhi, &nlevels);
1110 bAv = opt2bSet("-a", NFILE, fnm);
1112 /* Write the RMSD's to file */
1115 sprintf(buf, "%s", whatxvgname[ewhat]);
1119 sprintf(buf, "%s with frame %g %s ago", whatxvgname[ewhat],
1120 time[prev*freq]-time[0], output_env_get_time_label(oenv));
1122 fp = xvgropen(opt2fn("-o", NFILE, fnm), buf, output_env_get_xvgr_tlabel(oenv),
1123 whatxvglabel[ewhat], oenv);
1124 if (output_env_get_print_xvgr_codes(oenv))
1126 fprintf(fp, "@ subtitle \"%s%s after %s%s%s\"\n",
1127 (nrms == 1) ? "" : "of ", gn_rms[0], fitgraphlabel[efit],
1128 bFit ? " to " : "", bFit ? gn_fit : "");
1132 xvgr_legend(fp, nrms, (const char**)gn_rms, oenv);
1134 for (i = 0; (i < teller); i++)
1136 if (bSplit && i > 0 &&
1137 abs(time[bPrev ? freq*i : i]/output_env_get_time_factor(oenv)) < 1e-5)
1141 fprintf(fp, "%12.7f", time[bPrev ? freq*i : i]);
1142 for (j = 0; (j < nrms); j++)
1144 fprintf(fp, " %12.7f", rls[j][i]);
1147 rlstot += rls[j][i];
1156 /* Write the mirror RMSD's to file */
1157 sprintf(buf, "%s with Mirror", whatxvgname[ewhat]);
1158 sprintf(buf2, "Mirror %s", whatxvglabel[ewhat]);
1159 fp = xvgropen(opt2fn("-mir", NFILE, fnm), buf, output_env_get_xvgr_tlabel(oenv),
1163 if (output_env_get_print_xvgr_codes(oenv))
1165 fprintf(fp, "@ subtitle \"of %s after lsq fit to mirror of %s\"\n",
1171 if (output_env_get_print_xvgr_codes(oenv))
1173 fprintf(fp, "@ subtitle \"after lsq fit to mirror %s\"\n", gn_fit);
1175 xvgr_legend(fp, nrms, (const char**)gn_rms, oenv);
1177 for (i = 0; (i < teller); i++)
1179 if (bSplit && i > 0 && abs(time[i]) < 1e-5)
1183 fprintf(fp, "%12.7f", time[i]);
1184 for (j = 0; (j < nrms); j++)
1186 fprintf(fp, " %12.7f", rlsm[j][i]);
1195 sprintf(buf, "Average %s", whatxvgname[ewhat]);
1196 sprintf(buf2, "Average %s", whatxvglabel[ewhat]);
1197 fp = xvgropen(opt2fn("-a", NFILE, fnm), buf, "Residue", buf2, oenv);
1198 for (j = 0; (j < nrms); j++)
1200 fprintf(fp, "%10d %10g\n", j, rlstot/teller);
1207 fp = xvgropen("aver.xvg", gn_rms[0], "Residue", whatxvglabel[ewhat], oenv);
1208 for (j = 0; (j < irms[0]); j++)
1210 fprintf(fp, "%10d %10g\n", j, rlsnorm[j]/teller);
1214 do_view(oenv, opt2fn_null("-a", NFILE, fnm), "-graphtype bar");
1215 do_view(oenv, opt2fn("-o", NFILE, fnm), NULL);
1216 do_view(oenv, opt2fn_null("-mir", NFILE, fnm), NULL);
1217 do_view(oenv, opt2fn_null("-m", NFILE, fnm), NULL);
1218 do_view(oenv, opt2fn_null("-bm", NFILE, fnm), NULL);
1219 do_view(oenv, opt2fn_null("-dist", NFILE, fnm), NULL);