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