2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
37 #include "pme_error.h"
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"
69 /* #define TAKETIME */
72 /* Enum for situations that can occur during log file parsing */
78 eParselogResetProblem,
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 */
109 /* Returns TRUE when atom is charged */
110 static gmx_bool is_charge(real charge)
112 return charge*charge > GMX_REAL_EPS;
116 /* calculate charge density */
117 static void calc_q2all(const gmx_mtop_t *mtop, /* molecular topology */
118 real *q2all, real *q2allnr)
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 */
126 fprintf(stderr, "\nCharge density:\n");
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 */
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++)
137 qi = molecule.atoms.atom[i].q;
138 /* Is this charge worth to be considered? */
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;
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);
160 /* Estimate the direct space part error of the SPME Ewald sum */
161 static real estimate_direct(
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 */
170 beta = info->ewald_beta[0];
171 r_coulomb = info->rcoulomb[0];
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);
176 return ONE_4PI_EPS0*e_dir;
181 /* the following 4 functions determine polynomials required for the reciprocal error estimate */
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 */
189 real nom = 0; /* nominator */
190 real denom = 0; /* denominator */
198 for (i = -SUMORDER; i < 0; i++)
202 nom += std::pow( tmp, -n );
205 for (i = SUMORDER; i > 0; i--)
209 nom += std::pow( tmp, -n );
214 denom = std::pow( tmp, -n )+nom;
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 */
226 real nom = 0; /* nominator */
227 real denom = 0; /* denominator */
235 for (i = -SUMORDER; i < 0; i++)
239 nom += std::pow( tmp, -2*n );
242 for (i = SUMORDER; i > 0; i--)
246 nom += std::pow( tmp, -2*n );
249 for (i = -SUMORDER; i < SUMORDER+1; i++)
253 denom += std::pow( tmp, -n );
255 tmp = eps_poly1(m, K, n);
256 return nom / denom / denom + tmp*tmp;
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 */
266 real nom = 0; /* nominator */
267 real denom = 0; /* denominator */
275 for (i = -SUMORDER; i < 0; i++)
279 nom += i * std::pow( tmp, -2*n );
282 for (i = SUMORDER; i > 0; i--)
286 nom += i * std::pow( tmp, -2*n );
289 for (i = -SUMORDER; i < SUMORDER+1; i++)
293 denom += std::pow( tmp, -n );
296 return 2.0 * M_PI * nom / denom / denom;
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 */
306 real nom = 0; /* nominator */
307 real denom = 0; /* denominator */
315 for (i = -SUMORDER; i < 0; i++)
319 nom += i * i * std::pow( tmp, -2*n );
322 for (i = SUMORDER; i > 0; i--)
326 nom += i * i * std::pow( tmp, -2*n );
329 for (i = -SUMORDER; i < SUMORDER+1; i++)
333 denom += std::pow( tmp, -n );
336 return 4.0 * M_PI * M_PI * nom / denom / denom;
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 */
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 */
361 rcoord = iprod(rboxv, x);
364 for (i = -SUMORDER; i < 0; i++)
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;
373 for (i = SUMORDER; i > 0; i--)
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;
383 tmp = 2.0 * M_PI * m / K;
387 return 2.0 * M_PI * nom / denom * K;
393 /* The following routine is just a copy from pme.c */
395 static void calc_recipbox(matrix box, matrix recipbox)
397 /* Save some time by assuming upper right part is zero */
399 real tmp = 1.0/(box[XX][XX]*box[YY][YY]*box[ZZ][ZZ]);
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;
413 /* Estimate the reciprocal space part error of the SPME Ewald sum. */
414 static real estimate_reciprocal(
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,
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 */
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 */
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 */
448 int *numbers = nullptr;
450 /* Index variables for parallel work distribution */
451 int startglobal, stopglobal;
452 int startlocal, stoplocal;
463 seed = static_cast<int>(gmx::makeRandomSeed());
465 fprintf(stderr, "Using random seed %d.\n", seed);
467 gmx::DefaultRandomEngine rng(seed);
468 gmx::UniformIntDistribution<int> dist(0, nr-1);
476 for (i = 0; i < nr; i++)
481 /* Calculate indices for work distribution */
482 startglobal = -info->nkx[0]/2;
483 stopglobal = info->nkx[0]/2;
484 xtot = stopglobal*2+1;
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)
492 stoplocal = stopglobal;
497 startlocal = startglobal;
498 stoplocal = stopglobal;
514 fprintf(stderr, "Calculating reciprocal error part 1 ...");
518 for (nx = startlocal; nx <= stoplocal; nx++)
520 svmul(nx, info->recipbox[XX], gridpx);
521 for (ny = -info->nky[0]/2; ny < info->nky[0]/2+1; ny++)
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++)
527 if (0 == nx && 0 == ny && 0 == nz)
531 svmul(nz, info->recipbox[ZZ], tmpvec);
532 rvec_add(gridpxy, tmpvec, 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;
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]);
543 tmp1 = eps_poly1(nx, info->nkx[0], info->pme_order[0]);
544 tmp2 = eps_poly1(ny, info->nky[0], info->pme_order[0]);
546 tmp += 2.0 * tmp1 * tmp2;
548 tmp1 = eps_poly1(nz, info->nkz[0], info->pme_order[0]);
549 tmp2 = eps_poly1(ny, info->nky[0], info->pme_order[0]);
551 tmp += 2.0 * tmp1 * tmp2;
553 tmp1 = eps_poly1(nz, info->nkz[0], info->pme_order[0]);
554 tmp2 = eps_poly1(nx, info->nkx[0], info->pme_order[0]);
556 tmp += 2.0 * tmp1 * tmp2;
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]);
564 e_rec1 += 32.0 * M_PI * M_PI * coeff * coeff * coeff2 * tmp * q2_all * q2_all / nr;
566 tmp1 = eps_poly3(nx, info->nkx[0], info->pme_order[0]);
567 tmp1 *= info->nkx[0];
568 tmp2 = iprod(gridp, info->recipbox[XX]);
572 tmp1 = eps_poly3(ny, info->nky[0], info->pme_order[0]);
573 tmp1 *= info->nky[0];
574 tmp2 = iprod(gridp, info->recipbox[YY]);
578 tmp1 = eps_poly3(nz, info->nkz[0], info->pme_order[0]);
579 tmp1 *= info->nkz[0];
580 tmp2 = iprod(gridp, info->recipbox[ZZ]);
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];
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];
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];
604 e_rec2 += 4.0 * coeff * coeff * tmp * q2_all * q2_all / nr;
610 fprintf(stderr, "\rCalculating reciprocal error part 1 ... %3.0f%%", 100.0*(nx-startlocal+1)/(x_per_core));
618 fprintf(stderr, "\n");
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);
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;
634 /* In this case we use all nr particle positions */
636 x_per_core = static_cast<int>(std::ceil(static_cast<real>(xtot) / cr->nnodes));
639 startlocal = x_per_core * cr->nodeid;
640 stoplocal = std::min(startlocal + x_per_core, xtot); /* min needed if xtot == nr */
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 */
650 for (i = 0; i < xtot; i++)
652 numbers[i] = dist(rng); // [0,nr-1]
655 /* Broadcast the random number array to the other nodes */
658 nblock_bc(cr, xtot, numbers);
661 if (bVerbose && MASTER(cr))
663 fprintf(stdout, "Using %d sample%s to approximate the self interaction error term",
664 xtot, xtot == 1 ? "" : "s");
667 fprintf(stdout, " (%d sample%s per rank)", x_per_core, x_per_core == 1 ? "" : "s");
669 fprintf(stdout, ".\n");
673 /* Return the number of positions used for the Monte Carlo algorithm */
676 for (i = startlocal; i < stoplocal; i++)
684 /* Randomly pick a charge */
689 /* Use all charges */
693 /* for(nx=startlocal; nx<=stoplocal; nx++)*/
694 for (nx = -info->nkx[0]/2; nx < info->nkx[0]/2+1; nx++)
696 svmul(nx, info->recipbox[XX], gridpx);
697 for (ny = -info->nky[0]/2; ny < info->nky[0]/2+1; ny++)
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++)
704 if (0 == nx && 0 == ny && 0 == nz)
709 svmul(nz, info->recipbox[ZZ], tmpvec);
710 rvec_add(gridpxy, tmpvec, gridp);
712 coeff = std::exp(-1.0 * M_PI * M_PI * tmp / info->ewald_beta[0] / info->ewald_beta[0] );
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]);
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);
731 e_rec3 += q[ci]*q[ci]*q[ci]*q[ci]*norm2(tmpvec2) / ( xtot * M_PI * info->volume * M_PI * info->volume);
734 fprintf(stderr, "\rCalculating reciprocal error part 2 ... %3.0f%%",
735 100.0*(i+1)/stoplocal);
742 fprintf(stderr, "\n");
749 t1 = MPI_Wtime() - t0;
750 fprintf(fp_out, "Recip. err. est. took : %lf s\n", t1);
758 fprintf(stderr, "Rank %3d: nx=[%3d...%3d] e_rec3=%e\n",
759 cr->nodeid, startlocal, stoplocal, e_rec3);
765 gmx_sum(1, &e_rec1, cr);
766 gmx_sum(1, &e_rec2, cr);
767 gmx_sum(1, &e_rec3, cr);
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 ;
774 e_rec = std::sqrt(e_rec1+e_rec2+e_rec3);
777 return ONE_4PI_EPS0 * e_rec;
781 /* Allocate memory for the inputinfo struct: */
782 static void create_info(t_inputinfo *info)
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);
800 /* Allocate and fill an array with coordinates and charges,
801 * returns the number of charges found
803 static int prepare_x_q(real *q[], rvec *x[], const gmx_mtop_t *mtop, const rvec x_orig[], t_commrec *cr)
805 int nq; /* number of charged particles */
810 snew(*q, mtop->natoms);
811 snew(*x, mtop->natoms);
814 for (const AtomProxy atomP : AtomRange(*mtop))
816 const t_atom &local = atomP.atom();
817 int i = atomP.globalAtomNumber();
818 if (is_charge(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];
827 /* Give back some unneeded memory */
831 /* Broadcast x and q in the parallel case */
834 /* Transfer the number of charges */
838 nblock_bc(cr, nq, *x);
839 nblock_bc(cr, nq, *q);
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)
850 read_tpx_state(fn_sim_tpr, ir, state, mtop);
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;
865 info->ewald_beta[0] = user_beta;
869 info->ewald_beta[0] = calc_ewaldcoeff_q(info->rcoulomb[0], info->ewald_rtol[0]);
872 /* Check if PME was chosen */
873 if (EEL_PME(ir->coulombtype) == FALSE)
875 gmx_fatal(FARGS, "Can only do optimizations for simulations with PME");
878 /* Check if rcoulomb == rlist, which is necessary for PME */
879 if (!(ir->rcoulomb == ir->rlist))
881 gmx_fatal(FARGS, "PME requires rcoulomb (%f) to be equal to rlist (%f).", ir->rcoulomb, ir->rlist);
886 /* Transfer what we need for parallelizing the reciprocal error estimate */
887 static void bcast_info(t_inputinfo *info, const t_commrec *cr)
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);
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.
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,
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 */
929 fprintf(fp_out, "\n--- PME ERROR ESTIMATE ---\n");
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);
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]);
953 bcast_info(info, cr);
957 /* Calculate direct space error */
958 info->e_dir[0] = estimate_direct(info);
960 /* Calculate reciprocal space error */
961 info->e_rec[0] = estimate_reciprocal(info, x, q, ncharges, fp_out, bVerbose,
962 seed, &nsamples, cr);
966 bcast_info(info, cr);
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);
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]);
985 fprintf(stderr, "Starting tuning ...\n");
987 edir = info->e_dir[0];
988 erec = info->e_rec[0];
990 beta0 = info->ewald_beta[0];
993 info->ewald_beta[0] += 0.1;
997 info->ewald_beta[0] -= 0.1;
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);
1005 bcast_info(info, cr);
1009 edir = info->e_dir[0];
1010 erec = info->e_rec[0];
1012 while (std::abs(derr/std::min(erec, edir)) > 1e-4)
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;
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);
1027 bcast_info(info, cr);
1030 edir = info->e_dir[0];
1031 erec = info->e_rec[0];
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);
1043 info->ewald_rtol[0] = std::erfc(info->rcoulomb[0]*info->ewald_beta[0]);
1047 /* Write some info to log file */
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]);
1065 int gmx_pme_error(int argc, char *argv[])
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]",
1078 real fs = 0.0; /* 0 indicates: not set by the user */
1079 real user_beta = -1.0;
1080 real fracself = 1.0;
1082 t_state state; /* The state from the tpr input file */
1083 gmx_mtop_t mtop; /* The topology from the tpr input file */
1086 unsigned long PCA_Flags;
1087 gmx_bool bTUNE = FALSE;
1088 gmx_bool bVerbose = FALSE;
1092 static t_filenm fnm[] = {
1093 { efTPR, "-s", nullptr, ffREAD },
1094 { efOUT, "-o", "error", ffWRITE },
1095 { efTPR, "-so", "tuned", ffOPTWR }
1098 gmx_output_env_t *oenv = nullptr;
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" }
1114 #define NFILE asize(fnm)
1116 cr = init_commrec();
1117 PCA_Flags = PCA_NOEXIT_ON_ARGS;
1119 if (!parse_common_args(&argc, argv, PCA_Flags,
1120 NFILE, fnm, asize(pa), pa, asize(desc), desc,
1129 bTUNE = opt2bSet("-so", NFILE, fnm);
1134 /* Allocate memory for the inputinfo struct: */
1136 info.fourier_sp[0] = fs;
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");
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;
1152 /* Check consistency if the user provided fourierspacing */
1153 if (fs > 0 && MASTER(cr))
1155 /* Recalculate the grid dimensions using fourierspacing from user input */
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]) )
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]);
1168 /* Estimate (S)PME force error */
1172 bcast_info(&info, cr);
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);
1180 /* Write out optimized tpr file if requested */
1181 if (opt2bSet("-so", NFILE, fnm) || bTUNE)
1183 ir.ewald_rtol = info.ewald_rtol[0];
1184 write_tpx_state(opt2fn("-so", NFILE, fnm), &ir, &state, &mtop);
1186 please_cite(fp, "Wang2010");