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
46 #include "gromacs/fileio/tpxio.h"
47 #include "gromacs/fileio/trxio.h"
56 static const double bohr = 0.529177249; /* conversion factor to compensate for VMD plugin conversion... */
58 static void mequit(void)
60 printf("Memory allocation error\n");
64 int gmx_spatial(int argc, char *argv[])
66 const char *desc[] = {
67 "[THISMODULE] calculates the spatial distribution function and ",
68 "outputs it in a form that can be read by VMD as Gaussian98 cube format. ",
69 "For a system of 32,000 atoms and a 50 ns trajectory, the SDF can be generated ",
70 "in about 30 minutes, with most of the time dedicated to the two runs through ",
71 "[TT]trjconv[tt] that are required to center everything properly. ",
72 "This also takes a whole bunch of space (3 copies of the [TT].xtc[tt] file). ",
73 "Still, the pictures are pretty and very informative when the fitted selection is properly made. ",
74 "3-4 atoms in a widely mobile group (like a free amino acid in solution) works ",
75 "well, or select the protein backbone in a stable folded structure to get the SDF ",
76 "of solvent and look at the time-averaged solvation shell. ",
77 "It is also possible using this program to generate the SDF based on some arbitrary ",
78 "Cartesian coordinate. To do that, simply omit the preliminary [gmx-trjconv] steps. \n",
80 "1. Use [gmx-make_ndx] to create a group containing the atoms around which you want the SDF \n",
81 "2. [TT]gmx trjconv -s a.tpr -f a.xtc -o b.xtc -boxcenter tric -ur compact -pbc none[tt] \n",
82 "3. [TT]gmx trjconv -s a.tpr -f b.xtc -o c.xtc -fit rot+trans[tt] \n",
83 "4. run [THISMODULE] on the [TT].xtc[tt] output of step #3. \n",
84 "5. Load [TT]grid.cube[tt] into VMD and view as an isosurface. \n",
85 "[BB]Note[bb] that systems such as micelles will require [TT]gmx trjconv -pbc cluster[tt] between steps 1 and 2\n",
87 "The SDF will be generated for a cube that contains all bins that have some non-zero occupancy. ",
88 "However, the preparatory [TT]-fit rot+trans[tt] option to [gmx-trjconv] implies that your system will be rotating ",
89 "and translating in space (in order that the selected group does not). Therefore the values that are ",
90 "returned will only be valid for some region around your central group/coordinate that has full overlap ",
91 "with system volume throughout the entire translated/rotated system over the course of the trajectory. ",
92 "It is up to the user to ensure that this is the case. \n",
94 "When the allocated memory is not large enough, a segmentation fault may occur. This is usually detected ",
95 "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)",
96 "option. However, the program does not detect all such events. If you encounter a segmentation fault, run it again ",
97 "with an increased [TT]-nab[tt] value. \n",
99 "To reduce the amount of space and time required, you can output only the coords ",
100 "that are going to be used in the first and subsequent run through [gmx-trjconv]. ",
101 "However, be sure to set the [TT]-nab[tt] option to a sufficiently high value since ",
102 "memory is allocated for cube bins based on the initial coordinates and the [TT]-nab[tt] ",
106 static gmx_bool bPBC = FALSE;
107 static gmx_bool bSHIFT = FALSE;
108 static int iIGNOREOUTER = -1; /*Positive values may help if the surface is spikey */
109 static gmx_bool bCUTDOWN = TRUE;
110 static real rBINWIDTH = 0.05; /* nm */
111 static gmx_bool bCALCDIV = TRUE;
115 { "-pbc", FALSE, etBOOL, {&bPBC},
116 "Use periodic boundary conditions for computing distances" },
117 { "-div", FALSE, etBOOL, {&bCALCDIV},
118 "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" },
119 { "-ign", FALSE, etINT, {&iIGNOREOUTER},
120 "Do not display this number of outer cubes (positive values may reduce boundary speckles; -1 ensures outer surface is visible)" },
121 /* { "-cut", bCUTDOWN, etBOOL, {&bCUTDOWN},*/
122 /* "Display a total cube that is of minimal size" }, */
123 { "-bin", FALSE, etREAL, {&rBINWIDTH},
124 "Width of the bins (nm)" },
125 { "-nab", FALSE, etINT, {&iNAB},
126 "Number of additional bins to ensure proper memory allocation" }
135 rvec *xtop, *shx[26];
138 int flags = TRX_READ_X;
142 char *grpnm, *grpnmp;
143 atom_id *index, *indexp;
147 long ***bin = (long ***)NULL;
150 long x, y, z, minx, miny, minz, maxx, maxy, maxz;
155 gmx_rmpbc_t gpbc = NULL;
158 { efTPS, NULL, NULL, ffREAD }, /* this is for the topology */
159 { efTRX, "-f", NULL, ffREAD }, /* and this for the trajectory */
160 { efNDX, NULL, NULL, ffOPTRD }
163 #define NFILE asize(fnm)
165 /* This is the routine responsible for adding default options,
166 * calling the X/motif interface, etc. */
167 if (!parse_common_args(&argc, argv, PCA_CAN_TIME | PCA_CAN_VIEW,
168 NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL, &oenv))
173 read_tps_conf(ftp2fn(efTPS, NFILE, fnm), title, &top, &ePBC, &xtop, NULL, box, TRUE);
176 atoms = &(top.atoms);
177 printf("Select group to generate SDF:\n");
178 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &nidx, &index, &grpnm);
179 printf("Select group to output coords (e.g. solute):\n");
180 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &nidxp, &indexp, &grpnmp);
182 /* The first time we read data is a little special */
183 natoms = read_first_frame(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &fr, flags);
185 /* Memory Allocation */
186 MINBIN[XX] = MAXBIN[XX] = fr.x[0][XX];
187 MINBIN[YY] = MAXBIN[YY] = fr.x[0][YY];
188 MINBIN[ZZ] = MAXBIN[ZZ] = fr.x[0][ZZ];
189 for (i = 1; i < top.atoms.nr; ++i)
191 if (fr.x[i][XX] < MINBIN[XX])
193 MINBIN[XX] = fr.x[i][XX];
195 if (fr.x[i][XX] > MAXBIN[XX])
197 MAXBIN[XX] = fr.x[i][XX];
199 if (fr.x[i][YY] < MINBIN[YY])
201 MINBIN[YY] = fr.x[i][YY];
203 if (fr.x[i][YY] > MAXBIN[YY])
205 MAXBIN[YY] = fr.x[i][YY];
207 if (fr.x[i][ZZ] < MINBIN[ZZ])
209 MINBIN[ZZ] = fr.x[i][ZZ];
211 if (fr.x[i][ZZ] > MAXBIN[ZZ])
213 MAXBIN[ZZ] = fr.x[i][ZZ];
216 for (i = ZZ; i >= XX; --i)
218 MAXBIN[i] = (ceil((MAXBIN[i]-MINBIN[i])/rBINWIDTH)+(double)iNAB)*rBINWIDTH+MINBIN[i];
219 MINBIN[i] -= (double)iNAB*rBINWIDTH;
220 nbin[i] = (long)ceil((MAXBIN[i]-MINBIN[i])/rBINWIDTH);
222 bin = (long ***)malloc(nbin[XX]*sizeof(long **));
227 for (i = 0; i < nbin[XX]; ++i)
229 bin[i] = (long **)malloc(nbin[YY]*sizeof(long *));
234 for (j = 0; j < nbin[YY]; ++j)
236 bin[i][j] = (long *)calloc(nbin[ZZ], sizeof(long));
243 copy_mat(box, box_pbc);
245 minx = miny = minz = 999;
246 maxx = maxy = maxz = 0;
250 gpbc = gmx_rmpbc_init(&top.idef, ePBC, natoms);
252 /* This is the main loop over frames */
255 /* Must init pbc every step because of pressure coupling */
257 copy_mat(box, box_pbc);
260 gmx_rmpbc_trxfr(gpbc, &fr);
261 set_pbc(&pbc, ePBC, box_pbc);
264 for (i = 0; i < nidx; i++)
266 if (fr.x[index[i]][XX] < MINBIN[XX] || fr.x[index[i]][XX] > MAXBIN[XX] ||
267 fr.x[index[i]][YY] < MINBIN[YY] || fr.x[index[i]][YY] > MAXBIN[YY] ||
268 fr.x[index[i]][ZZ] < MINBIN[ZZ] || fr.x[index[i]][ZZ] > MAXBIN[ZZ])
270 printf("There was an item outside of the allocated memory. Increase the value given with the -nab option.\n");
271 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]);
272 printf("Memory was required for [%f,%f,%f]\n", fr.x[index[i]][XX], fr.x[index[i]][YY], fr.x[index[i]][ZZ]);
275 x = (long)ceil((fr.x[index[i]][XX]-MINBIN[XX])/rBINWIDTH);
276 y = (long)ceil((fr.x[index[i]][YY]-MINBIN[YY])/rBINWIDTH);
277 z = (long)ceil((fr.x[index[i]][ZZ]-MINBIN[ZZ])/rBINWIDTH);
305 /* printf("%f\t%f\t%f\n",box[XX][XX],box[YY][YY],box[ZZ][ZZ]); */
308 while (read_next_frame(oenv, status, &fr));
312 gmx_rmpbc_done(gpbc);
317 minx = miny = minz = 0;
324 flp = ffopen("grid.cube", "w");
325 fprintf(flp, "Spatial Distribution Function\n");
326 fprintf(flp, "test\n");
327 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);
328 fprintf(flp, "%5ld%12.6f%12.6f%12.6f\n", maxx-minx+1-(2*iIGNOREOUTER), rBINWIDTH*10./bohr, 0., 0.);
329 fprintf(flp, "%5ld%12.6f%12.6f%12.6f\n", maxy-miny+1-(2*iIGNOREOUTER), 0., rBINWIDTH*10./bohr, 0.);
330 fprintf(flp, "%5ld%12.6f%12.6f%12.6f\n", maxz-minz+1-(2*iIGNOREOUTER), 0., 0., rBINWIDTH*10./bohr);
331 for (i = 0; i < nidxp; i++)
334 if (*(top.atoms.atomname[indexp[i]][0]) == 'C')
338 if (*(top.atoms.atomname[indexp[i]][0]) == 'N')
342 if (*(top.atoms.atomname[indexp[i]][0]) == 'O')
346 if (*(top.atoms.atomname[indexp[i]][0]) == 'H')
350 if (*(top.atoms.atomname[indexp[i]][0]) == 'S')
354 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);
358 for (k = 0; k < nbin[XX]; k++)
360 if (!(k < minx || k > maxx))
364 for (j = 0; j < nbin[YY]; j++)
366 if (!(j < miny || j > maxy))
370 for (i = 0; i < nbin[ZZ]; i++)
372 if (!(i < minz || i > maxz))
376 if (bin[k][j][i] != 0)
378 printf("A bin was not empty when it should have been empty. Programming error.\n");
379 printf("bin[%d][%d][%d] was = %ld\n", k, j, i, bin[k][j][i]);
388 for (k = 0; k < nbin[XX]; k++)
390 if (k < minx+iIGNOREOUTER || k > maxx-iIGNOREOUTER)
394 for (j = 0; j < nbin[YY]; j++)
396 if (j < miny+iIGNOREOUTER || j > maxy-iIGNOREOUTER)
400 for (i = 0; i < nbin[ZZ]; i++)
402 if (i < minz+iIGNOREOUTER || i > maxz-iIGNOREOUTER)
407 if (bin[k][j][i] > max)
411 if (bin[k][j][i] < min)
419 numcu = (maxx-minx+1-(2*iIGNOREOUTER))*(maxy-miny+1-(2*iIGNOREOUTER))*(maxz-minz+1-(2*iIGNOREOUTER));
422 norm = ((double)numcu*(double)numfr) / (double)tot;
429 for (k = 0; k < nbin[XX]; k++)
431 if (k < minx+iIGNOREOUTER || k > maxx-iIGNOREOUTER)
435 for (j = 0; j < nbin[YY]; j++)
437 if (j < miny+iIGNOREOUTER || j > maxy-iIGNOREOUTER)
441 for (i = 0; i < nbin[ZZ]; i++)
443 if (i < minz+iIGNOREOUTER || i > maxz-iIGNOREOUTER)
447 fprintf(flp, "%12.6f ", norm*(double)bin[k][j][i]/(double)numfr);
455 /* printf("x=%d to %d\n",minx,maxx); */
456 /* printf("y=%d to %d\n",miny,maxy); */
457 /* printf("z=%d to %d\n",minz,maxz); */
461 printf("Counts per frame in all %ld cubes divided by %le\n", numcu, 1.0/norm);
462 printf("Normalized data: average %le, min %le, max %le\n", 1.0, norm*(double)min/(double)numfr, norm*(double)max/(double)numfr);
466 printf("grid.cube contains counts per frame in all %ld cubes\n", numcu);
467 printf("Raw data: average %le, min %le, max %le\n", 1.0/norm, (double)min/(double)numfr, (double)max/(double)numfr);