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_double 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_double.h"
50 #include "kernelutil_x86_avx_128_fma_double.h"
53 * Gromacs nonbonded kernel: nb_kernel_ElecGB_VdwNone_GeomP1P1_VF_avx_128_fma_double
54 * Electrostatics interaction: GeneralizedBorn
55 * VdW interaction: None
56 * Geometry: Particle-Particle
57 * Calculate force/pot: PotentialAndForce
60 nb_kernel_ElecGB_VdwNone_GeomP1P1_VF_avx_128_fma_double
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 refer to j loop unrolling done with SSE double precision, e.g. for the two 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;
77 int j_coord_offsetA,j_coord_offsetB;
78 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
80 real *shiftvec,*fshift,*x,*f;
81 __m128d tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
83 __m128d ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
84 int vdwjidx0A,vdwjidx0B;
85 __m128d jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
86 __m128d dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
87 __m128d velec,felec,velecsum,facel,crf,krf,krf2;
90 __m128d vgb,fgb,vgbsum,dvdasum,gbscale,gbtabscale,isaprod,gbqqfactor,gbinvepsdiff,dvdaj,gbeps,twogbeps,dvdatmp;
91 __m128d minushalf = _mm_set1_pd(-0.5);
92 real *invsqrta,*dvda,*gbtab;
94 __m128i ifour = _mm_set1_epi32(4);
95 __m128d rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF,twovfeps;
97 __m128d dummy_mask,cutoff_mask;
98 __m128d signbit = gmx_mm_castsi128_pd( _mm_set_epi32(0x80000000,0x00000000,0x80000000,0x00000000) );
99 __m128d one = _mm_set1_pd(1.0);
100 __m128d two = _mm_set1_pd(2.0);
106 jindex = nlist->jindex;
108 shiftidx = nlist->shift;
110 shiftvec = fr->shift_vec[0];
111 fshift = fr->fshift[0];
112 facel = _mm_set1_pd(fr->epsfac);
113 charge = mdatoms->chargeA;
115 invsqrta = fr->invsqrta;
117 gbtabscale = _mm_set1_pd(fr->gbtab.scale);
118 gbtab = fr->gbtab.data;
119 gbinvepsdiff = _mm_set1_pd((1.0/fr->epsilon_r) - (1.0/fr->gb_epsilon_solvent));
121 /* Avoid stupid compiler warnings */
129 /* Start outer loop over neighborlists */
130 for(iidx=0; iidx<nri; iidx++)
132 /* Load shift vector for this list */
133 i_shift_offset = DIM*shiftidx[iidx];
135 /* Load limits for loop over neighbors */
136 j_index_start = jindex[iidx];
137 j_index_end = jindex[iidx+1];
139 /* Get outer coordinate index */
141 i_coord_offset = DIM*inr;
143 /* Load i particle coords and add shift vector */
144 gmx_mm_load_shift_and_1rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
146 fix0 = _mm_setzero_pd();
147 fiy0 = _mm_setzero_pd();
148 fiz0 = _mm_setzero_pd();
150 /* Load parameters for i particles */
151 iq0 = _mm_mul_pd(facel,_mm_load1_pd(charge+inr+0));
152 isai0 = _mm_load1_pd(invsqrta+inr+0);
154 /* Reset potential sums */
155 velecsum = _mm_setzero_pd();
156 vgbsum = _mm_setzero_pd();
157 dvdasum = _mm_setzero_pd();
159 /* Start inner kernel loop */
160 for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
163 /* Get j neighbor index, and coordinate index */
166 j_coord_offsetA = DIM*jnrA;
167 j_coord_offsetB = DIM*jnrB;
169 /* load j atom coordinates */
170 gmx_mm_load_1rvec_2ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
173 /* Calculate displacement vector */
174 dx00 = _mm_sub_pd(ix0,jx0);
175 dy00 = _mm_sub_pd(iy0,jy0);
176 dz00 = _mm_sub_pd(iz0,jz0);
178 /* Calculate squared distance and things based on it */
179 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
181 rinv00 = gmx_mm_invsqrt_pd(rsq00);
183 /* Load parameters for j particles */
184 jq0 = gmx_mm_load_2real_swizzle_pd(charge+jnrA+0,charge+jnrB+0);
185 isaj0 = gmx_mm_load_2real_swizzle_pd(invsqrta+jnrA+0,invsqrta+jnrB+0);
187 /**************************
188 * CALCULATE INTERACTIONS *
189 **************************/
191 r00 = _mm_mul_pd(rsq00,rinv00);
193 /* Compute parameters for interactions between i and j atoms */
194 qq00 = _mm_mul_pd(iq0,jq0);
196 /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
197 isaprod = _mm_mul_pd(isai0,isaj0);
198 gbqqfactor = _mm_xor_pd(signbit,_mm_mul_pd(qq00,_mm_mul_pd(isaprod,gbinvepsdiff)));
199 gbscale = _mm_mul_pd(isaprod,gbtabscale);
201 /* Calculate generalized born table index - this is a separate table from the normal one,
202 * but we use the same procedure by multiplying r with scale and truncating to integer.
204 rt = _mm_mul_pd(r00,gbscale);
205 gbitab = _mm_cvttpd_epi32(rt);
207 gbeps = _mm_frcz_pd(rt);
209 gbeps = _mm_sub_pd(rt,_mm_round_pd(rt, _MM_FROUND_FLOOR));
211 gbitab = _mm_slli_epi32(gbitab,2);
213 Y = _mm_load_pd( gbtab + _mm_extract_epi32(gbitab,0) );
214 F = _mm_load_pd( gbtab + _mm_extract_epi32(gbitab,1) );
215 GMX_MM_TRANSPOSE2_PD(Y,F);
216 G = _mm_load_pd( gbtab + _mm_extract_epi32(gbitab,0) +2);
217 H = _mm_load_pd( gbtab + _mm_extract_epi32(gbitab,1) +2);
218 GMX_MM_TRANSPOSE2_PD(G,H);
219 Fp = _mm_macc_pd(gbeps,_mm_macc_pd(gbeps,H,G),F);
220 VV = _mm_macc_pd(gbeps,Fp,Y);
221 vgb = _mm_mul_pd(gbqqfactor,VV);
223 twogbeps = _mm_add_pd(gbeps,gbeps);
224 FF = _mm_macc_pd(_mm_macc_pd(twogbeps,H,G),gbeps,Fp);
225 fgb = _mm_mul_pd(gbqqfactor,_mm_mul_pd(FF,gbscale));
226 dvdatmp = _mm_mul_pd(minushalf,_mm_macc_pd(fgb,r00,vgb));
227 dvdasum = _mm_add_pd(dvdasum,dvdatmp);
228 gmx_mm_increment_2real_swizzle_pd(dvda+jnrA,dvda+jnrB,_mm_mul_pd(dvdatmp,_mm_mul_pd(isaj0,isaj0)));
229 velec = _mm_mul_pd(qq00,rinv00);
230 felec = _mm_mul_pd(_mm_msub_pd(velec,rinv00,fgb),rinv00);
232 /* Update potential sum for this i atom from the interaction with this j atom. */
233 velecsum = _mm_add_pd(velecsum,velec);
234 vgbsum = _mm_add_pd(vgbsum,vgb);
238 /* Update vectorial force */
239 fix0 = _mm_macc_pd(dx00,fscal,fix0);
240 fiy0 = _mm_macc_pd(dy00,fscal,fiy0);
241 fiz0 = _mm_macc_pd(dz00,fscal,fiz0);
243 gmx_mm_decrement_1rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,
244 _mm_mul_pd(dx00,fscal),
245 _mm_mul_pd(dy00,fscal),
246 _mm_mul_pd(dz00,fscal));
248 /* Inner loop uses 61 flops */
255 j_coord_offsetA = DIM*jnrA;
257 /* load j atom coordinates */
258 gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
261 /* Calculate displacement vector */
262 dx00 = _mm_sub_pd(ix0,jx0);
263 dy00 = _mm_sub_pd(iy0,jy0);
264 dz00 = _mm_sub_pd(iz0,jz0);
266 /* Calculate squared distance and things based on it */
267 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
269 rinv00 = gmx_mm_invsqrt_pd(rsq00);
271 /* Load parameters for j particles */
272 jq0 = _mm_load_sd(charge+jnrA+0);
273 isaj0 = _mm_load_sd(invsqrta+jnrA+0);
275 /**************************
276 * CALCULATE INTERACTIONS *
277 **************************/
279 r00 = _mm_mul_pd(rsq00,rinv00);
281 /* Compute parameters for interactions between i and j atoms */
282 qq00 = _mm_mul_pd(iq0,jq0);
284 /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
285 isaprod = _mm_mul_pd(isai0,isaj0);
286 gbqqfactor = _mm_xor_pd(signbit,_mm_mul_pd(qq00,_mm_mul_pd(isaprod,gbinvepsdiff)));
287 gbscale = _mm_mul_pd(isaprod,gbtabscale);
289 /* Calculate generalized born table index - this is a separate table from the normal one,
290 * but we use the same procedure by multiplying r with scale and truncating to integer.
292 rt = _mm_mul_pd(r00,gbscale);
293 gbitab = _mm_cvttpd_epi32(rt);
295 gbeps = _mm_frcz_pd(rt);
297 gbeps = _mm_sub_pd(rt,_mm_round_pd(rt, _MM_FROUND_FLOOR));
299 gbitab = _mm_slli_epi32(gbitab,2);
301 Y = _mm_load_pd( gbtab + _mm_extract_epi32(gbitab,0) );
302 F = _mm_setzero_pd();
303 GMX_MM_TRANSPOSE2_PD(Y,F);
304 G = _mm_load_pd( gbtab + _mm_extract_epi32(gbitab,0) +2);
305 H = _mm_setzero_pd();
306 GMX_MM_TRANSPOSE2_PD(G,H);
307 Fp = _mm_macc_pd(gbeps,_mm_macc_pd(gbeps,H,G),F);
308 VV = _mm_macc_pd(gbeps,Fp,Y);
309 vgb = _mm_mul_pd(gbqqfactor,VV);
311 twogbeps = _mm_add_pd(gbeps,gbeps);
312 FF = _mm_macc_pd(_mm_macc_pd(twogbeps,H,G),gbeps,Fp);
313 fgb = _mm_mul_pd(gbqqfactor,_mm_mul_pd(FF,gbscale));
314 dvdatmp = _mm_mul_pd(minushalf,_mm_macc_pd(fgb,r00,vgb));
315 dvdatmp = _mm_unpacklo_pd(dvdatmp,_mm_setzero_pd());
316 dvdasum = _mm_add_pd(dvdasum,dvdatmp);
317 gmx_mm_increment_1real_pd(dvda+jnrA,_mm_mul_pd(dvdatmp,_mm_mul_pd(isaj0,isaj0)));
318 velec = _mm_mul_pd(qq00,rinv00);
319 felec = _mm_mul_pd(_mm_msub_pd(velec,rinv00,fgb),rinv00);
321 /* Update potential sum for this i atom from the interaction with this j atom. */
322 velec = _mm_unpacklo_pd(velec,_mm_setzero_pd());
323 velecsum = _mm_add_pd(velecsum,velec);
324 vgb = _mm_unpacklo_pd(vgb,_mm_setzero_pd());
325 vgbsum = _mm_add_pd(vgbsum,vgb);
329 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
331 /* Update vectorial force */
332 fix0 = _mm_macc_pd(dx00,fscal,fix0);
333 fiy0 = _mm_macc_pd(dy00,fscal,fiy0);
334 fiz0 = _mm_macc_pd(dz00,fscal,fiz0);
336 gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,
337 _mm_mul_pd(dx00,fscal),
338 _mm_mul_pd(dy00,fscal),
339 _mm_mul_pd(dz00,fscal));
341 /* Inner loop uses 61 flops */
344 /* End of innermost loop */
346 gmx_mm_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
347 f+i_coord_offset,fshift+i_shift_offset);
350 /* Update potential energies */
351 gmx_mm_update_1pot_pd(velecsum,kernel_data->energygrp_elec+ggid);
352 gmx_mm_update_1pot_pd(vgbsum,kernel_data->energygrp_polarization+ggid);
353 dvdasum = _mm_mul_pd(dvdasum, _mm_mul_pd(isai0,isai0));
354 gmx_mm_update_1pot_pd(dvdasum,dvda+inr);
356 /* Increment number of inner iterations */
357 inneriter += j_index_end - j_index_start;
359 /* Outer loop uses 9 flops */
362 /* Increment number of outer iterations */
365 /* Update outer/inner flops */
367 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VF,outeriter*9 + inneriter*61);
370 * Gromacs nonbonded kernel: nb_kernel_ElecGB_VdwNone_GeomP1P1_F_avx_128_fma_double
371 * Electrostatics interaction: GeneralizedBorn
372 * VdW interaction: None
373 * Geometry: Particle-Particle
374 * Calculate force/pot: Force
377 nb_kernel_ElecGB_VdwNone_GeomP1P1_F_avx_128_fma_double
378 (t_nblist * gmx_restrict nlist,
379 rvec * gmx_restrict xx,
380 rvec * gmx_restrict ff,
381 t_forcerec * gmx_restrict fr,
382 t_mdatoms * gmx_restrict mdatoms,
383 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
384 t_nrnb * gmx_restrict nrnb)
386 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
387 * just 0 for non-waters.
388 * Suffixes A,B refer to j loop unrolling done with SSE double precision, e.g. for the two different
389 * jnr indices corresponding to data put in the four positions in the SIMD register.
391 int i_shift_offset,i_coord_offset,outeriter,inneriter;
392 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
394 int j_coord_offsetA,j_coord_offsetB;
395 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
397 real *shiftvec,*fshift,*x,*f;
398 __m128d tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
400 __m128d ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
401 int vdwjidx0A,vdwjidx0B;
402 __m128d jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
403 __m128d dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
404 __m128d velec,felec,velecsum,facel,crf,krf,krf2;
407 __m128d vgb,fgb,vgbsum,dvdasum,gbscale,gbtabscale,isaprod,gbqqfactor,gbinvepsdiff,dvdaj,gbeps,twogbeps,dvdatmp;
408 __m128d minushalf = _mm_set1_pd(-0.5);
409 real *invsqrta,*dvda,*gbtab;
411 __m128i ifour = _mm_set1_epi32(4);
412 __m128d rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF,twovfeps;
414 __m128d dummy_mask,cutoff_mask;
415 __m128d signbit = gmx_mm_castsi128_pd( _mm_set_epi32(0x80000000,0x00000000,0x80000000,0x00000000) );
416 __m128d one = _mm_set1_pd(1.0);
417 __m128d two = _mm_set1_pd(2.0);
423 jindex = nlist->jindex;
425 shiftidx = nlist->shift;
427 shiftvec = fr->shift_vec[0];
428 fshift = fr->fshift[0];
429 facel = _mm_set1_pd(fr->epsfac);
430 charge = mdatoms->chargeA;
432 invsqrta = fr->invsqrta;
434 gbtabscale = _mm_set1_pd(fr->gbtab.scale);
435 gbtab = fr->gbtab.data;
436 gbinvepsdiff = _mm_set1_pd((1.0/fr->epsilon_r) - (1.0/fr->gb_epsilon_solvent));
438 /* Avoid stupid compiler warnings */
446 /* Start outer loop over neighborlists */
447 for(iidx=0; iidx<nri; iidx++)
449 /* Load shift vector for this list */
450 i_shift_offset = DIM*shiftidx[iidx];
452 /* Load limits for loop over neighbors */
453 j_index_start = jindex[iidx];
454 j_index_end = jindex[iidx+1];
456 /* Get outer coordinate index */
458 i_coord_offset = DIM*inr;
460 /* Load i particle coords and add shift vector */
461 gmx_mm_load_shift_and_1rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
463 fix0 = _mm_setzero_pd();
464 fiy0 = _mm_setzero_pd();
465 fiz0 = _mm_setzero_pd();
467 /* Load parameters for i particles */
468 iq0 = _mm_mul_pd(facel,_mm_load1_pd(charge+inr+0));
469 isai0 = _mm_load1_pd(invsqrta+inr+0);
471 dvdasum = _mm_setzero_pd();
473 /* Start inner kernel loop */
474 for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
477 /* Get j neighbor index, and coordinate index */
480 j_coord_offsetA = DIM*jnrA;
481 j_coord_offsetB = DIM*jnrB;
483 /* load j atom coordinates */
484 gmx_mm_load_1rvec_2ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
487 /* Calculate displacement vector */
488 dx00 = _mm_sub_pd(ix0,jx0);
489 dy00 = _mm_sub_pd(iy0,jy0);
490 dz00 = _mm_sub_pd(iz0,jz0);
492 /* Calculate squared distance and things based on it */
493 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
495 rinv00 = gmx_mm_invsqrt_pd(rsq00);
497 /* Load parameters for j particles */
498 jq0 = gmx_mm_load_2real_swizzle_pd(charge+jnrA+0,charge+jnrB+0);
499 isaj0 = gmx_mm_load_2real_swizzle_pd(invsqrta+jnrA+0,invsqrta+jnrB+0);
501 /**************************
502 * CALCULATE INTERACTIONS *
503 **************************/
505 r00 = _mm_mul_pd(rsq00,rinv00);
507 /* Compute parameters for interactions between i and j atoms */
508 qq00 = _mm_mul_pd(iq0,jq0);
510 /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
511 isaprod = _mm_mul_pd(isai0,isaj0);
512 gbqqfactor = _mm_xor_pd(signbit,_mm_mul_pd(qq00,_mm_mul_pd(isaprod,gbinvepsdiff)));
513 gbscale = _mm_mul_pd(isaprod,gbtabscale);
515 /* Calculate generalized born table index - this is a separate table from the normal one,
516 * but we use the same procedure by multiplying r with scale and truncating to integer.
518 rt = _mm_mul_pd(r00,gbscale);
519 gbitab = _mm_cvttpd_epi32(rt);
521 gbeps = _mm_frcz_pd(rt);
523 gbeps = _mm_sub_pd(rt,_mm_round_pd(rt, _MM_FROUND_FLOOR));
525 gbitab = _mm_slli_epi32(gbitab,2);
527 Y = _mm_load_pd( gbtab + _mm_extract_epi32(gbitab,0) );
528 F = _mm_load_pd( gbtab + _mm_extract_epi32(gbitab,1) );
529 GMX_MM_TRANSPOSE2_PD(Y,F);
530 G = _mm_load_pd( gbtab + _mm_extract_epi32(gbitab,0) +2);
531 H = _mm_load_pd( gbtab + _mm_extract_epi32(gbitab,1) +2);
532 GMX_MM_TRANSPOSE2_PD(G,H);
533 Fp = _mm_macc_pd(gbeps,_mm_macc_pd(gbeps,H,G),F);
534 VV = _mm_macc_pd(gbeps,Fp,Y);
535 vgb = _mm_mul_pd(gbqqfactor,VV);
537 twogbeps = _mm_add_pd(gbeps,gbeps);
538 FF = _mm_macc_pd(_mm_macc_pd(twogbeps,H,G),gbeps,Fp);
539 fgb = _mm_mul_pd(gbqqfactor,_mm_mul_pd(FF,gbscale));
540 dvdatmp = _mm_mul_pd(minushalf,_mm_macc_pd(fgb,r00,vgb));
541 dvdasum = _mm_add_pd(dvdasum,dvdatmp);
542 gmx_mm_increment_2real_swizzle_pd(dvda+jnrA,dvda+jnrB,_mm_mul_pd(dvdatmp,_mm_mul_pd(isaj0,isaj0)));
543 velec = _mm_mul_pd(qq00,rinv00);
544 felec = _mm_mul_pd(_mm_msub_pd(velec,rinv00,fgb),rinv00);
548 /* Update vectorial force */
549 fix0 = _mm_macc_pd(dx00,fscal,fix0);
550 fiy0 = _mm_macc_pd(dy00,fscal,fiy0);
551 fiz0 = _mm_macc_pd(dz00,fscal,fiz0);
553 gmx_mm_decrement_1rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,
554 _mm_mul_pd(dx00,fscal),
555 _mm_mul_pd(dy00,fscal),
556 _mm_mul_pd(dz00,fscal));
558 /* Inner loop uses 59 flops */
565 j_coord_offsetA = DIM*jnrA;
567 /* load j atom coordinates */
568 gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
571 /* Calculate displacement vector */
572 dx00 = _mm_sub_pd(ix0,jx0);
573 dy00 = _mm_sub_pd(iy0,jy0);
574 dz00 = _mm_sub_pd(iz0,jz0);
576 /* Calculate squared distance and things based on it */
577 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
579 rinv00 = gmx_mm_invsqrt_pd(rsq00);
581 /* Load parameters for j particles */
582 jq0 = _mm_load_sd(charge+jnrA+0);
583 isaj0 = _mm_load_sd(invsqrta+jnrA+0);
585 /**************************
586 * CALCULATE INTERACTIONS *
587 **************************/
589 r00 = _mm_mul_pd(rsq00,rinv00);
591 /* Compute parameters for interactions between i and j atoms */
592 qq00 = _mm_mul_pd(iq0,jq0);
594 /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
595 isaprod = _mm_mul_pd(isai0,isaj0);
596 gbqqfactor = _mm_xor_pd(signbit,_mm_mul_pd(qq00,_mm_mul_pd(isaprod,gbinvepsdiff)));
597 gbscale = _mm_mul_pd(isaprod,gbtabscale);
599 /* Calculate generalized born table index - this is a separate table from the normal one,
600 * but we use the same procedure by multiplying r with scale and truncating to integer.
602 rt = _mm_mul_pd(r00,gbscale);
603 gbitab = _mm_cvttpd_epi32(rt);
605 gbeps = _mm_frcz_pd(rt);
607 gbeps = _mm_sub_pd(rt,_mm_round_pd(rt, _MM_FROUND_FLOOR));
609 gbitab = _mm_slli_epi32(gbitab,2);
611 Y = _mm_load_pd( gbtab + _mm_extract_epi32(gbitab,0) );
612 F = _mm_setzero_pd();
613 GMX_MM_TRANSPOSE2_PD(Y,F);
614 G = _mm_load_pd( gbtab + _mm_extract_epi32(gbitab,0) +2);
615 H = _mm_setzero_pd();
616 GMX_MM_TRANSPOSE2_PD(G,H);
617 Fp = _mm_macc_pd(gbeps,_mm_macc_pd(gbeps,H,G),F);
618 VV = _mm_macc_pd(gbeps,Fp,Y);
619 vgb = _mm_mul_pd(gbqqfactor,VV);
621 twogbeps = _mm_add_pd(gbeps,gbeps);
622 FF = _mm_macc_pd(_mm_macc_pd(twogbeps,H,G),gbeps,Fp);
623 fgb = _mm_mul_pd(gbqqfactor,_mm_mul_pd(FF,gbscale));
624 dvdatmp = _mm_mul_pd(minushalf,_mm_macc_pd(fgb,r00,vgb));
625 dvdatmp = _mm_unpacklo_pd(dvdatmp,_mm_setzero_pd());
626 dvdasum = _mm_add_pd(dvdasum,dvdatmp);
627 gmx_mm_increment_1real_pd(dvda+jnrA,_mm_mul_pd(dvdatmp,_mm_mul_pd(isaj0,isaj0)));
628 velec = _mm_mul_pd(qq00,rinv00);
629 felec = _mm_mul_pd(_mm_msub_pd(velec,rinv00,fgb),rinv00);
633 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
635 /* Update vectorial force */
636 fix0 = _mm_macc_pd(dx00,fscal,fix0);
637 fiy0 = _mm_macc_pd(dy00,fscal,fiy0);
638 fiz0 = _mm_macc_pd(dz00,fscal,fiz0);
640 gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,
641 _mm_mul_pd(dx00,fscal),
642 _mm_mul_pd(dy00,fscal),
643 _mm_mul_pd(dz00,fscal));
645 /* Inner loop uses 59 flops */
648 /* End of innermost loop */
650 gmx_mm_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
651 f+i_coord_offset,fshift+i_shift_offset);
653 dvdasum = _mm_mul_pd(dvdasum, _mm_mul_pd(isai0,isai0));
654 gmx_mm_update_1pot_pd(dvdasum,dvda+inr);
656 /* Increment number of inner iterations */
657 inneriter += j_index_end - j_index_start;
659 /* Outer loop uses 7 flops */
662 /* Increment number of outer iterations */
665 /* Update outer/inner flops */
667 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_F,outeriter*7 + inneriter*59);