956f68d0c4415fc1bd73dc7b0451904120992870
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_c / nb_kernel_ElecEwSh_VdwBhamSh_GeomW3P1_c.c
1 /*
2  * Note: this file was generated by the Gromacs c kernel generator.
3  *
4  *                This source code is part of
5  *
6  *                 G   R   O   M   A   C   S
7  *
8  * Copyright (c) 2001-2012, The GROMACS Development Team
9  *
10  * Gromacs is a library for molecular simulation and trajectory analysis,
11  * written by Erik Lindahl, David van der Spoel, Berk Hess, and others - for
12  * a full list of developers and information, check out http://www.gromacs.org
13  *
14  * This program is free software; you can redistribute it and/or modify it under
15  * the terms of the GNU Lesser General Public License as published by the Free
16  * Software Foundation; either version 2 of the License, or (at your option) any
17  * later version.
18  *
19  * To help fund GROMACS development, we humbly ask that you cite
20  * the papers people have written on it - you can find them on the website.
21  */
22 #ifdef HAVE_CONFIG_H
23 #include <config.h>
24 #endif
25
26 #include <math.h>
27
28 #include "../nb_kernel.h"
29 #include "types/simple.h"
30 #include "vec.h"
31 #include "nrnb.h"
32
33 /*
34  * Gromacs nonbonded kernel:   nb_kernel_ElecEwSh_VdwBhamSh_GeomW3P1_VF_c
35  * Electrostatics interaction: Ewald
36  * VdW interaction:            Buckingham
37  * Geometry:                   Water3-Particle
38  * Calculate force/pot:        PotentialAndForce
39  */
40 void
41 nb_kernel_ElecEwSh_VdwBhamSh_GeomW3P1_VF_c
42                     (t_nblist * gmx_restrict                nlist,
43                      rvec * gmx_restrict                    xx,
44                      rvec * gmx_restrict                    ff,
45                      t_forcerec * gmx_restrict              fr,
46                      t_mdatoms * gmx_restrict               mdatoms,
47                      nb_kernel_data_t * gmx_restrict        kernel_data,
48                      t_nrnb * gmx_restrict                  nrnb)
49 {
50     int              i_shift_offset,i_coord_offset,j_coord_offset;
51     int              j_index_start,j_index_end;
52     int              nri,inr,ggid,iidx,jidx,jnr,outeriter,inneriter;
53     real             shX,shY,shZ,tx,ty,tz,fscal,rcutoff,rcutoff2;
54     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
55     real             *shiftvec,*fshift,*x,*f;
56     int              vdwioffset0;
57     real             ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
58     int              vdwioffset1;
59     real             ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
60     int              vdwioffset2;
61     real             ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
62     int              vdwjidx0;
63     real             jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
64     real             dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00,cexp1_00,cexp2_00;
65     real             dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10,cexp1_10,cexp2_10;
66     real             dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20,cexp1_20,cexp2_20;
67     real             velec,felec,velecsum,facel,crf,krf,krf2;
68     real             *charge;
69     int              nvdwtype;
70     real             rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,br,vvdwexp,sh_vdw_invrcut6;
71     int              *vdwtype;
72     real             *vdwparam;
73     int              ewitab;
74     real             ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace;
75     real             *ewtab;
76
77     x                = xx[0];
78     f                = ff[0];
79
80     nri              = nlist->nri;
81     iinr             = nlist->iinr;
82     jindex           = nlist->jindex;
83     jjnr             = nlist->jjnr;
84     shiftidx         = nlist->shift;
85     gid              = nlist->gid;
86     shiftvec         = fr->shift_vec[0];
87     fshift           = fr->fshift[0];
88     facel            = fr->epsfac;
89     charge           = mdatoms->chargeA;
90     nvdwtype         = fr->ntype;
91     vdwparam         = fr->nbfp;
92     vdwtype          = mdatoms->typeA;
93
94     sh_ewald         = fr->ic->sh_ewald;
95     ewtab            = fr->ic->tabq_coul_FDV0;
96     ewtabscale       = fr->ic->tabq_scale;
97     ewtabhalfspace   = 0.5/ewtabscale;
98
99     /* Setup water-specific parameters */
100     inr              = nlist->iinr[0];
101     iq0              = facel*charge[inr+0];
102     iq1              = facel*charge[inr+1];
103     iq2              = facel*charge[inr+2];
104     vdwioffset0      = 3*nvdwtype*vdwtype[inr+0];
105
106     /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
107     rcutoff          = fr->rcoulomb;
108     rcutoff2         = rcutoff*rcutoff;
109
110     sh_vdw_invrcut6  = fr->ic->sh_invrc6;
111     rvdw             = fr->rvdw;
112
113     outeriter        = 0;
114     inneriter        = 0;
115
116     /* Start outer loop over neighborlists */
117     for(iidx=0; iidx<nri; iidx++)
118     {
119         /* Load shift vector for this list */
120         i_shift_offset   = DIM*shiftidx[iidx];
121         shX              = shiftvec[i_shift_offset+XX];
122         shY              = shiftvec[i_shift_offset+YY];
123         shZ              = shiftvec[i_shift_offset+ZZ];
124
125         /* Load limits for loop over neighbors */
126         j_index_start    = jindex[iidx];
127         j_index_end      = jindex[iidx+1];
128
129         /* Get outer coordinate index */
130         inr              = iinr[iidx];
131         i_coord_offset   = DIM*inr;
132
133         /* Load i particle coords and add shift vector */
134         ix0              = shX + x[i_coord_offset+DIM*0+XX];
135         iy0              = shY + x[i_coord_offset+DIM*0+YY];
136         iz0              = shZ + x[i_coord_offset+DIM*0+ZZ];
137         ix1              = shX + x[i_coord_offset+DIM*1+XX];
138         iy1              = shY + x[i_coord_offset+DIM*1+YY];
139         iz1              = shZ + x[i_coord_offset+DIM*1+ZZ];
140         ix2              = shX + x[i_coord_offset+DIM*2+XX];
141         iy2              = shY + x[i_coord_offset+DIM*2+YY];
142         iz2              = shZ + x[i_coord_offset+DIM*2+ZZ];
143
144         fix0             = 0.0;
145         fiy0             = 0.0;
146         fiz0             = 0.0;
147         fix1             = 0.0;
148         fiy1             = 0.0;
149         fiz1             = 0.0;
150         fix2             = 0.0;
151         fiy2             = 0.0;
152         fiz2             = 0.0;
153
154         /* Reset potential sums */
155         velecsum         = 0.0;
156         vvdwsum          = 0.0;
157
158         /* Start inner kernel loop */
159         for(jidx=j_index_start; jidx<j_index_end; jidx++)
160         {
161             /* Get j neighbor index, and coordinate index */
162             jnr              = jjnr[jidx];
163             j_coord_offset   = DIM*jnr;
164
165             /* load j atom coordinates */
166             jx0              = x[j_coord_offset+DIM*0+XX];
167             jy0              = x[j_coord_offset+DIM*0+YY];
168             jz0              = x[j_coord_offset+DIM*0+ZZ];
169
170             /* Calculate displacement vector */
171             dx00             = ix0 - jx0;
172             dy00             = iy0 - jy0;
173             dz00             = iz0 - jz0;
174             dx10             = ix1 - jx0;
175             dy10             = iy1 - jy0;
176             dz10             = iz1 - jz0;
177             dx20             = ix2 - jx0;
178             dy20             = iy2 - jy0;
179             dz20             = iz2 - jz0;
180
181             /* Calculate squared distance and things based on it */
182             rsq00            = dx00*dx00+dy00*dy00+dz00*dz00;
183             rsq10            = dx10*dx10+dy10*dy10+dz10*dz10;
184             rsq20            = dx20*dx20+dy20*dy20+dz20*dz20;
185
186             rinv00           = gmx_invsqrt(rsq00);
187             rinv10           = gmx_invsqrt(rsq10);
188             rinv20           = gmx_invsqrt(rsq20);
189
190             rinvsq00         = rinv00*rinv00;
191             rinvsq10         = rinv10*rinv10;
192             rinvsq20         = rinv20*rinv20;
193
194             /* Load parameters for j particles */
195             jq0              = charge[jnr+0];
196             vdwjidx0         = 3*vdwtype[jnr+0];
197
198             /**************************
199              * CALCULATE INTERACTIONS *
200              **************************/
201
202             if (rsq00<rcutoff2)
203             {
204
205             r00              = rsq00*rinv00;
206
207             qq00             = iq0*jq0;
208             c6_00            = vdwparam[vdwioffset0+vdwjidx0];
209             cexp1_00         = vdwparam[vdwioffset0+vdwjidx0+1];
210             cexp2_00         = vdwparam[vdwioffset0+vdwjidx0+2];
211
212             /* EWALD ELECTROSTATICS */
213
214             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
215             ewrt             = r00*ewtabscale;
216             ewitab           = ewrt;
217             eweps            = ewrt-ewitab;
218             ewitab           = 4*ewitab;
219             felec            = ewtab[ewitab]+eweps*ewtab[ewitab+1];
220             velec            = qq00*((rinv00-sh_ewald)-(ewtab[ewitab+2]-ewtabhalfspace*eweps*(ewtab[ewitab]+felec)));
221             felec            = qq00*rinv00*(rinvsq00-felec);
222
223             /* BUCKINGHAM DISPERSION/REPULSION */
224             rinvsix          = rinvsq00*rinvsq00*rinvsq00;
225             vvdw6            = c6_00*rinvsix;
226             br               = cexp2_00*r00;
227             vvdwexp          = cexp1_00*exp(-br);
228             vvdw             = (vvdwexp-cexp1_00*exp(-cexp2_00*rvdw)) - (vvdw6 - c6_00*sh_vdw_invrcut6)*(1.0/6.0);
229             fvdw             = (br*vvdwexp-vvdw6)*rinvsq00;
230
231             /* Update potential sums from outer loop */
232             velecsum        += velec;
233             vvdwsum         += vvdw;
234
235             fscal            = felec+fvdw;
236
237             /* Calculate temporary vectorial force */
238             tx               = fscal*dx00;
239             ty               = fscal*dy00;
240             tz               = fscal*dz00;
241
242             /* Update vectorial force */
243             fix0            += tx;
244             fiy0            += ty;
245             fiz0            += tz;
246             f[j_coord_offset+DIM*0+XX] -= tx;
247             f[j_coord_offset+DIM*0+YY] -= ty;
248             f[j_coord_offset+DIM*0+ZZ] -= tz;
249
250             }
251
252             /**************************
253              * CALCULATE INTERACTIONS *
254              **************************/
255
256             if (rsq10<rcutoff2)
257             {
258
259             r10              = rsq10*rinv10;
260
261             qq10             = iq1*jq0;
262
263             /* EWALD ELECTROSTATICS */
264
265             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
266             ewrt             = r10*ewtabscale;
267             ewitab           = ewrt;
268             eweps            = ewrt-ewitab;
269             ewitab           = 4*ewitab;
270             felec            = ewtab[ewitab]+eweps*ewtab[ewitab+1];
271             velec            = qq10*((rinv10-sh_ewald)-(ewtab[ewitab+2]-ewtabhalfspace*eweps*(ewtab[ewitab]+felec)));
272             felec            = qq10*rinv10*(rinvsq10-felec);
273
274             /* Update potential sums from outer loop */
275             velecsum        += velec;
276
277             fscal            = felec;
278
279             /* Calculate temporary vectorial force */
280             tx               = fscal*dx10;
281             ty               = fscal*dy10;
282             tz               = fscal*dz10;
283
284             /* Update vectorial force */
285             fix1            += tx;
286             fiy1            += ty;
287             fiz1            += 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             }
293
294             /**************************
295              * CALCULATE INTERACTIONS *
296              **************************/
297
298             if (rsq20<rcutoff2)
299             {
300
301             r20              = rsq20*rinv20;
302
303             qq20             = iq2*jq0;
304
305             /* EWALD ELECTROSTATICS */
306
307             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
308             ewrt             = r20*ewtabscale;
309             ewitab           = ewrt;
310             eweps            = ewrt-ewitab;
311             ewitab           = 4*ewitab;
312             felec            = ewtab[ewitab]+eweps*ewtab[ewitab+1];
313             velec            = qq20*((rinv20-sh_ewald)-(ewtab[ewitab+2]-ewtabhalfspace*eweps*(ewtab[ewitab]+felec)));
314             felec            = qq20*rinv20*(rinvsq20-felec);
315
316             /* Update potential sums from outer loop */
317             velecsum        += velec;
318
319             fscal            = felec;
320
321             /* Calculate temporary vectorial force */
322             tx               = fscal*dx20;
323             ty               = fscal*dy20;
324             tz               = fscal*dz20;
325
326             /* Update vectorial force */
327             fix2            += tx;
328             fiy2            += ty;
329             fiz2            += tz;
330             f[j_coord_offset+DIM*0+XX] -= tx;
331             f[j_coord_offset+DIM*0+YY] -= ty;
332             f[j_coord_offset+DIM*0+ZZ] -= tz;
333
334             }
335
336             /* Inner loop uses 195 flops */
337         }
338         /* End of innermost loop */
339
340         tx = ty = tz = 0;
341         f[i_coord_offset+DIM*0+XX] += fix0;
342         f[i_coord_offset+DIM*0+YY] += fiy0;
343         f[i_coord_offset+DIM*0+ZZ] += fiz0;
344         tx                         += fix0;
345         ty                         += fiy0;
346         tz                         += fiz0;
347         f[i_coord_offset+DIM*1+XX] += fix1;
348         f[i_coord_offset+DIM*1+YY] += fiy1;
349         f[i_coord_offset+DIM*1+ZZ] += fiz1;
350         tx                         += fix1;
351         ty                         += fiy1;
352         tz                         += fiz1;
353         f[i_coord_offset+DIM*2+XX] += fix2;
354         f[i_coord_offset+DIM*2+YY] += fiy2;
355         f[i_coord_offset+DIM*2+ZZ] += fiz2;
356         tx                         += fix2;
357         ty                         += fiy2;
358         tz                         += fiz2;
359         fshift[i_shift_offset+XX]  += tx;
360         fshift[i_shift_offset+YY]  += ty;
361         fshift[i_shift_offset+ZZ]  += tz;
362
363         ggid                        = gid[iidx];
364         /* Update potential energies */
365         kernel_data->energygrp_elec[ggid] += velecsum;
366         kernel_data->energygrp_vdw[ggid] += vvdwsum;
367
368         /* Increment number of inner iterations */
369         inneriter                  += j_index_end - j_index_start;
370
371         /* Outer loop uses 32 flops */
372     }
373
374     /* Increment number of outer iterations */
375     outeriter        += nri;
376
377     /* Update outer/inner flops */
378
379     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W3_VF,outeriter*32 + inneriter*195);
380 }
381 /*
382  * Gromacs nonbonded kernel:   nb_kernel_ElecEwSh_VdwBhamSh_GeomW3P1_F_c
383  * Electrostatics interaction: Ewald
384  * VdW interaction:            Buckingham
385  * Geometry:                   Water3-Particle
386  * Calculate force/pot:        Force
387  */
388 void
389 nb_kernel_ElecEwSh_VdwBhamSh_GeomW3P1_F_c
390                     (t_nblist * gmx_restrict                nlist,
391                      rvec * gmx_restrict                    xx,
392                      rvec * gmx_restrict                    ff,
393                      t_forcerec * gmx_restrict              fr,
394                      t_mdatoms * gmx_restrict               mdatoms,
395                      nb_kernel_data_t * gmx_restrict        kernel_data,
396                      t_nrnb * gmx_restrict                  nrnb)
397 {
398     int              i_shift_offset,i_coord_offset,j_coord_offset;
399     int              j_index_start,j_index_end;
400     int              nri,inr,ggid,iidx,jidx,jnr,outeriter,inneriter;
401     real             shX,shY,shZ,tx,ty,tz,fscal,rcutoff,rcutoff2;
402     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
403     real             *shiftvec,*fshift,*x,*f;
404     int              vdwioffset0;
405     real             ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
406     int              vdwioffset1;
407     real             ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
408     int              vdwioffset2;
409     real             ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
410     int              vdwjidx0;
411     real             jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
412     real             dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00,cexp1_00,cexp2_00;
413     real             dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10,cexp1_10,cexp2_10;
414     real             dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20,cexp1_20,cexp2_20;
415     real             velec,felec,velecsum,facel,crf,krf,krf2;
416     real             *charge;
417     int              nvdwtype;
418     real             rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,br,vvdwexp,sh_vdw_invrcut6;
419     int              *vdwtype;
420     real             *vdwparam;
421     int              ewitab;
422     real             ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace;
423     real             *ewtab;
424
425     x                = xx[0];
426     f                = ff[0];
427
428     nri              = nlist->nri;
429     iinr             = nlist->iinr;
430     jindex           = nlist->jindex;
431     jjnr             = nlist->jjnr;
432     shiftidx         = nlist->shift;
433     gid              = nlist->gid;
434     shiftvec         = fr->shift_vec[0];
435     fshift           = fr->fshift[0];
436     facel            = fr->epsfac;
437     charge           = mdatoms->chargeA;
438     nvdwtype         = fr->ntype;
439     vdwparam         = fr->nbfp;
440     vdwtype          = mdatoms->typeA;
441
442     sh_ewald         = fr->ic->sh_ewald;
443     ewtab            = fr->ic->tabq_coul_F;
444     ewtabscale       = fr->ic->tabq_scale;
445     ewtabhalfspace   = 0.5/ewtabscale;
446
447     /* Setup water-specific parameters */
448     inr              = nlist->iinr[0];
449     iq0              = facel*charge[inr+0];
450     iq1              = facel*charge[inr+1];
451     iq2              = facel*charge[inr+2];
452     vdwioffset0      = 3*nvdwtype*vdwtype[inr+0];
453
454     /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
455     rcutoff          = fr->rcoulomb;
456     rcutoff2         = rcutoff*rcutoff;
457
458     sh_vdw_invrcut6  = fr->ic->sh_invrc6;
459     rvdw             = fr->rvdw;
460
461     outeriter        = 0;
462     inneriter        = 0;
463
464     /* Start outer loop over neighborlists */
465     for(iidx=0; iidx<nri; iidx++)
466     {
467         /* Load shift vector for this list */
468         i_shift_offset   = DIM*shiftidx[iidx];
469         shX              = shiftvec[i_shift_offset+XX];
470         shY              = shiftvec[i_shift_offset+YY];
471         shZ              = shiftvec[i_shift_offset+ZZ];
472
473         /* Load limits for loop over neighbors */
474         j_index_start    = jindex[iidx];
475         j_index_end      = jindex[iidx+1];
476
477         /* Get outer coordinate index */
478         inr              = iinr[iidx];
479         i_coord_offset   = DIM*inr;
480
481         /* Load i particle coords and add shift vector */
482         ix0              = shX + x[i_coord_offset+DIM*0+XX];
483         iy0              = shY + x[i_coord_offset+DIM*0+YY];
484         iz0              = shZ + x[i_coord_offset+DIM*0+ZZ];
485         ix1              = shX + x[i_coord_offset+DIM*1+XX];
486         iy1              = shY + x[i_coord_offset+DIM*1+YY];
487         iz1              = shZ + x[i_coord_offset+DIM*1+ZZ];
488         ix2              = shX + x[i_coord_offset+DIM*2+XX];
489         iy2              = shY + x[i_coord_offset+DIM*2+YY];
490         iz2              = shZ + x[i_coord_offset+DIM*2+ZZ];
491
492         fix0             = 0.0;
493         fiy0             = 0.0;
494         fiz0             = 0.0;
495         fix1             = 0.0;
496         fiy1             = 0.0;
497         fiz1             = 0.0;
498         fix2             = 0.0;
499         fiy2             = 0.0;
500         fiz2             = 0.0;
501
502         /* Start inner kernel loop */
503         for(jidx=j_index_start; jidx<j_index_end; jidx++)
504         {
505             /* Get j neighbor index, and coordinate index */
506             jnr              = jjnr[jidx];
507             j_coord_offset   = DIM*jnr;
508
509             /* load j atom coordinates */
510             jx0              = x[j_coord_offset+DIM*0+XX];
511             jy0              = x[j_coord_offset+DIM*0+YY];
512             jz0              = x[j_coord_offset+DIM*0+ZZ];
513
514             /* Calculate displacement vector */
515             dx00             = ix0 - jx0;
516             dy00             = iy0 - jy0;
517             dz00             = iz0 - jz0;
518             dx10             = ix1 - jx0;
519             dy10             = iy1 - jy0;
520             dz10             = iz1 - jz0;
521             dx20             = ix2 - jx0;
522             dy20             = iy2 - jy0;
523             dz20             = iz2 - jz0;
524
525             /* Calculate squared distance and things based on it */
526             rsq00            = dx00*dx00+dy00*dy00+dz00*dz00;
527             rsq10            = dx10*dx10+dy10*dy10+dz10*dz10;
528             rsq20            = dx20*dx20+dy20*dy20+dz20*dz20;
529
530             rinv00           = gmx_invsqrt(rsq00);
531             rinv10           = gmx_invsqrt(rsq10);
532             rinv20           = gmx_invsqrt(rsq20);
533
534             rinvsq00         = rinv00*rinv00;
535             rinvsq10         = rinv10*rinv10;
536             rinvsq20         = rinv20*rinv20;
537
538             /* Load parameters for j particles */
539             jq0              = charge[jnr+0];
540             vdwjidx0         = 3*vdwtype[jnr+0];
541
542             /**************************
543              * CALCULATE INTERACTIONS *
544              **************************/
545
546             if (rsq00<rcutoff2)
547             {
548
549             r00              = rsq00*rinv00;
550
551             qq00             = iq0*jq0;
552             c6_00            = vdwparam[vdwioffset0+vdwjidx0];
553             cexp1_00         = vdwparam[vdwioffset0+vdwjidx0+1];
554             cexp2_00         = vdwparam[vdwioffset0+vdwjidx0+2];
555
556             /* EWALD ELECTROSTATICS */
557
558             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
559             ewrt             = r00*ewtabscale;
560             ewitab           = ewrt;
561             eweps            = ewrt-ewitab;
562             felec            = (1.0-eweps)*ewtab[ewitab]+eweps*ewtab[ewitab+1];
563             felec            = qq00*rinv00*(rinvsq00-felec);
564
565             /* BUCKINGHAM DISPERSION/REPULSION */
566             rinvsix          = rinvsq00*rinvsq00*rinvsq00;
567             vvdw6            = c6_00*rinvsix;
568             br               = cexp2_00*r00;
569             vvdwexp          = cexp1_00*exp(-br);
570             fvdw             = (br*vvdwexp-vvdw6)*rinvsq00;
571
572             fscal            = felec+fvdw;
573
574             /* Calculate temporary vectorial force */
575             tx               = fscal*dx00;
576             ty               = fscal*dy00;
577             tz               = fscal*dz00;
578
579             /* Update vectorial force */
580             fix0            += tx;
581             fiy0            += ty;
582             fiz0            += tz;
583             f[j_coord_offset+DIM*0+XX] -= tx;
584             f[j_coord_offset+DIM*0+YY] -= ty;
585             f[j_coord_offset+DIM*0+ZZ] -= tz;
586
587             }
588
589             /**************************
590              * CALCULATE INTERACTIONS *
591              **************************/
592
593             if (rsq10<rcutoff2)
594             {
595
596             r10              = rsq10*rinv10;
597
598             qq10             = iq1*jq0;
599
600             /* EWALD ELECTROSTATICS */
601
602             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
603             ewrt             = r10*ewtabscale;
604             ewitab           = ewrt;
605             eweps            = ewrt-ewitab;
606             felec            = (1.0-eweps)*ewtab[ewitab]+eweps*ewtab[ewitab+1];
607             felec            = qq10*rinv10*(rinvsq10-felec);
608
609             fscal            = felec;
610
611             /* Calculate temporary vectorial force */
612             tx               = fscal*dx10;
613             ty               = fscal*dy10;
614             tz               = fscal*dz10;
615
616             /* Update vectorial force */
617             fix1            += tx;
618             fiy1            += ty;
619             fiz1            += tz;
620             f[j_coord_offset+DIM*0+XX] -= tx;
621             f[j_coord_offset+DIM*0+YY] -= ty;
622             f[j_coord_offset+DIM*0+ZZ] -= tz;
623
624             }
625
626             /**************************
627              * CALCULATE INTERACTIONS *
628              **************************/
629
630             if (rsq20<rcutoff2)
631             {
632
633             r20              = rsq20*rinv20;
634
635             qq20             = iq2*jq0;
636
637             /* EWALD ELECTROSTATICS */
638
639             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
640             ewrt             = r20*ewtabscale;
641             ewitab           = ewrt;
642             eweps            = ewrt-ewitab;
643             felec            = (1.0-eweps)*ewtab[ewitab]+eweps*ewtab[ewitab+1];
644             felec            = qq20*rinv20*(rinvsq20-felec);
645
646             fscal            = felec;
647
648             /* Calculate temporary vectorial force */
649             tx               = fscal*dx20;
650             ty               = fscal*dy20;
651             tz               = fscal*dz20;
652
653             /* Update vectorial force */
654             fix2            += tx;
655             fiy2            += ty;
656             fiz2            += tz;
657             f[j_coord_offset+DIM*0+XX] -= tx;
658             f[j_coord_offset+DIM*0+YY] -= ty;
659             f[j_coord_offset+DIM*0+ZZ] -= tz;
660
661             }
662
663             /* Inner loop uses 137 flops */
664         }
665         /* End of innermost loop */
666
667         tx = ty = tz = 0;
668         f[i_coord_offset+DIM*0+XX] += fix0;
669         f[i_coord_offset+DIM*0+YY] += fiy0;
670         f[i_coord_offset+DIM*0+ZZ] += fiz0;
671         tx                         += fix0;
672         ty                         += fiy0;
673         tz                         += fiz0;
674         f[i_coord_offset+DIM*1+XX] += fix1;
675         f[i_coord_offset+DIM*1+YY] += fiy1;
676         f[i_coord_offset+DIM*1+ZZ] += fiz1;
677         tx                         += fix1;
678         ty                         += fiy1;
679         tz                         += fiz1;
680         f[i_coord_offset+DIM*2+XX] += fix2;
681         f[i_coord_offset+DIM*2+YY] += fiy2;
682         f[i_coord_offset+DIM*2+ZZ] += fiz2;
683         tx                         += fix2;
684         ty                         += fiy2;
685         tz                         += fiz2;
686         fshift[i_shift_offset+XX]  += tx;
687         fshift[i_shift_offset+YY]  += ty;
688         fshift[i_shift_offset+ZZ]  += tz;
689
690         /* Increment number of inner iterations */
691         inneriter                  += j_index_end - j_index_start;
692
693         /* Outer loop uses 30 flops */
694     }
695
696     /* Increment number of outer iterations */
697     outeriter        += nri;
698
699     /* Update outer/inner flops */
700
701     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W3_F,outeriter*30 + inneriter*137);
702 }