2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
36 * Note: this file was generated by the GROMACS avx_128_fma_double kernel generator.
44 #include "../nb_kernel.h"
45 #include "types/simple.h"
49 #include "gmx_math_x86_avx_128_fma_double.h"
50 #include "kernelutil_x86_avx_128_fma_double.h"
53 * Gromacs nonbonded kernel: nb_kernel_ElecCSTab_VdwNone_GeomP1P1_VF_avx_128_fma_double
54 * Electrostatics interaction: CubicSplineTable
55 * VdW interaction: None
56 * Geometry: Particle-Particle
57 * Calculate force/pot: PotentialAndForce
60 nb_kernel_ElecCSTab_VdwNone_GeomP1P1_VF_avx_128_fma_double
61 (t_nblist * gmx_restrict nlist,
62 rvec * gmx_restrict xx,
63 rvec * gmx_restrict ff,
64 t_forcerec * gmx_restrict fr,
65 t_mdatoms * gmx_restrict mdatoms,
66 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
67 t_nrnb * gmx_restrict nrnb)
69 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
70 * just 0 for non-waters.
71 * Suffixes A,B refer to j loop unrolling done with SSE double precision, e.g. for the two different
72 * jnr indices corresponding to data put in the four positions in the SIMD register.
74 int i_shift_offset,i_coord_offset,outeriter,inneriter;
75 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
77 int j_coord_offsetA,j_coord_offsetB;
78 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
80 real *shiftvec,*fshift,*x,*f;
81 __m128d tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
83 __m128d ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
84 int vdwjidx0A,vdwjidx0B;
85 __m128d jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
86 __m128d dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
87 __m128d velec,felec,velecsum,facel,crf,krf,krf2;
90 __m128i ifour = _mm_set1_epi32(4);
91 __m128d rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF,twovfeps;
93 __m128d dummy_mask,cutoff_mask;
94 __m128d signbit = gmx_mm_castsi128_pd( _mm_set_epi32(0x80000000,0x00000000,0x80000000,0x00000000) );
95 __m128d one = _mm_set1_pd(1.0);
96 __m128d two = _mm_set1_pd(2.0);
102 jindex = nlist->jindex;
104 shiftidx = nlist->shift;
106 shiftvec = fr->shift_vec[0];
107 fshift = fr->fshift[0];
108 facel = _mm_set1_pd(fr->epsfac);
109 charge = mdatoms->chargeA;
111 vftab = kernel_data->table_elec->data;
112 vftabscale = _mm_set1_pd(kernel_data->table_elec->scale);
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));
146 /* Reset potential sums */
147 velecsum = _mm_setzero_pd();
149 /* Start inner kernel loop */
150 for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
153 /* Get j neighbor index, and coordinate index */
156 j_coord_offsetA = DIM*jnrA;
157 j_coord_offsetB = DIM*jnrB;
159 /* load j atom coordinates */
160 gmx_mm_load_1rvec_2ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
163 /* Calculate displacement vector */
164 dx00 = _mm_sub_pd(ix0,jx0);
165 dy00 = _mm_sub_pd(iy0,jy0);
166 dz00 = _mm_sub_pd(iz0,jz0);
168 /* Calculate squared distance and things based on it */
169 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
171 rinv00 = gmx_mm_invsqrt_pd(rsq00);
173 /* Load parameters for j particles */
174 jq0 = gmx_mm_load_2real_swizzle_pd(charge+jnrA+0,charge+jnrB+0);
176 /**************************
177 * CALCULATE INTERACTIONS *
178 **************************/
180 r00 = _mm_mul_pd(rsq00,rinv00);
182 /* Compute parameters for interactions between i and j atoms */
183 qq00 = _mm_mul_pd(iq0,jq0);
185 /* Calculate table index by multiplying r with table scale and truncate to integer */
186 rt = _mm_mul_pd(r00,vftabscale);
187 vfitab = _mm_cvttpd_epi32(rt);
189 vfeps = _mm_frcz_pd(rt);
191 vfeps = _mm_sub_pd(rt,_mm_round_pd(rt, _MM_FROUND_FLOOR));
193 twovfeps = _mm_add_pd(vfeps,vfeps);
194 vfitab = _mm_slli_epi32(vfitab,2);
196 /* CUBIC SPLINE TABLE ELECTROSTATICS */
197 Y = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) );
198 F = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,1) );
199 GMX_MM_TRANSPOSE2_PD(Y,F);
200 G = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) +2);
201 H = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,1) +2);
202 GMX_MM_TRANSPOSE2_PD(G,H);
203 Fp = _mm_macc_pd(vfeps,_mm_macc_pd(vfeps,H,G),F);
204 VV = _mm_macc_pd(vfeps,Fp,Y);
205 velec = _mm_mul_pd(qq00,VV);
206 FF = _mm_macc_pd(_mm_macc_pd(twovfeps,H,G),vfeps,Fp);
207 felec = _mm_xor_pd(signbit,_mm_mul_pd(_mm_mul_pd(qq00,FF),_mm_mul_pd(vftabscale,rinv00)));
209 /* Update potential sum for this i atom from the interaction with this j atom. */
210 velecsum = _mm_add_pd(velecsum,velec);
214 /* Update vectorial force */
215 fix0 = _mm_macc_pd(dx00,fscal,fix0);
216 fiy0 = _mm_macc_pd(dy00,fscal,fiy0);
217 fiz0 = _mm_macc_pd(dz00,fscal,fiz0);
219 gmx_mm_decrement_1rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,
220 _mm_mul_pd(dx00,fscal),
221 _mm_mul_pd(dy00,fscal),
222 _mm_mul_pd(dz00,fscal));
224 /* Inner loop uses 46 flops */
231 j_coord_offsetA = DIM*jnrA;
233 /* load j atom coordinates */
234 gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
237 /* Calculate displacement vector */
238 dx00 = _mm_sub_pd(ix0,jx0);
239 dy00 = _mm_sub_pd(iy0,jy0);
240 dz00 = _mm_sub_pd(iz0,jz0);
242 /* Calculate squared distance and things based on it */
243 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
245 rinv00 = gmx_mm_invsqrt_pd(rsq00);
247 /* Load parameters for j particles */
248 jq0 = _mm_load_sd(charge+jnrA+0);
250 /**************************
251 * CALCULATE INTERACTIONS *
252 **************************/
254 r00 = _mm_mul_pd(rsq00,rinv00);
256 /* Compute parameters for interactions between i and j atoms */
257 qq00 = _mm_mul_pd(iq0,jq0);
259 /* Calculate table index by multiplying r with table scale and truncate to integer */
260 rt = _mm_mul_pd(r00,vftabscale);
261 vfitab = _mm_cvttpd_epi32(rt);
263 vfeps = _mm_frcz_pd(rt);
265 vfeps = _mm_sub_pd(rt,_mm_round_pd(rt, _MM_FROUND_FLOOR));
267 twovfeps = _mm_add_pd(vfeps,vfeps);
268 vfitab = _mm_slli_epi32(vfitab,2);
270 /* CUBIC SPLINE TABLE ELECTROSTATICS */
271 Y = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) );
272 F = _mm_setzero_pd();
273 GMX_MM_TRANSPOSE2_PD(Y,F);
274 G = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) +2);
275 H = _mm_setzero_pd();
276 GMX_MM_TRANSPOSE2_PD(G,H);
277 Fp = _mm_macc_pd(vfeps,_mm_macc_pd(vfeps,H,G),F);
278 VV = _mm_macc_pd(vfeps,Fp,Y);
279 velec = _mm_mul_pd(qq00,VV);
280 FF = _mm_macc_pd(_mm_macc_pd(twovfeps,H,G),vfeps,Fp);
281 felec = _mm_xor_pd(signbit,_mm_mul_pd(_mm_mul_pd(qq00,FF),_mm_mul_pd(vftabscale,rinv00)));
283 /* Update potential sum for this i atom from the interaction with this j atom. */
284 velec = _mm_unpacklo_pd(velec,_mm_setzero_pd());
285 velecsum = _mm_add_pd(velecsum,velec);
289 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
291 /* Update vectorial force */
292 fix0 = _mm_macc_pd(dx00,fscal,fix0);
293 fiy0 = _mm_macc_pd(dy00,fscal,fiy0);
294 fiz0 = _mm_macc_pd(dz00,fscal,fiz0);
296 gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,
297 _mm_mul_pd(dx00,fscal),
298 _mm_mul_pd(dy00,fscal),
299 _mm_mul_pd(dz00,fscal));
301 /* Inner loop uses 46 flops */
304 /* End of innermost loop */
306 gmx_mm_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
307 f+i_coord_offset,fshift+i_shift_offset);
310 /* Update potential energies */
311 gmx_mm_update_1pot_pd(velecsum,kernel_data->energygrp_elec+ggid);
313 /* Increment number of inner iterations */
314 inneriter += j_index_end - j_index_start;
316 /* Outer loop uses 8 flops */
319 /* Increment number of outer iterations */
322 /* Update outer/inner flops */
324 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VF,outeriter*8 + inneriter*46);
327 * Gromacs nonbonded kernel: nb_kernel_ElecCSTab_VdwNone_GeomP1P1_F_avx_128_fma_double
328 * Electrostatics interaction: CubicSplineTable
329 * VdW interaction: None
330 * Geometry: Particle-Particle
331 * Calculate force/pot: Force
334 nb_kernel_ElecCSTab_VdwNone_GeomP1P1_F_avx_128_fma_double
335 (t_nblist * gmx_restrict nlist,
336 rvec * gmx_restrict xx,
337 rvec * gmx_restrict ff,
338 t_forcerec * gmx_restrict fr,
339 t_mdatoms * gmx_restrict mdatoms,
340 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
341 t_nrnb * gmx_restrict nrnb)
343 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
344 * just 0 for non-waters.
345 * Suffixes A,B refer to j loop unrolling done with SSE double precision, e.g. for the two different
346 * jnr indices corresponding to data put in the four positions in the SIMD register.
348 int i_shift_offset,i_coord_offset,outeriter,inneriter;
349 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
351 int j_coord_offsetA,j_coord_offsetB;
352 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
354 real *shiftvec,*fshift,*x,*f;
355 __m128d tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
357 __m128d ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
358 int vdwjidx0A,vdwjidx0B;
359 __m128d jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
360 __m128d dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
361 __m128d velec,felec,velecsum,facel,crf,krf,krf2;
364 __m128i ifour = _mm_set1_epi32(4);
365 __m128d rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF,twovfeps;
367 __m128d dummy_mask,cutoff_mask;
368 __m128d signbit = gmx_mm_castsi128_pd( _mm_set_epi32(0x80000000,0x00000000,0x80000000,0x00000000) );
369 __m128d one = _mm_set1_pd(1.0);
370 __m128d two = _mm_set1_pd(2.0);
376 jindex = nlist->jindex;
378 shiftidx = nlist->shift;
380 shiftvec = fr->shift_vec[0];
381 fshift = fr->fshift[0];
382 facel = _mm_set1_pd(fr->epsfac);
383 charge = mdatoms->chargeA;
385 vftab = kernel_data->table_elec->data;
386 vftabscale = _mm_set1_pd(kernel_data->table_elec->scale);
388 /* Avoid stupid compiler warnings */
396 /* Start outer loop over neighborlists */
397 for(iidx=0; iidx<nri; iidx++)
399 /* Load shift vector for this list */
400 i_shift_offset = DIM*shiftidx[iidx];
402 /* Load limits for loop over neighbors */
403 j_index_start = jindex[iidx];
404 j_index_end = jindex[iidx+1];
406 /* Get outer coordinate index */
408 i_coord_offset = DIM*inr;
410 /* Load i particle coords and add shift vector */
411 gmx_mm_load_shift_and_1rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
413 fix0 = _mm_setzero_pd();
414 fiy0 = _mm_setzero_pd();
415 fiz0 = _mm_setzero_pd();
417 /* Load parameters for i particles */
418 iq0 = _mm_mul_pd(facel,_mm_load1_pd(charge+inr+0));
420 /* Start inner kernel loop */
421 for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
424 /* Get j neighbor index, and coordinate index */
427 j_coord_offsetA = DIM*jnrA;
428 j_coord_offsetB = DIM*jnrB;
430 /* load j atom coordinates */
431 gmx_mm_load_1rvec_2ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
434 /* Calculate displacement vector */
435 dx00 = _mm_sub_pd(ix0,jx0);
436 dy00 = _mm_sub_pd(iy0,jy0);
437 dz00 = _mm_sub_pd(iz0,jz0);
439 /* Calculate squared distance and things based on it */
440 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
442 rinv00 = gmx_mm_invsqrt_pd(rsq00);
444 /* Load parameters for j particles */
445 jq0 = gmx_mm_load_2real_swizzle_pd(charge+jnrA+0,charge+jnrB+0);
447 /**************************
448 * CALCULATE INTERACTIONS *
449 **************************/
451 r00 = _mm_mul_pd(rsq00,rinv00);
453 /* Compute parameters for interactions between i and j atoms */
454 qq00 = _mm_mul_pd(iq0,jq0);
456 /* Calculate table index by multiplying r with table scale and truncate to integer */
457 rt = _mm_mul_pd(r00,vftabscale);
458 vfitab = _mm_cvttpd_epi32(rt);
460 vfeps = _mm_frcz_pd(rt);
462 vfeps = _mm_sub_pd(rt,_mm_round_pd(rt, _MM_FROUND_FLOOR));
464 twovfeps = _mm_add_pd(vfeps,vfeps);
465 vfitab = _mm_slli_epi32(vfitab,2);
467 /* CUBIC SPLINE TABLE ELECTROSTATICS */
468 Y = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) );
469 F = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,1) );
470 GMX_MM_TRANSPOSE2_PD(Y,F);
471 G = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) +2);
472 H = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,1) +2);
473 GMX_MM_TRANSPOSE2_PD(G,H);
474 Fp = _mm_macc_pd(vfeps,_mm_macc_pd(vfeps,H,G),F);
475 FF = _mm_macc_pd(_mm_macc_pd(twovfeps,H,G),vfeps,Fp);
476 felec = _mm_xor_pd(signbit,_mm_mul_pd(_mm_mul_pd(qq00,FF),_mm_mul_pd(vftabscale,rinv00)));
480 /* Update vectorial force */
481 fix0 = _mm_macc_pd(dx00,fscal,fix0);
482 fiy0 = _mm_macc_pd(dy00,fscal,fiy0);
483 fiz0 = _mm_macc_pd(dz00,fscal,fiz0);
485 gmx_mm_decrement_1rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,
486 _mm_mul_pd(dx00,fscal),
487 _mm_mul_pd(dy00,fscal),
488 _mm_mul_pd(dz00,fscal));
490 /* Inner loop uses 42 flops */
497 j_coord_offsetA = DIM*jnrA;
499 /* load j atom coordinates */
500 gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
503 /* Calculate displacement vector */
504 dx00 = _mm_sub_pd(ix0,jx0);
505 dy00 = _mm_sub_pd(iy0,jy0);
506 dz00 = _mm_sub_pd(iz0,jz0);
508 /* Calculate squared distance and things based on it */
509 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
511 rinv00 = gmx_mm_invsqrt_pd(rsq00);
513 /* Load parameters for j particles */
514 jq0 = _mm_load_sd(charge+jnrA+0);
516 /**************************
517 * CALCULATE INTERACTIONS *
518 **************************/
520 r00 = _mm_mul_pd(rsq00,rinv00);
522 /* Compute parameters for interactions between i and j atoms */
523 qq00 = _mm_mul_pd(iq0,jq0);
525 /* Calculate table index by multiplying r with table scale and truncate to integer */
526 rt = _mm_mul_pd(r00,vftabscale);
527 vfitab = _mm_cvttpd_epi32(rt);
529 vfeps = _mm_frcz_pd(rt);
531 vfeps = _mm_sub_pd(rt,_mm_round_pd(rt, _MM_FROUND_FLOOR));
533 twovfeps = _mm_add_pd(vfeps,vfeps);
534 vfitab = _mm_slli_epi32(vfitab,2);
536 /* CUBIC SPLINE TABLE ELECTROSTATICS */
537 Y = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) );
538 F = _mm_setzero_pd();
539 GMX_MM_TRANSPOSE2_PD(Y,F);
540 G = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) +2);
541 H = _mm_setzero_pd();
542 GMX_MM_TRANSPOSE2_PD(G,H);
543 Fp = _mm_macc_pd(vfeps,_mm_macc_pd(vfeps,H,G),F);
544 FF = _mm_macc_pd(_mm_macc_pd(twovfeps,H,G),vfeps,Fp);
545 felec = _mm_xor_pd(signbit,_mm_mul_pd(_mm_mul_pd(qq00,FF),_mm_mul_pd(vftabscale,rinv00)));
549 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
551 /* Update vectorial force */
552 fix0 = _mm_macc_pd(dx00,fscal,fix0);
553 fiy0 = _mm_macc_pd(dy00,fscal,fiy0);
554 fiz0 = _mm_macc_pd(dz00,fscal,fiz0);
556 gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,
557 _mm_mul_pd(dx00,fscal),
558 _mm_mul_pd(dy00,fscal),
559 _mm_mul_pd(dz00,fscal));
561 /* Inner loop uses 42 flops */
564 /* End of innermost loop */
566 gmx_mm_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
567 f+i_coord_offset,fshift+i_shift_offset);
569 /* Increment number of inner iterations */
570 inneriter += j_index_end - j_index_start;
572 /* Outer loop uses 7 flops */
575 /* Increment number of outer iterations */
578 /* Update outer/inner flops */
580 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_F,outeriter*7 + inneriter*42);