2 * Note: this file was generated by the Gromacs sse2_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_sse2_single.h"
34 #include "kernelutil_x86_sse2_single.h"
37 * Gromacs nonbonded kernel: nb_kernel_ElecRFCut_VdwLJSw_GeomP1P1_VF_sse2_single
38 * Electrostatics interaction: ReactionField
39 * VdW interaction: LennardJones
40 * Geometry: Particle-Particle
41 * Calculate force/pot: PotentialAndForce
44 nb_kernel_ElecRFCut_VdwLJSw_GeomP1P1_VF_sse2_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);
82 __m128 rswitch,swV3,swV4,swV5,swF2,swF3,swF4,d,d2,sw,dsw;
83 real rswitch_scalar,d_scalar;
84 __m128 dummy_mask,cutoff_mask;
85 __m128 signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
86 __m128 one = _mm_set1_ps(1.0);
87 __m128 two = _mm_set1_ps(2.0);
93 jindex = nlist->jindex;
95 shiftidx = nlist->shift;
97 shiftvec = fr->shift_vec[0];
98 fshift = fr->fshift[0];
99 facel = _mm_set1_ps(fr->epsfac);
100 charge = mdatoms->chargeA;
101 krf = _mm_set1_ps(fr->ic->k_rf);
102 krf2 = _mm_set1_ps(fr->ic->k_rf*2.0);
103 crf = _mm_set1_ps(fr->ic->c_rf);
104 nvdwtype = fr->ntype;
106 vdwtype = mdatoms->typeA;
108 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
109 rcutoff_scalar = fr->rcoulomb;
110 rcutoff = _mm_set1_ps(rcutoff_scalar);
111 rcutoff2 = _mm_mul_ps(rcutoff,rcutoff);
113 rswitch_scalar = fr->rvdw_switch;
114 rswitch = _mm_set1_ps(rswitch_scalar);
115 /* Setup switch parameters */
116 d_scalar = rcutoff_scalar-rswitch_scalar;
117 d = _mm_set1_ps(d_scalar);
118 swV3 = _mm_set1_ps(-10.0/(d_scalar*d_scalar*d_scalar));
119 swV4 = _mm_set1_ps( 15.0/(d_scalar*d_scalar*d_scalar*d_scalar));
120 swV5 = _mm_set1_ps( -6.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
121 swF2 = _mm_set1_ps(-30.0/(d_scalar*d_scalar*d_scalar));
122 swF3 = _mm_set1_ps( 60.0/(d_scalar*d_scalar*d_scalar*d_scalar));
123 swF4 = _mm_set1_ps(-30.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
125 /* Avoid stupid compiler warnings */
126 jnrA = jnrB = jnrC = jnrD = 0;
135 for(iidx=0;iidx<4*DIM;iidx++)
140 /* Start outer loop over neighborlists */
141 for(iidx=0; iidx<nri; iidx++)
143 /* Load shift vector for this list */
144 i_shift_offset = DIM*shiftidx[iidx];
146 /* Load limits for loop over neighbors */
147 j_index_start = jindex[iidx];
148 j_index_end = jindex[iidx+1];
150 /* Get outer coordinate index */
152 i_coord_offset = DIM*inr;
154 /* Load i particle coords and add shift vector */
155 gmx_mm_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
157 fix0 = _mm_setzero_ps();
158 fiy0 = _mm_setzero_ps();
159 fiz0 = _mm_setzero_ps();
161 /* Load parameters for i particles */
162 iq0 = _mm_mul_ps(facel,_mm_load1_ps(charge+inr+0));
163 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
165 /* Reset potential sums */
166 velecsum = _mm_setzero_ps();
167 vvdwsum = _mm_setzero_ps();
169 /* Start inner kernel loop */
170 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
173 /* Get j neighbor index, and coordinate index */
178 j_coord_offsetA = DIM*jnrA;
179 j_coord_offsetB = DIM*jnrB;
180 j_coord_offsetC = DIM*jnrC;
181 j_coord_offsetD = DIM*jnrD;
183 /* load j atom coordinates */
184 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
185 x+j_coord_offsetC,x+j_coord_offsetD,
188 /* Calculate displacement vector */
189 dx00 = _mm_sub_ps(ix0,jx0);
190 dy00 = _mm_sub_ps(iy0,jy0);
191 dz00 = _mm_sub_ps(iz0,jz0);
193 /* Calculate squared distance and things based on it */
194 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
196 rinv00 = gmx_mm_invsqrt_ps(rsq00);
198 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
200 /* Load parameters for j particles */
201 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
202 charge+jnrC+0,charge+jnrD+0);
203 vdwjidx0A = 2*vdwtype[jnrA+0];
204 vdwjidx0B = 2*vdwtype[jnrB+0];
205 vdwjidx0C = 2*vdwtype[jnrC+0];
206 vdwjidx0D = 2*vdwtype[jnrD+0];
208 /**************************
209 * CALCULATE INTERACTIONS *
210 **************************/
212 if (gmx_mm_any_lt(rsq00,rcutoff2))
215 r00 = _mm_mul_ps(rsq00,rinv00);
217 /* Compute parameters for interactions between i and j atoms */
218 qq00 = _mm_mul_ps(iq0,jq0);
219 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
220 vdwparam+vdwioffset0+vdwjidx0B,
221 vdwparam+vdwioffset0+vdwjidx0C,
222 vdwparam+vdwioffset0+vdwjidx0D,
225 /* REACTION-FIELD ELECTROSTATICS */
226 velec = _mm_mul_ps(qq00,_mm_sub_ps(_mm_add_ps(rinv00,_mm_mul_ps(krf,rsq00)),crf));
227 felec = _mm_mul_ps(qq00,_mm_sub_ps(_mm_mul_ps(rinv00,rinvsq00),krf2));
229 /* LENNARD-JONES DISPERSION/REPULSION */
231 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
232 vvdw6 = _mm_mul_ps(c6_00,rinvsix);
233 vvdw12 = _mm_mul_ps(c12_00,_mm_mul_ps(rinvsix,rinvsix));
234 vvdw = _mm_sub_ps( _mm_mul_ps(vvdw12,one_twelfth) , _mm_mul_ps(vvdw6,one_sixth) );
235 fvdw = _mm_mul_ps(_mm_sub_ps(vvdw12,vvdw6),rinvsq00);
237 d = _mm_sub_ps(r00,rswitch);
238 d = _mm_max_ps(d,_mm_setzero_ps());
239 d2 = _mm_mul_ps(d,d);
240 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)))))));
242 dsw = _mm_mul_ps(d2,_mm_add_ps(swF2,_mm_mul_ps(d,_mm_add_ps(swF3,_mm_mul_ps(d,swF4)))));
244 /* Evaluate switch function */
245 /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
246 fvdw = _mm_sub_ps( _mm_mul_ps(fvdw,sw) , _mm_mul_ps(rinv00,_mm_mul_ps(vvdw,dsw)) );
247 vvdw = _mm_mul_ps(vvdw,sw);
248 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
250 /* Update potential sum for this i atom from the interaction with this j atom. */
251 velec = _mm_and_ps(velec,cutoff_mask);
252 velecsum = _mm_add_ps(velecsum,velec);
253 vvdw = _mm_and_ps(vvdw,cutoff_mask);
254 vvdwsum = _mm_add_ps(vvdwsum,vvdw);
256 fscal = _mm_add_ps(felec,fvdw);
258 fscal = _mm_and_ps(fscal,cutoff_mask);
260 /* Calculate temporary vectorial force */
261 tx = _mm_mul_ps(fscal,dx00);
262 ty = _mm_mul_ps(fscal,dy00);
263 tz = _mm_mul_ps(fscal,dz00);
265 /* Update vectorial force */
266 fix0 = _mm_add_ps(fix0,tx);
267 fiy0 = _mm_add_ps(fiy0,ty);
268 fiz0 = _mm_add_ps(fiz0,tz);
270 fjptrA = f+j_coord_offsetA;
271 fjptrB = f+j_coord_offsetB;
272 fjptrC = f+j_coord_offsetC;
273 fjptrD = f+j_coord_offsetD;
274 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
278 /* Inner loop uses 70 flops */
284 /* Get j neighbor index, and coordinate index */
285 jnrlistA = jjnr[jidx];
286 jnrlistB = jjnr[jidx+1];
287 jnrlistC = jjnr[jidx+2];
288 jnrlistD = jjnr[jidx+3];
289 /* Sign of each element will be negative for non-real atoms.
290 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
291 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
293 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
294 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
295 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
296 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
297 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
298 j_coord_offsetA = DIM*jnrA;
299 j_coord_offsetB = DIM*jnrB;
300 j_coord_offsetC = DIM*jnrC;
301 j_coord_offsetD = DIM*jnrD;
303 /* load j atom coordinates */
304 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
305 x+j_coord_offsetC,x+j_coord_offsetD,
308 /* Calculate displacement vector */
309 dx00 = _mm_sub_ps(ix0,jx0);
310 dy00 = _mm_sub_ps(iy0,jy0);
311 dz00 = _mm_sub_ps(iz0,jz0);
313 /* Calculate squared distance and things based on it */
314 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
316 rinv00 = gmx_mm_invsqrt_ps(rsq00);
318 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
320 /* Load parameters for j particles */
321 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
322 charge+jnrC+0,charge+jnrD+0);
323 vdwjidx0A = 2*vdwtype[jnrA+0];
324 vdwjidx0B = 2*vdwtype[jnrB+0];
325 vdwjidx0C = 2*vdwtype[jnrC+0];
326 vdwjidx0D = 2*vdwtype[jnrD+0];
328 /**************************
329 * CALCULATE INTERACTIONS *
330 **************************/
332 if (gmx_mm_any_lt(rsq00,rcutoff2))
335 r00 = _mm_mul_ps(rsq00,rinv00);
336 r00 = _mm_andnot_ps(dummy_mask,r00);
338 /* Compute parameters for interactions between i and j atoms */
339 qq00 = _mm_mul_ps(iq0,jq0);
340 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
341 vdwparam+vdwioffset0+vdwjidx0B,
342 vdwparam+vdwioffset0+vdwjidx0C,
343 vdwparam+vdwioffset0+vdwjidx0D,
346 /* REACTION-FIELD ELECTROSTATICS */
347 velec = _mm_mul_ps(qq00,_mm_sub_ps(_mm_add_ps(rinv00,_mm_mul_ps(krf,rsq00)),crf));
348 felec = _mm_mul_ps(qq00,_mm_sub_ps(_mm_mul_ps(rinv00,rinvsq00),krf2));
350 /* LENNARD-JONES DISPERSION/REPULSION */
352 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
353 vvdw6 = _mm_mul_ps(c6_00,rinvsix);
354 vvdw12 = _mm_mul_ps(c12_00,_mm_mul_ps(rinvsix,rinvsix));
355 vvdw = _mm_sub_ps( _mm_mul_ps(vvdw12,one_twelfth) , _mm_mul_ps(vvdw6,one_sixth) );
356 fvdw = _mm_mul_ps(_mm_sub_ps(vvdw12,vvdw6),rinvsq00);
358 d = _mm_sub_ps(r00,rswitch);
359 d = _mm_max_ps(d,_mm_setzero_ps());
360 d2 = _mm_mul_ps(d,d);
361 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)))))));
363 dsw = _mm_mul_ps(d2,_mm_add_ps(swF2,_mm_mul_ps(d,_mm_add_ps(swF3,_mm_mul_ps(d,swF4)))));
365 /* Evaluate switch function */
366 /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
367 fvdw = _mm_sub_ps( _mm_mul_ps(fvdw,sw) , _mm_mul_ps(rinv00,_mm_mul_ps(vvdw,dsw)) );
368 vvdw = _mm_mul_ps(vvdw,sw);
369 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
371 /* Update potential sum for this i atom from the interaction with this j atom. */
372 velec = _mm_and_ps(velec,cutoff_mask);
373 velec = _mm_andnot_ps(dummy_mask,velec);
374 velecsum = _mm_add_ps(velecsum,velec);
375 vvdw = _mm_and_ps(vvdw,cutoff_mask);
376 vvdw = _mm_andnot_ps(dummy_mask,vvdw);
377 vvdwsum = _mm_add_ps(vvdwsum,vvdw);
379 fscal = _mm_add_ps(felec,fvdw);
381 fscal = _mm_and_ps(fscal,cutoff_mask);
383 fscal = _mm_andnot_ps(dummy_mask,fscal);
385 /* Calculate temporary vectorial force */
386 tx = _mm_mul_ps(fscal,dx00);
387 ty = _mm_mul_ps(fscal,dy00);
388 tz = _mm_mul_ps(fscal,dz00);
390 /* Update vectorial force */
391 fix0 = _mm_add_ps(fix0,tx);
392 fiy0 = _mm_add_ps(fiy0,ty);
393 fiz0 = _mm_add_ps(fiz0,tz);
395 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
396 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
397 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
398 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
399 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
403 /* Inner loop uses 71 flops */
406 /* End of innermost loop */
408 gmx_mm_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
409 f+i_coord_offset,fshift+i_shift_offset);
412 /* Update potential energies */
413 gmx_mm_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
414 gmx_mm_update_1pot_ps(vvdwsum,kernel_data->energygrp_vdw+ggid);
416 /* Increment number of inner iterations */
417 inneriter += j_index_end - j_index_start;
419 /* Outer loop uses 9 flops */
422 /* Increment number of outer iterations */
425 /* Update outer/inner flops */
427 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_VF,outeriter*9 + inneriter*71);
430 * Gromacs nonbonded kernel: nb_kernel_ElecRFCut_VdwLJSw_GeomP1P1_F_sse2_single
431 * Electrostatics interaction: ReactionField
432 * VdW interaction: LennardJones
433 * Geometry: Particle-Particle
434 * Calculate force/pot: Force
437 nb_kernel_ElecRFCut_VdwLJSw_GeomP1P1_F_sse2_single
438 (t_nblist * gmx_restrict nlist,
439 rvec * gmx_restrict xx,
440 rvec * gmx_restrict ff,
441 t_forcerec * gmx_restrict fr,
442 t_mdatoms * gmx_restrict mdatoms,
443 nb_kernel_data_t * gmx_restrict kernel_data,
444 t_nrnb * gmx_restrict nrnb)
446 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
447 * just 0 for non-waters.
448 * Suffixes A,B,C,D refer to j loop unrolling done with SSE, e.g. for the four different
449 * jnr indices corresponding to data put in the four positions in the SIMD register.
451 int i_shift_offset,i_coord_offset,outeriter,inneriter;
452 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
453 int jnrA,jnrB,jnrC,jnrD;
454 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
455 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
456 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
458 real *shiftvec,*fshift,*x,*f;
459 real *fjptrA,*fjptrB,*fjptrC,*fjptrD;
461 __m128 tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
463 __m128 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
464 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
465 __m128 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
466 __m128 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
467 __m128 velec,felec,velecsum,facel,crf,krf,krf2;
470 __m128 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
473 __m128 one_sixth = _mm_set1_ps(1.0/6.0);
474 __m128 one_twelfth = _mm_set1_ps(1.0/12.0);
475 __m128 rswitch,swV3,swV4,swV5,swF2,swF3,swF4,d,d2,sw,dsw;
476 real rswitch_scalar,d_scalar;
477 __m128 dummy_mask,cutoff_mask;
478 __m128 signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
479 __m128 one = _mm_set1_ps(1.0);
480 __m128 two = _mm_set1_ps(2.0);
486 jindex = nlist->jindex;
488 shiftidx = nlist->shift;
490 shiftvec = fr->shift_vec[0];
491 fshift = fr->fshift[0];
492 facel = _mm_set1_ps(fr->epsfac);
493 charge = mdatoms->chargeA;
494 krf = _mm_set1_ps(fr->ic->k_rf);
495 krf2 = _mm_set1_ps(fr->ic->k_rf*2.0);
496 crf = _mm_set1_ps(fr->ic->c_rf);
497 nvdwtype = fr->ntype;
499 vdwtype = mdatoms->typeA;
501 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
502 rcutoff_scalar = fr->rcoulomb;
503 rcutoff = _mm_set1_ps(rcutoff_scalar);
504 rcutoff2 = _mm_mul_ps(rcutoff,rcutoff);
506 rswitch_scalar = fr->rvdw_switch;
507 rswitch = _mm_set1_ps(rswitch_scalar);
508 /* Setup switch parameters */
509 d_scalar = rcutoff_scalar-rswitch_scalar;
510 d = _mm_set1_ps(d_scalar);
511 swV3 = _mm_set1_ps(-10.0/(d_scalar*d_scalar*d_scalar));
512 swV4 = _mm_set1_ps( 15.0/(d_scalar*d_scalar*d_scalar*d_scalar));
513 swV5 = _mm_set1_ps( -6.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
514 swF2 = _mm_set1_ps(-30.0/(d_scalar*d_scalar*d_scalar));
515 swF3 = _mm_set1_ps( 60.0/(d_scalar*d_scalar*d_scalar*d_scalar));
516 swF4 = _mm_set1_ps(-30.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
518 /* Avoid stupid compiler warnings */
519 jnrA = jnrB = jnrC = jnrD = 0;
528 for(iidx=0;iidx<4*DIM;iidx++)
533 /* Start outer loop over neighborlists */
534 for(iidx=0; iidx<nri; iidx++)
536 /* Load shift vector for this list */
537 i_shift_offset = DIM*shiftidx[iidx];
539 /* Load limits for loop over neighbors */
540 j_index_start = jindex[iidx];
541 j_index_end = jindex[iidx+1];
543 /* Get outer coordinate index */
545 i_coord_offset = DIM*inr;
547 /* Load i particle coords and add shift vector */
548 gmx_mm_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
550 fix0 = _mm_setzero_ps();
551 fiy0 = _mm_setzero_ps();
552 fiz0 = _mm_setzero_ps();
554 /* Load parameters for i particles */
555 iq0 = _mm_mul_ps(facel,_mm_load1_ps(charge+inr+0));
556 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
558 /* Start inner kernel loop */
559 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
562 /* Get j neighbor index, and coordinate index */
567 j_coord_offsetA = DIM*jnrA;
568 j_coord_offsetB = DIM*jnrB;
569 j_coord_offsetC = DIM*jnrC;
570 j_coord_offsetD = DIM*jnrD;
572 /* load j atom coordinates */
573 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
574 x+j_coord_offsetC,x+j_coord_offsetD,
577 /* Calculate displacement vector */
578 dx00 = _mm_sub_ps(ix0,jx0);
579 dy00 = _mm_sub_ps(iy0,jy0);
580 dz00 = _mm_sub_ps(iz0,jz0);
582 /* Calculate squared distance and things based on it */
583 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
585 rinv00 = gmx_mm_invsqrt_ps(rsq00);
587 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
589 /* Load parameters for j particles */
590 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
591 charge+jnrC+0,charge+jnrD+0);
592 vdwjidx0A = 2*vdwtype[jnrA+0];
593 vdwjidx0B = 2*vdwtype[jnrB+0];
594 vdwjidx0C = 2*vdwtype[jnrC+0];
595 vdwjidx0D = 2*vdwtype[jnrD+0];
597 /**************************
598 * CALCULATE INTERACTIONS *
599 **************************/
601 if (gmx_mm_any_lt(rsq00,rcutoff2))
604 r00 = _mm_mul_ps(rsq00,rinv00);
606 /* Compute parameters for interactions between i and j atoms */
607 qq00 = _mm_mul_ps(iq0,jq0);
608 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
609 vdwparam+vdwioffset0+vdwjidx0B,
610 vdwparam+vdwioffset0+vdwjidx0C,
611 vdwparam+vdwioffset0+vdwjidx0D,
614 /* REACTION-FIELD ELECTROSTATICS */
615 felec = _mm_mul_ps(qq00,_mm_sub_ps(_mm_mul_ps(rinv00,rinvsq00),krf2));
617 /* LENNARD-JONES DISPERSION/REPULSION */
619 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
620 vvdw6 = _mm_mul_ps(c6_00,rinvsix);
621 vvdw12 = _mm_mul_ps(c12_00,_mm_mul_ps(rinvsix,rinvsix));
622 vvdw = _mm_sub_ps( _mm_mul_ps(vvdw12,one_twelfth) , _mm_mul_ps(vvdw6,one_sixth) );
623 fvdw = _mm_mul_ps(_mm_sub_ps(vvdw12,vvdw6),rinvsq00);
625 d = _mm_sub_ps(r00,rswitch);
626 d = _mm_max_ps(d,_mm_setzero_ps());
627 d2 = _mm_mul_ps(d,d);
628 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)))))));
630 dsw = _mm_mul_ps(d2,_mm_add_ps(swF2,_mm_mul_ps(d,_mm_add_ps(swF3,_mm_mul_ps(d,swF4)))));
632 /* Evaluate switch function */
633 /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
634 fvdw = _mm_sub_ps( _mm_mul_ps(fvdw,sw) , _mm_mul_ps(rinv00,_mm_mul_ps(vvdw,dsw)) );
635 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
637 fscal = _mm_add_ps(felec,fvdw);
639 fscal = _mm_and_ps(fscal,cutoff_mask);
641 /* Calculate temporary vectorial force */
642 tx = _mm_mul_ps(fscal,dx00);
643 ty = _mm_mul_ps(fscal,dy00);
644 tz = _mm_mul_ps(fscal,dz00);
646 /* Update vectorial force */
647 fix0 = _mm_add_ps(fix0,tx);
648 fiy0 = _mm_add_ps(fiy0,ty);
649 fiz0 = _mm_add_ps(fiz0,tz);
651 fjptrA = f+j_coord_offsetA;
652 fjptrB = f+j_coord_offsetB;
653 fjptrC = f+j_coord_offsetC;
654 fjptrD = f+j_coord_offsetD;
655 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
659 /* Inner loop uses 61 flops */
665 /* Get j neighbor index, and coordinate index */
666 jnrlistA = jjnr[jidx];
667 jnrlistB = jjnr[jidx+1];
668 jnrlistC = jjnr[jidx+2];
669 jnrlistD = jjnr[jidx+3];
670 /* Sign of each element will be negative for non-real atoms.
671 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
672 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
674 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
675 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
676 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
677 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
678 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
679 j_coord_offsetA = DIM*jnrA;
680 j_coord_offsetB = DIM*jnrB;
681 j_coord_offsetC = DIM*jnrC;
682 j_coord_offsetD = DIM*jnrD;
684 /* load j atom coordinates */
685 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
686 x+j_coord_offsetC,x+j_coord_offsetD,
689 /* Calculate displacement vector */
690 dx00 = _mm_sub_ps(ix0,jx0);
691 dy00 = _mm_sub_ps(iy0,jy0);
692 dz00 = _mm_sub_ps(iz0,jz0);
694 /* Calculate squared distance and things based on it */
695 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
697 rinv00 = gmx_mm_invsqrt_ps(rsq00);
699 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
701 /* Load parameters for j particles */
702 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
703 charge+jnrC+0,charge+jnrD+0);
704 vdwjidx0A = 2*vdwtype[jnrA+0];
705 vdwjidx0B = 2*vdwtype[jnrB+0];
706 vdwjidx0C = 2*vdwtype[jnrC+0];
707 vdwjidx0D = 2*vdwtype[jnrD+0];
709 /**************************
710 * CALCULATE INTERACTIONS *
711 **************************/
713 if (gmx_mm_any_lt(rsq00,rcutoff2))
716 r00 = _mm_mul_ps(rsq00,rinv00);
717 r00 = _mm_andnot_ps(dummy_mask,r00);
719 /* Compute parameters for interactions between i and j atoms */
720 qq00 = _mm_mul_ps(iq0,jq0);
721 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
722 vdwparam+vdwioffset0+vdwjidx0B,
723 vdwparam+vdwioffset0+vdwjidx0C,
724 vdwparam+vdwioffset0+vdwjidx0D,
727 /* REACTION-FIELD ELECTROSTATICS */
728 felec = _mm_mul_ps(qq00,_mm_sub_ps(_mm_mul_ps(rinv00,rinvsq00),krf2));
730 /* LENNARD-JONES DISPERSION/REPULSION */
732 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
733 vvdw6 = _mm_mul_ps(c6_00,rinvsix);
734 vvdw12 = _mm_mul_ps(c12_00,_mm_mul_ps(rinvsix,rinvsix));
735 vvdw = _mm_sub_ps( _mm_mul_ps(vvdw12,one_twelfth) , _mm_mul_ps(vvdw6,one_sixth) );
736 fvdw = _mm_mul_ps(_mm_sub_ps(vvdw12,vvdw6),rinvsq00);
738 d = _mm_sub_ps(r00,rswitch);
739 d = _mm_max_ps(d,_mm_setzero_ps());
740 d2 = _mm_mul_ps(d,d);
741 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)))))));
743 dsw = _mm_mul_ps(d2,_mm_add_ps(swF2,_mm_mul_ps(d,_mm_add_ps(swF3,_mm_mul_ps(d,swF4)))));
745 /* Evaluate switch function */
746 /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
747 fvdw = _mm_sub_ps( _mm_mul_ps(fvdw,sw) , _mm_mul_ps(rinv00,_mm_mul_ps(vvdw,dsw)) );
748 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
750 fscal = _mm_add_ps(felec,fvdw);
752 fscal = _mm_and_ps(fscal,cutoff_mask);
754 fscal = _mm_andnot_ps(dummy_mask,fscal);
756 /* Calculate temporary vectorial force */
757 tx = _mm_mul_ps(fscal,dx00);
758 ty = _mm_mul_ps(fscal,dy00);
759 tz = _mm_mul_ps(fscal,dz00);
761 /* Update vectorial force */
762 fix0 = _mm_add_ps(fix0,tx);
763 fiy0 = _mm_add_ps(fiy0,ty);
764 fiz0 = _mm_add_ps(fiz0,tz);
766 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
767 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
768 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
769 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
770 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
774 /* Inner loop uses 62 flops */
777 /* End of innermost loop */
779 gmx_mm_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
780 f+i_coord_offset,fshift+i_shift_offset);
782 /* Increment number of inner iterations */
783 inneriter += j_index_end - j_index_start;
785 /* Outer loop uses 7 flops */
788 /* Increment number of outer iterations */
791 /* Update outer/inner flops */
793 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_F,outeriter*7 + inneriter*62);