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