2 * Note: this file was generated by the Gromacs avx_256_double 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_256_double.h"
34 #include "kernelutil_x86_avx_256_double.h"
37 * Gromacs nonbonded kernel: nb_kernel_ElecEwSh_VdwNone_GeomP1P1_VF_avx_256_double
38 * Electrostatics interaction: Ewald
39 * VdW interaction: None
40 * Geometry: Particle-Particle
41 * Calculate force/pot: PotentialAndForce
44 nb_kernel_ElecEwSh_VdwNone_GeomP1P1_VF_avx_256_double
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, 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 jnrlistE,jnrlistF,jnrlistG,jnrlistH;
63 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
64 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
66 real *shiftvec,*fshift,*x,*f;
67 real *fjptrA,*fjptrB,*fjptrC,*fjptrD;
69 __m256d tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
70 real * vdwioffsetptr0;
71 __m256d ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
72 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
73 __m256d jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
74 __m256d dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
75 __m256d velec,felec,velecsum,facel,crf,krf,krf2;
78 __m256d ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
79 __m256d beta,beta2,beta3,zeta2,pmecorrF,pmecorrV,rinv3;
81 __m256d dummy_mask,cutoff_mask;
82 __m128 tmpmask0,tmpmask1;
83 __m256d signbit = _mm256_castsi256_pd( _mm256_set1_epi32(0x80000000) );
84 __m256d one = _mm256_set1_pd(1.0);
85 __m256d two = _mm256_set1_pd(2.0);
91 jindex = nlist->jindex;
93 shiftidx = nlist->shift;
95 shiftvec = fr->shift_vec[0];
96 fshift = fr->fshift[0];
97 facel = _mm256_set1_pd(fr->epsfac);
98 charge = mdatoms->chargeA;
100 sh_ewald = _mm256_set1_pd(fr->ic->sh_ewald);
101 beta = _mm256_set1_pd(fr->ic->ewaldcoeff);
102 beta2 = _mm256_mul_pd(beta,beta);
103 beta3 = _mm256_mul_pd(beta,beta2);
105 ewtab = fr->ic->tabq_coul_FDV0;
106 ewtabscale = _mm256_set1_pd(fr->ic->tabq_scale);
107 ewtabhalfspace = _mm256_set1_pd(0.5/fr->ic->tabq_scale);
109 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
110 rcutoff_scalar = fr->rcoulomb;
111 rcutoff = _mm256_set1_pd(rcutoff_scalar);
112 rcutoff2 = _mm256_mul_pd(rcutoff,rcutoff);
114 /* Avoid stupid compiler warnings */
115 jnrA = jnrB = jnrC = jnrD = 0;
124 for(iidx=0;iidx<4*DIM;iidx++)
129 /* Start outer loop over neighborlists */
130 for(iidx=0; iidx<nri; iidx++)
132 /* Load shift vector for this list */
133 i_shift_offset = DIM*shiftidx[iidx];
135 /* Load limits for loop over neighbors */
136 j_index_start = jindex[iidx];
137 j_index_end = jindex[iidx+1];
139 /* Get outer coordinate index */
141 i_coord_offset = DIM*inr;
143 /* Load i particle coords and add shift vector */
144 gmx_mm256_load_shift_and_1rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
146 fix0 = _mm256_setzero_pd();
147 fiy0 = _mm256_setzero_pd();
148 fiz0 = _mm256_setzero_pd();
150 /* Load parameters for i particles */
151 iq0 = _mm256_mul_pd(facel,_mm256_set1_pd(charge[inr+0]));
153 /* Reset potential sums */
154 velecsum = _mm256_setzero_pd();
156 /* Start inner kernel loop */
157 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
160 /* Get j neighbor index, and coordinate index */
165 j_coord_offsetA = DIM*jnrA;
166 j_coord_offsetB = DIM*jnrB;
167 j_coord_offsetC = DIM*jnrC;
168 j_coord_offsetD = DIM*jnrD;
170 /* load j atom coordinates */
171 gmx_mm256_load_1rvec_4ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
172 x+j_coord_offsetC,x+j_coord_offsetD,
175 /* Calculate displacement vector */
176 dx00 = _mm256_sub_pd(ix0,jx0);
177 dy00 = _mm256_sub_pd(iy0,jy0);
178 dz00 = _mm256_sub_pd(iz0,jz0);
180 /* Calculate squared distance and things based on it */
181 rsq00 = gmx_mm256_calc_rsq_pd(dx00,dy00,dz00);
183 rinv00 = gmx_mm256_invsqrt_pd(rsq00);
185 rinvsq00 = _mm256_mul_pd(rinv00,rinv00);
187 /* Load parameters for j particles */
188 jq0 = gmx_mm256_load_4real_swizzle_pd(charge+jnrA+0,charge+jnrB+0,
189 charge+jnrC+0,charge+jnrD+0);
191 /**************************
192 * CALCULATE INTERACTIONS *
193 **************************/
195 if (gmx_mm256_any_lt(rsq00,rcutoff2))
198 r00 = _mm256_mul_pd(rsq00,rinv00);
200 /* Compute parameters for interactions between i and j atoms */
201 qq00 = _mm256_mul_pd(iq0,jq0);
203 /* EWALD ELECTROSTATICS */
205 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
206 ewrt = _mm256_mul_pd(r00,ewtabscale);
207 ewitab = _mm256_cvttpd_epi32(ewrt);
208 eweps = _mm256_sub_pd(ewrt,_mm256_round_pd(ewrt, _MM_FROUND_FLOOR));
209 ewitab = _mm_slli_epi32(ewitab,2);
210 ewtabF = _mm256_load_pd( ewtab + _mm_extract_epi32(ewitab,0) );
211 ewtabD = _mm256_load_pd( ewtab + _mm_extract_epi32(ewitab,1) );
212 ewtabV = _mm256_load_pd( ewtab + _mm_extract_epi32(ewitab,2) );
213 ewtabFn = _mm256_load_pd( ewtab + _mm_extract_epi32(ewitab,3) );
214 GMX_MM256_FULLTRANSPOSE4_PD(ewtabF,ewtabD,ewtabV,ewtabFn);
215 felec = _mm256_add_pd(ewtabF,_mm256_mul_pd(eweps,ewtabD));
216 velec = _mm256_sub_pd(ewtabV,_mm256_mul_pd(_mm256_mul_pd(ewtabhalfspace,eweps),_mm256_add_pd(ewtabF,felec)));
217 velec = _mm256_mul_pd(qq00,_mm256_sub_pd(_mm256_sub_pd(rinv00,sh_ewald),velec));
218 felec = _mm256_mul_pd(_mm256_mul_pd(qq00,rinv00),_mm256_sub_pd(rinvsq00,felec));
220 cutoff_mask = _mm256_cmp_pd(rsq00,rcutoff2,_CMP_LT_OQ);
222 /* Update potential sum for this i atom from the interaction with this j atom. */
223 velec = _mm256_and_pd(velec,cutoff_mask);
224 velecsum = _mm256_add_pd(velecsum,velec);
228 fscal = _mm256_and_pd(fscal,cutoff_mask);
230 /* Calculate temporary vectorial force */
231 tx = _mm256_mul_pd(fscal,dx00);
232 ty = _mm256_mul_pd(fscal,dy00);
233 tz = _mm256_mul_pd(fscal,dz00);
235 /* Update vectorial force */
236 fix0 = _mm256_add_pd(fix0,tx);
237 fiy0 = _mm256_add_pd(fiy0,ty);
238 fiz0 = _mm256_add_pd(fiz0,tz);
240 fjptrA = f+j_coord_offsetA;
241 fjptrB = f+j_coord_offsetB;
242 fjptrC = f+j_coord_offsetC;
243 fjptrD = f+j_coord_offsetD;
244 gmx_mm256_decrement_1rvec_4ptr_swizzle_pd(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
248 /* Inner loop uses 46 flops */
254 /* Get j neighbor index, and coordinate index */
255 jnrlistA = jjnr[jidx];
256 jnrlistB = jjnr[jidx+1];
257 jnrlistC = jjnr[jidx+2];
258 jnrlistD = jjnr[jidx+3];
259 /* Sign of each element will be negative for non-real atoms.
260 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
261 * so use it as val = _mm_andnot_pd(mask,val) to clear dummy entries.
263 tmpmask0 = gmx_mm_castsi128_pd(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
265 tmpmask1 = _mm_permute_ps(tmpmask0,_GMX_MM_PERMUTE(3,3,2,2));
266 tmpmask0 = _mm_permute_ps(tmpmask0,_GMX_MM_PERMUTE(1,1,0,0));
267 dummy_mask = _mm256_castps_pd(gmx_mm256_set_m128(tmpmask1,tmpmask0));
269 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
270 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
271 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
272 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
273 j_coord_offsetA = DIM*jnrA;
274 j_coord_offsetB = DIM*jnrB;
275 j_coord_offsetC = DIM*jnrC;
276 j_coord_offsetD = DIM*jnrD;
278 /* load j atom coordinates */
279 gmx_mm256_load_1rvec_4ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
280 x+j_coord_offsetC,x+j_coord_offsetD,
283 /* Calculate displacement vector */
284 dx00 = _mm256_sub_pd(ix0,jx0);
285 dy00 = _mm256_sub_pd(iy0,jy0);
286 dz00 = _mm256_sub_pd(iz0,jz0);
288 /* Calculate squared distance and things based on it */
289 rsq00 = gmx_mm256_calc_rsq_pd(dx00,dy00,dz00);
291 rinv00 = gmx_mm256_invsqrt_pd(rsq00);
293 rinvsq00 = _mm256_mul_pd(rinv00,rinv00);
295 /* Load parameters for j particles */
296 jq0 = gmx_mm256_load_4real_swizzle_pd(charge+jnrA+0,charge+jnrB+0,
297 charge+jnrC+0,charge+jnrD+0);
299 /**************************
300 * CALCULATE INTERACTIONS *
301 **************************/
303 if (gmx_mm256_any_lt(rsq00,rcutoff2))
306 r00 = _mm256_mul_pd(rsq00,rinv00);
307 r00 = _mm256_andnot_pd(dummy_mask,r00);
309 /* Compute parameters for interactions between i and j atoms */
310 qq00 = _mm256_mul_pd(iq0,jq0);
312 /* EWALD ELECTROSTATICS */
314 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
315 ewrt = _mm256_mul_pd(r00,ewtabscale);
316 ewitab = _mm256_cvttpd_epi32(ewrt);
317 eweps = _mm256_sub_pd(ewrt,_mm256_round_pd(ewrt, _MM_FROUND_FLOOR));
318 ewitab = _mm_slli_epi32(ewitab,2);
319 ewtabF = _mm256_load_pd( ewtab + _mm_extract_epi32(ewitab,0) );
320 ewtabD = _mm256_load_pd( ewtab + _mm_extract_epi32(ewitab,1) );
321 ewtabV = _mm256_load_pd( ewtab + _mm_extract_epi32(ewitab,2) );
322 ewtabFn = _mm256_load_pd( ewtab + _mm_extract_epi32(ewitab,3) );
323 GMX_MM256_FULLTRANSPOSE4_PD(ewtabF,ewtabD,ewtabV,ewtabFn);
324 felec = _mm256_add_pd(ewtabF,_mm256_mul_pd(eweps,ewtabD));
325 velec = _mm256_sub_pd(ewtabV,_mm256_mul_pd(_mm256_mul_pd(ewtabhalfspace,eweps),_mm256_add_pd(ewtabF,felec)));
326 velec = _mm256_mul_pd(qq00,_mm256_sub_pd(_mm256_sub_pd(rinv00,sh_ewald),velec));
327 felec = _mm256_mul_pd(_mm256_mul_pd(qq00,rinv00),_mm256_sub_pd(rinvsq00,felec));
329 cutoff_mask = _mm256_cmp_pd(rsq00,rcutoff2,_CMP_LT_OQ);
331 /* Update potential sum for this i atom from the interaction with this j atom. */
332 velec = _mm256_and_pd(velec,cutoff_mask);
333 velec = _mm256_andnot_pd(dummy_mask,velec);
334 velecsum = _mm256_add_pd(velecsum,velec);
338 fscal = _mm256_and_pd(fscal,cutoff_mask);
340 fscal = _mm256_andnot_pd(dummy_mask,fscal);
342 /* Calculate temporary vectorial force */
343 tx = _mm256_mul_pd(fscal,dx00);
344 ty = _mm256_mul_pd(fscal,dy00);
345 tz = _mm256_mul_pd(fscal,dz00);
347 /* Update vectorial force */
348 fix0 = _mm256_add_pd(fix0,tx);
349 fiy0 = _mm256_add_pd(fiy0,ty);
350 fiz0 = _mm256_add_pd(fiz0,tz);
352 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
353 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
354 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
355 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
356 gmx_mm256_decrement_1rvec_4ptr_swizzle_pd(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
360 /* Inner loop uses 47 flops */
363 /* End of innermost loop */
365 gmx_mm256_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
366 f+i_coord_offset,fshift+i_shift_offset);
369 /* Update potential energies */
370 gmx_mm256_update_1pot_pd(velecsum,kernel_data->energygrp_elec+ggid);
372 /* Increment number of inner iterations */
373 inneriter += j_index_end - j_index_start;
375 /* Outer loop uses 8 flops */
378 /* Increment number of outer iterations */
381 /* Update outer/inner flops */
383 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VF,outeriter*8 + inneriter*47);
386 * Gromacs nonbonded kernel: nb_kernel_ElecEwSh_VdwNone_GeomP1P1_F_avx_256_double
387 * Electrostatics interaction: Ewald
388 * VdW interaction: None
389 * Geometry: Particle-Particle
390 * Calculate force/pot: Force
393 nb_kernel_ElecEwSh_VdwNone_GeomP1P1_F_avx_256_double
394 (t_nblist * gmx_restrict nlist,
395 rvec * gmx_restrict xx,
396 rvec * gmx_restrict ff,
397 t_forcerec * gmx_restrict fr,
398 t_mdatoms * gmx_restrict mdatoms,
399 nb_kernel_data_t * gmx_restrict kernel_data,
400 t_nrnb * gmx_restrict nrnb)
402 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
403 * just 0 for non-waters.
404 * Suffixes A,B,C,D refer to j loop unrolling done with AVX, e.g. for the four different
405 * jnr indices corresponding to data put in the four positions in the SIMD register.
407 int i_shift_offset,i_coord_offset,outeriter,inneriter;
408 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
409 int jnrA,jnrB,jnrC,jnrD;
410 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
411 int jnrlistE,jnrlistF,jnrlistG,jnrlistH;
412 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
413 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
415 real *shiftvec,*fshift,*x,*f;
416 real *fjptrA,*fjptrB,*fjptrC,*fjptrD;
418 __m256d tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
419 real * vdwioffsetptr0;
420 __m256d ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
421 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
422 __m256d jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
423 __m256d dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
424 __m256d velec,felec,velecsum,facel,crf,krf,krf2;
427 __m256d ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
428 __m256d beta,beta2,beta3,zeta2,pmecorrF,pmecorrV,rinv3;
430 __m256d dummy_mask,cutoff_mask;
431 __m128 tmpmask0,tmpmask1;
432 __m256d signbit = _mm256_castsi256_pd( _mm256_set1_epi32(0x80000000) );
433 __m256d one = _mm256_set1_pd(1.0);
434 __m256d two = _mm256_set1_pd(2.0);
440 jindex = nlist->jindex;
442 shiftidx = nlist->shift;
444 shiftvec = fr->shift_vec[0];
445 fshift = fr->fshift[0];
446 facel = _mm256_set1_pd(fr->epsfac);
447 charge = mdatoms->chargeA;
449 sh_ewald = _mm256_set1_pd(fr->ic->sh_ewald);
450 beta = _mm256_set1_pd(fr->ic->ewaldcoeff);
451 beta2 = _mm256_mul_pd(beta,beta);
452 beta3 = _mm256_mul_pd(beta,beta2);
454 ewtab = fr->ic->tabq_coul_F;
455 ewtabscale = _mm256_set1_pd(fr->ic->tabq_scale);
456 ewtabhalfspace = _mm256_set1_pd(0.5/fr->ic->tabq_scale);
458 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
459 rcutoff_scalar = fr->rcoulomb;
460 rcutoff = _mm256_set1_pd(rcutoff_scalar);
461 rcutoff2 = _mm256_mul_pd(rcutoff,rcutoff);
463 /* Avoid stupid compiler warnings */
464 jnrA = jnrB = jnrC = jnrD = 0;
473 for(iidx=0;iidx<4*DIM;iidx++)
478 /* Start outer loop over neighborlists */
479 for(iidx=0; iidx<nri; iidx++)
481 /* Load shift vector for this list */
482 i_shift_offset = DIM*shiftidx[iidx];
484 /* Load limits for loop over neighbors */
485 j_index_start = jindex[iidx];
486 j_index_end = jindex[iidx+1];
488 /* Get outer coordinate index */
490 i_coord_offset = DIM*inr;
492 /* Load i particle coords and add shift vector */
493 gmx_mm256_load_shift_and_1rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
495 fix0 = _mm256_setzero_pd();
496 fiy0 = _mm256_setzero_pd();
497 fiz0 = _mm256_setzero_pd();
499 /* Load parameters for i particles */
500 iq0 = _mm256_mul_pd(facel,_mm256_set1_pd(charge[inr+0]));
502 /* Start inner kernel loop */
503 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
506 /* Get j neighbor index, and coordinate index */
511 j_coord_offsetA = DIM*jnrA;
512 j_coord_offsetB = DIM*jnrB;
513 j_coord_offsetC = DIM*jnrC;
514 j_coord_offsetD = DIM*jnrD;
516 /* load j atom coordinates */
517 gmx_mm256_load_1rvec_4ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
518 x+j_coord_offsetC,x+j_coord_offsetD,
521 /* Calculate displacement vector */
522 dx00 = _mm256_sub_pd(ix0,jx0);
523 dy00 = _mm256_sub_pd(iy0,jy0);
524 dz00 = _mm256_sub_pd(iz0,jz0);
526 /* Calculate squared distance and things based on it */
527 rsq00 = gmx_mm256_calc_rsq_pd(dx00,dy00,dz00);
529 rinv00 = gmx_mm256_invsqrt_pd(rsq00);
531 rinvsq00 = _mm256_mul_pd(rinv00,rinv00);
533 /* Load parameters for j particles */
534 jq0 = gmx_mm256_load_4real_swizzle_pd(charge+jnrA+0,charge+jnrB+0,
535 charge+jnrC+0,charge+jnrD+0);
537 /**************************
538 * CALCULATE INTERACTIONS *
539 **************************/
541 if (gmx_mm256_any_lt(rsq00,rcutoff2))
544 r00 = _mm256_mul_pd(rsq00,rinv00);
546 /* Compute parameters for interactions between i and j atoms */
547 qq00 = _mm256_mul_pd(iq0,jq0);
549 /* EWALD ELECTROSTATICS */
551 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
552 ewrt = _mm256_mul_pd(r00,ewtabscale);
553 ewitab = _mm256_cvttpd_epi32(ewrt);
554 eweps = _mm256_sub_pd(ewrt,_mm256_round_pd(ewrt, _MM_FROUND_FLOOR));
555 gmx_mm256_load_4pair_swizzle_pd(ewtab + _mm_extract_epi32(ewitab,0),ewtab + _mm_extract_epi32(ewitab,1),
556 ewtab + _mm_extract_epi32(ewitab,2),ewtab + _mm_extract_epi32(ewitab,3),
558 felec = _mm256_add_pd(_mm256_mul_pd( _mm256_sub_pd(one,eweps),ewtabF),_mm256_mul_pd(eweps,ewtabFn));
559 felec = _mm256_mul_pd(_mm256_mul_pd(qq00,rinv00),_mm256_sub_pd(rinvsq00,felec));
561 cutoff_mask = _mm256_cmp_pd(rsq00,rcutoff2,_CMP_LT_OQ);
565 fscal = _mm256_and_pd(fscal,cutoff_mask);
567 /* Calculate temporary vectorial force */
568 tx = _mm256_mul_pd(fscal,dx00);
569 ty = _mm256_mul_pd(fscal,dy00);
570 tz = _mm256_mul_pd(fscal,dz00);
572 /* Update vectorial force */
573 fix0 = _mm256_add_pd(fix0,tx);
574 fiy0 = _mm256_add_pd(fiy0,ty);
575 fiz0 = _mm256_add_pd(fiz0,tz);
577 fjptrA = f+j_coord_offsetA;
578 fjptrB = f+j_coord_offsetB;
579 fjptrC = f+j_coord_offsetC;
580 fjptrD = f+j_coord_offsetD;
581 gmx_mm256_decrement_1rvec_4ptr_swizzle_pd(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
585 /* Inner loop uses 39 flops */
591 /* Get j neighbor index, and coordinate index */
592 jnrlistA = jjnr[jidx];
593 jnrlistB = jjnr[jidx+1];
594 jnrlistC = jjnr[jidx+2];
595 jnrlistD = jjnr[jidx+3];
596 /* Sign of each element will be negative for non-real atoms.
597 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
598 * so use it as val = _mm_andnot_pd(mask,val) to clear dummy entries.
600 tmpmask0 = gmx_mm_castsi128_pd(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
602 tmpmask1 = _mm_permute_ps(tmpmask0,_GMX_MM_PERMUTE(3,3,2,2));
603 tmpmask0 = _mm_permute_ps(tmpmask0,_GMX_MM_PERMUTE(1,1,0,0));
604 dummy_mask = _mm256_castps_pd(gmx_mm256_set_m128(tmpmask1,tmpmask0));
606 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
607 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
608 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
609 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
610 j_coord_offsetA = DIM*jnrA;
611 j_coord_offsetB = DIM*jnrB;
612 j_coord_offsetC = DIM*jnrC;
613 j_coord_offsetD = DIM*jnrD;
615 /* load j atom coordinates */
616 gmx_mm256_load_1rvec_4ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
617 x+j_coord_offsetC,x+j_coord_offsetD,
620 /* Calculate displacement vector */
621 dx00 = _mm256_sub_pd(ix0,jx0);
622 dy00 = _mm256_sub_pd(iy0,jy0);
623 dz00 = _mm256_sub_pd(iz0,jz0);
625 /* Calculate squared distance and things based on it */
626 rsq00 = gmx_mm256_calc_rsq_pd(dx00,dy00,dz00);
628 rinv00 = gmx_mm256_invsqrt_pd(rsq00);
630 rinvsq00 = _mm256_mul_pd(rinv00,rinv00);
632 /* Load parameters for j particles */
633 jq0 = gmx_mm256_load_4real_swizzle_pd(charge+jnrA+0,charge+jnrB+0,
634 charge+jnrC+0,charge+jnrD+0);
636 /**************************
637 * CALCULATE INTERACTIONS *
638 **************************/
640 if (gmx_mm256_any_lt(rsq00,rcutoff2))
643 r00 = _mm256_mul_pd(rsq00,rinv00);
644 r00 = _mm256_andnot_pd(dummy_mask,r00);
646 /* Compute parameters for interactions between i and j atoms */
647 qq00 = _mm256_mul_pd(iq0,jq0);
649 /* EWALD ELECTROSTATICS */
651 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
652 ewrt = _mm256_mul_pd(r00,ewtabscale);
653 ewitab = _mm256_cvttpd_epi32(ewrt);
654 eweps = _mm256_sub_pd(ewrt,_mm256_round_pd(ewrt, _MM_FROUND_FLOOR));
655 gmx_mm256_load_4pair_swizzle_pd(ewtab + _mm_extract_epi32(ewitab,0),ewtab + _mm_extract_epi32(ewitab,1),
656 ewtab + _mm_extract_epi32(ewitab,2),ewtab + _mm_extract_epi32(ewitab,3),
658 felec = _mm256_add_pd(_mm256_mul_pd( _mm256_sub_pd(one,eweps),ewtabF),_mm256_mul_pd(eweps,ewtabFn));
659 felec = _mm256_mul_pd(_mm256_mul_pd(qq00,rinv00),_mm256_sub_pd(rinvsq00,felec));
661 cutoff_mask = _mm256_cmp_pd(rsq00,rcutoff2,_CMP_LT_OQ);
665 fscal = _mm256_and_pd(fscal,cutoff_mask);
667 fscal = _mm256_andnot_pd(dummy_mask,fscal);
669 /* Calculate temporary vectorial force */
670 tx = _mm256_mul_pd(fscal,dx00);
671 ty = _mm256_mul_pd(fscal,dy00);
672 tz = _mm256_mul_pd(fscal,dz00);
674 /* Update vectorial force */
675 fix0 = _mm256_add_pd(fix0,tx);
676 fiy0 = _mm256_add_pd(fiy0,ty);
677 fiz0 = _mm256_add_pd(fiz0,tz);
679 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
680 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
681 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
682 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
683 gmx_mm256_decrement_1rvec_4ptr_swizzle_pd(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
687 /* Inner loop uses 40 flops */
690 /* End of innermost loop */
692 gmx_mm256_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
693 f+i_coord_offset,fshift+i_shift_offset);
695 /* Increment number of inner iterations */
696 inneriter += j_index_end - j_index_start;
698 /* Outer loop uses 7 flops */
701 /* Increment number of outer iterations */
704 /* Update outer/inner flops */
706 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_F,outeriter*7 + inneriter*40);