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