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_ElecGB_VdwNone_GeomP1P1_VF_avx_256_single
38 * Electrostatics interaction: GeneralizedBorn
39 * VdW interaction: None
40 * Geometry: Particle-Particle
41 * Calculate force/pot: PotentialAndForce
44 nb_kernel_ElecGB_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 gbitab_lo,gbitab_hi;
81 __m256 vgb,fgb,vgbsum,dvdasum,gbscale,gbtabscale,isaprod,gbqqfactor,gbinvepsdiff,gbeps,dvdatmp;
82 __m256 minushalf = _mm256_set1_ps(-0.5);
83 real *invsqrta,*dvda,*gbtab;
85 __m128i vfitab_lo,vfitab_hi;
86 __m128i ifour = _mm_set1_epi32(4);
87 __m256 rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF;
89 __m256 dummy_mask,cutoff_mask;
90 __m256 signbit = _mm256_castsi256_ps( _mm256_set1_epi32(0x80000000) );
91 __m256 one = _mm256_set1_ps(1.0);
92 __m256 two = _mm256_set1_ps(2.0);
98 jindex = nlist->jindex;
100 shiftidx = nlist->shift;
102 shiftvec = fr->shift_vec[0];
103 fshift = fr->fshift[0];
104 facel = _mm256_set1_ps(fr->epsfac);
105 charge = mdatoms->chargeA;
107 invsqrta = fr->invsqrta;
109 gbtabscale = _mm256_set1_ps(fr->gbtab.scale);
110 gbtab = fr->gbtab.data;
111 gbinvepsdiff = _mm256_set1_ps((1.0/fr->epsilon_r) - (1.0/fr->gb_epsilon_solvent));
113 /* Avoid stupid compiler warnings */
114 jnrA = jnrB = jnrC = jnrD = jnrE = jnrF = jnrG = jnrH = 0;
127 for(iidx=0;iidx<4*DIM;iidx++)
132 /* Start outer loop over neighborlists */
133 for(iidx=0; iidx<nri; iidx++)
135 /* Load shift vector for this list */
136 i_shift_offset = DIM*shiftidx[iidx];
138 /* Load limits for loop over neighbors */
139 j_index_start = jindex[iidx];
140 j_index_end = jindex[iidx+1];
142 /* Get outer coordinate index */
144 i_coord_offset = DIM*inr;
146 /* Load i particle coords and add shift vector */
147 gmx_mm256_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
149 fix0 = _mm256_setzero_ps();
150 fiy0 = _mm256_setzero_ps();
151 fiz0 = _mm256_setzero_ps();
153 /* Load parameters for i particles */
154 iq0 = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+0]));
155 isai0 = _mm256_set1_ps(invsqrta[inr+0]);
157 /* Reset potential sums */
158 velecsum = _mm256_setzero_ps();
159 vgbsum = _mm256_setzero_ps();
160 dvdasum = _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 /* Load parameters for j particles */
202 jq0 = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
203 charge+jnrC+0,charge+jnrD+0,
204 charge+jnrE+0,charge+jnrF+0,
205 charge+jnrG+0,charge+jnrH+0);
206 isaj0 = gmx_mm256_load_8real_swizzle_ps(invsqrta+jnrA+0,invsqrta+jnrB+0,
207 invsqrta+jnrC+0,invsqrta+jnrD+0,
208 invsqrta+jnrE+0,invsqrta+jnrF+0,
209 invsqrta+jnrG+0,invsqrta+jnrH+0);
211 /**************************
212 * CALCULATE INTERACTIONS *
213 **************************/
215 r00 = _mm256_mul_ps(rsq00,rinv00);
217 /* Compute parameters for interactions between i and j atoms */
218 qq00 = _mm256_mul_ps(iq0,jq0);
220 /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
221 isaprod = _mm256_mul_ps(isai0,isaj0);
222 gbqqfactor = _mm256_xor_ps(signbit,_mm256_mul_ps(qq00,_mm256_mul_ps(isaprod,gbinvepsdiff)));
223 gbscale = _mm256_mul_ps(isaprod,gbtabscale);
225 /* Calculate generalized born table index - this is a separate table from the normal one,
226 * but we use the same procedure by multiplying r with scale and truncating to integer.
228 rt = _mm256_mul_ps(r00,gbscale);
229 gbitab = _mm256_cvttps_epi32(rt);
230 gbeps = _mm256_sub_ps(rt,_mm256_round_ps(rt, _MM_FROUND_FLOOR));
231 /* AVX1 does not support 256-bit integer operations, so now we go to 128-bit mode... */
232 gbitab_lo = _mm256_extractf128_si256(gbitab,0x0);
233 gbitab_hi = _mm256_extractf128_si256(gbitab,0x1);
234 gbitab_lo = _mm_slli_epi32(gbitab_lo,2);
235 gbitab_hi = _mm_slli_epi32(gbitab_hi,2);
236 Y = gmx_mm256_set_m128(_mm_load_ps(gbtab + _mm_extract_epi32(gbitab_hi,0)),
237 _mm_load_ps(gbtab + _mm_extract_epi32(gbitab_lo,0)));
238 F = gmx_mm256_set_m128(_mm_load_ps(gbtab + _mm_extract_epi32(gbitab_hi,1)),
239 _mm_load_ps(gbtab + _mm_extract_epi32(gbitab_lo,1)));
240 G = gmx_mm256_set_m128(_mm_load_ps(gbtab + _mm_extract_epi32(gbitab_hi,2)),
241 _mm_load_ps(gbtab + _mm_extract_epi32(gbitab_lo,2)));
242 H = gmx_mm256_set_m128(_mm_load_ps(gbtab + _mm_extract_epi32(gbitab_hi,3)),
243 _mm_load_ps(gbtab + _mm_extract_epi32(gbitab_lo,3)));
244 GMX_MM256_HALFTRANSPOSE4_PS(Y,F,G,H);
245 Heps = _mm256_mul_ps(gbeps,H);
246 Fp = _mm256_add_ps(F,_mm256_mul_ps(gbeps,_mm256_add_ps(G,Heps)));
247 VV = _mm256_add_ps(Y,_mm256_mul_ps(gbeps,Fp));
248 vgb = _mm256_mul_ps(gbqqfactor,VV);
250 FF = _mm256_add_ps(Fp,_mm256_mul_ps(gbeps,_mm256_add_ps(G,_mm256_add_ps(Heps,Heps))));
251 fgb = _mm256_mul_ps(gbqqfactor,_mm256_mul_ps(FF,gbscale));
252 dvdatmp = _mm256_mul_ps(minushalf,_mm256_add_ps(vgb,_mm256_mul_ps(fgb,r00)));
253 dvdasum = _mm256_add_ps(dvdasum,dvdatmp);
262 gmx_mm256_increment_8real_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,
263 _mm256_mul_ps(dvdatmp,_mm256_mul_ps(isaj0,isaj0)));
264 velec = _mm256_mul_ps(qq00,rinv00);
265 felec = _mm256_mul_ps(_mm256_sub_ps(_mm256_mul_ps(velec,rinv00),fgb),rinv00);
267 /* Update potential sum for this i atom from the interaction with this j atom. */
268 velecsum = _mm256_add_ps(velecsum,velec);
269 vgbsum = _mm256_add_ps(vgbsum,vgb);
273 /* Calculate temporary vectorial force */
274 tx = _mm256_mul_ps(fscal,dx00);
275 ty = _mm256_mul_ps(fscal,dy00);
276 tz = _mm256_mul_ps(fscal,dz00);
278 /* Update vectorial force */
279 fix0 = _mm256_add_ps(fix0,tx);
280 fiy0 = _mm256_add_ps(fiy0,ty);
281 fiz0 = _mm256_add_ps(fiz0,tz);
283 fjptrA = f+j_coord_offsetA;
284 fjptrB = f+j_coord_offsetB;
285 fjptrC = f+j_coord_offsetC;
286 fjptrD = f+j_coord_offsetD;
287 fjptrE = f+j_coord_offsetE;
288 fjptrF = f+j_coord_offsetF;
289 fjptrG = f+j_coord_offsetG;
290 fjptrH = f+j_coord_offsetH;
291 gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,tx,ty,tz);
293 /* Inner loop uses 57 flops */
299 /* Get j neighbor index, and coordinate index */
300 jnrlistA = jjnr[jidx];
301 jnrlistB = jjnr[jidx+1];
302 jnrlistC = jjnr[jidx+2];
303 jnrlistD = jjnr[jidx+3];
304 jnrlistE = jjnr[jidx+4];
305 jnrlistF = jjnr[jidx+5];
306 jnrlistG = jjnr[jidx+6];
307 jnrlistH = jjnr[jidx+7];
308 /* Sign of each element will be negative for non-real atoms.
309 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
310 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
312 dummy_mask = gmx_mm256_set_m128(gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx+4)),_mm_setzero_si128())),
313 gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128())));
315 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
316 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
317 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
318 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
319 jnrE = (jnrlistE>=0) ? jnrlistE : 0;
320 jnrF = (jnrlistF>=0) ? jnrlistF : 0;
321 jnrG = (jnrlistG>=0) ? jnrlistG : 0;
322 jnrH = (jnrlistH>=0) ? jnrlistH : 0;
323 j_coord_offsetA = DIM*jnrA;
324 j_coord_offsetB = DIM*jnrB;
325 j_coord_offsetC = DIM*jnrC;
326 j_coord_offsetD = DIM*jnrD;
327 j_coord_offsetE = DIM*jnrE;
328 j_coord_offsetF = DIM*jnrF;
329 j_coord_offsetG = DIM*jnrG;
330 j_coord_offsetH = DIM*jnrH;
332 /* load j atom coordinates */
333 gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
334 x+j_coord_offsetC,x+j_coord_offsetD,
335 x+j_coord_offsetE,x+j_coord_offsetF,
336 x+j_coord_offsetG,x+j_coord_offsetH,
339 /* Calculate displacement vector */
340 dx00 = _mm256_sub_ps(ix0,jx0);
341 dy00 = _mm256_sub_ps(iy0,jy0);
342 dz00 = _mm256_sub_ps(iz0,jz0);
344 /* Calculate squared distance and things based on it */
345 rsq00 = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
347 rinv00 = gmx_mm256_invsqrt_ps(rsq00);
349 /* Load parameters for j particles */
350 jq0 = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
351 charge+jnrC+0,charge+jnrD+0,
352 charge+jnrE+0,charge+jnrF+0,
353 charge+jnrG+0,charge+jnrH+0);
354 isaj0 = gmx_mm256_load_8real_swizzle_ps(invsqrta+jnrA+0,invsqrta+jnrB+0,
355 invsqrta+jnrC+0,invsqrta+jnrD+0,
356 invsqrta+jnrE+0,invsqrta+jnrF+0,
357 invsqrta+jnrG+0,invsqrta+jnrH+0);
359 /**************************
360 * CALCULATE INTERACTIONS *
361 **************************/
363 r00 = _mm256_mul_ps(rsq00,rinv00);
364 r00 = _mm256_andnot_ps(dummy_mask,r00);
366 /* Compute parameters for interactions between i and j atoms */
367 qq00 = _mm256_mul_ps(iq0,jq0);
369 /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
370 isaprod = _mm256_mul_ps(isai0,isaj0);
371 gbqqfactor = _mm256_xor_ps(signbit,_mm256_mul_ps(qq00,_mm256_mul_ps(isaprod,gbinvepsdiff)));
372 gbscale = _mm256_mul_ps(isaprod,gbtabscale);
374 /* Calculate generalized born table index - this is a separate table from the normal one,
375 * but we use the same procedure by multiplying r with scale and truncating to integer.
377 rt = _mm256_mul_ps(r00,gbscale);
378 gbitab = _mm256_cvttps_epi32(rt);
379 gbeps = _mm256_sub_ps(rt,_mm256_round_ps(rt, _MM_FROUND_FLOOR));
380 /* AVX1 does not support 256-bit integer operations, so now we go to 128-bit mode... */
381 gbitab_lo = _mm256_extractf128_si256(gbitab,0x0);
382 gbitab_hi = _mm256_extractf128_si256(gbitab,0x1);
383 gbitab_lo = _mm_slli_epi32(gbitab_lo,2);
384 gbitab_hi = _mm_slli_epi32(gbitab_hi,2);
385 Y = gmx_mm256_set_m128(_mm_load_ps(gbtab + _mm_extract_epi32(gbitab_hi,0)),
386 _mm_load_ps(gbtab + _mm_extract_epi32(gbitab_lo,0)));
387 F = gmx_mm256_set_m128(_mm_load_ps(gbtab + _mm_extract_epi32(gbitab_hi,1)),
388 _mm_load_ps(gbtab + _mm_extract_epi32(gbitab_lo,1)));
389 G = gmx_mm256_set_m128(_mm_load_ps(gbtab + _mm_extract_epi32(gbitab_hi,2)),
390 _mm_load_ps(gbtab + _mm_extract_epi32(gbitab_lo,2)));
391 H = gmx_mm256_set_m128(_mm_load_ps(gbtab + _mm_extract_epi32(gbitab_hi,3)),
392 _mm_load_ps(gbtab + _mm_extract_epi32(gbitab_lo,3)));
393 GMX_MM256_HALFTRANSPOSE4_PS(Y,F,G,H);
394 Heps = _mm256_mul_ps(gbeps,H);
395 Fp = _mm256_add_ps(F,_mm256_mul_ps(gbeps,_mm256_add_ps(G,Heps)));
396 VV = _mm256_add_ps(Y,_mm256_mul_ps(gbeps,Fp));
397 vgb = _mm256_mul_ps(gbqqfactor,VV);
399 FF = _mm256_add_ps(Fp,_mm256_mul_ps(gbeps,_mm256_add_ps(G,_mm256_add_ps(Heps,Heps))));
400 fgb = _mm256_mul_ps(gbqqfactor,_mm256_mul_ps(FF,gbscale));
401 dvdatmp = _mm256_mul_ps(minushalf,_mm256_add_ps(vgb,_mm256_mul_ps(fgb,r00)));
402 dvdatmp = _mm256_andnot_ps(dummy_mask,dvdatmp);
403 dvdasum = _mm256_add_ps(dvdasum,dvdatmp);
404 /* The pointers to scratch make sure that this code with compilers that take gmx_restrict seriously (e.g. icc 13) really can't screw things up. */
405 fjptrA = (jnrlistA>=0) ? dvda+jnrA : scratch;
406 fjptrB = (jnrlistB>=0) ? dvda+jnrB : scratch;
407 fjptrC = (jnrlistC>=0) ? dvda+jnrC : scratch;
408 fjptrD = (jnrlistD>=0) ? dvda+jnrD : scratch;
409 fjptrE = (jnrlistE>=0) ? dvda+jnrE : scratch;
410 fjptrF = (jnrlistF>=0) ? dvda+jnrF : scratch;
411 fjptrG = (jnrlistG>=0) ? dvda+jnrG : scratch;
412 fjptrH = (jnrlistH>=0) ? dvda+jnrH : scratch;
413 gmx_mm256_increment_8real_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,
414 _mm256_mul_ps(dvdatmp,_mm256_mul_ps(isaj0,isaj0)));
415 velec = _mm256_mul_ps(qq00,rinv00);
416 felec = _mm256_mul_ps(_mm256_sub_ps(_mm256_mul_ps(velec,rinv00),fgb),rinv00);
418 /* Update potential sum for this i atom from the interaction with this j atom. */
419 velec = _mm256_andnot_ps(dummy_mask,velec);
420 velecsum = _mm256_add_ps(velecsum,velec);
421 vgb = _mm256_andnot_ps(dummy_mask,vgb);
422 vgbsum = _mm256_add_ps(vgbsum,vgb);
426 fscal = _mm256_andnot_ps(dummy_mask,fscal);
428 /* Calculate temporary vectorial force */
429 tx = _mm256_mul_ps(fscal,dx00);
430 ty = _mm256_mul_ps(fscal,dy00);
431 tz = _mm256_mul_ps(fscal,dz00);
433 /* Update vectorial force */
434 fix0 = _mm256_add_ps(fix0,tx);
435 fiy0 = _mm256_add_ps(fiy0,ty);
436 fiz0 = _mm256_add_ps(fiz0,tz);
438 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
439 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
440 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
441 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
442 fjptrE = (jnrlistE>=0) ? f+j_coord_offsetE : scratch;
443 fjptrF = (jnrlistF>=0) ? f+j_coord_offsetF : scratch;
444 fjptrG = (jnrlistG>=0) ? f+j_coord_offsetG : scratch;
445 fjptrH = (jnrlistH>=0) ? f+j_coord_offsetH : scratch;
446 gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,tx,ty,tz);
448 /* Inner loop uses 58 flops */
451 /* End of innermost loop */
453 gmx_mm256_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
454 f+i_coord_offset,fshift+i_shift_offset);
457 /* Update potential energies */
458 gmx_mm256_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
459 gmx_mm256_update_1pot_ps(vgbsum,kernel_data->energygrp_polarization+ggid);
460 dvdasum = _mm256_mul_ps(dvdasum, _mm256_mul_ps(isai0,isai0));
461 gmx_mm256_update_1pot_ps(dvdasum,dvda+inr);
463 /* Increment number of inner iterations */
464 inneriter += j_index_end - j_index_start;
466 /* Outer loop uses 9 flops */
469 /* Increment number of outer iterations */
472 /* Update outer/inner flops */
474 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VF,outeriter*9 + inneriter*58);
477 * Gromacs nonbonded kernel: nb_kernel_ElecGB_VdwNone_GeomP1P1_F_avx_256_single
478 * Electrostatics interaction: GeneralizedBorn
479 * VdW interaction: None
480 * Geometry: Particle-Particle
481 * Calculate force/pot: Force
484 nb_kernel_ElecGB_VdwNone_GeomP1P1_F_avx_256_single
485 (t_nblist * gmx_restrict nlist,
486 rvec * gmx_restrict xx,
487 rvec * gmx_restrict ff,
488 t_forcerec * gmx_restrict fr,
489 t_mdatoms * gmx_restrict mdatoms,
490 nb_kernel_data_t * gmx_restrict kernel_data,
491 t_nrnb * gmx_restrict nrnb)
493 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
494 * just 0 for non-waters.
495 * Suffixes A,B,C,D,E,F,G,H refer to j loop unrolling done with AVX, e.g. for the eight different
496 * jnr indices corresponding to data put in the four positions in the SIMD register.
498 int i_shift_offset,i_coord_offset,outeriter,inneriter;
499 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
500 int jnrA,jnrB,jnrC,jnrD;
501 int jnrE,jnrF,jnrG,jnrH;
502 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
503 int jnrlistE,jnrlistF,jnrlistG,jnrlistH;
504 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
505 int j_coord_offsetE,j_coord_offsetF,j_coord_offsetG,j_coord_offsetH;
506 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
508 real *shiftvec,*fshift,*x,*f;
509 real *fjptrA,*fjptrB,*fjptrC,*fjptrD,*fjptrE,*fjptrF,*fjptrG,*fjptrH;
511 __m256 tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
512 real * vdwioffsetptr0;
513 __m256 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
514 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D,vdwjidx0E,vdwjidx0F,vdwjidx0G,vdwjidx0H;
515 __m256 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
516 __m256 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
517 __m256 velec,felec,velecsum,facel,crf,krf,krf2;
520 __m128i gbitab_lo,gbitab_hi;
521 __m256 vgb,fgb,vgbsum,dvdasum,gbscale,gbtabscale,isaprod,gbqqfactor,gbinvepsdiff,gbeps,dvdatmp;
522 __m256 minushalf = _mm256_set1_ps(-0.5);
523 real *invsqrta,*dvda,*gbtab;
525 __m128i vfitab_lo,vfitab_hi;
526 __m128i ifour = _mm_set1_epi32(4);
527 __m256 rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF;
529 __m256 dummy_mask,cutoff_mask;
530 __m256 signbit = _mm256_castsi256_ps( _mm256_set1_epi32(0x80000000) );
531 __m256 one = _mm256_set1_ps(1.0);
532 __m256 two = _mm256_set1_ps(2.0);
538 jindex = nlist->jindex;
540 shiftidx = nlist->shift;
542 shiftvec = fr->shift_vec[0];
543 fshift = fr->fshift[0];
544 facel = _mm256_set1_ps(fr->epsfac);
545 charge = mdatoms->chargeA;
547 invsqrta = fr->invsqrta;
549 gbtabscale = _mm256_set1_ps(fr->gbtab.scale);
550 gbtab = fr->gbtab.data;
551 gbinvepsdiff = _mm256_set1_ps((1.0/fr->epsilon_r) - (1.0/fr->gb_epsilon_solvent));
553 /* Avoid stupid compiler warnings */
554 jnrA = jnrB = jnrC = jnrD = jnrE = jnrF = jnrG = jnrH = 0;
567 for(iidx=0;iidx<4*DIM;iidx++)
572 /* Start outer loop over neighborlists */
573 for(iidx=0; iidx<nri; iidx++)
575 /* Load shift vector for this list */
576 i_shift_offset = DIM*shiftidx[iidx];
578 /* Load limits for loop over neighbors */
579 j_index_start = jindex[iidx];
580 j_index_end = jindex[iidx+1];
582 /* Get outer coordinate index */
584 i_coord_offset = DIM*inr;
586 /* Load i particle coords and add shift vector */
587 gmx_mm256_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
589 fix0 = _mm256_setzero_ps();
590 fiy0 = _mm256_setzero_ps();
591 fiz0 = _mm256_setzero_ps();
593 /* Load parameters for i particles */
594 iq0 = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+0]));
595 isai0 = _mm256_set1_ps(invsqrta[inr+0]);
597 dvdasum = _mm256_setzero_ps();
599 /* Start inner kernel loop */
600 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+7]>=0; jidx+=8)
603 /* Get j neighbor index, and coordinate index */
612 j_coord_offsetA = DIM*jnrA;
613 j_coord_offsetB = DIM*jnrB;
614 j_coord_offsetC = DIM*jnrC;
615 j_coord_offsetD = DIM*jnrD;
616 j_coord_offsetE = DIM*jnrE;
617 j_coord_offsetF = DIM*jnrF;
618 j_coord_offsetG = DIM*jnrG;
619 j_coord_offsetH = DIM*jnrH;
621 /* load j atom coordinates */
622 gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
623 x+j_coord_offsetC,x+j_coord_offsetD,
624 x+j_coord_offsetE,x+j_coord_offsetF,
625 x+j_coord_offsetG,x+j_coord_offsetH,
628 /* Calculate displacement vector */
629 dx00 = _mm256_sub_ps(ix0,jx0);
630 dy00 = _mm256_sub_ps(iy0,jy0);
631 dz00 = _mm256_sub_ps(iz0,jz0);
633 /* Calculate squared distance and things based on it */
634 rsq00 = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
636 rinv00 = gmx_mm256_invsqrt_ps(rsq00);
638 /* Load parameters for j particles */
639 jq0 = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
640 charge+jnrC+0,charge+jnrD+0,
641 charge+jnrE+0,charge+jnrF+0,
642 charge+jnrG+0,charge+jnrH+0);
643 isaj0 = gmx_mm256_load_8real_swizzle_ps(invsqrta+jnrA+0,invsqrta+jnrB+0,
644 invsqrta+jnrC+0,invsqrta+jnrD+0,
645 invsqrta+jnrE+0,invsqrta+jnrF+0,
646 invsqrta+jnrG+0,invsqrta+jnrH+0);
648 /**************************
649 * CALCULATE INTERACTIONS *
650 **************************/
652 r00 = _mm256_mul_ps(rsq00,rinv00);
654 /* Compute parameters for interactions between i and j atoms */
655 qq00 = _mm256_mul_ps(iq0,jq0);
657 /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
658 isaprod = _mm256_mul_ps(isai0,isaj0);
659 gbqqfactor = _mm256_xor_ps(signbit,_mm256_mul_ps(qq00,_mm256_mul_ps(isaprod,gbinvepsdiff)));
660 gbscale = _mm256_mul_ps(isaprod,gbtabscale);
662 /* Calculate generalized born table index - this is a separate table from the normal one,
663 * but we use the same procedure by multiplying r with scale and truncating to integer.
665 rt = _mm256_mul_ps(r00,gbscale);
666 gbitab = _mm256_cvttps_epi32(rt);
667 gbeps = _mm256_sub_ps(rt,_mm256_round_ps(rt, _MM_FROUND_FLOOR));
668 /* AVX1 does not support 256-bit integer operations, so now we go to 128-bit mode... */
669 gbitab_lo = _mm256_extractf128_si256(gbitab,0x0);
670 gbitab_hi = _mm256_extractf128_si256(gbitab,0x1);
671 gbitab_lo = _mm_slli_epi32(gbitab_lo,2);
672 gbitab_hi = _mm_slli_epi32(gbitab_hi,2);
673 Y = gmx_mm256_set_m128(_mm_load_ps(gbtab + _mm_extract_epi32(gbitab_hi,0)),
674 _mm_load_ps(gbtab + _mm_extract_epi32(gbitab_lo,0)));
675 F = gmx_mm256_set_m128(_mm_load_ps(gbtab + _mm_extract_epi32(gbitab_hi,1)),
676 _mm_load_ps(gbtab + _mm_extract_epi32(gbitab_lo,1)));
677 G = gmx_mm256_set_m128(_mm_load_ps(gbtab + _mm_extract_epi32(gbitab_hi,2)),
678 _mm_load_ps(gbtab + _mm_extract_epi32(gbitab_lo,2)));
679 H = gmx_mm256_set_m128(_mm_load_ps(gbtab + _mm_extract_epi32(gbitab_hi,3)),
680 _mm_load_ps(gbtab + _mm_extract_epi32(gbitab_lo,3)));
681 GMX_MM256_HALFTRANSPOSE4_PS(Y,F,G,H);
682 Heps = _mm256_mul_ps(gbeps,H);
683 Fp = _mm256_add_ps(F,_mm256_mul_ps(gbeps,_mm256_add_ps(G,Heps)));
684 VV = _mm256_add_ps(Y,_mm256_mul_ps(gbeps,Fp));
685 vgb = _mm256_mul_ps(gbqqfactor,VV);
687 FF = _mm256_add_ps(Fp,_mm256_mul_ps(gbeps,_mm256_add_ps(G,_mm256_add_ps(Heps,Heps))));
688 fgb = _mm256_mul_ps(gbqqfactor,_mm256_mul_ps(FF,gbscale));
689 dvdatmp = _mm256_mul_ps(minushalf,_mm256_add_ps(vgb,_mm256_mul_ps(fgb,r00)));
690 dvdasum = _mm256_add_ps(dvdasum,dvdatmp);
699 gmx_mm256_increment_8real_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,
700 _mm256_mul_ps(dvdatmp,_mm256_mul_ps(isaj0,isaj0)));
701 velec = _mm256_mul_ps(qq00,rinv00);
702 felec = _mm256_mul_ps(_mm256_sub_ps(_mm256_mul_ps(velec,rinv00),fgb),rinv00);
706 /* Calculate temporary vectorial force */
707 tx = _mm256_mul_ps(fscal,dx00);
708 ty = _mm256_mul_ps(fscal,dy00);
709 tz = _mm256_mul_ps(fscal,dz00);
711 /* Update vectorial force */
712 fix0 = _mm256_add_ps(fix0,tx);
713 fiy0 = _mm256_add_ps(fiy0,ty);
714 fiz0 = _mm256_add_ps(fiz0,tz);
716 fjptrA = f+j_coord_offsetA;
717 fjptrB = f+j_coord_offsetB;
718 fjptrC = f+j_coord_offsetC;
719 fjptrD = f+j_coord_offsetD;
720 fjptrE = f+j_coord_offsetE;
721 fjptrF = f+j_coord_offsetF;
722 fjptrG = f+j_coord_offsetG;
723 fjptrH = f+j_coord_offsetH;
724 gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,tx,ty,tz);
726 /* Inner loop uses 55 flops */
732 /* Get j neighbor index, and coordinate index */
733 jnrlistA = jjnr[jidx];
734 jnrlistB = jjnr[jidx+1];
735 jnrlistC = jjnr[jidx+2];
736 jnrlistD = jjnr[jidx+3];
737 jnrlistE = jjnr[jidx+4];
738 jnrlistF = jjnr[jidx+5];
739 jnrlistG = jjnr[jidx+6];
740 jnrlistH = jjnr[jidx+7];
741 /* Sign of each element will be negative for non-real atoms.
742 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
743 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
745 dummy_mask = gmx_mm256_set_m128(gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx+4)),_mm_setzero_si128())),
746 gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128())));
748 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
749 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
750 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
751 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
752 jnrE = (jnrlistE>=0) ? jnrlistE : 0;
753 jnrF = (jnrlistF>=0) ? jnrlistF : 0;
754 jnrG = (jnrlistG>=0) ? jnrlistG : 0;
755 jnrH = (jnrlistH>=0) ? jnrlistH : 0;
756 j_coord_offsetA = DIM*jnrA;
757 j_coord_offsetB = DIM*jnrB;
758 j_coord_offsetC = DIM*jnrC;
759 j_coord_offsetD = DIM*jnrD;
760 j_coord_offsetE = DIM*jnrE;
761 j_coord_offsetF = DIM*jnrF;
762 j_coord_offsetG = DIM*jnrG;
763 j_coord_offsetH = DIM*jnrH;
765 /* load j atom coordinates */
766 gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
767 x+j_coord_offsetC,x+j_coord_offsetD,
768 x+j_coord_offsetE,x+j_coord_offsetF,
769 x+j_coord_offsetG,x+j_coord_offsetH,
772 /* Calculate displacement vector */
773 dx00 = _mm256_sub_ps(ix0,jx0);
774 dy00 = _mm256_sub_ps(iy0,jy0);
775 dz00 = _mm256_sub_ps(iz0,jz0);
777 /* Calculate squared distance and things based on it */
778 rsq00 = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
780 rinv00 = gmx_mm256_invsqrt_ps(rsq00);
782 /* Load parameters for j particles */
783 jq0 = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
784 charge+jnrC+0,charge+jnrD+0,
785 charge+jnrE+0,charge+jnrF+0,
786 charge+jnrG+0,charge+jnrH+0);
787 isaj0 = gmx_mm256_load_8real_swizzle_ps(invsqrta+jnrA+0,invsqrta+jnrB+0,
788 invsqrta+jnrC+0,invsqrta+jnrD+0,
789 invsqrta+jnrE+0,invsqrta+jnrF+0,
790 invsqrta+jnrG+0,invsqrta+jnrH+0);
792 /**************************
793 * CALCULATE INTERACTIONS *
794 **************************/
796 r00 = _mm256_mul_ps(rsq00,rinv00);
797 r00 = _mm256_andnot_ps(dummy_mask,r00);
799 /* Compute parameters for interactions between i and j atoms */
800 qq00 = _mm256_mul_ps(iq0,jq0);
802 /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
803 isaprod = _mm256_mul_ps(isai0,isaj0);
804 gbqqfactor = _mm256_xor_ps(signbit,_mm256_mul_ps(qq00,_mm256_mul_ps(isaprod,gbinvepsdiff)));
805 gbscale = _mm256_mul_ps(isaprod,gbtabscale);
807 /* Calculate generalized born table index - this is a separate table from the normal one,
808 * but we use the same procedure by multiplying r with scale and truncating to integer.
810 rt = _mm256_mul_ps(r00,gbscale);
811 gbitab = _mm256_cvttps_epi32(rt);
812 gbeps = _mm256_sub_ps(rt,_mm256_round_ps(rt, _MM_FROUND_FLOOR));
813 /* AVX1 does not support 256-bit integer operations, so now we go to 128-bit mode... */
814 gbitab_lo = _mm256_extractf128_si256(gbitab,0x0);
815 gbitab_hi = _mm256_extractf128_si256(gbitab,0x1);
816 gbitab_lo = _mm_slli_epi32(gbitab_lo,2);
817 gbitab_hi = _mm_slli_epi32(gbitab_hi,2);
818 Y = gmx_mm256_set_m128(_mm_load_ps(gbtab + _mm_extract_epi32(gbitab_hi,0)),
819 _mm_load_ps(gbtab + _mm_extract_epi32(gbitab_lo,0)));
820 F = gmx_mm256_set_m128(_mm_load_ps(gbtab + _mm_extract_epi32(gbitab_hi,1)),
821 _mm_load_ps(gbtab + _mm_extract_epi32(gbitab_lo,1)));
822 G = gmx_mm256_set_m128(_mm_load_ps(gbtab + _mm_extract_epi32(gbitab_hi,2)),
823 _mm_load_ps(gbtab + _mm_extract_epi32(gbitab_lo,2)));
824 H = gmx_mm256_set_m128(_mm_load_ps(gbtab + _mm_extract_epi32(gbitab_hi,3)),
825 _mm_load_ps(gbtab + _mm_extract_epi32(gbitab_lo,3)));
826 GMX_MM256_HALFTRANSPOSE4_PS(Y,F,G,H);
827 Heps = _mm256_mul_ps(gbeps,H);
828 Fp = _mm256_add_ps(F,_mm256_mul_ps(gbeps,_mm256_add_ps(G,Heps)));
829 VV = _mm256_add_ps(Y,_mm256_mul_ps(gbeps,Fp));
830 vgb = _mm256_mul_ps(gbqqfactor,VV);
832 FF = _mm256_add_ps(Fp,_mm256_mul_ps(gbeps,_mm256_add_ps(G,_mm256_add_ps(Heps,Heps))));
833 fgb = _mm256_mul_ps(gbqqfactor,_mm256_mul_ps(FF,gbscale));
834 dvdatmp = _mm256_mul_ps(minushalf,_mm256_add_ps(vgb,_mm256_mul_ps(fgb,r00)));
835 dvdatmp = _mm256_andnot_ps(dummy_mask,dvdatmp);
836 dvdasum = _mm256_add_ps(dvdasum,dvdatmp);
837 /* The pointers to scratch make sure that this code with compilers that take gmx_restrict seriously (e.g. icc 13) really can't screw things up. */
838 fjptrA = (jnrlistA>=0) ? dvda+jnrA : scratch;
839 fjptrB = (jnrlistB>=0) ? dvda+jnrB : scratch;
840 fjptrC = (jnrlistC>=0) ? dvda+jnrC : scratch;
841 fjptrD = (jnrlistD>=0) ? dvda+jnrD : scratch;
842 fjptrE = (jnrlistE>=0) ? dvda+jnrE : scratch;
843 fjptrF = (jnrlistF>=0) ? dvda+jnrF : scratch;
844 fjptrG = (jnrlistG>=0) ? dvda+jnrG : scratch;
845 fjptrH = (jnrlistH>=0) ? dvda+jnrH : scratch;
846 gmx_mm256_increment_8real_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,
847 _mm256_mul_ps(dvdatmp,_mm256_mul_ps(isaj0,isaj0)));
848 velec = _mm256_mul_ps(qq00,rinv00);
849 felec = _mm256_mul_ps(_mm256_sub_ps(_mm256_mul_ps(velec,rinv00),fgb),rinv00);
853 fscal = _mm256_andnot_ps(dummy_mask,fscal);
855 /* Calculate temporary vectorial force */
856 tx = _mm256_mul_ps(fscal,dx00);
857 ty = _mm256_mul_ps(fscal,dy00);
858 tz = _mm256_mul_ps(fscal,dz00);
860 /* Update vectorial force */
861 fix0 = _mm256_add_ps(fix0,tx);
862 fiy0 = _mm256_add_ps(fiy0,ty);
863 fiz0 = _mm256_add_ps(fiz0,tz);
865 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
866 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
867 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
868 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
869 fjptrE = (jnrlistE>=0) ? f+j_coord_offsetE : scratch;
870 fjptrF = (jnrlistF>=0) ? f+j_coord_offsetF : scratch;
871 fjptrG = (jnrlistG>=0) ? f+j_coord_offsetG : scratch;
872 fjptrH = (jnrlistH>=0) ? f+j_coord_offsetH : scratch;
873 gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,tx,ty,tz);
875 /* Inner loop uses 56 flops */
878 /* End of innermost loop */
880 gmx_mm256_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
881 f+i_coord_offset,fshift+i_shift_offset);
883 dvdasum = _mm256_mul_ps(dvdasum, _mm256_mul_ps(isai0,isai0));
884 gmx_mm256_update_1pot_ps(dvdasum,dvda+inr);
886 /* Increment number of inner iterations */
887 inneriter += j_index_end - j_index_start;
889 /* Outer loop uses 7 flops */
892 /* Increment number of outer iterations */
895 /* Update outer/inner flops */
897 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_F,outeriter*7 + inneriter*56);