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