c7be83e493e82ba6f7ebcd73c53b3c38f75d6530
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_sse2_single / nb_kernel_ElecCoul_VdwNone_GeomW3P1_sse2_single.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 sse2_single 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 #include "gromacs/simd/math_x86_sse2_single.h"
48 #include "kernelutil_x86_sse2_single.h"
49
50 /*
51  * Gromacs nonbonded kernel:   nb_kernel_ElecCoul_VdwNone_GeomW3P1_VF_sse2_single
52  * Electrostatics interaction: Coulomb
53  * VdW interaction:            None
54  * Geometry:                   Water3-Particle
55  * Calculate force/pot:        PotentialAndForce
56  */
57 void
58 nb_kernel_ElecCoul_VdwNone_GeomW3P1_VF_sse2_single
59                     (t_nblist                    * gmx_restrict       nlist,
60                      rvec                        * gmx_restrict          xx,
61                      rvec                        * gmx_restrict          ff,
62                      t_forcerec                  * gmx_restrict          fr,
63                      t_mdatoms                   * gmx_restrict     mdatoms,
64                      nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
65                      t_nrnb                      * gmx_restrict        nrnb)
66 {
67     /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or 
68      * just 0 for non-waters.
69      * Suffixes A,B,C,D refer to j loop unrolling done with SSE, e.g. for the four different
70      * jnr indices corresponding to data put in the four positions in the SIMD register.
71      */
72     int              i_shift_offset,i_coord_offset,outeriter,inneriter;
73     int              j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
74     int              jnrA,jnrB,jnrC,jnrD;
75     int              jnrlistA,jnrlistB,jnrlistC,jnrlistD;
76     int              j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
77     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
78     real             rcutoff_scalar;
79     real             *shiftvec,*fshift,*x,*f;
80     real             *fjptrA,*fjptrB,*fjptrC,*fjptrD;
81     real             scratch[4*DIM];
82     __m128           tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
83     int              vdwioffset0;
84     __m128           ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
85     int              vdwioffset1;
86     __m128           ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
87     int              vdwioffset2;
88     __m128           ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
89     int              vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
90     __m128           jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
91     __m128           dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
92     __m128           dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
93     __m128           dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
94     __m128           velec,felec,velecsum,facel,crf,krf,krf2;
95     real             *charge;
96     __m128           dummy_mask,cutoff_mask;
97     __m128           signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
98     __m128           one     = _mm_set1_ps(1.0);
99     __m128           two     = _mm_set1_ps(2.0);
100     x                = xx[0];
101     f                = ff[0];
102
103     nri              = nlist->nri;
104     iinr             = nlist->iinr;
105     jindex           = nlist->jindex;
106     jjnr             = nlist->jjnr;
107     shiftidx         = nlist->shift;
108     gid              = nlist->gid;
109     shiftvec         = fr->shift_vec[0];
110     fshift           = fr->fshift[0];
111     facel            = _mm_set1_ps(fr->epsfac);
112     charge           = mdatoms->chargeA;
113
114     /* Setup water-specific parameters */
115     inr              = nlist->iinr[0];
116     iq0              = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+0]));
117     iq1              = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+1]));
118     iq2              = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+2]));
119
120     /* Avoid stupid compiler warnings */
121     jnrA = jnrB = jnrC = jnrD = 0;
122     j_coord_offsetA = 0;
123     j_coord_offsetB = 0;
124     j_coord_offsetC = 0;
125     j_coord_offsetD = 0;
126
127     outeriter        = 0;
128     inneriter        = 0;
129
130     for(iidx=0;iidx<4*DIM;iidx++)
131     {
132         scratch[iidx] = 0.0;
133     }  
134
135     /* Start outer loop over neighborlists */
136     for(iidx=0; iidx<nri; iidx++)
137     {
138         /* Load shift vector for this list */
139         i_shift_offset   = DIM*shiftidx[iidx];
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         gmx_mm_load_shift_and_3rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,
151                                                  &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2);
152         
153         fix0             = _mm_setzero_ps();
154         fiy0             = _mm_setzero_ps();
155         fiz0             = _mm_setzero_ps();
156         fix1             = _mm_setzero_ps();
157         fiy1             = _mm_setzero_ps();
158         fiz1             = _mm_setzero_ps();
159         fix2             = _mm_setzero_ps();
160         fiy2             = _mm_setzero_ps();
161         fiz2             = _mm_setzero_ps();
162
163         /* Reset potential sums */
164         velecsum         = _mm_setzero_ps();
165
166         /* Start inner kernel loop */
167         for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
168         {
169
170             /* Get j neighbor index, and coordinate index */
171             jnrA             = jjnr[jidx];
172             jnrB             = jjnr[jidx+1];
173             jnrC             = jjnr[jidx+2];
174             jnrD             = jjnr[jidx+3];
175             j_coord_offsetA  = DIM*jnrA;
176             j_coord_offsetB  = DIM*jnrB;
177             j_coord_offsetC  = DIM*jnrC;
178             j_coord_offsetD  = DIM*jnrD;
179
180             /* load j atom coordinates */
181             gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
182                                               x+j_coord_offsetC,x+j_coord_offsetD,
183                                               &jx0,&jy0,&jz0);
184
185             /* Calculate displacement vector */
186             dx00             = _mm_sub_ps(ix0,jx0);
187             dy00             = _mm_sub_ps(iy0,jy0);
188             dz00             = _mm_sub_ps(iz0,jz0);
189             dx10             = _mm_sub_ps(ix1,jx0);
190             dy10             = _mm_sub_ps(iy1,jy0);
191             dz10             = _mm_sub_ps(iz1,jz0);
192             dx20             = _mm_sub_ps(ix2,jx0);
193             dy20             = _mm_sub_ps(iy2,jy0);
194             dz20             = _mm_sub_ps(iz2,jz0);
195
196             /* Calculate squared distance and things based on it */
197             rsq00            = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
198             rsq10            = gmx_mm_calc_rsq_ps(dx10,dy10,dz10);
199             rsq20            = gmx_mm_calc_rsq_ps(dx20,dy20,dz20);
200
201             rinv00           = gmx_mm_invsqrt_ps(rsq00);
202             rinv10           = gmx_mm_invsqrt_ps(rsq10);
203             rinv20           = gmx_mm_invsqrt_ps(rsq20);
204
205             rinvsq00         = _mm_mul_ps(rinv00,rinv00);
206             rinvsq10         = _mm_mul_ps(rinv10,rinv10);
207             rinvsq20         = _mm_mul_ps(rinv20,rinv20);
208
209             /* Load parameters for j particles */
210             jq0              = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
211                                                               charge+jnrC+0,charge+jnrD+0);
212
213             fjx0             = _mm_setzero_ps();
214             fjy0             = _mm_setzero_ps();
215             fjz0             = _mm_setzero_ps();
216
217             /**************************
218              * CALCULATE INTERACTIONS *
219              **************************/
220
221             /* Compute parameters for interactions between i and j atoms */
222             qq00             = _mm_mul_ps(iq0,jq0);
223
224             /* COULOMB ELECTROSTATICS */
225             velec            = _mm_mul_ps(qq00,rinv00);
226             felec            = _mm_mul_ps(velec,rinvsq00);
227
228             /* Update potential sum for this i atom from the interaction with this j atom. */
229             velecsum         = _mm_add_ps(velecsum,velec);
230
231             fscal            = felec;
232
233             /* Calculate temporary vectorial force */
234             tx               = _mm_mul_ps(fscal,dx00);
235             ty               = _mm_mul_ps(fscal,dy00);
236             tz               = _mm_mul_ps(fscal,dz00);
237
238             /* Update vectorial force */
239             fix0             = _mm_add_ps(fix0,tx);
240             fiy0             = _mm_add_ps(fiy0,ty);
241             fiz0             = _mm_add_ps(fiz0,tz);
242
243             fjx0             = _mm_add_ps(fjx0,tx);
244             fjy0             = _mm_add_ps(fjy0,ty);
245             fjz0             = _mm_add_ps(fjz0,tz);
246             
247             /**************************
248              * CALCULATE INTERACTIONS *
249              **************************/
250
251             /* Compute parameters for interactions between i and j atoms */
252             qq10             = _mm_mul_ps(iq1,jq0);
253
254             /* COULOMB ELECTROSTATICS */
255             velec            = _mm_mul_ps(qq10,rinv10);
256             felec            = _mm_mul_ps(velec,rinvsq10);
257
258             /* Update potential sum for this i atom from the interaction with this j atom. */
259             velecsum         = _mm_add_ps(velecsum,velec);
260
261             fscal            = felec;
262
263             /* Calculate temporary vectorial force */
264             tx               = _mm_mul_ps(fscal,dx10);
265             ty               = _mm_mul_ps(fscal,dy10);
266             tz               = _mm_mul_ps(fscal,dz10);
267
268             /* Update vectorial force */
269             fix1             = _mm_add_ps(fix1,tx);
270             fiy1             = _mm_add_ps(fiy1,ty);
271             fiz1             = _mm_add_ps(fiz1,tz);
272
273             fjx0             = _mm_add_ps(fjx0,tx);
274             fjy0             = _mm_add_ps(fjy0,ty);
275             fjz0             = _mm_add_ps(fjz0,tz);
276             
277             /**************************
278              * CALCULATE INTERACTIONS *
279              **************************/
280
281             /* Compute parameters for interactions between i and j atoms */
282             qq20             = _mm_mul_ps(iq2,jq0);
283
284             /* COULOMB ELECTROSTATICS */
285             velec            = _mm_mul_ps(qq20,rinv20);
286             felec            = _mm_mul_ps(velec,rinvsq20);
287
288             /* Update potential sum for this i atom from the interaction with this j atom. */
289             velecsum         = _mm_add_ps(velecsum,velec);
290
291             fscal            = felec;
292
293             /* Calculate temporary vectorial force */
294             tx               = _mm_mul_ps(fscal,dx20);
295             ty               = _mm_mul_ps(fscal,dy20);
296             tz               = _mm_mul_ps(fscal,dz20);
297
298             /* Update vectorial force */
299             fix2             = _mm_add_ps(fix2,tx);
300             fiy2             = _mm_add_ps(fiy2,ty);
301             fiz2             = _mm_add_ps(fiz2,tz);
302
303             fjx0             = _mm_add_ps(fjx0,tx);
304             fjy0             = _mm_add_ps(fjy0,ty);
305             fjz0             = _mm_add_ps(fjz0,tz);
306             
307             fjptrA             = f+j_coord_offsetA;
308             fjptrB             = f+j_coord_offsetB;
309             fjptrC             = f+j_coord_offsetC;
310             fjptrD             = f+j_coord_offsetD;
311
312             gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjx0,fjy0,fjz0);
313
314             /* Inner loop uses 84 flops */
315         }
316
317         if(jidx<j_index_end)
318         {
319
320             /* Get j neighbor index, and coordinate index */
321             jnrlistA         = jjnr[jidx];
322             jnrlistB         = jjnr[jidx+1];
323             jnrlistC         = jjnr[jidx+2];
324             jnrlistD         = jjnr[jidx+3];
325             /* Sign of each element will be negative for non-real atoms.
326              * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
327              * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
328              */
329             dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
330             jnrA       = (jnrlistA>=0) ? jnrlistA : 0;
331             jnrB       = (jnrlistB>=0) ? jnrlistB : 0;
332             jnrC       = (jnrlistC>=0) ? jnrlistC : 0;
333             jnrD       = (jnrlistD>=0) ? jnrlistD : 0;
334             j_coord_offsetA  = DIM*jnrA;
335             j_coord_offsetB  = DIM*jnrB;
336             j_coord_offsetC  = DIM*jnrC;
337             j_coord_offsetD  = DIM*jnrD;
338
339             /* load j atom coordinates */
340             gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
341                                               x+j_coord_offsetC,x+j_coord_offsetD,
342                                               &jx0,&jy0,&jz0);
343
344             /* Calculate displacement vector */
345             dx00             = _mm_sub_ps(ix0,jx0);
346             dy00             = _mm_sub_ps(iy0,jy0);
347             dz00             = _mm_sub_ps(iz0,jz0);
348             dx10             = _mm_sub_ps(ix1,jx0);
349             dy10             = _mm_sub_ps(iy1,jy0);
350             dz10             = _mm_sub_ps(iz1,jz0);
351             dx20             = _mm_sub_ps(ix2,jx0);
352             dy20             = _mm_sub_ps(iy2,jy0);
353             dz20             = _mm_sub_ps(iz2,jz0);
354
355             /* Calculate squared distance and things based on it */
356             rsq00            = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
357             rsq10            = gmx_mm_calc_rsq_ps(dx10,dy10,dz10);
358             rsq20            = gmx_mm_calc_rsq_ps(dx20,dy20,dz20);
359
360             rinv00           = gmx_mm_invsqrt_ps(rsq00);
361             rinv10           = gmx_mm_invsqrt_ps(rsq10);
362             rinv20           = gmx_mm_invsqrt_ps(rsq20);
363
364             rinvsq00         = _mm_mul_ps(rinv00,rinv00);
365             rinvsq10         = _mm_mul_ps(rinv10,rinv10);
366             rinvsq20         = _mm_mul_ps(rinv20,rinv20);
367
368             /* Load parameters for j particles */
369             jq0              = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
370                                                               charge+jnrC+0,charge+jnrD+0);
371
372             fjx0             = _mm_setzero_ps();
373             fjy0             = _mm_setzero_ps();
374             fjz0             = _mm_setzero_ps();
375
376             /**************************
377              * CALCULATE INTERACTIONS *
378              **************************/
379
380             /* Compute parameters for interactions between i and j atoms */
381             qq00             = _mm_mul_ps(iq0,jq0);
382
383             /* COULOMB ELECTROSTATICS */
384             velec            = _mm_mul_ps(qq00,rinv00);
385             felec            = _mm_mul_ps(velec,rinvsq00);
386
387             /* Update potential sum for this i atom from the interaction with this j atom. */
388             velec            = _mm_andnot_ps(dummy_mask,velec);
389             velecsum         = _mm_add_ps(velecsum,velec);
390
391             fscal            = felec;
392
393             fscal            = _mm_andnot_ps(dummy_mask,fscal);
394
395             /* Calculate temporary vectorial force */
396             tx               = _mm_mul_ps(fscal,dx00);
397             ty               = _mm_mul_ps(fscal,dy00);
398             tz               = _mm_mul_ps(fscal,dz00);
399
400             /* Update vectorial force */
401             fix0             = _mm_add_ps(fix0,tx);
402             fiy0             = _mm_add_ps(fiy0,ty);
403             fiz0             = _mm_add_ps(fiz0,tz);
404
405             fjx0             = _mm_add_ps(fjx0,tx);
406             fjy0             = _mm_add_ps(fjy0,ty);
407             fjz0             = _mm_add_ps(fjz0,tz);
408             
409             /**************************
410              * CALCULATE INTERACTIONS *
411              **************************/
412
413             /* Compute parameters for interactions between i and j atoms */
414             qq10             = _mm_mul_ps(iq1,jq0);
415
416             /* COULOMB ELECTROSTATICS */
417             velec            = _mm_mul_ps(qq10,rinv10);
418             felec            = _mm_mul_ps(velec,rinvsq10);
419
420             /* Update potential sum for this i atom from the interaction with this j atom. */
421             velec            = _mm_andnot_ps(dummy_mask,velec);
422             velecsum         = _mm_add_ps(velecsum,velec);
423
424             fscal            = felec;
425
426             fscal            = _mm_andnot_ps(dummy_mask,fscal);
427
428             /* Calculate temporary vectorial force */
429             tx               = _mm_mul_ps(fscal,dx10);
430             ty               = _mm_mul_ps(fscal,dy10);
431             tz               = _mm_mul_ps(fscal,dz10);
432
433             /* Update vectorial force */
434             fix1             = _mm_add_ps(fix1,tx);
435             fiy1             = _mm_add_ps(fiy1,ty);
436             fiz1             = _mm_add_ps(fiz1,tz);
437
438             fjx0             = _mm_add_ps(fjx0,tx);
439             fjy0             = _mm_add_ps(fjy0,ty);
440             fjz0             = _mm_add_ps(fjz0,tz);
441             
442             /**************************
443              * CALCULATE INTERACTIONS *
444              **************************/
445
446             /* Compute parameters for interactions between i and j atoms */
447             qq20             = _mm_mul_ps(iq2,jq0);
448
449             /* COULOMB ELECTROSTATICS */
450             velec            = _mm_mul_ps(qq20,rinv20);
451             felec            = _mm_mul_ps(velec,rinvsq20);
452
453             /* Update potential sum for this i atom from the interaction with this j atom. */
454             velec            = _mm_andnot_ps(dummy_mask,velec);
455             velecsum         = _mm_add_ps(velecsum,velec);
456
457             fscal            = felec;
458
459             fscal            = _mm_andnot_ps(dummy_mask,fscal);
460
461             /* Calculate temporary vectorial force */
462             tx               = _mm_mul_ps(fscal,dx20);
463             ty               = _mm_mul_ps(fscal,dy20);
464             tz               = _mm_mul_ps(fscal,dz20);
465
466             /* Update vectorial force */
467             fix2             = _mm_add_ps(fix2,tx);
468             fiy2             = _mm_add_ps(fiy2,ty);
469             fiz2             = _mm_add_ps(fiz2,tz);
470
471             fjx0             = _mm_add_ps(fjx0,tx);
472             fjy0             = _mm_add_ps(fjy0,ty);
473             fjz0             = _mm_add_ps(fjz0,tz);
474             
475             fjptrA             = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
476             fjptrB             = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
477             fjptrC             = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
478             fjptrD             = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
479
480             gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjx0,fjy0,fjz0);
481
482             /* Inner loop uses 84 flops */
483         }
484
485         /* End of innermost loop */
486
487         gmx_mm_update_iforce_3atom_swizzle_ps(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,
488                                               f+i_coord_offset,fshift+i_shift_offset);
489
490         ggid                        = gid[iidx];
491         /* Update potential energies */
492         gmx_mm_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
493
494         /* Increment number of inner iterations */
495         inneriter                  += j_index_end - j_index_start;
496
497         /* Outer loop uses 19 flops */
498     }
499
500     /* Increment number of outer iterations */
501     outeriter        += nri;
502
503     /* Update outer/inner flops */
504
505     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_W3_VF,outeriter*19 + inneriter*84);
506 }
507 /*
508  * Gromacs nonbonded kernel:   nb_kernel_ElecCoul_VdwNone_GeomW3P1_F_sse2_single
509  * Electrostatics interaction: Coulomb
510  * VdW interaction:            None
511  * Geometry:                   Water3-Particle
512  * Calculate force/pot:        Force
513  */
514 void
515 nb_kernel_ElecCoul_VdwNone_GeomW3P1_F_sse2_single
516                     (t_nblist                    * gmx_restrict       nlist,
517                      rvec                        * gmx_restrict          xx,
518                      rvec                        * gmx_restrict          ff,
519                      t_forcerec                  * gmx_restrict          fr,
520                      t_mdatoms                   * gmx_restrict     mdatoms,
521                      nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
522                      t_nrnb                      * gmx_restrict        nrnb)
523 {
524     /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or 
525      * just 0 for non-waters.
526      * Suffixes A,B,C,D refer to j loop unrolling done with SSE, e.g. for the four different
527      * jnr indices corresponding to data put in the four positions in the SIMD register.
528      */
529     int              i_shift_offset,i_coord_offset,outeriter,inneriter;
530     int              j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
531     int              jnrA,jnrB,jnrC,jnrD;
532     int              jnrlistA,jnrlistB,jnrlistC,jnrlistD;
533     int              j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
534     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
535     real             rcutoff_scalar;
536     real             *shiftvec,*fshift,*x,*f;
537     real             *fjptrA,*fjptrB,*fjptrC,*fjptrD;
538     real             scratch[4*DIM];
539     __m128           tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
540     int              vdwioffset0;
541     __m128           ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
542     int              vdwioffset1;
543     __m128           ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
544     int              vdwioffset2;
545     __m128           ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
546     int              vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
547     __m128           jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
548     __m128           dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
549     __m128           dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
550     __m128           dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
551     __m128           velec,felec,velecsum,facel,crf,krf,krf2;
552     real             *charge;
553     __m128           dummy_mask,cutoff_mask;
554     __m128           signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
555     __m128           one     = _mm_set1_ps(1.0);
556     __m128           two     = _mm_set1_ps(2.0);
557     x                = xx[0];
558     f                = ff[0];
559
560     nri              = nlist->nri;
561     iinr             = nlist->iinr;
562     jindex           = nlist->jindex;
563     jjnr             = nlist->jjnr;
564     shiftidx         = nlist->shift;
565     gid              = nlist->gid;
566     shiftvec         = fr->shift_vec[0];
567     fshift           = fr->fshift[0];
568     facel            = _mm_set1_ps(fr->epsfac);
569     charge           = mdatoms->chargeA;
570
571     /* Setup water-specific parameters */
572     inr              = nlist->iinr[0];
573     iq0              = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+0]));
574     iq1              = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+1]));
575     iq2              = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+2]));
576
577     /* Avoid stupid compiler warnings */
578     jnrA = jnrB = jnrC = jnrD = 0;
579     j_coord_offsetA = 0;
580     j_coord_offsetB = 0;
581     j_coord_offsetC = 0;
582     j_coord_offsetD = 0;
583
584     outeriter        = 0;
585     inneriter        = 0;
586
587     for(iidx=0;iidx<4*DIM;iidx++)
588     {
589         scratch[iidx] = 0.0;
590     }  
591
592     /* Start outer loop over neighborlists */
593     for(iidx=0; iidx<nri; iidx++)
594     {
595         /* Load shift vector for this list */
596         i_shift_offset   = DIM*shiftidx[iidx];
597
598         /* Load limits for loop over neighbors */
599         j_index_start    = jindex[iidx];
600         j_index_end      = jindex[iidx+1];
601
602         /* Get outer coordinate index */
603         inr              = iinr[iidx];
604         i_coord_offset   = DIM*inr;
605
606         /* Load i particle coords and add shift vector */
607         gmx_mm_load_shift_and_3rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,
608                                                  &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2);
609         
610         fix0             = _mm_setzero_ps();
611         fiy0             = _mm_setzero_ps();
612         fiz0             = _mm_setzero_ps();
613         fix1             = _mm_setzero_ps();
614         fiy1             = _mm_setzero_ps();
615         fiz1             = _mm_setzero_ps();
616         fix2             = _mm_setzero_ps();
617         fiy2             = _mm_setzero_ps();
618         fiz2             = _mm_setzero_ps();
619
620         /* Start inner kernel loop */
621         for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
622         {
623
624             /* Get j neighbor index, and coordinate index */
625             jnrA             = jjnr[jidx];
626             jnrB             = jjnr[jidx+1];
627             jnrC             = jjnr[jidx+2];
628             jnrD             = jjnr[jidx+3];
629             j_coord_offsetA  = DIM*jnrA;
630             j_coord_offsetB  = DIM*jnrB;
631             j_coord_offsetC  = DIM*jnrC;
632             j_coord_offsetD  = DIM*jnrD;
633
634             /* load j atom coordinates */
635             gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
636                                               x+j_coord_offsetC,x+j_coord_offsetD,
637                                               &jx0,&jy0,&jz0);
638
639             /* Calculate displacement vector */
640             dx00             = _mm_sub_ps(ix0,jx0);
641             dy00             = _mm_sub_ps(iy0,jy0);
642             dz00             = _mm_sub_ps(iz0,jz0);
643             dx10             = _mm_sub_ps(ix1,jx0);
644             dy10             = _mm_sub_ps(iy1,jy0);
645             dz10             = _mm_sub_ps(iz1,jz0);
646             dx20             = _mm_sub_ps(ix2,jx0);
647             dy20             = _mm_sub_ps(iy2,jy0);
648             dz20             = _mm_sub_ps(iz2,jz0);
649
650             /* Calculate squared distance and things based on it */
651             rsq00            = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
652             rsq10            = gmx_mm_calc_rsq_ps(dx10,dy10,dz10);
653             rsq20            = gmx_mm_calc_rsq_ps(dx20,dy20,dz20);
654
655             rinv00           = gmx_mm_invsqrt_ps(rsq00);
656             rinv10           = gmx_mm_invsqrt_ps(rsq10);
657             rinv20           = gmx_mm_invsqrt_ps(rsq20);
658
659             rinvsq00         = _mm_mul_ps(rinv00,rinv00);
660             rinvsq10         = _mm_mul_ps(rinv10,rinv10);
661             rinvsq20         = _mm_mul_ps(rinv20,rinv20);
662
663             /* Load parameters for j particles */
664             jq0              = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
665                                                               charge+jnrC+0,charge+jnrD+0);
666
667             fjx0             = _mm_setzero_ps();
668             fjy0             = _mm_setzero_ps();
669             fjz0             = _mm_setzero_ps();
670
671             /**************************
672              * CALCULATE INTERACTIONS *
673              **************************/
674
675             /* Compute parameters for interactions between i and j atoms */
676             qq00             = _mm_mul_ps(iq0,jq0);
677
678             /* COULOMB ELECTROSTATICS */
679             velec            = _mm_mul_ps(qq00,rinv00);
680             felec            = _mm_mul_ps(velec,rinvsq00);
681
682             fscal            = felec;
683
684             /* Calculate temporary vectorial force */
685             tx               = _mm_mul_ps(fscal,dx00);
686             ty               = _mm_mul_ps(fscal,dy00);
687             tz               = _mm_mul_ps(fscal,dz00);
688
689             /* Update vectorial force */
690             fix0             = _mm_add_ps(fix0,tx);
691             fiy0             = _mm_add_ps(fiy0,ty);
692             fiz0             = _mm_add_ps(fiz0,tz);
693
694             fjx0             = _mm_add_ps(fjx0,tx);
695             fjy0             = _mm_add_ps(fjy0,ty);
696             fjz0             = _mm_add_ps(fjz0,tz);
697             
698             /**************************
699              * CALCULATE INTERACTIONS *
700              **************************/
701
702             /* Compute parameters for interactions between i and j atoms */
703             qq10             = _mm_mul_ps(iq1,jq0);
704
705             /* COULOMB ELECTROSTATICS */
706             velec            = _mm_mul_ps(qq10,rinv10);
707             felec            = _mm_mul_ps(velec,rinvsq10);
708
709             fscal            = felec;
710
711             /* Calculate temporary vectorial force */
712             tx               = _mm_mul_ps(fscal,dx10);
713             ty               = _mm_mul_ps(fscal,dy10);
714             tz               = _mm_mul_ps(fscal,dz10);
715
716             /* Update vectorial force */
717             fix1             = _mm_add_ps(fix1,tx);
718             fiy1             = _mm_add_ps(fiy1,ty);
719             fiz1             = _mm_add_ps(fiz1,tz);
720
721             fjx0             = _mm_add_ps(fjx0,tx);
722             fjy0             = _mm_add_ps(fjy0,ty);
723             fjz0             = _mm_add_ps(fjz0,tz);
724             
725             /**************************
726              * CALCULATE INTERACTIONS *
727              **************************/
728
729             /* Compute parameters for interactions between i and j atoms */
730             qq20             = _mm_mul_ps(iq2,jq0);
731
732             /* COULOMB ELECTROSTATICS */
733             velec            = _mm_mul_ps(qq20,rinv20);
734             felec            = _mm_mul_ps(velec,rinvsq20);
735
736             fscal            = felec;
737
738             /* Calculate temporary vectorial force */
739             tx               = _mm_mul_ps(fscal,dx20);
740             ty               = _mm_mul_ps(fscal,dy20);
741             tz               = _mm_mul_ps(fscal,dz20);
742
743             /* Update vectorial force */
744             fix2             = _mm_add_ps(fix2,tx);
745             fiy2             = _mm_add_ps(fiy2,ty);
746             fiz2             = _mm_add_ps(fiz2,tz);
747
748             fjx0             = _mm_add_ps(fjx0,tx);
749             fjy0             = _mm_add_ps(fjy0,ty);
750             fjz0             = _mm_add_ps(fjz0,tz);
751             
752             fjptrA             = f+j_coord_offsetA;
753             fjptrB             = f+j_coord_offsetB;
754             fjptrC             = f+j_coord_offsetC;
755             fjptrD             = f+j_coord_offsetD;
756
757             gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjx0,fjy0,fjz0);
758
759             /* Inner loop uses 81 flops */
760         }
761
762         if(jidx<j_index_end)
763         {
764
765             /* Get j neighbor index, and coordinate index */
766             jnrlistA         = jjnr[jidx];
767             jnrlistB         = jjnr[jidx+1];
768             jnrlistC         = jjnr[jidx+2];
769             jnrlistD         = jjnr[jidx+3];
770             /* Sign of each element will be negative for non-real atoms.
771              * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
772              * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
773              */
774             dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
775             jnrA       = (jnrlistA>=0) ? jnrlistA : 0;
776             jnrB       = (jnrlistB>=0) ? jnrlistB : 0;
777             jnrC       = (jnrlistC>=0) ? jnrlistC : 0;
778             jnrD       = (jnrlistD>=0) ? jnrlistD : 0;
779             j_coord_offsetA  = DIM*jnrA;
780             j_coord_offsetB  = DIM*jnrB;
781             j_coord_offsetC  = DIM*jnrC;
782             j_coord_offsetD  = DIM*jnrD;
783
784             /* load j atom coordinates */
785             gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
786                                               x+j_coord_offsetC,x+j_coord_offsetD,
787                                               &jx0,&jy0,&jz0);
788
789             /* Calculate displacement vector */
790             dx00             = _mm_sub_ps(ix0,jx0);
791             dy00             = _mm_sub_ps(iy0,jy0);
792             dz00             = _mm_sub_ps(iz0,jz0);
793             dx10             = _mm_sub_ps(ix1,jx0);
794             dy10             = _mm_sub_ps(iy1,jy0);
795             dz10             = _mm_sub_ps(iz1,jz0);
796             dx20             = _mm_sub_ps(ix2,jx0);
797             dy20             = _mm_sub_ps(iy2,jy0);
798             dz20             = _mm_sub_ps(iz2,jz0);
799
800             /* Calculate squared distance and things based on it */
801             rsq00            = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
802             rsq10            = gmx_mm_calc_rsq_ps(dx10,dy10,dz10);
803             rsq20            = gmx_mm_calc_rsq_ps(dx20,dy20,dz20);
804
805             rinv00           = gmx_mm_invsqrt_ps(rsq00);
806             rinv10           = gmx_mm_invsqrt_ps(rsq10);
807             rinv20           = gmx_mm_invsqrt_ps(rsq20);
808
809             rinvsq00         = _mm_mul_ps(rinv00,rinv00);
810             rinvsq10         = _mm_mul_ps(rinv10,rinv10);
811             rinvsq20         = _mm_mul_ps(rinv20,rinv20);
812
813             /* Load parameters for j particles */
814             jq0              = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
815                                                               charge+jnrC+0,charge+jnrD+0);
816
817             fjx0             = _mm_setzero_ps();
818             fjy0             = _mm_setzero_ps();
819             fjz0             = _mm_setzero_ps();
820
821             /**************************
822              * CALCULATE INTERACTIONS *
823              **************************/
824
825             /* Compute parameters for interactions between i and j atoms */
826             qq00             = _mm_mul_ps(iq0,jq0);
827
828             /* COULOMB ELECTROSTATICS */
829             velec            = _mm_mul_ps(qq00,rinv00);
830             felec            = _mm_mul_ps(velec,rinvsq00);
831
832             fscal            = felec;
833
834             fscal            = _mm_andnot_ps(dummy_mask,fscal);
835
836             /* Calculate temporary vectorial force */
837             tx               = _mm_mul_ps(fscal,dx00);
838             ty               = _mm_mul_ps(fscal,dy00);
839             tz               = _mm_mul_ps(fscal,dz00);
840
841             /* Update vectorial force */
842             fix0             = _mm_add_ps(fix0,tx);
843             fiy0             = _mm_add_ps(fiy0,ty);
844             fiz0             = _mm_add_ps(fiz0,tz);
845
846             fjx0             = _mm_add_ps(fjx0,tx);
847             fjy0             = _mm_add_ps(fjy0,ty);
848             fjz0             = _mm_add_ps(fjz0,tz);
849             
850             /**************************
851              * CALCULATE INTERACTIONS *
852              **************************/
853
854             /* Compute parameters for interactions between i and j atoms */
855             qq10             = _mm_mul_ps(iq1,jq0);
856
857             /* COULOMB ELECTROSTATICS */
858             velec            = _mm_mul_ps(qq10,rinv10);
859             felec            = _mm_mul_ps(velec,rinvsq10);
860
861             fscal            = felec;
862
863             fscal            = _mm_andnot_ps(dummy_mask,fscal);
864
865             /* Calculate temporary vectorial force */
866             tx               = _mm_mul_ps(fscal,dx10);
867             ty               = _mm_mul_ps(fscal,dy10);
868             tz               = _mm_mul_ps(fscal,dz10);
869
870             /* Update vectorial force */
871             fix1             = _mm_add_ps(fix1,tx);
872             fiy1             = _mm_add_ps(fiy1,ty);
873             fiz1             = _mm_add_ps(fiz1,tz);
874
875             fjx0             = _mm_add_ps(fjx0,tx);
876             fjy0             = _mm_add_ps(fjy0,ty);
877             fjz0             = _mm_add_ps(fjz0,tz);
878             
879             /**************************
880              * CALCULATE INTERACTIONS *
881              **************************/
882
883             /* Compute parameters for interactions between i and j atoms */
884             qq20             = _mm_mul_ps(iq2,jq0);
885
886             /* COULOMB ELECTROSTATICS */
887             velec            = _mm_mul_ps(qq20,rinv20);
888             felec            = _mm_mul_ps(velec,rinvsq20);
889
890             fscal            = felec;
891
892             fscal            = _mm_andnot_ps(dummy_mask,fscal);
893
894             /* Calculate temporary vectorial force */
895             tx               = _mm_mul_ps(fscal,dx20);
896             ty               = _mm_mul_ps(fscal,dy20);
897             tz               = _mm_mul_ps(fscal,dz20);
898
899             /* Update vectorial force */
900             fix2             = _mm_add_ps(fix2,tx);
901             fiy2             = _mm_add_ps(fiy2,ty);
902             fiz2             = _mm_add_ps(fiz2,tz);
903
904             fjx0             = _mm_add_ps(fjx0,tx);
905             fjy0             = _mm_add_ps(fjy0,ty);
906             fjz0             = _mm_add_ps(fjz0,tz);
907             
908             fjptrA             = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
909             fjptrB             = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
910             fjptrC             = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
911             fjptrD             = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
912
913             gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjx0,fjy0,fjz0);
914
915             /* Inner loop uses 81 flops */
916         }
917
918         /* End of innermost loop */
919
920         gmx_mm_update_iforce_3atom_swizzle_ps(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,
921                                               f+i_coord_offset,fshift+i_shift_offset);
922
923         /* Increment number of inner iterations */
924         inneriter                  += j_index_end - j_index_start;
925
926         /* Outer loop uses 18 flops */
927     }
928
929     /* Increment number of outer iterations */
930     outeriter        += nri;
931
932     /* Update outer/inner flops */
933
934     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_W3_F,outeriter*18 + inneriter*81);
935 }