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_128_fma_single 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_128_fma_single.h"
48 #include "kernelutil_x86_avx_128_fma_single.h"
51 * Gromacs nonbonded kernel: nb_kernel_ElecNone_VdwLJEwSh_GeomP1P1_VF_avx_128_fma_single
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_128_fma_single
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_128, 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 j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
77 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
79 real *shiftvec,*fshift,*x,*f;
80 real *fjptrA,*fjptrB,*fjptrC,*fjptrD;
82 __m128 fscal,rcutoff,rcutoff2,jidxall;
84 __m128 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
85 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
86 __m128 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
87 __m128 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
89 __m128 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
92 __m128 one_sixth = _mm_set1_ps(1.0/6.0);
93 __m128 one_twelfth = _mm_set1_ps(1.0/12.0);
96 __m128 ewclj,ewclj2,ewclj6,ewcljrsq,poly,exponent,f6A,f6B,sh_lj_ewald;
97 __m128 one_half = _mm_set1_ps(0.5);
98 __m128 minus_one = _mm_set1_ps(-1.0);
99 __m128 dummy_mask,cutoff_mask;
100 __m128 signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
101 __m128 one = _mm_set1_ps(1.0);
102 __m128 two = _mm_set1_ps(2.0);
108 jindex = nlist->jindex;
110 shiftidx = nlist->shift;
112 shiftvec = fr->shift_vec[0];
113 fshift = fr->fshift[0];
114 nvdwtype = fr->ntype;
116 vdwtype = mdatoms->typeA;
117 vdwgridparam = fr->ljpme_c6grid;
118 sh_lj_ewald = _mm_set1_ps(fr->ic->sh_lj_ewald);
119 ewclj = _mm_set1_ps(fr->ewaldcoeff_lj);
120 ewclj2 = _mm_mul_ps(minus_one,_mm_mul_ps(ewclj,ewclj));
122 rcutoff_scalar = fr->rvdw;
123 rcutoff = _mm_set1_ps(rcutoff_scalar);
124 rcutoff2 = _mm_mul_ps(rcutoff,rcutoff);
126 sh_vdw_invrcut6 = _mm_set1_ps(fr->ic->sh_invrc6);
127 rvdw = _mm_set1_ps(fr->rvdw);
129 /* Avoid stupid compiler warnings */
130 jnrA = jnrB = jnrC = jnrD = 0;
139 for(iidx=0;iidx<4*DIM;iidx++)
144 /* Start outer loop over neighborlists */
145 for(iidx=0; iidx<nri; iidx++)
147 /* Load shift vector for this list */
148 i_shift_offset = DIM*shiftidx[iidx];
150 /* Load limits for loop over neighbors */
151 j_index_start = jindex[iidx];
152 j_index_end = jindex[iidx+1];
154 /* Get outer coordinate index */
156 i_coord_offset = DIM*inr;
158 /* Load i particle coords and add shift vector */
159 gmx_mm_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
161 fix0 = _mm_setzero_ps();
162 fiy0 = _mm_setzero_ps();
163 fiz0 = _mm_setzero_ps();
165 /* Load parameters for i particles */
166 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
168 /* Reset potential sums */
169 vvdwsum = _mm_setzero_ps();
171 /* Start inner kernel loop */
172 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
175 /* Get j neighbor index, and coordinate index */
180 j_coord_offsetA = DIM*jnrA;
181 j_coord_offsetB = DIM*jnrB;
182 j_coord_offsetC = DIM*jnrC;
183 j_coord_offsetD = DIM*jnrD;
185 /* load j atom coordinates */
186 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
187 x+j_coord_offsetC,x+j_coord_offsetD,
190 /* Calculate displacement vector */
191 dx00 = _mm_sub_ps(ix0,jx0);
192 dy00 = _mm_sub_ps(iy0,jy0);
193 dz00 = _mm_sub_ps(iz0,jz0);
195 /* Calculate squared distance and things based on it */
196 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
198 rinv00 = gmx_mm_invsqrt_ps(rsq00);
200 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
202 /* Load parameters for j particles */
203 vdwjidx0A = 2*vdwtype[jnrA+0];
204 vdwjidx0B = 2*vdwtype[jnrB+0];
205 vdwjidx0C = 2*vdwtype[jnrC+0];
206 vdwjidx0D = 2*vdwtype[jnrD+0];
208 /**************************
209 * CALCULATE INTERACTIONS *
210 **************************/
212 if (gmx_mm_any_lt(rsq00,rcutoff2))
215 r00 = _mm_mul_ps(rsq00,rinv00);
217 /* Compute parameters for interactions between i and j atoms */
218 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
219 vdwparam+vdwioffset0+vdwjidx0B,
220 vdwparam+vdwioffset0+vdwjidx0C,
221 vdwparam+vdwioffset0+vdwjidx0D,
224 c6grid_00 = gmx_mm_load_4real_swizzle_ps(vdwgridparam+vdwioffset0+vdwjidx0A,
225 vdwgridparam+vdwioffset0+vdwjidx0B,
226 vdwgridparam+vdwioffset0+vdwjidx0C,
227 vdwgridparam+vdwioffset0+vdwjidx0D);
229 /* Analytical LJ-PME */
230 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
231 ewcljrsq = _mm_mul_ps(ewclj2,rsq00);
232 ewclj6 = _mm_mul_ps(ewclj2,_mm_mul_ps(ewclj2,ewclj2));
233 exponent = gmx_simd_exp_r(ewcljrsq);
234 /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
235 poly = _mm_mul_ps(exponent,_mm_macc_ps(_mm_mul_ps(ewcljrsq,ewcljrsq),one_half,_mm_sub_ps(one,ewcljrsq)));
236 /* vvdw6 = [C6 - C6grid * (1-poly)]/r6 */
237 vvdw6 = _mm_mul_ps(_mm_macc_ps(-c6grid_00,_mm_sub_ps(one,poly),c6_00),rinvsix);
238 vvdw12 = _mm_mul_ps(c12_00,_mm_mul_ps(rinvsix,rinvsix));
239 vvdw = _mm_msub_ps(_mm_nmacc_ps(c12_00,_mm_mul_ps(sh_vdw_invrcut6,sh_vdw_invrcut6),vvdw12),one_twelfth,
240 _mm_mul_ps(_mm_sub_ps(vvdw6,_mm_macc_ps(c6grid_00,sh_lj_ewald,_mm_mul_ps(c6_00,sh_vdw_invrcut6))),one_sixth));
241 /* fvdw = vvdw12/r - (vvdw6/r + (C6grid * exponent * beta^6)/r) */
242 fvdw = _mm_mul_ps(_mm_add_ps(vvdw12,_mm_msub_ps(_mm_mul_ps(c6grid_00,one_sixth),_mm_mul_ps(exponent,ewclj6),vvdw6)),rinvsq00);
244 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
246 /* Update potential sum for this i atom from the interaction with this j atom. */
247 vvdw = _mm_and_ps(vvdw,cutoff_mask);
248 vvdwsum = _mm_add_ps(vvdwsum,vvdw);
252 fscal = _mm_and_ps(fscal,cutoff_mask);
254 /* Update vectorial force */
255 fix0 = _mm_macc_ps(dx00,fscal,fix0);
256 fiy0 = _mm_macc_ps(dy00,fscal,fiy0);
257 fiz0 = _mm_macc_ps(dz00,fscal,fiz0);
259 fjptrA = f+j_coord_offsetA;
260 fjptrB = f+j_coord_offsetB;
261 fjptrC = f+j_coord_offsetC;
262 fjptrD = f+j_coord_offsetD;
263 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,
264 _mm_mul_ps(dx00,fscal),
265 _mm_mul_ps(dy00,fscal),
266 _mm_mul_ps(dz00,fscal));
270 /* Inner loop uses 59 flops */
276 /* Get j neighbor index, and coordinate index */
277 jnrlistA = jjnr[jidx];
278 jnrlistB = jjnr[jidx+1];
279 jnrlistC = jjnr[jidx+2];
280 jnrlistD = jjnr[jidx+3];
281 /* Sign of each element will be negative for non-real atoms.
282 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
283 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
285 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
286 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
287 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
288 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
289 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
290 j_coord_offsetA = DIM*jnrA;
291 j_coord_offsetB = DIM*jnrB;
292 j_coord_offsetC = DIM*jnrC;
293 j_coord_offsetD = DIM*jnrD;
295 /* load j atom coordinates */
296 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
297 x+j_coord_offsetC,x+j_coord_offsetD,
300 /* Calculate displacement vector */
301 dx00 = _mm_sub_ps(ix0,jx0);
302 dy00 = _mm_sub_ps(iy0,jy0);
303 dz00 = _mm_sub_ps(iz0,jz0);
305 /* Calculate squared distance and things based on it */
306 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
308 rinv00 = gmx_mm_invsqrt_ps(rsq00);
310 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
312 /* Load parameters for j particles */
313 vdwjidx0A = 2*vdwtype[jnrA+0];
314 vdwjidx0B = 2*vdwtype[jnrB+0];
315 vdwjidx0C = 2*vdwtype[jnrC+0];
316 vdwjidx0D = 2*vdwtype[jnrD+0];
318 /**************************
319 * CALCULATE INTERACTIONS *
320 **************************/
322 if (gmx_mm_any_lt(rsq00,rcutoff2))
325 r00 = _mm_mul_ps(rsq00,rinv00);
326 r00 = _mm_andnot_ps(dummy_mask,r00);
328 /* Compute parameters for interactions between i and j atoms */
329 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
330 vdwparam+vdwioffset0+vdwjidx0B,
331 vdwparam+vdwioffset0+vdwjidx0C,
332 vdwparam+vdwioffset0+vdwjidx0D,
335 c6grid_00 = gmx_mm_load_4real_swizzle_ps(vdwgridparam+vdwioffset0+vdwjidx0A,
336 vdwgridparam+vdwioffset0+vdwjidx0B,
337 vdwgridparam+vdwioffset0+vdwjidx0C,
338 vdwgridparam+vdwioffset0+vdwjidx0D);
340 /* Analytical LJ-PME */
341 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
342 ewcljrsq = _mm_mul_ps(ewclj2,rsq00);
343 ewclj6 = _mm_mul_ps(ewclj2,_mm_mul_ps(ewclj2,ewclj2));
344 exponent = gmx_simd_exp_r(ewcljrsq);
345 /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
346 poly = _mm_mul_ps(exponent,_mm_macc_ps(_mm_mul_ps(ewcljrsq,ewcljrsq),one_half,_mm_sub_ps(one,ewcljrsq)));
347 /* vvdw6 = [C6 - C6grid * (1-poly)]/r6 */
348 vvdw6 = _mm_mul_ps(_mm_macc_ps(-c6grid_00,_mm_sub_ps(one,poly),c6_00),rinvsix);
349 vvdw12 = _mm_mul_ps(c12_00,_mm_mul_ps(rinvsix,rinvsix));
350 vvdw = _mm_msub_ps(_mm_nmacc_ps(c12_00,_mm_mul_ps(sh_vdw_invrcut6,sh_vdw_invrcut6),vvdw12),one_twelfth,
351 _mm_mul_ps(_mm_sub_ps(vvdw6,_mm_macc_ps(c6grid_00,sh_lj_ewald,_mm_mul_ps(c6_00,sh_vdw_invrcut6))),one_sixth));
352 /* fvdw = vvdw12/r - (vvdw6/r + (C6grid * exponent * beta^6)/r) */
353 fvdw = _mm_mul_ps(_mm_add_ps(vvdw12,_mm_msub_ps(_mm_mul_ps(c6grid_00,one_sixth),_mm_mul_ps(exponent,ewclj6),vvdw6)),rinvsq00);
355 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
357 /* Update potential sum for this i atom from the interaction with this j atom. */
358 vvdw = _mm_and_ps(vvdw,cutoff_mask);
359 vvdw = _mm_andnot_ps(dummy_mask,vvdw);
360 vvdwsum = _mm_add_ps(vvdwsum,vvdw);
364 fscal = _mm_and_ps(fscal,cutoff_mask);
366 fscal = _mm_andnot_ps(dummy_mask,fscal);
368 /* Update vectorial force */
369 fix0 = _mm_macc_ps(dx00,fscal,fix0);
370 fiy0 = _mm_macc_ps(dy00,fscal,fiy0);
371 fiz0 = _mm_macc_ps(dz00,fscal,fiz0);
373 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
374 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
375 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
376 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
377 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,
378 _mm_mul_ps(dx00,fscal),
379 _mm_mul_ps(dy00,fscal),
380 _mm_mul_ps(dz00,fscal));
384 /* Inner loop uses 60 flops */
387 /* End of innermost loop */
389 gmx_mm_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
390 f+i_coord_offset,fshift+i_shift_offset);
393 /* Update potential energies */
394 gmx_mm_update_1pot_ps(vvdwsum,kernel_data->energygrp_vdw+ggid);
396 /* Increment number of inner iterations */
397 inneriter += j_index_end - j_index_start;
399 /* Outer loop uses 7 flops */
402 /* Increment number of outer iterations */
405 /* Update outer/inner flops */
407 inc_nrnb(nrnb,eNR_NBKERNEL_VDW_VF,outeriter*7 + inneriter*60);
410 * Gromacs nonbonded kernel: nb_kernel_ElecNone_VdwLJEwSh_GeomP1P1_F_avx_128_fma_single
411 * Electrostatics interaction: None
412 * VdW interaction: LJEwald
413 * Geometry: Particle-Particle
414 * Calculate force/pot: Force
417 nb_kernel_ElecNone_VdwLJEwSh_GeomP1P1_F_avx_128_fma_single
418 (t_nblist * gmx_restrict nlist,
419 rvec * gmx_restrict xx,
420 rvec * gmx_restrict ff,
421 t_forcerec * gmx_restrict fr,
422 t_mdatoms * gmx_restrict mdatoms,
423 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
424 t_nrnb * gmx_restrict nrnb)
426 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
427 * just 0 for non-waters.
428 * Suffixes A,B,C,D refer to j loop unrolling done with AVX_128, e.g. for the four different
429 * jnr indices corresponding to data put in the four positions in the SIMD register.
431 int i_shift_offset,i_coord_offset,outeriter,inneriter;
432 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
433 int jnrA,jnrB,jnrC,jnrD;
434 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
435 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
436 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
438 real *shiftvec,*fshift,*x,*f;
439 real *fjptrA,*fjptrB,*fjptrC,*fjptrD;
441 __m128 fscal,rcutoff,rcutoff2,jidxall;
443 __m128 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
444 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
445 __m128 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
446 __m128 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
448 __m128 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
451 __m128 one_sixth = _mm_set1_ps(1.0/6.0);
452 __m128 one_twelfth = _mm_set1_ps(1.0/12.0);
455 __m128 ewclj,ewclj2,ewclj6,ewcljrsq,poly,exponent,f6A,f6B,sh_lj_ewald;
456 __m128 one_half = _mm_set1_ps(0.5);
457 __m128 minus_one = _mm_set1_ps(-1.0);
458 __m128 dummy_mask,cutoff_mask;
459 __m128 signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
460 __m128 one = _mm_set1_ps(1.0);
461 __m128 two = _mm_set1_ps(2.0);
467 jindex = nlist->jindex;
469 shiftidx = nlist->shift;
471 shiftvec = fr->shift_vec[0];
472 fshift = fr->fshift[0];
473 nvdwtype = fr->ntype;
475 vdwtype = mdatoms->typeA;
476 vdwgridparam = fr->ljpme_c6grid;
477 sh_lj_ewald = _mm_set1_ps(fr->ic->sh_lj_ewald);
478 ewclj = _mm_set1_ps(fr->ewaldcoeff_lj);
479 ewclj2 = _mm_mul_ps(minus_one,_mm_mul_ps(ewclj,ewclj));
481 rcutoff_scalar = fr->rvdw;
482 rcutoff = _mm_set1_ps(rcutoff_scalar);
483 rcutoff2 = _mm_mul_ps(rcutoff,rcutoff);
485 sh_vdw_invrcut6 = _mm_set1_ps(fr->ic->sh_invrc6);
486 rvdw = _mm_set1_ps(fr->rvdw);
488 /* Avoid stupid compiler warnings */
489 jnrA = jnrB = jnrC = jnrD = 0;
498 for(iidx=0;iidx<4*DIM;iidx++)
503 /* Start outer loop over neighborlists */
504 for(iidx=0; iidx<nri; iidx++)
506 /* Load shift vector for this list */
507 i_shift_offset = DIM*shiftidx[iidx];
509 /* Load limits for loop over neighbors */
510 j_index_start = jindex[iidx];
511 j_index_end = jindex[iidx+1];
513 /* Get outer coordinate index */
515 i_coord_offset = DIM*inr;
517 /* Load i particle coords and add shift vector */
518 gmx_mm_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
520 fix0 = _mm_setzero_ps();
521 fiy0 = _mm_setzero_ps();
522 fiz0 = _mm_setzero_ps();
524 /* Load parameters for i particles */
525 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
527 /* Start inner kernel loop */
528 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
531 /* Get j neighbor index, and coordinate index */
536 j_coord_offsetA = DIM*jnrA;
537 j_coord_offsetB = DIM*jnrB;
538 j_coord_offsetC = DIM*jnrC;
539 j_coord_offsetD = DIM*jnrD;
541 /* load j atom coordinates */
542 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
543 x+j_coord_offsetC,x+j_coord_offsetD,
546 /* Calculate displacement vector */
547 dx00 = _mm_sub_ps(ix0,jx0);
548 dy00 = _mm_sub_ps(iy0,jy0);
549 dz00 = _mm_sub_ps(iz0,jz0);
551 /* Calculate squared distance and things based on it */
552 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
554 rinv00 = gmx_mm_invsqrt_ps(rsq00);
556 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
558 /* Load parameters for j particles */
559 vdwjidx0A = 2*vdwtype[jnrA+0];
560 vdwjidx0B = 2*vdwtype[jnrB+0];
561 vdwjidx0C = 2*vdwtype[jnrC+0];
562 vdwjidx0D = 2*vdwtype[jnrD+0];
564 /**************************
565 * CALCULATE INTERACTIONS *
566 **************************/
568 if (gmx_mm_any_lt(rsq00,rcutoff2))
571 r00 = _mm_mul_ps(rsq00,rinv00);
573 /* Compute parameters for interactions between i and j atoms */
574 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
575 vdwparam+vdwioffset0+vdwjidx0B,
576 vdwparam+vdwioffset0+vdwjidx0C,
577 vdwparam+vdwioffset0+vdwjidx0D,
580 c6grid_00 = gmx_mm_load_4real_swizzle_ps(vdwgridparam+vdwioffset0+vdwjidx0A,
581 vdwgridparam+vdwioffset0+vdwjidx0B,
582 vdwgridparam+vdwioffset0+vdwjidx0C,
583 vdwgridparam+vdwioffset0+vdwjidx0D);
585 /* Analytical LJ-PME */
586 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
587 ewcljrsq = _mm_mul_ps(ewclj2,rsq00);
588 ewclj6 = _mm_mul_ps(ewclj2,_mm_mul_ps(ewclj2,ewclj2));
589 exponent = gmx_simd_exp_r(ewcljrsq);
590 /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
591 poly = _mm_mul_ps(exponent,_mm_macc_ps(_mm_mul_ps(ewcljrsq,ewcljrsq),one_half,_mm_sub_ps(one,ewcljrsq)));
592 /* f6A = 6 * C6grid * (1 - poly) */
593 f6A = _mm_mul_ps(c6grid_00,_mm_sub_ps(one,poly));
594 /* f6B = C6grid * exponent * beta^6 */
595 f6B = _mm_mul_ps(_mm_mul_ps(c6grid_00,one_sixth),_mm_mul_ps(exponent,ewclj6));
596 /* fvdw = 12*C12/r13 - ((6*C6 - f6A)/r6 + f6B)/r */
597 fvdw = _mm_mul_ps(_mm_macc_ps(_mm_msub_ps(c12_00,rinvsix,_mm_sub_ps(c6_00,f6A)),rinvsix,f6B),rinvsq00);
599 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
603 fscal = _mm_and_ps(fscal,cutoff_mask);
605 /* Update vectorial force */
606 fix0 = _mm_macc_ps(dx00,fscal,fix0);
607 fiy0 = _mm_macc_ps(dy00,fscal,fiy0);
608 fiz0 = _mm_macc_ps(dz00,fscal,fiz0);
610 fjptrA = f+j_coord_offsetA;
611 fjptrB = f+j_coord_offsetB;
612 fjptrC = f+j_coord_offsetC;
613 fjptrD = f+j_coord_offsetD;
614 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,
615 _mm_mul_ps(dx00,fscal),
616 _mm_mul_ps(dy00,fscal),
617 _mm_mul_ps(dz00,fscal));
621 /* Inner loop uses 50 flops */
627 /* Get j neighbor index, and coordinate index */
628 jnrlistA = jjnr[jidx];
629 jnrlistB = jjnr[jidx+1];
630 jnrlistC = jjnr[jidx+2];
631 jnrlistD = jjnr[jidx+3];
632 /* Sign of each element will be negative for non-real atoms.
633 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
634 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
636 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
637 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
638 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
639 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
640 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
641 j_coord_offsetA = DIM*jnrA;
642 j_coord_offsetB = DIM*jnrB;
643 j_coord_offsetC = DIM*jnrC;
644 j_coord_offsetD = DIM*jnrD;
646 /* load j atom coordinates */
647 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
648 x+j_coord_offsetC,x+j_coord_offsetD,
651 /* Calculate displacement vector */
652 dx00 = _mm_sub_ps(ix0,jx0);
653 dy00 = _mm_sub_ps(iy0,jy0);
654 dz00 = _mm_sub_ps(iz0,jz0);
656 /* Calculate squared distance and things based on it */
657 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
659 rinv00 = gmx_mm_invsqrt_ps(rsq00);
661 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
663 /* Load parameters for j particles */
664 vdwjidx0A = 2*vdwtype[jnrA+0];
665 vdwjidx0B = 2*vdwtype[jnrB+0];
666 vdwjidx0C = 2*vdwtype[jnrC+0];
667 vdwjidx0D = 2*vdwtype[jnrD+0];
669 /**************************
670 * CALCULATE INTERACTIONS *
671 **************************/
673 if (gmx_mm_any_lt(rsq00,rcutoff2))
676 r00 = _mm_mul_ps(rsq00,rinv00);
677 r00 = _mm_andnot_ps(dummy_mask,r00);
679 /* Compute parameters for interactions between i and j atoms */
680 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
681 vdwparam+vdwioffset0+vdwjidx0B,
682 vdwparam+vdwioffset0+vdwjidx0C,
683 vdwparam+vdwioffset0+vdwjidx0D,
686 c6grid_00 = gmx_mm_load_4real_swizzle_ps(vdwgridparam+vdwioffset0+vdwjidx0A,
687 vdwgridparam+vdwioffset0+vdwjidx0B,
688 vdwgridparam+vdwioffset0+vdwjidx0C,
689 vdwgridparam+vdwioffset0+vdwjidx0D);
691 /* Analytical LJ-PME */
692 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
693 ewcljrsq = _mm_mul_ps(ewclj2,rsq00);
694 ewclj6 = _mm_mul_ps(ewclj2,_mm_mul_ps(ewclj2,ewclj2));
695 exponent = gmx_simd_exp_r(ewcljrsq);
696 /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
697 poly = _mm_mul_ps(exponent,_mm_macc_ps(_mm_mul_ps(ewcljrsq,ewcljrsq),one_half,_mm_sub_ps(one,ewcljrsq)));
698 /* f6A = 6 * C6grid * (1 - poly) */
699 f6A = _mm_mul_ps(c6grid_00,_mm_sub_ps(one,poly));
700 /* f6B = C6grid * exponent * beta^6 */
701 f6B = _mm_mul_ps(_mm_mul_ps(c6grid_00,one_sixth),_mm_mul_ps(exponent,ewclj6));
702 /* fvdw = 12*C12/r13 - ((6*C6 - f6A)/r6 + f6B)/r */
703 fvdw = _mm_mul_ps(_mm_macc_ps(_mm_msub_ps(c12_00,rinvsix,_mm_sub_ps(c6_00,f6A)),rinvsix,f6B),rinvsq00);
705 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
709 fscal = _mm_and_ps(fscal,cutoff_mask);
711 fscal = _mm_andnot_ps(dummy_mask,fscal);
713 /* Update vectorial force */
714 fix0 = _mm_macc_ps(dx00,fscal,fix0);
715 fiy0 = _mm_macc_ps(dy00,fscal,fiy0);
716 fiz0 = _mm_macc_ps(dz00,fscal,fiz0);
718 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
719 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
720 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
721 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
722 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,
723 _mm_mul_ps(dx00,fscal),
724 _mm_mul_ps(dy00,fscal),
725 _mm_mul_ps(dz00,fscal));
729 /* Inner loop uses 51 flops */
732 /* End of innermost loop */
734 gmx_mm_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
735 f+i_coord_offset,fshift+i_shift_offset);
737 /* Increment number of inner iterations */
738 inneriter += j_index_end - j_index_start;
740 /* Outer loop uses 6 flops */
743 /* Increment number of outer iterations */
746 /* Update outer/inner flops */
748 inc_nrnb(nrnb,eNR_NBKERNEL_VDW_F,outeriter*6 + inneriter*51);