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 sse2_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_sse2_single.h"
50 #include "kernelutil_x86_sse2_single.h"
53 * Gromacs nonbonded kernel: nb_kernel_ElecNone_VdwLJEwSh_GeomP1P1_VF_sse2_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_sse2_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,
225 c6grid_00 = gmx_mm_load_4real_swizzle_ps(vdwgridparam+vdwioffset0+vdwjidx0A,
226 vdwgridparam+vdwioffset0+vdwjidx0B,
227 vdwgridparam+vdwioffset0+vdwjidx0C,
228 vdwgridparam+vdwioffset0+vdwjidx0D);
230 /* Analytical LJ-PME */
231 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
232 ewcljrsq = _mm_mul_ps(ewclj2,rsq00);
233 ewclj6 = _mm_mul_ps(ewclj2,_mm_mul_ps(ewclj2,ewclj2));
234 exponent = gmx_simd_exp_r(ewcljrsq);
235 /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
236 poly = _mm_mul_ps(exponent,_mm_add_ps(_mm_sub_ps(one,ewcljrsq),_mm_mul_ps(_mm_mul_ps(ewcljrsq,ewcljrsq),one_half)));
237 /* vvdw6 = [C6 - C6grid * (1-poly)]/r6 */
238 vvdw6 = _mm_mul_ps(_mm_sub_ps(c6_00,_mm_mul_ps(c6grid_00,_mm_sub_ps(one,poly))),rinvsix);
239 vvdw12 = _mm_mul_ps(c12_00,_mm_mul_ps(rinvsix,rinvsix));
240 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) ,
241 _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));
242 /* fvdw = vvdw12/r - (vvdw6/r + (C6grid * exponent * beta^6)/r) */
243 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);
245 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
247 /* Update potential sum for this i atom from the interaction with this j atom. */
248 vvdw = _mm_and_ps(vvdw,cutoff_mask);
249 vvdwsum = _mm_add_ps(vvdwsum,vvdw);
253 fscal = _mm_and_ps(fscal,cutoff_mask);
255 /* Calculate temporary vectorial force */
256 tx = _mm_mul_ps(fscal,dx00);
257 ty = _mm_mul_ps(fscal,dy00);
258 tz = _mm_mul_ps(fscal,dz00);
260 /* Update vectorial force */
261 fix0 = _mm_add_ps(fix0,tx);
262 fiy0 = _mm_add_ps(fiy0,ty);
263 fiz0 = _mm_add_ps(fiz0,tz);
265 fjptrA = f+j_coord_offsetA;
266 fjptrB = f+j_coord_offsetB;
267 fjptrC = f+j_coord_offsetC;
268 fjptrD = f+j_coord_offsetD;
269 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
273 /* Inner loop uses 62 flops */
279 /* Get j neighbor index, and coordinate index */
280 jnrlistA = jjnr[jidx];
281 jnrlistB = jjnr[jidx+1];
282 jnrlistC = jjnr[jidx+2];
283 jnrlistD = jjnr[jidx+3];
284 /* Sign of each element will be negative for non-real atoms.
285 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
286 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
288 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
289 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
290 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
291 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
292 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
293 j_coord_offsetA = DIM*jnrA;
294 j_coord_offsetB = DIM*jnrB;
295 j_coord_offsetC = DIM*jnrC;
296 j_coord_offsetD = DIM*jnrD;
298 /* load j atom coordinates */
299 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
300 x+j_coord_offsetC,x+j_coord_offsetD,
303 /* Calculate displacement vector */
304 dx00 = _mm_sub_ps(ix0,jx0);
305 dy00 = _mm_sub_ps(iy0,jy0);
306 dz00 = _mm_sub_ps(iz0,jz0);
308 /* Calculate squared distance and things based on it */
309 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
311 rinv00 = gmx_mm_invsqrt_ps(rsq00);
313 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
315 /* Load parameters for j particles */
316 vdwjidx0A = 2*vdwtype[jnrA+0];
317 vdwjidx0B = 2*vdwtype[jnrB+0];
318 vdwjidx0C = 2*vdwtype[jnrC+0];
319 vdwjidx0D = 2*vdwtype[jnrD+0];
321 /**************************
322 * CALCULATE INTERACTIONS *
323 **************************/
325 if (gmx_mm_any_lt(rsq00,rcutoff2))
328 r00 = _mm_mul_ps(rsq00,rinv00);
329 r00 = _mm_andnot_ps(dummy_mask,r00);
331 /* Compute parameters for interactions between i and j atoms */
332 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
333 vdwparam+vdwioffset0+vdwjidx0B,
334 vdwparam+vdwioffset0+vdwjidx0C,
335 vdwparam+vdwioffset0+vdwjidx0D,
337 c6grid_00 = gmx_mm_load_4real_swizzle_ps(vdwgridparam+vdwioffset0+vdwjidx0A,
338 vdwgridparam+vdwioffset0+vdwjidx0B,
339 vdwgridparam+vdwioffset0+vdwjidx0C,
340 vdwgridparam+vdwioffset0+vdwjidx0D);
342 /* Analytical LJ-PME */
343 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
344 ewcljrsq = _mm_mul_ps(ewclj2,rsq00);
345 ewclj6 = _mm_mul_ps(ewclj2,_mm_mul_ps(ewclj2,ewclj2));
346 exponent = gmx_simd_exp_r(ewcljrsq);
347 /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
348 poly = _mm_mul_ps(exponent,_mm_add_ps(_mm_sub_ps(one,ewcljrsq),_mm_mul_ps(_mm_mul_ps(ewcljrsq,ewcljrsq),one_half)));
349 /* vvdw6 = [C6 - C6grid * (1-poly)]/r6 */
350 vvdw6 = _mm_mul_ps(_mm_sub_ps(c6_00,_mm_mul_ps(c6grid_00,_mm_sub_ps(one,poly))),rinvsix);
351 vvdw12 = _mm_mul_ps(c12_00,_mm_mul_ps(rinvsix,rinvsix));
352 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) ,
353 _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));
354 /* fvdw = vvdw12/r - (vvdw6/r + (C6grid * exponent * beta^6)/r) */
355 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);
357 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
359 /* Update potential sum for this i atom from the interaction with this j atom. */
360 vvdw = _mm_and_ps(vvdw,cutoff_mask);
361 vvdw = _mm_andnot_ps(dummy_mask,vvdw);
362 vvdwsum = _mm_add_ps(vvdwsum,vvdw);
366 fscal = _mm_and_ps(fscal,cutoff_mask);
368 fscal = _mm_andnot_ps(dummy_mask,fscal);
370 /* Calculate temporary vectorial force */
371 tx = _mm_mul_ps(fscal,dx00);
372 ty = _mm_mul_ps(fscal,dy00);
373 tz = _mm_mul_ps(fscal,dz00);
375 /* Update vectorial force */
376 fix0 = _mm_add_ps(fix0,tx);
377 fiy0 = _mm_add_ps(fiy0,ty);
378 fiz0 = _mm_add_ps(fiz0,tz);
380 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
381 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
382 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
383 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
384 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
388 /* Inner loop uses 63 flops */
391 /* End of innermost loop */
393 gmx_mm_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
394 f+i_coord_offset,fshift+i_shift_offset);
397 /* Update potential energies */
398 gmx_mm_update_1pot_ps(vvdwsum,kernel_data->energygrp_vdw+ggid);
400 /* Increment number of inner iterations */
401 inneriter += j_index_end - j_index_start;
403 /* Outer loop uses 7 flops */
406 /* Increment number of outer iterations */
409 /* Update outer/inner flops */
411 inc_nrnb(nrnb,eNR_NBKERNEL_VDW_VF,outeriter*7 + inneriter*63);
414 * Gromacs nonbonded kernel: nb_kernel_ElecNone_VdwLJEwSh_GeomP1P1_F_sse2_single
415 * Electrostatics interaction: None
416 * VdW interaction: LJEwald
417 * Geometry: Particle-Particle
418 * Calculate force/pot: Force
421 nb_kernel_ElecNone_VdwLJEwSh_GeomP1P1_F_sse2_single
422 (t_nblist * gmx_restrict nlist,
423 rvec * gmx_restrict xx,
424 rvec * gmx_restrict ff,
425 t_forcerec * gmx_restrict fr,
426 t_mdatoms * gmx_restrict mdatoms,
427 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
428 t_nrnb * gmx_restrict nrnb)
430 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
431 * just 0 for non-waters.
432 * Suffixes A,B,C,D refer to j loop unrolling done with SSE, e.g. for the four different
433 * jnr indices corresponding to data put in the four positions in the SIMD register.
435 int i_shift_offset,i_coord_offset,outeriter,inneriter;
436 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
437 int jnrA,jnrB,jnrC,jnrD;
438 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
439 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
440 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
442 real *shiftvec,*fshift,*x,*f;
443 real *fjptrA,*fjptrB,*fjptrC,*fjptrD;
445 __m128 tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
447 __m128 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
448 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
449 __m128 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
450 __m128 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
452 __m128 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
455 __m128 one_sixth = _mm_set1_ps(1.0/6.0);
456 __m128 one_twelfth = _mm_set1_ps(1.0/12.0);
458 __m128 ewclj,ewclj2,ewclj6,ewcljrsq,poly,exponent,f6A,f6B,sh_lj_ewald;
460 __m128 one_half = _mm_set1_ps(0.5);
461 __m128 minus_one = _mm_set1_ps(-1.0);
462 __m128 dummy_mask,cutoff_mask;
463 __m128 signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
464 __m128 one = _mm_set1_ps(1.0);
465 __m128 two = _mm_set1_ps(2.0);
471 jindex = nlist->jindex;
473 shiftidx = nlist->shift;
475 shiftvec = fr->shift_vec[0];
476 fshift = fr->fshift[0];
477 nvdwtype = fr->ntype;
479 vdwtype = mdatoms->typeA;
480 vdwgridparam = fr->ljpme_c6grid;
481 sh_lj_ewald = _mm_set1_ps(fr->ic->sh_lj_ewald);
482 ewclj = _mm_set1_ps(fr->ewaldcoeff_lj);
483 ewclj2 = _mm_mul_ps(minus_one,_mm_mul_ps(ewclj,ewclj));
485 rcutoff_scalar = fr->rvdw;
486 rcutoff = _mm_set1_ps(rcutoff_scalar);
487 rcutoff2 = _mm_mul_ps(rcutoff,rcutoff);
489 sh_vdw_invrcut6 = _mm_set1_ps(fr->ic->sh_invrc6);
490 rvdw = _mm_set1_ps(fr->rvdw);
492 /* Avoid stupid compiler warnings */
493 jnrA = jnrB = jnrC = jnrD = 0;
502 for(iidx=0;iidx<4*DIM;iidx++)
507 /* Start outer loop over neighborlists */
508 for(iidx=0; iidx<nri; iidx++)
510 /* Load shift vector for this list */
511 i_shift_offset = DIM*shiftidx[iidx];
513 /* Load limits for loop over neighbors */
514 j_index_start = jindex[iidx];
515 j_index_end = jindex[iidx+1];
517 /* Get outer coordinate index */
519 i_coord_offset = DIM*inr;
521 /* Load i particle coords and add shift vector */
522 gmx_mm_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
524 fix0 = _mm_setzero_ps();
525 fiy0 = _mm_setzero_ps();
526 fiz0 = _mm_setzero_ps();
528 /* Load parameters for i particles */
529 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
531 /* Start inner kernel loop */
532 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
535 /* Get j neighbor index, and coordinate index */
540 j_coord_offsetA = DIM*jnrA;
541 j_coord_offsetB = DIM*jnrB;
542 j_coord_offsetC = DIM*jnrC;
543 j_coord_offsetD = DIM*jnrD;
545 /* load j atom coordinates */
546 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
547 x+j_coord_offsetC,x+j_coord_offsetD,
550 /* Calculate displacement vector */
551 dx00 = _mm_sub_ps(ix0,jx0);
552 dy00 = _mm_sub_ps(iy0,jy0);
553 dz00 = _mm_sub_ps(iz0,jz0);
555 /* Calculate squared distance and things based on it */
556 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
558 rinv00 = gmx_mm_invsqrt_ps(rsq00);
560 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
562 /* Load parameters for j particles */
563 vdwjidx0A = 2*vdwtype[jnrA+0];
564 vdwjidx0B = 2*vdwtype[jnrB+0];
565 vdwjidx0C = 2*vdwtype[jnrC+0];
566 vdwjidx0D = 2*vdwtype[jnrD+0];
568 /**************************
569 * CALCULATE INTERACTIONS *
570 **************************/
572 if (gmx_mm_any_lt(rsq00,rcutoff2))
575 r00 = _mm_mul_ps(rsq00,rinv00);
577 /* Compute parameters for interactions between i and j atoms */
578 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
579 vdwparam+vdwioffset0+vdwjidx0B,
580 vdwparam+vdwioffset0+vdwjidx0C,
581 vdwparam+vdwioffset0+vdwjidx0D,
583 c6grid_00 = gmx_mm_load_4real_swizzle_ps(vdwgridparam+vdwioffset0+vdwjidx0A,
584 vdwgridparam+vdwioffset0+vdwjidx0B,
585 vdwgridparam+vdwioffset0+vdwjidx0C,
586 vdwgridparam+vdwioffset0+vdwjidx0D);
588 /* Analytical LJ-PME */
589 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
590 ewcljrsq = _mm_mul_ps(ewclj2,rsq00);
591 ewclj6 = _mm_mul_ps(ewclj2,_mm_mul_ps(ewclj2,ewclj2));
592 exponent = gmx_simd_exp_r(ewcljrsq);
593 /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
594 poly = _mm_mul_ps(exponent,_mm_add_ps(_mm_sub_ps(one,ewcljrsq),_mm_mul_ps(_mm_mul_ps(ewcljrsq,ewcljrsq),one_half)));
595 /* f6A = 6 * C6grid * (1 - poly) */
596 f6A = _mm_mul_ps(c6grid_00,_mm_sub_ps(one,poly));
597 /* f6B = C6grid * exponent * beta^6 */
598 f6B = _mm_mul_ps(_mm_mul_ps(c6grid_00,one_sixth),_mm_mul_ps(exponent,ewclj6));
599 /* fvdw = 12*C12/r13 - ((6*C6 - f6A)/r6 + f6B)/r */
600 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);
602 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
606 fscal = _mm_and_ps(fscal,cutoff_mask);
608 /* Calculate temporary vectorial force */
609 tx = _mm_mul_ps(fscal,dx00);
610 ty = _mm_mul_ps(fscal,dy00);
611 tz = _mm_mul_ps(fscal,dz00);
613 /* Update vectorial force */
614 fix0 = _mm_add_ps(fix0,tx);
615 fiy0 = _mm_add_ps(fiy0,ty);
616 fiz0 = _mm_add_ps(fiz0,tz);
618 fjptrA = f+j_coord_offsetA;
619 fjptrB = f+j_coord_offsetB;
620 fjptrC = f+j_coord_offsetC;
621 fjptrD = f+j_coord_offsetD;
622 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
626 /* Inner loop uses 49 flops */
632 /* Get j neighbor index, and coordinate index */
633 jnrlistA = jjnr[jidx];
634 jnrlistB = jjnr[jidx+1];
635 jnrlistC = jjnr[jidx+2];
636 jnrlistD = jjnr[jidx+3];
637 /* Sign of each element will be negative for non-real atoms.
638 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
639 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
641 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
642 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
643 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
644 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
645 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
646 j_coord_offsetA = DIM*jnrA;
647 j_coord_offsetB = DIM*jnrB;
648 j_coord_offsetC = DIM*jnrC;
649 j_coord_offsetD = DIM*jnrD;
651 /* load j atom coordinates */
652 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
653 x+j_coord_offsetC,x+j_coord_offsetD,
656 /* Calculate displacement vector */
657 dx00 = _mm_sub_ps(ix0,jx0);
658 dy00 = _mm_sub_ps(iy0,jy0);
659 dz00 = _mm_sub_ps(iz0,jz0);
661 /* Calculate squared distance and things based on it */
662 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
664 rinv00 = gmx_mm_invsqrt_ps(rsq00);
666 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
668 /* Load parameters for j particles */
669 vdwjidx0A = 2*vdwtype[jnrA+0];
670 vdwjidx0B = 2*vdwtype[jnrB+0];
671 vdwjidx0C = 2*vdwtype[jnrC+0];
672 vdwjidx0D = 2*vdwtype[jnrD+0];
674 /**************************
675 * CALCULATE INTERACTIONS *
676 **************************/
678 if (gmx_mm_any_lt(rsq00,rcutoff2))
681 r00 = _mm_mul_ps(rsq00,rinv00);
682 r00 = _mm_andnot_ps(dummy_mask,r00);
684 /* Compute parameters for interactions between i and j atoms */
685 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
686 vdwparam+vdwioffset0+vdwjidx0B,
687 vdwparam+vdwioffset0+vdwjidx0C,
688 vdwparam+vdwioffset0+vdwjidx0D,
690 c6grid_00 = gmx_mm_load_4real_swizzle_ps(vdwgridparam+vdwioffset0+vdwjidx0A,
691 vdwgridparam+vdwioffset0+vdwjidx0B,
692 vdwgridparam+vdwioffset0+vdwjidx0C,
693 vdwgridparam+vdwioffset0+vdwjidx0D);
695 /* Analytical LJ-PME */
696 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
697 ewcljrsq = _mm_mul_ps(ewclj2,rsq00);
698 ewclj6 = _mm_mul_ps(ewclj2,_mm_mul_ps(ewclj2,ewclj2));
699 exponent = gmx_simd_exp_r(ewcljrsq);
700 /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
701 poly = _mm_mul_ps(exponent,_mm_add_ps(_mm_sub_ps(one,ewcljrsq),_mm_mul_ps(_mm_mul_ps(ewcljrsq,ewcljrsq),one_half)));
702 /* f6A = 6 * C6grid * (1 - poly) */
703 f6A = _mm_mul_ps(c6grid_00,_mm_sub_ps(one,poly));
704 /* f6B = C6grid * exponent * beta^6 */
705 f6B = _mm_mul_ps(_mm_mul_ps(c6grid_00,one_sixth),_mm_mul_ps(exponent,ewclj6));
706 /* fvdw = 12*C12/r13 - ((6*C6 - f6A)/r6 + f6B)/r */
707 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);
709 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
713 fscal = _mm_and_ps(fscal,cutoff_mask);
715 fscal = _mm_andnot_ps(dummy_mask,fscal);
717 /* Calculate temporary vectorial force */
718 tx = _mm_mul_ps(fscal,dx00);
719 ty = _mm_mul_ps(fscal,dy00);
720 tz = _mm_mul_ps(fscal,dz00);
722 /* Update vectorial force */
723 fix0 = _mm_add_ps(fix0,tx);
724 fiy0 = _mm_add_ps(fiy0,ty);
725 fiz0 = _mm_add_ps(fiz0,tz);
727 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
728 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
729 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
730 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
731 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
735 /* Inner loop uses 50 flops */
738 /* End of innermost loop */
740 gmx_mm_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
741 f+i_coord_offset,fshift+i_shift_offset);
743 /* Increment number of inner iterations */
744 inneriter += j_index_end - j_index_start;
746 /* Outer loop uses 6 flops */
749 /* Increment number of outer iterations */
752 /* Update outer/inner flops */
754 inc_nrnb(nrnb,eNR_NBKERNEL_VDW_F,outeriter*6 + inneriter*50);