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.
47 #include "gromacs/math/functions.h"
48 #include "gromacs/math/utilities.h"
49 #include "gromacs/math/vec.h"
50 #include "gromacs/mdlib/gmx_omp_nthreads.h"
51 #include "gromacs/mdtypes/forcerec.h" // only for GET_CGINFO_*
52 #include "gromacs/mdtypes/md_enums.h"
53 #include "gromacs/nbnxm/nbnxm.h"
54 #include "gromacs/pbcutil/ishift.h"
55 #include "gromacs/simd/simd.h"
56 #include "gromacs/utility/exceptions.h"
57 #include "gromacs/utility/fatalerror.h"
58 #include "gromacs/utility/gmxomp.h"
59 #include "gromacs/utility/logger.h"
60 #include "gromacs/utility/strconvert.h"
61 #include "gromacs/utility/stringutil.h"
65 #include "nbnxm_geometry.h"
66 #include "nbnxm_gpu.h"
69 using namespace gmx; // TODO: Remove when this file is moved into gmx namespace
72 const char* enumValueToString(LJCombinationRule enumValue)
74 static constexpr gmx::EnumerationArray<LJCombinationRule, const char*> s_ljCombinationRuleNames = {
75 "Geometric", "Lorentz-Berthelot", "None"
77 return s_ljCombinationRuleNames[enumValue];
80 void nbnxn_atomdata_t::resizeCoordinateBuffer(int numAtoms)
84 x_.resize(numAtoms * xstride);
87 void nbnxn_atomdata_t::resizeForceBuffers()
89 /* Force buffers need padding up to a multiple of the buffer flag size */
90 const int paddedSize =
91 (numAtoms() + NBNXN_BUFFERFLAG_SIZE - 1) / NBNXN_BUFFERFLAG_SIZE * NBNXN_BUFFERFLAG_SIZE;
93 /* Should we let each thread allocate it's own data instead? */
94 for (nbnxn_atomdata_output_t& outBuffer : out)
96 outBuffer.f.resize(paddedSize * fstride);
100 /* Initializes an nbnxn_atomdata_output_t data structure */
101 nbnxn_atomdata_output_t::nbnxn_atomdata_output_t(Nbnxm::KernelType kernelType,
103 int simdEnergyBufferStride,
104 gmx::PinningPolicy pinningPolicy) :
105 f({}, { pinningPolicy }),
106 fshift({}, { pinningPolicy }),
107 Vvdw({}, { pinningPolicy }),
108 Vc({}, { pinningPolicy })
110 fshift.resize(gmx::c_numShiftVectors * DIM);
111 Vvdw.resize(numEnergyGroups * numEnergyGroups);
112 Vc.resize(numEnergyGroups * numEnergyGroups);
114 if (Nbnxm::kernelTypeIsSimd(kernelType))
116 int cj_size = Nbnxm::JClusterSizePerKernelType[kernelType];
118 numEnergyGroups * numEnergyGroups * simdEnergyBufferStride * (cj_size / 2) * cj_size;
119 VSvdw.resize(numElements);
120 VSc.resize(numElements);
124 static void copy_int_to_nbat_int(const int* a, int na, int na_round, const int* in, int fill, int* innb)
127 for (int i = 0; i < na; i++)
129 innb[j++] = in[a[i]];
131 /* Complete the partially filled last cell with fill */
132 for (int i = na; i < na_round; i++)
138 void copy_rvec_to_nbat_real(const int* a, int na, int na_round, const rvec* x, int nbatFormat, real* xnb, int a0)
140 /* We complete partially filled cells, can only be the last one in each
141 * column, with coordinates farAway. The actual coordinate value does
142 * not influence the results, since these filler particles do not interact.
143 * Clusters with normal atoms + fillers have a bounding box based only
144 * on the coordinates of the atoms. Clusters with only fillers have as
145 * the bounding box the coordinates of the first filler. Such clusters
146 * are not considered as i-entries, but they are considered as j-entries.
147 * So for performance it is better to have their bounding boxes far away,
148 * such that filler only clusters don't end up in the pair list.
150 const real farAway = -1000000;
152 if (nbatFormat == nbatXYZ)
155 int j = a0 * STRIDE_XYZ;
158 xnb[j++] = x[a[i]][XX];
159 xnb[j++] = x[a[i]][YY];
160 xnb[j++] = x[a[i]][ZZ];
162 /* Complete the partially filled last cell with farAway elements */
163 for (; i < na_round; i++)
170 else if (nbatFormat == nbatXYZQ)
173 int j = a0 * STRIDE_XYZQ;
176 xnb[j++] = x[a[i]][XX];
177 xnb[j++] = x[a[i]][YY];
178 xnb[j++] = x[a[i]][ZZ];
181 /* Complete the partially filled last cell with zeros */
182 for (; i < na_round; i++)
190 else if (nbatFormat == nbatX4)
193 int j = atom_to_x_index<c_packX4>(a0);
194 int c = a0 & (c_packX4 - 1);
197 xnb[j + XX * c_packX4] = x[a[i]][XX];
198 xnb[j + YY * c_packX4] = x[a[i]][YY];
199 xnb[j + ZZ * c_packX4] = x[a[i]][ZZ];
204 j += (DIM - 1) * c_packX4;
208 /* Complete the partially filled last cell with zeros */
209 for (; i < na_round; i++)
211 xnb[j + XX * c_packX4] = farAway;
212 xnb[j + YY * c_packX4] = farAway;
213 xnb[j + ZZ * c_packX4] = farAway;
218 j += (DIM - 1) * c_packX4;
223 else if (nbatFormat == nbatX8)
226 int j = atom_to_x_index<c_packX8>(a0);
227 int c = a0 & (c_packX8 - 1);
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;
258 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 subtract 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)
446 /* Initializes an nbnxn_atomdata_t::Params data structure */
447 static void nbnxn_atomdata_params_init(const gmx::MDLogger& mdlog,
448 nbnxn_atomdata_t::Params* params,
449 const Nbnxm::KernelType kernelType,
450 int enbnxninitcombrule,
452 ArrayRef<const real> nbfp,
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 const char* tolOverrideString = getenv("GMX_LJCOMB_TOL");
468 if (tolOverrideString != nullptr)
470 double tolOverride = std::strtod(tolOverrideString, nullptr);
473 bool bCombGeom = true;
476 /* Temporarily fill params->nbfp_comb with sigma and epsilon
477 * to check for the LB rule.
479 for (int i = 0; i < ntype; i++)
481 const real c6 = nbfp[(i * ntype + i) * 2] / 6.0;
482 const real c12 = nbfp[(i * ntype + i) * 2 + 1] / 12.0;
483 if (c6 > 0 && c12 > 0)
485 params->nbfp_comb[i * 2] = gmx::sixthroot(c12 / c6);
486 params->nbfp_comb[i * 2 + 1] = 0.25 * c6 * c6 / c12;
488 else if (c6 == 0 && c12 == 0)
490 params->nbfp_comb[i * 2] = 0;
491 params->nbfp_comb[i * 2 + 1] = 0;
495 /* Can not use LB rule with only dispersion or repulsion */
500 for (int i = 0; i < params->numTypes; i++)
502 for (int j = 0; j < params->numTypes; j++)
504 if (i < ntype && j < ntype)
506 /* fr->nbfp has been updated, so that array too now stores c6/c12 including
507 * the 6.0/12.0 prefactors to save 2 flops in the most common case (force-only).
509 real c6 = nbfp[(i * ntype + j) * 2];
510 real c12 = nbfp[(i * ntype + j) * 2 + 1];
512 params->nbfp[(i * params->numTypes + j) * 2] = c6;
513 params->nbfp[(i * params->numTypes + j) * 2 + 1] = c12;
515 /* Compare 6*C6 and 12*C12 for geometric cobination rule */
518 && gmx_within_tol(c6 * c6, nbfp[(i * ntype + i) * 2] * nbfp[(j * ntype + j) * 2], tol)
519 && gmx_within_tol(c12 * c12,
520 nbfp[(i * ntype + i) * 2 + 1] * nbfp[(j * ntype + j) * 2 + 1],
523 /* Compare C6 and C12 for Lorentz-Berthelot combination rule */
528 && ((c6 == 0 && c12 == 0
529 && (params->nbfp_comb[i * 2 + 1] == 0 || params->nbfp_comb[j * 2 + 1] == 0))
530 || (c6 > 0 && c12 > 0
531 && gmx_within_tol(gmx::sixthroot(c12 / c6),
532 0.5 * (params->nbfp_comb[i * 2] + params->nbfp_comb[j * 2]),
534 && gmx_within_tol(0.25 * c6 * c6 / c12,
535 std::sqrt(params->nbfp_comb[i * 2 + 1]
536 * params->nbfp_comb[j * 2 + 1]),
541 /* Add zero parameters for the additional dummy atom type */
542 params->nbfp[(i * params->numTypes + j) * 2] = 0;
543 params->nbfp[(i * params->numTypes + j) * 2 + 1] = 0;
550 "Combination rules: geometric %s Lorentz-Berthelot %s\n",
551 gmx::boolToString(bCombGeom),
552 gmx::boolToString(bCombLB));
555 const bool simple = Nbnxm::kernelTypeUsesSimplePairlist(kernelType);
557 switch (enbnxninitcombrule)
559 case enbnxninitcombruleDETECT:
560 /* We prefer the geometric combination rule,
561 * as that gives a slightly faster kernel than the LB rule.
565 params->ljCombinationRule = LJCombinationRule::Geometric;
569 params->ljCombinationRule = LJCombinationRule::LorentzBerthelot;
573 params->ljCombinationRule = LJCombinationRule::None;
575 params->nbfp_comb.clear();
580 if (params->ljCombinationRule == LJCombinationRule::None)
582 mesg = "Using full Lennard-Jones parameter combination matrix";
586 mesg = gmx::formatString("Using %s Lennard-Jones combination rule",
587 enumValueToString(params->ljCombinationRule));
589 GMX_LOG(mdlog.info).asParagraph().appendText(mesg);
592 case enbnxninitcombruleGEOM:
593 params->ljCombinationRule = LJCombinationRule::Geometric;
595 case enbnxninitcombruleLB:
596 params->ljCombinationRule = LJCombinationRule::LorentzBerthelot;
598 case enbnxninitcombruleNONE:
599 params->ljCombinationRule = LJCombinationRule::None;
601 params->nbfp_comb.clear();
603 default: gmx_incons("Unknown enbnxninitcombrule");
606 const bool bSIMD = Nbnxm::kernelTypeIsSimd(kernelType);
608 set_lj_parameter_data(params, bSIMD);
610 params->nenergrp = n_energygroups;
613 // We now check for energy groups already when starting mdrun
614 GMX_RELEASE_ASSERT(n_energygroups == 1, "GPU kernels do not support energy groups");
616 /* Temporary storage goes as #grp^3*simd_width^2/2, so limit to 64 */
617 if (params->nenergrp > 64)
619 gmx_fatal(FARGS, "With NxN kernels not more than 64 energy groups are supported\n");
621 params->neg_2log = 1;
622 while (params->nenergrp > (1 << params->neg_2log))
628 /* Initializes an nbnxn_atomdata_t data structure */
629 void nbnxn_atomdata_init(const gmx::MDLogger& mdlog,
630 nbnxn_atomdata_t* nbat,
631 const Nbnxm::KernelType kernelType,
632 int enbnxninitcombrule,
634 ArrayRef<const real> nbfp,
638 nbnxn_atomdata_params_init(
639 mdlog, &nbat->paramsDeprecated(), kernelType, enbnxninitcombrule, ntype, nbfp, n_energygroups);
641 const bool simple = Nbnxm::kernelTypeUsesSimplePairlist(kernelType);
642 const bool bSIMD = Nbnxm::kernelTypeIsSimd(kernelType);
648 int pack_x = std::max(c_nbnxnCpuIClusterSize, Nbnxm::JClusterSizePerKernelType[kernelType]);
651 case 4: nbat->XFormat = nbatX4; break;
652 case 8: nbat->XFormat = nbatX8; break;
653 default: gmx_incons("Unsupported packing width");
658 nbat->XFormat = nbatXYZ;
661 nbat->FFormat = nbat->XFormat;
665 nbat->XFormat = nbatXYZQ;
666 nbat->FFormat = nbatXYZ;
669 nbat->shift_vec.resize(gmx::c_numShiftVectors);
671 nbat->xstride = (nbat->XFormat == nbatXYZQ ? STRIDE_XYZQ : DIM);
672 nbat->fstride = (nbat->FFormat == nbatXYZQ ? STRIDE_XYZQ : DIM);
674 /* Initialize the output data structures */
675 for (int i = 0; i < nout; i++)
677 const auto& pinningPolicy = nbat->params().type.get_allocator().pinningPolicy();
678 nbat->out.emplace_back(
679 kernelType, nbat->params().nenergrp, 1 << nbat->params().neg_2log, pinningPolicy);
682 nbat->buffer_flags.clear();
685 template<int packSize>
686 static void copy_lj_to_nbat_lj_comb(gmx::ArrayRef<const real> ljparam_type, const int* type, int na, real* ljparam_at)
688 /* The LJ params follow the combination rule:
689 * copy the params for the type array to the atom array.
691 for (int is = 0; is < na; is += packSize)
693 for (int k = 0; k < packSize; k++)
696 ljparam_at[is * 2 + k] = ljparam_type[type[i] * 2];
697 ljparam_at[is * 2 + packSize + k] = ljparam_type[type[i] * 2 + 1];
702 /* Sets the atom type in nbnxn_atomdata_t */
703 static void nbnxn_atomdata_set_atomtypes(nbnxn_atomdata_t::Params* params,
704 const Nbnxm::GridSet& gridSet,
705 ArrayRef<const int> atomTypes)
707 params->type.resize(gridSet.numGridAtomsTotal());
709 for (const Nbnxm::Grid& grid : gridSet.grids())
711 /* Loop over all columns and copy and fill */
712 for (int i = 0; i < grid.numColumns(); i++)
714 const int numAtoms = grid.paddedNumAtomsInColumn(i);
715 const int atomOffset = grid.firstAtomInColumn(i);
717 copy_int_to_nbat_int(gridSet.atomIndices().data() + atomOffset,
718 grid.numAtomsInColumn(i),
721 params->numTypes - 1,
722 params->type.data() + atomOffset);
727 /* Sets the LJ combination rule parameters in nbnxn_atomdata_t */
728 static void nbnxn_atomdata_set_ljcombparams(nbnxn_atomdata_t::Params* params,
730 const Nbnxm::GridSet& gridSet)
732 params->lj_comb.resize(gridSet.numGridAtomsTotal() * 2);
734 if (params->ljCombinationRule != LJCombinationRule::None)
736 for (const Nbnxm::Grid& grid : gridSet.grids())
738 /* Loop over all columns and copy and fill */
739 for (int i = 0; i < grid.numColumns(); i++)
741 const int numAtoms = grid.paddedNumAtomsInColumn(i);
742 const int atomOffset = grid.firstAtomInColumn(i);
744 if (XFormat == nbatX4)
746 copy_lj_to_nbat_lj_comb<c_packX4>(params->nbfp_comb,
747 params->type.data() + atomOffset,
749 params->lj_comb.data() + atomOffset * 2);
751 else if (XFormat == nbatX8)
753 copy_lj_to_nbat_lj_comb<c_packX8>(params->nbfp_comb,
754 params->type.data() + atomOffset,
756 params->lj_comb.data() + atomOffset * 2);
758 else if (XFormat == nbatXYZQ)
760 copy_lj_to_nbat_lj_comb<1>(params->nbfp_comb,
761 params->type.data() + atomOffset,
763 params->lj_comb.data() + atomOffset * 2);
770 /* Sets the charges in nbnxn_atomdata_t *nbat */
771 static void nbnxn_atomdata_set_charges(nbnxn_atomdata_t* nbat,
772 const Nbnxm::GridSet& gridSet,
773 ArrayRef<const real> charges)
775 if (nbat->XFormat != nbatXYZQ)
777 nbat->paramsDeprecated().q.resize(nbat->numAtoms());
780 for (const Nbnxm::Grid& grid : gridSet.grids())
782 /* Loop over all columns and copy and fill */
783 for (int cxy = 0; cxy < grid.numColumns(); cxy++)
785 const int atomOffset = grid.firstAtomInColumn(cxy);
786 const int numAtoms = grid.numAtomsInColumn(cxy);
787 const int paddedNumAtoms = grid.paddedNumAtomsInColumn(cxy);
789 if (nbat->XFormat == nbatXYZQ)
791 real* q = nbat->x().data() + atomOffset * STRIDE_XYZQ + ZZ + 1;
792 for (int i = 0; i < numAtoms; i++)
794 *q = charges[gridSet.atomIndices()[atomOffset + i]];
797 /* Complete the partially filled last cell with zeros */
798 for (int i = numAtoms; i < paddedNumAtoms; i++)
806 real* q = nbat->paramsDeprecated().q.data() + atomOffset;
807 for (int i = 0; i < numAtoms; i++)
809 *q = charges[gridSet.atomIndices()[atomOffset + i]];
812 /* Complete the partially filled last cell with zeros */
813 for (int i = numAtoms; i < paddedNumAtoms; i++)
823 /* Set the charges of perturbed atoms in nbnxn_atomdata_t to 0.
824 * This is to automatically remove the RF/PME self term in the nbnxn kernels.
825 * Part of the zero interactions are still calculated in the normal kernels.
826 * All perturbed interactions are calculated in the free energy kernel,
827 * using the original charge and LJ data, not nbnxn_atomdata_t.
829 static void nbnxn_atomdata_mask_fep(nbnxn_atomdata_t* nbat, const Nbnxm::GridSet& gridSet)
831 nbnxn_atomdata_t::Params& params = nbat->paramsDeprecated();
833 const bool formatIsXYZQ = (nbat->XFormat == nbatXYZQ);
835 real* q = formatIsXYZQ ? (nbat->x().data() + ZZ + 1) : params.q.data();
836 int stride_q = formatIsXYZQ ? STRIDE_XYZQ : 1;
838 for (const Nbnxm::Grid& grid : gridSet.grids())
840 const int nsubc = (grid.geometry().isSimple) ? 1 : c_gpuNumClusterPerCell;
842 const int c_offset = grid.firstAtomInColumn(0);
844 /* Loop over all columns and copy and fill */
845 for (int c = 0; c < grid.numCells() * nsubc; c++)
847 /* Does this cluster contain perturbed particles? */
848 if (grid.clusterIsPerturbed(c))
850 const int numAtomsPerCluster = grid.geometry().numAtomsICluster;
851 for (int i = 0; i < numAtomsPerCluster; i++)
853 /* Is this a perturbed particle? */
854 if (grid.atomIsPerturbed(c, i))
856 int ind = c_offset + c * numAtomsPerCluster + i;
857 /* Set atom type and charge to non-interacting */
858 params.type[ind] = params.numTypes - 1;
859 q[ind * stride_q] = 0;
867 /* Copies the energy group indices to a reordered and packed array */
869 copy_egp_to_nbat_egps(const int* a, int na, int na_round, int na_c, int bit_shift, const int* in, int* innb)
872 for (; i < na; i += na_c)
874 /* Store na_c energy group numbers into one int */
876 for (int sa = 0; sa < na_c; sa++)
881 comb |= (GET_CGINFO_GID(in[at]) << (sa * bit_shift));
886 /* Complete the partially filled last cell with fill */
887 for (; i < na_round; i += na_c)
893 /* Set the energy group indices for atoms in nbnxn_atomdata_t */
894 static void nbnxn_atomdata_set_energygroups(nbnxn_atomdata_t::Params* params,
895 const Nbnxm::GridSet& gridSet,
896 ArrayRef<const int> atomInfo)
898 if (params->nenergrp == 1)
903 params->energrp.resize(gridSet.numGridAtomsTotal());
905 for (const Nbnxm::Grid& grid : gridSet.grids())
907 /* Loop over all columns and copy and fill */
908 for (int i = 0; i < grid.numColumns(); i++)
910 const int numAtoms = grid.paddedNumAtomsInColumn(i);
911 const int atomOffset = grid.firstAtomInColumn(i);
913 copy_egp_to_nbat_egps(gridSet.atomIndices().data() + atomOffset,
914 grid.numAtomsInColumn(i),
916 c_nbnxnCpuIClusterSize,
919 params->energrp.data() + grid.atomToCluster(atomOffset));
924 /* Sets all required atom parameter data in nbnxn_atomdata_t */
925 void nbnxn_atomdata_set(nbnxn_atomdata_t* nbat,
926 const Nbnxm::GridSet& gridSet,
927 ArrayRef<const int> atomTypes,
928 ArrayRef<const real> atomCharges,
929 ArrayRef<const int> atomInfo)
931 nbnxn_atomdata_t::Params& params = nbat->paramsDeprecated();
933 nbnxn_atomdata_set_atomtypes(¶ms, gridSet, atomTypes);
935 nbnxn_atomdata_set_charges(nbat, gridSet, atomCharges);
937 if (gridSet.haveFep())
939 nbnxn_atomdata_mask_fep(nbat, gridSet);
942 /* This must be done after masking types for FEP */
943 nbnxn_atomdata_set_ljcombparams(¶ms, nbat->XFormat, gridSet);
945 nbnxn_atomdata_set_energygroups(¶ms, gridSet, atomInfo);
948 /* Copies the shift vector array to nbnxn_atomdata_t */
949 void nbnxn_atomdata_copy_shiftvec(gmx_bool bDynamicBox, gmx::ArrayRef<gmx::RVec> shift_vec, nbnxn_atomdata_t* nbat)
951 nbat->bDynamicBox = bDynamicBox;
952 std::copy(shift_vec.begin(), shift_vec.end(), nbat->shift_vec.begin());
955 // This is slightly different from nbnxn_get_atom_range(...) at the end of the file
956 // TODO: Combine if possible
957 static void getAtomRanges(const Nbnxm::GridSet& gridSet,
958 const gmx::AtomLocality locality,
964 case gmx::AtomLocality::All:
966 *gridEnd = gridSet.grids().size();
968 case gmx::AtomLocality::Local:
972 case gmx::AtomLocality::NonLocal:
974 *gridEnd = gridSet.grids().size();
976 case gmx::AtomLocality::Count:
977 GMX_ASSERT(false, "Count is invalid locality specifier");
982 /* Copies (and reorders) the coordinates to nbnxn_atomdata_t */
983 void nbnxn_atomdata_copy_x_to_nbat_x(const Nbnxm::GridSet& gridSet,
984 const gmx::AtomLocality locality,
985 const rvec* coordinates,
986 nbnxn_atomdata_t* nbat)
991 getAtomRanges(gridSet, locality, &gridBegin, &gridEnd);
993 const int nth = gmx_omp_nthreads_get(ModuleMultiThread::Pairsearch);
994 #pragma omp parallel for num_threads(nth) schedule(static)
995 for (int th = 0; th < nth; th++)
999 for (int g = gridBegin; g < gridEnd; g++)
1001 const Nbnxm::Grid& grid = gridSet.grids()[g];
1002 const int numCellsXY = grid.numColumns();
1004 const int cxy0 = (numCellsXY * th + nth - 1) / nth;
1005 const int cxy1 = (numCellsXY * (th + 1) + nth - 1) / nth;
1007 for (int cxy = cxy0; cxy < cxy1; cxy++)
1009 const int na = grid.numAtomsInColumn(cxy);
1010 const int ash = grid.firstAtomInColumn(cxy);
1012 copy_rvec_to_nbat_real(gridSet.atomIndices().data() + ash,
1022 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1026 /* Copies (and reorders) the coordinates to nbnxn_atomdata_t on the GPU*/
1027 void nbnxn_atomdata_x_to_nbat_x_gpu(const Nbnxm::GridSet& gridSet,
1028 const gmx::AtomLocality locality,
1030 DeviceBuffer<RVec> d_x,
1031 GpuEventSynchronizer* xReadyOnDevice)
1033 GMX_ASSERT(xReadyOnDevice != nullptr, "Need a valid GpuEventSynchronizer object");
1037 getAtomRanges(gridSet, locality, &gridBegin, &gridEnd);
1039 for (int g = gridBegin; g < gridEnd; g++)
1041 nbnxn_gpu_x_to_nbat_x(gridSet.grids()[g],
1044 (g == gridBegin) ? xReadyOnDevice : nullptr, // Sync on first iteration only
1047 gridSet.numColumnsMax(),
1048 (g == gridEnd - 1));
1052 static void nbnxn_atomdata_clear_reals(gmx::ArrayRef<real> dest, int i0, int i1)
1054 for (int i = i0; i < i1; i++)
1060 gmx_unused static void nbnxn_atomdata_reduce_reals(real* gmx_restrict dest,
1062 const real** gmx_restrict src,
1069 /* The destination buffer contains data, add to it */
1070 for (int i = i0; i < i1; i++)
1072 for (int s = 0; s < nsrc; s++)
1074 dest[i] += src[s][i];
1080 /* The destination buffer is uninitialized, set it first */
1081 for (int i = i0; i < i1; i++)
1083 dest[i] = src[0][i];
1084 for (int s = 1; s < nsrc; s++)
1086 dest[i] += src[s][i];
1092 gmx_unused static void nbnxn_atomdata_reduce_reals_simd(real gmx_unused* gmx_restrict dest,
1093 gmx_bool gmx_unused bDestSet,
1094 const gmx_unused real** gmx_restrict src,
1095 int gmx_unused nsrc,
1100 /* The SIMD width here is actually independent of that in the kernels,
1101 * but we use the same width for simplicity (usually optimal anyhow).
1103 SimdReal dest_SSE, src_SSE;
1107 for (int i = i0; i < i1; i += GMX_SIMD_REAL_WIDTH)
1109 dest_SSE = load<SimdReal>(dest + i);
1110 for (int s = 0; s < nsrc; s++)
1112 src_SSE = load<SimdReal>(src[s] + i);
1113 dest_SSE = dest_SSE + src_SSE;
1115 store(dest + i, dest_SSE);
1120 for (int i = i0; i < i1; i += GMX_SIMD_REAL_WIDTH)
1122 dest_SSE = load<SimdReal>(src[0] + i);
1123 for (int s = 1; s < nsrc; s++)
1125 src_SSE = load<SimdReal>(src[s] + i);
1126 dest_SSE = dest_SSE + src_SSE;
1128 store(dest + i, dest_SSE);
1134 /* Add part of the force array(s) from nbnxn_atomdata_t to f
1136 * Note: Adding restrict to f makes this function 50% slower with gcc 7.3
1138 static void nbnxn_atomdata_add_nbat_f_to_f_part(const Nbnxm::GridSet& gridSet,
1139 const nbnxn_atomdata_t& nbat,
1140 const nbnxn_atomdata_output_t& out,
1145 gmx::ArrayRef<const int> cell = gridSet.cells();
1146 // Note: Using ArrayRef instead makes this code 25% slower with gcc 7.3
1147 const real* fnb = out.f.data();
1149 /* Loop over all columns and copy and fill */
1150 switch (nbat.FFormat)
1154 for (int a = a0; a < a1; a++)
1156 int i = cell[a] * nbat.fstride;
1159 f[a][YY] += fnb[i + 1];
1160 f[a][ZZ] += fnb[i + 2];
1164 for (int a = a0; a < a1; a++)
1166 int i = atom_to_x_index<c_packX4>(cell[a]);
1168 f[a][XX] += fnb[i + XX * c_packX4];
1169 f[a][YY] += fnb[i + YY * c_packX4];
1170 f[a][ZZ] += fnb[i + ZZ * c_packX4];
1174 for (int a = a0; a < a1; a++)
1176 int i = atom_to_x_index<c_packX8>(cell[a]);
1178 f[a][XX] += fnb[i + XX * c_packX8];
1179 f[a][YY] += fnb[i + YY * c_packX8];
1180 f[a][ZZ] += fnb[i + ZZ * c_packX8];
1183 default: gmx_incons("Unsupported nbnxn_atomdata_t format");
1187 static void nbnxn_atomdata_add_nbat_f_to_f_reduce(nbnxn_atomdata_t* nbat, int nth)
1189 #pragma omp parallel for num_threads(nth) schedule(static)
1190 for (int th = 0; th < nth; th++)
1194 const real* fptr[NBNXN_BUFFERFLAG_MAX_THREADS];
1196 gmx::ArrayRef<const gmx_bitmask_t> flags = nbat->buffer_flags;
1198 /* Calculate the cell-block range for our thread */
1199 const int b0 = (flags.size() * th) / nth;
1200 const int b1 = (flags.size() * (th + 1)) / nth;
1202 for (int b = b0; b < b1; b++)
1204 const int i0 = b * NBNXN_BUFFERFLAG_SIZE * nbat->fstride;
1205 const int i1 = (b + 1) * NBNXN_BUFFERFLAG_SIZE * nbat->fstride;
1208 for (gmx::index out = 1; out < gmx::ssize(nbat->out); out++)
1210 if (bitmask_is_set(flags[b], out))
1212 fptr[nfptr++] = nbat->out[out].f.data();
1218 nbnxn_atomdata_reduce_reals_simd
1220 nbnxn_atomdata_reduce_reals
1222 (nbat->out[0].f.data(), bitmask_is_set(flags[b], 0), fptr, nfptr, i0, i1);
1224 else if (!bitmask_is_set(flags[b], 0))
1226 nbnxn_atomdata_clear_reals(nbat->out[0].f, i0, i1);
1230 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1235 /* Add the force array(s) from nbnxn_atomdata_t to f */
1236 void reduceForces(nbnxn_atomdata_t* nbat, const gmx::AtomLocality locality, const Nbnxm::GridSet& gridSet, rvec* f)
1241 nbnxn_get_atom_range(locality, gridSet, &a0, &na);
1245 /* The are no atoms for this reduction, avoid some overhead */
1249 int nth = gmx_omp_nthreads_get(ModuleMultiThread::Nonbonded);
1251 if (nbat->out.size() > 1)
1253 if (locality != gmx::AtomLocality::All)
1255 gmx_incons("add_f_to_f called with nout>1 and locality!=eatAll");
1258 /* Reduce the force thread output buffers into buffer 0, before adding
1259 * them to the, differently ordered, "real" force buffer.
1261 nbnxn_atomdata_add_nbat_f_to_f_reduce(nbat, nth);
1263 #pragma omp parallel for num_threads(nth) schedule(static)
1264 for (int th = 0; th < nth; th++)
1268 nbnxn_atomdata_add_nbat_f_to_f_part(
1269 gridSet, *nbat, nbat->out[0], a0 + ((th + 0) * na) / nth, a0 + ((th + 1) * na) / nth, f);
1271 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1275 void nbnxn_atomdata_add_nbat_fshift_to_fshift(const nbnxn_atomdata_t& nbat, gmx::ArrayRef<gmx::RVec> fshift)
1277 gmx::ArrayRef<const nbnxn_atomdata_output_t> outputBuffers = nbat.out;
1279 for (int s = 0; s < gmx::c_numShiftVectors; s++)
1283 for (const nbnxn_atomdata_output_t& out : outputBuffers)
1285 sum[XX] += out.fshift[s * DIM + XX];
1286 sum[YY] += out.fshift[s * DIM + YY];
1287 sum[ZZ] += out.fshift[s * DIM + ZZ];
1293 void nbnxn_get_atom_range(const gmx::AtomLocality atomLocality,
1294 const Nbnxm::GridSet& gridSet,
1299 switch (atomLocality)
1301 case gmx::AtomLocality::All:
1303 *nAtoms = gridSet.numRealAtomsTotal();
1305 case gmx::AtomLocality::Local:
1307 *nAtoms = gridSet.numRealAtomsLocal();
1309 case gmx::AtomLocality::NonLocal:
1310 *atomStart = gridSet.numRealAtomsLocal();
1311 *nAtoms = gridSet.numRealAtomsTotal() - gridSet.numRealAtomsLocal();
1313 case gmx::AtomLocality::Count:
1314 GMX_ASSERT(false, "Count is invalid locality specifier");