93b7f50ab056797153b09ccf9361e87fc4956bbe
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_c / nb_kernel_ElecRFCut_VdwCSTab_GeomP1P1_c.c
1 /*
2  * Note: this file was generated by the Gromacs c kernel generator.
3  *
4  *                This source code is part of
5  *
6  *                 G   R   O   M   A   C   S
7  *
8  * Copyright (c) 2001-2012, The GROMACS Development Team
9  *
10  * Gromacs is a library for molecular simulation and trajectory analysis,
11  * written by Erik Lindahl, David van der Spoel, Berk Hess, and others - for
12  * a full list of developers and information, check out http://www.gromacs.org
13  *
14  * This program is free software; you can redistribute it and/or modify it under
15  * the terms of the GNU Lesser General Public License as published by the Free
16  * Software Foundation; either version 2 of the License, or (at your option) any
17  * later version.
18  *
19  * To help fund GROMACS development, we humbly ask that you cite
20  * the papers people have written on it - you can find them on the website.
21  */
22 #ifdef HAVE_CONFIG_H
23 #include <config.h>
24 #endif
25
26 #include <math.h>
27
28 #include "../nb_kernel.h"
29 #include "types/simple.h"
30 #include "vec.h"
31 #include "nrnb.h"
32
33 /*
34  * Gromacs nonbonded kernel:   nb_kernel_ElecRFCut_VdwCSTab_GeomP1P1_VF_c
35  * Electrostatics interaction: ReactionField
36  * VdW interaction:            CubicSplineTable
37  * Geometry:                   Particle-Particle
38  * Calculate force/pot:        PotentialAndForce
39  */
40 void
41 nb_kernel_ElecRFCut_VdwCSTab_GeomP1P1_VF_c
42                     (t_nblist * gmx_restrict                nlist,
43                      rvec * gmx_restrict                    xx,
44                      rvec * gmx_restrict                    ff,
45                      t_forcerec * gmx_restrict              fr,
46                      t_mdatoms * gmx_restrict               mdatoms,
47                      nb_kernel_data_t * gmx_restrict        kernel_data,
48                      t_nrnb * gmx_restrict                  nrnb)
49 {
50     int              i_shift_offset,i_coord_offset,j_coord_offset;
51     int              j_index_start,j_index_end;
52     int              nri,inr,ggid,iidx,jidx,jnr,outeriter,inneriter;
53     real             shX,shY,shZ,tx,ty,tz,fscal,rcutoff,rcutoff2;
54     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
55     real             *shiftvec,*fshift,*x,*f;
56     int              vdwioffset0;
57     real             ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
58     int              vdwjidx0;
59     real             jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
60     real             dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00,cexp1_00,cexp2_00;
61     real             velec,felec,velecsum,facel,crf,krf,krf2;
62     real             *charge;
63     int              nvdwtype;
64     real             rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,br,vvdwexp,sh_vdw_invrcut6;
65     int              *vdwtype;
66     real             *vdwparam;
67     int              vfitab;
68     real             rt,vfeps,vftabscale,Y,F,Geps,Heps2,Fp,VV,FF;
69     real             *vftab;
70
71     x                = xx[0];
72     f                = ff[0];
73
74     nri              = nlist->nri;
75     iinr             = nlist->iinr;
76     jindex           = nlist->jindex;
77     jjnr             = nlist->jjnr;
78     shiftidx         = nlist->shift;
79     gid              = nlist->gid;
80     shiftvec         = fr->shift_vec[0];
81     fshift           = fr->fshift[0];
82     facel            = fr->epsfac;
83     charge           = mdatoms->chargeA;
84     krf              = fr->ic->k_rf;
85     krf2             = krf*2.0;
86     crf              = fr->ic->c_rf;
87     nvdwtype         = fr->ntype;
88     vdwparam         = fr->nbfp;
89     vdwtype          = mdatoms->typeA;
90
91     vftab            = kernel_data->table_vdw->data;
92     vftabscale       = kernel_data->table_vdw->scale;
93
94     /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
95     rcutoff          = fr->rcoulomb;
96     rcutoff2         = rcutoff*rcutoff;
97
98     outeriter        = 0;
99     inneriter        = 0;
100
101     /* Start outer loop over neighborlists */
102     for(iidx=0; iidx<nri; iidx++)
103     {
104         /* Load shift vector for this list */
105         i_shift_offset   = DIM*shiftidx[iidx];
106         shX              = shiftvec[i_shift_offset+XX];
107         shY              = shiftvec[i_shift_offset+YY];
108         shZ              = shiftvec[i_shift_offset+ZZ];
109
110         /* Load limits for loop over neighbors */
111         j_index_start    = jindex[iidx];
112         j_index_end      = jindex[iidx+1];
113
114         /* Get outer coordinate index */
115         inr              = iinr[iidx];
116         i_coord_offset   = DIM*inr;
117
118         /* Load i particle coords and add shift vector */
119         ix0              = shX + x[i_coord_offset+DIM*0+XX];
120         iy0              = shY + x[i_coord_offset+DIM*0+YY];
121         iz0              = shZ + x[i_coord_offset+DIM*0+ZZ];
122
123         fix0             = 0.0;
124         fiy0             = 0.0;
125         fiz0             = 0.0;
126
127         /* Load parameters for i particles */
128         iq0              = facel*charge[inr+0];
129         vdwioffset0      = 2*nvdwtype*vdwtype[inr+0];
130
131         /* Reset potential sums */
132         velecsum         = 0.0;
133         vvdwsum          = 0.0;
134
135         /* Start inner kernel loop */
136         for(jidx=j_index_start; jidx<j_index_end; jidx++)
137         {
138             /* Get j neighbor index, and coordinate index */
139             jnr              = jjnr[jidx];
140             j_coord_offset   = DIM*jnr;
141
142             /* load j atom coordinates */
143             jx0              = x[j_coord_offset+DIM*0+XX];
144             jy0              = x[j_coord_offset+DIM*0+YY];
145             jz0              = x[j_coord_offset+DIM*0+ZZ];
146
147             /* Calculate displacement vector */
148             dx00             = ix0 - jx0;
149             dy00             = iy0 - jy0;
150             dz00             = iz0 - jz0;
151
152             /* Calculate squared distance and things based on it */
153             rsq00            = dx00*dx00+dy00*dy00+dz00*dz00;
154
155             rinv00           = gmx_invsqrt(rsq00);
156
157             rinvsq00         = rinv00*rinv00;
158
159             /* Load parameters for j particles */
160             jq0              = charge[jnr+0];
161             vdwjidx0         = 2*vdwtype[jnr+0];
162
163             /**************************
164              * CALCULATE INTERACTIONS *
165              **************************/
166
167             if (rsq00<rcutoff2)
168             {
169
170             r00              = rsq00*rinv00;
171
172             qq00             = iq0*jq0;
173             c6_00            = vdwparam[vdwioffset0+vdwjidx0];
174             c12_00           = vdwparam[vdwioffset0+vdwjidx0+1];
175
176             /* Calculate table index by multiplying r with table scale and truncate to integer */
177             rt               = r00*vftabscale;
178             vfitab           = rt;
179             vfeps            = rt-vfitab;
180             vfitab           = 2*4*vfitab;
181
182             /* REACTION-FIELD ELECTROSTATICS */
183             velec            = qq00*(rinv00+krf*rsq00-crf);
184             felec            = qq00*(rinv00*rinvsq00-krf2);
185
186             /* CUBIC SPLINE TABLE DISPERSION */
187             vfitab          += 0;
188             Y                = vftab[vfitab];
189             F                = vftab[vfitab+1];
190             Geps             = vfeps*vftab[vfitab+2];
191             Heps2            = vfeps*vfeps*vftab[vfitab+3];
192             Fp               = F+Geps+Heps2;
193             VV               = Y+vfeps*Fp;
194             vvdw6            = c6_00*VV;
195             FF               = Fp+Geps+2.0*Heps2;
196             fvdw6            = c6_00*FF;
197
198             /* CUBIC SPLINE TABLE REPULSION */
199             Y                = vftab[vfitab+4];
200             F                = vftab[vfitab+5];
201             Geps             = vfeps*vftab[vfitab+6];
202             Heps2            = vfeps*vfeps*vftab[vfitab+7];
203             Fp               = F+Geps+Heps2;
204             VV               = Y+vfeps*Fp;
205             vvdw12           = c12_00*VV;
206             FF               = Fp+Geps+2.0*Heps2;
207             fvdw12           = c12_00*FF;
208             vvdw             = vvdw12+vvdw6;
209             fvdw             = -(fvdw6+fvdw12)*vftabscale*rinv00;
210
211             /* Update potential sums from outer loop */
212             velecsum        += velec;
213             vvdwsum         += vvdw;
214
215             fscal            = felec+fvdw;
216
217             /* Calculate temporary vectorial force */
218             tx               = fscal*dx00;
219             ty               = fscal*dy00;
220             tz               = fscal*dz00;
221
222             /* Update vectorial force */
223             fix0            += tx;
224             fiy0            += ty;
225             fiz0            += tz;
226             f[j_coord_offset+DIM*0+XX] -= tx;
227             f[j_coord_offset+DIM*0+YY] -= ty;
228             f[j_coord_offset+DIM*0+ZZ] -= tz;
229
230             }
231
232             /* Inner loop uses 66 flops */
233         }
234         /* End of innermost loop */
235
236         tx = ty = tz = 0;
237         f[i_coord_offset+DIM*0+XX] += fix0;
238         f[i_coord_offset+DIM*0+YY] += fiy0;
239         f[i_coord_offset+DIM*0+ZZ] += fiz0;
240         tx                         += fix0;
241         ty                         += fiy0;
242         tz                         += fiz0;
243         fshift[i_shift_offset+XX]  += tx;
244         fshift[i_shift_offset+YY]  += ty;
245         fshift[i_shift_offset+ZZ]  += tz;
246
247         ggid                        = gid[iidx];
248         /* Update potential energies */
249         kernel_data->energygrp_elec[ggid] += velecsum;
250         kernel_data->energygrp_vdw[ggid] += vvdwsum;
251
252         /* Increment number of inner iterations */
253         inneriter                  += j_index_end - j_index_start;
254
255         /* Outer loop uses 15 flops */
256     }
257
258     /* Increment number of outer iterations */
259     outeriter        += nri;
260
261     /* Update outer/inner flops */
262
263     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_VF,outeriter*15 + inneriter*66);
264 }
265 /*
266  * Gromacs nonbonded kernel:   nb_kernel_ElecRFCut_VdwCSTab_GeomP1P1_F_c
267  * Electrostatics interaction: ReactionField
268  * VdW interaction:            CubicSplineTable
269  * Geometry:                   Particle-Particle
270  * Calculate force/pot:        Force
271  */
272 void
273 nb_kernel_ElecRFCut_VdwCSTab_GeomP1P1_F_c
274                     (t_nblist * gmx_restrict                nlist,
275                      rvec * gmx_restrict                    xx,
276                      rvec * gmx_restrict                    ff,
277                      t_forcerec * gmx_restrict              fr,
278                      t_mdatoms * gmx_restrict               mdatoms,
279                      nb_kernel_data_t * gmx_restrict        kernel_data,
280                      t_nrnb * gmx_restrict                  nrnb)
281 {
282     int              i_shift_offset,i_coord_offset,j_coord_offset;
283     int              j_index_start,j_index_end;
284     int              nri,inr,ggid,iidx,jidx,jnr,outeriter,inneriter;
285     real             shX,shY,shZ,tx,ty,tz,fscal,rcutoff,rcutoff2;
286     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
287     real             *shiftvec,*fshift,*x,*f;
288     int              vdwioffset0;
289     real             ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
290     int              vdwjidx0;
291     real             jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
292     real             dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00,cexp1_00,cexp2_00;
293     real             velec,felec,velecsum,facel,crf,krf,krf2;
294     real             *charge;
295     int              nvdwtype;
296     real             rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,br,vvdwexp,sh_vdw_invrcut6;
297     int              *vdwtype;
298     real             *vdwparam;
299     int              vfitab;
300     real             rt,vfeps,vftabscale,Y,F,Geps,Heps2,Fp,VV,FF;
301     real             *vftab;
302
303     x                = xx[0];
304     f                = ff[0];
305
306     nri              = nlist->nri;
307     iinr             = nlist->iinr;
308     jindex           = nlist->jindex;
309     jjnr             = nlist->jjnr;
310     shiftidx         = nlist->shift;
311     gid              = nlist->gid;
312     shiftvec         = fr->shift_vec[0];
313     fshift           = fr->fshift[0];
314     facel            = fr->epsfac;
315     charge           = mdatoms->chargeA;
316     krf              = fr->ic->k_rf;
317     krf2             = krf*2.0;
318     crf              = fr->ic->c_rf;
319     nvdwtype         = fr->ntype;
320     vdwparam         = fr->nbfp;
321     vdwtype          = mdatoms->typeA;
322
323     vftab            = kernel_data->table_vdw->data;
324     vftabscale       = kernel_data->table_vdw->scale;
325
326     /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
327     rcutoff          = fr->rcoulomb;
328     rcutoff2         = rcutoff*rcutoff;
329
330     outeriter        = 0;
331     inneriter        = 0;
332
333     /* Start outer loop over neighborlists */
334     for(iidx=0; iidx<nri; iidx++)
335     {
336         /* Load shift vector for this list */
337         i_shift_offset   = DIM*shiftidx[iidx];
338         shX              = shiftvec[i_shift_offset+XX];
339         shY              = shiftvec[i_shift_offset+YY];
340         shZ              = shiftvec[i_shift_offset+ZZ];
341
342         /* Load limits for loop over neighbors */
343         j_index_start    = jindex[iidx];
344         j_index_end      = jindex[iidx+1];
345
346         /* Get outer coordinate index */
347         inr              = iinr[iidx];
348         i_coord_offset   = DIM*inr;
349
350         /* Load i particle coords and add shift vector */
351         ix0              = shX + x[i_coord_offset+DIM*0+XX];
352         iy0              = shY + x[i_coord_offset+DIM*0+YY];
353         iz0              = shZ + x[i_coord_offset+DIM*0+ZZ];
354
355         fix0             = 0.0;
356         fiy0             = 0.0;
357         fiz0             = 0.0;
358
359         /* Load parameters for i particles */
360         iq0              = facel*charge[inr+0];
361         vdwioffset0      = 2*nvdwtype*vdwtype[inr+0];
362
363         /* Start inner kernel loop */
364         for(jidx=j_index_start; jidx<j_index_end; jidx++)
365         {
366             /* Get j neighbor index, and coordinate index */
367             jnr              = jjnr[jidx];
368             j_coord_offset   = DIM*jnr;
369
370             /* load j atom coordinates */
371             jx0              = x[j_coord_offset+DIM*0+XX];
372             jy0              = x[j_coord_offset+DIM*0+YY];
373             jz0              = x[j_coord_offset+DIM*0+ZZ];
374
375             /* Calculate displacement vector */
376             dx00             = ix0 - jx0;
377             dy00             = iy0 - jy0;
378             dz00             = iz0 - jz0;
379
380             /* Calculate squared distance and things based on it */
381             rsq00            = dx00*dx00+dy00*dy00+dz00*dz00;
382
383             rinv00           = gmx_invsqrt(rsq00);
384
385             rinvsq00         = rinv00*rinv00;
386
387             /* Load parameters for j particles */
388             jq0              = charge[jnr+0];
389             vdwjidx0         = 2*vdwtype[jnr+0];
390
391             /**************************
392              * CALCULATE INTERACTIONS *
393              **************************/
394
395             if (rsq00<rcutoff2)
396             {
397
398             r00              = rsq00*rinv00;
399
400             qq00             = iq0*jq0;
401             c6_00            = vdwparam[vdwioffset0+vdwjidx0];
402             c12_00           = vdwparam[vdwioffset0+vdwjidx0+1];
403
404             /* Calculate table index by multiplying r with table scale and truncate to integer */
405             rt               = r00*vftabscale;
406             vfitab           = rt;
407             vfeps            = rt-vfitab;
408             vfitab           = 2*4*vfitab;
409
410             /* REACTION-FIELD ELECTROSTATICS */
411             felec            = qq00*(rinv00*rinvsq00-krf2);
412
413             /* CUBIC SPLINE TABLE DISPERSION */
414             vfitab          += 0;
415             F                = vftab[vfitab+1];
416             Geps             = vfeps*vftab[vfitab+2];
417             Heps2            = vfeps*vfeps*vftab[vfitab+3];
418             Fp               = F+Geps+Heps2;
419             FF               = Fp+Geps+2.0*Heps2;
420             fvdw6            = c6_00*FF;
421
422             /* CUBIC SPLINE TABLE REPULSION */
423             F                = vftab[vfitab+5];
424             Geps             = vfeps*vftab[vfitab+6];
425             Heps2            = vfeps*vfeps*vftab[vfitab+7];
426             Fp               = F+Geps+Heps2;
427             FF               = Fp+Geps+2.0*Heps2;
428             fvdw12           = c12_00*FF;
429             fvdw             = -(fvdw6+fvdw12)*vftabscale*rinv00;
430
431             fscal            = felec+fvdw;
432
433             /* Calculate temporary vectorial force */
434             tx               = fscal*dx00;
435             ty               = fscal*dy00;
436             tz               = fscal*dz00;
437
438             /* Update vectorial force */
439             fix0            += tx;
440             fiy0            += ty;
441             fiz0            += tz;
442             f[j_coord_offset+DIM*0+XX] -= tx;
443             f[j_coord_offset+DIM*0+YY] -= ty;
444             f[j_coord_offset+DIM*0+ZZ] -= tz;
445
446             }
447
448             /* Inner loop uses 53 flops */
449         }
450         /* End of innermost loop */
451
452         tx = ty = tz = 0;
453         f[i_coord_offset+DIM*0+XX] += fix0;
454         f[i_coord_offset+DIM*0+YY] += fiy0;
455         f[i_coord_offset+DIM*0+ZZ] += fiz0;
456         tx                         += fix0;
457         ty                         += fiy0;
458         tz                         += fiz0;
459         fshift[i_shift_offset+XX]  += tx;
460         fshift[i_shift_offset+YY]  += ty;
461         fshift[i_shift_offset+ZZ]  += tz;
462
463         /* Increment number of inner iterations */
464         inneriter                  += j_index_end - j_index_start;
465
466         /* Outer loop uses 13 flops */
467     }
468
469     /* Increment number of outer iterations */
470     outeriter        += nri;
471
472     /* Update outer/inner flops */
473
474     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_F,outeriter*13 + inneriter*53);
475 }