2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014, 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.
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.
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.
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.
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.
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.
44 #include "gromacs/commandline/pargs.h"
46 #include "gromacs/utility/smalloc.h"
48 #include "gromacs/math/vec.h"
50 #include "gromacs/utility/futil.h"
53 #include "gromacs/fileio/tpxio.h"
54 #include "gromacs/fileio/trxio.h"
56 #include "gromacs/fileio/matio.h"
61 #include "gromacs/utility/fatalerror.h"
63 int gmx_densmap(int argc, char *argv[])
65 const char *desc[] = {
66 "[THISMODULE] computes 2D number-density maps.",
67 "It can make planar and axial-radial density maps.",
68 "The output [TT].xpm[tt] file can be visualized with for instance xv",
69 "and can be converted to postscript with [TT]xpm2ps[tt].",
70 "Optionally, output can be in text form to a [TT].dat[tt] file with [TT]-od[tt], instead of the usual [TT].xpm[tt] file with [TT]-o[tt].",
72 "The default analysis is a 2-D number-density map for a selected",
73 "group of atoms in the x-y plane.",
74 "The averaging direction can be changed with the option [TT]-aver[tt].",
75 "When [TT]-xmin[tt] and/or [TT]-xmax[tt] are set only atoms that are",
76 "within the limit(s) in the averaging direction are taken into account.",
77 "The grid spacing is set with the option [TT]-bin[tt].",
78 "When [TT]-n1[tt] or [TT]-n2[tt] is non-zero, the grid",
79 "size is set by this option.",
80 "Box size fluctuations are properly taken into account.",
82 "When options [TT]-amax[tt] and [TT]-rmax[tt] are set, an axial-radial",
83 "number-density map is made. Three groups should be supplied, the centers",
84 "of mass of the first two groups define the axis, the third defines the",
85 "analysis group. The axial direction goes from -amax to +amax, where",
86 "the center is defined as the midpoint between the centers of mass and",
87 "the positive direction goes from the first to the second center of mass.",
88 "The radial direction goes from 0 to rmax or from -rmax to +rmax",
89 "when the [TT]-mirror[tt] option has been set.",
91 "The normalization of the output is set with the [TT]-unit[tt] option.",
92 "The default produces a true number density. Unit [TT]nm-2[tt] leaves out",
93 "the normalization for the averaging or the angular direction.",
94 "Option [TT]count[tt] produces the count for each grid cell.",
95 "When you do not want the scale in the output to go",
96 "from zero to the maximum density, you can set the maximum",
97 "with the option [TT]-dmax[tt]."
99 static int n1 = 0, n2 = 0;
100 static real xmin = -1, xmax = -1, bin = 0.02, dmin = 0, dmax = 0, amax = 0, rmax = 0;
101 static gmx_bool bMirror = FALSE, bSums = FALSE;
102 static const char *eaver[] = { NULL, "z", "y", "x", NULL };
103 static const char *eunit[] = { NULL, "nm-3", "nm-2", "count", NULL };
106 { "-bin", FALSE, etREAL, {&bin},
108 { "-aver", FALSE, etENUM, {eaver},
109 "The direction to average over" },
110 { "-xmin", FALSE, etREAL, {&xmin},
111 "Minimum coordinate for averaging" },
112 { "-xmax", FALSE, etREAL, {&xmax},
113 "Maximum coordinate for averaging" },
114 { "-n1", FALSE, etINT, {&n1},
115 "Number of grid cells in the first direction" },
116 { "-n2", FALSE, etINT, {&n2},
117 "Number of grid cells in the second direction" },
118 { "-amax", FALSE, etREAL, {&amax},
119 "Maximum axial distance from the center"},
120 { "-rmax", FALSE, etREAL, {&rmax},
121 "Maximum radial distance" },
122 { "-mirror", FALSE, etBOOL, {&bMirror},
123 "Add the mirror image below the axial axis" },
124 { "-sums", FALSE, etBOOL, {&bSums},
125 "Print density sums (1D map) to stdout" },
126 { "-unit", FALSE, etENUM, {eunit},
127 "Unit for the output" },
128 { "-dmin", FALSE, etREAL, {&dmin},
129 "Minimum density in output"},
130 { "-dmax", FALSE, etREAL, {&dmax},
131 "Maximum density in output (0 means calculate it)"},
133 gmx_bool bXmin, bXmax, bRadial;
138 rvec *x, xcom[2], direction, center, dx;
142 int cav = 0, c1 = 0, c2 = 0, natoms;
143 char **grpname, title[256], buf[STRLEN];
145 int i, j, k, l, ngrps, anagrp, *gnx = NULL, nindex, nradial = 0, nfr, nmpower;
146 atom_id **ind = NULL, *index;
147 real **grid, maxgrid, m1, m2, box1, box2, *tickx, *tickz, invcellvol;
148 real invspa = 0, invspz = 0, axial, r, vol_old, vol, rowsum;
150 t_rgb rlo = {1, 1, 1}, rhi = {0, 0, 0};
152 const char *label[] = { "x (nm)", "y (nm)", "z (nm)" };
154 { efTRX, "-f", NULL, ffREAD },
155 { efTPS, NULL, NULL, ffOPTRD },
156 { efNDX, NULL, NULL, ffOPTRD },
157 { efDAT, "-od", "densmap", ffOPTWR },
158 { efXPM, "-o", "densmap", ffWRITE }
160 #define NFILE asize(fnm)
165 if (!parse_common_args(&argc, argv, PCA_CAN_TIME | PCA_CAN_VIEW | PCA_BE_NICE,
166 NFILE, fnm, npargs, pa, asize(desc), desc, 0, NULL, &oenv))
171 bXmin = opt2parg_bSet("-xmin", npargs, pa);
172 bXmax = opt2parg_bSet("-xmax", npargs, pa);
173 bRadial = (amax > 0 || rmax > 0);
176 if (amax <= 0 || rmax <= 0)
178 gmx_fatal(FARGS, "Both amax and rmax should be larger than zero");
182 if (strcmp(eunit[0], "nm-3") == 0)
187 else if (strcmp(eunit[0], "nm-2") == 0)
198 if (ftp2bSet(efTPS, NFILE, fnm) || !ftp2bSet(efNDX, NFILE, fnm))
200 read_tps_conf(ftp2fn(efTPS, NFILE, fnm), title, &top, &ePBC, &x, NULL, box,
206 fprintf(stderr, "\nSelect an analysis group\n");
212 "\nSelect two groups to define the axis and an analysis group\n");
215 snew(grpname, ngrps);
217 get_index(&top.atoms, ftp2fn_null(efNDX, NFILE, fnm), ngrps, gnx, ind, grpname);
219 nindex = gnx[anagrp];
223 if ((gnx[0] > 1 || gnx[1] > 1) && !ftp2bSet(efTPS, NFILE, fnm))
225 gmx_fatal(FARGS, "No run input file was supplied (option -s), this is required for the center of mass calculation");
231 case 'x': cav = XX; c1 = YY; c2 = ZZ; break;
232 case 'y': cav = YY; c1 = XX; c2 = ZZ; break;
233 case 'z': cav = ZZ; c1 = XX; c2 = YY; break;
236 natoms = read_first_x(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &t, &x, box);
242 n1 = (int)(box[c1][c1]/bin + 0.5);
246 n2 = (int)(box[c2][c2]/bin + 0.5);
251 n1 = (int)(2*amax/bin + 0.5);
252 nradial = (int)(rmax/bin + 0.5);
253 invspa = n1/(2*amax);
254 invspz = nradial/rmax;
266 for (i = 0; i < n1; i++)
283 invcellvol /= det(box);
285 else if (nmpower == -2)
287 invcellvol /= box[c1][c1]*box[c2][c2];
289 for (i = 0; i < nindex; i++)
292 if ((!bXmin || x[j][cav] >= xmin) &&
293 (!bXmax || x[j][cav] <= xmax))
295 m1 = x[j][c1]/box[c1][c1];
304 m2 = x[j][c2]/box[c2][c2];
313 grid[(int)(m1*n1)][(int)(m2*n2)] += invcellvol;
319 set_pbc(&pbc, ePBC, box);
320 for (i = 0; i < 2; i++)
324 /* One atom, just copy the coordinates */
325 copy_rvec(x[ind[i][0]], xcom[i]);
329 /* Calculate the center of mass */
332 for (j = 0; j < gnx[i]; j++)
335 m = top.atoms.atom[k].m;
336 for (l = 0; l < DIM; l++)
338 xcom[i][l] += m*x[k][l];
342 svmul(1/mtot, xcom[i], xcom[i]);
345 pbc_dx(&pbc, xcom[1], xcom[0], direction);
346 for (i = 0; i < DIM; i++)
348 center[i] = xcom[0][i] + 0.5*direction[i];
350 unitv(direction, direction);
351 for (i = 0; i < nindex; i++)
354 pbc_dx(&pbc, x[j], center, dx);
355 axial = iprod(dx, direction);
356 r = sqrt(norm2(dx) - axial*axial);
357 if (axial >= -amax && axial < amax && r < rmax)
363 grid[(int)((axial + amax)*invspa)][(int)(r*invspz)] += 1;
369 while (read_next_x(oenv, status, &t, x, box));
372 /* normalize gridpoints */
376 for (i = 0; i < n1; i++)
378 for (j = 0; j < n2; j++)
381 if (grid[i][j] > maxgrid)
383 maxgrid = grid[i][j];
390 for (i = 0; i < n1; i++)
393 for (j = 0; j < nradial; j++)
397 case -3: vol = M_PI*(j+1)*(j+1)/(invspz*invspz*invspa); break;
398 case -2: vol = (j+1)/(invspz*invspa); break;
399 default: vol = j+1; break;
409 grid[i][k] /= nfr*(vol - vol_old);
412 grid[i][nradial-1-j] = grid[i][k];
415 if (grid[i][k] > maxgrid)
417 maxgrid = grid[i][k];
422 fprintf(stdout, "\n The maximum density is %f %s\n", maxgrid, unit);
432 /* normalize box-axes */
435 for (i = 0; i <= n1; i++)
437 tickx[i] = i*box1/n1;
439 for (i = 0; i <= n2; i++)
441 tickz[i] = i*box2/n2;
446 for (i = 0; i <= n1; i++)
448 tickx[i] = i/invspa - amax;
452 for (i = 0; i <= n2; i++)
454 tickz[i] = i/invspz - rmax;
459 for (i = 0; i <= n2; i++)
468 for (i = 0; i < n1; ++i)
470 fprintf(stdout, "Density sums:\n");
472 for (j = 0; j < n2; ++j)
474 rowsum += grid[i][j];
476 fprintf(stdout, "%g\t", rowsum);
478 fprintf(stdout, "\n");
481 sprintf(buf, "%s number density", grpname[anagrp]);
482 if (!bRadial && (bXmin || bXmax))
486 sprintf(buf+strlen(buf), ", %c > %g nm", eaver[0][0], xmin);
490 sprintf(buf+strlen(buf), ", %c < %g nm", eaver[0][0], xmax);
494 sprintf(buf+strlen(buf), ", %c: %g - %g nm", eaver[0][0], xmin, xmax);
497 if (ftp2bSet(efDAT, NFILE, fnm))
499 fp = gmx_ffopen(ftp2fn(efDAT, NFILE, fnm), "w");
500 /*optional text form output: first row is tickz; first col is tickx */
502 for (j = 0; j < n2; ++j)
504 fprintf(fp, "%g\t", tickz[j]);
508 for (i = 0; i < n1; ++i)
510 fprintf(fp, "%g\t", tickx[i]);
511 for (j = 0; j < n2; ++j)
513 fprintf(fp, "%g\t", grid[i][j]);
521 fp = gmx_ffopen(ftp2fn(efXPM, NFILE, fnm), "w");
522 write_xpm(fp, MAT_SPATIAL_X | MAT_SPATIAL_Y, buf, unit,
523 bRadial ? "axial (nm)" : label[c1], bRadial ? "r (nm)" : label[c2],
524 n1, n2, tickx, tickz, grid, dmin, maxgrid, rlo, rhi, &nlev);
528 do_view(oenv, opt2fn("-o", NFILE, fnm), NULL);