Fujitsu Sparc64 acceleration and general fixes for non-x86 builds
[alexxy/gromacs.git] / src / gmxlib / nonbonded / nb_kernel_sparc64_hpc_ace_double / nb_kernel_ElecEwSw_VdwLJSw_GeomW4P1_sparc64_hpc_ace_double.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012, by the GROMACS development team, led by
5  * David van der Spoel, Berk Hess, Erik Lindahl, and including many
6  * others, as listed in the AUTHORS file in the top-level source
7  * 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 #ifdef HAVE_CONFIG_H
39 #include <config.h>
40 #endif
41
42 #include <math.h>
43
44 #include "../nb_kernel.h"
45 #include "types/simple.h"
46 #include "vec.h"
47 #include "nrnb.h"
48
49 #include "kernelutil_sparc64_hpc_ace_double.h"
50
51 /*
52  * Gromacs nonbonded kernel:   nb_kernel_ElecEwSw_VdwLJSw_GeomW4P1_VF_sparc64_hpc_ace_double
53  * Electrostatics interaction: Ewald
54  * VdW interaction:            LennardJones
55  * Geometry:                   Water4-Particle
56  * Calculate force/pot:        PotentialAndForce
57  */
58 void
59 nb_kernel_ElecEwSw_VdwLJSw_GeomW4P1_VF_sparc64_hpc_ace_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_restrict        kernel_data,
66                      t_nrnb * gmx_restrict                  nrnb)
67 {
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 double precision SIMD, e.g. for the two different
71      * jnr indices corresponding to data put in the four positions in the SIMD register.
72      */
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     _fjsp_v2r8       tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
81     int              vdwioffset0;
82     _fjsp_v2r8       ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
83     int              vdwioffset1;
84     _fjsp_v2r8       ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
85     int              vdwioffset2;
86     _fjsp_v2r8       ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
87     int              vdwioffset3;
88     _fjsp_v2r8       ix3,iy3,iz3,fix3,fiy3,fiz3,iq3,isai3;
89     int              vdwjidx0A,vdwjidx0B;
90     _fjsp_v2r8       jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
91     _fjsp_v2r8       dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
92     _fjsp_v2r8       dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
93     _fjsp_v2r8       dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
94     _fjsp_v2r8       dx30,dy30,dz30,rsq30,rinv30,rinvsq30,r30,qq30,c6_30,c12_30;
95     _fjsp_v2r8       velec,felec,velecsum,facel,crf,krf,krf2;
96     real             *charge;
97     int              nvdwtype;
98     _fjsp_v2r8       rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
99     int              *vdwtype;
100     real             *vdwparam;
101     _fjsp_v2r8       one_sixth   = gmx_fjsp_set1_v2r8(1.0/6.0);
102     _fjsp_v2r8       one_twelfth = gmx_fjsp_set1_v2r8(1.0/12.0);
103     _fjsp_v2r8       ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
104     real             *ewtab;
105     _fjsp_v2r8       rswitch,swV3,swV4,swV5,swF2,swF3,swF4,d,d2,sw,dsw;
106     real             rswitch_scalar,d_scalar;
107     _fjsp_v2r8       itab_tmp;
108     _fjsp_v2r8       dummy_mask,cutoff_mask;
109     _fjsp_v2r8       one     = gmx_fjsp_set1_v2r8(1.0);
110     _fjsp_v2r8       two     = gmx_fjsp_set1_v2r8(2.0);
111     union { _fjsp_v2r8 simd; long long int i[2]; } vfconv,gbconv,ewconv;
112
113     x                = xx[0];
114     f                = ff[0];
115
116     nri              = nlist->nri;
117     iinr             = nlist->iinr;
118     jindex           = nlist->jindex;
119     jjnr             = nlist->jjnr;
120     shiftidx         = nlist->shift;
121     gid              = nlist->gid;
122     shiftvec         = fr->shift_vec[0];
123     fshift           = fr->fshift[0];
124     facel            = gmx_fjsp_set1_v2r8(fr->epsfac);
125     charge           = mdatoms->chargeA;
126     nvdwtype         = fr->ntype;
127     vdwparam         = fr->nbfp;
128     vdwtype          = mdatoms->typeA;
129
130     sh_ewald         = gmx_fjsp_set1_v2r8(fr->ic->sh_ewald);
131     ewtab            = fr->ic->tabq_coul_FDV0;
132     ewtabscale       = gmx_fjsp_set1_v2r8(fr->ic->tabq_scale);
133     ewtabhalfspace   = gmx_fjsp_set1_v2r8(0.5/fr->ic->tabq_scale);
134
135     /* Setup water-specific parameters */
136     inr              = nlist->iinr[0];
137     iq1              = _fjsp_mul_v2r8(facel,gmx_fjsp_set1_v2r8(charge[inr+1]));
138     iq2              = _fjsp_mul_v2r8(facel,gmx_fjsp_set1_v2r8(charge[inr+2]));
139     iq3              = _fjsp_mul_v2r8(facel,gmx_fjsp_set1_v2r8(charge[inr+3]));
140     vdwioffset0      = 2*nvdwtype*vdwtype[inr+0];
141
142     /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
143     rcutoff_scalar   = fr->rcoulomb;
144     rcutoff          = gmx_fjsp_set1_v2r8(rcutoff_scalar);
145     rcutoff2         = _fjsp_mul_v2r8(rcutoff,rcutoff);
146
147     rswitch_scalar   = fr->rcoulomb_switch;
148     rswitch          = gmx_fjsp_set1_v2r8(rswitch_scalar);
149     /* Setup switch parameters */
150     d_scalar         = rcutoff_scalar-rswitch_scalar;
151     d                = gmx_fjsp_set1_v2r8(d_scalar);
152     swV3             = gmx_fjsp_set1_v2r8(-10.0/(d_scalar*d_scalar*d_scalar));
153     swV4             = gmx_fjsp_set1_v2r8( 15.0/(d_scalar*d_scalar*d_scalar*d_scalar));
154     swV5             = gmx_fjsp_set1_v2r8( -6.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
155     swF2             = gmx_fjsp_set1_v2r8(-30.0/(d_scalar*d_scalar*d_scalar));
156     swF3             = gmx_fjsp_set1_v2r8( 60.0/(d_scalar*d_scalar*d_scalar*d_scalar));
157     swF4             = gmx_fjsp_set1_v2r8(-30.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
158
159     /* Avoid stupid compiler warnings */
160     jnrA = jnrB = 0;
161     j_coord_offsetA = 0;
162     j_coord_offsetB = 0;
163
164     outeriter        = 0;
165     inneriter        = 0;
166
167     /* Start outer loop over neighborlists */
168     for(iidx=0; iidx<nri; iidx++)
169     {
170         /* Load shift vector for this list */
171         i_shift_offset   = DIM*shiftidx[iidx];
172
173         /* Load limits for loop over neighbors */
174         j_index_start    = jindex[iidx];
175         j_index_end      = jindex[iidx+1];
176
177         /* Get outer coordinate index */
178         inr              = iinr[iidx];
179         i_coord_offset   = DIM*inr;
180
181         /* Load i particle coords and add shift vector */
182         gmx_fjsp_load_shift_and_4rvec_broadcast_v2r8(shiftvec+i_shift_offset,x+i_coord_offset,
183                                                  &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
184
185         fix0             = _fjsp_setzero_v2r8();
186         fiy0             = _fjsp_setzero_v2r8();
187         fiz0             = _fjsp_setzero_v2r8();
188         fix1             = _fjsp_setzero_v2r8();
189         fiy1             = _fjsp_setzero_v2r8();
190         fiz1             = _fjsp_setzero_v2r8();
191         fix2             = _fjsp_setzero_v2r8();
192         fiy2             = _fjsp_setzero_v2r8();
193         fiz2             = _fjsp_setzero_v2r8();
194         fix3             = _fjsp_setzero_v2r8();
195         fiy3             = _fjsp_setzero_v2r8();
196         fiz3             = _fjsp_setzero_v2r8();
197
198         /* Reset potential sums */
199         velecsum         = _fjsp_setzero_v2r8();
200         vvdwsum          = _fjsp_setzero_v2r8();
201
202         /* Start inner kernel loop */
203         for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
204         {
205
206             /* Get j neighbor index, and coordinate index */
207             jnrA             = jjnr[jidx];
208             jnrB             = jjnr[jidx+1];
209             j_coord_offsetA  = DIM*jnrA;
210             j_coord_offsetB  = DIM*jnrB;
211
212             /* load j atom coordinates */
213             gmx_fjsp_load_1rvec_2ptr_swizzle_v2r8(x+j_coord_offsetA,x+j_coord_offsetB,
214                                               &jx0,&jy0,&jz0);
215
216             /* Calculate displacement vector */
217             dx00             = _fjsp_sub_v2r8(ix0,jx0);
218             dy00             = _fjsp_sub_v2r8(iy0,jy0);
219             dz00             = _fjsp_sub_v2r8(iz0,jz0);
220             dx10             = _fjsp_sub_v2r8(ix1,jx0);
221             dy10             = _fjsp_sub_v2r8(iy1,jy0);
222             dz10             = _fjsp_sub_v2r8(iz1,jz0);
223             dx20             = _fjsp_sub_v2r8(ix2,jx0);
224             dy20             = _fjsp_sub_v2r8(iy2,jy0);
225             dz20             = _fjsp_sub_v2r8(iz2,jz0);
226             dx30             = _fjsp_sub_v2r8(ix3,jx0);
227             dy30             = _fjsp_sub_v2r8(iy3,jy0);
228             dz30             = _fjsp_sub_v2r8(iz3,jz0);
229
230             /* Calculate squared distance and things based on it */
231             rsq00            = gmx_fjsp_calc_rsq_v2r8(dx00,dy00,dz00);
232             rsq10            = gmx_fjsp_calc_rsq_v2r8(dx10,dy10,dz10);
233             rsq20            = gmx_fjsp_calc_rsq_v2r8(dx20,dy20,dz20);
234             rsq30            = gmx_fjsp_calc_rsq_v2r8(dx30,dy30,dz30);
235
236             rinv00           = gmx_fjsp_invsqrt_v2r8(rsq00);
237             rinv10           = gmx_fjsp_invsqrt_v2r8(rsq10);
238             rinv20           = gmx_fjsp_invsqrt_v2r8(rsq20);
239             rinv30           = gmx_fjsp_invsqrt_v2r8(rsq30);
240
241             rinvsq00         = _fjsp_mul_v2r8(rinv00,rinv00);
242             rinvsq10         = _fjsp_mul_v2r8(rinv10,rinv10);
243             rinvsq20         = _fjsp_mul_v2r8(rinv20,rinv20);
244             rinvsq30         = _fjsp_mul_v2r8(rinv30,rinv30);
245
246             /* Load parameters for j particles */
247             jq0              = gmx_fjsp_load_2real_swizzle_v2r8(charge+jnrA+0,charge+jnrB+0);
248             vdwjidx0A        = 2*vdwtype[jnrA+0];
249             vdwjidx0B        = 2*vdwtype[jnrB+0];
250
251             fjx0             = _fjsp_setzero_v2r8();
252             fjy0             = _fjsp_setzero_v2r8();
253             fjz0             = _fjsp_setzero_v2r8();
254
255             /**************************
256              * CALCULATE INTERACTIONS *
257              **************************/
258
259             if (gmx_fjsp_any_lt_v2r8(rsq00,rcutoff2))
260             {
261
262             r00              = _fjsp_mul_v2r8(rsq00,rinv00);
263
264             /* Compute parameters for interactions between i and j atoms */
265             gmx_fjsp_load_2pair_swizzle_v2r8(vdwparam+vdwioffset0+vdwjidx0A,
266                                          vdwparam+vdwioffset0+vdwjidx0B,&c6_00,&c12_00);
267
268             /* LENNARD-JONES DISPERSION/REPULSION */
269
270             rinvsix          = _fjsp_mul_v2r8(_fjsp_mul_v2r8(rinvsq00,rinvsq00),rinvsq00);
271             vvdw6            = _fjsp_mul_v2r8(c6_00,rinvsix);
272             vvdw12           = _fjsp_mul_v2r8(c12_00,_fjsp_mul_v2r8(rinvsix,rinvsix));
273             vvdw             = _fjsp_msub_v2r8( vvdw12,one_twelfth, _fjsp_mul_v2r8(vvdw6,one_sixth) );
274             fvdw             = _fjsp_mul_v2r8(_fjsp_sub_v2r8(vvdw12,vvdw6),rinvsq00);
275
276             d                = _fjsp_sub_v2r8(r00,rswitch);
277             d                = _fjsp_max_v2r8(d,_fjsp_setzero_v2r8());
278             d2               = _fjsp_mul_v2r8(d,d);
279             sw               = _fjsp_add_v2r8(one,_fjsp_mul_v2r8(d2,_fjsp_mul_v2r8(d,_fjsp_madd_v2r8(d,_fjsp_madd_v2r8(d,swV5,swV4),swV3))));
280
281             dsw              = _fjsp_mul_v2r8(d2,_fjsp_madd_v2r8(d,_fjsp_madd_v2r8(d,swF4,swF3),swF2));
282
283             /* Evaluate switch function */
284             /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
285             fvdw             = _fjsp_msub_v2r8( fvdw,sw , _fjsp_mul_v2r8(rinv00,_fjsp_mul_v2r8(vvdw,dsw)) );
286             vvdw             = _fjsp_mul_v2r8(vvdw,sw);
287             cutoff_mask      = _fjsp_cmplt_v2r8(rsq00,rcutoff2);
288
289             /* Update potential sum for this i atom from the interaction with this j atom. */
290             vvdw             = _fjsp_and_v2r8(vvdw,cutoff_mask);
291             vvdwsum          = _fjsp_add_v2r8(vvdwsum,vvdw);
292
293             fscal            = fvdw;
294
295             fscal            = _fjsp_and_v2r8(fscal,cutoff_mask);
296
297             /* Update vectorial force */
298             fix0             = _fjsp_madd_v2r8(dx00,fscal,fix0);
299             fiy0             = _fjsp_madd_v2r8(dy00,fscal,fiy0);
300             fiz0             = _fjsp_madd_v2r8(dz00,fscal,fiz0);
301             
302             fjx0             = _fjsp_madd_v2r8(dx00,fscal,fjx0);
303             fjy0             = _fjsp_madd_v2r8(dy00,fscal,fjy0);
304             fjz0             = _fjsp_madd_v2r8(dz00,fscal,fjz0);
305
306             }
307
308             /**************************
309              * CALCULATE INTERACTIONS *
310              **************************/
311
312             if (gmx_fjsp_any_lt_v2r8(rsq10,rcutoff2))
313             {
314
315             r10              = _fjsp_mul_v2r8(rsq10,rinv10);
316
317             /* Compute parameters for interactions between i and j atoms */
318             qq10             = _fjsp_mul_v2r8(iq1,jq0);
319
320             /* EWALD ELECTROSTATICS */
321
322             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
323             ewrt             = _fjsp_mul_v2r8(r10,ewtabscale);
324             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
325             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
326             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
327
328             ewtabF           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[0] );
329             ewtabD           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[1] );
330             GMX_FJSP_TRANSPOSE2_V2R8(ewtabF,ewtabD);
331             ewtabV           = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[0] +2);
332             ewtabFn          = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[1] +2);
333             GMX_FJSP_TRANSPOSE2_V2R8(ewtabV,ewtabFn);
334             felec            = _fjsp_madd_v2r8(eweps,ewtabD,ewtabF);
335             velec            = _fjsp_nmsub_v2r8(_fjsp_mul_v2r8(ewtabhalfspace,eweps) ,_fjsp_add_v2r8(ewtabF,felec), ewtabV);
336             velec            = _fjsp_mul_v2r8(qq10,_fjsp_sub_v2r8(rinv10,velec));
337             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq10,rinv10),_fjsp_sub_v2r8(rinvsq10,felec));
338
339             d                = _fjsp_sub_v2r8(r10,rswitch);
340             d                = _fjsp_max_v2r8(d,_fjsp_setzero_v2r8());
341             d2               = _fjsp_mul_v2r8(d,d);
342             sw               = _fjsp_add_v2r8(one,_fjsp_mul_v2r8(d2,_fjsp_mul_v2r8(d,_fjsp_madd_v2r8(d,_fjsp_madd_v2r8(d,swV5,swV4),swV3))));
343
344             dsw              = _fjsp_mul_v2r8(d2,_fjsp_madd_v2r8(d,_fjsp_madd_v2r8(d,swF4,swF3),swF2));
345
346             /* Evaluate switch function */
347             /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
348             felec            = _fjsp_msub_v2r8( felec,sw , _fjsp_mul_v2r8(rinv10,_fjsp_mul_v2r8(velec,dsw)) );
349             velec            = _fjsp_mul_v2r8(velec,sw);
350             cutoff_mask      = _fjsp_cmplt_v2r8(rsq10,rcutoff2);
351
352             /* Update potential sum for this i atom from the interaction with this j atom. */
353             velec            = _fjsp_and_v2r8(velec,cutoff_mask);
354             velecsum         = _fjsp_add_v2r8(velecsum,velec);
355
356             fscal            = felec;
357
358             fscal            = _fjsp_and_v2r8(fscal,cutoff_mask);
359
360             /* Update vectorial force */
361             fix1             = _fjsp_madd_v2r8(dx10,fscal,fix1);
362             fiy1             = _fjsp_madd_v2r8(dy10,fscal,fiy1);
363             fiz1             = _fjsp_madd_v2r8(dz10,fscal,fiz1);
364             
365             fjx0             = _fjsp_madd_v2r8(dx10,fscal,fjx0);
366             fjy0             = _fjsp_madd_v2r8(dy10,fscal,fjy0);
367             fjz0             = _fjsp_madd_v2r8(dz10,fscal,fjz0);
368
369             }
370
371             /**************************
372              * CALCULATE INTERACTIONS *
373              **************************/
374
375             if (gmx_fjsp_any_lt_v2r8(rsq20,rcutoff2))
376             {
377
378             r20              = _fjsp_mul_v2r8(rsq20,rinv20);
379
380             /* Compute parameters for interactions between i and j atoms */
381             qq20             = _fjsp_mul_v2r8(iq2,jq0);
382
383             /* EWALD ELECTROSTATICS */
384
385             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
386             ewrt             = _fjsp_mul_v2r8(r20,ewtabscale);
387             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
388             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
389             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
390
391             ewtabF           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[0] );
392             ewtabD           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[1] );
393             GMX_FJSP_TRANSPOSE2_V2R8(ewtabF,ewtabD);
394             ewtabV           = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[0] +2);
395             ewtabFn          = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[1] +2);
396             GMX_FJSP_TRANSPOSE2_V2R8(ewtabV,ewtabFn);
397             felec            = _fjsp_madd_v2r8(eweps,ewtabD,ewtabF);
398             velec            = _fjsp_nmsub_v2r8(_fjsp_mul_v2r8(ewtabhalfspace,eweps) ,_fjsp_add_v2r8(ewtabF,felec), ewtabV);
399             velec            = _fjsp_mul_v2r8(qq20,_fjsp_sub_v2r8(rinv20,velec));
400             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq20,rinv20),_fjsp_sub_v2r8(rinvsq20,felec));
401
402             d                = _fjsp_sub_v2r8(r20,rswitch);
403             d                = _fjsp_max_v2r8(d,_fjsp_setzero_v2r8());
404             d2               = _fjsp_mul_v2r8(d,d);
405             sw               = _fjsp_add_v2r8(one,_fjsp_mul_v2r8(d2,_fjsp_mul_v2r8(d,_fjsp_madd_v2r8(d,_fjsp_madd_v2r8(d,swV5,swV4),swV3))));
406
407             dsw              = _fjsp_mul_v2r8(d2,_fjsp_madd_v2r8(d,_fjsp_madd_v2r8(d,swF4,swF3),swF2));
408
409             /* Evaluate switch function */
410             /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
411             felec            = _fjsp_msub_v2r8( felec,sw , _fjsp_mul_v2r8(rinv20,_fjsp_mul_v2r8(velec,dsw)) );
412             velec            = _fjsp_mul_v2r8(velec,sw);
413             cutoff_mask      = _fjsp_cmplt_v2r8(rsq20,rcutoff2);
414
415             /* Update potential sum for this i atom from the interaction with this j atom. */
416             velec            = _fjsp_and_v2r8(velec,cutoff_mask);
417             velecsum         = _fjsp_add_v2r8(velecsum,velec);
418
419             fscal            = felec;
420
421             fscal            = _fjsp_and_v2r8(fscal,cutoff_mask);
422
423             /* Update vectorial force */
424             fix2             = _fjsp_madd_v2r8(dx20,fscal,fix2);
425             fiy2             = _fjsp_madd_v2r8(dy20,fscal,fiy2);
426             fiz2             = _fjsp_madd_v2r8(dz20,fscal,fiz2);
427             
428             fjx0             = _fjsp_madd_v2r8(dx20,fscal,fjx0);
429             fjy0             = _fjsp_madd_v2r8(dy20,fscal,fjy0);
430             fjz0             = _fjsp_madd_v2r8(dz20,fscal,fjz0);
431
432             }
433
434             /**************************
435              * CALCULATE INTERACTIONS *
436              **************************/
437
438             if (gmx_fjsp_any_lt_v2r8(rsq30,rcutoff2))
439             {
440
441             r30              = _fjsp_mul_v2r8(rsq30,rinv30);
442
443             /* Compute parameters for interactions between i and j atoms */
444             qq30             = _fjsp_mul_v2r8(iq3,jq0);
445
446             /* EWALD ELECTROSTATICS */
447
448             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
449             ewrt             = _fjsp_mul_v2r8(r30,ewtabscale);
450             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
451             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
452             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
453
454             ewtabF           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[0] );
455             ewtabD           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[1] );
456             GMX_FJSP_TRANSPOSE2_V2R8(ewtabF,ewtabD);
457             ewtabV           = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[0] +2);
458             ewtabFn          = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[1] +2);
459             GMX_FJSP_TRANSPOSE2_V2R8(ewtabV,ewtabFn);
460             felec            = _fjsp_madd_v2r8(eweps,ewtabD,ewtabF);
461             velec            = _fjsp_nmsub_v2r8(_fjsp_mul_v2r8(ewtabhalfspace,eweps) ,_fjsp_add_v2r8(ewtabF,felec), ewtabV);
462             velec            = _fjsp_mul_v2r8(qq30,_fjsp_sub_v2r8(rinv30,velec));
463             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq30,rinv30),_fjsp_sub_v2r8(rinvsq30,felec));
464
465             d                = _fjsp_sub_v2r8(r30,rswitch);
466             d                = _fjsp_max_v2r8(d,_fjsp_setzero_v2r8());
467             d2               = _fjsp_mul_v2r8(d,d);
468             sw               = _fjsp_add_v2r8(one,_fjsp_mul_v2r8(d2,_fjsp_mul_v2r8(d,_fjsp_madd_v2r8(d,_fjsp_madd_v2r8(d,swV5,swV4),swV3))));
469
470             dsw              = _fjsp_mul_v2r8(d2,_fjsp_madd_v2r8(d,_fjsp_madd_v2r8(d,swF4,swF3),swF2));
471
472             /* Evaluate switch function */
473             /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
474             felec            = _fjsp_msub_v2r8( felec,sw , _fjsp_mul_v2r8(rinv30,_fjsp_mul_v2r8(velec,dsw)) );
475             velec            = _fjsp_mul_v2r8(velec,sw);
476             cutoff_mask      = _fjsp_cmplt_v2r8(rsq30,rcutoff2);
477
478             /* Update potential sum for this i atom from the interaction with this j atom. */
479             velec            = _fjsp_and_v2r8(velec,cutoff_mask);
480             velecsum         = _fjsp_add_v2r8(velecsum,velec);
481
482             fscal            = felec;
483
484             fscal            = _fjsp_and_v2r8(fscal,cutoff_mask);
485
486             /* Update vectorial force */
487             fix3             = _fjsp_madd_v2r8(dx30,fscal,fix3);
488             fiy3             = _fjsp_madd_v2r8(dy30,fscal,fiy3);
489             fiz3             = _fjsp_madd_v2r8(dz30,fscal,fiz3);
490             
491             fjx0             = _fjsp_madd_v2r8(dx30,fscal,fjx0);
492             fjy0             = _fjsp_madd_v2r8(dy30,fscal,fjy0);
493             fjz0             = _fjsp_madd_v2r8(dz30,fscal,fjz0);
494
495             }
496
497             gmx_fjsp_decrement_1rvec_2ptr_swizzle_v2r8(f+j_coord_offsetA,f+j_coord_offsetB,fjx0,fjy0,fjz0);
498
499             /* Inner loop uses 269 flops */
500         }
501
502         if(jidx<j_index_end)
503         {
504
505             jnrA             = jjnr[jidx];
506             j_coord_offsetA  = DIM*jnrA;
507
508             /* load j atom coordinates */
509             gmx_fjsp_load_1rvec_1ptr_swizzle_v2r8(x+j_coord_offsetA,
510                                               &jx0,&jy0,&jz0);
511
512             /* Calculate displacement vector */
513             dx00             = _fjsp_sub_v2r8(ix0,jx0);
514             dy00             = _fjsp_sub_v2r8(iy0,jy0);
515             dz00             = _fjsp_sub_v2r8(iz0,jz0);
516             dx10             = _fjsp_sub_v2r8(ix1,jx0);
517             dy10             = _fjsp_sub_v2r8(iy1,jy0);
518             dz10             = _fjsp_sub_v2r8(iz1,jz0);
519             dx20             = _fjsp_sub_v2r8(ix2,jx0);
520             dy20             = _fjsp_sub_v2r8(iy2,jy0);
521             dz20             = _fjsp_sub_v2r8(iz2,jz0);
522             dx30             = _fjsp_sub_v2r8(ix3,jx0);
523             dy30             = _fjsp_sub_v2r8(iy3,jy0);
524             dz30             = _fjsp_sub_v2r8(iz3,jz0);
525
526             /* Calculate squared distance and things based on it */
527             rsq00            = gmx_fjsp_calc_rsq_v2r8(dx00,dy00,dz00);
528             rsq10            = gmx_fjsp_calc_rsq_v2r8(dx10,dy10,dz10);
529             rsq20            = gmx_fjsp_calc_rsq_v2r8(dx20,dy20,dz20);
530             rsq30            = gmx_fjsp_calc_rsq_v2r8(dx30,dy30,dz30);
531
532             rinv00           = gmx_fjsp_invsqrt_v2r8(rsq00);
533             rinv10           = gmx_fjsp_invsqrt_v2r8(rsq10);
534             rinv20           = gmx_fjsp_invsqrt_v2r8(rsq20);
535             rinv30           = gmx_fjsp_invsqrt_v2r8(rsq30);
536
537             rinvsq00         = _fjsp_mul_v2r8(rinv00,rinv00);
538             rinvsq10         = _fjsp_mul_v2r8(rinv10,rinv10);
539             rinvsq20         = _fjsp_mul_v2r8(rinv20,rinv20);
540             rinvsq30         = _fjsp_mul_v2r8(rinv30,rinv30);
541
542             /* Load parameters for j particles */
543             jq0              = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(),charge+jnrA+0);
544             vdwjidx0A        = 2*vdwtype[jnrA+0];
545
546             fjx0             = _fjsp_setzero_v2r8();
547             fjy0             = _fjsp_setzero_v2r8();
548             fjz0             = _fjsp_setzero_v2r8();
549
550             /**************************
551              * CALCULATE INTERACTIONS *
552              **************************/
553
554             if (gmx_fjsp_any_lt_v2r8(rsq00,rcutoff2))
555             {
556
557             r00              = _fjsp_mul_v2r8(rsq00,rinv00);
558
559             /* Compute parameters for interactions between i and j atoms */
560             gmx_fjsp_load_1pair_swizzle_v2r8(vdwparam+vdwioffset0+vdwjidx0A,&c6_00,&c12_00);
561
562             /* LENNARD-JONES DISPERSION/REPULSION */
563
564             rinvsix          = _fjsp_mul_v2r8(_fjsp_mul_v2r8(rinvsq00,rinvsq00),rinvsq00);
565             vvdw6            = _fjsp_mul_v2r8(c6_00,rinvsix);
566             vvdw12           = _fjsp_mul_v2r8(c12_00,_fjsp_mul_v2r8(rinvsix,rinvsix));
567             vvdw             = _fjsp_msub_v2r8( vvdw12,one_twelfth, _fjsp_mul_v2r8(vvdw6,one_sixth) );
568             fvdw             = _fjsp_mul_v2r8(_fjsp_sub_v2r8(vvdw12,vvdw6),rinvsq00);
569
570             d                = _fjsp_sub_v2r8(r00,rswitch);
571             d                = _fjsp_max_v2r8(d,_fjsp_setzero_v2r8());
572             d2               = _fjsp_mul_v2r8(d,d);
573             sw               = _fjsp_add_v2r8(one,_fjsp_mul_v2r8(d2,_fjsp_mul_v2r8(d,_fjsp_madd_v2r8(d,_fjsp_madd_v2r8(d,swV5,swV4),swV3))));
574
575             dsw              = _fjsp_mul_v2r8(d2,_fjsp_madd_v2r8(d,_fjsp_madd_v2r8(d,swF4,swF3),swF2));
576
577             /* Evaluate switch function */
578             /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
579             fvdw             = _fjsp_msub_v2r8( fvdw,sw , _fjsp_mul_v2r8(rinv00,_fjsp_mul_v2r8(vvdw,dsw)) );
580             vvdw             = _fjsp_mul_v2r8(vvdw,sw);
581             cutoff_mask      = _fjsp_cmplt_v2r8(rsq00,rcutoff2);
582
583             /* Update potential sum for this i atom from the interaction with this j atom. */
584             vvdw             = _fjsp_and_v2r8(vvdw,cutoff_mask);
585             vvdw             = _fjsp_unpacklo_v2r8(vvdw,_fjsp_setzero_v2r8());
586             vvdwsum          = _fjsp_add_v2r8(vvdwsum,vvdw);
587
588             fscal            = fvdw;
589
590             fscal            = _fjsp_and_v2r8(fscal,cutoff_mask);
591
592             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
593
594             /* Update vectorial force */
595             fix0             = _fjsp_madd_v2r8(dx00,fscal,fix0);
596             fiy0             = _fjsp_madd_v2r8(dy00,fscal,fiy0);
597             fiz0             = _fjsp_madd_v2r8(dz00,fscal,fiz0);
598             
599             fjx0             = _fjsp_madd_v2r8(dx00,fscal,fjx0);
600             fjy0             = _fjsp_madd_v2r8(dy00,fscal,fjy0);
601             fjz0             = _fjsp_madd_v2r8(dz00,fscal,fjz0);
602
603             }
604
605             /**************************
606              * CALCULATE INTERACTIONS *
607              **************************/
608
609             if (gmx_fjsp_any_lt_v2r8(rsq10,rcutoff2))
610             {
611
612             r10              = _fjsp_mul_v2r8(rsq10,rinv10);
613
614             /* Compute parameters for interactions between i and j atoms */
615             qq10             = _fjsp_mul_v2r8(iq1,jq0);
616
617             /* EWALD ELECTROSTATICS */
618
619             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
620             ewrt             = _fjsp_mul_v2r8(r10,ewtabscale);
621             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
622             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
623             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
624
625             ewtabF           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[0] );
626             ewtabD           = _fjsp_setzero_v2r8();
627             GMX_FJSP_TRANSPOSE2_V2R8(ewtabF,ewtabD);
628             ewtabV           = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[0] +2);
629             ewtabFn          = _fjsp_setzero_v2r8();
630             GMX_FJSP_TRANSPOSE2_V2R8(ewtabV,ewtabFn);
631             felec            = _fjsp_madd_v2r8(eweps,ewtabD,ewtabF);
632             velec            = _fjsp_nmsub_v2r8(_fjsp_mul_v2r8(ewtabhalfspace,eweps) ,_fjsp_add_v2r8(ewtabF,felec), ewtabV);
633             velec            = _fjsp_mul_v2r8(qq10,_fjsp_sub_v2r8(rinv10,velec));
634             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq10,rinv10),_fjsp_sub_v2r8(rinvsq10,felec));
635
636             d                = _fjsp_sub_v2r8(r10,rswitch);
637             d                = _fjsp_max_v2r8(d,_fjsp_setzero_v2r8());
638             d2               = _fjsp_mul_v2r8(d,d);
639             sw               = _fjsp_add_v2r8(one,_fjsp_mul_v2r8(d2,_fjsp_mul_v2r8(d,_fjsp_madd_v2r8(d,_fjsp_madd_v2r8(d,swV5,swV4),swV3))));
640
641             dsw              = _fjsp_mul_v2r8(d2,_fjsp_madd_v2r8(d,_fjsp_madd_v2r8(d,swF4,swF3),swF2));
642
643             /* Evaluate switch function */
644             /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
645             felec            = _fjsp_msub_v2r8( felec,sw , _fjsp_mul_v2r8(rinv10,_fjsp_mul_v2r8(velec,dsw)) );
646             velec            = _fjsp_mul_v2r8(velec,sw);
647             cutoff_mask      = _fjsp_cmplt_v2r8(rsq10,rcutoff2);
648
649             /* Update potential sum for this i atom from the interaction with this j atom. */
650             velec            = _fjsp_and_v2r8(velec,cutoff_mask);
651             velec            = _fjsp_unpacklo_v2r8(velec,_fjsp_setzero_v2r8());
652             velecsum         = _fjsp_add_v2r8(velecsum,velec);
653
654             fscal            = felec;
655
656             fscal            = _fjsp_and_v2r8(fscal,cutoff_mask);
657
658             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
659
660             /* Update vectorial force */
661             fix1             = _fjsp_madd_v2r8(dx10,fscal,fix1);
662             fiy1             = _fjsp_madd_v2r8(dy10,fscal,fiy1);
663             fiz1             = _fjsp_madd_v2r8(dz10,fscal,fiz1);
664             
665             fjx0             = _fjsp_madd_v2r8(dx10,fscal,fjx0);
666             fjy0             = _fjsp_madd_v2r8(dy10,fscal,fjy0);
667             fjz0             = _fjsp_madd_v2r8(dz10,fscal,fjz0);
668
669             }
670
671             /**************************
672              * CALCULATE INTERACTIONS *
673              **************************/
674
675             if (gmx_fjsp_any_lt_v2r8(rsq20,rcutoff2))
676             {
677
678             r20              = _fjsp_mul_v2r8(rsq20,rinv20);
679
680             /* Compute parameters for interactions between i and j atoms */
681             qq20             = _fjsp_mul_v2r8(iq2,jq0);
682
683             /* EWALD ELECTROSTATICS */
684
685             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
686             ewrt             = _fjsp_mul_v2r8(r20,ewtabscale);
687             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
688             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
689             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
690
691             ewtabF           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[0] );
692             ewtabD           = _fjsp_setzero_v2r8();
693             GMX_FJSP_TRANSPOSE2_V2R8(ewtabF,ewtabD);
694             ewtabV           = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[0] +2);
695             ewtabFn          = _fjsp_setzero_v2r8();
696             GMX_FJSP_TRANSPOSE2_V2R8(ewtabV,ewtabFn);
697             felec            = _fjsp_madd_v2r8(eweps,ewtabD,ewtabF);
698             velec            = _fjsp_nmsub_v2r8(_fjsp_mul_v2r8(ewtabhalfspace,eweps) ,_fjsp_add_v2r8(ewtabF,felec), ewtabV);
699             velec            = _fjsp_mul_v2r8(qq20,_fjsp_sub_v2r8(rinv20,velec));
700             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq20,rinv20),_fjsp_sub_v2r8(rinvsq20,felec));
701
702             d                = _fjsp_sub_v2r8(r20,rswitch);
703             d                = _fjsp_max_v2r8(d,_fjsp_setzero_v2r8());
704             d2               = _fjsp_mul_v2r8(d,d);
705             sw               = _fjsp_add_v2r8(one,_fjsp_mul_v2r8(d2,_fjsp_mul_v2r8(d,_fjsp_madd_v2r8(d,_fjsp_madd_v2r8(d,swV5,swV4),swV3))));
706
707             dsw              = _fjsp_mul_v2r8(d2,_fjsp_madd_v2r8(d,_fjsp_madd_v2r8(d,swF4,swF3),swF2));
708
709             /* Evaluate switch function */
710             /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
711             felec            = _fjsp_msub_v2r8( felec,sw , _fjsp_mul_v2r8(rinv20,_fjsp_mul_v2r8(velec,dsw)) );
712             velec            = _fjsp_mul_v2r8(velec,sw);
713             cutoff_mask      = _fjsp_cmplt_v2r8(rsq20,rcutoff2);
714
715             /* Update potential sum for this i atom from the interaction with this j atom. */
716             velec            = _fjsp_and_v2r8(velec,cutoff_mask);
717             velec            = _fjsp_unpacklo_v2r8(velec,_fjsp_setzero_v2r8());
718             velecsum         = _fjsp_add_v2r8(velecsum,velec);
719
720             fscal            = felec;
721
722             fscal            = _fjsp_and_v2r8(fscal,cutoff_mask);
723
724             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
725
726             /* Update vectorial force */
727             fix2             = _fjsp_madd_v2r8(dx20,fscal,fix2);
728             fiy2             = _fjsp_madd_v2r8(dy20,fscal,fiy2);
729             fiz2             = _fjsp_madd_v2r8(dz20,fscal,fiz2);
730             
731             fjx0             = _fjsp_madd_v2r8(dx20,fscal,fjx0);
732             fjy0             = _fjsp_madd_v2r8(dy20,fscal,fjy0);
733             fjz0             = _fjsp_madd_v2r8(dz20,fscal,fjz0);
734
735             }
736
737             /**************************
738              * CALCULATE INTERACTIONS *
739              **************************/
740
741             if (gmx_fjsp_any_lt_v2r8(rsq30,rcutoff2))
742             {
743
744             r30              = _fjsp_mul_v2r8(rsq30,rinv30);
745
746             /* Compute parameters for interactions between i and j atoms */
747             qq30             = _fjsp_mul_v2r8(iq3,jq0);
748
749             /* EWALD ELECTROSTATICS */
750
751             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
752             ewrt             = _fjsp_mul_v2r8(r30,ewtabscale);
753             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
754             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
755             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
756
757             ewtabF           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[0] );
758             ewtabD           = _fjsp_setzero_v2r8();
759             GMX_FJSP_TRANSPOSE2_V2R8(ewtabF,ewtabD);
760             ewtabV           = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[0] +2);
761             ewtabFn          = _fjsp_setzero_v2r8();
762             GMX_FJSP_TRANSPOSE2_V2R8(ewtabV,ewtabFn);
763             felec            = _fjsp_madd_v2r8(eweps,ewtabD,ewtabF);
764             velec            = _fjsp_nmsub_v2r8(_fjsp_mul_v2r8(ewtabhalfspace,eweps) ,_fjsp_add_v2r8(ewtabF,felec), ewtabV);
765             velec            = _fjsp_mul_v2r8(qq30,_fjsp_sub_v2r8(rinv30,velec));
766             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq30,rinv30),_fjsp_sub_v2r8(rinvsq30,felec));
767
768             d                = _fjsp_sub_v2r8(r30,rswitch);
769             d                = _fjsp_max_v2r8(d,_fjsp_setzero_v2r8());
770             d2               = _fjsp_mul_v2r8(d,d);
771             sw               = _fjsp_add_v2r8(one,_fjsp_mul_v2r8(d2,_fjsp_mul_v2r8(d,_fjsp_madd_v2r8(d,_fjsp_madd_v2r8(d,swV5,swV4),swV3))));
772
773             dsw              = _fjsp_mul_v2r8(d2,_fjsp_madd_v2r8(d,_fjsp_madd_v2r8(d,swF4,swF3),swF2));
774
775             /* Evaluate switch function */
776             /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
777             felec            = _fjsp_msub_v2r8( felec,sw , _fjsp_mul_v2r8(rinv30,_fjsp_mul_v2r8(velec,dsw)) );
778             velec            = _fjsp_mul_v2r8(velec,sw);
779             cutoff_mask      = _fjsp_cmplt_v2r8(rsq30,rcutoff2);
780
781             /* Update potential sum for this i atom from the interaction with this j atom. */
782             velec            = _fjsp_and_v2r8(velec,cutoff_mask);
783             velec            = _fjsp_unpacklo_v2r8(velec,_fjsp_setzero_v2r8());
784             velecsum         = _fjsp_add_v2r8(velecsum,velec);
785
786             fscal            = felec;
787
788             fscal            = _fjsp_and_v2r8(fscal,cutoff_mask);
789
790             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
791
792             /* Update vectorial force */
793             fix3             = _fjsp_madd_v2r8(dx30,fscal,fix3);
794             fiy3             = _fjsp_madd_v2r8(dy30,fscal,fiy3);
795             fiz3             = _fjsp_madd_v2r8(dz30,fscal,fiz3);
796             
797             fjx0             = _fjsp_madd_v2r8(dx30,fscal,fjx0);
798             fjy0             = _fjsp_madd_v2r8(dy30,fscal,fjy0);
799             fjz0             = _fjsp_madd_v2r8(dz30,fscal,fjz0);
800
801             }
802
803             gmx_fjsp_decrement_1rvec_1ptr_swizzle_v2r8(f+j_coord_offsetA,fjx0,fjy0,fjz0);
804
805             /* Inner loop uses 269 flops */
806         }
807
808         /* End of innermost loop */
809
810         gmx_fjsp_update_iforce_4atom_swizzle_v2r8(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
811                                               f+i_coord_offset,fshift+i_shift_offset);
812
813         ggid                        = gid[iidx];
814         /* Update potential energies */
815         gmx_fjsp_update_1pot_v2r8(velecsum,kernel_data->energygrp_elec+ggid);
816         gmx_fjsp_update_1pot_v2r8(vvdwsum,kernel_data->energygrp_vdw+ggid);
817
818         /* Increment number of inner iterations */
819         inneriter                  += j_index_end - j_index_start;
820
821         /* Outer loop uses 26 flops */
822     }
823
824     /* Increment number of outer iterations */
825     outeriter        += nri;
826
827     /* Update outer/inner flops */
828
829     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W4_VF,outeriter*26 + inneriter*269);
830 }
831 /*
832  * Gromacs nonbonded kernel:   nb_kernel_ElecEwSw_VdwLJSw_GeomW4P1_F_sparc64_hpc_ace_double
833  * Electrostatics interaction: Ewald
834  * VdW interaction:            LennardJones
835  * Geometry:                   Water4-Particle
836  * Calculate force/pot:        Force
837  */
838 void
839 nb_kernel_ElecEwSw_VdwLJSw_GeomW4P1_F_sparc64_hpc_ace_double
840                     (t_nblist * gmx_restrict                nlist,
841                      rvec * gmx_restrict                    xx,
842                      rvec * gmx_restrict                    ff,
843                      t_forcerec * gmx_restrict              fr,
844                      t_mdatoms * gmx_restrict               mdatoms,
845                      nb_kernel_data_t * gmx_restrict        kernel_data,
846                      t_nrnb * gmx_restrict                  nrnb)
847 {
848     /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
849      * just 0 for non-waters.
850      * Suffixes A,B refer to j loop unrolling done with double precision SIMD, e.g. for the two different
851      * jnr indices corresponding to data put in the four positions in the SIMD register.
852      */
853     int              i_shift_offset,i_coord_offset,outeriter,inneriter;
854     int              j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
855     int              jnrA,jnrB;
856     int              j_coord_offsetA,j_coord_offsetB;
857     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
858     real             rcutoff_scalar;
859     real             *shiftvec,*fshift,*x,*f;
860     _fjsp_v2r8       tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
861     int              vdwioffset0;
862     _fjsp_v2r8       ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
863     int              vdwioffset1;
864     _fjsp_v2r8       ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
865     int              vdwioffset2;
866     _fjsp_v2r8       ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
867     int              vdwioffset3;
868     _fjsp_v2r8       ix3,iy3,iz3,fix3,fiy3,fiz3,iq3,isai3;
869     int              vdwjidx0A,vdwjidx0B;
870     _fjsp_v2r8       jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
871     _fjsp_v2r8       dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
872     _fjsp_v2r8       dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
873     _fjsp_v2r8       dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
874     _fjsp_v2r8       dx30,dy30,dz30,rsq30,rinv30,rinvsq30,r30,qq30,c6_30,c12_30;
875     _fjsp_v2r8       velec,felec,velecsum,facel,crf,krf,krf2;
876     real             *charge;
877     int              nvdwtype;
878     _fjsp_v2r8       rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
879     int              *vdwtype;
880     real             *vdwparam;
881     _fjsp_v2r8       one_sixth   = gmx_fjsp_set1_v2r8(1.0/6.0);
882     _fjsp_v2r8       one_twelfth = gmx_fjsp_set1_v2r8(1.0/12.0);
883     _fjsp_v2r8       ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
884     real             *ewtab;
885     _fjsp_v2r8       rswitch,swV3,swV4,swV5,swF2,swF3,swF4,d,d2,sw,dsw;
886     real             rswitch_scalar,d_scalar;
887     _fjsp_v2r8       itab_tmp;
888     _fjsp_v2r8       dummy_mask,cutoff_mask;
889     _fjsp_v2r8       one     = gmx_fjsp_set1_v2r8(1.0);
890     _fjsp_v2r8       two     = gmx_fjsp_set1_v2r8(2.0);
891     union { _fjsp_v2r8 simd; long long int i[2]; } vfconv,gbconv,ewconv;
892
893     x                = xx[0];
894     f                = ff[0];
895
896     nri              = nlist->nri;
897     iinr             = nlist->iinr;
898     jindex           = nlist->jindex;
899     jjnr             = nlist->jjnr;
900     shiftidx         = nlist->shift;
901     gid              = nlist->gid;
902     shiftvec         = fr->shift_vec[0];
903     fshift           = fr->fshift[0];
904     facel            = gmx_fjsp_set1_v2r8(fr->epsfac);
905     charge           = mdatoms->chargeA;
906     nvdwtype         = fr->ntype;
907     vdwparam         = fr->nbfp;
908     vdwtype          = mdatoms->typeA;
909
910     sh_ewald         = gmx_fjsp_set1_v2r8(fr->ic->sh_ewald);
911     ewtab            = fr->ic->tabq_coul_FDV0;
912     ewtabscale       = gmx_fjsp_set1_v2r8(fr->ic->tabq_scale);
913     ewtabhalfspace   = gmx_fjsp_set1_v2r8(0.5/fr->ic->tabq_scale);
914
915     /* Setup water-specific parameters */
916     inr              = nlist->iinr[0];
917     iq1              = _fjsp_mul_v2r8(facel,gmx_fjsp_set1_v2r8(charge[inr+1]));
918     iq2              = _fjsp_mul_v2r8(facel,gmx_fjsp_set1_v2r8(charge[inr+2]));
919     iq3              = _fjsp_mul_v2r8(facel,gmx_fjsp_set1_v2r8(charge[inr+3]));
920     vdwioffset0      = 2*nvdwtype*vdwtype[inr+0];
921
922     /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
923     rcutoff_scalar   = fr->rcoulomb;
924     rcutoff          = gmx_fjsp_set1_v2r8(rcutoff_scalar);
925     rcutoff2         = _fjsp_mul_v2r8(rcutoff,rcutoff);
926
927     rswitch_scalar   = fr->rcoulomb_switch;
928     rswitch          = gmx_fjsp_set1_v2r8(rswitch_scalar);
929     /* Setup switch parameters */
930     d_scalar         = rcutoff_scalar-rswitch_scalar;
931     d                = gmx_fjsp_set1_v2r8(d_scalar);
932     swV3             = gmx_fjsp_set1_v2r8(-10.0/(d_scalar*d_scalar*d_scalar));
933     swV4             = gmx_fjsp_set1_v2r8( 15.0/(d_scalar*d_scalar*d_scalar*d_scalar));
934     swV5             = gmx_fjsp_set1_v2r8( -6.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
935     swF2             = gmx_fjsp_set1_v2r8(-30.0/(d_scalar*d_scalar*d_scalar));
936     swF3             = gmx_fjsp_set1_v2r8( 60.0/(d_scalar*d_scalar*d_scalar*d_scalar));
937     swF4             = gmx_fjsp_set1_v2r8(-30.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
938
939     /* Avoid stupid compiler warnings */
940     jnrA = jnrB = 0;
941     j_coord_offsetA = 0;
942     j_coord_offsetB = 0;
943
944     outeriter        = 0;
945     inneriter        = 0;
946
947     /* Start outer loop over neighborlists */
948     for(iidx=0; iidx<nri; iidx++)
949     {
950         /* Load shift vector for this list */
951         i_shift_offset   = DIM*shiftidx[iidx];
952
953         /* Load limits for loop over neighbors */
954         j_index_start    = jindex[iidx];
955         j_index_end      = jindex[iidx+1];
956
957         /* Get outer coordinate index */
958         inr              = iinr[iidx];
959         i_coord_offset   = DIM*inr;
960
961         /* Load i particle coords and add shift vector */
962         gmx_fjsp_load_shift_and_4rvec_broadcast_v2r8(shiftvec+i_shift_offset,x+i_coord_offset,
963                                                  &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
964
965         fix0             = _fjsp_setzero_v2r8();
966         fiy0             = _fjsp_setzero_v2r8();
967         fiz0             = _fjsp_setzero_v2r8();
968         fix1             = _fjsp_setzero_v2r8();
969         fiy1             = _fjsp_setzero_v2r8();
970         fiz1             = _fjsp_setzero_v2r8();
971         fix2             = _fjsp_setzero_v2r8();
972         fiy2             = _fjsp_setzero_v2r8();
973         fiz2             = _fjsp_setzero_v2r8();
974         fix3             = _fjsp_setzero_v2r8();
975         fiy3             = _fjsp_setzero_v2r8();
976         fiz3             = _fjsp_setzero_v2r8();
977
978         /* Start inner kernel loop */
979         for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
980         {
981
982             /* Get j neighbor index, and coordinate index */
983             jnrA             = jjnr[jidx];
984             jnrB             = jjnr[jidx+1];
985             j_coord_offsetA  = DIM*jnrA;
986             j_coord_offsetB  = DIM*jnrB;
987
988             /* load j atom coordinates */
989             gmx_fjsp_load_1rvec_2ptr_swizzle_v2r8(x+j_coord_offsetA,x+j_coord_offsetB,
990                                               &jx0,&jy0,&jz0);
991
992             /* Calculate displacement vector */
993             dx00             = _fjsp_sub_v2r8(ix0,jx0);
994             dy00             = _fjsp_sub_v2r8(iy0,jy0);
995             dz00             = _fjsp_sub_v2r8(iz0,jz0);
996             dx10             = _fjsp_sub_v2r8(ix1,jx0);
997             dy10             = _fjsp_sub_v2r8(iy1,jy0);
998             dz10             = _fjsp_sub_v2r8(iz1,jz0);
999             dx20             = _fjsp_sub_v2r8(ix2,jx0);
1000             dy20             = _fjsp_sub_v2r8(iy2,jy0);
1001             dz20             = _fjsp_sub_v2r8(iz2,jz0);
1002             dx30             = _fjsp_sub_v2r8(ix3,jx0);
1003             dy30             = _fjsp_sub_v2r8(iy3,jy0);
1004             dz30             = _fjsp_sub_v2r8(iz3,jz0);
1005
1006             /* Calculate squared distance and things based on it */
1007             rsq00            = gmx_fjsp_calc_rsq_v2r8(dx00,dy00,dz00);
1008             rsq10            = gmx_fjsp_calc_rsq_v2r8(dx10,dy10,dz10);
1009             rsq20            = gmx_fjsp_calc_rsq_v2r8(dx20,dy20,dz20);
1010             rsq30            = gmx_fjsp_calc_rsq_v2r8(dx30,dy30,dz30);
1011
1012             rinv00           = gmx_fjsp_invsqrt_v2r8(rsq00);
1013             rinv10           = gmx_fjsp_invsqrt_v2r8(rsq10);
1014             rinv20           = gmx_fjsp_invsqrt_v2r8(rsq20);
1015             rinv30           = gmx_fjsp_invsqrt_v2r8(rsq30);
1016
1017             rinvsq00         = _fjsp_mul_v2r8(rinv00,rinv00);
1018             rinvsq10         = _fjsp_mul_v2r8(rinv10,rinv10);
1019             rinvsq20         = _fjsp_mul_v2r8(rinv20,rinv20);
1020             rinvsq30         = _fjsp_mul_v2r8(rinv30,rinv30);
1021
1022             /* Load parameters for j particles */
1023             jq0              = gmx_fjsp_load_2real_swizzle_v2r8(charge+jnrA+0,charge+jnrB+0);
1024             vdwjidx0A        = 2*vdwtype[jnrA+0];
1025             vdwjidx0B        = 2*vdwtype[jnrB+0];
1026
1027             fjx0             = _fjsp_setzero_v2r8();
1028             fjy0             = _fjsp_setzero_v2r8();
1029             fjz0             = _fjsp_setzero_v2r8();
1030
1031             /**************************
1032              * CALCULATE INTERACTIONS *
1033              **************************/
1034
1035             if (gmx_fjsp_any_lt_v2r8(rsq00,rcutoff2))
1036             {
1037
1038             r00              = _fjsp_mul_v2r8(rsq00,rinv00);
1039
1040             /* Compute parameters for interactions between i and j atoms */
1041             gmx_fjsp_load_2pair_swizzle_v2r8(vdwparam+vdwioffset0+vdwjidx0A,
1042                                          vdwparam+vdwioffset0+vdwjidx0B,&c6_00,&c12_00);
1043
1044             /* LENNARD-JONES DISPERSION/REPULSION */
1045
1046             rinvsix          = _fjsp_mul_v2r8(_fjsp_mul_v2r8(rinvsq00,rinvsq00),rinvsq00);
1047             vvdw6            = _fjsp_mul_v2r8(c6_00,rinvsix);
1048             vvdw12           = _fjsp_mul_v2r8(c12_00,_fjsp_mul_v2r8(rinvsix,rinvsix));
1049             vvdw             = _fjsp_msub_v2r8( vvdw12,one_twelfth, _fjsp_mul_v2r8(vvdw6,one_sixth) );
1050             fvdw             = _fjsp_mul_v2r8(_fjsp_sub_v2r8(vvdw12,vvdw6),rinvsq00);
1051
1052             d                = _fjsp_sub_v2r8(r00,rswitch);
1053             d                = _fjsp_max_v2r8(d,_fjsp_setzero_v2r8());
1054             d2               = _fjsp_mul_v2r8(d,d);
1055             sw               = _fjsp_add_v2r8(one,_fjsp_mul_v2r8(d2,_fjsp_mul_v2r8(d,_fjsp_madd_v2r8(d,_fjsp_madd_v2r8(d,swV5,swV4),swV3))));
1056
1057             dsw              = _fjsp_mul_v2r8(d2,_fjsp_madd_v2r8(d,_fjsp_madd_v2r8(d,swF4,swF3),swF2));
1058
1059             /* Evaluate switch function */
1060             /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
1061             fvdw             = _fjsp_msub_v2r8( fvdw,sw , _fjsp_mul_v2r8(rinv00,_fjsp_mul_v2r8(vvdw,dsw)) );
1062             cutoff_mask      = _fjsp_cmplt_v2r8(rsq00,rcutoff2);
1063
1064             fscal            = fvdw;
1065
1066             fscal            = _fjsp_and_v2r8(fscal,cutoff_mask);
1067
1068             /* Update vectorial force */
1069             fix0             = _fjsp_madd_v2r8(dx00,fscal,fix0);
1070             fiy0             = _fjsp_madd_v2r8(dy00,fscal,fiy0);
1071             fiz0             = _fjsp_madd_v2r8(dz00,fscal,fiz0);
1072             
1073             fjx0             = _fjsp_madd_v2r8(dx00,fscal,fjx0);
1074             fjy0             = _fjsp_madd_v2r8(dy00,fscal,fjy0);
1075             fjz0             = _fjsp_madd_v2r8(dz00,fscal,fjz0);
1076
1077             }
1078
1079             /**************************
1080              * CALCULATE INTERACTIONS *
1081              **************************/
1082
1083             if (gmx_fjsp_any_lt_v2r8(rsq10,rcutoff2))
1084             {
1085
1086             r10              = _fjsp_mul_v2r8(rsq10,rinv10);
1087
1088             /* Compute parameters for interactions between i and j atoms */
1089             qq10             = _fjsp_mul_v2r8(iq1,jq0);
1090
1091             /* EWALD ELECTROSTATICS */
1092
1093             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
1094             ewrt             = _fjsp_mul_v2r8(r10,ewtabscale);
1095             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
1096             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
1097             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
1098
1099             ewtabF           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[0] );
1100             ewtabD           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[1] );
1101             GMX_FJSP_TRANSPOSE2_V2R8(ewtabF,ewtabD);
1102             ewtabV           = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[0] +2);
1103             ewtabFn          = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[1] +2);
1104             GMX_FJSP_TRANSPOSE2_V2R8(ewtabV,ewtabFn);
1105             felec            = _fjsp_madd_v2r8(eweps,ewtabD,ewtabF);
1106             velec            = _fjsp_nmsub_v2r8(_fjsp_mul_v2r8(ewtabhalfspace,eweps) ,_fjsp_add_v2r8(ewtabF,felec), ewtabV);
1107             velec            = _fjsp_mul_v2r8(qq10,_fjsp_sub_v2r8(rinv10,velec));
1108             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq10,rinv10),_fjsp_sub_v2r8(rinvsq10,felec));
1109
1110             d                = _fjsp_sub_v2r8(r10,rswitch);
1111             d                = _fjsp_max_v2r8(d,_fjsp_setzero_v2r8());
1112             d2               = _fjsp_mul_v2r8(d,d);
1113             sw               = _fjsp_add_v2r8(one,_fjsp_mul_v2r8(d2,_fjsp_mul_v2r8(d,_fjsp_madd_v2r8(d,_fjsp_madd_v2r8(d,swV5,swV4),swV3))));
1114
1115             dsw              = _fjsp_mul_v2r8(d2,_fjsp_madd_v2r8(d,_fjsp_madd_v2r8(d,swF4,swF3),swF2));
1116
1117             /* Evaluate switch function */
1118             /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
1119             felec            = _fjsp_msub_v2r8( felec,sw , _fjsp_mul_v2r8(rinv10,_fjsp_mul_v2r8(velec,dsw)) );
1120             cutoff_mask      = _fjsp_cmplt_v2r8(rsq10,rcutoff2);
1121
1122             fscal            = felec;
1123
1124             fscal            = _fjsp_and_v2r8(fscal,cutoff_mask);
1125
1126             /* Update vectorial force */
1127             fix1             = _fjsp_madd_v2r8(dx10,fscal,fix1);
1128             fiy1             = _fjsp_madd_v2r8(dy10,fscal,fiy1);
1129             fiz1             = _fjsp_madd_v2r8(dz10,fscal,fiz1);
1130             
1131             fjx0             = _fjsp_madd_v2r8(dx10,fscal,fjx0);
1132             fjy0             = _fjsp_madd_v2r8(dy10,fscal,fjy0);
1133             fjz0             = _fjsp_madd_v2r8(dz10,fscal,fjz0);
1134
1135             }
1136
1137             /**************************
1138              * CALCULATE INTERACTIONS *
1139              **************************/
1140
1141             if (gmx_fjsp_any_lt_v2r8(rsq20,rcutoff2))
1142             {
1143
1144             r20              = _fjsp_mul_v2r8(rsq20,rinv20);
1145
1146             /* Compute parameters for interactions between i and j atoms */
1147             qq20             = _fjsp_mul_v2r8(iq2,jq0);
1148
1149             /* EWALD ELECTROSTATICS */
1150
1151             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
1152             ewrt             = _fjsp_mul_v2r8(r20,ewtabscale);
1153             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
1154             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
1155             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
1156
1157             ewtabF           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[0] );
1158             ewtabD           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[1] );
1159             GMX_FJSP_TRANSPOSE2_V2R8(ewtabF,ewtabD);
1160             ewtabV           = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[0] +2);
1161             ewtabFn          = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[1] +2);
1162             GMX_FJSP_TRANSPOSE2_V2R8(ewtabV,ewtabFn);
1163             felec            = _fjsp_madd_v2r8(eweps,ewtabD,ewtabF);
1164             velec            = _fjsp_nmsub_v2r8(_fjsp_mul_v2r8(ewtabhalfspace,eweps) ,_fjsp_add_v2r8(ewtabF,felec), ewtabV);
1165             velec            = _fjsp_mul_v2r8(qq20,_fjsp_sub_v2r8(rinv20,velec));
1166             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq20,rinv20),_fjsp_sub_v2r8(rinvsq20,felec));
1167
1168             d                = _fjsp_sub_v2r8(r20,rswitch);
1169             d                = _fjsp_max_v2r8(d,_fjsp_setzero_v2r8());
1170             d2               = _fjsp_mul_v2r8(d,d);
1171             sw               = _fjsp_add_v2r8(one,_fjsp_mul_v2r8(d2,_fjsp_mul_v2r8(d,_fjsp_madd_v2r8(d,_fjsp_madd_v2r8(d,swV5,swV4),swV3))));
1172
1173             dsw              = _fjsp_mul_v2r8(d2,_fjsp_madd_v2r8(d,_fjsp_madd_v2r8(d,swF4,swF3),swF2));
1174
1175             /* Evaluate switch function */
1176             /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
1177             felec            = _fjsp_msub_v2r8( felec,sw , _fjsp_mul_v2r8(rinv20,_fjsp_mul_v2r8(velec,dsw)) );
1178             cutoff_mask      = _fjsp_cmplt_v2r8(rsq20,rcutoff2);
1179
1180             fscal            = felec;
1181
1182             fscal            = _fjsp_and_v2r8(fscal,cutoff_mask);
1183
1184             /* Update vectorial force */
1185             fix2             = _fjsp_madd_v2r8(dx20,fscal,fix2);
1186             fiy2             = _fjsp_madd_v2r8(dy20,fscal,fiy2);
1187             fiz2             = _fjsp_madd_v2r8(dz20,fscal,fiz2);
1188             
1189             fjx0             = _fjsp_madd_v2r8(dx20,fscal,fjx0);
1190             fjy0             = _fjsp_madd_v2r8(dy20,fscal,fjy0);
1191             fjz0             = _fjsp_madd_v2r8(dz20,fscal,fjz0);
1192
1193             }
1194
1195             /**************************
1196              * CALCULATE INTERACTIONS *
1197              **************************/
1198
1199             if (gmx_fjsp_any_lt_v2r8(rsq30,rcutoff2))
1200             {
1201
1202             r30              = _fjsp_mul_v2r8(rsq30,rinv30);
1203
1204             /* Compute parameters for interactions between i and j atoms */
1205             qq30             = _fjsp_mul_v2r8(iq3,jq0);
1206
1207             /* EWALD ELECTROSTATICS */
1208
1209             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
1210             ewrt             = _fjsp_mul_v2r8(r30,ewtabscale);
1211             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
1212             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
1213             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
1214
1215             ewtabF           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[0] );
1216             ewtabD           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[1] );
1217             GMX_FJSP_TRANSPOSE2_V2R8(ewtabF,ewtabD);
1218             ewtabV           = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[0] +2);
1219             ewtabFn          = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[1] +2);
1220             GMX_FJSP_TRANSPOSE2_V2R8(ewtabV,ewtabFn);
1221             felec            = _fjsp_madd_v2r8(eweps,ewtabD,ewtabF);
1222             velec            = _fjsp_nmsub_v2r8(_fjsp_mul_v2r8(ewtabhalfspace,eweps) ,_fjsp_add_v2r8(ewtabF,felec), ewtabV);
1223             velec            = _fjsp_mul_v2r8(qq30,_fjsp_sub_v2r8(rinv30,velec));
1224             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq30,rinv30),_fjsp_sub_v2r8(rinvsq30,felec));
1225
1226             d                = _fjsp_sub_v2r8(r30,rswitch);
1227             d                = _fjsp_max_v2r8(d,_fjsp_setzero_v2r8());
1228             d2               = _fjsp_mul_v2r8(d,d);
1229             sw               = _fjsp_add_v2r8(one,_fjsp_mul_v2r8(d2,_fjsp_mul_v2r8(d,_fjsp_madd_v2r8(d,_fjsp_madd_v2r8(d,swV5,swV4),swV3))));
1230
1231             dsw              = _fjsp_mul_v2r8(d2,_fjsp_madd_v2r8(d,_fjsp_madd_v2r8(d,swF4,swF3),swF2));
1232
1233             /* Evaluate switch function */
1234             /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
1235             felec            = _fjsp_msub_v2r8( felec,sw , _fjsp_mul_v2r8(rinv30,_fjsp_mul_v2r8(velec,dsw)) );
1236             cutoff_mask      = _fjsp_cmplt_v2r8(rsq30,rcutoff2);
1237
1238             fscal            = felec;
1239
1240             fscal            = _fjsp_and_v2r8(fscal,cutoff_mask);
1241
1242             /* Update vectorial force */
1243             fix3             = _fjsp_madd_v2r8(dx30,fscal,fix3);
1244             fiy3             = _fjsp_madd_v2r8(dy30,fscal,fiy3);
1245             fiz3             = _fjsp_madd_v2r8(dz30,fscal,fiz3);
1246             
1247             fjx0             = _fjsp_madd_v2r8(dx30,fscal,fjx0);
1248             fjy0             = _fjsp_madd_v2r8(dy30,fscal,fjy0);
1249             fjz0             = _fjsp_madd_v2r8(dz30,fscal,fjz0);
1250
1251             }
1252
1253             gmx_fjsp_decrement_1rvec_2ptr_swizzle_v2r8(f+j_coord_offsetA,f+j_coord_offsetB,fjx0,fjy0,fjz0);
1254
1255             /* Inner loop uses 257 flops */
1256         }
1257
1258         if(jidx<j_index_end)
1259         {
1260
1261             jnrA             = jjnr[jidx];
1262             j_coord_offsetA  = DIM*jnrA;
1263
1264             /* load j atom coordinates */
1265             gmx_fjsp_load_1rvec_1ptr_swizzle_v2r8(x+j_coord_offsetA,
1266                                               &jx0,&jy0,&jz0);
1267
1268             /* Calculate displacement vector */
1269             dx00             = _fjsp_sub_v2r8(ix0,jx0);
1270             dy00             = _fjsp_sub_v2r8(iy0,jy0);
1271             dz00             = _fjsp_sub_v2r8(iz0,jz0);
1272             dx10             = _fjsp_sub_v2r8(ix1,jx0);
1273             dy10             = _fjsp_sub_v2r8(iy1,jy0);
1274             dz10             = _fjsp_sub_v2r8(iz1,jz0);
1275             dx20             = _fjsp_sub_v2r8(ix2,jx0);
1276             dy20             = _fjsp_sub_v2r8(iy2,jy0);
1277             dz20             = _fjsp_sub_v2r8(iz2,jz0);
1278             dx30             = _fjsp_sub_v2r8(ix3,jx0);
1279             dy30             = _fjsp_sub_v2r8(iy3,jy0);
1280             dz30             = _fjsp_sub_v2r8(iz3,jz0);
1281
1282             /* Calculate squared distance and things based on it */
1283             rsq00            = gmx_fjsp_calc_rsq_v2r8(dx00,dy00,dz00);
1284             rsq10            = gmx_fjsp_calc_rsq_v2r8(dx10,dy10,dz10);
1285             rsq20            = gmx_fjsp_calc_rsq_v2r8(dx20,dy20,dz20);
1286             rsq30            = gmx_fjsp_calc_rsq_v2r8(dx30,dy30,dz30);
1287
1288             rinv00           = gmx_fjsp_invsqrt_v2r8(rsq00);
1289             rinv10           = gmx_fjsp_invsqrt_v2r8(rsq10);
1290             rinv20           = gmx_fjsp_invsqrt_v2r8(rsq20);
1291             rinv30           = gmx_fjsp_invsqrt_v2r8(rsq30);
1292
1293             rinvsq00         = _fjsp_mul_v2r8(rinv00,rinv00);
1294             rinvsq10         = _fjsp_mul_v2r8(rinv10,rinv10);
1295             rinvsq20         = _fjsp_mul_v2r8(rinv20,rinv20);
1296             rinvsq30         = _fjsp_mul_v2r8(rinv30,rinv30);
1297
1298             /* Load parameters for j particles */
1299             jq0              = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(),charge+jnrA+0);
1300             vdwjidx0A        = 2*vdwtype[jnrA+0];
1301
1302             fjx0             = _fjsp_setzero_v2r8();
1303             fjy0             = _fjsp_setzero_v2r8();
1304             fjz0             = _fjsp_setzero_v2r8();
1305
1306             /**************************
1307              * CALCULATE INTERACTIONS *
1308              **************************/
1309
1310             if (gmx_fjsp_any_lt_v2r8(rsq00,rcutoff2))
1311             {
1312
1313             r00              = _fjsp_mul_v2r8(rsq00,rinv00);
1314
1315             /* Compute parameters for interactions between i and j atoms */
1316             gmx_fjsp_load_1pair_swizzle_v2r8(vdwparam+vdwioffset0+vdwjidx0A,&c6_00,&c12_00);
1317
1318             /* LENNARD-JONES DISPERSION/REPULSION */
1319
1320             rinvsix          = _fjsp_mul_v2r8(_fjsp_mul_v2r8(rinvsq00,rinvsq00),rinvsq00);
1321             vvdw6            = _fjsp_mul_v2r8(c6_00,rinvsix);
1322             vvdw12           = _fjsp_mul_v2r8(c12_00,_fjsp_mul_v2r8(rinvsix,rinvsix));
1323             vvdw             = _fjsp_msub_v2r8( vvdw12,one_twelfth, _fjsp_mul_v2r8(vvdw6,one_sixth) );
1324             fvdw             = _fjsp_mul_v2r8(_fjsp_sub_v2r8(vvdw12,vvdw6),rinvsq00);
1325
1326             d                = _fjsp_sub_v2r8(r00,rswitch);
1327             d                = _fjsp_max_v2r8(d,_fjsp_setzero_v2r8());
1328             d2               = _fjsp_mul_v2r8(d,d);
1329             sw               = _fjsp_add_v2r8(one,_fjsp_mul_v2r8(d2,_fjsp_mul_v2r8(d,_fjsp_madd_v2r8(d,_fjsp_madd_v2r8(d,swV5,swV4),swV3))));
1330
1331             dsw              = _fjsp_mul_v2r8(d2,_fjsp_madd_v2r8(d,_fjsp_madd_v2r8(d,swF4,swF3),swF2));
1332
1333             /* Evaluate switch function */
1334             /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
1335             fvdw             = _fjsp_msub_v2r8( fvdw,sw , _fjsp_mul_v2r8(rinv00,_fjsp_mul_v2r8(vvdw,dsw)) );
1336             cutoff_mask      = _fjsp_cmplt_v2r8(rsq00,rcutoff2);
1337
1338             fscal            = fvdw;
1339
1340             fscal            = _fjsp_and_v2r8(fscal,cutoff_mask);
1341
1342             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
1343
1344             /* Update vectorial force */
1345             fix0             = _fjsp_madd_v2r8(dx00,fscal,fix0);
1346             fiy0             = _fjsp_madd_v2r8(dy00,fscal,fiy0);
1347             fiz0             = _fjsp_madd_v2r8(dz00,fscal,fiz0);
1348             
1349             fjx0             = _fjsp_madd_v2r8(dx00,fscal,fjx0);
1350             fjy0             = _fjsp_madd_v2r8(dy00,fscal,fjy0);
1351             fjz0             = _fjsp_madd_v2r8(dz00,fscal,fjz0);
1352
1353             }
1354
1355             /**************************
1356              * CALCULATE INTERACTIONS *
1357              **************************/
1358
1359             if (gmx_fjsp_any_lt_v2r8(rsq10,rcutoff2))
1360             {
1361
1362             r10              = _fjsp_mul_v2r8(rsq10,rinv10);
1363
1364             /* Compute parameters for interactions between i and j atoms */
1365             qq10             = _fjsp_mul_v2r8(iq1,jq0);
1366
1367             /* EWALD ELECTROSTATICS */
1368
1369             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
1370             ewrt             = _fjsp_mul_v2r8(r10,ewtabscale);
1371             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
1372             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
1373             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
1374
1375             ewtabF           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[0] );
1376             ewtabD           = _fjsp_setzero_v2r8();
1377             GMX_FJSP_TRANSPOSE2_V2R8(ewtabF,ewtabD);
1378             ewtabV           = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[0] +2);
1379             ewtabFn          = _fjsp_setzero_v2r8();
1380             GMX_FJSP_TRANSPOSE2_V2R8(ewtabV,ewtabFn);
1381             felec            = _fjsp_madd_v2r8(eweps,ewtabD,ewtabF);
1382             velec            = _fjsp_nmsub_v2r8(_fjsp_mul_v2r8(ewtabhalfspace,eweps) ,_fjsp_add_v2r8(ewtabF,felec), ewtabV);
1383             velec            = _fjsp_mul_v2r8(qq10,_fjsp_sub_v2r8(rinv10,velec));
1384             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq10,rinv10),_fjsp_sub_v2r8(rinvsq10,felec));
1385
1386             d                = _fjsp_sub_v2r8(r10,rswitch);
1387             d                = _fjsp_max_v2r8(d,_fjsp_setzero_v2r8());
1388             d2               = _fjsp_mul_v2r8(d,d);
1389             sw               = _fjsp_add_v2r8(one,_fjsp_mul_v2r8(d2,_fjsp_mul_v2r8(d,_fjsp_madd_v2r8(d,_fjsp_madd_v2r8(d,swV5,swV4),swV3))));
1390
1391             dsw              = _fjsp_mul_v2r8(d2,_fjsp_madd_v2r8(d,_fjsp_madd_v2r8(d,swF4,swF3),swF2));
1392
1393             /* Evaluate switch function */
1394             /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
1395             felec            = _fjsp_msub_v2r8( felec,sw , _fjsp_mul_v2r8(rinv10,_fjsp_mul_v2r8(velec,dsw)) );
1396             cutoff_mask      = _fjsp_cmplt_v2r8(rsq10,rcutoff2);
1397
1398             fscal            = felec;
1399
1400             fscal            = _fjsp_and_v2r8(fscal,cutoff_mask);
1401
1402             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
1403
1404             /* Update vectorial force */
1405             fix1             = _fjsp_madd_v2r8(dx10,fscal,fix1);
1406             fiy1             = _fjsp_madd_v2r8(dy10,fscal,fiy1);
1407             fiz1             = _fjsp_madd_v2r8(dz10,fscal,fiz1);
1408             
1409             fjx0             = _fjsp_madd_v2r8(dx10,fscal,fjx0);
1410             fjy0             = _fjsp_madd_v2r8(dy10,fscal,fjy0);
1411             fjz0             = _fjsp_madd_v2r8(dz10,fscal,fjz0);
1412
1413             }
1414
1415             /**************************
1416              * CALCULATE INTERACTIONS *
1417              **************************/
1418
1419             if (gmx_fjsp_any_lt_v2r8(rsq20,rcutoff2))
1420             {
1421
1422             r20              = _fjsp_mul_v2r8(rsq20,rinv20);
1423
1424             /* Compute parameters for interactions between i and j atoms */
1425             qq20             = _fjsp_mul_v2r8(iq2,jq0);
1426
1427             /* EWALD ELECTROSTATICS */
1428
1429             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
1430             ewrt             = _fjsp_mul_v2r8(r20,ewtabscale);
1431             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
1432             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
1433             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
1434
1435             ewtabF           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[0] );
1436             ewtabD           = _fjsp_setzero_v2r8();
1437             GMX_FJSP_TRANSPOSE2_V2R8(ewtabF,ewtabD);
1438             ewtabV           = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[0] +2);
1439             ewtabFn          = _fjsp_setzero_v2r8();
1440             GMX_FJSP_TRANSPOSE2_V2R8(ewtabV,ewtabFn);
1441             felec            = _fjsp_madd_v2r8(eweps,ewtabD,ewtabF);
1442             velec            = _fjsp_nmsub_v2r8(_fjsp_mul_v2r8(ewtabhalfspace,eweps) ,_fjsp_add_v2r8(ewtabF,felec), ewtabV);
1443             velec            = _fjsp_mul_v2r8(qq20,_fjsp_sub_v2r8(rinv20,velec));
1444             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq20,rinv20),_fjsp_sub_v2r8(rinvsq20,felec));
1445
1446             d                = _fjsp_sub_v2r8(r20,rswitch);
1447             d                = _fjsp_max_v2r8(d,_fjsp_setzero_v2r8());
1448             d2               = _fjsp_mul_v2r8(d,d);
1449             sw               = _fjsp_add_v2r8(one,_fjsp_mul_v2r8(d2,_fjsp_mul_v2r8(d,_fjsp_madd_v2r8(d,_fjsp_madd_v2r8(d,swV5,swV4),swV3))));
1450
1451             dsw              = _fjsp_mul_v2r8(d2,_fjsp_madd_v2r8(d,_fjsp_madd_v2r8(d,swF4,swF3),swF2));
1452
1453             /* Evaluate switch function */
1454             /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
1455             felec            = _fjsp_msub_v2r8( felec,sw , _fjsp_mul_v2r8(rinv20,_fjsp_mul_v2r8(velec,dsw)) );
1456             cutoff_mask      = _fjsp_cmplt_v2r8(rsq20,rcutoff2);
1457
1458             fscal            = felec;
1459
1460             fscal            = _fjsp_and_v2r8(fscal,cutoff_mask);
1461
1462             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
1463
1464             /* Update vectorial force */
1465             fix2             = _fjsp_madd_v2r8(dx20,fscal,fix2);
1466             fiy2             = _fjsp_madd_v2r8(dy20,fscal,fiy2);
1467             fiz2             = _fjsp_madd_v2r8(dz20,fscal,fiz2);
1468             
1469             fjx0             = _fjsp_madd_v2r8(dx20,fscal,fjx0);
1470             fjy0             = _fjsp_madd_v2r8(dy20,fscal,fjy0);
1471             fjz0             = _fjsp_madd_v2r8(dz20,fscal,fjz0);
1472
1473             }
1474
1475             /**************************
1476              * CALCULATE INTERACTIONS *
1477              **************************/
1478
1479             if (gmx_fjsp_any_lt_v2r8(rsq30,rcutoff2))
1480             {
1481
1482             r30              = _fjsp_mul_v2r8(rsq30,rinv30);
1483
1484             /* Compute parameters for interactions between i and j atoms */
1485             qq30             = _fjsp_mul_v2r8(iq3,jq0);
1486
1487             /* EWALD ELECTROSTATICS */
1488
1489             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
1490             ewrt             = _fjsp_mul_v2r8(r30,ewtabscale);
1491             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
1492             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
1493             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
1494
1495             ewtabF           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[0] );
1496             ewtabD           = _fjsp_setzero_v2r8();
1497             GMX_FJSP_TRANSPOSE2_V2R8(ewtabF,ewtabD);
1498             ewtabV           = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[0] +2);
1499             ewtabFn          = _fjsp_setzero_v2r8();
1500             GMX_FJSP_TRANSPOSE2_V2R8(ewtabV,ewtabFn);
1501             felec            = _fjsp_madd_v2r8(eweps,ewtabD,ewtabF);
1502             velec            = _fjsp_nmsub_v2r8(_fjsp_mul_v2r8(ewtabhalfspace,eweps) ,_fjsp_add_v2r8(ewtabF,felec), ewtabV);
1503             velec            = _fjsp_mul_v2r8(qq30,_fjsp_sub_v2r8(rinv30,velec));
1504             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq30,rinv30),_fjsp_sub_v2r8(rinvsq30,felec));
1505
1506             d                = _fjsp_sub_v2r8(r30,rswitch);
1507             d                = _fjsp_max_v2r8(d,_fjsp_setzero_v2r8());
1508             d2               = _fjsp_mul_v2r8(d,d);
1509             sw               = _fjsp_add_v2r8(one,_fjsp_mul_v2r8(d2,_fjsp_mul_v2r8(d,_fjsp_madd_v2r8(d,_fjsp_madd_v2r8(d,swV5,swV4),swV3))));
1510
1511             dsw              = _fjsp_mul_v2r8(d2,_fjsp_madd_v2r8(d,_fjsp_madd_v2r8(d,swF4,swF3),swF2));
1512
1513             /* Evaluate switch function */
1514             /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
1515             felec            = _fjsp_msub_v2r8( felec,sw , _fjsp_mul_v2r8(rinv30,_fjsp_mul_v2r8(velec,dsw)) );
1516             cutoff_mask      = _fjsp_cmplt_v2r8(rsq30,rcutoff2);
1517
1518             fscal            = felec;
1519
1520             fscal            = _fjsp_and_v2r8(fscal,cutoff_mask);
1521
1522             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
1523
1524             /* Update vectorial force */
1525             fix3             = _fjsp_madd_v2r8(dx30,fscal,fix3);
1526             fiy3             = _fjsp_madd_v2r8(dy30,fscal,fiy3);
1527             fiz3             = _fjsp_madd_v2r8(dz30,fscal,fiz3);
1528             
1529             fjx0             = _fjsp_madd_v2r8(dx30,fscal,fjx0);
1530             fjy0             = _fjsp_madd_v2r8(dy30,fscal,fjy0);
1531             fjz0             = _fjsp_madd_v2r8(dz30,fscal,fjz0);
1532
1533             }
1534
1535             gmx_fjsp_decrement_1rvec_1ptr_swizzle_v2r8(f+j_coord_offsetA,fjx0,fjy0,fjz0);
1536
1537             /* Inner loop uses 257 flops */
1538         }
1539
1540         /* End of innermost loop */
1541
1542         gmx_fjsp_update_iforce_4atom_swizzle_v2r8(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
1543                                               f+i_coord_offset,fshift+i_shift_offset);
1544
1545         /* Increment number of inner iterations */
1546         inneriter                  += j_index_end - j_index_start;
1547
1548         /* Outer loop uses 24 flops */
1549     }
1550
1551     /* Increment number of outer iterations */
1552     outeriter        += nri;
1553
1554     /* Update outer/inner flops */
1555
1556     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W4_F,outeriter*24 + inneriter*257);
1557 }