e49fbaec1744a097c1c457db73116a52050a68a8
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_c / nb_kernel_ElecRF_VdwCSTab_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_ElecRF_VdwCSTab_GeomW4P1_VF_c
49  * Electrostatics interaction: ReactionField
50  * VdW interaction:            CubicSplineTable
51  * Geometry:                   Water4-Particle
52  * Calculate force/pot:        PotentialAndForce
53  */
54 void
55 nb_kernel_ElecRF_VdwCSTab_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     int              vfitab;
91     real             rt,vfeps,vftabscale,Y,F,Geps,Heps2,Fp,VV,FF;
92     real             *vftab;
93
94     x                = xx[0];
95     f                = ff[0];
96
97     nri              = nlist->nri;
98     iinr             = nlist->iinr;
99     jindex           = nlist->jindex;
100     jjnr             = nlist->jjnr;
101     shiftidx         = nlist->shift;
102     gid              = nlist->gid;
103     shiftvec         = fr->shift_vec[0];
104     fshift           = fr->fshift[0];
105     facel            = fr->epsfac;
106     charge           = mdatoms->chargeA;
107     krf              = fr->ic->k_rf;
108     krf2             = krf*2.0;
109     crf              = fr->ic->c_rf;
110     nvdwtype         = fr->ntype;
111     vdwparam         = fr->nbfp;
112     vdwtype          = mdatoms->typeA;
113
114     vftab            = kernel_data->table_vdw->data;
115     vftabscale       = kernel_data->table_vdw->scale;
116
117     /* Setup water-specific parameters */
118     inr              = nlist->iinr[0];
119     iq1              = facel*charge[inr+1];
120     iq2              = facel*charge[inr+2];
121     iq3              = facel*charge[inr+3];
122     vdwioffset0      = 2*nvdwtype*vdwtype[inr+0];
123
124     outeriter        = 0;
125     inneriter        = 0;
126
127     /* Start outer loop over neighborlists */
128     for(iidx=0; iidx<nri; iidx++)
129     {
130         /* Load shift vector for this list */
131         i_shift_offset   = DIM*shiftidx[iidx];
132         shX              = shiftvec[i_shift_offset+XX];
133         shY              = shiftvec[i_shift_offset+YY];
134         shZ              = shiftvec[i_shift_offset+ZZ];
135
136         /* Load limits for loop over neighbors */
137         j_index_start    = jindex[iidx];
138         j_index_end      = jindex[iidx+1];
139
140         /* Get outer coordinate index */
141         inr              = iinr[iidx];
142         i_coord_offset   = DIM*inr;
143
144         /* Load i particle coords and add shift vector */
145         ix0              = shX + x[i_coord_offset+DIM*0+XX];
146         iy0              = shY + x[i_coord_offset+DIM*0+YY];
147         iz0              = shZ + x[i_coord_offset+DIM*0+ZZ];
148         ix1              = shX + x[i_coord_offset+DIM*1+XX];
149         iy1              = shY + x[i_coord_offset+DIM*1+YY];
150         iz1              = shZ + x[i_coord_offset+DIM*1+ZZ];
151         ix2              = shX + x[i_coord_offset+DIM*2+XX];
152         iy2              = shY + x[i_coord_offset+DIM*2+YY];
153         iz2              = shZ + x[i_coord_offset+DIM*2+ZZ];
154         ix3              = shX + x[i_coord_offset+DIM*3+XX];
155         iy3              = shY + x[i_coord_offset+DIM*3+YY];
156         iz3              = shZ + x[i_coord_offset+DIM*3+ZZ];
157
158         fix0             = 0.0;
159         fiy0             = 0.0;
160         fiz0             = 0.0;
161         fix1             = 0.0;
162         fiy1             = 0.0;
163         fiz1             = 0.0;
164         fix2             = 0.0;
165         fiy2             = 0.0;
166         fiz2             = 0.0;
167         fix3             = 0.0;
168         fiy3             = 0.0;
169         fiz3             = 0.0;
170
171         /* Reset potential sums */
172         velecsum         = 0.0;
173         vvdwsum          = 0.0;
174
175         /* Start inner kernel loop */
176         for(jidx=j_index_start; jidx<j_index_end; jidx++)
177         {
178             /* Get j neighbor index, and coordinate index */
179             jnr              = jjnr[jidx];
180             j_coord_offset   = DIM*jnr;
181
182             /* load j atom coordinates */
183             jx0              = x[j_coord_offset+DIM*0+XX];
184             jy0              = x[j_coord_offset+DIM*0+YY];
185             jz0              = x[j_coord_offset+DIM*0+ZZ];
186
187             /* Calculate displacement vector */
188             dx00             = ix0 - jx0;
189             dy00             = iy0 - jy0;
190             dz00             = iz0 - jz0;
191             dx10             = ix1 - jx0;
192             dy10             = iy1 - jy0;
193             dz10             = iz1 - jz0;
194             dx20             = ix2 - jx0;
195             dy20             = iy2 - jy0;
196             dz20             = iz2 - jz0;
197             dx30             = ix3 - jx0;
198             dy30             = iy3 - jy0;
199             dz30             = iz3 - jz0;
200
201             /* Calculate squared distance and things based on it */
202             rsq00            = dx00*dx00+dy00*dy00+dz00*dz00;
203             rsq10            = dx10*dx10+dy10*dy10+dz10*dz10;
204             rsq20            = dx20*dx20+dy20*dy20+dz20*dz20;
205             rsq30            = dx30*dx30+dy30*dy30+dz30*dz30;
206
207             rinv00           = gmx_invsqrt(rsq00);
208             rinv10           = gmx_invsqrt(rsq10);
209             rinv20           = gmx_invsqrt(rsq20);
210             rinv30           = gmx_invsqrt(rsq30);
211
212             rinvsq10         = rinv10*rinv10;
213             rinvsq20         = rinv20*rinv20;
214             rinvsq30         = rinv30*rinv30;
215
216             /* Load parameters for j particles */
217             jq0              = charge[jnr+0];
218             vdwjidx0         = 2*vdwtype[jnr+0];
219
220             /**************************
221              * CALCULATE INTERACTIONS *
222              **************************/
223
224             r00              = rsq00*rinv00;
225
226             c6_00            = vdwparam[vdwioffset0+vdwjidx0];
227             c12_00           = vdwparam[vdwioffset0+vdwjidx0+1];
228
229             /* Calculate table index by multiplying r with table scale and truncate to integer */
230             rt               = r00*vftabscale;
231             vfitab           = rt;
232             vfeps            = rt-vfitab;
233             vfitab           = 2*4*vfitab;
234
235             /* CUBIC SPLINE TABLE DISPERSION */
236             vfitab          += 0;
237             Y                = vftab[vfitab];
238             F                = vftab[vfitab+1];
239             Geps             = vfeps*vftab[vfitab+2];
240             Heps2            = vfeps*vfeps*vftab[vfitab+3];
241             Fp               = F+Geps+Heps2;
242             VV               = Y+vfeps*Fp;
243             vvdw6            = c6_00*VV;
244             FF               = Fp+Geps+2.0*Heps2;
245             fvdw6            = c6_00*FF;
246
247             /* CUBIC SPLINE TABLE REPULSION */
248             Y                = vftab[vfitab+4];
249             F                = vftab[vfitab+5];
250             Geps             = vfeps*vftab[vfitab+6];
251             Heps2            = vfeps*vfeps*vftab[vfitab+7];
252             Fp               = F+Geps+Heps2;
253             VV               = Y+vfeps*Fp;
254             vvdw12           = c12_00*VV;
255             FF               = Fp+Geps+2.0*Heps2;
256             fvdw12           = c12_00*FF;
257             vvdw             = vvdw12+vvdw6;
258             fvdw             = -(fvdw6+fvdw12)*vftabscale*rinv00;
259
260             /* Update potential sums from outer loop */
261             vvdwsum         += vvdw;
262
263             fscal            = fvdw;
264
265             /* Calculate temporary vectorial force */
266             tx               = fscal*dx00;
267             ty               = fscal*dy00;
268             tz               = fscal*dz00;
269
270             /* Update vectorial force */
271             fix0            += tx;
272             fiy0            += ty;
273             fiz0            += tz;
274             f[j_coord_offset+DIM*0+XX] -= tx;
275             f[j_coord_offset+DIM*0+YY] -= ty;
276             f[j_coord_offset+DIM*0+ZZ] -= tz;
277
278             /**************************
279              * CALCULATE INTERACTIONS *
280              **************************/
281
282             qq10             = iq1*jq0;
283
284             /* REACTION-FIELD ELECTROSTATICS */
285             velec            = qq10*(rinv10+krf*rsq10-crf);
286             felec            = qq10*(rinv10*rinvsq10-krf2);
287
288             /* Update potential sums from outer loop */
289             velecsum        += velec;
290
291             fscal            = felec;
292
293             /* Calculate temporary vectorial force */
294             tx               = fscal*dx10;
295             ty               = fscal*dy10;
296             tz               = fscal*dz10;
297
298             /* Update vectorial force */
299             fix1            += tx;
300             fiy1            += ty;
301             fiz1            += tz;
302             f[j_coord_offset+DIM*0+XX] -= tx;
303             f[j_coord_offset+DIM*0+YY] -= ty;
304             f[j_coord_offset+DIM*0+ZZ] -= tz;
305
306             /**************************
307              * CALCULATE INTERACTIONS *
308              **************************/
309
310             qq20             = iq2*jq0;
311
312             /* REACTION-FIELD ELECTROSTATICS */
313             velec            = qq20*(rinv20+krf*rsq20-crf);
314             felec            = qq20*(rinv20*rinvsq20-krf2);
315
316             /* Update potential sums from outer loop */
317             velecsum        += velec;
318
319             fscal            = felec;
320
321             /* Calculate temporary vectorial force */
322             tx               = fscal*dx20;
323             ty               = fscal*dy20;
324             tz               = fscal*dz20;
325
326             /* Update vectorial force */
327             fix2            += tx;
328             fiy2            += ty;
329             fiz2            += tz;
330             f[j_coord_offset+DIM*0+XX] -= tx;
331             f[j_coord_offset+DIM*0+YY] -= ty;
332             f[j_coord_offset+DIM*0+ZZ] -= tz;
333
334             /**************************
335              * CALCULATE INTERACTIONS *
336              **************************/
337
338             qq30             = iq3*jq0;
339
340             /* REACTION-FIELD ELECTROSTATICS */
341             velec            = qq30*(rinv30+krf*rsq30-crf);
342             felec            = qq30*(rinv30*rinvsq30-krf2);
343
344             /* Update potential sums from outer loop */
345             velecsum        += velec;
346
347             fscal            = felec;
348
349             /* Calculate temporary vectorial force */
350             tx               = fscal*dx30;
351             ty               = fscal*dy30;
352             tz               = fscal*dz30;
353
354             /* Update vectorial force */
355             fix3            += tx;
356             fiy3            += ty;
357             fiz3            += tz;
358             f[j_coord_offset+DIM*0+XX] -= tx;
359             f[j_coord_offset+DIM*0+YY] -= ty;
360             f[j_coord_offset+DIM*0+ZZ] -= tz;
361
362             /* Inner loop uses 151 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*151);
412 }
413 /*
414  * Gromacs nonbonded kernel:   nb_kernel_ElecRF_VdwCSTab_GeomW4P1_F_c
415  * Electrostatics interaction: ReactionField
416  * VdW interaction:            CubicSplineTable
417  * Geometry:                   Water4-Particle
418  * Calculate force/pot:        Force
419  */
420 void
421 nb_kernel_ElecRF_VdwCSTab_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     int              vfitab;
457     real             rt,vfeps,vftabscale,Y,F,Geps,Heps2,Fp,VV,FF;
458     real             *vftab;
459
460     x                = xx[0];
461     f                = ff[0];
462
463     nri              = nlist->nri;
464     iinr             = nlist->iinr;
465     jindex           = nlist->jindex;
466     jjnr             = nlist->jjnr;
467     shiftidx         = nlist->shift;
468     gid              = nlist->gid;
469     shiftvec         = fr->shift_vec[0];
470     fshift           = fr->fshift[0];
471     facel            = fr->epsfac;
472     charge           = mdatoms->chargeA;
473     krf              = fr->ic->k_rf;
474     krf2             = krf*2.0;
475     crf              = fr->ic->c_rf;
476     nvdwtype         = fr->ntype;
477     vdwparam         = fr->nbfp;
478     vdwtype          = mdatoms->typeA;
479
480     vftab            = kernel_data->table_vdw->data;
481     vftabscale       = kernel_data->table_vdw->scale;
482
483     /* Setup water-specific parameters */
484     inr              = nlist->iinr[0];
485     iq1              = facel*charge[inr+1];
486     iq2              = facel*charge[inr+2];
487     iq3              = facel*charge[inr+3];
488     vdwioffset0      = 2*nvdwtype*vdwtype[inr+0];
489
490     outeriter        = 0;
491     inneriter        = 0;
492
493     /* Start outer loop over neighborlists */
494     for(iidx=0; iidx<nri; iidx++)
495     {
496         /* Load shift vector for this list */
497         i_shift_offset   = DIM*shiftidx[iidx];
498         shX              = shiftvec[i_shift_offset+XX];
499         shY              = shiftvec[i_shift_offset+YY];
500         shZ              = shiftvec[i_shift_offset+ZZ];
501
502         /* Load limits for loop over neighbors */
503         j_index_start    = jindex[iidx];
504         j_index_end      = jindex[iidx+1];
505
506         /* Get outer coordinate index */
507         inr              = iinr[iidx];
508         i_coord_offset   = DIM*inr;
509
510         /* Load i particle coords and add shift vector */
511         ix0              = shX + x[i_coord_offset+DIM*0+XX];
512         iy0              = shY + x[i_coord_offset+DIM*0+YY];
513         iz0              = shZ + x[i_coord_offset+DIM*0+ZZ];
514         ix1              = shX + x[i_coord_offset+DIM*1+XX];
515         iy1              = shY + x[i_coord_offset+DIM*1+YY];
516         iz1              = shZ + x[i_coord_offset+DIM*1+ZZ];
517         ix2              = shX + x[i_coord_offset+DIM*2+XX];
518         iy2              = shY + x[i_coord_offset+DIM*2+YY];
519         iz2              = shZ + x[i_coord_offset+DIM*2+ZZ];
520         ix3              = shX + x[i_coord_offset+DIM*3+XX];
521         iy3              = shY + x[i_coord_offset+DIM*3+YY];
522         iz3              = shZ + x[i_coord_offset+DIM*3+ZZ];
523
524         fix0             = 0.0;
525         fiy0             = 0.0;
526         fiz0             = 0.0;
527         fix1             = 0.0;
528         fiy1             = 0.0;
529         fiz1             = 0.0;
530         fix2             = 0.0;
531         fiy2             = 0.0;
532         fiz2             = 0.0;
533         fix3             = 0.0;
534         fiy3             = 0.0;
535         fiz3             = 0.0;
536
537         /* Start inner kernel loop */
538         for(jidx=j_index_start; jidx<j_index_end; jidx++)
539         {
540             /* Get j neighbor index, and coordinate index */
541             jnr              = jjnr[jidx];
542             j_coord_offset   = DIM*jnr;
543
544             /* load j atom coordinates */
545             jx0              = x[j_coord_offset+DIM*0+XX];
546             jy0              = x[j_coord_offset+DIM*0+YY];
547             jz0              = x[j_coord_offset+DIM*0+ZZ];
548
549             /* Calculate displacement vector */
550             dx00             = ix0 - jx0;
551             dy00             = iy0 - jy0;
552             dz00             = iz0 - jz0;
553             dx10             = ix1 - jx0;
554             dy10             = iy1 - jy0;
555             dz10             = iz1 - jz0;
556             dx20             = ix2 - jx0;
557             dy20             = iy2 - jy0;
558             dz20             = iz2 - jz0;
559             dx30             = ix3 - jx0;
560             dy30             = iy3 - jy0;
561             dz30             = iz3 - jz0;
562
563             /* Calculate squared distance and things based on it */
564             rsq00            = dx00*dx00+dy00*dy00+dz00*dz00;
565             rsq10            = dx10*dx10+dy10*dy10+dz10*dz10;
566             rsq20            = dx20*dx20+dy20*dy20+dz20*dz20;
567             rsq30            = dx30*dx30+dy30*dy30+dz30*dz30;
568
569             rinv00           = gmx_invsqrt(rsq00);
570             rinv10           = gmx_invsqrt(rsq10);
571             rinv20           = gmx_invsqrt(rsq20);
572             rinv30           = gmx_invsqrt(rsq30);
573
574             rinvsq10         = rinv10*rinv10;
575             rinvsq20         = rinv20*rinv20;
576             rinvsq30         = rinv30*rinv30;
577
578             /* Load parameters for j particles */
579             jq0              = charge[jnr+0];
580             vdwjidx0         = 2*vdwtype[jnr+0];
581
582             /**************************
583              * CALCULATE INTERACTIONS *
584              **************************/
585
586             r00              = rsq00*rinv00;
587
588             c6_00            = vdwparam[vdwioffset0+vdwjidx0];
589             c12_00           = vdwparam[vdwioffset0+vdwjidx0+1];
590
591             /* Calculate table index by multiplying r with table scale and truncate to integer */
592             rt               = r00*vftabscale;
593             vfitab           = rt;
594             vfeps            = rt-vfitab;
595             vfitab           = 2*4*vfitab;
596
597             /* CUBIC SPLINE TABLE DISPERSION */
598             vfitab          += 0;
599             F                = vftab[vfitab+1];
600             Geps             = vfeps*vftab[vfitab+2];
601             Heps2            = vfeps*vfeps*vftab[vfitab+3];
602             Fp               = F+Geps+Heps2;
603             FF               = Fp+Geps+2.0*Heps2;
604             fvdw6            = c6_00*FF;
605
606             /* CUBIC SPLINE TABLE REPULSION */
607             F                = vftab[vfitab+5];
608             Geps             = vfeps*vftab[vfitab+6];
609             Heps2            = vfeps*vfeps*vftab[vfitab+7];
610             Fp               = F+Geps+Heps2;
611             FF               = Fp+Geps+2.0*Heps2;
612             fvdw12           = c12_00*FF;
613             fvdw             = -(fvdw6+fvdw12)*vftabscale*rinv00;
614
615             fscal            = fvdw;
616
617             /* Calculate temporary vectorial force */
618             tx               = fscal*dx00;
619             ty               = fscal*dy00;
620             tz               = fscal*dz00;
621
622             /* Update vectorial force */
623             fix0            += tx;
624             fiy0            += ty;
625             fiz0            += tz;
626             f[j_coord_offset+DIM*0+XX] -= tx;
627             f[j_coord_offset+DIM*0+YY] -= ty;
628             f[j_coord_offset+DIM*0+ZZ] -= tz;
629
630             /**************************
631              * CALCULATE INTERACTIONS *
632              **************************/
633
634             qq10             = iq1*jq0;
635
636             /* REACTION-FIELD ELECTROSTATICS */
637             felec            = qq10*(rinv10*rinvsq10-krf2);
638
639             fscal            = felec;
640
641             /* Calculate temporary vectorial force */
642             tx               = fscal*dx10;
643             ty               = fscal*dy10;
644             tz               = fscal*dz10;
645
646             /* Update vectorial force */
647             fix1            += tx;
648             fiy1            += ty;
649             fiz1            += tz;
650             f[j_coord_offset+DIM*0+XX] -= tx;
651             f[j_coord_offset+DIM*0+YY] -= ty;
652             f[j_coord_offset+DIM*0+ZZ] -= tz;
653
654             /**************************
655              * CALCULATE INTERACTIONS *
656              **************************/
657
658             qq20             = iq2*jq0;
659
660             /* REACTION-FIELD ELECTROSTATICS */
661             felec            = qq20*(rinv20*rinvsq20-krf2);
662
663             fscal            = felec;
664
665             /* Calculate temporary vectorial force */
666             tx               = fscal*dx20;
667             ty               = fscal*dy20;
668             tz               = fscal*dz20;
669
670             /* Update vectorial force */
671             fix2            += tx;
672             fiy2            += ty;
673             fiz2            += tz;
674             f[j_coord_offset+DIM*0+XX] -= tx;
675             f[j_coord_offset+DIM*0+YY] -= ty;
676             f[j_coord_offset+DIM*0+ZZ] -= tz;
677
678             /**************************
679              * CALCULATE INTERACTIONS *
680              **************************/
681
682             qq30             = iq3*jq0;
683
684             /* REACTION-FIELD ELECTROSTATICS */
685             felec            = qq30*(rinv30*rinvsq30-krf2);
686
687             fscal            = felec;
688
689             /* Calculate temporary vectorial force */
690             tx               = fscal*dx30;
691             ty               = fscal*dy30;
692             tz               = fscal*dz30;
693
694             /* Update vectorial force */
695             fix3            += tx;
696             fiy3            += ty;
697             fiz3            += tz;
698             f[j_coord_offset+DIM*0+XX] -= tx;
699             f[j_coord_offset+DIM*0+YY] -= ty;
700             f[j_coord_offset+DIM*0+ZZ] -= tz;
701
702             /* Inner loop uses 128 flops */
703         }
704         /* End of innermost loop */
705
706         tx = ty = tz = 0;
707         f[i_coord_offset+DIM*0+XX] += fix0;
708         f[i_coord_offset+DIM*0+YY] += fiy0;
709         f[i_coord_offset+DIM*0+ZZ] += fiz0;
710         tx                         += fix0;
711         ty                         += fiy0;
712         tz                         += fiz0;
713         f[i_coord_offset+DIM*1+XX] += fix1;
714         f[i_coord_offset+DIM*1+YY] += fiy1;
715         f[i_coord_offset+DIM*1+ZZ] += fiz1;
716         tx                         += fix1;
717         ty                         += fiy1;
718         tz                         += fiz1;
719         f[i_coord_offset+DIM*2+XX] += fix2;
720         f[i_coord_offset+DIM*2+YY] += fiy2;
721         f[i_coord_offset+DIM*2+ZZ] += fiz2;
722         tx                         += fix2;
723         ty                         += fiy2;
724         tz                         += fiz2;
725         f[i_coord_offset+DIM*3+XX] += fix3;
726         f[i_coord_offset+DIM*3+YY] += fiy3;
727         f[i_coord_offset+DIM*3+ZZ] += fiz3;
728         tx                         += fix3;
729         ty                         += fiy3;
730         tz                         += fiz3;
731         fshift[i_shift_offset+XX]  += tx;
732         fshift[i_shift_offset+YY]  += ty;
733         fshift[i_shift_offset+ZZ]  += tz;
734
735         /* Increment number of inner iterations */
736         inneriter                  += j_index_end - j_index_start;
737
738         /* Outer loop uses 39 flops */
739     }
740
741     /* Increment number of outer iterations */
742     outeriter        += nri;
743
744     /* Update outer/inner flops */
745
746     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W4_F,outeriter*39 + inneriter*128);
747 }