2 /* ## This file is part of the GROMACS molecular simulation package. */
4 /* ## Copyright (c) 2012, by the GROMACS development team, led by */
5 /* ## David van der Spoel, Berk Hess, Erik Lindahl, and including many */
6 /* ## others, as listed in the AUTHORS file in the top-level source */
7 /* ## directory and at http://www.gromacs.org. */
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. */
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. */
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. */
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. */
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. */
36 #error This file must be processed with the Gromacs pre-preprocessor
38 /* #if INCLUDE_HEADER */
45 #include "../nb_kernel.h"
46 #include "types/simple.h"
50 #include "gmx_math_x86_avx_256_single.h"
51 #include "kernelutil_x86_avx_256_single.h"
54 /* ## List of variables set by the generating script: */
56 /* ## Setttings that apply to the entire kernel: */
57 /* ## KERNEL_ELEC: String, choice for electrostatic interactions */
58 /* ## KERNEL_VDW: String, choice for van der Waals interactions */
59 /* ## KERNEL_NAME: String, name of this kernel */
60 /* ## KERNEL_VF: String telling if we calculate potential, force, or both */
61 /* ## GEOMETRY_I/GEOMETRY_J: String, name of each geometry, e.g. 'Water3' or '1Particle' */
63 /* ## Setttings that apply to particles in the outer (I) or inner (J) loops: */
64 /* ## PARTICLES_I[]/ Arrays with lists of i/j particles to use in kernel. It is */
65 /* ## PARTICLES_J[]: just [0] for particle geometry, but can be longer for water */
66 /* ## PARTICLES_ELEC_I[]/ Arrays with lists of i/j particle that have electrostatics */
67 /* ## PARTICLES_ELEC_J[]: interactions that should be calculated in this kernel. */
68 /* ## PARTICLES_VDW_I[]/ Arrays with the list of i/j particle that have VdW */
69 /* ## PARTICLES_VDW_J[]: interactions that should be calculated in this kernel. */
71 /* ## Setttings for pairs of interactions (e.g. 2nd i particle against 1st j particle) */
72 /* ## PAIRS_IJ[]: Array with (i,j) tuples of pairs for which interactions */
73 /* ## should be calculated in this kernel. Zero-charge particles */
74 /* ## do not have interactions with particles without vdw, and */
75 /* ## Vdw-only interactions are not evaluated in a no-vdw-kernel. */
76 /* ## INTERACTION_FLAGS[][]: 2D matrix, dimension e.g. 3*3 for water-water interactions. */
77 /* ## For each i-j pair, the element [I][J] is a list of strings */
78 /* ## defining properties/flags of this interaction. Examples */
79 /* ## include 'electrostatics'/'vdw' if that type of interaction */
80 /* ## should be evaluated, 'rsq'/'rinv'/'rinvsq' if those values */
81 /* ## are needed, and 'exactcutoff' or 'shift','switch' to */
82 /* ## decide if the force/potential should be modified. This way */
83 /* ## we only calculate values absolutely needed for each case. */
85 /* ## Calculate the size and offset for (merged/interleaved) table data */
88 * Gromacs nonbonded kernel: {KERNEL_NAME}
89 * Electrostatics interaction: {KERNEL_ELEC}
90 * VdW interaction: {KERNEL_VDW}
91 * Geometry: {GEOMETRY_I}-{GEOMETRY_J}
92 * Calculate force/pot: {KERNEL_VF}
96 (t_nblist * gmx_restrict nlist,
97 rvec * gmx_restrict xx,
98 rvec * gmx_restrict ff,
99 t_forcerec * gmx_restrict fr,
100 t_mdatoms * gmx_restrict mdatoms,
101 nb_kernel_data_t * gmx_restrict kernel_data,
102 t_nrnb * gmx_restrict nrnb)
104 /* ## Not all variables are used for all kernels, but any optimizing compiler fixes that, */
105 /* ## so there is no point in going to extremes to exclude variables that are not needed. */
106 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
107 * just 0 for non-waters.
108 * Suffixes A,B,C,D,E,F,G,H refer to j loop unrolling done with AVX, e.g. for the eight different
109 * jnr indices corresponding to data put in the four positions in the SIMD register.
111 int i_shift_offset,i_coord_offset,outeriter,inneriter;
112 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
113 int jnrA,jnrB,jnrC,jnrD;
114 int jnrE,jnrF,jnrG,jnrH;
115 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
116 int jnrlistE,jnrlistF,jnrlistG,jnrlistH;
117 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
118 int j_coord_offsetE,j_coord_offsetF,j_coord_offsetG,j_coord_offsetH;
119 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
121 real *shiftvec,*fshift,*x,*f;
122 real *fjptrA,*fjptrB,*fjptrC,*fjptrD,*fjptrE,*fjptrF,*fjptrG,*fjptrH;
124 __m256 tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
125 /* #for I in PARTICLES_I */
126 real * vdwioffsetptr{I};
127 __m256 ix{I},iy{I},iz{I},fix{I},fiy{I},fiz{I},iq{I},isai{I};
129 /* #for J in PARTICLES_J */
130 int vdwjidx{J}A,vdwjidx{J}B,vdwjidx{J}C,vdwjidx{J}D,vdwjidx{J}E,vdwjidx{J}F,vdwjidx{J}G,vdwjidx{J}H;
131 __m256 jx{J},jy{J},jz{J},fjx{J},fjy{J},fjz{J},jq{J},isaj{J};
133 /* #for I,J in PAIRS_IJ */
134 __m256 dx{I}{J},dy{I}{J},dz{I}{J},rsq{I}{J},rinv{I}{J},rinvsq{I}{J},r{I}{J},qq{I}{J},c6_{I}{J},c12_{I}{J};
136 /* #if KERNEL_ELEC != 'None' */
137 __m256 velec,felec,velecsum,facel,crf,krf,krf2;
140 /* #if 'GeneralizedBorn' in KERNEL_ELEC */
142 __m128i gbitab_lo,gbitab_hi;
143 __m256 vgb,fgb,vgbsum,dvdasum,gbscale,gbtabscale,isaprod,gbqqfactor,gbinvepsdiff,gbeps,dvdatmp;
144 __m256 minushalf = _mm256_set1_ps(-0.5);
145 real *invsqrta,*dvda,*gbtab;
147 /* #if KERNEL_VDW != 'None' */
149 __m256 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
152 __m256 one_sixth = _mm256_set1_ps(1.0/6.0);
153 __m256 one_twelfth = _mm256_set1_ps(1.0/12.0);
155 /* #if 'Table' in KERNEL_ELEC or 'GeneralizedBorn' in KERNEL_ELEC or 'Table' in KERNEL_VDW */
157 __m128i vfitab_lo,vfitab_hi;
158 __m128i ifour = _mm_set1_epi32(4);
159 __m256 rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF;
162 /* #if 'Ewald' in KERNEL_ELEC */
164 __m128i ewitab_lo,ewitab_hi;
165 __m256 ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
166 __m256 beta,beta2,beta3,zeta2,pmecorrF,pmecorrV,rinv3;
169 /* #if 'PotentialSwitch' in [KERNEL_MOD_ELEC,KERNEL_MOD_VDW] */
170 __m256 rswitch,swV3,swV4,swV5,swF2,swF3,swF4,d,d2,sw,dsw;
171 real rswitch_scalar,d_scalar;
173 __m256 dummy_mask,cutoff_mask;
174 __m256 signbit = _mm256_castsi256_ps( _mm256_set1_epi32(0x80000000) );
175 __m256 one = _mm256_set1_ps(1.0);
176 __m256 two = _mm256_set1_ps(2.0);
182 jindex = nlist->jindex;
184 shiftidx = nlist->shift;
186 shiftvec = fr->shift_vec[0];
187 fshift = fr->fshift[0];
188 /* #if KERNEL_ELEC != 'None' */
189 facel = _mm256_set1_ps(fr->epsfac);
190 charge = mdatoms->chargeA;
191 /* #if 'ReactionField' in KERNEL_ELEC */
192 krf = _mm256_set1_ps(fr->ic->k_rf);
193 krf2 = _mm256_set1_ps(fr->ic->k_rf*2.0);
194 crf = _mm256_set1_ps(fr->ic->c_rf);
197 /* #if KERNEL_VDW != 'None' */
198 nvdwtype = fr->ntype;
200 vdwtype = mdatoms->typeA;
203 /* #if 'Table' in KERNEL_ELEC and 'Table' in KERNEL_VDW */
204 vftab = kernel_data->table_elec_vdw->data;
205 vftabscale = _mm256_set1_ps(kernel_data->table_elec_vdw->scale);
206 /* #elif 'Table' in KERNEL_ELEC */
207 vftab = kernel_data->table_elec->data;
208 vftabscale = _mm256_set1_ps(kernel_data->table_elec->scale);
209 /* #elif 'Table' in KERNEL_VDW */
210 vftab = kernel_data->table_vdw->data;
211 vftabscale = _mm256_set1_ps(kernel_data->table_vdw->scale);
214 /* #if 'Ewald' in KERNEL_ELEC */
215 sh_ewald = _mm256_set1_ps(fr->ic->sh_ewald);
216 beta = _mm256_set1_ps(fr->ic->ewaldcoeff);
217 beta2 = _mm256_mul_ps(beta,beta);
218 beta3 = _mm256_mul_ps(beta,beta2);
220 /* #if KERNEL_VF=='Force' and KERNEL_MOD_ELEC!='PotentialSwitch' */
221 ewtab = fr->ic->tabq_coul_F;
222 ewtabscale = _mm256_set1_ps(fr->ic->tabq_scale);
223 ewtabhalfspace = _mm256_set1_ps(0.5/fr->ic->tabq_scale);
225 ewtab = fr->ic->tabq_coul_FDV0;
226 ewtabscale = _mm256_set1_ps(fr->ic->tabq_scale);
227 ewtabhalfspace = _mm256_set1_ps(0.5/fr->ic->tabq_scale);
231 /* #if KERNEL_ELEC=='GeneralizedBorn' */
232 invsqrta = fr->invsqrta;
234 gbtabscale = _mm256_set1_ps(fr->gbtab.scale);
235 gbtab = fr->gbtab.data;
236 gbinvepsdiff = _mm256_set1_ps((1.0/fr->epsilon_r) - (1.0/fr->gb_epsilon_solvent));
239 /* #if 'Water' in GEOMETRY_I */
240 /* Setup water-specific parameters */
241 inr = nlist->iinr[0];
242 /* #for I in PARTICLES_ELEC_I */
243 iq{I} = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+{I}]));
245 /* #for I in PARTICLES_VDW_I */
246 vdwioffsetptr{I} = vdwparam+2*nvdwtype*vdwtype[inr+{I}];
250 /* #if 'Water' in GEOMETRY_J */
251 /* #for J in PARTICLES_ELEC_J */
252 jq{J} = _mm256_set1_ps(charge[inr+{J}]);
254 /* #for J in PARTICLES_VDW_J */
255 vdwjidx{J}A = 2*vdwtype[inr+{J}];
257 /* #for I,J in PAIRS_IJ */
258 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] */
259 qq{I}{J} = _mm256_mul_ps(iq{I},jq{J});
261 /* #if 'vdw' in INTERACTION_FLAGS[I][J] */
262 c6_{I}{J} = _mm256_set1_ps(vdwioffsetptr{I}[vdwjidx{J}A]);
263 c12_{I}{J} = _mm256_set1_ps(vdwioffsetptr{I}[vdwjidx{J}A+1]);
268 /* #if KERNEL_MOD_ELEC!='None' or KERNEL_MOD_VDW!='None' */
269 /* #if KERNEL_ELEC!='None' */
270 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
271 rcutoff_scalar = fr->rcoulomb;
273 rcutoff_scalar = fr->rvdw;
275 rcutoff = _mm256_set1_ps(rcutoff_scalar);
276 rcutoff2 = _mm256_mul_ps(rcutoff,rcutoff);
279 /* #if KERNEL_MOD_VDW=='PotentialShift' */
280 sh_vdw_invrcut6 = _mm256_set1_ps(fr->ic->sh_invrc6);
281 rvdw = _mm256_set1_ps(fr->rvdw);
284 /* #if 'PotentialSwitch' in [KERNEL_MOD_ELEC,KERNEL_MOD_VDW] */
285 /* #if KERNEL_MOD_ELEC=='PotentialSwitch' */
286 rswitch_scalar = fr->rcoulomb_switch;
287 rswitch = _mm256_set1_ps(rswitch_scalar);
289 rswitch_scalar = fr->rvdw_switch;
290 rswitch = _mm256_set1_ps(rswitch_scalar);
292 /* Setup switch parameters */
293 d_scalar = rcutoff_scalar-rswitch_scalar;
294 d = _mm256_set1_ps(d_scalar);
295 swV3 = _mm256_set1_ps(-10.0/(d_scalar*d_scalar*d_scalar));
296 swV4 = _mm256_set1_ps( 15.0/(d_scalar*d_scalar*d_scalar*d_scalar));
297 swV5 = _mm256_set1_ps( -6.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
298 /* #if 'Force' in KERNEL_VF */
299 swF2 = _mm256_set1_ps(-30.0/(d_scalar*d_scalar*d_scalar));
300 swF3 = _mm256_set1_ps( 60.0/(d_scalar*d_scalar*d_scalar*d_scalar));
301 swF4 = _mm256_set1_ps(-30.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
305 /* Avoid stupid compiler warnings */
306 jnrA = jnrB = jnrC = jnrD = jnrE = jnrF = jnrG = jnrH = 0;
316 /* ## Keep track of the floating point operations we issue for reporting! */
317 /* #define OUTERFLOPS 0 */
321 for(iidx=0;iidx<4*DIM;iidx++)
326 /* Start outer loop over neighborlists */
327 for(iidx=0; iidx<nri; iidx++)
329 /* Load shift vector for this list */
330 i_shift_offset = DIM*shiftidx[iidx];
332 /* Load limits for loop over neighbors */
333 j_index_start = jindex[iidx];
334 j_index_end = jindex[iidx+1];
336 /* Get outer coordinate index */
338 i_coord_offset = DIM*inr;
340 /* Load i particle coords and add shift vector */
341 /* #if GEOMETRY_I == 'Particle' */
342 gmx_mm256_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
343 /* #elif GEOMETRY_I == 'Water3' */
344 gmx_mm256_load_shift_and_3rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,
345 &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2);
346 /* #elif GEOMETRY_I == 'Water4' */
347 /* #if 0 in PARTICLES_I */
348 gmx_mm256_load_shift_and_4rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,
349 &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
351 gmx_mm256_load_shift_and_3rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset+DIM,
352 &ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
356 /* #if 'Force' in KERNEL_VF */
357 /* #for I in PARTICLES_I */
358 fix{I} = _mm256_setzero_ps();
359 fiy{I} = _mm256_setzero_ps();
360 fiz{I} = _mm256_setzero_ps();
364 /* ## For water we already preloaded parameters at the start of the kernel */
365 /* #if not 'Water' in GEOMETRY_I */
366 /* Load parameters for i particles */
367 /* #for I in PARTICLES_ELEC_I */
368 iq{I} = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+{I}]));
369 /* #define OUTERFLOPS OUTERFLOPS+1 */
370 /* #if KERNEL_ELEC=='GeneralizedBorn' */
371 isai{I} = _mm256_set1_ps(invsqrta[inr+{I}]);
374 /* #for I in PARTICLES_VDW_I */
375 vdwioffsetptr{I} = vdwparam+2*nvdwtype*vdwtype[inr+{I}];
379 /* #if 'Potential' in KERNEL_VF */
380 /* Reset potential sums */
381 /* #if KERNEL_ELEC != 'None' */
382 velecsum = _mm256_setzero_ps();
384 /* #if 'GeneralizedBorn' in KERNEL_ELEC */
385 vgbsum = _mm256_setzero_ps();
387 /* #if KERNEL_VDW != 'None' */
388 vvdwsum = _mm256_setzero_ps();
391 /* #if 'GeneralizedBorn' in KERNEL_ELEC and 'Force' in KERNEL_VF */
392 dvdasum = _mm256_setzero_ps();
395 /* #for ROUND in ['Loop','Epilogue'] */
397 /* #if ROUND =='Loop' */
398 /* Start inner kernel loop */
399 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+7]>=0; jidx+=8)
401 /* ## First round is normal loop (next statement resets indentation) */
408 /* ## Second round is epilogue */
410 /* #define INNERFLOPS 0 */
412 /* Get j neighbor index, and coordinate index */
413 /* #if ROUND =='Loop' */
423 jnrlistA = jjnr[jidx];
424 jnrlistB = jjnr[jidx+1];
425 jnrlistC = jjnr[jidx+2];
426 jnrlistD = jjnr[jidx+3];
427 jnrlistE = jjnr[jidx+4];
428 jnrlistF = jjnr[jidx+5];
429 jnrlistG = jjnr[jidx+6];
430 jnrlistH = jjnr[jidx+7];
431 /* Sign of each element will be negative for non-real atoms.
432 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
433 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
435 dummy_mask = gmx_mm256_set_m128(gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx+4)),_mm_setzero_si128())),
436 gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128())));
438 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
439 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
440 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
441 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
442 jnrE = (jnrlistE>=0) ? jnrlistE : 0;
443 jnrF = (jnrlistF>=0) ? jnrlistF : 0;
444 jnrG = (jnrlistG>=0) ? jnrlistG : 0;
445 jnrH = (jnrlistH>=0) ? jnrlistH : 0;
447 j_coord_offsetA = DIM*jnrA;
448 j_coord_offsetB = DIM*jnrB;
449 j_coord_offsetC = DIM*jnrC;
450 j_coord_offsetD = DIM*jnrD;
451 j_coord_offsetE = DIM*jnrE;
452 j_coord_offsetF = DIM*jnrF;
453 j_coord_offsetG = DIM*jnrG;
454 j_coord_offsetH = DIM*jnrH;
456 /* load j atom coordinates */
457 /* #if GEOMETRY_J == 'Particle' */
458 gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
459 x+j_coord_offsetC,x+j_coord_offsetD,
460 x+j_coord_offsetE,x+j_coord_offsetF,
461 x+j_coord_offsetG,x+j_coord_offsetH,
463 /* #elif GEOMETRY_J == 'Water3' */
464 gmx_mm256_load_3rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
465 x+j_coord_offsetC,x+j_coord_offsetD,
466 x+j_coord_offsetE,x+j_coord_offsetF,
467 x+j_coord_offsetG,x+j_coord_offsetH,
468 &jx0,&jy0,&jz0,&jx1,&jy1,&jz1,&jx2,&jy2,&jz2);
469 /* #elif GEOMETRY_J == 'Water4' */
470 /* #if 0 in PARTICLES_J */
471 gmx_mm256_load_4rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
472 x+j_coord_offsetC,x+j_coord_offsetD,
473 x+j_coord_offsetE,x+j_coord_offsetF,
474 x+j_coord_offsetG,x+j_coord_offsetH,
475 &jx0,&jy0,&jz0,&jx1,&jy1,&jz1,&jx2,
476 &jy2,&jz2,&jx3,&jy3,&jz3);
478 gmx_mm256_load_3rvec_8ptr_swizzle_ps(x+j_coord_offsetA+DIM,x+j_coord_offsetB+DIM,
479 x+j_coord_offsetC+DIM,x+j_coord_offsetD+DIM,
480 x+j_coord_offsetE+DIM,x+j_coord_offsetF+DIM,
481 x+j_coord_offsetG+DIM,x+j_coord_offsetH+DIM,
482 &jx1,&jy1,&jz1,&jx2,&jy2,&jz2,&jx3,&jy3,&jz3);
486 /* Calculate displacement vector */
487 /* #for I,J in PAIRS_IJ */
488 dx{I}{J} = _mm256_sub_ps(ix{I},jx{J});
489 dy{I}{J} = _mm256_sub_ps(iy{I},jy{J});
490 dz{I}{J} = _mm256_sub_ps(iz{I},jz{J});
491 /* #define INNERFLOPS INNERFLOPS+3 */
494 /* Calculate squared distance and things based on it */
495 /* #for I,J in PAIRS_IJ */
496 rsq{I}{J} = gmx_mm256_calc_rsq_ps(dx{I}{J},dy{I}{J},dz{I}{J});
497 /* #define INNERFLOPS INNERFLOPS+5 */
500 /* #for I,J in PAIRS_IJ */
501 /* #if 'rinv' in INTERACTION_FLAGS[I][J] */
502 rinv{I}{J} = gmx_mm256_invsqrt_ps(rsq{I}{J});
503 /* #define INNERFLOPS INNERFLOPS+5 */
507 /* #for I,J in PAIRS_IJ */
508 /* #if 'rinvsq' in INTERACTION_FLAGS[I][J] */
509 /* # if 'rinv' not in INTERACTION_FLAGS[I][J] */
510 rinvsq{I}{J} = gmx_mm256_inv_ps(rsq{I}{J});
511 /* #define INNERFLOPS INNERFLOPS+4 */
513 rinvsq{I}{J} = _mm256_mul_ps(rinv{I}{J},rinv{I}{J});
514 /* #define INNERFLOPS INNERFLOPS+1 */
519 /* #if not 'Water' in GEOMETRY_J */
520 /* Load parameters for j particles */
521 /* #for J in PARTICLES_ELEC_J */
522 jq{J} = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+{J},charge+jnrB+{J},
523 charge+jnrC+{J},charge+jnrD+{J},
524 charge+jnrE+{J},charge+jnrF+{J},
525 charge+jnrG+{J},charge+jnrH+{J});
526 /* #if KERNEL_ELEC=='GeneralizedBorn' */
527 isaj{J} = gmx_mm256_load_8real_swizzle_ps(invsqrta+jnrA+{J},invsqrta+jnrB+{J},
528 invsqrta+jnrC+{J},invsqrta+jnrD+{J},
529 invsqrta+jnrE+{J},invsqrta+jnrF+{J},
530 invsqrta+jnrG+{J},invsqrta+jnrH+{J});
533 /* #for J in PARTICLES_VDW_J */
534 vdwjidx{J}A = 2*vdwtype[jnrA+{J}];
535 vdwjidx{J}B = 2*vdwtype[jnrB+{J}];
536 vdwjidx{J}C = 2*vdwtype[jnrC+{J}];
537 vdwjidx{J}D = 2*vdwtype[jnrD+{J}];
538 vdwjidx{J}E = 2*vdwtype[jnrE+{J}];
539 vdwjidx{J}F = 2*vdwtype[jnrF+{J}];
540 vdwjidx{J}G = 2*vdwtype[jnrG+{J}];
541 vdwjidx{J}H = 2*vdwtype[jnrH+{J}];
545 /* #if 'Force' in KERNEL_VF and not 'Particle' in GEOMETRY_I */
546 /* #for J in PARTICLES_J */
547 fjx{J} = _mm256_setzero_ps();
548 fjy{J} = _mm256_setzero_ps();
549 fjz{J} = _mm256_setzero_ps();
553 /* #for I,J in PAIRS_IJ */
555 /**************************
556 * CALCULATE INTERACTIONS *
557 **************************/
559 /* ## Note special check for TIP4P-TIP4P. Since we are cutting of all hydrogen interactions we also cut the LJ-only O-O interaction */
560 /* #if 'exactcutoff' in INTERACTION_FLAGS[I][J] or (GEOMETRY_I=='Water4' and GEOMETRY_J=='Water4' and 'exactcutoff' in INTERACTION_FLAGS[1][1]) */
561 /* ## We always calculate rinv/rinvsq above to enable pipelineing in compilers (performance tested on x86) */
562 if (gmx_mm256_any_lt(rsq{I}{J},rcutoff2))
564 /* #if 0 ## this and the next two lines is a hack to maintain auto-indentation in template file */
567 /* #define INNERFLOPS INNERFLOPS+1 */
570 /* #if 'r' in INTERACTION_FLAGS[I][J] */
571 r{I}{J} = _mm256_mul_ps(rsq{I}{J},rinv{I}{J});
572 /* #if ROUND == 'Epilogue' */
573 r{I}{J} = _mm256_andnot_ps(dummy_mask,r{I}{J});
574 /* #define INNERFLOPS INNERFLOPS+1 */
576 /* #define INNERFLOPS INNERFLOPS+1 */
579 /* ## For water geometries we already loaded parameters at the start of the kernel */
580 /* #if not 'Water' in GEOMETRY_J */
581 /* Compute parameters for interactions between i and j atoms */
582 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] */
583 qq{I}{J} = _mm256_mul_ps(iq{I},jq{J});
584 /* #define INNERFLOPS INNERFLOPS+1 */
586 /* #if 'vdw' in INTERACTION_FLAGS[I][J] */
587 gmx_mm256_load_8pair_swizzle_ps(vdwioffsetptr{I}+vdwjidx{J}A,
588 vdwioffsetptr{I}+vdwjidx{J}B,
589 vdwioffsetptr{I}+vdwjidx{J}C,
590 vdwioffsetptr{I}+vdwjidx{J}D,
591 vdwioffsetptr{I}+vdwjidx{J}E,
592 vdwioffsetptr{I}+vdwjidx{J}F,
593 vdwioffsetptr{I}+vdwjidx{J}G,
594 vdwioffsetptr{I}+vdwjidx{J}H,
595 &c6_{I}{J},&c12_{I}{J});
599 /* #if 'table' in INTERACTION_FLAGS[I][J] */
600 /* Calculate table index by multiplying r with table scale and truncate to integer */
601 rt = _mm256_mul_ps(r{I}{J},vftabscale);
602 vfitab = _mm256_cvttps_epi32(rt);
603 vfeps = _mm256_sub_ps(rt,_mm256_round_ps(rt, _MM_FROUND_FLOOR));
604 /* #define INNERFLOPS INNERFLOPS+4 */
605 /* AVX1 does not support 256-bit integer operations, so now we go to 128-bit mode... */
606 vfitab_lo = _mm256_extractf128_si256(vfitab,0x0);
607 vfitab_hi = _mm256_extractf128_si256(vfitab,0x1);
608 /* #if 'Table' in KERNEL_ELEC and 'Table' in KERNEL_VDW */
609 /* ## 3 tables, 4 bytes per point: multiply index by 12 */
610 vfitab_lo = _mm_slli_epi32(_mm_add_epi32(vfitab_lo,_mm_slli_epi32(vfitab_lo,1)),2);
611 vfitab_hi = _mm_slli_epi32(_mm_add_epi32(vfitab_hi,_mm_slli_epi32(vfitab_hi,1)),2);
612 /* #elif 'Table' in KERNEL_ELEC */
613 /* ## 1 table, 4 bytes per point: multiply index by 4 */
614 vfitab_lo = _mm_slli_epi32(vfitab_lo,2);
615 vfitab_hi = _mm_slli_epi32(vfitab_hi,2);
616 /* #elif 'Table' in KERNEL_VDW */
617 /* ## 2 tables, 4 bytes per point: multiply index by 8 */
618 vfitab_lo = _mm_slli_epi32(vfitab_lo,3);
619 vfitab_hi = _mm_slli_epi32(vfitab_hi,3);
623 /* ## ELECTROSTATIC INTERACTIONS */
624 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] */
626 /* #if KERNEL_ELEC=='Coulomb' */
628 /* COULOMB ELECTROSTATICS */
629 velec = _mm256_mul_ps(qq{I}{J},rinv{I}{J});
630 /* #define INNERFLOPS INNERFLOPS+1 */
631 /* #if 'Force' in KERNEL_VF */
632 felec = _mm256_mul_ps(velec,rinvsq{I}{J});
633 /* #define INNERFLOPS INNERFLOPS+1 */
636 /* #elif KERNEL_ELEC=='ReactionField' */
638 /* REACTION-FIELD ELECTROSTATICS */
639 /* #if 'Potential' in KERNEL_VF */
640 velec = _mm256_mul_ps(qq{I}{J},_mm256_sub_ps(_mm256_add_ps(rinv{I}{J},_mm256_mul_ps(krf,rsq{I}{J})),crf));
641 /* #define INNERFLOPS INNERFLOPS+4 */
643 /* #if 'Force' in KERNEL_VF */
644 felec = _mm256_mul_ps(qq{I}{J},_mm256_sub_ps(_mm256_mul_ps(rinv{I}{J},rinvsq{I}{J}),krf2));
645 /* #define INNERFLOPS INNERFLOPS+3 */
648 /* #elif KERNEL_ELEC=='GeneralizedBorn' */
650 /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
651 isaprod = _mm256_mul_ps(isai{I},isaj{J});
652 gbqqfactor = _mm256_xor_ps(signbit,_mm256_mul_ps(qq{I}{J},_mm256_mul_ps(isaprod,gbinvepsdiff)));
653 gbscale = _mm256_mul_ps(isaprod,gbtabscale);
654 /* #define INNERFLOPS INNERFLOPS+5 */
656 /* Calculate generalized born table index - this is a separate table from the normal one,
657 * but we use the same procedure by multiplying r with scale and truncating to integer.
659 rt = _mm256_mul_ps(r{I}{J},gbscale);
660 gbitab = _mm256_cvttps_epi32(rt);
661 gbeps = _mm256_sub_ps(rt,_mm256_round_ps(rt, _MM_FROUND_FLOOR));
662 /* AVX1 does not support 256-bit integer operations, so now we go to 128-bit mode... */
663 gbitab_lo = _mm256_extractf128_si256(gbitab,0x0);
664 gbitab_hi = _mm256_extractf128_si256(gbitab,0x1);
665 gbitab_lo = _mm_slli_epi32(gbitab_lo,2);
666 gbitab_hi = _mm_slli_epi32(gbitab_hi,2);
667 Y = gmx_mm256_set_m128(_mm_load_ps(gbtab + _mm_extract_epi32(gbitab_hi,0)),
668 _mm_load_ps(gbtab + _mm_extract_epi32(gbitab_lo,0)));
669 F = gmx_mm256_set_m128(_mm_load_ps(gbtab + _mm_extract_epi32(gbitab_hi,1)),
670 _mm_load_ps(gbtab + _mm_extract_epi32(gbitab_lo,1)));
671 G = gmx_mm256_set_m128(_mm_load_ps(gbtab + _mm_extract_epi32(gbitab_hi,2)),
672 _mm_load_ps(gbtab + _mm_extract_epi32(gbitab_lo,2)));
673 H = gmx_mm256_set_m128(_mm_load_ps(gbtab + _mm_extract_epi32(gbitab_hi,3)),
674 _mm_load_ps(gbtab + _mm_extract_epi32(gbitab_lo,3)));
675 GMX_MM256_HALFTRANSPOSE4_PS(Y,F,G,H);
676 Heps = _mm256_mul_ps(gbeps,H);
677 Fp = _mm256_add_ps(F,_mm256_mul_ps(gbeps,_mm256_add_ps(G,Heps)));
678 VV = _mm256_add_ps(Y,_mm256_mul_ps(gbeps,Fp));
679 vgb = _mm256_mul_ps(gbqqfactor,VV);
680 /* #define INNERFLOPS INNERFLOPS+10 */
682 /* #if 'Force' in KERNEL_VF */
683 FF = _mm256_add_ps(Fp,_mm256_mul_ps(gbeps,_mm256_add_ps(G,_mm256_add_ps(Heps,Heps))));
684 fgb = _mm256_mul_ps(gbqqfactor,_mm256_mul_ps(FF,gbscale));
685 dvdatmp = _mm256_mul_ps(minushalf,_mm256_add_ps(vgb,_mm256_mul_ps(fgb,r{I}{J})));
686 dvdasum = _mm256_add_ps(dvdasum,dvdatmp);
687 /* #if ROUND == 'Loop' */
697 /* The pointers to scratch make sure that this code with compilers that take gmx_restrict seriously (e.g. icc 13) really can't screw things up. */
698 fjptrA = (jnrlistA>=0) ? dvda+jnrA : scratch;
699 fjptrB = (jnrlistB>=0) ? dvda+jnrB : scratch;
700 fjptrC = (jnrlistC>=0) ? dvda+jnrC : scratch;
701 fjptrD = (jnrlistD>=0) ? dvda+jnrD : scratch;
702 fjptrE = (jnrlistE>=0) ? dvda+jnrE : scratch;
703 fjptrF = (jnrlistF>=0) ? dvda+jnrF : scratch;
704 fjptrG = (jnrlistG>=0) ? dvda+jnrG : scratch;
705 fjptrH = (jnrlistH>=0) ? dvda+jnrH : scratch;
707 gmx_mm256_increment_8real_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,
708 _mm256_mul_ps(dvdatmp,_mm256_mul_ps(isaj{J},isaj{J})));
709 /* #define INNERFLOPS INNERFLOPS+12 */
711 velec = _mm256_mul_ps(qq{I}{J},rinv{I}{J});
712 /* #define INNERFLOPS INNERFLOPS+1 */
713 /* #if 'Force' in KERNEL_VF */
714 felec = _mm256_mul_ps(_mm256_sub_ps(_mm256_mul_ps(velec,rinv{I}{J}),fgb),rinv{I}{J});
715 /* #define INNERFLOPS INNERFLOPS+3 */
718 /* #elif KERNEL_ELEC=='Ewald' */
719 /* EWALD ELECTROSTATICS */
721 /* Analytical PME correction */
722 zeta2 = _mm256_mul_ps(beta2,rsq{I}{J});
723 /* #if 'Force' in KERNEL_VF */
724 rinv3 = _mm256_mul_ps(rinvsq{I}{J},rinv{I}{J});
725 pmecorrF = gmx_mm256_pmecorrF_ps(zeta2);
726 felec = _mm256_add_ps( _mm256_mul_ps(pmecorrF,beta3), rinv3);
727 felec = _mm256_mul_ps(qq{I}{J},felec);
728 /* #define INNERFLOPS INNERFLOPS+31 */
730 /* #if 'Potential' in KERNEL_VF or KERNEL_MOD_ELEC=='PotentialSwitch' */
731 pmecorrV = gmx_mm256_pmecorrV_ps(zeta2);
732 pmecorrV = _mm256_mul_ps(pmecorrV,beta);
733 /* #define INNERFLOPS INNERFLOPS+27 */
734 /* #if KERNEL_MOD_ELEC=='PotentialShift' */
735 velec = _mm256_sub_ps(_mm256_sub_ps(rinv{I}{J},sh_ewald),pmecorrV);
736 /* #define INNERFLOPS INNERFLOPS+21 */
738 velec = _mm256_sub_ps(rinv{I}{J},pmecorrV);
740 velec = _mm256_mul_ps(qq{I}{J},velec);
743 /* #elif KERNEL_ELEC=='CubicSplineTable' */
745 /* CUBIC SPLINE TABLE ELECTROSTATICS */
746 Y = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,0)),
747 _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,0)));
748 F = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,1)),
749 _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,1)));
750 G = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,2)),
751 _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,2)));
752 H = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,3)),
753 _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,3)));
754 GMX_MM256_HALFTRANSPOSE4_PS(Y,F,G,H);
755 Heps = _mm256_mul_ps(vfeps,H);
756 Fp = _mm256_add_ps(F,_mm256_mul_ps(vfeps,_mm256_add_ps(G,Heps)));
757 /* #define INNERFLOPS INNERFLOPS+4 */
758 /* #if 'Potential' in KERNEL_VF */
759 VV = _mm256_add_ps(Y,_mm256_mul_ps(vfeps,Fp));
760 velec = _mm256_mul_ps(qq{I}{J},VV);
761 /* #define INNERFLOPS INNERFLOPS+3 */
763 /* #if 'Force' in KERNEL_VF */
764 FF = _mm256_add_ps(Fp,_mm256_mul_ps(vfeps,_mm256_add_ps(G,_mm256_add_ps(Heps,Heps))));
765 felec = _mm256_xor_ps(signbit,_mm256_mul_ps(_mm256_mul_ps(qq{I}{J},FF),_mm256_mul_ps(vftabscale,rinv{I}{J})));
766 /* #define INNERFLOPS INNERFLOPS+7 */
769 /* ## End of check for electrostatics interaction forms */
771 /* ## END OF ELECTROSTATIC INTERACTION CHECK FOR PAIR I-J */
773 /* #if 'vdw' in INTERACTION_FLAGS[I][J] */
775 /* #if KERNEL_VDW=='LennardJones' */
777 /* LENNARD-JONES DISPERSION/REPULSION */
779 rinvsix = _mm256_mul_ps(_mm256_mul_ps(rinvsq{I}{J},rinvsq{I}{J}),rinvsq{I}{J});
780 /* #define INNERFLOPS INNERFLOPS+2 */
781 /* #if 'Potential' in KERNEL_VF or KERNEL_MOD_VDW=='PotentialSwitch' */
782 vvdw6 = _mm256_mul_ps(c6_{I}{J},rinvsix);
783 vvdw12 = _mm256_mul_ps(c12_{I}{J},_mm256_mul_ps(rinvsix,rinvsix));
784 /* #define INNERFLOPS INNERFLOPS+3 */
785 /* #if KERNEL_MOD_VDW=='PotentialShift' */
786 vvdw = _mm256_sub_ps(_mm256_mul_ps( _mm256_sub_ps(vvdw12 , _mm256_mul_ps(c12_{I}{J},_mm256_mul_ps(sh_vdw_invrcut6,sh_vdw_invrcut6))), one_twelfth) ,
787 _mm256_mul_ps( _mm256_sub_ps(vvdw6,_mm256_mul_ps(c6_{I}{J},sh_vdw_invrcut6)),one_sixth));
788 /* #define INNERFLOPS INNERFLOPS+8 */
790 vvdw = _mm256_sub_ps( _mm256_mul_ps(vvdw12,one_twelfth) , _mm256_mul_ps(vvdw6,one_sixth) );
791 /* #define INNERFLOPS INNERFLOPS+3 */
793 /* ## Check for force inside potential check, i.e. this means we already did the potential part */
794 /* #if 'Force' in KERNEL_VF */
795 fvdw = _mm256_mul_ps(_mm256_sub_ps(vvdw12,vvdw6),rinvsq{I}{J});
796 /* #define INNERFLOPS INNERFLOPS+2 */
798 /* #elif KERNEL_VF=='Force' */
799 /* ## Force-only LennardJones makes it possible to save 1 flop (they do add up...) */
800 fvdw = _mm256_mul_ps(_mm256_sub_ps(_mm256_mul_ps(c12_{I}{J},rinvsix),c6_{I}{J}),_mm256_mul_ps(rinvsix,rinvsq{I}{J}));
801 /* #define INNERFLOPS INNERFLOPS+4 */
804 /* #elif KERNEL_VDW=='CubicSplineTable' */
806 /* CUBIC SPLINE TABLE DISPERSION */
807 /* #if 'Table' in KERNEL_ELEC */
808 vfitab_lo = _mm_add_epi32(vfitab_lo,ifour);
809 vfitab_hi = _mm_add_epi32(vfitab_hi,ifour);
811 Y = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,0)),
812 _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,0)));
813 F = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,1)),
814 _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,1)));
815 G = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,2)),
816 _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,2)));
817 H = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,3)),
818 _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,3)));
819 GMX_MM256_HALFTRANSPOSE4_PS(Y,F,G,H);
820 Heps = _mm256_mul_ps(vfeps,H);
821 Fp = _mm256_add_ps(F,_mm256_mul_ps(vfeps,_mm256_add_ps(G,Heps)));
822 /* #define INNERFLOPS INNERFLOPS+4 */
823 /* #if 'Potential' in KERNEL_VF */
824 VV = _mm256_add_ps(Y,_mm256_mul_ps(vfeps,Fp));
825 vvdw6 = _mm256_mul_ps(c6_{I}{J},VV);
826 /* #define INNERFLOPS INNERFLOPS+3 */
828 /* #if 'Force' in KERNEL_VF */
829 FF = _mm256_add_ps(Fp,_mm256_mul_ps(vfeps,_mm256_add_ps(G,_mm256_add_ps(Heps,Heps))));
830 fvdw6 = _mm256_mul_ps(c6_{I}{J},FF);
831 /* #define INNERFLOPS INNERFLOPS+4 */
834 /* CUBIC SPLINE TABLE REPULSION */
835 vfitab_lo = _mm_add_epi32(vfitab_lo,ifour);
836 vfitab_hi = _mm_add_epi32(vfitab_hi,ifour);
837 Y = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,0)),
838 _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,0)));
839 F = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,1)),
840 _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,1)));
841 G = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,2)),
842 _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,2)));
843 H = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,3)),
844 _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,3)));
845 GMX_MM256_HALFTRANSPOSE4_PS(Y,F,G,H);
846 Heps = _mm256_mul_ps(vfeps,H);
847 Fp = _mm256_add_ps(F,_mm256_mul_ps(vfeps,_mm256_add_ps(G,Heps)));
848 /* #define INNERFLOPS INNERFLOPS+4 */
849 /* #if 'Potential' in KERNEL_VF */
850 VV = _mm256_add_ps(Y,_mm256_mul_ps(vfeps,Fp));
851 vvdw12 = _mm256_mul_ps(c12_{I}{J},VV);
852 /* #define INNERFLOPS INNERFLOPS+3 */
854 /* #if 'Force' in KERNEL_VF */
855 FF = _mm256_add_ps(Fp,_mm256_mul_ps(vfeps,_mm256_add_ps(G,_mm256_add_ps(Heps,Heps))));
856 fvdw12 = _mm256_mul_ps(c12_{I}{J},FF);
857 /* #define INNERFLOPS INNERFLOPS+5 */
859 /* #if 'Potential' in KERNEL_VF */
860 vvdw = _mm256_add_ps(vvdw12,vvdw6);
861 /* #define INNERFLOPS INNERFLOPS+1 */
863 /* #if 'Force' in KERNEL_VF */
864 fvdw = _mm256_xor_ps(signbit,_mm256_mul_ps(_mm256_add_ps(fvdw6,fvdw12),_mm256_mul_ps(vftabscale,rinv{I}{J})));
865 /* #define INNERFLOPS INNERFLOPS+4 */
868 /* ## End of check for vdw interaction forms */
870 /* ## END OF VDW INTERACTION CHECK FOR PAIR I-J */
872 /* #if 'switch' in INTERACTION_FLAGS[I][J] */
873 d = _mm256_sub_ps(r{I}{J},rswitch);
874 d = _mm256_max_ps(d,_mm256_setzero_ps());
875 d2 = _mm256_mul_ps(d,d);
876 sw = _mm256_add_ps(one,_mm256_mul_ps(d2,_mm256_mul_ps(d,_mm256_add_ps(swV3,_mm256_mul_ps(d,_mm256_add_ps(swV4,_mm256_mul_ps(d,swV5)))))));
877 /* #define INNERFLOPS INNERFLOPS+10 */
879 /* #if 'Force' in KERNEL_VF */
880 dsw = _mm256_mul_ps(d2,_mm256_add_ps(swF2,_mm256_mul_ps(d,_mm256_add_ps(swF3,_mm256_mul_ps(d,swF4)))));
881 /* #define INNERFLOPS INNERFLOPS+5 */
884 /* Evaluate switch function */
885 /* #if 'Force' in KERNEL_VF */
886 /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
887 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] and KERNEL_MOD_ELEC=='PotentialSwitch' */
888 felec = _mm256_sub_ps( _mm256_mul_ps(felec,sw) , _mm256_mul_ps(rinv{I}{J},_mm256_mul_ps(velec,dsw)) );
889 /* #define INNERFLOPS INNERFLOPS+4 */
891 /* #if 'vdw' in INTERACTION_FLAGS[I][J] and KERNEL_MOD_VDW=='PotentialSwitch' */
892 fvdw = _mm256_sub_ps( _mm256_mul_ps(fvdw,sw) , _mm256_mul_ps(rinv{I}{J},_mm256_mul_ps(vvdw,dsw)) );
893 /* #define INNERFLOPS INNERFLOPS+4 */
896 /* #if 'Potential' in KERNEL_VF */
897 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] and KERNEL_MOD_ELEC=='PotentialSwitch' */
898 velec = _mm256_mul_ps(velec,sw);
899 /* #define INNERFLOPS INNERFLOPS+1 */
901 /* #if 'vdw' in INTERACTION_FLAGS[I][J] and KERNEL_MOD_VDW=='PotentialSwitch' */
902 vvdw = _mm256_mul_ps(vvdw,sw);
903 /* #define INNERFLOPS INNERFLOPS+1 */
907 /* ## Note special check for TIP4P-TIP4P. Since we are cutting of all hydrogen interactions we also cut the LJ-only O-O interaction */
908 /* #if 'exactcutoff' in INTERACTION_FLAGS[I][J] or (GEOMETRY_I=='Water4' and GEOMETRY_J=='Water4' and 'exactcutoff' in INTERACTION_FLAGS[1][1]) */
909 cutoff_mask = _mm256_cmp_ps(rsq{I}{J},rcutoff2,_CMP_LT_OQ);
910 /* #define INNERFLOPS INNERFLOPS+1 */
913 /* #if 'Potential' in KERNEL_VF */
914 /* Update potential sum for this i atom from the interaction with this j atom. */
915 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] */
916 /* #if 'exactcutoff' in INTERACTION_FLAGS[I][J] */
917 velec = _mm256_and_ps(velec,cutoff_mask);
918 /* #define INNERFLOPS INNERFLOPS+1 */
920 /* #if ROUND == 'Epilogue' */
921 velec = _mm256_andnot_ps(dummy_mask,velec);
923 velecsum = _mm256_add_ps(velecsum,velec);
924 /* #define INNERFLOPS INNERFLOPS+1 */
925 /* #if KERNEL_ELEC=='GeneralizedBorn' */
926 /* #if 'exactcutoff' in INTERACTION_FLAGS[I][J] */
927 vgb = _mm256_and_ps(vgb,cutoff_mask);
928 /* #define INNERFLOPS INNERFLOPS+1 */
930 /* #if ROUND == 'Epilogue' */
931 vgb = _mm256_andnot_ps(dummy_mask,vgb);
933 vgbsum = _mm256_add_ps(vgbsum,vgb);
934 /* #define INNERFLOPS INNERFLOPS+1 */
937 /* #if 'vdw' in INTERACTION_FLAGS[I][J] */
938 /* ## Note special check for TIP4P-TIP4P. Since we are cutting of all hydrogen interactions we also cut the LJ-only O-O interaction */
939 /* #if 'exactcutoff' in INTERACTION_FLAGS[I][J] or (GEOMETRY_I=='Water4' and GEOMETRY_J=='Water4' and 'exactcutoff' in INTERACTION_FLAGS[1][1]) */
940 vvdw = _mm256_and_ps(vvdw,cutoff_mask);
941 /* #define INNERFLOPS INNERFLOPS+1 */
943 /* #if ROUND == 'Epilogue' */
944 vvdw = _mm256_andnot_ps(dummy_mask,vvdw);
946 vvdwsum = _mm256_add_ps(vvdwsum,vvdw);
947 /* #define INNERFLOPS INNERFLOPS+1 */
951 /* #if 'Force' in KERNEL_VF */
953 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] and 'vdw' in INTERACTION_FLAGS[I][J] */
954 fscal = _mm256_add_ps(felec,fvdw);
955 /* #define INNERFLOPS INNERFLOPS+1 */
956 /* #elif 'electrostatics' in INTERACTION_FLAGS[I][J] */
958 /* #elif 'vdw' in INTERACTION_FLAGS[I][J] */
962 /* ## Note special check for TIP4P-TIP4P. Since we are cutting of all hydrogen interactions we also cut the LJ-only O-O interaction */
963 /* #if 'exactcutoff' in INTERACTION_FLAGS[I][J] or (GEOMETRY_I=='Water4' and GEOMETRY_J=='Water4' and 'exactcutoff' in INTERACTION_FLAGS[1][1]) */
964 fscal = _mm256_and_ps(fscal,cutoff_mask);
965 /* #define INNERFLOPS INNERFLOPS+1 */
968 /* #if ROUND == 'Epilogue' */
969 fscal = _mm256_andnot_ps(dummy_mask,fscal);
972 /* Calculate temporary vectorial force */
973 tx = _mm256_mul_ps(fscal,dx{I}{J});
974 ty = _mm256_mul_ps(fscal,dy{I}{J});
975 tz = _mm256_mul_ps(fscal,dz{I}{J});
977 /* Update vectorial force */
978 fix{I} = _mm256_add_ps(fix{I},tx);
979 fiy{I} = _mm256_add_ps(fiy{I},ty);
980 fiz{I} = _mm256_add_ps(fiz{I},tz);
981 /* #define INNERFLOPS INNERFLOPS+6 */
983 /* #if GEOMETRY_I == 'Particle' */
984 /* #if ROUND == 'Loop' */
985 fjptrA = f+j_coord_offsetA;
986 fjptrB = f+j_coord_offsetB;
987 fjptrC = f+j_coord_offsetC;
988 fjptrD = f+j_coord_offsetD;
989 fjptrE = f+j_coord_offsetE;
990 fjptrF = f+j_coord_offsetF;
991 fjptrG = f+j_coord_offsetG;
992 fjptrH = f+j_coord_offsetH;
994 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
995 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
996 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
997 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
998 fjptrE = (jnrlistE>=0) ? f+j_coord_offsetE : scratch;
999 fjptrF = (jnrlistF>=0) ? f+j_coord_offsetF : scratch;
1000 fjptrG = (jnrlistG>=0) ? f+j_coord_offsetG : scratch;
1001 fjptrH = (jnrlistH>=0) ? f+j_coord_offsetH : scratch;
1003 gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,tx,ty,tz);
1004 /* #define INNERFLOPS INNERFLOPS+3 */
1006 fjx{J} = _mm256_add_ps(fjx{J},tx);
1007 fjy{J} = _mm256_add_ps(fjy{J},ty);
1008 fjz{J} = _mm256_add_ps(fjz{J},tz);
1009 /* #define INNERFLOPS INNERFLOPS+3 */
1014 /* ## Note special check for TIP4P-TIP4P. Since we are cutting of all hydrogen interactions we also cut the LJ-only O-O interaction */
1015 /* #if 'exactcutoff' in INTERACTION_FLAGS[I][J] or (GEOMETRY_I=='Water4' and GEOMETRY_J=='Water4' and 'exactcutoff' in INTERACTION_FLAGS[1][1]) */
1016 /* #if 0 ## This and next two lines is a hack to maintain indentation in template file */
1021 /* ## End of check for the interaction being outside the cutoff */
1024 /* ## End of loop over i-j interaction pairs */
1026 /* #if GEOMETRY_I != 'Particle' */
1027 /* #if ROUND == 'Loop' */
1028 fjptrA = f+j_coord_offsetA;
1029 fjptrB = f+j_coord_offsetB;
1030 fjptrC = f+j_coord_offsetC;
1031 fjptrD = f+j_coord_offsetD;
1032 fjptrE = f+j_coord_offsetE;
1033 fjptrF = f+j_coord_offsetF;
1034 fjptrG = f+j_coord_offsetG;
1035 fjptrH = f+j_coord_offsetH;
1037 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
1038 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
1039 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
1040 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
1041 fjptrE = (jnrlistE>=0) ? f+j_coord_offsetE : scratch;
1042 fjptrF = (jnrlistF>=0) ? f+j_coord_offsetF : scratch;
1043 fjptrG = (jnrlistG>=0) ? f+j_coord_offsetG : scratch;
1044 fjptrH = (jnrlistH>=0) ? f+j_coord_offsetH : scratch;
1048 /* #if 'Water' in GEOMETRY_I and GEOMETRY_J == 'Particle' */
1049 gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,fjx0,fjy0,fjz0);
1050 /* #define INNERFLOPS INNERFLOPS+3 */
1051 /* #elif GEOMETRY_J == 'Water3' */
1052 gmx_mm256_decrement_3rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,
1053 fjx0,fjy0,fjz0,fjx1,fjy1,fjz1,fjx2,fjy2,fjz2);
1054 /* #define INNERFLOPS INNERFLOPS+9 */
1055 /* #elif GEOMETRY_J == 'Water4' */
1056 /* #if 0 in PARTICLES_J */
1057 gmx_mm256_decrement_4rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,
1058 fjx0,fjy0,fjz0,fjx1,fjy1,fjz1,
1059 fjx2,fjy2,fjz2,fjx3,fjy3,fjz3);
1060 /* #define INNERFLOPS INNERFLOPS+12 */
1062 gmx_mm256_decrement_3rvec_8ptr_swizzle_ps(fjptrA+DIM,fjptrB+DIM,fjptrC+DIM,fjptrD+DIM,
1063 fjptrE+DIM,fjptrF+DIM,fjptrG+DIM,fjptrH+DIM,
1064 fjx1,fjy1,fjz1,fjx2,fjy2,fjz2,fjx3,fjy3,fjz3);
1065 /* #define INNERFLOPS INNERFLOPS+9 */
1069 /* Inner loop uses {INNERFLOPS} flops */
1074 /* End of innermost loop */
1076 /* #if 'Force' in KERNEL_VF */
1077 /* #if GEOMETRY_I == 'Particle' */
1078 gmx_mm256_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
1079 f+i_coord_offset,fshift+i_shift_offset);
1080 /* #define OUTERFLOPS OUTERFLOPS+6 */
1081 /* #elif GEOMETRY_I == 'Water3' */
1082 gmx_mm256_update_iforce_3atom_swizzle_ps(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,
1083 f+i_coord_offset,fshift+i_shift_offset);
1084 /* #define OUTERFLOPS OUTERFLOPS+18 */
1085 /* #elif GEOMETRY_I == 'Water4' */
1086 /* #if 0 in PARTICLES_I */
1087 gmx_mm256_update_iforce_4atom_swizzle_ps(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
1088 f+i_coord_offset,fshift+i_shift_offset);
1089 /* #define OUTERFLOPS OUTERFLOPS+24 */
1091 gmx_mm256_update_iforce_3atom_swizzle_ps(fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
1092 f+i_coord_offset+DIM,fshift+i_shift_offset);
1093 /* #define OUTERFLOPS OUTERFLOPS+18 */
1098 /* #if 'Potential' in KERNEL_VF */
1100 /* Update potential energies */
1101 /* #if KERNEL_ELEC != 'None' */
1102 gmx_mm256_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
1103 /* #define OUTERFLOPS OUTERFLOPS+1 */
1105 /* #if 'GeneralizedBorn' in KERNEL_ELEC */
1106 gmx_mm256_update_1pot_ps(vgbsum,kernel_data->energygrp_polarization+ggid);
1107 /* #define OUTERFLOPS OUTERFLOPS+1 */
1109 /* #if KERNEL_VDW != 'None' */
1110 gmx_mm256_update_1pot_ps(vvdwsum,kernel_data->energygrp_vdw+ggid);
1111 /* #define OUTERFLOPS OUTERFLOPS+1 */
1114 /* #if 'GeneralizedBorn' in KERNEL_ELEC and 'Force' in KERNEL_VF */
1115 dvdasum = _mm256_mul_ps(dvdasum, _mm256_mul_ps(isai{I},isai{I}));
1116 gmx_mm256_update_1pot_ps(dvdasum,dvda+inr);
1119 /* Increment number of inner iterations */
1120 inneriter += j_index_end - j_index_start;
1122 /* Outer loop uses {OUTERFLOPS} flops */
1125 /* Increment number of outer iterations */
1128 /* Update outer/inner flops */
1129 /* ## NB: This is not important, it just affects the flopcount. However, since our preprocessor is */
1130 /* ## primitive and replaces aggressively even in strings inside these directives, we need to */
1131 /* ## assemble the main part of the name (containing KERNEL/ELEC/VDW) directly in the source. */
1132 /* #if GEOMETRY_I == 'Water3' */
1133 /* #define ISUFFIX '_W3' */
1134 /* #elif GEOMETRY_I == 'Water4' */
1135 /* #define ISUFFIX '_W4' */
1137 /* #define ISUFFIX '' */
1139 /* #if GEOMETRY_J == 'Water3' */
1140 /* #define JSUFFIX 'W3' */
1141 /* #elif GEOMETRY_J == 'Water4' */
1142 /* #define JSUFFIX 'W4' */
1144 /* #define JSUFFIX '' */
1146 /* #if 'PotentialAndForce' in KERNEL_VF */
1147 /* #define VFSUFFIX '_VF' */
1148 /* #elif 'Potential' in KERNEL_VF */
1149 /* #define VFSUFFIX '_V' */
1151 /* #define VFSUFFIX '_F' */
1154 /* #if KERNEL_ELEC != 'None' and KERNEL_VDW != 'None' */
1155 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW{ISUFFIX}{JSUFFIX}{VFSUFFIX},outeriter*{OUTERFLOPS} + inneriter*{INNERFLOPS});
1156 /* #elif KERNEL_ELEC != 'None' */
1157 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC{ISUFFIX}{JSUFFIX}{VFSUFFIX},outeriter*{OUTERFLOPS} + inneriter*{INNERFLOPS});
1159 inc_nrnb(nrnb,eNR_NBKERNEL_VDW{ISUFFIX}{JSUFFIX}{VFSUFFIX},outeriter*{OUTERFLOPS} + inneriter*{INNERFLOPS});