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, 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.
44 #include "gromacs/commandline/filenm.h"
45 #include "gromacs/commandline/pargs.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/gmx_ana.h"
51 #include "gromacs/math/functions.h"
52 #include "gromacs/math/vec.h"
53 #include "gromacs/pbcutil/pbc.h"
54 #include "gromacs/pbcutil/rmpbc.h"
55 #include "gromacs/topology/index.h"
56 #include "gromacs/topology/topology.h"
57 #include "gromacs/utility/arraysize.h"
58 #include "gromacs/utility/cstringutil.h"
59 #include "gromacs/utility/fatalerror.h"
60 #include "gromacs/utility/futil.h"
61 #include "gromacs/utility/smalloc.h"
66 static int *res_ndx(t_atoms
*atoms
)
75 snew(rndx
, atoms
->nr
);
76 r0
= atoms
->atom
[0].resind
;
77 for (i
= 0; (i
< atoms
->nr
); i
++)
79 rndx
[i
] = atoms
->atom
[i
].resind
-r0
;
85 static int *res_natm(t_atoms
*atoms
)
94 snew(natm
, atoms
->nres
);
95 r0
= atoms
->atom
[0].resind
;
97 for (i
= 0; (i
< atoms
->nres
); i
++)
99 while ((atoms
->atom
[j
].resind
)-r0
== i
)
109 static void calc_mat(int nres
, int natoms
, const int rndx
[],
110 rvec x
[], const int *index
,
111 real trunc
, real
**mdmat
, int **nmat
, int ePBC
, matrix box
)
113 int i
, j
, resi
, resj
;
118 set_pbc(&pbc
, ePBC
, box
);
119 trunc2
= gmx::square(trunc
);
120 for (resi
= 0; (resi
< nres
); resi
++)
122 for (resj
= 0; (resj
< nres
); resj
++)
124 mdmat
[resi
][resj
] = FARAWAY
;
127 for (i
= 0; (i
< natoms
); i
++)
130 for (j
= i
+1; (j
< natoms
); j
++)
133 pbc_dx(&pbc
, x
[index
[i
]], x
[index
[j
]], ddx
);
140 mdmat
[resi
][resj
] = std::min(r2
, mdmat
[resi
][resj
]);
144 for (resi
= 0; (resi
< nres
); resi
++)
146 mdmat
[resi
][resi
] = 0;
147 for (resj
= resi
+1; (resj
< nres
); resj
++)
149 r
= std::sqrt(mdmat
[resi
][resj
]);
150 mdmat
[resi
][resj
] = r
;
151 mdmat
[resj
][resi
] = r
;
156 static void tot_nmat(int nres
, int natoms
, int nframes
, int **nmat
,
157 int *tot_n
, real
*mean_n
)
161 for (i
= 0; (i
< nres
); i
++)
163 for (j
= 0; (j
< natoms
); j
++)
168 mean_n
[i
] += nmat
[i
][j
];
171 mean_n
[i
] /= nframes
;
175 int gmx_mdmat(int argc
, char *argv
[])
177 const char *desc
[] = {
178 "[THISMODULE] makes distance matrices consisting of the smallest distance",
179 "between residue pairs. With [TT]-frames[tt], these distance matrices can be",
180 "stored in order to see differences in tertiary structure as a",
181 "function of time. If you choose your options unwisely, this may generate",
182 "a large output file. By default, only an averaged matrix over the whole",
183 "trajectory is output.",
184 "Also a count of the number of different atomic contacts between",
185 "residues over the whole trajectory can be made.",
186 "The output can be processed with [gmx-xpm2ps] to make a PostScript (tm) plot."
188 static real truncate
= 1.5;
189 static int nlevels
= 40;
191 { "-t", FALSE
, etREAL
, {&truncate
},
193 { "-nlevels", FALSE
, etINT
, {&nlevels
},
194 "Discretize distance in this number of levels" }
197 { efTRX
, "-f", nullptr, ffREAD
},
198 { efTPS
, nullptr, nullptr, ffREAD
},
199 { efNDX
, nullptr, nullptr, ffOPTRD
},
200 { efXPM
, "-mean", "dm", ffWRITE
},
201 { efXPM
, "-frames", "dmf", ffOPTWR
},
202 { efXVG
, "-no", "num", ffOPTWR
},
204 #define NFILE asize(fnm)
206 FILE *out
= nullptr, *fp
;
213 int *rndx
, *natm
, prevres
, newres
;
215 int i
, j
, nres
, natoms
, nframes
, trxnat
;
217 gmx_bool bCalcN
, bFrames
;
222 real
**mdmat
, *resnr
, **totmdmat
;
223 int **nmat
, **totnmat
;
227 gmx_output_env_t
*oenv
;
228 gmx_rmpbc_t gpbc
= nullptr;
230 if (!parse_common_args(&argc
, argv
, PCA_CAN_TIME
, NFILE
, fnm
,
231 asize(pa
), pa
, asize(desc
), desc
, 0, nullptr, &oenv
))
236 fprintf(stderr
, "Will truncate at %f nm\n", truncate
);
237 bCalcN
= opt2bSet("-no", NFILE
, fnm
);
238 bFrames
= opt2bSet("-frames", NFILE
, fnm
);
241 fprintf(stderr
, "Will calculate number of different contacts\n");
244 read_tps_conf(ftp2fn(efTPS
, NFILE
, fnm
), &top
, &ePBC
, &x
, nullptr, box
, FALSE
);
246 fprintf(stderr
, "Select group for analysis\n");
247 get_index(&top
.atoms
, ftp2fn_null(efNDX
, NFILE
, fnm
), 1, &isize
, &index
, &grpname
);
250 snew(useatoms
.atom
, natoms
);
251 snew(useatoms
.atomname
, natoms
);
254 snew(useatoms
.resinfo
, natoms
);
256 prevres
= top
.atoms
.atom
[index
[0]].resind
;
258 for (i
= 0; (i
< isize
); i
++)
261 useatoms
.atomname
[i
] = top
.atoms
.atomname
[ii
];
262 if (top
.atoms
.atom
[ii
].resind
!= prevres
)
264 prevres
= top
.atoms
.atom
[ii
].resind
;
266 useatoms
.resinfo
[i
] = top
.atoms
.resinfo
[prevres
];
269 fprintf(debug
, "New residue: atom %5s %5s %6d, index entry %5d, newres %5d\n",
270 *(top
.atoms
.resinfo
[top
.atoms
.atom
[ii
].resind
].name
),
271 *(top
.atoms
.atomname
[ii
]),
275 useatoms
.atom
[i
].resind
= newres
;
277 useatoms
.nres
= newres
+1;
280 rndx
= res_ndx(&(useatoms
));
281 natm
= res_natm(&(useatoms
));
282 nres
= useatoms
.nres
;
283 fprintf(stderr
, "There are %d residues with %d atoms\n", nres
, natoms
);
291 for (i
= 0; (i
< nres
); i
++)
293 snew(mdmat
[i
], nres
);
294 snew(nmat
[i
], natoms
);
295 snew(totnmat
[i
], natoms
);
298 snew(totmdmat
, nres
);
299 for (i
= 0; (i
< nres
); i
++)
301 snew(totmdmat
[i
], nres
);
304 trxnat
= read_first_x(oenv
, &status
, ftp2fn(efTRX
, NFILE
, fnm
), &t
, &x
, box
);
308 rlo
.r
= 1.0; rlo
.g
= 1.0; rlo
.b
= 1.0;
309 rhi
.r
= 0.0; rhi
.g
= 0.0; rhi
.b
= 0.0;
311 gpbc
= gmx_rmpbc_init(&top
.idef
, ePBC
, trxnat
);
315 out
= opt2FILE("-frames", NFILE
, fnm
, "w");
319 gmx_rmpbc(gpbc
, trxnat
, box
, x
);
321 calc_mat(nres
, natoms
, rndx
, x
, index
, truncate
, mdmat
, nmat
, ePBC
, box
);
322 for (i
= 0; (i
< nres
); i
++)
324 for (j
= 0; (j
< natoms
); j
++)
332 for (i
= 0; (i
< nres
); i
++)
334 for (j
= 0; (j
< nres
); j
++)
336 totmdmat
[i
][j
] += mdmat
[i
][j
];
341 sprintf(label
, "t=%.0f ps", t
);
342 write_xpm(out
, 0, label
, "Distance (nm)", "Residue Index", "Residue Index",
343 nres
, nres
, resnr
, resnr
, mdmat
, 0, truncate
, rlo
, rhi
, &nlevels
);
346 while (read_next_x(oenv
, status
, &t
, x
, box
));
347 fprintf(stderr
, "\n");
349 gmx_rmpbc_done(gpbc
);
355 fprintf(stderr
, "Processed %d frames\n", nframes
);
357 for (i
= 0; (i
< nres
); i
++)
359 for (j
= 0; (j
< nres
); j
++)
361 totmdmat
[i
][j
] /= nframes
;
364 write_xpm(opt2FILE("-mean", NFILE
, fnm
, "w"), 0, "Mean smallest distance",
365 "Distance (nm)", "Residue Index", "Residue Index",
366 nres
, nres
, resnr
, resnr
, totmdmat
, 0, truncate
, rlo
, rhi
, &nlevels
);
373 for (i
= 0; i
< 5; i
++)
375 snew(legend
[i
], STRLEN
);
377 tot_nmat(nres
, natoms
, nframes
, totnmat
, tot_n
, mean_n
);
378 fp
= xvgropen(ftp2fn(efXVG
, NFILE
, fnm
),
379 "Increase in number of contacts", "Residue", "Ratio", oenv
);
380 sprintf(legend
[0], "Total/mean");
381 sprintf(legend
[1], "Total");
382 sprintf(legend
[2], "Mean");
383 sprintf(legend
[3], "# atoms");
384 sprintf(legend
[4], "Mean/# atoms");
385 xvgr_legend(fp
, 5, legend
, oenv
);
386 for (i
= 0; (i
< nres
); i
++)
394 ratio
= tot_n
[i
]/mean_n
[i
];
396 fprintf(fp
, "%3d %8.3f %3d %8.3f %3d %8.3f\n",
397 i
+1, ratio
, tot_n
[i
], mean_n
[i
], natm
[i
], mean_n
[i
]/natm
[i
]);