Fix gmx density for non-mass calculations
[gromacs.git] / src / gromacs / gmxana / gmx_density.cpp
blob29c8297c3cfd0734e03ef8ce6028f4a70362dddd
1 /*
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,2015,2017,2018, 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.
37 #include "gmxpre.h"
39 #include <cctype>
40 #include <cmath>
41 #include <cstdlib>
42 #include <cstring>
44 #include "gromacs/commandline/pargs.h"
45 #include "gromacs/commandline/viewit.h"
46 #include "gromacs/fileio/trxio.h"
47 #include "gromacs/fileio/xvgr.h"
48 #include "gromacs/gmxana/gmx_ana.h"
49 #include "gromacs/gmxana/gstat.h"
50 #include "gromacs/math/units.h"
51 #include "gromacs/math/vec.h"
52 #include "gromacs/pbcutil/pbc.h"
53 #include "gromacs/pbcutil/rmpbc.h"
54 #include "gromacs/topology/index.h"
55 #include "gromacs/topology/topology.h"
56 #include "gromacs/utility/arraysize.h"
57 #include "gromacs/utility/cstringutil.h"
58 #include "gromacs/utility/fatalerror.h"
59 #include "gromacs/utility/futil.h"
60 #include "gromacs/utility/gmxassert.h"
61 #include "gromacs/utility/smalloc.h"
63 typedef struct {
64 char *atomname;
65 int nr_el;
66 } t_electron;
68 /****************************************************************************/
69 /* This program calculates the partial density across the box. */
70 /* Peter Tieleman, Mei 1995 */
71 /****************************************************************************/
73 /* used for sorting the list */
74 static int compare(void *a, void *b)
76 t_electron *tmp1, *tmp2;
77 tmp1 = (t_electron *)a; tmp2 = (t_electron *)b;
79 return std::strcmp(tmp1->atomname, tmp2->atomname);
82 static int get_electrons(t_electron **eltab, const char *fn)
84 char buffer[256]; /* to read in a line */
85 char tempname[80]; /* buffer to hold name */
86 int tempnr;
88 FILE *in;
89 int nr; /* number of atomstypes to read */
90 int i;
92 if (!(in = gmx_ffopen(fn, "r")))
94 gmx_fatal(FARGS, "Couldn't open %s. Exiting.\n", fn);
97 if (nullptr == fgets(buffer, 255, in))
99 gmx_fatal(FARGS, "Error reading from file %s", fn);
102 if (sscanf(buffer, "%d", &nr) != 1)
104 gmx_fatal(FARGS, "Invalid number of atomtypes in datafile\n");
107 snew(*eltab, nr);
109 for (i = 0; i < nr; i++)
111 if (fgets(buffer, 255, in) == nullptr)
113 gmx_fatal(FARGS, "reading datafile. Check your datafile.\n");
115 if (sscanf(buffer, "%s = %d", tempname, &tempnr) != 2)
117 gmx_fatal(FARGS, "Invalid line in datafile at line %d\n", i+1);
119 (*eltab)[i].nr_el = tempnr;
120 (*eltab)[i].atomname = gmx_strdup(tempname);
122 gmx_ffclose(in);
124 /* sort the list */
125 fprintf(stderr, "Sorting list..\n");
126 qsort ((void*)*eltab, nr, sizeof(t_electron),
127 (int(*)(const void*, const void*))compare);
129 return nr;
132 static void center_coords(t_atoms *atoms, int *index_center, int ncenter,
133 matrix box, rvec x0[])
135 int i, k, m;
136 real tmass, mm;
137 rvec com, shift, box_center;
139 tmass = 0;
140 clear_rvec(com);
141 for (k = 0; (k < ncenter); k++)
143 i = index_center[k];
144 if (i >= atoms->nr)
146 gmx_fatal(FARGS, "Index %d refers to atom %d, which is larger than natoms (%d).",
147 k+1, i+1, atoms->nr);
149 mm = atoms->atom[i].m;
150 tmass += mm;
151 for (m = 0; (m < DIM); m++)
153 com[m] += mm*x0[i][m];
156 for (m = 0; (m < DIM); m++)
158 com[m] /= tmass;
160 calc_box_center(ecenterDEF, box, box_center);
161 rvec_sub(com, box_center, shift);
163 /* Important - while the center was calculated based on a group, we should move all atoms */
164 for (i = 0; (i < atoms->nr); i++)
166 rvec_dec(x0[i], shift);
170 static void calc_electron_density(const char *fn, int **index, int gnx[],
171 double ***slDensity, int *nslices, t_topology *top,
172 int ePBC,
173 int axis, int nr_grps, real *slWidth,
174 t_electron eltab[], int nr, gmx_bool bCenter,
175 int *index_center, int ncenter,
176 gmx_bool bRelative, const gmx_output_env_t *oenv)
178 rvec *x0; /* coordinates without pbc */
179 matrix box; /* box (3x3) */
180 double invvol;
181 int natoms; /* nr. atoms in trj */
182 t_trxstatus *status;
183 int i, n, /* loop indices */
184 nr_frames = 0, /* number of frames */
185 slice; /* current slice */
186 t_electron *found; /* found by bsearch */
187 t_electron sought; /* thingie thought by bsearch */
188 real boxSz, aveBox;
189 gmx_rmpbc_t gpbc = nullptr;
191 real t,
194 if (axis < 0 || axis >= DIM)
196 gmx_fatal(FARGS, "Invalid axes. Terminating\n");
199 if ((natoms = read_first_x(oenv, &status, fn, &t, &x0, box)) == 0)
201 gmx_fatal(FARGS, "Could not read coordinates from statusfile\n");
204 aveBox = 0;
206 if (!*nslices)
208 *nslices = static_cast<int>(box[axis][axis] * 10); /* default value */
209 fprintf(stderr, "\nDividing the box in %d slices\n", *nslices);
212 snew(*slDensity, nr_grps);
213 for (i = 0; i < nr_grps; i++)
215 snew((*slDensity)[i], *nslices);
218 gpbc = gmx_rmpbc_init(&top->idef, ePBC, top->atoms.nr);
219 /*********** Start processing trajectory ***********/
222 gmx_rmpbc(gpbc, natoms, box, x0);
224 /* Translate atoms so the com of the center-group is in the
225 * box geometrical center.
227 if (bCenter)
229 center_coords(&top->atoms, index_center, ncenter, box, x0);
232 invvol = *nslices/(box[XX][XX]*box[YY][YY]*box[ZZ][ZZ]);
234 if (bRelative)
236 *slWidth = 1.0/(*nslices);
237 boxSz = 1.0;
239 else
241 *slWidth = box[axis][axis]/(*nslices);
242 boxSz = box[axis][axis];
245 aveBox += box[axis][axis];
247 for (n = 0; n < nr_grps; n++)
249 for (i = 0; i < gnx[n]; i++) /* loop over all atoms in index file */
251 z = x0[index[n][i]][axis];
252 while (z < 0)
254 z += box[axis][axis];
256 while (z > box[axis][axis])
258 z -= box[axis][axis];
261 if (bRelative)
263 z = z/box[axis][axis];
266 /* determine which slice atom is in */
267 if (bCenter)
269 slice = static_cast<int>(std::floor( (z-(boxSz/2.0)) / (*slWidth) ) + *nslices/2);
271 else
273 slice = static_cast<int>(z / (*slWidth));
275 sought.nr_el = 0;
276 sought.atomname = gmx_strdup(*(top->atoms.atomname[index[n][i]]));
278 /* now find the number of electrons. This is not efficient. */
279 found = (t_electron *)
280 bsearch((const void *)&sought,
281 (const void *)eltab, nr, sizeof(t_electron),
282 (int(*)(const void*, const void*))compare);
284 if (found == nullptr)
286 fprintf(stderr, "Couldn't find %s. Add it to the .dat file\n",
287 *(top->atoms.atomname[index[n][i]]));
289 else
291 (*slDensity)[n][slice] += (found->nr_el -
292 top->atoms.atom[index[n][i]].q)*invvol;
294 free(sought.atomname);
297 nr_frames++;
299 while (read_next_x(oenv, status, &t, x0, box));
300 gmx_rmpbc_done(gpbc);
302 /*********** done with status file **********/
303 close_trx(status);
305 /* slDensity now contains the total number of electrons per slice, summed
306 over all frames. Now divide by nr_frames and volume of slice
309 fprintf(stderr, "\nRead %d frames from trajectory. Counting electrons\n",
310 nr_frames);
312 if (bRelative)
314 aveBox /= nr_frames;
315 *slWidth = aveBox/(*nslices);
318 for (n = 0; n < nr_grps; n++)
320 for (i = 0; i < *nslices; i++)
322 (*slDensity)[n][i] /= nr_frames;
326 sfree(x0); /* free memory used by coordinate array */
329 static void calc_density(const char *fn, int **index, int gnx[],
330 double ***slDensity, int *nslices, t_topology *top, int ePBC,
331 int axis, int nr_grps, real *slWidth, gmx_bool bCenter,
332 int *index_center, int ncenter,
333 gmx_bool bRelative, const gmx_output_env_t *oenv, const char **dens_opt)
335 rvec *x0; /* coordinates without pbc */
336 matrix box; /* box (3x3) */
337 double invvol;
338 int natoms; /* nr. atoms in trj */
339 t_trxstatus *status;
340 int i, n, /* loop indices */
341 nr_frames = 0, /* number of frames */
342 slice; /* current slice */
343 real t,
345 real boxSz, aveBox;
346 real *den_val; /* values from which the density is calculated */
347 gmx_rmpbc_t gpbc = nullptr;
349 if (axis < 0 || axis >= DIM)
351 gmx_fatal(FARGS, "Invalid axes. Terminating\n");
354 if ((natoms = read_first_x(oenv, &status, fn, &t, &x0, box)) == 0)
356 gmx_fatal(FARGS, "Could not read coordinates from statusfile\n");
359 aveBox = 0;
361 if (!*nslices)
363 *nslices = static_cast<int>(box[axis][axis] * 10); /* default value */
364 fprintf(stderr, "\nDividing the box in %d slices\n", *nslices);
367 snew(*slDensity, nr_grps);
368 for (i = 0; i < nr_grps; i++)
370 snew((*slDensity)[i], *nslices);
373 gpbc = gmx_rmpbc_init(&top->idef, ePBC, top->atoms.nr);
374 /*********** Start processing trajectory ***********/
376 snew(den_val, top->atoms.nr);
377 if (dens_opt[0][0] == 'n')
379 for (i = 0; (i < top->atoms.nr); i++)
381 den_val[i] = 1;
384 else if (dens_opt[0][0] == 'c')
386 for (i = 0; (i < top->atoms.nr); i++)
388 den_val[i] = top->atoms.atom[i].q;
391 else
393 for (i = 0; (i < top->atoms.nr); i++)
395 den_val[i] = top->atoms.atom[i].m;
401 gmx_rmpbc(gpbc, natoms, box, x0);
403 /* Translate atoms so the com of the center-group is in the
404 * box geometrical center.
406 if (bCenter)
408 center_coords(&top->atoms, index_center, ncenter, box, x0);
411 invvol = *nslices/(box[XX][XX]*box[YY][YY]*box[ZZ][ZZ]);
413 if (bRelative)
415 *slWidth = 1.0/(*nslices);
416 boxSz = 1.0;
418 else
420 *slWidth = box[axis][axis]/(*nslices);
421 boxSz = box[axis][axis];
424 aveBox += box[axis][axis];
426 for (n = 0; n < nr_grps; n++)
428 for (i = 0; i < gnx[n]; i++) /* loop over all atoms in index file */
430 z = x0[index[n][i]][axis];
431 while (z < 0)
433 z += box[axis][axis];
435 while (z > box[axis][axis])
437 z -= box[axis][axis];
440 if (bRelative)
442 z = z/box[axis][axis];
445 /* determine which slice atom is in */
446 if (bCenter)
448 slice = static_cast<int>(std::floor( (z-(boxSz/2.0)) / (*slWidth) ) + *nslices/2);
450 else
452 slice = static_cast<int>(std::floor(z / (*slWidth)));
455 /* Slice should already be 0<=slice<nslices, but we just make
456 * sure we are not hit by IEEE rounding errors since we do
457 * math operations after applying PBC above.
459 if (slice < 0)
461 slice += *nslices;
463 else if (slice >= *nslices)
465 slice -= *nslices;
468 (*slDensity)[n][slice] += den_val[index[n][i]]*invvol;
471 nr_frames++;
473 while (read_next_x(oenv, status, &t, x0, box));
474 gmx_rmpbc_done(gpbc);
476 /*********** done with status file **********/
477 close_trx(status);
479 /* slDensity now contains the total mass per slice, summed over all
480 frames. Now divide by nr_frames and volume of slice
483 fprintf(stderr, "\nRead %d frames from trajectory. Calculating density\n",
484 nr_frames);
486 if (bRelative)
488 aveBox /= nr_frames;
489 *slWidth = aveBox/(*nslices);
492 for (n = 0; n < nr_grps; n++)
494 for (i = 0; i < *nslices; i++)
496 (*slDensity)[n][i] /= nr_frames;
500 sfree(x0); /* free memory used by coordinate array */
501 sfree(den_val);
504 static void plot_density(double *slDensity[], const char *afile, int nslices,
505 int nr_grps, char *grpname[], real slWidth,
506 const char **dens_opt,
507 gmx_bool bCenter, gmx_bool bRelative, gmx_bool bSymmetrize,
508 const gmx_output_env_t *oenv)
510 FILE *den;
511 const char *title = nullptr;
512 const char *xlabel = nullptr;
513 const char *ylabel = nullptr;
514 int slice, n;
515 real ddd;
516 real axispos;
518 title = bSymmetrize ? "Symmetrized partial density" : "Partial density";
520 if (bCenter)
522 xlabel = bRelative ?
523 "Average relative position from center (nm)" :
524 "Relative position from center (nm)";
526 else
528 xlabel = bRelative ? "Average coordinate (nm)" : "Coordinate (nm)";
531 switch (dens_opt[0][0])
533 case 'm': ylabel = "Density (kg m\\S-3\\N)"; break;
534 case 'n': ylabel = "Number density (nm\\S-3\\N)"; break;
535 case 'c': ylabel = "Charge density (e nm\\S-3\\N)"; break;
536 case 'e': ylabel = "Electron density (e nm\\S-3\\N)"; break;
539 den = xvgropen(afile,
540 title, xlabel, ylabel, oenv);
542 xvgr_legend(den, nr_grps, (const char**)grpname, oenv);
544 for (slice = 0; (slice < nslices); slice++)
546 if (bCenter)
548 axispos = (slice - nslices/2.0 + 0.5) * slWidth;
550 else
552 axispos = (slice + 0.5) * slWidth;
554 fprintf(den, "%12g ", axispos);
555 for (n = 0; (n < nr_grps); n++)
557 if (bSymmetrize)
559 ddd = (slDensity[n][slice]+slDensity[n][nslices-slice-1])*0.5;
561 else
563 ddd = slDensity[n][slice];
565 if (dens_opt[0][0] == 'm')
567 fprintf(den, " %12g", ddd*AMU/(NANO*NANO*NANO));
569 else
571 fprintf(den, " %12g", ddd);
574 fprintf(den, "\n");
577 xvgrclose(den);
580 int gmx_density(int argc, char *argv[])
582 const char *desc[] = {
583 "[THISMODULE] computes partial densities across the box, using an index file.[PAR]",
584 "For the total density of NPT simulations, use [gmx-energy] instead.",
585 "[PAR]",
587 "Option [TT]-center[tt] performs the histogram binning relative to the center",
588 "of an arbitrary group, in absolute box coordinates. If you are calculating",
589 "profiles along the Z axis box dimension bZ, output would be from -bZ/2 to",
590 "bZ/2 if you center based on the entire system.",
591 "Note that this behaviour has changed in GROMACS 5.0; earlier versions",
592 "merely performed a static binning in (0,bZ) and shifted the output. Now",
593 "we compute the center for each frame and bin in (-bZ/2,bZ/2).[PAR]",
595 "Option [TT]-symm[tt] symmetrizes the output around the center. This will",
596 "automatically turn on [TT]-center[tt] too.",
598 "Option [TT]-relative[tt] performs the binning in relative instead of absolute",
599 "box coordinates, and scales the final output with the average box dimension",
600 "along the output axis. This can be used in combination with [TT]-center[tt].[PAR]",
602 "Densities are in kg/m^3, and number densities or electron densities can also be",
603 "calculated. For electron densities, a file describing the number of",
604 "electrons for each type of atom should be provided using [TT]-ei[tt].",
605 "It should look like::",
607 " 2",
608 " atomname = nrelectrons",
609 " atomname = nrelectrons",
611 "The first line contains the number of lines to read from the file.",
612 "There should be one line for each unique atom name in your system.",
613 "The number of electrons for each atom is modified by its atomic",
614 "partial charge.[PAR]",
616 "IMPORTANT CONSIDERATIONS FOR BILAYERS[PAR]",
617 "One of the most common usage scenarios is to calculate the density of various",
618 "groups across a lipid bilayer, typically with the z axis being the normal",
619 "direction. For short simulations, small systems, and fixed box sizes this",
620 "will work fine, but for the more general case lipid bilayers can be complicated.",
621 "The first problem that while both proteins and lipids have low volume",
622 "compressibility, lipids have quite high area compressiblity. This means the",
623 "shape of the box (thickness and area/lipid) will fluctuate substantially even",
624 "for a fully relaxed system. Since GROMACS places the box between the origin",
625 "and positive coordinates, this in turn means that a bilayer centered in the",
626 "box will move a bit up/down due to these fluctuations, and smear out your",
627 "profile. The easiest way to fix this (if you want pressure coupling) is",
628 "to use the [TT]-center[tt] option that calculates the density profile with",
629 "respect to the center of the box. Note that you can still center on the",
630 "bilayer part even if you have a complex non-symmetric system with a bilayer",
631 "and, say, membrane proteins - then our output will simply have more values",
632 "on one side of the (center) origin reference.[PAR]",
634 "Even the centered calculation will lead to some smearing out the output",
635 "profiles, as lipids themselves are compressed and expanded. In most cases",
636 "you probably want this (since it corresponds to macroscopic experiments),",
637 "but if you want to look at molecular details you can use the [TT]-relative[tt]",
638 "option to attempt to remove even more of the effects of volume fluctuations.[PAR]",
640 "Finally, large bilayers that are not subject to a surface tension will exhibit",
641 "undulatory fluctuations, where there are 'waves' forming in the system.",
642 "This is a fundamental property of the biological system, and if you are",
643 "comparing against experiments you likely want to include the undulation",
644 "smearing effect.",
648 gmx_output_env_t *oenv;
649 static const char *dens_opt[] =
650 { nullptr, "mass", "number", "charge", "electron", nullptr };
651 static int axis = 2; /* normal to memb. default z */
652 static const char *axtitle = "Z";
653 static int nslices = 50; /* nr of slices defined */
654 static int ngrps = 1; /* nr. of groups */
655 static gmx_bool bSymmetrize = FALSE;
656 static gmx_bool bCenter = FALSE;
657 static gmx_bool bRelative = FALSE;
659 t_pargs pa[] = {
660 { "-d", FALSE, etSTR, {&axtitle},
661 "Take the normal on the membrane in direction X, Y or Z." },
662 { "-sl", FALSE, etINT, {&nslices},
663 "Divide the box in this number of slices." },
664 { "-dens", FALSE, etENUM, {dens_opt},
665 "Density"},
666 { "-ng", FALSE, etINT, {&ngrps},
667 "Number of groups of which to compute densities." },
668 { "-center", FALSE, etBOOL, {&bCenter},
669 "Perform the binning relative to the center of the (changing) box. Useful for bilayers." },
670 { "-symm", FALSE, etBOOL, {&bSymmetrize},
671 "Symmetrize the density along the axis, with respect to the center. Useful for bilayers." },
672 { "-relative", FALSE, etBOOL, {&bRelative},
673 "Use relative coordinates for changing boxes and scale output by average dimensions." }
676 const char *bugs[] = {
677 "When calculating electron densities, atomnames are used instead of types. This is bad.",
680 double **density; /* density per slice */
681 real slWidth; /* width of one slice */
682 char *grpname_center; /* centering group name */
683 char **grpname; /* groupnames */
684 int nr_electrons; /* nr. electrons */
685 int ncenter; /* size of centering group */
686 int *ngx; /* sizes of groups */
687 t_electron *el_tab; /* tabel with nr. of electrons*/
688 t_topology *top; /* topology */
689 int ePBC;
690 int *index_center; /* index for centering group */
691 int **index; /* indices for all groups */
693 t_filenm fnm[] = { /* files for g_density */
694 { efTRX, "-f", nullptr, ffREAD },
695 { efNDX, nullptr, nullptr, ffOPTRD },
696 { efTPR, nullptr, nullptr, ffREAD },
697 { efDAT, "-ei", "electrons", ffOPTRD }, /* file with nr. of electrons */
698 { efXVG, "-o", "density", ffWRITE },
701 #define NFILE asize(fnm)
703 if (!parse_common_args(&argc, argv, PCA_CAN_VIEW | PCA_CAN_TIME,
704 NFILE, fnm, asize(pa), pa, asize(desc), desc, asize(bugs), bugs,
705 &oenv))
707 return 0;
710 GMX_RELEASE_ASSERT(dens_opt[0] != nullptr, "Option setting inconsistency; dens_opt[0] is NULL");
712 if (bSymmetrize && !bCenter)
714 fprintf(stderr, "Can not symmetrize without centering. Turning on -center\n");
715 bCenter = TRUE;
717 /* Calculate axis */
718 axis = toupper(axtitle[0]) - 'X';
720 top = read_top(ftp2fn(efTPR, NFILE, fnm), &ePBC); /* read topology file */
722 snew(grpname, ngrps);
723 snew(index, ngrps);
724 snew(ngx, ngrps);
726 if (bCenter)
728 fprintf(stderr,
729 "\nNote: that the center of mass is calculated inside the box without applying\n"
730 "any special periodicity. If necessary, it is your responsibility to first use\n"
731 "trjconv to make sure atoms in this group are placed in the right periodicity.\n\n"
732 "Select the group to center density profiles around:\n");
733 get_index(&top->atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &ncenter,
734 &index_center, &grpname_center);
736 else
738 ncenter = 0;
739 index_center = nullptr;
742 fprintf(stderr, "\nSelect %d group%s to calculate density for:\n", ngrps, (ngrps > 1) ? "s" : "");
743 get_index(&top->atoms, ftp2fn_null(efNDX, NFILE, fnm), ngrps, ngx, index, grpname);
745 if (dens_opt[0][0] == 'e')
747 nr_electrons = get_electrons(&el_tab, ftp2fn(efDAT, NFILE, fnm));
748 fprintf(stderr, "Read %d atomtypes from datafile\n", nr_electrons);
750 calc_electron_density(ftp2fn(efTRX, NFILE, fnm), index, ngx, &density,
751 &nslices, top, ePBC, axis, ngrps, &slWidth, el_tab,
752 nr_electrons, bCenter, index_center, ncenter,
753 bRelative, oenv);
755 else
757 calc_density(ftp2fn(efTRX, NFILE, fnm), index, ngx, &density, &nslices, top,
758 ePBC, axis, ngrps, &slWidth, bCenter, index_center, ncenter,
759 bRelative, oenv, dens_opt);
762 plot_density(density, opt2fn("-o", NFILE, fnm),
763 nslices, ngrps, grpname, slWidth, dens_opt,
764 bCenter, bRelative, bSymmetrize, oenv);
766 do_view(oenv, opt2fn("-o", NFILE, fnm), "-nxy"); /* view xvgr file */
767 return 0;