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 "gmx_math_x86_avx_128_fma_single.h"
50 #include "kernelutil_x86_avx_128_fma_single.h"
53 * Gromacs nonbonded kernel: nb_kernel_ElecRFCut_VdwLJSh_GeomW4P1_VF_avx_128_fma_single
54 * Electrostatics interaction: ReactionField
55 * VdW interaction: LennardJones
56 * Geometry: Water4-Particle
57 * Calculate force/pot: PotentialAndForce
60 nb_kernel_ElecRFCut_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);
107 __m128 dummy_mask,cutoff_mask;
108 __m128 signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
109 __m128 one = _mm_set1_ps(1.0);
110 __m128 two = _mm_set1_ps(2.0);
116 jindex = nlist->jindex;
118 shiftidx = nlist->shift;
120 shiftvec = fr->shift_vec[0];
121 fshift = fr->fshift[0];
122 facel = _mm_set1_ps(fr->epsfac);
123 charge = mdatoms->chargeA;
124 krf = _mm_set1_ps(fr->ic->k_rf);
125 krf2 = _mm_set1_ps(fr->ic->k_rf*2.0);
126 crf = _mm_set1_ps(fr->ic->c_rf);
127 nvdwtype = fr->ntype;
129 vdwtype = mdatoms->typeA;
131 /* Setup water-specific parameters */
132 inr = nlist->iinr[0];
133 iq1 = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+1]));
134 iq2 = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+2]));
135 iq3 = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+3]));
136 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
138 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
139 rcutoff_scalar = fr->rcoulomb;
140 rcutoff = _mm_set1_ps(rcutoff_scalar);
141 rcutoff2 = _mm_mul_ps(rcutoff,rcutoff);
143 sh_vdw_invrcut6 = _mm_set1_ps(fr->ic->sh_invrc6);
144 rvdw = _mm_set1_ps(fr->rvdw);
146 /* Avoid stupid compiler warnings */
147 jnrA = jnrB = jnrC = jnrD = 0;
156 for(iidx=0;iidx<4*DIM;iidx++)
161 /* Start outer loop over neighborlists */
162 for(iidx=0; iidx<nri; iidx++)
164 /* Load shift vector for this list */
165 i_shift_offset = DIM*shiftidx[iidx];
167 /* Load limits for loop over neighbors */
168 j_index_start = jindex[iidx];
169 j_index_end = jindex[iidx+1];
171 /* Get outer coordinate index */
173 i_coord_offset = DIM*inr;
175 /* Load i particle coords and add shift vector */
176 gmx_mm_load_shift_and_4rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,
177 &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
179 fix0 = _mm_setzero_ps();
180 fiy0 = _mm_setzero_ps();
181 fiz0 = _mm_setzero_ps();
182 fix1 = _mm_setzero_ps();
183 fiy1 = _mm_setzero_ps();
184 fiz1 = _mm_setzero_ps();
185 fix2 = _mm_setzero_ps();
186 fiy2 = _mm_setzero_ps();
187 fiz2 = _mm_setzero_ps();
188 fix3 = _mm_setzero_ps();
189 fiy3 = _mm_setzero_ps();
190 fiz3 = _mm_setzero_ps();
192 /* Reset potential sums */
193 velecsum = _mm_setzero_ps();
194 vvdwsum = _mm_setzero_ps();
196 /* Start inner kernel loop */
197 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
200 /* Get j neighbor index, and coordinate index */
205 j_coord_offsetA = DIM*jnrA;
206 j_coord_offsetB = DIM*jnrB;
207 j_coord_offsetC = DIM*jnrC;
208 j_coord_offsetD = DIM*jnrD;
210 /* load j atom coordinates */
211 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
212 x+j_coord_offsetC,x+j_coord_offsetD,
215 /* Calculate displacement vector */
216 dx00 = _mm_sub_ps(ix0,jx0);
217 dy00 = _mm_sub_ps(iy0,jy0);
218 dz00 = _mm_sub_ps(iz0,jz0);
219 dx10 = _mm_sub_ps(ix1,jx0);
220 dy10 = _mm_sub_ps(iy1,jy0);
221 dz10 = _mm_sub_ps(iz1,jz0);
222 dx20 = _mm_sub_ps(ix2,jx0);
223 dy20 = _mm_sub_ps(iy2,jy0);
224 dz20 = _mm_sub_ps(iz2,jz0);
225 dx30 = _mm_sub_ps(ix3,jx0);
226 dy30 = _mm_sub_ps(iy3,jy0);
227 dz30 = _mm_sub_ps(iz3,jz0);
229 /* Calculate squared distance and things based on it */
230 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
231 rsq10 = gmx_mm_calc_rsq_ps(dx10,dy10,dz10);
232 rsq20 = gmx_mm_calc_rsq_ps(dx20,dy20,dz20);
233 rsq30 = gmx_mm_calc_rsq_ps(dx30,dy30,dz30);
235 rinv10 = gmx_mm_invsqrt_ps(rsq10);
236 rinv20 = gmx_mm_invsqrt_ps(rsq20);
237 rinv30 = gmx_mm_invsqrt_ps(rsq30);
239 rinvsq00 = gmx_mm_inv_ps(rsq00);
240 rinvsq10 = _mm_mul_ps(rinv10,rinv10);
241 rinvsq20 = _mm_mul_ps(rinv20,rinv20);
242 rinvsq30 = _mm_mul_ps(rinv30,rinv30);
244 /* Load parameters for j particles */
245 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
246 charge+jnrC+0,charge+jnrD+0);
247 vdwjidx0A = 2*vdwtype[jnrA+0];
248 vdwjidx0B = 2*vdwtype[jnrB+0];
249 vdwjidx0C = 2*vdwtype[jnrC+0];
250 vdwjidx0D = 2*vdwtype[jnrD+0];
252 fjx0 = _mm_setzero_ps();
253 fjy0 = _mm_setzero_ps();
254 fjz0 = _mm_setzero_ps();
256 /**************************
257 * CALCULATE INTERACTIONS *
258 **************************/
260 if (gmx_mm_any_lt(rsq00,rcutoff2))
263 /* Compute parameters for interactions between i and j atoms */
264 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
265 vdwparam+vdwioffset0+vdwjidx0B,
266 vdwparam+vdwioffset0+vdwjidx0C,
267 vdwparam+vdwioffset0+vdwjidx0D,
270 /* LENNARD-JONES DISPERSION/REPULSION */
272 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
273 vvdw6 = _mm_mul_ps(c6_00,rinvsix);
274 vvdw12 = _mm_mul_ps(c12_00,_mm_mul_ps(rinvsix,rinvsix));
275 vvdw = _mm_msub_ps(_mm_nmacc_ps(c12_00,_mm_mul_ps(sh_vdw_invrcut6,sh_vdw_invrcut6),vvdw12),one_twelfth,
276 _mm_mul_ps( _mm_nmacc_ps(c6_00,sh_vdw_invrcut6,vvdw6),one_sixth));
277 fvdw = _mm_mul_ps(_mm_sub_ps(vvdw12,vvdw6),rinvsq00);
279 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
281 /* Update potential sum for this i atom from the interaction with this j atom. */
282 vvdw = _mm_and_ps(vvdw,cutoff_mask);
283 vvdwsum = _mm_add_ps(vvdwsum,vvdw);
287 fscal = _mm_and_ps(fscal,cutoff_mask);
289 /* Update vectorial force */
290 fix0 = _mm_macc_ps(dx00,fscal,fix0);
291 fiy0 = _mm_macc_ps(dy00,fscal,fiy0);
292 fiz0 = _mm_macc_ps(dz00,fscal,fiz0);
294 fjx0 = _mm_macc_ps(dx00,fscal,fjx0);
295 fjy0 = _mm_macc_ps(dy00,fscal,fjy0);
296 fjz0 = _mm_macc_ps(dz00,fscal,fjz0);
300 /**************************
301 * CALCULATE INTERACTIONS *
302 **************************/
304 if (gmx_mm_any_lt(rsq10,rcutoff2))
307 /* Compute parameters for interactions between i and j atoms */
308 qq10 = _mm_mul_ps(iq1,jq0);
310 /* REACTION-FIELD ELECTROSTATICS */
311 velec = _mm_mul_ps(qq10,_mm_sub_ps(_mm_macc_ps(krf,rsq10,rinv10),crf));
312 felec = _mm_mul_ps(qq10,_mm_msub_ps(rinv10,rinvsq10,krf2));
314 cutoff_mask = _mm_cmplt_ps(rsq10,rcutoff2);
316 /* Update potential sum for this i atom from the interaction with this j atom. */
317 velec = _mm_and_ps(velec,cutoff_mask);
318 velecsum = _mm_add_ps(velecsum,velec);
322 fscal = _mm_and_ps(fscal,cutoff_mask);
324 /* Update vectorial force */
325 fix1 = _mm_macc_ps(dx10,fscal,fix1);
326 fiy1 = _mm_macc_ps(dy10,fscal,fiy1);
327 fiz1 = _mm_macc_ps(dz10,fscal,fiz1);
329 fjx0 = _mm_macc_ps(dx10,fscal,fjx0);
330 fjy0 = _mm_macc_ps(dy10,fscal,fjy0);
331 fjz0 = _mm_macc_ps(dz10,fscal,fjz0);
335 /**************************
336 * CALCULATE INTERACTIONS *
337 **************************/
339 if (gmx_mm_any_lt(rsq20,rcutoff2))
342 /* Compute parameters for interactions between i and j atoms */
343 qq20 = _mm_mul_ps(iq2,jq0);
345 /* REACTION-FIELD ELECTROSTATICS */
346 velec = _mm_mul_ps(qq20,_mm_sub_ps(_mm_macc_ps(krf,rsq20,rinv20),crf));
347 felec = _mm_mul_ps(qq20,_mm_msub_ps(rinv20,rinvsq20,krf2));
349 cutoff_mask = _mm_cmplt_ps(rsq20,rcutoff2);
351 /* Update potential sum for this i atom from the interaction with this j atom. */
352 velec = _mm_and_ps(velec,cutoff_mask);
353 velecsum = _mm_add_ps(velecsum,velec);
357 fscal = _mm_and_ps(fscal,cutoff_mask);
359 /* Update vectorial force */
360 fix2 = _mm_macc_ps(dx20,fscal,fix2);
361 fiy2 = _mm_macc_ps(dy20,fscal,fiy2);
362 fiz2 = _mm_macc_ps(dz20,fscal,fiz2);
364 fjx0 = _mm_macc_ps(dx20,fscal,fjx0);
365 fjy0 = _mm_macc_ps(dy20,fscal,fjy0);
366 fjz0 = _mm_macc_ps(dz20,fscal,fjz0);
370 /**************************
371 * CALCULATE INTERACTIONS *
372 **************************/
374 if (gmx_mm_any_lt(rsq30,rcutoff2))
377 /* Compute parameters for interactions between i and j atoms */
378 qq30 = _mm_mul_ps(iq3,jq0);
380 /* REACTION-FIELD ELECTROSTATICS */
381 velec = _mm_mul_ps(qq30,_mm_sub_ps(_mm_macc_ps(krf,rsq30,rinv30),crf));
382 felec = _mm_mul_ps(qq30,_mm_msub_ps(rinv30,rinvsq30,krf2));
384 cutoff_mask = _mm_cmplt_ps(rsq30,rcutoff2);
386 /* Update potential sum for this i atom from the interaction with this j atom. */
387 velec = _mm_and_ps(velec,cutoff_mask);
388 velecsum = _mm_add_ps(velecsum,velec);
392 fscal = _mm_and_ps(fscal,cutoff_mask);
394 /* Update vectorial force */
395 fix3 = _mm_macc_ps(dx30,fscal,fix3);
396 fiy3 = _mm_macc_ps(dy30,fscal,fiy3);
397 fiz3 = _mm_macc_ps(dz30,fscal,fiz3);
399 fjx0 = _mm_macc_ps(dx30,fscal,fjx0);
400 fjy0 = _mm_macc_ps(dy30,fscal,fjy0);
401 fjz0 = _mm_macc_ps(dz30,fscal,fjz0);
405 fjptrA = f+j_coord_offsetA;
406 fjptrB = f+j_coord_offsetB;
407 fjptrC = f+j_coord_offsetC;
408 fjptrD = f+j_coord_offsetD;
410 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjx0,fjy0,fjz0);
412 /* Inner loop uses 161 flops */
418 /* Get j neighbor index, and coordinate index */
419 jnrlistA = jjnr[jidx];
420 jnrlistB = jjnr[jidx+1];
421 jnrlistC = jjnr[jidx+2];
422 jnrlistD = jjnr[jidx+3];
423 /* Sign of each element will be negative for non-real atoms.
424 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
425 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
427 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
428 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
429 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
430 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
431 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
432 j_coord_offsetA = DIM*jnrA;
433 j_coord_offsetB = DIM*jnrB;
434 j_coord_offsetC = DIM*jnrC;
435 j_coord_offsetD = DIM*jnrD;
437 /* load j atom coordinates */
438 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
439 x+j_coord_offsetC,x+j_coord_offsetD,
442 /* Calculate displacement vector */
443 dx00 = _mm_sub_ps(ix0,jx0);
444 dy00 = _mm_sub_ps(iy0,jy0);
445 dz00 = _mm_sub_ps(iz0,jz0);
446 dx10 = _mm_sub_ps(ix1,jx0);
447 dy10 = _mm_sub_ps(iy1,jy0);
448 dz10 = _mm_sub_ps(iz1,jz0);
449 dx20 = _mm_sub_ps(ix2,jx0);
450 dy20 = _mm_sub_ps(iy2,jy0);
451 dz20 = _mm_sub_ps(iz2,jz0);
452 dx30 = _mm_sub_ps(ix3,jx0);
453 dy30 = _mm_sub_ps(iy3,jy0);
454 dz30 = _mm_sub_ps(iz3,jz0);
456 /* Calculate squared distance and things based on it */
457 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
458 rsq10 = gmx_mm_calc_rsq_ps(dx10,dy10,dz10);
459 rsq20 = gmx_mm_calc_rsq_ps(dx20,dy20,dz20);
460 rsq30 = gmx_mm_calc_rsq_ps(dx30,dy30,dz30);
462 rinv10 = gmx_mm_invsqrt_ps(rsq10);
463 rinv20 = gmx_mm_invsqrt_ps(rsq20);
464 rinv30 = gmx_mm_invsqrt_ps(rsq30);
466 rinvsq00 = gmx_mm_inv_ps(rsq00);
467 rinvsq10 = _mm_mul_ps(rinv10,rinv10);
468 rinvsq20 = _mm_mul_ps(rinv20,rinv20);
469 rinvsq30 = _mm_mul_ps(rinv30,rinv30);
471 /* Load parameters for j particles */
472 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
473 charge+jnrC+0,charge+jnrD+0);
474 vdwjidx0A = 2*vdwtype[jnrA+0];
475 vdwjidx0B = 2*vdwtype[jnrB+0];
476 vdwjidx0C = 2*vdwtype[jnrC+0];
477 vdwjidx0D = 2*vdwtype[jnrD+0];
479 fjx0 = _mm_setzero_ps();
480 fjy0 = _mm_setzero_ps();
481 fjz0 = _mm_setzero_ps();
483 /**************************
484 * CALCULATE INTERACTIONS *
485 **************************/
487 if (gmx_mm_any_lt(rsq00,rcutoff2))
490 /* Compute parameters for interactions between i and j atoms */
491 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
492 vdwparam+vdwioffset0+vdwjidx0B,
493 vdwparam+vdwioffset0+vdwjidx0C,
494 vdwparam+vdwioffset0+vdwjidx0D,
497 /* LENNARD-JONES DISPERSION/REPULSION */
499 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
500 vvdw6 = _mm_mul_ps(c6_00,rinvsix);
501 vvdw12 = _mm_mul_ps(c12_00,_mm_mul_ps(rinvsix,rinvsix));
502 vvdw = _mm_msub_ps(_mm_nmacc_ps(c12_00,_mm_mul_ps(sh_vdw_invrcut6,sh_vdw_invrcut6),vvdw12),one_twelfth,
503 _mm_mul_ps( _mm_nmacc_ps(c6_00,sh_vdw_invrcut6,vvdw6),one_sixth));
504 fvdw = _mm_mul_ps(_mm_sub_ps(vvdw12,vvdw6),rinvsq00);
506 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
508 /* Update potential sum for this i atom from the interaction with this j atom. */
509 vvdw = _mm_and_ps(vvdw,cutoff_mask);
510 vvdw = _mm_andnot_ps(dummy_mask,vvdw);
511 vvdwsum = _mm_add_ps(vvdwsum,vvdw);
515 fscal = _mm_and_ps(fscal,cutoff_mask);
517 fscal = _mm_andnot_ps(dummy_mask,fscal);
519 /* Update vectorial force */
520 fix0 = _mm_macc_ps(dx00,fscal,fix0);
521 fiy0 = _mm_macc_ps(dy00,fscal,fiy0);
522 fiz0 = _mm_macc_ps(dz00,fscal,fiz0);
524 fjx0 = _mm_macc_ps(dx00,fscal,fjx0);
525 fjy0 = _mm_macc_ps(dy00,fscal,fjy0);
526 fjz0 = _mm_macc_ps(dz00,fscal,fjz0);
530 /**************************
531 * CALCULATE INTERACTIONS *
532 **************************/
534 if (gmx_mm_any_lt(rsq10,rcutoff2))
537 /* Compute parameters for interactions between i and j atoms */
538 qq10 = _mm_mul_ps(iq1,jq0);
540 /* REACTION-FIELD ELECTROSTATICS */
541 velec = _mm_mul_ps(qq10,_mm_sub_ps(_mm_macc_ps(krf,rsq10,rinv10),crf));
542 felec = _mm_mul_ps(qq10,_mm_msub_ps(rinv10,rinvsq10,krf2));
544 cutoff_mask = _mm_cmplt_ps(rsq10,rcutoff2);
546 /* Update potential sum for this i atom from the interaction with this j atom. */
547 velec = _mm_and_ps(velec,cutoff_mask);
548 velec = _mm_andnot_ps(dummy_mask,velec);
549 velecsum = _mm_add_ps(velecsum,velec);
553 fscal = _mm_and_ps(fscal,cutoff_mask);
555 fscal = _mm_andnot_ps(dummy_mask,fscal);
557 /* Update vectorial force */
558 fix1 = _mm_macc_ps(dx10,fscal,fix1);
559 fiy1 = _mm_macc_ps(dy10,fscal,fiy1);
560 fiz1 = _mm_macc_ps(dz10,fscal,fiz1);
562 fjx0 = _mm_macc_ps(dx10,fscal,fjx0);
563 fjy0 = _mm_macc_ps(dy10,fscal,fjy0);
564 fjz0 = _mm_macc_ps(dz10,fscal,fjz0);
568 /**************************
569 * CALCULATE INTERACTIONS *
570 **************************/
572 if (gmx_mm_any_lt(rsq20,rcutoff2))
575 /* Compute parameters for interactions between i and j atoms */
576 qq20 = _mm_mul_ps(iq2,jq0);
578 /* REACTION-FIELD ELECTROSTATICS */
579 velec = _mm_mul_ps(qq20,_mm_sub_ps(_mm_macc_ps(krf,rsq20,rinv20),crf));
580 felec = _mm_mul_ps(qq20,_mm_msub_ps(rinv20,rinvsq20,krf2));
582 cutoff_mask = _mm_cmplt_ps(rsq20,rcutoff2);
584 /* Update potential sum for this i atom from the interaction with this j atom. */
585 velec = _mm_and_ps(velec,cutoff_mask);
586 velec = _mm_andnot_ps(dummy_mask,velec);
587 velecsum = _mm_add_ps(velecsum,velec);
591 fscal = _mm_and_ps(fscal,cutoff_mask);
593 fscal = _mm_andnot_ps(dummy_mask,fscal);
595 /* Update vectorial force */
596 fix2 = _mm_macc_ps(dx20,fscal,fix2);
597 fiy2 = _mm_macc_ps(dy20,fscal,fiy2);
598 fiz2 = _mm_macc_ps(dz20,fscal,fiz2);
600 fjx0 = _mm_macc_ps(dx20,fscal,fjx0);
601 fjy0 = _mm_macc_ps(dy20,fscal,fjy0);
602 fjz0 = _mm_macc_ps(dz20,fscal,fjz0);
606 /**************************
607 * CALCULATE INTERACTIONS *
608 **************************/
610 if (gmx_mm_any_lt(rsq30,rcutoff2))
613 /* Compute parameters for interactions between i and j atoms */
614 qq30 = _mm_mul_ps(iq3,jq0);
616 /* REACTION-FIELD ELECTROSTATICS */
617 velec = _mm_mul_ps(qq30,_mm_sub_ps(_mm_macc_ps(krf,rsq30,rinv30),crf));
618 felec = _mm_mul_ps(qq30,_mm_msub_ps(rinv30,rinvsq30,krf2));
620 cutoff_mask = _mm_cmplt_ps(rsq30,rcutoff2);
622 /* Update potential sum for this i atom from the interaction with this j atom. */
623 velec = _mm_and_ps(velec,cutoff_mask);
624 velec = _mm_andnot_ps(dummy_mask,velec);
625 velecsum = _mm_add_ps(velecsum,velec);
629 fscal = _mm_and_ps(fscal,cutoff_mask);
631 fscal = _mm_andnot_ps(dummy_mask,fscal);
633 /* Update vectorial force */
634 fix3 = _mm_macc_ps(dx30,fscal,fix3);
635 fiy3 = _mm_macc_ps(dy30,fscal,fiy3);
636 fiz3 = _mm_macc_ps(dz30,fscal,fiz3);
638 fjx0 = _mm_macc_ps(dx30,fscal,fjx0);
639 fjy0 = _mm_macc_ps(dy30,fscal,fjy0);
640 fjz0 = _mm_macc_ps(dz30,fscal,fjz0);
644 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
645 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
646 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
647 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
649 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjx0,fjy0,fjz0);
651 /* Inner loop uses 161 flops */
654 /* End of innermost loop */
656 gmx_mm_update_iforce_4atom_swizzle_ps(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
657 f+i_coord_offset,fshift+i_shift_offset);
660 /* Update potential energies */
661 gmx_mm_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
662 gmx_mm_update_1pot_ps(vvdwsum,kernel_data->energygrp_vdw+ggid);
664 /* Increment number of inner iterations */
665 inneriter += j_index_end - j_index_start;
667 /* Outer loop uses 26 flops */
670 /* Increment number of outer iterations */
673 /* Update outer/inner flops */
675 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W4_VF,outeriter*26 + inneriter*161);
678 * Gromacs nonbonded kernel: nb_kernel_ElecRFCut_VdwLJSh_GeomW4P1_F_avx_128_fma_single
679 * Electrostatics interaction: ReactionField
680 * VdW interaction: LennardJones
681 * Geometry: Water4-Particle
682 * Calculate force/pot: Force
685 nb_kernel_ElecRFCut_VdwLJSh_GeomW4P1_F_avx_128_fma_single
686 (t_nblist * gmx_restrict nlist,
687 rvec * gmx_restrict xx,
688 rvec * gmx_restrict ff,
689 t_forcerec * gmx_restrict fr,
690 t_mdatoms * gmx_restrict mdatoms,
691 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
692 t_nrnb * gmx_restrict nrnb)
694 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
695 * just 0 for non-waters.
696 * Suffixes A,B,C,D refer to j loop unrolling done with AVX_128, e.g. for the four different
697 * jnr indices corresponding to data put in the four positions in the SIMD register.
699 int i_shift_offset,i_coord_offset,outeriter,inneriter;
700 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
701 int jnrA,jnrB,jnrC,jnrD;
702 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
703 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
704 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
706 real *shiftvec,*fshift,*x,*f;
707 real *fjptrA,*fjptrB,*fjptrC,*fjptrD;
709 __m128 fscal,rcutoff,rcutoff2,jidxall;
711 __m128 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
713 __m128 ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
715 __m128 ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
717 __m128 ix3,iy3,iz3,fix3,fiy3,fiz3,iq3,isai3;
718 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
719 __m128 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
720 __m128 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
721 __m128 dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
722 __m128 dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
723 __m128 dx30,dy30,dz30,rsq30,rinv30,rinvsq30,r30,qq30,c6_30,c12_30;
724 __m128 velec,felec,velecsum,facel,crf,krf,krf2;
727 __m128 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
730 __m128 one_sixth = _mm_set1_ps(1.0/6.0);
731 __m128 one_twelfth = _mm_set1_ps(1.0/12.0);
732 __m128 dummy_mask,cutoff_mask;
733 __m128 signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
734 __m128 one = _mm_set1_ps(1.0);
735 __m128 two = _mm_set1_ps(2.0);
741 jindex = nlist->jindex;
743 shiftidx = nlist->shift;
745 shiftvec = fr->shift_vec[0];
746 fshift = fr->fshift[0];
747 facel = _mm_set1_ps(fr->epsfac);
748 charge = mdatoms->chargeA;
749 krf = _mm_set1_ps(fr->ic->k_rf);
750 krf2 = _mm_set1_ps(fr->ic->k_rf*2.0);
751 crf = _mm_set1_ps(fr->ic->c_rf);
752 nvdwtype = fr->ntype;
754 vdwtype = mdatoms->typeA;
756 /* Setup water-specific parameters */
757 inr = nlist->iinr[0];
758 iq1 = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+1]));
759 iq2 = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+2]));
760 iq3 = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+3]));
761 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
763 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
764 rcutoff_scalar = fr->rcoulomb;
765 rcutoff = _mm_set1_ps(rcutoff_scalar);
766 rcutoff2 = _mm_mul_ps(rcutoff,rcutoff);
768 sh_vdw_invrcut6 = _mm_set1_ps(fr->ic->sh_invrc6);
769 rvdw = _mm_set1_ps(fr->rvdw);
771 /* Avoid stupid compiler warnings */
772 jnrA = jnrB = jnrC = jnrD = 0;
781 for(iidx=0;iidx<4*DIM;iidx++)
786 /* Start outer loop over neighborlists */
787 for(iidx=0; iidx<nri; iidx++)
789 /* Load shift vector for this list */
790 i_shift_offset = DIM*shiftidx[iidx];
792 /* Load limits for loop over neighbors */
793 j_index_start = jindex[iidx];
794 j_index_end = jindex[iidx+1];
796 /* Get outer coordinate index */
798 i_coord_offset = DIM*inr;
800 /* Load i particle coords and add shift vector */
801 gmx_mm_load_shift_and_4rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,
802 &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
804 fix0 = _mm_setzero_ps();
805 fiy0 = _mm_setzero_ps();
806 fiz0 = _mm_setzero_ps();
807 fix1 = _mm_setzero_ps();
808 fiy1 = _mm_setzero_ps();
809 fiz1 = _mm_setzero_ps();
810 fix2 = _mm_setzero_ps();
811 fiy2 = _mm_setzero_ps();
812 fiz2 = _mm_setzero_ps();
813 fix3 = _mm_setzero_ps();
814 fiy3 = _mm_setzero_ps();
815 fiz3 = _mm_setzero_ps();
817 /* Start inner kernel loop */
818 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
821 /* Get j neighbor index, and coordinate index */
826 j_coord_offsetA = DIM*jnrA;
827 j_coord_offsetB = DIM*jnrB;
828 j_coord_offsetC = DIM*jnrC;
829 j_coord_offsetD = DIM*jnrD;
831 /* load j atom coordinates */
832 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
833 x+j_coord_offsetC,x+j_coord_offsetD,
836 /* Calculate displacement vector */
837 dx00 = _mm_sub_ps(ix0,jx0);
838 dy00 = _mm_sub_ps(iy0,jy0);
839 dz00 = _mm_sub_ps(iz0,jz0);
840 dx10 = _mm_sub_ps(ix1,jx0);
841 dy10 = _mm_sub_ps(iy1,jy0);
842 dz10 = _mm_sub_ps(iz1,jz0);
843 dx20 = _mm_sub_ps(ix2,jx0);
844 dy20 = _mm_sub_ps(iy2,jy0);
845 dz20 = _mm_sub_ps(iz2,jz0);
846 dx30 = _mm_sub_ps(ix3,jx0);
847 dy30 = _mm_sub_ps(iy3,jy0);
848 dz30 = _mm_sub_ps(iz3,jz0);
850 /* Calculate squared distance and things based on it */
851 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
852 rsq10 = gmx_mm_calc_rsq_ps(dx10,dy10,dz10);
853 rsq20 = gmx_mm_calc_rsq_ps(dx20,dy20,dz20);
854 rsq30 = gmx_mm_calc_rsq_ps(dx30,dy30,dz30);
856 rinv10 = gmx_mm_invsqrt_ps(rsq10);
857 rinv20 = gmx_mm_invsqrt_ps(rsq20);
858 rinv30 = gmx_mm_invsqrt_ps(rsq30);
860 rinvsq00 = gmx_mm_inv_ps(rsq00);
861 rinvsq10 = _mm_mul_ps(rinv10,rinv10);
862 rinvsq20 = _mm_mul_ps(rinv20,rinv20);
863 rinvsq30 = _mm_mul_ps(rinv30,rinv30);
865 /* Load parameters for j particles */
866 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
867 charge+jnrC+0,charge+jnrD+0);
868 vdwjidx0A = 2*vdwtype[jnrA+0];
869 vdwjidx0B = 2*vdwtype[jnrB+0];
870 vdwjidx0C = 2*vdwtype[jnrC+0];
871 vdwjidx0D = 2*vdwtype[jnrD+0];
873 fjx0 = _mm_setzero_ps();
874 fjy0 = _mm_setzero_ps();
875 fjz0 = _mm_setzero_ps();
877 /**************************
878 * CALCULATE INTERACTIONS *
879 **************************/
881 if (gmx_mm_any_lt(rsq00,rcutoff2))
884 /* Compute parameters for interactions between i and j atoms */
885 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
886 vdwparam+vdwioffset0+vdwjidx0B,
887 vdwparam+vdwioffset0+vdwjidx0C,
888 vdwparam+vdwioffset0+vdwjidx0D,
891 /* LENNARD-JONES DISPERSION/REPULSION */
893 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
894 fvdw = _mm_mul_ps(_mm_msub_ps(c12_00,rinvsix,c6_00),_mm_mul_ps(rinvsix,rinvsq00));
896 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
900 fscal = _mm_and_ps(fscal,cutoff_mask);
902 /* Update vectorial force */
903 fix0 = _mm_macc_ps(dx00,fscal,fix0);
904 fiy0 = _mm_macc_ps(dy00,fscal,fiy0);
905 fiz0 = _mm_macc_ps(dz00,fscal,fiz0);
907 fjx0 = _mm_macc_ps(dx00,fscal,fjx0);
908 fjy0 = _mm_macc_ps(dy00,fscal,fjy0);
909 fjz0 = _mm_macc_ps(dz00,fscal,fjz0);
913 /**************************
914 * CALCULATE INTERACTIONS *
915 **************************/
917 if (gmx_mm_any_lt(rsq10,rcutoff2))
920 /* Compute parameters for interactions between i and j atoms */
921 qq10 = _mm_mul_ps(iq1,jq0);
923 /* REACTION-FIELD ELECTROSTATICS */
924 felec = _mm_mul_ps(qq10,_mm_msub_ps(rinv10,rinvsq10,krf2));
926 cutoff_mask = _mm_cmplt_ps(rsq10,rcutoff2);
930 fscal = _mm_and_ps(fscal,cutoff_mask);
932 /* Update vectorial force */
933 fix1 = _mm_macc_ps(dx10,fscal,fix1);
934 fiy1 = _mm_macc_ps(dy10,fscal,fiy1);
935 fiz1 = _mm_macc_ps(dz10,fscal,fiz1);
937 fjx0 = _mm_macc_ps(dx10,fscal,fjx0);
938 fjy0 = _mm_macc_ps(dy10,fscal,fjy0);
939 fjz0 = _mm_macc_ps(dz10,fscal,fjz0);
943 /**************************
944 * CALCULATE INTERACTIONS *
945 **************************/
947 if (gmx_mm_any_lt(rsq20,rcutoff2))
950 /* Compute parameters for interactions between i and j atoms */
951 qq20 = _mm_mul_ps(iq2,jq0);
953 /* REACTION-FIELD ELECTROSTATICS */
954 felec = _mm_mul_ps(qq20,_mm_msub_ps(rinv20,rinvsq20,krf2));
956 cutoff_mask = _mm_cmplt_ps(rsq20,rcutoff2);
960 fscal = _mm_and_ps(fscal,cutoff_mask);
962 /* Update vectorial force */
963 fix2 = _mm_macc_ps(dx20,fscal,fix2);
964 fiy2 = _mm_macc_ps(dy20,fscal,fiy2);
965 fiz2 = _mm_macc_ps(dz20,fscal,fiz2);
967 fjx0 = _mm_macc_ps(dx20,fscal,fjx0);
968 fjy0 = _mm_macc_ps(dy20,fscal,fjy0);
969 fjz0 = _mm_macc_ps(dz20,fscal,fjz0);
973 /**************************
974 * CALCULATE INTERACTIONS *
975 **************************/
977 if (gmx_mm_any_lt(rsq30,rcutoff2))
980 /* Compute parameters for interactions between i and j atoms */
981 qq30 = _mm_mul_ps(iq3,jq0);
983 /* REACTION-FIELD ELECTROSTATICS */
984 felec = _mm_mul_ps(qq30,_mm_msub_ps(rinv30,rinvsq30,krf2));
986 cutoff_mask = _mm_cmplt_ps(rsq30,rcutoff2);
990 fscal = _mm_and_ps(fscal,cutoff_mask);
992 /* Update vectorial force */
993 fix3 = _mm_macc_ps(dx30,fscal,fix3);
994 fiy3 = _mm_macc_ps(dy30,fscal,fiy3);
995 fiz3 = _mm_macc_ps(dz30,fscal,fiz3);
997 fjx0 = _mm_macc_ps(dx30,fscal,fjx0);
998 fjy0 = _mm_macc_ps(dy30,fscal,fjy0);
999 fjz0 = _mm_macc_ps(dz30,fscal,fjz0);
1003 fjptrA = f+j_coord_offsetA;
1004 fjptrB = f+j_coord_offsetB;
1005 fjptrC = f+j_coord_offsetC;
1006 fjptrD = f+j_coord_offsetD;
1008 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjx0,fjy0,fjz0);
1010 /* Inner loop uses 132 flops */
1013 if(jidx<j_index_end)
1016 /* Get j neighbor index, and coordinate index */
1017 jnrlistA = jjnr[jidx];
1018 jnrlistB = jjnr[jidx+1];
1019 jnrlistC = jjnr[jidx+2];
1020 jnrlistD = jjnr[jidx+3];
1021 /* Sign of each element will be negative for non-real atoms.
1022 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
1023 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
1025 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
1026 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
1027 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
1028 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
1029 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
1030 j_coord_offsetA = DIM*jnrA;
1031 j_coord_offsetB = DIM*jnrB;
1032 j_coord_offsetC = DIM*jnrC;
1033 j_coord_offsetD = DIM*jnrD;
1035 /* load j atom coordinates */
1036 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
1037 x+j_coord_offsetC,x+j_coord_offsetD,
1040 /* Calculate displacement vector */
1041 dx00 = _mm_sub_ps(ix0,jx0);
1042 dy00 = _mm_sub_ps(iy0,jy0);
1043 dz00 = _mm_sub_ps(iz0,jz0);
1044 dx10 = _mm_sub_ps(ix1,jx0);
1045 dy10 = _mm_sub_ps(iy1,jy0);
1046 dz10 = _mm_sub_ps(iz1,jz0);
1047 dx20 = _mm_sub_ps(ix2,jx0);
1048 dy20 = _mm_sub_ps(iy2,jy0);
1049 dz20 = _mm_sub_ps(iz2,jz0);
1050 dx30 = _mm_sub_ps(ix3,jx0);
1051 dy30 = _mm_sub_ps(iy3,jy0);
1052 dz30 = _mm_sub_ps(iz3,jz0);
1054 /* Calculate squared distance and things based on it */
1055 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
1056 rsq10 = gmx_mm_calc_rsq_ps(dx10,dy10,dz10);
1057 rsq20 = gmx_mm_calc_rsq_ps(dx20,dy20,dz20);
1058 rsq30 = gmx_mm_calc_rsq_ps(dx30,dy30,dz30);
1060 rinv10 = gmx_mm_invsqrt_ps(rsq10);
1061 rinv20 = gmx_mm_invsqrt_ps(rsq20);
1062 rinv30 = gmx_mm_invsqrt_ps(rsq30);
1064 rinvsq00 = gmx_mm_inv_ps(rsq00);
1065 rinvsq10 = _mm_mul_ps(rinv10,rinv10);
1066 rinvsq20 = _mm_mul_ps(rinv20,rinv20);
1067 rinvsq30 = _mm_mul_ps(rinv30,rinv30);
1069 /* Load parameters for j particles */
1070 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
1071 charge+jnrC+0,charge+jnrD+0);
1072 vdwjidx0A = 2*vdwtype[jnrA+0];
1073 vdwjidx0B = 2*vdwtype[jnrB+0];
1074 vdwjidx0C = 2*vdwtype[jnrC+0];
1075 vdwjidx0D = 2*vdwtype[jnrD+0];
1077 fjx0 = _mm_setzero_ps();
1078 fjy0 = _mm_setzero_ps();
1079 fjz0 = _mm_setzero_ps();
1081 /**************************
1082 * CALCULATE INTERACTIONS *
1083 **************************/
1085 if (gmx_mm_any_lt(rsq00,rcutoff2))
1088 /* Compute parameters for interactions between i and j atoms */
1089 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
1090 vdwparam+vdwioffset0+vdwjidx0B,
1091 vdwparam+vdwioffset0+vdwjidx0C,
1092 vdwparam+vdwioffset0+vdwjidx0D,
1095 /* LENNARD-JONES DISPERSION/REPULSION */
1097 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
1098 fvdw = _mm_mul_ps(_mm_msub_ps(c12_00,rinvsix,c6_00),_mm_mul_ps(rinvsix,rinvsq00));
1100 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
1104 fscal = _mm_and_ps(fscal,cutoff_mask);
1106 fscal = _mm_andnot_ps(dummy_mask,fscal);
1108 /* Update vectorial force */
1109 fix0 = _mm_macc_ps(dx00,fscal,fix0);
1110 fiy0 = _mm_macc_ps(dy00,fscal,fiy0);
1111 fiz0 = _mm_macc_ps(dz00,fscal,fiz0);
1113 fjx0 = _mm_macc_ps(dx00,fscal,fjx0);
1114 fjy0 = _mm_macc_ps(dy00,fscal,fjy0);
1115 fjz0 = _mm_macc_ps(dz00,fscal,fjz0);
1119 /**************************
1120 * CALCULATE INTERACTIONS *
1121 **************************/
1123 if (gmx_mm_any_lt(rsq10,rcutoff2))
1126 /* Compute parameters for interactions between i and j atoms */
1127 qq10 = _mm_mul_ps(iq1,jq0);
1129 /* REACTION-FIELD ELECTROSTATICS */
1130 felec = _mm_mul_ps(qq10,_mm_msub_ps(rinv10,rinvsq10,krf2));
1132 cutoff_mask = _mm_cmplt_ps(rsq10,rcutoff2);
1136 fscal = _mm_and_ps(fscal,cutoff_mask);
1138 fscal = _mm_andnot_ps(dummy_mask,fscal);
1140 /* Update vectorial force */
1141 fix1 = _mm_macc_ps(dx10,fscal,fix1);
1142 fiy1 = _mm_macc_ps(dy10,fscal,fiy1);
1143 fiz1 = _mm_macc_ps(dz10,fscal,fiz1);
1145 fjx0 = _mm_macc_ps(dx10,fscal,fjx0);
1146 fjy0 = _mm_macc_ps(dy10,fscal,fjy0);
1147 fjz0 = _mm_macc_ps(dz10,fscal,fjz0);
1151 /**************************
1152 * CALCULATE INTERACTIONS *
1153 **************************/
1155 if (gmx_mm_any_lt(rsq20,rcutoff2))
1158 /* Compute parameters for interactions between i and j atoms */
1159 qq20 = _mm_mul_ps(iq2,jq0);
1161 /* REACTION-FIELD ELECTROSTATICS */
1162 felec = _mm_mul_ps(qq20,_mm_msub_ps(rinv20,rinvsq20,krf2));
1164 cutoff_mask = _mm_cmplt_ps(rsq20,rcutoff2);
1168 fscal = _mm_and_ps(fscal,cutoff_mask);
1170 fscal = _mm_andnot_ps(dummy_mask,fscal);
1172 /* Update vectorial force */
1173 fix2 = _mm_macc_ps(dx20,fscal,fix2);
1174 fiy2 = _mm_macc_ps(dy20,fscal,fiy2);
1175 fiz2 = _mm_macc_ps(dz20,fscal,fiz2);
1177 fjx0 = _mm_macc_ps(dx20,fscal,fjx0);
1178 fjy0 = _mm_macc_ps(dy20,fscal,fjy0);
1179 fjz0 = _mm_macc_ps(dz20,fscal,fjz0);
1183 /**************************
1184 * CALCULATE INTERACTIONS *
1185 **************************/
1187 if (gmx_mm_any_lt(rsq30,rcutoff2))
1190 /* Compute parameters for interactions between i and j atoms */
1191 qq30 = _mm_mul_ps(iq3,jq0);
1193 /* REACTION-FIELD ELECTROSTATICS */
1194 felec = _mm_mul_ps(qq30,_mm_msub_ps(rinv30,rinvsq30,krf2));
1196 cutoff_mask = _mm_cmplt_ps(rsq30,rcutoff2);
1200 fscal = _mm_and_ps(fscal,cutoff_mask);
1202 fscal = _mm_andnot_ps(dummy_mask,fscal);
1204 /* Update vectorial force */
1205 fix3 = _mm_macc_ps(dx30,fscal,fix3);
1206 fiy3 = _mm_macc_ps(dy30,fscal,fiy3);
1207 fiz3 = _mm_macc_ps(dz30,fscal,fiz3);
1209 fjx0 = _mm_macc_ps(dx30,fscal,fjx0);
1210 fjy0 = _mm_macc_ps(dy30,fscal,fjy0);
1211 fjz0 = _mm_macc_ps(dz30,fscal,fjz0);
1215 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
1216 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
1217 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
1218 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
1220 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjx0,fjy0,fjz0);
1222 /* Inner loop uses 132 flops */
1225 /* End of innermost loop */
1227 gmx_mm_update_iforce_4atom_swizzle_ps(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
1228 f+i_coord_offset,fshift+i_shift_offset);
1230 /* Increment number of inner iterations */
1231 inneriter += j_index_end - j_index_start;
1233 /* Outer loop uses 24 flops */
1236 /* Increment number of outer iterations */
1239 /* Update outer/inner flops */
1241 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W4_F,outeriter*24 + inneriter*132);