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