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