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.
42 #include "../nb_kernel.h"
43 #include "gromacs/legacyheaders/types/simple.h"
44 #include "gromacs/math/vec.h"
45 #include "gromacs/legacyheaders/nrnb.h"
47 #include "gromacs/simd/math_x86_avx_128_fma_double.h"
48 #include "kernelutil_x86_avx_128_fma_double.h"
51 * Gromacs nonbonded kernel: nb_kernel_ElecEwSh_VdwLJEwSh_GeomP1P1_VF_avx_128_fma_double
52 * Electrostatics interaction: Ewald
53 * VdW interaction: LJEwald
54 * Geometry: Particle-Particle
55 * Calculate force/pot: PotentialAndForce
58 nb_kernel_ElecEwSh_VdwLJEwSh_GeomP1P1_VF_avx_128_fma_double
59 (t_nblist * gmx_restrict nlist,
60 rvec * gmx_restrict xx,
61 rvec * gmx_restrict ff,
62 t_forcerec * gmx_restrict fr,
63 t_mdatoms * gmx_restrict mdatoms,
64 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
65 t_nrnb * gmx_restrict nrnb)
67 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
68 * just 0 for non-waters.
69 * Suffixes A,B refer to j loop unrolling done with SSE double precision, e.g. for the two different
70 * jnr indices corresponding to data put in the four positions in the SIMD register.
72 int i_shift_offset,i_coord_offset,outeriter,inneriter;
73 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
75 int j_coord_offsetA,j_coord_offsetB;
76 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
78 real *shiftvec,*fshift,*x,*f;
79 __m128d tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
81 __m128d ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
82 int vdwjidx0A,vdwjidx0B;
83 __m128d jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
84 __m128d dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
85 __m128d velec,felec,velecsum,facel,crf,krf,krf2;
88 __m128d rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
91 __m128d one_sixth = _mm_set1_pd(1.0/6.0);
92 __m128d one_twelfth = _mm_set1_pd(1.0/12.0);
95 __m128d ewclj,ewclj2,ewclj6,ewcljrsq,poly,exponent,f6A,f6B,sh_lj_ewald;
96 __m128d one_half = _mm_set1_pd(0.5);
97 __m128d minus_one = _mm_set1_pd(-1.0);
99 __m128d ewtabscale,eweps,twoeweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
101 __m128d dummy_mask,cutoff_mask;
102 __m128d signbit = gmx_mm_castsi128_pd( _mm_set_epi32(0x80000000,0x00000000,0x80000000,0x00000000) );
103 __m128d one = _mm_set1_pd(1.0);
104 __m128d two = _mm_set1_pd(2.0);
110 jindex = nlist->jindex;
112 shiftidx = nlist->shift;
114 shiftvec = fr->shift_vec[0];
115 fshift = fr->fshift[0];
116 facel = _mm_set1_pd(fr->epsfac);
117 charge = mdatoms->chargeA;
118 nvdwtype = fr->ntype;
120 vdwtype = mdatoms->typeA;
121 vdwgridparam = fr->ljpme_c6grid;
122 sh_lj_ewald = _mm_set1_pd(fr->ic->sh_lj_ewald);
123 ewclj = _mm_set1_pd(fr->ewaldcoeff_lj);
124 ewclj2 = _mm_mul_pd(minus_one,_mm_mul_pd(ewclj,ewclj));
126 sh_ewald = _mm_set1_pd(fr->ic->sh_ewald);
127 ewtab = fr->ic->tabq_coul_FDV0;
128 ewtabscale = _mm_set1_pd(fr->ic->tabq_scale);
129 ewtabhalfspace = _mm_set1_pd(0.5/fr->ic->tabq_scale);
131 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
132 rcutoff_scalar = fr->rcoulomb;
133 rcutoff = _mm_set1_pd(rcutoff_scalar);
134 rcutoff2 = _mm_mul_pd(rcutoff,rcutoff);
136 sh_vdw_invrcut6 = _mm_set1_pd(fr->ic->sh_invrc6);
137 rvdw = _mm_set1_pd(fr->rvdw);
139 /* Avoid stupid compiler warnings */
147 /* Start outer loop over neighborlists */
148 for(iidx=0; iidx<nri; iidx++)
150 /* Load shift vector for this list */
151 i_shift_offset = DIM*shiftidx[iidx];
153 /* Load limits for loop over neighbors */
154 j_index_start = jindex[iidx];
155 j_index_end = jindex[iidx+1];
157 /* Get outer coordinate index */
159 i_coord_offset = DIM*inr;
161 /* Load i particle coords and add shift vector */
162 gmx_mm_load_shift_and_1rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
164 fix0 = _mm_setzero_pd();
165 fiy0 = _mm_setzero_pd();
166 fiz0 = _mm_setzero_pd();
168 /* Load parameters for i particles */
169 iq0 = _mm_mul_pd(facel,_mm_load1_pd(charge+inr+0));
170 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
172 /* Reset potential sums */
173 velecsum = _mm_setzero_pd();
174 vvdwsum = _mm_setzero_pd();
176 /* Start inner kernel loop */
177 for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
180 /* Get j neighbor index, and coordinate index */
183 j_coord_offsetA = DIM*jnrA;
184 j_coord_offsetB = DIM*jnrB;
186 /* load j atom coordinates */
187 gmx_mm_load_1rvec_2ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
190 /* Calculate displacement vector */
191 dx00 = _mm_sub_pd(ix0,jx0);
192 dy00 = _mm_sub_pd(iy0,jy0);
193 dz00 = _mm_sub_pd(iz0,jz0);
195 /* Calculate squared distance and things based on it */
196 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
198 rinv00 = gmx_mm_invsqrt_pd(rsq00);
200 rinvsq00 = _mm_mul_pd(rinv00,rinv00);
202 /* Load parameters for j particles */
203 jq0 = gmx_mm_load_2real_swizzle_pd(charge+jnrA+0,charge+jnrB+0);
204 vdwjidx0A = 2*vdwtype[jnrA+0];
205 vdwjidx0B = 2*vdwtype[jnrB+0];
207 /**************************
208 * CALCULATE INTERACTIONS *
209 **************************/
211 if (gmx_mm_any_lt(rsq00,rcutoff2))
214 r00 = _mm_mul_pd(rsq00,rinv00);
216 /* Compute parameters for interactions between i and j atoms */
217 qq00 = _mm_mul_pd(iq0,jq0);
218 gmx_mm_load_2pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,
219 vdwparam+vdwioffset0+vdwjidx0B,&c6_00,&c12_00);
220 c6grid_00 = gmx_mm_load_2real_swizzle_pd(vdwgridparam+vdwioffset0+vdwjidx0A,
221 vdwgridparam+vdwioffset0+vdwjidx0B);
223 /* EWALD ELECTROSTATICS */
225 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
226 ewrt = _mm_mul_pd(r00,ewtabscale);
227 ewitab = _mm_cvttpd_epi32(ewrt);
229 eweps = _mm_frcz_pd(ewrt);
231 eweps = _mm_sub_pd(ewrt,_mm_round_pd(ewrt, _MM_FROUND_FLOOR));
233 twoeweps = _mm_add_pd(eweps,eweps);
234 ewitab = _mm_slli_epi32(ewitab,2);
235 ewtabF = _mm_load_pd( ewtab + _mm_extract_epi32(ewitab,0) );
236 ewtabD = _mm_load_pd( ewtab + _mm_extract_epi32(ewitab,1) );
237 GMX_MM_TRANSPOSE2_PD(ewtabF,ewtabD);
238 ewtabV = _mm_load_sd( ewtab + _mm_extract_epi32(ewitab,0) +2);
239 ewtabFn = _mm_load_sd( ewtab + _mm_extract_epi32(ewitab,1) +2);
240 GMX_MM_TRANSPOSE2_PD(ewtabV,ewtabFn);
241 felec = _mm_macc_pd(eweps,ewtabD,ewtabF);
242 velec = _mm_nmacc_pd(_mm_mul_pd(ewtabhalfspace,eweps) ,_mm_add_pd(ewtabF,felec), ewtabV);
243 velec = _mm_mul_pd(qq00,_mm_sub_pd(_mm_sub_pd(rinv00,sh_ewald),velec));
244 felec = _mm_mul_pd(_mm_mul_pd(qq00,rinv00),_mm_sub_pd(rinvsq00,felec));
246 /* Analytical LJ-PME */
247 rinvsix = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
248 ewcljrsq = _mm_mul_pd(ewclj2,rsq00);
249 ewclj6 = _mm_mul_pd(ewclj2,_mm_mul_pd(ewclj2,ewclj2));
250 exponent = gmx_simd_exp_d(ewcljrsq);
251 /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
252 poly = _mm_mul_pd(exponent,_mm_macc_pd(_mm_mul_pd(ewcljrsq,ewcljrsq),one_half,_mm_sub_pd(one,ewcljrsq)));
253 /* vvdw6 = [C6 - C6grid * (1-poly)]/r6 */
254 vvdw6 = _mm_mul_pd(_mm_macc_pd(-c6grid_00,_mm_sub_pd(one,poly),c6_00),rinvsix);
255 vvdw12 = _mm_mul_pd(c12_00,_mm_mul_pd(rinvsix,rinvsix));
256 vvdw = _mm_msub_pd(_mm_nmacc_pd(c12_00,_mm_mul_pd(sh_vdw_invrcut6,sh_vdw_invrcut6),vvdw12),one_twelfth,
257 _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));
258 /* fvdw = vvdw12/r - (vvdw6/r + (C6grid * exponent * beta^6)/r) */
259 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);
261 cutoff_mask = _mm_cmplt_pd(rsq00,rcutoff2);
263 /* Update potential sum for this i atom from the interaction with this j atom. */
264 velec = _mm_and_pd(velec,cutoff_mask);
265 velecsum = _mm_add_pd(velecsum,velec);
266 vvdw = _mm_and_pd(vvdw,cutoff_mask);
267 vvdwsum = _mm_add_pd(vvdwsum,vvdw);
269 fscal = _mm_add_pd(felec,fvdw);
271 fscal = _mm_and_pd(fscal,cutoff_mask);
273 /* Update vectorial force */
274 fix0 = _mm_macc_pd(dx00,fscal,fix0);
275 fiy0 = _mm_macc_pd(dy00,fscal,fiy0);
276 fiz0 = _mm_macc_pd(dz00,fscal,fiz0);
278 gmx_mm_decrement_1rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,
279 _mm_mul_pd(dx00,fscal),
280 _mm_mul_pd(dy00,fscal),
281 _mm_mul_pd(dz00,fscal));
285 /* Inner loop uses 78 flops */
292 j_coord_offsetA = DIM*jnrA;
294 /* load j atom coordinates */
295 gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
298 /* Calculate displacement vector */
299 dx00 = _mm_sub_pd(ix0,jx0);
300 dy00 = _mm_sub_pd(iy0,jy0);
301 dz00 = _mm_sub_pd(iz0,jz0);
303 /* Calculate squared distance and things based on it */
304 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
306 rinv00 = gmx_mm_invsqrt_pd(rsq00);
308 rinvsq00 = _mm_mul_pd(rinv00,rinv00);
310 /* Load parameters for j particles */
311 jq0 = _mm_load_sd(charge+jnrA+0);
312 vdwjidx0A = 2*vdwtype[jnrA+0];
314 /**************************
315 * CALCULATE INTERACTIONS *
316 **************************/
318 if (gmx_mm_any_lt(rsq00,rcutoff2))
321 r00 = _mm_mul_pd(rsq00,rinv00);
323 /* Compute parameters for interactions between i and j atoms */
324 qq00 = _mm_mul_pd(iq0,jq0);
325 gmx_mm_load_1pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,&c6_00,&c12_00);
326 c6grid_00 = gmx_mm_load_1real_pd(vdwgridparam+vdwioffset0+vdwjidx0A);
328 /* EWALD ELECTROSTATICS */
330 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
331 ewrt = _mm_mul_pd(r00,ewtabscale);
332 ewitab = _mm_cvttpd_epi32(ewrt);
334 eweps = _mm_frcz_pd(ewrt);
336 eweps = _mm_sub_pd(ewrt,_mm_round_pd(ewrt, _MM_FROUND_FLOOR));
338 twoeweps = _mm_add_pd(eweps,eweps);
339 ewitab = _mm_slli_epi32(ewitab,2);
340 ewtabF = _mm_load_pd( ewtab + _mm_extract_epi32(ewitab,0) );
341 ewtabD = _mm_setzero_pd();
342 GMX_MM_TRANSPOSE2_PD(ewtabF,ewtabD);
343 ewtabV = _mm_load_sd( ewtab + _mm_extract_epi32(ewitab,0) +2);
344 ewtabFn = _mm_setzero_pd();
345 GMX_MM_TRANSPOSE2_PD(ewtabV,ewtabFn);
346 felec = _mm_macc_pd(eweps,ewtabD,ewtabF);
347 velec = _mm_nmacc_pd(_mm_mul_pd(ewtabhalfspace,eweps) ,_mm_add_pd(ewtabF,felec), ewtabV);
348 velec = _mm_mul_pd(qq00,_mm_sub_pd(_mm_sub_pd(rinv00,sh_ewald),velec));
349 felec = _mm_mul_pd(_mm_mul_pd(qq00,rinv00),_mm_sub_pd(rinvsq00,felec));
351 /* Analytical LJ-PME */
352 rinvsix = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
353 ewcljrsq = _mm_mul_pd(ewclj2,rsq00);
354 ewclj6 = _mm_mul_pd(ewclj2,_mm_mul_pd(ewclj2,ewclj2));
355 exponent = gmx_simd_exp_d(ewcljrsq);
356 /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
357 poly = _mm_mul_pd(exponent,_mm_macc_pd(_mm_mul_pd(ewcljrsq,ewcljrsq),one_half,_mm_sub_pd(one,ewcljrsq)));
358 /* vvdw6 = [C6 - C6grid * (1-poly)]/r6 */
359 vvdw6 = _mm_mul_pd(_mm_macc_pd(-c6grid_00,_mm_sub_pd(one,poly),c6_00),rinvsix);
360 vvdw12 = _mm_mul_pd(c12_00,_mm_mul_pd(rinvsix,rinvsix));
361 vvdw = _mm_msub_pd(_mm_nmacc_pd(c12_00,_mm_mul_pd(sh_vdw_invrcut6,sh_vdw_invrcut6),vvdw12),one_twelfth,
362 _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));
363 /* fvdw = vvdw12/r - (vvdw6/r + (C6grid * exponent * beta^6)/r) */
364 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);
366 cutoff_mask = _mm_cmplt_pd(rsq00,rcutoff2);
368 /* Update potential sum for this i atom from the interaction with this j atom. */
369 velec = _mm_and_pd(velec,cutoff_mask);
370 velec = _mm_unpacklo_pd(velec,_mm_setzero_pd());
371 velecsum = _mm_add_pd(velecsum,velec);
372 vvdw = _mm_and_pd(vvdw,cutoff_mask);
373 vvdw = _mm_unpacklo_pd(vvdw,_mm_setzero_pd());
374 vvdwsum = _mm_add_pd(vvdwsum,vvdw);
376 fscal = _mm_add_pd(felec,fvdw);
378 fscal = _mm_and_pd(fscal,cutoff_mask);
380 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
382 /* Update vectorial force */
383 fix0 = _mm_macc_pd(dx00,fscal,fix0);
384 fiy0 = _mm_macc_pd(dy00,fscal,fiy0);
385 fiz0 = _mm_macc_pd(dz00,fscal,fiz0);
387 gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,
388 _mm_mul_pd(dx00,fscal),
389 _mm_mul_pd(dy00,fscal),
390 _mm_mul_pd(dz00,fscal));
394 /* Inner loop uses 78 flops */
397 /* End of innermost loop */
399 gmx_mm_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
400 f+i_coord_offset,fshift+i_shift_offset);
403 /* Update potential energies */
404 gmx_mm_update_1pot_pd(velecsum,kernel_data->energygrp_elec+ggid);
405 gmx_mm_update_1pot_pd(vvdwsum,kernel_data->energygrp_vdw+ggid);
407 /* Increment number of inner iterations */
408 inneriter += j_index_end - j_index_start;
410 /* Outer loop uses 9 flops */
413 /* Increment number of outer iterations */
416 /* Update outer/inner flops */
418 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_VF,outeriter*9 + inneriter*78);
421 * Gromacs nonbonded kernel: nb_kernel_ElecEwSh_VdwLJEwSh_GeomP1P1_F_avx_128_fma_double
422 * Electrostatics interaction: Ewald
423 * VdW interaction: LJEwald
424 * Geometry: Particle-Particle
425 * Calculate force/pot: Force
428 nb_kernel_ElecEwSh_VdwLJEwSh_GeomP1P1_F_avx_128_fma_double
429 (t_nblist * gmx_restrict nlist,
430 rvec * gmx_restrict xx,
431 rvec * gmx_restrict ff,
432 t_forcerec * gmx_restrict fr,
433 t_mdatoms * gmx_restrict mdatoms,
434 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
435 t_nrnb * gmx_restrict nrnb)
437 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
438 * just 0 for non-waters.
439 * Suffixes A,B refer to j loop unrolling done with SSE double precision, e.g. for the two different
440 * jnr indices corresponding to data put in the four positions in the SIMD register.
442 int i_shift_offset,i_coord_offset,outeriter,inneriter;
443 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
445 int j_coord_offsetA,j_coord_offsetB;
446 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
448 real *shiftvec,*fshift,*x,*f;
449 __m128d tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
451 __m128d ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
452 int vdwjidx0A,vdwjidx0B;
453 __m128d jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
454 __m128d dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
455 __m128d velec,felec,velecsum,facel,crf,krf,krf2;
458 __m128d rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
461 __m128d one_sixth = _mm_set1_pd(1.0/6.0);
462 __m128d one_twelfth = _mm_set1_pd(1.0/12.0);
465 __m128d ewclj,ewclj2,ewclj6,ewcljrsq,poly,exponent,f6A,f6B,sh_lj_ewald;
466 __m128d one_half = _mm_set1_pd(0.5);
467 __m128d minus_one = _mm_set1_pd(-1.0);
469 __m128d ewtabscale,eweps,twoeweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
471 __m128d dummy_mask,cutoff_mask;
472 __m128d signbit = gmx_mm_castsi128_pd( _mm_set_epi32(0x80000000,0x00000000,0x80000000,0x00000000) );
473 __m128d one = _mm_set1_pd(1.0);
474 __m128d two = _mm_set1_pd(2.0);
480 jindex = nlist->jindex;
482 shiftidx = nlist->shift;
484 shiftvec = fr->shift_vec[0];
485 fshift = fr->fshift[0];
486 facel = _mm_set1_pd(fr->epsfac);
487 charge = mdatoms->chargeA;
488 nvdwtype = fr->ntype;
490 vdwtype = mdatoms->typeA;
491 vdwgridparam = fr->ljpme_c6grid;
492 sh_lj_ewald = _mm_set1_pd(fr->ic->sh_lj_ewald);
493 ewclj = _mm_set1_pd(fr->ewaldcoeff_lj);
494 ewclj2 = _mm_mul_pd(minus_one,_mm_mul_pd(ewclj,ewclj));
496 sh_ewald = _mm_set1_pd(fr->ic->sh_ewald);
497 ewtab = fr->ic->tabq_coul_F;
498 ewtabscale = _mm_set1_pd(fr->ic->tabq_scale);
499 ewtabhalfspace = _mm_set1_pd(0.5/fr->ic->tabq_scale);
501 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
502 rcutoff_scalar = fr->rcoulomb;
503 rcutoff = _mm_set1_pd(rcutoff_scalar);
504 rcutoff2 = _mm_mul_pd(rcutoff,rcutoff);
506 sh_vdw_invrcut6 = _mm_set1_pd(fr->ic->sh_invrc6);
507 rvdw = _mm_set1_pd(fr->rvdw);
509 /* Avoid stupid compiler warnings */
517 /* Start outer loop over neighborlists */
518 for(iidx=0; iidx<nri; iidx++)
520 /* Load shift vector for this list */
521 i_shift_offset = DIM*shiftidx[iidx];
523 /* Load limits for loop over neighbors */
524 j_index_start = jindex[iidx];
525 j_index_end = jindex[iidx+1];
527 /* Get outer coordinate index */
529 i_coord_offset = DIM*inr;
531 /* Load i particle coords and add shift vector */
532 gmx_mm_load_shift_and_1rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
534 fix0 = _mm_setzero_pd();
535 fiy0 = _mm_setzero_pd();
536 fiz0 = _mm_setzero_pd();
538 /* Load parameters for i particles */
539 iq0 = _mm_mul_pd(facel,_mm_load1_pd(charge+inr+0));
540 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
542 /* Start inner kernel loop */
543 for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
546 /* Get j neighbor index, and coordinate index */
549 j_coord_offsetA = DIM*jnrA;
550 j_coord_offsetB = DIM*jnrB;
552 /* load j atom coordinates */
553 gmx_mm_load_1rvec_2ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
556 /* Calculate displacement vector */
557 dx00 = _mm_sub_pd(ix0,jx0);
558 dy00 = _mm_sub_pd(iy0,jy0);
559 dz00 = _mm_sub_pd(iz0,jz0);
561 /* Calculate squared distance and things based on it */
562 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
564 rinv00 = gmx_mm_invsqrt_pd(rsq00);
566 rinvsq00 = _mm_mul_pd(rinv00,rinv00);
568 /* Load parameters for j particles */
569 jq0 = gmx_mm_load_2real_swizzle_pd(charge+jnrA+0,charge+jnrB+0);
570 vdwjidx0A = 2*vdwtype[jnrA+0];
571 vdwjidx0B = 2*vdwtype[jnrB+0];
573 /**************************
574 * CALCULATE INTERACTIONS *
575 **************************/
577 if (gmx_mm_any_lt(rsq00,rcutoff2))
580 r00 = _mm_mul_pd(rsq00,rinv00);
582 /* Compute parameters for interactions between i and j atoms */
583 qq00 = _mm_mul_pd(iq0,jq0);
584 gmx_mm_load_2pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,
585 vdwparam+vdwioffset0+vdwjidx0B,&c6_00,&c12_00);
586 c6grid_00 = gmx_mm_load_2real_swizzle_pd(vdwgridparam+vdwioffset0+vdwjidx0A,
587 vdwgridparam+vdwioffset0+vdwjidx0B);
589 /* EWALD ELECTROSTATICS */
591 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
592 ewrt = _mm_mul_pd(r00,ewtabscale);
593 ewitab = _mm_cvttpd_epi32(ewrt);
595 eweps = _mm_frcz_pd(ewrt);
597 eweps = _mm_sub_pd(ewrt,_mm_round_pd(ewrt, _MM_FROUND_FLOOR));
599 twoeweps = _mm_add_pd(eweps,eweps);
600 gmx_mm_load_2pair_swizzle_pd(ewtab+_mm_extract_epi32(ewitab,0),ewtab+_mm_extract_epi32(ewitab,1),
602 felec = _mm_macc_pd(eweps,ewtabFn,_mm_mul_pd( _mm_sub_pd(one,eweps),ewtabF));
603 felec = _mm_mul_pd(_mm_mul_pd(qq00,rinv00),_mm_sub_pd(rinvsq00,felec));
605 /* Analytical LJ-PME */
606 rinvsix = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
607 ewcljrsq = _mm_mul_pd(ewclj2,rsq00);
608 ewclj6 = _mm_mul_pd(ewclj2,_mm_mul_pd(ewclj2,ewclj2));
609 exponent = gmx_simd_exp_d(ewcljrsq);
610 /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
611 poly = _mm_mul_pd(exponent,_mm_macc_pd(_mm_mul_pd(ewcljrsq,ewcljrsq),one_half,_mm_sub_pd(one,ewcljrsq)));
612 /* f6A = 6 * C6grid * (1 - poly) */
613 f6A = _mm_mul_pd(c6grid_00,_mm_sub_pd(one,poly));
614 /* f6B = C6grid * exponent * beta^6 */
615 f6B = _mm_mul_pd(_mm_mul_pd(c6grid_00,one_sixth),_mm_mul_pd(exponent,ewclj6));
616 /* fvdw = 12*C12/r13 - ((6*C6 - f6A)/r6 + f6B)/r */
617 fvdw = _mm_mul_pd(_mm_macc_pd(_mm_msub_pd(c12_00,rinvsix,_mm_sub_pd(c6_00,f6A)),rinvsix,f6B),rinvsq00);
619 cutoff_mask = _mm_cmplt_pd(rsq00,rcutoff2);
621 fscal = _mm_add_pd(felec,fvdw);
623 fscal = _mm_and_pd(fscal,cutoff_mask);
625 /* Update vectorial force */
626 fix0 = _mm_macc_pd(dx00,fscal,fix0);
627 fiy0 = _mm_macc_pd(dy00,fscal,fiy0);
628 fiz0 = _mm_macc_pd(dz00,fscal,fiz0);
630 gmx_mm_decrement_1rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,
631 _mm_mul_pd(dx00,fscal),
632 _mm_mul_pd(dy00,fscal),
633 _mm_mul_pd(dz00,fscal));
637 /* Inner loop uses 63 flops */
644 j_coord_offsetA = DIM*jnrA;
646 /* load j atom coordinates */
647 gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
650 /* Calculate displacement vector */
651 dx00 = _mm_sub_pd(ix0,jx0);
652 dy00 = _mm_sub_pd(iy0,jy0);
653 dz00 = _mm_sub_pd(iz0,jz0);
655 /* Calculate squared distance and things based on it */
656 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
658 rinv00 = gmx_mm_invsqrt_pd(rsq00);
660 rinvsq00 = _mm_mul_pd(rinv00,rinv00);
662 /* Load parameters for j particles */
663 jq0 = _mm_load_sd(charge+jnrA+0);
664 vdwjidx0A = 2*vdwtype[jnrA+0];
666 /**************************
667 * CALCULATE INTERACTIONS *
668 **************************/
670 if (gmx_mm_any_lt(rsq00,rcutoff2))
673 r00 = _mm_mul_pd(rsq00,rinv00);
675 /* Compute parameters for interactions between i and j atoms */
676 qq00 = _mm_mul_pd(iq0,jq0);
677 gmx_mm_load_1pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,&c6_00,&c12_00);
678 c6grid_00 = gmx_mm_load_1real_pd(vdwgridparam+vdwioffset0+vdwjidx0A);
680 /* EWALD ELECTROSTATICS */
682 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
683 ewrt = _mm_mul_pd(r00,ewtabscale);
684 ewitab = _mm_cvttpd_epi32(ewrt);
686 eweps = _mm_frcz_pd(ewrt);
688 eweps = _mm_sub_pd(ewrt,_mm_round_pd(ewrt, _MM_FROUND_FLOOR));
690 twoeweps = _mm_add_pd(eweps,eweps);
691 gmx_mm_load_1pair_swizzle_pd(ewtab+_mm_extract_epi32(ewitab,0),&ewtabF,&ewtabFn);
692 felec = _mm_macc_pd(eweps,ewtabFn,_mm_mul_pd( _mm_sub_pd(one,eweps),ewtabF));
693 felec = _mm_mul_pd(_mm_mul_pd(qq00,rinv00),_mm_sub_pd(rinvsq00,felec));
695 /* Analytical LJ-PME */
696 rinvsix = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
697 ewcljrsq = _mm_mul_pd(ewclj2,rsq00);
698 ewclj6 = _mm_mul_pd(ewclj2,_mm_mul_pd(ewclj2,ewclj2));
699 exponent = gmx_simd_exp_d(ewcljrsq);
700 /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
701 poly = _mm_mul_pd(exponent,_mm_macc_pd(_mm_mul_pd(ewcljrsq,ewcljrsq),one_half,_mm_sub_pd(one,ewcljrsq)));
702 /* f6A = 6 * C6grid * (1 - poly) */
703 f6A = _mm_mul_pd(c6grid_00,_mm_sub_pd(one,poly));
704 /* f6B = C6grid * exponent * beta^6 */
705 f6B = _mm_mul_pd(_mm_mul_pd(c6grid_00,one_sixth),_mm_mul_pd(exponent,ewclj6));
706 /* fvdw = 12*C12/r13 - ((6*C6 - f6A)/r6 + f6B)/r */
707 fvdw = _mm_mul_pd(_mm_macc_pd(_mm_msub_pd(c12_00,rinvsix,_mm_sub_pd(c6_00,f6A)),rinvsix,f6B),rinvsq00);
709 cutoff_mask = _mm_cmplt_pd(rsq00,rcutoff2);
711 fscal = _mm_add_pd(felec,fvdw);
713 fscal = _mm_and_pd(fscal,cutoff_mask);
715 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
717 /* Update vectorial force */
718 fix0 = _mm_macc_pd(dx00,fscal,fix0);
719 fiy0 = _mm_macc_pd(dy00,fscal,fiy0);
720 fiz0 = _mm_macc_pd(dz00,fscal,fiz0);
722 gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,
723 _mm_mul_pd(dx00,fscal),
724 _mm_mul_pd(dy00,fscal),
725 _mm_mul_pd(dz00,fscal));
729 /* Inner loop uses 63 flops */
732 /* End of innermost loop */
734 gmx_mm_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
735 f+i_coord_offset,fshift+i_shift_offset);
737 /* Increment number of inner iterations */
738 inneriter += j_index_end - j_index_start;
740 /* Outer loop uses 7 flops */
743 /* Increment number of outer iterations */
746 /* Update outer/inner flops */
748 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_F,outeriter*7 + inneriter*63);