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