d77e5d67eb85d6badcefafdeb11af58f9259b7b8
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_sparc64_hpc_ace_double / nb_kernel_ElecEw_VdwLJEw_GeomW4P1_sparc64_hpc_ace_double.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35 /*
36  * Note: this file was generated by the GROMACS sparc64_hpc_ace_double kernel generator.
37  */
38 #include "config.h"
39
40 #include <math.h>
41
42 #include "../nb_kernel.h"
43 #include "types/simple.h"
44 #include "gromacs/math/vec.h"
45 #include "nrnb.h"
46
47 #include "kernelutil_sparc64_hpc_ace_double.h"
48
49 /*
50  * Gromacs nonbonded kernel:   nb_kernel_ElecEw_VdwLJEw_GeomW4P1_VF_sparc64_hpc_ace_double
51  * Electrostatics interaction: Ewald
52  * VdW interaction:            LJEwald
53  * Geometry:                   Water4-Particle
54  * Calculate force/pot:        PotentialAndForce
55  */
56 void
57 nb_kernel_ElecEw_VdwLJEw_GeomW4P1_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)
65 {
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.
70      */
71     int              i_shift_offset,i_coord_offset,outeriter,inneriter;
72     int              j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
73     int              jnrA,jnrB;
74     int              j_coord_offsetA,j_coord_offsetB;
75     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
76     real             rcutoff_scalar;
77     real             *shiftvec,*fshift,*x,*f;
78     _fjsp_v2r8       tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
79     int              vdwioffset0;
80     _fjsp_v2r8       ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
81     int              vdwioffset1;
82     _fjsp_v2r8       ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
83     int              vdwioffset2;
84     _fjsp_v2r8       ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
85     int              vdwioffset3;
86     _fjsp_v2r8       ix3,iy3,iz3,fix3,fiy3,fiz3,iq3,isai3;
87     int              vdwjidx0A,vdwjidx0B;
88     _fjsp_v2r8       jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
89     _fjsp_v2r8       dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
90     _fjsp_v2r8       dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
91     _fjsp_v2r8       dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
92     _fjsp_v2r8       dx30,dy30,dz30,rsq30,rinv30,rinvsq30,r30,qq30,c6_30,c12_30;
93     _fjsp_v2r8       velec,felec,velecsum,facel,crf,krf,krf2;
94     real             *charge;
95     int              nvdwtype;
96     _fjsp_v2r8       rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
97     int              *vdwtype;
98     real             *vdwparam;
99     _fjsp_v2r8       one_sixth   = gmx_fjsp_set1_v2r8(1.0/6.0);
100     _fjsp_v2r8       one_twelfth = gmx_fjsp_set1_v2r8(1.0/12.0);
101     _fjsp_v2r8           c6grid_00;
102     _fjsp_v2r8           c6grid_10;
103     _fjsp_v2r8           c6grid_20;
104     _fjsp_v2r8           c6grid_30;
105     real                 *vdwgridparam;
106     _fjsp_v2r8           ewclj,ewclj2,ewclj6,ewcljrsq,poly,exponent,f6A,f6B,sh_lj_ewald;
107     _fjsp_v2r8           one_half = gmx_fjsp_set1_v2r8(0.5);
108     _fjsp_v2r8           minus_one = gmx_fjsp_set1_v2r8(-1.0);
109     _fjsp_v2r8       ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
110     real             *ewtab;
111     _fjsp_v2r8       itab_tmp;
112     _fjsp_v2r8       dummy_mask,cutoff_mask;
113     _fjsp_v2r8       one     = gmx_fjsp_set1_v2r8(1.0);
114     _fjsp_v2r8       two     = gmx_fjsp_set1_v2r8(2.0);
115     union { _fjsp_v2r8 simd; long long int i[2]; } vfconv,gbconv,ewconv;
116
117     x                = xx[0];
118     f                = ff[0];
119
120     nri              = nlist->nri;
121     iinr             = nlist->iinr;
122     jindex           = nlist->jindex;
123     jjnr             = nlist->jjnr;
124     shiftidx         = nlist->shift;
125     gid              = nlist->gid;
126     shiftvec         = fr->shift_vec[0];
127     fshift           = fr->fshift[0];
128     facel            = gmx_fjsp_set1_v2r8(fr->epsfac);
129     charge           = mdatoms->chargeA;
130     nvdwtype         = fr->ntype;
131     vdwparam         = fr->nbfp;
132     vdwtype          = mdatoms->typeA;
133     vdwgridparam     = fr->ljpme_c6grid;
134     sh_lj_ewald      = gmx_fjsp_set1_v2r8(fr->ic->sh_lj_ewald);
135     ewclj            = gmx_fjsp_set1_v2r8(fr->ewaldcoeff_lj);
136     ewclj2           = _fjsp_mul_v2r8(minus_one,_fjsp_mul_v2r8(ewclj,ewclj));
137
138     sh_ewald         = gmx_fjsp_set1_v2r8(fr->ic->sh_ewald);
139     ewtab            = fr->ic->tabq_coul_FDV0;
140     ewtabscale       = gmx_fjsp_set1_v2r8(fr->ic->tabq_scale);
141     ewtabhalfspace   = gmx_fjsp_set1_v2r8(0.5/fr->ic->tabq_scale);
142
143     /* Setup water-specific parameters */
144     inr              = nlist->iinr[0];
145     iq1              = _fjsp_mul_v2r8(facel,gmx_fjsp_set1_v2r8(charge[inr+1]));
146     iq2              = _fjsp_mul_v2r8(facel,gmx_fjsp_set1_v2r8(charge[inr+2]));
147     iq3              = _fjsp_mul_v2r8(facel,gmx_fjsp_set1_v2r8(charge[inr+3]));
148     vdwioffset0      = 2*nvdwtype*vdwtype[inr+0];
149
150     /* Avoid stupid compiler warnings */
151     jnrA = jnrB = 0;
152     j_coord_offsetA = 0;
153     j_coord_offsetB = 0;
154
155     outeriter        = 0;
156     inneriter        = 0;
157
158     /* Start outer loop over neighborlists */
159     for(iidx=0; iidx<nri; iidx++)
160     {
161         /* Load shift vector for this list */
162         i_shift_offset   = DIM*shiftidx[iidx];
163
164         /* Load limits for loop over neighbors */
165         j_index_start    = jindex[iidx];
166         j_index_end      = jindex[iidx+1];
167
168         /* Get outer coordinate index */
169         inr              = iinr[iidx];
170         i_coord_offset   = DIM*inr;
171
172         /* Load i particle coords and add shift vector */
173         gmx_fjsp_load_shift_and_4rvec_broadcast_v2r8(shiftvec+i_shift_offset,x+i_coord_offset,
174                                                  &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
175
176         fix0             = _fjsp_setzero_v2r8();
177         fiy0             = _fjsp_setzero_v2r8();
178         fiz0             = _fjsp_setzero_v2r8();
179         fix1             = _fjsp_setzero_v2r8();
180         fiy1             = _fjsp_setzero_v2r8();
181         fiz1             = _fjsp_setzero_v2r8();
182         fix2             = _fjsp_setzero_v2r8();
183         fiy2             = _fjsp_setzero_v2r8();
184         fiz2             = _fjsp_setzero_v2r8();
185         fix3             = _fjsp_setzero_v2r8();
186         fiy3             = _fjsp_setzero_v2r8();
187         fiz3             = _fjsp_setzero_v2r8();
188
189         /* Reset potential sums */
190         velecsum         = _fjsp_setzero_v2r8();
191         vvdwsum          = _fjsp_setzero_v2r8();
192
193         /* Start inner kernel loop */
194         for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
195         {
196
197             /* Get j neighbor index, and coordinate index */
198             jnrA             = jjnr[jidx];
199             jnrB             = jjnr[jidx+1];
200             j_coord_offsetA  = DIM*jnrA;
201             j_coord_offsetB  = DIM*jnrB;
202
203             /* load j atom coordinates */
204             gmx_fjsp_load_1rvec_2ptr_swizzle_v2r8(x+j_coord_offsetA,x+j_coord_offsetB,
205                                               &jx0,&jy0,&jz0);
206
207             /* Calculate displacement vector */
208             dx00             = _fjsp_sub_v2r8(ix0,jx0);
209             dy00             = _fjsp_sub_v2r8(iy0,jy0);
210             dz00             = _fjsp_sub_v2r8(iz0,jz0);
211             dx10             = _fjsp_sub_v2r8(ix1,jx0);
212             dy10             = _fjsp_sub_v2r8(iy1,jy0);
213             dz10             = _fjsp_sub_v2r8(iz1,jz0);
214             dx20             = _fjsp_sub_v2r8(ix2,jx0);
215             dy20             = _fjsp_sub_v2r8(iy2,jy0);
216             dz20             = _fjsp_sub_v2r8(iz2,jz0);
217             dx30             = _fjsp_sub_v2r8(ix3,jx0);
218             dy30             = _fjsp_sub_v2r8(iy3,jy0);
219             dz30             = _fjsp_sub_v2r8(iz3,jz0);
220
221             /* Calculate squared distance and things based on it */
222             rsq00            = gmx_fjsp_calc_rsq_v2r8(dx00,dy00,dz00);
223             rsq10            = gmx_fjsp_calc_rsq_v2r8(dx10,dy10,dz10);
224             rsq20            = gmx_fjsp_calc_rsq_v2r8(dx20,dy20,dz20);
225             rsq30            = gmx_fjsp_calc_rsq_v2r8(dx30,dy30,dz30);
226
227             rinv00           = gmx_fjsp_invsqrt_v2r8(rsq00);
228             rinv10           = gmx_fjsp_invsqrt_v2r8(rsq10);
229             rinv20           = gmx_fjsp_invsqrt_v2r8(rsq20);
230             rinv30           = gmx_fjsp_invsqrt_v2r8(rsq30);
231
232             rinvsq00         = _fjsp_mul_v2r8(rinv00,rinv00);
233             rinvsq10         = _fjsp_mul_v2r8(rinv10,rinv10);
234             rinvsq20         = _fjsp_mul_v2r8(rinv20,rinv20);
235             rinvsq30         = _fjsp_mul_v2r8(rinv30,rinv30);
236
237             /* Load parameters for j particles */
238             jq0              = gmx_fjsp_load_2real_swizzle_v2r8(charge+jnrA+0,charge+jnrB+0);
239             vdwjidx0A        = 2*vdwtype[jnrA+0];
240             vdwjidx0B        = 2*vdwtype[jnrB+0];
241
242             fjx0             = _fjsp_setzero_v2r8();
243             fjy0             = _fjsp_setzero_v2r8();
244             fjz0             = _fjsp_setzero_v2r8();
245
246             /**************************
247              * CALCULATE INTERACTIONS *
248              **************************/
249
250             r00              = _fjsp_mul_v2r8(rsq00,rinv00);
251
252             /* Compute parameters for interactions between i and j atoms */
253             gmx_fjsp_load_2pair_swizzle_v2r8(vdwparam+vdwioffset0+vdwjidx0A,
254                                          vdwparam+vdwioffset0+vdwjidx0B,&c6_00,&c12_00);
255
256             c6grid_00       = gmx_fjsp_load_2real_swizzle_v2r8(vdwgridparam+vdwioffset0+vdwjidx0A,
257                                                                    vdwgridparam+vdwioffset0+vdwjidx0B);
258
259             /* Analytical LJ-PME */
260             rinvsix          = _fjsp_mul_v2r8(_fjsp_mul_v2r8(rinvsq00,rinvsq00),rinvsq00);
261             ewcljrsq         = _fjsp_mul_v2r8(ewclj2,rsq00);
262             ewclj6           = _fjsp_mul_v2r8(ewclj2,_fjsp_mul_v2r8(ewclj2,ewclj2));
263             exponent         = gmx_simd_exp_d(-ewcljrsq);
264             /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
265             poly             = _fjsp_mul_v2r8(exponent,_fjsp_madd_v2r8(_fjsp_mul_v2r8(ewcljrsq,ewcljrsq),one_half,_fjsp_sub_v2r8(one,ewcljrsq)));
266             /* vvdw6 = [C6 - C6grid * (1-poly)]/r6 */
267             vvdw6            = _fjsp_mul_v2r8(_fjsp_madd_v2r8(-c6grid_00,_fjsp_sub_v2r8(one,poly),c6_00),rinvsix);
268             vvdw12           = _fjsp_mul_v2r8(c12_00,_fjsp_mul_v2r8(rinvsix,rinvsix));
269             vvdw             = _fjsp_msub_v2r8(vvdw12,one_twelfth,_fjsp_mul_v2r8(vvdw6,one_sixth));         
270             /* fvdw = vvdw12/r - (vvdw6/r + (C6grid * exponent * beta^6)/r) */
271             fvdw             = _fjsp_mul_v2r8(_fjsp_add_v2r8(vvdw12,_fjsp_msub_v2r8(_fjsp_mul_v2r8(c6grid_00,one_sixth),_fjsp_mul_v2r8(exponent,ewclj6),vvdw6)),rinvsq00);
272
273             /* Update potential sum for this i atom from the interaction with this j atom. */
274             vvdwsum          = _fjsp_add_v2r8(vvdwsum,vvdw);
275
276             fscal            = fvdw;
277
278             /* Update vectorial force */
279             fix0             = _fjsp_madd_v2r8(dx00,fscal,fix0);
280             fiy0             = _fjsp_madd_v2r8(dy00,fscal,fiy0);
281             fiz0             = _fjsp_madd_v2r8(dz00,fscal,fiz0);
282             
283             fjx0             = _fjsp_madd_v2r8(dx00,fscal,fjx0);
284             fjy0             = _fjsp_madd_v2r8(dy00,fscal,fjy0);
285             fjz0             = _fjsp_madd_v2r8(dz00,fscal,fjz0);
286
287             /**************************
288              * CALCULATE INTERACTIONS *
289              **************************/
290
291             r10              = _fjsp_mul_v2r8(rsq10,rinv10);
292
293             /* Compute parameters for interactions between i and j atoms */
294             qq10             = _fjsp_mul_v2r8(iq1,jq0);
295
296             /* EWALD ELECTROSTATICS */
297
298             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
299             ewrt             = _fjsp_mul_v2r8(r10,ewtabscale);
300             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
301             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
302             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
303
304             ewtabF           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[0] );
305             ewtabD           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[1] );
306             GMX_FJSP_TRANSPOSE2_V2R8(ewtabF,ewtabD);
307             ewtabV           = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[0] +2);
308             ewtabFn          = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[1] +2);
309             GMX_FJSP_TRANSPOSE2_V2R8(ewtabV,ewtabFn);
310             felec            = _fjsp_madd_v2r8(eweps,ewtabD,ewtabF);
311             velec            = _fjsp_nmsub_v2r8(_fjsp_mul_v2r8(ewtabhalfspace,eweps) ,_fjsp_add_v2r8(ewtabF,felec), ewtabV);
312             velec            = _fjsp_mul_v2r8(qq10,_fjsp_sub_v2r8(rinv10,velec));
313             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq10,rinv10),_fjsp_sub_v2r8(rinvsq10,felec));
314
315             /* Update potential sum for this i atom from the interaction with this j atom. */
316             velecsum         = _fjsp_add_v2r8(velecsum,velec);
317
318             fscal            = felec;
319
320             /* Update vectorial force */
321             fix1             = _fjsp_madd_v2r8(dx10,fscal,fix1);
322             fiy1             = _fjsp_madd_v2r8(dy10,fscal,fiy1);
323             fiz1             = _fjsp_madd_v2r8(dz10,fscal,fiz1);
324             
325             fjx0             = _fjsp_madd_v2r8(dx10,fscal,fjx0);
326             fjy0             = _fjsp_madd_v2r8(dy10,fscal,fjy0);
327             fjz0             = _fjsp_madd_v2r8(dz10,fscal,fjz0);
328
329             /**************************
330              * CALCULATE INTERACTIONS *
331              **************************/
332
333             r20              = _fjsp_mul_v2r8(rsq20,rinv20);
334
335             /* Compute parameters for interactions between i and j atoms */
336             qq20             = _fjsp_mul_v2r8(iq2,jq0);
337
338             /* EWALD ELECTROSTATICS */
339
340             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
341             ewrt             = _fjsp_mul_v2r8(r20,ewtabscale);
342             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
343             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
344             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
345
346             ewtabF           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[0] );
347             ewtabD           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[1] );
348             GMX_FJSP_TRANSPOSE2_V2R8(ewtabF,ewtabD);
349             ewtabV           = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[0] +2);
350             ewtabFn          = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[1] +2);
351             GMX_FJSP_TRANSPOSE2_V2R8(ewtabV,ewtabFn);
352             felec            = _fjsp_madd_v2r8(eweps,ewtabD,ewtabF);
353             velec            = _fjsp_nmsub_v2r8(_fjsp_mul_v2r8(ewtabhalfspace,eweps) ,_fjsp_add_v2r8(ewtabF,felec), ewtabV);
354             velec            = _fjsp_mul_v2r8(qq20,_fjsp_sub_v2r8(rinv20,velec));
355             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq20,rinv20),_fjsp_sub_v2r8(rinvsq20,felec));
356
357             /* Update potential sum for this i atom from the interaction with this j atom. */
358             velecsum         = _fjsp_add_v2r8(velecsum,velec);
359
360             fscal            = felec;
361
362             /* Update vectorial force */
363             fix2             = _fjsp_madd_v2r8(dx20,fscal,fix2);
364             fiy2             = _fjsp_madd_v2r8(dy20,fscal,fiy2);
365             fiz2             = _fjsp_madd_v2r8(dz20,fscal,fiz2);
366             
367             fjx0             = _fjsp_madd_v2r8(dx20,fscal,fjx0);
368             fjy0             = _fjsp_madd_v2r8(dy20,fscal,fjy0);
369             fjz0             = _fjsp_madd_v2r8(dz20,fscal,fjz0);
370
371             /**************************
372              * CALCULATE INTERACTIONS *
373              **************************/
374
375             r30              = _fjsp_mul_v2r8(rsq30,rinv30);
376
377             /* Compute parameters for interactions between i and j atoms */
378             qq30             = _fjsp_mul_v2r8(iq3,jq0);
379
380             /* EWALD ELECTROSTATICS */
381
382             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
383             ewrt             = _fjsp_mul_v2r8(r30,ewtabscale);
384             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
385             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
386             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
387
388             ewtabF           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[0] );
389             ewtabD           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[1] );
390             GMX_FJSP_TRANSPOSE2_V2R8(ewtabF,ewtabD);
391             ewtabV           = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[0] +2);
392             ewtabFn          = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[1] +2);
393             GMX_FJSP_TRANSPOSE2_V2R8(ewtabV,ewtabFn);
394             felec            = _fjsp_madd_v2r8(eweps,ewtabD,ewtabF);
395             velec            = _fjsp_nmsub_v2r8(_fjsp_mul_v2r8(ewtabhalfspace,eweps) ,_fjsp_add_v2r8(ewtabF,felec), ewtabV);
396             velec            = _fjsp_mul_v2r8(qq30,_fjsp_sub_v2r8(rinv30,velec));
397             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq30,rinv30),_fjsp_sub_v2r8(rinvsq30,felec));
398
399             /* Update potential sum for this i atom from the interaction with this j atom. */
400             velecsum         = _fjsp_add_v2r8(velecsum,velec);
401
402             fscal            = felec;
403
404             /* Update vectorial force */
405             fix3             = _fjsp_madd_v2r8(dx30,fscal,fix3);
406             fiy3             = _fjsp_madd_v2r8(dy30,fscal,fiy3);
407             fiz3             = _fjsp_madd_v2r8(dz30,fscal,fiz3);
408             
409             fjx0             = _fjsp_madd_v2r8(dx30,fscal,fjx0);
410             fjy0             = _fjsp_madd_v2r8(dy30,fscal,fjy0);
411             fjz0             = _fjsp_madd_v2r8(dz30,fscal,fjz0);
412
413             gmx_fjsp_decrement_1rvec_2ptr_swizzle_v2r8(f+j_coord_offsetA,f+j_coord_offsetB,fjx0,fjy0,fjz0);
414
415             /* Inner loop uses 185 flops */
416         }
417
418         if(jidx<j_index_end)
419         {
420
421             jnrA             = jjnr[jidx];
422             j_coord_offsetA  = DIM*jnrA;
423
424             /* load j atom coordinates */
425             gmx_fjsp_load_1rvec_1ptr_swizzle_v2r8(x+j_coord_offsetA,
426                                               &jx0,&jy0,&jz0);
427
428             /* Calculate displacement vector */
429             dx00             = _fjsp_sub_v2r8(ix0,jx0);
430             dy00             = _fjsp_sub_v2r8(iy0,jy0);
431             dz00             = _fjsp_sub_v2r8(iz0,jz0);
432             dx10             = _fjsp_sub_v2r8(ix1,jx0);
433             dy10             = _fjsp_sub_v2r8(iy1,jy0);
434             dz10             = _fjsp_sub_v2r8(iz1,jz0);
435             dx20             = _fjsp_sub_v2r8(ix2,jx0);
436             dy20             = _fjsp_sub_v2r8(iy2,jy0);
437             dz20             = _fjsp_sub_v2r8(iz2,jz0);
438             dx30             = _fjsp_sub_v2r8(ix3,jx0);
439             dy30             = _fjsp_sub_v2r8(iy3,jy0);
440             dz30             = _fjsp_sub_v2r8(iz3,jz0);
441
442             /* Calculate squared distance and things based on it */
443             rsq00            = gmx_fjsp_calc_rsq_v2r8(dx00,dy00,dz00);
444             rsq10            = gmx_fjsp_calc_rsq_v2r8(dx10,dy10,dz10);
445             rsq20            = gmx_fjsp_calc_rsq_v2r8(dx20,dy20,dz20);
446             rsq30            = gmx_fjsp_calc_rsq_v2r8(dx30,dy30,dz30);
447
448             rinv00           = gmx_fjsp_invsqrt_v2r8(rsq00);
449             rinv10           = gmx_fjsp_invsqrt_v2r8(rsq10);
450             rinv20           = gmx_fjsp_invsqrt_v2r8(rsq20);
451             rinv30           = gmx_fjsp_invsqrt_v2r8(rsq30);
452
453             rinvsq00         = _fjsp_mul_v2r8(rinv00,rinv00);
454             rinvsq10         = _fjsp_mul_v2r8(rinv10,rinv10);
455             rinvsq20         = _fjsp_mul_v2r8(rinv20,rinv20);
456             rinvsq30         = _fjsp_mul_v2r8(rinv30,rinv30);
457
458             /* Load parameters for j particles */
459             jq0              = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(),charge+jnrA+0);
460             vdwjidx0A        = 2*vdwtype[jnrA+0];
461
462             fjx0             = _fjsp_setzero_v2r8();
463             fjy0             = _fjsp_setzero_v2r8();
464             fjz0             = _fjsp_setzero_v2r8();
465
466             /**************************
467              * CALCULATE INTERACTIONS *
468              **************************/
469
470             r00              = _fjsp_mul_v2r8(rsq00,rinv00);
471
472             /* Compute parameters for interactions between i and j atoms */
473             gmx_fjsp_load_1pair_swizzle_v2r8(vdwparam+vdwioffset0+vdwjidx0A,&c6_00,&c12_00);
474
475             c6grid_00       = gmx_fjsp_load_1real_swizzle_v2r8(vdwgridparam+vdwioffset0+vdwjidx0A);
476
477             /* Analytical LJ-PME */
478             rinvsix          = _fjsp_mul_v2r8(_fjsp_mul_v2r8(rinvsq00,rinvsq00),rinvsq00);
479             ewcljrsq         = _fjsp_mul_v2r8(ewclj2,rsq00);
480             ewclj6           = _fjsp_mul_v2r8(ewclj2,_fjsp_mul_v2r8(ewclj2,ewclj2));
481             exponent         = gmx_simd_exp_d(-ewcljrsq);
482             /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
483             poly             = _fjsp_mul_v2r8(exponent,_fjsp_madd_v2r8(_fjsp_mul_v2r8(ewcljrsq,ewcljrsq),one_half,_fjsp_sub_v2r8(one,ewcljrsq)));
484             /* vvdw6 = [C6 - C6grid * (1-poly)]/r6 */
485             vvdw6            = _fjsp_mul_v2r8(_fjsp_madd_v2r8(-c6grid_00,_fjsp_sub_v2r8(one,poly),c6_00),rinvsix);
486             vvdw12           = _fjsp_mul_v2r8(c12_00,_fjsp_mul_v2r8(rinvsix,rinvsix));
487             vvdw             = _fjsp_msub_v2r8(vvdw12,one_twelfth,_fjsp_mul_v2r8(vvdw6,one_sixth));         
488             /* fvdw = vvdw12/r - (vvdw6/r + (C6grid * exponent * beta^6)/r) */
489             fvdw             = _fjsp_mul_v2r8(_fjsp_add_v2r8(vvdw12,_fjsp_msub_v2r8(_fjsp_mul_v2r8(c6grid_00,one_sixth),_fjsp_mul_v2r8(exponent,ewclj6),vvdw6)),rinvsq00);
490
491             /* Update potential sum for this i atom from the interaction with this j atom. */
492             vvdw             = _fjsp_unpacklo_v2r8(vvdw,_fjsp_setzero_v2r8());
493             vvdwsum          = _fjsp_add_v2r8(vvdwsum,vvdw);
494
495             fscal            = fvdw;
496
497             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
498
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);
503             
504             fjx0             = _fjsp_madd_v2r8(dx00,fscal,fjx0);
505             fjy0             = _fjsp_madd_v2r8(dy00,fscal,fjy0);
506             fjz0             = _fjsp_madd_v2r8(dz00,fscal,fjz0);
507
508             /**************************
509              * CALCULATE INTERACTIONS *
510              **************************/
511
512             r10              = _fjsp_mul_v2r8(rsq10,rinv10);
513
514             /* Compute parameters for interactions between i and j atoms */
515             qq10             = _fjsp_mul_v2r8(iq1,jq0);
516
517             /* EWALD ELECTROSTATICS */
518
519             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
520             ewrt             = _fjsp_mul_v2r8(r10,ewtabscale);
521             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
522             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
523             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
524
525             ewtabF           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[0] );
526             ewtabD           = _fjsp_setzero_v2r8();
527             GMX_FJSP_TRANSPOSE2_V2R8(ewtabF,ewtabD);
528             ewtabV           = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[0] +2);
529             ewtabFn          = _fjsp_setzero_v2r8();
530             GMX_FJSP_TRANSPOSE2_V2R8(ewtabV,ewtabFn);
531             felec            = _fjsp_madd_v2r8(eweps,ewtabD,ewtabF);
532             velec            = _fjsp_nmsub_v2r8(_fjsp_mul_v2r8(ewtabhalfspace,eweps) ,_fjsp_add_v2r8(ewtabF,felec), ewtabV);
533             velec            = _fjsp_mul_v2r8(qq10,_fjsp_sub_v2r8(rinv10,velec));
534             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq10,rinv10),_fjsp_sub_v2r8(rinvsq10,felec));
535
536             /* Update potential sum for this i atom from the interaction with this j atom. */
537             velec            = _fjsp_unpacklo_v2r8(velec,_fjsp_setzero_v2r8());
538             velecsum         = _fjsp_add_v2r8(velecsum,velec);
539
540             fscal            = felec;
541
542             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
543
544             /* Update vectorial force */
545             fix1             = _fjsp_madd_v2r8(dx10,fscal,fix1);
546             fiy1             = _fjsp_madd_v2r8(dy10,fscal,fiy1);
547             fiz1             = _fjsp_madd_v2r8(dz10,fscal,fiz1);
548             
549             fjx0             = _fjsp_madd_v2r8(dx10,fscal,fjx0);
550             fjy0             = _fjsp_madd_v2r8(dy10,fscal,fjy0);
551             fjz0             = _fjsp_madd_v2r8(dz10,fscal,fjz0);
552
553             /**************************
554              * CALCULATE INTERACTIONS *
555              **************************/
556
557             r20              = _fjsp_mul_v2r8(rsq20,rinv20);
558
559             /* Compute parameters for interactions between i and j atoms */
560             qq20             = _fjsp_mul_v2r8(iq2,jq0);
561
562             /* EWALD ELECTROSTATICS */
563
564             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
565             ewrt             = _fjsp_mul_v2r8(r20,ewtabscale);
566             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
567             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
568             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
569
570             ewtabF           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[0] );
571             ewtabD           = _fjsp_setzero_v2r8();
572             GMX_FJSP_TRANSPOSE2_V2R8(ewtabF,ewtabD);
573             ewtabV           = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[0] +2);
574             ewtabFn          = _fjsp_setzero_v2r8();
575             GMX_FJSP_TRANSPOSE2_V2R8(ewtabV,ewtabFn);
576             felec            = _fjsp_madd_v2r8(eweps,ewtabD,ewtabF);
577             velec            = _fjsp_nmsub_v2r8(_fjsp_mul_v2r8(ewtabhalfspace,eweps) ,_fjsp_add_v2r8(ewtabF,felec), ewtabV);
578             velec            = _fjsp_mul_v2r8(qq20,_fjsp_sub_v2r8(rinv20,velec));
579             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq20,rinv20),_fjsp_sub_v2r8(rinvsq20,felec));
580
581             /* Update potential sum for this i atom from the interaction with this j atom. */
582             velec            = _fjsp_unpacklo_v2r8(velec,_fjsp_setzero_v2r8());
583             velecsum         = _fjsp_add_v2r8(velecsum,velec);
584
585             fscal            = felec;
586
587             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
588
589             /* Update vectorial force */
590             fix2             = _fjsp_madd_v2r8(dx20,fscal,fix2);
591             fiy2             = _fjsp_madd_v2r8(dy20,fscal,fiy2);
592             fiz2             = _fjsp_madd_v2r8(dz20,fscal,fiz2);
593             
594             fjx0             = _fjsp_madd_v2r8(dx20,fscal,fjx0);
595             fjy0             = _fjsp_madd_v2r8(dy20,fscal,fjy0);
596             fjz0             = _fjsp_madd_v2r8(dz20,fscal,fjz0);
597
598             /**************************
599              * CALCULATE INTERACTIONS *
600              **************************/
601
602             r30              = _fjsp_mul_v2r8(rsq30,rinv30);
603
604             /* Compute parameters for interactions between i and j atoms */
605             qq30             = _fjsp_mul_v2r8(iq3,jq0);
606
607             /* EWALD ELECTROSTATICS */
608
609             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
610             ewrt             = _fjsp_mul_v2r8(r30,ewtabscale);
611             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
612             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
613             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
614
615             ewtabF           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[0] );
616             ewtabD           = _fjsp_setzero_v2r8();
617             GMX_FJSP_TRANSPOSE2_V2R8(ewtabF,ewtabD);
618             ewtabV           = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[0] +2);
619             ewtabFn          = _fjsp_setzero_v2r8();
620             GMX_FJSP_TRANSPOSE2_V2R8(ewtabV,ewtabFn);
621             felec            = _fjsp_madd_v2r8(eweps,ewtabD,ewtabF);
622             velec            = _fjsp_nmsub_v2r8(_fjsp_mul_v2r8(ewtabhalfspace,eweps) ,_fjsp_add_v2r8(ewtabF,felec), ewtabV);
623             velec            = _fjsp_mul_v2r8(qq30,_fjsp_sub_v2r8(rinv30,velec));
624             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq30,rinv30),_fjsp_sub_v2r8(rinvsq30,felec));
625
626             /* Update potential sum for this i atom from the interaction with this j atom. */
627             velec            = _fjsp_unpacklo_v2r8(velec,_fjsp_setzero_v2r8());
628             velecsum         = _fjsp_add_v2r8(velecsum,velec);
629
630             fscal            = felec;
631
632             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
633
634             /* Update vectorial force */
635             fix3             = _fjsp_madd_v2r8(dx30,fscal,fix3);
636             fiy3             = _fjsp_madd_v2r8(dy30,fscal,fiy3);
637             fiz3             = _fjsp_madd_v2r8(dz30,fscal,fiz3);
638             
639             fjx0             = _fjsp_madd_v2r8(dx30,fscal,fjx0);
640             fjy0             = _fjsp_madd_v2r8(dy30,fscal,fjy0);
641             fjz0             = _fjsp_madd_v2r8(dz30,fscal,fjz0);
642
643             gmx_fjsp_decrement_1rvec_1ptr_swizzle_v2r8(f+j_coord_offsetA,fjx0,fjy0,fjz0);
644
645             /* Inner loop uses 185 flops */
646         }
647
648         /* End of innermost loop */
649
650         gmx_fjsp_update_iforce_4atom_swizzle_v2r8(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
651                                               f+i_coord_offset,fshift+i_shift_offset);
652
653         ggid                        = gid[iidx];
654         /* Update potential energies */
655         gmx_fjsp_update_1pot_v2r8(velecsum,kernel_data->energygrp_elec+ggid);
656         gmx_fjsp_update_1pot_v2r8(vvdwsum,kernel_data->energygrp_vdw+ggid);
657
658         /* Increment number of inner iterations */
659         inneriter                  += j_index_end - j_index_start;
660
661         /* Outer loop uses 26 flops */
662     }
663
664     /* Increment number of outer iterations */
665     outeriter        += nri;
666
667     /* Update outer/inner flops */
668
669     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W4_VF,outeriter*26 + inneriter*185);
670 }
671 /*
672  * Gromacs nonbonded kernel:   nb_kernel_ElecEw_VdwLJEw_GeomW4P1_F_sparc64_hpc_ace_double
673  * Electrostatics interaction: Ewald
674  * VdW interaction:            LJEwald
675  * Geometry:                   Water4-Particle
676  * Calculate force/pot:        Force
677  */
678 void
679 nb_kernel_ElecEw_VdwLJEw_GeomW4P1_F_sparc64_hpc_ace_double
680                     (t_nblist                    * gmx_restrict       nlist,
681                      rvec                        * gmx_restrict          xx,
682                      rvec                        * gmx_restrict          ff,
683                      t_forcerec                  * gmx_restrict          fr,
684                      t_mdatoms                   * gmx_restrict     mdatoms,
685                      nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
686                      t_nrnb                      * gmx_restrict        nrnb)
687 {
688     /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
689      * just 0 for non-waters.
690      * Suffixes A,B refer to j loop unrolling done with double precision SIMD, e.g. for the two different
691      * jnr indices corresponding to data put in the four positions in the SIMD register.
692      */
693     int              i_shift_offset,i_coord_offset,outeriter,inneriter;
694     int              j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
695     int              jnrA,jnrB;
696     int              j_coord_offsetA,j_coord_offsetB;
697     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
698     real             rcutoff_scalar;
699     real             *shiftvec,*fshift,*x,*f;
700     _fjsp_v2r8       tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
701     int              vdwioffset0;
702     _fjsp_v2r8       ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
703     int              vdwioffset1;
704     _fjsp_v2r8       ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
705     int              vdwioffset2;
706     _fjsp_v2r8       ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
707     int              vdwioffset3;
708     _fjsp_v2r8       ix3,iy3,iz3,fix3,fiy3,fiz3,iq3,isai3;
709     int              vdwjidx0A,vdwjidx0B;
710     _fjsp_v2r8       jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
711     _fjsp_v2r8       dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
712     _fjsp_v2r8       dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
713     _fjsp_v2r8       dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
714     _fjsp_v2r8       dx30,dy30,dz30,rsq30,rinv30,rinvsq30,r30,qq30,c6_30,c12_30;
715     _fjsp_v2r8       velec,felec,velecsum,facel,crf,krf,krf2;
716     real             *charge;
717     int              nvdwtype;
718     _fjsp_v2r8       rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
719     int              *vdwtype;
720     real             *vdwparam;
721     _fjsp_v2r8       one_sixth   = gmx_fjsp_set1_v2r8(1.0/6.0);
722     _fjsp_v2r8       one_twelfth = gmx_fjsp_set1_v2r8(1.0/12.0);
723     _fjsp_v2r8           c6grid_00;
724     _fjsp_v2r8           c6grid_10;
725     _fjsp_v2r8           c6grid_20;
726     _fjsp_v2r8           c6grid_30;
727     real                 *vdwgridparam;
728     _fjsp_v2r8           ewclj,ewclj2,ewclj6,ewcljrsq,poly,exponent,f6A,f6B,sh_lj_ewald;
729     _fjsp_v2r8           one_half = gmx_fjsp_set1_v2r8(0.5);
730     _fjsp_v2r8           minus_one = gmx_fjsp_set1_v2r8(-1.0);
731     _fjsp_v2r8       ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
732     real             *ewtab;
733     _fjsp_v2r8       itab_tmp;
734     _fjsp_v2r8       dummy_mask,cutoff_mask;
735     _fjsp_v2r8       one     = gmx_fjsp_set1_v2r8(1.0);
736     _fjsp_v2r8       two     = gmx_fjsp_set1_v2r8(2.0);
737     union { _fjsp_v2r8 simd; long long int i[2]; } vfconv,gbconv,ewconv;
738
739     x                = xx[0];
740     f                = ff[0];
741
742     nri              = nlist->nri;
743     iinr             = nlist->iinr;
744     jindex           = nlist->jindex;
745     jjnr             = nlist->jjnr;
746     shiftidx         = nlist->shift;
747     gid              = nlist->gid;
748     shiftvec         = fr->shift_vec[0];
749     fshift           = fr->fshift[0];
750     facel            = gmx_fjsp_set1_v2r8(fr->epsfac);
751     charge           = mdatoms->chargeA;
752     nvdwtype         = fr->ntype;
753     vdwparam         = fr->nbfp;
754     vdwtype          = mdatoms->typeA;
755     vdwgridparam     = fr->ljpme_c6grid;
756     sh_lj_ewald      = gmx_fjsp_set1_v2r8(fr->ic->sh_lj_ewald);
757     ewclj            = gmx_fjsp_set1_v2r8(fr->ewaldcoeff_lj);
758     ewclj2           = _fjsp_mul_v2r8(minus_one,_fjsp_mul_v2r8(ewclj,ewclj));
759
760     sh_ewald         = gmx_fjsp_set1_v2r8(fr->ic->sh_ewald);
761     ewtab            = fr->ic->tabq_coul_F;
762     ewtabscale       = gmx_fjsp_set1_v2r8(fr->ic->tabq_scale);
763     ewtabhalfspace   = gmx_fjsp_set1_v2r8(0.5/fr->ic->tabq_scale);
764
765     /* Setup water-specific parameters */
766     inr              = nlist->iinr[0];
767     iq1              = _fjsp_mul_v2r8(facel,gmx_fjsp_set1_v2r8(charge[inr+1]));
768     iq2              = _fjsp_mul_v2r8(facel,gmx_fjsp_set1_v2r8(charge[inr+2]));
769     iq3              = _fjsp_mul_v2r8(facel,gmx_fjsp_set1_v2r8(charge[inr+3]));
770     vdwioffset0      = 2*nvdwtype*vdwtype[inr+0];
771
772     /* Avoid stupid compiler warnings */
773     jnrA = jnrB = 0;
774     j_coord_offsetA = 0;
775     j_coord_offsetB = 0;
776
777     outeriter        = 0;
778     inneriter        = 0;
779
780     /* Start outer loop over neighborlists */
781     for(iidx=0; iidx<nri; iidx++)
782     {
783         /* Load shift vector for this list */
784         i_shift_offset   = DIM*shiftidx[iidx];
785
786         /* Load limits for loop over neighbors */
787         j_index_start    = jindex[iidx];
788         j_index_end      = jindex[iidx+1];
789
790         /* Get outer coordinate index */
791         inr              = iinr[iidx];
792         i_coord_offset   = DIM*inr;
793
794         /* Load i particle coords and add shift vector */
795         gmx_fjsp_load_shift_and_4rvec_broadcast_v2r8(shiftvec+i_shift_offset,x+i_coord_offset,
796                                                  &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
797
798         fix0             = _fjsp_setzero_v2r8();
799         fiy0             = _fjsp_setzero_v2r8();
800         fiz0             = _fjsp_setzero_v2r8();
801         fix1             = _fjsp_setzero_v2r8();
802         fiy1             = _fjsp_setzero_v2r8();
803         fiz1             = _fjsp_setzero_v2r8();
804         fix2             = _fjsp_setzero_v2r8();
805         fiy2             = _fjsp_setzero_v2r8();
806         fiz2             = _fjsp_setzero_v2r8();
807         fix3             = _fjsp_setzero_v2r8();
808         fiy3             = _fjsp_setzero_v2r8();
809         fiz3             = _fjsp_setzero_v2r8();
810
811         /* Start inner kernel loop */
812         for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
813         {
814
815             /* Get j neighbor index, and coordinate index */
816             jnrA             = jjnr[jidx];
817             jnrB             = jjnr[jidx+1];
818             j_coord_offsetA  = DIM*jnrA;
819             j_coord_offsetB  = DIM*jnrB;
820
821             /* load j atom coordinates */
822             gmx_fjsp_load_1rvec_2ptr_swizzle_v2r8(x+j_coord_offsetA,x+j_coord_offsetB,
823                                               &jx0,&jy0,&jz0);
824
825             /* Calculate displacement vector */
826             dx00             = _fjsp_sub_v2r8(ix0,jx0);
827             dy00             = _fjsp_sub_v2r8(iy0,jy0);
828             dz00             = _fjsp_sub_v2r8(iz0,jz0);
829             dx10             = _fjsp_sub_v2r8(ix1,jx0);
830             dy10             = _fjsp_sub_v2r8(iy1,jy0);
831             dz10             = _fjsp_sub_v2r8(iz1,jz0);
832             dx20             = _fjsp_sub_v2r8(ix2,jx0);
833             dy20             = _fjsp_sub_v2r8(iy2,jy0);
834             dz20             = _fjsp_sub_v2r8(iz2,jz0);
835             dx30             = _fjsp_sub_v2r8(ix3,jx0);
836             dy30             = _fjsp_sub_v2r8(iy3,jy0);
837             dz30             = _fjsp_sub_v2r8(iz3,jz0);
838
839             /* Calculate squared distance and things based on it */
840             rsq00            = gmx_fjsp_calc_rsq_v2r8(dx00,dy00,dz00);
841             rsq10            = gmx_fjsp_calc_rsq_v2r8(dx10,dy10,dz10);
842             rsq20            = gmx_fjsp_calc_rsq_v2r8(dx20,dy20,dz20);
843             rsq30            = gmx_fjsp_calc_rsq_v2r8(dx30,dy30,dz30);
844
845             rinv00           = gmx_fjsp_invsqrt_v2r8(rsq00);
846             rinv10           = gmx_fjsp_invsqrt_v2r8(rsq10);
847             rinv20           = gmx_fjsp_invsqrt_v2r8(rsq20);
848             rinv30           = gmx_fjsp_invsqrt_v2r8(rsq30);
849
850             rinvsq00         = _fjsp_mul_v2r8(rinv00,rinv00);
851             rinvsq10         = _fjsp_mul_v2r8(rinv10,rinv10);
852             rinvsq20         = _fjsp_mul_v2r8(rinv20,rinv20);
853             rinvsq30         = _fjsp_mul_v2r8(rinv30,rinv30);
854
855             /* Load parameters for j particles */
856             jq0              = gmx_fjsp_load_2real_swizzle_v2r8(charge+jnrA+0,charge+jnrB+0);
857             vdwjidx0A        = 2*vdwtype[jnrA+0];
858             vdwjidx0B        = 2*vdwtype[jnrB+0];
859
860             fjx0             = _fjsp_setzero_v2r8();
861             fjy0             = _fjsp_setzero_v2r8();
862             fjz0             = _fjsp_setzero_v2r8();
863
864             /**************************
865              * CALCULATE INTERACTIONS *
866              **************************/
867
868             r00              = _fjsp_mul_v2r8(rsq00,rinv00);
869
870             /* Compute parameters for interactions between i and j atoms */
871             gmx_fjsp_load_2pair_swizzle_v2r8(vdwparam+vdwioffset0+vdwjidx0A,
872                                          vdwparam+vdwioffset0+vdwjidx0B,&c6_00,&c12_00);
873
874             c6grid_00       = gmx_fjsp_load_2real_swizzle_v2r8(vdwgridparam+vdwioffset0+vdwjidx0A,
875                                                                    vdwgridparam+vdwioffset0+vdwjidx0B);
876
877             /* Analytical LJ-PME */
878             rinvsix          = _fjsp_mul_v2r8(_fjsp_mul_v2r8(rinvsq00,rinvsq00),rinvsq00);
879             ewcljrsq         = _fjsp_mul_v2r8(ewclj2,rsq00);
880             ewclj6           = _fjsp_mul_v2r8(ewclj2,_fjsp_mul_v2r8(ewclj2,ewclj2));
881             exponent         = gmx_simd_exp_d(-ewcljrsq);
882             /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
883             poly             = _fjsp_mul_v2r8(exponent,_fjsp_madd_v2r8(_fjsp_mul_v2r8(ewcljrsq,ewcljrsq),one_half,_fjsp_sub_v2r8(one,ewcljrsq)));
884             /* f6A = 6 * C6grid * (1 - poly) */
885             f6A              = _fjsp_mul_v2r8(c6grid_00,_fjsp_msub_v2r8(one,poly));
886             /* f6B = C6grid * exponent * beta^6 */
887             f6B              = _fjsp_mul_v2r8(_fjsp_mul_v2r8(c6grid_00,one_sixth),_fjsp_mul_v2r8(exponent,ewclj6));
888             /* fvdw = 12*C12/r13 - ((6*C6 - f6A)/r6 + f6B)/r */
889             fvdw              = _fjsp_mul_v2r8(_fjsp_madd_v2r8(_fjsp_msub_v2r8(c12_00,rinvsix,_fjsp_sub_v2r8(c6_00,f6A)),rinvsix,f6B),rinvsq00);
890
891             fscal            = fvdw;
892
893             /* Update vectorial force */
894             fix0             = _fjsp_madd_v2r8(dx00,fscal,fix0);
895             fiy0             = _fjsp_madd_v2r8(dy00,fscal,fiy0);
896             fiz0             = _fjsp_madd_v2r8(dz00,fscal,fiz0);
897             
898             fjx0             = _fjsp_madd_v2r8(dx00,fscal,fjx0);
899             fjy0             = _fjsp_madd_v2r8(dy00,fscal,fjy0);
900             fjz0             = _fjsp_madd_v2r8(dz00,fscal,fjz0);
901
902             /**************************
903              * CALCULATE INTERACTIONS *
904              **************************/
905
906             r10              = _fjsp_mul_v2r8(rsq10,rinv10);
907
908             /* Compute parameters for interactions between i and j atoms */
909             qq10             = _fjsp_mul_v2r8(iq1,jq0);
910
911             /* EWALD ELECTROSTATICS */
912
913             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
914             ewrt             = _fjsp_mul_v2r8(r10,ewtabscale);
915             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
916             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
917             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
918
919             gmx_fjsp_load_2pair_swizzle_v2r8(ewtab+ewconv.i[0],ewtab+ewconv.i[1],
920                                          &ewtabF,&ewtabFn);
921             felec            = _fjsp_madd_v2r8(eweps,ewtabFn,_fjsp_nmsub_v2r8(eweps,ewtabF,ewtabF));
922             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq10,rinv10),_fjsp_sub_v2r8(rinvsq10,felec));
923
924             fscal            = felec;
925
926             /* Update vectorial force */
927             fix1             = _fjsp_madd_v2r8(dx10,fscal,fix1);
928             fiy1             = _fjsp_madd_v2r8(dy10,fscal,fiy1);
929             fiz1             = _fjsp_madd_v2r8(dz10,fscal,fiz1);
930             
931             fjx0             = _fjsp_madd_v2r8(dx10,fscal,fjx0);
932             fjy0             = _fjsp_madd_v2r8(dy10,fscal,fjy0);
933             fjz0             = _fjsp_madd_v2r8(dz10,fscal,fjz0);
934
935             /**************************
936              * CALCULATE INTERACTIONS *
937              **************************/
938
939             r20              = _fjsp_mul_v2r8(rsq20,rinv20);
940
941             /* Compute parameters for interactions between i and j atoms */
942             qq20             = _fjsp_mul_v2r8(iq2,jq0);
943
944             /* EWALD ELECTROSTATICS */
945
946             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
947             ewrt             = _fjsp_mul_v2r8(r20,ewtabscale);
948             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
949             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
950             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
951
952             gmx_fjsp_load_2pair_swizzle_v2r8(ewtab+ewconv.i[0],ewtab+ewconv.i[1],
953                                          &ewtabF,&ewtabFn);
954             felec            = _fjsp_madd_v2r8(eweps,ewtabFn,_fjsp_nmsub_v2r8(eweps,ewtabF,ewtabF));
955             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq20,rinv20),_fjsp_sub_v2r8(rinvsq20,felec));
956
957             fscal            = felec;
958
959             /* Update vectorial force */
960             fix2             = _fjsp_madd_v2r8(dx20,fscal,fix2);
961             fiy2             = _fjsp_madd_v2r8(dy20,fscal,fiy2);
962             fiz2             = _fjsp_madd_v2r8(dz20,fscal,fiz2);
963             
964             fjx0             = _fjsp_madd_v2r8(dx20,fscal,fjx0);
965             fjy0             = _fjsp_madd_v2r8(dy20,fscal,fjy0);
966             fjz0             = _fjsp_madd_v2r8(dz20,fscal,fjz0);
967
968             /**************************
969              * CALCULATE INTERACTIONS *
970              **************************/
971
972             r30              = _fjsp_mul_v2r8(rsq30,rinv30);
973
974             /* Compute parameters for interactions between i and j atoms */
975             qq30             = _fjsp_mul_v2r8(iq3,jq0);
976
977             /* EWALD ELECTROSTATICS */
978
979             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
980             ewrt             = _fjsp_mul_v2r8(r30,ewtabscale);
981             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
982             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
983             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
984
985             gmx_fjsp_load_2pair_swizzle_v2r8(ewtab+ewconv.i[0],ewtab+ewconv.i[1],
986                                          &ewtabF,&ewtabFn);
987             felec            = _fjsp_madd_v2r8(eweps,ewtabFn,_fjsp_nmsub_v2r8(eweps,ewtabF,ewtabF));
988             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq30,rinv30),_fjsp_sub_v2r8(rinvsq30,felec));
989
990             fscal            = felec;
991
992             /* Update vectorial force */
993             fix3             = _fjsp_madd_v2r8(dx30,fscal,fix3);
994             fiy3             = _fjsp_madd_v2r8(dy30,fscal,fiy3);
995             fiz3             = _fjsp_madd_v2r8(dz30,fscal,fiz3);
996             
997             fjx0             = _fjsp_madd_v2r8(dx30,fscal,fjx0);
998             fjy0             = _fjsp_madd_v2r8(dy30,fscal,fjy0);
999             fjz0             = _fjsp_madd_v2r8(dz30,fscal,fjz0);
1000
1001             gmx_fjsp_decrement_1rvec_2ptr_swizzle_v2r8(f+j_coord_offsetA,f+j_coord_offsetB,fjx0,fjy0,fjz0);
1002
1003             /* Inner loop uses 168 flops */
1004         }
1005
1006         if(jidx<j_index_end)
1007         {
1008
1009             jnrA             = jjnr[jidx];
1010             j_coord_offsetA  = DIM*jnrA;
1011
1012             /* load j atom coordinates */
1013             gmx_fjsp_load_1rvec_1ptr_swizzle_v2r8(x+j_coord_offsetA,
1014                                               &jx0,&jy0,&jz0);
1015
1016             /* Calculate displacement vector */
1017             dx00             = _fjsp_sub_v2r8(ix0,jx0);
1018             dy00             = _fjsp_sub_v2r8(iy0,jy0);
1019             dz00             = _fjsp_sub_v2r8(iz0,jz0);
1020             dx10             = _fjsp_sub_v2r8(ix1,jx0);
1021             dy10             = _fjsp_sub_v2r8(iy1,jy0);
1022             dz10             = _fjsp_sub_v2r8(iz1,jz0);
1023             dx20             = _fjsp_sub_v2r8(ix2,jx0);
1024             dy20             = _fjsp_sub_v2r8(iy2,jy0);
1025             dz20             = _fjsp_sub_v2r8(iz2,jz0);
1026             dx30             = _fjsp_sub_v2r8(ix3,jx0);
1027             dy30             = _fjsp_sub_v2r8(iy3,jy0);
1028             dz30             = _fjsp_sub_v2r8(iz3,jz0);
1029
1030             /* Calculate squared distance and things based on it */
1031             rsq00            = gmx_fjsp_calc_rsq_v2r8(dx00,dy00,dz00);
1032             rsq10            = gmx_fjsp_calc_rsq_v2r8(dx10,dy10,dz10);
1033             rsq20            = gmx_fjsp_calc_rsq_v2r8(dx20,dy20,dz20);
1034             rsq30            = gmx_fjsp_calc_rsq_v2r8(dx30,dy30,dz30);
1035
1036             rinv00           = gmx_fjsp_invsqrt_v2r8(rsq00);
1037             rinv10           = gmx_fjsp_invsqrt_v2r8(rsq10);
1038             rinv20           = gmx_fjsp_invsqrt_v2r8(rsq20);
1039             rinv30           = gmx_fjsp_invsqrt_v2r8(rsq30);
1040
1041             rinvsq00         = _fjsp_mul_v2r8(rinv00,rinv00);
1042             rinvsq10         = _fjsp_mul_v2r8(rinv10,rinv10);
1043             rinvsq20         = _fjsp_mul_v2r8(rinv20,rinv20);
1044             rinvsq30         = _fjsp_mul_v2r8(rinv30,rinv30);
1045
1046             /* Load parameters for j particles */
1047             jq0              = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(),charge+jnrA+0);
1048             vdwjidx0A        = 2*vdwtype[jnrA+0];
1049
1050             fjx0             = _fjsp_setzero_v2r8();
1051             fjy0             = _fjsp_setzero_v2r8();
1052             fjz0             = _fjsp_setzero_v2r8();
1053
1054             /**************************
1055              * CALCULATE INTERACTIONS *
1056              **************************/
1057
1058             r00              = _fjsp_mul_v2r8(rsq00,rinv00);
1059
1060             /* Compute parameters for interactions between i and j atoms */
1061             gmx_fjsp_load_1pair_swizzle_v2r8(vdwparam+vdwioffset0+vdwjidx0A,&c6_00,&c12_00);
1062
1063             c6grid_00       = gmx_fjsp_load_1real_swizzle_v2r8(vdwgridparam+vdwioffset0+vdwjidx0A);
1064
1065             /* Analytical LJ-PME */
1066             rinvsix          = _fjsp_mul_v2r8(_fjsp_mul_v2r8(rinvsq00,rinvsq00),rinvsq00);
1067             ewcljrsq         = _fjsp_mul_v2r8(ewclj2,rsq00);
1068             ewclj6           = _fjsp_mul_v2r8(ewclj2,_fjsp_mul_v2r8(ewclj2,ewclj2));
1069             exponent         = gmx_simd_exp_d(-ewcljrsq);
1070             /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
1071             poly             = _fjsp_mul_v2r8(exponent,_fjsp_madd_v2r8(_fjsp_mul_v2r8(ewcljrsq,ewcljrsq),one_half,_fjsp_sub_v2r8(one,ewcljrsq)));
1072             /* f6A = 6 * C6grid * (1 - poly) */
1073             f6A              = _fjsp_mul_v2r8(c6grid_00,_fjsp_msub_v2r8(one,poly));
1074             /* f6B = C6grid * exponent * beta^6 */
1075             f6B              = _fjsp_mul_v2r8(_fjsp_mul_v2r8(c6grid_00,one_sixth),_fjsp_mul_v2r8(exponent,ewclj6));
1076             /* fvdw = 12*C12/r13 - ((6*C6 - f6A)/r6 + f6B)/r */
1077             fvdw              = _fjsp_mul_v2r8(_fjsp_madd_v2r8(_fjsp_msub_v2r8(c12_00,rinvsix,_fjsp_sub_v2r8(c6_00,f6A)),rinvsix,f6B),rinvsq00);
1078
1079             fscal            = fvdw;
1080
1081             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
1082
1083             /* Update vectorial force */
1084             fix0             = _fjsp_madd_v2r8(dx00,fscal,fix0);
1085             fiy0             = _fjsp_madd_v2r8(dy00,fscal,fiy0);
1086             fiz0             = _fjsp_madd_v2r8(dz00,fscal,fiz0);
1087             
1088             fjx0             = _fjsp_madd_v2r8(dx00,fscal,fjx0);
1089             fjy0             = _fjsp_madd_v2r8(dy00,fscal,fjy0);
1090             fjz0             = _fjsp_madd_v2r8(dz00,fscal,fjz0);
1091
1092             /**************************
1093              * CALCULATE INTERACTIONS *
1094              **************************/
1095
1096             r10              = _fjsp_mul_v2r8(rsq10,rinv10);
1097
1098             /* Compute parameters for interactions between i and j atoms */
1099             qq10             = _fjsp_mul_v2r8(iq1,jq0);
1100
1101             /* EWALD ELECTROSTATICS */
1102
1103             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
1104             ewrt             = _fjsp_mul_v2r8(r10,ewtabscale);
1105             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
1106             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
1107             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
1108
1109             gmx_fjsp_load_1pair_swizzle_v2r8(ewtab+ewconv.i[0],&ewtabF,&ewtabFn);
1110             felec            = _fjsp_madd_v2r8(eweps,ewtabFn,_fjsp_nmsub_v2r8(eweps,ewtabF,ewtabF));
1111             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq10,rinv10),_fjsp_sub_v2r8(rinvsq10,felec));
1112
1113             fscal            = felec;
1114
1115             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
1116
1117             /* Update vectorial force */
1118             fix1             = _fjsp_madd_v2r8(dx10,fscal,fix1);
1119             fiy1             = _fjsp_madd_v2r8(dy10,fscal,fiy1);
1120             fiz1             = _fjsp_madd_v2r8(dz10,fscal,fiz1);
1121             
1122             fjx0             = _fjsp_madd_v2r8(dx10,fscal,fjx0);
1123             fjy0             = _fjsp_madd_v2r8(dy10,fscal,fjy0);
1124             fjz0             = _fjsp_madd_v2r8(dz10,fscal,fjz0);
1125
1126             /**************************
1127              * CALCULATE INTERACTIONS *
1128              **************************/
1129
1130             r20              = _fjsp_mul_v2r8(rsq20,rinv20);
1131
1132             /* Compute parameters for interactions between i and j atoms */
1133             qq20             = _fjsp_mul_v2r8(iq2,jq0);
1134
1135             /* EWALD ELECTROSTATICS */
1136
1137             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
1138             ewrt             = _fjsp_mul_v2r8(r20,ewtabscale);
1139             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
1140             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
1141             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
1142
1143             gmx_fjsp_load_1pair_swizzle_v2r8(ewtab+ewconv.i[0],&ewtabF,&ewtabFn);
1144             felec            = _fjsp_madd_v2r8(eweps,ewtabFn,_fjsp_nmsub_v2r8(eweps,ewtabF,ewtabF));
1145             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq20,rinv20),_fjsp_sub_v2r8(rinvsq20,felec));
1146
1147             fscal            = felec;
1148
1149             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
1150
1151             /* Update vectorial force */
1152             fix2             = _fjsp_madd_v2r8(dx20,fscal,fix2);
1153             fiy2             = _fjsp_madd_v2r8(dy20,fscal,fiy2);
1154             fiz2             = _fjsp_madd_v2r8(dz20,fscal,fiz2);
1155             
1156             fjx0             = _fjsp_madd_v2r8(dx20,fscal,fjx0);
1157             fjy0             = _fjsp_madd_v2r8(dy20,fscal,fjy0);
1158             fjz0             = _fjsp_madd_v2r8(dz20,fscal,fjz0);
1159
1160             /**************************
1161              * CALCULATE INTERACTIONS *
1162              **************************/
1163
1164             r30              = _fjsp_mul_v2r8(rsq30,rinv30);
1165
1166             /* Compute parameters for interactions between i and j atoms */
1167             qq30             = _fjsp_mul_v2r8(iq3,jq0);
1168
1169             /* EWALD ELECTROSTATICS */
1170
1171             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
1172             ewrt             = _fjsp_mul_v2r8(r30,ewtabscale);
1173             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
1174             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
1175             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
1176
1177             gmx_fjsp_load_1pair_swizzle_v2r8(ewtab+ewconv.i[0],&ewtabF,&ewtabFn);
1178             felec            = _fjsp_madd_v2r8(eweps,ewtabFn,_fjsp_nmsub_v2r8(eweps,ewtabF,ewtabF));
1179             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq30,rinv30),_fjsp_sub_v2r8(rinvsq30,felec));
1180
1181             fscal            = felec;
1182
1183             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
1184
1185             /* Update vectorial force */
1186             fix3             = _fjsp_madd_v2r8(dx30,fscal,fix3);
1187             fiy3             = _fjsp_madd_v2r8(dy30,fscal,fiy3);
1188             fiz3             = _fjsp_madd_v2r8(dz30,fscal,fiz3);
1189             
1190             fjx0             = _fjsp_madd_v2r8(dx30,fscal,fjx0);
1191             fjy0             = _fjsp_madd_v2r8(dy30,fscal,fjy0);
1192             fjz0             = _fjsp_madd_v2r8(dz30,fscal,fjz0);
1193
1194             gmx_fjsp_decrement_1rvec_1ptr_swizzle_v2r8(f+j_coord_offsetA,fjx0,fjy0,fjz0);
1195
1196             /* Inner loop uses 168 flops */
1197         }
1198
1199         /* End of innermost loop */
1200
1201         gmx_fjsp_update_iforce_4atom_swizzle_v2r8(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
1202                                               f+i_coord_offset,fshift+i_shift_offset);
1203
1204         /* Increment number of inner iterations */
1205         inneriter                  += j_index_end - j_index_start;
1206
1207         /* Outer loop uses 24 flops */
1208     }
1209
1210     /* Increment number of outer iterations */
1211     outeriter        += nri;
1212
1213     /* Update outer/inner flops */
1214
1215     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W4_F,outeriter*24 + inneriter*168);
1216 }