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,2015,2016,2017, 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.
39 #include "long-range-correction.h"
43 #include "gromacs/math/functions.h"
44 #include "gromacs/math/units.h"
45 #include "gromacs/math/utilities.h"
46 #include "gromacs/math/vec.h"
47 #include "gromacs/mdtypes/commrec.h"
48 #include "gromacs/mdtypes/forcerec.h"
49 #include "gromacs/mdtypes/md_enums.h"
50 #include "gromacs/utility/fatalerror.h"
51 #include "gromacs/utility/gmxassert.h"
53 /* There's nothing special to do here if just masses are perturbed,
54 * but if either charge or type is perturbed then the implementation
55 * requires that B states are defined for both charge and type, and
56 * does not optimize for the cases where only one changes.
58 * The parameter vectors for B states are left undefined in atoms2md()
59 * when either FEP is inactive, or when there are no mass/charge/type
60 * perturbations. The parameter vectors for LJ-PME are likewise
61 * undefined when LJ-PME is not active. This works because
62 * bHaveChargeOrTypePerturbed handles the control flow. */
63 void ewald_LRcorrection(int numAtomsLocal
,
65 int numThreads
, int thread
,
67 real
*chargeA
, real
*chargeB
,
69 real
*sigmaA
, real
*sigmaB
,
70 real
*sigma3A
, real
*sigma3B
,
71 gmx_bool bHaveChargeOrTypePerturbed
,
72 gmx_bool calc_excl_corr
,
73 t_blocka
*excl
, rvec x
[],
74 matrix box
, rvec mu_tot
[],
75 int ewald_geometry
, real epsilon_surface
,
76 rvec
*f
, tensor vir_q
, tensor vir_lj
,
77 real
*Vcorr_q
, real
*Vcorr_lj
,
78 real lambda_q
, real lambda_lj
,
79 real
*dvdlambda_q
, real
*dvdlambda_lj
)
81 int numAtomsToBeCorrected
;
84 /* We need to correct all exclusion pairs (cutoff-scheme = group) */
85 numAtomsToBeCorrected
= excl
->nr
;
87 GMX_RELEASE_ASSERT(numAtomsToBeCorrected
>= numAtomsLocal
, "We might need to do self-corrections");
91 /* We need to correct only self interactions */
92 numAtomsToBeCorrected
= numAtomsLocal
;
94 int start
= (numAtomsToBeCorrected
* thread
)/numThreads
;
95 int end
= (numAtomsToBeCorrected
*(thread
+ 1))/numThreads
;
97 int i
, i1
, i2
, j
, k
, m
, iv
, jv
, q
;
99 double Vexcl_q
, dvdl_excl_q
, dvdl_excl_lj
; /* Necessary for precision */
102 real v
, vc
, qiA
, qiB
, dr2
, rinv
;
103 real Vself_q
[2], Vself_lj
[2], Vdipole
[2], rinv2
, ewc_q
= fr
->ic
->ewaldcoeff_q
, ewcdr
;
104 real ewc_lj
= fr
->ic
->ewaldcoeff_lj
, ewc_lj2
= ewc_lj
* ewc_lj
;
105 real c6Ai
= 0, c6Bi
= 0, c6A
= 0, c6B
= 0, ewcdr2
, ewcdr4
, c6L
= 0, rinv6
;
106 rvec df
, dx
, mutot
[2], dipcorrA
, dipcorrB
;
107 tensor dxdf_q
= {{0}}, dxdf_lj
= {{0}};
108 real vol
= box
[XX
][XX
]*box
[YY
][YY
]*box
[ZZ
][ZZ
];
109 real L1_q
, L1_lj
, dipole_coeff
, qqA
, qqB
, qqL
, vr0_q
, vr0_lj
= 0;
110 real chargecorr
[2] = { 0, 0 };
111 gmx_bool bMolPBC
= fr
->bMolPBC
;
112 gmx_bool bDoingLBRule
= (fr
->ljpme_combination_rule
== eljpmeLB
);
113 gmx_bool bNeedLongRangeCorrection
;
115 /* This routine can be made faster by using tables instead of analytical interactions
116 * However, that requires a thorough verification that they are correct in all cases.
119 bool vdwPme
= EVDW_PME(fr
->ic
->vdwtype
);
121 one_4pi_eps
= ONE_4PI_EPS0
/fr
->ic
->epsilon_r
;
122 vr0_q
= ewc_q
*M_2_SQRTPI
;
125 vr0_lj
= -gmx::power6(ewc_lj
)/6.0;
136 L1_lj
= 1.0-lambda_lj
;
137 /* Note that we have to transform back to gromacs units, since
138 * mu_tot contains the dipole in debye units (for output).
140 for (i
= 0; (i
< DIM
); i
++)
142 mutot
[0][i
] = mu_tot
[0][i
]*DEBYE2ENM
;
143 mutot
[1][i
] = mu_tot
[1][i
]*DEBYE2ENM
;
148 switch (ewald_geometry
)
151 if (epsilon_surface
!= 0)
154 2*M_PI
*ONE_4PI_EPS0
/((2*epsilon_surface
+ fr
->ic
->epsilon_r
)*vol
);
155 for (i
= 0; (i
< DIM
); i
++)
157 dipcorrA
[i
] = 2*dipole_coeff
*mutot
[0][i
];
158 dipcorrB
[i
] = 2*dipole_coeff
*mutot
[1][i
];
163 dipole_coeff
= 2*M_PI
*one_4pi_eps
/vol
;
164 dipcorrA
[ZZ
] = 2*dipole_coeff
*mutot
[0][ZZ
];
165 dipcorrB
[ZZ
] = 2*dipole_coeff
*mutot
[1][ZZ
];
166 for (int q
= 0; q
< (bHaveChargeOrTypePerturbed
? 2 : 1); q
++)
168 /* Avoid charge corrections with near-zero net charge */
169 if (fabs(fr
->qsum
[q
]) > 1e-4)
171 chargecorr
[q
] = 2*dipole_coeff
*fr
->qsum
[q
];
176 gmx_incons("Unsupported Ewald geometry");
181 fprintf(debug
, "dipcorr = %8.3f %8.3f %8.3f\n",
182 dipcorrA
[XX
], dipcorrA
[YY
], dipcorrA
[ZZ
]);
183 fprintf(debug
, "mutot = %8.3f %8.3f %8.3f\n",
184 mutot
[0][XX
], mutot
[0][YY
], mutot
[0][ZZ
]);
186 bNeedLongRangeCorrection
= (calc_excl_corr
|| dipole_coeff
!= 0);
187 if (bNeedLongRangeCorrection
&& !bHaveChargeOrTypePerturbed
)
189 for (i
= start
; (i
< end
); i
++)
191 /* Initiate local variables (for this i-particle) to 0 */
192 qiA
= chargeA
[i
]*one_4pi_eps
;
204 i2
= excl
->index
[i
+1];
206 /* Loop over excluded neighbours */
207 for (j
= i1
; (j
< i2
); j
++)
211 * First we must test whether k <> i, and then,
212 * because the exclusions are all listed twice i->k
213 * and k->i we must select just one of the two. As
214 * a minor optimization we only compute forces when
215 * the charges are non-zero.
219 qqA
= qiA
*chargeA
[k
];
225 c6A
*= gmx::power6(0.5*(sigmaA
[i
]+sigmaA
[k
]))*sigma3A
[k
];
228 if (qqA
!= 0.0 || c6A
!= 0.0)
230 rvec_sub(x
[i
], x
[k
], dx
);
233 /* Cheap pbc_dx, assume excluded pairs are at short distance. */
234 for (m
= DIM
-1; (m
>= 0); m
--)
236 if (dx
[m
] > 0.5*box
[m
][m
])
238 rvec_dec(dx
, box
[m
]);
240 else if (dx
[m
] < -0.5*box
[m
][m
])
242 rvec_inc(dx
, box
[m
]);
247 /* Distance between two excluded particles
248 * may be zero in the case of shells
252 rinv
= gmx::invsqrt(dr2
);
260 vc
= qqA
*std::erf(ewcdr
)*rinv
;
263 /* Relative accuracy at R_ERF_R_INACC of 3e-10 */
264 #define R_ERF_R_INACC 0.006
266 /* Relative accuracy at R_ERF_R_INACC of 2e-5 */
267 #define R_ERF_R_INACC 0.1
269 /* fscal is the scalar force pre-multiplied by rinv,
270 * to normalise the relative position vector dx */
271 if (ewcdr
> R_ERF_R_INACC
)
273 fscal
= rinv2
*(vc
- qqA
*ewc_q
*M_2_SQRTPI
*exp(-ewcdr
*ewcdr
));
277 /* Use a fourth order series expansion for small ewcdr */
278 fscal
= ewc_q
*ewc_q
*qqA
*vr0_q
*(2.0/3.0 - 0.4*ewcdr
*ewcdr
);
281 /* The force vector is obtained by multiplication with
282 * the relative position vector
284 svmul(fscal
, dx
, df
);
287 for (iv
= 0; (iv
< DIM
); iv
++)
289 for (jv
= 0; (jv
< DIM
); jv
++)
291 dxdf_q
[iv
][jv
] += dx
[iv
]*df
[jv
];
300 rinv6
= rinv2
*rinv2
*rinv2
;
301 ewcdr2
= ewc_lj2
*dr2
;
302 ewcdr4
= ewcdr2
*ewcdr2
;
303 /* We get the excluded long-range contribution from -C6*(1-g(r))
304 * g(r) is also defined in the manual under LJ-PME
306 vc
= -c6A
*rinv6
*(1.0 - exp(-ewcdr2
)*(1 + ewcdr2
+ 0.5*ewcdr4
));
308 /* The force is the derivative of the potential vc.
309 * fscal is the scalar force pre-multiplied by rinv,
310 * to normalise the relative position vector dx */
311 fscal
= 6.0*vc
*rinv2
+ c6A
*rinv6
*exp(-ewcdr2
)*ewc_lj2
*ewcdr4
;
313 /* The force vector is obtained by multiplication with
314 * the relative position vector
316 svmul(fscal
, dx
, df
);
319 for (iv
= 0; (iv
< DIM
); iv
++)
321 for (jv
= 0; (jv
< DIM
); jv
++)
323 dxdf_lj
[iv
][jv
] += dx
[iv
]*df
[jv
];
330 Vexcl_q
+= qqA
*vr0_q
;
331 Vexcl_lj
+= c6A
*vr0_lj
;
337 /* Dipole correction on force */
338 if (dipole_coeff
!= 0 && i
< numAtomsLocal
)
340 for (j
= 0; (j
< DIM
); j
++)
342 f
[i
][j
] -= dipcorrA
[j
]*chargeA
[i
];
344 if (chargecorr
[0] != 0)
346 f
[i
][ZZ
] += chargecorr
[0]*chargeA
[i
]*x
[i
][ZZ
];
351 else if (bNeedLongRangeCorrection
)
353 for (i
= start
; (i
< end
); i
++)
355 /* Initiate local variables (for this i-particle) to 0 */
356 qiA
= chargeA
[i
]*one_4pi_eps
;
357 qiB
= chargeB
[i
]*one_4pi_eps
;
371 i2
= excl
->index
[i
+1];
373 /* Loop over excluded neighbours */
374 for (j
= i1
; (j
< i2
); j
++)
379 qqA
= qiA
*chargeA
[k
];
380 qqB
= qiB
*chargeB
[k
];
387 c6A
*= gmx::power6(0.5*(sigmaA
[i
]+sigmaA
[k
]))*sigma3A
[k
];
388 c6B
*= gmx::power6(0.5*(sigmaB
[i
]+sigmaB
[k
]))*sigma3B
[k
];
391 if (qqA
!= 0.0 || qqB
!= 0.0 || c6A
!= 0.0 || c6B
!= 0.0)
395 qqL
= L1_q
*qqA
+ lambda_q
*qqB
;
398 c6L
= L1_lj
*c6A
+ lambda_lj
*c6B
;
400 rvec_sub(x
[i
], x
[k
], dx
);
403 /* Cheap pbc_dx, assume excluded pairs are at short distance. */
404 for (m
= DIM
-1; (m
>= 0); m
--)
406 if (dx
[m
] > 0.5*box
[m
][m
])
408 rvec_dec(dx
, box
[m
]);
410 else if (dx
[m
] < -0.5*box
[m
][m
])
412 rvec_inc(dx
, box
[m
]);
419 rinv
= gmx::invsqrt(dr2
);
421 if (qqA
!= 0.0 || qqB
!= 0.0)
426 v
= std::erf(ewc_q
*dr
)*rinv
;
429 /* fscal is the scalar force pre-multiplied by rinv,
430 * to normalise the relative position vector dx */
431 fscal
= rinv2
*(vc
-qqL
*ewc_q
*M_2_SQRTPI
*exp(-ewc_q
*ewc_q
*dr2
));
432 dvdl_excl_q
+= (qqB
- qqA
)*v
;
434 /* The force vector is obtained by multiplication with
435 * the relative position vector
437 svmul(fscal
, dx
, df
);
440 for (iv
= 0; (iv
< DIM
); iv
++)
442 for (jv
= 0; (jv
< DIM
); jv
++)
444 dxdf_q
[iv
][jv
] += dx
[iv
]*df
[jv
];
449 if ((c6A
!= 0.0 || c6B
!= 0.0) && vdwPme
)
451 rinv6
= rinv2
*rinv2
*rinv2
;
452 ewcdr2
= ewc_lj2
*dr2
;
453 ewcdr4
= ewcdr2
*ewcdr2
;
454 v
= -rinv6
*(1.0 - exp(-ewcdr2
)*(1 + ewcdr2
+ 0.5*ewcdr4
));
457 /* fscal is the scalar force pre-multiplied by rinv,
458 * to normalise the relative position vector dx */
459 fscal
= 6.0*vc
*rinv2
+ c6L
*rinv6
*exp(-ewcdr2
)*ewc_lj2
*ewcdr4
;
460 dvdl_excl_lj
+= (c6B
- c6A
)*v
;
462 /* The force vector is obtained by multiplication with
463 * the relative position vector
465 svmul(fscal
, dx
, df
);
468 for (iv
= 0; (iv
< DIM
); iv
++)
470 for (jv
= 0; (jv
< DIM
); jv
++)
472 dxdf_lj
[iv
][jv
] += dx
[iv
]*df
[jv
];
479 Vexcl_q
+= qqL
*vr0_q
;
480 dvdl_excl_q
+= (qqB
- qqA
)*vr0_q
;
481 Vexcl_lj
+= c6L
*vr0_lj
;
482 dvdl_excl_lj
+= (c6B
- c6A
)*vr0_lj
;
488 /* Dipole correction on force */
489 if (dipole_coeff
!= 0 && i
< numAtomsLocal
)
491 for (j
= 0; (j
< DIM
); j
++)
493 f
[i
][j
] -= L1_q
*dipcorrA
[j
]*chargeA
[i
]
494 + lambda_q
*dipcorrB
[j
]*chargeB
[i
];
496 if (chargecorr
[0] != 0 || chargecorr
[1] != 0)
498 f
[i
][ZZ
] += (L1_q
*chargecorr
[0]*chargeA
[i
]
499 + lambda_q
*chargecorr
[1])*x
[i
][ZZ
];
504 for (iv
= 0; (iv
< DIM
); iv
++)
506 for (jv
= 0; (jv
< DIM
); jv
++)
508 vir_q
[iv
][jv
] += 0.5*dxdf_q
[iv
][jv
];
509 vir_lj
[iv
][jv
] += 0.5*dxdf_lj
[iv
][jv
];
518 /* Global corrections only on master process */
519 if (MASTER(cr
) && thread
== 0)
521 for (q
= 0; q
< (bHaveChargeOrTypePerturbed
? 2 : 1); q
++)
525 /* Self-energy correction */
526 Vself_q
[q
] = ewc_q
*one_4pi_eps
*fr
->q2sum
[q
]*M_1_SQRTPI
;
529 Vself_lj
[q
] = fr
->c6sum
[q
]*0.5*vr0_lj
;
533 /* Apply surface and charged surface dipole correction:
534 * correction = dipole_coeff * ( (dipole)^2
535 * - qsum*sum_i q_i z_i^2 - qsum^2 * box_z^2 / 12 )
537 if (dipole_coeff
!= 0)
539 if (ewald_geometry
== eewg3D
)
541 Vdipole
[q
] = dipole_coeff
*iprod(mutot
[q
], mutot
[q
]);
543 else if (ewald_geometry
== eewg3DC
)
545 Vdipole
[q
] = dipole_coeff
*mutot
[q
][ZZ
]*mutot
[q
][ZZ
];
547 if (chargecorr
[q
] != 0)
549 /* Here we use a non thread-parallelized loop,
550 * because this is the only loop over atoms for
551 * energies and they need reduction (unlike forces).
552 * We could implement a reduction over threads,
553 * but this case is rarely used.
555 const real
*qPtr
= (q
== 0 ? chargeA
: chargeB
);
557 for (int i
= 0; i
< numAtomsLocal
; i
++)
559 sumQZ2
+= qPtr
[i
]*x
[i
][ZZ
]*x
[i
][ZZ
];
561 Vdipole
[q
] -= dipole_coeff
*fr
->qsum
[q
]*(sumQZ2
+ fr
->qsum
[q
]*box
[ZZ
][ZZ
]*box
[ZZ
][ZZ
]/12);
567 if (!bHaveChargeOrTypePerturbed
)
569 *Vcorr_q
= Vdipole
[0] - Vself_q
[0] - Vexcl_q
;
572 *Vcorr_lj
= -Vself_lj
[0] - Vexcl_lj
;
577 *Vcorr_q
= L1_q
*(Vdipole
[0] - Vself_q
[0])
578 + lambda_q
*(Vdipole
[1] - Vself_q
[1])
580 *dvdlambda_q
+= Vdipole
[1] - Vself_q
[1]
581 - (Vdipole
[0] - Vself_q
[0]) - dvdl_excl_q
;
584 *Vcorr_lj
= -(L1_lj
*Vself_lj
[0] + lambda_lj
*Vself_lj
[1]) - Vexcl_lj
;
585 *dvdlambda_lj
+= -Vself_lj
[1] + Vself_lj
[0] - dvdl_excl_lj
;
591 fprintf(debug
, "Long Range corrections for Ewald interactions:\n");
592 fprintf(debug
, "q2sum = %g, Vself_q=%g c6sum = %g, Vself_lj=%g\n",
593 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]);
594 fprintf(debug
, "Electrostatic Long Range correction: Vexcl=%g\n", Vexcl_q
);
595 fprintf(debug
, "Lennard-Jones Long Range correction: Vexcl=%g\n", Vexcl_lj
);
596 if (MASTER(cr
) && thread
== 0)
598 if (epsilon_surface
> 0 || ewald_geometry
== eewg3DC
)
600 fprintf(debug
, "Total dipole correction: Vdipole=%g\n",
601 L1_q
*Vdipole
[0]+lambda_q
*Vdipole
[1]);