Update instructions in containers.rst
[gromacs.git] / src / gromacs / gmxana / gmx_spol.cpp
blobd67045be907debdfe1c247d4f380ea00665b474d
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,2016,2017 by the GROMACS development team.
7 * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source 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.
38 #include "gmxpre.h"
40 #include <cmath>
42 #include "gromacs/commandline/pargs.h"
43 #include "gromacs/commandline/viewit.h"
44 #include "gromacs/fileio/tpxio.h"
45 #include "gromacs/fileio/trxio.h"
46 #include "gromacs/fileio/xvgr.h"
47 #include "gromacs/gmxana/gmx_ana.h"
48 #include "gromacs/gmxana/gstat.h"
49 #include "gromacs/math/functions.h"
50 #include "gromacs/math/units.h"
51 #include "gromacs/math/vec.h"
52 #include "gromacs/mdtypes/inputrec.h"
53 #include "gromacs/pbcutil/pbc.h"
54 #include "gromacs/pbcutil/rmpbc.h"
55 #include "gromacs/topology/index.h"
56 #include "gromacs/topology/topology.h"
57 #include "gromacs/utility/arraysize.h"
58 #include "gromacs/utility/cstringutil.h"
59 #include "gromacs/utility/fatalerror.h"
60 #include "gromacs/utility/smalloc.h"
62 static void calc_com_pbc(int nrefat, const t_topology* top, rvec x[], t_pbc* pbc, const int index[], rvec xref, PbcType pbcType)
64 const real tol = 1e-4;
65 gmx_bool bChanged;
66 int m, j, ai, iter;
67 real mass, mtot;
68 rvec dx, xtest;
70 /* First simple calculation */
71 clear_rvec(xref);
72 mtot = 0;
73 for (m = 0; (m < nrefat); m++)
75 ai = index[m];
76 mass = top->atoms.atom[ai].m;
77 for (j = 0; (j < DIM); j++)
79 xref[j] += mass * x[ai][j];
81 mtot += mass;
83 svmul(1 / mtot, xref, xref);
84 /* Now check if any atom is more than half the box from the COM */
85 if (pbcType != PbcType::No)
87 iter = 0;
90 bChanged = FALSE;
91 for (m = 0; (m < nrefat); m++)
93 ai = index[m];
94 mass = top->atoms.atom[ai].m / mtot;
95 pbc_dx(pbc, x[ai], xref, dx);
96 rvec_add(xref, dx, xtest);
97 for (j = 0; (j < DIM); j++)
99 if (std::abs(xtest[j] - x[ai][j]) > tol)
101 /* Here we have used the wrong image for contributing to the COM */
102 xref[j] += mass * (xtest[j] - x[ai][j]);
103 x[ai][j] = xtest[j];
104 bChanged = TRUE;
108 if (bChanged)
110 printf("COM: %8.3f %8.3f %8.3f iter = %d\n", xref[XX], xref[YY], xref[ZZ], iter);
112 iter++;
113 } while (bChanged);
117 static void spol_atom2molindex(int* n, int* index, const t_block* mols)
119 int nmol, i, j, m;
121 nmol = 0;
122 i = 0;
123 while (i < *n)
125 m = 0;
126 while (m < mols->nr && index[i] != mols->index[m])
128 m++;
130 if (m == mols->nr)
132 gmx_fatal(FARGS, "index[%d]=%d does not correspond to the first atom of a molecule",
133 i + 1, index[i] + 1);
135 for (j = mols->index[m]; j < mols->index[m + 1]; j++)
137 if (i >= *n || index[i] != j)
139 gmx_fatal(FARGS, "The index group is not a set of whole molecules");
141 i++;
143 /* Modify the index in place */
144 index[nmol++] = m;
146 printf("There are %d molecules in the selection\n", nmol);
148 *n = nmol;
151 int gmx_spol(int argc, char* argv[])
153 t_topology* top;
154 t_atom* atom;
155 t_trxstatus* status;
156 int nrefat, natoms, nf, ntot;
157 real t;
158 rvec * x, xref, trial, dx = { 0 }, dip, dir;
159 matrix box;
161 FILE* fp;
162 int * isize, nrefgrp;
163 int ** index, *molindex;
164 char** grpname;
165 real rmin2, rmax2, rcut, rcut2, rdx2 = 0, rtry2, qav, q, dip2, invbw;
166 int nbin, i, m, mol, a0, a1, a, d;
167 double sdip, sdip2, sinp, sdinp, nmol;
168 int* hist;
169 t_pbc pbc;
170 gmx_rmpbc_t gpbc = nullptr;
173 const char* desc[] = {
174 "[THISMODULE] analyzes dipoles around a solute; it is especially useful",
175 "for polarizable water. A group of reference atoms, or a center",
176 "of mass reference (option [TT]-com[tt]) and a group of solvent",
177 "atoms is required. The program splits the group of solvent atoms",
178 "into molecules. For each solvent molecule the distance to the",
179 "closest atom in reference group or to the COM is determined.",
180 "A cumulative distribution of these distances is plotted.",
181 "For each distance between [TT]-rmin[tt] and [TT]-rmax[tt]",
182 "the inner product of the distance vector",
183 "and the dipole of the solvent molecule is determined.",
184 "For solvent molecules with net charge (ions), the net charge of the ion",
185 "is subtracted evenly from all atoms in the selection of each ion.",
186 "The average of these dipole components is printed.",
187 "The same is done for the polarization, where the average dipole is",
188 "subtracted from the instantaneous dipole. The magnitude of the average",
189 "dipole is set with the option [TT]-dip[tt], the direction is defined",
190 "by the vector from the first atom in the selected solvent group",
191 "to the midpoint between the second and the third atom."
194 gmx_output_env_t* oenv;
195 static gmx_bool bCom = FALSE;
196 static int srefat = 1;
197 static real rmin = 0.0, rmax = 0.32, refdip = 0, bw = 0.01;
198 t_pargs pa[] = {
199 { "-com", FALSE, etBOOL, { &bCom }, "Use the center of mass as the reference position" },
200 { "-refat", FALSE, etINT, { &srefat }, "The reference atom of the solvent molecule" },
201 { "-rmin", FALSE, etREAL, { &rmin }, "Maximum distance (nm)" },
202 { "-rmax", FALSE, etREAL, { &rmax }, "Maximum distance (nm)" },
203 { "-dip", FALSE, etREAL, { &refdip }, "The average dipole (D)" },
204 { "-bw", FALSE, etREAL, { &bw }, "The bin width" }
207 t_filenm fnm[] = { { efTRX, nullptr, nullptr, ffREAD },
208 { efTPR, nullptr, nullptr, ffREAD },
209 { efNDX, nullptr, nullptr, ffOPTRD },
210 { efXVG, nullptr, "scdist", ffWRITE } };
211 #define NFILE asize(fnm)
213 if (!parse_common_args(&argc, argv, PCA_CAN_TIME | PCA_CAN_VIEW, NFILE, fnm, asize(pa), pa,
214 asize(desc), desc, 0, nullptr, &oenv))
216 return 0;
219 snew(top, 1);
220 // TODO: Only pbcType is used, not the full inputrec.
221 t_inputrec irInstance;
222 t_inputrec* ir = &irInstance;
223 read_tpx_top(ftp2fn(efTPR, NFILE, fnm), ir, box, &natoms, nullptr, nullptr, top);
225 /* get index groups */
226 printf("Select a group of reference particles and a solvent group:\n");
227 snew(grpname, 2);
228 snew(index, 2);
229 snew(isize, 2);
230 get_index(&top->atoms, ftp2fn_null(efNDX, NFILE, fnm), 2, isize, index, grpname);
232 if (bCom)
234 nrefgrp = 1;
235 nrefat = isize[0];
237 else
239 nrefgrp = isize[0];
240 nrefat = 1;
243 spol_atom2molindex(&(isize[1]), index[1], &(top->mols));
244 srefat--;
246 /* initialize reading trajectory: */
247 natoms = read_first_x(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &t, &x, box);
249 rcut = 0.99 * std::sqrt(max_cutoff2(ir->pbcType, box));
250 if (rcut == 0)
252 rcut = 10 * rmax;
254 rcut2 = gmx::square(rcut);
255 invbw = 1 / bw;
256 nbin = static_cast<int>(rcut * invbw) + 2;
257 snew(hist, nbin);
259 rmin2 = gmx::square(rmin);
260 rmax2 = gmx::square(rmax);
262 nf = 0;
263 ntot = 0;
264 sdip = 0;
265 sdip2 = 0;
266 sinp = 0;
267 sdinp = 0;
269 molindex = top->mols.index;
270 atom = top->atoms.atom;
272 gpbc = gmx_rmpbc_init(&top->idef, ir->pbcType, natoms);
274 /* start analysis of trajectory */
277 /* make molecules whole again */
278 gmx_rmpbc(gpbc, natoms, box, x);
280 set_pbc(&pbc, ir->pbcType, box);
281 if (bCom)
283 calc_com_pbc(nrefat, top, x, &pbc, index[0], xref, ir->pbcType);
286 for (m = 0; m < isize[1]; m++)
288 mol = index[1][m];
289 a0 = molindex[mol];
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);
298 rdx2 = rtry2;
301 if (rdx2 < rcut2)
303 hist[static_cast<int>(std::sqrt(rdx2) * invbw) + 1]++;
305 if (rdx2 >= rmin2 && rdx2 < rmax2)
307 unitv(dx, dx);
308 clear_rvec(dip);
309 qav = 0;
310 for (a = a0; a < a1; a++)
312 qav += atom[a].q;
314 qav /= (a1 - a0);
315 for (a = a0; a < a1; a++)
317 q = atom[a].q - qav;
318 for (d = 0; d < DIM; d++)
320 dip[d] += q * x[a][d];
323 for (d = 0; d < DIM; d++)
325 dir[d] = -x[a0][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];
334 unitv(dir, dir);
336 svmul(ENM2DEBYE, dip, dip);
337 dip2 = norm2(dip);
338 sdip += std::sqrt(dip2);
339 sdip2 += dip2;
340 for (d = 0; d < DIM; d++)
342 sinp += dx[d] * dip[d];
343 sdinp += dx[d] * (dip[d] - refdip * dir[d]);
346 ntot++;
349 nf++;
351 } while (read_next_x(oenv, status, &t, x, box));
353 gmx_rmpbc_done(gpbc);
355 /* clean up */
356 sfree(x);
357 close_trx(status);
359 fprintf(stderr, "Average number of molecules within %g nm is %.1f\n", rmax,
360 static_cast<real>(ntot) / nf);
361 if (ntot > 0)
363 sdip /= ntot;
364 sdip2 /= ntot;
365 sinp /= ntot;
366 sdinp /= ntot;
367 fprintf(stderr, "Average dipole: %f (D), std.dev. %f\n", sdip,
368 std::sqrt(sdip2 - gmx::square(sdip)));
369 fprintf(stderr, "Average radial component of the dipole: %f (D)\n", sinp);
370 fprintf(stderr, "Average radial component of the polarization: %f (D)\n", sdinp);
373 fp = xvgropen(opt2fn("-o", NFILE, fnm), "Cumulative solvent distribution", "r (nm)", "molecules", oenv);
374 nmol = 0;
375 for (i = 0; i <= nbin; i++)
377 nmol += hist[i];
378 fprintf(fp, "%g %g\n", i * bw, nmol / nf);
380 xvgrclose(fp);
382 do_view(oenv, opt2fn("-o", NFILE, fnm), nullptr);
384 return 0;