Fix segmentation fault in minimize
[gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_avx_128_fma_single / nb_kernel_ElecNone_VdwLJEwSh_GeomP1P1_avx_128_fma_single.cpp
blobfa98f2d4b2cae8ac7a191c69d3fad2b6d75afd9d
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015,2017,2018, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
36 * Note: this file was generated by the GROMACS avx_128_fma_single kernel generator.
38 #include "gmxpre.h"
40 #include "config.h"
42 #include <math.h>
44 #include "../nb_kernel.h"
45 #include "gromacs/gmxlib/nrnb.h"
47 #include "kernelutil_x86_avx_128_fma_single.h"
50 * Gromacs nonbonded kernel: nb_kernel_ElecNone_VdwLJEwSh_GeomP1P1_VF_avx_128_fma_single
51 * Electrostatics interaction: None
52 * VdW interaction: LJEwald
53 * Geometry: Particle-Particle
54 * Calculate force/pot: PotentialAndForce
56 void
57 nb_kernel_ElecNone_VdwLJEwSh_GeomP1P1_VF_avx_128_fma_single
58 (t_nblist * gmx_restrict nlist,
59 rvec * gmx_restrict xx,
60 rvec * gmx_restrict ff,
61 struct t_forcerec * gmx_restrict fr,
62 t_mdatoms * gmx_restrict mdatoms,
63 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
64 t_nrnb * gmx_restrict nrnb)
66 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
67 * just 0 for non-waters.
68 * Suffixes A,B,C,D refer to j loop unrolling done with AVX_128, e.g. for the four different
69 * jnr indices corresponding to data put in the four positions in the SIMD register.
71 int i_shift_offset,i_coord_offset,outeriter,inneriter;
72 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
73 int jnrA,jnrB,jnrC,jnrD;
74 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
75 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
76 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
77 real rcutoff_scalar;
78 real *shiftvec,*fshift,*x,*f;
79 real *fjptrA,*fjptrB,*fjptrC,*fjptrD;
80 real scratch[4*DIM];
81 __m128 fscal,rcutoff,rcutoff2,jidxall;
82 int vdwioffset0;
83 __m128 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
84 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
85 __m128 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
86 __m128 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
87 int nvdwtype;
88 __m128 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
89 int *vdwtype;
90 real *vdwparam;
91 __m128 one_sixth = _mm_set1_ps(1.0/6.0);
92 __m128 one_twelfth = _mm_set1_ps(1.0/12.0);
93 __m128 c6grid_00;
94 real *vdwgridparam;
95 __m128 ewclj,ewclj2,ewclj6,ewcljrsq,poly,exponent,f6A,f6B,sh_lj_ewald;
96 __m128 one_half = _mm_set1_ps(0.5);
97 __m128 minus_one = _mm_set1_ps(-1.0);
98 __m128 dummy_mask,cutoff_mask;
99 __m128 signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
100 __m128 one = _mm_set1_ps(1.0);
101 __m128 two = _mm_set1_ps(2.0);
102 x = xx[0];
103 f = ff[0];
105 nri = nlist->nri;
106 iinr = nlist->iinr;
107 jindex = nlist->jindex;
108 jjnr = nlist->jjnr;
109 shiftidx = nlist->shift;
110 gid = nlist->gid;
111 shiftvec = fr->shift_vec[0];
112 fshift = fr->fshift[0];
113 nvdwtype = fr->ntype;
114 vdwparam = fr->nbfp;
115 vdwtype = mdatoms->typeA;
116 vdwgridparam = fr->ljpme_c6grid;
117 sh_lj_ewald = _mm_set1_ps(fr->ic->sh_lj_ewald);
118 ewclj = _mm_set1_ps(fr->ic->ewaldcoeff_lj);
119 ewclj2 = _mm_mul_ps(minus_one,_mm_mul_ps(ewclj,ewclj));
121 rcutoff_scalar = fr->ic->rvdw;
122 rcutoff = _mm_set1_ps(rcutoff_scalar);
123 rcutoff2 = _mm_mul_ps(rcutoff,rcutoff);
125 sh_vdw_invrcut6 = _mm_set1_ps(fr->ic->sh_invrc6);
126 rvdw = _mm_set1_ps(fr->ic->rvdw);
128 /* Avoid stupid compiler warnings */
129 jnrA = jnrB = jnrC = jnrD = 0;
130 j_coord_offsetA = 0;
131 j_coord_offsetB = 0;
132 j_coord_offsetC = 0;
133 j_coord_offsetD = 0;
135 outeriter = 0;
136 inneriter = 0;
138 for(iidx=0;iidx<4*DIM;iidx++)
140 scratch[iidx] = 0.0;
143 /* Start outer loop over neighborlists */
144 for(iidx=0; iidx<nri; iidx++)
146 /* Load shift vector for this list */
147 i_shift_offset = DIM*shiftidx[iidx];
149 /* Load limits for loop over neighbors */
150 j_index_start = jindex[iidx];
151 j_index_end = jindex[iidx+1];
153 /* Get outer coordinate index */
154 inr = iinr[iidx];
155 i_coord_offset = DIM*inr;
157 /* Load i particle coords and add shift vector */
158 gmx_mm_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
160 fix0 = _mm_setzero_ps();
161 fiy0 = _mm_setzero_ps();
162 fiz0 = _mm_setzero_ps();
164 /* Load parameters for i particles */
165 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
167 /* Reset potential sums */
168 vvdwsum = _mm_setzero_ps();
170 /* Start inner kernel loop */
171 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
174 /* Get j neighbor index, and coordinate index */
175 jnrA = jjnr[jidx];
176 jnrB = jjnr[jidx+1];
177 jnrC = jjnr[jidx+2];
178 jnrD = jjnr[jidx+3];
179 j_coord_offsetA = DIM*jnrA;
180 j_coord_offsetB = DIM*jnrB;
181 j_coord_offsetC = DIM*jnrC;
182 j_coord_offsetD = DIM*jnrD;
184 /* load j atom coordinates */
185 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
186 x+j_coord_offsetC,x+j_coord_offsetD,
187 &jx0,&jy0,&jz0);
189 /* Calculate displacement vector */
190 dx00 = _mm_sub_ps(ix0,jx0);
191 dy00 = _mm_sub_ps(iy0,jy0);
192 dz00 = _mm_sub_ps(iz0,jz0);
194 /* Calculate squared distance and things based on it */
195 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
197 rinv00 = avx128fma_invsqrt_f(rsq00);
199 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
201 /* Load parameters for j particles */
202 vdwjidx0A = 2*vdwtype[jnrA+0];
203 vdwjidx0B = 2*vdwtype[jnrB+0];
204 vdwjidx0C = 2*vdwtype[jnrC+0];
205 vdwjidx0D = 2*vdwtype[jnrD+0];
207 /**************************
208 * CALCULATE INTERACTIONS *
209 **************************/
211 if (gmx_mm_any_lt(rsq00,rcutoff2))
214 r00 = _mm_mul_ps(rsq00,rinv00);
216 /* Compute parameters for interactions between i and j atoms */
217 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
218 vdwparam+vdwioffset0+vdwjidx0B,
219 vdwparam+vdwioffset0+vdwjidx0C,
220 vdwparam+vdwioffset0+vdwjidx0D,
221 &c6_00,&c12_00);
223 c6grid_00 = gmx_mm_load_4real_swizzle_ps(vdwgridparam+vdwioffset0+vdwjidx0A,
224 vdwgridparam+vdwioffset0+vdwjidx0B,
225 vdwgridparam+vdwioffset0+vdwjidx0C,
226 vdwgridparam+vdwioffset0+vdwjidx0D);
228 /* Analytical LJ-PME */
229 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
230 ewcljrsq = _mm_mul_ps(ewclj2,rsq00);
231 ewclj6 = _mm_mul_ps(ewclj2,_mm_mul_ps(ewclj2,ewclj2));
232 exponent = avx128fma_exp_f(ewcljrsq);
233 /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
234 poly = _mm_mul_ps(exponent,_mm_macc_ps(_mm_mul_ps(ewcljrsq,ewcljrsq),one_half,_mm_sub_ps(one,ewcljrsq)));
235 /* vvdw6 = [C6 - C6grid * (1-poly)]/r6 */
236 vvdw6 = _mm_mul_ps(_mm_macc_ps(-c6grid_00,_mm_sub_ps(one,poly),c6_00),rinvsix);
237 vvdw12 = _mm_mul_ps(c12_00,_mm_mul_ps(rinvsix,rinvsix));
238 vvdw = _mm_msub_ps(_mm_nmacc_ps(c12_00,_mm_mul_ps(sh_vdw_invrcut6,sh_vdw_invrcut6),vvdw12),one_twelfth,
239 _mm_mul_ps(_mm_sub_ps(vvdw6,_mm_macc_ps(c6grid_00,sh_lj_ewald,_mm_mul_ps(c6_00,sh_vdw_invrcut6))),one_sixth));
240 /* fvdw = vvdw12/r - (vvdw6/r + (C6grid * exponent * beta^6)/r) */
241 fvdw = _mm_mul_ps(_mm_add_ps(vvdw12,_mm_msub_ps(_mm_mul_ps(c6grid_00,one_sixth),_mm_mul_ps(exponent,ewclj6),vvdw6)),rinvsq00);
243 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
245 /* Update potential sum for this i atom from the interaction with this j atom. */
246 vvdw = _mm_and_ps(vvdw,cutoff_mask);
247 vvdwsum = _mm_add_ps(vvdwsum,vvdw);
249 fscal = fvdw;
251 fscal = _mm_and_ps(fscal,cutoff_mask);
253 /* Update vectorial force */
254 fix0 = _mm_macc_ps(dx00,fscal,fix0);
255 fiy0 = _mm_macc_ps(dy00,fscal,fiy0);
256 fiz0 = _mm_macc_ps(dz00,fscal,fiz0);
258 fjptrA = f+j_coord_offsetA;
259 fjptrB = f+j_coord_offsetB;
260 fjptrC = f+j_coord_offsetC;
261 fjptrD = f+j_coord_offsetD;
262 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,
263 _mm_mul_ps(dx00,fscal),
264 _mm_mul_ps(dy00,fscal),
265 _mm_mul_ps(dz00,fscal));
269 /* Inner loop uses 59 flops */
272 if(jidx<j_index_end)
275 /* Get j neighbor index, and coordinate index */
276 jnrlistA = jjnr[jidx];
277 jnrlistB = jjnr[jidx+1];
278 jnrlistC = jjnr[jidx+2];
279 jnrlistD = jjnr[jidx+3];
280 /* Sign of each element will be negative for non-real atoms.
281 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
282 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
284 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
285 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
286 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
287 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
288 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
289 j_coord_offsetA = DIM*jnrA;
290 j_coord_offsetB = DIM*jnrB;
291 j_coord_offsetC = DIM*jnrC;
292 j_coord_offsetD = DIM*jnrD;
294 /* load j atom coordinates */
295 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
296 x+j_coord_offsetC,x+j_coord_offsetD,
297 &jx0,&jy0,&jz0);
299 /* Calculate displacement vector */
300 dx00 = _mm_sub_ps(ix0,jx0);
301 dy00 = _mm_sub_ps(iy0,jy0);
302 dz00 = _mm_sub_ps(iz0,jz0);
304 /* Calculate squared distance and things based on it */
305 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
307 rinv00 = avx128fma_invsqrt_f(rsq00);
309 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
311 /* Load parameters for j particles */
312 vdwjidx0A = 2*vdwtype[jnrA+0];
313 vdwjidx0B = 2*vdwtype[jnrB+0];
314 vdwjidx0C = 2*vdwtype[jnrC+0];
315 vdwjidx0D = 2*vdwtype[jnrD+0];
317 /**************************
318 * CALCULATE INTERACTIONS *
319 **************************/
321 if (gmx_mm_any_lt(rsq00,rcutoff2))
324 r00 = _mm_mul_ps(rsq00,rinv00);
325 r00 = _mm_andnot_ps(dummy_mask,r00);
327 /* Compute parameters for interactions between i and j atoms */
328 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
329 vdwparam+vdwioffset0+vdwjidx0B,
330 vdwparam+vdwioffset0+vdwjidx0C,
331 vdwparam+vdwioffset0+vdwjidx0D,
332 &c6_00,&c12_00);
334 c6grid_00 = gmx_mm_load_4real_swizzle_ps(vdwgridparam+vdwioffset0+vdwjidx0A,
335 vdwgridparam+vdwioffset0+vdwjidx0B,
336 vdwgridparam+vdwioffset0+vdwjidx0C,
337 vdwgridparam+vdwioffset0+vdwjidx0D);
339 /* Analytical LJ-PME */
340 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
341 ewcljrsq = _mm_mul_ps(ewclj2,rsq00);
342 ewclj6 = _mm_mul_ps(ewclj2,_mm_mul_ps(ewclj2,ewclj2));
343 exponent = avx128fma_exp_f(ewcljrsq);
344 /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
345 poly = _mm_mul_ps(exponent,_mm_macc_ps(_mm_mul_ps(ewcljrsq,ewcljrsq),one_half,_mm_sub_ps(one,ewcljrsq)));
346 /* vvdw6 = [C6 - C6grid * (1-poly)]/r6 */
347 vvdw6 = _mm_mul_ps(_mm_macc_ps(-c6grid_00,_mm_sub_ps(one,poly),c6_00),rinvsix);
348 vvdw12 = _mm_mul_ps(c12_00,_mm_mul_ps(rinvsix,rinvsix));
349 vvdw = _mm_msub_ps(_mm_nmacc_ps(c12_00,_mm_mul_ps(sh_vdw_invrcut6,sh_vdw_invrcut6),vvdw12),one_twelfth,
350 _mm_mul_ps(_mm_sub_ps(vvdw6,_mm_macc_ps(c6grid_00,sh_lj_ewald,_mm_mul_ps(c6_00,sh_vdw_invrcut6))),one_sixth));
351 /* fvdw = vvdw12/r - (vvdw6/r + (C6grid * exponent * beta^6)/r) */
352 fvdw = _mm_mul_ps(_mm_add_ps(vvdw12,_mm_msub_ps(_mm_mul_ps(c6grid_00,one_sixth),_mm_mul_ps(exponent,ewclj6),vvdw6)),rinvsq00);
354 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
356 /* Update potential sum for this i atom from the interaction with this j atom. */
357 vvdw = _mm_and_ps(vvdw,cutoff_mask);
358 vvdw = _mm_andnot_ps(dummy_mask,vvdw);
359 vvdwsum = _mm_add_ps(vvdwsum,vvdw);
361 fscal = fvdw;
363 fscal = _mm_and_ps(fscal,cutoff_mask);
365 fscal = _mm_andnot_ps(dummy_mask,fscal);
367 /* Update vectorial force */
368 fix0 = _mm_macc_ps(dx00,fscal,fix0);
369 fiy0 = _mm_macc_ps(dy00,fscal,fiy0);
370 fiz0 = _mm_macc_ps(dz00,fscal,fiz0);
372 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
373 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
374 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
375 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
376 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,
377 _mm_mul_ps(dx00,fscal),
378 _mm_mul_ps(dy00,fscal),
379 _mm_mul_ps(dz00,fscal));
383 /* Inner loop uses 60 flops */
386 /* End of innermost loop */
388 gmx_mm_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
389 f+i_coord_offset,fshift+i_shift_offset);
391 ggid = gid[iidx];
392 /* Update potential energies */
393 gmx_mm_update_1pot_ps(vvdwsum,kernel_data->energygrp_vdw+ggid);
395 /* Increment number of inner iterations */
396 inneriter += j_index_end - j_index_start;
398 /* Outer loop uses 7 flops */
401 /* Increment number of outer iterations */
402 outeriter += nri;
404 /* Update outer/inner flops */
406 inc_nrnb(nrnb,eNR_NBKERNEL_VDW_VF,outeriter*7 + inneriter*60);
409 * Gromacs nonbonded kernel: nb_kernel_ElecNone_VdwLJEwSh_GeomP1P1_F_avx_128_fma_single
410 * Electrostatics interaction: None
411 * VdW interaction: LJEwald
412 * Geometry: Particle-Particle
413 * Calculate force/pot: Force
415 void
416 nb_kernel_ElecNone_VdwLJEwSh_GeomP1P1_F_avx_128_fma_single
417 (t_nblist * gmx_restrict nlist,
418 rvec * gmx_restrict xx,
419 rvec * gmx_restrict ff,
420 struct t_forcerec * gmx_restrict fr,
421 t_mdatoms * gmx_restrict mdatoms,
422 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
423 t_nrnb * gmx_restrict nrnb)
425 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
426 * just 0 for non-waters.
427 * Suffixes A,B,C,D refer to j loop unrolling done with AVX_128, e.g. for the four different
428 * jnr indices corresponding to data put in the four positions in the SIMD register.
430 int i_shift_offset,i_coord_offset,outeriter,inneriter;
431 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
432 int jnrA,jnrB,jnrC,jnrD;
433 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
434 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
435 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
436 real rcutoff_scalar;
437 real *shiftvec,*fshift,*x,*f;
438 real *fjptrA,*fjptrB,*fjptrC,*fjptrD;
439 real scratch[4*DIM];
440 __m128 fscal,rcutoff,rcutoff2,jidxall;
441 int vdwioffset0;
442 __m128 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
443 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
444 __m128 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
445 __m128 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
446 int nvdwtype;
447 __m128 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
448 int *vdwtype;
449 real *vdwparam;
450 __m128 one_sixth = _mm_set1_ps(1.0/6.0);
451 __m128 one_twelfth = _mm_set1_ps(1.0/12.0);
452 __m128 c6grid_00;
453 real *vdwgridparam;
454 __m128 ewclj,ewclj2,ewclj6,ewcljrsq,poly,exponent,f6A,f6B,sh_lj_ewald;
455 __m128 one_half = _mm_set1_ps(0.5);
456 __m128 minus_one = _mm_set1_ps(-1.0);
457 __m128 dummy_mask,cutoff_mask;
458 __m128 signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
459 __m128 one = _mm_set1_ps(1.0);
460 __m128 two = _mm_set1_ps(2.0);
461 x = xx[0];
462 f = ff[0];
464 nri = nlist->nri;
465 iinr = nlist->iinr;
466 jindex = nlist->jindex;
467 jjnr = nlist->jjnr;
468 shiftidx = nlist->shift;
469 gid = nlist->gid;
470 shiftvec = fr->shift_vec[0];
471 fshift = fr->fshift[0];
472 nvdwtype = fr->ntype;
473 vdwparam = fr->nbfp;
474 vdwtype = mdatoms->typeA;
475 vdwgridparam = fr->ljpme_c6grid;
476 sh_lj_ewald = _mm_set1_ps(fr->ic->sh_lj_ewald);
477 ewclj = _mm_set1_ps(fr->ic->ewaldcoeff_lj);
478 ewclj2 = _mm_mul_ps(minus_one,_mm_mul_ps(ewclj,ewclj));
480 rcutoff_scalar = fr->ic->rvdw;
481 rcutoff = _mm_set1_ps(rcutoff_scalar);
482 rcutoff2 = _mm_mul_ps(rcutoff,rcutoff);
484 sh_vdw_invrcut6 = _mm_set1_ps(fr->ic->sh_invrc6);
485 rvdw = _mm_set1_ps(fr->ic->rvdw);
487 /* Avoid stupid compiler warnings */
488 jnrA = jnrB = jnrC = jnrD = 0;
489 j_coord_offsetA = 0;
490 j_coord_offsetB = 0;
491 j_coord_offsetC = 0;
492 j_coord_offsetD = 0;
494 outeriter = 0;
495 inneriter = 0;
497 for(iidx=0;iidx<4*DIM;iidx++)
499 scratch[iidx] = 0.0;
502 /* Start outer loop over neighborlists */
503 for(iidx=0; iidx<nri; iidx++)
505 /* Load shift vector for this list */
506 i_shift_offset = DIM*shiftidx[iidx];
508 /* Load limits for loop over neighbors */
509 j_index_start = jindex[iidx];
510 j_index_end = jindex[iidx+1];
512 /* Get outer coordinate index */
513 inr = iinr[iidx];
514 i_coord_offset = DIM*inr;
516 /* Load i particle coords and add shift vector */
517 gmx_mm_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
519 fix0 = _mm_setzero_ps();
520 fiy0 = _mm_setzero_ps();
521 fiz0 = _mm_setzero_ps();
523 /* Load parameters for i particles */
524 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
526 /* Start inner kernel loop */
527 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
530 /* Get j neighbor index, and coordinate index */
531 jnrA = jjnr[jidx];
532 jnrB = jjnr[jidx+1];
533 jnrC = jjnr[jidx+2];
534 jnrD = jjnr[jidx+3];
535 j_coord_offsetA = DIM*jnrA;
536 j_coord_offsetB = DIM*jnrB;
537 j_coord_offsetC = DIM*jnrC;
538 j_coord_offsetD = DIM*jnrD;
540 /* load j atom coordinates */
541 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
542 x+j_coord_offsetC,x+j_coord_offsetD,
543 &jx0,&jy0,&jz0);
545 /* Calculate displacement vector */
546 dx00 = _mm_sub_ps(ix0,jx0);
547 dy00 = _mm_sub_ps(iy0,jy0);
548 dz00 = _mm_sub_ps(iz0,jz0);
550 /* Calculate squared distance and things based on it */
551 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
553 rinv00 = avx128fma_invsqrt_f(rsq00);
555 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
557 /* Load parameters for j particles */
558 vdwjidx0A = 2*vdwtype[jnrA+0];
559 vdwjidx0B = 2*vdwtype[jnrB+0];
560 vdwjidx0C = 2*vdwtype[jnrC+0];
561 vdwjidx0D = 2*vdwtype[jnrD+0];
563 /**************************
564 * CALCULATE INTERACTIONS *
565 **************************/
567 if (gmx_mm_any_lt(rsq00,rcutoff2))
570 r00 = _mm_mul_ps(rsq00,rinv00);
572 /* Compute parameters for interactions between i and j atoms */
573 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
574 vdwparam+vdwioffset0+vdwjidx0B,
575 vdwparam+vdwioffset0+vdwjidx0C,
576 vdwparam+vdwioffset0+vdwjidx0D,
577 &c6_00,&c12_00);
579 c6grid_00 = gmx_mm_load_4real_swizzle_ps(vdwgridparam+vdwioffset0+vdwjidx0A,
580 vdwgridparam+vdwioffset0+vdwjidx0B,
581 vdwgridparam+vdwioffset0+vdwjidx0C,
582 vdwgridparam+vdwioffset0+vdwjidx0D);
584 /* Analytical LJ-PME */
585 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
586 ewcljrsq = _mm_mul_ps(ewclj2,rsq00);
587 ewclj6 = _mm_mul_ps(ewclj2,_mm_mul_ps(ewclj2,ewclj2));
588 exponent = avx128fma_exp_f(ewcljrsq);
589 /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
590 poly = _mm_mul_ps(exponent,_mm_macc_ps(_mm_mul_ps(ewcljrsq,ewcljrsq),one_half,_mm_sub_ps(one,ewcljrsq)));
591 /* f6A = 6 * C6grid * (1 - poly) */
592 f6A = _mm_mul_ps(c6grid_00,_mm_sub_ps(one,poly));
593 /* f6B = C6grid * exponent * beta^6 */
594 f6B = _mm_mul_ps(_mm_mul_ps(c6grid_00,one_sixth),_mm_mul_ps(exponent,ewclj6));
595 /* fvdw = 12*C12/r13 - ((6*C6 - f6A)/r6 + f6B)/r */
596 fvdw = _mm_mul_ps(_mm_macc_ps(_mm_msub_ps(c12_00,rinvsix,_mm_sub_ps(c6_00,f6A)),rinvsix,f6B),rinvsq00);
598 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
600 fscal = fvdw;
602 fscal = _mm_and_ps(fscal,cutoff_mask);
604 /* Update vectorial force */
605 fix0 = _mm_macc_ps(dx00,fscal,fix0);
606 fiy0 = _mm_macc_ps(dy00,fscal,fiy0);
607 fiz0 = _mm_macc_ps(dz00,fscal,fiz0);
609 fjptrA = f+j_coord_offsetA;
610 fjptrB = f+j_coord_offsetB;
611 fjptrC = f+j_coord_offsetC;
612 fjptrD = f+j_coord_offsetD;
613 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,
614 _mm_mul_ps(dx00,fscal),
615 _mm_mul_ps(dy00,fscal),
616 _mm_mul_ps(dz00,fscal));
620 /* Inner loop uses 50 flops */
623 if(jidx<j_index_end)
626 /* Get j neighbor index, and coordinate index */
627 jnrlistA = jjnr[jidx];
628 jnrlistB = jjnr[jidx+1];
629 jnrlistC = jjnr[jidx+2];
630 jnrlistD = jjnr[jidx+3];
631 /* Sign of each element will be negative for non-real atoms.
632 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
633 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
635 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
636 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
637 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
638 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
639 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
640 j_coord_offsetA = DIM*jnrA;
641 j_coord_offsetB = DIM*jnrB;
642 j_coord_offsetC = DIM*jnrC;
643 j_coord_offsetD = DIM*jnrD;
645 /* load j atom coordinates */
646 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
647 x+j_coord_offsetC,x+j_coord_offsetD,
648 &jx0,&jy0,&jz0);
650 /* Calculate displacement vector */
651 dx00 = _mm_sub_ps(ix0,jx0);
652 dy00 = _mm_sub_ps(iy0,jy0);
653 dz00 = _mm_sub_ps(iz0,jz0);
655 /* Calculate squared distance and things based on it */
656 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
658 rinv00 = avx128fma_invsqrt_f(rsq00);
660 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
662 /* Load parameters for j particles */
663 vdwjidx0A = 2*vdwtype[jnrA+0];
664 vdwjidx0B = 2*vdwtype[jnrB+0];
665 vdwjidx0C = 2*vdwtype[jnrC+0];
666 vdwjidx0D = 2*vdwtype[jnrD+0];
668 /**************************
669 * CALCULATE INTERACTIONS *
670 **************************/
672 if (gmx_mm_any_lt(rsq00,rcutoff2))
675 r00 = _mm_mul_ps(rsq00,rinv00);
676 r00 = _mm_andnot_ps(dummy_mask,r00);
678 /* Compute parameters for interactions between i and j atoms */
679 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
680 vdwparam+vdwioffset0+vdwjidx0B,
681 vdwparam+vdwioffset0+vdwjidx0C,
682 vdwparam+vdwioffset0+vdwjidx0D,
683 &c6_00,&c12_00);
685 c6grid_00 = gmx_mm_load_4real_swizzle_ps(vdwgridparam+vdwioffset0+vdwjidx0A,
686 vdwgridparam+vdwioffset0+vdwjidx0B,
687 vdwgridparam+vdwioffset0+vdwjidx0C,
688 vdwgridparam+vdwioffset0+vdwjidx0D);
690 /* Analytical LJ-PME */
691 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
692 ewcljrsq = _mm_mul_ps(ewclj2,rsq00);
693 ewclj6 = _mm_mul_ps(ewclj2,_mm_mul_ps(ewclj2,ewclj2));
694 exponent = avx128fma_exp_f(ewcljrsq);
695 /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
696 poly = _mm_mul_ps(exponent,_mm_macc_ps(_mm_mul_ps(ewcljrsq,ewcljrsq),one_half,_mm_sub_ps(one,ewcljrsq)));
697 /* f6A = 6 * C6grid * (1 - poly) */
698 f6A = _mm_mul_ps(c6grid_00,_mm_sub_ps(one,poly));
699 /* f6B = C6grid * exponent * beta^6 */
700 f6B = _mm_mul_ps(_mm_mul_ps(c6grid_00,one_sixth),_mm_mul_ps(exponent,ewclj6));
701 /* fvdw = 12*C12/r13 - ((6*C6 - f6A)/r6 + f6B)/r */
702 fvdw = _mm_mul_ps(_mm_macc_ps(_mm_msub_ps(c12_00,rinvsix,_mm_sub_ps(c6_00,f6A)),rinvsix,f6B),rinvsq00);
704 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
706 fscal = fvdw;
708 fscal = _mm_and_ps(fscal,cutoff_mask);
710 fscal = _mm_andnot_ps(dummy_mask,fscal);
712 /* Update vectorial force */
713 fix0 = _mm_macc_ps(dx00,fscal,fix0);
714 fiy0 = _mm_macc_ps(dy00,fscal,fiy0);
715 fiz0 = _mm_macc_ps(dz00,fscal,fiz0);
717 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
718 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
719 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
720 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
721 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,
722 _mm_mul_ps(dx00,fscal),
723 _mm_mul_ps(dy00,fscal),
724 _mm_mul_ps(dz00,fscal));
728 /* Inner loop uses 51 flops */
731 /* End of innermost loop */
733 gmx_mm_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
734 f+i_coord_offset,fshift+i_shift_offset);
736 /* Increment number of inner iterations */
737 inneriter += j_index_end - j_index_start;
739 /* Outer loop uses 6 flops */
742 /* Increment number of outer iterations */
743 outeriter += nri;
745 /* Update outer/inner flops */
747 inc_nrnb(nrnb,eNR_NBKERNEL_VDW_F,outeriter*6 + inneriter*51);