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