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