190902b76ea9e2de51487aafb62a099694d9a3e0
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_c / nb_kernel_ElecCoul_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 #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_ElecCoul_VdwBham_GeomP1P1_VF_c
51  * Electrostatics interaction: Coulomb
52  * VdW interaction:            Buckingham
53  * Geometry:                   Particle-Particle
54  * Calculate force/pot:        PotentialAndForce
55  */
56 void
57 nb_kernel_ElecCoul_VdwBham_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     real             velec,felec,velecsum,facel,crf,krf,krf2;
78     real             *charge;
79     int              nvdwtype;
80     real             rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,br,vvdwexp,sh_vdw_invrcut6;
81     int              *vdwtype;
82     real             *vdwparam;
83
84     x                = xx[0];
85     f                = ff[0];
86
87     nri              = nlist->nri;
88     iinr             = nlist->iinr;
89     jindex           = nlist->jindex;
90     jjnr             = nlist->jjnr;
91     shiftidx         = nlist->shift;
92     gid              = nlist->gid;
93     shiftvec         = fr->shift_vec[0];
94     fshift           = fr->fshift[0];
95     facel            = fr->epsfac;
96     charge           = mdatoms->chargeA;
97     nvdwtype         = fr->ntype;
98     vdwparam         = fr->nbfp;
99     vdwtype          = mdatoms->typeA;
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         iq0              = facel*charge[inr+0];
132         vdwioffset0      = 3*nvdwtype*vdwtype[inr+0];
133
134         /* Reset potential sums */
135         velecsum         = 0.0;
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             rinv00           = gmx_invsqrt(rsq00);
159
160             rinvsq00         = rinv00*rinv00;
161
162             /* Load parameters for j particles */
163             jq0              = charge[jnr+0];
164             vdwjidx0         = 3*vdwtype[jnr+0];
165
166             /**************************
167              * CALCULATE INTERACTIONS *
168              **************************/
169
170             r00              = rsq00*rinv00;
171
172             qq00             = iq0*jq0;
173             c6_00            = vdwparam[vdwioffset0+vdwjidx0];
174             cexp1_00         = vdwparam[vdwioffset0+vdwjidx0+1];
175             cexp2_00         = vdwparam[vdwioffset0+vdwjidx0+2];
176
177             /* COULOMB ELECTROSTATICS */
178             velec            = qq00*rinv00;
179             felec            = velec*rinvsq00;
180
181             /* BUCKINGHAM DISPERSION/REPULSION */
182             rinvsix          = rinvsq00*rinvsq00*rinvsq00;
183             vvdw6            = c6_00*rinvsix;
184             br               = cexp2_00*r00;
185             vvdwexp          = cexp1_00*exp(-br);
186             vvdw             = vvdwexp - vvdw6*(1.0/6.0);
187             fvdw             = (br*vvdwexp-vvdw6)*rinvsq00;
188
189             /* Update potential sums from outer loop */
190             velecsum        += velec;
191             vvdwsum         += vvdw;
192
193             fscal            = felec+fvdw;
194
195             /* Calculate temporary vectorial force */
196             tx               = fscal*dx00;
197             ty               = fscal*dy00;
198             tz               = fscal*dz00;
199
200             /* Update vectorial force */
201             fix0            += tx;
202             fiy0            += ty;
203             fiz0            += tz;
204             f[j_coord_offset+DIM*0+XX] -= tx;
205             f[j_coord_offset+DIM*0+YY] -= ty;
206             f[j_coord_offset+DIM*0+ZZ] -= tz;
207
208             /* Inner loop uses 67 flops */
209         }
210         /* End of innermost loop */
211
212         tx = ty = tz = 0;
213         f[i_coord_offset+DIM*0+XX] += fix0;
214         f[i_coord_offset+DIM*0+YY] += fiy0;
215         f[i_coord_offset+DIM*0+ZZ] += fiz0;
216         tx                         += fix0;
217         ty                         += fiy0;
218         tz                         += fiz0;
219         fshift[i_shift_offset+XX]  += tx;
220         fshift[i_shift_offset+YY]  += ty;
221         fshift[i_shift_offset+ZZ]  += tz;
222
223         ggid                        = gid[iidx];
224         /* Update potential energies */
225         kernel_data->energygrp_elec[ggid] += velecsum;
226         kernel_data->energygrp_vdw[ggid] += vvdwsum;
227
228         /* Increment number of inner iterations */
229         inneriter                  += j_index_end - j_index_start;
230
231         /* Outer loop uses 15 flops */
232     }
233
234     /* Increment number of outer iterations */
235     outeriter        += nri;
236
237     /* Update outer/inner flops */
238
239     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_VF,outeriter*15 + inneriter*67);
240 }
241 /*
242  * Gromacs nonbonded kernel:   nb_kernel_ElecCoul_VdwBham_GeomP1P1_F_c
243  * Electrostatics interaction: Coulomb
244  * VdW interaction:            Buckingham
245  * Geometry:                   Particle-Particle
246  * Calculate force/pot:        Force
247  */
248 void
249 nb_kernel_ElecCoul_VdwBham_GeomP1P1_F_c
250                     (t_nblist                    * gmx_restrict       nlist,
251                      rvec                        * gmx_restrict          xx,
252                      rvec                        * gmx_restrict          ff,
253                      t_forcerec                  * gmx_restrict          fr,
254                      t_mdatoms                   * gmx_restrict     mdatoms,
255                      nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
256                      t_nrnb                      * gmx_restrict        nrnb)
257 {
258     int              i_shift_offset,i_coord_offset,j_coord_offset;
259     int              j_index_start,j_index_end;
260     int              nri,inr,ggid,iidx,jidx,jnr,outeriter,inneriter;
261     real             shX,shY,shZ,tx,ty,tz,fscal,rcutoff,rcutoff2;
262     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
263     real             *shiftvec,*fshift,*x,*f;
264     int              vdwioffset0;
265     real             ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
266     int              vdwjidx0;
267     real             jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
268     real             dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00,cexp1_00,cexp2_00;
269     real             velec,felec,velecsum,facel,crf,krf,krf2;
270     real             *charge;
271     int              nvdwtype;
272     real             rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,br,vvdwexp,sh_vdw_invrcut6;
273     int              *vdwtype;
274     real             *vdwparam;
275
276     x                = xx[0];
277     f                = ff[0];
278
279     nri              = nlist->nri;
280     iinr             = nlist->iinr;
281     jindex           = nlist->jindex;
282     jjnr             = nlist->jjnr;
283     shiftidx         = nlist->shift;
284     gid              = nlist->gid;
285     shiftvec         = fr->shift_vec[0];
286     fshift           = fr->fshift[0];
287     facel            = fr->epsfac;
288     charge           = mdatoms->chargeA;
289     nvdwtype         = fr->ntype;
290     vdwparam         = fr->nbfp;
291     vdwtype          = mdatoms->typeA;
292
293     outeriter        = 0;
294     inneriter        = 0;
295
296     /* Start outer loop over neighborlists */
297     for(iidx=0; iidx<nri; iidx++)
298     {
299         /* Load shift vector for this list */
300         i_shift_offset   = DIM*shiftidx[iidx];
301         shX              = shiftvec[i_shift_offset+XX];
302         shY              = shiftvec[i_shift_offset+YY];
303         shZ              = shiftvec[i_shift_offset+ZZ];
304
305         /* Load limits for loop over neighbors */
306         j_index_start    = jindex[iidx];
307         j_index_end      = jindex[iidx+1];
308
309         /* Get outer coordinate index */
310         inr              = iinr[iidx];
311         i_coord_offset   = DIM*inr;
312
313         /* Load i particle coords and add shift vector */
314         ix0              = shX + x[i_coord_offset+DIM*0+XX];
315         iy0              = shY + x[i_coord_offset+DIM*0+YY];
316         iz0              = shZ + x[i_coord_offset+DIM*0+ZZ];
317
318         fix0             = 0.0;
319         fiy0             = 0.0;
320         fiz0             = 0.0;
321
322         /* Load parameters for i particles */
323         iq0              = facel*charge[inr+0];
324         vdwioffset0      = 3*nvdwtype*vdwtype[inr+0];
325
326         /* Start inner kernel loop */
327         for(jidx=j_index_start; jidx<j_index_end; jidx++)
328         {
329             /* Get j neighbor index, and coordinate index */
330             jnr              = jjnr[jidx];
331             j_coord_offset   = DIM*jnr;
332
333             /* load j atom coordinates */
334             jx0              = x[j_coord_offset+DIM*0+XX];
335             jy0              = x[j_coord_offset+DIM*0+YY];
336             jz0              = x[j_coord_offset+DIM*0+ZZ];
337
338             /* Calculate displacement vector */
339             dx00             = ix0 - jx0;
340             dy00             = iy0 - jy0;
341             dz00             = iz0 - jz0;
342
343             /* Calculate squared distance and things based on it */
344             rsq00            = dx00*dx00+dy00*dy00+dz00*dz00;
345
346             rinv00           = gmx_invsqrt(rsq00);
347
348             rinvsq00         = rinv00*rinv00;
349
350             /* Load parameters for j particles */
351             jq0              = charge[jnr+0];
352             vdwjidx0         = 3*vdwtype[jnr+0];
353
354             /**************************
355              * CALCULATE INTERACTIONS *
356              **************************/
357
358             r00              = rsq00*rinv00;
359
360             qq00             = iq0*jq0;
361             c6_00            = vdwparam[vdwioffset0+vdwjidx0];
362             cexp1_00         = vdwparam[vdwioffset0+vdwjidx0+1];
363             cexp2_00         = vdwparam[vdwioffset0+vdwjidx0+2];
364
365             /* COULOMB ELECTROSTATICS */
366             velec            = qq00*rinv00;
367             felec            = velec*rinvsq00;
368
369             /* BUCKINGHAM DISPERSION/REPULSION */
370             rinvsix          = rinvsq00*rinvsq00*rinvsq00;
371             vvdw6            = c6_00*rinvsix;
372             br               = cexp2_00*r00;
373             vvdwexp          = cexp1_00*exp(-br);
374             fvdw             = (br*vvdwexp-vvdw6)*rinvsq00;
375
376             fscal            = felec+fvdw;
377
378             /* Calculate temporary vectorial force */
379             tx               = fscal*dx00;
380             ty               = fscal*dy00;
381             tz               = fscal*dz00;
382
383             /* Update vectorial force */
384             fix0            += tx;
385             fiy0            += ty;
386             fiz0            += tz;
387             f[j_coord_offset+DIM*0+XX] -= tx;
388             f[j_coord_offset+DIM*0+YY] -= ty;
389             f[j_coord_offset+DIM*0+ZZ] -= tz;
390
391             /* Inner loop uses 63 flops */
392         }
393         /* End of innermost loop */
394
395         tx = ty = tz = 0;
396         f[i_coord_offset+DIM*0+XX] += fix0;
397         f[i_coord_offset+DIM*0+YY] += fiy0;
398         f[i_coord_offset+DIM*0+ZZ] += fiz0;
399         tx                         += fix0;
400         ty                         += fiy0;
401         tz                         += fiz0;
402         fshift[i_shift_offset+XX]  += tx;
403         fshift[i_shift_offset+YY]  += ty;
404         fshift[i_shift_offset+ZZ]  += tz;
405
406         /* Increment number of inner iterations */
407         inneriter                  += j_index_end - j_index_start;
408
409         /* Outer loop uses 13 flops */
410     }
411
412     /* Increment number of outer iterations */
413     outeriter        += nri;
414
415     /* Update outer/inner flops */
416
417     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_F,outeriter*13 + inneriter*63);
418 }