2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2009,2010,2011,2012,2013, 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).
40 * The grid implementation could still be optimized in several different ways:
41 * - Triclinic grid cells are not the most efficient shape, but make PBC
43 * - Precalculating the required PBC shift for a pair of cells outside the
44 * inner loop. After this is done, it should be quite straightforward to
45 * move to rectangular cells.
46 * - Pruning grid cells from the search list if they are completely outside
47 * the sphere that is being considered.
48 * - A better heuristic could be added for falling back to simple loops for a
49 * small number of reference particles.
50 * - A better heuristic for selecting the grid size.
51 * - A multi-level grid implementation could be used to be able to use small
52 * grids for short cutoffs with very inhomogeneous particle distributions
53 * without a memory cost.
55 * \author Teemu Murtola <teemu.murtola@gmail.com>
56 * \ingroup module_selection
58 #include "gromacs/selection/nbsearch.h"
65 #include "gromacs/legacyheaders/smalloc.h"
66 #include "gromacs/legacyheaders/typedefs.h"
67 #include "gromacs/legacyheaders/pbc.h"
68 #include "gromacs/legacyheaders/vec.h"
69 #include "gromacs/legacyheaders/thread_mpi/mutex.h"
71 #include "gromacs/selection/position.h"
72 #include "gromacs/utility/arrayref.h"
73 #include "gromacs/utility/gmxassert.h"
81 /********************************************************************
82 * Implementation class declarations
85 class AnalysisNeighborhoodSearchImpl
88 typedef AnalysisNeighborhoodPairSearch::ImplPointer
89 PairSearchImplPointer;
90 typedef std::vector<PairSearchImplPointer> PairSearchList;
91 typedef std::vector<std::vector<int> > CellList;
93 explicit AnalysisNeighborhoodSearchImpl(real cutoff);
94 ~AnalysisNeighborhoodSearchImpl();
97 * Initializes the search with a given box and reference positions.
99 * \param[in] mode Search mode to use.
100 * \param[in] pbc PBC information.
101 * \param[in] positions Set of reference positions.
103 void init(AnalysisNeighborhood::SearchMode mode,
105 const AnalysisNeighborhoodPositions &positions);
106 PairSearchImplPointer getPairSearch();
108 real cutoffSquared() const { return cutoff2_; }
109 bool usesGridSearch() const { return bGrid_; }
112 //! Calculates offsets to neighboring grid cells that should be considered.
113 void initGridCellNeighborList();
115 * Determines a suitable grid size and sets up the cells.
117 * \param[in] pbc Information about the box.
118 * \returns false if grid search is not suitable.
120 bool initGridCells(const t_pbc *pbc);
122 * Sets ua a search grid for a given box.
124 * \param[in] pbc Information about the box.
125 * \returns false if grid search is not suitable.
127 bool initGrid(const t_pbc *pbc);
129 * Maps a point into a grid cell.
131 * \param[in] x Point to map.
132 * \param[out] cell Indices of the grid cell in which \p x lies.
134 * \p x should be within the triclinic unit cell.
136 void mapPointToGridCell(const rvec x, ivec cell) const;
138 * Calculates linear index of a grid cell.
140 * \param[in] cell Cell indices.
141 * \returns Linear index of \p cell.
143 int getGridCellIndex(const ivec cell) const;
145 * Adds an index into a grid cell.
147 * \param[in] cell Cell into which \p i should be added.
148 * \param[in] i Index to add.
150 void addToGridCell(const ivec cell, int i);
152 //! Whether to try grid searching.
156 //! The cutoff squared.
159 //! Number of reference points for the current frame.
161 //! Reference point positions.
163 //! Reference position ids (NULL if not available).
168 //! Number of excluded reference positions for current test particle.
170 //! Exclusions for current test particle.
173 //! Whether grid searching is actually used for the current positions.
175 //! Array allocated for storing in-unit-cell reference positions.
177 //! Allocation count for xref_alloc.
179 //! false if the box is rectangular.
181 //! Box vectors of a single grid cell.
183 //! The reciprocal cell vectors as columns; the inverse of \p cellbox.
185 //! Number of cells along each dimension.
187 //! Data structure to hold the grid cell contents.
189 //! Number of neighboring cells to consider.
191 //! Offsets of the neighboring cells to consider.
193 //! Allocation count for \p gnboffs.
196 tMPI::mutex createPairSearchMutex_;
197 PairSearchList pairSearchList_;
199 friend class AnalysisNeighborhoodPairSearchImpl;
201 GMX_DISALLOW_COPY_AND_ASSIGN(AnalysisNeighborhoodSearchImpl);
204 class AnalysisNeighborhoodPairSearchImpl
207 explicit AnalysisNeighborhoodPairSearchImpl(const AnalysisNeighborhoodSearchImpl &search)
211 clear_ivec(testcell_);
215 //! Initializes a search to find reference positions neighboring \p x.
216 void startSearch(const AnalysisNeighborhoodPositions &positions);
217 //! Searches for the next neighbor.
218 template <class Action>
219 bool searchNext(Action action);
220 //! Initializes a pair representing the pair found by searchNext().
221 void initFoundPair(AnalysisNeighborhoodPair *pair) const;
222 //! Advances to the next test position, skipping any remaining pairs.
223 void nextTestPosition();
226 //! Clears the loop indices.
227 void reset(int testIndex);
228 //! Checks whether a reference positiong should be excluded.
229 bool isExcluded(int j);
231 //! Parent search object.
232 const AnalysisNeighborhoodSearchImpl &search_;
233 //! Reference to the test positions.
234 ConstArrayRef<rvec> testPositions_;
235 //! Index of the currently active test position in \p testPositions_.
237 //! Stores test position during a pair loop.
239 //! Stores the previous returned position during a pair loop.
241 //! Stores the current exclusion index during loops.
243 //! Stores the test particle cell index during loops.
245 //! Stores the current cell neighbor index during pair loops.
247 //! Stores the index within the current cell during pair loops.
250 GMX_DISALLOW_COPY_AND_ASSIGN(AnalysisNeighborhoodPairSearchImpl);
253 /********************************************************************
254 * AnalysisNeighborhoodSearchImpl
257 AnalysisNeighborhoodSearchImpl::AnalysisNeighborhoodSearchImpl(real cutoff)
263 cutoff_ = GMX_REAL_MAX;
266 cutoff2_ = sqr(cutoff_);
282 clear_mat(recipcell_);
283 clear_ivec(ncelldim_);
290 AnalysisNeighborhoodSearchImpl::~AnalysisNeighborhoodSearchImpl()
292 PairSearchList::const_iterator i;
293 for (i = pairSearchList_.begin(); i != pairSearchList_.end(); ++i)
295 GMX_RELEASE_ASSERT(i->unique(),
296 "Dangling AnalysisNeighborhoodPairSearch reference");
302 AnalysisNeighborhoodSearchImpl::PairSearchImplPointer
303 AnalysisNeighborhoodSearchImpl::getPairSearch()
305 tMPI::lock_guard<tMPI::mutex> lock(createPairSearchMutex_);
306 // TODO: Consider whether this needs to/can be faster, e.g., by keeping a
307 // separate pool of unused search objects.
308 PairSearchList::const_iterator i;
309 for (i = pairSearchList_.begin(); i != pairSearchList_.end(); ++i)
316 PairSearchImplPointer pairSearch(new AnalysisNeighborhoodPairSearchImpl(*this));
317 pairSearchList_.push_back(pairSearch);
321 void AnalysisNeighborhoodSearchImpl::initGridCellNeighborList()
323 int maxx, maxy, maxz;
326 /* Find the extent of the sphere in triclinic coordinates */
327 maxz = static_cast<int>(cutoff_ * recipcell_[ZZ][ZZ]) + 1;
328 rvnorm = sqrt(sqr(recipcell_[YY][YY]) + sqr(recipcell_[ZZ][YY]));
329 maxy = static_cast<int>(cutoff_ * rvnorm) + 1;
330 rvnorm = sqrt(sqr(recipcell_[XX][XX]) + sqr(recipcell_[YY][XX])
331 + sqr(recipcell_[ZZ][XX]));
332 maxx = static_cast<int>(cutoff_ * rvnorm) + 1;
334 /* Calculate the number of cells and reallocate if necessary */
335 ngridnb_ = (2 * maxx + 1) * (2 * maxy + 1) * (2 * maxz + 1);
336 if (gnboffs_nalloc_ < ngridnb_)
338 gnboffs_nalloc_ = ngridnb_;
339 srenew(gnboffs_, gnboffs_nalloc_);
342 /* Store the whole cube */
343 /* TODO: Prune off corners that are not needed */
345 for (int x = -maxx; x <= maxx; ++x)
347 for (int y = -maxy; y <= maxy; ++y)
349 for (int z = -maxz; z <= maxz; ++z)
360 bool AnalysisNeighborhoodSearchImpl::initGridCells(const t_pbc *pbc)
362 const real targetsize =
363 pow(pbc->box[XX][XX] * pbc->box[YY][YY] * pbc->box[ZZ][ZZ]
364 * 10 / nref_, static_cast<real>(1./3.));
367 for (int dd = 0; dd < DIM; ++dd)
369 ncelldim_[dd] = static_cast<int>(pbc->box[dd][dd] / targetsize);
370 cellCount *= ncelldim_[dd];
371 if (ncelldim_[dd] < 3)
376 // Never decrease the size of the cell vector to avoid reallocating
377 // memory for the nested vectors. The actual size of the vector is not
378 // used outside this function.
379 if (cells_.size() < static_cast<size_t>(cellCount))
381 cells_.resize(cellCount);
383 for (int ci = 0; ci < cellCount; ++ci)
390 bool AnalysisNeighborhoodSearchImpl::initGrid(const t_pbc *pbc)
392 /* TODO: This check could be improved. */
393 if (0.5*pbc->max_cutoff2 < cutoff2_)
398 if (!initGridCells(pbc))
403 bTric_ = TRICLINIC(pbc->box);
406 for (int dd = 0; dd < DIM; ++dd)
408 svmul(1.0 / ncelldim_[dd], pbc->box[dd], cellbox_[dd]);
410 m_inv_ur0(cellbox_, recipcell_);
414 for (int dd = 0; dd < DIM; ++dd)
416 cellbox_[dd][dd] = pbc->box[dd][dd] / ncelldim_[dd];
417 recipcell_[dd][dd] = 1.0 / cellbox_[dd][dd];
420 initGridCellNeighborList();
424 void AnalysisNeighborhoodSearchImpl::mapPointToGridCell(const rvec x,
431 tmvmul_ur0(recipcell_, x, tx);
432 for (int dd = 0; dd < DIM; ++dd)
434 cell[dd] = static_cast<int>(tx[dd]);
439 for (int dd = 0; dd < DIM; ++dd)
441 cell[dd] = static_cast<int>(x[dd] * recipcell_[dd][dd]);
446 int AnalysisNeighborhoodSearchImpl::getGridCellIndex(const ivec cell) const
448 GMX_ASSERT(cell[XX] >= 0 && cell[XX] < ncelldim_[XX],
449 "Grid cell X index out of range");
450 GMX_ASSERT(cell[YY] >= 0 && cell[YY] < ncelldim_[YY],
451 "Grid cell Y index out of range");
452 GMX_ASSERT(cell[ZZ] >= 0 && cell[ZZ] < ncelldim_[ZZ],
453 "Grid cell Z index out of range");
454 return cell[XX] + cell[YY] * ncelldim_[XX]
455 + cell[ZZ] * ncelldim_[XX] * ncelldim_[YY];
458 void AnalysisNeighborhoodSearchImpl::addToGridCell(const ivec cell, int i)
460 const int ci = getGridCellIndex(cell);
461 cells_[ci].push_back(i);
464 void AnalysisNeighborhoodSearchImpl::init(
465 AnalysisNeighborhood::SearchMode mode,
467 const AnalysisNeighborhoodPositions &positions)
469 GMX_RELEASE_ASSERT(positions.index_ == -1,
470 "Individual indexed positions not supported as reference");
471 pbc_ = const_cast<t_pbc *>(pbc);
472 nref_ = positions.count_;
473 // TODO: Consider whether it would be possible to support grid searching in
475 if (mode == AnalysisNeighborhood::eSearchMode_Simple
476 || pbc_ == NULL || pbc_->ePBC != epbcXYZ)
482 // TODO: Actually implement forcing eSearchMode_Grid
483 bGrid_ = initGrid(pbc_);
487 if (xref_nalloc_ < nref_)
489 srenew(xref_alloc_, nref_);
490 xref_nalloc_ = nref_;
494 for (int i = 0; i < nref_; ++i)
496 copy_rvec(positions.x_[i], xref_alloc_[i]);
498 put_atoms_in_triclinic_unitcell(ecenterTRIC, pbc_->box,
500 for (int i = 0; i < nref_; ++i)
504 mapPointToGridCell(xref_[i], refcell);
505 addToGridCell(refcell, i);
510 xref_ = positions.x_;
512 // TODO: Once exclusions are supported, this may need to be initialized.
518 * Sets the exclusions for the next neighborhood search.
520 * \param[in,out] d Neighborhood search data structure.
521 * \param[in] nexcl Number of reference positions to exclude from next
523 * \param[in] excl Indices of reference positions to exclude.
525 * The set exclusions remain in effect until the next call of this function.
528 gmx_ana_nbsearch_set_excl(gmx_ana_nbsearch_t *d, int nexcl, int excl[])
536 /********************************************************************
537 * AnalysisNeighborhoodPairSearchImpl
540 void AnalysisNeighborhoodPairSearchImpl::reset(int testIndex)
542 testIndex_ = testIndex;
543 if (testIndex_ >= 0 && testIndex_ < static_cast<int>(testPositions_.size()))
545 copy_rvec(testPositions_[testIndex_], xtest_);
548 put_atoms_in_triclinic_unitcell(ecenterTRIC, search_.pbc_->box,
550 search_.mapPointToGridCell(xtest_, testcell_);
559 void AnalysisNeighborhoodPairSearchImpl::nextTestPosition()
561 if (testIndex_ < static_cast<int>(testPositions_.size()))
568 bool AnalysisNeighborhoodPairSearchImpl::isExcluded(int j)
570 if (exclind_ < search_.nexcl_)
574 while (exclind_ < search_.nexcl_
575 && search_.excl_[exclind_] < search_.refid_[j])
579 if (exclind_ < search_.nexcl_
580 && search_.refid_[j] == search_.excl_[exclind_])
588 while (search_.bGrid_ && exclind_ < search_.nexcl_
589 && search_.excl_[exclind_] < j)
593 if (search_.excl_[exclind_] == j)
603 void AnalysisNeighborhoodPairSearchImpl::startSearch(
604 const AnalysisNeighborhoodPositions &positions)
606 if (positions.index_ < 0)
608 testPositions_ = ConstArrayRef<rvec>(positions.x_, positions.count_);
613 // Somewhat of a hack: setup the array such that only the last position
615 testPositions_ = ConstArrayRef<rvec>(positions.x_, positions.index_ + 1);
616 reset(positions.index_);
620 template <class Action>
621 bool AnalysisNeighborhoodPairSearchImpl::searchNext(Action action)
623 while (testIndex_ < static_cast<int>(testPositions_.size()))
628 int cai = prevcai_ + 1;
630 for (; nbi < search_.ngridnb_; ++nbi)
634 ivec_add(testcell_, search_.gnboffs_[nbi], cell);
635 cell[XX] = (cell[XX] + search_.ncelldim_[XX]) % search_.ncelldim_[XX];
636 cell[YY] = (cell[YY] + search_.ncelldim_[YY]) % search_.ncelldim_[YY];
637 cell[ZZ] = (cell[ZZ] + search_.ncelldim_[ZZ]) % search_.ncelldim_[ZZ];
639 const int ci = search_.getGridCellIndex(cell);
640 const int cellSize = static_cast<int>(search_.cells_[ci].size());
641 /* TODO: Calculate the required PBC shift outside the inner loop */
642 for (; cai < cellSize; ++cai)
644 const int i = search_.cells_[ci][cai];
650 pbc_dx_aiuc(search_.pbc_, xtest_, search_.xref_[i], dx);
651 const real r2 = norm2(dx);
652 if (r2 <= search_.cutoff2_)
669 for (int i = previ_ + 1; i < search_.nref_; ++i)
678 pbc_dx(search_.pbc_, xtest_, search_.xref_[i], dx);
682 rvec_sub(xtest_, search_.xref_[i], dx);
684 const real r2 = norm2(dx);
685 if (r2 <= search_.cutoff2_)
700 void AnalysisNeighborhoodPairSearchImpl::initFoundPair(
701 AnalysisNeighborhoodPair *pair) const
705 *pair = AnalysisNeighborhoodPair();
709 *pair = AnalysisNeighborhoodPair(previ_, testIndex_);
713 } // namespace internal
719 * Search action to find the next neighbor.
721 * Used as the action for AnalysisNeighborhoodPairSearchImpl::searchNext() to
722 * find the next neighbor.
724 * Simply breaks the loop on the first found neighbor.
726 bool withinAction(int /*i*/, real /*r2*/)
732 * Search action find the minimum distance.
734 * Used as the action for AnalysisNeighborhoodPairSearchImpl::searchNext() to
735 * find the nearest neighbor.
737 * With this action, AnalysisNeighborhoodPairSearchImpl::searchNext() always
738 * returns false, and the output is put into the variables passed by pointer
739 * into the constructor. If no neighbors are found, the output variables are
740 * not modified, i.e., the caller must initialize them.
746 * Initializes the action with given output locations.
748 * \param[out] closestPoint Index of the closest reference location.
749 * \param[out] minDist2 Minimum distance squared.
751 * The constructor call does not modify the pointed values, but only
752 * stores the pointers for later use.
753 * See the class description for additional semantics.
755 MindistAction(int *closestPoint, real *minDist2)
756 : closestPoint_(*closestPoint), minDist2_(*minDist2)
760 //! Processes a neighbor to find the nearest point.
761 bool operator()(int i, real r2)
775 GMX_DISALLOW_ASSIGN(MindistAction);
780 /********************************************************************
781 * AnalysisNeighborhood::Impl
784 class AnalysisNeighborhood::Impl
787 typedef AnalysisNeighborhoodSearch::ImplPointer SearchImplPointer;
788 typedef std::vector<SearchImplPointer> SearchList;
790 Impl() : cutoff_(0), mode_(eSearchMode_Automatic)
795 SearchList::const_iterator i;
796 for (i = searchList_.begin(); i != searchList_.end(); ++i)
798 GMX_RELEASE_ASSERT(i->unique(),
799 "Dangling AnalysisNeighborhoodSearch reference");
803 SearchImplPointer getSearch();
805 tMPI::mutex createSearchMutex_;
806 SearchList searchList_;
811 AnalysisNeighborhood::Impl::SearchImplPointer
812 AnalysisNeighborhood::Impl::getSearch()
814 tMPI::lock_guard<tMPI::mutex> lock(createSearchMutex_);
815 // TODO: Consider whether this needs to/can be faster, e.g., by keeping a
816 // separate pool of unused search objects.
817 SearchList::const_iterator i;
818 for (i = searchList_.begin(); i != searchList_.end(); ++i)
825 SearchImplPointer search(new internal::AnalysisNeighborhoodSearchImpl(cutoff_));
826 searchList_.push_back(search);
830 /********************************************************************
831 * AnalysisNeighborhood
834 AnalysisNeighborhood::AnalysisNeighborhood()
839 AnalysisNeighborhood::~AnalysisNeighborhood()
843 void AnalysisNeighborhood::setCutoff(real cutoff)
845 GMX_RELEASE_ASSERT(impl_->searchList_.empty(),
846 "Changing the cutoff after initSearch() not currently supported");
847 impl_->cutoff_ = cutoff;
850 void AnalysisNeighborhood::setMode(SearchMode mode)
855 AnalysisNeighborhood::SearchMode AnalysisNeighborhood::mode() const
860 AnalysisNeighborhoodSearch
861 AnalysisNeighborhood::initSearch(const t_pbc *pbc,
862 const AnalysisNeighborhoodPositions &positions)
864 Impl::SearchImplPointer search(impl_->getSearch());
865 search->init(mode(), pbc, positions);
866 return AnalysisNeighborhoodSearch(search);
869 /********************************************************************
870 * AnalysisNeighborhoodSearch
873 AnalysisNeighborhoodSearch::AnalysisNeighborhoodSearch()
877 AnalysisNeighborhoodSearch::AnalysisNeighborhoodSearch(const ImplPointer &impl)
882 void AnalysisNeighborhoodSearch::reset()
887 AnalysisNeighborhood::SearchMode AnalysisNeighborhoodSearch::mode() const
889 GMX_RELEASE_ASSERT(impl_, "Accessing an invalid search object");
890 return (impl_->usesGridSearch()
891 ? AnalysisNeighborhood::eSearchMode_Grid
892 : AnalysisNeighborhood::eSearchMode_Simple);
895 bool AnalysisNeighborhoodSearch::isWithin(
896 const AnalysisNeighborhoodPositions &positions) const
898 GMX_RELEASE_ASSERT(impl_, "Accessing an invalid search object");
899 internal::AnalysisNeighborhoodPairSearchImpl pairSearch(*impl_);
900 pairSearch.startSearch(positions);
901 return pairSearch.searchNext(&withinAction);
904 real AnalysisNeighborhoodSearch::minimumDistance(
905 const AnalysisNeighborhoodPositions &positions) const
907 GMX_RELEASE_ASSERT(impl_, "Accessing an invalid search object");
908 internal::AnalysisNeighborhoodPairSearchImpl pairSearch(*impl_);
909 pairSearch.startSearch(positions);
910 real minDist2 = impl_->cutoffSquared();
911 int closestPoint = -1;
912 MindistAction action(&closestPoint, &minDist2);
913 (void)pairSearch.searchNext(action);
914 return sqrt(minDist2);
917 AnalysisNeighborhoodPair
918 AnalysisNeighborhoodSearch::nearestPoint(
919 const AnalysisNeighborhoodPositions &positions) const
921 GMX_RELEASE_ASSERT(impl_, "Accessing an invalid search object");
922 internal::AnalysisNeighborhoodPairSearchImpl pairSearch(*impl_);
923 pairSearch.startSearch(positions);
924 real minDist2 = impl_->cutoffSquared();
925 int closestPoint = -1;
926 MindistAction action(&closestPoint, &minDist2);
927 (void)pairSearch.searchNext(action);
928 return AnalysisNeighborhoodPair(closestPoint, 0);
931 AnalysisNeighborhoodPairSearch
932 AnalysisNeighborhoodSearch::startPairSearch(
933 const AnalysisNeighborhoodPositions &positions) const
935 GMX_RELEASE_ASSERT(impl_, "Accessing an invalid search object");
936 Impl::PairSearchImplPointer pairSearch(impl_->getPairSearch());
937 pairSearch->startSearch(positions);
938 return AnalysisNeighborhoodPairSearch(pairSearch);
941 /********************************************************************
942 * AnalysisNeighborhoodPairSearch
945 AnalysisNeighborhoodPairSearch::AnalysisNeighborhoodPairSearch(
946 const ImplPointer &impl)
951 bool AnalysisNeighborhoodPairSearch::findNextPair(AnalysisNeighborhoodPair *pair)
953 bool bFound = impl_->searchNext(&withinAction);
954 impl_->initFoundPair(pair);
958 void AnalysisNeighborhoodPairSearch::skipRemainingPairsForTestPosition()
960 impl_->nextTestPosition();