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"
62 int *res_ndx(t_atoms *atoms)
71 snew(rndx, atoms->nr);
72 r0 = atoms->atom[0].resind;
73 for (i = 0; (i < atoms->nr); i++)
75 rndx[i] = atoms->atom[i].resind-r0;
81 int *res_natm(t_atoms *atoms)
90 snew(natm, atoms->nres);
91 r0 = atoms->atom[0].resind;
93 for (i = 0; (i < atoms->nres); i++)
95 while ((atoms->atom[j].resind)-r0 == i)
105 static void calc_mat(int nres, int natoms, int rndx[],
106 rvec x[], atom_id *index,
107 real trunc, real **mdmat, int **nmat, int ePBC, matrix box)
109 int i, j, resi, resj;
114 set_pbc(&pbc, ePBC, box);
116 for (resi = 0; (resi < nres); resi++)
118 for (resj = 0; (resj < nres); resj++)
120 mdmat[resi][resj] = FARAWAY;
123 for (i = 0; (i < natoms); i++)
126 for (j = i+1; (j < natoms); j++)
129 pbc_dx(&pbc, x[index[i]], x[index[j]], ddx);
136 mdmat[resi][resj] = min(r2, mdmat[resi][resj]);
140 for (resi = 0; (resi < nres); resi++)
142 mdmat[resi][resi] = 0;
143 for (resj = resi+1; (resj < nres); resj++)
145 r = sqrt(mdmat[resi][resj]);
146 mdmat[resi][resj] = r;
147 mdmat[resj][resi] = r;
152 static void tot_nmat(int nres, int natoms, int nframes, int **nmat,
153 int *tot_n, real *mean_n)
157 for (i = 0; (i < nres); i++)
159 for (j = 0; (j < natoms); j++)
164 mean_n[i] += nmat[i][j];
167 mean_n[i] /= nframes;
171 int gmx_mdmat(int argc, char *argv[])
173 const char *desc[] = {
174 "[TT]g_mdmat[tt] makes distance matrices consisting of the smallest distance",
175 "between residue pairs. With [TT]-frames[tt], these distance matrices can be",
176 "stored in order to see differences in tertiary structure as a",
177 "function of time. If you choose your options unwisely, this may generate",
178 "a large output file. By default, only an averaged matrix over the whole",
179 "trajectory is output.",
180 "Also a count of the number of different atomic contacts between",
181 "residues over the whole trajectory can be made.",
182 "The output can be processed with [TT]xpm2ps[tt] to make a PostScript (tm) plot."
184 static real truncate = 1.5;
185 static gmx_bool bAtom = FALSE;
186 static int nlevels = 40;
188 { "-t", FALSE, etREAL, {&truncate},
190 { "-nlevels", FALSE, etINT, {&nlevels},
191 "Discretize distance in this number of levels" }
194 { efTRX, "-f", NULL, ffREAD },
195 { efTPS, NULL, NULL, ffREAD },
196 { efNDX, NULL, NULL, ffOPTRD },
197 { efXPM, "-mean", "dm", ffWRITE },
198 { efXPM, "-frames", "dmf", ffOPTWR },
199 { efXVG, "-no", "num", ffOPTWR },
201 #define NFILE asize(fnm)
203 FILE *out = NULL, *fp;
210 int *rndx, *natm, prevres, newres;
212 int i, j, nres, natoms, nframes, it, trxnat;
215 gmx_bool bCalcN, bFrames;
217 char title[256], label[234];
220 real **mdmat, *resnr, **totmdmat;
221 int **nmat, **totnmat;
226 gmx_rmpbc_t gpbc = NULL;
228 parse_common_args(&argc, argv, PCA_CAN_TIME | PCA_BE_NICE, NFILE, fnm,
229 asize(pa), pa, asize(desc), desc, 0, NULL, &oenv);
231 fprintf(stderr, "Will truncate at %f nm\n", truncate);
232 bCalcN = opt2bSet("-no", NFILE, fnm);
233 bFrames = opt2bSet("-frames", NFILE, fnm);
236 fprintf(stderr, "Will calculate number of different contacts\n");
239 read_tps_conf(ftp2fn(efTPS, NFILE, fnm), title, &top, &ePBC, &x, NULL, box, FALSE);
241 fprintf(stderr, "Select group for analysis\n");
242 get_index(&top.atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &isize, &index, &grpname);
245 snew(useatoms.atom, natoms);
246 snew(useatoms.atomname, natoms);
249 snew(useatoms.resinfo, natoms);
251 prevres = top.atoms.atom[index[0]].resind;
253 for (i = 0; (i < isize); i++)
256 useatoms.atomname[i] = top.atoms.atomname[ii];
257 if (top.atoms.atom[ii].resind != prevres)
259 prevres = top.atoms.atom[ii].resind;
261 useatoms.resinfo[i] = top.atoms.resinfo[prevres];
264 fprintf(debug, "New residue: atom %5s %5s %6d, index entry %5d, newres %5d\n",
265 *(top.atoms.resinfo[top.atoms.atom[ii].resind].name),
266 *(top.atoms.atomname[ii]),
270 useatoms.atom[i].resind = newres;
272 useatoms.nres = newres+1;
275 rndx = res_ndx(&(useatoms));
276 natm = res_natm(&(useatoms));
277 nres = useatoms.nres;
278 fprintf(stderr, "There are %d residues with %d atoms\n", nres, natoms);
286 for (i = 0; (i < nres); i++)
288 snew(mdmat[i], nres);
289 snew(nmat[i], natoms);
290 snew(totnmat[i], natoms);
293 snew(totmdmat, nres);
294 for (i = 0; (i < nres); i++)
296 snew(totmdmat[i], nres);
299 trxnat = read_first_x(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &t, &x, box);
303 rlo.r = 1.0, rlo.g = 1.0, rlo.b = 1.0;
304 rhi.r = 0.0, rhi.g = 0.0, rhi.b = 0.0;
306 gpbc = gmx_rmpbc_init(&top.idef, ePBC, trxnat, box);
310 out = opt2FILE("-frames", NFILE, fnm, "w");
314 gmx_rmpbc(gpbc, trxnat, box, x);
316 calc_mat(nres, natoms, rndx, x, index, truncate, mdmat, nmat, ePBC, box);
317 for (i = 0; (i < nres); i++)
319 for (j = 0; (j < natoms); j++)
327 for (i = 0; (i < nres); i++)
329 for (j = 0; (j < nres); j++)
331 totmdmat[i][j] += mdmat[i][j];
336 sprintf(label, "t=%.0f ps", t);
337 write_xpm(out, 0, label, "Distance (nm)", "Residue Index", "Residue Index",
338 nres, nres, resnr, resnr, mdmat, 0, truncate, rlo, rhi, &nlevels);
341 while (read_next_x(oenv, status, &t, trxnat, x, box));
342 fprintf(stderr, "\n");
344 gmx_rmpbc_done(gpbc);
350 fprintf(stderr, "Processed %d frames\n", nframes);
352 for (i = 0; (i < nres); i++)
354 for (j = 0; (j < nres); j++)
356 totmdmat[i][j] /= nframes;
359 write_xpm(opt2FILE("-mean", NFILE, fnm, "w"), 0, "Mean smallest distance",
360 "Distance (nm)", "Residue Index", "Residue Index",
361 nres, nres, resnr, resnr, totmdmat, 0, truncate, rlo, rhi, &nlevels);
365 tot_nmat(nres, natoms, nframes, totnmat, tot_n, mean_n);
366 fp = xvgropen(ftp2fn(efXVG, NFILE, fnm),
367 "Increase in number of contacts", "Residue", "Ratio", oenv);
368 fprintf(fp, "@ legend on\n");
369 fprintf(fp, "@ legend box on\n");
370 fprintf(fp, "@ legend loctype view\n");
371 fprintf(fp, "@ legend 0.75, 0.8\n");
372 fprintf(fp, "@ legend string 0 \"Total/mean\"\n");
373 fprintf(fp, "@ legend string 1 \"Total\"\n");
374 fprintf(fp, "@ legend string 2 \"Mean\"\n");
375 fprintf(fp, "@ legend string 3 \"# atoms\"\n");
376 fprintf(fp, "@ legend string 4 \"Mean/# atoms\"\n");
377 fprintf(fp, "#%3s %8s %3s %8s %3s %8s\n",
378 "res", "ratio", "tot", "mean", "natm", "mean/atm");
379 for (i = 0; (i < nres); i++)
387 ratio = tot_n[i]/mean_n[i];
389 fprintf(fp, "%3d %8.3f %3d %8.3f %3d %8.3f\n",
390 i+1, ratio, tot_n[i], mean_n[i], natm[i], mean_n[i]/natm[i]);