Fujitsu Sparc64 acceleration and general fixes for non-x86 builds
[alexxy/gromacs.git] / src / gmxlib / nonbonded / nb_kernel_sparc64_hpc_ace_double / nb_kernel_ElecRF_VdwLJ_GeomW4W4_sparc64_hpc_ace_double.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012, by the GROMACS development team, led by
5  * David van der Spoel, Berk Hess, Erik Lindahl, and including many
6  * others, as listed in the AUTHORS file in the top-level source
7  * directory and at http://www.gromacs.org.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35 /*
36  * Note: this file was generated by the GROMACS sparc64_hpc_ace_double kernel generator.
37  */
38 #ifdef HAVE_CONFIG_H
39 #include <config.h>
40 #endif
41
42 #include <math.h>
43
44 #include "../nb_kernel.h"
45 #include "types/simple.h"
46 #include "vec.h"
47 #include "nrnb.h"
48
49 #include "kernelutil_sparc64_hpc_ace_double.h"
50
51 /*
52  * Gromacs nonbonded kernel:   nb_kernel_ElecRF_VdwLJ_GeomW4W4_VF_sparc64_hpc_ace_double
53  * Electrostatics interaction: ReactionField
54  * VdW interaction:            LennardJones
55  * Geometry:                   Water4-Water4
56  * Calculate force/pot:        PotentialAndForce
57  */
58 void
59 nb_kernel_ElecRF_VdwLJ_GeomW4W4_VF_sparc64_hpc_ace_double
60                     (t_nblist * gmx_restrict                nlist,
61                      rvec * gmx_restrict                    xx,
62                      rvec * gmx_restrict                    ff,
63                      t_forcerec * gmx_restrict              fr,
64                      t_mdatoms * gmx_restrict               mdatoms,
65                      nb_kernel_data_t * gmx_restrict        kernel_data,
66                      t_nrnb * gmx_restrict                  nrnb)
67 {
68     /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
69      * just 0 for non-waters.
70      * Suffixes A,B refer to j loop unrolling done with double precision SIMD, e.g. for the two different
71      * jnr indices corresponding to data put in the four positions in the SIMD register.
72      */
73     int              i_shift_offset,i_coord_offset,outeriter,inneriter;
74     int              j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
75     int              jnrA,jnrB;
76     int              j_coord_offsetA,j_coord_offsetB;
77     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
78     real             rcutoff_scalar;
79     real             *shiftvec,*fshift,*x,*f;
80     _fjsp_v2r8       tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
81     int              vdwioffset0;
82     _fjsp_v2r8       ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
83     int              vdwioffset1;
84     _fjsp_v2r8       ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
85     int              vdwioffset2;
86     _fjsp_v2r8       ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
87     int              vdwioffset3;
88     _fjsp_v2r8       ix3,iy3,iz3,fix3,fiy3,fiz3,iq3,isai3;
89     int              vdwjidx0A,vdwjidx0B;
90     _fjsp_v2r8       jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
91     int              vdwjidx1A,vdwjidx1B;
92     _fjsp_v2r8       jx1,jy1,jz1,fjx1,fjy1,fjz1,jq1,isaj1;
93     int              vdwjidx2A,vdwjidx2B;
94     _fjsp_v2r8       jx2,jy2,jz2,fjx2,fjy2,fjz2,jq2,isaj2;
95     int              vdwjidx3A,vdwjidx3B;
96     _fjsp_v2r8       jx3,jy3,jz3,fjx3,fjy3,fjz3,jq3,isaj3;
97     _fjsp_v2r8       dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
98     _fjsp_v2r8       dx11,dy11,dz11,rsq11,rinv11,rinvsq11,r11,qq11,c6_11,c12_11;
99     _fjsp_v2r8       dx12,dy12,dz12,rsq12,rinv12,rinvsq12,r12,qq12,c6_12,c12_12;
100     _fjsp_v2r8       dx13,dy13,dz13,rsq13,rinv13,rinvsq13,r13,qq13,c6_13,c12_13;
101     _fjsp_v2r8       dx21,dy21,dz21,rsq21,rinv21,rinvsq21,r21,qq21,c6_21,c12_21;
102     _fjsp_v2r8       dx22,dy22,dz22,rsq22,rinv22,rinvsq22,r22,qq22,c6_22,c12_22;
103     _fjsp_v2r8       dx23,dy23,dz23,rsq23,rinv23,rinvsq23,r23,qq23,c6_23,c12_23;
104     _fjsp_v2r8       dx31,dy31,dz31,rsq31,rinv31,rinvsq31,r31,qq31,c6_31,c12_31;
105     _fjsp_v2r8       dx32,dy32,dz32,rsq32,rinv32,rinvsq32,r32,qq32,c6_32,c12_32;
106     _fjsp_v2r8       dx33,dy33,dz33,rsq33,rinv33,rinvsq33,r33,qq33,c6_33,c12_33;
107     _fjsp_v2r8       velec,felec,velecsum,facel,crf,krf,krf2;
108     real             *charge;
109     int              nvdwtype;
110     _fjsp_v2r8       rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
111     int              *vdwtype;
112     real             *vdwparam;
113     _fjsp_v2r8       one_sixth   = gmx_fjsp_set1_v2r8(1.0/6.0);
114     _fjsp_v2r8       one_twelfth = gmx_fjsp_set1_v2r8(1.0/12.0);
115     _fjsp_v2r8       itab_tmp;
116     _fjsp_v2r8       dummy_mask,cutoff_mask;
117     _fjsp_v2r8       one     = gmx_fjsp_set1_v2r8(1.0);
118     _fjsp_v2r8       two     = gmx_fjsp_set1_v2r8(2.0);
119     union { _fjsp_v2r8 simd; long long int i[2]; } vfconv,gbconv,ewconv;
120
121     x                = xx[0];
122     f                = ff[0];
123
124     nri              = nlist->nri;
125     iinr             = nlist->iinr;
126     jindex           = nlist->jindex;
127     jjnr             = nlist->jjnr;
128     shiftidx         = nlist->shift;
129     gid              = nlist->gid;
130     shiftvec         = fr->shift_vec[0];
131     fshift           = fr->fshift[0];
132     facel            = gmx_fjsp_set1_v2r8(fr->epsfac);
133     charge           = mdatoms->chargeA;
134     krf              = gmx_fjsp_set1_v2r8(fr->ic->k_rf);
135     krf2             = gmx_fjsp_set1_v2r8(fr->ic->k_rf*2.0);
136     crf              = gmx_fjsp_set1_v2r8(fr->ic->c_rf);
137     nvdwtype         = fr->ntype;
138     vdwparam         = fr->nbfp;
139     vdwtype          = mdatoms->typeA;
140
141     /* Setup water-specific parameters */
142     inr              = nlist->iinr[0];
143     iq1              = _fjsp_mul_v2r8(facel,gmx_fjsp_set1_v2r8(charge[inr+1]));
144     iq2              = _fjsp_mul_v2r8(facel,gmx_fjsp_set1_v2r8(charge[inr+2]));
145     iq3              = _fjsp_mul_v2r8(facel,gmx_fjsp_set1_v2r8(charge[inr+3]));
146     vdwioffset0      = 2*nvdwtype*vdwtype[inr+0];
147
148     jq1              = gmx_fjsp_set1_v2r8(charge[inr+1]);
149     jq2              = gmx_fjsp_set1_v2r8(charge[inr+2]);
150     jq3              = gmx_fjsp_set1_v2r8(charge[inr+3]);
151     vdwjidx0A        = 2*vdwtype[inr+0];
152     c6_00            = gmx_fjsp_set1_v2r8(vdwparam[vdwioffset0+vdwjidx0A]);
153     c12_00           = gmx_fjsp_set1_v2r8(vdwparam[vdwioffset0+vdwjidx0A+1]);
154     qq11             = _fjsp_mul_v2r8(iq1,jq1);
155     qq12             = _fjsp_mul_v2r8(iq1,jq2);
156     qq13             = _fjsp_mul_v2r8(iq1,jq3);
157     qq21             = _fjsp_mul_v2r8(iq2,jq1);
158     qq22             = _fjsp_mul_v2r8(iq2,jq2);
159     qq23             = _fjsp_mul_v2r8(iq2,jq3);
160     qq31             = _fjsp_mul_v2r8(iq3,jq1);
161     qq32             = _fjsp_mul_v2r8(iq3,jq2);
162     qq33             = _fjsp_mul_v2r8(iq3,jq3);
163
164     /* Avoid stupid compiler warnings */
165     jnrA = jnrB = 0;
166     j_coord_offsetA = 0;
167     j_coord_offsetB = 0;
168
169     outeriter        = 0;
170     inneriter        = 0;
171
172     /* Start outer loop over neighborlists */
173     for(iidx=0; iidx<nri; iidx++)
174     {
175         /* Load shift vector for this list */
176         i_shift_offset   = DIM*shiftidx[iidx];
177
178         /* Load limits for loop over neighbors */
179         j_index_start    = jindex[iidx];
180         j_index_end      = jindex[iidx+1];
181
182         /* Get outer coordinate index */
183         inr              = iinr[iidx];
184         i_coord_offset   = DIM*inr;
185
186         /* Load i particle coords and add shift vector */
187         gmx_fjsp_load_shift_and_4rvec_broadcast_v2r8(shiftvec+i_shift_offset,x+i_coord_offset,
188                                                  &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
189
190         fix0             = _fjsp_setzero_v2r8();
191         fiy0             = _fjsp_setzero_v2r8();
192         fiz0             = _fjsp_setzero_v2r8();
193         fix1             = _fjsp_setzero_v2r8();
194         fiy1             = _fjsp_setzero_v2r8();
195         fiz1             = _fjsp_setzero_v2r8();
196         fix2             = _fjsp_setzero_v2r8();
197         fiy2             = _fjsp_setzero_v2r8();
198         fiz2             = _fjsp_setzero_v2r8();
199         fix3             = _fjsp_setzero_v2r8();
200         fiy3             = _fjsp_setzero_v2r8();
201         fiz3             = _fjsp_setzero_v2r8();
202
203         /* Reset potential sums */
204         velecsum         = _fjsp_setzero_v2r8();
205         vvdwsum          = _fjsp_setzero_v2r8();
206
207         /* Start inner kernel loop */
208         for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
209         {
210
211             /* Get j neighbor index, and coordinate index */
212             jnrA             = jjnr[jidx];
213             jnrB             = jjnr[jidx+1];
214             j_coord_offsetA  = DIM*jnrA;
215             j_coord_offsetB  = DIM*jnrB;
216
217             /* load j atom coordinates */
218             gmx_fjsp_load_4rvec_2ptr_swizzle_v2r8(x+j_coord_offsetA,x+j_coord_offsetB,
219                                               &jx0,&jy0,&jz0,&jx1,&jy1,&jz1,&jx2,
220                                               &jy2,&jz2,&jx3,&jy3,&jz3);
221
222             /* Calculate displacement vector */
223             dx00             = _fjsp_sub_v2r8(ix0,jx0);
224             dy00             = _fjsp_sub_v2r8(iy0,jy0);
225             dz00             = _fjsp_sub_v2r8(iz0,jz0);
226             dx11             = _fjsp_sub_v2r8(ix1,jx1);
227             dy11             = _fjsp_sub_v2r8(iy1,jy1);
228             dz11             = _fjsp_sub_v2r8(iz1,jz1);
229             dx12             = _fjsp_sub_v2r8(ix1,jx2);
230             dy12             = _fjsp_sub_v2r8(iy1,jy2);
231             dz12             = _fjsp_sub_v2r8(iz1,jz2);
232             dx13             = _fjsp_sub_v2r8(ix1,jx3);
233             dy13             = _fjsp_sub_v2r8(iy1,jy3);
234             dz13             = _fjsp_sub_v2r8(iz1,jz3);
235             dx21             = _fjsp_sub_v2r8(ix2,jx1);
236             dy21             = _fjsp_sub_v2r8(iy2,jy1);
237             dz21             = _fjsp_sub_v2r8(iz2,jz1);
238             dx22             = _fjsp_sub_v2r8(ix2,jx2);
239             dy22             = _fjsp_sub_v2r8(iy2,jy2);
240             dz22             = _fjsp_sub_v2r8(iz2,jz2);
241             dx23             = _fjsp_sub_v2r8(ix2,jx3);
242             dy23             = _fjsp_sub_v2r8(iy2,jy3);
243             dz23             = _fjsp_sub_v2r8(iz2,jz3);
244             dx31             = _fjsp_sub_v2r8(ix3,jx1);
245             dy31             = _fjsp_sub_v2r8(iy3,jy1);
246             dz31             = _fjsp_sub_v2r8(iz3,jz1);
247             dx32             = _fjsp_sub_v2r8(ix3,jx2);
248             dy32             = _fjsp_sub_v2r8(iy3,jy2);
249             dz32             = _fjsp_sub_v2r8(iz3,jz2);
250             dx33             = _fjsp_sub_v2r8(ix3,jx3);
251             dy33             = _fjsp_sub_v2r8(iy3,jy3);
252             dz33             = _fjsp_sub_v2r8(iz3,jz3);
253
254             /* Calculate squared distance and things based on it */
255             rsq00            = gmx_fjsp_calc_rsq_v2r8(dx00,dy00,dz00);
256             rsq11            = gmx_fjsp_calc_rsq_v2r8(dx11,dy11,dz11);
257             rsq12            = gmx_fjsp_calc_rsq_v2r8(dx12,dy12,dz12);
258             rsq13            = gmx_fjsp_calc_rsq_v2r8(dx13,dy13,dz13);
259             rsq21            = gmx_fjsp_calc_rsq_v2r8(dx21,dy21,dz21);
260             rsq22            = gmx_fjsp_calc_rsq_v2r8(dx22,dy22,dz22);
261             rsq23            = gmx_fjsp_calc_rsq_v2r8(dx23,dy23,dz23);
262             rsq31            = gmx_fjsp_calc_rsq_v2r8(dx31,dy31,dz31);
263             rsq32            = gmx_fjsp_calc_rsq_v2r8(dx32,dy32,dz32);
264             rsq33            = gmx_fjsp_calc_rsq_v2r8(dx33,dy33,dz33);
265
266             rinv11           = gmx_fjsp_invsqrt_v2r8(rsq11);
267             rinv12           = gmx_fjsp_invsqrt_v2r8(rsq12);
268             rinv13           = gmx_fjsp_invsqrt_v2r8(rsq13);
269             rinv21           = gmx_fjsp_invsqrt_v2r8(rsq21);
270             rinv22           = gmx_fjsp_invsqrt_v2r8(rsq22);
271             rinv23           = gmx_fjsp_invsqrt_v2r8(rsq23);
272             rinv31           = gmx_fjsp_invsqrt_v2r8(rsq31);
273             rinv32           = gmx_fjsp_invsqrt_v2r8(rsq32);
274             rinv33           = gmx_fjsp_invsqrt_v2r8(rsq33);
275
276             rinvsq00         = gmx_fjsp_inv_v2r8(rsq00);
277             rinvsq11         = _fjsp_mul_v2r8(rinv11,rinv11);
278             rinvsq12         = _fjsp_mul_v2r8(rinv12,rinv12);
279             rinvsq13         = _fjsp_mul_v2r8(rinv13,rinv13);
280             rinvsq21         = _fjsp_mul_v2r8(rinv21,rinv21);
281             rinvsq22         = _fjsp_mul_v2r8(rinv22,rinv22);
282             rinvsq23         = _fjsp_mul_v2r8(rinv23,rinv23);
283             rinvsq31         = _fjsp_mul_v2r8(rinv31,rinv31);
284             rinvsq32         = _fjsp_mul_v2r8(rinv32,rinv32);
285             rinvsq33         = _fjsp_mul_v2r8(rinv33,rinv33);
286
287             fjx0             = _fjsp_setzero_v2r8();
288             fjy0             = _fjsp_setzero_v2r8();
289             fjz0             = _fjsp_setzero_v2r8();
290             fjx1             = _fjsp_setzero_v2r8();
291             fjy1             = _fjsp_setzero_v2r8();
292             fjz1             = _fjsp_setzero_v2r8();
293             fjx2             = _fjsp_setzero_v2r8();
294             fjy2             = _fjsp_setzero_v2r8();
295             fjz2             = _fjsp_setzero_v2r8();
296             fjx3             = _fjsp_setzero_v2r8();
297             fjy3             = _fjsp_setzero_v2r8();
298             fjz3             = _fjsp_setzero_v2r8();
299
300             /**************************
301              * CALCULATE INTERACTIONS *
302              **************************/
303
304             /* LENNARD-JONES DISPERSION/REPULSION */
305
306             rinvsix          = _fjsp_mul_v2r8(_fjsp_mul_v2r8(rinvsq00,rinvsq00),rinvsq00);
307             vvdw6            = _fjsp_mul_v2r8(c6_00,rinvsix);
308             vvdw12           = _fjsp_mul_v2r8(c12_00,_fjsp_mul_v2r8(rinvsix,rinvsix));
309             vvdw             = _fjsp_msub_v2r8( vvdw12,one_twelfth, _fjsp_mul_v2r8(vvdw6,one_sixth) );
310             fvdw             = _fjsp_mul_v2r8(_fjsp_sub_v2r8(vvdw12,vvdw6),rinvsq00);
311
312             /* Update potential sum for this i atom from the interaction with this j atom. */
313             vvdwsum          = _fjsp_add_v2r8(vvdwsum,vvdw);
314
315             fscal            = fvdw;
316
317             /* Update vectorial force */
318             fix0             = _fjsp_madd_v2r8(dx00,fscal,fix0);
319             fiy0             = _fjsp_madd_v2r8(dy00,fscal,fiy0);
320             fiz0             = _fjsp_madd_v2r8(dz00,fscal,fiz0);
321             
322             fjx0             = _fjsp_madd_v2r8(dx00,fscal,fjx0);
323             fjy0             = _fjsp_madd_v2r8(dy00,fscal,fjy0);
324             fjz0             = _fjsp_madd_v2r8(dz00,fscal,fjz0);
325
326             /**************************
327              * CALCULATE INTERACTIONS *
328              **************************/
329
330             /* REACTION-FIELD ELECTROSTATICS */
331             velec            = _fjsp_mul_v2r8(qq11,_fjsp_sub_v2r8(_fjsp_madd_v2r8(krf,rsq11,rinv11),crf));
332             felec            = _fjsp_mul_v2r8(qq11,_fjsp_msub_v2r8(rinv11,rinvsq11,krf2));
333
334             /* Update potential sum for this i atom from the interaction with this j atom. */
335             velecsum         = _fjsp_add_v2r8(velecsum,velec);
336
337             fscal            = felec;
338
339             /* Update vectorial force */
340             fix1             = _fjsp_madd_v2r8(dx11,fscal,fix1);
341             fiy1             = _fjsp_madd_v2r8(dy11,fscal,fiy1);
342             fiz1             = _fjsp_madd_v2r8(dz11,fscal,fiz1);
343             
344             fjx1             = _fjsp_madd_v2r8(dx11,fscal,fjx1);
345             fjy1             = _fjsp_madd_v2r8(dy11,fscal,fjy1);
346             fjz1             = _fjsp_madd_v2r8(dz11,fscal,fjz1);
347
348             /**************************
349              * CALCULATE INTERACTIONS *
350              **************************/
351
352             /* REACTION-FIELD ELECTROSTATICS */
353             velec            = _fjsp_mul_v2r8(qq12,_fjsp_sub_v2r8(_fjsp_madd_v2r8(krf,rsq12,rinv12),crf));
354             felec            = _fjsp_mul_v2r8(qq12,_fjsp_msub_v2r8(rinv12,rinvsq12,krf2));
355
356             /* Update potential sum for this i atom from the interaction with this j atom. */
357             velecsum         = _fjsp_add_v2r8(velecsum,velec);
358
359             fscal            = felec;
360
361             /* Update vectorial force */
362             fix1             = _fjsp_madd_v2r8(dx12,fscal,fix1);
363             fiy1             = _fjsp_madd_v2r8(dy12,fscal,fiy1);
364             fiz1             = _fjsp_madd_v2r8(dz12,fscal,fiz1);
365             
366             fjx2             = _fjsp_madd_v2r8(dx12,fscal,fjx2);
367             fjy2             = _fjsp_madd_v2r8(dy12,fscal,fjy2);
368             fjz2             = _fjsp_madd_v2r8(dz12,fscal,fjz2);
369
370             /**************************
371              * CALCULATE INTERACTIONS *
372              **************************/
373
374             /* REACTION-FIELD ELECTROSTATICS */
375             velec            = _fjsp_mul_v2r8(qq13,_fjsp_sub_v2r8(_fjsp_madd_v2r8(krf,rsq13,rinv13),crf));
376             felec            = _fjsp_mul_v2r8(qq13,_fjsp_msub_v2r8(rinv13,rinvsq13,krf2));
377
378             /* Update potential sum for this i atom from the interaction with this j atom. */
379             velecsum         = _fjsp_add_v2r8(velecsum,velec);
380
381             fscal            = felec;
382
383             /* Update vectorial force */
384             fix1             = _fjsp_madd_v2r8(dx13,fscal,fix1);
385             fiy1             = _fjsp_madd_v2r8(dy13,fscal,fiy1);
386             fiz1             = _fjsp_madd_v2r8(dz13,fscal,fiz1);
387             
388             fjx3             = _fjsp_madd_v2r8(dx13,fscal,fjx3);
389             fjy3             = _fjsp_madd_v2r8(dy13,fscal,fjy3);
390             fjz3             = _fjsp_madd_v2r8(dz13,fscal,fjz3);
391
392             /**************************
393              * CALCULATE INTERACTIONS *
394              **************************/
395
396             /* REACTION-FIELD ELECTROSTATICS */
397             velec            = _fjsp_mul_v2r8(qq21,_fjsp_sub_v2r8(_fjsp_madd_v2r8(krf,rsq21,rinv21),crf));
398             felec            = _fjsp_mul_v2r8(qq21,_fjsp_msub_v2r8(rinv21,rinvsq21,krf2));
399
400             /* Update potential sum for this i atom from the interaction with this j atom. */
401             velecsum         = _fjsp_add_v2r8(velecsum,velec);
402
403             fscal            = felec;
404
405             /* Update vectorial force */
406             fix2             = _fjsp_madd_v2r8(dx21,fscal,fix2);
407             fiy2             = _fjsp_madd_v2r8(dy21,fscal,fiy2);
408             fiz2             = _fjsp_madd_v2r8(dz21,fscal,fiz2);
409             
410             fjx1             = _fjsp_madd_v2r8(dx21,fscal,fjx1);
411             fjy1             = _fjsp_madd_v2r8(dy21,fscal,fjy1);
412             fjz1             = _fjsp_madd_v2r8(dz21,fscal,fjz1);
413
414             /**************************
415              * CALCULATE INTERACTIONS *
416              **************************/
417
418             /* REACTION-FIELD ELECTROSTATICS */
419             velec            = _fjsp_mul_v2r8(qq22,_fjsp_sub_v2r8(_fjsp_madd_v2r8(krf,rsq22,rinv22),crf));
420             felec            = _fjsp_mul_v2r8(qq22,_fjsp_msub_v2r8(rinv22,rinvsq22,krf2));
421
422             /* Update potential sum for this i atom from the interaction with this j atom. */
423             velecsum         = _fjsp_add_v2r8(velecsum,velec);
424
425             fscal            = felec;
426
427             /* Update vectorial force */
428             fix2             = _fjsp_madd_v2r8(dx22,fscal,fix2);
429             fiy2             = _fjsp_madd_v2r8(dy22,fscal,fiy2);
430             fiz2             = _fjsp_madd_v2r8(dz22,fscal,fiz2);
431             
432             fjx2             = _fjsp_madd_v2r8(dx22,fscal,fjx2);
433             fjy2             = _fjsp_madd_v2r8(dy22,fscal,fjy2);
434             fjz2             = _fjsp_madd_v2r8(dz22,fscal,fjz2);
435
436             /**************************
437              * CALCULATE INTERACTIONS *
438              **************************/
439
440             /* REACTION-FIELD ELECTROSTATICS */
441             velec            = _fjsp_mul_v2r8(qq23,_fjsp_sub_v2r8(_fjsp_madd_v2r8(krf,rsq23,rinv23),crf));
442             felec            = _fjsp_mul_v2r8(qq23,_fjsp_msub_v2r8(rinv23,rinvsq23,krf2));
443
444             /* Update potential sum for this i atom from the interaction with this j atom. */
445             velecsum         = _fjsp_add_v2r8(velecsum,velec);
446
447             fscal            = felec;
448
449             /* Update vectorial force */
450             fix2             = _fjsp_madd_v2r8(dx23,fscal,fix2);
451             fiy2             = _fjsp_madd_v2r8(dy23,fscal,fiy2);
452             fiz2             = _fjsp_madd_v2r8(dz23,fscal,fiz2);
453             
454             fjx3             = _fjsp_madd_v2r8(dx23,fscal,fjx3);
455             fjy3             = _fjsp_madd_v2r8(dy23,fscal,fjy3);
456             fjz3             = _fjsp_madd_v2r8(dz23,fscal,fjz3);
457
458             /**************************
459              * CALCULATE INTERACTIONS *
460              **************************/
461
462             /* REACTION-FIELD ELECTROSTATICS */
463             velec            = _fjsp_mul_v2r8(qq31,_fjsp_sub_v2r8(_fjsp_madd_v2r8(krf,rsq31,rinv31),crf));
464             felec            = _fjsp_mul_v2r8(qq31,_fjsp_msub_v2r8(rinv31,rinvsq31,krf2));
465
466             /* Update potential sum for this i atom from the interaction with this j atom. */
467             velecsum         = _fjsp_add_v2r8(velecsum,velec);
468
469             fscal            = felec;
470
471             /* Update vectorial force */
472             fix3             = _fjsp_madd_v2r8(dx31,fscal,fix3);
473             fiy3             = _fjsp_madd_v2r8(dy31,fscal,fiy3);
474             fiz3             = _fjsp_madd_v2r8(dz31,fscal,fiz3);
475             
476             fjx1             = _fjsp_madd_v2r8(dx31,fscal,fjx1);
477             fjy1             = _fjsp_madd_v2r8(dy31,fscal,fjy1);
478             fjz1             = _fjsp_madd_v2r8(dz31,fscal,fjz1);
479
480             /**************************
481              * CALCULATE INTERACTIONS *
482              **************************/
483
484             /* REACTION-FIELD ELECTROSTATICS */
485             velec            = _fjsp_mul_v2r8(qq32,_fjsp_sub_v2r8(_fjsp_madd_v2r8(krf,rsq32,rinv32),crf));
486             felec            = _fjsp_mul_v2r8(qq32,_fjsp_msub_v2r8(rinv32,rinvsq32,krf2));
487
488             /* Update potential sum for this i atom from the interaction with this j atom. */
489             velecsum         = _fjsp_add_v2r8(velecsum,velec);
490
491             fscal            = felec;
492
493             /* Update vectorial force */
494             fix3             = _fjsp_madd_v2r8(dx32,fscal,fix3);
495             fiy3             = _fjsp_madd_v2r8(dy32,fscal,fiy3);
496             fiz3             = _fjsp_madd_v2r8(dz32,fscal,fiz3);
497             
498             fjx2             = _fjsp_madd_v2r8(dx32,fscal,fjx2);
499             fjy2             = _fjsp_madd_v2r8(dy32,fscal,fjy2);
500             fjz2             = _fjsp_madd_v2r8(dz32,fscal,fjz2);
501
502             /**************************
503              * CALCULATE INTERACTIONS *
504              **************************/
505
506             /* REACTION-FIELD ELECTROSTATICS */
507             velec            = _fjsp_mul_v2r8(qq33,_fjsp_sub_v2r8(_fjsp_madd_v2r8(krf,rsq33,rinv33),crf));
508             felec            = _fjsp_mul_v2r8(qq33,_fjsp_msub_v2r8(rinv33,rinvsq33,krf2));
509
510             /* Update potential sum for this i atom from the interaction with this j atom. */
511             velecsum         = _fjsp_add_v2r8(velecsum,velec);
512
513             fscal            = felec;
514
515             /* Update vectorial force */
516             fix3             = _fjsp_madd_v2r8(dx33,fscal,fix3);
517             fiy3             = _fjsp_madd_v2r8(dy33,fscal,fiy3);
518             fiz3             = _fjsp_madd_v2r8(dz33,fscal,fiz3);
519             
520             fjx3             = _fjsp_madd_v2r8(dx33,fscal,fjx3);
521             fjy3             = _fjsp_madd_v2r8(dy33,fscal,fjy3);
522             fjz3             = _fjsp_madd_v2r8(dz33,fscal,fjz3);
523
524             gmx_fjsp_decrement_4rvec_2ptr_swizzle_v2r8(f+j_coord_offsetA,f+j_coord_offsetB,fjx0,fjy0,fjz0,fjx1,fjy1,fjz1,fjx2,fjy2,fjz2,fjx3,fjy3,fjz3);
525
526             /* Inner loop uses 353 flops */
527         }
528
529         if(jidx<j_index_end)
530         {
531
532             jnrA             = jjnr[jidx];
533             j_coord_offsetA  = DIM*jnrA;
534
535             /* load j atom coordinates */
536             gmx_fjsp_load_4rvec_1ptr_swizzle_v2r8(x+j_coord_offsetA,
537                                               &jx0,&jy0,&jz0,&jx1,&jy1,&jz1,&jx2,
538                                               &jy2,&jz2,&jx3,&jy3,&jz3);
539
540             /* Calculate displacement vector */
541             dx00             = _fjsp_sub_v2r8(ix0,jx0);
542             dy00             = _fjsp_sub_v2r8(iy0,jy0);
543             dz00             = _fjsp_sub_v2r8(iz0,jz0);
544             dx11             = _fjsp_sub_v2r8(ix1,jx1);
545             dy11             = _fjsp_sub_v2r8(iy1,jy1);
546             dz11             = _fjsp_sub_v2r8(iz1,jz1);
547             dx12             = _fjsp_sub_v2r8(ix1,jx2);
548             dy12             = _fjsp_sub_v2r8(iy1,jy2);
549             dz12             = _fjsp_sub_v2r8(iz1,jz2);
550             dx13             = _fjsp_sub_v2r8(ix1,jx3);
551             dy13             = _fjsp_sub_v2r8(iy1,jy3);
552             dz13             = _fjsp_sub_v2r8(iz1,jz3);
553             dx21             = _fjsp_sub_v2r8(ix2,jx1);
554             dy21             = _fjsp_sub_v2r8(iy2,jy1);
555             dz21             = _fjsp_sub_v2r8(iz2,jz1);
556             dx22             = _fjsp_sub_v2r8(ix2,jx2);
557             dy22             = _fjsp_sub_v2r8(iy2,jy2);
558             dz22             = _fjsp_sub_v2r8(iz2,jz2);
559             dx23             = _fjsp_sub_v2r8(ix2,jx3);
560             dy23             = _fjsp_sub_v2r8(iy2,jy3);
561             dz23             = _fjsp_sub_v2r8(iz2,jz3);
562             dx31             = _fjsp_sub_v2r8(ix3,jx1);
563             dy31             = _fjsp_sub_v2r8(iy3,jy1);
564             dz31             = _fjsp_sub_v2r8(iz3,jz1);
565             dx32             = _fjsp_sub_v2r8(ix3,jx2);
566             dy32             = _fjsp_sub_v2r8(iy3,jy2);
567             dz32             = _fjsp_sub_v2r8(iz3,jz2);
568             dx33             = _fjsp_sub_v2r8(ix3,jx3);
569             dy33             = _fjsp_sub_v2r8(iy3,jy3);
570             dz33             = _fjsp_sub_v2r8(iz3,jz3);
571
572             /* Calculate squared distance and things based on it */
573             rsq00            = gmx_fjsp_calc_rsq_v2r8(dx00,dy00,dz00);
574             rsq11            = gmx_fjsp_calc_rsq_v2r8(dx11,dy11,dz11);
575             rsq12            = gmx_fjsp_calc_rsq_v2r8(dx12,dy12,dz12);
576             rsq13            = gmx_fjsp_calc_rsq_v2r8(dx13,dy13,dz13);
577             rsq21            = gmx_fjsp_calc_rsq_v2r8(dx21,dy21,dz21);
578             rsq22            = gmx_fjsp_calc_rsq_v2r8(dx22,dy22,dz22);
579             rsq23            = gmx_fjsp_calc_rsq_v2r8(dx23,dy23,dz23);
580             rsq31            = gmx_fjsp_calc_rsq_v2r8(dx31,dy31,dz31);
581             rsq32            = gmx_fjsp_calc_rsq_v2r8(dx32,dy32,dz32);
582             rsq33            = gmx_fjsp_calc_rsq_v2r8(dx33,dy33,dz33);
583
584             rinv11           = gmx_fjsp_invsqrt_v2r8(rsq11);
585             rinv12           = gmx_fjsp_invsqrt_v2r8(rsq12);
586             rinv13           = gmx_fjsp_invsqrt_v2r8(rsq13);
587             rinv21           = gmx_fjsp_invsqrt_v2r8(rsq21);
588             rinv22           = gmx_fjsp_invsqrt_v2r8(rsq22);
589             rinv23           = gmx_fjsp_invsqrt_v2r8(rsq23);
590             rinv31           = gmx_fjsp_invsqrt_v2r8(rsq31);
591             rinv32           = gmx_fjsp_invsqrt_v2r8(rsq32);
592             rinv33           = gmx_fjsp_invsqrt_v2r8(rsq33);
593
594             rinvsq00         = gmx_fjsp_inv_v2r8(rsq00);
595             rinvsq11         = _fjsp_mul_v2r8(rinv11,rinv11);
596             rinvsq12         = _fjsp_mul_v2r8(rinv12,rinv12);
597             rinvsq13         = _fjsp_mul_v2r8(rinv13,rinv13);
598             rinvsq21         = _fjsp_mul_v2r8(rinv21,rinv21);
599             rinvsq22         = _fjsp_mul_v2r8(rinv22,rinv22);
600             rinvsq23         = _fjsp_mul_v2r8(rinv23,rinv23);
601             rinvsq31         = _fjsp_mul_v2r8(rinv31,rinv31);
602             rinvsq32         = _fjsp_mul_v2r8(rinv32,rinv32);
603             rinvsq33         = _fjsp_mul_v2r8(rinv33,rinv33);
604
605             fjx0             = _fjsp_setzero_v2r8();
606             fjy0             = _fjsp_setzero_v2r8();
607             fjz0             = _fjsp_setzero_v2r8();
608             fjx1             = _fjsp_setzero_v2r8();
609             fjy1             = _fjsp_setzero_v2r8();
610             fjz1             = _fjsp_setzero_v2r8();
611             fjx2             = _fjsp_setzero_v2r8();
612             fjy2             = _fjsp_setzero_v2r8();
613             fjz2             = _fjsp_setzero_v2r8();
614             fjx3             = _fjsp_setzero_v2r8();
615             fjy3             = _fjsp_setzero_v2r8();
616             fjz3             = _fjsp_setzero_v2r8();
617
618             /**************************
619              * CALCULATE INTERACTIONS *
620              **************************/
621
622             /* LENNARD-JONES DISPERSION/REPULSION */
623
624             rinvsix          = _fjsp_mul_v2r8(_fjsp_mul_v2r8(rinvsq00,rinvsq00),rinvsq00);
625             vvdw6            = _fjsp_mul_v2r8(c6_00,rinvsix);
626             vvdw12           = _fjsp_mul_v2r8(c12_00,_fjsp_mul_v2r8(rinvsix,rinvsix));
627             vvdw             = _fjsp_msub_v2r8( vvdw12,one_twelfth, _fjsp_mul_v2r8(vvdw6,one_sixth) );
628             fvdw             = _fjsp_mul_v2r8(_fjsp_sub_v2r8(vvdw12,vvdw6),rinvsq00);
629
630             /* Update potential sum for this i atom from the interaction with this j atom. */
631             vvdw             = _fjsp_unpacklo_v2r8(vvdw,_fjsp_setzero_v2r8());
632             vvdwsum          = _fjsp_add_v2r8(vvdwsum,vvdw);
633
634             fscal            = fvdw;
635
636             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
637
638             /* Update vectorial force */
639             fix0             = _fjsp_madd_v2r8(dx00,fscal,fix0);
640             fiy0             = _fjsp_madd_v2r8(dy00,fscal,fiy0);
641             fiz0             = _fjsp_madd_v2r8(dz00,fscal,fiz0);
642             
643             fjx0             = _fjsp_madd_v2r8(dx00,fscal,fjx0);
644             fjy0             = _fjsp_madd_v2r8(dy00,fscal,fjy0);
645             fjz0             = _fjsp_madd_v2r8(dz00,fscal,fjz0);
646
647             /**************************
648              * CALCULATE INTERACTIONS *
649              **************************/
650
651             /* REACTION-FIELD ELECTROSTATICS */
652             velec            = _fjsp_mul_v2r8(qq11,_fjsp_sub_v2r8(_fjsp_madd_v2r8(krf,rsq11,rinv11),crf));
653             felec            = _fjsp_mul_v2r8(qq11,_fjsp_msub_v2r8(rinv11,rinvsq11,krf2));
654
655             /* Update potential sum for this i atom from the interaction with this j atom. */
656             velec            = _fjsp_unpacklo_v2r8(velec,_fjsp_setzero_v2r8());
657             velecsum         = _fjsp_add_v2r8(velecsum,velec);
658
659             fscal            = felec;
660
661             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
662
663             /* Update vectorial force */
664             fix1             = _fjsp_madd_v2r8(dx11,fscal,fix1);
665             fiy1             = _fjsp_madd_v2r8(dy11,fscal,fiy1);
666             fiz1             = _fjsp_madd_v2r8(dz11,fscal,fiz1);
667             
668             fjx1             = _fjsp_madd_v2r8(dx11,fscal,fjx1);
669             fjy1             = _fjsp_madd_v2r8(dy11,fscal,fjy1);
670             fjz1             = _fjsp_madd_v2r8(dz11,fscal,fjz1);
671
672             /**************************
673              * CALCULATE INTERACTIONS *
674              **************************/
675
676             /* REACTION-FIELD ELECTROSTATICS */
677             velec            = _fjsp_mul_v2r8(qq12,_fjsp_sub_v2r8(_fjsp_madd_v2r8(krf,rsq12,rinv12),crf));
678             felec            = _fjsp_mul_v2r8(qq12,_fjsp_msub_v2r8(rinv12,rinvsq12,krf2));
679
680             /* Update potential sum for this i atom from the interaction with this j atom. */
681             velec            = _fjsp_unpacklo_v2r8(velec,_fjsp_setzero_v2r8());
682             velecsum         = _fjsp_add_v2r8(velecsum,velec);
683
684             fscal            = felec;
685
686             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
687
688             /* Update vectorial force */
689             fix1             = _fjsp_madd_v2r8(dx12,fscal,fix1);
690             fiy1             = _fjsp_madd_v2r8(dy12,fscal,fiy1);
691             fiz1             = _fjsp_madd_v2r8(dz12,fscal,fiz1);
692             
693             fjx2             = _fjsp_madd_v2r8(dx12,fscal,fjx2);
694             fjy2             = _fjsp_madd_v2r8(dy12,fscal,fjy2);
695             fjz2             = _fjsp_madd_v2r8(dz12,fscal,fjz2);
696
697             /**************************
698              * CALCULATE INTERACTIONS *
699              **************************/
700
701             /* REACTION-FIELD ELECTROSTATICS */
702             velec            = _fjsp_mul_v2r8(qq13,_fjsp_sub_v2r8(_fjsp_madd_v2r8(krf,rsq13,rinv13),crf));
703             felec            = _fjsp_mul_v2r8(qq13,_fjsp_msub_v2r8(rinv13,rinvsq13,krf2));
704
705             /* Update potential sum for this i atom from the interaction with this j atom. */
706             velec            = _fjsp_unpacklo_v2r8(velec,_fjsp_setzero_v2r8());
707             velecsum         = _fjsp_add_v2r8(velecsum,velec);
708
709             fscal            = felec;
710
711             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
712
713             /* Update vectorial force */
714             fix1             = _fjsp_madd_v2r8(dx13,fscal,fix1);
715             fiy1             = _fjsp_madd_v2r8(dy13,fscal,fiy1);
716             fiz1             = _fjsp_madd_v2r8(dz13,fscal,fiz1);
717             
718             fjx3             = _fjsp_madd_v2r8(dx13,fscal,fjx3);
719             fjy3             = _fjsp_madd_v2r8(dy13,fscal,fjy3);
720             fjz3             = _fjsp_madd_v2r8(dz13,fscal,fjz3);
721
722             /**************************
723              * CALCULATE INTERACTIONS *
724              **************************/
725
726             /* REACTION-FIELD ELECTROSTATICS */
727             velec            = _fjsp_mul_v2r8(qq21,_fjsp_sub_v2r8(_fjsp_madd_v2r8(krf,rsq21,rinv21),crf));
728             felec            = _fjsp_mul_v2r8(qq21,_fjsp_msub_v2r8(rinv21,rinvsq21,krf2));
729
730             /* Update potential sum for this i atom from the interaction with this j atom. */
731             velec            = _fjsp_unpacklo_v2r8(velec,_fjsp_setzero_v2r8());
732             velecsum         = _fjsp_add_v2r8(velecsum,velec);
733
734             fscal            = felec;
735
736             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
737
738             /* Update vectorial force */
739             fix2             = _fjsp_madd_v2r8(dx21,fscal,fix2);
740             fiy2             = _fjsp_madd_v2r8(dy21,fscal,fiy2);
741             fiz2             = _fjsp_madd_v2r8(dz21,fscal,fiz2);
742             
743             fjx1             = _fjsp_madd_v2r8(dx21,fscal,fjx1);
744             fjy1             = _fjsp_madd_v2r8(dy21,fscal,fjy1);
745             fjz1             = _fjsp_madd_v2r8(dz21,fscal,fjz1);
746
747             /**************************
748              * CALCULATE INTERACTIONS *
749              **************************/
750
751             /* REACTION-FIELD ELECTROSTATICS */
752             velec            = _fjsp_mul_v2r8(qq22,_fjsp_sub_v2r8(_fjsp_madd_v2r8(krf,rsq22,rinv22),crf));
753             felec            = _fjsp_mul_v2r8(qq22,_fjsp_msub_v2r8(rinv22,rinvsq22,krf2));
754
755             /* Update potential sum for this i atom from the interaction with this j atom. */
756             velec            = _fjsp_unpacklo_v2r8(velec,_fjsp_setzero_v2r8());
757             velecsum         = _fjsp_add_v2r8(velecsum,velec);
758
759             fscal            = felec;
760
761             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
762
763             /* Update vectorial force */
764             fix2             = _fjsp_madd_v2r8(dx22,fscal,fix2);
765             fiy2             = _fjsp_madd_v2r8(dy22,fscal,fiy2);
766             fiz2             = _fjsp_madd_v2r8(dz22,fscal,fiz2);
767             
768             fjx2             = _fjsp_madd_v2r8(dx22,fscal,fjx2);
769             fjy2             = _fjsp_madd_v2r8(dy22,fscal,fjy2);
770             fjz2             = _fjsp_madd_v2r8(dz22,fscal,fjz2);
771
772             /**************************
773              * CALCULATE INTERACTIONS *
774              **************************/
775
776             /* REACTION-FIELD ELECTROSTATICS */
777             velec            = _fjsp_mul_v2r8(qq23,_fjsp_sub_v2r8(_fjsp_madd_v2r8(krf,rsq23,rinv23),crf));
778             felec            = _fjsp_mul_v2r8(qq23,_fjsp_msub_v2r8(rinv23,rinvsq23,krf2));
779
780             /* Update potential sum for this i atom from the interaction with this j atom. */
781             velec            = _fjsp_unpacklo_v2r8(velec,_fjsp_setzero_v2r8());
782             velecsum         = _fjsp_add_v2r8(velecsum,velec);
783
784             fscal            = felec;
785
786             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
787
788             /* Update vectorial force */
789             fix2             = _fjsp_madd_v2r8(dx23,fscal,fix2);
790             fiy2             = _fjsp_madd_v2r8(dy23,fscal,fiy2);
791             fiz2             = _fjsp_madd_v2r8(dz23,fscal,fiz2);
792             
793             fjx3             = _fjsp_madd_v2r8(dx23,fscal,fjx3);
794             fjy3             = _fjsp_madd_v2r8(dy23,fscal,fjy3);
795             fjz3             = _fjsp_madd_v2r8(dz23,fscal,fjz3);
796
797             /**************************
798              * CALCULATE INTERACTIONS *
799              **************************/
800
801             /* REACTION-FIELD ELECTROSTATICS */
802             velec            = _fjsp_mul_v2r8(qq31,_fjsp_sub_v2r8(_fjsp_madd_v2r8(krf,rsq31,rinv31),crf));
803             felec            = _fjsp_mul_v2r8(qq31,_fjsp_msub_v2r8(rinv31,rinvsq31,krf2));
804
805             /* Update potential sum for this i atom from the interaction with this j atom. */
806             velec            = _fjsp_unpacklo_v2r8(velec,_fjsp_setzero_v2r8());
807             velecsum         = _fjsp_add_v2r8(velecsum,velec);
808
809             fscal            = felec;
810
811             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
812
813             /* Update vectorial force */
814             fix3             = _fjsp_madd_v2r8(dx31,fscal,fix3);
815             fiy3             = _fjsp_madd_v2r8(dy31,fscal,fiy3);
816             fiz3             = _fjsp_madd_v2r8(dz31,fscal,fiz3);
817             
818             fjx1             = _fjsp_madd_v2r8(dx31,fscal,fjx1);
819             fjy1             = _fjsp_madd_v2r8(dy31,fscal,fjy1);
820             fjz1             = _fjsp_madd_v2r8(dz31,fscal,fjz1);
821
822             /**************************
823              * CALCULATE INTERACTIONS *
824              **************************/
825
826             /* REACTION-FIELD ELECTROSTATICS */
827             velec            = _fjsp_mul_v2r8(qq32,_fjsp_sub_v2r8(_fjsp_madd_v2r8(krf,rsq32,rinv32),crf));
828             felec            = _fjsp_mul_v2r8(qq32,_fjsp_msub_v2r8(rinv32,rinvsq32,krf2));
829
830             /* Update potential sum for this i atom from the interaction with this j atom. */
831             velec            = _fjsp_unpacklo_v2r8(velec,_fjsp_setzero_v2r8());
832             velecsum         = _fjsp_add_v2r8(velecsum,velec);
833
834             fscal            = felec;
835
836             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
837
838             /* Update vectorial force */
839             fix3             = _fjsp_madd_v2r8(dx32,fscal,fix3);
840             fiy3             = _fjsp_madd_v2r8(dy32,fscal,fiy3);
841             fiz3             = _fjsp_madd_v2r8(dz32,fscal,fiz3);
842             
843             fjx2             = _fjsp_madd_v2r8(dx32,fscal,fjx2);
844             fjy2             = _fjsp_madd_v2r8(dy32,fscal,fjy2);
845             fjz2             = _fjsp_madd_v2r8(dz32,fscal,fjz2);
846
847             /**************************
848              * CALCULATE INTERACTIONS *
849              **************************/
850
851             /* REACTION-FIELD ELECTROSTATICS */
852             velec            = _fjsp_mul_v2r8(qq33,_fjsp_sub_v2r8(_fjsp_madd_v2r8(krf,rsq33,rinv33),crf));
853             felec            = _fjsp_mul_v2r8(qq33,_fjsp_msub_v2r8(rinv33,rinvsq33,krf2));
854
855             /* Update potential sum for this i atom from the interaction with this j atom. */
856             velec            = _fjsp_unpacklo_v2r8(velec,_fjsp_setzero_v2r8());
857             velecsum         = _fjsp_add_v2r8(velecsum,velec);
858
859             fscal            = felec;
860
861             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
862
863             /* Update vectorial force */
864             fix3             = _fjsp_madd_v2r8(dx33,fscal,fix3);
865             fiy3             = _fjsp_madd_v2r8(dy33,fscal,fiy3);
866             fiz3             = _fjsp_madd_v2r8(dz33,fscal,fiz3);
867             
868             fjx3             = _fjsp_madd_v2r8(dx33,fscal,fjx3);
869             fjy3             = _fjsp_madd_v2r8(dy33,fscal,fjy3);
870             fjz3             = _fjsp_madd_v2r8(dz33,fscal,fjz3);
871
872             gmx_fjsp_decrement_4rvec_1ptr_swizzle_v2r8(f+j_coord_offsetA,fjx0,fjy0,fjz0,fjx1,fjy1,fjz1,fjx2,fjy2,fjz2,fjx3,fjy3,fjz3);
873
874             /* Inner loop uses 353 flops */
875         }
876
877         /* End of innermost loop */
878
879         gmx_fjsp_update_iforce_4atom_swizzle_v2r8(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
880                                               f+i_coord_offset,fshift+i_shift_offset);
881
882         ggid                        = gid[iidx];
883         /* Update potential energies */
884         gmx_fjsp_update_1pot_v2r8(velecsum,kernel_data->energygrp_elec+ggid);
885         gmx_fjsp_update_1pot_v2r8(vvdwsum,kernel_data->energygrp_vdw+ggid);
886
887         /* Increment number of inner iterations */
888         inneriter                  += j_index_end - j_index_start;
889
890         /* Outer loop uses 26 flops */
891     }
892
893     /* Increment number of outer iterations */
894     outeriter        += nri;
895
896     /* Update outer/inner flops */
897
898     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W4W4_VF,outeriter*26 + inneriter*353);
899 }
900 /*
901  * Gromacs nonbonded kernel:   nb_kernel_ElecRF_VdwLJ_GeomW4W4_F_sparc64_hpc_ace_double
902  * Electrostatics interaction: ReactionField
903  * VdW interaction:            LennardJones
904  * Geometry:                   Water4-Water4
905  * Calculate force/pot:        Force
906  */
907 void
908 nb_kernel_ElecRF_VdwLJ_GeomW4W4_F_sparc64_hpc_ace_double
909                     (t_nblist * gmx_restrict                nlist,
910                      rvec * gmx_restrict                    xx,
911                      rvec * gmx_restrict                    ff,
912                      t_forcerec * gmx_restrict              fr,
913                      t_mdatoms * gmx_restrict               mdatoms,
914                      nb_kernel_data_t * gmx_restrict        kernel_data,
915                      t_nrnb * gmx_restrict                  nrnb)
916 {
917     /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
918      * just 0 for non-waters.
919      * Suffixes A,B refer to j loop unrolling done with double precision SIMD, e.g. for the two different
920      * jnr indices corresponding to data put in the four positions in the SIMD register.
921      */
922     int              i_shift_offset,i_coord_offset,outeriter,inneriter;
923     int              j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
924     int              jnrA,jnrB;
925     int              j_coord_offsetA,j_coord_offsetB;
926     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
927     real             rcutoff_scalar;
928     real             *shiftvec,*fshift,*x,*f;
929     _fjsp_v2r8       tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
930     int              vdwioffset0;
931     _fjsp_v2r8       ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
932     int              vdwioffset1;
933     _fjsp_v2r8       ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
934     int              vdwioffset2;
935     _fjsp_v2r8       ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
936     int              vdwioffset3;
937     _fjsp_v2r8       ix3,iy3,iz3,fix3,fiy3,fiz3,iq3,isai3;
938     int              vdwjidx0A,vdwjidx0B;
939     _fjsp_v2r8       jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
940     int              vdwjidx1A,vdwjidx1B;
941     _fjsp_v2r8       jx1,jy1,jz1,fjx1,fjy1,fjz1,jq1,isaj1;
942     int              vdwjidx2A,vdwjidx2B;
943     _fjsp_v2r8       jx2,jy2,jz2,fjx2,fjy2,fjz2,jq2,isaj2;
944     int              vdwjidx3A,vdwjidx3B;
945     _fjsp_v2r8       jx3,jy3,jz3,fjx3,fjy3,fjz3,jq3,isaj3;
946     _fjsp_v2r8       dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
947     _fjsp_v2r8       dx11,dy11,dz11,rsq11,rinv11,rinvsq11,r11,qq11,c6_11,c12_11;
948     _fjsp_v2r8       dx12,dy12,dz12,rsq12,rinv12,rinvsq12,r12,qq12,c6_12,c12_12;
949     _fjsp_v2r8       dx13,dy13,dz13,rsq13,rinv13,rinvsq13,r13,qq13,c6_13,c12_13;
950     _fjsp_v2r8       dx21,dy21,dz21,rsq21,rinv21,rinvsq21,r21,qq21,c6_21,c12_21;
951     _fjsp_v2r8       dx22,dy22,dz22,rsq22,rinv22,rinvsq22,r22,qq22,c6_22,c12_22;
952     _fjsp_v2r8       dx23,dy23,dz23,rsq23,rinv23,rinvsq23,r23,qq23,c6_23,c12_23;
953     _fjsp_v2r8       dx31,dy31,dz31,rsq31,rinv31,rinvsq31,r31,qq31,c6_31,c12_31;
954     _fjsp_v2r8       dx32,dy32,dz32,rsq32,rinv32,rinvsq32,r32,qq32,c6_32,c12_32;
955     _fjsp_v2r8       dx33,dy33,dz33,rsq33,rinv33,rinvsq33,r33,qq33,c6_33,c12_33;
956     _fjsp_v2r8       velec,felec,velecsum,facel,crf,krf,krf2;
957     real             *charge;
958     int              nvdwtype;
959     _fjsp_v2r8       rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
960     int              *vdwtype;
961     real             *vdwparam;
962     _fjsp_v2r8       one_sixth   = gmx_fjsp_set1_v2r8(1.0/6.0);
963     _fjsp_v2r8       one_twelfth = gmx_fjsp_set1_v2r8(1.0/12.0);
964     _fjsp_v2r8       itab_tmp;
965     _fjsp_v2r8       dummy_mask,cutoff_mask;
966     _fjsp_v2r8       one     = gmx_fjsp_set1_v2r8(1.0);
967     _fjsp_v2r8       two     = gmx_fjsp_set1_v2r8(2.0);
968     union { _fjsp_v2r8 simd; long long int i[2]; } vfconv,gbconv,ewconv;
969
970     x                = xx[0];
971     f                = ff[0];
972
973     nri              = nlist->nri;
974     iinr             = nlist->iinr;
975     jindex           = nlist->jindex;
976     jjnr             = nlist->jjnr;
977     shiftidx         = nlist->shift;
978     gid              = nlist->gid;
979     shiftvec         = fr->shift_vec[0];
980     fshift           = fr->fshift[0];
981     facel            = gmx_fjsp_set1_v2r8(fr->epsfac);
982     charge           = mdatoms->chargeA;
983     krf              = gmx_fjsp_set1_v2r8(fr->ic->k_rf);
984     krf2             = gmx_fjsp_set1_v2r8(fr->ic->k_rf*2.0);
985     crf              = gmx_fjsp_set1_v2r8(fr->ic->c_rf);
986     nvdwtype         = fr->ntype;
987     vdwparam         = fr->nbfp;
988     vdwtype          = mdatoms->typeA;
989
990     /* Setup water-specific parameters */
991     inr              = nlist->iinr[0];
992     iq1              = _fjsp_mul_v2r8(facel,gmx_fjsp_set1_v2r8(charge[inr+1]));
993     iq2              = _fjsp_mul_v2r8(facel,gmx_fjsp_set1_v2r8(charge[inr+2]));
994     iq3              = _fjsp_mul_v2r8(facel,gmx_fjsp_set1_v2r8(charge[inr+3]));
995     vdwioffset0      = 2*nvdwtype*vdwtype[inr+0];
996
997     jq1              = gmx_fjsp_set1_v2r8(charge[inr+1]);
998     jq2              = gmx_fjsp_set1_v2r8(charge[inr+2]);
999     jq3              = gmx_fjsp_set1_v2r8(charge[inr+3]);
1000     vdwjidx0A        = 2*vdwtype[inr+0];
1001     c6_00            = gmx_fjsp_set1_v2r8(vdwparam[vdwioffset0+vdwjidx0A]);
1002     c12_00           = gmx_fjsp_set1_v2r8(vdwparam[vdwioffset0+vdwjidx0A+1]);
1003     qq11             = _fjsp_mul_v2r8(iq1,jq1);
1004     qq12             = _fjsp_mul_v2r8(iq1,jq2);
1005     qq13             = _fjsp_mul_v2r8(iq1,jq3);
1006     qq21             = _fjsp_mul_v2r8(iq2,jq1);
1007     qq22             = _fjsp_mul_v2r8(iq2,jq2);
1008     qq23             = _fjsp_mul_v2r8(iq2,jq3);
1009     qq31             = _fjsp_mul_v2r8(iq3,jq1);
1010     qq32             = _fjsp_mul_v2r8(iq3,jq2);
1011     qq33             = _fjsp_mul_v2r8(iq3,jq3);
1012
1013     /* Avoid stupid compiler warnings */
1014     jnrA = jnrB = 0;
1015     j_coord_offsetA = 0;
1016     j_coord_offsetB = 0;
1017
1018     outeriter        = 0;
1019     inneriter        = 0;
1020
1021     /* Start outer loop over neighborlists */
1022     for(iidx=0; iidx<nri; iidx++)
1023     {
1024         /* Load shift vector for this list */
1025         i_shift_offset   = DIM*shiftidx[iidx];
1026
1027         /* Load limits for loop over neighbors */
1028         j_index_start    = jindex[iidx];
1029         j_index_end      = jindex[iidx+1];
1030
1031         /* Get outer coordinate index */
1032         inr              = iinr[iidx];
1033         i_coord_offset   = DIM*inr;
1034
1035         /* Load i particle coords and add shift vector */
1036         gmx_fjsp_load_shift_and_4rvec_broadcast_v2r8(shiftvec+i_shift_offset,x+i_coord_offset,
1037                                                  &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
1038
1039         fix0             = _fjsp_setzero_v2r8();
1040         fiy0             = _fjsp_setzero_v2r8();
1041         fiz0             = _fjsp_setzero_v2r8();
1042         fix1             = _fjsp_setzero_v2r8();
1043         fiy1             = _fjsp_setzero_v2r8();
1044         fiz1             = _fjsp_setzero_v2r8();
1045         fix2             = _fjsp_setzero_v2r8();
1046         fiy2             = _fjsp_setzero_v2r8();
1047         fiz2             = _fjsp_setzero_v2r8();
1048         fix3             = _fjsp_setzero_v2r8();
1049         fiy3             = _fjsp_setzero_v2r8();
1050         fiz3             = _fjsp_setzero_v2r8();
1051
1052         /* Start inner kernel loop */
1053         for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
1054         {
1055
1056             /* Get j neighbor index, and coordinate index */
1057             jnrA             = jjnr[jidx];
1058             jnrB             = jjnr[jidx+1];
1059             j_coord_offsetA  = DIM*jnrA;
1060             j_coord_offsetB  = DIM*jnrB;
1061
1062             /* load j atom coordinates */
1063             gmx_fjsp_load_4rvec_2ptr_swizzle_v2r8(x+j_coord_offsetA,x+j_coord_offsetB,
1064                                               &jx0,&jy0,&jz0,&jx1,&jy1,&jz1,&jx2,
1065                                               &jy2,&jz2,&jx3,&jy3,&jz3);
1066
1067             /* Calculate displacement vector */
1068             dx00             = _fjsp_sub_v2r8(ix0,jx0);
1069             dy00             = _fjsp_sub_v2r8(iy0,jy0);
1070             dz00             = _fjsp_sub_v2r8(iz0,jz0);
1071             dx11             = _fjsp_sub_v2r8(ix1,jx1);
1072             dy11             = _fjsp_sub_v2r8(iy1,jy1);
1073             dz11             = _fjsp_sub_v2r8(iz1,jz1);
1074             dx12             = _fjsp_sub_v2r8(ix1,jx2);
1075             dy12             = _fjsp_sub_v2r8(iy1,jy2);
1076             dz12             = _fjsp_sub_v2r8(iz1,jz2);
1077             dx13             = _fjsp_sub_v2r8(ix1,jx3);
1078             dy13             = _fjsp_sub_v2r8(iy1,jy3);
1079             dz13             = _fjsp_sub_v2r8(iz1,jz3);
1080             dx21             = _fjsp_sub_v2r8(ix2,jx1);
1081             dy21             = _fjsp_sub_v2r8(iy2,jy1);
1082             dz21             = _fjsp_sub_v2r8(iz2,jz1);
1083             dx22             = _fjsp_sub_v2r8(ix2,jx2);
1084             dy22             = _fjsp_sub_v2r8(iy2,jy2);
1085             dz22             = _fjsp_sub_v2r8(iz2,jz2);
1086             dx23             = _fjsp_sub_v2r8(ix2,jx3);
1087             dy23             = _fjsp_sub_v2r8(iy2,jy3);
1088             dz23             = _fjsp_sub_v2r8(iz2,jz3);
1089             dx31             = _fjsp_sub_v2r8(ix3,jx1);
1090             dy31             = _fjsp_sub_v2r8(iy3,jy1);
1091             dz31             = _fjsp_sub_v2r8(iz3,jz1);
1092             dx32             = _fjsp_sub_v2r8(ix3,jx2);
1093             dy32             = _fjsp_sub_v2r8(iy3,jy2);
1094             dz32             = _fjsp_sub_v2r8(iz3,jz2);
1095             dx33             = _fjsp_sub_v2r8(ix3,jx3);
1096             dy33             = _fjsp_sub_v2r8(iy3,jy3);
1097             dz33             = _fjsp_sub_v2r8(iz3,jz3);
1098
1099             /* Calculate squared distance and things based on it */
1100             rsq00            = gmx_fjsp_calc_rsq_v2r8(dx00,dy00,dz00);
1101             rsq11            = gmx_fjsp_calc_rsq_v2r8(dx11,dy11,dz11);
1102             rsq12            = gmx_fjsp_calc_rsq_v2r8(dx12,dy12,dz12);
1103             rsq13            = gmx_fjsp_calc_rsq_v2r8(dx13,dy13,dz13);
1104             rsq21            = gmx_fjsp_calc_rsq_v2r8(dx21,dy21,dz21);
1105             rsq22            = gmx_fjsp_calc_rsq_v2r8(dx22,dy22,dz22);
1106             rsq23            = gmx_fjsp_calc_rsq_v2r8(dx23,dy23,dz23);
1107             rsq31            = gmx_fjsp_calc_rsq_v2r8(dx31,dy31,dz31);
1108             rsq32            = gmx_fjsp_calc_rsq_v2r8(dx32,dy32,dz32);
1109             rsq33            = gmx_fjsp_calc_rsq_v2r8(dx33,dy33,dz33);
1110
1111             rinv11           = gmx_fjsp_invsqrt_v2r8(rsq11);
1112             rinv12           = gmx_fjsp_invsqrt_v2r8(rsq12);
1113             rinv13           = gmx_fjsp_invsqrt_v2r8(rsq13);
1114             rinv21           = gmx_fjsp_invsqrt_v2r8(rsq21);
1115             rinv22           = gmx_fjsp_invsqrt_v2r8(rsq22);
1116             rinv23           = gmx_fjsp_invsqrt_v2r8(rsq23);
1117             rinv31           = gmx_fjsp_invsqrt_v2r8(rsq31);
1118             rinv32           = gmx_fjsp_invsqrt_v2r8(rsq32);
1119             rinv33           = gmx_fjsp_invsqrt_v2r8(rsq33);
1120
1121             rinvsq00         = gmx_fjsp_inv_v2r8(rsq00);
1122             rinvsq11         = _fjsp_mul_v2r8(rinv11,rinv11);
1123             rinvsq12         = _fjsp_mul_v2r8(rinv12,rinv12);
1124             rinvsq13         = _fjsp_mul_v2r8(rinv13,rinv13);
1125             rinvsq21         = _fjsp_mul_v2r8(rinv21,rinv21);
1126             rinvsq22         = _fjsp_mul_v2r8(rinv22,rinv22);
1127             rinvsq23         = _fjsp_mul_v2r8(rinv23,rinv23);
1128             rinvsq31         = _fjsp_mul_v2r8(rinv31,rinv31);
1129             rinvsq32         = _fjsp_mul_v2r8(rinv32,rinv32);
1130             rinvsq33         = _fjsp_mul_v2r8(rinv33,rinv33);
1131
1132             fjx0             = _fjsp_setzero_v2r8();
1133             fjy0             = _fjsp_setzero_v2r8();
1134             fjz0             = _fjsp_setzero_v2r8();
1135             fjx1             = _fjsp_setzero_v2r8();
1136             fjy1             = _fjsp_setzero_v2r8();
1137             fjz1             = _fjsp_setzero_v2r8();
1138             fjx2             = _fjsp_setzero_v2r8();
1139             fjy2             = _fjsp_setzero_v2r8();
1140             fjz2             = _fjsp_setzero_v2r8();
1141             fjx3             = _fjsp_setzero_v2r8();
1142             fjy3             = _fjsp_setzero_v2r8();
1143             fjz3             = _fjsp_setzero_v2r8();
1144
1145             /**************************
1146              * CALCULATE INTERACTIONS *
1147              **************************/
1148
1149             /* LENNARD-JONES DISPERSION/REPULSION */
1150
1151             rinvsix          = _fjsp_mul_v2r8(_fjsp_mul_v2r8(rinvsq00,rinvsq00),rinvsq00);
1152             fvdw             = _fjsp_mul_v2r8(_fjsp_msub_v2r8(c12_00,rinvsix,c6_00),_fjsp_mul_v2r8(rinvsix,rinvsq00));
1153
1154             fscal            = fvdw;
1155
1156             /* Update vectorial force */
1157             fix0             = _fjsp_madd_v2r8(dx00,fscal,fix0);
1158             fiy0             = _fjsp_madd_v2r8(dy00,fscal,fiy0);
1159             fiz0             = _fjsp_madd_v2r8(dz00,fscal,fiz0);
1160             
1161             fjx0             = _fjsp_madd_v2r8(dx00,fscal,fjx0);
1162             fjy0             = _fjsp_madd_v2r8(dy00,fscal,fjy0);
1163             fjz0             = _fjsp_madd_v2r8(dz00,fscal,fjz0);
1164
1165             /**************************
1166              * CALCULATE INTERACTIONS *
1167              **************************/
1168
1169             /* REACTION-FIELD ELECTROSTATICS */
1170             felec            = _fjsp_mul_v2r8(qq11,_fjsp_msub_v2r8(rinv11,rinvsq11,krf2));
1171
1172             fscal            = felec;
1173
1174             /* Update vectorial force */
1175             fix1             = _fjsp_madd_v2r8(dx11,fscal,fix1);
1176             fiy1             = _fjsp_madd_v2r8(dy11,fscal,fiy1);
1177             fiz1             = _fjsp_madd_v2r8(dz11,fscal,fiz1);
1178             
1179             fjx1             = _fjsp_madd_v2r8(dx11,fscal,fjx1);
1180             fjy1             = _fjsp_madd_v2r8(dy11,fscal,fjy1);
1181             fjz1             = _fjsp_madd_v2r8(dz11,fscal,fjz1);
1182
1183             /**************************
1184              * CALCULATE INTERACTIONS *
1185              **************************/
1186
1187             /* REACTION-FIELD ELECTROSTATICS */
1188             felec            = _fjsp_mul_v2r8(qq12,_fjsp_msub_v2r8(rinv12,rinvsq12,krf2));
1189
1190             fscal            = felec;
1191
1192             /* Update vectorial force */
1193             fix1             = _fjsp_madd_v2r8(dx12,fscal,fix1);
1194             fiy1             = _fjsp_madd_v2r8(dy12,fscal,fiy1);
1195             fiz1             = _fjsp_madd_v2r8(dz12,fscal,fiz1);
1196             
1197             fjx2             = _fjsp_madd_v2r8(dx12,fscal,fjx2);
1198             fjy2             = _fjsp_madd_v2r8(dy12,fscal,fjy2);
1199             fjz2             = _fjsp_madd_v2r8(dz12,fscal,fjz2);
1200
1201             /**************************
1202              * CALCULATE INTERACTIONS *
1203              **************************/
1204
1205             /* REACTION-FIELD ELECTROSTATICS */
1206             felec            = _fjsp_mul_v2r8(qq13,_fjsp_msub_v2r8(rinv13,rinvsq13,krf2));
1207
1208             fscal            = felec;
1209
1210             /* Update vectorial force */
1211             fix1             = _fjsp_madd_v2r8(dx13,fscal,fix1);
1212             fiy1             = _fjsp_madd_v2r8(dy13,fscal,fiy1);
1213             fiz1             = _fjsp_madd_v2r8(dz13,fscal,fiz1);
1214             
1215             fjx3             = _fjsp_madd_v2r8(dx13,fscal,fjx3);
1216             fjy3             = _fjsp_madd_v2r8(dy13,fscal,fjy3);
1217             fjz3             = _fjsp_madd_v2r8(dz13,fscal,fjz3);
1218
1219             /**************************
1220              * CALCULATE INTERACTIONS *
1221              **************************/
1222
1223             /* REACTION-FIELD ELECTROSTATICS */
1224             felec            = _fjsp_mul_v2r8(qq21,_fjsp_msub_v2r8(rinv21,rinvsq21,krf2));
1225
1226             fscal            = felec;
1227
1228             /* Update vectorial force */
1229             fix2             = _fjsp_madd_v2r8(dx21,fscal,fix2);
1230             fiy2             = _fjsp_madd_v2r8(dy21,fscal,fiy2);
1231             fiz2             = _fjsp_madd_v2r8(dz21,fscal,fiz2);
1232             
1233             fjx1             = _fjsp_madd_v2r8(dx21,fscal,fjx1);
1234             fjy1             = _fjsp_madd_v2r8(dy21,fscal,fjy1);
1235             fjz1             = _fjsp_madd_v2r8(dz21,fscal,fjz1);
1236
1237             /**************************
1238              * CALCULATE INTERACTIONS *
1239              **************************/
1240
1241             /* REACTION-FIELD ELECTROSTATICS */
1242             felec            = _fjsp_mul_v2r8(qq22,_fjsp_msub_v2r8(rinv22,rinvsq22,krf2));
1243
1244             fscal            = felec;
1245
1246             /* Update vectorial force */
1247             fix2             = _fjsp_madd_v2r8(dx22,fscal,fix2);
1248             fiy2             = _fjsp_madd_v2r8(dy22,fscal,fiy2);
1249             fiz2             = _fjsp_madd_v2r8(dz22,fscal,fiz2);
1250             
1251             fjx2             = _fjsp_madd_v2r8(dx22,fscal,fjx2);
1252             fjy2             = _fjsp_madd_v2r8(dy22,fscal,fjy2);
1253             fjz2             = _fjsp_madd_v2r8(dz22,fscal,fjz2);
1254
1255             /**************************
1256              * CALCULATE INTERACTIONS *
1257              **************************/
1258
1259             /* REACTION-FIELD ELECTROSTATICS */
1260             felec            = _fjsp_mul_v2r8(qq23,_fjsp_msub_v2r8(rinv23,rinvsq23,krf2));
1261
1262             fscal            = felec;
1263
1264             /* Update vectorial force */
1265             fix2             = _fjsp_madd_v2r8(dx23,fscal,fix2);
1266             fiy2             = _fjsp_madd_v2r8(dy23,fscal,fiy2);
1267             fiz2             = _fjsp_madd_v2r8(dz23,fscal,fiz2);
1268             
1269             fjx3             = _fjsp_madd_v2r8(dx23,fscal,fjx3);
1270             fjy3             = _fjsp_madd_v2r8(dy23,fscal,fjy3);
1271             fjz3             = _fjsp_madd_v2r8(dz23,fscal,fjz3);
1272
1273             /**************************
1274              * CALCULATE INTERACTIONS *
1275              **************************/
1276
1277             /* REACTION-FIELD ELECTROSTATICS */
1278             felec            = _fjsp_mul_v2r8(qq31,_fjsp_msub_v2r8(rinv31,rinvsq31,krf2));
1279
1280             fscal            = felec;
1281
1282             /* Update vectorial force */
1283             fix3             = _fjsp_madd_v2r8(dx31,fscal,fix3);
1284             fiy3             = _fjsp_madd_v2r8(dy31,fscal,fiy3);
1285             fiz3             = _fjsp_madd_v2r8(dz31,fscal,fiz3);
1286             
1287             fjx1             = _fjsp_madd_v2r8(dx31,fscal,fjx1);
1288             fjy1             = _fjsp_madd_v2r8(dy31,fscal,fjy1);
1289             fjz1             = _fjsp_madd_v2r8(dz31,fscal,fjz1);
1290
1291             /**************************
1292              * CALCULATE INTERACTIONS *
1293              **************************/
1294
1295             /* REACTION-FIELD ELECTROSTATICS */
1296             felec            = _fjsp_mul_v2r8(qq32,_fjsp_msub_v2r8(rinv32,rinvsq32,krf2));
1297
1298             fscal            = felec;
1299
1300             /* Update vectorial force */
1301             fix3             = _fjsp_madd_v2r8(dx32,fscal,fix3);
1302             fiy3             = _fjsp_madd_v2r8(dy32,fscal,fiy3);
1303             fiz3             = _fjsp_madd_v2r8(dz32,fscal,fiz3);
1304             
1305             fjx2             = _fjsp_madd_v2r8(dx32,fscal,fjx2);
1306             fjy2             = _fjsp_madd_v2r8(dy32,fscal,fjy2);
1307             fjz2             = _fjsp_madd_v2r8(dz32,fscal,fjz2);
1308
1309             /**************************
1310              * CALCULATE INTERACTIONS *
1311              **************************/
1312
1313             /* REACTION-FIELD ELECTROSTATICS */
1314             felec            = _fjsp_mul_v2r8(qq33,_fjsp_msub_v2r8(rinv33,rinvsq33,krf2));
1315
1316             fscal            = felec;
1317
1318             /* Update vectorial force */
1319             fix3             = _fjsp_madd_v2r8(dx33,fscal,fix3);
1320             fiy3             = _fjsp_madd_v2r8(dy33,fscal,fiy3);
1321             fiz3             = _fjsp_madd_v2r8(dz33,fscal,fiz3);
1322             
1323             fjx3             = _fjsp_madd_v2r8(dx33,fscal,fjx3);
1324             fjy3             = _fjsp_madd_v2r8(dy33,fscal,fjy3);
1325             fjz3             = _fjsp_madd_v2r8(dz33,fscal,fjz3);
1326
1327             gmx_fjsp_decrement_4rvec_2ptr_swizzle_v2r8(f+j_coord_offsetA,f+j_coord_offsetB,fjx0,fjy0,fjz0,fjx1,fjy1,fjz1,fjx2,fjy2,fjz2,fjx3,fjy3,fjz3);
1328
1329             /* Inner loop uses 303 flops */
1330         }
1331
1332         if(jidx<j_index_end)
1333         {
1334
1335             jnrA             = jjnr[jidx];
1336             j_coord_offsetA  = DIM*jnrA;
1337
1338             /* load j atom coordinates */
1339             gmx_fjsp_load_4rvec_1ptr_swizzle_v2r8(x+j_coord_offsetA,
1340                                               &jx0,&jy0,&jz0,&jx1,&jy1,&jz1,&jx2,
1341                                               &jy2,&jz2,&jx3,&jy3,&jz3);
1342
1343             /* Calculate displacement vector */
1344             dx00             = _fjsp_sub_v2r8(ix0,jx0);
1345             dy00             = _fjsp_sub_v2r8(iy0,jy0);
1346             dz00             = _fjsp_sub_v2r8(iz0,jz0);
1347             dx11             = _fjsp_sub_v2r8(ix1,jx1);
1348             dy11             = _fjsp_sub_v2r8(iy1,jy1);
1349             dz11             = _fjsp_sub_v2r8(iz1,jz1);
1350             dx12             = _fjsp_sub_v2r8(ix1,jx2);
1351             dy12             = _fjsp_sub_v2r8(iy1,jy2);
1352             dz12             = _fjsp_sub_v2r8(iz1,jz2);
1353             dx13             = _fjsp_sub_v2r8(ix1,jx3);
1354             dy13             = _fjsp_sub_v2r8(iy1,jy3);
1355             dz13             = _fjsp_sub_v2r8(iz1,jz3);
1356             dx21             = _fjsp_sub_v2r8(ix2,jx1);
1357             dy21             = _fjsp_sub_v2r8(iy2,jy1);
1358             dz21             = _fjsp_sub_v2r8(iz2,jz1);
1359             dx22             = _fjsp_sub_v2r8(ix2,jx2);
1360             dy22             = _fjsp_sub_v2r8(iy2,jy2);
1361             dz22             = _fjsp_sub_v2r8(iz2,jz2);
1362             dx23             = _fjsp_sub_v2r8(ix2,jx3);
1363             dy23             = _fjsp_sub_v2r8(iy2,jy3);
1364             dz23             = _fjsp_sub_v2r8(iz2,jz3);
1365             dx31             = _fjsp_sub_v2r8(ix3,jx1);
1366             dy31             = _fjsp_sub_v2r8(iy3,jy1);
1367             dz31             = _fjsp_sub_v2r8(iz3,jz1);
1368             dx32             = _fjsp_sub_v2r8(ix3,jx2);
1369             dy32             = _fjsp_sub_v2r8(iy3,jy2);
1370             dz32             = _fjsp_sub_v2r8(iz3,jz2);
1371             dx33             = _fjsp_sub_v2r8(ix3,jx3);
1372             dy33             = _fjsp_sub_v2r8(iy3,jy3);
1373             dz33             = _fjsp_sub_v2r8(iz3,jz3);
1374
1375             /* Calculate squared distance and things based on it */
1376             rsq00            = gmx_fjsp_calc_rsq_v2r8(dx00,dy00,dz00);
1377             rsq11            = gmx_fjsp_calc_rsq_v2r8(dx11,dy11,dz11);
1378             rsq12            = gmx_fjsp_calc_rsq_v2r8(dx12,dy12,dz12);
1379             rsq13            = gmx_fjsp_calc_rsq_v2r8(dx13,dy13,dz13);
1380             rsq21            = gmx_fjsp_calc_rsq_v2r8(dx21,dy21,dz21);
1381             rsq22            = gmx_fjsp_calc_rsq_v2r8(dx22,dy22,dz22);
1382             rsq23            = gmx_fjsp_calc_rsq_v2r8(dx23,dy23,dz23);
1383             rsq31            = gmx_fjsp_calc_rsq_v2r8(dx31,dy31,dz31);
1384             rsq32            = gmx_fjsp_calc_rsq_v2r8(dx32,dy32,dz32);
1385             rsq33            = gmx_fjsp_calc_rsq_v2r8(dx33,dy33,dz33);
1386
1387             rinv11           = gmx_fjsp_invsqrt_v2r8(rsq11);
1388             rinv12           = gmx_fjsp_invsqrt_v2r8(rsq12);
1389             rinv13           = gmx_fjsp_invsqrt_v2r8(rsq13);
1390             rinv21           = gmx_fjsp_invsqrt_v2r8(rsq21);
1391             rinv22           = gmx_fjsp_invsqrt_v2r8(rsq22);
1392             rinv23           = gmx_fjsp_invsqrt_v2r8(rsq23);
1393             rinv31           = gmx_fjsp_invsqrt_v2r8(rsq31);
1394             rinv32           = gmx_fjsp_invsqrt_v2r8(rsq32);
1395             rinv33           = gmx_fjsp_invsqrt_v2r8(rsq33);
1396
1397             rinvsq00         = gmx_fjsp_inv_v2r8(rsq00);
1398             rinvsq11         = _fjsp_mul_v2r8(rinv11,rinv11);
1399             rinvsq12         = _fjsp_mul_v2r8(rinv12,rinv12);
1400             rinvsq13         = _fjsp_mul_v2r8(rinv13,rinv13);
1401             rinvsq21         = _fjsp_mul_v2r8(rinv21,rinv21);
1402             rinvsq22         = _fjsp_mul_v2r8(rinv22,rinv22);
1403             rinvsq23         = _fjsp_mul_v2r8(rinv23,rinv23);
1404             rinvsq31         = _fjsp_mul_v2r8(rinv31,rinv31);
1405             rinvsq32         = _fjsp_mul_v2r8(rinv32,rinv32);
1406             rinvsq33         = _fjsp_mul_v2r8(rinv33,rinv33);
1407
1408             fjx0             = _fjsp_setzero_v2r8();
1409             fjy0             = _fjsp_setzero_v2r8();
1410             fjz0             = _fjsp_setzero_v2r8();
1411             fjx1             = _fjsp_setzero_v2r8();
1412             fjy1             = _fjsp_setzero_v2r8();
1413             fjz1             = _fjsp_setzero_v2r8();
1414             fjx2             = _fjsp_setzero_v2r8();
1415             fjy2             = _fjsp_setzero_v2r8();
1416             fjz2             = _fjsp_setzero_v2r8();
1417             fjx3             = _fjsp_setzero_v2r8();
1418             fjy3             = _fjsp_setzero_v2r8();
1419             fjz3             = _fjsp_setzero_v2r8();
1420
1421             /**************************
1422              * CALCULATE INTERACTIONS *
1423              **************************/
1424
1425             /* LENNARD-JONES DISPERSION/REPULSION */
1426
1427             rinvsix          = _fjsp_mul_v2r8(_fjsp_mul_v2r8(rinvsq00,rinvsq00),rinvsq00);
1428             fvdw             = _fjsp_mul_v2r8(_fjsp_msub_v2r8(c12_00,rinvsix,c6_00),_fjsp_mul_v2r8(rinvsix,rinvsq00));
1429
1430             fscal            = fvdw;
1431
1432             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
1433
1434             /* Update vectorial force */
1435             fix0             = _fjsp_madd_v2r8(dx00,fscal,fix0);
1436             fiy0             = _fjsp_madd_v2r8(dy00,fscal,fiy0);
1437             fiz0             = _fjsp_madd_v2r8(dz00,fscal,fiz0);
1438             
1439             fjx0             = _fjsp_madd_v2r8(dx00,fscal,fjx0);
1440             fjy0             = _fjsp_madd_v2r8(dy00,fscal,fjy0);
1441             fjz0             = _fjsp_madd_v2r8(dz00,fscal,fjz0);
1442
1443             /**************************
1444              * CALCULATE INTERACTIONS *
1445              **************************/
1446
1447             /* REACTION-FIELD ELECTROSTATICS */
1448             felec            = _fjsp_mul_v2r8(qq11,_fjsp_msub_v2r8(rinv11,rinvsq11,krf2));
1449
1450             fscal            = felec;
1451
1452             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
1453
1454             /* Update vectorial force */
1455             fix1             = _fjsp_madd_v2r8(dx11,fscal,fix1);
1456             fiy1             = _fjsp_madd_v2r8(dy11,fscal,fiy1);
1457             fiz1             = _fjsp_madd_v2r8(dz11,fscal,fiz1);
1458             
1459             fjx1             = _fjsp_madd_v2r8(dx11,fscal,fjx1);
1460             fjy1             = _fjsp_madd_v2r8(dy11,fscal,fjy1);
1461             fjz1             = _fjsp_madd_v2r8(dz11,fscal,fjz1);
1462
1463             /**************************
1464              * CALCULATE INTERACTIONS *
1465              **************************/
1466
1467             /* REACTION-FIELD ELECTROSTATICS */
1468             felec            = _fjsp_mul_v2r8(qq12,_fjsp_msub_v2r8(rinv12,rinvsq12,krf2));
1469
1470             fscal            = felec;
1471
1472             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
1473
1474             /* Update vectorial force */
1475             fix1             = _fjsp_madd_v2r8(dx12,fscal,fix1);
1476             fiy1             = _fjsp_madd_v2r8(dy12,fscal,fiy1);
1477             fiz1             = _fjsp_madd_v2r8(dz12,fscal,fiz1);
1478             
1479             fjx2             = _fjsp_madd_v2r8(dx12,fscal,fjx2);
1480             fjy2             = _fjsp_madd_v2r8(dy12,fscal,fjy2);
1481             fjz2             = _fjsp_madd_v2r8(dz12,fscal,fjz2);
1482
1483             /**************************
1484              * CALCULATE INTERACTIONS *
1485              **************************/
1486
1487             /* REACTION-FIELD ELECTROSTATICS */
1488             felec            = _fjsp_mul_v2r8(qq13,_fjsp_msub_v2r8(rinv13,rinvsq13,krf2));
1489
1490             fscal            = felec;
1491
1492             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
1493
1494             /* Update vectorial force */
1495             fix1             = _fjsp_madd_v2r8(dx13,fscal,fix1);
1496             fiy1             = _fjsp_madd_v2r8(dy13,fscal,fiy1);
1497             fiz1             = _fjsp_madd_v2r8(dz13,fscal,fiz1);
1498             
1499             fjx3             = _fjsp_madd_v2r8(dx13,fscal,fjx3);
1500             fjy3             = _fjsp_madd_v2r8(dy13,fscal,fjy3);
1501             fjz3             = _fjsp_madd_v2r8(dz13,fscal,fjz3);
1502
1503             /**************************
1504              * CALCULATE INTERACTIONS *
1505              **************************/
1506
1507             /* REACTION-FIELD ELECTROSTATICS */
1508             felec            = _fjsp_mul_v2r8(qq21,_fjsp_msub_v2r8(rinv21,rinvsq21,krf2));
1509
1510             fscal            = felec;
1511
1512             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
1513
1514             /* Update vectorial force */
1515             fix2             = _fjsp_madd_v2r8(dx21,fscal,fix2);
1516             fiy2             = _fjsp_madd_v2r8(dy21,fscal,fiy2);
1517             fiz2             = _fjsp_madd_v2r8(dz21,fscal,fiz2);
1518             
1519             fjx1             = _fjsp_madd_v2r8(dx21,fscal,fjx1);
1520             fjy1             = _fjsp_madd_v2r8(dy21,fscal,fjy1);
1521             fjz1             = _fjsp_madd_v2r8(dz21,fscal,fjz1);
1522
1523             /**************************
1524              * CALCULATE INTERACTIONS *
1525              **************************/
1526
1527             /* REACTION-FIELD ELECTROSTATICS */
1528             felec            = _fjsp_mul_v2r8(qq22,_fjsp_msub_v2r8(rinv22,rinvsq22,krf2));
1529
1530             fscal            = felec;
1531
1532             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
1533
1534             /* Update vectorial force */
1535             fix2             = _fjsp_madd_v2r8(dx22,fscal,fix2);
1536             fiy2             = _fjsp_madd_v2r8(dy22,fscal,fiy2);
1537             fiz2             = _fjsp_madd_v2r8(dz22,fscal,fiz2);
1538             
1539             fjx2             = _fjsp_madd_v2r8(dx22,fscal,fjx2);
1540             fjy2             = _fjsp_madd_v2r8(dy22,fscal,fjy2);
1541             fjz2             = _fjsp_madd_v2r8(dz22,fscal,fjz2);
1542
1543             /**************************
1544              * CALCULATE INTERACTIONS *
1545              **************************/
1546
1547             /* REACTION-FIELD ELECTROSTATICS */
1548             felec            = _fjsp_mul_v2r8(qq23,_fjsp_msub_v2r8(rinv23,rinvsq23,krf2));
1549
1550             fscal            = felec;
1551
1552             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
1553
1554             /* Update vectorial force */
1555             fix2             = _fjsp_madd_v2r8(dx23,fscal,fix2);
1556             fiy2             = _fjsp_madd_v2r8(dy23,fscal,fiy2);
1557             fiz2             = _fjsp_madd_v2r8(dz23,fscal,fiz2);
1558             
1559             fjx3             = _fjsp_madd_v2r8(dx23,fscal,fjx3);
1560             fjy3             = _fjsp_madd_v2r8(dy23,fscal,fjy3);
1561             fjz3             = _fjsp_madd_v2r8(dz23,fscal,fjz3);
1562
1563             /**************************
1564              * CALCULATE INTERACTIONS *
1565              **************************/
1566
1567             /* REACTION-FIELD ELECTROSTATICS */
1568             felec            = _fjsp_mul_v2r8(qq31,_fjsp_msub_v2r8(rinv31,rinvsq31,krf2));
1569
1570             fscal            = felec;
1571
1572             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
1573
1574             /* Update vectorial force */
1575             fix3             = _fjsp_madd_v2r8(dx31,fscal,fix3);
1576             fiy3             = _fjsp_madd_v2r8(dy31,fscal,fiy3);
1577             fiz3             = _fjsp_madd_v2r8(dz31,fscal,fiz3);
1578             
1579             fjx1             = _fjsp_madd_v2r8(dx31,fscal,fjx1);
1580             fjy1             = _fjsp_madd_v2r8(dy31,fscal,fjy1);
1581             fjz1             = _fjsp_madd_v2r8(dz31,fscal,fjz1);
1582
1583             /**************************
1584              * CALCULATE INTERACTIONS *
1585              **************************/
1586
1587             /* REACTION-FIELD ELECTROSTATICS */
1588             felec            = _fjsp_mul_v2r8(qq32,_fjsp_msub_v2r8(rinv32,rinvsq32,krf2));
1589
1590             fscal            = felec;
1591
1592             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
1593
1594             /* Update vectorial force */
1595             fix3             = _fjsp_madd_v2r8(dx32,fscal,fix3);
1596             fiy3             = _fjsp_madd_v2r8(dy32,fscal,fiy3);
1597             fiz3             = _fjsp_madd_v2r8(dz32,fscal,fiz3);
1598             
1599             fjx2             = _fjsp_madd_v2r8(dx32,fscal,fjx2);
1600             fjy2             = _fjsp_madd_v2r8(dy32,fscal,fjy2);
1601             fjz2             = _fjsp_madd_v2r8(dz32,fscal,fjz2);
1602
1603             /**************************
1604              * CALCULATE INTERACTIONS *
1605              **************************/
1606
1607             /* REACTION-FIELD ELECTROSTATICS */
1608             felec            = _fjsp_mul_v2r8(qq33,_fjsp_msub_v2r8(rinv33,rinvsq33,krf2));
1609
1610             fscal            = felec;
1611
1612             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
1613
1614             /* Update vectorial force */
1615             fix3             = _fjsp_madd_v2r8(dx33,fscal,fix3);
1616             fiy3             = _fjsp_madd_v2r8(dy33,fscal,fiy3);
1617             fiz3             = _fjsp_madd_v2r8(dz33,fscal,fiz3);
1618             
1619             fjx3             = _fjsp_madd_v2r8(dx33,fscal,fjx3);
1620             fjy3             = _fjsp_madd_v2r8(dy33,fscal,fjy3);
1621             fjz3             = _fjsp_madd_v2r8(dz33,fscal,fjz3);
1622
1623             gmx_fjsp_decrement_4rvec_1ptr_swizzle_v2r8(f+j_coord_offsetA,fjx0,fjy0,fjz0,fjx1,fjy1,fjz1,fjx2,fjy2,fjz2,fjx3,fjy3,fjz3);
1624
1625             /* Inner loop uses 303 flops */
1626         }
1627
1628         /* End of innermost loop */
1629
1630         gmx_fjsp_update_iforce_4atom_swizzle_v2r8(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
1631                                               f+i_coord_offset,fshift+i_shift_offset);
1632
1633         /* Increment number of inner iterations */
1634         inneriter                  += j_index_end - j_index_start;
1635
1636         /* Outer loop uses 24 flops */
1637     }
1638
1639     /* Increment number of outer iterations */
1640     outeriter        += nri;
1641
1642     /* Update outer/inner flops */
1643
1644     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W4W4_F,outeriter*24 + inneriter*303);
1645 }