3 * This file is part of the GROMACS molecular simulation package.
5 * Copyright (c) 2012,2013,2014, by the GROMACS development team, led by
6 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7 * and including many others, as listed in the AUTHORS file in the
8 * top-level source directory and at http://www.gromacs.org.
10 * GROMACS is free software; you can redistribute it and/or
11 * modify it under the terms of the GNU Lesser General Public License
12 * as published by the Free Software Foundation; either version 2.1
13 * of the License, or (at your option) any later version.
15 * GROMACS is distributed in the hope that it will be useful,
16 * but WITHOUT ANY WARRANTY; without even the implied warranty of
17 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
18 * Lesser General Public License for more details.
20 * You should have received a copy of the GNU Lesser General Public
21 * License along with GROMACS; if not, see
22 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
25 * If you want to redistribute modifications to GROMACS, please
26 * consider that scientific software is very special. Version
27 * control is crucial - bugs must be traceable. We will be happy to
28 * consider code for inclusion in the official distribution, but
29 * derived work must not be called official GROMACS. Details are found
30 * in the README & COPYING files - if they are missing, get the
31 * official version at http://www.gromacs.org.
33 * To help us fund GROMACS development, we humbly ask that you cite
34 * 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 */
43 #include "../nb_kernel.h"
44 #include "gromacs/legacyheaders/types/simple.h"
45 #include "gromacs/math/vec.h"
46 #include "gromacs/legacyheaders/nrnb.h"
48 #include "gromacs/simd/math_x86_avx_256_single.h"
49 #include "kernelutil_x86_avx_256_single.h"
52 /* ## List of variables set by the generating script: */
54 /* ## Setttings that apply to the entire kernel: */
55 /* ## KERNEL_ELEC: String, choice for electrostatic interactions */
56 /* ## KERNEL_VDW: String, choice for van der Waals interactions */
57 /* ## KERNEL_NAME: String, name of this kernel */
58 /* ## KERNEL_VF: String telling if we calculate potential, force, or both */
59 /* ## GEOMETRY_I/GEOMETRY_J: String, name of each geometry, e.g. 'Water3' or '1Particle' */
61 /* ## Setttings that apply to particles in the outer (I) or inner (J) loops: */
62 /* ## PARTICLES_I[]/ Arrays with lists of i/j particles to use in kernel. It is */
63 /* ## PARTICLES_J[]: just [0] for particle geometry, but can be longer for water */
64 /* ## PARTICLES_ELEC_I[]/ Arrays with lists of i/j particle that have electrostatics */
65 /* ## PARTICLES_ELEC_J[]: interactions that should be calculated in this kernel. */
66 /* ## PARTICLES_VDW_I[]/ Arrays with the list of i/j particle that have VdW */
67 /* ## PARTICLES_VDW_J[]: interactions that should be calculated in this kernel. */
69 /* ## Setttings for pairs of interactions (e.g. 2nd i particle against 1st j particle) */
70 /* ## PAIRS_IJ[]: Array with (i,j) tuples of pairs for which interactions */
71 /* ## should be calculated in this kernel. Zero-charge particles */
72 /* ## do not have interactions with particles without vdw, and */
73 /* ## Vdw-only interactions are not evaluated in a no-vdw-kernel. */
74 /* ## INTERACTION_FLAGS[][]: 2D matrix, dimension e.g. 3*3 for water-water interactions. */
75 /* ## For each i-j pair, the element [I][J] is a list of strings */
76 /* ## defining properties/flags of this interaction. Examples */
77 /* ## include 'electrostatics'/'vdw' if that type of interaction */
78 /* ## should be evaluated, 'rsq'/'rinv'/'rinvsq' if those values */
79 /* ## are needed, and 'exactcutoff' or 'shift','switch' to */
80 /* ## decide if the force/potential should be modified. This way */
81 /* ## we only calculate values absolutely needed for each case. */
83 /* ## Calculate the size and offset for (merged/interleaved) table data */
89 * Gromacs nonbonded kernel: {KERNEL_NAME}
90 * Electrostatics interaction: {KERNEL_ELEC}
91 * VdW interaction: {KERNEL_VDW}
92 * Geometry: {GEOMETRY_I}-{GEOMETRY_J}
93 * Calculate force/pot: {KERNEL_VF}
97 (t_nblist * gmx_restrict nlist,
98 rvec * gmx_restrict xx,
99 rvec * gmx_restrict ff,
100 t_forcerec * gmx_restrict fr,
101 t_mdatoms * gmx_restrict mdatoms,
102 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
103 t_nrnb * gmx_restrict nrnb)
105 /* ## Not all variables are used for all kernels, but any optimizing compiler fixes that, */
106 /* ## so there is no point in going to extremes to exclude variables that are not needed. */
107 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
108 * just 0 for non-waters.
109 * Suffixes A,B,C,D,E,F,G,H refer to j loop unrolling done with AVX, e.g. for the eight different
110 * jnr indices corresponding to data put in the four positions in the SIMD register.
112 int i_shift_offset,i_coord_offset,outeriter,inneriter;
113 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
114 int jnrA,jnrB,jnrC,jnrD;
115 int jnrE,jnrF,jnrG,jnrH;
116 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
117 int jnrlistE,jnrlistF,jnrlistG,jnrlistH;
118 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
119 int j_coord_offsetE,j_coord_offsetF,j_coord_offsetG,j_coord_offsetH;
120 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
122 real *shiftvec,*fshift,*x,*f;
123 real *fjptrA,*fjptrB,*fjptrC,*fjptrD,*fjptrE,*fjptrF,*fjptrG,*fjptrH;
125 __m256 tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
126 /* #for I in PARTICLES_I */
127 real * vdwioffsetptr{I};
128 /* #if 'LJEwald' in KERNEL_VDW */
129 real * vdwgridioffsetptr{I};
131 __m256 ix{I},iy{I},iz{I},fix{I},fiy{I},fiz{I},iq{I},isai{I};
133 /* #for J in PARTICLES_J */
134 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;
135 __m256 jx{J},jy{J},jz{J},fjx{J},fjy{J},fjz{J},jq{J},isaj{J};
137 /* #for I,J in PAIRS_IJ */
138 __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};
140 /* #if KERNEL_ELEC != 'None' */
141 __m256 velec,felec,velecsum,facel,crf,krf,krf2;
144 /* #if 'GeneralizedBorn' in KERNEL_ELEC */
146 __m128i gbitab_lo,gbitab_hi;
147 __m256 vgb,fgb,vgbsum,dvdasum,gbscale,gbtabscale,isaprod,gbqqfactor,gbinvepsdiff,gbeps,dvdatmp;
148 __m256 minushalf = _mm256_set1_ps(-0.5);
149 real *invsqrta,*dvda,*gbtab;
151 /* #if KERNEL_VDW != 'None' */
153 __m256 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
156 __m256 one_sixth = _mm256_set1_ps(1.0/6.0);
157 __m256 one_twelfth = _mm256_set1_ps(1.0/12.0);
159 /* #if 'Table' in KERNEL_ELEC or 'GeneralizedBorn' in KERNEL_ELEC or 'Table' in KERNEL_VDW */
161 __m128i vfitab_lo,vfitab_hi;
162 __m128i ifour = _mm_set1_epi32(4);
163 __m256 rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF;
166 /* #if 'LJEwald' in KERNEL_VDW */
167 /* #for I,J in PAIRS_IJ */
168 __m256 c6grid_{I}{J};
171 __m256 ewclj,ewclj2,ewclj6,ewcljrsq,poly,exponent,f6A,f6B,sh_lj_ewald;
172 __m256 one_half = _mm256_set1_ps(0.5);
173 __m256 minus_one = _mm256_set1_ps(-1.0);
175 /* #if 'Ewald' in KERNEL_ELEC */
177 __m128i ewitab_lo,ewitab_hi;
178 __m256 ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
179 __m256 beta,beta2,beta3,zeta2,pmecorrF,pmecorrV,rinv3;
182 /* #if 'PotentialSwitch' in [KERNEL_MOD_ELEC,KERNEL_MOD_VDW] */
183 __m256 rswitch,swV3,swV4,swV5,swF2,swF3,swF4,d,d2,sw,dsw;
184 real rswitch_scalar,d_scalar;
186 __m256 dummy_mask,cutoff_mask;
187 __m256 signbit = _mm256_castsi256_ps( _mm256_set1_epi32(0x80000000) );
188 __m256 one = _mm256_set1_ps(1.0);
189 __m256 two = _mm256_set1_ps(2.0);
195 jindex = nlist->jindex;
197 shiftidx = nlist->shift;
199 shiftvec = fr->shift_vec[0];
200 fshift = fr->fshift[0];
201 /* #if KERNEL_ELEC != 'None' */
202 facel = _mm256_set1_ps(fr->epsfac);
203 charge = mdatoms->chargeA;
204 /* #if 'ReactionField' in KERNEL_ELEC */
205 krf = _mm256_set1_ps(fr->ic->k_rf);
206 krf2 = _mm256_set1_ps(fr->ic->k_rf*2.0);
207 crf = _mm256_set1_ps(fr->ic->c_rf);
210 /* #if KERNEL_VDW != 'None' */
211 nvdwtype = fr->ntype;
213 vdwtype = mdatoms->typeA;
215 /* #if 'LJEwald' in KERNEL_VDW */
216 vdwgridparam = fr->ljpme_c6grid;
217 sh_lj_ewald = _mm256_set1_ps(fr->ic->sh_lj_ewald);
218 ewclj = _mm256_set1_ps(fr->ewaldcoeff_lj);
219 ewclj2 = _mm256_mul_ps(minus_one,_mm256_mul_ps(ewclj,ewclj));
222 /* #if 'Table' in KERNEL_ELEC and 'Table' in KERNEL_VDW */
223 vftab = kernel_data->table_elec_vdw->data;
224 vftabscale = _mm256_set1_ps(kernel_data->table_elec_vdw->scale);
225 /* #elif 'Table' in KERNEL_ELEC */
226 vftab = kernel_data->table_elec->data;
227 vftabscale = _mm256_set1_ps(kernel_data->table_elec->scale);
228 /* #elif 'Table' in KERNEL_VDW */
229 vftab = kernel_data->table_vdw->data;
230 vftabscale = _mm256_set1_ps(kernel_data->table_vdw->scale);
233 /* #if 'Ewald' in KERNEL_ELEC */
234 sh_ewald = _mm256_set1_ps(fr->ic->sh_ewald);
235 beta = _mm256_set1_ps(fr->ic->ewaldcoeff_q);
236 beta2 = _mm256_mul_ps(beta,beta);
237 beta3 = _mm256_mul_ps(beta,beta2);
239 /* #if KERNEL_VF=='Force' and KERNEL_MOD_ELEC!='PotentialSwitch' */
240 ewtab = fr->ic->tabq_coul_F;
241 ewtabscale = _mm256_set1_ps(fr->ic->tabq_scale);
242 ewtabhalfspace = _mm256_set1_ps(0.5/fr->ic->tabq_scale);
244 ewtab = fr->ic->tabq_coul_FDV0;
245 ewtabscale = _mm256_set1_ps(fr->ic->tabq_scale);
246 ewtabhalfspace = _mm256_set1_ps(0.5/fr->ic->tabq_scale);
250 /* #if KERNEL_ELEC=='GeneralizedBorn' */
251 invsqrta = fr->invsqrta;
253 gbtabscale = _mm256_set1_ps(fr->gbtab.scale);
254 gbtab = fr->gbtab.data;
255 gbinvepsdiff = _mm256_set1_ps((1.0/fr->epsilon_r) - (1.0/fr->gb_epsilon_solvent));
258 /* #if 'Water' in GEOMETRY_I */
259 /* Setup water-specific parameters */
260 inr = nlist->iinr[0];
261 /* #for I in PARTICLES_ELEC_I */
262 iq{I} = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+{I}]));
264 /* #for I in PARTICLES_VDW_I */
265 vdwioffsetptr{I} = vdwparam+2*nvdwtype*vdwtype[inr+{I}];
266 /* #if 'LJEwald' in KERNEL_VDW */
267 vdwgridioffsetptr{I} = vdwgridparam+2*nvdwtype*vdwtype[inr+{I}];
272 /* #if 'Water' in GEOMETRY_J */
273 /* #for J in PARTICLES_ELEC_J */
274 jq{J} = _mm256_set1_ps(charge[inr+{J}]);
276 /* #for J in PARTICLES_VDW_J */
277 vdwjidx{J}A = 2*vdwtype[inr+{J}];
279 /* #for I,J in PAIRS_IJ */
280 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] */
281 qq{I}{J} = _mm256_mul_ps(iq{I},jq{J});
283 /* #if 'vdw' in INTERACTION_FLAGS[I][J] */
284 /* #if 'LJEwald' in KERNEL_VDW */
285 c6_{I}{J} = _mm256_set1_ps(vdwioffsetptr{I}[vdwjidx{J}A]);
286 c12_{I}{J} = _mm256_set1_ps(vdwioffsetptr{I}[vdwjidx{J}A+1]);
287 c6grid_{I}{J} = _mm256_set1_ps(vdwgridioffsetptr{I}[vdwjidx{J}A]);
289 c6_{I}{J} = _mm256_set1_ps(vdwioffsetptr{I}[vdwjidx{J}A]);
290 c12_{I}{J} = _mm256_set1_ps(vdwioffsetptr{I}[vdwjidx{J}A+1]);
296 /* #if KERNEL_MOD_ELEC!='None' or KERNEL_MOD_VDW!='None' */
297 /* #if KERNEL_ELEC!='None' */
298 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
299 rcutoff_scalar = fr->rcoulomb;
301 rcutoff_scalar = fr->rvdw;
303 rcutoff = _mm256_set1_ps(rcutoff_scalar);
304 rcutoff2 = _mm256_mul_ps(rcutoff,rcutoff);
307 /* #if KERNEL_MOD_VDW=='PotentialShift' */
308 sh_vdw_invrcut6 = _mm256_set1_ps(fr->ic->sh_invrc6);
309 rvdw = _mm256_set1_ps(fr->rvdw);
312 /* #if 'PotentialSwitch' in [KERNEL_MOD_ELEC,KERNEL_MOD_VDW] */
313 /* #if KERNEL_MOD_ELEC=='PotentialSwitch' */
314 rswitch_scalar = fr->rcoulomb_switch;
315 rswitch = _mm256_set1_ps(rswitch_scalar);
317 rswitch_scalar = fr->rvdw_switch;
318 rswitch = _mm256_set1_ps(rswitch_scalar);
320 /* Setup switch parameters */
321 d_scalar = rcutoff_scalar-rswitch_scalar;
322 d = _mm256_set1_ps(d_scalar);
323 swV3 = _mm256_set1_ps(-10.0/(d_scalar*d_scalar*d_scalar));
324 swV4 = _mm256_set1_ps( 15.0/(d_scalar*d_scalar*d_scalar*d_scalar));
325 swV5 = _mm256_set1_ps( -6.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
326 /* #if 'Force' in KERNEL_VF */
327 swF2 = _mm256_set1_ps(-30.0/(d_scalar*d_scalar*d_scalar));
328 swF3 = _mm256_set1_ps( 60.0/(d_scalar*d_scalar*d_scalar*d_scalar));
329 swF4 = _mm256_set1_ps(-30.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
333 /* Avoid stupid compiler warnings */
334 jnrA = jnrB = jnrC = jnrD = jnrE = jnrF = jnrG = jnrH = 0;
344 /* ## Keep track of the floating point operations we issue for reporting! */
345 /* #define OUTERFLOPS 0 */
349 for(iidx=0;iidx<4*DIM;iidx++)
354 /* Start outer loop over neighborlists */
355 for(iidx=0; iidx<nri; iidx++)
357 /* Load shift vector for this list */
358 i_shift_offset = DIM*shiftidx[iidx];
360 /* Load limits for loop over neighbors */
361 j_index_start = jindex[iidx];
362 j_index_end = jindex[iidx+1];
364 /* Get outer coordinate index */
366 i_coord_offset = DIM*inr;
368 /* Load i particle coords and add shift vector */
369 /* #if GEOMETRY_I == 'Particle' */
370 gmx_mm256_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
371 /* #elif GEOMETRY_I == 'Water3' */
372 gmx_mm256_load_shift_and_3rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,
373 &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2);
374 /* #elif GEOMETRY_I == 'Water4' */
375 /* #if 0 in PARTICLES_I */
376 gmx_mm256_load_shift_and_4rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,
377 &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
379 gmx_mm256_load_shift_and_3rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset+DIM,
380 &ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
384 /* #if 'Force' in KERNEL_VF */
385 /* #for I in PARTICLES_I */
386 fix{I} = _mm256_setzero_ps();
387 fiy{I} = _mm256_setzero_ps();
388 fiz{I} = _mm256_setzero_ps();
392 /* ## For water we already preloaded parameters at the start of the kernel */
393 /* #if not 'Water' in GEOMETRY_I */
394 /* Load parameters for i particles */
395 /* #for I in PARTICLES_ELEC_I */
396 iq{I} = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+{I}]));
397 /* #define OUTERFLOPS OUTERFLOPS+1 */
398 /* #if KERNEL_ELEC=='GeneralizedBorn' */
399 isai{I} = _mm256_set1_ps(invsqrta[inr+{I}]);
402 /* #for I in PARTICLES_VDW_I */
403 vdwioffsetptr{I} = vdwparam+2*nvdwtype*vdwtype[inr+{I}];
404 /* #if 'LJEwald' in KERNEL_VDW */
405 vdwgridioffsetptr{I} = vdwgridparam+2*nvdwtype*vdwtype[inr+{I}];
410 /* #if 'Potential' in KERNEL_VF */
411 /* Reset potential sums */
412 /* #if KERNEL_ELEC != 'None' */
413 velecsum = _mm256_setzero_ps();
415 /* #if 'GeneralizedBorn' in KERNEL_ELEC */
416 vgbsum = _mm256_setzero_ps();
418 /* #if KERNEL_VDW != 'None' */
419 vvdwsum = _mm256_setzero_ps();
422 /* #if 'GeneralizedBorn' in KERNEL_ELEC and 'Force' in KERNEL_VF */
423 dvdasum = _mm256_setzero_ps();
426 /* #for ROUND in ['Loop','Epilogue'] */
428 /* #if ROUND =='Loop' */
429 /* Start inner kernel loop */
430 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+7]>=0; jidx+=8)
432 /* ## First round is normal loop (next statement resets indentation) */
439 /* ## Second round is epilogue */
441 /* #define INNERFLOPS 0 */
443 /* Get j neighbor index, and coordinate index */
444 /* #if ROUND =='Loop' */
454 jnrlistA = jjnr[jidx];
455 jnrlistB = jjnr[jidx+1];
456 jnrlistC = jjnr[jidx+2];
457 jnrlistD = jjnr[jidx+3];
458 jnrlistE = jjnr[jidx+4];
459 jnrlistF = jjnr[jidx+5];
460 jnrlistG = jjnr[jidx+6];
461 jnrlistH = jjnr[jidx+7];
462 /* Sign of each element will be negative for non-real atoms.
463 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
464 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
466 dummy_mask = gmx_mm256_set_m128(gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx+4)),_mm_setzero_si128())),
467 gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128())));
469 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
470 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
471 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
472 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
473 jnrE = (jnrlistE>=0) ? jnrlistE : 0;
474 jnrF = (jnrlistF>=0) ? jnrlistF : 0;
475 jnrG = (jnrlistG>=0) ? jnrlistG : 0;
476 jnrH = (jnrlistH>=0) ? jnrlistH : 0;
478 j_coord_offsetA = DIM*jnrA;
479 j_coord_offsetB = DIM*jnrB;
480 j_coord_offsetC = DIM*jnrC;
481 j_coord_offsetD = DIM*jnrD;
482 j_coord_offsetE = DIM*jnrE;
483 j_coord_offsetF = DIM*jnrF;
484 j_coord_offsetG = DIM*jnrG;
485 j_coord_offsetH = DIM*jnrH;
487 /* load j atom coordinates */
488 /* #if GEOMETRY_J == 'Particle' */
489 gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
490 x+j_coord_offsetC,x+j_coord_offsetD,
491 x+j_coord_offsetE,x+j_coord_offsetF,
492 x+j_coord_offsetG,x+j_coord_offsetH,
494 /* #elif GEOMETRY_J == 'Water3' */
495 gmx_mm256_load_3rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
496 x+j_coord_offsetC,x+j_coord_offsetD,
497 x+j_coord_offsetE,x+j_coord_offsetF,
498 x+j_coord_offsetG,x+j_coord_offsetH,
499 &jx0,&jy0,&jz0,&jx1,&jy1,&jz1,&jx2,&jy2,&jz2);
500 /* #elif GEOMETRY_J == 'Water4' */
501 /* #if 0 in PARTICLES_J */
502 gmx_mm256_load_4rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
503 x+j_coord_offsetC,x+j_coord_offsetD,
504 x+j_coord_offsetE,x+j_coord_offsetF,
505 x+j_coord_offsetG,x+j_coord_offsetH,
506 &jx0,&jy0,&jz0,&jx1,&jy1,&jz1,&jx2,
507 &jy2,&jz2,&jx3,&jy3,&jz3);
509 gmx_mm256_load_3rvec_8ptr_swizzle_ps(x+j_coord_offsetA+DIM,x+j_coord_offsetB+DIM,
510 x+j_coord_offsetC+DIM,x+j_coord_offsetD+DIM,
511 x+j_coord_offsetE+DIM,x+j_coord_offsetF+DIM,
512 x+j_coord_offsetG+DIM,x+j_coord_offsetH+DIM,
513 &jx1,&jy1,&jz1,&jx2,&jy2,&jz2,&jx3,&jy3,&jz3);
517 /* Calculate displacement vector */
518 /* #for I,J in PAIRS_IJ */
519 dx{I}{J} = _mm256_sub_ps(ix{I},jx{J});
520 dy{I}{J} = _mm256_sub_ps(iy{I},jy{J});
521 dz{I}{J} = _mm256_sub_ps(iz{I},jz{J});
522 /* #define INNERFLOPS INNERFLOPS+3 */
525 /* Calculate squared distance and things based on it */
526 /* #for I,J in PAIRS_IJ */
527 rsq{I}{J} = gmx_mm256_calc_rsq_ps(dx{I}{J},dy{I}{J},dz{I}{J});
528 /* #define INNERFLOPS INNERFLOPS+5 */
531 /* #for I,J in PAIRS_IJ */
532 /* #if 'rinv' in INTERACTION_FLAGS[I][J] */
533 rinv{I}{J} = gmx_mm256_invsqrt_ps(rsq{I}{J});
534 /* #define INNERFLOPS INNERFLOPS+5 */
538 /* #for I,J in PAIRS_IJ */
539 /* #if 'rinvsq' in INTERACTION_FLAGS[I][J] */
540 /* # if 'rinv' not in INTERACTION_FLAGS[I][J] */
541 rinvsq{I}{J} = gmx_mm256_inv_ps(rsq{I}{J});
542 /* #define INNERFLOPS INNERFLOPS+4 */
544 rinvsq{I}{J} = _mm256_mul_ps(rinv{I}{J},rinv{I}{J});
545 /* #define INNERFLOPS INNERFLOPS+1 */
550 /* #if not 'Water' in GEOMETRY_J */
551 /* Load parameters for j particles */
552 /* #for J in PARTICLES_ELEC_J */
553 jq{J} = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+{J},charge+jnrB+{J},
554 charge+jnrC+{J},charge+jnrD+{J},
555 charge+jnrE+{J},charge+jnrF+{J},
556 charge+jnrG+{J},charge+jnrH+{J});
557 /* #if KERNEL_ELEC=='GeneralizedBorn' */
558 isaj{J} = gmx_mm256_load_8real_swizzle_ps(invsqrta+jnrA+{J},invsqrta+jnrB+{J},
559 invsqrta+jnrC+{J},invsqrta+jnrD+{J},
560 invsqrta+jnrE+{J},invsqrta+jnrF+{J},
561 invsqrta+jnrG+{J},invsqrta+jnrH+{J});
564 /* #for J in PARTICLES_VDW_J */
565 vdwjidx{J}A = 2*vdwtype[jnrA+{J}];
566 vdwjidx{J}B = 2*vdwtype[jnrB+{J}];
567 vdwjidx{J}C = 2*vdwtype[jnrC+{J}];
568 vdwjidx{J}D = 2*vdwtype[jnrD+{J}];
569 vdwjidx{J}E = 2*vdwtype[jnrE+{J}];
570 vdwjidx{J}F = 2*vdwtype[jnrF+{J}];
571 vdwjidx{J}G = 2*vdwtype[jnrG+{J}];
572 vdwjidx{J}H = 2*vdwtype[jnrH+{J}];
576 /* #if 'Force' in KERNEL_VF and not 'Particle' in GEOMETRY_I */
577 /* #for J in PARTICLES_J */
578 fjx{J} = _mm256_setzero_ps();
579 fjy{J} = _mm256_setzero_ps();
580 fjz{J} = _mm256_setzero_ps();
584 /* #for I,J in PAIRS_IJ */
586 /**************************
587 * CALCULATE INTERACTIONS *
588 **************************/
590 /* ## Note special check for TIP4P-TIP4P. Since we are cutting of all hydrogen interactions we also cut the LJ-only O-O interaction */
591 /* #if 'exactcutoff' in INTERACTION_FLAGS[I][J] or (GEOMETRY_I=='Water4' and GEOMETRY_J=='Water4' and 'exactcutoff' in INTERACTION_FLAGS[1][1]) */
592 /* ## We always calculate rinv/rinvsq above to enable pipelineing in compilers (performance tested on x86) */
593 if (gmx_mm256_any_lt(rsq{I}{J},rcutoff2))
595 /* #if 0 ## this and the next two lines is a hack to maintain auto-indentation in template file */
598 /* #define INNERFLOPS INNERFLOPS+1 */
601 /* #if 'r' in INTERACTION_FLAGS[I][J] */
602 r{I}{J} = _mm256_mul_ps(rsq{I}{J},rinv{I}{J});
603 /* #if ROUND == 'Epilogue' */
604 r{I}{J} = _mm256_andnot_ps(dummy_mask,r{I}{J});
605 /* #define INNERFLOPS INNERFLOPS+1 */
607 /* #define INNERFLOPS INNERFLOPS+1 */
610 /* ## For water geometries we already loaded parameters at the start of the kernel */
611 /* #if not 'Water' in GEOMETRY_J */
612 /* Compute parameters for interactions between i and j atoms */
613 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] */
614 qq{I}{J} = _mm256_mul_ps(iq{I},jq{J});
615 /* #define INNERFLOPS INNERFLOPS+1 */
617 /* #if 'vdw' in INTERACTION_FLAGS[I][J] */
618 gmx_mm256_load_8pair_swizzle_ps(vdwioffsetptr{I}+vdwjidx{J}A,
619 vdwioffsetptr{I}+vdwjidx{J}B,
620 vdwioffsetptr{I}+vdwjidx{J}C,
621 vdwioffsetptr{I}+vdwjidx{J}D,
622 vdwioffsetptr{I}+vdwjidx{J}E,
623 vdwioffsetptr{I}+vdwjidx{J}F,
624 vdwioffsetptr{I}+vdwjidx{J}G,
625 vdwioffsetptr{I}+vdwjidx{J}H,
626 &c6_{I}{J},&c12_{I}{J});
628 /* #if 'LJEwald' in KERNEL_VDW */
629 c6grid_{I}{J} = gmx_mm256_load_8real_swizzle_ps(vdwgridioffsetptr{I}+vdwjidx{J}A,
630 vdwgridioffsetptr{I}+vdwjidx{J}B,
631 vdwgridioffsetptr{I}+vdwjidx{J}C,
632 vdwgridioffsetptr{I}+vdwjidx{J}D,
633 vdwgridioffsetptr{I}+vdwjidx{J}E,
634 vdwgridioffsetptr{I}+vdwjidx{J}F,
635 vdwgridioffsetptr{I}+vdwjidx{J}G,
636 vdwgridioffsetptr{I}+vdwjidx{J}H);
641 /* #if 'table' in INTERACTION_FLAGS[I][J] */
642 /* Calculate table index by multiplying r with table scale and truncate to integer */
643 rt = _mm256_mul_ps(r{I}{J},vftabscale);
644 vfitab = _mm256_cvttps_epi32(rt);
645 vfeps = _mm256_sub_ps(rt,_mm256_round_ps(rt, _MM_FROUND_FLOOR));
646 /* #define INNERFLOPS INNERFLOPS+4 */
647 /* AVX1 does not support 256-bit integer operations, so now we go to 128-bit mode... */
648 vfitab_lo = _mm256_extractf128_si256(vfitab,0x0);
649 vfitab_hi = _mm256_extractf128_si256(vfitab,0x1);
650 /* #if 'Table' in KERNEL_ELEC and 'Table' in KERNEL_VDW */
651 /* ## 3 tables, 4 bytes per point: multiply index by 12 */
652 vfitab_lo = _mm_slli_epi32(_mm_add_epi32(vfitab_lo,_mm_slli_epi32(vfitab_lo,1)),2);
653 vfitab_hi = _mm_slli_epi32(_mm_add_epi32(vfitab_hi,_mm_slli_epi32(vfitab_hi,1)),2);
654 /* #elif 'Table' in KERNEL_ELEC */
655 /* ## 1 table, 4 bytes per point: multiply index by 4 */
656 vfitab_lo = _mm_slli_epi32(vfitab_lo,2);
657 vfitab_hi = _mm_slli_epi32(vfitab_hi,2);
658 /* #elif 'Table' in KERNEL_VDW */
659 /* ## 2 tables, 4 bytes per point: multiply index by 8 */
660 vfitab_lo = _mm_slli_epi32(vfitab_lo,3);
661 vfitab_hi = _mm_slli_epi32(vfitab_hi,3);
665 /* ## ELECTROSTATIC INTERACTIONS */
666 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] */
668 /* #if KERNEL_ELEC=='Coulomb' */
670 /* COULOMB ELECTROSTATICS */
671 velec = _mm256_mul_ps(qq{I}{J},rinv{I}{J});
672 /* #define INNERFLOPS INNERFLOPS+1 */
673 /* #if 'Force' in KERNEL_VF */
674 felec = _mm256_mul_ps(velec,rinvsq{I}{J});
675 /* #define INNERFLOPS INNERFLOPS+1 */
678 /* #elif KERNEL_ELEC=='ReactionField' */
680 /* REACTION-FIELD ELECTROSTATICS */
681 /* #if 'Potential' in KERNEL_VF */
682 velec = _mm256_mul_ps(qq{I}{J},_mm256_sub_ps(_mm256_add_ps(rinv{I}{J},_mm256_mul_ps(krf,rsq{I}{J})),crf));
683 /* #define INNERFLOPS INNERFLOPS+4 */
685 /* #if 'Force' in KERNEL_VF */
686 felec = _mm256_mul_ps(qq{I}{J},_mm256_sub_ps(_mm256_mul_ps(rinv{I}{J},rinvsq{I}{J}),krf2));
687 /* #define INNERFLOPS INNERFLOPS+3 */
690 /* #elif KERNEL_ELEC=='GeneralizedBorn' */
692 /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
693 isaprod = _mm256_mul_ps(isai{I},isaj{J});
694 gbqqfactor = _mm256_xor_ps(signbit,_mm256_mul_ps(qq{I}{J},_mm256_mul_ps(isaprod,gbinvepsdiff)));
695 gbscale = _mm256_mul_ps(isaprod,gbtabscale);
696 /* #define INNERFLOPS INNERFLOPS+5 */
698 /* Calculate generalized born table index - this is a separate table from the normal one,
699 * but we use the same procedure by multiplying r with scale and truncating to integer.
701 rt = _mm256_mul_ps(r{I}{J},gbscale);
702 gbitab = _mm256_cvttps_epi32(rt);
703 gbeps = _mm256_sub_ps(rt,_mm256_round_ps(rt, _MM_FROUND_FLOOR));
704 /* AVX1 does not support 256-bit integer operations, so now we go to 128-bit mode... */
705 gbitab_lo = _mm256_extractf128_si256(gbitab,0x0);
706 gbitab_hi = _mm256_extractf128_si256(gbitab,0x1);
707 gbitab_lo = _mm_slli_epi32(gbitab_lo,2);
708 gbitab_hi = _mm_slli_epi32(gbitab_hi,2);
709 Y = gmx_mm256_set_m128(_mm_load_ps(gbtab + _mm_extract_epi32(gbitab_hi,0)),
710 _mm_load_ps(gbtab + _mm_extract_epi32(gbitab_lo,0)));
711 F = gmx_mm256_set_m128(_mm_load_ps(gbtab + _mm_extract_epi32(gbitab_hi,1)),
712 _mm_load_ps(gbtab + _mm_extract_epi32(gbitab_lo,1)));
713 G = gmx_mm256_set_m128(_mm_load_ps(gbtab + _mm_extract_epi32(gbitab_hi,2)),
714 _mm_load_ps(gbtab + _mm_extract_epi32(gbitab_lo,2)));
715 H = gmx_mm256_set_m128(_mm_load_ps(gbtab + _mm_extract_epi32(gbitab_hi,3)),
716 _mm_load_ps(gbtab + _mm_extract_epi32(gbitab_lo,3)));
717 GMX_MM256_HALFTRANSPOSE4_PS(Y,F,G,H);
718 Heps = _mm256_mul_ps(gbeps,H);
719 Fp = _mm256_add_ps(F,_mm256_mul_ps(gbeps,_mm256_add_ps(G,Heps)));
720 VV = _mm256_add_ps(Y,_mm256_mul_ps(gbeps,Fp));
721 vgb = _mm256_mul_ps(gbqqfactor,VV);
722 /* #define INNERFLOPS INNERFLOPS+10 */
724 /* #if 'Force' in KERNEL_VF */
725 FF = _mm256_add_ps(Fp,_mm256_mul_ps(gbeps,_mm256_add_ps(G,_mm256_add_ps(Heps,Heps))));
726 fgb = _mm256_mul_ps(gbqqfactor,_mm256_mul_ps(FF,gbscale));
727 dvdatmp = _mm256_mul_ps(minushalf,_mm256_add_ps(vgb,_mm256_mul_ps(fgb,r{I}{J})));
728 /* #if ROUND == 'Epilogue' */
729 dvdatmp = _mm256_andnot_ps(dummy_mask,dvdatmp);
731 dvdasum = _mm256_add_ps(dvdasum,dvdatmp);
732 /* #if ROUND == 'Loop' */
742 /* 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. */
743 fjptrA = (jnrlistA>=0) ? dvda+jnrA : scratch;
744 fjptrB = (jnrlistB>=0) ? dvda+jnrB : scratch;
745 fjptrC = (jnrlistC>=0) ? dvda+jnrC : scratch;
746 fjptrD = (jnrlistD>=0) ? dvda+jnrD : scratch;
747 fjptrE = (jnrlistE>=0) ? dvda+jnrE : scratch;
748 fjptrF = (jnrlistF>=0) ? dvda+jnrF : scratch;
749 fjptrG = (jnrlistG>=0) ? dvda+jnrG : scratch;
750 fjptrH = (jnrlistH>=0) ? dvda+jnrH : scratch;
752 gmx_mm256_increment_8real_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,
753 _mm256_mul_ps(dvdatmp,_mm256_mul_ps(isaj{J},isaj{J})));
754 /* #define INNERFLOPS INNERFLOPS+12 */
756 velec = _mm256_mul_ps(qq{I}{J},rinv{I}{J});
757 /* #define INNERFLOPS INNERFLOPS+1 */
758 /* #if 'Force' in KERNEL_VF */
759 felec = _mm256_mul_ps(_mm256_sub_ps(_mm256_mul_ps(velec,rinv{I}{J}),fgb),rinv{I}{J});
760 /* #define INNERFLOPS INNERFLOPS+3 */
763 /* #elif KERNEL_ELEC=='Ewald' */
764 /* EWALD ELECTROSTATICS */
766 /* Analytical PME correction */
767 zeta2 = _mm256_mul_ps(beta2,rsq{I}{J});
768 /* #if 'Force' in KERNEL_VF */
769 rinv3 = _mm256_mul_ps(rinvsq{I}{J},rinv{I}{J});
770 pmecorrF = gmx_mm256_pmecorrF_ps(zeta2);
771 felec = _mm256_add_ps( _mm256_mul_ps(pmecorrF,beta3), rinv3);
772 felec = _mm256_mul_ps(qq{I}{J},felec);
773 /* #define INNERFLOPS INNERFLOPS+31 */
775 /* #if 'Potential' in KERNEL_VF or KERNEL_MOD_ELEC=='PotentialSwitch' */
776 pmecorrV = gmx_mm256_pmecorrV_ps(zeta2);
777 pmecorrV = _mm256_mul_ps(pmecorrV,beta);
778 /* #define INNERFLOPS INNERFLOPS+27 */
779 /* #if KERNEL_MOD_ELEC=='PotentialShift' */
780 velec = _mm256_sub_ps(_mm256_sub_ps(rinv{I}{J},sh_ewald),pmecorrV);
781 /* #define INNERFLOPS INNERFLOPS+21 */
783 velec = _mm256_sub_ps(rinv{I}{J},pmecorrV);
785 velec = _mm256_mul_ps(qq{I}{J},velec);
788 /* #elif KERNEL_ELEC=='CubicSplineTable' */
790 /* CUBIC SPLINE TABLE ELECTROSTATICS */
791 Y = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,0)),
792 _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,0)));
793 F = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,1)),
794 _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,1)));
795 G = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,2)),
796 _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,2)));
797 H = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,3)),
798 _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,3)));
799 GMX_MM256_HALFTRANSPOSE4_PS(Y,F,G,H);
800 Heps = _mm256_mul_ps(vfeps,H);
801 Fp = _mm256_add_ps(F,_mm256_mul_ps(vfeps,_mm256_add_ps(G,Heps)));
802 /* #define INNERFLOPS INNERFLOPS+4 */
803 /* #if 'Potential' in KERNEL_VF */
804 VV = _mm256_add_ps(Y,_mm256_mul_ps(vfeps,Fp));
805 velec = _mm256_mul_ps(qq{I}{J},VV);
806 /* #define INNERFLOPS INNERFLOPS+3 */
808 /* #if 'Force' in KERNEL_VF */
809 FF = _mm256_add_ps(Fp,_mm256_mul_ps(vfeps,_mm256_add_ps(G,_mm256_add_ps(Heps,Heps))));
810 felec = _mm256_xor_ps(signbit,_mm256_mul_ps(_mm256_mul_ps(qq{I}{J},FF),_mm256_mul_ps(vftabscale,rinv{I}{J})));
811 /* #define INNERFLOPS INNERFLOPS+7 */
814 /* ## End of check for electrostatics interaction forms */
816 /* ## END OF ELECTROSTATIC INTERACTION CHECK FOR PAIR I-J */
818 /* #if 'vdw' in INTERACTION_FLAGS[I][J] */
820 /* #if KERNEL_VDW=='LennardJones' */
822 /* LENNARD-JONES DISPERSION/REPULSION */
824 rinvsix = _mm256_mul_ps(_mm256_mul_ps(rinvsq{I}{J},rinvsq{I}{J}),rinvsq{I}{J});
825 /* #define INNERFLOPS INNERFLOPS+2 */
826 /* #if 'Potential' in KERNEL_VF or KERNEL_MOD_VDW=='PotentialSwitch' */
827 vvdw6 = _mm256_mul_ps(c6_{I}{J},rinvsix);
828 vvdw12 = _mm256_mul_ps(c12_{I}{J},_mm256_mul_ps(rinvsix,rinvsix));
829 /* #define INNERFLOPS INNERFLOPS+3 */
830 /* #if KERNEL_MOD_VDW=='PotentialShift' */
831 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) ,
832 _mm256_mul_ps( _mm256_sub_ps(vvdw6,_mm256_mul_ps(c6_{I}{J},sh_vdw_invrcut6)),one_sixth));
833 /* #define INNERFLOPS INNERFLOPS+8 */
835 vvdw = _mm256_sub_ps( _mm256_mul_ps(vvdw12,one_twelfth) , _mm256_mul_ps(vvdw6,one_sixth) );
836 /* #define INNERFLOPS INNERFLOPS+3 */
838 /* ## Check for force inside potential check, i.e. this means we already did the potential part */
839 /* #if 'Force' in KERNEL_VF */
840 fvdw = _mm256_mul_ps(_mm256_sub_ps(vvdw12,vvdw6),rinvsq{I}{J});
841 /* #define INNERFLOPS INNERFLOPS+2 */
843 /* #elif KERNEL_VF=='Force' */
844 /* ## Force-only LennardJones makes it possible to save 1 flop (they do add up...) */
845 fvdw = _mm256_mul_ps(_mm256_sub_ps(_mm256_mul_ps(c12_{I}{J},rinvsix),c6_{I}{J}),_mm256_mul_ps(rinvsix,rinvsq{I}{J}));
846 /* #define INNERFLOPS INNERFLOPS+4 */
849 /* #elif KERNEL_VDW=='CubicSplineTable' */
851 /* CUBIC SPLINE TABLE DISPERSION */
852 /* #if 'Table' in KERNEL_ELEC */
853 vfitab_lo = _mm_add_epi32(vfitab_lo,ifour);
854 vfitab_hi = _mm_add_epi32(vfitab_hi,ifour);
856 Y = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,0)),
857 _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,0)));
858 F = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,1)),
859 _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,1)));
860 G = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,2)),
861 _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,2)));
862 H = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,3)),
863 _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,3)));
864 GMX_MM256_HALFTRANSPOSE4_PS(Y,F,G,H);
865 Heps = _mm256_mul_ps(vfeps,H);
866 Fp = _mm256_add_ps(F,_mm256_mul_ps(vfeps,_mm256_add_ps(G,Heps)));
867 /* #define INNERFLOPS INNERFLOPS+4 */
868 /* #if 'Potential' in KERNEL_VF */
869 VV = _mm256_add_ps(Y,_mm256_mul_ps(vfeps,Fp));
870 vvdw6 = _mm256_mul_ps(c6_{I}{J},VV);
871 /* #define INNERFLOPS INNERFLOPS+3 */
873 /* #if 'Force' in KERNEL_VF */
874 FF = _mm256_add_ps(Fp,_mm256_mul_ps(vfeps,_mm256_add_ps(G,_mm256_add_ps(Heps,Heps))));
875 fvdw6 = _mm256_mul_ps(c6_{I}{J},FF);
876 /* #define INNERFLOPS INNERFLOPS+4 */
879 /* CUBIC SPLINE TABLE REPULSION */
880 vfitab_lo = _mm_add_epi32(vfitab_lo,ifour);
881 vfitab_hi = _mm_add_epi32(vfitab_hi,ifour);
882 Y = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,0)),
883 _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,0)));
884 F = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,1)),
885 _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,1)));
886 G = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,2)),
887 _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,2)));
888 H = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,3)),
889 _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,3)));
890 GMX_MM256_HALFTRANSPOSE4_PS(Y,F,G,H);
891 Heps = _mm256_mul_ps(vfeps,H);
892 Fp = _mm256_add_ps(F,_mm256_mul_ps(vfeps,_mm256_add_ps(G,Heps)));
893 /* #define INNERFLOPS INNERFLOPS+4 */
894 /* #if 'Potential' in KERNEL_VF */
895 VV = _mm256_add_ps(Y,_mm256_mul_ps(vfeps,Fp));
896 vvdw12 = _mm256_mul_ps(c12_{I}{J},VV);
897 /* #define INNERFLOPS INNERFLOPS+3 */
899 /* #if 'Force' in KERNEL_VF */
900 FF = _mm256_add_ps(Fp,_mm256_mul_ps(vfeps,_mm256_add_ps(G,_mm256_add_ps(Heps,Heps))));
901 fvdw12 = _mm256_mul_ps(c12_{I}{J},FF);
902 /* #define INNERFLOPS INNERFLOPS+5 */
904 /* #if 'Potential' in KERNEL_VF */
905 vvdw = _mm256_add_ps(vvdw12,vvdw6);
906 /* #define INNERFLOPS INNERFLOPS+1 */
908 /* #if 'Force' in KERNEL_VF */
909 fvdw = _mm256_xor_ps(signbit,_mm256_mul_ps(_mm256_add_ps(fvdw6,fvdw12),_mm256_mul_ps(vftabscale,rinv{I}{J})));
910 /* #define INNERFLOPS INNERFLOPS+4 */
913 /* #elif KERNEL_VDW=='LJEwald' */
915 /* Analytical LJ-PME */
916 rinvsix = _mm256_mul_ps(_mm256_mul_ps(rinvsq{I}{J},rinvsq{I}{J}),rinvsq{I}{J});
917 ewcljrsq = _mm256_mul_ps(ewclj2,rsq{I}{J});
918 ewclj6 = _mm256_mul_ps(ewclj2,_mm256_mul_ps(ewclj2,ewclj2));
919 exponent = gmx_simd_exp_r(ewcljrsq);
920 /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
921 poly = _mm256_mul_ps(exponent,_mm256_add_ps(_mm256_sub_ps(one,ewcljrsq),_mm256_mul_ps(_mm256_mul_ps(ewcljrsq,ewcljrsq),one_half)));
922 /* #define INNERFLOPS INNERFLOPS+11 */
923 /* #if 'Potential' in KERNEL_VF or KERNEL_MOD_VDW=='PotentialSwitch' */
924 /* vvdw6 = [C6 - C6grid * (1-poly)]/r6 */
925 vvdw6 = _mm256_mul_ps(_mm256_sub_ps(c6_{I}{J},_mm256_mul_ps(c6grid_{I}{J},_mm256_sub_ps(one,poly))),rinvsix);
926 vvdw12 = _mm256_mul_ps(c12_{I}{J},_mm256_mul_ps(rinvsix,rinvsix));
927 /* #define INNERFLOPS INNERFLOPS+6 */
928 /* #if KERNEL_MOD_VDW=='PotentialShift' */
929 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) ,
930 _mm256_mul_ps( _mm256_sub_ps(vvdw6,_mm256_add_ps(_mm256_mul_ps(c6_{I}{J},sh_vdw_invrcut6),_mm256_mul_ps(c6grid_{I}{J},sh_lj_ewald))),one_sixth));
931 /* #define INNERFLOPS INNERFLOPS+10 */
933 vvdw = _mm256_sub_ps(_mm256_mul_ps(vvdw12,one_twelfth),_mm256_mul_ps(vvdw6,one_sixth));
934 /* #define INNERFLOPS INNERFLOPS+3 */
936 /* ## Check for force inside potential check, i.e. this means we already did the potential part */
937 /* #if 'Force' in KERNEL_VF */
938 /* fvdw = vvdw12/r - (vvdw6/r + (C6grid * exponent * beta^6)/r) */
939 fvdw = _mm256_mul_ps(_mm256_sub_ps(vvdw12,_mm256_sub_ps(vvdw6,_mm256_mul_ps(_mm256_mul_ps(c6grid_{I}{J},one_sixth),_mm256_mul_ps(exponent,ewclj6)))),rinvsq{I}{J});
940 /* #define INNERFLOPS INNERFLOPS+6 */
942 /* #elif KERNEL_VF=='Force' */
943 /* f6A = 6 * C6grid * (1 - poly) */
944 f6A = _mm256_mul_ps(c6grid_{I}{J},_mm256_sub_ps(one,poly));
945 /* f6B = C6grid * exponent * beta^6 */
946 f6B = _mm256_mul_ps(_mm256_mul_ps(c6grid_{I}{J},one_sixth),_mm256_mul_ps(exponent,ewclj6));
947 /* fvdw = 12*C12/r13 - ((6*C6 - f6A)/r6 + f6B)/r */
948 fvdw = _mm256_mul_ps(_mm256_add_ps(_mm256_mul_ps(_mm256_sub_ps(_mm256_mul_ps(c12_{I}{J},rinvsix),_mm256_sub_ps(c6_{I}{J},f6A)),rinvsix),f6B),rinvsq{I}{J});
949 /* #define INNERFLOPS INNERFLOPS+11 */
952 /* ## End of check for vdw interaction forms */
954 /* ## END OF VDW INTERACTION CHECK FOR PAIR I-J */
956 /* #if 'switch' in INTERACTION_FLAGS[I][J] */
957 d = _mm256_sub_ps(r{I}{J},rswitch);
958 d = _mm256_max_ps(d,_mm256_setzero_ps());
959 d2 = _mm256_mul_ps(d,d);
960 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)))))));
961 /* #define INNERFLOPS INNERFLOPS+10 */
963 /* #if 'Force' in KERNEL_VF */
964 dsw = _mm256_mul_ps(d2,_mm256_add_ps(swF2,_mm256_mul_ps(d,_mm256_add_ps(swF3,_mm256_mul_ps(d,swF4)))));
965 /* #define INNERFLOPS INNERFLOPS+5 */
968 /* Evaluate switch function */
969 /* #if 'Force' in KERNEL_VF */
970 /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
971 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] and KERNEL_MOD_ELEC=='PotentialSwitch' */
972 felec = _mm256_sub_ps( _mm256_mul_ps(felec,sw) , _mm256_mul_ps(rinv{I}{J},_mm256_mul_ps(velec,dsw)) );
973 /* #define INNERFLOPS INNERFLOPS+4 */
975 /* #if 'vdw' in INTERACTION_FLAGS[I][J] and KERNEL_MOD_VDW=='PotentialSwitch' */
976 fvdw = _mm256_sub_ps( _mm256_mul_ps(fvdw,sw) , _mm256_mul_ps(rinv{I}{J},_mm256_mul_ps(vvdw,dsw)) );
977 /* #define INNERFLOPS INNERFLOPS+4 */
980 /* #if 'Potential' in KERNEL_VF */
981 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] and KERNEL_MOD_ELEC=='PotentialSwitch' */
982 velec = _mm256_mul_ps(velec,sw);
983 /* #define INNERFLOPS INNERFLOPS+1 */
985 /* #if 'vdw' in INTERACTION_FLAGS[I][J] and KERNEL_MOD_VDW=='PotentialSwitch' */
986 vvdw = _mm256_mul_ps(vvdw,sw);
987 /* #define INNERFLOPS INNERFLOPS+1 */
991 /* ## Note special check for TIP4P-TIP4P. Since we are cutting of all hydrogen interactions we also cut the LJ-only O-O interaction */
992 /* #if 'exactcutoff' in INTERACTION_FLAGS[I][J] or (GEOMETRY_I=='Water4' and GEOMETRY_J=='Water4' and 'exactcutoff' in INTERACTION_FLAGS[1][1]) */
993 cutoff_mask = _mm256_cmp_ps(rsq{I}{J},rcutoff2,_CMP_LT_OQ);
994 /* #define INNERFLOPS INNERFLOPS+1 */
997 /* #if 'Potential' in KERNEL_VF */
998 /* Update potential sum for this i atom from the interaction with this j atom. */
999 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] */
1000 /* #if 'exactcutoff' in INTERACTION_FLAGS[I][J] */
1001 velec = _mm256_and_ps(velec,cutoff_mask);
1002 /* #define INNERFLOPS INNERFLOPS+1 */
1004 /* #if ROUND == 'Epilogue' */
1005 velec = _mm256_andnot_ps(dummy_mask,velec);
1007 velecsum = _mm256_add_ps(velecsum,velec);
1008 /* #define INNERFLOPS INNERFLOPS+1 */
1009 /* #if KERNEL_ELEC=='GeneralizedBorn' */
1010 /* #if 'exactcutoff' in INTERACTION_FLAGS[I][J] */
1011 vgb = _mm256_and_ps(vgb,cutoff_mask);
1012 /* #define INNERFLOPS INNERFLOPS+1 */
1014 /* #if ROUND == 'Epilogue' */
1015 vgb = _mm256_andnot_ps(dummy_mask,vgb);
1017 vgbsum = _mm256_add_ps(vgbsum,vgb);
1018 /* #define INNERFLOPS INNERFLOPS+1 */
1021 /* #if 'vdw' in INTERACTION_FLAGS[I][J] */
1022 /* ## Note special check for TIP4P-TIP4P. Since we are cutting of all hydrogen interactions we also cut the LJ-only O-O interaction */
1023 /* #if 'exactcutoff' in INTERACTION_FLAGS[I][J] or (GEOMETRY_I=='Water4' and GEOMETRY_J=='Water4' and 'exactcutoff' in INTERACTION_FLAGS[1][1]) */
1024 vvdw = _mm256_and_ps(vvdw,cutoff_mask);
1025 /* #define INNERFLOPS INNERFLOPS+1 */
1027 /* #if ROUND == 'Epilogue' */
1028 vvdw = _mm256_andnot_ps(dummy_mask,vvdw);
1030 vvdwsum = _mm256_add_ps(vvdwsum,vvdw);
1031 /* #define INNERFLOPS INNERFLOPS+1 */
1035 /* #if 'Force' in KERNEL_VF */
1037 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] and 'vdw' in INTERACTION_FLAGS[I][J] */
1038 fscal = _mm256_add_ps(felec,fvdw);
1039 /* #define INNERFLOPS INNERFLOPS+1 */
1040 /* #elif 'electrostatics' in INTERACTION_FLAGS[I][J] */
1042 /* #elif 'vdw' in INTERACTION_FLAGS[I][J] */
1046 /* ## Note special check for TIP4P-TIP4P. Since we are cutting of all hydrogen interactions we also cut the LJ-only O-O interaction */
1047 /* #if 'exactcutoff' in INTERACTION_FLAGS[I][J] or (GEOMETRY_I=='Water4' and GEOMETRY_J=='Water4' and 'exactcutoff' in INTERACTION_FLAGS[1][1]) */
1048 fscal = _mm256_and_ps(fscal,cutoff_mask);
1049 /* #define INNERFLOPS INNERFLOPS+1 */
1052 /* #if ROUND == 'Epilogue' */
1053 fscal = _mm256_andnot_ps(dummy_mask,fscal);
1056 /* Calculate temporary vectorial force */
1057 tx = _mm256_mul_ps(fscal,dx{I}{J});
1058 ty = _mm256_mul_ps(fscal,dy{I}{J});
1059 tz = _mm256_mul_ps(fscal,dz{I}{J});
1061 /* Update vectorial force */
1062 fix{I} = _mm256_add_ps(fix{I},tx);
1063 fiy{I} = _mm256_add_ps(fiy{I},ty);
1064 fiz{I} = _mm256_add_ps(fiz{I},tz);
1065 /* #define INNERFLOPS INNERFLOPS+6 */
1067 /* #if GEOMETRY_I == 'Particle' */
1068 /* #if ROUND == 'Loop' */
1069 fjptrA = f+j_coord_offsetA;
1070 fjptrB = f+j_coord_offsetB;
1071 fjptrC = f+j_coord_offsetC;
1072 fjptrD = f+j_coord_offsetD;
1073 fjptrE = f+j_coord_offsetE;
1074 fjptrF = f+j_coord_offsetF;
1075 fjptrG = f+j_coord_offsetG;
1076 fjptrH = f+j_coord_offsetH;
1078 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
1079 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
1080 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
1081 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
1082 fjptrE = (jnrlistE>=0) ? f+j_coord_offsetE : scratch;
1083 fjptrF = (jnrlistF>=0) ? f+j_coord_offsetF : scratch;
1084 fjptrG = (jnrlistG>=0) ? f+j_coord_offsetG : scratch;
1085 fjptrH = (jnrlistH>=0) ? f+j_coord_offsetH : scratch;
1087 gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,tx,ty,tz);
1088 /* #define INNERFLOPS INNERFLOPS+3 */
1090 fjx{J} = _mm256_add_ps(fjx{J},tx);
1091 fjy{J} = _mm256_add_ps(fjy{J},ty);
1092 fjz{J} = _mm256_add_ps(fjz{J},tz);
1093 /* #define INNERFLOPS INNERFLOPS+3 */
1098 /* ## Note special check for TIP4P-TIP4P. Since we are cutting of all hydrogen interactions we also cut the LJ-only O-O interaction */
1099 /* #if 'exactcutoff' in INTERACTION_FLAGS[I][J] or (GEOMETRY_I=='Water4' and GEOMETRY_J=='Water4' and 'exactcutoff' in INTERACTION_FLAGS[1][1]) */
1100 /* #if 0 ## This and next two lines is a hack to maintain indentation in template file */
1105 /* ## End of check for the interaction being outside the cutoff */
1108 /* ## End of loop over i-j interaction pairs */
1110 /* #if GEOMETRY_I != 'Particle' */
1111 /* #if ROUND == 'Loop' */
1112 fjptrA = f+j_coord_offsetA;
1113 fjptrB = f+j_coord_offsetB;
1114 fjptrC = f+j_coord_offsetC;
1115 fjptrD = f+j_coord_offsetD;
1116 fjptrE = f+j_coord_offsetE;
1117 fjptrF = f+j_coord_offsetF;
1118 fjptrG = f+j_coord_offsetG;
1119 fjptrH = f+j_coord_offsetH;
1121 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
1122 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
1123 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
1124 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
1125 fjptrE = (jnrlistE>=0) ? f+j_coord_offsetE : scratch;
1126 fjptrF = (jnrlistF>=0) ? f+j_coord_offsetF : scratch;
1127 fjptrG = (jnrlistG>=0) ? f+j_coord_offsetG : scratch;
1128 fjptrH = (jnrlistH>=0) ? f+j_coord_offsetH : scratch;
1132 /* #if 'Water' in GEOMETRY_I and GEOMETRY_J == 'Particle' */
1133 gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,fjx0,fjy0,fjz0);
1134 /* #define INNERFLOPS INNERFLOPS+3 */
1135 /* #elif GEOMETRY_J == 'Water3' */
1136 gmx_mm256_decrement_3rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,
1137 fjx0,fjy0,fjz0,fjx1,fjy1,fjz1,fjx2,fjy2,fjz2);
1138 /* #define INNERFLOPS INNERFLOPS+9 */
1139 /* #elif GEOMETRY_J == 'Water4' */
1140 /* #if 0 in PARTICLES_J */
1141 gmx_mm256_decrement_4rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,
1142 fjx0,fjy0,fjz0,fjx1,fjy1,fjz1,
1143 fjx2,fjy2,fjz2,fjx3,fjy3,fjz3);
1144 /* #define INNERFLOPS INNERFLOPS+12 */
1146 gmx_mm256_decrement_3rvec_8ptr_swizzle_ps(fjptrA+DIM,fjptrB+DIM,fjptrC+DIM,fjptrD+DIM,
1147 fjptrE+DIM,fjptrF+DIM,fjptrG+DIM,fjptrH+DIM,
1148 fjx1,fjy1,fjz1,fjx2,fjy2,fjz2,fjx3,fjy3,fjz3);
1149 /* #define INNERFLOPS INNERFLOPS+9 */
1153 /* Inner loop uses {INNERFLOPS} flops */
1158 /* End of innermost loop */
1160 /* #if 'Force' in KERNEL_VF */
1161 /* #if GEOMETRY_I == 'Particle' */
1162 gmx_mm256_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
1163 f+i_coord_offset,fshift+i_shift_offset);
1164 /* #define OUTERFLOPS OUTERFLOPS+6 */
1165 /* #elif GEOMETRY_I == 'Water3' */
1166 gmx_mm256_update_iforce_3atom_swizzle_ps(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,
1167 f+i_coord_offset,fshift+i_shift_offset);
1168 /* #define OUTERFLOPS OUTERFLOPS+18 */
1169 /* #elif GEOMETRY_I == 'Water4' */
1170 /* #if 0 in PARTICLES_I */
1171 gmx_mm256_update_iforce_4atom_swizzle_ps(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
1172 f+i_coord_offset,fshift+i_shift_offset);
1173 /* #define OUTERFLOPS OUTERFLOPS+24 */
1175 gmx_mm256_update_iforce_3atom_swizzle_ps(fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
1176 f+i_coord_offset+DIM,fshift+i_shift_offset);
1177 /* #define OUTERFLOPS OUTERFLOPS+18 */
1182 /* #if 'Potential' in KERNEL_VF */
1184 /* Update potential energies */
1185 /* #if KERNEL_ELEC != 'None' */
1186 gmx_mm256_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
1187 /* #define OUTERFLOPS OUTERFLOPS+1 */
1189 /* #if 'GeneralizedBorn' in KERNEL_ELEC */
1190 gmx_mm256_update_1pot_ps(vgbsum,kernel_data->energygrp_polarization+ggid);
1191 /* #define OUTERFLOPS OUTERFLOPS+1 */
1193 /* #if KERNEL_VDW != 'None' */
1194 gmx_mm256_update_1pot_ps(vvdwsum,kernel_data->energygrp_vdw+ggid);
1195 /* #define OUTERFLOPS OUTERFLOPS+1 */
1198 /* #if 'GeneralizedBorn' in KERNEL_ELEC and 'Force' in KERNEL_VF */
1199 dvdasum = _mm256_mul_ps(dvdasum, _mm256_mul_ps(isai{I},isai{I}));
1200 gmx_mm256_update_1pot_ps(dvdasum,dvda+inr);
1203 /* Increment number of inner iterations */
1204 inneriter += j_index_end - j_index_start;
1206 /* Outer loop uses {OUTERFLOPS} flops */
1209 /* Increment number of outer iterations */
1212 /* Update outer/inner flops */
1213 /* ## NB: This is not important, it just affects the flopcount. However, since our preprocessor is */
1214 /* ## primitive and replaces aggressively even in strings inside these directives, we need to */
1215 /* ## assemble the main part of the name (containing KERNEL/ELEC/VDW) directly in the source. */
1216 /* #if GEOMETRY_I == 'Water3' */
1217 /* #define ISUFFIX '_W3' */
1218 /* #elif GEOMETRY_I == 'Water4' */
1219 /* #define ISUFFIX '_W4' */
1221 /* #define ISUFFIX '' */
1223 /* #if GEOMETRY_J == 'Water3' */
1224 /* #define JSUFFIX 'W3' */
1225 /* #elif GEOMETRY_J == 'Water4' */
1226 /* #define JSUFFIX 'W4' */
1228 /* #define JSUFFIX '' */
1230 /* #if 'PotentialAndForce' in KERNEL_VF */
1231 /* #define VFSUFFIX '_VF' */
1232 /* #elif 'Potential' in KERNEL_VF */
1233 /* #define VFSUFFIX '_V' */
1235 /* #define VFSUFFIX '_F' */
1238 /* #if KERNEL_ELEC != 'None' and KERNEL_VDW != 'None' */
1239 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW{ISUFFIX}{JSUFFIX}{VFSUFFIX},outeriter*{OUTERFLOPS} + inneriter*{INNERFLOPS});
1240 /* #elif KERNEL_ELEC != 'None' */
1241 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC{ISUFFIX}{JSUFFIX}{VFSUFFIX},outeriter*{OUTERFLOPS} + inneriter*{INNERFLOPS});
1243 inc_nrnb(nrnb,eNR_NBKERNEL_VDW{ISUFFIX}{JSUFFIX}{VFSUFFIX},outeriter*{OUTERFLOPS} + inneriter*{INNERFLOPS});