2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
36 * Note: this file was generated by the GROMACS sse4_1_single kernel generator.
42 #include "../nb_kernel.h"
43 #include "gromacs/legacyheaders/types/simple.h"
44 #include "gromacs/math/vec.h"
45 #include "gromacs/legacyheaders/nrnb.h"
47 #include "gromacs/simd/math_x86_sse4_1_single.h"
48 #include "kernelutil_x86_sse4_1_single.h"
51 * Gromacs nonbonded kernel: nb_kernel_ElecEwSw_VdwNone_GeomP1P1_VF_sse4_1_single
52 * Electrostatics interaction: Ewald
53 * VdW interaction: None
54 * Geometry: Particle-Particle
55 * Calculate force/pot: PotentialAndForce
58 nb_kernel_ElecEwSw_VdwNone_GeomP1P1_VF_sse4_1_single
59 (t_nblist * gmx_restrict nlist,
60 rvec * gmx_restrict xx,
61 rvec * gmx_restrict ff,
62 t_forcerec * gmx_restrict fr,
63 t_mdatoms * gmx_restrict mdatoms,
64 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
65 t_nrnb * gmx_restrict nrnb)
67 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
68 * just 0 for non-waters.
69 * Suffixes A,B,C,D refer to j loop unrolling done with SSE, e.g. for the four different
70 * jnr indices corresponding to data put in the four positions in the SIMD register.
72 int i_shift_offset,i_coord_offset,outeriter,inneriter;
73 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
74 int jnrA,jnrB,jnrC,jnrD;
75 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
76 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
77 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
79 real *shiftvec,*fshift,*x,*f;
80 real *fjptrA,*fjptrB,*fjptrC,*fjptrD;
82 __m128 tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
84 __m128 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
85 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
86 __m128 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
87 __m128 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
88 __m128 velec,felec,velecsum,facel,crf,krf,krf2;
91 __m128 ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
93 __m128 rswitch,swV3,swV4,swV5,swF2,swF3,swF4,d,d2,sw,dsw;
94 real rswitch_scalar,d_scalar;
95 __m128 dummy_mask,cutoff_mask;
96 __m128 signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
97 __m128 one = _mm_set1_ps(1.0);
98 __m128 two = _mm_set1_ps(2.0);
104 jindex = nlist->jindex;
106 shiftidx = nlist->shift;
108 shiftvec = fr->shift_vec[0];
109 fshift = fr->fshift[0];
110 facel = _mm_set1_ps(fr->epsfac);
111 charge = mdatoms->chargeA;
113 sh_ewald = _mm_set1_ps(fr->ic->sh_ewald);
114 ewtab = fr->ic->tabq_coul_FDV0;
115 ewtabscale = _mm_set1_ps(fr->ic->tabq_scale);
116 ewtabhalfspace = _mm_set1_ps(0.5/fr->ic->tabq_scale);
118 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
119 rcutoff_scalar = fr->rcoulomb;
120 rcutoff = _mm_set1_ps(rcutoff_scalar);
121 rcutoff2 = _mm_mul_ps(rcutoff,rcutoff);
123 rswitch_scalar = fr->rcoulomb_switch;
124 rswitch = _mm_set1_ps(rswitch_scalar);
125 /* Setup switch parameters */
126 d_scalar = rcutoff_scalar-rswitch_scalar;
127 d = _mm_set1_ps(d_scalar);
128 swV3 = _mm_set1_ps(-10.0/(d_scalar*d_scalar*d_scalar));
129 swV4 = _mm_set1_ps( 15.0/(d_scalar*d_scalar*d_scalar*d_scalar));
130 swV5 = _mm_set1_ps( -6.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
131 swF2 = _mm_set1_ps(-30.0/(d_scalar*d_scalar*d_scalar));
132 swF3 = _mm_set1_ps( 60.0/(d_scalar*d_scalar*d_scalar*d_scalar));
133 swF4 = _mm_set1_ps(-30.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
135 /* Avoid stupid compiler warnings */
136 jnrA = jnrB = jnrC = jnrD = 0;
145 for(iidx=0;iidx<4*DIM;iidx++)
150 /* Start outer loop over neighborlists */
151 for(iidx=0; iidx<nri; iidx++)
153 /* Load shift vector for this list */
154 i_shift_offset = DIM*shiftidx[iidx];
156 /* Load limits for loop over neighbors */
157 j_index_start = jindex[iidx];
158 j_index_end = jindex[iidx+1];
160 /* Get outer coordinate index */
162 i_coord_offset = DIM*inr;
164 /* Load i particle coords and add shift vector */
165 gmx_mm_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
167 fix0 = _mm_setzero_ps();
168 fiy0 = _mm_setzero_ps();
169 fiz0 = _mm_setzero_ps();
171 /* Load parameters for i particles */
172 iq0 = _mm_mul_ps(facel,_mm_load1_ps(charge+inr+0));
174 /* Reset potential sums */
175 velecsum = _mm_setzero_ps();
177 /* Start inner kernel loop */
178 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
181 /* Get j neighbor index, and coordinate index */
186 j_coord_offsetA = DIM*jnrA;
187 j_coord_offsetB = DIM*jnrB;
188 j_coord_offsetC = DIM*jnrC;
189 j_coord_offsetD = DIM*jnrD;
191 /* load j atom coordinates */
192 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
193 x+j_coord_offsetC,x+j_coord_offsetD,
196 /* Calculate displacement vector */
197 dx00 = _mm_sub_ps(ix0,jx0);
198 dy00 = _mm_sub_ps(iy0,jy0);
199 dz00 = _mm_sub_ps(iz0,jz0);
201 /* Calculate squared distance and things based on it */
202 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
204 rinv00 = gmx_mm_invsqrt_ps(rsq00);
206 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
208 /* Load parameters for j particles */
209 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
210 charge+jnrC+0,charge+jnrD+0);
212 /**************************
213 * CALCULATE INTERACTIONS *
214 **************************/
216 if (gmx_mm_any_lt(rsq00,rcutoff2))
219 r00 = _mm_mul_ps(rsq00,rinv00);
221 /* Compute parameters for interactions between i and j atoms */
222 qq00 = _mm_mul_ps(iq0,jq0);
224 /* EWALD ELECTROSTATICS */
226 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
227 ewrt = _mm_mul_ps(r00,ewtabscale);
228 ewitab = _mm_cvttps_epi32(ewrt);
229 eweps = _mm_sub_ps(ewrt,_mm_round_ps(ewrt, _MM_FROUND_FLOOR));
230 ewitab = _mm_slli_epi32(ewitab,2);
231 ewtabF = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,0) );
232 ewtabD = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,1) );
233 ewtabV = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,2) );
234 ewtabFn = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,3) );
235 _MM_TRANSPOSE4_PS(ewtabF,ewtabD,ewtabV,ewtabFn);
236 felec = _mm_add_ps(ewtabF,_mm_mul_ps(eweps,ewtabD));
237 velec = _mm_sub_ps(ewtabV,_mm_mul_ps(_mm_mul_ps(ewtabhalfspace,eweps),_mm_add_ps(ewtabF,felec)));
238 velec = _mm_mul_ps(qq00,_mm_sub_ps(rinv00,velec));
239 felec = _mm_mul_ps(_mm_mul_ps(qq00,rinv00),_mm_sub_ps(rinvsq00,felec));
241 d = _mm_sub_ps(r00,rswitch);
242 d = _mm_max_ps(d,_mm_setzero_ps());
243 d2 = _mm_mul_ps(d,d);
244 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)))))));
246 dsw = _mm_mul_ps(d2,_mm_add_ps(swF2,_mm_mul_ps(d,_mm_add_ps(swF3,_mm_mul_ps(d,swF4)))));
248 /* Evaluate switch function */
249 /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
250 felec = _mm_sub_ps( _mm_mul_ps(felec,sw) , _mm_mul_ps(rinv00,_mm_mul_ps(velec,dsw)) );
251 velec = _mm_mul_ps(velec,sw);
252 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
254 /* Update potential sum for this i atom from the interaction with this j atom. */
255 velec = _mm_and_ps(velec,cutoff_mask);
256 velecsum = _mm_add_ps(velecsum,velec);
260 fscal = _mm_and_ps(fscal,cutoff_mask);
262 /* Calculate temporary vectorial force */
263 tx = _mm_mul_ps(fscal,dx00);
264 ty = _mm_mul_ps(fscal,dy00);
265 tz = _mm_mul_ps(fscal,dz00);
267 /* Update vectorial force */
268 fix0 = _mm_add_ps(fix0,tx);
269 fiy0 = _mm_add_ps(fiy0,ty);
270 fiz0 = _mm_add_ps(fiz0,tz);
272 fjptrA = f+j_coord_offsetA;
273 fjptrB = f+j_coord_offsetB;
274 fjptrC = f+j_coord_offsetC;
275 fjptrD = f+j_coord_offsetD;
276 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
280 /* Inner loop uses 65 flops */
286 /* Get j neighbor index, and coordinate index */
287 jnrlistA = jjnr[jidx];
288 jnrlistB = jjnr[jidx+1];
289 jnrlistC = jjnr[jidx+2];
290 jnrlistD = jjnr[jidx+3];
291 /* Sign of each element will be negative for non-real atoms.
292 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
293 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
295 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
296 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
297 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
298 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
299 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
300 j_coord_offsetA = DIM*jnrA;
301 j_coord_offsetB = DIM*jnrB;
302 j_coord_offsetC = DIM*jnrC;
303 j_coord_offsetD = DIM*jnrD;
305 /* load j atom coordinates */
306 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
307 x+j_coord_offsetC,x+j_coord_offsetD,
310 /* Calculate displacement vector */
311 dx00 = _mm_sub_ps(ix0,jx0);
312 dy00 = _mm_sub_ps(iy0,jy0);
313 dz00 = _mm_sub_ps(iz0,jz0);
315 /* Calculate squared distance and things based on it */
316 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
318 rinv00 = gmx_mm_invsqrt_ps(rsq00);
320 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
322 /* Load parameters for j particles */
323 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
324 charge+jnrC+0,charge+jnrD+0);
326 /**************************
327 * CALCULATE INTERACTIONS *
328 **************************/
330 if (gmx_mm_any_lt(rsq00,rcutoff2))
333 r00 = _mm_mul_ps(rsq00,rinv00);
334 r00 = _mm_andnot_ps(dummy_mask,r00);
336 /* Compute parameters for interactions between i and j atoms */
337 qq00 = _mm_mul_ps(iq0,jq0);
339 /* EWALD ELECTROSTATICS */
341 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
342 ewrt = _mm_mul_ps(r00,ewtabscale);
343 ewitab = _mm_cvttps_epi32(ewrt);
344 eweps = _mm_sub_ps(ewrt,_mm_round_ps(ewrt, _MM_FROUND_FLOOR));
345 ewitab = _mm_slli_epi32(ewitab,2);
346 ewtabF = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,0) );
347 ewtabD = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,1) );
348 ewtabV = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,2) );
349 ewtabFn = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,3) );
350 _MM_TRANSPOSE4_PS(ewtabF,ewtabD,ewtabV,ewtabFn);
351 felec = _mm_add_ps(ewtabF,_mm_mul_ps(eweps,ewtabD));
352 velec = _mm_sub_ps(ewtabV,_mm_mul_ps(_mm_mul_ps(ewtabhalfspace,eweps),_mm_add_ps(ewtabF,felec)));
353 velec = _mm_mul_ps(qq00,_mm_sub_ps(rinv00,velec));
354 felec = _mm_mul_ps(_mm_mul_ps(qq00,rinv00),_mm_sub_ps(rinvsq00,felec));
356 d = _mm_sub_ps(r00,rswitch);
357 d = _mm_max_ps(d,_mm_setzero_ps());
358 d2 = _mm_mul_ps(d,d);
359 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)))))));
361 dsw = _mm_mul_ps(d2,_mm_add_ps(swF2,_mm_mul_ps(d,_mm_add_ps(swF3,_mm_mul_ps(d,swF4)))));
363 /* Evaluate switch function */
364 /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
365 felec = _mm_sub_ps( _mm_mul_ps(felec,sw) , _mm_mul_ps(rinv00,_mm_mul_ps(velec,dsw)) );
366 velec = _mm_mul_ps(velec,sw);
367 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
369 /* Update potential sum for this i atom from the interaction with this j atom. */
370 velec = _mm_and_ps(velec,cutoff_mask);
371 velec = _mm_andnot_ps(dummy_mask,velec);
372 velecsum = _mm_add_ps(velecsum,velec);
376 fscal = _mm_and_ps(fscal,cutoff_mask);
378 fscal = _mm_andnot_ps(dummy_mask,fscal);
380 /* Calculate temporary vectorial force */
381 tx = _mm_mul_ps(fscal,dx00);
382 ty = _mm_mul_ps(fscal,dy00);
383 tz = _mm_mul_ps(fscal,dz00);
385 /* Update vectorial force */
386 fix0 = _mm_add_ps(fix0,tx);
387 fiy0 = _mm_add_ps(fiy0,ty);
388 fiz0 = _mm_add_ps(fiz0,tz);
390 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
391 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
392 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
393 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
394 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
398 /* Inner loop uses 66 flops */
401 /* End of innermost loop */
403 gmx_mm_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
404 f+i_coord_offset,fshift+i_shift_offset);
407 /* Update potential energies */
408 gmx_mm_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
410 /* Increment number of inner iterations */
411 inneriter += j_index_end - j_index_start;
413 /* Outer loop uses 8 flops */
416 /* Increment number of outer iterations */
419 /* Update outer/inner flops */
421 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VF,outeriter*8 + inneriter*66);
424 * Gromacs nonbonded kernel: nb_kernel_ElecEwSw_VdwNone_GeomP1P1_F_sse4_1_single
425 * Electrostatics interaction: Ewald
426 * VdW interaction: None
427 * Geometry: Particle-Particle
428 * Calculate force/pot: Force
431 nb_kernel_ElecEwSw_VdwNone_GeomP1P1_F_sse4_1_single
432 (t_nblist * gmx_restrict nlist,
433 rvec * gmx_restrict xx,
434 rvec * gmx_restrict ff,
435 t_forcerec * gmx_restrict fr,
436 t_mdatoms * gmx_restrict mdatoms,
437 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
438 t_nrnb * gmx_restrict nrnb)
440 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
441 * just 0 for non-waters.
442 * Suffixes A,B,C,D refer to j loop unrolling done with SSE, e.g. for the four different
443 * jnr indices corresponding to data put in the four positions in the SIMD register.
445 int i_shift_offset,i_coord_offset,outeriter,inneriter;
446 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
447 int jnrA,jnrB,jnrC,jnrD;
448 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
449 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
450 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
452 real *shiftvec,*fshift,*x,*f;
453 real *fjptrA,*fjptrB,*fjptrC,*fjptrD;
455 __m128 tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
457 __m128 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
458 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
459 __m128 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
460 __m128 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
461 __m128 velec,felec,velecsum,facel,crf,krf,krf2;
464 __m128 ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
466 __m128 rswitch,swV3,swV4,swV5,swF2,swF3,swF4,d,d2,sw,dsw;
467 real rswitch_scalar,d_scalar;
468 __m128 dummy_mask,cutoff_mask;
469 __m128 signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
470 __m128 one = _mm_set1_ps(1.0);
471 __m128 two = _mm_set1_ps(2.0);
477 jindex = nlist->jindex;
479 shiftidx = nlist->shift;
481 shiftvec = fr->shift_vec[0];
482 fshift = fr->fshift[0];
483 facel = _mm_set1_ps(fr->epsfac);
484 charge = mdatoms->chargeA;
486 sh_ewald = _mm_set1_ps(fr->ic->sh_ewald);
487 ewtab = fr->ic->tabq_coul_FDV0;
488 ewtabscale = _mm_set1_ps(fr->ic->tabq_scale);
489 ewtabhalfspace = _mm_set1_ps(0.5/fr->ic->tabq_scale);
491 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
492 rcutoff_scalar = fr->rcoulomb;
493 rcutoff = _mm_set1_ps(rcutoff_scalar);
494 rcutoff2 = _mm_mul_ps(rcutoff,rcutoff);
496 rswitch_scalar = fr->rcoulomb_switch;
497 rswitch = _mm_set1_ps(rswitch_scalar);
498 /* Setup switch parameters */
499 d_scalar = rcutoff_scalar-rswitch_scalar;
500 d = _mm_set1_ps(d_scalar);
501 swV3 = _mm_set1_ps(-10.0/(d_scalar*d_scalar*d_scalar));
502 swV4 = _mm_set1_ps( 15.0/(d_scalar*d_scalar*d_scalar*d_scalar));
503 swV5 = _mm_set1_ps( -6.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
504 swF2 = _mm_set1_ps(-30.0/(d_scalar*d_scalar*d_scalar));
505 swF3 = _mm_set1_ps( 60.0/(d_scalar*d_scalar*d_scalar*d_scalar));
506 swF4 = _mm_set1_ps(-30.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
508 /* Avoid stupid compiler warnings */
509 jnrA = jnrB = jnrC = jnrD = 0;
518 for(iidx=0;iidx<4*DIM;iidx++)
523 /* Start outer loop over neighborlists */
524 for(iidx=0; iidx<nri; iidx++)
526 /* Load shift vector for this list */
527 i_shift_offset = DIM*shiftidx[iidx];
529 /* Load limits for loop over neighbors */
530 j_index_start = jindex[iidx];
531 j_index_end = jindex[iidx+1];
533 /* Get outer coordinate index */
535 i_coord_offset = DIM*inr;
537 /* Load i particle coords and add shift vector */
538 gmx_mm_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
540 fix0 = _mm_setzero_ps();
541 fiy0 = _mm_setzero_ps();
542 fiz0 = _mm_setzero_ps();
544 /* Load parameters for i particles */
545 iq0 = _mm_mul_ps(facel,_mm_load1_ps(charge+inr+0));
547 /* Start inner kernel loop */
548 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
551 /* Get j neighbor index, and coordinate index */
556 j_coord_offsetA = DIM*jnrA;
557 j_coord_offsetB = DIM*jnrB;
558 j_coord_offsetC = DIM*jnrC;
559 j_coord_offsetD = DIM*jnrD;
561 /* load j atom coordinates */
562 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
563 x+j_coord_offsetC,x+j_coord_offsetD,
566 /* Calculate displacement vector */
567 dx00 = _mm_sub_ps(ix0,jx0);
568 dy00 = _mm_sub_ps(iy0,jy0);
569 dz00 = _mm_sub_ps(iz0,jz0);
571 /* Calculate squared distance and things based on it */
572 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
574 rinv00 = gmx_mm_invsqrt_ps(rsq00);
576 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
578 /* Load parameters for j particles */
579 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
580 charge+jnrC+0,charge+jnrD+0);
582 /**************************
583 * CALCULATE INTERACTIONS *
584 **************************/
586 if (gmx_mm_any_lt(rsq00,rcutoff2))
589 r00 = _mm_mul_ps(rsq00,rinv00);
591 /* Compute parameters for interactions between i and j atoms */
592 qq00 = _mm_mul_ps(iq0,jq0);
594 /* EWALD ELECTROSTATICS */
596 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
597 ewrt = _mm_mul_ps(r00,ewtabscale);
598 ewitab = _mm_cvttps_epi32(ewrt);
599 eweps = _mm_sub_ps(ewrt,_mm_round_ps(ewrt, _MM_FROUND_FLOOR));
600 ewitab = _mm_slli_epi32(ewitab,2);
601 ewtabF = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,0) );
602 ewtabD = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,1) );
603 ewtabV = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,2) );
604 ewtabFn = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,3) );
605 _MM_TRANSPOSE4_PS(ewtabF,ewtabD,ewtabV,ewtabFn);
606 felec = _mm_add_ps(ewtabF,_mm_mul_ps(eweps,ewtabD));
607 velec = _mm_sub_ps(ewtabV,_mm_mul_ps(_mm_mul_ps(ewtabhalfspace,eweps),_mm_add_ps(ewtabF,felec)));
608 velec = _mm_mul_ps(qq00,_mm_sub_ps(rinv00,velec));
609 felec = _mm_mul_ps(_mm_mul_ps(qq00,rinv00),_mm_sub_ps(rinvsq00,felec));
611 d = _mm_sub_ps(r00,rswitch);
612 d = _mm_max_ps(d,_mm_setzero_ps());
613 d2 = _mm_mul_ps(d,d);
614 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)))))));
616 dsw = _mm_mul_ps(d2,_mm_add_ps(swF2,_mm_mul_ps(d,_mm_add_ps(swF3,_mm_mul_ps(d,swF4)))));
618 /* Evaluate switch function */
619 /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
620 felec = _mm_sub_ps( _mm_mul_ps(felec,sw) , _mm_mul_ps(rinv00,_mm_mul_ps(velec,dsw)) );
621 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
625 fscal = _mm_and_ps(fscal,cutoff_mask);
627 /* Calculate temporary vectorial force */
628 tx = _mm_mul_ps(fscal,dx00);
629 ty = _mm_mul_ps(fscal,dy00);
630 tz = _mm_mul_ps(fscal,dz00);
632 /* Update vectorial force */
633 fix0 = _mm_add_ps(fix0,tx);
634 fiy0 = _mm_add_ps(fiy0,ty);
635 fiz0 = _mm_add_ps(fiz0,tz);
637 fjptrA = f+j_coord_offsetA;
638 fjptrB = f+j_coord_offsetB;
639 fjptrC = f+j_coord_offsetC;
640 fjptrD = f+j_coord_offsetD;
641 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
645 /* Inner loop uses 62 flops */
651 /* Get j neighbor index, and coordinate index */
652 jnrlistA = jjnr[jidx];
653 jnrlistB = jjnr[jidx+1];
654 jnrlistC = jjnr[jidx+2];
655 jnrlistD = jjnr[jidx+3];
656 /* Sign of each element will be negative for non-real atoms.
657 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
658 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
660 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
661 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
662 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
663 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
664 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
665 j_coord_offsetA = DIM*jnrA;
666 j_coord_offsetB = DIM*jnrB;
667 j_coord_offsetC = DIM*jnrC;
668 j_coord_offsetD = DIM*jnrD;
670 /* load j atom coordinates */
671 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
672 x+j_coord_offsetC,x+j_coord_offsetD,
675 /* Calculate displacement vector */
676 dx00 = _mm_sub_ps(ix0,jx0);
677 dy00 = _mm_sub_ps(iy0,jy0);
678 dz00 = _mm_sub_ps(iz0,jz0);
680 /* Calculate squared distance and things based on it */
681 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
683 rinv00 = gmx_mm_invsqrt_ps(rsq00);
685 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
687 /* Load parameters for j particles */
688 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
689 charge+jnrC+0,charge+jnrD+0);
691 /**************************
692 * CALCULATE INTERACTIONS *
693 **************************/
695 if (gmx_mm_any_lt(rsq00,rcutoff2))
698 r00 = _mm_mul_ps(rsq00,rinv00);
699 r00 = _mm_andnot_ps(dummy_mask,r00);
701 /* Compute parameters for interactions between i and j atoms */
702 qq00 = _mm_mul_ps(iq0,jq0);
704 /* EWALD ELECTROSTATICS */
706 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
707 ewrt = _mm_mul_ps(r00,ewtabscale);
708 ewitab = _mm_cvttps_epi32(ewrt);
709 eweps = _mm_sub_ps(ewrt,_mm_round_ps(ewrt, _MM_FROUND_FLOOR));
710 ewitab = _mm_slli_epi32(ewitab,2);
711 ewtabF = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,0) );
712 ewtabD = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,1) );
713 ewtabV = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,2) );
714 ewtabFn = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,3) );
715 _MM_TRANSPOSE4_PS(ewtabF,ewtabD,ewtabV,ewtabFn);
716 felec = _mm_add_ps(ewtabF,_mm_mul_ps(eweps,ewtabD));
717 velec = _mm_sub_ps(ewtabV,_mm_mul_ps(_mm_mul_ps(ewtabhalfspace,eweps),_mm_add_ps(ewtabF,felec)));
718 velec = _mm_mul_ps(qq00,_mm_sub_ps(rinv00,velec));
719 felec = _mm_mul_ps(_mm_mul_ps(qq00,rinv00),_mm_sub_ps(rinvsq00,felec));
721 d = _mm_sub_ps(r00,rswitch);
722 d = _mm_max_ps(d,_mm_setzero_ps());
723 d2 = _mm_mul_ps(d,d);
724 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)))))));
726 dsw = _mm_mul_ps(d2,_mm_add_ps(swF2,_mm_mul_ps(d,_mm_add_ps(swF3,_mm_mul_ps(d,swF4)))));
728 /* Evaluate switch function */
729 /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
730 felec = _mm_sub_ps( _mm_mul_ps(felec,sw) , _mm_mul_ps(rinv00,_mm_mul_ps(velec,dsw)) );
731 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
735 fscal = _mm_and_ps(fscal,cutoff_mask);
737 fscal = _mm_andnot_ps(dummy_mask,fscal);
739 /* Calculate temporary vectorial force */
740 tx = _mm_mul_ps(fscal,dx00);
741 ty = _mm_mul_ps(fscal,dy00);
742 tz = _mm_mul_ps(fscal,dz00);
744 /* Update vectorial force */
745 fix0 = _mm_add_ps(fix0,tx);
746 fiy0 = _mm_add_ps(fiy0,ty);
747 fiz0 = _mm_add_ps(fiz0,tz);
749 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
750 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
751 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
752 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
753 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
757 /* Inner loop uses 63 flops */
760 /* End of innermost loop */
762 gmx_mm_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
763 f+i_coord_offset,fshift+i_shift_offset);
765 /* Increment number of inner iterations */
766 inneriter += j_index_end - j_index_start;
768 /* Outer loop uses 7 flops */
771 /* Increment number of outer iterations */
774 /* Update outer/inner flops */
776 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_F,outeriter*7 + inneriter*63);