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