68627d6cf6b89df3ab3d96027978f670d49c05f3
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_c / nb_kernel_ElecCSTab_VdwCSTab_GeomW4W4_c.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2013,2014, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35 /*
36  * Note: this file was generated by the GROMACS c kernel generator.
37  */
38 #ifdef HAVE_CONFIG_H
39 #include <config.h>
40 #endif
41
42 #include <math.h>
43
44 #include "../nb_kernel.h"
45 #include "types/simple.h"
46 #include "gromacs/math/vec.h"
47 #include "nrnb.h"
48
49 /*
50  * Gromacs nonbonded kernel:   nb_kernel_ElecCSTab_VdwCSTab_GeomW4W4_VF_c
51  * Electrostatics interaction: CubicSplineTable
52  * VdW interaction:            CubicSplineTable
53  * Geometry:                   Water4-Water4
54  * Calculate force/pot:        PotentialAndForce
55  */
56 void
57 nb_kernel_ElecCSTab_VdwCSTab_GeomW4W4_VF_c
58                     (t_nblist                    * gmx_restrict       nlist,
59                      rvec                        * gmx_restrict          xx,
60                      rvec                        * gmx_restrict          ff,
61                      t_forcerec                  * gmx_restrict          fr,
62                      t_mdatoms                   * gmx_restrict     mdatoms,
63                      nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
64                      t_nrnb                      * gmx_restrict        nrnb)
65 {
66     int              i_shift_offset,i_coord_offset,j_coord_offset;
67     int              j_index_start,j_index_end;
68     int              nri,inr,ggid,iidx,jidx,jnr,outeriter,inneriter;
69     real             shX,shY,shZ,tx,ty,tz,fscal,rcutoff,rcutoff2;
70     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
71     real             *shiftvec,*fshift,*x,*f;
72     int              vdwioffset0;
73     real             ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
74     int              vdwioffset1;
75     real             ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
76     int              vdwioffset2;
77     real             ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
78     int              vdwioffset3;
79     real             ix3,iy3,iz3,fix3,fiy3,fiz3,iq3,isai3;
80     int              vdwjidx0;
81     real             jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
82     int              vdwjidx1;
83     real             jx1,jy1,jz1,fjx1,fjy1,fjz1,jq1,isaj1;
84     int              vdwjidx2;
85     real             jx2,jy2,jz2,fjx2,fjy2,fjz2,jq2,isaj2;
86     int              vdwjidx3;
87     real             jx3,jy3,jz3,fjx3,fjy3,fjz3,jq3,isaj3;
88     real             dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00,cexp1_00,cexp2_00;
89     real             dx11,dy11,dz11,rsq11,rinv11,rinvsq11,r11,qq11,c6_11,c12_11,cexp1_11,cexp2_11;
90     real             dx12,dy12,dz12,rsq12,rinv12,rinvsq12,r12,qq12,c6_12,c12_12,cexp1_12,cexp2_12;
91     real             dx13,dy13,dz13,rsq13,rinv13,rinvsq13,r13,qq13,c6_13,c12_13,cexp1_13,cexp2_13;
92     real             dx21,dy21,dz21,rsq21,rinv21,rinvsq21,r21,qq21,c6_21,c12_21,cexp1_21,cexp2_21;
93     real             dx22,dy22,dz22,rsq22,rinv22,rinvsq22,r22,qq22,c6_22,c12_22,cexp1_22,cexp2_22;
94     real             dx23,dy23,dz23,rsq23,rinv23,rinvsq23,r23,qq23,c6_23,c12_23,cexp1_23,cexp2_23;
95     real             dx31,dy31,dz31,rsq31,rinv31,rinvsq31,r31,qq31,c6_31,c12_31,cexp1_31,cexp2_31;
96     real             dx32,dy32,dz32,rsq32,rinv32,rinvsq32,r32,qq32,c6_32,c12_32,cexp1_32,cexp2_32;
97     real             dx33,dy33,dz33,rsq33,rinv33,rinvsq33,r33,qq33,c6_33,c12_33,cexp1_33,cexp2_33;
98     real             velec,felec,velecsum,facel,crf,krf,krf2;
99     real             *charge;
100     int              nvdwtype;
101     real             rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,br,vvdwexp,sh_vdw_invrcut6;
102     int              *vdwtype;
103     real             *vdwparam;
104     int              vfitab;
105     real             rt,vfeps,vftabscale,Y,F,Geps,Heps2,Fp,VV,FF;
106     real             *vftab;
107
108     x                = xx[0];
109     f                = ff[0];
110
111     nri              = nlist->nri;
112     iinr             = nlist->iinr;
113     jindex           = nlist->jindex;
114     jjnr             = nlist->jjnr;
115     shiftidx         = nlist->shift;
116     gid              = nlist->gid;
117     shiftvec         = fr->shift_vec[0];
118     fshift           = fr->fshift[0];
119     facel            = fr->epsfac;
120     charge           = mdatoms->chargeA;
121     nvdwtype         = fr->ntype;
122     vdwparam         = fr->nbfp;
123     vdwtype          = mdatoms->typeA;
124
125     vftab            = kernel_data->table_elec_vdw->data;
126     vftabscale       = kernel_data->table_elec_vdw->scale;
127
128     /* Setup water-specific parameters */
129     inr              = nlist->iinr[0];
130     iq1              = facel*charge[inr+1];
131     iq2              = facel*charge[inr+2];
132     iq3              = facel*charge[inr+3];
133     vdwioffset0      = 2*nvdwtype*vdwtype[inr+0];
134
135     jq1              = charge[inr+1];
136     jq2              = charge[inr+2];
137     jq3              = charge[inr+3];
138     vdwjidx0         = 2*vdwtype[inr+0];
139     c6_00            = vdwparam[vdwioffset0+vdwjidx0];
140     c12_00           = vdwparam[vdwioffset0+vdwjidx0+1];
141     qq11             = iq1*jq1;
142     qq12             = iq1*jq2;
143     qq13             = iq1*jq3;
144     qq21             = iq2*jq1;
145     qq22             = iq2*jq2;
146     qq23             = iq2*jq3;
147     qq31             = iq3*jq1;
148     qq32             = iq3*jq2;
149     qq33             = iq3*jq3;
150
151     outeriter        = 0;
152     inneriter        = 0;
153
154     /* Start outer loop over neighborlists */
155     for(iidx=0; iidx<nri; iidx++)
156     {
157         /* Load shift vector for this list */
158         i_shift_offset   = DIM*shiftidx[iidx];
159         shX              = shiftvec[i_shift_offset+XX];
160         shY              = shiftvec[i_shift_offset+YY];
161         shZ              = shiftvec[i_shift_offset+ZZ];
162
163         /* Load limits for loop over neighbors */
164         j_index_start    = jindex[iidx];
165         j_index_end      = jindex[iidx+1];
166
167         /* Get outer coordinate index */
168         inr              = iinr[iidx];
169         i_coord_offset   = DIM*inr;
170
171         /* Load i particle coords and add shift vector */
172         ix0              = shX + x[i_coord_offset+DIM*0+XX];
173         iy0              = shY + x[i_coord_offset+DIM*0+YY];
174         iz0              = shZ + x[i_coord_offset+DIM*0+ZZ];
175         ix1              = shX + x[i_coord_offset+DIM*1+XX];
176         iy1              = shY + x[i_coord_offset+DIM*1+YY];
177         iz1              = shZ + x[i_coord_offset+DIM*1+ZZ];
178         ix2              = shX + x[i_coord_offset+DIM*2+XX];
179         iy2              = shY + x[i_coord_offset+DIM*2+YY];
180         iz2              = shZ + x[i_coord_offset+DIM*2+ZZ];
181         ix3              = shX + x[i_coord_offset+DIM*3+XX];
182         iy3              = shY + x[i_coord_offset+DIM*3+YY];
183         iz3              = shZ + x[i_coord_offset+DIM*3+ZZ];
184
185         fix0             = 0.0;
186         fiy0             = 0.0;
187         fiz0             = 0.0;
188         fix1             = 0.0;
189         fiy1             = 0.0;
190         fiz1             = 0.0;
191         fix2             = 0.0;
192         fiy2             = 0.0;
193         fiz2             = 0.0;
194         fix3             = 0.0;
195         fiy3             = 0.0;
196         fiz3             = 0.0;
197
198         /* Reset potential sums */
199         velecsum         = 0.0;
200         vvdwsum          = 0.0;
201
202         /* Start inner kernel loop */
203         for(jidx=j_index_start; jidx<j_index_end; jidx++)
204         {
205             /* Get j neighbor index, and coordinate index */
206             jnr              = jjnr[jidx];
207             j_coord_offset   = DIM*jnr;
208
209             /* load j atom coordinates */
210             jx0              = x[j_coord_offset+DIM*0+XX];
211             jy0              = x[j_coord_offset+DIM*0+YY];
212             jz0              = x[j_coord_offset+DIM*0+ZZ];
213             jx1              = x[j_coord_offset+DIM*1+XX];
214             jy1              = x[j_coord_offset+DIM*1+YY];
215             jz1              = x[j_coord_offset+DIM*1+ZZ];
216             jx2              = x[j_coord_offset+DIM*2+XX];
217             jy2              = x[j_coord_offset+DIM*2+YY];
218             jz2              = x[j_coord_offset+DIM*2+ZZ];
219             jx3              = x[j_coord_offset+DIM*3+XX];
220             jy3              = x[j_coord_offset+DIM*3+YY];
221             jz3              = x[j_coord_offset+DIM*3+ZZ];
222
223             /* Calculate displacement vector */
224             dx00             = ix0 - jx0;
225             dy00             = iy0 - jy0;
226             dz00             = iz0 - jz0;
227             dx11             = ix1 - jx1;
228             dy11             = iy1 - jy1;
229             dz11             = iz1 - jz1;
230             dx12             = ix1 - jx2;
231             dy12             = iy1 - jy2;
232             dz12             = iz1 - jz2;
233             dx13             = ix1 - jx3;
234             dy13             = iy1 - jy3;
235             dz13             = iz1 - jz3;
236             dx21             = ix2 - jx1;
237             dy21             = iy2 - jy1;
238             dz21             = iz2 - jz1;
239             dx22             = ix2 - jx2;
240             dy22             = iy2 - jy2;
241             dz22             = iz2 - jz2;
242             dx23             = ix2 - jx3;
243             dy23             = iy2 - jy3;
244             dz23             = iz2 - jz3;
245             dx31             = ix3 - jx1;
246             dy31             = iy3 - jy1;
247             dz31             = iz3 - jz1;
248             dx32             = ix3 - jx2;
249             dy32             = iy3 - jy2;
250             dz32             = iz3 - jz2;
251             dx33             = ix3 - jx3;
252             dy33             = iy3 - jy3;
253             dz33             = iz3 - jz3;
254
255             /* Calculate squared distance and things based on it */
256             rsq00            = dx00*dx00+dy00*dy00+dz00*dz00;
257             rsq11            = dx11*dx11+dy11*dy11+dz11*dz11;
258             rsq12            = dx12*dx12+dy12*dy12+dz12*dz12;
259             rsq13            = dx13*dx13+dy13*dy13+dz13*dz13;
260             rsq21            = dx21*dx21+dy21*dy21+dz21*dz21;
261             rsq22            = dx22*dx22+dy22*dy22+dz22*dz22;
262             rsq23            = dx23*dx23+dy23*dy23+dz23*dz23;
263             rsq31            = dx31*dx31+dy31*dy31+dz31*dz31;
264             rsq32            = dx32*dx32+dy32*dy32+dz32*dz32;
265             rsq33            = dx33*dx33+dy33*dy33+dz33*dz33;
266
267             rinv00           = gmx_invsqrt(rsq00);
268             rinv11           = gmx_invsqrt(rsq11);
269             rinv12           = gmx_invsqrt(rsq12);
270             rinv13           = gmx_invsqrt(rsq13);
271             rinv21           = gmx_invsqrt(rsq21);
272             rinv22           = gmx_invsqrt(rsq22);
273             rinv23           = gmx_invsqrt(rsq23);
274             rinv31           = gmx_invsqrt(rsq31);
275             rinv32           = gmx_invsqrt(rsq32);
276             rinv33           = gmx_invsqrt(rsq33);
277
278             /**************************
279              * CALCULATE INTERACTIONS *
280              **************************/
281
282             r00              = rsq00*rinv00;
283
284             /* Calculate table index by multiplying r with table scale and truncate to integer */
285             rt               = r00*vftabscale;
286             vfitab           = rt;
287             vfeps            = rt-vfitab;
288             vfitab           = 3*4*vfitab;
289
290             /* CUBIC SPLINE TABLE DISPERSION */
291             vfitab          += 4;
292             Y                = vftab[vfitab];
293             F                = vftab[vfitab+1];
294             Geps             = vfeps*vftab[vfitab+2];
295             Heps2            = vfeps*vfeps*vftab[vfitab+3];
296             Fp               = F+Geps+Heps2;
297             VV               = Y+vfeps*Fp;
298             vvdw6            = c6_00*VV;
299             FF               = Fp+Geps+2.0*Heps2;
300             fvdw6            = c6_00*FF;
301
302             /* CUBIC SPLINE TABLE REPULSION */
303             Y                = vftab[vfitab+4];
304             F                = vftab[vfitab+5];
305             Geps             = vfeps*vftab[vfitab+6];
306             Heps2            = vfeps*vfeps*vftab[vfitab+7];
307             Fp               = F+Geps+Heps2;
308             VV               = Y+vfeps*Fp;
309             vvdw12           = c12_00*VV;
310             FF               = Fp+Geps+2.0*Heps2;
311             fvdw12           = c12_00*FF;
312             vvdw             = vvdw12+vvdw6;
313             fvdw             = -(fvdw6+fvdw12)*vftabscale*rinv00;
314
315             /* Update potential sums from outer loop */
316             vvdwsum         += vvdw;
317
318             fscal            = fvdw;
319
320             /* Calculate temporary vectorial force */
321             tx               = fscal*dx00;
322             ty               = fscal*dy00;
323             tz               = fscal*dz00;
324
325             /* Update vectorial force */
326             fix0            += tx;
327             fiy0            += ty;
328             fiz0            += tz;
329             f[j_coord_offset+DIM*0+XX] -= tx;
330             f[j_coord_offset+DIM*0+YY] -= ty;
331             f[j_coord_offset+DIM*0+ZZ] -= tz;
332
333             /**************************
334              * CALCULATE INTERACTIONS *
335              **************************/
336
337             r11              = rsq11*rinv11;
338
339             /* Calculate table index by multiplying r with table scale and truncate to integer */
340             rt               = r11*vftabscale;
341             vfitab           = rt;
342             vfeps            = rt-vfitab;
343             vfitab           = 3*4*vfitab;
344
345             /* CUBIC SPLINE TABLE ELECTROSTATICS */
346             Y                = vftab[vfitab];
347             F                = vftab[vfitab+1];
348             Geps             = vfeps*vftab[vfitab+2];
349             Heps2            = vfeps*vfeps*vftab[vfitab+3];
350             Fp               = F+Geps+Heps2;
351             VV               = Y+vfeps*Fp;
352             velec            = qq11*VV;
353             FF               = Fp+Geps+2.0*Heps2;
354             felec            = -qq11*FF*vftabscale*rinv11;
355
356             /* Update potential sums from outer loop */
357             velecsum        += velec;
358
359             fscal            = felec;
360
361             /* Calculate temporary vectorial force */
362             tx               = fscal*dx11;
363             ty               = fscal*dy11;
364             tz               = fscal*dz11;
365
366             /* Update vectorial force */
367             fix1            += tx;
368             fiy1            += ty;
369             fiz1            += tz;
370             f[j_coord_offset+DIM*1+XX] -= tx;
371             f[j_coord_offset+DIM*1+YY] -= ty;
372             f[j_coord_offset+DIM*1+ZZ] -= tz;
373
374             /**************************
375              * CALCULATE INTERACTIONS *
376              **************************/
377
378             r12              = rsq12*rinv12;
379
380             /* Calculate table index by multiplying r with table scale and truncate to integer */
381             rt               = r12*vftabscale;
382             vfitab           = rt;
383             vfeps            = rt-vfitab;
384             vfitab           = 3*4*vfitab;
385
386             /* CUBIC SPLINE TABLE ELECTROSTATICS */
387             Y                = vftab[vfitab];
388             F                = vftab[vfitab+1];
389             Geps             = vfeps*vftab[vfitab+2];
390             Heps2            = vfeps*vfeps*vftab[vfitab+3];
391             Fp               = F+Geps+Heps2;
392             VV               = Y+vfeps*Fp;
393             velec            = qq12*VV;
394             FF               = Fp+Geps+2.0*Heps2;
395             felec            = -qq12*FF*vftabscale*rinv12;
396
397             /* Update potential sums from outer loop */
398             velecsum        += velec;
399
400             fscal            = felec;
401
402             /* Calculate temporary vectorial force */
403             tx               = fscal*dx12;
404             ty               = fscal*dy12;
405             tz               = fscal*dz12;
406
407             /* Update vectorial force */
408             fix1            += tx;
409             fiy1            += ty;
410             fiz1            += tz;
411             f[j_coord_offset+DIM*2+XX] -= tx;
412             f[j_coord_offset+DIM*2+YY] -= ty;
413             f[j_coord_offset+DIM*2+ZZ] -= tz;
414
415             /**************************
416              * CALCULATE INTERACTIONS *
417              **************************/
418
419             r13              = rsq13*rinv13;
420
421             /* Calculate table index by multiplying r with table scale and truncate to integer */
422             rt               = r13*vftabscale;
423             vfitab           = rt;
424             vfeps            = rt-vfitab;
425             vfitab           = 3*4*vfitab;
426
427             /* CUBIC SPLINE TABLE ELECTROSTATICS */
428             Y                = vftab[vfitab];
429             F                = vftab[vfitab+1];
430             Geps             = vfeps*vftab[vfitab+2];
431             Heps2            = vfeps*vfeps*vftab[vfitab+3];
432             Fp               = F+Geps+Heps2;
433             VV               = Y+vfeps*Fp;
434             velec            = qq13*VV;
435             FF               = Fp+Geps+2.0*Heps2;
436             felec            = -qq13*FF*vftabscale*rinv13;
437
438             /* Update potential sums from outer loop */
439             velecsum        += velec;
440
441             fscal            = felec;
442
443             /* Calculate temporary vectorial force */
444             tx               = fscal*dx13;
445             ty               = fscal*dy13;
446             tz               = fscal*dz13;
447
448             /* Update vectorial force */
449             fix1            += tx;
450             fiy1            += ty;
451             fiz1            += tz;
452             f[j_coord_offset+DIM*3+XX] -= tx;
453             f[j_coord_offset+DIM*3+YY] -= ty;
454             f[j_coord_offset+DIM*3+ZZ] -= tz;
455
456             /**************************
457              * CALCULATE INTERACTIONS *
458              **************************/
459
460             r21              = rsq21*rinv21;
461
462             /* Calculate table index by multiplying r with table scale and truncate to integer */
463             rt               = r21*vftabscale;
464             vfitab           = rt;
465             vfeps            = rt-vfitab;
466             vfitab           = 3*4*vfitab;
467
468             /* CUBIC SPLINE TABLE ELECTROSTATICS */
469             Y                = vftab[vfitab];
470             F                = vftab[vfitab+1];
471             Geps             = vfeps*vftab[vfitab+2];
472             Heps2            = vfeps*vfeps*vftab[vfitab+3];
473             Fp               = F+Geps+Heps2;
474             VV               = Y+vfeps*Fp;
475             velec            = qq21*VV;
476             FF               = Fp+Geps+2.0*Heps2;
477             felec            = -qq21*FF*vftabscale*rinv21;
478
479             /* Update potential sums from outer loop */
480             velecsum        += velec;
481
482             fscal            = felec;
483
484             /* Calculate temporary vectorial force */
485             tx               = fscal*dx21;
486             ty               = fscal*dy21;
487             tz               = fscal*dz21;
488
489             /* Update vectorial force */
490             fix2            += tx;
491             fiy2            += ty;
492             fiz2            += tz;
493             f[j_coord_offset+DIM*1+XX] -= tx;
494             f[j_coord_offset+DIM*1+YY] -= ty;
495             f[j_coord_offset+DIM*1+ZZ] -= tz;
496
497             /**************************
498              * CALCULATE INTERACTIONS *
499              **************************/
500
501             r22              = rsq22*rinv22;
502
503             /* Calculate table index by multiplying r with table scale and truncate to integer */
504             rt               = r22*vftabscale;
505             vfitab           = rt;
506             vfeps            = rt-vfitab;
507             vfitab           = 3*4*vfitab;
508
509             /* CUBIC SPLINE TABLE ELECTROSTATICS */
510             Y                = vftab[vfitab];
511             F                = vftab[vfitab+1];
512             Geps             = vfeps*vftab[vfitab+2];
513             Heps2            = vfeps*vfeps*vftab[vfitab+3];
514             Fp               = F+Geps+Heps2;
515             VV               = Y+vfeps*Fp;
516             velec            = qq22*VV;
517             FF               = Fp+Geps+2.0*Heps2;
518             felec            = -qq22*FF*vftabscale*rinv22;
519
520             /* Update potential sums from outer loop */
521             velecsum        += velec;
522
523             fscal            = felec;
524
525             /* Calculate temporary vectorial force */
526             tx               = fscal*dx22;
527             ty               = fscal*dy22;
528             tz               = fscal*dz22;
529
530             /* Update vectorial force */
531             fix2            += tx;
532             fiy2            += ty;
533             fiz2            += tz;
534             f[j_coord_offset+DIM*2+XX] -= tx;
535             f[j_coord_offset+DIM*2+YY] -= ty;
536             f[j_coord_offset+DIM*2+ZZ] -= tz;
537
538             /**************************
539              * CALCULATE INTERACTIONS *
540              **************************/
541
542             r23              = rsq23*rinv23;
543
544             /* Calculate table index by multiplying r with table scale and truncate to integer */
545             rt               = r23*vftabscale;
546             vfitab           = rt;
547             vfeps            = rt-vfitab;
548             vfitab           = 3*4*vfitab;
549
550             /* CUBIC SPLINE TABLE ELECTROSTATICS */
551             Y                = vftab[vfitab];
552             F                = vftab[vfitab+1];
553             Geps             = vfeps*vftab[vfitab+2];
554             Heps2            = vfeps*vfeps*vftab[vfitab+3];
555             Fp               = F+Geps+Heps2;
556             VV               = Y+vfeps*Fp;
557             velec            = qq23*VV;
558             FF               = Fp+Geps+2.0*Heps2;
559             felec            = -qq23*FF*vftabscale*rinv23;
560
561             /* Update potential sums from outer loop */
562             velecsum        += velec;
563
564             fscal            = felec;
565
566             /* Calculate temporary vectorial force */
567             tx               = fscal*dx23;
568             ty               = fscal*dy23;
569             tz               = fscal*dz23;
570
571             /* Update vectorial force */
572             fix2            += tx;
573             fiy2            += ty;
574             fiz2            += tz;
575             f[j_coord_offset+DIM*3+XX] -= tx;
576             f[j_coord_offset+DIM*3+YY] -= ty;
577             f[j_coord_offset+DIM*3+ZZ] -= tz;
578
579             /**************************
580              * CALCULATE INTERACTIONS *
581              **************************/
582
583             r31              = rsq31*rinv31;
584
585             /* Calculate table index by multiplying r with table scale and truncate to integer */
586             rt               = r31*vftabscale;
587             vfitab           = rt;
588             vfeps            = rt-vfitab;
589             vfitab           = 3*4*vfitab;
590
591             /* CUBIC SPLINE TABLE ELECTROSTATICS */
592             Y                = vftab[vfitab];
593             F                = vftab[vfitab+1];
594             Geps             = vfeps*vftab[vfitab+2];
595             Heps2            = vfeps*vfeps*vftab[vfitab+3];
596             Fp               = F+Geps+Heps2;
597             VV               = Y+vfeps*Fp;
598             velec            = qq31*VV;
599             FF               = Fp+Geps+2.0*Heps2;
600             felec            = -qq31*FF*vftabscale*rinv31;
601
602             /* Update potential sums from outer loop */
603             velecsum        += velec;
604
605             fscal            = felec;
606
607             /* Calculate temporary vectorial force */
608             tx               = fscal*dx31;
609             ty               = fscal*dy31;
610             tz               = fscal*dz31;
611
612             /* Update vectorial force */
613             fix3            += tx;
614             fiy3            += ty;
615             fiz3            += tz;
616             f[j_coord_offset+DIM*1+XX] -= tx;
617             f[j_coord_offset+DIM*1+YY] -= ty;
618             f[j_coord_offset+DIM*1+ZZ] -= tz;
619
620             /**************************
621              * CALCULATE INTERACTIONS *
622              **************************/
623
624             r32              = rsq32*rinv32;
625
626             /* Calculate table index by multiplying r with table scale and truncate to integer */
627             rt               = r32*vftabscale;
628             vfitab           = rt;
629             vfeps            = rt-vfitab;
630             vfitab           = 3*4*vfitab;
631
632             /* CUBIC SPLINE TABLE ELECTROSTATICS */
633             Y                = vftab[vfitab];
634             F                = vftab[vfitab+1];
635             Geps             = vfeps*vftab[vfitab+2];
636             Heps2            = vfeps*vfeps*vftab[vfitab+3];
637             Fp               = F+Geps+Heps2;
638             VV               = Y+vfeps*Fp;
639             velec            = qq32*VV;
640             FF               = Fp+Geps+2.0*Heps2;
641             felec            = -qq32*FF*vftabscale*rinv32;
642
643             /* Update potential sums from outer loop */
644             velecsum        += velec;
645
646             fscal            = felec;
647
648             /* Calculate temporary vectorial force */
649             tx               = fscal*dx32;
650             ty               = fscal*dy32;
651             tz               = fscal*dz32;
652
653             /* Update vectorial force */
654             fix3            += tx;
655             fiy3            += ty;
656             fiz3            += tz;
657             f[j_coord_offset+DIM*2+XX] -= tx;
658             f[j_coord_offset+DIM*2+YY] -= ty;
659             f[j_coord_offset+DIM*2+ZZ] -= tz;
660
661             /**************************
662              * CALCULATE INTERACTIONS *
663              **************************/
664
665             r33              = rsq33*rinv33;
666
667             /* Calculate table index by multiplying r with table scale and truncate to integer */
668             rt               = r33*vftabscale;
669             vfitab           = rt;
670             vfeps            = rt-vfitab;
671             vfitab           = 3*4*vfitab;
672
673             /* CUBIC SPLINE TABLE ELECTROSTATICS */
674             Y                = vftab[vfitab];
675             F                = vftab[vfitab+1];
676             Geps             = vfeps*vftab[vfitab+2];
677             Heps2            = vfeps*vfeps*vftab[vfitab+3];
678             Fp               = F+Geps+Heps2;
679             VV               = Y+vfeps*Fp;
680             velec            = qq33*VV;
681             FF               = Fp+Geps+2.0*Heps2;
682             felec            = -qq33*FF*vftabscale*rinv33;
683
684             /* Update potential sums from outer loop */
685             velecsum        += velec;
686
687             fscal            = felec;
688
689             /* Calculate temporary vectorial force */
690             tx               = fscal*dx33;
691             ty               = fscal*dy33;
692             tz               = fscal*dz33;
693
694             /* Update vectorial force */
695             fix3            += tx;
696             fiy3            += ty;
697             fiz3            += tz;
698             f[j_coord_offset+DIM*3+XX] -= tx;
699             f[j_coord_offset+DIM*3+YY] -= ty;
700             f[j_coord_offset+DIM*3+ZZ] -= tz;
701
702             /* Inner loop uses 424 flops */
703         }
704         /* End of innermost loop */
705
706         tx = ty = tz = 0;
707         f[i_coord_offset+DIM*0+XX] += fix0;
708         f[i_coord_offset+DIM*0+YY] += fiy0;
709         f[i_coord_offset+DIM*0+ZZ] += fiz0;
710         tx                         += fix0;
711         ty                         += fiy0;
712         tz                         += fiz0;
713         f[i_coord_offset+DIM*1+XX] += fix1;
714         f[i_coord_offset+DIM*1+YY] += fiy1;
715         f[i_coord_offset+DIM*1+ZZ] += fiz1;
716         tx                         += fix1;
717         ty                         += fiy1;
718         tz                         += fiz1;
719         f[i_coord_offset+DIM*2+XX] += fix2;
720         f[i_coord_offset+DIM*2+YY] += fiy2;
721         f[i_coord_offset+DIM*2+ZZ] += fiz2;
722         tx                         += fix2;
723         ty                         += fiy2;
724         tz                         += fiz2;
725         f[i_coord_offset+DIM*3+XX] += fix3;
726         f[i_coord_offset+DIM*3+YY] += fiy3;
727         f[i_coord_offset+DIM*3+ZZ] += fiz3;
728         tx                         += fix3;
729         ty                         += fiy3;
730         tz                         += fiz3;
731         fshift[i_shift_offset+XX]  += tx;
732         fshift[i_shift_offset+YY]  += ty;
733         fshift[i_shift_offset+ZZ]  += tz;
734
735         ggid                        = gid[iidx];
736         /* Update potential energies */
737         kernel_data->energygrp_elec[ggid] += velecsum;
738         kernel_data->energygrp_vdw[ggid] += vvdwsum;
739
740         /* Increment number of inner iterations */
741         inneriter                  += j_index_end - j_index_start;
742
743         /* Outer loop uses 41 flops */
744     }
745
746     /* Increment number of outer iterations */
747     outeriter        += nri;
748
749     /* Update outer/inner flops */
750
751     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W4W4_VF,outeriter*41 + inneriter*424);
752 }
753 /*
754  * Gromacs nonbonded kernel:   nb_kernel_ElecCSTab_VdwCSTab_GeomW4W4_F_c
755  * Electrostatics interaction: CubicSplineTable
756  * VdW interaction:            CubicSplineTable
757  * Geometry:                   Water4-Water4
758  * Calculate force/pot:        Force
759  */
760 void
761 nb_kernel_ElecCSTab_VdwCSTab_GeomW4W4_F_c
762                     (t_nblist                    * gmx_restrict       nlist,
763                      rvec                        * gmx_restrict          xx,
764                      rvec                        * gmx_restrict          ff,
765                      t_forcerec                  * gmx_restrict          fr,
766                      t_mdatoms                   * gmx_restrict     mdatoms,
767                      nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
768                      t_nrnb                      * gmx_restrict        nrnb)
769 {
770     int              i_shift_offset,i_coord_offset,j_coord_offset;
771     int              j_index_start,j_index_end;
772     int              nri,inr,ggid,iidx,jidx,jnr,outeriter,inneriter;
773     real             shX,shY,shZ,tx,ty,tz,fscal,rcutoff,rcutoff2;
774     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
775     real             *shiftvec,*fshift,*x,*f;
776     int              vdwioffset0;
777     real             ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
778     int              vdwioffset1;
779     real             ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
780     int              vdwioffset2;
781     real             ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
782     int              vdwioffset3;
783     real             ix3,iy3,iz3,fix3,fiy3,fiz3,iq3,isai3;
784     int              vdwjidx0;
785     real             jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
786     int              vdwjidx1;
787     real             jx1,jy1,jz1,fjx1,fjy1,fjz1,jq1,isaj1;
788     int              vdwjidx2;
789     real             jx2,jy2,jz2,fjx2,fjy2,fjz2,jq2,isaj2;
790     int              vdwjidx3;
791     real             jx3,jy3,jz3,fjx3,fjy3,fjz3,jq3,isaj3;
792     real             dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00,cexp1_00,cexp2_00;
793     real             dx11,dy11,dz11,rsq11,rinv11,rinvsq11,r11,qq11,c6_11,c12_11,cexp1_11,cexp2_11;
794     real             dx12,dy12,dz12,rsq12,rinv12,rinvsq12,r12,qq12,c6_12,c12_12,cexp1_12,cexp2_12;
795     real             dx13,dy13,dz13,rsq13,rinv13,rinvsq13,r13,qq13,c6_13,c12_13,cexp1_13,cexp2_13;
796     real             dx21,dy21,dz21,rsq21,rinv21,rinvsq21,r21,qq21,c6_21,c12_21,cexp1_21,cexp2_21;
797     real             dx22,dy22,dz22,rsq22,rinv22,rinvsq22,r22,qq22,c6_22,c12_22,cexp1_22,cexp2_22;
798     real             dx23,dy23,dz23,rsq23,rinv23,rinvsq23,r23,qq23,c6_23,c12_23,cexp1_23,cexp2_23;
799     real             dx31,dy31,dz31,rsq31,rinv31,rinvsq31,r31,qq31,c6_31,c12_31,cexp1_31,cexp2_31;
800     real             dx32,dy32,dz32,rsq32,rinv32,rinvsq32,r32,qq32,c6_32,c12_32,cexp1_32,cexp2_32;
801     real             dx33,dy33,dz33,rsq33,rinv33,rinvsq33,r33,qq33,c6_33,c12_33,cexp1_33,cexp2_33;
802     real             velec,felec,velecsum,facel,crf,krf,krf2;
803     real             *charge;
804     int              nvdwtype;
805     real             rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,br,vvdwexp,sh_vdw_invrcut6;
806     int              *vdwtype;
807     real             *vdwparam;
808     int              vfitab;
809     real             rt,vfeps,vftabscale,Y,F,Geps,Heps2,Fp,VV,FF;
810     real             *vftab;
811
812     x                = xx[0];
813     f                = ff[0];
814
815     nri              = nlist->nri;
816     iinr             = nlist->iinr;
817     jindex           = nlist->jindex;
818     jjnr             = nlist->jjnr;
819     shiftidx         = nlist->shift;
820     gid              = nlist->gid;
821     shiftvec         = fr->shift_vec[0];
822     fshift           = fr->fshift[0];
823     facel            = fr->epsfac;
824     charge           = mdatoms->chargeA;
825     nvdwtype         = fr->ntype;
826     vdwparam         = fr->nbfp;
827     vdwtype          = mdatoms->typeA;
828
829     vftab            = kernel_data->table_elec_vdw->data;
830     vftabscale       = kernel_data->table_elec_vdw->scale;
831
832     /* Setup water-specific parameters */
833     inr              = nlist->iinr[0];
834     iq1              = facel*charge[inr+1];
835     iq2              = facel*charge[inr+2];
836     iq3              = facel*charge[inr+3];
837     vdwioffset0      = 2*nvdwtype*vdwtype[inr+0];
838
839     jq1              = charge[inr+1];
840     jq2              = charge[inr+2];
841     jq3              = charge[inr+3];
842     vdwjidx0         = 2*vdwtype[inr+0];
843     c6_00            = vdwparam[vdwioffset0+vdwjidx0];
844     c12_00           = vdwparam[vdwioffset0+vdwjidx0+1];
845     qq11             = iq1*jq1;
846     qq12             = iq1*jq2;
847     qq13             = iq1*jq3;
848     qq21             = iq2*jq1;
849     qq22             = iq2*jq2;
850     qq23             = iq2*jq3;
851     qq31             = iq3*jq1;
852     qq32             = iq3*jq2;
853     qq33             = iq3*jq3;
854
855     outeriter        = 0;
856     inneriter        = 0;
857
858     /* Start outer loop over neighborlists */
859     for(iidx=0; iidx<nri; iidx++)
860     {
861         /* Load shift vector for this list */
862         i_shift_offset   = DIM*shiftidx[iidx];
863         shX              = shiftvec[i_shift_offset+XX];
864         shY              = shiftvec[i_shift_offset+YY];
865         shZ              = shiftvec[i_shift_offset+ZZ];
866
867         /* Load limits for loop over neighbors */
868         j_index_start    = jindex[iidx];
869         j_index_end      = jindex[iidx+1];
870
871         /* Get outer coordinate index */
872         inr              = iinr[iidx];
873         i_coord_offset   = DIM*inr;
874
875         /* Load i particle coords and add shift vector */
876         ix0              = shX + x[i_coord_offset+DIM*0+XX];
877         iy0              = shY + x[i_coord_offset+DIM*0+YY];
878         iz0              = shZ + x[i_coord_offset+DIM*0+ZZ];
879         ix1              = shX + x[i_coord_offset+DIM*1+XX];
880         iy1              = shY + x[i_coord_offset+DIM*1+YY];
881         iz1              = shZ + x[i_coord_offset+DIM*1+ZZ];
882         ix2              = shX + x[i_coord_offset+DIM*2+XX];
883         iy2              = shY + x[i_coord_offset+DIM*2+YY];
884         iz2              = shZ + x[i_coord_offset+DIM*2+ZZ];
885         ix3              = shX + x[i_coord_offset+DIM*3+XX];
886         iy3              = shY + x[i_coord_offset+DIM*3+YY];
887         iz3              = shZ + x[i_coord_offset+DIM*3+ZZ];
888
889         fix0             = 0.0;
890         fiy0             = 0.0;
891         fiz0             = 0.0;
892         fix1             = 0.0;
893         fiy1             = 0.0;
894         fiz1             = 0.0;
895         fix2             = 0.0;
896         fiy2             = 0.0;
897         fiz2             = 0.0;
898         fix3             = 0.0;
899         fiy3             = 0.0;
900         fiz3             = 0.0;
901
902         /* Start inner kernel loop */
903         for(jidx=j_index_start; jidx<j_index_end; jidx++)
904         {
905             /* Get j neighbor index, and coordinate index */
906             jnr              = jjnr[jidx];
907             j_coord_offset   = DIM*jnr;
908
909             /* load j atom coordinates */
910             jx0              = x[j_coord_offset+DIM*0+XX];
911             jy0              = x[j_coord_offset+DIM*0+YY];
912             jz0              = x[j_coord_offset+DIM*0+ZZ];
913             jx1              = x[j_coord_offset+DIM*1+XX];
914             jy1              = x[j_coord_offset+DIM*1+YY];
915             jz1              = x[j_coord_offset+DIM*1+ZZ];
916             jx2              = x[j_coord_offset+DIM*2+XX];
917             jy2              = x[j_coord_offset+DIM*2+YY];
918             jz2              = x[j_coord_offset+DIM*2+ZZ];
919             jx3              = x[j_coord_offset+DIM*3+XX];
920             jy3              = x[j_coord_offset+DIM*3+YY];
921             jz3              = x[j_coord_offset+DIM*3+ZZ];
922
923             /* Calculate displacement vector */
924             dx00             = ix0 - jx0;
925             dy00             = iy0 - jy0;
926             dz00             = iz0 - jz0;
927             dx11             = ix1 - jx1;
928             dy11             = iy1 - jy1;
929             dz11             = iz1 - jz1;
930             dx12             = ix1 - jx2;
931             dy12             = iy1 - jy2;
932             dz12             = iz1 - jz2;
933             dx13             = ix1 - jx3;
934             dy13             = iy1 - jy3;
935             dz13             = iz1 - jz3;
936             dx21             = ix2 - jx1;
937             dy21             = iy2 - jy1;
938             dz21             = iz2 - jz1;
939             dx22             = ix2 - jx2;
940             dy22             = iy2 - jy2;
941             dz22             = iz2 - jz2;
942             dx23             = ix2 - jx3;
943             dy23             = iy2 - jy3;
944             dz23             = iz2 - jz3;
945             dx31             = ix3 - jx1;
946             dy31             = iy3 - jy1;
947             dz31             = iz3 - jz1;
948             dx32             = ix3 - jx2;
949             dy32             = iy3 - jy2;
950             dz32             = iz3 - jz2;
951             dx33             = ix3 - jx3;
952             dy33             = iy3 - jy3;
953             dz33             = iz3 - jz3;
954
955             /* Calculate squared distance and things based on it */
956             rsq00            = dx00*dx00+dy00*dy00+dz00*dz00;
957             rsq11            = dx11*dx11+dy11*dy11+dz11*dz11;
958             rsq12            = dx12*dx12+dy12*dy12+dz12*dz12;
959             rsq13            = dx13*dx13+dy13*dy13+dz13*dz13;
960             rsq21            = dx21*dx21+dy21*dy21+dz21*dz21;
961             rsq22            = dx22*dx22+dy22*dy22+dz22*dz22;
962             rsq23            = dx23*dx23+dy23*dy23+dz23*dz23;
963             rsq31            = dx31*dx31+dy31*dy31+dz31*dz31;
964             rsq32            = dx32*dx32+dy32*dy32+dz32*dz32;
965             rsq33            = dx33*dx33+dy33*dy33+dz33*dz33;
966
967             rinv00           = gmx_invsqrt(rsq00);
968             rinv11           = gmx_invsqrt(rsq11);
969             rinv12           = gmx_invsqrt(rsq12);
970             rinv13           = gmx_invsqrt(rsq13);
971             rinv21           = gmx_invsqrt(rsq21);
972             rinv22           = gmx_invsqrt(rsq22);
973             rinv23           = gmx_invsqrt(rsq23);
974             rinv31           = gmx_invsqrt(rsq31);
975             rinv32           = gmx_invsqrt(rsq32);
976             rinv33           = gmx_invsqrt(rsq33);
977
978             /**************************
979              * CALCULATE INTERACTIONS *
980              **************************/
981
982             r00              = rsq00*rinv00;
983
984             /* Calculate table index by multiplying r with table scale and truncate to integer */
985             rt               = r00*vftabscale;
986             vfitab           = rt;
987             vfeps            = rt-vfitab;
988             vfitab           = 3*4*vfitab;
989
990             /* CUBIC SPLINE TABLE DISPERSION */
991             vfitab          += 4;
992             F                = vftab[vfitab+1];
993             Geps             = vfeps*vftab[vfitab+2];
994             Heps2            = vfeps*vfeps*vftab[vfitab+3];
995             Fp               = F+Geps+Heps2;
996             FF               = Fp+Geps+2.0*Heps2;
997             fvdw6            = c6_00*FF;
998
999             /* CUBIC SPLINE TABLE REPULSION */
1000             F                = vftab[vfitab+5];
1001             Geps             = vfeps*vftab[vfitab+6];
1002             Heps2            = vfeps*vfeps*vftab[vfitab+7];
1003             Fp               = F+Geps+Heps2;
1004             FF               = Fp+Geps+2.0*Heps2;
1005             fvdw12           = c12_00*FF;
1006             fvdw             = -(fvdw6+fvdw12)*vftabscale*rinv00;
1007
1008             fscal            = fvdw;
1009
1010             /* Calculate temporary vectorial force */
1011             tx               = fscal*dx00;
1012             ty               = fscal*dy00;
1013             tz               = fscal*dz00;
1014
1015             /* Update vectorial force */
1016             fix0            += tx;
1017             fiy0            += ty;
1018             fiz0            += tz;
1019             f[j_coord_offset+DIM*0+XX] -= tx;
1020             f[j_coord_offset+DIM*0+YY] -= ty;
1021             f[j_coord_offset+DIM*0+ZZ] -= tz;
1022
1023             /**************************
1024              * CALCULATE INTERACTIONS *
1025              **************************/
1026
1027             r11              = rsq11*rinv11;
1028
1029             /* Calculate table index by multiplying r with table scale and truncate to integer */
1030             rt               = r11*vftabscale;
1031             vfitab           = rt;
1032             vfeps            = rt-vfitab;
1033             vfitab           = 3*4*vfitab;
1034
1035             /* CUBIC SPLINE TABLE ELECTROSTATICS */
1036             F                = vftab[vfitab+1];
1037             Geps             = vfeps*vftab[vfitab+2];
1038             Heps2            = vfeps*vfeps*vftab[vfitab+3];
1039             Fp               = F+Geps+Heps2;
1040             FF               = Fp+Geps+2.0*Heps2;
1041             felec            = -qq11*FF*vftabscale*rinv11;
1042
1043             fscal            = felec;
1044
1045             /* Calculate temporary vectorial force */
1046             tx               = fscal*dx11;
1047             ty               = fscal*dy11;
1048             tz               = fscal*dz11;
1049
1050             /* Update vectorial force */
1051             fix1            += tx;
1052             fiy1            += ty;
1053             fiz1            += tz;
1054             f[j_coord_offset+DIM*1+XX] -= tx;
1055             f[j_coord_offset+DIM*1+YY] -= ty;
1056             f[j_coord_offset+DIM*1+ZZ] -= tz;
1057
1058             /**************************
1059              * CALCULATE INTERACTIONS *
1060              **************************/
1061
1062             r12              = rsq12*rinv12;
1063
1064             /* Calculate table index by multiplying r with table scale and truncate to integer */
1065             rt               = r12*vftabscale;
1066             vfitab           = rt;
1067             vfeps            = rt-vfitab;
1068             vfitab           = 3*4*vfitab;
1069
1070             /* CUBIC SPLINE TABLE ELECTROSTATICS */
1071             F                = vftab[vfitab+1];
1072             Geps             = vfeps*vftab[vfitab+2];
1073             Heps2            = vfeps*vfeps*vftab[vfitab+3];
1074             Fp               = F+Geps+Heps2;
1075             FF               = Fp+Geps+2.0*Heps2;
1076             felec            = -qq12*FF*vftabscale*rinv12;
1077
1078             fscal            = felec;
1079
1080             /* Calculate temporary vectorial force */
1081             tx               = fscal*dx12;
1082             ty               = fscal*dy12;
1083             tz               = fscal*dz12;
1084
1085             /* Update vectorial force */
1086             fix1            += tx;
1087             fiy1            += ty;
1088             fiz1            += tz;
1089             f[j_coord_offset+DIM*2+XX] -= tx;
1090             f[j_coord_offset+DIM*2+YY] -= ty;
1091             f[j_coord_offset+DIM*2+ZZ] -= tz;
1092
1093             /**************************
1094              * CALCULATE INTERACTIONS *
1095              **************************/
1096
1097             r13              = rsq13*rinv13;
1098
1099             /* Calculate table index by multiplying r with table scale and truncate to integer */
1100             rt               = r13*vftabscale;
1101             vfitab           = rt;
1102             vfeps            = rt-vfitab;
1103             vfitab           = 3*4*vfitab;
1104
1105             /* CUBIC SPLINE TABLE ELECTROSTATICS */
1106             F                = vftab[vfitab+1];
1107             Geps             = vfeps*vftab[vfitab+2];
1108             Heps2            = vfeps*vfeps*vftab[vfitab+3];
1109             Fp               = F+Geps+Heps2;
1110             FF               = Fp+Geps+2.0*Heps2;
1111             felec            = -qq13*FF*vftabscale*rinv13;
1112
1113             fscal            = felec;
1114
1115             /* Calculate temporary vectorial force */
1116             tx               = fscal*dx13;
1117             ty               = fscal*dy13;
1118             tz               = fscal*dz13;
1119
1120             /* Update vectorial force */
1121             fix1            += tx;
1122             fiy1            += ty;
1123             fiz1            += tz;
1124             f[j_coord_offset+DIM*3+XX] -= tx;
1125             f[j_coord_offset+DIM*3+YY] -= ty;
1126             f[j_coord_offset+DIM*3+ZZ] -= tz;
1127
1128             /**************************
1129              * CALCULATE INTERACTIONS *
1130              **************************/
1131
1132             r21              = rsq21*rinv21;
1133
1134             /* Calculate table index by multiplying r with table scale and truncate to integer */
1135             rt               = r21*vftabscale;
1136             vfitab           = rt;
1137             vfeps            = rt-vfitab;
1138             vfitab           = 3*4*vfitab;
1139
1140             /* CUBIC SPLINE TABLE ELECTROSTATICS */
1141             F                = vftab[vfitab+1];
1142             Geps             = vfeps*vftab[vfitab+2];
1143             Heps2            = vfeps*vfeps*vftab[vfitab+3];
1144             Fp               = F+Geps+Heps2;
1145             FF               = Fp+Geps+2.0*Heps2;
1146             felec            = -qq21*FF*vftabscale*rinv21;
1147
1148             fscal            = felec;
1149
1150             /* Calculate temporary vectorial force */
1151             tx               = fscal*dx21;
1152             ty               = fscal*dy21;
1153             tz               = fscal*dz21;
1154
1155             /* Update vectorial force */
1156             fix2            += tx;
1157             fiy2            += ty;
1158             fiz2            += tz;
1159             f[j_coord_offset+DIM*1+XX] -= tx;
1160             f[j_coord_offset+DIM*1+YY] -= ty;
1161             f[j_coord_offset+DIM*1+ZZ] -= tz;
1162
1163             /**************************
1164              * CALCULATE INTERACTIONS *
1165              **************************/
1166
1167             r22              = rsq22*rinv22;
1168
1169             /* Calculate table index by multiplying r with table scale and truncate to integer */
1170             rt               = r22*vftabscale;
1171             vfitab           = rt;
1172             vfeps            = rt-vfitab;
1173             vfitab           = 3*4*vfitab;
1174
1175             /* CUBIC SPLINE TABLE ELECTROSTATICS */
1176             F                = vftab[vfitab+1];
1177             Geps             = vfeps*vftab[vfitab+2];
1178             Heps2            = vfeps*vfeps*vftab[vfitab+3];
1179             Fp               = F+Geps+Heps2;
1180             FF               = Fp+Geps+2.0*Heps2;
1181             felec            = -qq22*FF*vftabscale*rinv22;
1182
1183             fscal            = felec;
1184
1185             /* Calculate temporary vectorial force */
1186             tx               = fscal*dx22;
1187             ty               = fscal*dy22;
1188             tz               = fscal*dz22;
1189
1190             /* Update vectorial force */
1191             fix2            += tx;
1192             fiy2            += ty;
1193             fiz2            += tz;
1194             f[j_coord_offset+DIM*2+XX] -= tx;
1195             f[j_coord_offset+DIM*2+YY] -= ty;
1196             f[j_coord_offset+DIM*2+ZZ] -= tz;
1197
1198             /**************************
1199              * CALCULATE INTERACTIONS *
1200              **************************/
1201
1202             r23              = rsq23*rinv23;
1203
1204             /* Calculate table index by multiplying r with table scale and truncate to integer */
1205             rt               = r23*vftabscale;
1206             vfitab           = rt;
1207             vfeps            = rt-vfitab;
1208             vfitab           = 3*4*vfitab;
1209
1210             /* CUBIC SPLINE TABLE ELECTROSTATICS */
1211             F                = vftab[vfitab+1];
1212             Geps             = vfeps*vftab[vfitab+2];
1213             Heps2            = vfeps*vfeps*vftab[vfitab+3];
1214             Fp               = F+Geps+Heps2;
1215             FF               = Fp+Geps+2.0*Heps2;
1216             felec            = -qq23*FF*vftabscale*rinv23;
1217
1218             fscal            = felec;
1219
1220             /* Calculate temporary vectorial force */
1221             tx               = fscal*dx23;
1222             ty               = fscal*dy23;
1223             tz               = fscal*dz23;
1224
1225             /* Update vectorial force */
1226             fix2            += tx;
1227             fiy2            += ty;
1228             fiz2            += tz;
1229             f[j_coord_offset+DIM*3+XX] -= tx;
1230             f[j_coord_offset+DIM*3+YY] -= ty;
1231             f[j_coord_offset+DIM*3+ZZ] -= tz;
1232
1233             /**************************
1234              * CALCULATE INTERACTIONS *
1235              **************************/
1236
1237             r31              = rsq31*rinv31;
1238
1239             /* Calculate table index by multiplying r with table scale and truncate to integer */
1240             rt               = r31*vftabscale;
1241             vfitab           = rt;
1242             vfeps            = rt-vfitab;
1243             vfitab           = 3*4*vfitab;
1244
1245             /* CUBIC SPLINE TABLE ELECTROSTATICS */
1246             F                = vftab[vfitab+1];
1247             Geps             = vfeps*vftab[vfitab+2];
1248             Heps2            = vfeps*vfeps*vftab[vfitab+3];
1249             Fp               = F+Geps+Heps2;
1250             FF               = Fp+Geps+2.0*Heps2;
1251             felec            = -qq31*FF*vftabscale*rinv31;
1252
1253             fscal            = felec;
1254
1255             /* Calculate temporary vectorial force */
1256             tx               = fscal*dx31;
1257             ty               = fscal*dy31;
1258             tz               = fscal*dz31;
1259
1260             /* Update vectorial force */
1261             fix3            += tx;
1262             fiy3            += ty;
1263             fiz3            += tz;
1264             f[j_coord_offset+DIM*1+XX] -= tx;
1265             f[j_coord_offset+DIM*1+YY] -= ty;
1266             f[j_coord_offset+DIM*1+ZZ] -= tz;
1267
1268             /**************************
1269              * CALCULATE INTERACTIONS *
1270              **************************/
1271
1272             r32              = rsq32*rinv32;
1273
1274             /* Calculate table index by multiplying r with table scale and truncate to integer */
1275             rt               = r32*vftabscale;
1276             vfitab           = rt;
1277             vfeps            = rt-vfitab;
1278             vfitab           = 3*4*vfitab;
1279
1280             /* CUBIC SPLINE TABLE ELECTROSTATICS */
1281             F                = vftab[vfitab+1];
1282             Geps             = vfeps*vftab[vfitab+2];
1283             Heps2            = vfeps*vfeps*vftab[vfitab+3];
1284             Fp               = F+Geps+Heps2;
1285             FF               = Fp+Geps+2.0*Heps2;
1286             felec            = -qq32*FF*vftabscale*rinv32;
1287
1288             fscal            = felec;
1289
1290             /* Calculate temporary vectorial force */
1291             tx               = fscal*dx32;
1292             ty               = fscal*dy32;
1293             tz               = fscal*dz32;
1294
1295             /* Update vectorial force */
1296             fix3            += tx;
1297             fiy3            += ty;
1298             fiz3            += tz;
1299             f[j_coord_offset+DIM*2+XX] -= tx;
1300             f[j_coord_offset+DIM*2+YY] -= ty;
1301             f[j_coord_offset+DIM*2+ZZ] -= tz;
1302
1303             /**************************
1304              * CALCULATE INTERACTIONS *
1305              **************************/
1306
1307             r33              = rsq33*rinv33;
1308
1309             /* Calculate table index by multiplying r with table scale and truncate to integer */
1310             rt               = r33*vftabscale;
1311             vfitab           = rt;
1312             vfeps            = rt-vfitab;
1313             vfitab           = 3*4*vfitab;
1314
1315             /* CUBIC SPLINE TABLE ELECTROSTATICS */
1316             F                = vftab[vfitab+1];
1317             Geps             = vfeps*vftab[vfitab+2];
1318             Heps2            = vfeps*vfeps*vftab[vfitab+3];
1319             Fp               = F+Geps+Heps2;
1320             FF               = Fp+Geps+2.0*Heps2;
1321             felec            = -qq33*FF*vftabscale*rinv33;
1322
1323             fscal            = felec;
1324
1325             /* Calculate temporary vectorial force */
1326             tx               = fscal*dx33;
1327             ty               = fscal*dy33;
1328             tz               = fscal*dz33;
1329
1330             /* Update vectorial force */
1331             fix3            += tx;
1332             fiy3            += ty;
1333             fiz3            += tz;
1334             f[j_coord_offset+DIM*3+XX] -= tx;
1335             f[j_coord_offset+DIM*3+YY] -= ty;
1336             f[j_coord_offset+DIM*3+ZZ] -= tz;
1337
1338             /* Inner loop uses 380 flops */
1339         }
1340         /* End of innermost loop */
1341
1342         tx = ty = tz = 0;
1343         f[i_coord_offset+DIM*0+XX] += fix0;
1344         f[i_coord_offset+DIM*0+YY] += fiy0;
1345         f[i_coord_offset+DIM*0+ZZ] += fiz0;
1346         tx                         += fix0;
1347         ty                         += fiy0;
1348         tz                         += fiz0;
1349         f[i_coord_offset+DIM*1+XX] += fix1;
1350         f[i_coord_offset+DIM*1+YY] += fiy1;
1351         f[i_coord_offset+DIM*1+ZZ] += fiz1;
1352         tx                         += fix1;
1353         ty                         += fiy1;
1354         tz                         += fiz1;
1355         f[i_coord_offset+DIM*2+XX] += fix2;
1356         f[i_coord_offset+DIM*2+YY] += fiy2;
1357         f[i_coord_offset+DIM*2+ZZ] += fiz2;
1358         tx                         += fix2;
1359         ty                         += fiy2;
1360         tz                         += fiz2;
1361         f[i_coord_offset+DIM*3+XX] += fix3;
1362         f[i_coord_offset+DIM*3+YY] += fiy3;
1363         f[i_coord_offset+DIM*3+ZZ] += fiz3;
1364         tx                         += fix3;
1365         ty                         += fiy3;
1366         tz                         += fiz3;
1367         fshift[i_shift_offset+XX]  += tx;
1368         fshift[i_shift_offset+YY]  += ty;
1369         fshift[i_shift_offset+ZZ]  += tz;
1370
1371         /* Increment number of inner iterations */
1372         inneriter                  += j_index_end - j_index_start;
1373
1374         /* Outer loop uses 39 flops */
1375     }
1376
1377     /* Increment number of outer iterations */
1378     outeriter        += nri;
1379
1380     /* Update outer/inner flops */
1381
1382     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W4W4_F,outeriter*39 + inneriter*380);
1383 }