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