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