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_VdwCSTab_GeomP1P1_VF_sse2_double
38 * Electrostatics interaction: GeneralizedBorn
39 * VdW interaction: CubicSplineTable
40 * Geometry: Particle-Particle
41 * Calculate force/pot: PotentialAndForce
44 nb_kernel_ElecGB_VdwCSTab_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 vftab = kernel_data->table_vdw->data;
109 vftabscale = _mm_set1_pd(kernel_data->table_vdw->scale);
111 invsqrta = fr->invsqrta;
113 gbtabscale = _mm_set1_pd(fr->gbtab.scale);
114 gbtab = fr->gbtab.data;
115 gbinvepsdiff = _mm_set1_pd((1.0/fr->epsilon_r) - (1.0/fr->gb_epsilon_solvent));
117 /* Avoid stupid compiler warnings */
125 /* Start outer loop over neighborlists */
126 for(iidx=0; iidx<nri; iidx++)
128 /* Load shift vector for this list */
129 i_shift_offset = DIM*shiftidx[iidx];
131 /* Load limits for loop over neighbors */
132 j_index_start = jindex[iidx];
133 j_index_end = jindex[iidx+1];
135 /* Get outer coordinate index */
137 i_coord_offset = DIM*inr;
139 /* Load i particle coords and add shift vector */
140 gmx_mm_load_shift_and_1rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
142 fix0 = _mm_setzero_pd();
143 fiy0 = _mm_setzero_pd();
144 fiz0 = _mm_setzero_pd();
146 /* Load parameters for i particles */
147 iq0 = _mm_mul_pd(facel,_mm_load1_pd(charge+inr+0));
148 isai0 = _mm_load1_pd(invsqrta+inr+0);
149 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
151 /* Reset potential sums */
152 velecsum = _mm_setzero_pd();
153 vgbsum = _mm_setzero_pd();
154 vvdwsum = _mm_setzero_pd();
155 dvdasum = _mm_setzero_pd();
157 /* Start inner kernel loop */
158 for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
161 /* Get j neighbor index, and coordinate index */
164 j_coord_offsetA = DIM*jnrA;
165 j_coord_offsetB = DIM*jnrB;
167 /* load j atom coordinates */
168 gmx_mm_load_1rvec_2ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
171 /* Calculate displacement vector */
172 dx00 = _mm_sub_pd(ix0,jx0);
173 dy00 = _mm_sub_pd(iy0,jy0);
174 dz00 = _mm_sub_pd(iz0,jz0);
176 /* Calculate squared distance and things based on it */
177 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
179 rinv00 = gmx_mm_invsqrt_pd(rsq00);
181 /* Load parameters for j particles */
182 jq0 = gmx_mm_load_2real_swizzle_pd(charge+jnrA+0,charge+jnrB+0);
183 isaj0 = gmx_mm_load_2real_swizzle_pd(invsqrta+jnrA+0,invsqrta+jnrB+0);
184 vdwjidx0A = 2*vdwtype[jnrA+0];
185 vdwjidx0B = 2*vdwtype[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);
195 gmx_mm_load_2pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,
196 vdwparam+vdwioffset0+vdwjidx0B,&c6_00,&c12_00);
198 /* Calculate table index by multiplying r with table scale and truncate to integer */
199 rt = _mm_mul_pd(r00,vftabscale);
200 vfitab = _mm_cvttpd_epi32(rt);
201 vfeps = _mm_sub_pd(rt,_mm_cvtepi32_pd(vfitab));
202 vfitab = _mm_slli_epi32(vfitab,3);
204 /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
205 isaprod = _mm_mul_pd(isai0,isaj0);
206 gbqqfactor = _mm_xor_pd(signbit,_mm_mul_pd(qq00,_mm_mul_pd(isaprod,gbinvepsdiff)));
207 gbscale = _mm_mul_pd(isaprod,gbtabscale);
209 /* Calculate generalized born table index - this is a separate table from the normal one,
210 * but we use the same procedure by multiplying r with scale and truncating to integer.
212 rt = _mm_mul_pd(r00,gbscale);
213 gbitab = _mm_cvttpd_epi32(rt);
214 gbeps = _mm_sub_pd(rt,_mm_cvtepi32_pd(gbitab));
215 gbitab = _mm_slli_epi32(gbitab,2);
217 Y = _mm_load_pd( gbtab + gmx_mm_extract_epi32(gbitab,0) );
218 F = _mm_load_pd( gbtab + gmx_mm_extract_epi32(gbitab,1) );
219 GMX_MM_TRANSPOSE2_PD(Y,F);
220 G = _mm_load_pd( gbtab + gmx_mm_extract_epi32(gbitab,0) +2);
221 H = _mm_load_pd( gbtab + gmx_mm_extract_epi32(gbitab,1) +2);
222 GMX_MM_TRANSPOSE2_PD(G,H);
223 Heps = _mm_mul_pd(gbeps,H);
224 Fp = _mm_add_pd(F,_mm_mul_pd(gbeps,_mm_add_pd(G,Heps)));
225 VV = _mm_add_pd(Y,_mm_mul_pd(gbeps,Fp));
226 vgb = _mm_mul_pd(gbqqfactor,VV);
228 FF = _mm_add_pd(Fp,_mm_mul_pd(gbeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
229 fgb = _mm_mul_pd(gbqqfactor,_mm_mul_pd(FF,gbscale));
230 dvdatmp = _mm_mul_pd(minushalf,_mm_add_pd(vgb,_mm_mul_pd(fgb,r00)));
231 dvdasum = _mm_add_pd(dvdasum,dvdatmp);
232 gmx_mm_increment_2real_swizzle_pd(dvda+jnrA,dvda+jnrB,_mm_mul_pd(dvdatmp,_mm_mul_pd(isaj0,isaj0)));
233 velec = _mm_mul_pd(qq00,rinv00);
234 felec = _mm_mul_pd(_mm_sub_pd(_mm_mul_pd(velec,rinv00),fgb),rinv00);
236 /* CUBIC SPLINE TABLE DISPERSION */
237 Y = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) );
238 F = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,1) );
239 GMX_MM_TRANSPOSE2_PD(Y,F);
240 G = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) +2);
241 H = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,1) +2);
242 GMX_MM_TRANSPOSE2_PD(G,H);
243 Heps = _mm_mul_pd(vfeps,H);
244 Fp = _mm_add_pd(F,_mm_mul_pd(vfeps,_mm_add_pd(G,Heps)));
245 VV = _mm_add_pd(Y,_mm_mul_pd(vfeps,Fp));
246 vvdw6 = _mm_mul_pd(c6_00,VV);
247 FF = _mm_add_pd(Fp,_mm_mul_pd(vfeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
248 fvdw6 = _mm_mul_pd(c6_00,FF);
250 /* CUBIC SPLINE TABLE REPULSION */
251 vfitab = _mm_add_epi32(vfitab,ifour);
252 Y = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) );
253 F = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,1) );
254 GMX_MM_TRANSPOSE2_PD(Y,F);
255 G = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) +2);
256 H = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,1) +2);
257 GMX_MM_TRANSPOSE2_PD(G,H);
258 Heps = _mm_mul_pd(vfeps,H);
259 Fp = _mm_add_pd(F,_mm_mul_pd(vfeps,_mm_add_pd(G,Heps)));
260 VV = _mm_add_pd(Y,_mm_mul_pd(vfeps,Fp));
261 vvdw12 = _mm_mul_pd(c12_00,VV);
262 FF = _mm_add_pd(Fp,_mm_mul_pd(vfeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
263 fvdw12 = _mm_mul_pd(c12_00,FF);
264 vvdw = _mm_add_pd(vvdw12,vvdw6);
265 fvdw = _mm_xor_pd(signbit,_mm_mul_pd(_mm_add_pd(fvdw6,fvdw12),_mm_mul_pd(vftabscale,rinv00)));
267 /* Update potential sum for this i atom from the interaction with this j atom. */
268 velecsum = _mm_add_pd(velecsum,velec);
269 vgbsum = _mm_add_pd(vgbsum,vgb);
270 vvdwsum = _mm_add_pd(vvdwsum,vvdw);
272 fscal = _mm_add_pd(felec,fvdw);
274 /* Calculate temporary vectorial force */
275 tx = _mm_mul_pd(fscal,dx00);
276 ty = _mm_mul_pd(fscal,dy00);
277 tz = _mm_mul_pd(fscal,dz00);
279 /* Update vectorial force */
280 fix0 = _mm_add_pd(fix0,tx);
281 fiy0 = _mm_add_pd(fiy0,ty);
282 fiz0 = _mm_add_pd(fiz0,tz);
284 gmx_mm_decrement_1rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,tx,ty,tz);
286 /* Inner loop uses 92 flops */
293 j_coord_offsetA = DIM*jnrA;
295 /* load j atom coordinates */
296 gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
299 /* Calculate displacement vector */
300 dx00 = _mm_sub_pd(ix0,jx0);
301 dy00 = _mm_sub_pd(iy0,jy0);
302 dz00 = _mm_sub_pd(iz0,jz0);
304 /* Calculate squared distance and things based on it */
305 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
307 rinv00 = gmx_mm_invsqrt_pd(rsq00);
309 /* Load parameters for j particles */
310 jq0 = _mm_load_sd(charge+jnrA+0);
311 isaj0 = _mm_load_sd(invsqrta+jnrA+0);
312 vdwjidx0A = 2*vdwtype[jnrA+0];
314 /**************************
315 * CALCULATE INTERACTIONS *
316 **************************/
318 r00 = _mm_mul_pd(rsq00,rinv00);
320 /* Compute parameters for interactions between i and j atoms */
321 qq00 = _mm_mul_pd(iq0,jq0);
322 gmx_mm_load_1pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,&c6_00,&c12_00);
324 /* Calculate table index by multiplying r with table scale and truncate to integer */
325 rt = _mm_mul_pd(r00,vftabscale);
326 vfitab = _mm_cvttpd_epi32(rt);
327 vfeps = _mm_sub_pd(rt,_mm_cvtepi32_pd(vfitab));
328 vfitab = _mm_slli_epi32(vfitab,3);
330 /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
331 isaprod = _mm_mul_pd(isai0,isaj0);
332 gbqqfactor = _mm_xor_pd(signbit,_mm_mul_pd(qq00,_mm_mul_pd(isaprod,gbinvepsdiff)));
333 gbscale = _mm_mul_pd(isaprod,gbtabscale);
335 /* Calculate generalized born table index - this is a separate table from the normal one,
336 * but we use the same procedure by multiplying r with scale and truncating to integer.
338 rt = _mm_mul_pd(r00,gbscale);
339 gbitab = _mm_cvttpd_epi32(rt);
340 gbeps = _mm_sub_pd(rt,_mm_cvtepi32_pd(gbitab));
341 gbitab = _mm_slli_epi32(gbitab,2);
343 Y = _mm_load_pd( gbtab + gmx_mm_extract_epi32(gbitab,0) );
344 F = _mm_setzero_pd();
345 GMX_MM_TRANSPOSE2_PD(Y,F);
346 G = _mm_load_pd( gbtab + gmx_mm_extract_epi32(gbitab,0) +2);
347 H = _mm_setzero_pd();
348 GMX_MM_TRANSPOSE2_PD(G,H);
349 Heps = _mm_mul_pd(gbeps,H);
350 Fp = _mm_add_pd(F,_mm_mul_pd(gbeps,_mm_add_pd(G,Heps)));
351 VV = _mm_add_pd(Y,_mm_mul_pd(gbeps,Fp));
352 vgb = _mm_mul_pd(gbqqfactor,VV);
354 FF = _mm_add_pd(Fp,_mm_mul_pd(gbeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
355 fgb = _mm_mul_pd(gbqqfactor,_mm_mul_pd(FF,gbscale));
356 dvdatmp = _mm_mul_pd(minushalf,_mm_add_pd(vgb,_mm_mul_pd(fgb,r00)));
357 dvdasum = _mm_add_pd(dvdasum,dvdatmp);
358 gmx_mm_increment_1real_pd(dvda+jnrA,_mm_mul_pd(dvdatmp,_mm_mul_pd(isaj0,isaj0)));
359 velec = _mm_mul_pd(qq00,rinv00);
360 felec = _mm_mul_pd(_mm_sub_pd(_mm_mul_pd(velec,rinv00),fgb),rinv00);
362 /* CUBIC SPLINE TABLE DISPERSION */
363 Y = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) );
364 F = _mm_setzero_pd();
365 GMX_MM_TRANSPOSE2_PD(Y,F);
366 G = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) +2);
367 H = _mm_setzero_pd();
368 GMX_MM_TRANSPOSE2_PD(G,H);
369 Heps = _mm_mul_pd(vfeps,H);
370 Fp = _mm_add_pd(F,_mm_mul_pd(vfeps,_mm_add_pd(G,Heps)));
371 VV = _mm_add_pd(Y,_mm_mul_pd(vfeps,Fp));
372 vvdw6 = _mm_mul_pd(c6_00,VV);
373 FF = _mm_add_pd(Fp,_mm_mul_pd(vfeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
374 fvdw6 = _mm_mul_pd(c6_00,FF);
376 /* CUBIC SPLINE TABLE REPULSION */
377 vfitab = _mm_add_epi32(vfitab,ifour);
378 Y = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) );
379 F = _mm_setzero_pd();
380 GMX_MM_TRANSPOSE2_PD(Y,F);
381 G = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) +2);
382 H = _mm_setzero_pd();
383 GMX_MM_TRANSPOSE2_PD(G,H);
384 Heps = _mm_mul_pd(vfeps,H);
385 Fp = _mm_add_pd(F,_mm_mul_pd(vfeps,_mm_add_pd(G,Heps)));
386 VV = _mm_add_pd(Y,_mm_mul_pd(vfeps,Fp));
387 vvdw12 = _mm_mul_pd(c12_00,VV);
388 FF = _mm_add_pd(Fp,_mm_mul_pd(vfeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
389 fvdw12 = _mm_mul_pd(c12_00,FF);
390 vvdw = _mm_add_pd(vvdw12,vvdw6);
391 fvdw = _mm_xor_pd(signbit,_mm_mul_pd(_mm_add_pd(fvdw6,fvdw12),_mm_mul_pd(vftabscale,rinv00)));
393 /* Update potential sum for this i atom from the interaction with this j atom. */
394 velec = _mm_unpacklo_pd(velec,_mm_setzero_pd());
395 velecsum = _mm_add_pd(velecsum,velec);
396 vgb = _mm_unpacklo_pd(vgb,_mm_setzero_pd());
397 vgbsum = _mm_add_pd(vgbsum,vgb);
398 vvdw = _mm_unpacklo_pd(vvdw,_mm_setzero_pd());
399 vvdwsum = _mm_add_pd(vvdwsum,vvdw);
401 fscal = _mm_add_pd(felec,fvdw);
403 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
405 /* Calculate temporary vectorial force */
406 tx = _mm_mul_pd(fscal,dx00);
407 ty = _mm_mul_pd(fscal,dy00);
408 tz = _mm_mul_pd(fscal,dz00);
410 /* Update vectorial force */
411 fix0 = _mm_add_pd(fix0,tx);
412 fiy0 = _mm_add_pd(fiy0,ty);
413 fiz0 = _mm_add_pd(fiz0,tz);
415 gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,tx,ty,tz);
417 /* Inner loop uses 92 flops */
420 /* End of innermost loop */
422 gmx_mm_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
423 f+i_coord_offset,fshift+i_shift_offset);
426 /* Update potential energies */
427 gmx_mm_update_1pot_pd(velecsum,kernel_data->energygrp_elec+ggid);
428 gmx_mm_update_1pot_pd(vgbsum,kernel_data->energygrp_polarization+ggid);
429 gmx_mm_update_1pot_pd(vvdwsum,kernel_data->energygrp_vdw+ggid);
430 dvdasum = _mm_mul_pd(dvdasum, _mm_mul_pd(isai0,isai0));
431 gmx_mm_update_1pot_pd(dvdasum,dvda+inr);
433 /* Increment number of inner iterations */
434 inneriter += j_index_end - j_index_start;
436 /* Outer loop uses 10 flops */
439 /* Increment number of outer iterations */
442 /* Update outer/inner flops */
444 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_VF,outeriter*10 + inneriter*92);
447 * Gromacs nonbonded kernel: nb_kernel_ElecGB_VdwCSTab_GeomP1P1_F_sse2_double
448 * Electrostatics interaction: GeneralizedBorn
449 * VdW interaction: CubicSplineTable
450 * Geometry: Particle-Particle
451 * Calculate force/pot: Force
454 nb_kernel_ElecGB_VdwCSTab_GeomP1P1_F_sse2_double
455 (t_nblist * gmx_restrict nlist,
456 rvec * gmx_restrict xx,
457 rvec * gmx_restrict ff,
458 t_forcerec * gmx_restrict fr,
459 t_mdatoms * gmx_restrict mdatoms,
460 nb_kernel_data_t * gmx_restrict kernel_data,
461 t_nrnb * gmx_restrict nrnb)
463 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
464 * just 0 for non-waters.
465 * Suffixes A,B refer to j loop unrolling done with SSE double precision, e.g. for the two different
466 * jnr indices corresponding to data put in the four positions in the SIMD register.
468 int i_shift_offset,i_coord_offset,outeriter,inneriter;
469 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
471 int j_coord_offsetA,j_coord_offsetB;
472 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
474 real *shiftvec,*fshift,*x,*f;
475 __m128d tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
477 __m128d ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
478 int vdwjidx0A,vdwjidx0B;
479 __m128d jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
480 __m128d dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
481 __m128d velec,felec,velecsum,facel,crf,krf,krf2;
484 __m128d vgb,fgb,vgbsum,dvdasum,gbscale,gbtabscale,isaprod,gbqqfactor,gbinvepsdiff,dvdaj,gbeps,dvdatmp;
485 __m128d minushalf = _mm_set1_pd(-0.5);
486 real *invsqrta,*dvda,*gbtab;
488 __m128d rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
491 __m128d one_sixth = _mm_set1_pd(1.0/6.0);
492 __m128d one_twelfth = _mm_set1_pd(1.0/12.0);
494 __m128i ifour = _mm_set1_epi32(4);
495 __m128d rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF;
497 __m128d dummy_mask,cutoff_mask;
498 __m128d signbit = gmx_mm_castsi128_pd( _mm_set_epi32(0x80000000,0x00000000,0x80000000,0x00000000) );
499 __m128d one = _mm_set1_pd(1.0);
500 __m128d two = _mm_set1_pd(2.0);
506 jindex = nlist->jindex;
508 shiftidx = nlist->shift;
510 shiftvec = fr->shift_vec[0];
511 fshift = fr->fshift[0];
512 facel = _mm_set1_pd(fr->epsfac);
513 charge = mdatoms->chargeA;
514 nvdwtype = fr->ntype;
516 vdwtype = mdatoms->typeA;
518 vftab = kernel_data->table_vdw->data;
519 vftabscale = _mm_set1_pd(kernel_data->table_vdw->scale);
521 invsqrta = fr->invsqrta;
523 gbtabscale = _mm_set1_pd(fr->gbtab.scale);
524 gbtab = fr->gbtab.data;
525 gbinvepsdiff = _mm_set1_pd((1.0/fr->epsilon_r) - (1.0/fr->gb_epsilon_solvent));
527 /* Avoid stupid compiler warnings */
535 /* Start outer loop over neighborlists */
536 for(iidx=0; iidx<nri; iidx++)
538 /* Load shift vector for this list */
539 i_shift_offset = DIM*shiftidx[iidx];
541 /* Load limits for loop over neighbors */
542 j_index_start = jindex[iidx];
543 j_index_end = jindex[iidx+1];
545 /* Get outer coordinate index */
547 i_coord_offset = DIM*inr;
549 /* Load i particle coords and add shift vector */
550 gmx_mm_load_shift_and_1rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
552 fix0 = _mm_setzero_pd();
553 fiy0 = _mm_setzero_pd();
554 fiz0 = _mm_setzero_pd();
556 /* Load parameters for i particles */
557 iq0 = _mm_mul_pd(facel,_mm_load1_pd(charge+inr+0));
558 isai0 = _mm_load1_pd(invsqrta+inr+0);
559 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
561 dvdasum = _mm_setzero_pd();
563 /* Start inner kernel loop */
564 for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
567 /* Get j neighbor index, and coordinate index */
570 j_coord_offsetA = DIM*jnrA;
571 j_coord_offsetB = DIM*jnrB;
573 /* load j atom coordinates */
574 gmx_mm_load_1rvec_2ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
577 /* Calculate displacement vector */
578 dx00 = _mm_sub_pd(ix0,jx0);
579 dy00 = _mm_sub_pd(iy0,jy0);
580 dz00 = _mm_sub_pd(iz0,jz0);
582 /* Calculate squared distance and things based on it */
583 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
585 rinv00 = gmx_mm_invsqrt_pd(rsq00);
587 /* Load parameters for j particles */
588 jq0 = gmx_mm_load_2real_swizzle_pd(charge+jnrA+0,charge+jnrB+0);
589 isaj0 = gmx_mm_load_2real_swizzle_pd(invsqrta+jnrA+0,invsqrta+jnrB+0);
590 vdwjidx0A = 2*vdwtype[jnrA+0];
591 vdwjidx0B = 2*vdwtype[jnrB+0];
593 /**************************
594 * CALCULATE INTERACTIONS *
595 **************************/
597 r00 = _mm_mul_pd(rsq00,rinv00);
599 /* Compute parameters for interactions between i and j atoms */
600 qq00 = _mm_mul_pd(iq0,jq0);
601 gmx_mm_load_2pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,
602 vdwparam+vdwioffset0+vdwjidx0B,&c6_00,&c12_00);
604 /* Calculate table index by multiplying r with table scale and truncate to integer */
605 rt = _mm_mul_pd(r00,vftabscale);
606 vfitab = _mm_cvttpd_epi32(rt);
607 vfeps = _mm_sub_pd(rt,_mm_cvtepi32_pd(vfitab));
608 vfitab = _mm_slli_epi32(vfitab,3);
610 /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
611 isaprod = _mm_mul_pd(isai0,isaj0);
612 gbqqfactor = _mm_xor_pd(signbit,_mm_mul_pd(qq00,_mm_mul_pd(isaprod,gbinvepsdiff)));
613 gbscale = _mm_mul_pd(isaprod,gbtabscale);
615 /* Calculate generalized born table index - this is a separate table from the normal one,
616 * but we use the same procedure by multiplying r with scale and truncating to integer.
618 rt = _mm_mul_pd(r00,gbscale);
619 gbitab = _mm_cvttpd_epi32(rt);
620 gbeps = _mm_sub_pd(rt,_mm_cvtepi32_pd(gbitab));
621 gbitab = _mm_slli_epi32(gbitab,2);
623 Y = _mm_load_pd( gbtab + gmx_mm_extract_epi32(gbitab,0) );
624 F = _mm_load_pd( gbtab + gmx_mm_extract_epi32(gbitab,1) );
625 GMX_MM_TRANSPOSE2_PD(Y,F);
626 G = _mm_load_pd( gbtab + gmx_mm_extract_epi32(gbitab,0) +2);
627 H = _mm_load_pd( gbtab + gmx_mm_extract_epi32(gbitab,1) +2);
628 GMX_MM_TRANSPOSE2_PD(G,H);
629 Heps = _mm_mul_pd(gbeps,H);
630 Fp = _mm_add_pd(F,_mm_mul_pd(gbeps,_mm_add_pd(G,Heps)));
631 VV = _mm_add_pd(Y,_mm_mul_pd(gbeps,Fp));
632 vgb = _mm_mul_pd(gbqqfactor,VV);
634 FF = _mm_add_pd(Fp,_mm_mul_pd(gbeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
635 fgb = _mm_mul_pd(gbqqfactor,_mm_mul_pd(FF,gbscale));
636 dvdatmp = _mm_mul_pd(minushalf,_mm_add_pd(vgb,_mm_mul_pd(fgb,r00)));
637 dvdasum = _mm_add_pd(dvdasum,dvdatmp);
638 gmx_mm_increment_2real_swizzle_pd(dvda+jnrA,dvda+jnrB,_mm_mul_pd(dvdatmp,_mm_mul_pd(isaj0,isaj0)));
639 velec = _mm_mul_pd(qq00,rinv00);
640 felec = _mm_mul_pd(_mm_sub_pd(_mm_mul_pd(velec,rinv00),fgb),rinv00);
642 /* CUBIC SPLINE TABLE DISPERSION */
643 Y = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) );
644 F = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,1) );
645 GMX_MM_TRANSPOSE2_PD(Y,F);
646 G = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) +2);
647 H = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,1) +2);
648 GMX_MM_TRANSPOSE2_PD(G,H);
649 Heps = _mm_mul_pd(vfeps,H);
650 Fp = _mm_add_pd(F,_mm_mul_pd(vfeps,_mm_add_pd(G,Heps)));
651 FF = _mm_add_pd(Fp,_mm_mul_pd(vfeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
652 fvdw6 = _mm_mul_pd(c6_00,FF);
654 /* CUBIC SPLINE TABLE REPULSION */
655 vfitab = _mm_add_epi32(vfitab,ifour);
656 Y = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) );
657 F = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,1) );
658 GMX_MM_TRANSPOSE2_PD(Y,F);
659 G = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) +2);
660 H = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,1) +2);
661 GMX_MM_TRANSPOSE2_PD(G,H);
662 Heps = _mm_mul_pd(vfeps,H);
663 Fp = _mm_add_pd(F,_mm_mul_pd(vfeps,_mm_add_pd(G,Heps)));
664 FF = _mm_add_pd(Fp,_mm_mul_pd(vfeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
665 fvdw12 = _mm_mul_pd(c12_00,FF);
666 fvdw = _mm_xor_pd(signbit,_mm_mul_pd(_mm_add_pd(fvdw6,fvdw12),_mm_mul_pd(vftabscale,rinv00)));
668 fscal = _mm_add_pd(felec,fvdw);
670 /* Calculate temporary vectorial force */
671 tx = _mm_mul_pd(fscal,dx00);
672 ty = _mm_mul_pd(fscal,dy00);
673 tz = _mm_mul_pd(fscal,dz00);
675 /* Update vectorial force */
676 fix0 = _mm_add_pd(fix0,tx);
677 fiy0 = _mm_add_pd(fiy0,ty);
678 fiz0 = _mm_add_pd(fiz0,tz);
680 gmx_mm_decrement_1rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,tx,ty,tz);
682 /* Inner loop uses 82 flops */
689 j_coord_offsetA = DIM*jnrA;
691 /* load j atom coordinates */
692 gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
695 /* Calculate displacement vector */
696 dx00 = _mm_sub_pd(ix0,jx0);
697 dy00 = _mm_sub_pd(iy0,jy0);
698 dz00 = _mm_sub_pd(iz0,jz0);
700 /* Calculate squared distance and things based on it */
701 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
703 rinv00 = gmx_mm_invsqrt_pd(rsq00);
705 /* Load parameters for j particles */
706 jq0 = _mm_load_sd(charge+jnrA+0);
707 isaj0 = _mm_load_sd(invsqrta+jnrA+0);
708 vdwjidx0A = 2*vdwtype[jnrA+0];
710 /**************************
711 * CALCULATE INTERACTIONS *
712 **************************/
714 r00 = _mm_mul_pd(rsq00,rinv00);
716 /* Compute parameters for interactions between i and j atoms */
717 qq00 = _mm_mul_pd(iq0,jq0);
718 gmx_mm_load_1pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,&c6_00,&c12_00);
720 /* Calculate table index by multiplying r with table scale and truncate to integer */
721 rt = _mm_mul_pd(r00,vftabscale);
722 vfitab = _mm_cvttpd_epi32(rt);
723 vfeps = _mm_sub_pd(rt,_mm_cvtepi32_pd(vfitab));
724 vfitab = _mm_slli_epi32(vfitab,3);
726 /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
727 isaprod = _mm_mul_pd(isai0,isaj0);
728 gbqqfactor = _mm_xor_pd(signbit,_mm_mul_pd(qq00,_mm_mul_pd(isaprod,gbinvepsdiff)));
729 gbscale = _mm_mul_pd(isaprod,gbtabscale);
731 /* Calculate generalized born table index - this is a separate table from the normal one,
732 * but we use the same procedure by multiplying r with scale and truncating to integer.
734 rt = _mm_mul_pd(r00,gbscale);
735 gbitab = _mm_cvttpd_epi32(rt);
736 gbeps = _mm_sub_pd(rt,_mm_cvtepi32_pd(gbitab));
737 gbitab = _mm_slli_epi32(gbitab,2);
739 Y = _mm_load_pd( gbtab + gmx_mm_extract_epi32(gbitab,0) );
740 F = _mm_setzero_pd();
741 GMX_MM_TRANSPOSE2_PD(Y,F);
742 G = _mm_load_pd( gbtab + gmx_mm_extract_epi32(gbitab,0) +2);
743 H = _mm_setzero_pd();
744 GMX_MM_TRANSPOSE2_PD(G,H);
745 Heps = _mm_mul_pd(gbeps,H);
746 Fp = _mm_add_pd(F,_mm_mul_pd(gbeps,_mm_add_pd(G,Heps)));
747 VV = _mm_add_pd(Y,_mm_mul_pd(gbeps,Fp));
748 vgb = _mm_mul_pd(gbqqfactor,VV);
750 FF = _mm_add_pd(Fp,_mm_mul_pd(gbeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
751 fgb = _mm_mul_pd(gbqqfactor,_mm_mul_pd(FF,gbscale));
752 dvdatmp = _mm_mul_pd(minushalf,_mm_add_pd(vgb,_mm_mul_pd(fgb,r00)));
753 dvdasum = _mm_add_pd(dvdasum,dvdatmp);
754 gmx_mm_increment_1real_pd(dvda+jnrA,_mm_mul_pd(dvdatmp,_mm_mul_pd(isaj0,isaj0)));
755 velec = _mm_mul_pd(qq00,rinv00);
756 felec = _mm_mul_pd(_mm_sub_pd(_mm_mul_pd(velec,rinv00),fgb),rinv00);
758 /* CUBIC SPLINE TABLE DISPERSION */
759 Y = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) );
760 F = _mm_setzero_pd();
761 GMX_MM_TRANSPOSE2_PD(Y,F);
762 G = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) +2);
763 H = _mm_setzero_pd();
764 GMX_MM_TRANSPOSE2_PD(G,H);
765 Heps = _mm_mul_pd(vfeps,H);
766 Fp = _mm_add_pd(F,_mm_mul_pd(vfeps,_mm_add_pd(G,Heps)));
767 FF = _mm_add_pd(Fp,_mm_mul_pd(vfeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
768 fvdw6 = _mm_mul_pd(c6_00,FF);
770 /* CUBIC SPLINE TABLE REPULSION */
771 vfitab = _mm_add_epi32(vfitab,ifour);
772 Y = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) );
773 F = _mm_setzero_pd();
774 GMX_MM_TRANSPOSE2_PD(Y,F);
775 G = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) +2);
776 H = _mm_setzero_pd();
777 GMX_MM_TRANSPOSE2_PD(G,H);
778 Heps = _mm_mul_pd(vfeps,H);
779 Fp = _mm_add_pd(F,_mm_mul_pd(vfeps,_mm_add_pd(G,Heps)));
780 FF = _mm_add_pd(Fp,_mm_mul_pd(vfeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
781 fvdw12 = _mm_mul_pd(c12_00,FF);
782 fvdw = _mm_xor_pd(signbit,_mm_mul_pd(_mm_add_pd(fvdw6,fvdw12),_mm_mul_pd(vftabscale,rinv00)));
784 fscal = _mm_add_pd(felec,fvdw);
786 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
788 /* Calculate temporary vectorial force */
789 tx = _mm_mul_pd(fscal,dx00);
790 ty = _mm_mul_pd(fscal,dy00);
791 tz = _mm_mul_pd(fscal,dz00);
793 /* Update vectorial force */
794 fix0 = _mm_add_pd(fix0,tx);
795 fiy0 = _mm_add_pd(fiy0,ty);
796 fiz0 = _mm_add_pd(fiz0,tz);
798 gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,tx,ty,tz);
800 /* Inner loop uses 82 flops */
803 /* End of innermost loop */
805 gmx_mm_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
806 f+i_coord_offset,fshift+i_shift_offset);
808 dvdasum = _mm_mul_pd(dvdasum, _mm_mul_pd(isai0,isai0));
809 gmx_mm_update_1pot_pd(dvdasum,dvda+inr);
811 /* Increment number of inner iterations */
812 inneriter += j_index_end - j_index_start;
814 /* Outer loop uses 7 flops */
817 /* Increment number of outer iterations */
820 /* Update outer/inner flops */
822 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_F,outeriter*7 + inneriter*82);