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