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