ffce96fc8695635dbf2c90208e72ae167b1b9165
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_avx_128_fma_single / nb_kernel_ElecCSTab_VdwNone_GeomP1P1_avx_128_fma_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 avx_128_fma_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_avx_128_fma_single.h"
48 #include "kernelutil_x86_avx_128_fma_single.h"
49
50 /*
51  * Gromacs nonbonded kernel:   nb_kernel_ElecCSTab_VdwNone_GeomP1P1_VF_avx_128_fma_single
52  * Electrostatics interaction: CubicSplineTable
53  * VdW interaction:            None
54  * Geometry:                   Particle-Particle
55  * Calculate force/pot:        PotentialAndForce
56  */
57 void
58 nb_kernel_ElecCSTab_VdwNone_GeomP1P1_VF_avx_128_fma_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 AVX_128, 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           fscal,rcutoff,rcutoff2,jidxall;
83     int              vdwioffset0;
84     __m128           ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
85     int              vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
86     __m128           jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
87     __m128           dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
88     __m128           velec,felec,velecsum,facel,crf,krf,krf2;
89     real             *charge;
90     __m128i          vfitab;
91     __m128i          ifour       = _mm_set1_epi32(4);
92     __m128           rt,vfeps,twovfeps,vftabscale,Y,F,G,H,Fp,VV,FF;
93     real             *vftab;
94     __m128           dummy_mask,cutoff_mask;
95     __m128           signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
96     __m128           one     = _mm_set1_ps(1.0);
97     __m128           two     = _mm_set1_ps(2.0);
98     x                = xx[0];
99     f                = ff[0];
100
101     nri              = nlist->nri;
102     iinr             = nlist->iinr;
103     jindex           = nlist->jindex;
104     jjnr             = nlist->jjnr;
105     shiftidx         = nlist->shift;
106     gid              = nlist->gid;
107     shiftvec         = fr->shift_vec[0];
108     fshift           = fr->fshift[0];
109     facel            = _mm_set1_ps(fr->epsfac);
110     charge           = mdatoms->chargeA;
111
112     vftab            = kernel_data->table_elec->data;
113     vftabscale       = _mm_set1_ps(kernel_data->table_elec->scale);
114
115     /* Avoid stupid compiler warnings */
116     jnrA = jnrB = jnrC = jnrD = 0;
117     j_coord_offsetA = 0;
118     j_coord_offsetB = 0;
119     j_coord_offsetC = 0;
120     j_coord_offsetD = 0;
121
122     outeriter        = 0;
123     inneriter        = 0;
124
125     for(iidx=0;iidx<4*DIM;iidx++)
126     {
127         scratch[iidx] = 0.0;
128     }
129
130     /* Start outer loop over neighborlists */
131     for(iidx=0; iidx<nri; iidx++)
132     {
133         /* Load shift vector for this list */
134         i_shift_offset   = DIM*shiftidx[iidx];
135
136         /* Load limits for loop over neighbors */
137         j_index_start    = jindex[iidx];
138         j_index_end      = jindex[iidx+1];
139
140         /* Get outer coordinate index */
141         inr              = iinr[iidx];
142         i_coord_offset   = DIM*inr;
143
144         /* Load i particle coords and add shift vector */
145         gmx_mm_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
146
147         fix0             = _mm_setzero_ps();
148         fiy0             = _mm_setzero_ps();
149         fiz0             = _mm_setzero_ps();
150
151         /* Load parameters for i particles */
152         iq0              = _mm_mul_ps(facel,_mm_load1_ps(charge+inr+0));
153
154         /* Reset potential sums */
155         velecsum         = _mm_setzero_ps();
156
157         /* Start inner kernel loop */
158         for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
159         {
160
161             /* Get j neighbor index, and coordinate index */
162             jnrA             = jjnr[jidx];
163             jnrB             = jjnr[jidx+1];
164             jnrC             = jjnr[jidx+2];
165             jnrD             = jjnr[jidx+3];
166             j_coord_offsetA  = DIM*jnrA;
167             j_coord_offsetB  = DIM*jnrB;
168             j_coord_offsetC  = DIM*jnrC;
169             j_coord_offsetD  = DIM*jnrD;
170
171             /* load j atom coordinates */
172             gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
173                                               x+j_coord_offsetC,x+j_coord_offsetD,
174                                               &jx0,&jy0,&jz0);
175
176             /* Calculate displacement vector */
177             dx00             = _mm_sub_ps(ix0,jx0);
178             dy00             = _mm_sub_ps(iy0,jy0);
179             dz00             = _mm_sub_ps(iz0,jz0);
180
181             /* Calculate squared distance and things based on it */
182             rsq00            = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
183
184             rinv00           = gmx_mm_invsqrt_ps(rsq00);
185
186             /* Load parameters for j particles */
187             jq0              = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
188                                                               charge+jnrC+0,charge+jnrD+0);
189
190             /**************************
191              * CALCULATE INTERACTIONS *
192              **************************/
193
194             r00              = _mm_mul_ps(rsq00,rinv00);
195
196             /* Compute parameters for interactions between i and j atoms */
197             qq00             = _mm_mul_ps(iq0,jq0);
198
199             /* Calculate table index by multiplying r with table scale and truncate to integer */
200             rt               = _mm_mul_ps(r00,vftabscale);
201             vfitab           = _mm_cvttps_epi32(rt);
202 #ifdef __XOP__
203             vfeps            = _mm_frcz_ps(rt);
204 #else
205             vfeps            = _mm_sub_ps(rt,_mm_round_ps(rt, _MM_FROUND_FLOOR));
206 #endif
207             twovfeps         = _mm_add_ps(vfeps,vfeps);
208             vfitab           = _mm_slli_epi32(vfitab,2);
209
210             /* CUBIC SPLINE TABLE ELECTROSTATICS */
211             Y                = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,0) );
212             F                = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,1) );
213             G                = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,2) );
214             H                = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,3) );
215             _MM_TRANSPOSE4_PS(Y,F,G,H);
216             Fp               = _mm_macc_ps(vfeps,_mm_macc_ps(H,vfeps,G),F);
217             VV               = _mm_macc_ps(vfeps,Fp,Y);
218             velec            = _mm_mul_ps(qq00,VV);
219             FF               = _mm_macc_ps(vfeps,_mm_macc_ps(twovfeps,H,G),Fp);
220             felec            = _mm_xor_ps(signbit,_mm_mul_ps(_mm_mul_ps(qq00,FF),_mm_mul_ps(vftabscale,rinv00)));
221
222             /* Update potential sum for this i atom from the interaction with this j atom. */
223             velecsum         = _mm_add_ps(velecsum,velec);
224
225             fscal            = felec;
226
227              /* Update vectorial force */
228             fix0             = _mm_macc_ps(dx00,fscal,fix0);
229             fiy0             = _mm_macc_ps(dy00,fscal,fiy0);
230             fiz0             = _mm_macc_ps(dz00,fscal,fiz0);
231
232             fjptrA             = f+j_coord_offsetA;
233             fjptrB             = f+j_coord_offsetB;
234             fjptrC             = f+j_coord_offsetC;
235             fjptrD             = f+j_coord_offsetD;
236             gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,
237                                                    _mm_mul_ps(dx00,fscal),
238                                                    _mm_mul_ps(dy00,fscal),
239                                                    _mm_mul_ps(dz00,fscal));
240
241             /* Inner loop uses 46 flops */
242         }
243
244         if(jidx<j_index_end)
245         {
246
247             /* Get j neighbor index, and coordinate index */
248             jnrlistA         = jjnr[jidx];
249             jnrlistB         = jjnr[jidx+1];
250             jnrlistC         = jjnr[jidx+2];
251             jnrlistD         = jjnr[jidx+3];
252             /* Sign of each element will be negative for non-real atoms.
253              * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
254              * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
255              */
256             dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
257             jnrA       = (jnrlistA>=0) ? jnrlistA : 0;
258             jnrB       = (jnrlistB>=0) ? jnrlistB : 0;
259             jnrC       = (jnrlistC>=0) ? jnrlistC : 0;
260             jnrD       = (jnrlistD>=0) ? jnrlistD : 0;
261             j_coord_offsetA  = DIM*jnrA;
262             j_coord_offsetB  = DIM*jnrB;
263             j_coord_offsetC  = DIM*jnrC;
264             j_coord_offsetD  = DIM*jnrD;
265
266             /* load j atom coordinates */
267             gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
268                                               x+j_coord_offsetC,x+j_coord_offsetD,
269                                               &jx0,&jy0,&jz0);
270
271             /* Calculate displacement vector */
272             dx00             = _mm_sub_ps(ix0,jx0);
273             dy00             = _mm_sub_ps(iy0,jy0);
274             dz00             = _mm_sub_ps(iz0,jz0);
275
276             /* Calculate squared distance and things based on it */
277             rsq00            = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
278
279             rinv00           = gmx_mm_invsqrt_ps(rsq00);
280
281             /* Load parameters for j particles */
282             jq0              = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
283                                                               charge+jnrC+0,charge+jnrD+0);
284
285             /**************************
286              * CALCULATE INTERACTIONS *
287              **************************/
288
289             r00              = _mm_mul_ps(rsq00,rinv00);
290             r00              = _mm_andnot_ps(dummy_mask,r00);
291
292             /* Compute parameters for interactions between i and j atoms */
293             qq00             = _mm_mul_ps(iq0,jq0);
294
295             /* Calculate table index by multiplying r with table scale and truncate to integer */
296             rt               = _mm_mul_ps(r00,vftabscale);
297             vfitab           = _mm_cvttps_epi32(rt);
298 #ifdef __XOP__
299             vfeps            = _mm_frcz_ps(rt);
300 #else
301             vfeps            = _mm_sub_ps(rt,_mm_round_ps(rt, _MM_FROUND_FLOOR));
302 #endif
303             twovfeps         = _mm_add_ps(vfeps,vfeps);
304             vfitab           = _mm_slli_epi32(vfitab,2);
305
306             /* CUBIC SPLINE TABLE ELECTROSTATICS */
307             Y                = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,0) );
308             F                = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,1) );
309             G                = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,2) );
310             H                = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,3) );
311             _MM_TRANSPOSE4_PS(Y,F,G,H);
312             Fp               = _mm_macc_ps(vfeps,_mm_macc_ps(H,vfeps,G),F);
313             VV               = _mm_macc_ps(vfeps,Fp,Y);
314             velec            = _mm_mul_ps(qq00,VV);
315             FF               = _mm_macc_ps(vfeps,_mm_macc_ps(twovfeps,H,G),Fp);
316             felec            = _mm_xor_ps(signbit,_mm_mul_ps(_mm_mul_ps(qq00,FF),_mm_mul_ps(vftabscale,rinv00)));
317
318             /* Update potential sum for this i atom from the interaction with this j atom. */
319             velec            = _mm_andnot_ps(dummy_mask,velec);
320             velecsum         = _mm_add_ps(velecsum,velec);
321
322             fscal            = felec;
323
324             fscal            = _mm_andnot_ps(dummy_mask,fscal);
325
326              /* Update vectorial force */
327             fix0             = _mm_macc_ps(dx00,fscal,fix0);
328             fiy0             = _mm_macc_ps(dy00,fscal,fiy0);
329             fiz0             = _mm_macc_ps(dz00,fscal,fiz0);
330
331             fjptrA             = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
332             fjptrB             = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
333             fjptrC             = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
334             fjptrD             = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
335             gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,
336                                                    _mm_mul_ps(dx00,fscal),
337                                                    _mm_mul_ps(dy00,fscal),
338                                                    _mm_mul_ps(dz00,fscal));
339
340             /* Inner loop uses 47 flops */
341         }
342
343         /* End of innermost loop */
344
345         gmx_mm_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
346                                               f+i_coord_offset,fshift+i_shift_offset);
347
348         ggid                        = gid[iidx];
349         /* Update potential energies */
350         gmx_mm_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
351
352         /* Increment number of inner iterations */
353         inneriter                  += j_index_end - j_index_start;
354
355         /* Outer loop uses 8 flops */
356     }
357
358     /* Increment number of outer iterations */
359     outeriter        += nri;
360
361     /* Update outer/inner flops */
362
363     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VF,outeriter*8 + inneriter*47);
364 }
365 /*
366  * Gromacs nonbonded kernel:   nb_kernel_ElecCSTab_VdwNone_GeomP1P1_F_avx_128_fma_single
367  * Electrostatics interaction: CubicSplineTable
368  * VdW interaction:            None
369  * Geometry:                   Particle-Particle
370  * Calculate force/pot:        Force
371  */
372 void
373 nb_kernel_ElecCSTab_VdwNone_GeomP1P1_F_avx_128_fma_single
374                     (t_nblist                    * gmx_restrict       nlist,
375                      rvec                        * gmx_restrict          xx,
376                      rvec                        * gmx_restrict          ff,
377                      t_forcerec                  * gmx_restrict          fr,
378                      t_mdatoms                   * gmx_restrict     mdatoms,
379                      nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
380                      t_nrnb                      * gmx_restrict        nrnb)
381 {
382     /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
383      * just 0 for non-waters.
384      * Suffixes A,B,C,D refer to j loop unrolling done with AVX_128, e.g. for the four different
385      * jnr indices corresponding to data put in the four positions in the SIMD register.
386      */
387     int              i_shift_offset,i_coord_offset,outeriter,inneriter;
388     int              j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
389     int              jnrA,jnrB,jnrC,jnrD;
390     int              jnrlistA,jnrlistB,jnrlistC,jnrlistD;
391     int              j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
392     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
393     real             rcutoff_scalar;
394     real             *shiftvec,*fshift,*x,*f;
395     real             *fjptrA,*fjptrB,*fjptrC,*fjptrD;
396     real             scratch[4*DIM];
397     __m128           fscal,rcutoff,rcutoff2,jidxall;
398     int              vdwioffset0;
399     __m128           ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
400     int              vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
401     __m128           jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
402     __m128           dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
403     __m128           velec,felec,velecsum,facel,crf,krf,krf2;
404     real             *charge;
405     __m128i          vfitab;
406     __m128i          ifour       = _mm_set1_epi32(4);
407     __m128           rt,vfeps,twovfeps,vftabscale,Y,F,G,H,Fp,VV,FF;
408     real             *vftab;
409     __m128           dummy_mask,cutoff_mask;
410     __m128           signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
411     __m128           one     = _mm_set1_ps(1.0);
412     __m128           two     = _mm_set1_ps(2.0);
413     x                = xx[0];
414     f                = ff[0];
415
416     nri              = nlist->nri;
417     iinr             = nlist->iinr;
418     jindex           = nlist->jindex;
419     jjnr             = nlist->jjnr;
420     shiftidx         = nlist->shift;
421     gid              = nlist->gid;
422     shiftvec         = fr->shift_vec[0];
423     fshift           = fr->fshift[0];
424     facel            = _mm_set1_ps(fr->epsfac);
425     charge           = mdatoms->chargeA;
426
427     vftab            = kernel_data->table_elec->data;
428     vftabscale       = _mm_set1_ps(kernel_data->table_elec->scale);
429
430     /* Avoid stupid compiler warnings */
431     jnrA = jnrB = jnrC = jnrD = 0;
432     j_coord_offsetA = 0;
433     j_coord_offsetB = 0;
434     j_coord_offsetC = 0;
435     j_coord_offsetD = 0;
436
437     outeriter        = 0;
438     inneriter        = 0;
439
440     for(iidx=0;iidx<4*DIM;iidx++)
441     {
442         scratch[iidx] = 0.0;
443     }
444
445     /* Start outer loop over neighborlists */
446     for(iidx=0; iidx<nri; iidx++)
447     {
448         /* Load shift vector for this list */
449         i_shift_offset   = DIM*shiftidx[iidx];
450
451         /* Load limits for loop over neighbors */
452         j_index_start    = jindex[iidx];
453         j_index_end      = jindex[iidx+1];
454
455         /* Get outer coordinate index */
456         inr              = iinr[iidx];
457         i_coord_offset   = DIM*inr;
458
459         /* Load i particle coords and add shift vector */
460         gmx_mm_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
461
462         fix0             = _mm_setzero_ps();
463         fiy0             = _mm_setzero_ps();
464         fiz0             = _mm_setzero_ps();
465
466         /* Load parameters for i particles */
467         iq0              = _mm_mul_ps(facel,_mm_load1_ps(charge+inr+0));
468
469         /* Start inner kernel loop */
470         for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
471         {
472
473             /* Get j neighbor index, and coordinate index */
474             jnrA             = jjnr[jidx];
475             jnrB             = jjnr[jidx+1];
476             jnrC             = jjnr[jidx+2];
477             jnrD             = jjnr[jidx+3];
478             j_coord_offsetA  = DIM*jnrA;
479             j_coord_offsetB  = DIM*jnrB;
480             j_coord_offsetC  = DIM*jnrC;
481             j_coord_offsetD  = DIM*jnrD;
482
483             /* load j atom coordinates */
484             gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
485                                               x+j_coord_offsetC,x+j_coord_offsetD,
486                                               &jx0,&jy0,&jz0);
487
488             /* Calculate displacement vector */
489             dx00             = _mm_sub_ps(ix0,jx0);
490             dy00             = _mm_sub_ps(iy0,jy0);
491             dz00             = _mm_sub_ps(iz0,jz0);
492
493             /* Calculate squared distance and things based on it */
494             rsq00            = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
495
496             rinv00           = gmx_mm_invsqrt_ps(rsq00);
497
498             /* Load parameters for j particles */
499             jq0              = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
500                                                               charge+jnrC+0,charge+jnrD+0);
501
502             /**************************
503              * CALCULATE INTERACTIONS *
504              **************************/
505
506             r00              = _mm_mul_ps(rsq00,rinv00);
507
508             /* Compute parameters for interactions between i and j atoms */
509             qq00             = _mm_mul_ps(iq0,jq0);
510
511             /* Calculate table index by multiplying r with table scale and truncate to integer */
512             rt               = _mm_mul_ps(r00,vftabscale);
513             vfitab           = _mm_cvttps_epi32(rt);
514 #ifdef __XOP__
515             vfeps            = _mm_frcz_ps(rt);
516 #else
517             vfeps            = _mm_sub_ps(rt,_mm_round_ps(rt, _MM_FROUND_FLOOR));
518 #endif
519             twovfeps         = _mm_add_ps(vfeps,vfeps);
520             vfitab           = _mm_slli_epi32(vfitab,2);
521
522             /* CUBIC SPLINE TABLE ELECTROSTATICS */
523             Y                = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,0) );
524             F                = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,1) );
525             G                = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,2) );
526             H                = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,3) );
527             _MM_TRANSPOSE4_PS(Y,F,G,H);
528             Fp               = _mm_macc_ps(vfeps,_mm_macc_ps(H,vfeps,G),F);
529             FF               = _mm_macc_ps(vfeps,_mm_macc_ps(twovfeps,H,G),Fp);
530             felec            = _mm_xor_ps(signbit,_mm_mul_ps(_mm_mul_ps(qq00,FF),_mm_mul_ps(vftabscale,rinv00)));
531
532             fscal            = felec;
533
534              /* Update vectorial force */
535             fix0             = _mm_macc_ps(dx00,fscal,fix0);
536             fiy0             = _mm_macc_ps(dy00,fscal,fiy0);
537             fiz0             = _mm_macc_ps(dz00,fscal,fiz0);
538
539             fjptrA             = f+j_coord_offsetA;
540             fjptrB             = f+j_coord_offsetB;
541             fjptrC             = f+j_coord_offsetC;
542             fjptrD             = f+j_coord_offsetD;
543             gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,
544                                                    _mm_mul_ps(dx00,fscal),
545                                                    _mm_mul_ps(dy00,fscal),
546                                                    _mm_mul_ps(dz00,fscal));
547
548             /* Inner loop uses 42 flops */
549         }
550
551         if(jidx<j_index_end)
552         {
553
554             /* Get j neighbor index, and coordinate index */
555             jnrlistA         = jjnr[jidx];
556             jnrlistB         = jjnr[jidx+1];
557             jnrlistC         = jjnr[jidx+2];
558             jnrlistD         = jjnr[jidx+3];
559             /* Sign of each element will be negative for non-real atoms.
560              * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
561              * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
562              */
563             dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
564             jnrA       = (jnrlistA>=0) ? jnrlistA : 0;
565             jnrB       = (jnrlistB>=0) ? jnrlistB : 0;
566             jnrC       = (jnrlistC>=0) ? jnrlistC : 0;
567             jnrD       = (jnrlistD>=0) ? jnrlistD : 0;
568             j_coord_offsetA  = DIM*jnrA;
569             j_coord_offsetB  = DIM*jnrB;
570             j_coord_offsetC  = DIM*jnrC;
571             j_coord_offsetD  = DIM*jnrD;
572
573             /* load j atom coordinates */
574             gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
575                                               x+j_coord_offsetC,x+j_coord_offsetD,
576                                               &jx0,&jy0,&jz0);
577
578             /* Calculate displacement vector */
579             dx00             = _mm_sub_ps(ix0,jx0);
580             dy00             = _mm_sub_ps(iy0,jy0);
581             dz00             = _mm_sub_ps(iz0,jz0);
582
583             /* Calculate squared distance and things based on it */
584             rsq00            = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
585
586             rinv00           = gmx_mm_invsqrt_ps(rsq00);
587
588             /* Load parameters for j particles */
589             jq0              = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
590                                                               charge+jnrC+0,charge+jnrD+0);
591
592             /**************************
593              * CALCULATE INTERACTIONS *
594              **************************/
595
596             r00              = _mm_mul_ps(rsq00,rinv00);
597             r00              = _mm_andnot_ps(dummy_mask,r00);
598
599             /* Compute parameters for interactions between i and j atoms */
600             qq00             = _mm_mul_ps(iq0,jq0);
601
602             /* Calculate table index by multiplying r with table scale and truncate to integer */
603             rt               = _mm_mul_ps(r00,vftabscale);
604             vfitab           = _mm_cvttps_epi32(rt);
605 #ifdef __XOP__
606             vfeps            = _mm_frcz_ps(rt);
607 #else
608             vfeps            = _mm_sub_ps(rt,_mm_round_ps(rt, _MM_FROUND_FLOOR));
609 #endif
610             twovfeps         = _mm_add_ps(vfeps,vfeps);
611             vfitab           = _mm_slli_epi32(vfitab,2);
612
613             /* CUBIC SPLINE TABLE ELECTROSTATICS */
614             Y                = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,0) );
615             F                = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,1) );
616             G                = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,2) );
617             H                = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,3) );
618             _MM_TRANSPOSE4_PS(Y,F,G,H);
619             Fp               = _mm_macc_ps(vfeps,_mm_macc_ps(H,vfeps,G),F);
620             FF               = _mm_macc_ps(vfeps,_mm_macc_ps(twovfeps,H,G),Fp);
621             felec            = _mm_xor_ps(signbit,_mm_mul_ps(_mm_mul_ps(qq00,FF),_mm_mul_ps(vftabscale,rinv00)));
622
623             fscal            = felec;
624
625             fscal            = _mm_andnot_ps(dummy_mask,fscal);
626
627              /* Update vectorial force */
628             fix0             = _mm_macc_ps(dx00,fscal,fix0);
629             fiy0             = _mm_macc_ps(dy00,fscal,fiy0);
630             fiz0             = _mm_macc_ps(dz00,fscal,fiz0);
631
632             fjptrA             = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
633             fjptrB             = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
634             fjptrC             = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
635             fjptrD             = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
636             gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,
637                                                    _mm_mul_ps(dx00,fscal),
638                                                    _mm_mul_ps(dy00,fscal),
639                                                    _mm_mul_ps(dz00,fscal));
640
641             /* Inner loop uses 43 flops */
642         }
643
644         /* End of innermost loop */
645
646         gmx_mm_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
647                                               f+i_coord_offset,fshift+i_shift_offset);
648
649         /* Increment number of inner iterations */
650         inneriter                  += j_index_end - j_index_start;
651
652         /* Outer loop uses 7 flops */
653     }
654
655     /* Increment number of outer iterations */
656     outeriter        += nri;
657
658     /* Update outer/inner flops */
659
660     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_F,outeriter*7 + inneriter*43);
661 }