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_ElecRFCut_VdwLJSh_GeomW4P1_VF_avx_128_fma_double
38 * Electrostatics interaction: ReactionField
39 * VdW interaction: LennardJones
40 * Geometry: Water4-Particle
41 * Calculate force/pot: PotentialAndForce
44 nb_kernel_ElecRFCut_VdwLJSh_GeomW4P1_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;
69 __m128d ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
71 __m128d ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
73 __m128d ix3,iy3,iz3,fix3,fiy3,fiz3,iq3,isai3;
74 int vdwjidx0A,vdwjidx0B;
75 __m128d jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
76 __m128d dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
77 __m128d dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
78 __m128d dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
79 __m128d dx30,dy30,dz30,rsq30,rinv30,rinvsq30,r30,qq30,c6_30,c12_30;
80 __m128d velec,felec,velecsum,facel,crf,krf,krf2;
83 __m128d rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
86 __m128d one_sixth = _mm_set1_pd(1.0/6.0);
87 __m128d one_twelfth = _mm_set1_pd(1.0/12.0);
88 __m128d dummy_mask,cutoff_mask;
89 __m128d signbit = gmx_mm_castsi128_pd( _mm_set_epi32(0x80000000,0x00000000,0x80000000,0x00000000) );
90 __m128d one = _mm_set1_pd(1.0);
91 __m128d two = _mm_set1_pd(2.0);
97 jindex = nlist->jindex;
99 shiftidx = nlist->shift;
101 shiftvec = fr->shift_vec[0];
102 fshift = fr->fshift[0];
103 facel = _mm_set1_pd(fr->epsfac);
104 charge = mdatoms->chargeA;
105 krf = _mm_set1_pd(fr->ic->k_rf);
106 krf2 = _mm_set1_pd(fr->ic->k_rf*2.0);
107 crf = _mm_set1_pd(fr->ic->c_rf);
108 nvdwtype = fr->ntype;
110 vdwtype = mdatoms->typeA;
112 /* Setup water-specific parameters */
113 inr = nlist->iinr[0];
114 iq1 = _mm_mul_pd(facel,_mm_set1_pd(charge[inr+1]));
115 iq2 = _mm_mul_pd(facel,_mm_set1_pd(charge[inr+2]));
116 iq3 = _mm_mul_pd(facel,_mm_set1_pd(charge[inr+3]));
117 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
119 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
120 rcutoff_scalar = fr->rcoulomb;
121 rcutoff = _mm_set1_pd(rcutoff_scalar);
122 rcutoff2 = _mm_mul_pd(rcutoff,rcutoff);
124 sh_vdw_invrcut6 = _mm_set1_pd(fr->ic->sh_invrc6);
125 rvdw = _mm_set1_pd(fr->rvdw);
127 /* Avoid stupid compiler warnings */
135 /* Start outer loop over neighborlists */
136 for(iidx=0; iidx<nri; iidx++)
138 /* Load shift vector for this list */
139 i_shift_offset = DIM*shiftidx[iidx];
141 /* Load limits for loop over neighbors */
142 j_index_start = jindex[iidx];
143 j_index_end = jindex[iidx+1];
145 /* Get outer coordinate index */
147 i_coord_offset = DIM*inr;
149 /* Load i particle coords and add shift vector */
150 gmx_mm_load_shift_and_4rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,
151 &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
153 fix0 = _mm_setzero_pd();
154 fiy0 = _mm_setzero_pd();
155 fiz0 = _mm_setzero_pd();
156 fix1 = _mm_setzero_pd();
157 fiy1 = _mm_setzero_pd();
158 fiz1 = _mm_setzero_pd();
159 fix2 = _mm_setzero_pd();
160 fiy2 = _mm_setzero_pd();
161 fiz2 = _mm_setzero_pd();
162 fix3 = _mm_setzero_pd();
163 fiy3 = _mm_setzero_pd();
164 fiz3 = _mm_setzero_pd();
166 /* Reset potential sums */
167 velecsum = _mm_setzero_pd();
168 vvdwsum = _mm_setzero_pd();
170 /* Start inner kernel loop */
171 for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
174 /* Get j neighbor index, and coordinate index */
177 j_coord_offsetA = DIM*jnrA;
178 j_coord_offsetB = DIM*jnrB;
180 /* load j atom coordinates */
181 gmx_mm_load_1rvec_2ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
184 /* Calculate displacement vector */
185 dx00 = _mm_sub_pd(ix0,jx0);
186 dy00 = _mm_sub_pd(iy0,jy0);
187 dz00 = _mm_sub_pd(iz0,jz0);
188 dx10 = _mm_sub_pd(ix1,jx0);
189 dy10 = _mm_sub_pd(iy1,jy0);
190 dz10 = _mm_sub_pd(iz1,jz0);
191 dx20 = _mm_sub_pd(ix2,jx0);
192 dy20 = _mm_sub_pd(iy2,jy0);
193 dz20 = _mm_sub_pd(iz2,jz0);
194 dx30 = _mm_sub_pd(ix3,jx0);
195 dy30 = _mm_sub_pd(iy3,jy0);
196 dz30 = _mm_sub_pd(iz3,jz0);
198 /* Calculate squared distance and things based on it */
199 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
200 rsq10 = gmx_mm_calc_rsq_pd(dx10,dy10,dz10);
201 rsq20 = gmx_mm_calc_rsq_pd(dx20,dy20,dz20);
202 rsq30 = gmx_mm_calc_rsq_pd(dx30,dy30,dz30);
204 rinv10 = gmx_mm_invsqrt_pd(rsq10);
205 rinv20 = gmx_mm_invsqrt_pd(rsq20);
206 rinv30 = gmx_mm_invsqrt_pd(rsq30);
208 rinvsq00 = gmx_mm_inv_pd(rsq00);
209 rinvsq10 = _mm_mul_pd(rinv10,rinv10);
210 rinvsq20 = _mm_mul_pd(rinv20,rinv20);
211 rinvsq30 = _mm_mul_pd(rinv30,rinv30);
213 /* Load parameters for j particles */
214 jq0 = gmx_mm_load_2real_swizzle_pd(charge+jnrA+0,charge+jnrB+0);
215 vdwjidx0A = 2*vdwtype[jnrA+0];
216 vdwjidx0B = 2*vdwtype[jnrB+0];
218 fjx0 = _mm_setzero_pd();
219 fjy0 = _mm_setzero_pd();
220 fjz0 = _mm_setzero_pd();
222 /**************************
223 * CALCULATE INTERACTIONS *
224 **************************/
226 if (gmx_mm_any_lt(rsq00,rcutoff2))
229 /* Compute parameters for interactions between i and j atoms */
230 gmx_mm_load_2pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,
231 vdwparam+vdwioffset0+vdwjidx0B,&c6_00,&c12_00);
233 /* LENNARD-JONES DISPERSION/REPULSION */
235 rinvsix = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
236 vvdw6 = _mm_mul_pd(c6_00,rinvsix);
237 vvdw12 = _mm_mul_pd(c12_00,_mm_mul_pd(rinvsix,rinvsix));
238 vvdw = _mm_msub_pd(_mm_nmacc_pd(c12_00,_mm_mul_pd(sh_vdw_invrcut6,sh_vdw_invrcut6),vvdw12),one_twelfth,
239 _mm_mul_pd(_mm_nmacc_pd( c6_00,sh_vdw_invrcut6,vvdw6),one_sixth));
240 fvdw = _mm_mul_pd(_mm_sub_pd(vvdw12,vvdw6),rinvsq00);
242 cutoff_mask = _mm_cmplt_pd(rsq00,rcutoff2);
244 /* Update potential sum for this i atom from the interaction with this j atom. */
245 vvdw = _mm_and_pd(vvdw,cutoff_mask);
246 vvdwsum = _mm_add_pd(vvdwsum,vvdw);
250 fscal = _mm_and_pd(fscal,cutoff_mask);
252 /* Update vectorial force */
253 fix0 = _mm_macc_pd(dx00,fscal,fix0);
254 fiy0 = _mm_macc_pd(dy00,fscal,fiy0);
255 fiz0 = _mm_macc_pd(dz00,fscal,fiz0);
257 fjx0 = _mm_macc_pd(dx00,fscal,fjx0);
258 fjy0 = _mm_macc_pd(dy00,fscal,fjy0);
259 fjz0 = _mm_macc_pd(dz00,fscal,fjz0);
263 /**************************
264 * CALCULATE INTERACTIONS *
265 **************************/
267 if (gmx_mm_any_lt(rsq10,rcutoff2))
270 /* Compute parameters for interactions between i and j atoms */
271 qq10 = _mm_mul_pd(iq1,jq0);
273 /* REACTION-FIELD ELECTROSTATICS */
274 velec = _mm_mul_pd(qq10,_mm_sub_pd(_mm_macc_pd(krf,rsq10,rinv10),crf));
275 felec = _mm_mul_pd(qq10,_mm_msub_pd(rinv10,rinvsq10,krf2));
277 cutoff_mask = _mm_cmplt_pd(rsq10,rcutoff2);
279 /* Update potential sum for this i atom from the interaction with this j atom. */
280 velec = _mm_and_pd(velec,cutoff_mask);
281 velecsum = _mm_add_pd(velecsum,velec);
285 fscal = _mm_and_pd(fscal,cutoff_mask);
287 /* Update vectorial force */
288 fix1 = _mm_macc_pd(dx10,fscal,fix1);
289 fiy1 = _mm_macc_pd(dy10,fscal,fiy1);
290 fiz1 = _mm_macc_pd(dz10,fscal,fiz1);
292 fjx0 = _mm_macc_pd(dx10,fscal,fjx0);
293 fjy0 = _mm_macc_pd(dy10,fscal,fjy0);
294 fjz0 = _mm_macc_pd(dz10,fscal,fjz0);
298 /**************************
299 * CALCULATE INTERACTIONS *
300 **************************/
302 if (gmx_mm_any_lt(rsq20,rcutoff2))
305 /* Compute parameters for interactions between i and j atoms */
306 qq20 = _mm_mul_pd(iq2,jq0);
308 /* REACTION-FIELD ELECTROSTATICS */
309 velec = _mm_mul_pd(qq20,_mm_sub_pd(_mm_macc_pd(krf,rsq20,rinv20),crf));
310 felec = _mm_mul_pd(qq20,_mm_msub_pd(rinv20,rinvsq20,krf2));
312 cutoff_mask = _mm_cmplt_pd(rsq20,rcutoff2);
314 /* Update potential sum for this i atom from the interaction with this j atom. */
315 velec = _mm_and_pd(velec,cutoff_mask);
316 velecsum = _mm_add_pd(velecsum,velec);
320 fscal = _mm_and_pd(fscal,cutoff_mask);
322 /* Update vectorial force */
323 fix2 = _mm_macc_pd(dx20,fscal,fix2);
324 fiy2 = _mm_macc_pd(dy20,fscal,fiy2);
325 fiz2 = _mm_macc_pd(dz20,fscal,fiz2);
327 fjx0 = _mm_macc_pd(dx20,fscal,fjx0);
328 fjy0 = _mm_macc_pd(dy20,fscal,fjy0);
329 fjz0 = _mm_macc_pd(dz20,fscal,fjz0);
333 /**************************
334 * CALCULATE INTERACTIONS *
335 **************************/
337 if (gmx_mm_any_lt(rsq30,rcutoff2))
340 /* Compute parameters for interactions between i and j atoms */
341 qq30 = _mm_mul_pd(iq3,jq0);
343 /* REACTION-FIELD ELECTROSTATICS */
344 velec = _mm_mul_pd(qq30,_mm_sub_pd(_mm_macc_pd(krf,rsq30,rinv30),crf));
345 felec = _mm_mul_pd(qq30,_mm_msub_pd(rinv30,rinvsq30,krf2));
347 cutoff_mask = _mm_cmplt_pd(rsq30,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 velecsum = _mm_add_pd(velecsum,velec);
355 fscal = _mm_and_pd(fscal,cutoff_mask);
357 /* Update vectorial force */
358 fix3 = _mm_macc_pd(dx30,fscal,fix3);
359 fiy3 = _mm_macc_pd(dy30,fscal,fiy3);
360 fiz3 = _mm_macc_pd(dz30,fscal,fiz3);
362 fjx0 = _mm_macc_pd(dx30,fscal,fjx0);
363 fjy0 = _mm_macc_pd(dy30,fscal,fjy0);
364 fjz0 = _mm_macc_pd(dz30,fscal,fjz0);
368 gmx_mm_decrement_1rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,fjx0,fjy0,fjz0);
370 /* Inner loop uses 164 flops */
377 j_coord_offsetA = DIM*jnrA;
379 /* load j atom coordinates */
380 gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
383 /* Calculate displacement vector */
384 dx00 = _mm_sub_pd(ix0,jx0);
385 dy00 = _mm_sub_pd(iy0,jy0);
386 dz00 = _mm_sub_pd(iz0,jz0);
387 dx10 = _mm_sub_pd(ix1,jx0);
388 dy10 = _mm_sub_pd(iy1,jy0);
389 dz10 = _mm_sub_pd(iz1,jz0);
390 dx20 = _mm_sub_pd(ix2,jx0);
391 dy20 = _mm_sub_pd(iy2,jy0);
392 dz20 = _mm_sub_pd(iz2,jz0);
393 dx30 = _mm_sub_pd(ix3,jx0);
394 dy30 = _mm_sub_pd(iy3,jy0);
395 dz30 = _mm_sub_pd(iz3,jz0);
397 /* Calculate squared distance and things based on it */
398 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
399 rsq10 = gmx_mm_calc_rsq_pd(dx10,dy10,dz10);
400 rsq20 = gmx_mm_calc_rsq_pd(dx20,dy20,dz20);
401 rsq30 = gmx_mm_calc_rsq_pd(dx30,dy30,dz30);
403 rinv10 = gmx_mm_invsqrt_pd(rsq10);
404 rinv20 = gmx_mm_invsqrt_pd(rsq20);
405 rinv30 = gmx_mm_invsqrt_pd(rsq30);
407 rinvsq00 = gmx_mm_inv_pd(rsq00);
408 rinvsq10 = _mm_mul_pd(rinv10,rinv10);
409 rinvsq20 = _mm_mul_pd(rinv20,rinv20);
410 rinvsq30 = _mm_mul_pd(rinv30,rinv30);
412 /* Load parameters for j particles */
413 jq0 = _mm_load_sd(charge+jnrA+0);
414 vdwjidx0A = 2*vdwtype[jnrA+0];
416 fjx0 = _mm_setzero_pd();
417 fjy0 = _mm_setzero_pd();
418 fjz0 = _mm_setzero_pd();
420 /**************************
421 * CALCULATE INTERACTIONS *
422 **************************/
424 if (gmx_mm_any_lt(rsq00,rcutoff2))
427 /* Compute parameters for interactions between i and j atoms */
428 gmx_mm_load_1pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,&c6_00,&c12_00);
430 /* LENNARD-JONES DISPERSION/REPULSION */
432 rinvsix = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
433 vvdw6 = _mm_mul_pd(c6_00,rinvsix);
434 vvdw12 = _mm_mul_pd(c12_00,_mm_mul_pd(rinvsix,rinvsix));
435 vvdw = _mm_msub_pd(_mm_nmacc_pd(c12_00,_mm_mul_pd(sh_vdw_invrcut6,sh_vdw_invrcut6),vvdw12),one_twelfth,
436 _mm_mul_pd(_mm_nmacc_pd( c6_00,sh_vdw_invrcut6,vvdw6),one_sixth));
437 fvdw = _mm_mul_pd(_mm_sub_pd(vvdw12,vvdw6),rinvsq00);
439 cutoff_mask = _mm_cmplt_pd(rsq00,rcutoff2);
441 /* Update potential sum for this i atom from the interaction with this j atom. */
442 vvdw = _mm_and_pd(vvdw,cutoff_mask);
443 vvdw = _mm_unpacklo_pd(vvdw,_mm_setzero_pd());
444 vvdwsum = _mm_add_pd(vvdwsum,vvdw);
448 fscal = _mm_and_pd(fscal,cutoff_mask);
450 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
452 /* Update vectorial force */
453 fix0 = _mm_macc_pd(dx00,fscal,fix0);
454 fiy0 = _mm_macc_pd(dy00,fscal,fiy0);
455 fiz0 = _mm_macc_pd(dz00,fscal,fiz0);
457 fjx0 = _mm_macc_pd(dx00,fscal,fjx0);
458 fjy0 = _mm_macc_pd(dy00,fscal,fjy0);
459 fjz0 = _mm_macc_pd(dz00,fscal,fjz0);
463 /**************************
464 * CALCULATE INTERACTIONS *
465 **************************/
467 if (gmx_mm_any_lt(rsq10,rcutoff2))
470 /* Compute parameters for interactions between i and j atoms */
471 qq10 = _mm_mul_pd(iq1,jq0);
473 /* REACTION-FIELD ELECTROSTATICS */
474 velec = _mm_mul_pd(qq10,_mm_sub_pd(_mm_macc_pd(krf,rsq10,rinv10),crf));
475 felec = _mm_mul_pd(qq10,_mm_msub_pd(rinv10,rinvsq10,krf2));
477 cutoff_mask = _mm_cmplt_pd(rsq10,rcutoff2);
479 /* Update potential sum for this i atom from the interaction with this j atom. */
480 velec = _mm_and_pd(velec,cutoff_mask);
481 velec = _mm_unpacklo_pd(velec,_mm_setzero_pd());
482 velecsum = _mm_add_pd(velecsum,velec);
486 fscal = _mm_and_pd(fscal,cutoff_mask);
488 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
490 /* Update vectorial force */
491 fix1 = _mm_macc_pd(dx10,fscal,fix1);
492 fiy1 = _mm_macc_pd(dy10,fscal,fiy1);
493 fiz1 = _mm_macc_pd(dz10,fscal,fiz1);
495 fjx0 = _mm_macc_pd(dx10,fscal,fjx0);
496 fjy0 = _mm_macc_pd(dy10,fscal,fjy0);
497 fjz0 = _mm_macc_pd(dz10,fscal,fjz0);
501 /**************************
502 * CALCULATE INTERACTIONS *
503 **************************/
505 if (gmx_mm_any_lt(rsq20,rcutoff2))
508 /* Compute parameters for interactions between i and j atoms */
509 qq20 = _mm_mul_pd(iq2,jq0);
511 /* REACTION-FIELD ELECTROSTATICS */
512 velec = _mm_mul_pd(qq20,_mm_sub_pd(_mm_macc_pd(krf,rsq20,rinv20),crf));
513 felec = _mm_mul_pd(qq20,_mm_msub_pd(rinv20,rinvsq20,krf2));
515 cutoff_mask = _mm_cmplt_pd(rsq20,rcutoff2);
517 /* Update potential sum for this i atom from the interaction with this j atom. */
518 velec = _mm_and_pd(velec,cutoff_mask);
519 velec = _mm_unpacklo_pd(velec,_mm_setzero_pd());
520 velecsum = _mm_add_pd(velecsum,velec);
524 fscal = _mm_and_pd(fscal,cutoff_mask);
526 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
528 /* Update vectorial force */
529 fix2 = _mm_macc_pd(dx20,fscal,fix2);
530 fiy2 = _mm_macc_pd(dy20,fscal,fiy2);
531 fiz2 = _mm_macc_pd(dz20,fscal,fiz2);
533 fjx0 = _mm_macc_pd(dx20,fscal,fjx0);
534 fjy0 = _mm_macc_pd(dy20,fscal,fjy0);
535 fjz0 = _mm_macc_pd(dz20,fscal,fjz0);
539 /**************************
540 * CALCULATE INTERACTIONS *
541 **************************/
543 if (gmx_mm_any_lt(rsq30,rcutoff2))
546 /* Compute parameters for interactions between i and j atoms */
547 qq30 = _mm_mul_pd(iq3,jq0);
549 /* REACTION-FIELD ELECTROSTATICS */
550 velec = _mm_mul_pd(qq30,_mm_sub_pd(_mm_macc_pd(krf,rsq30,rinv30),crf));
551 felec = _mm_mul_pd(qq30,_mm_msub_pd(rinv30,rinvsq30,krf2));
553 cutoff_mask = _mm_cmplt_pd(rsq30,rcutoff2);
555 /* Update potential sum for this i atom from the interaction with this j atom. */
556 velec = _mm_and_pd(velec,cutoff_mask);
557 velec = _mm_unpacklo_pd(velec,_mm_setzero_pd());
558 velecsum = _mm_add_pd(velecsum,velec);
562 fscal = _mm_and_pd(fscal,cutoff_mask);
564 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
566 /* Update vectorial force */
567 fix3 = _mm_macc_pd(dx30,fscal,fix3);
568 fiy3 = _mm_macc_pd(dy30,fscal,fiy3);
569 fiz3 = _mm_macc_pd(dz30,fscal,fiz3);
571 fjx0 = _mm_macc_pd(dx30,fscal,fjx0);
572 fjy0 = _mm_macc_pd(dy30,fscal,fjy0);
573 fjz0 = _mm_macc_pd(dz30,fscal,fjz0);
577 gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,fjx0,fjy0,fjz0);
579 /* Inner loop uses 164 flops */
582 /* End of innermost loop */
584 gmx_mm_update_iforce_4atom_swizzle_pd(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
585 f+i_coord_offset,fshift+i_shift_offset);
588 /* Update potential energies */
589 gmx_mm_update_1pot_pd(velecsum,kernel_data->energygrp_elec+ggid);
590 gmx_mm_update_1pot_pd(vvdwsum,kernel_data->energygrp_vdw+ggid);
592 /* Increment number of inner iterations */
593 inneriter += j_index_end - j_index_start;
595 /* Outer loop uses 26 flops */
598 /* Increment number of outer iterations */
601 /* Update outer/inner flops */
603 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W4_VF,outeriter*26 + inneriter*164);
606 * Gromacs nonbonded kernel: nb_kernel_ElecRFCut_VdwLJSh_GeomW4P1_F_avx_128_fma_double
607 * Electrostatics interaction: ReactionField
608 * VdW interaction: LennardJones
609 * Geometry: Water4-Particle
610 * Calculate force/pot: Force
613 nb_kernel_ElecRFCut_VdwLJSh_GeomW4P1_F_avx_128_fma_double
614 (t_nblist * gmx_restrict nlist,
615 rvec * gmx_restrict xx,
616 rvec * gmx_restrict ff,
617 t_forcerec * gmx_restrict fr,
618 t_mdatoms * gmx_restrict mdatoms,
619 nb_kernel_data_t * gmx_restrict kernel_data,
620 t_nrnb * gmx_restrict nrnb)
622 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
623 * just 0 for non-waters.
624 * Suffixes A,B refer to j loop unrolling done with SSE double precision, e.g. for the two different
625 * jnr indices corresponding to data put in the four positions in the SIMD register.
627 int i_shift_offset,i_coord_offset,outeriter,inneriter;
628 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
630 int j_coord_offsetA,j_coord_offsetB;
631 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
633 real *shiftvec,*fshift,*x,*f;
634 __m128d tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
636 __m128d ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
638 __m128d ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
640 __m128d ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
642 __m128d ix3,iy3,iz3,fix3,fiy3,fiz3,iq3,isai3;
643 int vdwjidx0A,vdwjidx0B;
644 __m128d jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
645 __m128d dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
646 __m128d dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
647 __m128d dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
648 __m128d dx30,dy30,dz30,rsq30,rinv30,rinvsq30,r30,qq30,c6_30,c12_30;
649 __m128d velec,felec,velecsum,facel,crf,krf,krf2;
652 __m128d rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
655 __m128d one_sixth = _mm_set1_pd(1.0/6.0);
656 __m128d one_twelfth = _mm_set1_pd(1.0/12.0);
657 __m128d dummy_mask,cutoff_mask;
658 __m128d signbit = gmx_mm_castsi128_pd( _mm_set_epi32(0x80000000,0x00000000,0x80000000,0x00000000) );
659 __m128d one = _mm_set1_pd(1.0);
660 __m128d two = _mm_set1_pd(2.0);
666 jindex = nlist->jindex;
668 shiftidx = nlist->shift;
670 shiftvec = fr->shift_vec[0];
671 fshift = fr->fshift[0];
672 facel = _mm_set1_pd(fr->epsfac);
673 charge = mdatoms->chargeA;
674 krf = _mm_set1_pd(fr->ic->k_rf);
675 krf2 = _mm_set1_pd(fr->ic->k_rf*2.0);
676 crf = _mm_set1_pd(fr->ic->c_rf);
677 nvdwtype = fr->ntype;
679 vdwtype = mdatoms->typeA;
681 /* Setup water-specific parameters */
682 inr = nlist->iinr[0];
683 iq1 = _mm_mul_pd(facel,_mm_set1_pd(charge[inr+1]));
684 iq2 = _mm_mul_pd(facel,_mm_set1_pd(charge[inr+2]));
685 iq3 = _mm_mul_pd(facel,_mm_set1_pd(charge[inr+3]));
686 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
688 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
689 rcutoff_scalar = fr->rcoulomb;
690 rcutoff = _mm_set1_pd(rcutoff_scalar);
691 rcutoff2 = _mm_mul_pd(rcutoff,rcutoff);
693 sh_vdw_invrcut6 = _mm_set1_pd(fr->ic->sh_invrc6);
694 rvdw = _mm_set1_pd(fr->rvdw);
696 /* Avoid stupid compiler warnings */
704 /* Start outer loop over neighborlists */
705 for(iidx=0; iidx<nri; iidx++)
707 /* Load shift vector for this list */
708 i_shift_offset = DIM*shiftidx[iidx];
710 /* Load limits for loop over neighbors */
711 j_index_start = jindex[iidx];
712 j_index_end = jindex[iidx+1];
714 /* Get outer coordinate index */
716 i_coord_offset = DIM*inr;
718 /* Load i particle coords and add shift vector */
719 gmx_mm_load_shift_and_4rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,
720 &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
722 fix0 = _mm_setzero_pd();
723 fiy0 = _mm_setzero_pd();
724 fiz0 = _mm_setzero_pd();
725 fix1 = _mm_setzero_pd();
726 fiy1 = _mm_setzero_pd();
727 fiz1 = _mm_setzero_pd();
728 fix2 = _mm_setzero_pd();
729 fiy2 = _mm_setzero_pd();
730 fiz2 = _mm_setzero_pd();
731 fix3 = _mm_setzero_pd();
732 fiy3 = _mm_setzero_pd();
733 fiz3 = _mm_setzero_pd();
735 /* Start inner kernel loop */
736 for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
739 /* Get j neighbor index, and coordinate index */
742 j_coord_offsetA = DIM*jnrA;
743 j_coord_offsetB = DIM*jnrB;
745 /* load j atom coordinates */
746 gmx_mm_load_1rvec_2ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
749 /* Calculate displacement vector */
750 dx00 = _mm_sub_pd(ix0,jx0);
751 dy00 = _mm_sub_pd(iy0,jy0);
752 dz00 = _mm_sub_pd(iz0,jz0);
753 dx10 = _mm_sub_pd(ix1,jx0);
754 dy10 = _mm_sub_pd(iy1,jy0);
755 dz10 = _mm_sub_pd(iz1,jz0);
756 dx20 = _mm_sub_pd(ix2,jx0);
757 dy20 = _mm_sub_pd(iy2,jy0);
758 dz20 = _mm_sub_pd(iz2,jz0);
759 dx30 = _mm_sub_pd(ix3,jx0);
760 dy30 = _mm_sub_pd(iy3,jy0);
761 dz30 = _mm_sub_pd(iz3,jz0);
763 /* Calculate squared distance and things based on it */
764 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
765 rsq10 = gmx_mm_calc_rsq_pd(dx10,dy10,dz10);
766 rsq20 = gmx_mm_calc_rsq_pd(dx20,dy20,dz20);
767 rsq30 = gmx_mm_calc_rsq_pd(dx30,dy30,dz30);
769 rinv10 = gmx_mm_invsqrt_pd(rsq10);
770 rinv20 = gmx_mm_invsqrt_pd(rsq20);
771 rinv30 = gmx_mm_invsqrt_pd(rsq30);
773 rinvsq00 = gmx_mm_inv_pd(rsq00);
774 rinvsq10 = _mm_mul_pd(rinv10,rinv10);
775 rinvsq20 = _mm_mul_pd(rinv20,rinv20);
776 rinvsq30 = _mm_mul_pd(rinv30,rinv30);
778 /* Load parameters for j particles */
779 jq0 = gmx_mm_load_2real_swizzle_pd(charge+jnrA+0,charge+jnrB+0);
780 vdwjidx0A = 2*vdwtype[jnrA+0];
781 vdwjidx0B = 2*vdwtype[jnrB+0];
783 fjx0 = _mm_setzero_pd();
784 fjy0 = _mm_setzero_pd();
785 fjz0 = _mm_setzero_pd();
787 /**************************
788 * CALCULATE INTERACTIONS *
789 **************************/
791 if (gmx_mm_any_lt(rsq00,rcutoff2))
794 /* Compute parameters for interactions between i and j atoms */
795 gmx_mm_load_2pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,
796 vdwparam+vdwioffset0+vdwjidx0B,&c6_00,&c12_00);
798 /* LENNARD-JONES DISPERSION/REPULSION */
800 rinvsix = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
801 fvdw = _mm_mul_pd(_mm_msub_pd(c12_00,rinvsix,c6_00),_mm_mul_pd(rinvsix,rinvsq00));
803 cutoff_mask = _mm_cmplt_pd(rsq00,rcutoff2);
807 fscal = _mm_and_pd(fscal,cutoff_mask);
809 /* Update vectorial force */
810 fix0 = _mm_macc_pd(dx00,fscal,fix0);
811 fiy0 = _mm_macc_pd(dy00,fscal,fiy0);
812 fiz0 = _mm_macc_pd(dz00,fscal,fiz0);
814 fjx0 = _mm_macc_pd(dx00,fscal,fjx0);
815 fjy0 = _mm_macc_pd(dy00,fscal,fjy0);
816 fjz0 = _mm_macc_pd(dz00,fscal,fjz0);
820 /**************************
821 * CALCULATE INTERACTIONS *
822 **************************/
824 if (gmx_mm_any_lt(rsq10,rcutoff2))
827 /* Compute parameters for interactions between i and j atoms */
828 qq10 = _mm_mul_pd(iq1,jq0);
830 /* REACTION-FIELD ELECTROSTATICS */
831 felec = _mm_mul_pd(qq10,_mm_msub_pd(rinv10,rinvsq10,krf2));
833 cutoff_mask = _mm_cmplt_pd(rsq10,rcutoff2);
837 fscal = _mm_and_pd(fscal,cutoff_mask);
839 /* Update vectorial force */
840 fix1 = _mm_macc_pd(dx10,fscal,fix1);
841 fiy1 = _mm_macc_pd(dy10,fscal,fiy1);
842 fiz1 = _mm_macc_pd(dz10,fscal,fiz1);
844 fjx0 = _mm_macc_pd(dx10,fscal,fjx0);
845 fjy0 = _mm_macc_pd(dy10,fscal,fjy0);
846 fjz0 = _mm_macc_pd(dz10,fscal,fjz0);
850 /**************************
851 * CALCULATE INTERACTIONS *
852 **************************/
854 if (gmx_mm_any_lt(rsq20,rcutoff2))
857 /* Compute parameters for interactions between i and j atoms */
858 qq20 = _mm_mul_pd(iq2,jq0);
860 /* REACTION-FIELD ELECTROSTATICS */
861 felec = _mm_mul_pd(qq20,_mm_msub_pd(rinv20,rinvsq20,krf2));
863 cutoff_mask = _mm_cmplt_pd(rsq20,rcutoff2);
867 fscal = _mm_and_pd(fscal,cutoff_mask);
869 /* Update vectorial force */
870 fix2 = _mm_macc_pd(dx20,fscal,fix2);
871 fiy2 = _mm_macc_pd(dy20,fscal,fiy2);
872 fiz2 = _mm_macc_pd(dz20,fscal,fiz2);
874 fjx0 = _mm_macc_pd(dx20,fscal,fjx0);
875 fjy0 = _mm_macc_pd(dy20,fscal,fjy0);
876 fjz0 = _mm_macc_pd(dz20,fscal,fjz0);
880 /**************************
881 * CALCULATE INTERACTIONS *
882 **************************/
884 if (gmx_mm_any_lt(rsq30,rcutoff2))
887 /* Compute parameters for interactions between i and j atoms */
888 qq30 = _mm_mul_pd(iq3,jq0);
890 /* REACTION-FIELD ELECTROSTATICS */
891 felec = _mm_mul_pd(qq30,_mm_msub_pd(rinv30,rinvsq30,krf2));
893 cutoff_mask = _mm_cmplt_pd(rsq30,rcutoff2);
897 fscal = _mm_and_pd(fscal,cutoff_mask);
899 /* Update vectorial force */
900 fix3 = _mm_macc_pd(dx30,fscal,fix3);
901 fiy3 = _mm_macc_pd(dy30,fscal,fiy3);
902 fiz3 = _mm_macc_pd(dz30,fscal,fiz3);
904 fjx0 = _mm_macc_pd(dx30,fscal,fjx0);
905 fjy0 = _mm_macc_pd(dy30,fscal,fjy0);
906 fjz0 = _mm_macc_pd(dz30,fscal,fjz0);
910 gmx_mm_decrement_1rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,fjx0,fjy0,fjz0);
912 /* Inner loop uses 135 flops */
919 j_coord_offsetA = DIM*jnrA;
921 /* load j atom coordinates */
922 gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
925 /* Calculate displacement vector */
926 dx00 = _mm_sub_pd(ix0,jx0);
927 dy00 = _mm_sub_pd(iy0,jy0);
928 dz00 = _mm_sub_pd(iz0,jz0);
929 dx10 = _mm_sub_pd(ix1,jx0);
930 dy10 = _mm_sub_pd(iy1,jy0);
931 dz10 = _mm_sub_pd(iz1,jz0);
932 dx20 = _mm_sub_pd(ix2,jx0);
933 dy20 = _mm_sub_pd(iy2,jy0);
934 dz20 = _mm_sub_pd(iz2,jz0);
935 dx30 = _mm_sub_pd(ix3,jx0);
936 dy30 = _mm_sub_pd(iy3,jy0);
937 dz30 = _mm_sub_pd(iz3,jz0);
939 /* Calculate squared distance and things based on it */
940 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
941 rsq10 = gmx_mm_calc_rsq_pd(dx10,dy10,dz10);
942 rsq20 = gmx_mm_calc_rsq_pd(dx20,dy20,dz20);
943 rsq30 = gmx_mm_calc_rsq_pd(dx30,dy30,dz30);
945 rinv10 = gmx_mm_invsqrt_pd(rsq10);
946 rinv20 = gmx_mm_invsqrt_pd(rsq20);
947 rinv30 = gmx_mm_invsqrt_pd(rsq30);
949 rinvsq00 = gmx_mm_inv_pd(rsq00);
950 rinvsq10 = _mm_mul_pd(rinv10,rinv10);
951 rinvsq20 = _mm_mul_pd(rinv20,rinv20);
952 rinvsq30 = _mm_mul_pd(rinv30,rinv30);
954 /* Load parameters for j particles */
955 jq0 = _mm_load_sd(charge+jnrA+0);
956 vdwjidx0A = 2*vdwtype[jnrA+0];
958 fjx0 = _mm_setzero_pd();
959 fjy0 = _mm_setzero_pd();
960 fjz0 = _mm_setzero_pd();
962 /**************************
963 * CALCULATE INTERACTIONS *
964 **************************/
966 if (gmx_mm_any_lt(rsq00,rcutoff2))
969 /* Compute parameters for interactions between i and j atoms */
970 gmx_mm_load_1pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,&c6_00,&c12_00);
972 /* LENNARD-JONES DISPERSION/REPULSION */
974 rinvsix = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
975 fvdw = _mm_mul_pd(_mm_msub_pd(c12_00,rinvsix,c6_00),_mm_mul_pd(rinvsix,rinvsq00));
977 cutoff_mask = _mm_cmplt_pd(rsq00,rcutoff2);
981 fscal = _mm_and_pd(fscal,cutoff_mask);
983 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
985 /* Update vectorial force */
986 fix0 = _mm_macc_pd(dx00,fscal,fix0);
987 fiy0 = _mm_macc_pd(dy00,fscal,fiy0);
988 fiz0 = _mm_macc_pd(dz00,fscal,fiz0);
990 fjx0 = _mm_macc_pd(dx00,fscal,fjx0);
991 fjy0 = _mm_macc_pd(dy00,fscal,fjy0);
992 fjz0 = _mm_macc_pd(dz00,fscal,fjz0);
996 /**************************
997 * CALCULATE INTERACTIONS *
998 **************************/
1000 if (gmx_mm_any_lt(rsq10,rcutoff2))
1003 /* Compute parameters for interactions between i and j atoms */
1004 qq10 = _mm_mul_pd(iq1,jq0);
1006 /* REACTION-FIELD ELECTROSTATICS */
1007 felec = _mm_mul_pd(qq10,_mm_msub_pd(rinv10,rinvsq10,krf2));
1009 cutoff_mask = _mm_cmplt_pd(rsq10,rcutoff2);
1013 fscal = _mm_and_pd(fscal,cutoff_mask);
1015 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
1017 /* Update vectorial force */
1018 fix1 = _mm_macc_pd(dx10,fscal,fix1);
1019 fiy1 = _mm_macc_pd(dy10,fscal,fiy1);
1020 fiz1 = _mm_macc_pd(dz10,fscal,fiz1);
1022 fjx0 = _mm_macc_pd(dx10,fscal,fjx0);
1023 fjy0 = _mm_macc_pd(dy10,fscal,fjy0);
1024 fjz0 = _mm_macc_pd(dz10,fscal,fjz0);
1028 /**************************
1029 * CALCULATE INTERACTIONS *
1030 **************************/
1032 if (gmx_mm_any_lt(rsq20,rcutoff2))
1035 /* Compute parameters for interactions between i and j atoms */
1036 qq20 = _mm_mul_pd(iq2,jq0);
1038 /* REACTION-FIELD ELECTROSTATICS */
1039 felec = _mm_mul_pd(qq20,_mm_msub_pd(rinv20,rinvsq20,krf2));
1041 cutoff_mask = _mm_cmplt_pd(rsq20,rcutoff2);
1045 fscal = _mm_and_pd(fscal,cutoff_mask);
1047 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
1049 /* Update vectorial force */
1050 fix2 = _mm_macc_pd(dx20,fscal,fix2);
1051 fiy2 = _mm_macc_pd(dy20,fscal,fiy2);
1052 fiz2 = _mm_macc_pd(dz20,fscal,fiz2);
1054 fjx0 = _mm_macc_pd(dx20,fscal,fjx0);
1055 fjy0 = _mm_macc_pd(dy20,fscal,fjy0);
1056 fjz0 = _mm_macc_pd(dz20,fscal,fjz0);
1060 /**************************
1061 * CALCULATE INTERACTIONS *
1062 **************************/
1064 if (gmx_mm_any_lt(rsq30,rcutoff2))
1067 /* Compute parameters for interactions between i and j atoms */
1068 qq30 = _mm_mul_pd(iq3,jq0);
1070 /* REACTION-FIELD ELECTROSTATICS */
1071 felec = _mm_mul_pd(qq30,_mm_msub_pd(rinv30,rinvsq30,krf2));
1073 cutoff_mask = _mm_cmplt_pd(rsq30,rcutoff2);
1077 fscal = _mm_and_pd(fscal,cutoff_mask);
1079 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
1081 /* Update vectorial force */
1082 fix3 = _mm_macc_pd(dx30,fscal,fix3);
1083 fiy3 = _mm_macc_pd(dy30,fscal,fiy3);
1084 fiz3 = _mm_macc_pd(dz30,fscal,fiz3);
1086 fjx0 = _mm_macc_pd(dx30,fscal,fjx0);
1087 fjy0 = _mm_macc_pd(dy30,fscal,fjy0);
1088 fjz0 = _mm_macc_pd(dz30,fscal,fjz0);
1092 gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,fjx0,fjy0,fjz0);
1094 /* Inner loop uses 135 flops */
1097 /* End of innermost loop */
1099 gmx_mm_update_iforce_4atom_swizzle_pd(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
1100 f+i_coord_offset,fshift+i_shift_offset);
1102 /* Increment number of inner iterations */
1103 inneriter += j_index_end - j_index_start;
1105 /* Outer loop uses 24 flops */
1108 /* Increment number of outer iterations */
1111 /* Update outer/inner flops */
1113 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W4_F,outeriter*24 + inneriter*135);