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