Fujitsu Sparc64 acceleration and general fixes for non-x86 builds
[alexxy/gromacs.git] / src / gmxlib / nonbonded / nb_kernel_sparc64_hpc_ace_double / nb_kernel_ElecEw_VdwNone_GeomW4W4_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_ElecEw_VdwNone_GeomW4W4_VF_sparc64_hpc_ace_double
53  * Electrostatics interaction: Ewald
54  * VdW interaction:            None
55  * Geometry:                   Water4-Water4
56  * Calculate force/pot:        PotentialAndForce
57  */
58 void
59 nb_kernel_ElecEw_VdwNone_GeomW4W4_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              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              vdwjidx1A,vdwjidx1B;
88     _fjsp_v2r8       jx1,jy1,jz1,fjx1,fjy1,fjz1,jq1,isaj1;
89     int              vdwjidx2A,vdwjidx2B;
90     _fjsp_v2r8       jx2,jy2,jz2,fjx2,fjy2,fjz2,jq2,isaj2;
91     int              vdwjidx3A,vdwjidx3B;
92     _fjsp_v2r8       jx3,jy3,jz3,fjx3,fjy3,fjz3,jq3,isaj3;
93     _fjsp_v2r8       dx11,dy11,dz11,rsq11,rinv11,rinvsq11,r11,qq11,c6_11,c12_11;
94     _fjsp_v2r8       dx12,dy12,dz12,rsq12,rinv12,rinvsq12,r12,qq12,c6_12,c12_12;
95     _fjsp_v2r8       dx13,dy13,dz13,rsq13,rinv13,rinvsq13,r13,qq13,c6_13,c12_13;
96     _fjsp_v2r8       dx21,dy21,dz21,rsq21,rinv21,rinvsq21,r21,qq21,c6_21,c12_21;
97     _fjsp_v2r8       dx22,dy22,dz22,rsq22,rinv22,rinvsq22,r22,qq22,c6_22,c12_22;
98     _fjsp_v2r8       dx23,dy23,dz23,rsq23,rinv23,rinvsq23,r23,qq23,c6_23,c12_23;
99     _fjsp_v2r8       dx31,dy31,dz31,rsq31,rinv31,rinvsq31,r31,qq31,c6_31,c12_31;
100     _fjsp_v2r8       dx32,dy32,dz32,rsq32,rinv32,rinvsq32,r32,qq32,c6_32,c12_32;
101     _fjsp_v2r8       dx33,dy33,dz33,rsq33,rinv33,rinvsq33,r33,qq33,c6_33,c12_33;
102     _fjsp_v2r8       velec,felec,velecsum,facel,crf,krf,krf2;
103     real             *charge;
104     _fjsp_v2r8       ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
105     real             *ewtab;
106     _fjsp_v2r8       itab_tmp;
107     _fjsp_v2r8       dummy_mask,cutoff_mask;
108     _fjsp_v2r8       one     = gmx_fjsp_set1_v2r8(1.0);
109     _fjsp_v2r8       two     = gmx_fjsp_set1_v2r8(2.0);
110     union { _fjsp_v2r8 simd; long long int i[2]; } vfconv,gbconv,ewconv;
111
112     x                = xx[0];
113     f                = ff[0];
114
115     nri              = nlist->nri;
116     iinr             = nlist->iinr;
117     jindex           = nlist->jindex;
118     jjnr             = nlist->jjnr;
119     shiftidx         = nlist->shift;
120     gid              = nlist->gid;
121     shiftvec         = fr->shift_vec[0];
122     fshift           = fr->fshift[0];
123     facel            = gmx_fjsp_set1_v2r8(fr->epsfac);
124     charge           = mdatoms->chargeA;
125
126     sh_ewald         = gmx_fjsp_set1_v2r8(fr->ic->sh_ewald);
127     ewtab            = fr->ic->tabq_coul_FDV0;
128     ewtabscale       = gmx_fjsp_set1_v2r8(fr->ic->tabq_scale);
129     ewtabhalfspace   = gmx_fjsp_set1_v2r8(0.5/fr->ic->tabq_scale);
130
131     /* Setup water-specific parameters */
132     inr              = nlist->iinr[0];
133     iq1              = _fjsp_mul_v2r8(facel,gmx_fjsp_set1_v2r8(charge[inr+1]));
134     iq2              = _fjsp_mul_v2r8(facel,gmx_fjsp_set1_v2r8(charge[inr+2]));
135     iq3              = _fjsp_mul_v2r8(facel,gmx_fjsp_set1_v2r8(charge[inr+3]));
136
137     jq1              = gmx_fjsp_set1_v2r8(charge[inr+1]);
138     jq2              = gmx_fjsp_set1_v2r8(charge[inr+2]);
139     jq3              = gmx_fjsp_set1_v2r8(charge[inr+3]);
140     qq11             = _fjsp_mul_v2r8(iq1,jq1);
141     qq12             = _fjsp_mul_v2r8(iq1,jq2);
142     qq13             = _fjsp_mul_v2r8(iq1,jq3);
143     qq21             = _fjsp_mul_v2r8(iq2,jq1);
144     qq22             = _fjsp_mul_v2r8(iq2,jq2);
145     qq23             = _fjsp_mul_v2r8(iq2,jq3);
146     qq31             = _fjsp_mul_v2r8(iq3,jq1);
147     qq32             = _fjsp_mul_v2r8(iq3,jq2);
148     qq33             = _fjsp_mul_v2r8(iq3,jq3);
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_3rvec_broadcast_v2r8(shiftvec+i_shift_offset,x+i_coord_offset+DIM,
174                                                  &ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
175
176         fix1             = _fjsp_setzero_v2r8();
177         fiy1             = _fjsp_setzero_v2r8();
178         fiz1             = _fjsp_setzero_v2r8();
179         fix2             = _fjsp_setzero_v2r8();
180         fiy2             = _fjsp_setzero_v2r8();
181         fiz2             = _fjsp_setzero_v2r8();
182         fix3             = _fjsp_setzero_v2r8();
183         fiy3             = _fjsp_setzero_v2r8();
184         fiz3             = _fjsp_setzero_v2r8();
185
186         /* Reset potential sums */
187         velecsum         = _fjsp_setzero_v2r8();
188
189         /* Start inner kernel loop */
190         for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
191         {
192
193             /* Get j neighbor index, and coordinate index */
194             jnrA             = jjnr[jidx];
195             jnrB             = jjnr[jidx+1];
196             j_coord_offsetA  = DIM*jnrA;
197             j_coord_offsetB  = DIM*jnrB;
198
199             /* load j atom coordinates */
200             gmx_fjsp_load_3rvec_2ptr_swizzle_v2r8(x+j_coord_offsetA+DIM,x+j_coord_offsetB+DIM,
201                                               &jx1,&jy1,&jz1,&jx2,&jy2,&jz2,&jx3,&jy3,&jz3);
202
203             /* Calculate displacement vector */
204             dx11             = _fjsp_sub_v2r8(ix1,jx1);
205             dy11             = _fjsp_sub_v2r8(iy1,jy1);
206             dz11             = _fjsp_sub_v2r8(iz1,jz1);
207             dx12             = _fjsp_sub_v2r8(ix1,jx2);
208             dy12             = _fjsp_sub_v2r8(iy1,jy2);
209             dz12             = _fjsp_sub_v2r8(iz1,jz2);
210             dx13             = _fjsp_sub_v2r8(ix1,jx3);
211             dy13             = _fjsp_sub_v2r8(iy1,jy3);
212             dz13             = _fjsp_sub_v2r8(iz1,jz3);
213             dx21             = _fjsp_sub_v2r8(ix2,jx1);
214             dy21             = _fjsp_sub_v2r8(iy2,jy1);
215             dz21             = _fjsp_sub_v2r8(iz2,jz1);
216             dx22             = _fjsp_sub_v2r8(ix2,jx2);
217             dy22             = _fjsp_sub_v2r8(iy2,jy2);
218             dz22             = _fjsp_sub_v2r8(iz2,jz2);
219             dx23             = _fjsp_sub_v2r8(ix2,jx3);
220             dy23             = _fjsp_sub_v2r8(iy2,jy3);
221             dz23             = _fjsp_sub_v2r8(iz2,jz3);
222             dx31             = _fjsp_sub_v2r8(ix3,jx1);
223             dy31             = _fjsp_sub_v2r8(iy3,jy1);
224             dz31             = _fjsp_sub_v2r8(iz3,jz1);
225             dx32             = _fjsp_sub_v2r8(ix3,jx2);
226             dy32             = _fjsp_sub_v2r8(iy3,jy2);
227             dz32             = _fjsp_sub_v2r8(iz3,jz2);
228             dx33             = _fjsp_sub_v2r8(ix3,jx3);
229             dy33             = _fjsp_sub_v2r8(iy3,jy3);
230             dz33             = _fjsp_sub_v2r8(iz3,jz3);
231
232             /* Calculate squared distance and things based on it */
233             rsq11            = gmx_fjsp_calc_rsq_v2r8(dx11,dy11,dz11);
234             rsq12            = gmx_fjsp_calc_rsq_v2r8(dx12,dy12,dz12);
235             rsq13            = gmx_fjsp_calc_rsq_v2r8(dx13,dy13,dz13);
236             rsq21            = gmx_fjsp_calc_rsq_v2r8(dx21,dy21,dz21);
237             rsq22            = gmx_fjsp_calc_rsq_v2r8(dx22,dy22,dz22);
238             rsq23            = gmx_fjsp_calc_rsq_v2r8(dx23,dy23,dz23);
239             rsq31            = gmx_fjsp_calc_rsq_v2r8(dx31,dy31,dz31);
240             rsq32            = gmx_fjsp_calc_rsq_v2r8(dx32,dy32,dz32);
241             rsq33            = gmx_fjsp_calc_rsq_v2r8(dx33,dy33,dz33);
242
243             rinv11           = gmx_fjsp_invsqrt_v2r8(rsq11);
244             rinv12           = gmx_fjsp_invsqrt_v2r8(rsq12);
245             rinv13           = gmx_fjsp_invsqrt_v2r8(rsq13);
246             rinv21           = gmx_fjsp_invsqrt_v2r8(rsq21);
247             rinv22           = gmx_fjsp_invsqrt_v2r8(rsq22);
248             rinv23           = gmx_fjsp_invsqrt_v2r8(rsq23);
249             rinv31           = gmx_fjsp_invsqrt_v2r8(rsq31);
250             rinv32           = gmx_fjsp_invsqrt_v2r8(rsq32);
251             rinv33           = gmx_fjsp_invsqrt_v2r8(rsq33);
252
253             rinvsq11         = _fjsp_mul_v2r8(rinv11,rinv11);
254             rinvsq12         = _fjsp_mul_v2r8(rinv12,rinv12);
255             rinvsq13         = _fjsp_mul_v2r8(rinv13,rinv13);
256             rinvsq21         = _fjsp_mul_v2r8(rinv21,rinv21);
257             rinvsq22         = _fjsp_mul_v2r8(rinv22,rinv22);
258             rinvsq23         = _fjsp_mul_v2r8(rinv23,rinv23);
259             rinvsq31         = _fjsp_mul_v2r8(rinv31,rinv31);
260             rinvsq32         = _fjsp_mul_v2r8(rinv32,rinv32);
261             rinvsq33         = _fjsp_mul_v2r8(rinv33,rinv33);
262
263             fjx1             = _fjsp_setzero_v2r8();
264             fjy1             = _fjsp_setzero_v2r8();
265             fjz1             = _fjsp_setzero_v2r8();
266             fjx2             = _fjsp_setzero_v2r8();
267             fjy2             = _fjsp_setzero_v2r8();
268             fjz2             = _fjsp_setzero_v2r8();
269             fjx3             = _fjsp_setzero_v2r8();
270             fjy3             = _fjsp_setzero_v2r8();
271             fjz3             = _fjsp_setzero_v2r8();
272
273             /**************************
274              * CALCULATE INTERACTIONS *
275              **************************/
276
277             r11              = _fjsp_mul_v2r8(rsq11,rinv11);
278
279             /* EWALD ELECTROSTATICS */
280
281             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
282             ewrt             = _fjsp_mul_v2r8(r11,ewtabscale);
283             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
284             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
285             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
286
287             ewtabF           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[0] );
288             ewtabD           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[1] );
289             GMX_FJSP_TRANSPOSE2_V2R8(ewtabF,ewtabD);
290             ewtabV           = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[0] +2);
291             ewtabFn          = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[1] +2);
292             GMX_FJSP_TRANSPOSE2_V2R8(ewtabV,ewtabFn);
293             felec            = _fjsp_madd_v2r8(eweps,ewtabD,ewtabF);
294             velec            = _fjsp_nmsub_v2r8(_fjsp_mul_v2r8(ewtabhalfspace,eweps) ,_fjsp_add_v2r8(ewtabF,felec), ewtabV);
295             velec            = _fjsp_mul_v2r8(qq11,_fjsp_sub_v2r8(rinv11,velec));
296             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq11,rinv11),_fjsp_sub_v2r8(rinvsq11,felec));
297
298             /* Update potential sum for this i atom from the interaction with this j atom. */
299             velecsum         = _fjsp_add_v2r8(velecsum,velec);
300
301             fscal            = felec;
302
303             /* Update vectorial force */
304             fix1             = _fjsp_madd_v2r8(dx11,fscal,fix1);
305             fiy1             = _fjsp_madd_v2r8(dy11,fscal,fiy1);
306             fiz1             = _fjsp_madd_v2r8(dz11,fscal,fiz1);
307             
308             fjx1             = _fjsp_madd_v2r8(dx11,fscal,fjx1);
309             fjy1             = _fjsp_madd_v2r8(dy11,fscal,fjy1);
310             fjz1             = _fjsp_madd_v2r8(dz11,fscal,fjz1);
311
312             /**************************
313              * CALCULATE INTERACTIONS *
314              **************************/
315
316             r12              = _fjsp_mul_v2r8(rsq12,rinv12);
317
318             /* EWALD ELECTROSTATICS */
319
320             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
321             ewrt             = _fjsp_mul_v2r8(r12,ewtabscale);
322             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
323             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
324             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
325
326             ewtabF           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[0] );
327             ewtabD           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[1] );
328             GMX_FJSP_TRANSPOSE2_V2R8(ewtabF,ewtabD);
329             ewtabV           = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[0] +2);
330             ewtabFn          = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[1] +2);
331             GMX_FJSP_TRANSPOSE2_V2R8(ewtabV,ewtabFn);
332             felec            = _fjsp_madd_v2r8(eweps,ewtabD,ewtabF);
333             velec            = _fjsp_nmsub_v2r8(_fjsp_mul_v2r8(ewtabhalfspace,eweps) ,_fjsp_add_v2r8(ewtabF,felec), ewtabV);
334             velec            = _fjsp_mul_v2r8(qq12,_fjsp_sub_v2r8(rinv12,velec));
335             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq12,rinv12),_fjsp_sub_v2r8(rinvsq12,felec));
336
337             /* Update potential sum for this i atom from the interaction with this j atom. */
338             velecsum         = _fjsp_add_v2r8(velecsum,velec);
339
340             fscal            = felec;
341
342             /* Update vectorial force */
343             fix1             = _fjsp_madd_v2r8(dx12,fscal,fix1);
344             fiy1             = _fjsp_madd_v2r8(dy12,fscal,fiy1);
345             fiz1             = _fjsp_madd_v2r8(dz12,fscal,fiz1);
346             
347             fjx2             = _fjsp_madd_v2r8(dx12,fscal,fjx2);
348             fjy2             = _fjsp_madd_v2r8(dy12,fscal,fjy2);
349             fjz2             = _fjsp_madd_v2r8(dz12,fscal,fjz2);
350
351             /**************************
352              * CALCULATE INTERACTIONS *
353              **************************/
354
355             r13              = _fjsp_mul_v2r8(rsq13,rinv13);
356
357             /* EWALD ELECTROSTATICS */
358
359             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
360             ewrt             = _fjsp_mul_v2r8(r13,ewtabscale);
361             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
362             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
363             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
364
365             ewtabF           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[0] );
366             ewtabD           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[1] );
367             GMX_FJSP_TRANSPOSE2_V2R8(ewtabF,ewtabD);
368             ewtabV           = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[0] +2);
369             ewtabFn          = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[1] +2);
370             GMX_FJSP_TRANSPOSE2_V2R8(ewtabV,ewtabFn);
371             felec            = _fjsp_madd_v2r8(eweps,ewtabD,ewtabF);
372             velec            = _fjsp_nmsub_v2r8(_fjsp_mul_v2r8(ewtabhalfspace,eweps) ,_fjsp_add_v2r8(ewtabF,felec), ewtabV);
373             velec            = _fjsp_mul_v2r8(qq13,_fjsp_sub_v2r8(rinv13,velec));
374             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq13,rinv13),_fjsp_sub_v2r8(rinvsq13,felec));
375
376             /* Update potential sum for this i atom from the interaction with this j atom. */
377             velecsum         = _fjsp_add_v2r8(velecsum,velec);
378
379             fscal            = felec;
380
381             /* Update vectorial force */
382             fix1             = _fjsp_madd_v2r8(dx13,fscal,fix1);
383             fiy1             = _fjsp_madd_v2r8(dy13,fscal,fiy1);
384             fiz1             = _fjsp_madd_v2r8(dz13,fscal,fiz1);
385             
386             fjx3             = _fjsp_madd_v2r8(dx13,fscal,fjx3);
387             fjy3             = _fjsp_madd_v2r8(dy13,fscal,fjy3);
388             fjz3             = _fjsp_madd_v2r8(dz13,fscal,fjz3);
389
390             /**************************
391              * CALCULATE INTERACTIONS *
392              **************************/
393
394             r21              = _fjsp_mul_v2r8(rsq21,rinv21);
395
396             /* EWALD ELECTROSTATICS */
397
398             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
399             ewrt             = _fjsp_mul_v2r8(r21,ewtabscale);
400             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
401             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
402             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
403
404             ewtabF           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[0] );
405             ewtabD           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[1] );
406             GMX_FJSP_TRANSPOSE2_V2R8(ewtabF,ewtabD);
407             ewtabV           = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[0] +2);
408             ewtabFn          = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[1] +2);
409             GMX_FJSP_TRANSPOSE2_V2R8(ewtabV,ewtabFn);
410             felec            = _fjsp_madd_v2r8(eweps,ewtabD,ewtabF);
411             velec            = _fjsp_nmsub_v2r8(_fjsp_mul_v2r8(ewtabhalfspace,eweps) ,_fjsp_add_v2r8(ewtabF,felec), ewtabV);
412             velec            = _fjsp_mul_v2r8(qq21,_fjsp_sub_v2r8(rinv21,velec));
413             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq21,rinv21),_fjsp_sub_v2r8(rinvsq21,felec));
414
415             /* Update potential sum for this i atom from the interaction with this j atom. */
416             velecsum         = _fjsp_add_v2r8(velecsum,velec);
417
418             fscal            = felec;
419
420             /* Update vectorial force */
421             fix2             = _fjsp_madd_v2r8(dx21,fscal,fix2);
422             fiy2             = _fjsp_madd_v2r8(dy21,fscal,fiy2);
423             fiz2             = _fjsp_madd_v2r8(dz21,fscal,fiz2);
424             
425             fjx1             = _fjsp_madd_v2r8(dx21,fscal,fjx1);
426             fjy1             = _fjsp_madd_v2r8(dy21,fscal,fjy1);
427             fjz1             = _fjsp_madd_v2r8(dz21,fscal,fjz1);
428
429             /**************************
430              * CALCULATE INTERACTIONS *
431              **************************/
432
433             r22              = _fjsp_mul_v2r8(rsq22,rinv22);
434
435             /* EWALD ELECTROSTATICS */
436
437             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
438             ewrt             = _fjsp_mul_v2r8(r22,ewtabscale);
439             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
440             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
441             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
442
443             ewtabF           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[0] );
444             ewtabD           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[1] );
445             GMX_FJSP_TRANSPOSE2_V2R8(ewtabF,ewtabD);
446             ewtabV           = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[0] +2);
447             ewtabFn          = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[1] +2);
448             GMX_FJSP_TRANSPOSE2_V2R8(ewtabV,ewtabFn);
449             felec            = _fjsp_madd_v2r8(eweps,ewtabD,ewtabF);
450             velec            = _fjsp_nmsub_v2r8(_fjsp_mul_v2r8(ewtabhalfspace,eweps) ,_fjsp_add_v2r8(ewtabF,felec), ewtabV);
451             velec            = _fjsp_mul_v2r8(qq22,_fjsp_sub_v2r8(rinv22,velec));
452             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq22,rinv22),_fjsp_sub_v2r8(rinvsq22,felec));
453
454             /* Update potential sum for this i atom from the interaction with this j atom. */
455             velecsum         = _fjsp_add_v2r8(velecsum,velec);
456
457             fscal            = felec;
458
459             /* Update vectorial force */
460             fix2             = _fjsp_madd_v2r8(dx22,fscal,fix2);
461             fiy2             = _fjsp_madd_v2r8(dy22,fscal,fiy2);
462             fiz2             = _fjsp_madd_v2r8(dz22,fscal,fiz2);
463             
464             fjx2             = _fjsp_madd_v2r8(dx22,fscal,fjx2);
465             fjy2             = _fjsp_madd_v2r8(dy22,fscal,fjy2);
466             fjz2             = _fjsp_madd_v2r8(dz22,fscal,fjz2);
467
468             /**************************
469              * CALCULATE INTERACTIONS *
470              **************************/
471
472             r23              = _fjsp_mul_v2r8(rsq23,rinv23);
473
474             /* EWALD ELECTROSTATICS */
475
476             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
477             ewrt             = _fjsp_mul_v2r8(r23,ewtabscale);
478             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
479             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
480             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
481
482             ewtabF           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[0] );
483             ewtabD           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[1] );
484             GMX_FJSP_TRANSPOSE2_V2R8(ewtabF,ewtabD);
485             ewtabV           = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[0] +2);
486             ewtabFn          = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[1] +2);
487             GMX_FJSP_TRANSPOSE2_V2R8(ewtabV,ewtabFn);
488             felec            = _fjsp_madd_v2r8(eweps,ewtabD,ewtabF);
489             velec            = _fjsp_nmsub_v2r8(_fjsp_mul_v2r8(ewtabhalfspace,eweps) ,_fjsp_add_v2r8(ewtabF,felec), ewtabV);
490             velec            = _fjsp_mul_v2r8(qq23,_fjsp_sub_v2r8(rinv23,velec));
491             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq23,rinv23),_fjsp_sub_v2r8(rinvsq23,felec));
492
493             /* Update potential sum for this i atom from the interaction with this j atom. */
494             velecsum         = _fjsp_add_v2r8(velecsum,velec);
495
496             fscal            = felec;
497
498             /* Update vectorial force */
499             fix2             = _fjsp_madd_v2r8(dx23,fscal,fix2);
500             fiy2             = _fjsp_madd_v2r8(dy23,fscal,fiy2);
501             fiz2             = _fjsp_madd_v2r8(dz23,fscal,fiz2);
502             
503             fjx3             = _fjsp_madd_v2r8(dx23,fscal,fjx3);
504             fjy3             = _fjsp_madd_v2r8(dy23,fscal,fjy3);
505             fjz3             = _fjsp_madd_v2r8(dz23,fscal,fjz3);
506
507             /**************************
508              * CALCULATE INTERACTIONS *
509              **************************/
510
511             r31              = _fjsp_mul_v2r8(rsq31,rinv31);
512
513             /* EWALD ELECTROSTATICS */
514
515             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
516             ewrt             = _fjsp_mul_v2r8(r31,ewtabscale);
517             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
518             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
519             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
520
521             ewtabF           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[0] );
522             ewtabD           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[1] );
523             GMX_FJSP_TRANSPOSE2_V2R8(ewtabF,ewtabD);
524             ewtabV           = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[0] +2);
525             ewtabFn          = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[1] +2);
526             GMX_FJSP_TRANSPOSE2_V2R8(ewtabV,ewtabFn);
527             felec            = _fjsp_madd_v2r8(eweps,ewtabD,ewtabF);
528             velec            = _fjsp_nmsub_v2r8(_fjsp_mul_v2r8(ewtabhalfspace,eweps) ,_fjsp_add_v2r8(ewtabF,felec), ewtabV);
529             velec            = _fjsp_mul_v2r8(qq31,_fjsp_sub_v2r8(rinv31,velec));
530             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq31,rinv31),_fjsp_sub_v2r8(rinvsq31,felec));
531
532             /* Update potential sum for this i atom from the interaction with this j atom. */
533             velecsum         = _fjsp_add_v2r8(velecsum,velec);
534
535             fscal            = felec;
536
537             /* Update vectorial force */
538             fix3             = _fjsp_madd_v2r8(dx31,fscal,fix3);
539             fiy3             = _fjsp_madd_v2r8(dy31,fscal,fiy3);
540             fiz3             = _fjsp_madd_v2r8(dz31,fscal,fiz3);
541             
542             fjx1             = _fjsp_madd_v2r8(dx31,fscal,fjx1);
543             fjy1             = _fjsp_madd_v2r8(dy31,fscal,fjy1);
544             fjz1             = _fjsp_madd_v2r8(dz31,fscal,fjz1);
545
546             /**************************
547              * CALCULATE INTERACTIONS *
548              **************************/
549
550             r32              = _fjsp_mul_v2r8(rsq32,rinv32);
551
552             /* EWALD ELECTROSTATICS */
553
554             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
555             ewrt             = _fjsp_mul_v2r8(r32,ewtabscale);
556             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
557             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
558             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
559
560             ewtabF           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[0] );
561             ewtabD           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[1] );
562             GMX_FJSP_TRANSPOSE2_V2R8(ewtabF,ewtabD);
563             ewtabV           = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[0] +2);
564             ewtabFn          = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[1] +2);
565             GMX_FJSP_TRANSPOSE2_V2R8(ewtabV,ewtabFn);
566             felec            = _fjsp_madd_v2r8(eweps,ewtabD,ewtabF);
567             velec            = _fjsp_nmsub_v2r8(_fjsp_mul_v2r8(ewtabhalfspace,eweps) ,_fjsp_add_v2r8(ewtabF,felec), ewtabV);
568             velec            = _fjsp_mul_v2r8(qq32,_fjsp_sub_v2r8(rinv32,velec));
569             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq32,rinv32),_fjsp_sub_v2r8(rinvsq32,felec));
570
571             /* Update potential sum for this i atom from the interaction with this j atom. */
572             velecsum         = _fjsp_add_v2r8(velecsum,velec);
573
574             fscal            = felec;
575
576             /* Update vectorial force */
577             fix3             = _fjsp_madd_v2r8(dx32,fscal,fix3);
578             fiy3             = _fjsp_madd_v2r8(dy32,fscal,fiy3);
579             fiz3             = _fjsp_madd_v2r8(dz32,fscal,fiz3);
580             
581             fjx2             = _fjsp_madd_v2r8(dx32,fscal,fjx2);
582             fjy2             = _fjsp_madd_v2r8(dy32,fscal,fjy2);
583             fjz2             = _fjsp_madd_v2r8(dz32,fscal,fjz2);
584
585             /**************************
586              * CALCULATE INTERACTIONS *
587              **************************/
588
589             r33              = _fjsp_mul_v2r8(rsq33,rinv33);
590
591             /* EWALD ELECTROSTATICS */
592
593             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
594             ewrt             = _fjsp_mul_v2r8(r33,ewtabscale);
595             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
596             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
597             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
598
599             ewtabF           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[0] );
600             ewtabD           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[1] );
601             GMX_FJSP_TRANSPOSE2_V2R8(ewtabF,ewtabD);
602             ewtabV           = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[0] +2);
603             ewtabFn          = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[1] +2);
604             GMX_FJSP_TRANSPOSE2_V2R8(ewtabV,ewtabFn);
605             felec            = _fjsp_madd_v2r8(eweps,ewtabD,ewtabF);
606             velec            = _fjsp_nmsub_v2r8(_fjsp_mul_v2r8(ewtabhalfspace,eweps) ,_fjsp_add_v2r8(ewtabF,felec), ewtabV);
607             velec            = _fjsp_mul_v2r8(qq33,_fjsp_sub_v2r8(rinv33,velec));
608             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq33,rinv33),_fjsp_sub_v2r8(rinvsq33,felec));
609
610             /* Update potential sum for this i atom from the interaction with this j atom. */
611             velecsum         = _fjsp_add_v2r8(velecsum,velec);
612
613             fscal            = felec;
614
615             /* Update vectorial force */
616             fix3             = _fjsp_madd_v2r8(dx33,fscal,fix3);
617             fiy3             = _fjsp_madd_v2r8(dy33,fscal,fiy3);
618             fiz3             = _fjsp_madd_v2r8(dz33,fscal,fiz3);
619             
620             fjx3             = _fjsp_madd_v2r8(dx33,fscal,fjx3);
621             fjy3             = _fjsp_madd_v2r8(dy33,fscal,fjy3);
622             fjz3             = _fjsp_madd_v2r8(dz33,fscal,fjz3);
623
624             gmx_fjsp_decrement_3rvec_2ptr_swizzle_v2r8(f+j_coord_offsetA+DIM,f+j_coord_offsetB+DIM,fjx1,fjy1,fjz1,fjx2,fjy2,fjz2,fjx3,fjy3,fjz3);
625
626             /* Inner loop uses 396 flops */
627         }
628
629         if(jidx<j_index_end)
630         {
631
632             jnrA             = jjnr[jidx];
633             j_coord_offsetA  = DIM*jnrA;
634
635             /* load j atom coordinates */
636             gmx_fjsp_load_3rvec_1ptr_swizzle_v2r8(x+j_coord_offsetA+DIM,
637                                               &jx1,&jy1,&jz1,&jx2,&jy2,&jz2,&jx3,&jy3,&jz3);
638
639             /* Calculate displacement vector */
640             dx11             = _fjsp_sub_v2r8(ix1,jx1);
641             dy11             = _fjsp_sub_v2r8(iy1,jy1);
642             dz11             = _fjsp_sub_v2r8(iz1,jz1);
643             dx12             = _fjsp_sub_v2r8(ix1,jx2);
644             dy12             = _fjsp_sub_v2r8(iy1,jy2);
645             dz12             = _fjsp_sub_v2r8(iz1,jz2);
646             dx13             = _fjsp_sub_v2r8(ix1,jx3);
647             dy13             = _fjsp_sub_v2r8(iy1,jy3);
648             dz13             = _fjsp_sub_v2r8(iz1,jz3);
649             dx21             = _fjsp_sub_v2r8(ix2,jx1);
650             dy21             = _fjsp_sub_v2r8(iy2,jy1);
651             dz21             = _fjsp_sub_v2r8(iz2,jz1);
652             dx22             = _fjsp_sub_v2r8(ix2,jx2);
653             dy22             = _fjsp_sub_v2r8(iy2,jy2);
654             dz22             = _fjsp_sub_v2r8(iz2,jz2);
655             dx23             = _fjsp_sub_v2r8(ix2,jx3);
656             dy23             = _fjsp_sub_v2r8(iy2,jy3);
657             dz23             = _fjsp_sub_v2r8(iz2,jz3);
658             dx31             = _fjsp_sub_v2r8(ix3,jx1);
659             dy31             = _fjsp_sub_v2r8(iy3,jy1);
660             dz31             = _fjsp_sub_v2r8(iz3,jz1);
661             dx32             = _fjsp_sub_v2r8(ix3,jx2);
662             dy32             = _fjsp_sub_v2r8(iy3,jy2);
663             dz32             = _fjsp_sub_v2r8(iz3,jz2);
664             dx33             = _fjsp_sub_v2r8(ix3,jx3);
665             dy33             = _fjsp_sub_v2r8(iy3,jy3);
666             dz33             = _fjsp_sub_v2r8(iz3,jz3);
667
668             /* Calculate squared distance and things based on it */
669             rsq11            = gmx_fjsp_calc_rsq_v2r8(dx11,dy11,dz11);
670             rsq12            = gmx_fjsp_calc_rsq_v2r8(dx12,dy12,dz12);
671             rsq13            = gmx_fjsp_calc_rsq_v2r8(dx13,dy13,dz13);
672             rsq21            = gmx_fjsp_calc_rsq_v2r8(dx21,dy21,dz21);
673             rsq22            = gmx_fjsp_calc_rsq_v2r8(dx22,dy22,dz22);
674             rsq23            = gmx_fjsp_calc_rsq_v2r8(dx23,dy23,dz23);
675             rsq31            = gmx_fjsp_calc_rsq_v2r8(dx31,dy31,dz31);
676             rsq32            = gmx_fjsp_calc_rsq_v2r8(dx32,dy32,dz32);
677             rsq33            = gmx_fjsp_calc_rsq_v2r8(dx33,dy33,dz33);
678
679             rinv11           = gmx_fjsp_invsqrt_v2r8(rsq11);
680             rinv12           = gmx_fjsp_invsqrt_v2r8(rsq12);
681             rinv13           = gmx_fjsp_invsqrt_v2r8(rsq13);
682             rinv21           = gmx_fjsp_invsqrt_v2r8(rsq21);
683             rinv22           = gmx_fjsp_invsqrt_v2r8(rsq22);
684             rinv23           = gmx_fjsp_invsqrt_v2r8(rsq23);
685             rinv31           = gmx_fjsp_invsqrt_v2r8(rsq31);
686             rinv32           = gmx_fjsp_invsqrt_v2r8(rsq32);
687             rinv33           = gmx_fjsp_invsqrt_v2r8(rsq33);
688
689             rinvsq11         = _fjsp_mul_v2r8(rinv11,rinv11);
690             rinvsq12         = _fjsp_mul_v2r8(rinv12,rinv12);
691             rinvsq13         = _fjsp_mul_v2r8(rinv13,rinv13);
692             rinvsq21         = _fjsp_mul_v2r8(rinv21,rinv21);
693             rinvsq22         = _fjsp_mul_v2r8(rinv22,rinv22);
694             rinvsq23         = _fjsp_mul_v2r8(rinv23,rinv23);
695             rinvsq31         = _fjsp_mul_v2r8(rinv31,rinv31);
696             rinvsq32         = _fjsp_mul_v2r8(rinv32,rinv32);
697             rinvsq33         = _fjsp_mul_v2r8(rinv33,rinv33);
698
699             fjx1             = _fjsp_setzero_v2r8();
700             fjy1             = _fjsp_setzero_v2r8();
701             fjz1             = _fjsp_setzero_v2r8();
702             fjx2             = _fjsp_setzero_v2r8();
703             fjy2             = _fjsp_setzero_v2r8();
704             fjz2             = _fjsp_setzero_v2r8();
705             fjx3             = _fjsp_setzero_v2r8();
706             fjy3             = _fjsp_setzero_v2r8();
707             fjz3             = _fjsp_setzero_v2r8();
708
709             /**************************
710              * CALCULATE INTERACTIONS *
711              **************************/
712
713             r11              = _fjsp_mul_v2r8(rsq11,rinv11);
714
715             /* EWALD ELECTROSTATICS */
716
717             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
718             ewrt             = _fjsp_mul_v2r8(r11,ewtabscale);
719             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
720             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
721             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
722
723             ewtabF           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[0] );
724             ewtabD           = _fjsp_setzero_v2r8();
725             GMX_FJSP_TRANSPOSE2_V2R8(ewtabF,ewtabD);
726             ewtabV           = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[0] +2);
727             ewtabFn          = _fjsp_setzero_v2r8();
728             GMX_FJSP_TRANSPOSE2_V2R8(ewtabV,ewtabFn);
729             felec            = _fjsp_madd_v2r8(eweps,ewtabD,ewtabF);
730             velec            = _fjsp_nmsub_v2r8(_fjsp_mul_v2r8(ewtabhalfspace,eweps) ,_fjsp_add_v2r8(ewtabF,felec), ewtabV);
731             velec            = _fjsp_mul_v2r8(qq11,_fjsp_sub_v2r8(rinv11,velec));
732             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq11,rinv11),_fjsp_sub_v2r8(rinvsq11,felec));
733
734             /* Update potential sum for this i atom from the interaction with this j atom. */
735             velec            = _fjsp_unpacklo_v2r8(velec,_fjsp_setzero_v2r8());
736             velecsum         = _fjsp_add_v2r8(velecsum,velec);
737
738             fscal            = felec;
739
740             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
741
742             /* Update vectorial force */
743             fix1             = _fjsp_madd_v2r8(dx11,fscal,fix1);
744             fiy1             = _fjsp_madd_v2r8(dy11,fscal,fiy1);
745             fiz1             = _fjsp_madd_v2r8(dz11,fscal,fiz1);
746             
747             fjx1             = _fjsp_madd_v2r8(dx11,fscal,fjx1);
748             fjy1             = _fjsp_madd_v2r8(dy11,fscal,fjy1);
749             fjz1             = _fjsp_madd_v2r8(dz11,fscal,fjz1);
750
751             /**************************
752              * CALCULATE INTERACTIONS *
753              **************************/
754
755             r12              = _fjsp_mul_v2r8(rsq12,rinv12);
756
757             /* EWALD ELECTROSTATICS */
758
759             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
760             ewrt             = _fjsp_mul_v2r8(r12,ewtabscale);
761             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
762             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
763             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
764
765             ewtabF           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[0] );
766             ewtabD           = _fjsp_setzero_v2r8();
767             GMX_FJSP_TRANSPOSE2_V2R8(ewtabF,ewtabD);
768             ewtabV           = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[0] +2);
769             ewtabFn          = _fjsp_setzero_v2r8();
770             GMX_FJSP_TRANSPOSE2_V2R8(ewtabV,ewtabFn);
771             felec            = _fjsp_madd_v2r8(eweps,ewtabD,ewtabF);
772             velec            = _fjsp_nmsub_v2r8(_fjsp_mul_v2r8(ewtabhalfspace,eweps) ,_fjsp_add_v2r8(ewtabF,felec), ewtabV);
773             velec            = _fjsp_mul_v2r8(qq12,_fjsp_sub_v2r8(rinv12,velec));
774             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq12,rinv12),_fjsp_sub_v2r8(rinvsq12,felec));
775
776             /* Update potential sum for this i atom from the interaction with this j atom. */
777             velec            = _fjsp_unpacklo_v2r8(velec,_fjsp_setzero_v2r8());
778             velecsum         = _fjsp_add_v2r8(velecsum,velec);
779
780             fscal            = felec;
781
782             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
783
784             /* Update vectorial force */
785             fix1             = _fjsp_madd_v2r8(dx12,fscal,fix1);
786             fiy1             = _fjsp_madd_v2r8(dy12,fscal,fiy1);
787             fiz1             = _fjsp_madd_v2r8(dz12,fscal,fiz1);
788             
789             fjx2             = _fjsp_madd_v2r8(dx12,fscal,fjx2);
790             fjy2             = _fjsp_madd_v2r8(dy12,fscal,fjy2);
791             fjz2             = _fjsp_madd_v2r8(dz12,fscal,fjz2);
792
793             /**************************
794              * CALCULATE INTERACTIONS *
795              **************************/
796
797             r13              = _fjsp_mul_v2r8(rsq13,rinv13);
798
799             /* EWALD ELECTROSTATICS */
800
801             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
802             ewrt             = _fjsp_mul_v2r8(r13,ewtabscale);
803             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
804             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
805             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
806
807             ewtabF           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[0] );
808             ewtabD           = _fjsp_setzero_v2r8();
809             GMX_FJSP_TRANSPOSE2_V2R8(ewtabF,ewtabD);
810             ewtabV           = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[0] +2);
811             ewtabFn          = _fjsp_setzero_v2r8();
812             GMX_FJSP_TRANSPOSE2_V2R8(ewtabV,ewtabFn);
813             felec            = _fjsp_madd_v2r8(eweps,ewtabD,ewtabF);
814             velec            = _fjsp_nmsub_v2r8(_fjsp_mul_v2r8(ewtabhalfspace,eweps) ,_fjsp_add_v2r8(ewtabF,felec), ewtabV);
815             velec            = _fjsp_mul_v2r8(qq13,_fjsp_sub_v2r8(rinv13,velec));
816             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq13,rinv13),_fjsp_sub_v2r8(rinvsq13,felec));
817
818             /* Update potential sum for this i atom from the interaction with this j atom. */
819             velec            = _fjsp_unpacklo_v2r8(velec,_fjsp_setzero_v2r8());
820             velecsum         = _fjsp_add_v2r8(velecsum,velec);
821
822             fscal            = felec;
823
824             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
825
826             /* Update vectorial force */
827             fix1             = _fjsp_madd_v2r8(dx13,fscal,fix1);
828             fiy1             = _fjsp_madd_v2r8(dy13,fscal,fiy1);
829             fiz1             = _fjsp_madd_v2r8(dz13,fscal,fiz1);
830             
831             fjx3             = _fjsp_madd_v2r8(dx13,fscal,fjx3);
832             fjy3             = _fjsp_madd_v2r8(dy13,fscal,fjy3);
833             fjz3             = _fjsp_madd_v2r8(dz13,fscal,fjz3);
834
835             /**************************
836              * CALCULATE INTERACTIONS *
837              **************************/
838
839             r21              = _fjsp_mul_v2r8(rsq21,rinv21);
840
841             /* EWALD ELECTROSTATICS */
842
843             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
844             ewrt             = _fjsp_mul_v2r8(r21,ewtabscale);
845             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
846             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
847             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
848
849             ewtabF           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[0] );
850             ewtabD           = _fjsp_setzero_v2r8();
851             GMX_FJSP_TRANSPOSE2_V2R8(ewtabF,ewtabD);
852             ewtabV           = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[0] +2);
853             ewtabFn          = _fjsp_setzero_v2r8();
854             GMX_FJSP_TRANSPOSE2_V2R8(ewtabV,ewtabFn);
855             felec            = _fjsp_madd_v2r8(eweps,ewtabD,ewtabF);
856             velec            = _fjsp_nmsub_v2r8(_fjsp_mul_v2r8(ewtabhalfspace,eweps) ,_fjsp_add_v2r8(ewtabF,felec), ewtabV);
857             velec            = _fjsp_mul_v2r8(qq21,_fjsp_sub_v2r8(rinv21,velec));
858             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq21,rinv21),_fjsp_sub_v2r8(rinvsq21,felec));
859
860             /* Update potential sum for this i atom from the interaction with this j atom. */
861             velec            = _fjsp_unpacklo_v2r8(velec,_fjsp_setzero_v2r8());
862             velecsum         = _fjsp_add_v2r8(velecsum,velec);
863
864             fscal            = felec;
865
866             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
867
868             /* Update vectorial force */
869             fix2             = _fjsp_madd_v2r8(dx21,fscal,fix2);
870             fiy2             = _fjsp_madd_v2r8(dy21,fscal,fiy2);
871             fiz2             = _fjsp_madd_v2r8(dz21,fscal,fiz2);
872             
873             fjx1             = _fjsp_madd_v2r8(dx21,fscal,fjx1);
874             fjy1             = _fjsp_madd_v2r8(dy21,fscal,fjy1);
875             fjz1             = _fjsp_madd_v2r8(dz21,fscal,fjz1);
876
877             /**************************
878              * CALCULATE INTERACTIONS *
879              **************************/
880
881             r22              = _fjsp_mul_v2r8(rsq22,rinv22);
882
883             /* EWALD ELECTROSTATICS */
884
885             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
886             ewrt             = _fjsp_mul_v2r8(r22,ewtabscale);
887             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
888             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
889             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
890
891             ewtabF           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[0] );
892             ewtabD           = _fjsp_setzero_v2r8();
893             GMX_FJSP_TRANSPOSE2_V2R8(ewtabF,ewtabD);
894             ewtabV           = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[0] +2);
895             ewtabFn          = _fjsp_setzero_v2r8();
896             GMX_FJSP_TRANSPOSE2_V2R8(ewtabV,ewtabFn);
897             felec            = _fjsp_madd_v2r8(eweps,ewtabD,ewtabF);
898             velec            = _fjsp_nmsub_v2r8(_fjsp_mul_v2r8(ewtabhalfspace,eweps) ,_fjsp_add_v2r8(ewtabF,felec), ewtabV);
899             velec            = _fjsp_mul_v2r8(qq22,_fjsp_sub_v2r8(rinv22,velec));
900             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq22,rinv22),_fjsp_sub_v2r8(rinvsq22,felec));
901
902             /* Update potential sum for this i atom from the interaction with this j atom. */
903             velec            = _fjsp_unpacklo_v2r8(velec,_fjsp_setzero_v2r8());
904             velecsum         = _fjsp_add_v2r8(velecsum,velec);
905
906             fscal            = felec;
907
908             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
909
910             /* Update vectorial force */
911             fix2             = _fjsp_madd_v2r8(dx22,fscal,fix2);
912             fiy2             = _fjsp_madd_v2r8(dy22,fscal,fiy2);
913             fiz2             = _fjsp_madd_v2r8(dz22,fscal,fiz2);
914             
915             fjx2             = _fjsp_madd_v2r8(dx22,fscal,fjx2);
916             fjy2             = _fjsp_madd_v2r8(dy22,fscal,fjy2);
917             fjz2             = _fjsp_madd_v2r8(dz22,fscal,fjz2);
918
919             /**************************
920              * CALCULATE INTERACTIONS *
921              **************************/
922
923             r23              = _fjsp_mul_v2r8(rsq23,rinv23);
924
925             /* EWALD ELECTROSTATICS */
926
927             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
928             ewrt             = _fjsp_mul_v2r8(r23,ewtabscale);
929             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
930             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
931             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
932
933             ewtabF           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[0] );
934             ewtabD           = _fjsp_setzero_v2r8();
935             GMX_FJSP_TRANSPOSE2_V2R8(ewtabF,ewtabD);
936             ewtabV           = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[0] +2);
937             ewtabFn          = _fjsp_setzero_v2r8();
938             GMX_FJSP_TRANSPOSE2_V2R8(ewtabV,ewtabFn);
939             felec            = _fjsp_madd_v2r8(eweps,ewtabD,ewtabF);
940             velec            = _fjsp_nmsub_v2r8(_fjsp_mul_v2r8(ewtabhalfspace,eweps) ,_fjsp_add_v2r8(ewtabF,felec), ewtabV);
941             velec            = _fjsp_mul_v2r8(qq23,_fjsp_sub_v2r8(rinv23,velec));
942             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq23,rinv23),_fjsp_sub_v2r8(rinvsq23,felec));
943
944             /* Update potential sum for this i atom from the interaction with this j atom. */
945             velec            = _fjsp_unpacklo_v2r8(velec,_fjsp_setzero_v2r8());
946             velecsum         = _fjsp_add_v2r8(velecsum,velec);
947
948             fscal            = felec;
949
950             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
951
952             /* Update vectorial force */
953             fix2             = _fjsp_madd_v2r8(dx23,fscal,fix2);
954             fiy2             = _fjsp_madd_v2r8(dy23,fscal,fiy2);
955             fiz2             = _fjsp_madd_v2r8(dz23,fscal,fiz2);
956             
957             fjx3             = _fjsp_madd_v2r8(dx23,fscal,fjx3);
958             fjy3             = _fjsp_madd_v2r8(dy23,fscal,fjy3);
959             fjz3             = _fjsp_madd_v2r8(dz23,fscal,fjz3);
960
961             /**************************
962              * CALCULATE INTERACTIONS *
963              **************************/
964
965             r31              = _fjsp_mul_v2r8(rsq31,rinv31);
966
967             /* EWALD ELECTROSTATICS */
968
969             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
970             ewrt             = _fjsp_mul_v2r8(r31,ewtabscale);
971             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
972             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
973             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
974
975             ewtabF           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[0] );
976             ewtabD           = _fjsp_setzero_v2r8();
977             GMX_FJSP_TRANSPOSE2_V2R8(ewtabF,ewtabD);
978             ewtabV           = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[0] +2);
979             ewtabFn          = _fjsp_setzero_v2r8();
980             GMX_FJSP_TRANSPOSE2_V2R8(ewtabV,ewtabFn);
981             felec            = _fjsp_madd_v2r8(eweps,ewtabD,ewtabF);
982             velec            = _fjsp_nmsub_v2r8(_fjsp_mul_v2r8(ewtabhalfspace,eweps) ,_fjsp_add_v2r8(ewtabF,felec), ewtabV);
983             velec            = _fjsp_mul_v2r8(qq31,_fjsp_sub_v2r8(rinv31,velec));
984             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq31,rinv31),_fjsp_sub_v2r8(rinvsq31,felec));
985
986             /* Update potential sum for this i atom from the interaction with this j atom. */
987             velec            = _fjsp_unpacklo_v2r8(velec,_fjsp_setzero_v2r8());
988             velecsum         = _fjsp_add_v2r8(velecsum,velec);
989
990             fscal            = felec;
991
992             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
993
994             /* Update vectorial force */
995             fix3             = _fjsp_madd_v2r8(dx31,fscal,fix3);
996             fiy3             = _fjsp_madd_v2r8(dy31,fscal,fiy3);
997             fiz3             = _fjsp_madd_v2r8(dz31,fscal,fiz3);
998             
999             fjx1             = _fjsp_madd_v2r8(dx31,fscal,fjx1);
1000             fjy1             = _fjsp_madd_v2r8(dy31,fscal,fjy1);
1001             fjz1             = _fjsp_madd_v2r8(dz31,fscal,fjz1);
1002
1003             /**************************
1004              * CALCULATE INTERACTIONS *
1005              **************************/
1006
1007             r32              = _fjsp_mul_v2r8(rsq32,rinv32);
1008
1009             /* EWALD ELECTROSTATICS */
1010
1011             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
1012             ewrt             = _fjsp_mul_v2r8(r32,ewtabscale);
1013             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
1014             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
1015             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
1016
1017             ewtabF           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[0] );
1018             ewtabD           = _fjsp_setzero_v2r8();
1019             GMX_FJSP_TRANSPOSE2_V2R8(ewtabF,ewtabD);
1020             ewtabV           = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[0] +2);
1021             ewtabFn          = _fjsp_setzero_v2r8();
1022             GMX_FJSP_TRANSPOSE2_V2R8(ewtabV,ewtabFn);
1023             felec            = _fjsp_madd_v2r8(eweps,ewtabD,ewtabF);
1024             velec            = _fjsp_nmsub_v2r8(_fjsp_mul_v2r8(ewtabhalfspace,eweps) ,_fjsp_add_v2r8(ewtabF,felec), ewtabV);
1025             velec            = _fjsp_mul_v2r8(qq32,_fjsp_sub_v2r8(rinv32,velec));
1026             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq32,rinv32),_fjsp_sub_v2r8(rinvsq32,felec));
1027
1028             /* Update potential sum for this i atom from the interaction with this j atom. */
1029             velec            = _fjsp_unpacklo_v2r8(velec,_fjsp_setzero_v2r8());
1030             velecsum         = _fjsp_add_v2r8(velecsum,velec);
1031
1032             fscal            = felec;
1033
1034             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
1035
1036             /* Update vectorial force */
1037             fix3             = _fjsp_madd_v2r8(dx32,fscal,fix3);
1038             fiy3             = _fjsp_madd_v2r8(dy32,fscal,fiy3);
1039             fiz3             = _fjsp_madd_v2r8(dz32,fscal,fiz3);
1040             
1041             fjx2             = _fjsp_madd_v2r8(dx32,fscal,fjx2);
1042             fjy2             = _fjsp_madd_v2r8(dy32,fscal,fjy2);
1043             fjz2             = _fjsp_madd_v2r8(dz32,fscal,fjz2);
1044
1045             /**************************
1046              * CALCULATE INTERACTIONS *
1047              **************************/
1048
1049             r33              = _fjsp_mul_v2r8(rsq33,rinv33);
1050
1051             /* EWALD ELECTROSTATICS */
1052
1053             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
1054             ewrt             = _fjsp_mul_v2r8(r33,ewtabscale);
1055             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
1056             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
1057             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
1058
1059             ewtabF           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[0] );
1060             ewtabD           = _fjsp_setzero_v2r8();
1061             GMX_FJSP_TRANSPOSE2_V2R8(ewtabF,ewtabD);
1062             ewtabV           = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[0] +2);
1063             ewtabFn          = _fjsp_setzero_v2r8();
1064             GMX_FJSP_TRANSPOSE2_V2R8(ewtabV,ewtabFn);
1065             felec            = _fjsp_madd_v2r8(eweps,ewtabD,ewtabF);
1066             velec            = _fjsp_nmsub_v2r8(_fjsp_mul_v2r8(ewtabhalfspace,eweps) ,_fjsp_add_v2r8(ewtabF,felec), ewtabV);
1067             velec            = _fjsp_mul_v2r8(qq33,_fjsp_sub_v2r8(rinv33,velec));
1068             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq33,rinv33),_fjsp_sub_v2r8(rinvsq33,felec));
1069
1070             /* Update potential sum for this i atom from the interaction with this j atom. */
1071             velec            = _fjsp_unpacklo_v2r8(velec,_fjsp_setzero_v2r8());
1072             velecsum         = _fjsp_add_v2r8(velecsum,velec);
1073
1074             fscal            = felec;
1075
1076             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
1077
1078             /* Update vectorial force */
1079             fix3             = _fjsp_madd_v2r8(dx33,fscal,fix3);
1080             fiy3             = _fjsp_madd_v2r8(dy33,fscal,fiy3);
1081             fiz3             = _fjsp_madd_v2r8(dz33,fscal,fiz3);
1082             
1083             fjx3             = _fjsp_madd_v2r8(dx33,fscal,fjx3);
1084             fjy3             = _fjsp_madd_v2r8(dy33,fscal,fjy3);
1085             fjz3             = _fjsp_madd_v2r8(dz33,fscal,fjz3);
1086
1087             gmx_fjsp_decrement_3rvec_1ptr_swizzle_v2r8(f+j_coord_offsetA+DIM,fjx1,fjy1,fjz1,fjx2,fjy2,fjz2,fjx3,fjy3,fjz3);
1088
1089             /* Inner loop uses 396 flops */
1090         }
1091
1092         /* End of innermost loop */
1093
1094         gmx_fjsp_update_iforce_3atom_swizzle_v2r8(fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
1095                                               f+i_coord_offset+DIM,fshift+i_shift_offset);
1096
1097         ggid                        = gid[iidx];
1098         /* Update potential energies */
1099         gmx_fjsp_update_1pot_v2r8(velecsum,kernel_data->energygrp_elec+ggid);
1100
1101         /* Increment number of inner iterations */
1102         inneriter                  += j_index_end - j_index_start;
1103
1104         /* Outer loop uses 19 flops */
1105     }
1106
1107     /* Increment number of outer iterations */
1108     outeriter        += nri;
1109
1110     /* Update outer/inner flops */
1111
1112     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_W4W4_VF,outeriter*19 + inneriter*396);
1113 }
1114 /*
1115  * Gromacs nonbonded kernel:   nb_kernel_ElecEw_VdwNone_GeomW4W4_F_sparc64_hpc_ace_double
1116  * Electrostatics interaction: Ewald
1117  * VdW interaction:            None
1118  * Geometry:                   Water4-Water4
1119  * Calculate force/pot:        Force
1120  */
1121 void
1122 nb_kernel_ElecEw_VdwNone_GeomW4W4_F_sparc64_hpc_ace_double
1123                     (t_nblist * gmx_restrict                nlist,
1124                      rvec * gmx_restrict                    xx,
1125                      rvec * gmx_restrict                    ff,
1126                      t_forcerec * gmx_restrict              fr,
1127                      t_mdatoms * gmx_restrict               mdatoms,
1128                      nb_kernel_data_t * gmx_restrict        kernel_data,
1129                      t_nrnb * gmx_restrict                  nrnb)
1130 {
1131     /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
1132      * just 0 for non-waters.
1133      * Suffixes A,B refer to j loop unrolling done with double precision SIMD, e.g. for the two different
1134      * jnr indices corresponding to data put in the four positions in the SIMD register.
1135      */
1136     int              i_shift_offset,i_coord_offset,outeriter,inneriter;
1137     int              j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
1138     int              jnrA,jnrB;
1139     int              j_coord_offsetA,j_coord_offsetB;
1140     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
1141     real             rcutoff_scalar;
1142     real             *shiftvec,*fshift,*x,*f;
1143     _fjsp_v2r8       tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
1144     int              vdwioffset1;
1145     _fjsp_v2r8       ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
1146     int              vdwioffset2;
1147     _fjsp_v2r8       ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
1148     int              vdwioffset3;
1149     _fjsp_v2r8       ix3,iy3,iz3,fix3,fiy3,fiz3,iq3,isai3;
1150     int              vdwjidx1A,vdwjidx1B;
1151     _fjsp_v2r8       jx1,jy1,jz1,fjx1,fjy1,fjz1,jq1,isaj1;
1152     int              vdwjidx2A,vdwjidx2B;
1153     _fjsp_v2r8       jx2,jy2,jz2,fjx2,fjy2,fjz2,jq2,isaj2;
1154     int              vdwjidx3A,vdwjidx3B;
1155     _fjsp_v2r8       jx3,jy3,jz3,fjx3,fjy3,fjz3,jq3,isaj3;
1156     _fjsp_v2r8       dx11,dy11,dz11,rsq11,rinv11,rinvsq11,r11,qq11,c6_11,c12_11;
1157     _fjsp_v2r8       dx12,dy12,dz12,rsq12,rinv12,rinvsq12,r12,qq12,c6_12,c12_12;
1158     _fjsp_v2r8       dx13,dy13,dz13,rsq13,rinv13,rinvsq13,r13,qq13,c6_13,c12_13;
1159     _fjsp_v2r8       dx21,dy21,dz21,rsq21,rinv21,rinvsq21,r21,qq21,c6_21,c12_21;
1160     _fjsp_v2r8       dx22,dy22,dz22,rsq22,rinv22,rinvsq22,r22,qq22,c6_22,c12_22;
1161     _fjsp_v2r8       dx23,dy23,dz23,rsq23,rinv23,rinvsq23,r23,qq23,c6_23,c12_23;
1162     _fjsp_v2r8       dx31,dy31,dz31,rsq31,rinv31,rinvsq31,r31,qq31,c6_31,c12_31;
1163     _fjsp_v2r8       dx32,dy32,dz32,rsq32,rinv32,rinvsq32,r32,qq32,c6_32,c12_32;
1164     _fjsp_v2r8       dx33,dy33,dz33,rsq33,rinv33,rinvsq33,r33,qq33,c6_33,c12_33;
1165     _fjsp_v2r8       velec,felec,velecsum,facel,crf,krf,krf2;
1166     real             *charge;
1167     _fjsp_v2r8       ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
1168     real             *ewtab;
1169     _fjsp_v2r8       itab_tmp;
1170     _fjsp_v2r8       dummy_mask,cutoff_mask;
1171     _fjsp_v2r8       one     = gmx_fjsp_set1_v2r8(1.0);
1172     _fjsp_v2r8       two     = gmx_fjsp_set1_v2r8(2.0);
1173     union { _fjsp_v2r8 simd; long long int i[2]; } vfconv,gbconv,ewconv;
1174
1175     x                = xx[0];
1176     f                = ff[0];
1177
1178     nri              = nlist->nri;
1179     iinr             = nlist->iinr;
1180     jindex           = nlist->jindex;
1181     jjnr             = nlist->jjnr;
1182     shiftidx         = nlist->shift;
1183     gid              = nlist->gid;
1184     shiftvec         = fr->shift_vec[0];
1185     fshift           = fr->fshift[0];
1186     facel            = gmx_fjsp_set1_v2r8(fr->epsfac);
1187     charge           = mdatoms->chargeA;
1188
1189     sh_ewald         = gmx_fjsp_set1_v2r8(fr->ic->sh_ewald);
1190     ewtab            = fr->ic->tabq_coul_F;
1191     ewtabscale       = gmx_fjsp_set1_v2r8(fr->ic->tabq_scale);
1192     ewtabhalfspace   = gmx_fjsp_set1_v2r8(0.5/fr->ic->tabq_scale);
1193
1194     /* Setup water-specific parameters */
1195     inr              = nlist->iinr[0];
1196     iq1              = _fjsp_mul_v2r8(facel,gmx_fjsp_set1_v2r8(charge[inr+1]));
1197     iq2              = _fjsp_mul_v2r8(facel,gmx_fjsp_set1_v2r8(charge[inr+2]));
1198     iq3              = _fjsp_mul_v2r8(facel,gmx_fjsp_set1_v2r8(charge[inr+3]));
1199
1200     jq1              = gmx_fjsp_set1_v2r8(charge[inr+1]);
1201     jq2              = gmx_fjsp_set1_v2r8(charge[inr+2]);
1202     jq3              = gmx_fjsp_set1_v2r8(charge[inr+3]);
1203     qq11             = _fjsp_mul_v2r8(iq1,jq1);
1204     qq12             = _fjsp_mul_v2r8(iq1,jq2);
1205     qq13             = _fjsp_mul_v2r8(iq1,jq3);
1206     qq21             = _fjsp_mul_v2r8(iq2,jq1);
1207     qq22             = _fjsp_mul_v2r8(iq2,jq2);
1208     qq23             = _fjsp_mul_v2r8(iq2,jq3);
1209     qq31             = _fjsp_mul_v2r8(iq3,jq1);
1210     qq32             = _fjsp_mul_v2r8(iq3,jq2);
1211     qq33             = _fjsp_mul_v2r8(iq3,jq3);
1212
1213     /* Avoid stupid compiler warnings */
1214     jnrA = jnrB = 0;
1215     j_coord_offsetA = 0;
1216     j_coord_offsetB = 0;
1217
1218     outeriter        = 0;
1219     inneriter        = 0;
1220
1221     /* Start outer loop over neighborlists */
1222     for(iidx=0; iidx<nri; iidx++)
1223     {
1224         /* Load shift vector for this list */
1225         i_shift_offset   = DIM*shiftidx[iidx];
1226
1227         /* Load limits for loop over neighbors */
1228         j_index_start    = jindex[iidx];
1229         j_index_end      = jindex[iidx+1];
1230
1231         /* Get outer coordinate index */
1232         inr              = iinr[iidx];
1233         i_coord_offset   = DIM*inr;
1234
1235         /* Load i particle coords and add shift vector */
1236         gmx_fjsp_load_shift_and_3rvec_broadcast_v2r8(shiftvec+i_shift_offset,x+i_coord_offset+DIM,
1237                                                  &ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
1238
1239         fix1             = _fjsp_setzero_v2r8();
1240         fiy1             = _fjsp_setzero_v2r8();
1241         fiz1             = _fjsp_setzero_v2r8();
1242         fix2             = _fjsp_setzero_v2r8();
1243         fiy2             = _fjsp_setzero_v2r8();
1244         fiz2             = _fjsp_setzero_v2r8();
1245         fix3             = _fjsp_setzero_v2r8();
1246         fiy3             = _fjsp_setzero_v2r8();
1247         fiz3             = _fjsp_setzero_v2r8();
1248
1249         /* Start inner kernel loop */
1250         for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
1251         {
1252
1253             /* Get j neighbor index, and coordinate index */
1254             jnrA             = jjnr[jidx];
1255             jnrB             = jjnr[jidx+1];
1256             j_coord_offsetA  = DIM*jnrA;
1257             j_coord_offsetB  = DIM*jnrB;
1258
1259             /* load j atom coordinates */
1260             gmx_fjsp_load_3rvec_2ptr_swizzle_v2r8(x+j_coord_offsetA+DIM,x+j_coord_offsetB+DIM,
1261                                               &jx1,&jy1,&jz1,&jx2,&jy2,&jz2,&jx3,&jy3,&jz3);
1262
1263             /* Calculate displacement vector */
1264             dx11             = _fjsp_sub_v2r8(ix1,jx1);
1265             dy11             = _fjsp_sub_v2r8(iy1,jy1);
1266             dz11             = _fjsp_sub_v2r8(iz1,jz1);
1267             dx12             = _fjsp_sub_v2r8(ix1,jx2);
1268             dy12             = _fjsp_sub_v2r8(iy1,jy2);
1269             dz12             = _fjsp_sub_v2r8(iz1,jz2);
1270             dx13             = _fjsp_sub_v2r8(ix1,jx3);
1271             dy13             = _fjsp_sub_v2r8(iy1,jy3);
1272             dz13             = _fjsp_sub_v2r8(iz1,jz3);
1273             dx21             = _fjsp_sub_v2r8(ix2,jx1);
1274             dy21             = _fjsp_sub_v2r8(iy2,jy1);
1275             dz21             = _fjsp_sub_v2r8(iz2,jz1);
1276             dx22             = _fjsp_sub_v2r8(ix2,jx2);
1277             dy22             = _fjsp_sub_v2r8(iy2,jy2);
1278             dz22             = _fjsp_sub_v2r8(iz2,jz2);
1279             dx23             = _fjsp_sub_v2r8(ix2,jx3);
1280             dy23             = _fjsp_sub_v2r8(iy2,jy3);
1281             dz23             = _fjsp_sub_v2r8(iz2,jz3);
1282             dx31             = _fjsp_sub_v2r8(ix3,jx1);
1283             dy31             = _fjsp_sub_v2r8(iy3,jy1);
1284             dz31             = _fjsp_sub_v2r8(iz3,jz1);
1285             dx32             = _fjsp_sub_v2r8(ix3,jx2);
1286             dy32             = _fjsp_sub_v2r8(iy3,jy2);
1287             dz32             = _fjsp_sub_v2r8(iz3,jz2);
1288             dx33             = _fjsp_sub_v2r8(ix3,jx3);
1289             dy33             = _fjsp_sub_v2r8(iy3,jy3);
1290             dz33             = _fjsp_sub_v2r8(iz3,jz3);
1291
1292             /* Calculate squared distance and things based on it */
1293             rsq11            = gmx_fjsp_calc_rsq_v2r8(dx11,dy11,dz11);
1294             rsq12            = gmx_fjsp_calc_rsq_v2r8(dx12,dy12,dz12);
1295             rsq13            = gmx_fjsp_calc_rsq_v2r8(dx13,dy13,dz13);
1296             rsq21            = gmx_fjsp_calc_rsq_v2r8(dx21,dy21,dz21);
1297             rsq22            = gmx_fjsp_calc_rsq_v2r8(dx22,dy22,dz22);
1298             rsq23            = gmx_fjsp_calc_rsq_v2r8(dx23,dy23,dz23);
1299             rsq31            = gmx_fjsp_calc_rsq_v2r8(dx31,dy31,dz31);
1300             rsq32            = gmx_fjsp_calc_rsq_v2r8(dx32,dy32,dz32);
1301             rsq33            = gmx_fjsp_calc_rsq_v2r8(dx33,dy33,dz33);
1302
1303             rinv11           = gmx_fjsp_invsqrt_v2r8(rsq11);
1304             rinv12           = gmx_fjsp_invsqrt_v2r8(rsq12);
1305             rinv13           = gmx_fjsp_invsqrt_v2r8(rsq13);
1306             rinv21           = gmx_fjsp_invsqrt_v2r8(rsq21);
1307             rinv22           = gmx_fjsp_invsqrt_v2r8(rsq22);
1308             rinv23           = gmx_fjsp_invsqrt_v2r8(rsq23);
1309             rinv31           = gmx_fjsp_invsqrt_v2r8(rsq31);
1310             rinv32           = gmx_fjsp_invsqrt_v2r8(rsq32);
1311             rinv33           = gmx_fjsp_invsqrt_v2r8(rsq33);
1312
1313             rinvsq11         = _fjsp_mul_v2r8(rinv11,rinv11);
1314             rinvsq12         = _fjsp_mul_v2r8(rinv12,rinv12);
1315             rinvsq13         = _fjsp_mul_v2r8(rinv13,rinv13);
1316             rinvsq21         = _fjsp_mul_v2r8(rinv21,rinv21);
1317             rinvsq22         = _fjsp_mul_v2r8(rinv22,rinv22);
1318             rinvsq23         = _fjsp_mul_v2r8(rinv23,rinv23);
1319             rinvsq31         = _fjsp_mul_v2r8(rinv31,rinv31);
1320             rinvsq32         = _fjsp_mul_v2r8(rinv32,rinv32);
1321             rinvsq33         = _fjsp_mul_v2r8(rinv33,rinv33);
1322
1323             fjx1             = _fjsp_setzero_v2r8();
1324             fjy1             = _fjsp_setzero_v2r8();
1325             fjz1             = _fjsp_setzero_v2r8();
1326             fjx2             = _fjsp_setzero_v2r8();
1327             fjy2             = _fjsp_setzero_v2r8();
1328             fjz2             = _fjsp_setzero_v2r8();
1329             fjx3             = _fjsp_setzero_v2r8();
1330             fjy3             = _fjsp_setzero_v2r8();
1331             fjz3             = _fjsp_setzero_v2r8();
1332
1333             /**************************
1334              * CALCULATE INTERACTIONS *
1335              **************************/
1336
1337             r11              = _fjsp_mul_v2r8(rsq11,rinv11);
1338
1339             /* EWALD ELECTROSTATICS */
1340
1341             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
1342             ewrt             = _fjsp_mul_v2r8(r11,ewtabscale);
1343             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
1344             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
1345             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
1346
1347             gmx_fjsp_load_2pair_swizzle_v2r8(ewtab+ewconv.i[0],ewtab+ewconv.i[1],
1348                                          &ewtabF,&ewtabFn);
1349             felec            = _fjsp_madd_v2r8(eweps,ewtabFn,_fjsp_nmsub_v2r8(eweps,ewtabF,ewtabF));
1350             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq11,rinv11),_fjsp_sub_v2r8(rinvsq11,felec));
1351
1352             fscal            = felec;
1353
1354             /* Update vectorial force */
1355             fix1             = _fjsp_madd_v2r8(dx11,fscal,fix1);
1356             fiy1             = _fjsp_madd_v2r8(dy11,fscal,fiy1);
1357             fiz1             = _fjsp_madd_v2r8(dz11,fscal,fiz1);
1358             
1359             fjx1             = _fjsp_madd_v2r8(dx11,fscal,fjx1);
1360             fjy1             = _fjsp_madd_v2r8(dy11,fscal,fjy1);
1361             fjz1             = _fjsp_madd_v2r8(dz11,fscal,fjz1);
1362
1363             /**************************
1364              * CALCULATE INTERACTIONS *
1365              **************************/
1366
1367             r12              = _fjsp_mul_v2r8(rsq12,rinv12);
1368
1369             /* EWALD ELECTROSTATICS */
1370
1371             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
1372             ewrt             = _fjsp_mul_v2r8(r12,ewtabscale);
1373             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
1374             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
1375             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
1376
1377             gmx_fjsp_load_2pair_swizzle_v2r8(ewtab+ewconv.i[0],ewtab+ewconv.i[1],
1378                                          &ewtabF,&ewtabFn);
1379             felec            = _fjsp_madd_v2r8(eweps,ewtabFn,_fjsp_nmsub_v2r8(eweps,ewtabF,ewtabF));
1380             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq12,rinv12),_fjsp_sub_v2r8(rinvsq12,felec));
1381
1382             fscal            = felec;
1383
1384             /* Update vectorial force */
1385             fix1             = _fjsp_madd_v2r8(dx12,fscal,fix1);
1386             fiy1             = _fjsp_madd_v2r8(dy12,fscal,fiy1);
1387             fiz1             = _fjsp_madd_v2r8(dz12,fscal,fiz1);
1388             
1389             fjx2             = _fjsp_madd_v2r8(dx12,fscal,fjx2);
1390             fjy2             = _fjsp_madd_v2r8(dy12,fscal,fjy2);
1391             fjz2             = _fjsp_madd_v2r8(dz12,fscal,fjz2);
1392
1393             /**************************
1394              * CALCULATE INTERACTIONS *
1395              **************************/
1396
1397             r13              = _fjsp_mul_v2r8(rsq13,rinv13);
1398
1399             /* EWALD ELECTROSTATICS */
1400
1401             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
1402             ewrt             = _fjsp_mul_v2r8(r13,ewtabscale);
1403             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
1404             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
1405             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
1406
1407             gmx_fjsp_load_2pair_swizzle_v2r8(ewtab+ewconv.i[0],ewtab+ewconv.i[1],
1408                                          &ewtabF,&ewtabFn);
1409             felec            = _fjsp_madd_v2r8(eweps,ewtabFn,_fjsp_nmsub_v2r8(eweps,ewtabF,ewtabF));
1410             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq13,rinv13),_fjsp_sub_v2r8(rinvsq13,felec));
1411
1412             fscal            = felec;
1413
1414             /* Update vectorial force */
1415             fix1             = _fjsp_madd_v2r8(dx13,fscal,fix1);
1416             fiy1             = _fjsp_madd_v2r8(dy13,fscal,fiy1);
1417             fiz1             = _fjsp_madd_v2r8(dz13,fscal,fiz1);
1418             
1419             fjx3             = _fjsp_madd_v2r8(dx13,fscal,fjx3);
1420             fjy3             = _fjsp_madd_v2r8(dy13,fscal,fjy3);
1421             fjz3             = _fjsp_madd_v2r8(dz13,fscal,fjz3);
1422
1423             /**************************
1424              * CALCULATE INTERACTIONS *
1425              **************************/
1426
1427             r21              = _fjsp_mul_v2r8(rsq21,rinv21);
1428
1429             /* EWALD ELECTROSTATICS */
1430
1431             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
1432             ewrt             = _fjsp_mul_v2r8(r21,ewtabscale);
1433             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
1434             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
1435             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
1436
1437             gmx_fjsp_load_2pair_swizzle_v2r8(ewtab+ewconv.i[0],ewtab+ewconv.i[1],
1438                                          &ewtabF,&ewtabFn);
1439             felec            = _fjsp_madd_v2r8(eweps,ewtabFn,_fjsp_nmsub_v2r8(eweps,ewtabF,ewtabF));
1440             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq21,rinv21),_fjsp_sub_v2r8(rinvsq21,felec));
1441
1442             fscal            = felec;
1443
1444             /* Update vectorial force */
1445             fix2             = _fjsp_madd_v2r8(dx21,fscal,fix2);
1446             fiy2             = _fjsp_madd_v2r8(dy21,fscal,fiy2);
1447             fiz2             = _fjsp_madd_v2r8(dz21,fscal,fiz2);
1448             
1449             fjx1             = _fjsp_madd_v2r8(dx21,fscal,fjx1);
1450             fjy1             = _fjsp_madd_v2r8(dy21,fscal,fjy1);
1451             fjz1             = _fjsp_madd_v2r8(dz21,fscal,fjz1);
1452
1453             /**************************
1454              * CALCULATE INTERACTIONS *
1455              **************************/
1456
1457             r22              = _fjsp_mul_v2r8(rsq22,rinv22);
1458
1459             /* EWALD ELECTROSTATICS */
1460
1461             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
1462             ewrt             = _fjsp_mul_v2r8(r22,ewtabscale);
1463             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
1464             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
1465             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
1466
1467             gmx_fjsp_load_2pair_swizzle_v2r8(ewtab+ewconv.i[0],ewtab+ewconv.i[1],
1468                                          &ewtabF,&ewtabFn);
1469             felec            = _fjsp_madd_v2r8(eweps,ewtabFn,_fjsp_nmsub_v2r8(eweps,ewtabF,ewtabF));
1470             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq22,rinv22),_fjsp_sub_v2r8(rinvsq22,felec));
1471
1472             fscal            = felec;
1473
1474             /* Update vectorial force */
1475             fix2             = _fjsp_madd_v2r8(dx22,fscal,fix2);
1476             fiy2             = _fjsp_madd_v2r8(dy22,fscal,fiy2);
1477             fiz2             = _fjsp_madd_v2r8(dz22,fscal,fiz2);
1478             
1479             fjx2             = _fjsp_madd_v2r8(dx22,fscal,fjx2);
1480             fjy2             = _fjsp_madd_v2r8(dy22,fscal,fjy2);
1481             fjz2             = _fjsp_madd_v2r8(dz22,fscal,fjz2);
1482
1483             /**************************
1484              * CALCULATE INTERACTIONS *
1485              **************************/
1486
1487             r23              = _fjsp_mul_v2r8(rsq23,rinv23);
1488
1489             /* EWALD ELECTROSTATICS */
1490
1491             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
1492             ewrt             = _fjsp_mul_v2r8(r23,ewtabscale);
1493             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
1494             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
1495             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
1496
1497             gmx_fjsp_load_2pair_swizzle_v2r8(ewtab+ewconv.i[0],ewtab+ewconv.i[1],
1498                                          &ewtabF,&ewtabFn);
1499             felec            = _fjsp_madd_v2r8(eweps,ewtabFn,_fjsp_nmsub_v2r8(eweps,ewtabF,ewtabF));
1500             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq23,rinv23),_fjsp_sub_v2r8(rinvsq23,felec));
1501
1502             fscal            = felec;
1503
1504             /* Update vectorial force */
1505             fix2             = _fjsp_madd_v2r8(dx23,fscal,fix2);
1506             fiy2             = _fjsp_madd_v2r8(dy23,fscal,fiy2);
1507             fiz2             = _fjsp_madd_v2r8(dz23,fscal,fiz2);
1508             
1509             fjx3             = _fjsp_madd_v2r8(dx23,fscal,fjx3);
1510             fjy3             = _fjsp_madd_v2r8(dy23,fscal,fjy3);
1511             fjz3             = _fjsp_madd_v2r8(dz23,fscal,fjz3);
1512
1513             /**************************
1514              * CALCULATE INTERACTIONS *
1515              **************************/
1516
1517             r31              = _fjsp_mul_v2r8(rsq31,rinv31);
1518
1519             /* EWALD ELECTROSTATICS */
1520
1521             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
1522             ewrt             = _fjsp_mul_v2r8(r31,ewtabscale);
1523             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
1524             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
1525             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
1526
1527             gmx_fjsp_load_2pair_swizzle_v2r8(ewtab+ewconv.i[0],ewtab+ewconv.i[1],
1528                                          &ewtabF,&ewtabFn);
1529             felec            = _fjsp_madd_v2r8(eweps,ewtabFn,_fjsp_nmsub_v2r8(eweps,ewtabF,ewtabF));
1530             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq31,rinv31),_fjsp_sub_v2r8(rinvsq31,felec));
1531
1532             fscal            = felec;
1533
1534             /* Update vectorial force */
1535             fix3             = _fjsp_madd_v2r8(dx31,fscal,fix3);
1536             fiy3             = _fjsp_madd_v2r8(dy31,fscal,fiy3);
1537             fiz3             = _fjsp_madd_v2r8(dz31,fscal,fiz3);
1538             
1539             fjx1             = _fjsp_madd_v2r8(dx31,fscal,fjx1);
1540             fjy1             = _fjsp_madd_v2r8(dy31,fscal,fjy1);
1541             fjz1             = _fjsp_madd_v2r8(dz31,fscal,fjz1);
1542
1543             /**************************
1544              * CALCULATE INTERACTIONS *
1545              **************************/
1546
1547             r32              = _fjsp_mul_v2r8(rsq32,rinv32);
1548
1549             /* EWALD ELECTROSTATICS */
1550
1551             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
1552             ewrt             = _fjsp_mul_v2r8(r32,ewtabscale);
1553             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
1554             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
1555             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
1556
1557             gmx_fjsp_load_2pair_swizzle_v2r8(ewtab+ewconv.i[0],ewtab+ewconv.i[1],
1558                                          &ewtabF,&ewtabFn);
1559             felec            = _fjsp_madd_v2r8(eweps,ewtabFn,_fjsp_nmsub_v2r8(eweps,ewtabF,ewtabF));
1560             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq32,rinv32),_fjsp_sub_v2r8(rinvsq32,felec));
1561
1562             fscal            = felec;
1563
1564             /* Update vectorial force */
1565             fix3             = _fjsp_madd_v2r8(dx32,fscal,fix3);
1566             fiy3             = _fjsp_madd_v2r8(dy32,fscal,fiy3);
1567             fiz3             = _fjsp_madd_v2r8(dz32,fscal,fiz3);
1568             
1569             fjx2             = _fjsp_madd_v2r8(dx32,fscal,fjx2);
1570             fjy2             = _fjsp_madd_v2r8(dy32,fscal,fjy2);
1571             fjz2             = _fjsp_madd_v2r8(dz32,fscal,fjz2);
1572
1573             /**************************
1574              * CALCULATE INTERACTIONS *
1575              **************************/
1576
1577             r33              = _fjsp_mul_v2r8(rsq33,rinv33);
1578
1579             /* EWALD ELECTROSTATICS */
1580
1581             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
1582             ewrt             = _fjsp_mul_v2r8(r33,ewtabscale);
1583             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
1584             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
1585             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
1586
1587             gmx_fjsp_load_2pair_swizzle_v2r8(ewtab+ewconv.i[0],ewtab+ewconv.i[1],
1588                                          &ewtabF,&ewtabFn);
1589             felec            = _fjsp_madd_v2r8(eweps,ewtabFn,_fjsp_nmsub_v2r8(eweps,ewtabF,ewtabF));
1590             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq33,rinv33),_fjsp_sub_v2r8(rinvsq33,felec));
1591
1592             fscal            = felec;
1593
1594             /* Update vectorial force */
1595             fix3             = _fjsp_madd_v2r8(dx33,fscal,fix3);
1596             fiy3             = _fjsp_madd_v2r8(dy33,fscal,fiy3);
1597             fiz3             = _fjsp_madd_v2r8(dz33,fscal,fiz3);
1598             
1599             fjx3             = _fjsp_madd_v2r8(dx33,fscal,fjx3);
1600             fjy3             = _fjsp_madd_v2r8(dy33,fscal,fjy3);
1601             fjz3             = _fjsp_madd_v2r8(dz33,fscal,fjz3);
1602
1603             gmx_fjsp_decrement_3rvec_2ptr_swizzle_v2r8(f+j_coord_offsetA+DIM,f+j_coord_offsetB+DIM,fjx1,fjy1,fjz1,fjx2,fjy2,fjz2,fjx3,fjy3,fjz3);
1604
1605             /* Inner loop uses 351 flops */
1606         }
1607
1608         if(jidx<j_index_end)
1609         {
1610
1611             jnrA             = jjnr[jidx];
1612             j_coord_offsetA  = DIM*jnrA;
1613
1614             /* load j atom coordinates */
1615             gmx_fjsp_load_3rvec_1ptr_swizzle_v2r8(x+j_coord_offsetA+DIM,
1616                                               &jx1,&jy1,&jz1,&jx2,&jy2,&jz2,&jx3,&jy3,&jz3);
1617
1618             /* Calculate displacement vector */
1619             dx11             = _fjsp_sub_v2r8(ix1,jx1);
1620             dy11             = _fjsp_sub_v2r8(iy1,jy1);
1621             dz11             = _fjsp_sub_v2r8(iz1,jz1);
1622             dx12             = _fjsp_sub_v2r8(ix1,jx2);
1623             dy12             = _fjsp_sub_v2r8(iy1,jy2);
1624             dz12             = _fjsp_sub_v2r8(iz1,jz2);
1625             dx13             = _fjsp_sub_v2r8(ix1,jx3);
1626             dy13             = _fjsp_sub_v2r8(iy1,jy3);
1627             dz13             = _fjsp_sub_v2r8(iz1,jz3);
1628             dx21             = _fjsp_sub_v2r8(ix2,jx1);
1629             dy21             = _fjsp_sub_v2r8(iy2,jy1);
1630             dz21             = _fjsp_sub_v2r8(iz2,jz1);
1631             dx22             = _fjsp_sub_v2r8(ix2,jx2);
1632             dy22             = _fjsp_sub_v2r8(iy2,jy2);
1633             dz22             = _fjsp_sub_v2r8(iz2,jz2);
1634             dx23             = _fjsp_sub_v2r8(ix2,jx3);
1635             dy23             = _fjsp_sub_v2r8(iy2,jy3);
1636             dz23             = _fjsp_sub_v2r8(iz2,jz3);
1637             dx31             = _fjsp_sub_v2r8(ix3,jx1);
1638             dy31             = _fjsp_sub_v2r8(iy3,jy1);
1639             dz31             = _fjsp_sub_v2r8(iz3,jz1);
1640             dx32             = _fjsp_sub_v2r8(ix3,jx2);
1641             dy32             = _fjsp_sub_v2r8(iy3,jy2);
1642             dz32             = _fjsp_sub_v2r8(iz3,jz2);
1643             dx33             = _fjsp_sub_v2r8(ix3,jx3);
1644             dy33             = _fjsp_sub_v2r8(iy3,jy3);
1645             dz33             = _fjsp_sub_v2r8(iz3,jz3);
1646
1647             /* Calculate squared distance and things based on it */
1648             rsq11            = gmx_fjsp_calc_rsq_v2r8(dx11,dy11,dz11);
1649             rsq12            = gmx_fjsp_calc_rsq_v2r8(dx12,dy12,dz12);
1650             rsq13            = gmx_fjsp_calc_rsq_v2r8(dx13,dy13,dz13);
1651             rsq21            = gmx_fjsp_calc_rsq_v2r8(dx21,dy21,dz21);
1652             rsq22            = gmx_fjsp_calc_rsq_v2r8(dx22,dy22,dz22);
1653             rsq23            = gmx_fjsp_calc_rsq_v2r8(dx23,dy23,dz23);
1654             rsq31            = gmx_fjsp_calc_rsq_v2r8(dx31,dy31,dz31);
1655             rsq32            = gmx_fjsp_calc_rsq_v2r8(dx32,dy32,dz32);
1656             rsq33            = gmx_fjsp_calc_rsq_v2r8(dx33,dy33,dz33);
1657
1658             rinv11           = gmx_fjsp_invsqrt_v2r8(rsq11);
1659             rinv12           = gmx_fjsp_invsqrt_v2r8(rsq12);
1660             rinv13           = gmx_fjsp_invsqrt_v2r8(rsq13);
1661             rinv21           = gmx_fjsp_invsqrt_v2r8(rsq21);
1662             rinv22           = gmx_fjsp_invsqrt_v2r8(rsq22);
1663             rinv23           = gmx_fjsp_invsqrt_v2r8(rsq23);
1664             rinv31           = gmx_fjsp_invsqrt_v2r8(rsq31);
1665             rinv32           = gmx_fjsp_invsqrt_v2r8(rsq32);
1666             rinv33           = gmx_fjsp_invsqrt_v2r8(rsq33);
1667
1668             rinvsq11         = _fjsp_mul_v2r8(rinv11,rinv11);
1669             rinvsq12         = _fjsp_mul_v2r8(rinv12,rinv12);
1670             rinvsq13         = _fjsp_mul_v2r8(rinv13,rinv13);
1671             rinvsq21         = _fjsp_mul_v2r8(rinv21,rinv21);
1672             rinvsq22         = _fjsp_mul_v2r8(rinv22,rinv22);
1673             rinvsq23         = _fjsp_mul_v2r8(rinv23,rinv23);
1674             rinvsq31         = _fjsp_mul_v2r8(rinv31,rinv31);
1675             rinvsq32         = _fjsp_mul_v2r8(rinv32,rinv32);
1676             rinvsq33         = _fjsp_mul_v2r8(rinv33,rinv33);
1677
1678             fjx1             = _fjsp_setzero_v2r8();
1679             fjy1             = _fjsp_setzero_v2r8();
1680             fjz1             = _fjsp_setzero_v2r8();
1681             fjx2             = _fjsp_setzero_v2r8();
1682             fjy2             = _fjsp_setzero_v2r8();
1683             fjz2             = _fjsp_setzero_v2r8();
1684             fjx3             = _fjsp_setzero_v2r8();
1685             fjy3             = _fjsp_setzero_v2r8();
1686             fjz3             = _fjsp_setzero_v2r8();
1687
1688             /**************************
1689              * CALCULATE INTERACTIONS *
1690              **************************/
1691
1692             r11              = _fjsp_mul_v2r8(rsq11,rinv11);
1693
1694             /* EWALD ELECTROSTATICS */
1695
1696             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
1697             ewrt             = _fjsp_mul_v2r8(r11,ewtabscale);
1698             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
1699             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
1700             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
1701
1702             gmx_fjsp_load_1pair_swizzle_v2r8(ewtab+ewconv.i[0],&ewtabF,&ewtabFn);
1703             felec            = _fjsp_madd_v2r8(eweps,ewtabFn,_fjsp_nmsub_v2r8(eweps,ewtabF,ewtabF));
1704             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq11,rinv11),_fjsp_sub_v2r8(rinvsq11,felec));
1705
1706             fscal            = felec;
1707
1708             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
1709
1710             /* Update vectorial force */
1711             fix1             = _fjsp_madd_v2r8(dx11,fscal,fix1);
1712             fiy1             = _fjsp_madd_v2r8(dy11,fscal,fiy1);
1713             fiz1             = _fjsp_madd_v2r8(dz11,fscal,fiz1);
1714             
1715             fjx1             = _fjsp_madd_v2r8(dx11,fscal,fjx1);
1716             fjy1             = _fjsp_madd_v2r8(dy11,fscal,fjy1);
1717             fjz1             = _fjsp_madd_v2r8(dz11,fscal,fjz1);
1718
1719             /**************************
1720              * CALCULATE INTERACTIONS *
1721              **************************/
1722
1723             r12              = _fjsp_mul_v2r8(rsq12,rinv12);
1724
1725             /* EWALD ELECTROSTATICS */
1726
1727             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
1728             ewrt             = _fjsp_mul_v2r8(r12,ewtabscale);
1729             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
1730             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
1731             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
1732
1733             gmx_fjsp_load_1pair_swizzle_v2r8(ewtab+ewconv.i[0],&ewtabF,&ewtabFn);
1734             felec            = _fjsp_madd_v2r8(eweps,ewtabFn,_fjsp_nmsub_v2r8(eweps,ewtabF,ewtabF));
1735             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq12,rinv12),_fjsp_sub_v2r8(rinvsq12,felec));
1736
1737             fscal            = felec;
1738
1739             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
1740
1741             /* Update vectorial force */
1742             fix1             = _fjsp_madd_v2r8(dx12,fscal,fix1);
1743             fiy1             = _fjsp_madd_v2r8(dy12,fscal,fiy1);
1744             fiz1             = _fjsp_madd_v2r8(dz12,fscal,fiz1);
1745             
1746             fjx2             = _fjsp_madd_v2r8(dx12,fscal,fjx2);
1747             fjy2             = _fjsp_madd_v2r8(dy12,fscal,fjy2);
1748             fjz2             = _fjsp_madd_v2r8(dz12,fscal,fjz2);
1749
1750             /**************************
1751              * CALCULATE INTERACTIONS *
1752              **************************/
1753
1754             r13              = _fjsp_mul_v2r8(rsq13,rinv13);
1755
1756             /* EWALD ELECTROSTATICS */
1757
1758             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
1759             ewrt             = _fjsp_mul_v2r8(r13,ewtabscale);
1760             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
1761             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
1762             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
1763
1764             gmx_fjsp_load_1pair_swizzle_v2r8(ewtab+ewconv.i[0],&ewtabF,&ewtabFn);
1765             felec            = _fjsp_madd_v2r8(eweps,ewtabFn,_fjsp_nmsub_v2r8(eweps,ewtabF,ewtabF));
1766             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq13,rinv13),_fjsp_sub_v2r8(rinvsq13,felec));
1767
1768             fscal            = felec;
1769
1770             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
1771
1772             /* Update vectorial force */
1773             fix1             = _fjsp_madd_v2r8(dx13,fscal,fix1);
1774             fiy1             = _fjsp_madd_v2r8(dy13,fscal,fiy1);
1775             fiz1             = _fjsp_madd_v2r8(dz13,fscal,fiz1);
1776             
1777             fjx3             = _fjsp_madd_v2r8(dx13,fscal,fjx3);
1778             fjy3             = _fjsp_madd_v2r8(dy13,fscal,fjy3);
1779             fjz3             = _fjsp_madd_v2r8(dz13,fscal,fjz3);
1780
1781             /**************************
1782              * CALCULATE INTERACTIONS *
1783              **************************/
1784
1785             r21              = _fjsp_mul_v2r8(rsq21,rinv21);
1786
1787             /* EWALD ELECTROSTATICS */
1788
1789             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
1790             ewrt             = _fjsp_mul_v2r8(r21,ewtabscale);
1791             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
1792             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
1793             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
1794
1795             gmx_fjsp_load_1pair_swizzle_v2r8(ewtab+ewconv.i[0],&ewtabF,&ewtabFn);
1796             felec            = _fjsp_madd_v2r8(eweps,ewtabFn,_fjsp_nmsub_v2r8(eweps,ewtabF,ewtabF));
1797             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq21,rinv21),_fjsp_sub_v2r8(rinvsq21,felec));
1798
1799             fscal            = felec;
1800
1801             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
1802
1803             /* Update vectorial force */
1804             fix2             = _fjsp_madd_v2r8(dx21,fscal,fix2);
1805             fiy2             = _fjsp_madd_v2r8(dy21,fscal,fiy2);
1806             fiz2             = _fjsp_madd_v2r8(dz21,fscal,fiz2);
1807             
1808             fjx1             = _fjsp_madd_v2r8(dx21,fscal,fjx1);
1809             fjy1             = _fjsp_madd_v2r8(dy21,fscal,fjy1);
1810             fjz1             = _fjsp_madd_v2r8(dz21,fscal,fjz1);
1811
1812             /**************************
1813              * CALCULATE INTERACTIONS *
1814              **************************/
1815
1816             r22              = _fjsp_mul_v2r8(rsq22,rinv22);
1817
1818             /* EWALD ELECTROSTATICS */
1819
1820             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
1821             ewrt             = _fjsp_mul_v2r8(r22,ewtabscale);
1822             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
1823             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
1824             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
1825
1826             gmx_fjsp_load_1pair_swizzle_v2r8(ewtab+ewconv.i[0],&ewtabF,&ewtabFn);
1827             felec            = _fjsp_madd_v2r8(eweps,ewtabFn,_fjsp_nmsub_v2r8(eweps,ewtabF,ewtabF));
1828             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq22,rinv22),_fjsp_sub_v2r8(rinvsq22,felec));
1829
1830             fscal            = felec;
1831
1832             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
1833
1834             /* Update vectorial force */
1835             fix2             = _fjsp_madd_v2r8(dx22,fscal,fix2);
1836             fiy2             = _fjsp_madd_v2r8(dy22,fscal,fiy2);
1837             fiz2             = _fjsp_madd_v2r8(dz22,fscal,fiz2);
1838             
1839             fjx2             = _fjsp_madd_v2r8(dx22,fscal,fjx2);
1840             fjy2             = _fjsp_madd_v2r8(dy22,fscal,fjy2);
1841             fjz2             = _fjsp_madd_v2r8(dz22,fscal,fjz2);
1842
1843             /**************************
1844              * CALCULATE INTERACTIONS *
1845              **************************/
1846
1847             r23              = _fjsp_mul_v2r8(rsq23,rinv23);
1848
1849             /* EWALD ELECTROSTATICS */
1850
1851             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
1852             ewrt             = _fjsp_mul_v2r8(r23,ewtabscale);
1853             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
1854             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
1855             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
1856
1857             gmx_fjsp_load_1pair_swizzle_v2r8(ewtab+ewconv.i[0],&ewtabF,&ewtabFn);
1858             felec            = _fjsp_madd_v2r8(eweps,ewtabFn,_fjsp_nmsub_v2r8(eweps,ewtabF,ewtabF));
1859             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq23,rinv23),_fjsp_sub_v2r8(rinvsq23,felec));
1860
1861             fscal            = felec;
1862
1863             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
1864
1865             /* Update vectorial force */
1866             fix2             = _fjsp_madd_v2r8(dx23,fscal,fix2);
1867             fiy2             = _fjsp_madd_v2r8(dy23,fscal,fiy2);
1868             fiz2             = _fjsp_madd_v2r8(dz23,fscal,fiz2);
1869             
1870             fjx3             = _fjsp_madd_v2r8(dx23,fscal,fjx3);
1871             fjy3             = _fjsp_madd_v2r8(dy23,fscal,fjy3);
1872             fjz3             = _fjsp_madd_v2r8(dz23,fscal,fjz3);
1873
1874             /**************************
1875              * CALCULATE INTERACTIONS *
1876              **************************/
1877
1878             r31              = _fjsp_mul_v2r8(rsq31,rinv31);
1879
1880             /* EWALD ELECTROSTATICS */
1881
1882             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
1883             ewrt             = _fjsp_mul_v2r8(r31,ewtabscale);
1884             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
1885             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
1886             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
1887
1888             gmx_fjsp_load_1pair_swizzle_v2r8(ewtab+ewconv.i[0],&ewtabF,&ewtabFn);
1889             felec            = _fjsp_madd_v2r8(eweps,ewtabFn,_fjsp_nmsub_v2r8(eweps,ewtabF,ewtabF));
1890             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq31,rinv31),_fjsp_sub_v2r8(rinvsq31,felec));
1891
1892             fscal            = felec;
1893
1894             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
1895
1896             /* Update vectorial force */
1897             fix3             = _fjsp_madd_v2r8(dx31,fscal,fix3);
1898             fiy3             = _fjsp_madd_v2r8(dy31,fscal,fiy3);
1899             fiz3             = _fjsp_madd_v2r8(dz31,fscal,fiz3);
1900             
1901             fjx1             = _fjsp_madd_v2r8(dx31,fscal,fjx1);
1902             fjy1             = _fjsp_madd_v2r8(dy31,fscal,fjy1);
1903             fjz1             = _fjsp_madd_v2r8(dz31,fscal,fjz1);
1904
1905             /**************************
1906              * CALCULATE INTERACTIONS *
1907              **************************/
1908
1909             r32              = _fjsp_mul_v2r8(rsq32,rinv32);
1910
1911             /* EWALD ELECTROSTATICS */
1912
1913             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
1914             ewrt             = _fjsp_mul_v2r8(r32,ewtabscale);
1915             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
1916             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
1917             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
1918
1919             gmx_fjsp_load_1pair_swizzle_v2r8(ewtab+ewconv.i[0],&ewtabF,&ewtabFn);
1920             felec            = _fjsp_madd_v2r8(eweps,ewtabFn,_fjsp_nmsub_v2r8(eweps,ewtabF,ewtabF));
1921             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq32,rinv32),_fjsp_sub_v2r8(rinvsq32,felec));
1922
1923             fscal            = felec;
1924
1925             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
1926
1927             /* Update vectorial force */
1928             fix3             = _fjsp_madd_v2r8(dx32,fscal,fix3);
1929             fiy3             = _fjsp_madd_v2r8(dy32,fscal,fiy3);
1930             fiz3             = _fjsp_madd_v2r8(dz32,fscal,fiz3);
1931             
1932             fjx2             = _fjsp_madd_v2r8(dx32,fscal,fjx2);
1933             fjy2             = _fjsp_madd_v2r8(dy32,fscal,fjy2);
1934             fjz2             = _fjsp_madd_v2r8(dz32,fscal,fjz2);
1935
1936             /**************************
1937              * CALCULATE INTERACTIONS *
1938              **************************/
1939
1940             r33              = _fjsp_mul_v2r8(rsq33,rinv33);
1941
1942             /* EWALD ELECTROSTATICS */
1943
1944             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
1945             ewrt             = _fjsp_mul_v2r8(r33,ewtabscale);
1946             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
1947             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
1948             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
1949
1950             gmx_fjsp_load_1pair_swizzle_v2r8(ewtab+ewconv.i[0],&ewtabF,&ewtabFn);
1951             felec            = _fjsp_madd_v2r8(eweps,ewtabFn,_fjsp_nmsub_v2r8(eweps,ewtabF,ewtabF));
1952             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq33,rinv33),_fjsp_sub_v2r8(rinvsq33,felec));
1953
1954             fscal            = felec;
1955
1956             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
1957
1958             /* Update vectorial force */
1959             fix3             = _fjsp_madd_v2r8(dx33,fscal,fix3);
1960             fiy3             = _fjsp_madd_v2r8(dy33,fscal,fiy3);
1961             fiz3             = _fjsp_madd_v2r8(dz33,fscal,fiz3);
1962             
1963             fjx3             = _fjsp_madd_v2r8(dx33,fscal,fjx3);
1964             fjy3             = _fjsp_madd_v2r8(dy33,fscal,fjy3);
1965             fjz3             = _fjsp_madd_v2r8(dz33,fscal,fjz3);
1966
1967             gmx_fjsp_decrement_3rvec_1ptr_swizzle_v2r8(f+j_coord_offsetA+DIM,fjx1,fjy1,fjz1,fjx2,fjy2,fjz2,fjx3,fjy3,fjz3);
1968
1969             /* Inner loop uses 351 flops */
1970         }
1971
1972         /* End of innermost loop */
1973
1974         gmx_fjsp_update_iforce_3atom_swizzle_v2r8(fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
1975                                               f+i_coord_offset+DIM,fshift+i_shift_offset);
1976
1977         /* Increment number of inner iterations */
1978         inneriter                  += j_index_end - j_index_start;
1979
1980         /* Outer loop uses 18 flops */
1981     }
1982
1983     /* Increment number of outer iterations */
1984     outeriter        += nri;
1985
1986     /* Update outer/inner flops */
1987
1988     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_W4W4_F,outeriter*18 + inneriter*351);
1989 }