2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015,2016,2017,2018,2019, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
39 * Implements the Grid class.
41 * \author Berk Hess <hess@kth.se>
42 * \ingroup module_nbnxm
54 #include "gromacs/domdec/domdec_struct.h"
55 #include "gromacs/math/utilities.h"
56 #include "gromacs/math/vec.h"
57 #include "gromacs/mdlib/gmx_omp_nthreads.h"
58 #include "gromacs/mdlib/updategroupscog.h"
59 #include "gromacs/mdtypes/forcerec.h" // only for GET_CGINFO_*
60 #include "gromacs/nbnxm/atomdata.h"
61 #include "gromacs/nbnxm/nbnxm.h"
62 #include "gromacs/nbnxm/nbnxm_geometry.h"
63 #include "gromacs/nbnxm/pairlist.h"
64 #include "gromacs/nbnxm/pairlistset.h"
65 #include "gromacs/simd/simd.h"
66 #include "gromacs/simd/vector_operations.h"
67 #include "gromacs/utility/exceptions.h"
68 #include "gromacs/utility/smalloc.h"
72 struct gmx_domdec_zones_t;
77 Grid::Geometry::Geometry(const PairlistType pairlistType) :
78 isSimple(pairlistType != PairlistType::Hierarchical8x8),
79 numAtomsICluster(IClusterSizePerListType[pairlistType]),
80 numAtomsJCluster(JClusterSizePerListType[pairlistType]),
81 numAtomsPerCell((isSimple ? 1 : c_gpuNumClusterPerCell)*numAtomsICluster),
82 numAtomsICluster2Log(get_2log(numAtomsICluster))
86 Grid::Grid(const PairlistType pairlistType) :
87 geometry_(pairlistType)
91 /*! \brief Returns the atom density (> 0) of a rectangular grid */
92 static real gridAtomDensity(int numAtoms,
93 const rvec lowerCorner,
94 const rvec upperCorner)
100 /* To avoid zero density we use a minimum of 1 atom */
104 rvec_sub(upperCorner, lowerCorner, size);
106 return numAtoms/(size[XX]*size[YY]*size[ZZ]);
109 void Grid::setDimensions(const nbnxn_search *nbs,
112 const rvec lowerCorner,
113 const rvec upperCorner,
115 const real maxAtomGroupRadius)
117 /* For the home zone we compute the density when not set (=-1) or when =0 */
118 if (ddZone == 0 && atomDensity <= 0)
120 atomDensity = gridAtomDensity(numAtoms, lowerCorner, upperCorner);
123 dimensions_.atomDensity = atomDensity;
124 dimensions_.maxAtomGroupRadius = maxAtomGroupRadius;
127 rvec_sub(upperCorner, lowerCorner, size);
129 if (numAtoms > geometry_.numAtomsPerCell)
131 GMX_ASSERT(atomDensity > 0, "With one or more atoms, the density should be positive");
133 /* target cell length */
136 if (geometry_.isSimple)
138 /* To minimize the zero interactions, we should make
139 * the largest of the i/j cell cubic.
141 int numAtomsInCell = std::max(geometry_.numAtomsICluster,
142 geometry_.numAtomsJCluster);
144 /* Approximately cubic cells */
145 real tlen = std::cbrt(numAtomsInCell/atomDensity);
151 /* Approximately cubic sub cells */
152 real tlen = std::cbrt(geometry_.numAtomsICluster/atomDensity);
153 tlen_x = tlen*c_gpuNumClusterPerCellX;
154 tlen_y = tlen*c_gpuNumClusterPerCellY;
156 /* We round ncx and ncy down, because we get less cell pairs
157 * in the nbsist when the fixed cell dimensions (x,y) are
158 * larger than the variable one (z) than the other way around.
160 dimensions_.numCells[XX] = std::max(1, static_cast<int>(size[XX]/tlen_x));
161 dimensions_.numCells[YY] = std::max(1, static_cast<int>(size[YY]/tlen_y));
165 dimensions_.numCells[XX] = 1;
166 dimensions_.numCells[YY] = 1;
169 for (int d = 0; d < DIM - 1; d++)
171 dimensions_.cellSize[d] = size[d]/dimensions_.numCells[d];
172 dimensions_.invCellSize[d] = 1/dimensions_.cellSize[d];
177 /* This is a non-home zone, add an extra row of cells
178 * for particles communicated for bonded interactions.
179 * These can be beyond the cut-off. It doesn't matter where
180 * they end up on the grid, but for performance it's better
181 * if they don't end up in cells that can be within cut-off range.
183 dimensions_.numCells[XX]++;
184 dimensions_.numCells[YY]++;
187 /* We need one additional cell entry for particles moved by DD */
188 cxy_na_.resize(numColumns() + 1);
189 cxy_ind_.resize(numColumns() + 2);
191 for (nbnxn_search_work_t &work : nbs->work)
193 work.cxy_na.resize(numColumns() + 1);
196 /* Worst case scenario of 1 atom in each last cell */
198 if (geometry_.numAtomsJCluster <= geometry_.numAtomsICluster)
200 maxNumCells = numAtoms/geometry_.numAtomsPerCell + numColumns();
204 maxNumCells = numAtoms/geometry_.numAtomsPerCell + numColumns()*geometry_.numAtomsJCluster/geometry_.numAtomsICluster;
207 if (!geometry_.isSimple)
209 numClusters_.resize(maxNumCells);
211 bbcz_.resize(maxNumCells*NNBSBB_D);
213 /* This resize also zeros the contents, this avoid possible
214 * floating exceptions in SIMD with the unused bb elements.
216 if (geometry_.isSimple)
218 bb_.resize(maxNumCells);
223 pbb_.resize(maxNumCells*c_gpuNumClusterPerCell/STRIDE_PBB*NNBSBB_XXXX);
225 bb_.resize(maxNumCells*c_gpuNumClusterPerCell);
229 if (geometry_.numAtomsJCluster == geometry_.numAtomsICluster)
235 GMX_ASSERT(geometry_.isSimple, "Only CPU lists should have different i/j cluster sizes");
237 bbjStorage_.resize(maxNumCells*geometry_.numAtomsICluster/geometry_.numAtomsJCluster);
241 flags_.resize(maxNumCells);
244 fep_.resize(maxNumCells*geometry_.numAtomsPerCell/geometry_.numAtomsICluster);
247 copy_rvec(lowerCorner, dimensions_.lowerCorner);
248 copy_rvec(upperCorner, dimensions_.upperCorner);
249 copy_rvec(size, dimensions_.gridSize);
252 /* We need to sort paricles in grid columns on z-coordinate.
253 * As particle are very often distributed homogeneously, we use a sorting
254 * algorithm similar to pigeonhole sort. We multiply the z-coordinate
255 * by a factor, cast to an int and try to store in that hole. If the hole
256 * is full, we move this or another particle. A second pass is needed to make
257 * contiguous elements. SORT_GRID_OVERSIZE is the ratio of holes to particles.
258 * 4 is the optimal value for homogeneous particle distribution and allows
259 * for an O(#particles) sort up till distributions were all particles are
260 * concentrated in 1/4 of the space. No NlogN fallback is implemented,
261 * as it can be expensive to detect imhomogeneous particle distributions.
263 /*! \brief Ratio of grid cells to atoms */
264 static constexpr int c_sortGridRatio = 4;
265 /*! \brief Maximum ratio of holes used, in the worst case all particles
266 * end up in the last hole and we need num. atoms extra holes at the end.
268 static constexpr int c_sortGridMaxSizeFactor = c_sortGridRatio + 1;
270 /*! \brief Sorts particle index a on coordinates x along dim.
272 * Backwards tells if we want decreasing iso increasing coordinates.
273 * h0 is the minimum of the coordinate range.
274 * invh is the 1/length of the sorting range.
275 * n_per_h (>=n) is the expected average number of particles per 1/invh
276 * sort is the sorting work array.
277 * sort should have a size of at least n_per_h*c_sortGridRatio + n,
278 * or easier, allocate at least n*c_sortGridMaxSizeFactor elements.
280 static void sort_atoms(int dim, gmx_bool Backwards,
281 int gmx_unused dd_zone,
282 bool gmx_unused relevantAtomsAreWithinGridBounds,
284 gmx::ArrayRef<const gmx::RVec> x,
285 real h0, real invh, int n_per_h,
286 gmx::ArrayRef<int> sort)
294 GMX_ASSERT(n <= n_per_h, "We require n <= n_per_h");
296 /* Transform the inverse range height into the inverse hole height */
297 invh *= n_per_h*c_sortGridRatio;
299 /* Set nsort to the maximum possible number of holes used.
300 * In worst case all n elements end up in the last bin.
302 int nsort = n_per_h*c_sortGridRatio + n;
304 /* Determine the index range used, so we can limit it for the second pass */
305 int zi_min = INT_MAX;
308 /* Sort the particles using a simple index sort */
309 for (int i = 0; i < n; i++)
311 /* The cast takes care of float-point rounding effects below zero.
312 * This code assumes particles are less than 1/SORT_GRID_OVERSIZE
313 * times the box height out of the box.
315 int zi = static_cast<int>((x[a[i]][dim] - h0)*invh);
318 /* As we can have rounding effect, we use > iso >= here */
319 if (relevantAtomsAreWithinGridBounds &&
320 (zi < 0 || (dd_zone == 0 && zi > n_per_h*c_sortGridRatio)))
322 gmx_fatal(FARGS, "(int)((x[%d][%c]=%f - %f)*%f) = %d, not in 0 - %d*%d\n",
323 a[i], 'x'+dim, x[a[i]][dim], h0, invh, zi,
324 n_per_h, c_sortGridRatio);
332 /* In a non-local domain, particles communcated for bonded interactions
333 * can be far beyond the grid size, which is set by the non-bonded
334 * cut-off distance. We sort such particles into the last cell.
336 if (zi > n_per_h*c_sortGridRatio)
338 zi = n_per_h*c_sortGridRatio;
341 /* Ideally this particle should go in sort cell zi,
342 * but that might already be in use,
343 * in that case find the first empty cell higher up
348 zi_min = std::min(zi_min, zi);
349 zi_max = std::max(zi_max, zi);
353 /* We have multiple atoms in the same sorting slot.
354 * Sort on real z for minimal bounding box size.
355 * There is an extra check for identical z to ensure
356 * well-defined output order, independent of input order
357 * to ensure binary reproducibility after restarts.
359 while (sort[zi] >= 0 && ( x[a[i]][dim] > x[sort[zi]][dim] ||
360 (x[a[i]][dim] == x[sort[zi]][dim] &&
368 /* Shift all elements by one slot until we find an empty slot */
371 while (sort[zim] >= 0)
379 zi_max = std::max(zi_max, zim);
382 zi_max = std::max(zi_max, zi);
389 for (int zi = 0; zi < nsort; zi++)
400 for (int zi = zi_max; zi >= zi_min; zi--)
411 gmx_incons("Lost particles while sorting");
416 //! Returns double up to one least significant float bit smaller than x
417 static double R2F_D(const float x)
419 return static_cast<float>(x >= 0 ? (1 - GMX_FLOAT_EPS)*x : (1 + GMX_FLOAT_EPS)*x);
421 //! Returns double up to one least significant float bit larger than x
422 static double R2F_U(const float x)
424 return static_cast<float>(x >= 0 ? (1 + GMX_FLOAT_EPS)*x : (1 - GMX_FLOAT_EPS)*x);
428 static float R2F_D(const float x)
433 static float R2F_U(const float x)
439 //! Computes the bounding box for na coordinates in order x,y,z, bb order xyz0
440 static void calc_bounding_box(int na, int stride, const real *x,
444 real xl, xh, yl, yh, zl, zh;
454 for (int j = 1; j < na; j++)
456 xl = std::min(xl, x[i+XX]);
457 xh = std::max(xh, x[i+XX]);
458 yl = std::min(yl, x[i+YY]);
459 yh = std::max(yh, x[i+YY]);
460 zl = std::min(zl, x[i+ZZ]);
461 zh = std::max(zh, x[i+ZZ]);
464 /* Note: possible double to float conversion here */
465 bb->lower.x = R2F_D(xl);
466 bb->lower.y = R2F_D(yl);
467 bb->lower.z = R2F_D(zl);
468 bb->upper.x = R2F_U(xh);
469 bb->upper.y = R2F_U(yh);
470 bb->upper.z = R2F_U(zh);
473 /*! \brief Computes the bounding box for na packed coordinates, bb order xyz0 */
474 static void calc_bounding_box_x_x4(int na, const real *x,
477 real xl, xh, yl, yh, zl, zh;
485 for (int j = 1; j < na; j++)
487 xl = std::min(xl, x[j+XX*c_packX4]);
488 xh = std::max(xh, x[j+XX*c_packX4]);
489 yl = std::min(yl, x[j+YY*c_packX4]);
490 yh = std::max(yh, x[j+YY*c_packX4]);
491 zl = std::min(zl, x[j+ZZ*c_packX4]);
492 zh = std::max(zh, x[j+ZZ*c_packX4]);
494 /* Note: possible double to float conversion here */
495 bb->lower.x = R2F_D(xl);
496 bb->lower.y = R2F_D(yl);
497 bb->lower.z = R2F_D(zl);
498 bb->upper.x = R2F_U(xh);
499 bb->upper.y = R2F_U(yh);
500 bb->upper.z = R2F_U(zh);
503 /*! \brief Computes the bounding box for na coordinates, bb order xyz0 */
504 static void calc_bounding_box_x_x8(int na, const real *x,
507 real xl, xh, yl, yh, zl, zh;
515 for (int j = 1; j < na; j++)
517 xl = std::min(xl, x[j+XX*c_packX8]);
518 xh = std::max(xh, x[j+XX*c_packX8]);
519 yl = std::min(yl, x[j+YY*c_packX8]);
520 yh = std::max(yh, x[j+YY*c_packX8]);
521 zl = std::min(zl, x[j+ZZ*c_packX8]);
522 zh = std::max(zh, x[j+ZZ*c_packX8]);
524 /* Note: possible double to float conversion here */
525 bb->lower.x = R2F_D(xl);
526 bb->lower.y = R2F_D(yl);
527 bb->lower.z = R2F_D(zl);
528 bb->upper.x = R2F_U(xh);
529 bb->upper.y = R2F_U(yh);
530 bb->upper.z = R2F_U(zh);
533 /*! \brief Computes the bounding box for na packed coordinates, bb order xyz0 */
534 gmx_unused static void calc_bounding_box_x_x4_halves(int na, const real *x,
538 // TODO: During SIMDv2 transition only some archs use namespace (remove when done)
541 calc_bounding_box_x_x4(std::min(na, 2), x, bbj);
545 calc_bounding_box_x_x4(std::min(na-2, 2), x+(c_packX4>>1), bbj+1);
549 /* Set the "empty" bounding box to the same as the first one,
550 * so we don't need to treat special cases in the rest of the code.
552 #if NBNXN_SEARCH_BB_SIMD4
553 store4(bbj[1].lower.ptr(), load4(bbj[0].lower.ptr()));
554 store4(bbj[1].upper.ptr(), load4(bbj[0].upper.ptr()));
560 #if NBNXN_SEARCH_BB_SIMD4
561 store4(bb->lower.ptr(),
562 min(load4(bbj[0].lower.ptr()), load4(bbj[1].lower.ptr())));
563 store4(bb->upper.ptr(),
564 max(load4(bbj[0].upper.ptr()), load4(bbj[1].upper.ptr())));
567 bb->lower = BoundingBox::Corner::min(bbj[0].lower,
569 bb->upper = BoundingBox::Corner::max(bbj[0].upper,
575 #if NBNXN_SEARCH_BB_SIMD4
577 /*! \brief Computes the bounding box for na coordinates in order xyz, bb order xxxxyyyyzzzz */
578 static void calc_bounding_box_xxxx(int na, int stride, const real *x, float *bb)
581 real xl, xh, yl, yh, zl, zh;
591 for (int j = 1; j < na; j++)
593 xl = std::min(xl, x[i+XX]);
594 xh = std::max(xh, x[i+XX]);
595 yl = std::min(yl, x[i+YY]);
596 yh = std::max(yh, x[i+YY]);
597 zl = std::min(zl, x[i+ZZ]);
598 zh = std::max(zh, x[i+ZZ]);
601 /* Note: possible double to float conversion here */
602 bb[0*STRIDE_PBB] = R2F_D(xl);
603 bb[1*STRIDE_PBB] = R2F_D(yl);
604 bb[2*STRIDE_PBB] = R2F_D(zl);
605 bb[3*STRIDE_PBB] = R2F_U(xh);
606 bb[4*STRIDE_PBB] = R2F_U(yh);
607 bb[5*STRIDE_PBB] = R2F_U(zh);
610 #endif /* NBNXN_SEARCH_BB_SIMD4 */
612 #if NBNXN_SEARCH_SIMD4_FLOAT_X_BB
614 /*! \brief Computes the bounding box for na coordinates in order xyz?, bb order xyz0 */
615 static void calc_bounding_box_simd4(int na, const float *x,
618 // TODO: During SIMDv2 transition only some archs use namespace (remove when done)
621 static_assert(sizeof(BoundingBox::Corner) == GMX_SIMD4_WIDTH*sizeof(float),
622 "The Corner struct should hold exactly 4 floats");
624 Simd4Float bb_0_S = load4(x);
625 Simd4Float bb_1_S = bb_0_S;
627 for (int i = 1; i < na; i++)
629 Simd4Float x_S = load4(x + i*GMX_SIMD4_WIDTH);
630 bb_0_S = min(bb_0_S, x_S);
631 bb_1_S = max(bb_1_S, x_S);
634 store4(bb->lower.ptr(), bb_0_S);
635 store4(bb->upper.ptr(), bb_1_S);
638 /*! \brief Computes the bounding box for na coordinates in order xyz?, bb order xxxxyyyyzzzz */
639 static void calc_bounding_box_xxxx_simd4(int na, const float *x,
640 BoundingBox *bb_work_aligned,
643 calc_bounding_box_simd4(na, x, bb_work_aligned);
645 bb[0*STRIDE_PBB] = bb_work_aligned->lower.x;
646 bb[1*STRIDE_PBB] = bb_work_aligned->lower.y;
647 bb[2*STRIDE_PBB] = bb_work_aligned->lower.z;
648 bb[3*STRIDE_PBB] = bb_work_aligned->upper.x;
649 bb[4*STRIDE_PBB] = bb_work_aligned->upper.y;
650 bb[5*STRIDE_PBB] = bb_work_aligned->upper.z;
653 #endif /* NBNXN_SEARCH_SIMD4_FLOAT_X_BB */
656 /*! \brief Combines pairs of consecutive bounding boxes */
657 static void combine_bounding_box_pairs(const Grid &grid,
658 gmx::ArrayRef<const BoundingBox> bb,
659 gmx::ArrayRef<BoundingBox> bbj)
661 // TODO: During SIMDv2 transition only some archs use namespace (remove when done)
664 for (int i = 0; i < grid.numColumns(); i++)
666 /* Starting bb in a column is expected to be 2-aligned */
667 const int sc2 = grid.firstCellInColumn(i) >> 1;
668 /* For odd numbers skip the last bb here */
669 const int nc2 = (grid.numAtomsInColumn(i) + 3) >> (2 + 1);
671 for (c2 = sc2; c2 < sc2 + nc2; c2++)
673 #if NBNXN_SEARCH_BB_SIMD4
674 Simd4Float min_S, max_S;
676 min_S = min(load4(bb[c2*2+0].lower.ptr()),
677 load4(bb[c2*2+1].lower.ptr()));
678 max_S = max(load4(bb[c2*2+0].upper.ptr()),
679 load4(bb[c2*2+1].upper.ptr()));
680 store4(bbj[c2].lower.ptr(), min_S);
681 store4(bbj[c2].upper.ptr(), max_S);
683 bbj[c2].lower = BoundingBox::Corner::min(bb[c2*2 + 0].lower,
685 bbj[c2].upper = BoundingBox::Corner::max(bb[c2*2 + 0].upper,
689 if (((grid.numAtomsInColumn(i) + 3) >> 2) & 1)
691 /* The bb count in this column is odd: duplicate the last bb */
692 bbj[c2].lower = bb[c2*2].lower;
693 bbj[c2].upper = bb[c2*2].upper;
699 /*! \brief Prints the average bb size, used for debug output */
700 static void print_bbsizes_simple(FILE *fp,
704 for (int c = 0; c < grid.numCells(); c++)
706 const BoundingBox &bb = grid.iBoundingBoxes()[c];
707 ba[XX] += bb.upper.x - bb.lower.x;
708 ba[YY] += bb.upper.y - bb.lower.y;
709 ba[ZZ] += bb.upper.z - bb.lower.z;
711 dsvmul(1.0/grid.numCells(), ba, ba);
713 const Grid::Dimensions &dims = grid.dimensions();
714 real avgCellSizeZ = (dims.atomDensity > 0 ? grid.geometry().numAtomsICluster/(dims.atomDensity*dims.cellSize[XX]*dims.cellSize[YY]) : 0.0);
716 fprintf(fp, "ns bb: grid %4.2f %4.2f %4.2f abs %4.2f %4.2f %4.2f rel %4.2f %4.2f %4.2f\n",
720 ba[XX], ba[YY], ba[ZZ],
721 ba[XX]*dims.invCellSize[XX],
722 ba[YY]*dims.invCellSize[YY],
723 dims.atomDensity > 0 ? ba[ZZ]/avgCellSizeZ : 0.0);
726 /*! \brief Prints the average bb size, used for debug output */
727 static void print_bbsizes_supersub(FILE *fp,
735 for (int c = 0; c < grid.numCells(); c++)
738 for (int s = 0; s < grid.numClustersPerCell()[c]; s += STRIDE_PBB)
740 int cs_w = (c*c_gpuNumClusterPerCell + s)/STRIDE_PBB;
741 for (int i = 0; i < STRIDE_PBB; i++)
743 for (int d = 0; d < DIM; d++)
746 grid.packedBoundingBoxes()[cs_w*NNBSBB_XXXX + (DIM + d)*STRIDE_PBB + i] -
747 grid.packedBoundingBoxes()[cs_w*NNBSBB_XXXX + d *STRIDE_PBB + i];
752 for (int s = 0; s < grid.numClustersPerCell()[c]; s++)
754 const BoundingBox &bb = grid.iBoundingBoxes()[c*c_gpuNumClusterPerCell + s];
755 ba[XX] += bb.upper.x - bb.lower.x;
756 ba[YY] += bb.upper.y - bb.lower.y;
757 ba[ZZ] += bb.upper.z - bb.lower.z;
760 ns += grid.numClustersPerCell()[c];
762 dsvmul(1.0/ns, ba, ba);
764 const Grid::Dimensions &dims = grid.dimensions();
765 const real avgClusterSizeZ =
766 (dims.atomDensity > 0 ? grid.geometry().numAtomsPerCell/(dims.atomDensity*dims.cellSize[XX]*dims.cellSize[YY]*c_gpuNumClusterPerCellZ) : 0.0);
768 fprintf(fp, "ns bb: grid %4.2f %4.2f %4.2f abs %4.2f %4.2f %4.2f rel %4.2f %4.2f %4.2f\n",
769 dims.cellSize[XX]/c_gpuNumClusterPerCellX,
770 dims.cellSize[YY]/c_gpuNumClusterPerCellY,
772 ba[XX], ba[YY], ba[ZZ],
773 ba[XX]*c_gpuNumClusterPerCellX*dims.invCellSize[XX],
774 ba[YY]*c_gpuNumClusterPerCellY*dims.invCellSize[YY],
775 dims.atomDensity > 0 ? ba[ZZ]/avgClusterSizeZ : 0.0);
778 /*!\brief Set non-bonded interaction flags for the current cluster.
780 * Sorts atoms on LJ coefficients: !=0 first, ==0 at the end.
782 static void sort_cluster_on_flag(int numAtomsInCluster,
786 gmx::ArrayRef<int> order,
789 constexpr int c_maxNumAtomsInCluster = 8;
790 int sort1[c_maxNumAtomsInCluster];
791 int sort2[c_maxNumAtomsInCluster];
793 GMX_ASSERT(numAtomsInCluster <= c_maxNumAtomsInCluster, "Need to increase c_maxNumAtomsInCluster to support larger clusters");
798 for (int s = atomStart; s < atomEnd; s += numAtomsInCluster)
800 /* Make lists for this (sub-)cell on atoms with and without LJ */
803 gmx_bool haveQ = FALSE;
805 for (int a = s; a < std::min(s + numAtomsInCluster, atomEnd); a++)
807 haveQ = haveQ || GET_CGINFO_HAS_Q(atinfo[order[a]]);
809 if (GET_CGINFO_HAS_VDW(atinfo[order[a]]))
811 sort1[n1++] = order[a];
816 sort2[n2++] = order[a];
820 /* If we don't have atoms with LJ, there's nothing to sort */
823 *flags |= NBNXN_CI_DO_LJ(subc);
825 if (2*n1 <= numAtomsInCluster)
827 /* Only sort when strictly necessary. Ordering particles
828 * Ordering particles can lead to less accurate summation
829 * due to rounding, both for LJ and Coulomb interactions.
831 if (2*(a_lj_max - s) >= numAtomsInCluster)
833 for (int i = 0; i < n1; i++)
835 order[atomStart + i] = sort1[i];
837 for (int j = 0; j < n2; j++)
839 order[atomStart + n1 + j] = sort2[j];
843 *flags |= NBNXN_CI_HALF_LJ(subc);
848 *flags |= NBNXN_CI_DO_COUL(subc);
854 /*! \brief Fill a pair search cell with atoms.
856 * Potentially sorts atoms and sets the interaction flags.
858 void Grid::fillCell(nbnxn_search *nbs,
859 nbnxn_atomdata_t *nbat,
863 gmx::ArrayRef<const gmx::RVec> x,
864 BoundingBox gmx_unused *bb_work_aligned)
866 const int numAtoms = atomEnd - atomStart;
868 if (geometry_.isSimple)
870 /* Note that non-local grids are already sorted.
871 * Then sort_cluster_on_flag will only set the flags and the sorting
872 * will not affect the atom order.
874 sort_cluster_on_flag(geometry_.numAtomsICluster, atomStart, atomEnd, atinfo, nbs->a,
875 flags_.data() + atomToCluster(atomStart) - cellOffset_);
880 /* Set the fep flag for perturbed atoms in this (sub-)cell */
882 /* The grid-local cluster/(sub-)cell index */
883 int cell = atomToCluster(atomStart) - cellOffset_*(geometry_.isSimple ? 1 : c_gpuNumClusterPerCell);
885 for (int at = atomStart; at < atomEnd; at++)
887 if (nbs->a[at] >= 0 && GET_CGINFO_FEP(atinfo[nbs->a[at]]))
889 fep_[cell] |= (1 << (at - atomStart));
894 /* Now we have sorted the atoms, set the cell indices */
895 for (int at = atomStart; at < atomEnd; at++)
897 nbs->cell[nbs->a[at]] = at;
900 copy_rvec_to_nbat_real(nbs->a.data() + atomStart, numAtoms, geometry_.numAtomsICluster,
901 as_rvec_array(x.data()),
902 nbat->XFormat, nbat->x().data(), atomStart);
904 if (nbat->XFormat == nbatX4)
906 /* Store the bounding boxes as xyz.xyz. */
907 size_t offset = atomToCluster(atomStart - cellOffset_*geometry_.numAtomsICluster);
908 BoundingBox *bb_ptr = bb_.data() + offset;
910 #if GMX_SIMD && GMX_SIMD_REAL_WIDTH == 2
911 if (2*geometry_.numAtomsJCluster == geometry_.numAtomsICluster)
913 calc_bounding_box_x_x4_halves(numAtoms, nbat->x().data() + atom_to_x_index<c_packX4>(atomStart), bb_ptr,
914 bbj_.data() + offset*2);
919 calc_bounding_box_x_x4(numAtoms, nbat->x().data() + atom_to_x_index<c_packX4>(atomStart), bb_ptr);
922 else if (nbat->XFormat == nbatX8)
924 /* Store the bounding boxes as xyz.xyz. */
925 size_t offset = atomToCluster(atomStart - cellOffset_*geometry_.numAtomsICluster);
926 BoundingBox *bb_ptr = bb_.data() + offset;
928 calc_bounding_box_x_x8(numAtoms, nbat->x().data() + atom_to_x_index<c_packX8>(atomStart), bb_ptr);
931 else if (!geometry_.isSimple)
933 /* Store the bounding boxes in a format convenient
934 * for SIMD4 calculations: xxxxyyyyzzzz...
938 ((atomStart - cellOffset_*geometry_.numAtomsPerCell) >> (geometry_.numAtomsICluster2Log + STRIDE_PBB_2LOG))*NNBSBB_XXXX +
939 (((atomStart - cellOffset_*geometry_.numAtomsPerCell) >> geometry_.numAtomsICluster2Log) & (STRIDE_PBB - 1));
941 #if NBNXN_SEARCH_SIMD4_FLOAT_X_BB
942 if (nbat->XFormat == nbatXYZQ)
944 calc_bounding_box_xxxx_simd4(numAtoms, nbat->x().data() + atomStart*nbat->xstride,
945 bb_work_aligned, pbb_ptr);
950 calc_bounding_box_xxxx(numAtoms, nbat->xstride, nbat->x().data() + atomStart*nbat->xstride,
955 fprintf(debug, "cell %4d bb %5.2f %5.2f %5.2f %5.2f %5.2f %5.2f\n",
956 atomToCluster(atomStart),
957 pbb_ptr[0*STRIDE_PBB], pbb_ptr[3*STRIDE_PBB],
958 pbb_ptr[1*STRIDE_PBB], pbb_ptr[4*STRIDE_PBB],
959 pbb_ptr[2*STRIDE_PBB], pbb_ptr[5*STRIDE_PBB]);
965 /* Store the bounding boxes as xyz.xyz. */
966 BoundingBox *bb_ptr = bb_.data() + atomToCluster(atomStart - cellOffset_*geometry_.numAtomsPerCell);
968 calc_bounding_box(numAtoms, nbat->xstride, nbat->x().data() + atomStart*nbat->xstride,
973 int bbo = atomToCluster(atomStart - cellOffset_*geometry_.numAtomsPerCell);
974 fprintf(debug, "cell %4d bb %5.2f %5.2f %5.2f %5.2f %5.2f %5.2f\n",
975 atomToCluster(atomStart),
986 void Grid::sortColumnsCpuGeometry(nbnxn_search *nbs,
988 int atomStart, int atomEnd,
990 gmx::ArrayRef<const gmx::RVec> x,
991 nbnxn_atomdata_t *nbat,
992 int cxy_start, int cxy_end,
993 gmx::ArrayRef<int> sort_work)
997 fprintf(debug, "cell_offset %d sorting columns %d - %d, atoms %d - %d\n",
998 cellOffset_, cxy_start, cxy_end, atomStart, atomEnd);
1001 const bool relevantAtomsAreWithinGridBounds = (dimensions_.maxAtomGroupRadius == 0);
1003 const int numAtomsPerCell = geometry_.numAtomsPerCell;
1005 /* Sort the atoms within each x,y column in 3 dimensions */
1006 for (int cxy = cxy_start; cxy < cxy_end; cxy++)
1008 const int numAtoms = numAtomsInColumn(cxy);
1009 const int numCellsZ = cxy_ind_[cxy + 1] - cxy_ind_[cxy];
1010 const int atomOffset = firstAtomInColumn(cxy);
1012 /* Sort the atoms within each x,y column on z coordinate */
1013 sort_atoms(ZZ, FALSE, dd_zone,
1014 relevantAtomsAreWithinGridBounds,
1015 nbs->a.data() + atomOffset, numAtoms, x,
1016 dimensions_.lowerCorner[ZZ],
1017 1.0/dimensions_.gridSize[ZZ], numCellsZ*numAtomsPerCell,
1020 /* Fill the ncz cells in this column */
1021 const int firstCell = firstCellInColumn(cxy);
1022 int cellFilled = firstCell;
1023 for (int cellZ = 0; cellZ < numCellsZ; cellZ++)
1025 const int cell = firstCell + cellZ;
1027 const int atomOffsetCell = atomOffset + cellZ*numAtomsPerCell;
1028 const int numAtomsCell = std::min(numAtomsPerCell, numAtoms - (atomOffsetCell - atomOffset));
1031 atomOffsetCell, atomOffsetCell + numAtomsCell, atinfo, x,
1034 /* This copy to bbcz is not really necessary.
1035 * But it allows to use the same grid search code
1036 * for the simple and supersub cell setups.
1038 if (numAtomsCell > 0)
1042 bbcz_[cell*NNBSBB_D ] = bb_[cellFilled].lower.z;
1043 bbcz_[cell*NNBSBB_D + 1] = bb_[cellFilled].upper.z;
1046 /* Set the unused atom indices to -1 */
1047 for (int ind = numAtoms; ind < numCellsZ*numAtomsPerCell; ind++)
1049 nbs->a[atomOffset + ind] = -1;
1054 /* Spatially sort the atoms within one grid column */
1055 void Grid::sortColumnsGpuGeometry(nbnxn_search *nbs,
1057 int atomStart, int atomEnd,
1059 gmx::ArrayRef<const gmx::RVec> x,
1060 nbnxn_atomdata_t *nbat,
1061 int cxy_start, int cxy_end,
1062 gmx::ArrayRef<int> sort_work)
1064 BoundingBox bb_work_array[2];
1065 BoundingBox *bb_work_aligned = reinterpret_cast<BoundingBox *>((reinterpret_cast<std::size_t>(bb_work_array + 1)) & (~(static_cast<std::size_t>(15))));
1069 fprintf(debug, "cell_offset %d sorting columns %d - %d, atoms %d - %d\n",
1070 cellOffset_, cxy_start, cxy_end, atomStart, atomEnd);
1073 const bool relevantAtomsAreWithinGridBounds = (dimensions_.maxAtomGroupRadius == 0);
1075 const int numAtomsPerCell = geometry_.numAtomsPerCell;
1077 const int subdiv_x = geometry_.numAtomsICluster;
1078 const int subdiv_y = c_gpuNumClusterPerCellX*subdiv_x;
1079 const int subdiv_z = c_gpuNumClusterPerCellY*subdiv_y;
1081 /* Sort the atoms within each x,y column in 3 dimensions.
1082 * Loop over all columns on the x/y grid.
1084 for (int cxy = cxy_start; cxy < cxy_end; cxy++)
1086 const int gridX = cxy/dimensions_.numCells[YY];
1087 const int gridY = cxy - gridX*dimensions_.numCells[YY];
1089 const int numAtomsInColumn = cxy_na_[cxy];
1090 const int numCellsInColumn = cxy_ind_[cxy + 1] - cxy_ind_[cxy];
1091 const int atomOffset = firstAtomInColumn(cxy);
1093 /* Sort the atoms within each x,y column on z coordinate */
1094 sort_atoms(ZZ, FALSE, dd_zone,
1095 relevantAtomsAreWithinGridBounds,
1096 nbs->a.data() + atomOffset, numAtomsInColumn, x,
1097 dimensions_.lowerCorner[ZZ],
1098 1.0/dimensions_.gridSize[ZZ], numCellsInColumn*numAtomsPerCell,
1101 /* This loop goes over the cells and clusters along z at once */
1102 for (int sub_z = 0; sub_z < numCellsInColumn*c_gpuNumClusterPerCellZ; sub_z++)
1104 const int atomOffsetZ = atomOffset + sub_z*subdiv_z;
1105 const int numAtomsZ = std::min(subdiv_z, numAtomsInColumn - (atomOffsetZ - atomOffset));
1107 /* We have already sorted on z */
1109 if (sub_z % c_gpuNumClusterPerCellZ == 0)
1111 cz = sub_z/c_gpuNumClusterPerCellZ;
1112 const int cell = cxy_ind_[cxy] + cz;
1114 /* The number of atoms in this cell/super-cluster */
1115 const int numAtoms = std::min(numAtomsPerCell, numAtomsInColumn - (atomOffsetZ - atomOffset));
1117 numClusters_[cell] = std::min(c_gpuNumClusterPerCell,
1118 (numAtoms + geometry_.numAtomsICluster - 1)/ geometry_.numAtomsICluster);
1120 /* Store the z-boundaries of the bounding box of the cell */
1121 bbcz_[cell*NNBSBB_D ] = x[nbs->a[atomOffsetZ]][ZZ];
1122 bbcz_[cell*NNBSBB_D+1] = x[nbs->a[atomOffsetZ + numAtoms - 1]][ZZ];
1125 if (c_gpuNumClusterPerCellY > 1)
1127 /* Sort the atoms along y */
1128 sort_atoms(YY, (sub_z & 1) != 0, dd_zone,
1129 relevantAtomsAreWithinGridBounds,
1130 nbs->a.data() + atomOffsetZ, numAtomsZ, x,
1131 dimensions_.lowerCorner[YY] + gridY*dimensions_.cellSize[YY],
1132 dimensions_.invCellSize[YY], subdiv_z,
1136 for (int sub_y = 0; sub_y < c_gpuNumClusterPerCellY; sub_y++)
1138 const int atomOffsetY = atomOffsetZ + sub_y*subdiv_y;
1139 const int numAtomsY = std::min(subdiv_y, numAtomsInColumn - (atomOffsetY - atomOffset));
1141 if (c_gpuNumClusterPerCellX > 1)
1143 /* Sort the atoms along x */
1144 sort_atoms(XX, ((cz*c_gpuNumClusterPerCellY + sub_y) & 1) != 0, dd_zone,
1145 relevantAtomsAreWithinGridBounds,
1146 nbs->a.data() + atomOffsetY, numAtomsY, x,
1147 dimensions_.lowerCorner[XX] + gridX*dimensions_.cellSize[XX],
1148 dimensions_.invCellSize[XX], subdiv_y,
1152 for (int sub_x = 0; sub_x < c_gpuNumClusterPerCellX; sub_x++)
1154 const int atomOffsetX = atomOffsetY + sub_x*subdiv_x;
1155 const int numAtomsX = std::min(subdiv_x, numAtomsInColumn - (atomOffsetX - atomOffset));
1158 atomOffsetX, atomOffsetX + numAtomsX, atinfo, x,
1164 /* Set the unused atom indices to -1 */
1165 for (int ind = numAtomsInColumn; ind < numCellsInColumn*numAtomsPerCell; ind++)
1167 nbs->a[atomOffset + ind] = -1;
1172 /*! \brief Sets the cell index in the cell array for atom \p atomIndex and increments the atom count for the grid column */
1173 static void setCellAndAtomCount(gmx::ArrayRef<int> cell,
1175 gmx::ArrayRef<int> cxy_na,
1178 cell[atomIndex] = cellIndex;
1179 cxy_na[cellIndex] += 1;
1182 /*! \brief Determine in which grid column atoms should go */
1183 static void calc_column_indices(const Grid::Dimensions &gridDims,
1184 const gmx::UpdateGroupsCog *updateGroupsCog,
1185 int atomStart, int atomEnd,
1186 gmx::ArrayRef<const gmx::RVec> x,
1187 int dd_zone, const int *move,
1188 int thread, int nthread,
1189 gmx::ArrayRef<int> cell,
1190 gmx::ArrayRef<int> cxy_na)
1192 const int numColumns = gridDims.numCells[XX]*gridDims.numCells[YY];
1194 /* We add one extra cell for particles which moved during DD */
1195 for (int i = 0; i < numColumns; i++)
1200 int taskAtomStart = atomStart + static_cast<int>((thread + 0)*(atomEnd - atomStart))/nthread;
1201 int taskAtomEnd = atomStart + static_cast<int>((thread + 1)*(atomEnd - atomStart))/nthread;
1206 for (int i = taskAtomStart; i < taskAtomEnd; i++)
1208 if (move == nullptr || move[i] >= 0)
1211 const gmx::RVec &coord = (updateGroupsCog ? updateGroupsCog->cogForAtom(i) : x[i]);
1213 /* We need to be careful with rounding,
1214 * particles might be a few bits outside the local zone.
1215 * The int cast takes care of the lower bound,
1216 * we will explicitly take care of the upper bound.
1218 int cx = static_cast<int>((coord[XX] - gridDims.lowerCorner[XX])*gridDims.invCellSize[XX]);
1219 int cy = static_cast<int>((coord[YY] - gridDims.lowerCorner[YY])*gridDims.invCellSize[YY]);
1222 if (cx < 0 || cx > gridDims.numCells[XX] ||
1223 cy < 0 || cy > gridDims.numCells[YY])
1226 "grid cell cx %d cy %d out of range (max %d %d)\n"
1227 "atom %f %f %f, grid->c0 %f %f",
1229 gridDims.numCells[XX],
1230 gridDims.numCells[YY],
1231 x[i][XX], x[i][YY], x[i][ZZ],
1232 gridDims.lowerCorner[XX],
1233 gridDims.lowerCorner[YY]);
1236 /* Take care of potential rouding issues */
1237 cx = std::min(cx, gridDims.numCells[XX] - 1);
1238 cy = std::min(cy, gridDims.numCells[YY] - 1);
1240 /* For the moment cell will contain only the, grid local,
1241 * x and y indices, not z.
1243 setCellAndAtomCount(cell, cx*gridDims.numCells[YY] + cy, cxy_na, i);
1247 /* Put this moved particle after the end of the grid,
1248 * so we can process it later without using conditionals.
1250 setCellAndAtomCount(cell, numColumns, cxy_na, i);
1257 for (int i = taskAtomStart; i < taskAtomEnd; i++)
1259 int cx = static_cast<int>((x[i][XX] - gridDims.lowerCorner[XX])*gridDims.invCellSize[XX]);
1260 int cy = static_cast<int>((x[i][YY] - gridDims.lowerCorner[YY])*gridDims.invCellSize[YY]);
1262 /* For non-home zones there could be particles outside
1263 * the non-bonded cut-off range, which have been communicated
1264 * for bonded interactions only. For the result it doesn't
1265 * matter where these end up on the grid. For performance
1266 * we put them in an extra row at the border.
1268 cx = std::max(cx, 0);
1269 cx = std::min(cx, gridDims.numCells[XX] - 1);
1270 cy = std::max(cy, 0);
1271 cy = std::min(cy, gridDims.numCells[YY] - 1);
1273 /* For the moment cell will contain only the, grid local,
1274 * x and y indices, not z.
1276 setCellAndAtomCount(cell, cx*gridDims.numCells[YY] + cy, cxy_na, i);
1281 /*! \brief Resizes grid and atom data which depend on the number of cells */
1282 static void resizeForNumberOfCells(const int numNbnxnAtoms,
1283 const int numAtomsMoved,
1285 nbnxn_atomdata_t *nbat)
1287 /* Note: nbs->cell was already resized before */
1289 /* To avoid conditionals we store the moved particles at the end of a,
1290 * make sure we have enough space.
1292 nbs->a.resize(numNbnxnAtoms + numAtomsMoved);
1294 /* Make space in nbat for storing the atom coordinates */
1295 nbat->resizeCoordinateBuffer(numNbnxnAtoms);
1298 void Grid::calcCellIndices(nbnxn_search *nbs,
1301 const gmx::UpdateGroupsCog *updateGroupsCog,
1305 gmx::ArrayRef<const gmx::RVec> x,
1308 nbnxn_atomdata_t *nbat)
1310 cellOffset_ = cellOffset;
1312 /* First compute all grid/column indices and store them in nbs->cell */
1313 nbs->cell.resize(atomEnd);
1315 const int nthread = gmx_omp_nthreads_get(emntPairsearch);
1317 #pragma omp parallel for num_threads(nthread) schedule(static)
1318 for (int thread = 0; thread < nthread; thread++)
1322 calc_column_indices(dimensions_,
1324 atomStart, atomEnd, x,
1325 ddZone, move, thread, nthread,
1326 nbs->cell, nbs->work[thread].cxy_na);
1328 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1331 const int numAtomsPerCell = geometry_.numAtomsPerCell;
1333 /* Make the cell index as a function of x and y */
1337 for (int i = 0; i < numColumns() + 1; i++)
1339 /* We set ncz_max at the beginning of the loop iso at the end
1340 * to skip i=grid->ncx*grid->numCells[YY] which are moved particles
1341 * that do not need to be ordered on the grid.
1347 int cxy_na_i = nbs->work[0].cxy_na[i];
1348 for (int thread = 1; thread < nthread; thread++)
1350 cxy_na_i += nbs->work[thread].cxy_na[i];
1352 ncz = (cxy_na_i + numAtomsPerCell - 1)/numAtomsPerCell;
1353 if (nbat->XFormat == nbatX8)
1355 /* Make the number of cell a multiple of 2 */
1356 ncz = (ncz + 1) & ~1;
1358 cxy_ind_[i+1] = cxy_ind_[i] + ncz;
1359 /* Clear cxy_na_, so we can reuse the array below */
1362 numCellsTotal_ = cxy_ind_[numColumns()] - cxy_ind_[0];
1364 resizeForNumberOfCells(atomIndexEnd(), numAtomsMoved, nbs, nbat);
1368 fprintf(debug, "ns na_sc %d na_c %d super-cells: %d x %d y %d z %.1f maxz %d\n",
1369 numAtomsPerCell, geometry_.numAtomsICluster, numCellsTotal_,
1370 dimensions_.numCells[XX], dimensions_.numCells[YY],
1371 numCellsTotal_/(static_cast<double>(numColumns())),
1376 for (int cy = 0; cy < dimensions_.numCells[YY]; cy++)
1378 for (int cx = 0; cx < dimensions_.numCells[XX]; cx++)
1380 fprintf(debug, " %2d", cxy_ind_[i + 1] - cxy_ind_[i]);
1383 fprintf(debug, "\n");
1388 /* Make sure the work array for sorting is large enough */
1389 const int worstCaseSortBufferSize = ncz_max*numAtomsPerCell*c_sortGridMaxSizeFactor;
1390 if (worstCaseSortBufferSize > gmx::index(nbs->work[0].sortBuffer.size()))
1392 for (nbnxn_search_work_t &work : nbs->work)
1394 /* Elements not in use should be -1 */
1395 work.sortBuffer.resize(worstCaseSortBufferSize, -1);
1399 /* Now we know the dimensions we can fill the grid.
1400 * This is the first, unsorted fill. We sort the columns after this.
1402 for (int i = atomStart; i < atomEnd; i++)
1404 /* At this point nbs->cell contains the local grid x,y indices */
1405 int cxy = nbs->cell[i];
1406 nbs->a[firstAtomInColumn(cxy) + cxy_na_[cxy]++] = i;
1411 /* Set the cell indices for the moved particles */
1412 int n0 = numCellsTotal_*numAtomsPerCell;
1413 int n1 = numCellsTotal_*numAtomsPerCell + cxy_na_[numColumns()];
1416 for (int i = n0; i < n1; i++)
1418 nbs->cell[nbs->a[i]] = i;
1423 /* Sort the super-cell columns along z into the sub-cells. */
1424 #pragma omp parallel for num_threads(nthread) schedule(static)
1425 for (int thread = 0; thread < nthread; thread++)
1429 int columnStart = ((thread + 0)*numColumns())/nthread;
1430 int columnEnd = ((thread + 1)*numColumns())/nthread;
1431 if (geometry_.isSimple)
1433 sortColumnsCpuGeometry(nbs, ddZone, atomStart, atomEnd, atinfo, x, nbat,
1434 columnStart, columnEnd,
1435 nbs->work[thread].sortBuffer);
1439 sortColumnsGpuGeometry(nbs, ddZone, atomStart, atomEnd, atinfo, x, nbat,
1440 columnStart, columnEnd,
1441 nbs->work[thread].sortBuffer);
1444 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1447 if (geometry_.isSimple && nbat->XFormat == nbatX8)
1449 combine_bounding_box_pairs(*this, bb_, bbj_);
1452 if (!geometry_.isSimple)
1454 numClustersTotal_ = 0;
1455 for (int i = 0; i < numCellsTotal_; i++)
1457 numClustersTotal_ += numClusters_[i];
1463 if (geometry_.isSimple)
1465 print_bbsizes_simple(debug, *this);
1469 fprintf(debug, "ns non-zero sub-cells: %d average atoms %.2f\n",
1471 (atomEnd - atomStart)/static_cast<double>(numClustersTotal_));
1473 print_bbsizes_supersub(debug, *this);
1478 } // namespace Nbnxm
1480 // TODO: Move the functions below to nbnxn.cpp
1482 /* Sets up a grid and puts the atoms on the grid.
1483 * This function only operates on one domain of the domain decompostion.
1484 * Note that without domain decomposition there is only one domain.
1486 void nbnxn_put_on_grid(nonbonded_verlet_t *nbv,
1489 const rvec lowerCorner,
1490 const rvec upperCorner,
1491 const gmx::UpdateGroupsCog *updateGroupsCog,
1496 gmx::ArrayRef<const gmx::RVec> x,
1500 nbnxn_search *nbs = nbv->nbs.get();
1501 Nbnxm::Grid &grid = nbs->grid[ddZone];
1503 nbs_cycle_start(&nbs->cc[enbsCCgrid]);
1512 const Nbnxm::Grid &previousGrid = nbs->grid[ddZone - 1];
1513 cellOffset = previousGrid.atomIndexEnd()/previousGrid.geometry().numAtomsPerCell;
1516 const int n = atomEnd - atomStart;
1518 real maxAtomGroupRadius;
1521 copy_mat(box, nbs->box);
1523 nbs->natoms_local = atomEnd - numAtomsMoved;
1524 /* We assume that nbnxn_put_on_grid is called first
1525 * for the local atoms (ddZone=0).
1527 nbs->natoms_nonlocal = atomEnd - numAtomsMoved;
1529 maxAtomGroupRadius = (updateGroupsCog ? updateGroupsCog->maxUpdateGroupRadius() : 0);
1533 fprintf(debug, "natoms_local = %5d atom_density = %5.1f\n",
1534 nbs->natoms_local, atomDensity);
1539 const Nbnxm::Grid::Dimensions &dimsGrid0 = nbs->grid[0].dimensions();
1540 atomDensity = dimsGrid0.atomDensity;
1541 maxAtomGroupRadius = dimsGrid0.maxAtomGroupRadius;
1543 nbs->natoms_nonlocal = std::max(nbs->natoms_nonlocal, atomEnd);
1546 /* We always use the home zone (grid[0]) for setting the cell size,
1547 * since determining densities for non-local zones is difficult.
1549 grid.setDimensions(nbs,
1550 ddZone, n - numAtomsMoved,
1551 lowerCorner, upperCorner,
1553 maxAtomGroupRadius);
1555 nbnxn_atomdata_t *nbat = nbv->nbat.get();
1557 grid.calcCellIndices(nbs, ddZone, cellOffset, updateGroupsCog, atomStart, atomEnd, atinfo, x, numAtomsMoved, move, nbat);
1561 nbat->natoms_local = nbat->numAtoms();
1563 if (ddZone == gmx::ssize(nbs->grid) - 1)
1565 /* We are done setting up all grids, we can resize the force buffers */
1566 nbat->resizeForceBuffers();
1569 nbs_cycle_stop(&nbs->cc[enbsCCgrid]);
1572 /* Calls nbnxn_put_on_grid for all non-local domains */
1573 void nbnxn_put_on_grid_nonlocal(nonbonded_verlet_t *nbv,
1574 const struct gmx_domdec_zones_t *zones,
1576 gmx::ArrayRef<const gmx::RVec> x)
1578 for (int zone = 1; zone < zones->n; zone++)
1581 for (int d = 0; d < DIM; d++)
1583 c0[d] = zones->size[zone].bb_x0[d];
1584 c1[d] = zones->size[zone].bb_x1[d];
1587 nbnxn_put_on_grid(nbv, nullptr,
1590 zones->cg_range[zone],
1591 zones->cg_range[zone+1],
1599 void nbnxn_get_ncells(const nbnxn_search *nbs, int *ncx, int *ncy)
1601 *ncx = nbs->grid[0].dimensions().numCells[XX];
1602 *ncy = nbs->grid[0].dimensions().numCells[YY];
1605 gmx::ArrayRef<const int> nbnxn_get_atomorder(const nbnxn_search *nbs)
1607 /* Return the atom order for the home cell (index 0) */
1608 const Nbnxm::Grid &grid = nbs->grid[0];
1610 const int numIndices = grid.atomIndexEnd() - grid.firstAtomInColumn(0);
1612 return gmx::constArrayRefFromArray(nbs->a.data(), numIndices);
1615 void nbnxn_set_atomorder(nbnxn_search *nbs)
1617 /* Set the atom order for the home cell (index 0) */
1618 const Nbnxm::Grid &grid = nbs->grid[0];
1621 for (int cxy = 0; cxy < grid.numColumns(); cxy++)
1623 const int numAtoms = grid.numAtomsInColumn(cxy);
1624 int cellIndex = grid.firstCellInColumn(cxy)*grid.geometry().numAtomsPerCell;
1625 for (int i = 0; i < numAtoms; i++)
1627 nbs->a[cellIndex] = atomIndex;
1628 nbs->cell[atomIndex] = cellIndex;
1635 gmx::ArrayRef<const int> nbnxn_get_gridindices(const nbnxn_search *nbs)