4f14c165da729916b3e84779558417778d20d073
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_avx_128_fma_single / nb_kernel_ElecRFCut_VdwNone_GeomW3P1_avx_128_fma_single.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 avx_128_fma_single kernel generator.
37  */
38 #ifdef HAVE_CONFIG_H
39 #include <config.h>
40 #endif
41
42 #include <math.h>
43
44 #include "../nb_kernel.h"
45 #include "types/simple.h"
46 #include "gromacs/math/vec.h"
47 #include "nrnb.h"
48
49 #include "gromacs/simd/math_x86_avx_128_fma_single.h"
50 #include "kernelutil_x86_avx_128_fma_single.h"
51
52 /*
53  * Gromacs nonbonded kernel:   nb_kernel_ElecRFCut_VdwNone_GeomW3P1_VF_avx_128_fma_single
54  * Electrostatics interaction: ReactionField
55  * VdW interaction:            None
56  * Geometry:                   Water3-Particle
57  * Calculate force/pot:        PotentialAndForce
58  */
59 void
60 nb_kernel_ElecRFCut_VdwNone_GeomW3P1_VF_avx_128_fma_single
61                     (t_nblist                    * gmx_restrict       nlist,
62                      rvec                        * gmx_restrict          xx,
63                      rvec                        * gmx_restrict          ff,
64                      t_forcerec                  * gmx_restrict          fr,
65                      t_mdatoms                   * gmx_restrict     mdatoms,
66                      nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
67                      t_nrnb                      * gmx_restrict        nrnb)
68 {
69     /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
70      * just 0 for non-waters.
71      * Suffixes A,B,C,D refer to j loop unrolling done with AVX_128, e.g. for the four different
72      * jnr indices corresponding to data put in the four positions in the SIMD register.
73      */
74     int              i_shift_offset,i_coord_offset,outeriter,inneriter;
75     int              j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
76     int              jnrA,jnrB,jnrC,jnrD;
77     int              jnrlistA,jnrlistB,jnrlistC,jnrlistD;
78     int              j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
79     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
80     real             rcutoff_scalar;
81     real             *shiftvec,*fshift,*x,*f;
82     real             *fjptrA,*fjptrB,*fjptrC,*fjptrD;
83     real             scratch[4*DIM];
84     __m128           fscal,rcutoff,rcutoff2,jidxall;
85     int              vdwioffset0;
86     __m128           ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
87     int              vdwioffset1;
88     __m128           ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
89     int              vdwioffset2;
90     __m128           ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
91     int              vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
92     __m128           jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
93     __m128           dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
94     __m128           dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
95     __m128           dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
96     __m128           velec,felec,velecsum,facel,crf,krf,krf2;
97     real             *charge;
98     __m128           dummy_mask,cutoff_mask;
99     __m128           signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
100     __m128           one     = _mm_set1_ps(1.0);
101     __m128           two     = _mm_set1_ps(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_ps(fr->epsfac);
114     charge           = mdatoms->chargeA;
115     krf              = _mm_set1_ps(fr->ic->k_rf);
116     krf2             = _mm_set1_ps(fr->ic->k_rf*2.0);
117     crf              = _mm_set1_ps(fr->ic->c_rf);
118
119     /* Setup water-specific parameters */
120     inr              = nlist->iinr[0];
121     iq0              = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+0]));
122     iq1              = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+1]));
123     iq2              = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+2]));
124
125     /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
126     rcutoff_scalar   = fr->rcoulomb;
127     rcutoff          = _mm_set1_ps(rcutoff_scalar);
128     rcutoff2         = _mm_mul_ps(rcutoff,rcutoff);
129
130     /* Avoid stupid compiler warnings */
131     jnrA = jnrB = jnrC = jnrD = 0;
132     j_coord_offsetA = 0;
133     j_coord_offsetB = 0;
134     j_coord_offsetC = 0;
135     j_coord_offsetD = 0;
136
137     outeriter        = 0;
138     inneriter        = 0;
139
140     for(iidx=0;iidx<4*DIM;iidx++)
141     {
142         scratch[iidx] = 0.0;
143     }
144
145     /* Start outer loop over neighborlists */
146     for(iidx=0; iidx<nri; iidx++)
147     {
148         /* Load shift vector for this list */
149         i_shift_offset   = DIM*shiftidx[iidx];
150
151         /* Load limits for loop over neighbors */
152         j_index_start    = jindex[iidx];
153         j_index_end      = jindex[iidx+1];
154
155         /* Get outer coordinate index */
156         inr              = iinr[iidx];
157         i_coord_offset   = DIM*inr;
158
159         /* Load i particle coords and add shift vector */
160         gmx_mm_load_shift_and_3rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,
161                                                  &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2);
162
163         fix0             = _mm_setzero_ps();
164         fiy0             = _mm_setzero_ps();
165         fiz0             = _mm_setzero_ps();
166         fix1             = _mm_setzero_ps();
167         fiy1             = _mm_setzero_ps();
168         fiz1             = _mm_setzero_ps();
169         fix2             = _mm_setzero_ps();
170         fiy2             = _mm_setzero_ps();
171         fiz2             = _mm_setzero_ps();
172
173         /* Reset potential sums */
174         velecsum         = _mm_setzero_ps();
175
176         /* Start inner kernel loop */
177         for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
178         {
179
180             /* Get j neighbor index, and coordinate index */
181             jnrA             = jjnr[jidx];
182             jnrB             = jjnr[jidx+1];
183             jnrC             = jjnr[jidx+2];
184             jnrD             = jjnr[jidx+3];
185             j_coord_offsetA  = DIM*jnrA;
186             j_coord_offsetB  = DIM*jnrB;
187             j_coord_offsetC  = DIM*jnrC;
188             j_coord_offsetD  = DIM*jnrD;
189
190             /* load j atom coordinates */
191             gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
192                                               x+j_coord_offsetC,x+j_coord_offsetD,
193                                               &jx0,&jy0,&jz0);
194
195             /* Calculate displacement vector */
196             dx00             = _mm_sub_ps(ix0,jx0);
197             dy00             = _mm_sub_ps(iy0,jy0);
198             dz00             = _mm_sub_ps(iz0,jz0);
199             dx10             = _mm_sub_ps(ix1,jx0);
200             dy10             = _mm_sub_ps(iy1,jy0);
201             dz10             = _mm_sub_ps(iz1,jz0);
202             dx20             = _mm_sub_ps(ix2,jx0);
203             dy20             = _mm_sub_ps(iy2,jy0);
204             dz20             = _mm_sub_ps(iz2,jz0);
205
206             /* Calculate squared distance and things based on it */
207             rsq00            = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
208             rsq10            = gmx_mm_calc_rsq_ps(dx10,dy10,dz10);
209             rsq20            = gmx_mm_calc_rsq_ps(dx20,dy20,dz20);
210
211             rinv00           = gmx_mm_invsqrt_ps(rsq00);
212             rinv10           = gmx_mm_invsqrt_ps(rsq10);
213             rinv20           = gmx_mm_invsqrt_ps(rsq20);
214
215             rinvsq00         = _mm_mul_ps(rinv00,rinv00);
216             rinvsq10         = _mm_mul_ps(rinv10,rinv10);
217             rinvsq20         = _mm_mul_ps(rinv20,rinv20);
218
219             /* Load parameters for j particles */
220             jq0              = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
221                                                               charge+jnrC+0,charge+jnrD+0);
222
223             fjx0             = _mm_setzero_ps();
224             fjy0             = _mm_setzero_ps();
225             fjz0             = _mm_setzero_ps();
226
227             /**************************
228              * CALCULATE INTERACTIONS *
229              **************************/
230
231             if (gmx_mm_any_lt(rsq00,rcutoff2))
232             {
233
234             /* Compute parameters for interactions between i and j atoms */
235             qq00             = _mm_mul_ps(iq0,jq0);
236
237             /* REACTION-FIELD ELECTROSTATICS */
238             velec            = _mm_mul_ps(qq00,_mm_sub_ps(_mm_macc_ps(krf,rsq00,rinv00),crf));
239             felec            = _mm_mul_ps(qq00,_mm_msub_ps(rinv00,rinvsq00,krf2));
240
241             cutoff_mask      = _mm_cmplt_ps(rsq00,rcutoff2);
242
243             /* Update potential sum for this i atom from the interaction with this j atom. */
244             velec            = _mm_and_ps(velec,cutoff_mask);
245             velecsum         = _mm_add_ps(velecsum,velec);
246
247             fscal            = felec;
248
249             fscal            = _mm_and_ps(fscal,cutoff_mask);
250
251              /* Update vectorial force */
252             fix0             = _mm_macc_ps(dx00,fscal,fix0);
253             fiy0             = _mm_macc_ps(dy00,fscal,fiy0);
254             fiz0             = _mm_macc_ps(dz00,fscal,fiz0);
255
256             fjx0             = _mm_macc_ps(dx00,fscal,fjx0);
257             fjy0             = _mm_macc_ps(dy00,fscal,fjy0);
258             fjz0             = _mm_macc_ps(dz00,fscal,fjz0);
259
260             }
261
262             /**************************
263              * CALCULATE INTERACTIONS *
264              **************************/
265
266             if (gmx_mm_any_lt(rsq10,rcutoff2))
267             {
268
269             /* Compute parameters for interactions between i and j atoms */
270             qq10             = _mm_mul_ps(iq1,jq0);
271
272             /* REACTION-FIELD ELECTROSTATICS */
273             velec            = _mm_mul_ps(qq10,_mm_sub_ps(_mm_macc_ps(krf,rsq10,rinv10),crf));
274             felec            = _mm_mul_ps(qq10,_mm_msub_ps(rinv10,rinvsq10,krf2));
275
276             cutoff_mask      = _mm_cmplt_ps(rsq10,rcutoff2);
277
278             /* Update potential sum for this i atom from the interaction with this j atom. */
279             velec            = _mm_and_ps(velec,cutoff_mask);
280             velecsum         = _mm_add_ps(velecsum,velec);
281
282             fscal            = felec;
283
284             fscal            = _mm_and_ps(fscal,cutoff_mask);
285
286              /* Update vectorial force */
287             fix1             = _mm_macc_ps(dx10,fscal,fix1);
288             fiy1             = _mm_macc_ps(dy10,fscal,fiy1);
289             fiz1             = _mm_macc_ps(dz10,fscal,fiz1);
290
291             fjx0             = _mm_macc_ps(dx10,fscal,fjx0);
292             fjy0             = _mm_macc_ps(dy10,fscal,fjy0);
293             fjz0             = _mm_macc_ps(dz10,fscal,fjz0);
294
295             }
296
297             /**************************
298              * CALCULATE INTERACTIONS *
299              **************************/
300
301             if (gmx_mm_any_lt(rsq20,rcutoff2))
302             {
303
304             /* Compute parameters for interactions between i and j atoms */
305             qq20             = _mm_mul_ps(iq2,jq0);
306
307             /* REACTION-FIELD ELECTROSTATICS */
308             velec            = _mm_mul_ps(qq20,_mm_sub_ps(_mm_macc_ps(krf,rsq20,rinv20),crf));
309             felec            = _mm_mul_ps(qq20,_mm_msub_ps(rinv20,rinvsq20,krf2));
310
311             cutoff_mask      = _mm_cmplt_ps(rsq20,rcutoff2);
312
313             /* Update potential sum for this i atom from the interaction with this j atom. */
314             velec            = _mm_and_ps(velec,cutoff_mask);
315             velecsum         = _mm_add_ps(velecsum,velec);
316
317             fscal            = felec;
318
319             fscal            = _mm_and_ps(fscal,cutoff_mask);
320
321              /* Update vectorial force */
322             fix2             = _mm_macc_ps(dx20,fscal,fix2);
323             fiy2             = _mm_macc_ps(dy20,fscal,fiy2);
324             fiz2             = _mm_macc_ps(dz20,fscal,fiz2);
325
326             fjx0             = _mm_macc_ps(dx20,fscal,fjx0);
327             fjy0             = _mm_macc_ps(dy20,fscal,fjy0);
328             fjz0             = _mm_macc_ps(dz20,fscal,fjz0);
329
330             }
331
332             fjptrA             = f+j_coord_offsetA;
333             fjptrB             = f+j_coord_offsetB;
334             fjptrC             = f+j_coord_offsetC;
335             fjptrD             = f+j_coord_offsetD;
336
337             gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjx0,fjy0,fjz0);
338
339             /* Inner loop uses 117 flops */
340         }
341
342         if(jidx<j_index_end)
343         {
344
345             /* Get j neighbor index, and coordinate index */
346             jnrlistA         = jjnr[jidx];
347             jnrlistB         = jjnr[jidx+1];
348             jnrlistC         = jjnr[jidx+2];
349             jnrlistD         = jjnr[jidx+3];
350             /* Sign of each element will be negative for non-real atoms.
351              * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
352              * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
353              */
354             dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
355             jnrA       = (jnrlistA>=0) ? jnrlistA : 0;
356             jnrB       = (jnrlistB>=0) ? jnrlistB : 0;
357             jnrC       = (jnrlistC>=0) ? jnrlistC : 0;
358             jnrD       = (jnrlistD>=0) ? jnrlistD : 0;
359             j_coord_offsetA  = DIM*jnrA;
360             j_coord_offsetB  = DIM*jnrB;
361             j_coord_offsetC  = DIM*jnrC;
362             j_coord_offsetD  = DIM*jnrD;
363
364             /* load j atom coordinates */
365             gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
366                                               x+j_coord_offsetC,x+j_coord_offsetD,
367                                               &jx0,&jy0,&jz0);
368
369             /* Calculate displacement vector */
370             dx00             = _mm_sub_ps(ix0,jx0);
371             dy00             = _mm_sub_ps(iy0,jy0);
372             dz00             = _mm_sub_ps(iz0,jz0);
373             dx10             = _mm_sub_ps(ix1,jx0);
374             dy10             = _mm_sub_ps(iy1,jy0);
375             dz10             = _mm_sub_ps(iz1,jz0);
376             dx20             = _mm_sub_ps(ix2,jx0);
377             dy20             = _mm_sub_ps(iy2,jy0);
378             dz20             = _mm_sub_ps(iz2,jz0);
379
380             /* Calculate squared distance and things based on it */
381             rsq00            = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
382             rsq10            = gmx_mm_calc_rsq_ps(dx10,dy10,dz10);
383             rsq20            = gmx_mm_calc_rsq_ps(dx20,dy20,dz20);
384
385             rinv00           = gmx_mm_invsqrt_ps(rsq00);
386             rinv10           = gmx_mm_invsqrt_ps(rsq10);
387             rinv20           = gmx_mm_invsqrt_ps(rsq20);
388
389             rinvsq00         = _mm_mul_ps(rinv00,rinv00);
390             rinvsq10         = _mm_mul_ps(rinv10,rinv10);
391             rinvsq20         = _mm_mul_ps(rinv20,rinv20);
392
393             /* Load parameters for j particles */
394             jq0              = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
395                                                               charge+jnrC+0,charge+jnrD+0);
396
397             fjx0             = _mm_setzero_ps();
398             fjy0             = _mm_setzero_ps();
399             fjz0             = _mm_setzero_ps();
400
401             /**************************
402              * CALCULATE INTERACTIONS *
403              **************************/
404
405             if (gmx_mm_any_lt(rsq00,rcutoff2))
406             {
407
408             /* Compute parameters for interactions between i and j atoms */
409             qq00             = _mm_mul_ps(iq0,jq0);
410
411             /* REACTION-FIELD ELECTROSTATICS */
412             velec            = _mm_mul_ps(qq00,_mm_sub_ps(_mm_macc_ps(krf,rsq00,rinv00),crf));
413             felec            = _mm_mul_ps(qq00,_mm_msub_ps(rinv00,rinvsq00,krf2));
414
415             cutoff_mask      = _mm_cmplt_ps(rsq00,rcutoff2);
416
417             /* Update potential sum for this i atom from the interaction with this j atom. */
418             velec            = _mm_and_ps(velec,cutoff_mask);
419             velec            = _mm_andnot_ps(dummy_mask,velec);
420             velecsum         = _mm_add_ps(velecsum,velec);
421
422             fscal            = felec;
423
424             fscal            = _mm_and_ps(fscal,cutoff_mask);
425
426             fscal            = _mm_andnot_ps(dummy_mask,fscal);
427
428              /* Update vectorial force */
429             fix0             = _mm_macc_ps(dx00,fscal,fix0);
430             fiy0             = _mm_macc_ps(dy00,fscal,fiy0);
431             fiz0             = _mm_macc_ps(dz00,fscal,fiz0);
432
433             fjx0             = _mm_macc_ps(dx00,fscal,fjx0);
434             fjy0             = _mm_macc_ps(dy00,fscal,fjy0);
435             fjz0             = _mm_macc_ps(dz00,fscal,fjz0);
436
437             }
438
439             /**************************
440              * CALCULATE INTERACTIONS *
441              **************************/
442
443             if (gmx_mm_any_lt(rsq10,rcutoff2))
444             {
445
446             /* Compute parameters for interactions between i and j atoms */
447             qq10             = _mm_mul_ps(iq1,jq0);
448
449             /* REACTION-FIELD ELECTROSTATICS */
450             velec            = _mm_mul_ps(qq10,_mm_sub_ps(_mm_macc_ps(krf,rsq10,rinv10),crf));
451             felec            = _mm_mul_ps(qq10,_mm_msub_ps(rinv10,rinvsq10,krf2));
452
453             cutoff_mask      = _mm_cmplt_ps(rsq10,rcutoff2);
454
455             /* Update potential sum for this i atom from the interaction with this j atom. */
456             velec            = _mm_and_ps(velec,cutoff_mask);
457             velec            = _mm_andnot_ps(dummy_mask,velec);
458             velecsum         = _mm_add_ps(velecsum,velec);
459
460             fscal            = felec;
461
462             fscal            = _mm_and_ps(fscal,cutoff_mask);
463
464             fscal            = _mm_andnot_ps(dummy_mask,fscal);
465
466              /* Update vectorial force */
467             fix1             = _mm_macc_ps(dx10,fscal,fix1);
468             fiy1             = _mm_macc_ps(dy10,fscal,fiy1);
469             fiz1             = _mm_macc_ps(dz10,fscal,fiz1);
470
471             fjx0             = _mm_macc_ps(dx10,fscal,fjx0);
472             fjy0             = _mm_macc_ps(dy10,fscal,fjy0);
473             fjz0             = _mm_macc_ps(dz10,fscal,fjz0);
474
475             }
476
477             /**************************
478              * CALCULATE INTERACTIONS *
479              **************************/
480
481             if (gmx_mm_any_lt(rsq20,rcutoff2))
482             {
483
484             /* Compute parameters for interactions between i and j atoms */
485             qq20             = _mm_mul_ps(iq2,jq0);
486
487             /* REACTION-FIELD ELECTROSTATICS */
488             velec            = _mm_mul_ps(qq20,_mm_sub_ps(_mm_macc_ps(krf,rsq20,rinv20),crf));
489             felec            = _mm_mul_ps(qq20,_mm_msub_ps(rinv20,rinvsq20,krf2));
490
491             cutoff_mask      = _mm_cmplt_ps(rsq20,rcutoff2);
492
493             /* Update potential sum for this i atom from the interaction with this j atom. */
494             velec            = _mm_and_ps(velec,cutoff_mask);
495             velec            = _mm_andnot_ps(dummy_mask,velec);
496             velecsum         = _mm_add_ps(velecsum,velec);
497
498             fscal            = felec;
499
500             fscal            = _mm_and_ps(fscal,cutoff_mask);
501
502             fscal            = _mm_andnot_ps(dummy_mask,fscal);
503
504              /* Update vectorial force */
505             fix2             = _mm_macc_ps(dx20,fscal,fix2);
506             fiy2             = _mm_macc_ps(dy20,fscal,fiy2);
507             fiz2             = _mm_macc_ps(dz20,fscal,fiz2);
508
509             fjx0             = _mm_macc_ps(dx20,fscal,fjx0);
510             fjy0             = _mm_macc_ps(dy20,fscal,fjy0);
511             fjz0             = _mm_macc_ps(dz20,fscal,fjz0);
512
513             }
514
515             fjptrA             = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
516             fjptrB             = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
517             fjptrC             = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
518             fjptrD             = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
519
520             gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjx0,fjy0,fjz0);
521
522             /* Inner loop uses 117 flops */
523         }
524
525         /* End of innermost loop */
526
527         gmx_mm_update_iforce_3atom_swizzle_ps(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,
528                                               f+i_coord_offset,fshift+i_shift_offset);
529
530         ggid                        = gid[iidx];
531         /* Update potential energies */
532         gmx_mm_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
533
534         /* Increment number of inner iterations */
535         inneriter                  += j_index_end - j_index_start;
536
537         /* Outer loop uses 19 flops */
538     }
539
540     /* Increment number of outer iterations */
541     outeriter        += nri;
542
543     /* Update outer/inner flops */
544
545     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_W3_VF,outeriter*19 + inneriter*117);
546 }
547 /*
548  * Gromacs nonbonded kernel:   nb_kernel_ElecRFCut_VdwNone_GeomW3P1_F_avx_128_fma_single
549  * Electrostatics interaction: ReactionField
550  * VdW interaction:            None
551  * Geometry:                   Water3-Particle
552  * Calculate force/pot:        Force
553  */
554 void
555 nb_kernel_ElecRFCut_VdwNone_GeomW3P1_F_avx_128_fma_single
556                     (t_nblist                    * gmx_restrict       nlist,
557                      rvec                        * gmx_restrict          xx,
558                      rvec                        * gmx_restrict          ff,
559                      t_forcerec                  * gmx_restrict          fr,
560                      t_mdatoms                   * gmx_restrict     mdatoms,
561                      nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
562                      t_nrnb                      * gmx_restrict        nrnb)
563 {
564     /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
565      * just 0 for non-waters.
566      * Suffixes A,B,C,D refer to j loop unrolling done with AVX_128, e.g. for the four different
567      * jnr indices corresponding to data put in the four positions in the SIMD register.
568      */
569     int              i_shift_offset,i_coord_offset,outeriter,inneriter;
570     int              j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
571     int              jnrA,jnrB,jnrC,jnrD;
572     int              jnrlistA,jnrlistB,jnrlistC,jnrlistD;
573     int              j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
574     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
575     real             rcutoff_scalar;
576     real             *shiftvec,*fshift,*x,*f;
577     real             *fjptrA,*fjptrB,*fjptrC,*fjptrD;
578     real             scratch[4*DIM];
579     __m128           fscal,rcutoff,rcutoff2,jidxall;
580     int              vdwioffset0;
581     __m128           ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
582     int              vdwioffset1;
583     __m128           ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
584     int              vdwioffset2;
585     __m128           ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
586     int              vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
587     __m128           jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
588     __m128           dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
589     __m128           dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
590     __m128           dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
591     __m128           velec,felec,velecsum,facel,crf,krf,krf2;
592     real             *charge;
593     __m128           dummy_mask,cutoff_mask;
594     __m128           signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
595     __m128           one     = _mm_set1_ps(1.0);
596     __m128           two     = _mm_set1_ps(2.0);
597     x                = xx[0];
598     f                = ff[0];
599
600     nri              = nlist->nri;
601     iinr             = nlist->iinr;
602     jindex           = nlist->jindex;
603     jjnr             = nlist->jjnr;
604     shiftidx         = nlist->shift;
605     gid              = nlist->gid;
606     shiftvec         = fr->shift_vec[0];
607     fshift           = fr->fshift[0];
608     facel            = _mm_set1_ps(fr->epsfac);
609     charge           = mdatoms->chargeA;
610     krf              = _mm_set1_ps(fr->ic->k_rf);
611     krf2             = _mm_set1_ps(fr->ic->k_rf*2.0);
612     crf              = _mm_set1_ps(fr->ic->c_rf);
613
614     /* Setup water-specific parameters */
615     inr              = nlist->iinr[0];
616     iq0              = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+0]));
617     iq1              = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+1]));
618     iq2              = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+2]));
619
620     /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
621     rcutoff_scalar   = fr->rcoulomb;
622     rcutoff          = _mm_set1_ps(rcutoff_scalar);
623     rcutoff2         = _mm_mul_ps(rcutoff,rcutoff);
624
625     /* Avoid stupid compiler warnings */
626     jnrA = jnrB = jnrC = jnrD = 0;
627     j_coord_offsetA = 0;
628     j_coord_offsetB = 0;
629     j_coord_offsetC = 0;
630     j_coord_offsetD = 0;
631
632     outeriter        = 0;
633     inneriter        = 0;
634
635     for(iidx=0;iidx<4*DIM;iidx++)
636     {
637         scratch[iidx] = 0.0;
638     }
639
640     /* Start outer loop over neighborlists */
641     for(iidx=0; iidx<nri; iidx++)
642     {
643         /* Load shift vector for this list */
644         i_shift_offset   = DIM*shiftidx[iidx];
645
646         /* Load limits for loop over neighbors */
647         j_index_start    = jindex[iidx];
648         j_index_end      = jindex[iidx+1];
649
650         /* Get outer coordinate index */
651         inr              = iinr[iidx];
652         i_coord_offset   = DIM*inr;
653
654         /* Load i particle coords and add shift vector */
655         gmx_mm_load_shift_and_3rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,
656                                                  &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2);
657
658         fix0             = _mm_setzero_ps();
659         fiy0             = _mm_setzero_ps();
660         fiz0             = _mm_setzero_ps();
661         fix1             = _mm_setzero_ps();
662         fiy1             = _mm_setzero_ps();
663         fiz1             = _mm_setzero_ps();
664         fix2             = _mm_setzero_ps();
665         fiy2             = _mm_setzero_ps();
666         fiz2             = _mm_setzero_ps();
667
668         /* Start inner kernel loop */
669         for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
670         {
671
672             /* Get j neighbor index, and coordinate index */
673             jnrA             = jjnr[jidx];
674             jnrB             = jjnr[jidx+1];
675             jnrC             = jjnr[jidx+2];
676             jnrD             = jjnr[jidx+3];
677             j_coord_offsetA  = DIM*jnrA;
678             j_coord_offsetB  = DIM*jnrB;
679             j_coord_offsetC  = DIM*jnrC;
680             j_coord_offsetD  = DIM*jnrD;
681
682             /* load j atom coordinates */
683             gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
684                                               x+j_coord_offsetC,x+j_coord_offsetD,
685                                               &jx0,&jy0,&jz0);
686
687             /* Calculate displacement vector */
688             dx00             = _mm_sub_ps(ix0,jx0);
689             dy00             = _mm_sub_ps(iy0,jy0);
690             dz00             = _mm_sub_ps(iz0,jz0);
691             dx10             = _mm_sub_ps(ix1,jx0);
692             dy10             = _mm_sub_ps(iy1,jy0);
693             dz10             = _mm_sub_ps(iz1,jz0);
694             dx20             = _mm_sub_ps(ix2,jx0);
695             dy20             = _mm_sub_ps(iy2,jy0);
696             dz20             = _mm_sub_ps(iz2,jz0);
697
698             /* Calculate squared distance and things based on it */
699             rsq00            = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
700             rsq10            = gmx_mm_calc_rsq_ps(dx10,dy10,dz10);
701             rsq20            = gmx_mm_calc_rsq_ps(dx20,dy20,dz20);
702
703             rinv00           = gmx_mm_invsqrt_ps(rsq00);
704             rinv10           = gmx_mm_invsqrt_ps(rsq10);
705             rinv20           = gmx_mm_invsqrt_ps(rsq20);
706
707             rinvsq00         = _mm_mul_ps(rinv00,rinv00);
708             rinvsq10         = _mm_mul_ps(rinv10,rinv10);
709             rinvsq20         = _mm_mul_ps(rinv20,rinv20);
710
711             /* Load parameters for j particles */
712             jq0              = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
713                                                               charge+jnrC+0,charge+jnrD+0);
714
715             fjx0             = _mm_setzero_ps();
716             fjy0             = _mm_setzero_ps();
717             fjz0             = _mm_setzero_ps();
718
719             /**************************
720              * CALCULATE INTERACTIONS *
721              **************************/
722
723             if (gmx_mm_any_lt(rsq00,rcutoff2))
724             {
725
726             /* Compute parameters for interactions between i and j atoms */
727             qq00             = _mm_mul_ps(iq0,jq0);
728
729             /* REACTION-FIELD ELECTROSTATICS */
730             felec            = _mm_mul_ps(qq00,_mm_msub_ps(rinv00,rinvsq00,krf2));
731
732             cutoff_mask      = _mm_cmplt_ps(rsq00,rcutoff2);
733
734             fscal            = felec;
735
736             fscal            = _mm_and_ps(fscal,cutoff_mask);
737
738              /* Update vectorial force */
739             fix0             = _mm_macc_ps(dx00,fscal,fix0);
740             fiy0             = _mm_macc_ps(dy00,fscal,fiy0);
741             fiz0             = _mm_macc_ps(dz00,fscal,fiz0);
742
743             fjx0             = _mm_macc_ps(dx00,fscal,fjx0);
744             fjy0             = _mm_macc_ps(dy00,fscal,fjy0);
745             fjz0             = _mm_macc_ps(dz00,fscal,fjz0);
746
747             }
748
749             /**************************
750              * CALCULATE INTERACTIONS *
751              **************************/
752
753             if (gmx_mm_any_lt(rsq10,rcutoff2))
754             {
755
756             /* Compute parameters for interactions between i and j atoms */
757             qq10             = _mm_mul_ps(iq1,jq0);
758
759             /* REACTION-FIELD ELECTROSTATICS */
760             felec            = _mm_mul_ps(qq10,_mm_msub_ps(rinv10,rinvsq10,krf2));
761
762             cutoff_mask      = _mm_cmplt_ps(rsq10,rcutoff2);
763
764             fscal            = felec;
765
766             fscal            = _mm_and_ps(fscal,cutoff_mask);
767
768              /* Update vectorial force */
769             fix1             = _mm_macc_ps(dx10,fscal,fix1);
770             fiy1             = _mm_macc_ps(dy10,fscal,fiy1);
771             fiz1             = _mm_macc_ps(dz10,fscal,fiz1);
772
773             fjx0             = _mm_macc_ps(dx10,fscal,fjx0);
774             fjy0             = _mm_macc_ps(dy10,fscal,fjy0);
775             fjz0             = _mm_macc_ps(dz10,fscal,fjz0);
776
777             }
778
779             /**************************
780              * CALCULATE INTERACTIONS *
781              **************************/
782
783             if (gmx_mm_any_lt(rsq20,rcutoff2))
784             {
785
786             /* Compute parameters for interactions between i and j atoms */
787             qq20             = _mm_mul_ps(iq2,jq0);
788
789             /* REACTION-FIELD ELECTROSTATICS */
790             felec            = _mm_mul_ps(qq20,_mm_msub_ps(rinv20,rinvsq20,krf2));
791
792             cutoff_mask      = _mm_cmplt_ps(rsq20,rcutoff2);
793
794             fscal            = felec;
795
796             fscal            = _mm_and_ps(fscal,cutoff_mask);
797
798              /* Update vectorial force */
799             fix2             = _mm_macc_ps(dx20,fscal,fix2);
800             fiy2             = _mm_macc_ps(dy20,fscal,fiy2);
801             fiz2             = _mm_macc_ps(dz20,fscal,fiz2);
802
803             fjx0             = _mm_macc_ps(dx20,fscal,fjx0);
804             fjy0             = _mm_macc_ps(dy20,fscal,fjy0);
805             fjz0             = _mm_macc_ps(dz20,fscal,fjz0);
806
807             }
808
809             fjptrA             = f+j_coord_offsetA;
810             fjptrB             = f+j_coord_offsetB;
811             fjptrC             = f+j_coord_offsetC;
812             fjptrD             = f+j_coord_offsetD;
813
814             gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjx0,fjy0,fjz0);
815
816             /* Inner loop uses 99 flops */
817         }
818
819         if(jidx<j_index_end)
820         {
821
822             /* Get j neighbor index, and coordinate index */
823             jnrlistA         = jjnr[jidx];
824             jnrlistB         = jjnr[jidx+1];
825             jnrlistC         = jjnr[jidx+2];
826             jnrlistD         = jjnr[jidx+3];
827             /* Sign of each element will be negative for non-real atoms.
828              * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
829              * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
830              */
831             dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
832             jnrA       = (jnrlistA>=0) ? jnrlistA : 0;
833             jnrB       = (jnrlistB>=0) ? jnrlistB : 0;
834             jnrC       = (jnrlistC>=0) ? jnrlistC : 0;
835             jnrD       = (jnrlistD>=0) ? jnrlistD : 0;
836             j_coord_offsetA  = DIM*jnrA;
837             j_coord_offsetB  = DIM*jnrB;
838             j_coord_offsetC  = DIM*jnrC;
839             j_coord_offsetD  = DIM*jnrD;
840
841             /* load j atom coordinates */
842             gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
843                                               x+j_coord_offsetC,x+j_coord_offsetD,
844                                               &jx0,&jy0,&jz0);
845
846             /* Calculate displacement vector */
847             dx00             = _mm_sub_ps(ix0,jx0);
848             dy00             = _mm_sub_ps(iy0,jy0);
849             dz00             = _mm_sub_ps(iz0,jz0);
850             dx10             = _mm_sub_ps(ix1,jx0);
851             dy10             = _mm_sub_ps(iy1,jy0);
852             dz10             = _mm_sub_ps(iz1,jz0);
853             dx20             = _mm_sub_ps(ix2,jx0);
854             dy20             = _mm_sub_ps(iy2,jy0);
855             dz20             = _mm_sub_ps(iz2,jz0);
856
857             /* Calculate squared distance and things based on it */
858             rsq00            = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
859             rsq10            = gmx_mm_calc_rsq_ps(dx10,dy10,dz10);
860             rsq20            = gmx_mm_calc_rsq_ps(dx20,dy20,dz20);
861
862             rinv00           = gmx_mm_invsqrt_ps(rsq00);
863             rinv10           = gmx_mm_invsqrt_ps(rsq10);
864             rinv20           = gmx_mm_invsqrt_ps(rsq20);
865
866             rinvsq00         = _mm_mul_ps(rinv00,rinv00);
867             rinvsq10         = _mm_mul_ps(rinv10,rinv10);
868             rinvsq20         = _mm_mul_ps(rinv20,rinv20);
869
870             /* Load parameters for j particles */
871             jq0              = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
872                                                               charge+jnrC+0,charge+jnrD+0);
873
874             fjx0             = _mm_setzero_ps();
875             fjy0             = _mm_setzero_ps();
876             fjz0             = _mm_setzero_ps();
877
878             /**************************
879              * CALCULATE INTERACTIONS *
880              **************************/
881
882             if (gmx_mm_any_lt(rsq00,rcutoff2))
883             {
884
885             /* Compute parameters for interactions between i and j atoms */
886             qq00             = _mm_mul_ps(iq0,jq0);
887
888             /* REACTION-FIELD ELECTROSTATICS */
889             felec            = _mm_mul_ps(qq00,_mm_msub_ps(rinv00,rinvsq00,krf2));
890
891             cutoff_mask      = _mm_cmplt_ps(rsq00,rcutoff2);
892
893             fscal            = felec;
894
895             fscal            = _mm_and_ps(fscal,cutoff_mask);
896
897             fscal            = _mm_andnot_ps(dummy_mask,fscal);
898
899              /* Update vectorial force */
900             fix0             = _mm_macc_ps(dx00,fscal,fix0);
901             fiy0             = _mm_macc_ps(dy00,fscal,fiy0);
902             fiz0             = _mm_macc_ps(dz00,fscal,fiz0);
903
904             fjx0             = _mm_macc_ps(dx00,fscal,fjx0);
905             fjy0             = _mm_macc_ps(dy00,fscal,fjy0);
906             fjz0             = _mm_macc_ps(dz00,fscal,fjz0);
907
908             }
909
910             /**************************
911              * CALCULATE INTERACTIONS *
912              **************************/
913
914             if (gmx_mm_any_lt(rsq10,rcutoff2))
915             {
916
917             /* Compute parameters for interactions between i and j atoms */
918             qq10             = _mm_mul_ps(iq1,jq0);
919
920             /* REACTION-FIELD ELECTROSTATICS */
921             felec            = _mm_mul_ps(qq10,_mm_msub_ps(rinv10,rinvsq10,krf2));
922
923             cutoff_mask      = _mm_cmplt_ps(rsq10,rcutoff2);
924
925             fscal            = felec;
926
927             fscal            = _mm_and_ps(fscal,cutoff_mask);
928
929             fscal            = _mm_andnot_ps(dummy_mask,fscal);
930
931              /* Update vectorial force */
932             fix1             = _mm_macc_ps(dx10,fscal,fix1);
933             fiy1             = _mm_macc_ps(dy10,fscal,fiy1);
934             fiz1             = _mm_macc_ps(dz10,fscal,fiz1);
935
936             fjx0             = _mm_macc_ps(dx10,fscal,fjx0);
937             fjy0             = _mm_macc_ps(dy10,fscal,fjy0);
938             fjz0             = _mm_macc_ps(dz10,fscal,fjz0);
939
940             }
941
942             /**************************
943              * CALCULATE INTERACTIONS *
944              **************************/
945
946             if (gmx_mm_any_lt(rsq20,rcutoff2))
947             {
948
949             /* Compute parameters for interactions between i and j atoms */
950             qq20             = _mm_mul_ps(iq2,jq0);
951
952             /* REACTION-FIELD ELECTROSTATICS */
953             felec            = _mm_mul_ps(qq20,_mm_msub_ps(rinv20,rinvsq20,krf2));
954
955             cutoff_mask      = _mm_cmplt_ps(rsq20,rcutoff2);
956
957             fscal            = felec;
958
959             fscal            = _mm_and_ps(fscal,cutoff_mask);
960
961             fscal            = _mm_andnot_ps(dummy_mask,fscal);
962
963              /* Update vectorial force */
964             fix2             = _mm_macc_ps(dx20,fscal,fix2);
965             fiy2             = _mm_macc_ps(dy20,fscal,fiy2);
966             fiz2             = _mm_macc_ps(dz20,fscal,fiz2);
967
968             fjx0             = _mm_macc_ps(dx20,fscal,fjx0);
969             fjy0             = _mm_macc_ps(dy20,fscal,fjy0);
970             fjz0             = _mm_macc_ps(dz20,fscal,fjz0);
971
972             }
973
974             fjptrA             = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
975             fjptrB             = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
976             fjptrC             = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
977             fjptrD             = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
978
979             gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjx0,fjy0,fjz0);
980
981             /* Inner loop uses 99 flops */
982         }
983
984         /* End of innermost loop */
985
986         gmx_mm_update_iforce_3atom_swizzle_ps(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,
987                                               f+i_coord_offset,fshift+i_shift_offset);
988
989         /* Increment number of inner iterations */
990         inneriter                  += j_index_end - j_index_start;
991
992         /* Outer loop uses 18 flops */
993     }
994
995     /* Increment number of outer iterations */
996     outeriter        += nri;
997
998     /* Update outer/inner flops */
999
1000     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_W3_F,outeriter*18 + inneriter*99);
1001 }