95f47d4824a3265aa1dd4d7cfe63a770a1adc143
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_c / nb_kernel_ElecEw_VdwBham_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_ElecEw_VdwBham_GeomW3P1_VF_c
49  * Electrostatics interaction: Ewald
50  * VdW interaction:            Buckingham
51  * Geometry:                   Water3-Particle
52  * Calculate force/pot:        PotentialAndForce
53  */
54 void
55 nb_kernel_ElecEw_VdwBham_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              ewitab;
88     real             ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace;
89     real             *ewtab;
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     nvdwtype         = fr->ntype;
105     vdwparam         = fr->nbfp;
106     vdwtype          = mdatoms->typeA;
107
108     sh_ewald         = fr->ic->sh_ewald;
109     ewtab            = fr->ic->tabq_coul_FDV0;
110     ewtabscale       = fr->ic->tabq_scale;
111     ewtabhalfspace   = 0.5/ewtabscale;
112
113     /* Setup water-specific parameters */
114     inr              = nlist->iinr[0];
115     iq0              = facel*charge[inr+0];
116     iq1              = facel*charge[inr+1];
117     iq2              = facel*charge[inr+2];
118     vdwioffset0      = 3*nvdwtype*vdwtype[inr+0];
119
120     outeriter        = 0;
121     inneriter        = 0;
122
123     /* Start outer loop over neighborlists */
124     for(iidx=0; iidx<nri; iidx++)
125     {
126         /* Load shift vector for this list */
127         i_shift_offset   = DIM*shiftidx[iidx];
128         shX              = shiftvec[i_shift_offset+XX];
129         shY              = shiftvec[i_shift_offset+YY];
130         shZ              = shiftvec[i_shift_offset+ZZ];
131
132         /* Load limits for loop over neighbors */
133         j_index_start    = jindex[iidx];
134         j_index_end      = jindex[iidx+1];
135
136         /* Get outer coordinate index */
137         inr              = iinr[iidx];
138         i_coord_offset   = DIM*inr;
139
140         /* Load i particle coords and add shift vector */
141         ix0              = shX + x[i_coord_offset+DIM*0+XX];
142         iy0              = shY + x[i_coord_offset+DIM*0+YY];
143         iz0              = shZ + x[i_coord_offset+DIM*0+ZZ];
144         ix1              = shX + x[i_coord_offset+DIM*1+XX];
145         iy1              = shY + x[i_coord_offset+DIM*1+YY];
146         iz1              = shZ + x[i_coord_offset+DIM*1+ZZ];
147         ix2              = shX + x[i_coord_offset+DIM*2+XX];
148         iy2              = shY + x[i_coord_offset+DIM*2+YY];
149         iz2              = shZ + x[i_coord_offset+DIM*2+ZZ];
150
151         fix0             = 0.0;
152         fiy0             = 0.0;
153         fiz0             = 0.0;
154         fix1             = 0.0;
155         fiy1             = 0.0;
156         fiz1             = 0.0;
157         fix2             = 0.0;
158         fiy2             = 0.0;
159         fiz2             = 0.0;
160
161         /* Reset potential sums */
162         velecsum         = 0.0;
163         vvdwsum          = 0.0;
164
165         /* Start inner kernel loop */
166         for(jidx=j_index_start; jidx<j_index_end; jidx++)
167         {
168             /* Get j neighbor index, and coordinate index */
169             jnr              = jjnr[jidx];
170             j_coord_offset   = DIM*jnr;
171
172             /* load j atom coordinates */
173             jx0              = x[j_coord_offset+DIM*0+XX];
174             jy0              = x[j_coord_offset+DIM*0+YY];
175             jz0              = x[j_coord_offset+DIM*0+ZZ];
176
177             /* Calculate displacement vector */
178             dx00             = ix0 - jx0;
179             dy00             = iy0 - jy0;
180             dz00             = iz0 - jz0;
181             dx10             = ix1 - jx0;
182             dy10             = iy1 - jy0;
183             dz10             = iz1 - jz0;
184             dx20             = ix2 - jx0;
185             dy20             = iy2 - jy0;
186             dz20             = iz2 - jz0;
187
188             /* Calculate squared distance and things based on it */
189             rsq00            = dx00*dx00+dy00*dy00+dz00*dz00;
190             rsq10            = dx10*dx10+dy10*dy10+dz10*dz10;
191             rsq20            = dx20*dx20+dy20*dy20+dz20*dz20;
192
193             rinv00           = gmx_invsqrt(rsq00);
194             rinv10           = gmx_invsqrt(rsq10);
195             rinv20           = gmx_invsqrt(rsq20);
196
197             rinvsq00         = rinv00*rinv00;
198             rinvsq10         = rinv10*rinv10;
199             rinvsq20         = rinv20*rinv20;
200
201             /* Load parameters for j particles */
202             jq0              = charge[jnr+0];
203             vdwjidx0         = 3*vdwtype[jnr+0];
204
205             /**************************
206              * CALCULATE INTERACTIONS *
207              **************************/
208
209             r00              = rsq00*rinv00;
210
211             qq00             = iq0*jq0;
212             c6_00            = vdwparam[vdwioffset0+vdwjidx0];
213             cexp1_00         = vdwparam[vdwioffset0+vdwjidx0+1];
214             cexp2_00         = vdwparam[vdwioffset0+vdwjidx0+2];
215
216             /* EWALD ELECTROSTATICS */
217
218             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
219             ewrt             = r00*ewtabscale;
220             ewitab           = ewrt;
221             eweps            = ewrt-ewitab;
222             ewitab           = 4*ewitab;
223             felec            = ewtab[ewitab]+eweps*ewtab[ewitab+1];
224             velec            = qq00*(rinv00-(ewtab[ewitab+2]-ewtabhalfspace*eweps*(ewtab[ewitab]+felec)));
225             felec            = qq00*rinv00*(rinvsq00-felec);
226
227             /* BUCKINGHAM DISPERSION/REPULSION */
228             rinvsix          = rinvsq00*rinvsq00*rinvsq00;
229             vvdw6            = c6_00*rinvsix;
230             br               = cexp2_00*r00;
231             vvdwexp          = cexp1_00*exp(-br);
232             vvdw             = vvdwexp - vvdw6*(1.0/6.0);
233             fvdw             = (br*vvdwexp-vvdw6)*rinvsq00;
234
235             /* Update potential sums from outer loop */
236             velecsum        += velec;
237             vvdwsum         += vvdw;
238
239             fscal            = felec+fvdw;
240
241             /* Calculate temporary vectorial force */
242             tx               = fscal*dx00;
243             ty               = fscal*dy00;
244             tz               = fscal*dz00;
245
246             /* Update vectorial force */
247             fix0            += tx;
248             fiy0            += ty;
249             fiz0            += tz;
250             f[j_coord_offset+DIM*0+XX] -= tx;
251             f[j_coord_offset+DIM*0+YY] -= ty;
252             f[j_coord_offset+DIM*0+ZZ] -= tz;
253
254             /**************************
255              * CALCULATE INTERACTIONS *
256              **************************/
257
258             r10              = rsq10*rinv10;
259
260             qq10             = iq1*jq0;
261
262             /* EWALD ELECTROSTATICS */
263
264             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
265             ewrt             = r10*ewtabscale;
266             ewitab           = ewrt;
267             eweps            = ewrt-ewitab;
268             ewitab           = 4*ewitab;
269             felec            = ewtab[ewitab]+eweps*ewtab[ewitab+1];
270             velec            = qq10*(rinv10-(ewtab[ewitab+2]-ewtabhalfspace*eweps*(ewtab[ewitab]+felec)));
271             felec            = qq10*rinv10*(rinvsq10-felec);
272
273             /* Update potential sums from outer loop */
274             velecsum        += velec;
275
276             fscal            = felec;
277
278             /* Calculate temporary vectorial force */
279             tx               = fscal*dx10;
280             ty               = fscal*dy10;
281             tz               = fscal*dz10;
282
283             /* Update vectorial force */
284             fix1            += tx;
285             fiy1            += ty;
286             fiz1            += tz;
287             f[j_coord_offset+DIM*0+XX] -= tx;
288             f[j_coord_offset+DIM*0+YY] -= ty;
289             f[j_coord_offset+DIM*0+ZZ] -= tz;
290
291             /**************************
292              * CALCULATE INTERACTIONS *
293              **************************/
294
295             r20              = rsq20*rinv20;
296
297             qq20             = iq2*jq0;
298
299             /* EWALD ELECTROSTATICS */
300
301             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
302             ewrt             = r20*ewtabscale;
303             ewitab           = ewrt;
304             eweps            = ewrt-ewitab;
305             ewitab           = 4*ewitab;
306             felec            = ewtab[ewitab]+eweps*ewtab[ewitab+1];
307             velec            = qq20*(rinv20-(ewtab[ewitab+2]-ewtabhalfspace*eweps*(ewtab[ewitab]+felec)));
308             felec            = qq20*rinv20*(rinvsq20-felec);
309
310             /* Update potential sums from outer loop */
311             velecsum        += velec;
312
313             fscal            = felec;
314
315             /* Calculate temporary vectorial force */
316             tx               = fscal*dx20;
317             ty               = fscal*dy20;
318             tz               = fscal*dz20;
319
320             /* Update vectorial force */
321             fix2            += tx;
322             fiy2            += ty;
323             fiz2            += tz;
324             f[j_coord_offset+DIM*0+XX] -= tx;
325             f[j_coord_offset+DIM*0+YY] -= ty;
326             f[j_coord_offset+DIM*0+ZZ] -= tz;
327
328             /* Inner loop uses 161 flops */
329         }
330         /* End of innermost loop */
331
332         tx = ty = tz = 0;
333         f[i_coord_offset+DIM*0+XX] += fix0;
334         f[i_coord_offset+DIM*0+YY] += fiy0;
335         f[i_coord_offset+DIM*0+ZZ] += fiz0;
336         tx                         += fix0;
337         ty                         += fiy0;
338         tz                         += fiz0;
339         f[i_coord_offset+DIM*1+XX] += fix1;
340         f[i_coord_offset+DIM*1+YY] += fiy1;
341         f[i_coord_offset+DIM*1+ZZ] += fiz1;
342         tx                         += fix1;
343         ty                         += fiy1;
344         tz                         += fiz1;
345         f[i_coord_offset+DIM*2+XX] += fix2;
346         f[i_coord_offset+DIM*2+YY] += fiy2;
347         f[i_coord_offset+DIM*2+ZZ] += fiz2;
348         tx                         += fix2;
349         ty                         += fiy2;
350         tz                         += fiz2;
351         fshift[i_shift_offset+XX]  += tx;
352         fshift[i_shift_offset+YY]  += ty;
353         fshift[i_shift_offset+ZZ]  += tz;
354
355         ggid                        = gid[iidx];
356         /* Update potential energies */
357         kernel_data->energygrp_elec[ggid] += velecsum;
358         kernel_data->energygrp_vdw[ggid] += vvdwsum;
359
360         /* Increment number of inner iterations */
361         inneriter                  += j_index_end - j_index_start;
362
363         /* Outer loop uses 32 flops */
364     }
365
366     /* Increment number of outer iterations */
367     outeriter        += nri;
368
369     /* Update outer/inner flops */
370
371     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W3_VF,outeriter*32 + inneriter*161);
372 }
373 /*
374  * Gromacs nonbonded kernel:   nb_kernel_ElecEw_VdwBham_GeomW3P1_F_c
375  * Electrostatics interaction: Ewald
376  * VdW interaction:            Buckingham
377  * Geometry:                   Water3-Particle
378  * Calculate force/pot:        Force
379  */
380 void
381 nb_kernel_ElecEw_VdwBham_GeomW3P1_F_c
382                     (t_nblist                    * gmx_restrict       nlist,
383                      rvec                        * gmx_restrict          xx,
384                      rvec                        * gmx_restrict          ff,
385                      t_forcerec                  * gmx_restrict          fr,
386                      t_mdatoms                   * gmx_restrict     mdatoms,
387                      nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
388                      t_nrnb                      * gmx_restrict        nrnb)
389 {
390     int              i_shift_offset,i_coord_offset,j_coord_offset;
391     int              j_index_start,j_index_end;
392     int              nri,inr,ggid,iidx,jidx,jnr,outeriter,inneriter;
393     real             shX,shY,shZ,tx,ty,tz,fscal,rcutoff,rcutoff2;
394     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
395     real             *shiftvec,*fshift,*x,*f;
396     int              vdwioffset0;
397     real             ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
398     int              vdwioffset1;
399     real             ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
400     int              vdwioffset2;
401     real             ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
402     int              vdwjidx0;
403     real             jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
404     real             dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00,cexp1_00,cexp2_00;
405     real             dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10,cexp1_10,cexp2_10;
406     real             dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20,cexp1_20,cexp2_20;
407     real             velec,felec,velecsum,facel,crf,krf,krf2;
408     real             *charge;
409     int              nvdwtype;
410     real             rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,br,vvdwexp,sh_vdw_invrcut6;
411     int              *vdwtype;
412     real             *vdwparam;
413     int              ewitab;
414     real             ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace;
415     real             *ewtab;
416
417     x                = xx[0];
418     f                = ff[0];
419
420     nri              = nlist->nri;
421     iinr             = nlist->iinr;
422     jindex           = nlist->jindex;
423     jjnr             = nlist->jjnr;
424     shiftidx         = nlist->shift;
425     gid              = nlist->gid;
426     shiftvec         = fr->shift_vec[0];
427     fshift           = fr->fshift[0];
428     facel            = fr->epsfac;
429     charge           = mdatoms->chargeA;
430     nvdwtype         = fr->ntype;
431     vdwparam         = fr->nbfp;
432     vdwtype          = mdatoms->typeA;
433
434     sh_ewald         = fr->ic->sh_ewald;
435     ewtab            = fr->ic->tabq_coul_F;
436     ewtabscale       = fr->ic->tabq_scale;
437     ewtabhalfspace   = 0.5/ewtabscale;
438
439     /* Setup water-specific parameters */
440     inr              = nlist->iinr[0];
441     iq0              = facel*charge[inr+0];
442     iq1              = facel*charge[inr+1];
443     iq2              = facel*charge[inr+2];
444     vdwioffset0      = 3*nvdwtype*vdwtype[inr+0];
445
446     outeriter        = 0;
447     inneriter        = 0;
448
449     /* Start outer loop over neighborlists */
450     for(iidx=0; iidx<nri; iidx++)
451     {
452         /* Load shift vector for this list */
453         i_shift_offset   = DIM*shiftidx[iidx];
454         shX              = shiftvec[i_shift_offset+XX];
455         shY              = shiftvec[i_shift_offset+YY];
456         shZ              = shiftvec[i_shift_offset+ZZ];
457
458         /* Load limits for loop over neighbors */
459         j_index_start    = jindex[iidx];
460         j_index_end      = jindex[iidx+1];
461
462         /* Get outer coordinate index */
463         inr              = iinr[iidx];
464         i_coord_offset   = DIM*inr;
465
466         /* Load i particle coords and add shift vector */
467         ix0              = shX + x[i_coord_offset+DIM*0+XX];
468         iy0              = shY + x[i_coord_offset+DIM*0+YY];
469         iz0              = shZ + x[i_coord_offset+DIM*0+ZZ];
470         ix1              = shX + x[i_coord_offset+DIM*1+XX];
471         iy1              = shY + x[i_coord_offset+DIM*1+YY];
472         iz1              = shZ + x[i_coord_offset+DIM*1+ZZ];
473         ix2              = shX + x[i_coord_offset+DIM*2+XX];
474         iy2              = shY + x[i_coord_offset+DIM*2+YY];
475         iz2              = shZ + x[i_coord_offset+DIM*2+ZZ];
476
477         fix0             = 0.0;
478         fiy0             = 0.0;
479         fiz0             = 0.0;
480         fix1             = 0.0;
481         fiy1             = 0.0;
482         fiz1             = 0.0;
483         fix2             = 0.0;
484         fiy2             = 0.0;
485         fiz2             = 0.0;
486
487         /* Start inner kernel loop */
488         for(jidx=j_index_start; jidx<j_index_end; jidx++)
489         {
490             /* Get j neighbor index, and coordinate index */
491             jnr              = jjnr[jidx];
492             j_coord_offset   = DIM*jnr;
493
494             /* load j atom coordinates */
495             jx0              = x[j_coord_offset+DIM*0+XX];
496             jy0              = x[j_coord_offset+DIM*0+YY];
497             jz0              = x[j_coord_offset+DIM*0+ZZ];
498
499             /* Calculate displacement vector */
500             dx00             = ix0 - jx0;
501             dy00             = iy0 - jy0;
502             dz00             = iz0 - jz0;
503             dx10             = ix1 - jx0;
504             dy10             = iy1 - jy0;
505             dz10             = iz1 - jz0;
506             dx20             = ix2 - jx0;
507             dy20             = iy2 - jy0;
508             dz20             = iz2 - jz0;
509
510             /* Calculate squared distance and things based on it */
511             rsq00            = dx00*dx00+dy00*dy00+dz00*dz00;
512             rsq10            = dx10*dx10+dy10*dy10+dz10*dz10;
513             rsq20            = dx20*dx20+dy20*dy20+dz20*dz20;
514
515             rinv00           = gmx_invsqrt(rsq00);
516             rinv10           = gmx_invsqrt(rsq10);
517             rinv20           = gmx_invsqrt(rsq20);
518
519             rinvsq00         = rinv00*rinv00;
520             rinvsq10         = rinv10*rinv10;
521             rinvsq20         = rinv20*rinv20;
522
523             /* Load parameters for j particles */
524             jq0              = charge[jnr+0];
525             vdwjidx0         = 3*vdwtype[jnr+0];
526
527             /**************************
528              * CALCULATE INTERACTIONS *
529              **************************/
530
531             r00              = rsq00*rinv00;
532
533             qq00             = iq0*jq0;
534             c6_00            = vdwparam[vdwioffset0+vdwjidx0];
535             cexp1_00         = vdwparam[vdwioffset0+vdwjidx0+1];
536             cexp2_00         = vdwparam[vdwioffset0+vdwjidx0+2];
537
538             /* EWALD ELECTROSTATICS */
539
540             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
541             ewrt             = r00*ewtabscale;
542             ewitab           = ewrt;
543             eweps            = ewrt-ewitab;
544             felec            = (1.0-eweps)*ewtab[ewitab]+eweps*ewtab[ewitab+1];
545             felec            = qq00*rinv00*(rinvsq00-felec);
546
547             /* BUCKINGHAM DISPERSION/REPULSION */
548             rinvsix          = rinvsq00*rinvsq00*rinvsq00;
549             vvdw6            = c6_00*rinvsix;
550             br               = cexp2_00*r00;
551             vvdwexp          = cexp1_00*exp(-br);
552             fvdw             = (br*vvdwexp-vvdw6)*rinvsq00;
553
554             fscal            = felec+fvdw;
555
556             /* Calculate temporary vectorial force */
557             tx               = fscal*dx00;
558             ty               = fscal*dy00;
559             tz               = fscal*dz00;
560
561             /* Update vectorial force */
562             fix0            += tx;
563             fiy0            += ty;
564             fiz0            += tz;
565             f[j_coord_offset+DIM*0+XX] -= tx;
566             f[j_coord_offset+DIM*0+YY] -= ty;
567             f[j_coord_offset+DIM*0+ZZ] -= tz;
568
569             /**************************
570              * CALCULATE INTERACTIONS *
571              **************************/
572
573             r10              = rsq10*rinv10;
574
575             qq10             = iq1*jq0;
576
577             /* EWALD ELECTROSTATICS */
578
579             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
580             ewrt             = r10*ewtabscale;
581             ewitab           = ewrt;
582             eweps            = ewrt-ewitab;
583             felec            = (1.0-eweps)*ewtab[ewitab]+eweps*ewtab[ewitab+1];
584             felec            = qq10*rinv10*(rinvsq10-felec);
585
586             fscal            = felec;
587
588             /* Calculate temporary vectorial force */
589             tx               = fscal*dx10;
590             ty               = fscal*dy10;
591             tz               = fscal*dz10;
592
593             /* Update vectorial force */
594             fix1            += tx;
595             fiy1            += ty;
596             fiz1            += tz;
597             f[j_coord_offset+DIM*0+XX] -= tx;
598             f[j_coord_offset+DIM*0+YY] -= ty;
599             f[j_coord_offset+DIM*0+ZZ] -= tz;
600
601             /**************************
602              * CALCULATE INTERACTIONS *
603              **************************/
604
605             r20              = rsq20*rinv20;
606
607             qq20             = iq2*jq0;
608
609             /* EWALD ELECTROSTATICS */
610
611             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
612             ewrt             = r20*ewtabscale;
613             ewitab           = ewrt;
614             eweps            = ewrt-ewitab;
615             felec            = (1.0-eweps)*ewtab[ewitab]+eweps*ewtab[ewitab+1];
616             felec            = qq20*rinv20*(rinvsq20-felec);
617
618             fscal            = felec;
619
620             /* Calculate temporary vectorial force */
621             tx               = fscal*dx20;
622             ty               = fscal*dy20;
623             tz               = fscal*dz20;
624
625             /* Update vectorial force */
626             fix2            += tx;
627             fiy2            += ty;
628             fiz2            += 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             /* Inner loop uses 137 flops */
634         }
635         /* End of innermost loop */
636
637         tx = ty = tz = 0;
638         f[i_coord_offset+DIM*0+XX] += fix0;
639         f[i_coord_offset+DIM*0+YY] += fiy0;
640         f[i_coord_offset+DIM*0+ZZ] += fiz0;
641         tx                         += fix0;
642         ty                         += fiy0;
643         tz                         += fiz0;
644         f[i_coord_offset+DIM*1+XX] += fix1;
645         f[i_coord_offset+DIM*1+YY] += fiy1;
646         f[i_coord_offset+DIM*1+ZZ] += fiz1;
647         tx                         += fix1;
648         ty                         += fiy1;
649         tz                         += fiz1;
650         f[i_coord_offset+DIM*2+XX] += fix2;
651         f[i_coord_offset+DIM*2+YY] += fiy2;
652         f[i_coord_offset+DIM*2+ZZ] += fiz2;
653         tx                         += fix2;
654         ty                         += fiy2;
655         tz                         += fiz2;
656         fshift[i_shift_offset+XX]  += tx;
657         fshift[i_shift_offset+YY]  += ty;
658         fshift[i_shift_offset+ZZ]  += tz;
659
660         /* Increment number of inner iterations */
661         inneriter                  += j_index_end - j_index_start;
662
663         /* Outer loop uses 30 flops */
664     }
665
666     /* Increment number of outer iterations */
667     outeriter        += nri;
668
669     /* Update outer/inner flops */
670
671     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W3_F,outeriter*30 + inneriter*137);
672 }