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
46 #include "gromacs/fileio/filenm.h"
48 #include "gromacs/fileio/futil.h"
49 #include "gmx_fatal.h"
51 #include "gromacs/fileio/matio.h"
54 #include "gromacs/fileio/tpxio.h"
55 #include "gromacs/fileio/trxio.h"
63 int *res_ndx(t_atoms *atoms)
72 snew(rndx, atoms->nr);
73 r0 = atoms->atom[0].resind;
74 for (i = 0; (i < atoms->nr); i++)
76 rndx[i] = atoms->atom[i].resind-r0;
82 int *res_natm(t_atoms *atoms)
91 snew(natm, atoms->nres);
92 r0 = atoms->atom[0].resind;
94 for (i = 0; (i < atoms->nres); i++)
96 while ((atoms->atom[j].resind)-r0 == i)
106 static void calc_mat(int nres, int natoms, int rndx[],
107 rvec x[], atom_id *index,
108 real trunc, real **mdmat, int **nmat, int ePBC, matrix box)
110 int i, j, resi, resj;
115 set_pbc(&pbc, ePBC, box);
117 for (resi = 0; (resi < nres); resi++)
119 for (resj = 0; (resj < nres); resj++)
121 mdmat[resi][resj] = FARAWAY;
124 for (i = 0; (i < natoms); i++)
127 for (j = i+1; (j < natoms); j++)
130 pbc_dx(&pbc, x[index[i]], x[index[j]], ddx);
137 mdmat[resi][resj] = min(r2, mdmat[resi][resj]);
141 for (resi = 0; (resi < nres); resi++)
143 mdmat[resi][resi] = 0;
144 for (resj = resi+1; (resj < nres); resj++)
146 r = sqrt(mdmat[resi][resj]);
147 mdmat[resi][resj] = r;
148 mdmat[resj][resi] = r;
153 static void tot_nmat(int nres, int natoms, int nframes, int **nmat,
154 int *tot_n, real *mean_n)
158 for (i = 0; (i < nres); i++)
160 for (j = 0; (j < natoms); j++)
165 mean_n[i] += nmat[i][j];
168 mean_n[i] /= nframes;
172 int gmx_mdmat(int argc, char *argv[])
174 const char *desc[] = {
175 "[TT]g_mdmat[tt] makes distance matrices consisting of the smallest distance",
176 "between residue pairs. With [TT]-frames[tt], these distance matrices can be",
177 "stored in order to see differences in tertiary structure as a",
178 "function of time. If you choose your options unwisely, this may generate",
179 "a large output file. By default, only an averaged matrix over the whole",
180 "trajectory is output.",
181 "Also a count of the number of different atomic contacts between",
182 "residues over the whole trajectory can be made.",
183 "The output can be processed with [TT]xpm2ps[tt] to make a PostScript (tm) plot."
185 static real truncate = 1.5;
186 static gmx_bool bAtom = FALSE;
187 static int nlevels = 40;
189 { "-t", FALSE, etREAL, {&truncate},
191 { "-nlevels", FALSE, etINT, {&nlevels},
192 "Discretize distance in this number of levels" }
195 { efTRX, "-f", NULL, ffREAD },
196 { efTPS, NULL, NULL, ffREAD },
197 { efNDX, NULL, NULL, ffOPTRD },
198 { efXPM, "-mean", "dm", ffWRITE },
199 { efXPM, "-frames", "dmf", ffOPTWR },
200 { efXVG, "-no", "num", ffOPTWR },
202 #define NFILE asize(fnm)
204 FILE *out = NULL, *fp;
211 int *rndx, *natm, prevres, newres;
213 int i, j, nres, natoms, nframes, it, trxnat;
216 gmx_bool bCalcN, bFrames;
218 char title[256], label[234];
221 real **mdmat, *resnr, **totmdmat;
222 int **nmat, **totnmat;
227 gmx_rmpbc_t gpbc = NULL;
229 if (!parse_common_args(&argc, argv, PCA_CAN_TIME | PCA_BE_NICE, NFILE, fnm,
230 asize(pa), pa, asize(desc), desc, 0, NULL, &oenv))
235 fprintf(stderr, "Will truncate at %f nm\n", truncate);
236 bCalcN = opt2bSet("-no", NFILE, fnm);
237 bFrames = opt2bSet("-frames", NFILE, fnm);
240 fprintf(stderr, "Will calculate number of different contacts\n");
243 read_tps_conf(ftp2fn(efTPS, NFILE, fnm), title, &top, &ePBC, &x, NULL, box, FALSE);
245 fprintf(stderr, "Select group for analysis\n");
246 get_index(&top.atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &isize, &index, &grpname);
249 snew(useatoms.atom, natoms);
250 snew(useatoms.atomname, natoms);
253 snew(useatoms.resinfo, natoms);
255 prevres = top.atoms.atom[index[0]].resind;
257 for (i = 0; (i < isize); i++)
260 useatoms.atomname[i] = top.atoms.atomname[ii];
261 if (top.atoms.atom[ii].resind != prevres)
263 prevres = top.atoms.atom[ii].resind;
265 useatoms.resinfo[i] = top.atoms.resinfo[prevres];
268 fprintf(debug, "New residue: atom %5s %5s %6d, index entry %5d, newres %5d\n",
269 *(top.atoms.resinfo[top.atoms.atom[ii].resind].name),
270 *(top.atoms.atomname[ii]),
274 useatoms.atom[i].resind = newres;
276 useatoms.nres = newres+1;
279 rndx = res_ndx(&(useatoms));
280 natm = res_natm(&(useatoms));
281 nres = useatoms.nres;
282 fprintf(stderr, "There are %d residues with %d atoms\n", nres, natoms);
290 for (i = 0; (i < nres); i++)
292 snew(mdmat[i], nres);
293 snew(nmat[i], natoms);
294 snew(totnmat[i], natoms);
297 snew(totmdmat, nres);
298 for (i = 0; (i < nres); i++)
300 snew(totmdmat[i], nres);
303 trxnat = read_first_x(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &t, &x, box);
307 rlo.r = 1.0, rlo.g = 1.0, rlo.b = 1.0;
308 rhi.r = 0.0, rhi.g = 0.0, rhi.b = 0.0;
310 gpbc = gmx_rmpbc_init(&top.idef, ePBC, trxnat);
314 out = opt2FILE("-frames", NFILE, fnm, "w");
318 gmx_rmpbc(gpbc, trxnat, box, x);
320 calc_mat(nres, natoms, rndx, x, index, truncate, mdmat, nmat, ePBC, box);
321 for (i = 0; (i < nres); i++)
323 for (j = 0; (j < natoms); j++)
331 for (i = 0; (i < nres); i++)
333 for (j = 0; (j < nres); j++)
335 totmdmat[i][j] += mdmat[i][j];
340 sprintf(label, "t=%.0f ps", t);
341 write_xpm(out, 0, label, "Distance (nm)", "Residue Index", "Residue Index",
342 nres, nres, resnr, resnr, mdmat, 0, truncate, rlo, rhi, &nlevels);
345 while (read_next_x(oenv, status, &t, x, box));
346 fprintf(stderr, "\n");
348 gmx_rmpbc_done(gpbc);
354 fprintf(stderr, "Processed %d frames\n", nframes);
356 for (i = 0; (i < nres); i++)
358 for (j = 0; (j < nres); j++)
360 totmdmat[i][j] /= nframes;
363 write_xpm(opt2FILE("-mean", NFILE, fnm, "w"), 0, "Mean smallest distance",
364 "Distance (nm)", "Residue Index", "Residue Index",
365 nres, nres, resnr, resnr, totmdmat, 0, truncate, rlo, rhi, &nlevels);
369 tot_nmat(nres, natoms, nframes, totnmat, tot_n, mean_n);
370 fp = xvgropen(ftp2fn(efXVG, NFILE, fnm),
371 "Increase in number of contacts", "Residue", "Ratio", oenv);
372 fprintf(fp, "@ legend on\n");
373 fprintf(fp, "@ legend box on\n");
374 fprintf(fp, "@ legend loctype view\n");
375 fprintf(fp, "@ legend 0.75, 0.8\n");
376 fprintf(fp, "@ legend string 0 \"Total/mean\"\n");
377 fprintf(fp, "@ legend string 1 \"Total\"\n");
378 fprintf(fp, "@ legend string 2 \"Mean\"\n");
379 fprintf(fp, "@ legend string 3 \"# atoms\"\n");
380 fprintf(fp, "@ legend string 4 \"Mean/# atoms\"\n");
381 fprintf(fp, "#%3s %8s %3s %8s %3s %8s\n",
382 "res", "ratio", "tot", "mean", "natm", "mean/atm");
383 for (i = 0; (i < nres); i++)
391 ratio = tot_n[i]/mean_n[i];
393 fprintf(fp, "%3d %8.3f %3d %8.3f %3d %8.3f\n",
394 i+1, ratio, tot_n[i], mean_n[i], natm[i], mean_n[i]/natm[i]);