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