Removed legacyheaders/typedefs.h
[gromacs.git] / src / gromacs / gmxana / gmx_rms.cpp
blob603ccae674ffb25d030aa32799b177b9cc371a58
1 /*
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,2015, 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.
37 #include "gmxpre.h"
39 #include <cmath>
40 #include <cstdlib>
42 #include <algorithm>
44 #include "gromacs/commandline/pargs.h"
45 #include "gromacs/commandline/viewit.h"
46 #include "gromacs/fileio/confio.h"
47 #include "gromacs/fileio/matio.h"
48 #include "gromacs/fileio/trxio.h"
49 #include "gromacs/fileio/xvgr.h"
50 #include "gromacs/gmxana/cmat.h"
51 #include "gromacs/gmxana/gmx_ana.h"
52 #include "gromacs/gmxana/princ.h"
53 #include "gromacs/legacyheaders/copyrite.h"
54 #include "gromacs/legacyheaders/types/ifunc.h"
55 #include "gromacs/math/do_fit.h"
56 #include "gromacs/math/vec.h"
57 #include "gromacs/pbcutil/rmpbc.h"
58 #include "gromacs/topology/index.h"
59 #include "gromacs/topology/topology.h"
60 #include "gromacs/utility/arraysize.h"
61 #include "gromacs/utility/fatalerror.h"
62 #include "gromacs/utility/futil.h"
63 #include "gromacs/utility/smalloc.h"
65 static void norm_princ(t_atoms *atoms, int isize, atom_id *index, int natoms,
66 rvec *x)
68 int i, m;
69 rvec princ, vec;
71 /* equalize principal components: */
72 /* orient principal axes, get principal components */
73 orient_princ(atoms, isize, index, natoms, x, NULL, princ);
75 /* calc our own principal components */
76 clear_rvec(vec);
77 for (m = 0; m < DIM; m++)
79 for (i = 0; i < isize; i++)
81 vec[m] += sqr(x[index[i]][m]);
83 vec[m] = std::sqrt(vec[m] / isize);
84 /* calculate scaling constants */
85 vec[m] = 1.0 / (std::sqrt(3.0) * vec[m]);
88 /* scale coordinates */
89 for (i = 0; i < natoms; i++)
91 for (m = 0; m < DIM; m++)
93 x[i][m] *= vec[m];
98 int gmx_rms(int argc, char *argv[])
100 const char *desc[] =
102 "[THISMODULE] compares two structures by computing the root mean square",
103 "deviation (RMSD), the size-independent [GRK]rho[grk] similarity parameter",
104 "([TT]rho[tt]) or the scaled [GRK]rho[grk] ([TT]rhosc[tt]), ",
105 "see Maiorov & Crippen, Proteins [BB]22[bb], 273 (1995).",
106 "This is selected by [TT]-what[tt].[PAR]"
108 "Each structure from a trajectory ([TT]-f[tt]) is compared to a",
109 "reference structure. The reference structure",
110 "is taken from the structure file ([TT]-s[tt]).[PAR]",
112 "With option [TT]-mir[tt] also a comparison with the mirror image of",
113 "the reference structure is calculated.",
114 "This is useful as a reference for 'significant' values, see",
115 "Maiorov & Crippen, Proteins [BB]22[bb], 273 (1995).[PAR]",
117 "Option [TT]-prev[tt] produces the comparison with a previous frame",
118 "the specified number of frames ago.[PAR]",
120 "Option [TT]-m[tt] produces a matrix in [REF].xpm[ref] format of",
121 "comparison values of each structure in the trajectory with respect to",
122 "each other structure. This file can be visualized with for instance",
123 "[TT]xv[tt] and can be converted to postscript with [gmx-xpm2ps].[PAR]",
125 "Option [TT]-fit[tt] controls the least-squares fitting of",
126 "the structures on top of each other: complete fit (rotation and",
127 "translation), translation only, or no fitting at all.[PAR]",
129 "Option [TT]-mw[tt] controls whether mass weighting is done or not.",
130 "If you select the option (default) and ",
131 "supply a valid [REF].tpr[ref] file masses will be taken from there, ",
132 "otherwise the masses will be deduced from the [TT]atommass.dat[tt] file in",
133 "[TT]GMXLIB[tt]. This is fine for proteins, but not",
134 "necessarily for other molecules. A default mass of 12.011 amu (carbon)",
135 "is assigned to unknown atoms. You can check whether this happend by",
136 "turning on the [TT]-debug[tt] flag and inspecting the log file.[PAR]",
138 "With [TT]-f2[tt], the 'other structures' are taken from a second",
139 "trajectory, this generates a comparison matrix of one trajectory",
140 "versus the other.[PAR]",
142 "Option [TT]-bin[tt] does a binary dump of the comparison matrix.[PAR]",
144 "Option [TT]-bm[tt] produces a matrix of average bond angle deviations",
145 "analogously to the [TT]-m[tt] option. Only bonds between atoms in the",
146 "comparison group are considered."
148 static gmx_bool bPBC = TRUE, bFitAll = TRUE, bSplit = FALSE;
149 static gmx_bool bDeltaLog = FALSE;
150 static int prev = 0, freq = 1, freq2 = 1, nlevels = 80, avl = 0;
151 static real rmsd_user_max = -1, rmsd_user_min = -1, bond_user_max = -1,
152 bond_user_min = -1, delta_maxy = 0.0;
153 /* strings and things for selecting difference method */
154 enum
156 ewSel, ewRMSD, ewRho, ewRhoSc, ewNR
158 int ewhat;
159 const char *what[ewNR + 1] =
160 { NULL, "rmsd", "rho", "rhosc", NULL };
161 const char *whatname[ewNR] =
162 { NULL, "RMSD", "Rho", "Rho sc" };
163 const char *whatlabel[ewNR] =
164 { NULL, "RMSD (nm)", "Rho", "Rho sc" };
165 const char *whatxvgname[ewNR] =
166 { NULL, "RMSD", "\\8r\\4", "\\8r\\4\\ssc\\N" };
167 const char *whatxvglabel[ewNR] =
168 { NULL, "RMSD (nm)", "\\8r\\4", "\\8r\\4\\ssc\\N" };
169 /* strings and things for fitting methods */
170 enum
172 efSel, efFit, efReset, efNone, efNR
174 int efit;
175 const char *fit[efNR + 1] =
176 { NULL, "rot+trans", "translation", "none", NULL };
177 const char *fitgraphlabel[efNR + 1] =
178 { NULL, "lsq fit", "translational fit", "no fit" };
179 static int nrms = 1;
180 static gmx_bool bMassWeighted = TRUE;
181 t_pargs pa[] =
183 { "-what", FALSE, etENUM,
184 { what }, "Structural difference measure" },
185 { "-pbc", FALSE, etBOOL,
186 { &bPBC }, "PBC check" },
187 { "-fit", FALSE, etENUM,
188 { fit }, "Fit to reference structure" },
189 { "-prev", FALSE, etINT,
190 { &prev }, "Compare with previous frame" },
191 { "-split", FALSE, etBOOL,
192 { &bSplit }, "Split graph where time is zero" },
193 { "-fitall", FALSE, etBOOL,
194 { &bFitAll }, "HIDDENFit all pairs of structures in matrix" },
195 { "-skip", FALSE, etINT,
196 { &freq }, "Only write every nr-th frame to matrix" },
197 { "-skip2", FALSE, etINT,
198 { &freq2 }, "Only write every nr-th frame to matrix" },
199 { "-max", FALSE, etREAL,
200 { &rmsd_user_max }, "Maximum level in comparison matrix" },
201 { "-min", FALSE, etREAL,
202 { &rmsd_user_min }, "Minimum level in comparison matrix" },
203 { "-bmax", FALSE, etREAL,
204 { &bond_user_max }, "Maximum level in bond angle matrix" },
205 { "-bmin", FALSE, etREAL,
206 { &bond_user_min }, "Minimum level in bond angle matrix" },
207 { "-mw", FALSE, etBOOL,
208 { &bMassWeighted }, "Use mass weighting for superposition" },
209 { "-nlevels", FALSE, etINT,
210 { &nlevels }, "Number of levels in the matrices" },
211 { "-ng", FALSE, etINT,
212 { &nrms }, "Number of groups to compute RMS between" },
213 { "-dlog", FALSE, etBOOL,
214 { &bDeltaLog },
215 "HIDDENUse a log x-axis in the delta t matrix" },
216 { "-dmax", FALSE, etREAL,
217 { &delta_maxy }, "HIDDENMaximum level in delta matrix" },
218 { "-aver", FALSE, etINT,
219 { &avl },
220 "HIDDENAverage over this distance in the RMSD matrix" }
222 int natoms_trx, natoms_trx2, natoms;
223 int i, j, k, m, teller, teller2, tel_mat, tel_mat2;
224 #define NFRAME 5000
225 int maxframe = NFRAME, maxframe2 = NFRAME;
226 real t, *w_rls, *w_rms, *w_rls_m = NULL, *w_rms_m = NULL;
227 gmx_bool bNorm, bAv, bFreq2, bFile2, bMat, bBond, bDelta, bMirror, bMass;
228 gmx_bool bFit, bReset;
229 t_topology top;
230 int ePBC;
231 t_iatom *iatom = NULL;
233 matrix box = {{0}};
234 rvec *x, *xp, *xm = NULL, **mat_x = NULL, **mat_x2, *mat_x2_j = NULL, vec1,
235 vec2;
236 t_trxstatus *status;
237 char buf[256], buf2[256];
238 int ncons = 0;
239 FILE *fp;
240 real rlstot = 0, **rls, **rlsm = NULL, *time, *time2, *rlsnorm = NULL,
241 **rmsd_mat = NULL, **bond_mat = NULL, *axis, *axis2, *del_xaxis,
242 *del_yaxis, rmsd_max, rmsd_min, rmsd_avg, bond_max, bond_min, ang;
243 real **rmsdav_mat = NULL, av_tot, weight, weight_tot;
244 real **delta = NULL, delta_max, delta_scalex = 0, delta_scaley = 0,
245 *delta_tot;
246 int delta_xsize = 0, del_lev = 100, mx, my, abs_my;
247 gmx_bool bA1, bA2, bPrev, bTop, *bInMat = NULL;
248 int ifit, *irms, ibond = 0, *ind_bond1 = NULL, *ind_bond2 = NULL, n_ind_m =
250 atom_id *ind_fit, **ind_rms, *ind_m = NULL, *rev_ind_m = NULL, *ind_rms_m =
251 NULL;
252 char *gn_fit, **gn_rms;
253 t_rgb rlo, rhi;
254 gmx_output_env_t *oenv;
255 gmx_rmpbc_t gpbc = NULL;
257 t_filenm fnm[] =
259 { efTPS, NULL, NULL, ffREAD },
260 { efTRX, "-f", NULL, ffREAD },
261 { efTRX, "-f2", NULL, ffOPTRD },
262 { efNDX, NULL, NULL, ffOPTRD },
263 { efXVG, NULL, "rmsd", ffWRITE },
264 { efXVG, "-mir", "rmsdmir", ffOPTWR },
265 { efXVG, "-a", "avgrp", ffOPTWR },
266 { efXVG, "-dist", "rmsd-dist", ffOPTWR },
267 { efXPM, "-m", "rmsd", ffOPTWR },
268 { efDAT, "-bin", "rmsd", ffOPTWR },
269 { efXPM, "-bm", "bond", ffOPTWR }
271 #define NFILE asize(fnm)
273 if (!parse_common_args(&argc, argv, PCA_CAN_TIME | PCA_TIME_UNIT | PCA_CAN_VIEW,
274 NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL,
275 &oenv))
277 return 0;
279 /* parse enumerated options: */
280 ewhat = nenum(what);
281 if (ewhat == ewRho || ewhat == ewRhoSc)
283 please_cite(stdout, "Maiorov95");
285 efit = nenum(fit);
286 bFit = efit == efFit;
287 bReset = efit == efReset;
288 if (bFit)
290 bReset = TRUE; /* for fit, reset *must* be set */
292 else
294 bFitAll = FALSE;
297 /* mark active cmdline options */
298 bMirror = opt2bSet("-mir", NFILE, fnm); /* calc RMSD vs mirror of ref. */
299 bFile2 = opt2bSet("-f2", NFILE, fnm);
300 bMat = opt2bSet("-m", NFILE, fnm);
301 bBond = opt2bSet("-bm", NFILE, fnm);
302 bDelta = (delta_maxy > 0); /* calculate rmsd vs delta t matrix from *
303 * your RMSD matrix (hidden option */
304 bNorm = opt2bSet("-a", NFILE, fnm);
305 bFreq2 = opt2parg_bSet("-skip2", asize(pa), pa);
306 if (freq <= 0)
308 fprintf(stderr, "The number of frames to skip is <= 0. "
309 "Writing out all frames.\n\n");
310 freq = 1;
312 if (!bFreq2)
314 freq2 = freq;
316 else if (bFile2 && freq2 <= 0)
318 fprintf(stderr,
319 "The number of frames to skip in second trajectory is <= 0.\n"
320 " Writing out all frames.\n\n");
321 freq2 = 1;
324 bPrev = (prev > 0);
325 if (bPrev)
327 fprintf(stderr, "WARNING: using option -prev with large trajectories will\n"
328 " require a lot of memory and could lead to crashes\n");
329 prev = abs(prev);
330 if (freq != 1)
332 fprintf(stderr, "WARNING: option -skip also applies to -prev\n");
336 if (bFile2 && !bMat && !bBond)
338 fprintf(
339 stderr,
340 "WARNING: second trajectory (-f2) useless when not calculating matrix (-m/-bm),\n"
341 " will not read from %s\n", opt2fn("-f2", NFILE,
342 fnm));
343 bFile2 = FALSE;
346 if (bDelta)
348 bMat = TRUE;
349 if (bFile2)
351 fprintf(stderr,
352 "WARNING: second trajectory (-f2) useless when making delta matrix,\n"
353 " will not read from %s\n", opt2fn("-f2",
354 NFILE, fnm));
355 bFile2 = FALSE;
359 bTop = read_tps_conf(ftp2fn(efTPS, NFILE, fnm), &top, &ePBC, &xp,
360 NULL, box, TRUE);
361 snew(w_rls, top.atoms.nr);
362 snew(w_rms, top.atoms.nr);
364 if (!bTop && bBond)
366 fprintf(stderr,
367 "WARNING: Need a run input file for bond angle matrix,\n"
368 " will not calculate bond angle matrix.\n");
369 bBond = FALSE;
372 if (bReset)
374 fprintf(stderr, "Select group for %s fit\n", bFit ? "least squares"
375 : "translational");
376 get_index(&(top.atoms), ftp2fn_null(efNDX, NFILE, fnm), 1, &ifit,
377 &ind_fit, &gn_fit);
379 else
381 ifit = 0;
384 if (bReset)
386 if (bFit && ifit < 3)
388 gmx_fatal(FARGS, "Need >= 3 points to fit!\n" );
391 bMass = FALSE;
392 for (i = 0; i < ifit; i++)
394 if (bMassWeighted)
396 w_rls[ind_fit[i]] = top.atoms.atom[ind_fit[i]].m;
398 else
400 w_rls[ind_fit[i]] = 1;
402 bMass = bMass || (top.atoms.atom[ind_fit[i]].m != 0);
404 if (!bMass)
406 fprintf(stderr, "All masses in the fit group are 0, using masses of 1\n");
407 for (i = 0; i < ifit; i++)
409 w_rls[ind_fit[i]] = 1;
414 if (bMat || bBond)
416 nrms = 1;
419 snew(gn_rms, nrms);
420 snew(ind_rms, nrms);
421 snew(irms, nrms);
423 fprintf(stderr, "Select group%s for %s calculation\n",
424 (nrms > 1) ? "s" : "", whatname[ewhat]);
425 get_index(&(top.atoms), ftp2fn_null(efNDX, NFILE, fnm),
426 nrms, irms, ind_rms, gn_rms);
428 if (bNorm)
430 snew(rlsnorm, irms[0]);
432 snew(rls, nrms);
433 for (j = 0; j < nrms; j++)
435 snew(rls[j], maxframe);
437 if (bMirror)
439 snew(rlsm, nrms);
440 for (j = 0; j < nrms; j++)
442 snew(rlsm[j], maxframe);
445 snew(time, maxframe);
446 for (j = 0; j < nrms; j++)
448 bMass = FALSE;
449 for (i = 0; i < irms[j]; i++)
451 if (bMassWeighted)
453 w_rms[ind_rms[j][i]] = top.atoms.atom[ind_rms[j][i]].m;
455 else
457 w_rms[ind_rms[j][i]] = 1.0;
459 bMass = bMass || (top.atoms.atom[ind_rms[j][i]].m != 0);
461 if (!bMass)
463 fprintf(stderr, "All masses in group %d are 0, using masses of 1\n", j);
464 for (i = 0; i < irms[j]; i++)
466 w_rms[ind_rms[j][i]] = 1;
470 /* Prepare reference frame */
471 if (bPBC)
473 gpbc = gmx_rmpbc_init(&top.idef, ePBC, top.atoms.nr);
474 gmx_rmpbc(gpbc, top.atoms.nr, box, xp);
476 if (bReset)
478 reset_x(ifit, ind_fit, top.atoms.nr, NULL, xp, w_rls);
480 if (bMirror)
482 /* generate reference structure mirror image: */
483 snew(xm, top.atoms.nr);
484 for (i = 0; i < top.atoms.nr; i++)
486 copy_rvec(xp[i], xm[i]);
487 xm[i][XX] = -xm[i][XX];
490 if (ewhat == ewRhoSc)
492 norm_princ(&top.atoms, ifit, ind_fit, top.atoms.nr, xp);
495 /* read first frame */
496 natoms_trx = read_first_x(oenv, &status, opt2fn("-f", NFILE, fnm), &t, &x, box);
497 if (natoms_trx != top.atoms.nr)
499 fprintf(stderr,
500 "\nWARNING: topology has %d atoms, whereas trajectory has %d\n",
501 top.atoms.nr, natoms_trx);
503 natoms = std::min(top.atoms.nr, natoms_trx);
504 if (bMat || bBond || bPrev)
506 snew(mat_x, NFRAME);
508 if (bPrev)
510 /* With -prev we use all atoms for simplicity */
511 n_ind_m = natoms;
513 else
515 /* Check which atoms we need (fit/rms) */
516 snew(bInMat, natoms);
517 for (i = 0; i < ifit; i++)
519 bInMat[ind_fit[i]] = TRUE;
521 n_ind_m = ifit;
522 for (i = 0; i < irms[0]; i++)
524 if (!bInMat[ind_rms[0][i]])
526 bInMat[ind_rms[0][i]] = TRUE;
527 n_ind_m++;
531 /* Make an index of needed atoms */
532 snew(ind_m, n_ind_m);
533 snew(rev_ind_m, natoms);
534 j = 0;
535 for (i = 0; i < natoms; i++)
537 if (bPrev || bInMat[i])
539 ind_m[j] = i;
540 rev_ind_m[i] = j;
541 j++;
544 snew(w_rls_m, n_ind_m);
545 snew(ind_rms_m, irms[0]);
546 snew(w_rms_m, n_ind_m);
547 for (i = 0; i < ifit; i++)
549 w_rls_m[rev_ind_m[ind_fit[i]]] = w_rls[ind_fit[i]];
551 for (i = 0; i < irms[0]; i++)
553 ind_rms_m[i] = rev_ind_m[ind_rms[0][i]];
554 w_rms_m[ind_rms_m[i]] = w_rms[ind_rms[0][i]];
556 sfree(bInMat);
559 if (bBond)
561 ncons = 0;
562 for (k = 0; k < F_NRE; k++)
564 if (IS_CHEMBOND(k))
566 ncons += top.idef.il[k].nr/3;
569 fprintf(stderr, "Found %d bonds in topology\n", ncons);
570 snew(ind_bond1, ncons);
571 snew(ind_bond2, ncons);
572 ibond = 0;
573 for (k = 0; k < F_NRE; k++)
575 if (IS_CHEMBOND(k))
577 iatom = top.idef.il[k].iatoms;
578 ncons = top.idef.il[k].nr/3;
579 for (i = 0; i < ncons; i++)
581 bA1 = FALSE;
582 bA2 = FALSE;
583 for (j = 0; j < irms[0]; j++)
585 if (iatom[3*i+1] == ind_rms[0][j])
587 bA1 = TRUE;
589 if (iatom[3*i+2] == ind_rms[0][j])
591 bA2 = TRUE;
594 if (bA1 && bA2)
596 ind_bond1[ibond] = rev_ind_m[iatom[3*i+1]];
597 ind_bond2[ibond] = rev_ind_m[iatom[3*i+2]];
598 ibond++;
603 fprintf(stderr, "Using %d bonds for bond angle matrix\n", ibond);
604 if (ibond == 0)
606 gmx_fatal(FARGS, "0 bonds found");
610 /* start looping over frames: */
611 tel_mat = 0;
612 teller = 0;
615 if (bPBC)
617 gmx_rmpbc(gpbc, natoms, box, x);
620 if (bReset)
622 reset_x(ifit, ind_fit, natoms, NULL, x, w_rls);
624 if (ewhat == ewRhoSc)
626 norm_princ(&top.atoms, ifit, ind_fit, natoms, x);
629 if (bFit)
631 /*do the least squares fit to original structure*/
632 do_fit(natoms, w_rls, xp, x);
635 if (teller % freq == 0)
637 /* keep frame for matrix calculation */
638 if (bMat || bBond || bPrev)
640 if (tel_mat >= NFRAME)
642 srenew(mat_x, tel_mat+1);
644 snew(mat_x[tel_mat], n_ind_m);
645 for (i = 0; i < n_ind_m; i++)
647 copy_rvec(x[ind_m[i]], mat_x[tel_mat][i]);
650 tel_mat++;
653 /*calculate energy of root_least_squares*/
654 if (bPrev)
656 j = tel_mat-prev-1;
657 if (j < 0)
659 j = 0;
661 for (i = 0; i < n_ind_m; i++)
663 copy_rvec(mat_x[j][i], xp[ind_m[i]]);
665 if (bReset)
667 reset_x(ifit, ind_fit, natoms, NULL, xp, w_rls);
669 if (bFit)
671 do_fit(natoms, w_rls, x, xp);
674 for (j = 0; (j < nrms); j++)
676 rls[j][teller] =
677 calc_similar_ind(ewhat != ewRMSD, irms[j], ind_rms[j], w_rms, x, xp);
679 if (bNorm)
681 for (j = 0; (j < irms[0]); j++)
683 rlsnorm[j] +=
684 calc_similar_ind(ewhat != ewRMSD, 1, &(ind_rms[0][j]), w_rms, x, xp);
688 if (bMirror)
690 if (bFit)
692 /*do the least squares fit to mirror of original structure*/
693 do_fit(natoms, w_rls, xm, x);
696 for (j = 0; j < nrms; j++)
698 rlsm[j][teller] =
699 calc_similar_ind(ewhat != ewRMSD, irms[j], ind_rms[j], w_rms, x, xm);
702 time[teller] = output_env_conv_time(oenv, t);
704 teller++;
705 if (teller >= maxframe)
707 maxframe += NFRAME;
708 srenew(time, maxframe);
709 for (j = 0; (j < nrms); j++)
711 srenew(rls[j], maxframe);
713 if (bMirror)
715 for (j = 0; (j < nrms); j++)
717 srenew(rlsm[j], maxframe);
722 while (read_next_x(oenv, status, &t, x, box));
723 close_trj(status);
725 if (bFile2)
727 snew(time2, maxframe2);
729 fprintf(stderr, "\nWill read second trajectory file\n");
730 snew(mat_x2, NFRAME);
731 natoms_trx2 =
732 read_first_x(oenv, &status, opt2fn("-f2", NFILE, fnm), &t, &x, box);
733 if (natoms_trx2 != natoms_trx)
735 gmx_fatal(FARGS,
736 "Second trajectory (%d atoms) does not match the first one"
737 " (%d atoms)", natoms_trx2, natoms_trx);
739 tel_mat2 = 0;
740 teller2 = 0;
743 if (bPBC)
745 gmx_rmpbc(gpbc, natoms, box, x);
748 if (bReset)
750 reset_x(ifit, ind_fit, natoms, NULL, x, w_rls);
752 if (ewhat == ewRhoSc)
754 norm_princ(&top.atoms, ifit, ind_fit, natoms, x);
757 if (bFit)
759 /*do the least squares fit to original structure*/
760 do_fit(natoms, w_rls, xp, x);
763 if (teller2 % freq2 == 0)
765 /* keep frame for matrix calculation */
766 if (bMat)
768 if (tel_mat2 >= NFRAME)
770 srenew(mat_x2, tel_mat2+1);
772 snew(mat_x2[tel_mat2], n_ind_m);
773 for (i = 0; i < n_ind_m; i++)
775 copy_rvec(x[ind_m[i]], mat_x2[tel_mat2][i]);
778 tel_mat2++;
781 time2[teller2] = output_env_conv_time(oenv, t);
783 teller2++;
784 if (teller2 >= maxframe2)
786 maxframe2 += NFRAME;
787 srenew(time2, maxframe2);
790 while (read_next_x(oenv, status, &t, x, box));
791 close_trj(status);
793 else
795 mat_x2 = mat_x;
796 time2 = time;
797 tel_mat2 = tel_mat;
798 freq2 = freq;
800 gmx_rmpbc_done(gpbc);
802 if (bMat || bBond)
804 /* calculate RMS matrix */
805 fprintf(stderr, "\n");
806 if (bMat)
808 fprintf(stderr, "Building %s matrix, %dx%d elements\n",
809 whatname[ewhat], tel_mat, tel_mat2);
810 snew(rmsd_mat, tel_mat);
812 if (bBond)
814 fprintf(stderr, "Building bond angle matrix, %dx%d elements\n",
815 tel_mat, tel_mat2);
816 snew(bond_mat, tel_mat);
818 snew(axis, tel_mat);
819 snew(axis2, tel_mat2);
820 rmsd_max = 0;
821 if (bFile2)
823 rmsd_min = 1e10;
825 else
827 rmsd_min = 0;
829 rmsd_avg = 0;
830 bond_max = 0;
831 bond_min = 1e10;
832 for (j = 0; j < tel_mat2; j++)
834 axis2[j] = time2[freq2*j];
836 if (bDelta)
838 if (bDeltaLog)
840 delta_scalex = 8.0/std::log(2.0);
841 delta_xsize = static_cast<int>(std::log(static_cast<real>(tel_mat/2))*delta_scalex+0.5)+1;
843 else
845 delta_xsize = tel_mat/2;
847 delta_scaley = 1.0/delta_maxy;
848 snew(delta, delta_xsize);
849 for (j = 0; j < delta_xsize; j++)
851 snew(delta[j], del_lev+1);
853 if (avl > 0)
855 snew(rmsdav_mat, tel_mat);
856 for (j = 0; j < tel_mat; j++)
858 snew(rmsdav_mat[j], tel_mat);
863 if (bFitAll)
865 snew(mat_x2_j, natoms);
867 for (i = 0; i < tel_mat; i++)
869 axis[i] = time[freq*i];
870 fprintf(stderr, "\r element %5d; time %5.2f ", i, axis[i]);
871 if (bMat)
873 snew(rmsd_mat[i], tel_mat2);
875 if (bBond)
877 snew(bond_mat[i], tel_mat2);
879 for (j = 0; j < tel_mat2; j++)
881 if (bFitAll)
883 for (k = 0; k < n_ind_m; k++)
885 copy_rvec(mat_x2[j][k], mat_x2_j[k]);
887 do_fit(n_ind_m, w_rls_m, mat_x[i], mat_x2_j);
889 else
891 mat_x2_j = mat_x2[j];
893 if (bMat)
895 if (bFile2 || (i < j))
897 rmsd_mat[i][j] =
898 calc_similar_ind(ewhat != ewRMSD, irms[0], ind_rms_m,
899 w_rms_m, mat_x[i], mat_x2_j);
900 if (rmsd_mat[i][j] > rmsd_max)
902 rmsd_max = rmsd_mat[i][j];
904 if (rmsd_mat[i][j] < rmsd_min)
906 rmsd_min = rmsd_mat[i][j];
908 rmsd_avg += rmsd_mat[i][j];
910 else
912 rmsd_mat[i][j] = rmsd_mat[j][i];
915 if (bBond)
917 if (bFile2 || (i <= j))
919 ang = 0.0;
920 for (m = 0; m < ibond; m++)
922 rvec_sub(mat_x[i][ind_bond1[m]], mat_x[i][ind_bond2[m]], vec1);
923 rvec_sub(mat_x2_j[ind_bond1[m]], mat_x2_j[ind_bond2[m]], vec2);
924 ang += std::acos(cos_angle(vec1, vec2));
926 bond_mat[i][j] = ang*180.0/(M_PI*ibond);
927 if (bond_mat[i][j] > bond_max)
929 bond_max = bond_mat[i][j];
931 if (bond_mat[i][j] < bond_min)
933 bond_min = bond_mat[i][j];
936 else
938 bond_mat[i][j] = bond_mat[j][i];
943 if (bFile2)
945 rmsd_avg /= tel_mat*tel_mat2;
947 else
949 rmsd_avg /= tel_mat*(tel_mat - 1)/2;
951 if (bMat && (avl > 0))
953 rmsd_max = 0.0;
954 rmsd_min = 0.0;
955 rmsd_avg = 0.0;
956 for (j = 0; j < tel_mat-1; j++)
958 for (i = j+1; i < tel_mat; i++)
960 av_tot = 0;
961 weight_tot = 0;
962 for (my = -avl; my <= avl; my++)
964 if ((j+my >= 0) && (j+my < tel_mat))
966 abs_my = std::abs(my);
967 for (mx = -avl; mx <= avl; mx++)
969 if ((i+mx >= 0) && (i+mx < tel_mat))
971 weight = avl+1.0-std::max(std::abs(mx), abs_my);
972 av_tot += weight*rmsd_mat[i+mx][j+my];
973 weight_tot += weight;
978 rmsdav_mat[i][j] = av_tot/weight_tot;
979 rmsdav_mat[j][i] = rmsdav_mat[i][j];
980 if (rmsdav_mat[i][j] > rmsd_max)
982 rmsd_max = rmsdav_mat[i][j];
986 rmsd_mat = rmsdav_mat;
989 if (bMat)
991 fprintf(stderr, "\n%s: Min %f, Max %f, Avg %f\n",
992 whatname[ewhat], rmsd_min, rmsd_max, rmsd_avg);
993 rlo.r = 1; rlo.g = 1; rlo.b = 1;
994 rhi.r = 0; rhi.g = 0; rhi.b = 0;
995 if (rmsd_user_max != -1)
997 rmsd_max = rmsd_user_max;
999 if (rmsd_user_min != -1)
1001 rmsd_min = rmsd_user_min;
1003 if ((rmsd_user_max != -1) || (rmsd_user_min != -1))
1005 fprintf(stderr, "Min and Max value set to resp. %f and %f\n",
1006 rmsd_min, rmsd_max);
1008 sprintf(buf, "%s %s matrix", gn_rms[0], whatname[ewhat]);
1009 write_xpm(opt2FILE("-m", NFILE, fnm, "w"), 0, buf, whatlabel[ewhat],
1010 output_env_get_time_label(oenv), output_env_get_time_label(oenv), tel_mat, tel_mat2,
1011 axis, axis2, rmsd_mat, rmsd_min, rmsd_max, rlo, rhi, &nlevels);
1012 /* Print the distribution of RMSD values */
1013 if (opt2bSet("-dist", NFILE, fnm))
1015 low_rmsd_dist(opt2fn("-dist", NFILE, fnm), rmsd_max, tel_mat, rmsd_mat, oenv);
1018 if (bDelta)
1020 snew(delta_tot, delta_xsize);
1021 for (j = 0; j < tel_mat-1; j++)
1023 for (i = j+1; i < tel_mat; i++)
1025 mx = i-j;
1026 if (mx < tel_mat/2)
1028 if (bDeltaLog)
1030 mx = static_cast<int>(std::log(static_cast<real>(mx))*delta_scalex+0.5);
1032 my = static_cast<int>(rmsd_mat[i][j]*delta_scaley*del_lev+0.5);
1033 delta_tot[mx] += 1.0;
1034 if ((rmsd_mat[i][j] >= 0) && (rmsd_mat[i][j] <= delta_maxy))
1036 delta[mx][my] += 1.0;
1041 delta_max = 0;
1042 for (i = 0; i < delta_xsize; i++)
1044 if (delta_tot[i] > 0.0)
1046 delta_tot[i] = 1.0/delta_tot[i];
1047 for (j = 0; j <= del_lev; j++)
1049 delta[i][j] *= delta_tot[i];
1050 if (delta[i][j] > delta_max)
1052 delta_max = delta[i][j];
1057 fprintf(stderr, "Maximum in delta matrix: %f\n", delta_max);
1058 snew(del_xaxis, delta_xsize);
1059 snew(del_yaxis, del_lev+1);
1060 for (i = 0; i < delta_xsize; i++)
1062 del_xaxis[i] = axis[i]-axis[0];
1064 for (i = 0; i < del_lev+1; i++)
1066 del_yaxis[i] = delta_maxy*i/del_lev;
1068 sprintf(buf, "%s %s vs. delta t", gn_rms[0], whatname[ewhat]);
1069 fp = gmx_ffopen("delta.xpm", "w");
1070 write_xpm(fp, 0, buf, "density", output_env_get_time_label(oenv), whatlabel[ewhat],
1071 delta_xsize, del_lev+1, del_xaxis, del_yaxis,
1072 delta, 0.0, delta_max, rlo, rhi, &nlevels);
1073 gmx_ffclose(fp);
1075 if (opt2bSet("-bin", NFILE, fnm))
1077 /* NB: File must be binary if we use fwrite */
1078 fp = ftp2FILE(efDAT, NFILE, fnm, "wb");
1079 for (i = 0; i < tel_mat; i++)
1081 if (static_cast<int>(fwrite(rmsd_mat[i], sizeof(**rmsd_mat), tel_mat2, fp)) != tel_mat2)
1083 gmx_fatal(FARGS, "Error writing to output file");
1086 gmx_ffclose(fp);
1089 if (bBond)
1091 fprintf(stderr, "\nMin. angle: %f, Max. angle: %f\n", bond_min, bond_max);
1092 if (bond_user_max != -1)
1094 bond_max = bond_user_max;
1096 if (bond_user_min != -1)
1098 bond_min = bond_user_min;
1100 if ((bond_user_max != -1) || (bond_user_min != -1))
1102 fprintf(stderr, "Bond angle Min and Max set to:\n"
1103 "Min. angle: %f, Max. angle: %f\n", bond_min, bond_max);
1105 rlo.r = 1; rlo.g = 1; rlo.b = 1;
1106 rhi.r = 0; rhi.g = 0; rhi.b = 0;
1107 sprintf(buf, "%s av. bond angle deviation", gn_rms[0]);
1108 write_xpm(opt2FILE("-bm", NFILE, fnm, "w"), 0, buf, "degrees",
1109 output_env_get_time_label(oenv), output_env_get_time_label(oenv), tel_mat, tel_mat2,
1110 axis, axis2, bond_mat, bond_min, bond_max, rlo, rhi, &nlevels);
1114 bAv = opt2bSet("-a", NFILE, fnm);
1116 /* Write the RMSD's to file */
1117 if (!bPrev)
1119 sprintf(buf, "%s", whatxvgname[ewhat]);
1121 else
1123 sprintf(buf, "%s with frame %g %s ago", whatxvgname[ewhat],
1124 time[prev*freq]-time[0], output_env_get_time_label(oenv));
1126 fp = xvgropen(opt2fn("-o", NFILE, fnm), buf, output_env_get_xvgr_tlabel(oenv),
1127 whatxvglabel[ewhat], oenv);
1128 if (output_env_get_print_xvgr_codes(oenv))
1130 fprintf(fp, "@ subtitle \"%s%s after %s%s%s\"\n",
1131 (nrms == 1) ? "" : "of ", gn_rms[0], fitgraphlabel[efit],
1132 bFit ? " to " : "", bFit ? gn_fit : "");
1134 if (nrms != 1)
1136 xvgr_legend(fp, nrms, (const char**)gn_rms, oenv);
1138 for (i = 0; (i < teller); i++)
1140 if (bSplit && i > 0 &&
1141 std::abs(time[bPrev ? freq*i : i]/output_env_get_time_factor(oenv)) < 1e-5)
1143 fprintf(fp, "%s\n", output_env_get_print_xvgr_codes(oenv) ? "&" : "");
1145 fprintf(fp, "%12.7f", time[bPrev ? freq*i : i]);
1146 for (j = 0; (j < nrms); j++)
1148 fprintf(fp, " %12.7f", rls[j][i]);
1149 if (bAv)
1151 rlstot += rls[j][i];
1154 fprintf(fp, "\n");
1156 xvgrclose(fp);
1158 if (bMirror)
1160 /* Write the mirror RMSD's to file */
1161 sprintf(buf, "%s with Mirror", whatxvgname[ewhat]);
1162 sprintf(buf2, "Mirror %s", whatxvglabel[ewhat]);
1163 fp = xvgropen(opt2fn("-mir", NFILE, fnm), buf, output_env_get_xvgr_tlabel(oenv),
1164 buf2, oenv);
1165 if (nrms == 1)
1167 if (output_env_get_print_xvgr_codes(oenv))
1169 fprintf(fp, "@ subtitle \"of %s after lsq fit to mirror of %s\"\n",
1170 gn_rms[0], bFit ? gn_fit : "");
1173 else
1175 if (output_env_get_print_xvgr_codes(oenv))
1177 fprintf(fp, "@ subtitle \"after lsq fit to mirror %s\"\n", bFit ? gn_fit : "");
1179 xvgr_legend(fp, nrms, (const char**)gn_rms, oenv);
1181 for (i = 0; (i < teller); i++)
1183 if (bSplit && i > 0 && std::abs(time[i]) < 1e-5)
1185 fprintf(fp, "%s\n", output_env_get_print_xvgr_codes(oenv) ? "&" : "");
1187 fprintf(fp, "%12.7f", time[i]);
1188 for (j = 0; (j < nrms); j++)
1190 fprintf(fp, " %12.7f", rlsm[j][i]);
1192 fprintf(fp, "\n");
1194 xvgrclose(fp);
1197 if (bAv)
1199 sprintf(buf, "Average %s", whatxvgname[ewhat]);
1200 sprintf(buf2, "Average %s", whatxvglabel[ewhat]);
1201 fp = xvgropen(opt2fn("-a", NFILE, fnm), buf, "Residue", buf2, oenv);
1202 for (j = 0; (j < nrms); j++)
1204 fprintf(fp, "%10d %10g\n", j, rlstot/teller);
1206 xvgrclose(fp);
1209 if (bNorm)
1211 fp = xvgropen("aver.xvg", gn_rms[0], "Residue", whatxvglabel[ewhat], oenv);
1212 for (j = 0; (j < irms[0]); j++)
1214 fprintf(fp, "%10d %10g\n", j, rlsnorm[j]/teller);
1216 xvgrclose(fp);
1218 do_view(oenv, opt2fn_null("-a", NFILE, fnm), "-graphtype bar");
1219 do_view(oenv, opt2fn("-o", NFILE, fnm), NULL);
1220 do_view(oenv, opt2fn_null("-mir", NFILE, fnm), NULL);
1221 do_view(oenv, opt2fn_null("-m", NFILE, fnm), NULL);
1222 do_view(oenv, opt2fn_null("-bm", NFILE, fnm), NULL);
1223 do_view(oenv, opt2fn_null("-dist", NFILE, fnm), NULL);
1225 return 0;