Remove all unnecessary HAVE_CONFIG_H
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_sse2_double / nb_kernel_ElecCoul_VdwNone_GeomW4P1_sse2_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 sse2_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 "gromacs/simd/math_x86_sse2_double.h"
48 #include "kernelutil_x86_sse2_double.h"
49
50 /*
51  * Gromacs nonbonded kernel:   nb_kernel_ElecCoul_VdwNone_GeomW4P1_VF_sse2_double
52  * Electrostatics interaction: Coulomb
53  * VdW interaction:            None
54  * Geometry:                   Water4-Particle
55  * Calculate force/pot:        PotentialAndForce
56  */
57 void
58 nb_kernel_ElecCoul_VdwNone_GeomW4P1_VF_sse2_double
59                     (t_nblist                    * gmx_restrict       nlist,
60                      rvec                        * gmx_restrict          xx,
61                      rvec                        * gmx_restrict          ff,
62                      t_forcerec                  * gmx_restrict          fr,
63                      t_mdatoms                   * gmx_restrict     mdatoms,
64                      nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
65                      t_nrnb                      * gmx_restrict        nrnb)
66 {
67     /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
68      * just 0 for non-waters.
69      * Suffixes A,B refer to j loop unrolling done with SSE double precision, e.g. for the two different
70      * jnr indices corresponding to data put in the four positions in the SIMD register.
71      */
72     int              i_shift_offset,i_coord_offset,outeriter,inneriter;
73     int              j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
74     int              jnrA,jnrB;
75     int              j_coord_offsetA,j_coord_offsetB;
76     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
77     real             rcutoff_scalar;
78     real             *shiftvec,*fshift,*x,*f;
79     __m128d          tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
80     int              vdwioffset1;
81     __m128d          ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
82     int              vdwioffset2;
83     __m128d          ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
84     int              vdwioffset3;
85     __m128d          ix3,iy3,iz3,fix3,fiy3,fiz3,iq3,isai3;
86     int              vdwjidx0A,vdwjidx0B;
87     __m128d          jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
88     __m128d          dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
89     __m128d          dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
90     __m128d          dx30,dy30,dz30,rsq30,rinv30,rinvsq30,r30,qq30,c6_30,c12_30;
91     __m128d          velec,felec,velecsum,facel,crf,krf,krf2;
92     real             *charge;
93     __m128d          dummy_mask,cutoff_mask;
94     __m128d          signbit   = gmx_mm_castsi128_pd( _mm_set_epi32(0x80000000,0x00000000,0x80000000,0x00000000) );
95     __m128d          one     = _mm_set1_pd(1.0);
96     __m128d          two     = _mm_set1_pd(2.0);
97     x                = xx[0];
98     f                = ff[0];
99
100     nri              = nlist->nri;
101     iinr             = nlist->iinr;
102     jindex           = nlist->jindex;
103     jjnr             = nlist->jjnr;
104     shiftidx         = nlist->shift;
105     gid              = nlist->gid;
106     shiftvec         = fr->shift_vec[0];
107     fshift           = fr->fshift[0];
108     facel            = _mm_set1_pd(fr->epsfac);
109     charge           = mdatoms->chargeA;
110
111     /* Setup water-specific parameters */
112     inr              = nlist->iinr[0];
113     iq1              = _mm_mul_pd(facel,_mm_set1_pd(charge[inr+1]));
114     iq2              = _mm_mul_pd(facel,_mm_set1_pd(charge[inr+2]));
115     iq3              = _mm_mul_pd(facel,_mm_set1_pd(charge[inr+3]));
116
117     /* Avoid stupid compiler warnings */
118     jnrA = jnrB = 0;
119     j_coord_offsetA = 0;
120     j_coord_offsetB = 0;
121
122     outeriter        = 0;
123     inneriter        = 0;
124
125     /* Start outer loop over neighborlists */
126     for(iidx=0; iidx<nri; iidx++)
127     {
128         /* Load shift vector for this list */
129         i_shift_offset   = DIM*shiftidx[iidx];
130
131         /* Load limits for loop over neighbors */
132         j_index_start    = jindex[iidx];
133         j_index_end      = jindex[iidx+1];
134
135         /* Get outer coordinate index */
136         inr              = iinr[iidx];
137         i_coord_offset   = DIM*inr;
138
139         /* Load i particle coords and add shift vector */
140         gmx_mm_load_shift_and_3rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset+DIM,
141                                                  &ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
142
143         fix1             = _mm_setzero_pd();
144         fiy1             = _mm_setzero_pd();
145         fiz1             = _mm_setzero_pd();
146         fix2             = _mm_setzero_pd();
147         fiy2             = _mm_setzero_pd();
148         fiz2             = _mm_setzero_pd();
149         fix3             = _mm_setzero_pd();
150         fiy3             = _mm_setzero_pd();
151         fiz3             = _mm_setzero_pd();
152
153         /* Reset potential sums */
154         velecsum         = _mm_setzero_pd();
155
156         /* Start inner kernel loop */
157         for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
158         {
159
160             /* Get j neighbor index, and coordinate index */
161             jnrA             = jjnr[jidx];
162             jnrB             = jjnr[jidx+1];
163             j_coord_offsetA  = DIM*jnrA;
164             j_coord_offsetB  = DIM*jnrB;
165
166             /* load j atom coordinates */
167             gmx_mm_load_1rvec_2ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
168                                               &jx0,&jy0,&jz0);
169
170             /* Calculate displacement vector */
171             dx10             = _mm_sub_pd(ix1,jx0);
172             dy10             = _mm_sub_pd(iy1,jy0);
173             dz10             = _mm_sub_pd(iz1,jz0);
174             dx20             = _mm_sub_pd(ix2,jx0);
175             dy20             = _mm_sub_pd(iy2,jy0);
176             dz20             = _mm_sub_pd(iz2,jz0);
177             dx30             = _mm_sub_pd(ix3,jx0);
178             dy30             = _mm_sub_pd(iy3,jy0);
179             dz30             = _mm_sub_pd(iz3,jz0);
180
181             /* Calculate squared distance and things based on it */
182             rsq10            = gmx_mm_calc_rsq_pd(dx10,dy10,dz10);
183             rsq20            = gmx_mm_calc_rsq_pd(dx20,dy20,dz20);
184             rsq30            = gmx_mm_calc_rsq_pd(dx30,dy30,dz30);
185
186             rinv10           = gmx_mm_invsqrt_pd(rsq10);
187             rinv20           = gmx_mm_invsqrt_pd(rsq20);
188             rinv30           = gmx_mm_invsqrt_pd(rsq30);
189
190             rinvsq10         = _mm_mul_pd(rinv10,rinv10);
191             rinvsq20         = _mm_mul_pd(rinv20,rinv20);
192             rinvsq30         = _mm_mul_pd(rinv30,rinv30);
193
194             /* Load parameters for j particles */
195             jq0              = gmx_mm_load_2real_swizzle_pd(charge+jnrA+0,charge+jnrB+0);
196
197             fjx0             = _mm_setzero_pd();
198             fjy0             = _mm_setzero_pd();
199             fjz0             = _mm_setzero_pd();
200
201             /**************************
202              * CALCULATE INTERACTIONS *
203              **************************/
204
205             /* Compute parameters for interactions between i and j atoms */
206             qq10             = _mm_mul_pd(iq1,jq0);
207
208             /* COULOMB ELECTROSTATICS */
209             velec            = _mm_mul_pd(qq10,rinv10);
210             felec            = _mm_mul_pd(velec,rinvsq10);
211
212             /* Update potential sum for this i atom from the interaction with this j atom. */
213             velecsum         = _mm_add_pd(velecsum,velec);
214
215             fscal            = felec;
216
217             /* Calculate temporary vectorial force */
218             tx               = _mm_mul_pd(fscal,dx10);
219             ty               = _mm_mul_pd(fscal,dy10);
220             tz               = _mm_mul_pd(fscal,dz10);
221
222             /* Update vectorial force */
223             fix1             = _mm_add_pd(fix1,tx);
224             fiy1             = _mm_add_pd(fiy1,ty);
225             fiz1             = _mm_add_pd(fiz1,tz);
226
227             fjx0             = _mm_add_pd(fjx0,tx);
228             fjy0             = _mm_add_pd(fjy0,ty);
229             fjz0             = _mm_add_pd(fjz0,tz);
230
231             /**************************
232              * CALCULATE INTERACTIONS *
233              **************************/
234
235             /* Compute parameters for interactions between i and j atoms */
236             qq20             = _mm_mul_pd(iq2,jq0);
237
238             /* COULOMB ELECTROSTATICS */
239             velec            = _mm_mul_pd(qq20,rinv20);
240             felec            = _mm_mul_pd(velec,rinvsq20);
241
242             /* Update potential sum for this i atom from the interaction with this j atom. */
243             velecsum         = _mm_add_pd(velecsum,velec);
244
245             fscal            = felec;
246
247             /* Calculate temporary vectorial force */
248             tx               = _mm_mul_pd(fscal,dx20);
249             ty               = _mm_mul_pd(fscal,dy20);
250             tz               = _mm_mul_pd(fscal,dz20);
251
252             /* Update vectorial force */
253             fix2             = _mm_add_pd(fix2,tx);
254             fiy2             = _mm_add_pd(fiy2,ty);
255             fiz2             = _mm_add_pd(fiz2,tz);
256
257             fjx0             = _mm_add_pd(fjx0,tx);
258             fjy0             = _mm_add_pd(fjy0,ty);
259             fjz0             = _mm_add_pd(fjz0,tz);
260
261             /**************************
262              * CALCULATE INTERACTIONS *
263              **************************/
264
265             /* Compute parameters for interactions between i and j atoms */
266             qq30             = _mm_mul_pd(iq3,jq0);
267
268             /* COULOMB ELECTROSTATICS */
269             velec            = _mm_mul_pd(qq30,rinv30);
270             felec            = _mm_mul_pd(velec,rinvsq30);
271
272             /* Update potential sum for this i atom from the interaction with this j atom. */
273             velecsum         = _mm_add_pd(velecsum,velec);
274
275             fscal            = felec;
276
277             /* Calculate temporary vectorial force */
278             tx               = _mm_mul_pd(fscal,dx30);
279             ty               = _mm_mul_pd(fscal,dy30);
280             tz               = _mm_mul_pd(fscal,dz30);
281
282             /* Update vectorial force */
283             fix3             = _mm_add_pd(fix3,tx);
284             fiy3             = _mm_add_pd(fiy3,ty);
285             fiz3             = _mm_add_pd(fiz3,tz);
286
287             fjx0             = _mm_add_pd(fjx0,tx);
288             fjy0             = _mm_add_pd(fjy0,ty);
289             fjz0             = _mm_add_pd(fjz0,tz);
290
291             gmx_mm_decrement_1rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,fjx0,fjy0,fjz0);
292
293             /* Inner loop uses 87 flops */
294         }
295
296         if(jidx<j_index_end)
297         {
298
299             jnrA             = jjnr[jidx];
300             j_coord_offsetA  = DIM*jnrA;
301
302             /* load j atom coordinates */
303             gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
304                                               &jx0,&jy0,&jz0);
305
306             /* Calculate displacement vector */
307             dx10             = _mm_sub_pd(ix1,jx0);
308             dy10             = _mm_sub_pd(iy1,jy0);
309             dz10             = _mm_sub_pd(iz1,jz0);
310             dx20             = _mm_sub_pd(ix2,jx0);
311             dy20             = _mm_sub_pd(iy2,jy0);
312             dz20             = _mm_sub_pd(iz2,jz0);
313             dx30             = _mm_sub_pd(ix3,jx0);
314             dy30             = _mm_sub_pd(iy3,jy0);
315             dz30             = _mm_sub_pd(iz3,jz0);
316
317             /* Calculate squared distance and things based on it */
318             rsq10            = gmx_mm_calc_rsq_pd(dx10,dy10,dz10);
319             rsq20            = gmx_mm_calc_rsq_pd(dx20,dy20,dz20);
320             rsq30            = gmx_mm_calc_rsq_pd(dx30,dy30,dz30);
321
322             rinv10           = gmx_mm_invsqrt_pd(rsq10);
323             rinv20           = gmx_mm_invsqrt_pd(rsq20);
324             rinv30           = gmx_mm_invsqrt_pd(rsq30);
325
326             rinvsq10         = _mm_mul_pd(rinv10,rinv10);
327             rinvsq20         = _mm_mul_pd(rinv20,rinv20);
328             rinvsq30         = _mm_mul_pd(rinv30,rinv30);
329
330             /* Load parameters for j particles */
331             jq0              = _mm_load_sd(charge+jnrA+0);
332
333             fjx0             = _mm_setzero_pd();
334             fjy0             = _mm_setzero_pd();
335             fjz0             = _mm_setzero_pd();
336
337             /**************************
338              * CALCULATE INTERACTIONS *
339              **************************/
340
341             /* Compute parameters for interactions between i and j atoms */
342             qq10             = _mm_mul_pd(iq1,jq0);
343
344             /* COULOMB ELECTROSTATICS */
345             velec            = _mm_mul_pd(qq10,rinv10);
346             felec            = _mm_mul_pd(velec,rinvsq10);
347
348             /* Update potential sum for this i atom from the interaction with this j atom. */
349             velec            = _mm_unpacklo_pd(velec,_mm_setzero_pd());
350             velecsum         = _mm_add_pd(velecsum,velec);
351
352             fscal            = felec;
353
354             fscal            = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
355
356             /* Calculate temporary vectorial force */
357             tx               = _mm_mul_pd(fscal,dx10);
358             ty               = _mm_mul_pd(fscal,dy10);
359             tz               = _mm_mul_pd(fscal,dz10);
360
361             /* Update vectorial force */
362             fix1             = _mm_add_pd(fix1,tx);
363             fiy1             = _mm_add_pd(fiy1,ty);
364             fiz1             = _mm_add_pd(fiz1,tz);
365
366             fjx0             = _mm_add_pd(fjx0,tx);
367             fjy0             = _mm_add_pd(fjy0,ty);
368             fjz0             = _mm_add_pd(fjz0,tz);
369
370             /**************************
371              * CALCULATE INTERACTIONS *
372              **************************/
373
374             /* Compute parameters for interactions between i and j atoms */
375             qq20             = _mm_mul_pd(iq2,jq0);
376
377             /* COULOMB ELECTROSTATICS */
378             velec            = _mm_mul_pd(qq20,rinv20);
379             felec            = _mm_mul_pd(velec,rinvsq20);
380
381             /* Update potential sum for this i atom from the interaction with this j atom. */
382             velec            = _mm_unpacklo_pd(velec,_mm_setzero_pd());
383             velecsum         = _mm_add_pd(velecsum,velec);
384
385             fscal            = felec;
386
387             fscal            = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
388
389             /* Calculate temporary vectorial force */
390             tx               = _mm_mul_pd(fscal,dx20);
391             ty               = _mm_mul_pd(fscal,dy20);
392             tz               = _mm_mul_pd(fscal,dz20);
393
394             /* Update vectorial force */
395             fix2             = _mm_add_pd(fix2,tx);
396             fiy2             = _mm_add_pd(fiy2,ty);
397             fiz2             = _mm_add_pd(fiz2,tz);
398
399             fjx0             = _mm_add_pd(fjx0,tx);
400             fjy0             = _mm_add_pd(fjy0,ty);
401             fjz0             = _mm_add_pd(fjz0,tz);
402
403             /**************************
404              * CALCULATE INTERACTIONS *
405              **************************/
406
407             /* Compute parameters for interactions between i and j atoms */
408             qq30             = _mm_mul_pd(iq3,jq0);
409
410             /* COULOMB ELECTROSTATICS */
411             velec            = _mm_mul_pd(qq30,rinv30);
412             felec            = _mm_mul_pd(velec,rinvsq30);
413
414             /* Update potential sum for this i atom from the interaction with this j atom. */
415             velec            = _mm_unpacklo_pd(velec,_mm_setzero_pd());
416             velecsum         = _mm_add_pd(velecsum,velec);
417
418             fscal            = felec;
419
420             fscal            = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
421
422             /* Calculate temporary vectorial force */
423             tx               = _mm_mul_pd(fscal,dx30);
424             ty               = _mm_mul_pd(fscal,dy30);
425             tz               = _mm_mul_pd(fscal,dz30);
426
427             /* Update vectorial force */
428             fix3             = _mm_add_pd(fix3,tx);
429             fiy3             = _mm_add_pd(fiy3,ty);
430             fiz3             = _mm_add_pd(fiz3,tz);
431
432             fjx0             = _mm_add_pd(fjx0,tx);
433             fjy0             = _mm_add_pd(fjy0,ty);
434             fjz0             = _mm_add_pd(fjz0,tz);
435
436             gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,fjx0,fjy0,fjz0);
437
438             /* Inner loop uses 87 flops */
439         }
440
441         /* End of innermost loop */
442
443         gmx_mm_update_iforce_3atom_swizzle_pd(fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
444                                               f+i_coord_offset+DIM,fshift+i_shift_offset);
445
446         ggid                        = gid[iidx];
447         /* Update potential energies */
448         gmx_mm_update_1pot_pd(velecsum,kernel_data->energygrp_elec+ggid);
449
450         /* Increment number of inner iterations */
451         inneriter                  += j_index_end - j_index_start;
452
453         /* Outer loop uses 19 flops */
454     }
455
456     /* Increment number of outer iterations */
457     outeriter        += nri;
458
459     /* Update outer/inner flops */
460
461     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_W4_VF,outeriter*19 + inneriter*87);
462 }
463 /*
464  * Gromacs nonbonded kernel:   nb_kernel_ElecCoul_VdwNone_GeomW4P1_F_sse2_double
465  * Electrostatics interaction: Coulomb
466  * VdW interaction:            None
467  * Geometry:                   Water4-Particle
468  * Calculate force/pot:        Force
469  */
470 void
471 nb_kernel_ElecCoul_VdwNone_GeomW4P1_F_sse2_double
472                     (t_nblist                    * gmx_restrict       nlist,
473                      rvec                        * gmx_restrict          xx,
474                      rvec                        * gmx_restrict          ff,
475                      t_forcerec                  * gmx_restrict          fr,
476                      t_mdatoms                   * gmx_restrict     mdatoms,
477                      nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
478                      t_nrnb                      * gmx_restrict        nrnb)
479 {
480     /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
481      * just 0 for non-waters.
482      * Suffixes A,B refer to j loop unrolling done with SSE double precision, e.g. for the two different
483      * jnr indices corresponding to data put in the four positions in the SIMD register.
484      */
485     int              i_shift_offset,i_coord_offset,outeriter,inneriter;
486     int              j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
487     int              jnrA,jnrB;
488     int              j_coord_offsetA,j_coord_offsetB;
489     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
490     real             rcutoff_scalar;
491     real             *shiftvec,*fshift,*x,*f;
492     __m128d          tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
493     int              vdwioffset1;
494     __m128d          ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
495     int              vdwioffset2;
496     __m128d          ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
497     int              vdwioffset3;
498     __m128d          ix3,iy3,iz3,fix3,fiy3,fiz3,iq3,isai3;
499     int              vdwjidx0A,vdwjidx0B;
500     __m128d          jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
501     __m128d          dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
502     __m128d          dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
503     __m128d          dx30,dy30,dz30,rsq30,rinv30,rinvsq30,r30,qq30,c6_30,c12_30;
504     __m128d          velec,felec,velecsum,facel,crf,krf,krf2;
505     real             *charge;
506     __m128d          dummy_mask,cutoff_mask;
507     __m128d          signbit   = gmx_mm_castsi128_pd( _mm_set_epi32(0x80000000,0x00000000,0x80000000,0x00000000) );
508     __m128d          one     = _mm_set1_pd(1.0);
509     __m128d          two     = _mm_set1_pd(2.0);
510     x                = xx[0];
511     f                = ff[0];
512
513     nri              = nlist->nri;
514     iinr             = nlist->iinr;
515     jindex           = nlist->jindex;
516     jjnr             = nlist->jjnr;
517     shiftidx         = nlist->shift;
518     gid              = nlist->gid;
519     shiftvec         = fr->shift_vec[0];
520     fshift           = fr->fshift[0];
521     facel            = _mm_set1_pd(fr->epsfac);
522     charge           = mdatoms->chargeA;
523
524     /* Setup water-specific parameters */
525     inr              = nlist->iinr[0];
526     iq1              = _mm_mul_pd(facel,_mm_set1_pd(charge[inr+1]));
527     iq2              = _mm_mul_pd(facel,_mm_set1_pd(charge[inr+2]));
528     iq3              = _mm_mul_pd(facel,_mm_set1_pd(charge[inr+3]));
529
530     /* Avoid stupid compiler warnings */
531     jnrA = jnrB = 0;
532     j_coord_offsetA = 0;
533     j_coord_offsetB = 0;
534
535     outeriter        = 0;
536     inneriter        = 0;
537
538     /* Start outer loop over neighborlists */
539     for(iidx=0; iidx<nri; iidx++)
540     {
541         /* Load shift vector for this list */
542         i_shift_offset   = DIM*shiftidx[iidx];
543
544         /* Load limits for loop over neighbors */
545         j_index_start    = jindex[iidx];
546         j_index_end      = jindex[iidx+1];
547
548         /* Get outer coordinate index */
549         inr              = iinr[iidx];
550         i_coord_offset   = DIM*inr;
551
552         /* Load i particle coords and add shift vector */
553         gmx_mm_load_shift_and_3rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset+DIM,
554                                                  &ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
555
556         fix1             = _mm_setzero_pd();
557         fiy1             = _mm_setzero_pd();
558         fiz1             = _mm_setzero_pd();
559         fix2             = _mm_setzero_pd();
560         fiy2             = _mm_setzero_pd();
561         fiz2             = _mm_setzero_pd();
562         fix3             = _mm_setzero_pd();
563         fiy3             = _mm_setzero_pd();
564         fiz3             = _mm_setzero_pd();
565
566         /* Start inner kernel loop */
567         for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
568         {
569
570             /* Get j neighbor index, and coordinate index */
571             jnrA             = jjnr[jidx];
572             jnrB             = jjnr[jidx+1];
573             j_coord_offsetA  = DIM*jnrA;
574             j_coord_offsetB  = DIM*jnrB;
575
576             /* load j atom coordinates */
577             gmx_mm_load_1rvec_2ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
578                                               &jx0,&jy0,&jz0);
579
580             /* Calculate displacement vector */
581             dx10             = _mm_sub_pd(ix1,jx0);
582             dy10             = _mm_sub_pd(iy1,jy0);
583             dz10             = _mm_sub_pd(iz1,jz0);
584             dx20             = _mm_sub_pd(ix2,jx0);
585             dy20             = _mm_sub_pd(iy2,jy0);
586             dz20             = _mm_sub_pd(iz2,jz0);
587             dx30             = _mm_sub_pd(ix3,jx0);
588             dy30             = _mm_sub_pd(iy3,jy0);
589             dz30             = _mm_sub_pd(iz3,jz0);
590
591             /* Calculate squared distance and things based on it */
592             rsq10            = gmx_mm_calc_rsq_pd(dx10,dy10,dz10);
593             rsq20            = gmx_mm_calc_rsq_pd(dx20,dy20,dz20);
594             rsq30            = gmx_mm_calc_rsq_pd(dx30,dy30,dz30);
595
596             rinv10           = gmx_mm_invsqrt_pd(rsq10);
597             rinv20           = gmx_mm_invsqrt_pd(rsq20);
598             rinv30           = gmx_mm_invsqrt_pd(rsq30);
599
600             rinvsq10         = _mm_mul_pd(rinv10,rinv10);
601             rinvsq20         = _mm_mul_pd(rinv20,rinv20);
602             rinvsq30         = _mm_mul_pd(rinv30,rinv30);
603
604             /* Load parameters for j particles */
605             jq0              = gmx_mm_load_2real_swizzle_pd(charge+jnrA+0,charge+jnrB+0);
606
607             fjx0             = _mm_setzero_pd();
608             fjy0             = _mm_setzero_pd();
609             fjz0             = _mm_setzero_pd();
610
611             /**************************
612              * CALCULATE INTERACTIONS *
613              **************************/
614
615             /* Compute parameters for interactions between i and j atoms */
616             qq10             = _mm_mul_pd(iq1,jq0);
617
618             /* COULOMB ELECTROSTATICS */
619             velec            = _mm_mul_pd(qq10,rinv10);
620             felec            = _mm_mul_pd(velec,rinvsq10);
621
622             fscal            = felec;
623
624             /* Calculate temporary vectorial force */
625             tx               = _mm_mul_pd(fscal,dx10);
626             ty               = _mm_mul_pd(fscal,dy10);
627             tz               = _mm_mul_pd(fscal,dz10);
628
629             /* Update vectorial force */
630             fix1             = _mm_add_pd(fix1,tx);
631             fiy1             = _mm_add_pd(fiy1,ty);
632             fiz1             = _mm_add_pd(fiz1,tz);
633
634             fjx0             = _mm_add_pd(fjx0,tx);
635             fjy0             = _mm_add_pd(fjy0,ty);
636             fjz0             = _mm_add_pd(fjz0,tz);
637
638             /**************************
639              * CALCULATE INTERACTIONS *
640              **************************/
641
642             /* Compute parameters for interactions between i and j atoms */
643             qq20             = _mm_mul_pd(iq2,jq0);
644
645             /* COULOMB ELECTROSTATICS */
646             velec            = _mm_mul_pd(qq20,rinv20);
647             felec            = _mm_mul_pd(velec,rinvsq20);
648
649             fscal            = felec;
650
651             /* Calculate temporary vectorial force */
652             tx               = _mm_mul_pd(fscal,dx20);
653             ty               = _mm_mul_pd(fscal,dy20);
654             tz               = _mm_mul_pd(fscal,dz20);
655
656             /* Update vectorial force */
657             fix2             = _mm_add_pd(fix2,tx);
658             fiy2             = _mm_add_pd(fiy2,ty);
659             fiz2             = _mm_add_pd(fiz2,tz);
660
661             fjx0             = _mm_add_pd(fjx0,tx);
662             fjy0             = _mm_add_pd(fjy0,ty);
663             fjz0             = _mm_add_pd(fjz0,tz);
664
665             /**************************
666              * CALCULATE INTERACTIONS *
667              **************************/
668
669             /* Compute parameters for interactions between i and j atoms */
670             qq30             = _mm_mul_pd(iq3,jq0);
671
672             /* COULOMB ELECTROSTATICS */
673             velec            = _mm_mul_pd(qq30,rinv30);
674             felec            = _mm_mul_pd(velec,rinvsq30);
675
676             fscal            = felec;
677
678             /* Calculate temporary vectorial force */
679             tx               = _mm_mul_pd(fscal,dx30);
680             ty               = _mm_mul_pd(fscal,dy30);
681             tz               = _mm_mul_pd(fscal,dz30);
682
683             /* Update vectorial force */
684             fix3             = _mm_add_pd(fix3,tx);
685             fiy3             = _mm_add_pd(fiy3,ty);
686             fiz3             = _mm_add_pd(fiz3,tz);
687
688             fjx0             = _mm_add_pd(fjx0,tx);
689             fjy0             = _mm_add_pd(fjy0,ty);
690             fjz0             = _mm_add_pd(fjz0,tz);
691
692             gmx_mm_decrement_1rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,fjx0,fjy0,fjz0);
693
694             /* Inner loop uses 84 flops */
695         }
696
697         if(jidx<j_index_end)
698         {
699
700             jnrA             = jjnr[jidx];
701             j_coord_offsetA  = DIM*jnrA;
702
703             /* load j atom coordinates */
704             gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
705                                               &jx0,&jy0,&jz0);
706
707             /* Calculate displacement vector */
708             dx10             = _mm_sub_pd(ix1,jx0);
709             dy10             = _mm_sub_pd(iy1,jy0);
710             dz10             = _mm_sub_pd(iz1,jz0);
711             dx20             = _mm_sub_pd(ix2,jx0);
712             dy20             = _mm_sub_pd(iy2,jy0);
713             dz20             = _mm_sub_pd(iz2,jz0);
714             dx30             = _mm_sub_pd(ix3,jx0);
715             dy30             = _mm_sub_pd(iy3,jy0);
716             dz30             = _mm_sub_pd(iz3,jz0);
717
718             /* Calculate squared distance and things based on it */
719             rsq10            = gmx_mm_calc_rsq_pd(dx10,dy10,dz10);
720             rsq20            = gmx_mm_calc_rsq_pd(dx20,dy20,dz20);
721             rsq30            = gmx_mm_calc_rsq_pd(dx30,dy30,dz30);
722
723             rinv10           = gmx_mm_invsqrt_pd(rsq10);
724             rinv20           = gmx_mm_invsqrt_pd(rsq20);
725             rinv30           = gmx_mm_invsqrt_pd(rsq30);
726
727             rinvsq10         = _mm_mul_pd(rinv10,rinv10);
728             rinvsq20         = _mm_mul_pd(rinv20,rinv20);
729             rinvsq30         = _mm_mul_pd(rinv30,rinv30);
730
731             /* Load parameters for j particles */
732             jq0              = _mm_load_sd(charge+jnrA+0);
733
734             fjx0             = _mm_setzero_pd();
735             fjy0             = _mm_setzero_pd();
736             fjz0             = _mm_setzero_pd();
737
738             /**************************
739              * CALCULATE INTERACTIONS *
740              **************************/
741
742             /* Compute parameters for interactions between i and j atoms */
743             qq10             = _mm_mul_pd(iq1,jq0);
744
745             /* COULOMB ELECTROSTATICS */
746             velec            = _mm_mul_pd(qq10,rinv10);
747             felec            = _mm_mul_pd(velec,rinvsq10);
748
749             fscal            = felec;
750
751             fscal            = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
752
753             /* Calculate temporary vectorial force */
754             tx               = _mm_mul_pd(fscal,dx10);
755             ty               = _mm_mul_pd(fscal,dy10);
756             tz               = _mm_mul_pd(fscal,dz10);
757
758             /* Update vectorial force */
759             fix1             = _mm_add_pd(fix1,tx);
760             fiy1             = _mm_add_pd(fiy1,ty);
761             fiz1             = _mm_add_pd(fiz1,tz);
762
763             fjx0             = _mm_add_pd(fjx0,tx);
764             fjy0             = _mm_add_pd(fjy0,ty);
765             fjz0             = _mm_add_pd(fjz0,tz);
766
767             /**************************
768              * CALCULATE INTERACTIONS *
769              **************************/
770
771             /* Compute parameters for interactions between i and j atoms */
772             qq20             = _mm_mul_pd(iq2,jq0);
773
774             /* COULOMB ELECTROSTATICS */
775             velec            = _mm_mul_pd(qq20,rinv20);
776             felec            = _mm_mul_pd(velec,rinvsq20);
777
778             fscal            = felec;
779
780             fscal            = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
781
782             /* Calculate temporary vectorial force */
783             tx               = _mm_mul_pd(fscal,dx20);
784             ty               = _mm_mul_pd(fscal,dy20);
785             tz               = _mm_mul_pd(fscal,dz20);
786
787             /* Update vectorial force */
788             fix2             = _mm_add_pd(fix2,tx);
789             fiy2             = _mm_add_pd(fiy2,ty);
790             fiz2             = _mm_add_pd(fiz2,tz);
791
792             fjx0             = _mm_add_pd(fjx0,tx);
793             fjy0             = _mm_add_pd(fjy0,ty);
794             fjz0             = _mm_add_pd(fjz0,tz);
795
796             /**************************
797              * CALCULATE INTERACTIONS *
798              **************************/
799
800             /* Compute parameters for interactions between i and j atoms */
801             qq30             = _mm_mul_pd(iq3,jq0);
802
803             /* COULOMB ELECTROSTATICS */
804             velec            = _mm_mul_pd(qq30,rinv30);
805             felec            = _mm_mul_pd(velec,rinvsq30);
806
807             fscal            = felec;
808
809             fscal            = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
810
811             /* Calculate temporary vectorial force */
812             tx               = _mm_mul_pd(fscal,dx30);
813             ty               = _mm_mul_pd(fscal,dy30);
814             tz               = _mm_mul_pd(fscal,dz30);
815
816             /* Update vectorial force */
817             fix3             = _mm_add_pd(fix3,tx);
818             fiy3             = _mm_add_pd(fiy3,ty);
819             fiz3             = _mm_add_pd(fiz3,tz);
820
821             fjx0             = _mm_add_pd(fjx0,tx);
822             fjy0             = _mm_add_pd(fjy0,ty);
823             fjz0             = _mm_add_pd(fjz0,tz);
824
825             gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,fjx0,fjy0,fjz0);
826
827             /* Inner loop uses 84 flops */
828         }
829
830         /* End of innermost loop */
831
832         gmx_mm_update_iforce_3atom_swizzle_pd(fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
833                                               f+i_coord_offset+DIM,fshift+i_shift_offset);
834
835         /* Increment number of inner iterations */
836         inneriter                  += j_index_end - j_index_start;
837
838         /* Outer loop uses 18 flops */
839     }
840
841     /* Increment number of outer iterations */
842     outeriter        += nri;
843
844     /* Update outer/inner flops */
845
846     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_W4_F,outeriter*18 + inneriter*84);
847 }