2 /* ## This file is part of the GROMACS molecular simulation package. */
4 /* ## Copyright (c) 2012, by the GROMACS development team, led by */
5 /* ## David van der Spoel, Berk Hess, Erik Lindahl, and including many */
6 /* ## others, as listed in the AUTHORS file in the top-level source */
7 /* ## directory and at http://www.gromacs.org. */
9 /* ## GROMACS is free software; you can redistribute it and/or */
10 /* ## modify it under the terms of the GNU Lesser General Public License */
11 /* ## as published by the Free Software Foundation; either version 2.1 */
12 /* ## of the License, or (at your option) any later version. */
14 /* ## GROMACS is distributed in the hope that it will be useful, */
15 /* ## but WITHOUT ANY WARRANTY; without even the implied warranty of */
16 /* ## MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU */
17 /* ## Lesser General Public License for more details. */
19 /* ## You should have received a copy of the GNU Lesser General Public */
20 /* ## License along with GROMACS; if not, see */
21 /* ## http://www.gnu.org/licenses, or write to the Free Software Foundation, */
22 /* ## Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA. */
24 /* ## If you want to redistribute modifications to GROMACS, please */
25 /* ## consider that scientific software is very special. Version */
26 /* ## control is crucial - bugs must be traceable. We will be happy to */
27 /* ## consider code for inclusion in the official distribution, but */
28 /* ## derived work must not be called official GROMACS. Details are found */
29 /* ## in the README & COPYING files - if they are missing, get the */
30 /* ## official version at http://www.gromacs.org. */
32 /* ## To help us fund GROMACS development, we humbly ask that you cite */
33 /* ## the research papers on the package. Check out http://www.gromacs.org. */
36 #error This file must be processed with the Gromacs pre-preprocessor
38 /* #if INCLUDE_HEADER */
45 #include "../nb_kernel.h"
46 #include "types/simple.h"
50 #include "gmx_math_x86_sse4_1_single.h"
51 #include "kernelutil_x86_sse4_1_single.h"
54 /* ## List of variables set by the generating script: */
56 /* ## Setttings that apply to the entire kernel: */
57 /* ## KERNEL_ELEC: String, choice for electrostatic interactions */
58 /* ## KERNEL_VDW: String, choice for van der Waals interactions */
59 /* ## KERNEL_NAME: String, name of this kernel */
60 /* ## KERNEL_VF: String telling if we calculate potential, force, or both */
61 /* ## GEOMETRY_I/GEOMETRY_J: String, name of each geometry, e.g. 'Water3' or '1Particle' */
63 /* ## Setttings that apply to particles in the outer (I) or inner (J) loops: */
64 /* ## PARTICLES_I[]/ Arrays with lists of i/j particles to use in kernel. It is */
65 /* ## PARTICLES_J[]: just [0] for particle geometry, but can be longer for water */
66 /* ## PARTICLES_ELEC_I[]/ Arrays with lists of i/j particle that have electrostatics */
67 /* ## PARTICLES_ELEC_J[]: interactions that should be calculated in this kernel. */
68 /* ## PARTICLES_VDW_I[]/ Arrays with the list of i/j particle that have VdW */
69 /* ## PARTICLES_VDW_J[]: interactions that should be calculated in this kernel. */
71 /* ## Setttings for pairs of interactions (e.g. 2nd i particle against 1st j particle) */
72 /* ## PAIRS_IJ[]: Array with (i,j) tuples of pairs for which interactions */
73 /* ## should be calculated in this kernel. Zero-charge particles */
74 /* ## do not have interactions with particles without vdw, and */
75 /* ## Vdw-only interactions are not evaluated in a no-vdw-kernel. */
76 /* ## INTERACTION_FLAGS[][]: 2D matrix, dimension e.g. 3*3 for water-water interactions. */
77 /* ## For each i-j pair, the element [I][J] is a list of strings */
78 /* ## defining properties/flags of this interaction. Examples */
79 /* ## include 'electrostatics'/'vdw' if that type of interaction */
80 /* ## should be evaluated, 'rsq'/'rinv'/'rinvsq' if those values */
81 /* ## are needed, and 'exactcutoff' or 'shift','switch' to */
82 /* ## decide if the force/potential should be modified. This way */
83 /* ## we only calculate values absolutely needed for each case. */
85 /* ## Calculate the size and offset for (merged/interleaved) table data */
88 * Gromacs nonbonded kernel: {KERNEL_NAME}
89 * Electrostatics interaction: {KERNEL_ELEC}
90 * VdW interaction: {KERNEL_VDW}
91 * Geometry: {GEOMETRY_I}-{GEOMETRY_J}
92 * Calculate force/pot: {KERNEL_VF}
96 (t_nblist * gmx_restrict nlist,
97 rvec * gmx_restrict xx,
98 rvec * gmx_restrict ff,
99 t_forcerec * gmx_restrict fr,
100 t_mdatoms * gmx_restrict mdatoms,
101 nb_kernel_data_t * gmx_restrict kernel_data,
102 t_nrnb * gmx_restrict nrnb)
104 /* ## Not all variables are used for all kernels, but any optimizing compiler fixes that, */
105 /* ## so there is no point in going to extremes to exclude variables that are not needed. */
106 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
107 * just 0 for non-waters.
108 * Suffixes A,B,C,D refer to j loop unrolling done with SSE, e.g. for the four different
109 * jnr indices corresponding to data put in the four positions in the SIMD register.
111 int i_shift_offset,i_coord_offset,outeriter,inneriter;
112 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
113 int jnrA,jnrB,jnrC,jnrD;
114 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
115 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
116 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
118 real *shiftvec,*fshift,*x,*f;
119 real *fjptrA,*fjptrB,*fjptrC,*fjptrD;
121 __m128 tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
122 /* #for I in PARTICLES_I */
124 __m128 ix{I},iy{I},iz{I},fix{I},fiy{I},fiz{I},iq{I},isai{I};
126 /* #for J in PARTICLES_J */
127 int vdwjidx{J}A,vdwjidx{J}B,vdwjidx{J}C,vdwjidx{J}D;
128 __m128 jx{J},jy{J},jz{J},fjx{J},fjy{J},fjz{J},jq{J},isaj{J};
130 /* #for I,J in PAIRS_IJ */
131 __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};
133 /* #if KERNEL_ELEC != 'None' */
134 __m128 velec,felec,velecsum,facel,crf,krf,krf2;
137 /* #if 'GeneralizedBorn' in KERNEL_ELEC */
139 __m128 vgb,fgb,vgbsum,dvdasum,gbscale,gbtabscale,isaprod,gbqqfactor,gbinvepsdiff,gbeps,dvdatmp;
140 __m128 minushalf = _mm_set1_ps(-0.5);
141 real *invsqrta,*dvda,*gbtab;
143 /* #if KERNEL_VDW != 'None' */
145 __m128 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
148 __m128 one_sixth = _mm_set1_ps(1.0/6.0);
149 __m128 one_twelfth = _mm_set1_ps(1.0/12.0);
151 /* #if 'Table' in KERNEL_ELEC or 'GeneralizedBorn' in KERNEL_ELEC or 'Table' in KERNEL_VDW */
153 __m128i ifour = _mm_set1_epi32(4);
154 __m128 rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF;
157 /* #if 'Ewald' in KERNEL_ELEC */
159 __m128 ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
162 /* #if 'PotentialSwitch' in [KERNEL_MOD_ELEC,KERNEL_MOD_VDW] */
163 __m128 rswitch,swV3,swV4,swV5,swF2,swF3,swF4,d,d2,sw,dsw;
164 real rswitch_scalar,d_scalar;
166 __m128 dummy_mask,cutoff_mask;
167 __m128 signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
168 __m128 one = _mm_set1_ps(1.0);
169 __m128 two = _mm_set1_ps(2.0);
175 jindex = nlist->jindex;
177 shiftidx = nlist->shift;
179 shiftvec = fr->shift_vec[0];
180 fshift = fr->fshift[0];
181 /* #if KERNEL_ELEC != 'None' */
182 facel = _mm_set1_ps(fr->epsfac);
183 charge = mdatoms->chargeA;
184 /* #if 'ReactionField' in KERNEL_ELEC */
185 krf = _mm_set1_ps(fr->ic->k_rf);
186 krf2 = _mm_set1_ps(fr->ic->k_rf*2.0);
187 crf = _mm_set1_ps(fr->ic->c_rf);
190 /* #if KERNEL_VDW != 'None' */
191 nvdwtype = fr->ntype;
193 vdwtype = mdatoms->typeA;
196 /* #if 'Table' in KERNEL_ELEC and 'Table' in KERNEL_VDW */
197 vftab = kernel_data->table_elec_vdw->data;
198 vftabscale = _mm_set1_ps(kernel_data->table_elec_vdw->scale);
199 /* #elif 'Table' in KERNEL_ELEC */
200 vftab = kernel_data->table_elec->data;
201 vftabscale = _mm_set1_ps(kernel_data->table_elec->scale);
202 /* #elif 'Table' in KERNEL_VDW */
203 vftab = kernel_data->table_vdw->data;
204 vftabscale = _mm_set1_ps(kernel_data->table_vdw->scale);
207 /* #if 'Ewald' in KERNEL_ELEC */
208 sh_ewald = _mm_set1_ps(fr->ic->sh_ewald);
209 /* #if KERNEL_VF=='Force' and KERNEL_MOD_ELEC!='PotentialSwitch' */
210 ewtab = fr->ic->tabq_coul_F;
211 ewtabscale = _mm_set1_ps(fr->ic->tabq_scale);
212 ewtabhalfspace = _mm_set1_ps(0.5/fr->ic->tabq_scale);
214 ewtab = fr->ic->tabq_coul_FDV0;
215 ewtabscale = _mm_set1_ps(fr->ic->tabq_scale);
216 ewtabhalfspace = _mm_set1_ps(0.5/fr->ic->tabq_scale);
220 /* #if KERNEL_ELEC=='GeneralizedBorn' */
221 invsqrta = fr->invsqrta;
223 gbtabscale = _mm_set1_ps(fr->gbtab.scale);
224 gbtab = fr->gbtab.data;
225 gbinvepsdiff = _mm_set1_ps((1.0/fr->epsilon_r) - (1.0/fr->gb_epsilon_solvent));
228 /* #if 'Water' in GEOMETRY_I */
229 /* Setup water-specific parameters */
230 inr = nlist->iinr[0];
231 /* #for I in PARTICLES_ELEC_I */
232 iq{I} = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+{I}]));
234 /* #for I in PARTICLES_VDW_I */
235 vdwioffset{I} = 2*nvdwtype*vdwtype[inr+{I}];
239 /* #if 'Water' in GEOMETRY_J */
240 /* #for J in PARTICLES_ELEC_J */
241 jq{J} = _mm_set1_ps(charge[inr+{J}]);
243 /* #for J in PARTICLES_VDW_J */
244 vdwjidx{J}A = 2*vdwtype[inr+{J}];
246 /* #for I,J in PAIRS_IJ */
247 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] */
248 qq{I}{J} = _mm_mul_ps(iq{I},jq{J});
250 /* #if 'vdw' in INTERACTION_FLAGS[I][J] */
251 c6_{I}{J} = _mm_set1_ps(vdwparam[vdwioffset{I}+vdwjidx{J}A]);
252 c12_{I}{J} = _mm_set1_ps(vdwparam[vdwioffset{I}+vdwjidx{J}A+1]);
257 /* #if KERNEL_MOD_ELEC!='None' or KERNEL_MOD_VDW!='None' */
258 /* #if KERNEL_ELEC!='None' */
259 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
260 rcutoff_scalar = fr->rcoulomb;
262 rcutoff_scalar = fr->rvdw;
264 rcutoff = _mm_set1_ps(rcutoff_scalar);
265 rcutoff2 = _mm_mul_ps(rcutoff,rcutoff);
268 /* #if KERNEL_MOD_VDW=='PotentialShift' */
269 sh_vdw_invrcut6 = _mm_set1_ps(fr->ic->sh_invrc6);
270 rvdw = _mm_set1_ps(fr->rvdw);
273 /* #if 'PotentialSwitch' in [KERNEL_MOD_ELEC,KERNEL_MOD_VDW] */
274 /* #if KERNEL_MOD_ELEC=='PotentialSwitch' */
275 rswitch_scalar = fr->rcoulomb_switch;
276 rswitch = _mm_set1_ps(rswitch_scalar);
278 rswitch_scalar = fr->rvdw_switch;
279 rswitch = _mm_set1_ps(rswitch_scalar);
281 /* Setup switch parameters */
282 d_scalar = rcutoff_scalar-rswitch_scalar;
283 d = _mm_set1_ps(d_scalar);
284 swV3 = _mm_set1_ps(-10.0/(d_scalar*d_scalar*d_scalar));
285 swV4 = _mm_set1_ps( 15.0/(d_scalar*d_scalar*d_scalar*d_scalar));
286 swV5 = _mm_set1_ps( -6.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
287 /* #if 'Force' in KERNEL_VF */
288 swF2 = _mm_set1_ps(-30.0/(d_scalar*d_scalar*d_scalar));
289 swF3 = _mm_set1_ps( 60.0/(d_scalar*d_scalar*d_scalar*d_scalar));
290 swF4 = _mm_set1_ps(-30.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
294 /* Avoid stupid compiler warnings */
295 jnrA = jnrB = jnrC = jnrD = 0;
301 /* ## Keep track of the floating point operations we issue for reporting! */
302 /* #define OUTERFLOPS 0 */
306 for(iidx=0;iidx<4*DIM;iidx++)
311 /* Start outer loop over neighborlists */
312 for(iidx=0; iidx<nri; iidx++)
314 /* Load shift vector for this list */
315 i_shift_offset = DIM*shiftidx[iidx];
317 /* Load limits for loop over neighbors */
318 j_index_start = jindex[iidx];
319 j_index_end = jindex[iidx+1];
321 /* Get outer coordinate index */
323 i_coord_offset = DIM*inr;
325 /* Load i particle coords and add shift vector */
326 /* #if GEOMETRY_I == 'Particle' */
327 gmx_mm_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
328 /* #elif GEOMETRY_I == 'Water3' */
329 gmx_mm_load_shift_and_3rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,
330 &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2);
331 /* #elif GEOMETRY_I == 'Water4' */
332 /* #if 0 in PARTICLES_I */
333 gmx_mm_load_shift_and_4rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,
334 &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
336 gmx_mm_load_shift_and_3rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset+DIM,
337 &ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
341 /* #if 'Force' in KERNEL_VF */
342 /* #for I in PARTICLES_I */
343 fix{I} = _mm_setzero_ps();
344 fiy{I} = _mm_setzero_ps();
345 fiz{I} = _mm_setzero_ps();
349 /* ## For water we already preloaded parameters at the start of the kernel */
350 /* #if not 'Water' in GEOMETRY_I */
351 /* Load parameters for i particles */
352 /* #for I in PARTICLES_ELEC_I */
353 iq{I} = _mm_mul_ps(facel,_mm_load1_ps(charge+inr+{I}));
354 /* #define OUTERFLOPS OUTERFLOPS+1 */
355 /* #if KERNEL_ELEC=='GeneralizedBorn' */
356 isai{I} = _mm_load1_ps(invsqrta+inr+{I});
359 /* #for I in PARTICLES_VDW_I */
360 vdwioffset{I} = 2*nvdwtype*vdwtype[inr+{I}];
364 /* #if 'Potential' in KERNEL_VF */
365 /* Reset potential sums */
366 /* #if KERNEL_ELEC != 'None' */
367 velecsum = _mm_setzero_ps();
369 /* #if 'GeneralizedBorn' in KERNEL_ELEC */
370 vgbsum = _mm_setzero_ps();
372 /* #if KERNEL_VDW != 'None' */
373 vvdwsum = _mm_setzero_ps();
376 /* #if 'GeneralizedBorn' in KERNEL_ELEC and 'Force' in KERNEL_VF */
377 dvdasum = _mm_setzero_ps();
380 /* #for ROUND in ['Loop','Epilogue'] */
382 /* #if ROUND =='Loop' */
383 /* Start inner kernel loop */
384 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
386 /* ## First round is normal loop (next statement resets indentation) */
393 /* ## Second round is epilogue */
395 /* #define INNERFLOPS 0 */
397 /* Get j neighbor index, and coordinate index */
398 /* #if ROUND =='Loop' */
404 jnrlistA = jjnr[jidx];
405 jnrlistB = jjnr[jidx+1];
406 jnrlistC = jjnr[jidx+2];
407 jnrlistD = jjnr[jidx+3];
408 /* Sign of each element will be negative for non-real atoms.
409 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
410 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
412 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
413 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
414 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
415 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
416 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
418 j_coord_offsetA = DIM*jnrA;
419 j_coord_offsetB = DIM*jnrB;
420 j_coord_offsetC = DIM*jnrC;
421 j_coord_offsetD = DIM*jnrD;
423 /* load j atom coordinates */
424 /* #if GEOMETRY_J == 'Particle' */
425 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
426 x+j_coord_offsetC,x+j_coord_offsetD,
428 /* #elif GEOMETRY_J == 'Water3' */
429 gmx_mm_load_3rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
430 x+j_coord_offsetC,x+j_coord_offsetD,
431 &jx0,&jy0,&jz0,&jx1,&jy1,&jz1,&jx2,&jy2,&jz2);
432 /* #elif GEOMETRY_J == 'Water4' */
433 /* #if 0 in PARTICLES_J */
434 gmx_mm_load_4rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
435 x+j_coord_offsetC,x+j_coord_offsetD,
436 &jx0,&jy0,&jz0,&jx1,&jy1,&jz1,&jx2,
437 &jy2,&jz2,&jx3,&jy3,&jz3);
439 gmx_mm_load_3rvec_4ptr_swizzle_ps(x+j_coord_offsetA+DIM,x+j_coord_offsetB+DIM,
440 x+j_coord_offsetC+DIM,x+j_coord_offsetD+DIM,
441 &jx1,&jy1,&jz1,&jx2,&jy2,&jz2,&jx3,&jy3,&jz3);
445 /* Calculate displacement vector */
446 /* #for I,J in PAIRS_IJ */
447 dx{I}{J} = _mm_sub_ps(ix{I},jx{J});
448 dy{I}{J} = _mm_sub_ps(iy{I},jy{J});
449 dz{I}{J} = _mm_sub_ps(iz{I},jz{J});
450 /* #define INNERFLOPS INNERFLOPS+3 */
453 /* Calculate squared distance and things based on it */
454 /* #for I,J in PAIRS_IJ */
455 rsq{I}{J} = gmx_mm_calc_rsq_ps(dx{I}{J},dy{I}{J},dz{I}{J});
456 /* #define INNERFLOPS INNERFLOPS+5 */
459 /* #for I,J in PAIRS_IJ */
460 /* #if 'rinv' in INTERACTION_FLAGS[I][J] */
461 rinv{I}{J} = gmx_mm_invsqrt_ps(rsq{I}{J});
462 /* #define INNERFLOPS INNERFLOPS+5 */
466 /* #for I,J in PAIRS_IJ */
467 /* #if 'rinvsq' in INTERACTION_FLAGS[I][J] */
468 /* # if 'rinv' not in INTERACTION_FLAGS[I][J] */
469 rinvsq{I}{J} = gmx_mm_inv_ps(rsq{I}{J});
470 /* #define INNERFLOPS INNERFLOPS+4 */
472 rinvsq{I}{J} = _mm_mul_ps(rinv{I}{J},rinv{I}{J});
473 /* #define INNERFLOPS INNERFLOPS+1 */
478 /* #if not 'Water' in GEOMETRY_J */
479 /* Load parameters for j particles */
480 /* #for J in PARTICLES_ELEC_J */
481 jq{J} = gmx_mm_load_4real_swizzle_ps(charge+jnrA+{J},charge+jnrB+{J},
482 charge+jnrC+{J},charge+jnrD+{J});
483 /* #if KERNEL_ELEC=='GeneralizedBorn' */
484 isaj{J} = gmx_mm_load_4real_swizzle_ps(invsqrta+jnrA+{J},invsqrta+jnrB+{J},
485 invsqrta+jnrC+{J},invsqrta+jnrD+{J});
488 /* #for J in PARTICLES_VDW_J */
489 vdwjidx{J}A = 2*vdwtype[jnrA+{J}];
490 vdwjidx{J}B = 2*vdwtype[jnrB+{J}];
491 vdwjidx{J}C = 2*vdwtype[jnrC+{J}];
492 vdwjidx{J}D = 2*vdwtype[jnrD+{J}];
496 /* #if 'Force' in KERNEL_VF and not 'Particle' in GEOMETRY_I */
497 /* #for J in PARTICLES_J */
498 fjx{J} = _mm_setzero_ps();
499 fjy{J} = _mm_setzero_ps();
500 fjz{J} = _mm_setzero_ps();
504 /* #for I,J in PAIRS_IJ */
506 /**************************
507 * CALCULATE INTERACTIONS *
508 **************************/
510 /* ## Note special check for TIP4P-TIP4P. Since we are cutting of all hydrogen interactions we also cut the LJ-only O-O interaction */
511 /* #if 'exactcutoff' in INTERACTION_FLAGS[I][J] or (GEOMETRY_I=='Water4' and GEOMETRY_J=='Water4' and 'exactcutoff' in INTERACTION_FLAGS[1][1]) */
512 /* ## We always calculate rinv/rinvsq above to enable pipelineing in compilers (performance tested on x86) */
513 if (gmx_mm_any_lt(rsq{I}{J},rcutoff2))
515 /* #if 0 ## this and the next two lines is a hack to maintain auto-indentation in template file */
518 /* #define INNERFLOPS INNERFLOPS+1 */
521 /* #if 'r' in INTERACTION_FLAGS[I][J] */
522 r{I}{J} = _mm_mul_ps(rsq{I}{J},rinv{I}{J});
523 /* #if ROUND == 'Epilogue' */
524 r{I}{J} = _mm_andnot_ps(dummy_mask,r{I}{J});
525 /* #define INNERFLOPS INNERFLOPS+1 */
527 /* #define INNERFLOPS INNERFLOPS+1 */
530 /* ## For water geometries we already loaded parameters at the start of the kernel */
531 /* #if not 'Water' in GEOMETRY_J */
532 /* Compute parameters for interactions between i and j atoms */
533 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] */
534 qq{I}{J} = _mm_mul_ps(iq{I},jq{J});
535 /* #define INNERFLOPS INNERFLOPS+1 */
537 /* #if 'vdw' in INTERACTION_FLAGS[I][J] */
538 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset{I}+vdwjidx{J}A,
539 vdwparam+vdwioffset{I}+vdwjidx{J}B,
540 vdwparam+vdwioffset{I}+vdwjidx{J}C,
541 vdwparam+vdwioffset{I}+vdwjidx{J}D,
542 &c6_{I}{J},&c12_{I}{J});
546 /* #if 'table' in INTERACTION_FLAGS[I][J] */
547 /* Calculate table index by multiplying r with table scale and truncate to integer */
548 rt = _mm_mul_ps(r{I}{J},vftabscale);
549 vfitab = _mm_cvttps_epi32(rt);
550 vfeps = _mm_sub_ps(rt,_mm_round_ps(rt, _MM_FROUND_FLOOR));
551 /* #define INNERFLOPS INNERFLOPS+4 */
552 /* #if 'Table' in KERNEL_ELEC and 'Table' in KERNEL_VDW */
553 /* ## 3 tables, 4 bytes per point: multiply index by 12 */
554 vfitab = _mm_slli_epi32(_mm_add_epi32(vfitab,_mm_slli_epi32(vfitab,1)),2);
555 /* #elif 'Table' in KERNEL_ELEC */
556 /* ## 1 table, 4 bytes per point: multiply index by 4 */
557 vfitab = _mm_slli_epi32(vfitab,2);
558 /* #elif 'Table' in KERNEL_VDW */
559 /* ## 2 tables, 4 bytes per point: multiply index by 8 */
560 vfitab = _mm_slli_epi32(vfitab,3);
564 /* ## ELECTROSTATIC INTERACTIONS */
565 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] */
567 /* #if KERNEL_ELEC=='Coulomb' */
569 /* COULOMB ELECTROSTATICS */
570 velec = _mm_mul_ps(qq{I}{J},rinv{I}{J});
571 /* #define INNERFLOPS INNERFLOPS+1 */
572 /* #if 'Force' in KERNEL_VF */
573 felec = _mm_mul_ps(velec,rinvsq{I}{J});
574 /* #define INNERFLOPS INNERFLOPS+2 */
577 /* #elif KERNEL_ELEC=='ReactionField' */
579 /* REACTION-FIELD ELECTROSTATICS */
580 /* #if 'Potential' in KERNEL_VF */
581 velec = _mm_mul_ps(qq{I}{J},_mm_sub_ps(_mm_add_ps(rinv{I}{J},_mm_mul_ps(krf,rsq{I}{J})),crf));
582 /* #define INNERFLOPS INNERFLOPS+4 */
584 /* #if 'Force' in KERNEL_VF */
585 felec = _mm_mul_ps(qq{I}{J},_mm_sub_ps(_mm_mul_ps(rinv{I}{J},rinvsq{I}{J}),krf2));
586 /* #define INNERFLOPS INNERFLOPS+3 */
589 /* #elif KERNEL_ELEC=='GeneralizedBorn' */
591 /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
592 isaprod = _mm_mul_ps(isai{I},isaj{J});
593 gbqqfactor = _mm_xor_ps(signbit,_mm_mul_ps(qq{I}{J},_mm_mul_ps(isaprod,gbinvepsdiff)));
594 gbscale = _mm_mul_ps(isaprod,gbtabscale);
595 /* #define INNERFLOPS INNERFLOPS+5 */
597 /* Calculate generalized born table index - this is a separate table from the normal one,
598 * but we use the same procedure by multiplying r with scale and truncating to integer.
600 rt = _mm_mul_ps(r{I}{J},gbscale);
601 gbitab = _mm_cvttps_epi32(rt);
602 gbeps = _mm_sub_ps(rt,_mm_round_ps(rt, _MM_FROUND_FLOOR));
603 gbitab = _mm_slli_epi32(gbitab,2);
604 Y = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,0) );
605 F = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,1) );
606 G = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,2) );
607 H = _mm_load_ps( gbtab + gmx_mm_extract_epi32(gbitab,3) );
608 _MM_TRANSPOSE4_PS(Y,F,G,H);
609 Heps = _mm_mul_ps(gbeps,H);
610 Fp = _mm_add_ps(F,_mm_mul_ps(gbeps,_mm_add_ps(G,Heps)));
611 VV = _mm_add_ps(Y,_mm_mul_ps(gbeps,Fp));
612 vgb = _mm_mul_ps(gbqqfactor,VV);
613 /* #define INNERFLOPS INNERFLOPS+10 */
615 /* #if 'Force' in KERNEL_VF */
616 FF = _mm_add_ps(Fp,_mm_mul_ps(gbeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
617 fgb = _mm_mul_ps(gbqqfactor,_mm_mul_ps(FF,gbscale));
618 dvdatmp = _mm_mul_ps(minushalf,_mm_add_ps(vgb,_mm_mul_ps(fgb,r{I}{J})));
619 /* #if ROUND == 'Epilogue' */
620 dvdatmp = _mm_andnot_ps(dummy_mask,dvdatmp);
622 dvdasum = _mm_add_ps(dvdasum,dvdatmp);
623 /* #if ROUND == 'Loop' */
629 /* 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. */
630 fjptrA = (jnrlistA>=0) ? dvda+jnrA : scratch;
631 fjptrB = (jnrlistB>=0) ? dvda+jnrB : scratch;
632 fjptrC = (jnrlistC>=0) ? dvda+jnrC : scratch;
633 fjptrD = (jnrlistD>=0) ? dvda+jnrD : scratch;
635 gmx_mm_increment_4real_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,_mm_mul_ps(dvdatmp,_mm_mul_ps(isaj{J},isaj{J})));
636 /* #define INNERFLOPS INNERFLOPS+13 */
638 velec = _mm_mul_ps(qq{I}{J},rinv{I}{J});
639 /* #define INNERFLOPS INNERFLOPS+1 */
640 /* #if 'Force' in KERNEL_VF */
641 felec = _mm_mul_ps(_mm_sub_ps(_mm_mul_ps(velec,rinv{I}{J}),fgb),rinv{I}{J});
642 /* #define INNERFLOPS INNERFLOPS+3 */
645 /* #elif KERNEL_ELEC=='Ewald' */
646 /* EWALD ELECTROSTATICS */
648 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
649 ewrt = _mm_mul_ps(r{I}{J},ewtabscale);
650 ewitab = _mm_cvttps_epi32(ewrt);
651 eweps = _mm_sub_ps(ewrt,_mm_round_ps(ewrt, _MM_FROUND_FLOOR));
652 /* #define INNERFLOPS INNERFLOPS+4 */
653 /* #if 'Potential' in KERNEL_VF or KERNEL_MOD_ELEC=='PotentialSwitch' */
654 ewitab = _mm_slli_epi32(ewitab,2);
655 ewtabF = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,0) );
656 ewtabD = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,1) );
657 ewtabV = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,2) );
658 ewtabFn = _mm_load_ps( ewtab + gmx_mm_extract_epi32(ewitab,3) );
659 _MM_TRANSPOSE4_PS(ewtabF,ewtabD,ewtabV,ewtabFn);
660 felec = _mm_add_ps(ewtabF,_mm_mul_ps(eweps,ewtabD));
661 /* #define INNERFLOPS INNERFLOPS+2 */
662 /* #if KERNEL_MOD_ELEC=='PotentialShift' */
663 velec = _mm_sub_ps(ewtabV,_mm_mul_ps(_mm_mul_ps(ewtabhalfspace,eweps),_mm_add_ps(ewtabF,felec)));
664 velec = _mm_mul_ps(qq{I}{J},_mm_sub_ps(_mm_sub_ps(rinv{I}{J},sh_ewald),velec));
665 /* #define INNERFLOPS INNERFLOPS+7 */
667 velec = _mm_sub_ps(ewtabV,_mm_mul_ps(_mm_mul_ps(ewtabhalfspace,eweps),_mm_add_ps(ewtabF,felec)));
668 velec = _mm_mul_ps(qq{I}{J},_mm_sub_ps(rinv{I}{J},velec));
669 /* #define INNERFLOPS INNERFLOPS+6 */
671 /* #if 'Force' in KERNEL_VF */
672 felec = _mm_mul_ps(_mm_mul_ps(qq{I}{J},rinv{I}{J}),_mm_sub_ps(rinvsq{I}{J},felec));
673 /* #define INNERFLOPS INNERFLOPS+3 */
675 /* #elif KERNEL_VF=='Force' */
676 gmx_mm_load_4pair_swizzle_ps(ewtab + gmx_mm_extract_epi32(ewitab,0),ewtab + gmx_mm_extract_epi32(ewitab,1),
677 ewtab + gmx_mm_extract_epi32(ewitab,2),ewtab + gmx_mm_extract_epi32(ewitab,3),
679 felec = _mm_add_ps(_mm_mul_ps( _mm_sub_ps(one,eweps),ewtabF),_mm_mul_ps(eweps,ewtabFn));
680 felec = _mm_mul_ps(_mm_mul_ps(qq{I}{J},rinv{I}{J}),_mm_sub_ps(rinvsq{I}{J},felec));
681 /* #define INNERFLOPS INNERFLOPS+7 */
684 /* #elif KERNEL_ELEC=='CubicSplineTable' */
686 /* CUBIC SPLINE TABLE ELECTROSTATICS */
687 Y = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,0) );
688 F = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,1) );
689 G = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,2) );
690 H = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,3) );
691 _MM_TRANSPOSE4_PS(Y,F,G,H);
692 Heps = _mm_mul_ps(vfeps,H);
693 Fp = _mm_add_ps(F,_mm_mul_ps(vfeps,_mm_add_ps(G,Heps)));
694 /* #define INNERFLOPS INNERFLOPS+4 */
695 /* #if 'Potential' in KERNEL_VF */
696 VV = _mm_add_ps(Y,_mm_mul_ps(vfeps,Fp));
697 velec = _mm_mul_ps(qq{I}{J},VV);
698 /* #define INNERFLOPS INNERFLOPS+3 */
700 /* #if 'Force' in KERNEL_VF */
701 FF = _mm_add_ps(Fp,_mm_mul_ps(vfeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
702 felec = _mm_xor_ps(signbit,_mm_mul_ps(_mm_mul_ps(qq{I}{J},FF),_mm_mul_ps(vftabscale,rinv{I}{J})));
703 /* #define INNERFLOPS INNERFLOPS+7 */
706 /* ## End of check for electrostatics interaction forms */
708 /* ## END OF ELECTROSTATIC INTERACTION CHECK FOR PAIR I-J */
710 /* #if 'vdw' in INTERACTION_FLAGS[I][J] */
712 /* #if KERNEL_VDW=='LennardJones' */
714 /* LENNARD-JONES DISPERSION/REPULSION */
716 rinvsix = _mm_mul_ps(_mm_mul_ps(rinvsq{I}{J},rinvsq{I}{J}),rinvsq{I}{J});
717 /* #define INNERFLOPS INNERFLOPS+2 */
718 /* #if 'Potential' in KERNEL_VF or KERNEL_MOD_VDW=='PotentialSwitch' */
719 vvdw6 = _mm_mul_ps(c6_{I}{J},rinvsix);
720 vvdw12 = _mm_mul_ps(c12_{I}{J},_mm_mul_ps(rinvsix,rinvsix));
721 /* #define INNERFLOPS INNERFLOPS+3 */
722 /* #if KERNEL_MOD_VDW=='PotentialShift' */
723 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) ,
724 _mm_mul_ps( _mm_sub_ps(vvdw6,_mm_mul_ps(c6_{I}{J},sh_vdw_invrcut6)),one_sixth));
725 /* #define INNERFLOPS INNERFLOPS+8 */
727 vvdw = _mm_sub_ps( _mm_mul_ps(vvdw12,one_twelfth) , _mm_mul_ps(vvdw6,one_sixth) );
728 /* #define INNERFLOPS INNERFLOPS+3 */
730 /* ## Check for force inside potential check, i.e. this means we already did the potential part */
731 /* #if 'Force' in KERNEL_VF */
732 fvdw = _mm_mul_ps(_mm_sub_ps(vvdw12,vvdw6),rinvsq{I}{J});
733 /* #define INNERFLOPS INNERFLOPS+2 */
735 /* #elif KERNEL_VF=='Force' */
736 /* ## Force-only LennardJones makes it possible to save 1 flop (they do add up...) */
737 fvdw = _mm_mul_ps(_mm_sub_ps(_mm_mul_ps(c12_{I}{J},rinvsix),c6_{I}{J}),_mm_mul_ps(rinvsix,rinvsq{I}{J}));
738 /* #define INNERFLOPS INNERFLOPS+4 */
741 /* #elif KERNEL_VDW=='CubicSplineTable' */
743 /* CUBIC SPLINE TABLE DISPERSION */
744 /* #if 'Table' in KERNEL_ELEC */
745 vfitab = _mm_add_epi32(vfitab,ifour);
747 Y = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,0) );
748 F = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,1) );
749 G = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,2) );
750 H = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,3) );
751 _MM_TRANSPOSE4_PS(Y,F,G,H);
752 Heps = _mm_mul_ps(vfeps,H);
753 Fp = _mm_add_ps(F,_mm_mul_ps(vfeps,_mm_add_ps(G,Heps)));
754 /* #define INNERFLOPS INNERFLOPS+4 */
755 /* #if 'Potential' in KERNEL_VF */
756 VV = _mm_add_ps(Y,_mm_mul_ps(vfeps,Fp));
757 vvdw6 = _mm_mul_ps(c6_{I}{J},VV);
758 /* #define INNERFLOPS INNERFLOPS+3 */
760 /* #if 'Force' in KERNEL_VF */
761 FF = _mm_add_ps(Fp,_mm_mul_ps(vfeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
762 fvdw6 = _mm_mul_ps(c6_{I}{J},FF);
763 /* #define INNERFLOPS INNERFLOPS+4 */
766 /* CUBIC SPLINE TABLE REPULSION */
767 vfitab = _mm_add_epi32(vfitab,ifour);
768 Y = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,0) );
769 F = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,1) );
770 G = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,2) );
771 H = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,3) );
772 _MM_TRANSPOSE4_PS(Y,F,G,H);
773 Heps = _mm_mul_ps(vfeps,H);
774 Fp = _mm_add_ps(F,_mm_mul_ps(vfeps,_mm_add_ps(G,Heps)));
775 /* #define INNERFLOPS INNERFLOPS+4 */
776 /* #if 'Potential' in KERNEL_VF */
777 VV = _mm_add_ps(Y,_mm_mul_ps(vfeps,Fp));
778 vvdw12 = _mm_mul_ps(c12_{I}{J},VV);
779 /* #define INNERFLOPS INNERFLOPS+3 */
781 /* #if 'Force' in KERNEL_VF */
782 FF = _mm_add_ps(Fp,_mm_mul_ps(vfeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
783 fvdw12 = _mm_mul_ps(c12_{I}{J},FF);
784 /* #define INNERFLOPS INNERFLOPS+5 */
786 /* #if 'Potential' in KERNEL_VF */
787 vvdw = _mm_add_ps(vvdw12,vvdw6);
788 /* #define INNERFLOPS INNERFLOPS+1 */
790 /* #if 'Force' in KERNEL_VF */
791 fvdw = _mm_xor_ps(signbit,_mm_mul_ps(_mm_add_ps(fvdw6,fvdw12),_mm_mul_ps(vftabscale,rinv{I}{J})));
792 /* #define INNERFLOPS INNERFLOPS+4 */
795 /* ## End of check for vdw interaction forms */
797 /* ## END OF VDW INTERACTION CHECK FOR PAIR I-J */
799 /* #if 'switch' in INTERACTION_FLAGS[I][J] */
800 d = _mm_sub_ps(r{I}{J},rswitch);
801 d = _mm_max_ps(d,_mm_setzero_ps());
802 d2 = _mm_mul_ps(d,d);
803 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)))))));
804 /* #define INNERFLOPS INNERFLOPS+10 */
806 /* #if 'Force' in KERNEL_VF */
807 dsw = _mm_mul_ps(d2,_mm_add_ps(swF2,_mm_mul_ps(d,_mm_add_ps(swF3,_mm_mul_ps(d,swF4)))));
808 /* #define INNERFLOPS INNERFLOPS+5 */
811 /* Evaluate switch function */
812 /* #if 'Force' in KERNEL_VF */
813 /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
814 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] and KERNEL_MOD_ELEC=='PotentialSwitch' */
815 felec = _mm_sub_ps( _mm_mul_ps(felec,sw) , _mm_mul_ps(rinv{I}{J},_mm_mul_ps(velec,dsw)) );
816 /* #define INNERFLOPS INNERFLOPS+4 */
818 /* #if 'vdw' in INTERACTION_FLAGS[I][J] and KERNEL_MOD_VDW=='PotentialSwitch' */
819 fvdw = _mm_sub_ps( _mm_mul_ps(fvdw,sw) , _mm_mul_ps(rinv{I}{J},_mm_mul_ps(vvdw,dsw)) );
820 /* #define INNERFLOPS INNERFLOPS+4 */
823 /* #if 'Potential' in KERNEL_VF */
824 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] and KERNEL_MOD_ELEC=='PotentialSwitch' */
825 velec = _mm_mul_ps(velec,sw);
826 /* #define INNERFLOPS INNERFLOPS+1 */
828 /* #if 'vdw' in INTERACTION_FLAGS[I][J] and KERNEL_MOD_VDW=='PotentialSwitch' */
829 vvdw = _mm_mul_ps(vvdw,sw);
830 /* #define INNERFLOPS INNERFLOPS+1 */
834 /* ## Note special check for TIP4P-TIP4P. Since we are cutting of all hydrogen interactions we also cut the LJ-only O-O interaction */
835 /* #if 'exactcutoff' in INTERACTION_FLAGS[I][J] or (GEOMETRY_I=='Water4' and GEOMETRY_J=='Water4' and 'exactcutoff' in INTERACTION_FLAGS[1][1]) */
836 cutoff_mask = _mm_cmplt_ps(rsq{I}{J},rcutoff2);
837 /* #define INNERFLOPS INNERFLOPS+1 */
840 /* #if 'Potential' in KERNEL_VF */
841 /* Update potential sum for this i atom from the interaction with this j atom. */
842 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] */
843 /* #if 'exactcutoff' in INTERACTION_FLAGS[I][J] */
844 velec = _mm_and_ps(velec,cutoff_mask);
845 /* #define INNERFLOPS INNERFLOPS+1 */
847 /* #if ROUND == 'Epilogue' */
848 velec = _mm_andnot_ps(dummy_mask,velec);
850 velecsum = _mm_add_ps(velecsum,velec);
851 /* #define INNERFLOPS INNERFLOPS+1 */
852 /* #if KERNEL_ELEC=='GeneralizedBorn' */
853 /* #if 'exactcutoff' in INTERACTION_FLAGS[I][J] */
854 vgb = _mm_and_ps(vgb,cutoff_mask);
855 /* #define INNERFLOPS INNERFLOPS+1 */
857 /* #if ROUND == 'Epilogue' */
858 vgb = _mm_andnot_ps(dummy_mask,vgb);
860 vgbsum = _mm_add_ps(vgbsum,vgb);
861 /* #define INNERFLOPS INNERFLOPS+1 */
864 /* #if 'vdw' in INTERACTION_FLAGS[I][J] */
865 /* ## Note special check for TIP4P-TIP4P. Since we are cutting of all hydrogen interactions we also cut the LJ-only O-O interaction */
866 /* #if 'exactcutoff' in INTERACTION_FLAGS[I][J] or (GEOMETRY_I=='Water4' and GEOMETRY_J=='Water4' and 'exactcutoff' in INTERACTION_FLAGS[1][1]) */
867 vvdw = _mm_and_ps(vvdw,cutoff_mask);
868 /* #define INNERFLOPS INNERFLOPS+1 */
870 /* #if ROUND == 'Epilogue' */
871 vvdw = _mm_andnot_ps(dummy_mask,vvdw);
873 vvdwsum = _mm_add_ps(vvdwsum,vvdw);
874 /* #define INNERFLOPS INNERFLOPS+1 */
878 /* #if 'Force' in KERNEL_VF */
880 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] and 'vdw' in INTERACTION_FLAGS[I][J] */
881 fscal = _mm_add_ps(felec,fvdw);
882 /* #define INNERFLOPS INNERFLOPS+1 */
883 /* #elif 'electrostatics' in INTERACTION_FLAGS[I][J] */
885 /* #elif 'vdw' in INTERACTION_FLAGS[I][J] */
889 /* ## Note special check for TIP4P-TIP4P. Since we are cutting of all hydrogen interactions we also cut the LJ-only O-O interaction */
890 /* #if 'exactcutoff' in INTERACTION_FLAGS[I][J] or (GEOMETRY_I=='Water4' and GEOMETRY_J=='Water4' and 'exactcutoff' in INTERACTION_FLAGS[1][1]) */
891 fscal = _mm_and_ps(fscal,cutoff_mask);
892 /* #define INNERFLOPS INNERFLOPS+1 */
895 /* #if ROUND == 'Epilogue' */
896 fscal = _mm_andnot_ps(dummy_mask,fscal);
899 /* Calculate temporary vectorial force */
900 tx = _mm_mul_ps(fscal,dx{I}{J});
901 ty = _mm_mul_ps(fscal,dy{I}{J});
902 tz = _mm_mul_ps(fscal,dz{I}{J});
904 /* Update vectorial force */
905 fix{I} = _mm_add_ps(fix{I},tx);
906 fiy{I} = _mm_add_ps(fiy{I},ty);
907 fiz{I} = _mm_add_ps(fiz{I},tz);
908 /* #define INNERFLOPS INNERFLOPS+6 */
910 /* #if GEOMETRY_I == 'Particle' */
911 /* #if ROUND == 'Loop' */
912 fjptrA = f+j_coord_offsetA;
913 fjptrB = f+j_coord_offsetB;
914 fjptrC = f+j_coord_offsetC;
915 fjptrD = f+j_coord_offsetD;
917 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
918 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
919 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
920 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
922 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
923 /* #define INNERFLOPS INNERFLOPS+3 */
925 fjx{J} = _mm_add_ps(fjx{J},tx);
926 fjy{J} = _mm_add_ps(fjy{J},ty);
927 fjz{J} = _mm_add_ps(fjz{J},tz);
928 /* #define INNERFLOPS INNERFLOPS+3 */
933 /* ## Note special check for TIP4P-TIP4P. Since we are cutting of all hydrogen interactions we also cut the LJ-only O-O interaction */
934 /* #if 'exactcutoff' in INTERACTION_FLAGS[I][J] or (GEOMETRY_I=='Water4' and GEOMETRY_J=='Water4' and 'exactcutoff' in INTERACTION_FLAGS[1][1]) */
935 /* #if 0 ## This and next two lines is a hack to maintain indentation in template file */
940 /* ## End of check for the interaction being outside the cutoff */
943 /* ## End of loop over i-j interaction pairs */
945 /* #if GEOMETRY_I != 'Particle' */
946 /* #if ROUND == 'Loop' */
947 fjptrA = f+j_coord_offsetA;
948 fjptrB = f+j_coord_offsetB;
949 fjptrC = f+j_coord_offsetC;
950 fjptrD = f+j_coord_offsetD;
952 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
953 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
954 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
955 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
959 /* #if 'Water' in GEOMETRY_I and GEOMETRY_J == 'Particle' */
960 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjx0,fjy0,fjz0);
961 /* #elif GEOMETRY_J == 'Water3' */
962 gmx_mm_decrement_3rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,
963 fjx0,fjy0,fjz0,fjx1,fjy1,fjz1,fjx2,fjy2,fjz2);
964 /* #define INNERFLOPS INNERFLOPS+9 */
965 /* #elif GEOMETRY_J == 'Water4' */
966 /* #if 0 in PARTICLES_J */
967 gmx_mm_decrement_4rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,
968 fjx0,fjy0,fjz0,fjx1,fjy1,fjz1,
969 fjx2,fjy2,fjz2,fjx3,fjy3,fjz3);
970 /* #define INNERFLOPS INNERFLOPS+12 */
972 gmx_mm_decrement_3rvec_4ptr_swizzle_ps(fjptrA+DIM,fjptrB+DIM,fjptrC+DIM,fjptrD+DIM,
973 fjx1,fjy1,fjz1,fjx2,fjy2,fjz2,fjx3,fjy3,fjz3);
974 /* #define INNERFLOPS INNERFLOPS+9 */
978 /* Inner loop uses {INNERFLOPS} flops */
983 /* End of innermost loop */
985 /* #if 'Force' in KERNEL_VF */
986 /* #if GEOMETRY_I == 'Particle' */
987 gmx_mm_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
988 f+i_coord_offset,fshift+i_shift_offset);
989 /* #define OUTERFLOPS OUTERFLOPS+6 */
990 /* #elif GEOMETRY_I == 'Water3' */
991 gmx_mm_update_iforce_3atom_swizzle_ps(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,
992 f+i_coord_offset,fshift+i_shift_offset);
993 /* #define OUTERFLOPS OUTERFLOPS+18 */
994 /* #elif GEOMETRY_I == 'Water4' */
995 /* #if 0 in PARTICLES_I */
996 gmx_mm_update_iforce_4atom_swizzle_ps(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
997 f+i_coord_offset,fshift+i_shift_offset);
998 /* #define OUTERFLOPS OUTERFLOPS+24 */
1000 gmx_mm_update_iforce_3atom_swizzle_ps(fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
1001 f+i_coord_offset+DIM,fshift+i_shift_offset);
1002 /* #define OUTERFLOPS OUTERFLOPS+18 */
1007 /* #if 'Potential' in KERNEL_VF */
1009 /* Update potential energies */
1010 /* #if KERNEL_ELEC != 'None' */
1011 gmx_mm_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
1012 /* #define OUTERFLOPS OUTERFLOPS+1 */
1014 /* #if 'GeneralizedBorn' in KERNEL_ELEC */
1015 gmx_mm_update_1pot_ps(vgbsum,kernel_data->energygrp_polarization+ggid);
1016 /* #define OUTERFLOPS OUTERFLOPS+1 */
1018 /* #if KERNEL_VDW != 'None' */
1019 gmx_mm_update_1pot_ps(vvdwsum,kernel_data->energygrp_vdw+ggid);
1020 /* #define OUTERFLOPS OUTERFLOPS+1 */
1023 /* #if 'GeneralizedBorn' in KERNEL_ELEC and 'Force' in KERNEL_VF */
1024 dvdasum = _mm_mul_ps(dvdasum, _mm_mul_ps(isai{I},isai{I}));
1025 gmx_mm_update_1pot_ps(dvdasum,dvda+inr);
1028 /* Increment number of inner iterations */
1029 inneriter += j_index_end - j_index_start;
1031 /* Outer loop uses {OUTERFLOPS} flops */
1034 /* Increment number of outer iterations */
1037 /* Update outer/inner flops */
1038 /* ## NB: This is not important, it just affects the flopcount. However, since our preprocessor is */
1039 /* ## primitive and replaces aggressively even in strings inside these directives, we need to */
1040 /* ## assemble the main part of the name (containing KERNEL/ELEC/VDW) directly in the source. */
1041 /* #if GEOMETRY_I == 'Water3' */
1042 /* #define ISUFFIX '_W3' */
1043 /* #elif GEOMETRY_I == 'Water4' */
1044 /* #define ISUFFIX '_W4' */
1046 /* #define ISUFFIX '' */
1048 /* #if GEOMETRY_J == 'Water3' */
1049 /* #define JSUFFIX 'W3' */
1050 /* #elif GEOMETRY_J == 'Water4' */
1051 /* #define JSUFFIX 'W4' */
1053 /* #define JSUFFIX '' */
1055 /* #if 'PotentialAndForce' in KERNEL_VF */
1056 /* #define VFSUFFIX '_VF' */
1057 /* #elif 'Potential' in KERNEL_VF */
1058 /* #define VFSUFFIX '_V' */
1060 /* #define VFSUFFIX '_F' */
1063 /* #if KERNEL_ELEC != 'None' and KERNEL_VDW != 'None' */
1064 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW{ISUFFIX}{JSUFFIX}{VFSUFFIX},outeriter*{OUTERFLOPS} + inneriter*{INNERFLOPS});
1065 /* #elif KERNEL_ELEC != 'None' */
1066 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC{ISUFFIX}{JSUFFIX}{VFSUFFIX},outeriter*{OUTERFLOPS} + inneriter*{INNERFLOPS});
1068 inc_nrnb(nrnb,eNR_NBKERNEL_VDW{ISUFFIX}{JSUFFIX}{VFSUFFIX},outeriter*{OUTERFLOPS} + inneriter*{INNERFLOPS});