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