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 */
45 #include "../nb_kernel.h"
46 #include "types/simple.h"
50 #include "gromacs/simd/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 */
91 * Gromacs nonbonded kernel: {KERNEL_NAME}
92 * Electrostatics interaction: {KERNEL_ELEC}
93 * VdW interaction: {KERNEL_VDW}
94 * Geometry: {GEOMETRY_I}-{GEOMETRY_J}
95 * Calculate force/pot: {KERNEL_VF}
99 (t_nblist * gmx_restrict nlist,
100 rvec * gmx_restrict xx,
101 rvec * gmx_restrict ff,
102 t_forcerec * gmx_restrict fr,
103 t_mdatoms * gmx_restrict mdatoms,
104 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
105 t_nrnb * gmx_restrict nrnb)
107 /* ## Not all variables are used for all kernels, but any optimizing compiler fixes that, */
108 /* ## so there is no point in going to extremes to exclude variables that are not needed. */
109 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
110 * just 0 for non-waters.
111 * Suffixes A,B,C,D,E,F,G,H refer to j loop unrolling done with AVX, e.g. for the eight different
112 * jnr indices corresponding to data put in the four positions in the SIMD register.
114 int i_shift_offset,i_coord_offset,outeriter,inneriter;
115 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
116 int jnrA,jnrB,jnrC,jnrD;
117 int jnrE,jnrF,jnrG,jnrH;
118 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
119 int jnrlistE,jnrlistF,jnrlistG,jnrlistH;
120 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
121 int j_coord_offsetE,j_coord_offsetF,j_coord_offsetG,j_coord_offsetH;
122 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
124 real *shiftvec,*fshift,*x,*f;
125 real *fjptrA,*fjptrB,*fjptrC,*fjptrD,*fjptrE,*fjptrF,*fjptrG,*fjptrH;
127 __m256 tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
128 /* #for I in PARTICLES_I */
129 real * vdwioffsetptr{I};
130 /* #if 'LJEwald' in KERNEL_VDW */
131 real * vdwgridioffsetptr{I};
133 __m256 ix{I},iy{I},iz{I},fix{I},fiy{I},fiz{I},iq{I},isai{I};
135 /* #for J in PARTICLES_J */
136 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;
137 __m256 jx{J},jy{J},jz{J},fjx{J},fjy{J},fjz{J},jq{J},isaj{J};
139 /* #for I,J in PAIRS_IJ */
140 __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};
142 /* #if KERNEL_ELEC != 'None' */
143 __m256 velec,felec,velecsum,facel,crf,krf,krf2;
146 /* #if 'GeneralizedBorn' in KERNEL_ELEC */
148 __m128i gbitab_lo,gbitab_hi;
149 __m256 vgb,fgb,vgbsum,dvdasum,gbscale,gbtabscale,isaprod,gbqqfactor,gbinvepsdiff,gbeps,dvdatmp;
150 __m256 minushalf = _mm256_set1_ps(-0.5);
151 real *invsqrta,*dvda,*gbtab;
153 /* #if KERNEL_VDW != 'None' */
155 __m256 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
158 __m256 one_sixth = _mm256_set1_ps(1.0/6.0);
159 __m256 one_twelfth = _mm256_set1_ps(1.0/12.0);
161 /* #if 'Table' in KERNEL_ELEC or 'GeneralizedBorn' in KERNEL_ELEC or 'Table' in KERNEL_VDW */
163 __m128i vfitab_lo,vfitab_hi;
164 __m128i ifour = _mm_set1_epi32(4);
165 __m256 rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF;
168 /* #if 'LJEwald' in KERNEL_VDW */
169 /* #for I,J in PAIRS_IJ */
170 __m256 c6grid_{I}{J};
173 __m256 ewclj,ewclj2,ewclj6,ewcljrsq,poly,exponent,f6A,f6B,sh_lj_ewald;
174 __m256 one_half = _mm256_set1_ps(0.5);
175 __m256 minus_one = _mm256_set1_ps(-1.0);
177 /* #if 'Ewald' in KERNEL_ELEC */
179 __m128i ewitab_lo,ewitab_hi;
180 __m256 ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
181 __m256 beta,beta2,beta3,zeta2,pmecorrF,pmecorrV,rinv3;
184 /* #if 'PotentialSwitch' in [KERNEL_MOD_ELEC,KERNEL_MOD_VDW] */
185 __m256 rswitch,swV3,swV4,swV5,swF2,swF3,swF4,d,d2,sw,dsw;
186 real rswitch_scalar,d_scalar;
188 __m256 dummy_mask,cutoff_mask;
189 __m256 signbit = _mm256_castsi256_ps( _mm256_set1_epi32(0x80000000) );
190 __m256 one = _mm256_set1_ps(1.0);
191 __m256 two = _mm256_set1_ps(2.0);
197 jindex = nlist->jindex;
199 shiftidx = nlist->shift;
201 shiftvec = fr->shift_vec[0];
202 fshift = fr->fshift[0];
203 /* #if KERNEL_ELEC != 'None' */
204 facel = _mm256_set1_ps(fr->epsfac);
205 charge = mdatoms->chargeA;
206 /* #if 'ReactionField' in KERNEL_ELEC */
207 krf = _mm256_set1_ps(fr->ic->k_rf);
208 krf2 = _mm256_set1_ps(fr->ic->k_rf*2.0);
209 crf = _mm256_set1_ps(fr->ic->c_rf);
212 /* #if KERNEL_VDW != 'None' */
213 nvdwtype = fr->ntype;
215 vdwtype = mdatoms->typeA;
217 /* #if 'LJEwald' in KERNEL_VDW */
218 vdwgridparam = fr->ljpme_c6grid;
219 sh_lj_ewald = _mm256_set1_ps(fr->ic->sh_lj_ewald);
220 ewclj = _mm256_set1_ps(fr->ewaldcoeff_lj);
221 ewclj2 = _mm256_mul_ps(minus_one,_mm256_mul_ps(ewclj,ewclj));
224 /* #if 'Table' in KERNEL_ELEC and 'Table' in KERNEL_VDW */
225 vftab = kernel_data->table_elec_vdw->data;
226 vftabscale = _mm256_set1_ps(kernel_data->table_elec_vdw->scale);
227 /* #elif 'Table' in KERNEL_ELEC */
228 vftab = kernel_data->table_elec->data;
229 vftabscale = _mm256_set1_ps(kernel_data->table_elec->scale);
230 /* #elif 'Table' in KERNEL_VDW */
231 vftab = kernel_data->table_vdw->data;
232 vftabscale = _mm256_set1_ps(kernel_data->table_vdw->scale);
235 /* #if 'Ewald' in KERNEL_ELEC */
236 sh_ewald = _mm256_set1_ps(fr->ic->sh_ewald);
237 beta = _mm256_set1_ps(fr->ic->ewaldcoeff_q);
238 beta2 = _mm256_mul_ps(beta,beta);
239 beta3 = _mm256_mul_ps(beta,beta2);
241 /* #if KERNEL_VF=='Force' and KERNEL_MOD_ELEC!='PotentialSwitch' */
242 ewtab = fr->ic->tabq_coul_F;
243 ewtabscale = _mm256_set1_ps(fr->ic->tabq_scale);
244 ewtabhalfspace = _mm256_set1_ps(0.5/fr->ic->tabq_scale);
246 ewtab = fr->ic->tabq_coul_FDV0;
247 ewtabscale = _mm256_set1_ps(fr->ic->tabq_scale);
248 ewtabhalfspace = _mm256_set1_ps(0.5/fr->ic->tabq_scale);
252 /* #if KERNEL_ELEC=='GeneralizedBorn' */
253 invsqrta = fr->invsqrta;
255 gbtabscale = _mm256_set1_ps(fr->gbtab.scale);
256 gbtab = fr->gbtab.data;
257 gbinvepsdiff = _mm256_set1_ps((1.0/fr->epsilon_r) - (1.0/fr->gb_epsilon_solvent));
260 /* #if 'Water' in GEOMETRY_I */
261 /* Setup water-specific parameters */
262 inr = nlist->iinr[0];
263 /* #for I in PARTICLES_ELEC_I */
264 iq{I} = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+{I}]));
266 /* #for I in PARTICLES_VDW_I */
267 vdwioffsetptr{I} = vdwparam+2*nvdwtype*vdwtype[inr+{I}];
268 /* #if 'LJEwald' in KERNEL_VDW */
269 vdwgridioffsetptr{I} = vdwgridparam+2*nvdwtype*vdwtype[inr+{I}];
274 /* #if 'Water' in GEOMETRY_J */
275 /* #for J in PARTICLES_ELEC_J */
276 jq{J} = _mm256_set1_ps(charge[inr+{J}]);
278 /* #for J in PARTICLES_VDW_J */
279 vdwjidx{J}A = 2*vdwtype[inr+{J}];
281 /* #for I,J in PAIRS_IJ */
282 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] */
283 qq{I}{J} = _mm256_mul_ps(iq{I},jq{J});
285 /* #if 'vdw' in INTERACTION_FLAGS[I][J] */
286 /* #if 'LJEwald' in KERNEL_VDW */
287 c6_{I}{J} = _mm256_set1_ps(vdwioffsetptr{I}[vdwjidx{J}A]);
288 c12_{I}{J} = _mm256_set1_ps(vdwioffsetptr{I}[vdwjidx{J}A+1]);
289 c6grid_{I}{J} = _mm256_set1_ps(vdwgridioffsetptr{I}[vdwjidx{J}A]);
291 c6_{I}{J} = _mm256_set1_ps(vdwioffsetptr{I}[vdwjidx{J}A]);
292 c12_{I}{J} = _mm256_set1_ps(vdwioffsetptr{I}[vdwjidx{J}A+1]);
298 /* #if KERNEL_MOD_ELEC!='None' or KERNEL_MOD_VDW!='None' */
299 /* #if KERNEL_ELEC!='None' */
300 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
301 rcutoff_scalar = fr->rcoulomb;
303 rcutoff_scalar = fr->rvdw;
305 rcutoff = _mm256_set1_ps(rcutoff_scalar);
306 rcutoff2 = _mm256_mul_ps(rcutoff,rcutoff);
309 /* #if KERNEL_MOD_VDW=='PotentialShift' */
310 sh_vdw_invrcut6 = _mm256_set1_ps(fr->ic->sh_invrc6);
311 rvdw = _mm256_set1_ps(fr->rvdw);
314 /* #if 'PotentialSwitch' in [KERNEL_MOD_ELEC,KERNEL_MOD_VDW] */
315 /* #if KERNEL_MOD_ELEC=='PotentialSwitch' */
316 rswitch_scalar = fr->rcoulomb_switch;
317 rswitch = _mm256_set1_ps(rswitch_scalar);
319 rswitch_scalar = fr->rvdw_switch;
320 rswitch = _mm256_set1_ps(rswitch_scalar);
322 /* Setup switch parameters */
323 d_scalar = rcutoff_scalar-rswitch_scalar;
324 d = _mm256_set1_ps(d_scalar);
325 swV3 = _mm256_set1_ps(-10.0/(d_scalar*d_scalar*d_scalar));
326 swV4 = _mm256_set1_ps( 15.0/(d_scalar*d_scalar*d_scalar*d_scalar));
327 swV5 = _mm256_set1_ps( -6.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
328 /* #if 'Force' in KERNEL_VF */
329 swF2 = _mm256_set1_ps(-30.0/(d_scalar*d_scalar*d_scalar));
330 swF3 = _mm256_set1_ps( 60.0/(d_scalar*d_scalar*d_scalar*d_scalar));
331 swF4 = _mm256_set1_ps(-30.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
335 /* Avoid stupid compiler warnings */
336 jnrA = jnrB = jnrC = jnrD = jnrE = jnrF = jnrG = jnrH = 0;
346 /* ## Keep track of the floating point operations we issue for reporting! */
347 /* #define OUTERFLOPS 0 */
351 for(iidx=0;iidx<4*DIM;iidx++)
356 /* Start outer loop over neighborlists */
357 for(iidx=0; iidx<nri; iidx++)
359 /* Load shift vector for this list */
360 i_shift_offset = DIM*shiftidx[iidx];
362 /* Load limits for loop over neighbors */
363 j_index_start = jindex[iidx];
364 j_index_end = jindex[iidx+1];
366 /* Get outer coordinate index */
368 i_coord_offset = DIM*inr;
370 /* Load i particle coords and add shift vector */
371 /* #if GEOMETRY_I == 'Particle' */
372 gmx_mm256_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
373 /* #elif GEOMETRY_I == 'Water3' */
374 gmx_mm256_load_shift_and_3rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,
375 &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2);
376 /* #elif GEOMETRY_I == 'Water4' */
377 /* #if 0 in PARTICLES_I */
378 gmx_mm256_load_shift_and_4rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,
379 &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
381 gmx_mm256_load_shift_and_3rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset+DIM,
382 &ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
386 /* #if 'Force' in KERNEL_VF */
387 /* #for I in PARTICLES_I */
388 fix{I} = _mm256_setzero_ps();
389 fiy{I} = _mm256_setzero_ps();
390 fiz{I} = _mm256_setzero_ps();
394 /* ## For water we already preloaded parameters at the start of the kernel */
395 /* #if not 'Water' in GEOMETRY_I */
396 /* Load parameters for i particles */
397 /* #for I in PARTICLES_ELEC_I */
398 iq{I} = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+{I}]));
399 /* #define OUTERFLOPS OUTERFLOPS+1 */
400 /* #if KERNEL_ELEC=='GeneralizedBorn' */
401 isai{I} = _mm256_set1_ps(invsqrta[inr+{I}]);
404 /* #for I in PARTICLES_VDW_I */
405 vdwioffsetptr{I} = vdwparam+2*nvdwtype*vdwtype[inr+{I}];
406 /* #if 'LJEwald' in KERNEL_VDW */
407 vdwgridioffsetptr{I} = vdwgridparam+2*nvdwtype*vdwtype[inr+{I}];
412 /* #if 'Potential' in KERNEL_VF */
413 /* Reset potential sums */
414 /* #if KERNEL_ELEC != 'None' */
415 velecsum = _mm256_setzero_ps();
417 /* #if 'GeneralizedBorn' in KERNEL_ELEC */
418 vgbsum = _mm256_setzero_ps();
420 /* #if KERNEL_VDW != 'None' */
421 vvdwsum = _mm256_setzero_ps();
424 /* #if 'GeneralizedBorn' in KERNEL_ELEC and 'Force' in KERNEL_VF */
425 dvdasum = _mm256_setzero_ps();
428 /* #for ROUND in ['Loop','Epilogue'] */
430 /* #if ROUND =='Loop' */
431 /* Start inner kernel loop */
432 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+7]>=0; jidx+=8)
434 /* ## First round is normal loop (next statement resets indentation) */
441 /* ## Second round is epilogue */
443 /* #define INNERFLOPS 0 */
445 /* Get j neighbor index, and coordinate index */
446 /* #if ROUND =='Loop' */
456 jnrlistA = jjnr[jidx];
457 jnrlistB = jjnr[jidx+1];
458 jnrlistC = jjnr[jidx+2];
459 jnrlistD = jjnr[jidx+3];
460 jnrlistE = jjnr[jidx+4];
461 jnrlistF = jjnr[jidx+5];
462 jnrlistG = jjnr[jidx+6];
463 jnrlistH = jjnr[jidx+7];
464 /* Sign of each element will be negative for non-real atoms.
465 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
466 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
468 dummy_mask = gmx_mm256_set_m128(gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx+4)),_mm_setzero_si128())),
469 gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128())));
471 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
472 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
473 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
474 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
475 jnrE = (jnrlistE>=0) ? jnrlistE : 0;
476 jnrF = (jnrlistF>=0) ? jnrlistF : 0;
477 jnrG = (jnrlistG>=0) ? jnrlistG : 0;
478 jnrH = (jnrlistH>=0) ? jnrlistH : 0;
480 j_coord_offsetA = DIM*jnrA;
481 j_coord_offsetB = DIM*jnrB;
482 j_coord_offsetC = DIM*jnrC;
483 j_coord_offsetD = DIM*jnrD;
484 j_coord_offsetE = DIM*jnrE;
485 j_coord_offsetF = DIM*jnrF;
486 j_coord_offsetG = DIM*jnrG;
487 j_coord_offsetH = DIM*jnrH;
489 /* load j atom coordinates */
490 /* #if GEOMETRY_J == 'Particle' */
491 gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
492 x+j_coord_offsetC,x+j_coord_offsetD,
493 x+j_coord_offsetE,x+j_coord_offsetF,
494 x+j_coord_offsetG,x+j_coord_offsetH,
496 /* #elif GEOMETRY_J == 'Water3' */
497 gmx_mm256_load_3rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
498 x+j_coord_offsetC,x+j_coord_offsetD,
499 x+j_coord_offsetE,x+j_coord_offsetF,
500 x+j_coord_offsetG,x+j_coord_offsetH,
501 &jx0,&jy0,&jz0,&jx1,&jy1,&jz1,&jx2,&jy2,&jz2);
502 /* #elif GEOMETRY_J == 'Water4' */
503 /* #if 0 in PARTICLES_J */
504 gmx_mm256_load_4rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
505 x+j_coord_offsetC,x+j_coord_offsetD,
506 x+j_coord_offsetE,x+j_coord_offsetF,
507 x+j_coord_offsetG,x+j_coord_offsetH,
508 &jx0,&jy0,&jz0,&jx1,&jy1,&jz1,&jx2,
509 &jy2,&jz2,&jx3,&jy3,&jz3);
511 gmx_mm256_load_3rvec_8ptr_swizzle_ps(x+j_coord_offsetA+DIM,x+j_coord_offsetB+DIM,
512 x+j_coord_offsetC+DIM,x+j_coord_offsetD+DIM,
513 x+j_coord_offsetE+DIM,x+j_coord_offsetF+DIM,
514 x+j_coord_offsetG+DIM,x+j_coord_offsetH+DIM,
515 &jx1,&jy1,&jz1,&jx2,&jy2,&jz2,&jx3,&jy3,&jz3);
519 /* Calculate displacement vector */
520 /* #for I,J in PAIRS_IJ */
521 dx{I}{J} = _mm256_sub_ps(ix{I},jx{J});
522 dy{I}{J} = _mm256_sub_ps(iy{I},jy{J});
523 dz{I}{J} = _mm256_sub_ps(iz{I},jz{J});
524 /* #define INNERFLOPS INNERFLOPS+3 */
527 /* Calculate squared distance and things based on it */
528 /* #for I,J in PAIRS_IJ */
529 rsq{I}{J} = gmx_mm256_calc_rsq_ps(dx{I}{J},dy{I}{J},dz{I}{J});
530 /* #define INNERFLOPS INNERFLOPS+5 */
533 /* #for I,J in PAIRS_IJ */
534 /* #if 'rinv' in INTERACTION_FLAGS[I][J] */
535 rinv{I}{J} = gmx_mm256_invsqrt_ps(rsq{I}{J});
536 /* #define INNERFLOPS INNERFLOPS+5 */
540 /* #for I,J in PAIRS_IJ */
541 /* #if 'rinvsq' in INTERACTION_FLAGS[I][J] */
542 /* # if 'rinv' not in INTERACTION_FLAGS[I][J] */
543 rinvsq{I}{J} = gmx_mm256_inv_ps(rsq{I}{J});
544 /* #define INNERFLOPS INNERFLOPS+4 */
546 rinvsq{I}{J} = _mm256_mul_ps(rinv{I}{J},rinv{I}{J});
547 /* #define INNERFLOPS INNERFLOPS+1 */
552 /* #if not 'Water' in GEOMETRY_J */
553 /* Load parameters for j particles */
554 /* #for J in PARTICLES_ELEC_J */
555 jq{J} = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+{J},charge+jnrB+{J},
556 charge+jnrC+{J},charge+jnrD+{J},
557 charge+jnrE+{J},charge+jnrF+{J},
558 charge+jnrG+{J},charge+jnrH+{J});
559 /* #if KERNEL_ELEC=='GeneralizedBorn' */
560 isaj{J} = gmx_mm256_load_8real_swizzle_ps(invsqrta+jnrA+{J},invsqrta+jnrB+{J},
561 invsqrta+jnrC+{J},invsqrta+jnrD+{J},
562 invsqrta+jnrE+{J},invsqrta+jnrF+{J},
563 invsqrta+jnrG+{J},invsqrta+jnrH+{J});
566 /* #for J in PARTICLES_VDW_J */
567 vdwjidx{J}A = 2*vdwtype[jnrA+{J}];
568 vdwjidx{J}B = 2*vdwtype[jnrB+{J}];
569 vdwjidx{J}C = 2*vdwtype[jnrC+{J}];
570 vdwjidx{J}D = 2*vdwtype[jnrD+{J}];
571 vdwjidx{J}E = 2*vdwtype[jnrE+{J}];
572 vdwjidx{J}F = 2*vdwtype[jnrF+{J}];
573 vdwjidx{J}G = 2*vdwtype[jnrG+{J}];
574 vdwjidx{J}H = 2*vdwtype[jnrH+{J}];
578 /* #if 'Force' in KERNEL_VF and not 'Particle' in GEOMETRY_I */
579 /* #for J in PARTICLES_J */
580 fjx{J} = _mm256_setzero_ps();
581 fjy{J} = _mm256_setzero_ps();
582 fjz{J} = _mm256_setzero_ps();
586 /* #for I,J in PAIRS_IJ */
588 /**************************
589 * CALCULATE INTERACTIONS *
590 **************************/
592 /* ## Note special check for TIP4P-TIP4P. Since we are cutting of all hydrogen interactions we also cut the LJ-only O-O interaction */
593 /* #if 'exactcutoff' in INTERACTION_FLAGS[I][J] or (GEOMETRY_I=='Water4' and GEOMETRY_J=='Water4' and 'exactcutoff' in INTERACTION_FLAGS[1][1]) */
594 /* ## We always calculate rinv/rinvsq above to enable pipelineing in compilers (performance tested on x86) */
595 if (gmx_mm256_any_lt(rsq{I}{J},rcutoff2))
597 /* #if 0 ## this and the next two lines is a hack to maintain auto-indentation in template file */
600 /* #define INNERFLOPS INNERFLOPS+1 */
603 /* #if 'r' in INTERACTION_FLAGS[I][J] */
604 r{I}{J} = _mm256_mul_ps(rsq{I}{J},rinv{I}{J});
605 /* #if ROUND == 'Epilogue' */
606 r{I}{J} = _mm256_andnot_ps(dummy_mask,r{I}{J});
607 /* #define INNERFLOPS INNERFLOPS+1 */
609 /* #define INNERFLOPS INNERFLOPS+1 */
612 /* ## For water geometries we already loaded parameters at the start of the kernel */
613 /* #if not 'Water' in GEOMETRY_J */
614 /* Compute parameters for interactions between i and j atoms */
615 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] */
616 qq{I}{J} = _mm256_mul_ps(iq{I},jq{J});
617 /* #define INNERFLOPS INNERFLOPS+1 */
619 /* #if 'vdw' in INTERACTION_FLAGS[I][J] */
620 gmx_mm256_load_8pair_swizzle_ps(vdwioffsetptr{I}+vdwjidx{J}A,
621 vdwioffsetptr{I}+vdwjidx{J}B,
622 vdwioffsetptr{I}+vdwjidx{J}C,
623 vdwioffsetptr{I}+vdwjidx{J}D,
624 vdwioffsetptr{I}+vdwjidx{J}E,
625 vdwioffsetptr{I}+vdwjidx{J}F,
626 vdwioffsetptr{I}+vdwjidx{J}G,
627 vdwioffsetptr{I}+vdwjidx{J}H,
628 &c6_{I}{J},&c12_{I}{J});
630 /* #if 'LJEwald' in KERNEL_VDW */
631 c6grid_{I}{J} = gmx_mm256_load_8real_swizzle_ps(vdwgridioffsetptr{I}+vdwjidx{J}A,
632 vdwgridioffsetptr{I}+vdwjidx{J}B,
633 vdwgridioffsetptr{I}+vdwjidx{J}C,
634 vdwgridioffsetptr{I}+vdwjidx{J}D,
635 vdwgridioffsetptr{I}+vdwjidx{J}E,
636 vdwgridioffsetptr{I}+vdwjidx{J}F,
637 vdwgridioffsetptr{I}+vdwjidx{J}G,
638 vdwgridioffsetptr{I}+vdwjidx{J}H);
643 /* #if 'table' in INTERACTION_FLAGS[I][J] */
644 /* Calculate table index by multiplying r with table scale and truncate to integer */
645 rt = _mm256_mul_ps(r{I}{J},vftabscale);
646 vfitab = _mm256_cvttps_epi32(rt);
647 vfeps = _mm256_sub_ps(rt,_mm256_round_ps(rt, _MM_FROUND_FLOOR));
648 /* #define INNERFLOPS INNERFLOPS+4 */
649 /* AVX1 does not support 256-bit integer operations, so now we go to 128-bit mode... */
650 vfitab_lo = _mm256_extractf128_si256(vfitab,0x0);
651 vfitab_hi = _mm256_extractf128_si256(vfitab,0x1);
652 /* #if 'Table' in KERNEL_ELEC and 'Table' in KERNEL_VDW */
653 /* ## 3 tables, 4 bytes per point: multiply index by 12 */
654 vfitab_lo = _mm_slli_epi32(_mm_add_epi32(vfitab_lo,_mm_slli_epi32(vfitab_lo,1)),2);
655 vfitab_hi = _mm_slli_epi32(_mm_add_epi32(vfitab_hi,_mm_slli_epi32(vfitab_hi,1)),2);
656 /* #elif 'Table' in KERNEL_ELEC */
657 /* ## 1 table, 4 bytes per point: multiply index by 4 */
658 vfitab_lo = _mm_slli_epi32(vfitab_lo,2);
659 vfitab_hi = _mm_slli_epi32(vfitab_hi,2);
660 /* #elif 'Table' in KERNEL_VDW */
661 /* ## 2 tables, 4 bytes per point: multiply index by 8 */
662 vfitab_lo = _mm_slli_epi32(vfitab_lo,3);
663 vfitab_hi = _mm_slli_epi32(vfitab_hi,3);
667 /* ## ELECTROSTATIC INTERACTIONS */
668 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] */
670 /* #if KERNEL_ELEC=='Coulomb' */
672 /* COULOMB ELECTROSTATICS */
673 velec = _mm256_mul_ps(qq{I}{J},rinv{I}{J});
674 /* #define INNERFLOPS INNERFLOPS+1 */
675 /* #if 'Force' in KERNEL_VF */
676 felec = _mm256_mul_ps(velec,rinvsq{I}{J});
677 /* #define INNERFLOPS INNERFLOPS+1 */
680 /* #elif KERNEL_ELEC=='ReactionField' */
682 /* REACTION-FIELD ELECTROSTATICS */
683 /* #if 'Potential' in KERNEL_VF */
684 velec = _mm256_mul_ps(qq{I}{J},_mm256_sub_ps(_mm256_add_ps(rinv{I}{J},_mm256_mul_ps(krf,rsq{I}{J})),crf));
685 /* #define INNERFLOPS INNERFLOPS+4 */
687 /* #if 'Force' in KERNEL_VF */
688 felec = _mm256_mul_ps(qq{I}{J},_mm256_sub_ps(_mm256_mul_ps(rinv{I}{J},rinvsq{I}{J}),krf2));
689 /* #define INNERFLOPS INNERFLOPS+3 */
692 /* #elif KERNEL_ELEC=='GeneralizedBorn' */
694 /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
695 isaprod = _mm256_mul_ps(isai{I},isaj{J});
696 gbqqfactor = _mm256_xor_ps(signbit,_mm256_mul_ps(qq{I}{J},_mm256_mul_ps(isaprod,gbinvepsdiff)));
697 gbscale = _mm256_mul_ps(isaprod,gbtabscale);
698 /* #define INNERFLOPS INNERFLOPS+5 */
700 /* Calculate generalized born table index - this is a separate table from the normal one,
701 * but we use the same procedure by multiplying r with scale and truncating to integer.
703 rt = _mm256_mul_ps(r{I}{J},gbscale);
704 gbitab = _mm256_cvttps_epi32(rt);
705 gbeps = _mm256_sub_ps(rt,_mm256_round_ps(rt, _MM_FROUND_FLOOR));
706 /* AVX1 does not support 256-bit integer operations, so now we go to 128-bit mode... */
707 gbitab_lo = _mm256_extractf128_si256(gbitab,0x0);
708 gbitab_hi = _mm256_extractf128_si256(gbitab,0x1);
709 gbitab_lo = _mm_slli_epi32(gbitab_lo,2);
710 gbitab_hi = _mm_slli_epi32(gbitab_hi,2);
711 Y = gmx_mm256_set_m128(_mm_load_ps(gbtab + _mm_extract_epi32(gbitab_hi,0)),
712 _mm_load_ps(gbtab + _mm_extract_epi32(gbitab_lo,0)));
713 F = gmx_mm256_set_m128(_mm_load_ps(gbtab + _mm_extract_epi32(gbitab_hi,1)),
714 _mm_load_ps(gbtab + _mm_extract_epi32(gbitab_lo,1)));
715 G = gmx_mm256_set_m128(_mm_load_ps(gbtab + _mm_extract_epi32(gbitab_hi,2)),
716 _mm_load_ps(gbtab + _mm_extract_epi32(gbitab_lo,2)));
717 H = gmx_mm256_set_m128(_mm_load_ps(gbtab + _mm_extract_epi32(gbitab_hi,3)),
718 _mm_load_ps(gbtab + _mm_extract_epi32(gbitab_lo,3)));
719 GMX_MM256_HALFTRANSPOSE4_PS(Y,F,G,H);
720 Heps = _mm256_mul_ps(gbeps,H);
721 Fp = _mm256_add_ps(F,_mm256_mul_ps(gbeps,_mm256_add_ps(G,Heps)));
722 VV = _mm256_add_ps(Y,_mm256_mul_ps(gbeps,Fp));
723 vgb = _mm256_mul_ps(gbqqfactor,VV);
724 /* #define INNERFLOPS INNERFLOPS+10 */
726 /* #if 'Force' in KERNEL_VF */
727 FF = _mm256_add_ps(Fp,_mm256_mul_ps(gbeps,_mm256_add_ps(G,_mm256_add_ps(Heps,Heps))));
728 fgb = _mm256_mul_ps(gbqqfactor,_mm256_mul_ps(FF,gbscale));
729 dvdatmp = _mm256_mul_ps(minushalf,_mm256_add_ps(vgb,_mm256_mul_ps(fgb,r{I}{J})));
730 /* #if ROUND == 'Epilogue' */
731 dvdatmp = _mm256_andnot_ps(dummy_mask,dvdatmp);
733 dvdasum = _mm256_add_ps(dvdasum,dvdatmp);
734 /* #if ROUND == 'Loop' */
744 /* 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. */
745 fjptrA = (jnrlistA>=0) ? dvda+jnrA : scratch;
746 fjptrB = (jnrlistB>=0) ? dvda+jnrB : scratch;
747 fjptrC = (jnrlistC>=0) ? dvda+jnrC : scratch;
748 fjptrD = (jnrlistD>=0) ? dvda+jnrD : scratch;
749 fjptrE = (jnrlistE>=0) ? dvda+jnrE : scratch;
750 fjptrF = (jnrlistF>=0) ? dvda+jnrF : scratch;
751 fjptrG = (jnrlistG>=0) ? dvda+jnrG : scratch;
752 fjptrH = (jnrlistH>=0) ? dvda+jnrH : scratch;
754 gmx_mm256_increment_8real_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,
755 _mm256_mul_ps(dvdatmp,_mm256_mul_ps(isaj{J},isaj{J})));
756 /* #define INNERFLOPS INNERFLOPS+12 */
758 velec = _mm256_mul_ps(qq{I}{J},rinv{I}{J});
759 /* #define INNERFLOPS INNERFLOPS+1 */
760 /* #if 'Force' in KERNEL_VF */
761 felec = _mm256_mul_ps(_mm256_sub_ps(_mm256_mul_ps(velec,rinv{I}{J}),fgb),rinv{I}{J});
762 /* #define INNERFLOPS INNERFLOPS+3 */
765 /* #elif KERNEL_ELEC=='Ewald' */
766 /* EWALD ELECTROSTATICS */
768 /* Analytical PME correction */
769 zeta2 = _mm256_mul_ps(beta2,rsq{I}{J});
770 /* #if 'Force' in KERNEL_VF */
771 rinv3 = _mm256_mul_ps(rinvsq{I}{J},rinv{I}{J});
772 pmecorrF = gmx_mm256_pmecorrF_ps(zeta2);
773 felec = _mm256_add_ps( _mm256_mul_ps(pmecorrF,beta3), rinv3);
774 felec = _mm256_mul_ps(qq{I}{J},felec);
775 /* #define INNERFLOPS INNERFLOPS+31 */
777 /* #if 'Potential' in KERNEL_VF or KERNEL_MOD_ELEC=='PotentialSwitch' */
778 pmecorrV = gmx_mm256_pmecorrV_ps(zeta2);
779 pmecorrV = _mm256_mul_ps(pmecorrV,beta);
780 /* #define INNERFLOPS INNERFLOPS+27 */
781 /* #if KERNEL_MOD_ELEC=='PotentialShift' */
782 velec = _mm256_sub_ps(_mm256_sub_ps(rinv{I}{J},sh_ewald),pmecorrV);
783 /* #define INNERFLOPS INNERFLOPS+21 */
785 velec = _mm256_sub_ps(rinv{I}{J},pmecorrV);
787 velec = _mm256_mul_ps(qq{I}{J},velec);
790 /* #elif KERNEL_ELEC=='CubicSplineTable' */
792 /* CUBIC SPLINE TABLE ELECTROSTATICS */
793 Y = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,0)),
794 _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,0)));
795 F = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,1)),
796 _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,1)));
797 G = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,2)),
798 _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,2)));
799 H = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,3)),
800 _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,3)));
801 GMX_MM256_HALFTRANSPOSE4_PS(Y,F,G,H);
802 Heps = _mm256_mul_ps(vfeps,H);
803 Fp = _mm256_add_ps(F,_mm256_mul_ps(vfeps,_mm256_add_ps(G,Heps)));
804 /* #define INNERFLOPS INNERFLOPS+4 */
805 /* #if 'Potential' in KERNEL_VF */
806 VV = _mm256_add_ps(Y,_mm256_mul_ps(vfeps,Fp));
807 velec = _mm256_mul_ps(qq{I}{J},VV);
808 /* #define INNERFLOPS INNERFLOPS+3 */
810 /* #if 'Force' in KERNEL_VF */
811 FF = _mm256_add_ps(Fp,_mm256_mul_ps(vfeps,_mm256_add_ps(G,_mm256_add_ps(Heps,Heps))));
812 felec = _mm256_xor_ps(signbit,_mm256_mul_ps(_mm256_mul_ps(qq{I}{J},FF),_mm256_mul_ps(vftabscale,rinv{I}{J})));
813 /* #define INNERFLOPS INNERFLOPS+7 */
816 /* ## End of check for electrostatics interaction forms */
818 /* ## END OF ELECTROSTATIC INTERACTION CHECK FOR PAIR I-J */
820 /* #if 'vdw' in INTERACTION_FLAGS[I][J] */
822 /* #if KERNEL_VDW=='LennardJones' */
824 /* LENNARD-JONES DISPERSION/REPULSION */
826 rinvsix = _mm256_mul_ps(_mm256_mul_ps(rinvsq{I}{J},rinvsq{I}{J}),rinvsq{I}{J});
827 /* #define INNERFLOPS INNERFLOPS+2 */
828 /* #if 'Potential' in KERNEL_VF or KERNEL_MOD_VDW=='PotentialSwitch' */
829 vvdw6 = _mm256_mul_ps(c6_{I}{J},rinvsix);
830 vvdw12 = _mm256_mul_ps(c12_{I}{J},_mm256_mul_ps(rinvsix,rinvsix));
831 /* #define INNERFLOPS INNERFLOPS+3 */
832 /* #if KERNEL_MOD_VDW=='PotentialShift' */
833 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) ,
834 _mm256_mul_ps( _mm256_sub_ps(vvdw6,_mm256_mul_ps(c6_{I}{J},sh_vdw_invrcut6)),one_sixth));
835 /* #define INNERFLOPS INNERFLOPS+8 */
837 vvdw = _mm256_sub_ps( _mm256_mul_ps(vvdw12,one_twelfth) , _mm256_mul_ps(vvdw6,one_sixth) );
838 /* #define INNERFLOPS INNERFLOPS+3 */
840 /* ## Check for force inside potential check, i.e. this means we already did the potential part */
841 /* #if 'Force' in KERNEL_VF */
842 fvdw = _mm256_mul_ps(_mm256_sub_ps(vvdw12,vvdw6),rinvsq{I}{J});
843 /* #define INNERFLOPS INNERFLOPS+2 */
845 /* #elif KERNEL_VF=='Force' */
846 /* ## Force-only LennardJones makes it possible to save 1 flop (they do add up...) */
847 fvdw = _mm256_mul_ps(_mm256_sub_ps(_mm256_mul_ps(c12_{I}{J},rinvsix),c6_{I}{J}),_mm256_mul_ps(rinvsix,rinvsq{I}{J}));
848 /* #define INNERFLOPS INNERFLOPS+4 */
851 /* #elif KERNEL_VDW=='CubicSplineTable' */
853 /* CUBIC SPLINE TABLE DISPERSION */
854 /* #if 'Table' in KERNEL_ELEC */
855 vfitab_lo = _mm_add_epi32(vfitab_lo,ifour);
856 vfitab_hi = _mm_add_epi32(vfitab_hi,ifour);
858 Y = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,0)),
859 _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,0)));
860 F = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,1)),
861 _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,1)));
862 G = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,2)),
863 _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,2)));
864 H = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,3)),
865 _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,3)));
866 GMX_MM256_HALFTRANSPOSE4_PS(Y,F,G,H);
867 Heps = _mm256_mul_ps(vfeps,H);
868 Fp = _mm256_add_ps(F,_mm256_mul_ps(vfeps,_mm256_add_ps(G,Heps)));
869 /* #define INNERFLOPS INNERFLOPS+4 */
870 /* #if 'Potential' in KERNEL_VF */
871 VV = _mm256_add_ps(Y,_mm256_mul_ps(vfeps,Fp));
872 vvdw6 = _mm256_mul_ps(c6_{I}{J},VV);
873 /* #define INNERFLOPS INNERFLOPS+3 */
875 /* #if 'Force' in KERNEL_VF */
876 FF = _mm256_add_ps(Fp,_mm256_mul_ps(vfeps,_mm256_add_ps(G,_mm256_add_ps(Heps,Heps))));
877 fvdw6 = _mm256_mul_ps(c6_{I}{J},FF);
878 /* #define INNERFLOPS INNERFLOPS+4 */
881 /* CUBIC SPLINE TABLE REPULSION */
882 vfitab_lo = _mm_add_epi32(vfitab_lo,ifour);
883 vfitab_hi = _mm_add_epi32(vfitab_hi,ifour);
884 Y = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,0)),
885 _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,0)));
886 F = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,1)),
887 _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,1)));
888 G = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,2)),
889 _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,2)));
890 H = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,3)),
891 _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,3)));
892 GMX_MM256_HALFTRANSPOSE4_PS(Y,F,G,H);
893 Heps = _mm256_mul_ps(vfeps,H);
894 Fp = _mm256_add_ps(F,_mm256_mul_ps(vfeps,_mm256_add_ps(G,Heps)));
895 /* #define INNERFLOPS INNERFLOPS+4 */
896 /* #if 'Potential' in KERNEL_VF */
897 VV = _mm256_add_ps(Y,_mm256_mul_ps(vfeps,Fp));
898 vvdw12 = _mm256_mul_ps(c12_{I}{J},VV);
899 /* #define INNERFLOPS INNERFLOPS+3 */
901 /* #if 'Force' in KERNEL_VF */
902 FF = _mm256_add_ps(Fp,_mm256_mul_ps(vfeps,_mm256_add_ps(G,_mm256_add_ps(Heps,Heps))));
903 fvdw12 = _mm256_mul_ps(c12_{I}{J},FF);
904 /* #define INNERFLOPS INNERFLOPS+5 */
906 /* #if 'Potential' in KERNEL_VF */
907 vvdw = _mm256_add_ps(vvdw12,vvdw6);
908 /* #define INNERFLOPS INNERFLOPS+1 */
910 /* #if 'Force' in KERNEL_VF */
911 fvdw = _mm256_xor_ps(signbit,_mm256_mul_ps(_mm256_add_ps(fvdw6,fvdw12),_mm256_mul_ps(vftabscale,rinv{I}{J})));
912 /* #define INNERFLOPS INNERFLOPS+4 */
915 /* #elif KERNEL_VDW=='LJEwald' */
917 /* Analytical LJ-PME */
918 rinvsix = _mm256_mul_ps(_mm256_mul_ps(rinvsq{I}{J},rinvsq{I}{J}),rinvsq{I}{J});
919 ewcljrsq = _mm256_mul_ps(ewclj2,rsq{I}{J});
920 ewclj6 = _mm256_mul_ps(ewclj2,_mm256_mul_ps(ewclj2,ewclj2));
921 exponent = gmx_simd_exp_r(ewcljrsq);
922 /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
923 poly = _mm256_mul_ps(exponent,_mm256_add_ps(_mm256_sub_ps(one,ewcljrsq),_mm256_mul_ps(_mm256_mul_ps(ewcljrsq,ewcljrsq),one_half)));
924 /* #define INNERFLOPS INNERFLOPS+11 */
925 /* #if 'Potential' in KERNEL_VF or KERNEL_MOD_VDW=='PotentialSwitch' */
926 /* vvdw6 = [C6 - C6grid * (1-poly)]/r6 */
927 vvdw6 = _mm256_mul_ps(_mm256_sub_ps(c6_{I}{J},_mm256_mul_ps(c6grid_{I}{J},_mm256_sub_ps(one,poly))),rinvsix);
928 vvdw12 = _mm256_mul_ps(c12_{I}{J},_mm256_mul_ps(rinvsix,rinvsix));
929 /* #define INNERFLOPS INNERFLOPS+6 */
930 /* #if KERNEL_MOD_VDW=='PotentialShift' */
931 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) ,
932 _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));
933 /* #define INNERFLOPS INNERFLOPS+10 */
935 vvdw = _mm256_sub_ps(_mm256_mul_ps(vvdw12,one_twelfth),_mm256_mul_ps(vvdw6,one_sixth));
936 /* #define INNERFLOPS INNERFLOPS+3 */
938 /* ## Check for force inside potential check, i.e. this means we already did the potential part */
939 /* #if 'Force' in KERNEL_VF */
940 /* fvdw = vvdw12/r - (vvdw6/r + (C6grid * exponent * beta^6)/r) */
941 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});
942 /* #define INNERFLOPS INNERFLOPS+6 */
944 /* #elif KERNEL_VF=='Force' */
945 /* f6A = 6 * C6grid * (1 - poly) */
946 f6A = _mm256_mul_ps(c6grid_{I}{J},_mm256_sub_ps(one,poly));
947 /* f6B = C6grid * exponent * beta^6 */
948 f6B = _mm256_mul_ps(_mm256_mul_ps(c6grid_{I}{J},one_sixth),_mm256_mul_ps(exponent,ewclj6));
949 /* fvdw = 12*C12/r13 - ((6*C6 - f6A)/r6 + f6B)/r */
950 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});
951 /* #define INNERFLOPS INNERFLOPS+11 */
954 /* ## End of check for vdw interaction forms */
956 /* ## END OF VDW INTERACTION CHECK FOR PAIR I-J */
958 /* #if 'switch' in INTERACTION_FLAGS[I][J] */
959 d = _mm256_sub_ps(r{I}{J},rswitch);
960 d = _mm256_max_ps(d,_mm256_setzero_ps());
961 d2 = _mm256_mul_ps(d,d);
962 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)))))));
963 /* #define INNERFLOPS INNERFLOPS+10 */
965 /* #if 'Force' in KERNEL_VF */
966 dsw = _mm256_mul_ps(d2,_mm256_add_ps(swF2,_mm256_mul_ps(d,_mm256_add_ps(swF3,_mm256_mul_ps(d,swF4)))));
967 /* #define INNERFLOPS INNERFLOPS+5 */
970 /* Evaluate switch function */
971 /* #if 'Force' in KERNEL_VF */
972 /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
973 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] and KERNEL_MOD_ELEC=='PotentialSwitch' */
974 felec = _mm256_sub_ps( _mm256_mul_ps(felec,sw) , _mm256_mul_ps(rinv{I}{J},_mm256_mul_ps(velec,dsw)) );
975 /* #define INNERFLOPS INNERFLOPS+4 */
977 /* #if 'vdw' in INTERACTION_FLAGS[I][J] and KERNEL_MOD_VDW=='PotentialSwitch' */
978 fvdw = _mm256_sub_ps( _mm256_mul_ps(fvdw,sw) , _mm256_mul_ps(rinv{I}{J},_mm256_mul_ps(vvdw,dsw)) );
979 /* #define INNERFLOPS INNERFLOPS+4 */
982 /* #if 'Potential' in KERNEL_VF */
983 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] and KERNEL_MOD_ELEC=='PotentialSwitch' */
984 velec = _mm256_mul_ps(velec,sw);
985 /* #define INNERFLOPS INNERFLOPS+1 */
987 /* #if 'vdw' in INTERACTION_FLAGS[I][J] and KERNEL_MOD_VDW=='PotentialSwitch' */
988 vvdw = _mm256_mul_ps(vvdw,sw);
989 /* #define INNERFLOPS INNERFLOPS+1 */
993 /* ## Note special check for TIP4P-TIP4P. Since we are cutting of all hydrogen interactions we also cut the LJ-only O-O interaction */
994 /* #if 'exactcutoff' in INTERACTION_FLAGS[I][J] or (GEOMETRY_I=='Water4' and GEOMETRY_J=='Water4' and 'exactcutoff' in INTERACTION_FLAGS[1][1]) */
995 cutoff_mask = _mm256_cmp_ps(rsq{I}{J},rcutoff2,_CMP_LT_OQ);
996 /* #define INNERFLOPS INNERFLOPS+1 */
999 /* #if 'Potential' in KERNEL_VF */
1000 /* Update potential sum for this i atom from the interaction with this j atom. */
1001 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] */
1002 /* #if 'exactcutoff' in INTERACTION_FLAGS[I][J] */
1003 velec = _mm256_and_ps(velec,cutoff_mask);
1004 /* #define INNERFLOPS INNERFLOPS+1 */
1006 /* #if ROUND == 'Epilogue' */
1007 velec = _mm256_andnot_ps(dummy_mask,velec);
1009 velecsum = _mm256_add_ps(velecsum,velec);
1010 /* #define INNERFLOPS INNERFLOPS+1 */
1011 /* #if KERNEL_ELEC=='GeneralizedBorn' */
1012 /* #if 'exactcutoff' in INTERACTION_FLAGS[I][J] */
1013 vgb = _mm256_and_ps(vgb,cutoff_mask);
1014 /* #define INNERFLOPS INNERFLOPS+1 */
1016 /* #if ROUND == 'Epilogue' */
1017 vgb = _mm256_andnot_ps(dummy_mask,vgb);
1019 vgbsum = _mm256_add_ps(vgbsum,vgb);
1020 /* #define INNERFLOPS INNERFLOPS+1 */
1023 /* #if 'vdw' in INTERACTION_FLAGS[I][J] */
1024 /* ## Note special check for TIP4P-TIP4P. Since we are cutting of all hydrogen interactions we also cut the LJ-only O-O interaction */
1025 /* #if 'exactcutoff' in INTERACTION_FLAGS[I][J] or (GEOMETRY_I=='Water4' and GEOMETRY_J=='Water4' and 'exactcutoff' in INTERACTION_FLAGS[1][1]) */
1026 vvdw = _mm256_and_ps(vvdw,cutoff_mask);
1027 /* #define INNERFLOPS INNERFLOPS+1 */
1029 /* #if ROUND == 'Epilogue' */
1030 vvdw = _mm256_andnot_ps(dummy_mask,vvdw);
1032 vvdwsum = _mm256_add_ps(vvdwsum,vvdw);
1033 /* #define INNERFLOPS INNERFLOPS+1 */
1037 /* #if 'Force' in KERNEL_VF */
1039 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] and 'vdw' in INTERACTION_FLAGS[I][J] */
1040 fscal = _mm256_add_ps(felec,fvdw);
1041 /* #define INNERFLOPS INNERFLOPS+1 */
1042 /* #elif 'electrostatics' in INTERACTION_FLAGS[I][J] */
1044 /* #elif 'vdw' in INTERACTION_FLAGS[I][J] */
1048 /* ## Note special check for TIP4P-TIP4P. Since we are cutting of all hydrogen interactions we also cut the LJ-only O-O interaction */
1049 /* #if 'exactcutoff' in INTERACTION_FLAGS[I][J] or (GEOMETRY_I=='Water4' and GEOMETRY_J=='Water4' and 'exactcutoff' in INTERACTION_FLAGS[1][1]) */
1050 fscal = _mm256_and_ps(fscal,cutoff_mask);
1051 /* #define INNERFLOPS INNERFLOPS+1 */
1054 /* #if ROUND == 'Epilogue' */
1055 fscal = _mm256_andnot_ps(dummy_mask,fscal);
1058 /* Calculate temporary vectorial force */
1059 tx = _mm256_mul_ps(fscal,dx{I}{J});
1060 ty = _mm256_mul_ps(fscal,dy{I}{J});
1061 tz = _mm256_mul_ps(fscal,dz{I}{J});
1063 /* Update vectorial force */
1064 fix{I} = _mm256_add_ps(fix{I},tx);
1065 fiy{I} = _mm256_add_ps(fiy{I},ty);
1066 fiz{I} = _mm256_add_ps(fiz{I},tz);
1067 /* #define INNERFLOPS INNERFLOPS+6 */
1069 /* #if GEOMETRY_I == 'Particle' */
1070 /* #if ROUND == 'Loop' */
1071 fjptrA = f+j_coord_offsetA;
1072 fjptrB = f+j_coord_offsetB;
1073 fjptrC = f+j_coord_offsetC;
1074 fjptrD = f+j_coord_offsetD;
1075 fjptrE = f+j_coord_offsetE;
1076 fjptrF = f+j_coord_offsetF;
1077 fjptrG = f+j_coord_offsetG;
1078 fjptrH = f+j_coord_offsetH;
1080 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
1081 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
1082 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
1083 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
1084 fjptrE = (jnrlistE>=0) ? f+j_coord_offsetE : scratch;
1085 fjptrF = (jnrlistF>=0) ? f+j_coord_offsetF : scratch;
1086 fjptrG = (jnrlistG>=0) ? f+j_coord_offsetG : scratch;
1087 fjptrH = (jnrlistH>=0) ? f+j_coord_offsetH : scratch;
1089 gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,tx,ty,tz);
1090 /* #define INNERFLOPS INNERFLOPS+3 */
1092 fjx{J} = _mm256_add_ps(fjx{J},tx);
1093 fjy{J} = _mm256_add_ps(fjy{J},ty);
1094 fjz{J} = _mm256_add_ps(fjz{J},tz);
1095 /* #define INNERFLOPS INNERFLOPS+3 */
1100 /* ## Note special check for TIP4P-TIP4P. Since we are cutting of all hydrogen interactions we also cut the LJ-only O-O interaction */
1101 /* #if 'exactcutoff' in INTERACTION_FLAGS[I][J] or (GEOMETRY_I=='Water4' and GEOMETRY_J=='Water4' and 'exactcutoff' in INTERACTION_FLAGS[1][1]) */
1102 /* #if 0 ## This and next two lines is a hack to maintain indentation in template file */
1107 /* ## End of check for the interaction being outside the cutoff */
1110 /* ## End of loop over i-j interaction pairs */
1112 /* #if GEOMETRY_I != 'Particle' */
1113 /* #if ROUND == 'Loop' */
1114 fjptrA = f+j_coord_offsetA;
1115 fjptrB = f+j_coord_offsetB;
1116 fjptrC = f+j_coord_offsetC;
1117 fjptrD = f+j_coord_offsetD;
1118 fjptrE = f+j_coord_offsetE;
1119 fjptrF = f+j_coord_offsetF;
1120 fjptrG = f+j_coord_offsetG;
1121 fjptrH = f+j_coord_offsetH;
1123 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
1124 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
1125 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
1126 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
1127 fjptrE = (jnrlistE>=0) ? f+j_coord_offsetE : scratch;
1128 fjptrF = (jnrlistF>=0) ? f+j_coord_offsetF : scratch;
1129 fjptrG = (jnrlistG>=0) ? f+j_coord_offsetG : scratch;
1130 fjptrH = (jnrlistH>=0) ? f+j_coord_offsetH : scratch;
1134 /* #if 'Water' in GEOMETRY_I and GEOMETRY_J == 'Particle' */
1135 gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,fjx0,fjy0,fjz0);
1136 /* #define INNERFLOPS INNERFLOPS+3 */
1137 /* #elif GEOMETRY_J == 'Water3' */
1138 gmx_mm256_decrement_3rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,
1139 fjx0,fjy0,fjz0,fjx1,fjy1,fjz1,fjx2,fjy2,fjz2);
1140 /* #define INNERFLOPS INNERFLOPS+9 */
1141 /* #elif GEOMETRY_J == 'Water4' */
1142 /* #if 0 in PARTICLES_J */
1143 gmx_mm256_decrement_4rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,
1144 fjx0,fjy0,fjz0,fjx1,fjy1,fjz1,
1145 fjx2,fjy2,fjz2,fjx3,fjy3,fjz3);
1146 /* #define INNERFLOPS INNERFLOPS+12 */
1148 gmx_mm256_decrement_3rvec_8ptr_swizzle_ps(fjptrA+DIM,fjptrB+DIM,fjptrC+DIM,fjptrD+DIM,
1149 fjptrE+DIM,fjptrF+DIM,fjptrG+DIM,fjptrH+DIM,
1150 fjx1,fjy1,fjz1,fjx2,fjy2,fjz2,fjx3,fjy3,fjz3);
1151 /* #define INNERFLOPS INNERFLOPS+9 */
1155 /* Inner loop uses {INNERFLOPS} flops */
1160 /* End of innermost loop */
1162 /* #if 'Force' in KERNEL_VF */
1163 /* #if GEOMETRY_I == 'Particle' */
1164 gmx_mm256_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
1165 f+i_coord_offset,fshift+i_shift_offset);
1166 /* #define OUTERFLOPS OUTERFLOPS+6 */
1167 /* #elif GEOMETRY_I == 'Water3' */
1168 gmx_mm256_update_iforce_3atom_swizzle_ps(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,
1169 f+i_coord_offset,fshift+i_shift_offset);
1170 /* #define OUTERFLOPS OUTERFLOPS+18 */
1171 /* #elif GEOMETRY_I == 'Water4' */
1172 /* #if 0 in PARTICLES_I */
1173 gmx_mm256_update_iforce_4atom_swizzle_ps(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
1174 f+i_coord_offset,fshift+i_shift_offset);
1175 /* #define OUTERFLOPS OUTERFLOPS+24 */
1177 gmx_mm256_update_iforce_3atom_swizzle_ps(fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
1178 f+i_coord_offset+DIM,fshift+i_shift_offset);
1179 /* #define OUTERFLOPS OUTERFLOPS+18 */
1184 /* #if 'Potential' in KERNEL_VF */
1186 /* Update potential energies */
1187 /* #if KERNEL_ELEC != 'None' */
1188 gmx_mm256_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
1189 /* #define OUTERFLOPS OUTERFLOPS+1 */
1191 /* #if 'GeneralizedBorn' in KERNEL_ELEC */
1192 gmx_mm256_update_1pot_ps(vgbsum,kernel_data->energygrp_polarization+ggid);
1193 /* #define OUTERFLOPS OUTERFLOPS+1 */
1195 /* #if KERNEL_VDW != 'None' */
1196 gmx_mm256_update_1pot_ps(vvdwsum,kernel_data->energygrp_vdw+ggid);
1197 /* #define OUTERFLOPS OUTERFLOPS+1 */
1200 /* #if 'GeneralizedBorn' in KERNEL_ELEC and 'Force' in KERNEL_VF */
1201 dvdasum = _mm256_mul_ps(dvdasum, _mm256_mul_ps(isai{I},isai{I}));
1202 gmx_mm256_update_1pot_ps(dvdasum,dvda+inr);
1205 /* Increment number of inner iterations */
1206 inneriter += j_index_end - j_index_start;
1208 /* Outer loop uses {OUTERFLOPS} flops */
1211 /* Increment number of outer iterations */
1214 /* Update outer/inner flops */
1215 /* ## NB: This is not important, it just affects the flopcount. However, since our preprocessor is */
1216 /* ## primitive and replaces aggressively even in strings inside these directives, we need to */
1217 /* ## assemble the main part of the name (containing KERNEL/ELEC/VDW) directly in the source. */
1218 /* #if GEOMETRY_I == 'Water3' */
1219 /* #define ISUFFIX '_W3' */
1220 /* #elif GEOMETRY_I == 'Water4' */
1221 /* #define ISUFFIX '_W4' */
1223 /* #define ISUFFIX '' */
1225 /* #if GEOMETRY_J == 'Water3' */
1226 /* #define JSUFFIX 'W3' */
1227 /* #elif GEOMETRY_J == 'Water4' */
1228 /* #define JSUFFIX 'W4' */
1230 /* #define JSUFFIX '' */
1232 /* #if 'PotentialAndForce' in KERNEL_VF */
1233 /* #define VFSUFFIX '_VF' */
1234 /* #elif 'Potential' in KERNEL_VF */
1235 /* #define VFSUFFIX '_V' */
1237 /* #define VFSUFFIX '_F' */
1240 /* #if KERNEL_ELEC != 'None' and KERNEL_VDW != 'None' */
1241 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW{ISUFFIX}{JSUFFIX}{VFSUFFIX},outeriter*{OUTERFLOPS} + inneriter*{INNERFLOPS});
1242 /* #elif KERNEL_ELEC != 'None' */
1243 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC{ISUFFIX}{JSUFFIX}{VFSUFFIX},outeriter*{OUTERFLOPS} + inneriter*{INNERFLOPS});
1245 inc_nrnb(nrnb,eNR_NBKERNEL_VDW{ISUFFIX}{JSUFFIX}{VFSUFFIX},outeriter*{OUTERFLOPS} + inneriter*{INNERFLOPS});