Compile nonbonded kernels as C++
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_c / nb_kernel_ElecNone_VdwLJ_GeomP1P1_c.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2013,2014.2015,2017,2018, 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 "gmxpre.h"
39
40 #include "config.h"
41
42 #include <math.h>
43
44 #include "../nb_kernel.h"
45 #include "gromacs/gmxlib/nrnb.h"
46
47 /*
48  * Gromacs nonbonded kernel:   nb_kernel_ElecNone_VdwLJ_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_VdwLJ_GeomP1P1_VF_c
56                     (t_nblist                    * gmx_restrict       nlist,
57                      rvec                        * gmx_restrict          xx,
58                      rvec                        * gmx_restrict          ff,
59                      struct 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      = 2*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             rinvsq00         = 1.0/rsq00;
151
152             /* Load parameters for j particles */
153             vdwjidx0         = 2*vdwtype[jnr+0];
154
155             /**************************
156              * CALCULATE INTERACTIONS *
157              **************************/
158
159             c6_00            = vdwparam[vdwioffset0+vdwjidx0];
160             c12_00           = vdwparam[vdwioffset0+vdwjidx0+1];
161
162             /* LENNARD-JONES DISPERSION/REPULSION */
163
164             rinvsix          = rinvsq00*rinvsq00*rinvsq00;
165             vvdw6            = c6_00*rinvsix;
166             vvdw12           = c12_00*rinvsix*rinvsix;
167             vvdw             = vvdw12*(1.0/12.0) - vvdw6*(1.0/6.0);
168             fvdw             = (vvdw12-vvdw6)*rinvsq00;
169
170             /* Update potential sums from outer loop */
171             vvdwsum         += vvdw;
172
173             fscal            = fvdw;
174
175             /* Calculate temporary vectorial force */
176             tx               = fscal*dx00;
177             ty               = fscal*dy00;
178             tz               = fscal*dz00;
179
180             /* Update vectorial force */
181             fix0            += tx;
182             fiy0            += ty;
183             fiz0            += tz;
184             f[j_coord_offset+DIM*0+XX] -= tx;
185             f[j_coord_offset+DIM*0+YY] -= ty;
186             f[j_coord_offset+DIM*0+ZZ] -= tz;
187
188             /* Inner loop uses 32 flops */
189         }
190         /* End of innermost loop */
191
192         tx = ty = tz = 0;
193         f[i_coord_offset+DIM*0+XX] += fix0;
194         f[i_coord_offset+DIM*0+YY] += fiy0;
195         f[i_coord_offset+DIM*0+ZZ] += fiz0;
196         tx                         += fix0;
197         ty                         += fiy0;
198         tz                         += fiz0;
199         fshift[i_shift_offset+XX]  += tx;
200         fshift[i_shift_offset+YY]  += ty;
201         fshift[i_shift_offset+ZZ]  += tz;
202
203         ggid                        = gid[iidx];
204         /* Update potential energies */
205         kernel_data->energygrp_vdw[ggid] += vvdwsum;
206
207         /* Increment number of inner iterations */
208         inneriter                  += j_index_end - j_index_start;
209
210         /* Outer loop uses 13 flops */
211     }
212
213     /* Increment number of outer iterations */
214     outeriter        += nri;
215
216     /* Update outer/inner flops */
217
218     inc_nrnb(nrnb,eNR_NBKERNEL_VDW_VF,outeriter*13 + inneriter*32);
219 }
220 /*
221  * Gromacs nonbonded kernel:   nb_kernel_ElecNone_VdwLJ_GeomP1P1_F_c
222  * Electrostatics interaction: None
223  * VdW interaction:            LennardJones
224  * Geometry:                   Particle-Particle
225  * Calculate force/pot:        Force
226  */
227 void
228 nb_kernel_ElecNone_VdwLJ_GeomP1P1_F_c
229                     (t_nblist                    * gmx_restrict       nlist,
230                      rvec                        * gmx_restrict          xx,
231                      rvec                        * gmx_restrict          ff,
232                      struct t_forcerec           * gmx_restrict          fr,
233                      t_mdatoms                   * gmx_restrict     mdatoms,
234                      nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
235                      t_nrnb                      * gmx_restrict        nrnb)
236 {
237     int              i_shift_offset,i_coord_offset,j_coord_offset;
238     int              j_index_start,j_index_end;
239     int              nri,inr,ggid,iidx,jidx,jnr,outeriter,inneriter;
240     real             shX,shY,shZ,tx,ty,tz,fscal,rcutoff,rcutoff2;
241     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
242     real             *shiftvec,*fshift,*x,*f;
243     int              vdwioffset0;
244     real             ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
245     int              vdwjidx0;
246     real             jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
247     real             dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00,cexp1_00,cexp2_00;
248     int              nvdwtype;
249     real             rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,br,vvdwexp,sh_vdw_invrcut6;
250     int              *vdwtype;
251     real             *vdwparam;
252
253     x                = xx[0];
254     f                = ff[0];
255
256     nri              = nlist->nri;
257     iinr             = nlist->iinr;
258     jindex           = nlist->jindex;
259     jjnr             = nlist->jjnr;
260     shiftidx         = nlist->shift;
261     gid              = nlist->gid;
262     shiftvec         = fr->shift_vec[0];
263     fshift           = fr->fshift[0];
264     nvdwtype         = fr->ntype;
265     vdwparam         = fr->nbfp;
266     vdwtype          = mdatoms->typeA;
267
268     outeriter        = 0;
269     inneriter        = 0;
270
271     /* Start outer loop over neighborlists */
272     for(iidx=0; iidx<nri; iidx++)
273     {
274         /* Load shift vector for this list */
275         i_shift_offset   = DIM*shiftidx[iidx];
276         shX              = shiftvec[i_shift_offset+XX];
277         shY              = shiftvec[i_shift_offset+YY];
278         shZ              = shiftvec[i_shift_offset+ZZ];
279
280         /* Load limits for loop over neighbors */
281         j_index_start    = jindex[iidx];
282         j_index_end      = jindex[iidx+1];
283
284         /* Get outer coordinate index */
285         inr              = iinr[iidx];
286         i_coord_offset   = DIM*inr;
287
288         /* Load i particle coords and add shift vector */
289         ix0              = shX + x[i_coord_offset+DIM*0+XX];
290         iy0              = shY + x[i_coord_offset+DIM*0+YY];
291         iz0              = shZ + x[i_coord_offset+DIM*0+ZZ];
292
293         fix0             = 0.0;
294         fiy0             = 0.0;
295         fiz0             = 0.0;
296
297         /* Load parameters for i particles */
298         vdwioffset0      = 2*nvdwtype*vdwtype[inr+0];
299
300         /* Start inner kernel loop */
301         for(jidx=j_index_start; jidx<j_index_end; jidx++)
302         {
303             /* Get j neighbor index, and coordinate index */
304             jnr              = jjnr[jidx];
305             j_coord_offset   = DIM*jnr;
306
307             /* load j atom coordinates */
308             jx0              = x[j_coord_offset+DIM*0+XX];
309             jy0              = x[j_coord_offset+DIM*0+YY];
310             jz0              = x[j_coord_offset+DIM*0+ZZ];
311
312             /* Calculate displacement vector */
313             dx00             = ix0 - jx0;
314             dy00             = iy0 - jy0;
315             dz00             = iz0 - jz0;
316
317             /* Calculate squared distance and things based on it */
318             rsq00            = dx00*dx00+dy00*dy00+dz00*dz00;
319
320             rinvsq00         = 1.0/rsq00;
321
322             /* Load parameters for j particles */
323             vdwjidx0         = 2*vdwtype[jnr+0];
324
325             /**************************
326              * CALCULATE INTERACTIONS *
327              **************************/
328
329             c6_00            = vdwparam[vdwioffset0+vdwjidx0];
330             c12_00           = vdwparam[vdwioffset0+vdwjidx0+1];
331
332             /* LENNARD-JONES DISPERSION/REPULSION */
333
334             rinvsix          = rinvsq00*rinvsq00*rinvsq00;
335             fvdw             = (c12_00*rinvsix-c6_00)*rinvsix*rinvsq00;
336
337             fscal            = fvdw;
338
339             /* Calculate temporary vectorial force */
340             tx               = fscal*dx00;
341             ty               = fscal*dy00;
342             tz               = fscal*dz00;
343
344             /* Update vectorial force */
345             fix0            += tx;
346             fiy0            += ty;
347             fiz0            += tz;
348             f[j_coord_offset+DIM*0+XX] -= tx;
349             f[j_coord_offset+DIM*0+YY] -= ty;
350             f[j_coord_offset+DIM*0+ZZ] -= tz;
351
352             /* Inner loop uses 27 flops */
353         }
354         /* End of innermost loop */
355
356         tx = ty = tz = 0;
357         f[i_coord_offset+DIM*0+XX] += fix0;
358         f[i_coord_offset+DIM*0+YY] += fiy0;
359         f[i_coord_offset+DIM*0+ZZ] += fiz0;
360         tx                         += fix0;
361         ty                         += fiy0;
362         tz                         += fiz0;
363         fshift[i_shift_offset+XX]  += tx;
364         fshift[i_shift_offset+YY]  += ty;
365         fshift[i_shift_offset+ZZ]  += tz;
366
367         /* Increment number of inner iterations */
368         inneriter                  += j_index_end - j_index_start;
369
370         /* Outer loop uses 12 flops */
371     }
372
373     /* Increment number of outer iterations */
374     outeriter        += nri;
375
376     /* Update outer/inner flops */
377
378     inc_nrnb(nrnb,eNR_NBKERNEL_VDW_F,outeriter*12 + inneriter*27);
379 }