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_single 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_sse2_single.h"
50 #include "kernelutil_x86_sse2_single.h"
53 * Gromacs nonbonded kernel: nb_kernel_ElecEwSh_VdwLJEwSh_GeomP1P1_VF_sse2_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_sse2_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 SSE, 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 tx,ty,tz,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);
99 __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,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
106 __m128 dummy_mask,cutoff_mask;
107 __m128 signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
108 __m128 one = _mm_set1_ps(1.0);
109 __m128 two = _mm_set1_ps(2.0);
115 jindex = nlist->jindex;
117 shiftidx = nlist->shift;
119 shiftvec = fr->shift_vec[0];
120 fshift = fr->fshift[0];
121 facel = _mm_set1_ps(fr->epsfac);
122 charge = mdatoms->chargeA;
123 nvdwtype = fr->ntype;
125 vdwtype = mdatoms->typeA;
126 vdwgridparam = fr->ljpme_c6grid;
127 sh_lj_ewald = _mm_set1_ps(fr->ic->sh_lj_ewald);
128 ewclj = _mm_set1_ps(fr->ewaldcoeff_lj);
129 ewclj2 = _mm_mul_ps(minus_one,_mm_mul_ps(ewclj,ewclj));
131 sh_ewald = _mm_set1_ps(fr->ic->sh_ewald);
132 ewtab = fr->ic->tabq_coul_FDV0;
133 ewtabscale = _mm_set1_ps(fr->ic->tabq_scale);
134 ewtabhalfspace = _mm_set1_ps(0.5/fr->ic->tabq_scale);
136 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
137 rcutoff_scalar = fr->rcoulomb;
138 rcutoff = _mm_set1_ps(rcutoff_scalar);
139 rcutoff2 = _mm_mul_ps(rcutoff,rcutoff);
141 sh_vdw_invrcut6 = _mm_set1_ps(fr->ic->sh_invrc6);
142 rvdw = _mm_set1_ps(fr->rvdw);
144 /* Avoid stupid compiler warnings */
145 jnrA = jnrB = jnrC = jnrD = 0;
154 for(iidx=0;iidx<4*DIM;iidx++)
159 /* Start outer loop over neighborlists */
160 for(iidx=0; iidx<nri; iidx++)
162 /* Load shift vector for this list */
163 i_shift_offset = DIM*shiftidx[iidx];
165 /* Load limits for loop over neighbors */
166 j_index_start = jindex[iidx];
167 j_index_end = jindex[iidx+1];
169 /* Get outer coordinate index */
171 i_coord_offset = DIM*inr;
173 /* Load i particle coords and add shift vector */
174 gmx_mm_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
176 fix0 = _mm_setzero_ps();
177 fiy0 = _mm_setzero_ps();
178 fiz0 = _mm_setzero_ps();
180 /* Load parameters for i particles */
181 iq0 = _mm_mul_ps(facel,_mm_load1_ps(charge+inr+0));
182 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
184 /* Reset potential sums */
185 velecsum = _mm_setzero_ps();
186 vvdwsum = _mm_setzero_ps();
188 /* Start inner kernel loop */
189 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
192 /* Get j neighbor index, and coordinate index */
197 j_coord_offsetA = DIM*jnrA;
198 j_coord_offsetB = DIM*jnrB;
199 j_coord_offsetC = DIM*jnrC;
200 j_coord_offsetD = DIM*jnrD;
202 /* load j atom coordinates */
203 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
204 x+j_coord_offsetC,x+j_coord_offsetD,
207 /* Calculate displacement vector */
208 dx00 = _mm_sub_ps(ix0,jx0);
209 dy00 = _mm_sub_ps(iy0,jy0);
210 dz00 = _mm_sub_ps(iz0,jz0);
212 /* Calculate squared distance and things based on it */
213 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
215 rinv00 = gmx_mm_invsqrt_ps(rsq00);
217 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
219 /* Load parameters for j particles */
220 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
221 charge+jnrC+0,charge+jnrD+0);
222 vdwjidx0A = 2*vdwtype[jnrA+0];
223 vdwjidx0B = 2*vdwtype[jnrB+0];
224 vdwjidx0C = 2*vdwtype[jnrC+0];
225 vdwjidx0D = 2*vdwtype[jnrD+0];
227 /**************************
228 * CALCULATE INTERACTIONS *
229 **************************/
231 if (gmx_mm_any_lt(rsq00,rcutoff2))
234 r00 = _mm_mul_ps(rsq00,rinv00);
236 /* Compute parameters for interactions between i and j atoms */
237 qq00 = _mm_mul_ps(iq0,jq0);
238 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
239 vdwparam+vdwioffset0+vdwjidx0B,
240 vdwparam+vdwioffset0+vdwjidx0C,
241 vdwparam+vdwioffset0+vdwjidx0D,
243 c6grid_00 = gmx_mm_load_4real_swizzle_ps(vdwgridparam+vdwioffset0+vdwjidx0A,
244 vdwgridparam+vdwioffset0+vdwjidx0B,
245 vdwgridparam+vdwioffset0+vdwjidx0C,
246 vdwgridparam+vdwioffset0+vdwjidx0D);
248 /* EWALD ELECTROSTATICS */
250 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
251 ewrt = _mm_mul_ps(r00,ewtabscale);
252 ewitab = _mm_cvttps_epi32(ewrt);
253 eweps = _mm_sub_ps(ewrt,_mm_cvtepi32_ps(ewitab));
254 ewitab = _mm_slli_epi32(ewitab,2);
255 ewtabF = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,0) );
256 ewtabD = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,1) );
257 ewtabV = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,2) );
258 ewtabFn = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,3) );
259 _MM_TRANSPOSE4_PS(ewtabF,ewtabD,ewtabV,ewtabFn);
260 felec = _mm_add_ps(ewtabF,_mm_mul_ps(eweps,ewtabD));
261 velec = _mm_sub_ps(ewtabV,_mm_mul_ps(_mm_mul_ps(ewtabhalfspace,eweps),_mm_add_ps(ewtabF,felec)));
262 velec = _mm_mul_ps(qq00,_mm_sub_ps(_mm_sub_ps(rinv00,sh_ewald),velec));
263 felec = _mm_mul_ps(_mm_mul_ps(qq00,rinv00),_mm_sub_ps(rinvsq00,felec));
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_add_ps(_mm_sub_ps(one,ewcljrsq),_mm_mul_ps(_mm_mul_ps(ewcljrsq,ewcljrsq),one_half)));
272 /* vvdw6 = [C6 - C6grid * (1-poly)]/r6 */
273 vvdw6 = _mm_mul_ps(_mm_sub_ps(c6_00,_mm_mul_ps(c6grid_00,_mm_sub_ps(one,poly))),rinvsix);
274 vvdw12 = _mm_mul_ps(c12_00,_mm_mul_ps(rinvsix,rinvsix));
275 vvdw = _mm_sub_ps(_mm_mul_ps( _mm_sub_ps(vvdw12 , _mm_mul_ps(c12_00,_mm_mul_ps(sh_vdw_invrcut6,sh_vdw_invrcut6))), one_twelfth) ,
276 _mm_mul_ps( _mm_sub_ps(vvdw6,_mm_add_ps(_mm_mul_ps(c6_00,sh_vdw_invrcut6),_mm_mul_ps(c6grid_00,sh_lj_ewald))),one_sixth));
277 /* fvdw = vvdw12/r - (vvdw6/r + (C6grid * exponent * beta^6)/r) */
278 fvdw = _mm_mul_ps(_mm_sub_ps(vvdw12,_mm_sub_ps(vvdw6,_mm_mul_ps(_mm_mul_ps(c6grid_00,one_sixth),_mm_mul_ps(exponent,ewclj6)))),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 /* Calculate temporary vectorial force */
293 tx = _mm_mul_ps(fscal,dx00);
294 ty = _mm_mul_ps(fscal,dy00);
295 tz = _mm_mul_ps(fscal,dz00);
297 /* Update vectorial force */
298 fix0 = _mm_add_ps(fix0,tx);
299 fiy0 = _mm_add_ps(fiy0,ty);
300 fiz0 = _mm_add_ps(fiz0,tz);
302 fjptrA = f+j_coord_offsetA;
303 fjptrB = f+j_coord_offsetB;
304 fjptrC = f+j_coord_offsetC;
305 fjptrD = f+j_coord_offsetD;
306 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
310 /* Inner loop uses 82 flops */
316 /* Get j neighbor index, and coordinate index */
317 jnrlistA = jjnr[jidx];
318 jnrlistB = jjnr[jidx+1];
319 jnrlistC = jjnr[jidx+2];
320 jnrlistD = jjnr[jidx+3];
321 /* Sign of each element will be negative for non-real atoms.
322 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
323 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
325 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
326 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
327 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
328 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
329 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
330 j_coord_offsetA = DIM*jnrA;
331 j_coord_offsetB = DIM*jnrB;
332 j_coord_offsetC = DIM*jnrC;
333 j_coord_offsetD = DIM*jnrD;
335 /* load j atom coordinates */
336 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
337 x+j_coord_offsetC,x+j_coord_offsetD,
340 /* Calculate displacement vector */
341 dx00 = _mm_sub_ps(ix0,jx0);
342 dy00 = _mm_sub_ps(iy0,jy0);
343 dz00 = _mm_sub_ps(iz0,jz0);
345 /* Calculate squared distance and things based on it */
346 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
348 rinv00 = gmx_mm_invsqrt_ps(rsq00);
350 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
352 /* Load parameters for j particles */
353 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
354 charge+jnrC+0,charge+jnrD+0);
355 vdwjidx0A = 2*vdwtype[jnrA+0];
356 vdwjidx0B = 2*vdwtype[jnrB+0];
357 vdwjidx0C = 2*vdwtype[jnrC+0];
358 vdwjidx0D = 2*vdwtype[jnrD+0];
360 /**************************
361 * CALCULATE INTERACTIONS *
362 **************************/
364 if (gmx_mm_any_lt(rsq00,rcutoff2))
367 r00 = _mm_mul_ps(rsq00,rinv00);
368 r00 = _mm_andnot_ps(dummy_mask,r00);
370 /* Compute parameters for interactions between i and j atoms */
371 qq00 = _mm_mul_ps(iq0,jq0);
372 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
373 vdwparam+vdwioffset0+vdwjidx0B,
374 vdwparam+vdwioffset0+vdwjidx0C,
375 vdwparam+vdwioffset0+vdwjidx0D,
377 c6grid_00 = gmx_mm_load_4real_swizzle_ps(vdwgridparam+vdwioffset0+vdwjidx0A,
378 vdwgridparam+vdwioffset0+vdwjidx0B,
379 vdwgridparam+vdwioffset0+vdwjidx0C,
380 vdwgridparam+vdwioffset0+vdwjidx0D);
382 /* EWALD ELECTROSTATICS */
384 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
385 ewrt = _mm_mul_ps(r00,ewtabscale);
386 ewitab = _mm_cvttps_epi32(ewrt);
387 eweps = _mm_sub_ps(ewrt,_mm_cvtepi32_ps(ewitab));
388 ewitab = _mm_slli_epi32(ewitab,2);
389 ewtabF = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,0) );
390 ewtabD = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,1) );
391 ewtabV = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,2) );
392 ewtabFn = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,3) );
393 _MM_TRANSPOSE4_PS(ewtabF,ewtabD,ewtabV,ewtabFn);
394 felec = _mm_add_ps(ewtabF,_mm_mul_ps(eweps,ewtabD));
395 velec = _mm_sub_ps(ewtabV,_mm_mul_ps(_mm_mul_ps(ewtabhalfspace,eweps),_mm_add_ps(ewtabF,felec)));
396 velec = _mm_mul_ps(qq00,_mm_sub_ps(_mm_sub_ps(rinv00,sh_ewald),velec));
397 felec = _mm_mul_ps(_mm_mul_ps(qq00,rinv00),_mm_sub_ps(rinvsq00,felec));
399 /* Analytical LJ-PME */
400 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
401 ewcljrsq = _mm_mul_ps(ewclj2,rsq00);
402 ewclj6 = _mm_mul_ps(ewclj2,_mm_mul_ps(ewclj2,ewclj2));
403 exponent = gmx_simd_exp_r(ewcljrsq);
404 /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
405 poly = _mm_mul_ps(exponent,_mm_add_ps(_mm_sub_ps(one,ewcljrsq),_mm_mul_ps(_mm_mul_ps(ewcljrsq,ewcljrsq),one_half)));
406 /* vvdw6 = [C6 - C6grid * (1-poly)]/r6 */
407 vvdw6 = _mm_mul_ps(_mm_sub_ps(c6_00,_mm_mul_ps(c6grid_00,_mm_sub_ps(one,poly))),rinvsix);
408 vvdw12 = _mm_mul_ps(c12_00,_mm_mul_ps(rinvsix,rinvsix));
409 vvdw = _mm_sub_ps(_mm_mul_ps( _mm_sub_ps(vvdw12 , _mm_mul_ps(c12_00,_mm_mul_ps(sh_vdw_invrcut6,sh_vdw_invrcut6))), one_twelfth) ,
410 _mm_mul_ps( _mm_sub_ps(vvdw6,_mm_add_ps(_mm_mul_ps(c6_00,sh_vdw_invrcut6),_mm_mul_ps(c6grid_00,sh_lj_ewald))),one_sixth));
411 /* fvdw = vvdw12/r - (vvdw6/r + (C6grid * exponent * beta^6)/r) */
412 fvdw = _mm_mul_ps(_mm_sub_ps(vvdw12,_mm_sub_ps(vvdw6,_mm_mul_ps(_mm_mul_ps(c6grid_00,one_sixth),_mm_mul_ps(exponent,ewclj6)))),rinvsq00);
414 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
416 /* Update potential sum for this i atom from the interaction with this j atom. */
417 velec = _mm_and_ps(velec,cutoff_mask);
418 velec = _mm_andnot_ps(dummy_mask,velec);
419 velecsum = _mm_add_ps(velecsum,velec);
420 vvdw = _mm_and_ps(vvdw,cutoff_mask);
421 vvdw = _mm_andnot_ps(dummy_mask,vvdw);
422 vvdwsum = _mm_add_ps(vvdwsum,vvdw);
424 fscal = _mm_add_ps(felec,fvdw);
426 fscal = _mm_and_ps(fscal,cutoff_mask);
428 fscal = _mm_andnot_ps(dummy_mask,fscal);
430 /* Calculate temporary vectorial force */
431 tx = _mm_mul_ps(fscal,dx00);
432 ty = _mm_mul_ps(fscal,dy00);
433 tz = _mm_mul_ps(fscal,dz00);
435 /* Update vectorial force */
436 fix0 = _mm_add_ps(fix0,tx);
437 fiy0 = _mm_add_ps(fiy0,ty);
438 fiz0 = _mm_add_ps(fiz0,tz);
440 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
441 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
442 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
443 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
444 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
448 /* Inner loop uses 83 flops */
451 /* End of innermost loop */
453 gmx_mm_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
454 f+i_coord_offset,fshift+i_shift_offset);
457 /* Update potential energies */
458 gmx_mm_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
459 gmx_mm_update_1pot_ps(vvdwsum,kernel_data->energygrp_vdw+ggid);
461 /* Increment number of inner iterations */
462 inneriter += j_index_end - j_index_start;
464 /* Outer loop uses 9 flops */
467 /* Increment number of outer iterations */
470 /* Update outer/inner flops */
472 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_VF,outeriter*9 + inneriter*83);
475 * Gromacs nonbonded kernel: nb_kernel_ElecEwSh_VdwLJEwSh_GeomP1P1_F_sse2_single
476 * Electrostatics interaction: Ewald
477 * VdW interaction: LJEwald
478 * Geometry: Particle-Particle
479 * Calculate force/pot: Force
482 nb_kernel_ElecEwSh_VdwLJEwSh_GeomP1P1_F_sse2_single
483 (t_nblist * gmx_restrict nlist,
484 rvec * gmx_restrict xx,
485 rvec * gmx_restrict ff,
486 t_forcerec * gmx_restrict fr,
487 t_mdatoms * gmx_restrict mdatoms,
488 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
489 t_nrnb * gmx_restrict nrnb)
491 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
492 * just 0 for non-waters.
493 * Suffixes A,B,C,D refer to j loop unrolling done with SSE, e.g. for the four different
494 * jnr indices corresponding to data put in the four positions in the SIMD register.
496 int i_shift_offset,i_coord_offset,outeriter,inneriter;
497 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
498 int jnrA,jnrB,jnrC,jnrD;
499 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
500 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
501 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
503 real *shiftvec,*fshift,*x,*f;
504 real *fjptrA,*fjptrB,*fjptrC,*fjptrD;
506 __m128 tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
508 __m128 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
509 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
510 __m128 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
511 __m128 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
512 __m128 velec,felec,velecsum,facel,crf,krf,krf2;
515 __m128 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
518 __m128 one_sixth = _mm_set1_ps(1.0/6.0);
519 __m128 one_twelfth = _mm_set1_ps(1.0/12.0);
521 __m128 ewclj,ewclj2,ewclj6,ewcljrsq,poly,exponent,f6A,f6B,sh_lj_ewald;
523 __m128 one_half = _mm_set1_ps(0.5);
524 __m128 minus_one = _mm_set1_ps(-1.0);
526 __m128 ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
528 __m128 dummy_mask,cutoff_mask;
529 __m128 signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
530 __m128 one = _mm_set1_ps(1.0);
531 __m128 two = _mm_set1_ps(2.0);
537 jindex = nlist->jindex;
539 shiftidx = nlist->shift;
541 shiftvec = fr->shift_vec[0];
542 fshift = fr->fshift[0];
543 facel = _mm_set1_ps(fr->epsfac);
544 charge = mdatoms->chargeA;
545 nvdwtype = fr->ntype;
547 vdwtype = mdatoms->typeA;
548 vdwgridparam = fr->ljpme_c6grid;
549 sh_lj_ewald = _mm_set1_ps(fr->ic->sh_lj_ewald);
550 ewclj = _mm_set1_ps(fr->ewaldcoeff_lj);
551 ewclj2 = _mm_mul_ps(minus_one,_mm_mul_ps(ewclj,ewclj));
553 sh_ewald = _mm_set1_ps(fr->ic->sh_ewald);
554 ewtab = fr->ic->tabq_coul_F;
555 ewtabscale = _mm_set1_ps(fr->ic->tabq_scale);
556 ewtabhalfspace = _mm_set1_ps(0.5/fr->ic->tabq_scale);
558 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
559 rcutoff_scalar = fr->rcoulomb;
560 rcutoff = _mm_set1_ps(rcutoff_scalar);
561 rcutoff2 = _mm_mul_ps(rcutoff,rcutoff);
563 sh_vdw_invrcut6 = _mm_set1_ps(fr->ic->sh_invrc6);
564 rvdw = _mm_set1_ps(fr->rvdw);
566 /* Avoid stupid compiler warnings */
567 jnrA = jnrB = jnrC = jnrD = 0;
576 for(iidx=0;iidx<4*DIM;iidx++)
581 /* Start outer loop over neighborlists */
582 for(iidx=0; iidx<nri; iidx++)
584 /* Load shift vector for this list */
585 i_shift_offset = DIM*shiftidx[iidx];
587 /* Load limits for loop over neighbors */
588 j_index_start = jindex[iidx];
589 j_index_end = jindex[iidx+1];
591 /* Get outer coordinate index */
593 i_coord_offset = DIM*inr;
595 /* Load i particle coords and add shift vector */
596 gmx_mm_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
598 fix0 = _mm_setzero_ps();
599 fiy0 = _mm_setzero_ps();
600 fiz0 = _mm_setzero_ps();
602 /* Load parameters for i particles */
603 iq0 = _mm_mul_ps(facel,_mm_load1_ps(charge+inr+0));
604 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
606 /* Start inner kernel loop */
607 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
610 /* Get j neighbor index, and coordinate index */
615 j_coord_offsetA = DIM*jnrA;
616 j_coord_offsetB = DIM*jnrB;
617 j_coord_offsetC = DIM*jnrC;
618 j_coord_offsetD = DIM*jnrD;
620 /* load j atom coordinates */
621 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
622 x+j_coord_offsetC,x+j_coord_offsetD,
625 /* Calculate displacement vector */
626 dx00 = _mm_sub_ps(ix0,jx0);
627 dy00 = _mm_sub_ps(iy0,jy0);
628 dz00 = _mm_sub_ps(iz0,jz0);
630 /* Calculate squared distance and things based on it */
631 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
633 rinv00 = gmx_mm_invsqrt_ps(rsq00);
635 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
637 /* Load parameters for j particles */
638 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
639 charge+jnrC+0,charge+jnrD+0);
640 vdwjidx0A = 2*vdwtype[jnrA+0];
641 vdwjidx0B = 2*vdwtype[jnrB+0];
642 vdwjidx0C = 2*vdwtype[jnrC+0];
643 vdwjidx0D = 2*vdwtype[jnrD+0];
645 /**************************
646 * CALCULATE INTERACTIONS *
647 **************************/
649 if (gmx_mm_any_lt(rsq00,rcutoff2))
652 r00 = _mm_mul_ps(rsq00,rinv00);
654 /* Compute parameters for interactions between i and j atoms */
655 qq00 = _mm_mul_ps(iq0,jq0);
656 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
657 vdwparam+vdwioffset0+vdwjidx0B,
658 vdwparam+vdwioffset0+vdwjidx0C,
659 vdwparam+vdwioffset0+vdwjidx0D,
661 c6grid_00 = gmx_mm_load_4real_swizzle_ps(vdwgridparam+vdwioffset0+vdwjidx0A,
662 vdwgridparam+vdwioffset0+vdwjidx0B,
663 vdwgridparam+vdwioffset0+vdwjidx0C,
664 vdwgridparam+vdwioffset0+vdwjidx0D);
666 /* EWALD ELECTROSTATICS */
668 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
669 ewrt = _mm_mul_ps(r00,ewtabscale);
670 ewitab = _mm_cvttps_epi32(ewrt);
671 eweps = _mm_sub_ps(ewrt,_mm_cvtepi32_ps(ewitab));
672 gmx_mm_load_4pair_swizzle_ps(ewtab+gmx_mm_extract_epi32(ewitab,0),ewtab+gmx_mm_extract_epi32(ewitab,1),
673 ewtab+gmx_mm_extract_epi32(ewitab,2),ewtab+gmx_mm_extract_epi32(ewitab,3),
675 felec = _mm_add_ps(_mm_mul_ps( _mm_sub_ps(one,eweps),ewtabF),_mm_mul_ps(eweps,ewtabFn));
676 felec = _mm_mul_ps(_mm_mul_ps(qq00,rinv00),_mm_sub_ps(rinvsq00,felec));
678 /* Analytical LJ-PME */
679 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
680 ewcljrsq = _mm_mul_ps(ewclj2,rsq00);
681 ewclj6 = _mm_mul_ps(ewclj2,_mm_mul_ps(ewclj2,ewclj2));
682 exponent = gmx_simd_exp_r(ewcljrsq);
683 /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
684 poly = _mm_mul_ps(exponent,_mm_add_ps(_mm_sub_ps(one,ewcljrsq),_mm_mul_ps(_mm_mul_ps(ewcljrsq,ewcljrsq),one_half)));
685 /* f6A = 6 * C6grid * (1 - poly) */
686 f6A = _mm_mul_ps(c6grid_00,_mm_sub_ps(one,poly));
687 /* f6B = C6grid * exponent * beta^6 */
688 f6B = _mm_mul_ps(_mm_mul_ps(c6grid_00,one_sixth),_mm_mul_ps(exponent,ewclj6));
689 /* fvdw = 12*C12/r13 - ((6*C6 - f6A)/r6 + f6B)/r */
690 fvdw = _mm_mul_ps(_mm_add_ps(_mm_mul_ps(_mm_sub_ps(_mm_mul_ps(c12_00,rinvsix),_mm_sub_ps(c6_00,f6A)),rinvsix),f6B),rinvsq00);
692 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
694 fscal = _mm_add_ps(felec,fvdw);
696 fscal = _mm_and_ps(fscal,cutoff_mask);
698 /* Calculate temporary vectorial force */
699 tx = _mm_mul_ps(fscal,dx00);
700 ty = _mm_mul_ps(fscal,dy00);
701 tz = _mm_mul_ps(fscal,dz00);
703 /* Update vectorial force */
704 fix0 = _mm_add_ps(fix0,tx);
705 fiy0 = _mm_add_ps(fiy0,ty);
706 fiz0 = _mm_add_ps(fiz0,tz);
708 fjptrA = f+j_coord_offsetA;
709 fjptrB = f+j_coord_offsetB;
710 fjptrC = f+j_coord_offsetC;
711 fjptrD = f+j_coord_offsetD;
712 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
716 /* Inner loop uses 62 flops */
722 /* Get j neighbor index, and coordinate index */
723 jnrlistA = jjnr[jidx];
724 jnrlistB = jjnr[jidx+1];
725 jnrlistC = jjnr[jidx+2];
726 jnrlistD = jjnr[jidx+3];
727 /* Sign of each element will be negative for non-real atoms.
728 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
729 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
731 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
732 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
733 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
734 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
735 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
736 j_coord_offsetA = DIM*jnrA;
737 j_coord_offsetB = DIM*jnrB;
738 j_coord_offsetC = DIM*jnrC;
739 j_coord_offsetD = DIM*jnrD;
741 /* load j atom coordinates */
742 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
743 x+j_coord_offsetC,x+j_coord_offsetD,
746 /* Calculate displacement vector */
747 dx00 = _mm_sub_ps(ix0,jx0);
748 dy00 = _mm_sub_ps(iy0,jy0);
749 dz00 = _mm_sub_ps(iz0,jz0);
751 /* Calculate squared distance and things based on it */
752 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
754 rinv00 = gmx_mm_invsqrt_ps(rsq00);
756 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
758 /* Load parameters for j particles */
759 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
760 charge+jnrC+0,charge+jnrD+0);
761 vdwjidx0A = 2*vdwtype[jnrA+0];
762 vdwjidx0B = 2*vdwtype[jnrB+0];
763 vdwjidx0C = 2*vdwtype[jnrC+0];
764 vdwjidx0D = 2*vdwtype[jnrD+0];
766 /**************************
767 * CALCULATE INTERACTIONS *
768 **************************/
770 if (gmx_mm_any_lt(rsq00,rcutoff2))
773 r00 = _mm_mul_ps(rsq00,rinv00);
774 r00 = _mm_andnot_ps(dummy_mask,r00);
776 /* Compute parameters for interactions between i and j atoms */
777 qq00 = _mm_mul_ps(iq0,jq0);
778 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
779 vdwparam+vdwioffset0+vdwjidx0B,
780 vdwparam+vdwioffset0+vdwjidx0C,
781 vdwparam+vdwioffset0+vdwjidx0D,
783 c6grid_00 = gmx_mm_load_4real_swizzle_ps(vdwgridparam+vdwioffset0+vdwjidx0A,
784 vdwgridparam+vdwioffset0+vdwjidx0B,
785 vdwgridparam+vdwioffset0+vdwjidx0C,
786 vdwgridparam+vdwioffset0+vdwjidx0D);
788 /* EWALD ELECTROSTATICS */
790 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
791 ewrt = _mm_mul_ps(r00,ewtabscale);
792 ewitab = _mm_cvttps_epi32(ewrt);
793 eweps = _mm_sub_ps(ewrt,_mm_cvtepi32_ps(ewitab));
794 gmx_mm_load_4pair_swizzle_ps(ewtab+gmx_mm_extract_epi32(ewitab,0),ewtab+gmx_mm_extract_epi32(ewitab,1),
795 ewtab+gmx_mm_extract_epi32(ewitab,2),ewtab+gmx_mm_extract_epi32(ewitab,3),
797 felec = _mm_add_ps(_mm_mul_ps( _mm_sub_ps(one,eweps),ewtabF),_mm_mul_ps(eweps,ewtabFn));
798 felec = _mm_mul_ps(_mm_mul_ps(qq00,rinv00),_mm_sub_ps(rinvsq00,felec));
800 /* Analytical LJ-PME */
801 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
802 ewcljrsq = _mm_mul_ps(ewclj2,rsq00);
803 ewclj6 = _mm_mul_ps(ewclj2,_mm_mul_ps(ewclj2,ewclj2));
804 exponent = gmx_simd_exp_r(ewcljrsq);
805 /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
806 poly = _mm_mul_ps(exponent,_mm_add_ps(_mm_sub_ps(one,ewcljrsq),_mm_mul_ps(_mm_mul_ps(ewcljrsq,ewcljrsq),one_half)));
807 /* f6A = 6 * C6grid * (1 - poly) */
808 f6A = _mm_mul_ps(c6grid_00,_mm_sub_ps(one,poly));
809 /* f6B = C6grid * exponent * beta^6 */
810 f6B = _mm_mul_ps(_mm_mul_ps(c6grid_00,one_sixth),_mm_mul_ps(exponent,ewclj6));
811 /* fvdw = 12*C12/r13 - ((6*C6 - f6A)/r6 + f6B)/r */
812 fvdw = _mm_mul_ps(_mm_add_ps(_mm_mul_ps(_mm_sub_ps(_mm_mul_ps(c12_00,rinvsix),_mm_sub_ps(c6_00,f6A)),rinvsix),f6B),rinvsq00);
814 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
816 fscal = _mm_add_ps(felec,fvdw);
818 fscal = _mm_and_ps(fscal,cutoff_mask);
820 fscal = _mm_andnot_ps(dummy_mask,fscal);
822 /* Calculate temporary vectorial force */
823 tx = _mm_mul_ps(fscal,dx00);
824 ty = _mm_mul_ps(fscal,dy00);
825 tz = _mm_mul_ps(fscal,dz00);
827 /* Update vectorial force */
828 fix0 = _mm_add_ps(fix0,tx);
829 fiy0 = _mm_add_ps(fiy0,ty);
830 fiz0 = _mm_add_ps(fiz0,tz);
832 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
833 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
834 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
835 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
836 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
840 /* Inner loop uses 63 flops */
843 /* End of innermost loop */
845 gmx_mm_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
846 f+i_coord_offset,fshift+i_shift_offset);
848 /* Increment number of inner iterations */
849 inneriter += j_index_end - j_index_start;
851 /* Outer loop uses 7 flops */
854 /* Increment number of outer iterations */
857 /* Update outer/inner flops */
859 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_F,outeriter*7 + inneriter*63);