2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014, 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_double kernel generator.
44 #include "../nb_kernel.h"
45 #include "gromacs/legacyheaders/types/simple.h"
46 #include "gromacs/math/vec.h"
47 #include "gromacs/legacyheaders/nrnb.h"
49 #include "gromacs/simd/math_x86_avx_256_double.h"
50 #include "kernelutil_x86_avx_256_double.h"
53 * Gromacs nonbonded kernel: nb_kernel_ElecNone_VdwLJ_GeomP1P1_VF_avx_256_double
54 * Electrostatics interaction: None
55 * VdW interaction: LennardJones
56 * Geometry: Particle-Particle
57 * Calculate force/pot: PotentialAndForce
60 nb_kernel_ElecNone_VdwLJ_GeomP1P1_VF_avx_256_double
61 (t_nblist * gmx_restrict nlist,
62 rvec * gmx_restrict xx,
63 rvec * gmx_restrict ff,
64 t_forcerec * gmx_restrict fr,
65 t_mdatoms * gmx_restrict mdatoms,
66 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
67 t_nrnb * gmx_restrict nrnb)
69 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
70 * just 0 for non-waters.
71 * Suffixes A,B,C,D refer to j loop unrolling done with AVX, e.g. for the four different
72 * jnr indices corresponding to data put in the four positions in the SIMD register.
74 int i_shift_offset,i_coord_offset,outeriter,inneriter;
75 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
76 int jnrA,jnrB,jnrC,jnrD;
77 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
78 int jnrlistE,jnrlistF,jnrlistG,jnrlistH;
79 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
80 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
82 real *shiftvec,*fshift,*x,*f;
83 real *fjptrA,*fjptrB,*fjptrC,*fjptrD;
85 __m256d tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
86 real * vdwioffsetptr0;
87 __m256d ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
88 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
89 __m256d jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
90 __m256d dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
92 __m256d rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
95 __m256d one_sixth = _mm256_set1_pd(1.0/6.0);
96 __m256d one_twelfth = _mm256_set1_pd(1.0/12.0);
97 __m256d dummy_mask,cutoff_mask;
98 __m128 tmpmask0,tmpmask1;
99 __m256d signbit = _mm256_castsi256_pd( _mm256_set1_epi32(0x80000000) );
100 __m256d one = _mm256_set1_pd(1.0);
101 __m256d two = _mm256_set1_pd(2.0);
107 jindex = nlist->jindex;
109 shiftidx = nlist->shift;
111 shiftvec = fr->shift_vec[0];
112 fshift = fr->fshift[0];
113 nvdwtype = fr->ntype;
115 vdwtype = mdatoms->typeA;
117 /* Avoid stupid compiler warnings */
118 jnrA = jnrB = jnrC = jnrD = 0;
127 for(iidx=0;iidx<4*DIM;iidx++)
132 /* Start outer loop over neighborlists */
133 for(iidx=0; iidx<nri; iidx++)
135 /* Load shift vector for this list */
136 i_shift_offset = DIM*shiftidx[iidx];
138 /* Load limits for loop over neighbors */
139 j_index_start = jindex[iidx];
140 j_index_end = jindex[iidx+1];
142 /* Get outer coordinate index */
144 i_coord_offset = DIM*inr;
146 /* Load i particle coords and add shift vector */
147 gmx_mm256_load_shift_and_1rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
149 fix0 = _mm256_setzero_pd();
150 fiy0 = _mm256_setzero_pd();
151 fiz0 = _mm256_setzero_pd();
153 /* Load parameters for i particles */
154 vdwioffsetptr0 = vdwparam+2*nvdwtype*vdwtype[inr+0];
156 /* Reset potential sums */
157 vvdwsum = _mm256_setzero_pd();
159 /* Start inner kernel loop */
160 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
163 /* Get j neighbor index, and coordinate index */
168 j_coord_offsetA = DIM*jnrA;
169 j_coord_offsetB = DIM*jnrB;
170 j_coord_offsetC = DIM*jnrC;
171 j_coord_offsetD = DIM*jnrD;
173 /* load j atom coordinates */
174 gmx_mm256_load_1rvec_4ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
175 x+j_coord_offsetC,x+j_coord_offsetD,
178 /* Calculate displacement vector */
179 dx00 = _mm256_sub_pd(ix0,jx0);
180 dy00 = _mm256_sub_pd(iy0,jy0);
181 dz00 = _mm256_sub_pd(iz0,jz0);
183 /* Calculate squared distance and things based on it */
184 rsq00 = gmx_mm256_calc_rsq_pd(dx00,dy00,dz00);
186 rinvsq00 = gmx_mm256_inv_pd(rsq00);
188 /* Load parameters for j particles */
189 vdwjidx0A = 2*vdwtype[jnrA+0];
190 vdwjidx0B = 2*vdwtype[jnrB+0];
191 vdwjidx0C = 2*vdwtype[jnrC+0];
192 vdwjidx0D = 2*vdwtype[jnrD+0];
194 /**************************
195 * CALCULATE INTERACTIONS *
196 **************************/
198 /* Compute parameters for interactions between i and j atoms */
199 gmx_mm256_load_4pair_swizzle_pd(vdwioffsetptr0+vdwjidx0A,
200 vdwioffsetptr0+vdwjidx0B,
201 vdwioffsetptr0+vdwjidx0C,
202 vdwioffsetptr0+vdwjidx0D,
205 /* LENNARD-JONES DISPERSION/REPULSION */
207 rinvsix = _mm256_mul_pd(_mm256_mul_pd(rinvsq00,rinvsq00),rinvsq00);
208 vvdw6 = _mm256_mul_pd(c6_00,rinvsix);
209 vvdw12 = _mm256_mul_pd(c12_00,_mm256_mul_pd(rinvsix,rinvsix));
210 vvdw = _mm256_sub_pd( _mm256_mul_pd(vvdw12,one_twelfth) , _mm256_mul_pd(vvdw6,one_sixth) );
211 fvdw = _mm256_mul_pd(_mm256_sub_pd(vvdw12,vvdw6),rinvsq00);
213 /* Update potential sum for this i atom from the interaction with this j atom. */
214 vvdwsum = _mm256_add_pd(vvdwsum,vvdw);
218 /* Calculate temporary vectorial force */
219 tx = _mm256_mul_pd(fscal,dx00);
220 ty = _mm256_mul_pd(fscal,dy00);
221 tz = _mm256_mul_pd(fscal,dz00);
223 /* Update vectorial force */
224 fix0 = _mm256_add_pd(fix0,tx);
225 fiy0 = _mm256_add_pd(fiy0,ty);
226 fiz0 = _mm256_add_pd(fiz0,tz);
228 fjptrA = f+j_coord_offsetA;
229 fjptrB = f+j_coord_offsetB;
230 fjptrC = f+j_coord_offsetC;
231 fjptrD = f+j_coord_offsetD;
232 gmx_mm256_decrement_1rvec_4ptr_swizzle_pd(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
234 /* Inner loop uses 32 flops */
240 /* Get j neighbor index, and coordinate index */
241 jnrlistA = jjnr[jidx];
242 jnrlistB = jjnr[jidx+1];
243 jnrlistC = jjnr[jidx+2];
244 jnrlistD = jjnr[jidx+3];
245 /* Sign of each element will be negative for non-real atoms.
246 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
247 * so use it as val = _mm_andnot_pd(mask,val) to clear dummy entries.
249 tmpmask0 = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
251 tmpmask1 = _mm_permute_ps(tmpmask0,_GMX_MM_PERMUTE(3,3,2,2));
252 tmpmask0 = _mm_permute_ps(tmpmask0,_GMX_MM_PERMUTE(1,1,0,0));
253 dummy_mask = _mm256_castps_pd(gmx_mm256_set_m128(tmpmask1,tmpmask0));
255 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
256 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
257 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
258 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
259 j_coord_offsetA = DIM*jnrA;
260 j_coord_offsetB = DIM*jnrB;
261 j_coord_offsetC = DIM*jnrC;
262 j_coord_offsetD = DIM*jnrD;
264 /* load j atom coordinates */
265 gmx_mm256_load_1rvec_4ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
266 x+j_coord_offsetC,x+j_coord_offsetD,
269 /* Calculate displacement vector */
270 dx00 = _mm256_sub_pd(ix0,jx0);
271 dy00 = _mm256_sub_pd(iy0,jy0);
272 dz00 = _mm256_sub_pd(iz0,jz0);
274 /* Calculate squared distance and things based on it */
275 rsq00 = gmx_mm256_calc_rsq_pd(dx00,dy00,dz00);
277 rinvsq00 = gmx_mm256_inv_pd(rsq00);
279 /* Load parameters for j particles */
280 vdwjidx0A = 2*vdwtype[jnrA+0];
281 vdwjidx0B = 2*vdwtype[jnrB+0];
282 vdwjidx0C = 2*vdwtype[jnrC+0];
283 vdwjidx0D = 2*vdwtype[jnrD+0];
285 /**************************
286 * CALCULATE INTERACTIONS *
287 **************************/
289 /* Compute parameters for interactions between i and j atoms */
290 gmx_mm256_load_4pair_swizzle_pd(vdwioffsetptr0+vdwjidx0A,
291 vdwioffsetptr0+vdwjidx0B,
292 vdwioffsetptr0+vdwjidx0C,
293 vdwioffsetptr0+vdwjidx0D,
296 /* LENNARD-JONES DISPERSION/REPULSION */
298 rinvsix = _mm256_mul_pd(_mm256_mul_pd(rinvsq00,rinvsq00),rinvsq00);
299 vvdw6 = _mm256_mul_pd(c6_00,rinvsix);
300 vvdw12 = _mm256_mul_pd(c12_00,_mm256_mul_pd(rinvsix,rinvsix));
301 vvdw = _mm256_sub_pd( _mm256_mul_pd(vvdw12,one_twelfth) , _mm256_mul_pd(vvdw6,one_sixth) );
302 fvdw = _mm256_mul_pd(_mm256_sub_pd(vvdw12,vvdw6),rinvsq00);
304 /* Update potential sum for this i atom from the interaction with this j atom. */
305 vvdw = _mm256_andnot_pd(dummy_mask,vvdw);
306 vvdwsum = _mm256_add_pd(vvdwsum,vvdw);
310 fscal = _mm256_andnot_pd(dummy_mask,fscal);
312 /* Calculate temporary vectorial force */
313 tx = _mm256_mul_pd(fscal,dx00);
314 ty = _mm256_mul_pd(fscal,dy00);
315 tz = _mm256_mul_pd(fscal,dz00);
317 /* Update vectorial force */
318 fix0 = _mm256_add_pd(fix0,tx);
319 fiy0 = _mm256_add_pd(fiy0,ty);
320 fiz0 = _mm256_add_pd(fiz0,tz);
322 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
323 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
324 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
325 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
326 gmx_mm256_decrement_1rvec_4ptr_swizzle_pd(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
328 /* Inner loop uses 32 flops */
331 /* End of innermost loop */
333 gmx_mm256_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
334 f+i_coord_offset,fshift+i_shift_offset);
337 /* Update potential energies */
338 gmx_mm256_update_1pot_pd(vvdwsum,kernel_data->energygrp_vdw+ggid);
340 /* Increment number of inner iterations */
341 inneriter += j_index_end - j_index_start;
343 /* Outer loop uses 7 flops */
346 /* Increment number of outer iterations */
349 /* Update outer/inner flops */
351 inc_nrnb(nrnb,eNR_NBKERNEL_VDW_VF,outeriter*7 + inneriter*32);
354 * Gromacs nonbonded kernel: nb_kernel_ElecNone_VdwLJ_GeomP1P1_F_avx_256_double
355 * Electrostatics interaction: None
356 * VdW interaction: LennardJones
357 * Geometry: Particle-Particle
358 * Calculate force/pot: Force
361 nb_kernel_ElecNone_VdwLJ_GeomP1P1_F_avx_256_double
362 (t_nblist * gmx_restrict nlist,
363 rvec * gmx_restrict xx,
364 rvec * gmx_restrict ff,
365 t_forcerec * gmx_restrict fr,
366 t_mdatoms * gmx_restrict mdatoms,
367 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
368 t_nrnb * gmx_restrict nrnb)
370 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
371 * just 0 for non-waters.
372 * Suffixes A,B,C,D refer to j loop unrolling done with AVX, e.g. for the four different
373 * jnr indices corresponding to data put in the four positions in the SIMD register.
375 int i_shift_offset,i_coord_offset,outeriter,inneriter;
376 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
377 int jnrA,jnrB,jnrC,jnrD;
378 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
379 int jnrlistE,jnrlistF,jnrlistG,jnrlistH;
380 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
381 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
383 real *shiftvec,*fshift,*x,*f;
384 real *fjptrA,*fjptrB,*fjptrC,*fjptrD;
386 __m256d tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
387 real * vdwioffsetptr0;
388 __m256d ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
389 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
390 __m256d jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
391 __m256d dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
393 __m256d rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
396 __m256d one_sixth = _mm256_set1_pd(1.0/6.0);
397 __m256d one_twelfth = _mm256_set1_pd(1.0/12.0);
398 __m256d dummy_mask,cutoff_mask;
399 __m128 tmpmask0,tmpmask1;
400 __m256d signbit = _mm256_castsi256_pd( _mm256_set1_epi32(0x80000000) );
401 __m256d one = _mm256_set1_pd(1.0);
402 __m256d two = _mm256_set1_pd(2.0);
408 jindex = nlist->jindex;
410 shiftidx = nlist->shift;
412 shiftvec = fr->shift_vec[0];
413 fshift = fr->fshift[0];
414 nvdwtype = fr->ntype;
416 vdwtype = mdatoms->typeA;
418 /* Avoid stupid compiler warnings */
419 jnrA = jnrB = jnrC = jnrD = 0;
428 for(iidx=0;iidx<4*DIM;iidx++)
433 /* Start outer loop over neighborlists */
434 for(iidx=0; iidx<nri; iidx++)
436 /* Load shift vector for this list */
437 i_shift_offset = DIM*shiftidx[iidx];
439 /* Load limits for loop over neighbors */
440 j_index_start = jindex[iidx];
441 j_index_end = jindex[iidx+1];
443 /* Get outer coordinate index */
445 i_coord_offset = DIM*inr;
447 /* Load i particle coords and add shift vector */
448 gmx_mm256_load_shift_and_1rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
450 fix0 = _mm256_setzero_pd();
451 fiy0 = _mm256_setzero_pd();
452 fiz0 = _mm256_setzero_pd();
454 /* Load parameters for i particles */
455 vdwioffsetptr0 = vdwparam+2*nvdwtype*vdwtype[inr+0];
457 /* Start inner kernel loop */
458 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
461 /* Get j neighbor index, and coordinate index */
466 j_coord_offsetA = DIM*jnrA;
467 j_coord_offsetB = DIM*jnrB;
468 j_coord_offsetC = DIM*jnrC;
469 j_coord_offsetD = DIM*jnrD;
471 /* load j atom coordinates */
472 gmx_mm256_load_1rvec_4ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
473 x+j_coord_offsetC,x+j_coord_offsetD,
476 /* Calculate displacement vector */
477 dx00 = _mm256_sub_pd(ix0,jx0);
478 dy00 = _mm256_sub_pd(iy0,jy0);
479 dz00 = _mm256_sub_pd(iz0,jz0);
481 /* Calculate squared distance and things based on it */
482 rsq00 = gmx_mm256_calc_rsq_pd(dx00,dy00,dz00);
484 rinvsq00 = gmx_mm256_inv_pd(rsq00);
486 /* Load parameters for j particles */
487 vdwjidx0A = 2*vdwtype[jnrA+0];
488 vdwjidx0B = 2*vdwtype[jnrB+0];
489 vdwjidx0C = 2*vdwtype[jnrC+0];
490 vdwjidx0D = 2*vdwtype[jnrD+0];
492 /**************************
493 * CALCULATE INTERACTIONS *
494 **************************/
496 /* Compute parameters for interactions between i and j atoms */
497 gmx_mm256_load_4pair_swizzle_pd(vdwioffsetptr0+vdwjidx0A,
498 vdwioffsetptr0+vdwjidx0B,
499 vdwioffsetptr0+vdwjidx0C,
500 vdwioffsetptr0+vdwjidx0D,
503 /* LENNARD-JONES DISPERSION/REPULSION */
505 rinvsix = _mm256_mul_pd(_mm256_mul_pd(rinvsq00,rinvsq00),rinvsq00);
506 fvdw = _mm256_mul_pd(_mm256_sub_pd(_mm256_mul_pd(c12_00,rinvsix),c6_00),_mm256_mul_pd(rinvsix,rinvsq00));
510 /* Calculate temporary vectorial force */
511 tx = _mm256_mul_pd(fscal,dx00);
512 ty = _mm256_mul_pd(fscal,dy00);
513 tz = _mm256_mul_pd(fscal,dz00);
515 /* Update vectorial force */
516 fix0 = _mm256_add_pd(fix0,tx);
517 fiy0 = _mm256_add_pd(fiy0,ty);
518 fiz0 = _mm256_add_pd(fiz0,tz);
520 fjptrA = f+j_coord_offsetA;
521 fjptrB = f+j_coord_offsetB;
522 fjptrC = f+j_coord_offsetC;
523 fjptrD = f+j_coord_offsetD;
524 gmx_mm256_decrement_1rvec_4ptr_swizzle_pd(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
526 /* Inner loop uses 27 flops */
532 /* Get j neighbor index, and coordinate index */
533 jnrlistA = jjnr[jidx];
534 jnrlistB = jjnr[jidx+1];
535 jnrlistC = jjnr[jidx+2];
536 jnrlistD = jjnr[jidx+3];
537 /* Sign of each element will be negative for non-real atoms.
538 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
539 * so use it as val = _mm_andnot_pd(mask,val) to clear dummy entries.
541 tmpmask0 = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
543 tmpmask1 = _mm_permute_ps(tmpmask0,_GMX_MM_PERMUTE(3,3,2,2));
544 tmpmask0 = _mm_permute_ps(tmpmask0,_GMX_MM_PERMUTE(1,1,0,0));
545 dummy_mask = _mm256_castps_pd(gmx_mm256_set_m128(tmpmask1,tmpmask0));
547 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
548 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
549 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
550 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
551 j_coord_offsetA = DIM*jnrA;
552 j_coord_offsetB = DIM*jnrB;
553 j_coord_offsetC = DIM*jnrC;
554 j_coord_offsetD = DIM*jnrD;
556 /* load j atom coordinates */
557 gmx_mm256_load_1rvec_4ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
558 x+j_coord_offsetC,x+j_coord_offsetD,
561 /* Calculate displacement vector */
562 dx00 = _mm256_sub_pd(ix0,jx0);
563 dy00 = _mm256_sub_pd(iy0,jy0);
564 dz00 = _mm256_sub_pd(iz0,jz0);
566 /* Calculate squared distance and things based on it */
567 rsq00 = gmx_mm256_calc_rsq_pd(dx00,dy00,dz00);
569 rinvsq00 = gmx_mm256_inv_pd(rsq00);
571 /* Load parameters for j particles */
572 vdwjidx0A = 2*vdwtype[jnrA+0];
573 vdwjidx0B = 2*vdwtype[jnrB+0];
574 vdwjidx0C = 2*vdwtype[jnrC+0];
575 vdwjidx0D = 2*vdwtype[jnrD+0];
577 /**************************
578 * CALCULATE INTERACTIONS *
579 **************************/
581 /* Compute parameters for interactions between i and j atoms */
582 gmx_mm256_load_4pair_swizzle_pd(vdwioffsetptr0+vdwjidx0A,
583 vdwioffsetptr0+vdwjidx0B,
584 vdwioffsetptr0+vdwjidx0C,
585 vdwioffsetptr0+vdwjidx0D,
588 /* LENNARD-JONES DISPERSION/REPULSION */
590 rinvsix = _mm256_mul_pd(_mm256_mul_pd(rinvsq00,rinvsq00),rinvsq00);
591 fvdw = _mm256_mul_pd(_mm256_sub_pd(_mm256_mul_pd(c12_00,rinvsix),c6_00),_mm256_mul_pd(rinvsix,rinvsq00));
595 fscal = _mm256_andnot_pd(dummy_mask,fscal);
597 /* Calculate temporary vectorial force */
598 tx = _mm256_mul_pd(fscal,dx00);
599 ty = _mm256_mul_pd(fscal,dy00);
600 tz = _mm256_mul_pd(fscal,dz00);
602 /* Update vectorial force */
603 fix0 = _mm256_add_pd(fix0,tx);
604 fiy0 = _mm256_add_pd(fiy0,ty);
605 fiz0 = _mm256_add_pd(fiz0,tz);
607 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
608 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
609 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
610 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
611 gmx_mm256_decrement_1rvec_4ptr_swizzle_pd(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
613 /* Inner loop uses 27 flops */
616 /* End of innermost loop */
618 gmx_mm256_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
619 f+i_coord_offset,fshift+i_shift_offset);
621 /* Increment number of inner iterations */
622 inneriter += j_index_end - j_index_start;
624 /* Outer loop uses 6 flops */
627 /* Increment number of outer iterations */
630 /* Update outer/inner flops */
632 inc_nrnb(nrnb,eNR_NBKERNEL_VDW_F,outeriter*6 + inneriter*27);