Remove unused code detected by PGI compiler
[gromacs.git] / src / gromacs / gmxana / gmx_spatial.cpp
blobba4899e83ead484d4ef948678d835382ff84362e
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2007,2008,2009,2010,2011,2012,2013,2014,2015, 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.
35 #include "gmxpre.h"
37 #include <cmath>
38 #include <cstdlib>
40 #include "gromacs/commandline/pargs.h"
41 #include "gromacs/fileio/confio.h"
42 #include "gromacs/fileio/trx.h"
43 #include "gromacs/fileio/trxio.h"
44 #include "gromacs/gmxana/gmx_ana.h"
45 #include "gromacs/math/vec.h"
46 #include "gromacs/pbcutil/pbc.h"
47 #include "gromacs/pbcutil/rmpbc.h"
48 #include "gromacs/topology/index.h"
49 #include "gromacs/topology/topology.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 = 0.529177249; /* conversion factor to compensate for VMD plugin conversion... */
57 int gmx_spatial(int argc, char *argv[])
59 const char *desc[] = {
60 "[THISMODULE] calculates the spatial distribution function and",
61 "outputs it in a form that can be read by VMD as Gaussian98 cube format.",
62 "For a system of 32,000 atoms and a 50 ns trajectory, the SDF can be generated",
63 "in about 30 minutes, with most of the time dedicated to the two runs through",
64 "[TT]trjconv[tt] that are required to center everything properly.",
65 "This also takes a whole bunch of space (3 copies of the trajectory file).",
66 "Still, the pictures are pretty and very informative when the fitted selection is properly made.",
67 "3-4 atoms in a widely mobile group (like a free amino acid in solution) works",
68 "well, or select the protein backbone in a stable folded structure to get the SDF",
69 "of solvent and look at the time-averaged solvation shell.",
70 "It is also possible using this program to generate the SDF based on some arbitrary",
71 "Cartesian coordinate. To do that, simply omit the preliminary [gmx-trjconv] steps.",
72 "",
73 "Usage:",
74 "",
75 "1. Use [gmx-make_ndx] to create a group containing the atoms around which you want the SDF",
76 "2. [TT]gmx trjconv -s a.tpr -f a.tng -o b.tng -boxcenter tric -ur compact -pbc none[tt]",
77 "3. [TT]gmx trjconv -s a.tpr -f b.tng -o c.tng -fit rot+trans[tt]",
78 "4. run [THISMODULE] on the [TT]c.tng[tt] output of step #3.",
79 "5. Load [TT]grid.cube[tt] into VMD and view as an isosurface.",
80 "",
81 "[BB]Note[bb] that systems such as micelles will require [TT]gmx trjconv -pbc cluster[tt] between steps 1 and 2.",
82 "",
83 "Warnings",
84 "^^^^^^^^",
85 "",
86 "The SDF will be generated for a cube that contains all bins that have some non-zero occupancy.",
87 "However, the preparatory [TT]-fit rot+trans[tt] option to [gmx-trjconv] implies that your system will be rotating",
88 "and translating in space (in order that the selected group does not). Therefore the values that are",
89 "returned will only be valid for some region around your central group/coordinate that has full overlap",
90 "with system volume throughout the entire translated/rotated system over the course of the trajectory.",
91 "It is up to the user to ensure that this is the case.",
92 "",
93 "Risky options",
94 "^^^^^^^^^^^^^",
95 "",
96 "To reduce the amount of space and time required, you can output only the coords",
97 "that are going to be used in the first and subsequent run through [gmx-trjconv].",
98 "However, be sure to set the [TT]-nab[tt] option to a sufficiently high value since",
99 "memory is allocated for cube bins based on the initial coordinates and the [TT]-nab[tt]",
100 "option value."
102 const char *bugs[] = {
103 "When the allocated memory is not large enough, a segmentation fault may occur. This is usually detected "
104 "and the program is halted prior to the fault while displaying a warning message suggesting the use of the [TT]-nab[tt] (Number of Additional Bins) "
105 "option. However, the program does not detect all such events. If you encounter a segmentation fault, run it again "
106 "with an increased [TT]-nab[tt] value."
109 static gmx_bool bPBC = FALSE;
110 static int iIGNOREOUTER = -1; /*Positive values may help if the surface is spikey */
111 static gmx_bool bCUTDOWN = TRUE;
112 static real rBINWIDTH = 0.05; /* nm */
113 static gmx_bool bCALCDIV = TRUE;
114 static int iNAB = 4;
116 t_pargs pa[] = {
117 { "-pbc", FALSE, etBOOL, {&bPBC},
118 "Use periodic boundary conditions for computing distances" },
119 { "-div", FALSE, etBOOL, {&bCALCDIV},
120 "Calculate and apply the divisor for bin occupancies based on atoms/minimal cube size. Set as TRUE for visualization and as FALSE ([TT]-nodiv[tt]) to get accurate counts per frame" },
121 { "-ign", FALSE, etINT, {&iIGNOREOUTER},
122 "Do not display this number of outer cubes (positive values may reduce boundary speckles; -1 ensures outer surface is visible)" },
123 /* { "-cut", bCUTDOWN, etBOOL, {&bCUTDOWN},*/
124 /* "Display a total cube that is of minimal size" }, */
125 { "-bin", FALSE, etREAL, {&rBINWIDTH},
126 "Width of the bins (nm)" },
127 { "-nab", FALSE, etINT, {&iNAB},
128 "Number of additional bins to ensure proper memory allocation" }
131 double MINBIN[3];
132 double MAXBIN[3];
133 t_topology top;
134 int ePBC;
135 t_trxframe fr;
136 rvec *xtop;
137 matrix box, box_pbc;
138 t_trxstatus *status;
139 int flags = TRX_READ_X;
140 t_pbc pbc;
141 t_atoms *atoms;
142 int natoms;
143 char *grpnm, *grpnmp;
144 int *index, *indexp;
145 int i, nidx, nidxp;
146 int v;
147 int j, k;
148 int ***bin = NULL;
149 int nbin[3];
150 FILE *flp;
151 int x, y, z, minx, miny, minz, maxx, maxy, maxz;
152 int numfr, numcu;
153 int tot, maxval, minval;
154 double norm;
155 gmx_output_env_t *oenv;
156 gmx_rmpbc_t gpbc = NULL;
158 t_filenm fnm[] = {
159 { efTPS, NULL, NULL, ffREAD }, /* this is for the topology */
160 { efTRX, "-f", NULL, ffREAD }, /* and this for the trajectory */
161 { efNDX, NULL, NULL, ffOPTRD }
164 #define NFILE asize(fnm)
166 /* This is the routine responsible for adding default options,
167 * calling the X/motif interface, etc. */
168 if (!parse_common_args(&argc, argv, PCA_CAN_TIME | PCA_CAN_VIEW,
169 NFILE, fnm, asize(pa), pa, asize(desc), desc, asize(bugs), bugs, &oenv))
171 return 0;
174 read_tps_conf(ftp2fn(efTPS, NFILE, fnm), &top, &ePBC, &xtop, NULL, box, TRUE);
175 sfree(xtop);
177 atoms = &(top.atoms);
178 printf("Select group to generate SDF:\n");
179 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &nidx, &index, &grpnm);
180 printf("Select group to output coords (e.g. solute):\n");
181 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &nidxp, &indexp, &grpnmp);
183 /* The first time we read data is a little special */
184 natoms = read_first_frame(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &fr, flags);
186 /* Memory Allocation */
187 MINBIN[XX] = MAXBIN[XX] = fr.x[0][XX];
188 MINBIN[YY] = MAXBIN[YY] = fr.x[0][YY];
189 MINBIN[ZZ] = MAXBIN[ZZ] = fr.x[0][ZZ];
190 for (i = 1; i < top.atoms.nr; ++i)
192 if (fr.x[i][XX] < MINBIN[XX])
194 MINBIN[XX] = fr.x[i][XX];
196 if (fr.x[i][XX] > MAXBIN[XX])
198 MAXBIN[XX] = fr.x[i][XX];
200 if (fr.x[i][YY] < MINBIN[YY])
202 MINBIN[YY] = fr.x[i][YY];
204 if (fr.x[i][YY] > MAXBIN[YY])
206 MAXBIN[YY] = fr.x[i][YY];
208 if (fr.x[i][ZZ] < MINBIN[ZZ])
210 MINBIN[ZZ] = fr.x[i][ZZ];
212 if (fr.x[i][ZZ] > MAXBIN[ZZ])
214 MAXBIN[ZZ] = fr.x[i][ZZ];
217 for (i = ZZ; i >= XX; --i)
219 MAXBIN[i] = (std::ceil((MAXBIN[i]-MINBIN[i])/rBINWIDTH)+iNAB)*rBINWIDTH+MINBIN[i];
220 MINBIN[i] -= iNAB*rBINWIDTH;
221 nbin[i] = static_cast<int>(std::ceil((MAXBIN[i]-MINBIN[i])/rBINWIDTH));
223 snew(bin, nbin[XX]);
224 for (i = 0; i < nbin[XX]; ++i)
226 snew(bin[i], nbin[YY]);
227 for (j = 0; j < nbin[YY]; ++j)
229 snew(bin[i][j], nbin[ZZ]);
232 copy_mat(box, box_pbc);
233 numfr = 0;
234 minx = miny = minz = 999;
235 maxx = maxy = maxz = 0;
237 if (bPBC)
239 gpbc = gmx_rmpbc_init(&top.idef, ePBC, natoms);
241 /* This is the main loop over frames */
244 /* Must init pbc every step because of pressure coupling */
246 copy_mat(box, box_pbc);
247 if (bPBC)
249 gmx_rmpbc_trxfr(gpbc, &fr);
250 set_pbc(&pbc, ePBC, box_pbc);
253 for (i = 0; i < nidx; i++)
255 if (fr.x[index[i]][XX] < MINBIN[XX] || fr.x[index[i]][XX] > MAXBIN[XX] ||
256 fr.x[index[i]][YY] < MINBIN[YY] || fr.x[index[i]][YY] > MAXBIN[YY] ||
257 fr.x[index[i]][ZZ] < MINBIN[ZZ] || fr.x[index[i]][ZZ] > MAXBIN[ZZ])
259 printf("There was an item outside of the allocated memory. Increase the value given with the -nab option.\n");
260 printf("Memory was allocated for [%f,%f,%f]\tto\t[%f,%f,%f]\n", MINBIN[XX], MINBIN[YY], MINBIN[ZZ], MAXBIN[XX], MAXBIN[YY], MAXBIN[ZZ]);
261 printf("Memory was required for [%f,%f,%f]\n", fr.x[index[i]][XX], fr.x[index[i]][YY], fr.x[index[i]][ZZ]);
262 exit(1);
264 x = static_cast<int>(std::ceil((fr.x[index[i]][XX]-MINBIN[XX])/rBINWIDTH));
265 y = static_cast<int>(std::ceil((fr.x[index[i]][YY]-MINBIN[YY])/rBINWIDTH));
266 z = static_cast<int>(std::ceil((fr.x[index[i]][ZZ]-MINBIN[ZZ])/rBINWIDTH));
267 ++bin[x][y][z];
268 if (x < minx)
270 minx = x;
272 if (x > maxx)
274 maxx = x;
276 if (y < miny)
278 miny = y;
280 if (y > maxy)
282 maxy = y;
284 if (z < minz)
286 minz = z;
288 if (z > maxz)
290 maxz = z;
293 numfr++;
294 /* printf("%f\t%f\t%f\n",box[XX][XX],box[YY][YY],box[ZZ][ZZ]); */
297 while (read_next_frame(oenv, status, &fr));
299 if (bPBC)
301 gmx_rmpbc_done(gpbc);
304 if (!bCUTDOWN)
306 minx = miny = minz = 0;
307 maxx = nbin[XX];
308 maxy = nbin[YY];
309 maxz = nbin[ZZ];
312 /* OUTPUT */
313 flp = gmx_ffopen("grid.cube", "w");
314 fprintf(flp, "Spatial Distribution Function\n");
315 fprintf(flp, "test\n");
316 fprintf(flp, "%5d%12.6f%12.6f%12.6f\n", nidxp, (MINBIN[XX]+(minx+iIGNOREOUTER)*rBINWIDTH)*10./bohr, (MINBIN[YY]+(miny+iIGNOREOUTER)*rBINWIDTH)*10./bohr, (MINBIN[ZZ]+(minz+iIGNOREOUTER)*rBINWIDTH)*10./bohr);
317 fprintf(flp, "%5d%12.6f%12.6f%12.6f\n", maxx-minx+1-(2*iIGNOREOUTER), rBINWIDTH*10./bohr, 0., 0.);
318 fprintf(flp, "%5d%12.6f%12.6f%12.6f\n", maxy-miny+1-(2*iIGNOREOUTER), 0., rBINWIDTH*10./bohr, 0.);
319 fprintf(flp, "%5d%12.6f%12.6f%12.6f\n", maxz-minz+1-(2*iIGNOREOUTER), 0., 0., rBINWIDTH*10./bohr);
320 for (i = 0; i < nidxp; i++)
322 v = 2;
323 if (*(top.atoms.atomname[indexp[i]][0]) == 'C')
325 v = 6;
327 if (*(top.atoms.atomname[indexp[i]][0]) == 'N')
329 v = 7;
331 if (*(top.atoms.atomname[indexp[i]][0]) == 'O')
333 v = 8;
335 if (*(top.atoms.atomname[indexp[i]][0]) == 'H')
337 v = 1;
339 if (*(top.atoms.atomname[indexp[i]][0]) == 'S')
341 v = 16;
343 fprintf(flp, "%5d%12.6f%12.6f%12.6f%12.6f\n", v, 0., fr.x[indexp[i]][XX]*10.0/bohr, fr.x[indexp[i]][YY]*10.0/bohr, fr.x[indexp[i]][ZZ]*10.0/bohr);
346 tot = 0;
347 for (k = 0; k < nbin[XX]; k++)
349 if (!(k < minx || k > maxx))
351 continue;
353 for (j = 0; j < nbin[YY]; j++)
355 if (!(j < miny || j > maxy))
357 continue;
359 for (i = 0; i < nbin[ZZ]; i++)
361 if (!(i < minz || i > maxz))
363 continue;
365 if (bin[k][j][i] != 0)
367 printf("A bin was not empty when it should have been empty. Programming error.\n");
368 printf("bin[%d][%d][%d] was = %d\n", k, j, i, bin[k][j][i]);
369 exit(1);
375 minval = 999;
376 maxval = 0;
377 for (k = 0; k < nbin[XX]; k++)
379 if (k < minx+iIGNOREOUTER || k > maxx-iIGNOREOUTER)
381 continue;
383 for (j = 0; j < nbin[YY]; j++)
385 if (j < miny+iIGNOREOUTER || j > maxy-iIGNOREOUTER)
387 continue;
389 for (i = 0; i < nbin[ZZ]; i++)
391 if (i < minz+iIGNOREOUTER || i > maxz-iIGNOREOUTER)
393 continue;
395 tot += bin[k][j][i];
396 if (bin[k][j][i] > maxval)
398 maxval = bin[k][j][i];
400 if (bin[k][j][i] < minval)
402 minval = bin[k][j][i];
408 numcu = (maxx-minx+1-(2*iIGNOREOUTER))*(maxy-miny+1-(2*iIGNOREOUTER))*(maxz-minz+1-(2*iIGNOREOUTER));
409 if (bCALCDIV)
411 norm = static_cast<double>(numcu*numfr)/tot;
413 else
415 norm = 1.0;
418 for (k = 0; k < nbin[XX]; k++)
420 if (k < minx+iIGNOREOUTER || k > maxx-iIGNOREOUTER)
422 continue;
424 for (j = 0; j < nbin[YY]; j++)
426 if (j < miny+iIGNOREOUTER || j > maxy-iIGNOREOUTER)
428 continue;
430 for (i = 0; i < nbin[ZZ]; i++)
432 if (i < minz+iIGNOREOUTER || i > maxz-iIGNOREOUTER)
434 continue;
436 fprintf(flp, "%12.6f ", static_cast<double>(norm*bin[k][j][i])/numfr);
438 fprintf(flp, "\n");
440 fprintf(flp, "\n");
442 gmx_ffclose(flp);
444 if (bCALCDIV)
446 printf("Counts per frame in all %d cubes divided by %le\n", numcu, 1.0/norm);
447 printf("Normalized data: average %le, min %le, max %le\n", 1.0, minval*norm/numfr, maxval*norm/numfr);
449 else
451 printf("grid.cube contains counts per frame in all %d cubes\n", numcu);
452 printf("Raw data: average %le, min %le, max %le\n", 1.0/norm, static_cast<double>(minval)/numfr, static_cast<double>(maxval)/numfr);
455 return 0;