Merge release-4-6 into master
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_c / nb_kernel_ElecGB_VdwLJ_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_ElecGB_VdwLJ_GeomP1P1_VF_c
35  * Electrostatics interaction: GeneralizedBorn
36  * VdW interaction:            LennardJones
37  * Geometry:                   Particle-Particle
38  * Calculate force/pot:        PotentialAndForce
39  */
40 void
41 nb_kernel_ElecGB_VdwLJ_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              gbitab;
64     real             vgb,fgb,vgbsum,dvdasum,gbscale,gbtabscale,isaprod,gbqqfactor,gbinvepsdiff,dvdaj,gbeps,dvdatmp;
65     real             *invsqrta,*dvda,*gbtab;
66     int              nvdwtype;
67     real             rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,br,vvdwexp,sh_vdw_invrcut6;
68     int              *vdwtype;
69     real             *vdwparam;
70     int              vfitab;
71     real             rt,vfeps,vftabscale,Y,F,Geps,Heps2,Fp,VV,FF;
72     real             *vftab;
73
74     x                = xx[0];
75     f                = ff[0];
76
77     nri              = nlist->nri;
78     iinr             = nlist->iinr;
79     jindex           = nlist->jindex;
80     jjnr             = nlist->jjnr;
81     shiftidx         = nlist->shift;
82     gid              = nlist->gid;
83     shiftvec         = fr->shift_vec[0];
84     fshift           = fr->fshift[0];
85     facel            = fr->epsfac;
86     charge           = mdatoms->chargeA;
87     nvdwtype         = fr->ntype;
88     vdwparam         = fr->nbfp;
89     vdwtype          = mdatoms->typeA;
90
91     invsqrta         = fr->invsqrta;
92     dvda             = fr->dvda;
93     gbtabscale       = fr->gbtab.scale;
94     gbtab            = fr->gbtab.data;
95     gbinvepsdiff     = (1.0/fr->epsilon_r) - (1.0/fr->gb_epsilon_solvent);
96
97     outeriter        = 0;
98     inneriter        = 0;
99
100     /* Start outer loop over neighborlists */
101     for(iidx=0; iidx<nri; iidx++)
102     {
103         /* Load shift vector for this list */
104         i_shift_offset   = DIM*shiftidx[iidx];
105         shX              = shiftvec[i_shift_offset+XX];
106         shY              = shiftvec[i_shift_offset+YY];
107         shZ              = shiftvec[i_shift_offset+ZZ];
108
109         /* Load limits for loop over neighbors */
110         j_index_start    = jindex[iidx];
111         j_index_end      = jindex[iidx+1];
112
113         /* Get outer coordinate index */
114         inr              = iinr[iidx];
115         i_coord_offset   = DIM*inr;
116
117         /* Load i particle coords and add shift vector */
118         ix0              = shX + x[i_coord_offset+DIM*0+XX];
119         iy0              = shY + x[i_coord_offset+DIM*0+YY];
120         iz0              = shZ + x[i_coord_offset+DIM*0+ZZ];
121
122         fix0             = 0.0;
123         fiy0             = 0.0;
124         fiz0             = 0.0;
125
126         /* Load parameters for i particles */
127         iq0              = facel*charge[inr+0];
128         isai0            = invsqrta[inr+0];
129         vdwioffset0      = 2*nvdwtype*vdwtype[inr+0];
130
131         /* Reset potential sums */
132         velecsum         = 0.0;
133         vgbsum           = 0.0;
134         vvdwsum          = 0.0;
135         dvdasum          = 0.0;
136
137         /* Start inner kernel loop */
138         for(jidx=j_index_start; jidx<j_index_end; jidx++)
139         {
140             /* Get j neighbor index, and coordinate index */
141             jnr              = jjnr[jidx];
142             j_coord_offset   = DIM*jnr;
143
144             /* load j atom coordinates */
145             jx0              = x[j_coord_offset+DIM*0+XX];
146             jy0              = x[j_coord_offset+DIM*0+YY];
147             jz0              = x[j_coord_offset+DIM*0+ZZ];
148
149             /* Calculate displacement vector */
150             dx00             = ix0 - jx0;
151             dy00             = iy0 - jy0;
152             dz00             = iz0 - jz0;
153
154             /* Calculate squared distance and things based on it */
155             rsq00            = dx00*dx00+dy00*dy00+dz00*dz00;
156
157             rinv00           = gmx_invsqrt(rsq00);
158
159             rinvsq00         = rinv00*rinv00;
160
161             /* Load parameters for j particles */
162             jq0              = charge[jnr+0];
163             isaj0           = invsqrta[jnr+0];
164             vdwjidx0         = 2*vdwtype[jnr+0];
165
166             /**************************
167              * CALCULATE INTERACTIONS *
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             /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
177             isaprod          = isai0*isaj0;
178             gbqqfactor       = isaprod*(-qq00)*gbinvepsdiff;
179             gbscale          = isaprod*gbtabscale;
180             dvdaj            = dvda[jnr+0];
181
182             /* Calculate generalized born table index - this is a separate table from the normal one,
183              * but we use the same procedure by multiplying r with scale and truncating to integer.
184              */
185             rt               = r00*gbscale;
186             gbitab           = rt;
187             gbeps            = rt-gbitab;
188             gbitab           = 4*gbitab;
189
190             Y                = gbtab[gbitab];
191             F                = gbtab[gbitab+1];
192             Geps             = gbeps*gbtab[gbitab+2];
193             Heps2            = gbeps*gbeps*gbtab[gbitab+3];
194             Fp               = F+Geps+Heps2;
195             VV               = Y+gbeps*Fp;
196             vgb              = gbqqfactor*VV;
197
198             FF               = Fp+Geps+2.0*Heps2;
199             fgb              = gbqqfactor*FF*gbscale;
200             dvdatmp          = -0.5*(vgb+fgb*r00);
201             dvdasum          = dvdasum + dvdatmp;
202             dvda[jnr]        = dvdaj+dvdatmp*isaj0*isaj0;
203             velec            = qq00*rinv00;
204             felec            = (velec*rinv00-fgb)*rinv00;
205
206             /* LENNARD-JONES DISPERSION/REPULSION */
207
208             rinvsix          = rinvsq00*rinvsq00*rinvsq00;
209             vvdw6            = c6_00*rinvsix;
210             vvdw12           = c12_00*rinvsix*rinvsix;
211             vvdw             = vvdw12*(1.0/12.0) - vvdw6*(1.0/6.0);
212             fvdw             = (vvdw12-vvdw6)*rinvsq00;
213
214             /* Update potential sums from outer loop */
215             velecsum        += velec;
216             vgbsum          += vgb;
217             vvdwsum         += vvdw;
218
219             fscal            = felec+fvdw;
220
221             /* Calculate temporary vectorial force */
222             tx               = fscal*dx00;
223             ty               = fscal*dy00;
224             tz               = fscal*dz00;
225
226             /* Update vectorial force */
227             fix0            += tx;
228             fiy0            += ty;
229             fiz0            += tz;
230             f[j_coord_offset+DIM*0+XX] -= tx;
231             f[j_coord_offset+DIM*0+YY] -= ty;
232             f[j_coord_offset+DIM*0+ZZ] -= tz;
233
234             /* Inner loop uses 71 flops */
235         }
236         /* End of innermost loop */
237
238         tx = ty = tz = 0;
239         f[i_coord_offset+DIM*0+XX] += fix0;
240         f[i_coord_offset+DIM*0+YY] += fiy0;
241         f[i_coord_offset+DIM*0+ZZ] += fiz0;
242         tx                         += fix0;
243         ty                         += fiy0;
244         tz                         += fiz0;
245         fshift[i_shift_offset+XX]  += tx;
246         fshift[i_shift_offset+YY]  += ty;
247         fshift[i_shift_offset+ZZ]  += tz;
248
249         ggid                        = gid[iidx];
250         /* Update potential energies */
251         kernel_data->energygrp_elec[ggid] += velecsum;
252         kernel_data->energygrp_polarization[ggid] += vgbsum;
253         kernel_data->energygrp_vdw[ggid] += vvdwsum;
254         dvda[inr]                   = dvda[inr] + dvdasum*isai0*isai0;
255
256         /* Increment number of inner iterations */
257         inneriter                  += j_index_end - j_index_start;
258
259         /* Outer loop uses 16 flops */
260     }
261
262     /* Increment number of outer iterations */
263     outeriter        += nri;
264
265     /* Update outer/inner flops */
266
267     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_VF,outeriter*16 + inneriter*71);
268 }
269 /*
270  * Gromacs nonbonded kernel:   nb_kernel_ElecGB_VdwLJ_GeomP1P1_F_c
271  * Electrostatics interaction: GeneralizedBorn
272  * VdW interaction:            LennardJones
273  * Geometry:                   Particle-Particle
274  * Calculate force/pot:        Force
275  */
276 void
277 nb_kernel_ElecGB_VdwLJ_GeomP1P1_F_c
278                     (t_nblist * gmx_restrict                nlist,
279                      rvec * gmx_restrict                    xx,
280                      rvec * gmx_restrict                    ff,
281                      t_forcerec * gmx_restrict              fr,
282                      t_mdatoms * gmx_restrict               mdatoms,
283                      nb_kernel_data_t * gmx_restrict        kernel_data,
284                      t_nrnb * gmx_restrict                  nrnb)
285 {
286     int              i_shift_offset,i_coord_offset,j_coord_offset;
287     int              j_index_start,j_index_end;
288     int              nri,inr,ggid,iidx,jidx,jnr,outeriter,inneriter;
289     real             shX,shY,shZ,tx,ty,tz,fscal,rcutoff,rcutoff2;
290     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
291     real             *shiftvec,*fshift,*x,*f;
292     int              vdwioffset0;
293     real             ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
294     int              vdwjidx0;
295     real             jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
296     real             dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00,cexp1_00,cexp2_00;
297     real             velec,felec,velecsum,facel,crf,krf,krf2;
298     real             *charge;
299     int              gbitab;
300     real             vgb,fgb,vgbsum,dvdasum,gbscale,gbtabscale,isaprod,gbqqfactor,gbinvepsdiff,dvdaj,gbeps,dvdatmp;
301     real             *invsqrta,*dvda,*gbtab;
302     int              nvdwtype;
303     real             rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,br,vvdwexp,sh_vdw_invrcut6;
304     int              *vdwtype;
305     real             *vdwparam;
306     int              vfitab;
307     real             rt,vfeps,vftabscale,Y,F,Geps,Heps2,Fp,VV,FF;
308     real             *vftab;
309
310     x                = xx[0];
311     f                = ff[0];
312
313     nri              = nlist->nri;
314     iinr             = nlist->iinr;
315     jindex           = nlist->jindex;
316     jjnr             = nlist->jjnr;
317     shiftidx         = nlist->shift;
318     gid              = nlist->gid;
319     shiftvec         = fr->shift_vec[0];
320     fshift           = fr->fshift[0];
321     facel            = fr->epsfac;
322     charge           = mdatoms->chargeA;
323     nvdwtype         = fr->ntype;
324     vdwparam         = fr->nbfp;
325     vdwtype          = mdatoms->typeA;
326
327     invsqrta         = fr->invsqrta;
328     dvda             = fr->dvda;
329     gbtabscale       = fr->gbtab.scale;
330     gbtab            = fr->gbtab.data;
331     gbinvepsdiff     = (1.0/fr->epsilon_r) - (1.0/fr->gb_epsilon_solvent);
332
333     outeriter        = 0;
334     inneriter        = 0;
335
336     /* Start outer loop over neighborlists */
337     for(iidx=0; iidx<nri; iidx++)
338     {
339         /* Load shift vector for this list */
340         i_shift_offset   = DIM*shiftidx[iidx];
341         shX              = shiftvec[i_shift_offset+XX];
342         shY              = shiftvec[i_shift_offset+YY];
343         shZ              = shiftvec[i_shift_offset+ZZ];
344
345         /* Load limits for loop over neighbors */
346         j_index_start    = jindex[iidx];
347         j_index_end      = jindex[iidx+1];
348
349         /* Get outer coordinate index */
350         inr              = iinr[iidx];
351         i_coord_offset   = DIM*inr;
352
353         /* Load i particle coords and add shift vector */
354         ix0              = shX + x[i_coord_offset+DIM*0+XX];
355         iy0              = shY + x[i_coord_offset+DIM*0+YY];
356         iz0              = shZ + x[i_coord_offset+DIM*0+ZZ];
357
358         fix0             = 0.0;
359         fiy0             = 0.0;
360         fiz0             = 0.0;
361
362         /* Load parameters for i particles */
363         iq0              = facel*charge[inr+0];
364         isai0            = invsqrta[inr+0];
365         vdwioffset0      = 2*nvdwtype*vdwtype[inr+0];
366
367         dvdasum          = 0.0;
368
369         /* Start inner kernel loop */
370         for(jidx=j_index_start; jidx<j_index_end; jidx++)
371         {
372             /* Get j neighbor index, and coordinate index */
373             jnr              = jjnr[jidx];
374             j_coord_offset   = DIM*jnr;
375
376             /* load j atom coordinates */
377             jx0              = x[j_coord_offset+DIM*0+XX];
378             jy0              = x[j_coord_offset+DIM*0+YY];
379             jz0              = x[j_coord_offset+DIM*0+ZZ];
380
381             /* Calculate displacement vector */
382             dx00             = ix0 - jx0;
383             dy00             = iy0 - jy0;
384             dz00             = iz0 - jz0;
385
386             /* Calculate squared distance and things based on it */
387             rsq00            = dx00*dx00+dy00*dy00+dz00*dz00;
388
389             rinv00           = gmx_invsqrt(rsq00);
390
391             rinvsq00         = rinv00*rinv00;
392
393             /* Load parameters for j particles */
394             jq0              = charge[jnr+0];
395             isaj0           = invsqrta[jnr+0];
396             vdwjidx0         = 2*vdwtype[jnr+0];
397
398             /**************************
399              * CALCULATE INTERACTIONS *
400              **************************/
401
402             r00              = rsq00*rinv00;
403
404             qq00             = iq0*jq0;
405             c6_00            = vdwparam[vdwioffset0+vdwjidx0];
406             c12_00           = vdwparam[vdwioffset0+vdwjidx0+1];
407
408             /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
409             isaprod          = isai0*isaj0;
410             gbqqfactor       = isaprod*(-qq00)*gbinvepsdiff;
411             gbscale          = isaprod*gbtabscale;
412             dvdaj            = dvda[jnr+0];
413
414             /* Calculate generalized born table index - this is a separate table from the normal one,
415              * but we use the same procedure by multiplying r with scale and truncating to integer.
416              */
417             rt               = r00*gbscale;
418             gbitab           = rt;
419             gbeps            = rt-gbitab;
420             gbitab           = 4*gbitab;
421
422             Y                = gbtab[gbitab];
423             F                = gbtab[gbitab+1];
424             Geps             = gbeps*gbtab[gbitab+2];
425             Heps2            = gbeps*gbeps*gbtab[gbitab+3];
426             Fp               = F+Geps+Heps2;
427             VV               = Y+gbeps*Fp;
428             vgb              = gbqqfactor*VV;
429
430             FF               = Fp+Geps+2.0*Heps2;
431             fgb              = gbqqfactor*FF*gbscale;
432             dvdatmp          = -0.5*(vgb+fgb*r00);
433             dvdasum          = dvdasum + dvdatmp;
434             dvda[jnr]        = dvdaj+dvdatmp*isaj0*isaj0;
435             velec            = qq00*rinv00;
436             felec            = (velec*rinv00-fgb)*rinv00;
437
438             /* LENNARD-JONES DISPERSION/REPULSION */
439
440             rinvsix          = rinvsq00*rinvsq00*rinvsq00;
441             fvdw             = (c12_00*rinvsix-c6_00)*rinvsix*rinvsq00;
442
443             fscal            = felec+fvdw;
444
445             /* Calculate temporary vectorial force */
446             tx               = fscal*dx00;
447             ty               = fscal*dy00;
448             tz               = fscal*dz00;
449
450             /* Update vectorial force */
451             fix0            += tx;
452             fiy0            += ty;
453             fiz0            += tz;
454             f[j_coord_offset+DIM*0+XX] -= tx;
455             f[j_coord_offset+DIM*0+YY] -= ty;
456             f[j_coord_offset+DIM*0+ZZ] -= tz;
457
458             /* Inner loop uses 64 flops */
459         }
460         /* End of innermost loop */
461
462         tx = ty = tz = 0;
463         f[i_coord_offset+DIM*0+XX] += fix0;
464         f[i_coord_offset+DIM*0+YY] += fiy0;
465         f[i_coord_offset+DIM*0+ZZ] += fiz0;
466         tx                         += fix0;
467         ty                         += fiy0;
468         tz                         += fiz0;
469         fshift[i_shift_offset+XX]  += tx;
470         fshift[i_shift_offset+YY]  += ty;
471         fshift[i_shift_offset+ZZ]  += tz;
472
473         dvda[inr]                   = dvda[inr] + dvdasum*isai0*isai0;
474
475         /* Increment number of inner iterations */
476         inneriter                  += j_index_end - j_index_start;
477
478         /* Outer loop uses 13 flops */
479     }
480
481     /* Increment number of outer iterations */
482     outeriter        += nri;
483
484     /* Update outer/inner flops */
485
486     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_F,outeriter*13 + inneriter*64);
487 }