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 "gromacs/math/functions.h"
49 #include "gromacs/math/utilities.h"
50 #include "gromacs/math/vec.h"
51 #include "gromacs/mdlib/gmx_omp_nthreads.h"
52 #include "gromacs/mdtypes/forcerec.h" // only for GET_CGINFO_*
53 #include "gromacs/mdtypes/md_enums.h"
54 #include "gromacs/nbnxm/nbnxm.h"
55 #include "gromacs/pbcutil/ishift.h"
56 #include "gromacs/simd/simd.h"
57 #include "gromacs/utility/exceptions.h"
58 #include "gromacs/utility/fatalerror.h"
59 #include "gromacs/utility/gmxomp.h"
60 #include "gromacs/utility/logger.h"
61 #include "gromacs/utility/smalloc.h"
62 #include "gromacs/utility/strconvert.h"
63 #include "gromacs/utility/stringutil.h"
67 #include "nbnxm_geometry.h"
68 #include "nbnxm_gpu.h"
71 using namespace gmx; // TODO: Remove when this file is moved into gmx namespace
74 const char* enumValueToString(LJCombinationRule enumValue)
76 static constexpr gmx::EnumerationArray<LJCombinationRule, const char*> s_ljCombinationRuleNames = {
77 "Geometric", "Lorentz-Berthelot", "None"
79 return s_ljCombinationRuleNames[enumValue];
82 void nbnxn_atomdata_t::resizeCoordinateBuffer(int numAtoms)
86 x_.resize(numAtoms * xstride);
89 void nbnxn_atomdata_t::resizeForceBuffers()
91 /* Force buffers need padding up to a multiple of the buffer flag size */
92 const int paddedSize =
93 (numAtoms() + NBNXN_BUFFERFLAG_SIZE - 1) / NBNXN_BUFFERFLAG_SIZE * NBNXN_BUFFERFLAG_SIZE;
95 /* Should we let each thread allocate it's own data instead? */
96 for (nbnxn_atomdata_output_t& outBuffer : out)
98 outBuffer.f.resize(paddedSize * fstride);
102 /* Initializes an nbnxn_atomdata_output_t data structure */
103 nbnxn_atomdata_output_t::nbnxn_atomdata_output_t(Nbnxm::KernelType kernelType,
105 int simdEnergyBufferStride,
106 gmx::PinningPolicy pinningPolicy) :
107 f({}, { pinningPolicy }),
108 fshift({}, { pinningPolicy }),
109 Vvdw({}, { pinningPolicy }),
110 Vc({}, { pinningPolicy })
112 fshift.resize(SHIFTS * DIM);
113 Vvdw.resize(numEnergyGroups * numEnergyGroups);
114 Vc.resize(numEnergyGroups * numEnergyGroups);
116 if (Nbnxm::kernelTypeIsSimd(kernelType))
118 int cj_size = Nbnxm::JClusterSizePerKernelType[kernelType];
120 numEnergyGroups * numEnergyGroups * simdEnergyBufferStride * (cj_size / 2) * cj_size;
121 VSvdw.resize(numElements);
122 VSc.resize(numElements);
126 static void copy_int_to_nbat_int(const int* a, int na, int na_round, const int* in, int fill, int* innb)
131 for (i = 0; i < na; i++)
133 innb[j++] = in[a[i]];
135 /* Complete the partially filled last cell with fill */
136 for (; i < na_round; i++)
142 void copy_rvec_to_nbat_real(const int* a, int na, int na_round, const rvec* x, int nbatFormat, real* xnb, int a0)
144 /* We complete partially filled cells, can only be the last one in each
145 * column, with coordinates farAway. The actual coordinate value does
146 * not influence the results, since these filler particles do not interact.
147 * Clusters with normal atoms + fillers have a bounding box based only
148 * on the coordinates of the atoms. Clusters with only fillers have as
149 * the bounding box the coordinates of the first filler. Such clusters
150 * are not considered as i-entries, but they are considered as j-entries.
151 * So for performance it is better to have their bounding boxes far away,
152 * such that filler only clusters don't end up in the pair list.
154 const real farAway = -1000000;
162 for (i = 0; i < na; i++)
164 xnb[j++] = x[a[i]][XX];
165 xnb[j++] = x[a[i]][YY];
166 xnb[j++] = x[a[i]][ZZ];
168 /* Complete the partially filled last cell with farAway elements */
169 for (; i < na_round; i++)
177 j = a0 * STRIDE_XYZQ;
178 for (i = 0; i < na; i++)
180 xnb[j++] = x[a[i]][XX];
181 xnb[j++] = x[a[i]][YY];
182 xnb[j++] = x[a[i]][ZZ];
185 /* Complete the partially filled last cell with zeros */
186 for (; i < na_round; i++)
195 j = atom_to_x_index<c_packX4>(a0);
196 c = a0 & (c_packX4 - 1);
197 for (i = 0; i < na; i++)
199 xnb[j + XX * c_packX4] = x[a[i]][XX];
200 xnb[j + YY * c_packX4] = x[a[i]][YY];
201 xnb[j + ZZ * c_packX4] = x[a[i]][ZZ];
206 j += (DIM - 1) * c_packX4;
210 /* Complete the partially filled last cell with zeros */
211 for (; i < na_round; i++)
213 xnb[j + XX * c_packX4] = farAway;
214 xnb[j + YY * c_packX4] = farAway;
215 xnb[j + ZZ * c_packX4] = farAway;
220 j += (DIM - 1) * c_packX4;
226 j = atom_to_x_index<c_packX8>(a0);
227 c = a0 & (c_packX8 - 1);
228 for (i = 0; i < na; i++)
230 xnb[j + XX * c_packX8] = x[a[i]][XX];
231 xnb[j + YY * c_packX8] = x[a[i]][YY];
232 xnb[j + ZZ * c_packX8] = x[a[i]][ZZ];
237 j += (DIM - 1) * c_packX8;
241 /* Complete the partially filled last cell with zeros */
242 for (; i < na_round; i++)
244 xnb[j + XX * c_packX8] = farAway;
245 xnb[j + YY * c_packX8] = farAway;
246 xnb[j + ZZ * c_packX8] = farAway;
251 j += (DIM - 1) * c_packX8;
256 default: gmx_incons("Unsupported nbnxn_atomdata_t format");
260 /* Stores the LJ parameter data in a format convenient for different kernels */
261 static void set_lj_parameter_data(nbnxn_atomdata_t::Params* params, gmx_bool bSIMD)
263 int nt = params->numTypes;
268 /* nbfp_aligned stores two parameters using the stride most suitable
269 * for the present SIMD architecture, as specified by the constant
270 * c_simdBestPairAlignment from the SIMD header.
271 * There's a slight inefficiency in allocating and initializing nbfp_aligned
272 * when it might not be used, but introducing the conditional code is not
275 params->nbfp_aligned.resize(nt * nt * c_simdBestPairAlignment);
277 for (int i = 0; i < nt; i++)
279 for (int j = 0; j < nt; j++)
281 params->nbfp_aligned[(i * nt + j) * c_simdBestPairAlignment + 0] =
282 params->nbfp[(i * nt + j) * 2 + 0];
283 params->nbfp_aligned[(i * nt + j) * c_simdBestPairAlignment + 1] =
284 params->nbfp[(i * nt + j) * 2 + 1];
285 if (c_simdBestPairAlignment > 2)
287 params->nbfp_aligned[(i * nt + j) * c_simdBestPairAlignment + 2] = 0;
288 params->nbfp_aligned[(i * nt + j) * c_simdBestPairAlignment + 3] = 0;
295 /* We use combination rule data for SIMD combination rule kernels
296 * and with LJ-PME kernels. We then only need parameters per atom type,
297 * not per pair of atom types.
299 params->nbfp_comb.resize(nt * 2);
300 switch (params->ljCombinationRule)
302 case LJCombinationRule::Geometric:
303 for (int i = 0; i < nt; i++)
305 /* Store the sqrt of the diagonal from the nbfp matrix */
306 params->nbfp_comb[i * 2] = std::sqrt(params->nbfp[(i * nt + i) * 2]);
307 params->nbfp_comb[i * 2 + 1] = std::sqrt(params->nbfp[(i * nt + i) * 2 + 1]);
310 case LJCombinationRule::LorentzBerthelot:
311 for (int i = 0; i < nt; i++)
313 /* Get 6*C6 and 12*C12 from the diagonal of the nbfp matrix */
314 const real c6 = params->nbfp[(i * nt + i) * 2];
315 const real c12 = params->nbfp[(i * nt + i) * 2 + 1];
316 if (c6 > 0 && c12 > 0)
318 /* We store 0.5*2^1/6*sigma and sqrt(4*3*eps),
319 * so we get 6*C6 and 12*C12 after combining.
321 params->nbfp_comb[i * 2] = 0.5 * gmx::sixthroot(c12 / c6);
322 params->nbfp_comb[i * 2 + 1] = std::sqrt(c6 * c6 / c12);
326 params->nbfp_comb[i * 2] = 0;
327 params->nbfp_comb[i * 2 + 1] = 0;
331 case LJCombinationRule::None:
332 /* We always store the full matrix (see code above) */
334 default: gmx_incons("Unknown combination rule");
338 nbnxn_atomdata_t::SimdMasks::SimdMasks()
341 constexpr int simd_width = GMX_SIMD_REAL_WIDTH;
342 /* Set the diagonal cluster pair exclusion mask setup data.
343 * In the kernel we check 0 < j - i to generate the masks.
344 * Here we store j - i for generating the mask for the first i,
345 * we substract 0.5 to avoid rounding issues.
346 * In the kernel we can subtract 1 to generate the subsequent mask.
348 const int simd_4xn_diag_size = std::max(c_nbnxnCpuIClusterSize, simd_width);
349 diagonal_4xn_j_minus_i.resize(simd_4xn_diag_size);
350 for (int j = 0; j < simd_4xn_diag_size; j++)
352 diagonal_4xn_j_minus_i[j] = j - 0.5;
355 diagonal_2xnn_j_minus_i.resize(simd_width);
356 for (int j = 0; j < simd_width / 2; j++)
358 /* The j-cluster size is half the SIMD width */
359 diagonal_2xnn_j_minus_i[j] = j - 0.5;
360 /* The next half of the SIMD width is for i + 1 */
361 diagonal_2xnn_j_minus_i[simd_width / 2 + j] = j - 1 - 0.5;
364 /* We use up to 32 bits for exclusion masking.
365 * The same masks are used for the 4xN and 2x(N+N) kernels.
366 * The masks are read either into integer SIMD registers or into
367 * real SIMD registers (together with a cast).
368 * In single precision this means the real and integer SIMD registers
371 const int simd_excl_size = c_nbnxnCpuIClusterSize * simd_width;
372 # if GMX_DOUBLE && !GMX_SIMD_HAVE_INT32_LOGICAL
373 exclusion_filter64.resize(simd_excl_size);
375 exclusion_filter.resize(simd_excl_size);
378 for (int j = 0; j < simd_excl_size; j++)
380 /* Set the consecutive bits for masking pair exclusions */
381 # if GMX_DOUBLE && !GMX_SIMD_HAVE_INT32_LOGICAL
382 exclusion_filter64[j] = (1U << j);
384 exclusion_filter[j] = (1U << j);
388 if (!GMX_SIMD_HAVE_LOGICAL && !GMX_SIMD_HAVE_INT32_LOGICAL) // NOLINT(misc-redundant-expression)
390 // If the SIMD implementation has no bitwise logical operation support
391 // whatsoever we cannot use the normal masking. Instead,
392 // we generate a vector of all 2^4 possible ways an i atom
393 // interacts with its 4 j atoms. Each array entry contains
394 // GMX_SIMD_REAL_WIDTH values that are read with a single aligned SIMD load.
395 // Since there is no logical value representation in this case, we use
396 // any nonzero value to indicate 'true', while zero mean 'false'.
397 // This can then be converted to a SIMD boolean internally in the SIMD
398 // module by comparing to zero.
399 // Each array entry encodes how this i atom will interact with the 4 j atoms.
400 // Matching code exists in set_ci_top_excls() to generate indices into this array.
401 // Those indices are used in the kernels.
403 const int simd_excl_size = c_nbnxnCpuIClusterSize * c_nbnxnCpuIClusterSize;
404 const real simdFalse = 0.0;
405 const real simdTrue = 1.0;
407 interaction_array.resize(simd_excl_size * GMX_SIMD_REAL_WIDTH);
408 for (int j = 0; j < simd_excl_size; j++)
410 const int index = j * GMX_SIMD_REAL_WIDTH;
411 for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
413 interaction_array[index + i] = (j & (1 << i)) ? simdTrue : simdFalse;
420 nbnxn_atomdata_t::Params::Params(gmx::PinningPolicy pinningPolicy) :
422 nbfp({}, { pinningPolicy }),
423 nbfp_comb({}, { pinningPolicy }),
424 type({}, { pinningPolicy }),
425 lj_comb({}, { pinningPolicy }),
426 q({}, { pinningPolicy }),
429 energrp({}, { pinningPolicy })
433 nbnxn_atomdata_t::nbnxn_atomdata_t(gmx::PinningPolicy pinningPolicy) :
434 params_(pinningPolicy),
437 shift_vec({}, { pinningPolicy }),
438 x_({}, { pinningPolicy }),
440 bUseBufferFlags(FALSE)
444 /* Initializes an nbnxn_atomdata_t::Params data structure */
445 static void nbnxn_atomdata_params_init(const gmx::MDLogger& mdlog,
446 nbnxn_atomdata_t::Params* params,
447 const Nbnxm::KernelType kernelType,
448 int enbnxninitcombrule,
450 ArrayRef<const real> nbfp,
455 gmx_bool simple, bCombGeom, bCombLB, bSIMD;
459 fprintf(debug, "There are %d atom types in the system, adding one for nbnxn_atomdata_t\n", ntype);
461 params->numTypes = ntype + 1;
462 params->nbfp.resize(params->numTypes * params->numTypes * 2);
463 params->nbfp_comb.resize(params->numTypes * 2);
465 /* A tolerance of 1e-5 seems reasonable for (possibly hand-typed)
466 * force-field floating point parameters.
469 ptr = getenv("GMX_LJCOMB_TOL");
474 sscanf(ptr, "%lf", &dbl);
480 /* Temporarily fill params->nbfp_comb with sigma and epsilon
481 * to check for the LB rule.
483 for (int i = 0; i < ntype; i++)
485 c6 = nbfp[(i * ntype + i) * 2] / 6.0;
486 c12 = nbfp[(i * ntype + i) * 2 + 1] / 12.0;
487 if (c6 > 0 && c12 > 0)
489 params->nbfp_comb[i * 2] = gmx::sixthroot(c12 / c6);
490 params->nbfp_comb[i * 2 + 1] = 0.25 * c6 * c6 / c12;
492 else if (c6 == 0 && c12 == 0)
494 params->nbfp_comb[i * 2] = 0;
495 params->nbfp_comb[i * 2 + 1] = 0;
499 /* Can not use LB rule with only dispersion or repulsion */
504 for (int i = 0; i < params->numTypes; i++)
506 for (int j = 0; j < params->numTypes; j++)
508 if (i < ntype && j < ntype)
510 /* fr->nbfp has been updated, so that array too now stores c6/c12 including
511 * the 6.0/12.0 prefactors to save 2 flops in the most common case (force-only).
513 c6 = nbfp[(i * ntype + j) * 2];
514 c12 = nbfp[(i * ntype + j) * 2 + 1];
515 params->nbfp[(i * params->numTypes + j) * 2] = c6;
516 params->nbfp[(i * params->numTypes + j) * 2 + 1] = c12;
518 /* Compare 6*C6 and 12*C12 for geometric cobination rule */
521 && gmx_within_tol(c6 * c6, nbfp[(i * ntype + i) * 2] * nbfp[(j * ntype + j) * 2], tol)
522 && gmx_within_tol(c12 * c12,
523 nbfp[(i * ntype + i) * 2 + 1] * nbfp[(j * ntype + j) * 2 + 1],
526 /* Compare C6 and C12 for Lorentz-Berthelot combination rule */
531 && ((c6 == 0 && c12 == 0
532 && (params->nbfp_comb[i * 2 + 1] == 0 || params->nbfp_comb[j * 2 + 1] == 0))
533 || (c6 > 0 && c12 > 0
534 && gmx_within_tol(gmx::sixthroot(c12 / c6),
535 0.5 * (params->nbfp_comb[i * 2] + params->nbfp_comb[j * 2]),
537 && gmx_within_tol(0.25 * c6 * c6 / c12,
538 std::sqrt(params->nbfp_comb[i * 2 + 1]
539 * params->nbfp_comb[j * 2 + 1]),
544 /* Add zero parameters for the additional dummy atom type */
545 params->nbfp[(i * params->numTypes + j) * 2] = 0;
546 params->nbfp[(i * params->numTypes + j) * 2 + 1] = 0;
553 "Combination rules: geometric %s Lorentz-Berthelot %s\n",
554 gmx::boolToString(bCombGeom),
555 gmx::boolToString(bCombLB));
558 simple = Nbnxm::kernelTypeUsesSimplePairlist(kernelType);
560 switch (enbnxninitcombrule)
562 case enbnxninitcombruleDETECT:
563 /* We prefer the geometic combination rule,
564 * as that gives a slightly faster kernel than the LB rule.
568 params->ljCombinationRule = LJCombinationRule::Geometric;
572 params->ljCombinationRule = LJCombinationRule::LorentzBerthelot;
576 params->ljCombinationRule = LJCombinationRule::None;
578 params->nbfp_comb.clear();
583 if (params->ljCombinationRule == LJCombinationRule::None)
585 mesg = "Using full Lennard-Jones parameter combination matrix";
589 mesg = gmx::formatString("Using %s Lennard-Jones combination rule",
590 enumValueToString(params->ljCombinationRule));
592 GMX_LOG(mdlog.info).asParagraph().appendText(mesg);
595 case enbnxninitcombruleGEOM:
596 params->ljCombinationRule = LJCombinationRule::Geometric;
598 case enbnxninitcombruleLB:
599 params->ljCombinationRule = LJCombinationRule::LorentzBerthelot;
601 case enbnxninitcombruleNONE:
602 params->ljCombinationRule = LJCombinationRule::None;
604 params->nbfp_comb.clear();
606 default: gmx_incons("Unknown enbnxninitcombrule");
609 bSIMD = Nbnxm::kernelTypeIsSimd(kernelType);
611 set_lj_parameter_data(params, bSIMD);
613 params->nenergrp = n_energygroups;
616 // We now check for energy groups already when starting mdrun
617 GMX_RELEASE_ASSERT(n_energygroups == 1, "GPU kernels do not support energy groups");
619 /* Temporary storage goes as #grp^3*simd_width^2/2, so limit to 64 */
620 if (params->nenergrp > 64)
622 gmx_fatal(FARGS, "With NxN kernels not more than 64 energy groups are supported\n");
624 params->neg_2log = 1;
625 while (params->nenergrp > (1 << params->neg_2log))
631 /* Initializes an nbnxn_atomdata_t data structure */
632 void nbnxn_atomdata_init(const gmx::MDLogger& mdlog,
633 nbnxn_atomdata_t* nbat,
634 const Nbnxm::KernelType kernelType,
635 int enbnxninitcombrule,
637 ArrayRef<const real> nbfp,
641 nbnxn_atomdata_params_init(
642 mdlog, &nbat->paramsDeprecated(), kernelType, enbnxninitcombrule, ntype, nbfp, n_energygroups);
644 const bool simple = Nbnxm::kernelTypeUsesSimplePairlist(kernelType);
645 const bool bSIMD = Nbnxm::kernelTypeIsSimd(kernelType);
653 pack_x = std::max(c_nbnxnCpuIClusterSize, Nbnxm::JClusterSizePerKernelType[kernelType]);
656 case 4: nbat->XFormat = nbatX4; break;
657 case 8: nbat->XFormat = nbatX8; break;
658 default: gmx_incons("Unsupported packing width");
663 nbat->XFormat = nbatXYZ;
666 nbat->FFormat = nbat->XFormat;
670 nbat->XFormat = nbatXYZQ;
671 nbat->FFormat = nbatXYZ;
674 nbat->shift_vec.resize(SHIFTS);
676 nbat->xstride = (nbat->XFormat == nbatXYZQ ? STRIDE_XYZQ : DIM);
677 nbat->fstride = (nbat->FFormat == nbatXYZQ ? STRIDE_XYZQ : DIM);
679 /* Initialize the output data structures */
680 for (int i = 0; i < nout; i++)
682 const auto& pinningPolicy = nbat->params().type.get_allocator().pinningPolicy();
683 nbat->out.emplace_back(
684 kernelType, nbat->params().nenergrp, 1 << nbat->params().neg_2log, pinningPolicy);
687 nbat->buffer_flags.clear();
690 template<int packSize>
691 static void copy_lj_to_nbat_lj_comb(gmx::ArrayRef<const real> ljparam_type, const int* type, int na, real* ljparam_at)
693 /* The LJ params follow the combination rule:
694 * copy the params for the type array to the atom array.
696 for (int is = 0; is < na; is += packSize)
698 for (int k = 0; k < packSize; k++)
701 ljparam_at[is * 2 + k] = ljparam_type[type[i] * 2];
702 ljparam_at[is * 2 + packSize + k] = ljparam_type[type[i] * 2 + 1];
707 /* Sets the atom type in nbnxn_atomdata_t */
708 static void nbnxn_atomdata_set_atomtypes(nbnxn_atomdata_t::Params* params,
709 const Nbnxm::GridSet& gridSet,
710 ArrayRef<const int> atomTypes)
712 params->type.resize(gridSet.numGridAtomsTotal());
714 for (const Nbnxm::Grid& grid : gridSet.grids())
716 /* Loop over all columns and copy and fill */
717 for (int i = 0; i < grid.numColumns(); i++)
719 const int numAtoms = grid.paddedNumAtomsInColumn(i);
720 const int atomOffset = grid.firstAtomInColumn(i);
722 copy_int_to_nbat_int(gridSet.atomIndices().data() + atomOffset,
723 grid.numAtomsInColumn(i),
726 params->numTypes - 1,
727 params->type.data() + atomOffset);
732 /* Sets the LJ combination rule parameters in nbnxn_atomdata_t */
733 static void nbnxn_atomdata_set_ljcombparams(nbnxn_atomdata_t::Params* params,
735 const Nbnxm::GridSet& gridSet)
737 params->lj_comb.resize(gridSet.numGridAtomsTotal() * 2);
739 if (params->ljCombinationRule != LJCombinationRule::None)
741 for (const Nbnxm::Grid& grid : gridSet.grids())
743 /* Loop over all columns and copy and fill */
744 for (int i = 0; i < grid.numColumns(); i++)
746 const int numAtoms = grid.paddedNumAtomsInColumn(i);
747 const int atomOffset = grid.firstAtomInColumn(i);
749 if (XFormat == nbatX4)
751 copy_lj_to_nbat_lj_comb<c_packX4>(params->nbfp_comb,
752 params->type.data() + atomOffset,
754 params->lj_comb.data() + atomOffset * 2);
756 else if (XFormat == nbatX8)
758 copy_lj_to_nbat_lj_comb<c_packX8>(params->nbfp_comb,
759 params->type.data() + atomOffset,
761 params->lj_comb.data() + atomOffset * 2);
763 else if (XFormat == nbatXYZQ)
765 copy_lj_to_nbat_lj_comb<1>(params->nbfp_comb,
766 params->type.data() + atomOffset,
768 params->lj_comb.data() + atomOffset * 2);
775 /* Sets the charges in nbnxn_atomdata_t *nbat */
776 static void nbnxn_atomdata_set_charges(nbnxn_atomdata_t* nbat,
777 const Nbnxm::GridSet& gridSet,
778 ArrayRef<const real> charges)
780 if (nbat->XFormat != nbatXYZQ)
782 nbat->paramsDeprecated().q.resize(nbat->numAtoms());
785 for (const Nbnxm::Grid& grid : gridSet.grids())
787 /* Loop over all columns and copy and fill */
788 for (int cxy = 0; cxy < grid.numColumns(); cxy++)
790 const int atomOffset = grid.firstAtomInColumn(cxy);
791 const int numAtoms = grid.numAtomsInColumn(cxy);
792 const int paddedNumAtoms = grid.paddedNumAtomsInColumn(cxy);
794 if (nbat->XFormat == nbatXYZQ)
796 real* q = nbat->x().data() + atomOffset * STRIDE_XYZQ + ZZ + 1;
798 for (i = 0; i < numAtoms; i++)
800 *q = charges[gridSet.atomIndices()[atomOffset + i]];
803 /* Complete the partially filled last cell with zeros */
804 for (; i < paddedNumAtoms; i++)
812 real* q = nbat->paramsDeprecated().q.data() + atomOffset;
814 for (i = 0; i < numAtoms; i++)
816 *q = charges[gridSet.atomIndices()[atomOffset + i]];
819 /* Complete the partially filled last cell with zeros */
820 for (; i < paddedNumAtoms; i++)
830 /* Set the charges of perturbed atoms in nbnxn_atomdata_t to 0.
831 * This is to automatically remove the RF/PME self term in the nbnxn kernels.
832 * Part of the zero interactions are still calculated in the normal kernels.
833 * All perturbed interactions are calculated in the free energy kernel,
834 * using the original charge and LJ data, not nbnxn_atomdata_t.
836 static void nbnxn_atomdata_mask_fep(nbnxn_atomdata_t* nbat, const Nbnxm::GridSet& gridSet)
838 nbnxn_atomdata_t::Params& params = nbat->paramsDeprecated();
842 if (nbat->XFormat == nbatXYZQ)
844 q = nbat->x().data() + ZZ + 1;
845 stride_q = STRIDE_XYZQ;
853 for (const Nbnxm::Grid& grid : gridSet.grids())
856 if (grid.geometry().isSimple)
862 nsubc = c_gpuNumClusterPerCell;
865 int c_offset = grid.firstAtomInColumn(0);
867 /* Loop over all columns and copy and fill */
868 for (int c = 0; c < grid.numCells() * nsubc; c++)
870 /* Does this cluster contain perturbed particles? */
871 if (grid.clusterIsPerturbed(c))
873 const int numAtomsPerCluster = grid.geometry().numAtomsICluster;
874 for (int i = 0; i < numAtomsPerCluster; i++)
876 /* Is this a perturbed particle? */
877 if (grid.atomIsPerturbed(c, i))
879 int ind = c_offset + c * numAtomsPerCluster + i;
880 /* Set atom type and charge to non-interacting */
881 params.type[ind] = params.numTypes - 1;
882 q[ind * stride_q] = 0;
890 /* Copies the energy group indices to a reordered and packed array */
892 copy_egp_to_nbat_egps(const int* a, int na, int na_round, int na_c, int bit_shift, const int* in, int* innb)
898 for (i = 0; i < na; i += na_c)
900 /* Store na_c energy group numbers into one int */
902 for (int sa = 0; sa < na_c; sa++)
907 comb |= (GET_CGINFO_GID(in[at]) << (sa * bit_shift));
912 /* Complete the partially filled last cell with fill */
913 for (; i < na_round; i += na_c)
919 /* Set the energy group indices for atoms in nbnxn_atomdata_t */
920 static void nbnxn_atomdata_set_energygroups(nbnxn_atomdata_t::Params* params,
921 const Nbnxm::GridSet& gridSet,
922 ArrayRef<const int> atomInfo)
924 if (params->nenergrp == 1)
929 params->energrp.resize(gridSet.numGridAtomsTotal());
931 for (const Nbnxm::Grid& grid : gridSet.grids())
933 /* Loop over all columns and copy and fill */
934 for (int i = 0; i < grid.numColumns(); i++)
936 const int numAtoms = grid.paddedNumAtomsInColumn(i);
937 const int atomOffset = grid.firstAtomInColumn(i);
939 copy_egp_to_nbat_egps(gridSet.atomIndices().data() + atomOffset,
940 grid.numAtomsInColumn(i),
942 c_nbnxnCpuIClusterSize,
945 params->energrp.data() + grid.atomToCluster(atomOffset));
950 /* Sets all required atom parameter data in nbnxn_atomdata_t */
951 void nbnxn_atomdata_set(nbnxn_atomdata_t* nbat,
952 const Nbnxm::GridSet& gridSet,
953 ArrayRef<const int> atomTypes,
954 ArrayRef<const real> atomCharges,
955 ArrayRef<const int> atomInfo)
957 nbnxn_atomdata_t::Params& params = nbat->paramsDeprecated();
959 nbnxn_atomdata_set_atomtypes(¶ms, gridSet, atomTypes);
961 nbnxn_atomdata_set_charges(nbat, gridSet, atomCharges);
963 if (gridSet.haveFep())
965 nbnxn_atomdata_mask_fep(nbat, gridSet);
968 /* This must be done after masking types for FEP */
969 nbnxn_atomdata_set_ljcombparams(¶ms, nbat->XFormat, gridSet);
971 nbnxn_atomdata_set_energygroups(¶ms, gridSet, atomInfo);
974 /* Copies the shift vector array to nbnxn_atomdata_t */
975 void nbnxn_atomdata_copy_shiftvec(gmx_bool bDynamicBox, rvec* shift_vec, nbnxn_atomdata_t* nbat)
979 nbat->bDynamicBox = bDynamicBox;
980 for (i = 0; i < SHIFTS; i++)
982 copy_rvec(shift_vec[i], nbat->shift_vec[i]);
986 // This is slightly different from nbnxn_get_atom_range(...) at the end of the file
987 // TODO: Combine if possible
988 static void getAtomRanges(const Nbnxm::GridSet& gridSet,
989 const gmx::AtomLocality locality,
995 case gmx::AtomLocality::All:
997 *gridEnd = gridSet.grids().size();
999 case gmx::AtomLocality::Local:
1003 case gmx::AtomLocality::NonLocal:
1005 *gridEnd = gridSet.grids().size();
1007 case gmx::AtomLocality::Count:
1008 GMX_ASSERT(false, "Count is invalid locality specifier");
1013 /* Copies (and reorders) the coordinates to nbnxn_atomdata_t */
1014 void nbnxn_atomdata_copy_x_to_nbat_x(const Nbnxm::GridSet& gridSet,
1015 const gmx::AtomLocality locality,
1017 const rvec* coordinates,
1018 nbnxn_atomdata_t* nbat)
1023 getAtomRanges(gridSet, locality, &gridBegin, &gridEnd);
1027 nbat->natoms_local = gridSet.grids()[0].atomIndexEnd();
1030 const int nth = gmx_omp_nthreads_get(emntPairsearch);
1031 #pragma omp parallel for num_threads(nth) schedule(static)
1032 for (int th = 0; th < nth; th++)
1036 for (int g = gridBegin; g < gridEnd; g++)
1038 const Nbnxm::Grid& grid = gridSet.grids()[g];
1039 const int numCellsXY = grid.numColumns();
1041 const int cxy0 = (numCellsXY * th + nth - 1) / nth;
1042 const int cxy1 = (numCellsXY * (th + 1) + nth - 1) / nth;
1044 for (int cxy = cxy0; cxy < cxy1; cxy++)
1046 const int na = grid.numAtomsInColumn(cxy);
1047 const int ash = grid.firstAtomInColumn(cxy);
1050 if (g == 0 && fillLocal)
1052 na_fill = grid.paddedNumAtomsInColumn(cxy);
1056 /* We fill only the real particle locations.
1057 * We assume the filling entries at the end have been
1058 * properly set before during pair-list generation.
1062 copy_rvec_to_nbat_real(gridSet.atomIndices().data() + ash,
1072 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1076 /* Copies (and reorders) the coordinates to nbnxn_atomdata_t on the GPU*/
1077 void nbnxn_atomdata_x_to_nbat_x_gpu(const Nbnxm::GridSet& gridSet,
1078 const gmx::AtomLocality locality,
1081 DeviceBuffer<RVec> d_x,
1082 GpuEventSynchronizer* xReadyOnDevice)
1087 getAtomRanges(gridSet, locality, &gridBegin, &gridEnd);
1089 for (int g = gridBegin; g < gridEnd; g++)
1091 nbnxn_gpu_x_to_nbat_x(gridSet.grids()[g],
1092 fillLocal && g == 0,
1098 gridSet.numColumnsMax());
1102 static void nbnxn_atomdata_clear_reals(gmx::ArrayRef<real> dest, int i0, int i1)
1104 for (int i = i0; i < i1; i++)
1110 gmx_unused static void nbnxn_atomdata_reduce_reals(real* gmx_restrict dest,
1112 const real** gmx_restrict src,
1119 /* The destination buffer contains data, add to it */
1120 for (int i = i0; i < i1; i++)
1122 for (int s = 0; s < nsrc; s++)
1124 dest[i] += src[s][i];
1130 /* The destination buffer is unitialized, set it first */
1131 for (int i = i0; i < i1; i++)
1133 dest[i] = src[0][i];
1134 for (int s = 1; s < nsrc; s++)
1136 dest[i] += src[s][i];
1142 gmx_unused static void nbnxn_atomdata_reduce_reals_simd(real gmx_unused* gmx_restrict dest,
1143 gmx_bool gmx_unused bDestSet,
1144 const gmx_unused real** gmx_restrict src,
1145 int gmx_unused nsrc,
1150 /* The SIMD width here is actually independent of that in the kernels,
1151 * but we use the same width for simplicity (usually optimal anyhow).
1153 SimdReal dest_SSE, src_SSE;
1157 for (int i = i0; i < i1; i += GMX_SIMD_REAL_WIDTH)
1159 dest_SSE = load<SimdReal>(dest + i);
1160 for (int s = 0; s < nsrc; s++)
1162 src_SSE = load<SimdReal>(src[s] + i);
1163 dest_SSE = dest_SSE + src_SSE;
1165 store(dest + i, dest_SSE);
1170 for (int i = i0; i < i1; i += GMX_SIMD_REAL_WIDTH)
1172 dest_SSE = load<SimdReal>(src[0] + i);
1173 for (int s = 1; s < nsrc; s++)
1175 src_SSE = load<SimdReal>(src[s] + i);
1176 dest_SSE = dest_SSE + src_SSE;
1178 store(dest + i, dest_SSE);
1184 /* Add part of the force array(s) from nbnxn_atomdata_t to f
1186 * Note: Adding restrict to f makes this function 50% slower with gcc 7.3
1188 static void nbnxn_atomdata_add_nbat_f_to_f_part(const Nbnxm::GridSet& gridSet,
1189 const nbnxn_atomdata_t& nbat,
1190 const nbnxn_atomdata_output_t& out,
1195 gmx::ArrayRef<const int> cell = gridSet.cells();
1196 // Note: Using ArrayRef instead makes this code 25% slower with gcc 7.3
1197 const real* fnb = out.f.data();
1199 /* Loop over all columns and copy and fill */
1200 switch (nbat.FFormat)
1204 for (int a = a0; a < a1; a++)
1206 int i = cell[a] * nbat.fstride;
1209 f[a][YY] += fnb[i + 1];
1210 f[a][ZZ] += fnb[i + 2];
1214 for (int a = a0; a < a1; a++)
1216 int i = atom_to_x_index<c_packX4>(cell[a]);
1218 f[a][XX] += fnb[i + XX * c_packX4];
1219 f[a][YY] += fnb[i + YY * c_packX4];
1220 f[a][ZZ] += fnb[i + ZZ * c_packX4];
1224 for (int a = a0; a < a1; a++)
1226 int i = atom_to_x_index<c_packX8>(cell[a]);
1228 f[a][XX] += fnb[i + XX * c_packX8];
1229 f[a][YY] += fnb[i + YY * c_packX8];
1230 f[a][ZZ] += fnb[i + ZZ * c_packX8];
1233 default: gmx_incons("Unsupported nbnxn_atomdata_t format");
1237 static inline unsigned char reverse_bits(unsigned char b)
1239 /* http://graphics.stanford.edu/~seander/bithacks.html#ReverseByteWith64BitsDiv */
1240 return (b * 0x0202020202ULL & 0x010884422010ULL) % 1023;
1243 static void nbnxn_atomdata_add_nbat_f_to_f_reduce(nbnxn_atomdata_t* nbat, int nth)
1245 #pragma omp parallel for num_threads(nth) schedule(static)
1246 for (int th = 0; th < nth; th++)
1251 const real* fptr[NBNXN_BUFFERFLAG_MAX_THREADS];
1253 gmx::ArrayRef<const gmx_bitmask_t> flags = nbat->buffer_flags;
1255 /* Calculate the cell-block range for our thread */
1256 int b0 = (flags.size() * th) / nth;
1257 int b1 = (flags.size() * (th + 1)) / nth;
1259 for (int b = b0; b < b1; b++)
1261 int i0 = b * NBNXN_BUFFERFLAG_SIZE * nbat->fstride;
1262 int i1 = (b + 1) * NBNXN_BUFFERFLAG_SIZE * nbat->fstride;
1265 for (gmx::index out = 1; out < gmx::ssize(nbat->out); out++)
1267 if (bitmask_is_set(flags[b], out))
1269 fptr[nfptr++] = nbat->out[out].f.data();
1275 nbnxn_atomdata_reduce_reals_simd
1277 nbnxn_atomdata_reduce_reals
1279 (nbat->out[0].f.data(), bitmask_is_set(flags[b], 0), fptr, nfptr, i0, i1);
1281 else if (!bitmask_is_set(flags[b], 0))
1283 nbnxn_atomdata_clear_reals(nbat->out[0].f, i0, i1);
1287 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1292 /* Add the force array(s) from nbnxn_atomdata_t to f */
1293 void reduceForces(nbnxn_atomdata_t* nbat, const gmx::AtomLocality locality, const Nbnxm::GridSet& gridSet, rvec* f)
1298 nbnxn_get_atom_range(locality, gridSet, &a0, &na);
1302 /* The are no atoms for this reduction, avoid some overhead */
1306 int nth = gmx_omp_nthreads_get(emntNonbonded);
1308 if (nbat->out.size() > 1)
1310 if (locality != gmx::AtomLocality::All)
1312 gmx_incons("add_f_to_f called with nout>1 and locality!=eatAll");
1315 /* Reduce the force thread output buffers into buffer 0, before adding
1316 * them to the, differently ordered, "real" force buffer.
1318 nbnxn_atomdata_add_nbat_f_to_f_reduce(nbat, nth);
1320 #pragma omp parallel for num_threads(nth) schedule(static)
1321 for (int th = 0; th < nth; th++)
1325 nbnxn_atomdata_add_nbat_f_to_f_part(
1326 gridSet, *nbat, nbat->out[0], a0 + ((th + 0) * na) / nth, a0 + ((th + 1) * na) / nth, f);
1328 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1332 void nbnxn_atomdata_add_nbat_fshift_to_fshift(const nbnxn_atomdata_t& nbat, gmx::ArrayRef<gmx::RVec> fshift)
1334 gmx::ArrayRef<const nbnxn_atomdata_output_t> outputBuffers = nbat.out;
1336 for (int s = 0; s < SHIFTS; s++)
1340 for (const nbnxn_atomdata_output_t& out : outputBuffers)
1342 sum[XX] += out.fshift[s * DIM + XX];
1343 sum[YY] += out.fshift[s * DIM + YY];
1344 sum[ZZ] += out.fshift[s * DIM + ZZ];
1350 void nbnxn_get_atom_range(const gmx::AtomLocality atomLocality,
1351 const Nbnxm::GridSet& gridSet,
1356 switch (atomLocality)
1358 case gmx::AtomLocality::All:
1360 *nAtoms = gridSet.numRealAtomsTotal();
1362 case gmx::AtomLocality::Local:
1364 *nAtoms = gridSet.numRealAtomsLocal();
1366 case gmx::AtomLocality::NonLocal:
1367 *atomStart = gridSet.numRealAtomsLocal();
1368 *nAtoms = gridSet.numRealAtomsTotal() - gridSet.numRealAtomsLocal();
1370 case gmx::AtomLocality::Count:
1371 GMX_ASSERT(false, "Count is invalid locality specifier");