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