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_single kernel generator.
44 #include "../nb_kernel.h"
45 #include "types/simple.h"
49 #include "gromacs/simd/math_x86_avx_128_fma_single.h"
50 #include "kernelutil_x86_avx_128_fma_single.h"
53 * Gromacs nonbonded kernel: nb_kernel_ElecEwSh_VdwLJSh_GeomW4P1_VF_avx_128_fma_single
54 * Electrostatics interaction: Ewald
55 * VdW interaction: LennardJones
56 * Geometry: Water4-Particle
57 * Calculate force/pot: PotentialAndForce
60 nb_kernel_ElecEwSh_VdwLJSh_GeomW4P1_VF_avx_128_fma_single
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,C,D refer to j loop unrolling done with AVX_128, e.g. for the four 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;
76 int jnrA,jnrB,jnrC,jnrD;
77 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
78 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
79 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
81 real *shiftvec,*fshift,*x,*f;
82 real *fjptrA,*fjptrB,*fjptrC,*fjptrD;
84 __m128 fscal,rcutoff,rcutoff2,jidxall;
86 __m128 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
88 __m128 ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
90 __m128 ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
92 __m128 ix3,iy3,iz3,fix3,fiy3,fiz3,iq3,isai3;
93 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
94 __m128 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
95 __m128 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
96 __m128 dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
97 __m128 dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
98 __m128 dx30,dy30,dz30,rsq30,rinv30,rinvsq30,r30,qq30,c6_30,c12_30;
99 __m128 velec,felec,velecsum,facel,crf,krf,krf2;
102 __m128 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
105 __m128 one_sixth = _mm_set1_ps(1.0/6.0);
106 __m128 one_twelfth = _mm_set1_ps(1.0/12.0);
108 __m128 ewtabscale,eweps,twoeweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
109 __m128 beta,beta2,beta3,zeta2,pmecorrF,pmecorrV,rinv3;
111 __m128 dummy_mask,cutoff_mask;
112 __m128 signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
113 __m128 one = _mm_set1_ps(1.0);
114 __m128 two = _mm_set1_ps(2.0);
120 jindex = nlist->jindex;
122 shiftidx = nlist->shift;
124 shiftvec = fr->shift_vec[0];
125 fshift = fr->fshift[0];
126 facel = _mm_set1_ps(fr->epsfac);
127 charge = mdatoms->chargeA;
128 nvdwtype = fr->ntype;
130 vdwtype = mdatoms->typeA;
132 sh_ewald = _mm_set1_ps(fr->ic->sh_ewald);
133 beta = _mm_set1_ps(fr->ic->ewaldcoeff_q);
134 beta2 = _mm_mul_ps(beta,beta);
135 beta3 = _mm_mul_ps(beta,beta2);
136 ewtab = fr->ic->tabq_coul_FDV0;
137 ewtabscale = _mm_set1_ps(fr->ic->tabq_scale);
138 ewtabhalfspace = _mm_set1_ps(0.5/fr->ic->tabq_scale);
140 /* Setup water-specific parameters */
141 inr = nlist->iinr[0];
142 iq1 = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+1]));
143 iq2 = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+2]));
144 iq3 = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+3]));
145 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
147 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
148 rcutoff_scalar = fr->rcoulomb;
149 rcutoff = _mm_set1_ps(rcutoff_scalar);
150 rcutoff2 = _mm_mul_ps(rcutoff,rcutoff);
152 sh_vdw_invrcut6 = _mm_set1_ps(fr->ic->sh_invrc6);
153 rvdw = _mm_set1_ps(fr->rvdw);
155 /* Avoid stupid compiler warnings */
156 jnrA = jnrB = jnrC = jnrD = 0;
165 for(iidx=0;iidx<4*DIM;iidx++)
170 /* Start outer loop over neighborlists */
171 for(iidx=0; iidx<nri; iidx++)
173 /* Load shift vector for this list */
174 i_shift_offset = DIM*shiftidx[iidx];
176 /* Load limits for loop over neighbors */
177 j_index_start = jindex[iidx];
178 j_index_end = jindex[iidx+1];
180 /* Get outer coordinate index */
182 i_coord_offset = DIM*inr;
184 /* Load i particle coords and add shift vector */
185 gmx_mm_load_shift_and_4rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,
186 &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
188 fix0 = _mm_setzero_ps();
189 fiy0 = _mm_setzero_ps();
190 fiz0 = _mm_setzero_ps();
191 fix1 = _mm_setzero_ps();
192 fiy1 = _mm_setzero_ps();
193 fiz1 = _mm_setzero_ps();
194 fix2 = _mm_setzero_ps();
195 fiy2 = _mm_setzero_ps();
196 fiz2 = _mm_setzero_ps();
197 fix3 = _mm_setzero_ps();
198 fiy3 = _mm_setzero_ps();
199 fiz3 = _mm_setzero_ps();
201 /* Reset potential sums */
202 velecsum = _mm_setzero_ps();
203 vvdwsum = _mm_setzero_ps();
205 /* Start inner kernel loop */
206 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
209 /* Get j neighbor index, and coordinate index */
214 j_coord_offsetA = DIM*jnrA;
215 j_coord_offsetB = DIM*jnrB;
216 j_coord_offsetC = DIM*jnrC;
217 j_coord_offsetD = DIM*jnrD;
219 /* load j atom coordinates */
220 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
221 x+j_coord_offsetC,x+j_coord_offsetD,
224 /* Calculate displacement vector */
225 dx00 = _mm_sub_ps(ix0,jx0);
226 dy00 = _mm_sub_ps(iy0,jy0);
227 dz00 = _mm_sub_ps(iz0,jz0);
228 dx10 = _mm_sub_ps(ix1,jx0);
229 dy10 = _mm_sub_ps(iy1,jy0);
230 dz10 = _mm_sub_ps(iz1,jz0);
231 dx20 = _mm_sub_ps(ix2,jx0);
232 dy20 = _mm_sub_ps(iy2,jy0);
233 dz20 = _mm_sub_ps(iz2,jz0);
234 dx30 = _mm_sub_ps(ix3,jx0);
235 dy30 = _mm_sub_ps(iy3,jy0);
236 dz30 = _mm_sub_ps(iz3,jz0);
238 /* Calculate squared distance and things based on it */
239 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
240 rsq10 = gmx_mm_calc_rsq_ps(dx10,dy10,dz10);
241 rsq20 = gmx_mm_calc_rsq_ps(dx20,dy20,dz20);
242 rsq30 = gmx_mm_calc_rsq_ps(dx30,dy30,dz30);
244 rinv10 = gmx_mm_invsqrt_ps(rsq10);
245 rinv20 = gmx_mm_invsqrt_ps(rsq20);
246 rinv30 = gmx_mm_invsqrt_ps(rsq30);
248 rinvsq00 = gmx_mm_inv_ps(rsq00);
249 rinvsq10 = _mm_mul_ps(rinv10,rinv10);
250 rinvsq20 = _mm_mul_ps(rinv20,rinv20);
251 rinvsq30 = _mm_mul_ps(rinv30,rinv30);
253 /* Load parameters for j particles */
254 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
255 charge+jnrC+0,charge+jnrD+0);
256 vdwjidx0A = 2*vdwtype[jnrA+0];
257 vdwjidx0B = 2*vdwtype[jnrB+0];
258 vdwjidx0C = 2*vdwtype[jnrC+0];
259 vdwjidx0D = 2*vdwtype[jnrD+0];
261 fjx0 = _mm_setzero_ps();
262 fjy0 = _mm_setzero_ps();
263 fjz0 = _mm_setzero_ps();
265 /**************************
266 * CALCULATE INTERACTIONS *
267 **************************/
269 if (gmx_mm_any_lt(rsq00,rcutoff2))
272 /* Compute parameters for interactions between i and j atoms */
273 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
274 vdwparam+vdwioffset0+vdwjidx0B,
275 vdwparam+vdwioffset0+vdwjidx0C,
276 vdwparam+vdwioffset0+vdwjidx0D,
279 /* LENNARD-JONES DISPERSION/REPULSION */
281 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
282 vvdw6 = _mm_mul_ps(c6_00,rinvsix);
283 vvdw12 = _mm_mul_ps(c12_00,_mm_mul_ps(rinvsix,rinvsix));
284 vvdw = _mm_msub_ps(_mm_nmacc_ps(c12_00,_mm_mul_ps(sh_vdw_invrcut6,sh_vdw_invrcut6),vvdw12),one_twelfth,
285 _mm_mul_ps( _mm_nmacc_ps(c6_00,sh_vdw_invrcut6,vvdw6),one_sixth));
286 fvdw = _mm_mul_ps(_mm_sub_ps(vvdw12,vvdw6),rinvsq00);
288 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
290 /* Update potential sum for this i atom from the interaction with this j atom. */
291 vvdw = _mm_and_ps(vvdw,cutoff_mask);
292 vvdwsum = _mm_add_ps(vvdwsum,vvdw);
296 fscal = _mm_and_ps(fscal,cutoff_mask);
298 /* Update vectorial force */
299 fix0 = _mm_macc_ps(dx00,fscal,fix0);
300 fiy0 = _mm_macc_ps(dy00,fscal,fiy0);
301 fiz0 = _mm_macc_ps(dz00,fscal,fiz0);
303 fjx0 = _mm_macc_ps(dx00,fscal,fjx0);
304 fjy0 = _mm_macc_ps(dy00,fscal,fjy0);
305 fjz0 = _mm_macc_ps(dz00,fscal,fjz0);
309 /**************************
310 * CALCULATE INTERACTIONS *
311 **************************/
313 if (gmx_mm_any_lt(rsq10,rcutoff2))
316 r10 = _mm_mul_ps(rsq10,rinv10);
318 /* Compute parameters for interactions between i and j atoms */
319 qq10 = _mm_mul_ps(iq1,jq0);
321 /* EWALD ELECTROSTATICS */
323 /* Analytical PME correction */
324 zeta2 = _mm_mul_ps(beta2,rsq10);
325 rinv3 = _mm_mul_ps(rinvsq10,rinv10);
326 pmecorrF = gmx_mm_pmecorrF_ps(zeta2);
327 felec = _mm_macc_ps(pmecorrF,beta3,rinv3);
328 felec = _mm_mul_ps(qq10,felec);
329 pmecorrV = gmx_mm_pmecorrV_ps(zeta2);
330 velec = _mm_nmacc_ps(pmecorrV,beta,_mm_sub_ps(rinv10,sh_ewald));
331 velec = _mm_mul_ps(qq10,velec);
333 cutoff_mask = _mm_cmplt_ps(rsq10,rcutoff2);
335 /* Update potential sum for this i atom from the interaction with this j atom. */
336 velec = _mm_and_ps(velec,cutoff_mask);
337 velecsum = _mm_add_ps(velecsum,velec);
341 fscal = _mm_and_ps(fscal,cutoff_mask);
343 /* Update vectorial force */
344 fix1 = _mm_macc_ps(dx10,fscal,fix1);
345 fiy1 = _mm_macc_ps(dy10,fscal,fiy1);
346 fiz1 = _mm_macc_ps(dz10,fscal,fiz1);
348 fjx0 = _mm_macc_ps(dx10,fscal,fjx0);
349 fjy0 = _mm_macc_ps(dy10,fscal,fjy0);
350 fjz0 = _mm_macc_ps(dz10,fscal,fjz0);
354 /**************************
355 * CALCULATE INTERACTIONS *
356 **************************/
358 if (gmx_mm_any_lt(rsq20,rcutoff2))
361 r20 = _mm_mul_ps(rsq20,rinv20);
363 /* Compute parameters for interactions between i and j atoms */
364 qq20 = _mm_mul_ps(iq2,jq0);
366 /* EWALD ELECTROSTATICS */
368 /* Analytical PME correction */
369 zeta2 = _mm_mul_ps(beta2,rsq20);
370 rinv3 = _mm_mul_ps(rinvsq20,rinv20);
371 pmecorrF = gmx_mm_pmecorrF_ps(zeta2);
372 felec = _mm_macc_ps(pmecorrF,beta3,rinv3);
373 felec = _mm_mul_ps(qq20,felec);
374 pmecorrV = gmx_mm_pmecorrV_ps(zeta2);
375 velec = _mm_nmacc_ps(pmecorrV,beta,_mm_sub_ps(rinv20,sh_ewald));
376 velec = _mm_mul_ps(qq20,velec);
378 cutoff_mask = _mm_cmplt_ps(rsq20,rcutoff2);
380 /* Update potential sum for this i atom from the interaction with this j atom. */
381 velec = _mm_and_ps(velec,cutoff_mask);
382 velecsum = _mm_add_ps(velecsum,velec);
386 fscal = _mm_and_ps(fscal,cutoff_mask);
388 /* Update vectorial force */
389 fix2 = _mm_macc_ps(dx20,fscal,fix2);
390 fiy2 = _mm_macc_ps(dy20,fscal,fiy2);
391 fiz2 = _mm_macc_ps(dz20,fscal,fiz2);
393 fjx0 = _mm_macc_ps(dx20,fscal,fjx0);
394 fjy0 = _mm_macc_ps(dy20,fscal,fjy0);
395 fjz0 = _mm_macc_ps(dz20,fscal,fjz0);
399 /**************************
400 * CALCULATE INTERACTIONS *
401 **************************/
403 if (gmx_mm_any_lt(rsq30,rcutoff2))
406 r30 = _mm_mul_ps(rsq30,rinv30);
408 /* Compute parameters for interactions between i and j atoms */
409 qq30 = _mm_mul_ps(iq3,jq0);
411 /* EWALD ELECTROSTATICS */
413 /* Analytical PME correction */
414 zeta2 = _mm_mul_ps(beta2,rsq30);
415 rinv3 = _mm_mul_ps(rinvsq30,rinv30);
416 pmecorrF = gmx_mm_pmecorrF_ps(zeta2);
417 felec = _mm_macc_ps(pmecorrF,beta3,rinv3);
418 felec = _mm_mul_ps(qq30,felec);
419 pmecorrV = gmx_mm_pmecorrV_ps(zeta2);
420 velec = _mm_nmacc_ps(pmecorrV,beta,_mm_sub_ps(rinv30,sh_ewald));
421 velec = _mm_mul_ps(qq30,velec);
423 cutoff_mask = _mm_cmplt_ps(rsq30,rcutoff2);
425 /* Update potential sum for this i atom from the interaction with this j atom. */
426 velec = _mm_and_ps(velec,cutoff_mask);
427 velecsum = _mm_add_ps(velecsum,velec);
431 fscal = _mm_and_ps(fscal,cutoff_mask);
433 /* Update vectorial force */
434 fix3 = _mm_macc_ps(dx30,fscal,fix3);
435 fiy3 = _mm_macc_ps(dy30,fscal,fiy3);
436 fiz3 = _mm_macc_ps(dz30,fscal,fiz3);
438 fjx0 = _mm_macc_ps(dx30,fscal,fjx0);
439 fjy0 = _mm_macc_ps(dy30,fscal,fjy0);
440 fjz0 = _mm_macc_ps(dz30,fscal,fjz0);
444 fjptrA = f+j_coord_offsetA;
445 fjptrB = f+j_coord_offsetB;
446 fjptrC = f+j_coord_offsetC;
447 fjptrD = f+j_coord_offsetD;
449 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjx0,fjy0,fjz0);
451 /* Inner loop uses 143 flops */
457 /* Get j neighbor index, and coordinate index */
458 jnrlistA = jjnr[jidx];
459 jnrlistB = jjnr[jidx+1];
460 jnrlistC = jjnr[jidx+2];
461 jnrlistD = jjnr[jidx+3];
462 /* Sign of each element will be negative for non-real atoms.
463 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
464 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
466 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
467 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
468 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
469 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
470 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
471 j_coord_offsetA = DIM*jnrA;
472 j_coord_offsetB = DIM*jnrB;
473 j_coord_offsetC = DIM*jnrC;
474 j_coord_offsetD = DIM*jnrD;
476 /* load j atom coordinates */
477 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
478 x+j_coord_offsetC,x+j_coord_offsetD,
481 /* Calculate displacement vector */
482 dx00 = _mm_sub_ps(ix0,jx0);
483 dy00 = _mm_sub_ps(iy0,jy0);
484 dz00 = _mm_sub_ps(iz0,jz0);
485 dx10 = _mm_sub_ps(ix1,jx0);
486 dy10 = _mm_sub_ps(iy1,jy0);
487 dz10 = _mm_sub_ps(iz1,jz0);
488 dx20 = _mm_sub_ps(ix2,jx0);
489 dy20 = _mm_sub_ps(iy2,jy0);
490 dz20 = _mm_sub_ps(iz2,jz0);
491 dx30 = _mm_sub_ps(ix3,jx0);
492 dy30 = _mm_sub_ps(iy3,jy0);
493 dz30 = _mm_sub_ps(iz3,jz0);
495 /* Calculate squared distance and things based on it */
496 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
497 rsq10 = gmx_mm_calc_rsq_ps(dx10,dy10,dz10);
498 rsq20 = gmx_mm_calc_rsq_ps(dx20,dy20,dz20);
499 rsq30 = gmx_mm_calc_rsq_ps(dx30,dy30,dz30);
501 rinv10 = gmx_mm_invsqrt_ps(rsq10);
502 rinv20 = gmx_mm_invsqrt_ps(rsq20);
503 rinv30 = gmx_mm_invsqrt_ps(rsq30);
505 rinvsq00 = gmx_mm_inv_ps(rsq00);
506 rinvsq10 = _mm_mul_ps(rinv10,rinv10);
507 rinvsq20 = _mm_mul_ps(rinv20,rinv20);
508 rinvsq30 = _mm_mul_ps(rinv30,rinv30);
510 /* Load parameters for j particles */
511 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
512 charge+jnrC+0,charge+jnrD+0);
513 vdwjidx0A = 2*vdwtype[jnrA+0];
514 vdwjidx0B = 2*vdwtype[jnrB+0];
515 vdwjidx0C = 2*vdwtype[jnrC+0];
516 vdwjidx0D = 2*vdwtype[jnrD+0];
518 fjx0 = _mm_setzero_ps();
519 fjy0 = _mm_setzero_ps();
520 fjz0 = _mm_setzero_ps();
522 /**************************
523 * CALCULATE INTERACTIONS *
524 **************************/
526 if (gmx_mm_any_lt(rsq00,rcutoff2))
529 /* Compute parameters for interactions between i and j atoms */
530 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
531 vdwparam+vdwioffset0+vdwjidx0B,
532 vdwparam+vdwioffset0+vdwjidx0C,
533 vdwparam+vdwioffset0+vdwjidx0D,
536 /* LENNARD-JONES DISPERSION/REPULSION */
538 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
539 vvdw6 = _mm_mul_ps(c6_00,rinvsix);
540 vvdw12 = _mm_mul_ps(c12_00,_mm_mul_ps(rinvsix,rinvsix));
541 vvdw = _mm_msub_ps(_mm_nmacc_ps(c12_00,_mm_mul_ps(sh_vdw_invrcut6,sh_vdw_invrcut6),vvdw12),one_twelfth,
542 _mm_mul_ps( _mm_nmacc_ps(c6_00,sh_vdw_invrcut6,vvdw6),one_sixth));
543 fvdw = _mm_mul_ps(_mm_sub_ps(vvdw12,vvdw6),rinvsq00);
545 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
547 /* Update potential sum for this i atom from the interaction with this j atom. */
548 vvdw = _mm_and_ps(vvdw,cutoff_mask);
549 vvdw = _mm_andnot_ps(dummy_mask,vvdw);
550 vvdwsum = _mm_add_ps(vvdwsum,vvdw);
554 fscal = _mm_and_ps(fscal,cutoff_mask);
556 fscal = _mm_andnot_ps(dummy_mask,fscal);
558 /* Update vectorial force */
559 fix0 = _mm_macc_ps(dx00,fscal,fix0);
560 fiy0 = _mm_macc_ps(dy00,fscal,fiy0);
561 fiz0 = _mm_macc_ps(dz00,fscal,fiz0);
563 fjx0 = _mm_macc_ps(dx00,fscal,fjx0);
564 fjy0 = _mm_macc_ps(dy00,fscal,fjy0);
565 fjz0 = _mm_macc_ps(dz00,fscal,fjz0);
569 /**************************
570 * CALCULATE INTERACTIONS *
571 **************************/
573 if (gmx_mm_any_lt(rsq10,rcutoff2))
576 r10 = _mm_mul_ps(rsq10,rinv10);
577 r10 = _mm_andnot_ps(dummy_mask,r10);
579 /* Compute parameters for interactions between i and j atoms */
580 qq10 = _mm_mul_ps(iq1,jq0);
582 /* EWALD ELECTROSTATICS */
584 /* Analytical PME correction */
585 zeta2 = _mm_mul_ps(beta2,rsq10);
586 rinv3 = _mm_mul_ps(rinvsq10,rinv10);
587 pmecorrF = gmx_mm_pmecorrF_ps(zeta2);
588 felec = _mm_macc_ps(pmecorrF,beta3,rinv3);
589 felec = _mm_mul_ps(qq10,felec);
590 pmecorrV = gmx_mm_pmecorrV_ps(zeta2);
591 velec = _mm_nmacc_ps(pmecorrV,beta,_mm_sub_ps(rinv10,sh_ewald));
592 velec = _mm_mul_ps(qq10,velec);
594 cutoff_mask = _mm_cmplt_ps(rsq10,rcutoff2);
596 /* Update potential sum for this i atom from the interaction with this j atom. */
597 velec = _mm_and_ps(velec,cutoff_mask);
598 velec = _mm_andnot_ps(dummy_mask,velec);
599 velecsum = _mm_add_ps(velecsum,velec);
603 fscal = _mm_and_ps(fscal,cutoff_mask);
605 fscal = _mm_andnot_ps(dummy_mask,fscal);
607 /* Update vectorial force */
608 fix1 = _mm_macc_ps(dx10,fscal,fix1);
609 fiy1 = _mm_macc_ps(dy10,fscal,fiy1);
610 fiz1 = _mm_macc_ps(dz10,fscal,fiz1);
612 fjx0 = _mm_macc_ps(dx10,fscal,fjx0);
613 fjy0 = _mm_macc_ps(dy10,fscal,fjy0);
614 fjz0 = _mm_macc_ps(dz10,fscal,fjz0);
618 /**************************
619 * CALCULATE INTERACTIONS *
620 **************************/
622 if (gmx_mm_any_lt(rsq20,rcutoff2))
625 r20 = _mm_mul_ps(rsq20,rinv20);
626 r20 = _mm_andnot_ps(dummy_mask,r20);
628 /* Compute parameters for interactions between i and j atoms */
629 qq20 = _mm_mul_ps(iq2,jq0);
631 /* EWALD ELECTROSTATICS */
633 /* Analytical PME correction */
634 zeta2 = _mm_mul_ps(beta2,rsq20);
635 rinv3 = _mm_mul_ps(rinvsq20,rinv20);
636 pmecorrF = gmx_mm_pmecorrF_ps(zeta2);
637 felec = _mm_macc_ps(pmecorrF,beta3,rinv3);
638 felec = _mm_mul_ps(qq20,felec);
639 pmecorrV = gmx_mm_pmecorrV_ps(zeta2);
640 velec = _mm_nmacc_ps(pmecorrV,beta,_mm_sub_ps(rinv20,sh_ewald));
641 velec = _mm_mul_ps(qq20,velec);
643 cutoff_mask = _mm_cmplt_ps(rsq20,rcutoff2);
645 /* Update potential sum for this i atom from the interaction with this j atom. */
646 velec = _mm_and_ps(velec,cutoff_mask);
647 velec = _mm_andnot_ps(dummy_mask,velec);
648 velecsum = _mm_add_ps(velecsum,velec);
652 fscal = _mm_and_ps(fscal,cutoff_mask);
654 fscal = _mm_andnot_ps(dummy_mask,fscal);
656 /* Update vectorial force */
657 fix2 = _mm_macc_ps(dx20,fscal,fix2);
658 fiy2 = _mm_macc_ps(dy20,fscal,fiy2);
659 fiz2 = _mm_macc_ps(dz20,fscal,fiz2);
661 fjx0 = _mm_macc_ps(dx20,fscal,fjx0);
662 fjy0 = _mm_macc_ps(dy20,fscal,fjy0);
663 fjz0 = _mm_macc_ps(dz20,fscal,fjz0);
667 /**************************
668 * CALCULATE INTERACTIONS *
669 **************************/
671 if (gmx_mm_any_lt(rsq30,rcutoff2))
674 r30 = _mm_mul_ps(rsq30,rinv30);
675 r30 = _mm_andnot_ps(dummy_mask,r30);
677 /* Compute parameters for interactions between i and j atoms */
678 qq30 = _mm_mul_ps(iq3,jq0);
680 /* EWALD ELECTROSTATICS */
682 /* Analytical PME correction */
683 zeta2 = _mm_mul_ps(beta2,rsq30);
684 rinv3 = _mm_mul_ps(rinvsq30,rinv30);
685 pmecorrF = gmx_mm_pmecorrF_ps(zeta2);
686 felec = _mm_macc_ps(pmecorrF,beta3,rinv3);
687 felec = _mm_mul_ps(qq30,felec);
688 pmecorrV = gmx_mm_pmecorrV_ps(zeta2);
689 velec = _mm_nmacc_ps(pmecorrV,beta,_mm_sub_ps(rinv30,sh_ewald));
690 velec = _mm_mul_ps(qq30,velec);
692 cutoff_mask = _mm_cmplt_ps(rsq30,rcutoff2);
694 /* Update potential sum for this i atom from the interaction with this j atom. */
695 velec = _mm_and_ps(velec,cutoff_mask);
696 velec = _mm_andnot_ps(dummy_mask,velec);
697 velecsum = _mm_add_ps(velecsum,velec);
701 fscal = _mm_and_ps(fscal,cutoff_mask);
703 fscal = _mm_andnot_ps(dummy_mask,fscal);
705 /* Update vectorial force */
706 fix3 = _mm_macc_ps(dx30,fscal,fix3);
707 fiy3 = _mm_macc_ps(dy30,fscal,fiy3);
708 fiz3 = _mm_macc_ps(dz30,fscal,fiz3);
710 fjx0 = _mm_macc_ps(dx30,fscal,fjx0);
711 fjy0 = _mm_macc_ps(dy30,fscal,fjy0);
712 fjz0 = _mm_macc_ps(dz30,fscal,fjz0);
716 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
717 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
718 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
719 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
721 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjx0,fjy0,fjz0);
723 /* Inner loop uses 146 flops */
726 /* End of innermost loop */
728 gmx_mm_update_iforce_4atom_swizzle_ps(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
729 f+i_coord_offset,fshift+i_shift_offset);
732 /* Update potential energies */
733 gmx_mm_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
734 gmx_mm_update_1pot_ps(vvdwsum,kernel_data->energygrp_vdw+ggid);
736 /* Increment number of inner iterations */
737 inneriter += j_index_end - j_index_start;
739 /* Outer loop uses 26 flops */
742 /* Increment number of outer iterations */
745 /* Update outer/inner flops */
747 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W4_VF,outeriter*26 + inneriter*146);
750 * Gromacs nonbonded kernel: nb_kernel_ElecEwSh_VdwLJSh_GeomW4P1_F_avx_128_fma_single
751 * Electrostatics interaction: Ewald
752 * VdW interaction: LennardJones
753 * Geometry: Water4-Particle
754 * Calculate force/pot: Force
757 nb_kernel_ElecEwSh_VdwLJSh_GeomW4P1_F_avx_128_fma_single
758 (t_nblist * gmx_restrict nlist,
759 rvec * gmx_restrict xx,
760 rvec * gmx_restrict ff,
761 t_forcerec * gmx_restrict fr,
762 t_mdatoms * gmx_restrict mdatoms,
763 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
764 t_nrnb * gmx_restrict nrnb)
766 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
767 * just 0 for non-waters.
768 * Suffixes A,B,C,D refer to j loop unrolling done with AVX_128, e.g. for the four different
769 * jnr indices corresponding to data put in the four positions in the SIMD register.
771 int i_shift_offset,i_coord_offset,outeriter,inneriter;
772 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
773 int jnrA,jnrB,jnrC,jnrD;
774 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
775 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
776 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
778 real *shiftvec,*fshift,*x,*f;
779 real *fjptrA,*fjptrB,*fjptrC,*fjptrD;
781 __m128 fscal,rcutoff,rcutoff2,jidxall;
783 __m128 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
785 __m128 ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
787 __m128 ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
789 __m128 ix3,iy3,iz3,fix3,fiy3,fiz3,iq3,isai3;
790 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
791 __m128 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
792 __m128 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
793 __m128 dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
794 __m128 dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
795 __m128 dx30,dy30,dz30,rsq30,rinv30,rinvsq30,r30,qq30,c6_30,c12_30;
796 __m128 velec,felec,velecsum,facel,crf,krf,krf2;
799 __m128 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
802 __m128 one_sixth = _mm_set1_ps(1.0/6.0);
803 __m128 one_twelfth = _mm_set1_ps(1.0/12.0);
805 __m128 ewtabscale,eweps,twoeweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
806 __m128 beta,beta2,beta3,zeta2,pmecorrF,pmecorrV,rinv3;
808 __m128 dummy_mask,cutoff_mask;
809 __m128 signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
810 __m128 one = _mm_set1_ps(1.0);
811 __m128 two = _mm_set1_ps(2.0);
817 jindex = nlist->jindex;
819 shiftidx = nlist->shift;
821 shiftvec = fr->shift_vec[0];
822 fshift = fr->fshift[0];
823 facel = _mm_set1_ps(fr->epsfac);
824 charge = mdatoms->chargeA;
825 nvdwtype = fr->ntype;
827 vdwtype = mdatoms->typeA;
829 sh_ewald = _mm_set1_ps(fr->ic->sh_ewald);
830 beta = _mm_set1_ps(fr->ic->ewaldcoeff_q);
831 beta2 = _mm_mul_ps(beta,beta);
832 beta3 = _mm_mul_ps(beta,beta2);
833 ewtab = fr->ic->tabq_coul_F;
834 ewtabscale = _mm_set1_ps(fr->ic->tabq_scale);
835 ewtabhalfspace = _mm_set1_ps(0.5/fr->ic->tabq_scale);
837 /* Setup water-specific parameters */
838 inr = nlist->iinr[0];
839 iq1 = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+1]));
840 iq2 = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+2]));
841 iq3 = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+3]));
842 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
844 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
845 rcutoff_scalar = fr->rcoulomb;
846 rcutoff = _mm_set1_ps(rcutoff_scalar);
847 rcutoff2 = _mm_mul_ps(rcutoff,rcutoff);
849 sh_vdw_invrcut6 = _mm_set1_ps(fr->ic->sh_invrc6);
850 rvdw = _mm_set1_ps(fr->rvdw);
852 /* Avoid stupid compiler warnings */
853 jnrA = jnrB = jnrC = jnrD = 0;
862 for(iidx=0;iidx<4*DIM;iidx++)
867 /* Start outer loop over neighborlists */
868 for(iidx=0; iidx<nri; iidx++)
870 /* Load shift vector for this list */
871 i_shift_offset = DIM*shiftidx[iidx];
873 /* Load limits for loop over neighbors */
874 j_index_start = jindex[iidx];
875 j_index_end = jindex[iidx+1];
877 /* Get outer coordinate index */
879 i_coord_offset = DIM*inr;
881 /* Load i particle coords and add shift vector */
882 gmx_mm_load_shift_and_4rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,
883 &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
885 fix0 = _mm_setzero_ps();
886 fiy0 = _mm_setzero_ps();
887 fiz0 = _mm_setzero_ps();
888 fix1 = _mm_setzero_ps();
889 fiy1 = _mm_setzero_ps();
890 fiz1 = _mm_setzero_ps();
891 fix2 = _mm_setzero_ps();
892 fiy2 = _mm_setzero_ps();
893 fiz2 = _mm_setzero_ps();
894 fix3 = _mm_setzero_ps();
895 fiy3 = _mm_setzero_ps();
896 fiz3 = _mm_setzero_ps();
898 /* Start inner kernel loop */
899 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
902 /* Get j neighbor index, and coordinate index */
907 j_coord_offsetA = DIM*jnrA;
908 j_coord_offsetB = DIM*jnrB;
909 j_coord_offsetC = DIM*jnrC;
910 j_coord_offsetD = DIM*jnrD;
912 /* load j atom coordinates */
913 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
914 x+j_coord_offsetC,x+j_coord_offsetD,
917 /* Calculate displacement vector */
918 dx00 = _mm_sub_ps(ix0,jx0);
919 dy00 = _mm_sub_ps(iy0,jy0);
920 dz00 = _mm_sub_ps(iz0,jz0);
921 dx10 = _mm_sub_ps(ix1,jx0);
922 dy10 = _mm_sub_ps(iy1,jy0);
923 dz10 = _mm_sub_ps(iz1,jz0);
924 dx20 = _mm_sub_ps(ix2,jx0);
925 dy20 = _mm_sub_ps(iy2,jy0);
926 dz20 = _mm_sub_ps(iz2,jz0);
927 dx30 = _mm_sub_ps(ix3,jx0);
928 dy30 = _mm_sub_ps(iy3,jy0);
929 dz30 = _mm_sub_ps(iz3,jz0);
931 /* Calculate squared distance and things based on it */
932 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
933 rsq10 = gmx_mm_calc_rsq_ps(dx10,dy10,dz10);
934 rsq20 = gmx_mm_calc_rsq_ps(dx20,dy20,dz20);
935 rsq30 = gmx_mm_calc_rsq_ps(dx30,dy30,dz30);
937 rinv10 = gmx_mm_invsqrt_ps(rsq10);
938 rinv20 = gmx_mm_invsqrt_ps(rsq20);
939 rinv30 = gmx_mm_invsqrt_ps(rsq30);
941 rinvsq00 = gmx_mm_inv_ps(rsq00);
942 rinvsq10 = _mm_mul_ps(rinv10,rinv10);
943 rinvsq20 = _mm_mul_ps(rinv20,rinv20);
944 rinvsq30 = _mm_mul_ps(rinv30,rinv30);
946 /* Load parameters for j particles */
947 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
948 charge+jnrC+0,charge+jnrD+0);
949 vdwjidx0A = 2*vdwtype[jnrA+0];
950 vdwjidx0B = 2*vdwtype[jnrB+0];
951 vdwjidx0C = 2*vdwtype[jnrC+0];
952 vdwjidx0D = 2*vdwtype[jnrD+0];
954 fjx0 = _mm_setzero_ps();
955 fjy0 = _mm_setzero_ps();
956 fjz0 = _mm_setzero_ps();
958 /**************************
959 * CALCULATE INTERACTIONS *
960 **************************/
962 if (gmx_mm_any_lt(rsq00,rcutoff2))
965 /* Compute parameters for interactions between i and j atoms */
966 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
967 vdwparam+vdwioffset0+vdwjidx0B,
968 vdwparam+vdwioffset0+vdwjidx0C,
969 vdwparam+vdwioffset0+vdwjidx0D,
972 /* LENNARD-JONES DISPERSION/REPULSION */
974 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
975 fvdw = _mm_mul_ps(_mm_msub_ps(c12_00,rinvsix,c6_00),_mm_mul_ps(rinvsix,rinvsq00));
977 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
981 fscal = _mm_and_ps(fscal,cutoff_mask);
983 /* Update vectorial force */
984 fix0 = _mm_macc_ps(dx00,fscal,fix0);
985 fiy0 = _mm_macc_ps(dy00,fscal,fiy0);
986 fiz0 = _mm_macc_ps(dz00,fscal,fiz0);
988 fjx0 = _mm_macc_ps(dx00,fscal,fjx0);
989 fjy0 = _mm_macc_ps(dy00,fscal,fjy0);
990 fjz0 = _mm_macc_ps(dz00,fscal,fjz0);
994 /**************************
995 * CALCULATE INTERACTIONS *
996 **************************/
998 if (gmx_mm_any_lt(rsq10,rcutoff2))
1001 r10 = _mm_mul_ps(rsq10,rinv10);
1003 /* Compute parameters for interactions between i and j atoms */
1004 qq10 = _mm_mul_ps(iq1,jq0);
1006 /* EWALD ELECTROSTATICS */
1008 /* Analytical PME correction */
1009 zeta2 = _mm_mul_ps(beta2,rsq10);
1010 rinv3 = _mm_mul_ps(rinvsq10,rinv10);
1011 pmecorrF = gmx_mm_pmecorrF_ps(zeta2);
1012 felec = _mm_macc_ps(pmecorrF,beta3,rinv3);
1013 felec = _mm_mul_ps(qq10,felec);
1015 cutoff_mask = _mm_cmplt_ps(rsq10,rcutoff2);
1019 fscal = _mm_and_ps(fscal,cutoff_mask);
1021 /* Update vectorial force */
1022 fix1 = _mm_macc_ps(dx10,fscal,fix1);
1023 fiy1 = _mm_macc_ps(dy10,fscal,fiy1);
1024 fiz1 = _mm_macc_ps(dz10,fscal,fiz1);
1026 fjx0 = _mm_macc_ps(dx10,fscal,fjx0);
1027 fjy0 = _mm_macc_ps(dy10,fscal,fjy0);
1028 fjz0 = _mm_macc_ps(dz10,fscal,fjz0);
1032 /**************************
1033 * CALCULATE INTERACTIONS *
1034 **************************/
1036 if (gmx_mm_any_lt(rsq20,rcutoff2))
1039 r20 = _mm_mul_ps(rsq20,rinv20);
1041 /* Compute parameters for interactions between i and j atoms */
1042 qq20 = _mm_mul_ps(iq2,jq0);
1044 /* EWALD ELECTROSTATICS */
1046 /* Analytical PME correction */
1047 zeta2 = _mm_mul_ps(beta2,rsq20);
1048 rinv3 = _mm_mul_ps(rinvsq20,rinv20);
1049 pmecorrF = gmx_mm_pmecorrF_ps(zeta2);
1050 felec = _mm_macc_ps(pmecorrF,beta3,rinv3);
1051 felec = _mm_mul_ps(qq20,felec);
1053 cutoff_mask = _mm_cmplt_ps(rsq20,rcutoff2);
1057 fscal = _mm_and_ps(fscal,cutoff_mask);
1059 /* Update vectorial force */
1060 fix2 = _mm_macc_ps(dx20,fscal,fix2);
1061 fiy2 = _mm_macc_ps(dy20,fscal,fiy2);
1062 fiz2 = _mm_macc_ps(dz20,fscal,fiz2);
1064 fjx0 = _mm_macc_ps(dx20,fscal,fjx0);
1065 fjy0 = _mm_macc_ps(dy20,fscal,fjy0);
1066 fjz0 = _mm_macc_ps(dz20,fscal,fjz0);
1070 /**************************
1071 * CALCULATE INTERACTIONS *
1072 **************************/
1074 if (gmx_mm_any_lt(rsq30,rcutoff2))
1077 r30 = _mm_mul_ps(rsq30,rinv30);
1079 /* Compute parameters for interactions between i and j atoms */
1080 qq30 = _mm_mul_ps(iq3,jq0);
1082 /* EWALD ELECTROSTATICS */
1084 /* Analytical PME correction */
1085 zeta2 = _mm_mul_ps(beta2,rsq30);
1086 rinv3 = _mm_mul_ps(rinvsq30,rinv30);
1087 pmecorrF = gmx_mm_pmecorrF_ps(zeta2);
1088 felec = _mm_macc_ps(pmecorrF,beta3,rinv3);
1089 felec = _mm_mul_ps(qq30,felec);
1091 cutoff_mask = _mm_cmplt_ps(rsq30,rcutoff2);
1095 fscal = _mm_and_ps(fscal,cutoff_mask);
1097 /* Update vectorial force */
1098 fix3 = _mm_macc_ps(dx30,fscal,fix3);
1099 fiy3 = _mm_macc_ps(dy30,fscal,fiy3);
1100 fiz3 = _mm_macc_ps(dz30,fscal,fiz3);
1102 fjx0 = _mm_macc_ps(dx30,fscal,fjx0);
1103 fjy0 = _mm_macc_ps(dy30,fscal,fjy0);
1104 fjz0 = _mm_macc_ps(dz30,fscal,fjz0);
1108 fjptrA = f+j_coord_offsetA;
1109 fjptrB = f+j_coord_offsetB;
1110 fjptrC = f+j_coord_offsetC;
1111 fjptrD = f+j_coord_offsetD;
1113 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjx0,fjy0,fjz0);
1115 /* Inner loop uses 126 flops */
1118 if(jidx<j_index_end)
1121 /* Get j neighbor index, and coordinate index */
1122 jnrlistA = jjnr[jidx];
1123 jnrlistB = jjnr[jidx+1];
1124 jnrlistC = jjnr[jidx+2];
1125 jnrlistD = jjnr[jidx+3];
1126 /* Sign of each element will be negative for non-real atoms.
1127 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
1128 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
1130 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
1131 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
1132 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
1133 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
1134 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
1135 j_coord_offsetA = DIM*jnrA;
1136 j_coord_offsetB = DIM*jnrB;
1137 j_coord_offsetC = DIM*jnrC;
1138 j_coord_offsetD = DIM*jnrD;
1140 /* load j atom coordinates */
1141 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
1142 x+j_coord_offsetC,x+j_coord_offsetD,
1145 /* Calculate displacement vector */
1146 dx00 = _mm_sub_ps(ix0,jx0);
1147 dy00 = _mm_sub_ps(iy0,jy0);
1148 dz00 = _mm_sub_ps(iz0,jz0);
1149 dx10 = _mm_sub_ps(ix1,jx0);
1150 dy10 = _mm_sub_ps(iy1,jy0);
1151 dz10 = _mm_sub_ps(iz1,jz0);
1152 dx20 = _mm_sub_ps(ix2,jx0);
1153 dy20 = _mm_sub_ps(iy2,jy0);
1154 dz20 = _mm_sub_ps(iz2,jz0);
1155 dx30 = _mm_sub_ps(ix3,jx0);
1156 dy30 = _mm_sub_ps(iy3,jy0);
1157 dz30 = _mm_sub_ps(iz3,jz0);
1159 /* Calculate squared distance and things based on it */
1160 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
1161 rsq10 = gmx_mm_calc_rsq_ps(dx10,dy10,dz10);
1162 rsq20 = gmx_mm_calc_rsq_ps(dx20,dy20,dz20);
1163 rsq30 = gmx_mm_calc_rsq_ps(dx30,dy30,dz30);
1165 rinv10 = gmx_mm_invsqrt_ps(rsq10);
1166 rinv20 = gmx_mm_invsqrt_ps(rsq20);
1167 rinv30 = gmx_mm_invsqrt_ps(rsq30);
1169 rinvsq00 = gmx_mm_inv_ps(rsq00);
1170 rinvsq10 = _mm_mul_ps(rinv10,rinv10);
1171 rinvsq20 = _mm_mul_ps(rinv20,rinv20);
1172 rinvsq30 = _mm_mul_ps(rinv30,rinv30);
1174 /* Load parameters for j particles */
1175 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
1176 charge+jnrC+0,charge+jnrD+0);
1177 vdwjidx0A = 2*vdwtype[jnrA+0];
1178 vdwjidx0B = 2*vdwtype[jnrB+0];
1179 vdwjidx0C = 2*vdwtype[jnrC+0];
1180 vdwjidx0D = 2*vdwtype[jnrD+0];
1182 fjx0 = _mm_setzero_ps();
1183 fjy0 = _mm_setzero_ps();
1184 fjz0 = _mm_setzero_ps();
1186 /**************************
1187 * CALCULATE INTERACTIONS *
1188 **************************/
1190 if (gmx_mm_any_lt(rsq00,rcutoff2))
1193 /* Compute parameters for interactions between i and j atoms */
1194 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
1195 vdwparam+vdwioffset0+vdwjidx0B,
1196 vdwparam+vdwioffset0+vdwjidx0C,
1197 vdwparam+vdwioffset0+vdwjidx0D,
1200 /* LENNARD-JONES DISPERSION/REPULSION */
1202 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
1203 fvdw = _mm_mul_ps(_mm_msub_ps(c12_00,rinvsix,c6_00),_mm_mul_ps(rinvsix,rinvsq00));
1205 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
1209 fscal = _mm_and_ps(fscal,cutoff_mask);
1211 fscal = _mm_andnot_ps(dummy_mask,fscal);
1213 /* Update vectorial force */
1214 fix0 = _mm_macc_ps(dx00,fscal,fix0);
1215 fiy0 = _mm_macc_ps(dy00,fscal,fiy0);
1216 fiz0 = _mm_macc_ps(dz00,fscal,fiz0);
1218 fjx0 = _mm_macc_ps(dx00,fscal,fjx0);
1219 fjy0 = _mm_macc_ps(dy00,fscal,fjy0);
1220 fjz0 = _mm_macc_ps(dz00,fscal,fjz0);
1224 /**************************
1225 * CALCULATE INTERACTIONS *
1226 **************************/
1228 if (gmx_mm_any_lt(rsq10,rcutoff2))
1231 r10 = _mm_mul_ps(rsq10,rinv10);
1232 r10 = _mm_andnot_ps(dummy_mask,r10);
1234 /* Compute parameters for interactions between i and j atoms */
1235 qq10 = _mm_mul_ps(iq1,jq0);
1237 /* EWALD ELECTROSTATICS */
1239 /* Analytical PME correction */
1240 zeta2 = _mm_mul_ps(beta2,rsq10);
1241 rinv3 = _mm_mul_ps(rinvsq10,rinv10);
1242 pmecorrF = gmx_mm_pmecorrF_ps(zeta2);
1243 felec = _mm_macc_ps(pmecorrF,beta3,rinv3);
1244 felec = _mm_mul_ps(qq10,felec);
1246 cutoff_mask = _mm_cmplt_ps(rsq10,rcutoff2);
1250 fscal = _mm_and_ps(fscal,cutoff_mask);
1252 fscal = _mm_andnot_ps(dummy_mask,fscal);
1254 /* Update vectorial force */
1255 fix1 = _mm_macc_ps(dx10,fscal,fix1);
1256 fiy1 = _mm_macc_ps(dy10,fscal,fiy1);
1257 fiz1 = _mm_macc_ps(dz10,fscal,fiz1);
1259 fjx0 = _mm_macc_ps(dx10,fscal,fjx0);
1260 fjy0 = _mm_macc_ps(dy10,fscal,fjy0);
1261 fjz0 = _mm_macc_ps(dz10,fscal,fjz0);
1265 /**************************
1266 * CALCULATE INTERACTIONS *
1267 **************************/
1269 if (gmx_mm_any_lt(rsq20,rcutoff2))
1272 r20 = _mm_mul_ps(rsq20,rinv20);
1273 r20 = _mm_andnot_ps(dummy_mask,r20);
1275 /* Compute parameters for interactions between i and j atoms */
1276 qq20 = _mm_mul_ps(iq2,jq0);
1278 /* EWALD ELECTROSTATICS */
1280 /* Analytical PME correction */
1281 zeta2 = _mm_mul_ps(beta2,rsq20);
1282 rinv3 = _mm_mul_ps(rinvsq20,rinv20);
1283 pmecorrF = gmx_mm_pmecorrF_ps(zeta2);
1284 felec = _mm_macc_ps(pmecorrF,beta3,rinv3);
1285 felec = _mm_mul_ps(qq20,felec);
1287 cutoff_mask = _mm_cmplt_ps(rsq20,rcutoff2);
1291 fscal = _mm_and_ps(fscal,cutoff_mask);
1293 fscal = _mm_andnot_ps(dummy_mask,fscal);
1295 /* Update vectorial force */
1296 fix2 = _mm_macc_ps(dx20,fscal,fix2);
1297 fiy2 = _mm_macc_ps(dy20,fscal,fiy2);
1298 fiz2 = _mm_macc_ps(dz20,fscal,fiz2);
1300 fjx0 = _mm_macc_ps(dx20,fscal,fjx0);
1301 fjy0 = _mm_macc_ps(dy20,fscal,fjy0);
1302 fjz0 = _mm_macc_ps(dz20,fscal,fjz0);
1306 /**************************
1307 * CALCULATE INTERACTIONS *
1308 **************************/
1310 if (gmx_mm_any_lt(rsq30,rcutoff2))
1313 r30 = _mm_mul_ps(rsq30,rinv30);
1314 r30 = _mm_andnot_ps(dummy_mask,r30);
1316 /* Compute parameters for interactions between i and j atoms */
1317 qq30 = _mm_mul_ps(iq3,jq0);
1319 /* EWALD ELECTROSTATICS */
1321 /* Analytical PME correction */
1322 zeta2 = _mm_mul_ps(beta2,rsq30);
1323 rinv3 = _mm_mul_ps(rinvsq30,rinv30);
1324 pmecorrF = gmx_mm_pmecorrF_ps(zeta2);
1325 felec = _mm_macc_ps(pmecorrF,beta3,rinv3);
1326 felec = _mm_mul_ps(qq30,felec);
1328 cutoff_mask = _mm_cmplt_ps(rsq30,rcutoff2);
1332 fscal = _mm_and_ps(fscal,cutoff_mask);
1334 fscal = _mm_andnot_ps(dummy_mask,fscal);
1336 /* Update vectorial force */
1337 fix3 = _mm_macc_ps(dx30,fscal,fix3);
1338 fiy3 = _mm_macc_ps(dy30,fscal,fiy3);
1339 fiz3 = _mm_macc_ps(dz30,fscal,fiz3);
1341 fjx0 = _mm_macc_ps(dx30,fscal,fjx0);
1342 fjy0 = _mm_macc_ps(dy30,fscal,fjy0);
1343 fjz0 = _mm_macc_ps(dz30,fscal,fjz0);
1347 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
1348 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
1349 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
1350 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
1352 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjx0,fjy0,fjz0);
1354 /* Inner loop uses 129 flops */
1357 /* End of innermost loop */
1359 gmx_mm_update_iforce_4atom_swizzle_ps(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
1360 f+i_coord_offset,fshift+i_shift_offset);
1362 /* Increment number of inner iterations */
1363 inneriter += j_index_end - j_index_start;
1365 /* Outer loop uses 24 flops */
1368 /* Increment number of outer iterations */
1371 /* Update outer/inner flops */
1373 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W4_F,outeriter*24 + inneriter*129);