2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015,2017, 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/gmxlib/nrnb.h"
47 #include "kernelutil_x86_avx_128_fma_double.h"
50 * Gromacs nonbonded kernel: nb_kernel_ElecGB_VdwNone_GeomP1P1_VF_avx_128_fma_double
51 * Electrostatics interaction: GeneralizedBorn
52 * VdW interaction: None
53 * Geometry: Particle-Particle
54 * Calculate force/pot: PotentialAndForce
57 nb_kernel_ElecGB_VdwNone_GeomP1P1_VF_avx_128_fma_double
58 (t_nblist * gmx_restrict nlist,
59 rvec * gmx_restrict xx,
60 rvec * gmx_restrict ff,
61 struct t_forcerec * gmx_restrict fr,
62 t_mdatoms * gmx_restrict mdatoms,
63 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
64 t_nrnb * gmx_restrict nrnb)
66 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
67 * just 0 for non-waters.
68 * Suffixes A,B refer to j loop unrolling done with SSE double precision, e.g. for the two different
69 * jnr indices corresponding to data put in the four positions in the SIMD register.
71 int i_shift_offset,i_coord_offset,outeriter,inneriter;
72 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
74 int j_coord_offsetA,j_coord_offsetB;
75 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
77 real *shiftvec,*fshift,*x,*f;
78 __m128d tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
80 __m128d ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
81 int vdwjidx0A,vdwjidx0B;
82 __m128d jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
83 __m128d dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
84 __m128d velec,felec,velecsum,facel,crf,krf,krf2;
87 __m128d vgb,fgb,vgbsum,dvdasum,gbscale,gbtabscale,isaprod,gbqqfactor,gbinvepsdiff,dvdaj,gbeps,twogbeps,dvdatmp;
88 __m128d minushalf = _mm_set1_pd(-0.5);
89 real *invsqrta,*dvda,*gbtab;
91 __m128i ifour = _mm_set1_epi32(4);
92 __m128d rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF,twovfeps;
94 __m128d dummy_mask,cutoff_mask;
95 __m128d signbit = gmx_mm_castsi128_pd( _mm_set_epi32(0x80000000,0x00000000,0x80000000,0x00000000) );
96 __m128d one = _mm_set1_pd(1.0);
97 __m128d two = _mm_set1_pd(2.0);
103 jindex = nlist->jindex;
105 shiftidx = nlist->shift;
107 shiftvec = fr->shift_vec[0];
108 fshift = fr->fshift[0];
109 facel = _mm_set1_pd(fr->ic->epsfac);
110 charge = mdatoms->chargeA;
112 invsqrta = fr->invsqrta;
114 gbtabscale = _mm_set1_pd(fr->gbtab->scale);
115 gbtab = fr->gbtab->data;
116 gbinvepsdiff = _mm_set1_pd((1.0/fr->ic->epsilon_r) - (1.0/fr->gb_epsilon_solvent));
118 /* Avoid stupid compiler warnings */
126 /* Start outer loop over neighborlists */
127 for(iidx=0; iidx<nri; iidx++)
129 /* Load shift vector for this list */
130 i_shift_offset = DIM*shiftidx[iidx];
132 /* Load limits for loop over neighbors */
133 j_index_start = jindex[iidx];
134 j_index_end = jindex[iidx+1];
136 /* Get outer coordinate index */
138 i_coord_offset = DIM*inr;
140 /* Load i particle coords and add shift vector */
141 gmx_mm_load_shift_and_1rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
143 fix0 = _mm_setzero_pd();
144 fiy0 = _mm_setzero_pd();
145 fiz0 = _mm_setzero_pd();
147 /* Load parameters for i particles */
148 iq0 = _mm_mul_pd(facel,_mm_load1_pd(charge+inr+0));
149 isai0 = _mm_load1_pd(invsqrta+inr+0);
151 /* Reset potential sums */
152 velecsum = _mm_setzero_pd();
153 vgbsum = _mm_setzero_pd();
154 dvdasum = _mm_setzero_pd();
156 /* Start inner kernel loop */
157 for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
160 /* Get j neighbor index, and coordinate index */
163 j_coord_offsetA = DIM*jnrA;
164 j_coord_offsetB = DIM*jnrB;
166 /* load j atom coordinates */
167 gmx_mm_load_1rvec_2ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
170 /* Calculate displacement vector */
171 dx00 = _mm_sub_pd(ix0,jx0);
172 dy00 = _mm_sub_pd(iy0,jy0);
173 dz00 = _mm_sub_pd(iz0,jz0);
175 /* Calculate squared distance and things based on it */
176 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
178 rinv00 = avx128fma_invsqrt_d(rsq00);
180 /* Load parameters for j particles */
181 jq0 = gmx_mm_load_2real_swizzle_pd(charge+jnrA+0,charge+jnrB+0);
182 isaj0 = gmx_mm_load_2real_swizzle_pd(invsqrta+jnrA+0,invsqrta+jnrB+0);
184 /**************************
185 * CALCULATE INTERACTIONS *
186 **************************/
188 r00 = _mm_mul_pd(rsq00,rinv00);
190 /* Compute parameters for interactions between i and j atoms */
191 qq00 = _mm_mul_pd(iq0,jq0);
193 /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
194 isaprod = _mm_mul_pd(isai0,isaj0);
195 gbqqfactor = _mm_xor_pd(signbit,_mm_mul_pd(qq00,_mm_mul_pd(isaprod,gbinvepsdiff)));
196 gbscale = _mm_mul_pd(isaprod,gbtabscale);
198 /* Calculate generalized born table index - this is a separate table from the normal one,
199 * but we use the same procedure by multiplying r with scale and truncating to integer.
201 rt = _mm_mul_pd(r00,gbscale);
202 gbitab = _mm_cvttpd_epi32(rt);
204 gbeps = _mm_frcz_pd(rt);
206 gbeps = _mm_sub_pd(rt,_mm_round_pd(rt, _MM_FROUND_FLOOR));
208 gbitab = _mm_slli_epi32(gbitab,2);
210 Y = _mm_load_pd( gbtab + _mm_extract_epi32(gbitab,0) );
211 F = _mm_load_pd( gbtab + _mm_extract_epi32(gbitab,1) );
212 GMX_MM_TRANSPOSE2_PD(Y,F);
213 G = _mm_load_pd( gbtab + _mm_extract_epi32(gbitab,0) +2);
214 H = _mm_load_pd( gbtab + _mm_extract_epi32(gbitab,1) +2);
215 GMX_MM_TRANSPOSE2_PD(G,H);
216 Fp = _mm_macc_pd(gbeps,_mm_macc_pd(gbeps,H,G),F);
217 VV = _mm_macc_pd(gbeps,Fp,Y);
218 vgb = _mm_mul_pd(gbqqfactor,VV);
220 twogbeps = _mm_add_pd(gbeps,gbeps);
221 FF = _mm_macc_pd(_mm_macc_pd(twogbeps,H,G),gbeps,Fp);
222 fgb = _mm_mul_pd(gbqqfactor,_mm_mul_pd(FF,gbscale));
223 dvdatmp = _mm_mul_pd(minushalf,_mm_macc_pd(fgb,r00,vgb));
224 dvdasum = _mm_add_pd(dvdasum,dvdatmp);
225 gmx_mm_increment_2real_swizzle_pd(dvda+jnrA,dvda+jnrB,_mm_mul_pd(dvdatmp,_mm_mul_pd(isaj0,isaj0)));
226 velec = _mm_mul_pd(qq00,rinv00);
227 felec = _mm_mul_pd(_mm_msub_pd(velec,rinv00,fgb),rinv00);
229 /* Update potential sum for this i atom from the interaction with this j atom. */
230 velecsum = _mm_add_pd(velecsum,velec);
231 vgbsum = _mm_add_pd(vgbsum,vgb);
235 /* Update vectorial force */
236 fix0 = _mm_macc_pd(dx00,fscal,fix0);
237 fiy0 = _mm_macc_pd(dy00,fscal,fiy0);
238 fiz0 = _mm_macc_pd(dz00,fscal,fiz0);
240 gmx_mm_decrement_1rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,
241 _mm_mul_pd(dx00,fscal),
242 _mm_mul_pd(dy00,fscal),
243 _mm_mul_pd(dz00,fscal));
245 /* Inner loop uses 61 flops */
252 j_coord_offsetA = DIM*jnrA;
254 /* load j atom coordinates */
255 gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
258 /* Calculate displacement vector */
259 dx00 = _mm_sub_pd(ix0,jx0);
260 dy00 = _mm_sub_pd(iy0,jy0);
261 dz00 = _mm_sub_pd(iz0,jz0);
263 /* Calculate squared distance and things based on it */
264 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
266 rinv00 = avx128fma_invsqrt_d(rsq00);
268 /* Load parameters for j particles */
269 jq0 = _mm_load_sd(charge+jnrA+0);
270 isaj0 = _mm_load_sd(invsqrta+jnrA+0);
272 /**************************
273 * CALCULATE INTERACTIONS *
274 **************************/
276 r00 = _mm_mul_pd(rsq00,rinv00);
278 /* Compute parameters for interactions between i and j atoms */
279 qq00 = _mm_mul_pd(iq0,jq0);
281 /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
282 isaprod = _mm_mul_pd(isai0,isaj0);
283 gbqqfactor = _mm_xor_pd(signbit,_mm_mul_pd(qq00,_mm_mul_pd(isaprod,gbinvepsdiff)));
284 gbscale = _mm_mul_pd(isaprod,gbtabscale);
286 /* Calculate generalized born table index - this is a separate table from the normal one,
287 * but we use the same procedure by multiplying r with scale and truncating to integer.
289 rt = _mm_mul_pd(r00,gbscale);
290 gbitab = _mm_cvttpd_epi32(rt);
292 gbeps = _mm_frcz_pd(rt);
294 gbeps = _mm_sub_pd(rt,_mm_round_pd(rt, _MM_FROUND_FLOOR));
296 gbitab = _mm_slli_epi32(gbitab,2);
298 Y = _mm_load_pd( gbtab + _mm_extract_epi32(gbitab,0) );
299 F = _mm_setzero_pd();
300 GMX_MM_TRANSPOSE2_PD(Y,F);
301 G = _mm_load_pd( gbtab + _mm_extract_epi32(gbitab,0) +2);
302 H = _mm_setzero_pd();
303 GMX_MM_TRANSPOSE2_PD(G,H);
304 Fp = _mm_macc_pd(gbeps,_mm_macc_pd(gbeps,H,G),F);
305 VV = _mm_macc_pd(gbeps,Fp,Y);
306 vgb = _mm_mul_pd(gbqqfactor,VV);
308 twogbeps = _mm_add_pd(gbeps,gbeps);
309 FF = _mm_macc_pd(_mm_macc_pd(twogbeps,H,G),gbeps,Fp);
310 fgb = _mm_mul_pd(gbqqfactor,_mm_mul_pd(FF,gbscale));
311 dvdatmp = _mm_mul_pd(minushalf,_mm_macc_pd(fgb,r00,vgb));
312 dvdatmp = _mm_unpacklo_pd(dvdatmp,_mm_setzero_pd());
313 dvdasum = _mm_add_pd(dvdasum,dvdatmp);
314 gmx_mm_increment_1real_pd(dvda+jnrA,_mm_mul_pd(dvdatmp,_mm_mul_pd(isaj0,isaj0)));
315 velec = _mm_mul_pd(qq00,rinv00);
316 felec = _mm_mul_pd(_mm_msub_pd(velec,rinv00,fgb),rinv00);
318 /* Update potential sum for this i atom from the interaction with this j atom. */
319 velec = _mm_unpacklo_pd(velec,_mm_setzero_pd());
320 velecsum = _mm_add_pd(velecsum,velec);
321 vgb = _mm_unpacklo_pd(vgb,_mm_setzero_pd());
322 vgbsum = _mm_add_pd(vgbsum,vgb);
326 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
328 /* Update vectorial force */
329 fix0 = _mm_macc_pd(dx00,fscal,fix0);
330 fiy0 = _mm_macc_pd(dy00,fscal,fiy0);
331 fiz0 = _mm_macc_pd(dz00,fscal,fiz0);
333 gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,
334 _mm_mul_pd(dx00,fscal),
335 _mm_mul_pd(dy00,fscal),
336 _mm_mul_pd(dz00,fscal));
338 /* Inner loop uses 61 flops */
341 /* End of innermost loop */
343 gmx_mm_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
344 f+i_coord_offset,fshift+i_shift_offset);
347 /* Update potential energies */
348 gmx_mm_update_1pot_pd(velecsum,kernel_data->energygrp_elec+ggid);
349 gmx_mm_update_1pot_pd(vgbsum,kernel_data->energygrp_polarization+ggid);
350 dvdasum = _mm_mul_pd(dvdasum, _mm_mul_pd(isai0,isai0));
351 gmx_mm_update_1pot_pd(dvdasum,dvda+inr);
353 /* Increment number of inner iterations */
354 inneriter += j_index_end - j_index_start;
356 /* Outer loop uses 9 flops */
359 /* Increment number of outer iterations */
362 /* Update outer/inner flops */
364 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VF,outeriter*9 + inneriter*61);
367 * Gromacs nonbonded kernel: nb_kernel_ElecGB_VdwNone_GeomP1P1_F_avx_128_fma_double
368 * Electrostatics interaction: GeneralizedBorn
369 * VdW interaction: None
370 * Geometry: Particle-Particle
371 * Calculate force/pot: Force
374 nb_kernel_ElecGB_VdwNone_GeomP1P1_F_avx_128_fma_double
375 (t_nblist * gmx_restrict nlist,
376 rvec * gmx_restrict xx,
377 rvec * gmx_restrict ff,
378 struct t_forcerec * gmx_restrict fr,
379 t_mdatoms * gmx_restrict mdatoms,
380 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
381 t_nrnb * gmx_restrict nrnb)
383 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
384 * just 0 for non-waters.
385 * Suffixes A,B refer to j loop unrolling done with SSE double precision, e.g. for the two different
386 * jnr indices corresponding to data put in the four positions in the SIMD register.
388 int i_shift_offset,i_coord_offset,outeriter,inneriter;
389 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
391 int j_coord_offsetA,j_coord_offsetB;
392 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
394 real *shiftvec,*fshift,*x,*f;
395 __m128d tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
397 __m128d ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
398 int vdwjidx0A,vdwjidx0B;
399 __m128d jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
400 __m128d dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
401 __m128d velec,felec,velecsum,facel,crf,krf,krf2;
404 __m128d vgb,fgb,vgbsum,dvdasum,gbscale,gbtabscale,isaprod,gbqqfactor,gbinvepsdiff,dvdaj,gbeps,twogbeps,dvdatmp;
405 __m128d minushalf = _mm_set1_pd(-0.5);
406 real *invsqrta,*dvda,*gbtab;
408 __m128i ifour = _mm_set1_epi32(4);
409 __m128d rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF,twovfeps;
411 __m128d dummy_mask,cutoff_mask;
412 __m128d signbit = gmx_mm_castsi128_pd( _mm_set_epi32(0x80000000,0x00000000,0x80000000,0x00000000) );
413 __m128d one = _mm_set1_pd(1.0);
414 __m128d two = _mm_set1_pd(2.0);
420 jindex = nlist->jindex;
422 shiftidx = nlist->shift;
424 shiftvec = fr->shift_vec[0];
425 fshift = fr->fshift[0];
426 facel = _mm_set1_pd(fr->ic->epsfac);
427 charge = mdatoms->chargeA;
429 invsqrta = fr->invsqrta;
431 gbtabscale = _mm_set1_pd(fr->gbtab->scale);
432 gbtab = fr->gbtab->data;
433 gbinvepsdiff = _mm_set1_pd((1.0/fr->ic->epsilon_r) - (1.0/fr->gb_epsilon_solvent));
435 /* Avoid stupid compiler warnings */
443 /* Start outer loop over neighborlists */
444 for(iidx=0; iidx<nri; iidx++)
446 /* Load shift vector for this list */
447 i_shift_offset = DIM*shiftidx[iidx];
449 /* Load limits for loop over neighbors */
450 j_index_start = jindex[iidx];
451 j_index_end = jindex[iidx+1];
453 /* Get outer coordinate index */
455 i_coord_offset = DIM*inr;
457 /* Load i particle coords and add shift vector */
458 gmx_mm_load_shift_and_1rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
460 fix0 = _mm_setzero_pd();
461 fiy0 = _mm_setzero_pd();
462 fiz0 = _mm_setzero_pd();
464 /* Load parameters for i particles */
465 iq0 = _mm_mul_pd(facel,_mm_load1_pd(charge+inr+0));
466 isai0 = _mm_load1_pd(invsqrta+inr+0);
468 dvdasum = _mm_setzero_pd();
470 /* Start inner kernel loop */
471 for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
474 /* Get j neighbor index, and coordinate index */
477 j_coord_offsetA = DIM*jnrA;
478 j_coord_offsetB = DIM*jnrB;
480 /* load j atom coordinates */
481 gmx_mm_load_1rvec_2ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
484 /* Calculate displacement vector */
485 dx00 = _mm_sub_pd(ix0,jx0);
486 dy00 = _mm_sub_pd(iy0,jy0);
487 dz00 = _mm_sub_pd(iz0,jz0);
489 /* Calculate squared distance and things based on it */
490 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
492 rinv00 = avx128fma_invsqrt_d(rsq00);
494 /* Load parameters for j particles */
495 jq0 = gmx_mm_load_2real_swizzle_pd(charge+jnrA+0,charge+jnrB+0);
496 isaj0 = gmx_mm_load_2real_swizzle_pd(invsqrta+jnrA+0,invsqrta+jnrB+0);
498 /**************************
499 * CALCULATE INTERACTIONS *
500 **************************/
502 r00 = _mm_mul_pd(rsq00,rinv00);
504 /* Compute parameters for interactions between i and j atoms */
505 qq00 = _mm_mul_pd(iq0,jq0);
507 /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
508 isaprod = _mm_mul_pd(isai0,isaj0);
509 gbqqfactor = _mm_xor_pd(signbit,_mm_mul_pd(qq00,_mm_mul_pd(isaprod,gbinvepsdiff)));
510 gbscale = _mm_mul_pd(isaprod,gbtabscale);
512 /* Calculate generalized born table index - this is a separate table from the normal one,
513 * but we use the same procedure by multiplying r with scale and truncating to integer.
515 rt = _mm_mul_pd(r00,gbscale);
516 gbitab = _mm_cvttpd_epi32(rt);
518 gbeps = _mm_frcz_pd(rt);
520 gbeps = _mm_sub_pd(rt,_mm_round_pd(rt, _MM_FROUND_FLOOR));
522 gbitab = _mm_slli_epi32(gbitab,2);
524 Y = _mm_load_pd( gbtab + _mm_extract_epi32(gbitab,0) );
525 F = _mm_load_pd( gbtab + _mm_extract_epi32(gbitab,1) );
526 GMX_MM_TRANSPOSE2_PD(Y,F);
527 G = _mm_load_pd( gbtab + _mm_extract_epi32(gbitab,0) +2);
528 H = _mm_load_pd( gbtab + _mm_extract_epi32(gbitab,1) +2);
529 GMX_MM_TRANSPOSE2_PD(G,H);
530 Fp = _mm_macc_pd(gbeps,_mm_macc_pd(gbeps,H,G),F);
531 VV = _mm_macc_pd(gbeps,Fp,Y);
532 vgb = _mm_mul_pd(gbqqfactor,VV);
534 twogbeps = _mm_add_pd(gbeps,gbeps);
535 FF = _mm_macc_pd(_mm_macc_pd(twogbeps,H,G),gbeps,Fp);
536 fgb = _mm_mul_pd(gbqqfactor,_mm_mul_pd(FF,gbscale));
537 dvdatmp = _mm_mul_pd(minushalf,_mm_macc_pd(fgb,r00,vgb));
538 dvdasum = _mm_add_pd(dvdasum,dvdatmp);
539 gmx_mm_increment_2real_swizzle_pd(dvda+jnrA,dvda+jnrB,_mm_mul_pd(dvdatmp,_mm_mul_pd(isaj0,isaj0)));
540 velec = _mm_mul_pd(qq00,rinv00);
541 felec = _mm_mul_pd(_mm_msub_pd(velec,rinv00,fgb),rinv00);
545 /* Update vectorial force */
546 fix0 = _mm_macc_pd(dx00,fscal,fix0);
547 fiy0 = _mm_macc_pd(dy00,fscal,fiy0);
548 fiz0 = _mm_macc_pd(dz00,fscal,fiz0);
550 gmx_mm_decrement_1rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,
551 _mm_mul_pd(dx00,fscal),
552 _mm_mul_pd(dy00,fscal),
553 _mm_mul_pd(dz00,fscal));
555 /* Inner loop uses 59 flops */
562 j_coord_offsetA = DIM*jnrA;
564 /* load j atom coordinates */
565 gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
568 /* Calculate displacement vector */
569 dx00 = _mm_sub_pd(ix0,jx0);
570 dy00 = _mm_sub_pd(iy0,jy0);
571 dz00 = _mm_sub_pd(iz0,jz0);
573 /* Calculate squared distance and things based on it */
574 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
576 rinv00 = avx128fma_invsqrt_d(rsq00);
578 /* Load parameters for j particles */
579 jq0 = _mm_load_sd(charge+jnrA+0);
580 isaj0 = _mm_load_sd(invsqrta+jnrA+0);
582 /**************************
583 * CALCULATE INTERACTIONS *
584 **************************/
586 r00 = _mm_mul_pd(rsq00,rinv00);
588 /* Compute parameters for interactions between i and j atoms */
589 qq00 = _mm_mul_pd(iq0,jq0);
591 /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
592 isaprod = _mm_mul_pd(isai0,isaj0);
593 gbqqfactor = _mm_xor_pd(signbit,_mm_mul_pd(qq00,_mm_mul_pd(isaprod,gbinvepsdiff)));
594 gbscale = _mm_mul_pd(isaprod,gbtabscale);
596 /* Calculate generalized born table index - this is a separate table from the normal one,
597 * but we use the same procedure by multiplying r with scale and truncating to integer.
599 rt = _mm_mul_pd(r00,gbscale);
600 gbitab = _mm_cvttpd_epi32(rt);
602 gbeps = _mm_frcz_pd(rt);
604 gbeps = _mm_sub_pd(rt,_mm_round_pd(rt, _MM_FROUND_FLOOR));
606 gbitab = _mm_slli_epi32(gbitab,2);
608 Y = _mm_load_pd( gbtab + _mm_extract_epi32(gbitab,0) );
609 F = _mm_setzero_pd();
610 GMX_MM_TRANSPOSE2_PD(Y,F);
611 G = _mm_load_pd( gbtab + _mm_extract_epi32(gbitab,0) +2);
612 H = _mm_setzero_pd();
613 GMX_MM_TRANSPOSE2_PD(G,H);
614 Fp = _mm_macc_pd(gbeps,_mm_macc_pd(gbeps,H,G),F);
615 VV = _mm_macc_pd(gbeps,Fp,Y);
616 vgb = _mm_mul_pd(gbqqfactor,VV);
618 twogbeps = _mm_add_pd(gbeps,gbeps);
619 FF = _mm_macc_pd(_mm_macc_pd(twogbeps,H,G),gbeps,Fp);
620 fgb = _mm_mul_pd(gbqqfactor,_mm_mul_pd(FF,gbscale));
621 dvdatmp = _mm_mul_pd(minushalf,_mm_macc_pd(fgb,r00,vgb));
622 dvdatmp = _mm_unpacklo_pd(dvdatmp,_mm_setzero_pd());
623 dvdasum = _mm_add_pd(dvdasum,dvdatmp);
624 gmx_mm_increment_1real_pd(dvda+jnrA,_mm_mul_pd(dvdatmp,_mm_mul_pd(isaj0,isaj0)));
625 velec = _mm_mul_pd(qq00,rinv00);
626 felec = _mm_mul_pd(_mm_msub_pd(velec,rinv00,fgb),rinv00);
630 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
632 /* Update vectorial force */
633 fix0 = _mm_macc_pd(dx00,fscal,fix0);
634 fiy0 = _mm_macc_pd(dy00,fscal,fiy0);
635 fiz0 = _mm_macc_pd(dz00,fscal,fiz0);
637 gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,
638 _mm_mul_pd(dx00,fscal),
639 _mm_mul_pd(dy00,fscal),
640 _mm_mul_pd(dz00,fscal));
642 /* Inner loop uses 59 flops */
645 /* End of innermost loop */
647 gmx_mm_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
648 f+i_coord_offset,fshift+i_shift_offset);
650 dvdasum = _mm_mul_pd(dvdasum, _mm_mul_pd(isai0,isai0));
651 gmx_mm_update_1pot_pd(dvdasum,dvda+inr);
653 /* Increment number of inner iterations */
654 inneriter += j_index_end - j_index_start;
656 /* Outer loop uses 7 flops */
659 /* Increment number of outer iterations */
662 /* Update outer/inner flops */
664 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_F,outeriter*7 + inneriter*59);