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