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 * check out http://www.gromacs.org for more information.
7 * Copyright (c) 2012,2013, by the GROMACS development team, led by
8 * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9 * others, as listed in the AUTHORS file in the top-level source
10 * directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
56 static void calc_com_pbc(int nrefat,t_topology *top,rvec x[],t_pbc *pbc,
57 atom_id index[],rvec xref,int ePBC,matrix box)
65 /* First simple calculation */
68 for(m=0; (m<nrefat); m++) {
70 mass = top->atoms.atom[ai].m;
71 for(j=0; (j<DIM); j++)
72 xref[j] += mass*x[ai][j];
75 svmul(1/mtot,xref,xref);
76 /* Now check if any atom is more than half the box from the COM */
77 if (ePBC != epbcNONE) {
81 for(m=0; (m<nrefat); m++) {
83 mass = top->atoms.atom[ai].m/mtot;
84 pbc_dx(pbc,x[ai],xref,dx);
85 rvec_add(xref,dx,xtest);
86 for(j=0; (j<DIM); j++)
87 if (fabs(xtest[j]-x[ai][j]) > tol) {
88 /* Here we have used the wrong image for contributing to the COM */
89 xref[j] += mass*(xtest[j]-x[ai][j]);
95 printf("COM: %8.3f %8.3f %8.3f iter = %d\n",xref[XX],xref[YY],xref[ZZ],iter);
101 void spol_atom2molindex(int *n,int *index,t_block *mols)
109 while(m < mols->nr && index[i] != mols->index[m])
112 gmx_fatal(FARGS,"index[%d]=%d does not correspond to the first atom of a molecule",i+1,index[i]+1);
113 for(j=mols->index[m]; j<mols->index[m+1]; j++) {
114 if (i >= *n || index[i] != j)
115 gmx_fatal(FARGS,"The index group is not a set of whole molecules");
118 /* Modify the index in place */
121 printf("There are %d molecules in the selection\n",nmol);
126 int gmx_spol(int argc,char *argv[])
133 int nrefat,natoms,nf,ntot;
135 rvec *xtop,*x,xref,trial,dx={0},dip,dir;
140 atom_id **index,*molindex;
142 real rmin2,rmax2,rcut,rcut2,rdx2=0,rtry2,qav,q,dip2,invbw;
143 int nbin,i,m,mol,a0,a1,a,d;
144 double sdip,sdip2,sinp,sdinp,nmol;
147 gmx_rmpbc_t gpbc=NULL;
150 const char *desc[] = {
151 "[TT]g_spol[tt] analyzes dipoles around a solute; it is especially useful",
152 "for polarizable water. A group of reference atoms, or a center",
153 "of mass reference (option [TT]-com[tt]) and a group of solvent",
154 "atoms is required. The program splits the group of solvent atoms",
155 "into molecules. For each solvent molecule the distance to the",
156 "closest atom in reference group or to the COM is determined.",
157 "A cumulative distribution of these distances is plotted.",
158 "For each distance between [TT]-rmin[tt] and [TT]-rmax[tt]",
159 "the inner product of the distance vector",
160 "and the dipole of the solvent molecule is determined.",
161 "For solvent molecules with net charge (ions), the net charge of the ion",
162 "is subtracted evenly from all atoms in the selection of each ion.",
163 "The average of these dipole components is printed.",
164 "The same is done for the polarization, where the average dipole is",
165 "subtracted from the instantaneous dipole. The magnitude of the average",
166 "dipole is set with the option [TT]-dip[tt], the direction is defined",
167 "by the vector from the first atom in the selected solvent group",
168 "to the midpoint between the second and the third atom."
172 static gmx_bool bCom = FALSE,bPBC = FALSE;
174 static real rmin=0.0,rmax=0.32,refdip=0,bw=0.01;
176 { "-com", FALSE, etBOOL, {&bCom},
177 "Use the center of mass as the reference postion" },
178 { "-refat", FALSE, etINT, {&srefat},
179 "The reference atom of the solvent molecule" },
180 { "-rmin", FALSE, etREAL, {&rmin}, "Maximum distance (nm)" },
181 { "-rmax", FALSE, etREAL, {&rmax}, "Maximum distance (nm)" },
182 { "-dip", FALSE, etREAL, {&refdip}, "The average dipole (D)" },
183 { "-bw", FALSE, etREAL, {&bw}, "The bin width" }
187 { efTRX, NULL, NULL, ffREAD },
188 { efTPX, NULL, NULL, ffREAD },
189 { efNDX, NULL, NULL, ffOPTRD },
190 { efXVG, NULL, "scdist.xvg", ffWRITE }
192 #define NFILE asize(fnm)
194 CopyRight(stderr,argv[0]);
195 parse_common_args(&argc,argv,PCA_CAN_TIME | PCA_CAN_VIEW | PCA_BE_NICE,
196 NFILE,fnm,asize(pa),pa,asize(desc),desc,0,NULL,&oenv);
200 read_tpx_top(ftp2fn(efTPX,NFILE,fnm),
201 ir,box,&natoms,NULL,NULL,NULL,top);
203 /* get index groups */
204 printf("Select a group of reference particles and a solvent group:\n");
208 get_index(&top->atoms,ftp2fn_null(efNDX,NFILE,fnm),2,isize,index,grpname);
218 spol_atom2molindex(&(isize[1]),index[1],&(top->mols));
221 /* initialize reading trajectory: */
222 natoms=read_first_x(oenv,&status,ftp2fn(efTRX,NFILE,fnm),&t,&x,box);
224 rcut = 0.99*sqrt(max_cutoff2(ir->ePBC,box));
229 nbin = (int)(rcut*invbw)+2;
242 molindex = top->mols.index;
243 atom = top->atoms.atom;
245 gpbc = gmx_rmpbc_init(&top->idef,ir->ePBC,natoms,box);
247 /* start analysis of trajectory */
249 /* make molecules whole again */
250 gmx_rmpbc(gpbc,natoms,box,x);
252 set_pbc(&pbc,ir->ePBC,box);
254 calc_com_pbc(nrefat,top,x,&pbc,index[0],xref,ir->ePBC,box);
256 for(m=0; m<isize[1]; m++) {
259 a1 = molindex[mol+1];
260 for(i=0; i<nrefgrp; i++) {
261 pbc_dx(&pbc,x[a0+srefat],bCom ? xref : x[index[0][i]],trial);
262 rtry2 = norm2(trial);
263 if (i==0 || rtry2 < rdx2) {
269 hist[(int)(sqrt(rdx2)*invbw)+1]++;
270 if (rdx2 >= rmin2 && rdx2 < rmax2) {
274 for(a=a0; a<a1; a++) {
278 for(a=a0; a<a1; a++) {
285 for(a=a0+1; a<a0+3; a++) {
287 dir[d] += 0.5*x[a][d];
291 svmul(ENM2DEBYE,dip,dip);
295 for(d=0; d<DIM; d++) {
296 sinp += dx[d]*dip[d];
297 sdinp += dx[d]*(dip[d] - refdip*dir[d]);
305 } while (read_next_x(oenv,status,&t,natoms,x,box));
307 gmx_rmpbc_done(gpbc);
313 fprintf(stderr,"Average number of molecules within %g nm is %.1f\n",
314 rmax,(real)ntot/(real)nf);
320 fprintf(stderr,"Average dipole: %f (D), std.dev. %f\n",
321 sdip,sqrt(sdip2-sqr(sdip)));
322 fprintf(stderr,"Average radial component of the dipole: %f (D)\n",
324 fprintf(stderr,"Average radial component of the polarization: %f (D)\n",
328 fp = xvgropen(opt2fn("-o",NFILE,fnm),
329 "Cumulative solvent distribution","r (nm)","molecules",oenv);
331 for(i=0; i<=nbin; i++) {
333 fprintf(fp,"%g %g\n",i*bw,nmol/nf);
337 do_view(oenv,opt2fn("-o",NFILE,fnm),NULL);