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_VdwLJ_GeomP1P1_VF_avx_128_fma_double
38 * Electrostatics interaction: GeneralizedBorn
39 * VdW interaction: LennardJones
40 * Geometry: Particle-Particle
41 * Calculate force/pot: PotentialAndForce
44 nb_kernel_ElecGB_VdwLJ_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 __m128d rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
81 __m128d one_sixth = _mm_set1_pd(1.0/6.0);
82 __m128d one_twelfth = _mm_set1_pd(1.0/12.0);
84 __m128i ifour = _mm_set1_epi32(4);
85 __m128d rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF,twovfeps;
87 __m128d dummy_mask,cutoff_mask;
88 __m128d signbit = gmx_mm_castsi128_pd( _mm_set_epi32(0x80000000,0x00000000,0x80000000,0x00000000) );
89 __m128d one = _mm_set1_pd(1.0);
90 __m128d two = _mm_set1_pd(2.0);
96 jindex = nlist->jindex;
98 shiftidx = nlist->shift;
100 shiftvec = fr->shift_vec[0];
101 fshift = fr->fshift[0];
102 facel = _mm_set1_pd(fr->epsfac);
103 charge = mdatoms->chargeA;
104 nvdwtype = fr->ntype;
106 vdwtype = mdatoms->typeA;
108 invsqrta = fr->invsqrta;
110 gbtabscale = _mm_set1_pd(fr->gbtab.scale);
111 gbtab = fr->gbtab.data;
112 gbinvepsdiff = _mm_set1_pd((1.0/fr->epsilon_r) - (1.0/fr->gb_epsilon_solvent));
114 /* Avoid stupid compiler warnings */
122 /* Start outer loop over neighborlists */
123 for(iidx=0; iidx<nri; iidx++)
125 /* Load shift vector for this list */
126 i_shift_offset = DIM*shiftidx[iidx];
128 /* Load limits for loop over neighbors */
129 j_index_start = jindex[iidx];
130 j_index_end = jindex[iidx+1];
132 /* Get outer coordinate index */
134 i_coord_offset = DIM*inr;
136 /* Load i particle coords and add shift vector */
137 gmx_mm_load_shift_and_1rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
139 fix0 = _mm_setzero_pd();
140 fiy0 = _mm_setzero_pd();
141 fiz0 = _mm_setzero_pd();
143 /* Load parameters for i particles */
144 iq0 = _mm_mul_pd(facel,_mm_load1_pd(charge+inr+0));
145 isai0 = _mm_load1_pd(invsqrta+inr+0);
146 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
148 /* Reset potential sums */
149 velecsum = _mm_setzero_pd();
150 vgbsum = _mm_setzero_pd();
151 vvdwsum = _mm_setzero_pd();
152 dvdasum = _mm_setzero_pd();
154 /* Start inner kernel loop */
155 for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
158 /* Get j neighbor index, and coordinate index */
161 j_coord_offsetA = DIM*jnrA;
162 j_coord_offsetB = DIM*jnrB;
164 /* load j atom coordinates */
165 gmx_mm_load_1rvec_2ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
168 /* Calculate displacement vector */
169 dx00 = _mm_sub_pd(ix0,jx0);
170 dy00 = _mm_sub_pd(iy0,jy0);
171 dz00 = _mm_sub_pd(iz0,jz0);
173 /* Calculate squared distance and things based on it */
174 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
176 rinv00 = gmx_mm_invsqrt_pd(rsq00);
178 rinvsq00 = _mm_mul_pd(rinv00,rinv00);
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);
183 vdwjidx0A = 2*vdwtype[jnrA+0];
184 vdwjidx0B = 2*vdwtype[jnrB+0];
186 /**************************
187 * CALCULATE INTERACTIONS *
188 **************************/
190 r00 = _mm_mul_pd(rsq00,rinv00);
192 /* Compute parameters for interactions between i and j atoms */
193 qq00 = _mm_mul_pd(iq0,jq0);
194 gmx_mm_load_2pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,
195 vdwparam+vdwioffset0+vdwjidx0B,&c6_00,&c12_00);
197 /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
198 isaprod = _mm_mul_pd(isai0,isaj0);
199 gbqqfactor = _mm_xor_pd(signbit,_mm_mul_pd(qq00,_mm_mul_pd(isaprod,gbinvepsdiff)));
200 gbscale = _mm_mul_pd(isaprod,gbtabscale);
202 /* Calculate generalized born table index - this is a separate table from the normal one,
203 * but we use the same procedure by multiplying r with scale and truncating to integer.
205 rt = _mm_mul_pd(r00,gbscale);
206 gbitab = _mm_cvttpd_epi32(rt);
208 gbeps = _mm_frcz_pd(rt);
210 gbeps = _mm_sub_pd(rt,_mm_round_pd(rt, _MM_FROUND_FLOOR));
212 gbitab = _mm_slli_epi32(gbitab,2);
214 Y = _mm_load_pd( gbtab + _mm_extract_epi32(gbitab,0) );
215 F = _mm_load_pd( gbtab + _mm_extract_epi32(gbitab,1) );
216 GMX_MM_TRANSPOSE2_PD(Y,F);
217 G = _mm_load_pd( gbtab + _mm_extract_epi32(gbitab,0) +2);
218 H = _mm_load_pd( gbtab + _mm_extract_epi32(gbitab,1) +2);
219 GMX_MM_TRANSPOSE2_PD(G,H);
220 Fp = _mm_macc_pd(gbeps,_mm_macc_pd(gbeps,H,G),F);
221 VV = _mm_macc_pd(gbeps,Fp,Y);
222 vgb = _mm_mul_pd(gbqqfactor,VV);
224 twogbeps = _mm_add_pd(gbeps,gbeps);
225 FF = _mm_macc_pd(_mm_macc_pd(twogbeps,H,G),gbeps,Fp);
226 fgb = _mm_mul_pd(gbqqfactor,_mm_mul_pd(FF,gbscale));
227 dvdatmp = _mm_mul_pd(minushalf,_mm_macc_pd(fgb,r00,vgb));
228 dvdasum = _mm_add_pd(dvdasum,dvdatmp);
229 gmx_mm_increment_2real_swizzle_pd(dvda+jnrA,dvda+jnrB,_mm_mul_pd(dvdatmp,_mm_mul_pd(isaj0,isaj0)));
230 velec = _mm_mul_pd(qq00,rinv00);
231 felec = _mm_mul_pd(_mm_msub_pd(velec,rinv00,fgb),rinv00);
233 /* LENNARD-JONES DISPERSION/REPULSION */
235 rinvsix = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
236 vvdw6 = _mm_mul_pd(c6_00,rinvsix);
237 vvdw12 = _mm_mul_pd(c12_00,_mm_mul_pd(rinvsix,rinvsix));
238 vvdw = _mm_msub_pd( vvdw12,one_twelfth, _mm_mul_pd(vvdw6,one_sixth) );
239 fvdw = _mm_mul_pd(_mm_sub_pd(vvdw12,vvdw6),rinvsq00);
241 /* Update potential sum for this i atom from the interaction with this j atom. */
242 velecsum = _mm_add_pd(velecsum,velec);
243 vgbsum = _mm_add_pd(vgbsum,vgb);
244 vvdwsum = _mm_add_pd(vvdwsum,vvdw);
246 fscal = _mm_add_pd(felec,fvdw);
248 /* Update vectorial force */
249 fix0 = _mm_macc_pd(dx00,fscal,fix0);
250 fiy0 = _mm_macc_pd(dy00,fscal,fiy0);
251 fiz0 = _mm_macc_pd(dz00,fscal,fiz0);
253 gmx_mm_decrement_1rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,
254 _mm_mul_pd(dx00,fscal),
255 _mm_mul_pd(dy00,fscal),
256 _mm_mul_pd(dz00,fscal));
258 /* Inner loop uses 74 flops */
265 j_coord_offsetA = DIM*jnrA;
267 /* load j atom coordinates */
268 gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
271 /* Calculate displacement vector */
272 dx00 = _mm_sub_pd(ix0,jx0);
273 dy00 = _mm_sub_pd(iy0,jy0);
274 dz00 = _mm_sub_pd(iz0,jz0);
276 /* Calculate squared distance and things based on it */
277 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
279 rinv00 = gmx_mm_invsqrt_pd(rsq00);
281 rinvsq00 = _mm_mul_pd(rinv00,rinv00);
283 /* Load parameters for j particles */
284 jq0 = _mm_load_sd(charge+jnrA+0);
285 isaj0 = _mm_load_sd(invsqrta+jnrA+0);
286 vdwjidx0A = 2*vdwtype[jnrA+0];
288 /**************************
289 * CALCULATE INTERACTIONS *
290 **************************/
292 r00 = _mm_mul_pd(rsq00,rinv00);
294 /* Compute parameters for interactions between i and j atoms */
295 qq00 = _mm_mul_pd(iq0,jq0);
296 gmx_mm_load_1pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,&c6_00,&c12_00);
298 /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
299 isaprod = _mm_mul_pd(isai0,isaj0);
300 gbqqfactor = _mm_xor_pd(signbit,_mm_mul_pd(qq00,_mm_mul_pd(isaprod,gbinvepsdiff)));
301 gbscale = _mm_mul_pd(isaprod,gbtabscale);
303 /* Calculate generalized born table index - this is a separate table from the normal one,
304 * but we use the same procedure by multiplying r with scale and truncating to integer.
306 rt = _mm_mul_pd(r00,gbscale);
307 gbitab = _mm_cvttpd_epi32(rt);
309 gbeps = _mm_frcz_pd(rt);
311 gbeps = _mm_sub_pd(rt,_mm_round_pd(rt, _MM_FROUND_FLOOR));
313 gbitab = _mm_slli_epi32(gbitab,2);
315 Y = _mm_load_pd( gbtab + _mm_extract_epi32(gbitab,0) );
316 F = _mm_setzero_pd();
317 GMX_MM_TRANSPOSE2_PD(Y,F);
318 G = _mm_load_pd( gbtab + _mm_extract_epi32(gbitab,0) +2);
319 H = _mm_setzero_pd();
320 GMX_MM_TRANSPOSE2_PD(G,H);
321 Fp = _mm_macc_pd(gbeps,_mm_macc_pd(gbeps,H,G),F);
322 VV = _mm_macc_pd(gbeps,Fp,Y);
323 vgb = _mm_mul_pd(gbqqfactor,VV);
325 twogbeps = _mm_add_pd(gbeps,gbeps);
326 FF = _mm_macc_pd(_mm_macc_pd(twogbeps,H,G),gbeps,Fp);
327 fgb = _mm_mul_pd(gbqqfactor,_mm_mul_pd(FF,gbscale));
328 dvdatmp = _mm_mul_pd(minushalf,_mm_macc_pd(fgb,r00,vgb));
329 dvdatmp = _mm_unpacklo_pd(dvdatmp,_mm_setzero_pd());
330 dvdasum = _mm_add_pd(dvdasum,dvdatmp);
331 gmx_mm_increment_1real_pd(dvda+jnrA,_mm_mul_pd(dvdatmp,_mm_mul_pd(isaj0,isaj0)));
332 velec = _mm_mul_pd(qq00,rinv00);
333 felec = _mm_mul_pd(_mm_msub_pd(velec,rinv00,fgb),rinv00);
335 /* LENNARD-JONES DISPERSION/REPULSION */
337 rinvsix = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
338 vvdw6 = _mm_mul_pd(c6_00,rinvsix);
339 vvdw12 = _mm_mul_pd(c12_00,_mm_mul_pd(rinvsix,rinvsix));
340 vvdw = _mm_msub_pd( vvdw12,one_twelfth, _mm_mul_pd(vvdw6,one_sixth) );
341 fvdw = _mm_mul_pd(_mm_sub_pd(vvdw12,vvdw6),rinvsq00);
343 /* Update potential sum for this i atom from the interaction with this j atom. */
344 velec = _mm_unpacklo_pd(velec,_mm_setzero_pd());
345 velecsum = _mm_add_pd(velecsum,velec);
346 vgb = _mm_unpacklo_pd(vgb,_mm_setzero_pd());
347 vgbsum = _mm_add_pd(vgbsum,vgb);
348 vvdw = _mm_unpacklo_pd(vvdw,_mm_setzero_pd());
349 vvdwsum = _mm_add_pd(vvdwsum,vvdw);
351 fscal = _mm_add_pd(felec,fvdw);
353 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
355 /* Update vectorial force */
356 fix0 = _mm_macc_pd(dx00,fscal,fix0);
357 fiy0 = _mm_macc_pd(dy00,fscal,fiy0);
358 fiz0 = _mm_macc_pd(dz00,fscal,fiz0);
360 gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,
361 _mm_mul_pd(dx00,fscal),
362 _mm_mul_pd(dy00,fscal),
363 _mm_mul_pd(dz00,fscal));
365 /* Inner loop uses 74 flops */
368 /* End of innermost loop */
370 gmx_mm_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
371 f+i_coord_offset,fshift+i_shift_offset);
374 /* Update potential energies */
375 gmx_mm_update_1pot_pd(velecsum,kernel_data->energygrp_elec+ggid);
376 gmx_mm_update_1pot_pd(vgbsum,kernel_data->energygrp_polarization+ggid);
377 gmx_mm_update_1pot_pd(vvdwsum,kernel_data->energygrp_vdw+ggid);
378 dvdasum = _mm_mul_pd(dvdasum, _mm_mul_pd(isai0,isai0));
379 gmx_mm_update_1pot_pd(dvdasum,dvda+inr);
381 /* Increment number of inner iterations */
382 inneriter += j_index_end - j_index_start;
384 /* Outer loop uses 10 flops */
387 /* Increment number of outer iterations */
390 /* Update outer/inner flops */
392 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_VF,outeriter*10 + inneriter*74);
395 * Gromacs nonbonded kernel: nb_kernel_ElecGB_VdwLJ_GeomP1P1_F_avx_128_fma_double
396 * Electrostatics interaction: GeneralizedBorn
397 * VdW interaction: LennardJones
398 * Geometry: Particle-Particle
399 * Calculate force/pot: Force
402 nb_kernel_ElecGB_VdwLJ_GeomP1P1_F_avx_128_fma_double
403 (t_nblist * gmx_restrict nlist,
404 rvec * gmx_restrict xx,
405 rvec * gmx_restrict ff,
406 t_forcerec * gmx_restrict fr,
407 t_mdatoms * gmx_restrict mdatoms,
408 nb_kernel_data_t * gmx_restrict kernel_data,
409 t_nrnb * gmx_restrict nrnb)
411 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
412 * just 0 for non-waters.
413 * Suffixes A,B refer to j loop unrolling done with SSE double precision, e.g. for the two different
414 * jnr indices corresponding to data put in the four positions in the SIMD register.
416 int i_shift_offset,i_coord_offset,outeriter,inneriter;
417 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
419 int j_coord_offsetA,j_coord_offsetB;
420 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
422 real *shiftvec,*fshift,*x,*f;
423 __m128d tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
425 __m128d ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
426 int vdwjidx0A,vdwjidx0B;
427 __m128d jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
428 __m128d dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
429 __m128d velec,felec,velecsum,facel,crf,krf,krf2;
432 __m128d vgb,fgb,vgbsum,dvdasum,gbscale,gbtabscale,isaprod,gbqqfactor,gbinvepsdiff,dvdaj,gbeps,twogbeps,dvdatmp;
433 __m128d minushalf = _mm_set1_pd(-0.5);
434 real *invsqrta,*dvda,*gbtab;
436 __m128d rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
439 __m128d one_sixth = _mm_set1_pd(1.0/6.0);
440 __m128d one_twelfth = _mm_set1_pd(1.0/12.0);
442 __m128i ifour = _mm_set1_epi32(4);
443 __m128d rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF,twovfeps;
445 __m128d dummy_mask,cutoff_mask;
446 __m128d signbit = gmx_mm_castsi128_pd( _mm_set_epi32(0x80000000,0x00000000,0x80000000,0x00000000) );
447 __m128d one = _mm_set1_pd(1.0);
448 __m128d two = _mm_set1_pd(2.0);
454 jindex = nlist->jindex;
456 shiftidx = nlist->shift;
458 shiftvec = fr->shift_vec[0];
459 fshift = fr->fshift[0];
460 facel = _mm_set1_pd(fr->epsfac);
461 charge = mdatoms->chargeA;
462 nvdwtype = fr->ntype;
464 vdwtype = mdatoms->typeA;
466 invsqrta = fr->invsqrta;
468 gbtabscale = _mm_set1_pd(fr->gbtab.scale);
469 gbtab = fr->gbtab.data;
470 gbinvepsdiff = _mm_set1_pd((1.0/fr->epsilon_r) - (1.0/fr->gb_epsilon_solvent));
472 /* Avoid stupid compiler warnings */
480 /* Start outer loop over neighborlists */
481 for(iidx=0; iidx<nri; iidx++)
483 /* Load shift vector for this list */
484 i_shift_offset = DIM*shiftidx[iidx];
486 /* Load limits for loop over neighbors */
487 j_index_start = jindex[iidx];
488 j_index_end = jindex[iidx+1];
490 /* Get outer coordinate index */
492 i_coord_offset = DIM*inr;
494 /* Load i particle coords and add shift vector */
495 gmx_mm_load_shift_and_1rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
497 fix0 = _mm_setzero_pd();
498 fiy0 = _mm_setzero_pd();
499 fiz0 = _mm_setzero_pd();
501 /* Load parameters for i particles */
502 iq0 = _mm_mul_pd(facel,_mm_load1_pd(charge+inr+0));
503 isai0 = _mm_load1_pd(invsqrta+inr+0);
504 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
506 dvdasum = _mm_setzero_pd();
508 /* Start inner kernel loop */
509 for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
512 /* Get j neighbor index, and coordinate index */
515 j_coord_offsetA = DIM*jnrA;
516 j_coord_offsetB = DIM*jnrB;
518 /* load j atom coordinates */
519 gmx_mm_load_1rvec_2ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
522 /* Calculate displacement vector */
523 dx00 = _mm_sub_pd(ix0,jx0);
524 dy00 = _mm_sub_pd(iy0,jy0);
525 dz00 = _mm_sub_pd(iz0,jz0);
527 /* Calculate squared distance and things based on it */
528 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
530 rinv00 = gmx_mm_invsqrt_pd(rsq00);
532 rinvsq00 = _mm_mul_pd(rinv00,rinv00);
534 /* Load parameters for j particles */
535 jq0 = gmx_mm_load_2real_swizzle_pd(charge+jnrA+0,charge+jnrB+0);
536 isaj0 = gmx_mm_load_2real_swizzle_pd(invsqrta+jnrA+0,invsqrta+jnrB+0);
537 vdwjidx0A = 2*vdwtype[jnrA+0];
538 vdwjidx0B = 2*vdwtype[jnrB+0];
540 /**************************
541 * CALCULATE INTERACTIONS *
542 **************************/
544 r00 = _mm_mul_pd(rsq00,rinv00);
546 /* Compute parameters for interactions between i and j atoms */
547 qq00 = _mm_mul_pd(iq0,jq0);
548 gmx_mm_load_2pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,
549 vdwparam+vdwioffset0+vdwjidx0B,&c6_00,&c12_00);
551 /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
552 isaprod = _mm_mul_pd(isai0,isaj0);
553 gbqqfactor = _mm_xor_pd(signbit,_mm_mul_pd(qq00,_mm_mul_pd(isaprod,gbinvepsdiff)));
554 gbscale = _mm_mul_pd(isaprod,gbtabscale);
556 /* Calculate generalized born table index - this is a separate table from the normal one,
557 * but we use the same procedure by multiplying r with scale and truncating to integer.
559 rt = _mm_mul_pd(r00,gbscale);
560 gbitab = _mm_cvttpd_epi32(rt);
562 gbeps = _mm_frcz_pd(rt);
564 gbeps = _mm_sub_pd(rt,_mm_round_pd(rt, _MM_FROUND_FLOOR));
566 gbitab = _mm_slli_epi32(gbitab,2);
568 Y = _mm_load_pd( gbtab + _mm_extract_epi32(gbitab,0) );
569 F = _mm_load_pd( gbtab + _mm_extract_epi32(gbitab,1) );
570 GMX_MM_TRANSPOSE2_PD(Y,F);
571 G = _mm_load_pd( gbtab + _mm_extract_epi32(gbitab,0) +2);
572 H = _mm_load_pd( gbtab + _mm_extract_epi32(gbitab,1) +2);
573 GMX_MM_TRANSPOSE2_PD(G,H);
574 Fp = _mm_macc_pd(gbeps,_mm_macc_pd(gbeps,H,G),F);
575 VV = _mm_macc_pd(gbeps,Fp,Y);
576 vgb = _mm_mul_pd(gbqqfactor,VV);
578 twogbeps = _mm_add_pd(gbeps,gbeps);
579 FF = _mm_macc_pd(_mm_macc_pd(twogbeps,H,G),gbeps,Fp);
580 fgb = _mm_mul_pd(gbqqfactor,_mm_mul_pd(FF,gbscale));
581 dvdatmp = _mm_mul_pd(minushalf,_mm_macc_pd(fgb,r00,vgb));
582 dvdasum = _mm_add_pd(dvdasum,dvdatmp);
583 gmx_mm_increment_2real_swizzle_pd(dvda+jnrA,dvda+jnrB,_mm_mul_pd(dvdatmp,_mm_mul_pd(isaj0,isaj0)));
584 velec = _mm_mul_pd(qq00,rinv00);
585 felec = _mm_mul_pd(_mm_msub_pd(velec,rinv00,fgb),rinv00);
587 /* LENNARD-JONES DISPERSION/REPULSION */
589 rinvsix = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
590 fvdw = _mm_mul_pd(_mm_msub_pd(c12_00,rinvsix,c6_00),_mm_mul_pd(rinvsix,rinvsq00));
592 fscal = _mm_add_pd(felec,fvdw);
594 /* Update vectorial force */
595 fix0 = _mm_macc_pd(dx00,fscal,fix0);
596 fiy0 = _mm_macc_pd(dy00,fscal,fiy0);
597 fiz0 = _mm_macc_pd(dz00,fscal,fiz0);
599 gmx_mm_decrement_1rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,
600 _mm_mul_pd(dx00,fscal),
601 _mm_mul_pd(dy00,fscal),
602 _mm_mul_pd(dz00,fscal));
604 /* Inner loop uses 67 flops */
611 j_coord_offsetA = DIM*jnrA;
613 /* load j atom coordinates */
614 gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
617 /* Calculate displacement vector */
618 dx00 = _mm_sub_pd(ix0,jx0);
619 dy00 = _mm_sub_pd(iy0,jy0);
620 dz00 = _mm_sub_pd(iz0,jz0);
622 /* Calculate squared distance and things based on it */
623 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
625 rinv00 = gmx_mm_invsqrt_pd(rsq00);
627 rinvsq00 = _mm_mul_pd(rinv00,rinv00);
629 /* Load parameters for j particles */
630 jq0 = _mm_load_sd(charge+jnrA+0);
631 isaj0 = _mm_load_sd(invsqrta+jnrA+0);
632 vdwjidx0A = 2*vdwtype[jnrA+0];
634 /**************************
635 * CALCULATE INTERACTIONS *
636 **************************/
638 r00 = _mm_mul_pd(rsq00,rinv00);
640 /* Compute parameters for interactions between i and j atoms */
641 qq00 = _mm_mul_pd(iq0,jq0);
642 gmx_mm_load_1pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,&c6_00,&c12_00);
644 /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
645 isaprod = _mm_mul_pd(isai0,isaj0);
646 gbqqfactor = _mm_xor_pd(signbit,_mm_mul_pd(qq00,_mm_mul_pd(isaprod,gbinvepsdiff)));
647 gbscale = _mm_mul_pd(isaprod,gbtabscale);
649 /* Calculate generalized born table index - this is a separate table from the normal one,
650 * but we use the same procedure by multiplying r with scale and truncating to integer.
652 rt = _mm_mul_pd(r00,gbscale);
653 gbitab = _mm_cvttpd_epi32(rt);
655 gbeps = _mm_frcz_pd(rt);
657 gbeps = _mm_sub_pd(rt,_mm_round_pd(rt, _MM_FROUND_FLOOR));
659 gbitab = _mm_slli_epi32(gbitab,2);
661 Y = _mm_load_pd( gbtab + _mm_extract_epi32(gbitab,0) );
662 F = _mm_setzero_pd();
663 GMX_MM_TRANSPOSE2_PD(Y,F);
664 G = _mm_load_pd( gbtab + _mm_extract_epi32(gbitab,0) +2);
665 H = _mm_setzero_pd();
666 GMX_MM_TRANSPOSE2_PD(G,H);
667 Fp = _mm_macc_pd(gbeps,_mm_macc_pd(gbeps,H,G),F);
668 VV = _mm_macc_pd(gbeps,Fp,Y);
669 vgb = _mm_mul_pd(gbqqfactor,VV);
671 twogbeps = _mm_add_pd(gbeps,gbeps);
672 FF = _mm_macc_pd(_mm_macc_pd(twogbeps,H,G),gbeps,Fp);
673 fgb = _mm_mul_pd(gbqqfactor,_mm_mul_pd(FF,gbscale));
674 dvdatmp = _mm_mul_pd(minushalf,_mm_macc_pd(fgb,r00,vgb));
675 dvdatmp = _mm_unpacklo_pd(dvdatmp,_mm_setzero_pd());
676 dvdasum = _mm_add_pd(dvdasum,dvdatmp);
677 gmx_mm_increment_1real_pd(dvda+jnrA,_mm_mul_pd(dvdatmp,_mm_mul_pd(isaj0,isaj0)));
678 velec = _mm_mul_pd(qq00,rinv00);
679 felec = _mm_mul_pd(_mm_msub_pd(velec,rinv00,fgb),rinv00);
681 /* LENNARD-JONES DISPERSION/REPULSION */
683 rinvsix = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
684 fvdw = _mm_mul_pd(_mm_msub_pd(c12_00,rinvsix,c6_00),_mm_mul_pd(rinvsix,rinvsq00));
686 fscal = _mm_add_pd(felec,fvdw);
688 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
690 /* Update vectorial force */
691 fix0 = _mm_macc_pd(dx00,fscal,fix0);
692 fiy0 = _mm_macc_pd(dy00,fscal,fiy0);
693 fiz0 = _mm_macc_pd(dz00,fscal,fiz0);
695 gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,
696 _mm_mul_pd(dx00,fscal),
697 _mm_mul_pd(dy00,fscal),
698 _mm_mul_pd(dz00,fscal));
700 /* Inner loop uses 67 flops */
703 /* End of innermost loop */
705 gmx_mm_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
706 f+i_coord_offset,fshift+i_shift_offset);
708 dvdasum = _mm_mul_pd(dvdasum, _mm_mul_pd(isai0,isai0));
709 gmx_mm_update_1pot_pd(dvdasum,dvda+inr);
711 /* Increment number of inner iterations */
712 inneriter += j_index_end - j_index_start;
714 /* Outer loop uses 7 flops */
717 /* Increment number of outer iterations */
720 /* Update outer/inner flops */
722 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_F,outeriter*7 + inneriter*67);