2 #error This file must be processed with the Gromacs pre-preprocessor
4 /* #if INCLUDE_HEADER */
11 #include "../nb_kernel.h"
12 #include "types/simple.h"
16 #include "gmx_math_x86_sse2_single.h"
17 #include "kernelutil_x86_sse2_single.h"
20 /* ## List of variables set by the generating script: */
22 /* ## Setttings that apply to the entire kernel: */
23 /* ## KERNEL_ELEC: String, choice for electrostatic interactions */
24 /* ## KERNEL_VDW: String, choice for van der Waals interactions */
25 /* ## KERNEL_NAME: String, name of this kernel */
26 /* ## KERNEL_VF: String telling if we calculate potential, force, or both */
27 /* ## GEOMETRY_I/GEOMETRY_J: String, name of each geometry, e.g. 'Water3' or '1Particle' */
29 /* ## Setttings that apply to particles in the outer (I) or inner (J) loops: */
30 /* ## PARTICLES_I[]/ Arrays with lists of i/j particles to use in kernel. It is */
31 /* ## PARTICLES_J[]: just [0] for particle geometry, but can be longer for water */
32 /* ## PARTICLES_ELEC_I[]/ Arrays with lists of i/j particle that have electrostatics */
33 /* ## PARTICLES_ELEC_J[]: interactions that should be calculated in this kernel. */
34 /* ## PARTICLES_VDW_I[]/ Arrays with the list of i/j particle that have VdW */
35 /* ## PARTICLES_VDW_J[]: interactions that should be calculated in this kernel. */
37 /* ## Setttings for pairs of interactions (e.g. 2nd i particle against 1st j particle) */
38 /* ## PAIRS_IJ[]: Array with (i,j) tuples of pairs for which interactions */
39 /* ## should be calculated in this kernel. Zero-charge particles */
40 /* ## do not have interactions with particles without vdw, and */
41 /* ## Vdw-only interactions are not evaluated in a no-vdw-kernel. */
42 /* ## INTERACTION_FLAGS[][]: 2D matrix, dimension e.g. 3*3 for water-water interactions. */
43 /* ## For each i-j pair, the element [I][J] is a list of strings */
44 /* ## defining properties/flags of this interaction. Examples */
45 /* ## include 'electrostatics'/'vdw' if that type of interaction */
46 /* ## should be evaluated, 'rsq'/'rinv'/'rinvsq' if those values */
47 /* ## are needed, and 'exactcutoff' or 'shift','switch' to */
48 /* ## decide if the force/potential should be modified. This way */
49 /* ## we only calculate values absolutely needed for each case. */
51 /* ## Calculate the size and offset for (merged/interleaved) table data */
54 * Gromacs nonbonded kernel: {KERNEL_NAME}
55 * Electrostatics interaction: {KERNEL_ELEC}
56 * VdW interaction: {KERNEL_VDW}
57 * Geometry: {GEOMETRY_I}-{GEOMETRY_J}
58 * Calculate force/pot: {KERNEL_VF}
62 (t_nblist * gmx_restrict nlist,
63 rvec * gmx_restrict xx,
64 rvec * gmx_restrict ff,
65 t_forcerec * gmx_restrict fr,
66 t_mdatoms * gmx_restrict mdatoms,
67 nb_kernel_data_t * gmx_restrict kernel_data,
68 t_nrnb * gmx_restrict nrnb)
70 /* ## Not all variables are used for all kernels, but any optimizing compiler fixes that, */
71 /* ## so there is no point in going to extremes to exclude variables that are not needed. */
72 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
73 * just 0 for non-waters.
74 * Suffixes A,B,C,D refer to j loop unrolling done with SSE, e.g. for the four different
75 * jnr indices corresponding to data put in the four positions in the SIMD register.
77 int i_shift_offset,i_coord_offset,outeriter,inneriter;
78 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
79 int jnrA,jnrB,jnrC,jnrD;
80 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
81 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
82 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
84 real *shiftvec,*fshift,*x,*f;
85 real *fjptrA,*fjptrB,*fjptrC,*fjptrD;
87 __m128 tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
88 /* #for I in PARTICLES_I */
90 __m128 ix{I},iy{I},iz{I},fix{I},fiy{I},fiz{I},iq{I},isai{I};
92 /* #for J in PARTICLES_J */
93 int vdwjidx{J}A,vdwjidx{J}B,vdwjidx{J}C,vdwjidx{J}D;
94 __m128 jx{J},jy{J},jz{J},fjx{J},fjy{J},fjz{J},jq{J},isaj{J};
96 /* #for I,J in PAIRS_IJ */
97 __m128 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};
99 /* #if KERNEL_ELEC != 'None' */
100 __m128 velec,felec,velecsum,facel,crf,krf,krf2;
103 /* #if 'GeneralizedBorn' in KERNEL_ELEC */
105 __m128 vgb,fgb,vgbsum,dvdasum,gbscale,gbtabscale,isaprod,gbqqfactor,gbinvepsdiff,gbeps,dvdatmp;
106 __m128 minushalf = _mm_set1_ps(-0.5);
107 real *invsqrta,*dvda,*gbtab;
109 /* #if KERNEL_VDW != 'None' */
111 __m128 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
114 __m128 one_sixth = _mm_set1_ps(1.0/6.0);
115 __m128 one_twelfth = _mm_set1_ps(1.0/12.0);
117 /* #if 'Table' in KERNEL_ELEC or 'GeneralizedBorn' in KERNEL_ELEC or 'Table' in KERNEL_VDW */
119 __m128i ifour = _mm_set1_epi32(4);
120 __m128 rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF;
123 /* #if 'Ewald' in KERNEL_ELEC */
125 __m128 ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
128 /* #if 'PotentialSwitch' in [KERNEL_MOD_ELEC,KERNEL_MOD_VDW] */
129 __m128 rswitch,swV3,swV4,swV5,swF2,swF3,swF4,d,d2,sw,dsw;
130 real rswitch_scalar,d_scalar;
132 __m128 dummy_mask,cutoff_mask;
133 __m128 signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
134 __m128 one = _mm_set1_ps(1.0);
135 __m128 two = _mm_set1_ps(2.0);
141 jindex = nlist->jindex;
143 shiftidx = nlist->shift;
145 shiftvec = fr->shift_vec[0];
146 fshift = fr->fshift[0];
147 /* #if KERNEL_ELEC != 'None' */
148 facel = _mm_set1_ps(fr->epsfac);
149 charge = mdatoms->chargeA;
150 /* #if 'ReactionField' in KERNEL_ELEC */
151 krf = _mm_set1_ps(fr->ic->k_rf);
152 krf2 = _mm_set1_ps(fr->ic->k_rf*2.0);
153 crf = _mm_set1_ps(fr->ic->c_rf);
156 /* #if KERNEL_VDW != 'None' */
157 nvdwtype = fr->ntype;
159 vdwtype = mdatoms->typeA;
162 /* #if 'Table' in KERNEL_ELEC and 'Table' in KERNEL_VDW */
163 vftab = kernel_data->table_elec_vdw->data;
164 vftabscale = _mm_set1_ps(kernel_data->table_elec_vdw->scale);
165 /* #elif 'Table' in KERNEL_ELEC */
166 vftab = kernel_data->table_elec->data;
167 vftabscale = _mm_set1_ps(kernel_data->table_elec->scale);
168 /* #elif 'Table' in KERNEL_VDW */
169 vftab = kernel_data->table_vdw->data;
170 vftabscale = _mm_set1_ps(kernel_data->table_vdw->scale);
173 /* #if 'Ewald' in KERNEL_ELEC */
174 sh_ewald = _mm_set1_ps(fr->ic->sh_ewald);
175 /* #if KERNEL_VF=='Force' and KERNEL_MOD_ELEC!='PotentialSwitch' */
176 ewtab = fr->ic->tabq_coul_F;
177 ewtabscale = _mm_set1_ps(fr->ic->tabq_scale);
178 ewtabhalfspace = _mm_set1_ps(0.5/fr->ic->tabq_scale);
180 ewtab = fr->ic->tabq_coul_FDV0;
181 ewtabscale = _mm_set1_ps(fr->ic->tabq_scale);
182 ewtabhalfspace = _mm_set1_ps(0.5/fr->ic->tabq_scale);
186 /* #if KERNEL_ELEC=='GeneralizedBorn' */
187 invsqrta = fr->invsqrta;
189 gbtabscale = _mm_set1_ps(fr->gbtab.scale);
190 gbtab = fr->gbtab.data;
191 gbinvepsdiff = _mm_set1_ps((1.0/fr->epsilon_r) - (1.0/fr->gb_epsilon_solvent));
194 /* #if 'Water' in GEOMETRY_I */
195 /* Setup water-specific parameters */
196 inr = nlist->iinr[0];
197 /* #for I in PARTICLES_ELEC_I */
198 iq{I} = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+{I}]));
200 /* #for I in PARTICLES_VDW_I */
201 vdwioffset{I} = 2*nvdwtype*vdwtype[inr+{I}];
205 /* #if 'Water' in GEOMETRY_J */
206 /* #for J in PARTICLES_ELEC_J */
207 jq{J} = _mm_set1_ps(charge[inr+{J}]);
209 /* #for J in PARTICLES_VDW_J */
210 vdwjidx{J}A = 2*vdwtype[inr+{J}];
212 /* #for I,J in PAIRS_IJ */
213 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] */
214 qq{I}{J} = _mm_mul_ps(iq{I},jq{J});
216 /* #if 'vdw' in INTERACTION_FLAGS[I][J] */
217 c6_{I}{J} = _mm_set1_ps(vdwparam[vdwioffset{I}+vdwjidx{J}A]);
218 c12_{I}{J} = _mm_set1_ps(vdwparam[vdwioffset{I}+vdwjidx{J}A+1]);
223 /* #if KERNEL_MOD_ELEC!='None' or KERNEL_MOD_VDW!='None' */
224 /* #if KERNEL_ELEC!='None' */
225 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
226 rcutoff_scalar = fr->rcoulomb;
228 rcutoff_scalar = fr->rvdw;
230 rcutoff = _mm_set1_ps(rcutoff_scalar);
231 rcutoff2 = _mm_mul_ps(rcutoff,rcutoff);
234 /* #if KERNEL_MOD_VDW=='PotentialShift' */
235 sh_vdw_invrcut6 = _mm_set1_ps(fr->ic->sh_invrc6);
236 rvdw = _mm_set1_ps(fr->rvdw);
239 /* #if 'PotentialSwitch' in [KERNEL_MOD_ELEC,KERNEL_MOD_VDW] */
240 /* #if KERNEL_MOD_ELEC=='PotentialSwitch' */
241 rswitch_scalar = fr->rcoulomb_switch;
242 rswitch = _mm_set1_ps(rswitch_scalar);
244 rswitch_scalar = fr->rvdw_switch;
245 rswitch = _mm_set1_ps(rswitch_scalar);
247 /* Setup switch parameters */
248 d_scalar = rcutoff_scalar-rswitch_scalar;
249 d = _mm_set1_ps(d_scalar);
250 swV3 = _mm_set1_ps(-10.0/(d_scalar*d_scalar*d_scalar));
251 swV4 = _mm_set1_ps( 15.0/(d_scalar*d_scalar*d_scalar*d_scalar));
252 swV5 = _mm_set1_ps( -6.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
253 /* #if 'Force' in KERNEL_VF */
254 swF2 = _mm_set1_ps(-30.0/(d_scalar*d_scalar*d_scalar));
255 swF3 = _mm_set1_ps( 60.0/(d_scalar*d_scalar*d_scalar*d_scalar));
256 swF4 = _mm_set1_ps(-30.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
260 /* Avoid stupid compiler warnings */
261 jnrA = jnrB = jnrC = jnrD = 0;
267 /* ## Keep track of the floating point operations we issue for reporting! */
268 /* #define OUTERFLOPS 0 */
272 for(iidx=0;iidx<4*DIM;iidx++)
277 /* Start outer loop over neighborlists */
278 for(iidx=0; iidx<nri; iidx++)
280 /* Load shift vector for this list */
281 i_shift_offset = DIM*shiftidx[iidx];
283 /* Load limits for loop over neighbors */
284 j_index_start = jindex[iidx];
285 j_index_end = jindex[iidx+1];
287 /* Get outer coordinate index */
289 i_coord_offset = DIM*inr;
291 /* Load i particle coords and add shift vector */
292 /* #if GEOMETRY_I == 'Particle' */
293 gmx_mm_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
294 /* #elif GEOMETRY_I == 'Water3' */
295 gmx_mm_load_shift_and_3rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,
296 &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2);
297 /* #elif GEOMETRY_I == 'Water4' */
298 /* #if 0 in PARTICLES_I */
299 gmx_mm_load_shift_and_4rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,
300 &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
302 gmx_mm_load_shift_and_3rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset+DIM,
303 &ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
307 /* #if 'Force' in KERNEL_VF */
308 /* #for I in PARTICLES_I */
309 fix{I} = _mm_setzero_ps();
310 fiy{I} = _mm_setzero_ps();
311 fiz{I} = _mm_setzero_ps();
315 /* ## For water we already preloaded parameters at the start of the kernel */
316 /* #if not 'Water' in GEOMETRY_I */
317 /* Load parameters for i particles */
318 /* #for I in PARTICLES_ELEC_I */
319 iq{I} = _mm_mul_ps(facel,_mm_load1_ps(charge+inr+{I}));
320 /* #define OUTERFLOPS OUTERFLOPS+1 */
321 /* #if KERNEL_ELEC=='GeneralizedBorn' */
322 isai{I} = _mm_load1_ps(invsqrta+inr+{I});
325 /* #for I in PARTICLES_VDW_I */
326 vdwioffset{I} = 2*nvdwtype*vdwtype[inr+{I}];
330 /* #if 'Potential' in KERNEL_VF */
331 /* Reset potential sums */
332 /* #if KERNEL_ELEC != 'None' */
333 velecsum = _mm_setzero_ps();
335 /* #if 'GeneralizedBorn' in KERNEL_ELEC */
336 vgbsum = _mm_setzero_ps();
338 /* #if KERNEL_VDW != 'None' */
339 vvdwsum = _mm_setzero_ps();
342 /* #if 'GeneralizedBorn' in KERNEL_ELEC and 'Force' in KERNEL_VF */
343 dvdasum = _mm_setzero_ps();
346 /* #for ROUND in ['Loop','Epilogue'] */
348 /* #if ROUND =='Loop' */
349 /* Start inner kernel loop */
350 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
352 /* ## First round is normal loop (next statement resets indentation) */
359 /* ## Second round is epilogue */
361 /* #define INNERFLOPS 0 */
363 /* Get j neighbor index, and coordinate index */
364 /* #if ROUND =='Loop' */
370 jnrlistA = jjnr[jidx];
371 jnrlistB = jjnr[jidx+1];
372 jnrlistC = jjnr[jidx+2];
373 jnrlistD = jjnr[jidx+3];
374 /* Sign of each element will be negative for non-real atoms.
375 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
376 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
378 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
379 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
380 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
381 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
382 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
384 j_coord_offsetA = DIM*jnrA;
385 j_coord_offsetB = DIM*jnrB;
386 j_coord_offsetC = DIM*jnrC;
387 j_coord_offsetD = DIM*jnrD;
389 /* load j atom coordinates */
390 /* #if GEOMETRY_J == 'Particle' */
391 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
392 x+j_coord_offsetC,x+j_coord_offsetD,
394 /* #elif GEOMETRY_J == 'Water3' */
395 gmx_mm_load_3rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
396 x+j_coord_offsetC,x+j_coord_offsetD,
397 &jx0,&jy0,&jz0,&jx1,&jy1,&jz1,&jx2,&jy2,&jz2);
398 /* #elif GEOMETRY_J == 'Water4' */
399 /* #if 0 in PARTICLES_J */
400 gmx_mm_load_4rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
401 x+j_coord_offsetC,x+j_coord_offsetD,
402 &jx0,&jy0,&jz0,&jx1,&jy1,&jz1,&jx2,
403 &jy2,&jz2,&jx3,&jy3,&jz3);
405 gmx_mm_load_3rvec_4ptr_swizzle_ps(x+j_coord_offsetA+DIM,x+j_coord_offsetB+DIM,
406 x+j_coord_offsetC+DIM,x+j_coord_offsetD+DIM,
407 &jx1,&jy1,&jz1,&jx2,&jy2,&jz2,&jx3,&jy3,&jz3);
411 /* Calculate displacement vector */
412 /* #for I,J in PAIRS_IJ */
413 dx{I}{J} = _mm_sub_ps(ix{I},jx{J});
414 dy{I}{J} = _mm_sub_ps(iy{I},jy{J});
415 dz{I}{J} = _mm_sub_ps(iz{I},jz{J});
416 /* #define INNERFLOPS INNERFLOPS+3 */
419 /* Calculate squared distance and things based on it */
420 /* #for I,J in PAIRS_IJ */
421 rsq{I}{J} = gmx_mm_calc_rsq_ps(dx{I}{J},dy{I}{J},dz{I}{J});
422 /* #define INNERFLOPS INNERFLOPS+5 */
425 /* #for I,J in PAIRS_IJ */
426 /* #if 'rinv' in INTERACTION_FLAGS[I][J] */
427 rinv{I}{J} = gmx_mm_invsqrt_ps(rsq{I}{J});
428 /* #define INNERFLOPS INNERFLOPS+5 */
432 /* #for I,J in PAIRS_IJ */
433 /* #if 'rinvsq' in INTERACTION_FLAGS[I][J] */
434 /* # if 'rinv' not in INTERACTION_FLAGS[I][J] */
435 rinvsq{I}{J} = gmx_mm_inv_ps(rsq{I}{J});
436 /* #define INNERFLOPS INNERFLOPS+4 */
438 rinvsq{I}{J} = _mm_mul_ps(rinv{I}{J},rinv{I}{J});
439 /* #define INNERFLOPS INNERFLOPS+1 */
444 /* #if not 'Water' in GEOMETRY_J */
445 /* Load parameters for j particles */
446 /* #for J in PARTICLES_ELEC_J */
447 jq{J} = gmx_mm_load_4real_swizzle_ps(charge+jnrA+{J},charge+jnrB+{J},
448 charge+jnrC+{J},charge+jnrD+{J});
449 /* #if KERNEL_ELEC=='GeneralizedBorn' */
450 isaj{J} = gmx_mm_load_4real_swizzle_ps(invsqrta+jnrA+{J},invsqrta+jnrB+{J},
451 invsqrta+jnrC+{J},invsqrta+jnrD+{J});
454 /* #for J in PARTICLES_VDW_J */
455 vdwjidx{J}A = 2*vdwtype[jnrA+{J}];
456 vdwjidx{J}B = 2*vdwtype[jnrB+{J}];
457 vdwjidx{J}C = 2*vdwtype[jnrC+{J}];
458 vdwjidx{J}D = 2*vdwtype[jnrD+{J}];
462 /* #if 'Force' in KERNEL_VF and not 'Particle' in GEOMETRY_I */
463 /* #for J in PARTICLES_J */
464 fjx{J} = _mm_setzero_ps();
465 fjy{J} = _mm_setzero_ps();
466 fjz{J} = _mm_setzero_ps();
470 /* #for I,J in PAIRS_IJ */
472 /**************************
473 * CALCULATE INTERACTIONS *
474 **************************/
476 /* ## Note special check for TIP4P-TIP4P. Since we are cutting of all hydrogen interactions we also cut the LJ-only O-O interaction */
477 /* #if 'exactcutoff' in INTERACTION_FLAGS[I][J] or (GEOMETRY_I=='Water4' and GEOMETRY_J=='Water4' and 'exactcutoff' in INTERACTION_FLAGS[1][1]) */
478 /* ## We always calculate rinv/rinvsq above to enable pipelineing in compilers (performance tested on x86) */
479 if (gmx_mm_any_lt(rsq{I}{J},rcutoff2))
481 /* #if 0 ## this and the next two lines is a hack to maintain auto-indentation in template file */
484 /* #define INNERFLOPS INNERFLOPS+1 */
487 /* #if 'r' in INTERACTION_FLAGS[I][J] */
488 r{I}{J} = _mm_mul_ps(rsq{I}{J},rinv{I}{J});
489 /* #if ROUND == 'Epilogue' */
490 r{I}{J} = _mm_andnot_ps(dummy_mask,r{I}{J});
491 /* #define INNERFLOPS INNERFLOPS+1 */
493 /* #define INNERFLOPS INNERFLOPS+1 */
496 /* ## For water geometries we already loaded parameters at the start of the kernel */
497 /* #if not 'Water' in GEOMETRY_J */
498 /* Compute parameters for interactions between i and j atoms */
499 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] */
500 qq{I}{J} = _mm_mul_ps(iq{I},jq{J});
501 /* #define INNERFLOPS INNERFLOPS+1 */
503 /* #if 'vdw' in INTERACTION_FLAGS[I][J] */
504 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset{I}+vdwjidx{J}A,
505 vdwparam+vdwioffset{I}+vdwjidx{J}B,
506 vdwparam+vdwioffset{I}+vdwjidx{J}C,
507 vdwparam+vdwioffset{I}+vdwjidx{J}D,
508 &c6_{I}{J},&c12_{I}{J});
512 /* #if 'table' in INTERACTION_FLAGS[I][J] */
513 /* Calculate table index by multiplying r with table scale and truncate to integer */
514 rt = _mm_mul_ps(r{I}{J},vftabscale);
515 vfitab = _mm_cvttps_epi32(rt);
516 vfeps = _mm_sub_ps(rt,_mm_cvtepi32_ps(vfitab));
517 /* #define INNERFLOPS INNERFLOPS+4 */
518 /* #if 'Table' in KERNEL_ELEC and 'Table' in KERNEL_VDW */
519 /* ## 3 tables, 4 bytes per point: multiply index by 12 */
520 vfitab = _mm_slli_epi32(_mm_add_epi32(vfitab,_mm_slli_epi32(vfitab,1)),2);
521 /* #elif 'Table' in KERNEL_ELEC */
522 /* ## 1 table, 4 bytes per point: multiply index by 4 */
523 vfitab = _mm_slli_epi32(vfitab,2);
524 /* #elif 'Table' in KERNEL_VDW */
525 /* ## 2 tables, 4 bytes per point: multiply index by 8 */
526 vfitab = _mm_slli_epi32(vfitab,3);
530 /* ## ELECTROSTATIC INTERACTIONS */
531 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] */
533 /* #if KERNEL_ELEC=='Coulomb' */
535 /* COULOMB ELECTROSTATICS */
536 velec = _mm_mul_ps(qq{I}{J},rinv{I}{J});
537 /* #define INNERFLOPS INNERFLOPS+1 */
538 /* #if 'Force' in KERNEL_VF */
539 felec = _mm_mul_ps(velec,rinvsq{I}{J});
540 /* #define INNERFLOPS INNERFLOPS+2 */
543 /* #elif KERNEL_ELEC=='ReactionField' */
545 /* REACTION-FIELD ELECTROSTATICS */
546 /* #if 'Potential' in KERNEL_VF */
547 velec = _mm_mul_ps(qq{I}{J},_mm_sub_ps(_mm_add_ps(rinv{I}{J},_mm_mul_ps(krf,rsq{I}{J})),crf));
548 /* #define INNERFLOPS INNERFLOPS+4 */
550 /* #if 'Force' in KERNEL_VF */
551 felec = _mm_mul_ps(qq{I}{J},_mm_sub_ps(_mm_mul_ps(rinv{I}{J},rinvsq{I}{J}),krf2));
552 /* #define INNERFLOPS INNERFLOPS+3 */
555 /* #elif KERNEL_ELEC=='GeneralizedBorn' */
557 /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
558 isaprod = _mm_mul_ps(isai{I},isaj{J});
559 gbqqfactor = _mm_xor_ps(signbit,_mm_mul_ps(qq{I}{J},_mm_mul_ps(isaprod,gbinvepsdiff)));
560 gbscale = _mm_mul_ps(isaprod,gbtabscale);
561 /* #define INNERFLOPS INNERFLOPS+5 */
563 /* Calculate generalized born table index - this is a separate table from the normal one,
564 * but we use the same procedure by multiplying r with scale and truncating to integer.
566 rt = _mm_mul_ps(r{I}{J},gbscale);
567 gbitab = _mm_cvttps_epi32(rt);
568 gbeps = _mm_sub_ps(rt,_mm_cvtepi32_ps(gbitab));
569 gbitab = _mm_slli_epi32(gbitab,2);
571 Y = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,0) );
572 F = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,1) );
573 G = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,2) );
574 H = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,3) );
575 _MM_TRANSPOSE4_PS(Y,F,G,H);
576 Heps = _mm_mul_ps(gbeps,H);
577 Fp = _mm_add_ps(F,_mm_mul_ps(gbeps,_mm_add_ps(G,Heps)));
578 VV = _mm_add_ps(Y,_mm_mul_ps(gbeps,Fp));
579 vgb = _mm_mul_ps(gbqqfactor,VV);
580 /* #define INNERFLOPS INNERFLOPS+10 */
582 /* #if 'Force' in KERNEL_VF */
583 FF = _mm_add_ps(Fp,_mm_mul_ps(gbeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
584 fgb = _mm_mul_ps(gbqqfactor,_mm_mul_ps(FF,gbscale));
585 dvdatmp = _mm_mul_ps(minushalf,_mm_add_ps(vgb,_mm_mul_ps(fgb,r{I}{J})));
586 dvdasum = _mm_add_ps(dvdasum,dvdatmp);
587 /* #if ROUND == 'Loop' */
593 /* 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. */
594 fjptrA = (jnrlistA>=0) ? dvda+jnrA : scratch;
595 fjptrB = (jnrlistB>=0) ? dvda+jnrB : scratch;
596 fjptrC = (jnrlistC>=0) ? dvda+jnrC : scratch;
597 fjptrD = (jnrlistD>=0) ? dvda+jnrD : scratch;
599 gmx_mm_increment_4real_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,_mm_mul_ps(dvdatmp,_mm_mul_ps(isaj{J},isaj{J})));
600 /* #define INNERFLOPS INNERFLOPS+13 */
602 velec = _mm_mul_ps(qq{I}{J},rinv{I}{J});
603 /* #define INNERFLOPS INNERFLOPS+1 */
604 /* #if 'Force' in KERNEL_VF */
605 felec = _mm_mul_ps(_mm_sub_ps(_mm_mul_ps(velec,rinv{I}{J}),fgb),rinv{I}{J});
606 /* #define INNERFLOPS INNERFLOPS+3 */
609 /* #elif KERNEL_ELEC=='Ewald' */
610 /* EWALD ELECTROSTATICS */
612 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
613 ewrt = _mm_mul_ps(r{I}{J},ewtabscale);
614 ewitab = _mm_cvttps_epi32(ewrt);
615 eweps = _mm_sub_ps(ewrt,_mm_cvtepi32_ps(ewitab));
616 /* #define INNERFLOPS INNERFLOPS+4 */
617 /* #if 'Potential' in KERNEL_VF or KERNEL_MOD_ELEC=='PotentialSwitch' */
618 ewitab = _mm_slli_epi32(ewitab,2);
619 ewtabF = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,0) );
620 ewtabD = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,1) );
621 ewtabV = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,2) );
622 ewtabFn = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,3) );
623 _MM_TRANSPOSE4_PS(ewtabF,ewtabD,ewtabV,ewtabFn);
624 felec = _mm_add_ps(ewtabF,_mm_mul_ps(eweps,ewtabD));
625 /* #define INNERFLOPS INNERFLOPS+2 */
626 /* #if KERNEL_MOD_ELEC=='PotentialShift' */
627 velec = _mm_sub_ps(ewtabV,_mm_mul_ps(_mm_mul_ps(ewtabhalfspace,eweps),_mm_add_ps(ewtabF,felec)));
628 velec = _mm_mul_ps(qq{I}{J},_mm_sub_ps(_mm_sub_ps(rinv{I}{J},sh_ewald),velec));
629 /* #define INNERFLOPS INNERFLOPS+7 */
631 velec = _mm_sub_ps(ewtabV,_mm_mul_ps(_mm_mul_ps(ewtabhalfspace,eweps),_mm_add_ps(ewtabF,felec)));
632 velec = _mm_mul_ps(qq{I}{J},_mm_sub_ps(rinv{I}{J},velec));
633 /* #define INNERFLOPS INNERFLOPS+6 */
635 /* #if 'Force' in KERNEL_VF */
636 felec = _mm_mul_ps(_mm_mul_ps(qq{I}{J},rinv{I}{J}),_mm_sub_ps(rinvsq{I}{J},felec));
637 /* #define INNERFLOPS INNERFLOPS+3 */
639 /* #elif KERNEL_VF=='Force' */
640 gmx_mm_load_4pair_swizzle_ps(ewtab+gmx_mm_extract_epi32(ewitab,0),ewtab+gmx_mm_extract_epi32(ewitab,1),
641 ewtab+gmx_mm_extract_epi32(ewitab,2),ewtab+gmx_mm_extract_epi32(ewitab,3),
643 felec = _mm_add_ps(_mm_mul_ps( _mm_sub_ps(one,eweps),ewtabF),_mm_mul_ps(eweps,ewtabFn));
644 felec = _mm_mul_ps(_mm_mul_ps(qq{I}{J},rinv{I}{J}),_mm_sub_ps(rinvsq{I}{J},felec));
645 /* #define INNERFLOPS INNERFLOPS+7 */
648 /* #elif KERNEL_ELEC=='CubicSplineTable' */
650 /* CUBIC SPLINE TABLE ELECTROSTATICS */
651 Y = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,0) );
652 F = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,1) );
653 G = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,2) );
654 H = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,3) );
655 _MM_TRANSPOSE4_PS(Y,F,G,H);
656 Heps = _mm_mul_ps(vfeps,H);
657 Fp = _mm_add_ps(F,_mm_mul_ps(vfeps,_mm_add_ps(G,Heps)));
658 /* #define INNERFLOPS INNERFLOPS+4 */
659 /* #if 'Potential' in KERNEL_VF */
660 VV = _mm_add_ps(Y,_mm_mul_ps(vfeps,Fp));
661 velec = _mm_mul_ps(qq{I}{J},VV);
662 /* #define INNERFLOPS INNERFLOPS+3 */
664 /* #if 'Force' in KERNEL_VF */
665 FF = _mm_add_ps(Fp,_mm_mul_ps(vfeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
666 felec = _mm_xor_ps(signbit,_mm_mul_ps(_mm_mul_ps(qq{I}{J},FF),_mm_mul_ps(vftabscale,rinv{I}{J})));
667 /* #define INNERFLOPS INNERFLOPS+7 */
670 /* ## End of check for electrostatics interaction forms */
672 /* ## END OF ELECTROSTATIC INTERACTION CHECK FOR PAIR I-J */
674 /* #if 'vdw' in INTERACTION_FLAGS[I][J] */
676 /* #if KERNEL_VDW=='LennardJones' */
678 /* LENNARD-JONES DISPERSION/REPULSION */
680 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq{I}{J},rinvsq{I}{J}),rinvsq{I}{J});
681 /* #define INNERFLOPS INNERFLOPS+2 */
682 /* #if 'Potential' in KERNEL_VF or KERNEL_MOD_VDW=='PotentialSwitch' */
683 vvdw6 = _mm_mul_ps(c6_{I}{J},rinvsix);
684 vvdw12 = _mm_mul_ps(c12_{I}{J},_mm_mul_ps(rinvsix,rinvsix));
685 /* #define INNERFLOPS INNERFLOPS+3 */
686 /* #if KERNEL_MOD_VDW=='PotentialShift' */
687 vvdw = _mm_sub_ps(_mm_mul_ps( _mm_sub_ps(vvdw12 , _mm_mul_ps(c12_{I}{J},_mm_mul_ps(sh_vdw_invrcut6,sh_vdw_invrcut6))), one_twelfth) ,
688 _mm_mul_ps( _mm_sub_ps(vvdw6,_mm_mul_ps(c6_{I}{J},sh_vdw_invrcut6)),one_sixth));
689 /* #define INNERFLOPS INNERFLOPS+8 */
691 vvdw = _mm_sub_ps( _mm_mul_ps(vvdw12,one_twelfth) , _mm_mul_ps(vvdw6,one_sixth) );
692 /* #define INNERFLOPS INNERFLOPS+3 */
694 /* ## Check for force inside potential check, i.e. this means we already did the potential part */
695 /* #if 'Force' in KERNEL_VF */
696 fvdw = _mm_mul_ps(_mm_sub_ps(vvdw12,vvdw6),rinvsq{I}{J});
697 /* #define INNERFLOPS INNERFLOPS+2 */
699 /* #elif KERNEL_VF=='Force' */
700 /* ## Force-only LennardJones makes it possible to save 1 flop (they do add up...) */
701 fvdw = _mm_mul_ps(_mm_sub_ps(_mm_mul_ps(c12_{I}{J},rinvsix),c6_{I}{J}),_mm_mul_ps(rinvsix,rinvsq{I}{J}));
702 /* #define INNERFLOPS INNERFLOPS+4 */
705 /* #elif KERNEL_VDW=='CubicSplineTable' */
707 /* CUBIC SPLINE TABLE DISPERSION */
708 /* #if 'Table' in KERNEL_ELEC */
709 vfitab = _mm_add_epi32(vfitab,ifour);
711 Y = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,0) );
712 F = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,1) );
713 G = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,2) );
714 H = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,3) );
715 _MM_TRANSPOSE4_PS(Y,F,G,H);
716 Heps = _mm_mul_ps(vfeps,H);
717 Fp = _mm_add_ps(F,_mm_mul_ps(vfeps,_mm_add_ps(G,Heps)));
718 /* #define INNERFLOPS INNERFLOPS+4 */
719 /* #if 'Potential' in KERNEL_VF */
720 VV = _mm_add_ps(Y,_mm_mul_ps(vfeps,Fp));
721 vvdw6 = _mm_mul_ps(c6_{I}{J},VV);
722 /* #define INNERFLOPS INNERFLOPS+3 */
724 /* #if 'Force' in KERNEL_VF */
725 FF = _mm_add_ps(Fp,_mm_mul_ps(vfeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
726 fvdw6 = _mm_mul_ps(c6_{I}{J},FF);
727 /* #define INNERFLOPS INNERFLOPS+4 */
730 /* CUBIC SPLINE TABLE REPULSION */
731 vfitab = _mm_add_epi32(vfitab,ifour);
732 Y = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,0) );
733 F = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,1) );
734 G = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,2) );
735 H = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,3) );
736 _MM_TRANSPOSE4_PS(Y,F,G,H);
737 Heps = _mm_mul_ps(vfeps,H);
738 Fp = _mm_add_ps(F,_mm_mul_ps(vfeps,_mm_add_ps(G,Heps)));
739 /* #define INNERFLOPS INNERFLOPS+4 */
740 /* #if 'Potential' in KERNEL_VF */
741 VV = _mm_add_ps(Y,_mm_mul_ps(vfeps,Fp));
742 vvdw12 = _mm_mul_ps(c12_{I}{J},VV);
743 /* #define INNERFLOPS INNERFLOPS+3 */
745 /* #if 'Force' in KERNEL_VF */
746 FF = _mm_add_ps(Fp,_mm_mul_ps(vfeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
747 fvdw12 = _mm_mul_ps(c12_{I}{J},FF);
748 /* #define INNERFLOPS INNERFLOPS+5 */
750 /* #if 'Potential' in KERNEL_VF */
751 vvdw = _mm_add_ps(vvdw12,vvdw6);
752 /* #define INNERFLOPS INNERFLOPS+1 */
754 /* #if 'Force' in KERNEL_VF */
755 fvdw = _mm_xor_ps(signbit,_mm_mul_ps(_mm_add_ps(fvdw6,fvdw12),_mm_mul_ps(vftabscale,rinv{I}{J})));
756 /* #define INNERFLOPS INNERFLOPS+4 */
759 /* ## End of check for vdw interaction forms */
761 /* ## END OF VDW INTERACTION CHECK FOR PAIR I-J */
763 /* #if 'switch' in INTERACTION_FLAGS[I][J] */
764 d = _mm_sub_ps(r{I}{J},rswitch);
765 d = _mm_max_ps(d,_mm_setzero_ps());
766 d2 = _mm_mul_ps(d,d);
767 sw = _mm_add_ps(one,_mm_mul_ps(d2,_mm_mul_ps(d,_mm_add_ps(swV3,_mm_mul_ps(d,_mm_add_ps(swV4,_mm_mul_ps(d,swV5)))))));
768 /* #define INNERFLOPS INNERFLOPS+10 */
770 /* #if 'Force' in KERNEL_VF */
771 dsw = _mm_mul_ps(d2,_mm_add_ps(swF2,_mm_mul_ps(d,_mm_add_ps(swF3,_mm_mul_ps(d,swF4)))));
772 /* #define INNERFLOPS INNERFLOPS+5 */
775 /* Evaluate switch function */
776 /* #if 'Force' in KERNEL_VF */
777 /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
778 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] and KERNEL_MOD_ELEC=='PotentialSwitch' */
779 felec = _mm_sub_ps( _mm_mul_ps(felec,sw) , _mm_mul_ps(rinv{I}{J},_mm_mul_ps(velec,dsw)) );
780 /* #define INNERFLOPS INNERFLOPS+4 */
782 /* #if 'vdw' in INTERACTION_FLAGS[I][J] and KERNEL_MOD_VDW=='PotentialSwitch' */
783 fvdw = _mm_sub_ps( _mm_mul_ps(fvdw,sw) , _mm_mul_ps(rinv{I}{J},_mm_mul_ps(vvdw,dsw)) );
784 /* #define INNERFLOPS INNERFLOPS+4 */
787 /* #if 'Potential' in KERNEL_VF */
788 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] and KERNEL_MOD_ELEC=='PotentialSwitch' */
789 velec = _mm_mul_ps(velec,sw);
790 /* #define INNERFLOPS INNERFLOPS+1 */
792 /* #if 'vdw' in INTERACTION_FLAGS[I][J] and KERNEL_MOD_VDW=='PotentialSwitch' */
793 vvdw = _mm_mul_ps(vvdw,sw);
794 /* #define INNERFLOPS INNERFLOPS+1 */
798 /* ## Note special check for TIP4P-TIP4P. Since we are cutting of all hydrogen interactions we also cut the LJ-only O-O interaction */
799 /* #if 'exactcutoff' in INTERACTION_FLAGS[I][J] or (GEOMETRY_I=='Water4' and GEOMETRY_J=='Water4' and 'exactcutoff' in INTERACTION_FLAGS[1][1]) */
800 cutoff_mask = _mm_cmplt_ps(rsq{I}{J},rcutoff2);
801 /* #define INNERFLOPS INNERFLOPS+1 */
804 /* #if 'Potential' in KERNEL_VF */
805 /* Update potential sum for this i atom from the interaction with this j atom. */
806 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] */
807 /* #if 'exactcutoff' in INTERACTION_FLAGS[I][J] */
808 velec = _mm_and_ps(velec,cutoff_mask);
809 /* #define INNERFLOPS INNERFLOPS+1 */
811 /* #if ROUND == 'Epilogue' */
812 velec = _mm_andnot_ps(dummy_mask,velec);
814 velecsum = _mm_add_ps(velecsum,velec);
815 /* #define INNERFLOPS INNERFLOPS+1 */
816 /* #if KERNEL_ELEC=='GeneralizedBorn' */
817 /* #if 'exactcutoff' in INTERACTION_FLAGS[I][J] */
818 vgb = _mm_and_ps(vgb,cutoff_mask);
819 /* #define INNERFLOPS INNERFLOPS+1 */
821 /* #if ROUND == 'Epilogue' */
822 vgb = _mm_andnot_ps(dummy_mask,vgb);
824 vgbsum = _mm_add_ps(vgbsum,vgb);
825 /* #define INNERFLOPS INNERFLOPS+1 */
828 /* #if 'vdw' in INTERACTION_FLAGS[I][J] */
829 /* ## Note special check for TIP4P-TIP4P. Since we are cutting of all hydrogen interactions we also cut the LJ-only O-O interaction */
830 /* #if 'exactcutoff' in INTERACTION_FLAGS[I][J] or (GEOMETRY_I=='Water4' and GEOMETRY_J=='Water4' and 'exactcutoff' in INTERACTION_FLAGS[1][1]) */
831 vvdw = _mm_and_ps(vvdw,cutoff_mask);
832 /* #define INNERFLOPS INNERFLOPS+1 */
834 /* #if ROUND == 'Epilogue' */
835 vvdw = _mm_andnot_ps(dummy_mask,vvdw);
837 vvdwsum = _mm_add_ps(vvdwsum,vvdw);
838 /* #define INNERFLOPS INNERFLOPS+1 */
842 /* #if 'Force' in KERNEL_VF */
844 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] and 'vdw' in INTERACTION_FLAGS[I][J] */
845 fscal = _mm_add_ps(felec,fvdw);
846 /* #define INNERFLOPS INNERFLOPS+1 */
847 /* #elif 'electrostatics' in INTERACTION_FLAGS[I][J] */
849 /* #elif 'vdw' in INTERACTION_FLAGS[I][J] */
853 /* ## Note special check for TIP4P-TIP4P. Since we are cutting of all hydrogen interactions we also cut the LJ-only O-O interaction */
854 /* #if 'exactcutoff' in INTERACTION_FLAGS[I][J] or (GEOMETRY_I=='Water4' and GEOMETRY_J=='Water4' and 'exactcutoff' in INTERACTION_FLAGS[1][1]) */
855 fscal = _mm_and_ps(fscal,cutoff_mask);
856 /* #define INNERFLOPS INNERFLOPS+1 */
859 /* #if ROUND == 'Epilogue' */
860 fscal = _mm_andnot_ps(dummy_mask,fscal);
863 /* Calculate temporary vectorial force */
864 tx = _mm_mul_ps(fscal,dx{I}{J});
865 ty = _mm_mul_ps(fscal,dy{I}{J});
866 tz = _mm_mul_ps(fscal,dz{I}{J});
868 /* Update vectorial force */
869 fix{I} = _mm_add_ps(fix{I},tx);
870 fiy{I} = _mm_add_ps(fiy{I},ty);
871 fiz{I} = _mm_add_ps(fiz{I},tz);
872 /* #define INNERFLOPS INNERFLOPS+6 */
874 /* #if GEOMETRY_I == 'Particle' */
875 /* #if ROUND == 'Loop' */
876 fjptrA = f+j_coord_offsetA;
877 fjptrB = f+j_coord_offsetB;
878 fjptrC = f+j_coord_offsetC;
879 fjptrD = f+j_coord_offsetD;
881 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
882 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
883 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
884 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
886 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
887 /* #define INNERFLOPS INNERFLOPS+3 */
889 fjx{J} = _mm_add_ps(fjx{J},tx);
890 fjy{J} = _mm_add_ps(fjy{J},ty);
891 fjz{J} = _mm_add_ps(fjz{J},tz);
892 /* #define INNERFLOPS INNERFLOPS+3 */
897 /* ## Note special check for TIP4P-TIP4P. Since we are cutting of all hydrogen interactions we also cut the LJ-only O-O interaction */
898 /* #if 'exactcutoff' in INTERACTION_FLAGS[I][J] or (GEOMETRY_I=='Water4' and GEOMETRY_J=='Water4' and 'exactcutoff' in INTERACTION_FLAGS[1][1]) */
899 /* #if 0 ## This and next two lines is a hack to maintain indentation in template file */
904 /* ## End of check for the interaction being outside the cutoff */
907 /* ## End of loop over i-j interaction pairs */
909 /* #if GEOMETRY_I != 'Particle' */
910 /* #if ROUND == 'Loop' */
911 fjptrA = f+j_coord_offsetA;
912 fjptrB = f+j_coord_offsetB;
913 fjptrC = f+j_coord_offsetC;
914 fjptrD = f+j_coord_offsetD;
916 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
917 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
918 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
919 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
923 /* #if 'Water' in GEOMETRY_I and GEOMETRY_J == 'Particle' */
924 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjx0,fjy0,fjz0);
925 /* #elif GEOMETRY_J == 'Water3' */
926 gmx_mm_decrement_3rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,
927 fjx0,fjy0,fjz0,fjx1,fjy1,fjz1,fjx2,fjy2,fjz2);
928 /* #define INNERFLOPS INNERFLOPS+9 */
929 /* #elif GEOMETRY_J == 'Water4' */
930 /* #if 0 in PARTICLES_J */
931 gmx_mm_decrement_4rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,
932 fjx0,fjy0,fjz0,fjx1,fjy1,fjz1,
933 fjx2,fjy2,fjz2,fjx3,fjy3,fjz3);
934 /* #define INNERFLOPS INNERFLOPS+12 */
936 gmx_mm_decrement_3rvec_4ptr_swizzle_ps(fjptrA+DIM,fjptrB+DIM,fjptrC+DIM,fjptrD+DIM,
937 fjx1,fjy1,fjz1,fjx2,fjy2,fjz2,fjx3,fjy3,fjz3);
938 /* #define INNERFLOPS INNERFLOPS+9 */
942 /* Inner loop uses {INNERFLOPS} flops */
947 /* End of innermost loop */
949 /* #if 'Force' in KERNEL_VF */
950 /* #if GEOMETRY_I == 'Particle' */
951 gmx_mm_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
952 f+i_coord_offset,fshift+i_shift_offset);
953 /* #define OUTERFLOPS OUTERFLOPS+6 */
954 /* #elif GEOMETRY_I == 'Water3' */
955 gmx_mm_update_iforce_3atom_swizzle_ps(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,
956 f+i_coord_offset,fshift+i_shift_offset);
957 /* #define OUTERFLOPS OUTERFLOPS+18 */
958 /* #elif GEOMETRY_I == 'Water4' */
959 /* #if 0 in PARTICLES_I */
960 gmx_mm_update_iforce_4atom_swizzle_ps(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
961 f+i_coord_offset,fshift+i_shift_offset);
962 /* #define OUTERFLOPS OUTERFLOPS+24 */
964 gmx_mm_update_iforce_3atom_swizzle_ps(fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
965 f+i_coord_offset+DIM,fshift+i_shift_offset);
966 /* #define OUTERFLOPS OUTERFLOPS+18 */
971 /* #if 'Potential' in KERNEL_VF */
973 /* Update potential energies */
974 /* #if KERNEL_ELEC != 'None' */
975 gmx_mm_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
976 /* #define OUTERFLOPS OUTERFLOPS+1 */
978 /* #if 'GeneralizedBorn' in KERNEL_ELEC */
979 gmx_mm_update_1pot_ps(vgbsum,kernel_data->energygrp_polarization+ggid);
980 /* #define OUTERFLOPS OUTERFLOPS+1 */
982 /* #if KERNEL_VDW != 'None' */
983 gmx_mm_update_1pot_ps(vvdwsum,kernel_data->energygrp_vdw+ggid);
984 /* #define OUTERFLOPS OUTERFLOPS+1 */
987 /* #if 'GeneralizedBorn' in KERNEL_ELEC and 'Force' in KERNEL_VF */
988 dvdasum = _mm_mul_ps(dvdasum, _mm_mul_ps(isai{I},isai{I}));
989 gmx_mm_update_1pot_ps(dvdasum,dvda+inr);
992 /* Increment number of inner iterations */
993 inneriter += j_index_end - j_index_start;
995 /* Outer loop uses {OUTERFLOPS} flops */
998 /* Increment number of outer iterations */
1001 /* Update outer/inner flops */
1002 /* ## NB: This is not important, it just affects the flopcount. However, since our preprocessor is */
1003 /* ## primitive and replaces aggressively even in strings inside these directives, we need to */
1004 /* ## assemble the main part of the name (containing KERNEL/ELEC/VDW) directly in the source. */
1005 /* #if GEOMETRY_I == 'Water3' */
1006 /* #define ISUFFIX '_W3' */
1007 /* #elif GEOMETRY_I == 'Water4' */
1008 /* #define ISUFFIX '_W4' */
1010 /* #define ISUFFIX '' */
1012 /* #if GEOMETRY_J == 'Water3' */
1013 /* #define JSUFFIX 'W3' */
1014 /* #elif GEOMETRY_J == 'Water4' */
1015 /* #define JSUFFIX 'W4' */
1017 /* #define JSUFFIX '' */
1019 /* #if 'PotentialAndForce' in KERNEL_VF */
1020 /* #define VFSUFFIX '_VF' */
1021 /* #elif 'Potential' in KERNEL_VF */
1022 /* #define VFSUFFIX '_V' */
1024 /* #define VFSUFFIX '_F' */
1027 /* #if KERNEL_ELEC != 'None' and KERNEL_VDW != 'None' */
1028 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW{ISUFFIX}{JSUFFIX}{VFSUFFIX},outeriter*{OUTERFLOPS} + inneriter*{INNERFLOPS});
1029 /* #elif KERNEL_ELEC != 'None' */
1030 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC{ISUFFIX}{JSUFFIX}{VFSUFFIX},outeriter*{OUTERFLOPS} + inneriter*{INNERFLOPS});
1032 inc_nrnb(nrnb,eNR_NBKERNEL_VDW{ISUFFIX}{JSUFFIX}{VFSUFFIX},outeriter*{OUTERFLOPS} + inneriter*{INNERFLOPS});