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_256_single kernel generator.
44 #include "../nb_kernel.h"
45 #include "types/simple.h"
49 #include "gmx_math_x86_avx_256_single.h"
50 #include "kernelutil_x86_avx_256_single.h"
53 * Gromacs nonbonded kernel: nb_kernel_ElecRF_VdwLJ_GeomW4W4_VF_avx_256_single
54 * Electrostatics interaction: ReactionField
55 * VdW interaction: LennardJones
56 * Geometry: Water4-Water4
57 * Calculate force/pot: PotentialAndForce
60 nb_kernel_ElecRF_VdwLJ_GeomW4W4_VF_avx_256_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,E,F,G,H refer to j loop unrolling done with AVX, e.g. for the eight 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 jnrE,jnrF,jnrG,jnrH;
78 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
79 int jnrlistE,jnrlistF,jnrlistG,jnrlistH;
80 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
81 int j_coord_offsetE,j_coord_offsetF,j_coord_offsetG,j_coord_offsetH;
82 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
84 real *shiftvec,*fshift,*x,*f;
85 real *fjptrA,*fjptrB,*fjptrC,*fjptrD,*fjptrE,*fjptrF,*fjptrG,*fjptrH;
87 __m256 tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
88 real * vdwioffsetptr0;
89 __m256 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
90 real * vdwioffsetptr1;
91 __m256 ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
92 real * vdwioffsetptr2;
93 __m256 ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
94 real * vdwioffsetptr3;
95 __m256 ix3,iy3,iz3,fix3,fiy3,fiz3,iq3,isai3;
96 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D,vdwjidx0E,vdwjidx0F,vdwjidx0G,vdwjidx0H;
97 __m256 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
98 int vdwjidx1A,vdwjidx1B,vdwjidx1C,vdwjidx1D,vdwjidx1E,vdwjidx1F,vdwjidx1G,vdwjidx1H;
99 __m256 jx1,jy1,jz1,fjx1,fjy1,fjz1,jq1,isaj1;
100 int vdwjidx2A,vdwjidx2B,vdwjidx2C,vdwjidx2D,vdwjidx2E,vdwjidx2F,vdwjidx2G,vdwjidx2H;
101 __m256 jx2,jy2,jz2,fjx2,fjy2,fjz2,jq2,isaj2;
102 int vdwjidx3A,vdwjidx3B,vdwjidx3C,vdwjidx3D,vdwjidx3E,vdwjidx3F,vdwjidx3G,vdwjidx3H;
103 __m256 jx3,jy3,jz3,fjx3,fjy3,fjz3,jq3,isaj3;
104 __m256 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
105 __m256 dx11,dy11,dz11,rsq11,rinv11,rinvsq11,r11,qq11,c6_11,c12_11;
106 __m256 dx12,dy12,dz12,rsq12,rinv12,rinvsq12,r12,qq12,c6_12,c12_12;
107 __m256 dx13,dy13,dz13,rsq13,rinv13,rinvsq13,r13,qq13,c6_13,c12_13;
108 __m256 dx21,dy21,dz21,rsq21,rinv21,rinvsq21,r21,qq21,c6_21,c12_21;
109 __m256 dx22,dy22,dz22,rsq22,rinv22,rinvsq22,r22,qq22,c6_22,c12_22;
110 __m256 dx23,dy23,dz23,rsq23,rinv23,rinvsq23,r23,qq23,c6_23,c12_23;
111 __m256 dx31,dy31,dz31,rsq31,rinv31,rinvsq31,r31,qq31,c6_31,c12_31;
112 __m256 dx32,dy32,dz32,rsq32,rinv32,rinvsq32,r32,qq32,c6_32,c12_32;
113 __m256 dx33,dy33,dz33,rsq33,rinv33,rinvsq33,r33,qq33,c6_33,c12_33;
114 __m256 velec,felec,velecsum,facel,crf,krf,krf2;
117 __m256 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
120 __m256 one_sixth = _mm256_set1_ps(1.0/6.0);
121 __m256 one_twelfth = _mm256_set1_ps(1.0/12.0);
122 __m256 dummy_mask,cutoff_mask;
123 __m256 signbit = _mm256_castsi256_ps( _mm256_set1_epi32(0x80000000) );
124 __m256 one = _mm256_set1_ps(1.0);
125 __m256 two = _mm256_set1_ps(2.0);
131 jindex = nlist->jindex;
133 shiftidx = nlist->shift;
135 shiftvec = fr->shift_vec[0];
136 fshift = fr->fshift[0];
137 facel = _mm256_set1_ps(fr->epsfac);
138 charge = mdatoms->chargeA;
139 krf = _mm256_set1_ps(fr->ic->k_rf);
140 krf2 = _mm256_set1_ps(fr->ic->k_rf*2.0);
141 crf = _mm256_set1_ps(fr->ic->c_rf);
142 nvdwtype = fr->ntype;
144 vdwtype = mdatoms->typeA;
146 /* Setup water-specific parameters */
147 inr = nlist->iinr[0];
148 iq1 = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+1]));
149 iq2 = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+2]));
150 iq3 = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+3]));
151 vdwioffsetptr0 = vdwparam+2*nvdwtype*vdwtype[inr+0];
153 jq1 = _mm256_set1_ps(charge[inr+1]);
154 jq2 = _mm256_set1_ps(charge[inr+2]);
155 jq3 = _mm256_set1_ps(charge[inr+3]);
156 vdwjidx0A = 2*vdwtype[inr+0];
157 c6_00 = _mm256_set1_ps(vdwioffsetptr0[vdwjidx0A]);
158 c12_00 = _mm256_set1_ps(vdwioffsetptr0[vdwjidx0A+1]);
159 qq11 = _mm256_mul_ps(iq1,jq1);
160 qq12 = _mm256_mul_ps(iq1,jq2);
161 qq13 = _mm256_mul_ps(iq1,jq3);
162 qq21 = _mm256_mul_ps(iq2,jq1);
163 qq22 = _mm256_mul_ps(iq2,jq2);
164 qq23 = _mm256_mul_ps(iq2,jq3);
165 qq31 = _mm256_mul_ps(iq3,jq1);
166 qq32 = _mm256_mul_ps(iq3,jq2);
167 qq33 = _mm256_mul_ps(iq3,jq3);
169 /* Avoid stupid compiler warnings */
170 jnrA = jnrB = jnrC = jnrD = jnrE = jnrF = jnrG = jnrH = 0;
183 for(iidx=0;iidx<4*DIM;iidx++)
188 /* Start outer loop over neighborlists */
189 for(iidx=0; iidx<nri; iidx++)
191 /* Load shift vector for this list */
192 i_shift_offset = DIM*shiftidx[iidx];
194 /* Load limits for loop over neighbors */
195 j_index_start = jindex[iidx];
196 j_index_end = jindex[iidx+1];
198 /* Get outer coordinate index */
200 i_coord_offset = DIM*inr;
202 /* Load i particle coords and add shift vector */
203 gmx_mm256_load_shift_and_4rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,
204 &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
206 fix0 = _mm256_setzero_ps();
207 fiy0 = _mm256_setzero_ps();
208 fiz0 = _mm256_setzero_ps();
209 fix1 = _mm256_setzero_ps();
210 fiy1 = _mm256_setzero_ps();
211 fiz1 = _mm256_setzero_ps();
212 fix2 = _mm256_setzero_ps();
213 fiy2 = _mm256_setzero_ps();
214 fiz2 = _mm256_setzero_ps();
215 fix3 = _mm256_setzero_ps();
216 fiy3 = _mm256_setzero_ps();
217 fiz3 = _mm256_setzero_ps();
219 /* Reset potential sums */
220 velecsum = _mm256_setzero_ps();
221 vvdwsum = _mm256_setzero_ps();
223 /* Start inner kernel loop */
224 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+7]>=0; jidx+=8)
227 /* Get j neighbor index, and coordinate index */
236 j_coord_offsetA = DIM*jnrA;
237 j_coord_offsetB = DIM*jnrB;
238 j_coord_offsetC = DIM*jnrC;
239 j_coord_offsetD = DIM*jnrD;
240 j_coord_offsetE = DIM*jnrE;
241 j_coord_offsetF = DIM*jnrF;
242 j_coord_offsetG = DIM*jnrG;
243 j_coord_offsetH = DIM*jnrH;
245 /* load j atom coordinates */
246 gmx_mm256_load_4rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
247 x+j_coord_offsetC,x+j_coord_offsetD,
248 x+j_coord_offsetE,x+j_coord_offsetF,
249 x+j_coord_offsetG,x+j_coord_offsetH,
250 &jx0,&jy0,&jz0,&jx1,&jy1,&jz1,&jx2,
251 &jy2,&jz2,&jx3,&jy3,&jz3);
253 /* Calculate displacement vector */
254 dx00 = _mm256_sub_ps(ix0,jx0);
255 dy00 = _mm256_sub_ps(iy0,jy0);
256 dz00 = _mm256_sub_ps(iz0,jz0);
257 dx11 = _mm256_sub_ps(ix1,jx1);
258 dy11 = _mm256_sub_ps(iy1,jy1);
259 dz11 = _mm256_sub_ps(iz1,jz1);
260 dx12 = _mm256_sub_ps(ix1,jx2);
261 dy12 = _mm256_sub_ps(iy1,jy2);
262 dz12 = _mm256_sub_ps(iz1,jz2);
263 dx13 = _mm256_sub_ps(ix1,jx3);
264 dy13 = _mm256_sub_ps(iy1,jy3);
265 dz13 = _mm256_sub_ps(iz1,jz3);
266 dx21 = _mm256_sub_ps(ix2,jx1);
267 dy21 = _mm256_sub_ps(iy2,jy1);
268 dz21 = _mm256_sub_ps(iz2,jz1);
269 dx22 = _mm256_sub_ps(ix2,jx2);
270 dy22 = _mm256_sub_ps(iy2,jy2);
271 dz22 = _mm256_sub_ps(iz2,jz2);
272 dx23 = _mm256_sub_ps(ix2,jx3);
273 dy23 = _mm256_sub_ps(iy2,jy3);
274 dz23 = _mm256_sub_ps(iz2,jz3);
275 dx31 = _mm256_sub_ps(ix3,jx1);
276 dy31 = _mm256_sub_ps(iy3,jy1);
277 dz31 = _mm256_sub_ps(iz3,jz1);
278 dx32 = _mm256_sub_ps(ix3,jx2);
279 dy32 = _mm256_sub_ps(iy3,jy2);
280 dz32 = _mm256_sub_ps(iz3,jz2);
281 dx33 = _mm256_sub_ps(ix3,jx3);
282 dy33 = _mm256_sub_ps(iy3,jy3);
283 dz33 = _mm256_sub_ps(iz3,jz3);
285 /* Calculate squared distance and things based on it */
286 rsq00 = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
287 rsq11 = gmx_mm256_calc_rsq_ps(dx11,dy11,dz11);
288 rsq12 = gmx_mm256_calc_rsq_ps(dx12,dy12,dz12);
289 rsq13 = gmx_mm256_calc_rsq_ps(dx13,dy13,dz13);
290 rsq21 = gmx_mm256_calc_rsq_ps(dx21,dy21,dz21);
291 rsq22 = gmx_mm256_calc_rsq_ps(dx22,dy22,dz22);
292 rsq23 = gmx_mm256_calc_rsq_ps(dx23,dy23,dz23);
293 rsq31 = gmx_mm256_calc_rsq_ps(dx31,dy31,dz31);
294 rsq32 = gmx_mm256_calc_rsq_ps(dx32,dy32,dz32);
295 rsq33 = gmx_mm256_calc_rsq_ps(dx33,dy33,dz33);
297 rinv11 = gmx_mm256_invsqrt_ps(rsq11);
298 rinv12 = gmx_mm256_invsqrt_ps(rsq12);
299 rinv13 = gmx_mm256_invsqrt_ps(rsq13);
300 rinv21 = gmx_mm256_invsqrt_ps(rsq21);
301 rinv22 = gmx_mm256_invsqrt_ps(rsq22);
302 rinv23 = gmx_mm256_invsqrt_ps(rsq23);
303 rinv31 = gmx_mm256_invsqrt_ps(rsq31);
304 rinv32 = gmx_mm256_invsqrt_ps(rsq32);
305 rinv33 = gmx_mm256_invsqrt_ps(rsq33);
307 rinvsq00 = gmx_mm256_inv_ps(rsq00);
308 rinvsq11 = _mm256_mul_ps(rinv11,rinv11);
309 rinvsq12 = _mm256_mul_ps(rinv12,rinv12);
310 rinvsq13 = _mm256_mul_ps(rinv13,rinv13);
311 rinvsq21 = _mm256_mul_ps(rinv21,rinv21);
312 rinvsq22 = _mm256_mul_ps(rinv22,rinv22);
313 rinvsq23 = _mm256_mul_ps(rinv23,rinv23);
314 rinvsq31 = _mm256_mul_ps(rinv31,rinv31);
315 rinvsq32 = _mm256_mul_ps(rinv32,rinv32);
316 rinvsq33 = _mm256_mul_ps(rinv33,rinv33);
318 fjx0 = _mm256_setzero_ps();
319 fjy0 = _mm256_setzero_ps();
320 fjz0 = _mm256_setzero_ps();
321 fjx1 = _mm256_setzero_ps();
322 fjy1 = _mm256_setzero_ps();
323 fjz1 = _mm256_setzero_ps();
324 fjx2 = _mm256_setzero_ps();
325 fjy2 = _mm256_setzero_ps();
326 fjz2 = _mm256_setzero_ps();
327 fjx3 = _mm256_setzero_ps();
328 fjy3 = _mm256_setzero_ps();
329 fjz3 = _mm256_setzero_ps();
331 /**************************
332 * CALCULATE INTERACTIONS *
333 **************************/
335 /* LENNARD-JONES DISPERSION/REPULSION */
337 rinvsix = _mm256_mul_ps(_mm256_mul_ps(rinvsq00,rinvsq00),rinvsq00);
338 vvdw6 = _mm256_mul_ps(c6_00,rinvsix);
339 vvdw12 = _mm256_mul_ps(c12_00,_mm256_mul_ps(rinvsix,rinvsix));
340 vvdw = _mm256_sub_ps( _mm256_mul_ps(vvdw12,one_twelfth) , _mm256_mul_ps(vvdw6,one_sixth) );
341 fvdw = _mm256_mul_ps(_mm256_sub_ps(vvdw12,vvdw6),rinvsq00);
343 /* Update potential sum for this i atom from the interaction with this j atom. */
344 vvdwsum = _mm256_add_ps(vvdwsum,vvdw);
348 /* Calculate temporary vectorial force */
349 tx = _mm256_mul_ps(fscal,dx00);
350 ty = _mm256_mul_ps(fscal,dy00);
351 tz = _mm256_mul_ps(fscal,dz00);
353 /* Update vectorial force */
354 fix0 = _mm256_add_ps(fix0,tx);
355 fiy0 = _mm256_add_ps(fiy0,ty);
356 fiz0 = _mm256_add_ps(fiz0,tz);
358 fjx0 = _mm256_add_ps(fjx0,tx);
359 fjy0 = _mm256_add_ps(fjy0,ty);
360 fjz0 = _mm256_add_ps(fjz0,tz);
362 /**************************
363 * CALCULATE INTERACTIONS *
364 **************************/
366 /* REACTION-FIELD ELECTROSTATICS */
367 velec = _mm256_mul_ps(qq11,_mm256_sub_ps(_mm256_add_ps(rinv11,_mm256_mul_ps(krf,rsq11)),crf));
368 felec = _mm256_mul_ps(qq11,_mm256_sub_ps(_mm256_mul_ps(rinv11,rinvsq11),krf2));
370 /* Update potential sum for this i atom from the interaction with this j atom. */
371 velecsum = _mm256_add_ps(velecsum,velec);
375 /* Calculate temporary vectorial force */
376 tx = _mm256_mul_ps(fscal,dx11);
377 ty = _mm256_mul_ps(fscal,dy11);
378 tz = _mm256_mul_ps(fscal,dz11);
380 /* Update vectorial force */
381 fix1 = _mm256_add_ps(fix1,tx);
382 fiy1 = _mm256_add_ps(fiy1,ty);
383 fiz1 = _mm256_add_ps(fiz1,tz);
385 fjx1 = _mm256_add_ps(fjx1,tx);
386 fjy1 = _mm256_add_ps(fjy1,ty);
387 fjz1 = _mm256_add_ps(fjz1,tz);
389 /**************************
390 * CALCULATE INTERACTIONS *
391 **************************/
393 /* REACTION-FIELD ELECTROSTATICS */
394 velec = _mm256_mul_ps(qq12,_mm256_sub_ps(_mm256_add_ps(rinv12,_mm256_mul_ps(krf,rsq12)),crf));
395 felec = _mm256_mul_ps(qq12,_mm256_sub_ps(_mm256_mul_ps(rinv12,rinvsq12),krf2));
397 /* Update potential sum for this i atom from the interaction with this j atom. */
398 velecsum = _mm256_add_ps(velecsum,velec);
402 /* Calculate temporary vectorial force */
403 tx = _mm256_mul_ps(fscal,dx12);
404 ty = _mm256_mul_ps(fscal,dy12);
405 tz = _mm256_mul_ps(fscal,dz12);
407 /* Update vectorial force */
408 fix1 = _mm256_add_ps(fix1,tx);
409 fiy1 = _mm256_add_ps(fiy1,ty);
410 fiz1 = _mm256_add_ps(fiz1,tz);
412 fjx2 = _mm256_add_ps(fjx2,tx);
413 fjy2 = _mm256_add_ps(fjy2,ty);
414 fjz2 = _mm256_add_ps(fjz2,tz);
416 /**************************
417 * CALCULATE INTERACTIONS *
418 **************************/
420 /* REACTION-FIELD ELECTROSTATICS */
421 velec = _mm256_mul_ps(qq13,_mm256_sub_ps(_mm256_add_ps(rinv13,_mm256_mul_ps(krf,rsq13)),crf));
422 felec = _mm256_mul_ps(qq13,_mm256_sub_ps(_mm256_mul_ps(rinv13,rinvsq13),krf2));
424 /* Update potential sum for this i atom from the interaction with this j atom. */
425 velecsum = _mm256_add_ps(velecsum,velec);
429 /* Calculate temporary vectorial force */
430 tx = _mm256_mul_ps(fscal,dx13);
431 ty = _mm256_mul_ps(fscal,dy13);
432 tz = _mm256_mul_ps(fscal,dz13);
434 /* Update vectorial force */
435 fix1 = _mm256_add_ps(fix1,tx);
436 fiy1 = _mm256_add_ps(fiy1,ty);
437 fiz1 = _mm256_add_ps(fiz1,tz);
439 fjx3 = _mm256_add_ps(fjx3,tx);
440 fjy3 = _mm256_add_ps(fjy3,ty);
441 fjz3 = _mm256_add_ps(fjz3,tz);
443 /**************************
444 * CALCULATE INTERACTIONS *
445 **************************/
447 /* REACTION-FIELD ELECTROSTATICS */
448 velec = _mm256_mul_ps(qq21,_mm256_sub_ps(_mm256_add_ps(rinv21,_mm256_mul_ps(krf,rsq21)),crf));
449 felec = _mm256_mul_ps(qq21,_mm256_sub_ps(_mm256_mul_ps(rinv21,rinvsq21),krf2));
451 /* Update potential sum for this i atom from the interaction with this j atom. */
452 velecsum = _mm256_add_ps(velecsum,velec);
456 /* Calculate temporary vectorial force */
457 tx = _mm256_mul_ps(fscal,dx21);
458 ty = _mm256_mul_ps(fscal,dy21);
459 tz = _mm256_mul_ps(fscal,dz21);
461 /* Update vectorial force */
462 fix2 = _mm256_add_ps(fix2,tx);
463 fiy2 = _mm256_add_ps(fiy2,ty);
464 fiz2 = _mm256_add_ps(fiz2,tz);
466 fjx1 = _mm256_add_ps(fjx1,tx);
467 fjy1 = _mm256_add_ps(fjy1,ty);
468 fjz1 = _mm256_add_ps(fjz1,tz);
470 /**************************
471 * CALCULATE INTERACTIONS *
472 **************************/
474 /* REACTION-FIELD ELECTROSTATICS */
475 velec = _mm256_mul_ps(qq22,_mm256_sub_ps(_mm256_add_ps(rinv22,_mm256_mul_ps(krf,rsq22)),crf));
476 felec = _mm256_mul_ps(qq22,_mm256_sub_ps(_mm256_mul_ps(rinv22,rinvsq22),krf2));
478 /* Update potential sum for this i atom from the interaction with this j atom. */
479 velecsum = _mm256_add_ps(velecsum,velec);
483 /* Calculate temporary vectorial force */
484 tx = _mm256_mul_ps(fscal,dx22);
485 ty = _mm256_mul_ps(fscal,dy22);
486 tz = _mm256_mul_ps(fscal,dz22);
488 /* Update vectorial force */
489 fix2 = _mm256_add_ps(fix2,tx);
490 fiy2 = _mm256_add_ps(fiy2,ty);
491 fiz2 = _mm256_add_ps(fiz2,tz);
493 fjx2 = _mm256_add_ps(fjx2,tx);
494 fjy2 = _mm256_add_ps(fjy2,ty);
495 fjz2 = _mm256_add_ps(fjz2,tz);
497 /**************************
498 * CALCULATE INTERACTIONS *
499 **************************/
501 /* REACTION-FIELD ELECTROSTATICS */
502 velec = _mm256_mul_ps(qq23,_mm256_sub_ps(_mm256_add_ps(rinv23,_mm256_mul_ps(krf,rsq23)),crf));
503 felec = _mm256_mul_ps(qq23,_mm256_sub_ps(_mm256_mul_ps(rinv23,rinvsq23),krf2));
505 /* Update potential sum for this i atom from the interaction with this j atom. */
506 velecsum = _mm256_add_ps(velecsum,velec);
510 /* Calculate temporary vectorial force */
511 tx = _mm256_mul_ps(fscal,dx23);
512 ty = _mm256_mul_ps(fscal,dy23);
513 tz = _mm256_mul_ps(fscal,dz23);
515 /* Update vectorial force */
516 fix2 = _mm256_add_ps(fix2,tx);
517 fiy2 = _mm256_add_ps(fiy2,ty);
518 fiz2 = _mm256_add_ps(fiz2,tz);
520 fjx3 = _mm256_add_ps(fjx3,tx);
521 fjy3 = _mm256_add_ps(fjy3,ty);
522 fjz3 = _mm256_add_ps(fjz3,tz);
524 /**************************
525 * CALCULATE INTERACTIONS *
526 **************************/
528 /* REACTION-FIELD ELECTROSTATICS */
529 velec = _mm256_mul_ps(qq31,_mm256_sub_ps(_mm256_add_ps(rinv31,_mm256_mul_ps(krf,rsq31)),crf));
530 felec = _mm256_mul_ps(qq31,_mm256_sub_ps(_mm256_mul_ps(rinv31,rinvsq31),krf2));
532 /* Update potential sum for this i atom from the interaction with this j atom. */
533 velecsum = _mm256_add_ps(velecsum,velec);
537 /* Calculate temporary vectorial force */
538 tx = _mm256_mul_ps(fscal,dx31);
539 ty = _mm256_mul_ps(fscal,dy31);
540 tz = _mm256_mul_ps(fscal,dz31);
542 /* Update vectorial force */
543 fix3 = _mm256_add_ps(fix3,tx);
544 fiy3 = _mm256_add_ps(fiy3,ty);
545 fiz3 = _mm256_add_ps(fiz3,tz);
547 fjx1 = _mm256_add_ps(fjx1,tx);
548 fjy1 = _mm256_add_ps(fjy1,ty);
549 fjz1 = _mm256_add_ps(fjz1,tz);
551 /**************************
552 * CALCULATE INTERACTIONS *
553 **************************/
555 /* REACTION-FIELD ELECTROSTATICS */
556 velec = _mm256_mul_ps(qq32,_mm256_sub_ps(_mm256_add_ps(rinv32,_mm256_mul_ps(krf,rsq32)),crf));
557 felec = _mm256_mul_ps(qq32,_mm256_sub_ps(_mm256_mul_ps(rinv32,rinvsq32),krf2));
559 /* Update potential sum for this i atom from the interaction with this j atom. */
560 velecsum = _mm256_add_ps(velecsum,velec);
564 /* Calculate temporary vectorial force */
565 tx = _mm256_mul_ps(fscal,dx32);
566 ty = _mm256_mul_ps(fscal,dy32);
567 tz = _mm256_mul_ps(fscal,dz32);
569 /* Update vectorial force */
570 fix3 = _mm256_add_ps(fix3,tx);
571 fiy3 = _mm256_add_ps(fiy3,ty);
572 fiz3 = _mm256_add_ps(fiz3,tz);
574 fjx2 = _mm256_add_ps(fjx2,tx);
575 fjy2 = _mm256_add_ps(fjy2,ty);
576 fjz2 = _mm256_add_ps(fjz2,tz);
578 /**************************
579 * CALCULATE INTERACTIONS *
580 **************************/
582 /* REACTION-FIELD ELECTROSTATICS */
583 velec = _mm256_mul_ps(qq33,_mm256_sub_ps(_mm256_add_ps(rinv33,_mm256_mul_ps(krf,rsq33)),crf));
584 felec = _mm256_mul_ps(qq33,_mm256_sub_ps(_mm256_mul_ps(rinv33,rinvsq33),krf2));
586 /* Update potential sum for this i atom from the interaction with this j atom. */
587 velecsum = _mm256_add_ps(velecsum,velec);
591 /* Calculate temporary vectorial force */
592 tx = _mm256_mul_ps(fscal,dx33);
593 ty = _mm256_mul_ps(fscal,dy33);
594 tz = _mm256_mul_ps(fscal,dz33);
596 /* Update vectorial force */
597 fix3 = _mm256_add_ps(fix3,tx);
598 fiy3 = _mm256_add_ps(fiy3,ty);
599 fiz3 = _mm256_add_ps(fiz3,tz);
601 fjx3 = _mm256_add_ps(fjx3,tx);
602 fjy3 = _mm256_add_ps(fjy3,ty);
603 fjz3 = _mm256_add_ps(fjz3,tz);
605 fjptrA = f+j_coord_offsetA;
606 fjptrB = f+j_coord_offsetB;
607 fjptrC = f+j_coord_offsetC;
608 fjptrD = f+j_coord_offsetD;
609 fjptrE = f+j_coord_offsetE;
610 fjptrF = f+j_coord_offsetF;
611 fjptrG = f+j_coord_offsetG;
612 fjptrH = f+j_coord_offsetH;
614 gmx_mm256_decrement_4rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,
615 fjx0,fjy0,fjz0,fjx1,fjy1,fjz1,
616 fjx2,fjy2,fjz2,fjx3,fjy3,fjz3);
618 /* Inner loop uses 323 flops */
624 /* Get j neighbor index, and coordinate index */
625 jnrlistA = jjnr[jidx];
626 jnrlistB = jjnr[jidx+1];
627 jnrlistC = jjnr[jidx+2];
628 jnrlistD = jjnr[jidx+3];
629 jnrlistE = jjnr[jidx+4];
630 jnrlistF = jjnr[jidx+5];
631 jnrlistG = jjnr[jidx+6];
632 jnrlistH = jjnr[jidx+7];
633 /* Sign of each element will be negative for non-real atoms.
634 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
635 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
637 dummy_mask = gmx_mm256_set_m128(gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx+4)),_mm_setzero_si128())),
638 gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128())));
640 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
641 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
642 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
643 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
644 jnrE = (jnrlistE>=0) ? jnrlistE : 0;
645 jnrF = (jnrlistF>=0) ? jnrlistF : 0;
646 jnrG = (jnrlistG>=0) ? jnrlistG : 0;
647 jnrH = (jnrlistH>=0) ? jnrlistH : 0;
648 j_coord_offsetA = DIM*jnrA;
649 j_coord_offsetB = DIM*jnrB;
650 j_coord_offsetC = DIM*jnrC;
651 j_coord_offsetD = DIM*jnrD;
652 j_coord_offsetE = DIM*jnrE;
653 j_coord_offsetF = DIM*jnrF;
654 j_coord_offsetG = DIM*jnrG;
655 j_coord_offsetH = DIM*jnrH;
657 /* load j atom coordinates */
658 gmx_mm256_load_4rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
659 x+j_coord_offsetC,x+j_coord_offsetD,
660 x+j_coord_offsetE,x+j_coord_offsetF,
661 x+j_coord_offsetG,x+j_coord_offsetH,
662 &jx0,&jy0,&jz0,&jx1,&jy1,&jz1,&jx2,
663 &jy2,&jz2,&jx3,&jy3,&jz3);
665 /* Calculate displacement vector */
666 dx00 = _mm256_sub_ps(ix0,jx0);
667 dy00 = _mm256_sub_ps(iy0,jy0);
668 dz00 = _mm256_sub_ps(iz0,jz0);
669 dx11 = _mm256_sub_ps(ix1,jx1);
670 dy11 = _mm256_sub_ps(iy1,jy1);
671 dz11 = _mm256_sub_ps(iz1,jz1);
672 dx12 = _mm256_sub_ps(ix1,jx2);
673 dy12 = _mm256_sub_ps(iy1,jy2);
674 dz12 = _mm256_sub_ps(iz1,jz2);
675 dx13 = _mm256_sub_ps(ix1,jx3);
676 dy13 = _mm256_sub_ps(iy1,jy3);
677 dz13 = _mm256_sub_ps(iz1,jz3);
678 dx21 = _mm256_sub_ps(ix2,jx1);
679 dy21 = _mm256_sub_ps(iy2,jy1);
680 dz21 = _mm256_sub_ps(iz2,jz1);
681 dx22 = _mm256_sub_ps(ix2,jx2);
682 dy22 = _mm256_sub_ps(iy2,jy2);
683 dz22 = _mm256_sub_ps(iz2,jz2);
684 dx23 = _mm256_sub_ps(ix2,jx3);
685 dy23 = _mm256_sub_ps(iy2,jy3);
686 dz23 = _mm256_sub_ps(iz2,jz3);
687 dx31 = _mm256_sub_ps(ix3,jx1);
688 dy31 = _mm256_sub_ps(iy3,jy1);
689 dz31 = _mm256_sub_ps(iz3,jz1);
690 dx32 = _mm256_sub_ps(ix3,jx2);
691 dy32 = _mm256_sub_ps(iy3,jy2);
692 dz32 = _mm256_sub_ps(iz3,jz2);
693 dx33 = _mm256_sub_ps(ix3,jx3);
694 dy33 = _mm256_sub_ps(iy3,jy3);
695 dz33 = _mm256_sub_ps(iz3,jz3);
697 /* Calculate squared distance and things based on it */
698 rsq00 = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
699 rsq11 = gmx_mm256_calc_rsq_ps(dx11,dy11,dz11);
700 rsq12 = gmx_mm256_calc_rsq_ps(dx12,dy12,dz12);
701 rsq13 = gmx_mm256_calc_rsq_ps(dx13,dy13,dz13);
702 rsq21 = gmx_mm256_calc_rsq_ps(dx21,dy21,dz21);
703 rsq22 = gmx_mm256_calc_rsq_ps(dx22,dy22,dz22);
704 rsq23 = gmx_mm256_calc_rsq_ps(dx23,dy23,dz23);
705 rsq31 = gmx_mm256_calc_rsq_ps(dx31,dy31,dz31);
706 rsq32 = gmx_mm256_calc_rsq_ps(dx32,dy32,dz32);
707 rsq33 = gmx_mm256_calc_rsq_ps(dx33,dy33,dz33);
709 rinv11 = gmx_mm256_invsqrt_ps(rsq11);
710 rinv12 = gmx_mm256_invsqrt_ps(rsq12);
711 rinv13 = gmx_mm256_invsqrt_ps(rsq13);
712 rinv21 = gmx_mm256_invsqrt_ps(rsq21);
713 rinv22 = gmx_mm256_invsqrt_ps(rsq22);
714 rinv23 = gmx_mm256_invsqrt_ps(rsq23);
715 rinv31 = gmx_mm256_invsqrt_ps(rsq31);
716 rinv32 = gmx_mm256_invsqrt_ps(rsq32);
717 rinv33 = gmx_mm256_invsqrt_ps(rsq33);
719 rinvsq00 = gmx_mm256_inv_ps(rsq00);
720 rinvsq11 = _mm256_mul_ps(rinv11,rinv11);
721 rinvsq12 = _mm256_mul_ps(rinv12,rinv12);
722 rinvsq13 = _mm256_mul_ps(rinv13,rinv13);
723 rinvsq21 = _mm256_mul_ps(rinv21,rinv21);
724 rinvsq22 = _mm256_mul_ps(rinv22,rinv22);
725 rinvsq23 = _mm256_mul_ps(rinv23,rinv23);
726 rinvsq31 = _mm256_mul_ps(rinv31,rinv31);
727 rinvsq32 = _mm256_mul_ps(rinv32,rinv32);
728 rinvsq33 = _mm256_mul_ps(rinv33,rinv33);
730 fjx0 = _mm256_setzero_ps();
731 fjy0 = _mm256_setzero_ps();
732 fjz0 = _mm256_setzero_ps();
733 fjx1 = _mm256_setzero_ps();
734 fjy1 = _mm256_setzero_ps();
735 fjz1 = _mm256_setzero_ps();
736 fjx2 = _mm256_setzero_ps();
737 fjy2 = _mm256_setzero_ps();
738 fjz2 = _mm256_setzero_ps();
739 fjx3 = _mm256_setzero_ps();
740 fjy3 = _mm256_setzero_ps();
741 fjz3 = _mm256_setzero_ps();
743 /**************************
744 * CALCULATE INTERACTIONS *
745 **************************/
747 /* LENNARD-JONES DISPERSION/REPULSION */
749 rinvsix = _mm256_mul_ps(_mm256_mul_ps(rinvsq00,rinvsq00),rinvsq00);
750 vvdw6 = _mm256_mul_ps(c6_00,rinvsix);
751 vvdw12 = _mm256_mul_ps(c12_00,_mm256_mul_ps(rinvsix,rinvsix));
752 vvdw = _mm256_sub_ps( _mm256_mul_ps(vvdw12,one_twelfth) , _mm256_mul_ps(vvdw6,one_sixth) );
753 fvdw = _mm256_mul_ps(_mm256_sub_ps(vvdw12,vvdw6),rinvsq00);
755 /* Update potential sum for this i atom from the interaction with this j atom. */
756 vvdw = _mm256_andnot_ps(dummy_mask,vvdw);
757 vvdwsum = _mm256_add_ps(vvdwsum,vvdw);
761 fscal = _mm256_andnot_ps(dummy_mask,fscal);
763 /* Calculate temporary vectorial force */
764 tx = _mm256_mul_ps(fscal,dx00);
765 ty = _mm256_mul_ps(fscal,dy00);
766 tz = _mm256_mul_ps(fscal,dz00);
768 /* Update vectorial force */
769 fix0 = _mm256_add_ps(fix0,tx);
770 fiy0 = _mm256_add_ps(fiy0,ty);
771 fiz0 = _mm256_add_ps(fiz0,tz);
773 fjx0 = _mm256_add_ps(fjx0,tx);
774 fjy0 = _mm256_add_ps(fjy0,ty);
775 fjz0 = _mm256_add_ps(fjz0,tz);
777 /**************************
778 * CALCULATE INTERACTIONS *
779 **************************/
781 /* REACTION-FIELD ELECTROSTATICS */
782 velec = _mm256_mul_ps(qq11,_mm256_sub_ps(_mm256_add_ps(rinv11,_mm256_mul_ps(krf,rsq11)),crf));
783 felec = _mm256_mul_ps(qq11,_mm256_sub_ps(_mm256_mul_ps(rinv11,rinvsq11),krf2));
785 /* Update potential sum for this i atom from the interaction with this j atom. */
786 velec = _mm256_andnot_ps(dummy_mask,velec);
787 velecsum = _mm256_add_ps(velecsum,velec);
791 fscal = _mm256_andnot_ps(dummy_mask,fscal);
793 /* Calculate temporary vectorial force */
794 tx = _mm256_mul_ps(fscal,dx11);
795 ty = _mm256_mul_ps(fscal,dy11);
796 tz = _mm256_mul_ps(fscal,dz11);
798 /* Update vectorial force */
799 fix1 = _mm256_add_ps(fix1,tx);
800 fiy1 = _mm256_add_ps(fiy1,ty);
801 fiz1 = _mm256_add_ps(fiz1,tz);
803 fjx1 = _mm256_add_ps(fjx1,tx);
804 fjy1 = _mm256_add_ps(fjy1,ty);
805 fjz1 = _mm256_add_ps(fjz1,tz);
807 /**************************
808 * CALCULATE INTERACTIONS *
809 **************************/
811 /* REACTION-FIELD ELECTROSTATICS */
812 velec = _mm256_mul_ps(qq12,_mm256_sub_ps(_mm256_add_ps(rinv12,_mm256_mul_ps(krf,rsq12)),crf));
813 felec = _mm256_mul_ps(qq12,_mm256_sub_ps(_mm256_mul_ps(rinv12,rinvsq12),krf2));
815 /* Update potential sum for this i atom from the interaction with this j atom. */
816 velec = _mm256_andnot_ps(dummy_mask,velec);
817 velecsum = _mm256_add_ps(velecsum,velec);
821 fscal = _mm256_andnot_ps(dummy_mask,fscal);
823 /* Calculate temporary vectorial force */
824 tx = _mm256_mul_ps(fscal,dx12);
825 ty = _mm256_mul_ps(fscal,dy12);
826 tz = _mm256_mul_ps(fscal,dz12);
828 /* Update vectorial force */
829 fix1 = _mm256_add_ps(fix1,tx);
830 fiy1 = _mm256_add_ps(fiy1,ty);
831 fiz1 = _mm256_add_ps(fiz1,tz);
833 fjx2 = _mm256_add_ps(fjx2,tx);
834 fjy2 = _mm256_add_ps(fjy2,ty);
835 fjz2 = _mm256_add_ps(fjz2,tz);
837 /**************************
838 * CALCULATE INTERACTIONS *
839 **************************/
841 /* REACTION-FIELD ELECTROSTATICS */
842 velec = _mm256_mul_ps(qq13,_mm256_sub_ps(_mm256_add_ps(rinv13,_mm256_mul_ps(krf,rsq13)),crf));
843 felec = _mm256_mul_ps(qq13,_mm256_sub_ps(_mm256_mul_ps(rinv13,rinvsq13),krf2));
845 /* Update potential sum for this i atom from the interaction with this j atom. */
846 velec = _mm256_andnot_ps(dummy_mask,velec);
847 velecsum = _mm256_add_ps(velecsum,velec);
851 fscal = _mm256_andnot_ps(dummy_mask,fscal);
853 /* Calculate temporary vectorial force */
854 tx = _mm256_mul_ps(fscal,dx13);
855 ty = _mm256_mul_ps(fscal,dy13);
856 tz = _mm256_mul_ps(fscal,dz13);
858 /* Update vectorial force */
859 fix1 = _mm256_add_ps(fix1,tx);
860 fiy1 = _mm256_add_ps(fiy1,ty);
861 fiz1 = _mm256_add_ps(fiz1,tz);
863 fjx3 = _mm256_add_ps(fjx3,tx);
864 fjy3 = _mm256_add_ps(fjy3,ty);
865 fjz3 = _mm256_add_ps(fjz3,tz);
867 /**************************
868 * CALCULATE INTERACTIONS *
869 **************************/
871 /* REACTION-FIELD ELECTROSTATICS */
872 velec = _mm256_mul_ps(qq21,_mm256_sub_ps(_mm256_add_ps(rinv21,_mm256_mul_ps(krf,rsq21)),crf));
873 felec = _mm256_mul_ps(qq21,_mm256_sub_ps(_mm256_mul_ps(rinv21,rinvsq21),krf2));
875 /* Update potential sum for this i atom from the interaction with this j atom. */
876 velec = _mm256_andnot_ps(dummy_mask,velec);
877 velecsum = _mm256_add_ps(velecsum,velec);
881 fscal = _mm256_andnot_ps(dummy_mask,fscal);
883 /* Calculate temporary vectorial force */
884 tx = _mm256_mul_ps(fscal,dx21);
885 ty = _mm256_mul_ps(fscal,dy21);
886 tz = _mm256_mul_ps(fscal,dz21);
888 /* Update vectorial force */
889 fix2 = _mm256_add_ps(fix2,tx);
890 fiy2 = _mm256_add_ps(fiy2,ty);
891 fiz2 = _mm256_add_ps(fiz2,tz);
893 fjx1 = _mm256_add_ps(fjx1,tx);
894 fjy1 = _mm256_add_ps(fjy1,ty);
895 fjz1 = _mm256_add_ps(fjz1,tz);
897 /**************************
898 * CALCULATE INTERACTIONS *
899 **************************/
901 /* REACTION-FIELD ELECTROSTATICS */
902 velec = _mm256_mul_ps(qq22,_mm256_sub_ps(_mm256_add_ps(rinv22,_mm256_mul_ps(krf,rsq22)),crf));
903 felec = _mm256_mul_ps(qq22,_mm256_sub_ps(_mm256_mul_ps(rinv22,rinvsq22),krf2));
905 /* Update potential sum for this i atom from the interaction with this j atom. */
906 velec = _mm256_andnot_ps(dummy_mask,velec);
907 velecsum = _mm256_add_ps(velecsum,velec);
911 fscal = _mm256_andnot_ps(dummy_mask,fscal);
913 /* Calculate temporary vectorial force */
914 tx = _mm256_mul_ps(fscal,dx22);
915 ty = _mm256_mul_ps(fscal,dy22);
916 tz = _mm256_mul_ps(fscal,dz22);
918 /* Update vectorial force */
919 fix2 = _mm256_add_ps(fix2,tx);
920 fiy2 = _mm256_add_ps(fiy2,ty);
921 fiz2 = _mm256_add_ps(fiz2,tz);
923 fjx2 = _mm256_add_ps(fjx2,tx);
924 fjy2 = _mm256_add_ps(fjy2,ty);
925 fjz2 = _mm256_add_ps(fjz2,tz);
927 /**************************
928 * CALCULATE INTERACTIONS *
929 **************************/
931 /* REACTION-FIELD ELECTROSTATICS */
932 velec = _mm256_mul_ps(qq23,_mm256_sub_ps(_mm256_add_ps(rinv23,_mm256_mul_ps(krf,rsq23)),crf));
933 felec = _mm256_mul_ps(qq23,_mm256_sub_ps(_mm256_mul_ps(rinv23,rinvsq23),krf2));
935 /* Update potential sum for this i atom from the interaction with this j atom. */
936 velec = _mm256_andnot_ps(dummy_mask,velec);
937 velecsum = _mm256_add_ps(velecsum,velec);
941 fscal = _mm256_andnot_ps(dummy_mask,fscal);
943 /* Calculate temporary vectorial force */
944 tx = _mm256_mul_ps(fscal,dx23);
945 ty = _mm256_mul_ps(fscal,dy23);
946 tz = _mm256_mul_ps(fscal,dz23);
948 /* Update vectorial force */
949 fix2 = _mm256_add_ps(fix2,tx);
950 fiy2 = _mm256_add_ps(fiy2,ty);
951 fiz2 = _mm256_add_ps(fiz2,tz);
953 fjx3 = _mm256_add_ps(fjx3,tx);
954 fjy3 = _mm256_add_ps(fjy3,ty);
955 fjz3 = _mm256_add_ps(fjz3,tz);
957 /**************************
958 * CALCULATE INTERACTIONS *
959 **************************/
961 /* REACTION-FIELD ELECTROSTATICS */
962 velec = _mm256_mul_ps(qq31,_mm256_sub_ps(_mm256_add_ps(rinv31,_mm256_mul_ps(krf,rsq31)),crf));
963 felec = _mm256_mul_ps(qq31,_mm256_sub_ps(_mm256_mul_ps(rinv31,rinvsq31),krf2));
965 /* Update potential sum for this i atom from the interaction with this j atom. */
966 velec = _mm256_andnot_ps(dummy_mask,velec);
967 velecsum = _mm256_add_ps(velecsum,velec);
971 fscal = _mm256_andnot_ps(dummy_mask,fscal);
973 /* Calculate temporary vectorial force */
974 tx = _mm256_mul_ps(fscal,dx31);
975 ty = _mm256_mul_ps(fscal,dy31);
976 tz = _mm256_mul_ps(fscal,dz31);
978 /* Update vectorial force */
979 fix3 = _mm256_add_ps(fix3,tx);
980 fiy3 = _mm256_add_ps(fiy3,ty);
981 fiz3 = _mm256_add_ps(fiz3,tz);
983 fjx1 = _mm256_add_ps(fjx1,tx);
984 fjy1 = _mm256_add_ps(fjy1,ty);
985 fjz1 = _mm256_add_ps(fjz1,tz);
987 /**************************
988 * CALCULATE INTERACTIONS *
989 **************************/
991 /* REACTION-FIELD ELECTROSTATICS */
992 velec = _mm256_mul_ps(qq32,_mm256_sub_ps(_mm256_add_ps(rinv32,_mm256_mul_ps(krf,rsq32)),crf));
993 felec = _mm256_mul_ps(qq32,_mm256_sub_ps(_mm256_mul_ps(rinv32,rinvsq32),krf2));
995 /* Update potential sum for this i atom from the interaction with this j atom. */
996 velec = _mm256_andnot_ps(dummy_mask,velec);
997 velecsum = _mm256_add_ps(velecsum,velec);
1001 fscal = _mm256_andnot_ps(dummy_mask,fscal);
1003 /* Calculate temporary vectorial force */
1004 tx = _mm256_mul_ps(fscal,dx32);
1005 ty = _mm256_mul_ps(fscal,dy32);
1006 tz = _mm256_mul_ps(fscal,dz32);
1008 /* Update vectorial force */
1009 fix3 = _mm256_add_ps(fix3,tx);
1010 fiy3 = _mm256_add_ps(fiy3,ty);
1011 fiz3 = _mm256_add_ps(fiz3,tz);
1013 fjx2 = _mm256_add_ps(fjx2,tx);
1014 fjy2 = _mm256_add_ps(fjy2,ty);
1015 fjz2 = _mm256_add_ps(fjz2,tz);
1017 /**************************
1018 * CALCULATE INTERACTIONS *
1019 **************************/
1021 /* REACTION-FIELD ELECTROSTATICS */
1022 velec = _mm256_mul_ps(qq33,_mm256_sub_ps(_mm256_add_ps(rinv33,_mm256_mul_ps(krf,rsq33)),crf));
1023 felec = _mm256_mul_ps(qq33,_mm256_sub_ps(_mm256_mul_ps(rinv33,rinvsq33),krf2));
1025 /* Update potential sum for this i atom from the interaction with this j atom. */
1026 velec = _mm256_andnot_ps(dummy_mask,velec);
1027 velecsum = _mm256_add_ps(velecsum,velec);
1031 fscal = _mm256_andnot_ps(dummy_mask,fscal);
1033 /* Calculate temporary vectorial force */
1034 tx = _mm256_mul_ps(fscal,dx33);
1035 ty = _mm256_mul_ps(fscal,dy33);
1036 tz = _mm256_mul_ps(fscal,dz33);
1038 /* Update vectorial force */
1039 fix3 = _mm256_add_ps(fix3,tx);
1040 fiy3 = _mm256_add_ps(fiy3,ty);
1041 fiz3 = _mm256_add_ps(fiz3,tz);
1043 fjx3 = _mm256_add_ps(fjx3,tx);
1044 fjy3 = _mm256_add_ps(fjy3,ty);
1045 fjz3 = _mm256_add_ps(fjz3,tz);
1047 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
1048 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
1049 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
1050 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
1051 fjptrE = (jnrlistE>=0) ? f+j_coord_offsetE : scratch;
1052 fjptrF = (jnrlistF>=0) ? f+j_coord_offsetF : scratch;
1053 fjptrG = (jnrlistG>=0) ? f+j_coord_offsetG : scratch;
1054 fjptrH = (jnrlistH>=0) ? f+j_coord_offsetH : scratch;
1056 gmx_mm256_decrement_4rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,
1057 fjx0,fjy0,fjz0,fjx1,fjy1,fjz1,
1058 fjx2,fjy2,fjz2,fjx3,fjy3,fjz3);
1060 /* Inner loop uses 323 flops */
1063 /* End of innermost loop */
1065 gmx_mm256_update_iforce_4atom_swizzle_ps(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
1066 f+i_coord_offset,fshift+i_shift_offset);
1069 /* Update potential energies */
1070 gmx_mm256_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
1071 gmx_mm256_update_1pot_ps(vvdwsum,kernel_data->energygrp_vdw+ggid);
1073 /* Increment number of inner iterations */
1074 inneriter += j_index_end - j_index_start;
1076 /* Outer loop uses 26 flops */
1079 /* Increment number of outer iterations */
1082 /* Update outer/inner flops */
1084 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W4W4_VF,outeriter*26 + inneriter*323);
1087 * Gromacs nonbonded kernel: nb_kernel_ElecRF_VdwLJ_GeomW4W4_F_avx_256_single
1088 * Electrostatics interaction: ReactionField
1089 * VdW interaction: LennardJones
1090 * Geometry: Water4-Water4
1091 * Calculate force/pot: Force
1094 nb_kernel_ElecRF_VdwLJ_GeomW4W4_F_avx_256_single
1095 (t_nblist * gmx_restrict nlist,
1096 rvec * gmx_restrict xx,
1097 rvec * gmx_restrict ff,
1098 t_forcerec * gmx_restrict fr,
1099 t_mdatoms * gmx_restrict mdatoms,
1100 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
1101 t_nrnb * gmx_restrict nrnb)
1103 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
1104 * just 0 for non-waters.
1105 * Suffixes A,B,C,D,E,F,G,H refer to j loop unrolling done with AVX, e.g. for the eight different
1106 * jnr indices corresponding to data put in the four positions in the SIMD register.
1108 int i_shift_offset,i_coord_offset,outeriter,inneriter;
1109 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
1110 int jnrA,jnrB,jnrC,jnrD;
1111 int jnrE,jnrF,jnrG,jnrH;
1112 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
1113 int jnrlistE,jnrlistF,jnrlistG,jnrlistH;
1114 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
1115 int j_coord_offsetE,j_coord_offsetF,j_coord_offsetG,j_coord_offsetH;
1116 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
1117 real rcutoff_scalar;
1118 real *shiftvec,*fshift,*x,*f;
1119 real *fjptrA,*fjptrB,*fjptrC,*fjptrD,*fjptrE,*fjptrF,*fjptrG,*fjptrH;
1120 real scratch[4*DIM];
1121 __m256 tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
1122 real * vdwioffsetptr0;
1123 __m256 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
1124 real * vdwioffsetptr1;
1125 __m256 ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
1126 real * vdwioffsetptr2;
1127 __m256 ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
1128 real * vdwioffsetptr3;
1129 __m256 ix3,iy3,iz3,fix3,fiy3,fiz3,iq3,isai3;
1130 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D,vdwjidx0E,vdwjidx0F,vdwjidx0G,vdwjidx0H;
1131 __m256 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
1132 int vdwjidx1A,vdwjidx1B,vdwjidx1C,vdwjidx1D,vdwjidx1E,vdwjidx1F,vdwjidx1G,vdwjidx1H;
1133 __m256 jx1,jy1,jz1,fjx1,fjy1,fjz1,jq1,isaj1;
1134 int vdwjidx2A,vdwjidx2B,vdwjidx2C,vdwjidx2D,vdwjidx2E,vdwjidx2F,vdwjidx2G,vdwjidx2H;
1135 __m256 jx2,jy2,jz2,fjx2,fjy2,fjz2,jq2,isaj2;
1136 int vdwjidx3A,vdwjidx3B,vdwjidx3C,vdwjidx3D,vdwjidx3E,vdwjidx3F,vdwjidx3G,vdwjidx3H;
1137 __m256 jx3,jy3,jz3,fjx3,fjy3,fjz3,jq3,isaj3;
1138 __m256 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
1139 __m256 dx11,dy11,dz11,rsq11,rinv11,rinvsq11,r11,qq11,c6_11,c12_11;
1140 __m256 dx12,dy12,dz12,rsq12,rinv12,rinvsq12,r12,qq12,c6_12,c12_12;
1141 __m256 dx13,dy13,dz13,rsq13,rinv13,rinvsq13,r13,qq13,c6_13,c12_13;
1142 __m256 dx21,dy21,dz21,rsq21,rinv21,rinvsq21,r21,qq21,c6_21,c12_21;
1143 __m256 dx22,dy22,dz22,rsq22,rinv22,rinvsq22,r22,qq22,c6_22,c12_22;
1144 __m256 dx23,dy23,dz23,rsq23,rinv23,rinvsq23,r23,qq23,c6_23,c12_23;
1145 __m256 dx31,dy31,dz31,rsq31,rinv31,rinvsq31,r31,qq31,c6_31,c12_31;
1146 __m256 dx32,dy32,dz32,rsq32,rinv32,rinvsq32,r32,qq32,c6_32,c12_32;
1147 __m256 dx33,dy33,dz33,rsq33,rinv33,rinvsq33,r33,qq33,c6_33,c12_33;
1148 __m256 velec,felec,velecsum,facel,crf,krf,krf2;
1151 __m256 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
1154 __m256 one_sixth = _mm256_set1_ps(1.0/6.0);
1155 __m256 one_twelfth = _mm256_set1_ps(1.0/12.0);
1156 __m256 dummy_mask,cutoff_mask;
1157 __m256 signbit = _mm256_castsi256_ps( _mm256_set1_epi32(0x80000000) );
1158 __m256 one = _mm256_set1_ps(1.0);
1159 __m256 two = _mm256_set1_ps(2.0);
1165 jindex = nlist->jindex;
1167 shiftidx = nlist->shift;
1169 shiftvec = fr->shift_vec[0];
1170 fshift = fr->fshift[0];
1171 facel = _mm256_set1_ps(fr->epsfac);
1172 charge = mdatoms->chargeA;
1173 krf = _mm256_set1_ps(fr->ic->k_rf);
1174 krf2 = _mm256_set1_ps(fr->ic->k_rf*2.0);
1175 crf = _mm256_set1_ps(fr->ic->c_rf);
1176 nvdwtype = fr->ntype;
1177 vdwparam = fr->nbfp;
1178 vdwtype = mdatoms->typeA;
1180 /* Setup water-specific parameters */
1181 inr = nlist->iinr[0];
1182 iq1 = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+1]));
1183 iq2 = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+2]));
1184 iq3 = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+3]));
1185 vdwioffsetptr0 = vdwparam+2*nvdwtype*vdwtype[inr+0];
1187 jq1 = _mm256_set1_ps(charge[inr+1]);
1188 jq2 = _mm256_set1_ps(charge[inr+2]);
1189 jq3 = _mm256_set1_ps(charge[inr+3]);
1190 vdwjidx0A = 2*vdwtype[inr+0];
1191 c6_00 = _mm256_set1_ps(vdwioffsetptr0[vdwjidx0A]);
1192 c12_00 = _mm256_set1_ps(vdwioffsetptr0[vdwjidx0A+1]);
1193 qq11 = _mm256_mul_ps(iq1,jq1);
1194 qq12 = _mm256_mul_ps(iq1,jq2);
1195 qq13 = _mm256_mul_ps(iq1,jq3);
1196 qq21 = _mm256_mul_ps(iq2,jq1);
1197 qq22 = _mm256_mul_ps(iq2,jq2);
1198 qq23 = _mm256_mul_ps(iq2,jq3);
1199 qq31 = _mm256_mul_ps(iq3,jq1);
1200 qq32 = _mm256_mul_ps(iq3,jq2);
1201 qq33 = _mm256_mul_ps(iq3,jq3);
1203 /* Avoid stupid compiler warnings */
1204 jnrA = jnrB = jnrC = jnrD = jnrE = jnrF = jnrG = jnrH = 0;
1205 j_coord_offsetA = 0;
1206 j_coord_offsetB = 0;
1207 j_coord_offsetC = 0;
1208 j_coord_offsetD = 0;
1209 j_coord_offsetE = 0;
1210 j_coord_offsetF = 0;
1211 j_coord_offsetG = 0;
1212 j_coord_offsetH = 0;
1217 for(iidx=0;iidx<4*DIM;iidx++)
1219 scratch[iidx] = 0.0;
1222 /* Start outer loop over neighborlists */
1223 for(iidx=0; iidx<nri; iidx++)
1225 /* Load shift vector for this list */
1226 i_shift_offset = DIM*shiftidx[iidx];
1228 /* Load limits for loop over neighbors */
1229 j_index_start = jindex[iidx];
1230 j_index_end = jindex[iidx+1];
1232 /* Get outer coordinate index */
1234 i_coord_offset = DIM*inr;
1236 /* Load i particle coords and add shift vector */
1237 gmx_mm256_load_shift_and_4rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,
1238 &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
1240 fix0 = _mm256_setzero_ps();
1241 fiy0 = _mm256_setzero_ps();
1242 fiz0 = _mm256_setzero_ps();
1243 fix1 = _mm256_setzero_ps();
1244 fiy1 = _mm256_setzero_ps();
1245 fiz1 = _mm256_setzero_ps();
1246 fix2 = _mm256_setzero_ps();
1247 fiy2 = _mm256_setzero_ps();
1248 fiz2 = _mm256_setzero_ps();
1249 fix3 = _mm256_setzero_ps();
1250 fiy3 = _mm256_setzero_ps();
1251 fiz3 = _mm256_setzero_ps();
1253 /* Start inner kernel loop */
1254 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+7]>=0; jidx+=8)
1257 /* Get j neighbor index, and coordinate index */
1259 jnrB = jjnr[jidx+1];
1260 jnrC = jjnr[jidx+2];
1261 jnrD = jjnr[jidx+3];
1262 jnrE = jjnr[jidx+4];
1263 jnrF = jjnr[jidx+5];
1264 jnrG = jjnr[jidx+6];
1265 jnrH = jjnr[jidx+7];
1266 j_coord_offsetA = DIM*jnrA;
1267 j_coord_offsetB = DIM*jnrB;
1268 j_coord_offsetC = DIM*jnrC;
1269 j_coord_offsetD = DIM*jnrD;
1270 j_coord_offsetE = DIM*jnrE;
1271 j_coord_offsetF = DIM*jnrF;
1272 j_coord_offsetG = DIM*jnrG;
1273 j_coord_offsetH = DIM*jnrH;
1275 /* load j atom coordinates */
1276 gmx_mm256_load_4rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
1277 x+j_coord_offsetC,x+j_coord_offsetD,
1278 x+j_coord_offsetE,x+j_coord_offsetF,
1279 x+j_coord_offsetG,x+j_coord_offsetH,
1280 &jx0,&jy0,&jz0,&jx1,&jy1,&jz1,&jx2,
1281 &jy2,&jz2,&jx3,&jy3,&jz3);
1283 /* Calculate displacement vector */
1284 dx00 = _mm256_sub_ps(ix0,jx0);
1285 dy00 = _mm256_sub_ps(iy0,jy0);
1286 dz00 = _mm256_sub_ps(iz0,jz0);
1287 dx11 = _mm256_sub_ps(ix1,jx1);
1288 dy11 = _mm256_sub_ps(iy1,jy1);
1289 dz11 = _mm256_sub_ps(iz1,jz1);
1290 dx12 = _mm256_sub_ps(ix1,jx2);
1291 dy12 = _mm256_sub_ps(iy1,jy2);
1292 dz12 = _mm256_sub_ps(iz1,jz2);
1293 dx13 = _mm256_sub_ps(ix1,jx3);
1294 dy13 = _mm256_sub_ps(iy1,jy3);
1295 dz13 = _mm256_sub_ps(iz1,jz3);
1296 dx21 = _mm256_sub_ps(ix2,jx1);
1297 dy21 = _mm256_sub_ps(iy2,jy1);
1298 dz21 = _mm256_sub_ps(iz2,jz1);
1299 dx22 = _mm256_sub_ps(ix2,jx2);
1300 dy22 = _mm256_sub_ps(iy2,jy2);
1301 dz22 = _mm256_sub_ps(iz2,jz2);
1302 dx23 = _mm256_sub_ps(ix2,jx3);
1303 dy23 = _mm256_sub_ps(iy2,jy3);
1304 dz23 = _mm256_sub_ps(iz2,jz3);
1305 dx31 = _mm256_sub_ps(ix3,jx1);
1306 dy31 = _mm256_sub_ps(iy3,jy1);
1307 dz31 = _mm256_sub_ps(iz3,jz1);
1308 dx32 = _mm256_sub_ps(ix3,jx2);
1309 dy32 = _mm256_sub_ps(iy3,jy2);
1310 dz32 = _mm256_sub_ps(iz3,jz2);
1311 dx33 = _mm256_sub_ps(ix3,jx3);
1312 dy33 = _mm256_sub_ps(iy3,jy3);
1313 dz33 = _mm256_sub_ps(iz3,jz3);
1315 /* Calculate squared distance and things based on it */
1316 rsq00 = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
1317 rsq11 = gmx_mm256_calc_rsq_ps(dx11,dy11,dz11);
1318 rsq12 = gmx_mm256_calc_rsq_ps(dx12,dy12,dz12);
1319 rsq13 = gmx_mm256_calc_rsq_ps(dx13,dy13,dz13);
1320 rsq21 = gmx_mm256_calc_rsq_ps(dx21,dy21,dz21);
1321 rsq22 = gmx_mm256_calc_rsq_ps(dx22,dy22,dz22);
1322 rsq23 = gmx_mm256_calc_rsq_ps(dx23,dy23,dz23);
1323 rsq31 = gmx_mm256_calc_rsq_ps(dx31,dy31,dz31);
1324 rsq32 = gmx_mm256_calc_rsq_ps(dx32,dy32,dz32);
1325 rsq33 = gmx_mm256_calc_rsq_ps(dx33,dy33,dz33);
1327 rinv11 = gmx_mm256_invsqrt_ps(rsq11);
1328 rinv12 = gmx_mm256_invsqrt_ps(rsq12);
1329 rinv13 = gmx_mm256_invsqrt_ps(rsq13);
1330 rinv21 = gmx_mm256_invsqrt_ps(rsq21);
1331 rinv22 = gmx_mm256_invsqrt_ps(rsq22);
1332 rinv23 = gmx_mm256_invsqrt_ps(rsq23);
1333 rinv31 = gmx_mm256_invsqrt_ps(rsq31);
1334 rinv32 = gmx_mm256_invsqrt_ps(rsq32);
1335 rinv33 = gmx_mm256_invsqrt_ps(rsq33);
1337 rinvsq00 = gmx_mm256_inv_ps(rsq00);
1338 rinvsq11 = _mm256_mul_ps(rinv11,rinv11);
1339 rinvsq12 = _mm256_mul_ps(rinv12,rinv12);
1340 rinvsq13 = _mm256_mul_ps(rinv13,rinv13);
1341 rinvsq21 = _mm256_mul_ps(rinv21,rinv21);
1342 rinvsq22 = _mm256_mul_ps(rinv22,rinv22);
1343 rinvsq23 = _mm256_mul_ps(rinv23,rinv23);
1344 rinvsq31 = _mm256_mul_ps(rinv31,rinv31);
1345 rinvsq32 = _mm256_mul_ps(rinv32,rinv32);
1346 rinvsq33 = _mm256_mul_ps(rinv33,rinv33);
1348 fjx0 = _mm256_setzero_ps();
1349 fjy0 = _mm256_setzero_ps();
1350 fjz0 = _mm256_setzero_ps();
1351 fjx1 = _mm256_setzero_ps();
1352 fjy1 = _mm256_setzero_ps();
1353 fjz1 = _mm256_setzero_ps();
1354 fjx2 = _mm256_setzero_ps();
1355 fjy2 = _mm256_setzero_ps();
1356 fjz2 = _mm256_setzero_ps();
1357 fjx3 = _mm256_setzero_ps();
1358 fjy3 = _mm256_setzero_ps();
1359 fjz3 = _mm256_setzero_ps();
1361 /**************************
1362 * CALCULATE INTERACTIONS *
1363 **************************/
1365 /* LENNARD-JONES DISPERSION/REPULSION */
1367 rinvsix = _mm256_mul_ps(_mm256_mul_ps(rinvsq00,rinvsq00),rinvsq00);
1368 fvdw = _mm256_mul_ps(_mm256_sub_ps(_mm256_mul_ps(c12_00,rinvsix),c6_00),_mm256_mul_ps(rinvsix,rinvsq00));
1372 /* Calculate temporary vectorial force */
1373 tx = _mm256_mul_ps(fscal,dx00);
1374 ty = _mm256_mul_ps(fscal,dy00);
1375 tz = _mm256_mul_ps(fscal,dz00);
1377 /* Update vectorial force */
1378 fix0 = _mm256_add_ps(fix0,tx);
1379 fiy0 = _mm256_add_ps(fiy0,ty);
1380 fiz0 = _mm256_add_ps(fiz0,tz);
1382 fjx0 = _mm256_add_ps(fjx0,tx);
1383 fjy0 = _mm256_add_ps(fjy0,ty);
1384 fjz0 = _mm256_add_ps(fjz0,tz);
1386 /**************************
1387 * CALCULATE INTERACTIONS *
1388 **************************/
1390 /* REACTION-FIELD ELECTROSTATICS */
1391 felec = _mm256_mul_ps(qq11,_mm256_sub_ps(_mm256_mul_ps(rinv11,rinvsq11),krf2));
1395 /* Calculate temporary vectorial force */
1396 tx = _mm256_mul_ps(fscal,dx11);
1397 ty = _mm256_mul_ps(fscal,dy11);
1398 tz = _mm256_mul_ps(fscal,dz11);
1400 /* Update vectorial force */
1401 fix1 = _mm256_add_ps(fix1,tx);
1402 fiy1 = _mm256_add_ps(fiy1,ty);
1403 fiz1 = _mm256_add_ps(fiz1,tz);
1405 fjx1 = _mm256_add_ps(fjx1,tx);
1406 fjy1 = _mm256_add_ps(fjy1,ty);
1407 fjz1 = _mm256_add_ps(fjz1,tz);
1409 /**************************
1410 * CALCULATE INTERACTIONS *
1411 **************************/
1413 /* REACTION-FIELD ELECTROSTATICS */
1414 felec = _mm256_mul_ps(qq12,_mm256_sub_ps(_mm256_mul_ps(rinv12,rinvsq12),krf2));
1418 /* Calculate temporary vectorial force */
1419 tx = _mm256_mul_ps(fscal,dx12);
1420 ty = _mm256_mul_ps(fscal,dy12);
1421 tz = _mm256_mul_ps(fscal,dz12);
1423 /* Update vectorial force */
1424 fix1 = _mm256_add_ps(fix1,tx);
1425 fiy1 = _mm256_add_ps(fiy1,ty);
1426 fiz1 = _mm256_add_ps(fiz1,tz);
1428 fjx2 = _mm256_add_ps(fjx2,tx);
1429 fjy2 = _mm256_add_ps(fjy2,ty);
1430 fjz2 = _mm256_add_ps(fjz2,tz);
1432 /**************************
1433 * CALCULATE INTERACTIONS *
1434 **************************/
1436 /* REACTION-FIELD ELECTROSTATICS */
1437 felec = _mm256_mul_ps(qq13,_mm256_sub_ps(_mm256_mul_ps(rinv13,rinvsq13),krf2));
1441 /* Calculate temporary vectorial force */
1442 tx = _mm256_mul_ps(fscal,dx13);
1443 ty = _mm256_mul_ps(fscal,dy13);
1444 tz = _mm256_mul_ps(fscal,dz13);
1446 /* Update vectorial force */
1447 fix1 = _mm256_add_ps(fix1,tx);
1448 fiy1 = _mm256_add_ps(fiy1,ty);
1449 fiz1 = _mm256_add_ps(fiz1,tz);
1451 fjx3 = _mm256_add_ps(fjx3,tx);
1452 fjy3 = _mm256_add_ps(fjy3,ty);
1453 fjz3 = _mm256_add_ps(fjz3,tz);
1455 /**************************
1456 * CALCULATE INTERACTIONS *
1457 **************************/
1459 /* REACTION-FIELD ELECTROSTATICS */
1460 felec = _mm256_mul_ps(qq21,_mm256_sub_ps(_mm256_mul_ps(rinv21,rinvsq21),krf2));
1464 /* Calculate temporary vectorial force */
1465 tx = _mm256_mul_ps(fscal,dx21);
1466 ty = _mm256_mul_ps(fscal,dy21);
1467 tz = _mm256_mul_ps(fscal,dz21);
1469 /* Update vectorial force */
1470 fix2 = _mm256_add_ps(fix2,tx);
1471 fiy2 = _mm256_add_ps(fiy2,ty);
1472 fiz2 = _mm256_add_ps(fiz2,tz);
1474 fjx1 = _mm256_add_ps(fjx1,tx);
1475 fjy1 = _mm256_add_ps(fjy1,ty);
1476 fjz1 = _mm256_add_ps(fjz1,tz);
1478 /**************************
1479 * CALCULATE INTERACTIONS *
1480 **************************/
1482 /* REACTION-FIELD ELECTROSTATICS */
1483 felec = _mm256_mul_ps(qq22,_mm256_sub_ps(_mm256_mul_ps(rinv22,rinvsq22),krf2));
1487 /* Calculate temporary vectorial force */
1488 tx = _mm256_mul_ps(fscal,dx22);
1489 ty = _mm256_mul_ps(fscal,dy22);
1490 tz = _mm256_mul_ps(fscal,dz22);
1492 /* Update vectorial force */
1493 fix2 = _mm256_add_ps(fix2,tx);
1494 fiy2 = _mm256_add_ps(fiy2,ty);
1495 fiz2 = _mm256_add_ps(fiz2,tz);
1497 fjx2 = _mm256_add_ps(fjx2,tx);
1498 fjy2 = _mm256_add_ps(fjy2,ty);
1499 fjz2 = _mm256_add_ps(fjz2,tz);
1501 /**************************
1502 * CALCULATE INTERACTIONS *
1503 **************************/
1505 /* REACTION-FIELD ELECTROSTATICS */
1506 felec = _mm256_mul_ps(qq23,_mm256_sub_ps(_mm256_mul_ps(rinv23,rinvsq23),krf2));
1510 /* Calculate temporary vectorial force */
1511 tx = _mm256_mul_ps(fscal,dx23);
1512 ty = _mm256_mul_ps(fscal,dy23);
1513 tz = _mm256_mul_ps(fscal,dz23);
1515 /* Update vectorial force */
1516 fix2 = _mm256_add_ps(fix2,tx);
1517 fiy2 = _mm256_add_ps(fiy2,ty);
1518 fiz2 = _mm256_add_ps(fiz2,tz);
1520 fjx3 = _mm256_add_ps(fjx3,tx);
1521 fjy3 = _mm256_add_ps(fjy3,ty);
1522 fjz3 = _mm256_add_ps(fjz3,tz);
1524 /**************************
1525 * CALCULATE INTERACTIONS *
1526 **************************/
1528 /* REACTION-FIELD ELECTROSTATICS */
1529 felec = _mm256_mul_ps(qq31,_mm256_sub_ps(_mm256_mul_ps(rinv31,rinvsq31),krf2));
1533 /* Calculate temporary vectorial force */
1534 tx = _mm256_mul_ps(fscal,dx31);
1535 ty = _mm256_mul_ps(fscal,dy31);
1536 tz = _mm256_mul_ps(fscal,dz31);
1538 /* Update vectorial force */
1539 fix3 = _mm256_add_ps(fix3,tx);
1540 fiy3 = _mm256_add_ps(fiy3,ty);
1541 fiz3 = _mm256_add_ps(fiz3,tz);
1543 fjx1 = _mm256_add_ps(fjx1,tx);
1544 fjy1 = _mm256_add_ps(fjy1,ty);
1545 fjz1 = _mm256_add_ps(fjz1,tz);
1547 /**************************
1548 * CALCULATE INTERACTIONS *
1549 **************************/
1551 /* REACTION-FIELD ELECTROSTATICS */
1552 felec = _mm256_mul_ps(qq32,_mm256_sub_ps(_mm256_mul_ps(rinv32,rinvsq32),krf2));
1556 /* Calculate temporary vectorial force */
1557 tx = _mm256_mul_ps(fscal,dx32);
1558 ty = _mm256_mul_ps(fscal,dy32);
1559 tz = _mm256_mul_ps(fscal,dz32);
1561 /* Update vectorial force */
1562 fix3 = _mm256_add_ps(fix3,tx);
1563 fiy3 = _mm256_add_ps(fiy3,ty);
1564 fiz3 = _mm256_add_ps(fiz3,tz);
1566 fjx2 = _mm256_add_ps(fjx2,tx);
1567 fjy2 = _mm256_add_ps(fjy2,ty);
1568 fjz2 = _mm256_add_ps(fjz2,tz);
1570 /**************************
1571 * CALCULATE INTERACTIONS *
1572 **************************/
1574 /* REACTION-FIELD ELECTROSTATICS */
1575 felec = _mm256_mul_ps(qq33,_mm256_sub_ps(_mm256_mul_ps(rinv33,rinvsq33),krf2));
1579 /* Calculate temporary vectorial force */
1580 tx = _mm256_mul_ps(fscal,dx33);
1581 ty = _mm256_mul_ps(fscal,dy33);
1582 tz = _mm256_mul_ps(fscal,dz33);
1584 /* Update vectorial force */
1585 fix3 = _mm256_add_ps(fix3,tx);
1586 fiy3 = _mm256_add_ps(fiy3,ty);
1587 fiz3 = _mm256_add_ps(fiz3,tz);
1589 fjx3 = _mm256_add_ps(fjx3,tx);
1590 fjy3 = _mm256_add_ps(fjy3,ty);
1591 fjz3 = _mm256_add_ps(fjz3,tz);
1593 fjptrA = f+j_coord_offsetA;
1594 fjptrB = f+j_coord_offsetB;
1595 fjptrC = f+j_coord_offsetC;
1596 fjptrD = f+j_coord_offsetD;
1597 fjptrE = f+j_coord_offsetE;
1598 fjptrF = f+j_coord_offsetF;
1599 fjptrG = f+j_coord_offsetG;
1600 fjptrH = f+j_coord_offsetH;
1602 gmx_mm256_decrement_4rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,
1603 fjx0,fjy0,fjz0,fjx1,fjy1,fjz1,
1604 fjx2,fjy2,fjz2,fjx3,fjy3,fjz3);
1606 /* Inner loop uses 273 flops */
1609 if(jidx<j_index_end)
1612 /* Get j neighbor index, and coordinate index */
1613 jnrlistA = jjnr[jidx];
1614 jnrlistB = jjnr[jidx+1];
1615 jnrlistC = jjnr[jidx+2];
1616 jnrlistD = jjnr[jidx+3];
1617 jnrlistE = jjnr[jidx+4];
1618 jnrlistF = jjnr[jidx+5];
1619 jnrlistG = jjnr[jidx+6];
1620 jnrlistH = jjnr[jidx+7];
1621 /* Sign of each element will be negative for non-real atoms.
1622 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
1623 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
1625 dummy_mask = gmx_mm256_set_m128(gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx+4)),_mm_setzero_si128())),
1626 gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128())));
1628 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
1629 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
1630 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
1631 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
1632 jnrE = (jnrlistE>=0) ? jnrlistE : 0;
1633 jnrF = (jnrlistF>=0) ? jnrlistF : 0;
1634 jnrG = (jnrlistG>=0) ? jnrlistG : 0;
1635 jnrH = (jnrlistH>=0) ? jnrlistH : 0;
1636 j_coord_offsetA = DIM*jnrA;
1637 j_coord_offsetB = DIM*jnrB;
1638 j_coord_offsetC = DIM*jnrC;
1639 j_coord_offsetD = DIM*jnrD;
1640 j_coord_offsetE = DIM*jnrE;
1641 j_coord_offsetF = DIM*jnrF;
1642 j_coord_offsetG = DIM*jnrG;
1643 j_coord_offsetH = DIM*jnrH;
1645 /* load j atom coordinates */
1646 gmx_mm256_load_4rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
1647 x+j_coord_offsetC,x+j_coord_offsetD,
1648 x+j_coord_offsetE,x+j_coord_offsetF,
1649 x+j_coord_offsetG,x+j_coord_offsetH,
1650 &jx0,&jy0,&jz0,&jx1,&jy1,&jz1,&jx2,
1651 &jy2,&jz2,&jx3,&jy3,&jz3);
1653 /* Calculate displacement vector */
1654 dx00 = _mm256_sub_ps(ix0,jx0);
1655 dy00 = _mm256_sub_ps(iy0,jy0);
1656 dz00 = _mm256_sub_ps(iz0,jz0);
1657 dx11 = _mm256_sub_ps(ix1,jx1);
1658 dy11 = _mm256_sub_ps(iy1,jy1);
1659 dz11 = _mm256_sub_ps(iz1,jz1);
1660 dx12 = _mm256_sub_ps(ix1,jx2);
1661 dy12 = _mm256_sub_ps(iy1,jy2);
1662 dz12 = _mm256_sub_ps(iz1,jz2);
1663 dx13 = _mm256_sub_ps(ix1,jx3);
1664 dy13 = _mm256_sub_ps(iy1,jy3);
1665 dz13 = _mm256_sub_ps(iz1,jz3);
1666 dx21 = _mm256_sub_ps(ix2,jx1);
1667 dy21 = _mm256_sub_ps(iy2,jy1);
1668 dz21 = _mm256_sub_ps(iz2,jz1);
1669 dx22 = _mm256_sub_ps(ix2,jx2);
1670 dy22 = _mm256_sub_ps(iy2,jy2);
1671 dz22 = _mm256_sub_ps(iz2,jz2);
1672 dx23 = _mm256_sub_ps(ix2,jx3);
1673 dy23 = _mm256_sub_ps(iy2,jy3);
1674 dz23 = _mm256_sub_ps(iz2,jz3);
1675 dx31 = _mm256_sub_ps(ix3,jx1);
1676 dy31 = _mm256_sub_ps(iy3,jy1);
1677 dz31 = _mm256_sub_ps(iz3,jz1);
1678 dx32 = _mm256_sub_ps(ix3,jx2);
1679 dy32 = _mm256_sub_ps(iy3,jy2);
1680 dz32 = _mm256_sub_ps(iz3,jz2);
1681 dx33 = _mm256_sub_ps(ix3,jx3);
1682 dy33 = _mm256_sub_ps(iy3,jy3);
1683 dz33 = _mm256_sub_ps(iz3,jz3);
1685 /* Calculate squared distance and things based on it */
1686 rsq00 = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
1687 rsq11 = gmx_mm256_calc_rsq_ps(dx11,dy11,dz11);
1688 rsq12 = gmx_mm256_calc_rsq_ps(dx12,dy12,dz12);
1689 rsq13 = gmx_mm256_calc_rsq_ps(dx13,dy13,dz13);
1690 rsq21 = gmx_mm256_calc_rsq_ps(dx21,dy21,dz21);
1691 rsq22 = gmx_mm256_calc_rsq_ps(dx22,dy22,dz22);
1692 rsq23 = gmx_mm256_calc_rsq_ps(dx23,dy23,dz23);
1693 rsq31 = gmx_mm256_calc_rsq_ps(dx31,dy31,dz31);
1694 rsq32 = gmx_mm256_calc_rsq_ps(dx32,dy32,dz32);
1695 rsq33 = gmx_mm256_calc_rsq_ps(dx33,dy33,dz33);
1697 rinv11 = gmx_mm256_invsqrt_ps(rsq11);
1698 rinv12 = gmx_mm256_invsqrt_ps(rsq12);
1699 rinv13 = gmx_mm256_invsqrt_ps(rsq13);
1700 rinv21 = gmx_mm256_invsqrt_ps(rsq21);
1701 rinv22 = gmx_mm256_invsqrt_ps(rsq22);
1702 rinv23 = gmx_mm256_invsqrt_ps(rsq23);
1703 rinv31 = gmx_mm256_invsqrt_ps(rsq31);
1704 rinv32 = gmx_mm256_invsqrt_ps(rsq32);
1705 rinv33 = gmx_mm256_invsqrt_ps(rsq33);
1707 rinvsq00 = gmx_mm256_inv_ps(rsq00);
1708 rinvsq11 = _mm256_mul_ps(rinv11,rinv11);
1709 rinvsq12 = _mm256_mul_ps(rinv12,rinv12);
1710 rinvsq13 = _mm256_mul_ps(rinv13,rinv13);
1711 rinvsq21 = _mm256_mul_ps(rinv21,rinv21);
1712 rinvsq22 = _mm256_mul_ps(rinv22,rinv22);
1713 rinvsq23 = _mm256_mul_ps(rinv23,rinv23);
1714 rinvsq31 = _mm256_mul_ps(rinv31,rinv31);
1715 rinvsq32 = _mm256_mul_ps(rinv32,rinv32);
1716 rinvsq33 = _mm256_mul_ps(rinv33,rinv33);
1718 fjx0 = _mm256_setzero_ps();
1719 fjy0 = _mm256_setzero_ps();
1720 fjz0 = _mm256_setzero_ps();
1721 fjx1 = _mm256_setzero_ps();
1722 fjy1 = _mm256_setzero_ps();
1723 fjz1 = _mm256_setzero_ps();
1724 fjx2 = _mm256_setzero_ps();
1725 fjy2 = _mm256_setzero_ps();
1726 fjz2 = _mm256_setzero_ps();
1727 fjx3 = _mm256_setzero_ps();
1728 fjy3 = _mm256_setzero_ps();
1729 fjz3 = _mm256_setzero_ps();
1731 /**************************
1732 * CALCULATE INTERACTIONS *
1733 **************************/
1735 /* LENNARD-JONES DISPERSION/REPULSION */
1737 rinvsix = _mm256_mul_ps(_mm256_mul_ps(rinvsq00,rinvsq00),rinvsq00);
1738 fvdw = _mm256_mul_ps(_mm256_sub_ps(_mm256_mul_ps(c12_00,rinvsix),c6_00),_mm256_mul_ps(rinvsix,rinvsq00));
1742 fscal = _mm256_andnot_ps(dummy_mask,fscal);
1744 /* Calculate temporary vectorial force */
1745 tx = _mm256_mul_ps(fscal,dx00);
1746 ty = _mm256_mul_ps(fscal,dy00);
1747 tz = _mm256_mul_ps(fscal,dz00);
1749 /* Update vectorial force */
1750 fix0 = _mm256_add_ps(fix0,tx);
1751 fiy0 = _mm256_add_ps(fiy0,ty);
1752 fiz0 = _mm256_add_ps(fiz0,tz);
1754 fjx0 = _mm256_add_ps(fjx0,tx);
1755 fjy0 = _mm256_add_ps(fjy0,ty);
1756 fjz0 = _mm256_add_ps(fjz0,tz);
1758 /**************************
1759 * CALCULATE INTERACTIONS *
1760 **************************/
1762 /* REACTION-FIELD ELECTROSTATICS */
1763 felec = _mm256_mul_ps(qq11,_mm256_sub_ps(_mm256_mul_ps(rinv11,rinvsq11),krf2));
1767 fscal = _mm256_andnot_ps(dummy_mask,fscal);
1769 /* Calculate temporary vectorial force */
1770 tx = _mm256_mul_ps(fscal,dx11);
1771 ty = _mm256_mul_ps(fscal,dy11);
1772 tz = _mm256_mul_ps(fscal,dz11);
1774 /* Update vectorial force */
1775 fix1 = _mm256_add_ps(fix1,tx);
1776 fiy1 = _mm256_add_ps(fiy1,ty);
1777 fiz1 = _mm256_add_ps(fiz1,tz);
1779 fjx1 = _mm256_add_ps(fjx1,tx);
1780 fjy1 = _mm256_add_ps(fjy1,ty);
1781 fjz1 = _mm256_add_ps(fjz1,tz);
1783 /**************************
1784 * CALCULATE INTERACTIONS *
1785 **************************/
1787 /* REACTION-FIELD ELECTROSTATICS */
1788 felec = _mm256_mul_ps(qq12,_mm256_sub_ps(_mm256_mul_ps(rinv12,rinvsq12),krf2));
1792 fscal = _mm256_andnot_ps(dummy_mask,fscal);
1794 /* Calculate temporary vectorial force */
1795 tx = _mm256_mul_ps(fscal,dx12);
1796 ty = _mm256_mul_ps(fscal,dy12);
1797 tz = _mm256_mul_ps(fscal,dz12);
1799 /* Update vectorial force */
1800 fix1 = _mm256_add_ps(fix1,tx);
1801 fiy1 = _mm256_add_ps(fiy1,ty);
1802 fiz1 = _mm256_add_ps(fiz1,tz);
1804 fjx2 = _mm256_add_ps(fjx2,tx);
1805 fjy2 = _mm256_add_ps(fjy2,ty);
1806 fjz2 = _mm256_add_ps(fjz2,tz);
1808 /**************************
1809 * CALCULATE INTERACTIONS *
1810 **************************/
1812 /* REACTION-FIELD ELECTROSTATICS */
1813 felec = _mm256_mul_ps(qq13,_mm256_sub_ps(_mm256_mul_ps(rinv13,rinvsq13),krf2));
1817 fscal = _mm256_andnot_ps(dummy_mask,fscal);
1819 /* Calculate temporary vectorial force */
1820 tx = _mm256_mul_ps(fscal,dx13);
1821 ty = _mm256_mul_ps(fscal,dy13);
1822 tz = _mm256_mul_ps(fscal,dz13);
1824 /* Update vectorial force */
1825 fix1 = _mm256_add_ps(fix1,tx);
1826 fiy1 = _mm256_add_ps(fiy1,ty);
1827 fiz1 = _mm256_add_ps(fiz1,tz);
1829 fjx3 = _mm256_add_ps(fjx3,tx);
1830 fjy3 = _mm256_add_ps(fjy3,ty);
1831 fjz3 = _mm256_add_ps(fjz3,tz);
1833 /**************************
1834 * CALCULATE INTERACTIONS *
1835 **************************/
1837 /* REACTION-FIELD ELECTROSTATICS */
1838 felec = _mm256_mul_ps(qq21,_mm256_sub_ps(_mm256_mul_ps(rinv21,rinvsq21),krf2));
1842 fscal = _mm256_andnot_ps(dummy_mask,fscal);
1844 /* Calculate temporary vectorial force */
1845 tx = _mm256_mul_ps(fscal,dx21);
1846 ty = _mm256_mul_ps(fscal,dy21);
1847 tz = _mm256_mul_ps(fscal,dz21);
1849 /* Update vectorial force */
1850 fix2 = _mm256_add_ps(fix2,tx);
1851 fiy2 = _mm256_add_ps(fiy2,ty);
1852 fiz2 = _mm256_add_ps(fiz2,tz);
1854 fjx1 = _mm256_add_ps(fjx1,tx);
1855 fjy1 = _mm256_add_ps(fjy1,ty);
1856 fjz1 = _mm256_add_ps(fjz1,tz);
1858 /**************************
1859 * CALCULATE INTERACTIONS *
1860 **************************/
1862 /* REACTION-FIELD ELECTROSTATICS */
1863 felec = _mm256_mul_ps(qq22,_mm256_sub_ps(_mm256_mul_ps(rinv22,rinvsq22),krf2));
1867 fscal = _mm256_andnot_ps(dummy_mask,fscal);
1869 /* Calculate temporary vectorial force */
1870 tx = _mm256_mul_ps(fscal,dx22);
1871 ty = _mm256_mul_ps(fscal,dy22);
1872 tz = _mm256_mul_ps(fscal,dz22);
1874 /* Update vectorial force */
1875 fix2 = _mm256_add_ps(fix2,tx);
1876 fiy2 = _mm256_add_ps(fiy2,ty);
1877 fiz2 = _mm256_add_ps(fiz2,tz);
1879 fjx2 = _mm256_add_ps(fjx2,tx);
1880 fjy2 = _mm256_add_ps(fjy2,ty);
1881 fjz2 = _mm256_add_ps(fjz2,tz);
1883 /**************************
1884 * CALCULATE INTERACTIONS *
1885 **************************/
1887 /* REACTION-FIELD ELECTROSTATICS */
1888 felec = _mm256_mul_ps(qq23,_mm256_sub_ps(_mm256_mul_ps(rinv23,rinvsq23),krf2));
1892 fscal = _mm256_andnot_ps(dummy_mask,fscal);
1894 /* Calculate temporary vectorial force */
1895 tx = _mm256_mul_ps(fscal,dx23);
1896 ty = _mm256_mul_ps(fscal,dy23);
1897 tz = _mm256_mul_ps(fscal,dz23);
1899 /* Update vectorial force */
1900 fix2 = _mm256_add_ps(fix2,tx);
1901 fiy2 = _mm256_add_ps(fiy2,ty);
1902 fiz2 = _mm256_add_ps(fiz2,tz);
1904 fjx3 = _mm256_add_ps(fjx3,tx);
1905 fjy3 = _mm256_add_ps(fjy3,ty);
1906 fjz3 = _mm256_add_ps(fjz3,tz);
1908 /**************************
1909 * CALCULATE INTERACTIONS *
1910 **************************/
1912 /* REACTION-FIELD ELECTROSTATICS */
1913 felec = _mm256_mul_ps(qq31,_mm256_sub_ps(_mm256_mul_ps(rinv31,rinvsq31),krf2));
1917 fscal = _mm256_andnot_ps(dummy_mask,fscal);
1919 /* Calculate temporary vectorial force */
1920 tx = _mm256_mul_ps(fscal,dx31);
1921 ty = _mm256_mul_ps(fscal,dy31);
1922 tz = _mm256_mul_ps(fscal,dz31);
1924 /* Update vectorial force */
1925 fix3 = _mm256_add_ps(fix3,tx);
1926 fiy3 = _mm256_add_ps(fiy3,ty);
1927 fiz3 = _mm256_add_ps(fiz3,tz);
1929 fjx1 = _mm256_add_ps(fjx1,tx);
1930 fjy1 = _mm256_add_ps(fjy1,ty);
1931 fjz1 = _mm256_add_ps(fjz1,tz);
1933 /**************************
1934 * CALCULATE INTERACTIONS *
1935 **************************/
1937 /* REACTION-FIELD ELECTROSTATICS */
1938 felec = _mm256_mul_ps(qq32,_mm256_sub_ps(_mm256_mul_ps(rinv32,rinvsq32),krf2));
1942 fscal = _mm256_andnot_ps(dummy_mask,fscal);
1944 /* Calculate temporary vectorial force */
1945 tx = _mm256_mul_ps(fscal,dx32);
1946 ty = _mm256_mul_ps(fscal,dy32);
1947 tz = _mm256_mul_ps(fscal,dz32);
1949 /* Update vectorial force */
1950 fix3 = _mm256_add_ps(fix3,tx);
1951 fiy3 = _mm256_add_ps(fiy3,ty);
1952 fiz3 = _mm256_add_ps(fiz3,tz);
1954 fjx2 = _mm256_add_ps(fjx2,tx);
1955 fjy2 = _mm256_add_ps(fjy2,ty);
1956 fjz2 = _mm256_add_ps(fjz2,tz);
1958 /**************************
1959 * CALCULATE INTERACTIONS *
1960 **************************/
1962 /* REACTION-FIELD ELECTROSTATICS */
1963 felec = _mm256_mul_ps(qq33,_mm256_sub_ps(_mm256_mul_ps(rinv33,rinvsq33),krf2));
1967 fscal = _mm256_andnot_ps(dummy_mask,fscal);
1969 /* Calculate temporary vectorial force */
1970 tx = _mm256_mul_ps(fscal,dx33);
1971 ty = _mm256_mul_ps(fscal,dy33);
1972 tz = _mm256_mul_ps(fscal,dz33);
1974 /* Update vectorial force */
1975 fix3 = _mm256_add_ps(fix3,tx);
1976 fiy3 = _mm256_add_ps(fiy3,ty);
1977 fiz3 = _mm256_add_ps(fiz3,tz);
1979 fjx3 = _mm256_add_ps(fjx3,tx);
1980 fjy3 = _mm256_add_ps(fjy3,ty);
1981 fjz3 = _mm256_add_ps(fjz3,tz);
1983 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
1984 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
1985 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
1986 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
1987 fjptrE = (jnrlistE>=0) ? f+j_coord_offsetE : scratch;
1988 fjptrF = (jnrlistF>=0) ? f+j_coord_offsetF : scratch;
1989 fjptrG = (jnrlistG>=0) ? f+j_coord_offsetG : scratch;
1990 fjptrH = (jnrlistH>=0) ? f+j_coord_offsetH : scratch;
1992 gmx_mm256_decrement_4rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,
1993 fjx0,fjy0,fjz0,fjx1,fjy1,fjz1,
1994 fjx2,fjy2,fjz2,fjx3,fjy3,fjz3);
1996 /* Inner loop uses 273 flops */
1999 /* End of innermost loop */
2001 gmx_mm256_update_iforce_4atom_swizzle_ps(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
2002 f+i_coord_offset,fshift+i_shift_offset);
2004 /* Increment number of inner iterations */
2005 inneriter += j_index_end - j_index_start;
2007 /* Outer loop uses 24 flops */
2010 /* Increment number of outer iterations */
2013 /* Update outer/inner flops */
2015 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W4W4_F,outeriter*24 + inneriter*273);