Merge release-4-6 into master
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_avx_256_double / nb_kernel_ElecEw_VdwNone_GeomP1P1_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_ElecEw_VdwNone_GeomP1P1_VF_avx_256_double
38  * Electrostatics interaction: Ewald
39  * VdW interaction:            None
40  * Geometry:                   Particle-Particle
41  * Calculate force/pot:        PotentialAndForce
42  */
43 void
44 nb_kernel_ElecEw_VdwNone_GeomP1P1_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     int              vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
73     __m256d          jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
74     __m256d          dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
75     __m256d          velec,felec,velecsum,facel,crf,krf,krf2;
76     real             *charge;
77     __m128i          ewitab;
78     __m256d          ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
79     __m256d          beta,beta2,beta3,zeta2,pmecorrF,pmecorrV,rinv3;
80     real             *ewtab;
81     __m256d          dummy_mask,cutoff_mask;
82     __m128           tmpmask0,tmpmask1;
83     __m256d          signbit = _mm256_castsi256_pd( _mm256_set1_epi32(0x80000000) );
84     __m256d          one     = _mm256_set1_pd(1.0);
85     __m256d          two     = _mm256_set1_pd(2.0);
86     x                = xx[0];
87     f                = ff[0];
88
89     nri              = nlist->nri;
90     iinr             = nlist->iinr;
91     jindex           = nlist->jindex;
92     jjnr             = nlist->jjnr;
93     shiftidx         = nlist->shift;
94     gid              = nlist->gid;
95     shiftvec         = fr->shift_vec[0];
96     fshift           = fr->fshift[0];
97     facel            = _mm256_set1_pd(fr->epsfac);
98     charge           = mdatoms->chargeA;
99
100     sh_ewald         = _mm256_set1_pd(fr->ic->sh_ewald);
101     beta             = _mm256_set1_pd(fr->ic->ewaldcoeff);
102     beta2            = _mm256_mul_pd(beta,beta);
103     beta3            = _mm256_mul_pd(beta,beta2);
104
105     ewtab            = fr->ic->tabq_coul_FDV0;
106     ewtabscale       = _mm256_set1_pd(fr->ic->tabq_scale);
107     ewtabhalfspace   = _mm256_set1_pd(0.5/fr->ic->tabq_scale);
108
109     /* Avoid stupid compiler warnings */
110     jnrA = jnrB = jnrC = jnrD = 0;
111     j_coord_offsetA = 0;
112     j_coord_offsetB = 0;
113     j_coord_offsetC = 0;
114     j_coord_offsetD = 0;
115
116     outeriter        = 0;
117     inneriter        = 0;
118
119     for(iidx=0;iidx<4*DIM;iidx++)
120     {
121         scratch[iidx] = 0.0;
122     }
123
124     /* Start outer loop over neighborlists */
125     for(iidx=0; iidx<nri; iidx++)
126     {
127         /* Load shift vector for this list */
128         i_shift_offset   = DIM*shiftidx[iidx];
129
130         /* Load limits for loop over neighbors */
131         j_index_start    = jindex[iidx];
132         j_index_end      = jindex[iidx+1];
133
134         /* Get outer coordinate index */
135         inr              = iinr[iidx];
136         i_coord_offset   = DIM*inr;
137
138         /* Load i particle coords and add shift vector */
139         gmx_mm256_load_shift_and_1rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
140
141         fix0             = _mm256_setzero_pd();
142         fiy0             = _mm256_setzero_pd();
143         fiz0             = _mm256_setzero_pd();
144
145         /* Load parameters for i particles */
146         iq0              = _mm256_mul_pd(facel,_mm256_set1_pd(charge[inr+0]));
147
148         /* Reset potential sums */
149         velecsum         = _mm256_setzero_pd();
150
151         /* Start inner kernel loop */
152         for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
153         {
154
155             /* Get j neighbor index, and coordinate index */
156             jnrA             = jjnr[jidx];
157             jnrB             = jjnr[jidx+1];
158             jnrC             = jjnr[jidx+2];
159             jnrD             = jjnr[jidx+3];
160             j_coord_offsetA  = DIM*jnrA;
161             j_coord_offsetB  = DIM*jnrB;
162             j_coord_offsetC  = DIM*jnrC;
163             j_coord_offsetD  = DIM*jnrD;
164
165             /* load j atom coordinates */
166             gmx_mm256_load_1rvec_4ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
167                                                  x+j_coord_offsetC,x+j_coord_offsetD,
168                                                  &jx0,&jy0,&jz0);
169
170             /* Calculate displacement vector */
171             dx00             = _mm256_sub_pd(ix0,jx0);
172             dy00             = _mm256_sub_pd(iy0,jy0);
173             dz00             = _mm256_sub_pd(iz0,jz0);
174
175             /* Calculate squared distance and things based on it */
176             rsq00            = gmx_mm256_calc_rsq_pd(dx00,dy00,dz00);
177
178             rinv00           = gmx_mm256_invsqrt_pd(rsq00);
179
180             rinvsq00         = _mm256_mul_pd(rinv00,rinv00);
181
182             /* Load parameters for j particles */
183             jq0              = gmx_mm256_load_4real_swizzle_pd(charge+jnrA+0,charge+jnrB+0,
184                                                                  charge+jnrC+0,charge+jnrD+0);
185
186             /**************************
187              * CALCULATE INTERACTIONS *
188              **************************/
189
190             r00              = _mm256_mul_pd(rsq00,rinv00);
191
192             /* Compute parameters for interactions between i and j atoms */
193             qq00             = _mm256_mul_pd(iq0,jq0);
194
195             /* EWALD ELECTROSTATICS */
196
197             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
198             ewrt             = _mm256_mul_pd(r00,ewtabscale);
199             ewitab           = _mm256_cvttpd_epi32(ewrt);
200             eweps            = _mm256_sub_pd(ewrt,_mm256_round_pd(ewrt, _MM_FROUND_FLOOR));
201             ewitab           = _mm_slli_epi32(ewitab,2);
202             ewtabF           = _mm256_load_pd( ewtab + _mm_extract_epi32(ewitab,0) );
203             ewtabD           = _mm256_load_pd( ewtab + _mm_extract_epi32(ewitab,1) );
204             ewtabV           = _mm256_load_pd( ewtab + _mm_extract_epi32(ewitab,2) );
205             ewtabFn          = _mm256_load_pd( ewtab + _mm_extract_epi32(ewitab,3) );
206             GMX_MM256_FULLTRANSPOSE4_PD(ewtabF,ewtabD,ewtabV,ewtabFn);
207             felec            = _mm256_add_pd(ewtabF,_mm256_mul_pd(eweps,ewtabD));
208             velec            = _mm256_sub_pd(ewtabV,_mm256_mul_pd(_mm256_mul_pd(ewtabhalfspace,eweps),_mm256_add_pd(ewtabF,felec)));
209             velec            = _mm256_mul_pd(qq00,_mm256_sub_pd(rinv00,velec));
210             felec            = _mm256_mul_pd(_mm256_mul_pd(qq00,rinv00),_mm256_sub_pd(rinvsq00,felec));
211
212             /* Update potential sum for this i atom from the interaction with this j atom. */
213             velecsum         = _mm256_add_pd(velecsum,velec);
214
215             fscal            = felec;
216
217             /* Calculate temporary vectorial force */
218             tx               = _mm256_mul_pd(fscal,dx00);
219             ty               = _mm256_mul_pd(fscal,dy00);
220             tz               = _mm256_mul_pd(fscal,dz00);
221
222             /* Update vectorial force */
223             fix0             = _mm256_add_pd(fix0,tx);
224             fiy0             = _mm256_add_pd(fiy0,ty);
225             fiz0             = _mm256_add_pd(fiz0,tz);
226
227             fjptrA             = f+j_coord_offsetA;
228             fjptrB             = f+j_coord_offsetB;
229             fjptrC             = f+j_coord_offsetC;
230             fjptrD             = f+j_coord_offsetD;
231             gmx_mm256_decrement_1rvec_4ptr_swizzle_pd(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
232
233             /* Inner loop uses 41 flops */
234         }
235
236         if(jidx<j_index_end)
237         {
238
239             /* Get j neighbor index, and coordinate index */
240             jnrlistA         = jjnr[jidx];
241             jnrlistB         = jjnr[jidx+1];
242             jnrlistC         = jjnr[jidx+2];
243             jnrlistD         = jjnr[jidx+3];
244             /* Sign of each element will be negative for non-real atoms.
245              * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
246              * so use it as val = _mm_andnot_pd(mask,val) to clear dummy entries.
247              */
248             tmpmask0 = gmx_mm_castsi128_pd(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
249
250             tmpmask1 = _mm_permute_ps(tmpmask0,_GMX_MM_PERMUTE(3,3,2,2));
251             tmpmask0 = _mm_permute_ps(tmpmask0,_GMX_MM_PERMUTE(1,1,0,0));
252             dummy_mask = _mm256_castps_pd(gmx_mm256_set_m128(tmpmask1,tmpmask0));
253
254             jnrA       = (jnrlistA>=0) ? jnrlistA : 0;
255             jnrB       = (jnrlistB>=0) ? jnrlistB : 0;
256             jnrC       = (jnrlistC>=0) ? jnrlistC : 0;
257             jnrD       = (jnrlistD>=0) ? jnrlistD : 0;
258             j_coord_offsetA  = DIM*jnrA;
259             j_coord_offsetB  = DIM*jnrB;
260             j_coord_offsetC  = DIM*jnrC;
261             j_coord_offsetD  = DIM*jnrD;
262
263             /* load j atom coordinates */
264             gmx_mm256_load_1rvec_4ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
265                                                  x+j_coord_offsetC,x+j_coord_offsetD,
266                                                  &jx0,&jy0,&jz0);
267
268             /* Calculate displacement vector */
269             dx00             = _mm256_sub_pd(ix0,jx0);
270             dy00             = _mm256_sub_pd(iy0,jy0);
271             dz00             = _mm256_sub_pd(iz0,jz0);
272
273             /* Calculate squared distance and things based on it */
274             rsq00            = gmx_mm256_calc_rsq_pd(dx00,dy00,dz00);
275
276             rinv00           = gmx_mm256_invsqrt_pd(rsq00);
277
278             rinvsq00         = _mm256_mul_pd(rinv00,rinv00);
279
280             /* Load parameters for j particles */
281             jq0              = gmx_mm256_load_4real_swizzle_pd(charge+jnrA+0,charge+jnrB+0,
282                                                                  charge+jnrC+0,charge+jnrD+0);
283
284             /**************************
285              * CALCULATE INTERACTIONS *
286              **************************/
287
288             r00              = _mm256_mul_pd(rsq00,rinv00);
289             r00              = _mm256_andnot_pd(dummy_mask,r00);
290
291             /* Compute parameters for interactions between i and j atoms */
292             qq00             = _mm256_mul_pd(iq0,jq0);
293
294             /* EWALD ELECTROSTATICS */
295
296             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
297             ewrt             = _mm256_mul_pd(r00,ewtabscale);
298             ewitab           = _mm256_cvttpd_epi32(ewrt);
299             eweps            = _mm256_sub_pd(ewrt,_mm256_round_pd(ewrt, _MM_FROUND_FLOOR));
300             ewitab           = _mm_slli_epi32(ewitab,2);
301             ewtabF           = _mm256_load_pd( ewtab + _mm_extract_epi32(ewitab,0) );
302             ewtabD           = _mm256_load_pd( ewtab + _mm_extract_epi32(ewitab,1) );
303             ewtabV           = _mm256_load_pd( ewtab + _mm_extract_epi32(ewitab,2) );
304             ewtabFn          = _mm256_load_pd( ewtab + _mm_extract_epi32(ewitab,3) );
305             GMX_MM256_FULLTRANSPOSE4_PD(ewtabF,ewtabD,ewtabV,ewtabFn);
306             felec            = _mm256_add_pd(ewtabF,_mm256_mul_pd(eweps,ewtabD));
307             velec            = _mm256_sub_pd(ewtabV,_mm256_mul_pd(_mm256_mul_pd(ewtabhalfspace,eweps),_mm256_add_pd(ewtabF,felec)));
308             velec            = _mm256_mul_pd(qq00,_mm256_sub_pd(rinv00,velec));
309             felec            = _mm256_mul_pd(_mm256_mul_pd(qq00,rinv00),_mm256_sub_pd(rinvsq00,felec));
310
311             /* Update potential sum for this i atom from the interaction with this j atom. */
312             velec            = _mm256_andnot_pd(dummy_mask,velec);
313             velecsum         = _mm256_add_pd(velecsum,velec);
314
315             fscal            = felec;
316
317             fscal            = _mm256_andnot_pd(dummy_mask,fscal);
318
319             /* Calculate temporary vectorial force */
320             tx               = _mm256_mul_pd(fscal,dx00);
321             ty               = _mm256_mul_pd(fscal,dy00);
322             tz               = _mm256_mul_pd(fscal,dz00);
323
324             /* Update vectorial force */
325             fix0             = _mm256_add_pd(fix0,tx);
326             fiy0             = _mm256_add_pd(fiy0,ty);
327             fiz0             = _mm256_add_pd(fiz0,tz);
328
329             fjptrA             = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
330             fjptrB             = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
331             fjptrC             = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
332             fjptrD             = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
333             gmx_mm256_decrement_1rvec_4ptr_swizzle_pd(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
334
335             /* Inner loop uses 42 flops */
336         }
337
338         /* End of innermost loop */
339
340         gmx_mm256_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
341                                                  f+i_coord_offset,fshift+i_shift_offset);
342
343         ggid                        = gid[iidx];
344         /* Update potential energies */
345         gmx_mm256_update_1pot_pd(velecsum,kernel_data->energygrp_elec+ggid);
346
347         /* Increment number of inner iterations */
348         inneriter                  += j_index_end - j_index_start;
349
350         /* Outer loop uses 8 flops */
351     }
352
353     /* Increment number of outer iterations */
354     outeriter        += nri;
355
356     /* Update outer/inner flops */
357
358     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VF,outeriter*8 + inneriter*42);
359 }
360 /*
361  * Gromacs nonbonded kernel:   nb_kernel_ElecEw_VdwNone_GeomP1P1_F_avx_256_double
362  * Electrostatics interaction: Ewald
363  * VdW interaction:            None
364  * Geometry:                   Particle-Particle
365  * Calculate force/pot:        Force
366  */
367 void
368 nb_kernel_ElecEw_VdwNone_GeomP1P1_F_avx_256_double
369                     (t_nblist * gmx_restrict                nlist,
370                      rvec * gmx_restrict                    xx,
371                      rvec * gmx_restrict                    ff,
372                      t_forcerec * gmx_restrict              fr,
373                      t_mdatoms * gmx_restrict               mdatoms,
374                      nb_kernel_data_t * gmx_restrict        kernel_data,
375                      t_nrnb * gmx_restrict                  nrnb)
376 {
377     /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or 
378      * just 0 for non-waters.
379      * Suffixes A,B,C,D refer to j loop unrolling done with AVX, e.g. for the four different
380      * jnr indices corresponding to data put in the four positions in the SIMD register.
381      */
382     int              i_shift_offset,i_coord_offset,outeriter,inneriter;
383     int              j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
384     int              jnrA,jnrB,jnrC,jnrD;
385     int              jnrlistA,jnrlistB,jnrlistC,jnrlistD;
386     int              jnrlistE,jnrlistF,jnrlistG,jnrlistH;
387     int              j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
388     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
389     real             rcutoff_scalar;
390     real             *shiftvec,*fshift,*x,*f;
391     real             *fjptrA,*fjptrB,*fjptrC,*fjptrD;
392     real             scratch[4*DIM];
393     __m256d          tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
394     real *           vdwioffsetptr0;
395     __m256d          ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
396     int              vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
397     __m256d          jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
398     __m256d          dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
399     __m256d          velec,felec,velecsum,facel,crf,krf,krf2;
400     real             *charge;
401     __m128i          ewitab;
402     __m256d          ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
403     __m256d          beta,beta2,beta3,zeta2,pmecorrF,pmecorrV,rinv3;
404     real             *ewtab;
405     __m256d          dummy_mask,cutoff_mask;
406     __m128           tmpmask0,tmpmask1;
407     __m256d          signbit = _mm256_castsi256_pd( _mm256_set1_epi32(0x80000000) );
408     __m256d          one     = _mm256_set1_pd(1.0);
409     __m256d          two     = _mm256_set1_pd(2.0);
410     x                = xx[0];
411     f                = ff[0];
412
413     nri              = nlist->nri;
414     iinr             = nlist->iinr;
415     jindex           = nlist->jindex;
416     jjnr             = nlist->jjnr;
417     shiftidx         = nlist->shift;
418     gid              = nlist->gid;
419     shiftvec         = fr->shift_vec[0];
420     fshift           = fr->fshift[0];
421     facel            = _mm256_set1_pd(fr->epsfac);
422     charge           = mdatoms->chargeA;
423
424     sh_ewald         = _mm256_set1_pd(fr->ic->sh_ewald);
425     beta             = _mm256_set1_pd(fr->ic->ewaldcoeff);
426     beta2            = _mm256_mul_pd(beta,beta);
427     beta3            = _mm256_mul_pd(beta,beta2);
428
429     ewtab            = fr->ic->tabq_coul_F;
430     ewtabscale       = _mm256_set1_pd(fr->ic->tabq_scale);
431     ewtabhalfspace   = _mm256_set1_pd(0.5/fr->ic->tabq_scale);
432
433     /* Avoid stupid compiler warnings */
434     jnrA = jnrB = jnrC = jnrD = 0;
435     j_coord_offsetA = 0;
436     j_coord_offsetB = 0;
437     j_coord_offsetC = 0;
438     j_coord_offsetD = 0;
439
440     outeriter        = 0;
441     inneriter        = 0;
442
443     for(iidx=0;iidx<4*DIM;iidx++)
444     {
445         scratch[iidx] = 0.0;
446     }
447
448     /* Start outer loop over neighborlists */
449     for(iidx=0; iidx<nri; iidx++)
450     {
451         /* Load shift vector for this list */
452         i_shift_offset   = DIM*shiftidx[iidx];
453
454         /* Load limits for loop over neighbors */
455         j_index_start    = jindex[iidx];
456         j_index_end      = jindex[iidx+1];
457
458         /* Get outer coordinate index */
459         inr              = iinr[iidx];
460         i_coord_offset   = DIM*inr;
461
462         /* Load i particle coords and add shift vector */
463         gmx_mm256_load_shift_and_1rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
464
465         fix0             = _mm256_setzero_pd();
466         fiy0             = _mm256_setzero_pd();
467         fiz0             = _mm256_setzero_pd();
468
469         /* Load parameters for i particles */
470         iq0              = _mm256_mul_pd(facel,_mm256_set1_pd(charge[inr+0]));
471
472         /* Start inner kernel loop */
473         for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
474         {
475
476             /* Get j neighbor index, and coordinate index */
477             jnrA             = jjnr[jidx];
478             jnrB             = jjnr[jidx+1];
479             jnrC             = jjnr[jidx+2];
480             jnrD             = jjnr[jidx+3];
481             j_coord_offsetA  = DIM*jnrA;
482             j_coord_offsetB  = DIM*jnrB;
483             j_coord_offsetC  = DIM*jnrC;
484             j_coord_offsetD  = DIM*jnrD;
485
486             /* load j atom coordinates */
487             gmx_mm256_load_1rvec_4ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
488                                                  x+j_coord_offsetC,x+j_coord_offsetD,
489                                                  &jx0,&jy0,&jz0);
490
491             /* Calculate displacement vector */
492             dx00             = _mm256_sub_pd(ix0,jx0);
493             dy00             = _mm256_sub_pd(iy0,jy0);
494             dz00             = _mm256_sub_pd(iz0,jz0);
495
496             /* Calculate squared distance and things based on it */
497             rsq00            = gmx_mm256_calc_rsq_pd(dx00,dy00,dz00);
498
499             rinv00           = gmx_mm256_invsqrt_pd(rsq00);
500
501             rinvsq00         = _mm256_mul_pd(rinv00,rinv00);
502
503             /* Load parameters for j particles */
504             jq0              = gmx_mm256_load_4real_swizzle_pd(charge+jnrA+0,charge+jnrB+0,
505                                                                  charge+jnrC+0,charge+jnrD+0);
506
507             /**************************
508              * CALCULATE INTERACTIONS *
509              **************************/
510
511             r00              = _mm256_mul_pd(rsq00,rinv00);
512
513             /* Compute parameters for interactions between i and j atoms */
514             qq00             = _mm256_mul_pd(iq0,jq0);
515
516             /* EWALD ELECTROSTATICS */
517
518             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
519             ewrt             = _mm256_mul_pd(r00,ewtabscale);
520             ewitab           = _mm256_cvttpd_epi32(ewrt);
521             eweps            = _mm256_sub_pd(ewrt,_mm256_round_pd(ewrt, _MM_FROUND_FLOOR));
522             gmx_mm256_load_4pair_swizzle_pd(ewtab + _mm_extract_epi32(ewitab,0),ewtab + _mm_extract_epi32(ewitab,1),
523                                             ewtab + _mm_extract_epi32(ewitab,2),ewtab + _mm_extract_epi32(ewitab,3),
524                                             &ewtabF,&ewtabFn);
525             felec            = _mm256_add_pd(_mm256_mul_pd( _mm256_sub_pd(one,eweps),ewtabF),_mm256_mul_pd(eweps,ewtabFn));
526             felec            = _mm256_mul_pd(_mm256_mul_pd(qq00,rinv00),_mm256_sub_pd(rinvsq00,felec));
527
528             fscal            = felec;
529
530             /* Calculate temporary vectorial force */
531             tx               = _mm256_mul_pd(fscal,dx00);
532             ty               = _mm256_mul_pd(fscal,dy00);
533             tz               = _mm256_mul_pd(fscal,dz00);
534
535             /* Update vectorial force */
536             fix0             = _mm256_add_pd(fix0,tx);
537             fiy0             = _mm256_add_pd(fiy0,ty);
538             fiz0             = _mm256_add_pd(fiz0,tz);
539
540             fjptrA             = f+j_coord_offsetA;
541             fjptrB             = f+j_coord_offsetB;
542             fjptrC             = f+j_coord_offsetC;
543             fjptrD             = f+j_coord_offsetD;
544             gmx_mm256_decrement_1rvec_4ptr_swizzle_pd(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
545
546             /* Inner loop uses 36 flops */
547         }
548
549         if(jidx<j_index_end)
550         {
551
552             /* Get j neighbor index, and coordinate index */
553             jnrlistA         = jjnr[jidx];
554             jnrlistB         = jjnr[jidx+1];
555             jnrlistC         = jjnr[jidx+2];
556             jnrlistD         = jjnr[jidx+3];
557             /* Sign of each element will be negative for non-real atoms.
558              * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
559              * so use it as val = _mm_andnot_pd(mask,val) to clear dummy entries.
560              */
561             tmpmask0 = gmx_mm_castsi128_pd(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
562
563             tmpmask1 = _mm_permute_ps(tmpmask0,_GMX_MM_PERMUTE(3,3,2,2));
564             tmpmask0 = _mm_permute_ps(tmpmask0,_GMX_MM_PERMUTE(1,1,0,0));
565             dummy_mask = _mm256_castps_pd(gmx_mm256_set_m128(tmpmask1,tmpmask0));
566
567             jnrA       = (jnrlistA>=0) ? jnrlistA : 0;
568             jnrB       = (jnrlistB>=0) ? jnrlistB : 0;
569             jnrC       = (jnrlistC>=0) ? jnrlistC : 0;
570             jnrD       = (jnrlistD>=0) ? jnrlistD : 0;
571             j_coord_offsetA  = DIM*jnrA;
572             j_coord_offsetB  = DIM*jnrB;
573             j_coord_offsetC  = DIM*jnrC;
574             j_coord_offsetD  = DIM*jnrD;
575
576             /* load j atom coordinates */
577             gmx_mm256_load_1rvec_4ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
578                                                  x+j_coord_offsetC,x+j_coord_offsetD,
579                                                  &jx0,&jy0,&jz0);
580
581             /* Calculate displacement vector */
582             dx00             = _mm256_sub_pd(ix0,jx0);
583             dy00             = _mm256_sub_pd(iy0,jy0);
584             dz00             = _mm256_sub_pd(iz0,jz0);
585
586             /* Calculate squared distance and things based on it */
587             rsq00            = gmx_mm256_calc_rsq_pd(dx00,dy00,dz00);
588
589             rinv00           = gmx_mm256_invsqrt_pd(rsq00);
590
591             rinvsq00         = _mm256_mul_pd(rinv00,rinv00);
592
593             /* Load parameters for j particles */
594             jq0              = gmx_mm256_load_4real_swizzle_pd(charge+jnrA+0,charge+jnrB+0,
595                                                                  charge+jnrC+0,charge+jnrD+0);
596
597             /**************************
598              * CALCULATE INTERACTIONS *
599              **************************/
600
601             r00              = _mm256_mul_pd(rsq00,rinv00);
602             r00              = _mm256_andnot_pd(dummy_mask,r00);
603
604             /* Compute parameters for interactions between i and j atoms */
605             qq00             = _mm256_mul_pd(iq0,jq0);
606
607             /* EWALD ELECTROSTATICS */
608
609             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
610             ewrt             = _mm256_mul_pd(r00,ewtabscale);
611             ewitab           = _mm256_cvttpd_epi32(ewrt);
612             eweps            = _mm256_sub_pd(ewrt,_mm256_round_pd(ewrt, _MM_FROUND_FLOOR));
613             gmx_mm256_load_4pair_swizzle_pd(ewtab + _mm_extract_epi32(ewitab,0),ewtab + _mm_extract_epi32(ewitab,1),
614                                             ewtab + _mm_extract_epi32(ewitab,2),ewtab + _mm_extract_epi32(ewitab,3),
615                                             &ewtabF,&ewtabFn);
616             felec            = _mm256_add_pd(_mm256_mul_pd( _mm256_sub_pd(one,eweps),ewtabF),_mm256_mul_pd(eweps,ewtabFn));
617             felec            = _mm256_mul_pd(_mm256_mul_pd(qq00,rinv00),_mm256_sub_pd(rinvsq00,felec));
618
619             fscal            = felec;
620
621             fscal            = _mm256_andnot_pd(dummy_mask,fscal);
622
623             /* Calculate temporary vectorial force */
624             tx               = _mm256_mul_pd(fscal,dx00);
625             ty               = _mm256_mul_pd(fscal,dy00);
626             tz               = _mm256_mul_pd(fscal,dz00);
627
628             /* Update vectorial force */
629             fix0             = _mm256_add_pd(fix0,tx);
630             fiy0             = _mm256_add_pd(fiy0,ty);
631             fiz0             = _mm256_add_pd(fiz0,tz);
632
633             fjptrA             = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
634             fjptrB             = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
635             fjptrC             = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
636             fjptrD             = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
637             gmx_mm256_decrement_1rvec_4ptr_swizzle_pd(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
638
639             /* Inner loop uses 37 flops */
640         }
641
642         /* End of innermost loop */
643
644         gmx_mm256_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
645                                                  f+i_coord_offset,fshift+i_shift_offset);
646
647         /* Increment number of inner iterations */
648         inneriter                  += j_index_end - j_index_start;
649
650         /* Outer loop uses 7 flops */
651     }
652
653     /* Increment number of outer iterations */
654     outeriter        += nri;
655
656     /* Update outer/inner flops */
657
658     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_F,outeriter*7 + inneriter*37);
659 }