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.
39 #include "gromacs/legacyheaders/macros.h"
41 #include "gromacs/legacyheaders/viewit.h"
42 #include "gromacs/pbcutil/pbc.h"
43 #include "gromacs/topology/index.h"
44 #include "gromacs/fileio/tpxio.h"
45 #include "gromacs/fileio/trxio.h"
46 #include "gromacs/math/units.h"
49 #include "gromacs/commandline/pargs.h"
50 #include "gromacs/fileio/xvgr.h"
51 #include "gromacs/math/vec.h"
52 #include "gromacs/pbcutil/rmpbc.h"
53 #include "gromacs/utility/fatalerror.h"
54 #include "gromacs/utility/smalloc.h"
56 static void calc_com_pbc(int nrefat, t_topology *top, rvec x[], t_pbc *pbc,
57 atom_id index[], rvec xref, int ePBC)
59 const real tol = 1e-4;
65 /* First simple calculation */
68 for (m = 0; (m < nrefat); m++)
71 mass = top->atoms.atom[ai].m;
72 for (j = 0; (j < DIM); j++)
74 xref[j] += mass*x[ai][j];
78 svmul(1/mtot, xref, xref);
79 /* Now check if any atom is more than half the box from the COM */
86 for (m = 0; (m < nrefat); m++)
89 mass = top->atoms.atom[ai].m/mtot;
90 pbc_dx(pbc, x[ai], xref, dx);
91 rvec_add(xref, dx, xtest);
92 for (j = 0; (j < DIM); j++)
94 if (fabs(xtest[j]-x[ai][j]) > tol)
96 /* Here we have used the wrong image for contributing to the COM */
97 xref[j] += mass*(xtest[j]-x[ai][j]);
105 printf("COM: %8.3f %8.3f %8.3f iter = %d\n", xref[XX], xref[YY], xref[ZZ], iter);
113 void spol_atom2molindex(int *n, int *index, t_block *mols)
122 while (m < mols->nr && index[i] != mols->index[m])
128 gmx_fatal(FARGS, "index[%d]=%d does not correspond to the first atom of a molecule", i+1, index[i]+1);
130 for (j = mols->index[m]; j < mols->index[m+1]; j++)
132 if (i >= *n || index[i] != j)
134 gmx_fatal(FARGS, "The index group is not a set of whole molecules");
138 /* Modify the index in place */
141 printf("There are %d molecules in the selection\n", nmol);
146 int gmx_spol(int argc, char *argv[])
153 int nrefat, natoms, nf, ntot;
155 rvec *xtop, *x, xref, trial, dx = {0}, dip, dir;
160 atom_id **index, *molindex;
162 real rmin2, rmax2, rcut, rcut2, rdx2 = 0, rtry2, qav, q, dip2, invbw;
163 int nbin, i, m, mol, a0, a1, a, d;
164 double sdip, sdip2, sinp, sdinp, nmol;
167 gmx_rmpbc_t gpbc = NULL;
170 const char *desc[] = {
171 "[THISMODULE] analyzes dipoles around a solute; it is especially useful",
172 "for polarizable water. A group of reference atoms, or a center",
173 "of mass reference (option [TT]-com[tt]) and a group of solvent",
174 "atoms is required. The program splits the group of solvent atoms",
175 "into molecules. For each solvent molecule the distance to the",
176 "closest atom in reference group or to the COM is determined.",
177 "A cumulative distribution of these distances is plotted.",
178 "For each distance between [TT]-rmin[tt] and [TT]-rmax[tt]",
179 "the inner product of the distance vector",
180 "and the dipole of the solvent molecule is determined.",
181 "For solvent molecules with net charge (ions), the net charge of the ion",
182 "is subtracted evenly from all atoms in the selection of each ion.",
183 "The average of these dipole components is printed.",
184 "The same is done for the polarization, where the average dipole is",
185 "subtracted from the instantaneous dipole. The magnitude of the average",
186 "dipole is set with the option [TT]-dip[tt], the direction is defined",
187 "by the vector from the first atom in the selected solvent group",
188 "to the midpoint between the second and the third atom."
192 static gmx_bool bCom = FALSE, bPBC = FALSE;
193 static int srefat = 1;
194 static real rmin = 0.0, rmax = 0.32, refdip = 0, bw = 0.01;
196 { "-com", FALSE, etBOOL, {&bCom},
197 "Use the center of mass as the reference postion" },
198 { "-refat", FALSE, etINT, {&srefat},
199 "The reference atom of the solvent molecule" },
200 { "-rmin", FALSE, etREAL, {&rmin}, "Maximum distance (nm)" },
201 { "-rmax", FALSE, etREAL, {&rmax}, "Maximum distance (nm)" },
202 { "-dip", FALSE, etREAL, {&refdip}, "The average dipole (D)" },
203 { "-bw", FALSE, etREAL, {&bw}, "The bin width" }
207 { efTRX, NULL, NULL, ffREAD },
208 { efTPX, NULL, NULL, ffREAD },
209 { efNDX, NULL, NULL, ffOPTRD },
210 { efXVG, NULL, "scdist", ffWRITE }
212 #define NFILE asize(fnm)
214 if (!parse_common_args(&argc, argv, PCA_CAN_TIME | PCA_CAN_VIEW,
215 NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL, &oenv))
222 read_tpx_top(ftp2fn(efTPX, NFILE, fnm),
223 ir, box, &natoms, NULL, NULL, NULL, top);
225 /* get index groups */
226 printf("Select a group of reference particles and a solvent group:\n");
230 get_index(&top->atoms, ftp2fn_null(efNDX, NFILE, fnm), 2, isize, index, grpname);
243 spol_atom2molindex(&(isize[1]), index[1], &(top->mols));
246 /* initialize reading trajectory: */
247 natoms = read_first_x(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &t, &x, box);
249 rcut = 0.99*sqrt(max_cutoff2(ir->ePBC, box));
256 nbin = (int)(rcut*invbw)+2;
269 molindex = top->mols.index;
270 atom = top->atoms.atom;
272 gpbc = gmx_rmpbc_init(&top->idef, ir->ePBC, natoms);
274 /* start analysis of trajectory */
277 /* make molecules whole again */
278 gmx_rmpbc(gpbc, natoms, box, x);
280 set_pbc(&pbc, ir->ePBC, box);
283 calc_com_pbc(nrefat, top, x, &pbc, index[0], xref, ir->ePBC);
286 for (m = 0; m < isize[1]; m++)
290 a1 = molindex[mol+1];
291 for (i = 0; i < nrefgrp; i++)
293 pbc_dx(&pbc, x[a0+srefat], bCom ? xref : x[index[0][i]], trial);
294 rtry2 = norm2(trial);
295 if (i == 0 || rtry2 < rdx2)
297 copy_rvec(trial, dx);
303 hist[(int)(sqrt(rdx2)*invbw)+1]++;
305 if (rdx2 >= rmin2 && rdx2 < rmax2)
310 for (a = a0; a < a1; a++)
315 for (a = a0; a < a1; a++)
318 for (d = 0; d < DIM; d++)
323 for (d = 0; d < DIM; d++)
327 for (a = a0+1; a < a0+3; a++)
329 for (d = 0; d < DIM; d++)
331 dir[d] += 0.5*x[a][d];
336 svmul(ENM2DEBYE, dip, dip);
340 for (d = 0; d < DIM; d++)
342 sinp += dx[d]*dip[d];
343 sdinp += dx[d]*(dip[d] - refdip*dir[d]);
352 while (read_next_x(oenv, status, &t, x, box));
354 gmx_rmpbc_done(gpbc);
360 fprintf(stderr, "Average number of molecules within %g nm is %.1f\n",
361 rmax, (real)ntot/(real)nf);
368 fprintf(stderr, "Average dipole: %f (D), std.dev. %f\n",
369 sdip, sqrt(sdip2-sqr(sdip)));
370 fprintf(stderr, "Average radial component of the dipole: %f (D)\n",
372 fprintf(stderr, "Average radial component of the polarization: %f (D)\n",
376 fp = xvgropen(opt2fn("-o", NFILE, fnm),
377 "Cumulative solvent distribution", "r (nm)", "molecules", oenv);
379 for (i = 0; i <= nbin; i++)
382 fprintf(fp, "%g %g\n", i*bw, nmol/nf);
386 do_view(oenv, opt2fn("-o", NFILE, fnm), NULL);