2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013, 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 avx_128_fma_single kernel generator.
44 #include "../nb_kernel.h"
45 #include "types/simple.h"
49 #include "gromacs/simd/math_x86_avx_128_fma_single.h"
50 #include "kernelutil_x86_avx_128_fma_single.h"
53 * Gromacs nonbonded kernel: nb_kernel_ElecCSTab_VdwNone_GeomP1P1_VF_avx_128_fma_single
54 * Electrostatics interaction: CubicSplineTable
55 * VdW interaction: None
56 * Geometry: Particle-Particle
57 * Calculate force/pot: PotentialAndForce
60 nb_kernel_ElecCSTab_VdwNone_GeomP1P1_VF_avx_128_fma_single
61 (t_nblist * gmx_restrict nlist,
62 rvec * gmx_restrict xx,
63 rvec * gmx_restrict ff,
64 t_forcerec * gmx_restrict fr,
65 t_mdatoms * gmx_restrict mdatoms,
66 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
67 t_nrnb * gmx_restrict nrnb)
69 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
70 * just 0 for non-waters.
71 * Suffixes A,B,C,D refer to j loop unrolling done with AVX_128, e.g. for the four different
72 * jnr indices corresponding to data put in the four positions in the SIMD register.
74 int i_shift_offset,i_coord_offset,outeriter,inneriter;
75 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
76 int jnrA,jnrB,jnrC,jnrD;
77 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
78 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
79 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
81 real *shiftvec,*fshift,*x,*f;
82 real *fjptrA,*fjptrB,*fjptrC,*fjptrD;
84 __m128 fscal,rcutoff,rcutoff2,jidxall;
86 __m128 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
87 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
88 __m128 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
89 __m128 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
90 __m128 velec,felec,velecsum,facel,crf,krf,krf2;
93 __m128i ifour = _mm_set1_epi32(4);
94 __m128 rt,vfeps,twovfeps,vftabscale,Y,F,G,H,Fp,VV,FF;
96 __m128 dummy_mask,cutoff_mask;
97 __m128 signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
98 __m128 one = _mm_set1_ps(1.0);
99 __m128 two = _mm_set1_ps(2.0);
105 jindex = nlist->jindex;
107 shiftidx = nlist->shift;
109 shiftvec = fr->shift_vec[0];
110 fshift = fr->fshift[0];
111 facel = _mm_set1_ps(fr->epsfac);
112 charge = mdatoms->chargeA;
114 vftab = kernel_data->table_elec->data;
115 vftabscale = _mm_set1_ps(kernel_data->table_elec->scale);
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));
156 /* Reset potential sums */
157 velecsum = _mm_setzero_ps();
159 /* Start inner kernel loop */
160 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
163 /* Get j neighbor index, and coordinate index */
168 j_coord_offsetA = DIM*jnrA;
169 j_coord_offsetB = DIM*jnrB;
170 j_coord_offsetC = DIM*jnrC;
171 j_coord_offsetD = DIM*jnrD;
173 /* load j atom coordinates */
174 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
175 x+j_coord_offsetC,x+j_coord_offsetD,
178 /* Calculate displacement vector */
179 dx00 = _mm_sub_ps(ix0,jx0);
180 dy00 = _mm_sub_ps(iy0,jy0);
181 dz00 = _mm_sub_ps(iz0,jz0);
183 /* Calculate squared distance and things based on it */
184 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
186 rinv00 = gmx_mm_invsqrt_ps(rsq00);
188 /* Load parameters for j particles */
189 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
190 charge+jnrC+0,charge+jnrD+0);
192 /**************************
193 * CALCULATE INTERACTIONS *
194 **************************/
196 r00 = _mm_mul_ps(rsq00,rinv00);
198 /* Compute parameters for interactions between i and j atoms */
199 qq00 = _mm_mul_ps(iq0,jq0);
201 /* Calculate table index by multiplying r with table scale and truncate to integer */
202 rt = _mm_mul_ps(r00,vftabscale);
203 vfitab = _mm_cvttps_epi32(rt);
205 vfeps = _mm_frcz_ps(rt);
207 vfeps = _mm_sub_ps(rt,_mm_round_ps(rt, _MM_FROUND_FLOOR));
209 twovfeps = _mm_add_ps(vfeps,vfeps);
210 vfitab = _mm_slli_epi32(vfitab,2);
212 /* CUBIC SPLINE TABLE ELECTROSTATICS */
213 Y = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,0) );
214 F = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,1) );
215 G = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,2) );
216 H = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,3) );
217 _MM_TRANSPOSE4_PS(Y,F,G,H);
218 Fp = _mm_macc_ps(vfeps,_mm_macc_ps(H,vfeps,G),F);
219 VV = _mm_macc_ps(vfeps,Fp,Y);
220 velec = _mm_mul_ps(qq00,VV);
221 FF = _mm_macc_ps(vfeps,_mm_macc_ps(twovfeps,H,G),Fp);
222 felec = _mm_xor_ps(signbit,_mm_mul_ps(_mm_mul_ps(qq00,FF),_mm_mul_ps(vftabscale,rinv00)));
224 /* Update potential sum for this i atom from the interaction with this j atom. */
225 velecsum = _mm_add_ps(velecsum,velec);
229 /* Update vectorial force */
230 fix0 = _mm_macc_ps(dx00,fscal,fix0);
231 fiy0 = _mm_macc_ps(dy00,fscal,fiy0);
232 fiz0 = _mm_macc_ps(dz00,fscal,fiz0);
234 fjptrA = f+j_coord_offsetA;
235 fjptrB = f+j_coord_offsetB;
236 fjptrC = f+j_coord_offsetC;
237 fjptrD = f+j_coord_offsetD;
238 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,
239 _mm_mul_ps(dx00,fscal),
240 _mm_mul_ps(dy00,fscal),
241 _mm_mul_ps(dz00,fscal));
243 /* Inner loop uses 46 flops */
249 /* Get j neighbor index, and coordinate index */
250 jnrlistA = jjnr[jidx];
251 jnrlistB = jjnr[jidx+1];
252 jnrlistC = jjnr[jidx+2];
253 jnrlistD = jjnr[jidx+3];
254 /* Sign of each element will be negative for non-real atoms.
255 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
256 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
258 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
259 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
260 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
261 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
262 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
263 j_coord_offsetA = DIM*jnrA;
264 j_coord_offsetB = DIM*jnrB;
265 j_coord_offsetC = DIM*jnrC;
266 j_coord_offsetD = DIM*jnrD;
268 /* load j atom coordinates */
269 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
270 x+j_coord_offsetC,x+j_coord_offsetD,
273 /* Calculate displacement vector */
274 dx00 = _mm_sub_ps(ix0,jx0);
275 dy00 = _mm_sub_ps(iy0,jy0);
276 dz00 = _mm_sub_ps(iz0,jz0);
278 /* Calculate squared distance and things based on it */
279 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
281 rinv00 = gmx_mm_invsqrt_ps(rsq00);
283 /* Load parameters for j particles */
284 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
285 charge+jnrC+0,charge+jnrD+0);
287 /**************************
288 * CALCULATE INTERACTIONS *
289 **************************/
291 r00 = _mm_mul_ps(rsq00,rinv00);
292 r00 = _mm_andnot_ps(dummy_mask,r00);
294 /* Compute parameters for interactions between i and j atoms */
295 qq00 = _mm_mul_ps(iq0,jq0);
297 /* Calculate table index by multiplying r with table scale and truncate to integer */
298 rt = _mm_mul_ps(r00,vftabscale);
299 vfitab = _mm_cvttps_epi32(rt);
301 vfeps = _mm_frcz_ps(rt);
303 vfeps = _mm_sub_ps(rt,_mm_round_ps(rt, _MM_FROUND_FLOOR));
305 twovfeps = _mm_add_ps(vfeps,vfeps);
306 vfitab = _mm_slli_epi32(vfitab,2);
308 /* CUBIC SPLINE TABLE ELECTROSTATICS */
309 Y = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,0) );
310 F = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,1) );
311 G = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,2) );
312 H = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,3) );
313 _MM_TRANSPOSE4_PS(Y,F,G,H);
314 Fp = _mm_macc_ps(vfeps,_mm_macc_ps(H,vfeps,G),F);
315 VV = _mm_macc_ps(vfeps,Fp,Y);
316 velec = _mm_mul_ps(qq00,VV);
317 FF = _mm_macc_ps(vfeps,_mm_macc_ps(twovfeps,H,G),Fp);
318 felec = _mm_xor_ps(signbit,_mm_mul_ps(_mm_mul_ps(qq00,FF),_mm_mul_ps(vftabscale,rinv00)));
320 /* Update potential sum for this i atom from the interaction with this j atom. */
321 velec = _mm_andnot_ps(dummy_mask,velec);
322 velecsum = _mm_add_ps(velecsum,velec);
326 fscal = _mm_andnot_ps(dummy_mask,fscal);
328 /* Update vectorial force */
329 fix0 = _mm_macc_ps(dx00,fscal,fix0);
330 fiy0 = _mm_macc_ps(dy00,fscal,fiy0);
331 fiz0 = _mm_macc_ps(dz00,fscal,fiz0);
333 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
334 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
335 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
336 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
337 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,
338 _mm_mul_ps(dx00,fscal),
339 _mm_mul_ps(dy00,fscal),
340 _mm_mul_ps(dz00,fscal));
342 /* Inner loop uses 47 flops */
345 /* End of innermost loop */
347 gmx_mm_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
348 f+i_coord_offset,fshift+i_shift_offset);
351 /* Update potential energies */
352 gmx_mm_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
354 /* Increment number of inner iterations */
355 inneriter += j_index_end - j_index_start;
357 /* Outer loop uses 8 flops */
360 /* Increment number of outer iterations */
363 /* Update outer/inner flops */
365 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VF,outeriter*8 + inneriter*47);
368 * Gromacs nonbonded kernel: nb_kernel_ElecCSTab_VdwNone_GeomP1P1_F_avx_128_fma_single
369 * Electrostatics interaction: CubicSplineTable
370 * VdW interaction: None
371 * Geometry: Particle-Particle
372 * Calculate force/pot: Force
375 nb_kernel_ElecCSTab_VdwNone_GeomP1P1_F_avx_128_fma_single
376 (t_nblist * gmx_restrict nlist,
377 rvec * gmx_restrict xx,
378 rvec * gmx_restrict ff,
379 t_forcerec * gmx_restrict fr,
380 t_mdatoms * gmx_restrict mdatoms,
381 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
382 t_nrnb * gmx_restrict nrnb)
384 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
385 * just 0 for non-waters.
386 * Suffixes A,B,C,D refer to j loop unrolling done with AVX_128, e.g. for the four different
387 * jnr indices corresponding to data put in the four positions in the SIMD register.
389 int i_shift_offset,i_coord_offset,outeriter,inneriter;
390 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
391 int jnrA,jnrB,jnrC,jnrD;
392 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
393 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
394 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
396 real *shiftvec,*fshift,*x,*f;
397 real *fjptrA,*fjptrB,*fjptrC,*fjptrD;
399 __m128 fscal,rcutoff,rcutoff2,jidxall;
401 __m128 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
402 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
403 __m128 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
404 __m128 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
405 __m128 velec,felec,velecsum,facel,crf,krf,krf2;
408 __m128i ifour = _mm_set1_epi32(4);
409 __m128 rt,vfeps,twovfeps,vftabscale,Y,F,G,H,Fp,VV,FF;
411 __m128 dummy_mask,cutoff_mask;
412 __m128 signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
413 __m128 one = _mm_set1_ps(1.0);
414 __m128 two = _mm_set1_ps(2.0);
420 jindex = nlist->jindex;
422 shiftidx = nlist->shift;
424 shiftvec = fr->shift_vec[0];
425 fshift = fr->fshift[0];
426 facel = _mm_set1_ps(fr->epsfac);
427 charge = mdatoms->chargeA;
429 vftab = kernel_data->table_elec->data;
430 vftabscale = _mm_set1_ps(kernel_data->table_elec->scale);
432 /* Avoid stupid compiler warnings */
433 jnrA = jnrB = jnrC = jnrD = 0;
442 for(iidx=0;iidx<4*DIM;iidx++)
447 /* Start outer loop over neighborlists */
448 for(iidx=0; iidx<nri; iidx++)
450 /* Load shift vector for this list */
451 i_shift_offset = DIM*shiftidx[iidx];
453 /* Load limits for loop over neighbors */
454 j_index_start = jindex[iidx];
455 j_index_end = jindex[iidx+1];
457 /* Get outer coordinate index */
459 i_coord_offset = DIM*inr;
461 /* Load i particle coords and add shift vector */
462 gmx_mm_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
464 fix0 = _mm_setzero_ps();
465 fiy0 = _mm_setzero_ps();
466 fiz0 = _mm_setzero_ps();
468 /* Load parameters for i particles */
469 iq0 = _mm_mul_ps(facel,_mm_load1_ps(charge+inr+0));
471 /* Start inner kernel loop */
472 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
475 /* Get j neighbor index, and coordinate index */
480 j_coord_offsetA = DIM*jnrA;
481 j_coord_offsetB = DIM*jnrB;
482 j_coord_offsetC = DIM*jnrC;
483 j_coord_offsetD = DIM*jnrD;
485 /* load j atom coordinates */
486 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
487 x+j_coord_offsetC,x+j_coord_offsetD,
490 /* Calculate displacement vector */
491 dx00 = _mm_sub_ps(ix0,jx0);
492 dy00 = _mm_sub_ps(iy0,jy0);
493 dz00 = _mm_sub_ps(iz0,jz0);
495 /* Calculate squared distance and things based on it */
496 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
498 rinv00 = gmx_mm_invsqrt_ps(rsq00);
500 /* Load parameters for j particles */
501 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
502 charge+jnrC+0,charge+jnrD+0);
504 /**************************
505 * CALCULATE INTERACTIONS *
506 **************************/
508 r00 = _mm_mul_ps(rsq00,rinv00);
510 /* Compute parameters for interactions between i and j atoms */
511 qq00 = _mm_mul_ps(iq0,jq0);
513 /* Calculate table index by multiplying r with table scale and truncate to integer */
514 rt = _mm_mul_ps(r00,vftabscale);
515 vfitab = _mm_cvttps_epi32(rt);
517 vfeps = _mm_frcz_ps(rt);
519 vfeps = _mm_sub_ps(rt,_mm_round_ps(rt, _MM_FROUND_FLOOR));
521 twovfeps = _mm_add_ps(vfeps,vfeps);
522 vfitab = _mm_slli_epi32(vfitab,2);
524 /* CUBIC SPLINE TABLE ELECTROSTATICS */
525 Y = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,0) );
526 F = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,1) );
527 G = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,2) );
528 H = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,3) );
529 _MM_TRANSPOSE4_PS(Y,F,G,H);
530 Fp = _mm_macc_ps(vfeps,_mm_macc_ps(H,vfeps,G),F);
531 FF = _mm_macc_ps(vfeps,_mm_macc_ps(twovfeps,H,G),Fp);
532 felec = _mm_xor_ps(signbit,_mm_mul_ps(_mm_mul_ps(qq00,FF),_mm_mul_ps(vftabscale,rinv00)));
536 /* Update vectorial force */
537 fix0 = _mm_macc_ps(dx00,fscal,fix0);
538 fiy0 = _mm_macc_ps(dy00,fscal,fiy0);
539 fiz0 = _mm_macc_ps(dz00,fscal,fiz0);
541 fjptrA = f+j_coord_offsetA;
542 fjptrB = f+j_coord_offsetB;
543 fjptrC = f+j_coord_offsetC;
544 fjptrD = f+j_coord_offsetD;
545 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,
546 _mm_mul_ps(dx00,fscal),
547 _mm_mul_ps(dy00,fscal),
548 _mm_mul_ps(dz00,fscal));
550 /* Inner loop uses 42 flops */
556 /* Get j neighbor index, and coordinate index */
557 jnrlistA = jjnr[jidx];
558 jnrlistB = jjnr[jidx+1];
559 jnrlistC = jjnr[jidx+2];
560 jnrlistD = jjnr[jidx+3];
561 /* Sign of each element will be negative for non-real atoms.
562 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
563 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
565 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
566 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
567 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
568 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
569 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
570 j_coord_offsetA = DIM*jnrA;
571 j_coord_offsetB = DIM*jnrB;
572 j_coord_offsetC = DIM*jnrC;
573 j_coord_offsetD = DIM*jnrD;
575 /* load j atom coordinates */
576 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
577 x+j_coord_offsetC,x+j_coord_offsetD,
580 /* Calculate displacement vector */
581 dx00 = _mm_sub_ps(ix0,jx0);
582 dy00 = _mm_sub_ps(iy0,jy0);
583 dz00 = _mm_sub_ps(iz0,jz0);
585 /* Calculate squared distance and things based on it */
586 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
588 rinv00 = gmx_mm_invsqrt_ps(rsq00);
590 /* Load parameters for j particles */
591 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
592 charge+jnrC+0,charge+jnrD+0);
594 /**************************
595 * CALCULATE INTERACTIONS *
596 **************************/
598 r00 = _mm_mul_ps(rsq00,rinv00);
599 r00 = _mm_andnot_ps(dummy_mask,r00);
601 /* Compute parameters for interactions between i and j atoms */
602 qq00 = _mm_mul_ps(iq0,jq0);
604 /* Calculate table index by multiplying r with table scale and truncate to integer */
605 rt = _mm_mul_ps(r00,vftabscale);
606 vfitab = _mm_cvttps_epi32(rt);
608 vfeps = _mm_frcz_ps(rt);
610 vfeps = _mm_sub_ps(rt,_mm_round_ps(rt, _MM_FROUND_FLOOR));
612 twovfeps = _mm_add_ps(vfeps,vfeps);
613 vfitab = _mm_slli_epi32(vfitab,2);
615 /* CUBIC SPLINE TABLE ELECTROSTATICS */
616 Y = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,0) );
617 F = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,1) );
618 G = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,2) );
619 H = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,3) );
620 _MM_TRANSPOSE4_PS(Y,F,G,H);
621 Fp = _mm_macc_ps(vfeps,_mm_macc_ps(H,vfeps,G),F);
622 FF = _mm_macc_ps(vfeps,_mm_macc_ps(twovfeps,H,G),Fp);
623 felec = _mm_xor_ps(signbit,_mm_mul_ps(_mm_mul_ps(qq00,FF),_mm_mul_ps(vftabscale,rinv00)));
627 fscal = _mm_andnot_ps(dummy_mask,fscal);
629 /* Update vectorial force */
630 fix0 = _mm_macc_ps(dx00,fscal,fix0);
631 fiy0 = _mm_macc_ps(dy00,fscal,fiy0);
632 fiz0 = _mm_macc_ps(dz00,fscal,fiz0);
634 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
635 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
636 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
637 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
638 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,
639 _mm_mul_ps(dx00,fscal),
640 _mm_mul_ps(dy00,fscal),
641 _mm_mul_ps(dz00,fscal));
643 /* Inner loop uses 43 flops */
646 /* End of innermost loop */
648 gmx_mm_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
649 f+i_coord_offset,fshift+i_shift_offset);
651 /* Increment number of inner iterations */
652 inneriter += j_index_end - j_index_start;
654 /* Outer loop uses 7 flops */
657 /* Increment number of outer iterations */
660 /* Update outer/inner flops */
662 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_F,outeriter*7 + inneriter*43);