2 * Note: this file was generated by the Gromacs avx_256_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_256_single.h"
34 #include "kernelutil_x86_avx_256_single.h"
37 * Gromacs nonbonded kernel: nb_kernel_ElecEw_VdwNone_GeomP1P1_VF_avx_256_single
38 * Electrostatics interaction: Ewald
39 * VdW interaction: None
40 * Geometry: Particle-Particle
41 * Calculate force/pot: PotentialAndForce
44 nb_kernel_ElecEw_VdwNone_GeomP1P1_VF_avx_256_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,E,F,G,H refer to j loop unrolling done with AVX, e.g. for the eight 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 jnrE,jnrF,jnrG,jnrH;
62 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
63 int jnrlistE,jnrlistF,jnrlistG,jnrlistH;
64 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
65 int j_coord_offsetE,j_coord_offsetF,j_coord_offsetG,j_coord_offsetH;
66 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
68 real *shiftvec,*fshift,*x,*f;
69 real *fjptrA,*fjptrB,*fjptrC,*fjptrD,*fjptrE,*fjptrF,*fjptrG,*fjptrH;
71 __m256 tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
72 real * vdwioffsetptr0;
73 __m256 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
74 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D,vdwjidx0E,vdwjidx0F,vdwjidx0G,vdwjidx0H;
75 __m256 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
76 __m256 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
77 __m256 velec,felec,velecsum,facel,crf,krf,krf2;
80 __m128i ewitab_lo,ewitab_hi;
81 __m256 ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
82 __m256 beta,beta2,beta3,zeta2,pmecorrF,pmecorrV,rinv3;
84 __m256 dummy_mask,cutoff_mask;
85 __m256 signbit = _mm256_castsi256_ps( _mm256_set1_epi32(0x80000000) );
86 __m256 one = _mm256_set1_ps(1.0);
87 __m256 two = _mm256_set1_ps(2.0);
93 jindex = nlist->jindex;
95 shiftidx = nlist->shift;
97 shiftvec = fr->shift_vec[0];
98 fshift = fr->fshift[0];
99 facel = _mm256_set1_ps(fr->epsfac);
100 charge = mdatoms->chargeA;
102 sh_ewald = _mm256_set1_ps(fr->ic->sh_ewald);
103 beta = _mm256_set1_ps(fr->ic->ewaldcoeff);
104 beta2 = _mm256_mul_ps(beta,beta);
105 beta3 = _mm256_mul_ps(beta,beta2);
107 ewtab = fr->ic->tabq_coul_FDV0;
108 ewtabscale = _mm256_set1_ps(fr->ic->tabq_scale);
109 ewtabhalfspace = _mm256_set1_ps(0.5/fr->ic->tabq_scale);
111 /* Avoid stupid compiler warnings */
112 jnrA = jnrB = jnrC = jnrD = jnrE = jnrF = jnrG = jnrH = 0;
125 for(iidx=0;iidx<4*DIM;iidx++)
130 /* Start outer loop over neighborlists */
131 for(iidx=0; iidx<nri; iidx++)
133 /* Load shift vector for this list */
134 i_shift_offset = DIM*shiftidx[iidx];
136 /* Load limits for loop over neighbors */
137 j_index_start = jindex[iidx];
138 j_index_end = jindex[iidx+1];
140 /* Get outer coordinate index */
142 i_coord_offset = DIM*inr;
144 /* Load i particle coords and add shift vector */
145 gmx_mm256_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
147 fix0 = _mm256_setzero_ps();
148 fiy0 = _mm256_setzero_ps();
149 fiz0 = _mm256_setzero_ps();
151 /* Load parameters for i particles */
152 iq0 = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+0]));
154 /* Reset potential sums */
155 velecsum = _mm256_setzero_ps();
157 /* Start inner kernel loop */
158 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+7]>=0; jidx+=8)
161 /* Get j neighbor index, and coordinate index */
170 j_coord_offsetA = DIM*jnrA;
171 j_coord_offsetB = DIM*jnrB;
172 j_coord_offsetC = DIM*jnrC;
173 j_coord_offsetD = DIM*jnrD;
174 j_coord_offsetE = DIM*jnrE;
175 j_coord_offsetF = DIM*jnrF;
176 j_coord_offsetG = DIM*jnrG;
177 j_coord_offsetH = DIM*jnrH;
179 /* load j atom coordinates */
180 gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
181 x+j_coord_offsetC,x+j_coord_offsetD,
182 x+j_coord_offsetE,x+j_coord_offsetF,
183 x+j_coord_offsetG,x+j_coord_offsetH,
186 /* Calculate displacement vector */
187 dx00 = _mm256_sub_ps(ix0,jx0);
188 dy00 = _mm256_sub_ps(iy0,jy0);
189 dz00 = _mm256_sub_ps(iz0,jz0);
191 /* Calculate squared distance and things based on it */
192 rsq00 = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
194 rinv00 = gmx_mm256_invsqrt_ps(rsq00);
196 rinvsq00 = _mm256_mul_ps(rinv00,rinv00);
198 /* Load parameters for j particles */
199 jq0 = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
200 charge+jnrC+0,charge+jnrD+0,
201 charge+jnrE+0,charge+jnrF+0,
202 charge+jnrG+0,charge+jnrH+0);
204 /**************************
205 * CALCULATE INTERACTIONS *
206 **************************/
208 r00 = _mm256_mul_ps(rsq00,rinv00);
210 /* Compute parameters for interactions between i and j atoms */
211 qq00 = _mm256_mul_ps(iq0,jq0);
213 /* EWALD ELECTROSTATICS */
215 /* Analytical PME correction */
216 zeta2 = _mm256_mul_ps(beta2,rsq00);
217 rinv3 = _mm256_mul_ps(rinvsq00,rinv00);
218 pmecorrF = gmx_mm256_pmecorrF_ps(zeta2);
219 felec = _mm256_add_ps( _mm256_mul_ps(pmecorrF,beta3), rinv3);
220 felec = _mm256_mul_ps(qq00,felec);
221 pmecorrV = gmx_mm256_pmecorrV_ps(zeta2);
222 pmecorrV = _mm256_mul_ps(pmecorrV,beta);
223 velec = _mm256_sub_ps(rinv00,pmecorrV);
224 velec = _mm256_mul_ps(qq00,velec);
226 /* Update potential sum for this i atom from the interaction with this j atom. */
227 velecsum = _mm256_add_ps(velecsum,velec);
231 /* Calculate temporary vectorial force */
232 tx = _mm256_mul_ps(fscal,dx00);
233 ty = _mm256_mul_ps(fscal,dy00);
234 tz = _mm256_mul_ps(fscal,dz00);
236 /* Update vectorial force */
237 fix0 = _mm256_add_ps(fix0,tx);
238 fiy0 = _mm256_add_ps(fiy0,ty);
239 fiz0 = _mm256_add_ps(fiz0,tz);
241 fjptrA = f+j_coord_offsetA;
242 fjptrB = f+j_coord_offsetB;
243 fjptrC = f+j_coord_offsetC;
244 fjptrD = f+j_coord_offsetD;
245 fjptrE = f+j_coord_offsetE;
246 fjptrF = f+j_coord_offsetF;
247 fjptrG = f+j_coord_offsetG;
248 fjptrH = f+j_coord_offsetH;
249 gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,tx,ty,tz);
251 /* Inner loop uses 84 flops */
257 /* Get j neighbor index, and coordinate index */
258 jnrlistA = jjnr[jidx];
259 jnrlistB = jjnr[jidx+1];
260 jnrlistC = jjnr[jidx+2];
261 jnrlistD = jjnr[jidx+3];
262 jnrlistE = jjnr[jidx+4];
263 jnrlistF = jjnr[jidx+5];
264 jnrlistG = jjnr[jidx+6];
265 jnrlistH = jjnr[jidx+7];
266 /* Sign of each element will be negative for non-real atoms.
267 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
268 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
270 dummy_mask = gmx_mm256_set_m128(gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx+4)),_mm_setzero_si128())),
271 gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128())));
273 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
274 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
275 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
276 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
277 jnrE = (jnrlistE>=0) ? jnrlistE : 0;
278 jnrF = (jnrlistF>=0) ? jnrlistF : 0;
279 jnrG = (jnrlistG>=0) ? jnrlistG : 0;
280 jnrH = (jnrlistH>=0) ? jnrlistH : 0;
281 j_coord_offsetA = DIM*jnrA;
282 j_coord_offsetB = DIM*jnrB;
283 j_coord_offsetC = DIM*jnrC;
284 j_coord_offsetD = DIM*jnrD;
285 j_coord_offsetE = DIM*jnrE;
286 j_coord_offsetF = DIM*jnrF;
287 j_coord_offsetG = DIM*jnrG;
288 j_coord_offsetH = DIM*jnrH;
290 /* load j atom coordinates */
291 gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
292 x+j_coord_offsetC,x+j_coord_offsetD,
293 x+j_coord_offsetE,x+j_coord_offsetF,
294 x+j_coord_offsetG,x+j_coord_offsetH,
297 /* Calculate displacement vector */
298 dx00 = _mm256_sub_ps(ix0,jx0);
299 dy00 = _mm256_sub_ps(iy0,jy0);
300 dz00 = _mm256_sub_ps(iz0,jz0);
302 /* Calculate squared distance and things based on it */
303 rsq00 = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
305 rinv00 = gmx_mm256_invsqrt_ps(rsq00);
307 rinvsq00 = _mm256_mul_ps(rinv00,rinv00);
309 /* Load parameters for j particles */
310 jq0 = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
311 charge+jnrC+0,charge+jnrD+0,
312 charge+jnrE+0,charge+jnrF+0,
313 charge+jnrG+0,charge+jnrH+0);
315 /**************************
316 * CALCULATE INTERACTIONS *
317 **************************/
319 r00 = _mm256_mul_ps(rsq00,rinv00);
320 r00 = _mm256_andnot_ps(dummy_mask,r00);
322 /* Compute parameters for interactions between i and j atoms */
323 qq00 = _mm256_mul_ps(iq0,jq0);
325 /* EWALD ELECTROSTATICS */
327 /* Analytical PME correction */
328 zeta2 = _mm256_mul_ps(beta2,rsq00);
329 rinv3 = _mm256_mul_ps(rinvsq00,rinv00);
330 pmecorrF = gmx_mm256_pmecorrF_ps(zeta2);
331 felec = _mm256_add_ps( _mm256_mul_ps(pmecorrF,beta3), rinv3);
332 felec = _mm256_mul_ps(qq00,felec);
333 pmecorrV = gmx_mm256_pmecorrV_ps(zeta2);
334 pmecorrV = _mm256_mul_ps(pmecorrV,beta);
335 velec = _mm256_sub_ps(rinv00,pmecorrV);
336 velec = _mm256_mul_ps(qq00,velec);
338 /* Update potential sum for this i atom from the interaction with this j atom. */
339 velec = _mm256_andnot_ps(dummy_mask,velec);
340 velecsum = _mm256_add_ps(velecsum,velec);
344 fscal = _mm256_andnot_ps(dummy_mask,fscal);
346 /* Calculate temporary vectorial force */
347 tx = _mm256_mul_ps(fscal,dx00);
348 ty = _mm256_mul_ps(fscal,dy00);
349 tz = _mm256_mul_ps(fscal,dz00);
351 /* Update vectorial force */
352 fix0 = _mm256_add_ps(fix0,tx);
353 fiy0 = _mm256_add_ps(fiy0,ty);
354 fiz0 = _mm256_add_ps(fiz0,tz);
356 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
357 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
358 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
359 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
360 fjptrE = (jnrlistE>=0) ? f+j_coord_offsetE : scratch;
361 fjptrF = (jnrlistF>=0) ? f+j_coord_offsetF : scratch;
362 fjptrG = (jnrlistG>=0) ? f+j_coord_offsetG : scratch;
363 fjptrH = (jnrlistH>=0) ? f+j_coord_offsetH : scratch;
364 gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,tx,ty,tz);
366 /* Inner loop uses 85 flops */
369 /* End of innermost loop */
371 gmx_mm256_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
372 f+i_coord_offset,fshift+i_shift_offset);
375 /* Update potential energies */
376 gmx_mm256_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
378 /* Increment number of inner iterations */
379 inneriter += j_index_end - j_index_start;
381 /* Outer loop uses 8 flops */
384 /* Increment number of outer iterations */
387 /* Update outer/inner flops */
389 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VF,outeriter*8 + inneriter*85);
392 * Gromacs nonbonded kernel: nb_kernel_ElecEw_VdwNone_GeomP1P1_F_avx_256_single
393 * Electrostatics interaction: Ewald
394 * VdW interaction: None
395 * Geometry: Particle-Particle
396 * Calculate force/pot: Force
399 nb_kernel_ElecEw_VdwNone_GeomP1P1_F_avx_256_single
400 (t_nblist * gmx_restrict nlist,
401 rvec * gmx_restrict xx,
402 rvec * gmx_restrict ff,
403 t_forcerec * gmx_restrict fr,
404 t_mdatoms * gmx_restrict mdatoms,
405 nb_kernel_data_t * gmx_restrict kernel_data,
406 t_nrnb * gmx_restrict nrnb)
408 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
409 * just 0 for non-waters.
410 * Suffixes A,B,C,D,E,F,G,H refer to j loop unrolling done with AVX, e.g. for the eight different
411 * jnr indices corresponding to data put in the four positions in the SIMD register.
413 int i_shift_offset,i_coord_offset,outeriter,inneriter;
414 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
415 int jnrA,jnrB,jnrC,jnrD;
416 int jnrE,jnrF,jnrG,jnrH;
417 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
418 int jnrlistE,jnrlistF,jnrlistG,jnrlistH;
419 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
420 int j_coord_offsetE,j_coord_offsetF,j_coord_offsetG,j_coord_offsetH;
421 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
423 real *shiftvec,*fshift,*x,*f;
424 real *fjptrA,*fjptrB,*fjptrC,*fjptrD,*fjptrE,*fjptrF,*fjptrG,*fjptrH;
426 __m256 tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
427 real * vdwioffsetptr0;
428 __m256 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
429 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D,vdwjidx0E,vdwjidx0F,vdwjidx0G,vdwjidx0H;
430 __m256 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
431 __m256 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
432 __m256 velec,felec,velecsum,facel,crf,krf,krf2;
435 __m128i ewitab_lo,ewitab_hi;
436 __m256 ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
437 __m256 beta,beta2,beta3,zeta2,pmecorrF,pmecorrV,rinv3;
439 __m256 dummy_mask,cutoff_mask;
440 __m256 signbit = _mm256_castsi256_ps( _mm256_set1_epi32(0x80000000) );
441 __m256 one = _mm256_set1_ps(1.0);
442 __m256 two = _mm256_set1_ps(2.0);
448 jindex = nlist->jindex;
450 shiftidx = nlist->shift;
452 shiftvec = fr->shift_vec[0];
453 fshift = fr->fshift[0];
454 facel = _mm256_set1_ps(fr->epsfac);
455 charge = mdatoms->chargeA;
457 sh_ewald = _mm256_set1_ps(fr->ic->sh_ewald);
458 beta = _mm256_set1_ps(fr->ic->ewaldcoeff);
459 beta2 = _mm256_mul_ps(beta,beta);
460 beta3 = _mm256_mul_ps(beta,beta2);
462 ewtab = fr->ic->tabq_coul_F;
463 ewtabscale = _mm256_set1_ps(fr->ic->tabq_scale);
464 ewtabhalfspace = _mm256_set1_ps(0.5/fr->ic->tabq_scale);
466 /* Avoid stupid compiler warnings */
467 jnrA = jnrB = jnrC = jnrD = jnrE = jnrF = jnrG = jnrH = 0;
480 for(iidx=0;iidx<4*DIM;iidx++)
485 /* Start outer loop over neighborlists */
486 for(iidx=0; iidx<nri; iidx++)
488 /* Load shift vector for this list */
489 i_shift_offset = DIM*shiftidx[iidx];
491 /* Load limits for loop over neighbors */
492 j_index_start = jindex[iidx];
493 j_index_end = jindex[iidx+1];
495 /* Get outer coordinate index */
497 i_coord_offset = DIM*inr;
499 /* Load i particle coords and add shift vector */
500 gmx_mm256_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
502 fix0 = _mm256_setzero_ps();
503 fiy0 = _mm256_setzero_ps();
504 fiz0 = _mm256_setzero_ps();
506 /* Load parameters for i particles */
507 iq0 = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+0]));
509 /* Start inner kernel loop */
510 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+7]>=0; jidx+=8)
513 /* Get j neighbor index, and coordinate index */
522 j_coord_offsetA = DIM*jnrA;
523 j_coord_offsetB = DIM*jnrB;
524 j_coord_offsetC = DIM*jnrC;
525 j_coord_offsetD = DIM*jnrD;
526 j_coord_offsetE = DIM*jnrE;
527 j_coord_offsetF = DIM*jnrF;
528 j_coord_offsetG = DIM*jnrG;
529 j_coord_offsetH = DIM*jnrH;
531 /* load j atom coordinates */
532 gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
533 x+j_coord_offsetC,x+j_coord_offsetD,
534 x+j_coord_offsetE,x+j_coord_offsetF,
535 x+j_coord_offsetG,x+j_coord_offsetH,
538 /* Calculate displacement vector */
539 dx00 = _mm256_sub_ps(ix0,jx0);
540 dy00 = _mm256_sub_ps(iy0,jy0);
541 dz00 = _mm256_sub_ps(iz0,jz0);
543 /* Calculate squared distance and things based on it */
544 rsq00 = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
546 rinv00 = gmx_mm256_invsqrt_ps(rsq00);
548 rinvsq00 = _mm256_mul_ps(rinv00,rinv00);
550 /* Load parameters for j particles */
551 jq0 = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
552 charge+jnrC+0,charge+jnrD+0,
553 charge+jnrE+0,charge+jnrF+0,
554 charge+jnrG+0,charge+jnrH+0);
556 /**************************
557 * CALCULATE INTERACTIONS *
558 **************************/
560 r00 = _mm256_mul_ps(rsq00,rinv00);
562 /* Compute parameters for interactions between i and j atoms */
563 qq00 = _mm256_mul_ps(iq0,jq0);
565 /* EWALD ELECTROSTATICS */
567 /* Analytical PME correction */
568 zeta2 = _mm256_mul_ps(beta2,rsq00);
569 rinv3 = _mm256_mul_ps(rinvsq00,rinv00);
570 pmecorrF = gmx_mm256_pmecorrF_ps(zeta2);
571 felec = _mm256_add_ps( _mm256_mul_ps(pmecorrF,beta3), rinv3);
572 felec = _mm256_mul_ps(qq00,felec);
576 /* Calculate temporary vectorial force */
577 tx = _mm256_mul_ps(fscal,dx00);
578 ty = _mm256_mul_ps(fscal,dy00);
579 tz = _mm256_mul_ps(fscal,dz00);
581 /* Update vectorial force */
582 fix0 = _mm256_add_ps(fix0,tx);
583 fiy0 = _mm256_add_ps(fiy0,ty);
584 fiz0 = _mm256_add_ps(fiz0,tz);
586 fjptrA = f+j_coord_offsetA;
587 fjptrB = f+j_coord_offsetB;
588 fjptrC = f+j_coord_offsetC;
589 fjptrD = f+j_coord_offsetD;
590 fjptrE = f+j_coord_offsetE;
591 fjptrF = f+j_coord_offsetF;
592 fjptrG = f+j_coord_offsetG;
593 fjptrH = f+j_coord_offsetH;
594 gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,tx,ty,tz);
596 /* Inner loop uses 56 flops */
602 /* Get j neighbor index, and coordinate index */
603 jnrlistA = jjnr[jidx];
604 jnrlistB = jjnr[jidx+1];
605 jnrlistC = jjnr[jidx+2];
606 jnrlistD = jjnr[jidx+3];
607 jnrlistE = jjnr[jidx+4];
608 jnrlistF = jjnr[jidx+5];
609 jnrlistG = jjnr[jidx+6];
610 jnrlistH = jjnr[jidx+7];
611 /* Sign of each element will be negative for non-real atoms.
612 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
613 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
615 dummy_mask = gmx_mm256_set_m128(gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx+4)),_mm_setzero_si128())),
616 gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128())));
618 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
619 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
620 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
621 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
622 jnrE = (jnrlistE>=0) ? jnrlistE : 0;
623 jnrF = (jnrlistF>=0) ? jnrlistF : 0;
624 jnrG = (jnrlistG>=0) ? jnrlistG : 0;
625 jnrH = (jnrlistH>=0) ? jnrlistH : 0;
626 j_coord_offsetA = DIM*jnrA;
627 j_coord_offsetB = DIM*jnrB;
628 j_coord_offsetC = DIM*jnrC;
629 j_coord_offsetD = DIM*jnrD;
630 j_coord_offsetE = DIM*jnrE;
631 j_coord_offsetF = DIM*jnrF;
632 j_coord_offsetG = DIM*jnrG;
633 j_coord_offsetH = DIM*jnrH;
635 /* load j atom coordinates */
636 gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
637 x+j_coord_offsetC,x+j_coord_offsetD,
638 x+j_coord_offsetE,x+j_coord_offsetF,
639 x+j_coord_offsetG,x+j_coord_offsetH,
642 /* Calculate displacement vector */
643 dx00 = _mm256_sub_ps(ix0,jx0);
644 dy00 = _mm256_sub_ps(iy0,jy0);
645 dz00 = _mm256_sub_ps(iz0,jz0);
647 /* Calculate squared distance and things based on it */
648 rsq00 = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
650 rinv00 = gmx_mm256_invsqrt_ps(rsq00);
652 rinvsq00 = _mm256_mul_ps(rinv00,rinv00);
654 /* Load parameters for j particles */
655 jq0 = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
656 charge+jnrC+0,charge+jnrD+0,
657 charge+jnrE+0,charge+jnrF+0,
658 charge+jnrG+0,charge+jnrH+0);
660 /**************************
661 * CALCULATE INTERACTIONS *
662 **************************/
664 r00 = _mm256_mul_ps(rsq00,rinv00);
665 r00 = _mm256_andnot_ps(dummy_mask,r00);
667 /* Compute parameters for interactions between i and j atoms */
668 qq00 = _mm256_mul_ps(iq0,jq0);
670 /* EWALD ELECTROSTATICS */
672 /* Analytical PME correction */
673 zeta2 = _mm256_mul_ps(beta2,rsq00);
674 rinv3 = _mm256_mul_ps(rinvsq00,rinv00);
675 pmecorrF = gmx_mm256_pmecorrF_ps(zeta2);
676 felec = _mm256_add_ps( _mm256_mul_ps(pmecorrF,beta3), rinv3);
677 felec = _mm256_mul_ps(qq00,felec);
681 fscal = _mm256_andnot_ps(dummy_mask,fscal);
683 /* Calculate temporary vectorial force */
684 tx = _mm256_mul_ps(fscal,dx00);
685 ty = _mm256_mul_ps(fscal,dy00);
686 tz = _mm256_mul_ps(fscal,dz00);
688 /* Update vectorial force */
689 fix0 = _mm256_add_ps(fix0,tx);
690 fiy0 = _mm256_add_ps(fiy0,ty);
691 fiz0 = _mm256_add_ps(fiz0,tz);
693 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
694 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
695 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
696 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
697 fjptrE = (jnrlistE>=0) ? f+j_coord_offsetE : scratch;
698 fjptrF = (jnrlistF>=0) ? f+j_coord_offsetF : scratch;
699 fjptrG = (jnrlistG>=0) ? f+j_coord_offsetG : scratch;
700 fjptrH = (jnrlistH>=0) ? f+j_coord_offsetH : scratch;
701 gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,tx,ty,tz);
703 /* Inner loop uses 57 flops */
706 /* End of innermost loop */
708 gmx_mm256_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
709 f+i_coord_offset,fshift+i_shift_offset);
711 /* Increment number of inner iterations */
712 inneriter += j_index_end - j_index_start;
714 /* Outer loop uses 7 flops */
717 /* Increment number of outer iterations */
720 /* Update outer/inner flops */
722 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_F,outeriter*7 + inneriter*57);