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 sparc64_hpc_ace_double 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 "kernelutil_sparc64_hpc_ace_double.h"
50 * Gromacs nonbonded kernel: nb_kernel_ElecEwSh_VdwNone_GeomP1P1_VF_sparc64_hpc_ace_double
51 * Electrostatics interaction: Ewald
52 * VdW interaction: None
53 * Geometry: Particle-Particle
54 * Calculate force/pot: PotentialAndForce
57 nb_kernel_ElecEwSh_VdwNone_GeomP1P1_VF_sparc64_hpc_ace_double
58 (t_nblist * gmx_restrict nlist,
59 rvec * gmx_restrict xx,
60 rvec * gmx_restrict ff,
61 t_forcerec * gmx_restrict fr,
62 t_mdatoms * gmx_restrict mdatoms,
63 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
64 t_nrnb * gmx_restrict nrnb)
66 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
67 * just 0 for non-waters.
68 * Suffixes A,B refer to j loop unrolling done with double precision SIMD, e.g. for the two different
69 * jnr indices corresponding to data put in the four positions in the SIMD register.
71 int i_shift_offset,i_coord_offset,outeriter,inneriter;
72 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
74 int j_coord_offsetA,j_coord_offsetB;
75 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
77 real *shiftvec,*fshift,*x,*f;
78 _fjsp_v2r8 tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
80 _fjsp_v2r8 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
81 int vdwjidx0A,vdwjidx0B;
82 _fjsp_v2r8 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
83 _fjsp_v2r8 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
84 _fjsp_v2r8 velec,felec,velecsum,facel,crf,krf,krf2;
86 _fjsp_v2r8 ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
89 _fjsp_v2r8 dummy_mask,cutoff_mask;
90 _fjsp_v2r8 one = gmx_fjsp_set1_v2r8(1.0);
91 _fjsp_v2r8 two = gmx_fjsp_set1_v2r8(2.0);
92 union { _fjsp_v2r8 simd; long long int i[2]; } vfconv,gbconv,ewconv;
99 jindex = nlist->jindex;
101 shiftidx = nlist->shift;
103 shiftvec = fr->shift_vec[0];
104 fshift = fr->fshift[0];
105 facel = gmx_fjsp_set1_v2r8(fr->epsfac);
106 charge = mdatoms->chargeA;
108 sh_ewald = gmx_fjsp_set1_v2r8(fr->ic->sh_ewald);
109 ewtab = fr->ic->tabq_coul_FDV0;
110 ewtabscale = gmx_fjsp_set1_v2r8(fr->ic->tabq_scale);
111 ewtabhalfspace = gmx_fjsp_set1_v2r8(0.5/fr->ic->tabq_scale);
113 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
114 rcutoff_scalar = fr->rcoulomb;
115 rcutoff = gmx_fjsp_set1_v2r8(rcutoff_scalar);
116 rcutoff2 = _fjsp_mul_v2r8(rcutoff,rcutoff);
118 /* Avoid stupid compiler warnings */
126 /* Start outer loop over neighborlists */
127 for(iidx=0; iidx<nri; iidx++)
129 /* Load shift vector for this list */
130 i_shift_offset = DIM*shiftidx[iidx];
132 /* Load limits for loop over neighbors */
133 j_index_start = jindex[iidx];
134 j_index_end = jindex[iidx+1];
136 /* Get outer coordinate index */
138 i_coord_offset = DIM*inr;
140 /* Load i particle coords and add shift vector */
141 gmx_fjsp_load_shift_and_1rvec_broadcast_v2r8(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
143 fix0 = _fjsp_setzero_v2r8();
144 fiy0 = _fjsp_setzero_v2r8();
145 fiz0 = _fjsp_setzero_v2r8();
147 /* Load parameters for i particles */
148 iq0 = _fjsp_mul_v2r8(facel,gmx_fjsp_load1_v2r8(charge+inr+0));
150 /* Reset potential sums */
151 velecsum = _fjsp_setzero_v2r8();
153 /* Start inner kernel loop */
154 for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
157 /* Get j neighbor index, and coordinate index */
160 j_coord_offsetA = DIM*jnrA;
161 j_coord_offsetB = DIM*jnrB;
163 /* load j atom coordinates */
164 gmx_fjsp_load_1rvec_2ptr_swizzle_v2r8(x+j_coord_offsetA,x+j_coord_offsetB,
167 /* Calculate displacement vector */
168 dx00 = _fjsp_sub_v2r8(ix0,jx0);
169 dy00 = _fjsp_sub_v2r8(iy0,jy0);
170 dz00 = _fjsp_sub_v2r8(iz0,jz0);
172 /* Calculate squared distance and things based on it */
173 rsq00 = gmx_fjsp_calc_rsq_v2r8(dx00,dy00,dz00);
175 rinv00 = gmx_fjsp_invsqrt_v2r8(rsq00);
177 rinvsq00 = _fjsp_mul_v2r8(rinv00,rinv00);
179 /* Load parameters for j particles */
180 jq0 = gmx_fjsp_load_2real_swizzle_v2r8(charge+jnrA+0,charge+jnrB+0);
182 /**************************
183 * CALCULATE INTERACTIONS *
184 **************************/
186 if (gmx_fjsp_any_lt_v2r8(rsq00,rcutoff2))
189 r00 = _fjsp_mul_v2r8(rsq00,rinv00);
191 /* Compute parameters for interactions between i and j atoms */
192 qq00 = _fjsp_mul_v2r8(iq0,jq0);
194 /* EWALD ELECTROSTATICS */
196 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
197 ewrt = _fjsp_mul_v2r8(r00,ewtabscale);
198 itab_tmp = _fjsp_dtox_v2r8(ewrt);
199 eweps = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
200 _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
202 ewtabF = _fjsp_load_v2r8( ewtab + 4*ewconv.i[0] );
203 ewtabD = _fjsp_load_v2r8( ewtab + 4*ewconv.i[1] );
204 GMX_FJSP_TRANSPOSE2_V2R8(ewtabF,ewtabD);
205 ewtabV = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[0] +2);
206 ewtabFn = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[1] +2);
207 GMX_FJSP_TRANSPOSE2_V2R8(ewtabV,ewtabFn);
208 felec = _fjsp_madd_v2r8(eweps,ewtabD,ewtabF);
209 velec = _fjsp_nmsub_v2r8(_fjsp_mul_v2r8(ewtabhalfspace,eweps) ,_fjsp_add_v2r8(ewtabF,felec), ewtabV);
210 velec = _fjsp_mul_v2r8(qq00,_fjsp_sub_v2r8(_fjsp_sub_v2r8(rinv00,sh_ewald),velec));
211 felec = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq00,rinv00),_fjsp_sub_v2r8(rinvsq00,felec));
213 cutoff_mask = _fjsp_cmplt_v2r8(rsq00,rcutoff2);
215 /* Update potential sum for this i atom from the interaction with this j atom. */
216 velec = _fjsp_and_v2r8(velec,cutoff_mask);
217 velecsum = _fjsp_add_v2r8(velecsum,velec);
221 fscal = _fjsp_and_v2r8(fscal,cutoff_mask);
223 /* Update vectorial force */
224 fix0 = _fjsp_madd_v2r8(dx00,fscal,fix0);
225 fiy0 = _fjsp_madd_v2r8(dy00,fscal,fiy0);
226 fiz0 = _fjsp_madd_v2r8(dz00,fscal,fiz0);
228 gmx_fjsp_decrement_fma_1rvec_2ptr_swizzle_v2r8(f+j_coord_offsetA,f+j_coord_offsetB,fscal,dx00,dy00,dz00);
232 /* Inner loop uses 49 flops */
239 j_coord_offsetA = DIM*jnrA;
241 /* load j atom coordinates */
242 gmx_fjsp_load_1rvec_1ptr_swizzle_v2r8(x+j_coord_offsetA,
245 /* Calculate displacement vector */
246 dx00 = _fjsp_sub_v2r8(ix0,jx0);
247 dy00 = _fjsp_sub_v2r8(iy0,jy0);
248 dz00 = _fjsp_sub_v2r8(iz0,jz0);
250 /* Calculate squared distance and things based on it */
251 rsq00 = gmx_fjsp_calc_rsq_v2r8(dx00,dy00,dz00);
253 rinv00 = gmx_fjsp_invsqrt_v2r8(rsq00);
255 rinvsq00 = _fjsp_mul_v2r8(rinv00,rinv00);
257 /* Load parameters for j particles */
258 jq0 = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(),charge+jnrA+0);
260 /**************************
261 * CALCULATE INTERACTIONS *
262 **************************/
264 if (gmx_fjsp_any_lt_v2r8(rsq00,rcutoff2))
267 r00 = _fjsp_mul_v2r8(rsq00,rinv00);
269 /* Compute parameters for interactions between i and j atoms */
270 qq00 = _fjsp_mul_v2r8(iq0,jq0);
272 /* EWALD ELECTROSTATICS */
274 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
275 ewrt = _fjsp_mul_v2r8(r00,ewtabscale);
276 itab_tmp = _fjsp_dtox_v2r8(ewrt);
277 eweps = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
278 _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
280 ewtabF = _fjsp_load_v2r8( ewtab + 4*ewconv.i[0] );
281 ewtabD = _fjsp_setzero_v2r8();
282 GMX_FJSP_TRANSPOSE2_V2R8(ewtabF,ewtabD);
283 ewtabV = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[0] +2);
284 ewtabFn = _fjsp_setzero_v2r8();
285 GMX_FJSP_TRANSPOSE2_V2R8(ewtabV,ewtabFn);
286 felec = _fjsp_madd_v2r8(eweps,ewtabD,ewtabF);
287 velec = _fjsp_nmsub_v2r8(_fjsp_mul_v2r8(ewtabhalfspace,eweps) ,_fjsp_add_v2r8(ewtabF,felec), ewtabV);
288 velec = _fjsp_mul_v2r8(qq00,_fjsp_sub_v2r8(_fjsp_sub_v2r8(rinv00,sh_ewald),velec));
289 felec = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq00,rinv00),_fjsp_sub_v2r8(rinvsq00,felec));
291 cutoff_mask = _fjsp_cmplt_v2r8(rsq00,rcutoff2);
293 /* Update potential sum for this i atom from the interaction with this j atom. */
294 velec = _fjsp_and_v2r8(velec,cutoff_mask);
295 velec = _fjsp_unpacklo_v2r8(velec,_fjsp_setzero_v2r8());
296 velecsum = _fjsp_add_v2r8(velecsum,velec);
300 fscal = _fjsp_and_v2r8(fscal,cutoff_mask);
302 fscal = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
304 /* Update vectorial force */
305 fix0 = _fjsp_madd_v2r8(dx00,fscal,fix0);
306 fiy0 = _fjsp_madd_v2r8(dy00,fscal,fiy0);
307 fiz0 = _fjsp_madd_v2r8(dz00,fscal,fiz0);
309 gmx_fjsp_decrement_fma_1rvec_1ptr_swizzle_v2r8(f+j_coord_offsetA,fscal,dx00,dy00,dz00);
313 /* Inner loop uses 49 flops */
316 /* End of innermost loop */
318 gmx_fjsp_update_iforce_1atom_swizzle_v2r8(fix0,fiy0,fiz0,
319 f+i_coord_offset,fshift+i_shift_offset);
322 /* Update potential energies */
323 gmx_fjsp_update_1pot_v2r8(velecsum,kernel_data->energygrp_elec+ggid);
325 /* Increment number of inner iterations */
326 inneriter += j_index_end - j_index_start;
328 /* Outer loop uses 8 flops */
331 /* Increment number of outer iterations */
334 /* Update outer/inner flops */
336 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VF,outeriter*8 + inneriter*49);
339 * Gromacs nonbonded kernel: nb_kernel_ElecEwSh_VdwNone_GeomP1P1_F_sparc64_hpc_ace_double
340 * Electrostatics interaction: Ewald
341 * VdW interaction: None
342 * Geometry: Particle-Particle
343 * Calculate force/pot: Force
346 nb_kernel_ElecEwSh_VdwNone_GeomP1P1_F_sparc64_hpc_ace_double
347 (t_nblist * gmx_restrict nlist,
348 rvec * gmx_restrict xx,
349 rvec * gmx_restrict ff,
350 t_forcerec * gmx_restrict fr,
351 t_mdatoms * gmx_restrict mdatoms,
352 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
353 t_nrnb * gmx_restrict nrnb)
355 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
356 * just 0 for non-waters.
357 * Suffixes A,B refer to j loop unrolling done with double precision SIMD, e.g. for the two different
358 * jnr indices corresponding to data put in the four positions in the SIMD register.
360 int i_shift_offset,i_coord_offset,outeriter,inneriter;
361 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
363 int j_coord_offsetA,j_coord_offsetB;
364 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
366 real *shiftvec,*fshift,*x,*f;
367 _fjsp_v2r8 tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
369 _fjsp_v2r8 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
370 int vdwjidx0A,vdwjidx0B;
371 _fjsp_v2r8 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
372 _fjsp_v2r8 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
373 _fjsp_v2r8 velec,felec,velecsum,facel,crf,krf,krf2;
375 _fjsp_v2r8 ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
378 _fjsp_v2r8 dummy_mask,cutoff_mask;
379 _fjsp_v2r8 one = gmx_fjsp_set1_v2r8(1.0);
380 _fjsp_v2r8 two = gmx_fjsp_set1_v2r8(2.0);
381 union { _fjsp_v2r8 simd; long long int i[2]; } vfconv,gbconv,ewconv;
388 jindex = nlist->jindex;
390 shiftidx = nlist->shift;
392 shiftvec = fr->shift_vec[0];
393 fshift = fr->fshift[0];
394 facel = gmx_fjsp_set1_v2r8(fr->epsfac);
395 charge = mdatoms->chargeA;
397 sh_ewald = gmx_fjsp_set1_v2r8(fr->ic->sh_ewald);
398 ewtab = fr->ic->tabq_coul_F;
399 ewtabscale = gmx_fjsp_set1_v2r8(fr->ic->tabq_scale);
400 ewtabhalfspace = gmx_fjsp_set1_v2r8(0.5/fr->ic->tabq_scale);
402 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
403 rcutoff_scalar = fr->rcoulomb;
404 rcutoff = gmx_fjsp_set1_v2r8(rcutoff_scalar);
405 rcutoff2 = _fjsp_mul_v2r8(rcutoff,rcutoff);
407 /* Avoid stupid compiler warnings */
415 /* Start outer loop over neighborlists */
416 for(iidx=0; iidx<nri; iidx++)
418 /* Load shift vector for this list */
419 i_shift_offset = DIM*shiftidx[iidx];
421 /* Load limits for loop over neighbors */
422 j_index_start = jindex[iidx];
423 j_index_end = jindex[iidx+1];
425 /* Get outer coordinate index */
427 i_coord_offset = DIM*inr;
429 /* Load i particle coords and add shift vector */
430 gmx_fjsp_load_shift_and_1rvec_broadcast_v2r8(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
432 fix0 = _fjsp_setzero_v2r8();
433 fiy0 = _fjsp_setzero_v2r8();
434 fiz0 = _fjsp_setzero_v2r8();
436 /* Load parameters for i particles */
437 iq0 = _fjsp_mul_v2r8(facel,gmx_fjsp_load1_v2r8(charge+inr+0));
439 /* Start inner kernel loop */
440 for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
443 /* Get j neighbor index, and coordinate index */
446 j_coord_offsetA = DIM*jnrA;
447 j_coord_offsetB = DIM*jnrB;
449 /* load j atom coordinates */
450 gmx_fjsp_load_1rvec_2ptr_swizzle_v2r8(x+j_coord_offsetA,x+j_coord_offsetB,
453 /* Calculate displacement vector */
454 dx00 = _fjsp_sub_v2r8(ix0,jx0);
455 dy00 = _fjsp_sub_v2r8(iy0,jy0);
456 dz00 = _fjsp_sub_v2r8(iz0,jz0);
458 /* Calculate squared distance and things based on it */
459 rsq00 = gmx_fjsp_calc_rsq_v2r8(dx00,dy00,dz00);
461 rinv00 = gmx_fjsp_invsqrt_v2r8(rsq00);
463 rinvsq00 = _fjsp_mul_v2r8(rinv00,rinv00);
465 /* Load parameters for j particles */
466 jq0 = gmx_fjsp_load_2real_swizzle_v2r8(charge+jnrA+0,charge+jnrB+0);
468 /**************************
469 * CALCULATE INTERACTIONS *
470 **************************/
472 if (gmx_fjsp_any_lt_v2r8(rsq00,rcutoff2))
475 r00 = _fjsp_mul_v2r8(rsq00,rinv00);
477 /* Compute parameters for interactions between i and j atoms */
478 qq00 = _fjsp_mul_v2r8(iq0,jq0);
480 /* EWALD ELECTROSTATICS */
482 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
483 ewrt = _fjsp_mul_v2r8(r00,ewtabscale);
484 itab_tmp = _fjsp_dtox_v2r8(ewrt);
485 eweps = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
486 _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
488 gmx_fjsp_load_2pair_swizzle_v2r8(ewtab+ewconv.i[0],ewtab+ewconv.i[1],
490 felec = _fjsp_madd_v2r8(eweps,ewtabFn,_fjsp_nmsub_v2r8(eweps,ewtabF,ewtabF));
491 felec = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq00,rinv00),_fjsp_sub_v2r8(rinvsq00,felec));
493 cutoff_mask = _fjsp_cmplt_v2r8(rsq00,rcutoff2);
497 fscal = _fjsp_and_v2r8(fscal,cutoff_mask);
499 /* Update vectorial force */
500 fix0 = _fjsp_madd_v2r8(dx00,fscal,fix0);
501 fiy0 = _fjsp_madd_v2r8(dy00,fscal,fiy0);
502 fiz0 = _fjsp_madd_v2r8(dz00,fscal,fiz0);
504 gmx_fjsp_decrement_fma_1rvec_2ptr_swizzle_v2r8(f+j_coord_offsetA,f+j_coord_offsetB,fscal,dx00,dy00,dz00);
508 /* Inner loop uses 42 flops */
515 j_coord_offsetA = DIM*jnrA;
517 /* load j atom coordinates */
518 gmx_fjsp_load_1rvec_1ptr_swizzle_v2r8(x+j_coord_offsetA,
521 /* Calculate displacement vector */
522 dx00 = _fjsp_sub_v2r8(ix0,jx0);
523 dy00 = _fjsp_sub_v2r8(iy0,jy0);
524 dz00 = _fjsp_sub_v2r8(iz0,jz0);
526 /* Calculate squared distance and things based on it */
527 rsq00 = gmx_fjsp_calc_rsq_v2r8(dx00,dy00,dz00);
529 rinv00 = gmx_fjsp_invsqrt_v2r8(rsq00);
531 rinvsq00 = _fjsp_mul_v2r8(rinv00,rinv00);
533 /* Load parameters for j particles */
534 jq0 = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(),charge+jnrA+0);
536 /**************************
537 * CALCULATE INTERACTIONS *
538 **************************/
540 if (gmx_fjsp_any_lt_v2r8(rsq00,rcutoff2))
543 r00 = _fjsp_mul_v2r8(rsq00,rinv00);
545 /* Compute parameters for interactions between i and j atoms */
546 qq00 = _fjsp_mul_v2r8(iq0,jq0);
548 /* EWALD ELECTROSTATICS */
550 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
551 ewrt = _fjsp_mul_v2r8(r00,ewtabscale);
552 itab_tmp = _fjsp_dtox_v2r8(ewrt);
553 eweps = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
554 _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
556 gmx_fjsp_load_1pair_swizzle_v2r8(ewtab+ewconv.i[0],&ewtabF,&ewtabFn);
557 felec = _fjsp_madd_v2r8(eweps,ewtabFn,_fjsp_nmsub_v2r8(eweps,ewtabF,ewtabF));
558 felec = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq00,rinv00),_fjsp_sub_v2r8(rinvsq00,felec));
560 cutoff_mask = _fjsp_cmplt_v2r8(rsq00,rcutoff2);
564 fscal = _fjsp_and_v2r8(fscal,cutoff_mask);
566 fscal = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
568 /* Update vectorial force */
569 fix0 = _fjsp_madd_v2r8(dx00,fscal,fix0);
570 fiy0 = _fjsp_madd_v2r8(dy00,fscal,fiy0);
571 fiz0 = _fjsp_madd_v2r8(dz00,fscal,fiz0);
573 gmx_fjsp_decrement_fma_1rvec_1ptr_swizzle_v2r8(f+j_coord_offsetA,fscal,dx00,dy00,dz00);
577 /* Inner loop uses 42 flops */
580 /* End of innermost loop */
582 gmx_fjsp_update_iforce_1atom_swizzle_v2r8(fix0,fiy0,fiz0,
583 f+i_coord_offset,fshift+i_shift_offset);
585 /* Increment number of inner iterations */
586 inneriter += j_index_end - j_index_start;
588 /* Outer loop uses 7 flops */
591 /* Increment number of outer iterations */
594 /* Update outer/inner flops */
596 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_F,outeriter*7 + inneriter*42);