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