ced88fcf03fe8f0814873f1012a2a1bcfa28c169
[alexxy/gromacs.git] / src / gromacs / mdlib / nbnxn_kernels / simd_2xnn / nbnxn_kernel_simd_2xnn_outer.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2013,2014,2015,2016,2017,2018,2019, 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
37 {
38     using namespace gmx;
39     const nbnxn_cj_t   *l_cj;
40     const real         *q;
41     const real         *shiftvec;
42     const real         *x;
43     real                facel;
44     int                 ci, ci_sh;
45     int                 ish, ish3;
46     gmx_bool            do_LJ, half_LJ, do_coul;
47     int                 cjind0, cjind1, cjind;
48
49 #ifdef ENERGY_GROUPS
50     int         Vstride_i;
51     int         egps_ishift, egps_imask;
52     int         egps_jshift, egps_jmask, egps_jstride;
53     int         egps_i;
54     real       *vvdwtp[UNROLLI];
55     real       *vctp[UNROLLI];
56 #endif
57
58     SimdReal  shX_S;
59     SimdReal  shY_S;
60     SimdReal  shZ_S;
61     SimdReal  ix_S0, iy_S0, iz_S0;
62     SimdReal  ix_S2, iy_S2, iz_S2;
63     SimdReal  fix_S0, fiy_S0, fiz_S0;
64     SimdReal  fix_S2, fiy_S2, fiz_S2;
65
66     SimdReal  diagonal_jmi_S;
67 #if UNROLLI == UNROLLJ
68     SimdBool  diagonal_mask_S0, diagonal_mask_S2;
69 #else
70     SimdBool  diagonal_mask0_S0, diagonal_mask0_S2;
71     SimdBool  diagonal_mask1_S0, diagonal_mask1_S2;
72 #endif
73
74     unsigned            *exclusion_filter;
75     SimdBitMask          filter_S0, filter_S2;
76
77     SimdReal             zero_S(0.0);
78
79     SimdReal             one_S(1.0);
80     SimdReal             iq_S0 = setZero();
81     SimdReal             iq_S2 = setZero();
82
83 #ifdef CALC_COUL_RF
84     SimdReal      mrc_3_S;
85 #ifdef CALC_ENERGIES
86     SimdReal      hrc_3_S, moh_rc_S;
87 #endif
88 #endif
89
90 #ifdef CALC_COUL_TAB
91     /* Coulomb table variables */
92     SimdReal          invtsp_S;
93     const real       *tab_coul_F;
94 #if defined CALC_ENERGIES && !defined TAB_FDV0
95     const real       *tab_coul_V;
96 #endif
97
98 #ifdef CALC_ENERGIES
99     SimdReal   mhalfsp_S;
100 #endif
101 #endif
102
103 #ifdef CALC_COUL_EWALD
104     SimdReal beta2_S, beta_S;
105 #endif
106
107 #if defined CALC_ENERGIES && (defined CALC_COUL_EWALD || defined CALC_COUL_TAB)
108     SimdReal  sh_ewald_S;
109 #endif
110
111 #if defined LJ_CUT && defined CALC_ENERGIES
112     SimdReal   p6_cpot_S, p12_cpot_S;
113 #endif
114 #ifdef LJ_POT_SWITCH
115     SimdReal   rswitch_S;
116     SimdReal   swV3_S, swV4_S, swV5_S;
117     SimdReal   swF2_S, swF3_S, swF4_S;
118 #endif
119 #ifdef LJ_FORCE_SWITCH
120     SimdReal   rswitch_S;
121     SimdReal   p6_fc2_S, p6_fc3_S;
122     SimdReal   p12_fc2_S, p12_fc3_S;
123 #ifdef CALC_ENERGIES
124     SimdReal   p6_vc3_S, p6_vc4_S;
125     SimdReal   p12_vc3_S, p12_vc4_S;
126     SimdReal   p6_6cpot_S, p12_12cpot_S;
127 #endif
128 #endif
129 #ifdef LJ_EWALD_GEOM
130     real              lj_ewaldcoeff2, lj_ewaldcoeff6_6;
131     SimdReal          mone_S, half_S, lje_c2_S, lje_c6_6_S;
132 #endif
133
134 #ifdef LJ_COMB_LB
135     const real       *ljc;
136
137     SimdReal          hsig_i_S0, seps_i_S0;
138     SimdReal          hsig_i_S2, seps_i_S2;
139 #else
140 #ifdef FIX_LJ_C
141     alignas(GMX_SIMD_ALIGNMENT) real  pvdw_c6[2*UNROLLI*UNROLLJ];
142     real  *pvdw_c12 = pvdw_c6 + UNROLLI*UNROLLJ;
143 #endif
144
145 #if defined LJ_COMB_GEOM || defined LJ_EWALD_GEOM
146     const real       *ljc;
147 #endif
148 #endif /* LJ_COMB_LB */
149
150     SimdReal  minRsq_S;
151     SimdReal  rc2_S;
152 #ifdef VDW_CUTOFF_CHECK
153     SimdReal  rcvdw2_S;
154 #endif
155
156     int ninner;
157
158 #ifdef COUNT_PAIRS
159     int npair = 0;
160 #endif
161
162 #if defined LJ_COMB_GEOM || defined LJ_COMB_LB || defined LJ_EWALD_GEOM
163     ljc = nbat->lj_comb;
164 #endif
165 #if !(defined LJ_COMB_GEOM || defined LJ_COMB_LB || defined FIX_LJ_C)
166     /* No combination rule used */
167     real      *nbfp_ptr = nbat->nbfp_aligned;
168     const int *type     = nbat->type;
169 #endif
170
171     /* Load j-i for the first i */
172     diagonal_jmi_S    = load<SimdReal>(nbat->simd_2xnn_diagonal_j_minus_i);
173     /* Generate all the diagonal masks as comparison results */
174 #if UNROLLI == UNROLLJ
175     diagonal_mask_S0  = (zero_S < diagonal_jmi_S);
176     diagonal_jmi_S    = diagonal_jmi_S - one_S;
177     diagonal_jmi_S    = diagonal_jmi_S - one_S;
178     diagonal_mask_S2  = (zero_S < diagonal_jmi_S);
179 #else
180 #if 2*UNROLLI == UNROLLJ
181     diagonal_mask0_S0 = (zero_S < diagonal_jmi_S);
182     diagonal_jmi_S    = diagonal_jmi_S - one_S;
183     diagonal_jmi_S    = diagonal_jmi_S - one_S;
184     diagonal_mask0_S2 = (zero_S < diagonal_jmi_S);
185     diagonal_jmi_S    = diagonal_jmi_S - one_S;
186     diagonal_jmi_S    = diagonal_jmi_S - one_S;
187     diagonal_mask1_S0 = (zero_S < diagonal_jmi_S);
188     diagonal_jmi_S    = diagonal_jmi_S - one_S;
189     diagonal_jmi_S    = diagonal_jmi_S - one_S;
190     diagonal_mask1_S2 = (zero_S < diagonal_jmi_S);
191 #endif
192 #endif
193
194     /* Load masks for topology exclusion masking. filter_stride is
195        static const, so the conditional will be optimized away. */
196 #if GMX_DOUBLE && !GMX_SIMD_HAVE_INT32_LOGICAL
197     exclusion_filter = nbat->simd_exclusion_filter64;
198 #else
199     exclusion_filter = nbat->simd_exclusion_filter;
200 #endif
201
202     /* Here we cast the exclusion filters from unsigned * to int * or real *.
203      * Since we only check bits, the actual value they represent does not
204      * matter, as long as both filter and mask data are treated the same way.
205      */
206 #if GMX_SIMD_HAVE_INT32_LOGICAL
207     filter_S0 = load<SimdBitMask>(reinterpret_cast<const int *>(exclusion_filter + 0*UNROLLJ));
208     filter_S2 = load<SimdBitMask>(reinterpret_cast<const int *>(exclusion_filter + 2*UNROLLJ));
209 #else
210     filter_S0 = load<SimdBitMask>(reinterpret_cast<const real *>(exclusion_filter + 0*UNROLLJ));
211     filter_S2 = load<SimdBitMask>(reinterpret_cast<const real *>(exclusion_filter + 2*UNROLLJ));
212 #endif
213
214 #ifdef CALC_COUL_RF
215     /* Reaction-field constants */
216     mrc_3_S  = SimdReal(-2*ic->k_rf);
217 #ifdef CALC_ENERGIES
218     hrc_3_S  = SimdReal(ic->k_rf);
219     moh_rc_S = SimdReal(-ic->c_rf);
220 #endif
221 #endif
222
223 #ifdef CALC_COUL_TAB
224
225     invtsp_S  = SimdReal(ic->tabq_scale);
226 #ifdef CALC_ENERGIES
227     mhalfsp_S = SimdReal(-0.5/ic->tabq_scale);
228 #endif
229
230 #ifdef TAB_FDV0
231     tab_coul_F = ic->tabq_coul_FDV0;
232 #else
233     tab_coul_F = ic->tabq_coul_F;
234 #ifdef CALC_ENERGIES
235     tab_coul_V = ic->tabq_coul_V;
236 #endif
237 #endif
238 #endif /* CALC_COUL_TAB */
239
240 #ifdef CALC_COUL_EWALD
241     beta2_S = SimdReal(ic->ewaldcoeff_q*ic->ewaldcoeff_q);
242     beta_S  = SimdReal(ic->ewaldcoeff_q);
243 #endif
244
245 #if (defined CALC_COUL_TAB || defined CALC_COUL_EWALD) && defined CALC_ENERGIES
246     sh_ewald_S = SimdReal(ic->sh_ewald);
247 #endif
248
249     /* LJ function constants */
250 #if defined CALC_ENERGIES || defined LJ_POT_SWITCH
251     SimdReal sixth_S      = SimdReal(1.0/6.0);
252     SimdReal twelveth_S   = SimdReal(1.0/12.0);
253 #endif
254
255 #if defined LJ_CUT && defined CALC_ENERGIES
256     /* We shift the potential by cpot, which can be zero */
257     p6_cpot_S    = SimdReal(ic->dispersion_shift.cpot);
258     p12_cpot_S   = SimdReal(ic->repulsion_shift.cpot);
259 #endif
260 #ifdef LJ_POT_SWITCH
261     rswitch_S = SimdReal(ic->rvdw_switch);
262     swV3_S    = SimdReal(ic->vdw_switch.c3);
263     swV4_S    = SimdReal(ic->vdw_switch.c4);
264     swV5_S    = SimdReal(ic->vdw_switch.c5);
265     swF2_S    = SimdReal(3*ic->vdw_switch.c3);
266     swF3_S    = SimdReal(4*ic->vdw_switch.c4);
267     swF4_S    = SimdReal(5*ic->vdw_switch.c5);
268 #endif
269 #ifdef LJ_FORCE_SWITCH
270     rswitch_S = SimdReal(ic->rvdw_switch);
271     p6_fc2_S  = SimdReal(ic->dispersion_shift.c2);
272     p6_fc3_S  = SimdReal(ic->dispersion_shift.c3);
273     p12_fc2_S = SimdReal(ic->repulsion_shift.c2);
274     p12_fc3_S = SimdReal(ic->repulsion_shift.c3);
275 #ifdef CALC_ENERGIES
276     {
277         SimdReal mthird_S  = SimdReal(-1.0/3.0);
278         SimdReal mfourth_S = SimdReal(-1.0/4.0);
279
280         p6_vc3_S     = mthird_S * p6_fc2_S;
281         p6_vc4_S     = mfourth_S * p6_fc3_S;
282         p6_6cpot_S   = SimdReal(ic->dispersion_shift.cpot/6);
283         p12_vc3_S    = mthird_S * p12_fc2_S;
284         p12_vc4_S    = mfourth_S * p12_fc3_S;
285         p12_12cpot_S = SimdReal(ic->repulsion_shift.cpot/12);
286     }
287 #endif
288 #endif
289 #ifdef LJ_EWALD_GEOM
290     mone_S           = SimdReal(-1.0);
291     half_S           = SimdReal(0.5);
292     lj_ewaldcoeff2   = ic->ewaldcoeff_lj*ic->ewaldcoeff_lj;
293     lj_ewaldcoeff6_6 = lj_ewaldcoeff2*lj_ewaldcoeff2*lj_ewaldcoeff2/6;
294     lje_c2_S         = SimdReal(lj_ewaldcoeff2);
295     lje_c6_6_S       = SimdReal(lj_ewaldcoeff6_6);
296 #ifdef CALC_ENERGIES
297     /* Determine the grid potential at the cut-off */
298     SimdReal lje_vc_S = SimdReal(ic->sh_lj_ewald);
299 #endif
300 #endif
301
302     /* The kernel either supports rcoulomb = rvdw or rcoulomb >= rvdw */
303     rc2_S    = SimdReal(ic->rcoulomb*ic->rcoulomb);
304 #ifdef VDW_CUTOFF_CHECK
305     rcvdw2_S = SimdReal(ic->rvdw*ic->rvdw);
306 #endif
307
308     minRsq_S            = SimdReal(NBNXN_MIN_RSQ);
309
310     q                   = nbat->q;
311     facel               = ic->epsfac;
312     shiftvec            = shift_vec[0];
313     x                   = nbat->x;
314
315 #ifdef FIX_LJ_C
316
317     for (jp = 0; jp < UNROLLJ; jp++)
318     {
319         pvdw_c6 [0*UNROLLJ+jp] = nbat->nbfp[0*2];
320         pvdw_c6 [1*UNROLLJ+jp] = nbat->nbfp[0*2];
321         pvdw_c6 [2*UNROLLJ+jp] = nbat->nbfp[0*2];
322         pvdw_c6 [3*UNROLLJ+jp] = nbat->nbfp[0*2];
323
324         pvdw_c12[0*UNROLLJ+jp] = nbat->nbfp[0*2+1];
325         pvdw_c12[1*UNROLLJ+jp] = nbat->nbfp[0*2+1];
326         pvdw_c12[2*UNROLLJ+jp] = nbat->nbfp[0*2+1];
327         pvdw_c12[3*UNROLLJ+jp] = nbat->nbfp[0*2+1];
328     }
329     SimdReal c6_S0  = load<SimdReal>(pvdw_c6 +0*UNROLLJ);
330     SimdReal c6_S1  = load<SimdReal>(pvdw_c6 +1*UNROLLJ);
331     SimdReal c6_S2  = load<SimdReal>(pvdw_c6 +2*UNROLLJ);
332     SimdReal c6_S3  = load<SimdReal>(pvdw_c6 +3*UNROLLJ);
333
334     SimdReal c12_S0 = load<SimdReal>(pvdw_c12+0*UNROLLJ);
335     SimdReal c12_S1 = load<SimdReal>(pvdw_c12+1*UNROLLJ);
336     SimdReal c12_S2 = load<SimdReal>(pvdw_c12+2*UNROLLJ);
337     SimdReal c12_S3 = load<SimdReal>(pvdw_c12+3*UNROLLJ);
338 #endif /* FIX_LJ_C */
339
340 #ifdef ENERGY_GROUPS
341     egps_ishift  = nbat->neg_2log;
342     egps_imask   = (1<<egps_ishift) - 1;
343     egps_jshift  = 2*nbat->neg_2log;
344     egps_jmask   = (1<<egps_jshift) - 1;
345     egps_jstride = (UNROLLJ>>1)*UNROLLJ;
346     /* Major division is over i-particle energy groups, determine the stride */
347     Vstride_i    = nbat->nenergrp*(1<<nbat->neg_2log)*egps_jstride;
348 #endif
349
350     l_cj = nbl->cj.data();
351
352     ninner = 0;
353     for (const nbnxn_ci_t &ciEntry : nbl->ci)
354     {
355         ish              = (ciEntry.shift & NBNXN_CI_SHIFT);
356         ish3             = ish*3;
357         cjind0           = ciEntry.cj_ind_start;
358         cjind1           = ciEntry.cj_ind_end;
359         ci               = ciEntry.ci;
360         ci_sh            = (ish == CENTRAL ? ci : -1);
361
362         shX_S = SimdReal(shiftvec[ish3]);
363         shY_S = SimdReal(shiftvec[ish3+1]);
364         shZ_S = SimdReal(shiftvec[ish3+2]);
365
366 #if UNROLLJ <= 4
367         int sci              = ci*STRIDE;
368         int scix             = sci*DIM;
369 #if defined LJ_COMB_LB || defined LJ_COMB_GEOM || defined LJ_EWALD_GEOM
370         int sci2             = sci*2;
371 #endif
372 #else
373         int sci              = (ci>>1)*STRIDE;
374         int scix             = sci*DIM + (ci & 1)*(STRIDE>>1);
375 #if defined LJ_COMB_LB || defined LJ_COMB_GEOM || defined LJ_EWALD_GEOM
376         int sci2             = sci*2 + (ci & 1)*(STRIDE>>1);
377 #endif
378         sci             += (ci & 1)*(STRIDE>>1);
379 #endif
380
381         /* We have 5 LJ/C combinations, but use only three inner loops,
382          * as the other combinations are unlikely and/or not much faster:
383          * inner half-LJ + C for half-LJ + C / no-LJ + C
384          * inner LJ + C      for full-LJ + C
385          * inner LJ          for full-LJ + no-C / half-LJ + no-C
386          */
387         do_LJ   = ((ciEntry.shift & NBNXN_CI_DO_LJ(0)) != 0);
388         do_coul = ((ciEntry.shift & NBNXN_CI_DO_COUL(0)) != 0);
389         half_LJ = (((ciEntry.shift & NBNXN_CI_HALF_LJ(0)) != 0) || !do_LJ) && do_coul;
390
391 #ifdef ENERGY_GROUPS
392         egps_i = nbat->energrp[ci];
393         {
394             int ia, egp_ia;
395
396             for (ia = 0; ia < UNROLLI; ia++)
397             {
398                 egp_ia     = (egps_i >> (ia*egps_ishift)) & egps_imask;
399                 vvdwtp[ia] = Vvdw + egp_ia*Vstride_i;
400                 vctp[ia]   = Vc   + egp_ia*Vstride_i;
401             }
402         }
403 #endif
404
405 #ifdef CALC_ENERGIES
406 #ifdef LJ_EWALD_GEOM
407         gmx_bool do_self = TRUE;
408 #else
409         gmx_bool do_self = do_coul;
410 #endif
411 #if UNROLLJ == 4
412         if (do_self && l_cj[ciEntry.cj_ind_start].cj == ci_sh)
413 #endif
414 #if UNROLLJ == 8
415         if (do_self && l_cj[ciEntry.cj_ind_start].cj == (ci_sh>>1))
416 #endif
417         {
418             if (do_coul)
419             {
420                 real Vc_sub_self;
421                 int  ia;
422
423 #ifdef CALC_COUL_RF
424                 Vc_sub_self = 0.5*ic->c_rf;
425 #endif
426 #ifdef CALC_COUL_TAB
427 #ifdef TAB_FDV0
428                 Vc_sub_self = 0.5*tab_coul_F[2];
429 #else
430                 Vc_sub_self = 0.5*tab_coul_V[0];
431 #endif
432 #endif
433 #ifdef CALC_COUL_EWALD
434                 /* beta/sqrt(pi) */
435                 Vc_sub_self = 0.5*ic->ewaldcoeff_q*M_2_SQRTPI;
436 #endif
437
438                 for (ia = 0; ia < UNROLLI; ia++)
439                 {
440                     real qi;
441
442                     qi = q[sci+ia];
443 #ifdef ENERGY_GROUPS
444                     vctp[ia][((egps_i>>(ia*egps_ishift)) & egps_imask)*egps_jstride]
445 #else
446                     Vc[0]
447 #endif
448                         -= facel*qi*qi*Vc_sub_self;
449                 }
450             }
451
452 #ifdef LJ_EWALD_GEOM
453             {
454                 int  ia;
455
456                 for (ia = 0; ia < UNROLLI; ia++)
457                 {
458                     real c6_i;
459
460                     c6_i = nbat->nbfp[nbat->type[sci+ia]*(nbat->ntype + 1)*2]/6;
461 #ifdef ENERGY_GROUPS
462                     vvdwtp[ia][((egps_i>>(ia*egps_ishift)) & egps_imask)*egps_jstride]
463 #else
464                     Vvdw[0]
465 #endif
466                         += 0.5*c6_i*lj_ewaldcoeff6_6;
467                 }
468             }
469 #endif      /* LJ_EWALD */
470         }
471 #endif
472
473         /* Load i atom data */
474         int sciy             = scix + STRIDE;
475         int sciz             = sciy + STRIDE;
476         ix_S0          = loadU1DualHsimd(x+scix);
477         ix_S2          = loadU1DualHsimd(x+scix+2);
478         iy_S0          = loadU1DualHsimd(x+sciy);
479         iy_S2          = loadU1DualHsimd(x+sciy+2);
480         iz_S0          = loadU1DualHsimd(x+sciz);
481         iz_S2          = loadU1DualHsimd(x+sciz+2);
482         ix_S0          = ix_S0 + shX_S;
483         ix_S2          = ix_S2 + shX_S;
484         iy_S0          = iy_S0 + shY_S;
485         iy_S2          = iy_S2 + shY_S;
486         iz_S0          = iz_S0 + shZ_S;
487         iz_S2          = iz_S2 + shZ_S;
488
489         if (do_coul)
490         {
491             SimdReal facel_S;
492
493             facel_S    = SimdReal(facel);
494
495             iq_S0      = loadU1DualHsimd(q+sci);
496             iq_S2      = loadU1DualHsimd(q+sci+2);
497             iq_S0      = facel_S * iq_S0;
498             iq_S2      = facel_S * iq_S2;
499         }
500
501 #ifdef LJ_COMB_LB
502         hsig_i_S0 = loadU1DualHsimd(ljc+sci2);
503         hsig_i_S2 = loadU1DualHsimd(ljc+sci2+2);
504         seps_i_S0 = loadU1DualHsimd(ljc+sci2+STRIDE);
505         seps_i_S2 = loadU1DualHsimd(ljc+sci2+STRIDE+2);
506 #else
507 #ifdef LJ_COMB_GEOM
508         SimdReal   c6s_S0, c12s_S0;
509         SimdReal   c6s_S2, c12s_S2;
510
511         c6s_S0 = loadU1DualHsimd(ljc+sci2);
512
513         if (!half_LJ)
514         {
515             c6s_S2 = loadU1DualHsimd(ljc+sci2+2);
516         }
517         c12s_S0 = loadU1DualHsimd(ljc+sci2+STRIDE);
518         if (!half_LJ)
519         {
520             c12s_S2 = loadU1DualHsimd(ljc+sci2+STRIDE+2);
521         }
522 #elif !defined LJ_COMB_LB && !defined FIX_LJ_C
523         const real *nbfp0     = nbfp_ptr + type[sci  ]*nbat->ntype*c_simdBestPairAlignment;
524         const real *nbfp1     = nbfp_ptr + type[sci+1]*nbat->ntype*c_simdBestPairAlignment;
525         const real *nbfp2     = nullptr, *nbfp3 = nullptr;
526         if (!half_LJ)
527         {
528             nbfp2 = nbfp_ptr + type[sci+2]*nbat->ntype*c_simdBestPairAlignment;
529             nbfp3 = nbfp_ptr + type[sci+3]*nbat->ntype*c_simdBestPairAlignment;
530         }
531 #endif
532 #endif
533 #ifdef LJ_EWALD_GEOM
534         /* We need the geometrically combined C6 for the PME grid correction */
535         SimdReal c6s_S0, c6s_S2;
536         c6s_S0 = loadU1DualHsimd(ljc+sci2);
537         if (!half_LJ)
538         {
539             c6s_S2 = loadU1DualHsimd(ljc+sci2+2);
540         }
541 #endif
542
543         /* Zero the potential energy for this list */
544 #ifdef CALC_ENERGIES
545         SimdReal Vvdwtot_S = setZero();
546         SimdReal vctot_S   = setZero();
547 #endif
548
549         /* Clear i atom forces */
550         fix_S0           = setZero();
551         fix_S2           = setZero();
552         fiy_S0           = setZero();
553         fiy_S2           = setZero();
554         fiz_S0           = setZero();
555         fiz_S2           = setZero();
556
557         cjind = cjind0;
558
559         /* Currently all kernels use (at least half) LJ */
560 #define CALC_LJ
561         if (half_LJ)
562         {
563             /* Coulomb: all i-atoms, LJ: first half i-atoms */
564 #define CALC_COULOMB
565 #define HALF_LJ
566 #define CHECK_EXCLS
567             while (cjind < cjind1 && nbl->cj[cjind].excl != NBNXN_INTERACTION_MASK_ALL)
568             {
569 #include "gromacs/mdlib/nbnxn_kernels/simd_2xnn/nbnxn_kernel_simd_2xnn_inner.h"
570                 cjind++;
571             }
572 #undef CHECK_EXCLS
573             for (; (cjind < cjind1); cjind++)
574             {
575 #include "gromacs/mdlib/nbnxn_kernels/simd_2xnn/nbnxn_kernel_simd_2xnn_inner.h"
576             }
577 #undef HALF_LJ
578 #undef CALC_COULOMB
579         }
580         else if (do_coul)
581         {
582             /* Coulomb: all i-atoms, LJ: all i-atoms */
583 #define CALC_COULOMB
584 #define CHECK_EXCLS
585             while (cjind < cjind1 && nbl->cj[cjind].excl != NBNXN_INTERACTION_MASK_ALL)
586             {
587 #include "gromacs/mdlib/nbnxn_kernels/simd_2xnn/nbnxn_kernel_simd_2xnn_inner.h"
588                 cjind++;
589             }
590 #undef CHECK_EXCLS
591             for (; (cjind < cjind1); cjind++)
592             {
593 #include "gromacs/mdlib/nbnxn_kernels/simd_2xnn/nbnxn_kernel_simd_2xnn_inner.h"
594             }
595 #undef CALC_COULOMB
596         }
597         else
598         {
599             /* Coulomb: none, LJ: all i-atoms */
600 #define CHECK_EXCLS
601             while (cjind < cjind1 && nbl->cj[cjind].excl != NBNXN_INTERACTION_MASK_ALL)
602             {
603 #include "gromacs/mdlib/nbnxn_kernels/simd_2xnn/nbnxn_kernel_simd_2xnn_inner.h"
604                 cjind++;
605             }
606 #undef CHECK_EXCLS
607             for (; (cjind < cjind1); cjind++)
608             {
609 #include "gromacs/mdlib/nbnxn_kernels/simd_2xnn/nbnxn_kernel_simd_2xnn_inner.h"
610             }
611         }
612 #undef CALC_LJ
613         ninner += cjind1 - cjind0;
614
615         /* Add accumulated i-forces to the force array */
616         real fShiftX = reduceIncr4ReturnSumHsimd(f+scix, fix_S0, fix_S2);
617         real fShiftY = reduceIncr4ReturnSumHsimd(f+sciy, fiy_S0, fiy_S2);
618         real fShiftZ = reduceIncr4ReturnSumHsimd(f+sciz, fiz_S0, fiz_S2);
619
620 #ifdef CALC_SHIFTFORCES
621         fshift[ish3+0] += fShiftX;
622         fshift[ish3+1] += fShiftY;
623         fshift[ish3+2] += fShiftZ;
624 #endif
625
626 #ifdef CALC_ENERGIES
627         if (do_coul)
628         {
629             *Vc += reduce(vctot_S);
630         }
631         *Vvdw += reduce(Vvdwtot_S);
632 #endif
633
634         /* Outer loop uses 6 flops/iteration */
635     }
636
637 #ifdef COUNT_PAIRS
638     printf("atom pairs %d\n", npair);
639 #endif
640 }