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