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_ElecCoul_VdwNone_GeomW4P1_VF_sse4_1_single
38 * Electrostatics interaction: Coulomb
39 * VdW interaction: None
40 * Geometry: Water4-Particle
41 * Calculate force/pot: PotentialAndForce
44 nb_kernel_ElecCoul_VdwNone_GeomW4P1_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 ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
72 __m128 ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
74 __m128 ix3,iy3,iz3,fix3,fiy3,fiz3,iq3,isai3;
75 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
76 __m128 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
77 __m128 dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
78 __m128 dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
79 __m128 dx30,dy30,dz30,rsq30,rinv30,rinvsq30,r30,qq30,c6_30,c12_30;
80 __m128 velec,felec,velecsum,facel,crf,krf,krf2;
82 __m128 dummy_mask,cutoff_mask;
83 __m128 signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
84 __m128 one = _mm_set1_ps(1.0);
85 __m128 two = _mm_set1_ps(2.0);
91 jindex = nlist->jindex;
93 shiftidx = nlist->shift;
95 shiftvec = fr->shift_vec[0];
96 fshift = fr->fshift[0];
97 facel = _mm_set1_ps(fr->epsfac);
98 charge = mdatoms->chargeA;
100 /* Setup water-specific parameters */
101 inr = nlist->iinr[0];
102 iq1 = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+1]));
103 iq2 = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+2]));
104 iq3 = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+3]));
106 /* Avoid stupid compiler warnings */
107 jnrA = jnrB = jnrC = jnrD = 0;
116 for(iidx=0;iidx<4*DIM;iidx++)
121 /* Start outer loop over neighborlists */
122 for(iidx=0; iidx<nri; iidx++)
124 /* Load shift vector for this list */
125 i_shift_offset = DIM*shiftidx[iidx];
127 /* Load limits for loop over neighbors */
128 j_index_start = jindex[iidx];
129 j_index_end = jindex[iidx+1];
131 /* Get outer coordinate index */
133 i_coord_offset = DIM*inr;
135 /* Load i particle coords and add shift vector */
136 gmx_mm_load_shift_and_3rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset+DIM,
137 &ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
139 fix1 = _mm_setzero_ps();
140 fiy1 = _mm_setzero_ps();
141 fiz1 = _mm_setzero_ps();
142 fix2 = _mm_setzero_ps();
143 fiy2 = _mm_setzero_ps();
144 fiz2 = _mm_setzero_ps();
145 fix3 = _mm_setzero_ps();
146 fiy3 = _mm_setzero_ps();
147 fiz3 = _mm_setzero_ps();
149 /* Reset potential sums */
150 velecsum = _mm_setzero_ps();
152 /* Start inner kernel loop */
153 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
156 /* Get j neighbor index, and coordinate index */
161 j_coord_offsetA = DIM*jnrA;
162 j_coord_offsetB = DIM*jnrB;
163 j_coord_offsetC = DIM*jnrC;
164 j_coord_offsetD = DIM*jnrD;
166 /* load j atom coordinates */
167 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
168 x+j_coord_offsetC,x+j_coord_offsetD,
171 /* Calculate displacement vector */
172 dx10 = _mm_sub_ps(ix1,jx0);
173 dy10 = _mm_sub_ps(iy1,jy0);
174 dz10 = _mm_sub_ps(iz1,jz0);
175 dx20 = _mm_sub_ps(ix2,jx0);
176 dy20 = _mm_sub_ps(iy2,jy0);
177 dz20 = _mm_sub_ps(iz2,jz0);
178 dx30 = _mm_sub_ps(ix3,jx0);
179 dy30 = _mm_sub_ps(iy3,jy0);
180 dz30 = _mm_sub_ps(iz3,jz0);
182 /* Calculate squared distance and things based on it */
183 rsq10 = gmx_mm_calc_rsq_ps(dx10,dy10,dz10);
184 rsq20 = gmx_mm_calc_rsq_ps(dx20,dy20,dz20);
185 rsq30 = gmx_mm_calc_rsq_ps(dx30,dy30,dz30);
187 rinv10 = gmx_mm_invsqrt_ps(rsq10);
188 rinv20 = gmx_mm_invsqrt_ps(rsq20);
189 rinv30 = gmx_mm_invsqrt_ps(rsq30);
191 rinvsq10 = _mm_mul_ps(rinv10,rinv10);
192 rinvsq20 = _mm_mul_ps(rinv20,rinv20);
193 rinvsq30 = _mm_mul_ps(rinv30,rinv30);
195 /* Load parameters for j particles */
196 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
197 charge+jnrC+0,charge+jnrD+0);
199 /**************************
200 * CALCULATE INTERACTIONS *
201 **************************/
203 /* Compute parameters for interactions between i and j atoms */
204 qq10 = _mm_mul_ps(iq1,jq0);
206 /* COULOMB ELECTROSTATICS */
207 velec = _mm_mul_ps(qq10,rinv10);
208 felec = _mm_mul_ps(velec,rinvsq10);
210 /* Update potential sum for this i atom from the interaction with this j atom. */
211 velecsum = _mm_add_ps(velecsum,velec);
215 /* Calculate temporary vectorial force */
216 tx = _mm_mul_ps(fscal,dx10);
217 ty = _mm_mul_ps(fscal,dy10);
218 tz = _mm_mul_ps(fscal,dz10);
220 /* Update vectorial force */
221 fix1 = _mm_add_ps(fix1,tx);
222 fiy1 = _mm_add_ps(fiy1,ty);
223 fiz1 = _mm_add_ps(fiz1,tz);
225 fjptrA = f+j_coord_offsetA;
226 fjptrB = f+j_coord_offsetB;
227 fjptrC = f+j_coord_offsetC;
228 fjptrD = f+j_coord_offsetD;
229 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
231 /**************************
232 * CALCULATE INTERACTIONS *
233 **************************/
235 /* Compute parameters for interactions between i and j atoms */
236 qq20 = _mm_mul_ps(iq2,jq0);
238 /* COULOMB ELECTROSTATICS */
239 velec = _mm_mul_ps(qq20,rinv20);
240 felec = _mm_mul_ps(velec,rinvsq20);
242 /* Update potential sum for this i atom from the interaction with this j atom. */
243 velecsum = _mm_add_ps(velecsum,velec);
247 /* Calculate temporary vectorial force */
248 tx = _mm_mul_ps(fscal,dx20);
249 ty = _mm_mul_ps(fscal,dy20);
250 tz = _mm_mul_ps(fscal,dz20);
252 /* Update vectorial force */
253 fix2 = _mm_add_ps(fix2,tx);
254 fiy2 = _mm_add_ps(fiy2,ty);
255 fiz2 = _mm_add_ps(fiz2,tz);
257 fjptrA = f+j_coord_offsetA;
258 fjptrB = f+j_coord_offsetB;
259 fjptrC = f+j_coord_offsetC;
260 fjptrD = f+j_coord_offsetD;
261 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
263 /**************************
264 * CALCULATE INTERACTIONS *
265 **************************/
267 /* Compute parameters for interactions between i and j atoms */
268 qq30 = _mm_mul_ps(iq3,jq0);
270 /* COULOMB ELECTROSTATICS */
271 velec = _mm_mul_ps(qq30,rinv30);
272 felec = _mm_mul_ps(velec,rinvsq30);
274 /* Update potential sum for this i atom from the interaction with this j atom. */
275 velecsum = _mm_add_ps(velecsum,velec);
279 /* Calculate temporary vectorial force */
280 tx = _mm_mul_ps(fscal,dx30);
281 ty = _mm_mul_ps(fscal,dy30);
282 tz = _mm_mul_ps(fscal,dz30);
284 /* Update vectorial force */
285 fix3 = _mm_add_ps(fix3,tx);
286 fiy3 = _mm_add_ps(fiy3,ty);
287 fiz3 = _mm_add_ps(fiz3,tz);
289 fjptrA = f+j_coord_offsetA;
290 fjptrB = f+j_coord_offsetB;
291 fjptrC = f+j_coord_offsetC;
292 fjptrD = f+j_coord_offsetD;
293 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
295 /* Inner loop uses 84 flops */
301 /* Get j neighbor index, and coordinate index */
302 jnrlistA = jjnr[jidx];
303 jnrlistB = jjnr[jidx+1];
304 jnrlistC = jjnr[jidx+2];
305 jnrlistD = jjnr[jidx+3];
306 /* Sign of each element will be negative for non-real atoms.
307 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
308 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
310 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
311 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
312 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
313 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
314 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
315 j_coord_offsetA = DIM*jnrA;
316 j_coord_offsetB = DIM*jnrB;
317 j_coord_offsetC = DIM*jnrC;
318 j_coord_offsetD = DIM*jnrD;
320 /* load j atom coordinates */
321 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
322 x+j_coord_offsetC,x+j_coord_offsetD,
325 /* Calculate displacement vector */
326 dx10 = _mm_sub_ps(ix1,jx0);
327 dy10 = _mm_sub_ps(iy1,jy0);
328 dz10 = _mm_sub_ps(iz1,jz0);
329 dx20 = _mm_sub_ps(ix2,jx0);
330 dy20 = _mm_sub_ps(iy2,jy0);
331 dz20 = _mm_sub_ps(iz2,jz0);
332 dx30 = _mm_sub_ps(ix3,jx0);
333 dy30 = _mm_sub_ps(iy3,jy0);
334 dz30 = _mm_sub_ps(iz3,jz0);
336 /* Calculate squared distance and things based on it */
337 rsq10 = gmx_mm_calc_rsq_ps(dx10,dy10,dz10);
338 rsq20 = gmx_mm_calc_rsq_ps(dx20,dy20,dz20);
339 rsq30 = gmx_mm_calc_rsq_ps(dx30,dy30,dz30);
341 rinv10 = gmx_mm_invsqrt_ps(rsq10);
342 rinv20 = gmx_mm_invsqrt_ps(rsq20);
343 rinv30 = gmx_mm_invsqrt_ps(rsq30);
345 rinvsq10 = _mm_mul_ps(rinv10,rinv10);
346 rinvsq20 = _mm_mul_ps(rinv20,rinv20);
347 rinvsq30 = _mm_mul_ps(rinv30,rinv30);
349 /* Load parameters for j particles */
350 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
351 charge+jnrC+0,charge+jnrD+0);
353 /**************************
354 * CALCULATE INTERACTIONS *
355 **************************/
357 /* Compute parameters for interactions between i and j atoms */
358 qq10 = _mm_mul_ps(iq1,jq0);
360 /* COULOMB ELECTROSTATICS */
361 velec = _mm_mul_ps(qq10,rinv10);
362 felec = _mm_mul_ps(velec,rinvsq10);
364 /* Update potential sum for this i atom from the interaction with this j atom. */
365 velec = _mm_andnot_ps(dummy_mask,velec);
366 velecsum = _mm_add_ps(velecsum,velec);
370 fscal = _mm_andnot_ps(dummy_mask,fscal);
372 /* Calculate temporary vectorial force */
373 tx = _mm_mul_ps(fscal,dx10);
374 ty = _mm_mul_ps(fscal,dy10);
375 tz = _mm_mul_ps(fscal,dz10);
377 /* Update vectorial force */
378 fix1 = _mm_add_ps(fix1,tx);
379 fiy1 = _mm_add_ps(fiy1,ty);
380 fiz1 = _mm_add_ps(fiz1,tz);
382 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
383 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
384 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
385 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
386 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
388 /**************************
389 * CALCULATE INTERACTIONS *
390 **************************/
392 /* Compute parameters for interactions between i and j atoms */
393 qq20 = _mm_mul_ps(iq2,jq0);
395 /* COULOMB ELECTROSTATICS */
396 velec = _mm_mul_ps(qq20,rinv20);
397 felec = _mm_mul_ps(velec,rinvsq20);
399 /* Update potential sum for this i atom from the interaction with this j atom. */
400 velec = _mm_andnot_ps(dummy_mask,velec);
401 velecsum = _mm_add_ps(velecsum,velec);
405 fscal = _mm_andnot_ps(dummy_mask,fscal);
407 /* Calculate temporary vectorial force */
408 tx = _mm_mul_ps(fscal,dx20);
409 ty = _mm_mul_ps(fscal,dy20);
410 tz = _mm_mul_ps(fscal,dz20);
412 /* Update vectorial force */
413 fix2 = _mm_add_ps(fix2,tx);
414 fiy2 = _mm_add_ps(fiy2,ty);
415 fiz2 = _mm_add_ps(fiz2,tz);
417 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
418 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
419 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
420 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
421 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
423 /**************************
424 * CALCULATE INTERACTIONS *
425 **************************/
427 /* Compute parameters for interactions between i and j atoms */
428 qq30 = _mm_mul_ps(iq3,jq0);
430 /* COULOMB ELECTROSTATICS */
431 velec = _mm_mul_ps(qq30,rinv30);
432 felec = _mm_mul_ps(velec,rinvsq30);
434 /* Update potential sum for this i atom from the interaction with this j atom. */
435 velec = _mm_andnot_ps(dummy_mask,velec);
436 velecsum = _mm_add_ps(velecsum,velec);
440 fscal = _mm_andnot_ps(dummy_mask,fscal);
442 /* Calculate temporary vectorial force */
443 tx = _mm_mul_ps(fscal,dx30);
444 ty = _mm_mul_ps(fscal,dy30);
445 tz = _mm_mul_ps(fscal,dz30);
447 /* Update vectorial force */
448 fix3 = _mm_add_ps(fix3,tx);
449 fiy3 = _mm_add_ps(fiy3,ty);
450 fiz3 = _mm_add_ps(fiz3,tz);
452 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
453 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
454 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
455 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
456 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
458 /* Inner loop uses 84 flops */
461 /* End of innermost loop */
463 gmx_mm_update_iforce_3atom_swizzle_ps(fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
464 f+i_coord_offset+DIM,fshift+i_shift_offset);
467 /* Update potential energies */
468 gmx_mm_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
470 /* Increment number of inner iterations */
471 inneriter += j_index_end - j_index_start;
473 /* Outer loop uses 19 flops */
476 /* Increment number of outer iterations */
479 /* Update outer/inner flops */
481 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_W4_VF,outeriter*19 + inneriter*84);
484 * Gromacs nonbonded kernel: nb_kernel_ElecCoul_VdwNone_GeomW4P1_F_sse4_1_single
485 * Electrostatics interaction: Coulomb
486 * VdW interaction: None
487 * Geometry: Water4-Particle
488 * Calculate force/pot: Force
491 nb_kernel_ElecCoul_VdwNone_GeomW4P1_F_sse4_1_single
492 (t_nblist * gmx_restrict nlist,
493 rvec * gmx_restrict xx,
494 rvec * gmx_restrict ff,
495 t_forcerec * gmx_restrict fr,
496 t_mdatoms * gmx_restrict mdatoms,
497 nb_kernel_data_t * gmx_restrict kernel_data,
498 t_nrnb * gmx_restrict nrnb)
500 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
501 * just 0 for non-waters.
502 * Suffixes A,B,C,D refer to j loop unrolling done with SSE, e.g. for the four different
503 * jnr indices corresponding to data put in the four positions in the SIMD register.
505 int i_shift_offset,i_coord_offset,outeriter,inneriter;
506 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
507 int jnrA,jnrB,jnrC,jnrD;
508 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
509 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
510 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
512 real *shiftvec,*fshift,*x,*f;
513 real *fjptrA,*fjptrB,*fjptrC,*fjptrD;
515 __m128 tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
517 __m128 ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
519 __m128 ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
521 __m128 ix3,iy3,iz3,fix3,fiy3,fiz3,iq3,isai3;
522 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
523 __m128 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
524 __m128 dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
525 __m128 dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
526 __m128 dx30,dy30,dz30,rsq30,rinv30,rinvsq30,r30,qq30,c6_30,c12_30;
527 __m128 velec,felec,velecsum,facel,crf,krf,krf2;
529 __m128 dummy_mask,cutoff_mask;
530 __m128 signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
531 __m128 one = _mm_set1_ps(1.0);
532 __m128 two = _mm_set1_ps(2.0);
538 jindex = nlist->jindex;
540 shiftidx = nlist->shift;
542 shiftvec = fr->shift_vec[0];
543 fshift = fr->fshift[0];
544 facel = _mm_set1_ps(fr->epsfac);
545 charge = mdatoms->chargeA;
547 /* Setup water-specific parameters */
548 inr = nlist->iinr[0];
549 iq1 = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+1]));
550 iq2 = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+2]));
551 iq3 = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+3]));
553 /* Avoid stupid compiler warnings */
554 jnrA = jnrB = jnrC = jnrD = 0;
563 for(iidx=0;iidx<4*DIM;iidx++)
568 /* Start outer loop over neighborlists */
569 for(iidx=0; iidx<nri; iidx++)
571 /* Load shift vector for this list */
572 i_shift_offset = DIM*shiftidx[iidx];
574 /* Load limits for loop over neighbors */
575 j_index_start = jindex[iidx];
576 j_index_end = jindex[iidx+1];
578 /* Get outer coordinate index */
580 i_coord_offset = DIM*inr;
582 /* Load i particle coords and add shift vector */
583 gmx_mm_load_shift_and_3rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset+DIM,
584 &ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
586 fix1 = _mm_setzero_ps();
587 fiy1 = _mm_setzero_ps();
588 fiz1 = _mm_setzero_ps();
589 fix2 = _mm_setzero_ps();
590 fiy2 = _mm_setzero_ps();
591 fiz2 = _mm_setzero_ps();
592 fix3 = _mm_setzero_ps();
593 fiy3 = _mm_setzero_ps();
594 fiz3 = _mm_setzero_ps();
596 /* Start inner kernel loop */
597 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
600 /* Get j neighbor index, and coordinate index */
605 j_coord_offsetA = DIM*jnrA;
606 j_coord_offsetB = DIM*jnrB;
607 j_coord_offsetC = DIM*jnrC;
608 j_coord_offsetD = DIM*jnrD;
610 /* load j atom coordinates */
611 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
612 x+j_coord_offsetC,x+j_coord_offsetD,
615 /* Calculate displacement vector */
616 dx10 = _mm_sub_ps(ix1,jx0);
617 dy10 = _mm_sub_ps(iy1,jy0);
618 dz10 = _mm_sub_ps(iz1,jz0);
619 dx20 = _mm_sub_ps(ix2,jx0);
620 dy20 = _mm_sub_ps(iy2,jy0);
621 dz20 = _mm_sub_ps(iz2,jz0);
622 dx30 = _mm_sub_ps(ix3,jx0);
623 dy30 = _mm_sub_ps(iy3,jy0);
624 dz30 = _mm_sub_ps(iz3,jz0);
626 /* Calculate squared distance and things based on it */
627 rsq10 = gmx_mm_calc_rsq_ps(dx10,dy10,dz10);
628 rsq20 = gmx_mm_calc_rsq_ps(dx20,dy20,dz20);
629 rsq30 = gmx_mm_calc_rsq_ps(dx30,dy30,dz30);
631 rinv10 = gmx_mm_invsqrt_ps(rsq10);
632 rinv20 = gmx_mm_invsqrt_ps(rsq20);
633 rinv30 = gmx_mm_invsqrt_ps(rsq30);
635 rinvsq10 = _mm_mul_ps(rinv10,rinv10);
636 rinvsq20 = _mm_mul_ps(rinv20,rinv20);
637 rinvsq30 = _mm_mul_ps(rinv30,rinv30);
639 /* Load parameters for j particles */
640 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
641 charge+jnrC+0,charge+jnrD+0);
643 /**************************
644 * CALCULATE INTERACTIONS *
645 **************************/
647 /* Compute parameters for interactions between i and j atoms */
648 qq10 = _mm_mul_ps(iq1,jq0);
650 /* COULOMB ELECTROSTATICS */
651 velec = _mm_mul_ps(qq10,rinv10);
652 felec = _mm_mul_ps(velec,rinvsq10);
656 /* Calculate temporary vectorial force */
657 tx = _mm_mul_ps(fscal,dx10);
658 ty = _mm_mul_ps(fscal,dy10);
659 tz = _mm_mul_ps(fscal,dz10);
661 /* Update vectorial force */
662 fix1 = _mm_add_ps(fix1,tx);
663 fiy1 = _mm_add_ps(fiy1,ty);
664 fiz1 = _mm_add_ps(fiz1,tz);
666 fjptrA = f+j_coord_offsetA;
667 fjptrB = f+j_coord_offsetB;
668 fjptrC = f+j_coord_offsetC;
669 fjptrD = f+j_coord_offsetD;
670 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
672 /**************************
673 * CALCULATE INTERACTIONS *
674 **************************/
676 /* Compute parameters for interactions between i and j atoms */
677 qq20 = _mm_mul_ps(iq2,jq0);
679 /* COULOMB ELECTROSTATICS */
680 velec = _mm_mul_ps(qq20,rinv20);
681 felec = _mm_mul_ps(velec,rinvsq20);
685 /* Calculate temporary vectorial force */
686 tx = _mm_mul_ps(fscal,dx20);
687 ty = _mm_mul_ps(fscal,dy20);
688 tz = _mm_mul_ps(fscal,dz20);
690 /* Update vectorial force */
691 fix2 = _mm_add_ps(fix2,tx);
692 fiy2 = _mm_add_ps(fiy2,ty);
693 fiz2 = _mm_add_ps(fiz2,tz);
695 fjptrA = f+j_coord_offsetA;
696 fjptrB = f+j_coord_offsetB;
697 fjptrC = f+j_coord_offsetC;
698 fjptrD = f+j_coord_offsetD;
699 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
701 /**************************
702 * CALCULATE INTERACTIONS *
703 **************************/
705 /* Compute parameters for interactions between i and j atoms */
706 qq30 = _mm_mul_ps(iq3,jq0);
708 /* COULOMB ELECTROSTATICS */
709 velec = _mm_mul_ps(qq30,rinv30);
710 felec = _mm_mul_ps(velec,rinvsq30);
714 /* Calculate temporary vectorial force */
715 tx = _mm_mul_ps(fscal,dx30);
716 ty = _mm_mul_ps(fscal,dy30);
717 tz = _mm_mul_ps(fscal,dz30);
719 /* Update vectorial force */
720 fix3 = _mm_add_ps(fix3,tx);
721 fiy3 = _mm_add_ps(fiy3,ty);
722 fiz3 = _mm_add_ps(fiz3,tz);
724 fjptrA = f+j_coord_offsetA;
725 fjptrB = f+j_coord_offsetB;
726 fjptrC = f+j_coord_offsetC;
727 fjptrD = f+j_coord_offsetD;
728 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
730 /* Inner loop uses 81 flops */
736 /* Get j neighbor index, and coordinate index */
737 jnrlistA = jjnr[jidx];
738 jnrlistB = jjnr[jidx+1];
739 jnrlistC = jjnr[jidx+2];
740 jnrlistD = jjnr[jidx+3];
741 /* Sign of each element will be negative for non-real atoms.
742 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
743 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
745 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
746 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
747 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
748 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
749 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
750 j_coord_offsetA = DIM*jnrA;
751 j_coord_offsetB = DIM*jnrB;
752 j_coord_offsetC = DIM*jnrC;
753 j_coord_offsetD = DIM*jnrD;
755 /* load j atom coordinates */
756 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
757 x+j_coord_offsetC,x+j_coord_offsetD,
760 /* Calculate displacement vector */
761 dx10 = _mm_sub_ps(ix1,jx0);
762 dy10 = _mm_sub_ps(iy1,jy0);
763 dz10 = _mm_sub_ps(iz1,jz0);
764 dx20 = _mm_sub_ps(ix2,jx0);
765 dy20 = _mm_sub_ps(iy2,jy0);
766 dz20 = _mm_sub_ps(iz2,jz0);
767 dx30 = _mm_sub_ps(ix3,jx0);
768 dy30 = _mm_sub_ps(iy3,jy0);
769 dz30 = _mm_sub_ps(iz3,jz0);
771 /* Calculate squared distance and things based on it */
772 rsq10 = gmx_mm_calc_rsq_ps(dx10,dy10,dz10);
773 rsq20 = gmx_mm_calc_rsq_ps(dx20,dy20,dz20);
774 rsq30 = gmx_mm_calc_rsq_ps(dx30,dy30,dz30);
776 rinv10 = gmx_mm_invsqrt_ps(rsq10);
777 rinv20 = gmx_mm_invsqrt_ps(rsq20);
778 rinv30 = gmx_mm_invsqrt_ps(rsq30);
780 rinvsq10 = _mm_mul_ps(rinv10,rinv10);
781 rinvsq20 = _mm_mul_ps(rinv20,rinv20);
782 rinvsq30 = _mm_mul_ps(rinv30,rinv30);
784 /* Load parameters for j particles */
785 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
786 charge+jnrC+0,charge+jnrD+0);
788 /**************************
789 * CALCULATE INTERACTIONS *
790 **************************/
792 /* Compute parameters for interactions between i and j atoms */
793 qq10 = _mm_mul_ps(iq1,jq0);
795 /* COULOMB ELECTROSTATICS */
796 velec = _mm_mul_ps(qq10,rinv10);
797 felec = _mm_mul_ps(velec,rinvsq10);
801 fscal = _mm_andnot_ps(dummy_mask,fscal);
803 /* Calculate temporary vectorial force */
804 tx = _mm_mul_ps(fscal,dx10);
805 ty = _mm_mul_ps(fscal,dy10);
806 tz = _mm_mul_ps(fscal,dz10);
808 /* Update vectorial force */
809 fix1 = _mm_add_ps(fix1,tx);
810 fiy1 = _mm_add_ps(fiy1,ty);
811 fiz1 = _mm_add_ps(fiz1,tz);
813 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
814 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
815 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
816 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
817 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
819 /**************************
820 * CALCULATE INTERACTIONS *
821 **************************/
823 /* Compute parameters for interactions between i and j atoms */
824 qq20 = _mm_mul_ps(iq2,jq0);
826 /* COULOMB ELECTROSTATICS */
827 velec = _mm_mul_ps(qq20,rinv20);
828 felec = _mm_mul_ps(velec,rinvsq20);
832 fscal = _mm_andnot_ps(dummy_mask,fscal);
834 /* Calculate temporary vectorial force */
835 tx = _mm_mul_ps(fscal,dx20);
836 ty = _mm_mul_ps(fscal,dy20);
837 tz = _mm_mul_ps(fscal,dz20);
839 /* Update vectorial force */
840 fix2 = _mm_add_ps(fix2,tx);
841 fiy2 = _mm_add_ps(fiy2,ty);
842 fiz2 = _mm_add_ps(fiz2,tz);
844 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
845 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
846 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
847 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
848 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
850 /**************************
851 * CALCULATE INTERACTIONS *
852 **************************/
854 /* Compute parameters for interactions between i and j atoms */
855 qq30 = _mm_mul_ps(iq3,jq0);
857 /* COULOMB ELECTROSTATICS */
858 velec = _mm_mul_ps(qq30,rinv30);
859 felec = _mm_mul_ps(velec,rinvsq30);
863 fscal = _mm_andnot_ps(dummy_mask,fscal);
865 /* Calculate temporary vectorial force */
866 tx = _mm_mul_ps(fscal,dx30);
867 ty = _mm_mul_ps(fscal,dy30);
868 tz = _mm_mul_ps(fscal,dz30);
870 /* Update vectorial force */
871 fix3 = _mm_add_ps(fix3,tx);
872 fiy3 = _mm_add_ps(fiy3,ty);
873 fiz3 = _mm_add_ps(fiz3,tz);
875 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
876 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
877 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
878 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
879 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
881 /* Inner loop uses 81 flops */
884 /* End of innermost loop */
886 gmx_mm_update_iforce_3atom_swizzle_ps(fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
887 f+i_coord_offset+DIM,fshift+i_shift_offset);
889 /* Increment number of inner iterations */
890 inneriter += j_index_end - j_index_start;
892 /* Outer loop uses 18 flops */
895 /* Increment number of outer iterations */
898 /* Update outer/inner flops */
900 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_W4_F,outeriter*18 + inneriter*81);