Introduce SimulatorBuilder
[gromacs.git] / src / gromacs / gmxana / gmx_spatial.cpp
blobbdcb6b2020bb30c5d8e09772e592912e2f22cd81
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2007,2008,2009,2010,2011,2012,2013,2014,2015,2017, 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.
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.
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.
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.
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.
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.
35 #include "gmxpre.h"
37 #include <cmath>
38 #include <cstdlib>
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/math/vec.h"
45 #include "gromacs/pbcutil/pbc.h"
46 #include "gromacs/pbcutil/rmpbc.h"
47 #include "gromacs/topology/index.h"
48 #include "gromacs/topology/topology.h"
49 #include "gromacs/trajectory/trajectoryframe.h"
50 #include "gromacs/utility/arraysize.h"
51 #include "gromacs/utility/cstringutil.h"
52 #include "gromacs/utility/futil.h"
53 #include "gromacs/utility/smalloc.h"
55 static const double bohr = 0.529177249; /* conversion factor to compensate for VMD plugin conversion... */
57 int gmx_spatial(int argc, char *argv[])
59 const char *desc[] = {
60 "[THISMODULE] calculates the spatial distribution function and",
61 "outputs it in a form that can be read by VMD as Gaussian98 cube format.",
62 "For a system of 32,000 atoms and a 50 ns trajectory, the SDF can be generated",
63 "in about 30 minutes, with most of the time dedicated to the two runs through",
64 "[TT]trjconv[tt] that are required to center everything properly.",
65 "This also takes a whole bunch of space (3 copies of the trajectory file).",
66 "Still, the pictures are pretty and very informative when the fitted selection is properly made.",
67 "3-4 atoms in a widely mobile group (like a free amino acid in solution) works",
68 "well, or select the protein backbone in a stable folded structure to get the SDF",
69 "of solvent and look at the time-averaged solvation shell.",
70 "It is also possible using this program to generate the SDF based on some arbitrary",
71 "Cartesian coordinate. To do that, simply omit the preliminary [gmx-trjconv] steps.",
72 "",
73 "Usage:",
74 "",
75 "1. Use [gmx-make_ndx] to create a group containing the atoms around which you want the SDF",
76 "2. [TT]gmx trjconv -s a.tpr -f a.tng -o b.tng -boxcenter tric -ur compact -pbc none[tt]",
77 "3. [TT]gmx trjconv -s a.tpr -f b.tng -o c.tng -fit rot+trans[tt]",
78 "4. run [THISMODULE] on the [TT]c.tng[tt] output of step #3.",
79 "5. Load [TT]grid.cube[tt] into VMD and view as an isosurface.",
80 "",
81 "[BB]Note[bb] that systems such as micelles will require [TT]gmx trjconv -pbc cluster[tt] between steps 1 and 2.",
82 "",
83 "Warnings",
84 "^^^^^^^^",
85 "",
86 "The SDF will be generated for a cube that contains all bins that have some non-zero occupancy.",
87 "However, the preparatory [TT]-fit rot+trans[tt] option to [gmx-trjconv] implies that your system will be rotating",
88 "and translating in space (in order that the selected group does not). Therefore the values that are",
89 "returned will only be valid for some region around your central group/coordinate that has full overlap",
90 "with system volume throughout the entire translated/rotated system over the course of the trajectory.",
91 "It is up to the user to ensure that this is the case.",
92 "",
93 "Risky options",
94 "^^^^^^^^^^^^^",
95 "",
96 "To reduce the amount of space and time required, you can output only the coords",
97 "that are going to be used in the first and subsequent run through [gmx-trjconv].",
98 "However, be sure to set the [TT]-nab[tt] option to a sufficiently high value since",
99 "memory is allocated for cube bins based on the initial coordinates and the [TT]-nab[tt]",
100 "option value."
102 const char *bugs[] = {
103 "When the allocated memory is not large enough, a segmentation fault may occur. This is usually detected "
104 "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) "
105 "option. However, the program does not detect all such events. If you encounter a segmentation fault, run it again "
106 "with an increased [TT]-nab[tt] value."
109 static gmx_bool bPBC = 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;
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" }
131 double MINBIN[3];
132 double MAXBIN[3];
133 t_topology top;
134 int ePBC;
135 t_trxframe fr;
136 rvec *xtop;
137 matrix box, box_pbc;
138 t_trxstatus *status;
139 int flags = TRX_READ_X;
140 t_pbc pbc;
141 t_atoms *atoms;
142 int natoms;
143 char *grpnm, *grpnmp;
144 int *index, *indexp;
145 int i, nidx, nidxp;
146 int v;
147 int j, k;
148 int ***bin = nullptr;
149 int nbin[3];
150 FILE *flp;
151 int x, y, z, minx, miny, minz, maxx, maxy, maxz;
152 int numfr, numcu;
153 int tot, maxval, minval;
154 double norm;
155 gmx_output_env_t *oenv;
156 gmx_rmpbc_t gpbc = nullptr;
158 t_filenm fnm[] = {
159 { efTPS, nullptr, nullptr, ffREAD }, /* this is for the topology */
160 { efTRX, "-f", nullptr, ffREAD }, /* and this for the trajectory */
161 { efNDX, nullptr, nullptr, ffOPTRD }
164 #define NFILE asize(fnm)
166 /* This is the routine responsible for adding default options,
167 * calling the X/motif interface, etc. */
168 if (!parse_common_args(&argc, argv, PCA_CAN_TIME | PCA_CAN_VIEW,
169 NFILE, fnm, asize(pa), pa, asize(desc), desc, asize(bugs), bugs, &oenv))
171 return 0;
174 read_tps_conf(ftp2fn(efTPS, NFILE, fnm), &top, &ePBC, &xtop, nullptr, box, TRUE);
175 sfree(xtop);
177 atoms = &(top.atoms);
178 printf("Select group to generate SDF:\n");
179 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &nidx, &index, &grpnm);
180 printf("Select group to output coords (e.g. solute):\n");
181 get_index(atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &nidxp, &indexp, &grpnmp);
183 /* The first time we read data is a little special */
184 read_first_frame(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &fr, flags);
185 natoms = fr.natoms;
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)
193 if (fr.x[i][XX] < MINBIN[XX])
195 MINBIN[XX] = fr.x[i][XX];
197 if (fr.x[i][XX] > MAXBIN[XX])
199 MAXBIN[XX] = fr.x[i][XX];
201 if (fr.x[i][YY] < MINBIN[YY])
203 MINBIN[YY] = fr.x[i][YY];
205 if (fr.x[i][YY] > MAXBIN[YY])
207 MAXBIN[YY] = fr.x[i][YY];
209 if (fr.x[i][ZZ] < MINBIN[ZZ])
211 MINBIN[ZZ] = fr.x[i][ZZ];
213 if (fr.x[i][ZZ] > MAXBIN[ZZ])
215 MAXBIN[ZZ] = fr.x[i][ZZ];
218 for (i = ZZ; i >= XX; --i)
220 MAXBIN[i] = (std::ceil((MAXBIN[i]-MINBIN[i])/rBINWIDTH)+iNAB)*rBINWIDTH+MINBIN[i];
221 MINBIN[i] -= iNAB*rBINWIDTH;
222 nbin[i] = static_cast<int>(std::ceil((MAXBIN[i]-MINBIN[i])/rBINWIDTH));
224 snew(bin, nbin[XX]);
225 for (i = 0; i < nbin[XX]; ++i)
227 snew(bin[i], nbin[YY]);
228 for (j = 0; j < nbin[YY]; ++j)
230 snew(bin[i][j], nbin[ZZ]);
233 copy_mat(box, box_pbc);
234 numfr = 0;
235 minx = miny = minz = 999;
236 maxx = maxy = maxz = 0;
238 if (bPBC)
240 gpbc = gmx_rmpbc_init(&top.idef, ePBC, natoms);
242 /* This is the main loop over frames */
245 /* Must init pbc every step because of pressure coupling */
247 copy_mat(box, box_pbc);
248 if (bPBC)
250 gmx_rmpbc_trxfr(gpbc, &fr);
251 set_pbc(&pbc, ePBC, box_pbc);
254 for (i = 0; i < nidx; i++)
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])
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);
265 x = static_cast<int>(std::ceil((fr.x[index[i]][XX]-MINBIN[XX])/rBINWIDTH));
266 y = static_cast<int>(std::ceil((fr.x[index[i]][YY]-MINBIN[YY])/rBINWIDTH));
267 z = static_cast<int>(std::ceil((fr.x[index[i]][ZZ]-MINBIN[ZZ])/rBINWIDTH));
268 ++bin[x][y][z];
269 if (x < minx)
271 minx = x;
273 if (x > maxx)
275 maxx = x;
277 if (y < miny)
279 miny = y;
281 if (y > maxy)
283 maxy = y;
285 if (z < minz)
287 minz = z;
289 if (z > maxz)
291 maxz = z;
294 numfr++;
295 /* printf("%f\t%f\t%f\n",box[XX][XX],box[YY][YY],box[ZZ][ZZ]); */
298 while (read_next_frame(oenv, status, &fr));
300 if (bPBC)
302 gmx_rmpbc_done(gpbc);
305 if (!bCUTDOWN)
307 minx = miny = minz = 0;
308 maxx = nbin[XX];
309 maxy = nbin[YY];
310 maxz = nbin[ZZ];
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, "%5d%12.6f%12.6f%12.6f\n", maxx-minx+1-(2*iIGNOREOUTER), rBINWIDTH*10./bohr, 0., 0.);
319 fprintf(flp, "%5d%12.6f%12.6f%12.6f\n", maxy-miny+1-(2*iIGNOREOUTER), 0., rBINWIDTH*10./bohr, 0.);
320 fprintf(flp, "%5d%12.6f%12.6f%12.6f\n", maxz-minz+1-(2*iIGNOREOUTER), 0., 0., rBINWIDTH*10./bohr);
321 for (i = 0; i < nidxp; i++)
323 v = 2;
324 if (*(top.atoms.atomname[indexp[i]][0]) == 'C')
326 v = 6;
328 if (*(top.atoms.atomname[indexp[i]][0]) == 'N')
330 v = 7;
332 if (*(top.atoms.atomname[indexp[i]][0]) == 'O')
334 v = 8;
336 if (*(top.atoms.atomname[indexp[i]][0]) == 'H')
338 v = 1;
340 if (*(top.atoms.atomname[indexp[i]][0]) == 'S')
342 v = 16;
344 fprintf(flp, "%5d%12.6f%12.6f%12.6f%12.6f\n", v, 0., fr.x[indexp[i]][XX]*10.0/bohr, fr.x[indexp[i]][YY]*10.0/bohr, fr.x[indexp[i]][ZZ]*10.0/bohr);
347 tot = 0;
348 for (k = 0; k < nbin[XX]; k++)
350 if (!(k < minx || k > maxx))
352 continue;
354 for (j = 0; j < nbin[YY]; j++)
356 if (!(j < miny || j > maxy))
358 continue;
360 for (i = 0; i < nbin[ZZ]; i++)
362 if (!(i < minz || i > maxz))
364 continue;
366 if (bin[k][j][i] != 0)
368 printf("A bin was not empty when it should have been empty. Programming error.\n");
369 printf("bin[%d][%d][%d] was = %d\n", k, j, i, bin[k][j][i]);
370 exit(1);
376 minval = 999;
377 maxval = 0;
378 for (k = 0; k < nbin[XX]; k++)
380 if (k < minx+iIGNOREOUTER || k > maxx-iIGNOREOUTER)
382 continue;
384 for (j = 0; j < nbin[YY]; j++)
386 if (j < miny+iIGNOREOUTER || j > maxy-iIGNOREOUTER)
388 continue;
390 for (i = 0; i < nbin[ZZ]; i++)
392 if (i < minz+iIGNOREOUTER || i > maxz-iIGNOREOUTER)
394 continue;
396 tot += bin[k][j][i];
397 if (bin[k][j][i] > maxval)
399 maxval = bin[k][j][i];
401 if (bin[k][j][i] < minval)
403 minval = bin[k][j][i];
409 numcu = (maxx-minx+1-(2*iIGNOREOUTER))*(maxy-miny+1-(2*iIGNOREOUTER))*(maxz-minz+1-(2*iIGNOREOUTER));
410 if (bCALCDIV)
412 norm = static_cast<double>(numcu*numfr)/tot;
414 else
416 norm = 1.0;
419 for (k = 0; k < nbin[XX]; k++)
421 if (k < minx+iIGNOREOUTER || k > maxx-iIGNOREOUTER)
423 continue;
425 for (j = 0; j < nbin[YY]; j++)
427 if (j < miny+iIGNOREOUTER || j > maxy-iIGNOREOUTER)
429 continue;
431 for (i = 0; i < nbin[ZZ]; i++)
433 if (i < minz+iIGNOREOUTER || i > maxz-iIGNOREOUTER)
435 continue;
437 fprintf(flp, "%12.6f ", static_cast<double>(norm*bin[k][j][i])/numfr);
439 fprintf(flp, "\n");
441 fprintf(flp, "\n");
443 gmx_ffclose(flp);
445 if (bCALCDIV)
447 printf("Counts per frame in all %d cubes divided by %le\n", numcu, 1.0/norm);
448 printf("Normalized data: average %le, min %le, max %le\n", 1.0, minval*norm/numfr, maxval*norm/numfr);
450 else
452 printf("grid.cube contains counts per frame in all %d cubes\n", numcu);
453 printf("Raw data: average %le, min %le, max %le\n", 1.0/norm, static_cast<double>(minval)/numfr, static_cast<double>(maxval)/numfr);
456 return 0;