2 * Note: this file was generated by the Gromacs avx_128_fma_double 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_avx_128_fma_double.h"
34 #include "kernelutil_x86_avx_128_fma_double.h"
37 * Gromacs nonbonded kernel: nb_kernel_ElecGB_VdwNone_GeomP1P1_VF_avx_128_fma_double
38 * Electrostatics interaction: GeneralizedBorn
39 * VdW interaction: None
40 * Geometry: Particle-Particle
41 * Calculate force/pot: PotentialAndForce
44 nb_kernel_ElecGB_VdwNone_GeomP1P1_VF_avx_128_fma_double
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 refer to j loop unrolling done with SSE double precision, e.g. for the two 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;
61 int j_coord_offsetA,j_coord_offsetB;
62 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
64 real *shiftvec,*fshift,*x,*f;
65 __m128d tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
67 __m128d ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
68 int vdwjidx0A,vdwjidx0B;
69 __m128d jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
70 __m128d dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
71 __m128d velec,felec,velecsum,facel,crf,krf,krf2;
74 __m128d vgb,fgb,vgbsum,dvdasum,gbscale,gbtabscale,isaprod,gbqqfactor,gbinvepsdiff,dvdaj,gbeps,twogbeps,dvdatmp;
75 __m128d minushalf = _mm_set1_pd(-0.5);
76 real *invsqrta,*dvda,*gbtab;
78 __m128i ifour = _mm_set1_epi32(4);
79 __m128d rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF,twovfeps;
81 __m128d dummy_mask,cutoff_mask;
82 __m128d signbit = gmx_mm_castsi128_pd( _mm_set_epi32(0x80000000,0x00000000,0x80000000,0x00000000) );
83 __m128d one = _mm_set1_pd(1.0);
84 __m128d two = _mm_set1_pd(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_pd(fr->epsfac);
97 charge = mdatoms->chargeA;
99 invsqrta = fr->invsqrta;
101 gbtabscale = _mm_set1_pd(fr->gbtab.scale);
102 gbtab = fr->gbtab.data;
103 gbinvepsdiff = _mm_set1_pd((1.0/fr->epsilon_r) - (1.0/fr->gb_epsilon_solvent));
105 /* Avoid stupid compiler warnings */
113 /* Start outer loop over neighborlists */
114 for(iidx=0; iidx<nri; iidx++)
116 /* Load shift vector for this list */
117 i_shift_offset = DIM*shiftidx[iidx];
119 /* Load limits for loop over neighbors */
120 j_index_start = jindex[iidx];
121 j_index_end = jindex[iidx+1];
123 /* Get outer coordinate index */
125 i_coord_offset = DIM*inr;
127 /* Load i particle coords and add shift vector */
128 gmx_mm_load_shift_and_1rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
130 fix0 = _mm_setzero_pd();
131 fiy0 = _mm_setzero_pd();
132 fiz0 = _mm_setzero_pd();
134 /* Load parameters for i particles */
135 iq0 = _mm_mul_pd(facel,_mm_load1_pd(charge+inr+0));
136 isai0 = _mm_load1_pd(invsqrta+inr+0);
138 /* Reset potential sums */
139 velecsum = _mm_setzero_pd();
140 vgbsum = _mm_setzero_pd();
141 dvdasum = _mm_setzero_pd();
143 /* Start inner kernel loop */
144 for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
147 /* Get j neighbor index, and coordinate index */
150 j_coord_offsetA = DIM*jnrA;
151 j_coord_offsetB = DIM*jnrB;
153 /* load j atom coordinates */
154 gmx_mm_load_1rvec_2ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
157 /* Calculate displacement vector */
158 dx00 = _mm_sub_pd(ix0,jx0);
159 dy00 = _mm_sub_pd(iy0,jy0);
160 dz00 = _mm_sub_pd(iz0,jz0);
162 /* Calculate squared distance and things based on it */
163 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
165 rinv00 = gmx_mm_invsqrt_pd(rsq00);
167 /* Load parameters for j particles */
168 jq0 = gmx_mm_load_2real_swizzle_pd(charge+jnrA+0,charge+jnrB+0);
169 isaj0 = gmx_mm_load_2real_swizzle_pd(invsqrta+jnrA+0,invsqrta+jnrB+0);
171 /**************************
172 * CALCULATE INTERACTIONS *
173 **************************/
175 r00 = _mm_mul_pd(rsq00,rinv00);
177 /* Compute parameters for interactions between i and j atoms */
178 qq00 = _mm_mul_pd(iq0,jq0);
180 /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
181 isaprod = _mm_mul_pd(isai0,isaj0);
182 gbqqfactor = _mm_xor_pd(signbit,_mm_mul_pd(qq00,_mm_mul_pd(isaprod,gbinvepsdiff)));
183 gbscale = _mm_mul_pd(isaprod,gbtabscale);
185 /* Calculate generalized born table index - this is a separate table from the normal one,
186 * but we use the same procedure by multiplying r with scale and truncating to integer.
188 rt = _mm_mul_pd(r00,gbscale);
189 gbitab = _mm_cvttpd_epi32(rt);
191 gbeps = _mm_frcz_pd(rt);
193 gbeps = _mm_sub_pd(rt,_mm_round_pd(rt, _MM_FROUND_FLOOR));
195 gbitab = _mm_slli_epi32(gbitab,2);
197 Y = _mm_load_pd( gbtab + _mm_extract_epi32(gbitab,0) );
198 F = _mm_load_pd( gbtab + _mm_extract_epi32(gbitab,1) );
199 GMX_MM_TRANSPOSE2_PD(Y,F);
200 G = _mm_load_pd( gbtab + _mm_extract_epi32(gbitab,0) +2);
201 H = _mm_load_pd( gbtab + _mm_extract_epi32(gbitab,1) +2);
202 GMX_MM_TRANSPOSE2_PD(G,H);
203 Fp = _mm_macc_pd(gbeps,_mm_macc_pd(gbeps,H,G),F);
204 VV = _mm_macc_pd(gbeps,Fp,Y);
205 vgb = _mm_mul_pd(gbqqfactor,VV);
207 twogbeps = _mm_add_pd(gbeps,gbeps);
208 FF = _mm_macc_pd(_mm_macc_pd(twogbeps,H,G),gbeps,Fp);
209 fgb = _mm_mul_pd(gbqqfactor,_mm_mul_pd(FF,gbscale));
210 dvdatmp = _mm_mul_pd(minushalf,_mm_macc_pd(fgb,r00,vgb));
211 dvdasum = _mm_add_pd(dvdasum,dvdatmp);
212 gmx_mm_increment_2real_swizzle_pd(dvda+jnrA,dvda+jnrB,_mm_mul_pd(dvdatmp,_mm_mul_pd(isaj0,isaj0)));
213 velec = _mm_mul_pd(qq00,rinv00);
214 felec = _mm_mul_pd(_mm_msub_pd(velec,rinv00,fgb),rinv00);
216 /* Update potential sum for this i atom from the interaction with this j atom. */
217 velecsum = _mm_add_pd(velecsum,velec);
218 vgbsum = _mm_add_pd(vgbsum,vgb);
222 /* Update vectorial force */
223 fix0 = _mm_macc_pd(dx00,fscal,fix0);
224 fiy0 = _mm_macc_pd(dy00,fscal,fiy0);
225 fiz0 = _mm_macc_pd(dz00,fscal,fiz0);
227 gmx_mm_decrement_1rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,
228 _mm_mul_pd(dx00,fscal),
229 _mm_mul_pd(dy00,fscal),
230 _mm_mul_pd(dz00,fscal));
232 /* Inner loop uses 61 flops */
239 j_coord_offsetA = DIM*jnrA;
241 /* load j atom coordinates */
242 gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
245 /* Calculate displacement vector */
246 dx00 = _mm_sub_pd(ix0,jx0);
247 dy00 = _mm_sub_pd(iy0,jy0);
248 dz00 = _mm_sub_pd(iz0,jz0);
250 /* Calculate squared distance and things based on it */
251 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
253 rinv00 = gmx_mm_invsqrt_pd(rsq00);
255 /* Load parameters for j particles */
256 jq0 = _mm_load_sd(charge+jnrA+0);
257 isaj0 = _mm_load_sd(invsqrta+jnrA+0);
259 /**************************
260 * CALCULATE INTERACTIONS *
261 **************************/
263 r00 = _mm_mul_pd(rsq00,rinv00);
265 /* Compute parameters for interactions between i and j atoms */
266 qq00 = _mm_mul_pd(iq0,jq0);
268 /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
269 isaprod = _mm_mul_pd(isai0,isaj0);
270 gbqqfactor = _mm_xor_pd(signbit,_mm_mul_pd(qq00,_mm_mul_pd(isaprod,gbinvepsdiff)));
271 gbscale = _mm_mul_pd(isaprod,gbtabscale);
273 /* Calculate generalized born table index - this is a separate table from the normal one,
274 * but we use the same procedure by multiplying r with scale and truncating to integer.
276 rt = _mm_mul_pd(r00,gbscale);
277 gbitab = _mm_cvttpd_epi32(rt);
279 gbeps = _mm_frcz_pd(rt);
281 gbeps = _mm_sub_pd(rt,_mm_round_pd(rt, _MM_FROUND_FLOOR));
283 gbitab = _mm_slli_epi32(gbitab,2);
285 Y = _mm_load_pd( gbtab + _mm_extract_epi32(gbitab,0) );
286 F = _mm_setzero_pd();
287 GMX_MM_TRANSPOSE2_PD(Y,F);
288 G = _mm_load_pd( gbtab + _mm_extract_epi32(gbitab,0) +2);
289 H = _mm_setzero_pd();
290 GMX_MM_TRANSPOSE2_PD(G,H);
291 Fp = _mm_macc_pd(gbeps,_mm_macc_pd(gbeps,H,G),F);
292 VV = _mm_macc_pd(gbeps,Fp,Y);
293 vgb = _mm_mul_pd(gbqqfactor,VV);
295 twogbeps = _mm_add_pd(gbeps,gbeps);
296 FF = _mm_macc_pd(_mm_macc_pd(twogbeps,H,G),gbeps,Fp);
297 fgb = _mm_mul_pd(gbqqfactor,_mm_mul_pd(FF,gbscale));
298 dvdatmp = _mm_mul_pd(minushalf,_mm_macc_pd(fgb,r00,vgb));
299 dvdasum = _mm_add_pd(dvdasum,dvdatmp);
300 gmx_mm_increment_1real_pd(dvda+jnrA,_mm_mul_pd(dvdatmp,_mm_mul_pd(isaj0,isaj0)));
301 velec = _mm_mul_pd(qq00,rinv00);
302 felec = _mm_mul_pd(_mm_msub_pd(velec,rinv00,fgb),rinv00);
304 /* Update potential sum for this i atom from the interaction with this j atom. */
305 velec = _mm_unpacklo_pd(velec,_mm_setzero_pd());
306 velecsum = _mm_add_pd(velecsum,velec);
307 vgb = _mm_unpacklo_pd(vgb,_mm_setzero_pd());
308 vgbsum = _mm_add_pd(vgbsum,vgb);
312 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
314 /* Update vectorial force */
315 fix0 = _mm_macc_pd(dx00,fscal,fix0);
316 fiy0 = _mm_macc_pd(dy00,fscal,fiy0);
317 fiz0 = _mm_macc_pd(dz00,fscal,fiz0);
319 gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,
320 _mm_mul_pd(dx00,fscal),
321 _mm_mul_pd(dy00,fscal),
322 _mm_mul_pd(dz00,fscal));
324 /* Inner loop uses 61 flops */
327 /* End of innermost loop */
329 gmx_mm_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
330 f+i_coord_offset,fshift+i_shift_offset);
333 /* Update potential energies */
334 gmx_mm_update_1pot_pd(velecsum,kernel_data->energygrp_elec+ggid);
335 gmx_mm_update_1pot_pd(vgbsum,kernel_data->energygrp_polarization+ggid);
336 dvdasum = _mm_mul_pd(dvdasum, _mm_mul_pd(isai0,isai0));
337 gmx_mm_update_1pot_pd(dvdasum,dvda+inr);
339 /* Increment number of inner iterations */
340 inneriter += j_index_end - j_index_start;
342 /* Outer loop uses 9 flops */
345 /* Increment number of outer iterations */
348 /* Update outer/inner flops */
350 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VF,outeriter*9 + inneriter*61);
353 * Gromacs nonbonded kernel: nb_kernel_ElecGB_VdwNone_GeomP1P1_F_avx_128_fma_double
354 * Electrostatics interaction: GeneralizedBorn
355 * VdW interaction: None
356 * Geometry: Particle-Particle
357 * Calculate force/pot: Force
360 nb_kernel_ElecGB_VdwNone_GeomP1P1_F_avx_128_fma_double
361 (t_nblist * gmx_restrict nlist,
362 rvec * gmx_restrict xx,
363 rvec * gmx_restrict ff,
364 t_forcerec * gmx_restrict fr,
365 t_mdatoms * gmx_restrict mdatoms,
366 nb_kernel_data_t * gmx_restrict kernel_data,
367 t_nrnb * gmx_restrict nrnb)
369 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
370 * just 0 for non-waters.
371 * Suffixes A,B refer to j loop unrolling done with SSE double precision, e.g. for the two different
372 * jnr indices corresponding to data put in the four positions in the SIMD register.
374 int i_shift_offset,i_coord_offset,outeriter,inneriter;
375 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
377 int j_coord_offsetA,j_coord_offsetB;
378 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
380 real *shiftvec,*fshift,*x,*f;
381 __m128d tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
383 __m128d ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
384 int vdwjidx0A,vdwjidx0B;
385 __m128d jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
386 __m128d dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
387 __m128d velec,felec,velecsum,facel,crf,krf,krf2;
390 __m128d vgb,fgb,vgbsum,dvdasum,gbscale,gbtabscale,isaprod,gbqqfactor,gbinvepsdiff,dvdaj,gbeps,twogbeps,dvdatmp;
391 __m128d minushalf = _mm_set1_pd(-0.5);
392 real *invsqrta,*dvda,*gbtab;
394 __m128i ifour = _mm_set1_epi32(4);
395 __m128d rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF,twovfeps;
397 __m128d dummy_mask,cutoff_mask;
398 __m128d signbit = gmx_mm_castsi128_pd( _mm_set_epi32(0x80000000,0x00000000,0x80000000,0x00000000) );
399 __m128d one = _mm_set1_pd(1.0);
400 __m128d two = _mm_set1_pd(2.0);
406 jindex = nlist->jindex;
408 shiftidx = nlist->shift;
410 shiftvec = fr->shift_vec[0];
411 fshift = fr->fshift[0];
412 facel = _mm_set1_pd(fr->epsfac);
413 charge = mdatoms->chargeA;
415 invsqrta = fr->invsqrta;
417 gbtabscale = _mm_set1_pd(fr->gbtab.scale);
418 gbtab = fr->gbtab.data;
419 gbinvepsdiff = _mm_set1_pd((1.0/fr->epsilon_r) - (1.0/fr->gb_epsilon_solvent));
421 /* Avoid stupid compiler warnings */
429 /* Start outer loop over neighborlists */
430 for(iidx=0; iidx<nri; iidx++)
432 /* Load shift vector for this list */
433 i_shift_offset = DIM*shiftidx[iidx];
435 /* Load limits for loop over neighbors */
436 j_index_start = jindex[iidx];
437 j_index_end = jindex[iidx+1];
439 /* Get outer coordinate index */
441 i_coord_offset = DIM*inr;
443 /* Load i particle coords and add shift vector */
444 gmx_mm_load_shift_and_1rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
446 fix0 = _mm_setzero_pd();
447 fiy0 = _mm_setzero_pd();
448 fiz0 = _mm_setzero_pd();
450 /* Load parameters for i particles */
451 iq0 = _mm_mul_pd(facel,_mm_load1_pd(charge+inr+0));
452 isai0 = _mm_load1_pd(invsqrta+inr+0);
454 dvdasum = _mm_setzero_pd();
456 /* Start inner kernel loop */
457 for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
460 /* Get j neighbor index, and coordinate index */
463 j_coord_offsetA = DIM*jnrA;
464 j_coord_offsetB = DIM*jnrB;
466 /* load j atom coordinates */
467 gmx_mm_load_1rvec_2ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
470 /* Calculate displacement vector */
471 dx00 = _mm_sub_pd(ix0,jx0);
472 dy00 = _mm_sub_pd(iy0,jy0);
473 dz00 = _mm_sub_pd(iz0,jz0);
475 /* Calculate squared distance and things based on it */
476 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
478 rinv00 = gmx_mm_invsqrt_pd(rsq00);
480 /* Load parameters for j particles */
481 jq0 = gmx_mm_load_2real_swizzle_pd(charge+jnrA+0,charge+jnrB+0);
482 isaj0 = gmx_mm_load_2real_swizzle_pd(invsqrta+jnrA+0,invsqrta+jnrB+0);
484 /**************************
485 * CALCULATE INTERACTIONS *
486 **************************/
488 r00 = _mm_mul_pd(rsq00,rinv00);
490 /* Compute parameters for interactions between i and j atoms */
491 qq00 = _mm_mul_pd(iq0,jq0);
493 /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
494 isaprod = _mm_mul_pd(isai0,isaj0);
495 gbqqfactor = _mm_xor_pd(signbit,_mm_mul_pd(qq00,_mm_mul_pd(isaprod,gbinvepsdiff)));
496 gbscale = _mm_mul_pd(isaprod,gbtabscale);
498 /* Calculate generalized born table index - this is a separate table from the normal one,
499 * but we use the same procedure by multiplying r with scale and truncating to integer.
501 rt = _mm_mul_pd(r00,gbscale);
502 gbitab = _mm_cvttpd_epi32(rt);
504 gbeps = _mm_frcz_pd(rt);
506 gbeps = _mm_sub_pd(rt,_mm_round_pd(rt, _MM_FROUND_FLOOR));
508 gbitab = _mm_slli_epi32(gbitab,2);
510 Y = _mm_load_pd( gbtab + _mm_extract_epi32(gbitab,0) );
511 F = _mm_load_pd( gbtab + _mm_extract_epi32(gbitab,1) );
512 GMX_MM_TRANSPOSE2_PD(Y,F);
513 G = _mm_load_pd( gbtab + _mm_extract_epi32(gbitab,0) +2);
514 H = _mm_load_pd( gbtab + _mm_extract_epi32(gbitab,1) +2);
515 GMX_MM_TRANSPOSE2_PD(G,H);
516 Fp = _mm_macc_pd(gbeps,_mm_macc_pd(gbeps,H,G),F);
517 VV = _mm_macc_pd(gbeps,Fp,Y);
518 vgb = _mm_mul_pd(gbqqfactor,VV);
520 twogbeps = _mm_add_pd(gbeps,gbeps);
521 FF = _mm_macc_pd(_mm_macc_pd(twogbeps,H,G),gbeps,Fp);
522 fgb = _mm_mul_pd(gbqqfactor,_mm_mul_pd(FF,gbscale));
523 dvdatmp = _mm_mul_pd(minushalf,_mm_macc_pd(fgb,r00,vgb));
524 dvdasum = _mm_add_pd(dvdasum,dvdatmp);
525 gmx_mm_increment_2real_swizzle_pd(dvda+jnrA,dvda+jnrB,_mm_mul_pd(dvdatmp,_mm_mul_pd(isaj0,isaj0)));
526 velec = _mm_mul_pd(qq00,rinv00);
527 felec = _mm_mul_pd(_mm_msub_pd(velec,rinv00,fgb),rinv00);
531 /* Update vectorial force */
532 fix0 = _mm_macc_pd(dx00,fscal,fix0);
533 fiy0 = _mm_macc_pd(dy00,fscal,fiy0);
534 fiz0 = _mm_macc_pd(dz00,fscal,fiz0);
536 gmx_mm_decrement_1rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,
537 _mm_mul_pd(dx00,fscal),
538 _mm_mul_pd(dy00,fscal),
539 _mm_mul_pd(dz00,fscal));
541 /* Inner loop uses 59 flops */
548 j_coord_offsetA = DIM*jnrA;
550 /* load j atom coordinates */
551 gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
554 /* Calculate displacement vector */
555 dx00 = _mm_sub_pd(ix0,jx0);
556 dy00 = _mm_sub_pd(iy0,jy0);
557 dz00 = _mm_sub_pd(iz0,jz0);
559 /* Calculate squared distance and things based on it */
560 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
562 rinv00 = gmx_mm_invsqrt_pd(rsq00);
564 /* Load parameters for j particles */
565 jq0 = _mm_load_sd(charge+jnrA+0);
566 isaj0 = _mm_load_sd(invsqrta+jnrA+0);
568 /**************************
569 * CALCULATE INTERACTIONS *
570 **************************/
572 r00 = _mm_mul_pd(rsq00,rinv00);
574 /* Compute parameters for interactions between i and j atoms */
575 qq00 = _mm_mul_pd(iq0,jq0);
577 /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
578 isaprod = _mm_mul_pd(isai0,isaj0);
579 gbqqfactor = _mm_xor_pd(signbit,_mm_mul_pd(qq00,_mm_mul_pd(isaprod,gbinvepsdiff)));
580 gbscale = _mm_mul_pd(isaprod,gbtabscale);
582 /* Calculate generalized born table index - this is a separate table from the normal one,
583 * but we use the same procedure by multiplying r with scale and truncating to integer.
585 rt = _mm_mul_pd(r00,gbscale);
586 gbitab = _mm_cvttpd_epi32(rt);
588 gbeps = _mm_frcz_pd(rt);
590 gbeps = _mm_sub_pd(rt,_mm_round_pd(rt, _MM_FROUND_FLOOR));
592 gbitab = _mm_slli_epi32(gbitab,2);
594 Y = _mm_load_pd( gbtab + _mm_extract_epi32(gbitab,0) );
595 F = _mm_setzero_pd();
596 GMX_MM_TRANSPOSE2_PD(Y,F);
597 G = _mm_load_pd( gbtab + _mm_extract_epi32(gbitab,0) +2);
598 H = _mm_setzero_pd();
599 GMX_MM_TRANSPOSE2_PD(G,H);
600 Fp = _mm_macc_pd(gbeps,_mm_macc_pd(gbeps,H,G),F);
601 VV = _mm_macc_pd(gbeps,Fp,Y);
602 vgb = _mm_mul_pd(gbqqfactor,VV);
604 twogbeps = _mm_add_pd(gbeps,gbeps);
605 FF = _mm_macc_pd(_mm_macc_pd(twogbeps,H,G),gbeps,Fp);
606 fgb = _mm_mul_pd(gbqqfactor,_mm_mul_pd(FF,gbscale));
607 dvdatmp = _mm_mul_pd(minushalf,_mm_macc_pd(fgb,r00,vgb));
608 dvdasum = _mm_add_pd(dvdasum,dvdatmp);
609 gmx_mm_increment_1real_pd(dvda+jnrA,_mm_mul_pd(dvdatmp,_mm_mul_pd(isaj0,isaj0)));
610 velec = _mm_mul_pd(qq00,rinv00);
611 felec = _mm_mul_pd(_mm_msub_pd(velec,rinv00,fgb),rinv00);
615 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
617 /* Update vectorial force */
618 fix0 = _mm_macc_pd(dx00,fscal,fix0);
619 fiy0 = _mm_macc_pd(dy00,fscal,fiy0);
620 fiz0 = _mm_macc_pd(dz00,fscal,fiz0);
622 gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,
623 _mm_mul_pd(dx00,fscal),
624 _mm_mul_pd(dy00,fscal),
625 _mm_mul_pd(dz00,fscal));
627 /* Inner loop uses 59 flops */
630 /* End of innermost loop */
632 gmx_mm_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
633 f+i_coord_offset,fshift+i_shift_offset);
635 dvdasum = _mm_mul_pd(dvdasum, _mm_mul_pd(isai0,isai0));
636 gmx_mm_update_1pot_pd(dvdasum,dvda+inr);
638 /* Increment number of inner iterations */
639 inneriter += j_index_end - j_index_start;
641 /* Outer loop uses 7 flops */
644 /* Increment number of outer iterations */
647 /* Update outer/inner flops */
649 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_F,outeriter*7 + inneriter*59);