8c97de1c7ea26d2c2e98187ac098bb5ef5185ca4
[alexxy/gromacs.git] / src / gromacs / mdlib / nbnxn_kernels / simd_2xnn / nbnxn_kernel_simd_2xnn_inner.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2013,2014, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35
36 /* This is the innermost loop contents for the 4 x N atom SIMD kernel.
37  * This flavor of the kernel duplicates the data for N j-particles in
38  * 2xN wide SIMD registers to do operate on 2 i-particles at once.
39  * This leads to 4/2=2 sets of most instructions. Therefore we call
40  * this kernel 2x(N+N) = 2xnn
41  *
42  * This 2xnn kernel is basically the 4xn equivalent with half the registers
43  * and instructions removed.
44  *
45  * An alternative would be to load to different cluster of N j-particles
46  * into SIMD registers, giving a 4x(N+N) kernel. This doubles the amount
47  * of instructions, which could lead to better scheduling. But we actually
48  * observed worse scheduling for the AVX-256 4x8 normal analytical PME
49  * kernel, which has a lower pair throughput than 2x(4+4) with gcc 4.7.
50  * It could be worth trying this option, but it takes some more effort.
51  * This 2xnn kernel is basically the 4xn equivalent with
52  */
53
54
55 /* When calculating RF or Ewald interactions we calculate the electrostatic/LJ
56  * forces on excluded atom pairs here in the non-bonded loops.
57  * But when energies and/or virial is required we calculate them
58  * separately to as then it is easier to separate the energy and virial
59  * contributions.
60  */
61 #if defined CHECK_EXCLS && (defined CALC_COULOMB || defined LJ_EWALD_GEOM)
62 #define EXCL_FORCES
63 #endif
64
65 /* Without exclusions and energies we only need to mask the cut-off,
66  * this can be faster with blendv.
67  */
68 #if !(defined CHECK_EXCLS || defined CALC_ENERGIES) && defined GMX_SIMD_HAVE_BLENDV && !defined COUNT_PAIRS
69 /* With RF and tabulated Coulomb we replace cmp+and with sub+blendv.
70  * With gcc this is slower, except for RF on Sandy Bridge.
71  * Tested with gcc 4.6.2, 4.6.3 and 4.7.1.
72  */
73 #if (defined CALC_COUL_RF || defined CALC_COUL_TAB) && (!defined __GNUC__ || (defined CALC_COUL_RF && defined GMX_SIMD_X86_AVX_256_OR_HIGHER))
74 #define NBNXN_CUTOFF_USE_BLENDV
75 #endif
76 /* With analytical Ewald we replace cmp+and+and with sub+blendv+blendv.
77  * This is only faster with icc on Sandy Bridge (PS kernel slower than gcc 4.7).
78  * Tested with icc 13.
79  */
80 #if defined CALC_COUL_EWALD && defined __INTEL_COMPILER && defined GMX_SIMD_X86_AVX_256_OR_HIGHER
81 #define NBNXN_CUTOFF_USE_BLENDV
82 #endif
83 #endif
84
85 {
86     int        cj, aj, ajx, ajy, ajz;
87
88 #ifdef ENERGY_GROUPS
89     /* Energy group indices for two atoms packed into one int */
90     int        egp_jj[UNROLLJ/2];
91 #endif
92
93 #ifdef CHECK_EXCLS
94     /* Interaction (non-exclusion) mask of all 1's or 0's */
95     gmx_simd_bool_t  interact_S0;
96     gmx_simd_bool_t  interact_S2;
97 #endif
98
99     gmx_simd_real_t  jx_S, jy_S, jz_S;
100     gmx_simd_real_t  dx_S0, dy_S0, dz_S0;
101     gmx_simd_real_t  dx_S2, dy_S2, dz_S2;
102     gmx_simd_real_t  tx_S0, ty_S0, tz_S0;
103     gmx_simd_real_t  tx_S2, ty_S2, tz_S2;
104     gmx_simd_real_t  rsq_S0, rinv_S0, rinvsq_S0;
105     gmx_simd_real_t  rsq_S2, rinv_S2, rinvsq_S2;
106 #ifndef NBNXN_CUTOFF_USE_BLENDV
107     /* wco: within cut-off, mask of all 1's or 0's */
108     gmx_simd_bool_t  wco_S0;
109     gmx_simd_bool_t  wco_S2;
110 #endif
111 #ifdef VDW_CUTOFF_CHECK
112     gmx_simd_bool_t  wco_vdw_S0;
113 #ifndef HALF_LJ
114     gmx_simd_bool_t  wco_vdw_S2;
115 #endif
116 #endif
117
118 #if (defined CALC_COULOMB && defined CALC_COUL_TAB) || defined LJ_FORCE_SWITCH || defined LJ_POT_SWITCH
119     gmx_simd_real_t r_S0;
120     gmx_simd_real_t r_S2;
121 #endif
122
123 #if defined LJ_FORCE_SWITCH || defined LJ_POT_SWITCH
124     gmx_simd_real_t  rsw_S0, rsw2_S0, rsw2_r_S0;
125 #ifndef HALF_LJ
126     gmx_simd_real_t  rsw_S2, rsw2_S2, rsw2_r_S2;
127 #endif
128 #endif
129
130 #ifdef CALC_COULOMB
131 #ifdef CHECK_EXCLS
132     /* 1/r masked with the interaction mask */
133     gmx_simd_real_t  rinv_ex_S0;
134     gmx_simd_real_t  rinv_ex_S2;
135 #endif
136     gmx_simd_real_t  jq_S;
137     gmx_simd_real_t  qq_S0;
138     gmx_simd_real_t  qq_S2;
139 #ifdef CALC_COUL_TAB
140     /* The force (PME mesh force) we need to subtract from 1/r^2 */
141     gmx_simd_real_t  fsub_S0;
142     gmx_simd_real_t  fsub_S2;
143 #endif
144 #ifdef CALC_COUL_EWALD
145     gmx_simd_real_t  brsq_S0, brsq_S2;
146     gmx_simd_real_t  ewcorr_S0, ewcorr_S2;
147 #endif
148
149     /* frcoul = (1/r - fsub)*r */
150     gmx_simd_real_t  frcoul_S0;
151     gmx_simd_real_t  frcoul_S2;
152 #ifdef CALC_COUL_TAB
153     /* For tables: r, rs=r/sp, rf=floor(rs), frac=rs-rf */
154     gmx_simd_real_t         rs_S0, rf_S0, frac_S0;
155     gmx_simd_real_t         rs_S2, rf_S2, frac_S2;
156     /* Table index: rs truncated to an int */
157     gmx_simd_int32_t        ti_S0, ti_S2;
158     /* Linear force table values */
159     gmx_simd_real_t         ctab0_S0, ctab1_S0;
160     gmx_simd_real_t         ctab0_S2, ctab1_S2;
161 #ifdef CALC_ENERGIES
162     /* Quadratic energy table value */
163     gmx_simd_real_t  ctabv_S0;
164     gmx_simd_real_t  ctabv_S2;
165 #endif
166 #endif
167 #if defined CALC_ENERGIES && (defined CALC_COUL_EWALD || defined CALC_COUL_TAB)
168     /* The potential (PME mesh) we need to subtract from 1/r */
169     gmx_simd_real_t  vc_sub_S0;
170     gmx_simd_real_t  vc_sub_S2;
171 #endif
172 #ifdef CALC_ENERGIES
173     /* Electrostatic potential */
174     gmx_simd_real_t  vcoul_S0;
175     gmx_simd_real_t  vcoul_S2;
176 #endif
177 #endif
178     /* The force times 1/r */
179     gmx_simd_real_t  fscal_S0;
180     gmx_simd_real_t  fscal_S2;
181
182 #ifdef CALC_LJ
183 #ifdef LJ_COMB_LB
184     /* LJ sigma_j/2 and sqrt(epsilon_j) */
185     gmx_simd_real_t  hsig_j_S, seps_j_S;
186     /* LJ sigma_ij and epsilon_ij */
187     gmx_simd_real_t  sig_S0, eps_S0;
188 #ifndef HALF_LJ
189     gmx_simd_real_t  sig_S2, eps_S2;
190 #endif
191 #ifdef CALC_ENERGIES
192     gmx_simd_real_t  sig2_S0, sig6_S0;
193 #ifndef HALF_LJ
194     gmx_simd_real_t  sig2_S2, sig6_S2;
195 #endif
196 #endif /* LJ_COMB_LB */
197 #endif /* CALC_LJ */
198
199 #ifdef LJ_COMB_GEOM
200     gmx_simd_real_t  c6s_j_S, c12s_j_S;
201 #endif
202
203 #if defined LJ_COMB_GEOM || defined LJ_COMB_LB || defined LJ_EWALD_GEOM
204     /* Index for loading LJ parameters, complicated when interleaving */
205     int         aj2;
206 #endif
207
208 #ifndef FIX_LJ_C
209     /* LJ C6 and C12 parameters, used with geometric comb. rule */
210     gmx_simd_real_t  c6_S0, c12_S0;
211 #ifndef HALF_LJ
212     gmx_simd_real_t  c6_S2, c12_S2;
213 #endif
214 #endif
215
216     /* Intermediate variables for LJ calculation */
217 #ifndef LJ_COMB_LB
218     gmx_simd_real_t  rinvsix_S0;
219 #ifndef HALF_LJ
220     gmx_simd_real_t  rinvsix_S2;
221 #endif
222 #endif
223 #ifdef LJ_COMB_LB
224     gmx_simd_real_t  sir_S0, sir2_S0, sir6_S0;
225 #ifndef HALF_LJ
226     gmx_simd_real_t  sir_S2, sir2_S2, sir6_S2;
227 #endif
228 #endif
229
230     gmx_simd_real_t  FrLJ6_S0, FrLJ12_S0, frLJ_S0;
231 #ifndef HALF_LJ
232     gmx_simd_real_t  FrLJ6_S2, FrLJ12_S2, frLJ_S2;
233 #endif
234 #if defined CALC_ENERGIES || defined LJ_POT_SWITCH
235     gmx_simd_real_t  VLJ6_S0, VLJ12_S0, VLJ_S0;
236 #ifndef HALF_LJ
237     gmx_simd_real_t  VLJ6_S2, VLJ12_S2, VLJ_S2;
238 #endif
239 #endif
240 #endif /* CALC_LJ */
241
242     gmx_mm_hpr fjx_S, fjy_S, fjz_S;
243
244     /* j-cluster index */
245     cj            = l_cj[cjind].cj;
246
247     /* Atom indices (of the first atom in the cluster) */
248     aj            = cj*UNROLLJ;
249 #if defined CALC_LJ && (defined LJ_COMB_GEOM || defined LJ_COMB_LB || defined LJ_EWALD_GEOM)
250 #if UNROLLJ == STRIDE
251     aj2           = aj*2;
252 #else
253     aj2           = (cj>>1)*2*STRIDE + (cj & 1)*UNROLLJ;
254 #endif
255 #endif
256 #if UNROLLJ == STRIDE
257     ajx           = aj*DIM;
258 #else
259     ajx           = (cj>>1)*DIM*STRIDE + (cj & 1)*UNROLLJ;
260 #endif
261     ajy           = ajx + STRIDE;
262     ajz           = ajy + STRIDE;
263
264 #ifdef CHECK_EXCLS
265     gmx_load_simd_2xnn_interactions(l_cj[cjind].excl,
266                                     filter_S0, filter_S2,
267                                     &interact_S0, &interact_S2);
268 #endif /* CHECK_EXCLS */
269
270     /* load j atom coordinates */
271     gmx_loaddh_pr(&jx_S, x+ajx);
272     gmx_loaddh_pr(&jy_S, x+ajy);
273     gmx_loaddh_pr(&jz_S, x+ajz);
274
275     /* Calculate distance */
276     dx_S0       = gmx_simd_sub_r(ix_S0, jx_S);
277     dy_S0       = gmx_simd_sub_r(iy_S0, jy_S);
278     dz_S0       = gmx_simd_sub_r(iz_S0, jz_S);
279     dx_S2       = gmx_simd_sub_r(ix_S2, jx_S);
280     dy_S2       = gmx_simd_sub_r(iy_S2, jy_S);
281     dz_S2       = gmx_simd_sub_r(iz_S2, jz_S);
282
283     /* rsq = dx*dx+dy*dy+dz*dz */
284     rsq_S0      = gmx_simd_calc_rsq_r(dx_S0, dy_S0, dz_S0);
285     rsq_S2      = gmx_simd_calc_rsq_r(dx_S2, dy_S2, dz_S2);
286
287 #ifndef NBNXN_CUTOFF_USE_BLENDV
288     wco_S0      = gmx_simd_cmplt_r(rsq_S0, rc2_S);
289     wco_S2      = gmx_simd_cmplt_r(rsq_S2, rc2_S);
290 #endif
291
292 #ifdef CHECK_EXCLS
293 #ifdef EXCL_FORCES
294     /* Only remove the (sub-)diagonal to avoid double counting */
295 #if UNROLLJ == UNROLLI
296     if (cj == ci_sh)
297     {
298         wco_S0  = gmx_simd_and_b(wco_S0, diagonal_mask_S0);
299         wco_S2  = gmx_simd_and_b(wco_S2, diagonal_mask_S2);
300     }
301 #else
302 #if UNROLLJ == 2*UNROLLI
303     if (cj*2 == ci_sh)
304     {
305         wco_S0  = gmx_simd_and_b(wco_S0, diagonal_mask0_S0);
306         wco_S2  = gmx_simd_and_b(wco_S2, diagonal_mask0_S2);
307     }
308     else if (cj*2 + 1 == ci_sh)
309     {
310         wco_S0  = gmx_simd_and_b(wco_S0, diagonal_mask1_S0);
311         wco_S2  = gmx_simd_and_b(wco_S2, diagonal_mask1_S2);
312     }
313 #else
314 #error "only UNROLLJ == UNROLLI*(1 or 2) currently supported in 2xnn kernels"
315 #endif
316 #endif
317 #else /* EXCL_FORCES */
318       /* No exclusion forces: remove all excluded atom pairs from the list */
319     wco_S0      = gmx_simd_and_b(wco_S0, interact_S0);
320     wco_S2      = gmx_simd_and_b(wco_S2, interact_S2);
321 #endif
322 #endif
323
324 #ifdef COUNT_PAIRS
325     {
326         int  i, j;
327         real tmpa[2*GMX_SIMD_REAL_WIDTH], *tmp;
328         tmp = gmx_simd_align_r(tmpa);
329         for (i = 0; i < UNROLLI; i += 2)
330         {
331             gmx_simd_store_r(tmp, i == 0 ? wco_S0 : wco_S2);
332             for (j = 0; j < 2*UNROLLJ; j++)
333             {
334                 if (!(tmp[j] == 0))
335                 {
336                     npair++;
337                 }
338             }
339         }
340     }
341 #endif
342
343 #ifdef CHECK_EXCLS
344     /* For excluded pairs add a small number to avoid r^-6 = NaN */
345     rsq_S0      = gmx_masknot_add_pr(interact_S0, rsq_S0, avoid_sing_S);
346     rsq_S2      = gmx_masknot_add_pr(interact_S2, rsq_S2, avoid_sing_S);
347 #endif
348
349     /* Calculate 1/r */
350     rinv_S0     = gmx_simd_invsqrt_r(rsq_S0);
351     rinv_S2     = gmx_simd_invsqrt_r(rsq_S2);
352
353 #ifdef CALC_COULOMB
354     /* Load parameters for j atom */
355     gmx_loaddh_pr(&jq_S, q+aj);
356     qq_S0       = gmx_simd_mul_r(iq_S0, jq_S);
357     qq_S2       = gmx_simd_mul_r(iq_S2, jq_S);
358 #endif
359
360 #ifdef CALC_LJ
361
362 #if !defined LJ_COMB_GEOM && !defined LJ_COMB_LB && !defined FIX_LJ_C
363     load_lj_pair_params2(nbfp0, nbfp1, type, aj, &c6_S0, &c12_S0);
364 #ifndef HALF_LJ
365     load_lj_pair_params2(nbfp2, nbfp3, type, aj, &c6_S2, &c12_S2);
366 #endif
367 #endif /* not defined any LJ rule */
368
369 #ifdef LJ_COMB_GEOM
370     gmx_loaddh_pr(&c6s_j_S,  ljc+aj2+0);
371     gmx_loaddh_pr(&c12s_j_S, ljc+aj2+STRIDE);
372     c6_S0       = gmx_simd_mul_r(c6s_S0, c6s_j_S );
373 #ifndef HALF_LJ
374     c6_S2       = gmx_simd_mul_r(c6s_S2, c6s_j_S );
375 #endif
376     c12_S0      = gmx_simd_mul_r(c12s_S0, c12s_j_S);
377 #ifndef HALF_LJ
378     c12_S2      = gmx_simd_mul_r(c12s_S2, c12s_j_S);
379 #endif
380 #endif /* LJ_COMB_GEOM */
381
382 #ifdef LJ_COMB_LB
383     gmx_loaddh_pr(&hsig_j_S, ljc+aj2+0);
384     gmx_loaddh_pr(&seps_j_S, ljc+aj2+STRIDE);
385
386     sig_S0      = gmx_simd_add_r(hsig_i_S0, hsig_j_S);
387     eps_S0      = gmx_simd_mul_r(seps_i_S0, seps_j_S);
388 #ifndef HALF_LJ
389     sig_S2      = gmx_simd_add_r(hsig_i_S2, hsig_j_S);
390     eps_S2      = gmx_simd_mul_r(seps_i_S2, seps_j_S);
391 #endif
392 #endif /* LJ_COMB_LB */
393
394 #endif /* CALC_LJ */
395
396 #ifndef NBNXN_CUTOFF_USE_BLENDV
397     rinv_S0     = gmx_simd_blendzero_r(rinv_S0, wco_S0);
398     rinv_S2     = gmx_simd_blendzero_r(rinv_S2, wco_S2);
399 #else
400     /* We only need to mask for the cut-off: blendv is faster */
401     rinv_S0     = gmx_simd_blendv_r(rinv_S0, zero_S, gmx_simd_sub_r(rc2_S, rsq_S0));
402     rinv_S2     = gmx_simd_blendv_r(rinv_S2, zero_S, gmx_simd_sub_r(rc2_S, rsq_S2));
403 #endif
404
405     rinvsq_S0   = gmx_simd_mul_r(rinv_S0, rinv_S0);
406     rinvsq_S2   = gmx_simd_mul_r(rinv_S2, rinv_S2);
407
408 #ifdef CALC_COULOMB
409     /* Note that here we calculate force*r, not the usual force/r.
410      * This allows avoiding masking the reaction-field contribution,
411      * as frcoul is later multiplied by rinvsq which has been
412      * masked with the cut-off check.
413      */
414
415 #ifdef EXCL_FORCES
416     /* Only add 1/r for non-excluded atom pairs */
417     rinv_ex_S0  = gmx_simd_blendzero_r(rinv_S0, interact_S0);
418     rinv_ex_S2  = gmx_simd_blendzero_r(rinv_S2, interact_S2);
419 #else
420     /* No exclusion forces, we always need 1/r */
421 #define     rinv_ex_S0    rinv_S0
422 #define     rinv_ex_S2    rinv_S2
423 #endif
424
425 #ifdef CALC_COUL_RF
426     /* Electrostatic interactions */
427     frcoul_S0   = gmx_simd_mul_r(qq_S0, gmx_simd_fmadd_r(rsq_S0, mrc_3_S, rinv_ex_S0));
428     frcoul_S2   = gmx_simd_mul_r(qq_S2, gmx_simd_fmadd_r(rsq_S2, mrc_3_S, rinv_ex_S2));
429
430 #ifdef CALC_ENERGIES
431     vcoul_S0    = gmx_simd_mul_r(qq_S0, gmx_simd_add_r(rinv_ex_S0, gmx_simd_add_r(gmx_simd_mul_r(rsq_S0, hrc_3_S), moh_rc_S)));
432     vcoul_S2    = gmx_simd_mul_r(qq_S2, gmx_simd_add_r(rinv_ex_S2, gmx_simd_add_r(gmx_simd_mul_r(rsq_S2, hrc_3_S), moh_rc_S)));
433 #endif
434 #endif
435
436 #ifdef CALC_COUL_EWALD
437     /* We need to mask (or limit) rsq for the cut-off,
438      * as large distances can cause an overflow in gmx_pmecorrF/V.
439      */
440 #ifndef NBNXN_CUTOFF_USE_BLENDV
441     brsq_S0     = gmx_simd_mul_r(beta2_S, gmx_simd_blendzero_r(rsq_S0, wco_S0));
442     brsq_S2     = gmx_simd_mul_r(beta2_S, gmx_simd_blendzero_r(rsq_S2, wco_S2));
443 #else
444     /* Strangely, putting mul on a separate line is slower (icc 13) */
445     brsq_S0     = gmx_simd_mul_r(beta2_S, gmx_simd_blendv_r(rsq_S0, zero_S, gmx_simd_sub_r(rc2_S, rsq_S0)));
446     brsq_S2     = gmx_simd_mul_r(beta2_S, gmx_simd_blendv_r(rsq_S2, zero_S, gmx_simd_sub_r(rc2_S, rsq_S2)));
447 #endif
448     ewcorr_S0   = gmx_simd_mul_r(gmx_simd_pmecorrF_r(brsq_S0), beta_S);
449     ewcorr_S2   = gmx_simd_mul_r(gmx_simd_pmecorrF_r(brsq_S2), beta_S);
450     frcoul_S0   = gmx_simd_mul_r(qq_S0, gmx_simd_fmadd_r(ewcorr_S0, brsq_S0, rinv_ex_S0));
451     frcoul_S2   = gmx_simd_mul_r(qq_S2, gmx_simd_fmadd_r(ewcorr_S2, brsq_S2, rinv_ex_S2));
452
453 #ifdef CALC_ENERGIES
454     vc_sub_S0   = gmx_simd_mul_r(gmx_simd_pmecorrV_r(brsq_S0), beta_S);
455     vc_sub_S2   = gmx_simd_mul_r(gmx_simd_pmecorrV_r(brsq_S2), beta_S);
456 #endif
457
458 #endif /* CALC_COUL_EWALD */
459
460 #ifdef CALC_COUL_TAB
461     /* Electrostatic interactions */
462     r_S0        = gmx_simd_mul_r(rsq_S0, rinv_S0);
463     r_S2        = gmx_simd_mul_r(rsq_S2, rinv_S2);
464     /* Convert r to scaled table units */
465     rs_S0       = gmx_simd_mul_r(r_S0, invtsp_S);
466     rs_S2       = gmx_simd_mul_r(r_S2, invtsp_S);
467     /* Truncate scaled r to an int */
468     ti_S0       = gmx_simd_cvtt_r2i(rs_S0);
469     ti_S2       = gmx_simd_cvtt_r2i(rs_S2);
470 #ifdef GMX_SIMD_HAVE_FLOOR
471     rf_S0       = gmx_simd_floor_r(rs_S0);
472     rf_S2       = gmx_simd_floor_r(rs_S2);
473 #else
474     rf_S0       = gmx_simd_cvt_i2r(ti_S0);
475     rf_S2       = gmx_simd_cvt_i2r(ti_S2);
476 #endif
477     frac_S0     = gmx_simd_sub_r(rs_S0, rf_S0);
478     frac_S2     = gmx_simd_sub_r(rs_S2, rf_S2);
479
480     /* Load and interpolate table forces and possibly energies.
481      * Force and energy can be combined in one table, stride 4: FDV0
482      * or in two separate tables with stride 1: F and V
483      * Currently single precision uses FDV0, double F and V.
484      */
485 #ifndef CALC_ENERGIES
486     load_table_f(tab_coul_F, ti_S0, ti0, &ctab0_S0, &ctab1_S0);
487     load_table_f(tab_coul_F, ti_S2, ti2, &ctab0_S2, &ctab1_S2);
488 #else
489 #ifdef TAB_FDV0
490     load_table_f_v(tab_coul_F, ti_S0, ti0, &ctab0_S0, &ctab1_S0, &ctabv_S0);
491     load_table_f_v(tab_coul_F, ti_S2, ti2, &ctab0_S2, &ctab1_S2, &ctabv_S2);
492 #else
493     load_table_f_v(tab_coul_F, tab_coul_V, ti_S0, ti0, &ctab0_S0, &ctab1_S0, &ctabv_S0);
494     load_table_f_v(tab_coul_F, tab_coul_V, ti_S2, ti2, &ctab0_S2, &ctab1_S2, &ctabv_S2);
495 #endif
496 #endif
497     fsub_S0     = gmx_simd_add_r(ctab0_S0, gmx_simd_mul_r(frac_S0, ctab1_S0));
498     fsub_S2     = gmx_simd_add_r(ctab0_S2, gmx_simd_mul_r(frac_S2, ctab1_S2));
499     frcoul_S0   = gmx_simd_mul_r(qq_S0, gmx_simd_sub_r(rinv_ex_S0, gmx_simd_mul_r(fsub_S0, r_S0)));
500     frcoul_S2   = gmx_simd_mul_r(qq_S2, gmx_simd_sub_r(rinv_ex_S2, gmx_simd_mul_r(fsub_S2, r_S2)));
501
502 #ifdef CALC_ENERGIES
503     vc_sub_S0   = gmx_simd_add_r(ctabv_S0, gmx_simd_mul_r(gmx_simd_mul_r(mhalfsp_S, frac_S0), gmx_simd_add_r(ctab0_S0, fsub_S0)));
504     vc_sub_S2   = gmx_simd_add_r(ctabv_S2, gmx_simd_mul_r(gmx_simd_mul_r(mhalfsp_S, frac_S2), gmx_simd_add_r(ctab0_S2, fsub_S2)));
505 #endif
506 #endif /* CALC_COUL_TAB */
507
508 #if defined CALC_ENERGIES && (defined CALC_COUL_EWALD || defined CALC_COUL_TAB)
509 #ifndef NO_SHIFT_EWALD
510     /* Add Ewald potential shift to vc_sub for convenience */
511 #ifdef CHECK_EXCLS
512     vc_sub_S0   = gmx_simd_add_r(vc_sub_S0, gmx_simd_blendzero_r(sh_ewald_S, interact_S0));
513     vc_sub_S2   = gmx_simd_add_r(vc_sub_S2, gmx_simd_blendzero_r(sh_ewald_S, interact_S2));
514 #else
515     vc_sub_S0   = gmx_simd_add_r(vc_sub_S0, sh_ewald_S);
516     vc_sub_S2   = gmx_simd_add_r(vc_sub_S2, sh_ewald_S);
517 #endif
518 #endif
519
520     vcoul_S0    = gmx_simd_mul_r(qq_S0, gmx_simd_sub_r(rinv_ex_S0, vc_sub_S0));
521     vcoul_S2    = gmx_simd_mul_r(qq_S2, gmx_simd_sub_r(rinv_ex_S2, vc_sub_S2));
522 #endif
523
524 #ifdef CALC_ENERGIES
525     /* Mask energy for cut-off and diagonal */
526     vcoul_S0    = gmx_simd_blendzero_r(vcoul_S0, wco_S0);
527     vcoul_S2    = gmx_simd_blendzero_r(vcoul_S2, wco_S2);
528 #endif
529
530 #endif /* CALC_COULOMB */
531
532 #ifdef CALC_LJ
533     /* Lennard-Jones interaction */
534
535 #ifdef VDW_CUTOFF_CHECK
536     wco_vdw_S0  = gmx_simd_cmplt_r(rsq_S0, rcvdw2_S);
537 #ifndef HALF_LJ
538     wco_vdw_S2  = gmx_simd_cmplt_r(rsq_S2, rcvdw2_S);
539 #endif
540 #else
541     /* Same cut-off for Coulomb and VdW, reuse the registers */
542 #define     wco_vdw_S0    wco_S0
543 #define     wco_vdw_S2    wco_S2
544 #endif
545
546 #ifndef LJ_COMB_LB
547     rinvsix_S0  = gmx_simd_mul_r(rinvsq_S0, gmx_simd_mul_r(rinvsq_S0, rinvsq_S0));
548 #ifdef EXCL_FORCES
549     rinvsix_S0  = gmx_simd_blendzero_r(rinvsix_S0, interact_S0);
550 #endif
551 #ifndef HALF_LJ
552     rinvsix_S2  = gmx_simd_mul_r(rinvsq_S2, gmx_simd_mul_r(rinvsq_S2, rinvsq_S2));
553 #ifdef EXCL_FORCES
554     rinvsix_S2  = gmx_simd_blendzero_r(rinvsix_S2, interact_S2);
555 #endif
556 #endif
557
558 #if defined LJ_CUT || defined LJ_POT_SWITCH
559     /* We have plain LJ or LJ-PME with simple C6/6 C12/12 coefficients */
560     FrLJ6_S0    = gmx_simd_mul_r(c6_S0, rinvsix_S0);
561 #ifndef HALF_LJ
562     FrLJ6_S2    = gmx_simd_mul_r(c6_S2, rinvsix_S2);
563 #endif
564     FrLJ12_S0   = gmx_simd_mul_r(c12_S0, gmx_simd_mul_r(rinvsix_S0, rinvsix_S0));
565 #ifndef HALF_LJ
566     FrLJ12_S2   = gmx_simd_mul_r(c12_S2, gmx_simd_mul_r(rinvsix_S2, rinvsix_S2));
567 #endif
568 #endif
569
570 #if defined LJ_FORCE_SWITCH || defined LJ_POT_SWITCH
571     /* We switch the LJ force */
572     r_S0        = gmx_simd_mul_r(rsq_S0, rinv_S0);
573     rsw_S0      = gmx_simd_max_r(gmx_simd_sub_r(r_S0, rswitch_S), zero_S);
574     rsw2_S0     = gmx_simd_mul_r(rsw_S0, rsw_S0);
575     rsw2_r_S0   = gmx_simd_mul_r(rsw2_S0, r_S0);
576 #ifndef HALF_LJ
577     r_S2        = gmx_simd_mul_r(rsq_S2, rinv_S2);
578     rsw_S2      = gmx_simd_max_r(gmx_simd_sub_r(r_S2, rswitch_S), zero_S);
579     rsw2_S2     = gmx_simd_mul_r(rsw_S2, rsw_S2);
580     rsw2_r_S2   = gmx_simd_mul_r(rsw2_S2, r_S2);
581 #endif
582 #endif
583
584 #ifdef LJ_FORCE_SWITCH
585
586 #define add_fr_switch(fr, rsw, rsw2_r, c2, c3) gmx_simd_fmadd_r(gmx_simd_fmadd_r(c3, rsw, c2), rsw2_r, fr)
587
588     FrLJ6_S0    = gmx_simd_mul_r(c6_S0, add_fr_switch(rinvsix_S0, rsw_S0, rsw2_r_S0, p6_fc2_S, p6_fc3_S));
589 #ifndef HALF_LJ
590     FrLJ6_S2    = gmx_simd_mul_r(c6_S2, add_fr_switch(rinvsix_S2, rsw_S2, rsw2_r_S2, p6_fc2_S, p6_fc3_S));
591 #endif
592     FrLJ12_S0   = gmx_simd_mul_r(c12_S0, add_fr_switch(gmx_simd_mul_r(rinvsix_S0, rinvsix_S0), rsw_S0, rsw2_r_S0, p12_fc2_S, p12_fc3_S));
593 #ifndef HALF_LJ
594     FrLJ12_S2   = gmx_simd_mul_r(c12_S2, add_fr_switch(gmx_simd_mul_r(rinvsix_S2, rinvsix_S2), rsw_S2, rsw2_r_S2, p12_fc2_S, p12_fc3_S));
595 #endif
596 #undef add_fr_switch
597 #endif /* LJ_FORCE_SWITCH */
598
599 #endif /* not LJ_COMB_LB */
600
601 #ifdef LJ_COMB_LB
602     sir_S0      = gmx_simd_mul_r(sig_S0, rinv_S0);
603 #ifndef HALF_LJ
604     sir_S2      = gmx_simd_mul_r(sig_S2, rinv_S2);
605 #endif
606     sir2_S0     = gmx_simd_mul_r(sir_S0, sir_S0);
607 #ifndef HALF_LJ
608     sir2_S2     = gmx_simd_mul_r(sir_S2, sir_S2);
609 #endif
610     sir6_S0     = gmx_simd_mul_r(sir2_S0, gmx_simd_mul_r(sir2_S0, sir2_S0));
611 #ifdef EXCL_FORCES
612     sir6_S0     = gmx_simd_blendzero_r(sir6_S0, interact_S0);
613 #endif
614 #ifndef HALF_LJ
615     sir6_S2     = gmx_simd_mul_r(sir2_S2, gmx_simd_mul_r(sir2_S2, sir2_S2));
616 #ifdef EXCL_FORCES
617     sir6_S2     = gmx_simd_blendzero_r(sir6_S2, interact_S2);
618 #endif
619 #endif
620 #ifdef VDW_CUTOFF_CHECK
621     sir6_S0     = gmx_simd_blendzero_r(sir6_S0, wco_vdw_S0);
622 #ifndef HALF_LJ
623     sir6_S2     = gmx_simd_blendzero_r(sir6_S2, wco_vdw_S2);
624 #endif
625 #endif
626     FrLJ6_S0    = gmx_simd_mul_r(eps_S0, sir6_S0);
627 #ifndef HALF_LJ
628     FrLJ6_S2    = gmx_simd_mul_r(eps_S2, sir6_S2);
629 #endif
630     FrLJ12_S0   = gmx_simd_mul_r(FrLJ6_S0, sir6_S0);
631 #ifndef HALF_LJ
632     FrLJ12_S2   = gmx_simd_mul_r(FrLJ6_S2, sir6_S2);
633 #endif
634 #if defined CALC_ENERGIES
635     /* We need C6 and C12 to calculate the LJ potential shift */
636     sig2_S0     = gmx_simd_mul_r(sig_S0, sig_S0);
637 #ifndef HALF_LJ
638     sig2_S2     = gmx_simd_mul_r(sig_S2, sig_S2);
639 #endif
640     sig6_S0     = gmx_simd_mul_r(sig2_S0, gmx_simd_mul_r(sig2_S0, sig2_S0));
641 #ifndef HALF_LJ
642     sig6_S2     = gmx_simd_mul_r(sig2_S2, gmx_simd_mul_r(sig2_S2, sig2_S2));
643 #endif
644     c6_S0       = gmx_simd_mul_r(eps_S0, sig6_S0);
645 #ifndef HALF_LJ
646     c6_S2       = gmx_simd_mul_r(eps_S2, sig6_S2);
647 #endif
648     c12_S0      = gmx_simd_mul_r(c6_S0, sig6_S0);
649 #ifndef HALF_LJ
650     c12_S2      = gmx_simd_mul_r(c6_S2, sig6_S2);
651 #endif
652 #endif
653 #endif /* LJ_COMB_LB */
654
655     /* Determine the total scalar LJ force*r */
656     frLJ_S0     = gmx_simd_sub_r(FrLJ12_S0, FrLJ6_S0);
657 #ifndef HALF_LJ
658     frLJ_S2     = gmx_simd_sub_r(FrLJ12_S2, FrLJ6_S2);
659 #endif
660
661 #if (defined LJ_CUT || defined LJ_FORCE_SWITCH) && defined CALC_ENERGIES
662
663 #ifdef LJ_CUT
664     /* Calculate the LJ energies, with constant potential shift */
665     VLJ6_S0     = gmx_simd_mul_r(sixth_S, gmx_simd_fmadd_r(c6_S0, p6_cpot_S, FrLJ6_S0));
666 #ifndef HALF_LJ
667     VLJ6_S2     = gmx_simd_mul_r(sixth_S, gmx_simd_fmadd_r(c6_S2, p6_cpot_S, FrLJ6_S2));
668 #endif
669     VLJ12_S0    = gmx_simd_mul_r(twelveth_S, gmx_simd_fmadd_r(c12_S0, p12_cpot_S, FrLJ12_S0));
670 #ifndef HALF_LJ
671     VLJ12_S2    = gmx_simd_mul_r(twelveth_S, gmx_simd_fmadd_r(c12_S2, p12_cpot_S, FrLJ12_S2));
672 #endif
673 #endif /* LJ_CUT */
674
675 #ifdef LJ_FORCE_SWITCH
676 #define v_fswitch_pr(rsw, rsw2, c0, c3, c4) gmx_simd_fmadd_r(gmx_simd_fmadd_r(c4, rsw, c3), gmx_simd_mul_r(rsw2, rsw), c0)
677
678     VLJ6_S0     = gmx_simd_mul_r(c6_S0, gmx_simd_fmadd_r(sixth_S, rinvsix_S0, v_fswitch_pr(rsw_S0, rsw2_S0, p6_6cpot_S, p6_vc3_S, p6_vc4_S)));
679 #ifndef HALF_LJ
680     VLJ6_S2     = gmx_simd_mul_r(c6_S2, gmx_simd_fmadd_r(sixth_S, rinvsix_S2, v_fswitch_pr(rsw_S2, rsw2_S2, p6_6cpot_S, p6_vc3_S, p6_vc4_S)));
681 #endif
682     VLJ12_S0    = gmx_simd_mul_r(c12_S0, gmx_simd_fmadd_r(twelveth_S, gmx_simd_mul_r(rinvsix_S0, rinvsix_S0), v_fswitch_pr(rsw_S0, rsw2_S0, p12_12cpot_S, p12_vc3_S, p12_vc4_S)));
683 #ifndef HALF_LJ
684     VLJ12_S2    = gmx_simd_mul_r(c12_S2, gmx_simd_fmadd_r(twelveth_S, gmx_simd_mul_r(rinvsix_S2, rinvsix_S2), v_fswitch_pr(rsw_S2, rsw2_S2, p12_12cpot_S, p12_vc3_S, p12_vc4_S)));
685 #endif
686 #undef v_fswitch_pr
687 #endif /* LJ_FORCE_SWITCH */
688
689     /* Add up the repulsion and dispersion */
690     VLJ_S0      = gmx_simd_sub_r(VLJ12_S0, VLJ6_S0);
691 #ifndef HALF_LJ
692     VLJ_S2      = gmx_simd_sub_r(VLJ12_S2, VLJ6_S2);
693 #endif
694
695 #endif /* (LJ_CUT || LJ_FORCE_SWITCH) && CALC_ENERGIES */
696
697 #ifdef LJ_POT_SWITCH
698     /* We always need the potential, since it is needed for the force */
699     VLJ_S0 = gmx_simd_fnmadd_r(sixth_S, FrLJ6_S0, gmx_simd_mul_r(twelveth_S, FrLJ12_S0));
700 #ifndef HALF_LJ
701     VLJ_S2 = gmx_simd_fnmadd_r(sixth_S, FrLJ6_S2, gmx_simd_mul_r(twelveth_S, FrLJ12_S2));
702 #endif
703
704     {
705         gmx_simd_real_t sw_S0, dsw_S0;
706 #ifndef HALF_LJ
707         gmx_simd_real_t sw_S2, dsw_S2;
708 #endif
709
710 #define switch_pr(rsw, rsw2, c3, c4, c5) gmx_simd_fmadd_r(gmx_simd_fmadd_r(gmx_simd_fmadd_r(c5, rsw, c4), rsw, c3), gmx_simd_mul_r(rsw2, rsw), one_S)
711 #define dswitch_pr(rsw, rsw2, c2, c3, c4) gmx_simd_mul_r(gmx_simd_fmadd_r(gmx_simd_fmadd_r(c4, rsw, c3), rsw, c2), rsw2)
712
713         sw_S0  = switch_pr(rsw_S0, rsw2_S0, swV3_S, swV4_S, swV5_S);
714         dsw_S0 = dswitch_pr(rsw_S0, rsw2_S0, swF2_S, swF3_S, swF4_S);
715 #ifndef HALF_LJ
716         sw_S2  = switch_pr(rsw_S2, rsw2_S2, swV3_S, swV4_S, swV5_S);
717         dsw_S2 = dswitch_pr(rsw_S2, rsw2_S2, swF2_S, swF3_S, swF4_S);
718 #endif
719         frLJ_S0 = gmx_simd_fnmadd_r(gmx_simd_mul_r(dsw_S0, VLJ_S0), r_S0, gmx_simd_mul_r(sw_S0, frLJ_S0));
720 #ifndef HALF_LJ
721         frLJ_S2 = gmx_simd_fnmadd_r(gmx_simd_mul_r(dsw_S2, VLJ_S2), r_S2, gmx_simd_mul_r(sw_S2, frLJ_S2));
722 #endif
723 #ifdef CALC_ENERGIES
724         VLJ_S0  = gmx_simd_mul_r(sw_S0, VLJ_S0);
725 #ifndef HALF_LJ
726         VLJ_S2  = gmx_simd_mul_r(sw_S2, VLJ_S2);
727 #endif
728 #endif
729
730 #undef switch_pr
731 #undef dswitch_pr
732     }
733 #endif /* LJ_POT_SWITCH */
734
735 #if defined CALC_ENERGIES && defined CHECK_EXCLS
736     /* The potential shift should be removed for excluded pairs */
737     VLJ_S0      = gmx_simd_blendzero_r(VLJ_S0, interact_S0);
738 #ifndef HALF_LJ
739     VLJ_S2      = gmx_simd_blendzero_r(VLJ_S2, interact_S2);
740 #endif
741 #endif
742
743 #ifdef LJ_EWALD_GEOM
744     {
745         gmx_simd_real_t c6s_j_S;
746         gmx_simd_real_t c6grid_S0, rinvsix_nm_S0, cr2_S0, expmcr2_S0, poly_S0;
747 #ifndef HALF_LJ
748         gmx_simd_real_t c6grid_S2, rinvsix_nm_S2, cr2_S2, expmcr2_S2, poly_S2;
749 #endif
750 #ifdef CALC_ENERGIES
751         gmx_simd_real_t sh_mask_S0;
752 #ifndef HALF_LJ
753         gmx_simd_real_t sh_mask_S2;
754 #endif
755 #endif
756
757         /* Determine C6 for the grid using the geometric combination rule */
758         gmx_loaddh_pr(&c6s_j_S,  ljc+aj2+0);
759         c6grid_S0       = gmx_simd_mul_r(c6s_S0, c6s_j_S);
760 #ifndef HALF_LJ
761         c6grid_S2       = gmx_simd_mul_r(c6s_S2, c6s_j_S);
762 #endif
763
764 #ifdef CHECK_EXCLS
765         /* Recalculate rinvsix without exclusion mask (compiler might optimize) */
766         rinvsix_nm_S0 = gmx_simd_mul_r(rinvsq_S0, gmx_simd_mul_r(rinvsq_S0, rinvsq_S0));
767 #ifndef HALF_LJ
768         rinvsix_nm_S2 = gmx_simd_mul_r(rinvsq_S2, gmx_simd_mul_r(rinvsq_S2, rinvsq_S2));
769 #endif
770 #else
771         /* We didn't use a mask, so we can copy */
772         rinvsix_nm_S0 = rinvsix_S0;
773 #ifndef HALF_LJ
774         rinvsix_nm_S2 = rinvsix_S2;
775 #endif
776 #endif
777
778         cr2_S0        = gmx_simd_mul_r(lje_c2_S, rsq_S0);
779 #ifndef HALF_LJ
780         cr2_S2        = gmx_simd_mul_r(lje_c2_S, rsq_S2);
781 #endif
782         expmcr2_S0    = gmx_simd_exp_r(gmx_simd_mul_r(mone_S, cr2_S0));
783 #ifndef HALF_LJ
784         expmcr2_S2    = gmx_simd_exp_r(gmx_simd_mul_r(mone_S, cr2_S2));
785 #endif
786
787         /* 1 + cr2 + 1/2*cr2^2 */
788         poly_S0       = gmx_simd_fmadd_r(gmx_simd_fmadd_r(half_S, cr2_S0, one_S), cr2_S0, one_S);
789 #ifndef HALF_LJ
790         poly_S2       = gmx_simd_fmadd_r(gmx_simd_fmadd_r(half_S, cr2_S2, one_S), cr2_S2, one_S);
791 #endif
792
793         /* We calculate LJ F*r = (6*C6)*(r^-6 - F_mesh/6), we use:
794          * r^-6*cexp*(1 + cr2 + cr2^2/2 + cr2^3/6) = cexp*(r^-6*poly + c^6/6)
795          */
796         frLJ_S0       = gmx_simd_fmadd_r(c6grid_S0, gmx_simd_fnmadd_r(expmcr2_S0, gmx_simd_fmadd_r(rinvsix_nm_S0, poly_S0, lje_c6_6_S), rinvsix_nm_S0), frLJ_S0);
797 #ifndef HALF_LJ
798         frLJ_S2       = gmx_simd_fmadd_r(c6grid_S2, gmx_simd_fnmadd_r(expmcr2_S2, gmx_simd_fmadd_r(rinvsix_nm_S2, poly_S2, lje_c6_6_S), rinvsix_nm_S2), frLJ_S2);
799 #endif
800
801 #ifdef CALC_ENERGIES
802 #ifdef CHECK_EXCLS
803         sh_mask_S0    = gmx_simd_blendzero_r(lje_vc_S, interact_S0);
804 #ifndef HALF_LJ
805         sh_mask_S2    = gmx_simd_blendzero_r(lje_vc_S, interact_S2);
806 #endif
807 #else
808         sh_mask_S0    = lje_vc_S;
809 #ifndef HALF_LJ
810         sh_mask_S2    = lje_vc_S;
811 #endif
812 #endif
813
814         VLJ_S0        = gmx_simd_fmadd_r(gmx_simd_mul_r(sixth_S, c6grid_S0), gmx_simd_fmadd_r(rinvsix_nm_S0, gmx_simd_fnmadd_r(expmcr2_S0, poly_S0, one_S), sh_mask_S0), VLJ_S0);
815 #ifndef HALF_LJ
816         VLJ_S2        = gmx_simd_fmadd_r(gmx_simd_mul_r(sixth_S, c6grid_S2), gmx_simd_fmadd_r(rinvsix_nm_S2, gmx_simd_fnmadd_r(expmcr2_S2, poly_S2, one_S), sh_mask_S2), VLJ_S2);
817 #endif
818 #endif /* CALC_ENERGIES */
819     }
820 #endif /* LJ_EWALD_GEOM */
821
822 #if defined VDW_CUTOFF_CHECK
823     /* frLJ is multiplied later by rinvsq, which is masked for the Coulomb
824      * cut-off, but if the VdW cut-off is shorter, we need to mask with that.
825      */
826     frLJ_S0     = gmx_simd_blendzero_r(frLJ_S0, wco_vdw_S0);
827 #ifndef HALF_LJ
828     frLJ_S2     = gmx_simd_blendzero_r(frLJ_S2, wco_vdw_S2);
829 #endif
830 #endif
831
832 #ifdef CALC_ENERGIES
833     /* The potential shift should be removed for pairs beyond cut-off */
834     VLJ_S0      = gmx_simd_blendzero_r(VLJ_S0, wco_vdw_S0);
835 #ifndef HALF_LJ
836     VLJ_S2      = gmx_simd_blendzero_r(VLJ_S2, wco_vdw_S2);
837 #endif
838 #endif
839
840 #endif /* CALC_LJ */
841
842 #ifdef CALC_ENERGIES
843 #ifdef ENERGY_GROUPS
844     /* Extract the group pair index per j pair.
845      * Energy groups are stored per i-cluster, so things get
846      * complicated when the i- and j-cluster size don't match.
847      */
848     {
849         int egps_j;
850 #if UNROLLJ == 2
851         egps_j    = nbat->energrp[cj>>1];
852         egp_jj[0] = ((egps_j >> ((cj & 1)*egps_jshift)) & egps_jmask)*egps_jstride;
853 #else
854         /* We assume UNROLLI <= UNROLLJ */
855         int jdi;
856         for (jdi = 0; jdi < UNROLLJ/UNROLLI; jdi++)
857         {
858             int jj;
859             egps_j = nbat->energrp[cj*(UNROLLJ/UNROLLI)+jdi];
860             for (jj = 0; jj < (UNROLLI/2); jj++)
861             {
862                 egp_jj[jdi*(UNROLLI/2)+jj] = ((egps_j >> (jj*egps_jshift)) & egps_jmask)*egps_jstride;
863             }
864         }
865 #endif
866     }
867 #endif
868
869 #ifdef CALC_COULOMB
870 #ifndef ENERGY_GROUPS
871     vctot_S      = gmx_simd_add_r(vctot_S, gmx_simd_add_r(vcoul_S0, vcoul_S2));
872 #else
873     add_ener_grp_halves(vcoul_S0, vctp[0], vctp[1], egp_jj);
874     add_ener_grp_halves(vcoul_S2, vctp[2], vctp[3], egp_jj);
875 #endif
876 #endif
877
878 #ifdef CALC_LJ
879 #ifndef ENERGY_GROUPS
880     Vvdwtot_S    = gmx_simd_add_r(Vvdwtot_S,
881 #ifndef HALF_LJ
882                                   gmx_simd_add_r(VLJ_S0, VLJ_S2)
883 #else
884                                   VLJ_S0
885 #endif
886                                   );
887 #else
888     add_ener_grp_halves(VLJ_S0, vvdwtp[0], vvdwtp[1], egp_jj);
889 #ifndef HALF_LJ
890     add_ener_grp_halves(VLJ_S2, vvdwtp[2], vvdwtp[3], egp_jj);
891 #endif
892 #endif
893 #endif /* CALC_LJ */
894 #endif /* CALC_ENERGIES */
895
896 #ifdef CALC_LJ
897 #ifdef CALC_COULOMB
898     fscal_S0    = gmx_simd_mul_r(rinvsq_S0, gmx_simd_add_r(frcoul_S0, frLJ_S0));
899 #else
900     fscal_S0    = gmx_simd_mul_r(rinvsq_S0, frLJ_S0);
901 #endif
902 #else
903     fscal_S0    = gmx_simd_mul_r(rinvsq_S0, frcoul_S0);
904 #endif /* CALC_LJ */
905 #if defined CALC_LJ && !defined HALF_LJ
906 #ifdef CALC_COULOMB
907     fscal_S2    = gmx_simd_mul_r(rinvsq_S2, gmx_simd_add_r(frcoul_S2, frLJ_S2));
908 #else
909     fscal_S2    = gmx_simd_mul_r(rinvsq_S2, frLJ_S2);
910 #endif
911 #else
912     /* Atom 2 and 3 don't have LJ, so only add Coulomb forces */
913     fscal_S2    = gmx_simd_mul_r(rinvsq_S2, frcoul_S2);
914 #endif
915
916     /* Calculate temporary vectorial force */
917     tx_S0       = gmx_simd_mul_r(fscal_S0, dx_S0);
918     tx_S2       = gmx_simd_mul_r(fscal_S2, dx_S2);
919     ty_S0       = gmx_simd_mul_r(fscal_S0, dy_S0);
920     ty_S2       = gmx_simd_mul_r(fscal_S2, dy_S2);
921     tz_S0       = gmx_simd_mul_r(fscal_S0, dz_S0);
922     tz_S2       = gmx_simd_mul_r(fscal_S2, dz_S2);
923
924     /* Increment i atom force */
925     fix_S0      = gmx_simd_add_r(fix_S0, tx_S0);
926     fix_S2      = gmx_simd_add_r(fix_S2, tx_S2);
927     fiy_S0      = gmx_simd_add_r(fiy_S0, ty_S0);
928     fiy_S2      = gmx_simd_add_r(fiy_S2, ty_S2);
929     fiz_S0      = gmx_simd_add_r(fiz_S0, tz_S0);
930     fiz_S2      = gmx_simd_add_r(fiz_S2, tz_S2);
931
932     /* Decrement j atom force */
933     gmx_load_hpr(&fjx_S, f+ajx);
934     gmx_load_hpr(&fjy_S, f+ajy);
935     gmx_load_hpr(&fjz_S, f+ajz);
936     gmx_store_hpr(f+ajx, gmx_sub_hpr(fjx_S, gmx_sum4_hpr(tx_S0, tx_S2)));
937     gmx_store_hpr(f+ajy, gmx_sub_hpr(fjy_S, gmx_sum4_hpr(ty_S0, ty_S2)));
938     gmx_store_hpr(f+ajz, gmx_sub_hpr(fjz_S, gmx_sum4_hpr(tz_S0, tz_S2)));
939 }
940
941 #undef  rinv_ex_S0
942 #undef  rinv_ex_S2
943
944 #undef  wco_vdw_S0
945 #undef  wco_vdw_S2
946
947 #undef  NBNXN_CUTOFF_USE_BLENDV
948
949 #undef  EXCL_FORCES