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