Merge release-5-0 into master
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_c / nb_kernel_ElecCoul_VdwNone_GeomW4P1_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_ElecCoul_VdwNone_GeomW4P1_VF_c
51  * Electrostatics interaction: Coulomb
52  * VdW interaction:            None
53  * Geometry:                   Water4-Particle
54  * Calculate force/pot:        PotentialAndForce
55  */
56 void
57 nb_kernel_ElecCoul_VdwNone_GeomW4P1_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              vdwioffset1;
73     real             ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
74     int              vdwioffset2;
75     real             ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
76     int              vdwioffset3;
77     real             ix3,iy3,iz3,fix3,fiy3,fiz3,iq3,isai3;
78     int              vdwjidx0;
79     real             jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
80     real             dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10,cexp1_10,cexp2_10;
81     real             dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20,cexp1_20,cexp2_20;
82     real             dx30,dy30,dz30,rsq30,rinv30,rinvsq30,r30,qq30,c6_30,c12_30,cexp1_30,cexp2_30;
83     real             velec,felec,velecsum,facel,crf,krf,krf2;
84     real             *charge;
85
86     x                = xx[0];
87     f                = ff[0];
88
89     nri              = nlist->nri;
90     iinr             = nlist->iinr;
91     jindex           = nlist->jindex;
92     jjnr             = nlist->jjnr;
93     shiftidx         = nlist->shift;
94     gid              = nlist->gid;
95     shiftvec         = fr->shift_vec[0];
96     fshift           = fr->fshift[0];
97     facel            = fr->epsfac;
98     charge           = mdatoms->chargeA;
99
100     /* Setup water-specific parameters */
101     inr              = nlist->iinr[0];
102     iq1              = facel*charge[inr+1];
103     iq2              = facel*charge[inr+2];
104     iq3              = facel*charge[inr+3];
105
106     outeriter        = 0;
107     inneriter        = 0;
108
109     /* Start outer loop over neighborlists */
110     for(iidx=0; iidx<nri; iidx++)
111     {
112         /* Load shift vector for this list */
113         i_shift_offset   = DIM*shiftidx[iidx];
114         shX              = shiftvec[i_shift_offset+XX];
115         shY              = shiftvec[i_shift_offset+YY];
116         shZ              = shiftvec[i_shift_offset+ZZ];
117
118         /* Load limits for loop over neighbors */
119         j_index_start    = jindex[iidx];
120         j_index_end      = jindex[iidx+1];
121
122         /* Get outer coordinate index */
123         inr              = iinr[iidx];
124         i_coord_offset   = DIM*inr;
125
126         /* Load i particle coords and add shift vector */
127         ix1              = shX + x[i_coord_offset+DIM*1+XX];
128         iy1              = shY + x[i_coord_offset+DIM*1+YY];
129         iz1              = shZ + x[i_coord_offset+DIM*1+ZZ];
130         ix2              = shX + x[i_coord_offset+DIM*2+XX];
131         iy2              = shY + x[i_coord_offset+DIM*2+YY];
132         iz2              = shZ + x[i_coord_offset+DIM*2+ZZ];
133         ix3              = shX + x[i_coord_offset+DIM*3+XX];
134         iy3              = shY + x[i_coord_offset+DIM*3+YY];
135         iz3              = shZ + x[i_coord_offset+DIM*3+ZZ];
136
137         fix1             = 0.0;
138         fiy1             = 0.0;
139         fiz1             = 0.0;
140         fix2             = 0.0;
141         fiy2             = 0.0;
142         fiz2             = 0.0;
143         fix3             = 0.0;
144         fiy3             = 0.0;
145         fiz3             = 0.0;
146
147         /* Reset potential sums */
148         velecsum         = 0.0;
149
150         /* Start inner kernel loop */
151         for(jidx=j_index_start; jidx<j_index_end; jidx++)
152         {
153             /* Get j neighbor index, and coordinate index */
154             jnr              = jjnr[jidx];
155             j_coord_offset   = DIM*jnr;
156
157             /* load j atom coordinates */
158             jx0              = x[j_coord_offset+DIM*0+XX];
159             jy0              = x[j_coord_offset+DIM*0+YY];
160             jz0              = x[j_coord_offset+DIM*0+ZZ];
161
162             /* Calculate displacement vector */
163             dx10             = ix1 - jx0;
164             dy10             = iy1 - jy0;
165             dz10             = iz1 - jz0;
166             dx20             = ix2 - jx0;
167             dy20             = iy2 - jy0;
168             dz20             = iz2 - jz0;
169             dx30             = ix3 - jx0;
170             dy30             = iy3 - jy0;
171             dz30             = iz3 - jz0;
172
173             /* Calculate squared distance and things based on it */
174             rsq10            = dx10*dx10+dy10*dy10+dz10*dz10;
175             rsq20            = dx20*dx20+dy20*dy20+dz20*dz20;
176             rsq30            = dx30*dx30+dy30*dy30+dz30*dz30;
177
178             rinv10           = gmx_invsqrt(rsq10);
179             rinv20           = gmx_invsqrt(rsq20);
180             rinv30           = gmx_invsqrt(rsq30);
181
182             rinvsq10         = rinv10*rinv10;
183             rinvsq20         = rinv20*rinv20;
184             rinvsq30         = rinv30*rinv30;
185
186             /* Load parameters for j particles */
187             jq0              = charge[jnr+0];
188
189             /**************************
190              * CALCULATE INTERACTIONS *
191              **************************/
192
193             qq10             = iq1*jq0;
194
195             /* COULOMB ELECTROSTATICS */
196             velec            = qq10*rinv10;
197             felec            = velec*rinvsq10;
198
199             /* Update potential sums from outer loop */
200             velecsum        += velec;
201
202             fscal            = felec;
203
204             /* Calculate temporary vectorial force */
205             tx               = fscal*dx10;
206             ty               = fscal*dy10;
207             tz               = fscal*dz10;
208
209             /* Update vectorial force */
210             fix1            += tx;
211             fiy1            += ty;
212             fiz1            += tz;
213             f[j_coord_offset+DIM*0+XX] -= tx;
214             f[j_coord_offset+DIM*0+YY] -= ty;
215             f[j_coord_offset+DIM*0+ZZ] -= tz;
216
217             /**************************
218              * CALCULATE INTERACTIONS *
219              **************************/
220
221             qq20             = iq2*jq0;
222
223             /* COULOMB ELECTROSTATICS */
224             velec            = qq20*rinv20;
225             felec            = velec*rinvsq20;
226
227             /* Update potential sums from outer loop */
228             velecsum        += velec;
229
230             fscal            = felec;
231
232             /* Calculate temporary vectorial force */
233             tx               = fscal*dx20;
234             ty               = fscal*dy20;
235             tz               = fscal*dz20;
236
237             /* Update vectorial force */
238             fix2            += tx;
239             fiy2            += ty;
240             fiz2            += tz;
241             f[j_coord_offset+DIM*0+XX] -= tx;
242             f[j_coord_offset+DIM*0+YY] -= ty;
243             f[j_coord_offset+DIM*0+ZZ] -= tz;
244
245             /**************************
246              * CALCULATE INTERACTIONS *
247              **************************/
248
249             qq30             = iq3*jq0;
250
251             /* COULOMB ELECTROSTATICS */
252             velec            = qq30*rinv30;
253             felec            = velec*rinvsq30;
254
255             /* Update potential sums from outer loop */
256             velecsum        += velec;
257
258             fscal            = felec;
259
260             /* Calculate temporary vectorial force */
261             tx               = fscal*dx30;
262             ty               = fscal*dy30;
263             tz               = fscal*dz30;
264
265             /* Update vectorial force */
266             fix3            += tx;
267             fiy3            += ty;
268             fiz3            += tz;
269             f[j_coord_offset+DIM*0+XX] -= tx;
270             f[j_coord_offset+DIM*0+YY] -= ty;
271             f[j_coord_offset+DIM*0+ZZ] -= tz;
272
273             /* Inner loop uses 84 flops */
274         }
275         /* End of innermost loop */
276
277         tx = ty = tz = 0;
278         f[i_coord_offset+DIM*1+XX] += fix1;
279         f[i_coord_offset+DIM*1+YY] += fiy1;
280         f[i_coord_offset+DIM*1+ZZ] += fiz1;
281         tx                         += fix1;
282         ty                         += fiy1;
283         tz                         += fiz1;
284         f[i_coord_offset+DIM*2+XX] += fix2;
285         f[i_coord_offset+DIM*2+YY] += fiy2;
286         f[i_coord_offset+DIM*2+ZZ] += fiz2;
287         tx                         += fix2;
288         ty                         += fiy2;
289         tz                         += fiz2;
290         f[i_coord_offset+DIM*3+XX] += fix3;
291         f[i_coord_offset+DIM*3+YY] += fiy3;
292         f[i_coord_offset+DIM*3+ZZ] += fiz3;
293         tx                         += fix3;
294         ty                         += fiy3;
295         tz                         += fiz3;
296         fshift[i_shift_offset+XX]  += tx;
297         fshift[i_shift_offset+YY]  += ty;
298         fshift[i_shift_offset+ZZ]  += tz;
299
300         ggid                        = gid[iidx];
301         /* Update potential energies */
302         kernel_data->energygrp_elec[ggid] += velecsum;
303
304         /* Increment number of inner iterations */
305         inneriter                  += j_index_end - j_index_start;
306
307         /* Outer loop uses 31 flops */
308     }
309
310     /* Increment number of outer iterations */
311     outeriter        += nri;
312
313     /* Update outer/inner flops */
314
315     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_W4_VF,outeriter*31 + inneriter*84);
316 }
317 /*
318  * Gromacs nonbonded kernel:   nb_kernel_ElecCoul_VdwNone_GeomW4P1_F_c
319  * Electrostatics interaction: Coulomb
320  * VdW interaction:            None
321  * Geometry:                   Water4-Particle
322  * Calculate force/pot:        Force
323  */
324 void
325 nb_kernel_ElecCoul_VdwNone_GeomW4P1_F_c
326                     (t_nblist                    * gmx_restrict       nlist,
327                      rvec                        * gmx_restrict          xx,
328                      rvec                        * gmx_restrict          ff,
329                      t_forcerec                  * gmx_restrict          fr,
330                      t_mdatoms                   * gmx_restrict     mdatoms,
331                      nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
332                      t_nrnb                      * gmx_restrict        nrnb)
333 {
334     int              i_shift_offset,i_coord_offset,j_coord_offset;
335     int              j_index_start,j_index_end;
336     int              nri,inr,ggid,iidx,jidx,jnr,outeriter,inneriter;
337     real             shX,shY,shZ,tx,ty,tz,fscal,rcutoff,rcutoff2;
338     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
339     real             *shiftvec,*fshift,*x,*f;
340     int              vdwioffset1;
341     real             ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
342     int              vdwioffset2;
343     real             ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
344     int              vdwioffset3;
345     real             ix3,iy3,iz3,fix3,fiy3,fiz3,iq3,isai3;
346     int              vdwjidx0;
347     real             jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
348     real             dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10,cexp1_10,cexp2_10;
349     real             dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20,cexp1_20,cexp2_20;
350     real             dx30,dy30,dz30,rsq30,rinv30,rinvsq30,r30,qq30,c6_30,c12_30,cexp1_30,cexp2_30;
351     real             velec,felec,velecsum,facel,crf,krf,krf2;
352     real             *charge;
353
354     x                = xx[0];
355     f                = ff[0];
356
357     nri              = nlist->nri;
358     iinr             = nlist->iinr;
359     jindex           = nlist->jindex;
360     jjnr             = nlist->jjnr;
361     shiftidx         = nlist->shift;
362     gid              = nlist->gid;
363     shiftvec         = fr->shift_vec[0];
364     fshift           = fr->fshift[0];
365     facel            = fr->epsfac;
366     charge           = mdatoms->chargeA;
367
368     /* Setup water-specific parameters */
369     inr              = nlist->iinr[0];
370     iq1              = facel*charge[inr+1];
371     iq2              = facel*charge[inr+2];
372     iq3              = facel*charge[inr+3];
373
374     outeriter        = 0;
375     inneriter        = 0;
376
377     /* Start outer loop over neighborlists */
378     for(iidx=0; iidx<nri; iidx++)
379     {
380         /* Load shift vector for this list */
381         i_shift_offset   = DIM*shiftidx[iidx];
382         shX              = shiftvec[i_shift_offset+XX];
383         shY              = shiftvec[i_shift_offset+YY];
384         shZ              = shiftvec[i_shift_offset+ZZ];
385
386         /* Load limits for loop over neighbors */
387         j_index_start    = jindex[iidx];
388         j_index_end      = jindex[iidx+1];
389
390         /* Get outer coordinate index */
391         inr              = iinr[iidx];
392         i_coord_offset   = DIM*inr;
393
394         /* Load i particle coords and add shift vector */
395         ix1              = shX + x[i_coord_offset+DIM*1+XX];
396         iy1              = shY + x[i_coord_offset+DIM*1+YY];
397         iz1              = shZ + x[i_coord_offset+DIM*1+ZZ];
398         ix2              = shX + x[i_coord_offset+DIM*2+XX];
399         iy2              = shY + x[i_coord_offset+DIM*2+YY];
400         iz2              = shZ + x[i_coord_offset+DIM*2+ZZ];
401         ix3              = shX + x[i_coord_offset+DIM*3+XX];
402         iy3              = shY + x[i_coord_offset+DIM*3+YY];
403         iz3              = shZ + x[i_coord_offset+DIM*3+ZZ];
404
405         fix1             = 0.0;
406         fiy1             = 0.0;
407         fiz1             = 0.0;
408         fix2             = 0.0;
409         fiy2             = 0.0;
410         fiz2             = 0.0;
411         fix3             = 0.0;
412         fiy3             = 0.0;
413         fiz3             = 0.0;
414
415         /* Start inner kernel loop */
416         for(jidx=j_index_start; jidx<j_index_end; jidx++)
417         {
418             /* Get j neighbor index, and coordinate index */
419             jnr              = jjnr[jidx];
420             j_coord_offset   = DIM*jnr;
421
422             /* load j atom coordinates */
423             jx0              = x[j_coord_offset+DIM*0+XX];
424             jy0              = x[j_coord_offset+DIM*0+YY];
425             jz0              = x[j_coord_offset+DIM*0+ZZ];
426
427             /* Calculate displacement vector */
428             dx10             = ix1 - jx0;
429             dy10             = iy1 - jy0;
430             dz10             = iz1 - jz0;
431             dx20             = ix2 - jx0;
432             dy20             = iy2 - jy0;
433             dz20             = iz2 - jz0;
434             dx30             = ix3 - jx0;
435             dy30             = iy3 - jy0;
436             dz30             = iz3 - jz0;
437
438             /* Calculate squared distance and things based on it */
439             rsq10            = dx10*dx10+dy10*dy10+dz10*dz10;
440             rsq20            = dx20*dx20+dy20*dy20+dz20*dz20;
441             rsq30            = dx30*dx30+dy30*dy30+dz30*dz30;
442
443             rinv10           = gmx_invsqrt(rsq10);
444             rinv20           = gmx_invsqrt(rsq20);
445             rinv30           = gmx_invsqrt(rsq30);
446
447             rinvsq10         = rinv10*rinv10;
448             rinvsq20         = rinv20*rinv20;
449             rinvsq30         = rinv30*rinv30;
450
451             /* Load parameters for j particles */
452             jq0              = charge[jnr+0];
453
454             /**************************
455              * CALCULATE INTERACTIONS *
456              **************************/
457
458             qq10             = iq1*jq0;
459
460             /* COULOMB ELECTROSTATICS */
461             velec            = qq10*rinv10;
462             felec            = velec*rinvsq10;
463
464             fscal            = felec;
465
466             /* Calculate temporary vectorial force */
467             tx               = fscal*dx10;
468             ty               = fscal*dy10;
469             tz               = fscal*dz10;
470
471             /* Update vectorial force */
472             fix1            += tx;
473             fiy1            += ty;
474             fiz1            += tz;
475             f[j_coord_offset+DIM*0+XX] -= tx;
476             f[j_coord_offset+DIM*0+YY] -= ty;
477             f[j_coord_offset+DIM*0+ZZ] -= tz;
478
479             /**************************
480              * CALCULATE INTERACTIONS *
481              **************************/
482
483             qq20             = iq2*jq0;
484
485             /* COULOMB ELECTROSTATICS */
486             velec            = qq20*rinv20;
487             felec            = velec*rinvsq20;
488
489             fscal            = felec;
490
491             /* Calculate temporary vectorial force */
492             tx               = fscal*dx20;
493             ty               = fscal*dy20;
494             tz               = fscal*dz20;
495
496             /* Update vectorial force */
497             fix2            += tx;
498             fiy2            += ty;
499             fiz2            += tz;
500             f[j_coord_offset+DIM*0+XX] -= tx;
501             f[j_coord_offset+DIM*0+YY] -= ty;
502             f[j_coord_offset+DIM*0+ZZ] -= tz;
503
504             /**************************
505              * CALCULATE INTERACTIONS *
506              **************************/
507
508             qq30             = iq3*jq0;
509
510             /* COULOMB ELECTROSTATICS */
511             velec            = qq30*rinv30;
512             felec            = velec*rinvsq30;
513
514             fscal            = felec;
515
516             /* Calculate temporary vectorial force */
517             tx               = fscal*dx30;
518             ty               = fscal*dy30;
519             tz               = fscal*dz30;
520
521             /* Update vectorial force */
522             fix3            += tx;
523             fiy3            += ty;
524             fiz3            += tz;
525             f[j_coord_offset+DIM*0+XX] -= tx;
526             f[j_coord_offset+DIM*0+YY] -= ty;
527             f[j_coord_offset+DIM*0+ZZ] -= tz;
528
529             /* Inner loop uses 81 flops */
530         }
531         /* End of innermost loop */
532
533         tx = ty = tz = 0;
534         f[i_coord_offset+DIM*1+XX] += fix1;
535         f[i_coord_offset+DIM*1+YY] += fiy1;
536         f[i_coord_offset+DIM*1+ZZ] += fiz1;
537         tx                         += fix1;
538         ty                         += fiy1;
539         tz                         += fiz1;
540         f[i_coord_offset+DIM*2+XX] += fix2;
541         f[i_coord_offset+DIM*2+YY] += fiy2;
542         f[i_coord_offset+DIM*2+ZZ] += fiz2;
543         tx                         += fix2;
544         ty                         += fiy2;
545         tz                         += fiz2;
546         f[i_coord_offset+DIM*3+XX] += fix3;
547         f[i_coord_offset+DIM*3+YY] += fiy3;
548         f[i_coord_offset+DIM*3+ZZ] += fiz3;
549         tx                         += fix3;
550         ty                         += fiy3;
551         tz                         += fiz3;
552         fshift[i_shift_offset+XX]  += tx;
553         fshift[i_shift_offset+YY]  += ty;
554         fshift[i_shift_offset+ZZ]  += tz;
555
556         /* Increment number of inner iterations */
557         inneriter                  += j_index_end - j_index_start;
558
559         /* Outer loop uses 30 flops */
560     }
561
562     /* Increment number of outer iterations */
563     outeriter        += nri;
564
565     /* Update outer/inner flops */
566
567     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_W4_F,outeriter*30 + inneriter*81);
568 }