2 * Note: this file was generated by the Gromacs sse2_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_sse2_single.h"
34 #include "kernelutil_x86_sse2_single.h"
37 * Gromacs nonbonded kernel: nb_kernel_ElecEw_VdwNone_GeomP1P1_VF_sse2_single
38 * Electrostatics interaction: Ewald
39 * VdW interaction: None
40 * Geometry: Particle-Particle
41 * Calculate force/pot: PotentialAndForce
44 nb_kernel_ElecEw_VdwNone_GeomP1P1_VF_sse2_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 j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
62 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
63 real shX,shY,shZ,rcutoff_scalar;
64 real *shiftvec,*fshift,*x,*f;
65 __m128 tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
67 __m128 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
68 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
69 __m128 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
70 __m128 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
71 __m128 velec,felec,velecsum,facel,crf,krf,krf2;
74 __m128 ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
76 __m128 dummy_mask,cutoff_mask;
77 __m128 signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
78 __m128 one = _mm_set1_ps(1.0);
79 __m128 two = _mm_set1_ps(2.0);
85 jindex = nlist->jindex;
87 shiftidx = nlist->shift;
89 shiftvec = fr->shift_vec[0];
90 fshift = fr->fshift[0];
91 facel = _mm_set1_ps(fr->epsfac);
92 charge = mdatoms->chargeA;
94 sh_ewald = _mm_set1_ps(fr->ic->sh_ewald);
95 ewtab = fr->ic->tabq_coul_FDV0;
96 ewtabscale = _mm_set1_ps(fr->ic->tabq_scale);
97 ewtabhalfspace = _mm_set1_ps(0.5/fr->ic->tabq_scale);
99 /* Avoid stupid compiler warnings */
100 jnrA = jnrB = jnrC = jnrD = 0;
109 /* Start outer loop over neighborlists */
110 for(iidx=0; iidx<nri; iidx++)
112 /* Load shift vector for this list */
113 i_shift_offset = DIM*shiftidx[iidx];
114 shX = shiftvec[i_shift_offset+XX];
115 shY = shiftvec[i_shift_offset+YY];
116 shZ = shiftvec[i_shift_offset+ZZ];
118 /* Load limits for loop over neighbors */
119 j_index_start = jindex[iidx];
120 j_index_end = jindex[iidx+1];
122 /* Get outer coordinate index */
124 i_coord_offset = DIM*inr;
126 /* Load i particle coords and add shift vector */
127 ix0 = _mm_set1_ps(shX + x[i_coord_offset+DIM*0+XX]);
128 iy0 = _mm_set1_ps(shY + x[i_coord_offset+DIM*0+YY]);
129 iz0 = _mm_set1_ps(shZ + x[i_coord_offset+DIM*0+ZZ]);
131 fix0 = _mm_setzero_ps();
132 fiy0 = _mm_setzero_ps();
133 fiz0 = _mm_setzero_ps();
135 /* Load parameters for i particles */
136 iq0 = _mm_mul_ps(facel,_mm_load1_ps(charge+inr+0));
138 /* Reset potential sums */
139 velecsum = _mm_setzero_ps();
141 /* Start inner kernel loop */
142 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
145 /* Get j neighbor index, and coordinate index */
151 j_coord_offsetA = DIM*jnrA;
152 j_coord_offsetB = DIM*jnrB;
153 j_coord_offsetC = DIM*jnrC;
154 j_coord_offsetD = DIM*jnrD;
156 /* load j atom coordinates */
157 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
158 x+j_coord_offsetC,x+j_coord_offsetD,
161 /* Calculate displacement vector */
162 dx00 = _mm_sub_ps(ix0,jx0);
163 dy00 = _mm_sub_ps(iy0,jy0);
164 dz00 = _mm_sub_ps(iz0,jz0);
166 /* Calculate squared distance and things based on it */
167 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
169 rinv00 = gmx_mm_invsqrt_ps(rsq00);
171 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
173 /* Load parameters for j particles */
174 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
175 charge+jnrC+0,charge+jnrD+0);
177 /**************************
178 * CALCULATE INTERACTIONS *
179 **************************/
181 r00 = _mm_mul_ps(rsq00,rinv00);
183 /* Compute parameters for interactions between i and j atoms */
184 qq00 = _mm_mul_ps(iq0,jq0);
186 /* EWALD ELECTROSTATICS */
188 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
189 ewrt = _mm_mul_ps(r00,ewtabscale);
190 ewitab = _mm_cvttps_epi32(ewrt);
191 eweps = _mm_sub_ps(ewrt,_mm_cvtepi32_ps(ewitab));
192 ewitab = _mm_slli_epi32(ewitab,2);
193 ewtabF = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,0) );
194 ewtabD = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,1) );
195 ewtabV = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,2) );
196 ewtabFn = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,3) );
197 _MM_TRANSPOSE4_PS(ewtabF,ewtabD,ewtabV,ewtabFn);
198 felec = _mm_add_ps(ewtabF,_mm_mul_ps(eweps,ewtabD));
199 velec = _mm_sub_ps(ewtabV,_mm_mul_ps(_mm_mul_ps(ewtabhalfspace,eweps),_mm_add_ps(ewtabF,felec)));
200 velec = _mm_mul_ps(qq00,_mm_sub_ps(rinv00,velec));
201 felec = _mm_mul_ps(_mm_mul_ps(qq00,rinv00),_mm_sub_ps(rinvsq00,felec));
203 /* Update potential sum for this i atom from the interaction with this j atom. */
204 velecsum = _mm_add_ps(velecsum,velec);
208 /* Calculate temporary vectorial force */
209 tx = _mm_mul_ps(fscal,dx00);
210 ty = _mm_mul_ps(fscal,dy00);
211 tz = _mm_mul_ps(fscal,dz00);
213 /* Update vectorial force */
214 fix0 = _mm_add_ps(fix0,tx);
215 fiy0 = _mm_add_ps(fiy0,ty);
216 fiz0 = _mm_add_ps(fiz0,tz);
218 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(f+j_coord_offsetA,f+j_coord_offsetB,
219 f+j_coord_offsetC,f+j_coord_offsetD,
222 /* Inner loop uses 41 flops */
228 /* Get j neighbor index, and coordinate index */
234 /* Sign of each element will be negative for non-real atoms.
235 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
236 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
238 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
239 jnrA = (jnrA>=0) ? jnrA : 0;
240 jnrB = (jnrB>=0) ? jnrB : 0;
241 jnrC = (jnrC>=0) ? jnrC : 0;
242 jnrD = (jnrD>=0) ? jnrD : 0;
244 j_coord_offsetA = DIM*jnrA;
245 j_coord_offsetB = DIM*jnrB;
246 j_coord_offsetC = DIM*jnrC;
247 j_coord_offsetD = DIM*jnrD;
249 /* load j atom coordinates */
250 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
251 x+j_coord_offsetC,x+j_coord_offsetD,
254 /* Calculate displacement vector */
255 dx00 = _mm_sub_ps(ix0,jx0);
256 dy00 = _mm_sub_ps(iy0,jy0);
257 dz00 = _mm_sub_ps(iz0,jz0);
259 /* Calculate squared distance and things based on it */
260 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
262 rinv00 = gmx_mm_invsqrt_ps(rsq00);
264 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
266 /* Load parameters for j particles */
267 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
268 charge+jnrC+0,charge+jnrD+0);
270 /**************************
271 * CALCULATE INTERACTIONS *
272 **************************/
274 r00 = _mm_mul_ps(rsq00,rinv00);
275 r00 = _mm_andnot_ps(dummy_mask,r00);
277 /* Compute parameters for interactions between i and j atoms */
278 qq00 = _mm_mul_ps(iq0,jq0);
280 /* EWALD ELECTROSTATICS */
282 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
283 ewrt = _mm_mul_ps(r00,ewtabscale);
284 ewitab = _mm_cvttps_epi32(ewrt);
285 eweps = _mm_sub_ps(ewrt,_mm_cvtepi32_ps(ewitab));
286 ewitab = _mm_slli_epi32(ewitab,2);
287 ewtabF = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,0) );
288 ewtabD = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,1) );
289 ewtabV = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,2) );
290 ewtabFn = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,3) );
291 _MM_TRANSPOSE4_PS(ewtabF,ewtabD,ewtabV,ewtabFn);
292 felec = _mm_add_ps(ewtabF,_mm_mul_ps(eweps,ewtabD));
293 velec = _mm_sub_ps(ewtabV,_mm_mul_ps(_mm_mul_ps(ewtabhalfspace,eweps),_mm_add_ps(ewtabF,felec)));
294 velec = _mm_mul_ps(qq00,_mm_sub_ps(rinv00,velec));
295 felec = _mm_mul_ps(_mm_mul_ps(qq00,rinv00),_mm_sub_ps(rinvsq00,felec));
297 /* Update potential sum for this i atom from the interaction with this j atom. */
298 velec = _mm_andnot_ps(dummy_mask,velec);
299 velecsum = _mm_add_ps(velecsum,velec);
303 fscal = _mm_andnot_ps(dummy_mask,fscal);
305 /* Calculate temporary vectorial force */
306 tx = _mm_mul_ps(fscal,dx00);
307 ty = _mm_mul_ps(fscal,dy00);
308 tz = _mm_mul_ps(fscal,dz00);
310 /* Update vectorial force */
311 fix0 = _mm_add_ps(fix0,tx);
312 fiy0 = _mm_add_ps(fiy0,ty);
313 fiz0 = _mm_add_ps(fiz0,tz);
315 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(f+j_coord_offsetA,f+j_coord_offsetB,
316 f+j_coord_offsetC,f+j_coord_offsetD,
319 /* Inner loop uses 42 flops */
322 /* End of innermost loop */
324 gmx_mm_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
325 f+i_coord_offset,fshift+i_shift_offset);
328 /* Update potential energies */
329 gmx_mm_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
331 /* Increment number of inner iterations */
332 inneriter += j_index_end - j_index_start;
334 /* Outer loop uses 11 flops */
337 /* Increment number of outer iterations */
340 /* Update outer/inner flops */
342 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VF,outeriter*11 + inneriter*42);
345 * Gromacs nonbonded kernel: nb_kernel_ElecEw_VdwNone_GeomP1P1_F_sse2_single
346 * Electrostatics interaction: Ewald
347 * VdW interaction: None
348 * Geometry: Particle-Particle
349 * Calculate force/pot: Force
352 nb_kernel_ElecEw_VdwNone_GeomP1P1_F_sse2_single
353 (t_nblist * gmx_restrict nlist,
354 rvec * gmx_restrict xx,
355 rvec * gmx_restrict ff,
356 t_forcerec * gmx_restrict fr,
357 t_mdatoms * gmx_restrict mdatoms,
358 nb_kernel_data_t * gmx_restrict kernel_data,
359 t_nrnb * gmx_restrict nrnb)
361 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
362 * just 0 for non-waters.
363 * Suffixes A,B,C,D refer to j loop unrolling done with SSE, e.g. for the four different
364 * jnr indices corresponding to data put in the four positions in the SIMD register.
366 int i_shift_offset,i_coord_offset,outeriter,inneriter;
367 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
368 int jnrA,jnrB,jnrC,jnrD;
369 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
370 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
371 real shX,shY,shZ,rcutoff_scalar;
372 real *shiftvec,*fshift,*x,*f;
373 __m128 tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
375 __m128 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
376 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
377 __m128 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
378 __m128 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
379 __m128 velec,felec,velecsum,facel,crf,krf,krf2;
382 __m128 ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
384 __m128 dummy_mask,cutoff_mask;
385 __m128 signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
386 __m128 one = _mm_set1_ps(1.0);
387 __m128 two = _mm_set1_ps(2.0);
393 jindex = nlist->jindex;
395 shiftidx = nlist->shift;
397 shiftvec = fr->shift_vec[0];
398 fshift = fr->fshift[0];
399 facel = _mm_set1_ps(fr->epsfac);
400 charge = mdatoms->chargeA;
402 sh_ewald = _mm_set1_ps(fr->ic->sh_ewald);
403 ewtab = fr->ic->tabq_coul_F;
404 ewtabscale = _mm_set1_ps(fr->ic->tabq_scale);
405 ewtabhalfspace = _mm_set1_ps(0.5/fr->ic->tabq_scale);
407 /* Avoid stupid compiler warnings */
408 jnrA = jnrB = jnrC = jnrD = 0;
417 /* Start outer loop over neighborlists */
418 for(iidx=0; iidx<nri; iidx++)
420 /* Load shift vector for this list */
421 i_shift_offset = DIM*shiftidx[iidx];
422 shX = shiftvec[i_shift_offset+XX];
423 shY = shiftvec[i_shift_offset+YY];
424 shZ = shiftvec[i_shift_offset+ZZ];
426 /* Load limits for loop over neighbors */
427 j_index_start = jindex[iidx];
428 j_index_end = jindex[iidx+1];
430 /* Get outer coordinate index */
432 i_coord_offset = DIM*inr;
434 /* Load i particle coords and add shift vector */
435 ix0 = _mm_set1_ps(shX + x[i_coord_offset+DIM*0+XX]);
436 iy0 = _mm_set1_ps(shY + x[i_coord_offset+DIM*0+YY]);
437 iz0 = _mm_set1_ps(shZ + x[i_coord_offset+DIM*0+ZZ]);
439 fix0 = _mm_setzero_ps();
440 fiy0 = _mm_setzero_ps();
441 fiz0 = _mm_setzero_ps();
443 /* Load parameters for i particles */
444 iq0 = _mm_mul_ps(facel,_mm_load1_ps(charge+inr+0));
446 /* Start inner kernel loop */
447 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
450 /* Get j neighbor index, and coordinate index */
456 j_coord_offsetA = DIM*jnrA;
457 j_coord_offsetB = DIM*jnrB;
458 j_coord_offsetC = DIM*jnrC;
459 j_coord_offsetD = DIM*jnrD;
461 /* load j atom coordinates */
462 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
463 x+j_coord_offsetC,x+j_coord_offsetD,
466 /* Calculate displacement vector */
467 dx00 = _mm_sub_ps(ix0,jx0);
468 dy00 = _mm_sub_ps(iy0,jy0);
469 dz00 = _mm_sub_ps(iz0,jz0);
471 /* Calculate squared distance and things based on it */
472 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
474 rinv00 = gmx_mm_invsqrt_ps(rsq00);
476 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
478 /* Load parameters for j particles */
479 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
480 charge+jnrC+0,charge+jnrD+0);
482 /**************************
483 * CALCULATE INTERACTIONS *
484 **************************/
486 r00 = _mm_mul_ps(rsq00,rinv00);
488 /* Compute parameters for interactions between i and j atoms */
489 qq00 = _mm_mul_ps(iq0,jq0);
491 /* EWALD ELECTROSTATICS */
493 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
494 ewrt = _mm_mul_ps(r00,ewtabscale);
495 ewitab = _mm_cvttps_epi32(ewrt);
496 eweps = _mm_sub_ps(ewrt,_mm_cvtepi32_ps(ewitab));
497 gmx_mm_load_4pair_swizzle_ps(ewtab+gmx_mm_extract_epi32(ewitab,0),ewtab+gmx_mm_extract_epi32(ewitab,1),
498 ewtab+gmx_mm_extract_epi32(ewitab,2),ewtab+gmx_mm_extract_epi32(ewitab,3),
500 felec = _mm_add_ps(_mm_mul_ps( _mm_sub_ps(one,eweps),ewtabF),_mm_mul_ps(eweps,ewtabFn));
501 felec = _mm_mul_ps(_mm_mul_ps(qq00,rinv00),_mm_sub_ps(rinvsq00,felec));
505 /* Calculate temporary vectorial force */
506 tx = _mm_mul_ps(fscal,dx00);
507 ty = _mm_mul_ps(fscal,dy00);
508 tz = _mm_mul_ps(fscal,dz00);
510 /* Update vectorial force */
511 fix0 = _mm_add_ps(fix0,tx);
512 fiy0 = _mm_add_ps(fiy0,ty);
513 fiz0 = _mm_add_ps(fiz0,tz);
515 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(f+j_coord_offsetA,f+j_coord_offsetB,
516 f+j_coord_offsetC,f+j_coord_offsetD,
519 /* Inner loop uses 36 flops */
525 /* Get j neighbor index, and coordinate index */
531 /* Sign of each element will be negative for non-real atoms.
532 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
533 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
535 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
536 jnrA = (jnrA>=0) ? jnrA : 0;
537 jnrB = (jnrB>=0) ? jnrB : 0;
538 jnrC = (jnrC>=0) ? jnrC : 0;
539 jnrD = (jnrD>=0) ? jnrD : 0;
541 j_coord_offsetA = DIM*jnrA;
542 j_coord_offsetB = DIM*jnrB;
543 j_coord_offsetC = DIM*jnrC;
544 j_coord_offsetD = DIM*jnrD;
546 /* load j atom coordinates */
547 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
548 x+j_coord_offsetC,x+j_coord_offsetD,
551 /* Calculate displacement vector */
552 dx00 = _mm_sub_ps(ix0,jx0);
553 dy00 = _mm_sub_ps(iy0,jy0);
554 dz00 = _mm_sub_ps(iz0,jz0);
556 /* Calculate squared distance and things based on it */
557 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
559 rinv00 = gmx_mm_invsqrt_ps(rsq00);
561 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
563 /* Load parameters for j particles */
564 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
565 charge+jnrC+0,charge+jnrD+0);
567 /**************************
568 * CALCULATE INTERACTIONS *
569 **************************/
571 r00 = _mm_mul_ps(rsq00,rinv00);
572 r00 = _mm_andnot_ps(dummy_mask,r00);
574 /* Compute parameters for interactions between i and j atoms */
575 qq00 = _mm_mul_ps(iq0,jq0);
577 /* EWALD ELECTROSTATICS */
579 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
580 ewrt = _mm_mul_ps(r00,ewtabscale);
581 ewitab = _mm_cvttps_epi32(ewrt);
582 eweps = _mm_sub_ps(ewrt,_mm_cvtepi32_ps(ewitab));
583 gmx_mm_load_4pair_swizzle_ps(ewtab+gmx_mm_extract_epi32(ewitab,0),ewtab+gmx_mm_extract_epi32(ewitab,1),
584 ewtab+gmx_mm_extract_epi32(ewitab,2),ewtab+gmx_mm_extract_epi32(ewitab,3),
586 felec = _mm_add_ps(_mm_mul_ps( _mm_sub_ps(one,eweps),ewtabF),_mm_mul_ps(eweps,ewtabFn));
587 felec = _mm_mul_ps(_mm_mul_ps(qq00,rinv00),_mm_sub_ps(rinvsq00,felec));
591 fscal = _mm_andnot_ps(dummy_mask,fscal);
593 /* Calculate temporary vectorial force */
594 tx = _mm_mul_ps(fscal,dx00);
595 ty = _mm_mul_ps(fscal,dy00);
596 tz = _mm_mul_ps(fscal,dz00);
598 /* Update vectorial force */
599 fix0 = _mm_add_ps(fix0,tx);
600 fiy0 = _mm_add_ps(fiy0,ty);
601 fiz0 = _mm_add_ps(fiz0,tz);
603 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(f+j_coord_offsetA,f+j_coord_offsetB,
604 f+j_coord_offsetC,f+j_coord_offsetD,
607 /* Inner loop uses 37 flops */
610 /* End of innermost loop */
612 gmx_mm_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
613 f+i_coord_offset,fshift+i_shift_offset);
615 /* Increment number of inner iterations */
616 inneriter += j_index_end - j_index_start;
618 /* Outer loop uses 10 flops */
621 /* Increment number of outer iterations */
624 /* Update outer/inner flops */
626 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_F,outeriter*10 + inneriter*37);