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.
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_sse2_single.h"
48 #include "kernelutil_x86_sse2_single.h"
51 * Gromacs nonbonded kernel: nb_kernel_ElecNone_VdwCSTab_GeomP1P1_VF_sse2_single
52 * Electrostatics interaction: None
53 * VdW interaction: CubicSplineTable
54 * Geometry: Particle-Particle
55 * Calculate force/pot: PotentialAndForce
58 nb_kernel_ElecNone_VdwCSTab_GeomP1P1_VF_sse2_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 SSE, 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 tx,ty,tz,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);
95 __m128i ifour = _mm_set1_epi32(4);
96 __m128 rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF;
98 __m128 dummy_mask,cutoff_mask;
99 __m128 signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
100 __m128 one = _mm_set1_ps(1.0);
101 __m128 two = _mm_set1_ps(2.0);
107 jindex = nlist->jindex;
109 shiftidx = nlist->shift;
111 shiftvec = fr->shift_vec[0];
112 fshift = fr->fshift[0];
113 nvdwtype = fr->ntype;
115 vdwtype = mdatoms->typeA;
117 vftab = kernel_data->table_vdw->data;
118 vftabscale = _mm_set1_ps(kernel_data->table_vdw->scale);
120 /* Avoid stupid compiler warnings */
121 jnrA = jnrB = jnrC = jnrD = 0;
130 for(iidx=0;iidx<4*DIM;iidx++)
135 /* Start outer loop over neighborlists */
136 for(iidx=0; iidx<nri; iidx++)
138 /* Load shift vector for this list */
139 i_shift_offset = DIM*shiftidx[iidx];
141 /* Load limits for loop over neighbors */
142 j_index_start = jindex[iidx];
143 j_index_end = jindex[iidx+1];
145 /* Get outer coordinate index */
147 i_coord_offset = DIM*inr;
149 /* Load i particle coords and add shift vector */
150 gmx_mm_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
152 fix0 = _mm_setzero_ps();
153 fiy0 = _mm_setzero_ps();
154 fiz0 = _mm_setzero_ps();
156 /* Load parameters for i particles */
157 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
159 /* Reset potential sums */
160 vvdwsum = _mm_setzero_ps();
162 /* Start inner kernel loop */
163 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
166 /* Get j neighbor index, and coordinate index */
171 j_coord_offsetA = DIM*jnrA;
172 j_coord_offsetB = DIM*jnrB;
173 j_coord_offsetC = DIM*jnrC;
174 j_coord_offsetD = DIM*jnrD;
176 /* load j atom coordinates */
177 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
178 x+j_coord_offsetC,x+j_coord_offsetD,
181 /* Calculate displacement vector */
182 dx00 = _mm_sub_ps(ix0,jx0);
183 dy00 = _mm_sub_ps(iy0,jy0);
184 dz00 = _mm_sub_ps(iz0,jz0);
186 /* Calculate squared distance and things based on it */
187 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
189 rinv00 = gmx_mm_invsqrt_ps(rsq00);
191 /* Load parameters for j particles */
192 vdwjidx0A = 2*vdwtype[jnrA+0];
193 vdwjidx0B = 2*vdwtype[jnrB+0];
194 vdwjidx0C = 2*vdwtype[jnrC+0];
195 vdwjidx0D = 2*vdwtype[jnrD+0];
197 /**************************
198 * CALCULATE INTERACTIONS *
199 **************************/
201 r00 = _mm_mul_ps(rsq00,rinv00);
203 /* Compute parameters for interactions between i and j atoms */
204 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
205 vdwparam+vdwioffset0+vdwjidx0B,
206 vdwparam+vdwioffset0+vdwjidx0C,
207 vdwparam+vdwioffset0+vdwjidx0D,
210 /* Calculate table index by multiplying r with table scale and truncate to integer */
211 rt = _mm_mul_ps(r00,vftabscale);
212 vfitab = _mm_cvttps_epi32(rt);
213 vfeps = _mm_sub_ps(rt,_mm_cvtepi32_ps(vfitab));
214 vfitab = _mm_slli_epi32(vfitab,3);
216 /* CUBIC SPLINE TABLE DISPERSION */
217 Y = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,0) );
218 F = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,1) );
219 G = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,2) );
220 H = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,3) );
221 _MM_TRANSPOSE4_PS(Y,F,G,H);
222 Heps = _mm_mul_ps(vfeps,H);
223 Fp = _mm_add_ps(F,_mm_mul_ps(vfeps,_mm_add_ps(G,Heps)));
224 VV = _mm_add_ps(Y,_mm_mul_ps(vfeps,Fp));
225 vvdw6 = _mm_mul_ps(c6_00,VV);
226 FF = _mm_add_ps(Fp,_mm_mul_ps(vfeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
227 fvdw6 = _mm_mul_ps(c6_00,FF);
229 /* CUBIC SPLINE TABLE REPULSION */
230 vfitab = _mm_add_epi32(vfitab,ifour);
231 Y = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,0) );
232 F = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,1) );
233 G = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,2) );
234 H = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,3) );
235 _MM_TRANSPOSE4_PS(Y,F,G,H);
236 Heps = _mm_mul_ps(vfeps,H);
237 Fp = _mm_add_ps(F,_mm_mul_ps(vfeps,_mm_add_ps(G,Heps)));
238 VV = _mm_add_ps(Y,_mm_mul_ps(vfeps,Fp));
239 vvdw12 = _mm_mul_ps(c12_00,VV);
240 FF = _mm_add_ps(Fp,_mm_mul_ps(vfeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
241 fvdw12 = _mm_mul_ps(c12_00,FF);
242 vvdw = _mm_add_ps(vvdw12,vvdw6);
243 fvdw = _mm_xor_ps(signbit,_mm_mul_ps(_mm_add_ps(fvdw6,fvdw12),_mm_mul_ps(vftabscale,rinv00)));
245 /* Update potential sum for this i atom from the interaction with this j atom. */
246 vvdwsum = _mm_add_ps(vvdwsum,vvdw);
250 /* Calculate temporary vectorial force */
251 tx = _mm_mul_ps(fscal,dx00);
252 ty = _mm_mul_ps(fscal,dy00);
253 tz = _mm_mul_ps(fscal,dz00);
255 /* Update vectorial force */
256 fix0 = _mm_add_ps(fix0,tx);
257 fiy0 = _mm_add_ps(fiy0,ty);
258 fiz0 = _mm_add_ps(fiz0,tz);
260 fjptrA = f+j_coord_offsetA;
261 fjptrB = f+j_coord_offsetB;
262 fjptrC = f+j_coord_offsetC;
263 fjptrD = f+j_coord_offsetD;
264 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
266 /* Inner loop uses 56 flops */
272 /* Get j neighbor index, and coordinate index */
273 jnrlistA = jjnr[jidx];
274 jnrlistB = jjnr[jidx+1];
275 jnrlistC = jjnr[jidx+2];
276 jnrlistD = jjnr[jidx+3];
277 /* Sign of each element will be negative for non-real atoms.
278 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
279 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
281 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
282 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
283 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
284 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
285 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
286 j_coord_offsetA = DIM*jnrA;
287 j_coord_offsetB = DIM*jnrB;
288 j_coord_offsetC = DIM*jnrC;
289 j_coord_offsetD = DIM*jnrD;
291 /* load j atom coordinates */
292 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
293 x+j_coord_offsetC,x+j_coord_offsetD,
296 /* Calculate displacement vector */
297 dx00 = _mm_sub_ps(ix0,jx0);
298 dy00 = _mm_sub_ps(iy0,jy0);
299 dz00 = _mm_sub_ps(iz0,jz0);
301 /* Calculate squared distance and things based on it */
302 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
304 rinv00 = gmx_mm_invsqrt_ps(rsq00);
306 /* Load parameters for j particles */
307 vdwjidx0A = 2*vdwtype[jnrA+0];
308 vdwjidx0B = 2*vdwtype[jnrB+0];
309 vdwjidx0C = 2*vdwtype[jnrC+0];
310 vdwjidx0D = 2*vdwtype[jnrD+0];
312 /**************************
313 * CALCULATE INTERACTIONS *
314 **************************/
316 r00 = _mm_mul_ps(rsq00,rinv00);
317 r00 = _mm_andnot_ps(dummy_mask,r00);
319 /* Compute parameters for interactions between i and j atoms */
320 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
321 vdwparam+vdwioffset0+vdwjidx0B,
322 vdwparam+vdwioffset0+vdwjidx0C,
323 vdwparam+vdwioffset0+vdwjidx0D,
326 /* Calculate table index by multiplying r with table scale and truncate to integer */
327 rt = _mm_mul_ps(r00,vftabscale);
328 vfitab = _mm_cvttps_epi32(rt);
329 vfeps = _mm_sub_ps(rt,_mm_cvtepi32_ps(vfitab));
330 vfitab = _mm_slli_epi32(vfitab,3);
332 /* CUBIC SPLINE TABLE DISPERSION */
333 Y = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,0) );
334 F = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,1) );
335 G = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,2) );
336 H = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,3) );
337 _MM_TRANSPOSE4_PS(Y,F,G,H);
338 Heps = _mm_mul_ps(vfeps,H);
339 Fp = _mm_add_ps(F,_mm_mul_ps(vfeps,_mm_add_ps(G,Heps)));
340 VV = _mm_add_ps(Y,_mm_mul_ps(vfeps,Fp));
341 vvdw6 = _mm_mul_ps(c6_00,VV);
342 FF = _mm_add_ps(Fp,_mm_mul_ps(vfeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
343 fvdw6 = _mm_mul_ps(c6_00,FF);
345 /* CUBIC SPLINE TABLE REPULSION */
346 vfitab = _mm_add_epi32(vfitab,ifour);
347 Y = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,0) );
348 F = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,1) );
349 G = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,2) );
350 H = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,3) );
351 _MM_TRANSPOSE4_PS(Y,F,G,H);
352 Heps = _mm_mul_ps(vfeps,H);
353 Fp = _mm_add_ps(F,_mm_mul_ps(vfeps,_mm_add_ps(G,Heps)));
354 VV = _mm_add_ps(Y,_mm_mul_ps(vfeps,Fp));
355 vvdw12 = _mm_mul_ps(c12_00,VV);
356 FF = _mm_add_ps(Fp,_mm_mul_ps(vfeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
357 fvdw12 = _mm_mul_ps(c12_00,FF);
358 vvdw = _mm_add_ps(vvdw12,vvdw6);
359 fvdw = _mm_xor_ps(signbit,_mm_mul_ps(_mm_add_ps(fvdw6,fvdw12),_mm_mul_ps(vftabscale,rinv00)));
361 /* Update potential sum for this i atom from the interaction with this j atom. */
362 vvdw = _mm_andnot_ps(dummy_mask,vvdw);
363 vvdwsum = _mm_add_ps(vvdwsum,vvdw);
367 fscal = _mm_andnot_ps(dummy_mask,fscal);
369 /* Calculate temporary vectorial force */
370 tx = _mm_mul_ps(fscal,dx00);
371 ty = _mm_mul_ps(fscal,dy00);
372 tz = _mm_mul_ps(fscal,dz00);
374 /* Update vectorial force */
375 fix0 = _mm_add_ps(fix0,tx);
376 fiy0 = _mm_add_ps(fiy0,ty);
377 fiz0 = _mm_add_ps(fiz0,tz);
379 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
380 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
381 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
382 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
383 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
385 /* Inner loop uses 57 flops */
388 /* End of innermost loop */
390 gmx_mm_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
391 f+i_coord_offset,fshift+i_shift_offset);
394 /* Update potential energies */
395 gmx_mm_update_1pot_ps(vvdwsum,kernel_data->energygrp_vdw+ggid);
397 /* Increment number of inner iterations */
398 inneriter += j_index_end - j_index_start;
400 /* Outer loop uses 7 flops */
403 /* Increment number of outer iterations */
406 /* Update outer/inner flops */
408 inc_nrnb(nrnb,eNR_NBKERNEL_VDW_VF,outeriter*7 + inneriter*57);
411 * Gromacs nonbonded kernel: nb_kernel_ElecNone_VdwCSTab_GeomP1P1_F_sse2_single
412 * Electrostatics interaction: None
413 * VdW interaction: CubicSplineTable
414 * Geometry: Particle-Particle
415 * Calculate force/pot: Force
418 nb_kernel_ElecNone_VdwCSTab_GeomP1P1_F_sse2_single
419 (t_nblist * gmx_restrict nlist,
420 rvec * gmx_restrict xx,
421 rvec * gmx_restrict ff,
422 t_forcerec * gmx_restrict fr,
423 t_mdatoms * gmx_restrict mdatoms,
424 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
425 t_nrnb * gmx_restrict nrnb)
427 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
428 * just 0 for non-waters.
429 * Suffixes A,B,C,D refer to j loop unrolling done with SSE, e.g. for the four different
430 * jnr indices corresponding to data put in the four positions in the SIMD register.
432 int i_shift_offset,i_coord_offset,outeriter,inneriter;
433 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
434 int jnrA,jnrB,jnrC,jnrD;
435 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
436 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
437 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
439 real *shiftvec,*fshift,*x,*f;
440 real *fjptrA,*fjptrB,*fjptrC,*fjptrD;
442 __m128 tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
444 __m128 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
445 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
446 __m128 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
447 __m128 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
449 __m128 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
452 __m128 one_sixth = _mm_set1_ps(1.0/6.0);
453 __m128 one_twelfth = _mm_set1_ps(1.0/12.0);
455 __m128i ifour = _mm_set1_epi32(4);
456 __m128 rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF;
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;
477 vftab = kernel_data->table_vdw->data;
478 vftabscale = _mm_set1_ps(kernel_data->table_vdw->scale);
480 /* Avoid stupid compiler warnings */
481 jnrA = jnrB = jnrC = jnrD = 0;
490 for(iidx=0;iidx<4*DIM;iidx++)
495 /* Start outer loop over neighborlists */
496 for(iidx=0; iidx<nri; iidx++)
498 /* Load shift vector for this list */
499 i_shift_offset = DIM*shiftidx[iidx];
501 /* Load limits for loop over neighbors */
502 j_index_start = jindex[iidx];
503 j_index_end = jindex[iidx+1];
505 /* Get outer coordinate index */
507 i_coord_offset = DIM*inr;
509 /* Load i particle coords and add shift vector */
510 gmx_mm_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
512 fix0 = _mm_setzero_ps();
513 fiy0 = _mm_setzero_ps();
514 fiz0 = _mm_setzero_ps();
516 /* Load parameters for i particles */
517 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
519 /* Start inner kernel loop */
520 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
523 /* Get j neighbor index, and coordinate index */
528 j_coord_offsetA = DIM*jnrA;
529 j_coord_offsetB = DIM*jnrB;
530 j_coord_offsetC = DIM*jnrC;
531 j_coord_offsetD = DIM*jnrD;
533 /* load j atom coordinates */
534 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
535 x+j_coord_offsetC,x+j_coord_offsetD,
538 /* Calculate displacement vector */
539 dx00 = _mm_sub_ps(ix0,jx0);
540 dy00 = _mm_sub_ps(iy0,jy0);
541 dz00 = _mm_sub_ps(iz0,jz0);
543 /* Calculate squared distance and things based on it */
544 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
546 rinv00 = gmx_mm_invsqrt_ps(rsq00);
548 /* Load parameters for j particles */
549 vdwjidx0A = 2*vdwtype[jnrA+0];
550 vdwjidx0B = 2*vdwtype[jnrB+0];
551 vdwjidx0C = 2*vdwtype[jnrC+0];
552 vdwjidx0D = 2*vdwtype[jnrD+0];
554 /**************************
555 * CALCULATE INTERACTIONS *
556 **************************/
558 r00 = _mm_mul_ps(rsq00,rinv00);
560 /* Compute parameters for interactions between i and j atoms */
561 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
562 vdwparam+vdwioffset0+vdwjidx0B,
563 vdwparam+vdwioffset0+vdwjidx0C,
564 vdwparam+vdwioffset0+vdwjidx0D,
567 /* Calculate table index by multiplying r with table scale and truncate to integer */
568 rt = _mm_mul_ps(r00,vftabscale);
569 vfitab = _mm_cvttps_epi32(rt);
570 vfeps = _mm_sub_ps(rt,_mm_cvtepi32_ps(vfitab));
571 vfitab = _mm_slli_epi32(vfitab,3);
573 /* CUBIC SPLINE TABLE DISPERSION */
574 Y = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,0) );
575 F = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,1) );
576 G = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,2) );
577 H = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,3) );
578 _MM_TRANSPOSE4_PS(Y,F,G,H);
579 Heps = _mm_mul_ps(vfeps,H);
580 Fp = _mm_add_ps(F,_mm_mul_ps(vfeps,_mm_add_ps(G,Heps)));
581 FF = _mm_add_ps(Fp,_mm_mul_ps(vfeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
582 fvdw6 = _mm_mul_ps(c6_00,FF);
584 /* CUBIC SPLINE TABLE REPULSION */
585 vfitab = _mm_add_epi32(vfitab,ifour);
586 Y = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,0) );
587 F = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,1) );
588 G = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,2) );
589 H = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,3) );
590 _MM_TRANSPOSE4_PS(Y,F,G,H);
591 Heps = _mm_mul_ps(vfeps,H);
592 Fp = _mm_add_ps(F,_mm_mul_ps(vfeps,_mm_add_ps(G,Heps)));
593 FF = _mm_add_ps(Fp,_mm_mul_ps(vfeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
594 fvdw12 = _mm_mul_ps(c12_00,FF);
595 fvdw = _mm_xor_ps(signbit,_mm_mul_ps(_mm_add_ps(fvdw6,fvdw12),_mm_mul_ps(vftabscale,rinv00)));
599 /* Calculate temporary vectorial force */
600 tx = _mm_mul_ps(fscal,dx00);
601 ty = _mm_mul_ps(fscal,dy00);
602 tz = _mm_mul_ps(fscal,dz00);
604 /* Update vectorial force */
605 fix0 = _mm_add_ps(fix0,tx);
606 fiy0 = _mm_add_ps(fiy0,ty);
607 fiz0 = _mm_add_ps(fiz0,tz);
609 fjptrA = f+j_coord_offsetA;
610 fjptrB = f+j_coord_offsetB;
611 fjptrC = f+j_coord_offsetC;
612 fjptrD = f+j_coord_offsetD;
613 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
615 /* Inner loop uses 48 flops */
621 /* Get j neighbor index, and coordinate index */
622 jnrlistA = jjnr[jidx];
623 jnrlistB = jjnr[jidx+1];
624 jnrlistC = jjnr[jidx+2];
625 jnrlistD = jjnr[jidx+3];
626 /* Sign of each element will be negative for non-real atoms.
627 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
628 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
630 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
631 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
632 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
633 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
634 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
635 j_coord_offsetA = DIM*jnrA;
636 j_coord_offsetB = DIM*jnrB;
637 j_coord_offsetC = DIM*jnrC;
638 j_coord_offsetD = DIM*jnrD;
640 /* load j atom coordinates */
641 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
642 x+j_coord_offsetC,x+j_coord_offsetD,
645 /* Calculate displacement vector */
646 dx00 = _mm_sub_ps(ix0,jx0);
647 dy00 = _mm_sub_ps(iy0,jy0);
648 dz00 = _mm_sub_ps(iz0,jz0);
650 /* Calculate squared distance and things based on it */
651 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
653 rinv00 = gmx_mm_invsqrt_ps(rsq00);
655 /* Load parameters for j particles */
656 vdwjidx0A = 2*vdwtype[jnrA+0];
657 vdwjidx0B = 2*vdwtype[jnrB+0];
658 vdwjidx0C = 2*vdwtype[jnrC+0];
659 vdwjidx0D = 2*vdwtype[jnrD+0];
661 /**************************
662 * CALCULATE INTERACTIONS *
663 **************************/
665 r00 = _mm_mul_ps(rsq00,rinv00);
666 r00 = _mm_andnot_ps(dummy_mask,r00);
668 /* Compute parameters for interactions between i and j atoms */
669 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
670 vdwparam+vdwioffset0+vdwjidx0B,
671 vdwparam+vdwioffset0+vdwjidx0C,
672 vdwparam+vdwioffset0+vdwjidx0D,
675 /* Calculate table index by multiplying r with table scale and truncate to integer */
676 rt = _mm_mul_ps(r00,vftabscale);
677 vfitab = _mm_cvttps_epi32(rt);
678 vfeps = _mm_sub_ps(rt,_mm_cvtepi32_ps(vfitab));
679 vfitab = _mm_slli_epi32(vfitab,3);
681 /* CUBIC SPLINE TABLE DISPERSION */
682 Y = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,0) );
683 F = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,1) );
684 G = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,2) );
685 H = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,3) );
686 _MM_TRANSPOSE4_PS(Y,F,G,H);
687 Heps = _mm_mul_ps(vfeps,H);
688 Fp = _mm_add_ps(F,_mm_mul_ps(vfeps,_mm_add_ps(G,Heps)));
689 FF = _mm_add_ps(Fp,_mm_mul_ps(vfeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
690 fvdw6 = _mm_mul_ps(c6_00,FF);
692 /* CUBIC SPLINE TABLE REPULSION */
693 vfitab = _mm_add_epi32(vfitab,ifour);
694 Y = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,0) );
695 F = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,1) );
696 G = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,2) );
697 H = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,3) );
698 _MM_TRANSPOSE4_PS(Y,F,G,H);
699 Heps = _mm_mul_ps(vfeps,H);
700 Fp = _mm_add_ps(F,_mm_mul_ps(vfeps,_mm_add_ps(G,Heps)));
701 FF = _mm_add_ps(Fp,_mm_mul_ps(vfeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
702 fvdw12 = _mm_mul_ps(c12_00,FF);
703 fvdw = _mm_xor_ps(signbit,_mm_mul_ps(_mm_add_ps(fvdw6,fvdw12),_mm_mul_ps(vftabscale,rinv00)));
707 fscal = _mm_andnot_ps(dummy_mask,fscal);
709 /* Calculate temporary vectorial force */
710 tx = _mm_mul_ps(fscal,dx00);
711 ty = _mm_mul_ps(fscal,dy00);
712 tz = _mm_mul_ps(fscal,dz00);
714 /* Update vectorial force */
715 fix0 = _mm_add_ps(fix0,tx);
716 fiy0 = _mm_add_ps(fiy0,ty);
717 fiz0 = _mm_add_ps(fiz0,tz);
719 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
720 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
721 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
722 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
723 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
725 /* Inner loop uses 49 flops */
728 /* End of innermost loop */
730 gmx_mm_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
731 f+i_coord_offset,fshift+i_shift_offset);
733 /* Increment number of inner iterations */
734 inneriter += j_index_end - j_index_start;
736 /* Outer loop uses 6 flops */
739 /* Increment number of outer iterations */
742 /* Update outer/inner flops */
744 inc_nrnb(nrnb,eNR_NBKERNEL_VDW_F,outeriter*6 + inneriter*49);