Update instructions in containers.rst
[gromacs.git] / src / gromacs / gmxana / gmx_spatial.cpp
blob0ec43ffafda9f4a68fdc4bbe948b17d00c0755d8
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2007,2008,2009,2010,2011 by the GROMACS development team.
5 * Copyright (c) 2012,2013,2014,2015,2017 by the GROMACS development team.
6 * Copyright (c) 2019,2020, 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.
37 #include "gmxpre.h"
39 #include <cmath>
40 #include <cstdlib>
42 #include "gromacs/commandline/pargs.h"
43 #include "gromacs/fileio/confio.h"
44 #include "gromacs/fileio/trxio.h"
45 #include "gromacs/gmxana/gmx_ana.h"
46 #include "gromacs/math/vec.h"
47 #include "gromacs/pbcutil/pbc.h"
48 #include "gromacs/pbcutil/rmpbc.h"
49 #include "gromacs/topology/index.h"
50 #include "gromacs/topology/topology.h"
51 #include "gromacs/trajectory/trajectoryframe.h"
52 #include "gromacs/utility/arraysize.h"
53 #include "gromacs/utility/cstringutil.h"
54 #include "gromacs/utility/futil.h"
55 #include "gromacs/utility/smalloc.h"
57 static const double bohr =
58 0.529177249; /* conversion factor to compensate for VMD plugin conversion... */
60 int gmx_spatial(int argc, char* argv[])
62 const char* desc[] = {
63 "[THISMODULE] calculates the spatial distribution function and",
64 "outputs it in a form that can be read by VMD as Gaussian98 cube format.",
65 "For a system of 32,000 atoms and a 50 ns trajectory, the SDF can be generated",
66 "in about 30 minutes, with most of the time dedicated to the two runs through",
67 "[TT]trjconv[tt] that are required to center everything properly.",
68 "This also takes a whole bunch of space (3 copies of the trajectory file).",
69 "Still, the pictures are pretty and very informative when the fitted selection is ",
70 "properly ",
71 "made.",
72 "3-4 atoms in a widely mobile group (like a free amino acid in solution) works",
73 "well, or select the protein backbone in a stable folded structure to get the SDF",
74 "of solvent and look at the time-averaged solvation shell.",
75 "It is also possible using this program to generate the SDF based on some arbitrary",
76 "Cartesian coordinate. To do that, simply omit the preliminary [gmx-trjconv] steps.",
77 "",
78 "Usage:",
79 "",
80 "1. Use [gmx-make_ndx] to create a group containing the atoms around which you want the ",
81 "SDF",
82 "2. [TT]gmx trjconv -s a.tpr -f a.tng -o b.tng -boxcenter tric -ur compact -pbc none[tt]",
83 "3. [TT]gmx trjconv -s a.tpr -f b.tng -o c.tng -fit rot+trans[tt]",
84 "4. run [THISMODULE] on the [TT]c.tng[tt] output of step #3.",
85 "5. Load [TT]grid.cube[tt] into VMD and view as an isosurface.",
86 "",
87 "[BB]Note[bb] that systems such as micelles will require [TT]gmx trjconv -pbc cluster[tt] ",
88 "between steps 1 and 2.",
89 "",
90 "Warnings",
91 "^^^^^^^^",
92 "",
93 "The SDF will be generated for a cube that contains all bins that have some non-zero ",
94 "occupancy.",
95 "However, the preparatory [TT]-fit rot+trans[tt] option to [gmx-trjconv] implies that ",
96 "your system will be rotating",
97 "and translating in space (in order that the selected group does not). Therefore the ",
98 "values that are",
99 "returned will only be valid for some region around your central group/coordinate that ",
100 "has full overlap",
101 "with system volume throughout the entire translated/rotated system over the course of ",
102 "the trajectory.",
103 "It is up to the user to ensure that this is the case.",
105 "Risky options",
106 "^^^^^^^^^^^^^",
108 "To reduce the amount of space and time required, you can output only the coords",
109 "that are going to be used in the first and subsequent run through [gmx-trjconv].",
110 "However, be sure to set the [TT]-nab[tt] option to a sufficiently high value since",
111 "memory is allocated for cube bins based on the initial coordinates and the [TT]-nab[tt]",
112 "option value."
114 const char* bugs[] = {
115 "When the allocated memory is not large enough, a segmentation fault may occur. ",
116 "This is usually detected ",
117 "and the program is halted prior to the fault while displaying a warning message ",
118 "suggesting the use of the [TT]-nab[tt] (Number of Additional Bins) ",
119 "option. However, the program does not detect all such events. If you encounter a ",
120 "segmentation fault, run it again ",
121 "with an increased [TT]-nab[tt] value."
124 static gmx_bool bPBC = FALSE;
125 static int iIGNOREOUTER = -1; /*Positive values may help if the surface is spikey */
126 static gmx_bool bCUTDOWN = TRUE;
127 static real rBINWIDTH = 0.05; /* nm */
128 static gmx_bool bCALCDIV = TRUE;
129 static int iNAB = 4;
131 t_pargs pa[] = { { "-pbc",
132 FALSE,
133 etBOOL,
134 { &bPBC },
135 "Use periodic boundary conditions for computing distances" },
136 { "-div",
137 FALSE,
138 etBOOL,
139 { &bCALCDIV },
140 "Calculate and apply the divisor for bin occupancies based on atoms/minimal "
141 "cube size. Set as TRUE for visualization and as FALSE ([TT]-nodiv[tt]) to "
142 "get accurate counts per frame" },
143 { "-ign",
144 FALSE,
145 etINT,
146 { &iIGNOREOUTER },
147 "Do not display this number of outer cubes (positive values may reduce "
148 "boundary speckles; -1 ensures outer surface is visible)" },
149 /* { "-cut", bCUTDOWN, etBOOL, {&bCUTDOWN},*/
150 /* "Display a total cube that is of minimal size" }, */
151 { "-bin", FALSE, etREAL, { &rBINWIDTH }, "Width of the bins (nm)" },
152 { "-nab",
153 FALSE,
154 etINT,
155 { &iNAB },
156 "Number of additional bins to ensure proper memory allocation" } };
158 double MINBIN[3];
159 double MAXBIN[3];
160 t_topology top;
161 PbcType pbcType;
162 t_trxframe fr;
163 rvec* xtop;
164 matrix box, box_pbc;
165 t_trxstatus* status;
166 int flags = TRX_READ_X;
167 t_pbc pbc;
168 t_atoms* atoms;
169 int natoms;
170 char * grpnm, *grpnmp;
171 int * index, *indexp;
172 int i, nidx, nidxp;
173 int v;
174 int j, k;
175 int*** bin = nullptr;
176 int nbin[3];
177 FILE* flp;
178 int x, y, z, minx, miny, minz, maxx, maxy, maxz;
179 int numfr, numcu;
180 int tot, maxval, minval;
181 double norm;
182 gmx_output_env_t* oenv;
183 gmx_rmpbc_t gpbc = nullptr;
185 t_filenm fnm[] = { { efTPS, nullptr, nullptr, ffREAD }, /* this is for the topology */
186 { efTRX, "-f", nullptr, ffREAD }, /* and this for the trajectory */
187 { efNDX, nullptr, nullptr, ffOPTRD } };
189 #define NFILE asize(fnm)
191 /* This is the routine responsible for adding default options,
192 * calling the X/motif interface, etc. */
193 if (!parse_common_args(&argc, argv, PCA_CAN_TIME | PCA_CAN_VIEW, NFILE, fnm, asize(pa), pa,
194 asize(desc), desc, asize(bugs), bugs, &oenv))
196 return 0;
199 read_tps_conf(ftp2fn(efTPS, NFILE, fnm), &top, &pbcType, &xtop, nullptr, box, TRUE);
200 sfree(xtop);
202 atoms = &(top.atoms);
203 printf("Select group to generate SDF:\n");
204 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &nidx, &index, &grpnm);
205 printf("Select group to output coords (e.g. solute):\n");
206 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &nidxp, &indexp, &grpnmp);
208 /* The first time we read data is a little special */
209 read_first_frame(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &fr, flags);
210 natoms = fr.natoms;
212 /* Memory Allocation */
213 MINBIN[XX] = MAXBIN[XX] = fr.x[0][XX];
214 MINBIN[YY] = MAXBIN[YY] = fr.x[0][YY];
215 MINBIN[ZZ] = MAXBIN[ZZ] = fr.x[0][ZZ];
216 for (i = 1; i < top.atoms.nr; ++i)
218 if (fr.x[i][XX] < MINBIN[XX])
220 MINBIN[XX] = fr.x[i][XX];
222 if (fr.x[i][XX] > MAXBIN[XX])
224 MAXBIN[XX] = fr.x[i][XX];
226 if (fr.x[i][YY] < MINBIN[YY])
228 MINBIN[YY] = fr.x[i][YY];
230 if (fr.x[i][YY] > MAXBIN[YY])
232 MAXBIN[YY] = fr.x[i][YY];
234 if (fr.x[i][ZZ] < MINBIN[ZZ])
236 MINBIN[ZZ] = fr.x[i][ZZ];
238 if (fr.x[i][ZZ] > MAXBIN[ZZ])
240 MAXBIN[ZZ] = fr.x[i][ZZ];
243 for (i = ZZ; i >= XX; --i)
245 MAXBIN[i] = (std::ceil((MAXBIN[i] - MINBIN[i]) / rBINWIDTH) + iNAB) * rBINWIDTH + MINBIN[i];
246 MINBIN[i] -= iNAB * rBINWIDTH;
247 nbin[i] = static_cast<int>(std::ceil((MAXBIN[i] - MINBIN[i]) / rBINWIDTH));
249 snew(bin, nbin[XX]);
250 for (i = 0; i < nbin[XX]; ++i)
252 snew(bin[i], nbin[YY]);
253 for (j = 0; j < nbin[YY]; ++j)
255 snew(bin[i][j], nbin[ZZ]);
258 copy_mat(box, box_pbc);
259 numfr = 0;
260 minx = miny = minz = 999;
261 maxx = maxy = maxz = 0;
263 if (bPBC)
265 gpbc = gmx_rmpbc_init(&top.idef, pbcType, natoms);
267 /* This is the main loop over frames */
270 /* Must init pbc every step because of pressure coupling */
272 copy_mat(box, box_pbc);
273 if (bPBC)
275 gmx_rmpbc_trxfr(gpbc, &fr);
276 set_pbc(&pbc, pbcType, box_pbc);
279 for (i = 0; i < nidx; i++)
281 if (fr.x[index[i]][XX] < MINBIN[XX] || fr.x[index[i]][XX] > MAXBIN[XX]
282 || fr.x[index[i]][YY] < MINBIN[YY] || fr.x[index[i]][YY] > MAXBIN[YY]
283 || fr.x[index[i]][ZZ] < MINBIN[ZZ] || fr.x[index[i]][ZZ] > MAXBIN[ZZ])
285 printf("There was an item outside of the allocated memory. Increase the value "
286 "given with the -nab option.\n");
287 printf("Memory was allocated for [%f,%f,%f]\tto\t[%f,%f,%f]\n", MINBIN[XX],
288 MINBIN[YY], MINBIN[ZZ], MAXBIN[XX], MAXBIN[YY], MAXBIN[ZZ]);
289 printf("Memory was required for [%f,%f,%f]\n", fr.x[index[i]][XX],
290 fr.x[index[i]][YY], fr.x[index[i]][ZZ]);
291 exit(1);
293 x = static_cast<int>(std::ceil((fr.x[index[i]][XX] - MINBIN[XX]) / rBINWIDTH));
294 y = static_cast<int>(std::ceil((fr.x[index[i]][YY] - MINBIN[YY]) / rBINWIDTH));
295 z = static_cast<int>(std::ceil((fr.x[index[i]][ZZ] - MINBIN[ZZ]) / rBINWIDTH));
296 ++bin[x][y][z];
297 if (x < minx)
299 minx = x;
301 if (x > maxx)
303 maxx = x;
305 if (y < miny)
307 miny = y;
309 if (y > maxy)
311 maxy = y;
313 if (z < minz)
315 minz = z;
317 if (z > maxz)
319 maxz = z;
322 numfr++;
323 /* printf("%f\t%f\t%f\n",box[XX][XX],box[YY][YY],box[ZZ][ZZ]); */
325 } while (read_next_frame(oenv, status, &fr));
327 if (bPBC)
329 gmx_rmpbc_done(gpbc);
332 if (!bCUTDOWN)
334 minx = miny = minz = 0;
335 maxx = nbin[XX];
336 maxy = nbin[YY];
337 maxz = nbin[ZZ];
340 /* OUTPUT */
341 flp = gmx_ffopen("grid.cube", "w");
342 fprintf(flp, "Spatial Distribution Function\n");
343 fprintf(flp, "test\n");
344 fprintf(flp, "%5d%12.6f%12.6f%12.6f\n", nidxp,
345 (MINBIN[XX] + (minx + iIGNOREOUTER) * rBINWIDTH) * 10. / bohr,
346 (MINBIN[YY] + (miny + iIGNOREOUTER) * rBINWIDTH) * 10. / bohr,
347 (MINBIN[ZZ] + (minz + iIGNOREOUTER) * rBINWIDTH) * 10. / bohr);
348 fprintf(flp, "%5d%12.6f%12.6f%12.6f\n", maxx - minx + 1 - (2 * iIGNOREOUTER),
349 rBINWIDTH * 10. / bohr, 0., 0.);
350 fprintf(flp, "%5d%12.6f%12.6f%12.6f\n", maxy - miny + 1 - (2 * iIGNOREOUTER), 0.,
351 rBINWIDTH * 10. / bohr, 0.);
352 fprintf(flp, "%5d%12.6f%12.6f%12.6f\n", maxz - minz + 1 - (2 * iIGNOREOUTER), 0., 0.,
353 rBINWIDTH * 10. / bohr);
354 for (i = 0; i < nidxp; i++)
356 v = 2;
357 if (*(top.atoms.atomname[indexp[i]][0]) == 'C')
359 v = 6;
361 if (*(top.atoms.atomname[indexp[i]][0]) == 'N')
363 v = 7;
365 if (*(top.atoms.atomname[indexp[i]][0]) == 'O')
367 v = 8;
369 if (*(top.atoms.atomname[indexp[i]][0]) == 'H')
371 v = 1;
373 if (*(top.atoms.atomname[indexp[i]][0]) == 'S')
375 v = 16;
377 fprintf(flp, "%5d%12.6f%12.6f%12.6f%12.6f\n", v, 0., fr.x[indexp[i]][XX] * 10.0 / bohr,
378 fr.x[indexp[i]][YY] * 10.0 / bohr, fr.x[indexp[i]][ZZ] * 10.0 / bohr);
381 tot = 0;
382 for (k = 0; k < nbin[XX]; k++)
384 if (!(k < minx || k > maxx))
386 continue;
388 for (j = 0; j < nbin[YY]; j++)
390 if (!(j < miny || j > maxy))
392 continue;
394 for (i = 0; i < nbin[ZZ]; i++)
396 if (!(i < minz || i > maxz))
398 continue;
400 if (bin[k][j][i] != 0)
402 printf("A bin was not empty when it should have been empty. Programming "
403 "error.\n");
404 printf("bin[%d][%d][%d] was = %d\n", k, j, i, bin[k][j][i]);
405 exit(1);
411 minval = 999;
412 maxval = 0;
413 for (k = 0; k < nbin[XX]; k++)
415 if (k < minx + iIGNOREOUTER || k > maxx - iIGNOREOUTER)
417 continue;
419 for (j = 0; j < nbin[YY]; j++)
421 if (j < miny + iIGNOREOUTER || j > maxy - iIGNOREOUTER)
423 continue;
425 for (i = 0; i < nbin[ZZ]; i++)
427 if (i < minz + iIGNOREOUTER || i > maxz - iIGNOREOUTER)
429 continue;
431 tot += bin[k][j][i];
432 if (bin[k][j][i] > maxval)
434 maxval = bin[k][j][i];
436 if (bin[k][j][i] < minval)
438 minval = bin[k][j][i];
444 numcu = (maxx - minx + 1 - (2 * iIGNOREOUTER)) * (maxy - miny + 1 - (2 * iIGNOREOUTER))
445 * (maxz - minz + 1 - (2 * iIGNOREOUTER));
446 if (bCALCDIV)
448 norm = static_cast<double>(numcu * numfr) / tot;
450 else
452 norm = 1.0;
455 for (k = 0; k < nbin[XX]; k++)
457 if (k < minx + iIGNOREOUTER || k > maxx - iIGNOREOUTER)
459 continue;
461 for (j = 0; j < nbin[YY]; j++)
463 if (j < miny + iIGNOREOUTER || j > maxy - iIGNOREOUTER)
465 continue;
467 for (i = 0; i < nbin[ZZ]; i++)
469 if (i < minz + iIGNOREOUTER || i > maxz - iIGNOREOUTER)
471 continue;
473 fprintf(flp, "%12.6f ", static_cast<double>(norm * bin[k][j][i]) / numfr);
475 fprintf(flp, "\n");
477 fprintf(flp, "\n");
479 gmx_ffclose(flp);
481 if (bCALCDIV)
483 printf("Counts per frame in all %d cubes divided by %le\n", numcu, 1.0 / norm);
484 printf("Normalized data: average %le, min %le, max %le\n", 1.0, minval * norm / numfr,
485 maxval * norm / numfr);
487 else
489 printf("grid.cube contains counts per frame in all %d cubes\n", numcu);
490 printf("Raw data: average %le, min %le, max %le\n", 1.0 / norm,
491 static_cast<double>(minval) / numfr, static_cast<double>(maxval) / numfr);
494 return 0;