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
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 "[TT]g_spatial[tt] calculates the spatial distribution function and ",
68 "outputs it in a form that can be read by VMD as Gaussian98 cube format. ",
69 "This was developed from template.c (GROMACS-3.3). ",
70 "For a system of 32,000 atoms and a 50 ns trajectory, the SDF can be generated ",
71 "in about 30 minutes, with most of the time dedicated to the two runs through ",
72 "[TT]trjconv[tt] that are required to center everything properly. ",
73 "This also takes a whole bunch of space (3 copies of the [TT].xtc[tt] file). ",
74 "Still, the pictures are pretty and very informative when the fitted selection is properly made. ",
75 "3-4 atoms in a widely mobile group (like a free amino acid in solution) works ",
76 "well, or select the protein backbone in a stable folded structure to get the SDF ",
77 "of solvent and look at the time-averaged solvation shell. ",
78 "It is also possible using this program to generate the SDF based on some arbitrary ",
79 "Cartesian coordinate. To do that, simply omit the preliminary [TT]trjconv[tt] steps. \n",
81 "1. Use [TT]make_ndx[tt] to create a group containing the atoms around which you want the SDF \n",
82 "2. [TT]trjconv -s a.tpr -f a.xtc -o b.xtc -center tric -ur compact -pbc none[tt] \n",
83 "3. [TT]trjconv -s a.tpr -f b.xtc -o c.xtc -fit rot+trans[tt] \n",
84 "4. run [TT]g_spatial[tt] on the [TT].xtc[tt] output of step #3. \n",
85 "5. Load [TT]grid.cube[tt] into VMD and view as an isosurface. \n",
86 "[BB]Note[bb] that systems such as micelles will require [TT]trjconv -pbc cluster[tt] between steps 1 and 2\n",
88 "The SDF will be generated for a cube that contains all bins that have some non-zero occupancy. ",
89 "However, the preparatory [TT]-fit rot+trans[tt] option to [TT]trjconv[tt] implies that your system will be rotating ",
90 "and translating in space (in order that the selected group does not). Therefore the values that are ",
91 "returned will only be valid for some region around your central group/coordinate that has full overlap ",
92 "with system volume throughout the entire translated/rotated system over the course of the trajectory. ",
93 "It is up to the user to ensure that this is the case. \n",
95 "When the allocated memory is not large enough, a segmentation fault may occur. This is usually detected ",
96 "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)",
97 "option. However, the program does not detect all such events. If you encounter a segmentation fault, run it again ",
98 "with an increased [TT]-nab[tt] value. \n",
100 "To reduce the amount of space and time required, you can output only the coords ",
101 "that are going to be used in the first and subsequent run through [TT]trjconv[tt]. ",
102 "However, be sure to set the [TT]-nab[tt] option to a sufficiently high value since ",
103 "memory is allocated for cube bins based on the initial coordinates and the [TT]-nab[tt] ",
107 static gmx_bool bPBC = FALSE;
108 static gmx_bool bSHIFT = FALSE;
109 static int iIGNOREOUTER = -1; /*Positive values may help if the surface is spikey */
110 static gmx_bool bCUTDOWN = TRUE;
111 static real rBINWIDTH = 0.05; /* nm */
112 static gmx_bool bCALCDIV = TRUE;
116 { "-pbc", FALSE, etBOOL, {&bPBC},
117 "Use periodic boundary conditions for computing distances" },
118 { "-div", FALSE, etBOOL, {&bCALCDIV},
119 "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" },
120 { "-ign", FALSE, etINT, {&iIGNOREOUTER},
121 "Do not display this number of outer cubes (positive values may reduce boundary speckles; -1 ensures outer surface is visible)" },
122 /* { "-cut", bCUTDOWN, etBOOL, {&bCUTDOWN},*/
123 /* "Display a total cube that is of minimal size" }, */
124 { "-bin", FALSE, etREAL, {&rBINWIDTH},
125 "Width of the bins (nm)" },
126 { "-nab", FALSE, etINT, {&iNAB},
127 "Number of additional bins to ensure proper memory allocation" }
136 rvec *xtop, *shx[26];
139 int flags = TRX_READ_X;
143 char *grpnm, *grpnmp;
144 atom_id *index, *indexp;
148 long ***bin = (long ***)NULL;
151 long x, y, z, minx, miny, minz, maxx, maxy, maxz;
156 gmx_rmpbc_t gpbc = NULL;
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 parse_common_args(&argc, argv, PCA_CAN_TIME | PCA_CAN_VIEW,
169 NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL, &oenv);
171 read_tps_conf(ftp2fn(efTPS, NFILE, fnm), title, &top, &ePBC, &xtop, NULL, box, TRUE);
174 atoms = &(top.atoms);
175 printf("Select group to generate SDF:\n");
176 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &nidx, &index, &grpnm);
177 printf("Select group to output coords (e.g. solute):\n");
178 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &nidxp, &indexp, &grpnmp);
180 /* The first time we read data is a little special */
181 natoms = read_first_frame(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &fr, flags);
183 /* Memory Allocation */
184 MINBIN[XX] = MAXBIN[XX] = fr.x[0][XX];
185 MINBIN[YY] = MAXBIN[YY] = fr.x[0][YY];
186 MINBIN[ZZ] = MAXBIN[ZZ] = fr.x[0][ZZ];
187 for (i = 1; i < top.atoms.nr; ++i)
189 if (fr.x[i][XX] < MINBIN[XX])
191 MINBIN[XX] = fr.x[i][XX];
193 if (fr.x[i][XX] > MAXBIN[XX])
195 MAXBIN[XX] = fr.x[i][XX];
197 if (fr.x[i][YY] < MINBIN[YY])
199 MINBIN[YY] = fr.x[i][YY];
201 if (fr.x[i][YY] > MAXBIN[YY])
203 MAXBIN[YY] = fr.x[i][YY];
205 if (fr.x[i][ZZ] < MINBIN[ZZ])
207 MINBIN[ZZ] = fr.x[i][ZZ];
209 if (fr.x[i][ZZ] > MAXBIN[ZZ])
211 MAXBIN[ZZ] = fr.x[i][ZZ];
214 for (i = ZZ; i >= XX; --i)
216 MAXBIN[i] = (ceil((MAXBIN[i]-MINBIN[i])/rBINWIDTH)+(double)iNAB)*rBINWIDTH+MINBIN[i];
217 MINBIN[i] -= (double)iNAB*rBINWIDTH;
218 nbin[i] = (long)ceil((MAXBIN[i]-MINBIN[i])/rBINWIDTH);
220 bin = (long ***)malloc(nbin[XX]*sizeof(long **));
225 for (i = 0; i < nbin[XX]; ++i)
227 bin[i] = (long **)malloc(nbin[YY]*sizeof(long *));
232 for (j = 0; j < nbin[YY]; ++j)
234 bin[i][j] = (long *)calloc(nbin[ZZ], sizeof(long));
241 copy_mat(box, box_pbc);
243 minx = miny = minz = 999;
244 maxx = maxy = maxz = 0;
248 gpbc = gmx_rmpbc_init(&top.idef, ePBC, natoms, box);
250 /* This is the main loop over frames */
253 /* Must init pbc every step because of pressure coupling */
255 copy_mat(box, box_pbc);
258 gmx_rmpbc_trxfr(gpbc, &fr);
259 set_pbc(&pbc, ePBC, box_pbc);
262 for (i = 0; i < nidx; i++)
264 if (fr.x[index[i]][XX] < MINBIN[XX] || fr.x[index[i]][XX] > MAXBIN[XX] ||
265 fr.x[index[i]][YY] < MINBIN[YY] || fr.x[index[i]][YY] > MAXBIN[YY] ||
266 fr.x[index[i]][ZZ] < MINBIN[ZZ] || fr.x[index[i]][ZZ] > MAXBIN[ZZ])
268 printf("There was an item outside of the allocated memory. Increase the value given with the -nab option.\n");
269 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]);
270 printf("Memory was required for [%f,%f,%f]\n", fr.x[index[i]][XX], fr.x[index[i]][YY], fr.x[index[i]][ZZ]);
273 x = (long)ceil((fr.x[index[i]][XX]-MINBIN[XX])/rBINWIDTH);
274 y = (long)ceil((fr.x[index[i]][YY]-MINBIN[YY])/rBINWIDTH);
275 z = (long)ceil((fr.x[index[i]][ZZ]-MINBIN[ZZ])/rBINWIDTH);
303 /* printf("%f\t%f\t%f\n",box[XX][XX],box[YY][YY],box[ZZ][ZZ]); */
306 while (read_next_frame(oenv, status, &fr));
310 gmx_rmpbc_done(gpbc);
315 minx = miny = minz = 0;
322 flp = ffopen("grid.cube", "w");
323 fprintf(flp, "Spatial Distribution Function\n");
324 fprintf(flp, "test\n");
325 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);
326 fprintf(flp, "%5ld%12.6f%12.6f%12.6f\n", maxx-minx+1-(2*iIGNOREOUTER), rBINWIDTH*10./bohr, 0., 0.);
327 fprintf(flp, "%5ld%12.6f%12.6f%12.6f\n", maxy-miny+1-(2*iIGNOREOUTER), 0., rBINWIDTH*10./bohr, 0.);
328 fprintf(flp, "%5ld%12.6f%12.6f%12.6f\n", maxz-minz+1-(2*iIGNOREOUTER), 0., 0., rBINWIDTH*10./bohr);
329 for (i = 0; i < nidxp; i++)
332 if (*(top.atoms.atomname[indexp[i]][0]) == 'C')
336 if (*(top.atoms.atomname[indexp[i]][0]) == 'N')
340 if (*(top.atoms.atomname[indexp[i]][0]) == 'O')
344 if (*(top.atoms.atomname[indexp[i]][0]) == 'H')
348 if (*(top.atoms.atomname[indexp[i]][0]) == 'S')
352 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);
356 for (k = 0; k < nbin[XX]; k++)
358 if (!(k < minx || k > maxx))
362 for (j = 0; j < nbin[YY]; j++)
364 if (!(j < miny || j > maxy))
368 for (i = 0; i < nbin[ZZ]; i++)
370 if (!(i < minz || i > maxz))
374 if (bin[k][j][i] != 0)
376 printf("A bin was not empty when it should have been empty. Programming error.\n");
377 printf("bin[%d][%d][%d] was = %ld\n", k, j, i, bin[k][j][i]);
386 for (k = 0; k < nbin[XX]; k++)
388 if (k < minx+iIGNOREOUTER || k > maxx-iIGNOREOUTER)
392 for (j = 0; j < nbin[YY]; j++)
394 if (j < miny+iIGNOREOUTER || j > maxy-iIGNOREOUTER)
398 for (i = 0; i < nbin[ZZ]; i++)
400 if (i < minz+iIGNOREOUTER || i > maxz-iIGNOREOUTER)
405 if (bin[k][j][i] > max)
409 if (bin[k][j][i] < min)
417 numcu = (maxx-minx+1-(2*iIGNOREOUTER))*(maxy-miny+1-(2*iIGNOREOUTER))*(maxz-minz+1-(2*iIGNOREOUTER));
420 norm = ((double)numcu*(double)numfr) / (double)tot;
427 for (k = 0; k < nbin[XX]; k++)
429 if (k < minx+iIGNOREOUTER || k > maxx-iIGNOREOUTER)
433 for (j = 0; j < nbin[YY]; j++)
435 if (j < miny+iIGNOREOUTER || j > maxy-iIGNOREOUTER)
439 for (i = 0; i < nbin[ZZ]; i++)
441 if (i < minz+iIGNOREOUTER || i > maxz-iIGNOREOUTER)
445 fprintf(flp, "%12.6f ", norm*(double)bin[k][j][i]/(double)numfr);
453 /* printf("x=%d to %d\n",minx,maxx); */
454 /* printf("y=%d to %d\n",miny,maxy); */
455 /* printf("z=%d to %d\n",minz,maxz); */
459 printf("Counts per frame in all %ld cubes divided by %le\n", numcu, 1.0/norm);
460 printf("Normalized data: average %le, min %le, max %le\n", 1.0, norm*(double)min/(double)numfr, norm*(double)max/(double)numfr);
464 printf("grid.cube contains counts per frame in all %ld cubes\n", numcu);
465 printf("Raw data: average %le, min %le, max %le\n", 1.0/norm, (double)min/(double)numfr, (double)max/(double)numfr);