3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2004, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
33 * Green Red Orange Magenta Azure Cyan Skyblue
59 #define EPS0 8.85419E-12
60 #define ELC 1.60219E-19
62 /****************************************************************************/
63 /* This program calculates the electrostatic potential across the box by */
64 /* determining the charge density in slices of the box and integrating these*/
66 /* Peter Tieleman, April 1995 */
67 /* It now also calculates electrostatic potential in spherical micelles, */
68 /* using \frac{1}{r}\frac{d^2r\Psi}{r^2} = - \frac{\rho}{\epsilon_0} */
69 /* This probably sucks but it seems to work. */
70 /****************************************************************************/
72 static int ce = 0, cb = 0;
74 /* this routine integrates the array data and returns the resulting array */
75 /* routine uses simple trapezoid rule */
76 void p_integrate(double *result, double data[], int ndata, double slWidth)
83 fprintf(stderr, "Warning: nr of slices very small. This will result"
87 fprintf(stderr, "Integrating from slice %d to slice %d\n", cb, ndata-ce);
89 for (slice = cb; slice < (ndata-ce); slice++)
92 for (i = cb; i < slice; i++)
94 sum += slWidth * (data[i] + 0.5 * (data[i+1] - data[i]));
101 void calc_potential(const char *fn, atom_id **index, int gnx[],
102 double ***slPotential, double ***slCharge,
103 double ***slField, int *nslices,
104 t_topology *top, int ePBC,
105 int axis, int nr_grps, double *slWidth,
106 double fudge_z, gmx_bool bSpherical, gmx_bool bCorrect,
107 const output_env_t oenv)
109 rvec *x0; /* coordinates without pbc */
110 matrix box; /* box (3x3) */
111 int natoms; /* nr. atoms in trj */
113 int **slCount, /* nr. of atoms in one slice for a group */
114 i, j, n, /* loop indices */
117 nr_frames = 0, /* number of frames */
118 slice; /* current slice */
119 double slVolume; /* volume of slice for spherical averaging */
124 gmx_rmpbc_t gpbc = NULL;
138 gmx_fatal(FARGS, "Invalid axes. Terminating\n");
141 if ((natoms = read_first_x(oenv, &status, fn, &t, &x0, box)) == 0)
143 gmx_fatal(FARGS, "Could not read coordinates from statusfile\n");
148 *nslices = (int)(box[axis][axis] * 10); /* default value */
151 fprintf(stderr, "\nDividing the box in %d slices\n", *nslices);
153 snew(*slField, nr_grps);
154 snew(*slCharge, nr_grps);
155 snew(*slPotential, nr_grps);
157 for (i = 0; i < nr_grps; i++)
159 snew((*slField)[i], *nslices);
160 snew((*slCharge)[i], *nslices);
161 snew((*slPotential)[i], *nslices);
165 gpbc = gmx_rmpbc_init(&top->idef, ePBC, natoms, box);
167 /*********** Start processing trajectory ***********/
170 *slWidth = box[axis][axis]/(*nslices);
172 gmx_rmpbc(gpbc, natoms, box, x0);
174 /* calculate position of center of mass based on group 1 */
175 calc_xcm(x0, gnx[0], index[0], top->atoms.atom, xcm, FALSE);
178 for (n = 0; n < nr_grps; n++)
180 /* Check whether we actually have all positions of the requested index
181 * group in the trajectory file */
184 gmx_fatal(FARGS, "You selected a group with %d atoms, but only %d atoms\n"
185 "were found in the trajectory.\n", gnx[n], natoms);
187 for (i = 0; i < gnx[n]; i++) /* loop over all atoms in index file */
191 rvec_add(x0[index[n][i]], xcm, x0[index[n][i]]);
192 /* only distance from origin counts, not sign */
193 slice = norm(x0[index[n][i]])/(*slWidth);
195 /* this is a nice check for spherical groups but not for
196 all water in a cubic box since a lot will fall outside
198 if (slice > (*nslices))
200 fprintf(stderr,"Warning: slice = %d\n",slice);
203 (*slCharge)[n][slice] += top->atoms.atom[index[n][i]].q;
207 z = x0[index[n][i]][axis];
211 z += box[axis][axis];
213 if (z > box[axis][axis])
215 z -= box[axis][axis];
217 /* determine which slice atom is in */
218 slice = (z / (*slWidth));
219 (*slCharge)[n][slice] += top->atoms.atom[index[n][i]].q;
225 while (read_next_x(oenv, status, &t, natoms, x0, box));
227 gmx_rmpbc_done(gpbc);
229 /*********** done with status file **********/
232 /* slCharge now contains the total charge per slice, summed over all
233 frames. Now divide by nr_frames and integrate twice
239 fprintf(stderr, "\n\nRead %d frames from trajectory. Calculating potential"
240 "in spherical coordinates\n", nr_frames);
244 fprintf(stderr, "\n\nRead %d frames from trajectory. Calculating potential\n",
248 for (n = 0; n < nr_grps; n++)
250 for (i = 0; i < *nslices; i++)
254 /* charge per volume is now the summed charge, divided by the nr
255 of frames and by the volume of the slice it's in, 4pi r^2 dr
257 slVolume = 4*M_PI * sqr(i) * sqr(*slWidth) * *slWidth;
260 (*slCharge)[n][i] = 0;
264 (*slCharge)[n][i] = (*slCharge)[n][i] / (nr_frames * slVolume);
269 /* get charge per volume */
270 (*slCharge)[n][i] = (*slCharge)[n][i] * (*nslices) /
271 (nr_frames * box[axis][axis] * box[ax1][ax1] * box[ax2][ax2]);
274 /* Now we have charge densities */
277 if (bCorrect && !bSpherical)
279 for (n = 0; n < nr_grps; n++)
283 for (i = 0; i < *nslices; i++)
285 if (fabs((*slCharge)[n][i]) >= GMX_DOUBLE_MIN)
288 qsum += (*slCharge)[n][i];
292 for (i = 0; i < *nslices; i++)
294 if (fabs((*slCharge)[n][i]) >= GMX_DOUBLE_MIN)
296 (*slCharge)[n][i] -= qsum;
302 for (n = 0; n < nr_grps; n++)
304 /* integrate twice to get field and potential */
305 p_integrate((*slField)[n], (*slCharge)[n], *nslices, *slWidth);
309 if (bCorrect && !bSpherical)
311 for (n = 0; n < nr_grps; n++)
315 for (i = 0; i < *nslices; i++)
317 if (fabs((*slCharge)[n][i]) >= GMX_DOUBLE_MIN)
320 qsum += (*slField)[n][i];
324 for (i = 0; i < *nslices; i++)
326 if (fabs((*slCharge)[n][i]) >= GMX_DOUBLE_MIN)
328 (*slField)[n][i] -= qsum;
334 for (n = 0; n < nr_grps; n++)
336 p_integrate((*slPotential)[n], (*slField)[n], *nslices, *slWidth);
339 /* Now correct for eps0 and in spherical case for r*/
340 for (n = 0; n < nr_grps; n++)
342 for (i = 0; i < *nslices; i++)
346 (*slPotential)[n][i] = ELC * (*slPotential)[n][i] * -1.0E9 /
347 (EPS0 * i * (*slWidth));
348 (*slField)[n][i] = ELC * (*slField)[n][i] * 1E18 /
349 (EPS0 * i * (*slWidth));
353 (*slPotential)[n][i] = ELC * (*slPotential)[n][i] * -1.0E9 / EPS0;
354 (*slField)[n][i] = ELC * (*slField)[n][i] * 1E18 / EPS0;
359 sfree(x0); /* free memory used by coordinate array */
362 void plot_potential(double *potential[], double *charge[], double *field[],
363 const char *afile, const char *bfile, const char *cfile,
364 int nslices, int nr_grps, const char *grpname[], double slWidth,
365 const output_env_t oenv)
367 FILE *pot, /* xvgr file with potential */
368 *cha, /* xvgr file with charges */
369 *fie; /* xvgr files with fields */
370 char buf[256]; /* for xvgr title */
373 sprintf(buf, "Electrostatic Potential");
374 pot = xvgropen(afile, buf, "Box (nm)", "Potential (V)", oenv);
375 xvgr_legend(pot, nr_grps, grpname, oenv);
377 sprintf(buf, "Charge Distribution");
378 cha = xvgropen(bfile, buf, "Box (nm)", "Charge density (q/nm\\S3\\N)", oenv);
379 xvgr_legend(cha, nr_grps, grpname, oenv);
381 sprintf(buf, "Electric Field");
382 fie = xvgropen(cfile, buf, "Box (nm)", "Field (V/nm)", oenv);
383 xvgr_legend(fie, nr_grps, grpname, oenv);
385 for (slice = cb; slice < (nslices - ce); slice++)
387 fprintf(pot, "%20.16g ", slice * slWidth);
388 fprintf(cha, "%20.16g ", slice * slWidth);
389 fprintf(fie, "%20.16g ", slice * slWidth);
390 for (n = 0; n < nr_grps; n++)
392 fprintf(pot, " %20.16g", potential[n][slice]);
393 fprintf(fie, " %20.16g", field[n][slice]/1e9); /* convert to V/nm */
394 fprintf(cha, " %20.16g", charge[n][slice]);
406 int gmx_potential(int argc, char *argv[])
408 const char *desc[] = {
409 "[TT]g_potential[tt] computes the electrostatical potential across the box. The potential is",
410 "calculated by first summing the charges per slice and then integrating",
411 "twice of this charge distribution. Periodic boundaries are not taken",
412 "into account. Reference of potential is taken to be the left side of",
413 "the box. It is also possible to calculate the potential in spherical",
414 "coordinates as function of r by calculating a charge distribution in",
415 "spherical slices and twice integrating them. epsilon_r is taken as 1,",
416 "but 2 is more appropriate in many cases."
419 static int axis = 2; /* normal to memb. default z */
420 static const char *axtitle = "Z";
421 static int nslices = 10; /* nr of slices defined */
422 static int ngrps = 1;
423 static gmx_bool bSpherical = FALSE; /* default is bilayer types */
424 static real fudge_z = 0; /* translate coordinates */
425 static gmx_bool bCorrect = 0;
427 { "-d", FALSE, etSTR, {&axtitle},
428 "Take the normal on the membrane in direction X, Y or Z." },
429 { "-sl", FALSE, etINT, {&nslices},
430 "Calculate potential as function of boxlength, dividing the box"
431 " in this number of slices." },
432 { "-cb", FALSE, etINT, {&cb},
433 "Discard this number of first slices of box for integration" },
434 { "-ce", FALSE, etINT, {&ce},
435 "Discard this number of last slices of box for integration" },
436 { "-tz", FALSE, etREAL, {&fudge_z},
437 "Translate all coordinates by this distance in the direction of the box" },
438 { "-spherical", FALSE, etBOOL, {&bSpherical},
439 "Calculate spherical thingie" },
440 { "-ng", FALSE, etINT, {&ngrps},
441 "Number of groups to consider" },
442 { "-correct", FALSE, etBOOL, {&bCorrect},
443 "Assume net zero charge of groups to improve accuracy" }
445 const char *bugs[] = {
446 "Discarding slices for integration should not be necessary."
449 double **potential, /* potential per slice */
450 **charge, /* total charge per slice */
451 **field, /* field per slice */
452 slWidth; /* width of one slice */
453 char **grpname; /* groupnames */
454 int *ngx; /* sizes of groups */
455 t_topology *top; /* topology */
457 atom_id **index; /* indices for all groups */
458 t_filenm fnm[] = { /* files for g_order */
459 { efTRX, "-f", NULL, ffREAD }, /* trajectory file */
460 { efNDX, NULL, NULL, ffREAD }, /* index file */
461 { efTPX, NULL, NULL, ffREAD }, /* topology file */
462 { efXVG, "-o", "potential", ffWRITE }, /* xvgr output file */
463 { efXVG, "-oc", "charge", ffWRITE }, /* xvgr output file */
464 { efXVG, "-of", "field", ffWRITE }, /* xvgr output file */
467 #define NFILE asize(fnm)
469 parse_common_args(&argc, argv, PCA_CAN_VIEW | PCA_CAN_TIME | PCA_BE_NICE,
470 NFILE, fnm, asize(pa), pa, asize(desc), desc, asize(bugs), bugs,
474 axis = toupper(axtitle[0]) - 'X';
476 top = read_top(ftp2fn(efTPX, NFILE, fnm), &ePBC); /* read topology file */
478 snew(grpname, ngrps);
482 rd_index(ftp2fn(efNDX, NFILE, fnm), ngrps, ngx, index, grpname);
485 calc_potential(ftp2fn(efTRX, NFILE, fnm), index, ngx,
486 &potential, &charge, &field,
487 &nslices, top, ePBC, axis, ngrps, &slWidth, fudge_z,
488 bSpherical, bCorrect, oenv);
490 plot_potential(potential, charge, field, opt2fn("-o", NFILE, fnm),
491 opt2fn("-oc", NFILE, fnm), opt2fn("-of", NFILE, fnm),
492 nslices, ngrps, (const char**)grpname, slWidth, oenv);
494 do_view(oenv, opt2fn("-o", NFILE, fnm), NULL); /* view xvgr file */
495 do_view(oenv, opt2fn("-oc", NFILE, fnm), NULL); /* view xvgr file */
496 do_view(oenv, opt2fn("-of", NFILE, fnm), NULL); /* view xvgr file */