Removed simple.h from nb_kernel_sse2_XX
[gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_sse2_double / nb_kernel_ElecCoul_VdwLJ_GeomP1P1_sse2_double.c
blob245a0c584308abb386e444f1118329a5ae4dde48
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015, 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 sse2_double kernel generator.
38 #include "gmxpre.h"
40 #include "config.h"
42 #include <math.h>
44 #include "../nb_kernel.h"
45 #include "gromacs/math/vec.h"
46 #include "gromacs/legacyheaders/nrnb.h"
48 #include "gromacs/simd/math_x86_sse2_double.h"
49 #include "kernelutil_x86_sse2_double.h"
52 * Gromacs nonbonded kernel: nb_kernel_ElecCoul_VdwLJ_GeomP1P1_VF_sse2_double
53 * Electrostatics interaction: Coulomb
54 * VdW interaction: LennardJones
55 * Geometry: Particle-Particle
56 * Calculate force/pot: PotentialAndForce
58 void
59 nb_kernel_ElecCoul_VdwLJ_GeomP1P1_VF_sse2_double
60 (t_nblist * gmx_restrict nlist,
61 rvec * gmx_restrict xx,
62 rvec * gmx_restrict ff,
63 t_forcerec * gmx_restrict fr,
64 t_mdatoms * gmx_restrict mdatoms,
65 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
66 t_nrnb * gmx_restrict nrnb)
68 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
69 * just 0 for non-waters.
70 * Suffixes A,B refer to j loop unrolling done with SSE double precision, e.g. for the two different
71 * jnr indices corresponding to data put in the four positions in the SIMD register.
73 int i_shift_offset,i_coord_offset,outeriter,inneriter;
74 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
75 int jnrA,jnrB;
76 int j_coord_offsetA,j_coord_offsetB;
77 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
78 real rcutoff_scalar;
79 real *shiftvec,*fshift,*x,*f;
80 __m128d tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
81 int vdwioffset0;
82 __m128d ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
83 int vdwjidx0A,vdwjidx0B;
84 __m128d jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
85 __m128d dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
86 __m128d velec,felec,velecsum,facel,crf,krf,krf2;
87 real *charge;
88 int nvdwtype;
89 __m128d rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
90 int *vdwtype;
91 real *vdwparam;
92 __m128d one_sixth = _mm_set1_pd(1.0/6.0);
93 __m128d one_twelfth = _mm_set1_pd(1.0/12.0);
94 __m128d dummy_mask,cutoff_mask;
95 __m128d signbit = gmx_mm_castsi128_pd( _mm_set_epi32(0x80000000,0x00000000,0x80000000,0x00000000) );
96 __m128d one = _mm_set1_pd(1.0);
97 __m128d two = _mm_set1_pd(2.0);
98 x = xx[0];
99 f = ff[0];
101 nri = nlist->nri;
102 iinr = nlist->iinr;
103 jindex = nlist->jindex;
104 jjnr = nlist->jjnr;
105 shiftidx = nlist->shift;
106 gid = nlist->gid;
107 shiftvec = fr->shift_vec[0];
108 fshift = fr->fshift[0];
109 facel = _mm_set1_pd(fr->epsfac);
110 charge = mdatoms->chargeA;
111 nvdwtype = fr->ntype;
112 vdwparam = fr->nbfp;
113 vdwtype = mdatoms->typeA;
115 /* Avoid stupid compiler warnings */
116 jnrA = jnrB = 0;
117 j_coord_offsetA = 0;
118 j_coord_offsetB = 0;
120 outeriter = 0;
121 inneriter = 0;
123 /* Start outer loop over neighborlists */
124 for(iidx=0; iidx<nri; iidx++)
126 /* Load shift vector for this list */
127 i_shift_offset = DIM*shiftidx[iidx];
129 /* Load limits for loop over neighbors */
130 j_index_start = jindex[iidx];
131 j_index_end = jindex[iidx+1];
133 /* Get outer coordinate index */
134 inr = iinr[iidx];
135 i_coord_offset = DIM*inr;
137 /* Load i particle coords and add shift vector */
138 gmx_mm_load_shift_and_1rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
140 fix0 = _mm_setzero_pd();
141 fiy0 = _mm_setzero_pd();
142 fiz0 = _mm_setzero_pd();
144 /* Load parameters for i particles */
145 iq0 = _mm_mul_pd(facel,_mm_load1_pd(charge+inr+0));
146 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
148 /* Reset potential sums */
149 velecsum = _mm_setzero_pd();
150 vvdwsum = _mm_setzero_pd();
152 /* Start inner kernel loop */
153 for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
156 /* Get j neighbor index, and coordinate index */
157 jnrA = jjnr[jidx];
158 jnrB = jjnr[jidx+1];
159 j_coord_offsetA = DIM*jnrA;
160 j_coord_offsetB = DIM*jnrB;
162 /* load j atom coordinates */
163 gmx_mm_load_1rvec_2ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
164 &jx0,&jy0,&jz0);
166 /* Calculate displacement vector */
167 dx00 = _mm_sub_pd(ix0,jx0);
168 dy00 = _mm_sub_pd(iy0,jy0);
169 dz00 = _mm_sub_pd(iz0,jz0);
171 /* Calculate squared distance and things based on it */
172 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
174 rinv00 = gmx_mm_invsqrt_pd(rsq00);
176 rinvsq00 = _mm_mul_pd(rinv00,rinv00);
178 /* Load parameters for j particles */
179 jq0 = gmx_mm_load_2real_swizzle_pd(charge+jnrA+0,charge+jnrB+0);
180 vdwjidx0A = 2*vdwtype[jnrA+0];
181 vdwjidx0B = 2*vdwtype[jnrB+0];
183 /**************************
184 * CALCULATE INTERACTIONS *
185 **************************/
187 /* Compute parameters for interactions between i and j atoms */
188 qq00 = _mm_mul_pd(iq0,jq0);
189 gmx_mm_load_2pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,
190 vdwparam+vdwioffset0+vdwjidx0B,&c6_00,&c12_00);
192 /* COULOMB ELECTROSTATICS */
193 velec = _mm_mul_pd(qq00,rinv00);
194 felec = _mm_mul_pd(velec,rinvsq00);
196 /* LENNARD-JONES DISPERSION/REPULSION */
198 rinvsix = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
199 vvdw6 = _mm_mul_pd(c6_00,rinvsix);
200 vvdw12 = _mm_mul_pd(c12_00,_mm_mul_pd(rinvsix,rinvsix));
201 vvdw = _mm_sub_pd( _mm_mul_pd(vvdw12,one_twelfth) , _mm_mul_pd(vvdw6,one_sixth) );
202 fvdw = _mm_mul_pd(_mm_sub_pd(vvdw12,vvdw6),rinvsq00);
204 /* Update potential sum for this i atom from the interaction with this j atom. */
205 velecsum = _mm_add_pd(velecsum,velec);
206 vvdwsum = _mm_add_pd(vvdwsum,vvdw);
208 fscal = _mm_add_pd(felec,fvdw);
210 /* Calculate temporary vectorial force */
211 tx = _mm_mul_pd(fscal,dx00);
212 ty = _mm_mul_pd(fscal,dy00);
213 tz = _mm_mul_pd(fscal,dz00);
215 /* Update vectorial force */
216 fix0 = _mm_add_pd(fix0,tx);
217 fiy0 = _mm_add_pd(fiy0,ty);
218 fiz0 = _mm_add_pd(fiz0,tz);
220 gmx_mm_decrement_1rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,tx,ty,tz);
222 /* Inner loop uses 40 flops */
225 if(jidx<j_index_end)
228 jnrA = jjnr[jidx];
229 j_coord_offsetA = DIM*jnrA;
231 /* load j atom coordinates */
232 gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
233 &jx0,&jy0,&jz0);
235 /* Calculate displacement vector */
236 dx00 = _mm_sub_pd(ix0,jx0);
237 dy00 = _mm_sub_pd(iy0,jy0);
238 dz00 = _mm_sub_pd(iz0,jz0);
240 /* Calculate squared distance and things based on it */
241 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
243 rinv00 = gmx_mm_invsqrt_pd(rsq00);
245 rinvsq00 = _mm_mul_pd(rinv00,rinv00);
247 /* Load parameters for j particles */
248 jq0 = _mm_load_sd(charge+jnrA+0);
249 vdwjidx0A = 2*vdwtype[jnrA+0];
251 /**************************
252 * CALCULATE INTERACTIONS *
253 **************************/
255 /* Compute parameters for interactions between i and j atoms */
256 qq00 = _mm_mul_pd(iq0,jq0);
257 gmx_mm_load_1pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,&c6_00,&c12_00);
259 /* COULOMB ELECTROSTATICS */
260 velec = _mm_mul_pd(qq00,rinv00);
261 felec = _mm_mul_pd(velec,rinvsq00);
263 /* LENNARD-JONES DISPERSION/REPULSION */
265 rinvsix = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
266 vvdw6 = _mm_mul_pd(c6_00,rinvsix);
267 vvdw12 = _mm_mul_pd(c12_00,_mm_mul_pd(rinvsix,rinvsix));
268 vvdw = _mm_sub_pd( _mm_mul_pd(vvdw12,one_twelfth) , _mm_mul_pd(vvdw6,one_sixth) );
269 fvdw = _mm_mul_pd(_mm_sub_pd(vvdw12,vvdw6),rinvsq00);
271 /* Update potential sum for this i atom from the interaction with this j atom. */
272 velec = _mm_unpacklo_pd(velec,_mm_setzero_pd());
273 velecsum = _mm_add_pd(velecsum,velec);
274 vvdw = _mm_unpacklo_pd(vvdw,_mm_setzero_pd());
275 vvdwsum = _mm_add_pd(vvdwsum,vvdw);
277 fscal = _mm_add_pd(felec,fvdw);
279 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
281 /* Calculate temporary vectorial force */
282 tx = _mm_mul_pd(fscal,dx00);
283 ty = _mm_mul_pd(fscal,dy00);
284 tz = _mm_mul_pd(fscal,dz00);
286 /* Update vectorial force */
287 fix0 = _mm_add_pd(fix0,tx);
288 fiy0 = _mm_add_pd(fiy0,ty);
289 fiz0 = _mm_add_pd(fiz0,tz);
291 gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,tx,ty,tz);
293 /* Inner loop uses 40 flops */
296 /* End of innermost loop */
298 gmx_mm_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
299 f+i_coord_offset,fshift+i_shift_offset);
301 ggid = gid[iidx];
302 /* Update potential energies */
303 gmx_mm_update_1pot_pd(velecsum,kernel_data->energygrp_elec+ggid);
304 gmx_mm_update_1pot_pd(vvdwsum,kernel_data->energygrp_vdw+ggid);
306 /* Increment number of inner iterations */
307 inneriter += j_index_end - j_index_start;
309 /* Outer loop uses 9 flops */
312 /* Increment number of outer iterations */
313 outeriter += nri;
315 /* Update outer/inner flops */
317 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_VF,outeriter*9 + inneriter*40);
320 * Gromacs nonbonded kernel: nb_kernel_ElecCoul_VdwLJ_GeomP1P1_F_sse2_double
321 * Electrostatics interaction: Coulomb
322 * VdW interaction: LennardJones
323 * Geometry: Particle-Particle
324 * Calculate force/pot: Force
326 void
327 nb_kernel_ElecCoul_VdwLJ_GeomP1P1_F_sse2_double
328 (t_nblist * gmx_restrict nlist,
329 rvec * gmx_restrict xx,
330 rvec * gmx_restrict ff,
331 t_forcerec * gmx_restrict fr,
332 t_mdatoms * gmx_restrict mdatoms,
333 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
334 t_nrnb * gmx_restrict nrnb)
336 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
337 * just 0 for non-waters.
338 * Suffixes A,B refer to j loop unrolling done with SSE double precision, e.g. for the two different
339 * jnr indices corresponding to data put in the four positions in the SIMD register.
341 int i_shift_offset,i_coord_offset,outeriter,inneriter;
342 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
343 int jnrA,jnrB;
344 int j_coord_offsetA,j_coord_offsetB;
345 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
346 real rcutoff_scalar;
347 real *shiftvec,*fshift,*x,*f;
348 __m128d tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
349 int vdwioffset0;
350 __m128d ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
351 int vdwjidx0A,vdwjidx0B;
352 __m128d jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
353 __m128d dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
354 __m128d velec,felec,velecsum,facel,crf,krf,krf2;
355 real *charge;
356 int nvdwtype;
357 __m128d rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
358 int *vdwtype;
359 real *vdwparam;
360 __m128d one_sixth = _mm_set1_pd(1.0/6.0);
361 __m128d one_twelfth = _mm_set1_pd(1.0/12.0);
362 __m128d dummy_mask,cutoff_mask;
363 __m128d signbit = gmx_mm_castsi128_pd( _mm_set_epi32(0x80000000,0x00000000,0x80000000,0x00000000) );
364 __m128d one = _mm_set1_pd(1.0);
365 __m128d two = _mm_set1_pd(2.0);
366 x = xx[0];
367 f = ff[0];
369 nri = nlist->nri;
370 iinr = nlist->iinr;
371 jindex = nlist->jindex;
372 jjnr = nlist->jjnr;
373 shiftidx = nlist->shift;
374 gid = nlist->gid;
375 shiftvec = fr->shift_vec[0];
376 fshift = fr->fshift[0];
377 facel = _mm_set1_pd(fr->epsfac);
378 charge = mdatoms->chargeA;
379 nvdwtype = fr->ntype;
380 vdwparam = fr->nbfp;
381 vdwtype = mdatoms->typeA;
383 /* Avoid stupid compiler warnings */
384 jnrA = jnrB = 0;
385 j_coord_offsetA = 0;
386 j_coord_offsetB = 0;
388 outeriter = 0;
389 inneriter = 0;
391 /* Start outer loop over neighborlists */
392 for(iidx=0; iidx<nri; iidx++)
394 /* Load shift vector for this list */
395 i_shift_offset = DIM*shiftidx[iidx];
397 /* Load limits for loop over neighbors */
398 j_index_start = jindex[iidx];
399 j_index_end = jindex[iidx+1];
401 /* Get outer coordinate index */
402 inr = iinr[iidx];
403 i_coord_offset = DIM*inr;
405 /* Load i particle coords and add shift vector */
406 gmx_mm_load_shift_and_1rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
408 fix0 = _mm_setzero_pd();
409 fiy0 = _mm_setzero_pd();
410 fiz0 = _mm_setzero_pd();
412 /* Load parameters for i particles */
413 iq0 = _mm_mul_pd(facel,_mm_load1_pd(charge+inr+0));
414 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
416 /* Start inner kernel loop */
417 for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
420 /* Get j neighbor index, and coordinate index */
421 jnrA = jjnr[jidx];
422 jnrB = jjnr[jidx+1];
423 j_coord_offsetA = DIM*jnrA;
424 j_coord_offsetB = DIM*jnrB;
426 /* load j atom coordinates */
427 gmx_mm_load_1rvec_2ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
428 &jx0,&jy0,&jz0);
430 /* Calculate displacement vector */
431 dx00 = _mm_sub_pd(ix0,jx0);
432 dy00 = _mm_sub_pd(iy0,jy0);
433 dz00 = _mm_sub_pd(iz0,jz0);
435 /* Calculate squared distance and things based on it */
436 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
438 rinv00 = gmx_mm_invsqrt_pd(rsq00);
440 rinvsq00 = _mm_mul_pd(rinv00,rinv00);
442 /* Load parameters for j particles */
443 jq0 = gmx_mm_load_2real_swizzle_pd(charge+jnrA+0,charge+jnrB+0);
444 vdwjidx0A = 2*vdwtype[jnrA+0];
445 vdwjidx0B = 2*vdwtype[jnrB+0];
447 /**************************
448 * CALCULATE INTERACTIONS *
449 **************************/
451 /* Compute parameters for interactions between i and j atoms */
452 qq00 = _mm_mul_pd(iq0,jq0);
453 gmx_mm_load_2pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,
454 vdwparam+vdwioffset0+vdwjidx0B,&c6_00,&c12_00);
456 /* COULOMB ELECTROSTATICS */
457 velec = _mm_mul_pd(qq00,rinv00);
458 felec = _mm_mul_pd(velec,rinvsq00);
460 /* LENNARD-JONES DISPERSION/REPULSION */
462 rinvsix = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
463 fvdw = _mm_mul_pd(_mm_sub_pd(_mm_mul_pd(c12_00,rinvsix),c6_00),_mm_mul_pd(rinvsix,rinvsq00));
465 fscal = _mm_add_pd(felec,fvdw);
467 /* Calculate temporary vectorial force */
468 tx = _mm_mul_pd(fscal,dx00);
469 ty = _mm_mul_pd(fscal,dy00);
470 tz = _mm_mul_pd(fscal,dz00);
472 /* Update vectorial force */
473 fix0 = _mm_add_pd(fix0,tx);
474 fiy0 = _mm_add_pd(fiy0,ty);
475 fiz0 = _mm_add_pd(fiz0,tz);
477 gmx_mm_decrement_1rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,tx,ty,tz);
479 /* Inner loop uses 34 flops */
482 if(jidx<j_index_end)
485 jnrA = jjnr[jidx];
486 j_coord_offsetA = DIM*jnrA;
488 /* load j atom coordinates */
489 gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
490 &jx0,&jy0,&jz0);
492 /* Calculate displacement vector */
493 dx00 = _mm_sub_pd(ix0,jx0);
494 dy00 = _mm_sub_pd(iy0,jy0);
495 dz00 = _mm_sub_pd(iz0,jz0);
497 /* Calculate squared distance and things based on it */
498 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
500 rinv00 = gmx_mm_invsqrt_pd(rsq00);
502 rinvsq00 = _mm_mul_pd(rinv00,rinv00);
504 /* Load parameters for j particles */
505 jq0 = _mm_load_sd(charge+jnrA+0);
506 vdwjidx0A = 2*vdwtype[jnrA+0];
508 /**************************
509 * CALCULATE INTERACTIONS *
510 **************************/
512 /* Compute parameters for interactions between i and j atoms */
513 qq00 = _mm_mul_pd(iq0,jq0);
514 gmx_mm_load_1pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,&c6_00,&c12_00);
516 /* COULOMB ELECTROSTATICS */
517 velec = _mm_mul_pd(qq00,rinv00);
518 felec = _mm_mul_pd(velec,rinvsq00);
520 /* LENNARD-JONES DISPERSION/REPULSION */
522 rinvsix = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
523 fvdw = _mm_mul_pd(_mm_sub_pd(_mm_mul_pd(c12_00,rinvsix),c6_00),_mm_mul_pd(rinvsix,rinvsq00));
525 fscal = _mm_add_pd(felec,fvdw);
527 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
529 /* Calculate temporary vectorial force */
530 tx = _mm_mul_pd(fscal,dx00);
531 ty = _mm_mul_pd(fscal,dy00);
532 tz = _mm_mul_pd(fscal,dz00);
534 /* Update vectorial force */
535 fix0 = _mm_add_pd(fix0,tx);
536 fiy0 = _mm_add_pd(fiy0,ty);
537 fiz0 = _mm_add_pd(fiz0,tz);
539 gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,tx,ty,tz);
541 /* Inner loop uses 34 flops */
544 /* End of innermost loop */
546 gmx_mm_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
547 f+i_coord_offset,fshift+i_shift_offset);
549 /* Increment number of inner iterations */
550 inneriter += j_index_end - j_index_start;
552 /* Outer loop uses 7 flops */
555 /* Increment number of outer iterations */
556 outeriter += nri;
558 /* Update outer/inner flops */
560 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_F,outeriter*7 + inneriter*34);