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 j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
81 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
82 real shX,shY,shZ,rcutoff_scalar;
83 real *shiftvec,*fshift,*x,*f;
84 __m128 tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
85 /* #for I in PARTICLES_I */
87 __m128 ix{I},iy{I},iz{I},fix{I},fiy{I},fiz{I},iq{I},isai{I};
89 /* #for J in PARTICLES_J */
90 int vdwjidx{J}A,vdwjidx{J}B,vdwjidx{J}C,vdwjidx{J}D;
91 __m128 jx{J},jy{J},jz{J},fjx{J},fjy{J},fjz{J},jq{J},isaj{J};
93 /* #for I,J in PAIRS_IJ */
94 __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};
96 /* #if KERNEL_ELEC != 'None' */
97 __m128 velec,felec,velecsum,facel,crf,krf,krf2;
100 /* #if 'GeneralizedBorn' in KERNEL_ELEC */
102 __m128 vgb,fgb,vgbsum,dvdasum,gbscale,gbtabscale,isaprod,gbqqfactor,gbinvepsdiff,dvdaj,gbeps,dvdatmp;
103 __m128 minushalf = _mm_set1_ps(-0.5);
104 real *invsqrta,*dvda,*gbtab;
106 /* #if KERNEL_VDW != 'None' */
108 __m128 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
111 __m128 one_sixth = _mm_set1_ps(1.0/6.0);
112 __m128 one_twelfth = _mm_set1_ps(1.0/12.0);
114 /* #if 'Table' in KERNEL_ELEC or 'GeneralizedBorn' in KERNEL_ELEC or 'Table' in KERNEL_VDW */
116 __m128i ifour = _mm_set1_epi32(4);
117 __m128 rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF;
120 /* #if 'Ewald' in KERNEL_ELEC */
122 __m128 ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
125 /* #if 'PotentialSwitch' in [KERNEL_MOD_ELEC,KERNEL_MOD_VDW] */
126 __m128 rswitch,swV3,swV4,swV5,swF2,swF3,swF4,d,d2,sw,dsw;
127 real rswitch_scalar,d_scalar;
129 __m128 dummy_mask,cutoff_mask;
130 __m128 signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
131 __m128 one = _mm_set1_ps(1.0);
132 __m128 two = _mm_set1_ps(2.0);
138 jindex = nlist->jindex;
140 shiftidx = nlist->shift;
142 shiftvec = fr->shift_vec[0];
143 fshift = fr->fshift[0];
144 /* #if KERNEL_ELEC != 'None' */
145 facel = _mm_set1_ps(fr->epsfac);
146 charge = mdatoms->chargeA;
147 /* #if 'ReactionField' in KERNEL_ELEC */
148 krf = _mm_set1_ps(fr->ic->k_rf);
149 krf2 = _mm_set1_ps(fr->ic->k_rf*2.0);
150 crf = _mm_set1_ps(fr->ic->c_rf);
153 /* #if KERNEL_VDW != 'None' */
154 nvdwtype = fr->ntype;
156 vdwtype = mdatoms->typeA;
159 /* #if 'Table' in KERNEL_ELEC and 'Table' in KERNEL_VDW */
160 vftab = kernel_data->table_elec_vdw->data;
161 vftabscale = _mm_set1_ps(kernel_data->table_elec_vdw->scale);
162 /* #elif 'Table' in KERNEL_ELEC */
163 vftab = kernel_data->table_elec->data;
164 vftabscale = _mm_set1_ps(kernel_data->table_elec->scale);
165 /* #elif 'Table' in KERNEL_VDW */
166 vftab = kernel_data->table_vdw->data;
167 vftabscale = _mm_set1_ps(kernel_data->table_vdw->scale);
170 /* #if 'Ewald' in KERNEL_ELEC */
171 sh_ewald = _mm_set1_ps(fr->ic->sh_ewald);
172 /* #if KERNEL_VF=='Force' and KERNEL_MOD_ELEC!='PotentialSwitch' */
173 ewtab = fr->ic->tabq_coul_F;
174 ewtabscale = _mm_set1_ps(fr->ic->tabq_scale);
175 ewtabhalfspace = _mm_set1_ps(0.5/fr->ic->tabq_scale);
177 ewtab = fr->ic->tabq_coul_FDV0;
178 ewtabscale = _mm_set1_ps(fr->ic->tabq_scale);
179 ewtabhalfspace = _mm_set1_ps(0.5/fr->ic->tabq_scale);
183 /* #if KERNEL_ELEC=='GeneralizedBorn' */
184 invsqrta = fr->invsqrta;
186 gbtabscale = _mm_set1_ps(fr->gbtab.scale);
187 gbtab = fr->gbtab.data;
188 gbinvepsdiff = _mm_set1_ps((1.0/fr->epsilon_r) - (1.0/fr->gb_epsilon_solvent));
191 /* #if 'Water' in GEOMETRY_I */
192 /* Setup water-specific parameters */
193 inr = nlist->iinr[0];
194 /* #for I in PARTICLES_ELEC_I */
195 iq{I} = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+{I}]));
197 /* #for I in PARTICLES_VDW_I */
198 vdwioffset{I} = 2*nvdwtype*vdwtype[inr+{I}];
202 /* #if 'Water' in GEOMETRY_J */
203 /* #for J in PARTICLES_ELEC_J */
204 jq{J} = _mm_set1_ps(charge[inr+{J}]);
206 /* #for J in PARTICLES_VDW_J */
207 vdwjidx{J}A = 2*vdwtype[inr+{J}];
209 /* #for I,J in PAIRS_IJ */
210 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] */
211 qq{I}{J} = _mm_mul_ps(iq{I},jq{J});
213 /* #if 'vdw' in INTERACTION_FLAGS[I][J] */
214 c6_{I}{J} = _mm_set1_ps(vdwparam[vdwioffset{I}+vdwjidx{J}A]);
215 c12_{I}{J} = _mm_set1_ps(vdwparam[vdwioffset{I}+vdwjidx{J}A+1]);
220 /* #if KERNEL_MOD_ELEC!='None' or KERNEL_MOD_VDW!='None' */
221 /* #if KERNEL_ELEC!='None' */
222 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
223 rcutoff_scalar = fr->rcoulomb;
225 rcutoff_scalar = fr->rvdw;
227 rcutoff = _mm_set1_ps(rcutoff_scalar);
228 rcutoff2 = _mm_mul_ps(rcutoff,rcutoff);
231 /* #if KERNEL_MOD_VDW=='PotentialShift' */
232 sh_vdw_invrcut6 = _mm_set1_ps(fr->ic->sh_invrc6);
233 rvdw = _mm_set1_ps(fr->rvdw);
236 /* #if 'PotentialSwitch' in [KERNEL_MOD_ELEC,KERNEL_MOD_VDW] */
237 /* #if KERNEL_MOD_ELEC=='PotentialSwitch' */
238 rswitch_scalar = fr->rcoulomb_switch;
239 rswitch = _mm_set1_ps(rswitch_scalar);
241 rswitch_scalar = fr->rvdw_switch;
242 rswitch = _mm_set1_ps(rswitch_scalar);
244 /* Setup switch parameters */
245 d_scalar = rcutoff_scalar-rswitch_scalar;
246 d = _mm_set1_ps(d_scalar);
247 swV3 = _mm_set1_ps(-10.0/(d_scalar*d_scalar*d_scalar));
248 swV4 = _mm_set1_ps( 15.0/(d_scalar*d_scalar*d_scalar*d_scalar));
249 swV5 = _mm_set1_ps( -6.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
250 /* #if 'Force' in KERNEL_VF */
251 swF2 = _mm_set1_ps(-30.0/(d_scalar*d_scalar*d_scalar));
252 swF3 = _mm_set1_ps( 60.0/(d_scalar*d_scalar*d_scalar*d_scalar));
253 swF4 = _mm_set1_ps(-30.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
257 /* Avoid stupid compiler warnings */
258 jnrA = jnrB = jnrC = jnrD = 0;
264 /* ## Keep track of the floating point operations we issue for reporting! */
265 /* #define OUTERFLOPS 0 */
269 /* Start outer loop over neighborlists */
270 for(iidx=0; iidx<nri; iidx++)
272 /* Load shift vector for this list */
273 i_shift_offset = DIM*shiftidx[iidx];
274 shX = shiftvec[i_shift_offset+XX];
275 shY = shiftvec[i_shift_offset+YY];
276 shZ = shiftvec[i_shift_offset+ZZ];
278 /* Load limits for loop over neighbors */
279 j_index_start = jindex[iidx];
280 j_index_end = jindex[iidx+1];
282 /* Get outer coordinate index */
284 i_coord_offset = DIM*inr;
286 /* Load i particle coords and add shift vector */
287 /* ## Loop over i particles, but only include ones that we use - skip e.g. vdw-only sites for elec-only kernel */
288 /* #for I in PARTICLES_I */
289 ix{I} = _mm_set1_ps(shX + x[i_coord_offset+DIM*{I}+XX]);
290 iy{I} = _mm_set1_ps(shY + x[i_coord_offset+DIM*{I}+YY]);
291 iz{I} = _mm_set1_ps(shZ + x[i_coord_offset+DIM*{I}+ZZ]);
292 /* #define OUTERFLOPS OUTERFLOPS+3 */
295 /* #if 'Force' in KERNEL_VF */
296 /* #for I in PARTICLES_I */
297 fix{I} = _mm_setzero_ps();
298 fiy{I} = _mm_setzero_ps();
299 fiz{I} = _mm_setzero_ps();
303 /* ## For water we already preloaded parameters at the start of the kernel */
304 /* #if not 'Water' in GEOMETRY_I */
305 /* Load parameters for i particles */
306 /* #for I in PARTICLES_ELEC_I */
307 iq{I} = _mm_mul_ps(facel,_mm_load1_ps(charge+inr+{I}));
308 /* #define OUTERFLOPS OUTERFLOPS+1 */
309 /* #if KERNEL_ELEC=='GeneralizedBorn' */
310 isai{I} = _mm_load1_ps(invsqrta+inr+{I});
313 /* #for I in PARTICLES_VDW_I */
314 vdwioffset{I} = 2*nvdwtype*vdwtype[inr+{I}];
318 /* #if 'Potential' in KERNEL_VF */
319 /* Reset potential sums */
320 /* #if KERNEL_ELEC != 'None' */
321 velecsum = _mm_setzero_ps();
323 /* #if 'GeneralizedBorn' in KERNEL_ELEC */
324 vgbsum = _mm_setzero_ps();
326 /* #if KERNEL_VDW != 'None' */
327 vvdwsum = _mm_setzero_ps();
330 /* #if 'GeneralizedBorn' in KERNEL_ELEC and 'Force' in KERNEL_VF */
331 dvdasum = _mm_setzero_ps();
334 /* #for ROUND in ['Loop','Epilogue'] */
336 /* #if ROUND =='Loop' */
337 /* Start inner kernel loop */
338 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
340 /* ## First round is normal loop (next statement resets indentation) */
347 /* ## Second round is epilogue */
349 /* #define INNERFLOPS 0 */
351 /* Get j neighbor index, and coordinate index */
357 /* #if ROUND =='Epilogue' */
358 /* Sign of each element will be negative for non-real atoms.
359 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
360 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
362 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
363 jnrA = (jnrA>=0) ? jnrA : 0;
364 jnrB = (jnrB>=0) ? jnrB : 0;
365 jnrC = (jnrC>=0) ? jnrC : 0;
366 jnrD = (jnrD>=0) ? jnrD : 0;
369 j_coord_offsetA = DIM*jnrA;
370 j_coord_offsetB = DIM*jnrB;
371 j_coord_offsetC = DIM*jnrC;
372 j_coord_offsetD = DIM*jnrD;
374 /* load j atom coordinates */
375 /* #if GEOMETRY_J == 'Particle' */
376 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
377 x+j_coord_offsetC,x+j_coord_offsetD,
379 /* #elif GEOMETRY_J == 'Water3' */
380 gmx_mm_load_3rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
381 x+j_coord_offsetC,x+j_coord_offsetD,
382 &jx0,&jy0,&jz0,&jx1,&jy1,&jz1,&jx2,&jy2,&jz2);
383 /* #elif GEOMETRY_J == 'Water4' */
384 /* #if 0 in PARTICLES_J */
385 gmx_mm_load_4rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
386 x+j_coord_offsetC,x+j_coord_offsetD,
387 &jx0,&jy0,&jz0,&jx1,&jy1,&jz1,&jx2,
388 &jy2,&jz2,&jx3,&jy3,&jz3);
390 gmx_mm_load_3rvec_4ptr_swizzle_ps(x+j_coord_offsetA+DIM,x+j_coord_offsetB+DIM,
391 x+j_coord_offsetC+DIM,x+j_coord_offsetD+DIM,
392 &jx1,&jy1,&jz1,&jx2,&jy2,&jz2,&jx3,&jy3,&jz3);
396 /* Calculate displacement vector */
397 /* #for I,J in PAIRS_IJ */
398 dx{I}{J} = _mm_sub_ps(ix{I},jx{J});
399 dy{I}{J} = _mm_sub_ps(iy{I},jy{J});
400 dz{I}{J} = _mm_sub_ps(iz{I},jz{J});
401 /* #define INNERFLOPS INNERFLOPS+3 */
404 /* Calculate squared distance and things based on it */
405 /* #for I,J in PAIRS_IJ */
406 rsq{I}{J} = gmx_mm_calc_rsq_ps(dx{I}{J},dy{I}{J},dz{I}{J});
407 /* #define INNERFLOPS INNERFLOPS+5 */
410 /* #for I,J in PAIRS_IJ */
411 /* #if 'rinv' in INTERACTION_FLAGS[I][J] */
412 rinv{I}{J} = gmx_mm_invsqrt_ps(rsq{I}{J});
413 /* #define INNERFLOPS INNERFLOPS+5 */
417 /* #for I,J in PAIRS_IJ */
418 /* #if 'rinvsq' in INTERACTION_FLAGS[I][J] */
419 /* # if 'rinv' not in INTERACTION_FLAGS[I][J] */
420 rinvsq{I}{J} = gmx_mm_inv_ps(rsq{I}{J});
421 /* #define INNERFLOPS INNERFLOPS+4 */
423 rinvsq{I}{J} = _mm_mul_ps(rinv{I}{J},rinv{I}{J});
424 /* #define INNERFLOPS INNERFLOPS+1 */
429 /* #if not 'Water' in GEOMETRY_J */
430 /* Load parameters for j particles */
431 /* #for J in PARTICLES_ELEC_J */
432 jq{J} = gmx_mm_load_4real_swizzle_ps(charge+jnrA+{J},charge+jnrB+{J},
433 charge+jnrC+{J},charge+jnrD+{J});
434 /* #if KERNEL_ELEC=='GeneralizedBorn' */
435 isaj{J} = gmx_mm_load_4real_swizzle_ps(invsqrta+jnrA+{J},invsqrta+jnrB+{J},
436 invsqrta+jnrC+{J},invsqrta+jnrD+{J});
439 /* #for J in PARTICLES_VDW_J */
440 vdwjidx{J}A = 2*vdwtype[jnrA+{J}];
441 vdwjidx{J}B = 2*vdwtype[jnrB+{J}];
442 vdwjidx{J}C = 2*vdwtype[jnrC+{J}];
443 vdwjidx{J}D = 2*vdwtype[jnrD+{J}];
447 /* #if 'Force' in KERNEL_VF and not 'Particle' in GEOMETRY_J */
448 /* #for J in PARTICLES_J */
449 fjx{J} = _mm_setzero_ps();
450 fjy{J} = _mm_setzero_ps();
451 fjz{J} = _mm_setzero_ps();
455 /* #for I,J in PAIRS_IJ */
457 /**************************
458 * CALCULATE INTERACTIONS *
459 **************************/
461 /* #if 'exactcutoff' in INTERACTION_FLAGS[I][J] */
462 /* ## We always calculate rinv/rinvsq above to enable pipelineing in compilers (performance tested on x86) */
463 if (gmx_mm_any_lt(rsq{I}{J},rcutoff2))
465 /* #if 0 ## this and the next two lines is a hack to maintain auto-indentation in template file */
468 /* #define INNERFLOPS INNERFLOPS+1 */
471 /* #if 'r' in INTERACTION_FLAGS[I][J] */
472 r{I}{J} = _mm_mul_ps(rsq{I}{J},rinv{I}{J});
473 /* #if ROUND == 'Epilogue' */
474 r{I}{J} = _mm_andnot_ps(dummy_mask,r{I}{J});
475 /* #define INNERFLOPS INNERFLOPS+1 */
477 /* #define INNERFLOPS INNERFLOPS+1 */
480 /* ## For water geometries we already loaded parameters at the start of the kernel */
481 /* #if not 'Water' in GEOMETRY_J */
482 /* Compute parameters for interactions between i and j atoms */
483 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] */
484 qq{I}{J} = _mm_mul_ps(iq{I},jq{J});
485 /* #define INNERFLOPS INNERFLOPS+1 */
487 /* #if 'vdw' in INTERACTION_FLAGS[I][J] */
488 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset{I}+vdwjidx{J}A,
489 vdwparam+vdwioffset{I}+vdwjidx{J}B,
490 vdwparam+vdwioffset{I}+vdwjidx{J}C,
491 vdwparam+vdwioffset{I}+vdwjidx{J}D,
492 &c6_{I}{J},&c12_{I}{J});
496 /* #if 'table' in INTERACTION_FLAGS[I][J] */
497 /* Calculate table index by multiplying r with table scale and truncate to integer */
498 rt = _mm_mul_ps(r{I}{J},vftabscale);
499 vfitab = _mm_cvttps_epi32(rt);
500 vfeps = _mm_sub_ps(rt,_mm_cvtepi32_ps(vfitab));
501 /* #define INNERFLOPS INNERFLOPS+4 */
502 /* #if 'Table' in KERNEL_ELEC and 'Table' in KERNEL_VDW */
503 /* ## 3 tables, 4 bytes per point: multiply index by 12 */
504 vfitab = _mm_slli_epi32(_mm_add_epi32(vfitab,_mm_slli_epi32(vfitab,1)),2);
505 /* #elif 'Table' in KERNEL_ELEC */
506 /* ## 1 table, 4 bytes per point: multiply index by 4 */
507 vfitab = _mm_slli_epi32(vfitab,2);
508 /* #elif 'Table' in KERNEL_VDW */
509 /* ## 2 tables, 4 bytes per point: multiply index by 8 */
510 vfitab = _mm_slli_epi32(vfitab,3);
514 /* ## ELECTROSTATIC INTERACTIONS */
515 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] */
517 /* #if KERNEL_ELEC=='Coulomb' */
519 /* COULOMB ELECTROSTATICS */
520 velec = _mm_mul_ps(qq{I}{J},rinv{I}{J});
521 /* #define INNERFLOPS INNERFLOPS+1 */
522 /* #if 'Force' in KERNEL_VF */
523 felec = _mm_mul_ps(velec,rinvsq{I}{J});
524 /* #define INNERFLOPS INNERFLOPS+2 */
527 /* #elif KERNEL_ELEC=='ReactionField' */
529 /* REACTION-FIELD ELECTROSTATICS */
530 /* #if 'Potential' in KERNEL_VF */
531 velec = _mm_mul_ps(qq{I}{J},_mm_sub_ps(_mm_add_ps(rinv{I}{J},_mm_mul_ps(krf,rsq{I}{J})),crf));
532 /* #define INNERFLOPS INNERFLOPS+4 */
534 /* #if 'Force' in KERNEL_VF */
535 felec = _mm_mul_ps(qq{I}{J},_mm_sub_ps(_mm_mul_ps(rinv{I}{J},rinvsq{I}{J}),krf2));
536 /* #define INNERFLOPS INNERFLOPS+3 */
539 /* #elif KERNEL_ELEC=='GeneralizedBorn' */
541 /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
542 isaprod = _mm_mul_ps(isai{I},isaj{J});
543 gbqqfactor = _mm_xor_ps(signbit,_mm_mul_ps(qq{I}{J},_mm_mul_ps(isaprod,gbinvepsdiff)));
544 gbscale = _mm_mul_ps(isaprod,gbtabscale);
545 dvdaj = gmx_mm_load_4real_swizzle_ps(dvda+jnrA+{J},dvda+jnrB+{J},dvda+jnrC+{J},dvda+jnrD+{J});
546 /* #define INNERFLOPS INNERFLOPS+5 */
548 /* Calculate generalized born table index - this is a separate table from the normal one,
549 * but we use the same procedure by multiplying r with scale and truncating to integer.
551 rt = _mm_mul_ps(r{I}{J},gbscale);
552 gbitab = _mm_cvttps_epi32(rt);
553 gbeps = _mm_sub_ps(rt,_mm_cvtepi32_ps(gbitab));
554 gbitab = _mm_slli_epi32(gbitab,2);
556 Y = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,0) );
557 F = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,1) );
558 G = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,2) );
559 H = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,3) );
560 _MM_TRANSPOSE4_PS(Y,F,G,H);
561 Heps = _mm_mul_ps(gbeps,H);
562 Fp = _mm_add_ps(F,_mm_mul_ps(gbeps,_mm_add_ps(G,Heps)));
563 VV = _mm_add_ps(Y,_mm_mul_ps(gbeps,Fp));
564 vgb = _mm_mul_ps(gbqqfactor,VV);
565 /* #define INNERFLOPS INNERFLOPS+10 */
567 /* #if 'Force' in KERNEL_VF */
568 FF = _mm_add_ps(Fp,_mm_mul_ps(gbeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
569 fgb = _mm_mul_ps(gbqqfactor,_mm_mul_ps(FF,gbscale));
570 dvdatmp = _mm_mul_ps(minushalf,_mm_add_ps(vgb,_mm_mul_ps(fgb,r{I}{J})));
571 dvdasum = _mm_add_ps(dvdasum,dvdatmp);
572 gmx_mm_store_4real_swizzle_ps(dvda+jnrA,dvda+jnrB,dvda+jnrC,dvda+jnrD,
573 _mm_mul_ps(dvdatmp,_mm_mul_ps(isaj{J},isaj{J})));
574 /* #define INNERFLOPS INNERFLOPS+13 */
576 velec = _mm_mul_ps(qq{I}{J},rinv{I}{J});
577 /* #define INNERFLOPS INNERFLOPS+1 */
578 /* #if 'Force' in KERNEL_VF */
579 felec = _mm_mul_ps(_mm_sub_ps(_mm_mul_ps(velec,rinv{I}{J}),fgb),rinv{I}{J});
580 /* #define INNERFLOPS INNERFLOPS+3 */
583 /* #elif KERNEL_ELEC=='Ewald' */
584 /* EWALD ELECTROSTATICS */
586 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
587 ewrt = _mm_mul_ps(r{I}{J},ewtabscale);
588 ewitab = _mm_cvttps_epi32(ewrt);
589 eweps = _mm_sub_ps(ewrt,_mm_cvtepi32_ps(ewitab));
590 /* #define INNERFLOPS INNERFLOPS+4 */
591 /* #if 'Potential' in KERNEL_VF or KERNEL_MOD_ELEC=='PotentialSwitch' */
592 ewitab = _mm_slli_epi32(ewitab,2);
593 ewtabF = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,0) );
594 ewtabD = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,1) );
595 ewtabV = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,2) );
596 ewtabFn = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,3) );
597 _MM_TRANSPOSE4_PS(ewtabF,ewtabD,ewtabV,ewtabFn);
598 felec = _mm_add_ps(ewtabF,_mm_mul_ps(eweps,ewtabD));
599 /* #define INNERFLOPS INNERFLOPS+2 */
600 /* #if KERNEL_MOD_ELEC=='PotentialShift' */
601 velec = _mm_sub_ps(ewtabV,_mm_mul_ps(_mm_mul_ps(ewtabhalfspace,eweps),_mm_add_ps(ewtabF,felec)));
602 velec = _mm_mul_ps(qq{I}{J},_mm_sub_ps(_mm_sub_ps(rinv{I}{J},sh_ewald),velec));
603 /* #define INNERFLOPS INNERFLOPS+7 */
605 velec = _mm_sub_ps(ewtabV,_mm_mul_ps(_mm_mul_ps(ewtabhalfspace,eweps),_mm_add_ps(ewtabF,felec)));
606 velec = _mm_mul_ps(qq{I}{J},_mm_sub_ps(rinv{I}{J},velec));
607 /* #define INNERFLOPS INNERFLOPS+6 */
609 /* #if 'Force' in KERNEL_VF */
610 felec = _mm_mul_ps(_mm_mul_ps(qq{I}{J},rinv{I}{J}),_mm_sub_ps(rinvsq{I}{J},felec));
611 /* #define INNERFLOPS INNERFLOPS+3 */
613 /* #elif KERNEL_VF=='Force' */
614 gmx_mm_load_4pair_swizzle_ps(ewtab+gmx_mm_extract_epi32(ewitab,0),ewtab+gmx_mm_extract_epi32(ewitab,1),
615 ewtab+gmx_mm_extract_epi32(ewitab,2),ewtab+gmx_mm_extract_epi32(ewitab,3),
617 felec = _mm_add_ps(_mm_mul_ps( _mm_sub_ps(one,eweps),ewtabF),_mm_mul_ps(eweps,ewtabFn));
618 felec = _mm_mul_ps(_mm_mul_ps(qq{I}{J},rinv{I}{J}),_mm_sub_ps(rinvsq{I}{J},felec));
619 /* #define INNERFLOPS INNERFLOPS+7 */
622 /* #elif KERNEL_ELEC=='CubicSplineTable' */
624 /* CUBIC SPLINE TABLE ELECTROSTATICS */
625 Y = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,0) );
626 F = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,1) );
627 G = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,2) );
628 H = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,3) );
629 _MM_TRANSPOSE4_PS(Y,F,G,H);
630 Heps = _mm_mul_ps(vfeps,H);
631 Fp = _mm_add_ps(F,_mm_mul_ps(vfeps,_mm_add_ps(G,Heps)));
632 /* #define INNERFLOPS INNERFLOPS+4 */
633 /* #if 'Potential' in KERNEL_VF */
634 VV = _mm_add_ps(Y,_mm_mul_ps(vfeps,Fp));
635 velec = _mm_mul_ps(qq{I}{J},VV);
636 /* #define INNERFLOPS INNERFLOPS+3 */
638 /* #if 'Force' in KERNEL_VF */
639 FF = _mm_add_ps(Fp,_mm_mul_ps(vfeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
640 felec = _mm_xor_ps(signbit,_mm_mul_ps(_mm_mul_ps(qq{I}{J},FF),_mm_mul_ps(vftabscale,rinv{I}{J})));
641 /* #define INNERFLOPS INNERFLOPS+7 */
644 /* ## End of check for electrostatics interaction forms */
646 /* ## END OF ELECTROSTATIC INTERACTION CHECK FOR PAIR I-J */
648 /* #if 'vdw' in INTERACTION_FLAGS[I][J] */
650 /* #if KERNEL_VDW=='LennardJones' */
652 /* LENNARD-JONES DISPERSION/REPULSION */
654 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq{I}{J},rinvsq{I}{J}),rinvsq{I}{J});
655 /* #define INNERFLOPS INNERFLOPS+2 */
656 /* #if 'Potential' in KERNEL_VF or KERNEL_MOD_VDW=='PotentialSwitch' */
657 vvdw6 = _mm_mul_ps(c6_{I}{J},rinvsix);
658 vvdw12 = _mm_mul_ps(c12_{I}{J},_mm_mul_ps(rinvsix,rinvsix));
659 /* #define INNERFLOPS INNERFLOPS+3 */
660 /* #if KERNEL_MOD_VDW=='PotentialShift' */
661 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) ,
662 _mm_mul_ps( _mm_sub_ps(vvdw6,_mm_mul_ps(c6_{I}{J},sh_vdw_invrcut6)),one_sixth));
663 /* #define INNERFLOPS INNERFLOPS+8 */
665 vvdw = _mm_sub_ps( _mm_mul_ps(vvdw12,one_twelfth) , _mm_mul_ps(vvdw6,one_sixth) );
666 /* #define INNERFLOPS INNERFLOPS+3 */
668 /* ## Check for force inside potential check, i.e. this means we already did the potential part */
669 /* #if 'Force' in KERNEL_VF */
670 fvdw = _mm_mul_ps(_mm_sub_ps(vvdw12,vvdw6),rinvsq{I}{J});
671 /* #define INNERFLOPS INNERFLOPS+2 */
673 /* #elif KERNEL_VF=='Force' */
674 /* ## Force-only LennardJones makes it possible to save 1 flop (they do add up...) */
675 fvdw = _mm_mul_ps(_mm_sub_ps(_mm_mul_ps(c12_{I}{J},rinvsix),c6_{I}{J}),_mm_mul_ps(rinvsix,rinvsq{I}{J}));
676 /* #define INNERFLOPS INNERFLOPS+4 */
679 /* #elif KERNEL_VDW=='CubicSplineTable' */
681 /* CUBIC SPLINE TABLE DISPERSION */
682 /* #if 'Table' in KERNEL_ELEC */
683 vfitab = _mm_add_epi32(vfitab,ifour);
685 Y = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,0) );
686 F = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,1) );
687 G = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,2) );
688 H = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,3) );
689 _MM_TRANSPOSE4_PS(Y,F,G,H);
690 Heps = _mm_mul_ps(vfeps,H);
691 Fp = _mm_add_ps(F,_mm_mul_ps(vfeps,_mm_add_ps(G,Heps)));
692 /* #define INNERFLOPS INNERFLOPS+4 */
693 /* #if 'Potential' in KERNEL_VF */
694 VV = _mm_add_ps(Y,_mm_mul_ps(vfeps,Fp));
695 vvdw6 = _mm_mul_ps(c6_{I}{J},VV);
696 /* #define INNERFLOPS INNERFLOPS+3 */
698 /* #if 'Force' in KERNEL_VF */
699 FF = _mm_add_ps(Fp,_mm_mul_ps(vfeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
700 fvdw6 = _mm_mul_ps(c6_{I}{J},FF);
701 /* #define INNERFLOPS INNERFLOPS+4 */
704 /* CUBIC SPLINE TABLE REPULSION */
705 vfitab = _mm_add_epi32(vfitab,ifour);
706 Y = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,0) );
707 F = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,1) );
708 G = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,2) );
709 H = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,3) );
710 _MM_TRANSPOSE4_PS(Y,F,G,H);
711 Heps = _mm_mul_ps(vfeps,H);
712 Fp = _mm_add_ps(F,_mm_mul_ps(vfeps,_mm_add_ps(G,Heps)));
713 /* #define INNERFLOPS INNERFLOPS+4 */
714 /* #if 'Potential' in KERNEL_VF */
715 VV = _mm_add_ps(Y,_mm_mul_ps(vfeps,Fp));
716 vvdw12 = _mm_mul_ps(c12_{I}{J},VV);
717 /* #define INNERFLOPS INNERFLOPS+3 */
719 /* #if 'Force' in KERNEL_VF */
720 FF = _mm_add_ps(Fp,_mm_mul_ps(vfeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
721 fvdw12 = _mm_mul_ps(c12_{I}{J},FF);
722 /* #define INNERFLOPS INNERFLOPS+5 */
724 /* #if 'Potential' in KERNEL_VF */
725 vvdw = _mm_add_ps(vvdw12,vvdw6);
726 /* #define INNERFLOPS INNERFLOPS+1 */
728 /* #if 'Force' in KERNEL_VF */
729 fvdw = _mm_xor_ps(signbit,_mm_mul_ps(_mm_add_ps(fvdw6,fvdw12),_mm_mul_ps(vftabscale,rinv{I}{J})));
730 /* #define INNERFLOPS INNERFLOPS+4 */
733 /* ## End of check for vdw interaction forms */
735 /* ## END OF VDW INTERACTION CHECK FOR PAIR I-J */
737 /* #if 'switch' in INTERACTION_FLAGS[I][J] */
738 d = _mm_sub_ps(r{I}{J},rswitch);
739 d = _mm_max_ps(d,_mm_setzero_ps());
740 d2 = _mm_mul_ps(d,d);
741 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)))))));
742 /* #define INNERFLOPS INNERFLOPS+10 */
744 /* #if 'Force' in KERNEL_VF */
745 dsw = _mm_mul_ps(d2,_mm_add_ps(swF2,_mm_mul_ps(d,_mm_add_ps(swF3,_mm_mul_ps(d,swF4)))));
746 /* #define INNERFLOPS INNERFLOPS+5 */
749 /* Evaluate switch function */
750 /* #if 'Force' in KERNEL_VF */
751 /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
752 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] and KERNEL_MOD_ELEC=='PotentialSwitch' */
753 felec = _mm_sub_ps( _mm_mul_ps(felec,sw) , _mm_mul_ps(rinv{I}{J},_mm_mul_ps(velec,dsw)) );
754 /* #define INNERFLOPS INNERFLOPS+4 */
756 /* #if 'vdw' in INTERACTION_FLAGS[I][J] and KERNEL_MOD_VDW=='PotentialSwitch' */
757 fvdw = _mm_sub_ps( _mm_mul_ps(fvdw,sw) , _mm_mul_ps(rinv{I}{J},_mm_mul_ps(vvdw,dsw)) );
758 /* #define INNERFLOPS INNERFLOPS+4 */
761 /* #if 'Potential' in KERNEL_VF */
762 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] and KERNEL_MOD_ELEC=='PotentialSwitch' */
763 velec = _mm_mul_ps(velec,sw);
764 /* #define INNERFLOPS INNERFLOPS+1 */
766 /* #if 'vdw' in INTERACTION_FLAGS[I][J] and KERNEL_MOD_VDW=='PotentialSwitch' */
767 vvdw = _mm_mul_ps(vvdw,sw);
768 /* #define INNERFLOPS INNERFLOPS+1 */
772 /* #if 'exactcutoff' in INTERACTION_FLAGS[I][J] */
773 cutoff_mask = _mm_cmplt_ps(rsq{I}{J},rcutoff2);
774 /* #define INNERFLOPS INNERFLOPS+1 */
777 /* #if 'Potential' in KERNEL_VF */
778 /* Update potential sum for this i atom from the interaction with this j atom. */
779 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] */
780 /* #if 'exactcutoff' in INTERACTION_FLAGS[I][J] */
781 velec = _mm_and_ps(velec,cutoff_mask);
782 /* #define INNERFLOPS INNERFLOPS+1 */
784 /* #if ROUND == 'Epilogue' */
785 velec = _mm_andnot_ps(dummy_mask,velec);
787 velecsum = _mm_add_ps(velecsum,velec);
788 /* #define INNERFLOPS INNERFLOPS+1 */
789 /* #if KERNEL_ELEC=='GeneralizedBorn' */
790 /* #if 'exactcutoff' in INTERACTION_FLAGS[I][J] */
791 vgb = _mm_and_ps(vgb,cutoff_mask);
792 /* #define INNERFLOPS INNERFLOPS+1 */
794 /* #if ROUND == 'Epilogue' */
795 vgb = _mm_andnot_ps(dummy_mask,vgb);
797 vgbsum = _mm_add_ps(vgbsum,vgb);
798 /* #define INNERFLOPS INNERFLOPS+1 */
801 /* #if 'vdw' in INTERACTION_FLAGS[I][J] */
802 /* #if 'exactcutoff' in INTERACTION_FLAGS[I][J] */
803 vvdw = _mm_and_ps(vvdw,cutoff_mask);
804 /* #define INNERFLOPS INNERFLOPS+1 */
806 /* #if ROUND == 'Epilogue' */
807 vvdw = _mm_andnot_ps(dummy_mask,vvdw);
809 vvdwsum = _mm_add_ps(vvdwsum,vvdw);
810 /* #define INNERFLOPS INNERFLOPS+1 */
814 /* #if 'Force' in KERNEL_VF */
816 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] and 'vdw' in INTERACTION_FLAGS[I][J] */
817 fscal = _mm_add_ps(felec,fvdw);
818 /* #define INNERFLOPS INNERFLOPS+1 */
819 /* #elif 'electrostatics' in INTERACTION_FLAGS[I][J] */
821 /* #elif 'vdw' in INTERACTION_FLAGS[I][J] */
825 /* #if 'exactcutoff' in INTERACTION_FLAGS[I][J] */
826 fscal = _mm_and_ps(fscal,cutoff_mask);
827 /* #define INNERFLOPS INNERFLOPS+1 */
830 /* #if ROUND == 'Epilogue' */
831 fscal = _mm_andnot_ps(dummy_mask,fscal);
834 /* Calculate temporary vectorial force */
835 tx = _mm_mul_ps(fscal,dx{I}{J});
836 ty = _mm_mul_ps(fscal,dy{I}{J});
837 tz = _mm_mul_ps(fscal,dz{I}{J});
839 /* Update vectorial force */
840 fix{I} = _mm_add_ps(fix{I},tx);
841 fiy{I} = _mm_add_ps(fiy{I},ty);
842 fiz{I} = _mm_add_ps(fiz{I},tz);
843 /* #define INNERFLOPS INNERFLOPS+6 */
845 /* #if GEOMETRY_J == 'Particle' */
846 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(f+j_coord_offsetA,f+j_coord_offsetB,
847 f+j_coord_offsetC,f+j_coord_offsetD,
849 /* #define INNERFLOPS INNERFLOPS+3 */
851 fjx{J} = _mm_add_ps(fjx{J},tx);
852 fjy{J} = _mm_add_ps(fjy{J},ty);
853 fjz{J} = _mm_add_ps(fjz{J},tz);
854 /* #define INNERFLOPS INNERFLOPS+3 */
859 /* #if 'exactcutoff' in INTERACTION_FLAGS[I][J] */
860 /* #if 0 ## This and next two lines is a hack to maintain indentation in template file */
865 /* ## End of check for the interaction being outside the cutoff */
868 /* ## End of loop over i-j interaction pairs */
870 /* #if GEOMETRY_J == 'Water3' */
871 gmx_mm_decrement_3rvec_4ptr_swizzle_ps(f+j_coord_offsetA,f+j_coord_offsetB,
872 f+j_coord_offsetC,f+j_coord_offsetD,
873 fjx0,fjy0,fjz0,fjx1,fjy1,fjz1,fjx2,fjy2,fjz2);
874 /* #define INNERFLOPS INNERFLOPS+9 */
875 /* #elif GEOMETRY_J == 'Water4' */
876 /* #if 0 in PARTICLES_J */
877 gmx_mm_decrement_4rvec_4ptr_swizzle_ps(f+j_coord_offsetA,f+j_coord_offsetB,
878 f+j_coord_offsetC,f+j_coord_offsetD,
879 fjx0,fjy0,fjz0,fjx1,fjy1,fjz1,
880 fjx2,fjy2,fjz2,fjx3,fjy3,fjz3);
881 /* #define INNERFLOPS INNERFLOPS+12 */
883 gmx_mm_decrement_3rvec_4ptr_swizzle_ps(f+j_coord_offsetA+DIM,f+j_coord_offsetB+DIM,
884 f+j_coord_offsetC+DIM,f+j_coord_offsetD+DIM,
885 fjx1,fjy1,fjz1,fjx2,fjy2,fjz2,fjx3,fjy3,fjz3);
886 /* #define INNERFLOPS INNERFLOPS+9 */
890 /* Inner loop uses {INNERFLOPS} flops */
895 /* End of innermost loop */
897 /* #if 'Force' in KERNEL_VF */
898 /* #if GEOMETRY_I == 'Particle' */
899 gmx_mm_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
900 f+i_coord_offset,fshift+i_shift_offset);
901 /* #define OUTERFLOPS OUTERFLOPS+6 */
902 /* #elif GEOMETRY_I == 'Water3' */
903 gmx_mm_update_iforce_3atom_swizzle_ps(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,
904 f+i_coord_offset,fshift+i_shift_offset);
905 /* #define OUTERFLOPS OUTERFLOPS+18 */
906 /* #elif GEOMETRY_I == 'Water4' */
907 /* #if 0 in PARTICLES_I */
908 gmx_mm_update_iforce_4atom_swizzle_ps(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
909 f+i_coord_offset,fshift+i_shift_offset);
910 /* #define OUTERFLOPS OUTERFLOPS+24 */
912 gmx_mm_update_iforce_3atom_swizzle_ps(fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
913 f+i_coord_offset+DIM,fshift+i_shift_offset);
914 /* #define OUTERFLOPS OUTERFLOPS+18 */
919 /* #if 'Potential' in KERNEL_VF */
921 /* Update potential energies */
922 /* #if KERNEL_ELEC != 'None' */
923 gmx_mm_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
924 /* #define OUTERFLOPS OUTERFLOPS+1 */
926 /* #if 'GeneralizedBorn' in KERNEL_ELEC */
927 gmx_mm_update_1pot_ps(vgbsum,kernel_data->energygrp_polarization+ggid);
928 /* #define OUTERFLOPS OUTERFLOPS+1 */
930 /* #if KERNEL_VDW != 'None' */
931 gmx_mm_update_1pot_ps(vvdwsum,kernel_data->energygrp_vdw+ggid);
932 /* #define OUTERFLOPS OUTERFLOPS+1 */
935 /* #if 'GeneralizedBorn' in KERNEL_ELEC and 'Force' in KERNEL_VF */
936 dvdasum = _mm_mul_ps(dvdasum, _mm_mul_ps(isai{I},isai{I}));
937 gmx_mm_update_1pot_ps(dvdasum,dvda+inr);
940 /* Increment number of inner iterations */
941 inneriter += j_index_end - j_index_start;
943 /* Outer loop uses {OUTERFLOPS} flops */
946 /* Increment number of outer iterations */
949 /* Update outer/inner flops */
950 /* ## NB: This is not important, it just affects the flopcount. However, since our preprocessor is */
951 /* ## primitive and replaces aggressively even in strings inside these directives, we need to */
952 /* ## assemble the main part of the name (containing KERNEL/ELEC/VDW) directly in the source. */
953 /* #if GEOMETRY_I == 'Water3' */
954 /* #define ISUFFIX '_W3' */
955 /* #elif GEOMETRY_I == 'Water4' */
956 /* #define ISUFFIX '_W4' */
958 /* #define ISUFFIX '' */
960 /* #if GEOMETRY_J == 'Water3' */
961 /* #define JSUFFIX 'W3' */
962 /* #elif GEOMETRY_J == 'Water4' */
963 /* #define JSUFFIX 'W4' */
965 /* #define JSUFFIX '' */
967 /* #if 'PotentialAndForce' in KERNEL_VF */
968 /* #define VFSUFFIX '_VF' */
969 /* #elif 'Potential' in KERNEL_VF */
970 /* #define VFSUFFIX '_V' */
972 /* #define VFSUFFIX '_F' */
975 /* #if KERNEL_ELEC != 'None' and KERNEL_VDW != 'None' */
976 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW{ISUFFIX}{JSUFFIX}{VFSUFFIX},outeriter*{OUTERFLOPS} + inneriter*{INNERFLOPS});
977 /* #elif KERNEL_ELEC != 'None' */
978 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC{ISUFFIX}{JSUFFIX}{VFSUFFIX},outeriter*{OUTERFLOPS} + inneriter*{INNERFLOPS});
980 inc_nrnb(nrnb,eNR_NBKERNEL_VDW{ISUFFIX}{JSUFFIX}{VFSUFFIX},outeriter*{OUTERFLOPS} + inneriter*{INNERFLOPS});