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_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_128_fma_double.h"
50 #include "kernelutil_x86_avx_128_fma_double.h"
53 * Gromacs nonbonded kernel: nb_kernel_ElecEwSh_VdwLJEwSh_GeomP1P1_VF_avx_128_fma_double
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_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 refer to j loop unrolling done with SSE double precision, e.g. for the two 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;
77 int j_coord_offsetA,j_coord_offsetB;
78 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
80 real *shiftvec,*fshift,*x,*f;
81 __m128d tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
83 __m128d ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
84 int vdwjidx0A,vdwjidx0B;
85 __m128d jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
86 __m128d dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
87 __m128d velec,felec,velecsum,facel,crf,krf,krf2;
90 __m128d rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
93 __m128d one_sixth = _mm_set1_pd(1.0/6.0);
94 __m128d one_twelfth = _mm_set1_pd(1.0/12.0);
97 __m128d ewclj,ewclj2,ewclj6,ewcljrsq,poly,exponent,f6A,f6B,sh_lj_ewald;
98 __m128d one_half = _mm_set1_pd(0.5);
99 __m128d minus_one = _mm_set1_pd(-1.0);
101 __m128d ewtabscale,eweps,twoeweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
103 __m128d dummy_mask,cutoff_mask;
104 __m128d signbit = gmx_mm_castsi128_pd( _mm_set_epi32(0x80000000,0x00000000,0x80000000,0x00000000) );
105 __m128d one = _mm_set1_pd(1.0);
106 __m128d two = _mm_set1_pd(2.0);
112 jindex = nlist->jindex;
114 shiftidx = nlist->shift;
116 shiftvec = fr->shift_vec[0];
117 fshift = fr->fshift[0];
118 facel = _mm_set1_pd(fr->epsfac);
119 charge = mdatoms->chargeA;
120 nvdwtype = fr->ntype;
122 vdwtype = mdatoms->typeA;
123 vdwgridparam = fr->ljpme_c6grid;
124 sh_lj_ewald = _mm_set1_pd(fr->ic->sh_lj_ewald);
125 ewclj = _mm_set1_pd(fr->ewaldcoeff_lj);
126 ewclj2 = _mm_mul_pd(minus_one,_mm_mul_pd(ewclj,ewclj));
128 sh_ewald = _mm_set1_pd(fr->ic->sh_ewald);
129 ewtab = fr->ic->tabq_coul_FDV0;
130 ewtabscale = _mm_set1_pd(fr->ic->tabq_scale);
131 ewtabhalfspace = _mm_set1_pd(0.5/fr->ic->tabq_scale);
133 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
134 rcutoff_scalar = fr->rcoulomb;
135 rcutoff = _mm_set1_pd(rcutoff_scalar);
136 rcutoff2 = _mm_mul_pd(rcutoff,rcutoff);
138 sh_vdw_invrcut6 = _mm_set1_pd(fr->ic->sh_invrc6);
139 rvdw = _mm_set1_pd(fr->rvdw);
141 /* Avoid stupid compiler warnings */
149 /* Start outer loop over neighborlists */
150 for(iidx=0; iidx<nri; iidx++)
152 /* Load shift vector for this list */
153 i_shift_offset = DIM*shiftidx[iidx];
155 /* Load limits for loop over neighbors */
156 j_index_start = jindex[iidx];
157 j_index_end = jindex[iidx+1];
159 /* Get outer coordinate index */
161 i_coord_offset = DIM*inr;
163 /* Load i particle coords and add shift vector */
164 gmx_mm_load_shift_and_1rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
166 fix0 = _mm_setzero_pd();
167 fiy0 = _mm_setzero_pd();
168 fiz0 = _mm_setzero_pd();
170 /* Load parameters for i particles */
171 iq0 = _mm_mul_pd(facel,_mm_load1_pd(charge+inr+0));
172 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
174 /* Reset potential sums */
175 velecsum = _mm_setzero_pd();
176 vvdwsum = _mm_setzero_pd();
178 /* Start inner kernel loop */
179 for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
182 /* Get j neighbor index, and coordinate index */
185 j_coord_offsetA = DIM*jnrA;
186 j_coord_offsetB = DIM*jnrB;
188 /* load j atom coordinates */
189 gmx_mm_load_1rvec_2ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
192 /* Calculate displacement vector */
193 dx00 = _mm_sub_pd(ix0,jx0);
194 dy00 = _mm_sub_pd(iy0,jy0);
195 dz00 = _mm_sub_pd(iz0,jz0);
197 /* Calculate squared distance and things based on it */
198 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
200 rinv00 = gmx_mm_invsqrt_pd(rsq00);
202 rinvsq00 = _mm_mul_pd(rinv00,rinv00);
204 /* Load parameters for j particles */
205 jq0 = gmx_mm_load_2real_swizzle_pd(charge+jnrA+0,charge+jnrB+0);
206 vdwjidx0A = 2*vdwtype[jnrA+0];
207 vdwjidx0B = 2*vdwtype[jnrB+0];
209 /**************************
210 * CALCULATE INTERACTIONS *
211 **************************/
213 if (gmx_mm_any_lt(rsq00,rcutoff2))
216 r00 = _mm_mul_pd(rsq00,rinv00);
218 /* Compute parameters for interactions between i and j atoms */
219 qq00 = _mm_mul_pd(iq0,jq0);
220 gmx_mm_load_2pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,
221 vdwparam+vdwioffset0+vdwjidx0B,&c6_00,&c12_00);
222 c6grid_00 = gmx_mm_load_2real_swizzle_pd(vdwgridparam+vdwioffset0+vdwjidx0A,
223 vdwgridparam+vdwioffset0+vdwjidx0B);
225 /* EWALD ELECTROSTATICS */
227 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
228 ewrt = _mm_mul_pd(r00,ewtabscale);
229 ewitab = _mm_cvttpd_epi32(ewrt);
231 eweps = _mm_frcz_pd(ewrt);
233 eweps = _mm_sub_pd(ewrt,_mm_round_pd(ewrt, _MM_FROUND_FLOOR));
235 twoeweps = _mm_add_pd(eweps,eweps);
236 ewitab = _mm_slli_epi32(ewitab,2);
237 ewtabF = _mm_load_pd( ewtab + _mm_extract_epi32(ewitab,0) );
238 ewtabD = _mm_load_pd( ewtab + _mm_extract_epi32(ewitab,1) );
239 GMX_MM_TRANSPOSE2_PD(ewtabF,ewtabD);
240 ewtabV = _mm_load_sd( ewtab + _mm_extract_epi32(ewitab,0) +2);
241 ewtabFn = _mm_load_sd( ewtab + _mm_extract_epi32(ewitab,1) +2);
242 GMX_MM_TRANSPOSE2_PD(ewtabV,ewtabFn);
243 felec = _mm_macc_pd(eweps,ewtabD,ewtabF);
244 velec = _mm_nmacc_pd(_mm_mul_pd(ewtabhalfspace,eweps) ,_mm_add_pd(ewtabF,felec), ewtabV);
245 velec = _mm_mul_pd(qq00,_mm_sub_pd(_mm_sub_pd(rinv00,sh_ewald),velec));
246 felec = _mm_mul_pd(_mm_mul_pd(qq00,rinv00),_mm_sub_pd(rinvsq00,felec));
248 /* Analytical LJ-PME */
249 rinvsix = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
250 ewcljrsq = _mm_mul_pd(ewclj2,rsq00);
251 ewclj6 = _mm_mul_pd(ewclj2,_mm_mul_pd(ewclj2,ewclj2));
252 exponent = gmx_simd_exp_d(ewcljrsq);
253 /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
254 poly = _mm_mul_pd(exponent,_mm_macc_pd(_mm_mul_pd(ewcljrsq,ewcljrsq),one_half,_mm_sub_pd(one,ewcljrsq)));
255 /* vvdw6 = [C6 - C6grid * (1-poly)]/r6 */
256 vvdw6 = _mm_mul_pd(_mm_macc_pd(-c6grid_00,_mm_sub_pd(one,poly),c6_00),rinvsix);
257 vvdw12 = _mm_mul_pd(c12_00,_mm_mul_pd(rinvsix,rinvsix));
258 vvdw = _mm_msub_pd(_mm_nmacc_pd(c12_00,_mm_mul_pd(sh_vdw_invrcut6,sh_vdw_invrcut6),vvdw12),one_twelfth,
259 _mm_mul_pd(_mm_sub_pd(vvdw6,_mm_macc_pd(c6grid_00,sh_lj_ewald,_mm_mul_pd(c6_00,sh_vdw_invrcut6))),one_sixth));
260 /* fvdw = vvdw12/r - (vvdw6/r + (C6grid * exponent * beta^6)/r) */
261 fvdw = _mm_mul_pd(_mm_add_pd(vvdw12,_mm_msub_pd(_mm_mul_pd(c6grid_00,one_sixth),_mm_mul_pd(exponent,ewclj6),vvdw6)),rinvsq00);
263 cutoff_mask = _mm_cmplt_pd(rsq00,rcutoff2);
265 /* Update potential sum for this i atom from the interaction with this j atom. */
266 velec = _mm_and_pd(velec,cutoff_mask);
267 velecsum = _mm_add_pd(velecsum,velec);
268 vvdw = _mm_and_pd(vvdw,cutoff_mask);
269 vvdwsum = _mm_add_pd(vvdwsum,vvdw);
271 fscal = _mm_add_pd(felec,fvdw);
273 fscal = _mm_and_pd(fscal,cutoff_mask);
275 /* Update vectorial force */
276 fix0 = _mm_macc_pd(dx00,fscal,fix0);
277 fiy0 = _mm_macc_pd(dy00,fscal,fiy0);
278 fiz0 = _mm_macc_pd(dz00,fscal,fiz0);
280 gmx_mm_decrement_1rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,
281 _mm_mul_pd(dx00,fscal),
282 _mm_mul_pd(dy00,fscal),
283 _mm_mul_pd(dz00,fscal));
287 /* Inner loop uses 78 flops */
294 j_coord_offsetA = DIM*jnrA;
296 /* load j atom coordinates */
297 gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
300 /* Calculate displacement vector */
301 dx00 = _mm_sub_pd(ix0,jx0);
302 dy00 = _mm_sub_pd(iy0,jy0);
303 dz00 = _mm_sub_pd(iz0,jz0);
305 /* Calculate squared distance and things based on it */
306 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
308 rinv00 = gmx_mm_invsqrt_pd(rsq00);
310 rinvsq00 = _mm_mul_pd(rinv00,rinv00);
312 /* Load parameters for j particles */
313 jq0 = _mm_load_sd(charge+jnrA+0);
314 vdwjidx0A = 2*vdwtype[jnrA+0];
316 /**************************
317 * CALCULATE INTERACTIONS *
318 **************************/
320 if (gmx_mm_any_lt(rsq00,rcutoff2))
323 r00 = _mm_mul_pd(rsq00,rinv00);
325 /* Compute parameters for interactions between i and j atoms */
326 qq00 = _mm_mul_pd(iq0,jq0);
327 gmx_mm_load_1pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,&c6_00,&c12_00);
328 c6grid_00 = gmx_mm_load_1real_pd(vdwgridparam+vdwioffset0+vdwjidx0A);
330 /* EWALD ELECTROSTATICS */
332 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
333 ewrt = _mm_mul_pd(r00,ewtabscale);
334 ewitab = _mm_cvttpd_epi32(ewrt);
336 eweps = _mm_frcz_pd(ewrt);
338 eweps = _mm_sub_pd(ewrt,_mm_round_pd(ewrt, _MM_FROUND_FLOOR));
340 twoeweps = _mm_add_pd(eweps,eweps);
341 ewitab = _mm_slli_epi32(ewitab,2);
342 ewtabF = _mm_load_pd( ewtab + _mm_extract_epi32(ewitab,0) );
343 ewtabD = _mm_setzero_pd();
344 GMX_MM_TRANSPOSE2_PD(ewtabF,ewtabD);
345 ewtabV = _mm_load_sd( ewtab + _mm_extract_epi32(ewitab,0) +2);
346 ewtabFn = _mm_setzero_pd();
347 GMX_MM_TRANSPOSE2_PD(ewtabV,ewtabFn);
348 felec = _mm_macc_pd(eweps,ewtabD,ewtabF);
349 velec = _mm_nmacc_pd(_mm_mul_pd(ewtabhalfspace,eweps) ,_mm_add_pd(ewtabF,felec), ewtabV);
350 velec = _mm_mul_pd(qq00,_mm_sub_pd(_mm_sub_pd(rinv00,sh_ewald),velec));
351 felec = _mm_mul_pd(_mm_mul_pd(qq00,rinv00),_mm_sub_pd(rinvsq00,felec));
353 /* Analytical LJ-PME */
354 rinvsix = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
355 ewcljrsq = _mm_mul_pd(ewclj2,rsq00);
356 ewclj6 = _mm_mul_pd(ewclj2,_mm_mul_pd(ewclj2,ewclj2));
357 exponent = gmx_simd_exp_d(ewcljrsq);
358 /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
359 poly = _mm_mul_pd(exponent,_mm_macc_pd(_mm_mul_pd(ewcljrsq,ewcljrsq),one_half,_mm_sub_pd(one,ewcljrsq)));
360 /* vvdw6 = [C6 - C6grid * (1-poly)]/r6 */
361 vvdw6 = _mm_mul_pd(_mm_macc_pd(-c6grid_00,_mm_sub_pd(one,poly),c6_00),rinvsix);
362 vvdw12 = _mm_mul_pd(c12_00,_mm_mul_pd(rinvsix,rinvsix));
363 vvdw = _mm_msub_pd(_mm_nmacc_pd(c12_00,_mm_mul_pd(sh_vdw_invrcut6,sh_vdw_invrcut6),vvdw12),one_twelfth,
364 _mm_mul_pd(_mm_sub_pd(vvdw6,_mm_macc_pd(c6grid_00,sh_lj_ewald,_mm_mul_pd(c6_00,sh_vdw_invrcut6))),one_sixth));
365 /* fvdw = vvdw12/r - (vvdw6/r + (C6grid * exponent * beta^6)/r) */
366 fvdw = _mm_mul_pd(_mm_add_pd(vvdw12,_mm_msub_pd(_mm_mul_pd(c6grid_00,one_sixth),_mm_mul_pd(exponent,ewclj6),vvdw6)),rinvsq00);
368 cutoff_mask = _mm_cmplt_pd(rsq00,rcutoff2);
370 /* Update potential sum for this i atom from the interaction with this j atom. */
371 velec = _mm_and_pd(velec,cutoff_mask);
372 velec = _mm_unpacklo_pd(velec,_mm_setzero_pd());
373 velecsum = _mm_add_pd(velecsum,velec);
374 vvdw = _mm_and_pd(vvdw,cutoff_mask);
375 vvdw = _mm_unpacklo_pd(vvdw,_mm_setzero_pd());
376 vvdwsum = _mm_add_pd(vvdwsum,vvdw);
378 fscal = _mm_add_pd(felec,fvdw);
380 fscal = _mm_and_pd(fscal,cutoff_mask);
382 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
384 /* Update vectorial force */
385 fix0 = _mm_macc_pd(dx00,fscal,fix0);
386 fiy0 = _mm_macc_pd(dy00,fscal,fiy0);
387 fiz0 = _mm_macc_pd(dz00,fscal,fiz0);
389 gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,
390 _mm_mul_pd(dx00,fscal),
391 _mm_mul_pd(dy00,fscal),
392 _mm_mul_pd(dz00,fscal));
396 /* Inner loop uses 78 flops */
399 /* End of innermost loop */
401 gmx_mm_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
402 f+i_coord_offset,fshift+i_shift_offset);
405 /* Update potential energies */
406 gmx_mm_update_1pot_pd(velecsum,kernel_data->energygrp_elec+ggid);
407 gmx_mm_update_1pot_pd(vvdwsum,kernel_data->energygrp_vdw+ggid);
409 /* Increment number of inner iterations */
410 inneriter += j_index_end - j_index_start;
412 /* Outer loop uses 9 flops */
415 /* Increment number of outer iterations */
418 /* Update outer/inner flops */
420 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_VF,outeriter*9 + inneriter*78);
423 * Gromacs nonbonded kernel: nb_kernel_ElecEwSh_VdwLJEwSh_GeomP1P1_F_avx_128_fma_double
424 * Electrostatics interaction: Ewald
425 * VdW interaction: LJEwald
426 * Geometry: Particle-Particle
427 * Calculate force/pot: Force
430 nb_kernel_ElecEwSh_VdwLJEwSh_GeomP1P1_F_avx_128_fma_double
431 (t_nblist * gmx_restrict nlist,
432 rvec * gmx_restrict xx,
433 rvec * gmx_restrict ff,
434 t_forcerec * gmx_restrict fr,
435 t_mdatoms * gmx_restrict mdatoms,
436 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
437 t_nrnb * gmx_restrict nrnb)
439 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
440 * just 0 for non-waters.
441 * Suffixes A,B refer to j loop unrolling done with SSE double precision, e.g. for the two different
442 * jnr indices corresponding to data put in the four positions in the SIMD register.
444 int i_shift_offset,i_coord_offset,outeriter,inneriter;
445 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
447 int j_coord_offsetA,j_coord_offsetB;
448 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
450 real *shiftvec,*fshift,*x,*f;
451 __m128d tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
453 __m128d ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
454 int vdwjidx0A,vdwjidx0B;
455 __m128d jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
456 __m128d dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
457 __m128d velec,felec,velecsum,facel,crf,krf,krf2;
460 __m128d rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
463 __m128d one_sixth = _mm_set1_pd(1.0/6.0);
464 __m128d one_twelfth = _mm_set1_pd(1.0/12.0);
467 __m128d ewclj,ewclj2,ewclj6,ewcljrsq,poly,exponent,f6A,f6B,sh_lj_ewald;
468 __m128d one_half = _mm_set1_pd(0.5);
469 __m128d minus_one = _mm_set1_pd(-1.0);
471 __m128d ewtabscale,eweps,twoeweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
473 __m128d dummy_mask,cutoff_mask;
474 __m128d signbit = gmx_mm_castsi128_pd( _mm_set_epi32(0x80000000,0x00000000,0x80000000,0x00000000) );
475 __m128d one = _mm_set1_pd(1.0);
476 __m128d two = _mm_set1_pd(2.0);
482 jindex = nlist->jindex;
484 shiftidx = nlist->shift;
486 shiftvec = fr->shift_vec[0];
487 fshift = fr->fshift[0];
488 facel = _mm_set1_pd(fr->epsfac);
489 charge = mdatoms->chargeA;
490 nvdwtype = fr->ntype;
492 vdwtype = mdatoms->typeA;
493 vdwgridparam = fr->ljpme_c6grid;
494 sh_lj_ewald = _mm_set1_pd(fr->ic->sh_lj_ewald);
495 ewclj = _mm_set1_pd(fr->ewaldcoeff_lj);
496 ewclj2 = _mm_mul_pd(minus_one,_mm_mul_pd(ewclj,ewclj));
498 sh_ewald = _mm_set1_pd(fr->ic->sh_ewald);
499 ewtab = fr->ic->tabq_coul_F;
500 ewtabscale = _mm_set1_pd(fr->ic->tabq_scale);
501 ewtabhalfspace = _mm_set1_pd(0.5/fr->ic->tabq_scale);
503 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
504 rcutoff_scalar = fr->rcoulomb;
505 rcutoff = _mm_set1_pd(rcutoff_scalar);
506 rcutoff2 = _mm_mul_pd(rcutoff,rcutoff);
508 sh_vdw_invrcut6 = _mm_set1_pd(fr->ic->sh_invrc6);
509 rvdw = _mm_set1_pd(fr->rvdw);
511 /* Avoid stupid compiler warnings */
519 /* Start outer loop over neighborlists */
520 for(iidx=0; iidx<nri; iidx++)
522 /* Load shift vector for this list */
523 i_shift_offset = DIM*shiftidx[iidx];
525 /* Load limits for loop over neighbors */
526 j_index_start = jindex[iidx];
527 j_index_end = jindex[iidx+1];
529 /* Get outer coordinate index */
531 i_coord_offset = DIM*inr;
533 /* Load i particle coords and add shift vector */
534 gmx_mm_load_shift_and_1rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
536 fix0 = _mm_setzero_pd();
537 fiy0 = _mm_setzero_pd();
538 fiz0 = _mm_setzero_pd();
540 /* Load parameters for i particles */
541 iq0 = _mm_mul_pd(facel,_mm_load1_pd(charge+inr+0));
542 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
544 /* Start inner kernel loop */
545 for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
548 /* Get j neighbor index, and coordinate index */
551 j_coord_offsetA = DIM*jnrA;
552 j_coord_offsetB = DIM*jnrB;
554 /* load j atom coordinates */
555 gmx_mm_load_1rvec_2ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
558 /* Calculate displacement vector */
559 dx00 = _mm_sub_pd(ix0,jx0);
560 dy00 = _mm_sub_pd(iy0,jy0);
561 dz00 = _mm_sub_pd(iz0,jz0);
563 /* Calculate squared distance and things based on it */
564 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
566 rinv00 = gmx_mm_invsqrt_pd(rsq00);
568 rinvsq00 = _mm_mul_pd(rinv00,rinv00);
570 /* Load parameters for j particles */
571 jq0 = gmx_mm_load_2real_swizzle_pd(charge+jnrA+0,charge+jnrB+0);
572 vdwjidx0A = 2*vdwtype[jnrA+0];
573 vdwjidx0B = 2*vdwtype[jnrB+0];
575 /**************************
576 * CALCULATE INTERACTIONS *
577 **************************/
579 if (gmx_mm_any_lt(rsq00,rcutoff2))
582 r00 = _mm_mul_pd(rsq00,rinv00);
584 /* Compute parameters for interactions between i and j atoms */
585 qq00 = _mm_mul_pd(iq0,jq0);
586 gmx_mm_load_2pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,
587 vdwparam+vdwioffset0+vdwjidx0B,&c6_00,&c12_00);
588 c6grid_00 = gmx_mm_load_2real_swizzle_pd(vdwgridparam+vdwioffset0+vdwjidx0A,
589 vdwgridparam+vdwioffset0+vdwjidx0B);
591 /* EWALD ELECTROSTATICS */
593 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
594 ewrt = _mm_mul_pd(r00,ewtabscale);
595 ewitab = _mm_cvttpd_epi32(ewrt);
597 eweps = _mm_frcz_pd(ewrt);
599 eweps = _mm_sub_pd(ewrt,_mm_round_pd(ewrt, _MM_FROUND_FLOOR));
601 twoeweps = _mm_add_pd(eweps,eweps);
602 gmx_mm_load_2pair_swizzle_pd(ewtab+_mm_extract_epi32(ewitab,0),ewtab+_mm_extract_epi32(ewitab,1),
604 felec = _mm_macc_pd(eweps,ewtabFn,_mm_mul_pd( _mm_sub_pd(one,eweps),ewtabF));
605 felec = _mm_mul_pd(_mm_mul_pd(qq00,rinv00),_mm_sub_pd(rinvsq00,felec));
607 /* Analytical LJ-PME */
608 rinvsix = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
609 ewcljrsq = _mm_mul_pd(ewclj2,rsq00);
610 ewclj6 = _mm_mul_pd(ewclj2,_mm_mul_pd(ewclj2,ewclj2));
611 exponent = gmx_simd_exp_d(ewcljrsq);
612 /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
613 poly = _mm_mul_pd(exponent,_mm_macc_pd(_mm_mul_pd(ewcljrsq,ewcljrsq),one_half,_mm_sub_pd(one,ewcljrsq)));
614 /* f6A = 6 * C6grid * (1 - poly) */
615 f6A = _mm_mul_pd(c6grid_00,_mm_sub_pd(one,poly));
616 /* f6B = C6grid * exponent * beta^6 */
617 f6B = _mm_mul_pd(_mm_mul_pd(c6grid_00,one_sixth),_mm_mul_pd(exponent,ewclj6));
618 /* fvdw = 12*C12/r13 - ((6*C6 - f6A)/r6 + f6B)/r */
619 fvdw = _mm_mul_pd(_mm_macc_pd(_mm_msub_pd(c12_00,rinvsix,_mm_sub_pd(c6_00,f6A)),rinvsix,f6B),rinvsq00);
621 cutoff_mask = _mm_cmplt_pd(rsq00,rcutoff2);
623 fscal = _mm_add_pd(felec,fvdw);
625 fscal = _mm_and_pd(fscal,cutoff_mask);
627 /* Update vectorial force */
628 fix0 = _mm_macc_pd(dx00,fscal,fix0);
629 fiy0 = _mm_macc_pd(dy00,fscal,fiy0);
630 fiz0 = _mm_macc_pd(dz00,fscal,fiz0);
632 gmx_mm_decrement_1rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,
633 _mm_mul_pd(dx00,fscal),
634 _mm_mul_pd(dy00,fscal),
635 _mm_mul_pd(dz00,fscal));
639 /* Inner loop uses 63 flops */
646 j_coord_offsetA = DIM*jnrA;
648 /* load j atom coordinates */
649 gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
652 /* Calculate displacement vector */
653 dx00 = _mm_sub_pd(ix0,jx0);
654 dy00 = _mm_sub_pd(iy0,jy0);
655 dz00 = _mm_sub_pd(iz0,jz0);
657 /* Calculate squared distance and things based on it */
658 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
660 rinv00 = gmx_mm_invsqrt_pd(rsq00);
662 rinvsq00 = _mm_mul_pd(rinv00,rinv00);
664 /* Load parameters for j particles */
665 jq0 = _mm_load_sd(charge+jnrA+0);
666 vdwjidx0A = 2*vdwtype[jnrA+0];
668 /**************************
669 * CALCULATE INTERACTIONS *
670 **************************/
672 if (gmx_mm_any_lt(rsq00,rcutoff2))
675 r00 = _mm_mul_pd(rsq00,rinv00);
677 /* Compute parameters for interactions between i and j atoms */
678 qq00 = _mm_mul_pd(iq0,jq0);
679 gmx_mm_load_1pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,&c6_00,&c12_00);
680 c6grid_00 = gmx_mm_load_1real_pd(vdwgridparam+vdwioffset0+vdwjidx0A);
682 /* EWALD ELECTROSTATICS */
684 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
685 ewrt = _mm_mul_pd(r00,ewtabscale);
686 ewitab = _mm_cvttpd_epi32(ewrt);
688 eweps = _mm_frcz_pd(ewrt);
690 eweps = _mm_sub_pd(ewrt,_mm_round_pd(ewrt, _MM_FROUND_FLOOR));
692 twoeweps = _mm_add_pd(eweps,eweps);
693 gmx_mm_load_1pair_swizzle_pd(ewtab+_mm_extract_epi32(ewitab,0),&ewtabF,&ewtabFn);
694 felec = _mm_macc_pd(eweps,ewtabFn,_mm_mul_pd( _mm_sub_pd(one,eweps),ewtabF));
695 felec = _mm_mul_pd(_mm_mul_pd(qq00,rinv00),_mm_sub_pd(rinvsq00,felec));
697 /* Analytical LJ-PME */
698 rinvsix = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
699 ewcljrsq = _mm_mul_pd(ewclj2,rsq00);
700 ewclj6 = _mm_mul_pd(ewclj2,_mm_mul_pd(ewclj2,ewclj2));
701 exponent = gmx_simd_exp_d(ewcljrsq);
702 /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
703 poly = _mm_mul_pd(exponent,_mm_macc_pd(_mm_mul_pd(ewcljrsq,ewcljrsq),one_half,_mm_sub_pd(one,ewcljrsq)));
704 /* f6A = 6 * C6grid * (1 - poly) */
705 f6A = _mm_mul_pd(c6grid_00,_mm_sub_pd(one,poly));
706 /* f6B = C6grid * exponent * beta^6 */
707 f6B = _mm_mul_pd(_mm_mul_pd(c6grid_00,one_sixth),_mm_mul_pd(exponent,ewclj6));
708 /* fvdw = 12*C12/r13 - ((6*C6 - f6A)/r6 + f6B)/r */
709 fvdw = _mm_mul_pd(_mm_macc_pd(_mm_msub_pd(c12_00,rinvsix,_mm_sub_pd(c6_00,f6A)),rinvsix,f6B),rinvsq00);
711 cutoff_mask = _mm_cmplt_pd(rsq00,rcutoff2);
713 fscal = _mm_add_pd(felec,fvdw);
715 fscal = _mm_and_pd(fscal,cutoff_mask);
717 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
719 /* Update vectorial force */
720 fix0 = _mm_macc_pd(dx00,fscal,fix0);
721 fiy0 = _mm_macc_pd(dy00,fscal,fiy0);
722 fiz0 = _mm_macc_pd(dz00,fscal,fiz0);
724 gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,
725 _mm_mul_pd(dx00,fscal),
726 _mm_mul_pd(dy00,fscal),
727 _mm_mul_pd(dz00,fscal));
731 /* Inner loop uses 63 flops */
734 /* End of innermost loop */
736 gmx_mm_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
737 f+i_coord_offset,fshift+i_shift_offset);
739 /* Increment number of inner iterations */
740 inneriter += j_index_end - j_index_start;
742 /* Outer loop uses 7 flops */
745 /* Increment number of outer iterations */
748 /* Update outer/inner flops */
750 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_F,outeriter*7 + inneriter*63);