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_VdwLJEw_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_VdwLJEw_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 /* Avoid stupid compiler warnings */
126 jnrA = jnrB = jnrC = jnrD = 0;
135 for(iidx=0;iidx<4*DIM;iidx++)
140 /* Start outer loop over neighborlists */
141 for(iidx=0; iidx<nri; iidx++)
143 /* Load shift vector for this list */
144 i_shift_offset = DIM*shiftidx[iidx];
146 /* Load limits for loop over neighbors */
147 j_index_start = jindex[iidx];
148 j_index_end = jindex[iidx+1];
150 /* Get outer coordinate index */
152 i_coord_offset = DIM*inr;
154 /* Load i particle coords and add shift vector */
155 gmx_mm256_load_shift_and_1rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
157 fix0 = _mm256_setzero_pd();
158 fiy0 = _mm256_setzero_pd();
159 fiz0 = _mm256_setzero_pd();
161 /* Load parameters for i particles */
162 vdwioffsetptr0 = vdwparam+2*nvdwtype*vdwtype[inr+0];
163 vdwgridioffsetptr0 = vdwgridparam+2*nvdwtype*vdwtype[inr+0];
165 /* Reset potential sums */
166 vvdwsum = _mm256_setzero_pd();
168 /* Start inner kernel loop */
169 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
172 /* Get j neighbor index, and coordinate index */
177 j_coord_offsetA = DIM*jnrA;
178 j_coord_offsetB = DIM*jnrB;
179 j_coord_offsetC = DIM*jnrC;
180 j_coord_offsetD = DIM*jnrD;
182 /* load j atom coordinates */
183 gmx_mm256_load_1rvec_4ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
184 x+j_coord_offsetC,x+j_coord_offsetD,
187 /* Calculate displacement vector */
188 dx00 = _mm256_sub_pd(ix0,jx0);
189 dy00 = _mm256_sub_pd(iy0,jy0);
190 dz00 = _mm256_sub_pd(iz0,jz0);
192 /* Calculate squared distance and things based on it */
193 rsq00 = gmx_mm256_calc_rsq_pd(dx00,dy00,dz00);
195 rinv00 = gmx_mm256_invsqrt_pd(rsq00);
197 rinvsq00 = _mm256_mul_pd(rinv00,rinv00);
199 /* Load parameters for j particles */
200 vdwjidx0A = 2*vdwtype[jnrA+0];
201 vdwjidx0B = 2*vdwtype[jnrB+0];
202 vdwjidx0C = 2*vdwtype[jnrC+0];
203 vdwjidx0D = 2*vdwtype[jnrD+0];
205 /**************************
206 * CALCULATE INTERACTIONS *
207 **************************/
209 r00 = _mm256_mul_pd(rsq00,rinv00);
211 /* Compute parameters for interactions between i and j atoms */
212 gmx_mm256_load_4pair_swizzle_pd(vdwioffsetptr0+vdwjidx0A,
213 vdwioffsetptr0+vdwjidx0B,
214 vdwioffsetptr0+vdwjidx0C,
215 vdwioffsetptr0+vdwjidx0D,
218 c6grid_00 = gmx_mm256_load_4real_swizzle_pd(vdwgridioffsetptr0+vdwjidx0A,
219 vdwgridioffsetptr0+vdwjidx0B,
220 vdwgridioffsetptr0+vdwjidx0C,
221 vdwgridioffsetptr0+vdwjidx0D);
223 /* Analytical LJ-PME */
224 rinvsix = _mm256_mul_pd(_mm256_mul_pd(rinvsq00,rinvsq00),rinvsq00);
225 ewcljrsq = _mm256_mul_pd(ewclj2,rsq00);
226 ewclj6 = _mm256_mul_pd(ewclj2,_mm256_mul_pd(ewclj2,ewclj2));
227 exponent = gmx_simd_exp_d(ewcljrsq);
228 /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
229 poly = _mm256_mul_pd(exponent,_mm256_add_pd(_mm256_sub_pd(one,ewcljrsq),_mm256_mul_pd(_mm256_mul_pd(ewcljrsq,ewcljrsq),one_half)));
230 /* vvdw6 = [C6 - C6grid * (1-poly)]/r6 */
231 vvdw6 = _mm256_mul_pd(_mm256_sub_pd(c6_00,_mm256_mul_pd(c6grid_00,_mm256_sub_pd(one,poly))),rinvsix);
232 vvdw12 = _mm256_mul_pd(c12_00,_mm256_mul_pd(rinvsix,rinvsix));
233 vvdw = _mm256_sub_pd(_mm256_mul_pd(vvdw12,one_twelfth),_mm256_mul_pd(vvdw6,one_sixth));
234 /* fvdw = vvdw12/r - (vvdw6/r + (C6grid * exponent * beta^6)/r) */
235 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);
237 /* Update potential sum for this i atom from the interaction with this j atom. */
238 vvdwsum = _mm256_add_pd(vvdwsum,vvdw);
242 /* Calculate temporary vectorial force */
243 tx = _mm256_mul_pd(fscal,dx00);
244 ty = _mm256_mul_pd(fscal,dy00);
245 tz = _mm256_mul_pd(fscal,dz00);
247 /* Update vectorial force */
248 fix0 = _mm256_add_pd(fix0,tx);
249 fiy0 = _mm256_add_pd(fiy0,ty);
250 fiz0 = _mm256_add_pd(fiz0,tz);
252 fjptrA = f+j_coord_offsetA;
253 fjptrB = f+j_coord_offsetB;
254 fjptrC = f+j_coord_offsetC;
255 fjptrD = f+j_coord_offsetD;
256 gmx_mm256_decrement_1rvec_4ptr_swizzle_pd(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
258 /* Inner loop uses 54 flops */
264 /* Get j neighbor index, and coordinate index */
265 jnrlistA = jjnr[jidx];
266 jnrlistB = jjnr[jidx+1];
267 jnrlistC = jjnr[jidx+2];
268 jnrlistD = jjnr[jidx+3];
269 /* Sign of each element will be negative for non-real atoms.
270 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
271 * so use it as val = _mm_andnot_pd(mask,val) to clear dummy entries.
273 tmpmask0 = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
275 tmpmask1 = _mm_permute_ps(tmpmask0,_GMX_MM_PERMUTE(3,3,2,2));
276 tmpmask0 = _mm_permute_ps(tmpmask0,_GMX_MM_PERMUTE(1,1,0,0));
277 dummy_mask = _mm256_castps_pd(gmx_mm256_set_m128(tmpmask1,tmpmask0));
279 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
280 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
281 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
282 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
283 j_coord_offsetA = DIM*jnrA;
284 j_coord_offsetB = DIM*jnrB;
285 j_coord_offsetC = DIM*jnrC;
286 j_coord_offsetD = DIM*jnrD;
288 /* load j atom coordinates */
289 gmx_mm256_load_1rvec_4ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
290 x+j_coord_offsetC,x+j_coord_offsetD,
293 /* Calculate displacement vector */
294 dx00 = _mm256_sub_pd(ix0,jx0);
295 dy00 = _mm256_sub_pd(iy0,jy0);
296 dz00 = _mm256_sub_pd(iz0,jz0);
298 /* Calculate squared distance and things based on it */
299 rsq00 = gmx_mm256_calc_rsq_pd(dx00,dy00,dz00);
301 rinv00 = gmx_mm256_invsqrt_pd(rsq00);
303 rinvsq00 = _mm256_mul_pd(rinv00,rinv00);
305 /* Load parameters for j particles */
306 vdwjidx0A = 2*vdwtype[jnrA+0];
307 vdwjidx0B = 2*vdwtype[jnrB+0];
308 vdwjidx0C = 2*vdwtype[jnrC+0];
309 vdwjidx0D = 2*vdwtype[jnrD+0];
311 /**************************
312 * CALCULATE INTERACTIONS *
313 **************************/
315 r00 = _mm256_mul_pd(rsq00,rinv00);
316 r00 = _mm256_andnot_pd(dummy_mask,r00);
318 /* Compute parameters for interactions between i and j atoms */
319 gmx_mm256_load_4pair_swizzle_pd(vdwioffsetptr0+vdwjidx0A,
320 vdwioffsetptr0+vdwjidx0B,
321 vdwioffsetptr0+vdwjidx0C,
322 vdwioffsetptr0+vdwjidx0D,
325 c6grid_00 = gmx_mm256_load_4real_swizzle_pd(vdwgridioffsetptr0+vdwjidx0A,
326 vdwgridioffsetptr0+vdwjidx0B,
327 vdwgridioffsetptr0+vdwjidx0C,
328 vdwgridioffsetptr0+vdwjidx0D);
330 /* Analytical LJ-PME */
331 rinvsix = _mm256_mul_pd(_mm256_mul_pd(rinvsq00,rinvsq00),rinvsq00);
332 ewcljrsq = _mm256_mul_pd(ewclj2,rsq00);
333 ewclj6 = _mm256_mul_pd(ewclj2,_mm256_mul_pd(ewclj2,ewclj2));
334 exponent = gmx_simd_exp_d(ewcljrsq);
335 /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
336 poly = _mm256_mul_pd(exponent,_mm256_add_pd(_mm256_sub_pd(one,ewcljrsq),_mm256_mul_pd(_mm256_mul_pd(ewcljrsq,ewcljrsq),one_half)));
337 /* vvdw6 = [C6 - C6grid * (1-poly)]/r6 */
338 vvdw6 = _mm256_mul_pd(_mm256_sub_pd(c6_00,_mm256_mul_pd(c6grid_00,_mm256_sub_pd(one,poly))),rinvsix);
339 vvdw12 = _mm256_mul_pd(c12_00,_mm256_mul_pd(rinvsix,rinvsix));
340 vvdw = _mm256_sub_pd(_mm256_mul_pd(vvdw12,one_twelfth),_mm256_mul_pd(vvdw6,one_sixth));
341 /* fvdw = vvdw12/r - (vvdw6/r + (C6grid * exponent * beta^6)/r) */
342 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);
344 /* Update potential sum for this i atom from the interaction with this j atom. */
345 vvdw = _mm256_andnot_pd(dummy_mask,vvdw);
346 vvdwsum = _mm256_add_pd(vvdwsum,vvdw);
350 fscal = _mm256_andnot_pd(dummy_mask,fscal);
352 /* Calculate temporary vectorial force */
353 tx = _mm256_mul_pd(fscal,dx00);
354 ty = _mm256_mul_pd(fscal,dy00);
355 tz = _mm256_mul_pd(fscal,dz00);
357 /* Update vectorial force */
358 fix0 = _mm256_add_pd(fix0,tx);
359 fiy0 = _mm256_add_pd(fiy0,ty);
360 fiz0 = _mm256_add_pd(fiz0,tz);
362 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
363 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
364 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
365 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
366 gmx_mm256_decrement_1rvec_4ptr_swizzle_pd(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
368 /* Inner loop uses 55 flops */
371 /* End of innermost loop */
373 gmx_mm256_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
374 f+i_coord_offset,fshift+i_shift_offset);
377 /* Update potential energies */
378 gmx_mm256_update_1pot_pd(vvdwsum,kernel_data->energygrp_vdw+ggid);
380 /* Increment number of inner iterations */
381 inneriter += j_index_end - j_index_start;
383 /* Outer loop uses 7 flops */
386 /* Increment number of outer iterations */
389 /* Update outer/inner flops */
391 inc_nrnb(nrnb,eNR_NBKERNEL_VDW_VF,outeriter*7 + inneriter*55);
394 * Gromacs nonbonded kernel: nb_kernel_ElecNone_VdwLJEw_GeomP1P1_F_avx_256_double
395 * Electrostatics interaction: None
396 * VdW interaction: LJEwald
397 * Geometry: Particle-Particle
398 * Calculate force/pot: Force
401 nb_kernel_ElecNone_VdwLJEw_GeomP1P1_F_avx_256_double
402 (t_nblist * gmx_restrict nlist,
403 rvec * gmx_restrict xx,
404 rvec * gmx_restrict ff,
405 t_forcerec * gmx_restrict fr,
406 t_mdatoms * gmx_restrict mdatoms,
407 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
408 t_nrnb * gmx_restrict nrnb)
410 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
411 * just 0 for non-waters.
412 * Suffixes A,B,C,D refer to j loop unrolling done with AVX, e.g. for the four different
413 * jnr indices corresponding to data put in the four positions in the SIMD register.
415 int i_shift_offset,i_coord_offset,outeriter,inneriter;
416 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
417 int jnrA,jnrB,jnrC,jnrD;
418 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
419 int jnrlistE,jnrlistF,jnrlistG,jnrlistH;
420 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
421 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
423 real *shiftvec,*fshift,*x,*f;
424 real *fjptrA,*fjptrB,*fjptrC,*fjptrD;
426 __m256d tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
427 real * vdwioffsetptr0;
428 real * vdwgridioffsetptr0;
429 __m256d ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
430 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
431 __m256d jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
432 __m256d dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
434 __m256d rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
437 __m256d one_sixth = _mm256_set1_pd(1.0/6.0);
438 __m256d one_twelfth = _mm256_set1_pd(1.0/12.0);
441 __m256d ewclj,ewclj2,ewclj6,ewcljrsq,poly,exponent,f6A,f6B,sh_lj_ewald;
442 __m256d one_half = _mm256_set1_pd(0.5);
443 __m256d minus_one = _mm256_set1_pd(-1.0);
444 __m256d dummy_mask,cutoff_mask;
445 __m128 tmpmask0,tmpmask1;
446 __m256d signbit = _mm256_castsi256_pd( _mm256_set1_epi32(0x80000000) );
447 __m256d one = _mm256_set1_pd(1.0);
448 __m256d two = _mm256_set1_pd(2.0);
454 jindex = nlist->jindex;
456 shiftidx = nlist->shift;
458 shiftvec = fr->shift_vec[0];
459 fshift = fr->fshift[0];
460 nvdwtype = fr->ntype;
462 vdwtype = mdatoms->typeA;
463 vdwgridparam = fr->ljpme_c6grid;
464 sh_lj_ewald = _mm256_set1_pd(fr->ic->sh_lj_ewald);
465 ewclj = _mm256_set1_pd(fr->ewaldcoeff_lj);
466 ewclj2 = _mm256_mul_pd(minus_one,_mm256_mul_pd(ewclj,ewclj));
468 /* Avoid stupid compiler warnings */
469 jnrA = jnrB = jnrC = jnrD = 0;
478 for(iidx=0;iidx<4*DIM;iidx++)
483 /* Start outer loop over neighborlists */
484 for(iidx=0; iidx<nri; iidx++)
486 /* Load shift vector for this list */
487 i_shift_offset = DIM*shiftidx[iidx];
489 /* Load limits for loop over neighbors */
490 j_index_start = jindex[iidx];
491 j_index_end = jindex[iidx+1];
493 /* Get outer coordinate index */
495 i_coord_offset = DIM*inr;
497 /* Load i particle coords and add shift vector */
498 gmx_mm256_load_shift_and_1rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
500 fix0 = _mm256_setzero_pd();
501 fiy0 = _mm256_setzero_pd();
502 fiz0 = _mm256_setzero_pd();
504 /* Load parameters for i particles */
505 vdwioffsetptr0 = vdwparam+2*nvdwtype*vdwtype[inr+0];
506 vdwgridioffsetptr0 = vdwgridparam+2*nvdwtype*vdwtype[inr+0];
508 /* Start inner kernel loop */
509 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
512 /* Get j neighbor index, and coordinate index */
517 j_coord_offsetA = DIM*jnrA;
518 j_coord_offsetB = DIM*jnrB;
519 j_coord_offsetC = DIM*jnrC;
520 j_coord_offsetD = DIM*jnrD;
522 /* load j atom coordinates */
523 gmx_mm256_load_1rvec_4ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
524 x+j_coord_offsetC,x+j_coord_offsetD,
527 /* Calculate displacement vector */
528 dx00 = _mm256_sub_pd(ix0,jx0);
529 dy00 = _mm256_sub_pd(iy0,jy0);
530 dz00 = _mm256_sub_pd(iz0,jz0);
532 /* Calculate squared distance and things based on it */
533 rsq00 = gmx_mm256_calc_rsq_pd(dx00,dy00,dz00);
535 rinv00 = gmx_mm256_invsqrt_pd(rsq00);
537 rinvsq00 = _mm256_mul_pd(rinv00,rinv00);
539 /* Load parameters for j particles */
540 vdwjidx0A = 2*vdwtype[jnrA+0];
541 vdwjidx0B = 2*vdwtype[jnrB+0];
542 vdwjidx0C = 2*vdwtype[jnrC+0];
543 vdwjidx0D = 2*vdwtype[jnrD+0];
545 /**************************
546 * CALCULATE INTERACTIONS *
547 **************************/
549 r00 = _mm256_mul_pd(rsq00,rinv00);
551 /* Compute parameters for interactions between i and j atoms */
552 gmx_mm256_load_4pair_swizzle_pd(vdwioffsetptr0+vdwjidx0A,
553 vdwioffsetptr0+vdwjidx0B,
554 vdwioffsetptr0+vdwjidx0C,
555 vdwioffsetptr0+vdwjidx0D,
558 c6grid_00 = gmx_mm256_load_4real_swizzle_pd(vdwgridioffsetptr0+vdwjidx0A,
559 vdwgridioffsetptr0+vdwjidx0B,
560 vdwgridioffsetptr0+vdwjidx0C,
561 vdwgridioffsetptr0+vdwjidx0D);
563 /* Analytical LJ-PME */
564 rinvsix = _mm256_mul_pd(_mm256_mul_pd(rinvsq00,rinvsq00),rinvsq00);
565 ewcljrsq = _mm256_mul_pd(ewclj2,rsq00);
566 ewclj6 = _mm256_mul_pd(ewclj2,_mm256_mul_pd(ewclj2,ewclj2));
567 exponent = gmx_simd_exp_d(ewcljrsq);
568 /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
569 poly = _mm256_mul_pd(exponent,_mm256_add_pd(_mm256_sub_pd(one,ewcljrsq),_mm256_mul_pd(_mm256_mul_pd(ewcljrsq,ewcljrsq),one_half)));
570 /* f6A = 6 * C6grid * (1 - poly) */
571 f6A = _mm256_mul_pd(c6grid_00,_mm256_sub_pd(one,poly));
572 /* f6B = C6grid * exponent * beta^6 */
573 f6B = _mm256_mul_pd(_mm256_mul_pd(c6grid_00,one_sixth),_mm256_mul_pd(exponent,ewclj6));
574 /* fvdw = 12*C12/r13 - ((6*C6 - f6A)/r6 + f6B)/r */
575 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);
579 /* Calculate temporary vectorial force */
580 tx = _mm256_mul_pd(fscal,dx00);
581 ty = _mm256_mul_pd(fscal,dy00);
582 tz = _mm256_mul_pd(fscal,dz00);
584 /* Update vectorial force */
585 fix0 = _mm256_add_pd(fix0,tx);
586 fiy0 = _mm256_add_pd(fiy0,ty);
587 fiz0 = _mm256_add_pd(fiz0,tz);
589 fjptrA = f+j_coord_offsetA;
590 fjptrB = f+j_coord_offsetB;
591 fjptrC = f+j_coord_offsetC;
592 fjptrD = f+j_coord_offsetD;
593 gmx_mm256_decrement_1rvec_4ptr_swizzle_pd(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
595 /* Inner loop uses 46 flops */
601 /* Get j neighbor index, and coordinate index */
602 jnrlistA = jjnr[jidx];
603 jnrlistB = jjnr[jidx+1];
604 jnrlistC = jjnr[jidx+2];
605 jnrlistD = jjnr[jidx+3];
606 /* Sign of each element will be negative for non-real atoms.
607 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
608 * so use it as val = _mm_andnot_pd(mask,val) to clear dummy entries.
610 tmpmask0 = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
612 tmpmask1 = _mm_permute_ps(tmpmask0,_GMX_MM_PERMUTE(3,3,2,2));
613 tmpmask0 = _mm_permute_ps(tmpmask0,_GMX_MM_PERMUTE(1,1,0,0));
614 dummy_mask = _mm256_castps_pd(gmx_mm256_set_m128(tmpmask1,tmpmask0));
616 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
617 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
618 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
619 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
620 j_coord_offsetA = DIM*jnrA;
621 j_coord_offsetB = DIM*jnrB;
622 j_coord_offsetC = DIM*jnrC;
623 j_coord_offsetD = DIM*jnrD;
625 /* load j atom coordinates */
626 gmx_mm256_load_1rvec_4ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
627 x+j_coord_offsetC,x+j_coord_offsetD,
630 /* Calculate displacement vector */
631 dx00 = _mm256_sub_pd(ix0,jx0);
632 dy00 = _mm256_sub_pd(iy0,jy0);
633 dz00 = _mm256_sub_pd(iz0,jz0);
635 /* Calculate squared distance and things based on it */
636 rsq00 = gmx_mm256_calc_rsq_pd(dx00,dy00,dz00);
638 rinv00 = gmx_mm256_invsqrt_pd(rsq00);
640 rinvsq00 = _mm256_mul_pd(rinv00,rinv00);
642 /* Load parameters for j particles */
643 vdwjidx0A = 2*vdwtype[jnrA+0];
644 vdwjidx0B = 2*vdwtype[jnrB+0];
645 vdwjidx0C = 2*vdwtype[jnrC+0];
646 vdwjidx0D = 2*vdwtype[jnrD+0];
648 /**************************
649 * CALCULATE INTERACTIONS *
650 **************************/
652 r00 = _mm256_mul_pd(rsq00,rinv00);
653 r00 = _mm256_andnot_pd(dummy_mask,r00);
655 /* Compute parameters for interactions between i and j atoms */
656 gmx_mm256_load_4pair_swizzle_pd(vdwioffsetptr0+vdwjidx0A,
657 vdwioffsetptr0+vdwjidx0B,
658 vdwioffsetptr0+vdwjidx0C,
659 vdwioffsetptr0+vdwjidx0D,
662 c6grid_00 = gmx_mm256_load_4real_swizzle_pd(vdwgridioffsetptr0+vdwjidx0A,
663 vdwgridioffsetptr0+vdwjidx0B,
664 vdwgridioffsetptr0+vdwjidx0C,
665 vdwgridioffsetptr0+vdwjidx0D);
667 /* Analytical LJ-PME */
668 rinvsix = _mm256_mul_pd(_mm256_mul_pd(rinvsq00,rinvsq00),rinvsq00);
669 ewcljrsq = _mm256_mul_pd(ewclj2,rsq00);
670 ewclj6 = _mm256_mul_pd(ewclj2,_mm256_mul_pd(ewclj2,ewclj2));
671 exponent = gmx_simd_exp_d(ewcljrsq);
672 /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
673 poly = _mm256_mul_pd(exponent,_mm256_add_pd(_mm256_sub_pd(one,ewcljrsq),_mm256_mul_pd(_mm256_mul_pd(ewcljrsq,ewcljrsq),one_half)));
674 /* f6A = 6 * C6grid * (1 - poly) */
675 f6A = _mm256_mul_pd(c6grid_00,_mm256_sub_pd(one,poly));
676 /* f6B = C6grid * exponent * beta^6 */
677 f6B = _mm256_mul_pd(_mm256_mul_pd(c6grid_00,one_sixth),_mm256_mul_pd(exponent,ewclj6));
678 /* fvdw = 12*C12/r13 - ((6*C6 - f6A)/r6 + f6B)/r */
679 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);
683 fscal = _mm256_andnot_pd(dummy_mask,fscal);
685 /* Calculate temporary vectorial force */
686 tx = _mm256_mul_pd(fscal,dx00);
687 ty = _mm256_mul_pd(fscal,dy00);
688 tz = _mm256_mul_pd(fscal,dz00);
690 /* Update vectorial force */
691 fix0 = _mm256_add_pd(fix0,tx);
692 fiy0 = _mm256_add_pd(fiy0,ty);
693 fiz0 = _mm256_add_pd(fiz0,tz);
695 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
696 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
697 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
698 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
699 gmx_mm256_decrement_1rvec_4ptr_swizzle_pd(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
701 /* Inner loop uses 47 flops */
704 /* End of innermost loop */
706 gmx_mm256_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
707 f+i_coord_offset,fshift+i_shift_offset);
709 /* Increment number of inner iterations */
710 inneriter += j_index_end - j_index_start;
712 /* Outer loop uses 6 flops */
715 /* Increment number of outer iterations */
718 /* Update outer/inner flops */
720 inc_nrnb(nrnb,eNR_NBKERNEL_VDW_F,outeriter*6 + inneriter*47);