Remove all unnecessary HAVE_CONFIG_H
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_sparc64_hpc_ace_double / nb_kernel_ElecCSTab_VdwCSTab_GeomP1P1_sparc64_hpc_ace_double.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2013,2014, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35 /*
36  * Note: this file was generated by the GROMACS sparc64_hpc_ace_double kernel generator.
37  */
38 #include "config.h"
39
40 #include <math.h>
41
42 #include "../nb_kernel.h"
43 #include "types/simple.h"
44 #include "gromacs/math/vec.h"
45 #include "nrnb.h"
46
47 #include "kernelutil_sparc64_hpc_ace_double.h"
48
49 /*
50  * Gromacs nonbonded kernel:   nb_kernel_ElecCSTab_VdwCSTab_GeomP1P1_VF_sparc64_hpc_ace_double
51  * Electrostatics interaction: CubicSplineTable
52  * VdW interaction:            CubicSplineTable
53  * Geometry:                   Particle-Particle
54  * Calculate force/pot:        PotentialAndForce
55  */
56 void
57 nb_kernel_ElecCSTab_VdwCSTab_GeomP1P1_VF_sparc64_hpc_ace_double
58                     (t_nblist                    * gmx_restrict       nlist,
59                      rvec                        * gmx_restrict          xx,
60                      rvec                        * gmx_restrict          ff,
61                      t_forcerec                  * gmx_restrict          fr,
62                      t_mdatoms                   * gmx_restrict     mdatoms,
63                      nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
64                      t_nrnb                      * gmx_restrict        nrnb)
65 {
66     /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
67      * just 0 for non-waters.
68      * Suffixes A,B refer to j loop unrolling done with double precision SIMD, e.g. for the two different
69      * jnr indices corresponding to data put in the four positions in the SIMD register.
70      */
71     int              i_shift_offset,i_coord_offset,outeriter,inneriter;
72     int              j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
73     int              jnrA,jnrB;
74     int              j_coord_offsetA,j_coord_offsetB;
75     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
76     real             rcutoff_scalar;
77     real             *shiftvec,*fshift,*x,*f;
78     _fjsp_v2r8       tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
79     int              vdwioffset0;
80     _fjsp_v2r8       ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
81     int              vdwjidx0A,vdwjidx0B;
82     _fjsp_v2r8       jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
83     _fjsp_v2r8       dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
84     _fjsp_v2r8       velec,felec,velecsum,facel,crf,krf,krf2;
85     real             *charge;
86     int              nvdwtype;
87     _fjsp_v2r8       rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
88     int              *vdwtype;
89     real             *vdwparam;
90     _fjsp_v2r8       one_sixth   = gmx_fjsp_set1_v2r8(1.0/6.0);
91     _fjsp_v2r8       one_twelfth = gmx_fjsp_set1_v2r8(1.0/12.0);
92     _fjsp_v2r8       rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF,twovfeps;
93     real             *vftab;
94     _fjsp_v2r8       itab_tmp;
95     _fjsp_v2r8       dummy_mask,cutoff_mask;
96     _fjsp_v2r8       one     = gmx_fjsp_set1_v2r8(1.0);
97     _fjsp_v2r8       two     = gmx_fjsp_set1_v2r8(2.0);
98     union { _fjsp_v2r8 simd; long long int i[2]; } vfconv,gbconv,ewconv;
99
100     x                = xx[0];
101     f                = ff[0];
102
103     nri              = nlist->nri;
104     iinr             = nlist->iinr;
105     jindex           = nlist->jindex;
106     jjnr             = nlist->jjnr;
107     shiftidx         = nlist->shift;
108     gid              = nlist->gid;
109     shiftvec         = fr->shift_vec[0];
110     fshift           = fr->fshift[0];
111     facel            = gmx_fjsp_set1_v2r8(fr->epsfac);
112     charge           = mdatoms->chargeA;
113     nvdwtype         = fr->ntype;
114     vdwparam         = fr->nbfp;
115     vdwtype          = mdatoms->typeA;
116
117     vftab            = kernel_data->table_elec_vdw->data;
118     vftabscale       = gmx_fjsp_set1_v2r8(kernel_data->table_elec_vdw->scale);
119
120     /* Avoid stupid compiler warnings */
121     jnrA = jnrB = 0;
122     j_coord_offsetA = 0;
123     j_coord_offsetB = 0;
124
125     outeriter        = 0;
126     inneriter        = 0;
127
128     /* Start outer loop over neighborlists */
129     for(iidx=0; iidx<nri; iidx++)
130     {
131         /* Load shift vector for this list */
132         i_shift_offset   = DIM*shiftidx[iidx];
133
134         /* Load limits for loop over neighbors */
135         j_index_start    = jindex[iidx];
136         j_index_end      = jindex[iidx+1];
137
138         /* Get outer coordinate index */
139         inr              = iinr[iidx];
140         i_coord_offset   = DIM*inr;
141
142         /* Load i particle coords and add shift vector */
143         gmx_fjsp_load_shift_and_1rvec_broadcast_v2r8(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
144
145         fix0             = _fjsp_setzero_v2r8();
146         fiy0             = _fjsp_setzero_v2r8();
147         fiz0             = _fjsp_setzero_v2r8();
148
149         /* Load parameters for i particles */
150         iq0              = _fjsp_mul_v2r8(facel,gmx_fjsp_load1_v2r8(charge+inr+0));
151         vdwioffset0      = 2*nvdwtype*vdwtype[inr+0];
152
153         /* Reset potential sums */
154         velecsum         = _fjsp_setzero_v2r8();
155         vvdwsum          = _fjsp_setzero_v2r8();
156
157         /* Start inner kernel loop */
158         for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
159         {
160
161             /* Get j neighbor index, and coordinate index */
162             jnrA             = jjnr[jidx];
163             jnrB             = jjnr[jidx+1];
164             j_coord_offsetA  = DIM*jnrA;
165             j_coord_offsetB  = DIM*jnrB;
166
167             /* load j atom coordinates */
168             gmx_fjsp_load_1rvec_2ptr_swizzle_v2r8(x+j_coord_offsetA,x+j_coord_offsetB,
169                                               &jx0,&jy0,&jz0);
170
171             /* Calculate displacement vector */
172             dx00             = _fjsp_sub_v2r8(ix0,jx0);
173             dy00             = _fjsp_sub_v2r8(iy0,jy0);
174             dz00             = _fjsp_sub_v2r8(iz0,jz0);
175
176             /* Calculate squared distance and things based on it */
177             rsq00            = gmx_fjsp_calc_rsq_v2r8(dx00,dy00,dz00);
178
179             rinv00           = gmx_fjsp_invsqrt_v2r8(rsq00);
180
181             /* Load parameters for j particles */
182             jq0              = gmx_fjsp_load_2real_swizzle_v2r8(charge+jnrA+0,charge+jnrB+0);
183             vdwjidx0A        = 2*vdwtype[jnrA+0];
184             vdwjidx0B        = 2*vdwtype[jnrB+0];
185
186             /**************************
187              * CALCULATE INTERACTIONS *
188              **************************/
189
190             r00              = _fjsp_mul_v2r8(rsq00,rinv00);
191
192             /* Compute parameters for interactions between i and j atoms */
193             qq00             = _fjsp_mul_v2r8(iq0,jq0);
194             gmx_fjsp_load_2pair_swizzle_v2r8(vdwparam+vdwioffset0+vdwjidx0A,
195                                          vdwparam+vdwioffset0+vdwjidx0B,&c6_00,&c12_00);
196
197             /* Calculate table index by multiplying r with table scale and truncate to integer */
198             rt               = _fjsp_mul_v2r8(r00,vftabscale);
199             itab_tmp         = _fjsp_dtox_v2r8(rt);
200             vfeps            = _fjsp_sub_v2r8(rt, _fjsp_xtod_v2r8(itab_tmp));
201             twovfeps         = _fjsp_add_v2r8(vfeps,vfeps);
202             _fjsp_store_v2r8(&vfconv.simd,itab_tmp);
203
204             vfconv.i[0]     *= 12;
205             vfconv.i[1]     *= 12;
206
207             /* CUBIC SPLINE TABLE ELECTROSTATICS */
208             Y                = _fjsp_load_v2r8( vftab + vfconv.i[0] );
209             F                = _fjsp_load_v2r8( vftab + vfconv.i[1] );
210             GMX_FJSP_TRANSPOSE2_V2R8(Y,F);
211             G                = _fjsp_load_v2r8( vftab + vfconv.i[0] +2);
212             H                = _fjsp_load_v2r8( vftab + vfconv.i[1] +2);
213             GMX_FJSP_TRANSPOSE2_V2R8(G,H);
214             Fp               = _fjsp_madd_v2r8(vfeps,_fjsp_madd_v2r8(vfeps,H,G),F);
215             VV               = _fjsp_madd_v2r8(vfeps,Fp,Y);
216             velec            = _fjsp_mul_v2r8(qq00,VV);
217             FF               = _fjsp_madd_v2r8(_fjsp_madd_v2r8(twovfeps,H,G),vfeps,Fp);
218             felec            = _fjsp_neg_v2r8(_fjsp_mul_v2r8(_fjsp_mul_v2r8(qq00,FF),_fjsp_mul_v2r8(vftabscale,rinv00)));
219
220             /* CUBIC SPLINE TABLE DISPERSION */
221             vfconv.i[0]       += 4;
222             vfconv.i[1]       += 4;
223             Y                = _fjsp_load_v2r8( vftab + vfconv.i[0] );
224             F                = _fjsp_load_v2r8( vftab + vfconv.i[1] );
225             GMX_FJSP_TRANSPOSE2_V2R8(Y,F);
226             G                = _fjsp_load_v2r8( vftab + vfconv.i[0] + 2 );
227             H                = _fjsp_load_v2r8( vftab + vfconv.i[1] + 2 );
228             GMX_FJSP_TRANSPOSE2_V2R8(G,H);
229             Fp               = _fjsp_madd_v2r8(vfeps,_fjsp_madd_v2r8(H,vfeps,G),F);
230             VV               = _fjsp_madd_v2r8(vfeps,Fp,Y);
231             vvdw6            = _fjsp_mul_v2r8(c6_00,VV);
232             FF               = _fjsp_madd_v2r8(vfeps,_fjsp_madd_v2r8(twovfeps,H,G),Fp);
233             fvdw6            = _fjsp_mul_v2r8(c6_00,FF);
234
235             /* CUBIC SPLINE TABLE REPULSION */
236             Y                = _fjsp_load_v2r8( vftab + vfconv.i[0] + 4 );
237             F                = _fjsp_load_v2r8( vftab + vfconv.i[1] + 4 );
238             GMX_FJSP_TRANSPOSE2_V2R8(Y,F);
239             G                = _fjsp_load_v2r8( vftab + vfconv.i[0] + 6 );
240             H                = _fjsp_load_v2r8( vftab + vfconv.i[1] + 6 );
241             GMX_FJSP_TRANSPOSE2_V2R8(G,H);
242             Fp               = _fjsp_madd_v2r8(vfeps,_fjsp_madd_v2r8(H,vfeps,G),F);
243             VV               = _fjsp_madd_v2r8(vfeps,Fp,Y);
244             vvdw12           = _fjsp_mul_v2r8(c12_00,VV);
245             FF               = _fjsp_madd_v2r8(vfeps,_fjsp_madd_v2r8(twovfeps,H,G),Fp);
246             fvdw12           = _fjsp_mul_v2r8(c12_00,FF);
247             vvdw             = _fjsp_add_v2r8(vvdw12,vvdw6);
248             fvdw             = _fjsp_neg_v2r8(_fjsp_mul_v2r8(_fjsp_add_v2r8(fvdw6,fvdw12),_fjsp_mul_v2r8(vftabscale,rinv00)));
249
250             /* Update potential sum for this i atom from the interaction with this j atom. */
251             velecsum         = _fjsp_add_v2r8(velecsum,velec);
252             vvdwsum          = _fjsp_add_v2r8(vvdwsum,vvdw);
253
254             fscal            = _fjsp_add_v2r8(felec,fvdw);
255
256             /* Update vectorial force */
257             fix0             = _fjsp_madd_v2r8(dx00,fscal,fix0);
258             fiy0             = _fjsp_madd_v2r8(dy00,fscal,fiy0);
259             fiz0             = _fjsp_madd_v2r8(dz00,fscal,fiz0);
260             
261             gmx_fjsp_decrement_fma_1rvec_2ptr_swizzle_v2r8(f+j_coord_offsetA,f+j_coord_offsetB,fscal,dx00,dy00,dz00);
262
263             /* Inner loop uses 76 flops */
264         }
265
266         if(jidx<j_index_end)
267         {
268
269             jnrA             = jjnr[jidx];
270             j_coord_offsetA  = DIM*jnrA;
271
272             /* load j atom coordinates */
273             gmx_fjsp_load_1rvec_1ptr_swizzle_v2r8(x+j_coord_offsetA,
274                                               &jx0,&jy0,&jz0);
275
276             /* Calculate displacement vector */
277             dx00             = _fjsp_sub_v2r8(ix0,jx0);
278             dy00             = _fjsp_sub_v2r8(iy0,jy0);
279             dz00             = _fjsp_sub_v2r8(iz0,jz0);
280
281             /* Calculate squared distance and things based on it */
282             rsq00            = gmx_fjsp_calc_rsq_v2r8(dx00,dy00,dz00);
283
284             rinv00           = gmx_fjsp_invsqrt_v2r8(rsq00);
285
286             /* Load parameters for j particles */
287             jq0              = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(),charge+jnrA+0);
288             vdwjidx0A        = 2*vdwtype[jnrA+0];
289
290             /**************************
291              * CALCULATE INTERACTIONS *
292              **************************/
293
294             r00              = _fjsp_mul_v2r8(rsq00,rinv00);
295
296             /* Compute parameters for interactions between i and j atoms */
297             qq00             = _fjsp_mul_v2r8(iq0,jq0);
298             gmx_fjsp_load_1pair_swizzle_v2r8(vdwparam+vdwioffset0+vdwjidx0A,&c6_00,&c12_00);
299
300             /* Calculate table index by multiplying r with table scale and truncate to integer */
301             rt               = _fjsp_mul_v2r8(r00,vftabscale);
302             itab_tmp         = _fjsp_dtox_v2r8(rt);
303             vfeps            = _fjsp_sub_v2r8(rt, _fjsp_xtod_v2r8(itab_tmp));
304             twovfeps         = _fjsp_add_v2r8(vfeps,vfeps);
305             _fjsp_store_v2r8(&vfconv.simd,itab_tmp);
306
307             vfconv.i[0]     *= 12;
308             vfconv.i[1]     *= 12;
309
310             /* CUBIC SPLINE TABLE ELECTROSTATICS */
311             Y                = _fjsp_load_v2r8( vftab + vfconv.i[0] );
312             F                = _fjsp_setzero_v2r8();
313             GMX_FJSP_TRANSPOSE2_V2R8(Y,F);
314             G                = _fjsp_load_v2r8( vftab + vfconv.i[0] +2);
315             H                = _fjsp_setzero_v2r8();
316             GMX_FJSP_TRANSPOSE2_V2R8(G,H);
317             Fp               = _fjsp_madd_v2r8(vfeps,_fjsp_madd_v2r8(vfeps,H,G),F);
318             VV               = _fjsp_madd_v2r8(vfeps,Fp,Y);
319             velec            = _fjsp_mul_v2r8(qq00,VV);
320             FF               = _fjsp_madd_v2r8(_fjsp_madd_v2r8(twovfeps,H,G),vfeps,Fp);
321             felec            = _fjsp_neg_v2r8(_fjsp_mul_v2r8(_fjsp_mul_v2r8(qq00,FF),_fjsp_mul_v2r8(vftabscale,rinv00)));
322
323             /* CUBIC SPLINE TABLE DISPERSION */
324             vfconv.i[0]       += 4;
325             vfconv.i[1]       += 4;
326             Y                = _fjsp_load_v2r8( vftab + vfconv.i[0] );
327             F                = _fjsp_setzero_v2r8();
328             GMX_FJSP_TRANSPOSE2_V2R8(Y,F);
329             G                = _fjsp_load_v2r8( vftab + vfconv.i[0] + 2 );
330             H                = _fjsp_setzero_v2r8();
331             GMX_FJSP_TRANSPOSE2_V2R8(G,H);
332             Fp               = _fjsp_madd_v2r8(vfeps,_fjsp_madd_v2r8(H,vfeps,G),F);
333             VV               = _fjsp_madd_v2r8(vfeps,Fp,Y);
334             vvdw6            = _fjsp_mul_v2r8(c6_00,VV);
335             FF               = _fjsp_madd_v2r8(vfeps,_fjsp_madd_v2r8(twovfeps,H,G),Fp);
336             fvdw6            = _fjsp_mul_v2r8(c6_00,FF);
337
338             /* CUBIC SPLINE TABLE REPULSION */
339             Y                = _fjsp_load_v2r8( vftab + vfconv.i[0] + 4 );
340             F                = _fjsp_setzero_v2r8();
341             GMX_FJSP_TRANSPOSE2_V2R8(Y,F);
342             G                = _fjsp_load_v2r8( vftab + vfconv.i[0] + 6 );
343             H                = _fjsp_setzero_v2r8();
344             GMX_FJSP_TRANSPOSE2_V2R8(G,H);
345             Fp               = _fjsp_madd_v2r8(vfeps,_fjsp_madd_v2r8(H,vfeps,G),F);
346             VV               = _fjsp_madd_v2r8(vfeps,Fp,Y);
347             vvdw12           = _fjsp_mul_v2r8(c12_00,VV);
348             FF               = _fjsp_madd_v2r8(vfeps,_fjsp_madd_v2r8(twovfeps,H,G),Fp);
349             fvdw12           = _fjsp_mul_v2r8(c12_00,FF);
350             vvdw             = _fjsp_add_v2r8(vvdw12,vvdw6);
351             fvdw             = _fjsp_neg_v2r8(_fjsp_mul_v2r8(_fjsp_add_v2r8(fvdw6,fvdw12),_fjsp_mul_v2r8(vftabscale,rinv00)));
352
353             /* Update potential sum for this i atom from the interaction with this j atom. */
354             velec            = _fjsp_unpacklo_v2r8(velec,_fjsp_setzero_v2r8());
355             velecsum         = _fjsp_add_v2r8(velecsum,velec);
356             vvdw             = _fjsp_unpacklo_v2r8(vvdw,_fjsp_setzero_v2r8());
357             vvdwsum          = _fjsp_add_v2r8(vvdwsum,vvdw);
358
359             fscal            = _fjsp_add_v2r8(felec,fvdw);
360
361             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
362
363             /* Update vectorial force */
364             fix0             = _fjsp_madd_v2r8(dx00,fscal,fix0);
365             fiy0             = _fjsp_madd_v2r8(dy00,fscal,fiy0);
366             fiz0             = _fjsp_madd_v2r8(dz00,fscal,fiz0);
367             
368             gmx_fjsp_decrement_fma_1rvec_1ptr_swizzle_v2r8(f+j_coord_offsetA,fscal,dx00,dy00,dz00);
369
370             /* Inner loop uses 76 flops */
371         }
372
373         /* End of innermost loop */
374
375         gmx_fjsp_update_iforce_1atom_swizzle_v2r8(fix0,fiy0,fiz0,
376                                               f+i_coord_offset,fshift+i_shift_offset);
377
378         ggid                        = gid[iidx];
379         /* Update potential energies */
380         gmx_fjsp_update_1pot_v2r8(velecsum,kernel_data->energygrp_elec+ggid);
381         gmx_fjsp_update_1pot_v2r8(vvdwsum,kernel_data->energygrp_vdw+ggid);
382
383         /* Increment number of inner iterations */
384         inneriter                  += j_index_end - j_index_start;
385
386         /* Outer loop uses 9 flops */
387     }
388
389     /* Increment number of outer iterations */
390     outeriter        += nri;
391
392     /* Update outer/inner flops */
393
394     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_VF,outeriter*9 + inneriter*76);
395 }
396 /*
397  * Gromacs nonbonded kernel:   nb_kernel_ElecCSTab_VdwCSTab_GeomP1P1_F_sparc64_hpc_ace_double
398  * Electrostatics interaction: CubicSplineTable
399  * VdW interaction:            CubicSplineTable
400  * Geometry:                   Particle-Particle
401  * Calculate force/pot:        Force
402  */
403 void
404 nb_kernel_ElecCSTab_VdwCSTab_GeomP1P1_F_sparc64_hpc_ace_double
405                     (t_nblist                    * gmx_restrict       nlist,
406                      rvec                        * gmx_restrict          xx,
407                      rvec                        * gmx_restrict          ff,
408                      t_forcerec                  * gmx_restrict          fr,
409                      t_mdatoms                   * gmx_restrict     mdatoms,
410                      nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
411                      t_nrnb                      * gmx_restrict        nrnb)
412 {
413     /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
414      * just 0 for non-waters.
415      * Suffixes A,B refer to j loop unrolling done with double precision SIMD, e.g. for the two different
416      * jnr indices corresponding to data put in the four positions in the SIMD register.
417      */
418     int              i_shift_offset,i_coord_offset,outeriter,inneriter;
419     int              j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
420     int              jnrA,jnrB;
421     int              j_coord_offsetA,j_coord_offsetB;
422     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
423     real             rcutoff_scalar;
424     real             *shiftvec,*fshift,*x,*f;
425     _fjsp_v2r8       tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
426     int              vdwioffset0;
427     _fjsp_v2r8       ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
428     int              vdwjidx0A,vdwjidx0B;
429     _fjsp_v2r8       jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
430     _fjsp_v2r8       dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
431     _fjsp_v2r8       velec,felec,velecsum,facel,crf,krf,krf2;
432     real             *charge;
433     int              nvdwtype;
434     _fjsp_v2r8       rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
435     int              *vdwtype;
436     real             *vdwparam;
437     _fjsp_v2r8       one_sixth   = gmx_fjsp_set1_v2r8(1.0/6.0);
438     _fjsp_v2r8       one_twelfth = gmx_fjsp_set1_v2r8(1.0/12.0);
439     _fjsp_v2r8       rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF,twovfeps;
440     real             *vftab;
441     _fjsp_v2r8       itab_tmp;
442     _fjsp_v2r8       dummy_mask,cutoff_mask;
443     _fjsp_v2r8       one     = gmx_fjsp_set1_v2r8(1.0);
444     _fjsp_v2r8       two     = gmx_fjsp_set1_v2r8(2.0);
445     union { _fjsp_v2r8 simd; long long int i[2]; } vfconv,gbconv,ewconv;
446
447     x                = xx[0];
448     f                = ff[0];
449
450     nri              = nlist->nri;
451     iinr             = nlist->iinr;
452     jindex           = nlist->jindex;
453     jjnr             = nlist->jjnr;
454     shiftidx         = nlist->shift;
455     gid              = nlist->gid;
456     shiftvec         = fr->shift_vec[0];
457     fshift           = fr->fshift[0];
458     facel            = gmx_fjsp_set1_v2r8(fr->epsfac);
459     charge           = mdatoms->chargeA;
460     nvdwtype         = fr->ntype;
461     vdwparam         = fr->nbfp;
462     vdwtype          = mdatoms->typeA;
463
464     vftab            = kernel_data->table_elec_vdw->data;
465     vftabscale       = gmx_fjsp_set1_v2r8(kernel_data->table_elec_vdw->scale);
466
467     /* Avoid stupid compiler warnings */
468     jnrA = jnrB = 0;
469     j_coord_offsetA = 0;
470     j_coord_offsetB = 0;
471
472     outeriter        = 0;
473     inneriter        = 0;
474
475     /* Start outer loop over neighborlists */
476     for(iidx=0; iidx<nri; iidx++)
477     {
478         /* Load shift vector for this list */
479         i_shift_offset   = DIM*shiftidx[iidx];
480
481         /* Load limits for loop over neighbors */
482         j_index_start    = jindex[iidx];
483         j_index_end      = jindex[iidx+1];
484
485         /* Get outer coordinate index */
486         inr              = iinr[iidx];
487         i_coord_offset   = DIM*inr;
488
489         /* Load i particle coords and add shift vector */
490         gmx_fjsp_load_shift_and_1rvec_broadcast_v2r8(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
491
492         fix0             = _fjsp_setzero_v2r8();
493         fiy0             = _fjsp_setzero_v2r8();
494         fiz0             = _fjsp_setzero_v2r8();
495
496         /* Load parameters for i particles */
497         iq0              = _fjsp_mul_v2r8(facel,gmx_fjsp_load1_v2r8(charge+inr+0));
498         vdwioffset0      = 2*nvdwtype*vdwtype[inr+0];
499
500         /* Start inner kernel loop */
501         for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
502         {
503
504             /* Get j neighbor index, and coordinate index */
505             jnrA             = jjnr[jidx];
506             jnrB             = jjnr[jidx+1];
507             j_coord_offsetA  = DIM*jnrA;
508             j_coord_offsetB  = DIM*jnrB;
509
510             /* load j atom coordinates */
511             gmx_fjsp_load_1rvec_2ptr_swizzle_v2r8(x+j_coord_offsetA,x+j_coord_offsetB,
512                                               &jx0,&jy0,&jz0);
513
514             /* Calculate displacement vector */
515             dx00             = _fjsp_sub_v2r8(ix0,jx0);
516             dy00             = _fjsp_sub_v2r8(iy0,jy0);
517             dz00             = _fjsp_sub_v2r8(iz0,jz0);
518
519             /* Calculate squared distance and things based on it */
520             rsq00            = gmx_fjsp_calc_rsq_v2r8(dx00,dy00,dz00);
521
522             rinv00           = gmx_fjsp_invsqrt_v2r8(rsq00);
523
524             /* Load parameters for j particles */
525             jq0              = gmx_fjsp_load_2real_swizzle_v2r8(charge+jnrA+0,charge+jnrB+0);
526             vdwjidx0A        = 2*vdwtype[jnrA+0];
527             vdwjidx0B        = 2*vdwtype[jnrB+0];
528
529             /**************************
530              * CALCULATE INTERACTIONS *
531              **************************/
532
533             r00              = _fjsp_mul_v2r8(rsq00,rinv00);
534
535             /* Compute parameters for interactions between i and j atoms */
536             qq00             = _fjsp_mul_v2r8(iq0,jq0);
537             gmx_fjsp_load_2pair_swizzle_v2r8(vdwparam+vdwioffset0+vdwjidx0A,
538                                          vdwparam+vdwioffset0+vdwjidx0B,&c6_00,&c12_00);
539
540             /* Calculate table index by multiplying r with table scale and truncate to integer */
541             rt               = _fjsp_mul_v2r8(r00,vftabscale);
542             itab_tmp         = _fjsp_dtox_v2r8(rt);
543             vfeps            = _fjsp_sub_v2r8(rt, _fjsp_xtod_v2r8(itab_tmp));
544             twovfeps         = _fjsp_add_v2r8(vfeps,vfeps);
545             _fjsp_store_v2r8(&vfconv.simd,itab_tmp);
546
547             vfconv.i[0]     *= 12;
548             vfconv.i[1]     *= 12;
549
550             /* CUBIC SPLINE TABLE ELECTROSTATICS */
551             Y                = _fjsp_load_v2r8( vftab + vfconv.i[0] );
552             F                = _fjsp_load_v2r8( vftab + vfconv.i[1] );
553             GMX_FJSP_TRANSPOSE2_V2R8(Y,F);
554             G                = _fjsp_load_v2r8( vftab + vfconv.i[0] +2);
555             H                = _fjsp_load_v2r8( vftab + vfconv.i[1] +2);
556             GMX_FJSP_TRANSPOSE2_V2R8(G,H);
557             Fp               = _fjsp_madd_v2r8(vfeps,_fjsp_madd_v2r8(vfeps,H,G),F);
558             FF               = _fjsp_madd_v2r8(_fjsp_madd_v2r8(twovfeps,H,G),vfeps,Fp);
559             felec            = _fjsp_neg_v2r8(_fjsp_mul_v2r8(_fjsp_mul_v2r8(qq00,FF),_fjsp_mul_v2r8(vftabscale,rinv00)));
560
561             /* CUBIC SPLINE TABLE DISPERSION */
562             vfconv.i[0]       += 4;
563             vfconv.i[1]       += 4;
564             Y                = _fjsp_load_v2r8( vftab + vfconv.i[0] );
565             F                = _fjsp_load_v2r8( vftab + vfconv.i[1] );
566             GMX_FJSP_TRANSPOSE2_V2R8(Y,F);
567             G                = _fjsp_load_v2r8( vftab + vfconv.i[0] + 2 );
568             H                = _fjsp_load_v2r8( vftab + vfconv.i[1] + 2 );
569             GMX_FJSP_TRANSPOSE2_V2R8(G,H);
570             Fp               = _fjsp_madd_v2r8(vfeps,_fjsp_madd_v2r8(H,vfeps,G),F);
571             FF               = _fjsp_madd_v2r8(vfeps,_fjsp_madd_v2r8(twovfeps,H,G),Fp);
572             fvdw6            = _fjsp_mul_v2r8(c6_00,FF);
573
574             /* CUBIC SPLINE TABLE REPULSION */
575             Y                = _fjsp_load_v2r8( vftab + vfconv.i[0] + 4 );
576             F                = _fjsp_load_v2r8( vftab + vfconv.i[1] + 4 );
577             GMX_FJSP_TRANSPOSE2_V2R8(Y,F);
578             G                = _fjsp_load_v2r8( vftab + vfconv.i[0] + 6 );
579             H                = _fjsp_load_v2r8( vftab + vfconv.i[1] + 6 );
580             GMX_FJSP_TRANSPOSE2_V2R8(G,H);
581             Fp               = _fjsp_madd_v2r8(vfeps,_fjsp_madd_v2r8(H,vfeps,G),F);
582             FF               = _fjsp_madd_v2r8(vfeps,_fjsp_madd_v2r8(twovfeps,H,G),Fp);
583             fvdw12           = _fjsp_mul_v2r8(c12_00,FF);
584             fvdw             = _fjsp_neg_v2r8(_fjsp_mul_v2r8(_fjsp_add_v2r8(fvdw6,fvdw12),_fjsp_mul_v2r8(vftabscale,rinv00)));
585
586             fscal            = _fjsp_add_v2r8(felec,fvdw);
587
588             /* Update vectorial force */
589             fix0             = _fjsp_madd_v2r8(dx00,fscal,fix0);
590             fiy0             = _fjsp_madd_v2r8(dy00,fscal,fiy0);
591             fiz0             = _fjsp_madd_v2r8(dz00,fscal,fiz0);
592             
593             gmx_fjsp_decrement_fma_1rvec_2ptr_swizzle_v2r8(f+j_coord_offsetA,f+j_coord_offsetB,fscal,dx00,dy00,dz00);
594
595             /* Inner loop uses 64 flops */
596         }
597
598         if(jidx<j_index_end)
599         {
600
601             jnrA             = jjnr[jidx];
602             j_coord_offsetA  = DIM*jnrA;
603
604             /* load j atom coordinates */
605             gmx_fjsp_load_1rvec_1ptr_swizzle_v2r8(x+j_coord_offsetA,
606                                               &jx0,&jy0,&jz0);
607
608             /* Calculate displacement vector */
609             dx00             = _fjsp_sub_v2r8(ix0,jx0);
610             dy00             = _fjsp_sub_v2r8(iy0,jy0);
611             dz00             = _fjsp_sub_v2r8(iz0,jz0);
612
613             /* Calculate squared distance and things based on it */
614             rsq00            = gmx_fjsp_calc_rsq_v2r8(dx00,dy00,dz00);
615
616             rinv00           = gmx_fjsp_invsqrt_v2r8(rsq00);
617
618             /* Load parameters for j particles */
619             jq0              = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(),charge+jnrA+0);
620             vdwjidx0A        = 2*vdwtype[jnrA+0];
621
622             /**************************
623              * CALCULATE INTERACTIONS *
624              **************************/
625
626             r00              = _fjsp_mul_v2r8(rsq00,rinv00);
627
628             /* Compute parameters for interactions between i and j atoms */
629             qq00             = _fjsp_mul_v2r8(iq0,jq0);
630             gmx_fjsp_load_1pair_swizzle_v2r8(vdwparam+vdwioffset0+vdwjidx0A,&c6_00,&c12_00);
631
632             /* Calculate table index by multiplying r with table scale and truncate to integer */
633             rt               = _fjsp_mul_v2r8(r00,vftabscale);
634             itab_tmp         = _fjsp_dtox_v2r8(rt);
635             vfeps            = _fjsp_sub_v2r8(rt, _fjsp_xtod_v2r8(itab_tmp));
636             twovfeps         = _fjsp_add_v2r8(vfeps,vfeps);
637             _fjsp_store_v2r8(&vfconv.simd,itab_tmp);
638
639             vfconv.i[0]     *= 12;
640             vfconv.i[1]     *= 12;
641
642             /* CUBIC SPLINE TABLE ELECTROSTATICS */
643             Y                = _fjsp_load_v2r8( vftab + vfconv.i[0] );
644             F                = _fjsp_setzero_v2r8();
645             GMX_FJSP_TRANSPOSE2_V2R8(Y,F);
646             G                = _fjsp_load_v2r8( vftab + vfconv.i[0] +2);
647             H                = _fjsp_setzero_v2r8();
648             GMX_FJSP_TRANSPOSE2_V2R8(G,H);
649             Fp               = _fjsp_madd_v2r8(vfeps,_fjsp_madd_v2r8(vfeps,H,G),F);
650             FF               = _fjsp_madd_v2r8(_fjsp_madd_v2r8(twovfeps,H,G),vfeps,Fp);
651             felec            = _fjsp_neg_v2r8(_fjsp_mul_v2r8(_fjsp_mul_v2r8(qq00,FF),_fjsp_mul_v2r8(vftabscale,rinv00)));
652
653             /* CUBIC SPLINE TABLE DISPERSION */
654             vfconv.i[0]       += 4;
655             vfconv.i[1]       += 4;
656             Y                = _fjsp_load_v2r8( vftab + vfconv.i[0] );
657             F                = _fjsp_setzero_v2r8();
658             GMX_FJSP_TRANSPOSE2_V2R8(Y,F);
659             G                = _fjsp_load_v2r8( vftab + vfconv.i[0] + 2 );
660             H                = _fjsp_setzero_v2r8();
661             GMX_FJSP_TRANSPOSE2_V2R8(G,H);
662             Fp               = _fjsp_madd_v2r8(vfeps,_fjsp_madd_v2r8(H,vfeps,G),F);
663             FF               = _fjsp_madd_v2r8(vfeps,_fjsp_madd_v2r8(twovfeps,H,G),Fp);
664             fvdw6            = _fjsp_mul_v2r8(c6_00,FF);
665
666             /* CUBIC SPLINE TABLE REPULSION */
667             Y                = _fjsp_load_v2r8( vftab + vfconv.i[0] + 4 );
668             F                = _fjsp_setzero_v2r8();
669             GMX_FJSP_TRANSPOSE2_V2R8(Y,F);
670             G                = _fjsp_load_v2r8( vftab + vfconv.i[0] + 6 );
671             H                = _fjsp_setzero_v2r8();
672             GMX_FJSP_TRANSPOSE2_V2R8(G,H);
673             Fp               = _fjsp_madd_v2r8(vfeps,_fjsp_madd_v2r8(H,vfeps,G),F);
674             FF               = _fjsp_madd_v2r8(vfeps,_fjsp_madd_v2r8(twovfeps,H,G),Fp);
675             fvdw12           = _fjsp_mul_v2r8(c12_00,FF);
676             fvdw             = _fjsp_neg_v2r8(_fjsp_mul_v2r8(_fjsp_add_v2r8(fvdw6,fvdw12),_fjsp_mul_v2r8(vftabscale,rinv00)));
677
678             fscal            = _fjsp_add_v2r8(felec,fvdw);
679
680             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
681
682             /* Update vectorial force */
683             fix0             = _fjsp_madd_v2r8(dx00,fscal,fix0);
684             fiy0             = _fjsp_madd_v2r8(dy00,fscal,fiy0);
685             fiz0             = _fjsp_madd_v2r8(dz00,fscal,fiz0);
686             
687             gmx_fjsp_decrement_fma_1rvec_1ptr_swizzle_v2r8(f+j_coord_offsetA,fscal,dx00,dy00,dz00);
688
689             /* Inner loop uses 64 flops */
690         }
691
692         /* End of innermost loop */
693
694         gmx_fjsp_update_iforce_1atom_swizzle_v2r8(fix0,fiy0,fiz0,
695                                               f+i_coord_offset,fshift+i_shift_offset);
696
697         /* Increment number of inner iterations */
698         inneriter                  += j_index_end - j_index_start;
699
700         /* Outer loop uses 7 flops */
701     }
702
703     /* Increment number of outer iterations */
704     outeriter        += nri;
705
706     /* Update outer/inner flops */
707
708     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_F,outeriter*7 + inneriter*64);
709 }