Introduce gmxpre.h for truly global definitions
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_c / nb_kernel_ElecCoul_VdwNone_GeomW3W3_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_GeomW3W3_VF_c
51  * Electrostatics interaction: Coulomb
52  * VdW interaction:            None
53  * Geometry:                   Water3-Water3
54  * Calculate force/pot:        PotentialAndForce
55  */
56 void
57 nb_kernel_ElecCoul_VdwNone_GeomW3W3_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     int              vdwjidx1;
81     real             jx1,jy1,jz1,fjx1,fjy1,fjz1,jq1,isaj1;
82     int              vdwjidx2;
83     real             jx2,jy2,jz2,fjx2,fjy2,fjz2,jq2,isaj2;
84     real             dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00,cexp1_00,cexp2_00;
85     real             dx01,dy01,dz01,rsq01,rinv01,rinvsq01,r01,qq01,c6_01,c12_01,cexp1_01,cexp2_01;
86     real             dx02,dy02,dz02,rsq02,rinv02,rinvsq02,r02,qq02,c6_02,c12_02,cexp1_02,cexp2_02;
87     real             dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10,cexp1_10,cexp2_10;
88     real             dx11,dy11,dz11,rsq11,rinv11,rinvsq11,r11,qq11,c6_11,c12_11,cexp1_11,cexp2_11;
89     real             dx12,dy12,dz12,rsq12,rinv12,rinvsq12,r12,qq12,c6_12,c12_12,cexp1_12,cexp2_12;
90     real             dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20,cexp1_20,cexp2_20;
91     real             dx21,dy21,dz21,rsq21,rinv21,rinvsq21,r21,qq21,c6_21,c12_21,cexp1_21,cexp2_21;
92     real             dx22,dy22,dz22,rsq22,rinv22,rinvsq22,r22,qq22,c6_22,c12_22,cexp1_22,cexp2_22;
93     real             velec,felec,velecsum,facel,crf,krf,krf2;
94     real             *charge;
95
96     x                = xx[0];
97     f                = ff[0];
98
99     nri              = nlist->nri;
100     iinr             = nlist->iinr;
101     jindex           = nlist->jindex;
102     jjnr             = nlist->jjnr;
103     shiftidx         = nlist->shift;
104     gid              = nlist->gid;
105     shiftvec         = fr->shift_vec[0];
106     fshift           = fr->fshift[0];
107     facel            = fr->epsfac;
108     charge           = mdatoms->chargeA;
109
110     /* Setup water-specific parameters */
111     inr              = nlist->iinr[0];
112     iq0              = facel*charge[inr+0];
113     iq1              = facel*charge[inr+1];
114     iq2              = facel*charge[inr+2];
115
116     jq0              = charge[inr+0];
117     jq1              = charge[inr+1];
118     jq2              = charge[inr+2];
119     qq00             = iq0*jq0;
120     qq01             = iq0*jq1;
121     qq02             = iq0*jq2;
122     qq10             = iq1*jq0;
123     qq11             = iq1*jq1;
124     qq12             = iq1*jq2;
125     qq20             = iq2*jq0;
126     qq21             = iq2*jq1;
127     qq22             = iq2*jq2;
128
129     outeriter        = 0;
130     inneriter        = 0;
131
132     /* Start outer loop over neighborlists */
133     for(iidx=0; iidx<nri; iidx++)
134     {
135         /* Load shift vector for this list */
136         i_shift_offset   = DIM*shiftidx[iidx];
137         shX              = shiftvec[i_shift_offset+XX];
138         shY              = shiftvec[i_shift_offset+YY];
139         shZ              = shiftvec[i_shift_offset+ZZ];
140
141         /* Load limits for loop over neighbors */
142         j_index_start    = jindex[iidx];
143         j_index_end      = jindex[iidx+1];
144
145         /* Get outer coordinate index */
146         inr              = iinr[iidx];
147         i_coord_offset   = DIM*inr;
148
149         /* Load i particle coords and add shift vector */
150         ix0              = shX + x[i_coord_offset+DIM*0+XX];
151         iy0              = shY + x[i_coord_offset+DIM*0+YY];
152         iz0              = shZ + x[i_coord_offset+DIM*0+ZZ];
153         ix1              = shX + x[i_coord_offset+DIM*1+XX];
154         iy1              = shY + x[i_coord_offset+DIM*1+YY];
155         iz1              = shZ + x[i_coord_offset+DIM*1+ZZ];
156         ix2              = shX + x[i_coord_offset+DIM*2+XX];
157         iy2              = shY + x[i_coord_offset+DIM*2+YY];
158         iz2              = shZ + x[i_coord_offset+DIM*2+ZZ];
159
160         fix0             = 0.0;
161         fiy0             = 0.0;
162         fiz0             = 0.0;
163         fix1             = 0.0;
164         fiy1             = 0.0;
165         fiz1             = 0.0;
166         fix2             = 0.0;
167         fiy2             = 0.0;
168         fiz2             = 0.0;
169
170         /* Reset potential sums */
171         velecsum         = 0.0;
172
173         /* Start inner kernel loop */
174         for(jidx=j_index_start; jidx<j_index_end; jidx++)
175         {
176             /* Get j neighbor index, and coordinate index */
177             jnr              = jjnr[jidx];
178             j_coord_offset   = DIM*jnr;
179
180             /* load j atom coordinates */
181             jx0              = x[j_coord_offset+DIM*0+XX];
182             jy0              = x[j_coord_offset+DIM*0+YY];
183             jz0              = x[j_coord_offset+DIM*0+ZZ];
184             jx1              = x[j_coord_offset+DIM*1+XX];
185             jy1              = x[j_coord_offset+DIM*1+YY];
186             jz1              = x[j_coord_offset+DIM*1+ZZ];
187             jx2              = x[j_coord_offset+DIM*2+XX];
188             jy2              = x[j_coord_offset+DIM*2+YY];
189             jz2              = x[j_coord_offset+DIM*2+ZZ];
190
191             /* Calculate displacement vector */
192             dx00             = ix0 - jx0;
193             dy00             = iy0 - jy0;
194             dz00             = iz0 - jz0;
195             dx01             = ix0 - jx1;
196             dy01             = iy0 - jy1;
197             dz01             = iz0 - jz1;
198             dx02             = ix0 - jx2;
199             dy02             = iy0 - jy2;
200             dz02             = iz0 - jz2;
201             dx10             = ix1 - jx0;
202             dy10             = iy1 - jy0;
203             dz10             = iz1 - jz0;
204             dx11             = ix1 - jx1;
205             dy11             = iy1 - jy1;
206             dz11             = iz1 - jz1;
207             dx12             = ix1 - jx2;
208             dy12             = iy1 - jy2;
209             dz12             = iz1 - jz2;
210             dx20             = ix2 - jx0;
211             dy20             = iy2 - jy0;
212             dz20             = iz2 - jz0;
213             dx21             = ix2 - jx1;
214             dy21             = iy2 - jy1;
215             dz21             = iz2 - jz1;
216             dx22             = ix2 - jx2;
217             dy22             = iy2 - jy2;
218             dz22             = iz2 - jz2;
219
220             /* Calculate squared distance and things based on it */
221             rsq00            = dx00*dx00+dy00*dy00+dz00*dz00;
222             rsq01            = dx01*dx01+dy01*dy01+dz01*dz01;
223             rsq02            = dx02*dx02+dy02*dy02+dz02*dz02;
224             rsq10            = dx10*dx10+dy10*dy10+dz10*dz10;
225             rsq11            = dx11*dx11+dy11*dy11+dz11*dz11;
226             rsq12            = dx12*dx12+dy12*dy12+dz12*dz12;
227             rsq20            = dx20*dx20+dy20*dy20+dz20*dz20;
228             rsq21            = dx21*dx21+dy21*dy21+dz21*dz21;
229             rsq22            = dx22*dx22+dy22*dy22+dz22*dz22;
230
231             rinv00           = gmx_invsqrt(rsq00);
232             rinv01           = gmx_invsqrt(rsq01);
233             rinv02           = gmx_invsqrt(rsq02);
234             rinv10           = gmx_invsqrt(rsq10);
235             rinv11           = gmx_invsqrt(rsq11);
236             rinv12           = gmx_invsqrt(rsq12);
237             rinv20           = gmx_invsqrt(rsq20);
238             rinv21           = gmx_invsqrt(rsq21);
239             rinv22           = gmx_invsqrt(rsq22);
240
241             rinvsq00         = rinv00*rinv00;
242             rinvsq01         = rinv01*rinv01;
243             rinvsq02         = rinv02*rinv02;
244             rinvsq10         = rinv10*rinv10;
245             rinvsq11         = rinv11*rinv11;
246             rinvsq12         = rinv12*rinv12;
247             rinvsq20         = rinv20*rinv20;
248             rinvsq21         = rinv21*rinv21;
249             rinvsq22         = rinv22*rinv22;
250
251             /**************************
252              * CALCULATE INTERACTIONS *
253              **************************/
254
255             /* COULOMB ELECTROSTATICS */
256             velec            = qq00*rinv00;
257             felec            = velec*rinvsq00;
258
259             /* Update potential sums from outer loop */
260             velecsum        += velec;
261
262             fscal            = felec;
263
264             /* Calculate temporary vectorial force */
265             tx               = fscal*dx00;
266             ty               = fscal*dy00;
267             tz               = fscal*dz00;
268
269             /* Update vectorial force */
270             fix0            += tx;
271             fiy0            += ty;
272             fiz0            += tz;
273             f[j_coord_offset+DIM*0+XX] -= tx;
274             f[j_coord_offset+DIM*0+YY] -= ty;
275             f[j_coord_offset+DIM*0+ZZ] -= tz;
276
277             /**************************
278              * CALCULATE INTERACTIONS *
279              **************************/
280
281             /* COULOMB ELECTROSTATICS */
282             velec            = qq01*rinv01;
283             felec            = velec*rinvsq01;
284
285             /* Update potential sums from outer loop */
286             velecsum        += velec;
287
288             fscal            = felec;
289
290             /* Calculate temporary vectorial force */
291             tx               = fscal*dx01;
292             ty               = fscal*dy01;
293             tz               = fscal*dz01;
294
295             /* Update vectorial force */
296             fix0            += tx;
297             fiy0            += ty;
298             fiz0            += tz;
299             f[j_coord_offset+DIM*1+XX] -= tx;
300             f[j_coord_offset+DIM*1+YY] -= ty;
301             f[j_coord_offset+DIM*1+ZZ] -= tz;
302
303             /**************************
304              * CALCULATE INTERACTIONS *
305              **************************/
306
307             /* COULOMB ELECTROSTATICS */
308             velec            = qq02*rinv02;
309             felec            = velec*rinvsq02;
310
311             /* Update potential sums from outer loop */
312             velecsum        += velec;
313
314             fscal            = felec;
315
316             /* Calculate temporary vectorial force */
317             tx               = fscal*dx02;
318             ty               = fscal*dy02;
319             tz               = fscal*dz02;
320
321             /* Update vectorial force */
322             fix0            += tx;
323             fiy0            += ty;
324             fiz0            += tz;
325             f[j_coord_offset+DIM*2+XX] -= tx;
326             f[j_coord_offset+DIM*2+YY] -= ty;
327             f[j_coord_offset+DIM*2+ZZ] -= tz;
328
329             /**************************
330              * CALCULATE INTERACTIONS *
331              **************************/
332
333             /* COULOMB ELECTROSTATICS */
334             velec            = qq10*rinv10;
335             felec            = velec*rinvsq10;
336
337             /* Update potential sums from outer loop */
338             velecsum        += velec;
339
340             fscal            = felec;
341
342             /* Calculate temporary vectorial force */
343             tx               = fscal*dx10;
344             ty               = fscal*dy10;
345             tz               = fscal*dz10;
346
347             /* Update vectorial force */
348             fix1            += tx;
349             fiy1            += ty;
350             fiz1            += tz;
351             f[j_coord_offset+DIM*0+XX] -= tx;
352             f[j_coord_offset+DIM*0+YY] -= ty;
353             f[j_coord_offset+DIM*0+ZZ] -= tz;
354
355             /**************************
356              * CALCULATE INTERACTIONS *
357              **************************/
358
359             /* COULOMB ELECTROSTATICS */
360             velec            = qq11*rinv11;
361             felec            = velec*rinvsq11;
362
363             /* Update potential sums from outer loop */
364             velecsum        += velec;
365
366             fscal            = felec;
367
368             /* Calculate temporary vectorial force */
369             tx               = fscal*dx11;
370             ty               = fscal*dy11;
371             tz               = fscal*dz11;
372
373             /* Update vectorial force */
374             fix1            += tx;
375             fiy1            += ty;
376             fiz1            += tz;
377             f[j_coord_offset+DIM*1+XX] -= tx;
378             f[j_coord_offset+DIM*1+YY] -= ty;
379             f[j_coord_offset+DIM*1+ZZ] -= tz;
380
381             /**************************
382              * CALCULATE INTERACTIONS *
383              **************************/
384
385             /* COULOMB ELECTROSTATICS */
386             velec            = qq12*rinv12;
387             felec            = velec*rinvsq12;
388
389             /* Update potential sums from outer loop */
390             velecsum        += velec;
391
392             fscal            = felec;
393
394             /* Calculate temporary vectorial force */
395             tx               = fscal*dx12;
396             ty               = fscal*dy12;
397             tz               = fscal*dz12;
398
399             /* Update vectorial force */
400             fix1            += tx;
401             fiy1            += ty;
402             fiz1            += tz;
403             f[j_coord_offset+DIM*2+XX] -= tx;
404             f[j_coord_offset+DIM*2+YY] -= ty;
405             f[j_coord_offset+DIM*2+ZZ] -= tz;
406
407             /**************************
408              * CALCULATE INTERACTIONS *
409              **************************/
410
411             /* COULOMB ELECTROSTATICS */
412             velec            = qq20*rinv20;
413             felec            = velec*rinvsq20;
414
415             /* Update potential sums from outer loop */
416             velecsum        += velec;
417
418             fscal            = felec;
419
420             /* Calculate temporary vectorial force */
421             tx               = fscal*dx20;
422             ty               = fscal*dy20;
423             tz               = fscal*dz20;
424
425             /* Update vectorial force */
426             fix2            += tx;
427             fiy2            += ty;
428             fiz2            += tz;
429             f[j_coord_offset+DIM*0+XX] -= tx;
430             f[j_coord_offset+DIM*0+YY] -= ty;
431             f[j_coord_offset+DIM*0+ZZ] -= tz;
432
433             /**************************
434              * CALCULATE INTERACTIONS *
435              **************************/
436
437             /* COULOMB ELECTROSTATICS */
438             velec            = qq21*rinv21;
439             felec            = velec*rinvsq21;
440
441             /* Update potential sums from outer loop */
442             velecsum        += velec;
443
444             fscal            = felec;
445
446             /* Calculate temporary vectorial force */
447             tx               = fscal*dx21;
448             ty               = fscal*dy21;
449             tz               = fscal*dz21;
450
451             /* Update vectorial force */
452             fix2            += tx;
453             fiy2            += ty;
454             fiz2            += tz;
455             f[j_coord_offset+DIM*1+XX] -= tx;
456             f[j_coord_offset+DIM*1+YY] -= ty;
457             f[j_coord_offset+DIM*1+ZZ] -= tz;
458
459             /**************************
460              * CALCULATE INTERACTIONS *
461              **************************/
462
463             /* COULOMB ELECTROSTATICS */
464             velec            = qq22*rinv22;
465             felec            = velec*rinvsq22;
466
467             /* Update potential sums from outer loop */
468             velecsum        += velec;
469
470             fscal            = felec;
471
472             /* Calculate temporary vectorial force */
473             tx               = fscal*dx22;
474             ty               = fscal*dy22;
475             tz               = fscal*dz22;
476
477             /* Update vectorial force */
478             fix2            += tx;
479             fiy2            += ty;
480             fiz2            += tz;
481             f[j_coord_offset+DIM*2+XX] -= tx;
482             f[j_coord_offset+DIM*2+YY] -= ty;
483             f[j_coord_offset+DIM*2+ZZ] -= tz;
484
485             /* Inner loop uses 243 flops */
486         }
487         /* End of innermost loop */
488
489         tx = ty = tz = 0;
490         f[i_coord_offset+DIM*0+XX] += fix0;
491         f[i_coord_offset+DIM*0+YY] += fiy0;
492         f[i_coord_offset+DIM*0+ZZ] += fiz0;
493         tx                         += fix0;
494         ty                         += fiy0;
495         tz                         += fiz0;
496         f[i_coord_offset+DIM*1+XX] += fix1;
497         f[i_coord_offset+DIM*1+YY] += fiy1;
498         f[i_coord_offset+DIM*1+ZZ] += fiz1;
499         tx                         += fix1;
500         ty                         += fiy1;
501         tz                         += fiz1;
502         f[i_coord_offset+DIM*2+XX] += fix2;
503         f[i_coord_offset+DIM*2+YY] += fiy2;
504         f[i_coord_offset+DIM*2+ZZ] += fiz2;
505         tx                         += fix2;
506         ty                         += fiy2;
507         tz                         += fiz2;
508         fshift[i_shift_offset+XX]  += tx;
509         fshift[i_shift_offset+YY]  += ty;
510         fshift[i_shift_offset+ZZ]  += tz;
511
512         ggid                        = gid[iidx];
513         /* Update potential energies */
514         kernel_data->energygrp_elec[ggid] += velecsum;
515
516         /* Increment number of inner iterations */
517         inneriter                  += j_index_end - j_index_start;
518
519         /* Outer loop uses 31 flops */
520     }
521
522     /* Increment number of outer iterations */
523     outeriter        += nri;
524
525     /* Update outer/inner flops */
526
527     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_W3W3_VF,outeriter*31 + inneriter*243);
528 }
529 /*
530  * Gromacs nonbonded kernel:   nb_kernel_ElecCoul_VdwNone_GeomW3W3_F_c
531  * Electrostatics interaction: Coulomb
532  * VdW interaction:            None
533  * Geometry:                   Water3-Water3
534  * Calculate force/pot:        Force
535  */
536 void
537 nb_kernel_ElecCoul_VdwNone_GeomW3W3_F_c
538                     (t_nblist                    * gmx_restrict       nlist,
539                      rvec                        * gmx_restrict          xx,
540                      rvec                        * gmx_restrict          ff,
541                      t_forcerec                  * gmx_restrict          fr,
542                      t_mdatoms                   * gmx_restrict     mdatoms,
543                      nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
544                      t_nrnb                      * gmx_restrict        nrnb)
545 {
546     int              i_shift_offset,i_coord_offset,j_coord_offset;
547     int              j_index_start,j_index_end;
548     int              nri,inr,ggid,iidx,jidx,jnr,outeriter,inneriter;
549     real             shX,shY,shZ,tx,ty,tz,fscal,rcutoff,rcutoff2;
550     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
551     real             *shiftvec,*fshift,*x,*f;
552     int              vdwioffset0;
553     real             ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
554     int              vdwioffset1;
555     real             ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
556     int              vdwioffset2;
557     real             ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
558     int              vdwjidx0;
559     real             jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
560     int              vdwjidx1;
561     real             jx1,jy1,jz1,fjx1,fjy1,fjz1,jq1,isaj1;
562     int              vdwjidx2;
563     real             jx2,jy2,jz2,fjx2,fjy2,fjz2,jq2,isaj2;
564     real             dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00,cexp1_00,cexp2_00;
565     real             dx01,dy01,dz01,rsq01,rinv01,rinvsq01,r01,qq01,c6_01,c12_01,cexp1_01,cexp2_01;
566     real             dx02,dy02,dz02,rsq02,rinv02,rinvsq02,r02,qq02,c6_02,c12_02,cexp1_02,cexp2_02;
567     real             dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10,cexp1_10,cexp2_10;
568     real             dx11,dy11,dz11,rsq11,rinv11,rinvsq11,r11,qq11,c6_11,c12_11,cexp1_11,cexp2_11;
569     real             dx12,dy12,dz12,rsq12,rinv12,rinvsq12,r12,qq12,c6_12,c12_12,cexp1_12,cexp2_12;
570     real             dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20,cexp1_20,cexp2_20;
571     real             dx21,dy21,dz21,rsq21,rinv21,rinvsq21,r21,qq21,c6_21,c12_21,cexp1_21,cexp2_21;
572     real             dx22,dy22,dz22,rsq22,rinv22,rinvsq22,r22,qq22,c6_22,c12_22,cexp1_22,cexp2_22;
573     real             velec,felec,velecsum,facel,crf,krf,krf2;
574     real             *charge;
575
576     x                = xx[0];
577     f                = ff[0];
578
579     nri              = nlist->nri;
580     iinr             = nlist->iinr;
581     jindex           = nlist->jindex;
582     jjnr             = nlist->jjnr;
583     shiftidx         = nlist->shift;
584     gid              = nlist->gid;
585     shiftvec         = fr->shift_vec[0];
586     fshift           = fr->fshift[0];
587     facel            = fr->epsfac;
588     charge           = mdatoms->chargeA;
589
590     /* Setup water-specific parameters */
591     inr              = nlist->iinr[0];
592     iq0              = facel*charge[inr+0];
593     iq1              = facel*charge[inr+1];
594     iq2              = facel*charge[inr+2];
595
596     jq0              = charge[inr+0];
597     jq1              = charge[inr+1];
598     jq2              = charge[inr+2];
599     qq00             = iq0*jq0;
600     qq01             = iq0*jq1;
601     qq02             = iq0*jq2;
602     qq10             = iq1*jq0;
603     qq11             = iq1*jq1;
604     qq12             = iq1*jq2;
605     qq20             = iq2*jq0;
606     qq21             = iq2*jq1;
607     qq22             = iq2*jq2;
608
609     outeriter        = 0;
610     inneriter        = 0;
611
612     /* Start outer loop over neighborlists */
613     for(iidx=0; iidx<nri; iidx++)
614     {
615         /* Load shift vector for this list */
616         i_shift_offset   = DIM*shiftidx[iidx];
617         shX              = shiftvec[i_shift_offset+XX];
618         shY              = shiftvec[i_shift_offset+YY];
619         shZ              = shiftvec[i_shift_offset+ZZ];
620
621         /* Load limits for loop over neighbors */
622         j_index_start    = jindex[iidx];
623         j_index_end      = jindex[iidx+1];
624
625         /* Get outer coordinate index */
626         inr              = iinr[iidx];
627         i_coord_offset   = DIM*inr;
628
629         /* Load i particle coords and add shift vector */
630         ix0              = shX + x[i_coord_offset+DIM*0+XX];
631         iy0              = shY + x[i_coord_offset+DIM*0+YY];
632         iz0              = shZ + x[i_coord_offset+DIM*0+ZZ];
633         ix1              = shX + x[i_coord_offset+DIM*1+XX];
634         iy1              = shY + x[i_coord_offset+DIM*1+YY];
635         iz1              = shZ + x[i_coord_offset+DIM*1+ZZ];
636         ix2              = shX + x[i_coord_offset+DIM*2+XX];
637         iy2              = shY + x[i_coord_offset+DIM*2+YY];
638         iz2              = shZ + x[i_coord_offset+DIM*2+ZZ];
639
640         fix0             = 0.0;
641         fiy0             = 0.0;
642         fiz0             = 0.0;
643         fix1             = 0.0;
644         fiy1             = 0.0;
645         fiz1             = 0.0;
646         fix2             = 0.0;
647         fiy2             = 0.0;
648         fiz2             = 0.0;
649
650         /* Start inner kernel loop */
651         for(jidx=j_index_start; jidx<j_index_end; jidx++)
652         {
653             /* Get j neighbor index, and coordinate index */
654             jnr              = jjnr[jidx];
655             j_coord_offset   = DIM*jnr;
656
657             /* load j atom coordinates */
658             jx0              = x[j_coord_offset+DIM*0+XX];
659             jy0              = x[j_coord_offset+DIM*0+YY];
660             jz0              = x[j_coord_offset+DIM*0+ZZ];
661             jx1              = x[j_coord_offset+DIM*1+XX];
662             jy1              = x[j_coord_offset+DIM*1+YY];
663             jz1              = x[j_coord_offset+DIM*1+ZZ];
664             jx2              = x[j_coord_offset+DIM*2+XX];
665             jy2              = x[j_coord_offset+DIM*2+YY];
666             jz2              = x[j_coord_offset+DIM*2+ZZ];
667
668             /* Calculate displacement vector */
669             dx00             = ix0 - jx0;
670             dy00             = iy0 - jy0;
671             dz00             = iz0 - jz0;
672             dx01             = ix0 - jx1;
673             dy01             = iy0 - jy1;
674             dz01             = iz0 - jz1;
675             dx02             = ix0 - jx2;
676             dy02             = iy0 - jy2;
677             dz02             = iz0 - jz2;
678             dx10             = ix1 - jx0;
679             dy10             = iy1 - jy0;
680             dz10             = iz1 - jz0;
681             dx11             = ix1 - jx1;
682             dy11             = iy1 - jy1;
683             dz11             = iz1 - jz1;
684             dx12             = ix1 - jx2;
685             dy12             = iy1 - jy2;
686             dz12             = iz1 - jz2;
687             dx20             = ix2 - jx0;
688             dy20             = iy2 - jy0;
689             dz20             = iz2 - jz0;
690             dx21             = ix2 - jx1;
691             dy21             = iy2 - jy1;
692             dz21             = iz2 - jz1;
693             dx22             = ix2 - jx2;
694             dy22             = iy2 - jy2;
695             dz22             = iz2 - jz2;
696
697             /* Calculate squared distance and things based on it */
698             rsq00            = dx00*dx00+dy00*dy00+dz00*dz00;
699             rsq01            = dx01*dx01+dy01*dy01+dz01*dz01;
700             rsq02            = dx02*dx02+dy02*dy02+dz02*dz02;
701             rsq10            = dx10*dx10+dy10*dy10+dz10*dz10;
702             rsq11            = dx11*dx11+dy11*dy11+dz11*dz11;
703             rsq12            = dx12*dx12+dy12*dy12+dz12*dz12;
704             rsq20            = dx20*dx20+dy20*dy20+dz20*dz20;
705             rsq21            = dx21*dx21+dy21*dy21+dz21*dz21;
706             rsq22            = dx22*dx22+dy22*dy22+dz22*dz22;
707
708             rinv00           = gmx_invsqrt(rsq00);
709             rinv01           = gmx_invsqrt(rsq01);
710             rinv02           = gmx_invsqrt(rsq02);
711             rinv10           = gmx_invsqrt(rsq10);
712             rinv11           = gmx_invsqrt(rsq11);
713             rinv12           = gmx_invsqrt(rsq12);
714             rinv20           = gmx_invsqrt(rsq20);
715             rinv21           = gmx_invsqrt(rsq21);
716             rinv22           = gmx_invsqrt(rsq22);
717
718             rinvsq00         = rinv00*rinv00;
719             rinvsq01         = rinv01*rinv01;
720             rinvsq02         = rinv02*rinv02;
721             rinvsq10         = rinv10*rinv10;
722             rinvsq11         = rinv11*rinv11;
723             rinvsq12         = rinv12*rinv12;
724             rinvsq20         = rinv20*rinv20;
725             rinvsq21         = rinv21*rinv21;
726             rinvsq22         = rinv22*rinv22;
727
728             /**************************
729              * CALCULATE INTERACTIONS *
730              **************************/
731
732             /* COULOMB ELECTROSTATICS */
733             velec            = qq00*rinv00;
734             felec            = velec*rinvsq00;
735
736             fscal            = felec;
737
738             /* Calculate temporary vectorial force */
739             tx               = fscal*dx00;
740             ty               = fscal*dy00;
741             tz               = fscal*dz00;
742
743             /* Update vectorial force */
744             fix0            += tx;
745             fiy0            += ty;
746             fiz0            += tz;
747             f[j_coord_offset+DIM*0+XX] -= tx;
748             f[j_coord_offset+DIM*0+YY] -= ty;
749             f[j_coord_offset+DIM*0+ZZ] -= tz;
750
751             /**************************
752              * CALCULATE INTERACTIONS *
753              **************************/
754
755             /* COULOMB ELECTROSTATICS */
756             velec            = qq01*rinv01;
757             felec            = velec*rinvsq01;
758
759             fscal            = felec;
760
761             /* Calculate temporary vectorial force */
762             tx               = fscal*dx01;
763             ty               = fscal*dy01;
764             tz               = fscal*dz01;
765
766             /* Update vectorial force */
767             fix0            += tx;
768             fiy0            += ty;
769             fiz0            += tz;
770             f[j_coord_offset+DIM*1+XX] -= tx;
771             f[j_coord_offset+DIM*1+YY] -= ty;
772             f[j_coord_offset+DIM*1+ZZ] -= tz;
773
774             /**************************
775              * CALCULATE INTERACTIONS *
776              **************************/
777
778             /* COULOMB ELECTROSTATICS */
779             velec            = qq02*rinv02;
780             felec            = velec*rinvsq02;
781
782             fscal            = felec;
783
784             /* Calculate temporary vectorial force */
785             tx               = fscal*dx02;
786             ty               = fscal*dy02;
787             tz               = fscal*dz02;
788
789             /* Update vectorial force */
790             fix0            += tx;
791             fiy0            += ty;
792             fiz0            += tz;
793             f[j_coord_offset+DIM*2+XX] -= tx;
794             f[j_coord_offset+DIM*2+YY] -= ty;
795             f[j_coord_offset+DIM*2+ZZ] -= tz;
796
797             /**************************
798              * CALCULATE INTERACTIONS *
799              **************************/
800
801             /* COULOMB ELECTROSTATICS */
802             velec            = qq10*rinv10;
803             felec            = velec*rinvsq10;
804
805             fscal            = felec;
806
807             /* Calculate temporary vectorial force */
808             tx               = fscal*dx10;
809             ty               = fscal*dy10;
810             tz               = fscal*dz10;
811
812             /* Update vectorial force */
813             fix1            += tx;
814             fiy1            += ty;
815             fiz1            += tz;
816             f[j_coord_offset+DIM*0+XX] -= tx;
817             f[j_coord_offset+DIM*0+YY] -= ty;
818             f[j_coord_offset+DIM*0+ZZ] -= tz;
819
820             /**************************
821              * CALCULATE INTERACTIONS *
822              **************************/
823
824             /* COULOMB ELECTROSTATICS */
825             velec            = qq11*rinv11;
826             felec            = velec*rinvsq11;
827
828             fscal            = felec;
829
830             /* Calculate temporary vectorial force */
831             tx               = fscal*dx11;
832             ty               = fscal*dy11;
833             tz               = fscal*dz11;
834
835             /* Update vectorial force */
836             fix1            += tx;
837             fiy1            += ty;
838             fiz1            += tz;
839             f[j_coord_offset+DIM*1+XX] -= tx;
840             f[j_coord_offset+DIM*1+YY] -= ty;
841             f[j_coord_offset+DIM*1+ZZ] -= tz;
842
843             /**************************
844              * CALCULATE INTERACTIONS *
845              **************************/
846
847             /* COULOMB ELECTROSTATICS */
848             velec            = qq12*rinv12;
849             felec            = velec*rinvsq12;
850
851             fscal            = felec;
852
853             /* Calculate temporary vectorial force */
854             tx               = fscal*dx12;
855             ty               = fscal*dy12;
856             tz               = fscal*dz12;
857
858             /* Update vectorial force */
859             fix1            += tx;
860             fiy1            += ty;
861             fiz1            += tz;
862             f[j_coord_offset+DIM*2+XX] -= tx;
863             f[j_coord_offset+DIM*2+YY] -= ty;
864             f[j_coord_offset+DIM*2+ZZ] -= tz;
865
866             /**************************
867              * CALCULATE INTERACTIONS *
868              **************************/
869
870             /* COULOMB ELECTROSTATICS */
871             velec            = qq20*rinv20;
872             felec            = velec*rinvsq20;
873
874             fscal            = felec;
875
876             /* Calculate temporary vectorial force */
877             tx               = fscal*dx20;
878             ty               = fscal*dy20;
879             tz               = fscal*dz20;
880
881             /* Update vectorial force */
882             fix2            += tx;
883             fiy2            += ty;
884             fiz2            += tz;
885             f[j_coord_offset+DIM*0+XX] -= tx;
886             f[j_coord_offset+DIM*0+YY] -= ty;
887             f[j_coord_offset+DIM*0+ZZ] -= tz;
888
889             /**************************
890              * CALCULATE INTERACTIONS *
891              **************************/
892
893             /* COULOMB ELECTROSTATICS */
894             velec            = qq21*rinv21;
895             felec            = velec*rinvsq21;
896
897             fscal            = felec;
898
899             /* Calculate temporary vectorial force */
900             tx               = fscal*dx21;
901             ty               = fscal*dy21;
902             tz               = fscal*dz21;
903
904             /* Update vectorial force */
905             fix2            += tx;
906             fiy2            += ty;
907             fiz2            += tz;
908             f[j_coord_offset+DIM*1+XX] -= tx;
909             f[j_coord_offset+DIM*1+YY] -= ty;
910             f[j_coord_offset+DIM*1+ZZ] -= tz;
911
912             /**************************
913              * CALCULATE INTERACTIONS *
914              **************************/
915
916             /* COULOMB ELECTROSTATICS */
917             velec            = qq22*rinv22;
918             felec            = velec*rinvsq22;
919
920             fscal            = felec;
921
922             /* Calculate temporary vectorial force */
923             tx               = fscal*dx22;
924             ty               = fscal*dy22;
925             tz               = fscal*dz22;
926
927             /* Update vectorial force */
928             fix2            += tx;
929             fiy2            += ty;
930             fiz2            += tz;
931             f[j_coord_offset+DIM*2+XX] -= tx;
932             f[j_coord_offset+DIM*2+YY] -= ty;
933             f[j_coord_offset+DIM*2+ZZ] -= tz;
934
935             /* Inner loop uses 234 flops */
936         }
937         /* End of innermost loop */
938
939         tx = ty = tz = 0;
940         f[i_coord_offset+DIM*0+XX] += fix0;
941         f[i_coord_offset+DIM*0+YY] += fiy0;
942         f[i_coord_offset+DIM*0+ZZ] += fiz0;
943         tx                         += fix0;
944         ty                         += fiy0;
945         tz                         += fiz0;
946         f[i_coord_offset+DIM*1+XX] += fix1;
947         f[i_coord_offset+DIM*1+YY] += fiy1;
948         f[i_coord_offset+DIM*1+ZZ] += fiz1;
949         tx                         += fix1;
950         ty                         += fiy1;
951         tz                         += fiz1;
952         f[i_coord_offset+DIM*2+XX] += fix2;
953         f[i_coord_offset+DIM*2+YY] += fiy2;
954         f[i_coord_offset+DIM*2+ZZ] += fiz2;
955         tx                         += fix2;
956         ty                         += fiy2;
957         tz                         += fiz2;
958         fshift[i_shift_offset+XX]  += tx;
959         fshift[i_shift_offset+YY]  += ty;
960         fshift[i_shift_offset+ZZ]  += tz;
961
962         /* Increment number of inner iterations */
963         inneriter                  += j_index_end - j_index_start;
964
965         /* Outer loop uses 30 flops */
966     }
967
968     /* Increment number of outer iterations */
969     outeriter        += nri;
970
971     /* Update outer/inner flops */
972
973     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_W3W3_F,outeriter*30 + inneriter*234);
974 }