2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012-2018, The GROMACS development team.
5 * Copyright (c) 2019,2020,2021, by the GROMACS development team, led by
6 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7 * and including many others, as listed in the AUTHORS file in the
8 * top-level source directory and at http://www.gromacs.org.
10 * GROMACS is free software; you can redistribute it and/or
11 * modify it under the terms of the GNU Lesser General Public License
12 * as published by the Free Software Foundation; either version 2.1
13 * of the License, or (at your option) any later version.
15 * GROMACS is distributed in the hope that it will be useful,
16 * but WITHOUT ANY WARRANTY; without even the implied warranty of
17 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
18 * Lesser General Public License for more details.
20 * You should have received a copy of the GNU Lesser General Public
21 * License along with GROMACS; if not, see
22 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
25 * If you want to redistribute modifications to GROMACS, please
26 * consider that scientific software is very special. Version
27 * control is crucial - bugs must be traceable. We will be happy to
28 * consider code for inclusion in the official distribution, but
29 * derived work must not be called official GROMACS. Details are found
30 * in the README & COPYING files - if they are missing, get the
31 * official version at http://www.gromacs.org.
33 * To help us fund GROMACS development, we humbly ask that you cite
34 * the research papers on the package. Check out http://www.gromacs.org.
48 #include "thread_mpi/atomic.h"
50 #include "gromacs/math/functions.h"
51 #include "gromacs/math/utilities.h"
52 #include "gromacs/math/vec.h"
53 #include "gromacs/mdlib/gmx_omp_nthreads.h"
54 #include "gromacs/mdtypes/forcerec.h" // only for GET_CGINFO_*
55 #include "gromacs/mdtypes/md_enums.h"
56 #include "gromacs/nbnxm/nbnxm.h"
57 #include "gromacs/pbcutil/ishift.h"
58 #include "gromacs/simd/simd.h"
59 #include "gromacs/utility/exceptions.h"
60 #include "gromacs/utility/fatalerror.h"
61 #include "gromacs/utility/gmxomp.h"
62 #include "gromacs/utility/logger.h"
63 #include "gromacs/utility/smalloc.h"
64 #include "gromacs/utility/strconvert.h"
65 #include "gromacs/utility/stringutil.h"
69 #include "nbnxm_geometry.h"
70 #include "nbnxm_gpu.h"
73 using namespace gmx; // TODO: Remove when this file is moved into gmx namespace
75 const char* const c_ljcrNames[ljcrNR + 1] = { "none", "geometric", "Lorentz-Berthelot", nullptr };
77 void nbnxn_atomdata_t::resizeCoordinateBuffer(int numAtoms)
81 x_.resize(numAtoms * xstride);
84 void nbnxn_atomdata_t::resizeForceBuffers()
86 /* Force buffers need padding up to a multiple of the buffer flag size */
87 const int paddedSize =
88 (numAtoms() + NBNXN_BUFFERFLAG_SIZE - 1) / NBNXN_BUFFERFLAG_SIZE * NBNXN_BUFFERFLAG_SIZE;
90 /* Should we let each thread allocate it's own data instead? */
91 for (nbnxn_atomdata_output_t& outBuffer : out)
93 outBuffer.f.resize(paddedSize * fstride);
97 /* Initializes an nbnxn_atomdata_output_t data structure */
98 nbnxn_atomdata_output_t::nbnxn_atomdata_output_t(Nbnxm::KernelType kernelType,
100 int simdEnergyBufferStride,
101 gmx::PinningPolicy pinningPolicy) :
102 f({}, { pinningPolicy }),
103 fshift({}, { pinningPolicy }),
104 Vvdw({}, { pinningPolicy }),
105 Vc({}, { pinningPolicy })
107 fshift.resize(SHIFTS * DIM);
108 Vvdw.resize(numEnergyGroups * numEnergyGroups);
109 Vc.resize(numEnergyGroups * numEnergyGroups);
111 if (Nbnxm::kernelTypeIsSimd(kernelType))
113 int cj_size = Nbnxm::JClusterSizePerKernelType[kernelType];
115 numEnergyGroups * numEnergyGroups * simdEnergyBufferStride * (cj_size / 2) * cj_size;
116 VSvdw.resize(numElements);
117 VSc.resize(numElements);
121 static void copy_int_to_nbat_int(const int* a, int na, int na_round, const int* in, int fill, int* innb)
126 for (i = 0; i < na; i++)
128 innb[j++] = in[a[i]];
130 /* Complete the partially filled last cell with fill */
131 for (; i < na_round; i++)
137 void copy_rvec_to_nbat_real(const int* a, int na, int na_round, const rvec* x, int nbatFormat, real* xnb, int a0)
139 /* We complete partially filled cells, can only be the last one in each
140 * column, with coordinates farAway. The actual coordinate value does
141 * not influence the results, since these filler particles do not interact.
142 * Clusters with normal atoms + fillers have a bounding box based only
143 * on the coordinates of the atoms. Clusters with only fillers have as
144 * the bounding box the coordinates of the first filler. Such clusters
145 * are not considered as i-entries, but they are considered as j-entries.
146 * So for performance it is better to have their bounding boxes far away,
147 * such that filler only clusters don't end up in the pair list.
149 const real farAway = -1000000;
157 for (i = 0; i < na; i++)
159 xnb[j++] = x[a[i]][XX];
160 xnb[j++] = x[a[i]][YY];
161 xnb[j++] = x[a[i]][ZZ];
163 /* Complete the partially filled last cell with farAway elements */
164 for (; i < na_round; i++)
172 j = a0 * STRIDE_XYZQ;
173 for (i = 0; i < na; i++)
175 xnb[j++] = x[a[i]][XX];
176 xnb[j++] = x[a[i]][YY];
177 xnb[j++] = x[a[i]][ZZ];
180 /* Complete the partially filled last cell with zeros */
181 for (; i < na_round; i++)
190 j = atom_to_x_index<c_packX4>(a0);
191 c = a0 & (c_packX4 - 1);
192 for (i = 0; i < na; i++)
194 xnb[j + XX * c_packX4] = x[a[i]][XX];
195 xnb[j + YY * c_packX4] = x[a[i]][YY];
196 xnb[j + ZZ * c_packX4] = x[a[i]][ZZ];
201 j += (DIM - 1) * c_packX4;
205 /* Complete the partially filled last cell with zeros */
206 for (; i < na_round; i++)
208 xnb[j + XX * c_packX4] = farAway;
209 xnb[j + YY * c_packX4] = farAway;
210 xnb[j + ZZ * c_packX4] = farAway;
215 j += (DIM - 1) * c_packX4;
221 j = atom_to_x_index<c_packX8>(a0);
222 c = a0 & (c_packX8 - 1);
223 for (i = 0; i < na; i++)
225 xnb[j + XX * c_packX8] = x[a[i]][XX];
226 xnb[j + YY * c_packX8] = x[a[i]][YY];
227 xnb[j + ZZ * c_packX8] = x[a[i]][ZZ];
232 j += (DIM - 1) * c_packX8;
236 /* Complete the partially filled last cell with zeros */
237 for (; i < na_round; i++)
239 xnb[j + XX * c_packX8] = farAway;
240 xnb[j + YY * c_packX8] = farAway;
241 xnb[j + ZZ * c_packX8] = farAway;
246 j += (DIM - 1) * c_packX8;
251 default: gmx_incons("Unsupported nbnxn_atomdata_t format");
255 /* Stores the LJ parameter data in a format convenient for different kernels */
256 static void set_lj_parameter_data(nbnxn_atomdata_t::Params* params, gmx_bool bSIMD)
258 int nt = params->numTypes;
263 /* nbfp_aligned stores two parameters using the stride most suitable
264 * for the present SIMD architecture, as specified by the constant
265 * c_simdBestPairAlignment from the SIMD header.
266 * There's a slight inefficiency in allocating and initializing nbfp_aligned
267 * when it might not be used, but introducing the conditional code is not
270 params->nbfp_aligned.resize(nt * nt * c_simdBestPairAlignment);
272 for (int i = 0; i < nt; i++)
274 for (int j = 0; j < nt; j++)
276 params->nbfp_aligned[(i * nt + j) * c_simdBestPairAlignment + 0] =
277 params->nbfp[(i * nt + j) * 2 + 0];
278 params->nbfp_aligned[(i * nt + j) * c_simdBestPairAlignment + 1] =
279 params->nbfp[(i * nt + j) * 2 + 1];
280 if (c_simdBestPairAlignment > 2)
282 params->nbfp_aligned[(i * nt + j) * c_simdBestPairAlignment + 2] = 0;
283 params->nbfp_aligned[(i * nt + j) * c_simdBestPairAlignment + 3] = 0;
290 /* We use combination rule data for SIMD combination rule kernels
291 * and with LJ-PME kernels. We then only need parameters per atom type,
292 * not per pair of atom types.
294 params->nbfp_comb.resize(nt * 2);
295 switch (params->comb_rule)
298 params->comb_rule = ljcrGEOM;
300 for (int i = 0; i < nt; i++)
302 /* Store the sqrt of the diagonal from the nbfp matrix */
303 params->nbfp_comb[i * 2] = std::sqrt(params->nbfp[(i * nt + i) * 2]);
304 params->nbfp_comb[i * 2 + 1] = std::sqrt(params->nbfp[(i * nt + i) * 2 + 1]);
308 for (int i = 0; i < nt; i++)
310 /* Get 6*C6 and 12*C12 from the diagonal of the nbfp matrix */
311 const real c6 = params->nbfp[(i * nt + i) * 2];
312 const real c12 = params->nbfp[(i * nt + i) * 2 + 1];
313 if (c6 > 0 && c12 > 0)
315 /* We store 0.5*2^1/6*sigma and sqrt(4*3*eps),
316 * so we get 6*C6 and 12*C12 after combining.
318 params->nbfp_comb[i * 2] = 0.5 * gmx::sixthroot(c12 / c6);
319 params->nbfp_comb[i * 2 + 1] = std::sqrt(c6 * c6 / c12);
323 params->nbfp_comb[i * 2] = 0;
324 params->nbfp_comb[i * 2 + 1] = 0;
329 /* We always store the full matrix (see code above) */
331 default: gmx_incons("Unknown combination rule");
335 nbnxn_atomdata_t::SimdMasks::SimdMasks()
338 constexpr int simd_width = GMX_SIMD_REAL_WIDTH;
339 /* Set the diagonal cluster pair exclusion mask setup data.
340 * In the kernel we check 0 < j - i to generate the masks.
341 * Here we store j - i for generating the mask for the first i,
342 * we substract 0.5 to avoid rounding issues.
343 * In the kernel we can subtract 1 to generate the subsequent mask.
345 const int simd_4xn_diag_size = std::max(c_nbnxnCpuIClusterSize, simd_width);
346 diagonal_4xn_j_minus_i.resize(simd_4xn_diag_size);
347 for (int j = 0; j < simd_4xn_diag_size; j++)
349 diagonal_4xn_j_minus_i[j] = j - 0.5;
352 diagonal_2xnn_j_minus_i.resize(simd_width);
353 for (int j = 0; j < simd_width / 2; j++)
355 /* The j-cluster size is half the SIMD width */
356 diagonal_2xnn_j_minus_i[j] = j - 0.5;
357 /* The next half of the SIMD width is for i + 1 */
358 diagonal_2xnn_j_minus_i[simd_width / 2 + j] = j - 1 - 0.5;
361 /* We use up to 32 bits for exclusion masking.
362 * The same masks are used for the 4xN and 2x(N+N) kernels.
363 * The masks are read either into integer SIMD registers or into
364 * real SIMD registers (together with a cast).
365 * In single precision this means the real and integer SIMD registers
368 const int simd_excl_size = c_nbnxnCpuIClusterSize * simd_width;
369 # if GMX_DOUBLE && !GMX_SIMD_HAVE_INT32_LOGICAL
370 exclusion_filter64.resize(simd_excl_size);
372 exclusion_filter.resize(simd_excl_size);
375 for (int j = 0; j < simd_excl_size; j++)
377 /* Set the consecutive bits for masking pair exclusions */
378 # if GMX_DOUBLE && !GMX_SIMD_HAVE_INT32_LOGICAL
379 exclusion_filter64[j] = (1U << j);
381 exclusion_filter[j] = (1U << j);
385 if (!GMX_SIMD_HAVE_LOGICAL && !GMX_SIMD_HAVE_INT32_LOGICAL) // NOLINT(misc-redundant-expression)
387 // If the SIMD implementation has no bitwise logical operation support
388 // whatsoever we cannot use the normal masking. Instead,
389 // we generate a vector of all 2^4 possible ways an i atom
390 // interacts with its 4 j atoms. Each array entry contains
391 // GMX_SIMD_REAL_WIDTH values that are read with a single aligned SIMD load.
392 // Since there is no logical value representation in this case, we use
393 // any nonzero value to indicate 'true', while zero mean 'false'.
394 // This can then be converted to a SIMD boolean internally in the SIMD
395 // module by comparing to zero.
396 // Each array entry encodes how this i atom will interact with the 4 j atoms.
397 // Matching code exists in set_ci_top_excls() to generate indices into this array.
398 // Those indices are used in the kernels.
400 const int simd_excl_size = c_nbnxnCpuIClusterSize * c_nbnxnCpuIClusterSize;
401 const real simdFalse = 0.0;
402 const real simdTrue = 1.0;
404 interaction_array.resize(simd_excl_size * GMX_SIMD_REAL_WIDTH);
405 for (int j = 0; j < simd_excl_size; j++)
407 const int index = j * GMX_SIMD_REAL_WIDTH;
408 for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
410 interaction_array[index + i] = (j & (1 << i)) ? simdTrue : simdFalse;
417 nbnxn_atomdata_t::Params::Params(gmx::PinningPolicy pinningPolicy) :
419 nbfp({}, { pinningPolicy }),
420 nbfp_comb({}, { pinningPolicy }),
421 type({}, { pinningPolicy }),
422 lj_comb({}, { pinningPolicy }),
423 q({}, { pinningPolicy }),
426 energrp({}, { pinningPolicy })
430 nbnxn_atomdata_t::nbnxn_atomdata_t(gmx::PinningPolicy pinningPolicy) :
431 params_(pinningPolicy),
434 shift_vec({}, { pinningPolicy }),
435 x_({}, { pinningPolicy }),
437 bUseBufferFlags(FALSE),
438 bUseTreeReduce(FALSE)
442 /* Initializes an nbnxn_atomdata_t::Params data structure */
443 static void nbnxn_atomdata_params_init(const gmx::MDLogger& mdlog,
444 nbnxn_atomdata_t::Params* params,
445 const Nbnxm::KernelType kernelType,
446 int enbnxninitcombrule,
448 ArrayRef<const real> nbfp,
453 gmx_bool simple, bCombGeom, bCombLB, bSIMD;
457 fprintf(debug, "There are %d atom types in the system, adding one for nbnxn_atomdata_t\n", ntype);
459 params->numTypes = ntype + 1;
460 params->nbfp.resize(params->numTypes * params->numTypes * 2);
461 params->nbfp_comb.resize(params->numTypes * 2);
463 /* A tolerance of 1e-5 seems reasonable for (possibly hand-typed)
464 * force-field floating point parameters.
467 ptr = getenv("GMX_LJCOMB_TOL");
472 sscanf(ptr, "%lf", &dbl);
478 /* Temporarily fill params->nbfp_comb with sigma and epsilon
479 * to check for the LB rule.
481 for (int i = 0; i < ntype; i++)
483 c6 = nbfp[(i * ntype + i) * 2] / 6.0;
484 c12 = nbfp[(i * ntype + i) * 2 + 1] / 12.0;
485 if (c6 > 0 && c12 > 0)
487 params->nbfp_comb[i * 2] = gmx::sixthroot(c12 / c6);
488 params->nbfp_comb[i * 2 + 1] = 0.25 * c6 * c6 / c12;
490 else if (c6 == 0 && c12 == 0)
492 params->nbfp_comb[i * 2] = 0;
493 params->nbfp_comb[i * 2 + 1] = 0;
497 /* Can not use LB rule with only dispersion or repulsion */
502 for (int i = 0; i < params->numTypes; i++)
504 for (int j = 0; j < params->numTypes; j++)
506 if (i < ntype && j < ntype)
508 /* fr->nbfp has been updated, so that array too now stores c6/c12 including
509 * the 6.0/12.0 prefactors to save 2 flops in the most common case (force-only).
511 c6 = nbfp[(i * ntype + j) * 2];
512 c12 = nbfp[(i * ntype + j) * 2 + 1];
513 params->nbfp[(i * params->numTypes + j) * 2] = c6;
514 params->nbfp[(i * params->numTypes + j) * 2 + 1] = c12;
516 /* Compare 6*C6 and 12*C12 for geometric cobination rule */
519 && gmx_within_tol(c6 * c6, nbfp[(i * ntype + i) * 2] * nbfp[(j * ntype + j) * 2], tol)
520 && gmx_within_tol(c12 * c12,
521 nbfp[(i * ntype + i) * 2 + 1] * nbfp[(j * ntype + j) * 2 + 1],
524 /* Compare C6 and C12 for Lorentz-Berthelot combination rule */
529 && ((c6 == 0 && c12 == 0
530 && (params->nbfp_comb[i * 2 + 1] == 0 || params->nbfp_comb[j * 2 + 1] == 0))
531 || (c6 > 0 && c12 > 0
532 && gmx_within_tol(gmx::sixthroot(c12 / c6),
533 0.5 * (params->nbfp_comb[i * 2] + params->nbfp_comb[j * 2]),
535 && gmx_within_tol(0.25 * c6 * c6 / c12,
536 std::sqrt(params->nbfp_comb[i * 2 + 1]
537 * params->nbfp_comb[j * 2 + 1]),
542 /* Add zero parameters for the additional dummy atom type */
543 params->nbfp[(i * params->numTypes + j) * 2] = 0;
544 params->nbfp[(i * params->numTypes + j) * 2 + 1] = 0;
551 "Combination rules: geometric %s Lorentz-Berthelot %s\n",
552 gmx::boolToString(bCombGeom),
553 gmx::boolToString(bCombLB));
556 simple = Nbnxm::kernelTypeUsesSimplePairlist(kernelType);
558 switch (enbnxninitcombrule)
560 case enbnxninitcombruleDETECT:
561 /* We prefer the geometic combination rule,
562 * as that gives a slightly faster kernel than the LB rule.
566 params->comb_rule = ljcrGEOM;
570 params->comb_rule = ljcrLB;
574 params->comb_rule = ljcrNONE;
576 params->nbfp_comb.clear();
581 if (params->comb_rule == ljcrNONE)
583 mesg = "Using full Lennard-Jones parameter combination matrix";
587 mesg = gmx::formatString("Using %s Lennard-Jones combination rule",
588 enum_name(params->comb_rule, ljcrNR, c_ljcrNames));
590 GMX_LOG(mdlog.info).asParagraph().appendText(mesg);
593 case enbnxninitcombruleGEOM: params->comb_rule = ljcrGEOM; break;
594 case enbnxninitcombruleLB: params->comb_rule = ljcrLB; break;
595 case enbnxninitcombruleNONE:
596 params->comb_rule = ljcrNONE;
598 params->nbfp_comb.clear();
600 default: gmx_incons("Unknown enbnxninitcombrule");
603 bSIMD = Nbnxm::kernelTypeIsSimd(kernelType);
605 set_lj_parameter_data(params, bSIMD);
607 params->nenergrp = n_energygroups;
610 // We now check for energy groups already when starting mdrun
611 GMX_RELEASE_ASSERT(n_energygroups == 1, "GPU kernels do not support energy groups");
613 /* Temporary storage goes as #grp^3*simd_width^2/2, so limit to 64 */
614 if (params->nenergrp > 64)
616 gmx_fatal(FARGS, "With NxN kernels not more than 64 energy groups are supported\n");
618 params->neg_2log = 1;
619 while (params->nenergrp > (1 << params->neg_2log))
625 /* Initializes an nbnxn_atomdata_t data structure */
626 void nbnxn_atomdata_init(const gmx::MDLogger& mdlog,
627 nbnxn_atomdata_t* nbat,
628 const Nbnxm::KernelType kernelType,
629 int enbnxninitcombrule,
631 ArrayRef<const real> nbfp,
635 nbnxn_atomdata_params_init(
636 mdlog, &nbat->paramsDeprecated(), kernelType, enbnxninitcombrule, ntype, nbfp, n_energygroups);
638 const bool simple = Nbnxm::kernelTypeUsesSimplePairlist(kernelType);
639 const bool bSIMD = Nbnxm::kernelTypeIsSimd(kernelType);
647 pack_x = std::max(c_nbnxnCpuIClusterSize, Nbnxm::JClusterSizePerKernelType[kernelType]);
650 case 4: nbat->XFormat = nbatX4; break;
651 case 8: nbat->XFormat = nbatX8; break;
652 default: gmx_incons("Unsupported packing width");
657 nbat->XFormat = nbatXYZ;
660 nbat->FFormat = nbat->XFormat;
664 nbat->XFormat = nbatXYZQ;
665 nbat->FFormat = nbatXYZ;
668 nbat->shift_vec.resize(SHIFTS);
670 nbat->xstride = (nbat->XFormat == nbatXYZQ ? STRIDE_XYZQ : DIM);
671 nbat->fstride = (nbat->FFormat == nbatXYZQ ? STRIDE_XYZQ : DIM);
673 /* Initialize the output data structures */
674 for (int i = 0; i < nout; i++)
676 const auto& pinningPolicy = nbat->params().type.get_allocator().pinningPolicy();
677 nbat->out.emplace_back(
678 kernelType, nbat->params().nenergrp, 1 << nbat->params().neg_2log, pinningPolicy);
681 nbat->buffer_flags.clear();
683 const int nth = gmx_omp_nthreads_get(emntNonbonded);
685 const char* ptr = getenv("GMX_USE_TREEREDUCE");
688 nbat->bUseTreeReduce = (strtol(ptr, nullptr, 10) != 0);
692 nbat->bUseTreeReduce = false;
694 if (nbat->bUseTreeReduce)
696 GMX_LOG(mdlog.info).asParagraph().appendText("Using tree force reduction");
698 nbat->syncStep = new tMPI_Atomic[nth];
702 template<int packSize>
703 static void copy_lj_to_nbat_lj_comb(gmx::ArrayRef<const real> ljparam_type, const int* type, int na, real* ljparam_at)
705 /* The LJ params follow the combination rule:
706 * copy the params for the type array to the atom array.
708 for (int is = 0; is < na; is += packSize)
710 for (int k = 0; k < packSize; k++)
713 ljparam_at[is * 2 + k] = ljparam_type[type[i] * 2];
714 ljparam_at[is * 2 + packSize + k] = ljparam_type[type[i] * 2 + 1];
719 /* Sets the atom type in nbnxn_atomdata_t */
720 static void nbnxn_atomdata_set_atomtypes(nbnxn_atomdata_t::Params* params,
721 const Nbnxm::GridSet& gridSet,
722 ArrayRef<const int> atomTypes)
724 params->type.resize(gridSet.numGridAtomsTotal());
726 for (const Nbnxm::Grid& grid : gridSet.grids())
728 /* Loop over all columns and copy and fill */
729 for (int i = 0; i < grid.numColumns(); i++)
731 const int numAtoms = grid.paddedNumAtomsInColumn(i);
732 const int atomOffset = grid.firstAtomInColumn(i);
734 copy_int_to_nbat_int(gridSet.atomIndices().data() + atomOffset,
735 grid.numAtomsInColumn(i),
738 params->numTypes - 1,
739 params->type.data() + atomOffset);
744 /* Sets the LJ combination rule parameters in nbnxn_atomdata_t */
745 static void nbnxn_atomdata_set_ljcombparams(nbnxn_atomdata_t::Params* params,
747 const Nbnxm::GridSet& gridSet)
749 params->lj_comb.resize(gridSet.numGridAtomsTotal() * 2);
751 if (params->comb_rule != ljcrNONE)
753 for (const Nbnxm::Grid& grid : gridSet.grids())
755 /* Loop over all columns and copy and fill */
756 for (int i = 0; i < grid.numColumns(); i++)
758 const int numAtoms = grid.paddedNumAtomsInColumn(i);
759 const int atomOffset = grid.firstAtomInColumn(i);
761 if (XFormat == nbatX4)
763 copy_lj_to_nbat_lj_comb<c_packX4>(params->nbfp_comb,
764 params->type.data() + atomOffset,
766 params->lj_comb.data() + atomOffset * 2);
768 else if (XFormat == nbatX8)
770 copy_lj_to_nbat_lj_comb<c_packX8>(params->nbfp_comb,
771 params->type.data() + atomOffset,
773 params->lj_comb.data() + atomOffset * 2);
775 else if (XFormat == nbatXYZQ)
777 copy_lj_to_nbat_lj_comb<1>(params->nbfp_comb,
778 params->type.data() + atomOffset,
780 params->lj_comb.data() + atomOffset * 2);
787 /* Sets the charges in nbnxn_atomdata_t *nbat */
788 static void nbnxn_atomdata_set_charges(nbnxn_atomdata_t* nbat,
789 const Nbnxm::GridSet& gridSet,
790 ArrayRef<const real> charges)
792 if (nbat->XFormat != nbatXYZQ)
794 nbat->paramsDeprecated().q.resize(nbat->numAtoms());
797 for (const Nbnxm::Grid& grid : gridSet.grids())
799 /* Loop over all columns and copy and fill */
800 for (int cxy = 0; cxy < grid.numColumns(); cxy++)
802 const int atomOffset = grid.firstAtomInColumn(cxy);
803 const int numAtoms = grid.numAtomsInColumn(cxy);
804 const int paddedNumAtoms = grid.paddedNumAtomsInColumn(cxy);
806 if (nbat->XFormat == nbatXYZQ)
808 real* q = nbat->x().data() + atomOffset * STRIDE_XYZQ + ZZ + 1;
810 for (i = 0; i < numAtoms; i++)
812 *q = charges[gridSet.atomIndices()[atomOffset + i]];
815 /* Complete the partially filled last cell with zeros */
816 for (; i < paddedNumAtoms; i++)
824 real* q = nbat->paramsDeprecated().q.data() + atomOffset;
826 for (i = 0; i < numAtoms; i++)
828 *q = charges[gridSet.atomIndices()[atomOffset + i]];
831 /* Complete the partially filled last cell with zeros */
832 for (; i < paddedNumAtoms; i++)
842 /* Set the charges of perturbed atoms in nbnxn_atomdata_t to 0.
843 * This is to automatically remove the RF/PME self term in the nbnxn kernels.
844 * Part of the zero interactions are still calculated in the normal kernels.
845 * All perturbed interactions are calculated in the free energy kernel,
846 * using the original charge and LJ data, not nbnxn_atomdata_t.
848 static void nbnxn_atomdata_mask_fep(nbnxn_atomdata_t* nbat, const Nbnxm::GridSet& gridSet)
850 nbnxn_atomdata_t::Params& params = nbat->paramsDeprecated();
854 if (nbat->XFormat == nbatXYZQ)
856 q = nbat->x().data() + ZZ + 1;
857 stride_q = STRIDE_XYZQ;
865 for (const Nbnxm::Grid& grid : gridSet.grids())
868 if (grid.geometry().isSimple)
874 nsubc = c_gpuNumClusterPerCell;
877 int c_offset = grid.firstAtomInColumn(0);
879 /* Loop over all columns and copy and fill */
880 for (int c = 0; c < grid.numCells() * nsubc; c++)
882 /* Does this cluster contain perturbed particles? */
883 if (grid.clusterIsPerturbed(c))
885 const int numAtomsPerCluster = grid.geometry().numAtomsICluster;
886 for (int i = 0; i < numAtomsPerCluster; i++)
888 /* Is this a perturbed particle? */
889 if (grid.atomIsPerturbed(c, i))
891 int ind = c_offset + c * numAtomsPerCluster + i;
892 /* Set atom type and charge to non-interacting */
893 params.type[ind] = params.numTypes - 1;
894 q[ind * stride_q] = 0;
902 /* Copies the energy group indices to a reordered and packed array */
904 copy_egp_to_nbat_egps(const int* a, int na, int na_round, int na_c, int bit_shift, const int* in, int* innb)
910 for (i = 0; i < na; i += na_c)
912 /* Store na_c energy group numbers into one int */
914 for (int sa = 0; sa < na_c; sa++)
919 comb |= (GET_CGINFO_GID(in[at]) << (sa * bit_shift));
924 /* Complete the partially filled last cell with fill */
925 for (; i < na_round; i += na_c)
931 /* Set the energy group indices for atoms in nbnxn_atomdata_t */
932 static void nbnxn_atomdata_set_energygroups(nbnxn_atomdata_t::Params* params,
933 const Nbnxm::GridSet& gridSet,
934 ArrayRef<const int> atomInfo)
936 if (params->nenergrp == 1)
941 params->energrp.resize(gridSet.numGridAtomsTotal());
943 for (const Nbnxm::Grid& grid : gridSet.grids())
945 /* Loop over all columns and copy and fill */
946 for (int i = 0; i < grid.numColumns(); i++)
948 const int numAtoms = grid.paddedNumAtomsInColumn(i);
949 const int atomOffset = grid.firstAtomInColumn(i);
951 copy_egp_to_nbat_egps(gridSet.atomIndices().data() + atomOffset,
952 grid.numAtomsInColumn(i),
954 c_nbnxnCpuIClusterSize,
957 params->energrp.data() + grid.atomToCluster(atomOffset));
962 /* Sets all required atom parameter data in nbnxn_atomdata_t */
963 void nbnxn_atomdata_set(nbnxn_atomdata_t* nbat,
964 const Nbnxm::GridSet& gridSet,
965 ArrayRef<const int> atomTypes,
966 ArrayRef<const real> atomCharges,
967 ArrayRef<const int> atomInfo)
969 nbnxn_atomdata_t::Params& params = nbat->paramsDeprecated();
971 nbnxn_atomdata_set_atomtypes(¶ms, gridSet, atomTypes);
973 nbnxn_atomdata_set_charges(nbat, gridSet, atomCharges);
975 if (gridSet.haveFep())
977 nbnxn_atomdata_mask_fep(nbat, gridSet);
980 /* This must be done after masking types for FEP */
981 nbnxn_atomdata_set_ljcombparams(¶ms, nbat->XFormat, gridSet);
983 nbnxn_atomdata_set_energygroups(¶ms, gridSet, atomInfo);
986 /* Copies the shift vector array to nbnxn_atomdata_t */
987 void nbnxn_atomdata_copy_shiftvec(gmx_bool bDynamicBox, rvec* shift_vec, nbnxn_atomdata_t* nbat)
991 nbat->bDynamicBox = bDynamicBox;
992 for (i = 0; i < SHIFTS; i++)
994 copy_rvec(shift_vec[i], nbat->shift_vec[i]);
998 // This is slightly different from nbnxn_get_atom_range(...) at the end of the file
999 // TODO: Combine if possible
1000 static void getAtomRanges(const Nbnxm::GridSet& gridSet,
1001 const gmx::AtomLocality locality,
1007 case gmx::AtomLocality::All:
1009 *gridEnd = gridSet.grids().size();
1011 case gmx::AtomLocality::Local:
1015 case gmx::AtomLocality::NonLocal:
1017 *gridEnd = gridSet.grids().size();
1019 case gmx::AtomLocality::Count:
1020 GMX_ASSERT(false, "Count is invalid locality specifier");
1025 /* Copies (and reorders) the coordinates to nbnxn_atomdata_t */
1026 void nbnxn_atomdata_copy_x_to_nbat_x(const Nbnxm::GridSet& gridSet,
1027 const gmx::AtomLocality locality,
1029 const rvec* coordinates,
1030 nbnxn_atomdata_t* nbat)
1035 getAtomRanges(gridSet, locality, &gridBegin, &gridEnd);
1039 nbat->natoms_local = gridSet.grids()[0].atomIndexEnd();
1042 const int nth = gmx_omp_nthreads_get(emntPairsearch);
1043 #pragma omp parallel for num_threads(nth) schedule(static)
1044 for (int th = 0; th < nth; th++)
1048 for (int g = gridBegin; g < gridEnd; g++)
1050 const Nbnxm::Grid& grid = gridSet.grids()[g];
1051 const int numCellsXY = grid.numColumns();
1053 const int cxy0 = (numCellsXY * th + nth - 1) / nth;
1054 const int cxy1 = (numCellsXY * (th + 1) + nth - 1) / nth;
1056 for (int cxy = cxy0; cxy < cxy1; cxy++)
1058 const int na = grid.numAtomsInColumn(cxy);
1059 const int ash = grid.firstAtomInColumn(cxy);
1062 if (g == 0 && fillLocal)
1064 na_fill = grid.paddedNumAtomsInColumn(cxy);
1068 /* We fill only the real particle locations.
1069 * We assume the filling entries at the end have been
1070 * properly set before during pair-list generation.
1074 copy_rvec_to_nbat_real(gridSet.atomIndices().data() + ash,
1084 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1088 /* Copies (and reorders) the coordinates to nbnxn_atomdata_t on the GPU*/
1089 void nbnxn_atomdata_x_to_nbat_x_gpu(const Nbnxm::GridSet& gridSet,
1090 const gmx::AtomLocality locality,
1093 DeviceBuffer<RVec> d_x,
1094 GpuEventSynchronizer* xReadyOnDevice)
1099 getAtomRanges(gridSet, locality, &gridBegin, &gridEnd);
1101 for (int g = gridBegin; g < gridEnd; g++)
1103 nbnxn_gpu_x_to_nbat_x(gridSet.grids()[g],
1104 fillLocal && g == 0,
1110 gridSet.numColumnsMax());
1114 static void nbnxn_atomdata_clear_reals(gmx::ArrayRef<real> dest, int i0, int i1)
1116 for (int i = i0; i < i1; i++)
1122 gmx_unused static void nbnxn_atomdata_reduce_reals(real* gmx_restrict dest,
1124 const real** gmx_restrict src,
1131 /* The destination buffer contains data, add to it */
1132 for (int i = i0; i < i1; i++)
1134 for (int s = 0; s < nsrc; s++)
1136 dest[i] += src[s][i];
1142 /* The destination buffer is unitialized, set it first */
1143 for (int i = i0; i < i1; i++)
1145 dest[i] = src[0][i];
1146 for (int s = 1; s < nsrc; s++)
1148 dest[i] += src[s][i];
1154 gmx_unused static void nbnxn_atomdata_reduce_reals_simd(real gmx_unused* gmx_restrict dest,
1155 gmx_bool gmx_unused bDestSet,
1156 const gmx_unused real** gmx_restrict src,
1157 int gmx_unused nsrc,
1162 /* The SIMD width here is actually independent of that in the kernels,
1163 * but we use the same width for simplicity (usually optimal anyhow).
1165 SimdReal dest_SSE, src_SSE;
1169 for (int i = i0; i < i1; i += GMX_SIMD_REAL_WIDTH)
1171 dest_SSE = load<SimdReal>(dest + i);
1172 for (int s = 0; s < nsrc; s++)
1174 src_SSE = load<SimdReal>(src[s] + i);
1175 dest_SSE = dest_SSE + src_SSE;
1177 store(dest + i, dest_SSE);
1182 for (int i = i0; i < i1; i += GMX_SIMD_REAL_WIDTH)
1184 dest_SSE = load<SimdReal>(src[0] + i);
1185 for (int s = 1; s < nsrc; s++)
1187 src_SSE = load<SimdReal>(src[s] + i);
1188 dest_SSE = dest_SSE + src_SSE;
1190 store(dest + i, dest_SSE);
1196 /* Add part of the force array(s) from nbnxn_atomdata_t to f
1198 * Note: Adding restrict to f makes this function 50% slower with gcc 7.3
1200 static void nbnxn_atomdata_add_nbat_f_to_f_part(const Nbnxm::GridSet& gridSet,
1201 const nbnxn_atomdata_t& nbat,
1202 const nbnxn_atomdata_output_t& out,
1207 gmx::ArrayRef<const int> cell = gridSet.cells();
1208 // Note: Using ArrayRef instead makes this code 25% slower with gcc 7.3
1209 const real* fnb = out.f.data();
1211 /* Loop over all columns and copy and fill */
1212 switch (nbat.FFormat)
1216 for (int a = a0; a < a1; a++)
1218 int i = cell[a] * nbat.fstride;
1221 f[a][YY] += fnb[i + 1];
1222 f[a][ZZ] += fnb[i + 2];
1226 for (int a = a0; a < a1; a++)
1228 int i = atom_to_x_index<c_packX4>(cell[a]);
1230 f[a][XX] += fnb[i + XX * c_packX4];
1231 f[a][YY] += fnb[i + YY * c_packX4];
1232 f[a][ZZ] += fnb[i + ZZ * c_packX4];
1236 for (int a = a0; a < a1; a++)
1238 int i = atom_to_x_index<c_packX8>(cell[a]);
1240 f[a][XX] += fnb[i + XX * c_packX8];
1241 f[a][YY] += fnb[i + YY * c_packX8];
1242 f[a][ZZ] += fnb[i + ZZ * c_packX8];
1245 default: gmx_incons("Unsupported nbnxn_atomdata_t format");
1249 static inline unsigned char reverse_bits(unsigned char b)
1251 /* http://graphics.stanford.edu/~seander/bithacks.html#ReverseByteWith64BitsDiv */
1252 return (b * 0x0202020202ULL & 0x010884422010ULL) % 1023;
1255 static void nbnxn_atomdata_add_nbat_f_to_f_treereduce(nbnxn_atomdata_t* nbat, int nth)
1257 gmx::ArrayRef<const gmx_bitmask_t> flags = nbat->buffer_flags;
1259 int next_pow2 = 1 << (gmx::log2I(nth - 1) + 1);
1261 const int numOutputBuffers = nbat->out.size();
1262 GMX_ASSERT(numOutputBuffers == nth,
1263 "tree-reduce currently only works for numOutputBuffers==nth");
1265 memset(nbat->syncStep, 0, sizeof(*(nbat->syncStep)) * nth);
1267 #pragma omp parallel num_threads(nth)
1275 th = gmx_omp_get_thread_num();
1277 for (group_size = 2; group_size < 2 * next_pow2; group_size *= 2)
1279 int index[2], group_pos, partner_pos, wu;
1280 int partner_th = th ^ (group_size / 2);
1285 /* wait on partner thread - replaces full barrier */
1286 int sync_th, sync_group_size;
1288 # if defined(__clang__) && __clang_major__ >= 8
1289 // Suppress warnings that the use of memory_barrier may be excessive
1290 // Only exists beginning with clang-8
1291 # pragma clang diagnostic push
1292 # pragma clang diagnostic ignored "-Watomic-implicit-seq-cst"
1295 tMPI_Atomic_memory_barrier(); /* guarantee data is saved before marking work as done */
1296 tMPI_Atomic_set(&(nbat->syncStep[th]), group_size / 2); /* mark previous step as completed */
1298 /* find thread to sync with. Equal to partner_th unless nth is not a power of two. */
1299 for (sync_th = partner_th, sync_group_size = group_size;
1300 sync_th >= nth && sync_group_size > 2;
1301 sync_group_size /= 2)
1303 sync_th &= ~(sync_group_size / 4);
1305 if (sync_th < nth) /* otherwise nothing to sync index[1] will be >=nout */
1307 /* wait on the thread which computed input data in previous step */
1308 while (tMPI_Atomic_get(static_cast<volatile tMPI_Atomic_t*>(&(nbat->syncStep[sync_th])))
1313 /* guarantee that no later load happens before wait loop is finisehd */
1314 tMPI_Atomic_memory_barrier();
1316 # if defined(__clang__) && __clang_major__ >= 8
1317 # pragma clang diagnostic pop
1319 #else /* TMPI_ATOMICS */
1320 # pragma omp barrier
1324 /* Calculate buffers to sum (result goes into first buffer) */
1325 group_pos = th % group_size;
1326 index[0] = th - group_pos;
1327 index[1] = index[0] + group_size / 2;
1329 /* If no second buffer, nothing to do */
1330 if (index[1] >= numOutputBuffers && group_size > 2)
1335 #if NBNXN_BUFFERFLAG_MAX_THREADS > 256
1336 # error reverse_bits assumes max 256 threads
1338 /* Position is permuted so that one of the 2 vectors being added was computed on the same thread in the previous step.
1339 This improves locality and enables to sync with just a single thread between steps (=the levels in the btree).
1340 The permutation which allows this corresponds to reversing the bits of the group position.
1342 group_pos = reverse_bits(group_pos) / (256 / group_size);
1344 partner_pos = group_pos ^ 1;
1346 /* loop over two work-units (own and partner) */
1347 for (wu = 0; wu < 2; wu++)
1351 if (partner_th < nth)
1353 break; /* partner exists we don't have to do his work */
1357 group_pos = partner_pos;
1361 /* Calculate the cell-block range for our thread */
1362 b0 = (flags.size() * group_pos) / group_size;
1363 b1 = (flags.size() * (group_pos + 1)) / group_size;
1365 for (b = b0; b < b1; b++)
1367 i0 = b * NBNXN_BUFFERFLAG_SIZE * nbat->fstride;
1368 i1 = (b + 1) * NBNXN_BUFFERFLAG_SIZE * nbat->fstride;
1370 if (bitmask_is_set(flags[b], index[1]) || group_size > 2)
1372 const real* fIndex1 = nbat->out[index[1]].f.data();
1374 nbnxn_atomdata_reduce_reals_simd
1376 nbnxn_atomdata_reduce_reals
1378 (nbat->out[index[0]].f.data(),
1379 bitmask_is_set(flags[b], index[0]) || group_size > 2,
1385 else if (!bitmask_is_set(flags[b], index[0]))
1387 nbnxn_atomdata_clear_reals(nbat->out[index[0]].f, i0, i1);
1393 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1398 static void nbnxn_atomdata_add_nbat_f_to_f_stdreduce(nbnxn_atomdata_t* nbat, int nth)
1400 #pragma omp parallel for num_threads(nth) schedule(static)
1401 for (int th = 0; th < nth; th++)
1406 const real* fptr[NBNXN_BUFFERFLAG_MAX_THREADS];
1408 gmx::ArrayRef<const gmx_bitmask_t> flags = nbat->buffer_flags;
1410 /* Calculate the cell-block range for our thread */
1411 int b0 = (flags.size() * th) / nth;
1412 int b1 = (flags.size() * (th + 1)) / nth;
1414 for (int b = b0; b < b1; b++)
1416 int i0 = b * NBNXN_BUFFERFLAG_SIZE * nbat->fstride;
1417 int i1 = (b + 1) * NBNXN_BUFFERFLAG_SIZE * nbat->fstride;
1420 for (gmx::index out = 1; out < gmx::ssize(nbat->out); out++)
1422 if (bitmask_is_set(flags[b], out))
1424 fptr[nfptr++] = nbat->out[out].f.data();
1430 nbnxn_atomdata_reduce_reals_simd
1432 nbnxn_atomdata_reduce_reals
1434 (nbat->out[0].f.data(), bitmask_is_set(flags[b], 0), fptr, nfptr, i0, i1);
1436 else if (!bitmask_is_set(flags[b], 0))
1438 nbnxn_atomdata_clear_reals(nbat->out[0].f, i0, i1);
1442 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1447 /* Add the force array(s) from nbnxn_atomdata_t to f */
1448 void reduceForces(nbnxn_atomdata_t* nbat, const gmx::AtomLocality locality, const Nbnxm::GridSet& gridSet, rvec* f)
1453 nbnxn_get_atom_range(locality, gridSet, &a0, &na);
1457 /* The are no atoms for this reduction, avoid some overhead */
1461 int nth = gmx_omp_nthreads_get(emntNonbonded);
1463 if (nbat->out.size() > 1)
1465 if (locality != gmx::AtomLocality::All)
1467 gmx_incons("add_f_to_f called with nout>1 and locality!=eatAll");
1470 /* Reduce the force thread output buffers into buffer 0, before adding
1471 * them to the, differently ordered, "real" force buffer.
1473 if (nbat->bUseTreeReduce)
1475 nbnxn_atomdata_add_nbat_f_to_f_treereduce(nbat, nth);
1479 nbnxn_atomdata_add_nbat_f_to_f_stdreduce(nbat, nth);
1482 #pragma omp parallel for num_threads(nth) schedule(static)
1483 for (int th = 0; th < nth; th++)
1487 nbnxn_atomdata_add_nbat_f_to_f_part(
1488 gridSet, *nbat, nbat->out[0], a0 + ((th + 0) * na) / nth, a0 + ((th + 1) * na) / nth, f);
1490 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1494 void nbnxn_atomdata_add_nbat_fshift_to_fshift(const nbnxn_atomdata_t& nbat, gmx::ArrayRef<gmx::RVec> fshift)
1496 gmx::ArrayRef<const nbnxn_atomdata_output_t> outputBuffers = nbat.out;
1498 for (int s = 0; s < SHIFTS; s++)
1502 for (const nbnxn_atomdata_output_t& out : outputBuffers)
1504 sum[XX] += out.fshift[s * DIM + XX];
1505 sum[YY] += out.fshift[s * DIM + YY];
1506 sum[ZZ] += out.fshift[s * DIM + ZZ];
1512 void nbnxn_get_atom_range(const gmx::AtomLocality atomLocality,
1513 const Nbnxm::GridSet& gridSet,
1518 switch (atomLocality)
1520 case gmx::AtomLocality::All:
1522 *nAtoms = gridSet.numRealAtomsTotal();
1524 case gmx::AtomLocality::Local:
1526 *nAtoms = gridSet.numRealAtomsLocal();
1528 case gmx::AtomLocality::NonLocal:
1529 *atomStart = gridSet.numRealAtomsLocal();
1530 *nAtoms = gridSet.numRealAtomsTotal() - gridSet.numRealAtomsLocal();
1532 case gmx::AtomLocality::Count:
1533 GMX_ASSERT(false, "Count is invalid locality specifier");