2 * Note: this file was generated by the Gromacs avx_256_double kernel generator.
4 * This source code is part of
8 * Copyright (c) 2001-2012, The GROMACS Development Team
10 * Gromacs is a library for molecular simulation and trajectory analysis,
11 * written by Erik Lindahl, David van der Spoel, Berk Hess, and others - for
12 * a full list of developers and information, check out http://www.gromacs.org
14 * This program is free software; you can redistribute it and/or modify it under
15 * the terms of the GNU Lesser General Public License as published by the Free
16 * Software Foundation; either version 2 of the License, or (at your option) any
19 * To help fund GROMACS development, we humbly ask that you cite
20 * the papers people have written on it - you can find them on the website.
28 #include "../nb_kernel.h"
29 #include "types/simple.h"
33 #include "gmx_math_x86_avx_256_double.h"
34 #include "kernelutil_x86_avx_256_double.h"
37 * Gromacs nonbonded kernel: nb_kernel_ElecRFCut_VdwCSTab_GeomP1P1_VF_avx_256_double
38 * Electrostatics interaction: ReactionField
39 * VdW interaction: CubicSplineTable
40 * Geometry: Particle-Particle
41 * Calculate force/pot: PotentialAndForce
44 nb_kernel_ElecRFCut_VdwCSTab_GeomP1P1_VF_avx_256_double
45 (t_nblist * gmx_restrict nlist,
46 rvec * gmx_restrict xx,
47 rvec * gmx_restrict ff,
48 t_forcerec * gmx_restrict fr,
49 t_mdatoms * gmx_restrict mdatoms,
50 nb_kernel_data_t * gmx_restrict kernel_data,
51 t_nrnb * gmx_restrict nrnb)
53 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
54 * just 0 for non-waters.
55 * Suffixes A,B,C,D refer to j loop unrolling done with AVX, e.g. for the four different
56 * jnr indices corresponding to data put in the four positions in the SIMD register.
58 int i_shift_offset,i_coord_offset,outeriter,inneriter;
59 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
60 int jnrA,jnrB,jnrC,jnrD;
61 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
62 int jnrlistE,jnrlistF,jnrlistG,jnrlistH;
63 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
64 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
66 real *shiftvec,*fshift,*x,*f;
67 real *fjptrA,*fjptrB,*fjptrC,*fjptrD;
69 __m256d tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
70 real * vdwioffsetptr0;
71 __m256d ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
72 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
73 __m256d jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
74 __m256d dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
75 __m256d velec,felec,velecsum,facel,crf,krf,krf2;
78 __m256d rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
81 __m256d one_sixth = _mm256_set1_pd(1.0/6.0);
82 __m256d one_twelfth = _mm256_set1_pd(1.0/12.0);
84 __m128i ifour = _mm_set1_epi32(4);
85 __m256d rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF;
87 __m256d dummy_mask,cutoff_mask;
88 __m128 tmpmask0,tmpmask1;
89 __m256d signbit = _mm256_castsi256_pd( _mm256_set1_epi32(0x80000000) );
90 __m256d one = _mm256_set1_pd(1.0);
91 __m256d two = _mm256_set1_pd(2.0);
97 jindex = nlist->jindex;
99 shiftidx = nlist->shift;
101 shiftvec = fr->shift_vec[0];
102 fshift = fr->fshift[0];
103 facel = _mm256_set1_pd(fr->epsfac);
104 charge = mdatoms->chargeA;
105 krf = _mm256_set1_pd(fr->ic->k_rf);
106 krf2 = _mm256_set1_pd(fr->ic->k_rf*2.0);
107 crf = _mm256_set1_pd(fr->ic->c_rf);
108 nvdwtype = fr->ntype;
110 vdwtype = mdatoms->typeA;
112 vftab = kernel_data->table_vdw->data;
113 vftabscale = _mm256_set1_pd(kernel_data->table_vdw->scale);
115 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
116 rcutoff_scalar = fr->rcoulomb;
117 rcutoff = _mm256_set1_pd(rcutoff_scalar);
118 rcutoff2 = _mm256_mul_pd(rcutoff,rcutoff);
120 /* Avoid stupid compiler warnings */
121 jnrA = jnrB = jnrC = jnrD = 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_pd(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
152 fix0 = _mm256_setzero_pd();
153 fiy0 = _mm256_setzero_pd();
154 fiz0 = _mm256_setzero_pd();
156 /* Load parameters for i particles */
157 iq0 = _mm256_mul_pd(facel,_mm256_set1_pd(charge[inr+0]));
158 vdwioffsetptr0 = vdwparam+2*nvdwtype*vdwtype[inr+0];
160 /* Reset potential sums */
161 velecsum = _mm256_setzero_pd();
162 vvdwsum = _mm256_setzero_pd();
164 /* Start inner kernel loop */
165 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
168 /* Get j neighbor index, and coordinate index */
173 j_coord_offsetA = DIM*jnrA;
174 j_coord_offsetB = DIM*jnrB;
175 j_coord_offsetC = DIM*jnrC;
176 j_coord_offsetD = DIM*jnrD;
178 /* load j atom coordinates */
179 gmx_mm256_load_1rvec_4ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
180 x+j_coord_offsetC,x+j_coord_offsetD,
183 /* Calculate displacement vector */
184 dx00 = _mm256_sub_pd(ix0,jx0);
185 dy00 = _mm256_sub_pd(iy0,jy0);
186 dz00 = _mm256_sub_pd(iz0,jz0);
188 /* Calculate squared distance and things based on it */
189 rsq00 = gmx_mm256_calc_rsq_pd(dx00,dy00,dz00);
191 rinv00 = gmx_mm256_invsqrt_pd(rsq00);
193 rinvsq00 = _mm256_mul_pd(rinv00,rinv00);
195 /* Load parameters for j particles */
196 jq0 = gmx_mm256_load_4real_swizzle_pd(charge+jnrA+0,charge+jnrB+0,
197 charge+jnrC+0,charge+jnrD+0);
198 vdwjidx0A = 2*vdwtype[jnrA+0];
199 vdwjidx0B = 2*vdwtype[jnrB+0];
200 vdwjidx0C = 2*vdwtype[jnrC+0];
201 vdwjidx0D = 2*vdwtype[jnrD+0];
203 /**************************
204 * CALCULATE INTERACTIONS *
205 **************************/
207 if (gmx_mm256_any_lt(rsq00,rcutoff2))
210 r00 = _mm256_mul_pd(rsq00,rinv00);
212 /* Compute parameters for interactions between i and j atoms */
213 qq00 = _mm256_mul_pd(iq0,jq0);
214 gmx_mm256_load_4pair_swizzle_pd(vdwioffsetptr0+vdwjidx0A,
215 vdwioffsetptr0+vdwjidx0B,
216 vdwioffsetptr0+vdwjidx0C,
217 vdwioffsetptr0+vdwjidx0D,
220 /* Calculate table index by multiplying r with table scale and truncate to integer */
221 rt = _mm256_mul_pd(r00,vftabscale);
222 vfitab = _mm256_cvttpd_epi32(rt);
223 vfeps = _mm256_sub_pd(rt,_mm256_round_pd(rt, _MM_FROUND_FLOOR));
224 vfitab = _mm_slli_epi32(vfitab,3);
226 /* REACTION-FIELD ELECTROSTATICS */
227 velec = _mm256_mul_pd(qq00,_mm256_sub_pd(_mm256_add_pd(rinv00,_mm256_mul_pd(krf,rsq00)),crf));
228 felec = _mm256_mul_pd(qq00,_mm256_sub_pd(_mm256_mul_pd(rinv00,rinvsq00),krf2));
230 /* CUBIC SPLINE TABLE DISPERSION */
231 Y = _mm256_load_pd( vftab + _mm_extract_epi32(vfitab,0) );
232 F = _mm256_load_pd( vftab + _mm_extract_epi32(vfitab,1) );
233 G = _mm256_load_pd( vftab + _mm_extract_epi32(vfitab,2) );
234 H = _mm256_load_pd( vftab + _mm_extract_epi32(vfitab,3) );
235 GMX_MM256_FULLTRANSPOSE4_PD(Y,F,G,H);
236 Heps = _mm256_mul_pd(vfeps,H);
237 Fp = _mm256_add_pd(F,_mm256_mul_pd(vfeps,_mm256_add_pd(G,Heps)));
238 VV = _mm256_add_pd(Y,_mm256_mul_pd(vfeps,Fp));
239 vvdw6 = _mm256_mul_pd(c6_00,VV);
240 FF = _mm256_add_pd(Fp,_mm256_mul_pd(vfeps,_mm256_add_pd(G,_mm256_add_pd(Heps,Heps))));
241 fvdw6 = _mm256_mul_pd(c6_00,FF);
243 /* CUBIC SPLINE TABLE REPULSION */
244 vfitab = _mm_add_epi32(vfitab,ifour);
245 Y = _mm256_load_pd( vftab + _mm_extract_epi32(vfitab,0) );
246 F = _mm256_load_pd( vftab + _mm_extract_epi32(vfitab,1) );
247 G = _mm256_load_pd( vftab + _mm_extract_epi32(vfitab,2) );
248 H = _mm256_load_pd( vftab + _mm_extract_epi32(vfitab,3) );
249 GMX_MM256_FULLTRANSPOSE4_PD(Y,F,G,H);
250 Heps = _mm256_mul_pd(vfeps,H);
251 Fp = _mm256_add_pd(F,_mm256_mul_pd(vfeps,_mm256_add_pd(G,Heps)));
252 VV = _mm256_add_pd(Y,_mm256_mul_pd(vfeps,Fp));
253 vvdw12 = _mm256_mul_pd(c12_00,VV);
254 FF = _mm256_add_pd(Fp,_mm256_mul_pd(vfeps,_mm256_add_pd(G,_mm256_add_pd(Heps,Heps))));
255 fvdw12 = _mm256_mul_pd(c12_00,FF);
256 vvdw = _mm256_add_pd(vvdw12,vvdw6);
257 fvdw = _mm256_xor_pd(signbit,_mm256_mul_pd(_mm256_add_pd(fvdw6,fvdw12),_mm256_mul_pd(vftabscale,rinv00)));
259 cutoff_mask = _mm256_cmp_pd(rsq00,rcutoff2,_CMP_LT_OQ);
261 /* Update potential sum for this i atom from the interaction with this j atom. */
262 velec = _mm256_and_pd(velec,cutoff_mask);
263 velecsum = _mm256_add_pd(velecsum,velec);
264 vvdw = _mm256_and_pd(vvdw,cutoff_mask);
265 vvdwsum = _mm256_add_pd(vvdwsum,vvdw);
267 fscal = _mm256_add_pd(felec,fvdw);
269 fscal = _mm256_and_pd(fscal,cutoff_mask);
271 /* Calculate temporary vectorial force */
272 tx = _mm256_mul_pd(fscal,dx00);
273 ty = _mm256_mul_pd(fscal,dy00);
274 tz = _mm256_mul_pd(fscal,dz00);
276 /* Update vectorial force */
277 fix0 = _mm256_add_pd(fix0,tx);
278 fiy0 = _mm256_add_pd(fiy0,ty);
279 fiz0 = _mm256_add_pd(fiz0,tz);
281 fjptrA = f+j_coord_offsetA;
282 fjptrB = f+j_coord_offsetB;
283 fjptrC = f+j_coord_offsetC;
284 fjptrD = f+j_coord_offsetD;
285 gmx_mm256_decrement_1rvec_4ptr_swizzle_pd(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
289 /* Inner loop uses 72 flops */
295 /* Get j neighbor index, and coordinate index */
296 jnrlistA = jjnr[jidx];
297 jnrlistB = jjnr[jidx+1];
298 jnrlistC = jjnr[jidx+2];
299 jnrlistD = jjnr[jidx+3];
300 /* Sign of each element will be negative for non-real atoms.
301 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
302 * so use it as val = _mm_andnot_pd(mask,val) to clear dummy entries.
304 tmpmask0 = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
306 tmpmask1 = _mm_permute_ps(tmpmask0,_GMX_MM_PERMUTE(3,3,2,2));
307 tmpmask0 = _mm_permute_ps(tmpmask0,_GMX_MM_PERMUTE(1,1,0,0));
308 dummy_mask = _mm256_castps_pd(gmx_mm256_set_m128(tmpmask1,tmpmask0));
310 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
311 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
312 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
313 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
314 j_coord_offsetA = DIM*jnrA;
315 j_coord_offsetB = DIM*jnrB;
316 j_coord_offsetC = DIM*jnrC;
317 j_coord_offsetD = DIM*jnrD;
319 /* load j atom coordinates */
320 gmx_mm256_load_1rvec_4ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
321 x+j_coord_offsetC,x+j_coord_offsetD,
324 /* Calculate displacement vector */
325 dx00 = _mm256_sub_pd(ix0,jx0);
326 dy00 = _mm256_sub_pd(iy0,jy0);
327 dz00 = _mm256_sub_pd(iz0,jz0);
329 /* Calculate squared distance and things based on it */
330 rsq00 = gmx_mm256_calc_rsq_pd(dx00,dy00,dz00);
332 rinv00 = gmx_mm256_invsqrt_pd(rsq00);
334 rinvsq00 = _mm256_mul_pd(rinv00,rinv00);
336 /* Load parameters for j particles */
337 jq0 = gmx_mm256_load_4real_swizzle_pd(charge+jnrA+0,charge+jnrB+0,
338 charge+jnrC+0,charge+jnrD+0);
339 vdwjidx0A = 2*vdwtype[jnrA+0];
340 vdwjidx0B = 2*vdwtype[jnrB+0];
341 vdwjidx0C = 2*vdwtype[jnrC+0];
342 vdwjidx0D = 2*vdwtype[jnrD+0];
344 /**************************
345 * CALCULATE INTERACTIONS *
346 **************************/
348 if (gmx_mm256_any_lt(rsq00,rcutoff2))
351 r00 = _mm256_mul_pd(rsq00,rinv00);
352 r00 = _mm256_andnot_pd(dummy_mask,r00);
354 /* Compute parameters for interactions between i and j atoms */
355 qq00 = _mm256_mul_pd(iq0,jq0);
356 gmx_mm256_load_4pair_swizzle_pd(vdwioffsetptr0+vdwjidx0A,
357 vdwioffsetptr0+vdwjidx0B,
358 vdwioffsetptr0+vdwjidx0C,
359 vdwioffsetptr0+vdwjidx0D,
362 /* Calculate table index by multiplying r with table scale and truncate to integer */
363 rt = _mm256_mul_pd(r00,vftabscale);
364 vfitab = _mm256_cvttpd_epi32(rt);
365 vfeps = _mm256_sub_pd(rt,_mm256_round_pd(rt, _MM_FROUND_FLOOR));
366 vfitab = _mm_slli_epi32(vfitab,3);
368 /* REACTION-FIELD ELECTROSTATICS */
369 velec = _mm256_mul_pd(qq00,_mm256_sub_pd(_mm256_add_pd(rinv00,_mm256_mul_pd(krf,rsq00)),crf));
370 felec = _mm256_mul_pd(qq00,_mm256_sub_pd(_mm256_mul_pd(rinv00,rinvsq00),krf2));
372 /* CUBIC SPLINE TABLE DISPERSION */
373 Y = _mm256_load_pd( vftab + _mm_extract_epi32(vfitab,0) );
374 F = _mm256_load_pd( vftab + _mm_extract_epi32(vfitab,1) );
375 G = _mm256_load_pd( vftab + _mm_extract_epi32(vfitab,2) );
376 H = _mm256_load_pd( vftab + _mm_extract_epi32(vfitab,3) );
377 GMX_MM256_FULLTRANSPOSE4_PD(Y,F,G,H);
378 Heps = _mm256_mul_pd(vfeps,H);
379 Fp = _mm256_add_pd(F,_mm256_mul_pd(vfeps,_mm256_add_pd(G,Heps)));
380 VV = _mm256_add_pd(Y,_mm256_mul_pd(vfeps,Fp));
381 vvdw6 = _mm256_mul_pd(c6_00,VV);
382 FF = _mm256_add_pd(Fp,_mm256_mul_pd(vfeps,_mm256_add_pd(G,_mm256_add_pd(Heps,Heps))));
383 fvdw6 = _mm256_mul_pd(c6_00,FF);
385 /* CUBIC SPLINE TABLE REPULSION */
386 vfitab = _mm_add_epi32(vfitab,ifour);
387 Y = _mm256_load_pd( vftab + _mm_extract_epi32(vfitab,0) );
388 F = _mm256_load_pd( vftab + _mm_extract_epi32(vfitab,1) );
389 G = _mm256_load_pd( vftab + _mm_extract_epi32(vfitab,2) );
390 H = _mm256_load_pd( vftab + _mm_extract_epi32(vfitab,3) );
391 GMX_MM256_FULLTRANSPOSE4_PD(Y,F,G,H);
392 Heps = _mm256_mul_pd(vfeps,H);
393 Fp = _mm256_add_pd(F,_mm256_mul_pd(vfeps,_mm256_add_pd(G,Heps)));
394 VV = _mm256_add_pd(Y,_mm256_mul_pd(vfeps,Fp));
395 vvdw12 = _mm256_mul_pd(c12_00,VV);
396 FF = _mm256_add_pd(Fp,_mm256_mul_pd(vfeps,_mm256_add_pd(G,_mm256_add_pd(Heps,Heps))));
397 fvdw12 = _mm256_mul_pd(c12_00,FF);
398 vvdw = _mm256_add_pd(vvdw12,vvdw6);
399 fvdw = _mm256_xor_pd(signbit,_mm256_mul_pd(_mm256_add_pd(fvdw6,fvdw12),_mm256_mul_pd(vftabscale,rinv00)));
401 cutoff_mask = _mm256_cmp_pd(rsq00,rcutoff2,_CMP_LT_OQ);
403 /* Update potential sum for this i atom from the interaction with this j atom. */
404 velec = _mm256_and_pd(velec,cutoff_mask);
405 velec = _mm256_andnot_pd(dummy_mask,velec);
406 velecsum = _mm256_add_pd(velecsum,velec);
407 vvdw = _mm256_and_pd(vvdw,cutoff_mask);
408 vvdw = _mm256_andnot_pd(dummy_mask,vvdw);
409 vvdwsum = _mm256_add_pd(vvdwsum,vvdw);
411 fscal = _mm256_add_pd(felec,fvdw);
413 fscal = _mm256_and_pd(fscal,cutoff_mask);
415 fscal = _mm256_andnot_pd(dummy_mask,fscal);
417 /* Calculate temporary vectorial force */
418 tx = _mm256_mul_pd(fscal,dx00);
419 ty = _mm256_mul_pd(fscal,dy00);
420 tz = _mm256_mul_pd(fscal,dz00);
422 /* Update vectorial force */
423 fix0 = _mm256_add_pd(fix0,tx);
424 fiy0 = _mm256_add_pd(fiy0,ty);
425 fiz0 = _mm256_add_pd(fiz0,tz);
427 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
428 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
429 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
430 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
431 gmx_mm256_decrement_1rvec_4ptr_swizzle_pd(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
435 /* Inner loop uses 73 flops */
438 /* End of innermost loop */
440 gmx_mm256_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
441 f+i_coord_offset,fshift+i_shift_offset);
444 /* Update potential energies */
445 gmx_mm256_update_1pot_pd(velecsum,kernel_data->energygrp_elec+ggid);
446 gmx_mm256_update_1pot_pd(vvdwsum,kernel_data->energygrp_vdw+ggid);
448 /* Increment number of inner iterations */
449 inneriter += j_index_end - j_index_start;
451 /* Outer loop uses 9 flops */
454 /* Increment number of outer iterations */
457 /* Update outer/inner flops */
459 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_VF,outeriter*9 + inneriter*73);
462 * Gromacs nonbonded kernel: nb_kernel_ElecRFCut_VdwCSTab_GeomP1P1_F_avx_256_double
463 * Electrostatics interaction: ReactionField
464 * VdW interaction: CubicSplineTable
465 * Geometry: Particle-Particle
466 * Calculate force/pot: Force
469 nb_kernel_ElecRFCut_VdwCSTab_GeomP1P1_F_avx_256_double
470 (t_nblist * gmx_restrict nlist,
471 rvec * gmx_restrict xx,
472 rvec * gmx_restrict ff,
473 t_forcerec * gmx_restrict fr,
474 t_mdatoms * gmx_restrict mdatoms,
475 nb_kernel_data_t * gmx_restrict kernel_data,
476 t_nrnb * gmx_restrict nrnb)
478 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
479 * just 0 for non-waters.
480 * Suffixes A,B,C,D refer to j loop unrolling done with AVX, e.g. for the four different
481 * jnr indices corresponding to data put in the four positions in the SIMD register.
483 int i_shift_offset,i_coord_offset,outeriter,inneriter;
484 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
485 int jnrA,jnrB,jnrC,jnrD;
486 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
487 int jnrlistE,jnrlistF,jnrlistG,jnrlistH;
488 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
489 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
491 real *shiftvec,*fshift,*x,*f;
492 real *fjptrA,*fjptrB,*fjptrC,*fjptrD;
494 __m256d tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
495 real * vdwioffsetptr0;
496 __m256d ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
497 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
498 __m256d jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
499 __m256d dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
500 __m256d velec,felec,velecsum,facel,crf,krf,krf2;
503 __m256d rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
506 __m256d one_sixth = _mm256_set1_pd(1.0/6.0);
507 __m256d one_twelfth = _mm256_set1_pd(1.0/12.0);
509 __m128i ifour = _mm_set1_epi32(4);
510 __m256d rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF;
512 __m256d dummy_mask,cutoff_mask;
513 __m128 tmpmask0,tmpmask1;
514 __m256d signbit = _mm256_castsi256_pd( _mm256_set1_epi32(0x80000000) );
515 __m256d one = _mm256_set1_pd(1.0);
516 __m256d two = _mm256_set1_pd(2.0);
522 jindex = nlist->jindex;
524 shiftidx = nlist->shift;
526 shiftvec = fr->shift_vec[0];
527 fshift = fr->fshift[0];
528 facel = _mm256_set1_pd(fr->epsfac);
529 charge = mdatoms->chargeA;
530 krf = _mm256_set1_pd(fr->ic->k_rf);
531 krf2 = _mm256_set1_pd(fr->ic->k_rf*2.0);
532 crf = _mm256_set1_pd(fr->ic->c_rf);
533 nvdwtype = fr->ntype;
535 vdwtype = mdatoms->typeA;
537 vftab = kernel_data->table_vdw->data;
538 vftabscale = _mm256_set1_pd(kernel_data->table_vdw->scale);
540 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
541 rcutoff_scalar = fr->rcoulomb;
542 rcutoff = _mm256_set1_pd(rcutoff_scalar);
543 rcutoff2 = _mm256_mul_pd(rcutoff,rcutoff);
545 /* Avoid stupid compiler warnings */
546 jnrA = jnrB = jnrC = jnrD = 0;
555 for(iidx=0;iidx<4*DIM;iidx++)
560 /* Start outer loop over neighborlists */
561 for(iidx=0; iidx<nri; iidx++)
563 /* Load shift vector for this list */
564 i_shift_offset = DIM*shiftidx[iidx];
566 /* Load limits for loop over neighbors */
567 j_index_start = jindex[iidx];
568 j_index_end = jindex[iidx+1];
570 /* Get outer coordinate index */
572 i_coord_offset = DIM*inr;
574 /* Load i particle coords and add shift vector */
575 gmx_mm256_load_shift_and_1rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
577 fix0 = _mm256_setzero_pd();
578 fiy0 = _mm256_setzero_pd();
579 fiz0 = _mm256_setzero_pd();
581 /* Load parameters for i particles */
582 iq0 = _mm256_mul_pd(facel,_mm256_set1_pd(charge[inr+0]));
583 vdwioffsetptr0 = vdwparam+2*nvdwtype*vdwtype[inr+0];
585 /* Start inner kernel loop */
586 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
589 /* Get j neighbor index, and coordinate index */
594 j_coord_offsetA = DIM*jnrA;
595 j_coord_offsetB = DIM*jnrB;
596 j_coord_offsetC = DIM*jnrC;
597 j_coord_offsetD = DIM*jnrD;
599 /* load j atom coordinates */
600 gmx_mm256_load_1rvec_4ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
601 x+j_coord_offsetC,x+j_coord_offsetD,
604 /* Calculate displacement vector */
605 dx00 = _mm256_sub_pd(ix0,jx0);
606 dy00 = _mm256_sub_pd(iy0,jy0);
607 dz00 = _mm256_sub_pd(iz0,jz0);
609 /* Calculate squared distance and things based on it */
610 rsq00 = gmx_mm256_calc_rsq_pd(dx00,dy00,dz00);
612 rinv00 = gmx_mm256_invsqrt_pd(rsq00);
614 rinvsq00 = _mm256_mul_pd(rinv00,rinv00);
616 /* Load parameters for j particles */
617 jq0 = gmx_mm256_load_4real_swizzle_pd(charge+jnrA+0,charge+jnrB+0,
618 charge+jnrC+0,charge+jnrD+0);
619 vdwjidx0A = 2*vdwtype[jnrA+0];
620 vdwjidx0B = 2*vdwtype[jnrB+0];
621 vdwjidx0C = 2*vdwtype[jnrC+0];
622 vdwjidx0D = 2*vdwtype[jnrD+0];
624 /**************************
625 * CALCULATE INTERACTIONS *
626 **************************/
628 if (gmx_mm256_any_lt(rsq00,rcutoff2))
631 r00 = _mm256_mul_pd(rsq00,rinv00);
633 /* Compute parameters for interactions between i and j atoms */
634 qq00 = _mm256_mul_pd(iq0,jq0);
635 gmx_mm256_load_4pair_swizzle_pd(vdwioffsetptr0+vdwjidx0A,
636 vdwioffsetptr0+vdwjidx0B,
637 vdwioffsetptr0+vdwjidx0C,
638 vdwioffsetptr0+vdwjidx0D,
641 /* Calculate table index by multiplying r with table scale and truncate to integer */
642 rt = _mm256_mul_pd(r00,vftabscale);
643 vfitab = _mm256_cvttpd_epi32(rt);
644 vfeps = _mm256_sub_pd(rt,_mm256_round_pd(rt, _MM_FROUND_FLOOR));
645 vfitab = _mm_slli_epi32(vfitab,3);
647 /* REACTION-FIELD ELECTROSTATICS */
648 felec = _mm256_mul_pd(qq00,_mm256_sub_pd(_mm256_mul_pd(rinv00,rinvsq00),krf2));
650 /* CUBIC SPLINE TABLE DISPERSION */
651 Y = _mm256_load_pd( vftab + _mm_extract_epi32(vfitab,0) );
652 F = _mm256_load_pd( vftab + _mm_extract_epi32(vfitab,1) );
653 G = _mm256_load_pd( vftab + _mm_extract_epi32(vfitab,2) );
654 H = _mm256_load_pd( vftab + _mm_extract_epi32(vfitab,3) );
655 GMX_MM256_FULLTRANSPOSE4_PD(Y,F,G,H);
656 Heps = _mm256_mul_pd(vfeps,H);
657 Fp = _mm256_add_pd(F,_mm256_mul_pd(vfeps,_mm256_add_pd(G,Heps)));
658 FF = _mm256_add_pd(Fp,_mm256_mul_pd(vfeps,_mm256_add_pd(G,_mm256_add_pd(Heps,Heps))));
659 fvdw6 = _mm256_mul_pd(c6_00,FF);
661 /* CUBIC SPLINE TABLE REPULSION */
662 vfitab = _mm_add_epi32(vfitab,ifour);
663 Y = _mm256_load_pd( vftab + _mm_extract_epi32(vfitab,0) );
664 F = _mm256_load_pd( vftab + _mm_extract_epi32(vfitab,1) );
665 G = _mm256_load_pd( vftab + _mm_extract_epi32(vfitab,2) );
666 H = _mm256_load_pd( vftab + _mm_extract_epi32(vfitab,3) );
667 GMX_MM256_FULLTRANSPOSE4_PD(Y,F,G,H);
668 Heps = _mm256_mul_pd(vfeps,H);
669 Fp = _mm256_add_pd(F,_mm256_mul_pd(vfeps,_mm256_add_pd(G,Heps)));
670 FF = _mm256_add_pd(Fp,_mm256_mul_pd(vfeps,_mm256_add_pd(G,_mm256_add_pd(Heps,Heps))));
671 fvdw12 = _mm256_mul_pd(c12_00,FF);
672 fvdw = _mm256_xor_pd(signbit,_mm256_mul_pd(_mm256_add_pd(fvdw6,fvdw12),_mm256_mul_pd(vftabscale,rinv00)));
674 cutoff_mask = _mm256_cmp_pd(rsq00,rcutoff2,_CMP_LT_OQ);
676 fscal = _mm256_add_pd(felec,fvdw);
678 fscal = _mm256_and_pd(fscal,cutoff_mask);
680 /* Calculate temporary vectorial force */
681 tx = _mm256_mul_pd(fscal,dx00);
682 ty = _mm256_mul_pd(fscal,dy00);
683 tz = _mm256_mul_pd(fscal,dz00);
685 /* Update vectorial force */
686 fix0 = _mm256_add_pd(fix0,tx);
687 fiy0 = _mm256_add_pd(fiy0,ty);
688 fiz0 = _mm256_add_pd(fiz0,tz);
690 fjptrA = f+j_coord_offsetA;
691 fjptrB = f+j_coord_offsetB;
692 fjptrC = f+j_coord_offsetC;
693 fjptrD = f+j_coord_offsetD;
694 gmx_mm256_decrement_1rvec_4ptr_swizzle_pd(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
698 /* Inner loop uses 57 flops */
704 /* Get j neighbor index, and coordinate index */
705 jnrlistA = jjnr[jidx];
706 jnrlistB = jjnr[jidx+1];
707 jnrlistC = jjnr[jidx+2];
708 jnrlistD = jjnr[jidx+3];
709 /* Sign of each element will be negative for non-real atoms.
710 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
711 * so use it as val = _mm_andnot_pd(mask,val) to clear dummy entries.
713 tmpmask0 = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
715 tmpmask1 = _mm_permute_ps(tmpmask0,_GMX_MM_PERMUTE(3,3,2,2));
716 tmpmask0 = _mm_permute_ps(tmpmask0,_GMX_MM_PERMUTE(1,1,0,0));
717 dummy_mask = _mm256_castps_pd(gmx_mm256_set_m128(tmpmask1,tmpmask0));
719 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
720 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
721 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
722 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
723 j_coord_offsetA = DIM*jnrA;
724 j_coord_offsetB = DIM*jnrB;
725 j_coord_offsetC = DIM*jnrC;
726 j_coord_offsetD = DIM*jnrD;
728 /* load j atom coordinates */
729 gmx_mm256_load_1rvec_4ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
730 x+j_coord_offsetC,x+j_coord_offsetD,
733 /* Calculate displacement vector */
734 dx00 = _mm256_sub_pd(ix0,jx0);
735 dy00 = _mm256_sub_pd(iy0,jy0);
736 dz00 = _mm256_sub_pd(iz0,jz0);
738 /* Calculate squared distance and things based on it */
739 rsq00 = gmx_mm256_calc_rsq_pd(dx00,dy00,dz00);
741 rinv00 = gmx_mm256_invsqrt_pd(rsq00);
743 rinvsq00 = _mm256_mul_pd(rinv00,rinv00);
745 /* Load parameters for j particles */
746 jq0 = gmx_mm256_load_4real_swizzle_pd(charge+jnrA+0,charge+jnrB+0,
747 charge+jnrC+0,charge+jnrD+0);
748 vdwjidx0A = 2*vdwtype[jnrA+0];
749 vdwjidx0B = 2*vdwtype[jnrB+0];
750 vdwjidx0C = 2*vdwtype[jnrC+0];
751 vdwjidx0D = 2*vdwtype[jnrD+0];
753 /**************************
754 * CALCULATE INTERACTIONS *
755 **************************/
757 if (gmx_mm256_any_lt(rsq00,rcutoff2))
760 r00 = _mm256_mul_pd(rsq00,rinv00);
761 r00 = _mm256_andnot_pd(dummy_mask,r00);
763 /* Compute parameters for interactions between i and j atoms */
764 qq00 = _mm256_mul_pd(iq0,jq0);
765 gmx_mm256_load_4pair_swizzle_pd(vdwioffsetptr0+vdwjidx0A,
766 vdwioffsetptr0+vdwjidx0B,
767 vdwioffsetptr0+vdwjidx0C,
768 vdwioffsetptr0+vdwjidx0D,
771 /* Calculate table index by multiplying r with table scale and truncate to integer */
772 rt = _mm256_mul_pd(r00,vftabscale);
773 vfitab = _mm256_cvttpd_epi32(rt);
774 vfeps = _mm256_sub_pd(rt,_mm256_round_pd(rt, _MM_FROUND_FLOOR));
775 vfitab = _mm_slli_epi32(vfitab,3);
777 /* REACTION-FIELD ELECTROSTATICS */
778 felec = _mm256_mul_pd(qq00,_mm256_sub_pd(_mm256_mul_pd(rinv00,rinvsq00),krf2));
780 /* CUBIC SPLINE TABLE DISPERSION */
781 Y = _mm256_load_pd( vftab + _mm_extract_epi32(vfitab,0) );
782 F = _mm256_load_pd( vftab + _mm_extract_epi32(vfitab,1) );
783 G = _mm256_load_pd( vftab + _mm_extract_epi32(vfitab,2) );
784 H = _mm256_load_pd( vftab + _mm_extract_epi32(vfitab,3) );
785 GMX_MM256_FULLTRANSPOSE4_PD(Y,F,G,H);
786 Heps = _mm256_mul_pd(vfeps,H);
787 Fp = _mm256_add_pd(F,_mm256_mul_pd(vfeps,_mm256_add_pd(G,Heps)));
788 FF = _mm256_add_pd(Fp,_mm256_mul_pd(vfeps,_mm256_add_pd(G,_mm256_add_pd(Heps,Heps))));
789 fvdw6 = _mm256_mul_pd(c6_00,FF);
791 /* CUBIC SPLINE TABLE REPULSION */
792 vfitab = _mm_add_epi32(vfitab,ifour);
793 Y = _mm256_load_pd( vftab + _mm_extract_epi32(vfitab,0) );
794 F = _mm256_load_pd( vftab + _mm_extract_epi32(vfitab,1) );
795 G = _mm256_load_pd( vftab + _mm_extract_epi32(vfitab,2) );
796 H = _mm256_load_pd( vftab + _mm_extract_epi32(vfitab,3) );
797 GMX_MM256_FULLTRANSPOSE4_PD(Y,F,G,H);
798 Heps = _mm256_mul_pd(vfeps,H);
799 Fp = _mm256_add_pd(F,_mm256_mul_pd(vfeps,_mm256_add_pd(G,Heps)));
800 FF = _mm256_add_pd(Fp,_mm256_mul_pd(vfeps,_mm256_add_pd(G,_mm256_add_pd(Heps,Heps))));
801 fvdw12 = _mm256_mul_pd(c12_00,FF);
802 fvdw = _mm256_xor_pd(signbit,_mm256_mul_pd(_mm256_add_pd(fvdw6,fvdw12),_mm256_mul_pd(vftabscale,rinv00)));
804 cutoff_mask = _mm256_cmp_pd(rsq00,rcutoff2,_CMP_LT_OQ);
806 fscal = _mm256_add_pd(felec,fvdw);
808 fscal = _mm256_and_pd(fscal,cutoff_mask);
810 fscal = _mm256_andnot_pd(dummy_mask,fscal);
812 /* Calculate temporary vectorial force */
813 tx = _mm256_mul_pd(fscal,dx00);
814 ty = _mm256_mul_pd(fscal,dy00);
815 tz = _mm256_mul_pd(fscal,dz00);
817 /* Update vectorial force */
818 fix0 = _mm256_add_pd(fix0,tx);
819 fiy0 = _mm256_add_pd(fiy0,ty);
820 fiz0 = _mm256_add_pd(fiz0,tz);
822 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
823 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
824 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
825 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
826 gmx_mm256_decrement_1rvec_4ptr_swizzle_pd(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
830 /* Inner loop uses 58 flops */
833 /* End of innermost loop */
835 gmx_mm256_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
836 f+i_coord_offset,fshift+i_shift_offset);
838 /* Increment number of inner iterations */
839 inneriter += j_index_end - j_index_start;
841 /* Outer loop uses 7 flops */
844 /* Increment number of outer iterations */
847 /* Update outer/inner flops */
849 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_F,outeriter*7 + inneriter*58);