2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
36 * Note: this file was generated by the GROMACS avx_256_single kernel generator.
44 #include "../nb_kernel.h"
45 #include "types/simple.h"
49 #include "gromacs/simd/math_x86_avx_256_single.h"
50 #include "kernelutil_x86_avx_256_single.h"
53 * Gromacs nonbonded kernel: nb_kernel_ElecRF_VdwNone_GeomP1P1_VF_avx_256_single
54 * Electrostatics interaction: ReactionField
55 * VdW interaction: None
56 * Geometry: Particle-Particle
57 * Calculate force/pot: PotentialAndForce
60 nb_kernel_ElecRF_VdwNone_GeomP1P1_VF_avx_256_single
61 (t_nblist * gmx_restrict nlist,
62 rvec * gmx_restrict xx,
63 rvec * gmx_restrict ff,
64 t_forcerec * gmx_restrict fr,
65 t_mdatoms * gmx_restrict mdatoms,
66 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
67 t_nrnb * gmx_restrict nrnb)
69 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
70 * just 0 for non-waters.
71 * Suffixes A,B,C,D,E,F,G,H refer to j loop unrolling done with AVX, e.g. for the eight different
72 * jnr indices corresponding to data put in the four positions in the SIMD register.
74 int i_shift_offset,i_coord_offset,outeriter,inneriter;
75 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
76 int jnrA,jnrB,jnrC,jnrD;
77 int jnrE,jnrF,jnrG,jnrH;
78 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
79 int jnrlistE,jnrlistF,jnrlistG,jnrlistH;
80 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
81 int j_coord_offsetE,j_coord_offsetF,j_coord_offsetG,j_coord_offsetH;
82 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
84 real *shiftvec,*fshift,*x,*f;
85 real *fjptrA,*fjptrB,*fjptrC,*fjptrD,*fjptrE,*fjptrF,*fjptrG,*fjptrH;
87 __m256 tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
88 real * vdwioffsetptr0;
89 __m256 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
90 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D,vdwjidx0E,vdwjidx0F,vdwjidx0G,vdwjidx0H;
91 __m256 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
92 __m256 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
93 __m256 velec,felec,velecsum,facel,crf,krf,krf2;
95 __m256 dummy_mask,cutoff_mask;
96 __m256 signbit = _mm256_castsi256_ps( _mm256_set1_epi32(0x80000000) );
97 __m256 one = _mm256_set1_ps(1.0);
98 __m256 two = _mm256_set1_ps(2.0);
104 jindex = nlist->jindex;
106 shiftidx = nlist->shift;
108 shiftvec = fr->shift_vec[0];
109 fshift = fr->fshift[0];
110 facel = _mm256_set1_ps(fr->epsfac);
111 charge = mdatoms->chargeA;
112 krf = _mm256_set1_ps(fr->ic->k_rf);
113 krf2 = _mm256_set1_ps(fr->ic->k_rf*2.0);
114 crf = _mm256_set1_ps(fr->ic->c_rf);
116 /* Avoid stupid compiler warnings */
117 jnrA = jnrB = jnrC = jnrD = jnrE = jnrF = jnrG = jnrH = 0;
130 for(iidx=0;iidx<4*DIM;iidx++)
135 /* Start outer loop over neighborlists */
136 for(iidx=0; iidx<nri; iidx++)
138 /* Load shift vector for this list */
139 i_shift_offset = DIM*shiftidx[iidx];
141 /* Load limits for loop over neighbors */
142 j_index_start = jindex[iidx];
143 j_index_end = jindex[iidx+1];
145 /* Get outer coordinate index */
147 i_coord_offset = DIM*inr;
149 /* Load i particle coords and add shift vector */
150 gmx_mm256_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
152 fix0 = _mm256_setzero_ps();
153 fiy0 = _mm256_setzero_ps();
154 fiz0 = _mm256_setzero_ps();
156 /* Load parameters for i particles */
157 iq0 = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+0]));
159 /* Reset potential sums */
160 velecsum = _mm256_setzero_ps();
162 /* Start inner kernel loop */
163 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+7]>=0; jidx+=8)
166 /* Get j neighbor index, and coordinate index */
175 j_coord_offsetA = DIM*jnrA;
176 j_coord_offsetB = DIM*jnrB;
177 j_coord_offsetC = DIM*jnrC;
178 j_coord_offsetD = DIM*jnrD;
179 j_coord_offsetE = DIM*jnrE;
180 j_coord_offsetF = DIM*jnrF;
181 j_coord_offsetG = DIM*jnrG;
182 j_coord_offsetH = DIM*jnrH;
184 /* load j atom coordinates */
185 gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
186 x+j_coord_offsetC,x+j_coord_offsetD,
187 x+j_coord_offsetE,x+j_coord_offsetF,
188 x+j_coord_offsetG,x+j_coord_offsetH,
191 /* Calculate displacement vector */
192 dx00 = _mm256_sub_ps(ix0,jx0);
193 dy00 = _mm256_sub_ps(iy0,jy0);
194 dz00 = _mm256_sub_ps(iz0,jz0);
196 /* Calculate squared distance and things based on it */
197 rsq00 = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
199 rinv00 = gmx_mm256_invsqrt_ps(rsq00);
201 rinvsq00 = _mm256_mul_ps(rinv00,rinv00);
203 /* Load parameters for j particles */
204 jq0 = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
205 charge+jnrC+0,charge+jnrD+0,
206 charge+jnrE+0,charge+jnrF+0,
207 charge+jnrG+0,charge+jnrH+0);
209 /**************************
210 * CALCULATE INTERACTIONS *
211 **************************/
213 /* Compute parameters for interactions between i and j atoms */
214 qq00 = _mm256_mul_ps(iq0,jq0);
216 /* REACTION-FIELD ELECTROSTATICS */
217 velec = _mm256_mul_ps(qq00,_mm256_sub_ps(_mm256_add_ps(rinv00,_mm256_mul_ps(krf,rsq00)),crf));
218 felec = _mm256_mul_ps(qq00,_mm256_sub_ps(_mm256_mul_ps(rinv00,rinvsq00),krf2));
220 /* Update potential sum for this i atom from the interaction with this j atom. */
221 velecsum = _mm256_add_ps(velecsum,velec);
225 /* Calculate temporary vectorial force */
226 tx = _mm256_mul_ps(fscal,dx00);
227 ty = _mm256_mul_ps(fscal,dy00);
228 tz = _mm256_mul_ps(fscal,dz00);
230 /* Update vectorial force */
231 fix0 = _mm256_add_ps(fix0,tx);
232 fiy0 = _mm256_add_ps(fiy0,ty);
233 fiz0 = _mm256_add_ps(fiz0,tz);
235 fjptrA = f+j_coord_offsetA;
236 fjptrB = f+j_coord_offsetB;
237 fjptrC = f+j_coord_offsetC;
238 fjptrD = f+j_coord_offsetD;
239 fjptrE = f+j_coord_offsetE;
240 fjptrF = f+j_coord_offsetF;
241 fjptrG = f+j_coord_offsetG;
242 fjptrH = f+j_coord_offsetH;
243 gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,tx,ty,tz);
245 /* Inner loop uses 32 flops */
251 /* Get j neighbor index, and coordinate index */
252 jnrlistA = jjnr[jidx];
253 jnrlistB = jjnr[jidx+1];
254 jnrlistC = jjnr[jidx+2];
255 jnrlistD = jjnr[jidx+3];
256 jnrlistE = jjnr[jidx+4];
257 jnrlistF = jjnr[jidx+5];
258 jnrlistG = jjnr[jidx+6];
259 jnrlistH = jjnr[jidx+7];
260 /* Sign of each element will be negative for non-real atoms.
261 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
262 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
264 dummy_mask = gmx_mm256_set_m128(gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx+4)),_mm_setzero_si128())),
265 gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128())));
267 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
268 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
269 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
270 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
271 jnrE = (jnrlistE>=0) ? jnrlistE : 0;
272 jnrF = (jnrlistF>=0) ? jnrlistF : 0;
273 jnrG = (jnrlistG>=0) ? jnrlistG : 0;
274 jnrH = (jnrlistH>=0) ? jnrlistH : 0;
275 j_coord_offsetA = DIM*jnrA;
276 j_coord_offsetB = DIM*jnrB;
277 j_coord_offsetC = DIM*jnrC;
278 j_coord_offsetD = DIM*jnrD;
279 j_coord_offsetE = DIM*jnrE;
280 j_coord_offsetF = DIM*jnrF;
281 j_coord_offsetG = DIM*jnrG;
282 j_coord_offsetH = DIM*jnrH;
284 /* load j atom coordinates */
285 gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
286 x+j_coord_offsetC,x+j_coord_offsetD,
287 x+j_coord_offsetE,x+j_coord_offsetF,
288 x+j_coord_offsetG,x+j_coord_offsetH,
291 /* Calculate displacement vector */
292 dx00 = _mm256_sub_ps(ix0,jx0);
293 dy00 = _mm256_sub_ps(iy0,jy0);
294 dz00 = _mm256_sub_ps(iz0,jz0);
296 /* Calculate squared distance and things based on it */
297 rsq00 = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
299 rinv00 = gmx_mm256_invsqrt_ps(rsq00);
301 rinvsq00 = _mm256_mul_ps(rinv00,rinv00);
303 /* Load parameters for j particles */
304 jq0 = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
305 charge+jnrC+0,charge+jnrD+0,
306 charge+jnrE+0,charge+jnrF+0,
307 charge+jnrG+0,charge+jnrH+0);
309 /**************************
310 * CALCULATE INTERACTIONS *
311 **************************/
313 /* Compute parameters for interactions between i and j atoms */
314 qq00 = _mm256_mul_ps(iq0,jq0);
316 /* REACTION-FIELD ELECTROSTATICS */
317 velec = _mm256_mul_ps(qq00,_mm256_sub_ps(_mm256_add_ps(rinv00,_mm256_mul_ps(krf,rsq00)),crf));
318 felec = _mm256_mul_ps(qq00,_mm256_sub_ps(_mm256_mul_ps(rinv00,rinvsq00),krf2));
320 /* Update potential sum for this i atom from the interaction with this j atom. */
321 velec = _mm256_andnot_ps(dummy_mask,velec);
322 velecsum = _mm256_add_ps(velecsum,velec);
326 fscal = _mm256_andnot_ps(dummy_mask,fscal);
328 /* Calculate temporary vectorial force */
329 tx = _mm256_mul_ps(fscal,dx00);
330 ty = _mm256_mul_ps(fscal,dy00);
331 tz = _mm256_mul_ps(fscal,dz00);
333 /* Update vectorial force */
334 fix0 = _mm256_add_ps(fix0,tx);
335 fiy0 = _mm256_add_ps(fiy0,ty);
336 fiz0 = _mm256_add_ps(fiz0,tz);
338 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
339 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
340 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
341 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
342 fjptrE = (jnrlistE>=0) ? f+j_coord_offsetE : scratch;
343 fjptrF = (jnrlistF>=0) ? f+j_coord_offsetF : scratch;
344 fjptrG = (jnrlistG>=0) ? f+j_coord_offsetG : scratch;
345 fjptrH = (jnrlistH>=0) ? f+j_coord_offsetH : scratch;
346 gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,tx,ty,tz);
348 /* Inner loop uses 32 flops */
351 /* End of innermost loop */
353 gmx_mm256_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
354 f+i_coord_offset,fshift+i_shift_offset);
357 /* Update potential energies */
358 gmx_mm256_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
360 /* Increment number of inner iterations */
361 inneriter += j_index_end - j_index_start;
363 /* Outer loop uses 8 flops */
366 /* Increment number of outer iterations */
369 /* Update outer/inner flops */
371 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VF,outeriter*8 + inneriter*32);
374 * Gromacs nonbonded kernel: nb_kernel_ElecRF_VdwNone_GeomP1P1_F_avx_256_single
375 * Electrostatics interaction: ReactionField
376 * VdW interaction: None
377 * Geometry: Particle-Particle
378 * Calculate force/pot: Force
381 nb_kernel_ElecRF_VdwNone_GeomP1P1_F_avx_256_single
382 (t_nblist * gmx_restrict nlist,
383 rvec * gmx_restrict xx,
384 rvec * gmx_restrict ff,
385 t_forcerec * gmx_restrict fr,
386 t_mdatoms * gmx_restrict mdatoms,
387 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
388 t_nrnb * gmx_restrict nrnb)
390 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
391 * just 0 for non-waters.
392 * Suffixes A,B,C,D,E,F,G,H refer to j loop unrolling done with AVX, e.g. for the eight different
393 * jnr indices corresponding to data put in the four positions in the SIMD register.
395 int i_shift_offset,i_coord_offset,outeriter,inneriter;
396 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
397 int jnrA,jnrB,jnrC,jnrD;
398 int jnrE,jnrF,jnrG,jnrH;
399 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
400 int jnrlistE,jnrlistF,jnrlistG,jnrlistH;
401 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
402 int j_coord_offsetE,j_coord_offsetF,j_coord_offsetG,j_coord_offsetH;
403 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
405 real *shiftvec,*fshift,*x,*f;
406 real *fjptrA,*fjptrB,*fjptrC,*fjptrD,*fjptrE,*fjptrF,*fjptrG,*fjptrH;
408 __m256 tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
409 real * vdwioffsetptr0;
410 __m256 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
411 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D,vdwjidx0E,vdwjidx0F,vdwjidx0G,vdwjidx0H;
412 __m256 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
413 __m256 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
414 __m256 velec,felec,velecsum,facel,crf,krf,krf2;
416 __m256 dummy_mask,cutoff_mask;
417 __m256 signbit = _mm256_castsi256_ps( _mm256_set1_epi32(0x80000000) );
418 __m256 one = _mm256_set1_ps(1.0);
419 __m256 two = _mm256_set1_ps(2.0);
425 jindex = nlist->jindex;
427 shiftidx = nlist->shift;
429 shiftvec = fr->shift_vec[0];
430 fshift = fr->fshift[0];
431 facel = _mm256_set1_ps(fr->epsfac);
432 charge = mdatoms->chargeA;
433 krf = _mm256_set1_ps(fr->ic->k_rf);
434 krf2 = _mm256_set1_ps(fr->ic->k_rf*2.0);
435 crf = _mm256_set1_ps(fr->ic->c_rf);
437 /* Avoid stupid compiler warnings */
438 jnrA = jnrB = jnrC = jnrD = jnrE = jnrF = jnrG = jnrH = 0;
451 for(iidx=0;iidx<4*DIM;iidx++)
456 /* Start outer loop over neighborlists */
457 for(iidx=0; iidx<nri; iidx++)
459 /* Load shift vector for this list */
460 i_shift_offset = DIM*shiftidx[iidx];
462 /* Load limits for loop over neighbors */
463 j_index_start = jindex[iidx];
464 j_index_end = jindex[iidx+1];
466 /* Get outer coordinate index */
468 i_coord_offset = DIM*inr;
470 /* Load i particle coords and add shift vector */
471 gmx_mm256_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
473 fix0 = _mm256_setzero_ps();
474 fiy0 = _mm256_setzero_ps();
475 fiz0 = _mm256_setzero_ps();
477 /* Load parameters for i particles */
478 iq0 = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+0]));
480 /* Start inner kernel loop */
481 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+7]>=0; jidx+=8)
484 /* Get j neighbor index, and coordinate index */
493 j_coord_offsetA = DIM*jnrA;
494 j_coord_offsetB = DIM*jnrB;
495 j_coord_offsetC = DIM*jnrC;
496 j_coord_offsetD = DIM*jnrD;
497 j_coord_offsetE = DIM*jnrE;
498 j_coord_offsetF = DIM*jnrF;
499 j_coord_offsetG = DIM*jnrG;
500 j_coord_offsetH = DIM*jnrH;
502 /* load j atom coordinates */
503 gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
504 x+j_coord_offsetC,x+j_coord_offsetD,
505 x+j_coord_offsetE,x+j_coord_offsetF,
506 x+j_coord_offsetG,x+j_coord_offsetH,
509 /* Calculate displacement vector */
510 dx00 = _mm256_sub_ps(ix0,jx0);
511 dy00 = _mm256_sub_ps(iy0,jy0);
512 dz00 = _mm256_sub_ps(iz0,jz0);
514 /* Calculate squared distance and things based on it */
515 rsq00 = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
517 rinv00 = gmx_mm256_invsqrt_ps(rsq00);
519 rinvsq00 = _mm256_mul_ps(rinv00,rinv00);
521 /* Load parameters for j particles */
522 jq0 = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
523 charge+jnrC+0,charge+jnrD+0,
524 charge+jnrE+0,charge+jnrF+0,
525 charge+jnrG+0,charge+jnrH+0);
527 /**************************
528 * CALCULATE INTERACTIONS *
529 **************************/
531 /* Compute parameters for interactions between i and j atoms */
532 qq00 = _mm256_mul_ps(iq0,jq0);
534 /* REACTION-FIELD ELECTROSTATICS */
535 felec = _mm256_mul_ps(qq00,_mm256_sub_ps(_mm256_mul_ps(rinv00,rinvsq00),krf2));
539 /* Calculate temporary vectorial force */
540 tx = _mm256_mul_ps(fscal,dx00);
541 ty = _mm256_mul_ps(fscal,dy00);
542 tz = _mm256_mul_ps(fscal,dz00);
544 /* Update vectorial force */
545 fix0 = _mm256_add_ps(fix0,tx);
546 fiy0 = _mm256_add_ps(fiy0,ty);
547 fiz0 = _mm256_add_ps(fiz0,tz);
549 fjptrA = f+j_coord_offsetA;
550 fjptrB = f+j_coord_offsetB;
551 fjptrC = f+j_coord_offsetC;
552 fjptrD = f+j_coord_offsetD;
553 fjptrE = f+j_coord_offsetE;
554 fjptrF = f+j_coord_offsetF;
555 fjptrG = f+j_coord_offsetG;
556 fjptrH = f+j_coord_offsetH;
557 gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,tx,ty,tz);
559 /* Inner loop uses 27 flops */
565 /* Get j neighbor index, and coordinate index */
566 jnrlistA = jjnr[jidx];
567 jnrlistB = jjnr[jidx+1];
568 jnrlistC = jjnr[jidx+2];
569 jnrlistD = jjnr[jidx+3];
570 jnrlistE = jjnr[jidx+4];
571 jnrlistF = jjnr[jidx+5];
572 jnrlistG = jjnr[jidx+6];
573 jnrlistH = jjnr[jidx+7];
574 /* Sign of each element will be negative for non-real atoms.
575 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
576 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
578 dummy_mask = gmx_mm256_set_m128(gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx+4)),_mm_setzero_si128())),
579 gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128())));
581 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
582 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
583 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
584 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
585 jnrE = (jnrlistE>=0) ? jnrlistE : 0;
586 jnrF = (jnrlistF>=0) ? jnrlistF : 0;
587 jnrG = (jnrlistG>=0) ? jnrlistG : 0;
588 jnrH = (jnrlistH>=0) ? jnrlistH : 0;
589 j_coord_offsetA = DIM*jnrA;
590 j_coord_offsetB = DIM*jnrB;
591 j_coord_offsetC = DIM*jnrC;
592 j_coord_offsetD = DIM*jnrD;
593 j_coord_offsetE = DIM*jnrE;
594 j_coord_offsetF = DIM*jnrF;
595 j_coord_offsetG = DIM*jnrG;
596 j_coord_offsetH = DIM*jnrH;
598 /* load j atom coordinates */
599 gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
600 x+j_coord_offsetC,x+j_coord_offsetD,
601 x+j_coord_offsetE,x+j_coord_offsetF,
602 x+j_coord_offsetG,x+j_coord_offsetH,
605 /* Calculate displacement vector */
606 dx00 = _mm256_sub_ps(ix0,jx0);
607 dy00 = _mm256_sub_ps(iy0,jy0);
608 dz00 = _mm256_sub_ps(iz0,jz0);
610 /* Calculate squared distance and things based on it */
611 rsq00 = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
613 rinv00 = gmx_mm256_invsqrt_ps(rsq00);
615 rinvsq00 = _mm256_mul_ps(rinv00,rinv00);
617 /* Load parameters for j particles */
618 jq0 = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
619 charge+jnrC+0,charge+jnrD+0,
620 charge+jnrE+0,charge+jnrF+0,
621 charge+jnrG+0,charge+jnrH+0);
623 /**************************
624 * CALCULATE INTERACTIONS *
625 **************************/
627 /* Compute parameters for interactions between i and j atoms */
628 qq00 = _mm256_mul_ps(iq0,jq0);
630 /* REACTION-FIELD ELECTROSTATICS */
631 felec = _mm256_mul_ps(qq00,_mm256_sub_ps(_mm256_mul_ps(rinv00,rinvsq00),krf2));
635 fscal = _mm256_andnot_ps(dummy_mask,fscal);
637 /* Calculate temporary vectorial force */
638 tx = _mm256_mul_ps(fscal,dx00);
639 ty = _mm256_mul_ps(fscal,dy00);
640 tz = _mm256_mul_ps(fscal,dz00);
642 /* Update vectorial force */
643 fix0 = _mm256_add_ps(fix0,tx);
644 fiy0 = _mm256_add_ps(fiy0,ty);
645 fiz0 = _mm256_add_ps(fiz0,tz);
647 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
648 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
649 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
650 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
651 fjptrE = (jnrlistE>=0) ? f+j_coord_offsetE : scratch;
652 fjptrF = (jnrlistF>=0) ? f+j_coord_offsetF : scratch;
653 fjptrG = (jnrlistG>=0) ? f+j_coord_offsetG : scratch;
654 fjptrH = (jnrlistH>=0) ? f+j_coord_offsetH : scratch;
655 gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,tx,ty,tz);
657 /* Inner loop uses 27 flops */
660 /* End of innermost loop */
662 gmx_mm256_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
663 f+i_coord_offset,fshift+i_shift_offset);
665 /* Increment number of inner iterations */
666 inneriter += j_index_end - j_index_start;
668 /* Outer loop uses 7 flops */
671 /* Increment number of outer iterations */
674 /* Update outer/inner flops */
676 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_F,outeriter*7 + inneriter*27);