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