2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2007,2008,2009,2010,2011,2012,2013,2014,2015,2017,2019, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
40 #include "gromacs/commandline/pargs.h"
41 #include "gromacs/fileio/confio.h"
42 #include "gromacs/fileio/trxio.h"
43 #include "gromacs/gmxana/gmx_ana.h"
44 #include "gromacs/math/vec.h"
45 #include "gromacs/pbcutil/pbc.h"
46 #include "gromacs/pbcutil/rmpbc.h"
47 #include "gromacs/topology/index.h"
48 #include "gromacs/topology/topology.h"
49 #include "gromacs/trajectory/trajectoryframe.h"
50 #include "gromacs/utility/arraysize.h"
51 #include "gromacs/utility/cstringutil.h"
52 #include "gromacs/utility/futil.h"
53 #include "gromacs/utility/smalloc.h"
55 static const double bohr =
56 0.529177249; /* conversion factor to compensate for VMD plugin conversion... */
58 int gmx_spatial(int argc, char* argv[])
60 const char* desc[] = {
61 "[THISMODULE] calculates the spatial distribution function and",
62 "outputs it in a form that can be read by VMD as Gaussian98 cube format.",
63 "For a system of 32,000 atoms and a 50 ns trajectory, the SDF can be generated",
64 "in about 30 minutes, with most of the time dedicated to the two runs through",
65 "[TT]trjconv[tt] that are required to center everything properly.",
66 "This also takes a whole bunch of space (3 copies of the trajectory file).",
67 "Still, the pictures are pretty and very informative when the fitted selection is ",
70 "3-4 atoms in a widely mobile group (like a free amino acid in solution) works",
71 "well, or select the protein backbone in a stable folded structure to get the SDF",
72 "of solvent and look at the time-averaged solvation shell.",
73 "It is also possible using this program to generate the SDF based on some arbitrary",
74 "Cartesian coordinate. To do that, simply omit the preliminary [gmx-trjconv] steps.",
78 "1. Use [gmx-make_ndx] to create a group containing the atoms around which you want the ",
80 "2. [TT]gmx trjconv -s a.tpr -f a.tng -o b.tng -boxcenter tric -ur compact -pbc none[tt]",
81 "3. [TT]gmx trjconv -s a.tpr -f b.tng -o c.tng -fit rot+trans[tt]",
82 "4. run [THISMODULE] on the [TT]c.tng[tt] output of step #3.",
83 "5. Load [TT]grid.cube[tt] into VMD and view as an isosurface.",
85 "[BB]Note[bb] that systems such as micelles will require [TT]gmx trjconv -pbc cluster[tt] ",
86 "between steps 1 and 2.",
91 "The SDF will be generated for a cube that contains all bins that have some non-zero ",
93 "However, the preparatory [TT]-fit rot+trans[tt] option to [gmx-trjconv] implies that ",
94 "your system will be rotating",
95 "and translating in space (in order that the selected group does not). Therefore the ",
97 "returned will only be valid for some region around your central group/coordinate that ",
99 "with system volume throughout the entire translated/rotated system over the course of ",
101 "It is up to the user to ensure that this is the case.",
106 "To reduce the amount of space and time required, you can output only the coords",
107 "that are going to be used in the first and subsequent run through [gmx-trjconv].",
108 "However, be sure to set the [TT]-nab[tt] option to a sufficiently high value since",
109 "memory is allocated for cube bins based on the initial coordinates and the [TT]-nab[tt]",
112 const char* bugs[] = {
113 "When the allocated memory is not large enough, a segmentation fault may occur. ",
114 "This is usually detected ",
115 "and the program is halted prior to the fault while displaying a warning message ",
116 "suggesting the use of the [TT]-nab[tt] (Number of Additional Bins) ",
117 "option. However, the program does not detect all such events. If you encounter a ",
118 "segmentation fault, run it again ",
119 "with an increased [TT]-nab[tt] value."
122 static gmx_bool bPBC = FALSE;
123 static int iIGNOREOUTER = -1; /*Positive values may help if the surface is spikey */
124 static gmx_bool bCUTDOWN = TRUE;
125 static real rBINWIDTH = 0.05; /* nm */
126 static gmx_bool bCALCDIV = TRUE;
129 t_pargs pa[] = { { "-pbc",
133 "Use periodic boundary conditions for computing distances" },
138 "Calculate and apply the divisor for bin occupancies based on atoms/minimal "
139 "cube size. Set as TRUE for visualization and as FALSE ([TT]-nodiv[tt]) to "
140 "get accurate counts per frame" },
145 "Do not display this number of outer cubes (positive values may reduce "
146 "boundary speckles; -1 ensures outer surface is visible)" },
147 /* { "-cut", bCUTDOWN, etBOOL, {&bCUTDOWN},*/
148 /* "Display a total cube that is of minimal size" }, */
149 { "-bin", FALSE, etREAL, { &rBINWIDTH }, "Width of the bins (nm)" },
154 "Number of additional bins to ensure proper memory allocation" } };
164 int flags = TRX_READ_X;
168 char * grpnm, *grpnmp;
169 int * index, *indexp;
173 int*** bin = nullptr;
176 int x, y, z, minx, miny, minz, maxx, maxy, maxz;
178 int tot, maxval, minval;
180 gmx_output_env_t* oenv;
181 gmx_rmpbc_t gpbc = nullptr;
183 t_filenm fnm[] = { { efTPS, nullptr, nullptr, ffREAD }, /* this is for the topology */
184 { efTRX, "-f", nullptr, ffREAD }, /* and this for the trajectory */
185 { efNDX, nullptr, nullptr, ffOPTRD } };
187 #define NFILE asize(fnm)
189 /* This is the routine responsible for adding default options,
190 * calling the X/motif interface, etc. */
191 if (!parse_common_args(&argc, argv, PCA_CAN_TIME | PCA_CAN_VIEW, NFILE, fnm, asize(pa), pa,
192 asize(desc), desc, asize(bugs), bugs, &oenv))
197 read_tps_conf(ftp2fn(efTPS, NFILE, fnm), &top, &ePBC, &xtop, nullptr, box, TRUE);
200 atoms = &(top.atoms);
201 printf("Select group to generate SDF:\n");
202 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &nidx, &index, &grpnm);
203 printf("Select group to output coords (e.g. solute):\n");
204 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &nidxp, &indexp, &grpnmp);
206 /* The first time we read data is a little special */
207 read_first_frame(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &fr, flags);
210 /* Memory Allocation */
211 MINBIN[XX] = MAXBIN[XX] = fr.x[0][XX];
212 MINBIN[YY] = MAXBIN[YY] = fr.x[0][YY];
213 MINBIN[ZZ] = MAXBIN[ZZ] = fr.x[0][ZZ];
214 for (i = 1; i < top.atoms.nr; ++i)
216 if (fr.x[i][XX] < MINBIN[XX])
218 MINBIN[XX] = fr.x[i][XX];
220 if (fr.x[i][XX] > MAXBIN[XX])
222 MAXBIN[XX] = fr.x[i][XX];
224 if (fr.x[i][YY] < MINBIN[YY])
226 MINBIN[YY] = fr.x[i][YY];
228 if (fr.x[i][YY] > MAXBIN[YY])
230 MAXBIN[YY] = fr.x[i][YY];
232 if (fr.x[i][ZZ] < MINBIN[ZZ])
234 MINBIN[ZZ] = fr.x[i][ZZ];
236 if (fr.x[i][ZZ] > MAXBIN[ZZ])
238 MAXBIN[ZZ] = fr.x[i][ZZ];
241 for (i = ZZ; i >= XX; --i)
243 MAXBIN[i] = (std::ceil((MAXBIN[i] - MINBIN[i]) / rBINWIDTH) + iNAB) * rBINWIDTH + MINBIN[i];
244 MINBIN[i] -= iNAB * rBINWIDTH;
245 nbin[i] = static_cast<int>(std::ceil((MAXBIN[i] - MINBIN[i]) / rBINWIDTH));
248 for (i = 0; i < nbin[XX]; ++i)
250 snew(bin[i], nbin[YY]);
251 for (j = 0; j < nbin[YY]; ++j)
253 snew(bin[i][j], nbin[ZZ]);
256 copy_mat(box, box_pbc);
258 minx = miny = minz = 999;
259 maxx = maxy = maxz = 0;
263 gpbc = gmx_rmpbc_init(&top.idef, ePBC, natoms);
265 /* This is the main loop over frames */
268 /* Must init pbc every step because of pressure coupling */
270 copy_mat(box, box_pbc);
273 gmx_rmpbc_trxfr(gpbc, &fr);
274 set_pbc(&pbc, ePBC, box_pbc);
277 for (i = 0; i < nidx; i++)
279 if (fr.x[index[i]][XX] < MINBIN[XX] || fr.x[index[i]][XX] > MAXBIN[XX]
280 || fr.x[index[i]][YY] < MINBIN[YY] || fr.x[index[i]][YY] > MAXBIN[YY]
281 || fr.x[index[i]][ZZ] < MINBIN[ZZ] || fr.x[index[i]][ZZ] > MAXBIN[ZZ])
283 printf("There was an item outside of the allocated memory. Increase the value "
284 "given with the -nab option.\n");
285 printf("Memory was allocated for [%f,%f,%f]\tto\t[%f,%f,%f]\n", MINBIN[XX],
286 MINBIN[YY], MINBIN[ZZ], MAXBIN[XX], MAXBIN[YY], MAXBIN[ZZ]);
287 printf("Memory was required for [%f,%f,%f]\n", fr.x[index[i]][XX],
288 fr.x[index[i]][YY], fr.x[index[i]][ZZ]);
291 x = static_cast<int>(std::ceil((fr.x[index[i]][XX] - MINBIN[XX]) / rBINWIDTH));
292 y = static_cast<int>(std::ceil((fr.x[index[i]][YY] - MINBIN[YY]) / rBINWIDTH));
293 z = static_cast<int>(std::ceil((fr.x[index[i]][ZZ] - MINBIN[ZZ]) / rBINWIDTH));
321 /* printf("%f\t%f\t%f\n",box[XX][XX],box[YY][YY],box[ZZ][ZZ]); */
323 } while (read_next_frame(oenv, status, &fr));
327 gmx_rmpbc_done(gpbc);
332 minx = miny = minz = 0;
339 flp = gmx_ffopen("grid.cube", "w");
340 fprintf(flp, "Spatial Distribution Function\n");
341 fprintf(flp, "test\n");
342 fprintf(flp, "%5d%12.6f%12.6f%12.6f\n", nidxp,
343 (MINBIN[XX] + (minx + iIGNOREOUTER) * rBINWIDTH) * 10. / bohr,
344 (MINBIN[YY] + (miny + iIGNOREOUTER) * rBINWIDTH) * 10. / bohr,
345 (MINBIN[ZZ] + (minz + iIGNOREOUTER) * rBINWIDTH) * 10. / bohr);
346 fprintf(flp, "%5d%12.6f%12.6f%12.6f\n", maxx - minx + 1 - (2 * iIGNOREOUTER),
347 rBINWIDTH * 10. / bohr, 0., 0.);
348 fprintf(flp, "%5d%12.6f%12.6f%12.6f\n", maxy - miny + 1 - (2 * iIGNOREOUTER), 0.,
349 rBINWIDTH * 10. / bohr, 0.);
350 fprintf(flp, "%5d%12.6f%12.6f%12.6f\n", maxz - minz + 1 - (2 * iIGNOREOUTER), 0., 0.,
351 rBINWIDTH * 10. / bohr);
352 for (i = 0; i < nidxp; i++)
355 if (*(top.atoms.atomname[indexp[i]][0]) == 'C')
359 if (*(top.atoms.atomname[indexp[i]][0]) == 'N')
363 if (*(top.atoms.atomname[indexp[i]][0]) == 'O')
367 if (*(top.atoms.atomname[indexp[i]][0]) == 'H')
371 if (*(top.atoms.atomname[indexp[i]][0]) == 'S')
375 fprintf(flp, "%5d%12.6f%12.6f%12.6f%12.6f\n", v, 0., fr.x[indexp[i]][XX] * 10.0 / bohr,
376 fr.x[indexp[i]][YY] * 10.0 / bohr, fr.x[indexp[i]][ZZ] * 10.0 / bohr);
380 for (k = 0; k < nbin[XX]; k++)
382 if (!(k < minx || k > maxx))
386 for (j = 0; j < nbin[YY]; j++)
388 if (!(j < miny || j > maxy))
392 for (i = 0; i < nbin[ZZ]; i++)
394 if (!(i < minz || i > maxz))
398 if (bin[k][j][i] != 0)
400 printf("A bin was not empty when it should have been empty. Programming "
402 printf("bin[%d][%d][%d] was = %d\n", k, j, i, bin[k][j][i]);
411 for (k = 0; k < nbin[XX]; k++)
413 if (k < minx + iIGNOREOUTER || k > maxx - iIGNOREOUTER)
417 for (j = 0; j < nbin[YY]; j++)
419 if (j < miny + iIGNOREOUTER || j > maxy - iIGNOREOUTER)
423 for (i = 0; i < nbin[ZZ]; i++)
425 if (i < minz + iIGNOREOUTER || i > maxz - iIGNOREOUTER)
430 if (bin[k][j][i] > maxval)
432 maxval = bin[k][j][i];
434 if (bin[k][j][i] < minval)
436 minval = bin[k][j][i];
442 numcu = (maxx - minx + 1 - (2 * iIGNOREOUTER)) * (maxy - miny + 1 - (2 * iIGNOREOUTER))
443 * (maxz - minz + 1 - (2 * iIGNOREOUTER));
446 norm = static_cast<double>(numcu * numfr) / tot;
453 for (k = 0; k < nbin[XX]; k++)
455 if (k < minx + iIGNOREOUTER || k > maxx - iIGNOREOUTER)
459 for (j = 0; j < nbin[YY]; j++)
461 if (j < miny + iIGNOREOUTER || j > maxy - iIGNOREOUTER)
465 for (i = 0; i < nbin[ZZ]; i++)
467 if (i < minz + iIGNOREOUTER || i > maxz - iIGNOREOUTER)
471 fprintf(flp, "%12.6f ", static_cast<double>(norm * bin[k][j][i]) / numfr);
481 printf("Counts per frame in all %d cubes divided by %le\n", numcu, 1.0 / norm);
482 printf("Normalized data: average %le, min %le, max %le\n", 1.0, minval * norm / numfr,
483 maxval * norm / numfr);
487 printf("grid.cube contains counts per frame in all %d cubes\n", numcu);
488 printf("Raw data: average %le, min %le, max %le\n", 1.0 / norm,
489 static_cast<double>(minval) / numfr, static_cast<double>(maxval) / numfr);