2f45b53bafe0acd349b949ad2bf60e603a0e1220
[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 LJ_EWALD_GEOM) && defined GMX_SIMD_HAVE_BLENDV
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     aj2           = aj*2;
251 #endif
252     ajx           = aj*DIM;
253     ajy           = ajx + STRIDE;
254     ajz           = ajy + STRIDE;
255
256 #ifdef CHECK_EXCLS
257     gmx_load_simd_2xnn_interactions(l_cj[cjind].excl,
258                                     filter_S0, filter_S2,
259                                     &interact_S0, &interact_S2);
260 #endif /* CHECK_EXCLS */
261
262     /* load j atom coordinates */
263     gmx_loaddh_pr(&jx_S, x+ajx);
264     gmx_loaddh_pr(&jy_S, x+ajy);
265     gmx_loaddh_pr(&jz_S, x+ajz);
266
267     /* Calculate distance */
268     dx_S0       = gmx_simd_sub_r(ix_S0, jx_S);
269     dy_S0       = gmx_simd_sub_r(iy_S0, jy_S);
270     dz_S0       = gmx_simd_sub_r(iz_S0, jz_S);
271     dx_S2       = gmx_simd_sub_r(ix_S2, jx_S);
272     dy_S2       = gmx_simd_sub_r(iy_S2, jy_S);
273     dz_S2       = gmx_simd_sub_r(iz_S2, jz_S);
274
275     /* rsq = dx*dx+dy*dy+dz*dz */
276     rsq_S0      = gmx_simd_calc_rsq_r(dx_S0, dy_S0, dz_S0);
277     rsq_S2      = gmx_simd_calc_rsq_r(dx_S2, dy_S2, dz_S2);
278
279 #ifndef NBNXN_CUTOFF_USE_BLENDV
280     wco_S0      = gmx_simd_cmplt_r(rsq_S0, rc2_S);
281     wco_S2      = gmx_simd_cmplt_r(rsq_S2, rc2_S);
282 #endif
283
284 #ifdef CHECK_EXCLS
285 #ifdef EXCL_FORCES
286     /* Only remove the (sub-)diagonal to avoid double counting */
287 #if UNROLLJ == UNROLLI
288     if (cj == ci_sh)
289     {
290         wco_S0  = gmx_simd_and_b(wco_S0, diagonal_mask_S0);
291         wco_S2  = gmx_simd_and_b(wco_S2, diagonal_mask_S2);
292     }
293 #else
294 #if UNROLLJ == 2*UNROLLI
295     if (cj*2 == ci_sh)
296     {
297         wco_S0  = gmx_simd_and_b(wco_S0, diagonal_mask0_S0);
298         wco_S2  = gmx_simd_and_b(wco_S2, diagonal_mask0_S2);
299     }
300     else if (cj*2 + 1 == ci_sh)
301     {
302         wco_S0  = gmx_simd_and_b(wco_S0, diagonal_mask1_S0);
303         wco_S2  = gmx_simd_and_b(wco_S2, diagonal_mask1_S2);
304     }
305 #else
306 #error "only UNROLLJ == UNROLLI*(1 or 2) currently supported in 2xnn kernels"
307 #endif
308 #endif
309 #else /* EXCL_FORCES */
310       /* No exclusion forces: remove all excluded atom pairs from the list */
311     wco_S0      = gmx_simd_and_b(wco_S0, interact_S0);
312     wco_S2      = gmx_simd_and_b(wco_S2, interact_S2);
313 #endif
314 #endif
315
316 #ifdef COUNT_PAIRS
317     {
318         int  i, j;
319         real tmpa[2*GMX_SIMD_REAL_WIDTH], *tmp;
320         tmp = gmx_simd_align_r(tmpa);
321         for (i = 0; i < UNROLLI; i += 2)
322         {
323             gmx_simd_store_r(tmp, gmx_simd_sub_r(rc2_S, i == 0 ? rsq_S0 : rsq_S2));
324             for (j = 0; j < 2*UNROLLJ; j++)
325             {
326                 if (tmp[j] >= 0)
327                 {
328                     npair++;
329                 }
330             }
331         }
332     }
333 #endif
334
335 #ifdef CHECK_EXCLS
336     /* For excluded pairs add a small number to avoid r^-6 = NaN */
337     rsq_S0      = gmx_simd_add_r(rsq_S0, gmx_simd_blendnotzero_r(avoid_sing_S, interact_S0));
338     rsq_S2      = gmx_simd_add_r(rsq_S2, gmx_simd_blendnotzero_r(avoid_sing_S, interact_S2));
339 #endif
340
341     /* Calculate 1/r */
342     rinv_S0     = gmx_simd_invsqrt_r(rsq_S0);
343     rinv_S2     = gmx_simd_invsqrt_r(rsq_S2);
344
345 #ifdef CALC_COULOMB
346     /* Load parameters for j atom */
347     gmx_loaddh_pr(&jq_S, q+aj);
348     qq_S0       = gmx_simd_mul_r(iq_S0, jq_S);
349     qq_S2       = gmx_simd_mul_r(iq_S2, jq_S);
350 #endif
351
352 #ifdef CALC_LJ
353
354 #if !defined LJ_COMB_GEOM && !defined LJ_COMB_LB && !defined FIX_LJ_C
355     load_lj_pair_params2(nbfp0, nbfp1, type, aj, &c6_S0, &c12_S0);
356 #ifndef HALF_LJ
357     load_lj_pair_params2(nbfp2, nbfp3, type, aj, &c6_S2, &c12_S2);
358 #endif
359 #endif /* not defined any LJ rule */
360
361 #ifdef LJ_COMB_GEOM
362     gmx_loaddh_pr(&c6s_j_S,  ljc+aj2+0);
363     gmx_loaddh_pr(&c12s_j_S, ljc+aj2+STRIDE);
364     c6_S0       = gmx_simd_mul_r(c6s_S0, c6s_j_S );
365 #ifndef HALF_LJ
366     c6_S2       = gmx_simd_mul_r(c6s_S2, c6s_j_S );
367 #endif
368     c12_S0      = gmx_simd_mul_r(c12s_S0, c12s_j_S);
369 #ifndef HALF_LJ
370     c12_S2      = gmx_simd_mul_r(c12s_S2, c12s_j_S);
371 #endif
372 #endif /* LJ_COMB_GEOM */
373
374 #ifdef LJ_COMB_LB
375     gmx_loaddh_pr(&hsig_j_S, ljc+aj2+0);
376     gmx_loaddh_pr(&seps_j_S, ljc+aj2+STRIDE);
377
378     sig_S0      = gmx_simd_add_r(hsig_i_S0, hsig_j_S);
379     eps_S0      = gmx_simd_mul_r(seps_i_S0, seps_j_S);
380 #ifndef HALF_LJ
381     sig_S2      = gmx_simd_add_r(hsig_i_S2, hsig_j_S);
382     eps_S2      = gmx_simd_mul_r(seps_i_S2, seps_j_S);
383 #endif
384 #endif /* LJ_COMB_LB */
385
386 #endif /* CALC_LJ */
387
388 #ifndef NBNXN_CUTOFF_USE_BLENDV
389     rinv_S0     = gmx_simd_blendzero_r(rinv_S0, wco_S0);
390     rinv_S2     = gmx_simd_blendzero_r(rinv_S2, wco_S2);
391 #else
392     /* This needs to be modified: It makes assumptions about the internal storage
393      * of the SIMD representation, in particular that the blendv instruction always
394      * selects based on the sign bit. If the performance is really critical, it
395      * should be turned into a function that is platform-specific.
396      */
397     /* We only need to mask for the cut-off: blendv is faster */
398     rinv_S0     = gmx_simd_blendv_r(rinv_S0, zero_S, gmx_simd_sub_r(rc2_S, rsq_S0));
399     rinv_S2     = gmx_simd_blendv_r(rinv_S2, zero_S, gmx_simd_sub_r(rc2_S, rsq_S2));
400 #endif
401
402     rinvsq_S0   = gmx_simd_mul_r(rinv_S0, rinv_S0);
403     rinvsq_S2   = gmx_simd_mul_r(rinv_S2, rinv_S2);
404
405 #ifdef CALC_COULOMB
406     /* Note that here we calculate force*r, not the usual force/r.
407      * This allows avoiding masking the reaction-field contribution,
408      * as frcoul is later multiplied by rinvsq which has been
409      * masked with the cut-off check.
410      */
411
412 #ifdef EXCL_FORCES
413     /* Only add 1/r for non-excluded atom pairs */
414     rinv_ex_S0  = gmx_simd_blendzero_r(rinv_S0, interact_S0);
415     rinv_ex_S2  = gmx_simd_blendzero_r(rinv_S2, interact_S2);
416 #else
417     /* No exclusion forces, we always need 1/r */
418 #define     rinv_ex_S0    rinv_S0
419 #define     rinv_ex_S2    rinv_S2
420 #endif
421
422 #ifdef CALC_COUL_RF
423     /* Electrostatic interactions */
424     frcoul_S0   = gmx_simd_mul_r(qq_S0, gmx_simd_fmadd_r(rsq_S0, mrc_3_S, rinv_ex_S0));
425     frcoul_S2   = gmx_simd_mul_r(qq_S2, gmx_simd_fmadd_r(rsq_S2, mrc_3_S, rinv_ex_S2));
426
427 #ifdef CALC_ENERGIES
428     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)));
429     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)));
430 #endif
431 #endif
432
433 #ifdef CALC_COUL_EWALD
434     /* We need to mask (or limit) rsq for the cut-off,
435      * as large distances can cause an overflow in gmx_pmecorrF/V.
436      */
437 #ifndef NBNXN_CUTOFF_USE_BLENDV
438     brsq_S0     = gmx_simd_mul_r(beta2_S, gmx_simd_blendzero_r(rsq_S0, wco_S0));
439     brsq_S2     = gmx_simd_mul_r(beta2_S, gmx_simd_blendzero_r(rsq_S2, wco_S2));
440 #else
441     /* Strangely, putting mul on a separate line is slower (icc 13) */
442     brsq_S0     = gmx_simd_mul_r(beta2_S, gmx_simd_blendv_r(rsq_S0, zero_S, gmx_simd_sub_r(rc2_S, rsq_S0)));
443     brsq_S2     = gmx_simd_mul_r(beta2_S, gmx_simd_blendv_r(rsq_S2, zero_S, gmx_simd_sub_r(rc2_S, rsq_S2)));
444 #endif
445     ewcorr_S0   = gmx_simd_mul_r(gmx_simd_pmecorrF_r(brsq_S0), beta_S);
446     ewcorr_S2   = gmx_simd_mul_r(gmx_simd_pmecorrF_r(brsq_S2), beta_S);
447     frcoul_S0   = gmx_simd_mul_r(qq_S0, gmx_simd_fmadd_r(ewcorr_S0, brsq_S0, rinv_ex_S0));
448     frcoul_S2   = gmx_simd_mul_r(qq_S2, gmx_simd_fmadd_r(ewcorr_S2, brsq_S2, rinv_ex_S2));
449
450 #ifdef CALC_ENERGIES
451     vc_sub_S0   = gmx_simd_mul_r(gmx_simd_pmecorrV_r(brsq_S0), beta_S);
452     vc_sub_S2   = gmx_simd_mul_r(gmx_simd_pmecorrV_r(brsq_S2), beta_S);
453 #endif
454
455 #endif /* CALC_COUL_EWALD */
456
457 #ifdef CALC_COUL_TAB
458     /* Electrostatic interactions */
459     r_S0        = gmx_simd_mul_r(rsq_S0, rinv_S0);
460     r_S2        = gmx_simd_mul_r(rsq_S2, rinv_S2);
461     /* Convert r to scaled table units */
462     rs_S0       = gmx_simd_mul_r(r_S0, invtsp_S);
463     rs_S2       = gmx_simd_mul_r(r_S2, invtsp_S);
464     /* Truncate scaled r to an int */
465     ti_S0       = gmx_simd_cvtt_r2i(rs_S0);
466     ti_S2       = gmx_simd_cvtt_r2i(rs_S2);
467 #ifdef GMX_SIMD_HAVE_TRUNC
468     rf_S0       = gmx_simd_trunc_r(rs_S0);
469     rf_S2       = gmx_simd_trunc_r(rs_S2);
470 #else
471     rf_S0       = gmx_simd_cvt_i2r(ti_S0);
472     rf_S2       = gmx_simd_cvt_i2r(ti_S2);
473 #endif
474     frac_S0     = gmx_simd_sub_r(rs_S0, rf_S0);
475     frac_S2     = gmx_simd_sub_r(rs_S2, rf_S2);
476
477     /* Load and interpolate table forces and possibly energies.
478      * Force and energy can be combined in one table, stride 4: FDV0
479      * or in two separate tables with stride 1: F and V
480      * Currently single precision uses FDV0, double F and V.
481      */
482 #ifndef CALC_ENERGIES
483     load_table_f(tab_coul_F, ti_S0, ti0, &ctab0_S0, &ctab1_S0);
484     load_table_f(tab_coul_F, ti_S2, ti2, &ctab0_S2, &ctab1_S2);
485 #else
486 #ifdef TAB_FDV0
487     load_table_f_v(tab_coul_F, ti_S0, ti0, &ctab0_S0, &ctab1_S0, &ctabv_S0);
488     load_table_f_v(tab_coul_F, ti_S2, ti2, &ctab0_S2, &ctab1_S2, &ctabv_S2);
489 #else
490     load_table_f_v(tab_coul_F, tab_coul_V, ti_S0, ti0, &ctab0_S0, &ctab1_S0, &ctabv_S0);
491     load_table_f_v(tab_coul_F, tab_coul_V, ti_S2, ti2, &ctab0_S2, &ctab1_S2, &ctabv_S2);
492 #endif
493 #endif
494     fsub_S0     = gmx_simd_add_r(ctab0_S0, gmx_simd_mul_r(frac_S0, ctab1_S0));
495     fsub_S2     = gmx_simd_add_r(ctab0_S2, gmx_simd_mul_r(frac_S2, ctab1_S2));
496     frcoul_S0   = gmx_simd_mul_r(qq_S0, gmx_simd_sub_r(rinv_ex_S0, gmx_simd_mul_r(fsub_S0, r_S0)));
497     frcoul_S2   = gmx_simd_mul_r(qq_S2, gmx_simd_sub_r(rinv_ex_S2, gmx_simd_mul_r(fsub_S2, r_S2)));
498
499 #ifdef CALC_ENERGIES
500     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)));
501     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)));
502 #endif
503 #endif /* CALC_COUL_TAB */
504
505 #if defined CALC_ENERGIES && (defined CALC_COUL_EWALD || defined CALC_COUL_TAB)
506 #ifndef NO_SHIFT_EWALD
507     /* Add Ewald potential shift to vc_sub for convenience */
508 #ifdef CHECK_EXCLS
509     vc_sub_S0   = gmx_simd_add_r(vc_sub_S0, gmx_simd_blendzero_r(sh_ewald_S, interact_S0));
510     vc_sub_S2   = gmx_simd_add_r(vc_sub_S2, gmx_simd_blendzero_r(sh_ewald_S, interact_S2));
511 #else
512     vc_sub_S0   = gmx_simd_add_r(vc_sub_S0, sh_ewald_S);
513     vc_sub_S2   = gmx_simd_add_r(vc_sub_S2, sh_ewald_S);
514 #endif
515 #endif
516
517     vcoul_S0    = gmx_simd_mul_r(qq_S0, gmx_simd_sub_r(rinv_ex_S0, vc_sub_S0));
518     vcoul_S2    = gmx_simd_mul_r(qq_S2, gmx_simd_sub_r(rinv_ex_S2, vc_sub_S2));
519 #endif
520
521 #ifdef CALC_ENERGIES
522     /* Mask energy for cut-off and diagonal */
523     vcoul_S0    = gmx_simd_blendzero_r(vcoul_S0, wco_S0);
524     vcoul_S2    = gmx_simd_blendzero_r(vcoul_S2, wco_S2);
525 #endif
526
527 #endif /* CALC_COULOMB */
528
529 #ifdef CALC_LJ
530     /* Lennard-Jones interaction */
531
532 #ifdef VDW_CUTOFF_CHECK
533     wco_vdw_S0  = gmx_simd_cmplt_r(rsq_S0, rcvdw2_S);
534 #ifndef HALF_LJ
535     wco_vdw_S2  = gmx_simd_cmplt_r(rsq_S2, rcvdw2_S);
536 #endif
537 #else
538     /* Same cut-off for Coulomb and VdW, reuse the registers */
539 #define     wco_vdw_S0    wco_S0
540 #define     wco_vdw_S2    wco_S2
541 #endif
542
543 #ifndef LJ_COMB_LB
544     rinvsix_S0  = gmx_simd_mul_r(rinvsq_S0, gmx_simd_mul_r(rinvsq_S0, rinvsq_S0));
545 #ifdef EXCL_FORCES
546     rinvsix_S0  = gmx_simd_blendzero_r(rinvsix_S0, interact_S0);
547 #endif
548 #ifndef HALF_LJ
549     rinvsix_S2  = gmx_simd_mul_r(rinvsq_S2, gmx_simd_mul_r(rinvsq_S2, rinvsq_S2));
550 #ifdef EXCL_FORCES
551     rinvsix_S2  = gmx_simd_blendzero_r(rinvsix_S2, interact_S2);
552 #endif
553 #endif
554
555 #if defined LJ_CUT || defined LJ_POT_SWITCH
556     /* We have plain LJ or LJ-PME with simple C6/6 C12/12 coefficients */
557     FrLJ6_S0    = gmx_simd_mul_r(c6_S0, rinvsix_S0);
558 #ifndef HALF_LJ
559     FrLJ6_S2    = gmx_simd_mul_r(c6_S2, rinvsix_S2);
560 #endif
561     FrLJ12_S0   = gmx_simd_mul_r(c12_S0, gmx_simd_mul_r(rinvsix_S0, rinvsix_S0));
562 #ifndef HALF_LJ
563     FrLJ12_S2   = gmx_simd_mul_r(c12_S2, gmx_simd_mul_r(rinvsix_S2, rinvsix_S2));
564 #endif
565 #endif
566
567 #if defined LJ_FORCE_SWITCH || defined LJ_POT_SWITCH
568     /* We switch the LJ force */
569     r_S0        = gmx_simd_mul_r(rsq_S0, rinv_S0);
570     rsw_S0      = gmx_simd_max_r(gmx_simd_sub_r(r_S0, rswitch_S), zero_S);
571     rsw2_S0     = gmx_simd_mul_r(rsw_S0, rsw_S0);
572     rsw2_r_S0   = gmx_simd_mul_r(rsw2_S0, r_S0);
573 #ifndef HALF_LJ
574     r_S2        = gmx_simd_mul_r(rsq_S2, rinv_S2);
575     rsw_S2      = gmx_simd_max_r(gmx_simd_sub_r(r_S2, rswitch_S), zero_S);
576     rsw2_S2     = gmx_simd_mul_r(rsw_S2, rsw_S2);
577     rsw2_r_S2   = gmx_simd_mul_r(rsw2_S2, r_S2);
578 #endif
579 #endif
580
581 #ifdef LJ_FORCE_SWITCH
582
583 #define add_fr_switch(fr, rsw, rsw2_r, c2, c3) gmx_simd_fmadd_r(gmx_simd_fmadd_r(c3, rsw, c2), rsw2_r, fr)
584
585     FrLJ6_S0    = gmx_simd_mul_r(c6_S0, add_fr_switch(rinvsix_S0, rsw_S0, rsw2_r_S0, p6_fc2_S, p6_fc3_S));
586 #ifndef HALF_LJ
587     FrLJ6_S2    = gmx_simd_mul_r(c6_S2, add_fr_switch(rinvsix_S2, rsw_S2, rsw2_r_S2, p6_fc2_S, p6_fc3_S));
588 #endif
589     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));
590 #ifndef HALF_LJ
591     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));
592 #endif
593 #undef add_fr_switch
594 #endif /* LJ_FORCE_SWITCH */
595
596 #endif /* not LJ_COMB_LB */
597
598 #ifdef LJ_COMB_LB
599     sir_S0      = gmx_simd_mul_r(sig_S0, rinv_S0);
600 #ifndef HALF_LJ
601     sir_S2      = gmx_simd_mul_r(sig_S2, rinv_S2);
602 #endif
603     sir2_S0     = gmx_simd_mul_r(sir_S0, sir_S0);
604 #ifndef HALF_LJ
605     sir2_S2     = gmx_simd_mul_r(sir_S2, sir_S2);
606 #endif
607     sir6_S0     = gmx_simd_mul_r(sir2_S0, gmx_simd_mul_r(sir2_S0, sir2_S0));
608 #ifdef EXCL_FORCES
609     sir6_S0     = gmx_simd_blendzero_r(sir6_S0, interact_S0);
610 #endif
611 #ifndef HALF_LJ
612     sir6_S2     = gmx_simd_mul_r(sir2_S2, gmx_simd_mul_r(sir2_S2, sir2_S2));
613 #ifdef EXCL_FORCES
614     sir6_S2     = gmx_simd_blendzero_r(sir6_S2, interact_S2);
615 #endif
616 #endif
617 #ifdef VDW_CUTOFF_CHECK
618     sir6_S0     = gmx_simd_blendzero_r(sir6_S0, wco_vdw_S0);
619 #ifndef HALF_LJ
620     sir6_S2     = gmx_simd_blendzero_r(sir6_S2, wco_vdw_S2);
621 #endif
622 #endif
623     FrLJ6_S0    = gmx_simd_mul_r(eps_S0, sir6_S0);
624 #ifndef HALF_LJ
625     FrLJ6_S2    = gmx_simd_mul_r(eps_S2, sir6_S2);
626 #endif
627     FrLJ12_S0   = gmx_simd_mul_r(FrLJ6_S0, sir6_S0);
628 #ifndef HALF_LJ
629     FrLJ12_S2   = gmx_simd_mul_r(FrLJ6_S2, sir6_S2);
630 #endif
631 #if defined CALC_ENERGIES
632     /* We need C6 and C12 to calculate the LJ potential shift */
633     sig2_S0     = gmx_simd_mul_r(sig_S0, sig_S0);
634 #ifndef HALF_LJ
635     sig2_S2     = gmx_simd_mul_r(sig_S2, sig_S2);
636 #endif
637     sig6_S0     = gmx_simd_mul_r(sig2_S0, gmx_simd_mul_r(sig2_S0, sig2_S0));
638 #ifndef HALF_LJ
639     sig6_S2     = gmx_simd_mul_r(sig2_S2, gmx_simd_mul_r(sig2_S2, sig2_S2));
640 #endif
641     c6_S0       = gmx_simd_mul_r(eps_S0, sig6_S0);
642 #ifndef HALF_LJ
643     c6_S2       = gmx_simd_mul_r(eps_S2, sig6_S2);
644 #endif
645     c12_S0      = gmx_simd_mul_r(c6_S0, sig6_S0);
646 #ifndef HALF_LJ
647     c12_S2      = gmx_simd_mul_r(c6_S2, sig6_S2);
648 #endif
649 #endif
650 #endif /* LJ_COMB_LB */
651
652     /* Determine the total scalar LJ force*r */
653     frLJ_S0     = gmx_simd_sub_r(FrLJ12_S0, FrLJ6_S0);
654 #ifndef HALF_LJ
655     frLJ_S2     = gmx_simd_sub_r(FrLJ12_S2, FrLJ6_S2);
656 #endif
657
658 #if (defined LJ_CUT || defined LJ_FORCE_SWITCH) && defined CALC_ENERGIES
659
660 #ifdef LJ_CUT
661     /* Calculate the LJ energies, with constant potential shift */
662     VLJ6_S0     = gmx_simd_mul_r(sixth_S, gmx_simd_fmadd_r(c6_S0, p6_cpot_S, FrLJ6_S0));
663 #ifndef HALF_LJ
664     VLJ6_S2     = gmx_simd_mul_r(sixth_S, gmx_simd_fmadd_r(c6_S2, p6_cpot_S, FrLJ6_S2));
665 #endif
666     VLJ12_S0    = gmx_simd_mul_r(twelveth_S, gmx_simd_fmadd_r(c12_S0, p12_cpot_S, FrLJ12_S0));
667 #ifndef HALF_LJ
668     VLJ12_S2    = gmx_simd_mul_r(twelveth_S, gmx_simd_fmadd_r(c12_S2, p12_cpot_S, FrLJ12_S2));
669 #endif
670 #endif /* LJ_CUT */
671
672 #ifdef LJ_FORCE_SWITCH
673 #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)
674
675     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)));
676 #ifndef HALF_LJ
677     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)));
678 #endif
679     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)));
680 #ifndef HALF_LJ
681     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)));
682 #endif
683 #undef v_fswitch_pr
684 #endif /* LJ_FORCE_SWITCH */
685
686     /* Add up the repulsion and dispersion */
687     VLJ_S0      = gmx_simd_sub_r(VLJ12_S0, VLJ6_S0);
688 #ifndef HALF_LJ
689     VLJ_S2      = gmx_simd_sub_r(VLJ12_S2, VLJ6_S2);
690 #endif
691
692 #endif /* (LJ_CUT || LJ_FORCE_SWITCH) && CALC_ENERGIES */
693
694 #ifdef LJ_POT_SWITCH
695     /* We always need the potential, since it is needed for the force */
696     VLJ_S0 = gmx_simd_fnmadd_r(sixth_S, FrLJ6_S0, gmx_simd_mul_r(twelveth_S, FrLJ12_S0));
697 #ifndef HALF_LJ
698     VLJ_S2 = gmx_simd_fnmadd_r(sixth_S, FrLJ6_S2, gmx_simd_mul_r(twelveth_S, FrLJ12_S2));
699 #endif
700
701     {
702         gmx_simd_real_t sw_S0, dsw_S0;
703 #ifndef HALF_LJ
704         gmx_simd_real_t sw_S2, dsw_S2;
705 #endif
706
707 #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)
708 #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)
709
710         sw_S0  = switch_pr(rsw_S0, rsw2_S0, swV3_S, swV4_S, swV5_S);
711         dsw_S0 = dswitch_pr(rsw_S0, rsw2_S0, swF2_S, swF3_S, swF4_S);
712 #ifndef HALF_LJ
713         sw_S2  = switch_pr(rsw_S2, rsw2_S2, swV3_S, swV4_S, swV5_S);
714         dsw_S2 = dswitch_pr(rsw_S2, rsw2_S2, swF2_S, swF3_S, swF4_S);
715 #endif
716         frLJ_S0 = gmx_simd_fnmadd_r(gmx_simd_mul_r(dsw_S0, VLJ_S0), r_S0, gmx_simd_mul_r(sw_S0, frLJ_S0));
717 #ifndef HALF_LJ
718         frLJ_S2 = gmx_simd_fnmadd_r(gmx_simd_mul_r(dsw_S2, VLJ_S2), r_S2, gmx_simd_mul_r(sw_S2, frLJ_S2));
719 #endif
720 #ifdef CALC_ENERGIES
721         VLJ_S0  = gmx_simd_mul_r(sw_S0, VLJ_S0);
722 #ifndef HALF_LJ
723         VLJ_S2  = gmx_simd_mul_r(sw_S2, VLJ_S2);
724 #endif
725 #endif
726
727 #undef switch_pr
728 #undef dswitch_pr
729     }
730 #endif /* LJ_POT_SWITCH */
731
732 #if defined CALC_ENERGIES && defined CHECK_EXCLS
733     /* The potential shift should be removed for excluded pairs */
734     VLJ_S0      = gmx_simd_blendzero_r(VLJ_S0, interact_S0);
735 #ifndef HALF_LJ
736     VLJ_S2      = gmx_simd_blendzero_r(VLJ_S2, interact_S2);
737 #endif
738 #endif
739
740 #ifdef LJ_EWALD_GEOM
741     {
742         gmx_simd_real_t c6s_j_S;
743         gmx_simd_real_t c6grid_S0, rinvsix_nm_S0, cr2_S0, expmcr2_S0, poly_S0;
744 #ifndef HALF_LJ
745         gmx_simd_real_t c6grid_S2, rinvsix_nm_S2, cr2_S2, expmcr2_S2, poly_S2;
746 #endif
747 #ifdef CALC_ENERGIES
748         gmx_simd_real_t sh_mask_S0;
749 #ifndef HALF_LJ
750         gmx_simd_real_t sh_mask_S2;
751 #endif
752 #endif
753
754         /* Determine C6 for the grid using the geometric combination rule */
755         gmx_loaddh_pr(&c6s_j_S,  ljc+aj2+0);
756         c6grid_S0       = gmx_simd_mul_r(c6s_S0, c6s_j_S);
757 #ifndef HALF_LJ
758         c6grid_S2       = gmx_simd_mul_r(c6s_S2, c6s_j_S);
759 #endif
760
761 #ifdef CHECK_EXCLS
762         /* Recalculate rinvsix without exclusion mask (compiler might optimize) */
763         rinvsix_nm_S0 = gmx_simd_mul_r(rinvsq_S0, gmx_simd_mul_r(rinvsq_S0, rinvsq_S0));
764 #ifndef HALF_LJ
765         rinvsix_nm_S2 = gmx_simd_mul_r(rinvsq_S2, gmx_simd_mul_r(rinvsq_S2, rinvsq_S2));
766 #endif
767 #else
768         /* We didn't use a mask, so we can copy */
769         rinvsix_nm_S0 = rinvsix_S0;
770 #ifndef HALF_LJ
771         rinvsix_nm_S2 = rinvsix_S2;
772 #endif
773 #endif
774
775         /* Mask for the cut-off to avoid overflow of cr2^2 */
776         cr2_S0        = gmx_simd_mul_r(lje_c2_S, gmx_simd_blendzero_r(rsq_S0, wco_vdw_S0));
777 #ifndef HALF_LJ
778         cr2_S2        = gmx_simd_mul_r(lje_c2_S, gmx_simd_blendzero_r(rsq_S2, wco_vdw_S2));
779 #endif
780         expmcr2_S0    = gmx_simd_exp_r(gmx_simd_mul_r(mone_S, cr2_S0));
781 #ifndef HALF_LJ
782         expmcr2_S2    = gmx_simd_exp_r(gmx_simd_mul_r(mone_S, cr2_S2));
783 #endif
784
785         /* 1 + cr2 + 1/2*cr2^2 */
786         poly_S0       = gmx_simd_fmadd_r(gmx_simd_fmadd_r(half_S, cr2_S0, one_S), cr2_S0, one_S);
787 #ifndef HALF_LJ
788         poly_S2       = gmx_simd_fmadd_r(gmx_simd_fmadd_r(half_S, cr2_S2, one_S), cr2_S2, one_S);
789 #endif
790
791         /* We calculate LJ F*r = (6*C6)*(r^-6 - F_mesh/6), we use:
792          * r^-6*cexp*(1 + cr2 + cr2^2/2 + cr2^3/6) = cexp*(r^-6*poly + c^6/6)
793          */
794         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);
795 #ifndef HALF_LJ
796         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);
797 #endif
798
799 #ifdef CALC_ENERGIES
800 #ifdef CHECK_EXCLS
801         sh_mask_S0    = gmx_simd_blendzero_r(lje_vc_S, interact_S0);
802 #ifndef HALF_LJ
803         sh_mask_S2    = gmx_simd_blendzero_r(lje_vc_S, interact_S2);
804 #endif
805 #else
806         sh_mask_S0    = lje_vc_S;
807 #ifndef HALF_LJ
808         sh_mask_S2    = lje_vc_S;
809 #endif
810 #endif
811
812         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);
813 #ifndef HALF_LJ
814         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);
815 #endif
816 #endif /* CALC_ENERGIES */
817     }
818 #endif /* LJ_EWALD_GEOM */
819
820 #if defined VDW_CUTOFF_CHECK
821     /* frLJ is multiplied later by rinvsq, which is masked for the Coulomb
822      * cut-off, but if the VdW cut-off is shorter, we need to mask with that.
823      */
824     frLJ_S0     = gmx_simd_blendzero_r(frLJ_S0, wco_vdw_S0);
825 #ifndef HALF_LJ
826     frLJ_S2     = gmx_simd_blendzero_r(frLJ_S2, wco_vdw_S2);
827 #endif
828 #endif
829
830 #ifdef CALC_ENERGIES
831     /* The potential shift should be removed for pairs beyond cut-off */
832     VLJ_S0      = gmx_simd_blendzero_r(VLJ_S0, wco_vdw_S0);
833 #ifndef HALF_LJ
834     VLJ_S2      = gmx_simd_blendzero_r(VLJ_S2, wco_vdw_S2);
835 #endif
836 #endif
837
838 #endif /* CALC_LJ */
839
840 #ifdef CALC_ENERGIES
841 #ifdef ENERGY_GROUPS
842     /* Extract the group pair index per j pair.
843      * Energy groups are stored per i-cluster, so things get
844      * complicated when the i- and j-cluster size don't match.
845      */
846     {
847         int egps_j;
848 #if UNROLLJ == 2
849         egps_j    = nbat->energrp[cj>>1];
850         egp_jj[0] = ((egps_j >> ((cj & 1)*egps_jshift)) & egps_jmask)*egps_jstride;
851 #else
852         /* We assume UNROLLI <= UNROLLJ */
853         int jdi;
854         for (jdi = 0; jdi < UNROLLJ/UNROLLI; jdi++)
855         {
856             int jj;
857             egps_j = nbat->energrp[cj*(UNROLLJ/UNROLLI)+jdi];
858             for (jj = 0; jj < (UNROLLI/2); jj++)
859             {
860                 egp_jj[jdi*(UNROLLI/2)+jj] = ((egps_j >> (jj*egps_jshift)) & egps_jmask)*egps_jstride;
861             }
862         }
863 #endif
864     }
865 #endif
866
867 #ifdef CALC_COULOMB
868 #ifndef ENERGY_GROUPS
869     vctot_S      = gmx_simd_add_r(vctot_S, gmx_simd_add_r(vcoul_S0, vcoul_S2));
870 #else
871     add_ener_grp_halves(vcoul_S0, vctp[0], vctp[1], egp_jj);
872     add_ener_grp_halves(vcoul_S2, vctp[2], vctp[3], egp_jj);
873 #endif
874 #endif
875
876 #ifdef CALC_LJ
877 #ifndef ENERGY_GROUPS
878     Vvdwtot_S    = gmx_simd_add_r(Vvdwtot_S,
879 #ifndef HALF_LJ
880                                   gmx_simd_add_r(VLJ_S0, VLJ_S2)
881 #else
882                                   VLJ_S0
883 #endif
884                                   );
885 #else
886     add_ener_grp_halves(VLJ_S0, vvdwtp[0], vvdwtp[1], egp_jj);
887 #ifndef HALF_LJ
888     add_ener_grp_halves(VLJ_S2, vvdwtp[2], vvdwtp[3], egp_jj);
889 #endif
890 #endif
891 #endif /* CALC_LJ */
892 #endif /* CALC_ENERGIES */
893
894 #ifdef CALC_LJ
895 #ifdef CALC_COULOMB
896     fscal_S0    = gmx_simd_mul_r(rinvsq_S0, gmx_simd_add_r(frcoul_S0, frLJ_S0));
897 #else
898     fscal_S0    = gmx_simd_mul_r(rinvsq_S0, frLJ_S0);
899 #endif
900 #else
901     fscal_S0    = gmx_simd_mul_r(rinvsq_S0, frcoul_S0);
902 #endif /* CALC_LJ */
903 #if defined CALC_LJ && !defined HALF_LJ
904 #ifdef CALC_COULOMB
905     fscal_S2    = gmx_simd_mul_r(rinvsq_S2, gmx_simd_add_r(frcoul_S2, frLJ_S2));
906 #else
907     fscal_S2    = gmx_simd_mul_r(rinvsq_S2, frLJ_S2);
908 #endif
909 #else
910     /* Atom 2 and 3 don't have LJ, so only add Coulomb forces */
911     fscal_S2    = gmx_simd_mul_r(rinvsq_S2, frcoul_S2);
912 #endif
913
914     /* Calculate temporary vectorial force */
915     tx_S0       = gmx_simd_mul_r(fscal_S0, dx_S0);
916     tx_S2       = gmx_simd_mul_r(fscal_S2, dx_S2);
917     ty_S0       = gmx_simd_mul_r(fscal_S0, dy_S0);
918     ty_S2       = gmx_simd_mul_r(fscal_S2, dy_S2);
919     tz_S0       = gmx_simd_mul_r(fscal_S0, dz_S0);
920     tz_S2       = gmx_simd_mul_r(fscal_S2, dz_S2);
921
922     /* Increment i atom force */
923     fix_S0      = gmx_simd_add_r(fix_S0, tx_S0);
924     fix_S2      = gmx_simd_add_r(fix_S2, tx_S2);
925     fiy_S0      = gmx_simd_add_r(fiy_S0, ty_S0);
926     fiy_S2      = gmx_simd_add_r(fiy_S2, ty_S2);
927     fiz_S0      = gmx_simd_add_r(fiz_S0, tz_S0);
928     fiz_S2      = gmx_simd_add_r(fiz_S2, tz_S2);
929
930     /* Decrement j atom force */
931     gmx_load_hpr(&fjx_S, f+ajx);
932     gmx_load_hpr(&fjy_S, f+ajy);
933     gmx_load_hpr(&fjz_S, f+ajz);
934     gmx_store_hpr(f+ajx, gmx_sub_hpr(fjx_S, gmx_sum4_hpr(tx_S0, tx_S2)));
935     gmx_store_hpr(f+ajy, gmx_sub_hpr(fjy_S, gmx_sum4_hpr(ty_S0, ty_S2)));
936     gmx_store_hpr(f+ajz, gmx_sub_hpr(fjz_S, gmx_sum4_hpr(tz_S0, tz_S2)));
937 }
938
939 #undef  rinv_ex_S0
940 #undef  rinv_ex_S2
941
942 #undef  wco_vdw_S0
943 #undef  wco_vdw_S2
944
945 #undef  NBNXN_CUTOFF_USE_BLENDV
946
947 #undef  EXCL_FORCES