2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
43 #include "gromacs/math/utilities.h"
45 #include "types/commrec.h"
46 #include "gromacs/math/vec.h"
48 #include "gromacs/utility/smalloc.h"
49 #include "gromacs/math/units.h"
51 #include "gromacs/utility/futil.h"
55 real calc_ewaldcoeff_q(real rc, real dtol)
57 real x = 5, low, high;
66 while (gmx_erfc(x*rc) > dtol);
68 n = i+60; /* search tolerance is 2^-60 */
71 for (i = 0; i < n; i++)
74 if (gmx_erfc(x*rc) > dtol)
86 static real ewald_function_lj(real x, real rc)
88 real xrc, xrc2, xrc4, factor;
93 factor = exp(-xrc2)*(1 + xrc2 + xrc4/2.0);
95 factor = expf(-xrc2)*(1 + xrc2 + xrc4/2.0);
101 real calc_ewaldcoeff_lj(real rc, real dtol)
103 real x = 5, low, high;
111 while (ewald_function_lj(x, rc) > dtol);
113 n = i + 60; /* search tolerance is 2^-60 */
116 for (i = 0; i < n; ++i)
118 x = (low + high) / 2.0;
119 if (ewald_function_lj(x, rc) > dtol)
131 void ewald_LRcorrection(int start, int end,
132 t_commrec *cr, int thread, t_forcerec *fr,
133 real *chargeA, real *chargeB,
134 real *C6A, real *C6B,
135 real *sigmaA, real *sigmaB,
136 real *sigma3A, real *sigma3B,
137 gmx_bool calc_excl_corr,
138 t_blocka *excl, rvec x[],
139 matrix box, rvec mu_tot[],
140 int ewald_geometry, real epsilon_surface,
141 rvec *f, tensor vir_q, tensor vir_lj,
142 real *Vcorr_q, real *Vcorr_lj,
143 real lambda_q, real lambda_lj,
144 real *dvdlambda_q, real *dvdlambda_lj)
146 int i, i1, i2, j, k, m, iv, jv, q;
148 double Vexcl_q, dvdl_excl_q, dvdl_excl_lj; /* Necessary for precision */
151 real v, vc, qiA, qiB, dr2, rinv, enercorr;
152 real Vself_q[2], Vself_lj[2], Vdipole[2], rinv2, ewc_q = fr->ewaldcoeff_q, ewcdr;
153 real ewc_lj = fr->ewaldcoeff_lj, ewc_lj2 = ewc_lj * ewc_lj;
154 real c6Ai = 0, c6Bi = 0, c6A = 0, c6B = 0, ewcdr2, ewcdr4, c6L = 0, rinv6;
155 rvec df, dx, mutot[2], dipcorrA, dipcorrB;
156 tensor dxdf_q, dxdf_lj;
157 real vol = box[XX][XX]*box[YY][YY]*box[ZZ][ZZ];
158 real L1_q, L1_lj, dipole_coeff, qqA, qqB, qqL, vr0_q, vr0_lj = 0;
159 gmx_bool bFreeEnergy = (chargeB != NULL);
160 gmx_bool bMolPBC = fr->bMolPBC;
161 gmx_bool bDoingLBRule = (fr->ljpme_combination_rule == eljpmeLB);
163 /* This routine can be made faster by using tables instead of analytical interactions
164 * However, that requires a thorough verification that they are correct in all cases.
167 one_4pi_eps = ONE_4PI_EPS0/fr->epsilon_r;
168 vr0_q = ewc_q*M_2_SQRTPI;
169 if (EVDW_PME(fr->vdwtype))
171 vr0_lj = -pow(ewc_lj, 6)/6.0;
182 L1_lj = 1.0-lambda_lj;
183 /* Note that we have to transform back to gromacs units, since
184 * mu_tot contains the dipole in debye units (for output).
186 for (i = 0; (i < DIM); i++)
188 mutot[0][i] = mu_tot[0][i]*DEBYE2ENM;
189 mutot[1][i] = mu_tot[1][i]*DEBYE2ENM;
194 switch (ewald_geometry)
197 if (epsilon_surface != 0)
200 2*M_PI*ONE_4PI_EPS0/((2*epsilon_surface + fr->epsilon_r)*vol);
201 for (i = 0; (i < DIM); i++)
203 dipcorrA[i] = 2*dipole_coeff*mutot[0][i];
204 dipcorrB[i] = 2*dipole_coeff*mutot[1][i];
209 dipole_coeff = 2*M_PI*one_4pi_eps/vol;
210 dipcorrA[ZZ] = 2*dipole_coeff*mutot[0][ZZ];
211 dipcorrB[ZZ] = 2*dipole_coeff*mutot[1][ZZ];
214 gmx_incons("Unsupported Ewald geometry");
219 fprintf(debug, "dipcorr = %8.3f %8.3f %8.3f\n",
220 dipcorrA[XX], dipcorrA[YY], dipcorrA[ZZ]);
221 fprintf(debug, "mutot = %8.3f %8.3f %8.3f\n",
222 mutot[0][XX], mutot[0][YY], mutot[0][ZZ]);
225 if (EVDW_PME(fr->vdwtype))
229 if ((calc_excl_corr || dipole_coeff != 0) && !bFreeEnergy)
231 for (i = start; (i < end); i++)
233 /* Initiate local variables (for this i-particle) to 0 */
234 qiA = chargeA[i]*one_4pi_eps;
235 if (EVDW_PME(fr->vdwtype))
246 i2 = excl->index[i+1];
248 /* Loop over excluded neighbours */
249 for (j = i1; (j < i2); j++)
253 * First we must test whether k <> i, and then,
254 * because the exclusions are all listed twice i->k
255 * and k->i we must select just one of the two. As
256 * a minor optimization we only compute forces when
257 * the charges are non-zero.
261 qqA = qiA*chargeA[k];
262 if (EVDW_PME(fr->vdwtype))
267 c6A *= pow(0.5*(sigmaA[i]+sigmaA[k]), 6)*sigma3A[k];
270 if (qqA != 0.0 || c6A != 0.0)
275 rvec_sub(x[i], x[k], dx);
278 /* Cheap pbc_dx, assume excluded pairs are at short distance. */
279 for (m = DIM-1; (m >= 0); m--)
281 if (dx[m] > 0.5*box[m][m])
283 rvec_dec(dx, box[m]);
285 else if (dx[m] < -0.5*box[m][m])
287 rvec_inc(dx, box[m]);
292 /* Distance between two excluded particles
293 * may be zero in the case of shells
297 rinv = gmx_invsqrt(dr2);
305 vc = qqA*gmx_erf(ewcdr)*rinv;
308 /* Relative accuracy at R_ERF_R_INACC of 3e-10 */
309 #define R_ERF_R_INACC 0.006
311 /* Relative accuracy at R_ERF_R_INACC of 2e-5 */
312 #define R_ERF_R_INACC 0.1
314 /* fscal is the scalar force pre-multiplied by rinv,
315 * to normalise the relative position vector dx */
316 if (ewcdr > R_ERF_R_INACC)
318 fscal = rinv2*(vc - qqA*ewc_q*M_2_SQRTPI*exp(-ewcdr*ewcdr));
322 /* Use a fourth order series expansion for small ewcdr */
323 fscal = ewc_q*ewc_q*qqA*vr0_q*(2.0/3.0 - 0.4*ewcdr*ewcdr);
326 /* The force vector is obtained by multiplication with
327 * the relative position vector
329 svmul(fscal, dx, df);
332 for (iv = 0; (iv < DIM); iv++)
334 for (jv = 0; (jv < DIM); jv++)
336 dxdf_q[iv][jv] += dx[iv]*df[jv];
343 rinv6 = rinv2*rinv2*rinv2;
344 ewcdr2 = ewc_lj2*dr2;
345 ewcdr4 = ewcdr2*ewcdr2;
346 /* We get the excluded long-range contribution from -C6*(1-g(r))
347 * g(r) is also defined in the manual under LJ-PME
349 vc = -c6A*rinv6*(1.0 - exp(-ewcdr2)*(1 + ewcdr2 + 0.5*ewcdr4));
351 /* The force is the derivative of the potential vc.
352 * fscal is the scalar force pre-multiplied by rinv,
353 * to normalise the relative position vector dx */
354 fscal = 6.0*vc*rinv2 + c6A*rinv6*exp(-ewcdr2)*ewc_lj2*ewcdr4;
356 /* The force vector is obtained by multiplication with
357 * the relative position vector
359 svmul(fscal, dx, df);
362 for (iv = 0; (iv < DIM); iv++)
364 for (jv = 0; (jv < DIM); jv++)
366 dxdf_lj[iv][jv] += dx[iv]*df[jv];
373 Vexcl_q += qqA*vr0_q;
374 Vexcl_lj += c6A*vr0_lj;
380 /* Dipole correction on force */
381 if (dipole_coeff != 0)
383 for (j = 0; (j < DIM); j++)
385 f[i][j] -= dipcorrA[j]*chargeA[i];
390 else if (calc_excl_corr || dipole_coeff != 0)
392 for (i = start; (i < end); i++)
394 /* Initiate local variables (for this i-particle) to 0 */
395 qiA = chargeA[i]*one_4pi_eps;
396 qiB = chargeB[i]*one_4pi_eps;
397 if (EVDW_PME(fr->vdwtype))
410 i2 = excl->index[i+1];
412 /* Loop over excluded neighbours */
413 for (j = i1; (j < i2); j++)
418 qqA = qiA*chargeA[k];
419 qqB = qiB*chargeB[k];
420 if (EVDW_PME(fr->vdwtype))
426 c6A *= pow(0.5*(sigmaA[i]+sigmaA[k]), 6)*sigma3A[k];
427 c6B *= pow(0.5*(sigmaB[i]+sigmaB[k]), 6)*sigma3B[k];
430 if (qqA != 0.0 || qqB != 0.0 || c6A != 0.0 || c6B != 0.0)
435 qqL = L1_q*qqA + lambda_q*qqB;
436 if (EVDW_PME(fr->vdwtype))
438 c6L = L1_lj*c6A + lambda_lj*c6B;
440 rvec_sub(x[i], x[k], dx);
443 /* Cheap pbc_dx, assume excluded pairs are at short distance. */
444 for (m = DIM-1; (m >= 0); m--)
446 if (dx[m] > 0.5*box[m][m])
448 rvec_dec(dx, box[m]);
450 else if (dx[m] < -0.5*box[m][m])
452 rvec_inc(dx, box[m]);
459 rinv = gmx_invsqrt(dr2);
461 if (qqA != 0.0 || qqB != 0.0)
466 v = gmx_erf(ewc_q*dr)*rinv;
469 /* fscal is the scalar force pre-multiplied by rinv,
470 * to normalise the relative position vector dx */
471 fscal = rinv2*(vc-qqL*ewc_q*M_2_SQRTPI*exp(-ewc_q*ewc_q*dr2));
472 dvdl_excl_q += (qqB - qqA)*v;
474 /* The force vector is obtained by multiplication with
475 * the relative position vector
477 svmul(fscal, dx, df);
480 for (iv = 0; (iv < DIM); iv++)
482 for (jv = 0; (jv < DIM); jv++)
484 dxdf_q[iv][jv] += dx[iv]*df[jv];
489 if ((c6A != 0.0 || c6B != 0.0) && EVDW_PME(fr->vdwtype))
491 rinv6 = rinv2*rinv2*rinv2;
492 ewcdr2 = ewc_lj2*dr2;
493 ewcdr4 = ewcdr2*ewcdr2;
494 v = -rinv6*(1.0 - exp(-ewcdr2)*(1 + ewcdr2 + 0.5*ewcdr4));
497 /* fscal is the scalar force pre-multiplied by rinv,
498 * to normalise the relative position vector dx */
499 fscal = 6.0*vc*rinv2 + c6L*rinv6*exp(-ewcdr2)*ewc_lj2*ewcdr4;
500 dvdl_excl_lj += (c6B - c6A)*v;
502 /* The force vector is obtained by multiplication with
503 * the relative position vector
505 svmul(fscal, dx, df);
508 for (iv = 0; (iv < DIM); iv++)
510 for (jv = 0; (jv < DIM); jv++)
512 dxdf_lj[iv][jv] += dx[iv]*df[jv];
519 Vexcl_q += qqL*vr0_q;
520 dvdl_excl_q += (qqB - qqA)*vr0_q;
521 Vexcl_lj += c6L*vr0_lj;
522 dvdl_excl_lj += (c6B - c6A)*vr0_lj;
528 /* Dipole correction on force */
529 if (dipole_coeff != 0)
531 for (j = 0; (j < DIM); j++)
533 f[i][j] -= L1_q*dipcorrA[j]*chargeA[i]
534 + lambda_q*dipcorrB[j]*chargeB[i];
539 for (iv = 0; (iv < DIM); iv++)
541 for (jv = 0; (jv < DIM); jv++)
543 vir_q[iv][jv] += 0.5*dxdf_q[iv][jv];
544 vir_lj[iv][jv] += 0.5*dxdf_lj[iv][jv];
553 /* Global corrections only on master process */
554 if (MASTER(cr) && thread == 0)
556 for (q = 0; q < (bFreeEnergy ? 2 : 1); q++)
560 /* Self-energy correction */
561 Vself_q[q] = ewc_q*one_4pi_eps*fr->q2sum[q]*M_1_SQRTPI;
562 if (EVDW_PME(fr->vdwtype))
564 Vself_lj[q] = fr->c6sum[q]*0.5*vr0_lj;
568 /* Apply surface dipole correction:
569 * correction = dipole_coeff * (dipole)^2
571 if (dipole_coeff != 0)
573 if (ewald_geometry == eewg3D)
575 Vdipole[q] = dipole_coeff*iprod(mutot[q], mutot[q]);
577 else if (ewald_geometry == eewg3DC)
579 Vdipole[q] = dipole_coeff*mutot[q][ZZ]*mutot[q][ZZ];
586 *Vcorr_q = Vdipole[0] - Vself_q[0] - Vexcl_q;
587 if (EVDW_PME(fr->vdwtype))
589 *Vcorr_lj = -Vself_lj[0] - Vexcl_lj;
594 *Vcorr_q = L1_q*(Vdipole[0] - Vself_q[0])
595 + lambda_q*(Vdipole[1] - Vself_q[1])
597 *dvdlambda_q += Vdipole[1] - Vself_q[1]
598 - (Vdipole[0] - Vself_q[0]) - dvdl_excl_q;
599 if (EVDW_PME(fr->vdwtype))
601 *Vcorr_lj = -(L1_lj*Vself_lj[0] + lambda_lj*Vself_lj[1]) - Vexcl_lj;
602 *dvdlambda_lj += -Vself_lj[1] + Vself_lj[0] - dvdl_excl_lj;
608 fprintf(debug, "Long Range corrections for Ewald interactions:\n");
609 fprintf(debug, "start=%d,natoms=%d\n", start, end-start);
610 fprintf(debug, "q2sum = %g, Vself_q=%g c6sum = %g, Vself_lj=%g\n",
611 L1_q*fr->q2sum[0]+lambda_q*fr->q2sum[1], L1_q*Vself_q[0]+lambda_q*Vself_q[1], L1_lj*fr->c6sum[0]+lambda_lj*fr->c6sum[1], L1_lj*Vself_lj[0]+lambda_lj*Vself_lj[1]);
612 fprintf(debug, "Electrostatic Long Range correction: Vexcl=%g\n", Vexcl_q);
613 fprintf(debug, "Lennard-Jones Long Range correction: Vexcl=%g\n", Vexcl_lj);
614 if (MASTER(cr) && thread == 0)
616 if (epsilon_surface > 0 || ewald_geometry == eewg3DC)
618 fprintf(debug, "Total dipole correction: Vdipole=%g\n",
619 L1_q*Vdipole[0]+lambda_q*Vdipole[1]);
625 real ewald_charge_correction(t_commrec *cr, t_forcerec *fr, real lambda,
627 real *dvdlambda, tensor vir)
630 real vol, fac, qs2A, qs2B, vc, enercorr;
635 /* Apply charge correction */
636 vol = box[XX][XX]*box[YY][YY]*box[ZZ][ZZ];
638 fac = M_PI*ONE_4PI_EPS0/(fr->epsilon_r*2.0*vol*vol*sqr(fr->ewaldcoeff_q));
640 qs2A = fr->qsum[0]*fr->qsum[0];
641 qs2B = fr->qsum[1]*fr->qsum[1];
643 vc = (qs2A*(1 - lambda) + qs2B*lambda)*fac;
647 *dvdlambda += -vol*(qs2B - qs2A)*fac;
649 for (d = 0; d < DIM; d++)
656 fprintf(debug, "Total charge correction: Vcharge=%g\n", enercorr);