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