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