4eb86003d4e6c0b8ba8dc112c019fbbae7b65028
[alexxy/gromacs.git] / src / gromacs / gmxana / gmx_spatial.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2007,2008,2009,2010,2011,2012,2013,2014, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35 #include "gmxpre.h"
36
37 #include <math.h>
38 #include <stdlib.h>
39
40 #include "gromacs/commandline/pargs.h"
41 #include "gromacs/legacyheaders/typedefs.h"
42 #include "gromacs/utility/smalloc.h"
43 #include "gromacs/math/vec.h"
44 #include "gromacs/fileio/tpxio.h"
45 #include "gromacs/fileio/trxio.h"
46 #include "gromacs/topology/index.h"
47 #include "gromacs/pbcutil/pbc.h"
48 #include "gromacs/pbcutil/rmpbc.h"
49 #include "gmx_ana.h"
50 #include "gromacs/legacyheaders/macros.h"
51
52 static const double bohr = 0.529177249;  /* conversion factor to compensate for VMD plugin conversion... */
53
54 int gmx_spatial(int argc, char *argv[])
55 {
56     const char     *desc[] = {
57         "[THISMODULE] calculates the spatial distribution function and ",
58         "outputs it in a form that can be read by VMD as Gaussian98 cube format. ",
59         "For a system of 32,000 atoms and a 50 ns trajectory, the SDF can be generated ",
60         "in about 30 minutes, with most of the time dedicated to the two runs through ",
61         "[TT]trjconv[tt] that are required to center everything properly. ",
62         "This also takes a whole bunch of space (3 copies of the trajectory file). ",
63         "Still, the pictures are pretty and very informative when the fitted selection is properly made. ",
64         "3-4 atoms in a widely mobile group (like a free amino acid in solution) works ",
65         "well, or select the protein backbone in a stable folded structure to get the SDF ",
66         "of solvent and look at the time-averaged solvation shell. ",
67         "It is also possible using this program to generate the SDF based on some arbitrary ",
68         "Cartesian coordinate. To do that, simply omit the preliminary [gmx-trjconv] steps. \n",
69         "USAGE: \n",
70         "1. Use [gmx-make_ndx] to create a group containing the atoms around which you want the SDF \n",
71         "2. [TT]gmx trjconv -s a.tpr -f a.tng -o b.tng -boxcenter tric -ur compact -pbc none[tt] \n",
72         "3. [TT]gmx trjconv -s a.tpr -f b.tng -o c.tng -fit rot+trans[tt] \n",
73         "4. run [THISMODULE] on the [TT]c.tng[tt] output of step #3. \n",
74         "5. Load [TT]grid.cube[tt] into VMD and view as an isosurface. \n",
75         "[BB]Note[bb] that systems such as micelles will require [TT]gmx trjconv -pbc cluster[tt] between steps 1 and 2\n",
76         "WARNINGS:[BR]",
77         "The SDF will be generated for a cube that contains all bins that have some non-zero occupancy. ",
78         "However, the preparatory [TT]-fit rot+trans[tt] option to [gmx-trjconv] implies that your system will be rotating ",
79         "and translating in space (in order that the selected group does not). Therefore the values that are ",
80         "returned will only be valid for some region around your central group/coordinate that has full overlap ",
81         "with system volume throughout the entire translated/rotated system over the course of the trajectory. ",
82         "It is up to the user to ensure that this is the case. \n",
83         "BUGS:[BR]",
84         "When the allocated memory is not large enough, a segmentation fault may occur. This is usually detected ",
85         "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)",
86         "option. However, the program does not detect all such events. If you encounter a segmentation fault, run it again ",
87         "with an increased [TT]-nab[tt] value. \n",
88         "RISKY OPTIONS:[BR]",
89         "To reduce the amount of space and time required, you can output only the coords ",
90         "that are going to be used in the first and subsequent run through [gmx-trjconv]. ",
91         "However, be sure to set the [TT]-nab[tt] option to a sufficiently high value since ",
92         "memory is allocated for cube bins based on the initial coordinates and the [TT]-nab[tt] ",
93         "option value. \n"
94     };
95
96     static gmx_bool bPBC         = FALSE;
97     static gmx_bool bSHIFT       = FALSE;
98     static int      iIGNOREOUTER = -1;   /*Positive values may help if the surface is spikey */
99     static gmx_bool bCUTDOWN     = TRUE;
100     static real     rBINWIDTH    = 0.05; /* nm */
101     static gmx_bool bCALCDIV     = TRUE;
102     static int      iNAB         = 4;
103
104     t_pargs         pa[] = {
105         { "-pbc",      FALSE, etBOOL, {&bPBC},
106           "Use periodic boundary conditions for computing distances" },
107         { "-div",      FALSE, etBOOL, {&bCALCDIV},
108           "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" },
109         { "-ign",      FALSE, etINT, {&iIGNOREOUTER},
110           "Do not display this number of outer cubes (positive values may reduce boundary speckles; -1 ensures outer surface is visible)" },
111         /*    { "-cut",      bCUTDOWN, etBOOL, {&bCUTDOWN},*/
112         /*      "Display a total cube that is of minimal size" }, */
113         { "-bin",      FALSE, etREAL, {&rBINWIDTH},
114           "Width of the bins (nm)" },
115         { "-nab",      FALSE, etINT, {&iNAB},
116           "Number of additional bins to ensure proper memory allocation" }
117     };
118
119     double          MINBIN[3];
120     double          MAXBIN[3];
121     t_topology      top;
122     int             ePBC;
123     char            title[STRLEN];
124     t_trxframe      fr;
125     rvec           *xtop, *shx[26];
126     matrix          box, box_pbc;
127     t_trxstatus    *status;
128     int             flags = TRX_READ_X;
129     t_pbc           pbc;
130     t_atoms        *atoms;
131     int             natoms;
132     char           *grpnm, *grpnmp;
133     atom_id        *index, *indexp;
134     int             i, nidx, nidxp;
135     int             v;
136     int             j, k;
137     long         ***bin = (long ***)NULL;
138     long            nbin[3];
139     FILE           *flp;
140     long            x, y, z, minx, miny, minz, maxx, maxy, maxz;
141     long            numfr, numcu;
142     long            tot, max, min;
143     double          norm;
144     output_env_t    oenv;
145     gmx_rmpbc_t     gpbc = NULL;
146
147     t_filenm        fnm[] = {
148         { efTPS,  NULL,  NULL, ffREAD }, /* this is for the topology */
149         { efTRX, "-f", NULL, ffREAD },   /* and this for the trajectory */
150         { efNDX, NULL, NULL, ffOPTRD }
151     };
152
153 #define NFILE asize(fnm)
154
155     /* This is the routine responsible for adding default options,
156      * calling the X/motif interface, etc. */
157     if (!parse_common_args(&argc, argv, PCA_CAN_TIME | PCA_CAN_VIEW,
158                            NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL, &oenv))
159     {
160         return 0;
161     }
162
163     read_tps_conf(ftp2fn(efTPS, NFILE, fnm), title, &top, &ePBC, &xtop, NULL, box, TRUE);
164     sfree(xtop);
165
166     atoms = &(top.atoms);
167     printf("Select group to generate SDF:\n");
168     get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &nidx, &index, &grpnm);
169     printf("Select group to output coords (e.g. solute):\n");
170     get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &nidxp, &indexp, &grpnmp);
171
172     /* The first time we read data is a little special */
173     natoms = read_first_frame(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &fr, flags);
174
175     /* Memory Allocation */
176     MINBIN[XX] = MAXBIN[XX] = fr.x[0][XX];
177     MINBIN[YY] = MAXBIN[YY] = fr.x[0][YY];
178     MINBIN[ZZ] = MAXBIN[ZZ] = fr.x[0][ZZ];
179     for (i = 1; i < top.atoms.nr; ++i)
180     {
181         if (fr.x[i][XX] < MINBIN[XX])
182         {
183             MINBIN[XX] = fr.x[i][XX];
184         }
185         if (fr.x[i][XX] > MAXBIN[XX])
186         {
187             MAXBIN[XX] = fr.x[i][XX];
188         }
189         if (fr.x[i][YY] < MINBIN[YY])
190         {
191             MINBIN[YY] = fr.x[i][YY];
192         }
193         if (fr.x[i][YY] > MAXBIN[YY])
194         {
195             MAXBIN[YY] = fr.x[i][YY];
196         }
197         if (fr.x[i][ZZ] < MINBIN[ZZ])
198         {
199             MINBIN[ZZ] = fr.x[i][ZZ];
200         }
201         if (fr.x[i][ZZ] > MAXBIN[ZZ])
202         {
203             MAXBIN[ZZ] = fr.x[i][ZZ];
204         }
205     }
206     for (i = ZZ; i >= XX; --i)
207     {
208         MAXBIN[i]  = (ceil((MAXBIN[i]-MINBIN[i])/rBINWIDTH)+(double)iNAB)*rBINWIDTH+MINBIN[i];
209         MINBIN[i] -= (double)iNAB*rBINWIDTH;
210         nbin[i]    = (long)ceil((MAXBIN[i]-MINBIN[i])/rBINWIDTH);
211     }
212     snew(bin, nbin[XX]);
213     for (i = 0; i < nbin[XX]; ++i)
214     {
215         snew(bin[i], nbin[YY]);
216         for (j = 0; j < nbin[YY]; ++j)
217         {
218             snew(bin[i][j], nbin[ZZ]);
219         }
220     }
221     copy_mat(box, box_pbc);
222     numfr = 0;
223     minx  = miny = minz = 999;
224     maxx  = maxy = maxz = 0;
225
226     if (bPBC)
227     {
228         gpbc = gmx_rmpbc_init(&top.idef, ePBC, natoms);
229     }
230     /* This is the main loop over frames */
231     do
232     {
233         /* Must init pbc every step because of pressure coupling */
234
235         copy_mat(box, box_pbc);
236         if (bPBC)
237         {
238             gmx_rmpbc_trxfr(gpbc, &fr);
239             set_pbc(&pbc, ePBC, box_pbc);
240         }
241
242         for (i = 0; i < nidx; i++)
243         {
244             if (fr.x[index[i]][XX] < MINBIN[XX] || fr.x[index[i]][XX] > MAXBIN[XX] ||
245                 fr.x[index[i]][YY] < MINBIN[YY] || fr.x[index[i]][YY] > MAXBIN[YY] ||
246                 fr.x[index[i]][ZZ] < MINBIN[ZZ] || fr.x[index[i]][ZZ] > MAXBIN[ZZ])
247             {
248                 printf("There was an item outside of the allocated memory. Increase the value given with the -nab option.\n");
249                 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]);
250                 printf("Memory was required for [%f,%f,%f]\n", fr.x[index[i]][XX], fr.x[index[i]][YY], fr.x[index[i]][ZZ]);
251                 exit(1);
252             }
253             x = (long)ceil((fr.x[index[i]][XX]-MINBIN[XX])/rBINWIDTH);
254             y = (long)ceil((fr.x[index[i]][YY]-MINBIN[YY])/rBINWIDTH);
255             z = (long)ceil((fr.x[index[i]][ZZ]-MINBIN[ZZ])/rBINWIDTH);
256             ++bin[x][y][z];
257             if (x < minx)
258             {
259                 minx = x;
260             }
261             if (x > maxx)
262             {
263                 maxx = x;
264             }
265             if (y < miny)
266             {
267                 miny = y;
268             }
269             if (y > maxy)
270             {
271                 maxy = y;
272             }
273             if (z < minz)
274             {
275                 minz = z;
276             }
277             if (z > maxz)
278             {
279                 maxz = z;
280             }
281         }
282         numfr++;
283         /* printf("%f\t%f\t%f\n",box[XX][XX],box[YY][YY],box[ZZ][ZZ]); */
284
285     }
286     while (read_next_frame(oenv, status, &fr));
287
288     if (bPBC)
289     {
290         gmx_rmpbc_done(gpbc);
291     }
292
293     if (!bCUTDOWN)
294     {
295         minx = miny = minz = 0;
296         maxx = nbin[XX];
297         maxy = nbin[YY];
298         maxz = nbin[ZZ];
299     }
300
301     /* OUTPUT */
302     flp = gmx_ffopen("grid.cube", "w");
303     fprintf(flp, "Spatial Distribution Function\n");
304     fprintf(flp, "test\n");
305     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);
306     fprintf(flp, "%5ld%12.6f%12.6f%12.6f\n", maxx-minx+1-(2*iIGNOREOUTER), rBINWIDTH*10./bohr, 0., 0.);
307     fprintf(flp, "%5ld%12.6f%12.6f%12.6f\n", maxy-miny+1-(2*iIGNOREOUTER), 0., rBINWIDTH*10./bohr, 0.);
308     fprintf(flp, "%5ld%12.6f%12.6f%12.6f\n", maxz-minz+1-(2*iIGNOREOUTER), 0., 0., rBINWIDTH*10./bohr);
309     for (i = 0; i < nidxp; i++)
310     {
311         v = 2;
312         if (*(top.atoms.atomname[indexp[i]][0]) == 'C')
313         {
314             v = 6;
315         }
316         if (*(top.atoms.atomname[indexp[i]][0]) == 'N')
317         {
318             v = 7;
319         }
320         if (*(top.atoms.atomname[indexp[i]][0]) == 'O')
321         {
322             v = 8;
323         }
324         if (*(top.atoms.atomname[indexp[i]][0]) == 'H')
325         {
326             v = 1;
327         }
328         if (*(top.atoms.atomname[indexp[i]][0]) == 'S')
329         {
330             v = 16;
331         }
332         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);
333     }
334
335     tot = 0;
336     for (k = 0; k < nbin[XX]; k++)
337     {
338         if (!(k < minx || k > maxx))
339         {
340             continue;
341         }
342         for (j = 0; j < nbin[YY]; j++)
343         {
344             if (!(j < miny || j > maxy))
345             {
346                 continue;
347             }
348             for (i = 0; i < nbin[ZZ]; i++)
349             {
350                 if (!(i < minz || i > maxz))
351                 {
352                     continue;
353                 }
354                 if (bin[k][j][i] != 0)
355                 {
356                     printf("A bin was not empty when it should have been empty. Programming error.\n");
357                     printf("bin[%d][%d][%d] was = %ld\n", k, j, i, bin[k][j][i]);
358                     exit(1);
359                 }
360             }
361         }
362     }
363
364     min = 999;
365     max = 0;
366     for (k = 0; k < nbin[XX]; k++)
367     {
368         if (k < minx+iIGNOREOUTER || k > maxx-iIGNOREOUTER)
369         {
370             continue;
371         }
372         for (j = 0; j < nbin[YY]; j++)
373         {
374             if (j < miny+iIGNOREOUTER || j > maxy-iIGNOREOUTER)
375             {
376                 continue;
377             }
378             for (i = 0; i < nbin[ZZ]; i++)
379             {
380                 if (i < minz+iIGNOREOUTER || i > maxz-iIGNOREOUTER)
381                 {
382                     continue;
383                 }
384                 tot += bin[k][j][i];
385                 if (bin[k][j][i] > max)
386                 {
387                     max = bin[k][j][i];
388                 }
389                 if (bin[k][j][i] < min)
390                 {
391                     min = bin[k][j][i];
392                 }
393             }
394         }
395     }
396
397     numcu = (maxx-minx+1-(2*iIGNOREOUTER))*(maxy-miny+1-(2*iIGNOREOUTER))*(maxz-minz+1-(2*iIGNOREOUTER));
398     if (bCALCDIV)
399     {
400         norm = ((double)numcu*(double)numfr) / (double)tot;
401     }
402     else
403     {
404         norm = 1.0;
405     }
406
407     for (k = 0; k < nbin[XX]; k++)
408     {
409         if (k < minx+iIGNOREOUTER || k > maxx-iIGNOREOUTER)
410         {
411             continue;
412         }
413         for (j = 0; j < nbin[YY]; j++)
414         {
415             if (j < miny+iIGNOREOUTER || j > maxy-iIGNOREOUTER)
416             {
417                 continue;
418             }
419             for (i = 0; i < nbin[ZZ]; i++)
420             {
421                 if (i < minz+iIGNOREOUTER || i > maxz-iIGNOREOUTER)
422                 {
423                     continue;
424                 }
425                 fprintf(flp, "%12.6f ", norm*(double)bin[k][j][i]/(double)numfr);
426             }
427             fprintf(flp, "\n");
428         }
429         fprintf(flp, "\n");
430     }
431     gmx_ffclose(flp);
432
433     /* printf("x=%d to %d\n",minx,maxx); */
434     /* printf("y=%d to %d\n",miny,maxy); */
435     /* printf("z=%d to %d\n",minz,maxz); */
436
437     if (bCALCDIV)
438     {
439         printf("Counts per frame in all %ld cubes divided by %le\n", numcu, 1.0/norm);
440         printf("Normalized data: average %le, min %le, max %le\n", 1.0, norm*(double)min/(double)numfr, norm*(double)max/(double)numfr);
441     }
442     else
443     {
444         printf("grid.cube contains counts per frame in all %ld cubes\n", numcu);
445         printf("Raw data: average %le, min %le, max %le\n", 1.0/norm, (double)min/(double)numfr, (double)max/(double)numfr);
446     }
447
448     return 0;
449 }