Merge branch release-5-1 into release-2016
[alexxy/gromacs.git] / src / gromacs / ewald / long-range-correction.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team.
6  * Copyright (c) 2013,2014,2015,2016, 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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 #include "gmxpre.h"
38
39 #include "long-range-correction.h"
40
41 #include <cmath>
42
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"
52
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.
57  *
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,
64                         t_commrec *cr,
65                         int numThreads, int thread,
66                         t_forcerec *fr,
67                         real *chargeA, real *chargeB,
68                         real *C6A, real *C6B,
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)
80 {
81     int numAtomsToBeCorrected;
82     if (calc_excl_corr)
83     {
84         /* We need to correct all exclusion pairs (cutoff-scheme = group) */
85         numAtomsToBeCorrected = excl->nr;
86
87         GMX_RELEASE_ASSERT(numAtomsToBeCorrected >= numAtomsLocal, "We might need to do self-corrections");
88     }
89     else
90     {
91         /* We need to correct only self interactions */
92         numAtomsToBeCorrected = numAtomsLocal;
93     }
94     int         start =  (numAtomsToBeCorrected* thread     )/numThreads;
95     int         end   =  (numAtomsToBeCorrected*(thread + 1))/numThreads;
96
97     int         i, i1, i2, j, k, m, iv, jv, q;
98     int        *AA;
99     double      Vexcl_q, dvdl_excl_q, dvdl_excl_lj; /* Necessary for precision */
100     double      Vexcl_lj;
101     real        one_4pi_eps;
102     real        v, vc, qiA, qiB, dr2, rinv;
103     real        Vself_q[2], Vself_lj[2], Vdipole[2], rinv2, ewc_q = fr->ewaldcoeff_q, ewcdr;
104     real        ewc_lj = fr->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     gmx_bool    bMolPBC      = fr->bMolPBC;
111     gmx_bool    bDoingLBRule = (fr->ljpme_combination_rule == eljpmeLB);
112     gmx_bool    bNeedLongRangeCorrection;
113
114     /* This routine can be made faster by using tables instead of analytical interactions
115      * However, that requires a thorough verification that they are correct in all cases.
116      */
117
118     one_4pi_eps   = ONE_4PI_EPS0/fr->epsilon_r;
119     vr0_q         = ewc_q*M_2_SQRTPI;
120     if (EVDW_PME(fr->vdwtype))
121     {
122         vr0_lj    = -gmx::power6(ewc_lj)/6.0;
123     }
124
125     AA           = excl->a;
126     Vexcl_q      = 0;
127     Vexcl_lj     = 0;
128     dvdl_excl_q  = 0;
129     dvdl_excl_lj = 0;
130     Vdipole[0]   = 0;
131     Vdipole[1]   = 0;
132     L1_q         = 1.0-lambda_q;
133     L1_lj        = 1.0-lambda_lj;
134     /* Note that we have to transform back to gromacs units, since
135      * mu_tot contains the dipole in debye units (for output).
136      */
137     for (i = 0; (i < DIM); i++)
138     {
139         mutot[0][i] = mu_tot[0][i]*DEBYE2ENM;
140         mutot[1][i] = mu_tot[1][i]*DEBYE2ENM;
141         dipcorrA[i] = 0;
142         dipcorrB[i] = 0;
143     }
144     dipole_coeff = 0;
145     switch (ewald_geometry)
146     {
147         case eewg3D:
148             if (epsilon_surface != 0)
149             {
150                 dipole_coeff =
151                     2*M_PI*ONE_4PI_EPS0/((2*epsilon_surface + fr->epsilon_r)*vol);
152                 for (i = 0; (i < DIM); i++)
153                 {
154                     dipcorrA[i] = 2*dipole_coeff*mutot[0][i];
155                     dipcorrB[i] = 2*dipole_coeff*mutot[1][i];
156                 }
157             }
158             break;
159         case eewg3DC:
160             dipole_coeff = 2*M_PI*one_4pi_eps/vol;
161             dipcorrA[ZZ] = 2*dipole_coeff*mutot[0][ZZ];
162             dipcorrB[ZZ] = 2*dipole_coeff*mutot[1][ZZ];
163             break;
164         default:
165             gmx_incons("Unsupported Ewald geometry");
166             break;
167     }
168     if (debug)
169     {
170         fprintf(debug, "dipcorr = %8.3f  %8.3f  %8.3f\n",
171                 dipcorrA[XX], dipcorrA[YY], dipcorrA[ZZ]);
172         fprintf(debug, "mutot   = %8.3f  %8.3f  %8.3f\n",
173                 mutot[0][XX], mutot[0][YY], mutot[0][ZZ]);
174     }
175     bNeedLongRangeCorrection = (calc_excl_corr || dipole_coeff != 0);
176     if (bNeedLongRangeCorrection && !bHaveChargeOrTypePerturbed)
177     {
178         for (i = start; (i < end); i++)
179         {
180             /* Initiate local variables (for this i-particle) to 0 */
181             qiA = chargeA[i]*one_4pi_eps;
182             if (EVDW_PME(fr->vdwtype))
183             {
184                 c6Ai = C6A[i];
185                 if (bDoingLBRule)
186                 {
187                     c6Ai *= sigma3A[i];
188                 }
189             }
190             if (calc_excl_corr)
191             {
192                 i1  = excl->index[i];
193                 i2  = excl->index[i+1];
194
195                 /* Loop over excluded neighbours */
196                 for (j = i1; (j < i2); j++)
197                 {
198                     k = AA[j];
199                     /*
200                      * First we must test whether k <> i, and then,
201                      * because the exclusions are all listed twice i->k
202                      * and k->i we must select just one of the two.  As
203                      * a minor optimization we only compute forces when
204                      * the charges are non-zero.
205                      */
206                     if (k > i)
207                     {
208                         qqA = qiA*chargeA[k];
209                         if (EVDW_PME(fr->vdwtype))
210                         {
211                             c6A  = c6Ai * C6A[k];
212                             if (bDoingLBRule)
213                             {
214                                 c6A *= gmx::power6(0.5*(sigmaA[i]+sigmaA[k]))*sigma3A[k];
215                             }
216                         }
217                         if (qqA != 0.0 || c6A != 0.0)
218                         {
219                             rvec_sub(x[i], x[k], dx);
220                             if (bMolPBC)
221                             {
222                                 /* Cheap pbc_dx, assume excluded pairs are at short distance. */
223                                 for (m = DIM-1; (m >= 0); m--)
224                                 {
225                                     if (dx[m] > 0.5*box[m][m])
226                                     {
227                                         rvec_dec(dx, box[m]);
228                                     }
229                                     else if (dx[m] < -0.5*box[m][m])
230                                     {
231                                         rvec_inc(dx, box[m]);
232                                     }
233                                 }
234                             }
235                             dr2 = norm2(dx);
236                             /* Distance between two excluded particles
237                              * may be zero in the case of shells
238                              */
239                             if (dr2 != 0)
240                             {
241                                 rinv              = gmx::invsqrt(dr2);
242                                 rinv2             = rinv*rinv;
243                                 if (qqA != 0.0)
244                                 {
245                                     real dr, fscal;
246
247                                     dr       = 1.0/rinv;
248                                     ewcdr    = ewc_q*dr;
249                                     vc       = qqA*std::erf(ewcdr)*rinv;
250                                     Vexcl_q += vc;
251 #if GMX_DOUBLE
252                                     /* Relative accuracy at R_ERF_R_INACC of 3e-10 */
253 #define       R_ERF_R_INACC 0.006
254 #else
255                                     /* Relative accuracy at R_ERF_R_INACC of 2e-5 */
256 #define       R_ERF_R_INACC 0.1
257 #endif
258                                     /* fscal is the scalar force pre-multiplied by rinv,
259                                      * to normalise the relative position vector dx */
260                                     if (ewcdr > R_ERF_R_INACC)
261                                     {
262                                         fscal = rinv2*(vc - qqA*ewc_q*M_2_SQRTPI*exp(-ewcdr*ewcdr));
263                                     }
264                                     else
265                                     {
266                                         /* Use a fourth order series expansion for small ewcdr */
267                                         fscal = ewc_q*ewc_q*qqA*vr0_q*(2.0/3.0 - 0.4*ewcdr*ewcdr);
268                                     }
269
270                                     /* The force vector is obtained by multiplication with
271                                      * the relative position vector
272                                      */
273                                     svmul(fscal, dx, df);
274                                     rvec_inc(f[k], df);
275                                     rvec_dec(f[i], df);
276                                     for (iv = 0; (iv < DIM); iv++)
277                                     {
278                                         for (jv = 0; (jv < DIM); jv++)
279                                         {
280                                             dxdf_q[iv][jv] += dx[iv]*df[jv];
281                                         }
282                                     }
283                                 }
284
285                                 if (c6A != 0.0)
286                                 {
287                                     real fscal;
288
289                                     rinv6     = rinv2*rinv2*rinv2;
290                                     ewcdr2    = ewc_lj2*dr2;
291                                     ewcdr4    = ewcdr2*ewcdr2;
292                                     /* We get the excluded long-range contribution from -C6*(1-g(r))
293                                      * g(r) is also defined in the manual under LJ-PME
294                                      */
295                                     vc        = -c6A*rinv6*(1.0 - exp(-ewcdr2)*(1 + ewcdr2 + 0.5*ewcdr4));
296                                     Vexcl_lj += vc;
297                                     /* The force is the derivative of the potential vc.
298                                      * fscal is the scalar force pre-multiplied by rinv,
299                                      * to normalise the relative position vector dx */
300                                     fscal     = 6.0*vc*rinv2 + c6A*rinv6*exp(-ewcdr2)*ewc_lj2*ewcdr4;
301
302                                     /* The force vector is obtained by multiplication with
303                                      * the relative position vector
304                                      */
305                                     svmul(fscal, dx, df);
306                                     rvec_inc(f[k], df);
307                                     rvec_dec(f[i], df);
308                                     for (iv = 0; (iv < DIM); iv++)
309                                     {
310                                         for (jv = 0; (jv < DIM); jv++)
311                                         {
312                                             dxdf_lj[iv][jv] += dx[iv]*df[jv];
313                                         }
314                                     }
315                                 }
316                             }
317                             else
318                             {
319                                 Vexcl_q  += qqA*vr0_q;
320                                 Vexcl_lj += c6A*vr0_lj;
321                             }
322                         }
323                     }
324                 }
325             }
326             /* Dipole correction on force */
327             if (dipole_coeff != 0 && i < numAtomsLocal)
328             {
329                 for (j = 0; (j < DIM); j++)
330                 {
331                     f[i][j] -= dipcorrA[j]*chargeA[i];
332                 }
333             }
334         }
335     }
336     else if (bNeedLongRangeCorrection)
337     {
338         for (i = start; (i < end); i++)
339         {
340             /* Initiate local variables (for this i-particle) to 0 */
341             qiA = chargeA[i]*one_4pi_eps;
342             qiB = chargeB[i]*one_4pi_eps;
343             if (EVDW_PME(fr->vdwtype))
344             {
345                 c6Ai = C6A[i];
346                 c6Bi = C6B[i];
347                 if (bDoingLBRule)
348                 {
349                     c6Ai *= sigma3A[i];
350                     c6Bi *= sigma3B[i];
351                 }
352             }
353             if (calc_excl_corr)
354             {
355                 i1  = excl->index[i];
356                 i2  = excl->index[i+1];
357
358                 /* Loop over excluded neighbours */
359                 for (j = i1; (j < i2); j++)
360                 {
361                     k = AA[j];
362                     if (k > i)
363                     {
364                         qqA = qiA*chargeA[k];
365                         qqB = qiB*chargeB[k];
366                         if (EVDW_PME(fr->vdwtype))
367                         {
368                             c6A = c6Ai*C6A[k];
369                             c6B = c6Bi*C6B[k];
370                             if (bDoingLBRule)
371                             {
372                                 c6A *= gmx::power6(0.5*(sigmaA[i]+sigmaA[k]))*sigma3A[k];
373                                 c6B *= gmx::power6(0.5*(sigmaB[i]+sigmaB[k]))*sigma3B[k];
374                             }
375                         }
376                         if (qqA != 0.0 || qqB != 0.0 || c6A != 0.0 || c6B != 0.0)
377                         {
378                             real fscal;
379
380                             qqL   = L1_q*qqA + lambda_q*qqB;
381                             if (EVDW_PME(fr->vdwtype))
382                             {
383                                 c6L = L1_lj*c6A + lambda_lj*c6B;
384                             }
385                             rvec_sub(x[i], x[k], dx);
386                             if (bMolPBC)
387                             {
388                                 /* Cheap pbc_dx, assume excluded pairs are at short distance. */
389                                 for (m = DIM-1; (m >= 0); m--)
390                                 {
391                                     if (dx[m] > 0.5*box[m][m])
392                                     {
393                                         rvec_dec(dx, box[m]);
394                                     }
395                                     else if (dx[m] < -0.5*box[m][m])
396                                     {
397                                         rvec_inc(dx, box[m]);
398                                     }
399                                 }
400                             }
401                             dr2 = norm2(dx);
402                             if (dr2 != 0)
403                             {
404                                 rinv    = gmx::invsqrt(dr2);
405                                 rinv2   = rinv*rinv;
406                                 if (qqA != 0.0 || qqB != 0.0)
407                                 {
408                                     real dr;
409
410                                     dr           = 1.0/rinv;
411                                     v            = std::erf(ewc_q*dr)*rinv;
412                                     vc           = qqL*v;
413                                     Vexcl_q     += vc;
414                                     /* fscal is the scalar force pre-multiplied by rinv,
415                                      * to normalise the relative position vector dx */
416                                     fscal        = rinv2*(vc-qqL*ewc_q*M_2_SQRTPI*exp(-ewc_q*ewc_q*dr2));
417                                     dvdl_excl_q += (qqB - qqA)*v;
418
419                                     /* The force vector is obtained by multiplication with
420                                      * the relative position vector
421                                      */
422                                     svmul(fscal, dx, df);
423                                     rvec_inc(f[k], df);
424                                     rvec_dec(f[i], df);
425                                     for (iv = 0; (iv < DIM); iv++)
426                                     {
427                                         for (jv = 0; (jv < DIM); jv++)
428                                         {
429                                             dxdf_q[iv][jv] += dx[iv]*df[jv];
430                                         }
431                                     }
432                                 }
433
434                                 if ((c6A != 0.0 || c6B != 0.0) && EVDW_PME(fr->vdwtype))
435                                 {
436                                     rinv6         = rinv2*rinv2*rinv2;
437                                     ewcdr2        = ewc_lj2*dr2;
438                                     ewcdr4        = ewcdr2*ewcdr2;
439                                     v             = -rinv6*(1.0 - exp(-ewcdr2)*(1 + ewcdr2 + 0.5*ewcdr4));
440                                     vc            = c6L*v;
441                                     Vexcl_lj     += vc;
442                                     /* fscal is the scalar force pre-multiplied by rinv,
443                                      * to normalise the relative position vector dx */
444                                     fscal         = 6.0*vc*rinv2 + c6L*rinv6*exp(-ewcdr2)*ewc_lj2*ewcdr4;
445                                     dvdl_excl_lj += (c6B - c6A)*v;
446
447                                     /* The force vector is obtained by multiplication with
448                                      * the relative position vector
449                                      */
450                                     svmul(fscal, dx, df);
451                                     rvec_inc(f[k], df);
452                                     rvec_dec(f[i], df);
453                                     for (iv = 0; (iv < DIM); iv++)
454                                     {
455                                         for (jv = 0; (jv < DIM); jv++)
456                                         {
457                                             dxdf_lj[iv][jv] += dx[iv]*df[jv];
458                                         }
459                                     }
460                                 }
461                             }
462                             else
463                             {
464                                 Vexcl_q      += qqL*vr0_q;
465                                 dvdl_excl_q  += (qqB - qqA)*vr0_q;
466                                 Vexcl_lj     += c6L*vr0_lj;
467                                 dvdl_excl_lj += (c6B - c6A)*vr0_lj;
468                             }
469                         }
470                     }
471                 }
472             }
473             /* Dipole correction on force */
474             if (dipole_coeff != 0 && i < numAtomsLocal)
475             {
476                 for (j = 0; (j < DIM); j++)
477                 {
478                     f[i][j] -= L1_q*dipcorrA[j]*chargeA[i]
479                         + lambda_q*dipcorrB[j]*chargeB[i];
480                 }
481             }
482         }
483     }
484     for (iv = 0; (iv < DIM); iv++)
485     {
486         for (jv = 0; (jv < DIM); jv++)
487         {
488             vir_q[iv][jv]  += 0.5*dxdf_q[iv][jv];
489             vir_lj[iv][jv] += 0.5*dxdf_lj[iv][jv];
490         }
491     }
492
493     Vself_q[0]  = 0;
494     Vself_q[1]  = 0;
495     Vself_lj[0] = 0;
496     Vself_lj[1] = 0;
497
498     /* Global corrections only on master process */
499     if (MASTER(cr) && thread == 0)
500     {
501         for (q = 0; q < (bHaveChargeOrTypePerturbed ? 2 : 1); q++)
502         {
503             if (calc_excl_corr)
504             {
505                 /* Self-energy correction */
506                 Vself_q[q] = ewc_q*one_4pi_eps*fr->q2sum[q]*M_1_SQRTPI;
507                 if (EVDW_PME(fr->vdwtype))
508                 {
509                     Vself_lj[q] =  fr->c6sum[q]*0.5*vr0_lj;
510                 }
511             }
512
513             /* Apply surface dipole correction:
514              * correction = dipole_coeff * (dipole)^2
515              */
516             if (dipole_coeff != 0)
517             {
518                 if (ewald_geometry == eewg3D)
519                 {
520                     Vdipole[q] = dipole_coeff*iprod(mutot[q], mutot[q]);
521                 }
522                 else if (ewald_geometry == eewg3DC)
523                 {
524                     Vdipole[q] = dipole_coeff*mutot[q][ZZ]*mutot[q][ZZ];
525                 }
526             }
527         }
528     }
529     if (!bHaveChargeOrTypePerturbed)
530     {
531         *Vcorr_q = Vdipole[0] - Vself_q[0] - Vexcl_q;
532         if (EVDW_PME(fr->vdwtype))
533         {
534             *Vcorr_lj = -Vself_lj[0] - Vexcl_lj;
535         }
536     }
537     else
538     {
539         *Vcorr_q = L1_q*(Vdipole[0] - Vself_q[0])
540             + lambda_q*(Vdipole[1] - Vself_q[1])
541             - Vexcl_q;
542         *dvdlambda_q += Vdipole[1] - Vself_q[1]
543             - (Vdipole[0] - Vself_q[0]) - dvdl_excl_q;
544         if (EVDW_PME(fr->vdwtype))
545         {
546             *Vcorr_lj      = -(L1_lj*Vself_lj[0] + lambda_lj*Vself_lj[1]) - Vexcl_lj;
547             *dvdlambda_lj += -Vself_lj[1] + Vself_lj[0] - dvdl_excl_lj;
548         }
549     }
550
551     if (debug)
552     {
553         fprintf(debug, "Long Range corrections for Ewald interactions:\n");
554         fprintf(debug, "q2sum = %g, Vself_q=%g c6sum = %g, Vself_lj=%g\n",
555                 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]);
556         fprintf(debug, "Electrostatic Long Range correction: Vexcl=%g\n", Vexcl_q);
557         fprintf(debug, "Lennard-Jones Long Range correction: Vexcl=%g\n", Vexcl_lj);
558         if (MASTER(cr) && thread == 0)
559         {
560             if (epsilon_surface > 0 || ewald_geometry == eewg3DC)
561             {
562                 fprintf(debug, "Total dipole correction: Vdipole=%g\n",
563                         L1_q*Vdipole[0]+lambda_q*Vdipole[1]);
564             }
565         }
566     }
567 }