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_ElecGB_VdwLJ_GeomP1P1_VF_avx_256_single
54 * Electrostatics interaction: GeneralizedBorn
55 * VdW interaction: LennardJones
56 * Geometry: Particle-Particle
57 * Calculate force/pot: PotentialAndForce
60 nb_kernel_ElecGB_VdwLJ_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;
96 __m128i gbitab_lo,gbitab_hi;
97 __m256 vgb,fgb,vgbsum,dvdasum,gbscale,gbtabscale,isaprod,gbqqfactor,gbinvepsdiff,gbeps,dvdatmp;
98 __m256 minushalf = _mm256_set1_ps(-0.5);
99 real *invsqrta,*dvda,*gbtab;
101 __m256 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
104 __m256 one_sixth = _mm256_set1_ps(1.0/6.0);
105 __m256 one_twelfth = _mm256_set1_ps(1.0/12.0);
107 __m128i vfitab_lo,vfitab_hi;
108 __m128i ifour = _mm_set1_epi32(4);
109 __m256 rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF;
111 __m256 dummy_mask,cutoff_mask;
112 __m256 signbit = _mm256_castsi256_ps( _mm256_set1_epi32(0x80000000) );
113 __m256 one = _mm256_set1_ps(1.0);
114 __m256 two = _mm256_set1_ps(2.0);
120 jindex = nlist->jindex;
122 shiftidx = nlist->shift;
124 shiftvec = fr->shift_vec[0];
125 fshift = fr->fshift[0];
126 facel = _mm256_set1_ps(fr->epsfac);
127 charge = mdatoms->chargeA;
128 nvdwtype = fr->ntype;
130 vdwtype = mdatoms->typeA;
132 invsqrta = fr->invsqrta;
134 gbtabscale = _mm256_set1_ps(fr->gbtab.scale);
135 gbtab = fr->gbtab.data;
136 gbinvepsdiff = _mm256_set1_ps((1.0/fr->epsilon_r) - (1.0/fr->gb_epsilon_solvent));
138 /* Avoid stupid compiler warnings */
139 jnrA = jnrB = jnrC = jnrD = jnrE = jnrF = jnrG = jnrH = 0;
152 for(iidx=0;iidx<4*DIM;iidx++)
157 /* Start outer loop over neighborlists */
158 for(iidx=0; iidx<nri; iidx++)
160 /* Load shift vector for this list */
161 i_shift_offset = DIM*shiftidx[iidx];
163 /* Load limits for loop over neighbors */
164 j_index_start = jindex[iidx];
165 j_index_end = jindex[iidx+1];
167 /* Get outer coordinate index */
169 i_coord_offset = DIM*inr;
171 /* Load i particle coords and add shift vector */
172 gmx_mm256_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
174 fix0 = _mm256_setzero_ps();
175 fiy0 = _mm256_setzero_ps();
176 fiz0 = _mm256_setzero_ps();
178 /* Load parameters for i particles */
179 iq0 = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+0]));
180 isai0 = _mm256_set1_ps(invsqrta[inr+0]);
181 vdwioffsetptr0 = vdwparam+2*nvdwtype*vdwtype[inr+0];
183 /* Reset potential sums */
184 velecsum = _mm256_setzero_ps();
185 vgbsum = _mm256_setzero_ps();
186 vvdwsum = _mm256_setzero_ps();
187 dvdasum = _mm256_setzero_ps();
189 /* Start inner kernel loop */
190 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+7]>=0; jidx+=8)
193 /* Get j neighbor index, and coordinate index */
202 j_coord_offsetA = DIM*jnrA;
203 j_coord_offsetB = DIM*jnrB;
204 j_coord_offsetC = DIM*jnrC;
205 j_coord_offsetD = DIM*jnrD;
206 j_coord_offsetE = DIM*jnrE;
207 j_coord_offsetF = DIM*jnrF;
208 j_coord_offsetG = DIM*jnrG;
209 j_coord_offsetH = DIM*jnrH;
211 /* load j atom coordinates */
212 gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
213 x+j_coord_offsetC,x+j_coord_offsetD,
214 x+j_coord_offsetE,x+j_coord_offsetF,
215 x+j_coord_offsetG,x+j_coord_offsetH,
218 /* Calculate displacement vector */
219 dx00 = _mm256_sub_ps(ix0,jx0);
220 dy00 = _mm256_sub_ps(iy0,jy0);
221 dz00 = _mm256_sub_ps(iz0,jz0);
223 /* Calculate squared distance and things based on it */
224 rsq00 = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
226 rinv00 = gmx_mm256_invsqrt_ps(rsq00);
228 rinvsq00 = _mm256_mul_ps(rinv00,rinv00);
230 /* Load parameters for j particles */
231 jq0 = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
232 charge+jnrC+0,charge+jnrD+0,
233 charge+jnrE+0,charge+jnrF+0,
234 charge+jnrG+0,charge+jnrH+0);
235 isaj0 = gmx_mm256_load_8real_swizzle_ps(invsqrta+jnrA+0,invsqrta+jnrB+0,
236 invsqrta+jnrC+0,invsqrta+jnrD+0,
237 invsqrta+jnrE+0,invsqrta+jnrF+0,
238 invsqrta+jnrG+0,invsqrta+jnrH+0);
239 vdwjidx0A = 2*vdwtype[jnrA+0];
240 vdwjidx0B = 2*vdwtype[jnrB+0];
241 vdwjidx0C = 2*vdwtype[jnrC+0];
242 vdwjidx0D = 2*vdwtype[jnrD+0];
243 vdwjidx0E = 2*vdwtype[jnrE+0];
244 vdwjidx0F = 2*vdwtype[jnrF+0];
245 vdwjidx0G = 2*vdwtype[jnrG+0];
246 vdwjidx0H = 2*vdwtype[jnrH+0];
248 /**************************
249 * CALCULATE INTERACTIONS *
250 **************************/
252 r00 = _mm256_mul_ps(rsq00,rinv00);
254 /* Compute parameters for interactions between i and j atoms */
255 qq00 = _mm256_mul_ps(iq0,jq0);
256 gmx_mm256_load_8pair_swizzle_ps(vdwioffsetptr0+vdwjidx0A,
257 vdwioffsetptr0+vdwjidx0B,
258 vdwioffsetptr0+vdwjidx0C,
259 vdwioffsetptr0+vdwjidx0D,
260 vdwioffsetptr0+vdwjidx0E,
261 vdwioffsetptr0+vdwjidx0F,
262 vdwioffsetptr0+vdwjidx0G,
263 vdwioffsetptr0+vdwjidx0H,
266 /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
267 isaprod = _mm256_mul_ps(isai0,isaj0);
268 gbqqfactor = _mm256_xor_ps(signbit,_mm256_mul_ps(qq00,_mm256_mul_ps(isaprod,gbinvepsdiff)));
269 gbscale = _mm256_mul_ps(isaprod,gbtabscale);
271 /* Calculate generalized born table index - this is a separate table from the normal one,
272 * but we use the same procedure by multiplying r with scale and truncating to integer.
274 rt = _mm256_mul_ps(r00,gbscale);
275 gbitab = _mm256_cvttps_epi32(rt);
276 gbeps = _mm256_sub_ps(rt,_mm256_round_ps(rt, _MM_FROUND_FLOOR));
277 /* AVX1 does not support 256-bit integer operations, so now we go to 128-bit mode... */
278 gbitab_lo = _mm256_extractf128_si256(gbitab,0x0);
279 gbitab_hi = _mm256_extractf128_si256(gbitab,0x1);
280 gbitab_lo = _mm_slli_epi32(gbitab_lo,2);
281 gbitab_hi = _mm_slli_epi32(gbitab_hi,2);
282 Y = gmx_mm256_set_m128(_mm_load_ps(gbtab + _mm_extract_epi32(gbitab_hi,0)),
283 _mm_load_ps(gbtab + _mm_extract_epi32(gbitab_lo,0)));
284 F = gmx_mm256_set_m128(_mm_load_ps(gbtab + _mm_extract_epi32(gbitab_hi,1)),
285 _mm_load_ps(gbtab + _mm_extract_epi32(gbitab_lo,1)));
286 G = gmx_mm256_set_m128(_mm_load_ps(gbtab + _mm_extract_epi32(gbitab_hi,2)),
287 _mm_load_ps(gbtab + _mm_extract_epi32(gbitab_lo,2)));
288 H = gmx_mm256_set_m128(_mm_load_ps(gbtab + _mm_extract_epi32(gbitab_hi,3)),
289 _mm_load_ps(gbtab + _mm_extract_epi32(gbitab_lo,3)));
290 GMX_MM256_HALFTRANSPOSE4_PS(Y,F,G,H);
291 Heps = _mm256_mul_ps(gbeps,H);
292 Fp = _mm256_add_ps(F,_mm256_mul_ps(gbeps,_mm256_add_ps(G,Heps)));
293 VV = _mm256_add_ps(Y,_mm256_mul_ps(gbeps,Fp));
294 vgb = _mm256_mul_ps(gbqqfactor,VV);
296 FF = _mm256_add_ps(Fp,_mm256_mul_ps(gbeps,_mm256_add_ps(G,_mm256_add_ps(Heps,Heps))));
297 fgb = _mm256_mul_ps(gbqqfactor,_mm256_mul_ps(FF,gbscale));
298 dvdatmp = _mm256_mul_ps(minushalf,_mm256_add_ps(vgb,_mm256_mul_ps(fgb,r00)));
299 dvdasum = _mm256_add_ps(dvdasum,dvdatmp);
308 gmx_mm256_increment_8real_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,
309 _mm256_mul_ps(dvdatmp,_mm256_mul_ps(isaj0,isaj0)));
310 velec = _mm256_mul_ps(qq00,rinv00);
311 felec = _mm256_mul_ps(_mm256_sub_ps(_mm256_mul_ps(velec,rinv00),fgb),rinv00);
313 /* LENNARD-JONES DISPERSION/REPULSION */
315 rinvsix = _mm256_mul_ps(_mm256_mul_ps(rinvsq00,rinvsq00),rinvsq00);
316 vvdw6 = _mm256_mul_ps(c6_00,rinvsix);
317 vvdw12 = _mm256_mul_ps(c12_00,_mm256_mul_ps(rinvsix,rinvsix));
318 vvdw = _mm256_sub_ps( _mm256_mul_ps(vvdw12,one_twelfth) , _mm256_mul_ps(vvdw6,one_sixth) );
319 fvdw = _mm256_mul_ps(_mm256_sub_ps(vvdw12,vvdw6),rinvsq00);
321 /* Update potential sum for this i atom from the interaction with this j atom. */
322 velecsum = _mm256_add_ps(velecsum,velec);
323 vgbsum = _mm256_add_ps(vgbsum,vgb);
324 vvdwsum = _mm256_add_ps(vvdwsum,vvdw);
326 fscal = _mm256_add_ps(felec,fvdw);
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 = f+j_coord_offsetA;
339 fjptrB = f+j_coord_offsetB;
340 fjptrC = f+j_coord_offsetC;
341 fjptrD = f+j_coord_offsetD;
342 fjptrE = f+j_coord_offsetE;
343 fjptrF = f+j_coord_offsetF;
344 fjptrG = f+j_coord_offsetG;
345 fjptrH = f+j_coord_offsetH;
346 gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,tx,ty,tz);
348 /* Inner loop uses 70 flops */
354 /* Get j neighbor index, and coordinate index */
355 jnrlistA = jjnr[jidx];
356 jnrlistB = jjnr[jidx+1];
357 jnrlistC = jjnr[jidx+2];
358 jnrlistD = jjnr[jidx+3];
359 jnrlistE = jjnr[jidx+4];
360 jnrlistF = jjnr[jidx+5];
361 jnrlistG = jjnr[jidx+6];
362 jnrlistH = jjnr[jidx+7];
363 /* Sign of each element will be negative for non-real atoms.
364 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
365 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
367 dummy_mask = gmx_mm256_set_m128(gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx+4)),_mm_setzero_si128())),
368 gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128())));
370 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
371 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
372 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
373 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
374 jnrE = (jnrlistE>=0) ? jnrlistE : 0;
375 jnrF = (jnrlistF>=0) ? jnrlistF : 0;
376 jnrG = (jnrlistG>=0) ? jnrlistG : 0;
377 jnrH = (jnrlistH>=0) ? jnrlistH : 0;
378 j_coord_offsetA = DIM*jnrA;
379 j_coord_offsetB = DIM*jnrB;
380 j_coord_offsetC = DIM*jnrC;
381 j_coord_offsetD = DIM*jnrD;
382 j_coord_offsetE = DIM*jnrE;
383 j_coord_offsetF = DIM*jnrF;
384 j_coord_offsetG = DIM*jnrG;
385 j_coord_offsetH = DIM*jnrH;
387 /* load j atom coordinates */
388 gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
389 x+j_coord_offsetC,x+j_coord_offsetD,
390 x+j_coord_offsetE,x+j_coord_offsetF,
391 x+j_coord_offsetG,x+j_coord_offsetH,
394 /* Calculate displacement vector */
395 dx00 = _mm256_sub_ps(ix0,jx0);
396 dy00 = _mm256_sub_ps(iy0,jy0);
397 dz00 = _mm256_sub_ps(iz0,jz0);
399 /* Calculate squared distance and things based on it */
400 rsq00 = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
402 rinv00 = gmx_mm256_invsqrt_ps(rsq00);
404 rinvsq00 = _mm256_mul_ps(rinv00,rinv00);
406 /* Load parameters for j particles */
407 jq0 = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
408 charge+jnrC+0,charge+jnrD+0,
409 charge+jnrE+0,charge+jnrF+0,
410 charge+jnrG+0,charge+jnrH+0);
411 isaj0 = gmx_mm256_load_8real_swizzle_ps(invsqrta+jnrA+0,invsqrta+jnrB+0,
412 invsqrta+jnrC+0,invsqrta+jnrD+0,
413 invsqrta+jnrE+0,invsqrta+jnrF+0,
414 invsqrta+jnrG+0,invsqrta+jnrH+0);
415 vdwjidx0A = 2*vdwtype[jnrA+0];
416 vdwjidx0B = 2*vdwtype[jnrB+0];
417 vdwjidx0C = 2*vdwtype[jnrC+0];
418 vdwjidx0D = 2*vdwtype[jnrD+0];
419 vdwjidx0E = 2*vdwtype[jnrE+0];
420 vdwjidx0F = 2*vdwtype[jnrF+0];
421 vdwjidx0G = 2*vdwtype[jnrG+0];
422 vdwjidx0H = 2*vdwtype[jnrH+0];
424 /**************************
425 * CALCULATE INTERACTIONS *
426 **************************/
428 r00 = _mm256_mul_ps(rsq00,rinv00);
429 r00 = _mm256_andnot_ps(dummy_mask,r00);
431 /* Compute parameters for interactions between i and j atoms */
432 qq00 = _mm256_mul_ps(iq0,jq0);
433 gmx_mm256_load_8pair_swizzle_ps(vdwioffsetptr0+vdwjidx0A,
434 vdwioffsetptr0+vdwjidx0B,
435 vdwioffsetptr0+vdwjidx0C,
436 vdwioffsetptr0+vdwjidx0D,
437 vdwioffsetptr0+vdwjidx0E,
438 vdwioffsetptr0+vdwjidx0F,
439 vdwioffsetptr0+vdwjidx0G,
440 vdwioffsetptr0+vdwjidx0H,
443 /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
444 isaprod = _mm256_mul_ps(isai0,isaj0);
445 gbqqfactor = _mm256_xor_ps(signbit,_mm256_mul_ps(qq00,_mm256_mul_ps(isaprod,gbinvepsdiff)));
446 gbscale = _mm256_mul_ps(isaprod,gbtabscale);
448 /* Calculate generalized born table index - this is a separate table from the normal one,
449 * but we use the same procedure by multiplying r with scale and truncating to integer.
451 rt = _mm256_mul_ps(r00,gbscale);
452 gbitab = _mm256_cvttps_epi32(rt);
453 gbeps = _mm256_sub_ps(rt,_mm256_round_ps(rt, _MM_FROUND_FLOOR));
454 /* AVX1 does not support 256-bit integer operations, so now we go to 128-bit mode... */
455 gbitab_lo = _mm256_extractf128_si256(gbitab,0x0);
456 gbitab_hi = _mm256_extractf128_si256(gbitab,0x1);
457 gbitab_lo = _mm_slli_epi32(gbitab_lo,2);
458 gbitab_hi = _mm_slli_epi32(gbitab_hi,2);
459 Y = gmx_mm256_set_m128(_mm_load_ps(gbtab + _mm_extract_epi32(gbitab_hi,0)),
460 _mm_load_ps(gbtab + _mm_extract_epi32(gbitab_lo,0)));
461 F = gmx_mm256_set_m128(_mm_load_ps(gbtab + _mm_extract_epi32(gbitab_hi,1)),
462 _mm_load_ps(gbtab + _mm_extract_epi32(gbitab_lo,1)));
463 G = gmx_mm256_set_m128(_mm_load_ps(gbtab + _mm_extract_epi32(gbitab_hi,2)),
464 _mm_load_ps(gbtab + _mm_extract_epi32(gbitab_lo,2)));
465 H = gmx_mm256_set_m128(_mm_load_ps(gbtab + _mm_extract_epi32(gbitab_hi,3)),
466 _mm_load_ps(gbtab + _mm_extract_epi32(gbitab_lo,3)));
467 GMX_MM256_HALFTRANSPOSE4_PS(Y,F,G,H);
468 Heps = _mm256_mul_ps(gbeps,H);
469 Fp = _mm256_add_ps(F,_mm256_mul_ps(gbeps,_mm256_add_ps(G,Heps)));
470 VV = _mm256_add_ps(Y,_mm256_mul_ps(gbeps,Fp));
471 vgb = _mm256_mul_ps(gbqqfactor,VV);
473 FF = _mm256_add_ps(Fp,_mm256_mul_ps(gbeps,_mm256_add_ps(G,_mm256_add_ps(Heps,Heps))));
474 fgb = _mm256_mul_ps(gbqqfactor,_mm256_mul_ps(FF,gbscale));
475 dvdatmp = _mm256_mul_ps(minushalf,_mm256_add_ps(vgb,_mm256_mul_ps(fgb,r00)));
476 dvdatmp = _mm256_andnot_ps(dummy_mask,dvdatmp);
477 dvdasum = _mm256_add_ps(dvdasum,dvdatmp);
478 /* 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. */
479 fjptrA = (jnrlistA>=0) ? dvda+jnrA : scratch;
480 fjptrB = (jnrlistB>=0) ? dvda+jnrB : scratch;
481 fjptrC = (jnrlistC>=0) ? dvda+jnrC : scratch;
482 fjptrD = (jnrlistD>=0) ? dvda+jnrD : scratch;
483 fjptrE = (jnrlistE>=0) ? dvda+jnrE : scratch;
484 fjptrF = (jnrlistF>=0) ? dvda+jnrF : scratch;
485 fjptrG = (jnrlistG>=0) ? dvda+jnrG : scratch;
486 fjptrH = (jnrlistH>=0) ? dvda+jnrH : scratch;
487 gmx_mm256_increment_8real_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,
488 _mm256_mul_ps(dvdatmp,_mm256_mul_ps(isaj0,isaj0)));
489 velec = _mm256_mul_ps(qq00,rinv00);
490 felec = _mm256_mul_ps(_mm256_sub_ps(_mm256_mul_ps(velec,rinv00),fgb),rinv00);
492 /* LENNARD-JONES DISPERSION/REPULSION */
494 rinvsix = _mm256_mul_ps(_mm256_mul_ps(rinvsq00,rinvsq00),rinvsq00);
495 vvdw6 = _mm256_mul_ps(c6_00,rinvsix);
496 vvdw12 = _mm256_mul_ps(c12_00,_mm256_mul_ps(rinvsix,rinvsix));
497 vvdw = _mm256_sub_ps( _mm256_mul_ps(vvdw12,one_twelfth) , _mm256_mul_ps(vvdw6,one_sixth) );
498 fvdw = _mm256_mul_ps(_mm256_sub_ps(vvdw12,vvdw6),rinvsq00);
500 /* Update potential sum for this i atom from the interaction with this j atom. */
501 velec = _mm256_andnot_ps(dummy_mask,velec);
502 velecsum = _mm256_add_ps(velecsum,velec);
503 vgb = _mm256_andnot_ps(dummy_mask,vgb);
504 vgbsum = _mm256_add_ps(vgbsum,vgb);
505 vvdw = _mm256_andnot_ps(dummy_mask,vvdw);
506 vvdwsum = _mm256_add_ps(vvdwsum,vvdw);
508 fscal = _mm256_add_ps(felec,fvdw);
510 fscal = _mm256_andnot_ps(dummy_mask,fscal);
512 /* Calculate temporary vectorial force */
513 tx = _mm256_mul_ps(fscal,dx00);
514 ty = _mm256_mul_ps(fscal,dy00);
515 tz = _mm256_mul_ps(fscal,dz00);
517 /* Update vectorial force */
518 fix0 = _mm256_add_ps(fix0,tx);
519 fiy0 = _mm256_add_ps(fiy0,ty);
520 fiz0 = _mm256_add_ps(fiz0,tz);
522 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
523 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
524 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
525 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
526 fjptrE = (jnrlistE>=0) ? f+j_coord_offsetE : scratch;
527 fjptrF = (jnrlistF>=0) ? f+j_coord_offsetF : scratch;
528 fjptrG = (jnrlistG>=0) ? f+j_coord_offsetG : scratch;
529 fjptrH = (jnrlistH>=0) ? f+j_coord_offsetH : scratch;
530 gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,tx,ty,tz);
532 /* Inner loop uses 71 flops */
535 /* End of innermost loop */
537 gmx_mm256_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
538 f+i_coord_offset,fshift+i_shift_offset);
541 /* Update potential energies */
542 gmx_mm256_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
543 gmx_mm256_update_1pot_ps(vgbsum,kernel_data->energygrp_polarization+ggid);
544 gmx_mm256_update_1pot_ps(vvdwsum,kernel_data->energygrp_vdw+ggid);
545 dvdasum = _mm256_mul_ps(dvdasum, _mm256_mul_ps(isai0,isai0));
546 gmx_mm256_update_1pot_ps(dvdasum,dvda+inr);
548 /* Increment number of inner iterations */
549 inneriter += j_index_end - j_index_start;
551 /* Outer loop uses 10 flops */
554 /* Increment number of outer iterations */
557 /* Update outer/inner flops */
559 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_VF,outeriter*10 + inneriter*71);
562 * Gromacs nonbonded kernel: nb_kernel_ElecGB_VdwLJ_GeomP1P1_F_avx_256_single
563 * Electrostatics interaction: GeneralizedBorn
564 * VdW interaction: LennardJones
565 * Geometry: Particle-Particle
566 * Calculate force/pot: Force
569 nb_kernel_ElecGB_VdwLJ_GeomP1P1_F_avx_256_single
570 (t_nblist * gmx_restrict nlist,
571 rvec * gmx_restrict xx,
572 rvec * gmx_restrict ff,
573 t_forcerec * gmx_restrict fr,
574 t_mdatoms * gmx_restrict mdatoms,
575 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
576 t_nrnb * gmx_restrict nrnb)
578 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
579 * just 0 for non-waters.
580 * Suffixes A,B,C,D,E,F,G,H refer to j loop unrolling done with AVX, e.g. for the eight different
581 * jnr indices corresponding to data put in the four positions in the SIMD register.
583 int i_shift_offset,i_coord_offset,outeriter,inneriter;
584 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
585 int jnrA,jnrB,jnrC,jnrD;
586 int jnrE,jnrF,jnrG,jnrH;
587 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
588 int jnrlistE,jnrlistF,jnrlistG,jnrlistH;
589 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
590 int j_coord_offsetE,j_coord_offsetF,j_coord_offsetG,j_coord_offsetH;
591 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
593 real *shiftvec,*fshift,*x,*f;
594 real *fjptrA,*fjptrB,*fjptrC,*fjptrD,*fjptrE,*fjptrF,*fjptrG,*fjptrH;
596 __m256 tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
597 real * vdwioffsetptr0;
598 __m256 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
599 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D,vdwjidx0E,vdwjidx0F,vdwjidx0G,vdwjidx0H;
600 __m256 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
601 __m256 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
602 __m256 velec,felec,velecsum,facel,crf,krf,krf2;
605 __m128i gbitab_lo,gbitab_hi;
606 __m256 vgb,fgb,vgbsum,dvdasum,gbscale,gbtabscale,isaprod,gbqqfactor,gbinvepsdiff,gbeps,dvdatmp;
607 __m256 minushalf = _mm256_set1_ps(-0.5);
608 real *invsqrta,*dvda,*gbtab;
610 __m256 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
613 __m256 one_sixth = _mm256_set1_ps(1.0/6.0);
614 __m256 one_twelfth = _mm256_set1_ps(1.0/12.0);
616 __m128i vfitab_lo,vfitab_hi;
617 __m128i ifour = _mm_set1_epi32(4);
618 __m256 rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF;
620 __m256 dummy_mask,cutoff_mask;
621 __m256 signbit = _mm256_castsi256_ps( _mm256_set1_epi32(0x80000000) );
622 __m256 one = _mm256_set1_ps(1.0);
623 __m256 two = _mm256_set1_ps(2.0);
629 jindex = nlist->jindex;
631 shiftidx = nlist->shift;
633 shiftvec = fr->shift_vec[0];
634 fshift = fr->fshift[0];
635 facel = _mm256_set1_ps(fr->epsfac);
636 charge = mdatoms->chargeA;
637 nvdwtype = fr->ntype;
639 vdwtype = mdatoms->typeA;
641 invsqrta = fr->invsqrta;
643 gbtabscale = _mm256_set1_ps(fr->gbtab.scale);
644 gbtab = fr->gbtab.data;
645 gbinvepsdiff = _mm256_set1_ps((1.0/fr->epsilon_r) - (1.0/fr->gb_epsilon_solvent));
647 /* Avoid stupid compiler warnings */
648 jnrA = jnrB = jnrC = jnrD = jnrE = jnrF = jnrG = jnrH = 0;
661 for(iidx=0;iidx<4*DIM;iidx++)
666 /* Start outer loop over neighborlists */
667 for(iidx=0; iidx<nri; iidx++)
669 /* Load shift vector for this list */
670 i_shift_offset = DIM*shiftidx[iidx];
672 /* Load limits for loop over neighbors */
673 j_index_start = jindex[iidx];
674 j_index_end = jindex[iidx+1];
676 /* Get outer coordinate index */
678 i_coord_offset = DIM*inr;
680 /* Load i particle coords and add shift vector */
681 gmx_mm256_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
683 fix0 = _mm256_setzero_ps();
684 fiy0 = _mm256_setzero_ps();
685 fiz0 = _mm256_setzero_ps();
687 /* Load parameters for i particles */
688 iq0 = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+0]));
689 isai0 = _mm256_set1_ps(invsqrta[inr+0]);
690 vdwioffsetptr0 = vdwparam+2*nvdwtype*vdwtype[inr+0];
692 dvdasum = _mm256_setzero_ps();
694 /* Start inner kernel loop */
695 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+7]>=0; jidx+=8)
698 /* Get j neighbor index, and coordinate index */
707 j_coord_offsetA = DIM*jnrA;
708 j_coord_offsetB = DIM*jnrB;
709 j_coord_offsetC = DIM*jnrC;
710 j_coord_offsetD = DIM*jnrD;
711 j_coord_offsetE = DIM*jnrE;
712 j_coord_offsetF = DIM*jnrF;
713 j_coord_offsetG = DIM*jnrG;
714 j_coord_offsetH = DIM*jnrH;
716 /* load j atom coordinates */
717 gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
718 x+j_coord_offsetC,x+j_coord_offsetD,
719 x+j_coord_offsetE,x+j_coord_offsetF,
720 x+j_coord_offsetG,x+j_coord_offsetH,
723 /* Calculate displacement vector */
724 dx00 = _mm256_sub_ps(ix0,jx0);
725 dy00 = _mm256_sub_ps(iy0,jy0);
726 dz00 = _mm256_sub_ps(iz0,jz0);
728 /* Calculate squared distance and things based on it */
729 rsq00 = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
731 rinv00 = gmx_mm256_invsqrt_ps(rsq00);
733 rinvsq00 = _mm256_mul_ps(rinv00,rinv00);
735 /* Load parameters for j particles */
736 jq0 = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
737 charge+jnrC+0,charge+jnrD+0,
738 charge+jnrE+0,charge+jnrF+0,
739 charge+jnrG+0,charge+jnrH+0);
740 isaj0 = gmx_mm256_load_8real_swizzle_ps(invsqrta+jnrA+0,invsqrta+jnrB+0,
741 invsqrta+jnrC+0,invsqrta+jnrD+0,
742 invsqrta+jnrE+0,invsqrta+jnrF+0,
743 invsqrta+jnrG+0,invsqrta+jnrH+0);
744 vdwjidx0A = 2*vdwtype[jnrA+0];
745 vdwjidx0B = 2*vdwtype[jnrB+0];
746 vdwjidx0C = 2*vdwtype[jnrC+0];
747 vdwjidx0D = 2*vdwtype[jnrD+0];
748 vdwjidx0E = 2*vdwtype[jnrE+0];
749 vdwjidx0F = 2*vdwtype[jnrF+0];
750 vdwjidx0G = 2*vdwtype[jnrG+0];
751 vdwjidx0H = 2*vdwtype[jnrH+0];
753 /**************************
754 * CALCULATE INTERACTIONS *
755 **************************/
757 r00 = _mm256_mul_ps(rsq00,rinv00);
759 /* Compute parameters for interactions between i and j atoms */
760 qq00 = _mm256_mul_ps(iq0,jq0);
761 gmx_mm256_load_8pair_swizzle_ps(vdwioffsetptr0+vdwjidx0A,
762 vdwioffsetptr0+vdwjidx0B,
763 vdwioffsetptr0+vdwjidx0C,
764 vdwioffsetptr0+vdwjidx0D,
765 vdwioffsetptr0+vdwjidx0E,
766 vdwioffsetptr0+vdwjidx0F,
767 vdwioffsetptr0+vdwjidx0G,
768 vdwioffsetptr0+vdwjidx0H,
771 /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
772 isaprod = _mm256_mul_ps(isai0,isaj0);
773 gbqqfactor = _mm256_xor_ps(signbit,_mm256_mul_ps(qq00,_mm256_mul_ps(isaprod,gbinvepsdiff)));
774 gbscale = _mm256_mul_ps(isaprod,gbtabscale);
776 /* Calculate generalized born table index - this is a separate table from the normal one,
777 * but we use the same procedure by multiplying r with scale and truncating to integer.
779 rt = _mm256_mul_ps(r00,gbscale);
780 gbitab = _mm256_cvttps_epi32(rt);
781 gbeps = _mm256_sub_ps(rt,_mm256_round_ps(rt, _MM_FROUND_FLOOR));
782 /* AVX1 does not support 256-bit integer operations, so now we go to 128-bit mode... */
783 gbitab_lo = _mm256_extractf128_si256(gbitab,0x0);
784 gbitab_hi = _mm256_extractf128_si256(gbitab,0x1);
785 gbitab_lo = _mm_slli_epi32(gbitab_lo,2);
786 gbitab_hi = _mm_slli_epi32(gbitab_hi,2);
787 Y = gmx_mm256_set_m128(_mm_load_ps(gbtab + _mm_extract_epi32(gbitab_hi,0)),
788 _mm_load_ps(gbtab + _mm_extract_epi32(gbitab_lo,0)));
789 F = gmx_mm256_set_m128(_mm_load_ps(gbtab + _mm_extract_epi32(gbitab_hi,1)),
790 _mm_load_ps(gbtab + _mm_extract_epi32(gbitab_lo,1)));
791 G = gmx_mm256_set_m128(_mm_load_ps(gbtab + _mm_extract_epi32(gbitab_hi,2)),
792 _mm_load_ps(gbtab + _mm_extract_epi32(gbitab_lo,2)));
793 H = gmx_mm256_set_m128(_mm_load_ps(gbtab + _mm_extract_epi32(gbitab_hi,3)),
794 _mm_load_ps(gbtab + _mm_extract_epi32(gbitab_lo,3)));
795 GMX_MM256_HALFTRANSPOSE4_PS(Y,F,G,H);
796 Heps = _mm256_mul_ps(gbeps,H);
797 Fp = _mm256_add_ps(F,_mm256_mul_ps(gbeps,_mm256_add_ps(G,Heps)));
798 VV = _mm256_add_ps(Y,_mm256_mul_ps(gbeps,Fp));
799 vgb = _mm256_mul_ps(gbqqfactor,VV);
801 FF = _mm256_add_ps(Fp,_mm256_mul_ps(gbeps,_mm256_add_ps(G,_mm256_add_ps(Heps,Heps))));
802 fgb = _mm256_mul_ps(gbqqfactor,_mm256_mul_ps(FF,gbscale));
803 dvdatmp = _mm256_mul_ps(minushalf,_mm256_add_ps(vgb,_mm256_mul_ps(fgb,r00)));
804 dvdasum = _mm256_add_ps(dvdasum,dvdatmp);
813 gmx_mm256_increment_8real_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,
814 _mm256_mul_ps(dvdatmp,_mm256_mul_ps(isaj0,isaj0)));
815 velec = _mm256_mul_ps(qq00,rinv00);
816 felec = _mm256_mul_ps(_mm256_sub_ps(_mm256_mul_ps(velec,rinv00),fgb),rinv00);
818 /* LENNARD-JONES DISPERSION/REPULSION */
820 rinvsix = _mm256_mul_ps(_mm256_mul_ps(rinvsq00,rinvsq00),rinvsq00);
821 fvdw = _mm256_mul_ps(_mm256_sub_ps(_mm256_mul_ps(c12_00,rinvsix),c6_00),_mm256_mul_ps(rinvsix,rinvsq00));
823 fscal = _mm256_add_ps(felec,fvdw);
825 /* Calculate temporary vectorial force */
826 tx = _mm256_mul_ps(fscal,dx00);
827 ty = _mm256_mul_ps(fscal,dy00);
828 tz = _mm256_mul_ps(fscal,dz00);
830 /* Update vectorial force */
831 fix0 = _mm256_add_ps(fix0,tx);
832 fiy0 = _mm256_add_ps(fiy0,ty);
833 fiz0 = _mm256_add_ps(fiz0,tz);
835 fjptrA = f+j_coord_offsetA;
836 fjptrB = f+j_coord_offsetB;
837 fjptrC = f+j_coord_offsetC;
838 fjptrD = f+j_coord_offsetD;
839 fjptrE = f+j_coord_offsetE;
840 fjptrF = f+j_coord_offsetF;
841 fjptrG = f+j_coord_offsetG;
842 fjptrH = f+j_coord_offsetH;
843 gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,tx,ty,tz);
845 /* Inner loop uses 63 flops */
851 /* Get j neighbor index, and coordinate index */
852 jnrlistA = jjnr[jidx];
853 jnrlistB = jjnr[jidx+1];
854 jnrlistC = jjnr[jidx+2];
855 jnrlistD = jjnr[jidx+3];
856 jnrlistE = jjnr[jidx+4];
857 jnrlistF = jjnr[jidx+5];
858 jnrlistG = jjnr[jidx+6];
859 jnrlistH = jjnr[jidx+7];
860 /* Sign of each element will be negative for non-real atoms.
861 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
862 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
864 dummy_mask = gmx_mm256_set_m128(gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx+4)),_mm_setzero_si128())),
865 gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128())));
867 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
868 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
869 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
870 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
871 jnrE = (jnrlistE>=0) ? jnrlistE : 0;
872 jnrF = (jnrlistF>=0) ? jnrlistF : 0;
873 jnrG = (jnrlistG>=0) ? jnrlistG : 0;
874 jnrH = (jnrlistH>=0) ? jnrlistH : 0;
875 j_coord_offsetA = DIM*jnrA;
876 j_coord_offsetB = DIM*jnrB;
877 j_coord_offsetC = DIM*jnrC;
878 j_coord_offsetD = DIM*jnrD;
879 j_coord_offsetE = DIM*jnrE;
880 j_coord_offsetF = DIM*jnrF;
881 j_coord_offsetG = DIM*jnrG;
882 j_coord_offsetH = DIM*jnrH;
884 /* load j atom coordinates */
885 gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
886 x+j_coord_offsetC,x+j_coord_offsetD,
887 x+j_coord_offsetE,x+j_coord_offsetF,
888 x+j_coord_offsetG,x+j_coord_offsetH,
891 /* Calculate displacement vector */
892 dx00 = _mm256_sub_ps(ix0,jx0);
893 dy00 = _mm256_sub_ps(iy0,jy0);
894 dz00 = _mm256_sub_ps(iz0,jz0);
896 /* Calculate squared distance and things based on it */
897 rsq00 = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
899 rinv00 = gmx_mm256_invsqrt_ps(rsq00);
901 rinvsq00 = _mm256_mul_ps(rinv00,rinv00);
903 /* Load parameters for j particles */
904 jq0 = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
905 charge+jnrC+0,charge+jnrD+0,
906 charge+jnrE+0,charge+jnrF+0,
907 charge+jnrG+0,charge+jnrH+0);
908 isaj0 = gmx_mm256_load_8real_swizzle_ps(invsqrta+jnrA+0,invsqrta+jnrB+0,
909 invsqrta+jnrC+0,invsqrta+jnrD+0,
910 invsqrta+jnrE+0,invsqrta+jnrF+0,
911 invsqrta+jnrG+0,invsqrta+jnrH+0);
912 vdwjidx0A = 2*vdwtype[jnrA+0];
913 vdwjidx0B = 2*vdwtype[jnrB+0];
914 vdwjidx0C = 2*vdwtype[jnrC+0];
915 vdwjidx0D = 2*vdwtype[jnrD+0];
916 vdwjidx0E = 2*vdwtype[jnrE+0];
917 vdwjidx0F = 2*vdwtype[jnrF+0];
918 vdwjidx0G = 2*vdwtype[jnrG+0];
919 vdwjidx0H = 2*vdwtype[jnrH+0];
921 /**************************
922 * CALCULATE INTERACTIONS *
923 **************************/
925 r00 = _mm256_mul_ps(rsq00,rinv00);
926 r00 = _mm256_andnot_ps(dummy_mask,r00);
928 /* Compute parameters for interactions between i and j atoms */
929 qq00 = _mm256_mul_ps(iq0,jq0);
930 gmx_mm256_load_8pair_swizzle_ps(vdwioffsetptr0+vdwjidx0A,
931 vdwioffsetptr0+vdwjidx0B,
932 vdwioffsetptr0+vdwjidx0C,
933 vdwioffsetptr0+vdwjidx0D,
934 vdwioffsetptr0+vdwjidx0E,
935 vdwioffsetptr0+vdwjidx0F,
936 vdwioffsetptr0+vdwjidx0G,
937 vdwioffsetptr0+vdwjidx0H,
940 /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
941 isaprod = _mm256_mul_ps(isai0,isaj0);
942 gbqqfactor = _mm256_xor_ps(signbit,_mm256_mul_ps(qq00,_mm256_mul_ps(isaprod,gbinvepsdiff)));
943 gbscale = _mm256_mul_ps(isaprod,gbtabscale);
945 /* Calculate generalized born table index - this is a separate table from the normal one,
946 * but we use the same procedure by multiplying r with scale and truncating to integer.
948 rt = _mm256_mul_ps(r00,gbscale);
949 gbitab = _mm256_cvttps_epi32(rt);
950 gbeps = _mm256_sub_ps(rt,_mm256_round_ps(rt, _MM_FROUND_FLOOR));
951 /* AVX1 does not support 256-bit integer operations, so now we go to 128-bit mode... */
952 gbitab_lo = _mm256_extractf128_si256(gbitab,0x0);
953 gbitab_hi = _mm256_extractf128_si256(gbitab,0x1);
954 gbitab_lo = _mm_slli_epi32(gbitab_lo,2);
955 gbitab_hi = _mm_slli_epi32(gbitab_hi,2);
956 Y = gmx_mm256_set_m128(_mm_load_ps(gbtab + _mm_extract_epi32(gbitab_hi,0)),
957 _mm_load_ps(gbtab + _mm_extract_epi32(gbitab_lo,0)));
958 F = gmx_mm256_set_m128(_mm_load_ps(gbtab + _mm_extract_epi32(gbitab_hi,1)),
959 _mm_load_ps(gbtab + _mm_extract_epi32(gbitab_lo,1)));
960 G = gmx_mm256_set_m128(_mm_load_ps(gbtab + _mm_extract_epi32(gbitab_hi,2)),
961 _mm_load_ps(gbtab + _mm_extract_epi32(gbitab_lo,2)));
962 H = gmx_mm256_set_m128(_mm_load_ps(gbtab + _mm_extract_epi32(gbitab_hi,3)),
963 _mm_load_ps(gbtab + _mm_extract_epi32(gbitab_lo,3)));
964 GMX_MM256_HALFTRANSPOSE4_PS(Y,F,G,H);
965 Heps = _mm256_mul_ps(gbeps,H);
966 Fp = _mm256_add_ps(F,_mm256_mul_ps(gbeps,_mm256_add_ps(G,Heps)));
967 VV = _mm256_add_ps(Y,_mm256_mul_ps(gbeps,Fp));
968 vgb = _mm256_mul_ps(gbqqfactor,VV);
970 FF = _mm256_add_ps(Fp,_mm256_mul_ps(gbeps,_mm256_add_ps(G,_mm256_add_ps(Heps,Heps))));
971 fgb = _mm256_mul_ps(gbqqfactor,_mm256_mul_ps(FF,gbscale));
972 dvdatmp = _mm256_mul_ps(minushalf,_mm256_add_ps(vgb,_mm256_mul_ps(fgb,r00)));
973 dvdatmp = _mm256_andnot_ps(dummy_mask,dvdatmp);
974 dvdasum = _mm256_add_ps(dvdasum,dvdatmp);
975 /* 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. */
976 fjptrA = (jnrlistA>=0) ? dvda+jnrA : scratch;
977 fjptrB = (jnrlistB>=0) ? dvda+jnrB : scratch;
978 fjptrC = (jnrlistC>=0) ? dvda+jnrC : scratch;
979 fjptrD = (jnrlistD>=0) ? dvda+jnrD : scratch;
980 fjptrE = (jnrlistE>=0) ? dvda+jnrE : scratch;
981 fjptrF = (jnrlistF>=0) ? dvda+jnrF : scratch;
982 fjptrG = (jnrlistG>=0) ? dvda+jnrG : scratch;
983 fjptrH = (jnrlistH>=0) ? dvda+jnrH : scratch;
984 gmx_mm256_increment_8real_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,
985 _mm256_mul_ps(dvdatmp,_mm256_mul_ps(isaj0,isaj0)));
986 velec = _mm256_mul_ps(qq00,rinv00);
987 felec = _mm256_mul_ps(_mm256_sub_ps(_mm256_mul_ps(velec,rinv00),fgb),rinv00);
989 /* LENNARD-JONES DISPERSION/REPULSION */
991 rinvsix = _mm256_mul_ps(_mm256_mul_ps(rinvsq00,rinvsq00),rinvsq00);
992 fvdw = _mm256_mul_ps(_mm256_sub_ps(_mm256_mul_ps(c12_00,rinvsix),c6_00),_mm256_mul_ps(rinvsix,rinvsq00));
994 fscal = _mm256_add_ps(felec,fvdw);
996 fscal = _mm256_andnot_ps(dummy_mask,fscal);
998 /* Calculate temporary vectorial force */
999 tx = _mm256_mul_ps(fscal,dx00);
1000 ty = _mm256_mul_ps(fscal,dy00);
1001 tz = _mm256_mul_ps(fscal,dz00);
1003 /* Update vectorial force */
1004 fix0 = _mm256_add_ps(fix0,tx);
1005 fiy0 = _mm256_add_ps(fiy0,ty);
1006 fiz0 = _mm256_add_ps(fiz0,tz);
1008 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
1009 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
1010 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
1011 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
1012 fjptrE = (jnrlistE>=0) ? f+j_coord_offsetE : scratch;
1013 fjptrF = (jnrlistF>=0) ? f+j_coord_offsetF : scratch;
1014 fjptrG = (jnrlistG>=0) ? f+j_coord_offsetG : scratch;
1015 fjptrH = (jnrlistH>=0) ? f+j_coord_offsetH : scratch;
1016 gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,tx,ty,tz);
1018 /* Inner loop uses 64 flops */
1021 /* End of innermost loop */
1023 gmx_mm256_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
1024 f+i_coord_offset,fshift+i_shift_offset);
1026 dvdasum = _mm256_mul_ps(dvdasum, _mm256_mul_ps(isai0,isai0));
1027 gmx_mm256_update_1pot_ps(dvdasum,dvda+inr);
1029 /* Increment number of inner iterations */
1030 inneriter += j_index_end - j_index_start;
1032 /* Outer loop uses 7 flops */
1035 /* Increment number of outer iterations */
1038 /* Update outer/inner flops */
1040 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_F,outeriter*7 + inneriter*64);