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