Make PBC type enumeration into PbcType enum class
[alexxy/gromacs.git] / src / gromacs / gmxana / gmx_dipoles.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 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.
11  *
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.
16  *
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.
21  *
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.
26  *
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.
34  *
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.
37  */
38 #include "gmxpre.h"
39
40 #include <cmath>
41 #include <cstring>
42
43 #include <algorithm>
44
45 #include "gromacs/commandline/pargs.h"
46 #include "gromacs/commandline/viewit.h"
47 #include "gromacs/correlationfunctions/autocorr.h"
48 #include "gromacs/fileio/confio.h"
49 #include "gromacs/fileio/enxio.h"
50 #include "gromacs/fileio/matio.h"
51 #include "gromacs/fileio/tpxio.h"
52 #include "gromacs/fileio/trxio.h"
53 #include "gromacs/fileio/xvgr.h"
54 #include "gromacs/gmxana/gmx_ana.h"
55 #include "gromacs/linearalgebra/nrjac.h"
56 #include "gromacs/listed_forces/bonded.h"
57 #include "gromacs/math/functions.h"
58 #include "gromacs/math/units.h"
59 #include "gromacs/math/vec.h"
60 #include "gromacs/math/vecdump.h"
61 #include "gromacs/mdtypes/md_enums.h"
62 #include "gromacs/pbcutil/pbc.h"
63 #include "gromacs/pbcutil/rmpbc.h"
64 #include "gromacs/statistics/statistics.h"
65 #include "gromacs/topology/index.h"
66 #include "gromacs/topology/topology.h"
67 #include "gromacs/trajectory/energyframe.h"
68 #include "gromacs/utility/arraysize.h"
69 #include "gromacs/utility/binaryinformation.h"
70 #include "gromacs/utility/cstringutil.h"
71 #include "gromacs/utility/exceptions.h"
72 #include "gromacs/utility/fatalerror.h"
73 #include "gromacs/utility/futil.h"
74 #include "gromacs/utility/gmxassert.h"
75 #include "gromacs/utility/smalloc.h"
76
77 #define e2d(x) ENM2DEBYE*(x)
78 #define EANG2CM (E_CHARGE * 1.0e-10)    /* e Angstrom to Coulomb meter */
79 #define CM2D (SPEED_OF_LIGHT * 1.0e+24) /* Coulomb meter to Debye */
80
81 typedef struct
82 {
83     int      nelem;
84     real     spacing, radius;
85     real*    elem;
86     int*     count;
87     gmx_bool bPhi;
88     int      nx, ny;
89     real**   cmap;
90 } t_gkrbin;
91
92 static t_gkrbin* mk_gkrbin(real radius, real rcmax, gmx_bool bPhi, int ndegrees)
93 {
94     t_gkrbin* gb;
95     char*     ptr;
96     int       i;
97
98     snew(gb, 1);
99
100     if ((ptr = getenv("GMX_DIPOLE_SPACING")) != nullptr)
101     {
102         double bw   = strtod(ptr, nullptr);
103         gb->spacing = bw;
104     }
105     else
106     {
107         gb->spacing = 0.01; /* nm */
108     }
109     gb->nelem = 1 + static_cast<int>(radius / gb->spacing);
110     if (rcmax == 0)
111     {
112         gb->nx = gb->nelem;
113     }
114     else
115     {
116         gb->nx = 1 + static_cast<int>(rcmax / gb->spacing);
117     }
118     gb->radius = radius;
119     snew(gb->elem, gb->nelem);
120     snew(gb->count, gb->nelem);
121
122     snew(gb->cmap, gb->nx);
123     gb->ny = std::max(2, ndegrees);
124     for (i = 0; (i < gb->nx); i++)
125     {
126         snew(gb->cmap[i], gb->ny);
127     }
128     gb->bPhi = bPhi;
129
130     return gb;
131 }
132
133 static void done_gkrbin(t_gkrbin** gb)
134 {
135     sfree((*gb)->elem);
136     sfree((*gb)->count);
137     sfree((*gb));
138     *gb = nullptr;
139 }
140
141 static void add2gkr(t_gkrbin* gb, real r, real cosa, real phi)
142 {
143     int  cy, index = gmx::roundToInt(r / gb->spacing);
144     real alpha;
145
146     if (index < gb->nelem)
147     {
148         gb->elem[index] += cosa;
149         gb->count[index]++;
150     }
151     if (index < gb->nx)
152     {
153         alpha = std::acos(cosa);
154         if (gb->bPhi)
155         {
156             cy = static_cast<int>((M_PI + phi) * gb->ny / (2 * M_PI));
157         }
158         else
159         {
160             cy = static_cast<int>((alpha * gb->ny) / M_PI); /*((1+cosa)*0.5*(gb->ny));*/
161         }
162         cy = std::min(gb->ny - 1, std::max(0, cy));
163         if (debug)
164         {
165             fprintf(debug, "CY: %10f  %5d\n", alpha, cy);
166         }
167         gb->cmap[index][cy] += 1;
168     }
169 }
170
171 static void rvec2sprvec(rvec dipcart, rvec dipsp)
172 {
173     clear_rvec(dipsp);
174     dipsp[0] = std::sqrt(dipcart[XX] * dipcart[XX] + dipcart[YY] * dipcart[YY]
175                          + dipcart[ZZ] * dipcart[ZZ]); /* R */
176     dipsp[1] = std::atan2(dipcart[YY], dipcart[XX]);   /* Theta */
177     dipsp[2] = std::atan2(std::sqrt(dipcart[XX] * dipcart[XX] + dipcart[YY] * dipcart[YY]),
178                           dipcart[ZZ]); /* Phi */
179 }
180
181
182 static void do_gkr(t_gkrbin*     gb,
183                    int           ncos,
184                    int*          ngrp,
185                    int*          molindex[],
186                    const int     mindex[],
187                    rvec          x[],
188                    rvec          mu[],
189                    PbcType       pbcType,
190                    const matrix  box,
191                    const t_atom* atom,
192                    const int*    nAtom)
193 {
194     static rvec* xcm[2] = { nullptr, nullptr };
195     int          gi, gj, j0, j1, i, j, k, n, grp0, grp1;
196     real         qtot, q, cosa, rr, phi;
197     rvec         dx;
198     t_pbc        pbc;
199
200     for (n = 0; (n < ncos); n++)
201     {
202         if (!xcm[n])
203         {
204             snew(xcm[n], ngrp[n]);
205         }
206         for (i = 0; (i < ngrp[n]); i++)
207         {
208             /* Calculate center of mass of molecule */
209             gi = molindex[n][i];
210             j0 = mindex[gi];
211
212             if (nAtom[n] > 0)
213             {
214                 copy_rvec(x[j0 + nAtom[n] - 1], xcm[n][i]);
215             }
216             else
217             {
218                 j1 = mindex[gi + 1];
219                 clear_rvec(xcm[n][i]);
220                 qtot = 0;
221                 for (j = j0; j < j1; j++)
222                 {
223                     q = std::abs(atom[j].q);
224                     qtot += q;
225                     for (k = 0; k < DIM; k++)
226                     {
227                         xcm[n][i][k] += q * x[j][k];
228                     }
229                 }
230                 svmul(1 / qtot, xcm[n][i], xcm[n][i]);
231             }
232         }
233     }
234     set_pbc(&pbc, pbcType, box);
235     grp0 = 0;
236     grp1 = ncos - 1;
237     for (i = 0; i < ngrp[grp0]; i++)
238     {
239         gi = molindex[grp0][i];
240         j0 = (ncos == 2) ? 0 : i + 1;
241         for (j = j0; j < ngrp[grp1]; j++)
242         {
243             gj = molindex[grp1][j];
244             if ((iprod(mu[gi], mu[gi]) > 0) && (iprod(mu[gj], mu[gj]) > 0))
245             {
246                 /* Compute distance between molecules including PBC */
247                 pbc_dx(&pbc, xcm[grp0][i], xcm[grp1][j], dx);
248                 rr = norm(dx);
249
250                 if (gb->bPhi)
251                 {
252                     rvec xi, xj, xk, xl;
253                     rvec r_ij, r_kj, r_kl, mm, nn;
254                     int  t1, t2, t3;
255
256                     copy_rvec(xcm[grp0][i], xj);
257                     copy_rvec(xcm[grp1][j], xk);
258                     rvec_add(xj, mu[gi], xi);
259                     rvec_add(xk, mu[gj], xl);
260                     phi  = dih_angle(xi, xj, xk, xl, &pbc, r_ij, r_kj, r_kl, mm, nn, /* out */
261                                     &t1, &t2, &t3);
262                     cosa = std::cos(phi);
263                 }
264                 else
265                 {
266                     cosa = cos_angle(mu[gi], mu[gj]);
267                     phi  = 0;
268                 }
269                 if (debug || std::isnan(cosa))
270                 {
271                     fprintf(debug ? debug : stderr,
272                             "mu[%d] = %5.2f %5.2f %5.2f |mi| = %5.2f, mu[%d] = %5.2f %5.2f %5.2f "
273                             "|mj| = %5.2f rr = %5.2f cosa = %5.2f\n",
274                             gi, mu[gi][XX], mu[gi][YY], mu[gi][ZZ], norm(mu[gi]), gj, mu[gj][XX],
275                             mu[gj][YY], mu[gj][ZZ], norm(mu[gj]), rr, cosa);
276                 }
277
278                 add2gkr(gb, rr, cosa, phi);
279             }
280         }
281     }
282 }
283
284 static real normalize_cmap(t_gkrbin* gb)
285 {
286     int    i, j;
287     real   hi;
288     double vol;
289
290     hi = 0;
291     for (i = 0; (i < gb->nx); i++)
292     {
293         vol = 4 * M_PI * gmx::square(gb->spacing * i) * gb->spacing;
294         for (j = 0; (j < gb->ny); j++)
295         {
296             gb->cmap[i][j] /= vol;
297             hi = std::max(hi, gb->cmap[i][j]);
298         }
299     }
300     if (hi <= 0)
301     {
302         gmx_fatal(FARGS, "No data in the cmap");
303     }
304     return hi;
305 }
306
307 static void print_cmap(const char* cmap, t_gkrbin* gb, int* nlevels)
308 {
309     FILE* out;
310     int   i, j;
311     real  hi;
312
313     real *xaxis, *yaxis;
314     t_rgb rlo = { 1, 1, 1 };
315     t_rgb rhi = { 0, 0, 0 };
316
317     hi = normalize_cmap(gb);
318     snew(xaxis, gb->nx + 1);
319     for (i = 0; (i < gb->nx + 1); i++)
320     {
321         xaxis[i] = i * gb->spacing;
322     }
323     snew(yaxis, gb->ny);
324     for (j = 0; (j < gb->ny); j++)
325     {
326         if (gb->bPhi)
327         {
328             yaxis[j] = (360.0 * j) / (gb->ny - 1.0) - 180;
329         }
330         else
331         {
332             yaxis[j] = (180.0 * j) / (gb->ny - 1.0);
333         }
334         /*2.0*j/(gb->ny-1.0)-1.0;*/
335     }
336     out = gmx_ffopen(cmap, "w");
337     write_xpm(out, 0, "Dipole Orientation Distribution", "Fraction", "r (nm)", gb->bPhi ? "Phi" : "Alpha",
338               gb->nx, gb->ny, xaxis, yaxis, gb->cmap, 0, hi, rlo, rhi, nlevels);
339     gmx_ffclose(out);
340     sfree(xaxis);
341     sfree(yaxis);
342 }
343
344 static void print_gkrbin(const char* fn, t_gkrbin* gb, int ngrp, int nframes, real volume, const gmx_output_env_t* oenv)
345 {
346     /* We compute Gk(r), gOO and hOO according to
347      * Nymand & Linse, JCP 112 (2000) pp 6386-6395.
348      * In this implementation the angle between dipoles is stored
349      * rather than their inner product. This allows to take polarizible
350      * models into account. The RDF is calculated as well, almost for free!
351      */
352     FILE*       fp;
353     const char* leg[] = { "G\\sk\\N(r)", "< cos >", "h\\sOO\\N", "g\\sOO\\N", "Energy" };
354     int         i, last;
355     real        x0, x1, ggg, Gkr, vol_s, rho, gOO, hOO, cosav, ener;
356     double      fac;
357
358     fp = xvgropen(fn, "Distance dependent Gk", "r (nm)", "G\\sk\\N(r)", oenv);
359     xvgr_legend(fp, asize(leg), leg, oenv);
360
361     Gkr = 1; /* Self-dipole inproduct = 1 */
362     rho = ngrp / volume;
363
364     if (debug)
365     {
366         fprintf(debug, "Number density is %g molecules / nm^3\n", rho);
367         fprintf(debug, "ngrp = %d, nframes = %d\n", ngrp, nframes);
368     }
369
370     last = gb->nelem - 1;
371     while (last > 1 && gb->elem[last - 1] == 0)
372     {
373         last--;
374     }
375
376     /* Divide by dipole squared, by number of frames, by number of origins.
377      * Multiply by 2 because we only take half the matrix of interactions
378      * into account.
379      */
380     fac = 2.0 / (ngrp * nframes);
381
382     x0 = 0;
383     for (i = 0; i < last; i++)
384     {
385         /* Centre of the coordinate in the spherical layer */
386         x1 = x0 + gb->spacing;
387
388         /* Volume of the layer */
389         vol_s = (4.0 / 3.0) * M_PI * (x1 * x1 * x1 - x0 * x0 * x0);
390
391         /* gOO */
392         gOO = gb->count[i] * fac / (rho * vol_s);
393
394         /* Dipole correlation hOO, normalized by the relative number density, like
395          * in a Radial distribution function.
396          */
397         ggg = gb->elem[i] * fac;
398         hOO = 3.0 * ggg / (rho * vol_s);
399         Gkr += ggg;
400         if (gb->count[i])
401         {
402             cosav = gb->elem[i] / gb->count[i];
403         }
404         else
405         {
406             cosav = 0;
407         }
408         ener = -0.5 * cosav * ONE_4PI_EPS0 / (x1 * x1 * x1);
409
410         fprintf(fp, "%10.5e %12.5e %12.5e %12.5e %12.5e  %12.5e\n", x1, Gkr, cosav, hOO, gOO, ener);
411
412         /* Swap x0 and x1 */
413         x0 = x1;
414     }
415     xvgrclose(fp);
416 }
417
418 static gmx_bool
419 read_mu_from_enx(ener_file_t fmu, int Vol, const ivec iMu, rvec mu, real* vol, real* t, int nre, t_enxframe* fr)
420 {
421     int      i;
422     gmx_bool bCont;
423     char     buf[22];
424
425     bCont = do_enx(fmu, fr);
426     if (fr->nre != nre)
427     {
428         fprintf(stderr,
429                 "Something strange: expected %d entries in energy file at step %s\n(time %g) but "
430                 "found %d entries\n",
431                 nre, gmx_step_str(fr->step, buf), fr->t, fr->nre);
432     }
433
434     if (bCont)
435     {
436         if (Vol != -1) /* we've got Volume in the energy file */
437         {
438             *vol = fr->ener[Vol].e;
439         }
440         for (i = 0; i < DIM; i++)
441         {
442             mu[i] = fr->ener[iMu[i]].e;
443         }
444         *t = fr->t;
445     }
446
447     return bCont;
448 }
449
450 static void neutralize_mols(int n, const int* index, const t_block* mols, t_atom* atom)
451 {
452     double mtot, qtot;
453     int    ncharged, m, a0, a1, a;
454
455     ncharged = 0;
456     for (m = 0; m < n; m++)
457     {
458         a0   = mols->index[index[m]];
459         a1   = mols->index[index[m] + 1];
460         mtot = 0;
461         qtot = 0;
462         for (a = a0; a < a1; a++)
463         {
464             mtot += atom[a].m;
465             qtot += atom[a].q;
466         }
467         /* This check is only for the count print */
468         if (std::abs(qtot) > 0.01)
469         {
470             ncharged++;
471         }
472         if (mtot > 0)
473         {
474             /* Remove the net charge at the center of mass */
475             for (a = a0; a < a1; a++)
476             {
477                 atom[a].q -= qtot * atom[a].m / mtot;
478             }
479         }
480     }
481
482     if (ncharged)
483     {
484         printf("There are %d charged molecules in the selection,\n"
485                "will subtract their charge at their center of mass\n",
486                ncharged);
487     }
488 }
489
490 static void mol_dip(int k0, int k1, rvec x[], const t_atom atom[], rvec mu)
491 {
492     int  k, m;
493     real q;
494
495     clear_rvec(mu);
496     for (k = k0; (k < k1); k++)
497     {
498         q = e2d(atom[k].q);
499         for (m = 0; (m < DIM); m++)
500         {
501             mu[m] += q * x[k][m];
502         }
503     }
504 }
505
506 static void mol_quad(int k0, int k1, rvec x[], const t_atom atom[], rvec quad)
507 {
508     int      i, k, m, n, niter;
509     real     q, r2, mass, masstot;
510     rvec     com; /* center of mass */
511     rvec     r;   /* distance of atoms to center of mass */
512     double** inten;
513     double   dd[DIM], **ev;
514
515     snew(inten, DIM);
516     snew(ev, DIM);
517     for (i = 0; (i < DIM); i++)
518     {
519         snew(inten[i], DIM);
520         snew(ev[i], DIM);
521         dd[i] = 0.0;
522     }
523
524     /* Compute center of mass */
525     clear_rvec(com);
526     masstot = 0.0;
527     for (k = k0; (k < k1); k++)
528     {
529         mass = atom[k].m;
530         masstot += mass;
531         for (i = 0; (i < DIM); i++)
532         {
533             com[i] += mass * x[k][i];
534         }
535     }
536     svmul((1.0 / masstot), com, com);
537
538     /* We want traceless quadrupole moments, so let us calculate the complete
539      * quadrupole moment tensor and diagonalize this tensor to get
540      * the individual components on the diagonal.
541      */
542
543 #define delta(a, b) (((a) == (b)) ? 1.0 : 0.0)
544
545     for (m = 0; (m < DIM); m++)
546     {
547         for (n = 0; (n < DIM); n++)
548         {
549             inten[m][n] = 0;
550         }
551     }
552     for (k = k0; (k < k1); k++) /* loop over atoms in a molecule */
553     {
554         q = (atom[k].q) * 100.0;
555         rvec_sub(x[k], com, r);
556         r2 = iprod(r, r);
557         for (m = 0; (m < DIM); m++)
558         {
559             for (n = 0; (n < DIM); n++)
560             {
561                 inten[m][n] += 0.5 * q * (3.0 * r[m] * r[n] - r2 * delta(m, n)) * EANG2CM * CM2D;
562             }
563         }
564     }
565     if (debug)
566     {
567         for (i = 0; (i < DIM); i++)
568         {
569             fprintf(debug, "Q[%d] = %8.3f  %8.3f  %8.3f\n", i, inten[i][XX], inten[i][YY], inten[i][ZZ]);
570         }
571     }
572
573     /* We've got the quadrupole tensor, now diagonalize the sucker */
574     jacobi(inten, 3, dd, ev, &niter);
575
576     if (debug)
577     {
578         for (i = 0; (i < DIM); i++)
579         {
580             fprintf(debug, "ev[%d] = %8.3f  %8.3f  %8.3f\n", i, ev[i][XX], ev[i][YY], ev[i][ZZ]);
581         }
582         for (i = 0; (i < DIM); i++)
583         {
584             fprintf(debug, "Q'[%d] = %8.3f  %8.3f  %8.3f\n", i, inten[i][XX], inten[i][YY], inten[i][ZZ]);
585         }
586     }
587     /* Sort the eigenvalues, for water we know that the order is as follows:
588      *
589      * Q_yy, Q_zz, Q_xx
590      *
591      * At the moment I have no idea how this will work out for other molecules...
592      */
593
594     if (dd[1] > dd[0])
595     {
596         std::swap(dd[0], dd[1]);
597     }
598     if (dd[2] > dd[1])
599     {
600         std::swap(dd[1], dd[2]);
601     }
602     if (dd[1] > dd[0])
603     {
604         std::swap(dd[0], dd[1]);
605     }
606
607     for (m = 0; (m < DIM); m++)
608     {
609         quad[0] = dd[2]; /* yy */
610         quad[1] = dd[0]; /* zz */
611         quad[2] = dd[1]; /* xx */
612     }
613
614     if (debug)
615     {
616         pr_rvec(debug, 0, "Quadrupole", quad, DIM, TRUE);
617     }
618
619     /* clean-up */
620     for (i = 0; (i < DIM); i++)
621     {
622         sfree(inten[i]);
623         sfree(ev[i]);
624     }
625     sfree(inten);
626     sfree(ev);
627 }
628
629 /*
630  * Calculates epsilon according to M. Neumann, Mol. Phys. 50, 841 (1983)
631  */
632 static real calc_eps(double M_diff, double volume, double epsRF, double temp)
633 {
634     double eps, A, teller, noemer;
635     double eps_0 = 8.854187817e-12; /* epsilon_0 in C^2 J^-1 m^-1 */
636     double fac   = 1.112650021e-59; /* converts Debye^2 to C^2 m^2 */
637
638     A = M_diff * fac / (3 * eps_0 * volume * NANO * NANO * NANO * BOLTZMANN * temp);
639
640     if (epsRF == 0.0)
641     {
642         teller = 1 + A;
643         noemer = 1;
644     }
645     else
646     {
647         teller = 1 + (A * 2 * epsRF / (2 * epsRF + 1));
648         noemer = 1 - (A / (2 * epsRF + 1));
649     }
650     eps = teller / noemer;
651
652     return eps;
653 }
654
655 static void update_slab_dipoles(int k0, int k1, rvec x[], rvec mu, int idim, int nslice, rvec slab_dipole[], matrix box)
656 {
657     int  k;
658     real xdim = 0;
659
660     for (k = k0; (k < k1); k++)
661     {
662         xdim += x[k][idim];
663     }
664     xdim /= (k1 - k0);
665     k = (static_cast<int>(xdim * nslice / box[idim][idim] + nslice)) % nslice;
666     rvec_inc(slab_dipole[k], mu);
667 }
668
669 static void dump_slab_dipoles(const char*             fn,
670                               int                     idim,
671                               int                     nslice,
672                               rvec                    slab_dipole[],
673                               matrix                  box,
674                               int                     nframes,
675                               const gmx_output_env_t* oenv)
676 {
677     FILE*       fp;
678     char        buf[STRLEN];
679     int         i;
680     real        mutot;
681     const char* leg_dim[4] = { "\\f{12}m\\f{4}\\sX\\N", "\\f{12}m\\f{4}\\sY\\N",
682                                "\\f{12}m\\f{4}\\sZ\\N", "\\f{12}m\\f{4}\\stot\\N" };
683
684     sprintf(buf, "Box-%c (nm)", 'X' + idim);
685     fp = xvgropen(fn, "Average dipole moment per slab", buf, "\\f{12}m\\f{4} (D)", oenv);
686     xvgr_legend(fp, DIM, leg_dim, oenv);
687     for (i = 0; (i < nslice); i++)
688     {
689         mutot = norm(slab_dipole[i]) / nframes;
690         fprintf(fp, "%10.3f  %10.3f  %10.3f  %10.3f  %10.3f\n",
691                 ((i + 0.5) * box[idim][idim]) / nslice, slab_dipole[i][XX] / nframes,
692                 slab_dipole[i][YY] / nframes, slab_dipole[i][ZZ] / nframes, mutot);
693     }
694     xvgrclose(fp);
695     do_view(oenv, fn, "-autoscale xy -nxy");
696 }
697
698 static void compute_avercos(int n, rvec dip[], real* dd, rvec axis, gmx_bool bPairs)
699 {
700     int    i, j, k;
701     double dc, d, ddc1, ddc2, ddc3;
702     rvec   xxx = { 1, 0, 0 };
703     rvec   yyy = { 0, 1, 0 };
704     rvec   zzz = { 0, 0, 1 };
705
706     d    = 0;
707     ddc1 = ddc2 = ddc3 = 0;
708     for (i = k = 0; (i < n); i++)
709     {
710         ddc1 += std::abs(cos_angle(dip[i], xxx));
711         ddc2 += std::abs(cos_angle(dip[i], yyy));
712         ddc3 += std::abs(cos_angle(dip[i], zzz));
713         if (bPairs)
714         {
715             for (j = i + 1; (j < n); j++, k++)
716             {
717                 dc = cos_angle(dip[i], dip[j]);
718                 d += std::abs(dc);
719             }
720         }
721     }
722     *dd      = d / k;
723     axis[XX] = ddc1 / n;
724     axis[YY] = ddc2 / n;
725     axis[ZZ] = ddc3 / n;
726 }
727
728 static void do_dip(const t_topology*       top,
729                    PbcType                 pbcType,
730                    real                    volume,
731                    const char*             fn,
732                    const char*             out_mtot,
733                    const char*             out_eps,
734                    const char*             out_aver,
735                    const char*             dipdist,
736                    const char*             cosaver,
737                    const char*             fndip3d,
738                    const char*             fnadip,
739                    gmx_bool                bPairs,
740                    const char*             corrtype,
741                    const char*             corf,
742                    gmx_bool                bGkr,
743                    const char*             gkrfn,
744                    gmx_bool                bPhi,
745                    int*                    nlevels,
746                    int                     ndegrees,
747                    int                     ncos,
748                    const char*             cmap,
749                    real                    rcmax,
750                    gmx_bool                bQuad,
751                    gmx_bool                bMU,
752                    const char*             mufn,
753                    int*                    gnx,
754                    int*                    molindex[],
755                    real                    mu_max,
756                    real                    mu_aver,
757                    real                    epsilonRF,
758                    real                    temp,
759                    int*                    gkatom,
760                    int                     skip,
761                    gmx_bool                bSlab,
762                    int                     nslices,
763                    const char*             axtitle,
764                    const char*             slabfn,
765                    const gmx_output_env_t* oenv)
766 {
767     const char* leg_mtot[] = { "M\\sx \\N", "M\\sy \\N", "M\\sz \\N", "|M\\stot \\N|" };
768 #define NLEGMTOT asize(leg_mtot)
769     const char* leg_eps[] = { "epsilon", "G\\sk", "g\\sk" };
770 #define NLEGEPS asize(leg_eps)
771     const char* leg_aver[] = { "< |M|\\S2\\N >", "< |M| >\\S2\\N", "< |M|\\S2\\N > - < |M| >\\S2\\N",
772                                "< |M| >\\S2\\N / < |M|\\S2\\N >" };
773 #define NLEGAVER asize(leg_aver)
774     const char* leg_cosaver[] = { "\\f{4}<|cos\\f{12}q\\f{4}\\sij\\N|>", "RMSD cos",
775                                   "\\f{4}<|cos\\f{12}q\\f{4}\\siX\\N|>",
776                                   "\\f{4}<|cos\\f{12}q\\f{4}\\siY\\N|>",
777                                   "\\f{4}<|cos\\f{12}q\\f{4}\\siZ\\N|>" };
778 #define NLEGCOSAVER asize(leg_cosaver)
779     const char* leg_adip[] = { "<mu>", "Std. Dev.", "Error" };
780 #define NLEGADIP asize(leg_adip)
781
782     FILE *         outdd, *outmtot, *outaver, *outeps, *caver = nullptr;
783     FILE *         dip3d = nullptr, *adip = nullptr;
784     rvec *         x, *dipole = nullptr, mu_t, quad, *dipsp = nullptr;
785     t_gkrbin*      gkrbin = nullptr;
786     gmx_enxnm_t*   enm    = nullptr;
787     t_enxframe*    fr;
788     int            nframes = 1000, nre, timecheck = 0, ncolour = 0;
789     ener_file_t    fmu = nullptr;
790     int            i, n, m, natom = 0, gnx_tot, teller, tel3;
791     t_trxstatus*   status;
792     int *          dipole_bin, ndipbin, ibin, iVol, idim = -1;
793     unsigned long  mode;
794     real           rcut = 0, t, t0, t1, dt, dd, rms_cos;
795     rvec           dipaxis;
796     matrix         box;
797     gmx_bool       bCorr, bTotal, bCont;
798     double         M_diff = 0, epsilon, invtel, vol_aver;
799     double         mu_ave, mu_mol, M2_ave = 0, M_ave2 = 0, M_av[DIM], M_av2[DIM];
800     double         M[3], M2[3], M4[3], Gk = 0, g_k = 0;
801     gmx_stats_t *  Qlsq, mulsq, muframelsq = nullptr;
802     ivec           iMu;
803     real**         muall        = nullptr;
804     rvec*          slab_dipoles = nullptr;
805     const t_atom*  atom         = nullptr;
806     const t_block* mols         = nullptr;
807     gmx_rmpbc_t    gpbc         = nullptr;
808
809     gnx_tot = gnx[0];
810     if (ncos == 2)
811     {
812         gnx_tot += gnx[1];
813     }
814     GMX_RELEASE_ASSERT(ncos == 1 || ncos == 2, "Invalid number of groups used with -ncos");
815
816     vol_aver = 0.0;
817
818     iVol = -1;
819
820     /* initialize to a negative value so we can see that it wasn't picked up */
821     iMu[XX] = iMu[YY] = iMu[ZZ] = -1;
822     if (bMU)
823     {
824         fmu = open_enx(mufn, "r");
825         do_enxnms(fmu, &nre, &enm);
826
827         /* Determine the indexes of the energy grps we need */
828         for (i = 0; (i < nre); i++)
829         {
830             if (std::strstr(enm[i].name, "Volume"))
831             {
832                 iVol = i;
833             }
834             else if (std::strstr(enm[i].name, "Mu-X"))
835             {
836                 iMu[XX] = i;
837             }
838             else if (std::strstr(enm[i].name, "Mu-Y"))
839             {
840                 iMu[YY] = i;
841             }
842             else if (std::strstr(enm[i].name, "Mu-Z"))
843             {
844                 iMu[ZZ] = i;
845             }
846         }
847         if (iMu[XX] < 0 || iMu[YY] < 0 || iMu[ZZ] < 0)
848         {
849             gmx_fatal(FARGS, "No index for Mu-X, Mu-Y or Mu-Z energy group.");
850         }
851     }
852     else
853     {
854         atom = top->atoms.atom;
855         mols = &(top->mols);
856     }
857     if ((iVol == -1) && bMU)
858     {
859         printf("Using Volume from topology: %g nm^3\n", volume);
860     }
861
862     /* Correlation stuff */
863     bCorr  = (corrtype[0] != 'n');
864     bTotal = (corrtype[0] == 't');
865     if (bCorr)
866     {
867         if (bTotal)
868         {
869             snew(muall, 1);
870             snew(muall[0], nframes * DIM);
871         }
872         else
873         {
874             snew(muall, gnx[0]);
875             for (i = 0; (i < gnx[0]); i++)
876             {
877                 snew(muall[i], nframes * DIM);
878             }
879         }
880     }
881
882     /* Allocate array which contains for every molecule in a frame the
883      * dipole moment.
884      */
885     if (!bMU)
886     {
887         snew(dipole, gnx_tot);
888     }
889
890     /* Statistics */
891     snew(Qlsq, DIM);
892     for (i = 0; (i < DIM); i++)
893     {
894         Qlsq[i] = gmx_stats_init();
895     }
896     mulsq = gmx_stats_init();
897
898     /* Open all the files */
899     outmtot = xvgropen(out_mtot, "Total dipole moment of the simulation box vs. time", "Time (ps)",
900                        "Total Dipole Moment (Debye)", oenv);
901     outeps  = xvgropen(out_eps, "Epsilon and Kirkwood factors", "Time (ps)", "", oenv);
902     outaver = xvgropen(out_aver, "Total dipole moment", "Time (ps)", "D", oenv);
903     if (bSlab)
904     {
905         idim = axtitle[0] - 'X';
906         if ((idim < 0) || (idim >= DIM))
907         {
908             idim = axtitle[0] - 'x';
909         }
910         if ((idim < 0) || (idim >= DIM))
911         {
912             bSlab = FALSE;
913         }
914         if (nslices < 2)
915         {
916             bSlab = FALSE;
917         }
918         fprintf(stderr, "axtitle = %s, nslices = %d, idim = %d\n", axtitle, nslices, idim);
919         if (bSlab)
920         {
921             snew(slab_dipoles, nslices);
922             fprintf(stderr, "Doing slab analysis\n");
923         }
924     }
925
926     if (fnadip)
927     {
928         adip = xvgropen(fnadip, "Average molecular dipole", "Dipole (D)", "", oenv);
929         xvgr_legend(adip, NLEGADIP, leg_adip, oenv);
930     }
931     if (cosaver)
932     {
933         caver = xvgropen(cosaver, bPairs ? "Average pair orientation" : "Average absolute dipole orientation",
934                          "Time (ps)", "", oenv);
935         xvgr_legend(caver, NLEGCOSAVER, bPairs ? leg_cosaver : &(leg_cosaver[1]), oenv);
936     }
937
938     if (fndip3d)
939     {
940         snew(dipsp, gnx_tot);
941
942         /* we need a dummy file for gnuplot */
943         dip3d = gmx_ffopen("dummy.dat", "w");
944         fprintf(dip3d, "%f %f %f", 0.0, 0.0, 0.0);
945         gmx_ffclose(dip3d);
946
947         dip3d = gmx_ffopen(fndip3d, "w");
948         try
949         {
950             gmx::BinaryInformationSettings settings;
951             settings.generatedByHeader(true);
952             settings.linePrefix("# ");
953             gmx::printBinaryInformation(dip3d, output_env_get_program_context(oenv), settings);
954         }
955         GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
956     }
957
958     /* Write legends to all the files */
959     xvgr_legend(outmtot, NLEGMTOT, leg_mtot, oenv);
960     xvgr_legend(outaver, NLEGAVER, leg_aver, oenv);
961
962     if (bMU && (mu_aver == -1))
963     {
964         xvgr_legend(outeps, NLEGEPS - 2, leg_eps, oenv);
965     }
966     else
967     {
968         xvgr_legend(outeps, NLEGEPS, leg_eps, oenv);
969     }
970
971     snew(fr, 1);
972     clear_rvec(mu_t);
973     teller = 0;
974     /* Read the first frame from energy or traj file */
975     if (bMU)
976     {
977         do
978         {
979             bCont = read_mu_from_enx(fmu, iVol, iMu, mu_t, &volume, &t, nre, fr);
980             if (bCont)
981             {
982                 timecheck = check_times(t);
983                 if (timecheck < 0)
984                 {
985                     teller++;
986                 }
987                 if ((teller % 10) == 0)
988                 {
989                     fprintf(stderr, "\r Skipping Frame %6d, time: %8.3f", teller, t);
990                     fflush(stderr);
991                 }
992             }
993             else
994             {
995                 printf("End of %s reached\n", mufn);
996                 break;
997             }
998         } while (bCont && (timecheck < 0));
999     }
1000     else
1001     {
1002         natom = read_first_x(oenv, &status, fn, &t, &x, box);
1003     }
1004
1005     /* Calculate spacing for dipole bin (simple histogram) */
1006     ndipbin = 1 + static_cast<int>(mu_max / 0.01);
1007     snew(dipole_bin, ndipbin);
1008     mu_ave = 0.0;
1009     for (m = 0; (m < DIM); m++)
1010     {
1011         M[m] = M2[m] = M4[m] = 0.0;
1012     }
1013
1014     if (bGkr)
1015     {
1016         /* Use 0.7 iso 0.5 to account for pressure scaling */
1017         /*  rcut   = 0.7*sqrt(max_cutoff2(box)); */
1018         rcut = 0.7
1019                * std::sqrt(gmx::square(box[XX][XX]) + gmx::square(box[YY][YY]) + gmx::square(box[ZZ][ZZ]));
1020
1021         gkrbin = mk_gkrbin(rcut, rcmax, bPhi, ndegrees);
1022     }
1023     gpbc = gmx_rmpbc_init(&top->idef, pbcType, natom);
1024
1025     /* Start while loop over frames */
1026     t0     = t;
1027     teller = 0;
1028     do
1029     {
1030         if (bCorr && (teller >= nframes))
1031         {
1032             nframes += 1000;
1033             if (bTotal)
1034             {
1035                 srenew(muall[0], nframes * DIM);
1036             }
1037             else
1038             {
1039                 for (i = 0; (i < gnx_tot); i++)
1040                 {
1041                     srenew(muall[i], nframes * DIM);
1042                 }
1043             }
1044         }
1045         t1 = t;
1046
1047         muframelsq = gmx_stats_init();
1048
1049         /* Initialise */
1050         for (m = 0; (m < DIM); m++)
1051         {
1052             M_av2[m] = 0;
1053         }
1054
1055         if (bMU)
1056         {
1057             /* Copy rvec into double precision local variable */
1058             for (m = 0; (m < DIM); m++)
1059             {
1060                 M_av[m] = mu_t[m];
1061             }
1062         }
1063         else
1064         {
1065             /* Initialise */
1066             for (m = 0; (m < DIM); m++)
1067             {
1068                 M_av[m] = 0;
1069             }
1070
1071             gmx_rmpbc(gpbc, natom, box, x);
1072
1073             /* Begin loop of all molecules in frame */
1074             for (n = 0; (n < ncos); n++)
1075             {
1076                 for (i = 0; (i < gnx[n]); i++)
1077                 {
1078                     int ind0, ind1;
1079
1080                     ind0 = mols->index[molindex[n][i]];
1081                     ind1 = mols->index[molindex[n][i] + 1];
1082
1083                     mol_dip(ind0, ind1, x, atom, dipole[i]);
1084                     gmx_stats_add_point(mulsq, 0, norm(dipole[i]), 0, 0);
1085                     gmx_stats_add_point(muframelsq, 0, norm(dipole[i]), 0, 0);
1086                     if (bSlab)
1087                     {
1088                         update_slab_dipoles(ind0, ind1, x, dipole[i], idim, nslices, slab_dipoles, box);
1089                     }
1090                     if (bQuad)
1091                     {
1092                         mol_quad(ind0, ind1, x, atom, quad);
1093                         for (m = 0; (m < DIM); m++)
1094                         {
1095                             gmx_stats_add_point(Qlsq[m], 0, quad[m], 0, 0);
1096                         }
1097                     }
1098                     if (bCorr && !bTotal)
1099                     {
1100                         tel3                = DIM * teller;
1101                         muall[i][tel3 + XX] = dipole[i][XX];
1102                         muall[i][tel3 + YY] = dipole[i][YY];
1103                         muall[i][tel3 + ZZ] = dipole[i][ZZ];
1104                     }
1105                     mu_mol = 0.0;
1106                     for (m = 0; (m < DIM); m++)
1107                     {
1108                         M_av[m] += dipole[i][m];               /* M per frame */
1109                         mu_mol += dipole[i][m] * dipole[i][m]; /* calc. mu for distribution */
1110                     }
1111                     mu_mol = std::sqrt(mu_mol);
1112
1113                     mu_ave += mu_mol; /* calc. the average mu */
1114
1115                     /* Update the dipole distribution */
1116                     ibin = gmx::roundToInt(ndipbin * mu_mol / mu_max);
1117                     if (ibin < ndipbin)
1118                     {
1119                         dipole_bin[ibin]++;
1120                     }
1121
1122                     if (fndip3d)
1123                     {
1124                         rvec2sprvec(dipole[i], dipsp[i]);
1125
1126                         if (dipsp[i][YY] > -M_PI && dipsp[i][YY] < -0.5 * M_PI)
1127                         {
1128                             if (dipsp[i][ZZ] < 0.5 * M_PI)
1129                             {
1130                                 ncolour = 1;
1131                             }
1132                             else
1133                             {
1134                                 ncolour = 2;
1135                             }
1136                         }
1137                         else if (dipsp[i][YY] > -0.5 * M_PI && dipsp[i][YY] < 0.0 * M_PI)
1138                         {
1139                             if (dipsp[i][ZZ] < 0.5 * M_PI)
1140                             {
1141                                 ncolour = 3;
1142                             }
1143                             else
1144                             {
1145                                 ncolour = 4;
1146                             }
1147                         }
1148                         else if (dipsp[i][YY] > 0.0 && dipsp[i][YY] < 0.5 * M_PI)
1149                         {
1150                             if (dipsp[i][ZZ] < 0.5 * M_PI)
1151                             {
1152                                 ncolour = 5;
1153                             }
1154                             else
1155                             {
1156                                 ncolour = 6;
1157                             }
1158                         }
1159                         else if (dipsp[i][YY] > 0.5 * M_PI && dipsp[i][YY] < M_PI)
1160                         {
1161                             if (dipsp[i][ZZ] < 0.5 * M_PI)
1162                             {
1163                                 ncolour = 7;
1164                             }
1165                             else
1166                             {
1167                                 ncolour = 8;
1168                             }
1169                         }
1170                         if (dip3d)
1171                         {
1172                             fprintf(dip3d,
1173                                     "set arrow %d from %f, %f, %f to %f, %f, %f lt %d  # %d %d\n",
1174                                     i + 1, x[ind0][XX], x[ind0][YY], x[ind0][ZZ],
1175                                     x[ind0][XX] + dipole[i][XX] / 25, x[ind0][YY] + dipole[i][YY] / 25,
1176                                     x[ind0][ZZ] + dipole[i][ZZ] / 25, ncolour, ind0, i);
1177                         }
1178                     }
1179                 } /* End loop of all molecules in frame */
1180
1181                 if (dip3d)
1182                 {
1183                     fprintf(dip3d, "set title \"t = %4.3f\"\n", t);
1184                     fprintf(dip3d, "set xrange [0.0:%4.2f]\n", box[XX][XX]);
1185                     fprintf(dip3d, "set yrange [0.0:%4.2f]\n", box[YY][YY]);
1186                     fprintf(dip3d, "set zrange [0.0:%4.2f]\n\n", box[ZZ][ZZ]);
1187                     fprintf(dip3d, "splot 'dummy.dat' using 1:2:3 w vec\n");
1188                     fprintf(dip3d, "pause -1 'Hit return to continue'\n");
1189                 }
1190             }
1191         }
1192         /* Compute square of total dipole */
1193         for (m = 0; (m < DIM); m++)
1194         {
1195             M_av2[m] = M_av[m] * M_av[m];
1196         }
1197
1198         if (cosaver)
1199         {
1200             compute_avercos(gnx_tot, dipole, &dd, dipaxis, bPairs);
1201             rms_cos = std::sqrt(gmx::square(dipaxis[XX] - 0.5) + gmx::square(dipaxis[YY] - 0.5)
1202                                 + gmx::square(dipaxis[ZZ] - 0.5));
1203             if (bPairs)
1204             {
1205                 fprintf(caver, "%10.3e  %10.3e  %10.3e  %10.3e  %10.3e  %10.3e\n", t, dd, rms_cos,
1206                         dipaxis[XX], dipaxis[YY], dipaxis[ZZ]);
1207             }
1208             else
1209             {
1210                 fprintf(caver, "%10.3e  %10.3e  %10.3e  %10.3e  %10.3e\n", t, rms_cos, dipaxis[XX],
1211                         dipaxis[YY], dipaxis[ZZ]);
1212             }
1213         }
1214
1215         if (bGkr)
1216         {
1217             do_gkr(gkrbin, ncos, gnx, molindex, mols->index, x, dipole, pbcType, box, atom, gkatom);
1218         }
1219
1220         if (bTotal)
1221         {
1222             tel3                = DIM * teller;
1223             muall[0][tel3 + XX] = M_av[XX];
1224             muall[0][tel3 + YY] = M_av[YY];
1225             muall[0][tel3 + ZZ] = M_av[ZZ];
1226         }
1227
1228         /* Write to file the total dipole moment of the box, and its components
1229          * for this frame.
1230          */
1231         if ((skip == 0) || ((teller % skip) == 0))
1232         {
1233             fprintf(outmtot, "%10g  %12.8e %12.8e %12.8e %12.8e\n", t, M_av[XX], M_av[YY], M_av[ZZ],
1234                     std::sqrt(M_av2[XX] + M_av2[YY] + M_av2[ZZ]));
1235         }
1236
1237         for (m = 0; (m < DIM); m++)
1238         {
1239             M[m] += M_av[m];
1240             M2[m] += M_av2[m];
1241             M4[m] += gmx::square(M_av2[m]);
1242         }
1243         /* Increment loop counter */
1244         teller++;
1245
1246         /* Calculate for output the running averages */
1247         invtel = 1.0 / teller;
1248         M2_ave = (M2[XX] + M2[YY] + M2[ZZ]) * invtel;
1249         M_ave2 = invtel * (invtel * (M[XX] * M[XX] + M[YY] * M[YY] + M[ZZ] * M[ZZ]));
1250         M_diff = M2_ave - M_ave2;
1251
1252         /* Compute volume from box in traj, else we use the one from above */
1253         if (!bMU)
1254         {
1255             volume = det(box);
1256         }
1257         vol_aver += volume;
1258
1259         epsilon = calc_eps(M_diff, (vol_aver / teller), epsilonRF, temp);
1260
1261         /* Calculate running average for dipole */
1262         if (mu_ave != 0)
1263         {
1264             mu_aver = (mu_ave / gnx_tot) * invtel;
1265         }
1266
1267         if ((skip == 0) || ((teller % skip) == 0))
1268         {
1269             /* Write to file < |M|^2 >, |< M >|^2. And the difference between
1270              * the two. Here M is sum mu_i. Further write the finite system
1271              * Kirkwood G factor and epsilon.
1272              */
1273             fprintf(outaver, "%10g  %10.3e %10.3e %10.3e %10.3e\n", t, M2_ave, M_ave2, M_diff,
1274                     M_ave2 / M2_ave);
1275
1276             if (fnadip)
1277             {
1278                 real aver;
1279                 gmx_stats_get_average(muframelsq, &aver);
1280                 fprintf(adip, "%10g %f \n", t, aver);
1281             }
1282             /*if (dipole)
1283                printf("%f %f\n", norm(dipole[0]), norm(dipole[1]));
1284              */
1285             if (!bMU || (mu_aver != -1))
1286             {
1287                 /* Finite system Kirkwood G-factor */
1288                 Gk = M_diff / (gnx_tot * mu_aver * mu_aver);
1289                 /* Infinite system Kirkwood G-factor */
1290                 if (epsilonRF == 0.0)
1291                 {
1292                     g_k = ((2 * epsilon + 1) * Gk / (3 * epsilon));
1293                 }
1294                 else
1295                 {
1296                     g_k = ((2 * epsilonRF + epsilon) * (2 * epsilon + 1) * Gk
1297                            / (3 * epsilon * (2 * epsilonRF + 1)));
1298                 }
1299
1300                 fprintf(outeps, "%10g  %10.3e %10.3e %10.3e\n", t, epsilon, Gk, g_k);
1301             }
1302             else
1303             {
1304                 fprintf(outeps, "%10g  %12.8e\n", t, epsilon);
1305             }
1306         }
1307         gmx_stats_free(muframelsq);
1308
1309         if (bMU)
1310         {
1311             bCont = read_mu_from_enx(fmu, iVol, iMu, mu_t, &volume, &t, nre, fr);
1312         }
1313         else
1314         {
1315             bCont = read_next_x(oenv, status, &t, x, box);
1316         }
1317         timecheck = check_times(t);
1318     } while (bCont && (timecheck == 0));
1319
1320     gmx_rmpbc_done(gpbc);
1321
1322     if (!bMU)
1323     {
1324         close_trx(status);
1325     }
1326
1327     xvgrclose(outmtot);
1328     xvgrclose(outaver);
1329     xvgrclose(outeps);
1330
1331     if (fnadip)
1332     {
1333         xvgrclose(adip);
1334     }
1335
1336     if (cosaver)
1337     {
1338         xvgrclose(caver);
1339     }
1340
1341     if (dip3d)
1342     {
1343         fprintf(dip3d, "set xrange [0.0:%4.2f]\n", box[XX][XX]);
1344         fprintf(dip3d, "set yrange [0.0:%4.2f]\n", box[YY][YY]);
1345         fprintf(dip3d, "set zrange [0.0:%4.2f]\n\n", box[ZZ][ZZ]);
1346         fprintf(dip3d, "splot 'dummy.dat' using 1:2:3 w vec\n");
1347         fprintf(dip3d, "pause -1 'Hit return to continue'\n");
1348         gmx_ffclose(dip3d);
1349     }
1350
1351     if (bSlab)
1352     {
1353         dump_slab_dipoles(slabfn, idim, nslices, slab_dipoles, box, teller, oenv);
1354         sfree(slab_dipoles);
1355     }
1356
1357     vol_aver /= teller;
1358     printf("Average volume over run is %g\n", vol_aver);
1359     if (bGkr)
1360     {
1361         print_gkrbin(gkrfn, gkrbin, gnx[0], teller, vol_aver, oenv);
1362         print_cmap(cmap, gkrbin, nlevels);
1363     }
1364     /* Autocorrelation function */
1365     if (bCorr)
1366     {
1367         if (teller < 2)
1368         {
1369             printf("Not enough frames for autocorrelation\n");
1370         }
1371         else
1372         {
1373             dt = (t1 - t0) / (teller - 1);
1374             printf("t0 %g, t %g, teller %d\n", t0, t, teller);
1375
1376             mode = eacVector;
1377
1378             if (bTotal)
1379             {
1380                 do_autocorr(corf, oenv, "Autocorrelation Function of Total Dipole", teller, 1,
1381                             muall, dt, mode, TRUE);
1382             }
1383             else
1384             {
1385                 do_autocorr(corf, oenv, "Dipole Autocorrelation Function", teller, gnx_tot, muall,
1386                             dt, mode, std::strcmp(corrtype, "molsep") != 0);
1387             }
1388         }
1389     }
1390     if (!bMU)
1391     {
1392         real aver, sigma, error;
1393
1394         gmx_stats_get_ase(mulsq, &aver, &sigma, &error);
1395         printf("\nDipole moment (Debye)\n");
1396         printf("---------------------\n");
1397         printf("Average  = %8.4f  Std. Dev. = %8.4f  Error = %8.4f\n", aver, sigma, error);
1398         if (bQuad)
1399         {
1400             rvec a, s, e;
1401             for (m = 0; (m < DIM); m++)
1402             {
1403                 gmx_stats_get_ase(mulsq, &(a[m]), &(s[m]), &(e[m]));
1404             }
1405
1406             printf("\nQuadrupole moment (Debye-Ang)\n");
1407             printf("-----------------------------\n");
1408             printf("Averages  = %8.4f  %8.4f  %8.4f\n", a[XX], a[YY], a[ZZ]);
1409             printf("Std. Dev. = %8.4f  %8.4f  %8.4f\n", s[XX], s[YY], s[ZZ]);
1410             printf("Error     = %8.4f  %8.4f  %8.4f\n", e[XX], e[YY], e[ZZ]);
1411         }
1412         printf("\n");
1413     }
1414     printf("The following averages for the complete trajectory have been calculated:\n\n");
1415     printf(" Total < M_x > = %g Debye\n", M[XX] / teller);
1416     printf(" Total < M_y > = %g Debye\n", M[YY] / teller);
1417     printf(" Total < M_z > = %g Debye\n\n", M[ZZ] / teller);
1418
1419     printf(" Total < M_x^2 > = %g Debye^2\n", M2[XX] / teller);
1420     printf(" Total < M_y^2 > = %g Debye^2\n", M2[YY] / teller);
1421     printf(" Total < M_z^2 > = %g Debye^2\n\n", M2[ZZ] / teller);
1422
1423     printf(" Total < |M|^2 > = %g Debye^2\n", M2_ave);
1424     printf(" Total |< M >|^2 = %g Debye^2\n\n", M_ave2);
1425
1426     printf(" < |M|^2 > - |< M >|^2 = %g Debye^2\n\n", M_diff);
1427
1428     if (!bMU || (mu_aver != -1))
1429     {
1430         printf("Finite system Kirkwood g factor G_k = %g\n", Gk);
1431         printf("Infinite system Kirkwood g factor g_k = %g\n\n", g_k);
1432     }
1433     printf("Epsilon = %g\n", epsilon);
1434
1435     if (!bMU)
1436     {
1437         /* Write to file the dipole moment distibution during the simulation.
1438          */
1439         outdd = xvgropen(dipdist, "Dipole Moment Distribution", "mu (Debye)", "", oenv);
1440         for (i = 0; (i < ndipbin); i++)
1441         {
1442             fprintf(outdd, "%10g  %10f\n", (i * mu_max) / ndipbin, static_cast<real>(dipole_bin[i]) / teller);
1443         }
1444         xvgrclose(outdd);
1445         sfree(dipole_bin);
1446     }
1447     if (bGkr)
1448     {
1449         done_gkrbin(&gkrbin);
1450     }
1451 }
1452
1453 static void dipole_atom2molindex(int* n, int* index, const t_block* mols)
1454 {
1455     int nmol, i, j, m;
1456
1457     nmol = 0;
1458     i    = 0;
1459     while (i < *n)
1460     {
1461         m = 0;
1462         while (m < mols->nr && index[i] != mols->index[m])
1463         {
1464             m++;
1465         }
1466         if (m == mols->nr)
1467         {
1468             gmx_fatal(FARGS, "index[%d]=%d does not correspond to the first atom of a molecule",
1469                       i + 1, index[i] + 1);
1470         }
1471         for (j = mols->index[m]; j < mols->index[m + 1]; j++)
1472         {
1473             if (i >= *n || index[i] != j)
1474             {
1475                 gmx_fatal(FARGS, "The index group is not a set of whole molecules");
1476             }
1477             i++;
1478         }
1479         /* Modify the index in place */
1480         index[nmol++] = m;
1481     }
1482     printf("There are %d molecules in the selection\n", nmol);
1483
1484     *n = nmol;
1485 }
1486 int gmx_dipoles(int argc, char* argv[])
1487 {
1488     const char* desc[] = {
1489         "[THISMODULE] computes the total dipole plus fluctuations of a simulation",
1490         "system. From this you can compute e.g. the dielectric constant for",
1491         "low-dielectric media.",
1492         "For molecules with a net charge, the net charge is subtracted at",
1493         "center of mass of the molecule.[PAR]",
1494         "The file [TT]Mtot.xvg[tt] contains the total dipole moment of a frame, the",
1495         "components as well as the norm of the vector.",
1496         "The file [TT]aver.xvg[tt] contains [CHEVRON][MAG][GRK]mu[grk][mag]^2[chevron] and ",
1497         "[MAG][CHEVRON][GRK]mu[grk][chevron][mag]^2 during the",
1498         "simulation.",
1499         "The file [TT]dipdist.xvg[tt] contains the distribution of dipole moments during",
1500         "the simulation",
1501         "The value of [TT]-mumax[tt] is used as the highest value in the distribution graph.[PAR]",
1502         "Furthermore, the dipole autocorrelation function will be computed when",
1503         "option [TT]-corr[tt] is used. The output file name is given with the [TT]-c[tt]",
1504         "option.",
1505         "The correlation functions can be averaged over all molecules",
1506         "([TT]mol[tt]), plotted per molecule separately ([TT]molsep[tt])",
1507         "or it can be computed over the total dipole moment of the simulation box",
1508         "([TT]total[tt]).[PAR]",
1509         "Option [TT]-g[tt] produces a plot of the distance dependent Kirkwood",
1510         "G-factor, as well as the average cosine of the angle between the dipoles",
1511         "as a function of the distance. The plot also includes gOO and hOO",
1512         "according to Nymand & Linse, J. Chem. Phys. 112 (2000) pp 6386-6395. In the same plot, ",
1513         "we also include the energy per scale computed by taking the inner product of",
1514         "the dipoles divided by the distance to the third power.[PAR]",
1515         "[PAR]",
1516         "EXAMPLES[PAR]",
1517         "[TT]gmx dipoles -corr mol -P 1 -o dip_sqr -mu 2.273 -mumax 5.0[tt][PAR]",
1518         "This will calculate the autocorrelation function of the molecular",
1519         "dipoles using a first order Legendre polynomial of the angle of the",
1520         "dipole vector and itself a time t later. For this calculation 1001",
1521         "frames will be used. Further, the dielectric constant will be calculated",
1522         "using an [TT]-epsilonRF[tt] of infinity (default), temperature of 300 K (default) and",
1523         "an average dipole moment of the molecule of 2.273 (SPC). For the",
1524         "distribution function a maximum of 5.0 will be used."
1525     };
1526     real              mu_max = 5, mu_aver = -1, rcmax = 0;
1527     real              epsilonRF = 0.0, temp = 300;
1528     gmx_bool          bPairs = TRUE, bPhi = FALSE, bQuad = FALSE;
1529     const char*       corrtype[] = { nullptr, "none", "mol", "molsep", "total", nullptr };
1530     const char*       axtitle    = "Z";
1531     int               nslices    = 10; /* nr of slices defined       */
1532     int               skip = 0, nFA = 0, nFB = 0, ncos = 1;
1533     int               nlevels = 20, ndegrees = 90;
1534     gmx_output_env_t* oenv;
1535     t_pargs           pa[] = {
1536         { "-mu", FALSE, etREAL, { &mu_aver }, "dipole of a single molecule (in Debye)" },
1537         { "-mumax", FALSE, etREAL, { &mu_max }, "max dipole in Debye (for histogram)" },
1538         { "-epsilonRF",
1539           FALSE,
1540           etREAL,
1541           { &epsilonRF },
1542           "[GRK]epsilon[grk] of the reaction field used during the simulation, needed for "
1543           "dielectric constant calculation. WARNING: 0.0 means infinity (default)" },
1544         { "-skip",
1545           FALSE,
1546           etINT,
1547           { &skip },
1548           "Skip steps in the output (but not in the computations)" },
1549         { "-temp",
1550           FALSE,
1551           etREAL,
1552           { &temp },
1553           "Average temperature of the simulation (needed for dielectric constant calculation)" },
1554         { "-corr", FALSE, etENUM, { corrtype }, "Correlation function to calculate" },
1555         { "-pairs",
1556           FALSE,
1557           etBOOL,
1558           { &bPairs },
1559           "Calculate [MAG][COS][GRK]theta[grk][cos][mag] between all pairs of molecules. May be "
1560           "slow" },
1561         { "-quad", FALSE, etBOOL, { &bQuad }, "Take quadrupole into account" },
1562         { "-ncos",
1563           FALSE,
1564           etINT,
1565           { &ncos },
1566           "Must be 1 or 2. Determines whether the [CHEVRON][COS][GRK]theta[grk][cos][chevron] is "
1567           "computed between all molecules in one group, or between molecules in two different "
1568           "groups. This turns on the [TT]-g[tt] flag." },
1569         { "-axis",
1570           FALSE,
1571           etSTR,
1572           { &axtitle },
1573           "Take the normal on the computational box in direction X, Y or Z." },
1574         { "-sl", FALSE, etINT, { &nslices }, "Divide the box into this number of slices." },
1575         { "-gkratom",
1576           FALSE,
1577           etINT,
1578           { &nFA },
1579           "Use the n-th atom of a molecule (starting from 1) to calculate the distance between "
1580           "molecules rather than the center of charge (when 0) in the calculation of distance "
1581           "dependent Kirkwood factors" },
1582         { "-gkratom2",
1583           FALSE,
1584           etINT,
1585           { &nFB },
1586           "Same as previous option in case ncos = 2, i.e. dipole interaction between two groups of "
1587           "molecules" },
1588         { "-rcmax",
1589           FALSE,
1590           etREAL,
1591           { &rcmax },
1592           "Maximum distance to use in the dipole orientation distribution (with ncos == 2). If "
1593           "zero, a criterion based on the box length will be used." },
1594         { "-phi",
1595           FALSE,
1596           etBOOL,
1597           { &bPhi },
1598           "Plot the 'torsion angle' defined as the rotation of the two dipole vectors around the "
1599           "distance vector between the two molecules in the [REF].xpm[ref] file from the "
1600           "[TT]-cmap[tt] option. By default the cosine of the angle between the dipoles is "
1601           "plotted." },
1602         { "-nlevels", FALSE, etINT, { &nlevels }, "Number of colors in the cmap output" },
1603         { "-ndegrees",
1604           FALSE,
1605           etINT,
1606           { &ndegrees },
1607           "Number of divisions on the [IT]y[it]-axis in the cmap output (for 180 degrees)" }
1608     };
1609     int*     gnx;
1610     int      nFF[2];
1611     int**    grpindex;
1612     char**   grpname = nullptr;
1613     gmx_bool bGkr, bMU, bSlab;
1614     t_filenm fnm[] = { { efEDR, "-en", nullptr, ffOPTRD },    { efTRX, "-f", nullptr, ffREAD },
1615                        { efTPR, nullptr, nullptr, ffREAD },   { efNDX, nullptr, nullptr, ffOPTRD },
1616                        { efXVG, "-o", "Mtot", ffWRITE },      { efXVG, "-eps", "epsilon", ffWRITE },
1617                        { efXVG, "-a", "aver", ffWRITE },      { efXVG, "-d", "dipdist", ffWRITE },
1618                        { efXVG, "-c", "dipcorr", ffOPTWR },   { efXVG, "-g", "gkr", ffOPTWR },
1619                        { efXVG, "-adip", "adip", ffOPTWR },   { efXVG, "-dip3d", "dip3d", ffOPTWR },
1620                        { efXVG, "-cos", "cosaver", ffOPTWR }, { efXPM, "-cmap", "cmap", ffOPTWR },
1621                        { efXVG, "-slab", "slab", ffOPTWR } };
1622 #define NFILE asize(fnm)
1623     int         npargs;
1624     t_pargs*    ppa;
1625     t_topology* top;
1626     PbcType     pbcType;
1627     int         k, natoms;
1628     matrix      box;
1629
1630     npargs = asize(pa);
1631     ppa    = add_acf_pargs(&npargs, pa);
1632     if (!parse_common_args(&argc, argv, PCA_CAN_TIME | PCA_CAN_VIEW, NFILE, fnm, npargs, ppa,
1633                            asize(desc), desc, 0, nullptr, &oenv))
1634     {
1635         sfree(ppa);
1636         return 0;
1637     }
1638
1639     printf("Using %g as mu_max and %g as the dipole moment.\n", mu_max, mu_aver);
1640     if (epsilonRF == 0.0)
1641     {
1642         printf("WARNING: EpsilonRF = 0.0, this really means EpsilonRF = infinity\n");
1643     }
1644
1645     bMU = opt2bSet("-en", NFILE, fnm);
1646     if (bMU)
1647     {
1648         gmx_fatal(FARGS,
1649                   "Due to new ways of treating molecules in GROMACS the total dipole in the energy "
1650                   "file may be incorrect, because molecules can be split over periodic boundary "
1651                   "conditions before computing the dipole. Please use your trajectory file.");
1652     }
1653     bGkr = opt2bSet("-g", NFILE, fnm);
1654     if (opt2parg_bSet("-ncos", asize(pa), pa))
1655     {
1656         if ((ncos != 1) && (ncos != 2))
1657         {
1658             gmx_fatal(FARGS, "ncos has to be either 1 or 2");
1659         }
1660         bGkr = TRUE;
1661     }
1662     bSlab = (opt2bSet("-slab", NFILE, fnm) || opt2parg_bSet("-sl", asize(pa), pa)
1663              || opt2parg_bSet("-axis", asize(pa), pa));
1664     if (bMU)
1665     {
1666         if (bQuad)
1667         {
1668             printf("WARNING: Can not determine quadrupoles from energy file\n");
1669             bQuad = FALSE;
1670         }
1671         if (bGkr)
1672         {
1673             printf("WARNING: Can not determine Gk(r) from energy file\n");
1674             bGkr = FALSE;
1675             ncos = 1;
1676         }
1677         if (mu_aver == -1)
1678         {
1679             printf("WARNING: Can not calculate Gk and gk, since you did\n"
1680                    "         not enter a valid dipole for the molecules\n");
1681         }
1682     }
1683
1684     snew(top, 1);
1685     pbcType = read_tpx_top(ftp2fn(efTPR, NFILE, fnm), nullptr, box, &natoms, nullptr, nullptr, top);
1686
1687     snew(gnx, ncos);
1688     snew(grpname, ncos);
1689     snew(grpindex, ncos);
1690     get_index(&top->atoms, ftp2fn_null(efNDX, NFILE, fnm), ncos, gnx, grpindex, grpname);
1691     for (k = 0; (k < ncos); k++)
1692     {
1693         dipole_atom2molindex(&gnx[k], grpindex[k], &(top->mols));
1694         neutralize_mols(gnx[k], grpindex[k], &(top->mols), top->atoms.atom);
1695     }
1696     nFF[0] = nFA;
1697     nFF[1] = nFB;
1698     do_dip(top, pbcType, det(box), ftp2fn(efTRX, NFILE, fnm), opt2fn("-o", NFILE, fnm),
1699            opt2fn("-eps", NFILE, fnm), opt2fn("-a", NFILE, fnm), opt2fn("-d", NFILE, fnm),
1700            opt2fn_null("-cos", NFILE, fnm), opt2fn_null("-dip3d", NFILE, fnm),
1701            opt2fn_null("-adip", NFILE, fnm), bPairs, corrtype[0], opt2fn("-c", NFILE, fnm), bGkr,
1702            opt2fn("-g", NFILE, fnm), bPhi, &nlevels, ndegrees, ncos, opt2fn("-cmap", NFILE, fnm),
1703            rcmax, bQuad, bMU, opt2fn("-en", NFILE, fnm), gnx, grpindex, mu_max, mu_aver, epsilonRF,
1704            temp, nFF, skip, bSlab, nslices, axtitle, opt2fn("-slab", NFILE, fnm), oenv);
1705
1706     do_view(oenv, opt2fn("-o", NFILE, fnm), "-autoscale xy -nxy");
1707     do_view(oenv, opt2fn("-eps", NFILE, fnm), "-autoscale xy -nxy");
1708     do_view(oenv, opt2fn("-a", NFILE, fnm), "-autoscale xy -nxy");
1709     do_view(oenv, opt2fn("-d", NFILE, fnm), "-autoscale xy");
1710     do_view(oenv, opt2fn("-c", NFILE, fnm), "-autoscale xy");
1711
1712     return 0;
1713 }