Merge "Fix bug in selection subexpression handling."
[alexxy/gromacs.git] / src / gromacs / gmxana / gmx_spatial.c
1 /*
2  *
3  *                This source code is part of
4  *
5  *                 G   R   O   M   A   C   S
6  *
7  *          GROningen MAchine for Chemical Simulations
8  *
9  *                        VERSION 3.0
10  *
11  * Copyright (c) 1991-2001
12  * BIOSON Research Institute, Dept. of Biophysical Chemistry
13  * University of Groningen, The Netherlands
14  *
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.
19  *
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.
26  *
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.
29  *
30  * Do check out http://www.gromacs.org , or mail us at gromacs@gromacs.org .
31  *
32  * And Hey:
33  * Gyas ROwers Mature At Cryogenic Speed
34  */
35 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
38
39
40
41 #include "statutil.h"
42 #include "typedefs.h"
43 #include "smalloc.h"
44 #include "vec.h"
45 #include "copyrite.h"
46 #include "statutil.h"
47 #include "tpxio.h"
48 #include <math.h>
49 #include "index.h"
50 #include "pbc.h"
51 #include "rmpbc.h"
52 #include "gmx_ana.h"
53 #include "macros.h"
54
55
56 static const double bohr = 0.529177249;  /* conversion factor to compensate for VMD plugin conversion... */
57
58 static void mequit(void)
59 {
60     printf("Memory allocation error\n");
61     exit(1);
62 }
63
64 int gmx_spatial(int argc, char *argv[])
65 {
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",
80         "USAGE: \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",
87         "WARNINGS:[BR]",
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",
94         "BUGS:[BR]",
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",
99         "RISKY OPTIONS:[BR]",
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] ",
104         "option value. \n"
105     };
106
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;
113     static int      iNAB         = 4;
114
115     t_pargs         pa[] = {
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" }
128     };
129
130     double          MINBIN[3];
131     double          MAXBIN[3];
132     t_topology      top;
133     int             ePBC;
134     char            title[STRLEN];
135     t_trxframe      fr;
136     rvec           *xtop, *shx[26];
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     atom_id        *index, *indexp;
145     int             i, nidx, nidxp;
146     int             v;
147     int             j, k;
148     long         ***bin = (long ***)NULL;
149     long            nbin[3];
150     FILE           *flp;
151     long            x, y, z, minx, miny, minz, maxx, maxy, maxz;
152     long            numfr, numcu;
153     long            tot, max, min;
154     double          norm;
155     output_env_t    oenv;
156     gmx_rmpbc_t     gpbc = NULL;
157
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 }
162     };
163
164 #define NFILE asize(fnm)
165
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);
170
171     read_tps_conf(ftp2fn(efTPS, NFILE, fnm), title, &top, &ePBC, &xtop, NULL, box, TRUE);
172     sfree(xtop);
173
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);
179
180     /* The first time we read data is a little special */
181     natoms = read_first_frame(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &fr, flags);
182
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)
188     {
189         if (fr.x[i][XX] < MINBIN[XX])
190         {
191             MINBIN[XX] = fr.x[i][XX];
192         }
193         if (fr.x[i][XX] > MAXBIN[XX])
194         {
195             MAXBIN[XX] = fr.x[i][XX];
196         }
197         if (fr.x[i][YY] < MINBIN[YY])
198         {
199             MINBIN[YY] = fr.x[i][YY];
200         }
201         if (fr.x[i][YY] > MAXBIN[YY])
202         {
203             MAXBIN[YY] = fr.x[i][YY];
204         }
205         if (fr.x[i][ZZ] < MINBIN[ZZ])
206         {
207             MINBIN[ZZ] = fr.x[i][ZZ];
208         }
209         if (fr.x[i][ZZ] > MAXBIN[ZZ])
210         {
211             MAXBIN[ZZ] = fr.x[i][ZZ];
212         }
213     }
214     for (i = ZZ; i >= XX; --i)
215     {
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);
219     }
220     bin = (long ***)malloc(nbin[XX]*sizeof(long **));
221     if (!bin)
222     {
223         mequit();
224     }
225     for (i = 0; i < nbin[XX]; ++i)
226     {
227         bin[i] = (long **)malloc(nbin[YY]*sizeof(long *));
228         if (!bin[i])
229         {
230             mequit();
231         }
232         for (j = 0; j < nbin[YY]; ++j)
233         {
234             bin[i][j] = (long *)calloc(nbin[ZZ], sizeof(long));
235             if (!bin[i][j])
236             {
237                 mequit();
238             }
239         }
240     }
241     copy_mat(box, box_pbc);
242     numfr = 0;
243     minx  = miny = minz = 999;
244     maxx  = maxy = maxz = 0;
245
246     if (bPBC)
247     {
248         gpbc = gmx_rmpbc_init(&top.idef, ePBC, natoms, box);
249     }
250     /* This is the main loop over frames */
251     do
252     {
253         /* Must init pbc every step because of pressure coupling */
254
255         copy_mat(box, box_pbc);
256         if (bPBC)
257         {
258             gmx_rmpbc_trxfr(gpbc, &fr);
259             set_pbc(&pbc, ePBC, box_pbc);
260         }
261
262         for (i = 0; i < nidx; i++)
263         {
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])
267             {
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]);
271                 exit(1);
272             }
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);
276             ++bin[x][y][z];
277             if (x < minx)
278             {
279                 minx = x;
280             }
281             if (x > maxx)
282             {
283                 maxx = x;
284             }
285             if (y < miny)
286             {
287                 miny = y;
288             }
289             if (y > maxy)
290             {
291                 maxy = y;
292             }
293             if (z < minz)
294             {
295                 minz = z;
296             }
297             if (z > maxz)
298             {
299                 maxz = z;
300             }
301         }
302         numfr++;
303         /* printf("%f\t%f\t%f\n",box[XX][XX],box[YY][YY],box[ZZ][ZZ]); */
304
305     }
306     while (read_next_frame(oenv, status, &fr));
307
308     if (bPBC)
309     {
310         gmx_rmpbc_done(gpbc);
311     }
312
313     if (!bCUTDOWN)
314     {
315         minx = miny = minz = 0;
316         maxx = nbin[XX];
317         maxy = nbin[YY];
318         maxz = nbin[ZZ];
319     }
320
321     /* OUTPUT */
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++)
330     {
331         v = 2;
332         if (*(top.atoms.atomname[indexp[i]][0]) == 'C')
333         {
334             v = 6;
335         }
336         if (*(top.atoms.atomname[indexp[i]][0]) == 'N')
337         {
338             v = 7;
339         }
340         if (*(top.atoms.atomname[indexp[i]][0]) == 'O')
341         {
342             v = 8;
343         }
344         if (*(top.atoms.atomname[indexp[i]][0]) == 'H')
345         {
346             v = 1;
347         }
348         if (*(top.atoms.atomname[indexp[i]][0]) == 'S')
349         {
350             v = 16;
351         }
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);
353     }
354
355     tot = 0;
356     for (k = 0; k < nbin[XX]; k++)
357     {
358         if (!(k < minx || k > maxx))
359         {
360             continue;
361         }
362         for (j = 0; j < nbin[YY]; j++)
363         {
364             if (!(j < miny || j > maxy))
365             {
366                 continue;
367             }
368             for (i = 0; i < nbin[ZZ]; i++)
369             {
370                 if (!(i < minz || i > maxz))
371                 {
372                     continue;
373                 }
374                 if (bin[k][j][i] != 0)
375                 {
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]);
378                     exit(1);
379                 }
380             }
381         }
382     }
383
384     min = 999;
385     max = 0;
386     for (k = 0; k < nbin[XX]; k++)
387     {
388         if (k < minx+iIGNOREOUTER || k > maxx-iIGNOREOUTER)
389         {
390             continue;
391         }
392         for (j = 0; j < nbin[YY]; j++)
393         {
394             if (j < miny+iIGNOREOUTER || j > maxy-iIGNOREOUTER)
395             {
396                 continue;
397             }
398             for (i = 0; i < nbin[ZZ]; i++)
399             {
400                 if (i < minz+iIGNOREOUTER || i > maxz-iIGNOREOUTER)
401                 {
402                     continue;
403                 }
404                 tot += bin[k][j][i];
405                 if (bin[k][j][i] > max)
406                 {
407                     max = bin[k][j][i];
408                 }
409                 if (bin[k][j][i] < min)
410                 {
411                     min = bin[k][j][i];
412                 }
413             }
414         }
415     }
416
417     numcu = (maxx-minx+1-(2*iIGNOREOUTER))*(maxy-miny+1-(2*iIGNOREOUTER))*(maxz-minz+1-(2*iIGNOREOUTER));
418     if (bCALCDIV)
419     {
420         norm = ((double)numcu*(double)numfr) / (double)tot;
421     }
422     else
423     {
424         norm = 1.0;
425     }
426
427     for (k = 0; k < nbin[XX]; k++)
428     {
429         if (k < minx+iIGNOREOUTER || k > maxx-iIGNOREOUTER)
430         {
431             continue;
432         }
433         for (j = 0; j < nbin[YY]; j++)
434         {
435             if (j < miny+iIGNOREOUTER || j > maxy-iIGNOREOUTER)
436             {
437                 continue;
438             }
439             for (i = 0; i < nbin[ZZ]; i++)
440             {
441                 if (i < minz+iIGNOREOUTER || i > maxz-iIGNOREOUTER)
442                 {
443                     continue;
444                 }
445                 fprintf(flp, "%12.6f ", norm*(double)bin[k][j][i]/(double)numfr);
446             }
447             fprintf(flp, "\n");
448         }
449         fprintf(flp, "\n");
450     }
451     ffclose(flp);
452
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); */
456
457     if (bCALCDIV)
458     {
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);
461     }
462     else
463     {
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);
466     }
467
468     thanx(stderr);
469
470     return 0;
471 }