Python detection consolidation.
[alexxy/gromacs.git] / src / gromacs / tools / pme_error.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2010,2011,2012,2013,2014,2015,2016,2017,2018,2019,2019, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35 #include "gmxpre.h"
36
37 #include "pme_error.h"
38
39 #include "config.h"
40
41 #include <cmath>
42
43 #include <algorithm>
44
45 #include "gromacs/commandline/pargs.h"
46 #include "gromacs/ewald/ewald_utils.h"
47 #include "gromacs/ewald/pme.h"
48 #include "gromacs/fft/calcgrid.h"
49 #include "gromacs/fileio/checkpoint.h"
50 #include "gromacs/fileio/tpxio.h"
51 #include "gromacs/gmxlib/network.h"
52 #include "gromacs/math/functions.h"
53 #include "gromacs/math/units.h"
54 #include "gromacs/math/vec.h"
55 #include "gromacs/mdlib/broadcaststructs.h"
56 #include "gromacs/mdtypes/commrec.h"
57 #include "gromacs/mdtypes/inputrec.h"
58 #include "gromacs/mdtypes/md_enums.h"
59 #include "gromacs/mdtypes/state.h"
60 #include "gromacs/random/threefry.h"
61 #include "gromacs/random/uniformintdistribution.h"
62 #include "gromacs/topology/mtop_util.h"
63 #include "gromacs/topology/topology.h"
64 #include "gromacs/utility/arraysize.h"
65 #include "gromacs/utility/fatalerror.h"
66 #include "gromacs/utility/pleasecite.h"
67 #include "gromacs/utility/smalloc.h"
68
69 /* #define TAKETIME */
70 /* #define DEBUG  */
71
72 /* Enum for situations that can occur during log file parsing */
73 enum {
74     eParselogOK,
75     eParselogNotFound,
76     eParselogNoPerfData,
77     eParselogTerm,
78     eParselogResetProblem,
79     eParselogNr
80 };
81
82
83 typedef struct
84 {
85     int64_t         orig_sim_steps;  /* Number of steps to be done in the real simulation  */
86     int             n_entries;       /* Number of entries in arrays                        */
87     real            volume;          /* The volume of the box                              */
88     matrix          recipbox;        /* The reciprocal box                                 */
89     int             natoms;          /* The number of atoms in the MD system               */
90     real           *fac;             /* The scaling factor                                 */
91     real           *rcoulomb;        /* The coulomb radii [0...nr_inputfiles]              */
92     real           *rvdw;            /* The vdW radii                                      */
93     int            *nkx, *nky, *nkz; /* Number of k vectors in each spatial dimension      */
94     real           *fourier_sp;      /* Fourierspacing                                     */
95     real           *ewald_rtol;      /* Real space tolerance for Ewald, determines         */
96                                      /* the real/reciprocal space relative weight          */
97     real           *ewald_beta;      /* Splitting parameter [1/nm]                         */
98     real            fracself;        /* fraction of particles for SI error                 */
99     real            q2all;           /* sum ( q ^2 )                                       */
100     real            q2allnr;         /* nr of charges                                      */
101     int            *pme_order;       /* Interpolation order for PME (bsplines)             */
102     char          **fn_out;          /* Name of the output tpr file                        */
103     real           *e_dir;           /* Direct space part of PME error with these settings */
104     real           *e_rec;           /* Reciprocal space part of PME error                 */
105     gmx_bool        bTUNE;           /* flag for tuning */
106 } t_inputinfo;
107
108
109 /* Returns TRUE when atom is charged */
110 static gmx_bool is_charge(real charge)
111 {
112     return charge*charge > GMX_REAL_EPS;
113 }
114
115
116 /* calculate charge density */
117 static void calc_q2all(const gmx_mtop_t *mtop,   /* molecular topology */
118                        real *q2all, real *q2allnr)
119 {
120     real            q2_all = 0;  /* Sum of squared charges */
121     int             nrq_mol;     /* Number of charges in a single molecule */
122     int             nrq_all;     /* Total number of charges in the MD system */
123     real            qi, q2_mol;
124
125 #ifdef DEBUG
126     fprintf(stderr, "\nCharge density:\n");
127 #endif
128     q2_all  = 0.0;                                        /* total q squared */
129     nrq_all = 0;                                          /* total number of charges in the system */
130     for (const gmx_molblock_t &molblock : mtop->molblock) /* Loop over molecule types */
131     {
132         q2_mol   = 0.0;                                   /* q squared value of this molecule */
133         nrq_mol  = 0;                                     /* number of charges this molecule carries */
134         const gmx_moltype_t &molecule = mtop->moltype[molblock.type];
135         for (int i = 0; i < molecule.atoms.nr; i++)
136         {
137             qi = molecule.atoms.atom[i].q;
138             /* Is this charge worth to be considered? */
139             if (is_charge(qi))
140             {
141                 q2_mol += qi*qi;
142                 nrq_mol++;
143             }
144         }
145         /* Multiply with the number of molecules present of this type and add */
146         q2_all  += q2_mol*molblock.nmol;
147         nrq_all += nrq_mol*molblock.nmol;
148 #ifdef DEBUG
149         fprintf(stderr, "Molecule %2d (%5d atoms) q2_mol=%10.3e nr.mol.charges=%5d (%6dx)  q2_all=%10.3e  tot.charges=%d\n",
150                 imol, molecule.atoms.nr, q2_mol, nrq_mol, molblock.nmol, q2_all, nrq_all);
151 #endif
152     }
153
154     *q2all   = q2_all;
155     *q2allnr = nrq_all;
156
157 }
158
159
160 /* Estimate the direct space part error of the SPME Ewald sum */
161 static real estimate_direct(
162         t_inputinfo *info
163         )
164 {
165     real e_dir     = 0; /* Error estimate */
166     real beta      = 0; /* Splitting parameter (1/nm) */
167     real r_coulomb = 0; /* Cut-off in direct space */
168
169
170     beta      = info->ewald_beta[0];
171     r_coulomb = info->rcoulomb[0];
172
173     e_dir  = 2.0 * info->q2all * gmx::invsqrt( info->q2allnr  *  r_coulomb * info->volume );
174     e_dir *= std::exp (-beta*beta*r_coulomb*r_coulomb);
175
176     return ONE_4PI_EPS0*e_dir;
177 }
178
179 #define SUMORDER 6
180
181 /* the following 4 functions determine polynomials required for the reciprocal error estimate */
182
183 static inline real eps_poly1(
184         real m,          /* grid coordinate in certain direction */
185         real K,          /* grid size in corresponding direction */
186         real n)          /* spline interpolation order of the SPME */
187 {
188     int  i;
189     real nom   = 0; /* nominator */
190     real denom = 0; /* denominator */
191     real tmp   = 0;
192
193     if (m == 0.0)
194     {
195         return 0.0;
196     }
197
198     for (i = -SUMORDER; i < 0; i++)
199     {
200         tmp  = m / K + i;
201         tmp *= 2.0*M_PI;
202         nom += std::pow( tmp, -n );
203     }
204
205     for (i = SUMORDER; i > 0; i--)
206     {
207         tmp  = m / K + i;
208         tmp *= 2.0*M_PI;
209         nom += std::pow( tmp, -n );
210     }
211
212     tmp   = m / K;
213     tmp  *= 2.0*M_PI;
214     denom = std::pow( tmp, -n )+nom;
215
216     return -nom/denom;
217
218 }
219
220 static inline real eps_poly2(
221         real m,          /* grid coordinate in certain direction */
222         real K,          /* grid size in corresponding direction */
223         real n)          /* spline interpolation order of the SPME */
224 {
225     int  i;
226     real nom   = 0; /* nominator */
227     real denom = 0; /* denominator */
228     real tmp   = 0;
229
230     if (m == 0.0)
231     {
232         return 0.0;
233     }
234
235     for (i = -SUMORDER; i < 0; i++)
236     {
237         tmp  = m / K + i;
238         tmp *= 2.0*M_PI;
239         nom += std::pow( tmp, -2*n );
240     }
241
242     for (i = SUMORDER; i > 0; i--)
243     {
244         tmp  = m / K + i;
245         tmp *= 2.0*M_PI;
246         nom += std::pow( tmp, -2*n );
247     }
248
249     for (i = -SUMORDER; i < SUMORDER+1; i++)
250     {
251         tmp    = m / K + i;
252         tmp   *= 2.0*M_PI;
253         denom += std::pow( tmp, -n );
254     }
255     tmp = eps_poly1(m, K, n);
256     return nom / denom / denom + tmp*tmp;
257
258 }
259
260 static inline real eps_poly3(
261         real m,          /* grid coordinate in certain direction */
262         real K,          /* grid size in corresponding direction */
263         real n)          /* spline interpolation order of the SPME */
264 {
265     int  i;
266     real nom   = 0; /* nominator */
267     real denom = 0; /* denominator */
268     real tmp   = 0;
269
270     if (m == 0.0)
271     {
272         return 0.0;
273     }
274
275     for (i = -SUMORDER; i < 0; i++)
276     {
277         tmp  = m / K + i;
278         tmp *= 2.0*M_PI;
279         nom += i * std::pow( tmp, -2*n );
280     }
281
282     for (i = SUMORDER; i > 0; i--)
283     {
284         tmp  = m / K + i;
285         tmp *= 2.0*M_PI;
286         nom += i * std::pow( tmp, -2*n );
287     }
288
289     for (i = -SUMORDER; i < SUMORDER+1; i++)
290     {
291         tmp    = m / K + i;
292         tmp   *= 2.0*M_PI;
293         denom += std::pow( tmp, -n );
294     }
295
296     return 2.0 * M_PI * nom / denom / denom;
297
298 }
299
300 static inline real eps_poly4(
301         real m,          /* grid coordinate in certain direction */
302         real K,          /* grid size in corresponding direction */
303         real n)          /* spline interpolation order of the SPME */
304 {
305     int  i;
306     real nom   = 0; /* nominator */
307     real denom = 0; /* denominator */
308     real tmp   = 0;
309
310     if (m == 0.0)
311     {
312         return 0.0;
313     }
314
315     for (i = -SUMORDER; i < 0; i++)
316     {
317         tmp  = m / K + i;
318         tmp *= 2.0*M_PI;
319         nom += i * i * std::pow( tmp, -2*n );
320     }
321
322     for (i = SUMORDER; i > 0; i--)
323     {
324         tmp  = m / K + i;
325         tmp *= 2.0*M_PI;
326         nom += i * i * std::pow( tmp, -2*n );
327     }
328
329     for (i = -SUMORDER; i < SUMORDER+1; i++)
330     {
331         tmp    = m / K + i;
332         tmp   *= 2.0*M_PI;
333         denom += std::pow( tmp, -n );
334     }
335
336     return 4.0 * M_PI * M_PI * nom / denom / denom;
337
338 }
339
340 static inline real eps_self(
341         real m,       /* grid coordinate in certain direction */
342         real K,       /* grid size in corresponding direction */
343         rvec rboxv,   /* reciprocal box vector */
344         real n,       /* spline interpolation order of the SPME */
345         rvec x)       /* coordinate of charge */
346 {
347     int  i;
348     real tmp    = 0; /* temporary variables for computations */
349     real tmp1   = 0; /* temporary variables for computations */
350     real tmp2   = 0; /* temporary variables for computations */
351     real rcoord = 0; /* coordinate in certain reciprocal space direction */
352     real nom    = 0; /* nominator */
353     real denom  = 0; /* denominator */
354
355
356     if (m == 0.0)
357     {
358         return 0.0;
359     }
360
361     rcoord = iprod(rboxv, x);
362
363
364     for (i = -SUMORDER; i < 0; i++)
365     {
366         tmp    = -std::sin(2.0 * M_PI * i * K * rcoord);
367         tmp1   = 2.0 * M_PI * m / K + 2.0 * M_PI * i;
368         tmp2   = std::pow(tmp1, -n);
369         nom   += tmp * tmp2 * i;
370         denom += tmp2;
371     }
372
373     for (i = SUMORDER; i > 0; i--)
374     {
375         tmp    = -std::sin(2.0 * M_PI * i * K * rcoord);
376         tmp1   = 2.0 * M_PI * m / K + 2.0 * M_PI * i;
377         tmp2   = std::pow(tmp1, -n);
378         nom   += tmp * tmp2 * i;
379         denom += tmp2;
380     }
381
382
383     tmp    = 2.0 * M_PI * m / K;
384     tmp1   = pow(tmp, -n);
385     denom += tmp1;
386
387     return 2.0 * M_PI * nom / denom * K;
388
389 }
390
391 #undef SUMORDER
392
393 /* The following routine is just a copy from pme.c */
394
395 static void calc_recipbox(matrix box, matrix recipbox)
396 {
397     /* Save some time by assuming upper right part is zero */
398
399     real tmp = 1.0/(box[XX][XX]*box[YY][YY]*box[ZZ][ZZ]);
400
401     recipbox[XX][XX] = box[YY][YY]*box[ZZ][ZZ]*tmp;
402     recipbox[XX][YY] = 0;
403     recipbox[XX][ZZ] = 0;
404     recipbox[YY][XX] = -box[YY][XX]*box[ZZ][ZZ]*tmp;
405     recipbox[YY][YY] = box[XX][XX]*box[ZZ][ZZ]*tmp;
406     recipbox[YY][ZZ] = 0;
407     recipbox[ZZ][XX] = (box[YY][XX]*box[ZZ][YY]-box[YY][YY]*box[ZZ][XX])*tmp;
408     recipbox[ZZ][YY] = -box[ZZ][YY]*box[XX][XX]*tmp;
409     recipbox[ZZ][ZZ] = box[XX][XX]*box[YY][YY]*tmp;
410 }
411
412
413 /* Estimate the reciprocal space part error of the SPME Ewald sum. */
414 static real estimate_reciprocal(
415         t_inputinfo             *info,
416         rvec                     x[], /* array of particles */
417         const real               q[], /* array of charges */
418         int                      nr,  /* number of charges = size of the charge array */
419         FILE  gmx_unused        *fp_out,
420         gmx_bool                 bVerbose,
421         int                      seed,     /* The seed for the random number generator */
422         int                     *nsamples, /* Return the number of samples used if Monte Carlo
423                                             * algorithm is used for self energy error estimate */
424         t_commrec               *cr)
425 {
426     real      e_rec   = 0; /* reciprocal error estimate */
427     real      e_rec1  = 0; /* Error estimate term 1*/
428     real      e_rec2  = 0; /* Error estimate term 2*/
429     real      e_rec3  = 0; /* Error estimate term 3 */
430     real      e_rec3x = 0; /* part of Error estimate term 3 in x */
431     real      e_rec3y = 0; /* part of Error estimate term 3 in y */
432     real      e_rec3z = 0; /* part of Error estimate term 3 in z */
433     int       i, ci;
434     int       nx, ny, nz;  /* grid coordinates */
435     real      q2_all = 0;  /* sum of squared charges */
436     rvec      gridpx;      /* reciprocal grid point in x direction*/
437     rvec      gridpxy;     /* reciprocal grid point in x and y direction*/
438     rvec      gridp;       /* complete reciprocal grid point in 3 directions*/
439     rvec      tmpvec;      /* template to create points from basis vectors */
440     rvec      tmpvec2;     /* template to create points from basis vectors */
441     real      coeff  = 0;  /* variable to compute coefficients of the error estimate */
442     real      coeff2 = 0;  /* variable to compute coefficients of the error estimate */
443     real      tmp    = 0;  /* variables to compute different factors from vectors */
444     real      tmp1   = 0;
445     real      tmp2   = 0;
446     gmx_bool  bFraction;
447
448     int      *numbers = nullptr;
449
450     /* Index variables for parallel work distribution */
451     int startglobal, stopglobal;
452     int startlocal, stoplocal;
453     int x_per_core;
454     int xtot;
455
456 #ifdef TAKETIME
457     double t0 = 0.0;
458     double t1 = 0.0;
459 #endif
460
461     if (seed == 0)
462     {
463         seed = static_cast<int>(gmx::makeRandomSeed());
464     }
465     fprintf(stderr, "Using random seed %d.\n", seed);
466
467     gmx::DefaultRandomEngine           rng(seed);
468     gmx::UniformIntDistribution<int>   dist(0, nr-1);
469
470     clear_rvec(gridpx);
471     clear_rvec(gridpxy);
472     clear_rvec(gridp);
473     clear_rvec(tmpvec);
474     clear_rvec(tmpvec2);
475
476     for (i = 0; i < nr; i++)
477     {
478         q2_all += q[i]*q[i];
479     }
480
481     /* Calculate indices for work distribution */
482     startglobal = -info->nkx[0]/2;
483     stopglobal  = info->nkx[0]/2;
484     xtot        = stopglobal*2+1;
485     if (PAR(cr))
486     {
487         x_per_core = static_cast<int>(std::ceil(static_cast<real>(xtot) / cr->nnodes));
488         startlocal = startglobal + x_per_core*cr->nodeid;
489         stoplocal  = startlocal + x_per_core -1;
490         if (stoplocal > stopglobal)
491         {
492             stoplocal = stopglobal;
493         }
494     }
495     else
496     {
497         startlocal = startglobal;
498         stoplocal  = stopglobal;
499         x_per_core = xtot;
500     }
501
502 #if GMX_LIB_MPI
503 #ifdef TAKETIME
504     if (MASTER(cr))
505     {
506         t0 = MPI_Wtime();
507     }
508 #endif
509 #endif
510
511     if (MASTER(cr))
512     {
513
514         fprintf(stderr, "Calculating reciprocal error part 1 ...");
515
516     }
517
518     for (nx = startlocal; nx <= stoplocal; nx++)
519     {
520         svmul(nx, info->recipbox[XX], gridpx);
521         for (ny = -info->nky[0]/2; ny < info->nky[0]/2+1; ny++)
522         {
523             svmul(ny, info->recipbox[YY], tmpvec);
524             rvec_add(gridpx, tmpvec, gridpxy);
525             for (nz = -info->nkz[0]/2; nz < info->nkz[0]/2+1; nz++)
526             {
527                 if (0 == nx &&  0 == ny &&  0 == nz)
528                 {
529                     continue;
530                 }
531                 svmul(nz, info->recipbox[ZZ], tmpvec);
532                 rvec_add(gridpxy, tmpvec, gridp);
533                 tmp    = norm2(gridp);
534                 coeff  = std::exp(-1.0 * M_PI * M_PI * tmp / info->ewald_beta[0] / info->ewald_beta[0] );
535                 coeff /= 2.0 * M_PI * info->volume * tmp;
536                 coeff2 = tmp;
537
538
539                 tmp  = eps_poly2(nx, info->nkx[0], info->pme_order[0]);
540                 tmp += eps_poly2(ny, info->nkx[0], info->pme_order[0]);
541                 tmp += eps_poly2(nz, info->nkx[0], info->pme_order[0]);
542
543                 tmp1 = eps_poly1(nx, info->nkx[0], info->pme_order[0]);
544                 tmp2 = eps_poly1(ny, info->nky[0], info->pme_order[0]);
545
546                 tmp += 2.0 * tmp1 * tmp2;
547
548                 tmp1 = eps_poly1(nz, info->nkz[0], info->pme_order[0]);
549                 tmp2 = eps_poly1(ny, info->nky[0], info->pme_order[0]);
550
551                 tmp += 2.0 * tmp1 * tmp2;
552
553                 tmp1 = eps_poly1(nz, info->nkz[0], info->pme_order[0]);
554                 tmp2 = eps_poly1(nx, info->nkx[0], info->pme_order[0]);
555
556                 tmp += 2.0 * tmp1 * tmp2;
557
558                 tmp1  = eps_poly1(nx, info->nkx[0], info->pme_order[0]);
559                 tmp1 += eps_poly1(ny, info->nky[0], info->pme_order[0]);
560                 tmp1 += eps_poly1(nz, info->nkz[0], info->pme_order[0]);
561
562                 tmp += tmp1 * tmp1;
563
564                 e_rec1 += 32.0 * M_PI * M_PI * coeff * coeff * coeff2 * tmp  * q2_all * q2_all / nr;
565
566                 tmp1  = eps_poly3(nx, info->nkx[0], info->pme_order[0]);
567                 tmp1 *= info->nkx[0];
568                 tmp2  = iprod(gridp, info->recipbox[XX]);
569
570                 tmp = tmp1*tmp2;
571
572                 tmp1  = eps_poly3(ny, info->nky[0], info->pme_order[0]);
573                 tmp1 *= info->nky[0];
574                 tmp2  = iprod(gridp, info->recipbox[YY]);
575
576                 tmp += tmp1*tmp2;
577
578                 tmp1  = eps_poly3(nz, info->nkz[0], info->pme_order[0]);
579                 tmp1 *= info->nkz[0];
580                 tmp2  = iprod(gridp, info->recipbox[ZZ]);
581
582                 tmp += tmp1*tmp2;
583
584                 tmp *= 4.0 * M_PI;
585
586                 tmp1  = eps_poly4(nx, info->nkx[0], info->pme_order[0]);
587                 tmp1 *= norm2(info->recipbox[XX]);
588                 tmp1 *= info->nkx[0] * info->nkx[0];
589
590                 tmp += tmp1;
591
592                 tmp1  = eps_poly4(ny, info->nky[0], info->pme_order[0]);
593                 tmp1 *= norm2(info->recipbox[YY]);
594                 tmp1 *= info->nky[0] * info->nky[0];
595
596                 tmp += tmp1;
597
598                 tmp1  = eps_poly4(nz, info->nkz[0], info->pme_order[0]);
599                 tmp1 *= norm2(info->recipbox[ZZ]);
600                 tmp1 *= info->nkz[0] * info->nkz[0];
601
602                 tmp += tmp1;
603
604                 e_rec2 += 4.0 * coeff * coeff * tmp * q2_all * q2_all / nr;
605
606             }
607         }
608         if (MASTER(cr))
609         {
610             fprintf(stderr, "\rCalculating reciprocal error part 1 ... %3.0f%%", 100.0*(nx-startlocal+1)/(x_per_core));
611             fflush(stderr);
612         }
613
614     }
615
616     if (MASTER(cr))
617     {
618         fprintf(stderr, "\n");
619     }
620
621     /* Use just a fraction of all charges to estimate the self energy error term? */
622     bFraction =  (info->fracself > 0.0) && (info->fracself < 1.0);
623
624     if (bFraction)
625     {
626         /* Here xtot is the number of samples taken for the Monte Carlo calculation
627          * of the average of term IV of equation 35 in Wang2010. Round up to a
628          * number of samples that is divisible by the number of nodes */
629         x_per_core  = static_cast<int>(std::ceil(info->fracself * nr / cr->nnodes));
630         xtot        = x_per_core * cr->nnodes;
631     }
632     else
633     {
634         /* In this case we use all nr particle positions */
635         xtot       = nr;
636         x_per_core = static_cast<int>(std::ceil(static_cast<real>(xtot) / cr->nnodes));
637     }
638
639     startlocal = x_per_core *  cr->nodeid;
640     stoplocal  = std::min(startlocal + x_per_core, xtot);  /* min needed if xtot == nr */
641
642     if (bFraction)
643     {
644         /* Make shure we get identical results in serial and parallel. Therefore,
645          * take the sample indices from a single, global random number array that
646          * is constructed on the master node and that only depends on the seed */
647         snew(numbers, xtot);
648         if (MASTER(cr))
649         {
650             for (i = 0; i < xtot; i++)
651             {
652                 numbers[i] = dist(rng); // [0,nr-1]
653             }
654         }
655         /* Broadcast the random number array to the other nodes */
656         if (PAR(cr))
657         {
658             nblock_bc(cr, xtot, numbers);
659         }
660
661         if (bVerbose && MASTER(cr))
662         {
663             fprintf(stdout, "Using %d sample%s to approximate the self interaction error term",
664                     xtot, xtot == 1 ? "" : "s");
665             if (PAR(cr))
666             {
667                 fprintf(stdout, " (%d sample%s per rank)", x_per_core, x_per_core == 1 ? "" : "s");
668             }
669             fprintf(stdout, ".\n");
670         }
671     }
672
673     /* Return the number of positions used for the Monte Carlo algorithm */
674     *nsamples = xtot;
675
676     for (i = startlocal; i < stoplocal; i++)
677     {
678         e_rec3x = 0;
679         e_rec3y = 0;
680         e_rec3z = 0;
681
682         if (bFraction)
683         {
684             /* Randomly pick a charge */
685             ci = numbers[i];
686         }
687         else
688         {
689             /* Use all charges */
690             ci = i;
691         }
692
693         /* for(nx=startlocal; nx<=stoplocal; nx++)*/
694         for (nx = -info->nkx[0]/2; nx < info->nkx[0]/2+1; nx++)
695         {
696             svmul(nx, info->recipbox[XX], gridpx);
697             for (ny = -info->nky[0]/2; ny < info->nky[0]/2+1; ny++)
698             {
699                 svmul(ny, info->recipbox[YY], tmpvec);
700                 rvec_add(gridpx, tmpvec, gridpxy);
701                 for (nz = -info->nkz[0]/2; nz < info->nkz[0]/2+1; nz++)
702                 {
703
704                     if (0 == nx && 0 == ny && 0 == nz)
705                     {
706                         continue;
707                     }
708
709                     svmul(nz, info->recipbox[ZZ], tmpvec);
710                     rvec_add(gridpxy, tmpvec, gridp);
711                     tmp      = norm2(gridp);
712                     coeff    = std::exp(-1.0 * M_PI * M_PI * tmp / info->ewald_beta[0] / info->ewald_beta[0] );
713                     coeff   /= tmp;
714                     e_rec3x += coeff*eps_self(nx, info->nkx[0], info->recipbox[XX], info->pme_order[0], x[ci]);
715                     e_rec3y += coeff*eps_self(ny, info->nky[0], info->recipbox[YY], info->pme_order[0], x[ci]);
716                     e_rec3z += coeff*eps_self(nz, info->nkz[0], info->recipbox[ZZ], info->pme_order[0], x[ci]);
717
718                 }
719             }
720         }
721
722         clear_rvec(tmpvec2);
723
724         svmul(e_rec3x, info->recipbox[XX], tmpvec);
725         rvec_inc(tmpvec2, tmpvec);
726         svmul(e_rec3y, info->recipbox[YY], tmpvec);
727         rvec_inc(tmpvec2, tmpvec);
728         svmul(e_rec3z, info->recipbox[ZZ], tmpvec);
729         rvec_inc(tmpvec2, tmpvec);
730
731         e_rec3 += q[ci]*q[ci]*q[ci]*q[ci]*norm2(tmpvec2) / ( xtot * M_PI * info->volume * M_PI * info->volume);
732         if (MASTER(cr))
733         {
734             fprintf(stderr, "\rCalculating reciprocal error part 2 ... %3.0f%%",
735                     100.0*(i+1)/stoplocal);
736             fflush(stderr);
737         }
738     }
739
740     if (MASTER(cr))
741     {
742         fprintf(stderr, "\n");
743     }
744
745 #if GMX_LIB_MPI
746 #ifdef TAKETIME
747     if (MASTER(cr))
748     {
749         t1 = MPI_Wtime() - t0;
750         fprintf(fp_out, "Recip. err. est. took   : %lf s\n", t1);
751     }
752 #endif
753 #endif
754
755 #ifdef DEBUG
756     if (PAR(cr))
757     {
758         fprintf(stderr, "Rank %3d: nx=[%3d...%3d]  e_rec3=%e\n",
759                 cr->nodeid, startlocal, stoplocal, e_rec3);
760     }
761 #endif
762
763     if (PAR(cr))
764     {
765         gmx_sum(1, &e_rec1, cr);
766         gmx_sum(1, &e_rec2, cr);
767         gmx_sum(1, &e_rec3, cr);
768     }
769
770     /* e_rec1*=8.0 * q2_all / info->volume / info->volume / nr ;
771        e_rec2*=  q2_all / M_PI / M_PI / info->volume / info->volume / nr ;
772        e_rec3/= M_PI * M_PI * info->volume * info->volume * nr ;
773      */
774     e_rec = std::sqrt(e_rec1+e_rec2+e_rec3);
775
776
777     return ONE_4PI_EPS0 * e_rec;
778 }
779
780
781 /* Allocate memory for the inputinfo struct: */
782 static void create_info(t_inputinfo *info)
783 {
784     snew(info->fac, info->n_entries);
785     snew(info->rcoulomb, info->n_entries);
786     snew(info->rvdw, info->n_entries);
787     snew(info->nkx, info->n_entries);
788     snew(info->nky, info->n_entries);
789     snew(info->nkz, info->n_entries);
790     snew(info->fourier_sp, info->n_entries);
791     snew(info->ewald_rtol, info->n_entries);
792     snew(info->ewald_beta, info->n_entries);
793     snew(info->pme_order, info->n_entries);
794     snew(info->fn_out, info->n_entries);
795     snew(info->e_dir, info->n_entries);
796     snew(info->e_rec, info->n_entries);
797 }
798
799
800 /* Allocate and fill an array with coordinates and charges,
801  * returns the number of charges found
802  */
803 static int prepare_x_q(real *q[], rvec *x[], const gmx_mtop_t *mtop, const rvec x_orig[], t_commrec *cr)
804 {
805     int                     nq; /* number of charged particles */
806
807
808     if (MASTER(cr))
809     {
810         snew(*q, mtop->natoms);
811         snew(*x, mtop->natoms);
812         nq = 0;
813
814         for (const AtomProxy atomP : AtomRange(*mtop))
815         {
816             const t_atom &local = atomP.atom();
817             int           i     = atomP.globalAtomNumber();
818             if (is_charge(local.q))
819             {
820                 (*q)[nq]     = local.q;
821                 (*x)[nq][XX] = x_orig[i][XX];
822                 (*x)[nq][YY] = x_orig[i][YY];
823                 (*x)[nq][ZZ] = x_orig[i][ZZ];
824                 nq++;
825             }
826         }
827         /* Give back some unneeded memory */
828         srenew(*q, nq);
829         srenew(*x, nq);
830     }
831     /* Broadcast x and q in the parallel case */
832     if (PAR(cr))
833     {
834         /* Transfer the number of charges */
835         block_bc(cr, nq);
836         snew_bc(cr, *x, nq);
837         snew_bc(cr, *q, nq);
838         nblock_bc(cr, nq, *x);
839         nblock_bc(cr, nq, *q);
840     }
841
842     return nq;
843 }
844
845
846
847 /* Read in the tpr file and save information we need later in info */
848 static void read_tpr_file(const char *fn_sim_tpr, t_inputinfo *info, t_state *state, gmx_mtop_t *mtop, t_inputrec *ir, real user_beta, real fracself)
849 {
850     read_tpx_state(fn_sim_tpr, ir, state, mtop);
851
852     /* The values of the original tpr input file are save in the first
853      * place [0] of the arrays */
854     info->orig_sim_steps = ir->nsteps;
855     info->pme_order[0]   = ir->pme_order;
856     info->rcoulomb[0]    = ir->rcoulomb;
857     info->rvdw[0]        = ir->rvdw;
858     info->nkx[0]         = ir->nkx;
859     info->nky[0]         = ir->nky;
860     info->nkz[0]         = ir->nkz;
861     info->ewald_rtol[0]  = ir->ewald_rtol;
862     info->fracself       = fracself;
863     if (user_beta > 0)
864     {
865         info->ewald_beta[0] = user_beta;
866     }
867     else
868     {
869         info->ewald_beta[0]  = calc_ewaldcoeff_q(info->rcoulomb[0], info->ewald_rtol[0]);
870     }
871
872     /* Check if PME was chosen */
873     if (EEL_PME(ir->coulombtype) == FALSE)
874     {
875         gmx_fatal(FARGS, "Can only do optimizations for simulations with PME");
876     }
877
878     /* Check if rcoulomb == rlist, which is necessary for PME */
879     if (!(ir->rcoulomb == ir->rlist))
880     {
881         gmx_fatal(FARGS, "PME requires rcoulomb (%f) to be equal to rlist (%f).", ir->rcoulomb, ir->rlist);
882     }
883 }
884
885
886 /* Transfer what we need for parallelizing the reciprocal error estimate */
887 static void bcast_info(t_inputinfo *info, const t_commrec *cr)
888 {
889     nblock_bc(cr, info->n_entries, info->nkx);
890     nblock_bc(cr, info->n_entries, info->nky);
891     nblock_bc(cr, info->n_entries, info->nkz);
892     nblock_bc(cr, info->n_entries, info->ewald_beta);
893     nblock_bc(cr, info->n_entries, info->pme_order);
894     nblock_bc(cr, info->n_entries, info->e_dir);
895     nblock_bc(cr, info->n_entries, info->e_rec);
896     block_bc(cr, info->volume);
897     block_bc(cr, info->recipbox);
898     block_bc(cr, info->natoms);
899     block_bc(cr, info->fracself);
900     block_bc(cr, info->bTUNE);
901     block_bc(cr, info->q2all);
902     block_bc(cr, info->q2allnr);
903 }
904
905
906 /* Estimate the error of the SPME Ewald sum. This estimate is based upon
907  * a) a homogeneous distribution of the charges
908  * b) a total charge of zero.
909  */
910 static void estimate_PME_error(t_inputinfo *info, const t_state *state,
911                                const gmx_mtop_t *mtop, FILE *fp_out, gmx_bool bVerbose, unsigned int seed,
912                                t_commrec *cr)
913 {
914     rvec *x     = nullptr; /* The coordinates */
915     real *q     = nullptr; /* The charges     */
916     real  edir  = 0.0;     /* real space error */
917     real  erec  = 0.0;     /* reciprocal space error */
918     real  derr  = 0.0;     /* difference of real and reciprocal space error */
919     real  derr0 = 0.0;     /* difference of real and reciprocal space error */
920     real  beta  = 0.0;     /* splitting parameter beta */
921     real  beta0 = 0.0;     /* splitting parameter beta */
922     int   ncharges;        /* The number of atoms with charges */
923     int   nsamples;        /* The number of samples used for the calculation of the
924                             * self-energy error term */
925     int   i = 0;
926
927     if (MASTER(cr))
928     {
929         fprintf(fp_out, "\n--- PME ERROR ESTIMATE ---\n");
930     }
931
932     /* Prepare an x and q array with only the charged atoms */
933     ncharges = prepare_x_q(&q, &x, mtop, state->x.rvec_array(), cr);
934     if (MASTER(cr))
935     {
936         calc_q2all(mtop, &(info->q2all), &(info->q2allnr));
937         info->ewald_rtol[0] = std::erfc(info->rcoulomb[0]*info->ewald_beta[0]);
938         /* Write some info to log file */
939         fprintf(fp_out, "Box volume              : %g nm^3\n", info->volume);
940         fprintf(fp_out, "Number of charged atoms : %d (total atoms %d)\n", ncharges, info->natoms);
941         fprintf(fp_out, "Coulomb radius          : %g nm\n", info->rcoulomb[0]);
942         fprintf(fp_out, "Ewald_rtol              : %g\n", info->ewald_rtol[0]);
943         fprintf(fp_out, "Ewald parameter beta    : %g\n", info->ewald_beta[0]);
944         fprintf(fp_out, "Interpolation order     : %d\n", info->pme_order[0]);
945         fprintf(fp_out, "Fourier grid (nx,ny,nz) : %d x %d x %d\n",
946                 info->nkx[0], info->nky[0], info->nkz[0]);
947         fflush(fp_out);
948
949     }
950
951     if (PAR(cr))
952     {
953         bcast_info(info, cr);
954     }
955
956
957     /* Calculate direct space error */
958     info->e_dir[0] = estimate_direct(info);
959
960     /* Calculate reciprocal space error */
961     info->e_rec[0] = estimate_reciprocal(info, x, q, ncharges, fp_out, bVerbose,
962                                          seed, &nsamples, cr);
963
964     if (PAR(cr))
965     {
966         bcast_info(info, cr);
967     }
968
969     if (MASTER(cr))
970     {
971         fprintf(fp_out, "Direct space error est. : %10.3e kJ/(mol*nm)\n", info->e_dir[0]);
972         fprintf(fp_out, "Reciprocal sp. err. est.: %10.3e kJ/(mol*nm)\n", info->e_rec[0]);
973         fprintf(fp_out, "Self-energy error term was estimated using %d samples\n", nsamples);
974         fflush(fp_out);
975         fprintf(stderr, "Direct space error est. : %10.3e kJ/(mol*nm)\n", info->e_dir[0]);
976         fprintf(stderr, "Reciprocal sp. err. est.: %10.3e kJ/(mol*nm)\n", info->e_rec[0]);
977     }
978
979     i = 0;
980
981     if (info->bTUNE)
982     {
983         if (MASTER(cr))
984         {
985             fprintf(stderr, "Starting tuning ...\n");
986         }
987         edir  = info->e_dir[0];
988         erec  = info->e_rec[0];
989         derr0 = edir-erec;
990         beta0 = info->ewald_beta[0];
991         if (derr > 0.0)
992         {
993             info->ewald_beta[0] += 0.1;
994         }
995         else
996         {
997             info->ewald_beta[0] -= 0.1;
998         }
999         info->e_dir[0] = estimate_direct(info);
1000         info->e_rec[0] = estimate_reciprocal(info, x, q, ncharges, fp_out, bVerbose,
1001                                              seed, &nsamples, cr);
1002
1003         if (PAR(cr))
1004         {
1005             bcast_info(info, cr);
1006         }
1007
1008
1009         edir = info->e_dir[0];
1010         erec = info->e_rec[0];
1011         derr = edir-erec;
1012         while (std::abs(derr/std::min(erec, edir)) > 1e-4)
1013         {
1014
1015             beta                = info->ewald_beta[0];
1016             beta               -= derr*(info->ewald_beta[0]-beta0)/(derr-derr0);
1017             beta0               = info->ewald_beta[0];
1018             info->ewald_beta[0] = beta;
1019             derr0               = derr;
1020
1021             info->e_dir[0] = estimate_direct(info);
1022             info->e_rec[0] = estimate_reciprocal(info, x, q, ncharges, fp_out, bVerbose,
1023                                                  seed, &nsamples, cr);
1024
1025             if (PAR(cr))
1026             {
1027                 bcast_info(info, cr);
1028             }
1029
1030             edir = info->e_dir[0];
1031             erec = info->e_rec[0];
1032             derr = edir-erec;
1033
1034             if (MASTER(cr))
1035             {
1036                 i++;
1037                 fprintf(stderr, "difference between real and rec. space error (step %d): %g\n", i, std::abs(derr));
1038                 fprintf(stderr, "old beta: %f\n", beta0);
1039                 fprintf(stderr, "new beta: %f\n", beta);
1040             }
1041         }
1042
1043         info->ewald_rtol[0] = std::erfc(info->rcoulomb[0]*info->ewald_beta[0]);
1044
1045         if (MASTER(cr))
1046         {
1047             /* Write some info to log file */
1048             fflush(fp_out);
1049             fprintf(fp_out, "=========  After tuning ========\n");
1050             fprintf(fp_out, "Direct space error est. : %10.3e kJ/(mol*nm)\n", info->e_dir[0]);
1051             fprintf(fp_out, "Reciprocal sp. err. est.: %10.3e kJ/(mol*nm)\n", info->e_rec[0]);
1052             fprintf(stderr, "Direct space error est. : %10.3e kJ/(mol*nm)\n", info->e_dir[0]);
1053             fprintf(stderr, "Reciprocal sp. err. est.: %10.3e kJ/(mol*nm)\n", info->e_rec[0]);
1054             fprintf(fp_out, "Ewald_rtol              : %g\n", info->ewald_rtol[0]);
1055             fprintf(fp_out, "Ewald parameter beta    : %g\n", info->ewald_beta[0]);
1056             fflush(fp_out);
1057
1058         }
1059
1060     }
1061
1062 }
1063
1064
1065 int gmx_pme_error(int argc, char *argv[])
1066 {
1067     const char     *desc[] = {
1068         "[THISMODULE] estimates the error of the electrostatic forces",
1069         "if using the sPME algorithm. The flag [TT]-tune[tt] will determine",
1070         "the splitting parameter such that the error is equally",
1071         "distributed over the real and reciprocal space part.",
1072         "The part of the error that stems from self interaction of the particles ",
1073         "is computationally demanding. However, a good a approximation is to",
1074         "just use a fraction of the particles for this term which can be",
1075         "indicated by the flag [TT]-self[tt].[PAR]",
1076     };
1077
1078     real            fs        = 0.0; /* 0 indicates: not set by the user */
1079     real            user_beta = -1.0;
1080     real            fracself  = 1.0;
1081     t_inputinfo     info;
1082     t_state         state;        /* The state from the tpr input file */
1083     gmx_mtop_t      mtop;         /* The topology from the tpr input file */
1084     FILE           *fp = nullptr;
1085     t_commrec      *cr;
1086     unsigned long   PCA_Flags;
1087     gmx_bool        bTUNE    = FALSE;
1088     gmx_bool        bVerbose = FALSE;
1089     int             seed     = 0;
1090
1091
1092     static t_filenm   fnm[] = {
1093         { efTPR, "-s",     nullptr,    ffREAD },
1094         { efOUT, "-o",    "error",  ffWRITE },
1095         { efTPR, "-so",   "tuned",  ffOPTWR }
1096     };
1097
1098     gmx_output_env_t *oenv = nullptr;
1099
1100     t_pargs           pa[] = {
1101         { "-beta",     FALSE, etREAL, {&user_beta},
1102           "If positive, overwrite ewald_beta from [REF].tpr[ref] file with this value" },
1103         { "-tune",     FALSE, etBOOL, {&bTUNE},
1104           "Tune the splitting parameter such that the error is equally distributed between real and reciprocal space" },
1105         { "-self",     FALSE, etREAL, {&fracself},
1106           "If between 0.0 and 1.0, determine self interaction error from just this fraction of the charged particles" },
1107         { "-seed",     FALSE, etINT,  {&seed},
1108           "Random number seed used for Monte Carlo algorithm when [TT]-self[tt] is set to a value between 0.0 and 1.0" },
1109         { "-v",        FALSE, etBOOL, {&bVerbose},
1110           "Be loud and noisy" }
1111     };
1112
1113
1114 #define NFILE asize(fnm)
1115
1116     cr         = init_commrec();
1117     PCA_Flags  = PCA_NOEXIT_ON_ARGS;
1118
1119     if (!parse_common_args(&argc, argv, PCA_Flags,
1120                            NFILE, fnm, asize(pa), pa, asize(desc), desc,
1121                            0, nullptr, &oenv))
1122     {
1123         sfree(cr);
1124         return 0;
1125     }
1126
1127     if (!bTUNE)
1128     {
1129         bTUNE = opt2bSet("-so", NFILE, fnm);
1130     }
1131
1132     info.n_entries = 1;
1133
1134     /* Allocate memory for the inputinfo struct: */
1135     create_info(&info);
1136     info.fourier_sp[0] = fs;
1137
1138     t_inputrec ir;
1139     if (MASTER(cr))
1140     {
1141         read_tpr_file(opt2fn("-s", NFILE, fnm), &info, &state, &mtop, &ir, user_beta, fracself);
1142         /* Open logfile for reading */
1143         fp = fopen(opt2fn("-o", NFILE, fnm), "w");
1144
1145         /* Determine the volume of the simulation box */
1146         info.volume = det(state.box);
1147         calc_recipbox(state.box, info.recipbox);
1148         info.natoms = mtop.natoms;
1149         info.bTUNE  = bTUNE;
1150     }
1151
1152     /* Check consistency if the user provided fourierspacing */
1153     if (fs > 0 && MASTER(cr))
1154     {
1155         /* Recalculate the grid dimensions using fourierspacing from user input */
1156         info.nkx[0] = 0;
1157         info.nky[0] = 0;
1158         info.nkz[0] = 0;
1159         calcFftGrid(stdout, state.box, info.fourier_sp[0], minimalPmeGridSize(info.pme_order[0]),
1160                     &(info.nkx[0]), &(info.nky[0]), &(info.nkz[0]));
1161         if ( (ir.nkx != info.nkx[0]) || (ir.nky != info.nky[0]) || (ir.nkz != info.nkz[0]) )
1162         {
1163             gmx_fatal(FARGS, "Wrong fourierspacing %f nm, input file grid = %d x %d x %d, computed grid = %d x %d x %d",
1164                       fs, ir.nkx, ir.nky, ir.nkz, info.nkx[0], info.nky[0], info.nkz[0]);
1165         }
1166     }
1167
1168     /* Estimate (S)PME force error */
1169
1170     if (PAR(cr))
1171     {
1172         bcast_info(&info, cr);
1173     }
1174
1175     /* Get an error estimate of the input tpr file and do some tuning if requested */
1176     estimate_PME_error(&info, &state, &mtop, fp, bVerbose, seed, cr);
1177
1178     if (MASTER(cr))
1179     {
1180         /* Write out optimized tpr file if requested */
1181         if (opt2bSet("-so", NFILE, fnm) || bTUNE)
1182         {
1183             ir.ewald_rtol = info.ewald_rtol[0];
1184             write_tpx_state(opt2fn("-so", NFILE, fnm), &ir, &state, &mtop);
1185         }
1186         please_cite(fp, "Wang2010");
1187         fclose(fp);
1188     }
1189
1190     return 0;
1191 }