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