Merge release-5-0 into master
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_c / nb_kernel_ElecCSTab_VdwNone_GeomW4P1_c.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2013,2014, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35 /*
36  * Note: this file was generated by the GROMACS c kernel generator.
37  */
38 #include "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_ElecCSTab_VdwNone_GeomW4P1_VF_c
51  * Electrostatics interaction: CubicSplineTable
52  * VdW interaction:            None
53  * Geometry:                   Water4-Particle
54  * Calculate force/pot:        PotentialAndForce
55  */
56 void
57 nb_kernel_ElecCSTab_VdwNone_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              vdwioffset1;
73     real             ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
74     int              vdwioffset2;
75     real             ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
76     int              vdwioffset3;
77     real             ix3,iy3,iz3,fix3,fiy3,fiz3,iq3,isai3;
78     int              vdwjidx0;
79     real             jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
80     real             dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10,cexp1_10,cexp2_10;
81     real             dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20,cexp1_20,cexp2_20;
82     real             dx30,dy30,dz30,rsq30,rinv30,rinvsq30,r30,qq30,c6_30,c12_30,cexp1_30,cexp2_30;
83     real             velec,felec,velecsum,facel,crf,krf,krf2;
84     real             *charge;
85     int              vfitab;
86     real             rt,vfeps,vftabscale,Y,F,Geps,Heps2,Fp,VV,FF;
87     real             *vftab;
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     vftab            = kernel_data->table_elec->data;
104     vftabscale       = kernel_data->table_elec->scale;
105
106     /* Setup water-specific parameters */
107     inr              = nlist->iinr[0];
108     iq1              = facel*charge[inr+1];
109     iq2              = facel*charge[inr+2];
110     iq3              = facel*charge[inr+3];
111
112     outeriter        = 0;
113     inneriter        = 0;
114
115     /* Start outer loop over neighborlists */
116     for(iidx=0; iidx<nri; iidx++)
117     {
118         /* Load shift vector for this list */
119         i_shift_offset   = DIM*shiftidx[iidx];
120         shX              = shiftvec[i_shift_offset+XX];
121         shY              = shiftvec[i_shift_offset+YY];
122         shZ              = shiftvec[i_shift_offset+ZZ];
123
124         /* Load limits for loop over neighbors */
125         j_index_start    = jindex[iidx];
126         j_index_end      = jindex[iidx+1];
127
128         /* Get outer coordinate index */
129         inr              = iinr[iidx];
130         i_coord_offset   = DIM*inr;
131
132         /* Load i particle coords and add shift vector */
133         ix1              = shX + x[i_coord_offset+DIM*1+XX];
134         iy1              = shY + x[i_coord_offset+DIM*1+YY];
135         iz1              = shZ + x[i_coord_offset+DIM*1+ZZ];
136         ix2              = shX + x[i_coord_offset+DIM*2+XX];
137         iy2              = shY + x[i_coord_offset+DIM*2+YY];
138         iz2              = shZ + x[i_coord_offset+DIM*2+ZZ];
139         ix3              = shX + x[i_coord_offset+DIM*3+XX];
140         iy3              = shY + x[i_coord_offset+DIM*3+YY];
141         iz3              = shZ + x[i_coord_offset+DIM*3+ZZ];
142
143         fix1             = 0.0;
144         fiy1             = 0.0;
145         fiz1             = 0.0;
146         fix2             = 0.0;
147         fiy2             = 0.0;
148         fiz2             = 0.0;
149         fix3             = 0.0;
150         fiy3             = 0.0;
151         fiz3             = 0.0;
152
153         /* Reset potential sums */
154         velecsum         = 0.0;
155
156         /* Start inner kernel loop */
157         for(jidx=j_index_start; jidx<j_index_end; jidx++)
158         {
159             /* Get j neighbor index, and coordinate index */
160             jnr              = jjnr[jidx];
161             j_coord_offset   = DIM*jnr;
162
163             /* load j atom coordinates */
164             jx0              = x[j_coord_offset+DIM*0+XX];
165             jy0              = x[j_coord_offset+DIM*0+YY];
166             jz0              = x[j_coord_offset+DIM*0+ZZ];
167
168             /* Calculate displacement vector */
169             dx10             = ix1 - jx0;
170             dy10             = iy1 - jy0;
171             dz10             = iz1 - jz0;
172             dx20             = ix2 - jx0;
173             dy20             = iy2 - jy0;
174             dz20             = iz2 - jz0;
175             dx30             = ix3 - jx0;
176             dy30             = iy3 - jy0;
177             dz30             = iz3 - jz0;
178
179             /* Calculate squared distance and things based on it */
180             rsq10            = dx10*dx10+dy10*dy10+dz10*dz10;
181             rsq20            = dx20*dx20+dy20*dy20+dz20*dz20;
182             rsq30            = dx30*dx30+dy30*dy30+dz30*dz30;
183
184             rinv10           = gmx_invsqrt(rsq10);
185             rinv20           = gmx_invsqrt(rsq20);
186             rinv30           = gmx_invsqrt(rsq30);
187
188             /* Load parameters for j particles */
189             jq0              = charge[jnr+0];
190
191             /**************************
192              * CALCULATE INTERACTIONS *
193              **************************/
194
195             r10              = rsq10*rinv10;
196
197             qq10             = iq1*jq0;
198
199             /* Calculate table index by multiplying r with table scale and truncate to integer */
200             rt               = r10*vftabscale;
201             vfitab           = rt;
202             vfeps            = rt-vfitab;
203             vfitab           = 1*4*vfitab;
204
205             /* CUBIC SPLINE TABLE ELECTROSTATICS */
206             Y                = vftab[vfitab];
207             F                = vftab[vfitab+1];
208             Geps             = vfeps*vftab[vfitab+2];
209             Heps2            = vfeps*vfeps*vftab[vfitab+3];
210             Fp               = F+Geps+Heps2;
211             VV               = Y+vfeps*Fp;
212             velec            = qq10*VV;
213             FF               = Fp+Geps+2.0*Heps2;
214             felec            = -qq10*FF*vftabscale*rinv10;
215
216             /* Update potential sums from outer loop */
217             velecsum        += velec;
218
219             fscal            = felec;
220
221             /* Calculate temporary vectorial force */
222             tx               = fscal*dx10;
223             ty               = fscal*dy10;
224             tz               = fscal*dz10;
225
226             /* Update vectorial force */
227             fix1            += tx;
228             fiy1            += ty;
229             fiz1            += tz;
230             f[j_coord_offset+DIM*0+XX] -= tx;
231             f[j_coord_offset+DIM*0+YY] -= ty;
232             f[j_coord_offset+DIM*0+ZZ] -= tz;
233
234             /**************************
235              * CALCULATE INTERACTIONS *
236              **************************/
237
238             r20              = rsq20*rinv20;
239
240             qq20             = iq2*jq0;
241
242             /* Calculate table index by multiplying r with table scale and truncate to integer */
243             rt               = r20*vftabscale;
244             vfitab           = rt;
245             vfeps            = rt-vfitab;
246             vfitab           = 1*4*vfitab;
247
248             /* CUBIC SPLINE TABLE ELECTROSTATICS */
249             Y                = vftab[vfitab];
250             F                = vftab[vfitab+1];
251             Geps             = vfeps*vftab[vfitab+2];
252             Heps2            = vfeps*vfeps*vftab[vfitab+3];
253             Fp               = F+Geps+Heps2;
254             VV               = Y+vfeps*Fp;
255             velec            = qq20*VV;
256             FF               = Fp+Geps+2.0*Heps2;
257             felec            = -qq20*FF*vftabscale*rinv20;
258
259             /* Update potential sums from outer loop */
260             velecsum        += velec;
261
262             fscal            = felec;
263
264             /* Calculate temporary vectorial force */
265             tx               = fscal*dx20;
266             ty               = fscal*dy20;
267             tz               = fscal*dz20;
268
269             /* Update vectorial force */
270             fix2            += tx;
271             fiy2            += ty;
272             fiz2            += tz;
273             f[j_coord_offset+DIM*0+XX] -= tx;
274             f[j_coord_offset+DIM*0+YY] -= ty;
275             f[j_coord_offset+DIM*0+ZZ] -= tz;
276
277             /**************************
278              * CALCULATE INTERACTIONS *
279              **************************/
280
281             r30              = rsq30*rinv30;
282
283             qq30             = iq3*jq0;
284
285             /* Calculate table index by multiplying r with table scale and truncate to integer */
286             rt               = r30*vftabscale;
287             vfitab           = rt;
288             vfeps            = rt-vfitab;
289             vfitab           = 1*4*vfitab;
290
291             /* CUBIC SPLINE TABLE ELECTROSTATICS */
292             Y                = vftab[vfitab];
293             F                = vftab[vfitab+1];
294             Geps             = vfeps*vftab[vfitab+2];
295             Heps2            = vfeps*vfeps*vftab[vfitab+3];
296             Fp               = F+Geps+Heps2;
297             VV               = Y+vfeps*Fp;
298             velec            = qq30*VV;
299             FF               = Fp+Geps+2.0*Heps2;
300             felec            = -qq30*FF*vftabscale*rinv30;
301
302             /* Update potential sums from outer loop */
303             velecsum        += velec;
304
305             fscal            = felec;
306
307             /* Calculate temporary vectorial force */
308             tx               = fscal*dx30;
309             ty               = fscal*dy30;
310             tz               = fscal*dz30;
311
312             /* Update vectorial force */
313             fix3            += tx;
314             fiy3            += ty;
315             fiz3            += tz;
316             f[j_coord_offset+DIM*0+XX] -= tx;
317             f[j_coord_offset+DIM*0+YY] -= ty;
318             f[j_coord_offset+DIM*0+ZZ] -= tz;
319
320             /* Inner loop uses 126 flops */
321         }
322         /* End of innermost loop */
323
324         tx = ty = tz = 0;
325         f[i_coord_offset+DIM*1+XX] += fix1;
326         f[i_coord_offset+DIM*1+YY] += fiy1;
327         f[i_coord_offset+DIM*1+ZZ] += fiz1;
328         tx                         += fix1;
329         ty                         += fiy1;
330         tz                         += fiz1;
331         f[i_coord_offset+DIM*2+XX] += fix2;
332         f[i_coord_offset+DIM*2+YY] += fiy2;
333         f[i_coord_offset+DIM*2+ZZ] += fiz2;
334         tx                         += fix2;
335         ty                         += fiy2;
336         tz                         += fiz2;
337         f[i_coord_offset+DIM*3+XX] += fix3;
338         f[i_coord_offset+DIM*3+YY] += fiy3;
339         f[i_coord_offset+DIM*3+ZZ] += fiz3;
340         tx                         += fix3;
341         ty                         += fiy3;
342         tz                         += fiz3;
343         fshift[i_shift_offset+XX]  += tx;
344         fshift[i_shift_offset+YY]  += ty;
345         fshift[i_shift_offset+ZZ]  += tz;
346
347         ggid                        = gid[iidx];
348         /* Update potential energies */
349         kernel_data->energygrp_elec[ggid] += velecsum;
350
351         /* Increment number of inner iterations */
352         inneriter                  += j_index_end - j_index_start;
353
354         /* Outer loop uses 31 flops */
355     }
356
357     /* Increment number of outer iterations */
358     outeriter        += nri;
359
360     /* Update outer/inner flops */
361
362     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_W4_VF,outeriter*31 + inneriter*126);
363 }
364 /*
365  * Gromacs nonbonded kernel:   nb_kernel_ElecCSTab_VdwNone_GeomW4P1_F_c
366  * Electrostatics interaction: CubicSplineTable
367  * VdW interaction:            None
368  * Geometry:                   Water4-Particle
369  * Calculate force/pot:        Force
370  */
371 void
372 nb_kernel_ElecCSTab_VdwNone_GeomW4P1_F_c
373                     (t_nblist                    * gmx_restrict       nlist,
374                      rvec                        * gmx_restrict          xx,
375                      rvec                        * gmx_restrict          ff,
376                      t_forcerec                  * gmx_restrict          fr,
377                      t_mdatoms                   * gmx_restrict     mdatoms,
378                      nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
379                      t_nrnb                      * gmx_restrict        nrnb)
380 {
381     int              i_shift_offset,i_coord_offset,j_coord_offset;
382     int              j_index_start,j_index_end;
383     int              nri,inr,ggid,iidx,jidx,jnr,outeriter,inneriter;
384     real             shX,shY,shZ,tx,ty,tz,fscal,rcutoff,rcutoff2;
385     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
386     real             *shiftvec,*fshift,*x,*f;
387     int              vdwioffset1;
388     real             ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
389     int              vdwioffset2;
390     real             ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
391     int              vdwioffset3;
392     real             ix3,iy3,iz3,fix3,fiy3,fiz3,iq3,isai3;
393     int              vdwjidx0;
394     real             jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
395     real             dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10,cexp1_10,cexp2_10;
396     real             dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20,cexp1_20,cexp2_20;
397     real             dx30,dy30,dz30,rsq30,rinv30,rinvsq30,r30,qq30,c6_30,c12_30,cexp1_30,cexp2_30;
398     real             velec,felec,velecsum,facel,crf,krf,krf2;
399     real             *charge;
400     int              vfitab;
401     real             rt,vfeps,vftabscale,Y,F,Geps,Heps2,Fp,VV,FF;
402     real             *vftab;
403
404     x                = xx[0];
405     f                = ff[0];
406
407     nri              = nlist->nri;
408     iinr             = nlist->iinr;
409     jindex           = nlist->jindex;
410     jjnr             = nlist->jjnr;
411     shiftidx         = nlist->shift;
412     gid              = nlist->gid;
413     shiftvec         = fr->shift_vec[0];
414     fshift           = fr->fshift[0];
415     facel            = fr->epsfac;
416     charge           = mdatoms->chargeA;
417
418     vftab            = kernel_data->table_elec->data;
419     vftabscale       = kernel_data->table_elec->scale;
420
421     /* Setup water-specific parameters */
422     inr              = nlist->iinr[0];
423     iq1              = facel*charge[inr+1];
424     iq2              = facel*charge[inr+2];
425     iq3              = facel*charge[inr+3];
426
427     outeriter        = 0;
428     inneriter        = 0;
429
430     /* Start outer loop over neighborlists */
431     for(iidx=0; iidx<nri; iidx++)
432     {
433         /* Load shift vector for this list */
434         i_shift_offset   = DIM*shiftidx[iidx];
435         shX              = shiftvec[i_shift_offset+XX];
436         shY              = shiftvec[i_shift_offset+YY];
437         shZ              = shiftvec[i_shift_offset+ZZ];
438
439         /* Load limits for loop over neighbors */
440         j_index_start    = jindex[iidx];
441         j_index_end      = jindex[iidx+1];
442
443         /* Get outer coordinate index */
444         inr              = iinr[iidx];
445         i_coord_offset   = DIM*inr;
446
447         /* Load i particle coords and add shift vector */
448         ix1              = shX + x[i_coord_offset+DIM*1+XX];
449         iy1              = shY + x[i_coord_offset+DIM*1+YY];
450         iz1              = shZ + x[i_coord_offset+DIM*1+ZZ];
451         ix2              = shX + x[i_coord_offset+DIM*2+XX];
452         iy2              = shY + x[i_coord_offset+DIM*2+YY];
453         iz2              = shZ + x[i_coord_offset+DIM*2+ZZ];
454         ix3              = shX + x[i_coord_offset+DIM*3+XX];
455         iy3              = shY + x[i_coord_offset+DIM*3+YY];
456         iz3              = shZ + x[i_coord_offset+DIM*3+ZZ];
457
458         fix1             = 0.0;
459         fiy1             = 0.0;
460         fiz1             = 0.0;
461         fix2             = 0.0;
462         fiy2             = 0.0;
463         fiz2             = 0.0;
464         fix3             = 0.0;
465         fiy3             = 0.0;
466         fiz3             = 0.0;
467
468         /* Start inner kernel loop */
469         for(jidx=j_index_start; jidx<j_index_end; jidx++)
470         {
471             /* Get j neighbor index, and coordinate index */
472             jnr              = jjnr[jidx];
473             j_coord_offset   = DIM*jnr;
474
475             /* load j atom coordinates */
476             jx0              = x[j_coord_offset+DIM*0+XX];
477             jy0              = x[j_coord_offset+DIM*0+YY];
478             jz0              = x[j_coord_offset+DIM*0+ZZ];
479
480             /* Calculate displacement vector */
481             dx10             = ix1 - jx0;
482             dy10             = iy1 - jy0;
483             dz10             = iz1 - jz0;
484             dx20             = ix2 - jx0;
485             dy20             = iy2 - jy0;
486             dz20             = iz2 - jz0;
487             dx30             = ix3 - jx0;
488             dy30             = iy3 - jy0;
489             dz30             = iz3 - jz0;
490
491             /* Calculate squared distance and things based on it */
492             rsq10            = dx10*dx10+dy10*dy10+dz10*dz10;
493             rsq20            = dx20*dx20+dy20*dy20+dz20*dz20;
494             rsq30            = dx30*dx30+dy30*dy30+dz30*dz30;
495
496             rinv10           = gmx_invsqrt(rsq10);
497             rinv20           = gmx_invsqrt(rsq20);
498             rinv30           = gmx_invsqrt(rsq30);
499
500             /* Load parameters for j particles */
501             jq0              = charge[jnr+0];
502
503             /**************************
504              * CALCULATE INTERACTIONS *
505              **************************/
506
507             r10              = rsq10*rinv10;
508
509             qq10             = iq1*jq0;
510
511             /* Calculate table index by multiplying r with table scale and truncate to integer */
512             rt               = r10*vftabscale;
513             vfitab           = rt;
514             vfeps            = rt-vfitab;
515             vfitab           = 1*4*vfitab;
516
517             /* CUBIC SPLINE TABLE ELECTROSTATICS */
518             F                = vftab[vfitab+1];
519             Geps             = vfeps*vftab[vfitab+2];
520             Heps2            = vfeps*vfeps*vftab[vfitab+3];
521             Fp               = F+Geps+Heps2;
522             FF               = Fp+Geps+2.0*Heps2;
523             felec            = -qq10*FF*vftabscale*rinv10;
524
525             fscal            = felec;
526
527             /* Calculate temporary vectorial force */
528             tx               = fscal*dx10;
529             ty               = fscal*dy10;
530             tz               = fscal*dz10;
531
532             /* Update vectorial force */
533             fix1            += tx;
534             fiy1            += ty;
535             fiz1            += tz;
536             f[j_coord_offset+DIM*0+XX] -= tx;
537             f[j_coord_offset+DIM*0+YY] -= ty;
538             f[j_coord_offset+DIM*0+ZZ] -= tz;
539
540             /**************************
541              * CALCULATE INTERACTIONS *
542              **************************/
543
544             r20              = rsq20*rinv20;
545
546             qq20             = iq2*jq0;
547
548             /* Calculate table index by multiplying r with table scale and truncate to integer */
549             rt               = r20*vftabscale;
550             vfitab           = rt;
551             vfeps            = rt-vfitab;
552             vfitab           = 1*4*vfitab;
553
554             /* CUBIC SPLINE TABLE ELECTROSTATICS */
555             F                = vftab[vfitab+1];
556             Geps             = vfeps*vftab[vfitab+2];
557             Heps2            = vfeps*vfeps*vftab[vfitab+3];
558             Fp               = F+Geps+Heps2;
559             FF               = Fp+Geps+2.0*Heps2;
560             felec            = -qq20*FF*vftabscale*rinv20;
561
562             fscal            = felec;
563
564             /* Calculate temporary vectorial force */
565             tx               = fscal*dx20;
566             ty               = fscal*dy20;
567             tz               = fscal*dz20;
568
569             /* Update vectorial force */
570             fix2            += tx;
571             fiy2            += ty;
572             fiz2            += tz;
573             f[j_coord_offset+DIM*0+XX] -= tx;
574             f[j_coord_offset+DIM*0+YY] -= ty;
575             f[j_coord_offset+DIM*0+ZZ] -= tz;
576
577             /**************************
578              * CALCULATE INTERACTIONS *
579              **************************/
580
581             r30              = rsq30*rinv30;
582
583             qq30             = iq3*jq0;
584
585             /* Calculate table index by multiplying r with table scale and truncate to integer */
586             rt               = r30*vftabscale;
587             vfitab           = rt;
588             vfeps            = rt-vfitab;
589             vfitab           = 1*4*vfitab;
590
591             /* CUBIC SPLINE TABLE ELECTROSTATICS */
592             F                = vftab[vfitab+1];
593             Geps             = vfeps*vftab[vfitab+2];
594             Heps2            = vfeps*vfeps*vftab[vfitab+3];
595             Fp               = F+Geps+Heps2;
596             FF               = Fp+Geps+2.0*Heps2;
597             felec            = -qq30*FF*vftabscale*rinv30;
598
599             fscal            = felec;
600
601             /* Calculate temporary vectorial force */
602             tx               = fscal*dx30;
603             ty               = fscal*dy30;
604             tz               = fscal*dz30;
605
606             /* Update vectorial force */
607             fix3            += tx;
608             fiy3            += ty;
609             fiz3            += tz;
610             f[j_coord_offset+DIM*0+XX] -= tx;
611             f[j_coord_offset+DIM*0+YY] -= ty;
612             f[j_coord_offset+DIM*0+ZZ] -= tz;
613
614             /* Inner loop uses 114 flops */
615         }
616         /* End of innermost loop */
617
618         tx = ty = tz = 0;
619         f[i_coord_offset+DIM*1+XX] += fix1;
620         f[i_coord_offset+DIM*1+YY] += fiy1;
621         f[i_coord_offset+DIM*1+ZZ] += fiz1;
622         tx                         += fix1;
623         ty                         += fiy1;
624         tz                         += fiz1;
625         f[i_coord_offset+DIM*2+XX] += fix2;
626         f[i_coord_offset+DIM*2+YY] += fiy2;
627         f[i_coord_offset+DIM*2+ZZ] += fiz2;
628         tx                         += fix2;
629         ty                         += fiy2;
630         tz                         += fiz2;
631         f[i_coord_offset+DIM*3+XX] += fix3;
632         f[i_coord_offset+DIM*3+YY] += fiy3;
633         f[i_coord_offset+DIM*3+ZZ] += fiz3;
634         tx                         += fix3;
635         ty                         += fiy3;
636         tz                         += fiz3;
637         fshift[i_shift_offset+XX]  += tx;
638         fshift[i_shift_offset+YY]  += ty;
639         fshift[i_shift_offset+ZZ]  += tz;
640
641         /* Increment number of inner iterations */
642         inneriter                  += j_index_end - j_index_start;
643
644         /* Outer loop uses 30 flops */
645     }
646
647     /* Increment number of outer iterations */
648     outeriter        += nri;
649
650     /* Update outer/inner flops */
651
652     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_W4_F,outeriter*30 + inneriter*114);
653 }