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