2 * Note: this file was generated by the Gromacs sse4_1_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_sse4_1_double.h"
34 #include "kernelutil_x86_sse4_1_double.h"
37 * Gromacs nonbonded kernel: nb_kernel_ElecRFCut_VdwCSTab_GeomP1P1_VF_sse4_1_double
38 * Electrostatics interaction: ReactionField
39 * VdW interaction: CubicSplineTable
40 * Geometry: Particle-Particle
41 * Calculate force/pot: PotentialAndForce
44 nb_kernel_ElecRFCut_VdwCSTab_GeomP1P1_VF_sse4_1_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 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
77 __m128d one_sixth = _mm_set1_pd(1.0/6.0);
78 __m128d one_twelfth = _mm_set1_pd(1.0/12.0);
80 __m128i ifour = _mm_set1_epi32(4);
81 __m128d rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF;
83 __m128d dummy_mask,cutoff_mask;
84 __m128d signbit = gmx_mm_castsi128_pd( _mm_set_epi32(0x80000000,0x00000000,0x80000000,0x00000000) );
85 __m128d one = _mm_set1_pd(1.0);
86 __m128d two = _mm_set1_pd(2.0);
92 jindex = nlist->jindex;
94 shiftidx = nlist->shift;
96 shiftvec = fr->shift_vec[0];
97 fshift = fr->fshift[0];
98 facel = _mm_set1_pd(fr->epsfac);
99 charge = mdatoms->chargeA;
100 krf = _mm_set1_pd(fr->ic->k_rf);
101 krf2 = _mm_set1_pd(fr->ic->k_rf*2.0);
102 crf = _mm_set1_pd(fr->ic->c_rf);
103 nvdwtype = fr->ntype;
105 vdwtype = mdatoms->typeA;
107 vftab = kernel_data->table_vdw->data;
108 vftabscale = _mm_set1_pd(kernel_data->table_vdw->scale);
110 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
111 rcutoff_scalar = fr->rcoulomb;
112 rcutoff = _mm_set1_pd(rcutoff_scalar);
113 rcutoff2 = _mm_mul_pd(rcutoff,rcutoff);
115 /* Avoid stupid compiler warnings */
123 /* Start outer loop over neighborlists */
124 for(iidx=0; iidx<nri; iidx++)
126 /* Load shift vector for this list */
127 i_shift_offset = DIM*shiftidx[iidx];
129 /* Load limits for loop over neighbors */
130 j_index_start = jindex[iidx];
131 j_index_end = jindex[iidx+1];
133 /* Get outer coordinate index */
135 i_coord_offset = DIM*inr;
137 /* Load i particle coords and add shift vector */
138 gmx_mm_load_shift_and_1rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
140 fix0 = _mm_setzero_pd();
141 fiy0 = _mm_setzero_pd();
142 fiz0 = _mm_setzero_pd();
144 /* Load parameters for i particles */
145 iq0 = _mm_mul_pd(facel,_mm_load1_pd(charge+inr+0));
146 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
148 /* Reset potential sums */
149 velecsum = _mm_setzero_pd();
150 vvdwsum = _mm_setzero_pd();
152 /* Start inner kernel loop */
153 for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
156 /* Get j neighbor index, and coordinate index */
159 j_coord_offsetA = DIM*jnrA;
160 j_coord_offsetB = DIM*jnrB;
162 /* load j atom coordinates */
163 gmx_mm_load_1rvec_2ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
166 /* Calculate displacement vector */
167 dx00 = _mm_sub_pd(ix0,jx0);
168 dy00 = _mm_sub_pd(iy0,jy0);
169 dz00 = _mm_sub_pd(iz0,jz0);
171 /* Calculate squared distance and things based on it */
172 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
174 rinv00 = gmx_mm_invsqrt_pd(rsq00);
176 rinvsq00 = _mm_mul_pd(rinv00,rinv00);
178 /* Load parameters for j particles */
179 jq0 = gmx_mm_load_2real_swizzle_pd(charge+jnrA+0,charge+jnrB+0);
180 vdwjidx0A = 2*vdwtype[jnrA+0];
181 vdwjidx0B = 2*vdwtype[jnrB+0];
183 /**************************
184 * CALCULATE INTERACTIONS *
185 **************************/
187 if (gmx_mm_any_lt(rsq00,rcutoff2))
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 /* Calculate table index by multiplying r with table scale and truncate to integer */
198 rt = _mm_mul_pd(r00,vftabscale);
199 vfitab = _mm_cvttpd_epi32(rt);
200 vfeps = _mm_sub_pd(rt,_mm_round_pd(rt, _MM_FROUND_FLOOR));
201 vfitab = _mm_slli_epi32(vfitab,3);
203 /* REACTION-FIELD ELECTROSTATICS */
204 velec = _mm_mul_pd(qq00,_mm_sub_pd(_mm_add_pd(rinv00,_mm_mul_pd(krf,rsq00)),crf));
205 felec = _mm_mul_pd(qq00,_mm_sub_pd(_mm_mul_pd(rinv00,rinvsq00),krf2));
207 /* CUBIC SPLINE TABLE DISPERSION */
208 Y = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) );
209 F = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,1) );
210 GMX_MM_TRANSPOSE2_PD(Y,F);
211 G = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) +2);
212 H = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,1) +2);
213 GMX_MM_TRANSPOSE2_PD(G,H);
214 Heps = _mm_mul_pd(vfeps,H);
215 Fp = _mm_add_pd(F,_mm_mul_pd(vfeps,_mm_add_pd(G,Heps)));
216 VV = _mm_add_pd(Y,_mm_mul_pd(vfeps,Fp));
217 vvdw6 = _mm_mul_pd(c6_00,VV);
218 FF = _mm_add_pd(Fp,_mm_mul_pd(vfeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
219 fvdw6 = _mm_mul_pd(c6_00,FF);
221 /* CUBIC SPLINE TABLE REPULSION */
222 vfitab = _mm_add_epi32(vfitab,ifour);
223 Y = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) );
224 F = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,1) );
225 GMX_MM_TRANSPOSE2_PD(Y,F);
226 G = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) +2);
227 H = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,1) +2);
228 GMX_MM_TRANSPOSE2_PD(G,H);
229 Heps = _mm_mul_pd(vfeps,H);
230 Fp = _mm_add_pd(F,_mm_mul_pd(vfeps,_mm_add_pd(G,Heps)));
231 VV = _mm_add_pd(Y,_mm_mul_pd(vfeps,Fp));
232 vvdw12 = _mm_mul_pd(c12_00,VV);
233 FF = _mm_add_pd(Fp,_mm_mul_pd(vfeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
234 fvdw12 = _mm_mul_pd(c12_00,FF);
235 vvdw = _mm_add_pd(vvdw12,vvdw6);
236 fvdw = _mm_xor_pd(signbit,_mm_mul_pd(_mm_add_pd(fvdw6,fvdw12),_mm_mul_pd(vftabscale,rinv00)));
238 cutoff_mask = _mm_cmplt_pd(rsq00,rcutoff2);
240 /* Update potential sum for this i atom from the interaction with this j atom. */
241 velec = _mm_and_pd(velec,cutoff_mask);
242 velecsum = _mm_add_pd(velecsum,velec);
243 vvdw = _mm_and_pd(vvdw,cutoff_mask);
244 vvdwsum = _mm_add_pd(vvdwsum,vvdw);
246 fscal = _mm_add_pd(felec,fvdw);
248 fscal = _mm_and_pd(fscal,cutoff_mask);
250 /* Calculate temporary vectorial force */
251 tx = _mm_mul_pd(fscal,dx00);
252 ty = _mm_mul_pd(fscal,dy00);
253 tz = _mm_mul_pd(fscal,dz00);
255 /* Update vectorial force */
256 fix0 = _mm_add_pd(fix0,tx);
257 fiy0 = _mm_add_pd(fiy0,ty);
258 fiz0 = _mm_add_pd(fiz0,tz);
260 gmx_mm_decrement_1rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,tx,ty,tz);
264 /* Inner loop uses 72 flops */
271 j_coord_offsetA = DIM*jnrA;
273 /* load j atom coordinates */
274 gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
277 /* Calculate displacement vector */
278 dx00 = _mm_sub_pd(ix0,jx0);
279 dy00 = _mm_sub_pd(iy0,jy0);
280 dz00 = _mm_sub_pd(iz0,jz0);
282 /* Calculate squared distance and things based on it */
283 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
285 rinv00 = gmx_mm_invsqrt_pd(rsq00);
287 rinvsq00 = _mm_mul_pd(rinv00,rinv00);
289 /* Load parameters for j particles */
290 jq0 = _mm_load_sd(charge+jnrA+0);
291 vdwjidx0A = 2*vdwtype[jnrA+0];
293 /**************************
294 * CALCULATE INTERACTIONS *
295 **************************/
297 if (gmx_mm_any_lt(rsq00,rcutoff2))
300 r00 = _mm_mul_pd(rsq00,rinv00);
302 /* Compute parameters for interactions between i and j atoms */
303 qq00 = _mm_mul_pd(iq0,jq0);
304 gmx_mm_load_1pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,&c6_00,&c12_00);
306 /* Calculate table index by multiplying r with table scale and truncate to integer */
307 rt = _mm_mul_pd(r00,vftabscale);
308 vfitab = _mm_cvttpd_epi32(rt);
309 vfeps = _mm_sub_pd(rt,_mm_round_pd(rt, _MM_FROUND_FLOOR));
310 vfitab = _mm_slli_epi32(vfitab,3);
312 /* REACTION-FIELD ELECTROSTATICS */
313 velec = _mm_mul_pd(qq00,_mm_sub_pd(_mm_add_pd(rinv00,_mm_mul_pd(krf,rsq00)),crf));
314 felec = _mm_mul_pd(qq00,_mm_sub_pd(_mm_mul_pd(rinv00,rinvsq00),krf2));
316 /* CUBIC SPLINE TABLE DISPERSION */
317 Y = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) );
318 F = _mm_setzero_pd();
319 GMX_MM_TRANSPOSE2_PD(Y,F);
320 G = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) +2);
321 H = _mm_setzero_pd();
322 GMX_MM_TRANSPOSE2_PD(G,H);
323 Heps = _mm_mul_pd(vfeps,H);
324 Fp = _mm_add_pd(F,_mm_mul_pd(vfeps,_mm_add_pd(G,Heps)));
325 VV = _mm_add_pd(Y,_mm_mul_pd(vfeps,Fp));
326 vvdw6 = _mm_mul_pd(c6_00,VV);
327 FF = _mm_add_pd(Fp,_mm_mul_pd(vfeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
328 fvdw6 = _mm_mul_pd(c6_00,FF);
330 /* CUBIC SPLINE TABLE REPULSION */
331 vfitab = _mm_add_epi32(vfitab,ifour);
332 Y = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) );
333 F = _mm_setzero_pd();
334 GMX_MM_TRANSPOSE2_PD(Y,F);
335 G = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) +2);
336 H = _mm_setzero_pd();
337 GMX_MM_TRANSPOSE2_PD(G,H);
338 Heps = _mm_mul_pd(vfeps,H);
339 Fp = _mm_add_pd(F,_mm_mul_pd(vfeps,_mm_add_pd(G,Heps)));
340 VV = _mm_add_pd(Y,_mm_mul_pd(vfeps,Fp));
341 vvdw12 = _mm_mul_pd(c12_00,VV);
342 FF = _mm_add_pd(Fp,_mm_mul_pd(vfeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
343 fvdw12 = _mm_mul_pd(c12_00,FF);
344 vvdw = _mm_add_pd(vvdw12,vvdw6);
345 fvdw = _mm_xor_pd(signbit,_mm_mul_pd(_mm_add_pd(fvdw6,fvdw12),_mm_mul_pd(vftabscale,rinv00)));
347 cutoff_mask = _mm_cmplt_pd(rsq00,rcutoff2);
349 /* Update potential sum for this i atom from the interaction with this j atom. */
350 velec = _mm_and_pd(velec,cutoff_mask);
351 velec = _mm_unpacklo_pd(velec,_mm_setzero_pd());
352 velecsum = _mm_add_pd(velecsum,velec);
353 vvdw = _mm_and_pd(vvdw,cutoff_mask);
354 vvdw = _mm_unpacklo_pd(vvdw,_mm_setzero_pd());
355 vvdwsum = _mm_add_pd(vvdwsum,vvdw);
357 fscal = _mm_add_pd(felec,fvdw);
359 fscal = _mm_and_pd(fscal,cutoff_mask);
361 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
363 /* Calculate temporary vectorial force */
364 tx = _mm_mul_pd(fscal,dx00);
365 ty = _mm_mul_pd(fscal,dy00);
366 tz = _mm_mul_pd(fscal,dz00);
368 /* Update vectorial force */
369 fix0 = _mm_add_pd(fix0,tx);
370 fiy0 = _mm_add_pd(fiy0,ty);
371 fiz0 = _mm_add_pd(fiz0,tz);
373 gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,tx,ty,tz);
377 /* Inner loop uses 72 flops */
380 /* End of innermost loop */
382 gmx_mm_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
383 f+i_coord_offset,fshift+i_shift_offset);
386 /* Update potential energies */
387 gmx_mm_update_1pot_pd(velecsum,kernel_data->energygrp_elec+ggid);
388 gmx_mm_update_1pot_pd(vvdwsum,kernel_data->energygrp_vdw+ggid);
390 /* Increment number of inner iterations */
391 inneriter += j_index_end - j_index_start;
393 /* Outer loop uses 9 flops */
396 /* Increment number of outer iterations */
399 /* Update outer/inner flops */
401 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_VF,outeriter*9 + inneriter*72);
404 * Gromacs nonbonded kernel: nb_kernel_ElecRFCut_VdwCSTab_GeomP1P1_F_sse4_1_double
405 * Electrostatics interaction: ReactionField
406 * VdW interaction: CubicSplineTable
407 * Geometry: Particle-Particle
408 * Calculate force/pot: Force
411 nb_kernel_ElecRFCut_VdwCSTab_GeomP1P1_F_sse4_1_double
412 (t_nblist * gmx_restrict nlist,
413 rvec * gmx_restrict xx,
414 rvec * gmx_restrict ff,
415 t_forcerec * gmx_restrict fr,
416 t_mdatoms * gmx_restrict mdatoms,
417 nb_kernel_data_t * gmx_restrict kernel_data,
418 t_nrnb * gmx_restrict nrnb)
420 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
421 * just 0 for non-waters.
422 * Suffixes A,B refer to j loop unrolling done with SSE double precision, e.g. for the two different
423 * jnr indices corresponding to data put in the four positions in the SIMD register.
425 int i_shift_offset,i_coord_offset,outeriter,inneriter;
426 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
428 int j_coord_offsetA,j_coord_offsetB;
429 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
431 real *shiftvec,*fshift,*x,*f;
432 __m128d tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
434 __m128d ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
435 int vdwjidx0A,vdwjidx0B;
436 __m128d jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
437 __m128d dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
438 __m128d velec,felec,velecsum,facel,crf,krf,krf2;
441 __m128d rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
444 __m128d one_sixth = _mm_set1_pd(1.0/6.0);
445 __m128d one_twelfth = _mm_set1_pd(1.0/12.0);
447 __m128i ifour = _mm_set1_epi32(4);
448 __m128d rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF;
450 __m128d dummy_mask,cutoff_mask;
451 __m128d signbit = gmx_mm_castsi128_pd( _mm_set_epi32(0x80000000,0x00000000,0x80000000,0x00000000) );
452 __m128d one = _mm_set1_pd(1.0);
453 __m128d two = _mm_set1_pd(2.0);
459 jindex = nlist->jindex;
461 shiftidx = nlist->shift;
463 shiftvec = fr->shift_vec[0];
464 fshift = fr->fshift[0];
465 facel = _mm_set1_pd(fr->epsfac);
466 charge = mdatoms->chargeA;
467 krf = _mm_set1_pd(fr->ic->k_rf);
468 krf2 = _mm_set1_pd(fr->ic->k_rf*2.0);
469 crf = _mm_set1_pd(fr->ic->c_rf);
470 nvdwtype = fr->ntype;
472 vdwtype = mdatoms->typeA;
474 vftab = kernel_data->table_vdw->data;
475 vftabscale = _mm_set1_pd(kernel_data->table_vdw->scale);
477 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
478 rcutoff_scalar = fr->rcoulomb;
479 rcutoff = _mm_set1_pd(rcutoff_scalar);
480 rcutoff2 = _mm_mul_pd(rcutoff,rcutoff);
482 /* Avoid stupid compiler warnings */
490 /* Start outer loop over neighborlists */
491 for(iidx=0; iidx<nri; iidx++)
493 /* Load shift vector for this list */
494 i_shift_offset = DIM*shiftidx[iidx];
496 /* Load limits for loop over neighbors */
497 j_index_start = jindex[iidx];
498 j_index_end = jindex[iidx+1];
500 /* Get outer coordinate index */
502 i_coord_offset = DIM*inr;
504 /* Load i particle coords and add shift vector */
505 gmx_mm_load_shift_and_1rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
507 fix0 = _mm_setzero_pd();
508 fiy0 = _mm_setzero_pd();
509 fiz0 = _mm_setzero_pd();
511 /* Load parameters for i particles */
512 iq0 = _mm_mul_pd(facel,_mm_load1_pd(charge+inr+0));
513 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
515 /* Start inner kernel loop */
516 for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
519 /* Get j neighbor index, and coordinate index */
522 j_coord_offsetA = DIM*jnrA;
523 j_coord_offsetB = DIM*jnrB;
525 /* load j atom coordinates */
526 gmx_mm_load_1rvec_2ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
529 /* Calculate displacement vector */
530 dx00 = _mm_sub_pd(ix0,jx0);
531 dy00 = _mm_sub_pd(iy0,jy0);
532 dz00 = _mm_sub_pd(iz0,jz0);
534 /* Calculate squared distance and things based on it */
535 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
537 rinv00 = gmx_mm_invsqrt_pd(rsq00);
539 rinvsq00 = _mm_mul_pd(rinv00,rinv00);
541 /* Load parameters for j particles */
542 jq0 = gmx_mm_load_2real_swizzle_pd(charge+jnrA+0,charge+jnrB+0);
543 vdwjidx0A = 2*vdwtype[jnrA+0];
544 vdwjidx0B = 2*vdwtype[jnrB+0];
546 /**************************
547 * CALCULATE INTERACTIONS *
548 **************************/
550 if (gmx_mm_any_lt(rsq00,rcutoff2))
553 r00 = _mm_mul_pd(rsq00,rinv00);
555 /* Compute parameters for interactions between i and j atoms */
556 qq00 = _mm_mul_pd(iq0,jq0);
557 gmx_mm_load_2pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,
558 vdwparam+vdwioffset0+vdwjidx0B,&c6_00,&c12_00);
560 /* Calculate table index by multiplying r with table scale and truncate to integer */
561 rt = _mm_mul_pd(r00,vftabscale);
562 vfitab = _mm_cvttpd_epi32(rt);
563 vfeps = _mm_sub_pd(rt,_mm_round_pd(rt, _MM_FROUND_FLOOR));
564 vfitab = _mm_slli_epi32(vfitab,3);
566 /* REACTION-FIELD ELECTROSTATICS */
567 felec = _mm_mul_pd(qq00,_mm_sub_pd(_mm_mul_pd(rinv00,rinvsq00),krf2));
569 /* CUBIC SPLINE TABLE DISPERSION */
570 Y = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) );
571 F = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,1) );
572 GMX_MM_TRANSPOSE2_PD(Y,F);
573 G = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) +2);
574 H = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,1) +2);
575 GMX_MM_TRANSPOSE2_PD(G,H);
576 Heps = _mm_mul_pd(vfeps,H);
577 Fp = _mm_add_pd(F,_mm_mul_pd(vfeps,_mm_add_pd(G,Heps)));
578 FF = _mm_add_pd(Fp,_mm_mul_pd(vfeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
579 fvdw6 = _mm_mul_pd(c6_00,FF);
581 /* CUBIC SPLINE TABLE REPULSION */
582 vfitab = _mm_add_epi32(vfitab,ifour);
583 Y = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) );
584 F = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,1) );
585 GMX_MM_TRANSPOSE2_PD(Y,F);
586 G = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) +2);
587 H = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,1) +2);
588 GMX_MM_TRANSPOSE2_PD(G,H);
589 Heps = _mm_mul_pd(vfeps,H);
590 Fp = _mm_add_pd(F,_mm_mul_pd(vfeps,_mm_add_pd(G,Heps)));
591 FF = _mm_add_pd(Fp,_mm_mul_pd(vfeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
592 fvdw12 = _mm_mul_pd(c12_00,FF);
593 fvdw = _mm_xor_pd(signbit,_mm_mul_pd(_mm_add_pd(fvdw6,fvdw12),_mm_mul_pd(vftabscale,rinv00)));
595 cutoff_mask = _mm_cmplt_pd(rsq00,rcutoff2);
597 fscal = _mm_add_pd(felec,fvdw);
599 fscal = _mm_and_pd(fscal,cutoff_mask);
601 /* Calculate temporary vectorial force */
602 tx = _mm_mul_pd(fscal,dx00);
603 ty = _mm_mul_pd(fscal,dy00);
604 tz = _mm_mul_pd(fscal,dz00);
606 /* Update vectorial force */
607 fix0 = _mm_add_pd(fix0,tx);
608 fiy0 = _mm_add_pd(fiy0,ty);
609 fiz0 = _mm_add_pd(fiz0,tz);
611 gmx_mm_decrement_1rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,tx,ty,tz);
615 /* Inner loop uses 57 flops */
622 j_coord_offsetA = DIM*jnrA;
624 /* load j atom coordinates */
625 gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
628 /* Calculate displacement vector */
629 dx00 = _mm_sub_pd(ix0,jx0);
630 dy00 = _mm_sub_pd(iy0,jy0);
631 dz00 = _mm_sub_pd(iz0,jz0);
633 /* Calculate squared distance and things based on it */
634 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
636 rinv00 = gmx_mm_invsqrt_pd(rsq00);
638 rinvsq00 = _mm_mul_pd(rinv00,rinv00);
640 /* Load parameters for j particles */
641 jq0 = _mm_load_sd(charge+jnrA+0);
642 vdwjidx0A = 2*vdwtype[jnrA+0];
644 /**************************
645 * CALCULATE INTERACTIONS *
646 **************************/
648 if (gmx_mm_any_lt(rsq00,rcutoff2))
651 r00 = _mm_mul_pd(rsq00,rinv00);
653 /* Compute parameters for interactions between i and j atoms */
654 qq00 = _mm_mul_pd(iq0,jq0);
655 gmx_mm_load_1pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,&c6_00,&c12_00);
657 /* Calculate table index by multiplying r with table scale and truncate to integer */
658 rt = _mm_mul_pd(r00,vftabscale);
659 vfitab = _mm_cvttpd_epi32(rt);
660 vfeps = _mm_sub_pd(rt,_mm_round_pd(rt, _MM_FROUND_FLOOR));
661 vfitab = _mm_slli_epi32(vfitab,3);
663 /* REACTION-FIELD ELECTROSTATICS */
664 felec = _mm_mul_pd(qq00,_mm_sub_pd(_mm_mul_pd(rinv00,rinvsq00),krf2));
666 /* CUBIC SPLINE TABLE DISPERSION */
667 Y = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) );
668 F = _mm_setzero_pd();
669 GMX_MM_TRANSPOSE2_PD(Y,F);
670 G = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) +2);
671 H = _mm_setzero_pd();
672 GMX_MM_TRANSPOSE2_PD(G,H);
673 Heps = _mm_mul_pd(vfeps,H);
674 Fp = _mm_add_pd(F,_mm_mul_pd(vfeps,_mm_add_pd(G,Heps)));
675 FF = _mm_add_pd(Fp,_mm_mul_pd(vfeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
676 fvdw6 = _mm_mul_pd(c6_00,FF);
678 /* CUBIC SPLINE TABLE REPULSION */
679 vfitab = _mm_add_epi32(vfitab,ifour);
680 Y = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) );
681 F = _mm_setzero_pd();
682 GMX_MM_TRANSPOSE2_PD(Y,F);
683 G = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) +2);
684 H = _mm_setzero_pd();
685 GMX_MM_TRANSPOSE2_PD(G,H);
686 Heps = _mm_mul_pd(vfeps,H);
687 Fp = _mm_add_pd(F,_mm_mul_pd(vfeps,_mm_add_pd(G,Heps)));
688 FF = _mm_add_pd(Fp,_mm_mul_pd(vfeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
689 fvdw12 = _mm_mul_pd(c12_00,FF);
690 fvdw = _mm_xor_pd(signbit,_mm_mul_pd(_mm_add_pd(fvdw6,fvdw12),_mm_mul_pd(vftabscale,rinv00)));
692 cutoff_mask = _mm_cmplt_pd(rsq00,rcutoff2);
694 fscal = _mm_add_pd(felec,fvdw);
696 fscal = _mm_and_pd(fscal,cutoff_mask);
698 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
700 /* Calculate temporary vectorial force */
701 tx = _mm_mul_pd(fscal,dx00);
702 ty = _mm_mul_pd(fscal,dy00);
703 tz = _mm_mul_pd(fscal,dz00);
705 /* Update vectorial force */
706 fix0 = _mm_add_pd(fix0,tx);
707 fiy0 = _mm_add_pd(fiy0,ty);
708 fiz0 = _mm_add_pd(fiz0,tz);
710 gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,tx,ty,tz);
714 /* Inner loop uses 57 flops */
717 /* End of innermost loop */
719 gmx_mm_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
720 f+i_coord_offset,fshift+i_shift_offset);
722 /* Increment number of inner iterations */
723 inneriter += j_index_end - j_index_start;
725 /* Outer loop uses 7 flops */
728 /* Increment number of outer iterations */
731 /* Update outer/inner flops */
733 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_F,outeriter*7 + inneriter*57);