2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014, 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 "gromacs/legacyheaders/types/simple.h"
46 #include "gromacs/math/vec.h"
47 #include "gromacs/legacyheaders/nrnb.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_ElecGB_VdwCSTab_GeomP1P1_VF_avx_128_fma_single
54 * Electrostatics interaction: GeneralizedBorn
55 * VdW interaction: CubicSplineTable
56 * Geometry: Particle-Particle
57 * Calculate force/pot: PotentialAndForce
60 nb_kernel_ElecGB_VdwCSTab_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 __m128 vgb,fgb,vgbsum,dvdasum,gbscale,gbtabscale,isaprod,gbqqfactor,gbinvepsdiff,gbeps,twogbeps,dvdatmp;
94 __m128 minushalf = _mm_set1_ps(-0.5);
95 real *invsqrta,*dvda,*gbtab;
97 __m128 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
100 __m128 one_sixth = _mm_set1_ps(1.0/6.0);
101 __m128 one_twelfth = _mm_set1_ps(1.0/12.0);
103 __m128i ifour = _mm_set1_epi32(4);
104 __m128 rt,vfeps,twovfeps,vftabscale,Y,F,G,H,Fp,VV,FF;
106 __m128 dummy_mask,cutoff_mask;
107 __m128 signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
108 __m128 one = _mm_set1_ps(1.0);
109 __m128 two = _mm_set1_ps(2.0);
115 jindex = nlist->jindex;
117 shiftidx = nlist->shift;
119 shiftvec = fr->shift_vec[0];
120 fshift = fr->fshift[0];
121 facel = _mm_set1_ps(fr->epsfac);
122 charge = mdatoms->chargeA;
123 nvdwtype = fr->ntype;
125 vdwtype = mdatoms->typeA;
127 vftab = kernel_data->table_vdw->data;
128 vftabscale = _mm_set1_ps(kernel_data->table_vdw->scale);
130 invsqrta = fr->invsqrta;
132 gbtabscale = _mm_set1_ps(fr->gbtab.scale);
133 gbtab = fr->gbtab.data;
134 gbinvepsdiff = _mm_set1_ps((1.0/fr->epsilon_r) - (1.0/fr->gb_epsilon_solvent));
136 /* Avoid stupid compiler warnings */
137 jnrA = jnrB = jnrC = jnrD = 0;
146 for(iidx=0;iidx<4*DIM;iidx++)
151 /* Start outer loop over neighborlists */
152 for(iidx=0; iidx<nri; iidx++)
154 /* Load shift vector for this list */
155 i_shift_offset = DIM*shiftidx[iidx];
157 /* Load limits for loop over neighbors */
158 j_index_start = jindex[iidx];
159 j_index_end = jindex[iidx+1];
161 /* Get outer coordinate index */
163 i_coord_offset = DIM*inr;
165 /* Load i particle coords and add shift vector */
166 gmx_mm_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
168 fix0 = _mm_setzero_ps();
169 fiy0 = _mm_setzero_ps();
170 fiz0 = _mm_setzero_ps();
172 /* Load parameters for i particles */
173 iq0 = _mm_mul_ps(facel,_mm_load1_ps(charge+inr+0));
174 isai0 = _mm_load1_ps(invsqrta+inr+0);
175 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
177 /* Reset potential sums */
178 velecsum = _mm_setzero_ps();
179 vgbsum = _mm_setzero_ps();
180 vvdwsum = _mm_setzero_ps();
181 dvdasum = _mm_setzero_ps();
183 /* Start inner kernel loop */
184 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
187 /* Get j neighbor index, and coordinate index */
192 j_coord_offsetA = DIM*jnrA;
193 j_coord_offsetB = DIM*jnrB;
194 j_coord_offsetC = DIM*jnrC;
195 j_coord_offsetD = DIM*jnrD;
197 /* load j atom coordinates */
198 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
199 x+j_coord_offsetC,x+j_coord_offsetD,
202 /* Calculate displacement vector */
203 dx00 = _mm_sub_ps(ix0,jx0);
204 dy00 = _mm_sub_ps(iy0,jy0);
205 dz00 = _mm_sub_ps(iz0,jz0);
207 /* Calculate squared distance and things based on it */
208 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
210 rinv00 = gmx_mm_invsqrt_ps(rsq00);
212 /* Load parameters for j particles */
213 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
214 charge+jnrC+0,charge+jnrD+0);
215 isaj0 = gmx_mm_load_4real_swizzle_ps(invsqrta+jnrA+0,invsqrta+jnrB+0,
216 invsqrta+jnrC+0,invsqrta+jnrD+0);
217 vdwjidx0A = 2*vdwtype[jnrA+0];
218 vdwjidx0B = 2*vdwtype[jnrB+0];
219 vdwjidx0C = 2*vdwtype[jnrC+0];
220 vdwjidx0D = 2*vdwtype[jnrD+0];
222 /**************************
223 * CALCULATE INTERACTIONS *
224 **************************/
226 r00 = _mm_mul_ps(rsq00,rinv00);
228 /* Compute parameters for interactions between i and j atoms */
229 qq00 = _mm_mul_ps(iq0,jq0);
230 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
231 vdwparam+vdwioffset0+vdwjidx0B,
232 vdwparam+vdwioffset0+vdwjidx0C,
233 vdwparam+vdwioffset0+vdwjidx0D,
236 /* Calculate table index by multiplying r with table scale and truncate to integer */
237 rt = _mm_mul_ps(r00,vftabscale);
238 vfitab = _mm_cvttps_epi32(rt);
240 vfeps = _mm_frcz_ps(rt);
242 vfeps = _mm_sub_ps(rt,_mm_round_ps(rt, _MM_FROUND_FLOOR));
244 twovfeps = _mm_add_ps(vfeps,vfeps);
245 vfitab = _mm_slli_epi32(vfitab,3);
247 /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
248 isaprod = _mm_mul_ps(isai0,isaj0);
249 gbqqfactor = _mm_xor_ps(signbit,_mm_mul_ps(qq00,_mm_mul_ps(isaprod,gbinvepsdiff)));
250 gbscale = _mm_mul_ps(isaprod,gbtabscale);
252 /* Calculate generalized born table index - this is a separate table from the normal one,
253 * but we use the same procedure by multiplying r with scale and truncating to integer.
255 rt = _mm_mul_ps(r00,gbscale);
256 gbitab = _mm_cvttps_epi32(rt);
258 gbeps = _mm_frcz_ps(rt);
260 gbeps = _mm_sub_ps(rt,_mm_round_ps(rt, _MM_FROUND_FLOOR));
262 gbitab = _mm_slli_epi32(gbitab,2);
264 Y = _mm_load_ps( gbtab + _mm_extract_epi32(gbitab,0) );
265 F = _mm_load_ps( gbtab + _mm_extract_epi32(gbitab,1) );
266 G = _mm_load_ps( gbtab + _mm_extract_epi32(gbitab,2) );
267 H = _mm_load_ps( gbtab + _mm_extract_epi32(gbitab,3) );
268 _MM_TRANSPOSE4_PS(Y,F,G,H);
269 Fp = _mm_macc_ps(gbeps,_mm_macc_ps(gbeps,H,G),F);
270 VV = _mm_macc_ps(gbeps,Fp,Y);
271 vgb = _mm_mul_ps(gbqqfactor,VV);
273 twogbeps = _mm_add_ps(gbeps,gbeps);
274 FF = _mm_macc_ps(_mm_macc_ps(twogbeps,H,G),gbeps,Fp);
275 fgb = _mm_mul_ps(gbqqfactor,_mm_mul_ps(FF,gbscale));
276 dvdatmp = _mm_mul_ps(minushalf,_mm_macc_ps(fgb,r00,vgb));
277 dvdasum = _mm_add_ps(dvdasum,dvdatmp);
282 gmx_mm_increment_4real_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,_mm_mul_ps(dvdatmp,_mm_mul_ps(isaj0,isaj0)));
283 velec = _mm_mul_ps(qq00,rinv00);
284 felec = _mm_mul_ps(_mm_msub_ps(velec,rinv00,fgb),rinv00);
286 /* CUBIC SPLINE TABLE DISPERSION */
287 Y = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,0) );
288 F = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,1) );
289 G = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,2) );
290 H = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,3) );
291 _MM_TRANSPOSE4_PS(Y,F,G,H);
292 Fp = _mm_macc_ps(vfeps,_mm_macc_ps(H,vfeps,G),F);
293 VV = _mm_macc_ps(vfeps,Fp,Y);
294 vvdw6 = _mm_mul_ps(c6_00,VV);
295 FF = _mm_macc_ps(vfeps,_mm_macc_ps(twovfeps,H,G),Fp);
296 fvdw6 = _mm_mul_ps(c6_00,FF);
298 /* CUBIC SPLINE TABLE REPULSION */
299 vfitab = _mm_add_epi32(vfitab,ifour);
300 Y = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,0) );
301 F = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,1) );
302 G = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,2) );
303 H = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,3) );
304 _MM_TRANSPOSE4_PS(Y,F,G,H);
305 Fp = _mm_macc_ps(vfeps,_mm_macc_ps(H,vfeps,G),F);
306 VV = _mm_macc_ps(vfeps,Fp,Y);
307 vvdw12 = _mm_mul_ps(c12_00,VV);
308 FF = _mm_macc_ps(vfeps,_mm_macc_ps(twovfeps,H,G),Fp);
309 fvdw12 = _mm_mul_ps(c12_00,FF);
310 vvdw = _mm_add_ps(vvdw12,vvdw6);
311 fvdw = _mm_xor_ps(signbit,_mm_mul_ps(_mm_add_ps(fvdw6,fvdw12),_mm_mul_ps(vftabscale,rinv00)));
313 /* Update potential sum for this i atom from the interaction with this j atom. */
314 velecsum = _mm_add_ps(velecsum,velec);
315 vgbsum = _mm_add_ps(vgbsum,vgb);
316 vvdwsum = _mm_add_ps(vvdwsum,vvdw);
318 fscal = _mm_add_ps(felec,fvdw);
320 /* Update vectorial force */
321 fix0 = _mm_macc_ps(dx00,fscal,fix0);
322 fiy0 = _mm_macc_ps(dy00,fscal,fiy0);
323 fiz0 = _mm_macc_ps(dz00,fscal,fiz0);
325 fjptrA = f+j_coord_offsetA;
326 fjptrB = f+j_coord_offsetB;
327 fjptrC = f+j_coord_offsetC;
328 fjptrD = f+j_coord_offsetD;
329 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,
330 _mm_mul_ps(dx00,fscal),
331 _mm_mul_ps(dy00,fscal),
332 _mm_mul_ps(dz00,fscal));
334 /* Inner loop uses 95 flops */
340 /* Get j neighbor index, and coordinate index */
341 jnrlistA = jjnr[jidx];
342 jnrlistB = jjnr[jidx+1];
343 jnrlistC = jjnr[jidx+2];
344 jnrlistD = jjnr[jidx+3];
345 /* Sign of each element will be negative for non-real atoms.
346 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
347 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
349 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
350 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
351 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
352 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
353 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
354 j_coord_offsetA = DIM*jnrA;
355 j_coord_offsetB = DIM*jnrB;
356 j_coord_offsetC = DIM*jnrC;
357 j_coord_offsetD = DIM*jnrD;
359 /* load j atom coordinates */
360 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
361 x+j_coord_offsetC,x+j_coord_offsetD,
364 /* Calculate displacement vector */
365 dx00 = _mm_sub_ps(ix0,jx0);
366 dy00 = _mm_sub_ps(iy0,jy0);
367 dz00 = _mm_sub_ps(iz0,jz0);
369 /* Calculate squared distance and things based on it */
370 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
372 rinv00 = gmx_mm_invsqrt_ps(rsq00);
374 /* Load parameters for j particles */
375 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
376 charge+jnrC+0,charge+jnrD+0);
377 isaj0 = gmx_mm_load_4real_swizzle_ps(invsqrta+jnrA+0,invsqrta+jnrB+0,
378 invsqrta+jnrC+0,invsqrta+jnrD+0);
379 vdwjidx0A = 2*vdwtype[jnrA+0];
380 vdwjidx0B = 2*vdwtype[jnrB+0];
381 vdwjidx0C = 2*vdwtype[jnrC+0];
382 vdwjidx0D = 2*vdwtype[jnrD+0];
384 /**************************
385 * CALCULATE INTERACTIONS *
386 **************************/
388 r00 = _mm_mul_ps(rsq00,rinv00);
389 r00 = _mm_andnot_ps(dummy_mask,r00);
391 /* Compute parameters for interactions between i and j atoms */
392 qq00 = _mm_mul_ps(iq0,jq0);
393 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
394 vdwparam+vdwioffset0+vdwjidx0B,
395 vdwparam+vdwioffset0+vdwjidx0C,
396 vdwparam+vdwioffset0+vdwjidx0D,
399 /* Calculate table index by multiplying r with table scale and truncate to integer */
400 rt = _mm_mul_ps(r00,vftabscale);
401 vfitab = _mm_cvttps_epi32(rt);
403 vfeps = _mm_frcz_ps(rt);
405 vfeps = _mm_sub_ps(rt,_mm_round_ps(rt, _MM_FROUND_FLOOR));
407 twovfeps = _mm_add_ps(vfeps,vfeps);
408 vfitab = _mm_slli_epi32(vfitab,3);
410 /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
411 isaprod = _mm_mul_ps(isai0,isaj0);
412 gbqqfactor = _mm_xor_ps(signbit,_mm_mul_ps(qq00,_mm_mul_ps(isaprod,gbinvepsdiff)));
413 gbscale = _mm_mul_ps(isaprod,gbtabscale);
415 /* Calculate generalized born table index - this is a separate table from the normal one,
416 * but we use the same procedure by multiplying r with scale and truncating to integer.
418 rt = _mm_mul_ps(r00,gbscale);
419 gbitab = _mm_cvttps_epi32(rt);
421 gbeps = _mm_frcz_ps(rt);
423 gbeps = _mm_sub_ps(rt,_mm_round_ps(rt, _MM_FROUND_FLOOR));
425 gbitab = _mm_slli_epi32(gbitab,2);
427 Y = _mm_load_ps( gbtab + _mm_extract_epi32(gbitab,0) );
428 F = _mm_load_ps( gbtab + _mm_extract_epi32(gbitab,1) );
429 G = _mm_load_ps( gbtab + _mm_extract_epi32(gbitab,2) );
430 H = _mm_load_ps( gbtab + _mm_extract_epi32(gbitab,3) );
431 _MM_TRANSPOSE4_PS(Y,F,G,H);
432 Fp = _mm_macc_ps(gbeps,_mm_macc_ps(gbeps,H,G),F);
433 VV = _mm_macc_ps(gbeps,Fp,Y);
434 vgb = _mm_mul_ps(gbqqfactor,VV);
436 twogbeps = _mm_add_ps(gbeps,gbeps);
437 FF = _mm_macc_ps(_mm_macc_ps(twogbeps,H,G),gbeps,Fp);
438 fgb = _mm_mul_ps(gbqqfactor,_mm_mul_ps(FF,gbscale));
439 dvdatmp = _mm_mul_ps(minushalf,_mm_macc_ps(fgb,r00,vgb));
440 dvdatmp = _mm_andnot_ps(dummy_mask,dvdatmp);
441 dvdasum = _mm_add_ps(dvdasum,dvdatmp);
442 /* 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. */
443 fjptrA = (jnrlistA>=0) ? dvda+jnrA : scratch;
444 fjptrB = (jnrlistB>=0) ? dvda+jnrB : scratch;
445 fjptrC = (jnrlistC>=0) ? dvda+jnrC : scratch;
446 fjptrD = (jnrlistD>=0) ? dvda+jnrD : scratch;
447 gmx_mm_increment_4real_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,_mm_mul_ps(dvdatmp,_mm_mul_ps(isaj0,isaj0)));
448 velec = _mm_mul_ps(qq00,rinv00);
449 felec = _mm_mul_ps(_mm_msub_ps(velec,rinv00,fgb),rinv00);
451 /* CUBIC SPLINE TABLE DISPERSION */
452 Y = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,0) );
453 F = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,1) );
454 G = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,2) );
455 H = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,3) );
456 _MM_TRANSPOSE4_PS(Y,F,G,H);
457 Fp = _mm_macc_ps(vfeps,_mm_macc_ps(H,vfeps,G),F);
458 VV = _mm_macc_ps(vfeps,Fp,Y);
459 vvdw6 = _mm_mul_ps(c6_00,VV);
460 FF = _mm_macc_ps(vfeps,_mm_macc_ps(twovfeps,H,G),Fp);
461 fvdw6 = _mm_mul_ps(c6_00,FF);
463 /* CUBIC SPLINE TABLE REPULSION */
464 vfitab = _mm_add_epi32(vfitab,ifour);
465 Y = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,0) );
466 F = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,1) );
467 G = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,2) );
468 H = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,3) );
469 _MM_TRANSPOSE4_PS(Y,F,G,H);
470 Fp = _mm_macc_ps(vfeps,_mm_macc_ps(H,vfeps,G),F);
471 VV = _mm_macc_ps(vfeps,Fp,Y);
472 vvdw12 = _mm_mul_ps(c12_00,VV);
473 FF = _mm_macc_ps(vfeps,_mm_macc_ps(twovfeps,H,G),Fp);
474 fvdw12 = _mm_mul_ps(c12_00,FF);
475 vvdw = _mm_add_ps(vvdw12,vvdw6);
476 fvdw = _mm_xor_ps(signbit,_mm_mul_ps(_mm_add_ps(fvdw6,fvdw12),_mm_mul_ps(vftabscale,rinv00)));
478 /* Update potential sum for this i atom from the interaction with this j atom. */
479 velec = _mm_andnot_ps(dummy_mask,velec);
480 velecsum = _mm_add_ps(velecsum,velec);
481 vgb = _mm_andnot_ps(dummy_mask,vgb);
482 vgbsum = _mm_add_ps(vgbsum,vgb);
483 vvdw = _mm_andnot_ps(dummy_mask,vvdw);
484 vvdwsum = _mm_add_ps(vvdwsum,vvdw);
486 fscal = _mm_add_ps(felec,fvdw);
488 fscal = _mm_andnot_ps(dummy_mask,fscal);
490 /* Update vectorial force */
491 fix0 = _mm_macc_ps(dx00,fscal,fix0);
492 fiy0 = _mm_macc_ps(dy00,fscal,fiy0);
493 fiz0 = _mm_macc_ps(dz00,fscal,fiz0);
495 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
496 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
497 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
498 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
499 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,
500 _mm_mul_ps(dx00,fscal),
501 _mm_mul_ps(dy00,fscal),
502 _mm_mul_ps(dz00,fscal));
504 /* Inner loop uses 96 flops */
507 /* End of innermost loop */
509 gmx_mm_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
510 f+i_coord_offset,fshift+i_shift_offset);
513 /* Update potential energies */
514 gmx_mm_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
515 gmx_mm_update_1pot_ps(vgbsum,kernel_data->energygrp_polarization+ggid);
516 gmx_mm_update_1pot_ps(vvdwsum,kernel_data->energygrp_vdw+ggid);
517 dvdasum = _mm_mul_ps(dvdasum, _mm_mul_ps(isai0,isai0));
518 gmx_mm_update_1pot_ps(dvdasum,dvda+inr);
520 /* Increment number of inner iterations */
521 inneriter += j_index_end - j_index_start;
523 /* Outer loop uses 10 flops */
526 /* Increment number of outer iterations */
529 /* Update outer/inner flops */
531 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_VF,outeriter*10 + inneriter*96);
534 * Gromacs nonbonded kernel: nb_kernel_ElecGB_VdwCSTab_GeomP1P1_F_avx_128_fma_single
535 * Electrostatics interaction: GeneralizedBorn
536 * VdW interaction: CubicSplineTable
537 * Geometry: Particle-Particle
538 * Calculate force/pot: Force
541 nb_kernel_ElecGB_VdwCSTab_GeomP1P1_F_avx_128_fma_single
542 (t_nblist * gmx_restrict nlist,
543 rvec * gmx_restrict xx,
544 rvec * gmx_restrict ff,
545 t_forcerec * gmx_restrict fr,
546 t_mdatoms * gmx_restrict mdatoms,
547 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
548 t_nrnb * gmx_restrict nrnb)
550 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
551 * just 0 for non-waters.
552 * Suffixes A,B,C,D refer to j loop unrolling done with AVX_128, e.g. for the four different
553 * jnr indices corresponding to data put in the four positions in the SIMD register.
555 int i_shift_offset,i_coord_offset,outeriter,inneriter;
556 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
557 int jnrA,jnrB,jnrC,jnrD;
558 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
559 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
560 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
562 real *shiftvec,*fshift,*x,*f;
563 real *fjptrA,*fjptrB,*fjptrC,*fjptrD;
565 __m128 fscal,rcutoff,rcutoff2,jidxall;
567 __m128 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
568 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
569 __m128 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
570 __m128 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
571 __m128 velec,felec,velecsum,facel,crf,krf,krf2;
574 __m128 vgb,fgb,vgbsum,dvdasum,gbscale,gbtabscale,isaprod,gbqqfactor,gbinvepsdiff,gbeps,twogbeps,dvdatmp;
575 __m128 minushalf = _mm_set1_ps(-0.5);
576 real *invsqrta,*dvda,*gbtab;
578 __m128 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
581 __m128 one_sixth = _mm_set1_ps(1.0/6.0);
582 __m128 one_twelfth = _mm_set1_ps(1.0/12.0);
584 __m128i ifour = _mm_set1_epi32(4);
585 __m128 rt,vfeps,twovfeps,vftabscale,Y,F,G,H,Fp,VV,FF;
587 __m128 dummy_mask,cutoff_mask;
588 __m128 signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
589 __m128 one = _mm_set1_ps(1.0);
590 __m128 two = _mm_set1_ps(2.0);
596 jindex = nlist->jindex;
598 shiftidx = nlist->shift;
600 shiftvec = fr->shift_vec[0];
601 fshift = fr->fshift[0];
602 facel = _mm_set1_ps(fr->epsfac);
603 charge = mdatoms->chargeA;
604 nvdwtype = fr->ntype;
606 vdwtype = mdatoms->typeA;
608 vftab = kernel_data->table_vdw->data;
609 vftabscale = _mm_set1_ps(kernel_data->table_vdw->scale);
611 invsqrta = fr->invsqrta;
613 gbtabscale = _mm_set1_ps(fr->gbtab.scale);
614 gbtab = fr->gbtab.data;
615 gbinvepsdiff = _mm_set1_ps((1.0/fr->epsilon_r) - (1.0/fr->gb_epsilon_solvent));
617 /* Avoid stupid compiler warnings */
618 jnrA = jnrB = jnrC = jnrD = 0;
627 for(iidx=0;iidx<4*DIM;iidx++)
632 /* Start outer loop over neighborlists */
633 for(iidx=0; iidx<nri; iidx++)
635 /* Load shift vector for this list */
636 i_shift_offset = DIM*shiftidx[iidx];
638 /* Load limits for loop over neighbors */
639 j_index_start = jindex[iidx];
640 j_index_end = jindex[iidx+1];
642 /* Get outer coordinate index */
644 i_coord_offset = DIM*inr;
646 /* Load i particle coords and add shift vector */
647 gmx_mm_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
649 fix0 = _mm_setzero_ps();
650 fiy0 = _mm_setzero_ps();
651 fiz0 = _mm_setzero_ps();
653 /* Load parameters for i particles */
654 iq0 = _mm_mul_ps(facel,_mm_load1_ps(charge+inr+0));
655 isai0 = _mm_load1_ps(invsqrta+inr+0);
656 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
658 dvdasum = _mm_setzero_ps();
660 /* Start inner kernel loop */
661 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
664 /* Get j neighbor index, and coordinate index */
669 j_coord_offsetA = DIM*jnrA;
670 j_coord_offsetB = DIM*jnrB;
671 j_coord_offsetC = DIM*jnrC;
672 j_coord_offsetD = DIM*jnrD;
674 /* load j atom coordinates */
675 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
676 x+j_coord_offsetC,x+j_coord_offsetD,
679 /* Calculate displacement vector */
680 dx00 = _mm_sub_ps(ix0,jx0);
681 dy00 = _mm_sub_ps(iy0,jy0);
682 dz00 = _mm_sub_ps(iz0,jz0);
684 /* Calculate squared distance and things based on it */
685 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
687 rinv00 = gmx_mm_invsqrt_ps(rsq00);
689 /* Load parameters for j particles */
690 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
691 charge+jnrC+0,charge+jnrD+0);
692 isaj0 = gmx_mm_load_4real_swizzle_ps(invsqrta+jnrA+0,invsqrta+jnrB+0,
693 invsqrta+jnrC+0,invsqrta+jnrD+0);
694 vdwjidx0A = 2*vdwtype[jnrA+0];
695 vdwjidx0B = 2*vdwtype[jnrB+0];
696 vdwjidx0C = 2*vdwtype[jnrC+0];
697 vdwjidx0D = 2*vdwtype[jnrD+0];
699 /**************************
700 * CALCULATE INTERACTIONS *
701 **************************/
703 r00 = _mm_mul_ps(rsq00,rinv00);
705 /* Compute parameters for interactions between i and j atoms */
706 qq00 = _mm_mul_ps(iq0,jq0);
707 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
708 vdwparam+vdwioffset0+vdwjidx0B,
709 vdwparam+vdwioffset0+vdwjidx0C,
710 vdwparam+vdwioffset0+vdwjidx0D,
713 /* Calculate table index by multiplying r with table scale and truncate to integer */
714 rt = _mm_mul_ps(r00,vftabscale);
715 vfitab = _mm_cvttps_epi32(rt);
717 vfeps = _mm_frcz_ps(rt);
719 vfeps = _mm_sub_ps(rt,_mm_round_ps(rt, _MM_FROUND_FLOOR));
721 twovfeps = _mm_add_ps(vfeps,vfeps);
722 vfitab = _mm_slli_epi32(vfitab,3);
724 /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
725 isaprod = _mm_mul_ps(isai0,isaj0);
726 gbqqfactor = _mm_xor_ps(signbit,_mm_mul_ps(qq00,_mm_mul_ps(isaprod,gbinvepsdiff)));
727 gbscale = _mm_mul_ps(isaprod,gbtabscale);
729 /* Calculate generalized born table index - this is a separate table from the normal one,
730 * but we use the same procedure by multiplying r with scale and truncating to integer.
732 rt = _mm_mul_ps(r00,gbscale);
733 gbitab = _mm_cvttps_epi32(rt);
735 gbeps = _mm_frcz_ps(rt);
737 gbeps = _mm_sub_ps(rt,_mm_round_ps(rt, _MM_FROUND_FLOOR));
739 gbitab = _mm_slli_epi32(gbitab,2);
741 Y = _mm_load_ps( gbtab + _mm_extract_epi32(gbitab,0) );
742 F = _mm_load_ps( gbtab + _mm_extract_epi32(gbitab,1) );
743 G = _mm_load_ps( gbtab + _mm_extract_epi32(gbitab,2) );
744 H = _mm_load_ps( gbtab + _mm_extract_epi32(gbitab,3) );
745 _MM_TRANSPOSE4_PS(Y,F,G,H);
746 Fp = _mm_macc_ps(gbeps,_mm_macc_ps(gbeps,H,G),F);
747 VV = _mm_macc_ps(gbeps,Fp,Y);
748 vgb = _mm_mul_ps(gbqqfactor,VV);
750 twogbeps = _mm_add_ps(gbeps,gbeps);
751 FF = _mm_macc_ps(_mm_macc_ps(twogbeps,H,G),gbeps,Fp);
752 fgb = _mm_mul_ps(gbqqfactor,_mm_mul_ps(FF,gbscale));
753 dvdatmp = _mm_mul_ps(minushalf,_mm_macc_ps(fgb,r00,vgb));
754 dvdasum = _mm_add_ps(dvdasum,dvdatmp);
759 gmx_mm_increment_4real_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,_mm_mul_ps(dvdatmp,_mm_mul_ps(isaj0,isaj0)));
760 velec = _mm_mul_ps(qq00,rinv00);
761 felec = _mm_mul_ps(_mm_msub_ps(velec,rinv00,fgb),rinv00);
763 /* CUBIC SPLINE TABLE DISPERSION */
764 Y = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,0) );
765 F = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,1) );
766 G = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,2) );
767 H = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,3) );
768 _MM_TRANSPOSE4_PS(Y,F,G,H);
769 Fp = _mm_macc_ps(vfeps,_mm_macc_ps(H,vfeps,G),F);
770 FF = _mm_macc_ps(vfeps,_mm_macc_ps(twovfeps,H,G),Fp);
771 fvdw6 = _mm_mul_ps(c6_00,FF);
773 /* CUBIC SPLINE TABLE REPULSION */
774 vfitab = _mm_add_epi32(vfitab,ifour);
775 Y = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,0) );
776 F = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,1) );
777 G = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,2) );
778 H = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,3) );
779 _MM_TRANSPOSE4_PS(Y,F,G,H);
780 Fp = _mm_macc_ps(vfeps,_mm_macc_ps(H,vfeps,G),F);
781 FF = _mm_macc_ps(vfeps,_mm_macc_ps(twovfeps,H,G),Fp);
782 fvdw12 = _mm_mul_ps(c12_00,FF);
783 fvdw = _mm_xor_ps(signbit,_mm_mul_ps(_mm_add_ps(fvdw6,fvdw12),_mm_mul_ps(vftabscale,rinv00)));
785 fscal = _mm_add_ps(felec,fvdw);
787 /* Update vectorial force */
788 fix0 = _mm_macc_ps(dx00,fscal,fix0);
789 fiy0 = _mm_macc_ps(dy00,fscal,fiy0);
790 fiz0 = _mm_macc_ps(dz00,fscal,fiz0);
792 fjptrA = f+j_coord_offsetA;
793 fjptrB = f+j_coord_offsetB;
794 fjptrC = f+j_coord_offsetC;
795 fjptrD = f+j_coord_offsetD;
796 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,
797 _mm_mul_ps(dx00,fscal),
798 _mm_mul_ps(dy00,fscal),
799 _mm_mul_ps(dz00,fscal));
801 /* Inner loop uses 85 flops */
807 /* Get j neighbor index, and coordinate index */
808 jnrlistA = jjnr[jidx];
809 jnrlistB = jjnr[jidx+1];
810 jnrlistC = jjnr[jidx+2];
811 jnrlistD = jjnr[jidx+3];
812 /* Sign of each element will be negative for non-real atoms.
813 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
814 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
816 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
817 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
818 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
819 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
820 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
821 j_coord_offsetA = DIM*jnrA;
822 j_coord_offsetB = DIM*jnrB;
823 j_coord_offsetC = DIM*jnrC;
824 j_coord_offsetD = DIM*jnrD;
826 /* load j atom coordinates */
827 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
828 x+j_coord_offsetC,x+j_coord_offsetD,
831 /* Calculate displacement vector */
832 dx00 = _mm_sub_ps(ix0,jx0);
833 dy00 = _mm_sub_ps(iy0,jy0);
834 dz00 = _mm_sub_ps(iz0,jz0);
836 /* Calculate squared distance and things based on it */
837 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
839 rinv00 = gmx_mm_invsqrt_ps(rsq00);
841 /* Load parameters for j particles */
842 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
843 charge+jnrC+0,charge+jnrD+0);
844 isaj0 = gmx_mm_load_4real_swizzle_ps(invsqrta+jnrA+0,invsqrta+jnrB+0,
845 invsqrta+jnrC+0,invsqrta+jnrD+0);
846 vdwjidx0A = 2*vdwtype[jnrA+0];
847 vdwjidx0B = 2*vdwtype[jnrB+0];
848 vdwjidx0C = 2*vdwtype[jnrC+0];
849 vdwjidx0D = 2*vdwtype[jnrD+0];
851 /**************************
852 * CALCULATE INTERACTIONS *
853 **************************/
855 r00 = _mm_mul_ps(rsq00,rinv00);
856 r00 = _mm_andnot_ps(dummy_mask,r00);
858 /* Compute parameters for interactions between i and j atoms */
859 qq00 = _mm_mul_ps(iq0,jq0);
860 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
861 vdwparam+vdwioffset0+vdwjidx0B,
862 vdwparam+vdwioffset0+vdwjidx0C,
863 vdwparam+vdwioffset0+vdwjidx0D,
866 /* Calculate table index by multiplying r with table scale and truncate to integer */
867 rt = _mm_mul_ps(r00,vftabscale);
868 vfitab = _mm_cvttps_epi32(rt);
870 vfeps = _mm_frcz_ps(rt);
872 vfeps = _mm_sub_ps(rt,_mm_round_ps(rt, _MM_FROUND_FLOOR));
874 twovfeps = _mm_add_ps(vfeps,vfeps);
875 vfitab = _mm_slli_epi32(vfitab,3);
877 /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
878 isaprod = _mm_mul_ps(isai0,isaj0);
879 gbqqfactor = _mm_xor_ps(signbit,_mm_mul_ps(qq00,_mm_mul_ps(isaprod,gbinvepsdiff)));
880 gbscale = _mm_mul_ps(isaprod,gbtabscale);
882 /* Calculate generalized born table index - this is a separate table from the normal one,
883 * but we use the same procedure by multiplying r with scale and truncating to integer.
885 rt = _mm_mul_ps(r00,gbscale);
886 gbitab = _mm_cvttps_epi32(rt);
888 gbeps = _mm_frcz_ps(rt);
890 gbeps = _mm_sub_ps(rt,_mm_round_ps(rt, _MM_FROUND_FLOOR));
892 gbitab = _mm_slli_epi32(gbitab,2);
894 Y = _mm_load_ps( gbtab + _mm_extract_epi32(gbitab,0) );
895 F = _mm_load_ps( gbtab + _mm_extract_epi32(gbitab,1) );
896 G = _mm_load_ps( gbtab + _mm_extract_epi32(gbitab,2) );
897 H = _mm_load_ps( gbtab + _mm_extract_epi32(gbitab,3) );
898 _MM_TRANSPOSE4_PS(Y,F,G,H);
899 Fp = _mm_macc_ps(gbeps,_mm_macc_ps(gbeps,H,G),F);
900 VV = _mm_macc_ps(gbeps,Fp,Y);
901 vgb = _mm_mul_ps(gbqqfactor,VV);
903 twogbeps = _mm_add_ps(gbeps,gbeps);
904 FF = _mm_macc_ps(_mm_macc_ps(twogbeps,H,G),gbeps,Fp);
905 fgb = _mm_mul_ps(gbqqfactor,_mm_mul_ps(FF,gbscale));
906 dvdatmp = _mm_mul_ps(minushalf,_mm_macc_ps(fgb,r00,vgb));
907 dvdatmp = _mm_andnot_ps(dummy_mask,dvdatmp);
908 dvdasum = _mm_add_ps(dvdasum,dvdatmp);
909 /* 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. */
910 fjptrA = (jnrlistA>=0) ? dvda+jnrA : scratch;
911 fjptrB = (jnrlistB>=0) ? dvda+jnrB : scratch;
912 fjptrC = (jnrlistC>=0) ? dvda+jnrC : scratch;
913 fjptrD = (jnrlistD>=0) ? dvda+jnrD : scratch;
914 gmx_mm_increment_4real_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,_mm_mul_ps(dvdatmp,_mm_mul_ps(isaj0,isaj0)));
915 velec = _mm_mul_ps(qq00,rinv00);
916 felec = _mm_mul_ps(_mm_msub_ps(velec,rinv00,fgb),rinv00);
918 /* CUBIC SPLINE TABLE DISPERSION */
919 Y = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,0) );
920 F = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,1) );
921 G = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,2) );
922 H = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,3) );
923 _MM_TRANSPOSE4_PS(Y,F,G,H);
924 Fp = _mm_macc_ps(vfeps,_mm_macc_ps(H,vfeps,G),F);
925 FF = _mm_macc_ps(vfeps,_mm_macc_ps(twovfeps,H,G),Fp);
926 fvdw6 = _mm_mul_ps(c6_00,FF);
928 /* CUBIC SPLINE TABLE REPULSION */
929 vfitab = _mm_add_epi32(vfitab,ifour);
930 Y = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,0) );
931 F = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,1) );
932 G = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,2) );
933 H = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,3) );
934 _MM_TRANSPOSE4_PS(Y,F,G,H);
935 Fp = _mm_macc_ps(vfeps,_mm_macc_ps(H,vfeps,G),F);
936 FF = _mm_macc_ps(vfeps,_mm_macc_ps(twovfeps,H,G),Fp);
937 fvdw12 = _mm_mul_ps(c12_00,FF);
938 fvdw = _mm_xor_ps(signbit,_mm_mul_ps(_mm_add_ps(fvdw6,fvdw12),_mm_mul_ps(vftabscale,rinv00)));
940 fscal = _mm_add_ps(felec,fvdw);
942 fscal = _mm_andnot_ps(dummy_mask,fscal);
944 /* Update vectorial force */
945 fix0 = _mm_macc_ps(dx00,fscal,fix0);
946 fiy0 = _mm_macc_ps(dy00,fscal,fiy0);
947 fiz0 = _mm_macc_ps(dz00,fscal,fiz0);
949 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
950 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
951 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
952 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
953 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,
954 _mm_mul_ps(dx00,fscal),
955 _mm_mul_ps(dy00,fscal),
956 _mm_mul_ps(dz00,fscal));
958 /* Inner loop uses 86 flops */
961 /* End of innermost loop */
963 gmx_mm_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
964 f+i_coord_offset,fshift+i_shift_offset);
966 dvdasum = _mm_mul_ps(dvdasum, _mm_mul_ps(isai0,isai0));
967 gmx_mm_update_1pot_ps(dvdasum,dvda+inr);
969 /* Increment number of inner iterations */
970 inneriter += j_index_end - j_index_start;
972 /* Outer loop uses 7 flops */
975 /* Increment number of outer iterations */
978 /* Update outer/inner flops */
980 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_F,outeriter*7 + inneriter*86);