Removed simple.h from nb_kernel_sse2_XX
[gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_sse2_single / nb_kernel_ElecEwSh_VdwLJEwSh_GeomP1P1_sse2_single.c
blob0265422c6e8893af4b4830c091e79a9bd8270aba
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015, 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.
38 #include "gmxpre.h"
40 #include "config.h"
42 #include <math.h>
44 #include "../nb_kernel.h"
45 #include "gromacs/math/vec.h"
46 #include "gromacs/legacyheaders/nrnb.h"
48 #include "gromacs/simd/math_x86_sse2_single.h"
49 #include "kernelutil_x86_sse2_single.h"
52 * Gromacs nonbonded kernel: nb_kernel_ElecEwSh_VdwLJEwSh_GeomP1P1_VF_sse2_single
53 * Electrostatics interaction: Ewald
54 * VdW interaction: LJEwald
55 * Geometry: Particle-Particle
56 * Calculate force/pot: PotentialAndForce
58 void
59 nb_kernel_ElecEwSh_VdwLJEwSh_GeomP1P1_VF_sse2_single
60 (t_nblist * gmx_restrict nlist,
61 rvec * gmx_restrict xx,
62 rvec * gmx_restrict ff,
63 t_forcerec * gmx_restrict fr,
64 t_mdatoms * gmx_restrict mdatoms,
65 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
66 t_nrnb * gmx_restrict nrnb)
68 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
69 * just 0 for non-waters.
70 * Suffixes A,B,C,D refer to j loop unrolling done with SSE, e.g. for the four different
71 * jnr indices corresponding to data put in the four positions in the SIMD register.
73 int i_shift_offset,i_coord_offset,outeriter,inneriter;
74 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
75 int jnrA,jnrB,jnrC,jnrD;
76 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
77 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
78 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
79 real rcutoff_scalar;
80 real *shiftvec,*fshift,*x,*f;
81 real *fjptrA,*fjptrB,*fjptrC,*fjptrD;
82 real scratch[4*DIM];
83 __m128 tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
84 int vdwioffset0;
85 __m128 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
86 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
87 __m128 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
88 __m128 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
89 __m128 velec,felec,velecsum,facel,crf,krf,krf2;
90 real *charge;
91 int nvdwtype;
92 __m128 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
93 int *vdwtype;
94 real *vdwparam;
95 __m128 one_sixth = _mm_set1_ps(1.0/6.0);
96 __m128 one_twelfth = _mm_set1_ps(1.0/12.0);
97 __m128 c6grid_00;
98 __m128 ewclj,ewclj2,ewclj6,ewcljrsq,poly,exponent,f6A,f6B,sh_lj_ewald;
99 real *vdwgridparam;
100 __m128 one_half = _mm_set1_ps(0.5);
101 __m128 minus_one = _mm_set1_ps(-1.0);
102 __m128i ewitab;
103 __m128 ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
104 real *ewtab;
105 __m128 dummy_mask,cutoff_mask;
106 __m128 signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
107 __m128 one = _mm_set1_ps(1.0);
108 __m128 two = _mm_set1_ps(2.0);
109 x = xx[0];
110 f = ff[0];
112 nri = nlist->nri;
113 iinr = nlist->iinr;
114 jindex = nlist->jindex;
115 jjnr = nlist->jjnr;
116 shiftidx = nlist->shift;
117 gid = nlist->gid;
118 shiftvec = fr->shift_vec[0];
119 fshift = fr->fshift[0];
120 facel = _mm_set1_ps(fr->epsfac);
121 charge = mdatoms->chargeA;
122 nvdwtype = fr->ntype;
123 vdwparam = fr->nbfp;
124 vdwtype = mdatoms->typeA;
125 vdwgridparam = fr->ljpme_c6grid;
126 sh_lj_ewald = _mm_set1_ps(fr->ic->sh_lj_ewald);
127 ewclj = _mm_set1_ps(fr->ewaldcoeff_lj);
128 ewclj2 = _mm_mul_ps(minus_one,_mm_mul_ps(ewclj,ewclj));
130 sh_ewald = _mm_set1_ps(fr->ic->sh_ewald);
131 ewtab = fr->ic->tabq_coul_FDV0;
132 ewtabscale = _mm_set1_ps(fr->ic->tabq_scale);
133 ewtabhalfspace = _mm_set1_ps(0.5/fr->ic->tabq_scale);
135 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
136 rcutoff_scalar = fr->rcoulomb;
137 rcutoff = _mm_set1_ps(rcutoff_scalar);
138 rcutoff2 = _mm_mul_ps(rcutoff,rcutoff);
140 sh_vdw_invrcut6 = _mm_set1_ps(fr->ic->sh_invrc6);
141 rvdw = _mm_set1_ps(fr->rvdw);
143 /* Avoid stupid compiler warnings */
144 jnrA = jnrB = jnrC = jnrD = 0;
145 j_coord_offsetA = 0;
146 j_coord_offsetB = 0;
147 j_coord_offsetC = 0;
148 j_coord_offsetD = 0;
150 outeriter = 0;
151 inneriter = 0;
153 for(iidx=0;iidx<4*DIM;iidx++)
155 scratch[iidx] = 0.0;
158 /* Start outer loop over neighborlists */
159 for(iidx=0; iidx<nri; iidx++)
161 /* Load shift vector for this list */
162 i_shift_offset = DIM*shiftidx[iidx];
164 /* Load limits for loop over neighbors */
165 j_index_start = jindex[iidx];
166 j_index_end = jindex[iidx+1];
168 /* Get outer coordinate index */
169 inr = iinr[iidx];
170 i_coord_offset = DIM*inr;
172 /* Load i particle coords and add shift vector */
173 gmx_mm_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
175 fix0 = _mm_setzero_ps();
176 fiy0 = _mm_setzero_ps();
177 fiz0 = _mm_setzero_ps();
179 /* Load parameters for i particles */
180 iq0 = _mm_mul_ps(facel,_mm_load1_ps(charge+inr+0));
181 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
183 /* Reset potential sums */
184 velecsum = _mm_setzero_ps();
185 vvdwsum = _mm_setzero_ps();
187 /* Start inner kernel loop */
188 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
191 /* Get j neighbor index, and coordinate index */
192 jnrA = jjnr[jidx];
193 jnrB = jjnr[jidx+1];
194 jnrC = jjnr[jidx+2];
195 jnrD = jjnr[jidx+3];
196 j_coord_offsetA = DIM*jnrA;
197 j_coord_offsetB = DIM*jnrB;
198 j_coord_offsetC = DIM*jnrC;
199 j_coord_offsetD = DIM*jnrD;
201 /* load j atom coordinates */
202 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
203 x+j_coord_offsetC,x+j_coord_offsetD,
204 &jx0,&jy0,&jz0);
206 /* Calculate displacement vector */
207 dx00 = _mm_sub_ps(ix0,jx0);
208 dy00 = _mm_sub_ps(iy0,jy0);
209 dz00 = _mm_sub_ps(iz0,jz0);
211 /* Calculate squared distance and things based on it */
212 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
214 rinv00 = gmx_mm_invsqrt_ps(rsq00);
216 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
218 /* Load parameters for j particles */
219 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
220 charge+jnrC+0,charge+jnrD+0);
221 vdwjidx0A = 2*vdwtype[jnrA+0];
222 vdwjidx0B = 2*vdwtype[jnrB+0];
223 vdwjidx0C = 2*vdwtype[jnrC+0];
224 vdwjidx0D = 2*vdwtype[jnrD+0];
226 /**************************
227 * CALCULATE INTERACTIONS *
228 **************************/
230 if (gmx_mm_any_lt(rsq00,rcutoff2))
233 r00 = _mm_mul_ps(rsq00,rinv00);
235 /* Compute parameters for interactions between i and j atoms */
236 qq00 = _mm_mul_ps(iq0,jq0);
237 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
238 vdwparam+vdwioffset0+vdwjidx0B,
239 vdwparam+vdwioffset0+vdwjidx0C,
240 vdwparam+vdwioffset0+vdwjidx0D,
241 &c6_00,&c12_00);
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_cvtepi32_ps(ewitab));
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 */
312 if(jidx<j_index_end)
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,
337 &jx0,&jy0,&jz0);
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,
375 &c6_00,&c12_00);
376 c6grid_00 = gmx_mm_load_4real_swizzle_ps(vdwgridparam+vdwioffset0+vdwjidx0A,
377 vdwgridparam+vdwioffset0+vdwjidx0B,
378 vdwgridparam+vdwioffset0+vdwjidx0C,
379 vdwgridparam+vdwioffset0+vdwjidx0D);
381 /* EWALD ELECTROSTATICS */
383 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
384 ewrt = _mm_mul_ps(r00,ewtabscale);
385 ewitab = _mm_cvttps_epi32(ewrt);
386 eweps = _mm_sub_ps(ewrt,_mm_cvtepi32_ps(ewitab));
387 ewitab = _mm_slli_epi32(ewitab,2);
388 ewtabF = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,0) );
389 ewtabD = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,1) );
390 ewtabV = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,2) );
391 ewtabFn = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,3) );
392 _MM_TRANSPOSE4_PS(ewtabF,ewtabD,ewtabV,ewtabFn);
393 felec = _mm_add_ps(ewtabF,_mm_mul_ps(eweps,ewtabD));
394 velec = _mm_sub_ps(ewtabV,_mm_mul_ps(_mm_mul_ps(ewtabhalfspace,eweps),_mm_add_ps(ewtabF,felec)));
395 velec = _mm_mul_ps(qq00,_mm_sub_ps(_mm_sub_ps(rinv00,sh_ewald),velec));
396 felec = _mm_mul_ps(_mm_mul_ps(qq00,rinv00),_mm_sub_ps(rinvsq00,felec));
398 /* Analytical LJ-PME */
399 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
400 ewcljrsq = _mm_mul_ps(ewclj2,rsq00);
401 ewclj6 = _mm_mul_ps(ewclj2,_mm_mul_ps(ewclj2,ewclj2));
402 exponent = gmx_simd_exp_r(ewcljrsq);
403 /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
404 poly = _mm_mul_ps(exponent,_mm_add_ps(_mm_sub_ps(one,ewcljrsq),_mm_mul_ps(_mm_mul_ps(ewcljrsq,ewcljrsq),one_half)));
405 /* vvdw6 = [C6 - C6grid * (1-poly)]/r6 */
406 vvdw6 = _mm_mul_ps(_mm_sub_ps(c6_00,_mm_mul_ps(c6grid_00,_mm_sub_ps(one,poly))),rinvsix);
407 vvdw12 = _mm_mul_ps(c12_00,_mm_mul_ps(rinvsix,rinvsix));
408 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) ,
409 _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));
410 /* fvdw = vvdw12/r - (vvdw6/r + (C6grid * exponent * beta^6)/r) */
411 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);
413 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
415 /* Update potential sum for this i atom from the interaction with this j atom. */
416 velec = _mm_and_ps(velec,cutoff_mask);
417 velec = _mm_andnot_ps(dummy_mask,velec);
418 velecsum = _mm_add_ps(velecsum,velec);
419 vvdw = _mm_and_ps(vvdw,cutoff_mask);
420 vvdw = _mm_andnot_ps(dummy_mask,vvdw);
421 vvdwsum = _mm_add_ps(vvdwsum,vvdw);
423 fscal = _mm_add_ps(felec,fvdw);
425 fscal = _mm_and_ps(fscal,cutoff_mask);
427 fscal = _mm_andnot_ps(dummy_mask,fscal);
429 /* Calculate temporary vectorial force */
430 tx = _mm_mul_ps(fscal,dx00);
431 ty = _mm_mul_ps(fscal,dy00);
432 tz = _mm_mul_ps(fscal,dz00);
434 /* Update vectorial force */
435 fix0 = _mm_add_ps(fix0,tx);
436 fiy0 = _mm_add_ps(fiy0,ty);
437 fiz0 = _mm_add_ps(fiz0,tz);
439 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
440 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
441 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
442 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
443 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
447 /* Inner loop uses 83 flops */
450 /* End of innermost loop */
452 gmx_mm_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
453 f+i_coord_offset,fshift+i_shift_offset);
455 ggid = gid[iidx];
456 /* Update potential energies */
457 gmx_mm_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
458 gmx_mm_update_1pot_ps(vvdwsum,kernel_data->energygrp_vdw+ggid);
460 /* Increment number of inner iterations */
461 inneriter += j_index_end - j_index_start;
463 /* Outer loop uses 9 flops */
466 /* Increment number of outer iterations */
467 outeriter += nri;
469 /* Update outer/inner flops */
471 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_VF,outeriter*9 + inneriter*83);
474 * Gromacs nonbonded kernel: nb_kernel_ElecEwSh_VdwLJEwSh_GeomP1P1_F_sse2_single
475 * Electrostatics interaction: Ewald
476 * VdW interaction: LJEwald
477 * Geometry: Particle-Particle
478 * Calculate force/pot: Force
480 void
481 nb_kernel_ElecEwSh_VdwLJEwSh_GeomP1P1_F_sse2_single
482 (t_nblist * gmx_restrict nlist,
483 rvec * gmx_restrict xx,
484 rvec * gmx_restrict ff,
485 t_forcerec * gmx_restrict fr,
486 t_mdatoms * gmx_restrict mdatoms,
487 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
488 t_nrnb * gmx_restrict nrnb)
490 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
491 * just 0 for non-waters.
492 * Suffixes A,B,C,D refer to j loop unrolling done with SSE, e.g. for the four different
493 * jnr indices corresponding to data put in the four positions in the SIMD register.
495 int i_shift_offset,i_coord_offset,outeriter,inneriter;
496 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
497 int jnrA,jnrB,jnrC,jnrD;
498 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
499 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
500 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
501 real rcutoff_scalar;
502 real *shiftvec,*fshift,*x,*f;
503 real *fjptrA,*fjptrB,*fjptrC,*fjptrD;
504 real scratch[4*DIM];
505 __m128 tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
506 int vdwioffset0;
507 __m128 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
508 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
509 __m128 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
510 __m128 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
511 __m128 velec,felec,velecsum,facel,crf,krf,krf2;
512 real *charge;
513 int nvdwtype;
514 __m128 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
515 int *vdwtype;
516 real *vdwparam;
517 __m128 one_sixth = _mm_set1_ps(1.0/6.0);
518 __m128 one_twelfth = _mm_set1_ps(1.0/12.0);
519 __m128 c6grid_00;
520 __m128 ewclj,ewclj2,ewclj6,ewcljrsq,poly,exponent,f6A,f6B,sh_lj_ewald;
521 real *vdwgridparam;
522 __m128 one_half = _mm_set1_ps(0.5);
523 __m128 minus_one = _mm_set1_ps(-1.0);
524 __m128i ewitab;
525 __m128 ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
526 real *ewtab;
527 __m128 dummy_mask,cutoff_mask;
528 __m128 signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
529 __m128 one = _mm_set1_ps(1.0);
530 __m128 two = _mm_set1_ps(2.0);
531 x = xx[0];
532 f = ff[0];
534 nri = nlist->nri;
535 iinr = nlist->iinr;
536 jindex = nlist->jindex;
537 jjnr = nlist->jjnr;
538 shiftidx = nlist->shift;
539 gid = nlist->gid;
540 shiftvec = fr->shift_vec[0];
541 fshift = fr->fshift[0];
542 facel = _mm_set1_ps(fr->epsfac);
543 charge = mdatoms->chargeA;
544 nvdwtype = fr->ntype;
545 vdwparam = fr->nbfp;
546 vdwtype = mdatoms->typeA;
547 vdwgridparam = fr->ljpme_c6grid;
548 sh_lj_ewald = _mm_set1_ps(fr->ic->sh_lj_ewald);
549 ewclj = _mm_set1_ps(fr->ewaldcoeff_lj);
550 ewclj2 = _mm_mul_ps(minus_one,_mm_mul_ps(ewclj,ewclj));
552 sh_ewald = _mm_set1_ps(fr->ic->sh_ewald);
553 ewtab = fr->ic->tabq_coul_F;
554 ewtabscale = _mm_set1_ps(fr->ic->tabq_scale);
555 ewtabhalfspace = _mm_set1_ps(0.5/fr->ic->tabq_scale);
557 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
558 rcutoff_scalar = fr->rcoulomb;
559 rcutoff = _mm_set1_ps(rcutoff_scalar);
560 rcutoff2 = _mm_mul_ps(rcutoff,rcutoff);
562 sh_vdw_invrcut6 = _mm_set1_ps(fr->ic->sh_invrc6);
563 rvdw = _mm_set1_ps(fr->rvdw);
565 /* Avoid stupid compiler warnings */
566 jnrA = jnrB = jnrC = jnrD = 0;
567 j_coord_offsetA = 0;
568 j_coord_offsetB = 0;
569 j_coord_offsetC = 0;
570 j_coord_offsetD = 0;
572 outeriter = 0;
573 inneriter = 0;
575 for(iidx=0;iidx<4*DIM;iidx++)
577 scratch[iidx] = 0.0;
580 /* Start outer loop over neighborlists */
581 for(iidx=0; iidx<nri; iidx++)
583 /* Load shift vector for this list */
584 i_shift_offset = DIM*shiftidx[iidx];
586 /* Load limits for loop over neighbors */
587 j_index_start = jindex[iidx];
588 j_index_end = jindex[iidx+1];
590 /* Get outer coordinate index */
591 inr = iinr[iidx];
592 i_coord_offset = DIM*inr;
594 /* Load i particle coords and add shift vector */
595 gmx_mm_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
597 fix0 = _mm_setzero_ps();
598 fiy0 = _mm_setzero_ps();
599 fiz0 = _mm_setzero_ps();
601 /* Load parameters for i particles */
602 iq0 = _mm_mul_ps(facel,_mm_load1_ps(charge+inr+0));
603 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
605 /* Start inner kernel loop */
606 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
609 /* Get j neighbor index, and coordinate index */
610 jnrA = jjnr[jidx];
611 jnrB = jjnr[jidx+1];
612 jnrC = jjnr[jidx+2];
613 jnrD = jjnr[jidx+3];
614 j_coord_offsetA = DIM*jnrA;
615 j_coord_offsetB = DIM*jnrB;
616 j_coord_offsetC = DIM*jnrC;
617 j_coord_offsetD = DIM*jnrD;
619 /* load j atom coordinates */
620 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
621 x+j_coord_offsetC,x+j_coord_offsetD,
622 &jx0,&jy0,&jz0);
624 /* Calculate displacement vector */
625 dx00 = _mm_sub_ps(ix0,jx0);
626 dy00 = _mm_sub_ps(iy0,jy0);
627 dz00 = _mm_sub_ps(iz0,jz0);
629 /* Calculate squared distance and things based on it */
630 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
632 rinv00 = gmx_mm_invsqrt_ps(rsq00);
634 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
636 /* Load parameters for j particles */
637 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
638 charge+jnrC+0,charge+jnrD+0);
639 vdwjidx0A = 2*vdwtype[jnrA+0];
640 vdwjidx0B = 2*vdwtype[jnrB+0];
641 vdwjidx0C = 2*vdwtype[jnrC+0];
642 vdwjidx0D = 2*vdwtype[jnrD+0];
644 /**************************
645 * CALCULATE INTERACTIONS *
646 **************************/
648 if (gmx_mm_any_lt(rsq00,rcutoff2))
651 r00 = _mm_mul_ps(rsq00,rinv00);
653 /* Compute parameters for interactions between i and j atoms */
654 qq00 = _mm_mul_ps(iq0,jq0);
655 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
656 vdwparam+vdwioffset0+vdwjidx0B,
657 vdwparam+vdwioffset0+vdwjidx0C,
658 vdwparam+vdwioffset0+vdwjidx0D,
659 &c6_00,&c12_00);
660 c6grid_00 = gmx_mm_load_4real_swizzle_ps(vdwgridparam+vdwioffset0+vdwjidx0A,
661 vdwgridparam+vdwioffset0+vdwjidx0B,
662 vdwgridparam+vdwioffset0+vdwjidx0C,
663 vdwgridparam+vdwioffset0+vdwjidx0D);
665 /* EWALD ELECTROSTATICS */
667 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
668 ewrt = _mm_mul_ps(r00,ewtabscale);
669 ewitab = _mm_cvttps_epi32(ewrt);
670 eweps = _mm_sub_ps(ewrt,_mm_cvtepi32_ps(ewitab));
671 gmx_mm_load_4pair_swizzle_ps(ewtab+gmx_mm_extract_epi32(ewitab,0),ewtab+gmx_mm_extract_epi32(ewitab,1),
672 ewtab+gmx_mm_extract_epi32(ewitab,2),ewtab+gmx_mm_extract_epi32(ewitab,3),
673 &ewtabF,&ewtabFn);
674 felec = _mm_add_ps(_mm_mul_ps( _mm_sub_ps(one,eweps),ewtabF),_mm_mul_ps(eweps,ewtabFn));
675 felec = _mm_mul_ps(_mm_mul_ps(qq00,rinv00),_mm_sub_ps(rinvsq00,felec));
677 /* Analytical LJ-PME */
678 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
679 ewcljrsq = _mm_mul_ps(ewclj2,rsq00);
680 ewclj6 = _mm_mul_ps(ewclj2,_mm_mul_ps(ewclj2,ewclj2));
681 exponent = gmx_simd_exp_r(ewcljrsq);
682 /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
683 poly = _mm_mul_ps(exponent,_mm_add_ps(_mm_sub_ps(one,ewcljrsq),_mm_mul_ps(_mm_mul_ps(ewcljrsq,ewcljrsq),one_half)));
684 /* f6A = 6 * C6grid * (1 - poly) */
685 f6A = _mm_mul_ps(c6grid_00,_mm_sub_ps(one,poly));
686 /* f6B = C6grid * exponent * beta^6 */
687 f6B = _mm_mul_ps(_mm_mul_ps(c6grid_00,one_sixth),_mm_mul_ps(exponent,ewclj6));
688 /* fvdw = 12*C12/r13 - ((6*C6 - f6A)/r6 + f6B)/r */
689 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);
691 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
693 fscal = _mm_add_ps(felec,fvdw);
695 fscal = _mm_and_ps(fscal,cutoff_mask);
697 /* Calculate temporary vectorial force */
698 tx = _mm_mul_ps(fscal,dx00);
699 ty = _mm_mul_ps(fscal,dy00);
700 tz = _mm_mul_ps(fscal,dz00);
702 /* Update vectorial force */
703 fix0 = _mm_add_ps(fix0,tx);
704 fiy0 = _mm_add_ps(fiy0,ty);
705 fiz0 = _mm_add_ps(fiz0,tz);
707 fjptrA = f+j_coord_offsetA;
708 fjptrB = f+j_coord_offsetB;
709 fjptrC = f+j_coord_offsetC;
710 fjptrD = f+j_coord_offsetD;
711 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
715 /* Inner loop uses 62 flops */
718 if(jidx<j_index_end)
721 /* Get j neighbor index, and coordinate index */
722 jnrlistA = jjnr[jidx];
723 jnrlistB = jjnr[jidx+1];
724 jnrlistC = jjnr[jidx+2];
725 jnrlistD = jjnr[jidx+3];
726 /* Sign of each element will be negative for non-real atoms.
727 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
728 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
730 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
731 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
732 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
733 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
734 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
735 j_coord_offsetA = DIM*jnrA;
736 j_coord_offsetB = DIM*jnrB;
737 j_coord_offsetC = DIM*jnrC;
738 j_coord_offsetD = DIM*jnrD;
740 /* load j atom coordinates */
741 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
742 x+j_coord_offsetC,x+j_coord_offsetD,
743 &jx0,&jy0,&jz0);
745 /* Calculate displacement vector */
746 dx00 = _mm_sub_ps(ix0,jx0);
747 dy00 = _mm_sub_ps(iy0,jy0);
748 dz00 = _mm_sub_ps(iz0,jz0);
750 /* Calculate squared distance and things based on it */
751 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
753 rinv00 = gmx_mm_invsqrt_ps(rsq00);
755 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
757 /* Load parameters for j particles */
758 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
759 charge+jnrC+0,charge+jnrD+0);
760 vdwjidx0A = 2*vdwtype[jnrA+0];
761 vdwjidx0B = 2*vdwtype[jnrB+0];
762 vdwjidx0C = 2*vdwtype[jnrC+0];
763 vdwjidx0D = 2*vdwtype[jnrD+0];
765 /**************************
766 * CALCULATE INTERACTIONS *
767 **************************/
769 if (gmx_mm_any_lt(rsq00,rcutoff2))
772 r00 = _mm_mul_ps(rsq00,rinv00);
773 r00 = _mm_andnot_ps(dummy_mask,r00);
775 /* Compute parameters for interactions between i and j atoms */
776 qq00 = _mm_mul_ps(iq0,jq0);
777 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
778 vdwparam+vdwioffset0+vdwjidx0B,
779 vdwparam+vdwioffset0+vdwjidx0C,
780 vdwparam+vdwioffset0+vdwjidx0D,
781 &c6_00,&c12_00);
782 c6grid_00 = gmx_mm_load_4real_swizzle_ps(vdwgridparam+vdwioffset0+vdwjidx0A,
783 vdwgridparam+vdwioffset0+vdwjidx0B,
784 vdwgridparam+vdwioffset0+vdwjidx0C,
785 vdwgridparam+vdwioffset0+vdwjidx0D);
787 /* EWALD ELECTROSTATICS */
789 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
790 ewrt = _mm_mul_ps(r00,ewtabscale);
791 ewitab = _mm_cvttps_epi32(ewrt);
792 eweps = _mm_sub_ps(ewrt,_mm_cvtepi32_ps(ewitab));
793 gmx_mm_load_4pair_swizzle_ps(ewtab+gmx_mm_extract_epi32(ewitab,0),ewtab+gmx_mm_extract_epi32(ewitab,1),
794 ewtab+gmx_mm_extract_epi32(ewitab,2),ewtab+gmx_mm_extract_epi32(ewitab,3),
795 &ewtabF,&ewtabFn);
796 felec = _mm_add_ps(_mm_mul_ps( _mm_sub_ps(one,eweps),ewtabF),_mm_mul_ps(eweps,ewtabFn));
797 felec = _mm_mul_ps(_mm_mul_ps(qq00,rinv00),_mm_sub_ps(rinvsq00,felec));
799 /* Analytical LJ-PME */
800 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
801 ewcljrsq = _mm_mul_ps(ewclj2,rsq00);
802 ewclj6 = _mm_mul_ps(ewclj2,_mm_mul_ps(ewclj2,ewclj2));
803 exponent = gmx_simd_exp_r(ewcljrsq);
804 /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
805 poly = _mm_mul_ps(exponent,_mm_add_ps(_mm_sub_ps(one,ewcljrsq),_mm_mul_ps(_mm_mul_ps(ewcljrsq,ewcljrsq),one_half)));
806 /* f6A = 6 * C6grid * (1 - poly) */
807 f6A = _mm_mul_ps(c6grid_00,_mm_sub_ps(one,poly));
808 /* f6B = C6grid * exponent * beta^6 */
809 f6B = _mm_mul_ps(_mm_mul_ps(c6grid_00,one_sixth),_mm_mul_ps(exponent,ewclj6));
810 /* fvdw = 12*C12/r13 - ((6*C6 - f6A)/r6 + f6B)/r */
811 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);
813 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
815 fscal = _mm_add_ps(felec,fvdw);
817 fscal = _mm_and_ps(fscal,cutoff_mask);
819 fscal = _mm_andnot_ps(dummy_mask,fscal);
821 /* Calculate temporary vectorial force */
822 tx = _mm_mul_ps(fscal,dx00);
823 ty = _mm_mul_ps(fscal,dy00);
824 tz = _mm_mul_ps(fscal,dz00);
826 /* Update vectorial force */
827 fix0 = _mm_add_ps(fix0,tx);
828 fiy0 = _mm_add_ps(fiy0,ty);
829 fiz0 = _mm_add_ps(fiz0,tz);
831 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
832 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
833 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
834 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
835 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
839 /* Inner loop uses 63 flops */
842 /* End of innermost loop */
844 gmx_mm_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
845 f+i_coord_offset,fshift+i_shift_offset);
847 /* Increment number of inner iterations */
848 inneriter += j_index_end - j_index_start;
850 /* Outer loop uses 7 flops */
853 /* Increment number of outer iterations */
854 outeriter += nri;
856 /* Update outer/inner flops */
858 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_F,outeriter*7 + inneriter*63);