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