2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014, 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_double kernel generator.
42 #include "../nb_kernel.h"
43 #include "gromacs/legacyheaders/types/simple.h"
44 #include "gromacs/math/vec.h"
45 #include "gromacs/legacyheaders/nrnb.h"
47 #include "gromacs/simd/math_x86_avx_256_double.h"
48 #include "kernelutil_x86_avx_256_double.h"
51 * Gromacs nonbonded kernel: nb_kernel_ElecNone_VdwLJEwSh_GeomP1P1_VF_avx_256_double
52 * Electrostatics interaction: None
53 * VdW interaction: LJEwald
54 * Geometry: Particle-Particle
55 * Calculate force/pot: PotentialAndForce
58 nb_kernel_ElecNone_VdwLJEwSh_GeomP1P1_VF_avx_256_double
59 (t_nblist * gmx_restrict nlist,
60 rvec * gmx_restrict xx,
61 rvec * gmx_restrict ff,
62 t_forcerec * gmx_restrict fr,
63 t_mdatoms * gmx_restrict mdatoms,
64 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
65 t_nrnb * gmx_restrict nrnb)
67 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
68 * just 0 for non-waters.
69 * Suffixes A,B,C,D refer to j loop unrolling done with AVX, e.g. for the four different
70 * jnr indices corresponding to data put in the four positions in the SIMD register.
72 int i_shift_offset,i_coord_offset,outeriter,inneriter;
73 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
74 int jnrA,jnrB,jnrC,jnrD;
75 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
76 int jnrlistE,jnrlistF,jnrlistG,jnrlistH;
77 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
78 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
80 real *shiftvec,*fshift,*x,*f;
81 real *fjptrA,*fjptrB,*fjptrC,*fjptrD;
83 __m256d tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
84 real * vdwioffsetptr0;
85 real * vdwgridioffsetptr0;
86 __m256d ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
87 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
88 __m256d jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
89 __m256d dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
91 __m256d rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
94 __m256d one_sixth = _mm256_set1_pd(1.0/6.0);
95 __m256d one_twelfth = _mm256_set1_pd(1.0/12.0);
98 __m256d ewclj,ewclj2,ewclj6,ewcljrsq,poly,exponent,f6A,f6B,sh_lj_ewald;
99 __m256d one_half = _mm256_set1_pd(0.5);
100 __m256d minus_one = _mm256_set1_pd(-1.0);
101 __m256d dummy_mask,cutoff_mask;
102 __m128 tmpmask0,tmpmask1;
103 __m256d signbit = _mm256_castsi256_pd( _mm256_set1_epi32(0x80000000) );
104 __m256d one = _mm256_set1_pd(1.0);
105 __m256d two = _mm256_set1_pd(2.0);
111 jindex = nlist->jindex;
113 shiftidx = nlist->shift;
115 shiftvec = fr->shift_vec[0];
116 fshift = fr->fshift[0];
117 nvdwtype = fr->ntype;
119 vdwtype = mdatoms->typeA;
120 vdwgridparam = fr->ljpme_c6grid;
121 sh_lj_ewald = _mm256_set1_pd(fr->ic->sh_lj_ewald);
122 ewclj = _mm256_set1_pd(fr->ewaldcoeff_lj);
123 ewclj2 = _mm256_mul_pd(minus_one,_mm256_mul_pd(ewclj,ewclj));
125 rcutoff_scalar = fr->rvdw;
126 rcutoff = _mm256_set1_pd(rcutoff_scalar);
127 rcutoff2 = _mm256_mul_pd(rcutoff,rcutoff);
129 sh_vdw_invrcut6 = _mm256_set1_pd(fr->ic->sh_invrc6);
130 rvdw = _mm256_set1_pd(fr->rvdw);
132 /* Avoid stupid compiler warnings */
133 jnrA = jnrB = jnrC = jnrD = 0;
142 for(iidx=0;iidx<4*DIM;iidx++)
147 /* Start outer loop over neighborlists */
148 for(iidx=0; iidx<nri; iidx++)
150 /* Load shift vector for this list */
151 i_shift_offset = DIM*shiftidx[iidx];
153 /* Load limits for loop over neighbors */
154 j_index_start = jindex[iidx];
155 j_index_end = jindex[iidx+1];
157 /* Get outer coordinate index */
159 i_coord_offset = DIM*inr;
161 /* Load i particle coords and add shift vector */
162 gmx_mm256_load_shift_and_1rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
164 fix0 = _mm256_setzero_pd();
165 fiy0 = _mm256_setzero_pd();
166 fiz0 = _mm256_setzero_pd();
168 /* Load parameters for i particles */
169 vdwioffsetptr0 = vdwparam+2*nvdwtype*vdwtype[inr+0];
170 vdwgridioffsetptr0 = vdwgridparam+2*nvdwtype*vdwtype[inr+0];
172 /* Reset potential sums */
173 vvdwsum = _mm256_setzero_pd();
175 /* Start inner kernel loop */
176 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
179 /* Get j neighbor index, and coordinate index */
184 j_coord_offsetA = DIM*jnrA;
185 j_coord_offsetB = DIM*jnrB;
186 j_coord_offsetC = DIM*jnrC;
187 j_coord_offsetD = DIM*jnrD;
189 /* load j atom coordinates */
190 gmx_mm256_load_1rvec_4ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
191 x+j_coord_offsetC,x+j_coord_offsetD,
194 /* Calculate displacement vector */
195 dx00 = _mm256_sub_pd(ix0,jx0);
196 dy00 = _mm256_sub_pd(iy0,jy0);
197 dz00 = _mm256_sub_pd(iz0,jz0);
199 /* Calculate squared distance and things based on it */
200 rsq00 = gmx_mm256_calc_rsq_pd(dx00,dy00,dz00);
202 rinv00 = gmx_mm256_invsqrt_pd(rsq00);
204 rinvsq00 = _mm256_mul_pd(rinv00,rinv00);
206 /* Load parameters for j particles */
207 vdwjidx0A = 2*vdwtype[jnrA+0];
208 vdwjidx0B = 2*vdwtype[jnrB+0];
209 vdwjidx0C = 2*vdwtype[jnrC+0];
210 vdwjidx0D = 2*vdwtype[jnrD+0];
212 /**************************
213 * CALCULATE INTERACTIONS *
214 **************************/
216 if (gmx_mm256_any_lt(rsq00,rcutoff2))
219 r00 = _mm256_mul_pd(rsq00,rinv00);
221 /* Compute parameters for interactions between i and j atoms */
222 gmx_mm256_load_4pair_swizzle_pd(vdwioffsetptr0+vdwjidx0A,
223 vdwioffsetptr0+vdwjidx0B,
224 vdwioffsetptr0+vdwjidx0C,
225 vdwioffsetptr0+vdwjidx0D,
228 c6grid_00 = gmx_mm256_load_4real_swizzle_pd(vdwgridioffsetptr0+vdwjidx0A,
229 vdwgridioffsetptr0+vdwjidx0B,
230 vdwgridioffsetptr0+vdwjidx0C,
231 vdwgridioffsetptr0+vdwjidx0D);
233 /* Analytical LJ-PME */
234 rinvsix = _mm256_mul_pd(_mm256_mul_pd(rinvsq00,rinvsq00),rinvsq00);
235 ewcljrsq = _mm256_mul_pd(ewclj2,rsq00);
236 ewclj6 = _mm256_mul_pd(ewclj2,_mm256_mul_pd(ewclj2,ewclj2));
237 exponent = gmx_simd_exp_d(ewcljrsq);
238 /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
239 poly = _mm256_mul_pd(exponent,_mm256_add_pd(_mm256_sub_pd(one,ewcljrsq),_mm256_mul_pd(_mm256_mul_pd(ewcljrsq,ewcljrsq),one_half)));
240 /* vvdw6 = [C6 - C6grid * (1-poly)]/r6 */
241 vvdw6 = _mm256_mul_pd(_mm256_sub_pd(c6_00,_mm256_mul_pd(c6grid_00,_mm256_sub_pd(one,poly))),rinvsix);
242 vvdw12 = _mm256_mul_pd(c12_00,_mm256_mul_pd(rinvsix,rinvsix));
243 vvdw = _mm256_sub_pd(_mm256_mul_pd( _mm256_sub_pd(vvdw12 , _mm256_mul_pd(c12_00,_mm256_mul_pd(sh_vdw_invrcut6,sh_vdw_invrcut6))), one_twelfth) ,
244 _mm256_mul_pd( _mm256_sub_pd(vvdw6,_mm256_add_pd(_mm256_mul_pd(c6_00,sh_vdw_invrcut6),_mm256_mul_pd(c6grid_00,sh_lj_ewald))),one_sixth));
245 /* fvdw = vvdw12/r - (vvdw6/r + (C6grid * exponent * beta^6)/r) */
246 fvdw = _mm256_mul_pd(_mm256_sub_pd(vvdw12,_mm256_sub_pd(vvdw6,_mm256_mul_pd(_mm256_mul_pd(c6grid_00,one_sixth),_mm256_mul_pd(exponent,ewclj6)))),rinvsq00);
248 cutoff_mask = _mm256_cmp_pd(rsq00,rcutoff2,_CMP_LT_OQ);
250 /* Update potential sum for this i atom from the interaction with this j atom. */
251 vvdw = _mm256_and_pd(vvdw,cutoff_mask);
252 vvdwsum = _mm256_add_pd(vvdwsum,vvdw);
256 fscal = _mm256_and_pd(fscal,cutoff_mask);
258 /* Calculate temporary vectorial force */
259 tx = _mm256_mul_pd(fscal,dx00);
260 ty = _mm256_mul_pd(fscal,dy00);
261 tz = _mm256_mul_pd(fscal,dz00);
263 /* Update vectorial force */
264 fix0 = _mm256_add_pd(fix0,tx);
265 fiy0 = _mm256_add_pd(fiy0,ty);
266 fiz0 = _mm256_add_pd(fiz0,tz);
268 fjptrA = f+j_coord_offsetA;
269 fjptrB = f+j_coord_offsetB;
270 fjptrC = f+j_coord_offsetC;
271 fjptrD = f+j_coord_offsetD;
272 gmx_mm256_decrement_1rvec_4ptr_swizzle_pd(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
276 /* Inner loop uses 62 flops */
282 /* Get j neighbor index, and coordinate index */
283 jnrlistA = jjnr[jidx];
284 jnrlistB = jjnr[jidx+1];
285 jnrlistC = jjnr[jidx+2];
286 jnrlistD = jjnr[jidx+3];
287 /* Sign of each element will be negative for non-real atoms.
288 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
289 * so use it as val = _mm_andnot_pd(mask,val) to clear dummy entries.
291 tmpmask0 = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
293 tmpmask1 = _mm_permute_ps(tmpmask0,_GMX_MM_PERMUTE(3,3,2,2));
294 tmpmask0 = _mm_permute_ps(tmpmask0,_GMX_MM_PERMUTE(1,1,0,0));
295 dummy_mask = _mm256_castps_pd(gmx_mm256_set_m128(tmpmask1,tmpmask0));
297 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
298 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
299 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
300 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
301 j_coord_offsetA = DIM*jnrA;
302 j_coord_offsetB = DIM*jnrB;
303 j_coord_offsetC = DIM*jnrC;
304 j_coord_offsetD = DIM*jnrD;
306 /* load j atom coordinates */
307 gmx_mm256_load_1rvec_4ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
308 x+j_coord_offsetC,x+j_coord_offsetD,
311 /* Calculate displacement vector */
312 dx00 = _mm256_sub_pd(ix0,jx0);
313 dy00 = _mm256_sub_pd(iy0,jy0);
314 dz00 = _mm256_sub_pd(iz0,jz0);
316 /* Calculate squared distance and things based on it */
317 rsq00 = gmx_mm256_calc_rsq_pd(dx00,dy00,dz00);
319 rinv00 = gmx_mm256_invsqrt_pd(rsq00);
321 rinvsq00 = _mm256_mul_pd(rinv00,rinv00);
323 /* Load parameters for j particles */
324 vdwjidx0A = 2*vdwtype[jnrA+0];
325 vdwjidx0B = 2*vdwtype[jnrB+0];
326 vdwjidx0C = 2*vdwtype[jnrC+0];
327 vdwjidx0D = 2*vdwtype[jnrD+0];
329 /**************************
330 * CALCULATE INTERACTIONS *
331 **************************/
333 if (gmx_mm256_any_lt(rsq00,rcutoff2))
336 r00 = _mm256_mul_pd(rsq00,rinv00);
337 r00 = _mm256_andnot_pd(dummy_mask,r00);
339 /* Compute parameters for interactions between i and j atoms */
340 gmx_mm256_load_4pair_swizzle_pd(vdwioffsetptr0+vdwjidx0A,
341 vdwioffsetptr0+vdwjidx0B,
342 vdwioffsetptr0+vdwjidx0C,
343 vdwioffsetptr0+vdwjidx0D,
346 c6grid_00 = gmx_mm256_load_4real_swizzle_pd(vdwgridioffsetptr0+vdwjidx0A,
347 vdwgridioffsetptr0+vdwjidx0B,
348 vdwgridioffsetptr0+vdwjidx0C,
349 vdwgridioffsetptr0+vdwjidx0D);
351 /* Analytical LJ-PME */
352 rinvsix = _mm256_mul_pd(_mm256_mul_pd(rinvsq00,rinvsq00),rinvsq00);
353 ewcljrsq = _mm256_mul_pd(ewclj2,rsq00);
354 ewclj6 = _mm256_mul_pd(ewclj2,_mm256_mul_pd(ewclj2,ewclj2));
355 exponent = gmx_simd_exp_d(ewcljrsq);
356 /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
357 poly = _mm256_mul_pd(exponent,_mm256_add_pd(_mm256_sub_pd(one,ewcljrsq),_mm256_mul_pd(_mm256_mul_pd(ewcljrsq,ewcljrsq),one_half)));
358 /* vvdw6 = [C6 - C6grid * (1-poly)]/r6 */
359 vvdw6 = _mm256_mul_pd(_mm256_sub_pd(c6_00,_mm256_mul_pd(c6grid_00,_mm256_sub_pd(one,poly))),rinvsix);
360 vvdw12 = _mm256_mul_pd(c12_00,_mm256_mul_pd(rinvsix,rinvsix));
361 vvdw = _mm256_sub_pd(_mm256_mul_pd( _mm256_sub_pd(vvdw12 , _mm256_mul_pd(c12_00,_mm256_mul_pd(sh_vdw_invrcut6,sh_vdw_invrcut6))), one_twelfth) ,
362 _mm256_mul_pd( _mm256_sub_pd(vvdw6,_mm256_add_pd(_mm256_mul_pd(c6_00,sh_vdw_invrcut6),_mm256_mul_pd(c6grid_00,sh_lj_ewald))),one_sixth));
363 /* fvdw = vvdw12/r - (vvdw6/r + (C6grid * exponent * beta^6)/r) */
364 fvdw = _mm256_mul_pd(_mm256_sub_pd(vvdw12,_mm256_sub_pd(vvdw6,_mm256_mul_pd(_mm256_mul_pd(c6grid_00,one_sixth),_mm256_mul_pd(exponent,ewclj6)))),rinvsq00);
366 cutoff_mask = _mm256_cmp_pd(rsq00,rcutoff2,_CMP_LT_OQ);
368 /* Update potential sum for this i atom from the interaction with this j atom. */
369 vvdw = _mm256_and_pd(vvdw,cutoff_mask);
370 vvdw = _mm256_andnot_pd(dummy_mask,vvdw);
371 vvdwsum = _mm256_add_pd(vvdwsum,vvdw);
375 fscal = _mm256_and_pd(fscal,cutoff_mask);
377 fscal = _mm256_andnot_pd(dummy_mask,fscal);
379 /* Calculate temporary vectorial force */
380 tx = _mm256_mul_pd(fscal,dx00);
381 ty = _mm256_mul_pd(fscal,dy00);
382 tz = _mm256_mul_pd(fscal,dz00);
384 /* Update vectorial force */
385 fix0 = _mm256_add_pd(fix0,tx);
386 fiy0 = _mm256_add_pd(fiy0,ty);
387 fiz0 = _mm256_add_pd(fiz0,tz);
389 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
390 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
391 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
392 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
393 gmx_mm256_decrement_1rvec_4ptr_swizzle_pd(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
397 /* Inner loop uses 63 flops */
400 /* End of innermost loop */
402 gmx_mm256_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
403 f+i_coord_offset,fshift+i_shift_offset);
406 /* Update potential energies */
407 gmx_mm256_update_1pot_pd(vvdwsum,kernel_data->energygrp_vdw+ggid);
409 /* Increment number of inner iterations */
410 inneriter += j_index_end - j_index_start;
412 /* Outer loop uses 7 flops */
415 /* Increment number of outer iterations */
418 /* Update outer/inner flops */
420 inc_nrnb(nrnb,eNR_NBKERNEL_VDW_VF,outeriter*7 + inneriter*63);
423 * Gromacs nonbonded kernel: nb_kernel_ElecNone_VdwLJEwSh_GeomP1P1_F_avx_256_double
424 * Electrostatics interaction: None
425 * VdW interaction: LJEwald
426 * Geometry: Particle-Particle
427 * Calculate force/pot: Force
430 nb_kernel_ElecNone_VdwLJEwSh_GeomP1P1_F_avx_256_double
431 (t_nblist * gmx_restrict nlist,
432 rvec * gmx_restrict xx,
433 rvec * gmx_restrict ff,
434 t_forcerec * gmx_restrict fr,
435 t_mdatoms * gmx_restrict mdatoms,
436 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
437 t_nrnb * gmx_restrict nrnb)
439 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
440 * just 0 for non-waters.
441 * Suffixes A,B,C,D refer to j loop unrolling done with AVX, e.g. for the four different
442 * jnr indices corresponding to data put in the four positions in the SIMD register.
444 int i_shift_offset,i_coord_offset,outeriter,inneriter;
445 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
446 int jnrA,jnrB,jnrC,jnrD;
447 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
448 int jnrlistE,jnrlistF,jnrlistG,jnrlistH;
449 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
450 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
452 real *shiftvec,*fshift,*x,*f;
453 real *fjptrA,*fjptrB,*fjptrC,*fjptrD;
455 __m256d tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
456 real * vdwioffsetptr0;
457 real * vdwgridioffsetptr0;
458 __m256d ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
459 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
460 __m256d jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
461 __m256d dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
463 __m256d rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
466 __m256d one_sixth = _mm256_set1_pd(1.0/6.0);
467 __m256d one_twelfth = _mm256_set1_pd(1.0/12.0);
470 __m256d ewclj,ewclj2,ewclj6,ewcljrsq,poly,exponent,f6A,f6B,sh_lj_ewald;
471 __m256d one_half = _mm256_set1_pd(0.5);
472 __m256d minus_one = _mm256_set1_pd(-1.0);
473 __m256d dummy_mask,cutoff_mask;
474 __m128 tmpmask0,tmpmask1;
475 __m256d signbit = _mm256_castsi256_pd( _mm256_set1_epi32(0x80000000) );
476 __m256d one = _mm256_set1_pd(1.0);
477 __m256d two = _mm256_set1_pd(2.0);
483 jindex = nlist->jindex;
485 shiftidx = nlist->shift;
487 shiftvec = fr->shift_vec[0];
488 fshift = fr->fshift[0];
489 nvdwtype = fr->ntype;
491 vdwtype = mdatoms->typeA;
492 vdwgridparam = fr->ljpme_c6grid;
493 sh_lj_ewald = _mm256_set1_pd(fr->ic->sh_lj_ewald);
494 ewclj = _mm256_set1_pd(fr->ewaldcoeff_lj);
495 ewclj2 = _mm256_mul_pd(minus_one,_mm256_mul_pd(ewclj,ewclj));
497 rcutoff_scalar = fr->rvdw;
498 rcutoff = _mm256_set1_pd(rcutoff_scalar);
499 rcutoff2 = _mm256_mul_pd(rcutoff,rcutoff);
501 sh_vdw_invrcut6 = _mm256_set1_pd(fr->ic->sh_invrc6);
502 rvdw = _mm256_set1_pd(fr->rvdw);
504 /* Avoid stupid compiler warnings */
505 jnrA = jnrB = jnrC = jnrD = 0;
514 for(iidx=0;iidx<4*DIM;iidx++)
519 /* Start outer loop over neighborlists */
520 for(iidx=0; iidx<nri; iidx++)
522 /* Load shift vector for this list */
523 i_shift_offset = DIM*shiftidx[iidx];
525 /* Load limits for loop over neighbors */
526 j_index_start = jindex[iidx];
527 j_index_end = jindex[iidx+1];
529 /* Get outer coordinate index */
531 i_coord_offset = DIM*inr;
533 /* Load i particle coords and add shift vector */
534 gmx_mm256_load_shift_and_1rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
536 fix0 = _mm256_setzero_pd();
537 fiy0 = _mm256_setzero_pd();
538 fiz0 = _mm256_setzero_pd();
540 /* Load parameters for i particles */
541 vdwioffsetptr0 = vdwparam+2*nvdwtype*vdwtype[inr+0];
542 vdwgridioffsetptr0 = vdwgridparam+2*nvdwtype*vdwtype[inr+0];
544 /* Start inner kernel loop */
545 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
548 /* Get j neighbor index, and coordinate index */
553 j_coord_offsetA = DIM*jnrA;
554 j_coord_offsetB = DIM*jnrB;
555 j_coord_offsetC = DIM*jnrC;
556 j_coord_offsetD = DIM*jnrD;
558 /* load j atom coordinates */
559 gmx_mm256_load_1rvec_4ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
560 x+j_coord_offsetC,x+j_coord_offsetD,
563 /* Calculate displacement vector */
564 dx00 = _mm256_sub_pd(ix0,jx0);
565 dy00 = _mm256_sub_pd(iy0,jy0);
566 dz00 = _mm256_sub_pd(iz0,jz0);
568 /* Calculate squared distance and things based on it */
569 rsq00 = gmx_mm256_calc_rsq_pd(dx00,dy00,dz00);
571 rinv00 = gmx_mm256_invsqrt_pd(rsq00);
573 rinvsq00 = _mm256_mul_pd(rinv00,rinv00);
575 /* Load parameters for j particles */
576 vdwjidx0A = 2*vdwtype[jnrA+0];
577 vdwjidx0B = 2*vdwtype[jnrB+0];
578 vdwjidx0C = 2*vdwtype[jnrC+0];
579 vdwjidx0D = 2*vdwtype[jnrD+0];
581 /**************************
582 * CALCULATE INTERACTIONS *
583 **************************/
585 if (gmx_mm256_any_lt(rsq00,rcutoff2))
588 r00 = _mm256_mul_pd(rsq00,rinv00);
590 /* Compute parameters for interactions between i and j atoms */
591 gmx_mm256_load_4pair_swizzle_pd(vdwioffsetptr0+vdwjidx0A,
592 vdwioffsetptr0+vdwjidx0B,
593 vdwioffsetptr0+vdwjidx0C,
594 vdwioffsetptr0+vdwjidx0D,
597 c6grid_00 = gmx_mm256_load_4real_swizzle_pd(vdwgridioffsetptr0+vdwjidx0A,
598 vdwgridioffsetptr0+vdwjidx0B,
599 vdwgridioffsetptr0+vdwjidx0C,
600 vdwgridioffsetptr0+vdwjidx0D);
602 /* Analytical LJ-PME */
603 rinvsix = _mm256_mul_pd(_mm256_mul_pd(rinvsq00,rinvsq00),rinvsq00);
604 ewcljrsq = _mm256_mul_pd(ewclj2,rsq00);
605 ewclj6 = _mm256_mul_pd(ewclj2,_mm256_mul_pd(ewclj2,ewclj2));
606 exponent = gmx_simd_exp_d(ewcljrsq);
607 /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
608 poly = _mm256_mul_pd(exponent,_mm256_add_pd(_mm256_sub_pd(one,ewcljrsq),_mm256_mul_pd(_mm256_mul_pd(ewcljrsq,ewcljrsq),one_half)));
609 /* f6A = 6 * C6grid * (1 - poly) */
610 f6A = _mm256_mul_pd(c6grid_00,_mm256_sub_pd(one,poly));
611 /* f6B = C6grid * exponent * beta^6 */
612 f6B = _mm256_mul_pd(_mm256_mul_pd(c6grid_00,one_sixth),_mm256_mul_pd(exponent,ewclj6));
613 /* fvdw = 12*C12/r13 - ((6*C6 - f6A)/r6 + f6B)/r */
614 fvdw = _mm256_mul_pd(_mm256_add_pd(_mm256_mul_pd(_mm256_sub_pd(_mm256_mul_pd(c12_00,rinvsix),_mm256_sub_pd(c6_00,f6A)),rinvsix),f6B),rinvsq00);
616 cutoff_mask = _mm256_cmp_pd(rsq00,rcutoff2,_CMP_LT_OQ);
620 fscal = _mm256_and_pd(fscal,cutoff_mask);
622 /* Calculate temporary vectorial force */
623 tx = _mm256_mul_pd(fscal,dx00);
624 ty = _mm256_mul_pd(fscal,dy00);
625 tz = _mm256_mul_pd(fscal,dz00);
627 /* Update vectorial force */
628 fix0 = _mm256_add_pd(fix0,tx);
629 fiy0 = _mm256_add_pd(fiy0,ty);
630 fiz0 = _mm256_add_pd(fiz0,tz);
632 fjptrA = f+j_coord_offsetA;
633 fjptrB = f+j_coord_offsetB;
634 fjptrC = f+j_coord_offsetC;
635 fjptrD = f+j_coord_offsetD;
636 gmx_mm256_decrement_1rvec_4ptr_swizzle_pd(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
640 /* Inner loop uses 49 flops */
646 /* Get j neighbor index, and coordinate index */
647 jnrlistA = jjnr[jidx];
648 jnrlistB = jjnr[jidx+1];
649 jnrlistC = jjnr[jidx+2];
650 jnrlistD = jjnr[jidx+3];
651 /* Sign of each element will be negative for non-real atoms.
652 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
653 * so use it as val = _mm_andnot_pd(mask,val) to clear dummy entries.
655 tmpmask0 = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
657 tmpmask1 = _mm_permute_ps(tmpmask0,_GMX_MM_PERMUTE(3,3,2,2));
658 tmpmask0 = _mm_permute_ps(tmpmask0,_GMX_MM_PERMUTE(1,1,0,0));
659 dummy_mask = _mm256_castps_pd(gmx_mm256_set_m128(tmpmask1,tmpmask0));
661 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
662 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
663 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
664 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
665 j_coord_offsetA = DIM*jnrA;
666 j_coord_offsetB = DIM*jnrB;
667 j_coord_offsetC = DIM*jnrC;
668 j_coord_offsetD = DIM*jnrD;
670 /* load j atom coordinates */
671 gmx_mm256_load_1rvec_4ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
672 x+j_coord_offsetC,x+j_coord_offsetD,
675 /* Calculate displacement vector */
676 dx00 = _mm256_sub_pd(ix0,jx0);
677 dy00 = _mm256_sub_pd(iy0,jy0);
678 dz00 = _mm256_sub_pd(iz0,jz0);
680 /* Calculate squared distance and things based on it */
681 rsq00 = gmx_mm256_calc_rsq_pd(dx00,dy00,dz00);
683 rinv00 = gmx_mm256_invsqrt_pd(rsq00);
685 rinvsq00 = _mm256_mul_pd(rinv00,rinv00);
687 /* Load parameters for j particles */
688 vdwjidx0A = 2*vdwtype[jnrA+0];
689 vdwjidx0B = 2*vdwtype[jnrB+0];
690 vdwjidx0C = 2*vdwtype[jnrC+0];
691 vdwjidx0D = 2*vdwtype[jnrD+0];
693 /**************************
694 * CALCULATE INTERACTIONS *
695 **************************/
697 if (gmx_mm256_any_lt(rsq00,rcutoff2))
700 r00 = _mm256_mul_pd(rsq00,rinv00);
701 r00 = _mm256_andnot_pd(dummy_mask,r00);
703 /* Compute parameters for interactions between i and j atoms */
704 gmx_mm256_load_4pair_swizzle_pd(vdwioffsetptr0+vdwjidx0A,
705 vdwioffsetptr0+vdwjidx0B,
706 vdwioffsetptr0+vdwjidx0C,
707 vdwioffsetptr0+vdwjidx0D,
710 c6grid_00 = gmx_mm256_load_4real_swizzle_pd(vdwgridioffsetptr0+vdwjidx0A,
711 vdwgridioffsetptr0+vdwjidx0B,
712 vdwgridioffsetptr0+vdwjidx0C,
713 vdwgridioffsetptr0+vdwjidx0D);
715 /* Analytical LJ-PME */
716 rinvsix = _mm256_mul_pd(_mm256_mul_pd(rinvsq00,rinvsq00),rinvsq00);
717 ewcljrsq = _mm256_mul_pd(ewclj2,rsq00);
718 ewclj6 = _mm256_mul_pd(ewclj2,_mm256_mul_pd(ewclj2,ewclj2));
719 exponent = gmx_simd_exp_d(ewcljrsq);
720 /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
721 poly = _mm256_mul_pd(exponent,_mm256_add_pd(_mm256_sub_pd(one,ewcljrsq),_mm256_mul_pd(_mm256_mul_pd(ewcljrsq,ewcljrsq),one_half)));
722 /* f6A = 6 * C6grid * (1 - poly) */
723 f6A = _mm256_mul_pd(c6grid_00,_mm256_sub_pd(one,poly));
724 /* f6B = C6grid * exponent * beta^6 */
725 f6B = _mm256_mul_pd(_mm256_mul_pd(c6grid_00,one_sixth),_mm256_mul_pd(exponent,ewclj6));
726 /* fvdw = 12*C12/r13 - ((6*C6 - f6A)/r6 + f6B)/r */
727 fvdw = _mm256_mul_pd(_mm256_add_pd(_mm256_mul_pd(_mm256_sub_pd(_mm256_mul_pd(c12_00,rinvsix),_mm256_sub_pd(c6_00,f6A)),rinvsix),f6B),rinvsq00);
729 cutoff_mask = _mm256_cmp_pd(rsq00,rcutoff2,_CMP_LT_OQ);
733 fscal = _mm256_and_pd(fscal,cutoff_mask);
735 fscal = _mm256_andnot_pd(dummy_mask,fscal);
737 /* Calculate temporary vectorial force */
738 tx = _mm256_mul_pd(fscal,dx00);
739 ty = _mm256_mul_pd(fscal,dy00);
740 tz = _mm256_mul_pd(fscal,dz00);
742 /* Update vectorial force */
743 fix0 = _mm256_add_pd(fix0,tx);
744 fiy0 = _mm256_add_pd(fiy0,ty);
745 fiz0 = _mm256_add_pd(fiz0,tz);
747 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
748 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
749 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
750 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
751 gmx_mm256_decrement_1rvec_4ptr_swizzle_pd(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
755 /* Inner loop uses 50 flops */
758 /* End of innermost loop */
760 gmx_mm256_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
761 f+i_coord_offset,fshift+i_shift_offset);
763 /* Increment number of inner iterations */
764 inneriter += j_index_end - j_index_start;
766 /* Outer loop uses 6 flops */
769 /* Increment number of outer iterations */
772 /* Update outer/inner flops */
774 inc_nrnb(nrnb,eNR_NBKERNEL_VDW_F,outeriter*6 + inneriter*50);