Merge release-5-0 into master
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_free_energy.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 #ifdef HAVE_CONFIG_H
38 #include <config.h>
39 #endif
40
41 #include <math.h>
42
43 #include "gromacs/math/vec.h"
44 #include "typedefs.h"
45 #include "nonbonded.h"
46 #include "nb_kernel.h"
47 #include "nrnb.h"
48 #include "macros.h"
49 #include "nb_free_energy.h"
50
51 #include "gromacs/utility/fatalerror.h"
52
53 void
54 gmx_nb_free_energy_kernel(const t_nblist * gmx_restrict    nlist,
55                           rvec * gmx_restrict              xx,
56                           rvec * gmx_restrict              ff,
57                           t_forcerec * gmx_restrict        fr,
58                           const t_mdatoms * gmx_restrict   mdatoms,
59                           nb_kernel_data_t * gmx_restrict  kernel_data,
60                           t_nrnb * gmx_restrict            nrnb)
61 {
62
63 #define  STATE_A  0
64 #define  STATE_B  1
65 #define  NSTATES  2
66     int           i, j, n, ii, is3, ii3, k, nj0, nj1, jnr, j3, ggid;
67     real          shX, shY, shZ;
68     real          Fscal, FscalC[NSTATES], FscalV[NSTATES], tx, ty, tz;
69     real          Vcoul[NSTATES], Vvdw[NSTATES];
70     real          rinv6, r, rt, rtC, rtV;
71     real          iqA, iqB;
72     real          qq[NSTATES], vctot, krsq;
73     int           ntiA, ntiB, tj[NSTATES];
74     real          Vvdw6, Vvdw12, vvtot;
75     real          ix, iy, iz, fix, fiy, fiz;
76     real          dx, dy, dz, rsq, rinv;
77     real          c6[NSTATES], c12[NSTATES], c6grid;
78     real          LFC[NSTATES], LFV[NSTATES], DLF[NSTATES];
79     double        dvdl_coul, dvdl_vdw;
80     real          lfac_coul[NSTATES], dlfac_coul[NSTATES], lfac_vdw[NSTATES], dlfac_vdw[NSTATES];
81     real          sigma6[NSTATES], alpha_vdw_eff, alpha_coul_eff, sigma2_def, sigma2_min;
82     real          rp, rpm2, rC, rV, rinvC, rpinvC, rinvV, rpinvV;
83     real          sigma2[NSTATES], sigma_pow[NSTATES], sigma_powm2[NSTATES], rs, rs2;
84     int           do_tab, tab_elemsize;
85     int           n0, n1C, n1V, nnn;
86     real          Y, F, G, H, Fp, Geps, Heps2, epsC, eps2C, epsV, eps2V, VV, FF;
87     int           icoul, ivdw;
88     int           nri;
89     const int *   iinr;
90     const int *   jindex;
91     const int *   jjnr;
92     const int *   shift;
93     const int *   gid;
94     const int *   typeA;
95     const int *   typeB;
96     int           ntype;
97     const real *  shiftvec;
98     real          dvdl_part;
99     real *        fshift;
100     real          tabscale = 0;
101     const real *  VFtab    = NULL;
102     const real *  x;
103     real *        f;
104     real          facel, krf, crf;
105     const real *  chargeA;
106     const real *  chargeB;
107     real          sigma6_min, sigma6_def, lam_power, sc_power, sc_r_power;
108     real          alpha_coul, alpha_vdw, lambda_coul, lambda_vdw, ewc_lj;
109     real          ewcljrsq, ewclj, ewclj2, exponent, poly, vvdw_disp, vvdw_rep, sh_lj_ewald;
110     real          ewclj6;
111     const real *  nbfp, *nbfp_grid;
112     real *        dvdl;
113     real *        Vv;
114     real *        Vc;
115     gmx_bool      bDoForces, bDoShiftForces, bDoPotential;
116     real          rcoulomb, rvdw, sh_invrc6;
117     gmx_bool      bExactElecCutoff, bExactVdwCutoff, bExactCutoffAll;
118     gmx_bool      bEwald, bEwaldLJ;
119     real          rcutoff_max2;
120     const real *  tab_ewald_F_lj;
121     const real *  tab_ewald_V_lj;
122     real          d, d2, sw, dsw, rinvcorr;
123     real          elec_swV3, elec_swV4, elec_swV5, elec_swF2, elec_swF3, elec_swF4;
124     real          vdw_swV3, vdw_swV4, vdw_swV5, vdw_swF2, vdw_swF3, vdw_swF4;
125     gmx_bool      bConvertEwaldToCoulomb, bConvertLJEwaldToLJ6;
126     gmx_bool      bComputeVdwInteraction, bComputeElecInteraction;
127     const real *  ewtab;
128     int           ewitab;
129     real          ewrt, eweps, ewtabscale, ewtabhalfspace, sh_ewald;
130
131     sh_ewald            = fr->ic->sh_ewald;
132     ewtab               = fr->ic->tabq_coul_FDV0;
133     ewtabscale          = fr->ic->tabq_scale;
134     ewtabhalfspace      = 0.5/ewtabscale;
135     tab_ewald_F_lj      = fr->ic->tabq_vdw_F;
136     tab_ewald_V_lj      = fr->ic->tabq_vdw_V;
137
138     x                   = xx[0];
139     f                   = ff[0];
140
141     fshift              = fr->fshift[0];
142
143     nri                 = nlist->nri;
144     iinr                = nlist->iinr;
145     jindex              = nlist->jindex;
146     jjnr                = nlist->jjnr;
147     icoul               = nlist->ielec;
148     ivdw                = nlist->ivdw;
149     shift               = nlist->shift;
150     gid                 = nlist->gid;
151
152     shiftvec            = fr->shift_vec[0];
153     chargeA             = mdatoms->chargeA;
154     chargeB             = mdatoms->chargeB;
155     facel               = fr->epsfac;
156     krf                 = fr->k_rf;
157     crf                 = fr->c_rf;
158     ewc_lj              = fr->ewaldcoeff_lj;
159     Vc                  = kernel_data->energygrp_elec;
160     typeA               = mdatoms->typeA;
161     typeB               = mdatoms->typeB;
162     ntype               = fr->ntype;
163     nbfp                = fr->nbfp;
164     nbfp_grid           = fr->ljpme_c6grid;
165     Vv                  = kernel_data->energygrp_vdw;
166     lambda_coul         = kernel_data->lambda[efptCOUL];
167     lambda_vdw          = kernel_data->lambda[efptVDW];
168     dvdl                = kernel_data->dvdl;
169     alpha_coul          = fr->sc_alphacoul;
170     alpha_vdw           = fr->sc_alphavdw;
171     lam_power           = fr->sc_power;
172     sc_r_power          = fr->sc_r_power;
173     sigma6_def          = fr->sc_sigma6_def;
174     sigma6_min          = fr->sc_sigma6_min;
175     bDoForces           = kernel_data->flags & GMX_NONBONDED_DO_FORCE;
176     bDoShiftForces      = kernel_data->flags & GMX_NONBONDED_DO_SHIFTFORCE;
177     bDoPotential        = kernel_data->flags & GMX_NONBONDED_DO_POTENTIAL;
178
179     rcoulomb            = fr->rcoulomb;
180     rvdw                = fr->rvdw;
181     sh_invrc6           = fr->ic->sh_invrc6;
182     sh_lj_ewald         = fr->ic->sh_lj_ewald;
183     ewclj               = fr->ewaldcoeff_lj;
184     ewclj2              = ewclj*ewclj;
185     ewclj6              = ewclj2*ewclj2*ewclj2;
186
187     if (fr->coulomb_modifier == eintmodPOTSWITCH)
188     {
189         d               = fr->rcoulomb-fr->rcoulomb_switch;
190         elec_swV3       = -10.0/(d*d*d);
191         elec_swV4       =  15.0/(d*d*d*d);
192         elec_swV5       =  -6.0/(d*d*d*d*d);
193         elec_swF2       = -30.0/(d*d*d);
194         elec_swF3       =  60.0/(d*d*d*d);
195         elec_swF4       = -30.0/(d*d*d*d*d);
196     }
197     else
198     {
199         /* Avoid warnings from stupid compilers (looking at you, Clang!) */
200         elec_swV3 = elec_swV4 = elec_swV5 = elec_swF2 = elec_swF3 = elec_swF4 = 0.0;
201     }
202
203     if (fr->vdw_modifier == eintmodPOTSWITCH)
204     {
205         d               = fr->rvdw-fr->rvdw_switch;
206         vdw_swV3        = -10.0/(d*d*d);
207         vdw_swV4        =  15.0/(d*d*d*d);
208         vdw_swV5        =  -6.0/(d*d*d*d*d);
209         vdw_swF2        = -30.0/(d*d*d);
210         vdw_swF3        =  60.0/(d*d*d*d);
211         vdw_swF4        = -30.0/(d*d*d*d*d);
212     }
213     else
214     {
215         /* Avoid warnings from stupid compilers (looking at you, Clang!) */
216         vdw_swV3 = vdw_swV4 = vdw_swV5 = vdw_swF2 = vdw_swF3 = vdw_swF4 = 0.0;
217     }
218
219     if (fr->cutoff_scheme == ecutsVERLET)
220     {
221         const interaction_const_t *ic;
222
223         ic = fr->ic;
224         if (EVDW_PME(ic->vdwtype))
225         {
226             ivdw         = GMX_NBKERNEL_VDW_LJEWALD;
227         }
228         else
229         {
230             ivdw         = GMX_NBKERNEL_VDW_LENNARDJONES;
231         }
232
233         if (ic->eeltype == eelCUT || EEL_RF(ic->eeltype))
234         {
235             icoul        = GMX_NBKERNEL_ELEC_REACTIONFIELD;
236         }
237         else if (EEL_PME_EWALD(ic->eeltype))
238         {
239             icoul        = GMX_NBKERNEL_ELEC_EWALD;
240         }
241         else
242         {
243             gmx_incons("Unsupported eeltype with Verlet and free-energy");
244         }
245
246         bExactElecCutoff = TRUE;
247         bExactVdwCutoff  = TRUE;
248     }
249     else
250     {
251         bExactElecCutoff = (fr->coulomb_modifier != eintmodNONE) || fr->eeltype == eelRF_ZERO;
252         bExactVdwCutoff  = (fr->vdw_modifier != eintmodNONE);
253     }
254
255     bExactCutoffAll = (bExactElecCutoff && bExactVdwCutoff);
256     rcutoff_max2    = max(fr->rcoulomb, fr->rvdw);
257     rcutoff_max2    = rcutoff_max2*rcutoff_max2;
258
259     bEwald          = (icoul == GMX_NBKERNEL_ELEC_EWALD);
260     bEwaldLJ        = (ivdw == GMX_NBKERNEL_VDW_LJEWALD);
261
262     /* For Ewald/PME interactions we cannot easily apply the soft-core component to
263      * reciprocal space. When we use vanilla (not switch/shift) Ewald interactions, we
264      * can apply the small trick of subtracting the _reciprocal_ space contribution
265      * in this kernel, and instead apply the free energy interaction to the 1/r
266      * (standard coulomb) interaction.
267      *
268      * However, we cannot use this approach for switch-modified since we would then
269      * effectively end up evaluating a significantly different interaction here compared to the
270      * normal (non-free-energy) kernels, either by applying a cutoff at a different
271      * position than what the user requested, or by switching different
272      * things (1/r rather than short-range Ewald). For these settings, we just
273      * use the traditional short-range Ewald interaction in that case.
274      */
275     bConvertEwaldToCoulomb = (bEwald && (fr->coulomb_modifier != eintmodPOTSWITCH));
276     /* For now the below will always be true (since LJ-PME only works with Shift in Gromacs-5.0),
277      * but writing it this way means we stay in sync with coulomb, and it avoids future bugs.
278      */
279     bConvertLJEwaldToLJ6   = (bEwaldLJ && (fr->vdw_modifier   != eintmodPOTSWITCH));
280
281     /* We currently don't implement exclusion correction, needed with the Verlet cut-off scheme, without conversion */
282     if (fr->cutoff_scheme == ecutsVERLET &&
283         ((bEwald   && !bConvertEwaldToCoulomb) ||
284          (bEwaldLJ && !bConvertLJEwaldToLJ6)))
285     {
286         gmx_incons("Unimplemented non-bonded setup");
287     }
288
289     /* fix compiler warnings */
290     nj1   = 0;
291     n1C   = n1V   = 0;
292     epsC  = epsV  = 0;
293     eps2C = eps2V = 0;
294
295     dvdl_coul  = 0;
296     dvdl_vdw   = 0;
297
298     /* Lambda factor for state A, 1-lambda*/
299     LFC[STATE_A] = 1.0 - lambda_coul;
300     LFV[STATE_A] = 1.0 - lambda_vdw;
301
302     /* Lambda factor for state B, lambda*/
303     LFC[STATE_B] = lambda_coul;
304     LFV[STATE_B] = lambda_vdw;
305
306     /*derivative of the lambda factor for state A and B */
307     DLF[STATE_A] = -1;
308     DLF[STATE_B] = 1;
309
310     for (i = 0; i < NSTATES; i++)
311     {
312         lfac_coul[i]  = (lam_power == 2 ? (1-LFC[i])*(1-LFC[i]) : (1-LFC[i]));
313         dlfac_coul[i] = DLF[i]*lam_power/sc_r_power*(lam_power == 2 ? (1-LFC[i]) : 1);
314         lfac_vdw[i]   = (lam_power == 2 ? (1-LFV[i])*(1-LFV[i]) : (1-LFV[i]));
315         dlfac_vdw[i]  = DLF[i]*lam_power/sc_r_power*(lam_power == 2 ? (1-LFV[i]) : 1);
316     }
317     /* precalculate */
318     sigma2_def = pow(sigma6_def, 1.0/3.0);
319     sigma2_min = pow(sigma6_min, 1.0/3.0);
320
321     /* Ewald (not PME) table is special (icoul==enbcoulFEWALD) */
322
323     do_tab = (icoul == GMX_NBKERNEL_ELEC_CUBICSPLINETABLE ||
324               ivdw == GMX_NBKERNEL_VDW_CUBICSPLINETABLE);
325     if (do_tab)
326     {
327         tabscale         = kernel_data->table_elec_vdw->scale;
328         VFtab            = kernel_data->table_elec_vdw->data;
329         /* we always use the combined table here */
330         tab_elemsize     = 12;
331     }
332
333     for (n = 0; (n < nri); n++)
334     {
335         int npair_within_cutoff;
336
337         npair_within_cutoff = 0;
338
339         is3              = 3*shift[n];
340         shX              = shiftvec[is3];
341         shY              = shiftvec[is3+1];
342         shZ              = shiftvec[is3+2];
343         nj0              = jindex[n];
344         nj1              = jindex[n+1];
345         ii               = iinr[n];
346         ii3              = 3*ii;
347         ix               = shX + x[ii3+0];
348         iy               = shY + x[ii3+1];
349         iz               = shZ + x[ii3+2];
350         iqA              = facel*chargeA[ii];
351         iqB              = facel*chargeB[ii];
352         ntiA             = 2*ntype*typeA[ii];
353         ntiB             = 2*ntype*typeB[ii];
354         vctot            = 0;
355         vvtot            = 0;
356         fix              = 0;
357         fiy              = 0;
358         fiz              = 0;
359
360         for (k = nj0; (k < nj1); k++)
361         {
362             jnr              = jjnr[k];
363             j3               = 3*jnr;
364             dx               = ix - x[j3];
365             dy               = iy - x[j3+1];
366             dz               = iz - x[j3+2];
367             rsq              = dx*dx + dy*dy + dz*dz;
368
369             if (bExactCutoffAll && rsq >= rcutoff_max2)
370             {
371                 /* We save significant time by skipping all code below.
372                  * Note that with soft-core interactions, the actual cut-off
373                  * check might be different. But since the soft-core distance
374                  * is always larger than r, checking on r here is safe.
375                  */
376                 continue;
377             }
378             npair_within_cutoff++;
379
380             if (rsq > 0)
381             {
382                 rinv         = gmx_invsqrt(rsq);
383                 r            = rsq*rinv;
384             }
385             else
386             {
387                 /* The force at r=0 is zero, because of symmetry.
388                  * But note that the potential is in general non-zero,
389                  * since the soft-cored r will be non-zero.
390                  */
391                 rinv         = 0;
392                 r            = 0;
393             }
394
395             if (sc_r_power == 6.0)
396             {
397                 rpm2             = rsq*rsq;  /* r4 */
398                 rp               = rpm2*rsq; /* r6 */
399             }
400             else if (sc_r_power == 48.0)
401             {
402                 rp               = rsq*rsq*rsq; /* r6 */
403                 rp               = rp*rp;       /* r12 */
404                 rp               = rp*rp;       /* r24 */
405                 rp               = rp*rp;       /* r48 */
406                 rpm2             = rp/rsq;      /* r46 */
407             }
408             else
409             {
410                 rp             = pow(r, sc_r_power);  /* not currently supported as input, but can handle it */
411                 rpm2           = rp/rsq;
412             }
413
414             Fscal = 0;
415
416             qq[STATE_A]      = iqA*chargeA[jnr];
417             qq[STATE_B]      = iqB*chargeB[jnr];
418
419             tj[STATE_A]      = ntiA+2*typeA[jnr];
420             tj[STATE_B]      = ntiB+2*typeB[jnr];
421
422             if (nlist->excl_fep == NULL || nlist->excl_fep[k])
423             {
424                 c6[STATE_A]      = nbfp[tj[STATE_A]];
425                 c6[STATE_B]      = nbfp[tj[STATE_B]];
426
427                 for (i = 0; i < NSTATES; i++)
428                 {
429                     c12[i]             = nbfp[tj[i]+1];
430                     if ((c6[i] > 0) && (c12[i] > 0))
431                     {
432                         /* c12 is stored scaled with 12.0 and c6 is scaled with 6.0 - correct for this */
433                         sigma6[i]       = 0.5*c12[i]/c6[i];
434                         sigma2[i]       = pow(sigma6[i], 1.0/3.0);
435                         /* should be able to get rid of this ^^^ internal pow call eventually.  Will require agreement on
436                            what data to store externally.  Can't be fixed without larger scale changes, so not 4.6 */
437                         if (sigma6[i] < sigma6_min)   /* for disappearing coul and vdw with soft core at the same time */
438                         {
439                             sigma6[i] = sigma6_min;
440                             sigma2[i] = sigma2_min;
441                         }
442                     }
443                     else
444                     {
445                         sigma6[i]       = sigma6_def;
446                         sigma2[i]       = sigma2_def;
447                     }
448                     if (sc_r_power == 6.0)
449                     {
450                         sigma_pow[i]    = sigma6[i];
451                         sigma_powm2[i]  = sigma6[i]/sigma2[i];
452                     }
453                     else if (sc_r_power == 48.0)
454                     {
455                         sigma_pow[i]    = sigma6[i]*sigma6[i];       /* sigma^12 */
456                         sigma_pow[i]    = sigma_pow[i]*sigma_pow[i]; /* sigma^24 */
457                         sigma_pow[i]    = sigma_pow[i]*sigma_pow[i]; /* sigma^48 */
458                         sigma_powm2[i]  = sigma_pow[i]/sigma2[i];
459                     }
460                     else
461                     {    /* not really supported as input, but in here for testing the general case*/
462                         sigma_pow[i]    = pow(sigma2[i], sc_r_power/2);
463                         sigma_powm2[i]  = sigma_pow[i]/(sigma2[i]);
464                     }
465                 }
466
467                 /* only use softcore if one of the states has a zero endstate - softcore is for avoiding infinities!*/
468                 if ((c12[STATE_A] > 0) && (c12[STATE_B] > 0))
469                 {
470                     alpha_vdw_eff    = 0;
471                     alpha_coul_eff   = 0;
472                 }
473                 else
474                 {
475                     alpha_vdw_eff    = alpha_vdw;
476                     alpha_coul_eff   = alpha_coul;
477                 }
478
479                 for (i = 0; i < NSTATES; i++)
480                 {
481                     FscalC[i]    = 0;
482                     FscalV[i]    = 0;
483                     Vcoul[i]     = 0;
484                     Vvdw[i]      = 0;
485
486                     /* Only spend time on A or B state if it is non-zero */
487                     if ( (qq[i] != 0) || (c6[i] != 0) || (c12[i] != 0) )
488                     {
489                         /* this section has to be inside the loop because of the dependence on sigma_pow */
490                         rpinvC         = 1.0/(alpha_coul_eff*lfac_coul[i]*sigma_pow[i]+rp);
491                         rinvC          = pow(rpinvC, 1.0/sc_r_power);
492                         rC             = 1.0/rinvC;
493
494                         rpinvV         = 1.0/(alpha_vdw_eff*lfac_vdw[i]*sigma_pow[i]+rp);
495                         rinvV          = pow(rpinvV, 1.0/sc_r_power);
496                         rV             = 1.0/rinvV;
497
498                         if (do_tab)
499                         {
500                             rtC        = rC*tabscale;
501                             n0         = rtC;
502                             epsC       = rtC-n0;
503                             eps2C      = epsC*epsC;
504                             n1C        = tab_elemsize*n0;
505
506                             rtV        = rV*tabscale;
507                             n0         = rtV;
508                             epsV       = rtV-n0;
509                             eps2V      = epsV*epsV;
510                             n1V        = tab_elemsize*n0;
511                         }
512
513                         /* Only process the coulomb interactions if we have charges,
514                          * and if we either include all entries in the list (no cutoff
515                          * used in the kernel), or if we are within the cutoff.
516                          */
517                         bComputeElecInteraction = !bExactElecCutoff ||
518                             ( bConvertEwaldToCoulomb && r < rcoulomb) ||
519                             (!bConvertEwaldToCoulomb && rC < rcoulomb);
520
521                         if ( (qq[i] != 0) && bComputeElecInteraction)
522                         {
523                             switch (icoul)
524                             {
525                                 case GMX_NBKERNEL_ELEC_COULOMB:
526                                     /* simple cutoff */
527                                     Vcoul[i]   = qq[i]*rinvC;
528                                     FscalC[i]  = Vcoul[i];
529                                     /* The shift for the Coulomb potential is stored in
530                                      * the RF parameter c_rf, which is 0 without shift.
531                                      */
532                                     Vcoul[i]  -= qq[i]*fr->ic->c_rf;
533                                     break;
534
535                                 case GMX_NBKERNEL_ELEC_REACTIONFIELD:
536                                     /* reaction-field */
537                                     Vcoul[i]   = qq[i]*(rinvC + krf*rC*rC-crf);
538                                     FscalC[i]  = qq[i]*(rinvC - 2.0*krf*rC*rC);
539                                     break;
540
541                                 case GMX_NBKERNEL_ELEC_CUBICSPLINETABLE:
542                                     /* non-Ewald tabulated coulomb */
543                                     nnn        = n1C;
544                                     Y          = VFtab[nnn];
545                                     F          = VFtab[nnn+1];
546                                     Geps       = epsC*VFtab[nnn+2];
547                                     Heps2      = eps2C*VFtab[nnn+3];
548                                     Fp         = F+Geps+Heps2;
549                                     VV         = Y+epsC*Fp;
550                                     FF         = Fp+Geps+2.0*Heps2;
551                                     Vcoul[i]   = qq[i]*VV;
552                                     FscalC[i]  = -qq[i]*tabscale*FF*rC;
553                                     break;
554
555                                 case GMX_NBKERNEL_ELEC_GENERALIZEDBORN:
556                                     gmx_fatal(FARGS, "Free energy and GB not implemented.\n");
557                                     break;
558
559                                 case GMX_NBKERNEL_ELEC_EWALD:
560                                     if (bConvertEwaldToCoulomb)
561                                     {
562                                         /* Ewald FEP is done only on the 1/r part */
563                                         Vcoul[i]   = qq[i]*(rinvC-sh_ewald);
564                                         FscalC[i]  = qq[i]*rinvC;
565                                     }
566                                     else
567                                     {
568                                         ewrt      = rC*ewtabscale;
569                                         ewitab    = (int) ewrt;
570                                         eweps     = ewrt-ewitab;
571                                         ewitab    = 4*ewitab;
572                                         FscalC[i] = ewtab[ewitab]+eweps*ewtab[ewitab+1];
573                                         rinvcorr  = rinvC-sh_ewald;
574                                         Vcoul[i]  = qq[i]*(rinvcorr-(ewtab[ewitab+2]-ewtabhalfspace*eweps*(ewtab[ewitab]+FscalC[i])));
575                                         FscalC[i] = qq[i]*(rinvC-rC*FscalC[i]);
576                                     }
577                                     break;
578
579                                 case GMX_NBKERNEL_ELEC_NONE:
580                                     FscalC[i]  = 0.0;
581                                     Vcoul[i]   = 0.0;
582                                     break;
583
584                                 default:
585                                     gmx_incons("Invalid icoul in free energy kernel");
586                                     break;
587                             }
588
589                             if (fr->coulomb_modifier == eintmodPOTSWITCH)
590                             {
591                                 d                = rC-fr->rcoulomb_switch;
592                                 d                = (d > 0.0) ? d : 0.0;
593                                 d2               = d*d;
594                                 sw               = 1.0+d2*d*(elec_swV3+d*(elec_swV4+d*elec_swV5));
595                                 dsw              = d2*(elec_swF2+d*(elec_swF3+d*elec_swF4));
596
597                                 FscalC[i]        = FscalC[i]*sw - rC*Vcoul[i]*dsw;
598                                 Vcoul[i]        *= sw;
599
600                                 FscalC[i]        = (rC < rcoulomb) ? FscalC[i] : 0.0;
601                                 Vcoul[i]         = (rC < rcoulomb) ? Vcoul[i] : 0.0;
602                             }
603                         }
604
605                         /* Only process the VDW interactions if we have
606                          * some non-zero parameters, and if we either
607                          * include all entries in the list (no cutoff used
608                          * in the kernel), or if we are within the cutoff.
609                          */
610                         bComputeVdwInteraction = !bExactVdwCutoff ||
611                             ( bConvertLJEwaldToLJ6 && r < rvdw) ||
612                             (!bConvertLJEwaldToLJ6 && rV < rvdw);
613                         if ((c6[i] != 0 || c12[i] != 0) && bComputeVdwInteraction)
614                         {
615                             switch (ivdw)
616                             {
617                                 case GMX_NBKERNEL_VDW_LENNARDJONES:
618                                     /* cutoff LJ */
619                                     if (sc_r_power == 6.0)
620                                     {
621                                         rinv6            = rpinvV;
622                                     }
623                                     else
624                                     {
625                                         rinv6            = rinvV*rinvV;
626                                         rinv6            = rinv6*rinv6*rinv6;
627                                     }
628                                     Vvdw6            = c6[i]*rinv6;
629                                     Vvdw12           = c12[i]*rinv6*rinv6;
630
631                                     Vvdw[i]          = ( (Vvdw12 - c12[i]*sh_invrc6*sh_invrc6)*(1.0/12.0)
632                                                          - (Vvdw6 - c6[i]*sh_invrc6)*(1.0/6.0));
633                                     FscalV[i]        = Vvdw12 - Vvdw6;
634                                     break;
635
636                                 case GMX_NBKERNEL_VDW_BUCKINGHAM:
637                                     gmx_fatal(FARGS, "Buckingham free energy not supported.");
638                                     break;
639
640                                 case GMX_NBKERNEL_VDW_CUBICSPLINETABLE:
641                                     /* Table LJ */
642                                     nnn = n1V+4;
643                                     /* dispersion */
644                                     Y          = VFtab[nnn];
645                                     F          = VFtab[nnn+1];
646                                     Geps       = epsV*VFtab[nnn+2];
647                                     Heps2      = eps2V*VFtab[nnn+3];
648                                     Fp         = F+Geps+Heps2;
649                                     VV         = Y+epsV*Fp;
650                                     FF         = Fp+Geps+2.0*Heps2;
651                                     Vvdw[i]   += c6[i]*VV;
652                                     FscalV[i] -= c6[i]*tabscale*FF*rV;
653
654                                     /* repulsion */
655                                     Y          = VFtab[nnn+4];
656                                     F          = VFtab[nnn+5];
657                                     Geps       = epsV*VFtab[nnn+6];
658                                     Heps2      = eps2V*VFtab[nnn+7];
659                                     Fp         = F+Geps+Heps2;
660                                     VV         = Y+epsV*Fp;
661                                     FF         = Fp+Geps+2.0*Heps2;
662                                     Vvdw[i]   += c12[i]*VV;
663                                     FscalV[i] -= c12[i]*tabscale*FF*rV;
664                                     break;
665
666                                 case GMX_NBKERNEL_VDW_LJEWALD:
667                                     if (sc_r_power == 6.0)
668                                     {
669                                         rinv6            = rpinvV;
670                                     }
671                                     else
672                                     {
673                                         rinv6            = rinvV*rinvV;
674                                         rinv6            = rinv6*rinv6*rinv6;
675                                     }
676                                     c6grid           = nbfp_grid[tj[i]];
677
678                                     if (bConvertLJEwaldToLJ6)
679                                     {
680                                         /* cutoff LJ */
681                                         Vvdw6            = c6[i]*rinv6;
682                                         Vvdw12           = c12[i]*rinv6*rinv6;
683
684                                         Vvdw[i]          = ( (Vvdw12 - c12[i]*sh_invrc6*sh_invrc6)*(1.0/12.0)
685                                                              - (Vvdw6 - c6[i]*sh_invrc6 - c6grid*sh_lj_ewald)*(1.0/6.0));
686                                         FscalV[i]        = Vvdw12 - Vvdw6;
687                                     }
688                                     else
689                                     {
690                                         /* Normal LJ-PME */
691                                         ewcljrsq         = ewclj2*rV*rV;
692                                         exponent         = exp(-ewcljrsq);
693                                         poly             = exponent*(1.0 + ewcljrsq + ewcljrsq*ewcljrsq*0.5);
694                                         vvdw_disp        = (c6[i]-c6grid*(1.0-poly))*rinv6;
695                                         vvdw_rep         = c12[i]*rinv6*rinv6;
696                                         FscalV[i]        = vvdw_rep - vvdw_disp - c6grid*(1.0/6.0)*exponent*ewclj6;
697                                         Vvdw[i]          = (vvdw_rep - c12[i]*sh_invrc6*sh_invrc6)/12.0 - (vvdw_disp - c6[i]*sh_invrc6 - c6grid*sh_lj_ewald)/6.0;
698                                     }
699                                     break;
700
701                                 case GMX_NBKERNEL_VDW_NONE:
702                                     Vvdw[i]    = 0.0;
703                                     FscalV[i]  = 0.0;
704                                     break;
705
706                                 default:
707                                     gmx_incons("Invalid ivdw in free energy kernel");
708                                     break;
709                             }
710
711                             if (fr->vdw_modifier == eintmodPOTSWITCH)
712                             {
713                                 d                = rV-fr->rvdw_switch;
714                                 d                = (d > 0.0) ? d : 0.0;
715                                 d2               = d*d;
716                                 sw               = 1.0+d2*d*(vdw_swV3+d*(vdw_swV4+d*vdw_swV5));
717                                 dsw              = d2*(vdw_swF2+d*(vdw_swF3+d*vdw_swF4));
718
719                                 FscalV[i]        = FscalV[i]*sw - rV*Vvdw[i]*dsw;
720                                 Vvdw[i]         *= sw;
721
722                                 FscalV[i]  = (rV < rvdw) ? FscalV[i] : 0.0;
723                                 Vvdw[i]    = (rV < rvdw) ? Vvdw[i] : 0.0;
724                             }
725                         }
726
727                         /* FscalC (and FscalV) now contain: dV/drC * rC
728                          * Now we multiply by rC^-p, so it will be: dV/drC * rC^1-p
729                          * Further down we first multiply by r^p-2 and then by
730                          * the vector r, which in total gives: dV/drC * (r/rC)^1-p
731                          */
732                         FscalC[i] *= rpinvC;
733                         FscalV[i] *= rpinvV;
734                     }
735                 }
736
737                 /* Assemble A and B states */
738                 for (i = 0; i < NSTATES; i++)
739                 {
740                     vctot         += LFC[i]*Vcoul[i];
741                     vvtot         += LFV[i]*Vvdw[i];
742
743                     Fscal         += LFC[i]*FscalC[i]*rpm2;
744                     Fscal         += LFV[i]*FscalV[i]*rpm2;
745
746                     dvdl_coul     += Vcoul[i]*DLF[i] + LFC[i]*alpha_coul_eff*dlfac_coul[i]*FscalC[i]*sigma_pow[i];
747                     dvdl_vdw      += Vvdw[i]*DLF[i] + LFV[i]*alpha_vdw_eff*dlfac_vdw[i]*FscalV[i]*sigma_pow[i];
748                 }
749             }
750             else if (icoul == GMX_NBKERNEL_ELEC_REACTIONFIELD)
751             {
752                 /* For excluded pairs, which are only in this pair list when
753                  * using the Verlet scheme, we don't use soft-core.
754                  * The group scheme also doesn't soft-core for these.
755                  * As there is no singularity, there is no need for soft-core.
756                  */
757                 VV = krf*rsq - crf;
758                 FF = -2.0*krf;
759
760                 if (ii == jnr)
761                 {
762                     VV *= 0.5;
763                 }
764
765                 for (i = 0; i < NSTATES; i++)
766                 {
767                     vctot      += LFC[i]*qq[i]*VV;
768                     Fscal      += LFC[i]*qq[i]*FF;
769                     dvdl_coul  += DLF[i]*qq[i]*VV;
770                 }
771             }
772
773             if (bConvertEwaldToCoulomb && ( !bExactElecCutoff || r < rcoulomb ) )
774             {
775                 /* See comment in the preamble. When using Ewald interactions
776                  * (unless we use a switch modifier) we subtract the reciprocal-space
777                  * Ewald component here which made it possible to apply the free
778                  * energy interaction to 1/r (vanilla coulomb short-range part)
779                  * above. This gets us closer to the ideal case of applying
780                  * the softcore to the entire electrostatic interaction,
781                  * including the reciprocal-space component.
782                  */
783                 real v_lr, f_lr;
784
785                 ewrt      = r*ewtabscale;
786                 ewitab    = (int) ewrt;
787                 eweps     = ewrt-ewitab;
788                 ewitab    = 4*ewitab;
789                 f_lr      = ewtab[ewitab]+eweps*ewtab[ewitab+1];
790                 v_lr      = (ewtab[ewitab+2]-ewtabhalfspace*eweps*(ewtab[ewitab]+f_lr));
791                 f_lr     *= rinv;
792
793                 /* Note that any possible Ewald shift has already been applied in
794                  * the normal interaction part above.
795                  */
796
797                 if (ii == jnr)
798                 {
799                     /* If we get here, the i particle (ii) has itself (jnr)
800                      * in its neighborlist. This can only happen with the Verlet
801                      * scheme, and corresponds to a self-interaction that will
802                      * occur twice. Scale it down by 50% to only include it once.
803                      */
804                     v_lr *= 0.5;
805                 }
806
807                 for (i = 0; i < NSTATES; i++)
808                 {
809                     vctot      -= LFC[i]*qq[i]*v_lr;
810                     Fscal      -= LFC[i]*qq[i]*f_lr;
811                     dvdl_coul  -= (DLF[i]*qq[i])*v_lr;
812                 }
813             }
814
815             if (bConvertLJEwaldToLJ6 && (!bExactVdwCutoff || r < rvdw))
816             {
817                 /* See comment in the preamble. When using LJ-Ewald interactions
818                  * (unless we use a switch modifier) we subtract the reciprocal-space
819                  * Ewald component here which made it possible to apply the free
820                  * energy interaction to r^-6 (vanilla LJ6 short-range part)
821                  * above. This gets us closer to the ideal case of applying
822                  * the softcore to the entire VdW interaction,
823                  * including the reciprocal-space component.
824                  */
825                 /* We could also use the analytical form here
826                  * iso a table, but that can cause issues for
827                  * r close to 0 for non-interacting pairs.
828                  */
829                 real rs, frac, f_lr;
830                 int  ri;
831
832                 rs     = rsq*rinv*ewtabscale;
833                 ri     = (int)rs;
834                 frac   = rs - ri;
835                 f_lr   = (1 - frac)*tab_ewald_F_lj[ri] + frac*tab_ewald_F_lj[ri+1];
836                 /* TODO: Currently the Ewald LJ table does not contain
837                  * the factor 1/6, we should add this.
838                  */
839                 FF     = f_lr*rinv/6.0;
840                 VV     = (tab_ewald_V_lj[ri] - ewtabhalfspace*frac*(tab_ewald_F_lj[ri] + f_lr))/6.0;
841
842                 if (ii == jnr)
843                 {
844                     /* If we get here, the i particle (ii) has itself (jnr)
845                      * in its neighborlist. This can only happen with the Verlet
846                      * scheme, and corresponds to a self-interaction that will
847                      * occur twice. Scale it down by 50% to only include it once.
848                      */
849                     VV *= 0.5;
850                 }
851
852                 for (i = 0; i < NSTATES; i++)
853                 {
854                     c6grid      = nbfp_grid[tj[i]];
855                     vvtot      += LFV[i]*c6grid*VV;
856                     Fscal      += LFV[i]*c6grid*FF;
857                     dvdl_vdw   += (DLF[i]*c6grid)*VV;
858                 }
859             }
860
861             if (bDoForces)
862             {
863                 tx         = Fscal*dx;
864                 ty         = Fscal*dy;
865                 tz         = Fscal*dz;
866                 fix        = fix + tx;
867                 fiy        = fiy + ty;
868                 fiz        = fiz + tz;
869                 /* OpenMP atomics are expensive, but this kernels is also
870                  * expensive, so we can take this hit, instead of using
871                  * thread-local output buffers and extra reduction.
872                  */
873 #pragma omp atomic
874                 f[j3]     -= tx;
875 #pragma omp atomic
876                 f[j3+1]   -= ty;
877 #pragma omp atomic
878                 f[j3+2]   -= tz;
879             }
880         }
881
882         /* The atomics below are expensive with many OpenMP threads.
883          * Here unperturbed i-particles will usually only have a few
884          * (perturbed) j-particles in the list. Thus with a buffered list
885          * we can skip a significant number of i-reductions with a check.
886          */
887         if (npair_within_cutoff > 0)
888         {
889             if (bDoForces)
890             {
891 #pragma omp atomic
892                 f[ii3]        += fix;
893 #pragma omp atomic
894                 f[ii3+1]      += fiy;
895 #pragma omp atomic
896                 f[ii3+2]      += fiz;
897             }
898             if (bDoShiftForces)
899             {
900 #pragma omp atomic
901                 fshift[is3]   += fix;
902 #pragma omp atomic
903                 fshift[is3+1] += fiy;
904 #pragma omp atomic
905                 fshift[is3+2] += fiz;
906             }
907             if (bDoPotential)
908             {
909                 ggid               = gid[n];
910 #pragma omp atomic
911                 Vc[ggid]          += vctot;
912 #pragma omp atomic
913                 Vv[ggid]          += vvtot;
914             }
915         }
916     }
917
918 #pragma omp atomic
919     dvdl[efptCOUL]     += dvdl_coul;
920  #pragma omp atomic
921     dvdl[efptVDW]      += dvdl_vdw;
922
923     /* Estimate flops, average for free energy stuff:
924      * 12  flops per outer iteration
925      * 150 flops per inner iteration
926      */
927 #pragma omp atomic
928     inc_nrnb(nrnb, eNR_NBKERNEL_FREE_ENERGY, nlist->nri*12 + nlist->jindex[n]*150);
929 }
930
931 real
932 nb_free_energy_evaluate_single(real r2, real sc_r_power, real alpha_coul, real alpha_vdw,
933                                real tabscale, real *vftab,
934                                real qqA, real c6A, real c12A, real qqB, real c6B, real c12B,
935                                real LFC[2], real LFV[2], real DLF[2],
936                                real lfac_coul[2], real lfac_vdw[2], real dlfac_coul[2], real dlfac_vdw[2],
937                                real sigma6_def, real sigma6_min, real sigma2_def, real sigma2_min,
938                                real *velectot, real *vvdwtot, real *dvdl)
939 {
940     real       r, rp, rpm2, rtab, eps, eps2, Y, F, Geps, Heps2, Fp, VV, FF, fscal;
941     real       qq[2], c6[2], c12[2], sigma6[2], sigma2[2], sigma_pow[2], sigma_powm2[2];
942     real       alpha_coul_eff, alpha_vdw_eff, dvdl_coul, dvdl_vdw;
943     real       rpinv, r_coul, r_vdw, velecsum, vvdwsum;
944     real       fscal_vdw[2], fscal_elec[2];
945     real       velec[2], vvdw[2];
946     int        i, ntab;
947
948     qq[0]    = qqA;
949     qq[1]    = qqB;
950     c6[0]    = c6A;
951     c6[1]    = c6B;
952     c12[0]   = c12A;
953     c12[1]   = c12B;
954
955     if (sc_r_power == 6.0)
956     {
957         rpm2             = r2*r2;   /* r4 */
958         rp               = rpm2*r2; /* r6 */
959     }
960     else if (sc_r_power == 48.0)
961     {
962         rp               = r2*r2*r2; /* r6 */
963         rp               = rp*rp;    /* r12 */
964         rp               = rp*rp;    /* r24 */
965         rp               = rp*rp;    /* r48 */
966         rpm2             = rp/r2;    /* r46 */
967     }
968     else
969     {
970         rp             = pow(r2, 0.5*sc_r_power);  /* not currently supported as input, but can handle it */
971         rpm2           = rp/r2;
972     }
973
974     /* Loop over state A(0) and B(1) */
975     for (i = 0; i < 2; i++)
976     {
977         if ((c6[i] > 0) && (c12[i] > 0))
978         {
979             /* The c6 & c12 coefficients now contain the constants 6.0 and 12.0, respectively.
980              * Correct for this by multiplying with (1/12.0)/(1/6.0)=6.0/12.0=0.5.
981              */
982             sigma6[i]       = 0.5*c12[i]/c6[i];
983             sigma2[i]       = pow(0.5*c12[i]/c6[i], 1.0/3.0);
984             /* should be able to get rid of this ^^^ internal pow call eventually.  Will require agreement on
985                what data to store externally.  Can't be fixed without larger scale changes, so not 5.0 */
986             if (sigma6[i] < sigma6_min)   /* for disappearing coul and vdw with soft core at the same time */
987             {
988                 sigma6[i] = sigma6_min;
989                 sigma2[i] = sigma2_min;
990             }
991         }
992         else
993         {
994             sigma6[i]       = sigma6_def;
995             sigma2[i]       = sigma2_def;
996         }
997         if (sc_r_power == 6.0)
998         {
999             sigma_pow[i]    = sigma6[i];
1000             sigma_powm2[i]  = sigma6[i]/sigma2[i];
1001         }
1002         else if (sc_r_power == 48.0)
1003         {
1004             sigma_pow[i]    = sigma6[i]*sigma6[i];       /* sigma^12 */
1005             sigma_pow[i]    = sigma_pow[i]*sigma_pow[i]; /* sigma^24 */
1006             sigma_pow[i]    = sigma_pow[i]*sigma_pow[i]; /* sigma^48 */
1007             sigma_powm2[i]  = sigma_pow[i]/sigma2[i];
1008         }
1009         else
1010         {    /* not really supported as input, but in here for testing the general case*/
1011             sigma_pow[i]    = pow(sigma2[i], sc_r_power/2);
1012             sigma_powm2[i]  = sigma_pow[i]/(sigma2[i]);
1013         }
1014     }
1015
1016     /* only use softcore if one of the states has a zero endstate - softcore is for avoiding infinities!*/
1017     if ((c12[0] > 0) && (c12[1] > 0))
1018     {
1019         alpha_vdw_eff    = 0;
1020         alpha_coul_eff   = 0;
1021     }
1022     else
1023     {
1024         alpha_vdw_eff    = alpha_vdw;
1025         alpha_coul_eff   = alpha_coul;
1026     }
1027
1028     /* Loop over A and B states again */
1029     for (i = 0; i < 2; i++)
1030     {
1031         fscal_elec[i] = 0;
1032         fscal_vdw[i]  = 0;
1033         velec[i]      = 0;
1034         vvdw[i]       = 0;
1035
1036         /* Only spend time on A or B state if it is non-zero */
1037         if ( (qq[i] != 0) || (c6[i] != 0) || (c12[i] != 0) )
1038         {
1039             /* Coulomb */
1040             rpinv            = 1.0/(alpha_coul_eff*lfac_coul[i]*sigma_pow[i]+rp);
1041             r_coul           = pow(rpinv, -1.0/sc_r_power);
1042
1043             /* Electrostatics table lookup data */
1044             rtab             = r_coul*tabscale;
1045             ntab             = rtab;
1046             eps              = rtab-ntab;
1047             eps2             = eps*eps;
1048             ntab             = 12*ntab;
1049             /* Electrostatics */
1050             Y                = vftab[ntab];
1051             F                = vftab[ntab+1];
1052             Geps             = eps*vftab[ntab+2];
1053             Heps2            = eps2*vftab[ntab+3];
1054             Fp               = F+Geps+Heps2;
1055             VV               = Y+eps*Fp;
1056             FF               = Fp+Geps+2.0*Heps2;
1057             velec[i]         = qq[i]*VV;
1058             fscal_elec[i]    = -qq[i]*FF*r_coul*rpinv*tabscale;
1059
1060             /* Vdw */
1061             rpinv            = 1.0/(alpha_vdw_eff*lfac_vdw[i]*sigma_pow[i]+rp);
1062             r_vdw            = pow(rpinv, -1.0/sc_r_power);
1063             /* Vdw table lookup data */
1064             rtab             = r_vdw*tabscale;
1065             ntab             = rtab;
1066             eps              = rtab-ntab;
1067             eps2             = eps*eps;
1068             ntab             = 12*ntab;
1069             /* Dispersion */
1070             Y                = vftab[ntab+4];
1071             F                = vftab[ntab+5];
1072             Geps             = eps*vftab[ntab+6];
1073             Heps2            = eps2*vftab[ntab+7];
1074             Fp               = F+Geps+Heps2;
1075             VV               = Y+eps*Fp;
1076             FF               = Fp+Geps+2.0*Heps2;
1077             vvdw[i]          = c6[i]*VV;
1078             fscal_vdw[i]     = -c6[i]*FF;
1079
1080             /* Repulsion */
1081             Y                = vftab[ntab+8];
1082             F                = vftab[ntab+9];
1083             Geps             = eps*vftab[ntab+10];
1084             Heps2            = eps2*vftab[ntab+11];
1085             Fp               = F+Geps+Heps2;
1086             VV               = Y+eps*Fp;
1087             FF               = Fp+Geps+2.0*Heps2;
1088             vvdw[i]         += c12[i]*VV;
1089             fscal_vdw[i]    -= c12[i]*FF;
1090             fscal_vdw[i]    *= r_vdw*rpinv*tabscale;
1091         }
1092     }
1093     /* Now we have velec[i], vvdw[i], and fscal[i] for both states */
1094     /* Assemble A and B states */
1095     velecsum  = 0;
1096     vvdwsum   = 0;
1097     dvdl_coul = 0;
1098     dvdl_vdw  = 0;
1099     fscal     = 0;
1100     for (i = 0; i < 2; i++)
1101     {
1102         velecsum      += LFC[i]*velec[i];
1103         vvdwsum       += LFV[i]*vvdw[i];
1104
1105         fscal         += (LFC[i]*fscal_elec[i]+LFV[i]*fscal_vdw[i])*rpm2;
1106
1107         dvdl_coul     += velec[i]*DLF[i] + LFC[i]*alpha_coul_eff*dlfac_coul[i]*fscal_elec[i]*sigma_pow[i];
1108         dvdl_vdw      += vvdw[i]*DLF[i] + LFV[i]*alpha_vdw_eff*dlfac_vdw[i]*fscal_vdw[i]*sigma_pow[i];
1109     }
1110
1111     dvdl[efptCOUL]     += dvdl_coul;
1112     dvdl[efptVDW]      += dvdl_vdw;
1113
1114     *velectot           = velecsum;
1115     *vvdwtot            = vvdwsum;
1116
1117     return fscal;
1118 }