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 sse4_1_single kernel generator.
44 #include "../nb_kernel.h"
45 #include "gromacs/legacyheaders/types/simple.h"
46 #include "gromacs/math/vec.h"
47 #include "gromacs/legacyheaders/nrnb.h"
49 #include "gromacs/simd/math_x86_sse4_1_single.h"
50 #include "kernelutil_x86_sse4_1_single.h"
53 * Gromacs nonbonded kernel: nb_kernel_ElecNone_VdwLJEwSh_GeomP1P1_VF_sse4_1_single
54 * Electrostatics interaction: None
55 * VdW interaction: LJEwald
56 * Geometry: Particle-Particle
57 * Calculate force/pot: PotentialAndForce
60 nb_kernel_ElecNone_VdwLJEwSh_GeomP1P1_VF_sse4_1_single
61 (t_nblist * gmx_restrict nlist,
62 rvec * gmx_restrict xx,
63 rvec * gmx_restrict ff,
64 t_forcerec * gmx_restrict fr,
65 t_mdatoms * gmx_restrict mdatoms,
66 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
67 t_nrnb * gmx_restrict nrnb)
69 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
70 * just 0 for non-waters.
71 * Suffixes A,B,C,D refer to j loop unrolling done with SSE, e.g. for the four different
72 * jnr indices corresponding to data put in the four positions in the SIMD register.
74 int i_shift_offset,i_coord_offset,outeriter,inneriter;
75 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
76 int jnrA,jnrB,jnrC,jnrD;
77 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
78 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
79 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
81 real *shiftvec,*fshift,*x,*f;
82 real *fjptrA,*fjptrB,*fjptrC,*fjptrD;
84 __m128 tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
86 __m128 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
87 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
88 __m128 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
89 __m128 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
91 __m128 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
94 __m128 one_sixth = _mm_set1_ps(1.0/6.0);
95 __m128 one_twelfth = _mm_set1_ps(1.0/12.0);
97 __m128 ewclj,ewclj2,ewclj6,ewcljrsq,poly,exponent,f6A,f6B,sh_lj_ewald;
99 __m128 one_half = _mm_set1_ps(0.5);
100 __m128 minus_one = _mm_set1_ps(-1.0);
101 __m128 dummy_mask,cutoff_mask;
102 __m128 signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
103 __m128 one = _mm_set1_ps(1.0);
104 __m128 two = _mm_set1_ps(2.0);
110 jindex = nlist->jindex;
112 shiftidx = nlist->shift;
114 shiftvec = fr->shift_vec[0];
115 fshift = fr->fshift[0];
116 nvdwtype = fr->ntype;
118 vdwtype = mdatoms->typeA;
119 vdwgridparam = fr->ljpme_c6grid;
120 sh_lj_ewald = _mm_set1_ps(fr->ic->sh_lj_ewald);
121 ewclj = _mm_set1_ps(fr->ewaldcoeff_lj);
122 ewclj2 = _mm_mul_ps(minus_one,_mm_mul_ps(ewclj,ewclj));
124 rcutoff_scalar = fr->rvdw;
125 rcutoff = _mm_set1_ps(rcutoff_scalar);
126 rcutoff2 = _mm_mul_ps(rcutoff,rcutoff);
128 sh_vdw_invrcut6 = _mm_set1_ps(fr->ic->sh_invrc6);
129 rvdw = _mm_set1_ps(fr->rvdw);
131 /* Avoid stupid compiler warnings */
132 jnrA = jnrB = jnrC = jnrD = 0;
141 for(iidx=0;iidx<4*DIM;iidx++)
146 /* Start outer loop over neighborlists */
147 for(iidx=0; iidx<nri; iidx++)
149 /* Load shift vector for this list */
150 i_shift_offset = DIM*shiftidx[iidx];
152 /* Load limits for loop over neighbors */
153 j_index_start = jindex[iidx];
154 j_index_end = jindex[iidx+1];
156 /* Get outer coordinate index */
158 i_coord_offset = DIM*inr;
160 /* Load i particle coords and add shift vector */
161 gmx_mm_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
163 fix0 = _mm_setzero_ps();
164 fiy0 = _mm_setzero_ps();
165 fiz0 = _mm_setzero_ps();
167 /* Load parameters for i particles */
168 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
170 /* Reset potential sums */
171 vvdwsum = _mm_setzero_ps();
173 /* Start inner kernel loop */
174 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
177 /* Get j neighbor index, and coordinate index */
182 j_coord_offsetA = DIM*jnrA;
183 j_coord_offsetB = DIM*jnrB;
184 j_coord_offsetC = DIM*jnrC;
185 j_coord_offsetD = DIM*jnrD;
187 /* load j atom coordinates */
188 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
189 x+j_coord_offsetC,x+j_coord_offsetD,
192 /* Calculate displacement vector */
193 dx00 = _mm_sub_ps(ix0,jx0);
194 dy00 = _mm_sub_ps(iy0,jy0);
195 dz00 = _mm_sub_ps(iz0,jz0);
197 /* Calculate squared distance and things based on it */
198 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
200 rinv00 = gmx_mm_invsqrt_ps(rsq00);
202 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
204 /* Load parameters for j particles */
205 vdwjidx0A = 2*vdwtype[jnrA+0];
206 vdwjidx0B = 2*vdwtype[jnrB+0];
207 vdwjidx0C = 2*vdwtype[jnrC+0];
208 vdwjidx0D = 2*vdwtype[jnrD+0];
210 /**************************
211 * CALCULATE INTERACTIONS *
212 **************************/
214 if (gmx_mm_any_lt(rsq00,rcutoff2))
217 r00 = _mm_mul_ps(rsq00,rinv00);
219 /* Compute parameters for interactions between i and j atoms */
220 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
221 vdwparam+vdwioffset0+vdwjidx0B,
222 vdwparam+vdwioffset0+vdwjidx0C,
223 vdwparam+vdwioffset0+vdwjidx0D,
226 c6grid_00 = gmx_mm_load_4real_swizzle_ps(vdwgridparam+vdwioffset0+vdwjidx0A,
227 vdwgridparam+vdwioffset0+vdwjidx0B,
228 vdwgridparam+vdwioffset0+vdwjidx0C,
229 vdwgridparam+vdwioffset0+vdwjidx0D);
231 /* Analytical LJ-PME */
232 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
233 ewcljrsq = _mm_mul_ps(ewclj2,rsq00);
234 ewclj6 = _mm_mul_ps(ewclj2,_mm_mul_ps(ewclj2,ewclj2));
235 exponent = gmx_simd_exp_r(ewcljrsq);
236 /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
237 poly = _mm_mul_ps(exponent,_mm_add_ps(_mm_sub_ps(one,ewcljrsq),_mm_mul_ps(_mm_mul_ps(ewcljrsq,ewcljrsq),one_half)));
238 /* vvdw6 = [C6 - C6grid * (1-poly)]/r6 */
239 vvdw6 = _mm_mul_ps(_mm_sub_ps(c6_00,_mm_mul_ps(c6grid_00,_mm_sub_ps(one,poly))),rinvsix);
240 vvdw12 = _mm_mul_ps(c12_00,_mm_mul_ps(rinvsix,rinvsix));
241 vvdw = _mm_sub_ps(_mm_mul_ps( _mm_sub_ps(vvdw12 , _mm_mul_ps(c12_00,_mm_mul_ps(sh_vdw_invrcut6,sh_vdw_invrcut6))),one_twelfth),
242 _mm_mul_ps( _mm_sub_ps(vvdw6,_mm_add_ps(_mm_mul_ps(c6_00,sh_vdw_invrcut6),_mm_mul_ps(c6grid_00,sh_lj_ewald))),one_sixth));
243 /* fvdw = vvdw12/r - (vvdw6/r + (C6grid * exponent * beta^6)/r) */
244 fvdw = _mm_mul_ps(_mm_sub_ps(vvdw12,_mm_sub_ps(vvdw6,_mm_mul_ps(_mm_mul_ps(c6grid_00,one_sixth),_mm_mul_ps(exponent,ewclj6)))),rinvsq00);
246 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
248 /* Update potential sum for this i atom from the interaction with this j atom. */
249 vvdw = _mm_and_ps(vvdw,cutoff_mask);
250 vvdwsum = _mm_add_ps(vvdwsum,vvdw);
254 fscal = _mm_and_ps(fscal,cutoff_mask);
256 /* Calculate temporary vectorial force */
257 tx = _mm_mul_ps(fscal,dx00);
258 ty = _mm_mul_ps(fscal,dy00);
259 tz = _mm_mul_ps(fscal,dz00);
261 /* Update vectorial force */
262 fix0 = _mm_add_ps(fix0,tx);
263 fiy0 = _mm_add_ps(fiy0,ty);
264 fiz0 = _mm_add_ps(fiz0,tz);
266 fjptrA = f+j_coord_offsetA;
267 fjptrB = f+j_coord_offsetB;
268 fjptrC = f+j_coord_offsetC;
269 fjptrD = f+j_coord_offsetD;
270 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
274 /* Inner loop uses 62 flops */
280 /* Get j neighbor index, and coordinate index */
281 jnrlistA = jjnr[jidx];
282 jnrlistB = jjnr[jidx+1];
283 jnrlistC = jjnr[jidx+2];
284 jnrlistD = jjnr[jidx+3];
285 /* Sign of each element will be negative for non-real atoms.
286 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
287 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
289 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
290 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
291 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
292 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
293 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
294 j_coord_offsetA = DIM*jnrA;
295 j_coord_offsetB = DIM*jnrB;
296 j_coord_offsetC = DIM*jnrC;
297 j_coord_offsetD = DIM*jnrD;
299 /* load j atom coordinates */
300 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
301 x+j_coord_offsetC,x+j_coord_offsetD,
304 /* Calculate displacement vector */
305 dx00 = _mm_sub_ps(ix0,jx0);
306 dy00 = _mm_sub_ps(iy0,jy0);
307 dz00 = _mm_sub_ps(iz0,jz0);
309 /* Calculate squared distance and things based on it */
310 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
312 rinv00 = gmx_mm_invsqrt_ps(rsq00);
314 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
316 /* Load parameters for j particles */
317 vdwjidx0A = 2*vdwtype[jnrA+0];
318 vdwjidx0B = 2*vdwtype[jnrB+0];
319 vdwjidx0C = 2*vdwtype[jnrC+0];
320 vdwjidx0D = 2*vdwtype[jnrD+0];
322 /**************************
323 * CALCULATE INTERACTIONS *
324 **************************/
326 if (gmx_mm_any_lt(rsq00,rcutoff2))
329 r00 = _mm_mul_ps(rsq00,rinv00);
330 r00 = _mm_andnot_ps(dummy_mask,r00);
332 /* Compute parameters for interactions between i and j atoms */
333 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
334 vdwparam+vdwioffset0+vdwjidx0B,
335 vdwparam+vdwioffset0+vdwjidx0C,
336 vdwparam+vdwioffset0+vdwjidx0D,
339 c6grid_00 = gmx_mm_load_4real_swizzle_ps(vdwgridparam+vdwioffset0+vdwjidx0A,
340 vdwgridparam+vdwioffset0+vdwjidx0B,
341 vdwgridparam+vdwioffset0+vdwjidx0C,
342 vdwgridparam+vdwioffset0+vdwjidx0D);
344 /* Analytical LJ-PME */
345 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
346 ewcljrsq = _mm_mul_ps(ewclj2,rsq00);
347 ewclj6 = _mm_mul_ps(ewclj2,_mm_mul_ps(ewclj2,ewclj2));
348 exponent = gmx_simd_exp_r(ewcljrsq);
349 /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
350 poly = _mm_mul_ps(exponent,_mm_add_ps(_mm_sub_ps(one,ewcljrsq),_mm_mul_ps(_mm_mul_ps(ewcljrsq,ewcljrsq),one_half)));
351 /* vvdw6 = [C6 - C6grid * (1-poly)]/r6 */
352 vvdw6 = _mm_mul_ps(_mm_sub_ps(c6_00,_mm_mul_ps(c6grid_00,_mm_sub_ps(one,poly))),rinvsix);
353 vvdw12 = _mm_mul_ps(c12_00,_mm_mul_ps(rinvsix,rinvsix));
354 vvdw = _mm_sub_ps(_mm_mul_ps( _mm_sub_ps(vvdw12 , _mm_mul_ps(c12_00,_mm_mul_ps(sh_vdw_invrcut6,sh_vdw_invrcut6))),one_twelfth),
355 _mm_mul_ps( _mm_sub_ps(vvdw6,_mm_add_ps(_mm_mul_ps(c6_00,sh_vdw_invrcut6),_mm_mul_ps(c6grid_00,sh_lj_ewald))),one_sixth));
356 /* fvdw = vvdw12/r - (vvdw6/r + (C6grid * exponent * beta^6)/r) */
357 fvdw = _mm_mul_ps(_mm_sub_ps(vvdw12,_mm_sub_ps(vvdw6,_mm_mul_ps(_mm_mul_ps(c6grid_00,one_sixth),_mm_mul_ps(exponent,ewclj6)))),rinvsq00);
359 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
361 /* Update potential sum for this i atom from the interaction with this j atom. */
362 vvdw = _mm_and_ps(vvdw,cutoff_mask);
363 vvdw = _mm_andnot_ps(dummy_mask,vvdw);
364 vvdwsum = _mm_add_ps(vvdwsum,vvdw);
368 fscal = _mm_and_ps(fscal,cutoff_mask);
370 fscal = _mm_andnot_ps(dummy_mask,fscal);
372 /* Calculate temporary vectorial force */
373 tx = _mm_mul_ps(fscal,dx00);
374 ty = _mm_mul_ps(fscal,dy00);
375 tz = _mm_mul_ps(fscal,dz00);
377 /* Update vectorial force */
378 fix0 = _mm_add_ps(fix0,tx);
379 fiy0 = _mm_add_ps(fiy0,ty);
380 fiz0 = _mm_add_ps(fiz0,tz);
382 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
383 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
384 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
385 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
386 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
390 /* Inner loop uses 63 flops */
393 /* End of innermost loop */
395 gmx_mm_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
396 f+i_coord_offset,fshift+i_shift_offset);
399 /* Update potential energies */
400 gmx_mm_update_1pot_ps(vvdwsum,kernel_data->energygrp_vdw+ggid);
402 /* Increment number of inner iterations */
403 inneriter += j_index_end - j_index_start;
405 /* Outer loop uses 7 flops */
408 /* Increment number of outer iterations */
411 /* Update outer/inner flops */
413 inc_nrnb(nrnb,eNR_NBKERNEL_VDW_VF,outeriter*7 + inneriter*63);
416 * Gromacs nonbonded kernel: nb_kernel_ElecNone_VdwLJEwSh_GeomP1P1_F_sse4_1_single
417 * Electrostatics interaction: None
418 * VdW interaction: LJEwald
419 * Geometry: Particle-Particle
420 * Calculate force/pot: Force
423 nb_kernel_ElecNone_VdwLJEwSh_GeomP1P1_F_sse4_1_single
424 (t_nblist * gmx_restrict nlist,
425 rvec * gmx_restrict xx,
426 rvec * gmx_restrict ff,
427 t_forcerec * gmx_restrict fr,
428 t_mdatoms * gmx_restrict mdatoms,
429 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
430 t_nrnb * gmx_restrict nrnb)
432 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
433 * just 0 for non-waters.
434 * Suffixes A,B,C,D refer to j loop unrolling done with SSE, e.g. for the four different
435 * jnr indices corresponding to data put in the four positions in the SIMD register.
437 int i_shift_offset,i_coord_offset,outeriter,inneriter;
438 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
439 int jnrA,jnrB,jnrC,jnrD;
440 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
441 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
442 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
444 real *shiftvec,*fshift,*x,*f;
445 real *fjptrA,*fjptrB,*fjptrC,*fjptrD;
447 __m128 tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
449 __m128 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
450 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
451 __m128 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
452 __m128 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
454 __m128 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
457 __m128 one_sixth = _mm_set1_ps(1.0/6.0);
458 __m128 one_twelfth = _mm_set1_ps(1.0/12.0);
460 __m128 ewclj,ewclj2,ewclj6,ewcljrsq,poly,exponent,f6A,f6B,sh_lj_ewald;
462 __m128 one_half = _mm_set1_ps(0.5);
463 __m128 minus_one = _mm_set1_ps(-1.0);
464 __m128 dummy_mask,cutoff_mask;
465 __m128 signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
466 __m128 one = _mm_set1_ps(1.0);
467 __m128 two = _mm_set1_ps(2.0);
473 jindex = nlist->jindex;
475 shiftidx = nlist->shift;
477 shiftvec = fr->shift_vec[0];
478 fshift = fr->fshift[0];
479 nvdwtype = fr->ntype;
481 vdwtype = mdatoms->typeA;
482 vdwgridparam = fr->ljpme_c6grid;
483 sh_lj_ewald = _mm_set1_ps(fr->ic->sh_lj_ewald);
484 ewclj = _mm_set1_ps(fr->ewaldcoeff_lj);
485 ewclj2 = _mm_mul_ps(minus_one,_mm_mul_ps(ewclj,ewclj));
487 rcutoff_scalar = fr->rvdw;
488 rcutoff = _mm_set1_ps(rcutoff_scalar);
489 rcutoff2 = _mm_mul_ps(rcutoff,rcutoff);
491 sh_vdw_invrcut6 = _mm_set1_ps(fr->ic->sh_invrc6);
492 rvdw = _mm_set1_ps(fr->rvdw);
494 /* Avoid stupid compiler warnings */
495 jnrA = jnrB = jnrC = jnrD = 0;
504 for(iidx=0;iidx<4*DIM;iidx++)
509 /* Start outer loop over neighborlists */
510 for(iidx=0; iidx<nri; iidx++)
512 /* Load shift vector for this list */
513 i_shift_offset = DIM*shiftidx[iidx];
515 /* Load limits for loop over neighbors */
516 j_index_start = jindex[iidx];
517 j_index_end = jindex[iidx+1];
519 /* Get outer coordinate index */
521 i_coord_offset = DIM*inr;
523 /* Load i particle coords and add shift vector */
524 gmx_mm_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
526 fix0 = _mm_setzero_ps();
527 fiy0 = _mm_setzero_ps();
528 fiz0 = _mm_setzero_ps();
530 /* Load parameters for i particles */
531 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
533 /* Start inner kernel loop */
534 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
537 /* Get j neighbor index, and coordinate index */
542 j_coord_offsetA = DIM*jnrA;
543 j_coord_offsetB = DIM*jnrB;
544 j_coord_offsetC = DIM*jnrC;
545 j_coord_offsetD = DIM*jnrD;
547 /* load j atom coordinates */
548 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
549 x+j_coord_offsetC,x+j_coord_offsetD,
552 /* Calculate displacement vector */
553 dx00 = _mm_sub_ps(ix0,jx0);
554 dy00 = _mm_sub_ps(iy0,jy0);
555 dz00 = _mm_sub_ps(iz0,jz0);
557 /* Calculate squared distance and things based on it */
558 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
560 rinv00 = gmx_mm_invsqrt_ps(rsq00);
562 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
564 /* Load parameters for j particles */
565 vdwjidx0A = 2*vdwtype[jnrA+0];
566 vdwjidx0B = 2*vdwtype[jnrB+0];
567 vdwjidx0C = 2*vdwtype[jnrC+0];
568 vdwjidx0D = 2*vdwtype[jnrD+0];
570 /**************************
571 * CALCULATE INTERACTIONS *
572 **************************/
574 if (gmx_mm_any_lt(rsq00,rcutoff2))
577 r00 = _mm_mul_ps(rsq00,rinv00);
579 /* Compute parameters for interactions between i and j atoms */
580 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
581 vdwparam+vdwioffset0+vdwjidx0B,
582 vdwparam+vdwioffset0+vdwjidx0C,
583 vdwparam+vdwioffset0+vdwjidx0D,
586 c6grid_00 = gmx_mm_load_4real_swizzle_ps(vdwgridparam+vdwioffset0+vdwjidx0A,
587 vdwgridparam+vdwioffset0+vdwjidx0B,
588 vdwgridparam+vdwioffset0+vdwjidx0C,
589 vdwgridparam+vdwioffset0+vdwjidx0D);
591 /* Analytical LJ-PME */
592 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
593 ewcljrsq = _mm_mul_ps(ewclj2,rsq00);
594 ewclj6 = _mm_mul_ps(ewclj2,_mm_mul_ps(ewclj2,ewclj2));
595 exponent = gmx_simd_exp_r(ewcljrsq);
596 /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
597 poly = _mm_mul_ps(exponent,_mm_add_ps(_mm_sub_ps(one,ewcljrsq),_mm_mul_ps(_mm_mul_ps(ewcljrsq,ewcljrsq),one_half)));
598 /* f6A = 6 * C6grid * (1 - poly) */
599 f6A = _mm_mul_ps(c6grid_00,_mm_sub_ps(one,poly));
600 /* f6B = C6grid * exponent * beta^6 */
601 f6B = _mm_mul_ps(_mm_mul_ps(c6grid_00,one_sixth),_mm_mul_ps(exponent,ewclj6));
602 /* fvdw = 12*C12/r13 - ((6*C6 - f6A)/r6 + f6B)/r */
603 fvdw = _mm_mul_ps(_mm_add_ps(_mm_mul_ps(_mm_sub_ps(_mm_mul_ps(c12_00,rinvsix),_mm_sub_ps(c6_00,f6A)),rinvsix),f6B),rinvsq00);
605 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
609 fscal = _mm_and_ps(fscal,cutoff_mask);
611 /* Calculate temporary vectorial force */
612 tx = _mm_mul_ps(fscal,dx00);
613 ty = _mm_mul_ps(fscal,dy00);
614 tz = _mm_mul_ps(fscal,dz00);
616 /* Update vectorial force */
617 fix0 = _mm_add_ps(fix0,tx);
618 fiy0 = _mm_add_ps(fiy0,ty);
619 fiz0 = _mm_add_ps(fiz0,tz);
621 fjptrA = f+j_coord_offsetA;
622 fjptrB = f+j_coord_offsetB;
623 fjptrC = f+j_coord_offsetC;
624 fjptrD = f+j_coord_offsetD;
625 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
629 /* Inner loop uses 49 flops */
635 /* Get j neighbor index, and coordinate index */
636 jnrlistA = jjnr[jidx];
637 jnrlistB = jjnr[jidx+1];
638 jnrlistC = jjnr[jidx+2];
639 jnrlistD = jjnr[jidx+3];
640 /* Sign of each element will be negative for non-real atoms.
641 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
642 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
644 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
645 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
646 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
647 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
648 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
649 j_coord_offsetA = DIM*jnrA;
650 j_coord_offsetB = DIM*jnrB;
651 j_coord_offsetC = DIM*jnrC;
652 j_coord_offsetD = DIM*jnrD;
654 /* load j atom coordinates */
655 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
656 x+j_coord_offsetC,x+j_coord_offsetD,
659 /* Calculate displacement vector */
660 dx00 = _mm_sub_ps(ix0,jx0);
661 dy00 = _mm_sub_ps(iy0,jy0);
662 dz00 = _mm_sub_ps(iz0,jz0);
664 /* Calculate squared distance and things based on it */
665 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
667 rinv00 = gmx_mm_invsqrt_ps(rsq00);
669 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
671 /* Load parameters for j particles */
672 vdwjidx0A = 2*vdwtype[jnrA+0];
673 vdwjidx0B = 2*vdwtype[jnrB+0];
674 vdwjidx0C = 2*vdwtype[jnrC+0];
675 vdwjidx0D = 2*vdwtype[jnrD+0];
677 /**************************
678 * CALCULATE INTERACTIONS *
679 **************************/
681 if (gmx_mm_any_lt(rsq00,rcutoff2))
684 r00 = _mm_mul_ps(rsq00,rinv00);
685 r00 = _mm_andnot_ps(dummy_mask,r00);
687 /* Compute parameters for interactions between i and j atoms */
688 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
689 vdwparam+vdwioffset0+vdwjidx0B,
690 vdwparam+vdwioffset0+vdwjidx0C,
691 vdwparam+vdwioffset0+vdwjidx0D,
694 c6grid_00 = gmx_mm_load_4real_swizzle_ps(vdwgridparam+vdwioffset0+vdwjidx0A,
695 vdwgridparam+vdwioffset0+vdwjidx0B,
696 vdwgridparam+vdwioffset0+vdwjidx0C,
697 vdwgridparam+vdwioffset0+vdwjidx0D);
699 /* Analytical LJ-PME */
700 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
701 ewcljrsq = _mm_mul_ps(ewclj2,rsq00);
702 ewclj6 = _mm_mul_ps(ewclj2,_mm_mul_ps(ewclj2,ewclj2));
703 exponent = gmx_simd_exp_r(ewcljrsq);
704 /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
705 poly = _mm_mul_ps(exponent,_mm_add_ps(_mm_sub_ps(one,ewcljrsq),_mm_mul_ps(_mm_mul_ps(ewcljrsq,ewcljrsq),one_half)));
706 /* f6A = 6 * C6grid * (1 - poly) */
707 f6A = _mm_mul_ps(c6grid_00,_mm_sub_ps(one,poly));
708 /* f6B = C6grid * exponent * beta^6 */
709 f6B = _mm_mul_ps(_mm_mul_ps(c6grid_00,one_sixth),_mm_mul_ps(exponent,ewclj6));
710 /* fvdw = 12*C12/r13 - ((6*C6 - f6A)/r6 + f6B)/r */
711 fvdw = _mm_mul_ps(_mm_add_ps(_mm_mul_ps(_mm_sub_ps(_mm_mul_ps(c12_00,rinvsix),_mm_sub_ps(c6_00,f6A)),rinvsix),f6B),rinvsq00);
713 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
717 fscal = _mm_and_ps(fscal,cutoff_mask);
719 fscal = _mm_andnot_ps(dummy_mask,fscal);
721 /* Calculate temporary vectorial force */
722 tx = _mm_mul_ps(fscal,dx00);
723 ty = _mm_mul_ps(fscal,dy00);
724 tz = _mm_mul_ps(fscal,dz00);
726 /* Update vectorial force */
727 fix0 = _mm_add_ps(fix0,tx);
728 fiy0 = _mm_add_ps(fiy0,ty);
729 fiz0 = _mm_add_ps(fiz0,tz);
731 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
732 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
733 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
734 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
735 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
739 /* Inner loop uses 50 flops */
742 /* End of innermost loop */
744 gmx_mm_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
745 f+i_coord_offset,fshift+i_shift_offset);
747 /* Increment number of inner iterations */
748 inneriter += j_index_end - j_index_start;
750 /* Outer loop uses 6 flops */
753 /* Increment number of outer iterations */
756 /* Update outer/inner flops */
758 inc_nrnb(nrnb,eNR_NBKERNEL_VDW_F,outeriter*6 + inneriter*50);