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