Merge release-5-0 into release-5-1
[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,2015, 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/fileio/tpxio.h"
42 #include "gromacs/fileio/trxio.h"
43 #include "gromacs/gmxana/gmx_ana.h"
44 #include "gromacs/legacyheaders/macros.h"
45 #include "gromacs/legacyheaders/typedefs.h"
46 #include "gromacs/math/vec.h"
47 #include "gromacs/pbcutil/pbc.h"
48 #include "gromacs/pbcutil/rmpbc.h"
49 #include "gromacs/topology/index.h"
50 #include "gromacs/utility/smalloc.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.",
69         "",
70         "Usage:",
71         "",
72         "1. Use [gmx-make_ndx] to create a group containing the atoms around which you want the SDF",
73         "2. [TT]gmx trjconv -s a.tpr -f a.tng -o b.tng -boxcenter tric -ur compact -pbc none[tt]",
74         "3. [TT]gmx trjconv -s a.tpr -f b.tng -o c.tng -fit rot+trans[tt]",
75         "4. run [THISMODULE] on the [TT]c.tng[tt] output of step #3.",
76         "5. Load [TT]grid.cube[tt] into VMD and view as an isosurface.",
77         "",
78         "[BB]Note[bb] that systems such as micelles will require [TT]gmx trjconv -pbc cluster[tt] between steps 1 and 2.",
79         "",
80         "Warnings",
81         "^^^^^^^^",
82         "",
83         "The SDF will be generated for a cube that contains all bins that have some non-zero occupancy.",
84         "However, the preparatory [TT]-fit rot+trans[tt] option to [gmx-trjconv] implies that your system will be rotating",
85         "and translating in space (in order that the selected group does not). Therefore the values that are",
86         "returned will only be valid for some region around your central group/coordinate that has full overlap",
87         "with system volume throughout the entire translated/rotated system over the course of the trajectory.",
88         "It is up to the user to ensure that this is the case.",
89         "",
90         "Risky options",
91         "^^^^^^^^^^^^^",
92         "",
93         "To reduce the amount of space and time required, you can output only the coords",
94         "that are going to be used in the first and subsequent run through [gmx-trjconv].",
95         "However, be sure to set the [TT]-nab[tt] option to a sufficiently high value since",
96         "memory is allocated for cube bins based on the initial coordinates and the [TT]-nab[tt]",
97         "option value."
98     };
99     const char     *bugs[] = {
100         "When the allocated memory is not large enough, a segmentation fault may occur. This is usually detected "
101         "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) "
102         "option. However, the program does not detect all such events. If you encounter a segmentation fault, run it again "
103         "with an increased [TT]-nab[tt] value."
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     if (!parse_common_args(&argc, argv, PCA_CAN_TIME | PCA_CAN_VIEW,
168                            NFILE, fnm, asize(pa), pa, asize(desc), desc, asize(bugs), bugs, &oenv))
169     {
170         return 0;
171     }
172
173     read_tps_conf(ftp2fn(efTPS, NFILE, fnm), title, &top, &ePBC, &xtop, NULL, box, TRUE);
174     sfree(xtop);
175
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);
181
182     /* The first time we read data is a little special */
183     natoms = read_first_frame(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &fr, flags);
184
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)
190     {
191         if (fr.x[i][XX] < MINBIN[XX])
192         {
193             MINBIN[XX] = fr.x[i][XX];
194         }
195         if (fr.x[i][XX] > MAXBIN[XX])
196         {
197             MAXBIN[XX] = fr.x[i][XX];
198         }
199         if (fr.x[i][YY] < MINBIN[YY])
200         {
201             MINBIN[YY] = fr.x[i][YY];
202         }
203         if (fr.x[i][YY] > MAXBIN[YY])
204         {
205             MAXBIN[YY] = fr.x[i][YY];
206         }
207         if (fr.x[i][ZZ] < MINBIN[ZZ])
208         {
209             MINBIN[ZZ] = fr.x[i][ZZ];
210         }
211         if (fr.x[i][ZZ] > MAXBIN[ZZ])
212         {
213             MAXBIN[ZZ] = fr.x[i][ZZ];
214         }
215     }
216     for (i = ZZ; i >= XX; --i)
217     {
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);
221     }
222     snew(bin, nbin[XX]);
223     for (i = 0; i < nbin[XX]; ++i)
224     {
225         snew(bin[i], nbin[YY]);
226         for (j = 0; j < nbin[YY]; ++j)
227         {
228             snew(bin[i][j], nbin[ZZ]);
229         }
230     }
231     copy_mat(box, box_pbc);
232     numfr = 0;
233     minx  = miny = minz = 999;
234     maxx  = maxy = maxz = 0;
235
236     if (bPBC)
237     {
238         gpbc = gmx_rmpbc_init(&top.idef, ePBC, natoms);
239     }
240     /* This is the main loop over frames */
241     do
242     {
243         /* Must init pbc every step because of pressure coupling */
244
245         copy_mat(box, box_pbc);
246         if (bPBC)
247         {
248             gmx_rmpbc_trxfr(gpbc, &fr);
249             set_pbc(&pbc, ePBC, box_pbc);
250         }
251
252         for (i = 0; i < nidx; i++)
253         {
254             if (fr.x[index[i]][XX] < MINBIN[XX] || fr.x[index[i]][XX] > MAXBIN[XX] ||
255                 fr.x[index[i]][YY] < MINBIN[YY] || fr.x[index[i]][YY] > MAXBIN[YY] ||
256                 fr.x[index[i]][ZZ] < MINBIN[ZZ] || fr.x[index[i]][ZZ] > MAXBIN[ZZ])
257             {
258                 printf("There was an item outside of the allocated memory. Increase the value given with the -nab option.\n");
259                 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]);
260                 printf("Memory was required for [%f,%f,%f]\n", fr.x[index[i]][XX], fr.x[index[i]][YY], fr.x[index[i]][ZZ]);
261                 exit(1);
262             }
263             x = (long)ceil((fr.x[index[i]][XX]-MINBIN[XX])/rBINWIDTH);
264             y = (long)ceil((fr.x[index[i]][YY]-MINBIN[YY])/rBINWIDTH);
265             z = (long)ceil((fr.x[index[i]][ZZ]-MINBIN[ZZ])/rBINWIDTH);
266             ++bin[x][y][z];
267             if (x < minx)
268             {
269                 minx = x;
270             }
271             if (x > maxx)
272             {
273                 maxx = x;
274             }
275             if (y < miny)
276             {
277                 miny = y;
278             }
279             if (y > maxy)
280             {
281                 maxy = y;
282             }
283             if (z < minz)
284             {
285                 minz = z;
286             }
287             if (z > maxz)
288             {
289                 maxz = z;
290             }
291         }
292         numfr++;
293         /* printf("%f\t%f\t%f\n",box[XX][XX],box[YY][YY],box[ZZ][ZZ]); */
294
295     }
296     while (read_next_frame(oenv, status, &fr));
297
298     if (bPBC)
299     {
300         gmx_rmpbc_done(gpbc);
301     }
302
303     if (!bCUTDOWN)
304     {
305         minx = miny = minz = 0;
306         maxx = nbin[XX];
307         maxy = nbin[YY];
308         maxz = nbin[ZZ];
309     }
310
311     /* OUTPUT */
312     flp = gmx_ffopen("grid.cube", "w");
313     fprintf(flp, "Spatial Distribution Function\n");
314     fprintf(flp, "test\n");
315     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);
316     fprintf(flp, "%5ld%12.6f%12.6f%12.6f\n", maxx-minx+1-(2*iIGNOREOUTER), rBINWIDTH*10./bohr, 0., 0.);
317     fprintf(flp, "%5ld%12.6f%12.6f%12.6f\n", maxy-miny+1-(2*iIGNOREOUTER), 0., rBINWIDTH*10./bohr, 0.);
318     fprintf(flp, "%5ld%12.6f%12.6f%12.6f\n", maxz-minz+1-(2*iIGNOREOUTER), 0., 0., rBINWIDTH*10./bohr);
319     for (i = 0; i < nidxp; i++)
320     {
321         v = 2;
322         if (*(top.atoms.atomname[indexp[i]][0]) == 'C')
323         {
324             v = 6;
325         }
326         if (*(top.atoms.atomname[indexp[i]][0]) == 'N')
327         {
328             v = 7;
329         }
330         if (*(top.atoms.atomname[indexp[i]][0]) == 'O')
331         {
332             v = 8;
333         }
334         if (*(top.atoms.atomname[indexp[i]][0]) == 'H')
335         {
336             v = 1;
337         }
338         if (*(top.atoms.atomname[indexp[i]][0]) == 'S')
339         {
340             v = 16;
341         }
342         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);
343     }
344
345     tot = 0;
346     for (k = 0; k < nbin[XX]; k++)
347     {
348         if (!(k < minx || k > maxx))
349         {
350             continue;
351         }
352         for (j = 0; j < nbin[YY]; j++)
353         {
354             if (!(j < miny || j > maxy))
355             {
356                 continue;
357             }
358             for (i = 0; i < nbin[ZZ]; i++)
359             {
360                 if (!(i < minz || i > maxz))
361                 {
362                     continue;
363                 }
364                 if (bin[k][j][i] != 0)
365                 {
366                     printf("A bin was not empty when it should have been empty. Programming error.\n");
367                     printf("bin[%d][%d][%d] was = %ld\n", k, j, i, bin[k][j][i]);
368                     exit(1);
369                 }
370             }
371         }
372     }
373
374     min = 999;
375     max = 0;
376     for (k = 0; k < nbin[XX]; k++)
377     {
378         if (k < minx+iIGNOREOUTER || k > maxx-iIGNOREOUTER)
379         {
380             continue;
381         }
382         for (j = 0; j < nbin[YY]; j++)
383         {
384             if (j < miny+iIGNOREOUTER || j > maxy-iIGNOREOUTER)
385             {
386                 continue;
387             }
388             for (i = 0; i < nbin[ZZ]; i++)
389             {
390                 if (i < minz+iIGNOREOUTER || i > maxz-iIGNOREOUTER)
391                 {
392                     continue;
393                 }
394                 tot += bin[k][j][i];
395                 if (bin[k][j][i] > max)
396                 {
397                     max = bin[k][j][i];
398                 }
399                 if (bin[k][j][i] < min)
400                 {
401                     min = bin[k][j][i];
402                 }
403             }
404         }
405     }
406
407     numcu = (maxx-minx+1-(2*iIGNOREOUTER))*(maxy-miny+1-(2*iIGNOREOUTER))*(maxz-minz+1-(2*iIGNOREOUTER));
408     if (bCALCDIV)
409     {
410         norm = ((double)numcu*(double)numfr) / (double)tot;
411     }
412     else
413     {
414         norm = 1.0;
415     }
416
417     for (k = 0; k < nbin[XX]; k++)
418     {
419         if (k < minx+iIGNOREOUTER || k > maxx-iIGNOREOUTER)
420         {
421             continue;
422         }
423         for (j = 0; j < nbin[YY]; j++)
424         {
425             if (j < miny+iIGNOREOUTER || j > maxy-iIGNOREOUTER)
426             {
427                 continue;
428             }
429             for (i = 0; i < nbin[ZZ]; i++)
430             {
431                 if (i < minz+iIGNOREOUTER || i > maxz-iIGNOREOUTER)
432                 {
433                     continue;
434                 }
435                 fprintf(flp, "%12.6f ", norm*(double)bin[k][j][i]/(double)numfr);
436             }
437             fprintf(flp, "\n");
438         }
439         fprintf(flp, "\n");
440     }
441     gmx_ffclose(flp);
442
443     /* printf("x=%d to %d\n",minx,maxx); */
444     /* printf("y=%d to %d\n",miny,maxy); */
445     /* printf("z=%d to %d\n",minz,maxz); */
446
447     if (bCALCDIV)
448     {
449         printf("Counts per frame in all %ld cubes divided by %le\n", numcu, 1.0/norm);
450         printf("Normalized data: average %le, min %le, max %le\n", 1.0, norm*(double)min/(double)numfr, norm*(double)max/(double)numfr);
451     }
452     else
453     {
454         printf("grid.cube contains counts per frame in all %ld cubes\n", numcu);
455         printf("Raw data: average %le, min %le, max %le\n", 1.0/norm, (double)min/(double)numfr, (double)max/(double)numfr);
456     }
457
458     return 0;
459 }