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 sse2_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_sse2_double.h"
48 #include "kernelutil_x86_sse2_double.h"
51 * Gromacs nonbonded kernel: nb_kernel_ElecEwSh_VdwLJEwSh_GeomP1P1_VF_sse2_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_sse2_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);
94 __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,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);
221 c6grid_00 = gmx_mm_load_2real_swizzle_pd(vdwgridparam+vdwioffset0+vdwjidx0A,
222 vdwgridparam+vdwioffset0+vdwjidx0B);
224 /* EWALD ELECTROSTATICS */
226 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
227 ewrt = _mm_mul_pd(r00,ewtabscale);
228 ewitab = _mm_cvttpd_epi32(ewrt);
229 eweps = _mm_sub_pd(ewrt,_mm_cvtepi32_pd(ewitab));
230 ewitab = _mm_slli_epi32(ewitab,2);
231 ewtabF = _mm_load_pd( ewtab + gmx_mm_extract_epi32(ewitab,0) );
232 ewtabD = _mm_load_pd( ewtab + gmx_mm_extract_epi32(ewitab,1) );
233 GMX_MM_TRANSPOSE2_PD(ewtabF,ewtabD);
234 ewtabV = _mm_load_sd( ewtab + gmx_mm_extract_epi32(ewitab,0) +2);
235 ewtabFn = _mm_load_sd( ewtab + gmx_mm_extract_epi32(ewitab,1) +2);
236 GMX_MM_TRANSPOSE2_PD(ewtabV,ewtabFn);
237 felec = _mm_add_pd(ewtabF,_mm_mul_pd(eweps,ewtabD));
238 velec = _mm_sub_pd(ewtabV,_mm_mul_pd(_mm_mul_pd(ewtabhalfspace,eweps),_mm_add_pd(ewtabF,felec)));
239 velec = _mm_mul_pd(qq00,_mm_sub_pd(_mm_sub_pd(rinv00,sh_ewald),velec));
240 felec = _mm_mul_pd(_mm_mul_pd(qq00,rinv00),_mm_sub_pd(rinvsq00,felec));
242 /* Analytical LJ-PME */
243 rinvsix = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
244 ewcljrsq = _mm_mul_pd(ewclj2,rsq00);
245 ewclj6 = _mm_mul_pd(ewclj2,_mm_mul_pd(ewclj2,ewclj2));
246 exponent = gmx_simd_exp_d(ewcljrsq);
247 /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
248 poly = _mm_mul_pd(exponent,_mm_add_pd(_mm_sub_pd(one,ewcljrsq),_mm_mul_pd(_mm_mul_pd(ewcljrsq,ewcljrsq),one_half)));
249 /* vvdw6 = [C6 - C6grid * (1-poly)]/r6 */
250 vvdw6 = _mm_mul_pd(_mm_sub_pd(c6_00,_mm_mul_pd(c6grid_00,_mm_sub_pd(one,poly))),rinvsix);
251 vvdw12 = _mm_mul_pd(c12_00,_mm_mul_pd(rinvsix,rinvsix));
252 vvdw = _mm_sub_pd(_mm_mul_pd( _mm_sub_pd(vvdw12 , _mm_mul_pd(c12_00,_mm_mul_pd(sh_vdw_invrcut6,sh_vdw_invrcut6))),one_twelfth),
253 _mm_mul_pd( _mm_sub_pd(vvdw6,_mm_add_pd(_mm_mul_pd(c6_00,sh_vdw_invrcut6),_mm_mul_pd(c6grid_00,sh_lj_ewald))),one_sixth));
254 /* fvdw = vvdw12/r - (vvdw6/r + (C6grid * exponent * beta^6)/r) */
255 fvdw = _mm_mul_pd(_mm_sub_pd(vvdw12,_mm_sub_pd(vvdw6,_mm_mul_pd(_mm_mul_pd(c6grid_00,one_sixth),_mm_mul_pd(exponent,ewclj6)))),rinvsq00);
257 cutoff_mask = _mm_cmplt_pd(rsq00,rcutoff2);
259 /* Update potential sum for this i atom from the interaction with this j atom. */
260 velec = _mm_and_pd(velec,cutoff_mask);
261 velecsum = _mm_add_pd(velecsum,velec);
262 vvdw = _mm_and_pd(vvdw,cutoff_mask);
263 vvdwsum = _mm_add_pd(vvdwsum,vvdw);
265 fscal = _mm_add_pd(felec,fvdw);
267 fscal = _mm_and_pd(fscal,cutoff_mask);
269 /* Calculate temporary vectorial force */
270 tx = _mm_mul_pd(fscal,dx00);
271 ty = _mm_mul_pd(fscal,dy00);
272 tz = _mm_mul_pd(fscal,dz00);
274 /* Update vectorial force */
275 fix0 = _mm_add_pd(fix0,tx);
276 fiy0 = _mm_add_pd(fiy0,ty);
277 fiz0 = _mm_add_pd(fiz0,tz);
279 gmx_mm_decrement_1rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,tx,ty,tz);
283 /* Inner loop uses 82 flops */
290 j_coord_offsetA = DIM*jnrA;
292 /* load j atom coordinates */
293 gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
296 /* Calculate displacement vector */
297 dx00 = _mm_sub_pd(ix0,jx0);
298 dy00 = _mm_sub_pd(iy0,jy0);
299 dz00 = _mm_sub_pd(iz0,jz0);
301 /* Calculate squared distance and things based on it */
302 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
304 rinv00 = gmx_mm_invsqrt_pd(rsq00);
306 rinvsq00 = _mm_mul_pd(rinv00,rinv00);
308 /* Load parameters for j particles */
309 jq0 = _mm_load_sd(charge+jnrA+0);
310 vdwjidx0A = 2*vdwtype[jnrA+0];
312 /**************************
313 * CALCULATE INTERACTIONS *
314 **************************/
316 if (gmx_mm_any_lt(rsq00,rcutoff2))
319 r00 = _mm_mul_pd(rsq00,rinv00);
321 /* Compute parameters for interactions between i and j atoms */
322 qq00 = _mm_mul_pd(iq0,jq0);
323 gmx_mm_load_1pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,&c6_00,&c12_00);
325 c6grid_00 = gmx_mm_load_1real_pd(vdwgridparam+vdwioffset0+vdwjidx0A);
327 /* EWALD ELECTROSTATICS */
329 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
330 ewrt = _mm_mul_pd(r00,ewtabscale);
331 ewitab = _mm_cvttpd_epi32(ewrt);
332 eweps = _mm_sub_pd(ewrt,_mm_cvtepi32_pd(ewitab));
333 ewitab = _mm_slli_epi32(ewitab,2);
334 ewtabF = _mm_load_pd( ewtab + gmx_mm_extract_epi32(ewitab,0) );
335 ewtabD = _mm_setzero_pd();
336 GMX_MM_TRANSPOSE2_PD(ewtabF,ewtabD);
337 ewtabV = _mm_load_sd( ewtab + gmx_mm_extract_epi32(ewitab,0) +2);
338 ewtabFn = _mm_setzero_pd();
339 GMX_MM_TRANSPOSE2_PD(ewtabV,ewtabFn);
340 felec = _mm_add_pd(ewtabF,_mm_mul_pd(eweps,ewtabD));
341 velec = _mm_sub_pd(ewtabV,_mm_mul_pd(_mm_mul_pd(ewtabhalfspace,eweps),_mm_add_pd(ewtabF,felec)));
342 velec = _mm_mul_pd(qq00,_mm_sub_pd(_mm_sub_pd(rinv00,sh_ewald),velec));
343 felec = _mm_mul_pd(_mm_mul_pd(qq00,rinv00),_mm_sub_pd(rinvsq00,felec));
345 /* Analytical LJ-PME */
346 rinvsix = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
347 ewcljrsq = _mm_mul_pd(ewclj2,rsq00);
348 ewclj6 = _mm_mul_pd(ewclj2,_mm_mul_pd(ewclj2,ewclj2));
349 exponent = gmx_simd_exp_d(ewcljrsq);
350 /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
351 poly = _mm_mul_pd(exponent,_mm_add_pd(_mm_sub_pd(one,ewcljrsq),_mm_mul_pd(_mm_mul_pd(ewcljrsq,ewcljrsq),one_half)));
352 /* vvdw6 = [C6 - C6grid * (1-poly)]/r6 */
353 vvdw6 = _mm_mul_pd(_mm_sub_pd(c6_00,_mm_mul_pd(c6grid_00,_mm_sub_pd(one,poly))),rinvsix);
354 vvdw12 = _mm_mul_pd(c12_00,_mm_mul_pd(rinvsix,rinvsix));
355 vvdw = _mm_sub_pd(_mm_mul_pd( _mm_sub_pd(vvdw12 , _mm_mul_pd(c12_00,_mm_mul_pd(sh_vdw_invrcut6,sh_vdw_invrcut6))),one_twelfth),
356 _mm_mul_pd( _mm_sub_pd(vvdw6,_mm_add_pd(_mm_mul_pd(c6_00,sh_vdw_invrcut6),_mm_mul_pd(c6grid_00,sh_lj_ewald))),one_sixth));
357 /* fvdw = vvdw12/r - (vvdw6/r + (C6grid * exponent * beta^6)/r) */
358 fvdw = _mm_mul_pd(_mm_sub_pd(vvdw12,_mm_sub_pd(vvdw6,_mm_mul_pd(_mm_mul_pd(c6grid_00,one_sixth),_mm_mul_pd(exponent,ewclj6)))),rinvsq00);
360 cutoff_mask = _mm_cmplt_pd(rsq00,rcutoff2);
362 /* Update potential sum for this i atom from the interaction with this j atom. */
363 velec = _mm_and_pd(velec,cutoff_mask);
364 velec = _mm_unpacklo_pd(velec,_mm_setzero_pd());
365 velecsum = _mm_add_pd(velecsum,velec);
366 vvdw = _mm_and_pd(vvdw,cutoff_mask);
367 vvdw = _mm_unpacklo_pd(vvdw,_mm_setzero_pd());
368 vvdwsum = _mm_add_pd(vvdwsum,vvdw);
370 fscal = _mm_add_pd(felec,fvdw);
372 fscal = _mm_and_pd(fscal,cutoff_mask);
374 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
376 /* Calculate temporary vectorial force */
377 tx = _mm_mul_pd(fscal,dx00);
378 ty = _mm_mul_pd(fscal,dy00);
379 tz = _mm_mul_pd(fscal,dz00);
381 /* Update vectorial force */
382 fix0 = _mm_add_pd(fix0,tx);
383 fiy0 = _mm_add_pd(fiy0,ty);
384 fiz0 = _mm_add_pd(fiz0,tz);
386 gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,tx,ty,tz);
390 /* Inner loop uses 82 flops */
393 /* End of innermost loop */
395 gmx_mm_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
396 f+i_coord_offset,fshift+i_shift_offset);
399 /* Update potential energies */
400 gmx_mm_update_1pot_pd(velecsum,kernel_data->energygrp_elec+ggid);
401 gmx_mm_update_1pot_pd(vvdwsum,kernel_data->energygrp_vdw+ggid);
403 /* Increment number of inner iterations */
404 inneriter += j_index_end - j_index_start;
406 /* Outer loop uses 9 flops */
409 /* Increment number of outer iterations */
412 /* Update outer/inner flops */
414 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_VF,outeriter*9 + inneriter*82);
417 * Gromacs nonbonded kernel: nb_kernel_ElecEwSh_VdwLJEwSh_GeomP1P1_F_sse2_double
418 * Electrostatics interaction: Ewald
419 * VdW interaction: LJEwald
420 * Geometry: Particle-Particle
421 * Calculate force/pot: Force
424 nb_kernel_ElecEwSh_VdwLJEwSh_GeomP1P1_F_sse2_double
425 (t_nblist * gmx_restrict nlist,
426 rvec * gmx_restrict xx,
427 rvec * gmx_restrict ff,
428 t_forcerec * gmx_restrict fr,
429 t_mdatoms * gmx_restrict mdatoms,
430 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
431 t_nrnb * gmx_restrict nrnb)
433 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
434 * just 0 for non-waters.
435 * Suffixes A,B refer to j loop unrolling done with SSE double precision, e.g. for the two different
436 * jnr indices corresponding to data put in the four positions in the SIMD register.
438 int i_shift_offset,i_coord_offset,outeriter,inneriter;
439 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
441 int j_coord_offsetA,j_coord_offsetB;
442 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
444 real *shiftvec,*fshift,*x,*f;
445 __m128d tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
447 __m128d ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
448 int vdwjidx0A,vdwjidx0B;
449 __m128d jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
450 __m128d dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
451 __m128d velec,felec,velecsum,facel,crf,krf,krf2;
454 __m128d rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
457 __m128d one_sixth = _mm_set1_pd(1.0/6.0);
458 __m128d one_twelfth = _mm_set1_pd(1.0/12.0);
460 __m128d ewclj,ewclj2,ewclj6,ewcljrsq,poly,exponent,f6A,f6B,sh_lj_ewald;
462 __m128d one_half = _mm_set1_pd(0.5);
463 __m128d minus_one = _mm_set1_pd(-1.0);
465 __m128d ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
467 __m128d dummy_mask,cutoff_mask;
468 __m128d signbit = gmx_mm_castsi128_pd( _mm_set_epi32(0x80000000,0x00000000,0x80000000,0x00000000) );
469 __m128d one = _mm_set1_pd(1.0);
470 __m128d two = _mm_set1_pd(2.0);
476 jindex = nlist->jindex;
478 shiftidx = nlist->shift;
480 shiftvec = fr->shift_vec[0];
481 fshift = fr->fshift[0];
482 facel = _mm_set1_pd(fr->epsfac);
483 charge = mdatoms->chargeA;
484 nvdwtype = fr->ntype;
486 vdwtype = mdatoms->typeA;
487 vdwgridparam = fr->ljpme_c6grid;
488 sh_lj_ewald = _mm_set1_pd(fr->ic->sh_lj_ewald);
489 ewclj = _mm_set1_pd(fr->ewaldcoeff_lj);
490 ewclj2 = _mm_mul_pd(minus_one,_mm_mul_pd(ewclj,ewclj));
492 sh_ewald = _mm_set1_pd(fr->ic->sh_ewald);
493 ewtab = fr->ic->tabq_coul_F;
494 ewtabscale = _mm_set1_pd(fr->ic->tabq_scale);
495 ewtabhalfspace = _mm_set1_pd(0.5/fr->ic->tabq_scale);
497 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
498 rcutoff_scalar = fr->rcoulomb;
499 rcutoff = _mm_set1_pd(rcutoff_scalar);
500 rcutoff2 = _mm_mul_pd(rcutoff,rcutoff);
502 sh_vdw_invrcut6 = _mm_set1_pd(fr->ic->sh_invrc6);
503 rvdw = _mm_set1_pd(fr->rvdw);
505 /* Avoid stupid compiler warnings */
513 /* Start outer loop over neighborlists */
514 for(iidx=0; iidx<nri; iidx++)
516 /* Load shift vector for this list */
517 i_shift_offset = DIM*shiftidx[iidx];
519 /* Load limits for loop over neighbors */
520 j_index_start = jindex[iidx];
521 j_index_end = jindex[iidx+1];
523 /* Get outer coordinate index */
525 i_coord_offset = DIM*inr;
527 /* Load i particle coords and add shift vector */
528 gmx_mm_load_shift_and_1rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
530 fix0 = _mm_setzero_pd();
531 fiy0 = _mm_setzero_pd();
532 fiz0 = _mm_setzero_pd();
534 /* Load parameters for i particles */
535 iq0 = _mm_mul_pd(facel,_mm_load1_pd(charge+inr+0));
536 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
538 /* Start inner kernel loop */
539 for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
542 /* Get j neighbor index, and coordinate index */
545 j_coord_offsetA = DIM*jnrA;
546 j_coord_offsetB = DIM*jnrB;
548 /* load j atom coordinates */
549 gmx_mm_load_1rvec_2ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
552 /* Calculate displacement vector */
553 dx00 = _mm_sub_pd(ix0,jx0);
554 dy00 = _mm_sub_pd(iy0,jy0);
555 dz00 = _mm_sub_pd(iz0,jz0);
557 /* Calculate squared distance and things based on it */
558 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
560 rinv00 = gmx_mm_invsqrt_pd(rsq00);
562 rinvsq00 = _mm_mul_pd(rinv00,rinv00);
564 /* Load parameters for j particles */
565 jq0 = gmx_mm_load_2real_swizzle_pd(charge+jnrA+0,charge+jnrB+0);
566 vdwjidx0A = 2*vdwtype[jnrA+0];
567 vdwjidx0B = 2*vdwtype[jnrB+0];
569 /**************************
570 * CALCULATE INTERACTIONS *
571 **************************/
573 if (gmx_mm_any_lt(rsq00,rcutoff2))
576 r00 = _mm_mul_pd(rsq00,rinv00);
578 /* Compute parameters for interactions between i and j atoms */
579 qq00 = _mm_mul_pd(iq0,jq0);
580 gmx_mm_load_2pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,
581 vdwparam+vdwioffset0+vdwjidx0B,&c6_00,&c12_00);
583 c6grid_00 = gmx_mm_load_2real_swizzle_pd(vdwgridparam+vdwioffset0+vdwjidx0A,
584 vdwgridparam+vdwioffset0+vdwjidx0B);
586 /* EWALD ELECTROSTATICS */
588 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
589 ewrt = _mm_mul_pd(r00,ewtabscale);
590 ewitab = _mm_cvttpd_epi32(ewrt);
591 eweps = _mm_sub_pd(ewrt,_mm_cvtepi32_pd(ewitab));
592 gmx_mm_load_2pair_swizzle_pd(ewtab+gmx_mm_extract_epi32(ewitab,0),ewtab+gmx_mm_extract_epi32(ewitab,1),
594 felec = _mm_add_pd(_mm_mul_pd( _mm_sub_pd(one,eweps),ewtabF),_mm_mul_pd(eweps,ewtabFn));
595 felec = _mm_mul_pd(_mm_mul_pd(qq00,rinv00),_mm_sub_pd(rinvsq00,felec));
597 /* Analytical LJ-PME */
598 rinvsix = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
599 ewcljrsq = _mm_mul_pd(ewclj2,rsq00);
600 ewclj6 = _mm_mul_pd(ewclj2,_mm_mul_pd(ewclj2,ewclj2));
601 exponent = gmx_simd_exp_d(ewcljrsq);
602 /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
603 poly = _mm_mul_pd(exponent,_mm_add_pd(_mm_sub_pd(one,ewcljrsq),_mm_mul_pd(_mm_mul_pd(ewcljrsq,ewcljrsq),one_half)));
604 /* f6A = 6 * C6grid * (1 - poly) */
605 f6A = _mm_mul_pd(c6grid_00,_mm_sub_pd(one,poly));
606 /* f6B = C6grid * exponent * beta^6 */
607 f6B = _mm_mul_pd(_mm_mul_pd(c6grid_00,one_sixth),_mm_mul_pd(exponent,ewclj6));
608 /* fvdw = 12*C12/r13 - ((6*C6 - f6A)/r6 + f6B)/r */
609 fvdw = _mm_mul_pd(_mm_add_pd(_mm_mul_pd(_mm_sub_pd(_mm_mul_pd(c12_00,rinvsix),_mm_sub_pd(c6_00,f6A)),rinvsix),f6B),rinvsq00);
611 cutoff_mask = _mm_cmplt_pd(rsq00,rcutoff2);
613 fscal = _mm_add_pd(felec,fvdw);
615 fscal = _mm_and_pd(fscal,cutoff_mask);
617 /* Calculate temporary vectorial force */
618 tx = _mm_mul_pd(fscal,dx00);
619 ty = _mm_mul_pd(fscal,dy00);
620 tz = _mm_mul_pd(fscal,dz00);
622 /* Update vectorial force */
623 fix0 = _mm_add_pd(fix0,tx);
624 fiy0 = _mm_add_pd(fiy0,ty);
625 fiz0 = _mm_add_pd(fiz0,tz);
627 gmx_mm_decrement_1rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,tx,ty,tz);
631 /* Inner loop uses 62 flops */
638 j_coord_offsetA = DIM*jnrA;
640 /* load j atom coordinates */
641 gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
644 /* Calculate displacement vector */
645 dx00 = _mm_sub_pd(ix0,jx0);
646 dy00 = _mm_sub_pd(iy0,jy0);
647 dz00 = _mm_sub_pd(iz0,jz0);
649 /* Calculate squared distance and things based on it */
650 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
652 rinv00 = gmx_mm_invsqrt_pd(rsq00);
654 rinvsq00 = _mm_mul_pd(rinv00,rinv00);
656 /* Load parameters for j particles */
657 jq0 = _mm_load_sd(charge+jnrA+0);
658 vdwjidx0A = 2*vdwtype[jnrA+0];
660 /**************************
661 * CALCULATE INTERACTIONS *
662 **************************/
664 if (gmx_mm_any_lt(rsq00,rcutoff2))
667 r00 = _mm_mul_pd(rsq00,rinv00);
669 /* Compute parameters for interactions between i and j atoms */
670 qq00 = _mm_mul_pd(iq0,jq0);
671 gmx_mm_load_1pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,&c6_00,&c12_00);
673 c6grid_00 = gmx_mm_load_1real_pd(vdwgridparam+vdwioffset0+vdwjidx0A);
675 /* EWALD ELECTROSTATICS */
677 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
678 ewrt = _mm_mul_pd(r00,ewtabscale);
679 ewitab = _mm_cvttpd_epi32(ewrt);
680 eweps = _mm_sub_pd(ewrt,_mm_cvtepi32_pd(ewitab));
681 gmx_mm_load_1pair_swizzle_pd(ewtab+gmx_mm_extract_epi32(ewitab,0),&ewtabF,&ewtabFn);
682 felec = _mm_add_pd(_mm_mul_pd( _mm_sub_pd(one,eweps),ewtabF),_mm_mul_pd(eweps,ewtabFn));
683 felec = _mm_mul_pd(_mm_mul_pd(qq00,rinv00),_mm_sub_pd(rinvsq00,felec));
685 /* Analytical LJ-PME */
686 rinvsix = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
687 ewcljrsq = _mm_mul_pd(ewclj2,rsq00);
688 ewclj6 = _mm_mul_pd(ewclj2,_mm_mul_pd(ewclj2,ewclj2));
689 exponent = gmx_simd_exp_d(ewcljrsq);
690 /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
691 poly = _mm_mul_pd(exponent,_mm_add_pd(_mm_sub_pd(one,ewcljrsq),_mm_mul_pd(_mm_mul_pd(ewcljrsq,ewcljrsq),one_half)));
692 /* f6A = 6 * C6grid * (1 - poly) */
693 f6A = _mm_mul_pd(c6grid_00,_mm_sub_pd(one,poly));
694 /* f6B = C6grid * exponent * beta^6 */
695 f6B = _mm_mul_pd(_mm_mul_pd(c6grid_00,one_sixth),_mm_mul_pd(exponent,ewclj6));
696 /* fvdw = 12*C12/r13 - ((6*C6 - f6A)/r6 + f6B)/r */
697 fvdw = _mm_mul_pd(_mm_add_pd(_mm_mul_pd(_mm_sub_pd(_mm_mul_pd(c12_00,rinvsix),_mm_sub_pd(c6_00,f6A)),rinvsix),f6B),rinvsq00);
699 cutoff_mask = _mm_cmplt_pd(rsq00,rcutoff2);
701 fscal = _mm_add_pd(felec,fvdw);
703 fscal = _mm_and_pd(fscal,cutoff_mask);
705 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
707 /* Calculate temporary vectorial force */
708 tx = _mm_mul_pd(fscal,dx00);
709 ty = _mm_mul_pd(fscal,dy00);
710 tz = _mm_mul_pd(fscal,dz00);
712 /* Update vectorial force */
713 fix0 = _mm_add_pd(fix0,tx);
714 fiy0 = _mm_add_pd(fiy0,ty);
715 fiz0 = _mm_add_pd(fiz0,tz);
717 gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,tx,ty,tz);
721 /* Inner loop uses 62 flops */
724 /* End of innermost loop */
726 gmx_mm_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
727 f+i_coord_offset,fshift+i_shift_offset);
729 /* Increment number of inner iterations */
730 inneriter += j_index_end - j_index_start;
732 /* Outer loop uses 7 flops */
735 /* Increment number of outer iterations */
738 /* Update outer/inner flops */
740 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_F,outeriter*7 + inneriter*62);