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