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