2 * Note: this file was generated by the Gromacs sse2_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_sse2_double.h"
34 #include "kernelutil_x86_sse2_double.h"
37 * Gromacs nonbonded kernel: nb_kernel_ElecGB_VdwNone_GeomP1P1_VF_sse2_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_sse2_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,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;
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);
190 gbeps = _mm_sub_pd(rt,_mm_cvtepi32_pd(gbitab));
191 gbitab = _mm_slli_epi32(gbitab,2);
193 Y = _mm_load_pd( gbtab + gmx_mm_extract_epi32(gbitab,0) );
194 F = _mm_load_pd( gbtab + gmx_mm_extract_epi32(gbitab,1) );
195 GMX_MM_TRANSPOSE2_PD(Y,F);
196 G = _mm_load_pd( gbtab + gmx_mm_extract_epi32(gbitab,0) +2);
197 H = _mm_load_pd( gbtab + gmx_mm_extract_epi32(gbitab,1) +2);
198 GMX_MM_TRANSPOSE2_PD(G,H);
199 Heps = _mm_mul_pd(gbeps,H);
200 Fp = _mm_add_pd(F,_mm_mul_pd(gbeps,_mm_add_pd(G,Heps)));
201 VV = _mm_add_pd(Y,_mm_mul_pd(gbeps,Fp));
202 vgb = _mm_mul_pd(gbqqfactor,VV);
204 FF = _mm_add_pd(Fp,_mm_mul_pd(gbeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
205 fgb = _mm_mul_pd(gbqqfactor,_mm_mul_pd(FF,gbscale));
206 dvdatmp = _mm_mul_pd(minushalf,_mm_add_pd(vgb,_mm_mul_pd(fgb,r00)));
207 dvdasum = _mm_add_pd(dvdasum,dvdatmp);
208 gmx_mm_increment_2real_swizzle_pd(dvda+jnrA,dvda+jnrB,_mm_mul_pd(dvdatmp,_mm_mul_pd(isaj0,isaj0)));
209 velec = _mm_mul_pd(qq00,rinv00);
210 felec = _mm_mul_pd(_mm_sub_pd(_mm_mul_pd(velec,rinv00),fgb),rinv00);
212 /* Update potential sum for this i atom from the interaction with this j atom. */
213 velecsum = _mm_add_pd(velecsum,velec);
214 vgbsum = _mm_add_pd(vgbsum,vgb);
218 /* Calculate temporary vectorial force */
219 tx = _mm_mul_pd(fscal,dx00);
220 ty = _mm_mul_pd(fscal,dy00);
221 tz = _mm_mul_pd(fscal,dz00);
223 /* Update vectorial force */
224 fix0 = _mm_add_pd(fix0,tx);
225 fiy0 = _mm_add_pd(fiy0,ty);
226 fiz0 = _mm_add_pd(fiz0,tz);
228 gmx_mm_decrement_1rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,tx,ty,tz);
230 /* Inner loop uses 58 flops */
237 j_coord_offsetA = DIM*jnrA;
239 /* load j atom coordinates */
240 gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
243 /* Calculate displacement vector */
244 dx00 = _mm_sub_pd(ix0,jx0);
245 dy00 = _mm_sub_pd(iy0,jy0);
246 dz00 = _mm_sub_pd(iz0,jz0);
248 /* Calculate squared distance and things based on it */
249 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
251 rinv00 = gmx_mm_invsqrt_pd(rsq00);
253 /* Load parameters for j particles */
254 jq0 = _mm_load_sd(charge+jnrA+0);
255 isaj0 = _mm_load_sd(invsqrta+jnrA+0);
257 /**************************
258 * CALCULATE INTERACTIONS *
259 **************************/
261 r00 = _mm_mul_pd(rsq00,rinv00);
263 /* Compute parameters for interactions between i and j atoms */
264 qq00 = _mm_mul_pd(iq0,jq0);
266 /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
267 isaprod = _mm_mul_pd(isai0,isaj0);
268 gbqqfactor = _mm_xor_pd(signbit,_mm_mul_pd(qq00,_mm_mul_pd(isaprod,gbinvepsdiff)));
269 gbscale = _mm_mul_pd(isaprod,gbtabscale);
271 /* Calculate generalized born table index - this is a separate table from the normal one,
272 * but we use the same procedure by multiplying r with scale and truncating to integer.
274 rt = _mm_mul_pd(r00,gbscale);
275 gbitab = _mm_cvttpd_epi32(rt);
276 gbeps = _mm_sub_pd(rt,_mm_cvtepi32_pd(gbitab));
277 gbitab = _mm_slli_epi32(gbitab,2);
279 Y = _mm_load_pd( gbtab + gmx_mm_extract_epi32(gbitab,0) );
280 F = _mm_setzero_pd();
281 GMX_MM_TRANSPOSE2_PD(Y,F);
282 G = _mm_load_pd( gbtab + gmx_mm_extract_epi32(gbitab,0) +2);
283 H = _mm_setzero_pd();
284 GMX_MM_TRANSPOSE2_PD(G,H);
285 Heps = _mm_mul_pd(gbeps,H);
286 Fp = _mm_add_pd(F,_mm_mul_pd(gbeps,_mm_add_pd(G,Heps)));
287 VV = _mm_add_pd(Y,_mm_mul_pd(gbeps,Fp));
288 vgb = _mm_mul_pd(gbqqfactor,VV);
290 FF = _mm_add_pd(Fp,_mm_mul_pd(gbeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
291 fgb = _mm_mul_pd(gbqqfactor,_mm_mul_pd(FF,gbscale));
292 dvdatmp = _mm_mul_pd(minushalf,_mm_add_pd(vgb,_mm_mul_pd(fgb,r00)));
293 dvdasum = _mm_add_pd(dvdasum,dvdatmp);
294 gmx_mm_increment_1real_pd(dvda+jnrA,_mm_mul_pd(dvdatmp,_mm_mul_pd(isaj0,isaj0)));
295 velec = _mm_mul_pd(qq00,rinv00);
296 felec = _mm_mul_pd(_mm_sub_pd(_mm_mul_pd(velec,rinv00),fgb),rinv00);
298 /* Update potential sum for this i atom from the interaction with this j atom. */
299 velec = _mm_unpacklo_pd(velec,_mm_setzero_pd());
300 velecsum = _mm_add_pd(velecsum,velec);
301 vgb = _mm_unpacklo_pd(vgb,_mm_setzero_pd());
302 vgbsum = _mm_add_pd(vgbsum,vgb);
306 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
308 /* Calculate temporary vectorial force */
309 tx = _mm_mul_pd(fscal,dx00);
310 ty = _mm_mul_pd(fscal,dy00);
311 tz = _mm_mul_pd(fscal,dz00);
313 /* Update vectorial force */
314 fix0 = _mm_add_pd(fix0,tx);
315 fiy0 = _mm_add_pd(fiy0,ty);
316 fiz0 = _mm_add_pd(fiz0,tz);
318 gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,tx,ty,tz);
320 /* Inner loop uses 58 flops */
323 /* End of innermost loop */
325 gmx_mm_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
326 f+i_coord_offset,fshift+i_shift_offset);
329 /* Update potential energies */
330 gmx_mm_update_1pot_pd(velecsum,kernel_data->energygrp_elec+ggid);
331 gmx_mm_update_1pot_pd(vgbsum,kernel_data->energygrp_polarization+ggid);
332 dvdasum = _mm_mul_pd(dvdasum, _mm_mul_pd(isai0,isai0));
333 gmx_mm_update_1pot_pd(dvdasum,dvda+inr);
335 /* Increment number of inner iterations */
336 inneriter += j_index_end - j_index_start;
338 /* Outer loop uses 9 flops */
341 /* Increment number of outer iterations */
344 /* Update outer/inner flops */
346 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VF,outeriter*9 + inneriter*58);
349 * Gromacs nonbonded kernel: nb_kernel_ElecGB_VdwNone_GeomP1P1_F_sse2_double
350 * Electrostatics interaction: GeneralizedBorn
351 * VdW interaction: None
352 * Geometry: Particle-Particle
353 * Calculate force/pot: Force
356 nb_kernel_ElecGB_VdwNone_GeomP1P1_F_sse2_double
357 (t_nblist * gmx_restrict nlist,
358 rvec * gmx_restrict xx,
359 rvec * gmx_restrict ff,
360 t_forcerec * gmx_restrict fr,
361 t_mdatoms * gmx_restrict mdatoms,
362 nb_kernel_data_t * gmx_restrict kernel_data,
363 t_nrnb * gmx_restrict nrnb)
365 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
366 * just 0 for non-waters.
367 * Suffixes A,B refer to j loop unrolling done with SSE double precision, e.g. for the two different
368 * jnr indices corresponding to data put in the four positions in the SIMD register.
370 int i_shift_offset,i_coord_offset,outeriter,inneriter;
371 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
373 int j_coord_offsetA,j_coord_offsetB;
374 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
376 real *shiftvec,*fshift,*x,*f;
377 __m128d tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
379 __m128d ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
380 int vdwjidx0A,vdwjidx0B;
381 __m128d jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
382 __m128d dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
383 __m128d velec,felec,velecsum,facel,crf,krf,krf2;
386 __m128d vgb,fgb,vgbsum,dvdasum,gbscale,gbtabscale,isaprod,gbqqfactor,gbinvepsdiff,dvdaj,gbeps,dvdatmp;
387 __m128d minushalf = _mm_set1_pd(-0.5);
388 real *invsqrta,*dvda,*gbtab;
390 __m128i ifour = _mm_set1_epi32(4);
391 __m128d rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF;
393 __m128d dummy_mask,cutoff_mask;
394 __m128d signbit = gmx_mm_castsi128_pd( _mm_set_epi32(0x80000000,0x00000000,0x80000000,0x00000000) );
395 __m128d one = _mm_set1_pd(1.0);
396 __m128d two = _mm_set1_pd(2.0);
402 jindex = nlist->jindex;
404 shiftidx = nlist->shift;
406 shiftvec = fr->shift_vec[0];
407 fshift = fr->fshift[0];
408 facel = _mm_set1_pd(fr->epsfac);
409 charge = mdatoms->chargeA;
411 invsqrta = fr->invsqrta;
413 gbtabscale = _mm_set1_pd(fr->gbtab.scale);
414 gbtab = fr->gbtab.data;
415 gbinvepsdiff = _mm_set1_pd((1.0/fr->epsilon_r) - (1.0/fr->gb_epsilon_solvent));
417 /* Avoid stupid compiler warnings */
425 /* Start outer loop over neighborlists */
426 for(iidx=0; iidx<nri; iidx++)
428 /* Load shift vector for this list */
429 i_shift_offset = DIM*shiftidx[iidx];
431 /* Load limits for loop over neighbors */
432 j_index_start = jindex[iidx];
433 j_index_end = jindex[iidx+1];
435 /* Get outer coordinate index */
437 i_coord_offset = DIM*inr;
439 /* Load i particle coords and add shift vector */
440 gmx_mm_load_shift_and_1rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
442 fix0 = _mm_setzero_pd();
443 fiy0 = _mm_setzero_pd();
444 fiz0 = _mm_setzero_pd();
446 /* Load parameters for i particles */
447 iq0 = _mm_mul_pd(facel,_mm_load1_pd(charge+inr+0));
448 isai0 = _mm_load1_pd(invsqrta+inr+0);
450 dvdasum = _mm_setzero_pd();
452 /* Start inner kernel loop */
453 for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
456 /* Get j neighbor index, and coordinate index */
459 j_coord_offsetA = DIM*jnrA;
460 j_coord_offsetB = DIM*jnrB;
462 /* load j atom coordinates */
463 gmx_mm_load_1rvec_2ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
466 /* Calculate displacement vector */
467 dx00 = _mm_sub_pd(ix0,jx0);
468 dy00 = _mm_sub_pd(iy0,jy0);
469 dz00 = _mm_sub_pd(iz0,jz0);
471 /* Calculate squared distance and things based on it */
472 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
474 rinv00 = gmx_mm_invsqrt_pd(rsq00);
476 /* Load parameters for j particles */
477 jq0 = gmx_mm_load_2real_swizzle_pd(charge+jnrA+0,charge+jnrB+0);
478 isaj0 = gmx_mm_load_2real_swizzle_pd(invsqrta+jnrA+0,invsqrta+jnrB+0);
480 /**************************
481 * CALCULATE INTERACTIONS *
482 **************************/
484 r00 = _mm_mul_pd(rsq00,rinv00);
486 /* Compute parameters for interactions between i and j atoms */
487 qq00 = _mm_mul_pd(iq0,jq0);
489 /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
490 isaprod = _mm_mul_pd(isai0,isaj0);
491 gbqqfactor = _mm_xor_pd(signbit,_mm_mul_pd(qq00,_mm_mul_pd(isaprod,gbinvepsdiff)));
492 gbscale = _mm_mul_pd(isaprod,gbtabscale);
494 /* Calculate generalized born table index - this is a separate table from the normal one,
495 * but we use the same procedure by multiplying r with scale and truncating to integer.
497 rt = _mm_mul_pd(r00,gbscale);
498 gbitab = _mm_cvttpd_epi32(rt);
499 gbeps = _mm_sub_pd(rt,_mm_cvtepi32_pd(gbitab));
500 gbitab = _mm_slli_epi32(gbitab,2);
502 Y = _mm_load_pd( gbtab + gmx_mm_extract_epi32(gbitab,0) );
503 F = _mm_load_pd( gbtab + gmx_mm_extract_epi32(gbitab,1) );
504 GMX_MM_TRANSPOSE2_PD(Y,F);
505 G = _mm_load_pd( gbtab + gmx_mm_extract_epi32(gbitab,0) +2);
506 H = _mm_load_pd( gbtab + gmx_mm_extract_epi32(gbitab,1) +2);
507 GMX_MM_TRANSPOSE2_PD(G,H);
508 Heps = _mm_mul_pd(gbeps,H);
509 Fp = _mm_add_pd(F,_mm_mul_pd(gbeps,_mm_add_pd(G,Heps)));
510 VV = _mm_add_pd(Y,_mm_mul_pd(gbeps,Fp));
511 vgb = _mm_mul_pd(gbqqfactor,VV);
513 FF = _mm_add_pd(Fp,_mm_mul_pd(gbeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
514 fgb = _mm_mul_pd(gbqqfactor,_mm_mul_pd(FF,gbscale));
515 dvdatmp = _mm_mul_pd(minushalf,_mm_add_pd(vgb,_mm_mul_pd(fgb,r00)));
516 dvdasum = _mm_add_pd(dvdasum,dvdatmp);
517 gmx_mm_increment_2real_swizzle_pd(dvda+jnrA,dvda+jnrB,_mm_mul_pd(dvdatmp,_mm_mul_pd(isaj0,isaj0)));
518 velec = _mm_mul_pd(qq00,rinv00);
519 felec = _mm_mul_pd(_mm_sub_pd(_mm_mul_pd(velec,rinv00),fgb),rinv00);
523 /* Calculate temporary vectorial force */
524 tx = _mm_mul_pd(fscal,dx00);
525 ty = _mm_mul_pd(fscal,dy00);
526 tz = _mm_mul_pd(fscal,dz00);
528 /* Update vectorial force */
529 fix0 = _mm_add_pd(fix0,tx);
530 fiy0 = _mm_add_pd(fiy0,ty);
531 fiz0 = _mm_add_pd(fiz0,tz);
533 gmx_mm_decrement_1rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,tx,ty,tz);
535 /* Inner loop uses 56 flops */
542 j_coord_offsetA = DIM*jnrA;
544 /* load j atom coordinates */
545 gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
548 /* Calculate displacement vector */
549 dx00 = _mm_sub_pd(ix0,jx0);
550 dy00 = _mm_sub_pd(iy0,jy0);
551 dz00 = _mm_sub_pd(iz0,jz0);
553 /* Calculate squared distance and things based on it */
554 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
556 rinv00 = gmx_mm_invsqrt_pd(rsq00);
558 /* Load parameters for j particles */
559 jq0 = _mm_load_sd(charge+jnrA+0);
560 isaj0 = _mm_load_sd(invsqrta+jnrA+0);
562 /**************************
563 * CALCULATE INTERACTIONS *
564 **************************/
566 r00 = _mm_mul_pd(rsq00,rinv00);
568 /* Compute parameters for interactions between i and j atoms */
569 qq00 = _mm_mul_pd(iq0,jq0);
571 /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
572 isaprod = _mm_mul_pd(isai0,isaj0);
573 gbqqfactor = _mm_xor_pd(signbit,_mm_mul_pd(qq00,_mm_mul_pd(isaprod,gbinvepsdiff)));
574 gbscale = _mm_mul_pd(isaprod,gbtabscale);
576 /* Calculate generalized born table index - this is a separate table from the normal one,
577 * but we use the same procedure by multiplying r with scale and truncating to integer.
579 rt = _mm_mul_pd(r00,gbscale);
580 gbitab = _mm_cvttpd_epi32(rt);
581 gbeps = _mm_sub_pd(rt,_mm_cvtepi32_pd(gbitab));
582 gbitab = _mm_slli_epi32(gbitab,2);
584 Y = _mm_load_pd( gbtab + gmx_mm_extract_epi32(gbitab,0) );
585 F = _mm_setzero_pd();
586 GMX_MM_TRANSPOSE2_PD(Y,F);
587 G = _mm_load_pd( gbtab + gmx_mm_extract_epi32(gbitab,0) +2);
588 H = _mm_setzero_pd();
589 GMX_MM_TRANSPOSE2_PD(G,H);
590 Heps = _mm_mul_pd(gbeps,H);
591 Fp = _mm_add_pd(F,_mm_mul_pd(gbeps,_mm_add_pd(G,Heps)));
592 VV = _mm_add_pd(Y,_mm_mul_pd(gbeps,Fp));
593 vgb = _mm_mul_pd(gbqqfactor,VV);
595 FF = _mm_add_pd(Fp,_mm_mul_pd(gbeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
596 fgb = _mm_mul_pd(gbqqfactor,_mm_mul_pd(FF,gbscale));
597 dvdatmp = _mm_mul_pd(minushalf,_mm_add_pd(vgb,_mm_mul_pd(fgb,r00)));
598 dvdasum = _mm_add_pd(dvdasum,dvdatmp);
599 gmx_mm_increment_1real_pd(dvda+jnrA,_mm_mul_pd(dvdatmp,_mm_mul_pd(isaj0,isaj0)));
600 velec = _mm_mul_pd(qq00,rinv00);
601 felec = _mm_mul_pd(_mm_sub_pd(_mm_mul_pd(velec,rinv00),fgb),rinv00);
605 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
607 /* Calculate temporary vectorial force */
608 tx = _mm_mul_pd(fscal,dx00);
609 ty = _mm_mul_pd(fscal,dy00);
610 tz = _mm_mul_pd(fscal,dz00);
612 /* Update vectorial force */
613 fix0 = _mm_add_pd(fix0,tx);
614 fiy0 = _mm_add_pd(fiy0,ty);
615 fiz0 = _mm_add_pd(fiz0,tz);
617 gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,tx,ty,tz);
619 /* Inner loop uses 56 flops */
622 /* End of innermost loop */
624 gmx_mm_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
625 f+i_coord_offset,fshift+i_shift_offset);
627 dvdasum = _mm_mul_pd(dvdasum, _mm_mul_pd(isai0,isai0));
628 gmx_mm_update_1pot_pd(dvdasum,dvda+inr);
630 /* Increment number of inner iterations */
631 inneriter += j_index_end - j_index_start;
633 /* Outer loop uses 7 flops */
636 /* Increment number of outer iterations */
639 /* Update outer/inner flops */
641 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_F,outeriter*7 + inneriter*56);