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