2 * Note: this file was generated by the Gromacs avx_128_fma_single kernel generator.
4 * This source code is part of
8 * Copyright (c) 2001-2012, The GROMACS Development Team
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
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
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.
28 #include "../nb_kernel.h"
29 #include "types/simple.h"
33 #include "gmx_math_x86_avx_128_fma_single.h"
34 #include "kernelutil_x86_avx_128_fma_single.h"
37 * Gromacs nonbonded kernel: nb_kernel_ElecEwSh_VdwNone_GeomW4W4_VF_avx_128_fma_single
38 * Electrostatics interaction: Ewald
39 * VdW interaction: None
40 * Geometry: Water4-Water4
41 * Calculate force/pot: PotentialAndForce
44 nb_kernel_ElecEwSh_VdwNone_GeomW4W4_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)
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.
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;
65 real *shiftvec,*fshift,*x,*f;
66 real *fjptrA,*fjptrB,*fjptrC,*fjptrD;
68 __m128 fscal,rcutoff,rcutoff2,jidxall;
70 __m128 ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
72 __m128 ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
74 __m128 ix3,iy3,iz3,fix3,fiy3,fiz3,iq3,isai3;
75 int vdwjidx1A,vdwjidx1B,vdwjidx1C,vdwjidx1D;
76 __m128 jx1,jy1,jz1,fjx1,fjy1,fjz1,jq1,isaj1;
77 int vdwjidx2A,vdwjidx2B,vdwjidx2C,vdwjidx2D;
78 __m128 jx2,jy2,jz2,fjx2,fjy2,fjz2,jq2,isaj2;
79 int vdwjidx3A,vdwjidx3B,vdwjidx3C,vdwjidx3D;
80 __m128 jx3,jy3,jz3,fjx3,fjy3,fjz3,jq3,isaj3;
81 __m128 dx11,dy11,dz11,rsq11,rinv11,rinvsq11,r11,qq11,c6_11,c12_11;
82 __m128 dx12,dy12,dz12,rsq12,rinv12,rinvsq12,r12,qq12,c6_12,c12_12;
83 __m128 dx13,dy13,dz13,rsq13,rinv13,rinvsq13,r13,qq13,c6_13,c12_13;
84 __m128 dx21,dy21,dz21,rsq21,rinv21,rinvsq21,r21,qq21,c6_21,c12_21;
85 __m128 dx22,dy22,dz22,rsq22,rinv22,rinvsq22,r22,qq22,c6_22,c12_22;
86 __m128 dx23,dy23,dz23,rsq23,rinv23,rinvsq23,r23,qq23,c6_23,c12_23;
87 __m128 dx31,dy31,dz31,rsq31,rinv31,rinvsq31,r31,qq31,c6_31,c12_31;
88 __m128 dx32,dy32,dz32,rsq32,rinv32,rinvsq32,r32,qq32,c6_32,c12_32;
89 __m128 dx33,dy33,dz33,rsq33,rinv33,rinvsq33,r33,qq33,c6_33,c12_33;
90 __m128 velec,felec,velecsum,facel,crf,krf,krf2;
93 __m128 ewtabscale,eweps,twoeweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
94 __m128 beta,beta2,beta3,zeta2,pmecorrF,pmecorrV,rinv3;
96 __m128 dummy_mask,cutoff_mask;
97 __m128 signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
98 __m128 one = _mm_set1_ps(1.0);
99 __m128 two = _mm_set1_ps(2.0);
105 jindex = nlist->jindex;
107 shiftidx = nlist->shift;
109 shiftvec = fr->shift_vec[0];
110 fshift = fr->fshift[0];
111 facel = _mm_set1_ps(fr->epsfac);
112 charge = mdatoms->chargeA;
114 sh_ewald = _mm_set1_ps(fr->ic->sh_ewald);
115 beta = _mm_set1_ps(fr->ic->ewaldcoeff);
116 beta2 = _mm_mul_ps(beta,beta);
117 beta3 = _mm_mul_ps(beta,beta2);
118 ewtab = fr->ic->tabq_coul_FDV0;
119 ewtabscale = _mm_set1_ps(fr->ic->tabq_scale);
120 ewtabhalfspace = _mm_set1_ps(0.5/fr->ic->tabq_scale);
122 /* Setup water-specific parameters */
123 inr = nlist->iinr[0];
124 iq1 = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+1]));
125 iq2 = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+2]));
126 iq3 = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+3]));
128 jq1 = _mm_set1_ps(charge[inr+1]);
129 jq2 = _mm_set1_ps(charge[inr+2]);
130 jq3 = _mm_set1_ps(charge[inr+3]);
131 qq11 = _mm_mul_ps(iq1,jq1);
132 qq12 = _mm_mul_ps(iq1,jq2);
133 qq13 = _mm_mul_ps(iq1,jq3);
134 qq21 = _mm_mul_ps(iq2,jq1);
135 qq22 = _mm_mul_ps(iq2,jq2);
136 qq23 = _mm_mul_ps(iq2,jq3);
137 qq31 = _mm_mul_ps(iq3,jq1);
138 qq32 = _mm_mul_ps(iq3,jq2);
139 qq33 = _mm_mul_ps(iq3,jq3);
141 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
142 rcutoff_scalar = fr->rcoulomb;
143 rcutoff = _mm_set1_ps(rcutoff_scalar);
144 rcutoff2 = _mm_mul_ps(rcutoff,rcutoff);
146 /* Avoid stupid compiler warnings */
147 jnrA = jnrB = jnrC = jnrD = 0;
156 for(iidx=0;iidx<4*DIM;iidx++)
161 /* Start outer loop over neighborlists */
162 for(iidx=0; iidx<nri; iidx++)
164 /* Load shift vector for this list */
165 i_shift_offset = DIM*shiftidx[iidx];
167 /* Load limits for loop over neighbors */
168 j_index_start = jindex[iidx];
169 j_index_end = jindex[iidx+1];
171 /* Get outer coordinate index */
173 i_coord_offset = DIM*inr;
175 /* Load i particle coords and add shift vector */
176 gmx_mm_load_shift_and_3rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset+DIM,
177 &ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
179 fix1 = _mm_setzero_ps();
180 fiy1 = _mm_setzero_ps();
181 fiz1 = _mm_setzero_ps();
182 fix2 = _mm_setzero_ps();
183 fiy2 = _mm_setzero_ps();
184 fiz2 = _mm_setzero_ps();
185 fix3 = _mm_setzero_ps();
186 fiy3 = _mm_setzero_ps();
187 fiz3 = _mm_setzero_ps();
189 /* Reset potential sums */
190 velecsum = _mm_setzero_ps();
192 /* Start inner kernel loop */
193 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
196 /* Get j neighbor index, and coordinate index */
201 j_coord_offsetA = DIM*jnrA;
202 j_coord_offsetB = DIM*jnrB;
203 j_coord_offsetC = DIM*jnrC;
204 j_coord_offsetD = DIM*jnrD;
206 /* load j atom coordinates */
207 gmx_mm_load_3rvec_4ptr_swizzle_ps(x+j_coord_offsetA+DIM,x+j_coord_offsetB+DIM,
208 x+j_coord_offsetC+DIM,x+j_coord_offsetD+DIM,
209 &jx1,&jy1,&jz1,&jx2,&jy2,&jz2,&jx3,&jy3,&jz3);
211 /* Calculate displacement vector */
212 dx11 = _mm_sub_ps(ix1,jx1);
213 dy11 = _mm_sub_ps(iy1,jy1);
214 dz11 = _mm_sub_ps(iz1,jz1);
215 dx12 = _mm_sub_ps(ix1,jx2);
216 dy12 = _mm_sub_ps(iy1,jy2);
217 dz12 = _mm_sub_ps(iz1,jz2);
218 dx13 = _mm_sub_ps(ix1,jx3);
219 dy13 = _mm_sub_ps(iy1,jy3);
220 dz13 = _mm_sub_ps(iz1,jz3);
221 dx21 = _mm_sub_ps(ix2,jx1);
222 dy21 = _mm_sub_ps(iy2,jy1);
223 dz21 = _mm_sub_ps(iz2,jz1);
224 dx22 = _mm_sub_ps(ix2,jx2);
225 dy22 = _mm_sub_ps(iy2,jy2);
226 dz22 = _mm_sub_ps(iz2,jz2);
227 dx23 = _mm_sub_ps(ix2,jx3);
228 dy23 = _mm_sub_ps(iy2,jy3);
229 dz23 = _mm_sub_ps(iz2,jz3);
230 dx31 = _mm_sub_ps(ix3,jx1);
231 dy31 = _mm_sub_ps(iy3,jy1);
232 dz31 = _mm_sub_ps(iz3,jz1);
233 dx32 = _mm_sub_ps(ix3,jx2);
234 dy32 = _mm_sub_ps(iy3,jy2);
235 dz32 = _mm_sub_ps(iz3,jz2);
236 dx33 = _mm_sub_ps(ix3,jx3);
237 dy33 = _mm_sub_ps(iy3,jy3);
238 dz33 = _mm_sub_ps(iz3,jz3);
240 /* Calculate squared distance and things based on it */
241 rsq11 = gmx_mm_calc_rsq_ps(dx11,dy11,dz11);
242 rsq12 = gmx_mm_calc_rsq_ps(dx12,dy12,dz12);
243 rsq13 = gmx_mm_calc_rsq_ps(dx13,dy13,dz13);
244 rsq21 = gmx_mm_calc_rsq_ps(dx21,dy21,dz21);
245 rsq22 = gmx_mm_calc_rsq_ps(dx22,dy22,dz22);
246 rsq23 = gmx_mm_calc_rsq_ps(dx23,dy23,dz23);
247 rsq31 = gmx_mm_calc_rsq_ps(dx31,dy31,dz31);
248 rsq32 = gmx_mm_calc_rsq_ps(dx32,dy32,dz32);
249 rsq33 = gmx_mm_calc_rsq_ps(dx33,dy33,dz33);
251 rinv11 = gmx_mm_invsqrt_ps(rsq11);
252 rinv12 = gmx_mm_invsqrt_ps(rsq12);
253 rinv13 = gmx_mm_invsqrt_ps(rsq13);
254 rinv21 = gmx_mm_invsqrt_ps(rsq21);
255 rinv22 = gmx_mm_invsqrt_ps(rsq22);
256 rinv23 = gmx_mm_invsqrt_ps(rsq23);
257 rinv31 = gmx_mm_invsqrt_ps(rsq31);
258 rinv32 = gmx_mm_invsqrt_ps(rsq32);
259 rinv33 = gmx_mm_invsqrt_ps(rsq33);
261 rinvsq11 = _mm_mul_ps(rinv11,rinv11);
262 rinvsq12 = _mm_mul_ps(rinv12,rinv12);
263 rinvsq13 = _mm_mul_ps(rinv13,rinv13);
264 rinvsq21 = _mm_mul_ps(rinv21,rinv21);
265 rinvsq22 = _mm_mul_ps(rinv22,rinv22);
266 rinvsq23 = _mm_mul_ps(rinv23,rinv23);
267 rinvsq31 = _mm_mul_ps(rinv31,rinv31);
268 rinvsq32 = _mm_mul_ps(rinv32,rinv32);
269 rinvsq33 = _mm_mul_ps(rinv33,rinv33);
271 fjx1 = _mm_setzero_ps();
272 fjy1 = _mm_setzero_ps();
273 fjz1 = _mm_setzero_ps();
274 fjx2 = _mm_setzero_ps();
275 fjy2 = _mm_setzero_ps();
276 fjz2 = _mm_setzero_ps();
277 fjx3 = _mm_setzero_ps();
278 fjy3 = _mm_setzero_ps();
279 fjz3 = _mm_setzero_ps();
281 /**************************
282 * CALCULATE INTERACTIONS *
283 **************************/
285 if (gmx_mm_any_lt(rsq11,rcutoff2))
288 r11 = _mm_mul_ps(rsq11,rinv11);
290 /* EWALD ELECTROSTATICS */
292 /* Analytical PME correction */
293 zeta2 = _mm_mul_ps(beta2,rsq11);
294 rinv3 = _mm_mul_ps(rinvsq11,rinv11);
295 pmecorrF = gmx_mm_pmecorrF_ps(zeta2);
296 felec = _mm_macc_ps(pmecorrF,beta3,rinv3);
297 felec = _mm_mul_ps(qq11,felec);
298 pmecorrV = gmx_mm_pmecorrV_ps(zeta2);
299 velec = _mm_nmacc_ps(pmecorrV,beta,_mm_sub_ps(rinv11,sh_ewald));
300 velec = _mm_mul_ps(qq11,velec);
302 cutoff_mask = _mm_cmplt_ps(rsq11,rcutoff2);
304 /* Update potential sum for this i atom from the interaction with this j atom. */
305 velec = _mm_and_ps(velec,cutoff_mask);
306 velecsum = _mm_add_ps(velecsum,velec);
310 fscal = _mm_and_ps(fscal,cutoff_mask);
312 /* Update vectorial force */
313 fix1 = _mm_macc_ps(dx11,fscal,fix1);
314 fiy1 = _mm_macc_ps(dy11,fscal,fiy1);
315 fiz1 = _mm_macc_ps(dz11,fscal,fiz1);
317 fjx1 = _mm_macc_ps(dx11,fscal,fjx1);
318 fjy1 = _mm_macc_ps(dy11,fscal,fjy1);
319 fjz1 = _mm_macc_ps(dz11,fscal,fjz1);
323 /**************************
324 * CALCULATE INTERACTIONS *
325 **************************/
327 if (gmx_mm_any_lt(rsq12,rcutoff2))
330 r12 = _mm_mul_ps(rsq12,rinv12);
332 /* EWALD ELECTROSTATICS */
334 /* Analytical PME correction */
335 zeta2 = _mm_mul_ps(beta2,rsq12);
336 rinv3 = _mm_mul_ps(rinvsq12,rinv12);
337 pmecorrF = gmx_mm_pmecorrF_ps(zeta2);
338 felec = _mm_macc_ps(pmecorrF,beta3,rinv3);
339 felec = _mm_mul_ps(qq12,felec);
340 pmecorrV = gmx_mm_pmecorrV_ps(zeta2);
341 velec = _mm_nmacc_ps(pmecorrV,beta,_mm_sub_ps(rinv12,sh_ewald));
342 velec = _mm_mul_ps(qq12,velec);
344 cutoff_mask = _mm_cmplt_ps(rsq12,rcutoff2);
346 /* Update potential sum for this i atom from the interaction with this j atom. */
347 velec = _mm_and_ps(velec,cutoff_mask);
348 velecsum = _mm_add_ps(velecsum,velec);
352 fscal = _mm_and_ps(fscal,cutoff_mask);
354 /* Update vectorial force */
355 fix1 = _mm_macc_ps(dx12,fscal,fix1);
356 fiy1 = _mm_macc_ps(dy12,fscal,fiy1);
357 fiz1 = _mm_macc_ps(dz12,fscal,fiz1);
359 fjx2 = _mm_macc_ps(dx12,fscal,fjx2);
360 fjy2 = _mm_macc_ps(dy12,fscal,fjy2);
361 fjz2 = _mm_macc_ps(dz12,fscal,fjz2);
365 /**************************
366 * CALCULATE INTERACTIONS *
367 **************************/
369 if (gmx_mm_any_lt(rsq13,rcutoff2))
372 r13 = _mm_mul_ps(rsq13,rinv13);
374 /* EWALD ELECTROSTATICS */
376 /* Analytical PME correction */
377 zeta2 = _mm_mul_ps(beta2,rsq13);
378 rinv3 = _mm_mul_ps(rinvsq13,rinv13);
379 pmecorrF = gmx_mm_pmecorrF_ps(zeta2);
380 felec = _mm_macc_ps(pmecorrF,beta3,rinv3);
381 felec = _mm_mul_ps(qq13,felec);
382 pmecorrV = gmx_mm_pmecorrV_ps(zeta2);
383 velec = _mm_nmacc_ps(pmecorrV,beta,_mm_sub_ps(rinv13,sh_ewald));
384 velec = _mm_mul_ps(qq13,velec);
386 cutoff_mask = _mm_cmplt_ps(rsq13,rcutoff2);
388 /* Update potential sum for this i atom from the interaction with this j atom. */
389 velec = _mm_and_ps(velec,cutoff_mask);
390 velecsum = _mm_add_ps(velecsum,velec);
394 fscal = _mm_and_ps(fscal,cutoff_mask);
396 /* Update vectorial force */
397 fix1 = _mm_macc_ps(dx13,fscal,fix1);
398 fiy1 = _mm_macc_ps(dy13,fscal,fiy1);
399 fiz1 = _mm_macc_ps(dz13,fscal,fiz1);
401 fjx3 = _mm_macc_ps(dx13,fscal,fjx3);
402 fjy3 = _mm_macc_ps(dy13,fscal,fjy3);
403 fjz3 = _mm_macc_ps(dz13,fscal,fjz3);
407 /**************************
408 * CALCULATE INTERACTIONS *
409 **************************/
411 if (gmx_mm_any_lt(rsq21,rcutoff2))
414 r21 = _mm_mul_ps(rsq21,rinv21);
416 /* EWALD ELECTROSTATICS */
418 /* Analytical PME correction */
419 zeta2 = _mm_mul_ps(beta2,rsq21);
420 rinv3 = _mm_mul_ps(rinvsq21,rinv21);
421 pmecorrF = gmx_mm_pmecorrF_ps(zeta2);
422 felec = _mm_macc_ps(pmecorrF,beta3,rinv3);
423 felec = _mm_mul_ps(qq21,felec);
424 pmecorrV = gmx_mm_pmecorrV_ps(zeta2);
425 velec = _mm_nmacc_ps(pmecorrV,beta,_mm_sub_ps(rinv21,sh_ewald));
426 velec = _mm_mul_ps(qq21,velec);
428 cutoff_mask = _mm_cmplt_ps(rsq21,rcutoff2);
430 /* Update potential sum for this i atom from the interaction with this j atom. */
431 velec = _mm_and_ps(velec,cutoff_mask);
432 velecsum = _mm_add_ps(velecsum,velec);
436 fscal = _mm_and_ps(fscal,cutoff_mask);
438 /* Update vectorial force */
439 fix2 = _mm_macc_ps(dx21,fscal,fix2);
440 fiy2 = _mm_macc_ps(dy21,fscal,fiy2);
441 fiz2 = _mm_macc_ps(dz21,fscal,fiz2);
443 fjx1 = _mm_macc_ps(dx21,fscal,fjx1);
444 fjy1 = _mm_macc_ps(dy21,fscal,fjy1);
445 fjz1 = _mm_macc_ps(dz21,fscal,fjz1);
449 /**************************
450 * CALCULATE INTERACTIONS *
451 **************************/
453 if (gmx_mm_any_lt(rsq22,rcutoff2))
456 r22 = _mm_mul_ps(rsq22,rinv22);
458 /* EWALD ELECTROSTATICS */
460 /* Analytical PME correction */
461 zeta2 = _mm_mul_ps(beta2,rsq22);
462 rinv3 = _mm_mul_ps(rinvsq22,rinv22);
463 pmecorrF = gmx_mm_pmecorrF_ps(zeta2);
464 felec = _mm_macc_ps(pmecorrF,beta3,rinv3);
465 felec = _mm_mul_ps(qq22,felec);
466 pmecorrV = gmx_mm_pmecorrV_ps(zeta2);
467 velec = _mm_nmacc_ps(pmecorrV,beta,_mm_sub_ps(rinv22,sh_ewald));
468 velec = _mm_mul_ps(qq22,velec);
470 cutoff_mask = _mm_cmplt_ps(rsq22,rcutoff2);
472 /* Update potential sum for this i atom from the interaction with this j atom. */
473 velec = _mm_and_ps(velec,cutoff_mask);
474 velecsum = _mm_add_ps(velecsum,velec);
478 fscal = _mm_and_ps(fscal,cutoff_mask);
480 /* Update vectorial force */
481 fix2 = _mm_macc_ps(dx22,fscal,fix2);
482 fiy2 = _mm_macc_ps(dy22,fscal,fiy2);
483 fiz2 = _mm_macc_ps(dz22,fscal,fiz2);
485 fjx2 = _mm_macc_ps(dx22,fscal,fjx2);
486 fjy2 = _mm_macc_ps(dy22,fscal,fjy2);
487 fjz2 = _mm_macc_ps(dz22,fscal,fjz2);
491 /**************************
492 * CALCULATE INTERACTIONS *
493 **************************/
495 if (gmx_mm_any_lt(rsq23,rcutoff2))
498 r23 = _mm_mul_ps(rsq23,rinv23);
500 /* EWALD ELECTROSTATICS */
502 /* Analytical PME correction */
503 zeta2 = _mm_mul_ps(beta2,rsq23);
504 rinv3 = _mm_mul_ps(rinvsq23,rinv23);
505 pmecorrF = gmx_mm_pmecorrF_ps(zeta2);
506 felec = _mm_macc_ps(pmecorrF,beta3,rinv3);
507 felec = _mm_mul_ps(qq23,felec);
508 pmecorrV = gmx_mm_pmecorrV_ps(zeta2);
509 velec = _mm_nmacc_ps(pmecorrV,beta,_mm_sub_ps(rinv23,sh_ewald));
510 velec = _mm_mul_ps(qq23,velec);
512 cutoff_mask = _mm_cmplt_ps(rsq23,rcutoff2);
514 /* Update potential sum for this i atom from the interaction with this j atom. */
515 velec = _mm_and_ps(velec,cutoff_mask);
516 velecsum = _mm_add_ps(velecsum,velec);
520 fscal = _mm_and_ps(fscal,cutoff_mask);
522 /* Update vectorial force */
523 fix2 = _mm_macc_ps(dx23,fscal,fix2);
524 fiy2 = _mm_macc_ps(dy23,fscal,fiy2);
525 fiz2 = _mm_macc_ps(dz23,fscal,fiz2);
527 fjx3 = _mm_macc_ps(dx23,fscal,fjx3);
528 fjy3 = _mm_macc_ps(dy23,fscal,fjy3);
529 fjz3 = _mm_macc_ps(dz23,fscal,fjz3);
533 /**************************
534 * CALCULATE INTERACTIONS *
535 **************************/
537 if (gmx_mm_any_lt(rsq31,rcutoff2))
540 r31 = _mm_mul_ps(rsq31,rinv31);
542 /* EWALD ELECTROSTATICS */
544 /* Analytical PME correction */
545 zeta2 = _mm_mul_ps(beta2,rsq31);
546 rinv3 = _mm_mul_ps(rinvsq31,rinv31);
547 pmecorrF = gmx_mm_pmecorrF_ps(zeta2);
548 felec = _mm_macc_ps(pmecorrF,beta3,rinv3);
549 felec = _mm_mul_ps(qq31,felec);
550 pmecorrV = gmx_mm_pmecorrV_ps(zeta2);
551 velec = _mm_nmacc_ps(pmecorrV,beta,_mm_sub_ps(rinv31,sh_ewald));
552 velec = _mm_mul_ps(qq31,velec);
554 cutoff_mask = _mm_cmplt_ps(rsq31,rcutoff2);
556 /* Update potential sum for this i atom from the interaction with this j atom. */
557 velec = _mm_and_ps(velec,cutoff_mask);
558 velecsum = _mm_add_ps(velecsum,velec);
562 fscal = _mm_and_ps(fscal,cutoff_mask);
564 /* Update vectorial force */
565 fix3 = _mm_macc_ps(dx31,fscal,fix3);
566 fiy3 = _mm_macc_ps(dy31,fscal,fiy3);
567 fiz3 = _mm_macc_ps(dz31,fscal,fiz3);
569 fjx1 = _mm_macc_ps(dx31,fscal,fjx1);
570 fjy1 = _mm_macc_ps(dy31,fscal,fjy1);
571 fjz1 = _mm_macc_ps(dz31,fscal,fjz1);
575 /**************************
576 * CALCULATE INTERACTIONS *
577 **************************/
579 if (gmx_mm_any_lt(rsq32,rcutoff2))
582 r32 = _mm_mul_ps(rsq32,rinv32);
584 /* EWALD ELECTROSTATICS */
586 /* Analytical PME correction */
587 zeta2 = _mm_mul_ps(beta2,rsq32);
588 rinv3 = _mm_mul_ps(rinvsq32,rinv32);
589 pmecorrF = gmx_mm_pmecorrF_ps(zeta2);
590 felec = _mm_macc_ps(pmecorrF,beta3,rinv3);
591 felec = _mm_mul_ps(qq32,felec);
592 pmecorrV = gmx_mm_pmecorrV_ps(zeta2);
593 velec = _mm_nmacc_ps(pmecorrV,beta,_mm_sub_ps(rinv32,sh_ewald));
594 velec = _mm_mul_ps(qq32,velec);
596 cutoff_mask = _mm_cmplt_ps(rsq32,rcutoff2);
598 /* Update potential sum for this i atom from the interaction with this j atom. */
599 velec = _mm_and_ps(velec,cutoff_mask);
600 velecsum = _mm_add_ps(velecsum,velec);
604 fscal = _mm_and_ps(fscal,cutoff_mask);
606 /* Update vectorial force */
607 fix3 = _mm_macc_ps(dx32,fscal,fix3);
608 fiy3 = _mm_macc_ps(dy32,fscal,fiy3);
609 fiz3 = _mm_macc_ps(dz32,fscal,fiz3);
611 fjx2 = _mm_macc_ps(dx32,fscal,fjx2);
612 fjy2 = _mm_macc_ps(dy32,fscal,fjy2);
613 fjz2 = _mm_macc_ps(dz32,fscal,fjz2);
617 /**************************
618 * CALCULATE INTERACTIONS *
619 **************************/
621 if (gmx_mm_any_lt(rsq33,rcutoff2))
624 r33 = _mm_mul_ps(rsq33,rinv33);
626 /* EWALD ELECTROSTATICS */
628 /* Analytical PME correction */
629 zeta2 = _mm_mul_ps(beta2,rsq33);
630 rinv3 = _mm_mul_ps(rinvsq33,rinv33);
631 pmecorrF = gmx_mm_pmecorrF_ps(zeta2);
632 felec = _mm_macc_ps(pmecorrF,beta3,rinv3);
633 felec = _mm_mul_ps(qq33,felec);
634 pmecorrV = gmx_mm_pmecorrV_ps(zeta2);
635 velec = _mm_nmacc_ps(pmecorrV,beta,_mm_sub_ps(rinv33,sh_ewald));
636 velec = _mm_mul_ps(qq33,velec);
638 cutoff_mask = _mm_cmplt_ps(rsq33,rcutoff2);
640 /* Update potential sum for this i atom from the interaction with this j atom. */
641 velec = _mm_and_ps(velec,cutoff_mask);
642 velecsum = _mm_add_ps(velecsum,velec);
646 fscal = _mm_and_ps(fscal,cutoff_mask);
648 /* Update vectorial force */
649 fix3 = _mm_macc_ps(dx33,fscal,fix3);
650 fiy3 = _mm_macc_ps(dy33,fscal,fiy3);
651 fiz3 = _mm_macc_ps(dz33,fscal,fiz3);
653 fjx3 = _mm_macc_ps(dx33,fscal,fjx3);
654 fjy3 = _mm_macc_ps(dy33,fscal,fjy3);
655 fjz3 = _mm_macc_ps(dz33,fscal,fjz3);
659 fjptrA = f+j_coord_offsetA;
660 fjptrB = f+j_coord_offsetB;
661 fjptrC = f+j_coord_offsetC;
662 fjptrD = f+j_coord_offsetD;
664 gmx_mm_decrement_3rvec_4ptr_swizzle_ps(fjptrA+DIM,fjptrB+DIM,fjptrC+DIM,fjptrD+DIM,
665 fjx1,fjy1,fjz1,fjx2,fjy2,fjz2,fjx3,fjy3,fjz3);
667 /* Inner loop uses 297 flops */
673 /* Get j neighbor index, and coordinate index */
674 jnrlistA = jjnr[jidx];
675 jnrlistB = jjnr[jidx+1];
676 jnrlistC = jjnr[jidx+2];
677 jnrlistD = jjnr[jidx+3];
678 /* Sign of each element will be negative for non-real atoms.
679 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
680 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
682 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
683 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
684 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
685 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
686 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
687 j_coord_offsetA = DIM*jnrA;
688 j_coord_offsetB = DIM*jnrB;
689 j_coord_offsetC = DIM*jnrC;
690 j_coord_offsetD = DIM*jnrD;
692 /* load j atom coordinates */
693 gmx_mm_load_3rvec_4ptr_swizzle_ps(x+j_coord_offsetA+DIM,x+j_coord_offsetB+DIM,
694 x+j_coord_offsetC+DIM,x+j_coord_offsetD+DIM,
695 &jx1,&jy1,&jz1,&jx2,&jy2,&jz2,&jx3,&jy3,&jz3);
697 /* Calculate displacement vector */
698 dx11 = _mm_sub_ps(ix1,jx1);
699 dy11 = _mm_sub_ps(iy1,jy1);
700 dz11 = _mm_sub_ps(iz1,jz1);
701 dx12 = _mm_sub_ps(ix1,jx2);
702 dy12 = _mm_sub_ps(iy1,jy2);
703 dz12 = _mm_sub_ps(iz1,jz2);
704 dx13 = _mm_sub_ps(ix1,jx3);
705 dy13 = _mm_sub_ps(iy1,jy3);
706 dz13 = _mm_sub_ps(iz1,jz3);
707 dx21 = _mm_sub_ps(ix2,jx1);
708 dy21 = _mm_sub_ps(iy2,jy1);
709 dz21 = _mm_sub_ps(iz2,jz1);
710 dx22 = _mm_sub_ps(ix2,jx2);
711 dy22 = _mm_sub_ps(iy2,jy2);
712 dz22 = _mm_sub_ps(iz2,jz2);
713 dx23 = _mm_sub_ps(ix2,jx3);
714 dy23 = _mm_sub_ps(iy2,jy3);
715 dz23 = _mm_sub_ps(iz2,jz3);
716 dx31 = _mm_sub_ps(ix3,jx1);
717 dy31 = _mm_sub_ps(iy3,jy1);
718 dz31 = _mm_sub_ps(iz3,jz1);
719 dx32 = _mm_sub_ps(ix3,jx2);
720 dy32 = _mm_sub_ps(iy3,jy2);
721 dz32 = _mm_sub_ps(iz3,jz2);
722 dx33 = _mm_sub_ps(ix3,jx3);
723 dy33 = _mm_sub_ps(iy3,jy3);
724 dz33 = _mm_sub_ps(iz3,jz3);
726 /* Calculate squared distance and things based on it */
727 rsq11 = gmx_mm_calc_rsq_ps(dx11,dy11,dz11);
728 rsq12 = gmx_mm_calc_rsq_ps(dx12,dy12,dz12);
729 rsq13 = gmx_mm_calc_rsq_ps(dx13,dy13,dz13);
730 rsq21 = gmx_mm_calc_rsq_ps(dx21,dy21,dz21);
731 rsq22 = gmx_mm_calc_rsq_ps(dx22,dy22,dz22);
732 rsq23 = gmx_mm_calc_rsq_ps(dx23,dy23,dz23);
733 rsq31 = gmx_mm_calc_rsq_ps(dx31,dy31,dz31);
734 rsq32 = gmx_mm_calc_rsq_ps(dx32,dy32,dz32);
735 rsq33 = gmx_mm_calc_rsq_ps(dx33,dy33,dz33);
737 rinv11 = gmx_mm_invsqrt_ps(rsq11);
738 rinv12 = gmx_mm_invsqrt_ps(rsq12);
739 rinv13 = gmx_mm_invsqrt_ps(rsq13);
740 rinv21 = gmx_mm_invsqrt_ps(rsq21);
741 rinv22 = gmx_mm_invsqrt_ps(rsq22);
742 rinv23 = gmx_mm_invsqrt_ps(rsq23);
743 rinv31 = gmx_mm_invsqrt_ps(rsq31);
744 rinv32 = gmx_mm_invsqrt_ps(rsq32);
745 rinv33 = gmx_mm_invsqrt_ps(rsq33);
747 rinvsq11 = _mm_mul_ps(rinv11,rinv11);
748 rinvsq12 = _mm_mul_ps(rinv12,rinv12);
749 rinvsq13 = _mm_mul_ps(rinv13,rinv13);
750 rinvsq21 = _mm_mul_ps(rinv21,rinv21);
751 rinvsq22 = _mm_mul_ps(rinv22,rinv22);
752 rinvsq23 = _mm_mul_ps(rinv23,rinv23);
753 rinvsq31 = _mm_mul_ps(rinv31,rinv31);
754 rinvsq32 = _mm_mul_ps(rinv32,rinv32);
755 rinvsq33 = _mm_mul_ps(rinv33,rinv33);
757 fjx1 = _mm_setzero_ps();
758 fjy1 = _mm_setzero_ps();
759 fjz1 = _mm_setzero_ps();
760 fjx2 = _mm_setzero_ps();
761 fjy2 = _mm_setzero_ps();
762 fjz2 = _mm_setzero_ps();
763 fjx3 = _mm_setzero_ps();
764 fjy3 = _mm_setzero_ps();
765 fjz3 = _mm_setzero_ps();
767 /**************************
768 * CALCULATE INTERACTIONS *
769 **************************/
771 if (gmx_mm_any_lt(rsq11,rcutoff2))
774 r11 = _mm_mul_ps(rsq11,rinv11);
775 r11 = _mm_andnot_ps(dummy_mask,r11);
777 /* EWALD ELECTROSTATICS */
779 /* Analytical PME correction */
780 zeta2 = _mm_mul_ps(beta2,rsq11);
781 rinv3 = _mm_mul_ps(rinvsq11,rinv11);
782 pmecorrF = gmx_mm_pmecorrF_ps(zeta2);
783 felec = _mm_macc_ps(pmecorrF,beta3,rinv3);
784 felec = _mm_mul_ps(qq11,felec);
785 pmecorrV = gmx_mm_pmecorrV_ps(zeta2);
786 velec = _mm_nmacc_ps(pmecorrV,beta,_mm_sub_ps(rinv11,sh_ewald));
787 velec = _mm_mul_ps(qq11,velec);
789 cutoff_mask = _mm_cmplt_ps(rsq11,rcutoff2);
791 /* Update potential sum for this i atom from the interaction with this j atom. */
792 velec = _mm_and_ps(velec,cutoff_mask);
793 velec = _mm_andnot_ps(dummy_mask,velec);
794 velecsum = _mm_add_ps(velecsum,velec);
798 fscal = _mm_and_ps(fscal,cutoff_mask);
800 fscal = _mm_andnot_ps(dummy_mask,fscal);
802 /* Update vectorial force */
803 fix1 = _mm_macc_ps(dx11,fscal,fix1);
804 fiy1 = _mm_macc_ps(dy11,fscal,fiy1);
805 fiz1 = _mm_macc_ps(dz11,fscal,fiz1);
807 fjx1 = _mm_macc_ps(dx11,fscal,fjx1);
808 fjy1 = _mm_macc_ps(dy11,fscal,fjy1);
809 fjz1 = _mm_macc_ps(dz11,fscal,fjz1);
813 /**************************
814 * CALCULATE INTERACTIONS *
815 **************************/
817 if (gmx_mm_any_lt(rsq12,rcutoff2))
820 r12 = _mm_mul_ps(rsq12,rinv12);
821 r12 = _mm_andnot_ps(dummy_mask,r12);
823 /* EWALD ELECTROSTATICS */
825 /* Analytical PME correction */
826 zeta2 = _mm_mul_ps(beta2,rsq12);
827 rinv3 = _mm_mul_ps(rinvsq12,rinv12);
828 pmecorrF = gmx_mm_pmecorrF_ps(zeta2);
829 felec = _mm_macc_ps(pmecorrF,beta3,rinv3);
830 felec = _mm_mul_ps(qq12,felec);
831 pmecorrV = gmx_mm_pmecorrV_ps(zeta2);
832 velec = _mm_nmacc_ps(pmecorrV,beta,_mm_sub_ps(rinv12,sh_ewald));
833 velec = _mm_mul_ps(qq12,velec);
835 cutoff_mask = _mm_cmplt_ps(rsq12,rcutoff2);
837 /* Update potential sum for this i atom from the interaction with this j atom. */
838 velec = _mm_and_ps(velec,cutoff_mask);
839 velec = _mm_andnot_ps(dummy_mask,velec);
840 velecsum = _mm_add_ps(velecsum,velec);
844 fscal = _mm_and_ps(fscal,cutoff_mask);
846 fscal = _mm_andnot_ps(dummy_mask,fscal);
848 /* Update vectorial force */
849 fix1 = _mm_macc_ps(dx12,fscal,fix1);
850 fiy1 = _mm_macc_ps(dy12,fscal,fiy1);
851 fiz1 = _mm_macc_ps(dz12,fscal,fiz1);
853 fjx2 = _mm_macc_ps(dx12,fscal,fjx2);
854 fjy2 = _mm_macc_ps(dy12,fscal,fjy2);
855 fjz2 = _mm_macc_ps(dz12,fscal,fjz2);
859 /**************************
860 * CALCULATE INTERACTIONS *
861 **************************/
863 if (gmx_mm_any_lt(rsq13,rcutoff2))
866 r13 = _mm_mul_ps(rsq13,rinv13);
867 r13 = _mm_andnot_ps(dummy_mask,r13);
869 /* EWALD ELECTROSTATICS */
871 /* Analytical PME correction */
872 zeta2 = _mm_mul_ps(beta2,rsq13);
873 rinv3 = _mm_mul_ps(rinvsq13,rinv13);
874 pmecorrF = gmx_mm_pmecorrF_ps(zeta2);
875 felec = _mm_macc_ps(pmecorrF,beta3,rinv3);
876 felec = _mm_mul_ps(qq13,felec);
877 pmecorrV = gmx_mm_pmecorrV_ps(zeta2);
878 velec = _mm_nmacc_ps(pmecorrV,beta,_mm_sub_ps(rinv13,sh_ewald));
879 velec = _mm_mul_ps(qq13,velec);
881 cutoff_mask = _mm_cmplt_ps(rsq13,rcutoff2);
883 /* Update potential sum for this i atom from the interaction with this j atom. */
884 velec = _mm_and_ps(velec,cutoff_mask);
885 velec = _mm_andnot_ps(dummy_mask,velec);
886 velecsum = _mm_add_ps(velecsum,velec);
890 fscal = _mm_and_ps(fscal,cutoff_mask);
892 fscal = _mm_andnot_ps(dummy_mask,fscal);
894 /* Update vectorial force */
895 fix1 = _mm_macc_ps(dx13,fscal,fix1);
896 fiy1 = _mm_macc_ps(dy13,fscal,fiy1);
897 fiz1 = _mm_macc_ps(dz13,fscal,fiz1);
899 fjx3 = _mm_macc_ps(dx13,fscal,fjx3);
900 fjy3 = _mm_macc_ps(dy13,fscal,fjy3);
901 fjz3 = _mm_macc_ps(dz13,fscal,fjz3);
905 /**************************
906 * CALCULATE INTERACTIONS *
907 **************************/
909 if (gmx_mm_any_lt(rsq21,rcutoff2))
912 r21 = _mm_mul_ps(rsq21,rinv21);
913 r21 = _mm_andnot_ps(dummy_mask,r21);
915 /* EWALD ELECTROSTATICS */
917 /* Analytical PME correction */
918 zeta2 = _mm_mul_ps(beta2,rsq21);
919 rinv3 = _mm_mul_ps(rinvsq21,rinv21);
920 pmecorrF = gmx_mm_pmecorrF_ps(zeta2);
921 felec = _mm_macc_ps(pmecorrF,beta3,rinv3);
922 felec = _mm_mul_ps(qq21,felec);
923 pmecorrV = gmx_mm_pmecorrV_ps(zeta2);
924 velec = _mm_nmacc_ps(pmecorrV,beta,_mm_sub_ps(rinv21,sh_ewald));
925 velec = _mm_mul_ps(qq21,velec);
927 cutoff_mask = _mm_cmplt_ps(rsq21,rcutoff2);
929 /* Update potential sum for this i atom from the interaction with this j atom. */
930 velec = _mm_and_ps(velec,cutoff_mask);
931 velec = _mm_andnot_ps(dummy_mask,velec);
932 velecsum = _mm_add_ps(velecsum,velec);
936 fscal = _mm_and_ps(fscal,cutoff_mask);
938 fscal = _mm_andnot_ps(dummy_mask,fscal);
940 /* Update vectorial force */
941 fix2 = _mm_macc_ps(dx21,fscal,fix2);
942 fiy2 = _mm_macc_ps(dy21,fscal,fiy2);
943 fiz2 = _mm_macc_ps(dz21,fscal,fiz2);
945 fjx1 = _mm_macc_ps(dx21,fscal,fjx1);
946 fjy1 = _mm_macc_ps(dy21,fscal,fjy1);
947 fjz1 = _mm_macc_ps(dz21,fscal,fjz1);
951 /**************************
952 * CALCULATE INTERACTIONS *
953 **************************/
955 if (gmx_mm_any_lt(rsq22,rcutoff2))
958 r22 = _mm_mul_ps(rsq22,rinv22);
959 r22 = _mm_andnot_ps(dummy_mask,r22);
961 /* EWALD ELECTROSTATICS */
963 /* Analytical PME correction */
964 zeta2 = _mm_mul_ps(beta2,rsq22);
965 rinv3 = _mm_mul_ps(rinvsq22,rinv22);
966 pmecorrF = gmx_mm_pmecorrF_ps(zeta2);
967 felec = _mm_macc_ps(pmecorrF,beta3,rinv3);
968 felec = _mm_mul_ps(qq22,felec);
969 pmecorrV = gmx_mm_pmecorrV_ps(zeta2);
970 velec = _mm_nmacc_ps(pmecorrV,beta,_mm_sub_ps(rinv22,sh_ewald));
971 velec = _mm_mul_ps(qq22,velec);
973 cutoff_mask = _mm_cmplt_ps(rsq22,rcutoff2);
975 /* Update potential sum for this i atom from the interaction with this j atom. */
976 velec = _mm_and_ps(velec,cutoff_mask);
977 velec = _mm_andnot_ps(dummy_mask,velec);
978 velecsum = _mm_add_ps(velecsum,velec);
982 fscal = _mm_and_ps(fscal,cutoff_mask);
984 fscal = _mm_andnot_ps(dummy_mask,fscal);
986 /* Update vectorial force */
987 fix2 = _mm_macc_ps(dx22,fscal,fix2);
988 fiy2 = _mm_macc_ps(dy22,fscal,fiy2);
989 fiz2 = _mm_macc_ps(dz22,fscal,fiz2);
991 fjx2 = _mm_macc_ps(dx22,fscal,fjx2);
992 fjy2 = _mm_macc_ps(dy22,fscal,fjy2);
993 fjz2 = _mm_macc_ps(dz22,fscal,fjz2);
997 /**************************
998 * CALCULATE INTERACTIONS *
999 **************************/
1001 if (gmx_mm_any_lt(rsq23,rcutoff2))
1004 r23 = _mm_mul_ps(rsq23,rinv23);
1005 r23 = _mm_andnot_ps(dummy_mask,r23);
1007 /* EWALD ELECTROSTATICS */
1009 /* Analytical PME correction */
1010 zeta2 = _mm_mul_ps(beta2,rsq23);
1011 rinv3 = _mm_mul_ps(rinvsq23,rinv23);
1012 pmecorrF = gmx_mm_pmecorrF_ps(zeta2);
1013 felec = _mm_macc_ps(pmecorrF,beta3,rinv3);
1014 felec = _mm_mul_ps(qq23,felec);
1015 pmecorrV = gmx_mm_pmecorrV_ps(zeta2);
1016 velec = _mm_nmacc_ps(pmecorrV,beta,_mm_sub_ps(rinv23,sh_ewald));
1017 velec = _mm_mul_ps(qq23,velec);
1019 cutoff_mask = _mm_cmplt_ps(rsq23,rcutoff2);
1021 /* Update potential sum for this i atom from the interaction with this j atom. */
1022 velec = _mm_and_ps(velec,cutoff_mask);
1023 velec = _mm_andnot_ps(dummy_mask,velec);
1024 velecsum = _mm_add_ps(velecsum,velec);
1028 fscal = _mm_and_ps(fscal,cutoff_mask);
1030 fscal = _mm_andnot_ps(dummy_mask,fscal);
1032 /* Update vectorial force */
1033 fix2 = _mm_macc_ps(dx23,fscal,fix2);
1034 fiy2 = _mm_macc_ps(dy23,fscal,fiy2);
1035 fiz2 = _mm_macc_ps(dz23,fscal,fiz2);
1037 fjx3 = _mm_macc_ps(dx23,fscal,fjx3);
1038 fjy3 = _mm_macc_ps(dy23,fscal,fjy3);
1039 fjz3 = _mm_macc_ps(dz23,fscal,fjz3);
1043 /**************************
1044 * CALCULATE INTERACTIONS *
1045 **************************/
1047 if (gmx_mm_any_lt(rsq31,rcutoff2))
1050 r31 = _mm_mul_ps(rsq31,rinv31);
1051 r31 = _mm_andnot_ps(dummy_mask,r31);
1053 /* EWALD ELECTROSTATICS */
1055 /* Analytical PME correction */
1056 zeta2 = _mm_mul_ps(beta2,rsq31);
1057 rinv3 = _mm_mul_ps(rinvsq31,rinv31);
1058 pmecorrF = gmx_mm_pmecorrF_ps(zeta2);
1059 felec = _mm_macc_ps(pmecorrF,beta3,rinv3);
1060 felec = _mm_mul_ps(qq31,felec);
1061 pmecorrV = gmx_mm_pmecorrV_ps(zeta2);
1062 velec = _mm_nmacc_ps(pmecorrV,beta,_mm_sub_ps(rinv31,sh_ewald));
1063 velec = _mm_mul_ps(qq31,velec);
1065 cutoff_mask = _mm_cmplt_ps(rsq31,rcutoff2);
1067 /* Update potential sum for this i atom from the interaction with this j atom. */
1068 velec = _mm_and_ps(velec,cutoff_mask);
1069 velec = _mm_andnot_ps(dummy_mask,velec);
1070 velecsum = _mm_add_ps(velecsum,velec);
1074 fscal = _mm_and_ps(fscal,cutoff_mask);
1076 fscal = _mm_andnot_ps(dummy_mask,fscal);
1078 /* Update vectorial force */
1079 fix3 = _mm_macc_ps(dx31,fscal,fix3);
1080 fiy3 = _mm_macc_ps(dy31,fscal,fiy3);
1081 fiz3 = _mm_macc_ps(dz31,fscal,fiz3);
1083 fjx1 = _mm_macc_ps(dx31,fscal,fjx1);
1084 fjy1 = _mm_macc_ps(dy31,fscal,fjy1);
1085 fjz1 = _mm_macc_ps(dz31,fscal,fjz1);
1089 /**************************
1090 * CALCULATE INTERACTIONS *
1091 **************************/
1093 if (gmx_mm_any_lt(rsq32,rcutoff2))
1096 r32 = _mm_mul_ps(rsq32,rinv32);
1097 r32 = _mm_andnot_ps(dummy_mask,r32);
1099 /* EWALD ELECTROSTATICS */
1101 /* Analytical PME correction */
1102 zeta2 = _mm_mul_ps(beta2,rsq32);
1103 rinv3 = _mm_mul_ps(rinvsq32,rinv32);
1104 pmecorrF = gmx_mm_pmecorrF_ps(zeta2);
1105 felec = _mm_macc_ps(pmecorrF,beta3,rinv3);
1106 felec = _mm_mul_ps(qq32,felec);
1107 pmecorrV = gmx_mm_pmecorrV_ps(zeta2);
1108 velec = _mm_nmacc_ps(pmecorrV,beta,_mm_sub_ps(rinv32,sh_ewald));
1109 velec = _mm_mul_ps(qq32,velec);
1111 cutoff_mask = _mm_cmplt_ps(rsq32,rcutoff2);
1113 /* Update potential sum for this i atom from the interaction with this j atom. */
1114 velec = _mm_and_ps(velec,cutoff_mask);
1115 velec = _mm_andnot_ps(dummy_mask,velec);
1116 velecsum = _mm_add_ps(velecsum,velec);
1120 fscal = _mm_and_ps(fscal,cutoff_mask);
1122 fscal = _mm_andnot_ps(dummy_mask,fscal);
1124 /* Update vectorial force */
1125 fix3 = _mm_macc_ps(dx32,fscal,fix3);
1126 fiy3 = _mm_macc_ps(dy32,fscal,fiy3);
1127 fiz3 = _mm_macc_ps(dz32,fscal,fiz3);
1129 fjx2 = _mm_macc_ps(dx32,fscal,fjx2);
1130 fjy2 = _mm_macc_ps(dy32,fscal,fjy2);
1131 fjz2 = _mm_macc_ps(dz32,fscal,fjz2);
1135 /**************************
1136 * CALCULATE INTERACTIONS *
1137 **************************/
1139 if (gmx_mm_any_lt(rsq33,rcutoff2))
1142 r33 = _mm_mul_ps(rsq33,rinv33);
1143 r33 = _mm_andnot_ps(dummy_mask,r33);
1145 /* EWALD ELECTROSTATICS */
1147 /* Analytical PME correction */
1148 zeta2 = _mm_mul_ps(beta2,rsq33);
1149 rinv3 = _mm_mul_ps(rinvsq33,rinv33);
1150 pmecorrF = gmx_mm_pmecorrF_ps(zeta2);
1151 felec = _mm_macc_ps(pmecorrF,beta3,rinv3);
1152 felec = _mm_mul_ps(qq33,felec);
1153 pmecorrV = gmx_mm_pmecorrV_ps(zeta2);
1154 velec = _mm_nmacc_ps(pmecorrV,beta,_mm_sub_ps(rinv33,sh_ewald));
1155 velec = _mm_mul_ps(qq33,velec);
1157 cutoff_mask = _mm_cmplt_ps(rsq33,rcutoff2);
1159 /* Update potential sum for this i atom from the interaction with this j atom. */
1160 velec = _mm_and_ps(velec,cutoff_mask);
1161 velec = _mm_andnot_ps(dummy_mask,velec);
1162 velecsum = _mm_add_ps(velecsum,velec);
1166 fscal = _mm_and_ps(fscal,cutoff_mask);
1168 fscal = _mm_andnot_ps(dummy_mask,fscal);
1170 /* Update vectorial force */
1171 fix3 = _mm_macc_ps(dx33,fscal,fix3);
1172 fiy3 = _mm_macc_ps(dy33,fscal,fiy3);
1173 fiz3 = _mm_macc_ps(dz33,fscal,fiz3);
1175 fjx3 = _mm_macc_ps(dx33,fscal,fjx3);
1176 fjy3 = _mm_macc_ps(dy33,fscal,fjy3);
1177 fjz3 = _mm_macc_ps(dz33,fscal,fjz3);
1181 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
1182 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
1183 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
1184 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
1186 gmx_mm_decrement_3rvec_4ptr_swizzle_ps(fjptrA+DIM,fjptrB+DIM,fjptrC+DIM,fjptrD+DIM,
1187 fjx1,fjy1,fjz1,fjx2,fjy2,fjz2,fjx3,fjy3,fjz3);
1189 /* Inner loop uses 306 flops */
1192 /* End of innermost loop */
1194 gmx_mm_update_iforce_3atom_swizzle_ps(fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
1195 f+i_coord_offset+DIM,fshift+i_shift_offset);
1198 /* Update potential energies */
1199 gmx_mm_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
1201 /* Increment number of inner iterations */
1202 inneriter += j_index_end - j_index_start;
1204 /* Outer loop uses 19 flops */
1207 /* Increment number of outer iterations */
1210 /* Update outer/inner flops */
1212 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_W4W4_VF,outeriter*19 + inneriter*306);
1215 * Gromacs nonbonded kernel: nb_kernel_ElecEwSh_VdwNone_GeomW4W4_F_avx_128_fma_single
1216 * Electrostatics interaction: Ewald
1217 * VdW interaction: None
1218 * Geometry: Water4-Water4
1219 * Calculate force/pot: Force
1222 nb_kernel_ElecEwSh_VdwNone_GeomW4W4_F_avx_128_fma_single
1223 (t_nblist * gmx_restrict nlist,
1224 rvec * gmx_restrict xx,
1225 rvec * gmx_restrict ff,
1226 t_forcerec * gmx_restrict fr,
1227 t_mdatoms * gmx_restrict mdatoms,
1228 nb_kernel_data_t * gmx_restrict kernel_data,
1229 t_nrnb * gmx_restrict nrnb)
1231 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
1232 * just 0 for non-waters.
1233 * Suffixes A,B,C,D refer to j loop unrolling done with AVX_128, e.g. for the four different
1234 * jnr indices corresponding to data put in the four positions in the SIMD register.
1236 int i_shift_offset,i_coord_offset,outeriter,inneriter;
1237 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
1238 int jnrA,jnrB,jnrC,jnrD;
1239 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
1240 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
1241 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
1242 real rcutoff_scalar;
1243 real *shiftvec,*fshift,*x,*f;
1244 real *fjptrA,*fjptrB,*fjptrC,*fjptrD;
1245 real scratch[4*DIM];
1246 __m128 fscal,rcutoff,rcutoff2,jidxall;
1248 __m128 ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
1250 __m128 ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
1252 __m128 ix3,iy3,iz3,fix3,fiy3,fiz3,iq3,isai3;
1253 int vdwjidx1A,vdwjidx1B,vdwjidx1C,vdwjidx1D;
1254 __m128 jx1,jy1,jz1,fjx1,fjy1,fjz1,jq1,isaj1;
1255 int vdwjidx2A,vdwjidx2B,vdwjidx2C,vdwjidx2D;
1256 __m128 jx2,jy2,jz2,fjx2,fjy2,fjz2,jq2,isaj2;
1257 int vdwjidx3A,vdwjidx3B,vdwjidx3C,vdwjidx3D;
1258 __m128 jx3,jy3,jz3,fjx3,fjy3,fjz3,jq3,isaj3;
1259 __m128 dx11,dy11,dz11,rsq11,rinv11,rinvsq11,r11,qq11,c6_11,c12_11;
1260 __m128 dx12,dy12,dz12,rsq12,rinv12,rinvsq12,r12,qq12,c6_12,c12_12;
1261 __m128 dx13,dy13,dz13,rsq13,rinv13,rinvsq13,r13,qq13,c6_13,c12_13;
1262 __m128 dx21,dy21,dz21,rsq21,rinv21,rinvsq21,r21,qq21,c6_21,c12_21;
1263 __m128 dx22,dy22,dz22,rsq22,rinv22,rinvsq22,r22,qq22,c6_22,c12_22;
1264 __m128 dx23,dy23,dz23,rsq23,rinv23,rinvsq23,r23,qq23,c6_23,c12_23;
1265 __m128 dx31,dy31,dz31,rsq31,rinv31,rinvsq31,r31,qq31,c6_31,c12_31;
1266 __m128 dx32,dy32,dz32,rsq32,rinv32,rinvsq32,r32,qq32,c6_32,c12_32;
1267 __m128 dx33,dy33,dz33,rsq33,rinv33,rinvsq33,r33,qq33,c6_33,c12_33;
1268 __m128 velec,felec,velecsum,facel,crf,krf,krf2;
1271 __m128 ewtabscale,eweps,twoeweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
1272 __m128 beta,beta2,beta3,zeta2,pmecorrF,pmecorrV,rinv3;
1274 __m128 dummy_mask,cutoff_mask;
1275 __m128 signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
1276 __m128 one = _mm_set1_ps(1.0);
1277 __m128 two = _mm_set1_ps(2.0);
1283 jindex = nlist->jindex;
1285 shiftidx = nlist->shift;
1287 shiftvec = fr->shift_vec[0];
1288 fshift = fr->fshift[0];
1289 facel = _mm_set1_ps(fr->epsfac);
1290 charge = mdatoms->chargeA;
1292 sh_ewald = _mm_set1_ps(fr->ic->sh_ewald);
1293 beta = _mm_set1_ps(fr->ic->ewaldcoeff);
1294 beta2 = _mm_mul_ps(beta,beta);
1295 beta3 = _mm_mul_ps(beta,beta2);
1296 ewtab = fr->ic->tabq_coul_F;
1297 ewtabscale = _mm_set1_ps(fr->ic->tabq_scale);
1298 ewtabhalfspace = _mm_set1_ps(0.5/fr->ic->tabq_scale);
1300 /* Setup water-specific parameters */
1301 inr = nlist->iinr[0];
1302 iq1 = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+1]));
1303 iq2 = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+2]));
1304 iq3 = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+3]));
1306 jq1 = _mm_set1_ps(charge[inr+1]);
1307 jq2 = _mm_set1_ps(charge[inr+2]);
1308 jq3 = _mm_set1_ps(charge[inr+3]);
1309 qq11 = _mm_mul_ps(iq1,jq1);
1310 qq12 = _mm_mul_ps(iq1,jq2);
1311 qq13 = _mm_mul_ps(iq1,jq3);
1312 qq21 = _mm_mul_ps(iq2,jq1);
1313 qq22 = _mm_mul_ps(iq2,jq2);
1314 qq23 = _mm_mul_ps(iq2,jq3);
1315 qq31 = _mm_mul_ps(iq3,jq1);
1316 qq32 = _mm_mul_ps(iq3,jq2);
1317 qq33 = _mm_mul_ps(iq3,jq3);
1319 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
1320 rcutoff_scalar = fr->rcoulomb;
1321 rcutoff = _mm_set1_ps(rcutoff_scalar);
1322 rcutoff2 = _mm_mul_ps(rcutoff,rcutoff);
1324 /* Avoid stupid compiler warnings */
1325 jnrA = jnrB = jnrC = jnrD = 0;
1326 j_coord_offsetA = 0;
1327 j_coord_offsetB = 0;
1328 j_coord_offsetC = 0;
1329 j_coord_offsetD = 0;
1334 for(iidx=0;iidx<4*DIM;iidx++)
1336 scratch[iidx] = 0.0;
1339 /* Start outer loop over neighborlists */
1340 for(iidx=0; iidx<nri; iidx++)
1342 /* Load shift vector for this list */
1343 i_shift_offset = DIM*shiftidx[iidx];
1345 /* Load limits for loop over neighbors */
1346 j_index_start = jindex[iidx];
1347 j_index_end = jindex[iidx+1];
1349 /* Get outer coordinate index */
1351 i_coord_offset = DIM*inr;
1353 /* Load i particle coords and add shift vector */
1354 gmx_mm_load_shift_and_3rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset+DIM,
1355 &ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
1357 fix1 = _mm_setzero_ps();
1358 fiy1 = _mm_setzero_ps();
1359 fiz1 = _mm_setzero_ps();
1360 fix2 = _mm_setzero_ps();
1361 fiy2 = _mm_setzero_ps();
1362 fiz2 = _mm_setzero_ps();
1363 fix3 = _mm_setzero_ps();
1364 fiy3 = _mm_setzero_ps();
1365 fiz3 = _mm_setzero_ps();
1367 /* Start inner kernel loop */
1368 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
1371 /* Get j neighbor index, and coordinate index */
1373 jnrB = jjnr[jidx+1];
1374 jnrC = jjnr[jidx+2];
1375 jnrD = jjnr[jidx+3];
1376 j_coord_offsetA = DIM*jnrA;
1377 j_coord_offsetB = DIM*jnrB;
1378 j_coord_offsetC = DIM*jnrC;
1379 j_coord_offsetD = DIM*jnrD;
1381 /* load j atom coordinates */
1382 gmx_mm_load_3rvec_4ptr_swizzle_ps(x+j_coord_offsetA+DIM,x+j_coord_offsetB+DIM,
1383 x+j_coord_offsetC+DIM,x+j_coord_offsetD+DIM,
1384 &jx1,&jy1,&jz1,&jx2,&jy2,&jz2,&jx3,&jy3,&jz3);
1386 /* Calculate displacement vector */
1387 dx11 = _mm_sub_ps(ix1,jx1);
1388 dy11 = _mm_sub_ps(iy1,jy1);
1389 dz11 = _mm_sub_ps(iz1,jz1);
1390 dx12 = _mm_sub_ps(ix1,jx2);
1391 dy12 = _mm_sub_ps(iy1,jy2);
1392 dz12 = _mm_sub_ps(iz1,jz2);
1393 dx13 = _mm_sub_ps(ix1,jx3);
1394 dy13 = _mm_sub_ps(iy1,jy3);
1395 dz13 = _mm_sub_ps(iz1,jz3);
1396 dx21 = _mm_sub_ps(ix2,jx1);
1397 dy21 = _mm_sub_ps(iy2,jy1);
1398 dz21 = _mm_sub_ps(iz2,jz1);
1399 dx22 = _mm_sub_ps(ix2,jx2);
1400 dy22 = _mm_sub_ps(iy2,jy2);
1401 dz22 = _mm_sub_ps(iz2,jz2);
1402 dx23 = _mm_sub_ps(ix2,jx3);
1403 dy23 = _mm_sub_ps(iy2,jy3);
1404 dz23 = _mm_sub_ps(iz2,jz3);
1405 dx31 = _mm_sub_ps(ix3,jx1);
1406 dy31 = _mm_sub_ps(iy3,jy1);
1407 dz31 = _mm_sub_ps(iz3,jz1);
1408 dx32 = _mm_sub_ps(ix3,jx2);
1409 dy32 = _mm_sub_ps(iy3,jy2);
1410 dz32 = _mm_sub_ps(iz3,jz2);
1411 dx33 = _mm_sub_ps(ix3,jx3);
1412 dy33 = _mm_sub_ps(iy3,jy3);
1413 dz33 = _mm_sub_ps(iz3,jz3);
1415 /* Calculate squared distance and things based on it */
1416 rsq11 = gmx_mm_calc_rsq_ps(dx11,dy11,dz11);
1417 rsq12 = gmx_mm_calc_rsq_ps(dx12,dy12,dz12);
1418 rsq13 = gmx_mm_calc_rsq_ps(dx13,dy13,dz13);
1419 rsq21 = gmx_mm_calc_rsq_ps(dx21,dy21,dz21);
1420 rsq22 = gmx_mm_calc_rsq_ps(dx22,dy22,dz22);
1421 rsq23 = gmx_mm_calc_rsq_ps(dx23,dy23,dz23);
1422 rsq31 = gmx_mm_calc_rsq_ps(dx31,dy31,dz31);
1423 rsq32 = gmx_mm_calc_rsq_ps(dx32,dy32,dz32);
1424 rsq33 = gmx_mm_calc_rsq_ps(dx33,dy33,dz33);
1426 rinv11 = gmx_mm_invsqrt_ps(rsq11);
1427 rinv12 = gmx_mm_invsqrt_ps(rsq12);
1428 rinv13 = gmx_mm_invsqrt_ps(rsq13);
1429 rinv21 = gmx_mm_invsqrt_ps(rsq21);
1430 rinv22 = gmx_mm_invsqrt_ps(rsq22);
1431 rinv23 = gmx_mm_invsqrt_ps(rsq23);
1432 rinv31 = gmx_mm_invsqrt_ps(rsq31);
1433 rinv32 = gmx_mm_invsqrt_ps(rsq32);
1434 rinv33 = gmx_mm_invsqrt_ps(rsq33);
1436 rinvsq11 = _mm_mul_ps(rinv11,rinv11);
1437 rinvsq12 = _mm_mul_ps(rinv12,rinv12);
1438 rinvsq13 = _mm_mul_ps(rinv13,rinv13);
1439 rinvsq21 = _mm_mul_ps(rinv21,rinv21);
1440 rinvsq22 = _mm_mul_ps(rinv22,rinv22);
1441 rinvsq23 = _mm_mul_ps(rinv23,rinv23);
1442 rinvsq31 = _mm_mul_ps(rinv31,rinv31);
1443 rinvsq32 = _mm_mul_ps(rinv32,rinv32);
1444 rinvsq33 = _mm_mul_ps(rinv33,rinv33);
1446 fjx1 = _mm_setzero_ps();
1447 fjy1 = _mm_setzero_ps();
1448 fjz1 = _mm_setzero_ps();
1449 fjx2 = _mm_setzero_ps();
1450 fjy2 = _mm_setzero_ps();
1451 fjz2 = _mm_setzero_ps();
1452 fjx3 = _mm_setzero_ps();
1453 fjy3 = _mm_setzero_ps();
1454 fjz3 = _mm_setzero_ps();
1456 /**************************
1457 * CALCULATE INTERACTIONS *
1458 **************************/
1460 if (gmx_mm_any_lt(rsq11,rcutoff2))
1463 r11 = _mm_mul_ps(rsq11,rinv11);
1465 /* EWALD ELECTROSTATICS */
1467 /* Analytical PME correction */
1468 zeta2 = _mm_mul_ps(beta2,rsq11);
1469 rinv3 = _mm_mul_ps(rinvsq11,rinv11);
1470 pmecorrF = gmx_mm_pmecorrF_ps(zeta2);
1471 felec = _mm_macc_ps(pmecorrF,beta3,rinv3);
1472 felec = _mm_mul_ps(qq11,felec);
1474 cutoff_mask = _mm_cmplt_ps(rsq11,rcutoff2);
1478 fscal = _mm_and_ps(fscal,cutoff_mask);
1480 /* Update vectorial force */
1481 fix1 = _mm_macc_ps(dx11,fscal,fix1);
1482 fiy1 = _mm_macc_ps(dy11,fscal,fiy1);
1483 fiz1 = _mm_macc_ps(dz11,fscal,fiz1);
1485 fjx1 = _mm_macc_ps(dx11,fscal,fjx1);
1486 fjy1 = _mm_macc_ps(dy11,fscal,fjy1);
1487 fjz1 = _mm_macc_ps(dz11,fscal,fjz1);
1491 /**************************
1492 * CALCULATE INTERACTIONS *
1493 **************************/
1495 if (gmx_mm_any_lt(rsq12,rcutoff2))
1498 r12 = _mm_mul_ps(rsq12,rinv12);
1500 /* EWALD ELECTROSTATICS */
1502 /* Analytical PME correction */
1503 zeta2 = _mm_mul_ps(beta2,rsq12);
1504 rinv3 = _mm_mul_ps(rinvsq12,rinv12);
1505 pmecorrF = gmx_mm_pmecorrF_ps(zeta2);
1506 felec = _mm_macc_ps(pmecorrF,beta3,rinv3);
1507 felec = _mm_mul_ps(qq12,felec);
1509 cutoff_mask = _mm_cmplt_ps(rsq12,rcutoff2);
1513 fscal = _mm_and_ps(fscal,cutoff_mask);
1515 /* Update vectorial force */
1516 fix1 = _mm_macc_ps(dx12,fscal,fix1);
1517 fiy1 = _mm_macc_ps(dy12,fscal,fiy1);
1518 fiz1 = _mm_macc_ps(dz12,fscal,fiz1);
1520 fjx2 = _mm_macc_ps(dx12,fscal,fjx2);
1521 fjy2 = _mm_macc_ps(dy12,fscal,fjy2);
1522 fjz2 = _mm_macc_ps(dz12,fscal,fjz2);
1526 /**************************
1527 * CALCULATE INTERACTIONS *
1528 **************************/
1530 if (gmx_mm_any_lt(rsq13,rcutoff2))
1533 r13 = _mm_mul_ps(rsq13,rinv13);
1535 /* EWALD ELECTROSTATICS */
1537 /* Analytical PME correction */
1538 zeta2 = _mm_mul_ps(beta2,rsq13);
1539 rinv3 = _mm_mul_ps(rinvsq13,rinv13);
1540 pmecorrF = gmx_mm_pmecorrF_ps(zeta2);
1541 felec = _mm_macc_ps(pmecorrF,beta3,rinv3);
1542 felec = _mm_mul_ps(qq13,felec);
1544 cutoff_mask = _mm_cmplt_ps(rsq13,rcutoff2);
1548 fscal = _mm_and_ps(fscal,cutoff_mask);
1550 /* Update vectorial force */
1551 fix1 = _mm_macc_ps(dx13,fscal,fix1);
1552 fiy1 = _mm_macc_ps(dy13,fscal,fiy1);
1553 fiz1 = _mm_macc_ps(dz13,fscal,fiz1);
1555 fjx3 = _mm_macc_ps(dx13,fscal,fjx3);
1556 fjy3 = _mm_macc_ps(dy13,fscal,fjy3);
1557 fjz3 = _mm_macc_ps(dz13,fscal,fjz3);
1561 /**************************
1562 * CALCULATE INTERACTIONS *
1563 **************************/
1565 if (gmx_mm_any_lt(rsq21,rcutoff2))
1568 r21 = _mm_mul_ps(rsq21,rinv21);
1570 /* EWALD ELECTROSTATICS */
1572 /* Analytical PME correction */
1573 zeta2 = _mm_mul_ps(beta2,rsq21);
1574 rinv3 = _mm_mul_ps(rinvsq21,rinv21);
1575 pmecorrF = gmx_mm_pmecorrF_ps(zeta2);
1576 felec = _mm_macc_ps(pmecorrF,beta3,rinv3);
1577 felec = _mm_mul_ps(qq21,felec);
1579 cutoff_mask = _mm_cmplt_ps(rsq21,rcutoff2);
1583 fscal = _mm_and_ps(fscal,cutoff_mask);
1585 /* Update vectorial force */
1586 fix2 = _mm_macc_ps(dx21,fscal,fix2);
1587 fiy2 = _mm_macc_ps(dy21,fscal,fiy2);
1588 fiz2 = _mm_macc_ps(dz21,fscal,fiz2);
1590 fjx1 = _mm_macc_ps(dx21,fscal,fjx1);
1591 fjy1 = _mm_macc_ps(dy21,fscal,fjy1);
1592 fjz1 = _mm_macc_ps(dz21,fscal,fjz1);
1596 /**************************
1597 * CALCULATE INTERACTIONS *
1598 **************************/
1600 if (gmx_mm_any_lt(rsq22,rcutoff2))
1603 r22 = _mm_mul_ps(rsq22,rinv22);
1605 /* EWALD ELECTROSTATICS */
1607 /* Analytical PME correction */
1608 zeta2 = _mm_mul_ps(beta2,rsq22);
1609 rinv3 = _mm_mul_ps(rinvsq22,rinv22);
1610 pmecorrF = gmx_mm_pmecorrF_ps(zeta2);
1611 felec = _mm_macc_ps(pmecorrF,beta3,rinv3);
1612 felec = _mm_mul_ps(qq22,felec);
1614 cutoff_mask = _mm_cmplt_ps(rsq22,rcutoff2);
1618 fscal = _mm_and_ps(fscal,cutoff_mask);
1620 /* Update vectorial force */
1621 fix2 = _mm_macc_ps(dx22,fscal,fix2);
1622 fiy2 = _mm_macc_ps(dy22,fscal,fiy2);
1623 fiz2 = _mm_macc_ps(dz22,fscal,fiz2);
1625 fjx2 = _mm_macc_ps(dx22,fscal,fjx2);
1626 fjy2 = _mm_macc_ps(dy22,fscal,fjy2);
1627 fjz2 = _mm_macc_ps(dz22,fscal,fjz2);
1631 /**************************
1632 * CALCULATE INTERACTIONS *
1633 **************************/
1635 if (gmx_mm_any_lt(rsq23,rcutoff2))
1638 r23 = _mm_mul_ps(rsq23,rinv23);
1640 /* EWALD ELECTROSTATICS */
1642 /* Analytical PME correction */
1643 zeta2 = _mm_mul_ps(beta2,rsq23);
1644 rinv3 = _mm_mul_ps(rinvsq23,rinv23);
1645 pmecorrF = gmx_mm_pmecorrF_ps(zeta2);
1646 felec = _mm_macc_ps(pmecorrF,beta3,rinv3);
1647 felec = _mm_mul_ps(qq23,felec);
1649 cutoff_mask = _mm_cmplt_ps(rsq23,rcutoff2);
1653 fscal = _mm_and_ps(fscal,cutoff_mask);
1655 /* Update vectorial force */
1656 fix2 = _mm_macc_ps(dx23,fscal,fix2);
1657 fiy2 = _mm_macc_ps(dy23,fscal,fiy2);
1658 fiz2 = _mm_macc_ps(dz23,fscal,fiz2);
1660 fjx3 = _mm_macc_ps(dx23,fscal,fjx3);
1661 fjy3 = _mm_macc_ps(dy23,fscal,fjy3);
1662 fjz3 = _mm_macc_ps(dz23,fscal,fjz3);
1666 /**************************
1667 * CALCULATE INTERACTIONS *
1668 **************************/
1670 if (gmx_mm_any_lt(rsq31,rcutoff2))
1673 r31 = _mm_mul_ps(rsq31,rinv31);
1675 /* EWALD ELECTROSTATICS */
1677 /* Analytical PME correction */
1678 zeta2 = _mm_mul_ps(beta2,rsq31);
1679 rinv3 = _mm_mul_ps(rinvsq31,rinv31);
1680 pmecorrF = gmx_mm_pmecorrF_ps(zeta2);
1681 felec = _mm_macc_ps(pmecorrF,beta3,rinv3);
1682 felec = _mm_mul_ps(qq31,felec);
1684 cutoff_mask = _mm_cmplt_ps(rsq31,rcutoff2);
1688 fscal = _mm_and_ps(fscal,cutoff_mask);
1690 /* Update vectorial force */
1691 fix3 = _mm_macc_ps(dx31,fscal,fix3);
1692 fiy3 = _mm_macc_ps(dy31,fscal,fiy3);
1693 fiz3 = _mm_macc_ps(dz31,fscal,fiz3);
1695 fjx1 = _mm_macc_ps(dx31,fscal,fjx1);
1696 fjy1 = _mm_macc_ps(dy31,fscal,fjy1);
1697 fjz1 = _mm_macc_ps(dz31,fscal,fjz1);
1701 /**************************
1702 * CALCULATE INTERACTIONS *
1703 **************************/
1705 if (gmx_mm_any_lt(rsq32,rcutoff2))
1708 r32 = _mm_mul_ps(rsq32,rinv32);
1710 /* EWALD ELECTROSTATICS */
1712 /* Analytical PME correction */
1713 zeta2 = _mm_mul_ps(beta2,rsq32);
1714 rinv3 = _mm_mul_ps(rinvsq32,rinv32);
1715 pmecorrF = gmx_mm_pmecorrF_ps(zeta2);
1716 felec = _mm_macc_ps(pmecorrF,beta3,rinv3);
1717 felec = _mm_mul_ps(qq32,felec);
1719 cutoff_mask = _mm_cmplt_ps(rsq32,rcutoff2);
1723 fscal = _mm_and_ps(fscal,cutoff_mask);
1725 /* Update vectorial force */
1726 fix3 = _mm_macc_ps(dx32,fscal,fix3);
1727 fiy3 = _mm_macc_ps(dy32,fscal,fiy3);
1728 fiz3 = _mm_macc_ps(dz32,fscal,fiz3);
1730 fjx2 = _mm_macc_ps(dx32,fscal,fjx2);
1731 fjy2 = _mm_macc_ps(dy32,fscal,fjy2);
1732 fjz2 = _mm_macc_ps(dz32,fscal,fjz2);
1736 /**************************
1737 * CALCULATE INTERACTIONS *
1738 **************************/
1740 if (gmx_mm_any_lt(rsq33,rcutoff2))
1743 r33 = _mm_mul_ps(rsq33,rinv33);
1745 /* EWALD ELECTROSTATICS */
1747 /* Analytical PME correction */
1748 zeta2 = _mm_mul_ps(beta2,rsq33);
1749 rinv3 = _mm_mul_ps(rinvsq33,rinv33);
1750 pmecorrF = gmx_mm_pmecorrF_ps(zeta2);
1751 felec = _mm_macc_ps(pmecorrF,beta3,rinv3);
1752 felec = _mm_mul_ps(qq33,felec);
1754 cutoff_mask = _mm_cmplt_ps(rsq33,rcutoff2);
1758 fscal = _mm_and_ps(fscal,cutoff_mask);
1760 /* Update vectorial force */
1761 fix3 = _mm_macc_ps(dx33,fscal,fix3);
1762 fiy3 = _mm_macc_ps(dy33,fscal,fiy3);
1763 fiz3 = _mm_macc_ps(dz33,fscal,fiz3);
1765 fjx3 = _mm_macc_ps(dx33,fscal,fjx3);
1766 fjy3 = _mm_macc_ps(dy33,fscal,fjy3);
1767 fjz3 = _mm_macc_ps(dz33,fscal,fjz3);
1771 fjptrA = f+j_coord_offsetA;
1772 fjptrB = f+j_coord_offsetB;
1773 fjptrC = f+j_coord_offsetC;
1774 fjptrD = f+j_coord_offsetD;
1776 gmx_mm_decrement_3rvec_4ptr_swizzle_ps(fjptrA+DIM,fjptrB+DIM,fjptrC+DIM,fjptrD+DIM,
1777 fjx1,fjy1,fjz1,fjx2,fjy2,fjz2,fjx3,fjy3,fjz3);
1779 /* Inner loop uses 279 flops */
1782 if(jidx<j_index_end)
1785 /* Get j neighbor index, and coordinate index */
1786 jnrlistA = jjnr[jidx];
1787 jnrlistB = jjnr[jidx+1];
1788 jnrlistC = jjnr[jidx+2];
1789 jnrlistD = jjnr[jidx+3];
1790 /* Sign of each element will be negative for non-real atoms.
1791 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
1792 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
1794 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
1795 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
1796 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
1797 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
1798 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
1799 j_coord_offsetA = DIM*jnrA;
1800 j_coord_offsetB = DIM*jnrB;
1801 j_coord_offsetC = DIM*jnrC;
1802 j_coord_offsetD = DIM*jnrD;
1804 /* load j atom coordinates */
1805 gmx_mm_load_3rvec_4ptr_swizzle_ps(x+j_coord_offsetA+DIM,x+j_coord_offsetB+DIM,
1806 x+j_coord_offsetC+DIM,x+j_coord_offsetD+DIM,
1807 &jx1,&jy1,&jz1,&jx2,&jy2,&jz2,&jx3,&jy3,&jz3);
1809 /* Calculate displacement vector */
1810 dx11 = _mm_sub_ps(ix1,jx1);
1811 dy11 = _mm_sub_ps(iy1,jy1);
1812 dz11 = _mm_sub_ps(iz1,jz1);
1813 dx12 = _mm_sub_ps(ix1,jx2);
1814 dy12 = _mm_sub_ps(iy1,jy2);
1815 dz12 = _mm_sub_ps(iz1,jz2);
1816 dx13 = _mm_sub_ps(ix1,jx3);
1817 dy13 = _mm_sub_ps(iy1,jy3);
1818 dz13 = _mm_sub_ps(iz1,jz3);
1819 dx21 = _mm_sub_ps(ix2,jx1);
1820 dy21 = _mm_sub_ps(iy2,jy1);
1821 dz21 = _mm_sub_ps(iz2,jz1);
1822 dx22 = _mm_sub_ps(ix2,jx2);
1823 dy22 = _mm_sub_ps(iy2,jy2);
1824 dz22 = _mm_sub_ps(iz2,jz2);
1825 dx23 = _mm_sub_ps(ix2,jx3);
1826 dy23 = _mm_sub_ps(iy2,jy3);
1827 dz23 = _mm_sub_ps(iz2,jz3);
1828 dx31 = _mm_sub_ps(ix3,jx1);
1829 dy31 = _mm_sub_ps(iy3,jy1);
1830 dz31 = _mm_sub_ps(iz3,jz1);
1831 dx32 = _mm_sub_ps(ix3,jx2);
1832 dy32 = _mm_sub_ps(iy3,jy2);
1833 dz32 = _mm_sub_ps(iz3,jz2);
1834 dx33 = _mm_sub_ps(ix3,jx3);
1835 dy33 = _mm_sub_ps(iy3,jy3);
1836 dz33 = _mm_sub_ps(iz3,jz3);
1838 /* Calculate squared distance and things based on it */
1839 rsq11 = gmx_mm_calc_rsq_ps(dx11,dy11,dz11);
1840 rsq12 = gmx_mm_calc_rsq_ps(dx12,dy12,dz12);
1841 rsq13 = gmx_mm_calc_rsq_ps(dx13,dy13,dz13);
1842 rsq21 = gmx_mm_calc_rsq_ps(dx21,dy21,dz21);
1843 rsq22 = gmx_mm_calc_rsq_ps(dx22,dy22,dz22);
1844 rsq23 = gmx_mm_calc_rsq_ps(dx23,dy23,dz23);
1845 rsq31 = gmx_mm_calc_rsq_ps(dx31,dy31,dz31);
1846 rsq32 = gmx_mm_calc_rsq_ps(dx32,dy32,dz32);
1847 rsq33 = gmx_mm_calc_rsq_ps(dx33,dy33,dz33);
1849 rinv11 = gmx_mm_invsqrt_ps(rsq11);
1850 rinv12 = gmx_mm_invsqrt_ps(rsq12);
1851 rinv13 = gmx_mm_invsqrt_ps(rsq13);
1852 rinv21 = gmx_mm_invsqrt_ps(rsq21);
1853 rinv22 = gmx_mm_invsqrt_ps(rsq22);
1854 rinv23 = gmx_mm_invsqrt_ps(rsq23);
1855 rinv31 = gmx_mm_invsqrt_ps(rsq31);
1856 rinv32 = gmx_mm_invsqrt_ps(rsq32);
1857 rinv33 = gmx_mm_invsqrt_ps(rsq33);
1859 rinvsq11 = _mm_mul_ps(rinv11,rinv11);
1860 rinvsq12 = _mm_mul_ps(rinv12,rinv12);
1861 rinvsq13 = _mm_mul_ps(rinv13,rinv13);
1862 rinvsq21 = _mm_mul_ps(rinv21,rinv21);
1863 rinvsq22 = _mm_mul_ps(rinv22,rinv22);
1864 rinvsq23 = _mm_mul_ps(rinv23,rinv23);
1865 rinvsq31 = _mm_mul_ps(rinv31,rinv31);
1866 rinvsq32 = _mm_mul_ps(rinv32,rinv32);
1867 rinvsq33 = _mm_mul_ps(rinv33,rinv33);
1869 fjx1 = _mm_setzero_ps();
1870 fjy1 = _mm_setzero_ps();
1871 fjz1 = _mm_setzero_ps();
1872 fjx2 = _mm_setzero_ps();
1873 fjy2 = _mm_setzero_ps();
1874 fjz2 = _mm_setzero_ps();
1875 fjx3 = _mm_setzero_ps();
1876 fjy3 = _mm_setzero_ps();
1877 fjz3 = _mm_setzero_ps();
1879 /**************************
1880 * CALCULATE INTERACTIONS *
1881 **************************/
1883 if (gmx_mm_any_lt(rsq11,rcutoff2))
1886 r11 = _mm_mul_ps(rsq11,rinv11);
1887 r11 = _mm_andnot_ps(dummy_mask,r11);
1889 /* EWALD ELECTROSTATICS */
1891 /* Analytical PME correction */
1892 zeta2 = _mm_mul_ps(beta2,rsq11);
1893 rinv3 = _mm_mul_ps(rinvsq11,rinv11);
1894 pmecorrF = gmx_mm_pmecorrF_ps(zeta2);
1895 felec = _mm_macc_ps(pmecorrF,beta3,rinv3);
1896 felec = _mm_mul_ps(qq11,felec);
1898 cutoff_mask = _mm_cmplt_ps(rsq11,rcutoff2);
1902 fscal = _mm_and_ps(fscal,cutoff_mask);
1904 fscal = _mm_andnot_ps(dummy_mask,fscal);
1906 /* Update vectorial force */
1907 fix1 = _mm_macc_ps(dx11,fscal,fix1);
1908 fiy1 = _mm_macc_ps(dy11,fscal,fiy1);
1909 fiz1 = _mm_macc_ps(dz11,fscal,fiz1);
1911 fjx1 = _mm_macc_ps(dx11,fscal,fjx1);
1912 fjy1 = _mm_macc_ps(dy11,fscal,fjy1);
1913 fjz1 = _mm_macc_ps(dz11,fscal,fjz1);
1917 /**************************
1918 * CALCULATE INTERACTIONS *
1919 **************************/
1921 if (gmx_mm_any_lt(rsq12,rcutoff2))
1924 r12 = _mm_mul_ps(rsq12,rinv12);
1925 r12 = _mm_andnot_ps(dummy_mask,r12);
1927 /* EWALD ELECTROSTATICS */
1929 /* Analytical PME correction */
1930 zeta2 = _mm_mul_ps(beta2,rsq12);
1931 rinv3 = _mm_mul_ps(rinvsq12,rinv12);
1932 pmecorrF = gmx_mm_pmecorrF_ps(zeta2);
1933 felec = _mm_macc_ps(pmecorrF,beta3,rinv3);
1934 felec = _mm_mul_ps(qq12,felec);
1936 cutoff_mask = _mm_cmplt_ps(rsq12,rcutoff2);
1940 fscal = _mm_and_ps(fscal,cutoff_mask);
1942 fscal = _mm_andnot_ps(dummy_mask,fscal);
1944 /* Update vectorial force */
1945 fix1 = _mm_macc_ps(dx12,fscal,fix1);
1946 fiy1 = _mm_macc_ps(dy12,fscal,fiy1);
1947 fiz1 = _mm_macc_ps(dz12,fscal,fiz1);
1949 fjx2 = _mm_macc_ps(dx12,fscal,fjx2);
1950 fjy2 = _mm_macc_ps(dy12,fscal,fjy2);
1951 fjz2 = _mm_macc_ps(dz12,fscal,fjz2);
1955 /**************************
1956 * CALCULATE INTERACTIONS *
1957 **************************/
1959 if (gmx_mm_any_lt(rsq13,rcutoff2))
1962 r13 = _mm_mul_ps(rsq13,rinv13);
1963 r13 = _mm_andnot_ps(dummy_mask,r13);
1965 /* EWALD ELECTROSTATICS */
1967 /* Analytical PME correction */
1968 zeta2 = _mm_mul_ps(beta2,rsq13);
1969 rinv3 = _mm_mul_ps(rinvsq13,rinv13);
1970 pmecorrF = gmx_mm_pmecorrF_ps(zeta2);
1971 felec = _mm_macc_ps(pmecorrF,beta3,rinv3);
1972 felec = _mm_mul_ps(qq13,felec);
1974 cutoff_mask = _mm_cmplt_ps(rsq13,rcutoff2);
1978 fscal = _mm_and_ps(fscal,cutoff_mask);
1980 fscal = _mm_andnot_ps(dummy_mask,fscal);
1982 /* Update vectorial force */
1983 fix1 = _mm_macc_ps(dx13,fscal,fix1);
1984 fiy1 = _mm_macc_ps(dy13,fscal,fiy1);
1985 fiz1 = _mm_macc_ps(dz13,fscal,fiz1);
1987 fjx3 = _mm_macc_ps(dx13,fscal,fjx3);
1988 fjy3 = _mm_macc_ps(dy13,fscal,fjy3);
1989 fjz3 = _mm_macc_ps(dz13,fscal,fjz3);
1993 /**************************
1994 * CALCULATE INTERACTIONS *
1995 **************************/
1997 if (gmx_mm_any_lt(rsq21,rcutoff2))
2000 r21 = _mm_mul_ps(rsq21,rinv21);
2001 r21 = _mm_andnot_ps(dummy_mask,r21);
2003 /* EWALD ELECTROSTATICS */
2005 /* Analytical PME correction */
2006 zeta2 = _mm_mul_ps(beta2,rsq21);
2007 rinv3 = _mm_mul_ps(rinvsq21,rinv21);
2008 pmecorrF = gmx_mm_pmecorrF_ps(zeta2);
2009 felec = _mm_macc_ps(pmecorrF,beta3,rinv3);
2010 felec = _mm_mul_ps(qq21,felec);
2012 cutoff_mask = _mm_cmplt_ps(rsq21,rcutoff2);
2016 fscal = _mm_and_ps(fscal,cutoff_mask);
2018 fscal = _mm_andnot_ps(dummy_mask,fscal);
2020 /* Update vectorial force */
2021 fix2 = _mm_macc_ps(dx21,fscal,fix2);
2022 fiy2 = _mm_macc_ps(dy21,fscal,fiy2);
2023 fiz2 = _mm_macc_ps(dz21,fscal,fiz2);
2025 fjx1 = _mm_macc_ps(dx21,fscal,fjx1);
2026 fjy1 = _mm_macc_ps(dy21,fscal,fjy1);
2027 fjz1 = _mm_macc_ps(dz21,fscal,fjz1);
2031 /**************************
2032 * CALCULATE INTERACTIONS *
2033 **************************/
2035 if (gmx_mm_any_lt(rsq22,rcutoff2))
2038 r22 = _mm_mul_ps(rsq22,rinv22);
2039 r22 = _mm_andnot_ps(dummy_mask,r22);
2041 /* EWALD ELECTROSTATICS */
2043 /* Analytical PME correction */
2044 zeta2 = _mm_mul_ps(beta2,rsq22);
2045 rinv3 = _mm_mul_ps(rinvsq22,rinv22);
2046 pmecorrF = gmx_mm_pmecorrF_ps(zeta2);
2047 felec = _mm_macc_ps(pmecorrF,beta3,rinv3);
2048 felec = _mm_mul_ps(qq22,felec);
2050 cutoff_mask = _mm_cmplt_ps(rsq22,rcutoff2);
2054 fscal = _mm_and_ps(fscal,cutoff_mask);
2056 fscal = _mm_andnot_ps(dummy_mask,fscal);
2058 /* Update vectorial force */
2059 fix2 = _mm_macc_ps(dx22,fscal,fix2);
2060 fiy2 = _mm_macc_ps(dy22,fscal,fiy2);
2061 fiz2 = _mm_macc_ps(dz22,fscal,fiz2);
2063 fjx2 = _mm_macc_ps(dx22,fscal,fjx2);
2064 fjy2 = _mm_macc_ps(dy22,fscal,fjy2);
2065 fjz2 = _mm_macc_ps(dz22,fscal,fjz2);
2069 /**************************
2070 * CALCULATE INTERACTIONS *
2071 **************************/
2073 if (gmx_mm_any_lt(rsq23,rcutoff2))
2076 r23 = _mm_mul_ps(rsq23,rinv23);
2077 r23 = _mm_andnot_ps(dummy_mask,r23);
2079 /* EWALD ELECTROSTATICS */
2081 /* Analytical PME correction */
2082 zeta2 = _mm_mul_ps(beta2,rsq23);
2083 rinv3 = _mm_mul_ps(rinvsq23,rinv23);
2084 pmecorrF = gmx_mm_pmecorrF_ps(zeta2);
2085 felec = _mm_macc_ps(pmecorrF,beta3,rinv3);
2086 felec = _mm_mul_ps(qq23,felec);
2088 cutoff_mask = _mm_cmplt_ps(rsq23,rcutoff2);
2092 fscal = _mm_and_ps(fscal,cutoff_mask);
2094 fscal = _mm_andnot_ps(dummy_mask,fscal);
2096 /* Update vectorial force */
2097 fix2 = _mm_macc_ps(dx23,fscal,fix2);
2098 fiy2 = _mm_macc_ps(dy23,fscal,fiy2);
2099 fiz2 = _mm_macc_ps(dz23,fscal,fiz2);
2101 fjx3 = _mm_macc_ps(dx23,fscal,fjx3);
2102 fjy3 = _mm_macc_ps(dy23,fscal,fjy3);
2103 fjz3 = _mm_macc_ps(dz23,fscal,fjz3);
2107 /**************************
2108 * CALCULATE INTERACTIONS *
2109 **************************/
2111 if (gmx_mm_any_lt(rsq31,rcutoff2))
2114 r31 = _mm_mul_ps(rsq31,rinv31);
2115 r31 = _mm_andnot_ps(dummy_mask,r31);
2117 /* EWALD ELECTROSTATICS */
2119 /* Analytical PME correction */
2120 zeta2 = _mm_mul_ps(beta2,rsq31);
2121 rinv3 = _mm_mul_ps(rinvsq31,rinv31);
2122 pmecorrF = gmx_mm_pmecorrF_ps(zeta2);
2123 felec = _mm_macc_ps(pmecorrF,beta3,rinv3);
2124 felec = _mm_mul_ps(qq31,felec);
2126 cutoff_mask = _mm_cmplt_ps(rsq31,rcutoff2);
2130 fscal = _mm_and_ps(fscal,cutoff_mask);
2132 fscal = _mm_andnot_ps(dummy_mask,fscal);
2134 /* Update vectorial force */
2135 fix3 = _mm_macc_ps(dx31,fscal,fix3);
2136 fiy3 = _mm_macc_ps(dy31,fscal,fiy3);
2137 fiz3 = _mm_macc_ps(dz31,fscal,fiz3);
2139 fjx1 = _mm_macc_ps(dx31,fscal,fjx1);
2140 fjy1 = _mm_macc_ps(dy31,fscal,fjy1);
2141 fjz1 = _mm_macc_ps(dz31,fscal,fjz1);
2145 /**************************
2146 * CALCULATE INTERACTIONS *
2147 **************************/
2149 if (gmx_mm_any_lt(rsq32,rcutoff2))
2152 r32 = _mm_mul_ps(rsq32,rinv32);
2153 r32 = _mm_andnot_ps(dummy_mask,r32);
2155 /* EWALD ELECTROSTATICS */
2157 /* Analytical PME correction */
2158 zeta2 = _mm_mul_ps(beta2,rsq32);
2159 rinv3 = _mm_mul_ps(rinvsq32,rinv32);
2160 pmecorrF = gmx_mm_pmecorrF_ps(zeta2);
2161 felec = _mm_macc_ps(pmecorrF,beta3,rinv3);
2162 felec = _mm_mul_ps(qq32,felec);
2164 cutoff_mask = _mm_cmplt_ps(rsq32,rcutoff2);
2168 fscal = _mm_and_ps(fscal,cutoff_mask);
2170 fscal = _mm_andnot_ps(dummy_mask,fscal);
2172 /* Update vectorial force */
2173 fix3 = _mm_macc_ps(dx32,fscal,fix3);
2174 fiy3 = _mm_macc_ps(dy32,fscal,fiy3);
2175 fiz3 = _mm_macc_ps(dz32,fscal,fiz3);
2177 fjx2 = _mm_macc_ps(dx32,fscal,fjx2);
2178 fjy2 = _mm_macc_ps(dy32,fscal,fjy2);
2179 fjz2 = _mm_macc_ps(dz32,fscal,fjz2);
2183 /**************************
2184 * CALCULATE INTERACTIONS *
2185 **************************/
2187 if (gmx_mm_any_lt(rsq33,rcutoff2))
2190 r33 = _mm_mul_ps(rsq33,rinv33);
2191 r33 = _mm_andnot_ps(dummy_mask,r33);
2193 /* EWALD ELECTROSTATICS */
2195 /* Analytical PME correction */
2196 zeta2 = _mm_mul_ps(beta2,rsq33);
2197 rinv3 = _mm_mul_ps(rinvsq33,rinv33);
2198 pmecorrF = gmx_mm_pmecorrF_ps(zeta2);
2199 felec = _mm_macc_ps(pmecorrF,beta3,rinv3);
2200 felec = _mm_mul_ps(qq33,felec);
2202 cutoff_mask = _mm_cmplt_ps(rsq33,rcutoff2);
2206 fscal = _mm_and_ps(fscal,cutoff_mask);
2208 fscal = _mm_andnot_ps(dummy_mask,fscal);
2210 /* Update vectorial force */
2211 fix3 = _mm_macc_ps(dx33,fscal,fix3);
2212 fiy3 = _mm_macc_ps(dy33,fscal,fiy3);
2213 fiz3 = _mm_macc_ps(dz33,fscal,fiz3);
2215 fjx3 = _mm_macc_ps(dx33,fscal,fjx3);
2216 fjy3 = _mm_macc_ps(dy33,fscal,fjy3);
2217 fjz3 = _mm_macc_ps(dz33,fscal,fjz3);
2221 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
2222 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
2223 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
2224 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
2226 gmx_mm_decrement_3rvec_4ptr_swizzle_ps(fjptrA+DIM,fjptrB+DIM,fjptrC+DIM,fjptrD+DIM,
2227 fjx1,fjy1,fjz1,fjx2,fjy2,fjz2,fjx3,fjy3,fjz3);
2229 /* Inner loop uses 288 flops */
2232 /* End of innermost loop */
2234 gmx_mm_update_iforce_3atom_swizzle_ps(fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
2235 f+i_coord_offset+DIM,fshift+i_shift_offset);
2237 /* Increment number of inner iterations */
2238 inneriter += j_index_end - j_index_start;
2240 /* Outer loop uses 18 flops */
2243 /* Increment number of outer iterations */
2246 /* Update outer/inner flops */
2248 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_W4W4_F,outeriter*18 + inneriter*288);