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