2 * Note: this file was generated by the Gromacs sse4_1_single kernel generator.
4 * This source code is part of
8 * Copyright (c) 2001-2012, The GROMACS Development Team
10 * Gromacs is a library for molecular simulation and trajectory analysis,
11 * written by Erik Lindahl, David van der Spoel, Berk Hess, and others - for
12 * a full list of developers and information, check out http://www.gromacs.org
14 * This program is free software; you can redistribute it and/or modify it under
15 * the terms of the GNU Lesser General Public License as published by the Free
16 * Software Foundation; either version 2 of the License, or (at your option) any
19 * To help fund GROMACS development, we humbly ask that you cite
20 * the papers people have written on it - you can find them on the website.
28 #include "../nb_kernel.h"
29 #include "types/simple.h"
33 #include "gmx_math_x86_sse4_1_single.h"
34 #include "kernelutil_x86_sse4_1_single.h"
37 * Gromacs nonbonded kernel: nb_kernel_ElecEwSw_VdwLJSw_GeomP1P1_VF_sse4_1_single
38 * Electrostatics interaction: Ewald
39 * VdW interaction: LennardJones
40 * Geometry: Particle-Particle
41 * Calculate force/pot: PotentialAndForce
44 nb_kernel_ElecEwSw_VdwLJSw_GeomP1P1_VF_sse4_1_single
45 (t_nblist * gmx_restrict nlist,
46 rvec * gmx_restrict xx,
47 rvec * gmx_restrict ff,
48 t_forcerec * gmx_restrict fr,
49 t_mdatoms * gmx_restrict mdatoms,
50 nb_kernel_data_t * gmx_restrict kernel_data,
51 t_nrnb * gmx_restrict nrnb)
53 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
54 * just 0 for non-waters.
55 * Suffixes A,B,C,D refer to j loop unrolling done with SSE, e.g. for the four different
56 * jnr indices corresponding to data put in the four positions in the SIMD register.
58 int i_shift_offset,i_coord_offset,outeriter,inneriter;
59 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
60 int jnrA,jnrB,jnrC,jnrD;
61 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
62 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
63 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
65 real *shiftvec,*fshift,*x,*f;
66 real *fjptrA,*fjptrB,*fjptrC,*fjptrD;
68 __m128 tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
70 __m128 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
71 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
72 __m128 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
73 __m128 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
74 __m128 velec,felec,velecsum,facel,crf,krf,krf2;
77 __m128 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
80 __m128 one_sixth = _mm_set1_ps(1.0/6.0);
81 __m128 one_twelfth = _mm_set1_ps(1.0/12.0);
83 __m128 ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
85 __m128 rswitch,swV3,swV4,swV5,swF2,swF3,swF4,d,d2,sw,dsw;
86 real rswitch_scalar,d_scalar;
87 __m128 dummy_mask,cutoff_mask;
88 __m128 signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
89 __m128 one = _mm_set1_ps(1.0);
90 __m128 two = _mm_set1_ps(2.0);
96 jindex = nlist->jindex;
98 shiftidx = nlist->shift;
100 shiftvec = fr->shift_vec[0];
101 fshift = fr->fshift[0];
102 facel = _mm_set1_ps(fr->epsfac);
103 charge = mdatoms->chargeA;
104 nvdwtype = fr->ntype;
106 vdwtype = mdatoms->typeA;
108 sh_ewald = _mm_set1_ps(fr->ic->sh_ewald);
109 ewtab = fr->ic->tabq_coul_FDV0;
110 ewtabscale = _mm_set1_ps(fr->ic->tabq_scale);
111 ewtabhalfspace = _mm_set1_ps(0.5/fr->ic->tabq_scale);
113 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
114 rcutoff_scalar = fr->rcoulomb;
115 rcutoff = _mm_set1_ps(rcutoff_scalar);
116 rcutoff2 = _mm_mul_ps(rcutoff,rcutoff);
118 rswitch_scalar = fr->rcoulomb_switch;
119 rswitch = _mm_set1_ps(rswitch_scalar);
120 /* Setup switch parameters */
121 d_scalar = rcutoff_scalar-rswitch_scalar;
122 d = _mm_set1_ps(d_scalar);
123 swV3 = _mm_set1_ps(-10.0/(d_scalar*d_scalar*d_scalar));
124 swV4 = _mm_set1_ps( 15.0/(d_scalar*d_scalar*d_scalar*d_scalar));
125 swV5 = _mm_set1_ps( -6.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
126 swF2 = _mm_set1_ps(-30.0/(d_scalar*d_scalar*d_scalar));
127 swF3 = _mm_set1_ps( 60.0/(d_scalar*d_scalar*d_scalar*d_scalar));
128 swF4 = _mm_set1_ps(-30.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
130 /* Avoid stupid compiler warnings */
131 jnrA = jnrB = jnrC = jnrD = 0;
140 for(iidx=0;iidx<4*DIM;iidx++)
145 /* Start outer loop over neighborlists */
146 for(iidx=0; iidx<nri; iidx++)
148 /* Load shift vector for this list */
149 i_shift_offset = DIM*shiftidx[iidx];
151 /* Load limits for loop over neighbors */
152 j_index_start = jindex[iidx];
153 j_index_end = jindex[iidx+1];
155 /* Get outer coordinate index */
157 i_coord_offset = DIM*inr;
159 /* Load i particle coords and add shift vector */
160 gmx_mm_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
162 fix0 = _mm_setzero_ps();
163 fiy0 = _mm_setzero_ps();
164 fiz0 = _mm_setzero_ps();
166 /* Load parameters for i particles */
167 iq0 = _mm_mul_ps(facel,_mm_load1_ps(charge+inr+0));
168 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
170 /* Reset potential sums */
171 velecsum = _mm_setzero_ps();
172 vvdwsum = _mm_setzero_ps();
174 /* Start inner kernel loop */
175 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
178 /* Get j neighbor index, and coordinate index */
183 j_coord_offsetA = DIM*jnrA;
184 j_coord_offsetB = DIM*jnrB;
185 j_coord_offsetC = DIM*jnrC;
186 j_coord_offsetD = DIM*jnrD;
188 /* load j atom coordinates */
189 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
190 x+j_coord_offsetC,x+j_coord_offsetD,
193 /* Calculate displacement vector */
194 dx00 = _mm_sub_ps(ix0,jx0);
195 dy00 = _mm_sub_ps(iy0,jy0);
196 dz00 = _mm_sub_ps(iz0,jz0);
198 /* Calculate squared distance and things based on it */
199 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
201 rinv00 = gmx_mm_invsqrt_ps(rsq00);
203 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
205 /* Load parameters for j particles */
206 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
207 charge+jnrC+0,charge+jnrD+0);
208 vdwjidx0A = 2*vdwtype[jnrA+0];
209 vdwjidx0B = 2*vdwtype[jnrB+0];
210 vdwjidx0C = 2*vdwtype[jnrC+0];
211 vdwjidx0D = 2*vdwtype[jnrD+0];
213 /**************************
214 * CALCULATE INTERACTIONS *
215 **************************/
217 if (gmx_mm_any_lt(rsq00,rcutoff2))
220 r00 = _mm_mul_ps(rsq00,rinv00);
222 /* Compute parameters for interactions between i and j atoms */
223 qq00 = _mm_mul_ps(iq0,jq0);
224 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
225 vdwparam+vdwioffset0+vdwjidx0B,
226 vdwparam+vdwioffset0+vdwjidx0C,
227 vdwparam+vdwioffset0+vdwjidx0D,
230 /* EWALD ELECTROSTATICS */
232 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
233 ewrt = _mm_mul_ps(r00,ewtabscale);
234 ewitab = _mm_cvttps_epi32(ewrt);
235 eweps = _mm_sub_ps(ewrt,_mm_round_ps(ewrt, _MM_FROUND_FLOOR));
236 ewitab = _mm_slli_epi32(ewitab,2);
237 ewtabF = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,0) );
238 ewtabD = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,1) );
239 ewtabV = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,2) );
240 ewtabFn = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,3) );
241 _MM_TRANSPOSE4_PS(ewtabF,ewtabD,ewtabV,ewtabFn);
242 felec = _mm_add_ps(ewtabF,_mm_mul_ps(eweps,ewtabD));
243 velec = _mm_sub_ps(ewtabV,_mm_mul_ps(_mm_mul_ps(ewtabhalfspace,eweps),_mm_add_ps(ewtabF,felec)));
244 velec = _mm_mul_ps(qq00,_mm_sub_ps(rinv00,velec));
245 felec = _mm_mul_ps(_mm_mul_ps(qq00,rinv00),_mm_sub_ps(rinvsq00,felec));
247 /* LENNARD-JONES DISPERSION/REPULSION */
249 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
250 vvdw6 = _mm_mul_ps(c6_00,rinvsix);
251 vvdw12 = _mm_mul_ps(c12_00,_mm_mul_ps(rinvsix,rinvsix));
252 vvdw = _mm_sub_ps( _mm_mul_ps(vvdw12,one_twelfth) , _mm_mul_ps(vvdw6,one_sixth) );
253 fvdw = _mm_mul_ps(_mm_sub_ps(vvdw12,vvdw6),rinvsq00);
255 d = _mm_sub_ps(r00,rswitch);
256 d = _mm_max_ps(d,_mm_setzero_ps());
257 d2 = _mm_mul_ps(d,d);
258 sw = _mm_add_ps(one,_mm_mul_ps(d2,_mm_mul_ps(d,_mm_add_ps(swV3,_mm_mul_ps(d,_mm_add_ps(swV4,_mm_mul_ps(d,swV5)))))));
260 dsw = _mm_mul_ps(d2,_mm_add_ps(swF2,_mm_mul_ps(d,_mm_add_ps(swF3,_mm_mul_ps(d,swF4)))));
262 /* Evaluate switch function */
263 /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
264 felec = _mm_sub_ps( _mm_mul_ps(felec,sw) , _mm_mul_ps(rinv00,_mm_mul_ps(velec,dsw)) );
265 fvdw = _mm_sub_ps( _mm_mul_ps(fvdw,sw) , _mm_mul_ps(rinv00,_mm_mul_ps(vvdw,dsw)) );
266 velec = _mm_mul_ps(velec,sw);
267 vvdw = _mm_mul_ps(vvdw,sw);
268 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
270 /* Update potential sum for this i atom from the interaction with this j atom. */
271 velec = _mm_and_ps(velec,cutoff_mask);
272 velecsum = _mm_add_ps(velecsum,velec);
273 vvdw = _mm_and_ps(vvdw,cutoff_mask);
274 vvdwsum = _mm_add_ps(vvdwsum,vvdw);
276 fscal = _mm_add_ps(felec,fvdw);
278 fscal = _mm_and_ps(fscal,cutoff_mask);
280 /* Calculate temporary vectorial force */
281 tx = _mm_mul_ps(fscal,dx00);
282 ty = _mm_mul_ps(fscal,dy00);
283 tz = _mm_mul_ps(fscal,dz00);
285 /* Update vectorial force */
286 fix0 = _mm_add_ps(fix0,tx);
287 fiy0 = _mm_add_ps(fiy0,ty);
288 fiz0 = _mm_add_ps(fiz0,tz);
290 fjptrA = f+j_coord_offsetA;
291 fjptrB = f+j_coord_offsetB;
292 fjptrC = f+j_coord_offsetC;
293 fjptrD = f+j_coord_offsetD;
294 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
298 /* Inner loop uses 83 flops */
304 /* Get j neighbor index, and coordinate index */
305 jnrlistA = jjnr[jidx];
306 jnrlistB = jjnr[jidx+1];
307 jnrlistC = jjnr[jidx+2];
308 jnrlistD = jjnr[jidx+3];
309 /* Sign of each element will be negative for non-real atoms.
310 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
311 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
313 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
314 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
315 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
316 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
317 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
318 j_coord_offsetA = DIM*jnrA;
319 j_coord_offsetB = DIM*jnrB;
320 j_coord_offsetC = DIM*jnrC;
321 j_coord_offsetD = DIM*jnrD;
323 /* load j atom coordinates */
324 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
325 x+j_coord_offsetC,x+j_coord_offsetD,
328 /* Calculate displacement vector */
329 dx00 = _mm_sub_ps(ix0,jx0);
330 dy00 = _mm_sub_ps(iy0,jy0);
331 dz00 = _mm_sub_ps(iz0,jz0);
333 /* Calculate squared distance and things based on it */
334 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
336 rinv00 = gmx_mm_invsqrt_ps(rsq00);
338 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
340 /* Load parameters for j particles */
341 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
342 charge+jnrC+0,charge+jnrD+0);
343 vdwjidx0A = 2*vdwtype[jnrA+0];
344 vdwjidx0B = 2*vdwtype[jnrB+0];
345 vdwjidx0C = 2*vdwtype[jnrC+0];
346 vdwjidx0D = 2*vdwtype[jnrD+0];
348 /**************************
349 * CALCULATE INTERACTIONS *
350 **************************/
352 if (gmx_mm_any_lt(rsq00,rcutoff2))
355 r00 = _mm_mul_ps(rsq00,rinv00);
356 r00 = _mm_andnot_ps(dummy_mask,r00);
358 /* Compute parameters for interactions between i and j atoms */
359 qq00 = _mm_mul_ps(iq0,jq0);
360 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
361 vdwparam+vdwioffset0+vdwjidx0B,
362 vdwparam+vdwioffset0+vdwjidx0C,
363 vdwparam+vdwioffset0+vdwjidx0D,
366 /* EWALD ELECTROSTATICS */
368 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
369 ewrt = _mm_mul_ps(r00,ewtabscale);
370 ewitab = _mm_cvttps_epi32(ewrt);
371 eweps = _mm_sub_ps(ewrt,_mm_round_ps(ewrt, _MM_FROUND_FLOOR));
372 ewitab = _mm_slli_epi32(ewitab,2);
373 ewtabF = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,0) );
374 ewtabD = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,1) );
375 ewtabV = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,2) );
376 ewtabFn = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,3) );
377 _MM_TRANSPOSE4_PS(ewtabF,ewtabD,ewtabV,ewtabFn);
378 felec = _mm_add_ps(ewtabF,_mm_mul_ps(eweps,ewtabD));
379 velec = _mm_sub_ps(ewtabV,_mm_mul_ps(_mm_mul_ps(ewtabhalfspace,eweps),_mm_add_ps(ewtabF,felec)));
380 velec = _mm_mul_ps(qq00,_mm_sub_ps(rinv00,velec));
381 felec = _mm_mul_ps(_mm_mul_ps(qq00,rinv00),_mm_sub_ps(rinvsq00,felec));
383 /* LENNARD-JONES DISPERSION/REPULSION */
385 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
386 vvdw6 = _mm_mul_ps(c6_00,rinvsix);
387 vvdw12 = _mm_mul_ps(c12_00,_mm_mul_ps(rinvsix,rinvsix));
388 vvdw = _mm_sub_ps( _mm_mul_ps(vvdw12,one_twelfth) , _mm_mul_ps(vvdw6,one_sixth) );
389 fvdw = _mm_mul_ps(_mm_sub_ps(vvdw12,vvdw6),rinvsq00);
391 d = _mm_sub_ps(r00,rswitch);
392 d = _mm_max_ps(d,_mm_setzero_ps());
393 d2 = _mm_mul_ps(d,d);
394 sw = _mm_add_ps(one,_mm_mul_ps(d2,_mm_mul_ps(d,_mm_add_ps(swV3,_mm_mul_ps(d,_mm_add_ps(swV4,_mm_mul_ps(d,swV5)))))));
396 dsw = _mm_mul_ps(d2,_mm_add_ps(swF2,_mm_mul_ps(d,_mm_add_ps(swF3,_mm_mul_ps(d,swF4)))));
398 /* Evaluate switch function */
399 /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
400 felec = _mm_sub_ps( _mm_mul_ps(felec,sw) , _mm_mul_ps(rinv00,_mm_mul_ps(velec,dsw)) );
401 fvdw = _mm_sub_ps( _mm_mul_ps(fvdw,sw) , _mm_mul_ps(rinv00,_mm_mul_ps(vvdw,dsw)) );
402 velec = _mm_mul_ps(velec,sw);
403 vvdw = _mm_mul_ps(vvdw,sw);
404 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
406 /* Update potential sum for this i atom from the interaction with this j atom. */
407 velec = _mm_and_ps(velec,cutoff_mask);
408 velec = _mm_andnot_ps(dummy_mask,velec);
409 velecsum = _mm_add_ps(velecsum,velec);
410 vvdw = _mm_and_ps(vvdw,cutoff_mask);
411 vvdw = _mm_andnot_ps(dummy_mask,vvdw);
412 vvdwsum = _mm_add_ps(vvdwsum,vvdw);
414 fscal = _mm_add_ps(felec,fvdw);
416 fscal = _mm_and_ps(fscal,cutoff_mask);
418 fscal = _mm_andnot_ps(dummy_mask,fscal);
420 /* Calculate temporary vectorial force */
421 tx = _mm_mul_ps(fscal,dx00);
422 ty = _mm_mul_ps(fscal,dy00);
423 tz = _mm_mul_ps(fscal,dz00);
425 /* Update vectorial force */
426 fix0 = _mm_add_ps(fix0,tx);
427 fiy0 = _mm_add_ps(fiy0,ty);
428 fiz0 = _mm_add_ps(fiz0,tz);
430 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
431 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
432 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
433 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
434 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
438 /* Inner loop uses 84 flops */
441 /* End of innermost loop */
443 gmx_mm_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
444 f+i_coord_offset,fshift+i_shift_offset);
447 /* Update potential energies */
448 gmx_mm_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
449 gmx_mm_update_1pot_ps(vvdwsum,kernel_data->energygrp_vdw+ggid);
451 /* Increment number of inner iterations */
452 inneriter += j_index_end - j_index_start;
454 /* Outer loop uses 9 flops */
457 /* Increment number of outer iterations */
460 /* Update outer/inner flops */
462 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_VF,outeriter*9 + inneriter*84);
465 * Gromacs nonbonded kernel: nb_kernel_ElecEwSw_VdwLJSw_GeomP1P1_F_sse4_1_single
466 * Electrostatics interaction: Ewald
467 * VdW interaction: LennardJones
468 * Geometry: Particle-Particle
469 * Calculate force/pot: Force
472 nb_kernel_ElecEwSw_VdwLJSw_GeomP1P1_F_sse4_1_single
473 (t_nblist * gmx_restrict nlist,
474 rvec * gmx_restrict xx,
475 rvec * gmx_restrict ff,
476 t_forcerec * gmx_restrict fr,
477 t_mdatoms * gmx_restrict mdatoms,
478 nb_kernel_data_t * gmx_restrict kernel_data,
479 t_nrnb * gmx_restrict nrnb)
481 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
482 * just 0 for non-waters.
483 * Suffixes A,B,C,D refer to j loop unrolling done with SSE, e.g. for the four different
484 * jnr indices corresponding to data put in the four positions in the SIMD register.
486 int i_shift_offset,i_coord_offset,outeriter,inneriter;
487 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
488 int jnrA,jnrB,jnrC,jnrD;
489 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
490 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
491 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
493 real *shiftvec,*fshift,*x,*f;
494 real *fjptrA,*fjptrB,*fjptrC,*fjptrD;
496 __m128 tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
498 __m128 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
499 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
500 __m128 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
501 __m128 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
502 __m128 velec,felec,velecsum,facel,crf,krf,krf2;
505 __m128 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
508 __m128 one_sixth = _mm_set1_ps(1.0/6.0);
509 __m128 one_twelfth = _mm_set1_ps(1.0/12.0);
511 __m128 ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
513 __m128 rswitch,swV3,swV4,swV5,swF2,swF3,swF4,d,d2,sw,dsw;
514 real rswitch_scalar,d_scalar;
515 __m128 dummy_mask,cutoff_mask;
516 __m128 signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
517 __m128 one = _mm_set1_ps(1.0);
518 __m128 two = _mm_set1_ps(2.0);
524 jindex = nlist->jindex;
526 shiftidx = nlist->shift;
528 shiftvec = fr->shift_vec[0];
529 fshift = fr->fshift[0];
530 facel = _mm_set1_ps(fr->epsfac);
531 charge = mdatoms->chargeA;
532 nvdwtype = fr->ntype;
534 vdwtype = mdatoms->typeA;
536 sh_ewald = _mm_set1_ps(fr->ic->sh_ewald);
537 ewtab = fr->ic->tabq_coul_FDV0;
538 ewtabscale = _mm_set1_ps(fr->ic->tabq_scale);
539 ewtabhalfspace = _mm_set1_ps(0.5/fr->ic->tabq_scale);
541 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
542 rcutoff_scalar = fr->rcoulomb;
543 rcutoff = _mm_set1_ps(rcutoff_scalar);
544 rcutoff2 = _mm_mul_ps(rcutoff,rcutoff);
546 rswitch_scalar = fr->rcoulomb_switch;
547 rswitch = _mm_set1_ps(rswitch_scalar);
548 /* Setup switch parameters */
549 d_scalar = rcutoff_scalar-rswitch_scalar;
550 d = _mm_set1_ps(d_scalar);
551 swV3 = _mm_set1_ps(-10.0/(d_scalar*d_scalar*d_scalar));
552 swV4 = _mm_set1_ps( 15.0/(d_scalar*d_scalar*d_scalar*d_scalar));
553 swV5 = _mm_set1_ps( -6.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
554 swF2 = _mm_set1_ps(-30.0/(d_scalar*d_scalar*d_scalar));
555 swF3 = _mm_set1_ps( 60.0/(d_scalar*d_scalar*d_scalar*d_scalar));
556 swF4 = _mm_set1_ps(-30.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
558 /* Avoid stupid compiler warnings */
559 jnrA = jnrB = jnrC = jnrD = 0;
568 for(iidx=0;iidx<4*DIM;iidx++)
573 /* Start outer loop over neighborlists */
574 for(iidx=0; iidx<nri; iidx++)
576 /* Load shift vector for this list */
577 i_shift_offset = DIM*shiftidx[iidx];
579 /* Load limits for loop over neighbors */
580 j_index_start = jindex[iidx];
581 j_index_end = jindex[iidx+1];
583 /* Get outer coordinate index */
585 i_coord_offset = DIM*inr;
587 /* Load i particle coords and add shift vector */
588 gmx_mm_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
590 fix0 = _mm_setzero_ps();
591 fiy0 = _mm_setzero_ps();
592 fiz0 = _mm_setzero_ps();
594 /* Load parameters for i particles */
595 iq0 = _mm_mul_ps(facel,_mm_load1_ps(charge+inr+0));
596 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
598 /* Start inner kernel loop */
599 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
602 /* Get j neighbor index, and coordinate index */
607 j_coord_offsetA = DIM*jnrA;
608 j_coord_offsetB = DIM*jnrB;
609 j_coord_offsetC = DIM*jnrC;
610 j_coord_offsetD = DIM*jnrD;
612 /* load j atom coordinates */
613 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
614 x+j_coord_offsetC,x+j_coord_offsetD,
617 /* Calculate displacement vector */
618 dx00 = _mm_sub_ps(ix0,jx0);
619 dy00 = _mm_sub_ps(iy0,jy0);
620 dz00 = _mm_sub_ps(iz0,jz0);
622 /* Calculate squared distance and things based on it */
623 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
625 rinv00 = gmx_mm_invsqrt_ps(rsq00);
627 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
629 /* Load parameters for j particles */
630 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
631 charge+jnrC+0,charge+jnrD+0);
632 vdwjidx0A = 2*vdwtype[jnrA+0];
633 vdwjidx0B = 2*vdwtype[jnrB+0];
634 vdwjidx0C = 2*vdwtype[jnrC+0];
635 vdwjidx0D = 2*vdwtype[jnrD+0];
637 /**************************
638 * CALCULATE INTERACTIONS *
639 **************************/
641 if (gmx_mm_any_lt(rsq00,rcutoff2))
644 r00 = _mm_mul_ps(rsq00,rinv00);
646 /* Compute parameters for interactions between i and j atoms */
647 qq00 = _mm_mul_ps(iq0,jq0);
648 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
649 vdwparam+vdwioffset0+vdwjidx0B,
650 vdwparam+vdwioffset0+vdwjidx0C,
651 vdwparam+vdwioffset0+vdwjidx0D,
654 /* EWALD ELECTROSTATICS */
656 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
657 ewrt = _mm_mul_ps(r00,ewtabscale);
658 ewitab = _mm_cvttps_epi32(ewrt);
659 eweps = _mm_sub_ps(ewrt,_mm_round_ps(ewrt, _MM_FROUND_FLOOR));
660 ewitab = _mm_slli_epi32(ewitab,2);
661 ewtabF = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,0) );
662 ewtabD = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,1) );
663 ewtabV = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,2) );
664 ewtabFn = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,3) );
665 _MM_TRANSPOSE4_PS(ewtabF,ewtabD,ewtabV,ewtabFn);
666 felec = _mm_add_ps(ewtabF,_mm_mul_ps(eweps,ewtabD));
667 velec = _mm_sub_ps(ewtabV,_mm_mul_ps(_mm_mul_ps(ewtabhalfspace,eweps),_mm_add_ps(ewtabF,felec)));
668 velec = _mm_mul_ps(qq00,_mm_sub_ps(rinv00,velec));
669 felec = _mm_mul_ps(_mm_mul_ps(qq00,rinv00),_mm_sub_ps(rinvsq00,felec));
671 /* LENNARD-JONES DISPERSION/REPULSION */
673 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
674 vvdw6 = _mm_mul_ps(c6_00,rinvsix);
675 vvdw12 = _mm_mul_ps(c12_00,_mm_mul_ps(rinvsix,rinvsix));
676 vvdw = _mm_sub_ps( _mm_mul_ps(vvdw12,one_twelfth) , _mm_mul_ps(vvdw6,one_sixth) );
677 fvdw = _mm_mul_ps(_mm_sub_ps(vvdw12,vvdw6),rinvsq00);
679 d = _mm_sub_ps(r00,rswitch);
680 d = _mm_max_ps(d,_mm_setzero_ps());
681 d2 = _mm_mul_ps(d,d);
682 sw = _mm_add_ps(one,_mm_mul_ps(d2,_mm_mul_ps(d,_mm_add_ps(swV3,_mm_mul_ps(d,_mm_add_ps(swV4,_mm_mul_ps(d,swV5)))))));
684 dsw = _mm_mul_ps(d2,_mm_add_ps(swF2,_mm_mul_ps(d,_mm_add_ps(swF3,_mm_mul_ps(d,swF4)))));
686 /* Evaluate switch function */
687 /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
688 felec = _mm_sub_ps( _mm_mul_ps(felec,sw) , _mm_mul_ps(rinv00,_mm_mul_ps(velec,dsw)) );
689 fvdw = _mm_sub_ps( _mm_mul_ps(fvdw,sw) , _mm_mul_ps(rinv00,_mm_mul_ps(vvdw,dsw)) );
690 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
692 fscal = _mm_add_ps(felec,fvdw);
694 fscal = _mm_and_ps(fscal,cutoff_mask);
696 /* Calculate temporary vectorial force */
697 tx = _mm_mul_ps(fscal,dx00);
698 ty = _mm_mul_ps(fscal,dy00);
699 tz = _mm_mul_ps(fscal,dz00);
701 /* Update vectorial force */
702 fix0 = _mm_add_ps(fix0,tx);
703 fiy0 = _mm_add_ps(fiy0,ty);
704 fiz0 = _mm_add_ps(fiz0,tz);
706 fjptrA = f+j_coord_offsetA;
707 fjptrB = f+j_coord_offsetB;
708 fjptrC = f+j_coord_offsetC;
709 fjptrD = f+j_coord_offsetD;
710 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
714 /* Inner loop uses 77 flops */
720 /* Get j neighbor index, and coordinate index */
721 jnrlistA = jjnr[jidx];
722 jnrlistB = jjnr[jidx+1];
723 jnrlistC = jjnr[jidx+2];
724 jnrlistD = jjnr[jidx+3];
725 /* Sign of each element will be negative for non-real atoms.
726 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
727 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
729 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
730 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
731 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
732 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
733 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
734 j_coord_offsetA = DIM*jnrA;
735 j_coord_offsetB = DIM*jnrB;
736 j_coord_offsetC = DIM*jnrC;
737 j_coord_offsetD = DIM*jnrD;
739 /* load j atom coordinates */
740 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
741 x+j_coord_offsetC,x+j_coord_offsetD,
744 /* Calculate displacement vector */
745 dx00 = _mm_sub_ps(ix0,jx0);
746 dy00 = _mm_sub_ps(iy0,jy0);
747 dz00 = _mm_sub_ps(iz0,jz0);
749 /* Calculate squared distance and things based on it */
750 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
752 rinv00 = gmx_mm_invsqrt_ps(rsq00);
754 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
756 /* Load parameters for j particles */
757 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
758 charge+jnrC+0,charge+jnrD+0);
759 vdwjidx0A = 2*vdwtype[jnrA+0];
760 vdwjidx0B = 2*vdwtype[jnrB+0];
761 vdwjidx0C = 2*vdwtype[jnrC+0];
762 vdwjidx0D = 2*vdwtype[jnrD+0];
764 /**************************
765 * CALCULATE INTERACTIONS *
766 **************************/
768 if (gmx_mm_any_lt(rsq00,rcutoff2))
771 r00 = _mm_mul_ps(rsq00,rinv00);
772 r00 = _mm_andnot_ps(dummy_mask,r00);
774 /* Compute parameters for interactions between i and j atoms */
775 qq00 = _mm_mul_ps(iq0,jq0);
776 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
777 vdwparam+vdwioffset0+vdwjidx0B,
778 vdwparam+vdwioffset0+vdwjidx0C,
779 vdwparam+vdwioffset0+vdwjidx0D,
782 /* EWALD ELECTROSTATICS */
784 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
785 ewrt = _mm_mul_ps(r00,ewtabscale);
786 ewitab = _mm_cvttps_epi32(ewrt);
787 eweps = _mm_sub_ps(ewrt,_mm_round_ps(ewrt, _MM_FROUND_FLOOR));
788 ewitab = _mm_slli_epi32(ewitab,2);
789 ewtabF = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,0) );
790 ewtabD = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,1) );
791 ewtabV = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,2) );
792 ewtabFn = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,3) );
793 _MM_TRANSPOSE4_PS(ewtabF,ewtabD,ewtabV,ewtabFn);
794 felec = _mm_add_ps(ewtabF,_mm_mul_ps(eweps,ewtabD));
795 velec = _mm_sub_ps(ewtabV,_mm_mul_ps(_mm_mul_ps(ewtabhalfspace,eweps),_mm_add_ps(ewtabF,felec)));
796 velec = _mm_mul_ps(qq00,_mm_sub_ps(rinv00,velec));
797 felec = _mm_mul_ps(_mm_mul_ps(qq00,rinv00),_mm_sub_ps(rinvsq00,felec));
799 /* LENNARD-JONES DISPERSION/REPULSION */
801 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
802 vvdw6 = _mm_mul_ps(c6_00,rinvsix);
803 vvdw12 = _mm_mul_ps(c12_00,_mm_mul_ps(rinvsix,rinvsix));
804 vvdw = _mm_sub_ps( _mm_mul_ps(vvdw12,one_twelfth) , _mm_mul_ps(vvdw6,one_sixth) );
805 fvdw = _mm_mul_ps(_mm_sub_ps(vvdw12,vvdw6),rinvsq00);
807 d = _mm_sub_ps(r00,rswitch);
808 d = _mm_max_ps(d,_mm_setzero_ps());
809 d2 = _mm_mul_ps(d,d);
810 sw = _mm_add_ps(one,_mm_mul_ps(d2,_mm_mul_ps(d,_mm_add_ps(swV3,_mm_mul_ps(d,_mm_add_ps(swV4,_mm_mul_ps(d,swV5)))))));
812 dsw = _mm_mul_ps(d2,_mm_add_ps(swF2,_mm_mul_ps(d,_mm_add_ps(swF3,_mm_mul_ps(d,swF4)))));
814 /* Evaluate switch function */
815 /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
816 felec = _mm_sub_ps( _mm_mul_ps(felec,sw) , _mm_mul_ps(rinv00,_mm_mul_ps(velec,dsw)) );
817 fvdw = _mm_sub_ps( _mm_mul_ps(fvdw,sw) , _mm_mul_ps(rinv00,_mm_mul_ps(vvdw,dsw)) );
818 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
820 fscal = _mm_add_ps(felec,fvdw);
822 fscal = _mm_and_ps(fscal,cutoff_mask);
824 fscal = _mm_andnot_ps(dummy_mask,fscal);
826 /* Calculate temporary vectorial force */
827 tx = _mm_mul_ps(fscal,dx00);
828 ty = _mm_mul_ps(fscal,dy00);
829 tz = _mm_mul_ps(fscal,dz00);
831 /* Update vectorial force */
832 fix0 = _mm_add_ps(fix0,tx);
833 fiy0 = _mm_add_ps(fiy0,ty);
834 fiz0 = _mm_add_ps(fiz0,tz);
836 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
837 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
838 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
839 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
840 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
844 /* Inner loop uses 78 flops */
847 /* End of innermost loop */
849 gmx_mm_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
850 f+i_coord_offset,fshift+i_shift_offset);
852 /* Increment number of inner iterations */
853 inneriter += j_index_end - j_index_start;
855 /* Outer loop uses 7 flops */
858 /* Increment number of outer iterations */
861 /* Update outer/inner flops */
863 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_F,outeriter*7 + inneriter*78);