made errors during GPU detection non-fatal
[alexxy/gromacs.git] / src / gmxlib / nonbonded / nb_kernel_avx_128_fma_single / nb_kernel_ElecEwSh_VdwLJSh_GeomP1P1_avx_128_fma_single.c
1 /*
2  * Note: this file was generated by the Gromacs avx_128_fma_single kernel generator.
3  *
4  *                This source code is part of
5  *
6  *                 G   R   O   M   A   C   S
7  *
8  * Copyright (c) 2001-2012, The GROMACS Development Team
9  *
10  * Gromacs is a library for molecular simulation and trajectory analysis,
11  * written by Erik Lindahl, David van der Spoel, Berk Hess, and others - for
12  * a full list of developers and information, check out http://www.gromacs.org
13  *
14  * This program is free software; you can redistribute it and/or modify it under
15  * the terms of the GNU Lesser General Public License as published by the Free
16  * Software Foundation; either version 2 of the License, or (at your option) any
17  * later version.
18  *
19  * To help fund GROMACS development, we humbly ask that you cite
20  * the papers people have written on it - you can find them on the website.
21  */
22 #ifdef HAVE_CONFIG_H
23 #include <config.h>
24 #endif
25
26 #include <math.h>
27
28 #include "../nb_kernel.h"
29 #include "types/simple.h"
30 #include "vec.h"
31 #include "nrnb.h"
32
33 #include "gmx_math_x86_avx_128_fma_single.h"
34 #include "kernelutil_x86_avx_128_fma_single.h"
35
36 /*
37  * Gromacs nonbonded kernel:   nb_kernel_ElecEwSh_VdwLJSh_GeomP1P1_VF_avx_128_fma_single
38  * Electrostatics interaction: Ewald
39  * VdW interaction:            LennardJones
40  * Geometry:                   Particle-Particle
41  * Calculate force/pot:        PotentialAndForce
42  */
43 void
44 nb_kernel_ElecEwSh_VdwLJSh_GeomP1P1_VF_avx_128_fma_single
45                     (t_nblist * gmx_restrict                nlist,
46                      rvec * gmx_restrict                    xx,
47                      rvec * gmx_restrict                    ff,
48                      t_forcerec * gmx_restrict              fr,
49                      t_mdatoms * gmx_restrict               mdatoms,
50                      nb_kernel_data_t * gmx_restrict        kernel_data,
51                      t_nrnb * gmx_restrict                  nrnb)
52 {
53     /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
54      * just 0 for non-waters.
55      * Suffixes A,B,C,D refer to j loop unrolling done with AVX_128, e.g. for the four different
56      * jnr indices corresponding to data put in the four positions in the SIMD register.
57      */
58     int              i_shift_offset,i_coord_offset,outeriter,inneriter;
59     int              j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
60     int              jnrA,jnrB,jnrC,jnrD;
61     int              jnrlistA,jnrlistB,jnrlistC,jnrlistD;
62     int              j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
63     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
64     real             rcutoff_scalar;
65     real             *shiftvec,*fshift,*x,*f;
66     real             *fjptrA,*fjptrB,*fjptrC,*fjptrD;
67     real             scratch[4*DIM];
68     __m128           fscal,rcutoff,rcutoff2,jidxall;
69     int              vdwioffset0;
70     __m128           ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
71     int              vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
72     __m128           jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
73     __m128           dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
74     __m128           velec,felec,velecsum,facel,crf,krf,krf2;
75     real             *charge;
76     int              nvdwtype;
77     __m128           rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
78     int              *vdwtype;
79     real             *vdwparam;
80     __m128           one_sixth   = _mm_set1_ps(1.0/6.0);
81     __m128           one_twelfth = _mm_set1_ps(1.0/12.0);
82     __m128i          ewitab;
83     __m128           ewtabscale,eweps,twoeweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
84     __m128           beta,beta2,beta3,zeta2,pmecorrF,pmecorrV,rinv3;
85     real             *ewtab;
86     __m128           dummy_mask,cutoff_mask;
87     __m128           signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
88     __m128           one     = _mm_set1_ps(1.0);
89     __m128           two     = _mm_set1_ps(2.0);
90     x                = xx[0];
91     f                = ff[0];
92
93     nri              = nlist->nri;
94     iinr             = nlist->iinr;
95     jindex           = nlist->jindex;
96     jjnr             = nlist->jjnr;
97     shiftidx         = nlist->shift;
98     gid              = nlist->gid;
99     shiftvec         = fr->shift_vec[0];
100     fshift           = fr->fshift[0];
101     facel            = _mm_set1_ps(fr->epsfac);
102     charge           = mdatoms->chargeA;
103     nvdwtype         = fr->ntype;
104     vdwparam         = fr->nbfp;
105     vdwtype          = mdatoms->typeA;
106
107     sh_ewald         = _mm_set1_ps(fr->ic->sh_ewald);
108     beta             = _mm_set1_ps(fr->ic->ewaldcoeff);
109     beta2            = _mm_mul_ps(beta,beta);
110     beta3            = _mm_mul_ps(beta,beta2);
111     ewtab            = fr->ic->tabq_coul_FDV0;
112     ewtabscale       = _mm_set1_ps(fr->ic->tabq_scale);
113     ewtabhalfspace   = _mm_set1_ps(0.5/fr->ic->tabq_scale);
114
115     /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
116     rcutoff_scalar   = fr->rcoulomb;
117     rcutoff          = _mm_set1_ps(rcutoff_scalar);
118     rcutoff2         = _mm_mul_ps(rcutoff,rcutoff);
119
120     sh_vdw_invrcut6  = _mm_set1_ps(fr->ic->sh_invrc6);
121     rvdw             = _mm_set1_ps(fr->rvdw);
122
123     /* Avoid stupid compiler warnings */
124     jnrA = jnrB = jnrC = jnrD = 0;
125     j_coord_offsetA = 0;
126     j_coord_offsetB = 0;
127     j_coord_offsetC = 0;
128     j_coord_offsetD = 0;
129
130     outeriter        = 0;
131     inneriter        = 0;
132
133     for(iidx=0;iidx<4*DIM;iidx++)
134     {
135         scratch[iidx] = 0.0;
136     }
137
138     /* Start outer loop over neighborlists */
139     for(iidx=0; iidx<nri; iidx++)
140     {
141         /* Load shift vector for this list */
142         i_shift_offset   = DIM*shiftidx[iidx];
143
144         /* Load limits for loop over neighbors */
145         j_index_start    = jindex[iidx];
146         j_index_end      = jindex[iidx+1];
147
148         /* Get outer coordinate index */
149         inr              = iinr[iidx];
150         i_coord_offset   = DIM*inr;
151
152         /* Load i particle coords and add shift vector */
153         gmx_mm_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
154
155         fix0             = _mm_setzero_ps();
156         fiy0             = _mm_setzero_ps();
157         fiz0             = _mm_setzero_ps();
158
159         /* Load parameters for i particles */
160         iq0              = _mm_mul_ps(facel,_mm_load1_ps(charge+inr+0));
161         vdwioffset0      = 2*nvdwtype*vdwtype[inr+0];
162
163         /* Reset potential sums */
164         velecsum         = _mm_setzero_ps();
165         vvdwsum          = _mm_setzero_ps();
166
167         /* Start inner kernel loop */
168         for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
169         {
170
171             /* Get j neighbor index, and coordinate index */
172             jnrA             = jjnr[jidx];
173             jnrB             = jjnr[jidx+1];
174             jnrC             = jjnr[jidx+2];
175             jnrD             = jjnr[jidx+3];
176             j_coord_offsetA  = DIM*jnrA;
177             j_coord_offsetB  = DIM*jnrB;
178             j_coord_offsetC  = DIM*jnrC;
179             j_coord_offsetD  = DIM*jnrD;
180
181             /* load j atom coordinates */
182             gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
183                                               x+j_coord_offsetC,x+j_coord_offsetD,
184                                               &jx0,&jy0,&jz0);
185
186             /* Calculate displacement vector */
187             dx00             = _mm_sub_ps(ix0,jx0);
188             dy00             = _mm_sub_ps(iy0,jy0);
189             dz00             = _mm_sub_ps(iz0,jz0);
190
191             /* Calculate squared distance and things based on it */
192             rsq00            = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
193
194             rinv00           = gmx_mm_invsqrt_ps(rsq00);
195
196             rinvsq00         = _mm_mul_ps(rinv00,rinv00);
197
198             /* Load parameters for j particles */
199             jq0              = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
200                                                               charge+jnrC+0,charge+jnrD+0);
201             vdwjidx0A        = 2*vdwtype[jnrA+0];
202             vdwjidx0B        = 2*vdwtype[jnrB+0];
203             vdwjidx0C        = 2*vdwtype[jnrC+0];
204             vdwjidx0D        = 2*vdwtype[jnrD+0];
205
206             /**************************
207              * CALCULATE INTERACTIONS *
208              **************************/
209
210             if (gmx_mm_any_lt(rsq00,rcutoff2))
211             {
212
213             r00              = _mm_mul_ps(rsq00,rinv00);
214
215             /* Compute parameters for interactions between i and j atoms */
216             qq00             = _mm_mul_ps(iq0,jq0);
217             gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
218                                          vdwparam+vdwioffset0+vdwjidx0B,
219                                          vdwparam+vdwioffset0+vdwjidx0C,
220                                          vdwparam+vdwioffset0+vdwjidx0D,
221                                          &c6_00,&c12_00);
222
223             /* EWALD ELECTROSTATICS */
224
225             /* Analytical PME correction */
226             zeta2            = _mm_mul_ps(beta2,rsq00);
227             rinv3            = _mm_mul_ps(rinvsq00,rinv00);
228             pmecorrF         = gmx_mm_pmecorrF_ps(zeta2);
229             felec            = _mm_macc_ps(pmecorrF,beta3,rinv3);
230             felec            = _mm_mul_ps(qq00,felec);
231             pmecorrV         = gmx_mm_pmecorrV_ps(zeta2);
232             velec            = _mm_nmacc_ps(pmecorrV,beta,_mm_sub_ps(rinv00,sh_ewald));
233             velec            = _mm_mul_ps(qq00,velec);
234
235             /* LENNARD-JONES DISPERSION/REPULSION */
236
237             rinvsix          = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
238             vvdw6            = _mm_mul_ps(c6_00,rinvsix);
239             vvdw12           = _mm_mul_ps(c12_00,_mm_mul_ps(rinvsix,rinvsix));
240             vvdw             = _mm_msub_ps(_mm_nmacc_ps(c12_00,_mm_mul_ps(sh_vdw_invrcut6,sh_vdw_invrcut6),vvdw12),one_twelfth,
241                                           _mm_mul_ps( _mm_nmacc_ps(c6_00,sh_vdw_invrcut6,vvdw6),one_sixth));
242             fvdw             = _mm_mul_ps(_mm_sub_ps(vvdw12,vvdw6),rinvsq00);
243
244             cutoff_mask      = _mm_cmplt_ps(rsq00,rcutoff2);
245
246             /* Update potential sum for this i atom from the interaction with this j atom. */
247             velec            = _mm_and_ps(velec,cutoff_mask);
248             velecsum         = _mm_add_ps(velecsum,velec);
249             vvdw             = _mm_and_ps(vvdw,cutoff_mask);
250             vvdwsum          = _mm_add_ps(vvdwsum,vvdw);
251
252             fscal            = _mm_add_ps(felec,fvdw);
253
254             fscal            = _mm_and_ps(fscal,cutoff_mask);
255
256              /* Update vectorial force */
257             fix0             = _mm_macc_ps(dx00,fscal,fix0);
258             fiy0             = _mm_macc_ps(dy00,fscal,fiy0);
259             fiz0             = _mm_macc_ps(dz00,fscal,fiz0);
260
261             fjptrA             = f+j_coord_offsetA;
262             fjptrB             = f+j_coord_offsetB;
263             fjptrC             = f+j_coord_offsetC;
264             fjptrD             = f+j_coord_offsetD;
265             gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,
266                                                    _mm_mul_ps(dx00,fscal),
267                                                    _mm_mul_ps(dy00,fscal),
268                                                    _mm_mul_ps(dz00,fscal));
269
270             }
271
272             /* Inner loop uses 51 flops */
273         }
274
275         if(jidx<j_index_end)
276         {
277
278             /* Get j neighbor index, and coordinate index */
279             jnrlistA         = jjnr[jidx];
280             jnrlistB         = jjnr[jidx+1];
281             jnrlistC         = jjnr[jidx+2];
282             jnrlistD         = jjnr[jidx+3];
283             /* Sign of each element will be negative for non-real atoms.
284              * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
285              * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
286              */
287             dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
288             jnrA       = (jnrlistA>=0) ? jnrlistA : 0;
289             jnrB       = (jnrlistB>=0) ? jnrlistB : 0;
290             jnrC       = (jnrlistC>=0) ? jnrlistC : 0;
291             jnrD       = (jnrlistD>=0) ? jnrlistD : 0;
292             j_coord_offsetA  = DIM*jnrA;
293             j_coord_offsetB  = DIM*jnrB;
294             j_coord_offsetC  = DIM*jnrC;
295             j_coord_offsetD  = DIM*jnrD;
296
297             /* load j atom coordinates */
298             gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
299                                               x+j_coord_offsetC,x+j_coord_offsetD,
300                                               &jx0,&jy0,&jz0);
301
302             /* Calculate displacement vector */
303             dx00             = _mm_sub_ps(ix0,jx0);
304             dy00             = _mm_sub_ps(iy0,jy0);
305             dz00             = _mm_sub_ps(iz0,jz0);
306
307             /* Calculate squared distance and things based on it */
308             rsq00            = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
309
310             rinv00           = gmx_mm_invsqrt_ps(rsq00);
311
312             rinvsq00         = _mm_mul_ps(rinv00,rinv00);
313
314             /* Load parameters for j particles */
315             jq0              = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
316                                                               charge+jnrC+0,charge+jnrD+0);
317             vdwjidx0A        = 2*vdwtype[jnrA+0];
318             vdwjidx0B        = 2*vdwtype[jnrB+0];
319             vdwjidx0C        = 2*vdwtype[jnrC+0];
320             vdwjidx0D        = 2*vdwtype[jnrD+0];
321
322             /**************************
323              * CALCULATE INTERACTIONS *
324              **************************/
325
326             if (gmx_mm_any_lt(rsq00,rcutoff2))
327             {
328
329             r00              = _mm_mul_ps(rsq00,rinv00);
330             r00              = _mm_andnot_ps(dummy_mask,r00);
331
332             /* Compute parameters for interactions between i and j atoms */
333             qq00             = _mm_mul_ps(iq0,jq0);
334             gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
335                                          vdwparam+vdwioffset0+vdwjidx0B,
336                                          vdwparam+vdwioffset0+vdwjidx0C,
337                                          vdwparam+vdwioffset0+vdwjidx0D,
338                                          &c6_00,&c12_00);
339
340             /* EWALD ELECTROSTATICS */
341
342             /* Analytical PME correction */
343             zeta2            = _mm_mul_ps(beta2,rsq00);
344             rinv3            = _mm_mul_ps(rinvsq00,rinv00);
345             pmecorrF         = gmx_mm_pmecorrF_ps(zeta2);
346             felec            = _mm_macc_ps(pmecorrF,beta3,rinv3);
347             felec            = _mm_mul_ps(qq00,felec);
348             pmecorrV         = gmx_mm_pmecorrV_ps(zeta2);
349             velec            = _mm_nmacc_ps(pmecorrV,beta,_mm_sub_ps(rinv00,sh_ewald));
350             velec            = _mm_mul_ps(qq00,velec);
351
352             /* LENNARD-JONES DISPERSION/REPULSION */
353
354             rinvsix          = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
355             vvdw6            = _mm_mul_ps(c6_00,rinvsix);
356             vvdw12           = _mm_mul_ps(c12_00,_mm_mul_ps(rinvsix,rinvsix));
357             vvdw             = _mm_msub_ps(_mm_nmacc_ps(c12_00,_mm_mul_ps(sh_vdw_invrcut6,sh_vdw_invrcut6),vvdw12),one_twelfth,
358                                           _mm_mul_ps( _mm_nmacc_ps(c6_00,sh_vdw_invrcut6,vvdw6),one_sixth));
359             fvdw             = _mm_mul_ps(_mm_sub_ps(vvdw12,vvdw6),rinvsq00);
360
361             cutoff_mask      = _mm_cmplt_ps(rsq00,rcutoff2);
362
363             /* Update potential sum for this i atom from the interaction with this j atom. */
364             velec            = _mm_and_ps(velec,cutoff_mask);
365             velec            = _mm_andnot_ps(dummy_mask,velec);
366             velecsum         = _mm_add_ps(velecsum,velec);
367             vvdw             = _mm_and_ps(vvdw,cutoff_mask);
368             vvdw             = _mm_andnot_ps(dummy_mask,vvdw);
369             vvdwsum          = _mm_add_ps(vvdwsum,vvdw);
370
371             fscal            = _mm_add_ps(felec,fvdw);
372
373             fscal            = _mm_and_ps(fscal,cutoff_mask);
374
375             fscal            = _mm_andnot_ps(dummy_mask,fscal);
376
377              /* Update vectorial force */
378             fix0             = _mm_macc_ps(dx00,fscal,fix0);
379             fiy0             = _mm_macc_ps(dy00,fscal,fiy0);
380             fiz0             = _mm_macc_ps(dz00,fscal,fiz0);
381
382             fjptrA             = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
383             fjptrB             = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
384             fjptrC             = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
385             fjptrD             = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
386             gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,
387                                                    _mm_mul_ps(dx00,fscal),
388                                                    _mm_mul_ps(dy00,fscal),
389                                                    _mm_mul_ps(dz00,fscal));
390
391             }
392
393             /* Inner loop uses 52 flops */
394         }
395
396         /* End of innermost loop */
397
398         gmx_mm_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
399                                               f+i_coord_offset,fshift+i_shift_offset);
400
401         ggid                        = gid[iidx];
402         /* Update potential energies */
403         gmx_mm_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
404         gmx_mm_update_1pot_ps(vvdwsum,kernel_data->energygrp_vdw+ggid);
405
406         /* Increment number of inner iterations */
407         inneriter                  += j_index_end - j_index_start;
408
409         /* Outer loop uses 9 flops */
410     }
411
412     /* Increment number of outer iterations */
413     outeriter        += nri;
414
415     /* Update outer/inner flops */
416
417     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_VF,outeriter*9 + inneriter*52);
418 }
419 /*
420  * Gromacs nonbonded kernel:   nb_kernel_ElecEwSh_VdwLJSh_GeomP1P1_F_avx_128_fma_single
421  * Electrostatics interaction: Ewald
422  * VdW interaction:            LennardJones
423  * Geometry:                   Particle-Particle
424  * Calculate force/pot:        Force
425  */
426 void
427 nb_kernel_ElecEwSh_VdwLJSh_GeomP1P1_F_avx_128_fma_single
428                     (t_nblist * gmx_restrict                nlist,
429                      rvec * gmx_restrict                    xx,
430                      rvec * gmx_restrict                    ff,
431                      t_forcerec * gmx_restrict              fr,
432                      t_mdatoms * gmx_restrict               mdatoms,
433                      nb_kernel_data_t * gmx_restrict        kernel_data,
434                      t_nrnb * gmx_restrict                  nrnb)
435 {
436     /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
437      * just 0 for non-waters.
438      * Suffixes A,B,C,D refer to j loop unrolling done with AVX_128, e.g. for the four different
439      * jnr indices corresponding to data put in the four positions in the SIMD register.
440      */
441     int              i_shift_offset,i_coord_offset,outeriter,inneriter;
442     int              j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
443     int              jnrA,jnrB,jnrC,jnrD;
444     int              jnrlistA,jnrlistB,jnrlistC,jnrlistD;
445     int              j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
446     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
447     real             rcutoff_scalar;
448     real             *shiftvec,*fshift,*x,*f;
449     real             *fjptrA,*fjptrB,*fjptrC,*fjptrD;
450     real             scratch[4*DIM];
451     __m128           fscal,rcutoff,rcutoff2,jidxall;
452     int              vdwioffset0;
453     __m128           ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
454     int              vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
455     __m128           jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
456     __m128           dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
457     __m128           velec,felec,velecsum,facel,crf,krf,krf2;
458     real             *charge;
459     int              nvdwtype;
460     __m128           rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
461     int              *vdwtype;
462     real             *vdwparam;
463     __m128           one_sixth   = _mm_set1_ps(1.0/6.0);
464     __m128           one_twelfth = _mm_set1_ps(1.0/12.0);
465     __m128i          ewitab;
466     __m128           ewtabscale,eweps,twoeweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
467     __m128           beta,beta2,beta3,zeta2,pmecorrF,pmecorrV,rinv3;
468     real             *ewtab;
469     __m128           dummy_mask,cutoff_mask;
470     __m128           signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
471     __m128           one     = _mm_set1_ps(1.0);
472     __m128           two     = _mm_set1_ps(2.0);
473     x                = xx[0];
474     f                = ff[0];
475
476     nri              = nlist->nri;
477     iinr             = nlist->iinr;
478     jindex           = nlist->jindex;
479     jjnr             = nlist->jjnr;
480     shiftidx         = nlist->shift;
481     gid              = nlist->gid;
482     shiftvec         = fr->shift_vec[0];
483     fshift           = fr->fshift[0];
484     facel            = _mm_set1_ps(fr->epsfac);
485     charge           = mdatoms->chargeA;
486     nvdwtype         = fr->ntype;
487     vdwparam         = fr->nbfp;
488     vdwtype          = mdatoms->typeA;
489
490     sh_ewald         = _mm_set1_ps(fr->ic->sh_ewald);
491     beta             = _mm_set1_ps(fr->ic->ewaldcoeff);
492     beta2            = _mm_mul_ps(beta,beta);
493     beta3            = _mm_mul_ps(beta,beta2);
494     ewtab            = fr->ic->tabq_coul_F;
495     ewtabscale       = _mm_set1_ps(fr->ic->tabq_scale);
496     ewtabhalfspace   = _mm_set1_ps(0.5/fr->ic->tabq_scale);
497
498     /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
499     rcutoff_scalar   = fr->rcoulomb;
500     rcutoff          = _mm_set1_ps(rcutoff_scalar);
501     rcutoff2         = _mm_mul_ps(rcutoff,rcutoff);
502
503     sh_vdw_invrcut6  = _mm_set1_ps(fr->ic->sh_invrc6);
504     rvdw             = _mm_set1_ps(fr->rvdw);
505
506     /* Avoid stupid compiler warnings */
507     jnrA = jnrB = jnrC = jnrD = 0;
508     j_coord_offsetA = 0;
509     j_coord_offsetB = 0;
510     j_coord_offsetC = 0;
511     j_coord_offsetD = 0;
512
513     outeriter        = 0;
514     inneriter        = 0;
515
516     for(iidx=0;iidx<4*DIM;iidx++)
517     {
518         scratch[iidx] = 0.0;
519     }
520
521     /* Start outer loop over neighborlists */
522     for(iidx=0; iidx<nri; iidx++)
523     {
524         /* Load shift vector for this list */
525         i_shift_offset   = DIM*shiftidx[iidx];
526
527         /* Load limits for loop over neighbors */
528         j_index_start    = jindex[iidx];
529         j_index_end      = jindex[iidx+1];
530
531         /* Get outer coordinate index */
532         inr              = iinr[iidx];
533         i_coord_offset   = DIM*inr;
534
535         /* Load i particle coords and add shift vector */
536         gmx_mm_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
537
538         fix0             = _mm_setzero_ps();
539         fiy0             = _mm_setzero_ps();
540         fiz0             = _mm_setzero_ps();
541
542         /* Load parameters for i particles */
543         iq0              = _mm_mul_ps(facel,_mm_load1_ps(charge+inr+0));
544         vdwioffset0      = 2*nvdwtype*vdwtype[inr+0];
545
546         /* Start inner kernel loop */
547         for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
548         {
549
550             /* Get j neighbor index, and coordinate index */
551             jnrA             = jjnr[jidx];
552             jnrB             = jjnr[jidx+1];
553             jnrC             = jjnr[jidx+2];
554             jnrD             = jjnr[jidx+3];
555             j_coord_offsetA  = DIM*jnrA;
556             j_coord_offsetB  = DIM*jnrB;
557             j_coord_offsetC  = DIM*jnrC;
558             j_coord_offsetD  = DIM*jnrD;
559
560             /* load j atom coordinates */
561             gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
562                                               x+j_coord_offsetC,x+j_coord_offsetD,
563                                               &jx0,&jy0,&jz0);
564
565             /* Calculate displacement vector */
566             dx00             = _mm_sub_ps(ix0,jx0);
567             dy00             = _mm_sub_ps(iy0,jy0);
568             dz00             = _mm_sub_ps(iz0,jz0);
569
570             /* Calculate squared distance and things based on it */
571             rsq00            = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
572
573             rinv00           = gmx_mm_invsqrt_ps(rsq00);
574
575             rinvsq00         = _mm_mul_ps(rinv00,rinv00);
576
577             /* Load parameters for j particles */
578             jq0              = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
579                                                               charge+jnrC+0,charge+jnrD+0);
580             vdwjidx0A        = 2*vdwtype[jnrA+0];
581             vdwjidx0B        = 2*vdwtype[jnrB+0];
582             vdwjidx0C        = 2*vdwtype[jnrC+0];
583             vdwjidx0D        = 2*vdwtype[jnrD+0];
584
585             /**************************
586              * CALCULATE INTERACTIONS *
587              **************************/
588
589             if (gmx_mm_any_lt(rsq00,rcutoff2))
590             {
591
592             r00              = _mm_mul_ps(rsq00,rinv00);
593
594             /* Compute parameters for interactions between i and j atoms */
595             qq00             = _mm_mul_ps(iq0,jq0);
596             gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
597                                          vdwparam+vdwioffset0+vdwjidx0B,
598                                          vdwparam+vdwioffset0+vdwjidx0C,
599                                          vdwparam+vdwioffset0+vdwjidx0D,
600                                          &c6_00,&c12_00);
601
602             /* EWALD ELECTROSTATICS */
603
604             /* Analytical PME correction */
605             zeta2            = _mm_mul_ps(beta2,rsq00);
606             rinv3            = _mm_mul_ps(rinvsq00,rinv00);
607             pmecorrF         = gmx_mm_pmecorrF_ps(zeta2);
608             felec            = _mm_macc_ps(pmecorrF,beta3,rinv3);
609             felec            = _mm_mul_ps(qq00,felec);
610
611             /* LENNARD-JONES DISPERSION/REPULSION */
612
613             rinvsix          = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
614             fvdw             = _mm_mul_ps(_mm_msub_ps(c12_00,rinvsix,c6_00),_mm_mul_ps(rinvsix,rinvsq00));
615
616             cutoff_mask      = _mm_cmplt_ps(rsq00,rcutoff2);
617
618             fscal            = _mm_add_ps(felec,fvdw);
619
620             fscal            = _mm_and_ps(fscal,cutoff_mask);
621
622              /* Update vectorial force */
623             fix0             = _mm_macc_ps(dx00,fscal,fix0);
624             fiy0             = _mm_macc_ps(dy00,fscal,fiy0);
625             fiz0             = _mm_macc_ps(dz00,fscal,fiz0);
626
627             fjptrA             = f+j_coord_offsetA;
628             fjptrB             = f+j_coord_offsetB;
629             fjptrC             = f+j_coord_offsetC;
630             fjptrD             = f+j_coord_offsetD;
631             gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,
632                                                    _mm_mul_ps(dx00,fscal),
633                                                    _mm_mul_ps(dy00,fscal),
634                                                    _mm_mul_ps(dz00,fscal));
635
636             }
637
638             /* Inner loop uses 38 flops */
639         }
640
641         if(jidx<j_index_end)
642         {
643
644             /* Get j neighbor index, and coordinate index */
645             jnrlistA         = jjnr[jidx];
646             jnrlistB         = jjnr[jidx+1];
647             jnrlistC         = jjnr[jidx+2];
648             jnrlistD         = jjnr[jidx+3];
649             /* Sign of each element will be negative for non-real atoms.
650              * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
651              * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
652              */
653             dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
654             jnrA       = (jnrlistA>=0) ? jnrlistA : 0;
655             jnrB       = (jnrlistB>=0) ? jnrlistB : 0;
656             jnrC       = (jnrlistC>=0) ? jnrlistC : 0;
657             jnrD       = (jnrlistD>=0) ? jnrlistD : 0;
658             j_coord_offsetA  = DIM*jnrA;
659             j_coord_offsetB  = DIM*jnrB;
660             j_coord_offsetC  = DIM*jnrC;
661             j_coord_offsetD  = DIM*jnrD;
662
663             /* load j atom coordinates */
664             gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
665                                               x+j_coord_offsetC,x+j_coord_offsetD,
666                                               &jx0,&jy0,&jz0);
667
668             /* Calculate displacement vector */
669             dx00             = _mm_sub_ps(ix0,jx0);
670             dy00             = _mm_sub_ps(iy0,jy0);
671             dz00             = _mm_sub_ps(iz0,jz0);
672
673             /* Calculate squared distance and things based on it */
674             rsq00            = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
675
676             rinv00           = gmx_mm_invsqrt_ps(rsq00);
677
678             rinvsq00         = _mm_mul_ps(rinv00,rinv00);
679
680             /* Load parameters for j particles */
681             jq0              = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
682                                                               charge+jnrC+0,charge+jnrD+0);
683             vdwjidx0A        = 2*vdwtype[jnrA+0];
684             vdwjidx0B        = 2*vdwtype[jnrB+0];
685             vdwjidx0C        = 2*vdwtype[jnrC+0];
686             vdwjidx0D        = 2*vdwtype[jnrD+0];
687
688             /**************************
689              * CALCULATE INTERACTIONS *
690              **************************/
691
692             if (gmx_mm_any_lt(rsq00,rcutoff2))
693             {
694
695             r00              = _mm_mul_ps(rsq00,rinv00);
696             r00              = _mm_andnot_ps(dummy_mask,r00);
697
698             /* Compute parameters for interactions between i and j atoms */
699             qq00             = _mm_mul_ps(iq0,jq0);
700             gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
701                                          vdwparam+vdwioffset0+vdwjidx0B,
702                                          vdwparam+vdwioffset0+vdwjidx0C,
703                                          vdwparam+vdwioffset0+vdwjidx0D,
704                                          &c6_00,&c12_00);
705
706             /* EWALD ELECTROSTATICS */
707
708             /* Analytical PME correction */
709             zeta2            = _mm_mul_ps(beta2,rsq00);
710             rinv3            = _mm_mul_ps(rinvsq00,rinv00);
711             pmecorrF         = gmx_mm_pmecorrF_ps(zeta2);
712             felec            = _mm_macc_ps(pmecorrF,beta3,rinv3);
713             felec            = _mm_mul_ps(qq00,felec);
714
715             /* LENNARD-JONES DISPERSION/REPULSION */
716
717             rinvsix          = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
718             fvdw             = _mm_mul_ps(_mm_msub_ps(c12_00,rinvsix,c6_00),_mm_mul_ps(rinvsix,rinvsq00));
719
720             cutoff_mask      = _mm_cmplt_ps(rsq00,rcutoff2);
721
722             fscal            = _mm_add_ps(felec,fvdw);
723
724             fscal            = _mm_and_ps(fscal,cutoff_mask);
725
726             fscal            = _mm_andnot_ps(dummy_mask,fscal);
727
728              /* Update vectorial force */
729             fix0             = _mm_macc_ps(dx00,fscal,fix0);
730             fiy0             = _mm_macc_ps(dy00,fscal,fiy0);
731             fiz0             = _mm_macc_ps(dz00,fscal,fiz0);
732
733             fjptrA             = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
734             fjptrB             = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
735             fjptrC             = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
736             fjptrD             = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
737             gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,
738                                                    _mm_mul_ps(dx00,fscal),
739                                                    _mm_mul_ps(dy00,fscal),
740                                                    _mm_mul_ps(dz00,fscal));
741
742             }
743
744             /* Inner loop uses 39 flops */
745         }
746
747         /* End of innermost loop */
748
749         gmx_mm_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
750                                               f+i_coord_offset,fshift+i_shift_offset);
751
752         /* Increment number of inner iterations */
753         inneriter                  += j_index_end - j_index_start;
754
755         /* Outer loop uses 7 flops */
756     }
757
758     /* Increment number of outer iterations */
759     outeriter        += nri;
760
761     /* Update outer/inner flops */
762
763     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_F,outeriter*7 + inneriter*39);
764 }