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_VdwLJ_GeomP1P1_VF_sse2_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_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 __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;
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);
207 gbeps = _mm_sub_pd(rt,_mm_cvtepi32_pd(gbitab));
208 gbitab = _mm_slli_epi32(gbitab,2);
210 Y = _mm_load_pd( gbtab + gmx_mm_extract_epi32(gbitab,0) );
211 F = _mm_load_pd( gbtab + gmx_mm_extract_epi32(gbitab,1) );
212 GMX_MM_TRANSPOSE2_PD(Y,F);
213 G = _mm_load_pd( gbtab + gmx_mm_extract_epi32(gbitab,0) +2);
214 H = _mm_load_pd( gbtab + gmx_mm_extract_epi32(gbitab,1) +2);
215 GMX_MM_TRANSPOSE2_PD(G,H);
216 Heps = _mm_mul_pd(gbeps,H);
217 Fp = _mm_add_pd(F,_mm_mul_pd(gbeps,_mm_add_pd(G,Heps)));
218 VV = _mm_add_pd(Y,_mm_mul_pd(gbeps,Fp));
219 vgb = _mm_mul_pd(gbqqfactor,VV);
221 FF = _mm_add_pd(Fp,_mm_mul_pd(gbeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
222 fgb = _mm_mul_pd(gbqqfactor,_mm_mul_pd(FF,gbscale));
223 dvdatmp = _mm_mul_pd(minushalf,_mm_add_pd(vgb,_mm_mul_pd(fgb,r00)));
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_sub_pd(_mm_mul_pd(velec,rinv00),fgb),rinv00);
229 /* LENNARD-JONES DISPERSION/REPULSION */
231 rinvsix = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
232 vvdw6 = _mm_mul_pd(c6_00,rinvsix);
233 vvdw12 = _mm_mul_pd(c12_00,_mm_mul_pd(rinvsix,rinvsix));
234 vvdw = _mm_sub_pd( _mm_mul_pd(vvdw12,one_twelfth) , _mm_mul_pd(vvdw6,one_sixth) );
235 fvdw = _mm_mul_pd(_mm_sub_pd(vvdw12,vvdw6),rinvsq00);
237 /* Update potential sum for this i atom from the interaction with this j atom. */
238 velecsum = _mm_add_pd(velecsum,velec);
239 vgbsum = _mm_add_pd(vgbsum,vgb);
240 vvdwsum = _mm_add_pd(vvdwsum,vvdw);
242 fscal = _mm_add_pd(felec,fvdw);
244 /* Calculate temporary vectorial force */
245 tx = _mm_mul_pd(fscal,dx00);
246 ty = _mm_mul_pd(fscal,dy00);
247 tz = _mm_mul_pd(fscal,dz00);
249 /* Update vectorial force */
250 fix0 = _mm_add_pd(fix0,tx);
251 fiy0 = _mm_add_pd(fiy0,ty);
252 fiz0 = _mm_add_pd(fiz0,tz);
254 gmx_mm_decrement_1rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,tx,ty,tz);
256 /* Inner loop uses 71 flops */
263 j_coord_offsetA = DIM*jnrA;
265 /* load j atom coordinates */
266 gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
269 /* Calculate displacement vector */
270 dx00 = _mm_sub_pd(ix0,jx0);
271 dy00 = _mm_sub_pd(iy0,jy0);
272 dz00 = _mm_sub_pd(iz0,jz0);
274 /* Calculate squared distance and things based on it */
275 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
277 rinv00 = gmx_mm_invsqrt_pd(rsq00);
279 rinvsq00 = _mm_mul_pd(rinv00,rinv00);
281 /* Load parameters for j particles */
282 jq0 = _mm_load_sd(charge+jnrA+0);
283 isaj0 = _mm_load_sd(invsqrta+jnrA+0);
284 vdwjidx0A = 2*vdwtype[jnrA+0];
286 /**************************
287 * CALCULATE INTERACTIONS *
288 **************************/
290 r00 = _mm_mul_pd(rsq00,rinv00);
292 /* Compute parameters for interactions between i and j atoms */
293 qq00 = _mm_mul_pd(iq0,jq0);
294 gmx_mm_load_1pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,&c6_00,&c12_00);
296 /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
297 isaprod = _mm_mul_pd(isai0,isaj0);
298 gbqqfactor = _mm_xor_pd(signbit,_mm_mul_pd(qq00,_mm_mul_pd(isaprod,gbinvepsdiff)));
299 gbscale = _mm_mul_pd(isaprod,gbtabscale);
301 /* Calculate generalized born table index - this is a separate table from the normal one,
302 * but we use the same procedure by multiplying r with scale and truncating to integer.
304 rt = _mm_mul_pd(r00,gbscale);
305 gbitab = _mm_cvttpd_epi32(rt);
306 gbeps = _mm_sub_pd(rt,_mm_cvtepi32_pd(gbitab));
307 gbitab = _mm_slli_epi32(gbitab,2);
309 Y = _mm_load_pd( gbtab + gmx_mm_extract_epi32(gbitab,0) );
310 F = _mm_setzero_pd();
311 GMX_MM_TRANSPOSE2_PD(Y,F);
312 G = _mm_load_pd( gbtab + gmx_mm_extract_epi32(gbitab,0) +2);
313 H = _mm_setzero_pd();
314 GMX_MM_TRANSPOSE2_PD(G,H);
315 Heps = _mm_mul_pd(gbeps,H);
316 Fp = _mm_add_pd(F,_mm_mul_pd(gbeps,_mm_add_pd(G,Heps)));
317 VV = _mm_add_pd(Y,_mm_mul_pd(gbeps,Fp));
318 vgb = _mm_mul_pd(gbqqfactor,VV);
320 FF = _mm_add_pd(Fp,_mm_mul_pd(gbeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
321 fgb = _mm_mul_pd(gbqqfactor,_mm_mul_pd(FF,gbscale));
322 dvdatmp = _mm_mul_pd(minushalf,_mm_add_pd(vgb,_mm_mul_pd(fgb,r00)));
323 dvdatmp = _mm_unpacklo_pd(dvdatmp,_mm_setzero_pd());
324 dvdasum = _mm_add_pd(dvdasum,dvdatmp);
325 gmx_mm_increment_1real_pd(dvda+jnrA,_mm_mul_pd(dvdatmp,_mm_mul_pd(isaj0,isaj0)));
326 velec = _mm_mul_pd(qq00,rinv00);
327 felec = _mm_mul_pd(_mm_sub_pd(_mm_mul_pd(velec,rinv00),fgb),rinv00);
329 /* LENNARD-JONES DISPERSION/REPULSION */
331 rinvsix = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
332 vvdw6 = _mm_mul_pd(c6_00,rinvsix);
333 vvdw12 = _mm_mul_pd(c12_00,_mm_mul_pd(rinvsix,rinvsix));
334 vvdw = _mm_sub_pd( _mm_mul_pd(vvdw12,one_twelfth) , _mm_mul_pd(vvdw6,one_sixth) );
335 fvdw = _mm_mul_pd(_mm_sub_pd(vvdw12,vvdw6),rinvsq00);
337 /* Update potential sum for this i atom from the interaction with this j atom. */
338 velec = _mm_unpacklo_pd(velec,_mm_setzero_pd());
339 velecsum = _mm_add_pd(velecsum,velec);
340 vgb = _mm_unpacklo_pd(vgb,_mm_setzero_pd());
341 vgbsum = _mm_add_pd(vgbsum,vgb);
342 vvdw = _mm_unpacklo_pd(vvdw,_mm_setzero_pd());
343 vvdwsum = _mm_add_pd(vvdwsum,vvdw);
345 fscal = _mm_add_pd(felec,fvdw);
347 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
349 /* Calculate temporary vectorial force */
350 tx = _mm_mul_pd(fscal,dx00);
351 ty = _mm_mul_pd(fscal,dy00);
352 tz = _mm_mul_pd(fscal,dz00);
354 /* Update vectorial force */
355 fix0 = _mm_add_pd(fix0,tx);
356 fiy0 = _mm_add_pd(fiy0,ty);
357 fiz0 = _mm_add_pd(fiz0,tz);
359 gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,tx,ty,tz);
361 /* Inner loop uses 71 flops */
364 /* End of innermost loop */
366 gmx_mm_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
367 f+i_coord_offset,fshift+i_shift_offset);
370 /* Update potential energies */
371 gmx_mm_update_1pot_pd(velecsum,kernel_data->energygrp_elec+ggid);
372 gmx_mm_update_1pot_pd(vgbsum,kernel_data->energygrp_polarization+ggid);
373 gmx_mm_update_1pot_pd(vvdwsum,kernel_data->energygrp_vdw+ggid);
374 dvdasum = _mm_mul_pd(dvdasum, _mm_mul_pd(isai0,isai0));
375 gmx_mm_update_1pot_pd(dvdasum,dvda+inr);
377 /* Increment number of inner iterations */
378 inneriter += j_index_end - j_index_start;
380 /* Outer loop uses 10 flops */
383 /* Increment number of outer iterations */
386 /* Update outer/inner flops */
388 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_VF,outeriter*10 + inneriter*71);
391 * Gromacs nonbonded kernel: nb_kernel_ElecGB_VdwLJ_GeomP1P1_F_sse2_double
392 * Electrostatics interaction: GeneralizedBorn
393 * VdW interaction: LennardJones
394 * Geometry: Particle-Particle
395 * Calculate force/pot: Force
398 nb_kernel_ElecGB_VdwLJ_GeomP1P1_F_sse2_double
399 (t_nblist * gmx_restrict nlist,
400 rvec * gmx_restrict xx,
401 rvec * gmx_restrict ff,
402 t_forcerec * gmx_restrict fr,
403 t_mdatoms * gmx_restrict mdatoms,
404 nb_kernel_data_t * gmx_restrict kernel_data,
405 t_nrnb * gmx_restrict nrnb)
407 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
408 * just 0 for non-waters.
409 * Suffixes A,B refer to j loop unrolling done with SSE double precision, e.g. for the two different
410 * jnr indices corresponding to data put in the four positions in the SIMD register.
412 int i_shift_offset,i_coord_offset,outeriter,inneriter;
413 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
415 int j_coord_offsetA,j_coord_offsetB;
416 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
418 real *shiftvec,*fshift,*x,*f;
419 __m128d tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
421 __m128d ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
422 int vdwjidx0A,vdwjidx0B;
423 __m128d jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
424 __m128d dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
425 __m128d velec,felec,velecsum,facel,crf,krf,krf2;
428 __m128d vgb,fgb,vgbsum,dvdasum,gbscale,gbtabscale,isaprod,gbqqfactor,gbinvepsdiff,dvdaj,gbeps,dvdatmp;
429 __m128d minushalf = _mm_set1_pd(-0.5);
430 real *invsqrta,*dvda,*gbtab;
432 __m128d rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
435 __m128d one_sixth = _mm_set1_pd(1.0/6.0);
436 __m128d one_twelfth = _mm_set1_pd(1.0/12.0);
438 __m128i ifour = _mm_set1_epi32(4);
439 __m128d rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF;
441 __m128d dummy_mask,cutoff_mask;
442 __m128d signbit = gmx_mm_castsi128_pd( _mm_set_epi32(0x80000000,0x00000000,0x80000000,0x00000000) );
443 __m128d one = _mm_set1_pd(1.0);
444 __m128d two = _mm_set1_pd(2.0);
450 jindex = nlist->jindex;
452 shiftidx = nlist->shift;
454 shiftvec = fr->shift_vec[0];
455 fshift = fr->fshift[0];
456 facel = _mm_set1_pd(fr->epsfac);
457 charge = mdatoms->chargeA;
458 nvdwtype = fr->ntype;
460 vdwtype = mdatoms->typeA;
462 invsqrta = fr->invsqrta;
464 gbtabscale = _mm_set1_pd(fr->gbtab.scale);
465 gbtab = fr->gbtab.data;
466 gbinvepsdiff = _mm_set1_pd((1.0/fr->epsilon_r) - (1.0/fr->gb_epsilon_solvent));
468 /* Avoid stupid compiler warnings */
476 /* Start outer loop over neighborlists */
477 for(iidx=0; iidx<nri; iidx++)
479 /* Load shift vector for this list */
480 i_shift_offset = DIM*shiftidx[iidx];
482 /* Load limits for loop over neighbors */
483 j_index_start = jindex[iidx];
484 j_index_end = jindex[iidx+1];
486 /* Get outer coordinate index */
488 i_coord_offset = DIM*inr;
490 /* Load i particle coords and add shift vector */
491 gmx_mm_load_shift_and_1rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
493 fix0 = _mm_setzero_pd();
494 fiy0 = _mm_setzero_pd();
495 fiz0 = _mm_setzero_pd();
497 /* Load parameters for i particles */
498 iq0 = _mm_mul_pd(facel,_mm_load1_pd(charge+inr+0));
499 isai0 = _mm_load1_pd(invsqrta+inr+0);
500 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
502 dvdasum = _mm_setzero_pd();
504 /* Start inner kernel loop */
505 for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
508 /* Get j neighbor index, and coordinate index */
511 j_coord_offsetA = DIM*jnrA;
512 j_coord_offsetB = DIM*jnrB;
514 /* load j atom coordinates */
515 gmx_mm_load_1rvec_2ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
518 /* Calculate displacement vector */
519 dx00 = _mm_sub_pd(ix0,jx0);
520 dy00 = _mm_sub_pd(iy0,jy0);
521 dz00 = _mm_sub_pd(iz0,jz0);
523 /* Calculate squared distance and things based on it */
524 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
526 rinv00 = gmx_mm_invsqrt_pd(rsq00);
528 rinvsq00 = _mm_mul_pd(rinv00,rinv00);
530 /* Load parameters for j particles */
531 jq0 = gmx_mm_load_2real_swizzle_pd(charge+jnrA+0,charge+jnrB+0);
532 isaj0 = gmx_mm_load_2real_swizzle_pd(invsqrta+jnrA+0,invsqrta+jnrB+0);
533 vdwjidx0A = 2*vdwtype[jnrA+0];
534 vdwjidx0B = 2*vdwtype[jnrB+0];
536 /**************************
537 * CALCULATE INTERACTIONS *
538 **************************/
540 r00 = _mm_mul_pd(rsq00,rinv00);
542 /* Compute parameters for interactions between i and j atoms */
543 qq00 = _mm_mul_pd(iq0,jq0);
544 gmx_mm_load_2pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,
545 vdwparam+vdwioffset0+vdwjidx0B,&c6_00,&c12_00);
547 /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
548 isaprod = _mm_mul_pd(isai0,isaj0);
549 gbqqfactor = _mm_xor_pd(signbit,_mm_mul_pd(qq00,_mm_mul_pd(isaprod,gbinvepsdiff)));
550 gbscale = _mm_mul_pd(isaprod,gbtabscale);
552 /* Calculate generalized born table index - this is a separate table from the normal one,
553 * but we use the same procedure by multiplying r with scale and truncating to integer.
555 rt = _mm_mul_pd(r00,gbscale);
556 gbitab = _mm_cvttpd_epi32(rt);
557 gbeps = _mm_sub_pd(rt,_mm_cvtepi32_pd(gbitab));
558 gbitab = _mm_slli_epi32(gbitab,2);
560 Y = _mm_load_pd( gbtab + gmx_mm_extract_epi32(gbitab,0) );
561 F = _mm_load_pd( gbtab + gmx_mm_extract_epi32(gbitab,1) );
562 GMX_MM_TRANSPOSE2_PD(Y,F);
563 G = _mm_load_pd( gbtab + gmx_mm_extract_epi32(gbitab,0) +2);
564 H = _mm_load_pd( gbtab + gmx_mm_extract_epi32(gbitab,1) +2);
565 GMX_MM_TRANSPOSE2_PD(G,H);
566 Heps = _mm_mul_pd(gbeps,H);
567 Fp = _mm_add_pd(F,_mm_mul_pd(gbeps,_mm_add_pd(G,Heps)));
568 VV = _mm_add_pd(Y,_mm_mul_pd(gbeps,Fp));
569 vgb = _mm_mul_pd(gbqqfactor,VV);
571 FF = _mm_add_pd(Fp,_mm_mul_pd(gbeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
572 fgb = _mm_mul_pd(gbqqfactor,_mm_mul_pd(FF,gbscale));
573 dvdatmp = _mm_mul_pd(minushalf,_mm_add_pd(vgb,_mm_mul_pd(fgb,r00)));
574 dvdasum = _mm_add_pd(dvdasum,dvdatmp);
575 gmx_mm_increment_2real_swizzle_pd(dvda+jnrA,dvda+jnrB,_mm_mul_pd(dvdatmp,_mm_mul_pd(isaj0,isaj0)));
576 velec = _mm_mul_pd(qq00,rinv00);
577 felec = _mm_mul_pd(_mm_sub_pd(_mm_mul_pd(velec,rinv00),fgb),rinv00);
579 /* LENNARD-JONES DISPERSION/REPULSION */
581 rinvsix = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
582 fvdw = _mm_mul_pd(_mm_sub_pd(_mm_mul_pd(c12_00,rinvsix),c6_00),_mm_mul_pd(rinvsix,rinvsq00));
584 fscal = _mm_add_pd(felec,fvdw);
586 /* Calculate temporary vectorial force */
587 tx = _mm_mul_pd(fscal,dx00);
588 ty = _mm_mul_pd(fscal,dy00);
589 tz = _mm_mul_pd(fscal,dz00);
591 /* Update vectorial force */
592 fix0 = _mm_add_pd(fix0,tx);
593 fiy0 = _mm_add_pd(fiy0,ty);
594 fiz0 = _mm_add_pd(fiz0,tz);
596 gmx_mm_decrement_1rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,tx,ty,tz);
598 /* Inner loop uses 64 flops */
605 j_coord_offsetA = DIM*jnrA;
607 /* load j atom coordinates */
608 gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
611 /* Calculate displacement vector */
612 dx00 = _mm_sub_pd(ix0,jx0);
613 dy00 = _mm_sub_pd(iy0,jy0);
614 dz00 = _mm_sub_pd(iz0,jz0);
616 /* Calculate squared distance and things based on it */
617 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
619 rinv00 = gmx_mm_invsqrt_pd(rsq00);
621 rinvsq00 = _mm_mul_pd(rinv00,rinv00);
623 /* Load parameters for j particles */
624 jq0 = _mm_load_sd(charge+jnrA+0);
625 isaj0 = _mm_load_sd(invsqrta+jnrA+0);
626 vdwjidx0A = 2*vdwtype[jnrA+0];
628 /**************************
629 * CALCULATE INTERACTIONS *
630 **************************/
632 r00 = _mm_mul_pd(rsq00,rinv00);
634 /* Compute parameters for interactions between i and j atoms */
635 qq00 = _mm_mul_pd(iq0,jq0);
636 gmx_mm_load_1pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,&c6_00,&c12_00);
638 /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
639 isaprod = _mm_mul_pd(isai0,isaj0);
640 gbqqfactor = _mm_xor_pd(signbit,_mm_mul_pd(qq00,_mm_mul_pd(isaprod,gbinvepsdiff)));
641 gbscale = _mm_mul_pd(isaprod,gbtabscale);
643 /* Calculate generalized born table index - this is a separate table from the normal one,
644 * but we use the same procedure by multiplying r with scale and truncating to integer.
646 rt = _mm_mul_pd(r00,gbscale);
647 gbitab = _mm_cvttpd_epi32(rt);
648 gbeps = _mm_sub_pd(rt,_mm_cvtepi32_pd(gbitab));
649 gbitab = _mm_slli_epi32(gbitab,2);
651 Y = _mm_load_pd( gbtab + gmx_mm_extract_epi32(gbitab,0) );
652 F = _mm_setzero_pd();
653 GMX_MM_TRANSPOSE2_PD(Y,F);
654 G = _mm_load_pd( gbtab + gmx_mm_extract_epi32(gbitab,0) +2);
655 H = _mm_setzero_pd();
656 GMX_MM_TRANSPOSE2_PD(G,H);
657 Heps = _mm_mul_pd(gbeps,H);
658 Fp = _mm_add_pd(F,_mm_mul_pd(gbeps,_mm_add_pd(G,Heps)));
659 VV = _mm_add_pd(Y,_mm_mul_pd(gbeps,Fp));
660 vgb = _mm_mul_pd(gbqqfactor,VV);
662 FF = _mm_add_pd(Fp,_mm_mul_pd(gbeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
663 fgb = _mm_mul_pd(gbqqfactor,_mm_mul_pd(FF,gbscale));
664 dvdatmp = _mm_mul_pd(minushalf,_mm_add_pd(vgb,_mm_mul_pd(fgb,r00)));
665 dvdatmp = _mm_unpacklo_pd(dvdatmp,_mm_setzero_pd());
666 dvdasum = _mm_add_pd(dvdasum,dvdatmp);
667 gmx_mm_increment_1real_pd(dvda+jnrA,_mm_mul_pd(dvdatmp,_mm_mul_pd(isaj0,isaj0)));
668 velec = _mm_mul_pd(qq00,rinv00);
669 felec = _mm_mul_pd(_mm_sub_pd(_mm_mul_pd(velec,rinv00),fgb),rinv00);
671 /* LENNARD-JONES DISPERSION/REPULSION */
673 rinvsix = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
674 fvdw = _mm_mul_pd(_mm_sub_pd(_mm_mul_pd(c12_00,rinvsix),c6_00),_mm_mul_pd(rinvsix,rinvsq00));
676 fscal = _mm_add_pd(felec,fvdw);
678 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
680 /* Calculate temporary vectorial force */
681 tx = _mm_mul_pd(fscal,dx00);
682 ty = _mm_mul_pd(fscal,dy00);
683 tz = _mm_mul_pd(fscal,dz00);
685 /* Update vectorial force */
686 fix0 = _mm_add_pd(fix0,tx);
687 fiy0 = _mm_add_pd(fiy0,ty);
688 fiz0 = _mm_add_pd(fiz0,tz);
690 gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,tx,ty,tz);
692 /* Inner loop uses 64 flops */
695 /* End of innermost loop */
697 gmx_mm_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
698 f+i_coord_offset,fshift+i_shift_offset);
700 dvdasum = _mm_mul_pd(dvdasum, _mm_mul_pd(isai0,isai0));
701 gmx_mm_update_1pot_pd(dvdasum,dvda+inr);
703 /* Increment number of inner iterations */
704 inneriter += j_index_end - j_index_start;
706 /* Outer loop uses 7 flops */
709 /* Increment number of outer iterations */
712 /* Update outer/inner flops */
714 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_F,outeriter*7 + inneriter*64);