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"
49 #include "gromacs/pbcutil/pbc.h"
50 #include "gromacs/utility/futil.h"
51 #include "gromacs/topology/index.h"
53 #include "gromacs/fileio/tpxio.h"
54 #include "gromacs/fileio/trxio.h"
56 #include "gromacs/fileio/matio.h"
60 #include "gromacs/utility/fatalerror.h"
62 int gmx_densmap(int argc, char *argv[])
64 const char *desc[] = {
65 "[THISMODULE] computes 2D number-density maps.",
66 "It can make planar and axial-radial density maps.",
67 "The output [TT].xpm[tt] file can be visualized with for instance xv",
68 "and can be converted to postscript with [TT]xpm2ps[tt].",
69 "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].",
71 "The default analysis is a 2-D number-density map for a selected",
72 "group of atoms in the x-y plane.",
73 "The averaging direction can be changed with the option [TT]-aver[tt].",
74 "When [TT]-xmin[tt] and/or [TT]-xmax[tt] are set only atoms that are",
75 "within the limit(s) in the averaging direction are taken into account.",
76 "The grid spacing is set with the option [TT]-bin[tt].",
77 "When [TT]-n1[tt] or [TT]-n2[tt] is non-zero, the grid",
78 "size is set by this option.",
79 "Box size fluctuations are properly taken into account.",
81 "When options [TT]-amax[tt] and [TT]-rmax[tt] are set, an axial-radial",
82 "number-density map is made. Three groups should be supplied, the centers",
83 "of mass of the first two groups define the axis, the third defines the",
84 "analysis group. The axial direction goes from -amax to +amax, where",
85 "the center is defined as the midpoint between the centers of mass and",
86 "the positive direction goes from the first to the second center of mass.",
87 "The radial direction goes from 0 to rmax or from -rmax to +rmax",
88 "when the [TT]-mirror[tt] option has been set.",
90 "The normalization of the output is set with the [TT]-unit[tt] option.",
91 "The default produces a true number density. Unit [TT]nm-2[tt] leaves out",
92 "the normalization for the averaging or the angular direction.",
93 "Option [TT]count[tt] produces the count for each grid cell.",
94 "When you do not want the scale in the output to go",
95 "from zero to the maximum density, you can set the maximum",
96 "with the option [TT]-dmax[tt]."
98 static int n1 = 0, n2 = 0;
99 static real xmin = -1, xmax = -1, bin = 0.02, dmin = 0, dmax = 0, amax = 0, rmax = 0;
100 static gmx_bool bMirror = FALSE, bSums = FALSE;
101 static const char *eaver[] = { NULL, "z", "y", "x", NULL };
102 static const char *eunit[] = { NULL, "nm-3", "nm-2", "count", NULL };
105 { "-bin", FALSE, etREAL, {&bin},
107 { "-aver", FALSE, etENUM, {eaver},
108 "The direction to average over" },
109 { "-xmin", FALSE, etREAL, {&xmin},
110 "Minimum coordinate for averaging" },
111 { "-xmax", FALSE, etREAL, {&xmax},
112 "Maximum coordinate for averaging" },
113 { "-n1", FALSE, etINT, {&n1},
114 "Number of grid cells in the first direction" },
115 { "-n2", FALSE, etINT, {&n2},
116 "Number of grid cells in the second direction" },
117 { "-amax", FALSE, etREAL, {&amax},
118 "Maximum axial distance from the center"},
119 { "-rmax", FALSE, etREAL, {&rmax},
120 "Maximum radial distance" },
121 { "-mirror", FALSE, etBOOL, {&bMirror},
122 "Add the mirror image below the axial axis" },
123 { "-sums", FALSE, etBOOL, {&bSums},
124 "Print density sums (1D map) to stdout" },
125 { "-unit", FALSE, etENUM, {eunit},
126 "Unit for the output" },
127 { "-dmin", FALSE, etREAL, {&dmin},
128 "Minimum density in output"},
129 { "-dmax", FALSE, etREAL, {&dmax},
130 "Maximum density in output (0 means calculate it)"},
132 gmx_bool bXmin, bXmax, bRadial;
137 rvec *x, xcom[2], direction, center, dx;
141 int cav = 0, c1 = 0, c2 = 0, natoms;
142 char **grpname, title[256], buf[STRLEN];
144 int i, j, k, l, ngrps, anagrp, *gnx = NULL, nindex, nradial = 0, nfr, nmpower;
145 atom_id **ind = NULL, *index;
146 real **grid, maxgrid, m1, m2, box1, box2, *tickx, *tickz, invcellvol;
147 real invspa = 0, invspz = 0, axial, r, vol_old, vol, rowsum;
149 t_rgb rlo = {1, 1, 1}, rhi = {0, 0, 0};
151 const char *label[] = { "x (nm)", "y (nm)", "z (nm)" };
153 { efTRX, "-f", NULL, ffREAD },
154 { efTPS, NULL, NULL, ffOPTRD },
155 { efNDX, NULL, NULL, ffOPTRD },
156 { efDAT, "-od", "densmap", ffOPTWR },
157 { efXPM, "-o", "densmap", ffWRITE }
159 #define NFILE asize(fnm)
164 if (!parse_common_args(&argc, argv, PCA_CAN_TIME | PCA_CAN_VIEW | PCA_BE_NICE,
165 NFILE, fnm, npargs, pa, asize(desc), desc, 0, NULL, &oenv))
170 bXmin = opt2parg_bSet("-xmin", npargs, pa);
171 bXmax = opt2parg_bSet("-xmax", npargs, pa);
172 bRadial = (amax > 0 || rmax > 0);
175 if (amax <= 0 || rmax <= 0)
177 gmx_fatal(FARGS, "Both amax and rmax should be larger than zero");
181 if (strcmp(eunit[0], "nm-3") == 0)
186 else if (strcmp(eunit[0], "nm-2") == 0)
197 if (ftp2bSet(efTPS, NFILE, fnm) || !ftp2bSet(efNDX, NFILE, fnm))
199 read_tps_conf(ftp2fn(efTPS, NFILE, fnm), title, &top, &ePBC, &x, NULL, box,
205 fprintf(stderr, "\nSelect an analysis group\n");
211 "\nSelect two groups to define the axis and an analysis group\n");
214 snew(grpname, ngrps);
216 get_index(&top.atoms, ftp2fn_null(efNDX, NFILE, fnm), ngrps, gnx, ind, grpname);
218 nindex = gnx[anagrp];
222 if ((gnx[0] > 1 || gnx[1] > 1) && !ftp2bSet(efTPS, NFILE, fnm))
224 gmx_fatal(FARGS, "No run input file was supplied (option -s), this is required for the center of mass calculation");
230 case 'x': cav = XX; c1 = YY; c2 = ZZ; break;
231 case 'y': cav = YY; c1 = XX; c2 = ZZ; break;
232 case 'z': cav = ZZ; c1 = XX; c2 = YY; break;
235 natoms = read_first_x(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &t, &x, box);
241 n1 = (int)(box[c1][c1]/bin + 0.5);
245 n2 = (int)(box[c2][c2]/bin + 0.5);
250 n1 = (int)(2*amax/bin + 0.5);
251 nradial = (int)(rmax/bin + 0.5);
252 invspa = n1/(2*amax);
253 invspz = nradial/rmax;
265 for (i = 0; i < n1; i++)
282 invcellvol /= det(box);
284 else if (nmpower == -2)
286 invcellvol /= box[c1][c1]*box[c2][c2];
288 for (i = 0; i < nindex; i++)
291 if ((!bXmin || x[j][cav] >= xmin) &&
292 (!bXmax || x[j][cav] <= xmax))
294 m1 = x[j][c1]/box[c1][c1];
303 m2 = x[j][c2]/box[c2][c2];
312 grid[(int)(m1*n1)][(int)(m2*n2)] += invcellvol;
318 set_pbc(&pbc, ePBC, box);
319 for (i = 0; i < 2; i++)
323 /* One atom, just copy the coordinates */
324 copy_rvec(x[ind[i][0]], xcom[i]);
328 /* Calculate the center of mass */
331 for (j = 0; j < gnx[i]; j++)
334 m = top.atoms.atom[k].m;
335 for (l = 0; l < DIM; l++)
337 xcom[i][l] += m*x[k][l];
341 svmul(1/mtot, xcom[i], xcom[i]);
344 pbc_dx(&pbc, xcom[1], xcom[0], direction);
345 for (i = 0; i < DIM; i++)
347 center[i] = xcom[0][i] + 0.5*direction[i];
349 unitv(direction, direction);
350 for (i = 0; i < nindex; i++)
353 pbc_dx(&pbc, x[j], center, dx);
354 axial = iprod(dx, direction);
355 r = sqrt(norm2(dx) - axial*axial);
356 if (axial >= -amax && axial < amax && r < rmax)
362 grid[(int)((axial + amax)*invspa)][(int)(r*invspz)] += 1;
368 while (read_next_x(oenv, status, &t, x, box));
371 /* normalize gridpoints */
375 for (i = 0; i < n1; i++)
377 for (j = 0; j < n2; j++)
380 if (grid[i][j] > maxgrid)
382 maxgrid = grid[i][j];
389 for (i = 0; i < n1; i++)
392 for (j = 0; j < nradial; j++)
396 case -3: vol = M_PI*(j+1)*(j+1)/(invspz*invspz*invspa); break;
397 case -2: vol = (j+1)/(invspz*invspa); break;
398 default: vol = j+1; break;
408 grid[i][k] /= nfr*(vol - vol_old);
411 grid[i][nradial-1-j] = grid[i][k];
414 if (grid[i][k] > maxgrid)
416 maxgrid = grid[i][k];
421 fprintf(stdout, "\n The maximum density is %f %s\n", maxgrid, unit);
431 /* normalize box-axes */
434 for (i = 0; i <= n1; i++)
436 tickx[i] = i*box1/n1;
438 for (i = 0; i <= n2; i++)
440 tickz[i] = i*box2/n2;
445 for (i = 0; i <= n1; i++)
447 tickx[i] = i/invspa - amax;
451 for (i = 0; i <= n2; i++)
453 tickz[i] = i/invspz - rmax;
458 for (i = 0; i <= n2; i++)
467 for (i = 0; i < n1; ++i)
469 fprintf(stdout, "Density sums:\n");
471 for (j = 0; j < n2; ++j)
473 rowsum += grid[i][j];
475 fprintf(stdout, "%g\t", rowsum);
477 fprintf(stdout, "\n");
480 sprintf(buf, "%s number density", grpname[anagrp]);
481 if (!bRadial && (bXmin || bXmax))
485 sprintf(buf+strlen(buf), ", %c > %g nm", eaver[0][0], xmin);
489 sprintf(buf+strlen(buf), ", %c < %g nm", eaver[0][0], xmax);
493 sprintf(buf+strlen(buf), ", %c: %g - %g nm", eaver[0][0], xmin, xmax);
496 if (ftp2bSet(efDAT, NFILE, fnm))
498 fp = gmx_ffopen(ftp2fn(efDAT, NFILE, fnm), "w");
499 /*optional text form output: first row is tickz; first col is tickx */
501 for (j = 0; j < n2; ++j)
503 fprintf(fp, "%g\t", tickz[j]);
507 for (i = 0; i < n1; ++i)
509 fprintf(fp, "%g\t", tickx[i]);
510 for (j = 0; j < n2; ++j)
512 fprintf(fp, "%g\t", grid[i][j]);
520 fp = gmx_ffopen(ftp2fn(efXPM, NFILE, fnm), "w");
521 write_xpm(fp, MAT_SPATIAL_X | MAT_SPATIAL_Y, buf, unit,
522 bRadial ? "axial (nm)" : label[c1], bRadial ? "r (nm)" : label[c2],
523 n1, n2, tickx, tickz, grid, dmin, maxgrid, rlo, rhi, &nlev);
527 do_view(oenv, opt2fn("-o", NFILE, fnm), NULL);