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