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
76 const char* enumValueToString(LJCombinationRule enumValue)
78 static constexpr gmx::EnumerationArray<LJCombinationRule, const char*> s_ljCombinationRuleNames = {
79 "Geometric", "Lorentz-Berthelot", "None"
81 return s_ljCombinationRuleNames[enumValue];
84 void nbnxn_atomdata_t::resizeCoordinateBuffer(int numAtoms)
88 x_.resize(numAtoms * xstride);
91 void nbnxn_atomdata_t::resizeForceBuffers()
93 /* Force buffers need padding up to a multiple of the buffer flag size */
94 const int paddedSize =
95 (numAtoms() + NBNXN_BUFFERFLAG_SIZE - 1) / NBNXN_BUFFERFLAG_SIZE * NBNXN_BUFFERFLAG_SIZE;
97 /* Should we let each thread allocate it's own data instead? */
98 for (nbnxn_atomdata_output_t& outBuffer : out)
100 outBuffer.f.resize(paddedSize * fstride);
104 /* Initializes an nbnxn_atomdata_output_t data structure */
105 nbnxn_atomdata_output_t::nbnxn_atomdata_output_t(Nbnxm::KernelType kernelType,
107 int simdEnergyBufferStride,
108 gmx::PinningPolicy pinningPolicy) :
109 f({}, { pinningPolicy }),
110 fshift({}, { pinningPolicy }),
111 Vvdw({}, { pinningPolicy }),
112 Vc({}, { pinningPolicy })
114 fshift.resize(SHIFTS * DIM);
115 Vvdw.resize(numEnergyGroups * numEnergyGroups);
116 Vc.resize(numEnergyGroups * numEnergyGroups);
118 if (Nbnxm::kernelTypeIsSimd(kernelType))
120 int cj_size = Nbnxm::JClusterSizePerKernelType[kernelType];
122 numEnergyGroups * numEnergyGroups * simdEnergyBufferStride * (cj_size / 2) * cj_size;
123 VSvdw.resize(numElements);
124 VSc.resize(numElements);
128 static void copy_int_to_nbat_int(const int* a, int na, int na_round, const int* in, int fill, int* innb)
133 for (i = 0; i < na; i++)
135 innb[j++] = in[a[i]];
137 /* Complete the partially filled last cell with fill */
138 for (; i < na_round; i++)
144 void copy_rvec_to_nbat_real(const int* a, int na, int na_round, const rvec* x, int nbatFormat, real* xnb, int a0)
146 /* We complete partially filled cells, can only be the last one in each
147 * column, with coordinates farAway. The actual coordinate value does
148 * not influence the results, since these filler particles do not interact.
149 * Clusters with normal atoms + fillers have a bounding box based only
150 * on the coordinates of the atoms. Clusters with only fillers have as
151 * the bounding box the coordinates of the first filler. Such clusters
152 * are not considered as i-entries, but they are considered as j-entries.
153 * So for performance it is better to have their bounding boxes far away,
154 * such that filler only clusters don't end up in the pair list.
156 const real farAway = -1000000;
164 for (i = 0; i < na; i++)
166 xnb[j++] = x[a[i]][XX];
167 xnb[j++] = x[a[i]][YY];
168 xnb[j++] = x[a[i]][ZZ];
170 /* Complete the partially filled last cell with farAway elements */
171 for (; i < na_round; i++)
179 j = a0 * STRIDE_XYZQ;
180 for (i = 0; i < na; i++)
182 xnb[j++] = x[a[i]][XX];
183 xnb[j++] = x[a[i]][YY];
184 xnb[j++] = x[a[i]][ZZ];
187 /* Complete the partially filled last cell with zeros */
188 for (; i < na_round; i++)
197 j = atom_to_x_index<c_packX4>(a0);
198 c = a0 & (c_packX4 - 1);
199 for (i = 0; i < na; i++)
201 xnb[j + XX * c_packX4] = x[a[i]][XX];
202 xnb[j + YY * c_packX4] = x[a[i]][YY];
203 xnb[j + ZZ * c_packX4] = x[a[i]][ZZ];
208 j += (DIM - 1) * c_packX4;
212 /* Complete the partially filled last cell with zeros */
213 for (; i < na_round; i++)
215 xnb[j + XX * c_packX4] = farAway;
216 xnb[j + YY * c_packX4] = farAway;
217 xnb[j + ZZ * c_packX4] = farAway;
222 j += (DIM - 1) * c_packX4;
228 j = atom_to_x_index<c_packX8>(a0);
229 c = a0 & (c_packX8 - 1);
230 for (i = 0; i < na; i++)
232 xnb[j + XX * c_packX8] = x[a[i]][XX];
233 xnb[j + YY * c_packX8] = x[a[i]][YY];
234 xnb[j + ZZ * c_packX8] = x[a[i]][ZZ];
239 j += (DIM - 1) * c_packX8;
243 /* Complete the partially filled last cell with zeros */
244 for (; i < na_round; i++)
246 xnb[j + XX * c_packX8] = farAway;
247 xnb[j + YY * c_packX8] = farAway;
248 xnb[j + ZZ * c_packX8] = farAway;
253 j += (DIM - 1) * c_packX8;
258 default: gmx_incons("Unsupported nbnxn_atomdata_t format");
262 /* Stores the LJ parameter data in a format convenient for different kernels */
263 static void set_lj_parameter_data(nbnxn_atomdata_t::Params* params, gmx_bool bSIMD)
265 int nt = params->numTypes;
270 /* nbfp_aligned stores two parameters using the stride most suitable
271 * for the present SIMD architecture, as specified by the constant
272 * c_simdBestPairAlignment from the SIMD header.
273 * There's a slight inefficiency in allocating and initializing nbfp_aligned
274 * when it might not be used, but introducing the conditional code is not
277 params->nbfp_aligned.resize(nt * nt * c_simdBestPairAlignment);
279 for (int i = 0; i < nt; i++)
281 for (int j = 0; j < nt; j++)
283 params->nbfp_aligned[(i * nt + j) * c_simdBestPairAlignment + 0] =
284 params->nbfp[(i * nt + j) * 2 + 0];
285 params->nbfp_aligned[(i * nt + j) * c_simdBestPairAlignment + 1] =
286 params->nbfp[(i * nt + j) * 2 + 1];
287 if (c_simdBestPairAlignment > 2)
289 params->nbfp_aligned[(i * nt + j) * c_simdBestPairAlignment + 2] = 0;
290 params->nbfp_aligned[(i * nt + j) * c_simdBestPairAlignment + 3] = 0;
297 /* We use combination rule data for SIMD combination rule kernels
298 * and with LJ-PME kernels. We then only need parameters per atom type,
299 * not per pair of atom types.
301 params->nbfp_comb.resize(nt * 2);
302 switch (params->ljCombinationRule)
304 case LJCombinationRule::Geometric:
305 for (int i = 0; i < nt; i++)
307 /* Store the sqrt of the diagonal from the nbfp matrix */
308 params->nbfp_comb[i * 2] = std::sqrt(params->nbfp[(i * nt + i) * 2]);
309 params->nbfp_comb[i * 2 + 1] = std::sqrt(params->nbfp[(i * nt + i) * 2 + 1]);
312 case LJCombinationRule::LorentzBerthelot:
313 for (int i = 0; i < nt; i++)
315 /* Get 6*C6 and 12*C12 from the diagonal of the nbfp matrix */
316 const real c6 = params->nbfp[(i * nt + i) * 2];
317 const real c12 = params->nbfp[(i * nt + i) * 2 + 1];
318 if (c6 > 0 && c12 > 0)
320 /* We store 0.5*2^1/6*sigma and sqrt(4*3*eps),
321 * so we get 6*C6 and 12*C12 after combining.
323 params->nbfp_comb[i * 2] = 0.5 * gmx::sixthroot(c12 / c6);
324 params->nbfp_comb[i * 2 + 1] = std::sqrt(c6 * c6 / c12);
328 params->nbfp_comb[i * 2] = 0;
329 params->nbfp_comb[i * 2 + 1] = 0;
333 case LJCombinationRule::None:
334 /* We always store the full matrix (see code above) */
336 default: gmx_incons("Unknown combination rule");
340 nbnxn_atomdata_t::SimdMasks::SimdMasks()
343 constexpr int simd_width = GMX_SIMD_REAL_WIDTH;
344 /* Set the diagonal cluster pair exclusion mask setup data.
345 * In the kernel we check 0 < j - i to generate the masks.
346 * Here we store j - i for generating the mask for the first i,
347 * we substract 0.5 to avoid rounding issues.
348 * In the kernel we can subtract 1 to generate the subsequent mask.
350 const int simd_4xn_diag_size = std::max(c_nbnxnCpuIClusterSize, simd_width);
351 diagonal_4xn_j_minus_i.resize(simd_4xn_diag_size);
352 for (int j = 0; j < simd_4xn_diag_size; j++)
354 diagonal_4xn_j_minus_i[j] = j - 0.5;
357 diagonal_2xnn_j_minus_i.resize(simd_width);
358 for (int j = 0; j < simd_width / 2; j++)
360 /* The j-cluster size is half the SIMD width */
361 diagonal_2xnn_j_minus_i[j] = j - 0.5;
362 /* The next half of the SIMD width is for i + 1 */
363 diagonal_2xnn_j_minus_i[simd_width / 2 + j] = j - 1 - 0.5;
366 /* We use up to 32 bits for exclusion masking.
367 * The same masks are used for the 4xN and 2x(N+N) kernels.
368 * The masks are read either into integer SIMD registers or into
369 * real SIMD registers (together with a cast).
370 * In single precision this means the real and integer SIMD registers
373 const int simd_excl_size = c_nbnxnCpuIClusterSize * simd_width;
374 # if GMX_DOUBLE && !GMX_SIMD_HAVE_INT32_LOGICAL
375 exclusion_filter64.resize(simd_excl_size);
377 exclusion_filter.resize(simd_excl_size);
380 for (int j = 0; j < simd_excl_size; j++)
382 /* Set the consecutive bits for masking pair exclusions */
383 # if GMX_DOUBLE && !GMX_SIMD_HAVE_INT32_LOGICAL
384 exclusion_filter64[j] = (1U << j);
386 exclusion_filter[j] = (1U << j);
390 if (!GMX_SIMD_HAVE_LOGICAL && !GMX_SIMD_HAVE_INT32_LOGICAL) // NOLINT(misc-redundant-expression)
392 // If the SIMD implementation has no bitwise logical operation support
393 // whatsoever we cannot use the normal masking. Instead,
394 // we generate a vector of all 2^4 possible ways an i atom
395 // interacts with its 4 j atoms. Each array entry contains
396 // GMX_SIMD_REAL_WIDTH values that are read with a single aligned SIMD load.
397 // Since there is no logical value representation in this case, we use
398 // any nonzero value to indicate 'true', while zero mean 'false'.
399 // This can then be converted to a SIMD boolean internally in the SIMD
400 // module by comparing to zero.
401 // Each array entry encodes how this i atom will interact with the 4 j atoms.
402 // Matching code exists in set_ci_top_excls() to generate indices into this array.
403 // Those indices are used in the kernels.
405 const int simd_excl_size = c_nbnxnCpuIClusterSize * c_nbnxnCpuIClusterSize;
406 const real simdFalse = 0.0;
407 const real simdTrue = 1.0;
409 interaction_array.resize(simd_excl_size * GMX_SIMD_REAL_WIDTH);
410 for (int j = 0; j < simd_excl_size; j++)
412 const int index = j * GMX_SIMD_REAL_WIDTH;
413 for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
415 interaction_array[index + i] = (j & (1 << i)) ? simdTrue : simdFalse;
422 nbnxn_atomdata_t::Params::Params(gmx::PinningPolicy pinningPolicy) :
424 nbfp({}, { pinningPolicy }),
425 nbfp_comb({}, { pinningPolicy }),
426 type({}, { pinningPolicy }),
427 lj_comb({}, { pinningPolicy }),
428 q({}, { pinningPolicy }),
431 energrp({}, { pinningPolicy })
435 nbnxn_atomdata_t::nbnxn_atomdata_t(gmx::PinningPolicy pinningPolicy) :
436 params_(pinningPolicy),
439 shift_vec({}, { pinningPolicy }),
440 x_({}, { pinningPolicy }),
442 bUseBufferFlags(FALSE),
443 bUseTreeReduce(FALSE)
447 /* Initializes an nbnxn_atomdata_t::Params data structure */
448 static void nbnxn_atomdata_params_init(const gmx::MDLogger& mdlog,
449 nbnxn_atomdata_t::Params* params,
450 const Nbnxm::KernelType kernelType,
451 int enbnxninitcombrule,
453 ArrayRef<const real> nbfp,
458 gmx_bool simple, bCombGeom, bCombLB, bSIMD;
462 fprintf(debug, "There are %d atom types in the system, adding one for nbnxn_atomdata_t\n", ntype);
464 params->numTypes = ntype + 1;
465 params->nbfp.resize(params->numTypes * params->numTypes * 2);
466 params->nbfp_comb.resize(params->numTypes * 2);
468 /* A tolerance of 1e-5 seems reasonable for (possibly hand-typed)
469 * force-field floating point parameters.
472 ptr = getenv("GMX_LJCOMB_TOL");
477 sscanf(ptr, "%lf", &dbl);
483 /* Temporarily fill params->nbfp_comb with sigma and epsilon
484 * to check for the LB rule.
486 for (int i = 0; i < ntype; i++)
488 c6 = nbfp[(i * ntype + i) * 2] / 6.0;
489 c12 = nbfp[(i * ntype + i) * 2 + 1] / 12.0;
490 if (c6 > 0 && c12 > 0)
492 params->nbfp_comb[i * 2] = gmx::sixthroot(c12 / c6);
493 params->nbfp_comb[i * 2 + 1] = 0.25 * c6 * c6 / c12;
495 else if (c6 == 0 && c12 == 0)
497 params->nbfp_comb[i * 2] = 0;
498 params->nbfp_comb[i * 2 + 1] = 0;
502 /* Can not use LB rule with only dispersion or repulsion */
507 for (int i = 0; i < params->numTypes; i++)
509 for (int j = 0; j < params->numTypes; j++)
511 if (i < ntype && j < ntype)
513 /* fr->nbfp has been updated, so that array too now stores c6/c12 including
514 * the 6.0/12.0 prefactors to save 2 flops in the most common case (force-only).
516 c6 = nbfp[(i * ntype + j) * 2];
517 c12 = nbfp[(i * ntype + j) * 2 + 1];
518 params->nbfp[(i * params->numTypes + j) * 2] = c6;
519 params->nbfp[(i * params->numTypes + j) * 2 + 1] = c12;
521 /* Compare 6*C6 and 12*C12 for geometric cobination rule */
524 && gmx_within_tol(c6 * c6, nbfp[(i * ntype + i) * 2] * nbfp[(j * ntype + j) * 2], tol)
525 && gmx_within_tol(c12 * c12,
526 nbfp[(i * ntype + i) * 2 + 1] * nbfp[(j * ntype + j) * 2 + 1],
529 /* Compare C6 and C12 for Lorentz-Berthelot combination rule */
534 && ((c6 == 0 && c12 == 0
535 && (params->nbfp_comb[i * 2 + 1] == 0 || params->nbfp_comb[j * 2 + 1] == 0))
536 || (c6 > 0 && c12 > 0
537 && gmx_within_tol(gmx::sixthroot(c12 / c6),
538 0.5 * (params->nbfp_comb[i * 2] + params->nbfp_comb[j * 2]),
540 && gmx_within_tol(0.25 * c6 * c6 / c12,
541 std::sqrt(params->nbfp_comb[i * 2 + 1]
542 * params->nbfp_comb[j * 2 + 1]),
547 /* Add zero parameters for the additional dummy atom type */
548 params->nbfp[(i * params->numTypes + j) * 2] = 0;
549 params->nbfp[(i * params->numTypes + j) * 2 + 1] = 0;
556 "Combination rules: geometric %s Lorentz-Berthelot %s\n",
557 gmx::boolToString(bCombGeom),
558 gmx::boolToString(bCombLB));
561 simple = Nbnxm::kernelTypeUsesSimplePairlist(kernelType);
563 switch (enbnxninitcombrule)
565 case enbnxninitcombruleDETECT:
566 /* We prefer the geometic combination rule,
567 * as that gives a slightly faster kernel than the LB rule.
571 params->ljCombinationRule = LJCombinationRule::Geometric;
575 params->ljCombinationRule = LJCombinationRule::LorentzBerthelot;
579 params->ljCombinationRule = LJCombinationRule::None;
581 params->nbfp_comb.clear();
586 if (params->ljCombinationRule == LJCombinationRule::None)
588 mesg = "Using full Lennard-Jones parameter combination matrix";
592 mesg = gmx::formatString("Using %s Lennard-Jones combination rule",
593 enumValueToString(params->ljCombinationRule));
595 GMX_LOG(mdlog.info).asParagraph().appendText(mesg);
598 case enbnxninitcombruleGEOM:
599 params->ljCombinationRule = LJCombinationRule::Geometric;
601 case enbnxninitcombruleLB:
602 params->ljCombinationRule = LJCombinationRule::LorentzBerthelot;
604 case enbnxninitcombruleNONE:
605 params->ljCombinationRule = LJCombinationRule::None;
607 params->nbfp_comb.clear();
609 default: gmx_incons("Unknown enbnxninitcombrule");
612 bSIMD = Nbnxm::kernelTypeIsSimd(kernelType);
614 set_lj_parameter_data(params, bSIMD);
616 params->nenergrp = n_energygroups;
619 // We now check for energy groups already when starting mdrun
620 GMX_RELEASE_ASSERT(n_energygroups == 1, "GPU kernels do not support energy groups");
622 /* Temporary storage goes as #grp^3*simd_width^2/2, so limit to 64 */
623 if (params->nenergrp > 64)
625 gmx_fatal(FARGS, "With NxN kernels not more than 64 energy groups are supported\n");
627 params->neg_2log = 1;
628 while (params->nenergrp > (1 << params->neg_2log))
634 /* Initializes an nbnxn_atomdata_t data structure */
635 void nbnxn_atomdata_init(const gmx::MDLogger& mdlog,
636 nbnxn_atomdata_t* nbat,
637 const Nbnxm::KernelType kernelType,
638 int enbnxninitcombrule,
640 ArrayRef<const real> nbfp,
644 nbnxn_atomdata_params_init(
645 mdlog, &nbat->paramsDeprecated(), kernelType, enbnxninitcombrule, ntype, nbfp, n_energygroups);
647 const bool simple = Nbnxm::kernelTypeUsesSimplePairlist(kernelType);
648 const bool bSIMD = Nbnxm::kernelTypeIsSimd(kernelType);
656 pack_x = std::max(c_nbnxnCpuIClusterSize, Nbnxm::JClusterSizePerKernelType[kernelType]);
659 case 4: nbat->XFormat = nbatX4; break;
660 case 8: nbat->XFormat = nbatX8; break;
661 default: gmx_incons("Unsupported packing width");
666 nbat->XFormat = nbatXYZ;
669 nbat->FFormat = nbat->XFormat;
673 nbat->XFormat = nbatXYZQ;
674 nbat->FFormat = nbatXYZ;
677 nbat->shift_vec.resize(SHIFTS);
679 nbat->xstride = (nbat->XFormat == nbatXYZQ ? STRIDE_XYZQ : DIM);
680 nbat->fstride = (nbat->FFormat == nbatXYZQ ? STRIDE_XYZQ : DIM);
682 /* Initialize the output data structures */
683 for (int i = 0; i < nout; i++)
685 const auto& pinningPolicy = nbat->params().type.get_allocator().pinningPolicy();
686 nbat->out.emplace_back(
687 kernelType, nbat->params().nenergrp, 1 << nbat->params().neg_2log, pinningPolicy);
690 nbat->buffer_flags.clear();
692 const int nth = gmx_omp_nthreads_get(emntNonbonded);
694 const char* ptr = getenv("GMX_USE_TREEREDUCE");
697 nbat->bUseTreeReduce = (strtol(ptr, nullptr, 10) != 0);
701 nbat->bUseTreeReduce = false;
703 if (nbat->bUseTreeReduce)
705 GMX_LOG(mdlog.info).asParagraph().appendText("Using tree force reduction");
707 nbat->syncStep = new tMPI_Atomic[nth];
711 template<int packSize>
712 static void copy_lj_to_nbat_lj_comb(gmx::ArrayRef<const real> ljparam_type, const int* type, int na, real* ljparam_at)
714 /* The LJ params follow the combination rule:
715 * copy the params for the type array to the atom array.
717 for (int is = 0; is < na; is += packSize)
719 for (int k = 0; k < packSize; k++)
722 ljparam_at[is * 2 + k] = ljparam_type[type[i] * 2];
723 ljparam_at[is * 2 + packSize + k] = ljparam_type[type[i] * 2 + 1];
728 /* Sets the atom type in nbnxn_atomdata_t */
729 static void nbnxn_atomdata_set_atomtypes(nbnxn_atomdata_t::Params* params,
730 const Nbnxm::GridSet& gridSet,
731 ArrayRef<const int> atomTypes)
733 params->type.resize(gridSet.numGridAtomsTotal());
735 for (const Nbnxm::Grid& grid : gridSet.grids())
737 /* Loop over all columns and copy and fill */
738 for (int i = 0; i < grid.numColumns(); i++)
740 const int numAtoms = grid.paddedNumAtomsInColumn(i);
741 const int atomOffset = grid.firstAtomInColumn(i);
743 copy_int_to_nbat_int(gridSet.atomIndices().data() + atomOffset,
744 grid.numAtomsInColumn(i),
747 params->numTypes - 1,
748 params->type.data() + atomOffset);
753 /* Sets the LJ combination rule parameters in nbnxn_atomdata_t */
754 static void nbnxn_atomdata_set_ljcombparams(nbnxn_atomdata_t::Params* params,
756 const Nbnxm::GridSet& gridSet)
758 params->lj_comb.resize(gridSet.numGridAtomsTotal() * 2);
760 if (params->ljCombinationRule != LJCombinationRule::None)
762 for (const Nbnxm::Grid& grid : gridSet.grids())
764 /* Loop over all columns and copy and fill */
765 for (int i = 0; i < grid.numColumns(); i++)
767 const int numAtoms = grid.paddedNumAtomsInColumn(i);
768 const int atomOffset = grid.firstAtomInColumn(i);
770 if (XFormat == nbatX4)
772 copy_lj_to_nbat_lj_comb<c_packX4>(params->nbfp_comb,
773 params->type.data() + atomOffset,
775 params->lj_comb.data() + atomOffset * 2);
777 else if (XFormat == nbatX8)
779 copy_lj_to_nbat_lj_comb<c_packX8>(params->nbfp_comb,
780 params->type.data() + atomOffset,
782 params->lj_comb.data() + atomOffset * 2);
784 else if (XFormat == nbatXYZQ)
786 copy_lj_to_nbat_lj_comb<1>(params->nbfp_comb,
787 params->type.data() + atomOffset,
789 params->lj_comb.data() + atomOffset * 2);
796 /* Sets the charges in nbnxn_atomdata_t *nbat */
797 static void nbnxn_atomdata_set_charges(nbnxn_atomdata_t* nbat,
798 const Nbnxm::GridSet& gridSet,
799 ArrayRef<const real> charges)
801 if (nbat->XFormat != nbatXYZQ)
803 nbat->paramsDeprecated().q.resize(nbat->numAtoms());
806 for (const Nbnxm::Grid& grid : gridSet.grids())
808 /* Loop over all columns and copy and fill */
809 for (int cxy = 0; cxy < grid.numColumns(); cxy++)
811 const int atomOffset = grid.firstAtomInColumn(cxy);
812 const int numAtoms = grid.numAtomsInColumn(cxy);
813 const int paddedNumAtoms = grid.paddedNumAtomsInColumn(cxy);
815 if (nbat->XFormat == nbatXYZQ)
817 real* q = nbat->x().data() + atomOffset * STRIDE_XYZQ + ZZ + 1;
819 for (i = 0; i < numAtoms; i++)
821 *q = charges[gridSet.atomIndices()[atomOffset + i]];
824 /* Complete the partially filled last cell with zeros */
825 for (; i < paddedNumAtoms; i++)
833 real* q = nbat->paramsDeprecated().q.data() + atomOffset;
835 for (i = 0; i < numAtoms; i++)
837 *q = charges[gridSet.atomIndices()[atomOffset + i]];
840 /* Complete the partially filled last cell with zeros */
841 for (; i < paddedNumAtoms; i++)
851 /* Set the charges of perturbed atoms in nbnxn_atomdata_t to 0.
852 * This is to automatically remove the RF/PME self term in the nbnxn kernels.
853 * Part of the zero interactions are still calculated in the normal kernels.
854 * All perturbed interactions are calculated in the free energy kernel,
855 * using the original charge and LJ data, not nbnxn_atomdata_t.
857 static void nbnxn_atomdata_mask_fep(nbnxn_atomdata_t* nbat, const Nbnxm::GridSet& gridSet)
859 nbnxn_atomdata_t::Params& params = nbat->paramsDeprecated();
863 if (nbat->XFormat == nbatXYZQ)
865 q = nbat->x().data() + ZZ + 1;
866 stride_q = STRIDE_XYZQ;
874 for (const Nbnxm::Grid& grid : gridSet.grids())
877 if (grid.geometry().isSimple)
883 nsubc = c_gpuNumClusterPerCell;
886 int c_offset = grid.firstAtomInColumn(0);
888 /* Loop over all columns and copy and fill */
889 for (int c = 0; c < grid.numCells() * nsubc; c++)
891 /* Does this cluster contain perturbed particles? */
892 if (grid.clusterIsPerturbed(c))
894 const int numAtomsPerCluster = grid.geometry().numAtomsICluster;
895 for (int i = 0; i < numAtomsPerCluster; i++)
897 /* Is this a perturbed particle? */
898 if (grid.atomIsPerturbed(c, i))
900 int ind = c_offset + c * numAtomsPerCluster + i;
901 /* Set atom type and charge to non-interacting */
902 params.type[ind] = params.numTypes - 1;
903 q[ind * stride_q] = 0;
911 /* Copies the energy group indices to a reordered and packed array */
913 copy_egp_to_nbat_egps(const int* a, int na, int na_round, int na_c, int bit_shift, const int* in, int* innb)
919 for (i = 0; i < na; i += na_c)
921 /* Store na_c energy group numbers into one int */
923 for (int sa = 0; sa < na_c; sa++)
928 comb |= (GET_CGINFO_GID(in[at]) << (sa * bit_shift));
933 /* Complete the partially filled last cell with fill */
934 for (; i < na_round; i += na_c)
940 /* Set the energy group indices for atoms in nbnxn_atomdata_t */
941 static void nbnxn_atomdata_set_energygroups(nbnxn_atomdata_t::Params* params,
942 const Nbnxm::GridSet& gridSet,
943 ArrayRef<const int> atomInfo)
945 if (params->nenergrp == 1)
950 params->energrp.resize(gridSet.numGridAtomsTotal());
952 for (const Nbnxm::Grid& grid : gridSet.grids())
954 /* Loop over all columns and copy and fill */
955 for (int i = 0; i < grid.numColumns(); i++)
957 const int numAtoms = grid.paddedNumAtomsInColumn(i);
958 const int atomOffset = grid.firstAtomInColumn(i);
960 copy_egp_to_nbat_egps(gridSet.atomIndices().data() + atomOffset,
961 grid.numAtomsInColumn(i),
963 c_nbnxnCpuIClusterSize,
966 params->energrp.data() + grid.atomToCluster(atomOffset));
971 /* Sets all required atom parameter data in nbnxn_atomdata_t */
972 void nbnxn_atomdata_set(nbnxn_atomdata_t* nbat,
973 const Nbnxm::GridSet& gridSet,
974 ArrayRef<const int> atomTypes,
975 ArrayRef<const real> atomCharges,
976 ArrayRef<const int> atomInfo)
978 nbnxn_atomdata_t::Params& params = nbat->paramsDeprecated();
980 nbnxn_atomdata_set_atomtypes(¶ms, gridSet, atomTypes);
982 nbnxn_atomdata_set_charges(nbat, gridSet, atomCharges);
984 if (gridSet.haveFep())
986 nbnxn_atomdata_mask_fep(nbat, gridSet);
989 /* This must be done after masking types for FEP */
990 nbnxn_atomdata_set_ljcombparams(¶ms, nbat->XFormat, gridSet);
992 nbnxn_atomdata_set_energygroups(¶ms, gridSet, atomInfo);
995 /* Copies the shift vector array to nbnxn_atomdata_t */
996 void nbnxn_atomdata_copy_shiftvec(gmx_bool bDynamicBox, rvec* shift_vec, nbnxn_atomdata_t* nbat)
1000 nbat->bDynamicBox = bDynamicBox;
1001 for (i = 0; i < SHIFTS; i++)
1003 copy_rvec(shift_vec[i], nbat->shift_vec[i]);
1007 // This is slightly different from nbnxn_get_atom_range(...) at the end of the file
1008 // TODO: Combine if possible
1009 static void getAtomRanges(const Nbnxm::GridSet& gridSet,
1010 const gmx::AtomLocality locality,
1016 case gmx::AtomLocality::All:
1018 *gridEnd = gridSet.grids().size();
1020 case gmx::AtomLocality::Local:
1024 case gmx::AtomLocality::NonLocal:
1026 *gridEnd = gridSet.grids().size();
1028 case gmx::AtomLocality::Count:
1029 GMX_ASSERT(false, "Count is invalid locality specifier");
1034 /* Copies (and reorders) the coordinates to nbnxn_atomdata_t */
1035 void nbnxn_atomdata_copy_x_to_nbat_x(const Nbnxm::GridSet& gridSet,
1036 const gmx::AtomLocality locality,
1038 const rvec* coordinates,
1039 nbnxn_atomdata_t* nbat)
1044 getAtomRanges(gridSet, locality, &gridBegin, &gridEnd);
1048 nbat->natoms_local = gridSet.grids()[0].atomIndexEnd();
1051 const int nth = gmx_omp_nthreads_get(emntPairsearch);
1052 #pragma omp parallel for num_threads(nth) schedule(static)
1053 for (int th = 0; th < nth; th++)
1057 for (int g = gridBegin; g < gridEnd; g++)
1059 const Nbnxm::Grid& grid = gridSet.grids()[g];
1060 const int numCellsXY = grid.numColumns();
1062 const int cxy0 = (numCellsXY * th + nth - 1) / nth;
1063 const int cxy1 = (numCellsXY * (th + 1) + nth - 1) / nth;
1065 for (int cxy = cxy0; cxy < cxy1; cxy++)
1067 const int na = grid.numAtomsInColumn(cxy);
1068 const int ash = grid.firstAtomInColumn(cxy);
1071 if (g == 0 && fillLocal)
1073 na_fill = grid.paddedNumAtomsInColumn(cxy);
1077 /* We fill only the real particle locations.
1078 * We assume the filling entries at the end have been
1079 * properly set before during pair-list generation.
1083 copy_rvec_to_nbat_real(gridSet.atomIndices().data() + ash,
1093 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1097 /* Copies (and reorders) the coordinates to nbnxn_atomdata_t on the GPU*/
1098 void nbnxn_atomdata_x_to_nbat_x_gpu(const Nbnxm::GridSet& gridSet,
1099 const gmx::AtomLocality locality,
1102 DeviceBuffer<RVec> d_x,
1103 GpuEventSynchronizer* xReadyOnDevice)
1108 getAtomRanges(gridSet, locality, &gridBegin, &gridEnd);
1110 for (int g = gridBegin; g < gridEnd; g++)
1112 nbnxn_gpu_x_to_nbat_x(gridSet.grids()[g],
1113 fillLocal && g == 0,
1119 gridSet.numColumnsMax());
1123 static void nbnxn_atomdata_clear_reals(gmx::ArrayRef<real> dest, int i0, int i1)
1125 for (int i = i0; i < i1; i++)
1131 gmx_unused static void nbnxn_atomdata_reduce_reals(real* gmx_restrict dest,
1133 const real** gmx_restrict src,
1140 /* The destination buffer contains data, add to it */
1141 for (int i = i0; i < i1; i++)
1143 for (int s = 0; s < nsrc; s++)
1145 dest[i] += src[s][i];
1151 /* The destination buffer is unitialized, set it first */
1152 for (int i = i0; i < i1; i++)
1154 dest[i] = src[0][i];
1155 for (int s = 1; s < nsrc; s++)
1157 dest[i] += src[s][i];
1163 gmx_unused static void nbnxn_atomdata_reduce_reals_simd(real gmx_unused* gmx_restrict dest,
1164 gmx_bool gmx_unused bDestSet,
1165 const gmx_unused real** gmx_restrict src,
1166 int gmx_unused nsrc,
1171 /* The SIMD width here is actually independent of that in the kernels,
1172 * but we use the same width for simplicity (usually optimal anyhow).
1174 SimdReal dest_SSE, src_SSE;
1178 for (int i = i0; i < i1; i += GMX_SIMD_REAL_WIDTH)
1180 dest_SSE = load<SimdReal>(dest + i);
1181 for (int s = 0; s < nsrc; s++)
1183 src_SSE = load<SimdReal>(src[s] + i);
1184 dest_SSE = dest_SSE + src_SSE;
1186 store(dest + i, dest_SSE);
1191 for (int i = i0; i < i1; i += GMX_SIMD_REAL_WIDTH)
1193 dest_SSE = load<SimdReal>(src[0] + i);
1194 for (int s = 1; s < nsrc; s++)
1196 src_SSE = load<SimdReal>(src[s] + i);
1197 dest_SSE = dest_SSE + src_SSE;
1199 store(dest + i, dest_SSE);
1205 /* Add part of the force array(s) from nbnxn_atomdata_t to f
1207 * Note: Adding restrict to f makes this function 50% slower with gcc 7.3
1209 static void nbnxn_atomdata_add_nbat_f_to_f_part(const Nbnxm::GridSet& gridSet,
1210 const nbnxn_atomdata_t& nbat,
1211 const nbnxn_atomdata_output_t& out,
1216 gmx::ArrayRef<const int> cell = gridSet.cells();
1217 // Note: Using ArrayRef instead makes this code 25% slower with gcc 7.3
1218 const real* fnb = out.f.data();
1220 /* Loop over all columns and copy and fill */
1221 switch (nbat.FFormat)
1225 for (int a = a0; a < a1; a++)
1227 int i = cell[a] * nbat.fstride;
1230 f[a][YY] += fnb[i + 1];
1231 f[a][ZZ] += fnb[i + 2];
1235 for (int a = a0; a < a1; a++)
1237 int i = atom_to_x_index<c_packX4>(cell[a]);
1239 f[a][XX] += fnb[i + XX * c_packX4];
1240 f[a][YY] += fnb[i + YY * c_packX4];
1241 f[a][ZZ] += fnb[i + ZZ * c_packX4];
1245 for (int a = a0; a < a1; a++)
1247 int i = atom_to_x_index<c_packX8>(cell[a]);
1249 f[a][XX] += fnb[i + XX * c_packX8];
1250 f[a][YY] += fnb[i + YY * c_packX8];
1251 f[a][ZZ] += fnb[i + ZZ * c_packX8];
1254 default: gmx_incons("Unsupported nbnxn_atomdata_t format");
1258 static inline unsigned char reverse_bits(unsigned char b)
1260 /* http://graphics.stanford.edu/~seander/bithacks.html#ReverseByteWith64BitsDiv */
1261 return (b * 0x0202020202ULL & 0x010884422010ULL) % 1023;
1264 static void nbnxn_atomdata_add_nbat_f_to_f_treereduce(nbnxn_atomdata_t* nbat, int nth)
1266 gmx::ArrayRef<const gmx_bitmask_t> flags = nbat->buffer_flags;
1268 int next_pow2 = 1 << (gmx::log2I(nth - 1) + 1);
1270 const int numOutputBuffers = nbat->out.size();
1271 GMX_ASSERT(numOutputBuffers == nth,
1272 "tree-reduce currently only works for numOutputBuffers==nth");
1274 memset(nbat->syncStep, 0, sizeof(*(nbat->syncStep)) * nth);
1276 #pragma omp parallel num_threads(nth)
1284 th = gmx_omp_get_thread_num();
1286 for (group_size = 2; group_size < 2 * next_pow2; group_size *= 2)
1288 int index[2], group_pos, partner_pos, wu;
1289 int partner_th = th ^ (group_size / 2);
1294 /* wait on partner thread - replaces full barrier */
1295 int sync_th, sync_group_size;
1297 # if defined(__clang__) && __clang_major__ >= 8
1298 // Suppress warnings that the use of memory_barrier may be excessive
1299 // Only exists beginning with clang-8
1300 # pragma clang diagnostic push
1301 # pragma clang diagnostic ignored "-Watomic-implicit-seq-cst"
1304 tMPI_Atomic_memory_barrier(); /* guarantee data is saved before marking work as done */
1305 tMPI_Atomic_set(&(nbat->syncStep[th]), group_size / 2); /* mark previous step as completed */
1307 /* find thread to sync with. Equal to partner_th unless nth is not a power of two. */
1308 for (sync_th = partner_th, sync_group_size = group_size;
1309 sync_th >= nth && sync_group_size > 2;
1310 sync_group_size /= 2)
1312 sync_th &= ~(sync_group_size / 4);
1314 if (sync_th < nth) /* otherwise nothing to sync index[1] will be >=nout */
1316 /* wait on the thread which computed input data in previous step */
1317 while (tMPI_Atomic_get(static_cast<volatile tMPI_Atomic_t*>(&(nbat->syncStep[sync_th])))
1322 /* guarantee that no later load happens before wait loop is finisehd */
1323 tMPI_Atomic_memory_barrier();
1325 # if defined(__clang__) && __clang_major__ >= 8
1326 # pragma clang diagnostic pop
1328 #else /* TMPI_ATOMICS */
1329 # pragma omp barrier
1333 /* Calculate buffers to sum (result goes into first buffer) */
1334 group_pos = th % group_size;
1335 index[0] = th - group_pos;
1336 index[1] = index[0] + group_size / 2;
1338 /* If no second buffer, nothing to do */
1339 if (index[1] >= numOutputBuffers && group_size > 2)
1344 #if NBNXN_BUFFERFLAG_MAX_THREADS > 256
1345 # error reverse_bits assumes max 256 threads
1347 /* Position is permuted so that one of the 2 vectors being added was computed on the same thread in the previous step.
1348 This improves locality and enables to sync with just a single thread between steps (=the levels in the btree).
1349 The permutation which allows this corresponds to reversing the bits of the group position.
1351 group_pos = reverse_bits(group_pos) / (256 / group_size);
1353 partner_pos = group_pos ^ 1;
1355 /* loop over two work-units (own and partner) */
1356 for (wu = 0; wu < 2; wu++)
1360 if (partner_th < nth)
1362 break; /* partner exists we don't have to do his work */
1366 group_pos = partner_pos;
1370 /* Calculate the cell-block range for our thread */
1371 b0 = (flags.size() * group_pos) / group_size;
1372 b1 = (flags.size() * (group_pos + 1)) / group_size;
1374 for (b = b0; b < b1; b++)
1376 i0 = b * NBNXN_BUFFERFLAG_SIZE * nbat->fstride;
1377 i1 = (b + 1) * NBNXN_BUFFERFLAG_SIZE * nbat->fstride;
1379 if (bitmask_is_set(flags[b], index[1]) || group_size > 2)
1381 const real* fIndex1 = nbat->out[index[1]].f.data();
1383 nbnxn_atomdata_reduce_reals_simd
1385 nbnxn_atomdata_reduce_reals
1387 (nbat->out[index[0]].f.data(),
1388 bitmask_is_set(flags[b], index[0]) || group_size > 2,
1394 else if (!bitmask_is_set(flags[b], index[0]))
1396 nbnxn_atomdata_clear_reals(nbat->out[index[0]].f, i0, i1);
1402 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1407 static void nbnxn_atomdata_add_nbat_f_to_f_stdreduce(nbnxn_atomdata_t* nbat, int nth)
1409 #pragma omp parallel for num_threads(nth) schedule(static)
1410 for (int th = 0; th < nth; th++)
1415 const real* fptr[NBNXN_BUFFERFLAG_MAX_THREADS];
1417 gmx::ArrayRef<const gmx_bitmask_t> flags = nbat->buffer_flags;
1419 /* Calculate the cell-block range for our thread */
1420 int b0 = (flags.size() * th) / nth;
1421 int b1 = (flags.size() * (th + 1)) / nth;
1423 for (int b = b0; b < b1; b++)
1425 int i0 = b * NBNXN_BUFFERFLAG_SIZE * nbat->fstride;
1426 int i1 = (b + 1) * NBNXN_BUFFERFLAG_SIZE * nbat->fstride;
1429 for (gmx::index out = 1; out < gmx::ssize(nbat->out); out++)
1431 if (bitmask_is_set(flags[b], out))
1433 fptr[nfptr++] = nbat->out[out].f.data();
1439 nbnxn_atomdata_reduce_reals_simd
1441 nbnxn_atomdata_reduce_reals
1443 (nbat->out[0].f.data(), bitmask_is_set(flags[b], 0), fptr, nfptr, i0, i1);
1445 else if (!bitmask_is_set(flags[b], 0))
1447 nbnxn_atomdata_clear_reals(nbat->out[0].f, i0, i1);
1451 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1456 /* Add the force array(s) from nbnxn_atomdata_t to f */
1457 void reduceForces(nbnxn_atomdata_t* nbat, const gmx::AtomLocality locality, const Nbnxm::GridSet& gridSet, rvec* f)
1462 nbnxn_get_atom_range(locality, gridSet, &a0, &na);
1466 /* The are no atoms for this reduction, avoid some overhead */
1470 int nth = gmx_omp_nthreads_get(emntNonbonded);
1472 if (nbat->out.size() > 1)
1474 if (locality != gmx::AtomLocality::All)
1476 gmx_incons("add_f_to_f called with nout>1 and locality!=eatAll");
1479 /* Reduce the force thread output buffers into buffer 0, before adding
1480 * them to the, differently ordered, "real" force buffer.
1482 if (nbat->bUseTreeReduce)
1484 nbnxn_atomdata_add_nbat_f_to_f_treereduce(nbat, nth);
1488 nbnxn_atomdata_add_nbat_f_to_f_stdreduce(nbat, nth);
1491 #pragma omp parallel for num_threads(nth) schedule(static)
1492 for (int th = 0; th < nth; th++)
1496 nbnxn_atomdata_add_nbat_f_to_f_part(
1497 gridSet, *nbat, nbat->out[0], a0 + ((th + 0) * na) / nth, a0 + ((th + 1) * na) / nth, f);
1499 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1503 void nbnxn_atomdata_add_nbat_fshift_to_fshift(const nbnxn_atomdata_t& nbat, gmx::ArrayRef<gmx::RVec> fshift)
1505 gmx::ArrayRef<const nbnxn_atomdata_output_t> outputBuffers = nbat.out;
1507 for (int s = 0; s < SHIFTS; s++)
1511 for (const nbnxn_atomdata_output_t& out : outputBuffers)
1513 sum[XX] += out.fshift[s * DIM + XX];
1514 sum[YY] += out.fshift[s * DIM + YY];
1515 sum[ZZ] += out.fshift[s * DIM + ZZ];
1521 void nbnxn_get_atom_range(const gmx::AtomLocality atomLocality,
1522 const Nbnxm::GridSet& gridSet,
1527 switch (atomLocality)
1529 case gmx::AtomLocality::All:
1531 *nAtoms = gridSet.numRealAtomsTotal();
1533 case gmx::AtomLocality::Local:
1535 *nAtoms = gridSet.numRealAtomsLocal();
1537 case gmx::AtomLocality::NonLocal:
1538 *atomStart = gridSet.numRealAtomsLocal();
1539 *nAtoms = gridSet.numRealAtomsTotal() - gridSet.numRealAtomsLocal();
1541 case gmx::AtomLocality::Count:
1542 GMX_ASSERT(false, "Count is invalid locality specifier");