2f6a6e641cd2a4c2eff2c04f08b933738bb50830
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_c / nb_kernel_ElecNone_VdwLJEwSh_GeomP1P1_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 "types/simple.h"
44 #include "gromacs/math/vec.h"
45 #include "nrnb.h"
46
47 /*
48  * Gromacs nonbonded kernel:   nb_kernel_ElecNone_VdwLJEwSh_GeomP1P1_VF_c
49  * Electrostatics interaction: None
50  * VdW interaction:            LJEwald
51  * Geometry:                   Particle-Particle
52  * Calculate force/pot:        PotentialAndForce
53  */
54 void
55 nb_kernel_ElecNone_VdwLJEwSh_GeomP1P1_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              vdwjidx0;
73     real             jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
74     real             dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00,cexp1_00,cexp2_00;
75     int              nvdwtype;
76     real             rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,br,vvdwexp,sh_vdw_invrcut6;
77     int              *vdwtype;
78     real             *vdwparam;
79     real             c6grid_00;
80     real             ewclj,ewclj2,ewclj6,ewcljrsq,poly,exponent,sh_lj_ewald;
81     real             *vdwgridparam;
82
83     x                = xx[0];
84     f                = ff[0];
85
86     nri              = nlist->nri;
87     iinr             = nlist->iinr;
88     jindex           = nlist->jindex;
89     jjnr             = nlist->jjnr;
90     shiftidx         = nlist->shift;
91     gid              = nlist->gid;
92     shiftvec         = fr->shift_vec[0];
93     fshift           = fr->fshift[0];
94     nvdwtype         = fr->ntype;
95     vdwparam         = fr->nbfp;
96     vdwtype          = mdatoms->typeA;
97     vdwgridparam     = fr->ljpme_c6grid;
98     ewclj            = fr->ewaldcoeff_lj;
99     sh_lj_ewald      = fr->ic->sh_lj_ewald;
100     ewclj2           = ewclj*ewclj;
101     ewclj6           = ewclj2*ewclj2*ewclj2;
102
103     rcutoff          = fr->rvdw;
104     rcutoff2         = rcutoff*rcutoff;
105
106     sh_vdw_invrcut6  = fr->ic->sh_invrc6;
107     rvdw             = fr->rvdw;
108
109     outeriter        = 0;
110     inneriter        = 0;
111
112     /* Start outer loop over neighborlists */
113     for(iidx=0; iidx<nri; iidx++)
114     {
115         /* Load shift vector for this list */
116         i_shift_offset   = DIM*shiftidx[iidx];
117         shX              = shiftvec[i_shift_offset+XX];
118         shY              = shiftvec[i_shift_offset+YY];
119         shZ              = shiftvec[i_shift_offset+ZZ];
120
121         /* Load limits for loop over neighbors */
122         j_index_start    = jindex[iidx];
123         j_index_end      = jindex[iidx+1];
124
125         /* Get outer coordinate index */
126         inr              = iinr[iidx];
127         i_coord_offset   = DIM*inr;
128
129         /* Load i particle coords and add shift vector */
130         ix0              = shX + x[i_coord_offset+DIM*0+XX];
131         iy0              = shY + x[i_coord_offset+DIM*0+YY];
132         iz0              = shZ + x[i_coord_offset+DIM*0+ZZ];
133
134         fix0             = 0.0;
135         fiy0             = 0.0;
136         fiz0             = 0.0;
137
138         /* Load parameters for i particles */
139         vdwioffset0      = 2*nvdwtype*vdwtype[inr+0];
140
141         /* Reset potential sums */
142         vvdwsum          = 0.0;
143
144         /* Start inner kernel loop */
145         for(jidx=j_index_start; jidx<j_index_end; jidx++)
146         {
147             /* Get j neighbor index, and coordinate index */
148             jnr              = jjnr[jidx];
149             j_coord_offset   = DIM*jnr;
150
151             /* load j atom coordinates */
152             jx0              = x[j_coord_offset+DIM*0+XX];
153             jy0              = x[j_coord_offset+DIM*0+YY];
154             jz0              = x[j_coord_offset+DIM*0+ZZ];
155
156             /* Calculate displacement vector */
157             dx00             = ix0 - jx0;
158             dy00             = iy0 - jy0;
159             dz00             = iz0 - jz0;
160
161             /* Calculate squared distance and things based on it */
162             rsq00            = dx00*dx00+dy00*dy00+dz00*dz00;
163
164             rinv00           = gmx_invsqrt(rsq00);
165
166             rinvsq00         = rinv00*rinv00;
167
168             /* Load parameters for j particles */
169             vdwjidx0         = 2*vdwtype[jnr+0];
170
171             /**************************
172              * CALCULATE INTERACTIONS *
173              **************************/
174
175             if (rsq00<rcutoff2)
176             {
177
178             r00              = rsq00*rinv00;
179
180             c6_00            = vdwparam[vdwioffset0+vdwjidx0];
181             c12_00           = vdwparam[vdwioffset0+vdwjidx0+1];
182             c6grid_00        = vdwgridparam[vdwioffset0+vdwjidx0];
183
184             rinvsix          = rinvsq00*rinvsq00*rinvsq00;
185             ewcljrsq         = ewclj2*rsq00;
186             exponent         = exp(-ewcljrsq);
187             poly             = exponent*(1.0 + ewcljrsq + ewcljrsq*ewcljrsq*0.5);
188             vvdw6            = (c6_00-c6grid_00*(1.0-poly))*rinvsix;
189             vvdw12           = c12_00*rinvsix*rinvsix;
190             vvdw             = (vvdw12 - c12_00*sh_vdw_invrcut6*sh_vdw_invrcut6)*(1.0/12.0) - (vvdw6 - c6_00*sh_vdw_invrcut6 - c6grid_00*sh_lj_ewald)*(1.0/6.0);
191             fvdw             = (vvdw12 - vvdw6 - c6grid_00*(1.0/6.0)*exponent*ewclj6)*rinvsq00;
192
193             /* Update potential sums from outer loop */
194             vvdwsum         += vvdw;
195
196             fscal            = fvdw;
197
198             /* Calculate temporary vectorial force */
199             tx               = fscal*dx00;
200             ty               = fscal*dy00;
201             tz               = fscal*dz00;
202
203             /* Update vectorial force */
204             fix0            += tx;
205             fiy0            += ty;
206             fiz0            += tz;
207             f[j_coord_offset+DIM*0+XX] -= tx;
208             f[j_coord_offset+DIM*0+YY] -= ty;
209             f[j_coord_offset+DIM*0+ZZ] -= tz;
210
211             }
212
213             /* Inner loop uses 55 flops */
214         }
215         /* End of innermost loop */
216
217         tx = ty = tz = 0;
218         f[i_coord_offset+DIM*0+XX] += fix0;
219         f[i_coord_offset+DIM*0+YY] += fiy0;
220         f[i_coord_offset+DIM*0+ZZ] += fiz0;
221         tx                         += fix0;
222         ty                         += fiy0;
223         tz                         += fiz0;
224         fshift[i_shift_offset+XX]  += tx;
225         fshift[i_shift_offset+YY]  += ty;
226         fshift[i_shift_offset+ZZ]  += tz;
227
228         ggid                        = gid[iidx];
229         /* Update potential energies */
230         kernel_data->energygrp_vdw[ggid] += vvdwsum;
231
232         /* Increment number of inner iterations */
233         inneriter                  += j_index_end - j_index_start;
234
235         /* Outer loop uses 13 flops */
236     }
237
238     /* Increment number of outer iterations */
239     outeriter        += nri;
240
241     /* Update outer/inner flops */
242
243     inc_nrnb(nrnb,eNR_NBKERNEL_VDW_VF,outeriter*13 + inneriter*55);
244 }
245 /*
246  * Gromacs nonbonded kernel:   nb_kernel_ElecNone_VdwLJEwSh_GeomP1P1_F_c
247  * Electrostatics interaction: None
248  * VdW interaction:            LJEwald
249  * Geometry:                   Particle-Particle
250  * Calculate force/pot:        Force
251  */
252 void
253 nb_kernel_ElecNone_VdwLJEwSh_GeomP1P1_F_c
254                     (t_nblist                    * gmx_restrict       nlist,
255                      rvec                        * gmx_restrict          xx,
256                      rvec                        * gmx_restrict          ff,
257                      t_forcerec                  * gmx_restrict          fr,
258                      t_mdatoms                   * gmx_restrict     mdatoms,
259                      nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
260                      t_nrnb                      * gmx_restrict        nrnb)
261 {
262     int              i_shift_offset,i_coord_offset,j_coord_offset;
263     int              j_index_start,j_index_end;
264     int              nri,inr,ggid,iidx,jidx,jnr,outeriter,inneriter;
265     real             shX,shY,shZ,tx,ty,tz,fscal,rcutoff,rcutoff2;
266     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
267     real             *shiftvec,*fshift,*x,*f;
268     int              vdwioffset0;
269     real             ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
270     int              vdwjidx0;
271     real             jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
272     real             dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00,cexp1_00,cexp2_00;
273     int              nvdwtype;
274     real             rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,br,vvdwexp,sh_vdw_invrcut6;
275     int              *vdwtype;
276     real             *vdwparam;
277     real             c6grid_00;
278     real             ewclj,ewclj2,ewclj6,ewcljrsq,poly,exponent,sh_lj_ewald;
279     real             *vdwgridparam;
280
281     x                = xx[0];
282     f                = ff[0];
283
284     nri              = nlist->nri;
285     iinr             = nlist->iinr;
286     jindex           = nlist->jindex;
287     jjnr             = nlist->jjnr;
288     shiftidx         = nlist->shift;
289     gid              = nlist->gid;
290     shiftvec         = fr->shift_vec[0];
291     fshift           = fr->fshift[0];
292     nvdwtype         = fr->ntype;
293     vdwparam         = fr->nbfp;
294     vdwtype          = mdatoms->typeA;
295     vdwgridparam     = fr->ljpme_c6grid;
296     ewclj            = fr->ewaldcoeff_lj;
297     sh_lj_ewald      = fr->ic->sh_lj_ewald;
298     ewclj2           = ewclj*ewclj;
299     ewclj6           = ewclj2*ewclj2*ewclj2;
300
301     rcutoff          = fr->rvdw;
302     rcutoff2         = rcutoff*rcutoff;
303
304     sh_vdw_invrcut6  = fr->ic->sh_invrc6;
305     rvdw             = fr->rvdw;
306
307     outeriter        = 0;
308     inneriter        = 0;
309
310     /* Start outer loop over neighborlists */
311     for(iidx=0; iidx<nri; iidx++)
312     {
313         /* Load shift vector for this list */
314         i_shift_offset   = DIM*shiftidx[iidx];
315         shX              = shiftvec[i_shift_offset+XX];
316         shY              = shiftvec[i_shift_offset+YY];
317         shZ              = shiftvec[i_shift_offset+ZZ];
318
319         /* Load limits for loop over neighbors */
320         j_index_start    = jindex[iidx];
321         j_index_end      = jindex[iidx+1];
322
323         /* Get outer coordinate index */
324         inr              = iinr[iidx];
325         i_coord_offset   = DIM*inr;
326
327         /* Load i particle coords and add shift vector */
328         ix0              = shX + x[i_coord_offset+DIM*0+XX];
329         iy0              = shY + x[i_coord_offset+DIM*0+YY];
330         iz0              = shZ + x[i_coord_offset+DIM*0+ZZ];
331
332         fix0             = 0.0;
333         fiy0             = 0.0;
334         fiz0             = 0.0;
335
336         /* Load parameters for i particles */
337         vdwioffset0      = 2*nvdwtype*vdwtype[inr+0];
338
339         /* Start inner kernel loop */
340         for(jidx=j_index_start; jidx<j_index_end; jidx++)
341         {
342             /* Get j neighbor index, and coordinate index */
343             jnr              = jjnr[jidx];
344             j_coord_offset   = DIM*jnr;
345
346             /* load j atom coordinates */
347             jx0              = x[j_coord_offset+DIM*0+XX];
348             jy0              = x[j_coord_offset+DIM*0+YY];
349             jz0              = x[j_coord_offset+DIM*0+ZZ];
350
351             /* Calculate displacement vector */
352             dx00             = ix0 - jx0;
353             dy00             = iy0 - jy0;
354             dz00             = iz0 - jz0;
355
356             /* Calculate squared distance and things based on it */
357             rsq00            = dx00*dx00+dy00*dy00+dz00*dz00;
358
359             rinv00           = gmx_invsqrt(rsq00);
360
361             rinvsq00         = rinv00*rinv00;
362
363             /* Load parameters for j particles */
364             vdwjidx0         = 2*vdwtype[jnr+0];
365
366             /**************************
367              * CALCULATE INTERACTIONS *
368              **************************/
369
370             if (rsq00<rcutoff2)
371             {
372
373             r00              = rsq00*rinv00;
374
375             c6_00            = vdwparam[vdwioffset0+vdwjidx0];
376             c12_00           = vdwparam[vdwioffset0+vdwjidx0+1];
377             c6grid_00        = vdwgridparam[vdwioffset0+vdwjidx0];
378
379             rinvsix          = rinvsq00*rinvsq00*rinvsq00;
380             ewcljrsq         = ewclj2*rsq00;
381             exponent         = exp(-ewcljrsq);
382             poly             = exponent*(1.0 + ewcljrsq + ewcljrsq*ewcljrsq*0.5);
383             fvdw             = (((c12_00*rinvsix - c6_00 + c6grid_00*(1.0-poly))*rinvsix) - c6grid_00*(1.0/6.0)*exponent*ewclj6)*rinvsq00;
384
385             fscal            = fvdw;
386
387             /* Calculate temporary vectorial force */
388             tx               = fscal*dx00;
389             ty               = fscal*dy00;
390             tz               = fscal*dz00;
391
392             /* Update vectorial force */
393             fix0            += tx;
394             fiy0            += ty;
395             fiz0            += tz;
396             f[j_coord_offset+DIM*0+XX] -= tx;
397             f[j_coord_offset+DIM*0+YY] -= ty;
398             f[j_coord_offset+DIM*0+ZZ] -= tz;
399
400             }
401
402             /* Inner loop uses 44 flops */
403         }
404         /* End of innermost loop */
405
406         tx = ty = tz = 0;
407         f[i_coord_offset+DIM*0+XX] += fix0;
408         f[i_coord_offset+DIM*0+YY] += fiy0;
409         f[i_coord_offset+DIM*0+ZZ] += fiz0;
410         tx                         += fix0;
411         ty                         += fiy0;
412         tz                         += fiz0;
413         fshift[i_shift_offset+XX]  += tx;
414         fshift[i_shift_offset+YY]  += ty;
415         fshift[i_shift_offset+ZZ]  += tz;
416
417         /* Increment number of inner iterations */
418         inneriter                  += j_index_end - j_index_start;
419
420         /* Outer loop uses 12 flops */
421     }
422
423     /* Increment number of outer iterations */
424     outeriter        += nri;
425
426     /* Update outer/inner flops */
427
428     inc_nrnb(nrnb,eNR_NBKERNEL_VDW_F,outeriter*12 + inneriter*44);
429 }