Remove all unnecessary HAVE_CONFIG_H
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_c / nb_kernel_ElecEwSh_VdwBhamSh_GeomW4P1_c.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2013,2014, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35 /*
36  * Note: this file was generated by the GROMACS c kernel generator.
37  */
38 #include "config.h"
39
40 #include <math.h>
41
42 #include "../nb_kernel.h"
43 #include "types/simple.h"
44 #include "gromacs/math/vec.h"
45 #include "nrnb.h"
46
47 /*
48  * Gromacs nonbonded kernel:   nb_kernel_ElecEwSh_VdwBhamSh_GeomW4P1_VF_c
49  * Electrostatics interaction: Ewald
50  * VdW interaction:            Buckingham
51  * Geometry:                   Water4-Particle
52  * Calculate force/pot:        PotentialAndForce
53  */
54 void
55 nb_kernel_ElecEwSh_VdwBhamSh_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              ewitab;
91     real             ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace;
92     real             *ewtab;
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     sh_ewald         = fr->ic->sh_ewald;
112     ewtab            = fr->ic->tabq_coul_FDV0;
113     ewtabscale       = fr->ic->tabq_scale;
114     ewtabhalfspace   = 0.5/ewtabscale;
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      = 3*nvdwtype*vdwtype[inr+0];
122
123     /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
124     rcutoff          = fr->rcoulomb;
125     rcutoff2         = rcutoff*rcutoff;
126
127     sh_vdw_invrcut6  = fr->ic->sh_invrc6;
128     rvdw             = fr->rvdw;
129
130     outeriter        = 0;
131     inneriter        = 0;
132
133     /* Start outer loop over neighborlists */
134     for(iidx=0; iidx<nri; iidx++)
135     {
136         /* Load shift vector for this list */
137         i_shift_offset   = DIM*shiftidx[iidx];
138         shX              = shiftvec[i_shift_offset+XX];
139         shY              = shiftvec[i_shift_offset+YY];
140         shZ              = shiftvec[i_shift_offset+ZZ];
141
142         /* Load limits for loop over neighbors */
143         j_index_start    = jindex[iidx];
144         j_index_end      = jindex[iidx+1];
145
146         /* Get outer coordinate index */
147         inr              = iinr[iidx];
148         i_coord_offset   = DIM*inr;
149
150         /* Load i particle coords and add shift vector */
151         ix0              = shX + x[i_coord_offset+DIM*0+XX];
152         iy0              = shY + x[i_coord_offset+DIM*0+YY];
153         iz0              = shZ + x[i_coord_offset+DIM*0+ZZ];
154         ix1              = shX + x[i_coord_offset+DIM*1+XX];
155         iy1              = shY + x[i_coord_offset+DIM*1+YY];
156         iz1              = shZ + x[i_coord_offset+DIM*1+ZZ];
157         ix2              = shX + x[i_coord_offset+DIM*2+XX];
158         iy2              = shY + x[i_coord_offset+DIM*2+YY];
159         iz2              = shZ + x[i_coord_offset+DIM*2+ZZ];
160         ix3              = shX + x[i_coord_offset+DIM*3+XX];
161         iy3              = shY + x[i_coord_offset+DIM*3+YY];
162         iz3              = shZ + x[i_coord_offset+DIM*3+ZZ];
163
164         fix0             = 0.0;
165         fiy0             = 0.0;
166         fiz0             = 0.0;
167         fix1             = 0.0;
168         fiy1             = 0.0;
169         fiz1             = 0.0;
170         fix2             = 0.0;
171         fiy2             = 0.0;
172         fiz2             = 0.0;
173         fix3             = 0.0;
174         fiy3             = 0.0;
175         fiz3             = 0.0;
176
177         /* Reset potential sums */
178         velecsum         = 0.0;
179         vvdwsum          = 0.0;
180
181         /* Start inner kernel loop */
182         for(jidx=j_index_start; jidx<j_index_end; jidx++)
183         {
184             /* Get j neighbor index, and coordinate index */
185             jnr              = jjnr[jidx];
186             j_coord_offset   = DIM*jnr;
187
188             /* load j atom coordinates */
189             jx0              = x[j_coord_offset+DIM*0+XX];
190             jy0              = x[j_coord_offset+DIM*0+YY];
191             jz0              = x[j_coord_offset+DIM*0+ZZ];
192
193             /* Calculate displacement vector */
194             dx00             = ix0 - jx0;
195             dy00             = iy0 - jy0;
196             dz00             = iz0 - jz0;
197             dx10             = ix1 - jx0;
198             dy10             = iy1 - jy0;
199             dz10             = iz1 - jz0;
200             dx20             = ix2 - jx0;
201             dy20             = iy2 - jy0;
202             dz20             = iz2 - jz0;
203             dx30             = ix3 - jx0;
204             dy30             = iy3 - jy0;
205             dz30             = iz3 - jz0;
206
207             /* Calculate squared distance and things based on it */
208             rsq00            = dx00*dx00+dy00*dy00+dz00*dz00;
209             rsq10            = dx10*dx10+dy10*dy10+dz10*dz10;
210             rsq20            = dx20*dx20+dy20*dy20+dz20*dz20;
211             rsq30            = dx30*dx30+dy30*dy30+dz30*dz30;
212
213             rinv00           = gmx_invsqrt(rsq00);
214             rinv10           = gmx_invsqrt(rsq10);
215             rinv20           = gmx_invsqrt(rsq20);
216             rinv30           = gmx_invsqrt(rsq30);
217
218             rinvsq00         = rinv00*rinv00;
219             rinvsq10         = rinv10*rinv10;
220             rinvsq20         = rinv20*rinv20;
221             rinvsq30         = rinv30*rinv30;
222
223             /* Load parameters for j particles */
224             jq0              = charge[jnr+0];
225             vdwjidx0         = 3*vdwtype[jnr+0];
226
227             /**************************
228              * CALCULATE INTERACTIONS *
229              **************************/
230
231             if (rsq00<rcutoff2)
232             {
233
234             r00              = rsq00*rinv00;
235
236             c6_00            = vdwparam[vdwioffset0+vdwjidx0];
237             cexp1_00         = vdwparam[vdwioffset0+vdwjidx0+1];
238             cexp2_00         = vdwparam[vdwioffset0+vdwjidx0+2];
239
240             /* BUCKINGHAM DISPERSION/REPULSION */
241             rinvsix          = rinvsq00*rinvsq00*rinvsq00;
242             vvdw6            = c6_00*rinvsix;
243             br               = cexp2_00*r00;
244             vvdwexp          = cexp1_00*exp(-br);
245             vvdw             = (vvdwexp-cexp1_00*exp(-cexp2_00*rvdw)) - (vvdw6 - c6_00*sh_vdw_invrcut6)*(1.0/6.0);
246             fvdw             = (br*vvdwexp-vvdw6)*rinvsq00;
247
248             /* Update potential sums from outer loop */
249             vvdwsum         += vvdw;
250
251             fscal            = fvdw;
252
253             /* Calculate temporary vectorial force */
254             tx               = fscal*dx00;
255             ty               = fscal*dy00;
256             tz               = fscal*dz00;
257
258             /* Update vectorial force */
259             fix0            += tx;
260             fiy0            += ty;
261             fiz0            += tz;
262             f[j_coord_offset+DIM*0+XX] -= tx;
263             f[j_coord_offset+DIM*0+YY] -= ty;
264             f[j_coord_offset+DIM*0+ZZ] -= tz;
265
266             }
267
268             /**************************
269              * CALCULATE INTERACTIONS *
270              **************************/
271
272             if (rsq10<rcutoff2)
273             {
274
275             r10              = rsq10*rinv10;
276
277             qq10             = iq1*jq0;
278
279             /* EWALD ELECTROSTATICS */
280
281             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
282             ewrt             = r10*ewtabscale;
283             ewitab           = ewrt;
284             eweps            = ewrt-ewitab;
285             ewitab           = 4*ewitab;
286             felec            = ewtab[ewitab]+eweps*ewtab[ewitab+1];
287             velec            = qq10*((rinv10-sh_ewald)-(ewtab[ewitab+2]-ewtabhalfspace*eweps*(ewtab[ewitab]+felec)));
288             felec            = qq10*rinv10*(rinvsq10-felec);
289
290             /* Update potential sums from outer loop */
291             velecsum        += velec;
292
293             fscal            = felec;
294
295             /* Calculate temporary vectorial force */
296             tx               = fscal*dx10;
297             ty               = fscal*dy10;
298             tz               = fscal*dz10;
299
300             /* Update vectorial force */
301             fix1            += tx;
302             fiy1            += ty;
303             fiz1            += tz;
304             f[j_coord_offset+DIM*0+XX] -= tx;
305             f[j_coord_offset+DIM*0+YY] -= ty;
306             f[j_coord_offset+DIM*0+ZZ] -= tz;
307
308             }
309
310             /**************************
311              * CALCULATE INTERACTIONS *
312              **************************/
313
314             if (rsq20<rcutoff2)
315             {
316
317             r20              = rsq20*rinv20;
318
319             qq20             = iq2*jq0;
320
321             /* EWALD ELECTROSTATICS */
322
323             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
324             ewrt             = r20*ewtabscale;
325             ewitab           = ewrt;
326             eweps            = ewrt-ewitab;
327             ewitab           = 4*ewitab;
328             felec            = ewtab[ewitab]+eweps*ewtab[ewitab+1];
329             velec            = qq20*((rinv20-sh_ewald)-(ewtab[ewitab+2]-ewtabhalfspace*eweps*(ewtab[ewitab]+felec)));
330             felec            = qq20*rinv20*(rinvsq20-felec);
331
332             /* Update potential sums from outer loop */
333             velecsum        += velec;
334
335             fscal            = felec;
336
337             /* Calculate temporary vectorial force */
338             tx               = fscal*dx20;
339             ty               = fscal*dy20;
340             tz               = fscal*dz20;
341
342             /* Update vectorial force */
343             fix2            += tx;
344             fiy2            += ty;
345             fiz2            += tz;
346             f[j_coord_offset+DIM*0+XX] -= tx;
347             f[j_coord_offset+DIM*0+YY] -= ty;
348             f[j_coord_offset+DIM*0+ZZ] -= tz;
349
350             }
351
352             /**************************
353              * CALCULATE INTERACTIONS *
354              **************************/
355
356             if (rsq30<rcutoff2)
357             {
358
359             r30              = rsq30*rinv30;
360
361             qq30             = iq3*jq0;
362
363             /* EWALD ELECTROSTATICS */
364
365             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
366             ewrt             = r30*ewtabscale;
367             ewitab           = ewrt;
368             eweps            = ewrt-ewitab;
369             ewitab           = 4*ewitab;
370             felec            = ewtab[ewitab]+eweps*ewtab[ewitab+1];
371             velec            = qq30*((rinv30-sh_ewald)-(ewtab[ewitab+2]-ewtabhalfspace*eweps*(ewtab[ewitab]+felec)));
372             felec            = qq30*rinv30*(rinvsq30-felec);
373
374             /* Update potential sums from outer loop */
375             velecsum        += velec;
376
377             fscal            = felec;
378
379             /* Calculate temporary vectorial force */
380             tx               = fscal*dx30;
381             ty               = fscal*dy30;
382             tz               = fscal*dz30;
383
384             /* Update vectorial force */
385             fix3            += tx;
386             fiy3            += ty;
387             fiz3            += tz;
388             f[j_coord_offset+DIM*0+XX] -= tx;
389             f[j_coord_offset+DIM*0+YY] -= ty;
390             f[j_coord_offset+DIM*0+ZZ] -= tz;
391
392             }
393
394             /* Inner loop uses 218 flops */
395         }
396         /* End of innermost loop */
397
398         tx = ty = tz = 0;
399         f[i_coord_offset+DIM*0+XX] += fix0;
400         f[i_coord_offset+DIM*0+YY] += fiy0;
401         f[i_coord_offset+DIM*0+ZZ] += fiz0;
402         tx                         += fix0;
403         ty                         += fiy0;
404         tz                         += fiz0;
405         f[i_coord_offset+DIM*1+XX] += fix1;
406         f[i_coord_offset+DIM*1+YY] += fiy1;
407         f[i_coord_offset+DIM*1+ZZ] += fiz1;
408         tx                         += fix1;
409         ty                         += fiy1;
410         tz                         += fiz1;
411         f[i_coord_offset+DIM*2+XX] += fix2;
412         f[i_coord_offset+DIM*2+YY] += fiy2;
413         f[i_coord_offset+DIM*2+ZZ] += fiz2;
414         tx                         += fix2;
415         ty                         += fiy2;
416         tz                         += fiz2;
417         f[i_coord_offset+DIM*3+XX] += fix3;
418         f[i_coord_offset+DIM*3+YY] += fiy3;
419         f[i_coord_offset+DIM*3+ZZ] += fiz3;
420         tx                         += fix3;
421         ty                         += fiy3;
422         tz                         += fiz3;
423         fshift[i_shift_offset+XX]  += tx;
424         fshift[i_shift_offset+YY]  += ty;
425         fshift[i_shift_offset+ZZ]  += tz;
426
427         ggid                        = gid[iidx];
428         /* Update potential energies */
429         kernel_data->energygrp_elec[ggid] += velecsum;
430         kernel_data->energygrp_vdw[ggid] += vvdwsum;
431
432         /* Increment number of inner iterations */
433         inneriter                  += j_index_end - j_index_start;
434
435         /* Outer loop uses 41 flops */
436     }
437
438     /* Increment number of outer iterations */
439     outeriter        += nri;
440
441     /* Update outer/inner flops */
442
443     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W4_VF,outeriter*41 + inneriter*218);
444 }
445 /*
446  * Gromacs nonbonded kernel:   nb_kernel_ElecEwSh_VdwBhamSh_GeomW4P1_F_c
447  * Electrostatics interaction: Ewald
448  * VdW interaction:            Buckingham
449  * Geometry:                   Water4-Particle
450  * Calculate force/pot:        Force
451  */
452 void
453 nb_kernel_ElecEwSh_VdwBhamSh_GeomW4P1_F_c
454                     (t_nblist                    * gmx_restrict       nlist,
455                      rvec                        * gmx_restrict          xx,
456                      rvec                        * gmx_restrict          ff,
457                      t_forcerec                  * gmx_restrict          fr,
458                      t_mdatoms                   * gmx_restrict     mdatoms,
459                      nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
460                      t_nrnb                      * gmx_restrict        nrnb)
461 {
462     int              i_shift_offset,i_coord_offset,j_coord_offset;
463     int              j_index_start,j_index_end;
464     int              nri,inr,ggid,iidx,jidx,jnr,outeriter,inneriter;
465     real             shX,shY,shZ,tx,ty,tz,fscal,rcutoff,rcutoff2;
466     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
467     real             *shiftvec,*fshift,*x,*f;
468     int              vdwioffset0;
469     real             ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
470     int              vdwioffset1;
471     real             ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
472     int              vdwioffset2;
473     real             ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
474     int              vdwioffset3;
475     real             ix3,iy3,iz3,fix3,fiy3,fiz3,iq3,isai3;
476     int              vdwjidx0;
477     real             jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
478     real             dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00,cexp1_00,cexp2_00;
479     real             dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10,cexp1_10,cexp2_10;
480     real             dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20,cexp1_20,cexp2_20;
481     real             dx30,dy30,dz30,rsq30,rinv30,rinvsq30,r30,qq30,c6_30,c12_30,cexp1_30,cexp2_30;
482     real             velec,felec,velecsum,facel,crf,krf,krf2;
483     real             *charge;
484     int              nvdwtype;
485     real             rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,br,vvdwexp,sh_vdw_invrcut6;
486     int              *vdwtype;
487     real             *vdwparam;
488     int              ewitab;
489     real             ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace;
490     real             *ewtab;
491
492     x                = xx[0];
493     f                = ff[0];
494
495     nri              = nlist->nri;
496     iinr             = nlist->iinr;
497     jindex           = nlist->jindex;
498     jjnr             = nlist->jjnr;
499     shiftidx         = nlist->shift;
500     gid              = nlist->gid;
501     shiftvec         = fr->shift_vec[0];
502     fshift           = fr->fshift[0];
503     facel            = fr->epsfac;
504     charge           = mdatoms->chargeA;
505     nvdwtype         = fr->ntype;
506     vdwparam         = fr->nbfp;
507     vdwtype          = mdatoms->typeA;
508
509     sh_ewald         = fr->ic->sh_ewald;
510     ewtab            = fr->ic->tabq_coul_F;
511     ewtabscale       = fr->ic->tabq_scale;
512     ewtabhalfspace   = 0.5/ewtabscale;
513
514     /* Setup water-specific parameters */
515     inr              = nlist->iinr[0];
516     iq1              = facel*charge[inr+1];
517     iq2              = facel*charge[inr+2];
518     iq3              = facel*charge[inr+3];
519     vdwioffset0      = 3*nvdwtype*vdwtype[inr+0];
520
521     /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
522     rcutoff          = fr->rcoulomb;
523     rcutoff2         = rcutoff*rcutoff;
524
525     sh_vdw_invrcut6  = fr->ic->sh_invrc6;
526     rvdw             = fr->rvdw;
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         ix3              = shX + x[i_coord_offset+DIM*3+XX];
559         iy3              = shY + x[i_coord_offset+DIM*3+YY];
560         iz3              = shZ + x[i_coord_offset+DIM*3+ZZ];
561
562         fix0             = 0.0;
563         fiy0             = 0.0;
564         fiz0             = 0.0;
565         fix1             = 0.0;
566         fiy1             = 0.0;
567         fiz1             = 0.0;
568         fix2             = 0.0;
569         fiy2             = 0.0;
570         fiz2             = 0.0;
571         fix3             = 0.0;
572         fiy3             = 0.0;
573         fiz3             = 0.0;
574
575         /* Start inner kernel loop */
576         for(jidx=j_index_start; jidx<j_index_end; jidx++)
577         {
578             /* Get j neighbor index, and coordinate index */
579             jnr              = jjnr[jidx];
580             j_coord_offset   = DIM*jnr;
581
582             /* load j atom coordinates */
583             jx0              = x[j_coord_offset+DIM*0+XX];
584             jy0              = x[j_coord_offset+DIM*0+YY];
585             jz0              = x[j_coord_offset+DIM*0+ZZ];
586
587             /* Calculate displacement vector */
588             dx00             = ix0 - jx0;
589             dy00             = iy0 - jy0;
590             dz00             = iz0 - jz0;
591             dx10             = ix1 - jx0;
592             dy10             = iy1 - jy0;
593             dz10             = iz1 - jz0;
594             dx20             = ix2 - jx0;
595             dy20             = iy2 - jy0;
596             dz20             = iz2 - jz0;
597             dx30             = ix3 - jx0;
598             dy30             = iy3 - jy0;
599             dz30             = iz3 - jz0;
600
601             /* Calculate squared distance and things based on it */
602             rsq00            = dx00*dx00+dy00*dy00+dz00*dz00;
603             rsq10            = dx10*dx10+dy10*dy10+dz10*dz10;
604             rsq20            = dx20*dx20+dy20*dy20+dz20*dz20;
605             rsq30            = dx30*dx30+dy30*dy30+dz30*dz30;
606
607             rinv00           = gmx_invsqrt(rsq00);
608             rinv10           = gmx_invsqrt(rsq10);
609             rinv20           = gmx_invsqrt(rsq20);
610             rinv30           = gmx_invsqrt(rsq30);
611
612             rinvsq00         = rinv00*rinv00;
613             rinvsq10         = rinv10*rinv10;
614             rinvsq20         = rinv20*rinv20;
615             rinvsq30         = rinv30*rinv30;
616
617             /* Load parameters for j particles */
618             jq0              = charge[jnr+0];
619             vdwjidx0         = 3*vdwtype[jnr+0];
620
621             /**************************
622              * CALCULATE INTERACTIONS *
623              **************************/
624
625             if (rsq00<rcutoff2)
626             {
627
628             r00              = rsq00*rinv00;
629
630             c6_00            = vdwparam[vdwioffset0+vdwjidx0];
631             cexp1_00         = vdwparam[vdwioffset0+vdwjidx0+1];
632             cexp2_00         = vdwparam[vdwioffset0+vdwjidx0+2];
633
634             /* BUCKINGHAM DISPERSION/REPULSION */
635             rinvsix          = rinvsq00*rinvsq00*rinvsq00;
636             vvdw6            = c6_00*rinvsix;
637             br               = cexp2_00*r00;
638             vvdwexp          = cexp1_00*exp(-br);
639             fvdw             = (br*vvdwexp-vvdw6)*rinvsq00;
640
641             fscal            = fvdw;
642
643             /* Calculate temporary vectorial force */
644             tx               = fscal*dx00;
645             ty               = fscal*dy00;
646             tz               = fscal*dz00;
647
648             /* Update vectorial force */
649             fix0            += tx;
650             fiy0            += ty;
651             fiz0            += 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
658             /**************************
659              * CALCULATE INTERACTIONS *
660              **************************/
661
662             if (rsq10<rcutoff2)
663             {
664
665             r10              = rsq10*rinv10;
666
667             qq10             = iq1*jq0;
668
669             /* EWALD ELECTROSTATICS */
670
671             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
672             ewrt             = r10*ewtabscale;
673             ewitab           = ewrt;
674             eweps            = ewrt-ewitab;
675             felec            = (1.0-eweps)*ewtab[ewitab]+eweps*ewtab[ewitab+1];
676             felec            = qq10*rinv10*(rinvsq10-felec);
677
678             fscal            = felec;
679
680             /* Calculate temporary vectorial force */
681             tx               = fscal*dx10;
682             ty               = fscal*dy10;
683             tz               = fscal*dz10;
684
685             /* Update vectorial force */
686             fix1            += tx;
687             fiy1            += ty;
688             fiz1            += 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
695             /**************************
696              * CALCULATE INTERACTIONS *
697              **************************/
698
699             if (rsq20<rcutoff2)
700             {
701
702             r20              = rsq20*rinv20;
703
704             qq20             = iq2*jq0;
705
706             /* EWALD ELECTROSTATICS */
707
708             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
709             ewrt             = r20*ewtabscale;
710             ewitab           = ewrt;
711             eweps            = ewrt-ewitab;
712             felec            = (1.0-eweps)*ewtab[ewitab]+eweps*ewtab[ewitab+1];
713             felec            = qq20*rinv20*(rinvsq20-felec);
714
715             fscal            = felec;
716
717             /* Calculate temporary vectorial force */
718             tx               = fscal*dx20;
719             ty               = fscal*dy20;
720             tz               = fscal*dz20;
721
722             /* Update vectorial force */
723             fix2            += tx;
724             fiy2            += ty;
725             fiz2            += 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             }
731
732             /**************************
733              * CALCULATE INTERACTIONS *
734              **************************/
735
736             if (rsq30<rcutoff2)
737             {
738
739             r30              = rsq30*rinv30;
740
741             qq30             = iq3*jq0;
742
743             /* EWALD ELECTROSTATICS */
744
745             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
746             ewrt             = r30*ewtabscale;
747             ewitab           = ewrt;
748             eweps            = ewrt-ewitab;
749             felec            = (1.0-eweps)*ewtab[ewitab]+eweps*ewtab[ewitab+1];
750             felec            = qq30*rinv30*(rinvsq30-felec);
751
752             fscal            = felec;
753
754             /* Calculate temporary vectorial force */
755             tx               = fscal*dx30;
756             ty               = fscal*dy30;
757             tz               = fscal*dz30;
758
759             /* Update vectorial force */
760             fix3            += tx;
761             fiy3            += ty;
762             fiz3            += tz;
763             f[j_coord_offset+DIM*0+XX] -= tx;
764             f[j_coord_offset+DIM*0+YY] -= ty;
765             f[j_coord_offset+DIM*0+ZZ] -= tz;
766
767             }
768
769             /* Inner loop uses 160 flops */
770         }
771         /* End of innermost loop */
772
773         tx = ty = tz = 0;
774         f[i_coord_offset+DIM*0+XX] += fix0;
775         f[i_coord_offset+DIM*0+YY] += fiy0;
776         f[i_coord_offset+DIM*0+ZZ] += fiz0;
777         tx                         += fix0;
778         ty                         += fiy0;
779         tz                         += fiz0;
780         f[i_coord_offset+DIM*1+XX] += fix1;
781         f[i_coord_offset+DIM*1+YY] += fiy1;
782         f[i_coord_offset+DIM*1+ZZ] += fiz1;
783         tx                         += fix1;
784         ty                         += fiy1;
785         tz                         += fiz1;
786         f[i_coord_offset+DIM*2+XX] += fix2;
787         f[i_coord_offset+DIM*2+YY] += fiy2;
788         f[i_coord_offset+DIM*2+ZZ] += fiz2;
789         tx                         += fix2;
790         ty                         += fiy2;
791         tz                         += fiz2;
792         f[i_coord_offset+DIM*3+XX] += fix3;
793         f[i_coord_offset+DIM*3+YY] += fiy3;
794         f[i_coord_offset+DIM*3+ZZ] += fiz3;
795         tx                         += fix3;
796         ty                         += fiy3;
797         tz                         += fiz3;
798         fshift[i_shift_offset+XX]  += tx;
799         fshift[i_shift_offset+YY]  += ty;
800         fshift[i_shift_offset+ZZ]  += tz;
801
802         /* Increment number of inner iterations */
803         inneriter                  += j_index_end - j_index_start;
804
805         /* Outer loop uses 39 flops */
806     }
807
808     /* Increment number of outer iterations */
809     outeriter        += nri;
810
811     /* Update outer/inner flops */
812
813     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W4_F,outeriter*39 + inneriter*160);
814 }