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_ElecEwSw_VdwNone_GeomP1P1_VF_sse4_1_single
38 * Electrostatics interaction: Ewald
39 * VdW interaction: None
40 * Geometry: Particle-Particle
41 * Calculate force/pot: PotentialAndForce
44 nb_kernel_ElecEwSw_VdwNone_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 ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
79 __m128 rswitch,swV3,swV4,swV5,swF2,swF3,swF4,d,d2,sw,dsw;
80 real rswitch_scalar,d_scalar;
81 __m128 dummy_mask,cutoff_mask;
82 __m128 signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
83 __m128 one = _mm_set1_ps(1.0);
84 __m128 two = _mm_set1_ps(2.0);
90 jindex = nlist->jindex;
92 shiftidx = nlist->shift;
94 shiftvec = fr->shift_vec[0];
95 fshift = fr->fshift[0];
96 facel = _mm_set1_ps(fr->epsfac);
97 charge = mdatoms->chargeA;
99 sh_ewald = _mm_set1_ps(fr->ic->sh_ewald);
100 ewtab = fr->ic->tabq_coul_FDV0;
101 ewtabscale = _mm_set1_ps(fr->ic->tabq_scale);
102 ewtabhalfspace = _mm_set1_ps(0.5/fr->ic->tabq_scale);
104 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
105 rcutoff_scalar = fr->rcoulomb;
106 rcutoff = _mm_set1_ps(rcutoff_scalar);
107 rcutoff2 = _mm_mul_ps(rcutoff,rcutoff);
109 rswitch_scalar = fr->rcoulomb_switch;
110 rswitch = _mm_set1_ps(rswitch_scalar);
111 /* Setup switch parameters */
112 d_scalar = rcutoff_scalar-rswitch_scalar;
113 d = _mm_set1_ps(d_scalar);
114 swV3 = _mm_set1_ps(-10.0/(d_scalar*d_scalar*d_scalar));
115 swV4 = _mm_set1_ps( 15.0/(d_scalar*d_scalar*d_scalar*d_scalar));
116 swV5 = _mm_set1_ps( -6.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
117 swF2 = _mm_set1_ps(-30.0/(d_scalar*d_scalar*d_scalar));
118 swF3 = _mm_set1_ps( 60.0/(d_scalar*d_scalar*d_scalar*d_scalar));
119 swF4 = _mm_set1_ps(-30.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
121 /* Avoid stupid compiler warnings */
122 jnrA = jnrB = jnrC = jnrD = 0;
131 for(iidx=0;iidx<4*DIM;iidx++)
136 /* Start outer loop over neighborlists */
137 for(iidx=0; iidx<nri; iidx++)
139 /* Load shift vector for this list */
140 i_shift_offset = DIM*shiftidx[iidx];
142 /* Load limits for loop over neighbors */
143 j_index_start = jindex[iidx];
144 j_index_end = jindex[iidx+1];
146 /* Get outer coordinate index */
148 i_coord_offset = DIM*inr;
150 /* Load i particle coords and add shift vector */
151 gmx_mm_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
153 fix0 = _mm_setzero_ps();
154 fiy0 = _mm_setzero_ps();
155 fiz0 = _mm_setzero_ps();
157 /* Load parameters for i particles */
158 iq0 = _mm_mul_ps(facel,_mm_load1_ps(charge+inr+0));
160 /* Reset potential sums */
161 velecsum = _mm_setzero_ps();
163 /* Start inner kernel loop */
164 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
167 /* Get j neighbor index, and coordinate index */
172 j_coord_offsetA = DIM*jnrA;
173 j_coord_offsetB = DIM*jnrB;
174 j_coord_offsetC = DIM*jnrC;
175 j_coord_offsetD = DIM*jnrD;
177 /* load j atom coordinates */
178 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
179 x+j_coord_offsetC,x+j_coord_offsetD,
182 /* Calculate displacement vector */
183 dx00 = _mm_sub_ps(ix0,jx0);
184 dy00 = _mm_sub_ps(iy0,jy0);
185 dz00 = _mm_sub_ps(iz0,jz0);
187 /* Calculate squared distance and things based on it */
188 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
190 rinv00 = gmx_mm_invsqrt_ps(rsq00);
192 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
194 /* Load parameters for j particles */
195 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
196 charge+jnrC+0,charge+jnrD+0);
198 /**************************
199 * CALCULATE INTERACTIONS *
200 **************************/
202 if (gmx_mm_any_lt(rsq00,rcutoff2))
205 r00 = _mm_mul_ps(rsq00,rinv00);
207 /* Compute parameters for interactions between i and j atoms */
208 qq00 = _mm_mul_ps(iq0,jq0);
210 /* EWALD ELECTROSTATICS */
212 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
213 ewrt = _mm_mul_ps(r00,ewtabscale);
214 ewitab = _mm_cvttps_epi32(ewrt);
215 eweps = _mm_sub_ps(ewrt,_mm_round_ps(ewrt, _MM_FROUND_FLOOR));
216 ewitab = _mm_slli_epi32(ewitab,2);
217 ewtabF = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,0) );
218 ewtabD = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,1) );
219 ewtabV = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,2) );
220 ewtabFn = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,3) );
221 _MM_TRANSPOSE4_PS(ewtabF,ewtabD,ewtabV,ewtabFn);
222 felec = _mm_add_ps(ewtabF,_mm_mul_ps(eweps,ewtabD));
223 velec = _mm_sub_ps(ewtabV,_mm_mul_ps(_mm_mul_ps(ewtabhalfspace,eweps),_mm_add_ps(ewtabF,felec)));
224 velec = _mm_mul_ps(qq00,_mm_sub_ps(rinv00,velec));
225 felec = _mm_mul_ps(_mm_mul_ps(qq00,rinv00),_mm_sub_ps(rinvsq00,felec));
227 d = _mm_sub_ps(r00,rswitch);
228 d = _mm_max_ps(d,_mm_setzero_ps());
229 d2 = _mm_mul_ps(d,d);
230 sw = _mm_add_ps(one,_mm_mul_ps(d2,_mm_mul_ps(d,_mm_add_ps(swV3,_mm_mul_ps(d,_mm_add_ps(swV4,_mm_mul_ps(d,swV5)))))));
232 dsw = _mm_mul_ps(d2,_mm_add_ps(swF2,_mm_mul_ps(d,_mm_add_ps(swF3,_mm_mul_ps(d,swF4)))));
234 /* Evaluate switch function */
235 /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
236 felec = _mm_sub_ps( _mm_mul_ps(felec,sw) , _mm_mul_ps(rinv00,_mm_mul_ps(velec,dsw)) );
237 velec = _mm_mul_ps(velec,sw);
238 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
240 /* Update potential sum for this i atom from the interaction with this j atom. */
241 velec = _mm_and_ps(velec,cutoff_mask);
242 velecsum = _mm_add_ps(velecsum,velec);
246 fscal = _mm_and_ps(fscal,cutoff_mask);
248 /* Calculate temporary vectorial force */
249 tx = _mm_mul_ps(fscal,dx00);
250 ty = _mm_mul_ps(fscal,dy00);
251 tz = _mm_mul_ps(fscal,dz00);
253 /* Update vectorial force */
254 fix0 = _mm_add_ps(fix0,tx);
255 fiy0 = _mm_add_ps(fiy0,ty);
256 fiz0 = _mm_add_ps(fiz0,tz);
258 fjptrA = f+j_coord_offsetA;
259 fjptrB = f+j_coord_offsetB;
260 fjptrC = f+j_coord_offsetC;
261 fjptrD = f+j_coord_offsetD;
262 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
266 /* Inner loop uses 65 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 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
308 /* Load parameters for j particles */
309 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
310 charge+jnrC+0,charge+jnrD+0);
312 /**************************
313 * CALCULATE INTERACTIONS *
314 **************************/
316 if (gmx_mm_any_lt(rsq00,rcutoff2))
319 r00 = _mm_mul_ps(rsq00,rinv00);
320 r00 = _mm_andnot_ps(dummy_mask,r00);
322 /* Compute parameters for interactions between i and j atoms */
323 qq00 = _mm_mul_ps(iq0,jq0);
325 /* EWALD ELECTROSTATICS */
327 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
328 ewrt = _mm_mul_ps(r00,ewtabscale);
329 ewitab = _mm_cvttps_epi32(ewrt);
330 eweps = _mm_sub_ps(ewrt,_mm_round_ps(ewrt, _MM_FROUND_FLOOR));
331 ewitab = _mm_slli_epi32(ewitab,2);
332 ewtabF = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,0) );
333 ewtabD = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,1) );
334 ewtabV = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,2) );
335 ewtabFn = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,3) );
336 _MM_TRANSPOSE4_PS(ewtabF,ewtabD,ewtabV,ewtabFn);
337 felec = _mm_add_ps(ewtabF,_mm_mul_ps(eweps,ewtabD));
338 velec = _mm_sub_ps(ewtabV,_mm_mul_ps(_mm_mul_ps(ewtabhalfspace,eweps),_mm_add_ps(ewtabF,felec)));
339 velec = _mm_mul_ps(qq00,_mm_sub_ps(rinv00,velec));
340 felec = _mm_mul_ps(_mm_mul_ps(qq00,rinv00),_mm_sub_ps(rinvsq00,felec));
342 d = _mm_sub_ps(r00,rswitch);
343 d = _mm_max_ps(d,_mm_setzero_ps());
344 d2 = _mm_mul_ps(d,d);
345 sw = _mm_add_ps(one,_mm_mul_ps(d2,_mm_mul_ps(d,_mm_add_ps(swV3,_mm_mul_ps(d,_mm_add_ps(swV4,_mm_mul_ps(d,swV5)))))));
347 dsw = _mm_mul_ps(d2,_mm_add_ps(swF2,_mm_mul_ps(d,_mm_add_ps(swF3,_mm_mul_ps(d,swF4)))));
349 /* Evaluate switch function */
350 /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
351 felec = _mm_sub_ps( _mm_mul_ps(felec,sw) , _mm_mul_ps(rinv00,_mm_mul_ps(velec,dsw)) );
352 velec = _mm_mul_ps(velec,sw);
353 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
355 /* Update potential sum for this i atom from the interaction with this j atom. */
356 velec = _mm_and_ps(velec,cutoff_mask);
357 velec = _mm_andnot_ps(dummy_mask,velec);
358 velecsum = _mm_add_ps(velecsum,velec);
362 fscal = _mm_and_ps(fscal,cutoff_mask);
364 fscal = _mm_andnot_ps(dummy_mask,fscal);
366 /* Calculate temporary vectorial force */
367 tx = _mm_mul_ps(fscal,dx00);
368 ty = _mm_mul_ps(fscal,dy00);
369 tz = _mm_mul_ps(fscal,dz00);
371 /* Update vectorial force */
372 fix0 = _mm_add_ps(fix0,tx);
373 fiy0 = _mm_add_ps(fiy0,ty);
374 fiz0 = _mm_add_ps(fiz0,tz);
376 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
377 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
378 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
379 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
380 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
384 /* Inner loop uses 66 flops */
387 /* End of innermost loop */
389 gmx_mm_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
390 f+i_coord_offset,fshift+i_shift_offset);
393 /* Update potential energies */
394 gmx_mm_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
396 /* Increment number of inner iterations */
397 inneriter += j_index_end - j_index_start;
399 /* Outer loop uses 8 flops */
402 /* Increment number of outer iterations */
405 /* Update outer/inner flops */
407 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VF,outeriter*8 + inneriter*66);
410 * Gromacs nonbonded kernel: nb_kernel_ElecEwSw_VdwNone_GeomP1P1_F_sse4_1_single
411 * Electrostatics interaction: Ewald
412 * VdW interaction: None
413 * Geometry: Particle-Particle
414 * Calculate force/pot: Force
417 nb_kernel_ElecEwSw_VdwNone_GeomP1P1_F_sse4_1_single
418 (t_nblist * gmx_restrict nlist,
419 rvec * gmx_restrict xx,
420 rvec * gmx_restrict ff,
421 t_forcerec * gmx_restrict fr,
422 t_mdatoms * gmx_restrict mdatoms,
423 nb_kernel_data_t * gmx_restrict kernel_data,
424 t_nrnb * gmx_restrict nrnb)
426 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
427 * just 0 for non-waters.
428 * Suffixes A,B,C,D refer to j loop unrolling done with SSE, e.g. for the four different
429 * jnr indices corresponding to data put in the four positions in the SIMD register.
431 int i_shift_offset,i_coord_offset,outeriter,inneriter;
432 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
433 int jnrA,jnrB,jnrC,jnrD;
434 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
435 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
436 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
438 real *shiftvec,*fshift,*x,*f;
439 real *fjptrA,*fjptrB,*fjptrC,*fjptrD;
441 __m128 tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
443 __m128 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
444 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
445 __m128 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
446 __m128 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
447 __m128 velec,felec,velecsum,facel,crf,krf,krf2;
450 __m128 ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
452 __m128 rswitch,swV3,swV4,swV5,swF2,swF3,swF4,d,d2,sw,dsw;
453 real rswitch_scalar,d_scalar;
454 __m128 dummy_mask,cutoff_mask;
455 __m128 signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
456 __m128 one = _mm_set1_ps(1.0);
457 __m128 two = _mm_set1_ps(2.0);
463 jindex = nlist->jindex;
465 shiftidx = nlist->shift;
467 shiftvec = fr->shift_vec[0];
468 fshift = fr->fshift[0];
469 facel = _mm_set1_ps(fr->epsfac);
470 charge = mdatoms->chargeA;
472 sh_ewald = _mm_set1_ps(fr->ic->sh_ewald);
473 ewtab = fr->ic->tabq_coul_FDV0;
474 ewtabscale = _mm_set1_ps(fr->ic->tabq_scale);
475 ewtabhalfspace = _mm_set1_ps(0.5/fr->ic->tabq_scale);
477 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
478 rcutoff_scalar = fr->rcoulomb;
479 rcutoff = _mm_set1_ps(rcutoff_scalar);
480 rcutoff2 = _mm_mul_ps(rcutoff,rcutoff);
482 rswitch_scalar = fr->rcoulomb_switch;
483 rswitch = _mm_set1_ps(rswitch_scalar);
484 /* Setup switch parameters */
485 d_scalar = rcutoff_scalar-rswitch_scalar;
486 d = _mm_set1_ps(d_scalar);
487 swV3 = _mm_set1_ps(-10.0/(d_scalar*d_scalar*d_scalar));
488 swV4 = _mm_set1_ps( 15.0/(d_scalar*d_scalar*d_scalar*d_scalar));
489 swV5 = _mm_set1_ps( -6.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
490 swF2 = _mm_set1_ps(-30.0/(d_scalar*d_scalar*d_scalar));
491 swF3 = _mm_set1_ps( 60.0/(d_scalar*d_scalar*d_scalar*d_scalar));
492 swF4 = _mm_set1_ps(-30.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
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 iq0 = _mm_mul_ps(facel,_mm_load1_ps(charge+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 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
566 charge+jnrC+0,charge+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 qq00 = _mm_mul_ps(iq0,jq0);
580 /* EWALD ELECTROSTATICS */
582 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
583 ewrt = _mm_mul_ps(r00,ewtabscale);
584 ewitab = _mm_cvttps_epi32(ewrt);
585 eweps = _mm_sub_ps(ewrt,_mm_round_ps(ewrt, _MM_FROUND_FLOOR));
586 ewitab = _mm_slli_epi32(ewitab,2);
587 ewtabF = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,0) );
588 ewtabD = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,1) );
589 ewtabV = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,2) );
590 ewtabFn = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,3) );
591 _MM_TRANSPOSE4_PS(ewtabF,ewtabD,ewtabV,ewtabFn);
592 felec = _mm_add_ps(ewtabF,_mm_mul_ps(eweps,ewtabD));
593 velec = _mm_sub_ps(ewtabV,_mm_mul_ps(_mm_mul_ps(ewtabhalfspace,eweps),_mm_add_ps(ewtabF,felec)));
594 velec = _mm_mul_ps(qq00,_mm_sub_ps(rinv00,velec));
595 felec = _mm_mul_ps(_mm_mul_ps(qq00,rinv00),_mm_sub_ps(rinvsq00,felec));
597 d = _mm_sub_ps(r00,rswitch);
598 d = _mm_max_ps(d,_mm_setzero_ps());
599 d2 = _mm_mul_ps(d,d);
600 sw = _mm_add_ps(one,_mm_mul_ps(d2,_mm_mul_ps(d,_mm_add_ps(swV3,_mm_mul_ps(d,_mm_add_ps(swV4,_mm_mul_ps(d,swV5)))))));
602 dsw = _mm_mul_ps(d2,_mm_add_ps(swF2,_mm_mul_ps(d,_mm_add_ps(swF3,_mm_mul_ps(d,swF4)))));
604 /* Evaluate switch function */
605 /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
606 felec = _mm_sub_ps( _mm_mul_ps(felec,sw) , _mm_mul_ps(rinv00,_mm_mul_ps(velec,dsw)) );
607 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
611 fscal = _mm_and_ps(fscal,cutoff_mask);
613 /* Calculate temporary vectorial force */
614 tx = _mm_mul_ps(fscal,dx00);
615 ty = _mm_mul_ps(fscal,dy00);
616 tz = _mm_mul_ps(fscal,dz00);
618 /* Update vectorial force */
619 fix0 = _mm_add_ps(fix0,tx);
620 fiy0 = _mm_add_ps(fiy0,ty);
621 fiz0 = _mm_add_ps(fiz0,tz);
623 fjptrA = f+j_coord_offsetA;
624 fjptrB = f+j_coord_offsetB;
625 fjptrC = f+j_coord_offsetC;
626 fjptrD = f+j_coord_offsetD;
627 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
631 /* Inner loop uses 62 flops */
637 /* Get j neighbor index, and coordinate index */
638 jnrlistA = jjnr[jidx];
639 jnrlistB = jjnr[jidx+1];
640 jnrlistC = jjnr[jidx+2];
641 jnrlistD = jjnr[jidx+3];
642 /* Sign of each element will be negative for non-real atoms.
643 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
644 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
646 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
647 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
648 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
649 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
650 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
651 j_coord_offsetA = DIM*jnrA;
652 j_coord_offsetB = DIM*jnrB;
653 j_coord_offsetC = DIM*jnrC;
654 j_coord_offsetD = DIM*jnrD;
656 /* load j atom coordinates */
657 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
658 x+j_coord_offsetC,x+j_coord_offsetD,
661 /* Calculate displacement vector */
662 dx00 = _mm_sub_ps(ix0,jx0);
663 dy00 = _mm_sub_ps(iy0,jy0);
664 dz00 = _mm_sub_ps(iz0,jz0);
666 /* Calculate squared distance and things based on it */
667 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
669 rinv00 = gmx_mm_invsqrt_ps(rsq00);
671 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
673 /* Load parameters for j particles */
674 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
675 charge+jnrC+0,charge+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 qq00 = _mm_mul_ps(iq0,jq0);
690 /* EWALD ELECTROSTATICS */
692 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
693 ewrt = _mm_mul_ps(r00,ewtabscale);
694 ewitab = _mm_cvttps_epi32(ewrt);
695 eweps = _mm_sub_ps(ewrt,_mm_round_ps(ewrt, _MM_FROUND_FLOOR));
696 ewitab = _mm_slli_epi32(ewitab,2);
697 ewtabF = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,0) );
698 ewtabD = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,1) );
699 ewtabV = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,2) );
700 ewtabFn = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,3) );
701 _MM_TRANSPOSE4_PS(ewtabF,ewtabD,ewtabV,ewtabFn);
702 felec = _mm_add_ps(ewtabF,_mm_mul_ps(eweps,ewtabD));
703 velec = _mm_sub_ps(ewtabV,_mm_mul_ps(_mm_mul_ps(ewtabhalfspace,eweps),_mm_add_ps(ewtabF,felec)));
704 velec = _mm_mul_ps(qq00,_mm_sub_ps(rinv00,velec));
705 felec = _mm_mul_ps(_mm_mul_ps(qq00,rinv00),_mm_sub_ps(rinvsq00,felec));
707 d = _mm_sub_ps(r00,rswitch);
708 d = _mm_max_ps(d,_mm_setzero_ps());
709 d2 = _mm_mul_ps(d,d);
710 sw = _mm_add_ps(one,_mm_mul_ps(d2,_mm_mul_ps(d,_mm_add_ps(swV3,_mm_mul_ps(d,_mm_add_ps(swV4,_mm_mul_ps(d,swV5)))))));
712 dsw = _mm_mul_ps(d2,_mm_add_ps(swF2,_mm_mul_ps(d,_mm_add_ps(swF3,_mm_mul_ps(d,swF4)))));
714 /* Evaluate switch function */
715 /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
716 felec = _mm_sub_ps( _mm_mul_ps(felec,sw) , _mm_mul_ps(rinv00,_mm_mul_ps(velec,dsw)) );
717 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
721 fscal = _mm_and_ps(fscal,cutoff_mask);
723 fscal = _mm_andnot_ps(dummy_mask,fscal);
725 /* Calculate temporary vectorial force */
726 tx = _mm_mul_ps(fscal,dx00);
727 ty = _mm_mul_ps(fscal,dy00);
728 tz = _mm_mul_ps(fscal,dz00);
730 /* Update vectorial force */
731 fix0 = _mm_add_ps(fix0,tx);
732 fiy0 = _mm_add_ps(fiy0,ty);
733 fiz0 = _mm_add_ps(fiz0,tz);
735 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
736 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
737 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
738 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
739 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
743 /* Inner loop uses 63 flops */
746 /* End of innermost loop */
748 gmx_mm_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
749 f+i_coord_offset,fshift+i_shift_offset);
751 /* Increment number of inner iterations */
752 inneriter += j_index_end - j_index_start;
754 /* Outer loop uses 7 flops */
757 /* Increment number of outer iterations */
760 /* Update outer/inner flops */
762 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_F,outeriter*7 + inneriter*63);