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