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