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_128_fma_single kernel generator.
44 #include "../nb_kernel.h"
45 #include "types/simple.h"
49 #include "gromacs/simd/math_x86_avx_128_fma_single.h"
50 #include "kernelutil_x86_avx_128_fma_single.h"
53 * Gromacs nonbonded kernel: nb_kernel_ElecEwSh_VdwLJEwSh_GeomP1P1_VF_avx_128_fma_single
54 * Electrostatics interaction: Ewald
55 * VdW interaction: LJEwald
56 * Geometry: Particle-Particle
57 * Calculate force/pot: PotentialAndForce
60 nb_kernel_ElecEwSh_VdwLJEwSh_GeomP1P1_VF_avx_128_fma_single
61 (t_nblist * gmx_restrict nlist,
62 rvec * gmx_restrict xx,
63 rvec * gmx_restrict ff,
64 t_forcerec * gmx_restrict fr,
65 t_mdatoms * gmx_restrict mdatoms,
66 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
67 t_nrnb * gmx_restrict nrnb)
69 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
70 * just 0 for non-waters.
71 * Suffixes A,B,C,D refer to j loop unrolling done with AVX_128, e.g. for the four different
72 * jnr indices corresponding to data put in the four positions in the SIMD register.
74 int i_shift_offset,i_coord_offset,outeriter,inneriter;
75 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
76 int jnrA,jnrB,jnrC,jnrD;
77 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
78 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
79 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
81 real *shiftvec,*fshift,*x,*f;
82 real *fjptrA,*fjptrB,*fjptrC,*fjptrD;
84 __m128 fscal,rcutoff,rcutoff2,jidxall;
86 __m128 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
87 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
88 __m128 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
89 __m128 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
90 __m128 velec,felec,velecsum,facel,crf,krf,krf2;
93 __m128 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
96 __m128 one_sixth = _mm_set1_ps(1.0/6.0);
97 __m128 one_twelfth = _mm_set1_ps(1.0/12.0);
100 __m128 ewclj,ewclj2,ewclj6,ewcljrsq,poly,exponent,f6A,f6B,sh_lj_ewald;
101 __m128 one_half = _mm_set1_ps(0.5);
102 __m128 minus_one = _mm_set1_ps(-1.0);
104 __m128 ewtabscale,eweps,twoeweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
105 __m128 beta,beta2,beta3,zeta2,pmecorrF,pmecorrV,rinv3;
107 __m128 dummy_mask,cutoff_mask;
108 __m128 signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
109 __m128 one = _mm_set1_ps(1.0);
110 __m128 two = _mm_set1_ps(2.0);
116 jindex = nlist->jindex;
118 shiftidx = nlist->shift;
120 shiftvec = fr->shift_vec[0];
121 fshift = fr->fshift[0];
122 facel = _mm_set1_ps(fr->epsfac);
123 charge = mdatoms->chargeA;
124 nvdwtype = fr->ntype;
126 vdwtype = mdatoms->typeA;
127 vdwgridparam = fr->ljpme_c6grid;
128 sh_lj_ewald = _mm_set1_ps(fr->ic->sh_lj_ewald);
129 ewclj = _mm_set1_ps(fr->ewaldcoeff_lj);
130 ewclj2 = _mm_mul_ps(minus_one,_mm_mul_ps(ewclj,ewclj));
132 sh_ewald = _mm_set1_ps(fr->ic->sh_ewald);
133 beta = _mm_set1_ps(fr->ic->ewaldcoeff_q);
134 beta2 = _mm_mul_ps(beta,beta);
135 beta3 = _mm_mul_ps(beta,beta2);
136 ewtab = fr->ic->tabq_coul_FDV0;
137 ewtabscale = _mm_set1_ps(fr->ic->tabq_scale);
138 ewtabhalfspace = _mm_set1_ps(0.5/fr->ic->tabq_scale);
140 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
141 rcutoff_scalar = fr->rcoulomb;
142 rcutoff = _mm_set1_ps(rcutoff_scalar);
143 rcutoff2 = _mm_mul_ps(rcutoff,rcutoff);
145 sh_vdw_invrcut6 = _mm_set1_ps(fr->ic->sh_invrc6);
146 rvdw = _mm_set1_ps(fr->rvdw);
148 /* Avoid stupid compiler warnings */
149 jnrA = jnrB = jnrC = jnrD = 0;
158 for(iidx=0;iidx<4*DIM;iidx++)
163 /* Start outer loop over neighborlists */
164 for(iidx=0; iidx<nri; iidx++)
166 /* Load shift vector for this list */
167 i_shift_offset = DIM*shiftidx[iidx];
169 /* Load limits for loop over neighbors */
170 j_index_start = jindex[iidx];
171 j_index_end = jindex[iidx+1];
173 /* Get outer coordinate index */
175 i_coord_offset = DIM*inr;
177 /* Load i particle coords and add shift vector */
178 gmx_mm_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
180 fix0 = _mm_setzero_ps();
181 fiy0 = _mm_setzero_ps();
182 fiz0 = _mm_setzero_ps();
184 /* Load parameters for i particles */
185 iq0 = _mm_mul_ps(facel,_mm_load1_ps(charge+inr+0));
186 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
188 /* Reset potential sums */
189 velecsum = _mm_setzero_ps();
190 vvdwsum = _mm_setzero_ps();
192 /* Start inner kernel loop */
193 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
196 /* Get j neighbor index, and coordinate index */
201 j_coord_offsetA = DIM*jnrA;
202 j_coord_offsetB = DIM*jnrB;
203 j_coord_offsetC = DIM*jnrC;
204 j_coord_offsetD = DIM*jnrD;
206 /* load j atom coordinates */
207 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
208 x+j_coord_offsetC,x+j_coord_offsetD,
211 /* Calculate displacement vector */
212 dx00 = _mm_sub_ps(ix0,jx0);
213 dy00 = _mm_sub_ps(iy0,jy0);
214 dz00 = _mm_sub_ps(iz0,jz0);
216 /* Calculate squared distance and things based on it */
217 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
219 rinv00 = gmx_mm_invsqrt_ps(rsq00);
221 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
223 /* Load parameters for j particles */
224 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
225 charge+jnrC+0,charge+jnrD+0);
226 vdwjidx0A = 2*vdwtype[jnrA+0];
227 vdwjidx0B = 2*vdwtype[jnrB+0];
228 vdwjidx0C = 2*vdwtype[jnrC+0];
229 vdwjidx0D = 2*vdwtype[jnrD+0];
231 /**************************
232 * CALCULATE INTERACTIONS *
233 **************************/
235 if (gmx_mm_any_lt(rsq00,rcutoff2))
238 r00 = _mm_mul_ps(rsq00,rinv00);
240 /* Compute parameters for interactions between i and j atoms */
241 qq00 = _mm_mul_ps(iq0,jq0);
242 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
243 vdwparam+vdwioffset0+vdwjidx0B,
244 vdwparam+vdwioffset0+vdwjidx0C,
245 vdwparam+vdwioffset0+vdwjidx0D,
248 c6grid_00 = gmx_mm_load_4real_swizzle_ps(vdwgridparam+vdwioffset0+vdwjidx0A,
249 vdwgridparam+vdwioffset0+vdwjidx0B,
250 vdwgridparam+vdwioffset0+vdwjidx0C,
251 vdwgridparam+vdwioffset0+vdwjidx0D);
253 /* EWALD ELECTROSTATICS */
255 /* Analytical PME correction */
256 zeta2 = _mm_mul_ps(beta2,rsq00);
257 rinv3 = _mm_mul_ps(rinvsq00,rinv00);
258 pmecorrF = gmx_mm_pmecorrF_ps(zeta2);
259 felec = _mm_macc_ps(pmecorrF,beta3,rinv3);
260 felec = _mm_mul_ps(qq00,felec);
261 pmecorrV = gmx_mm_pmecorrV_ps(zeta2);
262 velec = _mm_nmacc_ps(pmecorrV,beta,_mm_sub_ps(rinv00,sh_ewald));
263 velec = _mm_mul_ps(qq00,velec);
265 /* Analytical LJ-PME */
266 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
267 ewcljrsq = _mm_mul_ps(ewclj2,rsq00);
268 ewclj6 = _mm_mul_ps(ewclj2,_mm_mul_ps(ewclj2,ewclj2));
269 exponent = gmx_simd_exp_r(ewcljrsq);
270 /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
271 poly = _mm_mul_ps(exponent,_mm_macc_ps(_mm_mul_ps(ewcljrsq,ewcljrsq),one_half,_mm_sub_ps(one,ewcljrsq)));
272 /* vvdw6 = [C6 - C6grid * (1-poly)]/r6 */
273 vvdw6 = _mm_mul_ps(_mm_macc_ps(-c6grid_00,_mm_sub_ps(one,poly),c6_00),rinvsix);
274 vvdw12 = _mm_mul_ps(c12_00,_mm_mul_ps(rinvsix,rinvsix));
275 vvdw = _mm_msub_ps(_mm_nmacc_ps(c12_00,_mm_mul_ps(sh_vdw_invrcut6,sh_vdw_invrcut6),vvdw12),one_twelfth,
276 _mm_mul_ps(_mm_sub_ps(vvdw6,_mm_macc_ps(c6grid_00,sh_lj_ewald,_mm_mul_ps(c6_00,sh_vdw_invrcut6))),one_sixth));
277 /* fvdw = vvdw12/r - (vvdw6/r + (C6grid * exponent * beta^6)/r) */
278 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);
280 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
282 /* Update potential sum for this i atom from the interaction with this j atom. */
283 velec = _mm_and_ps(velec,cutoff_mask);
284 velecsum = _mm_add_ps(velecsum,velec);
285 vvdw = _mm_and_ps(vvdw,cutoff_mask);
286 vvdwsum = _mm_add_ps(vvdwsum,vvdw);
288 fscal = _mm_add_ps(felec,fvdw);
290 fscal = _mm_and_ps(fscal,cutoff_mask);
292 /* Update vectorial force */
293 fix0 = _mm_macc_ps(dx00,fscal,fix0);
294 fiy0 = _mm_macc_ps(dy00,fscal,fiy0);
295 fiz0 = _mm_macc_ps(dz00,fscal,fiz0);
297 fjptrA = f+j_coord_offsetA;
298 fjptrB = f+j_coord_offsetB;
299 fjptrC = f+j_coord_offsetC;
300 fjptrD = f+j_coord_offsetD;
301 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,
302 _mm_mul_ps(dx00,fscal),
303 _mm_mul_ps(dy00,fscal),
304 _mm_mul_ps(dz00,fscal));
308 /* Inner loop uses 63 flops */
314 /* Get j neighbor index, and coordinate index */
315 jnrlistA = jjnr[jidx];
316 jnrlistB = jjnr[jidx+1];
317 jnrlistC = jjnr[jidx+2];
318 jnrlistD = jjnr[jidx+3];
319 /* Sign of each element will be negative for non-real atoms.
320 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
321 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
323 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
324 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
325 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
326 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
327 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
328 j_coord_offsetA = DIM*jnrA;
329 j_coord_offsetB = DIM*jnrB;
330 j_coord_offsetC = DIM*jnrC;
331 j_coord_offsetD = DIM*jnrD;
333 /* load j atom coordinates */
334 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
335 x+j_coord_offsetC,x+j_coord_offsetD,
338 /* Calculate displacement vector */
339 dx00 = _mm_sub_ps(ix0,jx0);
340 dy00 = _mm_sub_ps(iy0,jy0);
341 dz00 = _mm_sub_ps(iz0,jz0);
343 /* Calculate squared distance and things based on it */
344 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
346 rinv00 = gmx_mm_invsqrt_ps(rsq00);
348 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
350 /* Load parameters for j particles */
351 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
352 charge+jnrC+0,charge+jnrD+0);
353 vdwjidx0A = 2*vdwtype[jnrA+0];
354 vdwjidx0B = 2*vdwtype[jnrB+0];
355 vdwjidx0C = 2*vdwtype[jnrC+0];
356 vdwjidx0D = 2*vdwtype[jnrD+0];
358 /**************************
359 * CALCULATE INTERACTIONS *
360 **************************/
362 if (gmx_mm_any_lt(rsq00,rcutoff2))
365 r00 = _mm_mul_ps(rsq00,rinv00);
366 r00 = _mm_andnot_ps(dummy_mask,r00);
368 /* Compute parameters for interactions between i and j atoms */
369 qq00 = _mm_mul_ps(iq0,jq0);
370 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
371 vdwparam+vdwioffset0+vdwjidx0B,
372 vdwparam+vdwioffset0+vdwjidx0C,
373 vdwparam+vdwioffset0+vdwjidx0D,
376 c6grid_00 = gmx_mm_load_4real_swizzle_ps(vdwgridparam+vdwioffset0+vdwjidx0A,
377 vdwgridparam+vdwioffset0+vdwjidx0B,
378 vdwgridparam+vdwioffset0+vdwjidx0C,
379 vdwgridparam+vdwioffset0+vdwjidx0D);
381 /* EWALD ELECTROSTATICS */
383 /* Analytical PME correction */
384 zeta2 = _mm_mul_ps(beta2,rsq00);
385 rinv3 = _mm_mul_ps(rinvsq00,rinv00);
386 pmecorrF = gmx_mm_pmecorrF_ps(zeta2);
387 felec = _mm_macc_ps(pmecorrF,beta3,rinv3);
388 felec = _mm_mul_ps(qq00,felec);
389 pmecorrV = gmx_mm_pmecorrV_ps(zeta2);
390 velec = _mm_nmacc_ps(pmecorrV,beta,_mm_sub_ps(rinv00,sh_ewald));
391 velec = _mm_mul_ps(qq00,velec);
393 /* Analytical LJ-PME */
394 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
395 ewcljrsq = _mm_mul_ps(ewclj2,rsq00);
396 ewclj6 = _mm_mul_ps(ewclj2,_mm_mul_ps(ewclj2,ewclj2));
397 exponent = gmx_simd_exp_r(ewcljrsq);
398 /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
399 poly = _mm_mul_ps(exponent,_mm_macc_ps(_mm_mul_ps(ewcljrsq,ewcljrsq),one_half,_mm_sub_ps(one,ewcljrsq)));
400 /* vvdw6 = [C6 - C6grid * (1-poly)]/r6 */
401 vvdw6 = _mm_mul_ps(_mm_macc_ps(-c6grid_00,_mm_sub_ps(one,poly),c6_00),rinvsix);
402 vvdw12 = _mm_mul_ps(c12_00,_mm_mul_ps(rinvsix,rinvsix));
403 vvdw = _mm_msub_ps(_mm_nmacc_ps(c12_00,_mm_mul_ps(sh_vdw_invrcut6,sh_vdw_invrcut6),vvdw12),one_twelfth,
404 _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));
405 /* fvdw = vvdw12/r - (vvdw6/r + (C6grid * exponent * beta^6)/r) */
406 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);
408 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
410 /* Update potential sum for this i atom from the interaction with this j atom. */
411 velec = _mm_and_ps(velec,cutoff_mask);
412 velec = _mm_andnot_ps(dummy_mask,velec);
413 velecsum = _mm_add_ps(velecsum,velec);
414 vvdw = _mm_and_ps(vvdw,cutoff_mask);
415 vvdw = _mm_andnot_ps(dummy_mask,vvdw);
416 vvdwsum = _mm_add_ps(vvdwsum,vvdw);
418 fscal = _mm_add_ps(felec,fvdw);
420 fscal = _mm_and_ps(fscal,cutoff_mask);
422 fscal = _mm_andnot_ps(dummy_mask,fscal);
424 /* Update vectorial force */
425 fix0 = _mm_macc_ps(dx00,fscal,fix0);
426 fiy0 = _mm_macc_ps(dy00,fscal,fiy0);
427 fiz0 = _mm_macc_ps(dz00,fscal,fiz0);
429 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
430 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
431 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
432 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
433 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,
434 _mm_mul_ps(dx00,fscal),
435 _mm_mul_ps(dy00,fscal),
436 _mm_mul_ps(dz00,fscal));
440 /* Inner loop uses 64 flops */
443 /* End of innermost loop */
445 gmx_mm_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
446 f+i_coord_offset,fshift+i_shift_offset);
449 /* Update potential energies */
450 gmx_mm_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
451 gmx_mm_update_1pot_ps(vvdwsum,kernel_data->energygrp_vdw+ggid);
453 /* Increment number of inner iterations */
454 inneriter += j_index_end - j_index_start;
456 /* Outer loop uses 9 flops */
459 /* Increment number of outer iterations */
462 /* Update outer/inner flops */
464 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_VF,outeriter*9 + inneriter*64);
467 * Gromacs nonbonded kernel: nb_kernel_ElecEwSh_VdwLJEwSh_GeomP1P1_F_avx_128_fma_single
468 * Electrostatics interaction: Ewald
469 * VdW interaction: LJEwald
470 * Geometry: Particle-Particle
471 * Calculate force/pot: Force
474 nb_kernel_ElecEwSh_VdwLJEwSh_GeomP1P1_F_avx_128_fma_single
475 (t_nblist * gmx_restrict nlist,
476 rvec * gmx_restrict xx,
477 rvec * gmx_restrict ff,
478 t_forcerec * gmx_restrict fr,
479 t_mdatoms * gmx_restrict mdatoms,
480 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
481 t_nrnb * gmx_restrict nrnb)
483 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
484 * just 0 for non-waters.
485 * Suffixes A,B,C,D refer to j loop unrolling done with AVX_128, e.g. for the four different
486 * jnr indices corresponding to data put in the four positions in the SIMD register.
488 int i_shift_offset,i_coord_offset,outeriter,inneriter;
489 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
490 int jnrA,jnrB,jnrC,jnrD;
491 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
492 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
493 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
495 real *shiftvec,*fshift,*x,*f;
496 real *fjptrA,*fjptrB,*fjptrC,*fjptrD;
498 __m128 fscal,rcutoff,rcutoff2,jidxall;
500 __m128 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
501 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
502 __m128 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
503 __m128 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
504 __m128 velec,felec,velecsum,facel,crf,krf,krf2;
507 __m128 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
510 __m128 one_sixth = _mm_set1_ps(1.0/6.0);
511 __m128 one_twelfth = _mm_set1_ps(1.0/12.0);
514 __m128 ewclj,ewclj2,ewclj6,ewcljrsq,poly,exponent,f6A,f6B,sh_lj_ewald;
515 __m128 one_half = _mm_set1_ps(0.5);
516 __m128 minus_one = _mm_set1_ps(-1.0);
518 __m128 ewtabscale,eweps,twoeweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
519 __m128 beta,beta2,beta3,zeta2,pmecorrF,pmecorrV,rinv3;
521 __m128 dummy_mask,cutoff_mask;
522 __m128 signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
523 __m128 one = _mm_set1_ps(1.0);
524 __m128 two = _mm_set1_ps(2.0);
530 jindex = nlist->jindex;
532 shiftidx = nlist->shift;
534 shiftvec = fr->shift_vec[0];
535 fshift = fr->fshift[0];
536 facel = _mm_set1_ps(fr->epsfac);
537 charge = mdatoms->chargeA;
538 nvdwtype = fr->ntype;
540 vdwtype = mdatoms->typeA;
541 vdwgridparam = fr->ljpme_c6grid;
542 sh_lj_ewald = _mm_set1_ps(fr->ic->sh_lj_ewald);
543 ewclj = _mm_set1_ps(fr->ewaldcoeff_lj);
544 ewclj2 = _mm_mul_ps(minus_one,_mm_mul_ps(ewclj,ewclj));
546 sh_ewald = _mm_set1_ps(fr->ic->sh_ewald);
547 beta = _mm_set1_ps(fr->ic->ewaldcoeff_q);
548 beta2 = _mm_mul_ps(beta,beta);
549 beta3 = _mm_mul_ps(beta,beta2);
550 ewtab = fr->ic->tabq_coul_F;
551 ewtabscale = _mm_set1_ps(fr->ic->tabq_scale);
552 ewtabhalfspace = _mm_set1_ps(0.5/fr->ic->tabq_scale);
554 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
555 rcutoff_scalar = fr->rcoulomb;
556 rcutoff = _mm_set1_ps(rcutoff_scalar);
557 rcutoff2 = _mm_mul_ps(rcutoff,rcutoff);
559 sh_vdw_invrcut6 = _mm_set1_ps(fr->ic->sh_invrc6);
560 rvdw = _mm_set1_ps(fr->rvdw);
562 /* Avoid stupid compiler warnings */
563 jnrA = jnrB = jnrC = jnrD = 0;
572 for(iidx=0;iidx<4*DIM;iidx++)
577 /* Start outer loop over neighborlists */
578 for(iidx=0; iidx<nri; iidx++)
580 /* Load shift vector for this list */
581 i_shift_offset = DIM*shiftidx[iidx];
583 /* Load limits for loop over neighbors */
584 j_index_start = jindex[iidx];
585 j_index_end = jindex[iidx+1];
587 /* Get outer coordinate index */
589 i_coord_offset = DIM*inr;
591 /* Load i particle coords and add shift vector */
592 gmx_mm_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
594 fix0 = _mm_setzero_ps();
595 fiy0 = _mm_setzero_ps();
596 fiz0 = _mm_setzero_ps();
598 /* Load parameters for i particles */
599 iq0 = _mm_mul_ps(facel,_mm_load1_ps(charge+inr+0));
600 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
602 /* Start inner kernel loop */
603 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
606 /* Get j neighbor index, and coordinate index */
611 j_coord_offsetA = DIM*jnrA;
612 j_coord_offsetB = DIM*jnrB;
613 j_coord_offsetC = DIM*jnrC;
614 j_coord_offsetD = DIM*jnrD;
616 /* load j atom coordinates */
617 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
618 x+j_coord_offsetC,x+j_coord_offsetD,
621 /* Calculate displacement vector */
622 dx00 = _mm_sub_ps(ix0,jx0);
623 dy00 = _mm_sub_ps(iy0,jy0);
624 dz00 = _mm_sub_ps(iz0,jz0);
626 /* Calculate squared distance and things based on it */
627 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
629 rinv00 = gmx_mm_invsqrt_ps(rsq00);
631 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
633 /* Load parameters for j particles */
634 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
635 charge+jnrC+0,charge+jnrD+0);
636 vdwjidx0A = 2*vdwtype[jnrA+0];
637 vdwjidx0B = 2*vdwtype[jnrB+0];
638 vdwjidx0C = 2*vdwtype[jnrC+0];
639 vdwjidx0D = 2*vdwtype[jnrD+0];
641 /**************************
642 * CALCULATE INTERACTIONS *
643 **************************/
645 if (gmx_mm_any_lt(rsq00,rcutoff2))
648 r00 = _mm_mul_ps(rsq00,rinv00);
650 /* Compute parameters for interactions between i and j atoms */
651 qq00 = _mm_mul_ps(iq0,jq0);
652 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
653 vdwparam+vdwioffset0+vdwjidx0B,
654 vdwparam+vdwioffset0+vdwjidx0C,
655 vdwparam+vdwioffset0+vdwjidx0D,
658 c6grid_00 = gmx_mm_load_4real_swizzle_ps(vdwgridparam+vdwioffset0+vdwjidx0A,
659 vdwgridparam+vdwioffset0+vdwjidx0B,
660 vdwgridparam+vdwioffset0+vdwjidx0C,
661 vdwgridparam+vdwioffset0+vdwjidx0D);
663 /* EWALD ELECTROSTATICS */
665 /* Analytical PME correction */
666 zeta2 = _mm_mul_ps(beta2,rsq00);
667 rinv3 = _mm_mul_ps(rinvsq00,rinv00);
668 pmecorrF = gmx_mm_pmecorrF_ps(zeta2);
669 felec = _mm_macc_ps(pmecorrF,beta3,rinv3);
670 felec = _mm_mul_ps(qq00,felec);
672 /* Analytical LJ-PME */
673 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
674 ewcljrsq = _mm_mul_ps(ewclj2,rsq00);
675 ewclj6 = _mm_mul_ps(ewclj2,_mm_mul_ps(ewclj2,ewclj2));
676 exponent = gmx_simd_exp_r(ewcljrsq);
677 /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
678 poly = _mm_mul_ps(exponent,_mm_macc_ps(_mm_mul_ps(ewcljrsq,ewcljrsq),one_half,_mm_sub_ps(one,ewcljrsq)));
679 /* f6A = 6 * C6grid * (1 - poly) */
680 f6A = _mm_mul_ps(c6grid_00,_mm_sub_ps(one,poly));
681 /* f6B = C6grid * exponent * beta^6 */
682 f6B = _mm_mul_ps(_mm_mul_ps(c6grid_00,one_sixth),_mm_mul_ps(exponent,ewclj6));
683 /* fvdw = 12*C12/r13 - ((6*C6 - f6A)/r6 + f6B)/r */
684 fvdw = _mm_mul_ps(_mm_macc_ps(_mm_msub_ps(c12_00,rinvsix,_mm_sub_ps(c6_00,f6A)),rinvsix,f6B),rinvsq00);
686 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
688 fscal = _mm_add_ps(felec,fvdw);
690 fscal = _mm_and_ps(fscal,cutoff_mask);
692 /* Update vectorial force */
693 fix0 = _mm_macc_ps(dx00,fscal,fix0);
694 fiy0 = _mm_macc_ps(dy00,fscal,fiy0);
695 fiz0 = _mm_macc_ps(dz00,fscal,fiz0);
697 fjptrA = f+j_coord_offsetA;
698 fjptrB = f+j_coord_offsetB;
699 fjptrC = f+j_coord_offsetC;
700 fjptrD = f+j_coord_offsetD;
701 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,
702 _mm_mul_ps(dx00,fscal),
703 _mm_mul_ps(dy00,fscal),
704 _mm_mul_ps(dz00,fscal));
708 /* Inner loop uses 52 flops */
714 /* Get j neighbor index, and coordinate index */
715 jnrlistA = jjnr[jidx];
716 jnrlistB = jjnr[jidx+1];
717 jnrlistC = jjnr[jidx+2];
718 jnrlistD = jjnr[jidx+3];
719 /* Sign of each element will be negative for non-real atoms.
720 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
721 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
723 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
724 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
725 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
726 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
727 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
728 j_coord_offsetA = DIM*jnrA;
729 j_coord_offsetB = DIM*jnrB;
730 j_coord_offsetC = DIM*jnrC;
731 j_coord_offsetD = DIM*jnrD;
733 /* load j atom coordinates */
734 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
735 x+j_coord_offsetC,x+j_coord_offsetD,
738 /* Calculate displacement vector */
739 dx00 = _mm_sub_ps(ix0,jx0);
740 dy00 = _mm_sub_ps(iy0,jy0);
741 dz00 = _mm_sub_ps(iz0,jz0);
743 /* Calculate squared distance and things based on it */
744 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
746 rinv00 = gmx_mm_invsqrt_ps(rsq00);
748 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
750 /* Load parameters for j particles */
751 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
752 charge+jnrC+0,charge+jnrD+0);
753 vdwjidx0A = 2*vdwtype[jnrA+0];
754 vdwjidx0B = 2*vdwtype[jnrB+0];
755 vdwjidx0C = 2*vdwtype[jnrC+0];
756 vdwjidx0D = 2*vdwtype[jnrD+0];
758 /**************************
759 * CALCULATE INTERACTIONS *
760 **************************/
762 if (gmx_mm_any_lt(rsq00,rcutoff2))
765 r00 = _mm_mul_ps(rsq00,rinv00);
766 r00 = _mm_andnot_ps(dummy_mask,r00);
768 /* Compute parameters for interactions between i and j atoms */
769 qq00 = _mm_mul_ps(iq0,jq0);
770 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
771 vdwparam+vdwioffset0+vdwjidx0B,
772 vdwparam+vdwioffset0+vdwjidx0C,
773 vdwparam+vdwioffset0+vdwjidx0D,
776 c6grid_00 = gmx_mm_load_4real_swizzle_ps(vdwgridparam+vdwioffset0+vdwjidx0A,
777 vdwgridparam+vdwioffset0+vdwjidx0B,
778 vdwgridparam+vdwioffset0+vdwjidx0C,
779 vdwgridparam+vdwioffset0+vdwjidx0D);
781 /* EWALD ELECTROSTATICS */
783 /* Analytical PME correction */
784 zeta2 = _mm_mul_ps(beta2,rsq00);
785 rinv3 = _mm_mul_ps(rinvsq00,rinv00);
786 pmecorrF = gmx_mm_pmecorrF_ps(zeta2);
787 felec = _mm_macc_ps(pmecorrF,beta3,rinv3);
788 felec = _mm_mul_ps(qq00,felec);
790 /* Analytical LJ-PME */
791 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
792 ewcljrsq = _mm_mul_ps(ewclj2,rsq00);
793 ewclj6 = _mm_mul_ps(ewclj2,_mm_mul_ps(ewclj2,ewclj2));
794 exponent = gmx_simd_exp_r(ewcljrsq);
795 /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
796 poly = _mm_mul_ps(exponent,_mm_macc_ps(_mm_mul_ps(ewcljrsq,ewcljrsq),one_half,_mm_sub_ps(one,ewcljrsq)));
797 /* f6A = 6 * C6grid * (1 - poly) */
798 f6A = _mm_mul_ps(c6grid_00,_mm_sub_ps(one,poly));
799 /* f6B = C6grid * exponent * beta^6 */
800 f6B = _mm_mul_ps(_mm_mul_ps(c6grid_00,one_sixth),_mm_mul_ps(exponent,ewclj6));
801 /* fvdw = 12*C12/r13 - ((6*C6 - f6A)/r6 + f6B)/r */
802 fvdw = _mm_mul_ps(_mm_macc_ps(_mm_msub_ps(c12_00,rinvsix,_mm_sub_ps(c6_00,f6A)),rinvsix,f6B),rinvsq00);
804 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
806 fscal = _mm_add_ps(felec,fvdw);
808 fscal = _mm_and_ps(fscal,cutoff_mask);
810 fscal = _mm_andnot_ps(dummy_mask,fscal);
812 /* Update vectorial force */
813 fix0 = _mm_macc_ps(dx00,fscal,fix0);
814 fiy0 = _mm_macc_ps(dy00,fscal,fiy0);
815 fiz0 = _mm_macc_ps(dz00,fscal,fiz0);
817 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
818 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
819 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
820 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
821 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,
822 _mm_mul_ps(dx00,fscal),
823 _mm_mul_ps(dy00,fscal),
824 _mm_mul_ps(dz00,fscal));
828 /* Inner loop uses 53 flops */
831 /* End of innermost loop */
833 gmx_mm_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
834 f+i_coord_offset,fshift+i_shift_offset);
836 /* Increment number of inner iterations */
837 inneriter += j_index_end - j_index_start;
839 /* Outer loop uses 7 flops */
842 /* Increment number of outer iterations */
845 /* Update outer/inner flops */
847 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_F,outeriter*7 + inneriter*53);