664442fb632cb73c3ea48660c9d428ed70abd2f8
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_c / nb_kernel_ElecEwSh_VdwNone_GeomP1P1_c.c
1 /*
2  * Note: this file was generated by the Gromacs c kernel generator.
3  *
4  *                This source code is part of
5  *
6  *                 G   R   O   M   A   C   S
7  *
8  * Copyright (c) 2001-2012, The GROMACS Development Team
9  *
10  * Gromacs is a library for molecular simulation and trajectory analysis,
11  * written by Erik Lindahl, David van der Spoel, Berk Hess, and others - for
12  * a full list of developers and information, check out http://www.gromacs.org
13  *
14  * This program is free software; you can redistribute it and/or modify it under
15  * the terms of the GNU Lesser General Public License as published by the Free
16  * Software Foundation; either version 2 of the License, or (at your option) any
17  * later version.
18  *
19  * To help fund GROMACS development, we humbly ask that you cite
20  * the papers people have written on it - you can find them on the website.
21  */
22 #ifdef HAVE_CONFIG_H
23 #include <config.h>
24 #endif
25
26 #include <math.h>
27
28 #include "../nb_kernel.h"
29 #include "types/simple.h"
30 #include "vec.h"
31 #include "nrnb.h"
32
33 /*
34  * Gromacs nonbonded kernel:   nb_kernel_ElecEwSh_VdwNone_GeomP1P1_VF_c
35  * Electrostatics interaction: Ewald
36  * VdW interaction:            None
37  * Geometry:                   Particle-Particle
38  * Calculate force/pot:        PotentialAndForce
39  */
40 void
41 nb_kernel_ElecEwSh_VdwNone_GeomP1P1_VF_c
42                     (t_nblist * gmx_restrict                nlist,
43                      rvec * gmx_restrict                    xx,
44                      rvec * gmx_restrict                    ff,
45                      t_forcerec * gmx_restrict              fr,
46                      t_mdatoms * gmx_restrict               mdatoms,
47                      nb_kernel_data_t * gmx_restrict        kernel_data,
48                      t_nrnb * gmx_restrict                  nrnb)
49 {
50     int              i_shift_offset,i_coord_offset,j_coord_offset;
51     int              j_index_start,j_index_end;
52     int              nri,inr,ggid,iidx,jidx,jnr,outeriter,inneriter;
53     real             shX,shY,shZ,tx,ty,tz,fscal,rcutoff,rcutoff2;
54     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
55     real             *shiftvec,*fshift,*x,*f;
56     int              vdwioffset0;
57     real             ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
58     int              vdwjidx0;
59     real             jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
60     real             dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00,cexp1_00,cexp2_00;
61     real             velec,felec,velecsum,facel,crf,krf,krf2;
62     real             *charge;
63     int              ewitab;
64     real             ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace;
65     real             *ewtab;
66
67     x                = xx[0];
68     f                = ff[0];
69
70     nri              = nlist->nri;
71     iinr             = nlist->iinr;
72     jindex           = nlist->jindex;
73     jjnr             = nlist->jjnr;
74     shiftidx         = nlist->shift;
75     gid              = nlist->gid;
76     shiftvec         = fr->shift_vec[0];
77     fshift           = fr->fshift[0];
78     facel            = fr->epsfac;
79     charge           = mdatoms->chargeA;
80
81     sh_ewald         = fr->ic->sh_ewald;
82     ewtab            = fr->ic->tabq_coul_FDV0;
83     ewtabscale       = fr->ic->tabq_scale;
84     ewtabhalfspace   = 0.5/ewtabscale;
85
86     /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
87     rcutoff          = fr->rcoulomb;
88     rcutoff2         = rcutoff*rcutoff;
89
90     outeriter        = 0;
91     inneriter        = 0;
92
93     /* Start outer loop over neighborlists */
94     for(iidx=0; iidx<nri; iidx++)
95     {
96         /* Load shift vector for this list */
97         i_shift_offset   = DIM*shiftidx[iidx];
98         shX              = shiftvec[i_shift_offset+XX];
99         shY              = shiftvec[i_shift_offset+YY];
100         shZ              = shiftvec[i_shift_offset+ZZ];
101
102         /* Load limits for loop over neighbors */
103         j_index_start    = jindex[iidx];
104         j_index_end      = jindex[iidx+1];
105
106         /* Get outer coordinate index */
107         inr              = iinr[iidx];
108         i_coord_offset   = DIM*inr;
109
110         /* Load i particle coords and add shift vector */
111         ix0              = shX + x[i_coord_offset+DIM*0+XX];
112         iy0              = shY + x[i_coord_offset+DIM*0+YY];
113         iz0              = shZ + x[i_coord_offset+DIM*0+ZZ];
114
115         fix0             = 0.0;
116         fiy0             = 0.0;
117         fiz0             = 0.0;
118
119         /* Load parameters for i particles */
120         iq0              = facel*charge[inr+0];
121
122         /* Reset potential sums */
123         velecsum         = 0.0;
124
125         /* Start inner kernel loop */
126         for(jidx=j_index_start; jidx<j_index_end; jidx++)
127         {
128             /* Get j neighbor index, and coordinate index */
129             jnr              = jjnr[jidx];
130             j_coord_offset   = DIM*jnr;
131
132             /* load j atom coordinates */
133             jx0              = x[j_coord_offset+DIM*0+XX];
134             jy0              = x[j_coord_offset+DIM*0+YY];
135             jz0              = x[j_coord_offset+DIM*0+ZZ];
136
137             /* Calculate displacement vector */
138             dx00             = ix0 - jx0;
139             dy00             = iy0 - jy0;
140             dz00             = iz0 - jz0;
141
142             /* Calculate squared distance and things based on it */
143             rsq00            = dx00*dx00+dy00*dy00+dz00*dz00;
144
145             rinv00           = gmx_invsqrt(rsq00);
146
147             rinvsq00         = rinv00*rinv00;
148
149             /* Load parameters for j particles */
150             jq0              = charge[jnr+0];
151
152             /**************************
153              * CALCULATE INTERACTIONS *
154              **************************/
155
156             if (rsq00<rcutoff2)
157             {
158
159             r00              = rsq00*rinv00;
160
161             qq00             = iq0*jq0;
162
163             /* EWALD ELECTROSTATICS */
164
165             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
166             ewrt             = r00*ewtabscale;
167             ewitab           = ewrt;
168             eweps            = ewrt-ewitab;
169             ewitab           = 4*ewitab;
170             felec            = ewtab[ewitab]+eweps*ewtab[ewitab+1];
171             velec            = qq00*((rinv00-sh_ewald)-(ewtab[ewitab+2]-ewtabhalfspace*eweps*(ewtab[ewitab]+felec)));
172             felec            = qq00*rinv00*(rinvsq00-felec);
173
174             /* Update potential sums from outer loop */
175             velecsum        += velec;
176
177             fscal            = felec;
178
179             /* Calculate temporary vectorial force */
180             tx               = fscal*dx00;
181             ty               = fscal*dy00;
182             tz               = fscal*dz00;
183
184             /* Update vectorial force */
185             fix0            += tx;
186             fiy0            += ty;
187             fiz0            += tz;
188             f[j_coord_offset+DIM*0+XX] -= tx;
189             f[j_coord_offset+DIM*0+YY] -= ty;
190             f[j_coord_offset+DIM*0+ZZ] -= tz;
191
192             }
193
194             /* Inner loop uses 42 flops */
195         }
196         /* End of innermost loop */
197
198         tx = ty = tz = 0;
199         f[i_coord_offset+DIM*0+XX] += fix0;
200         f[i_coord_offset+DIM*0+YY] += fiy0;
201         f[i_coord_offset+DIM*0+ZZ] += fiz0;
202         tx                         += fix0;
203         ty                         += fiy0;
204         tz                         += fiz0;
205         fshift[i_shift_offset+XX]  += tx;
206         fshift[i_shift_offset+YY]  += ty;
207         fshift[i_shift_offset+ZZ]  += tz;
208
209         ggid                        = gid[iidx];
210         /* Update potential energies */
211         kernel_data->energygrp_elec[ggid] += velecsum;
212
213         /* Increment number of inner iterations */
214         inneriter                  += j_index_end - j_index_start;
215
216         /* Outer loop uses 14 flops */
217     }
218
219     /* Increment number of outer iterations */
220     outeriter        += nri;
221
222     /* Update outer/inner flops */
223
224     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VF,outeriter*14 + inneriter*42);
225 }
226 /*
227  * Gromacs nonbonded kernel:   nb_kernel_ElecEwSh_VdwNone_GeomP1P1_F_c
228  * Electrostatics interaction: Ewald
229  * VdW interaction:            None
230  * Geometry:                   Particle-Particle
231  * Calculate force/pot:        Force
232  */
233 void
234 nb_kernel_ElecEwSh_VdwNone_GeomP1P1_F_c
235                     (t_nblist * gmx_restrict                nlist,
236                      rvec * gmx_restrict                    xx,
237                      rvec * gmx_restrict                    ff,
238                      t_forcerec * gmx_restrict              fr,
239                      t_mdatoms * gmx_restrict               mdatoms,
240                      nb_kernel_data_t * gmx_restrict        kernel_data,
241                      t_nrnb * gmx_restrict                  nrnb)
242 {
243     int              i_shift_offset,i_coord_offset,j_coord_offset;
244     int              j_index_start,j_index_end;
245     int              nri,inr,ggid,iidx,jidx,jnr,outeriter,inneriter;
246     real             shX,shY,shZ,tx,ty,tz,fscal,rcutoff,rcutoff2;
247     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
248     real             *shiftvec,*fshift,*x,*f;
249     int              vdwioffset0;
250     real             ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
251     int              vdwjidx0;
252     real             jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
253     real             dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00,cexp1_00,cexp2_00;
254     real             velec,felec,velecsum,facel,crf,krf,krf2;
255     real             *charge;
256     int              ewitab;
257     real             ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace;
258     real             *ewtab;
259
260     x                = xx[0];
261     f                = ff[0];
262
263     nri              = nlist->nri;
264     iinr             = nlist->iinr;
265     jindex           = nlist->jindex;
266     jjnr             = nlist->jjnr;
267     shiftidx         = nlist->shift;
268     gid              = nlist->gid;
269     shiftvec         = fr->shift_vec[0];
270     fshift           = fr->fshift[0];
271     facel            = fr->epsfac;
272     charge           = mdatoms->chargeA;
273
274     sh_ewald         = fr->ic->sh_ewald;
275     ewtab            = fr->ic->tabq_coul_F;
276     ewtabscale       = fr->ic->tabq_scale;
277     ewtabhalfspace   = 0.5/ewtabscale;
278
279     /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
280     rcutoff          = fr->rcoulomb;
281     rcutoff2         = rcutoff*rcutoff;
282
283     outeriter        = 0;
284     inneriter        = 0;
285
286     /* Start outer loop over neighborlists */
287     for(iidx=0; iidx<nri; iidx++)
288     {
289         /* Load shift vector for this list */
290         i_shift_offset   = DIM*shiftidx[iidx];
291         shX              = shiftvec[i_shift_offset+XX];
292         shY              = shiftvec[i_shift_offset+YY];
293         shZ              = shiftvec[i_shift_offset+ZZ];
294
295         /* Load limits for loop over neighbors */
296         j_index_start    = jindex[iidx];
297         j_index_end      = jindex[iidx+1];
298
299         /* Get outer coordinate index */
300         inr              = iinr[iidx];
301         i_coord_offset   = DIM*inr;
302
303         /* Load i particle coords and add shift vector */
304         ix0              = shX + x[i_coord_offset+DIM*0+XX];
305         iy0              = shY + x[i_coord_offset+DIM*0+YY];
306         iz0              = shZ + x[i_coord_offset+DIM*0+ZZ];
307
308         fix0             = 0.0;
309         fiy0             = 0.0;
310         fiz0             = 0.0;
311
312         /* Load parameters for i particles */
313         iq0              = facel*charge[inr+0];
314
315         /* Start inner kernel loop */
316         for(jidx=j_index_start; jidx<j_index_end; jidx++)
317         {
318             /* Get j neighbor index, and coordinate index */
319             jnr              = jjnr[jidx];
320             j_coord_offset   = DIM*jnr;
321
322             /* load j atom coordinates */
323             jx0              = x[j_coord_offset+DIM*0+XX];
324             jy0              = x[j_coord_offset+DIM*0+YY];
325             jz0              = x[j_coord_offset+DIM*0+ZZ];
326
327             /* Calculate displacement vector */
328             dx00             = ix0 - jx0;
329             dy00             = iy0 - jy0;
330             dz00             = iz0 - jz0;
331
332             /* Calculate squared distance and things based on it */
333             rsq00            = dx00*dx00+dy00*dy00+dz00*dz00;
334
335             rinv00           = gmx_invsqrt(rsq00);
336
337             rinvsq00         = rinv00*rinv00;
338
339             /* Load parameters for j particles */
340             jq0              = charge[jnr+0];
341
342             /**************************
343              * CALCULATE INTERACTIONS *
344              **************************/
345
346             if (rsq00<rcutoff2)
347             {
348
349             r00              = rsq00*rinv00;
350
351             qq00             = iq0*jq0;
352
353             /* EWALD ELECTROSTATICS */
354
355             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
356             ewrt             = r00*ewtabscale;
357             ewitab           = ewrt;
358             eweps            = ewrt-ewitab;
359             felec            = (1.0-eweps)*ewtab[ewitab]+eweps*ewtab[ewitab+1];
360             felec            = qq00*rinv00*(rinvsq00-felec);
361
362             fscal            = felec;
363
364             /* Calculate temporary vectorial force */
365             tx               = fscal*dx00;
366             ty               = fscal*dy00;
367             tz               = fscal*dz00;
368
369             /* Update vectorial force */
370             fix0            += tx;
371             fiy0            += ty;
372             fiz0            += tz;
373             f[j_coord_offset+DIM*0+XX] -= tx;
374             f[j_coord_offset+DIM*0+YY] -= ty;
375             f[j_coord_offset+DIM*0+ZZ] -= tz;
376
377             }
378
379             /* Inner loop uses 34 flops */
380         }
381         /* End of innermost loop */
382
383         tx = ty = tz = 0;
384         f[i_coord_offset+DIM*0+XX] += fix0;
385         f[i_coord_offset+DIM*0+YY] += fiy0;
386         f[i_coord_offset+DIM*0+ZZ] += fiz0;
387         tx                         += fix0;
388         ty                         += fiy0;
389         tz                         += fiz0;
390         fshift[i_shift_offset+XX]  += tx;
391         fshift[i_shift_offset+YY]  += ty;
392         fshift[i_shift_offset+ZZ]  += tz;
393
394         /* Increment number of inner iterations */
395         inneriter                  += j_index_end - j_index_start;
396
397         /* Outer loop uses 13 flops */
398     }
399
400     /* Increment number of outer iterations */
401     outeriter        += nri;
402
403     /* Update outer/inner flops */
404
405     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_F,outeriter*13 + inneriter*34);
406 }