Remove all unnecessary HAVE_CONFIG_H
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_c / nb_kernel_ElecEwSh_VdwLJSh_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_ElecEwSh_VdwLJSh_GeomP1P1_VF_c
49  * Electrostatics interaction: Ewald
50  * VdW interaction:            LennardJones
51  * Geometry:                   Particle-Particle
52  * Calculate force/pot:        PotentialAndForce
53  */
54 void
55 nb_kernel_ElecEwSh_VdwLJSh_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     int              ewitab;
82     real             ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace;
83     real             *ewtab;
84
85     x                = xx[0];
86     f                = ff[0];
87
88     nri              = nlist->nri;
89     iinr             = nlist->iinr;
90     jindex           = nlist->jindex;
91     jjnr             = nlist->jjnr;
92     shiftidx         = nlist->shift;
93     gid              = nlist->gid;
94     shiftvec         = fr->shift_vec[0];
95     fshift           = fr->fshift[0];
96     facel            = fr->epsfac;
97     charge           = mdatoms->chargeA;
98     nvdwtype         = fr->ntype;
99     vdwparam         = fr->nbfp;
100     vdwtype          = mdatoms->typeA;
101
102     sh_ewald         = fr->ic->sh_ewald;
103     ewtab            = fr->ic->tabq_coul_FDV0;
104     ewtabscale       = fr->ic->tabq_scale;
105     ewtabhalfspace   = 0.5/ewtabscale;
106
107     /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
108     rcutoff          = fr->rcoulomb;
109     rcutoff2         = rcutoff*rcutoff;
110
111     sh_vdw_invrcut6  = fr->ic->sh_invrc6;
112     rvdw             = fr->rvdw;
113
114     outeriter        = 0;
115     inneriter        = 0;
116
117     /* Start outer loop over neighborlists */
118     for(iidx=0; iidx<nri; iidx++)
119     {
120         /* Load shift vector for this list */
121         i_shift_offset   = DIM*shiftidx[iidx];
122         shX              = shiftvec[i_shift_offset+XX];
123         shY              = shiftvec[i_shift_offset+YY];
124         shZ              = shiftvec[i_shift_offset+ZZ];
125
126         /* Load limits for loop over neighbors */
127         j_index_start    = jindex[iidx];
128         j_index_end      = jindex[iidx+1];
129
130         /* Get outer coordinate index */
131         inr              = iinr[iidx];
132         i_coord_offset   = DIM*inr;
133
134         /* Load i particle coords and add shift vector */
135         ix0              = shX + x[i_coord_offset+DIM*0+XX];
136         iy0              = shY + x[i_coord_offset+DIM*0+YY];
137         iz0              = shZ + x[i_coord_offset+DIM*0+ZZ];
138
139         fix0             = 0.0;
140         fiy0             = 0.0;
141         fiz0             = 0.0;
142
143         /* Load parameters for i particles */
144         iq0              = facel*charge[inr+0];
145         vdwioffset0      = 2*nvdwtype*vdwtype[inr+0];
146
147         /* Reset potential sums */
148         velecsum         = 0.0;
149         vvdwsum          = 0.0;
150
151         /* Start inner kernel loop */
152         for(jidx=j_index_start; jidx<j_index_end; jidx++)
153         {
154             /* Get j neighbor index, and coordinate index */
155             jnr              = jjnr[jidx];
156             j_coord_offset   = DIM*jnr;
157
158             /* load j atom coordinates */
159             jx0              = x[j_coord_offset+DIM*0+XX];
160             jy0              = x[j_coord_offset+DIM*0+YY];
161             jz0              = x[j_coord_offset+DIM*0+ZZ];
162
163             /* Calculate displacement vector */
164             dx00             = ix0 - jx0;
165             dy00             = iy0 - jy0;
166             dz00             = iz0 - jz0;
167
168             /* Calculate squared distance and things based on it */
169             rsq00            = dx00*dx00+dy00*dy00+dz00*dz00;
170
171             rinv00           = gmx_invsqrt(rsq00);
172
173             rinvsq00         = rinv00*rinv00;
174
175             /* Load parameters for j particles */
176             jq0              = charge[jnr+0];
177             vdwjidx0         = 2*vdwtype[jnr+0];
178
179             /**************************
180              * CALCULATE INTERACTIONS *
181              **************************/
182
183             if (rsq00<rcutoff2)
184             {
185
186             r00              = rsq00*rinv00;
187
188             qq00             = iq0*jq0;
189             c6_00            = vdwparam[vdwioffset0+vdwjidx0];
190             c12_00           = vdwparam[vdwioffset0+vdwjidx0+1];
191
192             /* EWALD ELECTROSTATICS */
193
194             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
195             ewrt             = r00*ewtabscale;
196             ewitab           = ewrt;
197             eweps            = ewrt-ewitab;
198             ewitab           = 4*ewitab;
199             felec            = ewtab[ewitab]+eweps*ewtab[ewitab+1];
200             velec            = qq00*((rinv00-sh_ewald)-(ewtab[ewitab+2]-ewtabhalfspace*eweps*(ewtab[ewitab]+felec)));
201             felec            = qq00*rinv00*(rinvsq00-felec);
202
203             /* LENNARD-JONES DISPERSION/REPULSION */
204
205             rinvsix          = rinvsq00*rinvsq00*rinvsq00;
206             vvdw6            = c6_00*rinvsix;
207             vvdw12           = c12_00*rinvsix*rinvsix;
208             vvdw             = (vvdw12 - c12_00*sh_vdw_invrcut6*sh_vdw_invrcut6)*(1.0/12.0) - (vvdw6 - c6_00*sh_vdw_invrcut6)*(1.0/6.0);
209             fvdw             = (vvdw12-vvdw6)*rinvsq00;
210
211             /* Update potential sums from outer loop */
212             velecsum        += velec;
213             vvdwsum         += vvdw;
214
215             fscal            = felec+fvdw;
216
217             /* Calculate temporary vectorial force */
218             tx               = fscal*dx00;
219             ty               = fscal*dy00;
220             tz               = fscal*dz00;
221
222             /* Update vectorial force */
223             fix0            += tx;
224             fiy0            += ty;
225             fiz0            += tz;
226             f[j_coord_offset+DIM*0+XX] -= tx;
227             f[j_coord_offset+DIM*0+YY] -= ty;
228             f[j_coord_offset+DIM*0+ZZ] -= tz;
229
230             }
231
232             /* Inner loop uses 59 flops */
233         }
234         /* End of innermost loop */
235
236         tx = ty = tz = 0;
237         f[i_coord_offset+DIM*0+XX] += fix0;
238         f[i_coord_offset+DIM*0+YY] += fiy0;
239         f[i_coord_offset+DIM*0+ZZ] += fiz0;
240         tx                         += fix0;
241         ty                         += fiy0;
242         tz                         += fiz0;
243         fshift[i_shift_offset+XX]  += tx;
244         fshift[i_shift_offset+YY]  += ty;
245         fshift[i_shift_offset+ZZ]  += tz;
246
247         ggid                        = gid[iidx];
248         /* Update potential energies */
249         kernel_data->energygrp_elec[ggid] += velecsum;
250         kernel_data->energygrp_vdw[ggid] += vvdwsum;
251
252         /* Increment number of inner iterations */
253         inneriter                  += j_index_end - j_index_start;
254
255         /* Outer loop uses 15 flops */
256     }
257
258     /* Increment number of outer iterations */
259     outeriter        += nri;
260
261     /* Update outer/inner flops */
262
263     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_VF,outeriter*15 + inneriter*59);
264 }
265 /*
266  * Gromacs nonbonded kernel:   nb_kernel_ElecEwSh_VdwLJSh_GeomP1P1_F_c
267  * Electrostatics interaction: Ewald
268  * VdW interaction:            LennardJones
269  * Geometry:                   Particle-Particle
270  * Calculate force/pot:        Force
271  */
272 void
273 nb_kernel_ElecEwSh_VdwLJSh_GeomP1P1_F_c
274                     (t_nblist                    * gmx_restrict       nlist,
275                      rvec                        * gmx_restrict          xx,
276                      rvec                        * gmx_restrict          ff,
277                      t_forcerec                  * gmx_restrict          fr,
278                      t_mdatoms                   * gmx_restrict     mdatoms,
279                      nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
280                      t_nrnb                      * gmx_restrict        nrnb)
281 {
282     int              i_shift_offset,i_coord_offset,j_coord_offset;
283     int              j_index_start,j_index_end;
284     int              nri,inr,ggid,iidx,jidx,jnr,outeriter,inneriter;
285     real             shX,shY,shZ,tx,ty,tz,fscal,rcutoff,rcutoff2;
286     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
287     real             *shiftvec,*fshift,*x,*f;
288     int              vdwioffset0;
289     real             ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
290     int              vdwjidx0;
291     real             jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
292     real             dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00,cexp1_00,cexp2_00;
293     real             velec,felec,velecsum,facel,crf,krf,krf2;
294     real             *charge;
295     int              nvdwtype;
296     real             rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,br,vvdwexp,sh_vdw_invrcut6;
297     int              *vdwtype;
298     real             *vdwparam;
299     int              ewitab;
300     real             ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace;
301     real             *ewtab;
302
303     x                = xx[0];
304     f                = ff[0];
305
306     nri              = nlist->nri;
307     iinr             = nlist->iinr;
308     jindex           = nlist->jindex;
309     jjnr             = nlist->jjnr;
310     shiftidx         = nlist->shift;
311     gid              = nlist->gid;
312     shiftvec         = fr->shift_vec[0];
313     fshift           = fr->fshift[0];
314     facel            = fr->epsfac;
315     charge           = mdatoms->chargeA;
316     nvdwtype         = fr->ntype;
317     vdwparam         = fr->nbfp;
318     vdwtype          = mdatoms->typeA;
319
320     sh_ewald         = fr->ic->sh_ewald;
321     ewtab            = fr->ic->tabq_coul_F;
322     ewtabscale       = fr->ic->tabq_scale;
323     ewtabhalfspace   = 0.5/ewtabscale;
324
325     /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
326     rcutoff          = fr->rcoulomb;
327     rcutoff2         = rcutoff*rcutoff;
328
329     sh_vdw_invrcut6  = fr->ic->sh_invrc6;
330     rvdw             = fr->rvdw;
331
332     outeriter        = 0;
333     inneriter        = 0;
334
335     /* Start outer loop over neighborlists */
336     for(iidx=0; iidx<nri; iidx++)
337     {
338         /* Load shift vector for this list */
339         i_shift_offset   = DIM*shiftidx[iidx];
340         shX              = shiftvec[i_shift_offset+XX];
341         shY              = shiftvec[i_shift_offset+YY];
342         shZ              = shiftvec[i_shift_offset+ZZ];
343
344         /* Load limits for loop over neighbors */
345         j_index_start    = jindex[iidx];
346         j_index_end      = jindex[iidx+1];
347
348         /* Get outer coordinate index */
349         inr              = iinr[iidx];
350         i_coord_offset   = DIM*inr;
351
352         /* Load i particle coords and add shift vector */
353         ix0              = shX + x[i_coord_offset+DIM*0+XX];
354         iy0              = shY + x[i_coord_offset+DIM*0+YY];
355         iz0              = shZ + x[i_coord_offset+DIM*0+ZZ];
356
357         fix0             = 0.0;
358         fiy0             = 0.0;
359         fiz0             = 0.0;
360
361         /* Load parameters for i particles */
362         iq0              = facel*charge[inr+0];
363         vdwioffset0      = 2*nvdwtype*vdwtype[inr+0];
364
365         /* Start inner kernel loop */
366         for(jidx=j_index_start; jidx<j_index_end; jidx++)
367         {
368             /* Get j neighbor index, and coordinate index */
369             jnr              = jjnr[jidx];
370             j_coord_offset   = DIM*jnr;
371
372             /* load j atom coordinates */
373             jx0              = x[j_coord_offset+DIM*0+XX];
374             jy0              = x[j_coord_offset+DIM*0+YY];
375             jz0              = x[j_coord_offset+DIM*0+ZZ];
376
377             /* Calculate displacement vector */
378             dx00             = ix0 - jx0;
379             dy00             = iy0 - jy0;
380             dz00             = iz0 - jz0;
381
382             /* Calculate squared distance and things based on it */
383             rsq00            = dx00*dx00+dy00*dy00+dz00*dz00;
384
385             rinv00           = gmx_invsqrt(rsq00);
386
387             rinvsq00         = rinv00*rinv00;
388
389             /* Load parameters for j particles */
390             jq0              = charge[jnr+0];
391             vdwjidx0         = 2*vdwtype[jnr+0];
392
393             /**************************
394              * CALCULATE INTERACTIONS *
395              **************************/
396
397             if (rsq00<rcutoff2)
398             {
399
400             r00              = rsq00*rinv00;
401
402             qq00             = iq0*jq0;
403             c6_00            = vdwparam[vdwioffset0+vdwjidx0];
404             c12_00           = vdwparam[vdwioffset0+vdwjidx0+1];
405
406             /* EWALD ELECTROSTATICS */
407
408             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
409             ewrt             = r00*ewtabscale;
410             ewitab           = ewrt;
411             eweps            = ewrt-ewitab;
412             felec            = (1.0-eweps)*ewtab[ewitab]+eweps*ewtab[ewitab+1];
413             felec            = qq00*rinv00*(rinvsq00-felec);
414
415             /* LENNARD-JONES DISPERSION/REPULSION */
416
417             rinvsix          = rinvsq00*rinvsq00*rinvsq00;
418             fvdw             = (c12_00*rinvsix-c6_00)*rinvsix*rinvsq00;
419
420             fscal            = felec+fvdw;
421
422             /* Calculate temporary vectorial force */
423             tx               = fscal*dx00;
424             ty               = fscal*dy00;
425             tz               = fscal*dz00;
426
427             /* Update vectorial force */
428             fix0            += tx;
429             fiy0            += ty;
430             fiz0            += tz;
431             f[j_coord_offset+DIM*0+XX] -= tx;
432             f[j_coord_offset+DIM*0+YY] -= ty;
433             f[j_coord_offset+DIM*0+ZZ] -= tz;
434
435             }
436
437             /* Inner loop uses 41 flops */
438         }
439         /* End of innermost loop */
440
441         tx = ty = tz = 0;
442         f[i_coord_offset+DIM*0+XX] += fix0;
443         f[i_coord_offset+DIM*0+YY] += fiy0;
444         f[i_coord_offset+DIM*0+ZZ] += fiz0;
445         tx                         += fix0;
446         ty                         += fiy0;
447         tz                         += fiz0;
448         fshift[i_shift_offset+XX]  += tx;
449         fshift[i_shift_offset+YY]  += ty;
450         fshift[i_shift_offset+ZZ]  += tz;
451
452         /* Increment number of inner iterations */
453         inneriter                  += j_index_end - j_index_start;
454
455         /* Outer loop uses 13 flops */
456     }
457
458     /* Increment number of outer iterations */
459     outeriter        += nri;
460
461     /* Update outer/inner flops */
462
463     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_F,outeriter*13 + inneriter*41);
464 }