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