clang-tidy: readability-non-const-parameter (2/2)
[alexxy/gromacs.git] / src / gromacs / gmxana / gmx_spol.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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,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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 #include "gmxpre.h"
38
39 #include <cmath>
40
41 #include "gromacs/commandline/pargs.h"
42 #include "gromacs/commandline/viewit.h"
43 #include "gromacs/fileio/tpxio.h"
44 #include "gromacs/fileio/trxio.h"
45 #include "gromacs/fileio/xvgr.h"
46 #include "gromacs/gmxana/gmx_ana.h"
47 #include "gromacs/gmxana/gstat.h"
48 #include "gromacs/math/functions.h"
49 #include "gromacs/math/units.h"
50 #include "gromacs/math/vec.h"
51 #include "gromacs/mdtypes/inputrec.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/smalloc.h"
60
61 static void calc_com_pbc(int nrefat, const t_topology *top, rvec x[], t_pbc *pbc,
62                          const int index[], rvec xref, int ePBC)
63 {
64     const real tol = 1e-4;
65     gmx_bool   bChanged;
66     int        m, j, ai, iter;
67     real       mass, mtot;
68     rvec       dx, xtest;
69
70     /* First simple calculation */
71     clear_rvec(xref);
72     mtot = 0;
73     for (m = 0; (m < nrefat); m++)
74     {
75         ai   = index[m];
76         mass = top->atoms.atom[ai].m;
77         for (j = 0; (j < DIM); j++)
78         {
79             xref[j] += mass*x[ai][j];
80         }
81         mtot += mass;
82     }
83     svmul(1/mtot, xref, xref);
84     /* Now check if any atom is more than half the box from the COM */
85     if (ePBC != epbcNONE)
86     {
87         iter = 0;
88         do
89         {
90             bChanged = FALSE;
91             for (m = 0; (m < nrefat); m++)
92             {
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++)
98                 {
99                     if (std::abs(xtest[j]-x[ai][j]) > tol)
100                     {
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;
105                     }
106                 }
107             }
108             if (bChanged)
109             {
110                 printf("COM: %8.3f  %8.3f  %8.3f  iter = %d\n", xref[XX], xref[YY], xref[ZZ], iter);
111             }
112             iter++;
113         }
114         while (bChanged);
115     }
116 }
117
118 static void spol_atom2molindex(int *n, int *index, const t_block *mols)
119 {
120     int nmol, i, j, m;
121
122     nmol = 0;
123     i    = 0;
124     while (i < *n)
125     {
126         m = 0;
127         while (m < mols->nr && index[i] != mols->index[m])
128         {
129             m++;
130         }
131         if (m == mols->nr)
132         {
133             gmx_fatal(FARGS, "index[%d]=%d does not correspond to the first atom of a molecule", i+1, index[i]+1);
134         }
135         for (j = mols->index[m]; j < mols->index[m+1]; j++)
136         {
137             if (i >= *n || index[i] != j)
138             {
139                 gmx_fatal(FARGS, "The index group is not a set of whole molecules");
140             }
141             i++;
142         }
143         /* Modify the index in place */
144         index[nmol++] = m;
145     }
146     printf("There are %d molecules in the selection\n", nmol);
147
148     *n = nmol;
149 }
150
151 int gmx_spol(int argc, char *argv[])
152 {
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;
160
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;
171
172
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."
192     };
193
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},
200           "Use the center of mass as the reference position" },
201         { "-refat",  FALSE, etINT, {&srefat},
202           "The reference atom of the solvent molecule" },
203         { "-rmin",  FALSE, etREAL, {&rmin}, "Maximum distance (nm)" },
204         { "-rmax",  FALSE, etREAL, {&rmax}, "Maximum distance (nm)" },
205         { "-dip",   FALSE, etREAL, {&refdip}, "The average dipole (D)" },
206         { "-bw",    FALSE, etREAL, {&bw}, "The bin width" }
207     };
208
209     t_filenm          fnm[] = {
210         { efTRX, nullptr,  nullptr,  ffREAD },
211         { efTPR, nullptr,  nullptr,  ffREAD },
212         { efNDX, nullptr,  nullptr,  ffOPTRD },
213         { efXVG, nullptr,  "scdist",  ffWRITE }
214     };
215 #define NFILE asize(fnm)
216
217     if (!parse_common_args(&argc, argv, PCA_CAN_TIME | PCA_CAN_VIEW,
218                            NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, nullptr, &oenv))
219     {
220         return 0;
221     }
222
223     snew(top, 1);
224     // TODO: Only ePBC is used, not the full inputrec.
225     t_inputrec  irInstance;
226     t_inputrec *ir = &irInstance;
227     read_tpx_top(ftp2fn(efTPR, NFILE, fnm),
228                  ir, box, &natoms, nullptr, nullptr, top);
229
230     /* get index groups */
231     printf("Select a group of reference particles and a solvent group:\n");
232     snew(grpname, 2);
233     snew(index, 2);
234     snew(isize, 2);
235     get_index(&top->atoms, ftp2fn_null(efNDX, NFILE, fnm), 2, isize, index, grpname);
236
237     if (bCom)
238     {
239         nrefgrp = 1;
240         nrefat  = isize[0];
241     }
242     else
243     {
244         nrefgrp = isize[0];
245         nrefat  = 1;
246     }
247
248     spol_atom2molindex(&(isize[1]), index[1], &(top->mols));
249     srefat--;
250
251     /* initialize reading trajectory:                         */
252     natoms = read_first_x(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &t, &x, box);
253
254     rcut  = 0.99*std::sqrt(max_cutoff2(ir->ePBC, box));
255     if (rcut == 0)
256     {
257         rcut = 10*rmax;
258     }
259     rcut2 = gmx::square(rcut);
260     invbw = 1/bw;
261     nbin  = static_cast<int>(rcut*invbw)+2;
262     snew(hist, nbin);
263
264     rmin2 = gmx::square(rmin);
265     rmax2 = gmx::square(rmax);
266
267     nf    = 0;
268     ntot  = 0;
269     sdip  = 0;
270     sdip2 = 0;
271     sinp  = 0;
272     sdinp = 0;
273
274     molindex = top->mols.index;
275     atom     = top->atoms.atom;
276
277     gpbc = gmx_rmpbc_init(&top->idef, ir->ePBC, natoms);
278
279     /* start analysis of trajectory */
280     do
281     {
282         /* make molecules whole again */
283         gmx_rmpbc(gpbc, natoms, box, x);
284
285         set_pbc(&pbc, ir->ePBC, box);
286         if (bCom)
287         {
288             calc_com_pbc(nrefat, top, x, &pbc, index[0], xref, ir->ePBC);
289         }
290
291         for (m = 0; m < isize[1]; m++)
292         {
293             mol = index[1][m];
294             a0  = molindex[mol];
295             a1  = molindex[mol+1];
296             for (i = 0; i < nrefgrp; i++)
297             {
298                 pbc_dx(&pbc, x[a0+srefat], bCom ? xref : x[index[0][i]], trial);
299                 rtry2 = norm2(trial);
300                 if (i == 0 || rtry2 < rdx2)
301                 {
302                     copy_rvec(trial, dx);
303                     rdx2 = rtry2;
304                 }
305             }
306             if (rdx2 < rcut2)
307             {
308                 hist[static_cast<int>(std::sqrt(rdx2)*invbw)+1]++;
309             }
310             if (rdx2 >= rmin2 && rdx2 < rmax2)
311             {
312                 unitv(dx, dx);
313                 clear_rvec(dip);
314                 qav = 0;
315                 for (a = a0; a < a1; a++)
316                 {
317                     qav += atom[a].q;
318                 }
319                 qav /= (a1 - a0);
320                 for (a = a0; a < a1; a++)
321                 {
322                     q = atom[a].q - qav;
323                     for (d = 0; d < DIM; d++)
324                     {
325                         dip[d] += q*x[a][d];
326                     }
327                 }
328                 for (d = 0; d < DIM; d++)
329                 {
330                     dir[d] = -x[a0][d];
331                 }
332                 for (a = a0+1; a < a0+3; a++)
333                 {
334                     for (d = 0; d < DIM; d++)
335                     {
336                         dir[d] += 0.5*x[a][d];
337                     }
338                 }
339                 unitv(dir, dir);
340
341                 svmul(ENM2DEBYE, dip, dip);
342                 dip2   = norm2(dip);
343                 sdip  += std::sqrt(dip2);
344                 sdip2 += dip2;
345                 for (d = 0; d < DIM; d++)
346                 {
347                     sinp  += dx[d]*dip[d];
348                     sdinp += dx[d]*(dip[d] - refdip*dir[d]);
349                 }
350
351                 ntot++;
352             }
353         }
354         nf++;
355
356     }
357     while (read_next_x(oenv, status, &t, x, box));
358
359     gmx_rmpbc_done(gpbc);
360
361     /* clean up */
362     sfree(x);
363     close_trx(status);
364
365     fprintf(stderr, "Average number of molecules within %g nm is %.1f\n",
366             rmax, static_cast<real>(ntot)/nf);
367     if (ntot > 0)
368     {
369         sdip  /= ntot;
370         sdip2 /= ntot;
371         sinp  /= ntot;
372         sdinp /= ntot;
373         fprintf(stderr, "Average dipole:                               %f (D), std.dev. %f\n",
374                 sdip, std::sqrt(sdip2-gmx::square(sdip)));
375         fprintf(stderr, "Average radial component of the dipole:       %f (D)\n",
376                 sinp);
377         fprintf(stderr, "Average radial component of the polarization: %f (D)\n",
378                 sdinp);
379     }
380
381     fp = xvgropen(opt2fn("-o", NFILE, fnm),
382                   "Cumulative solvent distribution", "r (nm)", "molecules", oenv);
383     nmol = 0;
384     for (i = 0; i <= nbin; i++)
385     {
386         nmol += hist[i];
387         fprintf(fp, "%g %g\n", i*bw, nmol/nf);
388     }
389     xvgrclose(fp);
390
391     do_view(oenv, opt2fn("-o", NFILE, fnm), nullptr);
392
393     return 0;
394 }