/*
* This file is part of the GROMACS molecular simulation package.
*
- * Copyright (c) 2009,2010,2011,2012, by the GROMACS development team, led by
- * David van der Spoel, Berk Hess, Erik Lindahl, and including many
- * others, as listed in the AUTHORS file in the top-level source
- * directory and at http://www.gromacs.org.
+ * Copyright (c) 2009,2010,2011,2012,2013,2014, by the GROMACS development team, led by
+ * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
+ * and including many others, as listed in the AUTHORS file in the
+ * top-level source directory and at http://www.gromacs.org.
*
* GROMACS is free software; you can redistribute it and/or
* modify it under the terms of the GNU Lesser General Public License
* To help us fund GROMACS development, we humbly ask that you cite
* the research papers on the package. Check out http://www.gromacs.org.
*/
-/*! \page nbsearch Neighborhood search routines
- *
- * Functions to find particles within a neighborhood of a set of particles
- * are defined in nbsearch.h.
- * The usage is simple: a data structure is allocated with
- * gmx_ana_nbsearch_create(), and the box shape and reference positions for a
- * frame are set using gmx_ana_nbsearch_init() or gmx_ana_nbsearch_pos_init().
- * Searches can then be performed with gmx_ana_nbsearch_is_within() and
- * gmx_ana_nbsearch_mindist(), or with versions that take the \c gmx_ana_pos_t
- * data structure.
- * When the data structure is no longer required, it can be freed with
- * gmx_ana_nbsearch_free().
- *
- * \internal
+/*! \internal \file
+ * \brief
+ * Implements neighborhood searching for analysis (from nbsearch.h).
*
* \todo
* The grid implementation could still be optimized in several different ways:
* - Triclinic grid cells are not the most efficient shape, but make PBC
* handling easier.
- * - Precalculating the required PBC shift for a pair of cells outside the
- * inner loop. After this is done, it should be quite straightforward to
- * move to rectangular cells.
* - Pruning grid cells from the search list if they are completely outside
* the sphere that is being considered.
* - A better heuristic could be added for falling back to simple loops for a
* - A multi-level grid implementation could be used to be able to use small
* grids for short cutoffs with very inhomogeneous particle distributions
* without a memory cost.
- */
-/*! \internal \file
- * \brief
- * Implements functions in nbsearch.h.
*
* \author Teemu Murtola <teemu.murtola@gmail.com>
* \ingroup module_selection
*/
-#ifdef HAVE_CONFIG_H
-#include "config.h"
-#endif
+#include "gmxpre.h"
+
+#include "nbsearch.h"
+
+#include <cmath>
+#include <cstring>
-#include <math.h>
+#include <algorithm>
+#include <vector>
-#include "gromacs/legacyheaders/smalloc.h"
-#include "gromacs/legacyheaders/typedefs.h"
-#include "gromacs/legacyheaders/pbc.h"
-#include "gromacs/legacyheaders/vec.h"
+#include "thread_mpi/mutex.h"
-#include "gromacs/selection/nbsearch.h"
+#include "gromacs/legacyheaders/names.h"
+#include "gromacs/math/vec.h"
+#include "gromacs/pbcutil/pbc.h"
#include "gromacs/selection/position.h"
+#include "gromacs/topology/block.h"
+#include "gromacs/utility/arrayref.h"
+#include "gromacs/utility/exceptions.h"
+#include "gromacs/utility/gmxassert.h"
+#include "gromacs/utility/smalloc.h"
+#include "gromacs/utility/stringutil.h"
-/*! \internal \brief
- * Data structure for neighborhood searches.
+namespace gmx
+{
+
+namespace internal
+{
+
+/********************************************************************
+ * Implementation class declarations
*/
-struct gmx_ana_nbsearch_t
-{
- /** The cutoff. */
- real cutoff;
- /** The cutoff squared. */
- real cutoff2;
- /** Maximum number of reference points. */
- int maxnref;
-
- /** Number of reference points for the current frame. */
- int nref;
- /** Reference point positions. */
- rvec *xref;
- /** Reference position ids (NULL if not available). */
- int *refid;
- /** PBC data. */
- t_pbc *pbc;
-
- /** Number of excluded reference positions for current test particle. */
- int nexcl;
- /** Exclusions for current test particle. */
- int *excl;
-
- /** Whether to try grid searching. */
- bool bTryGrid;
- /** Whether grid searching is actually used for the current positions. */
- bool bGrid;
- /** Array allocated for storing in-unit-cell reference positions. */
- rvec *xref_alloc;
- /** false if the box is rectangular. */
- bool bTric;
- /** Box vectors of a single grid cell. */
- matrix cellbox;
- /** The reciprocal cell vectors as columns; the inverse of \p cellbox. */
- matrix recipcell;
- /** Number of cells along each dimension. */
- ivec ncelldim;
- /** Total number of cells. */
- int ncells;
- /** Number of reference positions in each cell. */
- int *ncatoms;
- /** List of reference positions in each cell. */
- atom_id **catom;
- /** Allocation counts for each \p catom[i]. */
- int *catom_nalloc;
- /** Allocation count for the per-cell arrays. */
- int cells_nalloc;
- /** Number of neighboring cells to consider. */
- int ngridnb;
- /** Offsets of the neighboring cells to consider. */
- ivec *gnboffs;
- /** Allocation count for \p gnboffs. */
- int gnboffs_nalloc;
-
- /** Stores test position during a pair loop. */
- rvec xtest;
- /** Stores the previous returned position during a pair loop. */
- int previ;
- /** Stores the current exclusion index during loops. */
- int exclind;
- /** Stores the test particle cell index during loops. */
- ivec testcell;
- /** Stores the current cell neighbor index during pair loops. */
- int prevnbi;
- /** Stores the index within the current cell during pair loops. */
- int prevcai;
+
+class AnalysisNeighborhoodSearchImpl
+{
+ public:
+ typedef AnalysisNeighborhoodPairSearch::ImplPointer
+ PairSearchImplPointer;
+ typedef std::vector<PairSearchImplPointer> PairSearchList;
+ typedef std::vector<std::vector<int> > CellList;
+
+ explicit AnalysisNeighborhoodSearchImpl(real cutoff);
+ ~AnalysisNeighborhoodSearchImpl();
+
+ /*! \brief
+ * Initializes the search with a given box and reference positions.
+ *
+ * \param[in] mode Search mode to use.
+ * \param[in] bXY Whether to use 2D searching.
+ * \param[in] excls Exclusions.
+ * \param[in] pbc PBC information.
+ * \param[in] positions Set of reference positions.
+ */
+ void init(AnalysisNeighborhood::SearchMode mode,
+ bool bXY,
+ const t_blocka *excls,
+ const t_pbc *pbc,
+ const AnalysisNeighborhoodPositions &positions);
+ PairSearchImplPointer getPairSearch();
+
+ real cutoffSquared() const { return cutoff2_; }
+ bool usesGridSearch() const { return bGrid_; }
+
+ private:
+ //! Calculates offsets to neighboring grid cells that should be considered.
+ void initGridCellNeighborList();
+ /*! \brief
+ * Determines a suitable grid size and sets up the cells.
+ *
+ * \param[in] pbc Information about the box.
+ * \returns false if grid search is not suitable.
+ */
+ bool initGridCells(const t_pbc &pbc);
+ /*! \brief
+ * Sets ua a search grid for a given box.
+ *
+ * \param[in] pbc Information about the box.
+ * \returns false if grid search is not suitable.
+ */
+ bool initGrid(const t_pbc &pbc);
+ /*! \brief
+ * Maps a point into a grid cell.
+ *
+ * \param[in] x Point to map.
+ * \param[out] cell Indices of the grid cell in which \p x lies.
+ * \param[out] xout Coordinates to use
+ * (will be within the triclinic unit cell).
+ */
+ void mapPointToGridCell(const rvec x, ivec cell, rvec xout) const;
+ /*! \brief
+ * Calculates linear index of a grid cell.
+ *
+ * \param[in] cell Cell indices.
+ * \returns Linear index of \p cell.
+ */
+ int getGridCellIndex(const ivec cell) const;
+ /*! \brief
+ * Adds an index into a grid cell.
+ *
+ * \param[in] cell Cell into which \p i should be added.
+ * \param[in] i Index to add.
+ */
+ void addToGridCell(const ivec cell, int i);
+ /*! \brief
+ * Calculates the index of a neighboring grid cell.
+ *
+ * \param[in] sourceCell Location of the source cell.
+ * \param[in] index Index of the neighbor to calculate.
+ * \param[out] shift Shift to apply to get the periodic distance
+ * for distances between the cells.
+ * \returns Grid cell index of the neighboring cell.
+ */
+ int getNeighboringCell(const ivec sourceCell, int index, rvec shift) const;
+
+ //! Whether to try grid searching.
+ bool bTryGrid_;
+ //! The cutoff.
+ real cutoff_;
+ //! The cutoff squared.
+ real cutoff2_;
+ //! Whether to do searching in XY plane only.
+ bool bXY_;
+
+ //! Number of reference points for the current frame.
+ int nref_;
+ //! Reference point positions.
+ const rvec *xref_;
+ //! Reference position exclusion IDs.
+ const int *refExclusionIds_;
+ //! Exclusions.
+ const t_blocka *excls_;
+ //! PBC data.
+ t_pbc pbc_;
+
+ //! Whether grid searching is actually used for the current positions.
+ bool bGrid_;
+ //! Array allocated for storing in-unit-cell reference positions.
+ rvec *xref_alloc_;
+ //! Allocation count for xref_alloc.
+ int xref_nalloc_;
+ //! false if the box is rectangular.
+ bool bTric_;
+ //! Box vectors of a single grid cell.
+ matrix cellbox_;
+ //! The reciprocal cell vectors as columns; the inverse of \p cellbox.
+ matrix recipcell_;
+ //! Number of cells along each dimension.
+ ivec ncelldim_;
+ //! Data structure to hold the grid cell contents.
+ CellList cells_;
+ //! Number of neighboring cells to consider.
+ int ngridnb_;
+ //! Offsets of the neighboring cells to consider.
+ ivec *gnboffs_;
+ //! Allocation count for \p gnboffs.
+ int gnboffs_nalloc_;
+
+ tMPI::mutex createPairSearchMutex_;
+ PairSearchList pairSearchList_;
+
+ friend class AnalysisNeighborhoodPairSearchImpl;
+
+ GMX_DISALLOW_COPY_AND_ASSIGN(AnalysisNeighborhoodSearchImpl);
};
-/*!
- * \param[in] cutoff Cutoff distance for the search
- * (<=0 stands for no cutoff).
- * \param[in] maxn Maximum number of reference particles.
- * \returns Created neighborhood search data structure.
- */
-gmx_ana_nbsearch_t *
-gmx_ana_nbsearch_create(real cutoff, int maxn)
+class AnalysisNeighborhoodPairSearchImpl
{
- gmx_ana_nbsearch_t *d;
+ public:
+ explicit AnalysisNeighborhoodPairSearchImpl(const AnalysisNeighborhoodSearchImpl &search)
+ : search_(search)
+ {
+ testExclusionIds_ = NULL;
+ nexcl_ = 0;
+ excl_ = NULL;
+ clear_rvec(xtest_);
+ clear_ivec(testcell_);
+ reset(-1);
+ }
+
+ //! Initializes a search to find reference positions neighboring \p x.
+ void startSearch(const AnalysisNeighborhoodPositions &positions);
+ //! Searches for the next neighbor.
+ template <class Action>
+ bool searchNext(Action action);
+ //! Initializes a pair representing the pair found by searchNext().
+ void initFoundPair(AnalysisNeighborhoodPair *pair) const;
+ //! Advances to the next test position, skipping any remaining pairs.
+ void nextTestPosition();
+
+ private:
+ //! Clears the loop indices.
+ void reset(int testIndex);
+ //! Checks whether a reference positiong should be excluded.
+ bool isExcluded(int j);
- snew(d, 1);
- d->bTryGrid = true;
- if (cutoff <= 0)
+ //! Parent search object.
+ const AnalysisNeighborhoodSearchImpl &search_;
+ //! Reference to the test positions.
+ ConstArrayRef<rvec> testPositions_;
+ //! Reference to the test exclusion indices.
+ const int *testExclusionIds_;
+ //! Number of excluded reference positions for current test particle.
+ int nexcl_;
+ //! Exclusions for current test particle.
+ const int *excl_;
+ //! Index of the currently active test position in \p testPositions_.
+ int testIndex_;
+ //! Stores test position during a pair loop.
+ rvec xtest_;
+ //! Stores the previous returned position during a pair loop.
+ int previ_;
+ //! Stores the pair distance corresponding to previ_;
+ real prevr2_;
+ //! Stores the current exclusion index during loops.
+ int exclind_;
+ //! Stores the test particle cell index during loops.
+ ivec testcell_;
+ //! Stores the current cell neighbor index during pair loops.
+ int prevnbi_;
+ //! Stores the index within the current cell during pair loops.
+ int prevcai_;
+
+ GMX_DISALLOW_COPY_AND_ASSIGN(AnalysisNeighborhoodPairSearchImpl);
+};
+
+/********************************************************************
+ * AnalysisNeighborhoodSearchImpl
+ */
+
+AnalysisNeighborhoodSearchImpl::AnalysisNeighborhoodSearchImpl(real cutoff)
+{
+ bTryGrid_ = true;
+ cutoff_ = cutoff;
+ if (cutoff_ <= 0)
{
- cutoff = GMX_REAL_MAX;
- d->bTryGrid = false;
+ cutoff_ = GMX_REAL_MAX;
+ bTryGrid_ = false;
}
- d->cutoff = cutoff;
- d->cutoff2 = sqr(cutoff);
- d->maxnref = maxn;
+ cutoff2_ = sqr(cutoff_);
+ bXY_ = false;
- d->xref = NULL;
- d->nexcl = 0;
- d->exclind = 0;
+ nref_ = 0;
+ xref_ = NULL;
+ refExclusionIds_ = NULL;
+ std::memset(&pbc_, 0, sizeof(pbc_));
- d->xref_alloc = NULL;
- d->ncells = 0;
- d->ncatoms = NULL;
- d->catom = NULL;
- d->catom_nalloc = 0;
- d->cells_nalloc = 0;
+ bGrid_ = false;
- d->ngridnb = 0;
- d->gnboffs = NULL;
- d->gnboffs_nalloc = 0;
+ xref_alloc_ = NULL;
+ xref_nalloc_ = 0;
+ bTric_ = false;
+ clear_mat(cellbox_);
+ clear_mat(recipcell_);
+ clear_ivec(ncelldim_);
- return d;
+ ngridnb_ = 0;
+ gnboffs_ = NULL;
+ gnboffs_nalloc_ = 0;
}
-/*!
- * \param d Data structure to free.
- *
- * After the call, the pointer \p d is no longer valid.
- */
-void
-gmx_ana_nbsearch_free(gmx_ana_nbsearch_t *d)
+AnalysisNeighborhoodSearchImpl::~AnalysisNeighborhoodSearchImpl()
{
- sfree(d->xref_alloc);
- sfree(d->ncatoms);
- if (d->catom)
+ PairSearchList::const_iterator i;
+ for (i = pairSearchList_.begin(); i != pairSearchList_.end(); ++i)
{
- int ci;
+ GMX_RELEASE_ASSERT(i->unique(),
+ "Dangling AnalysisNeighborhoodPairSearch reference");
+ }
+ sfree(xref_alloc_);
+ sfree(gnboffs_);
+}
- for (ci = 0; ci < d->ncells; ++ci)
+AnalysisNeighborhoodSearchImpl::PairSearchImplPointer
+AnalysisNeighborhoodSearchImpl::getPairSearch()
+{
+ tMPI::lock_guard<tMPI::mutex> lock(createPairSearchMutex_);
+ // TODO: Consider whether this needs to/can be faster, e.g., by keeping a
+ // separate pool of unused search objects.
+ PairSearchList::const_iterator i;
+ for (i = pairSearchList_.begin(); i != pairSearchList_.end(); ++i)
+ {
+ if (i->unique())
{
- sfree(d->catom[ci]);
+ return *i;
}
- sfree(d->catom);
}
- sfree(d->catom_nalloc);
- sfree(d->gnboffs);
- sfree(d);
+ PairSearchImplPointer pairSearch(new AnalysisNeighborhoodPairSearchImpl(*this));
+ pairSearchList_.push_back(pairSearch);
+ return pairSearch;
}
-/*! \brief
- * Calculates offsets to neighboring grid cells that should be considered.
- *
- * \param[in,out] d Grid information.
- * \param[in] pbc Information about the box.
- */
-static void
-grid_init_cell_nblist(gmx_ana_nbsearch_t *d, t_pbc *pbc)
+void AnalysisNeighborhoodSearchImpl::initGridCellNeighborList()
{
int maxx, maxy, maxz;
- int x, y, z, i;
real rvnorm;
/* Find the extent of the sphere in triclinic coordinates */
- maxz = (int)(d->cutoff * d->recipcell[ZZ][ZZ]) + 1;
- rvnorm = sqrt(sqr(d->recipcell[YY][YY]) + sqr(d->recipcell[ZZ][YY]));
- maxy = (int)(d->cutoff * rvnorm) + 1;
- rvnorm = sqrt(sqr(d->recipcell[XX][XX]) + sqr(d->recipcell[YY][XX])
- + sqr(d->recipcell[ZZ][XX]));
- maxx = (int)(d->cutoff * rvnorm) + 1;
+ maxz = static_cast<int>(cutoff_ * recipcell_[ZZ][ZZ]) + 1;
+ rvnorm = sqrt(sqr(recipcell_[YY][YY]) + sqr(recipcell_[ZZ][YY]));
+ maxy = static_cast<int>(cutoff_ * rvnorm) + 1;
+ rvnorm = sqrt(sqr(recipcell_[XX][XX]) + sqr(recipcell_[YY][XX])
+ + sqr(recipcell_[ZZ][XX]));
+ maxx = static_cast<int>(cutoff_ * rvnorm) + 1;
/* Calculate the number of cells and reallocate if necessary */
- d->ngridnb = (2 * maxx + 1) * (2 * maxy + 1) * (2 * maxz + 1);
- if (d->gnboffs_nalloc < d->ngridnb)
+ ngridnb_ = (2 * maxx + 1) * (2 * maxy + 1) * (2 * maxz + 1);
+ if (gnboffs_nalloc_ < ngridnb_)
{
- d->gnboffs_nalloc = d->ngridnb;
- srenew(d->gnboffs, d->gnboffs_nalloc);
+ gnboffs_nalloc_ = ngridnb_;
+ srenew(gnboffs_, gnboffs_nalloc_);
}
/* Store the whole cube */
/* TODO: Prune off corners that are not needed */
- i = 0;
- for (x = -maxx; x <= maxx; ++x)
+ int i = 0;
+ for (int x = -maxx; x <= maxx; ++x)
{
- for (y = -maxy; y <= maxy; ++y)
+ for (int y = -maxy; y <= maxy; ++y)
{
- for (z = -maxz; z <= maxz; ++z)
+ for (int z = -maxz; z <= maxz; ++z)
{
- d->gnboffs[i][XX] = x;
- d->gnboffs[i][YY] = y;
- d->gnboffs[i][ZZ] = z;
+ gnboffs_[i][XX] = x;
+ gnboffs_[i][YY] = y;
+ gnboffs_[i][ZZ] = z;
++i;
}
}
}
}
-/*! \brief
- * Determines a suitable grid size.
- *
- * \param[in,out] d Grid information.
- * \param[in] pbc Information about the box.
- * \returns false if grid search is not suitable.
- */
-static bool
-grid_setup_cells(gmx_ana_nbsearch_t *d, t_pbc *pbc)
-{
- real targetsize;
- int dd;
-
-#ifdef HAVE_CBRT
- targetsize = cbrt(pbc->box[XX][XX] * pbc->box[YY][YY] * pbc->box[ZZ][ZZ]
- * 10 / d->nref);
-#else
- targetsize = pow(pbc->box[XX][XX] * pbc->box[YY][YY] * pbc->box[ZZ][ZZ]
- * 10 / d->nref, (real)(1./3.));
-#endif
-
- d->ncells = 1;
- for (dd = 0; dd < DIM; ++dd)
- {
- d->ncelldim[dd] = (int)(pbc->box[dd][dd] / targetsize);
- d->ncells *= d->ncelldim[dd];
- if (d->ncelldim[dd] < 3)
+bool AnalysisNeighborhoodSearchImpl::initGridCells(const t_pbc &pbc)
+{
+ const real targetsize =
+ pow(pbc.box[XX][XX] * pbc.box[YY][YY] * pbc.box[ZZ][ZZ]
+ * 10 / nref_, static_cast<real>(1./3.));
+
+ int cellCount = 1;
+ for (int dd = 0; dd < DIM; ++dd)
+ {
+ ncelldim_[dd] = static_cast<int>(pbc.box[dd][dd] / targetsize);
+ cellCount *= ncelldim_[dd];
+ if (ncelldim_[dd] < 3)
{
return false;
}
}
- /* Reallocate if necessary */
- if (d->cells_nalloc < d->ncells)
+ // Never decrease the size of the cell vector to avoid reallocating
+ // memory for the nested vectors. The actual size of the vector is not
+ // used outside this function.
+ if (cells_.size() < static_cast<size_t>(cellCount))
{
- int i;
-
- srenew(d->ncatoms, d->ncells);
- srenew(d->catom, d->ncells);
- srenew(d->catom_nalloc, d->ncells);
- for (i = d->cells_nalloc; i < d->ncells; ++i)
- {
- d->catom[i] = NULL;
- d->catom_nalloc[i] = 0;
- }
- d->cells_nalloc = d->ncells;
+ cells_.resize(cellCount);
+ }
+ for (int ci = 0; ci < cellCount; ++ci)
+ {
+ cells_[ci].clear();
}
return true;
}
-/*! \brief
- * Sets ua a search grid for a given box.
- *
- * \param[in,out] d Grid information.
- * \param[in] pbc Information about the box.
- * \returns false if grid search is not suitable.
- */
-static bool
-grid_set_box(gmx_ana_nbsearch_t *d, t_pbc *pbc)
+bool AnalysisNeighborhoodSearchImpl::initGrid(const t_pbc &pbc)
{
- int dd;
-
/* TODO: This check could be improved. */
- if (0.5*pbc->max_cutoff2 < d->cutoff2)
+ if (0.5*pbc.max_cutoff2 < cutoff2_)
{
return false;
}
- if (!grid_setup_cells(d, pbc))
+ if (!initGridCells(pbc))
{
return false;
}
- d->bTric = TRICLINIC(pbc->box);
- if (d->bTric)
+ bTric_ = TRICLINIC(pbc.box);
+ if (bTric_)
{
- for (dd = 0; dd < DIM; ++dd)
+ for (int dd = 0; dd < DIM; ++dd)
{
- svmul(1.0 / d->ncelldim[dd], pbc->box[dd], d->cellbox[dd]);
+ svmul(1.0 / ncelldim_[dd], pbc.box[dd], cellbox_[dd]);
}
- m_inv_ur0(d->cellbox, d->recipcell);
+ m_inv_ur0(cellbox_, recipcell_);
}
else
{
- for (dd = 0; dd < DIM; ++dd)
+ for (int dd = 0; dd < DIM; ++dd)
{
- d->cellbox[dd][dd] = pbc->box[dd][dd] / d->ncelldim[dd];
- d->recipcell[dd][dd] = 1 / d->cellbox[dd][dd];
+ cellbox_[dd][dd] = pbc.box[dd][dd] / ncelldim_[dd];
+ recipcell_[dd][dd] = 1.0 / cellbox_[dd][dd];
}
}
- grid_init_cell_nblist(d, pbc);
+ initGridCellNeighborList();
return true;
}
-/*! \brief
- * Maps a point into a grid cell.
- *
- * \param[in] d Grid information.
- * \param[in] x Point to map.
- * \param[out] cell Indices of the grid cell in which \p x lies.
- *
- * \p x should be in the triclinic unit cell.
- */
-static void
-grid_map_onto(gmx_ana_nbsearch_t *d, const rvec x, ivec cell)
+void AnalysisNeighborhoodSearchImpl::mapPointToGridCell(const rvec x,
+ ivec cell,
+ rvec xout) const
{
- int dd;
-
- if (d->bTric)
+ rvec xtmp;
+ copy_rvec(x, xtmp);
+ if (bTric_)
{
rvec tx;
-
- tmvmul_ur0(d->recipcell, x, tx);
- for (dd = 0; dd < DIM; ++dd)
+ tmvmul_ur0(recipcell_, xtmp, tx);
+ for (int dd = 0; dd < DIM; ++dd)
{
- cell[dd] = (int)tx[dd];
+ const int cellCount = ncelldim_[dd];
+ int cellIndex = static_cast<int>(floor(tx[dd]));
+ while (cellIndex < 0)
+ {
+ cellIndex += cellCount;
+ rvec_add(xtmp, pbc_.box[dd], xtmp);
+ }
+ while (cellIndex >= cellCount)
+ {
+ cellIndex -= cellCount;
+ rvec_sub(xtmp, pbc_.box[dd], xtmp);
+ }
+ cell[dd] = cellIndex;
}
}
else
{
- for (dd = 0; dd < DIM; ++dd)
+ for (int dd = 0; dd < DIM; ++dd)
{
- cell[dd] = (int)(x[dd] * d->recipcell[dd][dd]);
+ const int cellCount = ncelldim_[dd];
+ int cellIndex = static_cast<int>(floor(xtmp[dd] * recipcell_[dd][dd]));
+ while (cellIndex < 0)
+ {
+ cellIndex += cellCount;
+ xtmp[dd] += pbc_.box[dd][dd];
+ }
+ while (cellIndex >= cellCount)
+ {
+ cellIndex -= cellCount;
+ xtmp[dd] -= pbc_.box[dd][dd];
+ }
+ cell[dd] = cellIndex;
}
}
+ copy_rvec(xtmp, xout);
}
-/*! \brief
- * Calculates linear index of a grid cell.
- *
- * \param[in] d Grid information.
- * \param[in] cell Cell indices.
- * \returns Linear index of \p cell.
- */
-static int
-grid_index(gmx_ana_nbsearch_t *d, const ivec cell)
+int AnalysisNeighborhoodSearchImpl::getGridCellIndex(const ivec cell) const
{
- return cell[XX] + cell[YY] * d->ncelldim[XX]
- + cell[ZZ] * d->ncelldim[XX] * d->ncelldim[YY];
+ GMX_ASSERT(cell[XX] >= 0 && cell[XX] < ncelldim_[XX],
+ "Grid cell X index out of range");
+ GMX_ASSERT(cell[YY] >= 0 && cell[YY] < ncelldim_[YY],
+ "Grid cell Y index out of range");
+ GMX_ASSERT(cell[ZZ] >= 0 && cell[ZZ] < ncelldim_[ZZ],
+ "Grid cell Z index out of range");
+ return cell[XX] + cell[YY] * ncelldim_[XX]
+ + cell[ZZ] * ncelldim_[XX] * ncelldim_[YY];
}
-/*! \brief
- * Clears all grid cells.
- *
- * \param[in,out] d Grid information.
- */
-static void
-grid_clear_cells(gmx_ana_nbsearch_t *d)
+void AnalysisNeighborhoodSearchImpl::addToGridCell(const ivec cell, int i)
{
- int ci;
-
- for (ci = 0; ci < d->ncells; ++ci)
- {
- d->ncatoms[ci] = 0;
- }
+ const int ci = getGridCellIndex(cell);
+ cells_[ci].push_back(i);
}
-/*! \brief
- * Adds an index into a grid cell.
- *
- * \param[in,out] d Grid information.
- * \param[in] cell Cell into which \p i should be added.
- * \param[in] i Index to add.
- */
-static void
-grid_add_to_cell(gmx_ana_nbsearch_t *d, const ivec cell, int i)
+int AnalysisNeighborhoodSearchImpl::getNeighboringCell(
+ const ivec sourceCell, int index, rvec shift) const
{
- int ci = grid_index(d, cell);
+ ivec cell;
+ ivec_add(sourceCell, gnboffs_[index], cell);
- if (d->ncatoms[ci] == d->catom_nalloc[ci])
+ // TODO: Consider unifying with the similar shifting code in
+ // mapPointToGridCell().
+ clear_rvec(shift);
+ for (int d = 0; d < DIM; ++d)
{
- d->catom_nalloc[ci] += 10;
- srenew(d->catom[ci], d->catom_nalloc[ci]);
+ const int cellCount = ncelldim_[d];
+ if (cell[d] < 0)
+ {
+ cell[d] += cellCount;
+ rvec_add(shift, pbc_.box[d], shift);
+ }
+ else if (cell[d] >= cellCount)
+ {
+ cell[d] -= cellCount;
+ rvec_sub(shift, pbc_.box[d], shift);
+ }
}
- d->catom[ci][d->ncatoms[ci]++] = i;
+
+ return getGridCellIndex(cell);
}
-/*!
- * \param[in,out] d Neighborhood search data structure.
- * \param[in] pbc PBC information for the frame.
- * \param[in] n Number of reference positions for the frame.
- * \param[in] x \p n reference positions for the frame.
- *
- * Initializes the data structure \p d such that it can be used to search
- * for the neighbors of \p x.
- */
-void
-gmx_ana_nbsearch_init(gmx_ana_nbsearch_t *d, t_pbc *pbc, int n, const rvec x[])
+void AnalysisNeighborhoodSearchImpl::init(
+ AnalysisNeighborhood::SearchMode mode,
+ bool bXY,
+ const t_blocka *excls,
+ const t_pbc *pbc,
+ const AnalysisNeighborhoodPositions &positions)
{
- d->pbc = pbc;
- d->nref = n;
- if (!pbc)
+ GMX_RELEASE_ASSERT(positions.index_ == -1,
+ "Individual indexed positions not supported as reference");
+ bXY_ = bXY;
+ if (bXY_ && pbc->ePBC != epbcNONE)
{
- d->bGrid = false;
+ if (pbc->ePBC != epbcXY && pbc->ePBC != epbcXYZ)
+ {
+ std::string message =
+ formatString("Computations in the XY plane are not supported with PBC type '%s'",
+ EPBC(pbc->ePBC));
+ GMX_THROW(NotImplementedError(message));
+ }
+ if (std::fabs(pbc->box[ZZ][XX]) > GMX_REAL_EPS*pbc->box[ZZ][ZZ] ||
+ std::fabs(pbc->box[ZZ][YY]) > GMX_REAL_EPS*pbc->box[ZZ][ZZ])
+ {
+ GMX_THROW(NotImplementedError("Computations in the XY plane are not supported when the last box vector is not parallel to the Z axis"));
+ }
+ set_pbc(&pbc_, epbcXY, const_cast<rvec *>(pbc->box));
}
- else if (d->bTryGrid)
+ else if (pbc != NULL)
{
- d->bGrid = grid_set_box(d, pbc);
+ pbc_ = *pbc;
}
- if (d->bGrid)
+ else
{
- int i;
-
- if (!d->xref_alloc)
+ pbc_.ePBC = epbcNONE;
+ }
+ nref_ = positions.count_;
+ // TODO: Consider whether it would be possible to support grid searching in
+ // more cases.
+ if (mode == AnalysisNeighborhood::eSearchMode_Simple
+ || pbc_.ePBC != epbcXYZ)
+ {
+ bGrid_ = false;
+ }
+ else if (bTryGrid_)
+ {
+ // TODO: Actually implement forcing eSearchMode_Grid
+ bGrid_ = initGrid(pbc_);
+ }
+ if (bGrid_)
+ {
+ if (xref_nalloc_ < nref_)
{
- snew(d->xref_alloc, d->maxnref);
+ srenew(xref_alloc_, nref_);
+ xref_nalloc_ = nref_;
}
- d->xref = d->xref_alloc;
- grid_clear_cells(d);
+ xref_ = xref_alloc_;
- for (i = 0; i < n; ++i)
- {
- copy_rvec(x[i], d->xref[i]);
- }
- put_atoms_in_triclinic_unitcell(ecenterTRIC, pbc->box, n, d->xref);
- for (i = 0; i < n; ++i)
+ for (int i = 0; i < nref_; ++i)
{
ivec refcell;
-
- grid_map_onto(d, d->xref[i], refcell);
- grid_add_to_cell(d, refcell, i);
+ mapPointToGridCell(positions.x_[i], refcell, xref_alloc_[i]);
+ addToGridCell(refcell, i);
}
}
else
{
- // Won't be modified in this case, but when a grid is used,
- // xref _is_ modified, so it can't be const.
- d->xref = const_cast<rvec *>(x);
+ xref_ = positions.x_;
+ }
+ excls_ = excls;
+ refExclusionIds_ = NULL;
+ if (excls != NULL)
+ {
+ // TODO: Check that the IDs are ascending, or remove the limitation.
+ refExclusionIds_ = positions.exclusionIds_;
+ GMX_RELEASE_ASSERT(refExclusionIds_ != NULL,
+ "Exclusion IDs must be set for reference positions "
+ "when exclusions are enabled");
}
- d->refid = NULL;
}
-/*!
- * \param[in,out] d Neighborhood search data structure.
- * \param[in] pbc PBC information for the frame.
- * \param[in] p Reference positions for the frame.
- *
- * A convenience wrapper for gmx_ana_nbsearch_init().
+/********************************************************************
+ * AnalysisNeighborhoodPairSearchImpl
*/
-void
-gmx_ana_nbsearch_pos_init(gmx_ana_nbsearch_t *d, t_pbc *pbc, const gmx_ana_pos_t *p)
-{
- gmx_ana_nbsearch_init(d, pbc, p->nr, p->x);
- d->refid = (p->nr < d->maxnref ? p->m.refid : NULL);
-}
-/*!
- * \param[in,out] d Neighborhood search data structure.
- * \param[in] nexcl Number of reference positions to exclude from next
- * search.
- * \param[in] excl Indices of reference positions to exclude.
- *
- * The set exclusions remain in effect until the next call of this function.
- */
-void
-gmx_ana_nbsearch_set_excl(gmx_ana_nbsearch_t *d, int nexcl, int excl[])
+void AnalysisNeighborhoodPairSearchImpl::reset(int testIndex)
{
-
- d->nexcl = nexcl;
- d->excl = excl;
-}
-
-/*! \brief
- * Helper function to check whether a reference point should be excluded.
- */
-static bool
-is_excluded(gmx_ana_nbsearch_t *d, int j)
-{
- if (d->exclind < d->nexcl)
+ testIndex_ = testIndex;
+ if (testIndex_ >= 0 && testIndex_ < static_cast<int>(testPositions_.size()))
{
- if (d->refid)
+ copy_rvec(testPositions_[testIndex_], xtest_);
+ if (search_.bGrid_)
{
- while (d->exclind < d->nexcl &&d->refid[j] > d->excl[d->exclind])
- {
- ++d->exclind;
- }
- if (d->exclind < d->nexcl && d->refid[j] == d->excl[d->exclind])
- {
- ++d->exclind;
- return true;
- }
+ search_.mapPointToGridCell(testPositions_[testIndex], testcell_, xtest_);
}
else
{
- while (d->bGrid && d->exclind < d->nexcl && d->excl[d->exclind] < j)
+ copy_rvec(testPositions_[testIndex_], xtest_);
+ }
+ if (search_.excls_ != NULL)
+ {
+ const int exclIndex = testExclusionIds_[testIndex];
+ if (exclIndex < search_.excls_->nr)
{
- ++d->exclind;
+ const int startIndex = search_.excls_->index[exclIndex];
+ nexcl_ = search_.excls_->index[exclIndex + 1] - startIndex;
+ excl_ = &search_.excls_->a[startIndex];
}
- if (d->excl[d->exclind] == j)
+ else
{
- ++d->exclind;
- return true;
+ nexcl_ = 0;
+ excl_ = NULL;
}
}
}
+ previ_ = -1;
+ prevr2_ = 0.0;
+ exclind_ = 0;
+ prevnbi_ = 0;
+ prevcai_ = -1;
+}
+
+void AnalysisNeighborhoodPairSearchImpl::nextTestPosition()
+{
+ if (testIndex_ < static_cast<int>(testPositions_.size()))
+ {
+ ++testIndex_;
+ reset(testIndex_);
+ }
+}
+
+bool AnalysisNeighborhoodPairSearchImpl::isExcluded(int j)
+{
+ if (exclind_ < nexcl_)
+ {
+ const int refId = search_.refExclusionIds_[j];
+ while (exclind_ < nexcl_ && excl_[exclind_] < refId)
+ {
+ ++exclind_;
+ }
+ if (exclind_ < nexcl_ && refId == excl_[exclind_])
+ {
+ ++exclind_;
+ return true;
+ }
+ }
return false;
}
-/*! \brief
- * Initializes a grid search to find reference positions neighboring \p x.
- */
-static void
-grid_search_start(gmx_ana_nbsearch_t *d, const rvec x)
+void AnalysisNeighborhoodPairSearchImpl::startSearch(
+ const AnalysisNeighborhoodPositions &positions)
{
- copy_rvec(x, d->xtest);
- if (d->bGrid)
+ testExclusionIds_ = positions.exclusionIds_;
+ GMX_RELEASE_ASSERT(search_.excls_ == NULL || testExclusionIds_ != NULL,
+ "Exclusion IDs must be set when exclusions are enabled");
+ if (positions.index_ < 0)
{
- put_atoms_in_triclinic_unitcell(ecenterTRIC, d->pbc->box, 1, &d->xtest);
- grid_map_onto(d, d->xtest, d->testcell);
- d->prevnbi = 0;
- d->prevcai = -1;
+ testPositions_ = constArrayRefFromArray<rvec>(positions.x_, positions.count_);
+ reset(0);
}
else
{
- d->previ = -1;
+ // Somewhat of a hack: setup the array such that only the last position
+ // will be used.
+ testPositions_ = constArrayRefFromArray<rvec>(positions.x_, positions.index_ + 1);
+ reset(positions.index_);
}
- d->exclind = 0;
}
-/*! \brief
- * Does a grid search.
- */
-static bool
-grid_search(gmx_ana_nbsearch_t *d,
- bool (*action)(gmx_ana_nbsearch_t *d, int i, real r2))
+template <class Action>
+bool AnalysisNeighborhoodPairSearchImpl::searchNext(Action action)
{
- int i;
- rvec dx;
- real r2;
-
- if (d->bGrid)
+ while (testIndex_ < static_cast<int>(testPositions_.size()))
{
- int nbi, ci, cai;
+ if (search_.bGrid_)
+ {
+ GMX_RELEASE_ASSERT(!search_.bXY_, "Grid-based XY searches not implemented");
- nbi = d->prevnbi;
- cai = d->prevcai + 1;
+ int nbi = prevnbi_;
+ int cai = prevcai_ + 1;
- for (; nbi < d->ngridnb; ++nbi)
+ for (; nbi < search_.ngridnb_; ++nbi)
+ {
+ rvec shift;
+ const int ci = search_.getNeighboringCell(testcell_, nbi, shift);
+ const int cellSize = static_cast<int>(search_.cells_[ci].size());
+ for (; cai < cellSize; ++cai)
+ {
+ const int i = search_.cells_[ci][cai];
+ if (isExcluded(i))
+ {
+ continue;
+ }
+ rvec dx;
+ rvec_sub(xtest_, search_.xref_[i], dx);
+ rvec_add(dx, shift, dx);
+ const real r2 = norm2(dx);
+ if (r2 <= search_.cutoff2_)
+ {
+ if (action(i, r2))
+ {
+ prevnbi_ = nbi;
+ prevcai_ = cai;
+ previ_ = i;
+ prevr2_ = r2;
+ return true;
+ }
+ }
+ }
+ exclind_ = 0;
+ cai = 0;
+ }
+ }
+ else
{
- ivec cell;
-
- ivec_add(d->testcell, d->gnboffs[nbi], cell);
- /* TODO: Support for 2D and screw PBC */
- cell[XX] = (cell[XX] + d->ncelldim[XX]) % d->ncelldim[XX];
- cell[YY] = (cell[YY] + d->ncelldim[YY]) % d->ncelldim[YY];
- cell[ZZ] = (cell[ZZ] + d->ncelldim[ZZ]) % d->ncelldim[ZZ];
- ci = grid_index(d, cell);
- /* TODO: Calculate the required PBC shift outside the inner loop */
- for (; cai < d->ncatoms[ci]; ++cai)
+ for (int i = previ_ + 1; i < search_.nref_; ++i)
{
- i = d->catom[ci][cai];
- if (is_excluded(d, i))
+ if (isExcluded(i))
{
continue;
}
- pbc_dx_aiuc(d->pbc, d->xtest, d->xref[i], dx);
- r2 = norm2(dx);
- if (r2 <= d->cutoff2)
+ rvec dx;
+ if (search_.pbc_.ePBC != epbcNONE)
+ {
+ pbc_dx(&search_.pbc_, xtest_, search_.xref_[i], dx);
+ }
+ else
{
- if (action(d, i, r2))
+ rvec_sub(xtest_, search_.xref_[i], dx);
+ }
+ const real r2
+ = search_.bXY_
+ ? dx[XX]*dx[XX] + dx[YY]*dx[YY]
+ : norm2(dx);
+ if (r2 <= search_.cutoff2_)
+ {
+ if (action(i, r2))
{
- d->prevnbi = nbi;
- d->prevcai = cai;
- d->previ = i;
+ previ_ = i;
+ prevr2_ = r2;
return true;
}
}
}
- d->exclind = 0;
- cai = 0;
}
+ nextTestPosition();
+ }
+ return false;
+}
+
+void AnalysisNeighborhoodPairSearchImpl::initFoundPair(
+ AnalysisNeighborhoodPair *pair) const
+{
+ if (previ_ < 0)
+ {
+ *pair = AnalysisNeighborhoodPair();
}
else
{
- i = d->previ + 1;
- for (; i < d->nref; ++i)
- {
- if (is_excluded(d, i))
- {
- continue;
- }
- if (d->pbc)
- {
- pbc_dx(d->pbc, d->xtest, d->xref[i], dx);
- }
- else
- {
- rvec_sub(d->xtest, d->xref[i], dx);
- }
- r2 = norm2(dx);
- if (r2 <= d->cutoff2)
- {
- if (action(d, i, r2))
- {
- d->previ = i;
- return true;
- }
- }
- }
+ *pair = AnalysisNeighborhoodPair(previ_, testIndex_, prevr2_);
}
- return false;
}
+} // namespace internal
+
+namespace
+{
+
/*! \brief
- * Helper function to use with grid_search() to find the next neighbor.
+ * Search action to find the next neighbor.
+ *
+ * Used as the action for AnalysisNeighborhoodPairSearchImpl::searchNext() to
+ * find the next neighbor.
*
* Simply breaks the loop on the first found neighbor.
*/
-static bool
-within_action(gmx_ana_nbsearch_t *d, int i, real r2)
+bool withinAction(int /*i*/, real /*r2*/)
{
return true;
}
/*! \brief
- * Helper function to use with grid_search() to find the minimum distance.
+ * Search action find the minimum distance.
+ *
+ * Used as the action for AnalysisNeighborhoodPairSearchImpl::searchNext() to
+ * find the nearest neighbor.
+ *
+ * With this action, AnalysisNeighborhoodPairSearchImpl::searchNext() always
+ * returns false, and the output is put into the variables passed by pointer
+ * into the constructor. If no neighbors are found, the output variables are
+ * not modified, i.e., the caller must initialize them.
*/
-static bool
-mindist_action(gmx_ana_nbsearch_t *d, int i, real r2)
+class MindistAction
{
- d->cutoff2 = r2;
- return false;
-}
+ public:
+ /*! \brief
+ * Initializes the action with given output locations.
+ *
+ * \param[out] closestPoint Index of the closest reference location.
+ * \param[out] minDist2 Minimum distance squared.
+ *
+ * The constructor call does not modify the pointed values, but only
+ * stores the pointers for later use.
+ * See the class description for additional semantics.
+ */
+ MindistAction(int *closestPoint, real *minDist2)
+ : closestPoint_(*closestPoint), minDist2_(*minDist2)
+ {
+ }
+
+ //! Processes a neighbor to find the nearest point.
+ bool operator()(int i, real r2)
+ {
+ if (r2 < minDist2_)
+ {
+ closestPoint_ = i;
+ minDist2_ = r2;
+ }
+ return false;
+ }
+
+ private:
+ int &closestPoint_;
+ real &minDist2_;
+
+ GMX_DISALLOW_ASSIGN(MindistAction);
+};
-/*!
- * \param[in] d Neighborhood search data structure.
- * \param[in] x Test position.
- * \returns true if \p x is within the cutoff of any reference position,
- * false otherwise.
+} // namespace
+
+/********************************************************************
+ * AnalysisNeighborhood::Impl
*/
-bool
-gmx_ana_nbsearch_is_within(gmx_ana_nbsearch_t *d, const rvec x)
+
+class AnalysisNeighborhood::Impl
{
- grid_search_start(d, x);
- return grid_search(d, &within_action);
+ public:
+ typedef AnalysisNeighborhoodSearch::ImplPointer SearchImplPointer;
+ typedef std::vector<SearchImplPointer> SearchList;
+
+ Impl() : cutoff_(0), excls_(NULL), mode_(eSearchMode_Automatic), bXY_(false)
+ {
+ }
+ ~Impl()
+ {
+ SearchList::const_iterator i;
+ for (i = searchList_.begin(); i != searchList_.end(); ++i)
+ {
+ GMX_RELEASE_ASSERT(i->unique(),
+ "Dangling AnalysisNeighborhoodSearch reference");
+ }
+ }
+
+ SearchImplPointer getSearch();
+
+ tMPI::mutex createSearchMutex_;
+ SearchList searchList_;
+ real cutoff_;
+ const t_blocka *excls_;
+ SearchMode mode_;
+ bool bXY_;
+};
+
+AnalysisNeighborhood::Impl::SearchImplPointer
+AnalysisNeighborhood::Impl::getSearch()
+{
+ tMPI::lock_guard<tMPI::mutex> lock(createSearchMutex_);
+ // TODO: Consider whether this needs to/can be faster, e.g., by keeping a
+ // separate pool of unused search objects.
+ SearchList::const_iterator i;
+ for (i = searchList_.begin(); i != searchList_.end(); ++i)
+ {
+ if (i->unique())
+ {
+ return *i;
+ }
+ }
+ SearchImplPointer search(new internal::AnalysisNeighborhoodSearchImpl(cutoff_));
+ searchList_.push_back(search);
+ return search;
}
-/*!
- * \param[in] d Neighborhood search data structure.
- * \param[in] p Test positions.
- * \param[in] i Use the i'th position in \p p for testing.
- * \returns true if the test position is within the cutoff of any reference
- * position, false otherwise.
+/********************************************************************
+ * AnalysisNeighborhood
*/
-bool
-gmx_ana_nbsearch_pos_is_within(gmx_ana_nbsearch_t *d, const gmx_ana_pos_t *p, int i)
+
+AnalysisNeighborhood::AnalysisNeighborhood()
+ : impl_(new Impl)
{
- return gmx_ana_nbsearch_is_within(d, p->x[i]);
}
-/*!
- * \param[in] d Neighborhood search data structure.
- * \param[in] x Test position.
- * \returns The distance to the nearest reference position, or the cutoff
- * value if there are no reference positions within the cutoff.
- */
-real
-gmx_ana_nbsearch_mindist(gmx_ana_nbsearch_t *d, const rvec x)
+AnalysisNeighborhood::~AnalysisNeighborhood()
{
- real mind;
+}
- grid_search_start(d, x);
- grid_search(d, &mindist_action);
- mind = sqrt(d->cutoff2);
- d->cutoff2 = sqr(d->cutoff);
- return mind;
+void AnalysisNeighborhood::setCutoff(real cutoff)
+{
+ GMX_RELEASE_ASSERT(impl_->searchList_.empty(),
+ "Changing the cutoff after initSearch() not currently supported");
+ impl_->cutoff_ = cutoff;
}
-/*!
- * \param[in] d Neighborhood search data structure.
- * \param[in] p Test positions.
- * \param[in] i Use the i'th position in \p p for testing.
- * \returns The distance to the nearest reference position, or the cutoff
- * value if there are no reference positions within the cutoff.
- */
-real
-gmx_ana_nbsearch_pos_mindist(gmx_ana_nbsearch_t *d, const gmx_ana_pos_t *p, int i)
+void AnalysisNeighborhood::setXYMode(bool bXY)
{
- return gmx_ana_nbsearch_mindist(d, p->x[i]);
+ impl_->bXY_ = bXY;
}
-/*!
- * \param[in] d Neighborhood search data structure.
- * \param[in] x Test positions.
- * \param[out] jp Index of the reference position in the first pair.
- * \returns true if there are positions within the cutoff.
- */
-bool
-gmx_ana_nbsearch_first_within(gmx_ana_nbsearch_t *d, const rvec x, int *jp)
+void AnalysisNeighborhood::setTopologyExclusions(const t_blocka *excls)
+{
+ GMX_RELEASE_ASSERT(impl_->searchList_.empty(),
+ "Changing the exclusions after initSearch() not currently supported");
+ impl_->excls_ = excls;
+}
+
+void AnalysisNeighborhood::setMode(SearchMode mode)
+{
+ impl_->mode_ = mode;
+}
+
+AnalysisNeighborhood::SearchMode AnalysisNeighborhood::mode() const
+{
+ return impl_->mode_;
+}
+
+AnalysisNeighborhoodSearch
+AnalysisNeighborhood::initSearch(const t_pbc *pbc,
+ const AnalysisNeighborhoodPositions &positions)
{
- grid_search_start(d, x);
- return gmx_ana_nbsearch_next_within(d, jp);
+ Impl::SearchImplPointer search(impl_->getSearch());
+ search->init(mode(), impl_->bXY_, impl_->excls_, pbc, positions);
+ return AnalysisNeighborhoodSearch(search);
}
-/*!
- * \param[in] d Neighborhood search data structure.
- * \param[in] p Test positions.
- * \param[in] i Use the i'th position in \p p.
- * \param[out] jp Index of the reference position in the first pair.
- * \returns true if there are positions within the cutoff.
+/********************************************************************
+ * AnalysisNeighborhoodSearch
*/
-bool
-gmx_ana_nbsearch_pos_first_within(gmx_ana_nbsearch_t *d, const gmx_ana_pos_t *p,
- int i, int *jp)
+
+AnalysisNeighborhoodSearch::AnalysisNeighborhoodSearch()
{
- return gmx_ana_nbsearch_first_within(d, p->x[i], jp);
}
-/*!
- * \param[in] d Neighborhood search data structure.
- * \param[out] jp Index of the test position in the next pair.
- * \returns true if there are positions within the cutoff.
+AnalysisNeighborhoodSearch::AnalysisNeighborhoodSearch(const ImplPointer &impl)
+ : impl_(impl)
+{
+}
+
+void AnalysisNeighborhoodSearch::reset()
+{
+ impl_.reset();
+}
+
+AnalysisNeighborhood::SearchMode AnalysisNeighborhoodSearch::mode() const
+{
+ GMX_RELEASE_ASSERT(impl_, "Accessing an invalid search object");
+ return (impl_->usesGridSearch()
+ ? AnalysisNeighborhood::eSearchMode_Grid
+ : AnalysisNeighborhood::eSearchMode_Simple);
+}
+
+bool AnalysisNeighborhoodSearch::isWithin(
+ const AnalysisNeighborhoodPositions &positions) const
+{
+ GMX_RELEASE_ASSERT(impl_, "Accessing an invalid search object");
+ internal::AnalysisNeighborhoodPairSearchImpl pairSearch(*impl_);
+ pairSearch.startSearch(positions);
+ return pairSearch.searchNext(&withinAction);
+}
+
+real AnalysisNeighborhoodSearch::minimumDistance(
+ const AnalysisNeighborhoodPositions &positions) const
+{
+ GMX_RELEASE_ASSERT(impl_, "Accessing an invalid search object");
+ internal::AnalysisNeighborhoodPairSearchImpl pairSearch(*impl_);
+ pairSearch.startSearch(positions);
+ real minDist2 = impl_->cutoffSquared();
+ int closestPoint = -1;
+ MindistAction action(&closestPoint, &minDist2);
+ (void)pairSearch.searchNext(action);
+ return sqrt(minDist2);
+}
+
+AnalysisNeighborhoodPair
+AnalysisNeighborhoodSearch::nearestPoint(
+ const AnalysisNeighborhoodPositions &positions) const
+{
+ GMX_RELEASE_ASSERT(impl_, "Accessing an invalid search object");
+ internal::AnalysisNeighborhoodPairSearchImpl pairSearch(*impl_);
+ pairSearch.startSearch(positions);
+ real minDist2 = impl_->cutoffSquared();
+ int closestPoint = -1;
+ MindistAction action(&closestPoint, &minDist2);
+ (void)pairSearch.searchNext(action);
+ return AnalysisNeighborhoodPair(closestPoint, 0, minDist2);
+}
+
+AnalysisNeighborhoodPairSearch
+AnalysisNeighborhoodSearch::startPairSearch(
+ const AnalysisNeighborhoodPositions &positions) const
+{
+ GMX_RELEASE_ASSERT(impl_, "Accessing an invalid search object");
+ Impl::PairSearchImplPointer pairSearch(impl_->getPairSearch());
+ pairSearch->startSearch(positions);
+ return AnalysisNeighborhoodPairSearch(pairSearch);
+}
+
+/********************************************************************
+ * AnalysisNeighborhoodPairSearch
*/
-bool
-gmx_ana_nbsearch_next_within(gmx_ana_nbsearch_t *d, int *jp)
+
+AnalysisNeighborhoodPairSearch::AnalysisNeighborhoodPairSearch(
+ const ImplPointer &impl)
+ : impl_(impl)
{
- if (grid_search(d, &within_action))
- {
- *jp = d->previ;
- return true;
- }
- *jp = -1;
- return false;
}
+
+bool AnalysisNeighborhoodPairSearch::findNextPair(AnalysisNeighborhoodPair *pair)
+{
+ bool bFound = impl_->searchNext(&withinAction);
+ impl_->initFoundPair(pair);
+ return bFound;
+}
+
+void AnalysisNeighborhoodPairSearch::skipRemainingPairsForTestPosition()
+{
+ impl_->nextTestPosition();
+}
+
+} // namespace gmx