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_ElecCSTab_VdwNone_GeomP1P1_VF_avx_128_fma_double
38 * Electrostatics interaction: CubicSplineTable
39 * VdW interaction: None
40 * Geometry: Particle-Particle
41 * Calculate force/pot: PotentialAndForce
44 nb_kernel_ElecCSTab_VdwNone_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 __m128i ifour = _mm_set1_epi32(4);
75 __m128d rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF,twovfeps;
77 __m128d dummy_mask,cutoff_mask;
78 __m128d signbit = gmx_mm_castsi128_pd( _mm_set_epi32(0x80000000,0x00000000,0x80000000,0x00000000) );
79 __m128d one = _mm_set1_pd(1.0);
80 __m128d two = _mm_set1_pd(2.0);
86 jindex = nlist->jindex;
88 shiftidx = nlist->shift;
90 shiftvec = fr->shift_vec[0];
91 fshift = fr->fshift[0];
92 facel = _mm_set1_pd(fr->epsfac);
93 charge = mdatoms->chargeA;
95 vftab = kernel_data->table_elec->data;
96 vftabscale = _mm_set1_pd(kernel_data->table_elec->scale);
98 /* Avoid stupid compiler warnings */
106 /* Start outer loop over neighborlists */
107 for(iidx=0; iidx<nri; iidx++)
109 /* Load shift vector for this list */
110 i_shift_offset = DIM*shiftidx[iidx];
112 /* Load limits for loop over neighbors */
113 j_index_start = jindex[iidx];
114 j_index_end = jindex[iidx+1];
116 /* Get outer coordinate index */
118 i_coord_offset = DIM*inr;
120 /* Load i particle coords and add shift vector */
121 gmx_mm_load_shift_and_1rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
123 fix0 = _mm_setzero_pd();
124 fiy0 = _mm_setzero_pd();
125 fiz0 = _mm_setzero_pd();
127 /* Load parameters for i particles */
128 iq0 = _mm_mul_pd(facel,_mm_load1_pd(charge+inr+0));
130 /* Reset potential sums */
131 velecsum = _mm_setzero_pd();
133 /* Start inner kernel loop */
134 for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
137 /* Get j neighbor index, and coordinate index */
140 j_coord_offsetA = DIM*jnrA;
141 j_coord_offsetB = DIM*jnrB;
143 /* load j atom coordinates */
144 gmx_mm_load_1rvec_2ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
147 /* Calculate displacement vector */
148 dx00 = _mm_sub_pd(ix0,jx0);
149 dy00 = _mm_sub_pd(iy0,jy0);
150 dz00 = _mm_sub_pd(iz0,jz0);
152 /* Calculate squared distance and things based on it */
153 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
155 rinv00 = gmx_mm_invsqrt_pd(rsq00);
157 /* Load parameters for j particles */
158 jq0 = gmx_mm_load_2real_swizzle_pd(charge+jnrA+0,charge+jnrB+0);
160 /**************************
161 * CALCULATE INTERACTIONS *
162 **************************/
164 r00 = _mm_mul_pd(rsq00,rinv00);
166 /* Compute parameters for interactions between i and j atoms */
167 qq00 = _mm_mul_pd(iq0,jq0);
169 /* Calculate table index by multiplying r with table scale and truncate to integer */
170 rt = _mm_mul_pd(r00,vftabscale);
171 vfitab = _mm_cvttpd_epi32(rt);
173 vfeps = _mm_frcz_pd(rt);
175 vfeps = _mm_sub_pd(rt,_mm_round_pd(rt, _MM_FROUND_FLOOR));
177 twovfeps = _mm_add_pd(vfeps,vfeps);
178 vfitab = _mm_slli_epi32(vfitab,2);
180 /* CUBIC SPLINE TABLE ELECTROSTATICS */
181 Y = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) );
182 F = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,1) );
183 GMX_MM_TRANSPOSE2_PD(Y,F);
184 G = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) +2);
185 H = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,1) +2);
186 GMX_MM_TRANSPOSE2_PD(G,H);
187 Fp = _mm_macc_pd(vfeps,_mm_macc_pd(vfeps,H,G),F);
188 VV = _mm_macc_pd(vfeps,Fp,Y);
189 velec = _mm_mul_pd(qq00,VV);
190 FF = _mm_macc_pd(_mm_macc_pd(twovfeps,H,G),vfeps,Fp);
191 felec = _mm_xor_pd(signbit,_mm_mul_pd(_mm_mul_pd(qq00,FF),_mm_mul_pd(vftabscale,rinv00)));
193 /* Update potential sum for this i atom from the interaction with this j atom. */
194 velecsum = _mm_add_pd(velecsum,velec);
198 /* Update vectorial force */
199 fix0 = _mm_macc_pd(dx00,fscal,fix0);
200 fiy0 = _mm_macc_pd(dy00,fscal,fiy0);
201 fiz0 = _mm_macc_pd(dz00,fscal,fiz0);
203 gmx_mm_decrement_1rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,
204 _mm_mul_pd(dx00,fscal),
205 _mm_mul_pd(dy00,fscal),
206 _mm_mul_pd(dz00,fscal));
208 /* Inner loop uses 46 flops */
215 j_coord_offsetA = DIM*jnrA;
217 /* load j atom coordinates */
218 gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
221 /* Calculate displacement vector */
222 dx00 = _mm_sub_pd(ix0,jx0);
223 dy00 = _mm_sub_pd(iy0,jy0);
224 dz00 = _mm_sub_pd(iz0,jz0);
226 /* Calculate squared distance and things based on it */
227 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
229 rinv00 = gmx_mm_invsqrt_pd(rsq00);
231 /* Load parameters for j particles */
232 jq0 = _mm_load_sd(charge+jnrA+0);
234 /**************************
235 * CALCULATE INTERACTIONS *
236 **************************/
238 r00 = _mm_mul_pd(rsq00,rinv00);
240 /* Compute parameters for interactions between i and j atoms */
241 qq00 = _mm_mul_pd(iq0,jq0);
243 /* Calculate table index by multiplying r with table scale and truncate to integer */
244 rt = _mm_mul_pd(r00,vftabscale);
245 vfitab = _mm_cvttpd_epi32(rt);
247 vfeps = _mm_frcz_pd(rt);
249 vfeps = _mm_sub_pd(rt,_mm_round_pd(rt, _MM_FROUND_FLOOR));
251 twovfeps = _mm_add_pd(vfeps,vfeps);
252 vfitab = _mm_slli_epi32(vfitab,2);
254 /* CUBIC SPLINE TABLE ELECTROSTATICS */
255 Y = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) );
256 F = _mm_setzero_pd();
257 GMX_MM_TRANSPOSE2_PD(Y,F);
258 G = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) +2);
259 H = _mm_setzero_pd();
260 GMX_MM_TRANSPOSE2_PD(G,H);
261 Fp = _mm_macc_pd(vfeps,_mm_macc_pd(vfeps,H,G),F);
262 VV = _mm_macc_pd(vfeps,Fp,Y);
263 velec = _mm_mul_pd(qq00,VV);
264 FF = _mm_macc_pd(_mm_macc_pd(twovfeps,H,G),vfeps,Fp);
265 felec = _mm_xor_pd(signbit,_mm_mul_pd(_mm_mul_pd(qq00,FF),_mm_mul_pd(vftabscale,rinv00)));
267 /* Update potential sum for this i atom from the interaction with this j atom. */
268 velec = _mm_unpacklo_pd(velec,_mm_setzero_pd());
269 velecsum = _mm_add_pd(velecsum,velec);
273 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
275 /* Update vectorial force */
276 fix0 = _mm_macc_pd(dx00,fscal,fix0);
277 fiy0 = _mm_macc_pd(dy00,fscal,fiy0);
278 fiz0 = _mm_macc_pd(dz00,fscal,fiz0);
280 gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,
281 _mm_mul_pd(dx00,fscal),
282 _mm_mul_pd(dy00,fscal),
283 _mm_mul_pd(dz00,fscal));
285 /* Inner loop uses 46 flops */
288 /* End of innermost loop */
290 gmx_mm_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
291 f+i_coord_offset,fshift+i_shift_offset);
294 /* Update potential energies */
295 gmx_mm_update_1pot_pd(velecsum,kernel_data->energygrp_elec+ggid);
297 /* Increment number of inner iterations */
298 inneriter += j_index_end - j_index_start;
300 /* Outer loop uses 8 flops */
303 /* Increment number of outer iterations */
306 /* Update outer/inner flops */
308 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VF,outeriter*8 + inneriter*46);
311 * Gromacs nonbonded kernel: nb_kernel_ElecCSTab_VdwNone_GeomP1P1_F_avx_128_fma_double
312 * Electrostatics interaction: CubicSplineTable
313 * VdW interaction: None
314 * Geometry: Particle-Particle
315 * Calculate force/pot: Force
318 nb_kernel_ElecCSTab_VdwNone_GeomP1P1_F_avx_128_fma_double
319 (t_nblist * gmx_restrict nlist,
320 rvec * gmx_restrict xx,
321 rvec * gmx_restrict ff,
322 t_forcerec * gmx_restrict fr,
323 t_mdatoms * gmx_restrict mdatoms,
324 nb_kernel_data_t * gmx_restrict kernel_data,
325 t_nrnb * gmx_restrict nrnb)
327 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
328 * just 0 for non-waters.
329 * Suffixes A,B refer to j loop unrolling done with SSE double precision, e.g. for the two different
330 * jnr indices corresponding to data put in the four positions in the SIMD register.
332 int i_shift_offset,i_coord_offset,outeriter,inneriter;
333 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
335 int j_coord_offsetA,j_coord_offsetB;
336 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
338 real *shiftvec,*fshift,*x,*f;
339 __m128d tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
341 __m128d ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
342 int vdwjidx0A,vdwjidx0B;
343 __m128d jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
344 __m128d dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
345 __m128d velec,felec,velecsum,facel,crf,krf,krf2;
348 __m128i ifour = _mm_set1_epi32(4);
349 __m128d rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF,twovfeps;
351 __m128d dummy_mask,cutoff_mask;
352 __m128d signbit = gmx_mm_castsi128_pd( _mm_set_epi32(0x80000000,0x00000000,0x80000000,0x00000000) );
353 __m128d one = _mm_set1_pd(1.0);
354 __m128d two = _mm_set1_pd(2.0);
360 jindex = nlist->jindex;
362 shiftidx = nlist->shift;
364 shiftvec = fr->shift_vec[0];
365 fshift = fr->fshift[0];
366 facel = _mm_set1_pd(fr->epsfac);
367 charge = mdatoms->chargeA;
369 vftab = kernel_data->table_elec->data;
370 vftabscale = _mm_set1_pd(kernel_data->table_elec->scale);
372 /* Avoid stupid compiler warnings */
380 /* Start outer loop over neighborlists */
381 for(iidx=0; iidx<nri; iidx++)
383 /* Load shift vector for this list */
384 i_shift_offset = DIM*shiftidx[iidx];
386 /* Load limits for loop over neighbors */
387 j_index_start = jindex[iidx];
388 j_index_end = jindex[iidx+1];
390 /* Get outer coordinate index */
392 i_coord_offset = DIM*inr;
394 /* Load i particle coords and add shift vector */
395 gmx_mm_load_shift_and_1rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
397 fix0 = _mm_setzero_pd();
398 fiy0 = _mm_setzero_pd();
399 fiz0 = _mm_setzero_pd();
401 /* Load parameters for i particles */
402 iq0 = _mm_mul_pd(facel,_mm_load1_pd(charge+inr+0));
404 /* Start inner kernel loop */
405 for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
408 /* Get j neighbor index, and coordinate index */
411 j_coord_offsetA = DIM*jnrA;
412 j_coord_offsetB = DIM*jnrB;
414 /* load j atom coordinates */
415 gmx_mm_load_1rvec_2ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
418 /* Calculate displacement vector */
419 dx00 = _mm_sub_pd(ix0,jx0);
420 dy00 = _mm_sub_pd(iy0,jy0);
421 dz00 = _mm_sub_pd(iz0,jz0);
423 /* Calculate squared distance and things based on it */
424 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
426 rinv00 = gmx_mm_invsqrt_pd(rsq00);
428 /* Load parameters for j particles */
429 jq0 = gmx_mm_load_2real_swizzle_pd(charge+jnrA+0,charge+jnrB+0);
431 /**************************
432 * CALCULATE INTERACTIONS *
433 **************************/
435 r00 = _mm_mul_pd(rsq00,rinv00);
437 /* Compute parameters for interactions between i and j atoms */
438 qq00 = _mm_mul_pd(iq0,jq0);
440 /* Calculate table index by multiplying r with table scale and truncate to integer */
441 rt = _mm_mul_pd(r00,vftabscale);
442 vfitab = _mm_cvttpd_epi32(rt);
444 vfeps = _mm_frcz_pd(rt);
446 vfeps = _mm_sub_pd(rt,_mm_round_pd(rt, _MM_FROUND_FLOOR));
448 twovfeps = _mm_add_pd(vfeps,vfeps);
449 vfitab = _mm_slli_epi32(vfitab,2);
451 /* CUBIC SPLINE TABLE ELECTROSTATICS */
452 Y = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) );
453 F = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,1) );
454 GMX_MM_TRANSPOSE2_PD(Y,F);
455 G = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) +2);
456 H = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,1) +2);
457 GMX_MM_TRANSPOSE2_PD(G,H);
458 Fp = _mm_macc_pd(vfeps,_mm_macc_pd(vfeps,H,G),F);
459 FF = _mm_macc_pd(_mm_macc_pd(twovfeps,H,G),vfeps,Fp);
460 felec = _mm_xor_pd(signbit,_mm_mul_pd(_mm_mul_pd(qq00,FF),_mm_mul_pd(vftabscale,rinv00)));
464 /* Update vectorial force */
465 fix0 = _mm_macc_pd(dx00,fscal,fix0);
466 fiy0 = _mm_macc_pd(dy00,fscal,fiy0);
467 fiz0 = _mm_macc_pd(dz00,fscal,fiz0);
469 gmx_mm_decrement_1rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,
470 _mm_mul_pd(dx00,fscal),
471 _mm_mul_pd(dy00,fscal),
472 _mm_mul_pd(dz00,fscal));
474 /* Inner loop uses 42 flops */
481 j_coord_offsetA = DIM*jnrA;
483 /* load j atom coordinates */
484 gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
487 /* Calculate displacement vector */
488 dx00 = _mm_sub_pd(ix0,jx0);
489 dy00 = _mm_sub_pd(iy0,jy0);
490 dz00 = _mm_sub_pd(iz0,jz0);
492 /* Calculate squared distance and things based on it */
493 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
495 rinv00 = gmx_mm_invsqrt_pd(rsq00);
497 /* Load parameters for j particles */
498 jq0 = _mm_load_sd(charge+jnrA+0);
500 /**************************
501 * CALCULATE INTERACTIONS *
502 **************************/
504 r00 = _mm_mul_pd(rsq00,rinv00);
506 /* Compute parameters for interactions between i and j atoms */
507 qq00 = _mm_mul_pd(iq0,jq0);
509 /* Calculate table index by multiplying r with table scale and truncate to integer */
510 rt = _mm_mul_pd(r00,vftabscale);
511 vfitab = _mm_cvttpd_epi32(rt);
513 vfeps = _mm_frcz_pd(rt);
515 vfeps = _mm_sub_pd(rt,_mm_round_pd(rt, _MM_FROUND_FLOOR));
517 twovfeps = _mm_add_pd(vfeps,vfeps);
518 vfitab = _mm_slli_epi32(vfitab,2);
520 /* CUBIC SPLINE TABLE ELECTROSTATICS */
521 Y = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) );
522 F = _mm_setzero_pd();
523 GMX_MM_TRANSPOSE2_PD(Y,F);
524 G = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) +2);
525 H = _mm_setzero_pd();
526 GMX_MM_TRANSPOSE2_PD(G,H);
527 Fp = _mm_macc_pd(vfeps,_mm_macc_pd(vfeps,H,G),F);
528 FF = _mm_macc_pd(_mm_macc_pd(twovfeps,H,G),vfeps,Fp);
529 felec = _mm_xor_pd(signbit,_mm_mul_pd(_mm_mul_pd(qq00,FF),_mm_mul_pd(vftabscale,rinv00)));
533 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
535 /* Update vectorial force */
536 fix0 = _mm_macc_pd(dx00,fscal,fix0);
537 fiy0 = _mm_macc_pd(dy00,fscal,fiy0);
538 fiz0 = _mm_macc_pd(dz00,fscal,fiz0);
540 gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,
541 _mm_mul_pd(dx00,fscal),
542 _mm_mul_pd(dy00,fscal),
543 _mm_mul_pd(dz00,fscal));
545 /* Inner loop uses 42 flops */
548 /* End of innermost loop */
550 gmx_mm_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
551 f+i_coord_offset,fshift+i_shift_offset);
553 /* Increment number of inner iterations */
554 inneriter += j_index_end - j_index_start;
556 /* Outer loop uses 7 flops */
559 /* Increment number of outer iterations */
562 /* Update outer/inner flops */
564 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_F,outeriter*7 + inneriter*42);