Remove all unnecessary HAVE_CONFIG_H
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_sparc64_hpc_ace_double / nb_kernel_ElecEw_VdwNone_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_ElecEw_VdwNone_GeomP1P1_VF_sparc64_hpc_ace_double
51  * Electrostatics interaction: Ewald
52  * VdW interaction:            None
53  * Geometry:                   Particle-Particle
54  * Calculate force/pot:        PotentialAndForce
55  */
56 void
57 nb_kernel_ElecEw_VdwNone_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     _fjsp_v2r8       ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
87     real             *ewtab;
88     _fjsp_v2r8       itab_tmp;
89     _fjsp_v2r8       dummy_mask,cutoff_mask;
90     _fjsp_v2r8       one     = gmx_fjsp_set1_v2r8(1.0);
91     _fjsp_v2r8       two     = gmx_fjsp_set1_v2r8(2.0);
92     union { _fjsp_v2r8 simd; long long int i[2]; } vfconv,gbconv,ewconv;
93
94     x                = xx[0];
95     f                = ff[0];
96
97     nri              = nlist->nri;
98     iinr             = nlist->iinr;
99     jindex           = nlist->jindex;
100     jjnr             = nlist->jjnr;
101     shiftidx         = nlist->shift;
102     gid              = nlist->gid;
103     shiftvec         = fr->shift_vec[0];
104     fshift           = fr->fshift[0];
105     facel            = gmx_fjsp_set1_v2r8(fr->epsfac);
106     charge           = mdatoms->chargeA;
107
108     sh_ewald         = gmx_fjsp_set1_v2r8(fr->ic->sh_ewald);
109     ewtab            = fr->ic->tabq_coul_FDV0;
110     ewtabscale       = gmx_fjsp_set1_v2r8(fr->ic->tabq_scale);
111     ewtabhalfspace   = gmx_fjsp_set1_v2r8(0.5/fr->ic->tabq_scale);
112
113     /* Avoid stupid compiler warnings */
114     jnrA = jnrB = 0;
115     j_coord_offsetA = 0;
116     j_coord_offsetB = 0;
117
118     outeriter        = 0;
119     inneriter        = 0;
120
121     /* Start outer loop over neighborlists */
122     for(iidx=0; iidx<nri; iidx++)
123     {
124         /* Load shift vector for this list */
125         i_shift_offset   = DIM*shiftidx[iidx];
126
127         /* Load limits for loop over neighbors */
128         j_index_start    = jindex[iidx];
129         j_index_end      = jindex[iidx+1];
130
131         /* Get outer coordinate index */
132         inr              = iinr[iidx];
133         i_coord_offset   = DIM*inr;
134
135         /* Load i particle coords and add shift vector */
136         gmx_fjsp_load_shift_and_1rvec_broadcast_v2r8(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
137
138         fix0             = _fjsp_setzero_v2r8();
139         fiy0             = _fjsp_setzero_v2r8();
140         fiz0             = _fjsp_setzero_v2r8();
141
142         /* Load parameters for i particles */
143         iq0              = _fjsp_mul_v2r8(facel,gmx_fjsp_load1_v2r8(charge+inr+0));
144
145         /* Reset potential sums */
146         velecsum         = _fjsp_setzero_v2r8();
147
148         /* Start inner kernel loop */
149         for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
150         {
151
152             /* Get j neighbor index, and coordinate index */
153             jnrA             = jjnr[jidx];
154             jnrB             = jjnr[jidx+1];
155             j_coord_offsetA  = DIM*jnrA;
156             j_coord_offsetB  = DIM*jnrB;
157
158             /* load j atom coordinates */
159             gmx_fjsp_load_1rvec_2ptr_swizzle_v2r8(x+j_coord_offsetA,x+j_coord_offsetB,
160                                               &jx0,&jy0,&jz0);
161
162             /* Calculate displacement vector */
163             dx00             = _fjsp_sub_v2r8(ix0,jx0);
164             dy00             = _fjsp_sub_v2r8(iy0,jy0);
165             dz00             = _fjsp_sub_v2r8(iz0,jz0);
166
167             /* Calculate squared distance and things based on it */
168             rsq00            = gmx_fjsp_calc_rsq_v2r8(dx00,dy00,dz00);
169
170             rinv00           = gmx_fjsp_invsqrt_v2r8(rsq00);
171
172             rinvsq00         = _fjsp_mul_v2r8(rinv00,rinv00);
173
174             /* Load parameters for j particles */
175             jq0              = gmx_fjsp_load_2real_swizzle_v2r8(charge+jnrA+0,charge+jnrB+0);
176
177             /**************************
178              * CALCULATE INTERACTIONS *
179              **************************/
180
181             r00              = _fjsp_mul_v2r8(rsq00,rinv00);
182
183             /* Compute parameters for interactions between i and j atoms */
184             qq00             = _fjsp_mul_v2r8(iq0,jq0);
185
186             /* EWALD ELECTROSTATICS */
187
188             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
189             ewrt             = _fjsp_mul_v2r8(r00,ewtabscale);
190             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
191             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
192             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
193
194             ewtabF           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[0] );
195             ewtabD           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[1] );
196             GMX_FJSP_TRANSPOSE2_V2R8(ewtabF,ewtabD);
197             ewtabV           = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[0] +2);
198             ewtabFn          = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[1] +2);
199             GMX_FJSP_TRANSPOSE2_V2R8(ewtabV,ewtabFn);
200             felec            = _fjsp_madd_v2r8(eweps,ewtabD,ewtabF);
201             velec            = _fjsp_nmsub_v2r8(_fjsp_mul_v2r8(ewtabhalfspace,eweps) ,_fjsp_add_v2r8(ewtabF,felec), ewtabV);
202             velec            = _fjsp_mul_v2r8(qq00,_fjsp_sub_v2r8(rinv00,velec));
203             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq00,rinv00),_fjsp_sub_v2r8(rinvsq00,felec));
204
205             /* Update potential sum for this i atom from the interaction with this j atom. */
206             velecsum         = _fjsp_add_v2r8(velecsum,velec);
207
208             fscal            = felec;
209
210             /* Update vectorial force */
211             fix0             = _fjsp_madd_v2r8(dx00,fscal,fix0);
212             fiy0             = _fjsp_madd_v2r8(dy00,fscal,fiy0);
213             fiz0             = _fjsp_madd_v2r8(dz00,fscal,fiz0);
214             
215             gmx_fjsp_decrement_fma_1rvec_2ptr_swizzle_v2r8(f+j_coord_offsetA,f+j_coord_offsetB,fscal,dx00,dy00,dz00);
216
217             /* Inner loop uses 44 flops */
218         }
219
220         if(jidx<j_index_end)
221         {
222
223             jnrA             = jjnr[jidx];
224             j_coord_offsetA  = DIM*jnrA;
225
226             /* load j atom coordinates */
227             gmx_fjsp_load_1rvec_1ptr_swizzle_v2r8(x+j_coord_offsetA,
228                                               &jx0,&jy0,&jz0);
229
230             /* Calculate displacement vector */
231             dx00             = _fjsp_sub_v2r8(ix0,jx0);
232             dy00             = _fjsp_sub_v2r8(iy0,jy0);
233             dz00             = _fjsp_sub_v2r8(iz0,jz0);
234
235             /* Calculate squared distance and things based on it */
236             rsq00            = gmx_fjsp_calc_rsq_v2r8(dx00,dy00,dz00);
237
238             rinv00           = gmx_fjsp_invsqrt_v2r8(rsq00);
239
240             rinvsq00         = _fjsp_mul_v2r8(rinv00,rinv00);
241
242             /* Load parameters for j particles */
243             jq0              = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(),charge+jnrA+0);
244
245             /**************************
246              * CALCULATE INTERACTIONS *
247              **************************/
248
249             r00              = _fjsp_mul_v2r8(rsq00,rinv00);
250
251             /* Compute parameters for interactions between i and j atoms */
252             qq00             = _fjsp_mul_v2r8(iq0,jq0);
253
254             /* EWALD ELECTROSTATICS */
255
256             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
257             ewrt             = _fjsp_mul_v2r8(r00,ewtabscale);
258             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
259             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
260             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
261
262             ewtabF           = _fjsp_load_v2r8( ewtab + 4*ewconv.i[0] );
263             ewtabD           = _fjsp_setzero_v2r8();
264             GMX_FJSP_TRANSPOSE2_V2R8(ewtabF,ewtabD);
265             ewtabV           = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(), ewtab + 4*ewconv.i[0] +2);
266             ewtabFn          = _fjsp_setzero_v2r8();
267             GMX_FJSP_TRANSPOSE2_V2R8(ewtabV,ewtabFn);
268             felec            = _fjsp_madd_v2r8(eweps,ewtabD,ewtabF);
269             velec            = _fjsp_nmsub_v2r8(_fjsp_mul_v2r8(ewtabhalfspace,eweps) ,_fjsp_add_v2r8(ewtabF,felec), ewtabV);
270             velec            = _fjsp_mul_v2r8(qq00,_fjsp_sub_v2r8(rinv00,velec));
271             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq00,rinv00),_fjsp_sub_v2r8(rinvsq00,felec));
272
273             /* Update potential sum for this i atom from the interaction with this j atom. */
274             velec            = _fjsp_unpacklo_v2r8(velec,_fjsp_setzero_v2r8());
275             velecsum         = _fjsp_add_v2r8(velecsum,velec);
276
277             fscal            = felec;
278
279             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
280
281             /* Update vectorial force */
282             fix0             = _fjsp_madd_v2r8(dx00,fscal,fix0);
283             fiy0             = _fjsp_madd_v2r8(dy00,fscal,fiy0);
284             fiz0             = _fjsp_madd_v2r8(dz00,fscal,fiz0);
285             
286             gmx_fjsp_decrement_fma_1rvec_1ptr_swizzle_v2r8(f+j_coord_offsetA,fscal,dx00,dy00,dz00);
287
288             /* Inner loop uses 44 flops */
289         }
290
291         /* End of innermost loop */
292
293         gmx_fjsp_update_iforce_1atom_swizzle_v2r8(fix0,fiy0,fiz0,
294                                               f+i_coord_offset,fshift+i_shift_offset);
295
296         ggid                        = gid[iidx];
297         /* Update potential energies */
298         gmx_fjsp_update_1pot_v2r8(velecsum,kernel_data->energygrp_elec+ggid);
299
300         /* Increment number of inner iterations */
301         inneriter                  += j_index_end - j_index_start;
302
303         /* Outer loop uses 8 flops */
304     }
305
306     /* Increment number of outer iterations */
307     outeriter        += nri;
308
309     /* Update outer/inner flops */
310
311     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VF,outeriter*8 + inneriter*44);
312 }
313 /*
314  * Gromacs nonbonded kernel:   nb_kernel_ElecEw_VdwNone_GeomP1P1_F_sparc64_hpc_ace_double
315  * Electrostatics interaction: Ewald
316  * VdW interaction:            None
317  * Geometry:                   Particle-Particle
318  * Calculate force/pot:        Force
319  */
320 void
321 nb_kernel_ElecEw_VdwNone_GeomP1P1_F_sparc64_hpc_ace_double
322                     (t_nblist                    * gmx_restrict       nlist,
323                      rvec                        * gmx_restrict          xx,
324                      rvec                        * gmx_restrict          ff,
325                      t_forcerec                  * gmx_restrict          fr,
326                      t_mdatoms                   * gmx_restrict     mdatoms,
327                      nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
328                      t_nrnb                      * gmx_restrict        nrnb)
329 {
330     /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
331      * just 0 for non-waters.
332      * Suffixes A,B refer to j loop unrolling done with double precision SIMD, e.g. for the two different
333      * jnr indices corresponding to data put in the four positions in the SIMD register.
334      */
335     int              i_shift_offset,i_coord_offset,outeriter,inneriter;
336     int              j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
337     int              jnrA,jnrB;
338     int              j_coord_offsetA,j_coord_offsetB;
339     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
340     real             rcutoff_scalar;
341     real             *shiftvec,*fshift,*x,*f;
342     _fjsp_v2r8       tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
343     int              vdwioffset0;
344     _fjsp_v2r8       ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
345     int              vdwjidx0A,vdwjidx0B;
346     _fjsp_v2r8       jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
347     _fjsp_v2r8       dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
348     _fjsp_v2r8       velec,felec,velecsum,facel,crf,krf,krf2;
349     real             *charge;
350     _fjsp_v2r8       ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
351     real             *ewtab;
352     _fjsp_v2r8       itab_tmp;
353     _fjsp_v2r8       dummy_mask,cutoff_mask;
354     _fjsp_v2r8       one     = gmx_fjsp_set1_v2r8(1.0);
355     _fjsp_v2r8       two     = gmx_fjsp_set1_v2r8(2.0);
356     union { _fjsp_v2r8 simd; long long int i[2]; } vfconv,gbconv,ewconv;
357
358     x                = xx[0];
359     f                = ff[0];
360
361     nri              = nlist->nri;
362     iinr             = nlist->iinr;
363     jindex           = nlist->jindex;
364     jjnr             = nlist->jjnr;
365     shiftidx         = nlist->shift;
366     gid              = nlist->gid;
367     shiftvec         = fr->shift_vec[0];
368     fshift           = fr->fshift[0];
369     facel            = gmx_fjsp_set1_v2r8(fr->epsfac);
370     charge           = mdatoms->chargeA;
371
372     sh_ewald         = gmx_fjsp_set1_v2r8(fr->ic->sh_ewald);
373     ewtab            = fr->ic->tabq_coul_F;
374     ewtabscale       = gmx_fjsp_set1_v2r8(fr->ic->tabq_scale);
375     ewtabhalfspace   = gmx_fjsp_set1_v2r8(0.5/fr->ic->tabq_scale);
376
377     /* Avoid stupid compiler warnings */
378     jnrA = jnrB = 0;
379     j_coord_offsetA = 0;
380     j_coord_offsetB = 0;
381
382     outeriter        = 0;
383     inneriter        = 0;
384
385     /* Start outer loop over neighborlists */
386     for(iidx=0; iidx<nri; iidx++)
387     {
388         /* Load shift vector for this list */
389         i_shift_offset   = DIM*shiftidx[iidx];
390
391         /* Load limits for loop over neighbors */
392         j_index_start    = jindex[iidx];
393         j_index_end      = jindex[iidx+1];
394
395         /* Get outer coordinate index */
396         inr              = iinr[iidx];
397         i_coord_offset   = DIM*inr;
398
399         /* Load i particle coords and add shift vector */
400         gmx_fjsp_load_shift_and_1rvec_broadcast_v2r8(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
401
402         fix0             = _fjsp_setzero_v2r8();
403         fiy0             = _fjsp_setzero_v2r8();
404         fiz0             = _fjsp_setzero_v2r8();
405
406         /* Load parameters for i particles */
407         iq0              = _fjsp_mul_v2r8(facel,gmx_fjsp_load1_v2r8(charge+inr+0));
408
409         /* Start inner kernel loop */
410         for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
411         {
412
413             /* Get j neighbor index, and coordinate index */
414             jnrA             = jjnr[jidx];
415             jnrB             = jjnr[jidx+1];
416             j_coord_offsetA  = DIM*jnrA;
417             j_coord_offsetB  = DIM*jnrB;
418
419             /* load j atom coordinates */
420             gmx_fjsp_load_1rvec_2ptr_swizzle_v2r8(x+j_coord_offsetA,x+j_coord_offsetB,
421                                               &jx0,&jy0,&jz0);
422
423             /* Calculate displacement vector */
424             dx00             = _fjsp_sub_v2r8(ix0,jx0);
425             dy00             = _fjsp_sub_v2r8(iy0,jy0);
426             dz00             = _fjsp_sub_v2r8(iz0,jz0);
427
428             /* Calculate squared distance and things based on it */
429             rsq00            = gmx_fjsp_calc_rsq_v2r8(dx00,dy00,dz00);
430
431             rinv00           = gmx_fjsp_invsqrt_v2r8(rsq00);
432
433             rinvsq00         = _fjsp_mul_v2r8(rinv00,rinv00);
434
435             /* Load parameters for j particles */
436             jq0              = gmx_fjsp_load_2real_swizzle_v2r8(charge+jnrA+0,charge+jnrB+0);
437
438             /**************************
439              * CALCULATE INTERACTIONS *
440              **************************/
441
442             r00              = _fjsp_mul_v2r8(rsq00,rinv00);
443
444             /* Compute parameters for interactions between i and j atoms */
445             qq00             = _fjsp_mul_v2r8(iq0,jq0);
446
447             /* EWALD ELECTROSTATICS */
448
449             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
450             ewrt             = _fjsp_mul_v2r8(r00,ewtabscale);
451             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
452             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
453             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
454
455             gmx_fjsp_load_2pair_swizzle_v2r8(ewtab+ewconv.i[0],ewtab+ewconv.i[1],
456                                          &ewtabF,&ewtabFn);
457             felec            = _fjsp_madd_v2r8(eweps,ewtabFn,_fjsp_nmsub_v2r8(eweps,ewtabF,ewtabF));
458             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq00,rinv00),_fjsp_sub_v2r8(rinvsq00,felec));
459
460             fscal            = felec;
461
462             /* Update vectorial force */
463             fix0             = _fjsp_madd_v2r8(dx00,fscal,fix0);
464             fiy0             = _fjsp_madd_v2r8(dy00,fscal,fiy0);
465             fiz0             = _fjsp_madd_v2r8(dz00,fscal,fiz0);
466             
467             gmx_fjsp_decrement_fma_1rvec_2ptr_swizzle_v2r8(f+j_coord_offsetA,f+j_coord_offsetB,fscal,dx00,dy00,dz00);
468
469             /* Inner loop uses 39 flops */
470         }
471
472         if(jidx<j_index_end)
473         {
474
475             jnrA             = jjnr[jidx];
476             j_coord_offsetA  = DIM*jnrA;
477
478             /* load j atom coordinates */
479             gmx_fjsp_load_1rvec_1ptr_swizzle_v2r8(x+j_coord_offsetA,
480                                               &jx0,&jy0,&jz0);
481
482             /* Calculate displacement vector */
483             dx00             = _fjsp_sub_v2r8(ix0,jx0);
484             dy00             = _fjsp_sub_v2r8(iy0,jy0);
485             dz00             = _fjsp_sub_v2r8(iz0,jz0);
486
487             /* Calculate squared distance and things based on it */
488             rsq00            = gmx_fjsp_calc_rsq_v2r8(dx00,dy00,dz00);
489
490             rinv00           = gmx_fjsp_invsqrt_v2r8(rsq00);
491
492             rinvsq00         = _fjsp_mul_v2r8(rinv00,rinv00);
493
494             /* Load parameters for j particles */
495             jq0              = _fjsp_loadl_v2r8(_fjsp_setzero_v2r8(),charge+jnrA+0);
496
497             /**************************
498              * CALCULATE INTERACTIONS *
499              **************************/
500
501             r00              = _fjsp_mul_v2r8(rsq00,rinv00);
502
503             /* Compute parameters for interactions between i and j atoms */
504             qq00             = _fjsp_mul_v2r8(iq0,jq0);
505
506             /* EWALD ELECTROSTATICS */
507
508             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
509             ewrt             = _fjsp_mul_v2r8(r00,ewtabscale);
510             itab_tmp         = _fjsp_dtox_v2r8(ewrt);
511             eweps            = _fjsp_sub_v2r8(ewrt,_fjsp_xtod_v2r8(itab_tmp));
512             _fjsp_store_v2r8(&ewconv.simd,itab_tmp);
513
514             gmx_fjsp_load_1pair_swizzle_v2r8(ewtab+ewconv.i[0],&ewtabF,&ewtabFn);
515             felec            = _fjsp_madd_v2r8(eweps,ewtabFn,_fjsp_nmsub_v2r8(eweps,ewtabF,ewtabF));
516             felec            = _fjsp_mul_v2r8(_fjsp_mul_v2r8(qq00,rinv00),_fjsp_sub_v2r8(rinvsq00,felec));
517
518             fscal            = felec;
519
520             fscal            = _fjsp_unpacklo_v2r8(fscal,_fjsp_setzero_v2r8());
521
522             /* Update vectorial force */
523             fix0             = _fjsp_madd_v2r8(dx00,fscal,fix0);
524             fiy0             = _fjsp_madd_v2r8(dy00,fscal,fiy0);
525             fiz0             = _fjsp_madd_v2r8(dz00,fscal,fiz0);
526             
527             gmx_fjsp_decrement_fma_1rvec_1ptr_swizzle_v2r8(f+j_coord_offsetA,fscal,dx00,dy00,dz00);
528
529             /* Inner loop uses 39 flops */
530         }
531
532         /* End of innermost loop */
533
534         gmx_fjsp_update_iforce_1atom_swizzle_v2r8(fix0,fiy0,fiz0,
535                                               f+i_coord_offset,fshift+i_shift_offset);
536
537         /* Increment number of inner iterations */
538         inneriter                  += j_index_end - j_index_start;
539
540         /* Outer loop uses 7 flops */
541     }
542
543     /* Increment number of outer iterations */
544     outeriter        += nri;
545
546     /* Update outer/inner flops */
547
548     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_F,outeriter*7 + inneriter*39);
549 }