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,2017,2018 by the GROMACS development team.
7 * Copyright (c) 2019,2020, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
45 #include "gromacs/commandline/filenm.h"
46 #include "gromacs/commandline/pargs.h"
47 #include "gromacs/fileio/confio.h"
48 #include "gromacs/fileio/matio.h"
49 #include "gromacs/fileio/trxio.h"
50 #include "gromacs/fileio/xvgr.h"
51 #include "gromacs/gmxana/gmx_ana.h"
52 #include "gromacs/math/functions.h"
53 #include "gromacs/math/vec.h"
54 #include "gromacs/pbcutil/pbc.h"
55 #include "gromacs/pbcutil/rmpbc.h"
56 #include "gromacs/topology/index.h"
57 #include "gromacs/topology/topology.h"
58 #include "gromacs/utility/arraysize.h"
59 #include "gromacs/utility/cstringutil.h"
60 #include "gromacs/utility/fatalerror.h"
61 #include "gromacs/utility/futil.h"
62 #include "gromacs/utility/smalloc.h"
67 static int* res_ndx(t_atoms* atoms)
76 snew(rndx, atoms->nr);
77 r0 = atoms->atom[0].resind;
78 for (i = 0; (i < atoms->nr); i++)
80 rndx[i] = atoms->atom[i].resind - r0;
86 static int* res_natm(t_atoms* atoms)
95 snew(natm, atoms->nres);
96 r0 = atoms->atom[0].resind;
98 for (i = 0; (i < atoms->nres); i++)
100 while ((atoms->atom[j].resind) - r0 == i)
110 static void calc_mat(int nres,
121 int i, j, resi, resj;
126 set_pbc(&pbc, pbcType, box);
127 trunc2 = gmx::square(trunc);
128 for (resi = 0; (resi < nres); resi++)
130 for (resj = 0; (resj < nres); resj++)
132 mdmat[resi][resj] = FARAWAY;
135 for (i = 0; (i < natoms); i++)
138 for (j = i + 1; (j < natoms); j++)
141 pbc_dx(&pbc, x[index[i]], x[index[j]], ddx);
148 mdmat[resi][resj] = std::min(r2, mdmat[resi][resj]);
152 for (resi = 0; (resi < nres); resi++)
154 mdmat[resi][resi] = 0;
155 for (resj = resi + 1; (resj < nres); resj++)
157 r = std::sqrt(mdmat[resi][resj]);
158 mdmat[resi][resj] = r;
159 mdmat[resj][resi] = r;
164 static void tot_nmat(int nres, int natoms, int nframes, int** nmat, int* tot_n, real* mean_n)
168 for (i = 0; (i < nres); i++)
170 for (j = 0; (j < natoms); j++)
175 mean_n[i] += nmat[i][j];
178 mean_n[i] /= nframes;
182 int gmx_mdmat(int argc, char* argv[])
184 const char* desc[] = {
185 "[THISMODULE] makes distance matrices consisting of the smallest distance",
186 "between residue pairs. With [TT]-frames[tt], these distance matrices can be",
187 "stored in order to see differences in tertiary structure as a",
188 "function of time. If you choose your options unwisely, this may generate",
189 "a large output file. By default, only an averaged matrix over the whole",
190 "trajectory is output.",
191 "Also a count of the number of different atomic contacts between",
192 "residues over the whole trajectory can be made.",
193 "The output can be processed with [gmx-xpm2ps] to make a PostScript (tm) plot."
195 static real truncate = 1.5;
196 static int nlevels = 40;
198 { "-t", FALSE, etREAL, { &truncate }, "trunc distance" },
199 { "-nlevels", FALSE, etINT, { &nlevels }, "Discretize distance in this number of levels" }
202 { efTRX, "-f", nullptr, ffREAD }, { efTPS, nullptr, nullptr, ffREAD },
203 { efNDX, nullptr, nullptr, ffOPTRD }, { efXPM, "-mean", "dm", ffWRITE },
204 { efXPM, "-frames", "dmf", ffOPTWR }, { efXVG, "-no", "num", ffOPTWR },
206 #define NFILE asize(fnm)
208 FILE * out = nullptr, *fp;
215 int * rndx, *natm, prevres, newres;
217 int i, j, nres, natoms, nframes, trxnat;
219 gmx_bool bCalcN, bFrames;
224 real ** mdmat, *resnr, **totmdmat;
225 int ** nmat, **totnmat;
228 matrix box = { { 0 } };
229 gmx_output_env_t* oenv;
230 gmx_rmpbc_t gpbc = nullptr;
232 if (!parse_common_args(&argc, argv, PCA_CAN_TIME, NFILE, fnm, asize(pa), pa, asize(desc), desc,
238 fprintf(stderr, "Will truncate at %f nm\n", truncate);
239 bCalcN = opt2bSet("-no", NFILE, fnm);
240 bFrames = opt2bSet("-frames", NFILE, fnm);
243 fprintf(stderr, "Will calculate number of different contacts\n");
246 read_tps_conf(ftp2fn(efTPS, NFILE, fnm), &top, &pbcType, &x, nullptr, box, FALSE);
248 fprintf(stderr, "Select group for analysis\n");
249 get_index(&top.atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &isize, &index, &grpname);
252 snew(useatoms.atom, natoms);
253 snew(useatoms.atomname, natoms);
256 snew(useatoms.resinfo, natoms);
258 prevres = top.atoms.atom[index[0]].resind;
260 for (i = 0; (i < isize); i++)
263 useatoms.atomname[i] = top.atoms.atomname[ii];
264 if (top.atoms.atom[ii].resind != prevres)
266 prevres = top.atoms.atom[ii].resind;
268 useatoms.resinfo[i] = top.atoms.resinfo[prevres];
271 fprintf(debug, "New residue: atom %5s %5s %6d, index entry %5d, newres %5d\n",
272 *(top.atoms.resinfo[top.atoms.atom[ii].resind].name),
273 *(top.atoms.atomname[ii]), ii, i, newres);
276 useatoms.atom[i].resind = newres;
278 useatoms.nres = newres + 1;
281 rndx = res_ndx(&(useatoms));
282 natm = res_natm(&(useatoms));
283 nres = useatoms.nres;
284 fprintf(stderr, "There are %d residues with %d atoms\n", nres, natoms);
292 for (i = 0; (i < nres); i++)
294 snew(mdmat[i], nres);
295 snew(nmat[i], natoms);
296 snew(totnmat[i], natoms);
299 snew(totmdmat, nres);
300 for (i = 0; (i < nres); i++)
302 snew(totmdmat[i], nres);
305 trxnat = read_first_x(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &t, &x, box);
316 gpbc = gmx_rmpbc_init(&top.idef, pbcType, trxnat);
320 out = opt2FILE("-frames", NFILE, fnm, "w");
324 gmx_rmpbc(gpbc, trxnat, box, x);
326 calc_mat(nres, natoms, rndx, x, index, truncate, mdmat, nmat, pbcType, box);
327 for (i = 0; (i < nres); i++)
329 for (j = 0; (j < natoms); j++)
337 for (i = 0; (i < nres); i++)
339 for (j = 0; (j < nres); j++)
341 totmdmat[i][j] += mdmat[i][j];
346 sprintf(label, "t=%.0f ps", t);
347 write_xpm(out, 0, label, "Distance (nm)", "Residue Index", "Residue Index", nres, nres,
348 resnr, resnr, mdmat, 0, truncate, rlo, rhi, &nlevels);
350 } while (read_next_x(oenv, status, &t, x, box));
351 fprintf(stderr, "\n");
353 gmx_rmpbc_done(gpbc);
359 fprintf(stderr, "Processed %d frames\n", nframes);
361 for (i = 0; (i < nres); i++)
363 for (j = 0; (j < nres); j++)
365 totmdmat[i][j] /= nframes;
368 write_xpm(opt2FILE("-mean", NFILE, fnm, "w"), 0, "Mean smallest distance", "Distance (nm)",
369 "Residue Index", "Residue Index", nres, nres, resnr, resnr, totmdmat, 0, truncate,
377 for (i = 0; i < 5; i++)
379 snew(legend[i], STRLEN);
381 tot_nmat(nres, natoms, nframes, totnmat, tot_n, mean_n);
382 fp = xvgropen(ftp2fn(efXVG, NFILE, fnm), "Increase in number of contacts", "Residue",
384 sprintf(legend[0], "Total/mean");
385 sprintf(legend[1], "Total");
386 sprintf(legend[2], "Mean");
387 sprintf(legend[3], "# atoms");
388 sprintf(legend[4], "Mean/# atoms");
389 xvgr_legend(fp, 5, legend, oenv);
390 for (i = 0; (i < nres); i++)
398 ratio = tot_n[i] / mean_n[i];
400 fprintf(fp, "%3d %8.3f %3d %8.3f %3d %8.3f\n", i + 1, ratio, tot_n[i], mean_n[i],
401 natm[i], mean_n[i] / natm[i]);