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 "types/simple.h"
44 #include "gromacs/math/vec.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_ElecRFCut_VdwNone_GeomP1P1_VF_sse4_1_single
52 * Electrostatics interaction: ReactionField
53 * VdW interaction: None
54 * Geometry: Particle-Particle
55 * Calculate force/pot: PotentialAndForce
58 nb_kernel_ElecRFCut_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;
90 __m128 dummy_mask,cutoff_mask;
91 __m128 signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
92 __m128 one = _mm_set1_ps(1.0);
93 __m128 two = _mm_set1_ps(2.0);
99 jindex = nlist->jindex;
101 shiftidx = nlist->shift;
103 shiftvec = fr->shift_vec[0];
104 fshift = fr->fshift[0];
105 facel = _mm_set1_ps(fr->epsfac);
106 charge = mdatoms->chargeA;
107 krf = _mm_set1_ps(fr->ic->k_rf);
108 krf2 = _mm_set1_ps(fr->ic->k_rf*2.0);
109 crf = _mm_set1_ps(fr->ic->c_rf);
111 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
112 rcutoff_scalar = fr->rcoulomb;
113 rcutoff = _mm_set1_ps(rcutoff_scalar);
114 rcutoff2 = _mm_mul_ps(rcutoff,rcutoff);
116 /* Avoid stupid compiler warnings */
117 jnrA = jnrB = jnrC = jnrD = 0;
126 for(iidx=0;iidx<4*DIM;iidx++)
131 /* Start outer loop over neighborlists */
132 for(iidx=0; iidx<nri; iidx++)
134 /* Load shift vector for this list */
135 i_shift_offset = DIM*shiftidx[iidx];
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 gmx_mm_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
148 fix0 = _mm_setzero_ps();
149 fiy0 = _mm_setzero_ps();
150 fiz0 = _mm_setzero_ps();
152 /* Load parameters for i particles */
153 iq0 = _mm_mul_ps(facel,_mm_load1_ps(charge+inr+0));
155 /* Reset potential sums */
156 velecsum = _mm_setzero_ps();
158 /* Start inner kernel loop */
159 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
162 /* Get j neighbor index, and coordinate index */
167 j_coord_offsetA = DIM*jnrA;
168 j_coord_offsetB = DIM*jnrB;
169 j_coord_offsetC = DIM*jnrC;
170 j_coord_offsetD = DIM*jnrD;
172 /* load j atom coordinates */
173 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
174 x+j_coord_offsetC,x+j_coord_offsetD,
177 /* Calculate displacement vector */
178 dx00 = _mm_sub_ps(ix0,jx0);
179 dy00 = _mm_sub_ps(iy0,jy0);
180 dz00 = _mm_sub_ps(iz0,jz0);
182 /* Calculate squared distance and things based on it */
183 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
185 rinv00 = gmx_mm_invsqrt_ps(rsq00);
187 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
189 /* Load parameters for j particles */
190 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
191 charge+jnrC+0,charge+jnrD+0);
193 /**************************
194 * CALCULATE INTERACTIONS *
195 **************************/
197 if (gmx_mm_any_lt(rsq00,rcutoff2))
200 /* Compute parameters for interactions between i and j atoms */
201 qq00 = _mm_mul_ps(iq0,jq0);
203 /* REACTION-FIELD ELECTROSTATICS */
204 velec = _mm_mul_ps(qq00,_mm_sub_ps(_mm_add_ps(rinv00,_mm_mul_ps(krf,rsq00)),crf));
205 felec = _mm_mul_ps(qq00,_mm_sub_ps(_mm_mul_ps(rinv00,rinvsq00),krf2));
207 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
209 /* Update potential sum for this i atom from the interaction with this j atom. */
210 velec = _mm_and_ps(velec,cutoff_mask);
211 velecsum = _mm_add_ps(velecsum,velec);
215 fscal = _mm_and_ps(fscal,cutoff_mask);
217 /* Calculate temporary vectorial force */
218 tx = _mm_mul_ps(fscal,dx00);
219 ty = _mm_mul_ps(fscal,dy00);
220 tz = _mm_mul_ps(fscal,dz00);
222 /* Update vectorial force */
223 fix0 = _mm_add_ps(fix0,tx);
224 fiy0 = _mm_add_ps(fiy0,ty);
225 fiz0 = _mm_add_ps(fiz0,tz);
227 fjptrA = f+j_coord_offsetA;
228 fjptrB = f+j_coord_offsetB;
229 fjptrC = f+j_coord_offsetC;
230 fjptrD = f+j_coord_offsetD;
231 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
235 /* Inner loop uses 36 flops */
241 /* Get j neighbor index, and coordinate index */
242 jnrlistA = jjnr[jidx];
243 jnrlistB = jjnr[jidx+1];
244 jnrlistC = jjnr[jidx+2];
245 jnrlistD = jjnr[jidx+3];
246 /* Sign of each element will be negative for non-real atoms.
247 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
248 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
250 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
251 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
252 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
253 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
254 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
255 j_coord_offsetA = DIM*jnrA;
256 j_coord_offsetB = DIM*jnrB;
257 j_coord_offsetC = DIM*jnrC;
258 j_coord_offsetD = DIM*jnrD;
260 /* load j atom coordinates */
261 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
262 x+j_coord_offsetC,x+j_coord_offsetD,
265 /* Calculate displacement vector */
266 dx00 = _mm_sub_ps(ix0,jx0);
267 dy00 = _mm_sub_ps(iy0,jy0);
268 dz00 = _mm_sub_ps(iz0,jz0);
270 /* Calculate squared distance and things based on it */
271 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
273 rinv00 = gmx_mm_invsqrt_ps(rsq00);
275 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
277 /* Load parameters for j particles */
278 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
279 charge+jnrC+0,charge+jnrD+0);
281 /**************************
282 * CALCULATE INTERACTIONS *
283 **************************/
285 if (gmx_mm_any_lt(rsq00,rcutoff2))
288 /* Compute parameters for interactions between i and j atoms */
289 qq00 = _mm_mul_ps(iq0,jq0);
291 /* REACTION-FIELD ELECTROSTATICS */
292 velec = _mm_mul_ps(qq00,_mm_sub_ps(_mm_add_ps(rinv00,_mm_mul_ps(krf,rsq00)),crf));
293 felec = _mm_mul_ps(qq00,_mm_sub_ps(_mm_mul_ps(rinv00,rinvsq00),krf2));
295 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
297 /* Update potential sum for this i atom from the interaction with this j atom. */
298 velec = _mm_and_ps(velec,cutoff_mask);
299 velec = _mm_andnot_ps(dummy_mask,velec);
300 velecsum = _mm_add_ps(velecsum,velec);
304 fscal = _mm_and_ps(fscal,cutoff_mask);
306 fscal = _mm_andnot_ps(dummy_mask,fscal);
308 /* Calculate temporary vectorial force */
309 tx = _mm_mul_ps(fscal,dx00);
310 ty = _mm_mul_ps(fscal,dy00);
311 tz = _mm_mul_ps(fscal,dz00);
313 /* Update vectorial force */
314 fix0 = _mm_add_ps(fix0,tx);
315 fiy0 = _mm_add_ps(fiy0,ty);
316 fiz0 = _mm_add_ps(fiz0,tz);
318 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
319 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
320 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
321 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
322 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
326 /* Inner loop uses 36 flops */
329 /* End of innermost loop */
331 gmx_mm_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
332 f+i_coord_offset,fshift+i_shift_offset);
335 /* Update potential energies */
336 gmx_mm_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
338 /* Increment number of inner iterations */
339 inneriter += j_index_end - j_index_start;
341 /* Outer loop uses 8 flops */
344 /* Increment number of outer iterations */
347 /* Update outer/inner flops */
349 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VF,outeriter*8 + inneriter*36);
352 * Gromacs nonbonded kernel: nb_kernel_ElecRFCut_VdwNone_GeomP1P1_F_sse4_1_single
353 * Electrostatics interaction: ReactionField
354 * VdW interaction: None
355 * Geometry: Particle-Particle
356 * Calculate force/pot: Force
359 nb_kernel_ElecRFCut_VdwNone_GeomP1P1_F_sse4_1_single
360 (t_nblist * gmx_restrict nlist,
361 rvec * gmx_restrict xx,
362 rvec * gmx_restrict ff,
363 t_forcerec * gmx_restrict fr,
364 t_mdatoms * gmx_restrict mdatoms,
365 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
366 t_nrnb * gmx_restrict nrnb)
368 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
369 * just 0 for non-waters.
370 * Suffixes A,B,C,D refer to j loop unrolling done with SSE, e.g. for the four different
371 * jnr indices corresponding to data put in the four positions in the SIMD register.
373 int i_shift_offset,i_coord_offset,outeriter,inneriter;
374 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
375 int jnrA,jnrB,jnrC,jnrD;
376 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
377 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
378 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
380 real *shiftvec,*fshift,*x,*f;
381 real *fjptrA,*fjptrB,*fjptrC,*fjptrD;
383 __m128 tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
385 __m128 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
386 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
387 __m128 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
388 __m128 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
389 __m128 velec,felec,velecsum,facel,crf,krf,krf2;
391 __m128 dummy_mask,cutoff_mask;
392 __m128 signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
393 __m128 one = _mm_set1_ps(1.0);
394 __m128 two = _mm_set1_ps(2.0);
400 jindex = nlist->jindex;
402 shiftidx = nlist->shift;
404 shiftvec = fr->shift_vec[0];
405 fshift = fr->fshift[0];
406 facel = _mm_set1_ps(fr->epsfac);
407 charge = mdatoms->chargeA;
408 krf = _mm_set1_ps(fr->ic->k_rf);
409 krf2 = _mm_set1_ps(fr->ic->k_rf*2.0);
410 crf = _mm_set1_ps(fr->ic->c_rf);
412 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
413 rcutoff_scalar = fr->rcoulomb;
414 rcutoff = _mm_set1_ps(rcutoff_scalar);
415 rcutoff2 = _mm_mul_ps(rcutoff,rcutoff);
417 /* Avoid stupid compiler warnings */
418 jnrA = jnrB = jnrC = jnrD = 0;
427 for(iidx=0;iidx<4*DIM;iidx++)
432 /* Start outer loop over neighborlists */
433 for(iidx=0; iidx<nri; iidx++)
435 /* Load shift vector for this list */
436 i_shift_offset = DIM*shiftidx[iidx];
438 /* Load limits for loop over neighbors */
439 j_index_start = jindex[iidx];
440 j_index_end = jindex[iidx+1];
442 /* Get outer coordinate index */
444 i_coord_offset = DIM*inr;
446 /* Load i particle coords and add shift vector */
447 gmx_mm_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
449 fix0 = _mm_setzero_ps();
450 fiy0 = _mm_setzero_ps();
451 fiz0 = _mm_setzero_ps();
453 /* Load parameters for i particles */
454 iq0 = _mm_mul_ps(facel,_mm_load1_ps(charge+inr+0));
456 /* Start inner kernel loop */
457 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
460 /* Get j neighbor index, and coordinate index */
465 j_coord_offsetA = DIM*jnrA;
466 j_coord_offsetB = DIM*jnrB;
467 j_coord_offsetC = DIM*jnrC;
468 j_coord_offsetD = DIM*jnrD;
470 /* load j atom coordinates */
471 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
472 x+j_coord_offsetC,x+j_coord_offsetD,
475 /* Calculate displacement vector */
476 dx00 = _mm_sub_ps(ix0,jx0);
477 dy00 = _mm_sub_ps(iy0,jy0);
478 dz00 = _mm_sub_ps(iz0,jz0);
480 /* Calculate squared distance and things based on it */
481 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
483 rinv00 = gmx_mm_invsqrt_ps(rsq00);
485 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
487 /* Load parameters for j particles */
488 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
489 charge+jnrC+0,charge+jnrD+0);
491 /**************************
492 * CALCULATE INTERACTIONS *
493 **************************/
495 if (gmx_mm_any_lt(rsq00,rcutoff2))
498 /* Compute parameters for interactions between i and j atoms */
499 qq00 = _mm_mul_ps(iq0,jq0);
501 /* REACTION-FIELD ELECTROSTATICS */
502 felec = _mm_mul_ps(qq00,_mm_sub_ps(_mm_mul_ps(rinv00,rinvsq00),krf2));
504 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
508 fscal = _mm_and_ps(fscal,cutoff_mask);
510 /* Calculate temporary vectorial force */
511 tx = _mm_mul_ps(fscal,dx00);
512 ty = _mm_mul_ps(fscal,dy00);
513 tz = _mm_mul_ps(fscal,dz00);
515 /* Update vectorial force */
516 fix0 = _mm_add_ps(fix0,tx);
517 fiy0 = _mm_add_ps(fiy0,ty);
518 fiz0 = _mm_add_ps(fiz0,tz);
520 fjptrA = f+j_coord_offsetA;
521 fjptrB = f+j_coord_offsetB;
522 fjptrC = f+j_coord_offsetC;
523 fjptrD = f+j_coord_offsetD;
524 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
528 /* Inner loop uses 30 flops */
534 /* Get j neighbor index, and coordinate index */
535 jnrlistA = jjnr[jidx];
536 jnrlistB = jjnr[jidx+1];
537 jnrlistC = jjnr[jidx+2];
538 jnrlistD = jjnr[jidx+3];
539 /* Sign of each element will be negative for non-real atoms.
540 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
541 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
543 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
544 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
545 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
546 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
547 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
548 j_coord_offsetA = DIM*jnrA;
549 j_coord_offsetB = DIM*jnrB;
550 j_coord_offsetC = DIM*jnrC;
551 j_coord_offsetD = DIM*jnrD;
553 /* load j atom coordinates */
554 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
555 x+j_coord_offsetC,x+j_coord_offsetD,
558 /* Calculate displacement vector */
559 dx00 = _mm_sub_ps(ix0,jx0);
560 dy00 = _mm_sub_ps(iy0,jy0);
561 dz00 = _mm_sub_ps(iz0,jz0);
563 /* Calculate squared distance and things based on it */
564 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
566 rinv00 = gmx_mm_invsqrt_ps(rsq00);
568 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
570 /* Load parameters for j particles */
571 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
572 charge+jnrC+0,charge+jnrD+0);
574 /**************************
575 * CALCULATE INTERACTIONS *
576 **************************/
578 if (gmx_mm_any_lt(rsq00,rcutoff2))
581 /* Compute parameters for interactions between i and j atoms */
582 qq00 = _mm_mul_ps(iq0,jq0);
584 /* REACTION-FIELD ELECTROSTATICS */
585 felec = _mm_mul_ps(qq00,_mm_sub_ps(_mm_mul_ps(rinv00,rinvsq00),krf2));
587 cutoff_mask = _mm_cmplt_ps(rsq00,rcutoff2);
591 fscal = _mm_and_ps(fscal,cutoff_mask);
593 fscal = _mm_andnot_ps(dummy_mask,fscal);
595 /* Calculate temporary vectorial force */
596 tx = _mm_mul_ps(fscal,dx00);
597 ty = _mm_mul_ps(fscal,dy00);
598 tz = _mm_mul_ps(fscal,dz00);
600 /* Update vectorial force */
601 fix0 = _mm_add_ps(fix0,tx);
602 fiy0 = _mm_add_ps(fiy0,ty);
603 fiz0 = _mm_add_ps(fiz0,tz);
605 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
606 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
607 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
608 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
609 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
613 /* Inner loop uses 30 flops */
616 /* End of innermost loop */
618 gmx_mm_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
619 f+i_coord_offset,fshift+i_shift_offset);
621 /* Increment number of inner iterations */
622 inneriter += j_index_end - j_index_start;
624 /* Outer loop uses 7 flops */
627 /* Increment number of outer iterations */
630 /* Update outer/inner flops */
632 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_F,outeriter*7 + inneriter*30);