2 * Note: this file was generated by the Gromacs avx_256_double kernel generator.
4 * This source code is part of
8 * Copyright (c) 2001-2012, The GROMACS Development Team
10 * Gromacs is a library for molecular simulation and trajectory analysis,
11 * written by Erik Lindahl, David van der Spoel, Berk Hess, and others - for
12 * a full list of developers and information, check out http://www.gromacs.org
14 * This program is free software; you can redistribute it and/or modify it under
15 * the terms of the GNU Lesser General Public License as published by the Free
16 * Software Foundation; either version 2 of the License, or (at your option) any
19 * To help fund GROMACS development, we humbly ask that you cite
20 * the papers people have written on it - you can find them on the website.
28 #include "../nb_kernel.h"
29 #include "types/simple.h"
33 #include "gmx_math_x86_avx_256_double.h"
34 #include "kernelutil_x86_avx_256_double.h"
37 * Gromacs nonbonded kernel: nb_kernel_ElecRF_VdwLJ_GeomW3W3_VF_avx_256_double
38 * Electrostatics interaction: ReactionField
39 * VdW interaction: LennardJones
40 * Geometry: Water3-Water3
41 * Calculate force/pot: PotentialAndForce
44 nb_kernel_ElecRF_VdwLJ_GeomW3W3_VF_avx_256_double
45 (t_nblist * gmx_restrict nlist,
46 rvec * gmx_restrict xx,
47 rvec * gmx_restrict ff,
48 t_forcerec * gmx_restrict fr,
49 t_mdatoms * gmx_restrict mdatoms,
50 nb_kernel_data_t * gmx_restrict kernel_data,
51 t_nrnb * gmx_restrict nrnb)
53 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
54 * just 0 for non-waters.
55 * Suffixes A,B,C,D refer to j loop unrolling done with AVX, e.g. for the four different
56 * jnr indices corresponding to data put in the four positions in the SIMD register.
58 int i_shift_offset,i_coord_offset,outeriter,inneriter;
59 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
60 int jnrA,jnrB,jnrC,jnrD;
61 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
62 int jnrlistE,jnrlistF,jnrlistG,jnrlistH;
63 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
64 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
66 real *shiftvec,*fshift,*x,*f;
67 real *fjptrA,*fjptrB,*fjptrC,*fjptrD;
69 __m256d tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
70 real * vdwioffsetptr0;
71 __m256d ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
72 real * vdwioffsetptr1;
73 __m256d ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
74 real * vdwioffsetptr2;
75 __m256d ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
76 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
77 __m256d jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
78 int vdwjidx1A,vdwjidx1B,vdwjidx1C,vdwjidx1D;
79 __m256d jx1,jy1,jz1,fjx1,fjy1,fjz1,jq1,isaj1;
80 int vdwjidx2A,vdwjidx2B,vdwjidx2C,vdwjidx2D;
81 __m256d jx2,jy2,jz2,fjx2,fjy2,fjz2,jq2,isaj2;
82 __m256d dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
83 __m256d dx01,dy01,dz01,rsq01,rinv01,rinvsq01,r01,qq01,c6_01,c12_01;
84 __m256d dx02,dy02,dz02,rsq02,rinv02,rinvsq02,r02,qq02,c6_02,c12_02;
85 __m256d dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
86 __m256d dx11,dy11,dz11,rsq11,rinv11,rinvsq11,r11,qq11,c6_11,c12_11;
87 __m256d dx12,dy12,dz12,rsq12,rinv12,rinvsq12,r12,qq12,c6_12,c12_12;
88 __m256d dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
89 __m256d dx21,dy21,dz21,rsq21,rinv21,rinvsq21,r21,qq21,c6_21,c12_21;
90 __m256d dx22,dy22,dz22,rsq22,rinv22,rinvsq22,r22,qq22,c6_22,c12_22;
91 __m256d velec,felec,velecsum,facel,crf,krf,krf2;
94 __m256d rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
97 __m256d one_sixth = _mm256_set1_pd(1.0/6.0);
98 __m256d one_twelfth = _mm256_set1_pd(1.0/12.0);
99 __m256d dummy_mask,cutoff_mask;
100 __m128 tmpmask0,tmpmask1;
101 __m256d signbit = _mm256_castsi256_pd( _mm256_set1_epi32(0x80000000) );
102 __m256d one = _mm256_set1_pd(1.0);
103 __m256d two = _mm256_set1_pd(2.0);
109 jindex = nlist->jindex;
111 shiftidx = nlist->shift;
113 shiftvec = fr->shift_vec[0];
114 fshift = fr->fshift[0];
115 facel = _mm256_set1_pd(fr->epsfac);
116 charge = mdatoms->chargeA;
117 krf = _mm256_set1_pd(fr->ic->k_rf);
118 krf2 = _mm256_set1_pd(fr->ic->k_rf*2.0);
119 crf = _mm256_set1_pd(fr->ic->c_rf);
120 nvdwtype = fr->ntype;
122 vdwtype = mdatoms->typeA;
124 /* Setup water-specific parameters */
125 inr = nlist->iinr[0];
126 iq0 = _mm256_mul_pd(facel,_mm256_set1_pd(charge[inr+0]));
127 iq1 = _mm256_mul_pd(facel,_mm256_set1_pd(charge[inr+1]));
128 iq2 = _mm256_mul_pd(facel,_mm256_set1_pd(charge[inr+2]));
129 vdwioffsetptr0 = vdwparam+2*nvdwtype*vdwtype[inr+0];
131 jq0 = _mm256_set1_pd(charge[inr+0]);
132 jq1 = _mm256_set1_pd(charge[inr+1]);
133 jq2 = _mm256_set1_pd(charge[inr+2]);
134 vdwjidx0A = 2*vdwtype[inr+0];
135 qq00 = _mm256_mul_pd(iq0,jq0);
136 c6_00 = _mm256_set1_pd(vdwioffsetptr0[vdwjidx0A]);
137 c12_00 = _mm256_set1_pd(vdwioffsetptr0[vdwjidx0A+1]);
138 qq01 = _mm256_mul_pd(iq0,jq1);
139 qq02 = _mm256_mul_pd(iq0,jq2);
140 qq10 = _mm256_mul_pd(iq1,jq0);
141 qq11 = _mm256_mul_pd(iq1,jq1);
142 qq12 = _mm256_mul_pd(iq1,jq2);
143 qq20 = _mm256_mul_pd(iq2,jq0);
144 qq21 = _mm256_mul_pd(iq2,jq1);
145 qq22 = _mm256_mul_pd(iq2,jq2);
147 /* Avoid stupid compiler warnings */
148 jnrA = jnrB = jnrC = jnrD = 0;
157 for(iidx=0;iidx<4*DIM;iidx++)
162 /* Start outer loop over neighborlists */
163 for(iidx=0; iidx<nri; iidx++)
165 /* Load shift vector for this list */
166 i_shift_offset = DIM*shiftidx[iidx];
168 /* Load limits for loop over neighbors */
169 j_index_start = jindex[iidx];
170 j_index_end = jindex[iidx+1];
172 /* Get outer coordinate index */
174 i_coord_offset = DIM*inr;
176 /* Load i particle coords and add shift vector */
177 gmx_mm256_load_shift_and_3rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,
178 &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2);
180 fix0 = _mm256_setzero_pd();
181 fiy0 = _mm256_setzero_pd();
182 fiz0 = _mm256_setzero_pd();
183 fix1 = _mm256_setzero_pd();
184 fiy1 = _mm256_setzero_pd();
185 fiz1 = _mm256_setzero_pd();
186 fix2 = _mm256_setzero_pd();
187 fiy2 = _mm256_setzero_pd();
188 fiz2 = _mm256_setzero_pd();
190 /* Reset potential sums */
191 velecsum = _mm256_setzero_pd();
192 vvdwsum = _mm256_setzero_pd();
194 /* Start inner kernel loop */
195 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
198 /* Get j neighbor index, and coordinate index */
203 j_coord_offsetA = DIM*jnrA;
204 j_coord_offsetB = DIM*jnrB;
205 j_coord_offsetC = DIM*jnrC;
206 j_coord_offsetD = DIM*jnrD;
208 /* load j atom coordinates */
209 gmx_mm256_load_3rvec_4ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
210 x+j_coord_offsetC,x+j_coord_offsetD,
211 &jx0,&jy0,&jz0,&jx1,&jy1,&jz1,&jx2,&jy2,&jz2);
213 /* Calculate displacement vector */
214 dx00 = _mm256_sub_pd(ix0,jx0);
215 dy00 = _mm256_sub_pd(iy0,jy0);
216 dz00 = _mm256_sub_pd(iz0,jz0);
217 dx01 = _mm256_sub_pd(ix0,jx1);
218 dy01 = _mm256_sub_pd(iy0,jy1);
219 dz01 = _mm256_sub_pd(iz0,jz1);
220 dx02 = _mm256_sub_pd(ix0,jx2);
221 dy02 = _mm256_sub_pd(iy0,jy2);
222 dz02 = _mm256_sub_pd(iz0,jz2);
223 dx10 = _mm256_sub_pd(ix1,jx0);
224 dy10 = _mm256_sub_pd(iy1,jy0);
225 dz10 = _mm256_sub_pd(iz1,jz0);
226 dx11 = _mm256_sub_pd(ix1,jx1);
227 dy11 = _mm256_sub_pd(iy1,jy1);
228 dz11 = _mm256_sub_pd(iz1,jz1);
229 dx12 = _mm256_sub_pd(ix1,jx2);
230 dy12 = _mm256_sub_pd(iy1,jy2);
231 dz12 = _mm256_sub_pd(iz1,jz2);
232 dx20 = _mm256_sub_pd(ix2,jx0);
233 dy20 = _mm256_sub_pd(iy2,jy0);
234 dz20 = _mm256_sub_pd(iz2,jz0);
235 dx21 = _mm256_sub_pd(ix2,jx1);
236 dy21 = _mm256_sub_pd(iy2,jy1);
237 dz21 = _mm256_sub_pd(iz2,jz1);
238 dx22 = _mm256_sub_pd(ix2,jx2);
239 dy22 = _mm256_sub_pd(iy2,jy2);
240 dz22 = _mm256_sub_pd(iz2,jz2);
242 /* Calculate squared distance and things based on it */
243 rsq00 = gmx_mm256_calc_rsq_pd(dx00,dy00,dz00);
244 rsq01 = gmx_mm256_calc_rsq_pd(dx01,dy01,dz01);
245 rsq02 = gmx_mm256_calc_rsq_pd(dx02,dy02,dz02);
246 rsq10 = gmx_mm256_calc_rsq_pd(dx10,dy10,dz10);
247 rsq11 = gmx_mm256_calc_rsq_pd(dx11,dy11,dz11);
248 rsq12 = gmx_mm256_calc_rsq_pd(dx12,dy12,dz12);
249 rsq20 = gmx_mm256_calc_rsq_pd(dx20,dy20,dz20);
250 rsq21 = gmx_mm256_calc_rsq_pd(dx21,dy21,dz21);
251 rsq22 = gmx_mm256_calc_rsq_pd(dx22,dy22,dz22);
253 rinv00 = gmx_mm256_invsqrt_pd(rsq00);
254 rinv01 = gmx_mm256_invsqrt_pd(rsq01);
255 rinv02 = gmx_mm256_invsqrt_pd(rsq02);
256 rinv10 = gmx_mm256_invsqrt_pd(rsq10);
257 rinv11 = gmx_mm256_invsqrt_pd(rsq11);
258 rinv12 = gmx_mm256_invsqrt_pd(rsq12);
259 rinv20 = gmx_mm256_invsqrt_pd(rsq20);
260 rinv21 = gmx_mm256_invsqrt_pd(rsq21);
261 rinv22 = gmx_mm256_invsqrt_pd(rsq22);
263 rinvsq00 = _mm256_mul_pd(rinv00,rinv00);
264 rinvsq01 = _mm256_mul_pd(rinv01,rinv01);
265 rinvsq02 = _mm256_mul_pd(rinv02,rinv02);
266 rinvsq10 = _mm256_mul_pd(rinv10,rinv10);
267 rinvsq11 = _mm256_mul_pd(rinv11,rinv11);
268 rinvsq12 = _mm256_mul_pd(rinv12,rinv12);
269 rinvsq20 = _mm256_mul_pd(rinv20,rinv20);
270 rinvsq21 = _mm256_mul_pd(rinv21,rinv21);
271 rinvsq22 = _mm256_mul_pd(rinv22,rinv22);
273 fjx0 = _mm256_setzero_pd();
274 fjy0 = _mm256_setzero_pd();
275 fjz0 = _mm256_setzero_pd();
276 fjx1 = _mm256_setzero_pd();
277 fjy1 = _mm256_setzero_pd();
278 fjz1 = _mm256_setzero_pd();
279 fjx2 = _mm256_setzero_pd();
280 fjy2 = _mm256_setzero_pd();
281 fjz2 = _mm256_setzero_pd();
283 /**************************
284 * CALCULATE INTERACTIONS *
285 **************************/
287 /* REACTION-FIELD ELECTROSTATICS */
288 velec = _mm256_mul_pd(qq00,_mm256_sub_pd(_mm256_add_pd(rinv00,_mm256_mul_pd(krf,rsq00)),crf));
289 felec = _mm256_mul_pd(qq00,_mm256_sub_pd(_mm256_mul_pd(rinv00,rinvsq00),krf2));
291 /* LENNARD-JONES DISPERSION/REPULSION */
293 rinvsix = _mm256_mul_pd(_mm256_mul_pd(rinvsq00,rinvsq00),rinvsq00);
294 vvdw6 = _mm256_mul_pd(c6_00,rinvsix);
295 vvdw12 = _mm256_mul_pd(c12_00,_mm256_mul_pd(rinvsix,rinvsix));
296 vvdw = _mm256_sub_pd( _mm256_mul_pd(vvdw12,one_twelfth) , _mm256_mul_pd(vvdw6,one_sixth) );
297 fvdw = _mm256_mul_pd(_mm256_sub_pd(vvdw12,vvdw6),rinvsq00);
299 /* Update potential sum for this i atom from the interaction with this j atom. */
300 velecsum = _mm256_add_pd(velecsum,velec);
301 vvdwsum = _mm256_add_pd(vvdwsum,vvdw);
303 fscal = _mm256_add_pd(felec,fvdw);
305 /* Calculate temporary vectorial force */
306 tx = _mm256_mul_pd(fscal,dx00);
307 ty = _mm256_mul_pd(fscal,dy00);
308 tz = _mm256_mul_pd(fscal,dz00);
310 /* Update vectorial force */
311 fix0 = _mm256_add_pd(fix0,tx);
312 fiy0 = _mm256_add_pd(fiy0,ty);
313 fiz0 = _mm256_add_pd(fiz0,tz);
315 fjx0 = _mm256_add_pd(fjx0,tx);
316 fjy0 = _mm256_add_pd(fjy0,ty);
317 fjz0 = _mm256_add_pd(fjz0,tz);
319 /**************************
320 * CALCULATE INTERACTIONS *
321 **************************/
323 /* REACTION-FIELD ELECTROSTATICS */
324 velec = _mm256_mul_pd(qq01,_mm256_sub_pd(_mm256_add_pd(rinv01,_mm256_mul_pd(krf,rsq01)),crf));
325 felec = _mm256_mul_pd(qq01,_mm256_sub_pd(_mm256_mul_pd(rinv01,rinvsq01),krf2));
327 /* Update potential sum for this i atom from the interaction with this j atom. */
328 velecsum = _mm256_add_pd(velecsum,velec);
332 /* Calculate temporary vectorial force */
333 tx = _mm256_mul_pd(fscal,dx01);
334 ty = _mm256_mul_pd(fscal,dy01);
335 tz = _mm256_mul_pd(fscal,dz01);
337 /* Update vectorial force */
338 fix0 = _mm256_add_pd(fix0,tx);
339 fiy0 = _mm256_add_pd(fiy0,ty);
340 fiz0 = _mm256_add_pd(fiz0,tz);
342 fjx1 = _mm256_add_pd(fjx1,tx);
343 fjy1 = _mm256_add_pd(fjy1,ty);
344 fjz1 = _mm256_add_pd(fjz1,tz);
346 /**************************
347 * CALCULATE INTERACTIONS *
348 **************************/
350 /* REACTION-FIELD ELECTROSTATICS */
351 velec = _mm256_mul_pd(qq02,_mm256_sub_pd(_mm256_add_pd(rinv02,_mm256_mul_pd(krf,rsq02)),crf));
352 felec = _mm256_mul_pd(qq02,_mm256_sub_pd(_mm256_mul_pd(rinv02,rinvsq02),krf2));
354 /* Update potential sum for this i atom from the interaction with this j atom. */
355 velecsum = _mm256_add_pd(velecsum,velec);
359 /* Calculate temporary vectorial force */
360 tx = _mm256_mul_pd(fscal,dx02);
361 ty = _mm256_mul_pd(fscal,dy02);
362 tz = _mm256_mul_pd(fscal,dz02);
364 /* Update vectorial force */
365 fix0 = _mm256_add_pd(fix0,tx);
366 fiy0 = _mm256_add_pd(fiy0,ty);
367 fiz0 = _mm256_add_pd(fiz0,tz);
369 fjx2 = _mm256_add_pd(fjx2,tx);
370 fjy2 = _mm256_add_pd(fjy2,ty);
371 fjz2 = _mm256_add_pd(fjz2,tz);
373 /**************************
374 * CALCULATE INTERACTIONS *
375 **************************/
377 /* REACTION-FIELD ELECTROSTATICS */
378 velec = _mm256_mul_pd(qq10,_mm256_sub_pd(_mm256_add_pd(rinv10,_mm256_mul_pd(krf,rsq10)),crf));
379 felec = _mm256_mul_pd(qq10,_mm256_sub_pd(_mm256_mul_pd(rinv10,rinvsq10),krf2));
381 /* Update potential sum for this i atom from the interaction with this j atom. */
382 velecsum = _mm256_add_pd(velecsum,velec);
386 /* Calculate temporary vectorial force */
387 tx = _mm256_mul_pd(fscal,dx10);
388 ty = _mm256_mul_pd(fscal,dy10);
389 tz = _mm256_mul_pd(fscal,dz10);
391 /* Update vectorial force */
392 fix1 = _mm256_add_pd(fix1,tx);
393 fiy1 = _mm256_add_pd(fiy1,ty);
394 fiz1 = _mm256_add_pd(fiz1,tz);
396 fjx0 = _mm256_add_pd(fjx0,tx);
397 fjy0 = _mm256_add_pd(fjy0,ty);
398 fjz0 = _mm256_add_pd(fjz0,tz);
400 /**************************
401 * CALCULATE INTERACTIONS *
402 **************************/
404 /* REACTION-FIELD ELECTROSTATICS */
405 velec = _mm256_mul_pd(qq11,_mm256_sub_pd(_mm256_add_pd(rinv11,_mm256_mul_pd(krf,rsq11)),crf));
406 felec = _mm256_mul_pd(qq11,_mm256_sub_pd(_mm256_mul_pd(rinv11,rinvsq11),krf2));
408 /* Update potential sum for this i atom from the interaction with this j atom. */
409 velecsum = _mm256_add_pd(velecsum,velec);
413 /* Calculate temporary vectorial force */
414 tx = _mm256_mul_pd(fscal,dx11);
415 ty = _mm256_mul_pd(fscal,dy11);
416 tz = _mm256_mul_pd(fscal,dz11);
418 /* Update vectorial force */
419 fix1 = _mm256_add_pd(fix1,tx);
420 fiy1 = _mm256_add_pd(fiy1,ty);
421 fiz1 = _mm256_add_pd(fiz1,tz);
423 fjx1 = _mm256_add_pd(fjx1,tx);
424 fjy1 = _mm256_add_pd(fjy1,ty);
425 fjz1 = _mm256_add_pd(fjz1,tz);
427 /**************************
428 * CALCULATE INTERACTIONS *
429 **************************/
431 /* REACTION-FIELD ELECTROSTATICS */
432 velec = _mm256_mul_pd(qq12,_mm256_sub_pd(_mm256_add_pd(rinv12,_mm256_mul_pd(krf,rsq12)),crf));
433 felec = _mm256_mul_pd(qq12,_mm256_sub_pd(_mm256_mul_pd(rinv12,rinvsq12),krf2));
435 /* Update potential sum for this i atom from the interaction with this j atom. */
436 velecsum = _mm256_add_pd(velecsum,velec);
440 /* Calculate temporary vectorial force */
441 tx = _mm256_mul_pd(fscal,dx12);
442 ty = _mm256_mul_pd(fscal,dy12);
443 tz = _mm256_mul_pd(fscal,dz12);
445 /* Update vectorial force */
446 fix1 = _mm256_add_pd(fix1,tx);
447 fiy1 = _mm256_add_pd(fiy1,ty);
448 fiz1 = _mm256_add_pd(fiz1,tz);
450 fjx2 = _mm256_add_pd(fjx2,tx);
451 fjy2 = _mm256_add_pd(fjy2,ty);
452 fjz2 = _mm256_add_pd(fjz2,tz);
454 /**************************
455 * CALCULATE INTERACTIONS *
456 **************************/
458 /* REACTION-FIELD ELECTROSTATICS */
459 velec = _mm256_mul_pd(qq20,_mm256_sub_pd(_mm256_add_pd(rinv20,_mm256_mul_pd(krf,rsq20)),crf));
460 felec = _mm256_mul_pd(qq20,_mm256_sub_pd(_mm256_mul_pd(rinv20,rinvsq20),krf2));
462 /* Update potential sum for this i atom from the interaction with this j atom. */
463 velecsum = _mm256_add_pd(velecsum,velec);
467 /* Calculate temporary vectorial force */
468 tx = _mm256_mul_pd(fscal,dx20);
469 ty = _mm256_mul_pd(fscal,dy20);
470 tz = _mm256_mul_pd(fscal,dz20);
472 /* Update vectorial force */
473 fix2 = _mm256_add_pd(fix2,tx);
474 fiy2 = _mm256_add_pd(fiy2,ty);
475 fiz2 = _mm256_add_pd(fiz2,tz);
477 fjx0 = _mm256_add_pd(fjx0,tx);
478 fjy0 = _mm256_add_pd(fjy0,ty);
479 fjz0 = _mm256_add_pd(fjz0,tz);
481 /**************************
482 * CALCULATE INTERACTIONS *
483 **************************/
485 /* REACTION-FIELD ELECTROSTATICS */
486 velec = _mm256_mul_pd(qq21,_mm256_sub_pd(_mm256_add_pd(rinv21,_mm256_mul_pd(krf,rsq21)),crf));
487 felec = _mm256_mul_pd(qq21,_mm256_sub_pd(_mm256_mul_pd(rinv21,rinvsq21),krf2));
489 /* Update potential sum for this i atom from the interaction with this j atom. */
490 velecsum = _mm256_add_pd(velecsum,velec);
494 /* Calculate temporary vectorial force */
495 tx = _mm256_mul_pd(fscal,dx21);
496 ty = _mm256_mul_pd(fscal,dy21);
497 tz = _mm256_mul_pd(fscal,dz21);
499 /* Update vectorial force */
500 fix2 = _mm256_add_pd(fix2,tx);
501 fiy2 = _mm256_add_pd(fiy2,ty);
502 fiz2 = _mm256_add_pd(fiz2,tz);
504 fjx1 = _mm256_add_pd(fjx1,tx);
505 fjy1 = _mm256_add_pd(fjy1,ty);
506 fjz1 = _mm256_add_pd(fjz1,tz);
508 /**************************
509 * CALCULATE INTERACTIONS *
510 **************************/
512 /* REACTION-FIELD ELECTROSTATICS */
513 velec = _mm256_mul_pd(qq22,_mm256_sub_pd(_mm256_add_pd(rinv22,_mm256_mul_pd(krf,rsq22)),crf));
514 felec = _mm256_mul_pd(qq22,_mm256_sub_pd(_mm256_mul_pd(rinv22,rinvsq22),krf2));
516 /* Update potential sum for this i atom from the interaction with this j atom. */
517 velecsum = _mm256_add_pd(velecsum,velec);
521 /* Calculate temporary vectorial force */
522 tx = _mm256_mul_pd(fscal,dx22);
523 ty = _mm256_mul_pd(fscal,dy22);
524 tz = _mm256_mul_pd(fscal,dz22);
526 /* Update vectorial force */
527 fix2 = _mm256_add_pd(fix2,tx);
528 fiy2 = _mm256_add_pd(fiy2,ty);
529 fiz2 = _mm256_add_pd(fiz2,tz);
531 fjx2 = _mm256_add_pd(fjx2,tx);
532 fjy2 = _mm256_add_pd(fjy2,ty);
533 fjz2 = _mm256_add_pd(fjz2,tz);
535 fjptrA = f+j_coord_offsetA;
536 fjptrB = f+j_coord_offsetB;
537 fjptrC = f+j_coord_offsetC;
538 fjptrD = f+j_coord_offsetD;
540 gmx_mm256_decrement_3rvec_4ptr_swizzle_pd(fjptrA,fjptrB,fjptrC,fjptrD,
541 fjx0,fjy0,fjz0,fjx1,fjy1,fjz1,fjx2,fjy2,fjz2);
543 /* Inner loop uses 300 flops */
549 /* Get j neighbor index, and coordinate index */
550 jnrlistA = jjnr[jidx];
551 jnrlistB = jjnr[jidx+1];
552 jnrlistC = jjnr[jidx+2];
553 jnrlistD = jjnr[jidx+3];
554 /* Sign of each element will be negative for non-real atoms.
555 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
556 * so use it as val = _mm_andnot_pd(mask,val) to clear dummy entries.
558 tmpmask0 = gmx_mm_castsi128_pd(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
560 tmpmask1 = _mm_permute_ps(tmpmask0,_GMX_MM_PERMUTE(3,3,2,2));
561 tmpmask0 = _mm_permute_ps(tmpmask0,_GMX_MM_PERMUTE(1,1,0,0));
562 dummy_mask = _mm256_castps_pd(gmx_mm256_set_m128(tmpmask1,tmpmask0));
564 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
565 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
566 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
567 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
568 j_coord_offsetA = DIM*jnrA;
569 j_coord_offsetB = DIM*jnrB;
570 j_coord_offsetC = DIM*jnrC;
571 j_coord_offsetD = DIM*jnrD;
573 /* load j atom coordinates */
574 gmx_mm256_load_3rvec_4ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
575 x+j_coord_offsetC,x+j_coord_offsetD,
576 &jx0,&jy0,&jz0,&jx1,&jy1,&jz1,&jx2,&jy2,&jz2);
578 /* Calculate displacement vector */
579 dx00 = _mm256_sub_pd(ix0,jx0);
580 dy00 = _mm256_sub_pd(iy0,jy0);
581 dz00 = _mm256_sub_pd(iz0,jz0);
582 dx01 = _mm256_sub_pd(ix0,jx1);
583 dy01 = _mm256_sub_pd(iy0,jy1);
584 dz01 = _mm256_sub_pd(iz0,jz1);
585 dx02 = _mm256_sub_pd(ix0,jx2);
586 dy02 = _mm256_sub_pd(iy0,jy2);
587 dz02 = _mm256_sub_pd(iz0,jz2);
588 dx10 = _mm256_sub_pd(ix1,jx0);
589 dy10 = _mm256_sub_pd(iy1,jy0);
590 dz10 = _mm256_sub_pd(iz1,jz0);
591 dx11 = _mm256_sub_pd(ix1,jx1);
592 dy11 = _mm256_sub_pd(iy1,jy1);
593 dz11 = _mm256_sub_pd(iz1,jz1);
594 dx12 = _mm256_sub_pd(ix1,jx2);
595 dy12 = _mm256_sub_pd(iy1,jy2);
596 dz12 = _mm256_sub_pd(iz1,jz2);
597 dx20 = _mm256_sub_pd(ix2,jx0);
598 dy20 = _mm256_sub_pd(iy2,jy0);
599 dz20 = _mm256_sub_pd(iz2,jz0);
600 dx21 = _mm256_sub_pd(ix2,jx1);
601 dy21 = _mm256_sub_pd(iy2,jy1);
602 dz21 = _mm256_sub_pd(iz2,jz1);
603 dx22 = _mm256_sub_pd(ix2,jx2);
604 dy22 = _mm256_sub_pd(iy2,jy2);
605 dz22 = _mm256_sub_pd(iz2,jz2);
607 /* Calculate squared distance and things based on it */
608 rsq00 = gmx_mm256_calc_rsq_pd(dx00,dy00,dz00);
609 rsq01 = gmx_mm256_calc_rsq_pd(dx01,dy01,dz01);
610 rsq02 = gmx_mm256_calc_rsq_pd(dx02,dy02,dz02);
611 rsq10 = gmx_mm256_calc_rsq_pd(dx10,dy10,dz10);
612 rsq11 = gmx_mm256_calc_rsq_pd(dx11,dy11,dz11);
613 rsq12 = gmx_mm256_calc_rsq_pd(dx12,dy12,dz12);
614 rsq20 = gmx_mm256_calc_rsq_pd(dx20,dy20,dz20);
615 rsq21 = gmx_mm256_calc_rsq_pd(dx21,dy21,dz21);
616 rsq22 = gmx_mm256_calc_rsq_pd(dx22,dy22,dz22);
618 rinv00 = gmx_mm256_invsqrt_pd(rsq00);
619 rinv01 = gmx_mm256_invsqrt_pd(rsq01);
620 rinv02 = gmx_mm256_invsqrt_pd(rsq02);
621 rinv10 = gmx_mm256_invsqrt_pd(rsq10);
622 rinv11 = gmx_mm256_invsqrt_pd(rsq11);
623 rinv12 = gmx_mm256_invsqrt_pd(rsq12);
624 rinv20 = gmx_mm256_invsqrt_pd(rsq20);
625 rinv21 = gmx_mm256_invsqrt_pd(rsq21);
626 rinv22 = gmx_mm256_invsqrt_pd(rsq22);
628 rinvsq00 = _mm256_mul_pd(rinv00,rinv00);
629 rinvsq01 = _mm256_mul_pd(rinv01,rinv01);
630 rinvsq02 = _mm256_mul_pd(rinv02,rinv02);
631 rinvsq10 = _mm256_mul_pd(rinv10,rinv10);
632 rinvsq11 = _mm256_mul_pd(rinv11,rinv11);
633 rinvsq12 = _mm256_mul_pd(rinv12,rinv12);
634 rinvsq20 = _mm256_mul_pd(rinv20,rinv20);
635 rinvsq21 = _mm256_mul_pd(rinv21,rinv21);
636 rinvsq22 = _mm256_mul_pd(rinv22,rinv22);
638 fjx0 = _mm256_setzero_pd();
639 fjy0 = _mm256_setzero_pd();
640 fjz0 = _mm256_setzero_pd();
641 fjx1 = _mm256_setzero_pd();
642 fjy1 = _mm256_setzero_pd();
643 fjz1 = _mm256_setzero_pd();
644 fjx2 = _mm256_setzero_pd();
645 fjy2 = _mm256_setzero_pd();
646 fjz2 = _mm256_setzero_pd();
648 /**************************
649 * CALCULATE INTERACTIONS *
650 **************************/
652 /* REACTION-FIELD ELECTROSTATICS */
653 velec = _mm256_mul_pd(qq00,_mm256_sub_pd(_mm256_add_pd(rinv00,_mm256_mul_pd(krf,rsq00)),crf));
654 felec = _mm256_mul_pd(qq00,_mm256_sub_pd(_mm256_mul_pd(rinv00,rinvsq00),krf2));
656 /* LENNARD-JONES DISPERSION/REPULSION */
658 rinvsix = _mm256_mul_pd(_mm256_mul_pd(rinvsq00,rinvsq00),rinvsq00);
659 vvdw6 = _mm256_mul_pd(c6_00,rinvsix);
660 vvdw12 = _mm256_mul_pd(c12_00,_mm256_mul_pd(rinvsix,rinvsix));
661 vvdw = _mm256_sub_pd( _mm256_mul_pd(vvdw12,one_twelfth) , _mm256_mul_pd(vvdw6,one_sixth) );
662 fvdw = _mm256_mul_pd(_mm256_sub_pd(vvdw12,vvdw6),rinvsq00);
664 /* Update potential sum for this i atom from the interaction with this j atom. */
665 velec = _mm256_andnot_pd(dummy_mask,velec);
666 velecsum = _mm256_add_pd(velecsum,velec);
667 vvdw = _mm256_andnot_pd(dummy_mask,vvdw);
668 vvdwsum = _mm256_add_pd(vvdwsum,vvdw);
670 fscal = _mm256_add_pd(felec,fvdw);
672 fscal = _mm256_andnot_pd(dummy_mask,fscal);
674 /* Calculate temporary vectorial force */
675 tx = _mm256_mul_pd(fscal,dx00);
676 ty = _mm256_mul_pd(fscal,dy00);
677 tz = _mm256_mul_pd(fscal,dz00);
679 /* Update vectorial force */
680 fix0 = _mm256_add_pd(fix0,tx);
681 fiy0 = _mm256_add_pd(fiy0,ty);
682 fiz0 = _mm256_add_pd(fiz0,tz);
684 fjx0 = _mm256_add_pd(fjx0,tx);
685 fjy0 = _mm256_add_pd(fjy0,ty);
686 fjz0 = _mm256_add_pd(fjz0,tz);
688 /**************************
689 * CALCULATE INTERACTIONS *
690 **************************/
692 /* REACTION-FIELD ELECTROSTATICS */
693 velec = _mm256_mul_pd(qq01,_mm256_sub_pd(_mm256_add_pd(rinv01,_mm256_mul_pd(krf,rsq01)),crf));
694 felec = _mm256_mul_pd(qq01,_mm256_sub_pd(_mm256_mul_pd(rinv01,rinvsq01),krf2));
696 /* Update potential sum for this i atom from the interaction with this j atom. */
697 velec = _mm256_andnot_pd(dummy_mask,velec);
698 velecsum = _mm256_add_pd(velecsum,velec);
702 fscal = _mm256_andnot_pd(dummy_mask,fscal);
704 /* Calculate temporary vectorial force */
705 tx = _mm256_mul_pd(fscal,dx01);
706 ty = _mm256_mul_pd(fscal,dy01);
707 tz = _mm256_mul_pd(fscal,dz01);
709 /* Update vectorial force */
710 fix0 = _mm256_add_pd(fix0,tx);
711 fiy0 = _mm256_add_pd(fiy0,ty);
712 fiz0 = _mm256_add_pd(fiz0,tz);
714 fjx1 = _mm256_add_pd(fjx1,tx);
715 fjy1 = _mm256_add_pd(fjy1,ty);
716 fjz1 = _mm256_add_pd(fjz1,tz);
718 /**************************
719 * CALCULATE INTERACTIONS *
720 **************************/
722 /* REACTION-FIELD ELECTROSTATICS */
723 velec = _mm256_mul_pd(qq02,_mm256_sub_pd(_mm256_add_pd(rinv02,_mm256_mul_pd(krf,rsq02)),crf));
724 felec = _mm256_mul_pd(qq02,_mm256_sub_pd(_mm256_mul_pd(rinv02,rinvsq02),krf2));
726 /* Update potential sum for this i atom from the interaction with this j atom. */
727 velec = _mm256_andnot_pd(dummy_mask,velec);
728 velecsum = _mm256_add_pd(velecsum,velec);
732 fscal = _mm256_andnot_pd(dummy_mask,fscal);
734 /* Calculate temporary vectorial force */
735 tx = _mm256_mul_pd(fscal,dx02);
736 ty = _mm256_mul_pd(fscal,dy02);
737 tz = _mm256_mul_pd(fscal,dz02);
739 /* Update vectorial force */
740 fix0 = _mm256_add_pd(fix0,tx);
741 fiy0 = _mm256_add_pd(fiy0,ty);
742 fiz0 = _mm256_add_pd(fiz0,tz);
744 fjx2 = _mm256_add_pd(fjx2,tx);
745 fjy2 = _mm256_add_pd(fjy2,ty);
746 fjz2 = _mm256_add_pd(fjz2,tz);
748 /**************************
749 * CALCULATE INTERACTIONS *
750 **************************/
752 /* REACTION-FIELD ELECTROSTATICS */
753 velec = _mm256_mul_pd(qq10,_mm256_sub_pd(_mm256_add_pd(rinv10,_mm256_mul_pd(krf,rsq10)),crf));
754 felec = _mm256_mul_pd(qq10,_mm256_sub_pd(_mm256_mul_pd(rinv10,rinvsq10),krf2));
756 /* Update potential sum for this i atom from the interaction with this j atom. */
757 velec = _mm256_andnot_pd(dummy_mask,velec);
758 velecsum = _mm256_add_pd(velecsum,velec);
762 fscal = _mm256_andnot_pd(dummy_mask,fscal);
764 /* Calculate temporary vectorial force */
765 tx = _mm256_mul_pd(fscal,dx10);
766 ty = _mm256_mul_pd(fscal,dy10);
767 tz = _mm256_mul_pd(fscal,dz10);
769 /* Update vectorial force */
770 fix1 = _mm256_add_pd(fix1,tx);
771 fiy1 = _mm256_add_pd(fiy1,ty);
772 fiz1 = _mm256_add_pd(fiz1,tz);
774 fjx0 = _mm256_add_pd(fjx0,tx);
775 fjy0 = _mm256_add_pd(fjy0,ty);
776 fjz0 = _mm256_add_pd(fjz0,tz);
778 /**************************
779 * CALCULATE INTERACTIONS *
780 **************************/
782 /* REACTION-FIELD ELECTROSTATICS */
783 velec = _mm256_mul_pd(qq11,_mm256_sub_pd(_mm256_add_pd(rinv11,_mm256_mul_pd(krf,rsq11)),crf));
784 felec = _mm256_mul_pd(qq11,_mm256_sub_pd(_mm256_mul_pd(rinv11,rinvsq11),krf2));
786 /* Update potential sum for this i atom from the interaction with this j atom. */
787 velec = _mm256_andnot_pd(dummy_mask,velec);
788 velecsum = _mm256_add_pd(velecsum,velec);
792 fscal = _mm256_andnot_pd(dummy_mask,fscal);
794 /* Calculate temporary vectorial force */
795 tx = _mm256_mul_pd(fscal,dx11);
796 ty = _mm256_mul_pd(fscal,dy11);
797 tz = _mm256_mul_pd(fscal,dz11);
799 /* Update vectorial force */
800 fix1 = _mm256_add_pd(fix1,tx);
801 fiy1 = _mm256_add_pd(fiy1,ty);
802 fiz1 = _mm256_add_pd(fiz1,tz);
804 fjx1 = _mm256_add_pd(fjx1,tx);
805 fjy1 = _mm256_add_pd(fjy1,ty);
806 fjz1 = _mm256_add_pd(fjz1,tz);
808 /**************************
809 * CALCULATE INTERACTIONS *
810 **************************/
812 /* REACTION-FIELD ELECTROSTATICS */
813 velec = _mm256_mul_pd(qq12,_mm256_sub_pd(_mm256_add_pd(rinv12,_mm256_mul_pd(krf,rsq12)),crf));
814 felec = _mm256_mul_pd(qq12,_mm256_sub_pd(_mm256_mul_pd(rinv12,rinvsq12),krf2));
816 /* Update potential sum for this i atom from the interaction with this j atom. */
817 velec = _mm256_andnot_pd(dummy_mask,velec);
818 velecsum = _mm256_add_pd(velecsum,velec);
822 fscal = _mm256_andnot_pd(dummy_mask,fscal);
824 /* Calculate temporary vectorial force */
825 tx = _mm256_mul_pd(fscal,dx12);
826 ty = _mm256_mul_pd(fscal,dy12);
827 tz = _mm256_mul_pd(fscal,dz12);
829 /* Update vectorial force */
830 fix1 = _mm256_add_pd(fix1,tx);
831 fiy1 = _mm256_add_pd(fiy1,ty);
832 fiz1 = _mm256_add_pd(fiz1,tz);
834 fjx2 = _mm256_add_pd(fjx2,tx);
835 fjy2 = _mm256_add_pd(fjy2,ty);
836 fjz2 = _mm256_add_pd(fjz2,tz);
838 /**************************
839 * CALCULATE INTERACTIONS *
840 **************************/
842 /* REACTION-FIELD ELECTROSTATICS */
843 velec = _mm256_mul_pd(qq20,_mm256_sub_pd(_mm256_add_pd(rinv20,_mm256_mul_pd(krf,rsq20)),crf));
844 felec = _mm256_mul_pd(qq20,_mm256_sub_pd(_mm256_mul_pd(rinv20,rinvsq20),krf2));
846 /* Update potential sum for this i atom from the interaction with this j atom. */
847 velec = _mm256_andnot_pd(dummy_mask,velec);
848 velecsum = _mm256_add_pd(velecsum,velec);
852 fscal = _mm256_andnot_pd(dummy_mask,fscal);
854 /* Calculate temporary vectorial force */
855 tx = _mm256_mul_pd(fscal,dx20);
856 ty = _mm256_mul_pd(fscal,dy20);
857 tz = _mm256_mul_pd(fscal,dz20);
859 /* Update vectorial force */
860 fix2 = _mm256_add_pd(fix2,tx);
861 fiy2 = _mm256_add_pd(fiy2,ty);
862 fiz2 = _mm256_add_pd(fiz2,tz);
864 fjx0 = _mm256_add_pd(fjx0,tx);
865 fjy0 = _mm256_add_pd(fjy0,ty);
866 fjz0 = _mm256_add_pd(fjz0,tz);
868 /**************************
869 * CALCULATE INTERACTIONS *
870 **************************/
872 /* REACTION-FIELD ELECTROSTATICS */
873 velec = _mm256_mul_pd(qq21,_mm256_sub_pd(_mm256_add_pd(rinv21,_mm256_mul_pd(krf,rsq21)),crf));
874 felec = _mm256_mul_pd(qq21,_mm256_sub_pd(_mm256_mul_pd(rinv21,rinvsq21),krf2));
876 /* Update potential sum for this i atom from the interaction with this j atom. */
877 velec = _mm256_andnot_pd(dummy_mask,velec);
878 velecsum = _mm256_add_pd(velecsum,velec);
882 fscal = _mm256_andnot_pd(dummy_mask,fscal);
884 /* Calculate temporary vectorial force */
885 tx = _mm256_mul_pd(fscal,dx21);
886 ty = _mm256_mul_pd(fscal,dy21);
887 tz = _mm256_mul_pd(fscal,dz21);
889 /* Update vectorial force */
890 fix2 = _mm256_add_pd(fix2,tx);
891 fiy2 = _mm256_add_pd(fiy2,ty);
892 fiz2 = _mm256_add_pd(fiz2,tz);
894 fjx1 = _mm256_add_pd(fjx1,tx);
895 fjy1 = _mm256_add_pd(fjy1,ty);
896 fjz1 = _mm256_add_pd(fjz1,tz);
898 /**************************
899 * CALCULATE INTERACTIONS *
900 **************************/
902 /* REACTION-FIELD ELECTROSTATICS */
903 velec = _mm256_mul_pd(qq22,_mm256_sub_pd(_mm256_add_pd(rinv22,_mm256_mul_pd(krf,rsq22)),crf));
904 felec = _mm256_mul_pd(qq22,_mm256_sub_pd(_mm256_mul_pd(rinv22,rinvsq22),krf2));
906 /* Update potential sum for this i atom from the interaction with this j atom. */
907 velec = _mm256_andnot_pd(dummy_mask,velec);
908 velecsum = _mm256_add_pd(velecsum,velec);
912 fscal = _mm256_andnot_pd(dummy_mask,fscal);
914 /* Calculate temporary vectorial force */
915 tx = _mm256_mul_pd(fscal,dx22);
916 ty = _mm256_mul_pd(fscal,dy22);
917 tz = _mm256_mul_pd(fscal,dz22);
919 /* Update vectorial force */
920 fix2 = _mm256_add_pd(fix2,tx);
921 fiy2 = _mm256_add_pd(fiy2,ty);
922 fiz2 = _mm256_add_pd(fiz2,tz);
924 fjx2 = _mm256_add_pd(fjx2,tx);
925 fjy2 = _mm256_add_pd(fjy2,ty);
926 fjz2 = _mm256_add_pd(fjz2,tz);
928 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
929 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
930 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
931 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
933 gmx_mm256_decrement_3rvec_4ptr_swizzle_pd(fjptrA,fjptrB,fjptrC,fjptrD,
934 fjx0,fjy0,fjz0,fjx1,fjy1,fjz1,fjx2,fjy2,fjz2);
936 /* Inner loop uses 300 flops */
939 /* End of innermost loop */
941 gmx_mm256_update_iforce_3atom_swizzle_pd(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,
942 f+i_coord_offset,fshift+i_shift_offset);
945 /* Update potential energies */
946 gmx_mm256_update_1pot_pd(velecsum,kernel_data->energygrp_elec+ggid);
947 gmx_mm256_update_1pot_pd(vvdwsum,kernel_data->energygrp_vdw+ggid);
949 /* Increment number of inner iterations */
950 inneriter += j_index_end - j_index_start;
952 /* Outer loop uses 20 flops */
955 /* Increment number of outer iterations */
958 /* Update outer/inner flops */
960 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W3W3_VF,outeriter*20 + inneriter*300);
963 * Gromacs nonbonded kernel: nb_kernel_ElecRF_VdwLJ_GeomW3W3_F_avx_256_double
964 * Electrostatics interaction: ReactionField
965 * VdW interaction: LennardJones
966 * Geometry: Water3-Water3
967 * Calculate force/pot: Force
970 nb_kernel_ElecRF_VdwLJ_GeomW3W3_F_avx_256_double
971 (t_nblist * gmx_restrict nlist,
972 rvec * gmx_restrict xx,
973 rvec * gmx_restrict ff,
974 t_forcerec * gmx_restrict fr,
975 t_mdatoms * gmx_restrict mdatoms,
976 nb_kernel_data_t * gmx_restrict kernel_data,
977 t_nrnb * gmx_restrict nrnb)
979 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
980 * just 0 for non-waters.
981 * Suffixes A,B,C,D refer to j loop unrolling done with AVX, e.g. for the four different
982 * jnr indices corresponding to data put in the four positions in the SIMD register.
984 int i_shift_offset,i_coord_offset,outeriter,inneriter;
985 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
986 int jnrA,jnrB,jnrC,jnrD;
987 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
988 int jnrlistE,jnrlistF,jnrlistG,jnrlistH;
989 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
990 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
992 real *shiftvec,*fshift,*x,*f;
993 real *fjptrA,*fjptrB,*fjptrC,*fjptrD;
995 __m256d tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
996 real * vdwioffsetptr0;
997 __m256d ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
998 real * vdwioffsetptr1;
999 __m256d ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
1000 real * vdwioffsetptr2;
1001 __m256d ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
1002 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
1003 __m256d jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
1004 int vdwjidx1A,vdwjidx1B,vdwjidx1C,vdwjidx1D;
1005 __m256d jx1,jy1,jz1,fjx1,fjy1,fjz1,jq1,isaj1;
1006 int vdwjidx2A,vdwjidx2B,vdwjidx2C,vdwjidx2D;
1007 __m256d jx2,jy2,jz2,fjx2,fjy2,fjz2,jq2,isaj2;
1008 __m256d dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
1009 __m256d dx01,dy01,dz01,rsq01,rinv01,rinvsq01,r01,qq01,c6_01,c12_01;
1010 __m256d dx02,dy02,dz02,rsq02,rinv02,rinvsq02,r02,qq02,c6_02,c12_02;
1011 __m256d dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
1012 __m256d dx11,dy11,dz11,rsq11,rinv11,rinvsq11,r11,qq11,c6_11,c12_11;
1013 __m256d dx12,dy12,dz12,rsq12,rinv12,rinvsq12,r12,qq12,c6_12,c12_12;
1014 __m256d dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
1015 __m256d dx21,dy21,dz21,rsq21,rinv21,rinvsq21,r21,qq21,c6_21,c12_21;
1016 __m256d dx22,dy22,dz22,rsq22,rinv22,rinvsq22,r22,qq22,c6_22,c12_22;
1017 __m256d velec,felec,velecsum,facel,crf,krf,krf2;
1020 __m256d rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
1023 __m256d one_sixth = _mm256_set1_pd(1.0/6.0);
1024 __m256d one_twelfth = _mm256_set1_pd(1.0/12.0);
1025 __m256d dummy_mask,cutoff_mask;
1026 __m128 tmpmask0,tmpmask1;
1027 __m256d signbit = _mm256_castsi256_pd( _mm256_set1_epi32(0x80000000) );
1028 __m256d one = _mm256_set1_pd(1.0);
1029 __m256d two = _mm256_set1_pd(2.0);
1035 jindex = nlist->jindex;
1037 shiftidx = nlist->shift;
1039 shiftvec = fr->shift_vec[0];
1040 fshift = fr->fshift[0];
1041 facel = _mm256_set1_pd(fr->epsfac);
1042 charge = mdatoms->chargeA;
1043 krf = _mm256_set1_pd(fr->ic->k_rf);
1044 krf2 = _mm256_set1_pd(fr->ic->k_rf*2.0);
1045 crf = _mm256_set1_pd(fr->ic->c_rf);
1046 nvdwtype = fr->ntype;
1047 vdwparam = fr->nbfp;
1048 vdwtype = mdatoms->typeA;
1050 /* Setup water-specific parameters */
1051 inr = nlist->iinr[0];
1052 iq0 = _mm256_mul_pd(facel,_mm256_set1_pd(charge[inr+0]));
1053 iq1 = _mm256_mul_pd(facel,_mm256_set1_pd(charge[inr+1]));
1054 iq2 = _mm256_mul_pd(facel,_mm256_set1_pd(charge[inr+2]));
1055 vdwioffsetptr0 = vdwparam+2*nvdwtype*vdwtype[inr+0];
1057 jq0 = _mm256_set1_pd(charge[inr+0]);
1058 jq1 = _mm256_set1_pd(charge[inr+1]);
1059 jq2 = _mm256_set1_pd(charge[inr+2]);
1060 vdwjidx0A = 2*vdwtype[inr+0];
1061 qq00 = _mm256_mul_pd(iq0,jq0);
1062 c6_00 = _mm256_set1_pd(vdwioffsetptr0[vdwjidx0A]);
1063 c12_00 = _mm256_set1_pd(vdwioffsetptr0[vdwjidx0A+1]);
1064 qq01 = _mm256_mul_pd(iq0,jq1);
1065 qq02 = _mm256_mul_pd(iq0,jq2);
1066 qq10 = _mm256_mul_pd(iq1,jq0);
1067 qq11 = _mm256_mul_pd(iq1,jq1);
1068 qq12 = _mm256_mul_pd(iq1,jq2);
1069 qq20 = _mm256_mul_pd(iq2,jq0);
1070 qq21 = _mm256_mul_pd(iq2,jq1);
1071 qq22 = _mm256_mul_pd(iq2,jq2);
1073 /* Avoid stupid compiler warnings */
1074 jnrA = jnrB = jnrC = jnrD = 0;
1075 j_coord_offsetA = 0;
1076 j_coord_offsetB = 0;
1077 j_coord_offsetC = 0;
1078 j_coord_offsetD = 0;
1083 for(iidx=0;iidx<4*DIM;iidx++)
1085 scratch[iidx] = 0.0;
1088 /* Start outer loop over neighborlists */
1089 for(iidx=0; iidx<nri; iidx++)
1091 /* Load shift vector for this list */
1092 i_shift_offset = DIM*shiftidx[iidx];
1094 /* Load limits for loop over neighbors */
1095 j_index_start = jindex[iidx];
1096 j_index_end = jindex[iidx+1];
1098 /* Get outer coordinate index */
1100 i_coord_offset = DIM*inr;
1102 /* Load i particle coords and add shift vector */
1103 gmx_mm256_load_shift_and_3rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,
1104 &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2);
1106 fix0 = _mm256_setzero_pd();
1107 fiy0 = _mm256_setzero_pd();
1108 fiz0 = _mm256_setzero_pd();
1109 fix1 = _mm256_setzero_pd();
1110 fiy1 = _mm256_setzero_pd();
1111 fiz1 = _mm256_setzero_pd();
1112 fix2 = _mm256_setzero_pd();
1113 fiy2 = _mm256_setzero_pd();
1114 fiz2 = _mm256_setzero_pd();
1116 /* Start inner kernel loop */
1117 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
1120 /* Get j neighbor index, and coordinate index */
1122 jnrB = jjnr[jidx+1];
1123 jnrC = jjnr[jidx+2];
1124 jnrD = jjnr[jidx+3];
1125 j_coord_offsetA = DIM*jnrA;
1126 j_coord_offsetB = DIM*jnrB;
1127 j_coord_offsetC = DIM*jnrC;
1128 j_coord_offsetD = DIM*jnrD;
1130 /* load j atom coordinates */
1131 gmx_mm256_load_3rvec_4ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
1132 x+j_coord_offsetC,x+j_coord_offsetD,
1133 &jx0,&jy0,&jz0,&jx1,&jy1,&jz1,&jx2,&jy2,&jz2);
1135 /* Calculate displacement vector */
1136 dx00 = _mm256_sub_pd(ix0,jx0);
1137 dy00 = _mm256_sub_pd(iy0,jy0);
1138 dz00 = _mm256_sub_pd(iz0,jz0);
1139 dx01 = _mm256_sub_pd(ix0,jx1);
1140 dy01 = _mm256_sub_pd(iy0,jy1);
1141 dz01 = _mm256_sub_pd(iz0,jz1);
1142 dx02 = _mm256_sub_pd(ix0,jx2);
1143 dy02 = _mm256_sub_pd(iy0,jy2);
1144 dz02 = _mm256_sub_pd(iz0,jz2);
1145 dx10 = _mm256_sub_pd(ix1,jx0);
1146 dy10 = _mm256_sub_pd(iy1,jy0);
1147 dz10 = _mm256_sub_pd(iz1,jz0);
1148 dx11 = _mm256_sub_pd(ix1,jx1);
1149 dy11 = _mm256_sub_pd(iy1,jy1);
1150 dz11 = _mm256_sub_pd(iz1,jz1);
1151 dx12 = _mm256_sub_pd(ix1,jx2);
1152 dy12 = _mm256_sub_pd(iy1,jy2);
1153 dz12 = _mm256_sub_pd(iz1,jz2);
1154 dx20 = _mm256_sub_pd(ix2,jx0);
1155 dy20 = _mm256_sub_pd(iy2,jy0);
1156 dz20 = _mm256_sub_pd(iz2,jz0);
1157 dx21 = _mm256_sub_pd(ix2,jx1);
1158 dy21 = _mm256_sub_pd(iy2,jy1);
1159 dz21 = _mm256_sub_pd(iz2,jz1);
1160 dx22 = _mm256_sub_pd(ix2,jx2);
1161 dy22 = _mm256_sub_pd(iy2,jy2);
1162 dz22 = _mm256_sub_pd(iz2,jz2);
1164 /* Calculate squared distance and things based on it */
1165 rsq00 = gmx_mm256_calc_rsq_pd(dx00,dy00,dz00);
1166 rsq01 = gmx_mm256_calc_rsq_pd(dx01,dy01,dz01);
1167 rsq02 = gmx_mm256_calc_rsq_pd(dx02,dy02,dz02);
1168 rsq10 = gmx_mm256_calc_rsq_pd(dx10,dy10,dz10);
1169 rsq11 = gmx_mm256_calc_rsq_pd(dx11,dy11,dz11);
1170 rsq12 = gmx_mm256_calc_rsq_pd(dx12,dy12,dz12);
1171 rsq20 = gmx_mm256_calc_rsq_pd(dx20,dy20,dz20);
1172 rsq21 = gmx_mm256_calc_rsq_pd(dx21,dy21,dz21);
1173 rsq22 = gmx_mm256_calc_rsq_pd(dx22,dy22,dz22);
1175 rinv00 = gmx_mm256_invsqrt_pd(rsq00);
1176 rinv01 = gmx_mm256_invsqrt_pd(rsq01);
1177 rinv02 = gmx_mm256_invsqrt_pd(rsq02);
1178 rinv10 = gmx_mm256_invsqrt_pd(rsq10);
1179 rinv11 = gmx_mm256_invsqrt_pd(rsq11);
1180 rinv12 = gmx_mm256_invsqrt_pd(rsq12);
1181 rinv20 = gmx_mm256_invsqrt_pd(rsq20);
1182 rinv21 = gmx_mm256_invsqrt_pd(rsq21);
1183 rinv22 = gmx_mm256_invsqrt_pd(rsq22);
1185 rinvsq00 = _mm256_mul_pd(rinv00,rinv00);
1186 rinvsq01 = _mm256_mul_pd(rinv01,rinv01);
1187 rinvsq02 = _mm256_mul_pd(rinv02,rinv02);
1188 rinvsq10 = _mm256_mul_pd(rinv10,rinv10);
1189 rinvsq11 = _mm256_mul_pd(rinv11,rinv11);
1190 rinvsq12 = _mm256_mul_pd(rinv12,rinv12);
1191 rinvsq20 = _mm256_mul_pd(rinv20,rinv20);
1192 rinvsq21 = _mm256_mul_pd(rinv21,rinv21);
1193 rinvsq22 = _mm256_mul_pd(rinv22,rinv22);
1195 fjx0 = _mm256_setzero_pd();
1196 fjy0 = _mm256_setzero_pd();
1197 fjz0 = _mm256_setzero_pd();
1198 fjx1 = _mm256_setzero_pd();
1199 fjy1 = _mm256_setzero_pd();
1200 fjz1 = _mm256_setzero_pd();
1201 fjx2 = _mm256_setzero_pd();
1202 fjy2 = _mm256_setzero_pd();
1203 fjz2 = _mm256_setzero_pd();
1205 /**************************
1206 * CALCULATE INTERACTIONS *
1207 **************************/
1209 /* REACTION-FIELD ELECTROSTATICS */
1210 felec = _mm256_mul_pd(qq00,_mm256_sub_pd(_mm256_mul_pd(rinv00,rinvsq00),krf2));
1212 /* LENNARD-JONES DISPERSION/REPULSION */
1214 rinvsix = _mm256_mul_pd(_mm256_mul_pd(rinvsq00,rinvsq00),rinvsq00);
1215 fvdw = _mm256_mul_pd(_mm256_sub_pd(_mm256_mul_pd(c12_00,rinvsix),c6_00),_mm256_mul_pd(rinvsix,rinvsq00));
1217 fscal = _mm256_add_pd(felec,fvdw);
1219 /* Calculate temporary vectorial force */
1220 tx = _mm256_mul_pd(fscal,dx00);
1221 ty = _mm256_mul_pd(fscal,dy00);
1222 tz = _mm256_mul_pd(fscal,dz00);
1224 /* Update vectorial force */
1225 fix0 = _mm256_add_pd(fix0,tx);
1226 fiy0 = _mm256_add_pd(fiy0,ty);
1227 fiz0 = _mm256_add_pd(fiz0,tz);
1229 fjx0 = _mm256_add_pd(fjx0,tx);
1230 fjy0 = _mm256_add_pd(fjy0,ty);
1231 fjz0 = _mm256_add_pd(fjz0,tz);
1233 /**************************
1234 * CALCULATE INTERACTIONS *
1235 **************************/
1237 /* REACTION-FIELD ELECTROSTATICS */
1238 felec = _mm256_mul_pd(qq01,_mm256_sub_pd(_mm256_mul_pd(rinv01,rinvsq01),krf2));
1242 /* Calculate temporary vectorial force */
1243 tx = _mm256_mul_pd(fscal,dx01);
1244 ty = _mm256_mul_pd(fscal,dy01);
1245 tz = _mm256_mul_pd(fscal,dz01);
1247 /* Update vectorial force */
1248 fix0 = _mm256_add_pd(fix0,tx);
1249 fiy0 = _mm256_add_pd(fiy0,ty);
1250 fiz0 = _mm256_add_pd(fiz0,tz);
1252 fjx1 = _mm256_add_pd(fjx1,tx);
1253 fjy1 = _mm256_add_pd(fjy1,ty);
1254 fjz1 = _mm256_add_pd(fjz1,tz);
1256 /**************************
1257 * CALCULATE INTERACTIONS *
1258 **************************/
1260 /* REACTION-FIELD ELECTROSTATICS */
1261 felec = _mm256_mul_pd(qq02,_mm256_sub_pd(_mm256_mul_pd(rinv02,rinvsq02),krf2));
1265 /* Calculate temporary vectorial force */
1266 tx = _mm256_mul_pd(fscal,dx02);
1267 ty = _mm256_mul_pd(fscal,dy02);
1268 tz = _mm256_mul_pd(fscal,dz02);
1270 /* Update vectorial force */
1271 fix0 = _mm256_add_pd(fix0,tx);
1272 fiy0 = _mm256_add_pd(fiy0,ty);
1273 fiz0 = _mm256_add_pd(fiz0,tz);
1275 fjx2 = _mm256_add_pd(fjx2,tx);
1276 fjy2 = _mm256_add_pd(fjy2,ty);
1277 fjz2 = _mm256_add_pd(fjz2,tz);
1279 /**************************
1280 * CALCULATE INTERACTIONS *
1281 **************************/
1283 /* REACTION-FIELD ELECTROSTATICS */
1284 felec = _mm256_mul_pd(qq10,_mm256_sub_pd(_mm256_mul_pd(rinv10,rinvsq10),krf2));
1288 /* Calculate temporary vectorial force */
1289 tx = _mm256_mul_pd(fscal,dx10);
1290 ty = _mm256_mul_pd(fscal,dy10);
1291 tz = _mm256_mul_pd(fscal,dz10);
1293 /* Update vectorial force */
1294 fix1 = _mm256_add_pd(fix1,tx);
1295 fiy1 = _mm256_add_pd(fiy1,ty);
1296 fiz1 = _mm256_add_pd(fiz1,tz);
1298 fjx0 = _mm256_add_pd(fjx0,tx);
1299 fjy0 = _mm256_add_pd(fjy0,ty);
1300 fjz0 = _mm256_add_pd(fjz0,tz);
1302 /**************************
1303 * CALCULATE INTERACTIONS *
1304 **************************/
1306 /* REACTION-FIELD ELECTROSTATICS */
1307 felec = _mm256_mul_pd(qq11,_mm256_sub_pd(_mm256_mul_pd(rinv11,rinvsq11),krf2));
1311 /* Calculate temporary vectorial force */
1312 tx = _mm256_mul_pd(fscal,dx11);
1313 ty = _mm256_mul_pd(fscal,dy11);
1314 tz = _mm256_mul_pd(fscal,dz11);
1316 /* Update vectorial force */
1317 fix1 = _mm256_add_pd(fix1,tx);
1318 fiy1 = _mm256_add_pd(fiy1,ty);
1319 fiz1 = _mm256_add_pd(fiz1,tz);
1321 fjx1 = _mm256_add_pd(fjx1,tx);
1322 fjy1 = _mm256_add_pd(fjy1,ty);
1323 fjz1 = _mm256_add_pd(fjz1,tz);
1325 /**************************
1326 * CALCULATE INTERACTIONS *
1327 **************************/
1329 /* REACTION-FIELD ELECTROSTATICS */
1330 felec = _mm256_mul_pd(qq12,_mm256_sub_pd(_mm256_mul_pd(rinv12,rinvsq12),krf2));
1334 /* Calculate temporary vectorial force */
1335 tx = _mm256_mul_pd(fscal,dx12);
1336 ty = _mm256_mul_pd(fscal,dy12);
1337 tz = _mm256_mul_pd(fscal,dz12);
1339 /* Update vectorial force */
1340 fix1 = _mm256_add_pd(fix1,tx);
1341 fiy1 = _mm256_add_pd(fiy1,ty);
1342 fiz1 = _mm256_add_pd(fiz1,tz);
1344 fjx2 = _mm256_add_pd(fjx2,tx);
1345 fjy2 = _mm256_add_pd(fjy2,ty);
1346 fjz2 = _mm256_add_pd(fjz2,tz);
1348 /**************************
1349 * CALCULATE INTERACTIONS *
1350 **************************/
1352 /* REACTION-FIELD ELECTROSTATICS */
1353 felec = _mm256_mul_pd(qq20,_mm256_sub_pd(_mm256_mul_pd(rinv20,rinvsq20),krf2));
1357 /* Calculate temporary vectorial force */
1358 tx = _mm256_mul_pd(fscal,dx20);
1359 ty = _mm256_mul_pd(fscal,dy20);
1360 tz = _mm256_mul_pd(fscal,dz20);
1362 /* Update vectorial force */
1363 fix2 = _mm256_add_pd(fix2,tx);
1364 fiy2 = _mm256_add_pd(fiy2,ty);
1365 fiz2 = _mm256_add_pd(fiz2,tz);
1367 fjx0 = _mm256_add_pd(fjx0,tx);
1368 fjy0 = _mm256_add_pd(fjy0,ty);
1369 fjz0 = _mm256_add_pd(fjz0,tz);
1371 /**************************
1372 * CALCULATE INTERACTIONS *
1373 **************************/
1375 /* REACTION-FIELD ELECTROSTATICS */
1376 felec = _mm256_mul_pd(qq21,_mm256_sub_pd(_mm256_mul_pd(rinv21,rinvsq21),krf2));
1380 /* Calculate temporary vectorial force */
1381 tx = _mm256_mul_pd(fscal,dx21);
1382 ty = _mm256_mul_pd(fscal,dy21);
1383 tz = _mm256_mul_pd(fscal,dz21);
1385 /* Update vectorial force */
1386 fix2 = _mm256_add_pd(fix2,tx);
1387 fiy2 = _mm256_add_pd(fiy2,ty);
1388 fiz2 = _mm256_add_pd(fiz2,tz);
1390 fjx1 = _mm256_add_pd(fjx1,tx);
1391 fjy1 = _mm256_add_pd(fjy1,ty);
1392 fjz1 = _mm256_add_pd(fjz1,tz);
1394 /**************************
1395 * CALCULATE INTERACTIONS *
1396 **************************/
1398 /* REACTION-FIELD ELECTROSTATICS */
1399 felec = _mm256_mul_pd(qq22,_mm256_sub_pd(_mm256_mul_pd(rinv22,rinvsq22),krf2));
1403 /* Calculate temporary vectorial force */
1404 tx = _mm256_mul_pd(fscal,dx22);
1405 ty = _mm256_mul_pd(fscal,dy22);
1406 tz = _mm256_mul_pd(fscal,dz22);
1408 /* Update vectorial force */
1409 fix2 = _mm256_add_pd(fix2,tx);
1410 fiy2 = _mm256_add_pd(fiy2,ty);
1411 fiz2 = _mm256_add_pd(fiz2,tz);
1413 fjx2 = _mm256_add_pd(fjx2,tx);
1414 fjy2 = _mm256_add_pd(fjy2,ty);
1415 fjz2 = _mm256_add_pd(fjz2,tz);
1417 fjptrA = f+j_coord_offsetA;
1418 fjptrB = f+j_coord_offsetB;
1419 fjptrC = f+j_coord_offsetC;
1420 fjptrD = f+j_coord_offsetD;
1422 gmx_mm256_decrement_3rvec_4ptr_swizzle_pd(fjptrA,fjptrB,fjptrC,fjptrD,
1423 fjx0,fjy0,fjz0,fjx1,fjy1,fjz1,fjx2,fjy2,fjz2);
1425 /* Inner loop uses 250 flops */
1428 if(jidx<j_index_end)
1431 /* Get j neighbor index, and coordinate index */
1432 jnrlistA = jjnr[jidx];
1433 jnrlistB = jjnr[jidx+1];
1434 jnrlistC = jjnr[jidx+2];
1435 jnrlistD = jjnr[jidx+3];
1436 /* Sign of each element will be negative for non-real atoms.
1437 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
1438 * so use it as val = _mm_andnot_pd(mask,val) to clear dummy entries.
1440 tmpmask0 = gmx_mm_castsi128_pd(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
1442 tmpmask1 = _mm_permute_ps(tmpmask0,_GMX_MM_PERMUTE(3,3,2,2));
1443 tmpmask0 = _mm_permute_ps(tmpmask0,_GMX_MM_PERMUTE(1,1,0,0));
1444 dummy_mask = _mm256_castps_pd(gmx_mm256_set_m128(tmpmask1,tmpmask0));
1446 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
1447 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
1448 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
1449 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
1450 j_coord_offsetA = DIM*jnrA;
1451 j_coord_offsetB = DIM*jnrB;
1452 j_coord_offsetC = DIM*jnrC;
1453 j_coord_offsetD = DIM*jnrD;
1455 /* load j atom coordinates */
1456 gmx_mm256_load_3rvec_4ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
1457 x+j_coord_offsetC,x+j_coord_offsetD,
1458 &jx0,&jy0,&jz0,&jx1,&jy1,&jz1,&jx2,&jy2,&jz2);
1460 /* Calculate displacement vector */
1461 dx00 = _mm256_sub_pd(ix0,jx0);
1462 dy00 = _mm256_sub_pd(iy0,jy0);
1463 dz00 = _mm256_sub_pd(iz0,jz0);
1464 dx01 = _mm256_sub_pd(ix0,jx1);
1465 dy01 = _mm256_sub_pd(iy0,jy1);
1466 dz01 = _mm256_sub_pd(iz0,jz1);
1467 dx02 = _mm256_sub_pd(ix0,jx2);
1468 dy02 = _mm256_sub_pd(iy0,jy2);
1469 dz02 = _mm256_sub_pd(iz0,jz2);
1470 dx10 = _mm256_sub_pd(ix1,jx0);
1471 dy10 = _mm256_sub_pd(iy1,jy0);
1472 dz10 = _mm256_sub_pd(iz1,jz0);
1473 dx11 = _mm256_sub_pd(ix1,jx1);
1474 dy11 = _mm256_sub_pd(iy1,jy1);
1475 dz11 = _mm256_sub_pd(iz1,jz1);
1476 dx12 = _mm256_sub_pd(ix1,jx2);
1477 dy12 = _mm256_sub_pd(iy1,jy2);
1478 dz12 = _mm256_sub_pd(iz1,jz2);
1479 dx20 = _mm256_sub_pd(ix2,jx0);
1480 dy20 = _mm256_sub_pd(iy2,jy0);
1481 dz20 = _mm256_sub_pd(iz2,jz0);
1482 dx21 = _mm256_sub_pd(ix2,jx1);
1483 dy21 = _mm256_sub_pd(iy2,jy1);
1484 dz21 = _mm256_sub_pd(iz2,jz1);
1485 dx22 = _mm256_sub_pd(ix2,jx2);
1486 dy22 = _mm256_sub_pd(iy2,jy2);
1487 dz22 = _mm256_sub_pd(iz2,jz2);
1489 /* Calculate squared distance and things based on it */
1490 rsq00 = gmx_mm256_calc_rsq_pd(dx00,dy00,dz00);
1491 rsq01 = gmx_mm256_calc_rsq_pd(dx01,dy01,dz01);
1492 rsq02 = gmx_mm256_calc_rsq_pd(dx02,dy02,dz02);
1493 rsq10 = gmx_mm256_calc_rsq_pd(dx10,dy10,dz10);
1494 rsq11 = gmx_mm256_calc_rsq_pd(dx11,dy11,dz11);
1495 rsq12 = gmx_mm256_calc_rsq_pd(dx12,dy12,dz12);
1496 rsq20 = gmx_mm256_calc_rsq_pd(dx20,dy20,dz20);
1497 rsq21 = gmx_mm256_calc_rsq_pd(dx21,dy21,dz21);
1498 rsq22 = gmx_mm256_calc_rsq_pd(dx22,dy22,dz22);
1500 rinv00 = gmx_mm256_invsqrt_pd(rsq00);
1501 rinv01 = gmx_mm256_invsqrt_pd(rsq01);
1502 rinv02 = gmx_mm256_invsqrt_pd(rsq02);
1503 rinv10 = gmx_mm256_invsqrt_pd(rsq10);
1504 rinv11 = gmx_mm256_invsqrt_pd(rsq11);
1505 rinv12 = gmx_mm256_invsqrt_pd(rsq12);
1506 rinv20 = gmx_mm256_invsqrt_pd(rsq20);
1507 rinv21 = gmx_mm256_invsqrt_pd(rsq21);
1508 rinv22 = gmx_mm256_invsqrt_pd(rsq22);
1510 rinvsq00 = _mm256_mul_pd(rinv00,rinv00);
1511 rinvsq01 = _mm256_mul_pd(rinv01,rinv01);
1512 rinvsq02 = _mm256_mul_pd(rinv02,rinv02);
1513 rinvsq10 = _mm256_mul_pd(rinv10,rinv10);
1514 rinvsq11 = _mm256_mul_pd(rinv11,rinv11);
1515 rinvsq12 = _mm256_mul_pd(rinv12,rinv12);
1516 rinvsq20 = _mm256_mul_pd(rinv20,rinv20);
1517 rinvsq21 = _mm256_mul_pd(rinv21,rinv21);
1518 rinvsq22 = _mm256_mul_pd(rinv22,rinv22);
1520 fjx0 = _mm256_setzero_pd();
1521 fjy0 = _mm256_setzero_pd();
1522 fjz0 = _mm256_setzero_pd();
1523 fjx1 = _mm256_setzero_pd();
1524 fjy1 = _mm256_setzero_pd();
1525 fjz1 = _mm256_setzero_pd();
1526 fjx2 = _mm256_setzero_pd();
1527 fjy2 = _mm256_setzero_pd();
1528 fjz2 = _mm256_setzero_pd();
1530 /**************************
1531 * CALCULATE INTERACTIONS *
1532 **************************/
1534 /* REACTION-FIELD ELECTROSTATICS */
1535 felec = _mm256_mul_pd(qq00,_mm256_sub_pd(_mm256_mul_pd(rinv00,rinvsq00),krf2));
1537 /* LENNARD-JONES DISPERSION/REPULSION */
1539 rinvsix = _mm256_mul_pd(_mm256_mul_pd(rinvsq00,rinvsq00),rinvsq00);
1540 fvdw = _mm256_mul_pd(_mm256_sub_pd(_mm256_mul_pd(c12_00,rinvsix),c6_00),_mm256_mul_pd(rinvsix,rinvsq00));
1542 fscal = _mm256_add_pd(felec,fvdw);
1544 fscal = _mm256_andnot_pd(dummy_mask,fscal);
1546 /* Calculate temporary vectorial force */
1547 tx = _mm256_mul_pd(fscal,dx00);
1548 ty = _mm256_mul_pd(fscal,dy00);
1549 tz = _mm256_mul_pd(fscal,dz00);
1551 /* Update vectorial force */
1552 fix0 = _mm256_add_pd(fix0,tx);
1553 fiy0 = _mm256_add_pd(fiy0,ty);
1554 fiz0 = _mm256_add_pd(fiz0,tz);
1556 fjx0 = _mm256_add_pd(fjx0,tx);
1557 fjy0 = _mm256_add_pd(fjy0,ty);
1558 fjz0 = _mm256_add_pd(fjz0,tz);
1560 /**************************
1561 * CALCULATE INTERACTIONS *
1562 **************************/
1564 /* REACTION-FIELD ELECTROSTATICS */
1565 felec = _mm256_mul_pd(qq01,_mm256_sub_pd(_mm256_mul_pd(rinv01,rinvsq01),krf2));
1569 fscal = _mm256_andnot_pd(dummy_mask,fscal);
1571 /* Calculate temporary vectorial force */
1572 tx = _mm256_mul_pd(fscal,dx01);
1573 ty = _mm256_mul_pd(fscal,dy01);
1574 tz = _mm256_mul_pd(fscal,dz01);
1576 /* Update vectorial force */
1577 fix0 = _mm256_add_pd(fix0,tx);
1578 fiy0 = _mm256_add_pd(fiy0,ty);
1579 fiz0 = _mm256_add_pd(fiz0,tz);
1581 fjx1 = _mm256_add_pd(fjx1,tx);
1582 fjy1 = _mm256_add_pd(fjy1,ty);
1583 fjz1 = _mm256_add_pd(fjz1,tz);
1585 /**************************
1586 * CALCULATE INTERACTIONS *
1587 **************************/
1589 /* REACTION-FIELD ELECTROSTATICS */
1590 felec = _mm256_mul_pd(qq02,_mm256_sub_pd(_mm256_mul_pd(rinv02,rinvsq02),krf2));
1594 fscal = _mm256_andnot_pd(dummy_mask,fscal);
1596 /* Calculate temporary vectorial force */
1597 tx = _mm256_mul_pd(fscal,dx02);
1598 ty = _mm256_mul_pd(fscal,dy02);
1599 tz = _mm256_mul_pd(fscal,dz02);
1601 /* Update vectorial force */
1602 fix0 = _mm256_add_pd(fix0,tx);
1603 fiy0 = _mm256_add_pd(fiy0,ty);
1604 fiz0 = _mm256_add_pd(fiz0,tz);
1606 fjx2 = _mm256_add_pd(fjx2,tx);
1607 fjy2 = _mm256_add_pd(fjy2,ty);
1608 fjz2 = _mm256_add_pd(fjz2,tz);
1610 /**************************
1611 * CALCULATE INTERACTIONS *
1612 **************************/
1614 /* REACTION-FIELD ELECTROSTATICS */
1615 felec = _mm256_mul_pd(qq10,_mm256_sub_pd(_mm256_mul_pd(rinv10,rinvsq10),krf2));
1619 fscal = _mm256_andnot_pd(dummy_mask,fscal);
1621 /* Calculate temporary vectorial force */
1622 tx = _mm256_mul_pd(fscal,dx10);
1623 ty = _mm256_mul_pd(fscal,dy10);
1624 tz = _mm256_mul_pd(fscal,dz10);
1626 /* Update vectorial force */
1627 fix1 = _mm256_add_pd(fix1,tx);
1628 fiy1 = _mm256_add_pd(fiy1,ty);
1629 fiz1 = _mm256_add_pd(fiz1,tz);
1631 fjx0 = _mm256_add_pd(fjx0,tx);
1632 fjy0 = _mm256_add_pd(fjy0,ty);
1633 fjz0 = _mm256_add_pd(fjz0,tz);
1635 /**************************
1636 * CALCULATE INTERACTIONS *
1637 **************************/
1639 /* REACTION-FIELD ELECTROSTATICS */
1640 felec = _mm256_mul_pd(qq11,_mm256_sub_pd(_mm256_mul_pd(rinv11,rinvsq11),krf2));
1644 fscal = _mm256_andnot_pd(dummy_mask,fscal);
1646 /* Calculate temporary vectorial force */
1647 tx = _mm256_mul_pd(fscal,dx11);
1648 ty = _mm256_mul_pd(fscal,dy11);
1649 tz = _mm256_mul_pd(fscal,dz11);
1651 /* Update vectorial force */
1652 fix1 = _mm256_add_pd(fix1,tx);
1653 fiy1 = _mm256_add_pd(fiy1,ty);
1654 fiz1 = _mm256_add_pd(fiz1,tz);
1656 fjx1 = _mm256_add_pd(fjx1,tx);
1657 fjy1 = _mm256_add_pd(fjy1,ty);
1658 fjz1 = _mm256_add_pd(fjz1,tz);
1660 /**************************
1661 * CALCULATE INTERACTIONS *
1662 **************************/
1664 /* REACTION-FIELD ELECTROSTATICS */
1665 felec = _mm256_mul_pd(qq12,_mm256_sub_pd(_mm256_mul_pd(rinv12,rinvsq12),krf2));
1669 fscal = _mm256_andnot_pd(dummy_mask,fscal);
1671 /* Calculate temporary vectorial force */
1672 tx = _mm256_mul_pd(fscal,dx12);
1673 ty = _mm256_mul_pd(fscal,dy12);
1674 tz = _mm256_mul_pd(fscal,dz12);
1676 /* Update vectorial force */
1677 fix1 = _mm256_add_pd(fix1,tx);
1678 fiy1 = _mm256_add_pd(fiy1,ty);
1679 fiz1 = _mm256_add_pd(fiz1,tz);
1681 fjx2 = _mm256_add_pd(fjx2,tx);
1682 fjy2 = _mm256_add_pd(fjy2,ty);
1683 fjz2 = _mm256_add_pd(fjz2,tz);
1685 /**************************
1686 * CALCULATE INTERACTIONS *
1687 **************************/
1689 /* REACTION-FIELD ELECTROSTATICS */
1690 felec = _mm256_mul_pd(qq20,_mm256_sub_pd(_mm256_mul_pd(rinv20,rinvsq20),krf2));
1694 fscal = _mm256_andnot_pd(dummy_mask,fscal);
1696 /* Calculate temporary vectorial force */
1697 tx = _mm256_mul_pd(fscal,dx20);
1698 ty = _mm256_mul_pd(fscal,dy20);
1699 tz = _mm256_mul_pd(fscal,dz20);
1701 /* Update vectorial force */
1702 fix2 = _mm256_add_pd(fix2,tx);
1703 fiy2 = _mm256_add_pd(fiy2,ty);
1704 fiz2 = _mm256_add_pd(fiz2,tz);
1706 fjx0 = _mm256_add_pd(fjx0,tx);
1707 fjy0 = _mm256_add_pd(fjy0,ty);
1708 fjz0 = _mm256_add_pd(fjz0,tz);
1710 /**************************
1711 * CALCULATE INTERACTIONS *
1712 **************************/
1714 /* REACTION-FIELD ELECTROSTATICS */
1715 felec = _mm256_mul_pd(qq21,_mm256_sub_pd(_mm256_mul_pd(rinv21,rinvsq21),krf2));
1719 fscal = _mm256_andnot_pd(dummy_mask,fscal);
1721 /* Calculate temporary vectorial force */
1722 tx = _mm256_mul_pd(fscal,dx21);
1723 ty = _mm256_mul_pd(fscal,dy21);
1724 tz = _mm256_mul_pd(fscal,dz21);
1726 /* Update vectorial force */
1727 fix2 = _mm256_add_pd(fix2,tx);
1728 fiy2 = _mm256_add_pd(fiy2,ty);
1729 fiz2 = _mm256_add_pd(fiz2,tz);
1731 fjx1 = _mm256_add_pd(fjx1,tx);
1732 fjy1 = _mm256_add_pd(fjy1,ty);
1733 fjz1 = _mm256_add_pd(fjz1,tz);
1735 /**************************
1736 * CALCULATE INTERACTIONS *
1737 **************************/
1739 /* REACTION-FIELD ELECTROSTATICS */
1740 felec = _mm256_mul_pd(qq22,_mm256_sub_pd(_mm256_mul_pd(rinv22,rinvsq22),krf2));
1744 fscal = _mm256_andnot_pd(dummy_mask,fscal);
1746 /* Calculate temporary vectorial force */
1747 tx = _mm256_mul_pd(fscal,dx22);
1748 ty = _mm256_mul_pd(fscal,dy22);
1749 tz = _mm256_mul_pd(fscal,dz22);
1751 /* Update vectorial force */
1752 fix2 = _mm256_add_pd(fix2,tx);
1753 fiy2 = _mm256_add_pd(fiy2,ty);
1754 fiz2 = _mm256_add_pd(fiz2,tz);
1756 fjx2 = _mm256_add_pd(fjx2,tx);
1757 fjy2 = _mm256_add_pd(fjy2,ty);
1758 fjz2 = _mm256_add_pd(fjz2,tz);
1760 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
1761 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
1762 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
1763 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
1765 gmx_mm256_decrement_3rvec_4ptr_swizzle_pd(fjptrA,fjptrB,fjptrC,fjptrD,
1766 fjx0,fjy0,fjz0,fjx1,fjy1,fjz1,fjx2,fjy2,fjz2);
1768 /* Inner loop uses 250 flops */
1771 /* End of innermost loop */
1773 gmx_mm256_update_iforce_3atom_swizzle_pd(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,
1774 f+i_coord_offset,fshift+i_shift_offset);
1776 /* Increment number of inner iterations */
1777 inneriter += j_index_end - j_index_start;
1779 /* Outer loop uses 18 flops */
1782 /* Increment number of outer iterations */
1785 /* Update outer/inner flops */
1787 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W3W3_F,outeriter*18 + inneriter*250);