Use full path for legacyheaders
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_sparc64_hpc_ace_double / nb_kernel_ElecEw_VdwLJEw_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 "gromacs/legacyheaders/types/simple.h"
44 #include "gromacs/math/vec.h"
45 #include "gromacs/legacyheaders/nrnb.h"
46
47 #include "kernelutil_sparc64_hpc_ace_double.h"
48
49 /*
50  * Gromacs nonbonded kernel:   nb_kernel_ElecEw_VdwLJEw_GeomP1P1_VF_sparc64_hpc_ace_double
51  * Electrostatics interaction: Ewald
52  * VdW interaction:            LJEwald
53  * Geometry:                   Particle-Particle
54  * Calculate force/pot:        PotentialAndForce
55  */
56 void
57 nb_kernel_ElecEw_VdwLJEw_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           c6grid_00;
93     real                 *vdwgridparam;
94     _fjsp_v2r8           ewclj,ewclj2,ewclj6,ewcljrsq,poly,exponent,f6A,f6B,sh_lj_ewald;
95     _fjsp_v2r8           one_half = gmx_fjsp_set1_v2r8(0.5);
96     _fjsp_v2r8           minus_one = gmx_fjsp_set1_v2r8(-1.0);
97     _fjsp_v2r8       ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
98     real             *ewtab;
99     _fjsp_v2r8       itab_tmp;
100     _fjsp_v2r8       dummy_mask,cutoff_mask;
101     _fjsp_v2r8       one     = gmx_fjsp_set1_v2r8(1.0);
102     _fjsp_v2r8       two     = gmx_fjsp_set1_v2r8(2.0);
103     union { _fjsp_v2r8 simd; long long int i[2]; } vfconv,gbconv,ewconv;
104
105     x                = xx[0];
106     f                = ff[0];
107
108     nri              = nlist->nri;
109     iinr             = nlist->iinr;
110     jindex           = nlist->jindex;
111     jjnr             = nlist->jjnr;
112     shiftidx         = nlist->shift;
113     gid              = nlist->gid;
114     shiftvec         = fr->shift_vec[0];
115     fshift           = fr->fshift[0];
116     facel            = gmx_fjsp_set1_v2r8(fr->epsfac);
117     charge           = mdatoms->chargeA;
118     nvdwtype         = fr->ntype;
119     vdwparam         = fr->nbfp;
120     vdwtype          = mdatoms->typeA;
121     vdwgridparam     = fr->ljpme_c6grid;
122     sh_lj_ewald      = gmx_fjsp_set1_v2r8(fr->ic->sh_lj_ewald);
123     ewclj            = gmx_fjsp_set1_v2r8(fr->ewaldcoeff_lj);
124     ewclj2           = _fjsp_mul_v2r8(minus_one,_fjsp_mul_v2r8(ewclj,ewclj));
125
126     sh_ewald         = gmx_fjsp_set1_v2r8(fr->ic->sh_ewald);
127     ewtab            = fr->ic->tabq_coul_FDV0;
128     ewtabscale       = gmx_fjsp_set1_v2r8(fr->ic->tabq_scale);
129     ewtabhalfspace   = gmx_fjsp_set1_v2r8(0.5/fr->ic->tabq_scale);
130
131     /* Avoid stupid compiler warnings */
132     jnrA = jnrB = 0;
133     j_coord_offsetA = 0;
134     j_coord_offsetB = 0;
135
136     outeriter        = 0;
137     inneriter        = 0;
138
139     /* Start outer loop over neighborlists */
140     for(iidx=0; iidx<nri; iidx++)
141     {
142         /* Load shift vector for this list */
143         i_shift_offset   = DIM*shiftidx[iidx];
144
145         /* Load limits for loop over neighbors */
146         j_index_start    = jindex[iidx];
147         j_index_end      = jindex[iidx+1];
148
149         /* Get outer coordinate index */
150         inr              = iinr[iidx];
151         i_coord_offset   = DIM*inr;
152
153         /* Load i particle coords and add shift vector */
154         gmx_fjsp_load_shift_and_1rvec_broadcast_v2r8(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
155
156         fix0             = _fjsp_setzero_v2r8();
157         fiy0             = _fjsp_setzero_v2r8();
158         fiz0             = _fjsp_setzero_v2r8();
159
160         /* Load parameters for i particles */
161         iq0              = _fjsp_mul_v2r8(facel,gmx_fjsp_load1_v2r8(charge+inr+0));
162         vdwioffset0      = 2*nvdwtype*vdwtype[inr+0];
163
164         /* Reset potential sums */
165         velecsum         = _fjsp_setzero_v2r8();
166         vvdwsum          = _fjsp_setzero_v2r8();
167
168         /* Start inner kernel loop */
169         for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
170         {
171
172             /* Get j neighbor index, and coordinate index */
173             jnrA             = jjnr[jidx];
174             jnrB             = jjnr[jidx+1];
175             j_coord_offsetA  = DIM*jnrA;
176             j_coord_offsetB  = DIM*jnrB;
177
178             /* load j atom coordinates */
179             gmx_fjsp_load_1rvec_2ptr_swizzle_v2r8(x+j_coord_offsetA,x+j_coord_offsetB,
180                                               &jx0,&jy0,&jz0);
181
182             /* Calculate displacement vector */
183             dx00             = _fjsp_sub_v2r8(ix0,jx0);
184             dy00             = _fjsp_sub_v2r8(iy0,jy0);
185             dz00             = _fjsp_sub_v2r8(iz0,jz0);
186
187             /* Calculate squared distance and things based on it */
188             rsq00            = gmx_fjsp_calc_rsq_v2r8(dx00,dy00,dz00);
189
190             rinv00           = gmx_fjsp_invsqrt_v2r8(rsq00);
191
192             rinvsq00         = _fjsp_mul_v2r8(rinv00,rinv00);
193
194             /* Load parameters for j particles */
195             jq0              = gmx_fjsp_load_2real_swizzle_v2r8(charge+jnrA+0,charge+jnrB+0);
196             vdwjidx0A        = 2*vdwtype[jnrA+0];
197             vdwjidx0B        = 2*vdwtype[jnrB+0];
198
199             /**************************
200              * CALCULATE INTERACTIONS *
201              **************************/
202
203             r00              = _fjsp_mul_v2r8(rsq00,rinv00);
204
205             /* Compute parameters for interactions between i and j atoms */
206             qq00             = _fjsp_mul_v2r8(iq0,jq0);
207             gmx_fjsp_load_2pair_swizzle_v2r8(vdwparam+vdwioffset0+vdwjidx0A,
208                                          vdwparam+vdwioffset0+vdwjidx0B,&c6_00,&c12_00);
209
210             c6grid_00       = gmx_fjsp_load_2real_swizzle_v2r8(vdwgridparam+vdwioffset0+vdwjidx0A,
211                                                                    vdwgridparam+vdwioffset0+vdwjidx0B);
212
213             /* EWALD ELECTROSTATICS */
214
215             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
216             ewrt             = _fjsp_mul_v2r8(r00,ewtabscale);
217             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
218             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
219             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
220
221             ewtabF           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[0] );
222             ewtabD           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[1] );
223             GMX_FJSP_TRANSPOSE2_V2R8(ewtabF,ewtabD);
224             ewtabV           = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[0] +2);
225             ewtabFn          = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[1] +2);
226             GMX_FJSP_TRANSPOSE2_V2R8(ewtabV,ewtabFn);
227             felec            = _fjsp_madd_v2r8(eweps,ewtabD,ewtabF);
228             velec            = _fjsp_nmsub_v2r8(_fjsp_mul_v2r8(ewtabhalfspace,eweps) ,_fjsp_add_v2r8(ewtabF,felec), ewtabV);
229             velec            = _fjsp_mul_v2r8(qq00,_fjsp_sub_v2r8(rinv00,velec));
230             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq00,rinv00),_fjsp_sub_v2r8(rinvsq00,felec));
231
232             /* Analytical LJ-PME */
233             rinvsix          = _fjsp_mul_v2r8(_fjsp_mul_v2r8(rinvsq00,rinvsq00),rinvsq00);
234             ewcljrsq         = _fjsp_mul_v2r8(ewclj2,rsq00);
235             ewclj6           = _fjsp_mul_v2r8(ewclj2,_fjsp_mul_v2r8(ewclj2,ewclj2));
236             exponent         = gmx_simd_exp_d(-ewcljrsq);
237             /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
238             poly             = _fjsp_mul_v2r8(exponent,_fjsp_madd_v2r8(_fjsp_mul_v2r8(ewcljrsq,ewcljrsq),one_half,_fjsp_sub_v2r8(one,ewcljrsq)));
239             /* vvdw6 = [C6 - C6grid * (1-poly)]/r6 */
240             vvdw6            = _fjsp_mul_v2r8(_fjsp_madd_v2r8(-c6grid_00,_fjsp_sub_v2r8(one,poly),c6_00),rinvsix);
241             vvdw12           = _fjsp_mul_v2r8(c12_00,_fjsp_mul_v2r8(rinvsix,rinvsix));
242             vvdw             = _fjsp_msub_v2r8(vvdw12,one_twelfth,_fjsp_mul_v2r8(vvdw6,one_sixth));         
243             /* fvdw = vvdw12/r - (vvdw6/r + (C6grid * exponent * beta^6)/r) */
244             fvdw             = _fjsp_mul_v2r8(_fjsp_add_v2r8(vvdw12,_fjsp_msub_v2r8(_fjsp_mul_v2r8(c6grid_00,one_sixth),_fjsp_mul_v2r8(exponent,ewclj6),vvdw6)),rinvsq00);
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 68 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             c6grid_00       = gmx_fjsp_load_1real_swizzle_v2r8(vdwgridparam+vdwioffset0+vdwjidx0A);
299
300             /* EWALD ELECTROSTATICS */
301
302             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
303             ewrt             = _fjsp_mul_v2r8(r00,ewtabscale);
304             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
305             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
306             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
307
308             ewtabF           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[0] );
309             ewtabD           = _fjsp_setzero_v2r8();
310             GMX_FJSP_TRANSPOSE2_V2R8(ewtabF,ewtabD);
311             ewtabV           = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[0] +2);
312             ewtabFn          = _fjsp_setzero_v2r8();
313             GMX_FJSP_TRANSPOSE2_V2R8(ewtabV,ewtabFn);
314             felec            = _fjsp_madd_v2r8(eweps,ewtabD,ewtabF);
315             velec            = _fjsp_nmsub_v2r8(_fjsp_mul_v2r8(ewtabhalfspace,eweps) ,_fjsp_add_v2r8(ewtabF,felec), ewtabV);
316             velec            = _fjsp_mul_v2r8(qq00,_fjsp_sub_v2r8(rinv00,velec));
317             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq00,rinv00),_fjsp_sub_v2r8(rinvsq00,felec));
318
319             /* Analytical LJ-PME */
320             rinvsix          = _fjsp_mul_v2r8(_fjsp_mul_v2r8(rinvsq00,rinvsq00),rinvsq00);
321             ewcljrsq         = _fjsp_mul_v2r8(ewclj2,rsq00);
322             ewclj6           = _fjsp_mul_v2r8(ewclj2,_fjsp_mul_v2r8(ewclj2,ewclj2));
323             exponent         = gmx_simd_exp_d(-ewcljrsq);
324             /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
325             poly             = _fjsp_mul_v2r8(exponent,_fjsp_madd_v2r8(_fjsp_mul_v2r8(ewcljrsq,ewcljrsq),one_half,_fjsp_sub_v2r8(one,ewcljrsq)));
326             /* vvdw6 = [C6 - C6grid * (1-poly)]/r6 */
327             vvdw6            = _fjsp_mul_v2r8(_fjsp_madd_v2r8(-c6grid_00,_fjsp_sub_v2r8(one,poly),c6_00),rinvsix);
328             vvdw12           = _fjsp_mul_v2r8(c12_00,_fjsp_mul_v2r8(rinvsix,rinvsix));
329             vvdw             = _fjsp_msub_v2r8(vvdw12,one_twelfth,_fjsp_mul_v2r8(vvdw6,one_sixth));         
330             /* fvdw = vvdw12/r - (vvdw6/r + (C6grid * exponent * beta^6)/r) */
331             fvdw             = _fjsp_mul_v2r8(_fjsp_add_v2r8(vvdw12,_fjsp_msub_v2r8(_fjsp_mul_v2r8(c6grid_00,one_sixth),_fjsp_mul_v2r8(exponent,ewclj6),vvdw6)),rinvsq00);
332
333             /* Update potential sum for this i atom from the interaction with this j atom. */
334             velec            = _fjsp_unpacklo_v2r8(velec,_fjsp_setzero_v2r8());
335             velecsum         = _fjsp_add_v2r8(velecsum,velec);
336             vvdw             = _fjsp_unpacklo_v2r8(vvdw,_fjsp_setzero_v2r8());
337             vvdwsum          = _fjsp_add_v2r8(vvdwsum,vvdw);
338
339             fscal            = _fjsp_add_v2r8(felec,fvdw);
340
341             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
342
343             /* Update vectorial force */
344             fix0             = _fjsp_madd_v2r8(dx00,fscal,fix0);
345             fiy0             = _fjsp_madd_v2r8(dy00,fscal,fiy0);
346             fiz0             = _fjsp_madd_v2r8(dz00,fscal,fiz0);
347             
348             gmx_fjsp_decrement_fma_1rvec_1ptr_swizzle_v2r8(f+j_coord_offsetA,fscal,dx00,dy00,dz00);
349
350             /* Inner loop uses 68 flops */
351         }
352
353         /* End of innermost loop */
354
355         gmx_fjsp_update_iforce_1atom_swizzle_v2r8(fix0,fiy0,fiz0,
356                                               f+i_coord_offset,fshift+i_shift_offset);
357
358         ggid                        = gid[iidx];
359         /* Update potential energies */
360         gmx_fjsp_update_1pot_v2r8(velecsum,kernel_data->energygrp_elec+ggid);
361         gmx_fjsp_update_1pot_v2r8(vvdwsum,kernel_data->energygrp_vdw+ggid);
362
363         /* Increment number of inner iterations */
364         inneriter                  += j_index_end - j_index_start;
365
366         /* Outer loop uses 9 flops */
367     }
368
369     /* Increment number of outer iterations */
370     outeriter        += nri;
371
372     /* Update outer/inner flops */
373
374     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_VF,outeriter*9 + inneriter*68);
375 }
376 /*
377  * Gromacs nonbonded kernel:   nb_kernel_ElecEw_VdwLJEw_GeomP1P1_F_sparc64_hpc_ace_double
378  * Electrostatics interaction: Ewald
379  * VdW interaction:            LJEwald
380  * Geometry:                   Particle-Particle
381  * Calculate force/pot:        Force
382  */
383 void
384 nb_kernel_ElecEw_VdwLJEw_GeomP1P1_F_sparc64_hpc_ace_double
385                     (t_nblist                    * gmx_restrict       nlist,
386                      rvec                        * gmx_restrict          xx,
387                      rvec                        * gmx_restrict          ff,
388                      t_forcerec                  * gmx_restrict          fr,
389                      t_mdatoms                   * gmx_restrict     mdatoms,
390                      nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
391                      t_nrnb                      * gmx_restrict        nrnb)
392 {
393     /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
394      * just 0 for non-waters.
395      * Suffixes A,B refer to j loop unrolling done with double precision SIMD, e.g. for the two different
396      * jnr indices corresponding to data put in the four positions in the SIMD register.
397      */
398     int              i_shift_offset,i_coord_offset,outeriter,inneriter;
399     int              j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
400     int              jnrA,jnrB;
401     int              j_coord_offsetA,j_coord_offsetB;
402     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
403     real             rcutoff_scalar;
404     real             *shiftvec,*fshift,*x,*f;
405     _fjsp_v2r8       tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
406     int              vdwioffset0;
407     _fjsp_v2r8       ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
408     int              vdwjidx0A,vdwjidx0B;
409     _fjsp_v2r8       jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
410     _fjsp_v2r8       dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
411     _fjsp_v2r8       velec,felec,velecsum,facel,crf,krf,krf2;
412     real             *charge;
413     int              nvdwtype;
414     _fjsp_v2r8       rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
415     int              *vdwtype;
416     real             *vdwparam;
417     _fjsp_v2r8       one_sixth   = gmx_fjsp_set1_v2r8(1.0/6.0);
418     _fjsp_v2r8       one_twelfth = gmx_fjsp_set1_v2r8(1.0/12.0);
419     _fjsp_v2r8           c6grid_00;
420     real                 *vdwgridparam;
421     _fjsp_v2r8           ewclj,ewclj2,ewclj6,ewcljrsq,poly,exponent,f6A,f6B,sh_lj_ewald;
422     _fjsp_v2r8           one_half = gmx_fjsp_set1_v2r8(0.5);
423     _fjsp_v2r8           minus_one = gmx_fjsp_set1_v2r8(-1.0);
424     _fjsp_v2r8       ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
425     real             *ewtab;
426     _fjsp_v2r8       itab_tmp;
427     _fjsp_v2r8       dummy_mask,cutoff_mask;
428     _fjsp_v2r8       one     = gmx_fjsp_set1_v2r8(1.0);
429     _fjsp_v2r8       two     = gmx_fjsp_set1_v2r8(2.0);
430     union { _fjsp_v2r8 simd; long long int i[2]; } vfconv,gbconv,ewconv;
431
432     x                = xx[0];
433     f                = ff[0];
434
435     nri              = nlist->nri;
436     iinr             = nlist->iinr;
437     jindex           = nlist->jindex;
438     jjnr             = nlist->jjnr;
439     shiftidx         = nlist->shift;
440     gid              = nlist->gid;
441     shiftvec         = fr->shift_vec[0];
442     fshift           = fr->fshift[0];
443     facel            = gmx_fjsp_set1_v2r8(fr->epsfac);
444     charge           = mdatoms->chargeA;
445     nvdwtype         = fr->ntype;
446     vdwparam         = fr->nbfp;
447     vdwtype          = mdatoms->typeA;
448     vdwgridparam     = fr->ljpme_c6grid;
449     sh_lj_ewald      = gmx_fjsp_set1_v2r8(fr->ic->sh_lj_ewald);
450     ewclj            = gmx_fjsp_set1_v2r8(fr->ewaldcoeff_lj);
451     ewclj2           = _fjsp_mul_v2r8(minus_one,_fjsp_mul_v2r8(ewclj,ewclj));
452
453     sh_ewald         = gmx_fjsp_set1_v2r8(fr->ic->sh_ewald);
454     ewtab            = fr->ic->tabq_coul_F;
455     ewtabscale       = gmx_fjsp_set1_v2r8(fr->ic->tabq_scale);
456     ewtabhalfspace   = gmx_fjsp_set1_v2r8(0.5/fr->ic->tabq_scale);
457
458     /* Avoid stupid compiler warnings */
459     jnrA = jnrB = 0;
460     j_coord_offsetA = 0;
461     j_coord_offsetB = 0;
462
463     outeriter        = 0;
464     inneriter        = 0;
465
466     /* Start outer loop over neighborlists */
467     for(iidx=0; iidx<nri; iidx++)
468     {
469         /* Load shift vector for this list */
470         i_shift_offset   = DIM*shiftidx[iidx];
471
472         /* Load limits for loop over neighbors */
473         j_index_start    = jindex[iidx];
474         j_index_end      = jindex[iidx+1];
475
476         /* Get outer coordinate index */
477         inr              = iinr[iidx];
478         i_coord_offset   = DIM*inr;
479
480         /* Load i particle coords and add shift vector */
481         gmx_fjsp_load_shift_and_1rvec_broadcast_v2r8(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
482
483         fix0             = _fjsp_setzero_v2r8();
484         fiy0             = _fjsp_setzero_v2r8();
485         fiz0             = _fjsp_setzero_v2r8();
486
487         /* Load parameters for i particles */
488         iq0              = _fjsp_mul_v2r8(facel,gmx_fjsp_load1_v2r8(charge+inr+0));
489         vdwioffset0      = 2*nvdwtype*vdwtype[inr+0];
490
491         /* Start inner kernel loop */
492         for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
493         {
494
495             /* Get j neighbor index, and coordinate index */
496             jnrA             = jjnr[jidx];
497             jnrB             = jjnr[jidx+1];
498             j_coord_offsetA  = DIM*jnrA;
499             j_coord_offsetB  = DIM*jnrB;
500
501             /* load j atom coordinates */
502             gmx_fjsp_load_1rvec_2ptr_swizzle_v2r8(x+j_coord_offsetA,x+j_coord_offsetB,
503                                               &jx0,&jy0,&jz0);
504
505             /* Calculate displacement vector */
506             dx00             = _fjsp_sub_v2r8(ix0,jx0);
507             dy00             = _fjsp_sub_v2r8(iy0,jy0);
508             dz00             = _fjsp_sub_v2r8(iz0,jz0);
509
510             /* Calculate squared distance and things based on it */
511             rsq00            = gmx_fjsp_calc_rsq_v2r8(dx00,dy00,dz00);
512
513             rinv00           = gmx_fjsp_invsqrt_v2r8(rsq00);
514
515             rinvsq00         = _fjsp_mul_v2r8(rinv00,rinv00);
516
517             /* Load parameters for j particles */
518             jq0              = gmx_fjsp_load_2real_swizzle_v2r8(charge+jnrA+0,charge+jnrB+0);
519             vdwjidx0A        = 2*vdwtype[jnrA+0];
520             vdwjidx0B        = 2*vdwtype[jnrB+0];
521
522             /**************************
523              * CALCULATE INTERACTIONS *
524              **************************/
525
526             r00              = _fjsp_mul_v2r8(rsq00,rinv00);
527
528             /* Compute parameters for interactions between i and j atoms */
529             qq00             = _fjsp_mul_v2r8(iq0,jq0);
530             gmx_fjsp_load_2pair_swizzle_v2r8(vdwparam+vdwioffset0+vdwjidx0A,
531                                          vdwparam+vdwioffset0+vdwjidx0B,&c6_00,&c12_00);
532
533             c6grid_00       = gmx_fjsp_load_2real_swizzle_v2r8(vdwgridparam+vdwioffset0+vdwjidx0A,
534                                                                    vdwgridparam+vdwioffset0+vdwjidx0B);
535
536             /* EWALD ELECTROSTATICS */
537
538             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
539             ewrt             = _fjsp_mul_v2r8(r00,ewtabscale);
540             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
541             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
542             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
543
544             gmx_fjsp_load_2pair_swizzle_v2r8(ewtab+ewconv.i[0],ewtab+ewconv.i[1],
545                                          &ewtabF,&ewtabFn);
546             felec            = _fjsp_madd_v2r8(eweps,ewtabFn,_fjsp_nmsub_v2r8(eweps,ewtabF,ewtabF));
547             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq00,rinv00),_fjsp_sub_v2r8(rinvsq00,felec));
548
549             /* Analytical LJ-PME */
550             rinvsix          = _fjsp_mul_v2r8(_fjsp_mul_v2r8(rinvsq00,rinvsq00),rinvsq00);
551             ewcljrsq         = _fjsp_mul_v2r8(ewclj2,rsq00);
552             ewclj6           = _fjsp_mul_v2r8(ewclj2,_fjsp_mul_v2r8(ewclj2,ewclj2));
553             exponent         = gmx_simd_exp_d(-ewcljrsq);
554             /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
555             poly             = _fjsp_mul_v2r8(exponent,_fjsp_madd_v2r8(_fjsp_mul_v2r8(ewcljrsq,ewcljrsq),one_half,_fjsp_sub_v2r8(one,ewcljrsq)));
556             /* f6A = 6 * C6grid * (1 - poly) */
557             f6A              = _fjsp_mul_v2r8(c6grid_00,_fjsp_msub_v2r8(one,poly));
558             /* f6B = C6grid * exponent * beta^6 */
559             f6B              = _fjsp_mul_v2r8(_fjsp_mul_v2r8(c6grid_00,one_sixth),_fjsp_mul_v2r8(exponent,ewclj6));
560             /* fvdw = 12*C12/r13 - ((6*C6 - f6A)/r6 + f6B)/r */
561             fvdw              = _fjsp_mul_v2r8(_fjsp_madd_v2r8(_fjsp_msub_v2r8(c12_00,rinvsix,_fjsp_sub_v2r8(c6_00,f6A)),rinvsix,f6B),rinvsq00);
562
563             fscal            = _fjsp_add_v2r8(felec,fvdw);
564
565             /* Update vectorial force */
566             fix0             = _fjsp_madd_v2r8(dx00,fscal,fix0);
567             fiy0             = _fjsp_madd_v2r8(dy00,fscal,fiy0);
568             fiz0             = _fjsp_madd_v2r8(dz00,fscal,fiz0);
569             
570             gmx_fjsp_decrement_fma_1rvec_2ptr_swizzle_v2r8(f+j_coord_offsetA,f+j_coord_offsetB,fscal,dx00,dy00,dz00);
571
572             /* Inner loop uses 61 flops */
573         }
574
575         if(jidx<j_index_end)
576         {
577
578             jnrA             = jjnr[jidx];
579             j_coord_offsetA  = DIM*jnrA;
580
581             /* load j atom coordinates */
582             gmx_fjsp_load_1rvec_1ptr_swizzle_v2r8(x+j_coord_offsetA,
583                                               &jx0,&jy0,&jz0);
584
585             /* Calculate displacement vector */
586             dx00             = _fjsp_sub_v2r8(ix0,jx0);
587             dy00             = _fjsp_sub_v2r8(iy0,jy0);
588             dz00             = _fjsp_sub_v2r8(iz0,jz0);
589
590             /* Calculate squared distance and things based on it */
591             rsq00            = gmx_fjsp_calc_rsq_v2r8(dx00,dy00,dz00);
592
593             rinv00           = gmx_fjsp_invsqrt_v2r8(rsq00);
594
595             rinvsq00         = _fjsp_mul_v2r8(rinv00,rinv00);
596
597             /* Load parameters for j particles */
598             jq0              = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(),charge+jnrA+0);
599             vdwjidx0A        = 2*vdwtype[jnrA+0];
600
601             /**************************
602              * CALCULATE INTERACTIONS *
603              **************************/
604
605             r00              = _fjsp_mul_v2r8(rsq00,rinv00);
606
607             /* Compute parameters for interactions between i and j atoms */
608             qq00             = _fjsp_mul_v2r8(iq0,jq0);
609             gmx_fjsp_load_1pair_swizzle_v2r8(vdwparam+vdwioffset0+vdwjidx0A,&c6_00,&c12_00);
610
611             c6grid_00       = gmx_fjsp_load_1real_swizzle_v2r8(vdwgridparam+vdwioffset0+vdwjidx0A);
612
613             /* EWALD ELECTROSTATICS */
614
615             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
616             ewrt             = _fjsp_mul_v2r8(r00,ewtabscale);
617             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
618             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
619             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
620
621             gmx_fjsp_load_1pair_swizzle_v2r8(ewtab+ewconv.i[0],&ewtabF,&ewtabFn);
622             felec            = _fjsp_madd_v2r8(eweps,ewtabFn,_fjsp_nmsub_v2r8(eweps,ewtabF,ewtabF));
623             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq00,rinv00),_fjsp_sub_v2r8(rinvsq00,felec));
624
625             /* Analytical LJ-PME */
626             rinvsix          = _fjsp_mul_v2r8(_fjsp_mul_v2r8(rinvsq00,rinvsq00),rinvsq00);
627             ewcljrsq         = _fjsp_mul_v2r8(ewclj2,rsq00);
628             ewclj6           = _fjsp_mul_v2r8(ewclj2,_fjsp_mul_v2r8(ewclj2,ewclj2));
629             exponent         = gmx_simd_exp_d(-ewcljrsq);
630             /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
631             poly             = _fjsp_mul_v2r8(exponent,_fjsp_madd_v2r8(_fjsp_mul_v2r8(ewcljrsq,ewcljrsq),one_half,_fjsp_sub_v2r8(one,ewcljrsq)));
632             /* f6A = 6 * C6grid * (1 - poly) */
633             f6A              = _fjsp_mul_v2r8(c6grid_00,_fjsp_msub_v2r8(one,poly));
634             /* f6B = C6grid * exponent * beta^6 */
635             f6B              = _fjsp_mul_v2r8(_fjsp_mul_v2r8(c6grid_00,one_sixth),_fjsp_mul_v2r8(exponent,ewclj6));
636             /* fvdw = 12*C12/r13 - ((6*C6 - f6A)/r6 + f6B)/r */
637             fvdw              = _fjsp_mul_v2r8(_fjsp_madd_v2r8(_fjsp_msub_v2r8(c12_00,rinvsix,_fjsp_sub_v2r8(c6_00,f6A)),rinvsix,f6B),rinvsq00);
638
639             fscal            = _fjsp_add_v2r8(felec,fvdw);
640
641             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
642
643             /* Update vectorial force */
644             fix0             = _fjsp_madd_v2r8(dx00,fscal,fix0);
645             fiy0             = _fjsp_madd_v2r8(dy00,fscal,fiy0);
646             fiz0             = _fjsp_madd_v2r8(dz00,fscal,fiz0);
647             
648             gmx_fjsp_decrement_fma_1rvec_1ptr_swizzle_v2r8(f+j_coord_offsetA,fscal,dx00,dy00,dz00);
649
650             /* Inner loop uses 61 flops */
651         }
652
653         /* End of innermost loop */
654
655         gmx_fjsp_update_iforce_1atom_swizzle_v2r8(fix0,fiy0,fiz0,
656                                               f+i_coord_offset,fshift+i_shift_offset);
657
658         /* Increment number of inner iterations */
659         inneriter                  += j_index_end - j_index_start;
660
661         /* Outer loop uses 7 flops */
662     }
663
664     /* Increment number of outer iterations */
665     outeriter        += nri;
666
667     /* Update outer/inner flops */
668
669     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_F,outeriter*7 + inneriter*61);
670 }