Merge branch release-4-6 into master
[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 "statutil.h"
46 #include "tpxio.h"
47 #include <math.h>
48 #include "index.h"
49 #include "pbc.h"
50 #include "rmpbc.h"
51 #include "gmx_ana.h"
52 #include "macros.h"
53
54
55 static const double bohr = 0.529177249;  /* conversion factor to compensate for VMD plugin conversion... */
56
57 static void mequit(void)
58 {
59     printf("Memory allocation error\n");
60     exit(1);
61 }
62
63 int gmx_spatial(int argc, char *argv[])
64 {
65     const char     *desc[] = {
66         "[TT]g_spatial[tt] calculates the spatial distribution function and ",
67         "outputs it in a form that can be read by VMD as Gaussian98 cube format. ",
68         "This was developed from template.c (GROMACS-3.3). ",
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 [TT]trjconv[tt] steps. \n",
79         "USAGE: \n",
80         "1. Use [TT]make_ndx[tt] to create a group containing the atoms around which you want the SDF \n",
81         "2. [TT]trjconv -s a.tpr -f a.xtc -o b.xtc -center tric -ur compact -pbc none[tt] \n",
82         "3. [TT]trjconv -s a.tpr -f b.xtc -o c.xtc -fit rot+trans[tt] \n",
83         "4. run [TT]g_spatial[tt] 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]trjconv -pbc cluster[tt] between steps 1 and 2\n",
86         "WARNINGS:[BR]",
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 [TT]trjconv[tt] 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",
93         "BUGS:[BR]",
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",
98         "RISKY OPTIONS:[BR]",
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 [TT]trjconv[tt]. ",
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] ",
103         "option value. \n"
104     };
105
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;
112     static int      iNAB         = 4;
113
114     t_pargs         pa[] = {
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" }
127     };
128
129     double          MINBIN[3];
130     double          MAXBIN[3];
131     t_topology      top;
132     int             ePBC;
133     char            title[STRLEN];
134     t_trxframe      fr;
135     rvec           *xtop, *shx[26];
136     matrix          box, box_pbc;
137     t_trxstatus    *status;
138     int             flags = TRX_READ_X;
139     t_pbc           pbc;
140     t_atoms        *atoms;
141     int             natoms;
142     char           *grpnm, *grpnmp;
143     atom_id        *index, *indexp;
144     int             i, nidx, nidxp;
145     int             v;
146     int             j, k;
147     long         ***bin = (long ***)NULL;
148     long            nbin[3];
149     FILE           *flp;
150     long            x, y, z, minx, miny, minz, maxx, maxy, maxz;
151     long            numfr, numcu;
152     long            tot, max, min;
153     double          norm;
154     output_env_t    oenv;
155     gmx_rmpbc_t     gpbc = NULL;
156
157     t_filenm        fnm[] = {
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 }
161     };
162
163 #define NFILE asize(fnm)
164
165     /* This is the routine responsible for adding default options,
166      * calling the X/motif interface, etc. */
167     parse_common_args(&argc, argv, PCA_CAN_TIME | PCA_CAN_VIEW,
168                       NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL, &oenv);
169
170     read_tps_conf(ftp2fn(efTPS, NFILE, fnm), title, &top, &ePBC, &xtop, NULL, box, TRUE);
171     sfree(xtop);
172
173     atoms = &(top.atoms);
174     printf("Select group to generate SDF:\n");
175     get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &nidx, &index, &grpnm);
176     printf("Select group to output coords (e.g. solute):\n");
177     get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &nidxp, &indexp, &grpnmp);
178
179     /* The first time we read data is a little special */
180     natoms = read_first_frame(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &fr, flags);
181
182     /* Memory Allocation */
183     MINBIN[XX] = MAXBIN[XX] = fr.x[0][XX];
184     MINBIN[YY] = MAXBIN[YY] = fr.x[0][YY];
185     MINBIN[ZZ] = MAXBIN[ZZ] = fr.x[0][ZZ];
186     for (i = 1; i < top.atoms.nr; ++i)
187     {
188         if (fr.x[i][XX] < MINBIN[XX])
189         {
190             MINBIN[XX] = fr.x[i][XX];
191         }
192         if (fr.x[i][XX] > MAXBIN[XX])
193         {
194             MAXBIN[XX] = fr.x[i][XX];
195         }
196         if (fr.x[i][YY] < MINBIN[YY])
197         {
198             MINBIN[YY] = fr.x[i][YY];
199         }
200         if (fr.x[i][YY] > MAXBIN[YY])
201         {
202             MAXBIN[YY] = fr.x[i][YY];
203         }
204         if (fr.x[i][ZZ] < MINBIN[ZZ])
205         {
206             MINBIN[ZZ] = fr.x[i][ZZ];
207         }
208         if (fr.x[i][ZZ] > MAXBIN[ZZ])
209         {
210             MAXBIN[ZZ] = fr.x[i][ZZ];
211         }
212     }
213     for (i = ZZ; i >= XX; --i)
214     {
215         MAXBIN[i]  = (ceil((MAXBIN[i]-MINBIN[i])/rBINWIDTH)+(double)iNAB)*rBINWIDTH+MINBIN[i];
216         MINBIN[i] -= (double)iNAB*rBINWIDTH;
217         nbin[i]    = (long)ceil((MAXBIN[i]-MINBIN[i])/rBINWIDTH);
218     }
219     bin = (long ***)malloc(nbin[XX]*sizeof(long **));
220     if (!bin)
221     {
222         mequit();
223     }
224     for (i = 0; i < nbin[XX]; ++i)
225     {
226         bin[i] = (long **)malloc(nbin[YY]*sizeof(long *));
227         if (!bin[i])
228         {
229             mequit();
230         }
231         for (j = 0; j < nbin[YY]; ++j)
232         {
233             bin[i][j] = (long *)calloc(nbin[ZZ], sizeof(long));
234             if (!bin[i][j])
235             {
236                 mequit();
237             }
238         }
239     }
240     copy_mat(box, box_pbc);
241     numfr = 0;
242     minx  = miny = minz = 999;
243     maxx  = maxy = maxz = 0;
244
245     if (bPBC)
246     {
247         gpbc = gmx_rmpbc_init(&top.idef, ePBC, natoms);
248     }
249     /* This is the main loop over frames */
250     do
251     {
252         /* Must init pbc every step because of pressure coupling */
253
254         copy_mat(box, box_pbc);
255         if (bPBC)
256         {
257             gmx_rmpbc_trxfr(gpbc, &fr);
258             set_pbc(&pbc, ePBC, box_pbc);
259         }
260
261         for (i = 0; i < nidx; i++)
262         {
263             if (fr.x[index[i]][XX] < MINBIN[XX] || fr.x[index[i]][XX] > MAXBIN[XX] ||
264                 fr.x[index[i]][YY] < MINBIN[YY] || fr.x[index[i]][YY] > MAXBIN[YY] ||
265                 fr.x[index[i]][ZZ] < MINBIN[ZZ] || fr.x[index[i]][ZZ] > MAXBIN[ZZ])
266             {
267                 printf("There was an item outside of the allocated memory. Increase the value given with the -nab option.\n");
268                 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]);
269                 printf("Memory was required for [%f,%f,%f]\n", fr.x[index[i]][XX], fr.x[index[i]][YY], fr.x[index[i]][ZZ]);
270                 exit(1);
271             }
272             x = (long)ceil((fr.x[index[i]][XX]-MINBIN[XX])/rBINWIDTH);
273             y = (long)ceil((fr.x[index[i]][YY]-MINBIN[YY])/rBINWIDTH);
274             z = (long)ceil((fr.x[index[i]][ZZ]-MINBIN[ZZ])/rBINWIDTH);
275             ++bin[x][y][z];
276             if (x < minx)
277             {
278                 minx = x;
279             }
280             if (x > maxx)
281             {
282                 maxx = x;
283             }
284             if (y < miny)
285             {
286                 miny = y;
287             }
288             if (y > maxy)
289             {
290                 maxy = y;
291             }
292             if (z < minz)
293             {
294                 minz = z;
295             }
296             if (z > maxz)
297             {
298                 maxz = z;
299             }
300         }
301         numfr++;
302         /* printf("%f\t%f\t%f\n",box[XX][XX],box[YY][YY],box[ZZ][ZZ]); */
303
304     }
305     while (read_next_frame(oenv, status, &fr));
306
307     if (bPBC)
308     {
309         gmx_rmpbc_done(gpbc);
310     }
311
312     if (!bCUTDOWN)
313     {
314         minx = miny = minz = 0;
315         maxx = nbin[XX];
316         maxy = nbin[YY];
317         maxz = nbin[ZZ];
318     }
319
320     /* OUTPUT */
321     flp = ffopen("grid.cube", "w");
322     fprintf(flp, "Spatial Distribution Function\n");
323     fprintf(flp, "test\n");
324     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);
325     fprintf(flp, "%5ld%12.6f%12.6f%12.6f\n", maxx-minx+1-(2*iIGNOREOUTER), rBINWIDTH*10./bohr, 0., 0.);
326     fprintf(flp, "%5ld%12.6f%12.6f%12.6f\n", maxy-miny+1-(2*iIGNOREOUTER), 0., rBINWIDTH*10./bohr, 0.);
327     fprintf(flp, "%5ld%12.6f%12.6f%12.6f\n", maxz-minz+1-(2*iIGNOREOUTER), 0., 0., rBINWIDTH*10./bohr);
328     for (i = 0; i < nidxp; i++)
329     {
330         v = 2;
331         if (*(top.atoms.atomname[indexp[i]][0]) == 'C')
332         {
333             v = 6;
334         }
335         if (*(top.atoms.atomname[indexp[i]][0]) == 'N')
336         {
337             v = 7;
338         }
339         if (*(top.atoms.atomname[indexp[i]][0]) == 'O')
340         {
341             v = 8;
342         }
343         if (*(top.atoms.atomname[indexp[i]][0]) == 'H')
344         {
345             v = 1;
346         }
347         if (*(top.atoms.atomname[indexp[i]][0]) == 'S')
348         {
349             v = 16;
350         }
351         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);
352     }
353
354     tot = 0;
355     for (k = 0; k < nbin[XX]; k++)
356     {
357         if (!(k < minx || k > maxx))
358         {
359             continue;
360         }
361         for (j = 0; j < nbin[YY]; j++)
362         {
363             if (!(j < miny || j > maxy))
364             {
365                 continue;
366             }
367             for (i = 0; i < nbin[ZZ]; i++)
368             {
369                 if (!(i < minz || i > maxz))
370                 {
371                     continue;
372                 }
373                 if (bin[k][j][i] != 0)
374                 {
375                     printf("A bin was not empty when it should have been empty. Programming error.\n");
376                     printf("bin[%d][%d][%d] was = %ld\n", k, j, i, bin[k][j][i]);
377                     exit(1);
378                 }
379             }
380         }
381     }
382
383     min = 999;
384     max = 0;
385     for (k = 0; k < nbin[XX]; k++)
386     {
387         if (k < minx+iIGNOREOUTER || k > maxx-iIGNOREOUTER)
388         {
389             continue;
390         }
391         for (j = 0; j < nbin[YY]; j++)
392         {
393             if (j < miny+iIGNOREOUTER || j > maxy-iIGNOREOUTER)
394             {
395                 continue;
396             }
397             for (i = 0; i < nbin[ZZ]; i++)
398             {
399                 if (i < minz+iIGNOREOUTER || i > maxz-iIGNOREOUTER)
400                 {
401                     continue;
402                 }
403                 tot += bin[k][j][i];
404                 if (bin[k][j][i] > max)
405                 {
406                     max = bin[k][j][i];
407                 }
408                 if (bin[k][j][i] < min)
409                 {
410                     min = bin[k][j][i];
411                 }
412             }
413         }
414     }
415
416     numcu = (maxx-minx+1-(2*iIGNOREOUTER))*(maxy-miny+1-(2*iIGNOREOUTER))*(maxz-minz+1-(2*iIGNOREOUTER));
417     if (bCALCDIV)
418     {
419         norm = ((double)numcu*(double)numfr) / (double)tot;
420     }
421     else
422     {
423         norm = 1.0;
424     }
425
426     for (k = 0; k < nbin[XX]; k++)
427     {
428         if (k < minx+iIGNOREOUTER || k > maxx-iIGNOREOUTER)
429         {
430             continue;
431         }
432         for (j = 0; j < nbin[YY]; j++)
433         {
434             if (j < miny+iIGNOREOUTER || j > maxy-iIGNOREOUTER)
435             {
436                 continue;
437             }
438             for (i = 0; i < nbin[ZZ]; i++)
439             {
440                 if (i < minz+iIGNOREOUTER || i > maxz-iIGNOREOUTER)
441                 {
442                     continue;
443                 }
444                 fprintf(flp, "%12.6f ", norm*(double)bin[k][j][i]/(double)numfr);
445             }
446             fprintf(flp, "\n");
447         }
448         fprintf(flp, "\n");
449     }
450     ffclose(flp);
451
452     /* printf("x=%d to %d\n",minx,maxx); */
453     /* printf("y=%d to %d\n",miny,maxy); */
454     /* printf("z=%d to %d\n",minz,maxz); */
455
456     if (bCALCDIV)
457     {
458         printf("Counts per frame in all %ld cubes divided by %le\n", numcu, 1.0/norm);
459         printf("Normalized data: average %le, min %le, max %le\n", 1.0, norm*(double)min/(double)numfr, norm*(double)max/(double)numfr);
460     }
461     else
462     {
463         printf("grid.cube contains counts per frame in all %ld cubes\n", numcu);
464         printf("Raw data: average %le, min %le, max %le\n", 1.0/norm, (double)min/(double)numfr, (double)max/(double)numfr);
465     }
466
467     return 0;
468 }