Move vec.h to math/
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_sparc64_hpc_ace_double / nb_kernel_ElecRF_VdwCSTab_GeomW4P1_sparc64_hpc_ace_double.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2013,2014, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35 /*
36  * Note: this file was generated by the GROMACS sparc64_hpc_ace_double kernel generator.
37  */
38 #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 "gromacs/math/vec.h"
47 #include "nrnb.h"
48
49 #include "kernelutil_sparc64_hpc_ace_double.h"
50
51 /*
52  * Gromacs nonbonded kernel:   nb_kernel_ElecRF_VdwCSTab_GeomW4P1_VF_sparc64_hpc_ace_double
53  * Electrostatics interaction: ReactionField
54  * VdW interaction:            CubicSplineTable
55  * Geometry:                   Water4-Particle
56  * Calculate force/pot:        PotentialAndForce
57  */
58 void
59 nb_kernel_ElecRF_VdwCSTab_GeomW4P1_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_unused * 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     _fjsp_v2r8       dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
92     _fjsp_v2r8       dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
93     _fjsp_v2r8       dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
94     _fjsp_v2r8       dx30,dy30,dz30,rsq30,rinv30,rinvsq30,r30,qq30,c6_30,c12_30;
95     _fjsp_v2r8       velec,felec,velecsum,facel,crf,krf,krf2;
96     real             *charge;
97     int              nvdwtype;
98     _fjsp_v2r8       rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
99     int              *vdwtype;
100     real             *vdwparam;
101     _fjsp_v2r8       one_sixth   = gmx_fjsp_set1_v2r8(1.0/6.0);
102     _fjsp_v2r8       one_twelfth = gmx_fjsp_set1_v2r8(1.0/12.0);
103     _fjsp_v2r8       rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF,twovfeps;
104     real             *vftab;
105     _fjsp_v2r8       itab_tmp;
106     _fjsp_v2r8       dummy_mask,cutoff_mask;
107     _fjsp_v2r8       one     = gmx_fjsp_set1_v2r8(1.0);
108     _fjsp_v2r8       two     = gmx_fjsp_set1_v2r8(2.0);
109     union { _fjsp_v2r8 simd; long long int i[2]; } vfconv,gbconv,ewconv;
110
111     x                = xx[0];
112     f                = ff[0];
113
114     nri              = nlist->nri;
115     iinr             = nlist->iinr;
116     jindex           = nlist->jindex;
117     jjnr             = nlist->jjnr;
118     shiftidx         = nlist->shift;
119     gid              = nlist->gid;
120     shiftvec         = fr->shift_vec[0];
121     fshift           = fr->fshift[0];
122     facel            = gmx_fjsp_set1_v2r8(fr->epsfac);
123     charge           = mdatoms->chargeA;
124     krf              = gmx_fjsp_set1_v2r8(fr->ic->k_rf);
125     krf2             = gmx_fjsp_set1_v2r8(fr->ic->k_rf*2.0);
126     crf              = gmx_fjsp_set1_v2r8(fr->ic->c_rf);
127     nvdwtype         = fr->ntype;
128     vdwparam         = fr->nbfp;
129     vdwtype          = mdatoms->typeA;
130
131     vftab            = kernel_data->table_vdw->data;
132     vftabscale       = gmx_fjsp_set1_v2r8(kernel_data->table_vdw->scale);
133
134     /* Setup water-specific parameters */
135     inr              = nlist->iinr[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     iq3              = _fjsp_mul_v2r8(facel,gmx_fjsp_set1_v2r8(charge[inr+3]));
139     vdwioffset0      = 2*nvdwtype*vdwtype[inr+0];
140
141     /* Avoid stupid compiler warnings */
142     jnrA = jnrB = 0;
143     j_coord_offsetA = 0;
144     j_coord_offsetB = 0;
145
146     outeriter        = 0;
147     inneriter        = 0;
148
149     /* Start outer loop over neighborlists */
150     for(iidx=0; iidx<nri; iidx++)
151     {
152         /* Load shift vector for this list */
153         i_shift_offset   = DIM*shiftidx[iidx];
154
155         /* Load limits for loop over neighbors */
156         j_index_start    = jindex[iidx];
157         j_index_end      = jindex[iidx+1];
158
159         /* Get outer coordinate index */
160         inr              = iinr[iidx];
161         i_coord_offset   = DIM*inr;
162
163         /* Load i particle coords and add shift vector */
164         gmx_fjsp_load_shift_and_4rvec_broadcast_v2r8(shiftvec+i_shift_offset,x+i_coord_offset,
165                                                  &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
166
167         fix0             = _fjsp_setzero_v2r8();
168         fiy0             = _fjsp_setzero_v2r8();
169         fiz0             = _fjsp_setzero_v2r8();
170         fix1             = _fjsp_setzero_v2r8();
171         fiy1             = _fjsp_setzero_v2r8();
172         fiz1             = _fjsp_setzero_v2r8();
173         fix2             = _fjsp_setzero_v2r8();
174         fiy2             = _fjsp_setzero_v2r8();
175         fiz2             = _fjsp_setzero_v2r8();
176         fix3             = _fjsp_setzero_v2r8();
177         fiy3             = _fjsp_setzero_v2r8();
178         fiz3             = _fjsp_setzero_v2r8();
179
180         /* Reset potential sums */
181         velecsum         = _fjsp_setzero_v2r8();
182         vvdwsum          = _fjsp_setzero_v2r8();
183
184         /* Start inner kernel loop */
185         for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
186         {
187
188             /* Get j neighbor index, and coordinate index */
189             jnrA             = jjnr[jidx];
190             jnrB             = jjnr[jidx+1];
191             j_coord_offsetA  = DIM*jnrA;
192             j_coord_offsetB  = DIM*jnrB;
193
194             /* load j atom coordinates */
195             gmx_fjsp_load_1rvec_2ptr_swizzle_v2r8(x+j_coord_offsetA,x+j_coord_offsetB,
196                                               &jx0,&jy0,&jz0);
197
198             /* Calculate displacement vector */
199             dx00             = _fjsp_sub_v2r8(ix0,jx0);
200             dy00             = _fjsp_sub_v2r8(iy0,jy0);
201             dz00             = _fjsp_sub_v2r8(iz0,jz0);
202             dx10             = _fjsp_sub_v2r8(ix1,jx0);
203             dy10             = _fjsp_sub_v2r8(iy1,jy0);
204             dz10             = _fjsp_sub_v2r8(iz1,jz0);
205             dx20             = _fjsp_sub_v2r8(ix2,jx0);
206             dy20             = _fjsp_sub_v2r8(iy2,jy0);
207             dz20             = _fjsp_sub_v2r8(iz2,jz0);
208             dx30             = _fjsp_sub_v2r8(ix3,jx0);
209             dy30             = _fjsp_sub_v2r8(iy3,jy0);
210             dz30             = _fjsp_sub_v2r8(iz3,jz0);
211
212             /* Calculate squared distance and things based on it */
213             rsq00            = gmx_fjsp_calc_rsq_v2r8(dx00,dy00,dz00);
214             rsq10            = gmx_fjsp_calc_rsq_v2r8(dx10,dy10,dz10);
215             rsq20            = gmx_fjsp_calc_rsq_v2r8(dx20,dy20,dz20);
216             rsq30            = gmx_fjsp_calc_rsq_v2r8(dx30,dy30,dz30);
217
218             rinv00           = gmx_fjsp_invsqrt_v2r8(rsq00);
219             rinv10           = gmx_fjsp_invsqrt_v2r8(rsq10);
220             rinv20           = gmx_fjsp_invsqrt_v2r8(rsq20);
221             rinv30           = gmx_fjsp_invsqrt_v2r8(rsq30);
222
223             rinvsq10         = _fjsp_mul_v2r8(rinv10,rinv10);
224             rinvsq20         = _fjsp_mul_v2r8(rinv20,rinv20);
225             rinvsq30         = _fjsp_mul_v2r8(rinv30,rinv30);
226
227             /* Load parameters for j particles */
228             jq0              = gmx_fjsp_load_2real_swizzle_v2r8(charge+jnrA+0,charge+jnrB+0);
229             vdwjidx0A        = 2*vdwtype[jnrA+0];
230             vdwjidx0B        = 2*vdwtype[jnrB+0];
231
232             fjx0             = _fjsp_setzero_v2r8();
233             fjy0             = _fjsp_setzero_v2r8();
234             fjz0             = _fjsp_setzero_v2r8();
235
236             /**************************
237              * CALCULATE INTERACTIONS *
238              **************************/
239
240             r00              = _fjsp_mul_v2r8(rsq00,rinv00);
241
242             /* Compute parameters for interactions between i and j atoms */
243             gmx_fjsp_load_2pair_swizzle_v2r8(vdwparam+vdwioffset0+vdwjidx0A,
244                                          vdwparam+vdwioffset0+vdwjidx0B,&c6_00,&c12_00);
245
246             /* Calculate table index by multiplying r with table scale and truncate to integer */
247             rt               = _fjsp_mul_v2r8(r00,vftabscale);
248             itab_tmp         = _fjsp_dtox_v2r8(rt);
249             vfeps            = _fjsp_sub_v2r8(rt, _fjsp_xtod_v2r8(itab_tmp));
250             twovfeps         = _fjsp_add_v2r8(vfeps,vfeps);
251             _fjsp_store_v2r8(&vfconv.simd,itab_tmp);
252
253             vfconv.i[0]     *= 8;
254             vfconv.i[1]     *= 8;
255
256             /* CUBIC SPLINE TABLE DISPERSION */
257             Y                = _fjsp_load_v2r8( vftab + vfconv.i[0] );
258             F                = _fjsp_load_v2r8( vftab + vfconv.i[1] );
259             GMX_FJSP_TRANSPOSE2_V2R8(Y,F);
260             G                = _fjsp_load_v2r8( vftab + vfconv.i[0] + 2 );
261             H                = _fjsp_load_v2r8( vftab + vfconv.i[1] + 2 );
262             GMX_FJSP_TRANSPOSE2_V2R8(G,H);
263             Fp               = _fjsp_madd_v2r8(vfeps,_fjsp_madd_v2r8(H,vfeps,G),F);
264             VV               = _fjsp_madd_v2r8(vfeps,Fp,Y);
265             vvdw6            = _fjsp_mul_v2r8(c6_00,VV);
266             FF               = _fjsp_madd_v2r8(vfeps,_fjsp_madd_v2r8(twovfeps,H,G),Fp);
267             fvdw6            = _fjsp_mul_v2r8(c6_00,FF);
268
269             /* CUBIC SPLINE TABLE REPULSION */
270             Y                = _fjsp_load_v2r8( vftab + vfconv.i[0] + 4 );
271             F                = _fjsp_load_v2r8( vftab + vfconv.i[1] + 4 );
272             GMX_FJSP_TRANSPOSE2_V2R8(Y,F);
273             G                = _fjsp_load_v2r8( vftab + vfconv.i[0] + 6 );
274             H                = _fjsp_load_v2r8( vftab + vfconv.i[1] + 6 );
275             GMX_FJSP_TRANSPOSE2_V2R8(G,H);
276             Fp               = _fjsp_madd_v2r8(vfeps,_fjsp_madd_v2r8(H,vfeps,G),F);
277             VV               = _fjsp_madd_v2r8(vfeps,Fp,Y);
278             vvdw12           = _fjsp_mul_v2r8(c12_00,VV);
279             FF               = _fjsp_madd_v2r8(vfeps,_fjsp_madd_v2r8(twovfeps,H,G),Fp);
280             fvdw12           = _fjsp_mul_v2r8(c12_00,FF);
281             vvdw             = _fjsp_add_v2r8(vvdw12,vvdw6);
282             fvdw             = _fjsp_neg_v2r8(_fjsp_mul_v2r8(_fjsp_add_v2r8(fvdw6,fvdw12),_fjsp_mul_v2r8(vftabscale,rinv00)));
283
284             /* Update potential sum for this i atom from the interaction with this j atom. */
285             vvdwsum          = _fjsp_add_v2r8(vvdwsum,vvdw);
286
287             fscal            = fvdw;
288
289             /* Update vectorial force */
290             fix0             = _fjsp_madd_v2r8(dx00,fscal,fix0);
291             fiy0             = _fjsp_madd_v2r8(dy00,fscal,fiy0);
292             fiz0             = _fjsp_madd_v2r8(dz00,fscal,fiz0);
293             
294             fjx0             = _fjsp_madd_v2r8(dx00,fscal,fjx0);
295             fjy0             = _fjsp_madd_v2r8(dy00,fscal,fjy0);
296             fjz0             = _fjsp_madd_v2r8(dz00,fscal,fjz0);
297
298             /**************************
299              * CALCULATE INTERACTIONS *
300              **************************/
301
302             /* Compute parameters for interactions between i and j atoms */
303             qq10             = _fjsp_mul_v2r8(iq1,jq0);
304
305             /* REACTION-FIELD ELECTROSTATICS */
306             velec            = _fjsp_mul_v2r8(qq10,_fjsp_sub_v2r8(_fjsp_madd_v2r8(krf,rsq10,rinv10),crf));
307             felec            = _fjsp_mul_v2r8(qq10,_fjsp_msub_v2r8(rinv10,rinvsq10,krf2));
308
309             /* Update potential sum for this i atom from the interaction with this j atom. */
310             velecsum         = _fjsp_add_v2r8(velecsum,velec);
311
312             fscal            = felec;
313
314             /* Update vectorial force */
315             fix1             = _fjsp_madd_v2r8(dx10,fscal,fix1);
316             fiy1             = _fjsp_madd_v2r8(dy10,fscal,fiy1);
317             fiz1             = _fjsp_madd_v2r8(dz10,fscal,fiz1);
318             
319             fjx0             = _fjsp_madd_v2r8(dx10,fscal,fjx0);
320             fjy0             = _fjsp_madd_v2r8(dy10,fscal,fjy0);
321             fjz0             = _fjsp_madd_v2r8(dz10,fscal,fjz0);
322
323             /**************************
324              * CALCULATE INTERACTIONS *
325              **************************/
326
327             /* Compute parameters for interactions between i and j atoms */
328             qq20             = _fjsp_mul_v2r8(iq2,jq0);
329
330             /* REACTION-FIELD ELECTROSTATICS */
331             velec            = _fjsp_mul_v2r8(qq20,_fjsp_sub_v2r8(_fjsp_madd_v2r8(krf,rsq20,rinv20),crf));
332             felec            = _fjsp_mul_v2r8(qq20,_fjsp_msub_v2r8(rinv20,rinvsq20,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             fix2             = _fjsp_madd_v2r8(dx20,fscal,fix2);
341             fiy2             = _fjsp_madd_v2r8(dy20,fscal,fiy2);
342             fiz2             = _fjsp_madd_v2r8(dz20,fscal,fiz2);
343             
344             fjx0             = _fjsp_madd_v2r8(dx20,fscal,fjx0);
345             fjy0             = _fjsp_madd_v2r8(dy20,fscal,fjy0);
346             fjz0             = _fjsp_madd_v2r8(dz20,fscal,fjz0);
347
348             /**************************
349              * CALCULATE INTERACTIONS *
350              **************************/
351
352             /* Compute parameters for interactions between i and j atoms */
353             qq30             = _fjsp_mul_v2r8(iq3,jq0);
354
355             /* REACTION-FIELD ELECTROSTATICS */
356             velec            = _fjsp_mul_v2r8(qq30,_fjsp_sub_v2r8(_fjsp_madd_v2r8(krf,rsq30,rinv30),crf));
357             felec            = _fjsp_mul_v2r8(qq30,_fjsp_msub_v2r8(rinv30,rinvsq30,krf2));
358
359             /* Update potential sum for this i atom from the interaction with this j atom. */
360             velecsum         = _fjsp_add_v2r8(velecsum,velec);
361
362             fscal            = felec;
363
364             /* Update vectorial force */
365             fix3             = _fjsp_madd_v2r8(dx30,fscal,fix3);
366             fiy3             = _fjsp_madd_v2r8(dy30,fscal,fiy3);
367             fiz3             = _fjsp_madd_v2r8(dz30,fscal,fiz3);
368             
369             fjx0             = _fjsp_madd_v2r8(dx30,fscal,fjx0);
370             fjy0             = _fjsp_madd_v2r8(dy30,fscal,fjy0);
371             fjz0             = _fjsp_madd_v2r8(dz30,fscal,fjz0);
372
373             gmx_fjsp_decrement_1rvec_2ptr_swizzle_v2r8(f+j_coord_offsetA,f+j_coord_offsetB,fjx0,fjy0,fjz0);
374
375             /* Inner loop uses 167 flops */
376         }
377
378         if(jidx<j_index_end)
379         {
380
381             jnrA             = jjnr[jidx];
382             j_coord_offsetA  = DIM*jnrA;
383
384             /* load j atom coordinates */
385             gmx_fjsp_load_1rvec_1ptr_swizzle_v2r8(x+j_coord_offsetA,
386                                               &jx0,&jy0,&jz0);
387
388             /* Calculate displacement vector */
389             dx00             = _fjsp_sub_v2r8(ix0,jx0);
390             dy00             = _fjsp_sub_v2r8(iy0,jy0);
391             dz00             = _fjsp_sub_v2r8(iz0,jz0);
392             dx10             = _fjsp_sub_v2r8(ix1,jx0);
393             dy10             = _fjsp_sub_v2r8(iy1,jy0);
394             dz10             = _fjsp_sub_v2r8(iz1,jz0);
395             dx20             = _fjsp_sub_v2r8(ix2,jx0);
396             dy20             = _fjsp_sub_v2r8(iy2,jy0);
397             dz20             = _fjsp_sub_v2r8(iz2,jz0);
398             dx30             = _fjsp_sub_v2r8(ix3,jx0);
399             dy30             = _fjsp_sub_v2r8(iy3,jy0);
400             dz30             = _fjsp_sub_v2r8(iz3,jz0);
401
402             /* Calculate squared distance and things based on it */
403             rsq00            = gmx_fjsp_calc_rsq_v2r8(dx00,dy00,dz00);
404             rsq10            = gmx_fjsp_calc_rsq_v2r8(dx10,dy10,dz10);
405             rsq20            = gmx_fjsp_calc_rsq_v2r8(dx20,dy20,dz20);
406             rsq30            = gmx_fjsp_calc_rsq_v2r8(dx30,dy30,dz30);
407
408             rinv00           = gmx_fjsp_invsqrt_v2r8(rsq00);
409             rinv10           = gmx_fjsp_invsqrt_v2r8(rsq10);
410             rinv20           = gmx_fjsp_invsqrt_v2r8(rsq20);
411             rinv30           = gmx_fjsp_invsqrt_v2r8(rsq30);
412
413             rinvsq10         = _fjsp_mul_v2r8(rinv10,rinv10);
414             rinvsq20         = _fjsp_mul_v2r8(rinv20,rinv20);
415             rinvsq30         = _fjsp_mul_v2r8(rinv30,rinv30);
416
417             /* Load parameters for j particles */
418             jq0              = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(),charge+jnrA+0);
419             vdwjidx0A        = 2*vdwtype[jnrA+0];
420
421             fjx0             = _fjsp_setzero_v2r8();
422             fjy0             = _fjsp_setzero_v2r8();
423             fjz0             = _fjsp_setzero_v2r8();
424
425             /**************************
426              * CALCULATE INTERACTIONS *
427              **************************/
428
429             r00              = _fjsp_mul_v2r8(rsq00,rinv00);
430
431             /* Compute parameters for interactions between i and j atoms */
432             gmx_fjsp_load_1pair_swizzle_v2r8(vdwparam+vdwioffset0+vdwjidx0A,&c6_00,&c12_00);
433
434             /* Calculate table index by multiplying r with table scale and truncate to integer */
435             rt               = _fjsp_mul_v2r8(r00,vftabscale);
436             itab_tmp         = _fjsp_dtox_v2r8(rt);
437             vfeps            = _fjsp_sub_v2r8(rt, _fjsp_xtod_v2r8(itab_tmp));
438             twovfeps         = _fjsp_add_v2r8(vfeps,vfeps);
439             _fjsp_store_v2r8(&vfconv.simd,itab_tmp);
440
441             vfconv.i[0]     *= 8;
442             vfconv.i[1]     *= 8;
443
444             /* CUBIC SPLINE TABLE DISPERSION */
445             Y                = _fjsp_load_v2r8( vftab + vfconv.i[0] );
446             F                = _fjsp_setzero_v2r8();
447             GMX_FJSP_TRANSPOSE2_V2R8(Y,F);
448             G                = _fjsp_load_v2r8( vftab + vfconv.i[0] + 2 );
449             H                = _fjsp_setzero_v2r8();
450             GMX_FJSP_TRANSPOSE2_V2R8(G,H);
451             Fp               = _fjsp_madd_v2r8(vfeps,_fjsp_madd_v2r8(H,vfeps,G),F);
452             VV               = _fjsp_madd_v2r8(vfeps,Fp,Y);
453             vvdw6            = _fjsp_mul_v2r8(c6_00,VV);
454             FF               = _fjsp_madd_v2r8(vfeps,_fjsp_madd_v2r8(twovfeps,H,G),Fp);
455             fvdw6            = _fjsp_mul_v2r8(c6_00,FF);
456
457             /* CUBIC SPLINE TABLE REPULSION */
458             Y                = _fjsp_load_v2r8( vftab + vfconv.i[0] + 4 );
459             F                = _fjsp_setzero_v2r8();
460             GMX_FJSP_TRANSPOSE2_V2R8(Y,F);
461             G                = _fjsp_load_v2r8( vftab + vfconv.i[0] + 6 );
462             H                = _fjsp_setzero_v2r8();
463             GMX_FJSP_TRANSPOSE2_V2R8(G,H);
464             Fp               = _fjsp_madd_v2r8(vfeps,_fjsp_madd_v2r8(H,vfeps,G),F);
465             VV               = _fjsp_madd_v2r8(vfeps,Fp,Y);
466             vvdw12           = _fjsp_mul_v2r8(c12_00,VV);
467             FF               = _fjsp_madd_v2r8(vfeps,_fjsp_madd_v2r8(twovfeps,H,G),Fp);
468             fvdw12           = _fjsp_mul_v2r8(c12_00,FF);
469             vvdw             = _fjsp_add_v2r8(vvdw12,vvdw6);
470             fvdw             = _fjsp_neg_v2r8(_fjsp_mul_v2r8(_fjsp_add_v2r8(fvdw6,fvdw12),_fjsp_mul_v2r8(vftabscale,rinv00)));
471
472             /* Update potential sum for this i atom from the interaction with this j atom. */
473             vvdw             = _fjsp_unpacklo_v2r8(vvdw,_fjsp_setzero_v2r8());
474             vvdwsum          = _fjsp_add_v2r8(vvdwsum,vvdw);
475
476             fscal            = fvdw;
477
478             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
479
480             /* Update vectorial force */
481             fix0             = _fjsp_madd_v2r8(dx00,fscal,fix0);
482             fiy0             = _fjsp_madd_v2r8(dy00,fscal,fiy0);
483             fiz0             = _fjsp_madd_v2r8(dz00,fscal,fiz0);
484             
485             fjx0             = _fjsp_madd_v2r8(dx00,fscal,fjx0);
486             fjy0             = _fjsp_madd_v2r8(dy00,fscal,fjy0);
487             fjz0             = _fjsp_madd_v2r8(dz00,fscal,fjz0);
488
489             /**************************
490              * CALCULATE INTERACTIONS *
491              **************************/
492
493             /* Compute parameters for interactions between i and j atoms */
494             qq10             = _fjsp_mul_v2r8(iq1,jq0);
495
496             /* REACTION-FIELD ELECTROSTATICS */
497             velec            = _fjsp_mul_v2r8(qq10,_fjsp_sub_v2r8(_fjsp_madd_v2r8(krf,rsq10,rinv10),crf));
498             felec            = _fjsp_mul_v2r8(qq10,_fjsp_msub_v2r8(rinv10,rinvsq10,krf2));
499
500             /* Update potential sum for this i atom from the interaction with this j atom. */
501             velec            = _fjsp_unpacklo_v2r8(velec,_fjsp_setzero_v2r8());
502             velecsum         = _fjsp_add_v2r8(velecsum,velec);
503
504             fscal            = felec;
505
506             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
507
508             /* Update vectorial force */
509             fix1             = _fjsp_madd_v2r8(dx10,fscal,fix1);
510             fiy1             = _fjsp_madd_v2r8(dy10,fscal,fiy1);
511             fiz1             = _fjsp_madd_v2r8(dz10,fscal,fiz1);
512             
513             fjx0             = _fjsp_madd_v2r8(dx10,fscal,fjx0);
514             fjy0             = _fjsp_madd_v2r8(dy10,fscal,fjy0);
515             fjz0             = _fjsp_madd_v2r8(dz10,fscal,fjz0);
516
517             /**************************
518              * CALCULATE INTERACTIONS *
519              **************************/
520
521             /* Compute parameters for interactions between i and j atoms */
522             qq20             = _fjsp_mul_v2r8(iq2,jq0);
523
524             /* REACTION-FIELD ELECTROSTATICS */
525             velec            = _fjsp_mul_v2r8(qq20,_fjsp_sub_v2r8(_fjsp_madd_v2r8(krf,rsq20,rinv20),crf));
526             felec            = _fjsp_mul_v2r8(qq20,_fjsp_msub_v2r8(rinv20,rinvsq20,krf2));
527
528             /* Update potential sum for this i atom from the interaction with this j atom. */
529             velec            = _fjsp_unpacklo_v2r8(velec,_fjsp_setzero_v2r8());
530             velecsum         = _fjsp_add_v2r8(velecsum,velec);
531
532             fscal            = felec;
533
534             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
535
536             /* Update vectorial force */
537             fix2             = _fjsp_madd_v2r8(dx20,fscal,fix2);
538             fiy2             = _fjsp_madd_v2r8(dy20,fscal,fiy2);
539             fiz2             = _fjsp_madd_v2r8(dz20,fscal,fiz2);
540             
541             fjx0             = _fjsp_madd_v2r8(dx20,fscal,fjx0);
542             fjy0             = _fjsp_madd_v2r8(dy20,fscal,fjy0);
543             fjz0             = _fjsp_madd_v2r8(dz20,fscal,fjz0);
544
545             /**************************
546              * CALCULATE INTERACTIONS *
547              **************************/
548
549             /* Compute parameters for interactions between i and j atoms */
550             qq30             = _fjsp_mul_v2r8(iq3,jq0);
551
552             /* REACTION-FIELD ELECTROSTATICS */
553             velec            = _fjsp_mul_v2r8(qq30,_fjsp_sub_v2r8(_fjsp_madd_v2r8(krf,rsq30,rinv30),crf));
554             felec            = _fjsp_mul_v2r8(qq30,_fjsp_msub_v2r8(rinv30,rinvsq30,krf2));
555
556             /* Update potential sum for this i atom from the interaction with this j atom. */
557             velec            = _fjsp_unpacklo_v2r8(velec,_fjsp_setzero_v2r8());
558             velecsum         = _fjsp_add_v2r8(velecsum,velec);
559
560             fscal            = felec;
561
562             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
563
564             /* Update vectorial force */
565             fix3             = _fjsp_madd_v2r8(dx30,fscal,fix3);
566             fiy3             = _fjsp_madd_v2r8(dy30,fscal,fiy3);
567             fiz3             = _fjsp_madd_v2r8(dz30,fscal,fiz3);
568             
569             fjx0             = _fjsp_madd_v2r8(dx30,fscal,fjx0);
570             fjy0             = _fjsp_madd_v2r8(dy30,fscal,fjy0);
571             fjz0             = _fjsp_madd_v2r8(dz30,fscal,fjz0);
572
573             gmx_fjsp_decrement_1rvec_1ptr_swizzle_v2r8(f+j_coord_offsetA,fjx0,fjy0,fjz0);
574
575             /* Inner loop uses 167 flops */
576         }
577
578         /* End of innermost loop */
579
580         gmx_fjsp_update_iforce_4atom_swizzle_v2r8(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
581                                               f+i_coord_offset,fshift+i_shift_offset);
582
583         ggid                        = gid[iidx];
584         /* Update potential energies */
585         gmx_fjsp_update_1pot_v2r8(velecsum,kernel_data->energygrp_elec+ggid);
586         gmx_fjsp_update_1pot_v2r8(vvdwsum,kernel_data->energygrp_vdw+ggid);
587
588         /* Increment number of inner iterations */
589         inneriter                  += j_index_end - j_index_start;
590
591         /* Outer loop uses 26 flops */
592     }
593
594     /* Increment number of outer iterations */
595     outeriter        += nri;
596
597     /* Update outer/inner flops */
598
599     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W4_VF,outeriter*26 + inneriter*167);
600 }
601 /*
602  * Gromacs nonbonded kernel:   nb_kernel_ElecRF_VdwCSTab_GeomW4P1_F_sparc64_hpc_ace_double
603  * Electrostatics interaction: ReactionField
604  * VdW interaction:            CubicSplineTable
605  * Geometry:                   Water4-Particle
606  * Calculate force/pot:        Force
607  */
608 void
609 nb_kernel_ElecRF_VdwCSTab_GeomW4P1_F_sparc64_hpc_ace_double
610                     (t_nblist                    * gmx_restrict       nlist,
611                      rvec                        * gmx_restrict          xx,
612                      rvec                        * gmx_restrict          ff,
613                      t_forcerec                  * gmx_restrict          fr,
614                      t_mdatoms                   * gmx_restrict     mdatoms,
615                      nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
616                      t_nrnb                      * gmx_restrict        nrnb)
617 {
618     /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
619      * just 0 for non-waters.
620      * Suffixes A,B refer to j loop unrolling done with double precision SIMD, e.g. for the two different
621      * jnr indices corresponding to data put in the four positions in the SIMD register.
622      */
623     int              i_shift_offset,i_coord_offset,outeriter,inneriter;
624     int              j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
625     int              jnrA,jnrB;
626     int              j_coord_offsetA,j_coord_offsetB;
627     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
628     real             rcutoff_scalar;
629     real             *shiftvec,*fshift,*x,*f;
630     _fjsp_v2r8       tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
631     int              vdwioffset0;
632     _fjsp_v2r8       ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
633     int              vdwioffset1;
634     _fjsp_v2r8       ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
635     int              vdwioffset2;
636     _fjsp_v2r8       ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
637     int              vdwioffset3;
638     _fjsp_v2r8       ix3,iy3,iz3,fix3,fiy3,fiz3,iq3,isai3;
639     int              vdwjidx0A,vdwjidx0B;
640     _fjsp_v2r8       jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
641     _fjsp_v2r8       dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
642     _fjsp_v2r8       dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
643     _fjsp_v2r8       dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
644     _fjsp_v2r8       dx30,dy30,dz30,rsq30,rinv30,rinvsq30,r30,qq30,c6_30,c12_30;
645     _fjsp_v2r8       velec,felec,velecsum,facel,crf,krf,krf2;
646     real             *charge;
647     int              nvdwtype;
648     _fjsp_v2r8       rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
649     int              *vdwtype;
650     real             *vdwparam;
651     _fjsp_v2r8       one_sixth   = gmx_fjsp_set1_v2r8(1.0/6.0);
652     _fjsp_v2r8       one_twelfth = gmx_fjsp_set1_v2r8(1.0/12.0);
653     _fjsp_v2r8       rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF,twovfeps;
654     real             *vftab;
655     _fjsp_v2r8       itab_tmp;
656     _fjsp_v2r8       dummy_mask,cutoff_mask;
657     _fjsp_v2r8       one     = gmx_fjsp_set1_v2r8(1.0);
658     _fjsp_v2r8       two     = gmx_fjsp_set1_v2r8(2.0);
659     union { _fjsp_v2r8 simd; long long int i[2]; } vfconv,gbconv,ewconv;
660
661     x                = xx[0];
662     f                = ff[0];
663
664     nri              = nlist->nri;
665     iinr             = nlist->iinr;
666     jindex           = nlist->jindex;
667     jjnr             = nlist->jjnr;
668     shiftidx         = nlist->shift;
669     gid              = nlist->gid;
670     shiftvec         = fr->shift_vec[0];
671     fshift           = fr->fshift[0];
672     facel            = gmx_fjsp_set1_v2r8(fr->epsfac);
673     charge           = mdatoms->chargeA;
674     krf              = gmx_fjsp_set1_v2r8(fr->ic->k_rf);
675     krf2             = gmx_fjsp_set1_v2r8(fr->ic->k_rf*2.0);
676     crf              = gmx_fjsp_set1_v2r8(fr->ic->c_rf);
677     nvdwtype         = fr->ntype;
678     vdwparam         = fr->nbfp;
679     vdwtype          = mdatoms->typeA;
680
681     vftab            = kernel_data->table_vdw->data;
682     vftabscale       = gmx_fjsp_set1_v2r8(kernel_data->table_vdw->scale);
683
684     /* Setup water-specific parameters */
685     inr              = nlist->iinr[0];
686     iq1              = _fjsp_mul_v2r8(facel,gmx_fjsp_set1_v2r8(charge[inr+1]));
687     iq2              = _fjsp_mul_v2r8(facel,gmx_fjsp_set1_v2r8(charge[inr+2]));
688     iq3              = _fjsp_mul_v2r8(facel,gmx_fjsp_set1_v2r8(charge[inr+3]));
689     vdwioffset0      = 2*nvdwtype*vdwtype[inr+0];
690
691     /* Avoid stupid compiler warnings */
692     jnrA = jnrB = 0;
693     j_coord_offsetA = 0;
694     j_coord_offsetB = 0;
695
696     outeriter        = 0;
697     inneriter        = 0;
698
699     /* Start outer loop over neighborlists */
700     for(iidx=0; iidx<nri; iidx++)
701     {
702         /* Load shift vector for this list */
703         i_shift_offset   = DIM*shiftidx[iidx];
704
705         /* Load limits for loop over neighbors */
706         j_index_start    = jindex[iidx];
707         j_index_end      = jindex[iidx+1];
708
709         /* Get outer coordinate index */
710         inr              = iinr[iidx];
711         i_coord_offset   = DIM*inr;
712
713         /* Load i particle coords and add shift vector */
714         gmx_fjsp_load_shift_and_4rvec_broadcast_v2r8(shiftvec+i_shift_offset,x+i_coord_offset,
715                                                  &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
716
717         fix0             = _fjsp_setzero_v2r8();
718         fiy0             = _fjsp_setzero_v2r8();
719         fiz0             = _fjsp_setzero_v2r8();
720         fix1             = _fjsp_setzero_v2r8();
721         fiy1             = _fjsp_setzero_v2r8();
722         fiz1             = _fjsp_setzero_v2r8();
723         fix2             = _fjsp_setzero_v2r8();
724         fiy2             = _fjsp_setzero_v2r8();
725         fiz2             = _fjsp_setzero_v2r8();
726         fix3             = _fjsp_setzero_v2r8();
727         fiy3             = _fjsp_setzero_v2r8();
728         fiz3             = _fjsp_setzero_v2r8();
729
730         /* Start inner kernel loop */
731         for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
732         {
733
734             /* Get j neighbor index, and coordinate index */
735             jnrA             = jjnr[jidx];
736             jnrB             = jjnr[jidx+1];
737             j_coord_offsetA  = DIM*jnrA;
738             j_coord_offsetB  = DIM*jnrB;
739
740             /* load j atom coordinates */
741             gmx_fjsp_load_1rvec_2ptr_swizzle_v2r8(x+j_coord_offsetA,x+j_coord_offsetB,
742                                               &jx0,&jy0,&jz0);
743
744             /* Calculate displacement vector */
745             dx00             = _fjsp_sub_v2r8(ix0,jx0);
746             dy00             = _fjsp_sub_v2r8(iy0,jy0);
747             dz00             = _fjsp_sub_v2r8(iz0,jz0);
748             dx10             = _fjsp_sub_v2r8(ix1,jx0);
749             dy10             = _fjsp_sub_v2r8(iy1,jy0);
750             dz10             = _fjsp_sub_v2r8(iz1,jz0);
751             dx20             = _fjsp_sub_v2r8(ix2,jx0);
752             dy20             = _fjsp_sub_v2r8(iy2,jy0);
753             dz20             = _fjsp_sub_v2r8(iz2,jz0);
754             dx30             = _fjsp_sub_v2r8(ix3,jx0);
755             dy30             = _fjsp_sub_v2r8(iy3,jy0);
756             dz30             = _fjsp_sub_v2r8(iz3,jz0);
757
758             /* Calculate squared distance and things based on it */
759             rsq00            = gmx_fjsp_calc_rsq_v2r8(dx00,dy00,dz00);
760             rsq10            = gmx_fjsp_calc_rsq_v2r8(dx10,dy10,dz10);
761             rsq20            = gmx_fjsp_calc_rsq_v2r8(dx20,dy20,dz20);
762             rsq30            = gmx_fjsp_calc_rsq_v2r8(dx30,dy30,dz30);
763
764             rinv00           = gmx_fjsp_invsqrt_v2r8(rsq00);
765             rinv10           = gmx_fjsp_invsqrt_v2r8(rsq10);
766             rinv20           = gmx_fjsp_invsqrt_v2r8(rsq20);
767             rinv30           = gmx_fjsp_invsqrt_v2r8(rsq30);
768
769             rinvsq10         = _fjsp_mul_v2r8(rinv10,rinv10);
770             rinvsq20         = _fjsp_mul_v2r8(rinv20,rinv20);
771             rinvsq30         = _fjsp_mul_v2r8(rinv30,rinv30);
772
773             /* Load parameters for j particles */
774             jq0              = gmx_fjsp_load_2real_swizzle_v2r8(charge+jnrA+0,charge+jnrB+0);
775             vdwjidx0A        = 2*vdwtype[jnrA+0];
776             vdwjidx0B        = 2*vdwtype[jnrB+0];
777
778             fjx0             = _fjsp_setzero_v2r8();
779             fjy0             = _fjsp_setzero_v2r8();
780             fjz0             = _fjsp_setzero_v2r8();
781
782             /**************************
783              * CALCULATE INTERACTIONS *
784              **************************/
785
786             r00              = _fjsp_mul_v2r8(rsq00,rinv00);
787
788             /* Compute parameters for interactions between i and j atoms */
789             gmx_fjsp_load_2pair_swizzle_v2r8(vdwparam+vdwioffset0+vdwjidx0A,
790                                          vdwparam+vdwioffset0+vdwjidx0B,&c6_00,&c12_00);
791
792             /* Calculate table index by multiplying r with table scale and truncate to integer */
793             rt               = _fjsp_mul_v2r8(r00,vftabscale);
794             itab_tmp         = _fjsp_dtox_v2r8(rt);
795             vfeps            = _fjsp_sub_v2r8(rt, _fjsp_xtod_v2r8(itab_tmp));
796             twovfeps         = _fjsp_add_v2r8(vfeps,vfeps);
797             _fjsp_store_v2r8(&vfconv.simd,itab_tmp);
798
799             vfconv.i[0]     *= 8;
800             vfconv.i[1]     *= 8;
801
802             /* CUBIC SPLINE TABLE DISPERSION */
803             Y                = _fjsp_load_v2r8( vftab + vfconv.i[0] );
804             F                = _fjsp_load_v2r8( vftab + vfconv.i[1] );
805             GMX_FJSP_TRANSPOSE2_V2R8(Y,F);
806             G                = _fjsp_load_v2r8( vftab + vfconv.i[0] + 2 );
807             H                = _fjsp_load_v2r8( vftab + vfconv.i[1] + 2 );
808             GMX_FJSP_TRANSPOSE2_V2R8(G,H);
809             Fp               = _fjsp_madd_v2r8(vfeps,_fjsp_madd_v2r8(H,vfeps,G),F);
810             FF               = _fjsp_madd_v2r8(vfeps,_fjsp_madd_v2r8(twovfeps,H,G),Fp);
811             fvdw6            = _fjsp_mul_v2r8(c6_00,FF);
812
813             /* CUBIC SPLINE TABLE REPULSION */
814             Y                = _fjsp_load_v2r8( vftab + vfconv.i[0] + 4 );
815             F                = _fjsp_load_v2r8( vftab + vfconv.i[1] + 4 );
816             GMX_FJSP_TRANSPOSE2_V2R8(Y,F);
817             G                = _fjsp_load_v2r8( vftab + vfconv.i[0] + 6 );
818             H                = _fjsp_load_v2r8( vftab + vfconv.i[1] + 6 );
819             GMX_FJSP_TRANSPOSE2_V2R8(G,H);
820             Fp               = _fjsp_madd_v2r8(vfeps,_fjsp_madd_v2r8(H,vfeps,G),F);
821             FF               = _fjsp_madd_v2r8(vfeps,_fjsp_madd_v2r8(twovfeps,H,G),Fp);
822             fvdw12           = _fjsp_mul_v2r8(c12_00,FF);
823             fvdw             = _fjsp_neg_v2r8(_fjsp_mul_v2r8(_fjsp_add_v2r8(fvdw6,fvdw12),_fjsp_mul_v2r8(vftabscale,rinv00)));
824
825             fscal            = fvdw;
826
827             /* Update vectorial force */
828             fix0             = _fjsp_madd_v2r8(dx00,fscal,fix0);
829             fiy0             = _fjsp_madd_v2r8(dy00,fscal,fiy0);
830             fiz0             = _fjsp_madd_v2r8(dz00,fscal,fiz0);
831             
832             fjx0             = _fjsp_madd_v2r8(dx00,fscal,fjx0);
833             fjy0             = _fjsp_madd_v2r8(dy00,fscal,fjy0);
834             fjz0             = _fjsp_madd_v2r8(dz00,fscal,fjz0);
835
836             /**************************
837              * CALCULATE INTERACTIONS *
838              **************************/
839
840             /* Compute parameters for interactions between i and j atoms */
841             qq10             = _fjsp_mul_v2r8(iq1,jq0);
842
843             /* REACTION-FIELD ELECTROSTATICS */
844             felec            = _fjsp_mul_v2r8(qq10,_fjsp_msub_v2r8(rinv10,rinvsq10,krf2));
845
846             fscal            = felec;
847
848             /* Update vectorial force */
849             fix1             = _fjsp_madd_v2r8(dx10,fscal,fix1);
850             fiy1             = _fjsp_madd_v2r8(dy10,fscal,fiy1);
851             fiz1             = _fjsp_madd_v2r8(dz10,fscal,fiz1);
852             
853             fjx0             = _fjsp_madd_v2r8(dx10,fscal,fjx0);
854             fjy0             = _fjsp_madd_v2r8(dy10,fscal,fjy0);
855             fjz0             = _fjsp_madd_v2r8(dz10,fscal,fjz0);
856
857             /**************************
858              * CALCULATE INTERACTIONS *
859              **************************/
860
861             /* Compute parameters for interactions between i and j atoms */
862             qq20             = _fjsp_mul_v2r8(iq2,jq0);
863
864             /* REACTION-FIELD ELECTROSTATICS */
865             felec            = _fjsp_mul_v2r8(qq20,_fjsp_msub_v2r8(rinv20,rinvsq20,krf2));
866
867             fscal            = felec;
868
869             /* Update vectorial force */
870             fix2             = _fjsp_madd_v2r8(dx20,fscal,fix2);
871             fiy2             = _fjsp_madd_v2r8(dy20,fscal,fiy2);
872             fiz2             = _fjsp_madd_v2r8(dz20,fscal,fiz2);
873             
874             fjx0             = _fjsp_madd_v2r8(dx20,fscal,fjx0);
875             fjy0             = _fjsp_madd_v2r8(dy20,fscal,fjy0);
876             fjz0             = _fjsp_madd_v2r8(dz20,fscal,fjz0);
877
878             /**************************
879              * CALCULATE INTERACTIONS *
880              **************************/
881
882             /* Compute parameters for interactions between i and j atoms */
883             qq30             = _fjsp_mul_v2r8(iq3,jq0);
884
885             /* REACTION-FIELD ELECTROSTATICS */
886             felec            = _fjsp_mul_v2r8(qq30,_fjsp_msub_v2r8(rinv30,rinvsq30,krf2));
887
888             fscal            = felec;
889
890             /* Update vectorial force */
891             fix3             = _fjsp_madd_v2r8(dx30,fscal,fix3);
892             fiy3             = _fjsp_madd_v2r8(dy30,fscal,fiy3);
893             fiz3             = _fjsp_madd_v2r8(dz30,fscal,fiz3);
894             
895             fjx0             = _fjsp_madd_v2r8(dx30,fscal,fjx0);
896             fjy0             = _fjsp_madd_v2r8(dy30,fscal,fjy0);
897             fjz0             = _fjsp_madd_v2r8(dz30,fscal,fjz0);
898
899             gmx_fjsp_decrement_1rvec_2ptr_swizzle_v2r8(f+j_coord_offsetA,f+j_coord_offsetB,fjx0,fjy0,fjz0);
900
901             /* Inner loop uses 144 flops */
902         }
903
904         if(jidx<j_index_end)
905         {
906
907             jnrA             = jjnr[jidx];
908             j_coord_offsetA  = DIM*jnrA;
909
910             /* load j atom coordinates */
911             gmx_fjsp_load_1rvec_1ptr_swizzle_v2r8(x+j_coord_offsetA,
912                                               &jx0,&jy0,&jz0);
913
914             /* Calculate displacement vector */
915             dx00             = _fjsp_sub_v2r8(ix0,jx0);
916             dy00             = _fjsp_sub_v2r8(iy0,jy0);
917             dz00             = _fjsp_sub_v2r8(iz0,jz0);
918             dx10             = _fjsp_sub_v2r8(ix1,jx0);
919             dy10             = _fjsp_sub_v2r8(iy1,jy0);
920             dz10             = _fjsp_sub_v2r8(iz1,jz0);
921             dx20             = _fjsp_sub_v2r8(ix2,jx0);
922             dy20             = _fjsp_sub_v2r8(iy2,jy0);
923             dz20             = _fjsp_sub_v2r8(iz2,jz0);
924             dx30             = _fjsp_sub_v2r8(ix3,jx0);
925             dy30             = _fjsp_sub_v2r8(iy3,jy0);
926             dz30             = _fjsp_sub_v2r8(iz3,jz0);
927
928             /* Calculate squared distance and things based on it */
929             rsq00            = gmx_fjsp_calc_rsq_v2r8(dx00,dy00,dz00);
930             rsq10            = gmx_fjsp_calc_rsq_v2r8(dx10,dy10,dz10);
931             rsq20            = gmx_fjsp_calc_rsq_v2r8(dx20,dy20,dz20);
932             rsq30            = gmx_fjsp_calc_rsq_v2r8(dx30,dy30,dz30);
933
934             rinv00           = gmx_fjsp_invsqrt_v2r8(rsq00);
935             rinv10           = gmx_fjsp_invsqrt_v2r8(rsq10);
936             rinv20           = gmx_fjsp_invsqrt_v2r8(rsq20);
937             rinv30           = gmx_fjsp_invsqrt_v2r8(rsq30);
938
939             rinvsq10         = _fjsp_mul_v2r8(rinv10,rinv10);
940             rinvsq20         = _fjsp_mul_v2r8(rinv20,rinv20);
941             rinvsq30         = _fjsp_mul_v2r8(rinv30,rinv30);
942
943             /* Load parameters for j particles */
944             jq0              = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(),charge+jnrA+0);
945             vdwjidx0A        = 2*vdwtype[jnrA+0];
946
947             fjx0             = _fjsp_setzero_v2r8();
948             fjy0             = _fjsp_setzero_v2r8();
949             fjz0             = _fjsp_setzero_v2r8();
950
951             /**************************
952              * CALCULATE INTERACTIONS *
953              **************************/
954
955             r00              = _fjsp_mul_v2r8(rsq00,rinv00);
956
957             /* Compute parameters for interactions between i and j atoms */
958             gmx_fjsp_load_1pair_swizzle_v2r8(vdwparam+vdwioffset0+vdwjidx0A,&c6_00,&c12_00);
959
960             /* Calculate table index by multiplying r with table scale and truncate to integer */
961             rt               = _fjsp_mul_v2r8(r00,vftabscale);
962             itab_tmp         = _fjsp_dtox_v2r8(rt);
963             vfeps            = _fjsp_sub_v2r8(rt, _fjsp_xtod_v2r8(itab_tmp));
964             twovfeps         = _fjsp_add_v2r8(vfeps,vfeps);
965             _fjsp_store_v2r8(&vfconv.simd,itab_tmp);
966
967             vfconv.i[0]     *= 8;
968             vfconv.i[1]     *= 8;
969
970             /* CUBIC SPLINE TABLE DISPERSION */
971             Y                = _fjsp_load_v2r8( vftab + vfconv.i[0] );
972             F                = _fjsp_setzero_v2r8();
973             GMX_FJSP_TRANSPOSE2_V2R8(Y,F);
974             G                = _fjsp_load_v2r8( vftab + vfconv.i[0] + 2 );
975             H                = _fjsp_setzero_v2r8();
976             GMX_FJSP_TRANSPOSE2_V2R8(G,H);
977             Fp               = _fjsp_madd_v2r8(vfeps,_fjsp_madd_v2r8(H,vfeps,G),F);
978             FF               = _fjsp_madd_v2r8(vfeps,_fjsp_madd_v2r8(twovfeps,H,G),Fp);
979             fvdw6            = _fjsp_mul_v2r8(c6_00,FF);
980
981             /* CUBIC SPLINE TABLE REPULSION */
982             Y                = _fjsp_load_v2r8( vftab + vfconv.i[0] + 4 );
983             F                = _fjsp_setzero_v2r8();
984             GMX_FJSP_TRANSPOSE2_V2R8(Y,F);
985             G                = _fjsp_load_v2r8( vftab + vfconv.i[0] + 6 );
986             H                = _fjsp_setzero_v2r8();
987             GMX_FJSP_TRANSPOSE2_V2R8(G,H);
988             Fp               = _fjsp_madd_v2r8(vfeps,_fjsp_madd_v2r8(H,vfeps,G),F);
989             FF               = _fjsp_madd_v2r8(vfeps,_fjsp_madd_v2r8(twovfeps,H,G),Fp);
990             fvdw12           = _fjsp_mul_v2r8(c12_00,FF);
991             fvdw             = _fjsp_neg_v2r8(_fjsp_mul_v2r8(_fjsp_add_v2r8(fvdw6,fvdw12),_fjsp_mul_v2r8(vftabscale,rinv00)));
992
993             fscal            = fvdw;
994
995             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
996
997             /* Update vectorial force */
998             fix0             = _fjsp_madd_v2r8(dx00,fscal,fix0);
999             fiy0             = _fjsp_madd_v2r8(dy00,fscal,fiy0);
1000             fiz0             = _fjsp_madd_v2r8(dz00,fscal,fiz0);
1001             
1002             fjx0             = _fjsp_madd_v2r8(dx00,fscal,fjx0);
1003             fjy0             = _fjsp_madd_v2r8(dy00,fscal,fjy0);
1004             fjz0             = _fjsp_madd_v2r8(dz00,fscal,fjz0);
1005
1006             /**************************
1007              * CALCULATE INTERACTIONS *
1008              **************************/
1009
1010             /* Compute parameters for interactions between i and j atoms */
1011             qq10             = _fjsp_mul_v2r8(iq1,jq0);
1012
1013             /* REACTION-FIELD ELECTROSTATICS */
1014             felec            = _fjsp_mul_v2r8(qq10,_fjsp_msub_v2r8(rinv10,rinvsq10,krf2));
1015
1016             fscal            = felec;
1017
1018             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
1019
1020             /* Update vectorial force */
1021             fix1             = _fjsp_madd_v2r8(dx10,fscal,fix1);
1022             fiy1             = _fjsp_madd_v2r8(dy10,fscal,fiy1);
1023             fiz1             = _fjsp_madd_v2r8(dz10,fscal,fiz1);
1024             
1025             fjx0             = _fjsp_madd_v2r8(dx10,fscal,fjx0);
1026             fjy0             = _fjsp_madd_v2r8(dy10,fscal,fjy0);
1027             fjz0             = _fjsp_madd_v2r8(dz10,fscal,fjz0);
1028
1029             /**************************
1030              * CALCULATE INTERACTIONS *
1031              **************************/
1032
1033             /* Compute parameters for interactions between i and j atoms */
1034             qq20             = _fjsp_mul_v2r8(iq2,jq0);
1035
1036             /* REACTION-FIELD ELECTROSTATICS */
1037             felec            = _fjsp_mul_v2r8(qq20,_fjsp_msub_v2r8(rinv20,rinvsq20,krf2));
1038
1039             fscal            = felec;
1040
1041             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
1042
1043             /* Update vectorial force */
1044             fix2             = _fjsp_madd_v2r8(dx20,fscal,fix2);
1045             fiy2             = _fjsp_madd_v2r8(dy20,fscal,fiy2);
1046             fiz2             = _fjsp_madd_v2r8(dz20,fscal,fiz2);
1047             
1048             fjx0             = _fjsp_madd_v2r8(dx20,fscal,fjx0);
1049             fjy0             = _fjsp_madd_v2r8(dy20,fscal,fjy0);
1050             fjz0             = _fjsp_madd_v2r8(dz20,fscal,fjz0);
1051
1052             /**************************
1053              * CALCULATE INTERACTIONS *
1054              **************************/
1055
1056             /* Compute parameters for interactions between i and j atoms */
1057             qq30             = _fjsp_mul_v2r8(iq3,jq0);
1058
1059             /* REACTION-FIELD ELECTROSTATICS */
1060             felec            = _fjsp_mul_v2r8(qq30,_fjsp_msub_v2r8(rinv30,rinvsq30,krf2));
1061
1062             fscal            = felec;
1063
1064             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
1065
1066             /* Update vectorial force */
1067             fix3             = _fjsp_madd_v2r8(dx30,fscal,fix3);
1068             fiy3             = _fjsp_madd_v2r8(dy30,fscal,fiy3);
1069             fiz3             = _fjsp_madd_v2r8(dz30,fscal,fiz3);
1070             
1071             fjx0             = _fjsp_madd_v2r8(dx30,fscal,fjx0);
1072             fjy0             = _fjsp_madd_v2r8(dy30,fscal,fjy0);
1073             fjz0             = _fjsp_madd_v2r8(dz30,fscal,fjz0);
1074
1075             gmx_fjsp_decrement_1rvec_1ptr_swizzle_v2r8(f+j_coord_offsetA,fjx0,fjy0,fjz0);
1076
1077             /* Inner loop uses 144 flops */
1078         }
1079
1080         /* End of innermost loop */
1081
1082         gmx_fjsp_update_iforce_4atom_swizzle_v2r8(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
1083                                               f+i_coord_offset,fshift+i_shift_offset);
1084
1085         /* Increment number of inner iterations */
1086         inneriter                  += j_index_end - j_index_start;
1087
1088         /* Outer loop uses 24 flops */
1089     }
1090
1091     /* Increment number of outer iterations */
1092     outeriter        += nri;
1093
1094     /* Update outer/inner flops */
1095
1096     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W4_F,outeriter*24 + inneriter*144);
1097 }