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_ElecCSTab_VdwNone_GeomP1P1_VF_sse2_single
38 * Electrostatics interaction: CubicSplineTable
39 * VdW interaction: None
40 * Geometry: Particle-Particle
41 * Calculate force/pot: PotentialAndForce
44 nb_kernel_ElecCSTab_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 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 __m128i ifour = _mm_set1_epi32(4);
78 __m128 rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF;
80 __m128 dummy_mask,cutoff_mask;
81 __m128 signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
82 __m128 one = _mm_set1_ps(1.0);
83 __m128 two = _mm_set1_ps(2.0);
89 jindex = nlist->jindex;
91 shiftidx = nlist->shift;
93 shiftvec = fr->shift_vec[0];
94 fshift = fr->fshift[0];
95 facel = _mm_set1_ps(fr->epsfac);
96 charge = mdatoms->chargeA;
98 vftab = kernel_data->table_elec->data;
99 vftabscale = _mm_set1_ps(kernel_data->table_elec->scale);
101 /* Avoid stupid compiler warnings */
102 jnrA = jnrB = jnrC = jnrD = 0;
111 for(iidx=0;iidx<4*DIM;iidx++)
116 /* Start outer loop over neighborlists */
117 for(iidx=0; iidx<nri; iidx++)
119 /* Load shift vector for this list */
120 i_shift_offset = DIM*shiftidx[iidx];
122 /* Load limits for loop over neighbors */
123 j_index_start = jindex[iidx];
124 j_index_end = jindex[iidx+1];
126 /* Get outer coordinate index */
128 i_coord_offset = DIM*inr;
130 /* Load i particle coords and add shift vector */
131 gmx_mm_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
133 fix0 = _mm_setzero_ps();
134 fiy0 = _mm_setzero_ps();
135 fiz0 = _mm_setzero_ps();
137 /* Load parameters for i particles */
138 iq0 = _mm_mul_ps(facel,_mm_load1_ps(charge+inr+0));
140 /* Reset potential sums */
141 velecsum = _mm_setzero_ps();
143 /* Start inner kernel loop */
144 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
147 /* Get j neighbor index, and coordinate index */
152 j_coord_offsetA = DIM*jnrA;
153 j_coord_offsetB = DIM*jnrB;
154 j_coord_offsetC = DIM*jnrC;
155 j_coord_offsetD = DIM*jnrD;
157 /* load j atom coordinates */
158 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
159 x+j_coord_offsetC,x+j_coord_offsetD,
162 /* Calculate displacement vector */
163 dx00 = _mm_sub_ps(ix0,jx0);
164 dy00 = _mm_sub_ps(iy0,jy0);
165 dz00 = _mm_sub_ps(iz0,jz0);
167 /* Calculate squared distance and things based on it */
168 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
170 rinv00 = gmx_mm_invsqrt_ps(rsq00);
172 /* Load parameters for j particles */
173 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
174 charge+jnrC+0,charge+jnrD+0);
176 /**************************
177 * CALCULATE INTERACTIONS *
178 **************************/
180 r00 = _mm_mul_ps(rsq00,rinv00);
182 /* Compute parameters for interactions between i and j atoms */
183 qq00 = _mm_mul_ps(iq0,jq0);
185 /* Calculate table index by multiplying r with table scale and truncate to integer */
186 rt = _mm_mul_ps(r00,vftabscale);
187 vfitab = _mm_cvttps_epi32(rt);
188 vfeps = _mm_sub_ps(rt,_mm_cvtepi32_ps(vfitab));
189 vfitab = _mm_slli_epi32(vfitab,2);
191 /* CUBIC SPLINE TABLE ELECTROSTATICS */
192 Y = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,0) );
193 F = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,1) );
194 G = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,2) );
195 H = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,3) );
196 _MM_TRANSPOSE4_PS(Y,F,G,H);
197 Heps = _mm_mul_ps(vfeps,H);
198 Fp = _mm_add_ps(F,_mm_mul_ps(vfeps,_mm_add_ps(G,Heps)));
199 VV = _mm_add_ps(Y,_mm_mul_ps(vfeps,Fp));
200 velec = _mm_mul_ps(qq00,VV);
201 FF = _mm_add_ps(Fp,_mm_mul_ps(vfeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
202 felec = _mm_xor_ps(signbit,_mm_mul_ps(_mm_mul_ps(qq00,FF),_mm_mul_ps(vftabscale,rinv00)));
204 /* Update potential sum for this i atom from the interaction with this j atom. */
205 velecsum = _mm_add_ps(velecsum,velec);
209 /* Calculate temporary vectorial force */
210 tx = _mm_mul_ps(fscal,dx00);
211 ty = _mm_mul_ps(fscal,dy00);
212 tz = _mm_mul_ps(fscal,dz00);
214 /* Update vectorial force */
215 fix0 = _mm_add_ps(fix0,tx);
216 fiy0 = _mm_add_ps(fiy0,ty);
217 fiz0 = _mm_add_ps(fiz0,tz);
219 fjptrA = f+j_coord_offsetA;
220 fjptrB = f+j_coord_offsetB;
221 fjptrC = f+j_coord_offsetC;
222 fjptrD = f+j_coord_offsetD;
223 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
225 /* Inner loop uses 43 flops */
231 /* Get j neighbor index, and coordinate index */
232 jnrlistA = jjnr[jidx];
233 jnrlistB = jjnr[jidx+1];
234 jnrlistC = jjnr[jidx+2];
235 jnrlistD = jjnr[jidx+3];
236 /* Sign of each element will be negative for non-real atoms.
237 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
238 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
240 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
241 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
242 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
243 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
244 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
245 j_coord_offsetA = DIM*jnrA;
246 j_coord_offsetB = DIM*jnrB;
247 j_coord_offsetC = DIM*jnrC;
248 j_coord_offsetD = DIM*jnrD;
250 /* load j atom coordinates */
251 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
252 x+j_coord_offsetC,x+j_coord_offsetD,
255 /* Calculate displacement vector */
256 dx00 = _mm_sub_ps(ix0,jx0);
257 dy00 = _mm_sub_ps(iy0,jy0);
258 dz00 = _mm_sub_ps(iz0,jz0);
260 /* Calculate squared distance and things based on it */
261 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
263 rinv00 = gmx_mm_invsqrt_ps(rsq00);
265 /* Load parameters for j particles */
266 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
267 charge+jnrC+0,charge+jnrD+0);
269 /**************************
270 * CALCULATE INTERACTIONS *
271 **************************/
273 r00 = _mm_mul_ps(rsq00,rinv00);
274 r00 = _mm_andnot_ps(dummy_mask,r00);
276 /* Compute parameters for interactions between i and j atoms */
277 qq00 = _mm_mul_ps(iq0,jq0);
279 /* Calculate table index by multiplying r with table scale and truncate to integer */
280 rt = _mm_mul_ps(r00,vftabscale);
281 vfitab = _mm_cvttps_epi32(rt);
282 vfeps = _mm_sub_ps(rt,_mm_cvtepi32_ps(vfitab));
283 vfitab = _mm_slli_epi32(vfitab,2);
285 /* CUBIC SPLINE TABLE ELECTROSTATICS */
286 Y = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,0) );
287 F = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,1) );
288 G = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,2) );
289 H = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,3) );
290 _MM_TRANSPOSE4_PS(Y,F,G,H);
291 Heps = _mm_mul_ps(vfeps,H);
292 Fp = _mm_add_ps(F,_mm_mul_ps(vfeps,_mm_add_ps(G,Heps)));
293 VV = _mm_add_ps(Y,_mm_mul_ps(vfeps,Fp));
294 velec = _mm_mul_ps(qq00,VV);
295 FF = _mm_add_ps(Fp,_mm_mul_ps(vfeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
296 felec = _mm_xor_ps(signbit,_mm_mul_ps(_mm_mul_ps(qq00,FF),_mm_mul_ps(vftabscale,rinv00)));
298 /* Update potential sum for this i atom from the interaction with this j atom. */
299 velec = _mm_andnot_ps(dummy_mask,velec);
300 velecsum = _mm_add_ps(velecsum,velec);
304 fscal = _mm_andnot_ps(dummy_mask,fscal);
306 /* Calculate temporary vectorial force */
307 tx = _mm_mul_ps(fscal,dx00);
308 ty = _mm_mul_ps(fscal,dy00);
309 tz = _mm_mul_ps(fscal,dz00);
311 /* Update vectorial force */
312 fix0 = _mm_add_ps(fix0,tx);
313 fiy0 = _mm_add_ps(fiy0,ty);
314 fiz0 = _mm_add_ps(fiz0,tz);
316 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
317 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
318 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
319 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
320 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
322 /* Inner loop uses 44 flops */
325 /* End of innermost loop */
327 gmx_mm_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
328 f+i_coord_offset,fshift+i_shift_offset);
331 /* Update potential energies */
332 gmx_mm_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
334 /* Increment number of inner iterations */
335 inneriter += j_index_end - j_index_start;
337 /* Outer loop uses 8 flops */
340 /* Increment number of outer iterations */
343 /* Update outer/inner flops */
345 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VF,outeriter*8 + inneriter*44);
348 * Gromacs nonbonded kernel: nb_kernel_ElecCSTab_VdwNone_GeomP1P1_F_sse2_single
349 * Electrostatics interaction: CubicSplineTable
350 * VdW interaction: None
351 * Geometry: Particle-Particle
352 * Calculate force/pot: Force
355 nb_kernel_ElecCSTab_VdwNone_GeomP1P1_F_sse2_single
356 (t_nblist * gmx_restrict nlist,
357 rvec * gmx_restrict xx,
358 rvec * gmx_restrict ff,
359 t_forcerec * gmx_restrict fr,
360 t_mdatoms * gmx_restrict mdatoms,
361 nb_kernel_data_t * gmx_restrict kernel_data,
362 t_nrnb * gmx_restrict nrnb)
364 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
365 * just 0 for non-waters.
366 * Suffixes A,B,C,D refer to j loop unrolling done with SSE, e.g. for the four different
367 * jnr indices corresponding to data put in the four positions in the SIMD register.
369 int i_shift_offset,i_coord_offset,outeriter,inneriter;
370 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
371 int jnrA,jnrB,jnrC,jnrD;
372 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
373 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
374 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
376 real *shiftvec,*fshift,*x,*f;
377 real *fjptrA,*fjptrB,*fjptrC,*fjptrD;
379 __m128 tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
381 __m128 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
382 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
383 __m128 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
384 __m128 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
385 __m128 velec,felec,velecsum,facel,crf,krf,krf2;
388 __m128i ifour = _mm_set1_epi32(4);
389 __m128 rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF;
391 __m128 dummy_mask,cutoff_mask;
392 __m128 signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
393 __m128 one = _mm_set1_ps(1.0);
394 __m128 two = _mm_set1_ps(2.0);
400 jindex = nlist->jindex;
402 shiftidx = nlist->shift;
404 shiftvec = fr->shift_vec[0];
405 fshift = fr->fshift[0];
406 facel = _mm_set1_ps(fr->epsfac);
407 charge = mdatoms->chargeA;
409 vftab = kernel_data->table_elec->data;
410 vftabscale = _mm_set1_ps(kernel_data->table_elec->scale);
412 /* Avoid stupid compiler warnings */
413 jnrA = jnrB = jnrC = jnrD = 0;
422 for(iidx=0;iidx<4*DIM;iidx++)
427 /* Start outer loop over neighborlists */
428 for(iidx=0; iidx<nri; iidx++)
430 /* Load shift vector for this list */
431 i_shift_offset = DIM*shiftidx[iidx];
433 /* Load limits for loop over neighbors */
434 j_index_start = jindex[iidx];
435 j_index_end = jindex[iidx+1];
437 /* Get outer coordinate index */
439 i_coord_offset = DIM*inr;
441 /* Load i particle coords and add shift vector */
442 gmx_mm_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
444 fix0 = _mm_setzero_ps();
445 fiy0 = _mm_setzero_ps();
446 fiz0 = _mm_setzero_ps();
448 /* Load parameters for i particles */
449 iq0 = _mm_mul_ps(facel,_mm_load1_ps(charge+inr+0));
451 /* Start inner kernel loop */
452 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
455 /* Get j neighbor index, and coordinate index */
460 j_coord_offsetA = DIM*jnrA;
461 j_coord_offsetB = DIM*jnrB;
462 j_coord_offsetC = DIM*jnrC;
463 j_coord_offsetD = DIM*jnrD;
465 /* load j atom coordinates */
466 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
467 x+j_coord_offsetC,x+j_coord_offsetD,
470 /* Calculate displacement vector */
471 dx00 = _mm_sub_ps(ix0,jx0);
472 dy00 = _mm_sub_ps(iy0,jy0);
473 dz00 = _mm_sub_ps(iz0,jz0);
475 /* Calculate squared distance and things based on it */
476 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
478 rinv00 = gmx_mm_invsqrt_ps(rsq00);
480 /* Load parameters for j particles */
481 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
482 charge+jnrC+0,charge+jnrD+0);
484 /**************************
485 * CALCULATE INTERACTIONS *
486 **************************/
488 r00 = _mm_mul_ps(rsq00,rinv00);
490 /* Compute parameters for interactions between i and j atoms */
491 qq00 = _mm_mul_ps(iq0,jq0);
493 /* Calculate table index by multiplying r with table scale and truncate to integer */
494 rt = _mm_mul_ps(r00,vftabscale);
495 vfitab = _mm_cvttps_epi32(rt);
496 vfeps = _mm_sub_ps(rt,_mm_cvtepi32_ps(vfitab));
497 vfitab = _mm_slli_epi32(vfitab,2);
499 /* CUBIC SPLINE TABLE ELECTROSTATICS */
500 Y = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,0) );
501 F = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,1) );
502 G = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,2) );
503 H = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,3) );
504 _MM_TRANSPOSE4_PS(Y,F,G,H);
505 Heps = _mm_mul_ps(vfeps,H);
506 Fp = _mm_add_ps(F,_mm_mul_ps(vfeps,_mm_add_ps(G,Heps)));
507 FF = _mm_add_ps(Fp,_mm_mul_ps(vfeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
508 felec = _mm_xor_ps(signbit,_mm_mul_ps(_mm_mul_ps(qq00,FF),_mm_mul_ps(vftabscale,rinv00)));
512 /* Calculate temporary vectorial force */
513 tx = _mm_mul_ps(fscal,dx00);
514 ty = _mm_mul_ps(fscal,dy00);
515 tz = _mm_mul_ps(fscal,dz00);
517 /* Update vectorial force */
518 fix0 = _mm_add_ps(fix0,tx);
519 fiy0 = _mm_add_ps(fiy0,ty);
520 fiz0 = _mm_add_ps(fiz0,tz);
522 fjptrA = f+j_coord_offsetA;
523 fjptrB = f+j_coord_offsetB;
524 fjptrC = f+j_coord_offsetC;
525 fjptrD = f+j_coord_offsetD;
526 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
528 /* Inner loop uses 39 flops */
534 /* Get j neighbor index, and coordinate index */
535 jnrlistA = jjnr[jidx];
536 jnrlistB = jjnr[jidx+1];
537 jnrlistC = jjnr[jidx+2];
538 jnrlistD = jjnr[jidx+3];
539 /* Sign of each element will be negative for non-real atoms.
540 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
541 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
543 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
544 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
545 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
546 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
547 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
548 j_coord_offsetA = DIM*jnrA;
549 j_coord_offsetB = DIM*jnrB;
550 j_coord_offsetC = DIM*jnrC;
551 j_coord_offsetD = DIM*jnrD;
553 /* load j atom coordinates */
554 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
555 x+j_coord_offsetC,x+j_coord_offsetD,
558 /* Calculate displacement vector */
559 dx00 = _mm_sub_ps(ix0,jx0);
560 dy00 = _mm_sub_ps(iy0,jy0);
561 dz00 = _mm_sub_ps(iz0,jz0);
563 /* Calculate squared distance and things based on it */
564 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
566 rinv00 = gmx_mm_invsqrt_ps(rsq00);
568 /* Load parameters for j particles */
569 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
570 charge+jnrC+0,charge+jnrD+0);
572 /**************************
573 * CALCULATE INTERACTIONS *
574 **************************/
576 r00 = _mm_mul_ps(rsq00,rinv00);
577 r00 = _mm_andnot_ps(dummy_mask,r00);
579 /* Compute parameters for interactions between i and j atoms */
580 qq00 = _mm_mul_ps(iq0,jq0);
582 /* Calculate table index by multiplying r with table scale and truncate to integer */
583 rt = _mm_mul_ps(r00,vftabscale);
584 vfitab = _mm_cvttps_epi32(rt);
585 vfeps = _mm_sub_ps(rt,_mm_cvtepi32_ps(vfitab));
586 vfitab = _mm_slli_epi32(vfitab,2);
588 /* CUBIC SPLINE TABLE ELECTROSTATICS */
589 Y = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,0) );
590 F = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,1) );
591 G = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,2) );
592 H = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,3) );
593 _MM_TRANSPOSE4_PS(Y,F,G,H);
594 Heps = _mm_mul_ps(vfeps,H);
595 Fp = _mm_add_ps(F,_mm_mul_ps(vfeps,_mm_add_ps(G,Heps)));
596 FF = _mm_add_ps(Fp,_mm_mul_ps(vfeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
597 felec = _mm_xor_ps(signbit,_mm_mul_ps(_mm_mul_ps(qq00,FF),_mm_mul_ps(vftabscale,rinv00)));
601 fscal = _mm_andnot_ps(dummy_mask,fscal);
603 /* Calculate temporary vectorial force */
604 tx = _mm_mul_ps(fscal,dx00);
605 ty = _mm_mul_ps(fscal,dy00);
606 tz = _mm_mul_ps(fscal,dz00);
608 /* Update vectorial force */
609 fix0 = _mm_add_ps(fix0,tx);
610 fiy0 = _mm_add_ps(fiy0,ty);
611 fiz0 = _mm_add_ps(fiz0,tz);
613 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
614 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
615 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
616 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
617 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
619 /* Inner loop uses 40 flops */
622 /* End of innermost loop */
624 gmx_mm_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
625 f+i_coord_offset,fshift+i_shift_offset);
627 /* Increment number of inner iterations */
628 inneriter += j_index_end - j_index_start;
630 /* Outer loop uses 7 flops */
633 /* Increment number of outer iterations */
636 /* Update outer/inner flops */
638 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_F,outeriter*7 + inneriter*40);