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