3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
11 * Copyright (c) 1991-2001
12 * BIOSON Research Institute, Dept. of Biophysical Chemistry
13 * University of Groningen, The Netherlands
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * Do check out http://www.gromacs.org , or mail us at gromacs@gromacs.org .
33 * Gyas ROwers Mature At Cryogenic Speed
41 #include "gromacs/commandline/pargs.h"
45 #include "gromacs/fileio/tpxio.h"
46 #include "gromacs/fileio/trxio.h"
55 static const double bohr = 0.529177249; /* conversion factor to compensate for VMD plugin conversion... */
57 static void mequit(void)
59 printf("Memory allocation error\n");
63 int gmx_spatial(int argc, char *argv[])
65 const char *desc[] = {
66 "[THISMODULE] calculates the spatial distribution function and ",
67 "outputs it in a form that can be read by VMD as Gaussian98 cube format. ",
68 "For a system of 32,000 atoms and a 50 ns trajectory, the SDF can be generated ",
69 "in about 30 minutes, with most of the time dedicated to the two runs through ",
70 "[TT]trjconv[tt] that are required to center everything properly. ",
71 "This also takes a whole bunch of space (3 copies of the [TT].xtc[tt] file). ",
72 "Still, the pictures are pretty and very informative when the fitted selection is properly made. ",
73 "3-4 atoms in a widely mobile group (like a free amino acid in solution) works ",
74 "well, or select the protein backbone in a stable folded structure to get the SDF ",
75 "of solvent and look at the time-averaged solvation shell. ",
76 "It is also possible using this program to generate the SDF based on some arbitrary ",
77 "Cartesian coordinate. To do that, simply omit the preliminary [gmx-trjconv] steps. \n",
79 "1. Use [gmx-make_ndx] to create a group containing the atoms around which you want the SDF \n",
80 "2. [TT]gmx trjconv -s a.tpr -f a.xtc -o b.xtc -boxcenter tric -ur compact -pbc none[tt] \n",
81 "3. [TT]gmx trjconv -s a.tpr -f b.xtc -o c.xtc -fit rot+trans[tt] \n",
82 "4. run [THISMODULE] on the [TT].xtc[tt] output of step #3. \n",
83 "5. Load [TT]grid.cube[tt] into VMD and view as an isosurface. \n",
84 "[BB]Note[bb] that systems such as micelles will require [TT]gmx trjconv -pbc cluster[tt] between steps 1 and 2\n",
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. \n",
93 "When the allocated memory is not large enough, a segmentation fault may occur. This is usually detected ",
94 "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)",
95 "option. However, the program does not detect all such events. If you encounter a segmentation fault, run it again ",
96 "with an increased [TT]-nab[tt] value. \n",
98 "To reduce the amount of space and time required, you can output only the coords ",
99 "that are going to be used in the first and subsequent run through [gmx-trjconv]. ",
100 "However, be sure to set the [TT]-nab[tt] option to a sufficiently high value since ",
101 "memory is allocated for cube bins based on the initial coordinates and the [TT]-nab[tt] ",
105 static gmx_bool bPBC = FALSE;
106 static gmx_bool bSHIFT = FALSE;
107 static int iIGNOREOUTER = -1; /*Positive values may help if the surface is spikey */
108 static gmx_bool bCUTDOWN = TRUE;
109 static real rBINWIDTH = 0.05; /* nm */
110 static gmx_bool bCALCDIV = TRUE;
114 { "-pbc", FALSE, etBOOL, {&bPBC},
115 "Use periodic boundary conditions for computing distances" },
116 { "-div", FALSE, etBOOL, {&bCALCDIV},
117 "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" },
118 { "-ign", FALSE, etINT, {&iIGNOREOUTER},
119 "Do not display this number of outer cubes (positive values may reduce boundary speckles; -1 ensures outer surface is visible)" },
120 /* { "-cut", bCUTDOWN, etBOOL, {&bCUTDOWN},*/
121 /* "Display a total cube that is of minimal size" }, */
122 { "-bin", FALSE, etREAL, {&rBINWIDTH},
123 "Width of the bins (nm)" },
124 { "-nab", FALSE, etINT, {&iNAB},
125 "Number of additional bins to ensure proper memory allocation" }
134 rvec *xtop, *shx[26];
137 int flags = TRX_READ_X;
141 char *grpnm, *grpnmp;
142 atom_id *index, *indexp;
146 long ***bin = (long ***)NULL;
149 long x, y, z, minx, miny, minz, maxx, maxy, maxz;
154 gmx_rmpbc_t gpbc = NULL;
157 { efTPS, NULL, NULL, ffREAD }, /* this is for the topology */
158 { efTRX, "-f", NULL, ffREAD }, /* and this for the trajectory */
159 { efNDX, NULL, NULL, ffOPTRD }
162 #define NFILE asize(fnm)
164 /* This is the routine responsible for adding default options,
165 * calling the X/motif interface, etc. */
166 if (!parse_common_args(&argc, argv, PCA_CAN_TIME | PCA_CAN_VIEW,
167 NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL, &oenv))
172 read_tps_conf(ftp2fn(efTPS, NFILE, fnm), title, &top, &ePBC, &xtop, NULL, box, TRUE);
175 atoms = &(top.atoms);
176 printf("Select group to generate SDF:\n");
177 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &nidx, &index, &grpnm);
178 printf("Select group to output coords (e.g. solute):\n");
179 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &nidxp, &indexp, &grpnmp);
181 /* The first time we read data is a little special */
182 natoms = read_first_frame(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &fr, flags);
184 /* Memory Allocation */
185 MINBIN[XX] = MAXBIN[XX] = fr.x[0][XX];
186 MINBIN[YY] = MAXBIN[YY] = fr.x[0][YY];
187 MINBIN[ZZ] = MAXBIN[ZZ] = fr.x[0][ZZ];
188 for (i = 1; i < top.atoms.nr; ++i)
190 if (fr.x[i][XX] < MINBIN[XX])
192 MINBIN[XX] = fr.x[i][XX];
194 if (fr.x[i][XX] > MAXBIN[XX])
196 MAXBIN[XX] = fr.x[i][XX];
198 if (fr.x[i][YY] < MINBIN[YY])
200 MINBIN[YY] = fr.x[i][YY];
202 if (fr.x[i][YY] > MAXBIN[YY])
204 MAXBIN[YY] = fr.x[i][YY];
206 if (fr.x[i][ZZ] < MINBIN[ZZ])
208 MINBIN[ZZ] = fr.x[i][ZZ];
210 if (fr.x[i][ZZ] > MAXBIN[ZZ])
212 MAXBIN[ZZ] = fr.x[i][ZZ];
215 for (i = ZZ; i >= XX; --i)
217 MAXBIN[i] = (ceil((MAXBIN[i]-MINBIN[i])/rBINWIDTH)+(double)iNAB)*rBINWIDTH+MINBIN[i];
218 MINBIN[i] -= (double)iNAB*rBINWIDTH;
219 nbin[i] = (long)ceil((MAXBIN[i]-MINBIN[i])/rBINWIDTH);
221 bin = (long ***)malloc(nbin[XX]*sizeof(long **));
226 for (i = 0; i < nbin[XX]; ++i)
228 bin[i] = (long **)malloc(nbin[YY]*sizeof(long *));
233 for (j = 0; j < nbin[YY]; ++j)
235 bin[i][j] = (long *)calloc(nbin[ZZ], sizeof(long));
242 copy_mat(box, box_pbc);
244 minx = miny = minz = 999;
245 maxx = maxy = maxz = 0;
249 gpbc = gmx_rmpbc_init(&top.idef, ePBC, natoms);
251 /* This is the main loop over frames */
254 /* Must init pbc every step because of pressure coupling */
256 copy_mat(box, box_pbc);
259 gmx_rmpbc_trxfr(gpbc, &fr);
260 set_pbc(&pbc, ePBC, box_pbc);
263 for (i = 0; i < nidx; i++)
265 if (fr.x[index[i]][XX] < MINBIN[XX] || fr.x[index[i]][XX] > MAXBIN[XX] ||
266 fr.x[index[i]][YY] < MINBIN[YY] || fr.x[index[i]][YY] > MAXBIN[YY] ||
267 fr.x[index[i]][ZZ] < MINBIN[ZZ] || fr.x[index[i]][ZZ] > MAXBIN[ZZ])
269 printf("There was an item outside of the allocated memory. Increase the value given with the -nab option.\n");
270 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]);
271 printf("Memory was required for [%f,%f,%f]\n", fr.x[index[i]][XX], fr.x[index[i]][YY], fr.x[index[i]][ZZ]);
274 x = (long)ceil((fr.x[index[i]][XX]-MINBIN[XX])/rBINWIDTH);
275 y = (long)ceil((fr.x[index[i]][YY]-MINBIN[YY])/rBINWIDTH);
276 z = (long)ceil((fr.x[index[i]][ZZ]-MINBIN[ZZ])/rBINWIDTH);
304 /* printf("%f\t%f\t%f\n",box[XX][XX],box[YY][YY],box[ZZ][ZZ]); */
307 while (read_next_frame(oenv, status, &fr));
311 gmx_rmpbc_done(gpbc);
316 minx = miny = minz = 0;
323 flp = ffopen("grid.cube", "w");
324 fprintf(flp, "Spatial Distribution Function\n");
325 fprintf(flp, "test\n");
326 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);
327 fprintf(flp, "%5ld%12.6f%12.6f%12.6f\n", maxx-minx+1-(2*iIGNOREOUTER), rBINWIDTH*10./bohr, 0., 0.);
328 fprintf(flp, "%5ld%12.6f%12.6f%12.6f\n", maxy-miny+1-(2*iIGNOREOUTER), 0., rBINWIDTH*10./bohr, 0.);
329 fprintf(flp, "%5ld%12.6f%12.6f%12.6f\n", maxz-minz+1-(2*iIGNOREOUTER), 0., 0., rBINWIDTH*10./bohr);
330 for (i = 0; i < nidxp; i++)
333 if (*(top.atoms.atomname[indexp[i]][0]) == 'C')
337 if (*(top.atoms.atomname[indexp[i]][0]) == 'N')
341 if (*(top.atoms.atomname[indexp[i]][0]) == 'O')
345 if (*(top.atoms.atomname[indexp[i]][0]) == 'H')
349 if (*(top.atoms.atomname[indexp[i]][0]) == 'S')
353 fprintf(flp, "%5d%12.6f%12.6f%12.6f%12.6f\n", v, 0., (double)fr.x[indexp[i]][XX]*10./bohr, (double)fr.x[indexp[i]][YY]*10./bohr, (double)fr.x[indexp[i]][ZZ]*10./bohr);
357 for (k = 0; k < nbin[XX]; k++)
359 if (!(k < minx || k > maxx))
363 for (j = 0; j < nbin[YY]; j++)
365 if (!(j < miny || j > maxy))
369 for (i = 0; i < nbin[ZZ]; i++)
371 if (!(i < minz || i > maxz))
375 if (bin[k][j][i] != 0)
377 printf("A bin was not empty when it should have been empty. Programming error.\n");
378 printf("bin[%d][%d][%d] was = %ld\n", k, j, i, bin[k][j][i]);
387 for (k = 0; k < nbin[XX]; k++)
389 if (k < minx+iIGNOREOUTER || k > maxx-iIGNOREOUTER)
393 for (j = 0; j < nbin[YY]; j++)
395 if (j < miny+iIGNOREOUTER || j > maxy-iIGNOREOUTER)
399 for (i = 0; i < nbin[ZZ]; i++)
401 if (i < minz+iIGNOREOUTER || i > maxz-iIGNOREOUTER)
406 if (bin[k][j][i] > max)
410 if (bin[k][j][i] < min)
418 numcu = (maxx-minx+1-(2*iIGNOREOUTER))*(maxy-miny+1-(2*iIGNOREOUTER))*(maxz-minz+1-(2*iIGNOREOUTER));
421 norm = ((double)numcu*(double)numfr) / (double)tot;
428 for (k = 0; k < nbin[XX]; k++)
430 if (k < minx+iIGNOREOUTER || k > maxx-iIGNOREOUTER)
434 for (j = 0; j < nbin[YY]; j++)
436 if (j < miny+iIGNOREOUTER || j > maxy-iIGNOREOUTER)
440 for (i = 0; i < nbin[ZZ]; i++)
442 if (i < minz+iIGNOREOUTER || i > maxz-iIGNOREOUTER)
446 fprintf(flp, "%12.6f ", norm*(double)bin[k][j][i]/(double)numfr);
454 /* printf("x=%d to %d\n",minx,maxx); */
455 /* printf("y=%d to %d\n",miny,maxy); */
456 /* printf("z=%d to %d\n",minz,maxz); */
460 printf("Counts per frame in all %ld cubes divided by %le\n", numcu, 1.0/norm);
461 printf("Normalized data: average %le, min %le, max %le\n", 1.0, norm*(double)min/(double)numfr, norm*(double)max/(double)numfr);
465 printf("grid.cube contains counts per frame in all %ld cubes\n", numcu);
466 printf("Raw data: average %le, min %le, max %le\n", 1.0/norm, (double)min/(double)numfr, (double)max/(double)numfr);