2 * Note: this file was generated by the Gromacs sse4_1_single kernel generator.
4 * This source code is part of
8 * Copyright (c) 2001-2012, The GROMACS Development Team
10 * Gromacs is a library for molecular simulation and trajectory analysis,
11 * written by Erik Lindahl, David van der Spoel, Berk Hess, and others - for
12 * a full list of developers and information, check out http://www.gromacs.org
14 * This program is free software; you can redistribute it and/or modify it under
15 * the terms of the GNU Lesser General Public License as published by the Free
16 * Software Foundation; either version 2 of the License, or (at your option) any
19 * To help fund GROMACS development, we humbly ask that you cite
20 * the papers people have written on it - you can find them on the website.
28 #include "../nb_kernel.h"
29 #include "types/simple.h"
33 #include "gmx_math_x86_sse4_1_single.h"
34 #include "kernelutil_x86_sse4_1_single.h"
37 * Gromacs nonbonded kernel: nb_kernel_ElecGB_VdwLJ_GeomP1P1_VF_sse4_1_single
38 * Electrostatics interaction: GeneralizedBorn
39 * VdW interaction: LennardJones
40 * Geometry: Particle-Particle
41 * Calculate force/pot: PotentialAndForce
44 nb_kernel_ElecGB_VdwLJ_GeomP1P1_VF_sse4_1_single
45 (t_nblist * gmx_restrict nlist,
46 rvec * gmx_restrict xx,
47 rvec * gmx_restrict ff,
48 t_forcerec * gmx_restrict fr,
49 t_mdatoms * gmx_restrict mdatoms,
50 nb_kernel_data_t * gmx_restrict kernel_data,
51 t_nrnb * gmx_restrict nrnb)
53 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
54 * just 0 for non-waters.
55 * Suffixes A,B,C,D refer to j loop unrolling done with SSE, e.g. for the four different
56 * jnr indices corresponding to data put in the four positions in the SIMD register.
58 int i_shift_offset,i_coord_offset,outeriter,inneriter;
59 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
60 int jnrA,jnrB,jnrC,jnrD;
61 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
62 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
63 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
65 real *shiftvec,*fshift,*x,*f;
66 real *fjptrA,*fjptrB,*fjptrC,*fjptrD;
68 __m128 tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
70 __m128 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
71 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
72 __m128 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
73 __m128 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
74 __m128 velec,felec,velecsum,facel,crf,krf,krf2;
77 __m128 vgb,fgb,vgbsum,dvdasum,gbscale,gbtabscale,isaprod,gbqqfactor,gbinvepsdiff,gbeps,dvdatmp;
78 __m128 minushalf = _mm_set1_ps(-0.5);
79 real *invsqrta,*dvda,*gbtab;
81 __m128 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
84 __m128 one_sixth = _mm_set1_ps(1.0/6.0);
85 __m128 one_twelfth = _mm_set1_ps(1.0/12.0);
87 __m128i ifour = _mm_set1_epi32(4);
88 __m128 rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF;
90 __m128 dummy_mask,cutoff_mask;
91 __m128 signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
92 __m128 one = _mm_set1_ps(1.0);
93 __m128 two = _mm_set1_ps(2.0);
99 jindex = nlist->jindex;
101 shiftidx = nlist->shift;
103 shiftvec = fr->shift_vec[0];
104 fshift = fr->fshift[0];
105 facel = _mm_set1_ps(fr->epsfac);
106 charge = mdatoms->chargeA;
107 nvdwtype = fr->ntype;
109 vdwtype = mdatoms->typeA;
111 invsqrta = fr->invsqrta;
113 gbtabscale = _mm_set1_ps(fr->gbtab.scale);
114 gbtab = fr->gbtab.data;
115 gbinvepsdiff = _mm_set1_ps((1.0/fr->epsilon_r) - (1.0/fr->gb_epsilon_solvent));
117 /* Avoid stupid compiler warnings */
118 jnrA = jnrB = jnrC = jnrD = 0;
127 for(iidx=0;iidx<4*DIM;iidx++)
132 /* Start outer loop over neighborlists */
133 for(iidx=0; iidx<nri; iidx++)
135 /* Load shift vector for this list */
136 i_shift_offset = DIM*shiftidx[iidx];
138 /* Load limits for loop over neighbors */
139 j_index_start = jindex[iidx];
140 j_index_end = jindex[iidx+1];
142 /* Get outer coordinate index */
144 i_coord_offset = DIM*inr;
146 /* Load i particle coords and add shift vector */
147 gmx_mm_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
149 fix0 = _mm_setzero_ps();
150 fiy0 = _mm_setzero_ps();
151 fiz0 = _mm_setzero_ps();
153 /* Load parameters for i particles */
154 iq0 = _mm_mul_ps(facel,_mm_load1_ps(charge+inr+0));
155 isai0 = _mm_load1_ps(invsqrta+inr+0);
156 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
158 /* Reset potential sums */
159 velecsum = _mm_setzero_ps();
160 vgbsum = _mm_setzero_ps();
161 vvdwsum = _mm_setzero_ps();
162 dvdasum = _mm_setzero_ps();
164 /* Start inner kernel loop */
165 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
168 /* Get j neighbor index, and coordinate index */
173 j_coord_offsetA = DIM*jnrA;
174 j_coord_offsetB = DIM*jnrB;
175 j_coord_offsetC = DIM*jnrC;
176 j_coord_offsetD = DIM*jnrD;
178 /* load j atom coordinates */
179 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
180 x+j_coord_offsetC,x+j_coord_offsetD,
183 /* Calculate displacement vector */
184 dx00 = _mm_sub_ps(ix0,jx0);
185 dy00 = _mm_sub_ps(iy0,jy0);
186 dz00 = _mm_sub_ps(iz0,jz0);
188 /* Calculate squared distance and things based on it */
189 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
191 rinv00 = gmx_mm_invsqrt_ps(rsq00);
193 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
195 /* Load parameters for j particles */
196 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
197 charge+jnrC+0,charge+jnrD+0);
198 isaj0 = gmx_mm_load_4real_swizzle_ps(invsqrta+jnrA+0,invsqrta+jnrB+0,
199 invsqrta+jnrC+0,invsqrta+jnrD+0);
200 vdwjidx0A = 2*vdwtype[jnrA+0];
201 vdwjidx0B = 2*vdwtype[jnrB+0];
202 vdwjidx0C = 2*vdwtype[jnrC+0];
203 vdwjidx0D = 2*vdwtype[jnrD+0];
205 /**************************
206 * CALCULATE INTERACTIONS *
207 **************************/
209 r00 = _mm_mul_ps(rsq00,rinv00);
211 /* Compute parameters for interactions between i and j atoms */
212 qq00 = _mm_mul_ps(iq0,jq0);
213 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
214 vdwparam+vdwioffset0+vdwjidx0B,
215 vdwparam+vdwioffset0+vdwjidx0C,
216 vdwparam+vdwioffset0+vdwjidx0D,
219 /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
220 isaprod = _mm_mul_ps(isai0,isaj0);
221 gbqqfactor = _mm_xor_ps(signbit,_mm_mul_ps(qq00,_mm_mul_ps(isaprod,gbinvepsdiff)));
222 gbscale = _mm_mul_ps(isaprod,gbtabscale);
224 /* Calculate generalized born table index - this is a separate table from the normal one,
225 * but we use the same procedure by multiplying r with scale and truncating to integer.
227 rt = _mm_mul_ps(r00,gbscale);
228 gbitab = _mm_cvttps_epi32(rt);
229 gbeps = _mm_sub_ps(rt,_mm_round_ps(rt, _MM_FROUND_FLOOR));
230 gbitab = _mm_slli_epi32(gbitab,2);
231 Y = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,0) );
232 F = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,1) );
233 G = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,2) );
234 H = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,3) );
235 _MM_TRANSPOSE4_PS(Y,F,G,H);
236 Heps = _mm_mul_ps(gbeps,H);
237 Fp = _mm_add_ps(F,_mm_mul_ps(gbeps,_mm_add_ps(G,Heps)));
238 VV = _mm_add_ps(Y,_mm_mul_ps(gbeps,Fp));
239 vgb = _mm_mul_ps(gbqqfactor,VV);
241 FF = _mm_add_ps(Fp,_mm_mul_ps(gbeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
242 fgb = _mm_mul_ps(gbqqfactor,_mm_mul_ps(FF,gbscale));
243 dvdatmp = _mm_mul_ps(minushalf,_mm_add_ps(vgb,_mm_mul_ps(fgb,r00)));
244 dvdasum = _mm_add_ps(dvdasum,dvdatmp);
249 gmx_mm_increment_4real_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,_mm_mul_ps(dvdatmp,_mm_mul_ps(isaj0,isaj0)));
250 velec = _mm_mul_ps(qq00,rinv00);
251 felec = _mm_mul_ps(_mm_sub_ps(_mm_mul_ps(velec,rinv00),fgb),rinv00);
253 /* LENNARD-JONES DISPERSION/REPULSION */
255 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
256 vvdw6 = _mm_mul_ps(c6_00,rinvsix);
257 vvdw12 = _mm_mul_ps(c12_00,_mm_mul_ps(rinvsix,rinvsix));
258 vvdw = _mm_sub_ps( _mm_mul_ps(vvdw12,one_twelfth) , _mm_mul_ps(vvdw6,one_sixth) );
259 fvdw = _mm_mul_ps(_mm_sub_ps(vvdw12,vvdw6),rinvsq00);
261 /* Update potential sum for this i atom from the interaction with this j atom. */
262 velecsum = _mm_add_ps(velecsum,velec);
263 vgbsum = _mm_add_ps(vgbsum,vgb);
264 vvdwsum = _mm_add_ps(vvdwsum,vvdw);
266 fscal = _mm_add_ps(felec,fvdw);
268 /* Calculate temporary vectorial force */
269 tx = _mm_mul_ps(fscal,dx00);
270 ty = _mm_mul_ps(fscal,dy00);
271 tz = _mm_mul_ps(fscal,dz00);
273 /* Update vectorial force */
274 fix0 = _mm_add_ps(fix0,tx);
275 fiy0 = _mm_add_ps(fiy0,ty);
276 fiz0 = _mm_add_ps(fiz0,tz);
278 fjptrA = f+j_coord_offsetA;
279 fjptrB = f+j_coord_offsetB;
280 fjptrC = f+j_coord_offsetC;
281 fjptrD = f+j_coord_offsetD;
282 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
284 /* Inner loop uses 71 flops */
290 /* Get j neighbor index, and coordinate index */
291 jnrlistA = jjnr[jidx];
292 jnrlistB = jjnr[jidx+1];
293 jnrlistC = jjnr[jidx+2];
294 jnrlistD = jjnr[jidx+3];
295 /* Sign of each element will be negative for non-real atoms.
296 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
297 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
299 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
300 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
301 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
302 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
303 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
304 j_coord_offsetA = DIM*jnrA;
305 j_coord_offsetB = DIM*jnrB;
306 j_coord_offsetC = DIM*jnrC;
307 j_coord_offsetD = DIM*jnrD;
309 /* load j atom coordinates */
310 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
311 x+j_coord_offsetC,x+j_coord_offsetD,
314 /* Calculate displacement vector */
315 dx00 = _mm_sub_ps(ix0,jx0);
316 dy00 = _mm_sub_ps(iy0,jy0);
317 dz00 = _mm_sub_ps(iz0,jz0);
319 /* Calculate squared distance and things based on it */
320 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
322 rinv00 = gmx_mm_invsqrt_ps(rsq00);
324 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
326 /* Load parameters for j particles */
327 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
328 charge+jnrC+0,charge+jnrD+0);
329 isaj0 = gmx_mm_load_4real_swizzle_ps(invsqrta+jnrA+0,invsqrta+jnrB+0,
330 invsqrta+jnrC+0,invsqrta+jnrD+0);
331 vdwjidx0A = 2*vdwtype[jnrA+0];
332 vdwjidx0B = 2*vdwtype[jnrB+0];
333 vdwjidx0C = 2*vdwtype[jnrC+0];
334 vdwjidx0D = 2*vdwtype[jnrD+0];
336 /**************************
337 * CALCULATE INTERACTIONS *
338 **************************/
340 r00 = _mm_mul_ps(rsq00,rinv00);
341 r00 = _mm_andnot_ps(dummy_mask,r00);
343 /* Compute parameters for interactions between i and j atoms */
344 qq00 = _mm_mul_ps(iq0,jq0);
345 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
346 vdwparam+vdwioffset0+vdwjidx0B,
347 vdwparam+vdwioffset0+vdwjidx0C,
348 vdwparam+vdwioffset0+vdwjidx0D,
351 /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
352 isaprod = _mm_mul_ps(isai0,isaj0);
353 gbqqfactor = _mm_xor_ps(signbit,_mm_mul_ps(qq00,_mm_mul_ps(isaprod,gbinvepsdiff)));
354 gbscale = _mm_mul_ps(isaprod,gbtabscale);
356 /* Calculate generalized born table index - this is a separate table from the normal one,
357 * but we use the same procedure by multiplying r with scale and truncating to integer.
359 rt = _mm_mul_ps(r00,gbscale);
360 gbitab = _mm_cvttps_epi32(rt);
361 gbeps = _mm_sub_ps(rt,_mm_round_ps(rt, _MM_FROUND_FLOOR));
362 gbitab = _mm_slli_epi32(gbitab,2);
363 Y = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,0) );
364 F = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,1) );
365 G = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,2) );
366 H = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,3) );
367 _MM_TRANSPOSE4_PS(Y,F,G,H);
368 Heps = _mm_mul_ps(gbeps,H);
369 Fp = _mm_add_ps(F,_mm_mul_ps(gbeps,_mm_add_ps(G,Heps)));
370 VV = _mm_add_ps(Y,_mm_mul_ps(gbeps,Fp));
371 vgb = _mm_mul_ps(gbqqfactor,VV);
373 FF = _mm_add_ps(Fp,_mm_mul_ps(gbeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
374 fgb = _mm_mul_ps(gbqqfactor,_mm_mul_ps(FF,gbscale));
375 dvdatmp = _mm_mul_ps(minushalf,_mm_add_ps(vgb,_mm_mul_ps(fgb,r00)));
376 dvdatmp = _mm_andnot_ps(dummy_mask,dvdatmp);
377 dvdasum = _mm_add_ps(dvdasum,dvdatmp);
378 /* The pointers to scratch make sure that this code with compilers that take gmx_restrict seriously (e.g. icc 13) really can't screw things up. */
379 fjptrA = (jnrlistA>=0) ? dvda+jnrA : scratch;
380 fjptrB = (jnrlistB>=0) ? dvda+jnrB : scratch;
381 fjptrC = (jnrlistC>=0) ? dvda+jnrC : scratch;
382 fjptrD = (jnrlistD>=0) ? dvda+jnrD : scratch;
383 gmx_mm_increment_4real_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,_mm_mul_ps(dvdatmp,_mm_mul_ps(isaj0,isaj0)));
384 velec = _mm_mul_ps(qq00,rinv00);
385 felec = _mm_mul_ps(_mm_sub_ps(_mm_mul_ps(velec,rinv00),fgb),rinv00);
387 /* LENNARD-JONES DISPERSION/REPULSION */
389 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
390 vvdw6 = _mm_mul_ps(c6_00,rinvsix);
391 vvdw12 = _mm_mul_ps(c12_00,_mm_mul_ps(rinvsix,rinvsix));
392 vvdw = _mm_sub_ps( _mm_mul_ps(vvdw12,one_twelfth) , _mm_mul_ps(vvdw6,one_sixth) );
393 fvdw = _mm_mul_ps(_mm_sub_ps(vvdw12,vvdw6),rinvsq00);
395 /* Update potential sum for this i atom from the interaction with this j atom. */
396 velec = _mm_andnot_ps(dummy_mask,velec);
397 velecsum = _mm_add_ps(velecsum,velec);
398 vgb = _mm_andnot_ps(dummy_mask,vgb);
399 vgbsum = _mm_add_ps(vgbsum,vgb);
400 vvdw = _mm_andnot_ps(dummy_mask,vvdw);
401 vvdwsum = _mm_add_ps(vvdwsum,vvdw);
403 fscal = _mm_add_ps(felec,fvdw);
405 fscal = _mm_andnot_ps(dummy_mask,fscal);
407 /* Calculate temporary vectorial force */
408 tx = _mm_mul_ps(fscal,dx00);
409 ty = _mm_mul_ps(fscal,dy00);
410 tz = _mm_mul_ps(fscal,dz00);
412 /* Update vectorial force */
413 fix0 = _mm_add_ps(fix0,tx);
414 fiy0 = _mm_add_ps(fiy0,ty);
415 fiz0 = _mm_add_ps(fiz0,tz);
417 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
418 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
419 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
420 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
421 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
423 /* Inner loop uses 72 flops */
426 /* End of innermost loop */
428 gmx_mm_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
429 f+i_coord_offset,fshift+i_shift_offset);
432 /* Update potential energies */
433 gmx_mm_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
434 gmx_mm_update_1pot_ps(vgbsum,kernel_data->energygrp_polarization+ggid);
435 gmx_mm_update_1pot_ps(vvdwsum,kernel_data->energygrp_vdw+ggid);
436 dvdasum = _mm_mul_ps(dvdasum, _mm_mul_ps(isai0,isai0));
437 gmx_mm_update_1pot_ps(dvdasum,dvda+inr);
439 /* Increment number of inner iterations */
440 inneriter += j_index_end - j_index_start;
442 /* Outer loop uses 10 flops */
445 /* Increment number of outer iterations */
448 /* Update outer/inner flops */
450 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_VF,outeriter*10 + inneriter*72);
453 * Gromacs nonbonded kernel: nb_kernel_ElecGB_VdwLJ_GeomP1P1_F_sse4_1_single
454 * Electrostatics interaction: GeneralizedBorn
455 * VdW interaction: LennardJones
456 * Geometry: Particle-Particle
457 * Calculate force/pot: Force
460 nb_kernel_ElecGB_VdwLJ_GeomP1P1_F_sse4_1_single
461 (t_nblist * gmx_restrict nlist,
462 rvec * gmx_restrict xx,
463 rvec * gmx_restrict ff,
464 t_forcerec * gmx_restrict fr,
465 t_mdatoms * gmx_restrict mdatoms,
466 nb_kernel_data_t * gmx_restrict kernel_data,
467 t_nrnb * gmx_restrict nrnb)
469 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
470 * just 0 for non-waters.
471 * Suffixes A,B,C,D refer to j loop unrolling done with SSE, e.g. for the four different
472 * jnr indices corresponding to data put in the four positions in the SIMD register.
474 int i_shift_offset,i_coord_offset,outeriter,inneriter;
475 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
476 int jnrA,jnrB,jnrC,jnrD;
477 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
478 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
479 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
481 real *shiftvec,*fshift,*x,*f;
482 real *fjptrA,*fjptrB,*fjptrC,*fjptrD;
484 __m128 tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
486 __m128 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
487 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
488 __m128 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
489 __m128 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
490 __m128 velec,felec,velecsum,facel,crf,krf,krf2;
493 __m128 vgb,fgb,vgbsum,dvdasum,gbscale,gbtabscale,isaprod,gbqqfactor,gbinvepsdiff,gbeps,dvdatmp;
494 __m128 minushalf = _mm_set1_ps(-0.5);
495 real *invsqrta,*dvda,*gbtab;
497 __m128 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
500 __m128 one_sixth = _mm_set1_ps(1.0/6.0);
501 __m128 one_twelfth = _mm_set1_ps(1.0/12.0);
503 __m128i ifour = _mm_set1_epi32(4);
504 __m128 rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF;
506 __m128 dummy_mask,cutoff_mask;
507 __m128 signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
508 __m128 one = _mm_set1_ps(1.0);
509 __m128 two = _mm_set1_ps(2.0);
515 jindex = nlist->jindex;
517 shiftidx = nlist->shift;
519 shiftvec = fr->shift_vec[0];
520 fshift = fr->fshift[0];
521 facel = _mm_set1_ps(fr->epsfac);
522 charge = mdatoms->chargeA;
523 nvdwtype = fr->ntype;
525 vdwtype = mdatoms->typeA;
527 invsqrta = fr->invsqrta;
529 gbtabscale = _mm_set1_ps(fr->gbtab.scale);
530 gbtab = fr->gbtab.data;
531 gbinvepsdiff = _mm_set1_ps((1.0/fr->epsilon_r) - (1.0/fr->gb_epsilon_solvent));
533 /* Avoid stupid compiler warnings */
534 jnrA = jnrB = jnrC = jnrD = 0;
543 for(iidx=0;iidx<4*DIM;iidx++)
548 /* Start outer loop over neighborlists */
549 for(iidx=0; iidx<nri; iidx++)
551 /* Load shift vector for this list */
552 i_shift_offset = DIM*shiftidx[iidx];
554 /* Load limits for loop over neighbors */
555 j_index_start = jindex[iidx];
556 j_index_end = jindex[iidx+1];
558 /* Get outer coordinate index */
560 i_coord_offset = DIM*inr;
562 /* Load i particle coords and add shift vector */
563 gmx_mm_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
565 fix0 = _mm_setzero_ps();
566 fiy0 = _mm_setzero_ps();
567 fiz0 = _mm_setzero_ps();
569 /* Load parameters for i particles */
570 iq0 = _mm_mul_ps(facel,_mm_load1_ps(charge+inr+0));
571 isai0 = _mm_load1_ps(invsqrta+inr+0);
572 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
574 dvdasum = _mm_setzero_ps();
576 /* Start inner kernel loop */
577 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
580 /* Get j neighbor index, and coordinate index */
585 j_coord_offsetA = DIM*jnrA;
586 j_coord_offsetB = DIM*jnrB;
587 j_coord_offsetC = DIM*jnrC;
588 j_coord_offsetD = DIM*jnrD;
590 /* load j atom coordinates */
591 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
592 x+j_coord_offsetC,x+j_coord_offsetD,
595 /* Calculate displacement vector */
596 dx00 = _mm_sub_ps(ix0,jx0);
597 dy00 = _mm_sub_ps(iy0,jy0);
598 dz00 = _mm_sub_ps(iz0,jz0);
600 /* Calculate squared distance and things based on it */
601 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
603 rinv00 = gmx_mm_invsqrt_ps(rsq00);
605 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
607 /* Load parameters for j particles */
608 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
609 charge+jnrC+0,charge+jnrD+0);
610 isaj0 = gmx_mm_load_4real_swizzle_ps(invsqrta+jnrA+0,invsqrta+jnrB+0,
611 invsqrta+jnrC+0,invsqrta+jnrD+0);
612 vdwjidx0A = 2*vdwtype[jnrA+0];
613 vdwjidx0B = 2*vdwtype[jnrB+0];
614 vdwjidx0C = 2*vdwtype[jnrC+0];
615 vdwjidx0D = 2*vdwtype[jnrD+0];
617 /**************************
618 * CALCULATE INTERACTIONS *
619 **************************/
621 r00 = _mm_mul_ps(rsq00,rinv00);
623 /* Compute parameters for interactions between i and j atoms */
624 qq00 = _mm_mul_ps(iq0,jq0);
625 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
626 vdwparam+vdwioffset0+vdwjidx0B,
627 vdwparam+vdwioffset0+vdwjidx0C,
628 vdwparam+vdwioffset0+vdwjidx0D,
631 /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
632 isaprod = _mm_mul_ps(isai0,isaj0);
633 gbqqfactor = _mm_xor_ps(signbit,_mm_mul_ps(qq00,_mm_mul_ps(isaprod,gbinvepsdiff)));
634 gbscale = _mm_mul_ps(isaprod,gbtabscale);
636 /* Calculate generalized born table index - this is a separate table from the normal one,
637 * but we use the same procedure by multiplying r with scale and truncating to integer.
639 rt = _mm_mul_ps(r00,gbscale);
640 gbitab = _mm_cvttps_epi32(rt);
641 gbeps = _mm_sub_ps(rt,_mm_round_ps(rt, _MM_FROUND_FLOOR));
642 gbitab = _mm_slli_epi32(gbitab,2);
643 Y = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,0) );
644 F = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,1) );
645 G = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,2) );
646 H = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,3) );
647 _MM_TRANSPOSE4_PS(Y,F,G,H);
648 Heps = _mm_mul_ps(gbeps,H);
649 Fp = _mm_add_ps(F,_mm_mul_ps(gbeps,_mm_add_ps(G,Heps)));
650 VV = _mm_add_ps(Y,_mm_mul_ps(gbeps,Fp));
651 vgb = _mm_mul_ps(gbqqfactor,VV);
653 FF = _mm_add_ps(Fp,_mm_mul_ps(gbeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
654 fgb = _mm_mul_ps(gbqqfactor,_mm_mul_ps(FF,gbscale));
655 dvdatmp = _mm_mul_ps(minushalf,_mm_add_ps(vgb,_mm_mul_ps(fgb,r00)));
656 dvdasum = _mm_add_ps(dvdasum,dvdatmp);
661 gmx_mm_increment_4real_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,_mm_mul_ps(dvdatmp,_mm_mul_ps(isaj0,isaj0)));
662 velec = _mm_mul_ps(qq00,rinv00);
663 felec = _mm_mul_ps(_mm_sub_ps(_mm_mul_ps(velec,rinv00),fgb),rinv00);
665 /* LENNARD-JONES DISPERSION/REPULSION */
667 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
668 fvdw = _mm_mul_ps(_mm_sub_ps(_mm_mul_ps(c12_00,rinvsix),c6_00),_mm_mul_ps(rinvsix,rinvsq00));
670 fscal = _mm_add_ps(felec,fvdw);
672 /* Calculate temporary vectorial force */
673 tx = _mm_mul_ps(fscal,dx00);
674 ty = _mm_mul_ps(fscal,dy00);
675 tz = _mm_mul_ps(fscal,dz00);
677 /* Update vectorial force */
678 fix0 = _mm_add_ps(fix0,tx);
679 fiy0 = _mm_add_ps(fiy0,ty);
680 fiz0 = _mm_add_ps(fiz0,tz);
682 fjptrA = f+j_coord_offsetA;
683 fjptrB = f+j_coord_offsetB;
684 fjptrC = f+j_coord_offsetC;
685 fjptrD = f+j_coord_offsetD;
686 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
688 /* Inner loop uses 64 flops */
694 /* Get j neighbor index, and coordinate index */
695 jnrlistA = jjnr[jidx];
696 jnrlistB = jjnr[jidx+1];
697 jnrlistC = jjnr[jidx+2];
698 jnrlistD = jjnr[jidx+3];
699 /* Sign of each element will be negative for non-real atoms.
700 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
701 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
703 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
704 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
705 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
706 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
707 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
708 j_coord_offsetA = DIM*jnrA;
709 j_coord_offsetB = DIM*jnrB;
710 j_coord_offsetC = DIM*jnrC;
711 j_coord_offsetD = DIM*jnrD;
713 /* load j atom coordinates */
714 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
715 x+j_coord_offsetC,x+j_coord_offsetD,
718 /* Calculate displacement vector */
719 dx00 = _mm_sub_ps(ix0,jx0);
720 dy00 = _mm_sub_ps(iy0,jy0);
721 dz00 = _mm_sub_ps(iz0,jz0);
723 /* Calculate squared distance and things based on it */
724 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
726 rinv00 = gmx_mm_invsqrt_ps(rsq00);
728 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
730 /* Load parameters for j particles */
731 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
732 charge+jnrC+0,charge+jnrD+0);
733 isaj0 = gmx_mm_load_4real_swizzle_ps(invsqrta+jnrA+0,invsqrta+jnrB+0,
734 invsqrta+jnrC+0,invsqrta+jnrD+0);
735 vdwjidx0A = 2*vdwtype[jnrA+0];
736 vdwjidx0B = 2*vdwtype[jnrB+0];
737 vdwjidx0C = 2*vdwtype[jnrC+0];
738 vdwjidx0D = 2*vdwtype[jnrD+0];
740 /**************************
741 * CALCULATE INTERACTIONS *
742 **************************/
744 r00 = _mm_mul_ps(rsq00,rinv00);
745 r00 = _mm_andnot_ps(dummy_mask,r00);
747 /* Compute parameters for interactions between i and j atoms */
748 qq00 = _mm_mul_ps(iq0,jq0);
749 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
750 vdwparam+vdwioffset0+vdwjidx0B,
751 vdwparam+vdwioffset0+vdwjidx0C,
752 vdwparam+vdwioffset0+vdwjidx0D,
755 /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
756 isaprod = _mm_mul_ps(isai0,isaj0);
757 gbqqfactor = _mm_xor_ps(signbit,_mm_mul_ps(qq00,_mm_mul_ps(isaprod,gbinvepsdiff)));
758 gbscale = _mm_mul_ps(isaprod,gbtabscale);
760 /* Calculate generalized born table index - this is a separate table from the normal one,
761 * but we use the same procedure by multiplying r with scale and truncating to integer.
763 rt = _mm_mul_ps(r00,gbscale);
764 gbitab = _mm_cvttps_epi32(rt);
765 gbeps = _mm_sub_ps(rt,_mm_round_ps(rt, _MM_FROUND_FLOOR));
766 gbitab = _mm_slli_epi32(gbitab,2);
767 Y = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,0) );
768 F = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,1) );
769 G = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,2) );
770 H = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,3) );
771 _MM_TRANSPOSE4_PS(Y,F,G,H);
772 Heps = _mm_mul_ps(gbeps,H);
773 Fp = _mm_add_ps(F,_mm_mul_ps(gbeps,_mm_add_ps(G,Heps)));
774 VV = _mm_add_ps(Y,_mm_mul_ps(gbeps,Fp));
775 vgb = _mm_mul_ps(gbqqfactor,VV);
777 FF = _mm_add_ps(Fp,_mm_mul_ps(gbeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
778 fgb = _mm_mul_ps(gbqqfactor,_mm_mul_ps(FF,gbscale));
779 dvdatmp = _mm_mul_ps(minushalf,_mm_add_ps(vgb,_mm_mul_ps(fgb,r00)));
780 dvdatmp = _mm_andnot_ps(dummy_mask,dvdatmp);
781 dvdasum = _mm_add_ps(dvdasum,dvdatmp);
782 /* The pointers to scratch make sure that this code with compilers that take gmx_restrict seriously (e.g. icc 13) really can't screw things up. */
783 fjptrA = (jnrlistA>=0) ? dvda+jnrA : scratch;
784 fjptrB = (jnrlistB>=0) ? dvda+jnrB : scratch;
785 fjptrC = (jnrlistC>=0) ? dvda+jnrC : scratch;
786 fjptrD = (jnrlistD>=0) ? dvda+jnrD : scratch;
787 gmx_mm_increment_4real_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,_mm_mul_ps(dvdatmp,_mm_mul_ps(isaj0,isaj0)));
788 velec = _mm_mul_ps(qq00,rinv00);
789 felec = _mm_mul_ps(_mm_sub_ps(_mm_mul_ps(velec,rinv00),fgb),rinv00);
791 /* LENNARD-JONES DISPERSION/REPULSION */
793 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
794 fvdw = _mm_mul_ps(_mm_sub_ps(_mm_mul_ps(c12_00,rinvsix),c6_00),_mm_mul_ps(rinvsix,rinvsq00));
796 fscal = _mm_add_ps(felec,fvdw);
798 fscal = _mm_andnot_ps(dummy_mask,fscal);
800 /* Calculate temporary vectorial force */
801 tx = _mm_mul_ps(fscal,dx00);
802 ty = _mm_mul_ps(fscal,dy00);
803 tz = _mm_mul_ps(fscal,dz00);
805 /* Update vectorial force */
806 fix0 = _mm_add_ps(fix0,tx);
807 fiy0 = _mm_add_ps(fiy0,ty);
808 fiz0 = _mm_add_ps(fiz0,tz);
810 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
811 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
812 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
813 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
814 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
816 /* Inner loop uses 65 flops */
819 /* End of innermost loop */
821 gmx_mm_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
822 f+i_coord_offset,fshift+i_shift_offset);
824 dvdasum = _mm_mul_ps(dvdasum, _mm_mul_ps(isai0,isai0));
825 gmx_mm_update_1pot_ps(dvdasum,dvda+inr);
827 /* Increment number of inner iterations */
828 inneriter += j_index_end - j_index_start;
830 /* Outer loop uses 7 flops */
833 /* Increment number of outer iterations */
836 /* Update outer/inner flops */
838 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_F,outeriter*7 + inneriter*65);