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
49 #include "gmx_fatal.h"
61 static void norm_princ(t_atoms *atoms, int isize, atom_id *index, int natoms,
67 /* equalize principal components: */
68 /* orient principal axes, get principal components */
69 orient_princ(atoms, isize, index, natoms, x, NULL, princ);
71 /* calc our own principal components */
73 for (m = 0; m < DIM; m++)
75 for (i = 0; i < isize; i++)
77 vec[m] += sqr(x[index[i]][m]);
79 vec[m] = sqrt(vec[m] / isize);
80 /* calculate scaling constants */
81 vec[m] = 1 / (sqrt(3) * vec[m]);
84 /* scale coordinates */
85 for (i = 0; i < natoms; i++)
87 for (m = 0; m < DIM; m++)
94 int gmx_rms(int argc, char *argv[])
99 "[TT]g_rms[tt] 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 [TT]xpm2ps[tt].[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 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,
273 /* parse enumerated options: */
275 if (ewhat == ewRho || ewhat == ewRhoSc)
277 please_cite(stdout, "Maiorov95");
280 bFit = efit == efFit;
281 bReset = efit == efReset;
284 bReset = TRUE; /* for fit, reset *must* be set */
291 /* mark active cmdline options */
292 bMirror = opt2bSet("-mir", NFILE, fnm); /* calc RMSD vs mirror of ref. */
293 bFile2 = opt2bSet("-f2", NFILE, fnm);
294 bMat = opt2bSet("-m", NFILE, fnm);
295 bBond = opt2bSet("-bm", NFILE, fnm);
296 bDelta = (delta_maxy > 0); /* calculate rmsd vs delta t matrix from *
297 * your RMSD matrix (hidden option */
298 bNorm = opt2bSet("-a", NFILE, fnm);
299 bFreq2 = opt2parg_bSet("-skip2", asize(pa), pa);
302 fprintf(stderr, "The number of frames to skip is <= 0. "
303 "Writing out all frames.\n\n");
310 else if (bFile2 && freq2 <= 0)
313 "The number of frames to skip in second trajectory is <= 0.\n"
314 " Writing out all frames.\n\n");
324 fprintf(stderr, "WARNING: option -skip also applies to -prev\n");
328 if (bFile2 && !bMat && !bBond)
332 "WARNING: second trajectory (-f2) useless when not calculating matrix (-m/-bm),\n"
333 " will not read from %s\n", opt2fn("-f2", NFILE,
344 "WARNING: second trajectory (-f2) useless when making delta matrix,\n"
345 " will not read from %s\n", opt2fn("-f2",
351 bTop = read_tps_conf(ftp2fn(efTPS, NFILE, fnm), buf, &top, &ePBC, &xp,
353 snew(w_rls, top.atoms.nr);
354 snew(w_rms, top.atoms.nr);
359 "WARNING: Need a run input file for bond angle matrix,\n"
360 " will not calculate bond angle matrix.\n");
366 fprintf(stderr, "Select group for %s fit\n", bFit ? "least squares"
368 get_index(&(top.atoms), ftp2fn_null(efNDX, NFILE, fnm), 1, &ifit,
378 if (bFit && ifit < 3)
380 gmx_fatal(FARGS, "Need >= 3 points to fit!\n" );
384 for (i = 0; i < ifit; i++)
388 w_rls[ind_fit[i]] = top.atoms.atom[ind_fit[i]].m;
392 w_rls[ind_fit[i]] = 1;
394 bMass = bMass || (top.atoms.atom[ind_fit[i]].m != 0);
398 fprintf(stderr, "All masses in the fit group are 0, using masses of 1\n");
399 for (i = 0; i < ifit; i++)
401 w_rls[ind_fit[i]] = 1;
415 fprintf(stderr, "Select group%s for %s calculation\n",
416 (nrms > 1) ? "s" : "", whatname[ewhat]);
417 get_index(&(top.atoms), ftp2fn_null(efNDX, NFILE, fnm),
418 nrms, irms, ind_rms, gn_rms);
422 snew(rlsnorm, irms[0]);
425 for (j = 0; j < nrms; j++)
427 snew(rls[j], maxframe);
432 for (j = 0; j < nrms; j++)
434 snew(rlsm[j], maxframe);
437 snew(time, maxframe);
438 for (j = 0; j < nrms; j++)
441 for (i = 0; i < irms[j]; i++)
445 w_rms[ind_rms[j][i]] = top.atoms.atom[ind_rms[j][i]].m;
449 w_rms[ind_rms[j][i]] = 1.0;
451 bMass = bMass || (top.atoms.atom[ind_rms[j][i]].m != 0);
455 fprintf(stderr, "All masses in group %d are 0, using masses of 1\n", j);
456 for (i = 0; i < irms[j]; i++)
458 w_rms[ind_rms[j][i]] = 1;
462 /* Prepare reference frame */
465 gpbc = gmx_rmpbc_init(&top.idef, ePBC, top.atoms.nr);
466 gmx_rmpbc(gpbc, top.atoms.nr, box, xp);
470 reset_x(ifit, ind_fit, top.atoms.nr, NULL, xp, w_rls);
474 /* generate reference structure mirror image: */
475 snew(xm, top.atoms.nr);
476 for (i = 0; i < top.atoms.nr; i++)
478 copy_rvec(xp[i], xm[i]);
479 xm[i][XX] = -xm[i][XX];
482 if (ewhat == ewRhoSc)
484 norm_princ(&top.atoms, ifit, ind_fit, top.atoms.nr, xp);
487 /* read first frame */
488 natoms_trx = read_first_x(oenv, &status, opt2fn("-f", NFILE, fnm), &t, &x, box);
489 if (natoms_trx != top.atoms.nr)
492 "\nWARNING: topology has %d atoms, whereas trajectory has %d\n",
493 top.atoms.nr, natoms_trx);
495 natoms = min(top.atoms.nr, natoms_trx);
496 if (bMat || bBond || bPrev)
502 /* With -prev we use all atoms for simplicity */
507 /* Check which atoms we need (fit/rms) */
508 snew(bInMat, natoms);
509 for (i = 0; i < ifit; i++)
511 bInMat[ind_fit[i]] = TRUE;
514 for (i = 0; i < irms[0]; i++)
516 if (!bInMat[ind_rms[0][i]])
518 bInMat[ind_rms[0][i]] = TRUE;
523 /* Make an index of needed atoms */
524 snew(ind_m, n_ind_m);
525 snew(rev_ind_m, natoms);
527 for (i = 0; i < natoms; i++)
529 if (bPrev || bInMat[i])
536 snew(w_rls_m, n_ind_m);
537 snew(ind_rms_m, irms[0]);
538 snew(w_rms_m, n_ind_m);
539 for (i = 0; i < ifit; i++)
541 w_rls_m[rev_ind_m[ind_fit[i]]] = w_rls[ind_fit[i]];
543 for (i = 0; i < irms[0]; i++)
545 ind_rms_m[i] = rev_ind_m[ind_rms[0][i]];
546 w_rms_m[ind_rms_m[i]] = w_rms[ind_rms[0][i]];
554 for (k = 0; k < F_NRE; k++)
558 iatom = top.idef.il[k].iatoms;
559 ncons += top.idef.il[k].nr/3;
562 fprintf(stderr, "Found %d bonds in topology\n", ncons);
563 snew(ind_bond1, ncons);
564 snew(ind_bond2, ncons);
566 for (k = 0; k < F_NRE; k++)
570 iatom = top.idef.il[k].iatoms;
571 ncons = top.idef.il[k].nr/3;
572 for (i = 0; i < ncons; i++)
576 for (j = 0; j < irms[0]; j++)
578 if (iatom[3*i+1] == ind_rms[0][j])
582 if (iatom[3*i+2] == ind_rms[0][j])
589 ind_bond1[ibond] = rev_ind_m[iatom[3*i+1]];
590 ind_bond2[ibond] = rev_ind_m[iatom[3*i+2]];
596 fprintf(stderr, "Using %d bonds for bond angle matrix\n", ibond);
599 gmx_fatal(FARGS, "0 bonds found");
603 /* start looping over frames: */
610 gmx_rmpbc(gpbc, natoms, box, x);
615 reset_x(ifit, ind_fit, natoms, NULL, x, w_rls);
617 if (ewhat == ewRhoSc)
619 norm_princ(&top.atoms, ifit, ind_fit, natoms, x);
624 /*do the least squares fit to original structure*/
625 do_fit(natoms, w_rls, xp, x);
628 if (teller % freq == 0)
630 /* keep frame for matrix calculation */
631 if (bMat || bBond || bPrev)
633 if (tel_mat >= NFRAME)
635 srenew(mat_x, tel_mat+1);
637 snew(mat_x[tel_mat], n_ind_m);
638 for (i = 0; i < n_ind_m; i++)
640 copy_rvec(x[ind_m[i]], mat_x[tel_mat][i]);
646 /*calculate energy of root_least_squares*/
654 for (i = 0; i < n_ind_m; i++)
656 copy_rvec(mat_x[j][i], xp[ind_m[i]]);
660 reset_x(ifit, ind_fit, natoms, NULL, xp, w_rls);
664 do_fit(natoms, w_rls, x, xp);
667 for (j = 0; (j < nrms); j++)
670 calc_similar_ind(ewhat != ewRMSD, irms[j], ind_rms[j], w_rms, x, xp);
674 for (j = 0; (j < irms[0]); j++)
677 calc_similar_ind(ewhat != ewRMSD, 1, &(ind_rms[0][j]), w_rms, x, xp);
685 /*do the least squares fit to mirror of original structure*/
686 do_fit(natoms, w_rls, xm, x);
689 for (j = 0; j < nrms; j++)
692 calc_similar_ind(ewhat != ewRMSD, irms[j], ind_rms[j], w_rms, x, xm);
695 time[teller] = output_env_conv_time(oenv, t);
698 if (teller >= maxframe)
701 srenew(time, maxframe);
702 for (j = 0; (j < nrms); j++)
704 srenew(rls[j], maxframe);
708 for (j = 0; (j < nrms); j++)
710 srenew(rlsm[j], maxframe);
715 while (read_next_x(oenv, status, &t, x, box));
720 snew(time2, maxframe2);
722 fprintf(stderr, "\nWill read second trajectory file\n");
723 snew(mat_x2, NFRAME);
725 read_first_x(oenv, &status, opt2fn("-f2", NFILE, fnm), &t, &x, box);
726 if (natoms_trx2 != natoms_trx)
729 "Second trajectory (%d atoms) does not match the first one"
730 " (%d atoms)", natoms_trx2, natoms_trx);
738 gmx_rmpbc(gpbc, natoms, box, x);
743 reset_x(ifit, ind_fit, natoms, NULL, x, w_rls);
745 if (ewhat == ewRhoSc)
747 norm_princ(&top.atoms, ifit, ind_fit, natoms, x);
752 /*do the least squares fit to original structure*/
753 do_fit(natoms, w_rls, xp, x);
756 if (teller2 % freq2 == 0)
758 /* keep frame for matrix calculation */
761 if (tel_mat2 >= NFRAME)
763 srenew(mat_x2, tel_mat2+1);
765 snew(mat_x2[tel_mat2], n_ind_m);
766 for (i = 0; i < n_ind_m; i++)
768 copy_rvec(x[ind_m[i]], mat_x2[tel_mat2][i]);
774 time2[teller2] = output_env_conv_time(oenv, t);
777 if (teller2 >= maxframe2)
780 srenew(time2, maxframe2);
783 while (read_next_x(oenv, status, &t, x, box));
793 gmx_rmpbc_done(gpbc);
797 /* calculate RMS matrix */
798 fprintf(stderr, "\n");
801 fprintf(stderr, "Building %s matrix, %dx%d elements\n",
802 whatname[ewhat], tel_mat, tel_mat2);
803 snew(rmsd_mat, tel_mat);
807 fprintf(stderr, "Building bond angle matrix, %dx%d elements\n",
809 snew(bond_mat, tel_mat);
812 snew(axis2, tel_mat2);
825 for (j = 0; j < tel_mat2; j++)
827 axis2[j] = time2[freq2*j];
833 delta_scalex = 8.0/log(2.0);
834 delta_xsize = (int)(log(tel_mat/2)*delta_scalex+0.5)+1;
838 delta_xsize = tel_mat/2;
840 delta_scaley = 1.0/delta_maxy;
841 snew(delta, delta_xsize);
842 for (j = 0; j < delta_xsize; j++)
844 snew(delta[j], del_lev+1);
848 snew(rmsdav_mat, tel_mat);
849 for (j = 0; j < tel_mat; j++)
851 snew(rmsdav_mat[j], tel_mat);
858 snew(mat_x2_j, natoms);
860 for (i = 0; i < tel_mat; i++)
862 axis[i] = time[freq*i];
863 fprintf(stderr, "\r element %5d; time %5.2f ", i, axis[i]);
866 snew(rmsd_mat[i], tel_mat2);
870 snew(bond_mat[i], tel_mat2);
872 for (j = 0; j < tel_mat2; j++)
876 for (k = 0; k < n_ind_m; k++)
878 copy_rvec(mat_x2[j][k], mat_x2_j[k]);
880 do_fit(n_ind_m, w_rls_m, mat_x[i], mat_x2_j);
884 mat_x2_j = mat_x2[j];
888 if (bFile2 || (i < j))
891 calc_similar_ind(ewhat != ewRMSD, irms[0], ind_rms_m,
892 w_rms_m, mat_x[i], mat_x2_j);
893 if (rmsd_mat[i][j] > rmsd_max)
895 rmsd_max = rmsd_mat[i][j];
897 if (rmsd_mat[i][j] < rmsd_min)
899 rmsd_min = rmsd_mat[i][j];
901 rmsd_avg += rmsd_mat[i][j];
905 rmsd_mat[i][j] = rmsd_mat[j][i];
910 if (bFile2 || (i <= j))
913 for (m = 0; m < ibond; m++)
915 rvec_sub(mat_x[i][ind_bond1[m]], mat_x[i][ind_bond2[m]], vec1);
916 rvec_sub(mat_x2_j[ind_bond1[m]], mat_x2_j[ind_bond2[m]], vec2);
917 ang += acos(cos_angle(vec1, vec2));
919 bond_mat[i][j] = ang*180.0/(M_PI*ibond);
920 if (bond_mat[i][j] > bond_max)
922 bond_max = bond_mat[i][j];
924 if (bond_mat[i][j] < bond_min)
926 bond_min = bond_mat[i][j];
931 bond_mat[i][j] = bond_mat[j][i];
938 rmsd_avg /= tel_mat*tel_mat2;
942 rmsd_avg /= tel_mat*(tel_mat - 1)/2;
944 if (bMat && (avl > 0))
949 for (j = 0; j < tel_mat-1; j++)
951 for (i = j+1; i < tel_mat; i++)
955 for (my = -avl; my <= avl; my++)
957 if ((j+my >= 0) && (j+my < tel_mat))
960 for (mx = -avl; mx <= avl; mx++)
962 if ((i+mx >= 0) && (i+mx < tel_mat))
964 weight = (real)(avl+1-max(abs(mx), abs_my));
965 av_tot += weight*rmsd_mat[i+mx][j+my];
966 weight_tot += weight;
971 rmsdav_mat[i][j] = av_tot/weight_tot;
972 rmsdav_mat[j][i] = rmsdav_mat[i][j];
973 if (rmsdav_mat[i][j] > rmsd_max)
975 rmsd_max = rmsdav_mat[i][j];
979 rmsd_mat = rmsdav_mat;
984 fprintf(stderr, "\n%s: Min %f, Max %f, Avg %f\n",
985 whatname[ewhat], rmsd_min, rmsd_max, rmsd_avg);
986 rlo.r = 1; rlo.g = 1; rlo.b = 1;
987 rhi.r = 0; rhi.g = 0; rhi.b = 0;
988 if (rmsd_user_max != -1)
990 rmsd_max = rmsd_user_max;
992 if (rmsd_user_min != -1)
994 rmsd_min = rmsd_user_min;
996 if ((rmsd_user_max != -1) || (rmsd_user_min != -1))
998 fprintf(stderr, "Min and Max value set to resp. %f and %f\n",
1001 sprintf(buf, "%s %s matrix", gn_rms[0], whatname[ewhat]);
1002 write_xpm(opt2FILE("-m", NFILE, fnm, "w"), 0, buf, whatlabel[ewhat],
1003 output_env_get_time_label(oenv), output_env_get_time_label(oenv), tel_mat, tel_mat2,
1004 axis, axis2, rmsd_mat, rmsd_min, rmsd_max, rlo, rhi, &nlevels);
1005 /* Print the distribution of RMSD values */
1006 if (opt2bSet("-dist", NFILE, fnm))
1008 low_rmsd_dist(opt2fn("-dist", NFILE, fnm), rmsd_max, tel_mat, rmsd_mat, oenv);
1013 snew(delta_tot, delta_xsize);
1014 for (j = 0; j < tel_mat-1; j++)
1016 for (i = j+1; i < tel_mat; i++)
1023 mx = (int)(log(mx)*delta_scalex+0.5);
1025 my = (int)(rmsd_mat[i][j]*delta_scaley*del_lev+0.5);
1026 delta_tot[mx] += 1.0;
1027 if ((rmsd_mat[i][j] >= 0) && (rmsd_mat[i][j] <= delta_maxy))
1029 delta[mx][my] += 1.0;
1035 for (i = 0; i < delta_xsize; i++)
1037 if (delta_tot[i] > 0.0)
1039 delta_tot[i] = 1.0/delta_tot[i];
1040 for (j = 0; j <= del_lev; j++)
1042 delta[i][j] *= delta_tot[i];
1043 if (delta[i][j] > delta_max)
1045 delta_max = delta[i][j];
1050 fprintf(stderr, "Maximum in delta matrix: %f\n", delta_max);
1051 snew(del_xaxis, delta_xsize);
1052 snew(del_yaxis, del_lev+1);
1053 for (i = 0; i < delta_xsize; i++)
1055 del_xaxis[i] = axis[i]-axis[0];
1057 for (i = 0; i < del_lev+1; i++)
1059 del_yaxis[i] = delta_maxy*i/del_lev;
1061 sprintf(buf, "%s %s vs. delta t", gn_rms[0], whatname[ewhat]);
1062 fp = ffopen("delta.xpm", "w");
1063 write_xpm(fp, 0, buf, "density", output_env_get_time_label(oenv), whatlabel[ewhat],
1064 delta_xsize, del_lev+1, del_xaxis, del_yaxis,
1065 delta, 0.0, delta_max, rlo, rhi, &nlevels);
1068 if (opt2bSet("-bin", NFILE, fnm))
1070 /* NB: File must be binary if we use fwrite */
1071 fp = ftp2FILE(efDAT, NFILE, fnm, "wb");
1072 for (i = 0; i < tel_mat; i++)
1074 if (fwrite(rmsd_mat[i], sizeof(**rmsd_mat), tel_mat2, fp) != tel_mat2)
1076 gmx_fatal(FARGS, "Error writing to output file");
1084 fprintf(stderr, "\nMin. angle: %f, Max. angle: %f\n", bond_min, bond_max);
1085 if (bond_user_max != -1)
1087 bond_max = bond_user_max;
1089 if (bond_user_min != -1)
1091 bond_min = bond_user_min;
1093 if ((bond_user_max != -1) || (bond_user_min != -1))
1095 fprintf(stderr, "Bond angle Min and Max set to:\n"
1096 "Min. angle: %f, Max. angle: %f\n", bond_min, bond_max);
1098 rlo.r = 1; rlo.g = 1; rlo.b = 1;
1099 rhi.r = 0; rhi.g = 0; rhi.b = 0;
1100 sprintf(buf, "%s av. bond angle deviation", gn_rms[0]);
1101 write_xpm(opt2FILE("-bm", NFILE, fnm, "w"), 0, buf, "degrees",
1102 output_env_get_time_label(oenv), output_env_get_time_label(oenv), tel_mat, tel_mat2,
1103 axis, axis2, bond_mat, bond_min, bond_max, rlo, rhi, &nlevels);
1107 bAv = opt2bSet("-a", NFILE, fnm);
1109 /* Write the RMSD's to file */
1112 sprintf(buf, "%s", whatxvgname[ewhat]);
1116 sprintf(buf, "%s with frame %g %s ago", whatxvgname[ewhat],
1117 time[prev*freq]-time[0], output_env_get_time_label(oenv));
1119 fp = xvgropen(opt2fn("-o", NFILE, fnm), buf, output_env_get_xvgr_tlabel(oenv),
1120 whatxvglabel[ewhat], oenv);
1121 if (output_env_get_print_xvgr_codes(oenv))
1123 fprintf(fp, "@ subtitle \"%s%s after %s%s%s\"\n",
1124 (nrms == 1) ? "" : "of ", gn_rms[0], fitgraphlabel[efit],
1125 bFit ? " to " : "", bFit ? gn_fit : "");
1129 xvgr_legend(fp, nrms, (const char**)gn_rms, oenv);
1131 for (i = 0; (i < teller); i++)
1133 if (bSplit && i > 0 &&
1134 abs(time[bPrev ? freq*i : i]/output_env_get_time_factor(oenv)) < 1e-5)
1138 fprintf(fp, "%12.7f", time[bPrev ? freq*i : i]);
1139 for (j = 0; (j < nrms); j++)
1141 fprintf(fp, " %12.7f", rls[j][i]);
1144 rlstot += rls[j][i];
1153 /* Write the mirror RMSD's to file */
1154 sprintf(buf, "%s with Mirror", whatxvgname[ewhat]);
1155 sprintf(buf2, "Mirror %s", whatxvglabel[ewhat]);
1156 fp = xvgropen(opt2fn("-mir", NFILE, fnm), buf, output_env_get_xvgr_tlabel(oenv),
1160 if (output_env_get_print_xvgr_codes(oenv))
1162 fprintf(fp, "@ subtitle \"of %s after lsq fit to mirror of %s\"\n",
1168 if (output_env_get_print_xvgr_codes(oenv))
1170 fprintf(fp, "@ subtitle \"after lsq fit to mirror %s\"\n", gn_fit);
1172 xvgr_legend(fp, nrms, (const char**)gn_rms, oenv);
1174 for (i = 0; (i < teller); i++)
1176 if (bSplit && i > 0 && abs(time[i]) < 1e-5)
1180 fprintf(fp, "%12.7f", time[i]);
1181 for (j = 0; (j < nrms); j++)
1183 fprintf(fp, " %12.7f", rlsm[j][i]);
1192 sprintf(buf, "Average %s", whatxvgname[ewhat]);
1193 sprintf(buf2, "Average %s", whatxvglabel[ewhat]);
1194 fp = xvgropen(opt2fn("-a", NFILE, fnm), buf, "Residue", buf2, oenv);
1195 for (j = 0; (j < nrms); j++)
1197 fprintf(fp, "%10d %10g\n", j, rlstot/teller);
1204 fp = xvgropen("aver.xvg", gn_rms[0], "Residue", whatxvglabel[ewhat], oenv);
1205 for (j = 0; (j < irms[0]); j++)
1207 fprintf(fp, "%10d %10g\n", j, rlsnorm[j]/teller);
1211 do_view(oenv, opt2fn_null("-a", NFILE, fnm), "-graphtype bar");
1212 do_view(oenv, opt2fn("-o", NFILE, fnm), NULL);
1213 do_view(oenv, opt2fn_null("-mir", NFILE, fnm), NULL);
1214 do_view(oenv, opt2fn_null("-m", NFILE, fnm), NULL);
1215 do_view(oenv, opt2fn_null("-bm", NFILE, fnm), NULL);
1216 do_view(oenv, opt2fn_null("-dist", NFILE, fnm), NULL);