2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2009,2010,2011,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.
37 * Implements neighborhood searching for analysis (from nbsearch.h).
39 * High-level overview of the algorithm is at \ref page_analysisnbsearch.
42 * The grid implementation could still be optimized in several different ways:
43 * - A better heuristic for selecting the grid size or falling back to a
44 * simple all-pairs search.
45 * - A multi-level grid implementation could be used to be able to use small
46 * grids for short cutoffs with very inhomogeneous particle distributions
47 * without a memory cost.
49 * \author Teemu Murtola <teemu.murtola@gmail.com>
50 * \ingroup module_selection
62 #include "gromacs/math/functions.h"
63 #include "gromacs/math/vec.h"
64 #include "gromacs/pbcutil/pbc.h"
65 #include "gromacs/utility/arrayref.h"
66 #include "gromacs/utility/exceptions.h"
67 #include "gromacs/utility/gmxassert.h"
68 #include "gromacs/utility/listoflists.h"
69 #include "gromacs/utility/mutex.h"
70 #include "gromacs/utility/stringutil.h"
81 * Computes the bounding box for a set of positions.
83 * \param[in] posCount Number of positions in \p x.
84 * \param[in] x Positions to compute the bounding box for.
85 * \param[out] origin Origin of the bounding box.
86 * \param[out] size Size of the bounding box.
88 void computeBoundingBox(int posCount, const rvec x[], rvec origin, rvec size)
91 copy_rvec(x[0], origin);
92 copy_rvec(x[0], maxBound);
93 for (int i = 1; i < posCount; ++i)
95 for (int d = 0; d < DIM; ++d)
97 if (origin[d] > x[i][d])
101 if (maxBound[d] < x[i][d])
103 maxBound[d] = x[i][d];
107 rvec_sub(maxBound, origin, size);
115 /********************************************************************
116 * Implementation class declarations
119 class AnalysisNeighborhoodSearchImpl
122 typedef AnalysisNeighborhoodPairSearch::ImplPointer PairSearchImplPointer;
123 typedef std::vector<PairSearchImplPointer> PairSearchList;
124 typedef std::vector<std::vector<int>> CellList;
126 explicit AnalysisNeighborhoodSearchImpl(real cutoff);
127 ~AnalysisNeighborhoodSearchImpl();
130 * Initializes the search with a given box and reference positions.
132 * \param[in] mode Search mode to use.
133 * \param[in] bXY Whether to use 2D searching.
134 * \param[in] excls Exclusions.
135 * \param[in] pbc PBC information.
136 * \param[in] positions Set of reference positions.
138 void init(AnalysisNeighborhood::SearchMode mode,
140 const ListOfLists<int>* excls,
142 const AnalysisNeighborhoodPositions& positions);
143 PairSearchImplPointer getPairSearch();
145 real cutoffSquared() const { return cutoff2_; }
146 bool usesGridSearch() const { return bGrid_; }
150 * Determines a suitable grid size and sets up the cells.
152 * \param[in] box Box vectors (should not have zero vectors).
153 * \param[in] bSingleCell If `true`, the corresponding dimension will
154 * be forced to use a single cell.
155 * \param[in] posCount Number of positions that will be put on the
157 * \returns `false` if grid search is not suitable.
159 bool initGridCells(const matrix box, bool bSingleCell[DIM], int posCount);
161 * Sets ua a search grid for a given box.
163 * \param[in] pbc Information about the box.
164 * \param[in] posCount Number of positions in \p x.
165 * \param[in] x Reference positions that will be put on the grid.
166 * \param[in] bForce If `true`, grid searching will be used if at all
167 * possible, even if a simple search might give better performance.
168 * \returns `false` if grid search is not suitable.
170 bool initGrid(const t_pbc& pbc, int posCount, const rvec x[], bool bForce);
172 * Maps a point into a grid cell.
174 * \param[in] x Point to map.
175 * \param[out] cell Fractional cell coordinates of \p x on the grid.
176 * \param[out] xout Coordinates to use.
178 * \p xout will be within the rectangular unit cell in dimensions where
179 * the grid is periodic. For other dimensions, both \p xout and
180 * \p cell can be outside the grid/unit cell.
182 void mapPointToGridCell(const rvec x, rvec cell, rvec xout) const;
184 * Calculates linear index of a grid cell.
186 * \param[in] cell Cell indices (must be within the grid).
187 * \returns Linear index of \p cell.
189 int getGridCellIndex(const ivec cell) const;
191 * Calculates linear index of a grid cell from fractional coordinates.
193 * \param[in] cell Cell indices (must be within the grid).
194 * \returns Linear index of \p cell.
196 int getGridCellIndex(const rvec cell) const;
198 * Adds an index into a grid cell.
200 * \param[in] cell Fractional cell coordinates into which \p i should
202 * \param[in] i Index to add.
204 * \p cell should satisfy the conditions that \p mapPointToGridCell()
207 void addToGridCell(const rvec cell, int i);
209 * Initializes a cell pair loop for a dimension.
211 * \param[in] centerCell Fractional cell coordiates of the particle
212 * for which pairs are being searched.
213 * \param[in,out] cell Current/initial cell to loop over.
214 * \param[in,out] upperBound Last cell to loop over.
215 * \param[in] dim Dimension to initialize in this call.
217 * Initializes `cell[dim]` and `upperBound[dim]` for looping over
218 * neighbors of a particle at position given by \p centerCell.
219 * If 'dim != ZZ`, `cell[d]` (`d > dim`) set the plane/row of cells
220 * for which the loop is initialized. The loop should then go from
221 * `cell[dim]` until `upperBound[dim]`, inclusive.
222 * `cell[d]` with `d < dim` or `upperBound[d]` with `d != dim` are not
223 * modified by this function.
225 * `cell` and `upperBound` may be outside the grid for periodic
226 * dimensions and need to be shifted separately: to simplify the
227 * looping, the range is always (roughly) symmetric around the value in
230 void initCellRange(const rvec centerCell, ivec cell, ivec upperBound, int dim) const;
232 * Computes the extent of the cutoff sphere on a particular cell edge.
234 * \param[in] centerCell Fractional cell coordiates of the particle
235 * for which pairs are being searched.
236 * \param[in] cell Current cell (for dimensions `>dim`).
237 * \param[in] dim Dimension to compute in this call.
238 * \returns Fractional extent of the cutoff sphere when looping
239 * over cells in dimension `dim`, for `cell[d]` (`d > dim`).
241 * Input parameters are as for initCellRange(), except that if `cell`
242 * is over a periodic boundary from `centerCell`, triclinic shifts
243 * should have been applied to `centerCell` X/Y components.
245 real computeCutoffExtent(RVec centerCell, const ivec cell, int dim) const;
247 * Advances cell pair loop to the next cell.
249 * \param[in] centerCell Fractional cell coordiates of the particle
250 * for which pairs are being searched.
251 * \param[in,out] cell Current (in)/next (out) cell in the loop.
252 * \param[in,out] upperBound Last cell in the loop for each dimension.
254 bool nextCell(const rvec centerCell, ivec cell, ivec upperBound) const;
256 * Calculates the index and shift of a grid cell during looping.
258 * \param[in] cell Unshifted cell index.
259 * \param[out] shift Shift to apply to get the periodic distance
260 * for distances between the cells.
261 * \returns Grid cell index corresponding to `cell`.
263 int shiftCell(const ivec cell, rvec shift) const;
265 //! Whether to try grid searching.
269 //! The cutoff squared.
271 //! Whether to do searching in XY plane only.
274 //! Number of reference points for the current frame.
276 //! Reference point positions.
278 //! Reference position exclusion IDs.
279 const int* refExclusionIds_;
280 //! Reference position indices (NULL if no indices).
281 const int* refIndices_;
283 const ListOfLists<int>* excls_;
287 //! Whether grid searching is actually used for the current positions.
289 //! false if the box is rectangular.
291 //! Whether the grid is periodic in a dimension.
293 //! Array for storing in-unit-cell reference positions.
294 std::vector<RVec> xrefAlloc_;
295 //! Origin of the grid (zero for periodic dimensions).
297 //! Size of a single grid cell.
299 //! Inverse of \p cellSize_. Zero for dimensions where grid is not used.
302 * Shift in cell coordinates (for triclinic boxes) in X when crossing
303 * the Z periodic boundary.
307 * Shift in cell coordinates (for triclinic boxes) in Y when crossing
308 * the Z periodic boundary.
312 * Shift in cell coordinates (for triclinic boxes) in X when crossing
313 * the Y periodic boundary.
316 //! Number of cells along each dimension.
318 //! Data structure to hold the grid cell contents.
321 Mutex createPairSearchMutex_;
322 PairSearchList pairSearchList_;
324 friend class AnalysisNeighborhoodPairSearchImpl;
326 GMX_DISALLOW_COPY_AND_ASSIGN(AnalysisNeighborhoodSearchImpl);
329 class AnalysisNeighborhoodPairSearchImpl
332 explicit AnalysisNeighborhoodPairSearchImpl(const AnalysisNeighborhoodSearchImpl& search) :
335 selfSearchMode_ = false;
337 testPositions_ = nullptr;
338 testExclusionIds_ = nullptr;
339 testIndices_ = nullptr;
341 clear_rvec(testcell_);
342 clear_ivec(currCell_);
343 clear_ivec(cellBound_);
347 //! Initializes a search to find reference positions neighboring \p x.
348 void startSearch(const AnalysisNeighborhoodPositions& positions);
349 //! Initializes a search to find reference position pairs.
350 void startSelfSearch();
351 //! Searches for the next neighbor.
352 template<class Action>
353 bool searchNext(Action action);
354 //! Initializes a pair representing the pair found by searchNext().
355 void initFoundPair(AnalysisNeighborhoodPair* pair) const;
356 //! Advances to the next test position, skipping any remaining pairs.
357 void nextTestPosition();
360 //! Clears the loop indices.
361 void reset(int testIndex);
362 //! Checks whether a reference positiong should be excluded.
363 bool isExcluded(int j);
365 //! Parent search object.
366 const AnalysisNeighborhoodSearchImpl& search_;
367 //! Whether we are searching for ref-ref pairs.
368 bool selfSearchMode_;
369 //! Number of test positions.
371 //! Reference to the test positions.
372 const rvec* testPositions_;
373 //! Reference to the test exclusion indices.
374 const int* testExclusionIds_;
375 //! Reference to the test position indices.
376 const int* testIndices_;
377 //! Exclusions for current test particle.
378 ArrayRef<const int> excl_;
379 //! Index of the currently active test position in \p testPositions_.
381 //! Stores test position during a pair loop.
383 //! Stores the previous returned position during a pair loop.
385 //! Stores the pair distance corresponding to previ_.
387 //! Stores the shortest distance vector corresponding to previ_.
389 //! Stores the current exclusion index during loops.
391 //! Stores the fractional test particle cell location during loops.
393 //! Stores the cell index corresponding to testcell_.
395 //! Stores the current cell during pair loops.
397 //! Stores the current loop upper bounds for each dimension during pair loops.
399 //! Stores the index within the current cell during pair loops.
402 GMX_DISALLOW_COPY_AND_ASSIGN(AnalysisNeighborhoodPairSearchImpl);
405 /********************************************************************
406 * AnalysisNeighborhoodSearchImpl
409 AnalysisNeighborhoodSearchImpl::AnalysisNeighborhoodSearchImpl(real cutoff)
415 cutoff_ = cutoff2_ = GMX_REAL_MAX;
420 cutoff2_ = gmx::square(cutoff_);
425 refExclusionIds_ = nullptr;
426 refIndices_ = nullptr;
427 std::memset(&pbc_, 0, sizeof(pbc_));
431 bGridPBC_[XX] = true;
432 bGridPBC_[YY] = true;
433 bGridPBC_[ZZ] = true;
435 clear_rvec(gridOrigin_);
436 clear_rvec(cellSize_);
437 clear_rvec(invCellSize_);
438 clear_ivec(ncelldim_);
441 AnalysisNeighborhoodSearchImpl::~AnalysisNeighborhoodSearchImpl()
443 PairSearchList::const_iterator i;
444 for (i = pairSearchList_.begin(); i != pairSearchList_.end(); ++i)
446 GMX_RELEASE_ASSERT(i->unique(), "Dangling AnalysisNeighborhoodPairSearch reference");
450 AnalysisNeighborhoodSearchImpl::PairSearchImplPointer AnalysisNeighborhoodSearchImpl::getPairSearch()
452 lock_guard<Mutex> lock(createPairSearchMutex_);
453 // TODO: Consider whether this needs to/can be faster, e.g., by keeping a
454 // separate pool of unused search objects.
455 PairSearchList::const_iterator i;
456 for (i = pairSearchList_.begin(); i != pairSearchList_.end(); ++i)
463 PairSearchImplPointer pairSearch(new AnalysisNeighborhoodPairSearchImpl(*this));
464 pairSearchList_.push_back(pairSearch);
468 bool AnalysisNeighborhoodSearchImpl::initGridCells(const matrix box, bool bSingleCell[DIM], int posCount)
470 // Determine the size of cubes where there are on average 10 positions.
471 // The loop takes care of cases where some of the box edges are shorter
472 // than the desired cube size; in such cases, a single grid cell is
473 // used in these dimensions, and the cube size is determined only from the
474 // larger box vectors. Such boxes should be rare, but the bounding box
475 // approach can result in very flat boxes with certain types of selections
476 // (e.g., for interfacial systems or for small number of atoms).
477 real targetsize = 0.0;
478 int prevDimCount = 4;
483 for (int dd = 0; dd < DIM; ++dd)
485 const real boxSize = box[dd][dd];
486 if (boxSize < targetsize)
488 bSingleCell[dd] = true;
491 // TODO: Consider if a fallback would be possible/better.
504 if (dimCount == 0 || dimCount == prevDimCount)
508 targetsize = pow(volume * 10 / posCount, static_cast<real>(1. / dimCount));
509 prevDimCount = dimCount;
512 int totalCellCount = 1;
513 for (int dd = 0; dd < DIM; ++dd)
522 cellCount = std::max(1, static_cast<int>(box[dd][dd] / targetsize));
523 // TODO: If the cell count is one or two, it could be better to
524 // just fall back to bSingleCell[dd] = true.
525 if (bGridPBC_[dd] && cellCount < 3)
530 totalCellCount *= cellCount;
531 ncelldim_[dd] = cellCount;
533 if (totalCellCount <= 3)
537 // Never decrease the size of the cell vector to avoid reallocating
538 // memory for the nested vectors. The actual size of the vector is not
539 // used outside this function.
540 if (cells_.size() < static_cast<size_t>(totalCellCount))
542 cells_.resize(totalCellCount);
544 for (int ci = 0; ci < totalCellCount; ++ci)
551 bool AnalysisNeighborhoodSearchImpl::initGrid(const t_pbc& pbc, int posCount, const rvec x[], bool bForce)
558 // TODO: Use this again (can be useful when tuning initGridCells()),
559 // or remove throughout.
560 GMX_UNUSED_VALUE(bForce);
565 bGridPBC_[XX] = false;
566 bGridPBC_[YY] = false;
567 bGridPBC_[ZZ] = false;
570 bGridPBC_[XX] = true;
571 bGridPBC_[YY] = true;
572 bGridPBC_[ZZ] = false;
575 bGridPBC_[XX] = true;
576 bGridPBC_[YY] = true;
577 bGridPBC_[ZZ] = true;
580 // Grid searching not supported for now with screw.
584 bool bSingleCell[DIM] = { false, false, bXY_ };
586 copy_mat(pbc.box, box);
587 // TODO: In principle, we could use the bounding box for periodic
588 // dimensions as well if the bounding box is sufficiently far from the box
590 rvec origin, boundingBoxSize;
591 computeBoundingBox(posCount, x, origin, boundingBoxSize);
592 clear_rvec(gridOrigin_);
593 for (int dd = 0; dd < DIM; ++dd)
595 if (!bGridPBC_[dd] && !bSingleCell[dd])
597 gridOrigin_[dd] = origin[dd];
599 box[dd][dd] = boundingBoxSize[dd];
601 // TODO: In case the zero vector comes from the bounding box, this does
602 // not lead to a very efficient grid search, but that should be rare.
603 if (box[dd][dd] <= 0.0)
605 GMX_ASSERT(!bGridPBC_[dd], "Periodic box vector is zero");
606 bSingleCell[dd] = true;
612 if (!initGridCells(box, bSingleCell, posCount))
617 bTric_ = TRICLINIC(pbc.box);
618 for (int dd = 0; dd < DIM; ++dd)
620 cellSize_[dd] = box[dd][dd] / ncelldim_[dd];
623 invCellSize_[dd] = 0.0;
627 invCellSize_[dd] = 1.0 / cellSize_[dd];
628 // TODO: It could be better to avoid this when determining the cell
629 // size, but this can still remain here as a fallback to avoid
630 // incorrect results.
631 if (std::ceil(2 * cutoff_ * invCellSize_[dd]) >= ncelldim_[dd])
633 // Cutoff is too close to half the box size for grid searching
634 // (it is not possible to find a single shift for every pair of
642 cellShiftZY_ = box[ZZ][YY] * invCellSize_[YY];
643 cellShiftZX_ = box[ZZ][XX] * invCellSize_[XX];
644 cellShiftYX_ = box[YY][XX] * invCellSize_[XX];
649 void AnalysisNeighborhoodSearchImpl::mapPointToGridCell(const rvec x, rvec cell, rvec xout) const
652 rvec_sub(x, gridOrigin_, xtmp);
653 // The reverse order is necessary for triclinic cells: shifting in Z may
654 // modify also X and Y, and shifting in Y may modify X, so the mapping to
655 // a rectangular grid needs to be done in this order.
656 for (int dd = DIM - 1; dd >= 0; --dd)
658 real cellIndex = xtmp[dd] * invCellSize_[dd];
661 const real cellCount = ncelldim_[dd];
662 while (cellIndex < 0)
664 cellIndex += cellCount;
665 rvec_inc(xtmp, pbc_.box[dd]);
667 while (cellIndex >= cellCount)
669 cellIndex -= cellCount;
670 rvec_dec(xtmp, pbc_.box[dd]);
673 cell[dd] = cellIndex;
675 copy_rvec(xtmp, xout);
678 int AnalysisNeighborhoodSearchImpl::getGridCellIndex(const ivec cell) const
680 GMX_ASSERT(cell[XX] >= 0 && cell[XX] < ncelldim_[XX], "Grid cell X index out of range");
681 GMX_ASSERT(cell[YY] >= 0 && cell[YY] < ncelldim_[YY], "Grid cell Y index out of range");
682 GMX_ASSERT(cell[ZZ] >= 0 && cell[ZZ] < ncelldim_[ZZ], "Grid cell Z index out of range");
683 return cell[XX] + cell[YY] * ncelldim_[XX] + cell[ZZ] * ncelldim_[XX] * ncelldim_[YY];
686 int AnalysisNeighborhoodSearchImpl::getGridCellIndex(const rvec cell) const
689 for (int dd = 0; dd < DIM; ++dd)
691 int cellIndex = static_cast<int>(floor(cell[dd]));
694 const int cellCount = ncelldim_[dd];
699 else if (cellIndex >= cellCount)
701 cellIndex = cellCount - 1;
704 icell[dd] = cellIndex;
706 return getGridCellIndex(icell);
709 void AnalysisNeighborhoodSearchImpl::addToGridCell(const rvec cell, int i)
711 const int ci = getGridCellIndex(cell);
712 cells_[ci].push_back(i);
715 void AnalysisNeighborhoodSearchImpl::initCellRange(const rvec centerCell, ivec currCell, ivec upperBound, int dim) const
717 RVec shiftedCenter(centerCell);
718 // Shift the center to the cell coordinates of currCell, so that
719 // computeCutoffExtent() can assume simple rectangular grid.
724 if (currCell[ZZ] < 0)
726 shiftedCenter[XX] += cellShiftZX_;
728 else if (currCell[ZZ] >= ncelldim_[ZZ])
730 shiftedCenter[XX] -= cellShiftZX_;
732 if (currCell[YY] < 0)
734 shiftedCenter[XX] += cellShiftYX_;
736 else if (currCell[YY] >= ncelldim_[YY])
738 shiftedCenter[XX] -= cellShiftYX_;
741 if (dim == XX || dim == YY)
743 if (currCell[ZZ] < 0)
745 shiftedCenter[YY] += cellShiftZY_;
747 else if (currCell[ZZ] >= ncelldim_[ZZ])
749 shiftedCenter[YY] -= cellShiftZY_;
753 const real range = computeCutoffExtent(shiftedCenter, currCell, dim) * invCellSize_[dim];
754 real startOffset = shiftedCenter[dim] - range;
755 real endOffset = shiftedCenter[dim] + range;
756 // For non-periodic dimensions, clamp to the actual grid edges.
759 // If endOffset < 0 or startOffset > N, these may cause the whole
760 // test position/grid plane/grid row to be skipped.
765 const int cellCount = ncelldim_[dim];
766 if (endOffset > cellCount - 1)
768 endOffset = cellCount - 1;
771 currCell[dim] = static_cast<int>(floor(startOffset));
772 upperBound[dim] = static_cast<int>(floor(endOffset));
775 real AnalysisNeighborhoodSearchImpl::computeCutoffExtent(const RVec centerCell, const ivec cell, int dim) const
783 for (int d = dim + 1; d < DIM; ++d)
785 real dimDist = cell[d] - centerCell[d];
790 else if (dimDist <= 0)
794 dist2 += dimDist * dimDist * cellSize_[d] * cellSize_[d];
796 if (dist2 >= cutoff2_)
800 return std::sqrt(cutoff2_ - dist2);
803 bool AnalysisNeighborhoodSearchImpl::nextCell(const rvec centerCell, ivec cell, ivec upperBound) const
810 if (cell[dim] > upperBound[dim])
815 for (int d = dim - 1; d >= 0; --d)
817 initCellRange(centerCell, cell, upperBound, d);
818 if (cell[d] > upperBound[d])
829 int AnalysisNeighborhoodSearchImpl::shiftCell(const ivec cell, rvec shift) const
832 copy_ivec(cell, shiftedCell);
835 for (int d = 0; d < DIM; ++d)
837 const int cellCount = ncelldim_[d];
840 // A single shift may not be sufficient if the cell must be shifted
841 // in more than one dimension, although for each individual
842 // dimension it would be.
843 while (shiftedCell[d] < 0)
845 shiftedCell[d] += cellCount;
846 rvec_inc(shift, pbc_.box[d]);
848 while (shiftedCell[d] >= cellCount)
850 shiftedCell[d] -= cellCount;
851 rvec_dec(shift, pbc_.box[d]);
856 return getGridCellIndex(shiftedCell);
859 void AnalysisNeighborhoodSearchImpl::init(AnalysisNeighborhood::SearchMode mode,
861 const ListOfLists<int>* excls,
863 const AnalysisNeighborhoodPositions& positions)
865 GMX_RELEASE_ASSERT(positions.index_ == -1,
866 "Individual indexed positions not supported as reference");
868 if (bXY_ && pbc != nullptr && pbc->ePBC != epbcNONE)
870 if (pbc->ePBC != epbcXY && pbc->ePBC != epbcXYZ)
872 std::string message = formatString(
873 "Computations in the XY plane are not supported with PBC type '%s'",
874 epbc_names[pbc->ePBC]);
875 GMX_THROW(NotImplementedError(message));
877 if (pbc->ePBC == epbcXYZ
878 && (std::fabs(pbc->box[ZZ][XX]) > GMX_REAL_EPS * pbc->box[ZZ][ZZ]
879 || std::fabs(pbc->box[ZZ][YY]) > GMX_REAL_EPS * pbc->box[ZZ][ZZ]))
882 NotImplementedError("Computations in the XY plane are not supported when the "
883 "last box vector is not parallel to the Z axis"));
885 // Use a single grid cell in Z direction.
887 copy_mat(pbc->box, box);
889 set_pbc(&pbc_, epbcXY, box);
891 else if (pbc != nullptr)
897 pbc_.ePBC = epbcNONE;
900 nref_ = positions.count_;
901 if (mode == AnalysisNeighborhood::eSearchMode_Simple)
907 bGrid_ = initGrid(pbc_, positions.count_, positions.x_,
908 mode == AnalysisNeighborhood::eSearchMode_Grid);
910 refIndices_ = positions.indices_;
913 xrefAlloc_.resize(nref_);
914 xref_ = as_rvec_array(xrefAlloc_.data());
916 for (int i = 0; i < nref_; ++i)
918 const int ii = (refIndices_ != nullptr) ? refIndices_[i] : i;
920 mapPointToGridCell(positions.x_[ii], refcell, xrefAlloc_[i]);
921 addToGridCell(refcell, i);
924 else if (refIndices_ != nullptr)
926 xrefAlloc_.resize(nref_);
927 xref_ = as_rvec_array(xrefAlloc_.data());
928 for (int i = 0; i < nref_; ++i)
930 copy_rvec(positions.x_[refIndices_[i]], xrefAlloc_[i]);
935 xref_ = positions.x_;
938 refExclusionIds_ = nullptr;
939 if (excls != nullptr)
941 // TODO: Check that the IDs are ascending, or remove the limitation.
942 refExclusionIds_ = positions.exclusionIds_;
943 GMX_RELEASE_ASSERT(refExclusionIds_ != nullptr,
944 "Exclusion IDs must be set for reference positions "
945 "when exclusions are enabled");
949 /********************************************************************
950 * AnalysisNeighborhoodPairSearchImpl
953 void AnalysisNeighborhoodPairSearchImpl::reset(int testIndex)
955 testIndex_ = testIndex;
962 if (testIndex_ >= 0 && testIndex_ < testPosCount_)
964 const int index = (testIndices_ != nullptr ? testIndices_[testIndex] : testIndex);
967 search_.mapPointToGridCell(testPositions_[index], testcell_, xtest_);
968 search_.initCellRange(testcell_, currCell_, cellBound_, ZZ);
969 search_.initCellRange(testcell_, currCell_, cellBound_, YY);
970 search_.initCellRange(testcell_, currCell_, cellBound_, XX);
973 testCellIndex_ = search_.getGridCellIndex(testcell_);
978 copy_rvec(testPositions_[index], xtest_);
984 if (search_.excls_ != nullptr)
986 const int exclIndex = testExclusionIds_[index];
987 if (exclIndex < search_.excls_->ssize())
989 excl_ = (*search_.excls_)[exclIndex];
993 excl_ = ArrayRef<const int>();
999 void AnalysisNeighborhoodPairSearchImpl::nextTestPosition()
1001 if (testIndex_ < testPosCount_)
1008 bool AnalysisNeighborhoodPairSearchImpl::isExcluded(int j)
1010 const int nexcl = excl_.ssize();
1011 if (exclind_ < nexcl)
1013 const int index = (search_.refIndices_ != nullptr ? search_.refIndices_[j] : j);
1014 const int refId = search_.refExclusionIds_[index];
1015 while (exclind_ < nexcl && excl_[exclind_] < refId)
1019 if (exclind_ < nexcl && refId == excl_[exclind_])
1028 void AnalysisNeighborhoodPairSearchImpl::startSearch(const AnalysisNeighborhoodPositions& positions)
1030 selfSearchMode_ = false;
1031 testPosCount_ = positions.count_;
1032 testPositions_ = positions.x_;
1033 testExclusionIds_ = positions.exclusionIds_;
1034 testIndices_ = positions.indices_;
1035 GMX_RELEASE_ASSERT(search_.excls_ == nullptr || testExclusionIds_ != nullptr,
1036 "Exclusion IDs must be set when exclusions are enabled");
1037 if (positions.index_ < 0)
1043 // Somewhat of a hack: setup the array such that only the last position
1045 testPosCount_ = positions.index_ + 1;
1046 reset(positions.index_);
1050 void AnalysisNeighborhoodPairSearchImpl::startSelfSearch()
1052 selfSearchMode_ = true;
1053 testPosCount_ = search_.nref_;
1054 testPositions_ = search_.xref_;
1055 testExclusionIds_ = search_.refExclusionIds_;
1056 testIndices_ = search_.refIndices_;
1057 GMX_RELEASE_ASSERT(search_.excls_ == nullptr || testIndices_ == nullptr,
1058 "Exclusion IDs not implemented with indexed ref positions");
1062 template<class Action>
1063 bool AnalysisNeighborhoodPairSearchImpl::searchNext(Action action)
1065 while (testIndex_ < testPosCount_)
1069 int cai = prevcai_ + 1;
1074 const int ci = search_.shiftCell(currCell_, shift);
1075 if (selfSearchMode_ && ci > testCellIndex_)
1079 const int cellSize = ssize(search_.cells_[ci]);
1080 for (; cai < cellSize; ++cai)
1082 const int i = search_.cells_[ci][cai];
1083 if (selfSearchMode_ && ci == testCellIndex_ && i >= testIndex_)
1092 rvec_sub(search_.xref_[i], xtest_, dx);
1093 rvec_sub(dx, shift, dx);
1094 const real r2 = search_.bXY_ ? dx[XX] * dx[XX] + dx[YY] * dx[YY] : norm2(dx);
1095 if (r2 <= search_.cutoff2_)
1097 if (action(i, r2, dx))
1102 copy_rvec(dx, prevdx_);
1109 } while (search_.nextCell(testcell_, currCell_, cellBound_));
1113 for (int i = previ_ + 1; i < search_.nref_; ++i)
1120 if (search_.pbc_.ePBC != epbcNONE)
1122 pbc_dx(&search_.pbc_, search_.xref_[i], xtest_, dx);
1126 rvec_sub(search_.xref_[i], xtest_, dx);
1128 const real r2 = search_.bXY_ ? dx[XX] * dx[XX] + dx[YY] * dx[YY] : norm2(dx);
1129 if (r2 <= search_.cutoff2_)
1131 if (action(i, r2, dx))
1135 copy_rvec(dx, prevdx_);
1146 void AnalysisNeighborhoodPairSearchImpl::initFoundPair(AnalysisNeighborhoodPair* pair) const
1150 *pair = AnalysisNeighborhoodPair();
1154 *pair = AnalysisNeighborhoodPair(previ_, testIndex_, prevr2_, prevdx_);
1158 } // namespace internal
1164 * Search action to find the next neighbor.
1166 * Used as the action for AnalysisNeighborhoodPairSearchImpl::searchNext() to
1167 * find the next neighbor.
1169 * Simply breaks the loop on the first found neighbor.
1171 bool withinAction(int /*i*/, real /*r2*/, const rvec /* dx */)
1177 * Search action find the minimum distance.
1179 * Used as the action for AnalysisNeighborhoodPairSearchImpl::searchNext() to
1180 * find the nearest neighbor.
1182 * With this action, AnalysisNeighborhoodPairSearchImpl::searchNext() always
1183 * returns false, and the output is put into the variables passed by pointer
1184 * into the constructor. If no neighbors are found, the output variables are
1185 * not modified, i.e., the caller must initialize them.
1191 * Initializes the action with given output locations.
1193 * \param[out] closestPoint Index of the closest reference location.
1194 * \param[out] minDist2 Minimum distance squared.
1195 * \param[out] dx Shortest distance vector.
1197 * The constructor call does not modify the pointed values, but only
1198 * stores the pointers for later use.
1199 * See the class description for additional semantics.
1201 MindistAction(int* closestPoint, real* minDist2, rvec* dx) // NOLINT(readability-non-const-parameter)
1203 closestPoint_(*closestPoint),
1204 minDist2_(*minDist2),
1208 //! Copies the action.
1209 MindistAction(const MindistAction&) = default;
1211 //! Processes a neighbor to find the nearest point.
1212 bool operator()(int i, real r2, const rvec dx)
1228 GMX_DISALLOW_ASSIGN(MindistAction);
1233 /********************************************************************
1234 * AnalysisNeighborhood::Impl
1237 class AnalysisNeighborhood::Impl
1240 typedef AnalysisNeighborhoodSearch::ImplPointer SearchImplPointer;
1241 typedef std::vector<SearchImplPointer> SearchList;
1243 Impl() : cutoff_(0), excls_(nullptr), mode_(eSearchMode_Automatic), bXY_(false) {}
1246 SearchList::const_iterator i;
1247 for (i = searchList_.begin(); i != searchList_.end(); ++i)
1249 GMX_RELEASE_ASSERT(i->unique(), "Dangling AnalysisNeighborhoodSearch reference");
1253 SearchImplPointer getSearch();
1255 Mutex createSearchMutex_;
1256 SearchList searchList_;
1258 const ListOfLists<int>* excls_;
1263 AnalysisNeighborhood::Impl::SearchImplPointer AnalysisNeighborhood::Impl::getSearch()
1265 lock_guard<Mutex> lock(createSearchMutex_);
1266 // TODO: Consider whether this needs to/can be faster, e.g., by keeping a
1267 // separate pool of unused search objects.
1268 SearchList::const_iterator i;
1269 for (i = searchList_.begin(); i != searchList_.end(); ++i)
1276 SearchImplPointer search(new internal::AnalysisNeighborhoodSearchImpl(cutoff_));
1277 searchList_.push_back(search);
1281 /********************************************************************
1282 * AnalysisNeighborhood
1285 AnalysisNeighborhood::AnalysisNeighborhood() : impl_(new Impl) {}
1287 AnalysisNeighborhood::~AnalysisNeighborhood() {}
1289 void AnalysisNeighborhood::setCutoff(real cutoff)
1291 GMX_RELEASE_ASSERT(impl_->searchList_.empty(),
1292 "Changing the cutoff after initSearch() not currently supported");
1293 impl_->cutoff_ = cutoff;
1296 void AnalysisNeighborhood::setXYMode(bool bXY)
1301 void AnalysisNeighborhood::setTopologyExclusions(const ListOfLists<int>* excls)
1303 GMX_RELEASE_ASSERT(impl_->searchList_.empty(),
1304 "Changing the exclusions after initSearch() not currently supported");
1305 impl_->excls_ = excls;
1308 void AnalysisNeighborhood::setMode(SearchMode mode)
1310 impl_->mode_ = mode;
1313 AnalysisNeighborhood::SearchMode AnalysisNeighborhood::mode() const
1315 return impl_->mode_;
1318 AnalysisNeighborhoodSearch AnalysisNeighborhood::initSearch(const t_pbc* pbc,
1319 const AnalysisNeighborhoodPositions& positions)
1321 Impl::SearchImplPointer search(impl_->getSearch());
1322 search->init(mode(), impl_->bXY_, impl_->excls_, pbc, positions);
1323 return AnalysisNeighborhoodSearch(search);
1326 /********************************************************************
1327 * AnalysisNeighborhoodSearch
1330 AnalysisNeighborhoodSearch::AnalysisNeighborhoodSearch() {}
1332 AnalysisNeighborhoodSearch::AnalysisNeighborhoodSearch(const ImplPointer& impl) : impl_(impl) {}
1334 void AnalysisNeighborhoodSearch::reset()
1339 AnalysisNeighborhood::SearchMode AnalysisNeighborhoodSearch::mode() const
1341 GMX_RELEASE_ASSERT(impl_, "Accessing an invalid search object");
1342 return (impl_->usesGridSearch() ? AnalysisNeighborhood::eSearchMode_Grid
1343 : AnalysisNeighborhood::eSearchMode_Simple);
1346 bool AnalysisNeighborhoodSearch::isWithin(const AnalysisNeighborhoodPositions& positions) const
1348 GMX_RELEASE_ASSERT(impl_, "Accessing an invalid search object");
1349 internal::AnalysisNeighborhoodPairSearchImpl pairSearch(*impl_);
1350 pairSearch.startSearch(positions);
1351 return pairSearch.searchNext(&withinAction);
1354 real AnalysisNeighborhoodSearch::minimumDistance(const AnalysisNeighborhoodPositions& positions) const
1356 GMX_RELEASE_ASSERT(impl_, "Accessing an invalid search object");
1357 internal::AnalysisNeighborhoodPairSearchImpl pairSearch(*impl_);
1358 pairSearch.startSearch(positions);
1359 real minDist2 = impl_->cutoffSquared();
1360 int closestPoint = -1;
1361 rvec dx = { 0.0, 0.0, 0.0 };
1362 MindistAction action(&closestPoint, &minDist2, &dx);
1363 (void)pairSearch.searchNext(action);
1364 return std::sqrt(minDist2);
1367 AnalysisNeighborhoodPair AnalysisNeighborhoodSearch::nearestPoint(const AnalysisNeighborhoodPositions& positions) const
1369 GMX_RELEASE_ASSERT(impl_, "Accessing an invalid search object");
1370 internal::AnalysisNeighborhoodPairSearchImpl pairSearch(*impl_);
1371 pairSearch.startSearch(positions);
1372 real minDist2 = impl_->cutoffSquared();
1373 int closestPoint = -1;
1374 rvec dx = { 0.0, 0.0, 0.0 };
1375 MindistAction action(&closestPoint, &minDist2, &dx);
1376 (void)pairSearch.searchNext(action);
1377 return AnalysisNeighborhoodPair(closestPoint, 0, minDist2, dx);
1380 AnalysisNeighborhoodPairSearch AnalysisNeighborhoodSearch::startSelfPairSearch() const
1382 GMX_RELEASE_ASSERT(impl_, "Accessing an invalid search object");
1383 Impl::PairSearchImplPointer pairSearch(impl_->getPairSearch());
1384 pairSearch->startSelfSearch();
1385 return AnalysisNeighborhoodPairSearch(pairSearch);
1388 AnalysisNeighborhoodPairSearch
1389 AnalysisNeighborhoodSearch::startPairSearch(const AnalysisNeighborhoodPositions& positions) const
1391 GMX_RELEASE_ASSERT(impl_, "Accessing an invalid search object");
1392 Impl::PairSearchImplPointer pairSearch(impl_->getPairSearch());
1393 pairSearch->startSearch(positions);
1394 return AnalysisNeighborhoodPairSearch(pairSearch);
1397 /********************************************************************
1398 * AnalysisNeighborhoodPairSearch
1401 AnalysisNeighborhoodPairSearch::AnalysisNeighborhoodPairSearch(const ImplPointer& impl) :
1406 bool AnalysisNeighborhoodPairSearch::findNextPair(AnalysisNeighborhoodPair* pair)
1408 bool bFound = impl_->searchNext(&withinAction);
1409 impl_->initFoundPair(pair);
1413 void AnalysisNeighborhoodPairSearch::skipRemainingPairsForTestPosition()
1415 impl_->nextTestPosition();