Use full path for legacyheaders
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_c / nb_kernel_ElecRFCut_VdwCSTab_GeomW3P1_c.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2013,2014, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35 /*
36  * Note: this file was generated by the GROMACS c kernel generator.
37  */
38 #include "config.h"
39
40 #include <math.h>
41
42 #include "../nb_kernel.h"
43 #include "gromacs/legacyheaders/types/simple.h"
44 #include "gromacs/math/vec.h"
45 #include "gromacs/legacyheaders/nrnb.h"
46
47 /*
48  * Gromacs nonbonded kernel:   nb_kernel_ElecRFCut_VdwCSTab_GeomW3P1_VF_c
49  * Electrostatics interaction: ReactionField
50  * VdW interaction:            CubicSplineTable
51  * Geometry:                   Water3-Particle
52  * Calculate force/pot:        PotentialAndForce
53  */
54 void
55 nb_kernel_ElecRFCut_VdwCSTab_GeomW3P1_VF_c
56                     (t_nblist                    * gmx_restrict       nlist,
57                      rvec                        * gmx_restrict          xx,
58                      rvec                        * gmx_restrict          ff,
59                      t_forcerec                  * gmx_restrict          fr,
60                      t_mdatoms                   * gmx_restrict     mdatoms,
61                      nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
62                      t_nrnb                      * gmx_restrict        nrnb)
63 {
64     int              i_shift_offset,i_coord_offset,j_coord_offset;
65     int              j_index_start,j_index_end;
66     int              nri,inr,ggid,iidx,jidx,jnr,outeriter,inneriter;
67     real             shX,shY,shZ,tx,ty,tz,fscal,rcutoff,rcutoff2;
68     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
69     real             *shiftvec,*fshift,*x,*f;
70     int              vdwioffset0;
71     real             ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
72     int              vdwioffset1;
73     real             ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
74     int              vdwioffset2;
75     real             ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
76     int              vdwjidx0;
77     real             jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
78     real             dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00,cexp1_00,cexp2_00;
79     real             dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10,cexp1_10,cexp2_10;
80     real             dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20,cexp1_20,cexp2_20;
81     real             velec,felec,velecsum,facel,crf,krf,krf2;
82     real             *charge;
83     int              nvdwtype;
84     real             rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,br,vvdwexp,sh_vdw_invrcut6;
85     int              *vdwtype;
86     real             *vdwparam;
87     int              vfitab;
88     real             rt,vfeps,vftabscale,Y,F,Geps,Heps2,Fp,VV,FF;
89     real             *vftab;
90
91     x                = xx[0];
92     f                = ff[0];
93
94     nri              = nlist->nri;
95     iinr             = nlist->iinr;
96     jindex           = nlist->jindex;
97     jjnr             = nlist->jjnr;
98     shiftidx         = nlist->shift;
99     gid              = nlist->gid;
100     shiftvec         = fr->shift_vec[0];
101     fshift           = fr->fshift[0];
102     facel            = fr->epsfac;
103     charge           = mdatoms->chargeA;
104     krf              = fr->ic->k_rf;
105     krf2             = krf*2.0;
106     crf              = fr->ic->c_rf;
107     nvdwtype         = fr->ntype;
108     vdwparam         = fr->nbfp;
109     vdwtype          = mdatoms->typeA;
110
111     vftab            = kernel_data->table_vdw->data;
112     vftabscale       = kernel_data->table_vdw->scale;
113
114     /* Setup water-specific parameters */
115     inr              = nlist->iinr[0];
116     iq0              = facel*charge[inr+0];
117     iq1              = facel*charge[inr+1];
118     iq2              = facel*charge[inr+2];
119     vdwioffset0      = 2*nvdwtype*vdwtype[inr+0];
120
121     /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
122     rcutoff          = fr->rcoulomb;
123     rcutoff2         = rcutoff*rcutoff;
124
125     outeriter        = 0;
126     inneriter        = 0;
127
128     /* Start outer loop over neighborlists */
129     for(iidx=0; iidx<nri; iidx++)
130     {
131         /* Load shift vector for this list */
132         i_shift_offset   = DIM*shiftidx[iidx];
133         shX              = shiftvec[i_shift_offset+XX];
134         shY              = shiftvec[i_shift_offset+YY];
135         shZ              = shiftvec[i_shift_offset+ZZ];
136
137         /* Load limits for loop over neighbors */
138         j_index_start    = jindex[iidx];
139         j_index_end      = jindex[iidx+1];
140
141         /* Get outer coordinate index */
142         inr              = iinr[iidx];
143         i_coord_offset   = DIM*inr;
144
145         /* Load i particle coords and add shift vector */
146         ix0              = shX + x[i_coord_offset+DIM*0+XX];
147         iy0              = shY + x[i_coord_offset+DIM*0+YY];
148         iz0              = shZ + x[i_coord_offset+DIM*0+ZZ];
149         ix1              = shX + x[i_coord_offset+DIM*1+XX];
150         iy1              = shY + x[i_coord_offset+DIM*1+YY];
151         iz1              = shZ + x[i_coord_offset+DIM*1+ZZ];
152         ix2              = shX + x[i_coord_offset+DIM*2+XX];
153         iy2              = shY + x[i_coord_offset+DIM*2+YY];
154         iz2              = shZ + x[i_coord_offset+DIM*2+ZZ];
155
156         fix0             = 0.0;
157         fiy0             = 0.0;
158         fiz0             = 0.0;
159         fix1             = 0.0;
160         fiy1             = 0.0;
161         fiz1             = 0.0;
162         fix2             = 0.0;
163         fiy2             = 0.0;
164         fiz2             = 0.0;
165
166         /* Reset potential sums */
167         velecsum         = 0.0;
168         vvdwsum          = 0.0;
169
170         /* Start inner kernel loop */
171         for(jidx=j_index_start; jidx<j_index_end; jidx++)
172         {
173             /* Get j neighbor index, and coordinate index */
174             jnr              = jjnr[jidx];
175             j_coord_offset   = DIM*jnr;
176
177             /* load j atom coordinates */
178             jx0              = x[j_coord_offset+DIM*0+XX];
179             jy0              = x[j_coord_offset+DIM*0+YY];
180             jz0              = x[j_coord_offset+DIM*0+ZZ];
181
182             /* Calculate displacement vector */
183             dx00             = ix0 - jx0;
184             dy00             = iy0 - jy0;
185             dz00             = iz0 - jz0;
186             dx10             = ix1 - jx0;
187             dy10             = iy1 - jy0;
188             dz10             = iz1 - jz0;
189             dx20             = ix2 - jx0;
190             dy20             = iy2 - jy0;
191             dz20             = iz2 - jz0;
192
193             /* Calculate squared distance and things based on it */
194             rsq00            = dx00*dx00+dy00*dy00+dz00*dz00;
195             rsq10            = dx10*dx10+dy10*dy10+dz10*dz10;
196             rsq20            = dx20*dx20+dy20*dy20+dz20*dz20;
197
198             rinv00           = gmx_invsqrt(rsq00);
199             rinv10           = gmx_invsqrt(rsq10);
200             rinv20           = gmx_invsqrt(rsq20);
201
202             rinvsq00         = rinv00*rinv00;
203             rinvsq10         = rinv10*rinv10;
204             rinvsq20         = rinv20*rinv20;
205
206             /* Load parameters for j particles */
207             jq0              = charge[jnr+0];
208             vdwjidx0         = 2*vdwtype[jnr+0];
209
210             /**************************
211              * CALCULATE INTERACTIONS *
212              **************************/
213
214             if (rsq00<rcutoff2)
215             {
216
217             r00              = rsq00*rinv00;
218
219             qq00             = iq0*jq0;
220             c6_00            = vdwparam[vdwioffset0+vdwjidx0];
221             c12_00           = vdwparam[vdwioffset0+vdwjidx0+1];
222
223             /* Calculate table index by multiplying r with table scale and truncate to integer */
224             rt               = r00*vftabscale;
225             vfitab           = rt;
226             vfeps            = rt-vfitab;
227             vfitab           = 2*4*vfitab;
228
229             /* REACTION-FIELD ELECTROSTATICS */
230             velec            = qq00*(rinv00+krf*rsq00-crf);
231             felec            = qq00*(rinv00*rinvsq00-krf2);
232
233             /* CUBIC SPLINE TABLE DISPERSION */
234             vfitab          += 0;
235             Y                = vftab[vfitab];
236             F                = vftab[vfitab+1];
237             Geps             = vfeps*vftab[vfitab+2];
238             Heps2            = vfeps*vfeps*vftab[vfitab+3];
239             Fp               = F+Geps+Heps2;
240             VV               = Y+vfeps*Fp;
241             vvdw6            = c6_00*VV;
242             FF               = Fp+Geps+2.0*Heps2;
243             fvdw6            = c6_00*FF;
244
245             /* CUBIC SPLINE TABLE REPULSION */
246             Y                = vftab[vfitab+4];
247             F                = vftab[vfitab+5];
248             Geps             = vfeps*vftab[vfitab+6];
249             Heps2            = vfeps*vfeps*vftab[vfitab+7];
250             Fp               = F+Geps+Heps2;
251             VV               = Y+vfeps*Fp;
252             vvdw12           = c12_00*VV;
253             FF               = Fp+Geps+2.0*Heps2;
254             fvdw12           = c12_00*FF;
255             vvdw             = vvdw12+vvdw6;
256             fvdw             = -(fvdw6+fvdw12)*vftabscale*rinv00;
257
258             /* Update potential sums from outer loop */
259             velecsum        += velec;
260             vvdwsum         += vvdw;
261
262             fscal            = felec+fvdw;
263
264             /* Calculate temporary vectorial force */
265             tx               = fscal*dx00;
266             ty               = fscal*dy00;
267             tz               = fscal*dz00;
268
269             /* Update vectorial force */
270             fix0            += tx;
271             fiy0            += ty;
272             fiz0            += tz;
273             f[j_coord_offset+DIM*0+XX] -= tx;
274             f[j_coord_offset+DIM*0+YY] -= ty;
275             f[j_coord_offset+DIM*0+ZZ] -= tz;
276
277             }
278
279             /**************************
280              * CALCULATE INTERACTIONS *
281              **************************/
282
283             if (rsq10<rcutoff2)
284             {
285
286             qq10             = iq1*jq0;
287
288             /* REACTION-FIELD ELECTROSTATICS */
289             velec            = qq10*(rinv10+krf*rsq10-crf);
290             felec            = qq10*(rinv10*rinvsq10-krf2);
291
292             /* Update potential sums from outer loop */
293             velecsum        += velec;
294
295             fscal            = felec;
296
297             /* Calculate temporary vectorial force */
298             tx               = fscal*dx10;
299             ty               = fscal*dy10;
300             tz               = fscal*dz10;
301
302             /* Update vectorial force */
303             fix1            += tx;
304             fiy1            += ty;
305             fiz1            += tz;
306             f[j_coord_offset+DIM*0+XX] -= tx;
307             f[j_coord_offset+DIM*0+YY] -= ty;
308             f[j_coord_offset+DIM*0+ZZ] -= tz;
309
310             }
311
312             /**************************
313              * CALCULATE INTERACTIONS *
314              **************************/
315
316             if (rsq20<rcutoff2)
317             {
318
319             qq20             = iq2*jq0;
320
321             /* REACTION-FIELD ELECTROSTATICS */
322             velec            = qq20*(rinv20+krf*rsq20-crf);
323             felec            = qq20*(rinv20*rinvsq20-krf2);
324
325             /* Update potential sums from outer loop */
326             velecsum        += velec;
327
328             fscal            = felec;
329
330             /* Calculate temporary vectorial force */
331             tx               = fscal*dx20;
332             ty               = fscal*dy20;
333             tz               = fscal*dz20;
334
335             /* Update vectorial force */
336             fix2            += tx;
337             fiy2            += ty;
338             fiz2            += tz;
339             f[j_coord_offset+DIM*0+XX] -= tx;
340             f[j_coord_offset+DIM*0+YY] -= ty;
341             f[j_coord_offset+DIM*0+ZZ] -= tz;
342
343             }
344
345             /* Inner loop uses 130 flops */
346         }
347         /* End of innermost loop */
348
349         tx = ty = tz = 0;
350         f[i_coord_offset+DIM*0+XX] += fix0;
351         f[i_coord_offset+DIM*0+YY] += fiy0;
352         f[i_coord_offset+DIM*0+ZZ] += fiz0;
353         tx                         += fix0;
354         ty                         += fiy0;
355         tz                         += fiz0;
356         f[i_coord_offset+DIM*1+XX] += fix1;
357         f[i_coord_offset+DIM*1+YY] += fiy1;
358         f[i_coord_offset+DIM*1+ZZ] += fiz1;
359         tx                         += fix1;
360         ty                         += fiy1;
361         tz                         += fiz1;
362         f[i_coord_offset+DIM*2+XX] += fix2;
363         f[i_coord_offset+DIM*2+YY] += fiy2;
364         f[i_coord_offset+DIM*2+ZZ] += fiz2;
365         tx                         += fix2;
366         ty                         += fiy2;
367         tz                         += fiz2;
368         fshift[i_shift_offset+XX]  += tx;
369         fshift[i_shift_offset+YY]  += ty;
370         fshift[i_shift_offset+ZZ]  += tz;
371
372         ggid                        = gid[iidx];
373         /* Update potential energies */
374         kernel_data->energygrp_elec[ggid] += velecsum;
375         kernel_data->energygrp_vdw[ggid] += vvdwsum;
376
377         /* Increment number of inner iterations */
378         inneriter                  += j_index_end - j_index_start;
379
380         /* Outer loop uses 32 flops */
381     }
382
383     /* Increment number of outer iterations */
384     outeriter        += nri;
385
386     /* Update outer/inner flops */
387
388     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W3_VF,outeriter*32 + inneriter*130);
389 }
390 /*
391  * Gromacs nonbonded kernel:   nb_kernel_ElecRFCut_VdwCSTab_GeomW3P1_F_c
392  * Electrostatics interaction: ReactionField
393  * VdW interaction:            CubicSplineTable
394  * Geometry:                   Water3-Particle
395  * Calculate force/pot:        Force
396  */
397 void
398 nb_kernel_ElecRFCut_VdwCSTab_GeomW3P1_F_c
399                     (t_nblist                    * gmx_restrict       nlist,
400                      rvec                        * gmx_restrict          xx,
401                      rvec                        * gmx_restrict          ff,
402                      t_forcerec                  * gmx_restrict          fr,
403                      t_mdatoms                   * gmx_restrict     mdatoms,
404                      nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
405                      t_nrnb                      * gmx_restrict        nrnb)
406 {
407     int              i_shift_offset,i_coord_offset,j_coord_offset;
408     int              j_index_start,j_index_end;
409     int              nri,inr,ggid,iidx,jidx,jnr,outeriter,inneriter;
410     real             shX,shY,shZ,tx,ty,tz,fscal,rcutoff,rcutoff2;
411     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
412     real             *shiftvec,*fshift,*x,*f;
413     int              vdwioffset0;
414     real             ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
415     int              vdwioffset1;
416     real             ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
417     int              vdwioffset2;
418     real             ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
419     int              vdwjidx0;
420     real             jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
421     real             dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00,cexp1_00,cexp2_00;
422     real             dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10,cexp1_10,cexp2_10;
423     real             dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20,cexp1_20,cexp2_20;
424     real             velec,felec,velecsum,facel,crf,krf,krf2;
425     real             *charge;
426     int              nvdwtype;
427     real             rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,br,vvdwexp,sh_vdw_invrcut6;
428     int              *vdwtype;
429     real             *vdwparam;
430     int              vfitab;
431     real             rt,vfeps,vftabscale,Y,F,Geps,Heps2,Fp,VV,FF;
432     real             *vftab;
433
434     x                = xx[0];
435     f                = ff[0];
436
437     nri              = nlist->nri;
438     iinr             = nlist->iinr;
439     jindex           = nlist->jindex;
440     jjnr             = nlist->jjnr;
441     shiftidx         = nlist->shift;
442     gid              = nlist->gid;
443     shiftvec         = fr->shift_vec[0];
444     fshift           = fr->fshift[0];
445     facel            = fr->epsfac;
446     charge           = mdatoms->chargeA;
447     krf              = fr->ic->k_rf;
448     krf2             = krf*2.0;
449     crf              = fr->ic->c_rf;
450     nvdwtype         = fr->ntype;
451     vdwparam         = fr->nbfp;
452     vdwtype          = mdatoms->typeA;
453
454     vftab            = kernel_data->table_vdw->data;
455     vftabscale       = kernel_data->table_vdw->scale;
456
457     /* Setup water-specific parameters */
458     inr              = nlist->iinr[0];
459     iq0              = facel*charge[inr+0];
460     iq1              = facel*charge[inr+1];
461     iq2              = facel*charge[inr+2];
462     vdwioffset0      = 2*nvdwtype*vdwtype[inr+0];
463
464     /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
465     rcutoff          = fr->rcoulomb;
466     rcutoff2         = rcutoff*rcutoff;
467
468     outeriter        = 0;
469     inneriter        = 0;
470
471     /* Start outer loop over neighborlists */
472     for(iidx=0; iidx<nri; iidx++)
473     {
474         /* Load shift vector for this list */
475         i_shift_offset   = DIM*shiftidx[iidx];
476         shX              = shiftvec[i_shift_offset+XX];
477         shY              = shiftvec[i_shift_offset+YY];
478         shZ              = shiftvec[i_shift_offset+ZZ];
479
480         /* Load limits for loop over neighbors */
481         j_index_start    = jindex[iidx];
482         j_index_end      = jindex[iidx+1];
483
484         /* Get outer coordinate index */
485         inr              = iinr[iidx];
486         i_coord_offset   = DIM*inr;
487
488         /* Load i particle coords and add shift vector */
489         ix0              = shX + x[i_coord_offset+DIM*0+XX];
490         iy0              = shY + x[i_coord_offset+DIM*0+YY];
491         iz0              = shZ + x[i_coord_offset+DIM*0+ZZ];
492         ix1              = shX + x[i_coord_offset+DIM*1+XX];
493         iy1              = shY + x[i_coord_offset+DIM*1+YY];
494         iz1              = shZ + x[i_coord_offset+DIM*1+ZZ];
495         ix2              = shX + x[i_coord_offset+DIM*2+XX];
496         iy2              = shY + x[i_coord_offset+DIM*2+YY];
497         iz2              = shZ + x[i_coord_offset+DIM*2+ZZ];
498
499         fix0             = 0.0;
500         fiy0             = 0.0;
501         fiz0             = 0.0;
502         fix1             = 0.0;
503         fiy1             = 0.0;
504         fiz1             = 0.0;
505         fix2             = 0.0;
506         fiy2             = 0.0;
507         fiz2             = 0.0;
508
509         /* Start inner kernel loop */
510         for(jidx=j_index_start; jidx<j_index_end; jidx++)
511         {
512             /* Get j neighbor index, and coordinate index */
513             jnr              = jjnr[jidx];
514             j_coord_offset   = DIM*jnr;
515
516             /* load j atom coordinates */
517             jx0              = x[j_coord_offset+DIM*0+XX];
518             jy0              = x[j_coord_offset+DIM*0+YY];
519             jz0              = x[j_coord_offset+DIM*0+ZZ];
520
521             /* Calculate displacement vector */
522             dx00             = ix0 - jx0;
523             dy00             = iy0 - jy0;
524             dz00             = iz0 - jz0;
525             dx10             = ix1 - jx0;
526             dy10             = iy1 - jy0;
527             dz10             = iz1 - jz0;
528             dx20             = ix2 - jx0;
529             dy20             = iy2 - jy0;
530             dz20             = iz2 - jz0;
531
532             /* Calculate squared distance and things based on it */
533             rsq00            = dx00*dx00+dy00*dy00+dz00*dz00;
534             rsq10            = dx10*dx10+dy10*dy10+dz10*dz10;
535             rsq20            = dx20*dx20+dy20*dy20+dz20*dz20;
536
537             rinv00           = gmx_invsqrt(rsq00);
538             rinv10           = gmx_invsqrt(rsq10);
539             rinv20           = gmx_invsqrt(rsq20);
540
541             rinvsq00         = rinv00*rinv00;
542             rinvsq10         = rinv10*rinv10;
543             rinvsq20         = rinv20*rinv20;
544
545             /* Load parameters for j particles */
546             jq0              = charge[jnr+0];
547             vdwjidx0         = 2*vdwtype[jnr+0];
548
549             /**************************
550              * CALCULATE INTERACTIONS *
551              **************************/
552
553             if (rsq00<rcutoff2)
554             {
555
556             r00              = rsq00*rinv00;
557
558             qq00             = iq0*jq0;
559             c6_00            = vdwparam[vdwioffset0+vdwjidx0];
560             c12_00           = vdwparam[vdwioffset0+vdwjidx0+1];
561
562             /* Calculate table index by multiplying r with table scale and truncate to integer */
563             rt               = r00*vftabscale;
564             vfitab           = rt;
565             vfeps            = rt-vfitab;
566             vfitab           = 2*4*vfitab;
567
568             /* REACTION-FIELD ELECTROSTATICS */
569             felec            = qq00*(rinv00*rinvsq00-krf2);
570
571             /* CUBIC SPLINE TABLE DISPERSION */
572             vfitab          += 0;
573             F                = vftab[vfitab+1];
574             Geps             = vfeps*vftab[vfitab+2];
575             Heps2            = vfeps*vfeps*vftab[vfitab+3];
576             Fp               = F+Geps+Heps2;
577             FF               = Fp+Geps+2.0*Heps2;
578             fvdw6            = c6_00*FF;
579
580             /* CUBIC SPLINE TABLE REPULSION */
581             F                = vftab[vfitab+5];
582             Geps             = vfeps*vftab[vfitab+6];
583             Heps2            = vfeps*vfeps*vftab[vfitab+7];
584             Fp               = F+Geps+Heps2;
585             FF               = Fp+Geps+2.0*Heps2;
586             fvdw12           = c12_00*FF;
587             fvdw             = -(fvdw6+fvdw12)*vftabscale*rinv00;
588
589             fscal            = felec+fvdw;
590
591             /* Calculate temporary vectorial force */
592             tx               = fscal*dx00;
593             ty               = fscal*dy00;
594             tz               = fscal*dz00;
595
596             /* Update vectorial force */
597             fix0            += tx;
598             fiy0            += ty;
599             fiz0            += tz;
600             f[j_coord_offset+DIM*0+XX] -= tx;
601             f[j_coord_offset+DIM*0+YY] -= ty;
602             f[j_coord_offset+DIM*0+ZZ] -= tz;
603
604             }
605
606             /**************************
607              * CALCULATE INTERACTIONS *
608              **************************/
609
610             if (rsq10<rcutoff2)
611             {
612
613             qq10             = iq1*jq0;
614
615             /* REACTION-FIELD ELECTROSTATICS */
616             felec            = qq10*(rinv10*rinvsq10-krf2);
617
618             fscal            = felec;
619
620             /* Calculate temporary vectorial force */
621             tx               = fscal*dx10;
622             ty               = fscal*dy10;
623             tz               = fscal*dz10;
624
625             /* Update vectorial force */
626             fix1            += tx;
627             fiy1            += ty;
628             fiz1            += tz;
629             f[j_coord_offset+DIM*0+XX] -= tx;
630             f[j_coord_offset+DIM*0+YY] -= ty;
631             f[j_coord_offset+DIM*0+ZZ] -= tz;
632
633             }
634
635             /**************************
636              * CALCULATE INTERACTIONS *
637              **************************/
638
639             if (rsq20<rcutoff2)
640             {
641
642             qq20             = iq2*jq0;
643
644             /* REACTION-FIELD ELECTROSTATICS */
645             felec            = qq20*(rinv20*rinvsq20-krf2);
646
647             fscal            = felec;
648
649             /* Calculate temporary vectorial force */
650             tx               = fscal*dx20;
651             ty               = fscal*dy20;
652             tz               = fscal*dz20;
653
654             /* Update vectorial force */
655             fix2            += tx;
656             fiy2            += ty;
657             fiz2            += tz;
658             f[j_coord_offset+DIM*0+XX] -= tx;
659             f[j_coord_offset+DIM*0+YY] -= ty;
660             f[j_coord_offset+DIM*0+ZZ] -= tz;
661
662             }
663
664             /* Inner loop uses 107 flops */
665         }
666         /* End of innermost loop */
667
668         tx = ty = tz = 0;
669         f[i_coord_offset+DIM*0+XX] += fix0;
670         f[i_coord_offset+DIM*0+YY] += fiy0;
671         f[i_coord_offset+DIM*0+ZZ] += fiz0;
672         tx                         += fix0;
673         ty                         += fiy0;
674         tz                         += fiz0;
675         f[i_coord_offset+DIM*1+XX] += fix1;
676         f[i_coord_offset+DIM*1+YY] += fiy1;
677         f[i_coord_offset+DIM*1+ZZ] += fiz1;
678         tx                         += fix1;
679         ty                         += fiy1;
680         tz                         += fiz1;
681         f[i_coord_offset+DIM*2+XX] += fix2;
682         f[i_coord_offset+DIM*2+YY] += fiy2;
683         f[i_coord_offset+DIM*2+ZZ] += fiz2;
684         tx                         += fix2;
685         ty                         += fiy2;
686         tz                         += fiz2;
687         fshift[i_shift_offset+XX]  += tx;
688         fshift[i_shift_offset+YY]  += ty;
689         fshift[i_shift_offset+ZZ]  += tz;
690
691         /* Increment number of inner iterations */
692         inneriter                  += j_index_end - j_index_start;
693
694         /* Outer loop uses 30 flops */
695     }
696
697     /* Increment number of outer iterations */
698     outeriter        += nri;
699
700     /* Update outer/inner flops */
701
702     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W3_F,outeriter*30 + inneriter*107);
703 }