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