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_ElecEwSw_VdwNone_GeomP1P1_VF_sse2_single
38 * Electrostatics interaction: Ewald
39 * VdW interaction: None
40 * Geometry: Particle-Particle
41 * Calculate force/pot: PotentialAndForce
44 nb_kernel_ElecEwSw_VdwNone_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 j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
62 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
63 real shX,shY,shZ,rcutoff_scalar;
64 real *shiftvec,*fshift,*x,*f;
65 __m128 tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
67 __m128 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
68 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
69 __m128 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
70 __m128 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
71 __m128 velec,felec,velecsum,facel,crf,krf,krf2;
74 __m128 ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
76 __m128 rswitch,swV3,swV4,swV5,swF2,swF3,swF4,d,d2,sw,dsw;
77 real rswitch_scalar,d_scalar;
78 __m128 dummy_mask,cutoff_mask;
79 __m128 signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
80 __m128 one = _mm_set1_ps(1.0);
81 __m128 two = _mm_set1_ps(2.0);
87 jindex = nlist->jindex;
89 shiftidx = nlist->shift;
91 shiftvec = fr->shift_vec[0];
92 fshift = fr->fshift[0];
93 facel = _mm_set1_ps(fr->epsfac);
94 charge = mdatoms->chargeA;
96 sh_ewald = _mm_set1_ps(fr->ic->sh_ewald);
97 ewtab = fr->ic->tabq_coul_FDV0;
98 ewtabscale = _mm_set1_ps(fr->ic->tabq_scale);
99 ewtabhalfspace = _mm_set1_ps(0.5/fr->ic->tabq_scale);
101 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
102 rcutoff_scalar = fr->rcoulomb;
103 rcutoff = _mm_set1_ps(rcutoff_scalar);
104 rcutoff2 = _mm_mul_ps(rcutoff,rcutoff);
106 rswitch_scalar = fr->rcoulomb_switch;
107 rswitch = _mm_set1_ps(rswitch_scalar);
108 /* Setup switch parameters */
109 d_scalar = rcutoff_scalar-rswitch_scalar;
110 d = _mm_set1_ps(d_scalar);
111 swV3 = _mm_set1_ps(-10.0/(d_scalar*d_scalar*d_scalar));
112 swV4 = _mm_set1_ps( 15.0/(d_scalar*d_scalar*d_scalar*d_scalar));
113 swV5 = _mm_set1_ps( -6.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
114 swF2 = _mm_set1_ps(-30.0/(d_scalar*d_scalar*d_scalar));
115 swF3 = _mm_set1_ps( 60.0/(d_scalar*d_scalar*d_scalar*d_scalar));
116 swF4 = _mm_set1_ps(-30.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
118 /* Avoid stupid compiler warnings */
119 jnrA = jnrB = jnrC = jnrD = 0;
128 /* Start outer loop over neighborlists */
129 for(iidx=0; iidx<nri; iidx++)
131 /* Load shift vector for this list */
132 i_shift_offset = DIM*shiftidx[iidx];
133 shX = shiftvec[i_shift_offset+XX];
134 shY = shiftvec[i_shift_offset+YY];
135 shZ = shiftvec[i_shift_offset+ZZ];
137 /* Load limits for loop over neighbors */
138 j_index_start = jindex[iidx];
139 j_index_end = jindex[iidx+1];
141 /* Get outer coordinate index */
143 i_coord_offset = DIM*inr;
145 /* Load i particle coords and add shift vector */
146 ix0 = _mm_set1_ps(shX + x[i_coord_offset+DIM*0+XX]);
147 iy0 = _mm_set1_ps(shY + x[i_coord_offset+DIM*0+YY]);
148 iz0 = _mm_set1_ps(shZ + x[i_coord_offset+DIM*0+ZZ]);
150 fix0 = _mm_setzero_ps();
151 fiy0 = _mm_setzero_ps();
152 fiz0 = _mm_setzero_ps();
154 /* Load parameters for i particles */
155 iq0 = _mm_mul_ps(facel,_mm_load1_ps(charge+inr+0));
157 /* Reset potential sums */
158 velecsum = _mm_setzero_ps();
160 /* Start inner kernel loop */
161 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
164 /* Get j neighbor index, and coordinate index */
170 j_coord_offsetA = DIM*jnrA;
171 j_coord_offsetB = DIM*jnrB;
172 j_coord_offsetC = DIM*jnrC;
173 j_coord_offsetD = DIM*jnrD;
175 /* load j atom coordinates */
176 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
177 x+j_coord_offsetC,x+j_coord_offsetD,
180 /* Calculate displacement vector */
181 dx00 = _mm_sub_ps(ix0,jx0);
182 dy00 = _mm_sub_ps(iy0,jy0);
183 dz00 = _mm_sub_ps(iz0,jz0);
185 /* Calculate squared distance and things based on it */
186 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
188 rinv00 = gmx_mm_invsqrt_ps(rsq00);
190 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
192 /* Load parameters for j particles */
193 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
194 charge+jnrC+0,charge+jnrD+0);
196 /**************************
197 * CALCULATE INTERACTIONS *
198 **************************/
200 if (gmx_mm_any_lt(rsq00,rcutoff2))
203 r00 = _mm_mul_ps(rsq00,rinv00);
205 /* Compute parameters for interactions between i and j atoms */
206 qq00 = _mm_mul_ps(iq0,jq0);
208 /* EWALD ELECTROSTATICS */
210 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
211 ewrt = _mm_mul_ps(r00,ewtabscale);
212 ewitab = _mm_cvttps_epi32(ewrt);
213 eweps = _mm_sub_ps(ewrt,_mm_cvtepi32_ps(ewitab));
214 ewitab = _mm_slli_epi32(ewitab,2);
215 ewtabF = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,0) );
216 ewtabD = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,1) );
217 ewtabV = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,2) );
218 ewtabFn = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,3) );
219 _MM_TRANSPOSE4_PS(ewtabF,ewtabD,ewtabV,ewtabFn);
220 felec = _mm_add_ps(ewtabF,_mm_mul_ps(eweps,ewtabD));
221 velec = _mm_sub_ps(ewtabV,_mm_mul_ps(_mm_mul_ps(ewtabhalfspace,eweps),_mm_add_ps(ewtabF,felec)));
222 velec = _mm_mul_ps(qq00,_mm_sub_ps(rinv00,velec));
223 felec = _mm_mul_ps(_mm_mul_ps(qq00,rinv00),_mm_sub_ps(rinvsq00,felec));
225 d = _mm_sub_ps(r00,rswitch);
226 d = _mm_max_ps(d,_mm_setzero_ps());
227 d2 = _mm_mul_ps(d,d);
228 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)))))));
230 dsw = _mm_mul_ps(d2,_mm_add_ps(swF2,_mm_mul_ps(d,_mm_add_ps(swF3,_mm_mul_ps(d,swF4)))));
232 /* Evaluate switch function */
233 /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
234 felec = _mm_sub_ps( _mm_mul_ps(felec,sw) , _mm_mul_ps(rinv00,_mm_mul_ps(velec,dsw)) );
235 velec = _mm_mul_ps(velec,sw);
236 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
238 /* Update potential sum for this i atom from the interaction with this j atom. */
239 velec = _mm_and_ps(velec,cutoff_mask);
240 velecsum = _mm_add_ps(velecsum,velec);
244 fscal = _mm_and_ps(fscal,cutoff_mask);
246 /* Calculate temporary vectorial force */
247 tx = _mm_mul_ps(fscal,dx00);
248 ty = _mm_mul_ps(fscal,dy00);
249 tz = _mm_mul_ps(fscal,dz00);
251 /* Update vectorial force */
252 fix0 = _mm_add_ps(fix0,tx);
253 fiy0 = _mm_add_ps(fiy0,ty);
254 fiz0 = _mm_add_ps(fiz0,tz);
256 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(f+j_coord_offsetA,f+j_coord_offsetB,
257 f+j_coord_offsetC,f+j_coord_offsetD,
262 /* Inner loop uses 65 flops */
268 /* Get j neighbor index, and coordinate index */
274 /* Sign of each element will be negative for non-real atoms.
275 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
276 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
278 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
279 jnrA = (jnrA>=0) ? jnrA : 0;
280 jnrB = (jnrB>=0) ? jnrB : 0;
281 jnrC = (jnrC>=0) ? jnrC : 0;
282 jnrD = (jnrD>=0) ? jnrD : 0;
284 j_coord_offsetA = DIM*jnrA;
285 j_coord_offsetB = DIM*jnrB;
286 j_coord_offsetC = DIM*jnrC;
287 j_coord_offsetD = DIM*jnrD;
289 /* load j atom coordinates */
290 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
291 x+j_coord_offsetC,x+j_coord_offsetD,
294 /* Calculate displacement vector */
295 dx00 = _mm_sub_ps(ix0,jx0);
296 dy00 = _mm_sub_ps(iy0,jy0);
297 dz00 = _mm_sub_ps(iz0,jz0);
299 /* Calculate squared distance and things based on it */
300 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
302 rinv00 = gmx_mm_invsqrt_ps(rsq00);
304 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
306 /* Load parameters for j particles */
307 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
308 charge+jnrC+0,charge+jnrD+0);
310 /**************************
311 * CALCULATE INTERACTIONS *
312 **************************/
314 if (gmx_mm_any_lt(rsq00,rcutoff2))
317 r00 = _mm_mul_ps(rsq00,rinv00);
318 r00 = _mm_andnot_ps(dummy_mask,r00);
320 /* Compute parameters for interactions between i and j atoms */
321 qq00 = _mm_mul_ps(iq0,jq0);
323 /* EWALD ELECTROSTATICS */
325 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
326 ewrt = _mm_mul_ps(r00,ewtabscale);
327 ewitab = _mm_cvttps_epi32(ewrt);
328 eweps = _mm_sub_ps(ewrt,_mm_cvtepi32_ps(ewitab));
329 ewitab = _mm_slli_epi32(ewitab,2);
330 ewtabF = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,0) );
331 ewtabD = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,1) );
332 ewtabV = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,2) );
333 ewtabFn = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,3) );
334 _MM_TRANSPOSE4_PS(ewtabF,ewtabD,ewtabV,ewtabFn);
335 felec = _mm_add_ps(ewtabF,_mm_mul_ps(eweps,ewtabD));
336 velec = _mm_sub_ps(ewtabV,_mm_mul_ps(_mm_mul_ps(ewtabhalfspace,eweps),_mm_add_ps(ewtabF,felec)));
337 velec = _mm_mul_ps(qq00,_mm_sub_ps(rinv00,velec));
338 felec = _mm_mul_ps(_mm_mul_ps(qq00,rinv00),_mm_sub_ps(rinvsq00,felec));
340 d = _mm_sub_ps(r00,rswitch);
341 d = _mm_max_ps(d,_mm_setzero_ps());
342 d2 = _mm_mul_ps(d,d);
343 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)))))));
345 dsw = _mm_mul_ps(d2,_mm_add_ps(swF2,_mm_mul_ps(d,_mm_add_ps(swF3,_mm_mul_ps(d,swF4)))));
347 /* Evaluate switch function */
348 /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
349 felec = _mm_sub_ps( _mm_mul_ps(felec,sw) , _mm_mul_ps(rinv00,_mm_mul_ps(velec,dsw)) );
350 velec = _mm_mul_ps(velec,sw);
351 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
353 /* Update potential sum for this i atom from the interaction with this j atom. */
354 velec = _mm_and_ps(velec,cutoff_mask);
355 velec = _mm_andnot_ps(dummy_mask,velec);
356 velecsum = _mm_add_ps(velecsum,velec);
360 fscal = _mm_and_ps(fscal,cutoff_mask);
362 fscal = _mm_andnot_ps(dummy_mask,fscal);
364 /* Calculate temporary vectorial force */
365 tx = _mm_mul_ps(fscal,dx00);
366 ty = _mm_mul_ps(fscal,dy00);
367 tz = _mm_mul_ps(fscal,dz00);
369 /* Update vectorial force */
370 fix0 = _mm_add_ps(fix0,tx);
371 fiy0 = _mm_add_ps(fiy0,ty);
372 fiz0 = _mm_add_ps(fiz0,tz);
374 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(f+j_coord_offsetA,f+j_coord_offsetB,
375 f+j_coord_offsetC,f+j_coord_offsetD,
380 /* Inner loop uses 66 flops */
383 /* End of innermost loop */
385 gmx_mm_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
386 f+i_coord_offset,fshift+i_shift_offset);
389 /* Update potential energies */
390 gmx_mm_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
392 /* Increment number of inner iterations */
393 inneriter += j_index_end - j_index_start;
395 /* Outer loop uses 11 flops */
398 /* Increment number of outer iterations */
401 /* Update outer/inner flops */
403 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VF,outeriter*11 + inneriter*66);
406 * Gromacs nonbonded kernel: nb_kernel_ElecEwSw_VdwNone_GeomP1P1_F_sse2_single
407 * Electrostatics interaction: Ewald
408 * VdW interaction: None
409 * Geometry: Particle-Particle
410 * Calculate force/pot: Force
413 nb_kernel_ElecEwSw_VdwNone_GeomP1P1_F_sse2_single
414 (t_nblist * gmx_restrict nlist,
415 rvec * gmx_restrict xx,
416 rvec * gmx_restrict ff,
417 t_forcerec * gmx_restrict fr,
418 t_mdatoms * gmx_restrict mdatoms,
419 nb_kernel_data_t * gmx_restrict kernel_data,
420 t_nrnb * gmx_restrict nrnb)
422 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
423 * just 0 for non-waters.
424 * Suffixes A,B,C,D refer to j loop unrolling done with SSE, e.g. for the four different
425 * jnr indices corresponding to data put in the four positions in the SIMD register.
427 int i_shift_offset,i_coord_offset,outeriter,inneriter;
428 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
429 int jnrA,jnrB,jnrC,jnrD;
430 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
431 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
432 real shX,shY,shZ,rcutoff_scalar;
433 real *shiftvec,*fshift,*x,*f;
434 __m128 tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
436 __m128 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
437 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
438 __m128 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
439 __m128 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
440 __m128 velec,felec,velecsum,facel,crf,krf,krf2;
443 __m128 ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
445 __m128 rswitch,swV3,swV4,swV5,swF2,swF3,swF4,d,d2,sw,dsw;
446 real rswitch_scalar,d_scalar;
447 __m128 dummy_mask,cutoff_mask;
448 __m128 signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
449 __m128 one = _mm_set1_ps(1.0);
450 __m128 two = _mm_set1_ps(2.0);
456 jindex = nlist->jindex;
458 shiftidx = nlist->shift;
460 shiftvec = fr->shift_vec[0];
461 fshift = fr->fshift[0];
462 facel = _mm_set1_ps(fr->epsfac);
463 charge = mdatoms->chargeA;
465 sh_ewald = _mm_set1_ps(fr->ic->sh_ewald);
466 ewtab = fr->ic->tabq_coul_FDV0;
467 ewtabscale = _mm_set1_ps(fr->ic->tabq_scale);
468 ewtabhalfspace = _mm_set1_ps(0.5/fr->ic->tabq_scale);
470 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
471 rcutoff_scalar = fr->rcoulomb;
472 rcutoff = _mm_set1_ps(rcutoff_scalar);
473 rcutoff2 = _mm_mul_ps(rcutoff,rcutoff);
475 rswitch_scalar = fr->rcoulomb_switch;
476 rswitch = _mm_set1_ps(rswitch_scalar);
477 /* Setup switch parameters */
478 d_scalar = rcutoff_scalar-rswitch_scalar;
479 d = _mm_set1_ps(d_scalar);
480 swV3 = _mm_set1_ps(-10.0/(d_scalar*d_scalar*d_scalar));
481 swV4 = _mm_set1_ps( 15.0/(d_scalar*d_scalar*d_scalar*d_scalar));
482 swV5 = _mm_set1_ps( -6.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
483 swF2 = _mm_set1_ps(-30.0/(d_scalar*d_scalar*d_scalar));
484 swF3 = _mm_set1_ps( 60.0/(d_scalar*d_scalar*d_scalar*d_scalar));
485 swF4 = _mm_set1_ps(-30.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
487 /* Avoid stupid compiler warnings */
488 jnrA = jnrB = jnrC = jnrD = 0;
497 /* Start outer loop over neighborlists */
498 for(iidx=0; iidx<nri; iidx++)
500 /* Load shift vector for this list */
501 i_shift_offset = DIM*shiftidx[iidx];
502 shX = shiftvec[i_shift_offset+XX];
503 shY = shiftvec[i_shift_offset+YY];
504 shZ = shiftvec[i_shift_offset+ZZ];
506 /* Load limits for loop over neighbors */
507 j_index_start = jindex[iidx];
508 j_index_end = jindex[iidx+1];
510 /* Get outer coordinate index */
512 i_coord_offset = DIM*inr;
514 /* Load i particle coords and add shift vector */
515 ix0 = _mm_set1_ps(shX + x[i_coord_offset+DIM*0+XX]);
516 iy0 = _mm_set1_ps(shY + x[i_coord_offset+DIM*0+YY]);
517 iz0 = _mm_set1_ps(shZ + x[i_coord_offset+DIM*0+ZZ]);
519 fix0 = _mm_setzero_ps();
520 fiy0 = _mm_setzero_ps();
521 fiz0 = _mm_setzero_ps();
523 /* Load parameters for i particles */
524 iq0 = _mm_mul_ps(facel,_mm_load1_ps(charge+inr+0));
526 /* Start inner kernel loop */
527 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
530 /* Get j neighbor index, and coordinate index */
536 j_coord_offsetA = DIM*jnrA;
537 j_coord_offsetB = DIM*jnrB;
538 j_coord_offsetC = DIM*jnrC;
539 j_coord_offsetD = DIM*jnrD;
541 /* load j atom coordinates */
542 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
543 x+j_coord_offsetC,x+j_coord_offsetD,
546 /* Calculate displacement vector */
547 dx00 = _mm_sub_ps(ix0,jx0);
548 dy00 = _mm_sub_ps(iy0,jy0);
549 dz00 = _mm_sub_ps(iz0,jz0);
551 /* Calculate squared distance and things based on it */
552 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
554 rinv00 = gmx_mm_invsqrt_ps(rsq00);
556 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
558 /* Load parameters for j particles */
559 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
560 charge+jnrC+0,charge+jnrD+0);
562 /**************************
563 * CALCULATE INTERACTIONS *
564 **************************/
566 if (gmx_mm_any_lt(rsq00,rcutoff2))
569 r00 = _mm_mul_ps(rsq00,rinv00);
571 /* Compute parameters for interactions between i and j atoms */
572 qq00 = _mm_mul_ps(iq0,jq0);
574 /* EWALD ELECTROSTATICS */
576 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
577 ewrt = _mm_mul_ps(r00,ewtabscale);
578 ewitab = _mm_cvttps_epi32(ewrt);
579 eweps = _mm_sub_ps(ewrt,_mm_cvtepi32_ps(ewitab));
580 ewitab = _mm_slli_epi32(ewitab,2);
581 ewtabF = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,0) );
582 ewtabD = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,1) );
583 ewtabV = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,2) );
584 ewtabFn = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,3) );
585 _MM_TRANSPOSE4_PS(ewtabF,ewtabD,ewtabV,ewtabFn);
586 felec = _mm_add_ps(ewtabF,_mm_mul_ps(eweps,ewtabD));
587 velec = _mm_sub_ps(ewtabV,_mm_mul_ps(_mm_mul_ps(ewtabhalfspace,eweps),_mm_add_ps(ewtabF,felec)));
588 velec = _mm_mul_ps(qq00,_mm_sub_ps(rinv00,velec));
589 felec = _mm_mul_ps(_mm_mul_ps(qq00,rinv00),_mm_sub_ps(rinvsq00,felec));
591 d = _mm_sub_ps(r00,rswitch);
592 d = _mm_max_ps(d,_mm_setzero_ps());
593 d2 = _mm_mul_ps(d,d);
594 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)))))));
596 dsw = _mm_mul_ps(d2,_mm_add_ps(swF2,_mm_mul_ps(d,_mm_add_ps(swF3,_mm_mul_ps(d,swF4)))));
598 /* Evaluate switch function */
599 /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
600 felec = _mm_sub_ps( _mm_mul_ps(felec,sw) , _mm_mul_ps(rinv00,_mm_mul_ps(velec,dsw)) );
601 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
605 fscal = _mm_and_ps(fscal,cutoff_mask);
607 /* Calculate temporary vectorial force */
608 tx = _mm_mul_ps(fscal,dx00);
609 ty = _mm_mul_ps(fscal,dy00);
610 tz = _mm_mul_ps(fscal,dz00);
612 /* Update vectorial force */
613 fix0 = _mm_add_ps(fix0,tx);
614 fiy0 = _mm_add_ps(fiy0,ty);
615 fiz0 = _mm_add_ps(fiz0,tz);
617 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(f+j_coord_offsetA,f+j_coord_offsetB,
618 f+j_coord_offsetC,f+j_coord_offsetD,
623 /* Inner loop uses 62 flops */
629 /* Get j neighbor index, and coordinate index */
635 /* Sign of each element will be negative for non-real atoms.
636 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
637 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
639 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
640 jnrA = (jnrA>=0) ? jnrA : 0;
641 jnrB = (jnrB>=0) ? jnrB : 0;
642 jnrC = (jnrC>=0) ? jnrC : 0;
643 jnrD = (jnrD>=0) ? jnrD : 0;
645 j_coord_offsetA = DIM*jnrA;
646 j_coord_offsetB = DIM*jnrB;
647 j_coord_offsetC = DIM*jnrC;
648 j_coord_offsetD = DIM*jnrD;
650 /* load j atom coordinates */
651 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
652 x+j_coord_offsetC,x+j_coord_offsetD,
655 /* Calculate displacement vector */
656 dx00 = _mm_sub_ps(ix0,jx0);
657 dy00 = _mm_sub_ps(iy0,jy0);
658 dz00 = _mm_sub_ps(iz0,jz0);
660 /* Calculate squared distance and things based on it */
661 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
663 rinv00 = gmx_mm_invsqrt_ps(rsq00);
665 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
667 /* Load parameters for j particles */
668 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
669 charge+jnrC+0,charge+jnrD+0);
671 /**************************
672 * CALCULATE INTERACTIONS *
673 **************************/
675 if (gmx_mm_any_lt(rsq00,rcutoff2))
678 r00 = _mm_mul_ps(rsq00,rinv00);
679 r00 = _mm_andnot_ps(dummy_mask,r00);
681 /* Compute parameters for interactions between i and j atoms */
682 qq00 = _mm_mul_ps(iq0,jq0);
684 /* EWALD ELECTROSTATICS */
686 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
687 ewrt = _mm_mul_ps(r00,ewtabscale);
688 ewitab = _mm_cvttps_epi32(ewrt);
689 eweps = _mm_sub_ps(ewrt,_mm_cvtepi32_ps(ewitab));
690 ewitab = _mm_slli_epi32(ewitab,2);
691 ewtabF = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,0) );
692 ewtabD = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,1) );
693 ewtabV = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,2) );
694 ewtabFn = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,3) );
695 _MM_TRANSPOSE4_PS(ewtabF,ewtabD,ewtabV,ewtabFn);
696 felec = _mm_add_ps(ewtabF,_mm_mul_ps(eweps,ewtabD));
697 velec = _mm_sub_ps(ewtabV,_mm_mul_ps(_mm_mul_ps(ewtabhalfspace,eweps),_mm_add_ps(ewtabF,felec)));
698 velec = _mm_mul_ps(qq00,_mm_sub_ps(rinv00,velec));
699 felec = _mm_mul_ps(_mm_mul_ps(qq00,rinv00),_mm_sub_ps(rinvsq00,felec));
701 d = _mm_sub_ps(r00,rswitch);
702 d = _mm_max_ps(d,_mm_setzero_ps());
703 d2 = _mm_mul_ps(d,d);
704 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)))))));
706 dsw = _mm_mul_ps(d2,_mm_add_ps(swF2,_mm_mul_ps(d,_mm_add_ps(swF3,_mm_mul_ps(d,swF4)))));
708 /* Evaluate switch function */
709 /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
710 felec = _mm_sub_ps( _mm_mul_ps(felec,sw) , _mm_mul_ps(rinv00,_mm_mul_ps(velec,dsw)) );
711 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
715 fscal = _mm_and_ps(fscal,cutoff_mask);
717 fscal = _mm_andnot_ps(dummy_mask,fscal);
719 /* Calculate temporary vectorial force */
720 tx = _mm_mul_ps(fscal,dx00);
721 ty = _mm_mul_ps(fscal,dy00);
722 tz = _mm_mul_ps(fscal,dz00);
724 /* Update vectorial force */
725 fix0 = _mm_add_ps(fix0,tx);
726 fiy0 = _mm_add_ps(fiy0,ty);
727 fiz0 = _mm_add_ps(fiz0,tz);
729 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(f+j_coord_offsetA,f+j_coord_offsetB,
730 f+j_coord_offsetC,f+j_coord_offsetD,
735 /* Inner loop uses 63 flops */
738 /* End of innermost loop */
740 gmx_mm_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
741 f+i_coord_offset,fshift+i_shift_offset);
743 /* Increment number of inner iterations */
744 inneriter += j_index_end - j_index_start;
746 /* Outer loop uses 10 flops */
749 /* Increment number of outer iterations */
752 /* Update outer/inner flops */
754 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_F,outeriter*10 + inneriter*63);