eca058381a18d3d1ea803605f17d4ee5e08b8360
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_sparc64_hpc_ace_double / nb_kernel_ElecEwSh_VdwLJEwSh_GeomW3W3_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_ElecEwSh_VdwLJEwSh_GeomW3W3_VF_sparc64_hpc_ace_double
51  * Electrostatics interaction: Ewald
52  * VdW interaction:            LJEwald
53  * Geometry:                   Water3-Water3
54  * Calculate force/pot:        PotentialAndForce
55  */
56 void
57 nb_kernel_ElecEwSh_VdwLJEwSh_GeomW3W3_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              vdwjidx0A,vdwjidx0B;
86     _fjsp_v2r8       jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
87     int              vdwjidx1A,vdwjidx1B;
88     _fjsp_v2r8       jx1,jy1,jz1,fjx1,fjy1,fjz1,jq1,isaj1;
89     int              vdwjidx2A,vdwjidx2B;
90     _fjsp_v2r8       jx2,jy2,jz2,fjx2,fjy2,fjz2,jq2,isaj2;
91     _fjsp_v2r8       dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
92     _fjsp_v2r8       dx01,dy01,dz01,rsq01,rinv01,rinvsq01,r01,qq01,c6_01,c12_01;
93     _fjsp_v2r8       dx02,dy02,dz02,rsq02,rinv02,rinvsq02,r02,qq02,c6_02,c12_02;
94     _fjsp_v2r8       dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
95     _fjsp_v2r8       dx11,dy11,dz11,rsq11,rinv11,rinvsq11,r11,qq11,c6_11,c12_11;
96     _fjsp_v2r8       dx12,dy12,dz12,rsq12,rinv12,rinvsq12,r12,qq12,c6_12,c12_12;
97     _fjsp_v2r8       dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
98     _fjsp_v2r8       dx21,dy21,dz21,rsq21,rinv21,rinvsq21,r21,qq21,c6_21,c12_21;
99     _fjsp_v2r8       dx22,dy22,dz22,rsq22,rinv22,rinvsq22,r22,qq22,c6_22,c12_22;
100     _fjsp_v2r8       velec,felec,velecsum,facel,crf,krf,krf2;
101     real             *charge;
102     int              nvdwtype;
103     _fjsp_v2r8       rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
104     int              *vdwtype;
105     real             *vdwparam;
106     _fjsp_v2r8       one_sixth   = gmx_fjsp_set1_v2r8(1.0/6.0);
107     _fjsp_v2r8       one_twelfth = gmx_fjsp_set1_v2r8(1.0/12.0);
108     _fjsp_v2r8           c6grid_00;
109     _fjsp_v2r8           c6grid_01;
110     _fjsp_v2r8           c6grid_02;
111     _fjsp_v2r8           c6grid_10;
112     _fjsp_v2r8           c6grid_11;
113     _fjsp_v2r8           c6grid_12;
114     _fjsp_v2r8           c6grid_20;
115     _fjsp_v2r8           c6grid_21;
116     _fjsp_v2r8           c6grid_22;
117     real                 *vdwgridparam;
118     _fjsp_v2r8           ewclj,ewclj2,ewclj6,ewcljrsq,poly,exponent,f6A,f6B,sh_lj_ewald;
119     _fjsp_v2r8           one_half = gmx_fjsp_set1_v2r8(0.5);
120     _fjsp_v2r8           minus_one = gmx_fjsp_set1_v2r8(-1.0);
121     _fjsp_v2r8       ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
122     real             *ewtab;
123     _fjsp_v2r8       itab_tmp;
124     _fjsp_v2r8       dummy_mask,cutoff_mask;
125     _fjsp_v2r8       one     = gmx_fjsp_set1_v2r8(1.0);
126     _fjsp_v2r8       two     = gmx_fjsp_set1_v2r8(2.0);
127     union { _fjsp_v2r8 simd; long long int i[2]; } vfconv,gbconv,ewconv;
128
129     x                = xx[0];
130     f                = ff[0];
131
132     nri              = nlist->nri;
133     iinr             = nlist->iinr;
134     jindex           = nlist->jindex;
135     jjnr             = nlist->jjnr;
136     shiftidx         = nlist->shift;
137     gid              = nlist->gid;
138     shiftvec         = fr->shift_vec[0];
139     fshift           = fr->fshift[0];
140     facel            = gmx_fjsp_set1_v2r8(fr->epsfac);
141     charge           = mdatoms->chargeA;
142     nvdwtype         = fr->ntype;
143     vdwparam         = fr->nbfp;
144     vdwtype          = mdatoms->typeA;
145     vdwgridparam     = fr->ljpme_c6grid;
146     sh_lj_ewald      = gmx_fjsp_set1_v2r8(fr->ic->sh_lj_ewald);
147     ewclj            = gmx_fjsp_set1_v2r8(fr->ewaldcoeff_lj);
148     ewclj2           = _fjsp_mul_v2r8(minus_one,_fjsp_mul_v2r8(ewclj,ewclj));
149
150     sh_ewald         = gmx_fjsp_set1_v2r8(fr->ic->sh_ewald);
151     ewtab            = fr->ic->tabq_coul_FDV0;
152     ewtabscale       = gmx_fjsp_set1_v2r8(fr->ic->tabq_scale);
153     ewtabhalfspace   = gmx_fjsp_set1_v2r8(0.5/fr->ic->tabq_scale);
154
155     /* Setup water-specific parameters */
156     inr              = nlist->iinr[0];
157     iq0              = _fjsp_mul_v2r8(facel,gmx_fjsp_set1_v2r8(charge[inr+0]));
158     iq1              = _fjsp_mul_v2r8(facel,gmx_fjsp_set1_v2r8(charge[inr+1]));
159     iq2              = _fjsp_mul_v2r8(facel,gmx_fjsp_set1_v2r8(charge[inr+2]));
160     vdwioffset0      = 2*nvdwtype*vdwtype[inr+0];
161
162     jq0              = gmx_fjsp_set1_v2r8(charge[inr+0]);
163     jq1              = gmx_fjsp_set1_v2r8(charge[inr+1]);
164     jq2              = gmx_fjsp_set1_v2r8(charge[inr+2]);
165     vdwjidx0A        = 2*vdwtype[inr+0];
166     qq00             = _fjsp_mul_v2r8(iq0,jq0);
167     c6_00            = gmx_fjsp_set1_v2r8(vdwparam[vdwioffset0+vdwjidx0A]);
168     c12_00           = gmx_fjsp_set1_v2r8(vdwparam[vdwioffset0+vdwjidx0A+1]);
169     c6grid_00        = gmx_fjsp_set1_v2r8(vdwgridparam[vdwioffset0+vdwjidx0A]);
170     qq01             = _fjsp_mul_v2r8(iq0,jq1);
171     qq02             = _fjsp_mul_v2r8(iq0,jq2);
172     qq10             = _fjsp_mul_v2r8(iq1,jq0);
173     qq11             = _fjsp_mul_v2r8(iq1,jq1);
174     qq12             = _fjsp_mul_v2r8(iq1,jq2);
175     qq20             = _fjsp_mul_v2r8(iq2,jq0);
176     qq21             = _fjsp_mul_v2r8(iq2,jq1);
177     qq22             = _fjsp_mul_v2r8(iq2,jq2);
178
179     /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
180     rcutoff_scalar   = fr->rcoulomb;
181     rcutoff          = gmx_fjsp_set1_v2r8(rcutoff_scalar);
182     rcutoff2         = _fjsp_mul_v2r8(rcutoff,rcutoff);
183
184     sh_vdw_invrcut6  = gmx_fjsp_set1_v2r8(fr->ic->sh_invrc6);
185     rvdw             = gmx_fjsp_set1_v2r8(fr->rvdw);
186
187     /* Avoid stupid compiler warnings */
188     jnrA = jnrB = 0;
189     j_coord_offsetA = 0;
190     j_coord_offsetB = 0;
191
192     outeriter        = 0;
193     inneriter        = 0;
194
195     /* Start outer loop over neighborlists */
196     for(iidx=0; iidx<nri; iidx++)
197     {
198         /* Load shift vector for this list */
199         i_shift_offset   = DIM*shiftidx[iidx];
200
201         /* Load limits for loop over neighbors */
202         j_index_start    = jindex[iidx];
203         j_index_end      = jindex[iidx+1];
204
205         /* Get outer coordinate index */
206         inr              = iinr[iidx];
207         i_coord_offset   = DIM*inr;
208
209         /* Load i particle coords and add shift vector */
210         gmx_fjsp_load_shift_and_3rvec_broadcast_v2r8(shiftvec+i_shift_offset,x+i_coord_offset,
211                                                  &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2);
212
213         fix0             = _fjsp_setzero_v2r8();
214         fiy0             = _fjsp_setzero_v2r8();
215         fiz0             = _fjsp_setzero_v2r8();
216         fix1             = _fjsp_setzero_v2r8();
217         fiy1             = _fjsp_setzero_v2r8();
218         fiz1             = _fjsp_setzero_v2r8();
219         fix2             = _fjsp_setzero_v2r8();
220         fiy2             = _fjsp_setzero_v2r8();
221         fiz2             = _fjsp_setzero_v2r8();
222
223         /* Reset potential sums */
224         velecsum         = _fjsp_setzero_v2r8();
225         vvdwsum          = _fjsp_setzero_v2r8();
226
227         /* Start inner kernel loop */
228         for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
229         {
230
231             /* Get j neighbor index, and coordinate index */
232             jnrA             = jjnr[jidx];
233             jnrB             = jjnr[jidx+1];
234             j_coord_offsetA  = DIM*jnrA;
235             j_coord_offsetB  = DIM*jnrB;
236
237             /* load j atom coordinates */
238             gmx_fjsp_load_3rvec_2ptr_swizzle_v2r8(x+j_coord_offsetA,x+j_coord_offsetB,
239                                               &jx0,&jy0,&jz0,&jx1,&jy1,&jz1,&jx2,&jy2,&jz2);
240
241             /* Calculate displacement vector */
242             dx00             = _fjsp_sub_v2r8(ix0,jx0);
243             dy00             = _fjsp_sub_v2r8(iy0,jy0);
244             dz00             = _fjsp_sub_v2r8(iz0,jz0);
245             dx01             = _fjsp_sub_v2r8(ix0,jx1);
246             dy01             = _fjsp_sub_v2r8(iy0,jy1);
247             dz01             = _fjsp_sub_v2r8(iz0,jz1);
248             dx02             = _fjsp_sub_v2r8(ix0,jx2);
249             dy02             = _fjsp_sub_v2r8(iy0,jy2);
250             dz02             = _fjsp_sub_v2r8(iz0,jz2);
251             dx10             = _fjsp_sub_v2r8(ix1,jx0);
252             dy10             = _fjsp_sub_v2r8(iy1,jy0);
253             dz10             = _fjsp_sub_v2r8(iz1,jz0);
254             dx11             = _fjsp_sub_v2r8(ix1,jx1);
255             dy11             = _fjsp_sub_v2r8(iy1,jy1);
256             dz11             = _fjsp_sub_v2r8(iz1,jz1);
257             dx12             = _fjsp_sub_v2r8(ix1,jx2);
258             dy12             = _fjsp_sub_v2r8(iy1,jy2);
259             dz12             = _fjsp_sub_v2r8(iz1,jz2);
260             dx20             = _fjsp_sub_v2r8(ix2,jx0);
261             dy20             = _fjsp_sub_v2r8(iy2,jy0);
262             dz20             = _fjsp_sub_v2r8(iz2,jz0);
263             dx21             = _fjsp_sub_v2r8(ix2,jx1);
264             dy21             = _fjsp_sub_v2r8(iy2,jy1);
265             dz21             = _fjsp_sub_v2r8(iz2,jz1);
266             dx22             = _fjsp_sub_v2r8(ix2,jx2);
267             dy22             = _fjsp_sub_v2r8(iy2,jy2);
268             dz22             = _fjsp_sub_v2r8(iz2,jz2);
269
270             /* Calculate squared distance and things based on it */
271             rsq00            = gmx_fjsp_calc_rsq_v2r8(dx00,dy00,dz00);
272             rsq01            = gmx_fjsp_calc_rsq_v2r8(dx01,dy01,dz01);
273             rsq02            = gmx_fjsp_calc_rsq_v2r8(dx02,dy02,dz02);
274             rsq10            = gmx_fjsp_calc_rsq_v2r8(dx10,dy10,dz10);
275             rsq11            = gmx_fjsp_calc_rsq_v2r8(dx11,dy11,dz11);
276             rsq12            = gmx_fjsp_calc_rsq_v2r8(dx12,dy12,dz12);
277             rsq20            = gmx_fjsp_calc_rsq_v2r8(dx20,dy20,dz20);
278             rsq21            = gmx_fjsp_calc_rsq_v2r8(dx21,dy21,dz21);
279             rsq22            = gmx_fjsp_calc_rsq_v2r8(dx22,dy22,dz22);
280
281             rinv00           = gmx_fjsp_invsqrt_v2r8(rsq00);
282             rinv01           = gmx_fjsp_invsqrt_v2r8(rsq01);
283             rinv02           = gmx_fjsp_invsqrt_v2r8(rsq02);
284             rinv10           = gmx_fjsp_invsqrt_v2r8(rsq10);
285             rinv11           = gmx_fjsp_invsqrt_v2r8(rsq11);
286             rinv12           = gmx_fjsp_invsqrt_v2r8(rsq12);
287             rinv20           = gmx_fjsp_invsqrt_v2r8(rsq20);
288             rinv21           = gmx_fjsp_invsqrt_v2r8(rsq21);
289             rinv22           = gmx_fjsp_invsqrt_v2r8(rsq22);
290
291             rinvsq00         = _fjsp_mul_v2r8(rinv00,rinv00);
292             rinvsq01         = _fjsp_mul_v2r8(rinv01,rinv01);
293             rinvsq02         = _fjsp_mul_v2r8(rinv02,rinv02);
294             rinvsq10         = _fjsp_mul_v2r8(rinv10,rinv10);
295             rinvsq11         = _fjsp_mul_v2r8(rinv11,rinv11);
296             rinvsq12         = _fjsp_mul_v2r8(rinv12,rinv12);
297             rinvsq20         = _fjsp_mul_v2r8(rinv20,rinv20);
298             rinvsq21         = _fjsp_mul_v2r8(rinv21,rinv21);
299             rinvsq22         = _fjsp_mul_v2r8(rinv22,rinv22);
300
301             fjx0             = _fjsp_setzero_v2r8();
302             fjy0             = _fjsp_setzero_v2r8();
303             fjz0             = _fjsp_setzero_v2r8();
304             fjx1             = _fjsp_setzero_v2r8();
305             fjy1             = _fjsp_setzero_v2r8();
306             fjz1             = _fjsp_setzero_v2r8();
307             fjx2             = _fjsp_setzero_v2r8();
308             fjy2             = _fjsp_setzero_v2r8();
309             fjz2             = _fjsp_setzero_v2r8();
310
311             /**************************
312              * CALCULATE INTERACTIONS *
313              **************************/
314
315             if (gmx_fjsp_any_lt_v2r8(rsq00,rcutoff2))
316             {
317
318             r00              = _fjsp_mul_v2r8(rsq00,rinv00);
319
320             /* EWALD ELECTROSTATICS */
321
322             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
323             ewrt             = _fjsp_mul_v2r8(r00,ewtabscale);
324             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
325             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
326             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
327
328             ewtabF           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[0] );
329             ewtabD           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[1] );
330             GMX_FJSP_TRANSPOSE2_V2R8(ewtabF,ewtabD);
331             ewtabV           = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[0] +2);
332             ewtabFn          = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[1] +2);
333             GMX_FJSP_TRANSPOSE2_V2R8(ewtabV,ewtabFn);
334             felec            = _fjsp_madd_v2r8(eweps,ewtabD,ewtabF);
335             velec            = _fjsp_nmsub_v2r8(_fjsp_mul_v2r8(ewtabhalfspace,eweps) ,_fjsp_add_v2r8(ewtabF,felec), ewtabV);
336             velec            = _fjsp_mul_v2r8(qq00,_fjsp_sub_v2r8(_fjsp_sub_v2r8(rinv00,sh_ewald),velec));
337             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq00,rinv00),_fjsp_sub_v2r8(rinvsq00,felec));
338
339             /* Analytical LJ-PME */
340             rinvsix          = _fjsp_mul_v2r8(_fjsp_mul_v2r8(rinvsq00,rinvsq00),rinvsq00);
341             ewcljrsq         = _fjsp_mul_v2r8(ewclj2,rsq00);
342             ewclj6           = _fjsp_mul_v2r8(ewclj2,_fjsp_mul_v2r8(ewclj2,ewclj2));
343             exponent         = gmx_simd_exp_d(-ewcljrsq);
344             /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
345             poly             = _fjsp_mul_v2r8(exponent,_fjsp_madd_v2r8(_fjsp_mul_v2r8(ewcljrsq,ewcljrsq),one_half,_fjsp_sub_v2r8(one,ewcljrsq)));
346             /* vvdw6 = [C6 - C6grid * (1-poly)]/r6 */
347             vvdw6            = _fjsp_mul_v2r8(_fjsp_madd_v2r8(-c6grid_00,_fjsp_sub_v2r8(one,poly),c6_00),rinvsix);
348             vvdw12           = _fjsp_mul_v2r8(c12_00,_fjsp_mul_v2r8(rinvsix,rinvsix));
349             vvdw             = _fjsp_msub_v2r8(_fjsp_nmsub_v2r8(c12_00,_fjsp_mul_v2r8(sh_vdw_invrcut6,sh_vdw_invrcut6),vvdw12),one_twelfth,
350                                _fjsp_mul_v2r8(_fjsp_sub_v2r8(vvdw6,_fjsp_madd_v2r8(c6grid_00,sh_lj_ewald,_fjsp_mul_v2r8(c6_00,sh_vdw_invrcut6))),one_sixth));
351             /* fvdw = vvdw12/r - (vvdw6/r + (C6grid * exponent * beta^6)/r) */
352             fvdw             = _fjsp_mul_v2r8(_fjsp_add_v2r8(vvdw12,_fjsp_msub_v2r8(_fjsp_mul_v2r8(c6grid_00,one_sixth),_fjsp_mul_v2r8(exponent,ewclj6),vvdw6)),rinvsq00);
353
354             cutoff_mask      = _fjsp_cmplt_v2r8(rsq00,rcutoff2);
355
356             /* Update potential sum for this i atom from the interaction with this j atom. */
357             velec            = _fjsp_and_v2r8(velec,cutoff_mask);
358             velecsum         = _fjsp_add_v2r8(velecsum,velec);
359             vvdw             = _fjsp_and_v2r8(vvdw,cutoff_mask);
360             vvdwsum          = _fjsp_add_v2r8(vvdwsum,vvdw);
361
362             fscal            = _fjsp_add_v2r8(felec,fvdw);
363
364             fscal            = _fjsp_and_v2r8(fscal,cutoff_mask);
365
366             /* Update vectorial force */
367             fix0             = _fjsp_madd_v2r8(dx00,fscal,fix0);
368             fiy0             = _fjsp_madd_v2r8(dy00,fscal,fiy0);
369             fiz0             = _fjsp_madd_v2r8(dz00,fscal,fiz0);
370             
371             fjx0             = _fjsp_madd_v2r8(dx00,fscal,fjx0);
372             fjy0             = _fjsp_madd_v2r8(dy00,fscal,fjy0);
373             fjz0             = _fjsp_madd_v2r8(dz00,fscal,fjz0);
374
375             }
376
377             /**************************
378              * CALCULATE INTERACTIONS *
379              **************************/
380
381             if (gmx_fjsp_any_lt_v2r8(rsq01,rcutoff2))
382             {
383
384             r01              = _fjsp_mul_v2r8(rsq01,rinv01);
385
386             /* EWALD ELECTROSTATICS */
387
388             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
389             ewrt             = _fjsp_mul_v2r8(r01,ewtabscale);
390             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
391             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
392             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
393
394             ewtabF           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[0] );
395             ewtabD           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[1] );
396             GMX_FJSP_TRANSPOSE2_V2R8(ewtabF,ewtabD);
397             ewtabV           = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[0] +2);
398             ewtabFn          = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[1] +2);
399             GMX_FJSP_TRANSPOSE2_V2R8(ewtabV,ewtabFn);
400             felec            = _fjsp_madd_v2r8(eweps,ewtabD,ewtabF);
401             velec            = _fjsp_nmsub_v2r8(_fjsp_mul_v2r8(ewtabhalfspace,eweps) ,_fjsp_add_v2r8(ewtabF,felec), ewtabV);
402             velec            = _fjsp_mul_v2r8(qq01,_fjsp_sub_v2r8(_fjsp_sub_v2r8(rinv01,sh_ewald),velec));
403             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq01,rinv01),_fjsp_sub_v2r8(rinvsq01,felec));
404
405             cutoff_mask      = _fjsp_cmplt_v2r8(rsq01,rcutoff2);
406
407             /* Update potential sum for this i atom from the interaction with this j atom. */
408             velec            = _fjsp_and_v2r8(velec,cutoff_mask);
409             velecsum         = _fjsp_add_v2r8(velecsum,velec);
410
411             fscal            = felec;
412
413             fscal            = _fjsp_and_v2r8(fscal,cutoff_mask);
414
415             /* Update vectorial force */
416             fix0             = _fjsp_madd_v2r8(dx01,fscal,fix0);
417             fiy0             = _fjsp_madd_v2r8(dy01,fscal,fiy0);
418             fiz0             = _fjsp_madd_v2r8(dz01,fscal,fiz0);
419             
420             fjx1             = _fjsp_madd_v2r8(dx01,fscal,fjx1);
421             fjy1             = _fjsp_madd_v2r8(dy01,fscal,fjy1);
422             fjz1             = _fjsp_madd_v2r8(dz01,fscal,fjz1);
423
424             }
425
426             /**************************
427              * CALCULATE INTERACTIONS *
428              **************************/
429
430             if (gmx_fjsp_any_lt_v2r8(rsq02,rcutoff2))
431             {
432
433             r02              = _fjsp_mul_v2r8(rsq02,rinv02);
434
435             /* EWALD ELECTROSTATICS */
436
437             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
438             ewrt             = _fjsp_mul_v2r8(r02,ewtabscale);
439             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
440             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
441             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
442
443             ewtabF           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[0] );
444             ewtabD           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[1] );
445             GMX_FJSP_TRANSPOSE2_V2R8(ewtabF,ewtabD);
446             ewtabV           = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[0] +2);
447             ewtabFn          = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[1] +2);
448             GMX_FJSP_TRANSPOSE2_V2R8(ewtabV,ewtabFn);
449             felec            = _fjsp_madd_v2r8(eweps,ewtabD,ewtabF);
450             velec            = _fjsp_nmsub_v2r8(_fjsp_mul_v2r8(ewtabhalfspace,eweps) ,_fjsp_add_v2r8(ewtabF,felec), ewtabV);
451             velec            = _fjsp_mul_v2r8(qq02,_fjsp_sub_v2r8(_fjsp_sub_v2r8(rinv02,sh_ewald),velec));
452             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq02,rinv02),_fjsp_sub_v2r8(rinvsq02,felec));
453
454             cutoff_mask      = _fjsp_cmplt_v2r8(rsq02,rcutoff2);
455
456             /* Update potential sum for this i atom from the interaction with this j atom. */
457             velec            = _fjsp_and_v2r8(velec,cutoff_mask);
458             velecsum         = _fjsp_add_v2r8(velecsum,velec);
459
460             fscal            = felec;
461
462             fscal            = _fjsp_and_v2r8(fscal,cutoff_mask);
463
464             /* Update vectorial force */
465             fix0             = _fjsp_madd_v2r8(dx02,fscal,fix0);
466             fiy0             = _fjsp_madd_v2r8(dy02,fscal,fiy0);
467             fiz0             = _fjsp_madd_v2r8(dz02,fscal,fiz0);
468             
469             fjx2             = _fjsp_madd_v2r8(dx02,fscal,fjx2);
470             fjy2             = _fjsp_madd_v2r8(dy02,fscal,fjy2);
471             fjz2             = _fjsp_madd_v2r8(dz02,fscal,fjz2);
472
473             }
474
475             /**************************
476              * CALCULATE INTERACTIONS *
477              **************************/
478
479             if (gmx_fjsp_any_lt_v2r8(rsq10,rcutoff2))
480             {
481
482             r10              = _fjsp_mul_v2r8(rsq10,rinv10);
483
484             /* EWALD ELECTROSTATICS */
485
486             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
487             ewrt             = _fjsp_mul_v2r8(r10,ewtabscale);
488             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
489             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
490             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
491
492             ewtabF           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[0] );
493             ewtabD           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[1] );
494             GMX_FJSP_TRANSPOSE2_V2R8(ewtabF,ewtabD);
495             ewtabV           = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[0] +2);
496             ewtabFn          = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[1] +2);
497             GMX_FJSP_TRANSPOSE2_V2R8(ewtabV,ewtabFn);
498             felec            = _fjsp_madd_v2r8(eweps,ewtabD,ewtabF);
499             velec            = _fjsp_nmsub_v2r8(_fjsp_mul_v2r8(ewtabhalfspace,eweps) ,_fjsp_add_v2r8(ewtabF,felec), ewtabV);
500             velec            = _fjsp_mul_v2r8(qq10,_fjsp_sub_v2r8(_fjsp_sub_v2r8(rinv10,sh_ewald),velec));
501             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq10,rinv10),_fjsp_sub_v2r8(rinvsq10,felec));
502
503             cutoff_mask      = _fjsp_cmplt_v2r8(rsq10,rcutoff2);
504
505             /* Update potential sum for this i atom from the interaction with this j atom. */
506             velec            = _fjsp_and_v2r8(velec,cutoff_mask);
507             velecsum         = _fjsp_add_v2r8(velecsum,velec);
508
509             fscal            = felec;
510
511             fscal            = _fjsp_and_v2r8(fscal,cutoff_mask);
512
513             /* Update vectorial force */
514             fix1             = _fjsp_madd_v2r8(dx10,fscal,fix1);
515             fiy1             = _fjsp_madd_v2r8(dy10,fscal,fiy1);
516             fiz1             = _fjsp_madd_v2r8(dz10,fscal,fiz1);
517             
518             fjx0             = _fjsp_madd_v2r8(dx10,fscal,fjx0);
519             fjy0             = _fjsp_madd_v2r8(dy10,fscal,fjy0);
520             fjz0             = _fjsp_madd_v2r8(dz10,fscal,fjz0);
521
522             }
523
524             /**************************
525              * CALCULATE INTERACTIONS *
526              **************************/
527
528             if (gmx_fjsp_any_lt_v2r8(rsq11,rcutoff2))
529             {
530
531             r11              = _fjsp_mul_v2r8(rsq11,rinv11);
532
533             /* EWALD ELECTROSTATICS */
534
535             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
536             ewrt             = _fjsp_mul_v2r8(r11,ewtabscale);
537             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
538             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
539             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
540
541             ewtabF           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[0] );
542             ewtabD           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[1] );
543             GMX_FJSP_TRANSPOSE2_V2R8(ewtabF,ewtabD);
544             ewtabV           = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[0] +2);
545             ewtabFn          = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[1] +2);
546             GMX_FJSP_TRANSPOSE2_V2R8(ewtabV,ewtabFn);
547             felec            = _fjsp_madd_v2r8(eweps,ewtabD,ewtabF);
548             velec            = _fjsp_nmsub_v2r8(_fjsp_mul_v2r8(ewtabhalfspace,eweps) ,_fjsp_add_v2r8(ewtabF,felec), ewtabV);
549             velec            = _fjsp_mul_v2r8(qq11,_fjsp_sub_v2r8(_fjsp_sub_v2r8(rinv11,sh_ewald),velec));
550             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq11,rinv11),_fjsp_sub_v2r8(rinvsq11,felec));
551
552             cutoff_mask      = _fjsp_cmplt_v2r8(rsq11,rcutoff2);
553
554             /* Update potential sum for this i atom from the interaction with this j atom. */
555             velec            = _fjsp_and_v2r8(velec,cutoff_mask);
556             velecsum         = _fjsp_add_v2r8(velecsum,velec);
557
558             fscal            = felec;
559
560             fscal            = _fjsp_and_v2r8(fscal,cutoff_mask);
561
562             /* Update vectorial force */
563             fix1             = _fjsp_madd_v2r8(dx11,fscal,fix1);
564             fiy1             = _fjsp_madd_v2r8(dy11,fscal,fiy1);
565             fiz1             = _fjsp_madd_v2r8(dz11,fscal,fiz1);
566             
567             fjx1             = _fjsp_madd_v2r8(dx11,fscal,fjx1);
568             fjy1             = _fjsp_madd_v2r8(dy11,fscal,fjy1);
569             fjz1             = _fjsp_madd_v2r8(dz11,fscal,fjz1);
570
571             }
572
573             /**************************
574              * CALCULATE INTERACTIONS *
575              **************************/
576
577             if (gmx_fjsp_any_lt_v2r8(rsq12,rcutoff2))
578             {
579
580             r12              = _fjsp_mul_v2r8(rsq12,rinv12);
581
582             /* EWALD ELECTROSTATICS */
583
584             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
585             ewrt             = _fjsp_mul_v2r8(r12,ewtabscale);
586             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
587             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
588             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
589
590             ewtabF           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[0] );
591             ewtabD           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[1] );
592             GMX_FJSP_TRANSPOSE2_V2R8(ewtabF,ewtabD);
593             ewtabV           = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[0] +2);
594             ewtabFn          = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[1] +2);
595             GMX_FJSP_TRANSPOSE2_V2R8(ewtabV,ewtabFn);
596             felec            = _fjsp_madd_v2r8(eweps,ewtabD,ewtabF);
597             velec            = _fjsp_nmsub_v2r8(_fjsp_mul_v2r8(ewtabhalfspace,eweps) ,_fjsp_add_v2r8(ewtabF,felec), ewtabV);
598             velec            = _fjsp_mul_v2r8(qq12,_fjsp_sub_v2r8(_fjsp_sub_v2r8(rinv12,sh_ewald),velec));
599             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq12,rinv12),_fjsp_sub_v2r8(rinvsq12,felec));
600
601             cutoff_mask      = _fjsp_cmplt_v2r8(rsq12,rcutoff2);
602
603             /* Update potential sum for this i atom from the interaction with this j atom. */
604             velec            = _fjsp_and_v2r8(velec,cutoff_mask);
605             velecsum         = _fjsp_add_v2r8(velecsum,velec);
606
607             fscal            = felec;
608
609             fscal            = _fjsp_and_v2r8(fscal,cutoff_mask);
610
611             /* Update vectorial force */
612             fix1             = _fjsp_madd_v2r8(dx12,fscal,fix1);
613             fiy1             = _fjsp_madd_v2r8(dy12,fscal,fiy1);
614             fiz1             = _fjsp_madd_v2r8(dz12,fscal,fiz1);
615             
616             fjx2             = _fjsp_madd_v2r8(dx12,fscal,fjx2);
617             fjy2             = _fjsp_madd_v2r8(dy12,fscal,fjy2);
618             fjz2             = _fjsp_madd_v2r8(dz12,fscal,fjz2);
619
620             }
621
622             /**************************
623              * CALCULATE INTERACTIONS *
624              **************************/
625
626             if (gmx_fjsp_any_lt_v2r8(rsq20,rcutoff2))
627             {
628
629             r20              = _fjsp_mul_v2r8(rsq20,rinv20);
630
631             /* EWALD ELECTROSTATICS */
632
633             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
634             ewrt             = _fjsp_mul_v2r8(r20,ewtabscale);
635             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
636             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
637             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
638
639             ewtabF           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[0] );
640             ewtabD           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[1] );
641             GMX_FJSP_TRANSPOSE2_V2R8(ewtabF,ewtabD);
642             ewtabV           = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[0] +2);
643             ewtabFn          = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[1] +2);
644             GMX_FJSP_TRANSPOSE2_V2R8(ewtabV,ewtabFn);
645             felec            = _fjsp_madd_v2r8(eweps,ewtabD,ewtabF);
646             velec            = _fjsp_nmsub_v2r8(_fjsp_mul_v2r8(ewtabhalfspace,eweps) ,_fjsp_add_v2r8(ewtabF,felec), ewtabV);
647             velec            = _fjsp_mul_v2r8(qq20,_fjsp_sub_v2r8(_fjsp_sub_v2r8(rinv20,sh_ewald),velec));
648             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq20,rinv20),_fjsp_sub_v2r8(rinvsq20,felec));
649
650             cutoff_mask      = _fjsp_cmplt_v2r8(rsq20,rcutoff2);
651
652             /* Update potential sum for this i atom from the interaction with this j atom. */
653             velec            = _fjsp_and_v2r8(velec,cutoff_mask);
654             velecsum         = _fjsp_add_v2r8(velecsum,velec);
655
656             fscal            = felec;
657
658             fscal            = _fjsp_and_v2r8(fscal,cutoff_mask);
659
660             /* Update vectorial force */
661             fix2             = _fjsp_madd_v2r8(dx20,fscal,fix2);
662             fiy2             = _fjsp_madd_v2r8(dy20,fscal,fiy2);
663             fiz2             = _fjsp_madd_v2r8(dz20,fscal,fiz2);
664             
665             fjx0             = _fjsp_madd_v2r8(dx20,fscal,fjx0);
666             fjy0             = _fjsp_madd_v2r8(dy20,fscal,fjy0);
667             fjz0             = _fjsp_madd_v2r8(dz20,fscal,fjz0);
668
669             }
670
671             /**************************
672              * CALCULATE INTERACTIONS *
673              **************************/
674
675             if (gmx_fjsp_any_lt_v2r8(rsq21,rcutoff2))
676             {
677
678             r21              = _fjsp_mul_v2r8(rsq21,rinv21);
679
680             /* EWALD ELECTROSTATICS */
681
682             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
683             ewrt             = _fjsp_mul_v2r8(r21,ewtabscale);
684             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
685             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
686             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
687
688             ewtabF           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[0] );
689             ewtabD           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[1] );
690             GMX_FJSP_TRANSPOSE2_V2R8(ewtabF,ewtabD);
691             ewtabV           = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[0] +2);
692             ewtabFn          = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[1] +2);
693             GMX_FJSP_TRANSPOSE2_V2R8(ewtabV,ewtabFn);
694             felec            = _fjsp_madd_v2r8(eweps,ewtabD,ewtabF);
695             velec            = _fjsp_nmsub_v2r8(_fjsp_mul_v2r8(ewtabhalfspace,eweps) ,_fjsp_add_v2r8(ewtabF,felec), ewtabV);
696             velec            = _fjsp_mul_v2r8(qq21,_fjsp_sub_v2r8(_fjsp_sub_v2r8(rinv21,sh_ewald),velec));
697             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq21,rinv21),_fjsp_sub_v2r8(rinvsq21,felec));
698
699             cutoff_mask      = _fjsp_cmplt_v2r8(rsq21,rcutoff2);
700
701             /* Update potential sum for this i atom from the interaction with this j atom. */
702             velec            = _fjsp_and_v2r8(velec,cutoff_mask);
703             velecsum         = _fjsp_add_v2r8(velecsum,velec);
704
705             fscal            = felec;
706
707             fscal            = _fjsp_and_v2r8(fscal,cutoff_mask);
708
709             /* Update vectorial force */
710             fix2             = _fjsp_madd_v2r8(dx21,fscal,fix2);
711             fiy2             = _fjsp_madd_v2r8(dy21,fscal,fiy2);
712             fiz2             = _fjsp_madd_v2r8(dz21,fscal,fiz2);
713             
714             fjx1             = _fjsp_madd_v2r8(dx21,fscal,fjx1);
715             fjy1             = _fjsp_madd_v2r8(dy21,fscal,fjy1);
716             fjz1             = _fjsp_madd_v2r8(dz21,fscal,fjz1);
717
718             }
719
720             /**************************
721              * CALCULATE INTERACTIONS *
722              **************************/
723
724             if (gmx_fjsp_any_lt_v2r8(rsq22,rcutoff2))
725             {
726
727             r22              = _fjsp_mul_v2r8(rsq22,rinv22);
728
729             /* EWALD ELECTROSTATICS */
730
731             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
732             ewrt             = _fjsp_mul_v2r8(r22,ewtabscale);
733             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
734             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
735             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
736
737             ewtabF           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[0] );
738             ewtabD           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[1] );
739             GMX_FJSP_TRANSPOSE2_V2R8(ewtabF,ewtabD);
740             ewtabV           = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[0] +2);
741             ewtabFn          = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[1] +2);
742             GMX_FJSP_TRANSPOSE2_V2R8(ewtabV,ewtabFn);
743             felec            = _fjsp_madd_v2r8(eweps,ewtabD,ewtabF);
744             velec            = _fjsp_nmsub_v2r8(_fjsp_mul_v2r8(ewtabhalfspace,eweps) ,_fjsp_add_v2r8(ewtabF,felec), ewtabV);
745             velec            = _fjsp_mul_v2r8(qq22,_fjsp_sub_v2r8(_fjsp_sub_v2r8(rinv22,sh_ewald),velec));
746             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq22,rinv22),_fjsp_sub_v2r8(rinvsq22,felec));
747
748             cutoff_mask      = _fjsp_cmplt_v2r8(rsq22,rcutoff2);
749
750             /* Update potential sum for this i atom from the interaction with this j atom. */
751             velec            = _fjsp_and_v2r8(velec,cutoff_mask);
752             velecsum         = _fjsp_add_v2r8(velecsum,velec);
753
754             fscal            = felec;
755
756             fscal            = _fjsp_and_v2r8(fscal,cutoff_mask);
757
758             /* Update vectorial force */
759             fix2             = _fjsp_madd_v2r8(dx22,fscal,fix2);
760             fiy2             = _fjsp_madd_v2r8(dy22,fscal,fiy2);
761             fiz2             = _fjsp_madd_v2r8(dz22,fscal,fiz2);
762             
763             fjx2             = _fjsp_madd_v2r8(dx22,fscal,fjx2);
764             fjy2             = _fjsp_madd_v2r8(dy22,fscal,fjy2);
765             fjz2             = _fjsp_madd_v2r8(dz22,fscal,fjz2);
766
767             }
768
769             gmx_fjsp_decrement_3rvec_2ptr_swizzle_v2r8(f+j_coord_offsetA,f+j_coord_offsetB,fjx0,fjy0,fjz0,fjx1,fjy1,fjz1,fjx2,fjy2,fjz2);
770
771             /* Inner loop uses 471 flops */
772         }
773
774         if(jidx<j_index_end)
775         {
776
777             jnrA             = jjnr[jidx];
778             j_coord_offsetA  = DIM*jnrA;
779
780             /* load j atom coordinates */
781             gmx_fjsp_load_3rvec_1ptr_swizzle_v2r8(x+j_coord_offsetA,
782                                               &jx0,&jy0,&jz0,&jx1,&jy1,&jz1,&jx2,&jy2,&jz2);
783
784             /* Calculate displacement vector */
785             dx00             = _fjsp_sub_v2r8(ix0,jx0);
786             dy00             = _fjsp_sub_v2r8(iy0,jy0);
787             dz00             = _fjsp_sub_v2r8(iz0,jz0);
788             dx01             = _fjsp_sub_v2r8(ix0,jx1);
789             dy01             = _fjsp_sub_v2r8(iy0,jy1);
790             dz01             = _fjsp_sub_v2r8(iz0,jz1);
791             dx02             = _fjsp_sub_v2r8(ix0,jx2);
792             dy02             = _fjsp_sub_v2r8(iy0,jy2);
793             dz02             = _fjsp_sub_v2r8(iz0,jz2);
794             dx10             = _fjsp_sub_v2r8(ix1,jx0);
795             dy10             = _fjsp_sub_v2r8(iy1,jy0);
796             dz10             = _fjsp_sub_v2r8(iz1,jz0);
797             dx11             = _fjsp_sub_v2r8(ix1,jx1);
798             dy11             = _fjsp_sub_v2r8(iy1,jy1);
799             dz11             = _fjsp_sub_v2r8(iz1,jz1);
800             dx12             = _fjsp_sub_v2r8(ix1,jx2);
801             dy12             = _fjsp_sub_v2r8(iy1,jy2);
802             dz12             = _fjsp_sub_v2r8(iz1,jz2);
803             dx20             = _fjsp_sub_v2r8(ix2,jx0);
804             dy20             = _fjsp_sub_v2r8(iy2,jy0);
805             dz20             = _fjsp_sub_v2r8(iz2,jz0);
806             dx21             = _fjsp_sub_v2r8(ix2,jx1);
807             dy21             = _fjsp_sub_v2r8(iy2,jy1);
808             dz21             = _fjsp_sub_v2r8(iz2,jz1);
809             dx22             = _fjsp_sub_v2r8(ix2,jx2);
810             dy22             = _fjsp_sub_v2r8(iy2,jy2);
811             dz22             = _fjsp_sub_v2r8(iz2,jz2);
812
813             /* Calculate squared distance and things based on it */
814             rsq00            = gmx_fjsp_calc_rsq_v2r8(dx00,dy00,dz00);
815             rsq01            = gmx_fjsp_calc_rsq_v2r8(dx01,dy01,dz01);
816             rsq02            = gmx_fjsp_calc_rsq_v2r8(dx02,dy02,dz02);
817             rsq10            = gmx_fjsp_calc_rsq_v2r8(dx10,dy10,dz10);
818             rsq11            = gmx_fjsp_calc_rsq_v2r8(dx11,dy11,dz11);
819             rsq12            = gmx_fjsp_calc_rsq_v2r8(dx12,dy12,dz12);
820             rsq20            = gmx_fjsp_calc_rsq_v2r8(dx20,dy20,dz20);
821             rsq21            = gmx_fjsp_calc_rsq_v2r8(dx21,dy21,dz21);
822             rsq22            = gmx_fjsp_calc_rsq_v2r8(dx22,dy22,dz22);
823
824             rinv00           = gmx_fjsp_invsqrt_v2r8(rsq00);
825             rinv01           = gmx_fjsp_invsqrt_v2r8(rsq01);
826             rinv02           = gmx_fjsp_invsqrt_v2r8(rsq02);
827             rinv10           = gmx_fjsp_invsqrt_v2r8(rsq10);
828             rinv11           = gmx_fjsp_invsqrt_v2r8(rsq11);
829             rinv12           = gmx_fjsp_invsqrt_v2r8(rsq12);
830             rinv20           = gmx_fjsp_invsqrt_v2r8(rsq20);
831             rinv21           = gmx_fjsp_invsqrt_v2r8(rsq21);
832             rinv22           = gmx_fjsp_invsqrt_v2r8(rsq22);
833
834             rinvsq00         = _fjsp_mul_v2r8(rinv00,rinv00);
835             rinvsq01         = _fjsp_mul_v2r8(rinv01,rinv01);
836             rinvsq02         = _fjsp_mul_v2r8(rinv02,rinv02);
837             rinvsq10         = _fjsp_mul_v2r8(rinv10,rinv10);
838             rinvsq11         = _fjsp_mul_v2r8(rinv11,rinv11);
839             rinvsq12         = _fjsp_mul_v2r8(rinv12,rinv12);
840             rinvsq20         = _fjsp_mul_v2r8(rinv20,rinv20);
841             rinvsq21         = _fjsp_mul_v2r8(rinv21,rinv21);
842             rinvsq22         = _fjsp_mul_v2r8(rinv22,rinv22);
843
844             fjx0             = _fjsp_setzero_v2r8();
845             fjy0             = _fjsp_setzero_v2r8();
846             fjz0             = _fjsp_setzero_v2r8();
847             fjx1             = _fjsp_setzero_v2r8();
848             fjy1             = _fjsp_setzero_v2r8();
849             fjz1             = _fjsp_setzero_v2r8();
850             fjx2             = _fjsp_setzero_v2r8();
851             fjy2             = _fjsp_setzero_v2r8();
852             fjz2             = _fjsp_setzero_v2r8();
853
854             /**************************
855              * CALCULATE INTERACTIONS *
856              **************************/
857
858             if (gmx_fjsp_any_lt_v2r8(rsq00,rcutoff2))
859             {
860
861             r00              = _fjsp_mul_v2r8(rsq00,rinv00);
862
863             /* EWALD ELECTROSTATICS */
864
865             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
866             ewrt             = _fjsp_mul_v2r8(r00,ewtabscale);
867             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
868             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
869             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
870
871             ewtabF           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[0] );
872             ewtabD           = _fjsp_setzero_v2r8();
873             GMX_FJSP_TRANSPOSE2_V2R8(ewtabF,ewtabD);
874             ewtabV           = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[0] +2);
875             ewtabFn          = _fjsp_setzero_v2r8();
876             GMX_FJSP_TRANSPOSE2_V2R8(ewtabV,ewtabFn);
877             felec            = _fjsp_madd_v2r8(eweps,ewtabD,ewtabF);
878             velec            = _fjsp_nmsub_v2r8(_fjsp_mul_v2r8(ewtabhalfspace,eweps) ,_fjsp_add_v2r8(ewtabF,felec), ewtabV);
879             velec            = _fjsp_mul_v2r8(qq00,_fjsp_sub_v2r8(_fjsp_sub_v2r8(rinv00,sh_ewald),velec));
880             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq00,rinv00),_fjsp_sub_v2r8(rinvsq00,felec));
881
882             /* Analytical LJ-PME */
883             rinvsix          = _fjsp_mul_v2r8(_fjsp_mul_v2r8(rinvsq00,rinvsq00),rinvsq00);
884             ewcljrsq         = _fjsp_mul_v2r8(ewclj2,rsq00);
885             ewclj6           = _fjsp_mul_v2r8(ewclj2,_fjsp_mul_v2r8(ewclj2,ewclj2));
886             exponent         = gmx_simd_exp_d(-ewcljrsq);
887             /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
888             poly             = _fjsp_mul_v2r8(exponent,_fjsp_madd_v2r8(_fjsp_mul_v2r8(ewcljrsq,ewcljrsq),one_half,_fjsp_sub_v2r8(one,ewcljrsq)));
889             /* vvdw6 = [C6 - C6grid * (1-poly)]/r6 */
890             vvdw6            = _fjsp_mul_v2r8(_fjsp_madd_v2r8(-c6grid_00,_fjsp_sub_v2r8(one,poly),c6_00),rinvsix);
891             vvdw12           = _fjsp_mul_v2r8(c12_00,_fjsp_mul_v2r8(rinvsix,rinvsix));
892             vvdw             = _fjsp_msub_v2r8(_fjsp_nmsub_v2r8(c12_00,_fjsp_mul_v2r8(sh_vdw_invrcut6,sh_vdw_invrcut6),vvdw12),one_twelfth,
893                                _fjsp_mul_v2r8(_fjsp_sub_v2r8(vvdw6,_fjsp_madd_v2r8(c6grid_00,sh_lj_ewald,_fjsp_mul_v2r8(c6_00,sh_vdw_invrcut6))),one_sixth));
894             /* fvdw = vvdw12/r - (vvdw6/r + (C6grid * exponent * beta^6)/r) */
895             fvdw             = _fjsp_mul_v2r8(_fjsp_add_v2r8(vvdw12,_fjsp_msub_v2r8(_fjsp_mul_v2r8(c6grid_00,one_sixth),_fjsp_mul_v2r8(exponent,ewclj6),vvdw6)),rinvsq00);
896
897             cutoff_mask      = _fjsp_cmplt_v2r8(rsq00,rcutoff2);
898
899             /* Update potential sum for this i atom from the interaction with this j atom. */
900             velec            = _fjsp_and_v2r8(velec,cutoff_mask);
901             velec            = _fjsp_unpacklo_v2r8(velec,_fjsp_setzero_v2r8());
902             velecsum         = _fjsp_add_v2r8(velecsum,velec);
903             vvdw             = _fjsp_and_v2r8(vvdw,cutoff_mask);
904             vvdw             = _fjsp_unpacklo_v2r8(vvdw,_fjsp_setzero_v2r8());
905             vvdwsum          = _fjsp_add_v2r8(vvdwsum,vvdw);
906
907             fscal            = _fjsp_add_v2r8(felec,fvdw);
908
909             fscal            = _fjsp_and_v2r8(fscal,cutoff_mask);
910
911             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
912
913             /* Update vectorial force */
914             fix0             = _fjsp_madd_v2r8(dx00,fscal,fix0);
915             fiy0             = _fjsp_madd_v2r8(dy00,fscal,fiy0);
916             fiz0             = _fjsp_madd_v2r8(dz00,fscal,fiz0);
917             
918             fjx0             = _fjsp_madd_v2r8(dx00,fscal,fjx0);
919             fjy0             = _fjsp_madd_v2r8(dy00,fscal,fjy0);
920             fjz0             = _fjsp_madd_v2r8(dz00,fscal,fjz0);
921
922             }
923
924             /**************************
925              * CALCULATE INTERACTIONS *
926              **************************/
927
928             if (gmx_fjsp_any_lt_v2r8(rsq01,rcutoff2))
929             {
930
931             r01              = _fjsp_mul_v2r8(rsq01,rinv01);
932
933             /* EWALD ELECTROSTATICS */
934
935             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
936             ewrt             = _fjsp_mul_v2r8(r01,ewtabscale);
937             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
938             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
939             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
940
941             ewtabF           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[0] );
942             ewtabD           = _fjsp_setzero_v2r8();
943             GMX_FJSP_TRANSPOSE2_V2R8(ewtabF,ewtabD);
944             ewtabV           = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[0] +2);
945             ewtabFn          = _fjsp_setzero_v2r8();
946             GMX_FJSP_TRANSPOSE2_V2R8(ewtabV,ewtabFn);
947             felec            = _fjsp_madd_v2r8(eweps,ewtabD,ewtabF);
948             velec            = _fjsp_nmsub_v2r8(_fjsp_mul_v2r8(ewtabhalfspace,eweps) ,_fjsp_add_v2r8(ewtabF,felec), ewtabV);
949             velec            = _fjsp_mul_v2r8(qq01,_fjsp_sub_v2r8(_fjsp_sub_v2r8(rinv01,sh_ewald),velec));
950             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq01,rinv01),_fjsp_sub_v2r8(rinvsq01,felec));
951
952             cutoff_mask      = _fjsp_cmplt_v2r8(rsq01,rcutoff2);
953
954             /* Update potential sum for this i atom from the interaction with this j atom. */
955             velec            = _fjsp_and_v2r8(velec,cutoff_mask);
956             velec            = _fjsp_unpacklo_v2r8(velec,_fjsp_setzero_v2r8());
957             velecsum         = _fjsp_add_v2r8(velecsum,velec);
958
959             fscal            = felec;
960
961             fscal            = _fjsp_and_v2r8(fscal,cutoff_mask);
962
963             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
964
965             /* Update vectorial force */
966             fix0             = _fjsp_madd_v2r8(dx01,fscal,fix0);
967             fiy0             = _fjsp_madd_v2r8(dy01,fscal,fiy0);
968             fiz0             = _fjsp_madd_v2r8(dz01,fscal,fiz0);
969             
970             fjx1             = _fjsp_madd_v2r8(dx01,fscal,fjx1);
971             fjy1             = _fjsp_madd_v2r8(dy01,fscal,fjy1);
972             fjz1             = _fjsp_madd_v2r8(dz01,fscal,fjz1);
973
974             }
975
976             /**************************
977              * CALCULATE INTERACTIONS *
978              **************************/
979
980             if (gmx_fjsp_any_lt_v2r8(rsq02,rcutoff2))
981             {
982
983             r02              = _fjsp_mul_v2r8(rsq02,rinv02);
984
985             /* EWALD ELECTROSTATICS */
986
987             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
988             ewrt             = _fjsp_mul_v2r8(r02,ewtabscale);
989             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
990             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
991             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
992
993             ewtabF           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[0] );
994             ewtabD           = _fjsp_setzero_v2r8();
995             GMX_FJSP_TRANSPOSE2_V2R8(ewtabF,ewtabD);
996             ewtabV           = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[0] +2);
997             ewtabFn          = _fjsp_setzero_v2r8();
998             GMX_FJSP_TRANSPOSE2_V2R8(ewtabV,ewtabFn);
999             felec            = _fjsp_madd_v2r8(eweps,ewtabD,ewtabF);
1000             velec            = _fjsp_nmsub_v2r8(_fjsp_mul_v2r8(ewtabhalfspace,eweps) ,_fjsp_add_v2r8(ewtabF,felec), ewtabV);
1001             velec            = _fjsp_mul_v2r8(qq02,_fjsp_sub_v2r8(_fjsp_sub_v2r8(rinv02,sh_ewald),velec));
1002             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq02,rinv02),_fjsp_sub_v2r8(rinvsq02,felec));
1003
1004             cutoff_mask      = _fjsp_cmplt_v2r8(rsq02,rcutoff2);
1005
1006             /* Update potential sum for this i atom from the interaction with this j atom. */
1007             velec            = _fjsp_and_v2r8(velec,cutoff_mask);
1008             velec            = _fjsp_unpacklo_v2r8(velec,_fjsp_setzero_v2r8());
1009             velecsum         = _fjsp_add_v2r8(velecsum,velec);
1010
1011             fscal            = felec;
1012
1013             fscal            = _fjsp_and_v2r8(fscal,cutoff_mask);
1014
1015             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
1016
1017             /* Update vectorial force */
1018             fix0             = _fjsp_madd_v2r8(dx02,fscal,fix0);
1019             fiy0             = _fjsp_madd_v2r8(dy02,fscal,fiy0);
1020             fiz0             = _fjsp_madd_v2r8(dz02,fscal,fiz0);
1021             
1022             fjx2             = _fjsp_madd_v2r8(dx02,fscal,fjx2);
1023             fjy2             = _fjsp_madd_v2r8(dy02,fscal,fjy2);
1024             fjz2             = _fjsp_madd_v2r8(dz02,fscal,fjz2);
1025
1026             }
1027
1028             /**************************
1029              * CALCULATE INTERACTIONS *
1030              **************************/
1031
1032             if (gmx_fjsp_any_lt_v2r8(rsq10,rcutoff2))
1033             {
1034
1035             r10              = _fjsp_mul_v2r8(rsq10,rinv10);
1036
1037             /* EWALD ELECTROSTATICS */
1038
1039             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
1040             ewrt             = _fjsp_mul_v2r8(r10,ewtabscale);
1041             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
1042             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
1043             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
1044
1045             ewtabF           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[0] );
1046             ewtabD           = _fjsp_setzero_v2r8();
1047             GMX_FJSP_TRANSPOSE2_V2R8(ewtabF,ewtabD);
1048             ewtabV           = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[0] +2);
1049             ewtabFn          = _fjsp_setzero_v2r8();
1050             GMX_FJSP_TRANSPOSE2_V2R8(ewtabV,ewtabFn);
1051             felec            = _fjsp_madd_v2r8(eweps,ewtabD,ewtabF);
1052             velec            = _fjsp_nmsub_v2r8(_fjsp_mul_v2r8(ewtabhalfspace,eweps) ,_fjsp_add_v2r8(ewtabF,felec), ewtabV);
1053             velec            = _fjsp_mul_v2r8(qq10,_fjsp_sub_v2r8(_fjsp_sub_v2r8(rinv10,sh_ewald),velec));
1054             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq10,rinv10),_fjsp_sub_v2r8(rinvsq10,felec));
1055
1056             cutoff_mask      = _fjsp_cmplt_v2r8(rsq10,rcutoff2);
1057
1058             /* Update potential sum for this i atom from the interaction with this j atom. */
1059             velec            = _fjsp_and_v2r8(velec,cutoff_mask);
1060             velec            = _fjsp_unpacklo_v2r8(velec,_fjsp_setzero_v2r8());
1061             velecsum         = _fjsp_add_v2r8(velecsum,velec);
1062
1063             fscal            = felec;
1064
1065             fscal            = _fjsp_and_v2r8(fscal,cutoff_mask);
1066
1067             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
1068
1069             /* Update vectorial force */
1070             fix1             = _fjsp_madd_v2r8(dx10,fscal,fix1);
1071             fiy1             = _fjsp_madd_v2r8(dy10,fscal,fiy1);
1072             fiz1             = _fjsp_madd_v2r8(dz10,fscal,fiz1);
1073             
1074             fjx0             = _fjsp_madd_v2r8(dx10,fscal,fjx0);
1075             fjy0             = _fjsp_madd_v2r8(dy10,fscal,fjy0);
1076             fjz0             = _fjsp_madd_v2r8(dz10,fscal,fjz0);
1077
1078             }
1079
1080             /**************************
1081              * CALCULATE INTERACTIONS *
1082              **************************/
1083
1084             if (gmx_fjsp_any_lt_v2r8(rsq11,rcutoff2))
1085             {
1086
1087             r11              = _fjsp_mul_v2r8(rsq11,rinv11);
1088
1089             /* EWALD ELECTROSTATICS */
1090
1091             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
1092             ewrt             = _fjsp_mul_v2r8(r11,ewtabscale);
1093             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
1094             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
1095             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
1096
1097             ewtabF           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[0] );
1098             ewtabD           = _fjsp_setzero_v2r8();
1099             GMX_FJSP_TRANSPOSE2_V2R8(ewtabF,ewtabD);
1100             ewtabV           = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[0] +2);
1101             ewtabFn          = _fjsp_setzero_v2r8();
1102             GMX_FJSP_TRANSPOSE2_V2R8(ewtabV,ewtabFn);
1103             felec            = _fjsp_madd_v2r8(eweps,ewtabD,ewtabF);
1104             velec            = _fjsp_nmsub_v2r8(_fjsp_mul_v2r8(ewtabhalfspace,eweps) ,_fjsp_add_v2r8(ewtabF,felec), ewtabV);
1105             velec            = _fjsp_mul_v2r8(qq11,_fjsp_sub_v2r8(_fjsp_sub_v2r8(rinv11,sh_ewald),velec));
1106             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq11,rinv11),_fjsp_sub_v2r8(rinvsq11,felec));
1107
1108             cutoff_mask      = _fjsp_cmplt_v2r8(rsq11,rcutoff2);
1109
1110             /* Update potential sum for this i atom from the interaction with this j atom. */
1111             velec            = _fjsp_and_v2r8(velec,cutoff_mask);
1112             velec            = _fjsp_unpacklo_v2r8(velec,_fjsp_setzero_v2r8());
1113             velecsum         = _fjsp_add_v2r8(velecsum,velec);
1114
1115             fscal            = felec;
1116
1117             fscal            = _fjsp_and_v2r8(fscal,cutoff_mask);
1118
1119             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
1120
1121             /* Update vectorial force */
1122             fix1             = _fjsp_madd_v2r8(dx11,fscal,fix1);
1123             fiy1             = _fjsp_madd_v2r8(dy11,fscal,fiy1);
1124             fiz1             = _fjsp_madd_v2r8(dz11,fscal,fiz1);
1125             
1126             fjx1             = _fjsp_madd_v2r8(dx11,fscal,fjx1);
1127             fjy1             = _fjsp_madd_v2r8(dy11,fscal,fjy1);
1128             fjz1             = _fjsp_madd_v2r8(dz11,fscal,fjz1);
1129
1130             }
1131
1132             /**************************
1133              * CALCULATE INTERACTIONS *
1134              **************************/
1135
1136             if (gmx_fjsp_any_lt_v2r8(rsq12,rcutoff2))
1137             {
1138
1139             r12              = _fjsp_mul_v2r8(rsq12,rinv12);
1140
1141             /* EWALD ELECTROSTATICS */
1142
1143             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
1144             ewrt             = _fjsp_mul_v2r8(r12,ewtabscale);
1145             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
1146             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
1147             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
1148
1149             ewtabF           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[0] );
1150             ewtabD           = _fjsp_setzero_v2r8();
1151             GMX_FJSP_TRANSPOSE2_V2R8(ewtabF,ewtabD);
1152             ewtabV           = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[0] +2);
1153             ewtabFn          = _fjsp_setzero_v2r8();
1154             GMX_FJSP_TRANSPOSE2_V2R8(ewtabV,ewtabFn);
1155             felec            = _fjsp_madd_v2r8(eweps,ewtabD,ewtabF);
1156             velec            = _fjsp_nmsub_v2r8(_fjsp_mul_v2r8(ewtabhalfspace,eweps) ,_fjsp_add_v2r8(ewtabF,felec), ewtabV);
1157             velec            = _fjsp_mul_v2r8(qq12,_fjsp_sub_v2r8(_fjsp_sub_v2r8(rinv12,sh_ewald),velec));
1158             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq12,rinv12),_fjsp_sub_v2r8(rinvsq12,felec));
1159
1160             cutoff_mask      = _fjsp_cmplt_v2r8(rsq12,rcutoff2);
1161
1162             /* Update potential sum for this i atom from the interaction with this j atom. */
1163             velec            = _fjsp_and_v2r8(velec,cutoff_mask);
1164             velec            = _fjsp_unpacklo_v2r8(velec,_fjsp_setzero_v2r8());
1165             velecsum         = _fjsp_add_v2r8(velecsum,velec);
1166
1167             fscal            = felec;
1168
1169             fscal            = _fjsp_and_v2r8(fscal,cutoff_mask);
1170
1171             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
1172
1173             /* Update vectorial force */
1174             fix1             = _fjsp_madd_v2r8(dx12,fscal,fix1);
1175             fiy1             = _fjsp_madd_v2r8(dy12,fscal,fiy1);
1176             fiz1             = _fjsp_madd_v2r8(dz12,fscal,fiz1);
1177             
1178             fjx2             = _fjsp_madd_v2r8(dx12,fscal,fjx2);
1179             fjy2             = _fjsp_madd_v2r8(dy12,fscal,fjy2);
1180             fjz2             = _fjsp_madd_v2r8(dz12,fscal,fjz2);
1181
1182             }
1183
1184             /**************************
1185              * CALCULATE INTERACTIONS *
1186              **************************/
1187
1188             if (gmx_fjsp_any_lt_v2r8(rsq20,rcutoff2))
1189             {
1190
1191             r20              = _fjsp_mul_v2r8(rsq20,rinv20);
1192
1193             /* EWALD ELECTROSTATICS */
1194
1195             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
1196             ewrt             = _fjsp_mul_v2r8(r20,ewtabscale);
1197             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
1198             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
1199             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
1200
1201             ewtabF           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[0] );
1202             ewtabD           = _fjsp_setzero_v2r8();
1203             GMX_FJSP_TRANSPOSE2_V2R8(ewtabF,ewtabD);
1204             ewtabV           = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[0] +2);
1205             ewtabFn          = _fjsp_setzero_v2r8();
1206             GMX_FJSP_TRANSPOSE2_V2R8(ewtabV,ewtabFn);
1207             felec            = _fjsp_madd_v2r8(eweps,ewtabD,ewtabF);
1208             velec            = _fjsp_nmsub_v2r8(_fjsp_mul_v2r8(ewtabhalfspace,eweps) ,_fjsp_add_v2r8(ewtabF,felec), ewtabV);
1209             velec            = _fjsp_mul_v2r8(qq20,_fjsp_sub_v2r8(_fjsp_sub_v2r8(rinv20,sh_ewald),velec));
1210             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq20,rinv20),_fjsp_sub_v2r8(rinvsq20,felec));
1211
1212             cutoff_mask      = _fjsp_cmplt_v2r8(rsq20,rcutoff2);
1213
1214             /* Update potential sum for this i atom from the interaction with this j atom. */
1215             velec            = _fjsp_and_v2r8(velec,cutoff_mask);
1216             velec            = _fjsp_unpacklo_v2r8(velec,_fjsp_setzero_v2r8());
1217             velecsum         = _fjsp_add_v2r8(velecsum,velec);
1218
1219             fscal            = felec;
1220
1221             fscal            = _fjsp_and_v2r8(fscal,cutoff_mask);
1222
1223             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
1224
1225             /* Update vectorial force */
1226             fix2             = _fjsp_madd_v2r8(dx20,fscal,fix2);
1227             fiy2             = _fjsp_madd_v2r8(dy20,fscal,fiy2);
1228             fiz2             = _fjsp_madd_v2r8(dz20,fscal,fiz2);
1229             
1230             fjx0             = _fjsp_madd_v2r8(dx20,fscal,fjx0);
1231             fjy0             = _fjsp_madd_v2r8(dy20,fscal,fjy0);
1232             fjz0             = _fjsp_madd_v2r8(dz20,fscal,fjz0);
1233
1234             }
1235
1236             /**************************
1237              * CALCULATE INTERACTIONS *
1238              **************************/
1239
1240             if (gmx_fjsp_any_lt_v2r8(rsq21,rcutoff2))
1241             {
1242
1243             r21              = _fjsp_mul_v2r8(rsq21,rinv21);
1244
1245             /* EWALD ELECTROSTATICS */
1246
1247             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
1248             ewrt             = _fjsp_mul_v2r8(r21,ewtabscale);
1249             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
1250             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
1251             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
1252
1253             ewtabF           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[0] );
1254             ewtabD           = _fjsp_setzero_v2r8();
1255             GMX_FJSP_TRANSPOSE2_V2R8(ewtabF,ewtabD);
1256             ewtabV           = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[0] +2);
1257             ewtabFn          = _fjsp_setzero_v2r8();
1258             GMX_FJSP_TRANSPOSE2_V2R8(ewtabV,ewtabFn);
1259             felec            = _fjsp_madd_v2r8(eweps,ewtabD,ewtabF);
1260             velec            = _fjsp_nmsub_v2r8(_fjsp_mul_v2r8(ewtabhalfspace,eweps) ,_fjsp_add_v2r8(ewtabF,felec), ewtabV);
1261             velec            = _fjsp_mul_v2r8(qq21,_fjsp_sub_v2r8(_fjsp_sub_v2r8(rinv21,sh_ewald),velec));
1262             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq21,rinv21),_fjsp_sub_v2r8(rinvsq21,felec));
1263
1264             cutoff_mask      = _fjsp_cmplt_v2r8(rsq21,rcutoff2);
1265
1266             /* Update potential sum for this i atom from the interaction with this j atom. */
1267             velec            = _fjsp_and_v2r8(velec,cutoff_mask);
1268             velec            = _fjsp_unpacklo_v2r8(velec,_fjsp_setzero_v2r8());
1269             velecsum         = _fjsp_add_v2r8(velecsum,velec);
1270
1271             fscal            = felec;
1272
1273             fscal            = _fjsp_and_v2r8(fscal,cutoff_mask);
1274
1275             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
1276
1277             /* Update vectorial force */
1278             fix2             = _fjsp_madd_v2r8(dx21,fscal,fix2);
1279             fiy2             = _fjsp_madd_v2r8(dy21,fscal,fiy2);
1280             fiz2             = _fjsp_madd_v2r8(dz21,fscal,fiz2);
1281             
1282             fjx1             = _fjsp_madd_v2r8(dx21,fscal,fjx1);
1283             fjy1             = _fjsp_madd_v2r8(dy21,fscal,fjy1);
1284             fjz1             = _fjsp_madd_v2r8(dz21,fscal,fjz1);
1285
1286             }
1287
1288             /**************************
1289              * CALCULATE INTERACTIONS *
1290              **************************/
1291
1292             if (gmx_fjsp_any_lt_v2r8(rsq22,rcutoff2))
1293             {
1294
1295             r22              = _fjsp_mul_v2r8(rsq22,rinv22);
1296
1297             /* EWALD ELECTROSTATICS */
1298
1299             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
1300             ewrt             = _fjsp_mul_v2r8(r22,ewtabscale);
1301             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
1302             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
1303             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
1304
1305             ewtabF           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[0] );
1306             ewtabD           = _fjsp_setzero_v2r8();
1307             GMX_FJSP_TRANSPOSE2_V2R8(ewtabF,ewtabD);
1308             ewtabV           = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[0] +2);
1309             ewtabFn          = _fjsp_setzero_v2r8();
1310             GMX_FJSP_TRANSPOSE2_V2R8(ewtabV,ewtabFn);
1311             felec            = _fjsp_madd_v2r8(eweps,ewtabD,ewtabF);
1312             velec            = _fjsp_nmsub_v2r8(_fjsp_mul_v2r8(ewtabhalfspace,eweps) ,_fjsp_add_v2r8(ewtabF,felec), ewtabV);
1313             velec            = _fjsp_mul_v2r8(qq22,_fjsp_sub_v2r8(_fjsp_sub_v2r8(rinv22,sh_ewald),velec));
1314             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq22,rinv22),_fjsp_sub_v2r8(rinvsq22,felec));
1315
1316             cutoff_mask      = _fjsp_cmplt_v2r8(rsq22,rcutoff2);
1317
1318             /* Update potential sum for this i atom from the interaction with this j atom. */
1319             velec            = _fjsp_and_v2r8(velec,cutoff_mask);
1320             velec            = _fjsp_unpacklo_v2r8(velec,_fjsp_setzero_v2r8());
1321             velecsum         = _fjsp_add_v2r8(velecsum,velec);
1322
1323             fscal            = felec;
1324
1325             fscal            = _fjsp_and_v2r8(fscal,cutoff_mask);
1326
1327             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
1328
1329             /* Update vectorial force */
1330             fix2             = _fjsp_madd_v2r8(dx22,fscal,fix2);
1331             fiy2             = _fjsp_madd_v2r8(dy22,fscal,fiy2);
1332             fiz2             = _fjsp_madd_v2r8(dz22,fscal,fiz2);
1333             
1334             fjx2             = _fjsp_madd_v2r8(dx22,fscal,fjx2);
1335             fjy2             = _fjsp_madd_v2r8(dy22,fscal,fjy2);
1336             fjz2             = _fjsp_madd_v2r8(dz22,fscal,fjz2);
1337
1338             }
1339
1340             gmx_fjsp_decrement_3rvec_1ptr_swizzle_v2r8(f+j_coord_offsetA,fjx0,fjy0,fjz0,fjx1,fjy1,fjz1,fjx2,fjy2,fjz2);
1341
1342             /* Inner loop uses 471 flops */
1343         }
1344
1345         /* End of innermost loop */
1346
1347         gmx_fjsp_update_iforce_3atom_swizzle_v2r8(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,
1348                                               f+i_coord_offset,fshift+i_shift_offset);
1349
1350         ggid                        = gid[iidx];
1351         /* Update potential energies */
1352         gmx_fjsp_update_1pot_v2r8(velecsum,kernel_data->energygrp_elec+ggid);
1353         gmx_fjsp_update_1pot_v2r8(vvdwsum,kernel_data->energygrp_vdw+ggid);
1354
1355         /* Increment number of inner iterations */
1356         inneriter                  += j_index_end - j_index_start;
1357
1358         /* Outer loop uses 20 flops */
1359     }
1360
1361     /* Increment number of outer iterations */
1362     outeriter        += nri;
1363
1364     /* Update outer/inner flops */
1365
1366     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W3W3_VF,outeriter*20 + inneriter*471);
1367 }
1368 /*
1369  * Gromacs nonbonded kernel:   nb_kernel_ElecEwSh_VdwLJEwSh_GeomW3W3_F_sparc64_hpc_ace_double
1370  * Electrostatics interaction: Ewald
1371  * VdW interaction:            LJEwald
1372  * Geometry:                   Water3-Water3
1373  * Calculate force/pot:        Force
1374  */
1375 void
1376 nb_kernel_ElecEwSh_VdwLJEwSh_GeomW3W3_F_sparc64_hpc_ace_double
1377                     (t_nblist                    * gmx_restrict       nlist,
1378                      rvec                        * gmx_restrict          xx,
1379                      rvec                        * gmx_restrict          ff,
1380                      t_forcerec                  * gmx_restrict          fr,
1381                      t_mdatoms                   * gmx_restrict     mdatoms,
1382                      nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
1383                      t_nrnb                      * gmx_restrict        nrnb)
1384 {
1385     /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
1386      * just 0 for non-waters.
1387      * Suffixes A,B refer to j loop unrolling done with double precision SIMD, e.g. for the two different
1388      * jnr indices corresponding to data put in the four positions in the SIMD register.
1389      */
1390     int              i_shift_offset,i_coord_offset,outeriter,inneriter;
1391     int              j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
1392     int              jnrA,jnrB;
1393     int              j_coord_offsetA,j_coord_offsetB;
1394     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
1395     real             rcutoff_scalar;
1396     real             *shiftvec,*fshift,*x,*f;
1397     _fjsp_v2r8       tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
1398     int              vdwioffset0;
1399     _fjsp_v2r8       ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
1400     int              vdwioffset1;
1401     _fjsp_v2r8       ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
1402     int              vdwioffset2;
1403     _fjsp_v2r8       ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
1404     int              vdwjidx0A,vdwjidx0B;
1405     _fjsp_v2r8       jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
1406     int              vdwjidx1A,vdwjidx1B;
1407     _fjsp_v2r8       jx1,jy1,jz1,fjx1,fjy1,fjz1,jq1,isaj1;
1408     int              vdwjidx2A,vdwjidx2B;
1409     _fjsp_v2r8       jx2,jy2,jz2,fjx2,fjy2,fjz2,jq2,isaj2;
1410     _fjsp_v2r8       dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
1411     _fjsp_v2r8       dx01,dy01,dz01,rsq01,rinv01,rinvsq01,r01,qq01,c6_01,c12_01;
1412     _fjsp_v2r8       dx02,dy02,dz02,rsq02,rinv02,rinvsq02,r02,qq02,c6_02,c12_02;
1413     _fjsp_v2r8       dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
1414     _fjsp_v2r8       dx11,dy11,dz11,rsq11,rinv11,rinvsq11,r11,qq11,c6_11,c12_11;
1415     _fjsp_v2r8       dx12,dy12,dz12,rsq12,rinv12,rinvsq12,r12,qq12,c6_12,c12_12;
1416     _fjsp_v2r8       dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
1417     _fjsp_v2r8       dx21,dy21,dz21,rsq21,rinv21,rinvsq21,r21,qq21,c6_21,c12_21;
1418     _fjsp_v2r8       dx22,dy22,dz22,rsq22,rinv22,rinvsq22,r22,qq22,c6_22,c12_22;
1419     _fjsp_v2r8       velec,felec,velecsum,facel,crf,krf,krf2;
1420     real             *charge;
1421     int              nvdwtype;
1422     _fjsp_v2r8       rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
1423     int              *vdwtype;
1424     real             *vdwparam;
1425     _fjsp_v2r8       one_sixth   = gmx_fjsp_set1_v2r8(1.0/6.0);
1426     _fjsp_v2r8       one_twelfth = gmx_fjsp_set1_v2r8(1.0/12.0);
1427     _fjsp_v2r8           c6grid_00;
1428     _fjsp_v2r8           c6grid_01;
1429     _fjsp_v2r8           c6grid_02;
1430     _fjsp_v2r8           c6grid_10;
1431     _fjsp_v2r8           c6grid_11;
1432     _fjsp_v2r8           c6grid_12;
1433     _fjsp_v2r8           c6grid_20;
1434     _fjsp_v2r8           c6grid_21;
1435     _fjsp_v2r8           c6grid_22;
1436     real                 *vdwgridparam;
1437     _fjsp_v2r8           ewclj,ewclj2,ewclj6,ewcljrsq,poly,exponent,f6A,f6B,sh_lj_ewald;
1438     _fjsp_v2r8           one_half = gmx_fjsp_set1_v2r8(0.5);
1439     _fjsp_v2r8           minus_one = gmx_fjsp_set1_v2r8(-1.0);
1440     _fjsp_v2r8       ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
1441     real             *ewtab;
1442     _fjsp_v2r8       itab_tmp;
1443     _fjsp_v2r8       dummy_mask,cutoff_mask;
1444     _fjsp_v2r8       one     = gmx_fjsp_set1_v2r8(1.0);
1445     _fjsp_v2r8       two     = gmx_fjsp_set1_v2r8(2.0);
1446     union { _fjsp_v2r8 simd; long long int i[2]; } vfconv,gbconv,ewconv;
1447
1448     x                = xx[0];
1449     f                = ff[0];
1450
1451     nri              = nlist->nri;
1452     iinr             = nlist->iinr;
1453     jindex           = nlist->jindex;
1454     jjnr             = nlist->jjnr;
1455     shiftidx         = nlist->shift;
1456     gid              = nlist->gid;
1457     shiftvec         = fr->shift_vec[0];
1458     fshift           = fr->fshift[0];
1459     facel            = gmx_fjsp_set1_v2r8(fr->epsfac);
1460     charge           = mdatoms->chargeA;
1461     nvdwtype         = fr->ntype;
1462     vdwparam         = fr->nbfp;
1463     vdwtype          = mdatoms->typeA;
1464     vdwgridparam     = fr->ljpme_c6grid;
1465     sh_lj_ewald      = gmx_fjsp_set1_v2r8(fr->ic->sh_lj_ewald);
1466     ewclj            = gmx_fjsp_set1_v2r8(fr->ewaldcoeff_lj);
1467     ewclj2           = _fjsp_mul_v2r8(minus_one,_fjsp_mul_v2r8(ewclj,ewclj));
1468
1469     sh_ewald         = gmx_fjsp_set1_v2r8(fr->ic->sh_ewald);
1470     ewtab            = fr->ic->tabq_coul_F;
1471     ewtabscale       = gmx_fjsp_set1_v2r8(fr->ic->tabq_scale);
1472     ewtabhalfspace   = gmx_fjsp_set1_v2r8(0.5/fr->ic->tabq_scale);
1473
1474     /* Setup water-specific parameters */
1475     inr              = nlist->iinr[0];
1476     iq0              = _fjsp_mul_v2r8(facel,gmx_fjsp_set1_v2r8(charge[inr+0]));
1477     iq1              = _fjsp_mul_v2r8(facel,gmx_fjsp_set1_v2r8(charge[inr+1]));
1478     iq2              = _fjsp_mul_v2r8(facel,gmx_fjsp_set1_v2r8(charge[inr+2]));
1479     vdwioffset0      = 2*nvdwtype*vdwtype[inr+0];
1480
1481     jq0              = gmx_fjsp_set1_v2r8(charge[inr+0]);
1482     jq1              = gmx_fjsp_set1_v2r8(charge[inr+1]);
1483     jq2              = gmx_fjsp_set1_v2r8(charge[inr+2]);
1484     vdwjidx0A        = 2*vdwtype[inr+0];
1485     qq00             = _fjsp_mul_v2r8(iq0,jq0);
1486     c6_00            = gmx_fjsp_set1_v2r8(vdwparam[vdwioffset0+vdwjidx0A]);
1487     c12_00           = gmx_fjsp_set1_v2r8(vdwparam[vdwioffset0+vdwjidx0A+1]);
1488     c6grid_00        = gmx_fjsp_set1_v2r8(vdwgridparam[vdwioffset0+vdwjidx0A]);
1489     qq01             = _fjsp_mul_v2r8(iq0,jq1);
1490     qq02             = _fjsp_mul_v2r8(iq0,jq2);
1491     qq10             = _fjsp_mul_v2r8(iq1,jq0);
1492     qq11             = _fjsp_mul_v2r8(iq1,jq1);
1493     qq12             = _fjsp_mul_v2r8(iq1,jq2);
1494     qq20             = _fjsp_mul_v2r8(iq2,jq0);
1495     qq21             = _fjsp_mul_v2r8(iq2,jq1);
1496     qq22             = _fjsp_mul_v2r8(iq2,jq2);
1497
1498     /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
1499     rcutoff_scalar   = fr->rcoulomb;
1500     rcutoff          = gmx_fjsp_set1_v2r8(rcutoff_scalar);
1501     rcutoff2         = _fjsp_mul_v2r8(rcutoff,rcutoff);
1502
1503     sh_vdw_invrcut6  = gmx_fjsp_set1_v2r8(fr->ic->sh_invrc6);
1504     rvdw             = gmx_fjsp_set1_v2r8(fr->rvdw);
1505
1506     /* Avoid stupid compiler warnings */
1507     jnrA = jnrB = 0;
1508     j_coord_offsetA = 0;
1509     j_coord_offsetB = 0;
1510
1511     outeriter        = 0;
1512     inneriter        = 0;
1513
1514     /* Start outer loop over neighborlists */
1515     for(iidx=0; iidx<nri; iidx++)
1516     {
1517         /* Load shift vector for this list */
1518         i_shift_offset   = DIM*shiftidx[iidx];
1519
1520         /* Load limits for loop over neighbors */
1521         j_index_start    = jindex[iidx];
1522         j_index_end      = jindex[iidx+1];
1523
1524         /* Get outer coordinate index */
1525         inr              = iinr[iidx];
1526         i_coord_offset   = DIM*inr;
1527
1528         /* Load i particle coords and add shift vector */
1529         gmx_fjsp_load_shift_and_3rvec_broadcast_v2r8(shiftvec+i_shift_offset,x+i_coord_offset,
1530                                                  &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2);
1531
1532         fix0             = _fjsp_setzero_v2r8();
1533         fiy0             = _fjsp_setzero_v2r8();
1534         fiz0             = _fjsp_setzero_v2r8();
1535         fix1             = _fjsp_setzero_v2r8();
1536         fiy1             = _fjsp_setzero_v2r8();
1537         fiz1             = _fjsp_setzero_v2r8();
1538         fix2             = _fjsp_setzero_v2r8();
1539         fiy2             = _fjsp_setzero_v2r8();
1540         fiz2             = _fjsp_setzero_v2r8();
1541
1542         /* Start inner kernel loop */
1543         for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
1544         {
1545
1546             /* Get j neighbor index, and coordinate index */
1547             jnrA             = jjnr[jidx];
1548             jnrB             = jjnr[jidx+1];
1549             j_coord_offsetA  = DIM*jnrA;
1550             j_coord_offsetB  = DIM*jnrB;
1551
1552             /* load j atom coordinates */
1553             gmx_fjsp_load_3rvec_2ptr_swizzle_v2r8(x+j_coord_offsetA,x+j_coord_offsetB,
1554                                               &jx0,&jy0,&jz0,&jx1,&jy1,&jz1,&jx2,&jy2,&jz2);
1555
1556             /* Calculate displacement vector */
1557             dx00             = _fjsp_sub_v2r8(ix0,jx0);
1558             dy00             = _fjsp_sub_v2r8(iy0,jy0);
1559             dz00             = _fjsp_sub_v2r8(iz0,jz0);
1560             dx01             = _fjsp_sub_v2r8(ix0,jx1);
1561             dy01             = _fjsp_sub_v2r8(iy0,jy1);
1562             dz01             = _fjsp_sub_v2r8(iz0,jz1);
1563             dx02             = _fjsp_sub_v2r8(ix0,jx2);
1564             dy02             = _fjsp_sub_v2r8(iy0,jy2);
1565             dz02             = _fjsp_sub_v2r8(iz0,jz2);
1566             dx10             = _fjsp_sub_v2r8(ix1,jx0);
1567             dy10             = _fjsp_sub_v2r8(iy1,jy0);
1568             dz10             = _fjsp_sub_v2r8(iz1,jz0);
1569             dx11             = _fjsp_sub_v2r8(ix1,jx1);
1570             dy11             = _fjsp_sub_v2r8(iy1,jy1);
1571             dz11             = _fjsp_sub_v2r8(iz1,jz1);
1572             dx12             = _fjsp_sub_v2r8(ix1,jx2);
1573             dy12             = _fjsp_sub_v2r8(iy1,jy2);
1574             dz12             = _fjsp_sub_v2r8(iz1,jz2);
1575             dx20             = _fjsp_sub_v2r8(ix2,jx0);
1576             dy20             = _fjsp_sub_v2r8(iy2,jy0);
1577             dz20             = _fjsp_sub_v2r8(iz2,jz0);
1578             dx21             = _fjsp_sub_v2r8(ix2,jx1);
1579             dy21             = _fjsp_sub_v2r8(iy2,jy1);
1580             dz21             = _fjsp_sub_v2r8(iz2,jz1);
1581             dx22             = _fjsp_sub_v2r8(ix2,jx2);
1582             dy22             = _fjsp_sub_v2r8(iy2,jy2);
1583             dz22             = _fjsp_sub_v2r8(iz2,jz2);
1584
1585             /* Calculate squared distance and things based on it */
1586             rsq00            = gmx_fjsp_calc_rsq_v2r8(dx00,dy00,dz00);
1587             rsq01            = gmx_fjsp_calc_rsq_v2r8(dx01,dy01,dz01);
1588             rsq02            = gmx_fjsp_calc_rsq_v2r8(dx02,dy02,dz02);
1589             rsq10            = gmx_fjsp_calc_rsq_v2r8(dx10,dy10,dz10);
1590             rsq11            = gmx_fjsp_calc_rsq_v2r8(dx11,dy11,dz11);
1591             rsq12            = gmx_fjsp_calc_rsq_v2r8(dx12,dy12,dz12);
1592             rsq20            = gmx_fjsp_calc_rsq_v2r8(dx20,dy20,dz20);
1593             rsq21            = gmx_fjsp_calc_rsq_v2r8(dx21,dy21,dz21);
1594             rsq22            = gmx_fjsp_calc_rsq_v2r8(dx22,dy22,dz22);
1595
1596             rinv00           = gmx_fjsp_invsqrt_v2r8(rsq00);
1597             rinv01           = gmx_fjsp_invsqrt_v2r8(rsq01);
1598             rinv02           = gmx_fjsp_invsqrt_v2r8(rsq02);
1599             rinv10           = gmx_fjsp_invsqrt_v2r8(rsq10);
1600             rinv11           = gmx_fjsp_invsqrt_v2r8(rsq11);
1601             rinv12           = gmx_fjsp_invsqrt_v2r8(rsq12);
1602             rinv20           = gmx_fjsp_invsqrt_v2r8(rsq20);
1603             rinv21           = gmx_fjsp_invsqrt_v2r8(rsq21);
1604             rinv22           = gmx_fjsp_invsqrt_v2r8(rsq22);
1605
1606             rinvsq00         = _fjsp_mul_v2r8(rinv00,rinv00);
1607             rinvsq01         = _fjsp_mul_v2r8(rinv01,rinv01);
1608             rinvsq02         = _fjsp_mul_v2r8(rinv02,rinv02);
1609             rinvsq10         = _fjsp_mul_v2r8(rinv10,rinv10);
1610             rinvsq11         = _fjsp_mul_v2r8(rinv11,rinv11);
1611             rinvsq12         = _fjsp_mul_v2r8(rinv12,rinv12);
1612             rinvsq20         = _fjsp_mul_v2r8(rinv20,rinv20);
1613             rinvsq21         = _fjsp_mul_v2r8(rinv21,rinv21);
1614             rinvsq22         = _fjsp_mul_v2r8(rinv22,rinv22);
1615
1616             fjx0             = _fjsp_setzero_v2r8();
1617             fjy0             = _fjsp_setzero_v2r8();
1618             fjz0             = _fjsp_setzero_v2r8();
1619             fjx1             = _fjsp_setzero_v2r8();
1620             fjy1             = _fjsp_setzero_v2r8();
1621             fjz1             = _fjsp_setzero_v2r8();
1622             fjx2             = _fjsp_setzero_v2r8();
1623             fjy2             = _fjsp_setzero_v2r8();
1624             fjz2             = _fjsp_setzero_v2r8();
1625
1626             /**************************
1627              * CALCULATE INTERACTIONS *
1628              **************************/
1629
1630             if (gmx_fjsp_any_lt_v2r8(rsq00,rcutoff2))
1631             {
1632
1633             r00              = _fjsp_mul_v2r8(rsq00,rinv00);
1634
1635             /* EWALD ELECTROSTATICS */
1636
1637             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
1638             ewrt             = _fjsp_mul_v2r8(r00,ewtabscale);
1639             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
1640             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
1641             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
1642
1643             gmx_fjsp_load_2pair_swizzle_v2r8(ewtab+ewconv.i[0],ewtab+ewconv.i[1],
1644                                          &ewtabF,&ewtabFn);
1645             felec            = _fjsp_madd_v2r8(eweps,ewtabFn,_fjsp_nmsub_v2r8(eweps,ewtabF,ewtabF));
1646             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq00,rinv00),_fjsp_sub_v2r8(rinvsq00,felec));
1647
1648             /* Analytical LJ-PME */
1649             rinvsix          = _fjsp_mul_v2r8(_fjsp_mul_v2r8(rinvsq00,rinvsq00),rinvsq00);
1650             ewcljrsq         = _fjsp_mul_v2r8(ewclj2,rsq00);
1651             ewclj6           = _fjsp_mul_v2r8(ewclj2,_fjsp_mul_v2r8(ewclj2,ewclj2));
1652             exponent         = gmx_simd_exp_d(-ewcljrsq);
1653             /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
1654             poly             = _fjsp_mul_v2r8(exponent,_fjsp_madd_v2r8(_fjsp_mul_v2r8(ewcljrsq,ewcljrsq),one_half,_fjsp_sub_v2r8(one,ewcljrsq)));
1655             /* f6A = 6 * C6grid * (1 - poly) */
1656             f6A              = _fjsp_mul_v2r8(c6grid_00,_fjsp_msub_v2r8(one,poly));
1657             /* f6B = C6grid * exponent * beta^6 */
1658             f6B              = _fjsp_mul_v2r8(_fjsp_mul_v2r8(c6grid_00,one_sixth),_fjsp_mul_v2r8(exponent,ewclj6));
1659             /* fvdw = 12*C12/r13 - ((6*C6 - f6A)/r6 + f6B)/r */
1660             fvdw              = _fjsp_mul_v2r8(_fjsp_madd_v2r8(_fjsp_msub_v2r8(c12_00,rinvsix,_fjsp_sub_v2r8(c6_00,f6A)),rinvsix,f6B),rinvsq00);
1661
1662             cutoff_mask      = _fjsp_cmplt_v2r8(rsq00,rcutoff2);
1663
1664             fscal            = _fjsp_add_v2r8(felec,fvdw);
1665
1666             fscal            = _fjsp_and_v2r8(fscal,cutoff_mask);
1667
1668             /* Update vectorial force */
1669             fix0             = _fjsp_madd_v2r8(dx00,fscal,fix0);
1670             fiy0             = _fjsp_madd_v2r8(dy00,fscal,fiy0);
1671             fiz0             = _fjsp_madd_v2r8(dz00,fscal,fiz0);
1672             
1673             fjx0             = _fjsp_madd_v2r8(dx00,fscal,fjx0);
1674             fjy0             = _fjsp_madd_v2r8(dy00,fscal,fjy0);
1675             fjz0             = _fjsp_madd_v2r8(dz00,fscal,fjz0);
1676
1677             }
1678
1679             /**************************
1680              * CALCULATE INTERACTIONS *
1681              **************************/
1682
1683             if (gmx_fjsp_any_lt_v2r8(rsq01,rcutoff2))
1684             {
1685
1686             r01              = _fjsp_mul_v2r8(rsq01,rinv01);
1687
1688             /* EWALD ELECTROSTATICS */
1689
1690             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
1691             ewrt             = _fjsp_mul_v2r8(r01,ewtabscale);
1692             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
1693             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
1694             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
1695
1696             gmx_fjsp_load_2pair_swizzle_v2r8(ewtab+ewconv.i[0],ewtab+ewconv.i[1],
1697                                          &ewtabF,&ewtabFn);
1698             felec            = _fjsp_madd_v2r8(eweps,ewtabFn,_fjsp_nmsub_v2r8(eweps,ewtabF,ewtabF));
1699             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq01,rinv01),_fjsp_sub_v2r8(rinvsq01,felec));
1700
1701             cutoff_mask      = _fjsp_cmplt_v2r8(rsq01,rcutoff2);
1702
1703             fscal            = felec;
1704
1705             fscal            = _fjsp_and_v2r8(fscal,cutoff_mask);
1706
1707             /* Update vectorial force */
1708             fix0             = _fjsp_madd_v2r8(dx01,fscal,fix0);
1709             fiy0             = _fjsp_madd_v2r8(dy01,fscal,fiy0);
1710             fiz0             = _fjsp_madd_v2r8(dz01,fscal,fiz0);
1711             
1712             fjx1             = _fjsp_madd_v2r8(dx01,fscal,fjx1);
1713             fjy1             = _fjsp_madd_v2r8(dy01,fscal,fjy1);
1714             fjz1             = _fjsp_madd_v2r8(dz01,fscal,fjz1);
1715
1716             }
1717
1718             /**************************
1719              * CALCULATE INTERACTIONS *
1720              **************************/
1721
1722             if (gmx_fjsp_any_lt_v2r8(rsq02,rcutoff2))
1723             {
1724
1725             r02              = _fjsp_mul_v2r8(rsq02,rinv02);
1726
1727             /* EWALD ELECTROSTATICS */
1728
1729             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
1730             ewrt             = _fjsp_mul_v2r8(r02,ewtabscale);
1731             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
1732             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
1733             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
1734
1735             gmx_fjsp_load_2pair_swizzle_v2r8(ewtab+ewconv.i[0],ewtab+ewconv.i[1],
1736                                          &ewtabF,&ewtabFn);
1737             felec            = _fjsp_madd_v2r8(eweps,ewtabFn,_fjsp_nmsub_v2r8(eweps,ewtabF,ewtabF));
1738             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq02,rinv02),_fjsp_sub_v2r8(rinvsq02,felec));
1739
1740             cutoff_mask      = _fjsp_cmplt_v2r8(rsq02,rcutoff2);
1741
1742             fscal            = felec;
1743
1744             fscal            = _fjsp_and_v2r8(fscal,cutoff_mask);
1745
1746             /* Update vectorial force */
1747             fix0             = _fjsp_madd_v2r8(dx02,fscal,fix0);
1748             fiy0             = _fjsp_madd_v2r8(dy02,fscal,fiy0);
1749             fiz0             = _fjsp_madd_v2r8(dz02,fscal,fiz0);
1750             
1751             fjx2             = _fjsp_madd_v2r8(dx02,fscal,fjx2);
1752             fjy2             = _fjsp_madd_v2r8(dy02,fscal,fjy2);
1753             fjz2             = _fjsp_madd_v2r8(dz02,fscal,fjz2);
1754
1755             }
1756
1757             /**************************
1758              * CALCULATE INTERACTIONS *
1759              **************************/
1760
1761             if (gmx_fjsp_any_lt_v2r8(rsq10,rcutoff2))
1762             {
1763
1764             r10              = _fjsp_mul_v2r8(rsq10,rinv10);
1765
1766             /* EWALD ELECTROSTATICS */
1767
1768             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
1769             ewrt             = _fjsp_mul_v2r8(r10,ewtabscale);
1770             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
1771             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
1772             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
1773
1774             gmx_fjsp_load_2pair_swizzle_v2r8(ewtab+ewconv.i[0],ewtab+ewconv.i[1],
1775                                          &ewtabF,&ewtabFn);
1776             felec            = _fjsp_madd_v2r8(eweps,ewtabFn,_fjsp_nmsub_v2r8(eweps,ewtabF,ewtabF));
1777             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq10,rinv10),_fjsp_sub_v2r8(rinvsq10,felec));
1778
1779             cutoff_mask      = _fjsp_cmplt_v2r8(rsq10,rcutoff2);
1780
1781             fscal            = felec;
1782
1783             fscal            = _fjsp_and_v2r8(fscal,cutoff_mask);
1784
1785             /* Update vectorial force */
1786             fix1             = _fjsp_madd_v2r8(dx10,fscal,fix1);
1787             fiy1             = _fjsp_madd_v2r8(dy10,fscal,fiy1);
1788             fiz1             = _fjsp_madd_v2r8(dz10,fscal,fiz1);
1789             
1790             fjx0             = _fjsp_madd_v2r8(dx10,fscal,fjx0);
1791             fjy0             = _fjsp_madd_v2r8(dy10,fscal,fjy0);
1792             fjz0             = _fjsp_madd_v2r8(dz10,fscal,fjz0);
1793
1794             }
1795
1796             /**************************
1797              * CALCULATE INTERACTIONS *
1798              **************************/
1799
1800             if (gmx_fjsp_any_lt_v2r8(rsq11,rcutoff2))
1801             {
1802
1803             r11              = _fjsp_mul_v2r8(rsq11,rinv11);
1804
1805             /* EWALD ELECTROSTATICS */
1806
1807             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
1808             ewrt             = _fjsp_mul_v2r8(r11,ewtabscale);
1809             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
1810             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
1811             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
1812
1813             gmx_fjsp_load_2pair_swizzle_v2r8(ewtab+ewconv.i[0],ewtab+ewconv.i[1],
1814                                          &ewtabF,&ewtabFn);
1815             felec            = _fjsp_madd_v2r8(eweps,ewtabFn,_fjsp_nmsub_v2r8(eweps,ewtabF,ewtabF));
1816             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq11,rinv11),_fjsp_sub_v2r8(rinvsq11,felec));
1817
1818             cutoff_mask      = _fjsp_cmplt_v2r8(rsq11,rcutoff2);
1819
1820             fscal            = felec;
1821
1822             fscal            = _fjsp_and_v2r8(fscal,cutoff_mask);
1823
1824             /* Update vectorial force */
1825             fix1             = _fjsp_madd_v2r8(dx11,fscal,fix1);
1826             fiy1             = _fjsp_madd_v2r8(dy11,fscal,fiy1);
1827             fiz1             = _fjsp_madd_v2r8(dz11,fscal,fiz1);
1828             
1829             fjx1             = _fjsp_madd_v2r8(dx11,fscal,fjx1);
1830             fjy1             = _fjsp_madd_v2r8(dy11,fscal,fjy1);
1831             fjz1             = _fjsp_madd_v2r8(dz11,fscal,fjz1);
1832
1833             }
1834
1835             /**************************
1836              * CALCULATE INTERACTIONS *
1837              **************************/
1838
1839             if (gmx_fjsp_any_lt_v2r8(rsq12,rcutoff2))
1840             {
1841
1842             r12              = _fjsp_mul_v2r8(rsq12,rinv12);
1843
1844             /* EWALD ELECTROSTATICS */
1845
1846             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
1847             ewrt             = _fjsp_mul_v2r8(r12,ewtabscale);
1848             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
1849             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
1850             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
1851
1852             gmx_fjsp_load_2pair_swizzle_v2r8(ewtab+ewconv.i[0],ewtab+ewconv.i[1],
1853                                          &ewtabF,&ewtabFn);
1854             felec            = _fjsp_madd_v2r8(eweps,ewtabFn,_fjsp_nmsub_v2r8(eweps,ewtabF,ewtabF));
1855             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq12,rinv12),_fjsp_sub_v2r8(rinvsq12,felec));
1856
1857             cutoff_mask      = _fjsp_cmplt_v2r8(rsq12,rcutoff2);
1858
1859             fscal            = felec;
1860
1861             fscal            = _fjsp_and_v2r8(fscal,cutoff_mask);
1862
1863             /* Update vectorial force */
1864             fix1             = _fjsp_madd_v2r8(dx12,fscal,fix1);
1865             fiy1             = _fjsp_madd_v2r8(dy12,fscal,fiy1);
1866             fiz1             = _fjsp_madd_v2r8(dz12,fscal,fiz1);
1867             
1868             fjx2             = _fjsp_madd_v2r8(dx12,fscal,fjx2);
1869             fjy2             = _fjsp_madd_v2r8(dy12,fscal,fjy2);
1870             fjz2             = _fjsp_madd_v2r8(dz12,fscal,fjz2);
1871
1872             }
1873
1874             /**************************
1875              * CALCULATE INTERACTIONS *
1876              **************************/
1877
1878             if (gmx_fjsp_any_lt_v2r8(rsq20,rcutoff2))
1879             {
1880
1881             r20              = _fjsp_mul_v2r8(rsq20,rinv20);
1882
1883             /* EWALD ELECTROSTATICS */
1884
1885             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
1886             ewrt             = _fjsp_mul_v2r8(r20,ewtabscale);
1887             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
1888             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
1889             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
1890
1891             gmx_fjsp_load_2pair_swizzle_v2r8(ewtab+ewconv.i[0],ewtab+ewconv.i[1],
1892                                          &ewtabF,&ewtabFn);
1893             felec            = _fjsp_madd_v2r8(eweps,ewtabFn,_fjsp_nmsub_v2r8(eweps,ewtabF,ewtabF));
1894             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq20,rinv20),_fjsp_sub_v2r8(rinvsq20,felec));
1895
1896             cutoff_mask      = _fjsp_cmplt_v2r8(rsq20,rcutoff2);
1897
1898             fscal            = felec;
1899
1900             fscal            = _fjsp_and_v2r8(fscal,cutoff_mask);
1901
1902             /* Update vectorial force */
1903             fix2             = _fjsp_madd_v2r8(dx20,fscal,fix2);
1904             fiy2             = _fjsp_madd_v2r8(dy20,fscal,fiy2);
1905             fiz2             = _fjsp_madd_v2r8(dz20,fscal,fiz2);
1906             
1907             fjx0             = _fjsp_madd_v2r8(dx20,fscal,fjx0);
1908             fjy0             = _fjsp_madd_v2r8(dy20,fscal,fjy0);
1909             fjz0             = _fjsp_madd_v2r8(dz20,fscal,fjz0);
1910
1911             }
1912
1913             /**************************
1914              * CALCULATE INTERACTIONS *
1915              **************************/
1916
1917             if (gmx_fjsp_any_lt_v2r8(rsq21,rcutoff2))
1918             {
1919
1920             r21              = _fjsp_mul_v2r8(rsq21,rinv21);
1921
1922             /* EWALD ELECTROSTATICS */
1923
1924             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
1925             ewrt             = _fjsp_mul_v2r8(r21,ewtabscale);
1926             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
1927             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
1928             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
1929
1930             gmx_fjsp_load_2pair_swizzle_v2r8(ewtab+ewconv.i[0],ewtab+ewconv.i[1],
1931                                          &ewtabF,&ewtabFn);
1932             felec            = _fjsp_madd_v2r8(eweps,ewtabFn,_fjsp_nmsub_v2r8(eweps,ewtabF,ewtabF));
1933             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq21,rinv21),_fjsp_sub_v2r8(rinvsq21,felec));
1934
1935             cutoff_mask      = _fjsp_cmplt_v2r8(rsq21,rcutoff2);
1936
1937             fscal            = felec;
1938
1939             fscal            = _fjsp_and_v2r8(fscal,cutoff_mask);
1940
1941             /* Update vectorial force */
1942             fix2             = _fjsp_madd_v2r8(dx21,fscal,fix2);
1943             fiy2             = _fjsp_madd_v2r8(dy21,fscal,fiy2);
1944             fiz2             = _fjsp_madd_v2r8(dz21,fscal,fiz2);
1945             
1946             fjx1             = _fjsp_madd_v2r8(dx21,fscal,fjx1);
1947             fjy1             = _fjsp_madd_v2r8(dy21,fscal,fjy1);
1948             fjz1             = _fjsp_madd_v2r8(dz21,fscal,fjz1);
1949
1950             }
1951
1952             /**************************
1953              * CALCULATE INTERACTIONS *
1954              **************************/
1955
1956             if (gmx_fjsp_any_lt_v2r8(rsq22,rcutoff2))
1957             {
1958
1959             r22              = _fjsp_mul_v2r8(rsq22,rinv22);
1960
1961             /* EWALD ELECTROSTATICS */
1962
1963             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
1964             ewrt             = _fjsp_mul_v2r8(r22,ewtabscale);
1965             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
1966             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
1967             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
1968
1969             gmx_fjsp_load_2pair_swizzle_v2r8(ewtab+ewconv.i[0],ewtab+ewconv.i[1],
1970                                          &ewtabF,&ewtabFn);
1971             felec            = _fjsp_madd_v2r8(eweps,ewtabFn,_fjsp_nmsub_v2r8(eweps,ewtabF,ewtabF));
1972             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq22,rinv22),_fjsp_sub_v2r8(rinvsq22,felec));
1973
1974             cutoff_mask      = _fjsp_cmplt_v2r8(rsq22,rcutoff2);
1975
1976             fscal            = felec;
1977
1978             fscal            = _fjsp_and_v2r8(fscal,cutoff_mask);
1979
1980             /* Update vectorial force */
1981             fix2             = _fjsp_madd_v2r8(dx22,fscal,fix2);
1982             fiy2             = _fjsp_madd_v2r8(dy22,fscal,fiy2);
1983             fiz2             = _fjsp_madd_v2r8(dz22,fscal,fiz2);
1984             
1985             fjx2             = _fjsp_madd_v2r8(dx22,fscal,fjx2);
1986             fjy2             = _fjsp_madd_v2r8(dy22,fscal,fjy2);
1987             fjz2             = _fjsp_madd_v2r8(dz22,fscal,fjz2);
1988
1989             }
1990
1991             gmx_fjsp_decrement_3rvec_2ptr_swizzle_v2r8(f+j_coord_offsetA,f+j_coord_offsetB,fjx0,fjy0,fjz0,fjx1,fjy1,fjz1,fjx2,fjy2,fjz2);
1992
1993             /* Inner loop uses 400 flops */
1994         }
1995
1996         if(jidx<j_index_end)
1997         {
1998
1999             jnrA             = jjnr[jidx];
2000             j_coord_offsetA  = DIM*jnrA;
2001
2002             /* load j atom coordinates */
2003             gmx_fjsp_load_3rvec_1ptr_swizzle_v2r8(x+j_coord_offsetA,
2004                                               &jx0,&jy0,&jz0,&jx1,&jy1,&jz1,&jx2,&jy2,&jz2);
2005
2006             /* Calculate displacement vector */
2007             dx00             = _fjsp_sub_v2r8(ix0,jx0);
2008             dy00             = _fjsp_sub_v2r8(iy0,jy0);
2009             dz00             = _fjsp_sub_v2r8(iz0,jz0);
2010             dx01             = _fjsp_sub_v2r8(ix0,jx1);
2011             dy01             = _fjsp_sub_v2r8(iy0,jy1);
2012             dz01             = _fjsp_sub_v2r8(iz0,jz1);
2013             dx02             = _fjsp_sub_v2r8(ix0,jx2);
2014             dy02             = _fjsp_sub_v2r8(iy0,jy2);
2015             dz02             = _fjsp_sub_v2r8(iz0,jz2);
2016             dx10             = _fjsp_sub_v2r8(ix1,jx0);
2017             dy10             = _fjsp_sub_v2r8(iy1,jy0);
2018             dz10             = _fjsp_sub_v2r8(iz1,jz0);
2019             dx11             = _fjsp_sub_v2r8(ix1,jx1);
2020             dy11             = _fjsp_sub_v2r8(iy1,jy1);
2021             dz11             = _fjsp_sub_v2r8(iz1,jz1);
2022             dx12             = _fjsp_sub_v2r8(ix1,jx2);
2023             dy12             = _fjsp_sub_v2r8(iy1,jy2);
2024             dz12             = _fjsp_sub_v2r8(iz1,jz2);
2025             dx20             = _fjsp_sub_v2r8(ix2,jx0);
2026             dy20             = _fjsp_sub_v2r8(iy2,jy0);
2027             dz20             = _fjsp_sub_v2r8(iz2,jz0);
2028             dx21             = _fjsp_sub_v2r8(ix2,jx1);
2029             dy21             = _fjsp_sub_v2r8(iy2,jy1);
2030             dz21             = _fjsp_sub_v2r8(iz2,jz1);
2031             dx22             = _fjsp_sub_v2r8(ix2,jx2);
2032             dy22             = _fjsp_sub_v2r8(iy2,jy2);
2033             dz22             = _fjsp_sub_v2r8(iz2,jz2);
2034
2035             /* Calculate squared distance and things based on it */
2036             rsq00            = gmx_fjsp_calc_rsq_v2r8(dx00,dy00,dz00);
2037             rsq01            = gmx_fjsp_calc_rsq_v2r8(dx01,dy01,dz01);
2038             rsq02            = gmx_fjsp_calc_rsq_v2r8(dx02,dy02,dz02);
2039             rsq10            = gmx_fjsp_calc_rsq_v2r8(dx10,dy10,dz10);
2040             rsq11            = gmx_fjsp_calc_rsq_v2r8(dx11,dy11,dz11);
2041             rsq12            = gmx_fjsp_calc_rsq_v2r8(dx12,dy12,dz12);
2042             rsq20            = gmx_fjsp_calc_rsq_v2r8(dx20,dy20,dz20);
2043             rsq21            = gmx_fjsp_calc_rsq_v2r8(dx21,dy21,dz21);
2044             rsq22            = gmx_fjsp_calc_rsq_v2r8(dx22,dy22,dz22);
2045
2046             rinv00           = gmx_fjsp_invsqrt_v2r8(rsq00);
2047             rinv01           = gmx_fjsp_invsqrt_v2r8(rsq01);
2048             rinv02           = gmx_fjsp_invsqrt_v2r8(rsq02);
2049             rinv10           = gmx_fjsp_invsqrt_v2r8(rsq10);
2050             rinv11           = gmx_fjsp_invsqrt_v2r8(rsq11);
2051             rinv12           = gmx_fjsp_invsqrt_v2r8(rsq12);
2052             rinv20           = gmx_fjsp_invsqrt_v2r8(rsq20);
2053             rinv21           = gmx_fjsp_invsqrt_v2r8(rsq21);
2054             rinv22           = gmx_fjsp_invsqrt_v2r8(rsq22);
2055
2056             rinvsq00         = _fjsp_mul_v2r8(rinv00,rinv00);
2057             rinvsq01         = _fjsp_mul_v2r8(rinv01,rinv01);
2058             rinvsq02         = _fjsp_mul_v2r8(rinv02,rinv02);
2059             rinvsq10         = _fjsp_mul_v2r8(rinv10,rinv10);
2060             rinvsq11         = _fjsp_mul_v2r8(rinv11,rinv11);
2061             rinvsq12         = _fjsp_mul_v2r8(rinv12,rinv12);
2062             rinvsq20         = _fjsp_mul_v2r8(rinv20,rinv20);
2063             rinvsq21         = _fjsp_mul_v2r8(rinv21,rinv21);
2064             rinvsq22         = _fjsp_mul_v2r8(rinv22,rinv22);
2065
2066             fjx0             = _fjsp_setzero_v2r8();
2067             fjy0             = _fjsp_setzero_v2r8();
2068             fjz0             = _fjsp_setzero_v2r8();
2069             fjx1             = _fjsp_setzero_v2r8();
2070             fjy1             = _fjsp_setzero_v2r8();
2071             fjz1             = _fjsp_setzero_v2r8();
2072             fjx2             = _fjsp_setzero_v2r8();
2073             fjy2             = _fjsp_setzero_v2r8();
2074             fjz2             = _fjsp_setzero_v2r8();
2075
2076             /**************************
2077              * CALCULATE INTERACTIONS *
2078              **************************/
2079
2080             if (gmx_fjsp_any_lt_v2r8(rsq00,rcutoff2))
2081             {
2082
2083             r00              = _fjsp_mul_v2r8(rsq00,rinv00);
2084
2085             /* EWALD ELECTROSTATICS */
2086
2087             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
2088             ewrt             = _fjsp_mul_v2r8(r00,ewtabscale);
2089             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
2090             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
2091             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
2092
2093             gmx_fjsp_load_1pair_swizzle_v2r8(ewtab+ewconv.i[0],&ewtabF,&ewtabFn);
2094             felec            = _fjsp_madd_v2r8(eweps,ewtabFn,_fjsp_nmsub_v2r8(eweps,ewtabF,ewtabF));
2095             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq00,rinv00),_fjsp_sub_v2r8(rinvsq00,felec));
2096
2097             /* Analytical LJ-PME */
2098             rinvsix          = _fjsp_mul_v2r8(_fjsp_mul_v2r8(rinvsq00,rinvsq00),rinvsq00);
2099             ewcljrsq         = _fjsp_mul_v2r8(ewclj2,rsq00);
2100             ewclj6           = _fjsp_mul_v2r8(ewclj2,_fjsp_mul_v2r8(ewclj2,ewclj2));
2101             exponent         = gmx_simd_exp_d(-ewcljrsq);
2102             /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
2103             poly             = _fjsp_mul_v2r8(exponent,_fjsp_madd_v2r8(_fjsp_mul_v2r8(ewcljrsq,ewcljrsq),one_half,_fjsp_sub_v2r8(one,ewcljrsq)));
2104             /* f6A = 6 * C6grid * (1 - poly) */
2105             f6A              = _fjsp_mul_v2r8(c6grid_00,_fjsp_msub_v2r8(one,poly));
2106             /* f6B = C6grid * exponent * beta^6 */
2107             f6B              = _fjsp_mul_v2r8(_fjsp_mul_v2r8(c6grid_00,one_sixth),_fjsp_mul_v2r8(exponent,ewclj6));
2108             /* fvdw = 12*C12/r13 - ((6*C6 - f6A)/r6 + f6B)/r */
2109             fvdw              = _fjsp_mul_v2r8(_fjsp_madd_v2r8(_fjsp_msub_v2r8(c12_00,rinvsix,_fjsp_sub_v2r8(c6_00,f6A)),rinvsix,f6B),rinvsq00);
2110
2111             cutoff_mask      = _fjsp_cmplt_v2r8(rsq00,rcutoff2);
2112
2113             fscal            = _fjsp_add_v2r8(felec,fvdw);
2114
2115             fscal            = _fjsp_and_v2r8(fscal,cutoff_mask);
2116
2117             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
2118
2119             /* Update vectorial force */
2120             fix0             = _fjsp_madd_v2r8(dx00,fscal,fix0);
2121             fiy0             = _fjsp_madd_v2r8(dy00,fscal,fiy0);
2122             fiz0             = _fjsp_madd_v2r8(dz00,fscal,fiz0);
2123             
2124             fjx0             = _fjsp_madd_v2r8(dx00,fscal,fjx0);
2125             fjy0             = _fjsp_madd_v2r8(dy00,fscal,fjy0);
2126             fjz0             = _fjsp_madd_v2r8(dz00,fscal,fjz0);
2127
2128             }
2129
2130             /**************************
2131              * CALCULATE INTERACTIONS *
2132              **************************/
2133
2134             if (gmx_fjsp_any_lt_v2r8(rsq01,rcutoff2))
2135             {
2136
2137             r01              = _fjsp_mul_v2r8(rsq01,rinv01);
2138
2139             /* EWALD ELECTROSTATICS */
2140
2141             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
2142             ewrt             = _fjsp_mul_v2r8(r01,ewtabscale);
2143             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
2144             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
2145             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
2146
2147             gmx_fjsp_load_1pair_swizzle_v2r8(ewtab+ewconv.i[0],&ewtabF,&ewtabFn);
2148             felec            = _fjsp_madd_v2r8(eweps,ewtabFn,_fjsp_nmsub_v2r8(eweps,ewtabF,ewtabF));
2149             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq01,rinv01),_fjsp_sub_v2r8(rinvsq01,felec));
2150
2151             cutoff_mask      = _fjsp_cmplt_v2r8(rsq01,rcutoff2);
2152
2153             fscal            = felec;
2154
2155             fscal            = _fjsp_and_v2r8(fscal,cutoff_mask);
2156
2157             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
2158
2159             /* Update vectorial force */
2160             fix0             = _fjsp_madd_v2r8(dx01,fscal,fix0);
2161             fiy0             = _fjsp_madd_v2r8(dy01,fscal,fiy0);
2162             fiz0             = _fjsp_madd_v2r8(dz01,fscal,fiz0);
2163             
2164             fjx1             = _fjsp_madd_v2r8(dx01,fscal,fjx1);
2165             fjy1             = _fjsp_madd_v2r8(dy01,fscal,fjy1);
2166             fjz1             = _fjsp_madd_v2r8(dz01,fscal,fjz1);
2167
2168             }
2169
2170             /**************************
2171              * CALCULATE INTERACTIONS *
2172              **************************/
2173
2174             if (gmx_fjsp_any_lt_v2r8(rsq02,rcutoff2))
2175             {
2176
2177             r02              = _fjsp_mul_v2r8(rsq02,rinv02);
2178
2179             /* EWALD ELECTROSTATICS */
2180
2181             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
2182             ewrt             = _fjsp_mul_v2r8(r02,ewtabscale);
2183             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
2184             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
2185             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
2186
2187             gmx_fjsp_load_1pair_swizzle_v2r8(ewtab+ewconv.i[0],&ewtabF,&ewtabFn);
2188             felec            = _fjsp_madd_v2r8(eweps,ewtabFn,_fjsp_nmsub_v2r8(eweps,ewtabF,ewtabF));
2189             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq02,rinv02),_fjsp_sub_v2r8(rinvsq02,felec));
2190
2191             cutoff_mask      = _fjsp_cmplt_v2r8(rsq02,rcutoff2);
2192
2193             fscal            = felec;
2194
2195             fscal            = _fjsp_and_v2r8(fscal,cutoff_mask);
2196
2197             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
2198
2199             /* Update vectorial force */
2200             fix0             = _fjsp_madd_v2r8(dx02,fscal,fix0);
2201             fiy0             = _fjsp_madd_v2r8(dy02,fscal,fiy0);
2202             fiz0             = _fjsp_madd_v2r8(dz02,fscal,fiz0);
2203             
2204             fjx2             = _fjsp_madd_v2r8(dx02,fscal,fjx2);
2205             fjy2             = _fjsp_madd_v2r8(dy02,fscal,fjy2);
2206             fjz2             = _fjsp_madd_v2r8(dz02,fscal,fjz2);
2207
2208             }
2209
2210             /**************************
2211              * CALCULATE INTERACTIONS *
2212              **************************/
2213
2214             if (gmx_fjsp_any_lt_v2r8(rsq10,rcutoff2))
2215             {
2216
2217             r10              = _fjsp_mul_v2r8(rsq10,rinv10);
2218
2219             /* EWALD ELECTROSTATICS */
2220
2221             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
2222             ewrt             = _fjsp_mul_v2r8(r10,ewtabscale);
2223             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
2224             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
2225             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
2226
2227             gmx_fjsp_load_1pair_swizzle_v2r8(ewtab+ewconv.i[0],&ewtabF,&ewtabFn);
2228             felec            = _fjsp_madd_v2r8(eweps,ewtabFn,_fjsp_nmsub_v2r8(eweps,ewtabF,ewtabF));
2229             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq10,rinv10),_fjsp_sub_v2r8(rinvsq10,felec));
2230
2231             cutoff_mask      = _fjsp_cmplt_v2r8(rsq10,rcutoff2);
2232
2233             fscal            = felec;
2234
2235             fscal            = _fjsp_and_v2r8(fscal,cutoff_mask);
2236
2237             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
2238
2239             /* Update vectorial force */
2240             fix1             = _fjsp_madd_v2r8(dx10,fscal,fix1);
2241             fiy1             = _fjsp_madd_v2r8(dy10,fscal,fiy1);
2242             fiz1             = _fjsp_madd_v2r8(dz10,fscal,fiz1);
2243             
2244             fjx0             = _fjsp_madd_v2r8(dx10,fscal,fjx0);
2245             fjy0             = _fjsp_madd_v2r8(dy10,fscal,fjy0);
2246             fjz0             = _fjsp_madd_v2r8(dz10,fscal,fjz0);
2247
2248             }
2249
2250             /**************************
2251              * CALCULATE INTERACTIONS *
2252              **************************/
2253
2254             if (gmx_fjsp_any_lt_v2r8(rsq11,rcutoff2))
2255             {
2256
2257             r11              = _fjsp_mul_v2r8(rsq11,rinv11);
2258
2259             /* EWALD ELECTROSTATICS */
2260
2261             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
2262             ewrt             = _fjsp_mul_v2r8(r11,ewtabscale);
2263             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
2264             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
2265             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
2266
2267             gmx_fjsp_load_1pair_swizzle_v2r8(ewtab+ewconv.i[0],&ewtabF,&ewtabFn);
2268             felec            = _fjsp_madd_v2r8(eweps,ewtabFn,_fjsp_nmsub_v2r8(eweps,ewtabF,ewtabF));
2269             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq11,rinv11),_fjsp_sub_v2r8(rinvsq11,felec));
2270
2271             cutoff_mask      = _fjsp_cmplt_v2r8(rsq11,rcutoff2);
2272
2273             fscal            = felec;
2274
2275             fscal            = _fjsp_and_v2r8(fscal,cutoff_mask);
2276
2277             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
2278
2279             /* Update vectorial force */
2280             fix1             = _fjsp_madd_v2r8(dx11,fscal,fix1);
2281             fiy1             = _fjsp_madd_v2r8(dy11,fscal,fiy1);
2282             fiz1             = _fjsp_madd_v2r8(dz11,fscal,fiz1);
2283             
2284             fjx1             = _fjsp_madd_v2r8(dx11,fscal,fjx1);
2285             fjy1             = _fjsp_madd_v2r8(dy11,fscal,fjy1);
2286             fjz1             = _fjsp_madd_v2r8(dz11,fscal,fjz1);
2287
2288             }
2289
2290             /**************************
2291              * CALCULATE INTERACTIONS *
2292              **************************/
2293
2294             if (gmx_fjsp_any_lt_v2r8(rsq12,rcutoff2))
2295             {
2296
2297             r12              = _fjsp_mul_v2r8(rsq12,rinv12);
2298
2299             /* EWALD ELECTROSTATICS */
2300
2301             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
2302             ewrt             = _fjsp_mul_v2r8(r12,ewtabscale);
2303             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
2304             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
2305             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
2306
2307             gmx_fjsp_load_1pair_swizzle_v2r8(ewtab+ewconv.i[0],&ewtabF,&ewtabFn);
2308             felec            = _fjsp_madd_v2r8(eweps,ewtabFn,_fjsp_nmsub_v2r8(eweps,ewtabF,ewtabF));
2309             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq12,rinv12),_fjsp_sub_v2r8(rinvsq12,felec));
2310
2311             cutoff_mask      = _fjsp_cmplt_v2r8(rsq12,rcutoff2);
2312
2313             fscal            = felec;
2314
2315             fscal            = _fjsp_and_v2r8(fscal,cutoff_mask);
2316
2317             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
2318
2319             /* Update vectorial force */
2320             fix1             = _fjsp_madd_v2r8(dx12,fscal,fix1);
2321             fiy1             = _fjsp_madd_v2r8(dy12,fscal,fiy1);
2322             fiz1             = _fjsp_madd_v2r8(dz12,fscal,fiz1);
2323             
2324             fjx2             = _fjsp_madd_v2r8(dx12,fscal,fjx2);
2325             fjy2             = _fjsp_madd_v2r8(dy12,fscal,fjy2);
2326             fjz2             = _fjsp_madd_v2r8(dz12,fscal,fjz2);
2327
2328             }
2329
2330             /**************************
2331              * CALCULATE INTERACTIONS *
2332              **************************/
2333
2334             if (gmx_fjsp_any_lt_v2r8(rsq20,rcutoff2))
2335             {
2336
2337             r20              = _fjsp_mul_v2r8(rsq20,rinv20);
2338
2339             /* EWALD ELECTROSTATICS */
2340
2341             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
2342             ewrt             = _fjsp_mul_v2r8(r20,ewtabscale);
2343             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
2344             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
2345             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
2346
2347             gmx_fjsp_load_1pair_swizzle_v2r8(ewtab+ewconv.i[0],&ewtabF,&ewtabFn);
2348             felec            = _fjsp_madd_v2r8(eweps,ewtabFn,_fjsp_nmsub_v2r8(eweps,ewtabF,ewtabF));
2349             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq20,rinv20),_fjsp_sub_v2r8(rinvsq20,felec));
2350
2351             cutoff_mask      = _fjsp_cmplt_v2r8(rsq20,rcutoff2);
2352
2353             fscal            = felec;
2354
2355             fscal            = _fjsp_and_v2r8(fscal,cutoff_mask);
2356
2357             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
2358
2359             /* Update vectorial force */
2360             fix2             = _fjsp_madd_v2r8(dx20,fscal,fix2);
2361             fiy2             = _fjsp_madd_v2r8(dy20,fscal,fiy2);
2362             fiz2             = _fjsp_madd_v2r8(dz20,fscal,fiz2);
2363             
2364             fjx0             = _fjsp_madd_v2r8(dx20,fscal,fjx0);
2365             fjy0             = _fjsp_madd_v2r8(dy20,fscal,fjy0);
2366             fjz0             = _fjsp_madd_v2r8(dz20,fscal,fjz0);
2367
2368             }
2369
2370             /**************************
2371              * CALCULATE INTERACTIONS *
2372              **************************/
2373
2374             if (gmx_fjsp_any_lt_v2r8(rsq21,rcutoff2))
2375             {
2376
2377             r21              = _fjsp_mul_v2r8(rsq21,rinv21);
2378
2379             /* EWALD ELECTROSTATICS */
2380
2381             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
2382             ewrt             = _fjsp_mul_v2r8(r21,ewtabscale);
2383             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
2384             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
2385             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
2386
2387             gmx_fjsp_load_1pair_swizzle_v2r8(ewtab+ewconv.i[0],&ewtabF,&ewtabFn);
2388             felec            = _fjsp_madd_v2r8(eweps,ewtabFn,_fjsp_nmsub_v2r8(eweps,ewtabF,ewtabF));
2389             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq21,rinv21),_fjsp_sub_v2r8(rinvsq21,felec));
2390
2391             cutoff_mask      = _fjsp_cmplt_v2r8(rsq21,rcutoff2);
2392
2393             fscal            = felec;
2394
2395             fscal            = _fjsp_and_v2r8(fscal,cutoff_mask);
2396
2397             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
2398
2399             /* Update vectorial force */
2400             fix2             = _fjsp_madd_v2r8(dx21,fscal,fix2);
2401             fiy2             = _fjsp_madd_v2r8(dy21,fscal,fiy2);
2402             fiz2             = _fjsp_madd_v2r8(dz21,fscal,fiz2);
2403             
2404             fjx1             = _fjsp_madd_v2r8(dx21,fscal,fjx1);
2405             fjy1             = _fjsp_madd_v2r8(dy21,fscal,fjy1);
2406             fjz1             = _fjsp_madd_v2r8(dz21,fscal,fjz1);
2407
2408             }
2409
2410             /**************************
2411              * CALCULATE INTERACTIONS *
2412              **************************/
2413
2414             if (gmx_fjsp_any_lt_v2r8(rsq22,rcutoff2))
2415             {
2416
2417             r22              = _fjsp_mul_v2r8(rsq22,rinv22);
2418
2419             /* EWALD ELECTROSTATICS */
2420
2421             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
2422             ewrt             = _fjsp_mul_v2r8(r22,ewtabscale);
2423             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
2424             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
2425             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
2426
2427             gmx_fjsp_load_1pair_swizzle_v2r8(ewtab+ewconv.i[0],&ewtabF,&ewtabFn);
2428             felec            = _fjsp_madd_v2r8(eweps,ewtabFn,_fjsp_nmsub_v2r8(eweps,ewtabF,ewtabF));
2429             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq22,rinv22),_fjsp_sub_v2r8(rinvsq22,felec));
2430
2431             cutoff_mask      = _fjsp_cmplt_v2r8(rsq22,rcutoff2);
2432
2433             fscal            = felec;
2434
2435             fscal            = _fjsp_and_v2r8(fscal,cutoff_mask);
2436
2437             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
2438
2439             /* Update vectorial force */
2440             fix2             = _fjsp_madd_v2r8(dx22,fscal,fix2);
2441             fiy2             = _fjsp_madd_v2r8(dy22,fscal,fiy2);
2442             fiz2             = _fjsp_madd_v2r8(dz22,fscal,fiz2);
2443             
2444             fjx2             = _fjsp_madd_v2r8(dx22,fscal,fjx2);
2445             fjy2             = _fjsp_madd_v2r8(dy22,fscal,fjy2);
2446             fjz2             = _fjsp_madd_v2r8(dz22,fscal,fjz2);
2447
2448             }
2449
2450             gmx_fjsp_decrement_3rvec_1ptr_swizzle_v2r8(f+j_coord_offsetA,fjx0,fjy0,fjz0,fjx1,fjy1,fjz1,fjx2,fjy2,fjz2);
2451
2452             /* Inner loop uses 400 flops */
2453         }
2454
2455         /* End of innermost loop */
2456
2457         gmx_fjsp_update_iforce_3atom_swizzle_v2r8(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,
2458                                               f+i_coord_offset,fshift+i_shift_offset);
2459
2460         /* Increment number of inner iterations */
2461         inneriter                  += j_index_end - j_index_start;
2462
2463         /* Outer loop uses 18 flops */
2464     }
2465
2466     /* Increment number of outer iterations */
2467     outeriter        += nri;
2468
2469     /* Update outer/inner flops */
2470
2471     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W3W3_F,outeriter*18 + inneriter*400);
2472 }