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