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 sse4_1_single kernel generator.
42 #include "../nb_kernel.h"
43 #include "types/simple.h"
44 #include "gromacs/math/vec.h"
47 #include "gromacs/simd/math_x86_sse4_1_single.h"
48 #include "kernelutil_x86_sse4_1_single.h"
51 * Gromacs nonbonded kernel: nb_kernel_ElecEwSh_VdwLJEwSh_GeomP1P1_VF_sse4_1_single
52 * Electrostatics interaction: Ewald
53 * VdW interaction: LJEwald
54 * Geometry: Particle-Particle
55 * Calculate force/pot: PotentialAndForce
58 nb_kernel_ElecEwSh_VdwLJEwSh_GeomP1P1_VF_sse4_1_single
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,C,D refer to j loop unrolling done with SSE, e.g. for the four 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;
74 int jnrA,jnrB,jnrC,jnrD;
75 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
76 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
77 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
79 real *shiftvec,*fshift,*x,*f;
80 real *fjptrA,*fjptrB,*fjptrC,*fjptrD;
82 __m128 tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
84 __m128 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
85 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
86 __m128 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
87 __m128 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
88 __m128 velec,felec,velecsum,facel,crf,krf,krf2;
91 __m128 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
94 __m128 one_sixth = _mm_set1_ps(1.0/6.0);
95 __m128 one_twelfth = _mm_set1_ps(1.0/12.0);
97 __m128 ewclj,ewclj2,ewclj6,ewcljrsq,poly,exponent,f6A,f6B,sh_lj_ewald;
99 __m128 one_half = _mm_set1_ps(0.5);
100 __m128 minus_one = _mm_set1_ps(-1.0);
102 __m128 ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
104 __m128 dummy_mask,cutoff_mask;
105 __m128 signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
106 __m128 one = _mm_set1_ps(1.0);
107 __m128 two = _mm_set1_ps(2.0);
113 jindex = nlist->jindex;
115 shiftidx = nlist->shift;
117 shiftvec = fr->shift_vec[0];
118 fshift = fr->fshift[0];
119 facel = _mm_set1_ps(fr->epsfac);
120 charge = mdatoms->chargeA;
121 nvdwtype = fr->ntype;
123 vdwtype = mdatoms->typeA;
124 vdwgridparam = fr->ljpme_c6grid;
125 sh_lj_ewald = _mm_set1_ps(fr->ic->sh_lj_ewald);
126 ewclj = _mm_set1_ps(fr->ewaldcoeff_lj);
127 ewclj2 = _mm_mul_ps(minus_one,_mm_mul_ps(ewclj,ewclj));
129 sh_ewald = _mm_set1_ps(fr->ic->sh_ewald);
130 ewtab = fr->ic->tabq_coul_FDV0;
131 ewtabscale = _mm_set1_ps(fr->ic->tabq_scale);
132 ewtabhalfspace = _mm_set1_ps(0.5/fr->ic->tabq_scale);
134 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
135 rcutoff_scalar = fr->rcoulomb;
136 rcutoff = _mm_set1_ps(rcutoff_scalar);
137 rcutoff2 = _mm_mul_ps(rcutoff,rcutoff);
139 sh_vdw_invrcut6 = _mm_set1_ps(fr->ic->sh_invrc6);
140 rvdw = _mm_set1_ps(fr->rvdw);
142 /* Avoid stupid compiler warnings */
143 jnrA = jnrB = jnrC = jnrD = 0;
152 for(iidx=0;iidx<4*DIM;iidx++)
157 /* Start outer loop over neighborlists */
158 for(iidx=0; iidx<nri; iidx++)
160 /* Load shift vector for this list */
161 i_shift_offset = DIM*shiftidx[iidx];
163 /* Load limits for loop over neighbors */
164 j_index_start = jindex[iidx];
165 j_index_end = jindex[iidx+1];
167 /* Get outer coordinate index */
169 i_coord_offset = DIM*inr;
171 /* Load i particle coords and add shift vector */
172 gmx_mm_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
174 fix0 = _mm_setzero_ps();
175 fiy0 = _mm_setzero_ps();
176 fiz0 = _mm_setzero_ps();
178 /* Load parameters for i particles */
179 iq0 = _mm_mul_ps(facel,_mm_load1_ps(charge+inr+0));
180 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
182 /* Reset potential sums */
183 velecsum = _mm_setzero_ps();
184 vvdwsum = _mm_setzero_ps();
186 /* Start inner kernel loop */
187 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
190 /* Get j neighbor index, and coordinate index */
195 j_coord_offsetA = DIM*jnrA;
196 j_coord_offsetB = DIM*jnrB;
197 j_coord_offsetC = DIM*jnrC;
198 j_coord_offsetD = DIM*jnrD;
200 /* load j atom coordinates */
201 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
202 x+j_coord_offsetC,x+j_coord_offsetD,
205 /* Calculate displacement vector */
206 dx00 = _mm_sub_ps(ix0,jx0);
207 dy00 = _mm_sub_ps(iy0,jy0);
208 dz00 = _mm_sub_ps(iz0,jz0);
210 /* Calculate squared distance and things based on it */
211 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
213 rinv00 = gmx_mm_invsqrt_ps(rsq00);
215 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
217 /* Load parameters for j particles */
218 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
219 charge+jnrC+0,charge+jnrD+0);
220 vdwjidx0A = 2*vdwtype[jnrA+0];
221 vdwjidx0B = 2*vdwtype[jnrB+0];
222 vdwjidx0C = 2*vdwtype[jnrC+0];
223 vdwjidx0D = 2*vdwtype[jnrD+0];
225 /**************************
226 * CALCULATE INTERACTIONS *
227 **************************/
229 if (gmx_mm_any_lt(rsq00,rcutoff2))
232 r00 = _mm_mul_ps(rsq00,rinv00);
234 /* Compute parameters for interactions between i and j atoms */
235 qq00 = _mm_mul_ps(iq0,jq0);
236 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
237 vdwparam+vdwioffset0+vdwjidx0B,
238 vdwparam+vdwioffset0+vdwjidx0C,
239 vdwparam+vdwioffset0+vdwjidx0D,
242 c6grid_00 = gmx_mm_load_4real_swizzle_ps(vdwgridparam+vdwioffset0+vdwjidx0A,
243 vdwgridparam+vdwioffset0+vdwjidx0B,
244 vdwgridparam+vdwioffset0+vdwjidx0C,
245 vdwgridparam+vdwioffset0+vdwjidx0D);
247 /* EWALD ELECTROSTATICS */
249 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
250 ewrt = _mm_mul_ps(r00,ewtabscale);
251 ewitab = _mm_cvttps_epi32(ewrt);
252 eweps = _mm_sub_ps(ewrt,_mm_round_ps(ewrt, _MM_FROUND_FLOOR));
253 ewitab = _mm_slli_epi32(ewitab,2);
254 ewtabF = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,0) );
255 ewtabD = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,1) );
256 ewtabV = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,2) );
257 ewtabFn = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,3) );
258 _MM_TRANSPOSE4_PS(ewtabF,ewtabD,ewtabV,ewtabFn);
259 felec = _mm_add_ps(ewtabF,_mm_mul_ps(eweps,ewtabD));
260 velec = _mm_sub_ps(ewtabV,_mm_mul_ps(_mm_mul_ps(ewtabhalfspace,eweps),_mm_add_ps(ewtabF,felec)));
261 velec = _mm_mul_ps(qq00,_mm_sub_ps(_mm_sub_ps(rinv00,sh_ewald),velec));
262 felec = _mm_mul_ps(_mm_mul_ps(qq00,rinv00),_mm_sub_ps(rinvsq00,felec));
264 /* Analytical LJ-PME */
265 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
266 ewcljrsq = _mm_mul_ps(ewclj2,rsq00);
267 ewclj6 = _mm_mul_ps(ewclj2,_mm_mul_ps(ewclj2,ewclj2));
268 exponent = gmx_simd_exp_r(ewcljrsq);
269 /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
270 poly = _mm_mul_ps(exponent,_mm_add_ps(_mm_sub_ps(one,ewcljrsq),_mm_mul_ps(_mm_mul_ps(ewcljrsq,ewcljrsq),one_half)));
271 /* vvdw6 = [C6 - C6grid * (1-poly)]/r6 */
272 vvdw6 = _mm_mul_ps(_mm_sub_ps(c6_00,_mm_mul_ps(c6grid_00,_mm_sub_ps(one,poly))),rinvsix);
273 vvdw12 = _mm_mul_ps(c12_00,_mm_mul_ps(rinvsix,rinvsix));
274 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),
275 _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));
276 /* fvdw = vvdw12/r - (vvdw6/r + (C6grid * exponent * beta^6)/r) */
277 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);
279 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
281 /* Update potential sum for this i atom from the interaction with this j atom. */
282 velec = _mm_and_ps(velec,cutoff_mask);
283 velecsum = _mm_add_ps(velecsum,velec);
284 vvdw = _mm_and_ps(vvdw,cutoff_mask);
285 vvdwsum = _mm_add_ps(vvdwsum,vvdw);
287 fscal = _mm_add_ps(felec,fvdw);
289 fscal = _mm_and_ps(fscal,cutoff_mask);
291 /* Calculate temporary vectorial force */
292 tx = _mm_mul_ps(fscal,dx00);
293 ty = _mm_mul_ps(fscal,dy00);
294 tz = _mm_mul_ps(fscal,dz00);
296 /* Update vectorial force */
297 fix0 = _mm_add_ps(fix0,tx);
298 fiy0 = _mm_add_ps(fiy0,ty);
299 fiz0 = _mm_add_ps(fiz0,tz);
301 fjptrA = f+j_coord_offsetA;
302 fjptrB = f+j_coord_offsetB;
303 fjptrC = f+j_coord_offsetC;
304 fjptrD = f+j_coord_offsetD;
305 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
309 /* Inner loop uses 82 flops */
315 /* Get j neighbor index, and coordinate index */
316 jnrlistA = jjnr[jidx];
317 jnrlistB = jjnr[jidx+1];
318 jnrlistC = jjnr[jidx+2];
319 jnrlistD = jjnr[jidx+3];
320 /* Sign of each element will be negative for non-real atoms.
321 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
322 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
324 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
325 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
326 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
327 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
328 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
329 j_coord_offsetA = DIM*jnrA;
330 j_coord_offsetB = DIM*jnrB;
331 j_coord_offsetC = DIM*jnrC;
332 j_coord_offsetD = DIM*jnrD;
334 /* load j atom coordinates */
335 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
336 x+j_coord_offsetC,x+j_coord_offsetD,
339 /* Calculate displacement vector */
340 dx00 = _mm_sub_ps(ix0,jx0);
341 dy00 = _mm_sub_ps(iy0,jy0);
342 dz00 = _mm_sub_ps(iz0,jz0);
344 /* Calculate squared distance and things based on it */
345 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
347 rinv00 = gmx_mm_invsqrt_ps(rsq00);
349 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
351 /* Load parameters for j particles */
352 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
353 charge+jnrC+0,charge+jnrD+0);
354 vdwjidx0A = 2*vdwtype[jnrA+0];
355 vdwjidx0B = 2*vdwtype[jnrB+0];
356 vdwjidx0C = 2*vdwtype[jnrC+0];
357 vdwjidx0D = 2*vdwtype[jnrD+0];
359 /**************************
360 * CALCULATE INTERACTIONS *
361 **************************/
363 if (gmx_mm_any_lt(rsq00,rcutoff2))
366 r00 = _mm_mul_ps(rsq00,rinv00);
367 r00 = _mm_andnot_ps(dummy_mask,r00);
369 /* Compute parameters for interactions between i and j atoms */
370 qq00 = _mm_mul_ps(iq0,jq0);
371 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
372 vdwparam+vdwioffset0+vdwjidx0B,
373 vdwparam+vdwioffset0+vdwjidx0C,
374 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_round_ps(ewrt, _MM_FROUND_FLOOR));
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_sse4_1_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_sse4_1_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,
662 c6grid_00 = gmx_mm_load_4real_swizzle_ps(vdwgridparam+vdwioffset0+vdwjidx0A,
663 vdwgridparam+vdwioffset0+vdwjidx0B,
664 vdwgridparam+vdwioffset0+vdwjidx0C,
665 vdwgridparam+vdwioffset0+vdwjidx0D);
667 /* EWALD ELECTROSTATICS */
669 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
670 ewrt = _mm_mul_ps(r00,ewtabscale);
671 ewitab = _mm_cvttps_epi32(ewrt);
672 eweps = _mm_sub_ps(ewrt,_mm_round_ps(ewrt, _MM_FROUND_FLOOR));
673 gmx_mm_load_4pair_swizzle_ps(ewtab + gmx_mm_extract_epi32(ewitab,0),ewtab + gmx_mm_extract_epi32(ewitab,1),
674 ewtab + gmx_mm_extract_epi32(ewitab,2),ewtab + gmx_mm_extract_epi32(ewitab,3),
676 felec = _mm_add_ps(_mm_mul_ps( _mm_sub_ps(one,eweps),ewtabF),_mm_mul_ps(eweps,ewtabFn));
677 felec = _mm_mul_ps(_mm_mul_ps(qq00,rinv00),_mm_sub_ps(rinvsq00,felec));
679 /* Analytical LJ-PME */
680 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
681 ewcljrsq = _mm_mul_ps(ewclj2,rsq00);
682 ewclj6 = _mm_mul_ps(ewclj2,_mm_mul_ps(ewclj2,ewclj2));
683 exponent = gmx_simd_exp_r(ewcljrsq);
684 /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
685 poly = _mm_mul_ps(exponent,_mm_add_ps(_mm_sub_ps(one,ewcljrsq),_mm_mul_ps(_mm_mul_ps(ewcljrsq,ewcljrsq),one_half)));
686 /* f6A = 6 * C6grid * (1 - poly) */
687 f6A = _mm_mul_ps(c6grid_00,_mm_sub_ps(one,poly));
688 /* f6B = C6grid * exponent * beta^6 */
689 f6B = _mm_mul_ps(_mm_mul_ps(c6grid_00,one_sixth),_mm_mul_ps(exponent,ewclj6));
690 /* fvdw = 12*C12/r13 - ((6*C6 - f6A)/r6 + f6B)/r */
691 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);
693 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
695 fscal = _mm_add_ps(felec,fvdw);
697 fscal = _mm_and_ps(fscal,cutoff_mask);
699 /* Calculate temporary vectorial force */
700 tx = _mm_mul_ps(fscal,dx00);
701 ty = _mm_mul_ps(fscal,dy00);
702 tz = _mm_mul_ps(fscal,dz00);
704 /* Update vectorial force */
705 fix0 = _mm_add_ps(fix0,tx);
706 fiy0 = _mm_add_ps(fiy0,ty);
707 fiz0 = _mm_add_ps(fiz0,tz);
709 fjptrA = f+j_coord_offsetA;
710 fjptrB = f+j_coord_offsetB;
711 fjptrC = f+j_coord_offsetC;
712 fjptrD = f+j_coord_offsetD;
713 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
717 /* Inner loop uses 62 flops */
723 /* Get j neighbor index, and coordinate index */
724 jnrlistA = jjnr[jidx];
725 jnrlistB = jjnr[jidx+1];
726 jnrlistC = jjnr[jidx+2];
727 jnrlistD = jjnr[jidx+3];
728 /* Sign of each element will be negative for non-real atoms.
729 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
730 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
732 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
733 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
734 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
735 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
736 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
737 j_coord_offsetA = DIM*jnrA;
738 j_coord_offsetB = DIM*jnrB;
739 j_coord_offsetC = DIM*jnrC;
740 j_coord_offsetD = DIM*jnrD;
742 /* load j atom coordinates */
743 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
744 x+j_coord_offsetC,x+j_coord_offsetD,
747 /* Calculate displacement vector */
748 dx00 = _mm_sub_ps(ix0,jx0);
749 dy00 = _mm_sub_ps(iy0,jy0);
750 dz00 = _mm_sub_ps(iz0,jz0);
752 /* Calculate squared distance and things based on it */
753 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
755 rinv00 = gmx_mm_invsqrt_ps(rsq00);
757 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
759 /* Load parameters for j particles */
760 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
761 charge+jnrC+0,charge+jnrD+0);
762 vdwjidx0A = 2*vdwtype[jnrA+0];
763 vdwjidx0B = 2*vdwtype[jnrB+0];
764 vdwjidx0C = 2*vdwtype[jnrC+0];
765 vdwjidx0D = 2*vdwtype[jnrD+0];
767 /**************************
768 * CALCULATE INTERACTIONS *
769 **************************/
771 if (gmx_mm_any_lt(rsq00,rcutoff2))
774 r00 = _mm_mul_ps(rsq00,rinv00);
775 r00 = _mm_andnot_ps(dummy_mask,r00);
777 /* Compute parameters for interactions between i and j atoms */
778 qq00 = _mm_mul_ps(iq0,jq0);
779 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
780 vdwparam+vdwioffset0+vdwjidx0B,
781 vdwparam+vdwioffset0+vdwjidx0C,
782 vdwparam+vdwioffset0+vdwjidx0D,
785 c6grid_00 = gmx_mm_load_4real_swizzle_ps(vdwgridparam+vdwioffset0+vdwjidx0A,
786 vdwgridparam+vdwioffset0+vdwjidx0B,
787 vdwgridparam+vdwioffset0+vdwjidx0C,
788 vdwgridparam+vdwioffset0+vdwjidx0D);
790 /* EWALD ELECTROSTATICS */
792 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
793 ewrt = _mm_mul_ps(r00,ewtabscale);
794 ewitab = _mm_cvttps_epi32(ewrt);
795 eweps = _mm_sub_ps(ewrt,_mm_round_ps(ewrt, _MM_FROUND_FLOOR));
796 gmx_mm_load_4pair_swizzle_ps(ewtab + gmx_mm_extract_epi32(ewitab,0),ewtab + gmx_mm_extract_epi32(ewitab,1),
797 ewtab + gmx_mm_extract_epi32(ewitab,2),ewtab + gmx_mm_extract_epi32(ewitab,3),
799 felec = _mm_add_ps(_mm_mul_ps( _mm_sub_ps(one,eweps),ewtabF),_mm_mul_ps(eweps,ewtabFn));
800 felec = _mm_mul_ps(_mm_mul_ps(qq00,rinv00),_mm_sub_ps(rinvsq00,felec));
802 /* Analytical LJ-PME */
803 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
804 ewcljrsq = _mm_mul_ps(ewclj2,rsq00);
805 ewclj6 = _mm_mul_ps(ewclj2,_mm_mul_ps(ewclj2,ewclj2));
806 exponent = gmx_simd_exp_r(ewcljrsq);
807 /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
808 poly = _mm_mul_ps(exponent,_mm_add_ps(_mm_sub_ps(one,ewcljrsq),_mm_mul_ps(_mm_mul_ps(ewcljrsq,ewcljrsq),one_half)));
809 /* f6A = 6 * C6grid * (1 - poly) */
810 f6A = _mm_mul_ps(c6grid_00,_mm_sub_ps(one,poly));
811 /* f6B = C6grid * exponent * beta^6 */
812 f6B = _mm_mul_ps(_mm_mul_ps(c6grid_00,one_sixth),_mm_mul_ps(exponent,ewclj6));
813 /* fvdw = 12*C12/r13 - ((6*C6 - f6A)/r6 + f6B)/r */
814 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);
816 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
818 fscal = _mm_add_ps(felec,fvdw);
820 fscal = _mm_and_ps(fscal,cutoff_mask);
822 fscal = _mm_andnot_ps(dummy_mask,fscal);
824 /* Calculate temporary vectorial force */
825 tx = _mm_mul_ps(fscal,dx00);
826 ty = _mm_mul_ps(fscal,dy00);
827 tz = _mm_mul_ps(fscal,dz00);
829 /* Update vectorial force */
830 fix0 = _mm_add_ps(fix0,tx);
831 fiy0 = _mm_add_ps(fiy0,ty);
832 fiz0 = _mm_add_ps(fiz0,tz);
834 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
835 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
836 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
837 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
838 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
842 /* Inner loop uses 63 flops */
845 /* End of innermost loop */
847 gmx_mm_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
848 f+i_coord_offset,fshift+i_shift_offset);
850 /* Increment number of inner iterations */
851 inneriter += j_index_end - j_index_start;
853 /* Outer loop uses 7 flops */
856 /* Increment number of outer iterations */
859 /* Update outer/inner flops */
861 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_F,outeriter*7 + inneriter*63);