c6d605ecbd4779d218ee6765a249f346a1291b75
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_c / nb_kernel_ElecEwSh_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_ElecEwSh_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_ElecEwSh_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     /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
101     rcutoff          = fr->rcoulomb;
102     rcutoff2         = rcutoff*rcutoff;
103
104     outeriter        = 0;
105     inneriter        = 0;
106
107     /* Start outer loop over neighborlists */
108     for(iidx=0; iidx<nri; iidx++)
109     {
110         /* Load shift vector for this list */
111         i_shift_offset   = DIM*shiftidx[iidx];
112         shX              = shiftvec[i_shift_offset+XX];
113         shY              = shiftvec[i_shift_offset+YY];
114         shZ              = shiftvec[i_shift_offset+ZZ];
115
116         /* Load limits for loop over neighbors */
117         j_index_start    = jindex[iidx];
118         j_index_end      = jindex[iidx+1];
119
120         /* Get outer coordinate index */
121         inr              = iinr[iidx];
122         i_coord_offset   = DIM*inr;
123
124         /* Load i particle coords and add shift vector */
125         ix0              = shX + x[i_coord_offset+DIM*0+XX];
126         iy0              = shY + x[i_coord_offset+DIM*0+YY];
127         iz0              = shZ + x[i_coord_offset+DIM*0+ZZ];
128
129         fix0             = 0.0;
130         fiy0             = 0.0;
131         fiz0             = 0.0;
132
133         /* Load parameters for i particles */
134         iq0              = facel*charge[inr+0];
135
136         /* Reset potential sums */
137         velecsum         = 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
166             /**************************
167              * CALCULATE INTERACTIONS *
168              **************************/
169
170             if (rsq00<rcutoff2)
171             {
172
173             r00              = rsq00*rinv00;
174
175             qq00             = iq0*jq0;
176
177             /* EWALD ELECTROSTATICS */
178
179             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
180             ewrt             = r00*ewtabscale;
181             ewitab           = ewrt;
182             eweps            = ewrt-ewitab;
183             ewitab           = 4*ewitab;
184             felec            = ewtab[ewitab]+eweps*ewtab[ewitab+1];
185             velec            = qq00*((rinv00-sh_ewald)-(ewtab[ewitab+2]-ewtabhalfspace*eweps*(ewtab[ewitab]+felec)));
186             felec            = qq00*rinv00*(rinvsq00-felec);
187
188             /* Update potential sums from outer loop */
189             velecsum        += velec;
190
191             fscal            = felec;
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             }
207
208             /* Inner loop uses 42 flops */
209         }
210         /* End of innermost loop */
211
212         tx = ty = tz = 0;
213         f[i_coord_offset+DIM*0+XX] += fix0;
214         f[i_coord_offset+DIM*0+YY] += fiy0;
215         f[i_coord_offset+DIM*0+ZZ] += fiz0;
216         tx                         += fix0;
217         ty                         += fiy0;
218         tz                         += fiz0;
219         fshift[i_shift_offset+XX]  += tx;
220         fshift[i_shift_offset+YY]  += ty;
221         fshift[i_shift_offset+ZZ]  += tz;
222
223         ggid                        = gid[iidx];
224         /* Update potential energies */
225         kernel_data->energygrp_elec[ggid] += velecsum;
226
227         /* Increment number of inner iterations */
228         inneriter                  += j_index_end - j_index_start;
229
230         /* Outer loop uses 14 flops */
231     }
232
233     /* Increment number of outer iterations */
234     outeriter        += nri;
235
236     /* Update outer/inner flops */
237
238     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VF,outeriter*14 + inneriter*42);
239 }
240 /*
241  * Gromacs nonbonded kernel:   nb_kernel_ElecEwSh_VdwNone_GeomP1P1_F_c
242  * Electrostatics interaction: Ewald
243  * VdW interaction:            None
244  * Geometry:                   Particle-Particle
245  * Calculate force/pot:        Force
246  */
247 void
248 nb_kernel_ElecEwSh_VdwNone_GeomP1P1_F_c
249                     (t_nblist                    * gmx_restrict       nlist,
250                      rvec                        * gmx_restrict          xx,
251                      rvec                        * gmx_restrict          ff,
252                      t_forcerec                  * gmx_restrict          fr,
253                      t_mdatoms                   * gmx_restrict     mdatoms,
254                      nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
255                      t_nrnb                      * gmx_restrict        nrnb)
256 {
257     int              i_shift_offset,i_coord_offset,j_coord_offset;
258     int              j_index_start,j_index_end;
259     int              nri,inr,ggid,iidx,jidx,jnr,outeriter,inneriter;
260     real             shX,shY,shZ,tx,ty,tz,fscal,rcutoff,rcutoff2;
261     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
262     real             *shiftvec,*fshift,*x,*f;
263     int              vdwioffset0;
264     real             ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
265     int              vdwjidx0;
266     real             jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
267     real             dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00,cexp1_00,cexp2_00;
268     real             velec,felec,velecsum,facel,crf,krf,krf2;
269     real             *charge;
270     int              ewitab;
271     real             ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace;
272     real             *ewtab;
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
288     sh_ewald         = fr->ic->sh_ewald;
289     ewtab            = fr->ic->tabq_coul_F;
290     ewtabscale       = fr->ic->tabq_scale;
291     ewtabhalfspace   = 0.5/ewtabscale;
292
293     /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
294     rcutoff          = fr->rcoulomb;
295     rcutoff2         = rcutoff*rcutoff;
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
329         /* Start inner kernel loop */
330         for(jidx=j_index_start; jidx<j_index_end; jidx++)
331         {
332             /* Get j neighbor index, and coordinate index */
333             jnr              = jjnr[jidx];
334             j_coord_offset   = DIM*jnr;
335
336             /* load j atom coordinates */
337             jx0              = x[j_coord_offset+DIM*0+XX];
338             jy0              = x[j_coord_offset+DIM*0+YY];
339             jz0              = x[j_coord_offset+DIM*0+ZZ];
340
341             /* Calculate displacement vector */
342             dx00             = ix0 - jx0;
343             dy00             = iy0 - jy0;
344             dz00             = iz0 - jz0;
345
346             /* Calculate squared distance and things based on it */
347             rsq00            = dx00*dx00+dy00*dy00+dz00*dz00;
348
349             rinv00           = gmx_invsqrt(rsq00);
350
351             rinvsq00         = rinv00*rinv00;
352
353             /* Load parameters for j particles */
354             jq0              = charge[jnr+0];
355
356             /**************************
357              * CALCULATE INTERACTIONS *
358              **************************/
359
360             if (rsq00<rcutoff2)
361             {
362
363             r00              = rsq00*rinv00;
364
365             qq00             = iq0*jq0;
366
367             /* EWALD ELECTROSTATICS */
368
369             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
370             ewrt             = r00*ewtabscale;
371             ewitab           = ewrt;
372             eweps            = ewrt-ewitab;
373             felec            = (1.0-eweps)*ewtab[ewitab]+eweps*ewtab[ewitab+1];
374             felec            = qq00*rinv00*(rinvsq00-felec);
375
376             fscal            = felec;
377
378             /* Calculate temporary vectorial force */
379             tx               = fscal*dx00;
380             ty               = fscal*dy00;
381             tz               = fscal*dz00;
382
383             /* Update vectorial force */
384             fix0            += tx;
385             fiy0            += ty;
386             fiz0            += tz;
387             f[j_coord_offset+DIM*0+XX] -= tx;
388             f[j_coord_offset+DIM*0+YY] -= ty;
389             f[j_coord_offset+DIM*0+ZZ] -= tz;
390
391             }
392
393             /* Inner loop uses 34 flops */
394         }
395         /* End of innermost loop */
396
397         tx = ty = tz = 0;
398         f[i_coord_offset+DIM*0+XX] += fix0;
399         f[i_coord_offset+DIM*0+YY] += fiy0;
400         f[i_coord_offset+DIM*0+ZZ] += fiz0;
401         tx                         += fix0;
402         ty                         += fiy0;
403         tz                         += fiz0;
404         fshift[i_shift_offset+XX]  += tx;
405         fshift[i_shift_offset+YY]  += ty;
406         fshift[i_shift_offset+ZZ]  += tz;
407
408         /* Increment number of inner iterations */
409         inneriter                  += j_index_end - j_index_start;
410
411         /* Outer loop uses 13 flops */
412     }
413
414     /* Increment number of outer iterations */
415     outeriter        += nri;
416
417     /* Update outer/inner flops */
418
419     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_F,outeriter*13 + inneriter*34);
420 }