Proper C++ API for analysis neighborhood search.
[alexxy/gromacs.git] / src / gromacs / selection / nbsearch.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2009,2010,2011,2012,2013, by the GROMACS development team, led by
5  * David van der Spoel, Berk Hess, Erik Lindahl, and including many
6  * others, as listed in the AUTHORS file in the top-level source
7  * directory and at http://www.gromacs.org.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35 /*! \file
36  * \brief API for neighborhood searching.
37  *
38  * The API is documented in more detail on a separate page:
39  * \ref nbsearch
40  *
41  * The functions within this file can be used independently of the other parts
42  * of the library.
43  * The library also uses the functions internally.
44  *
45  * \author Teemu Murtola <teemu.murtola@gmail.com>
46  * \inpublicapi
47  * \ingroup module_selection
48  */
49 #ifndef GMX_SELECTION_NBSEARCH_H
50 #define GMX_SELECTION_NBSEARCH_H
51
52 #include <boost/shared_ptr.hpp>
53
54 #include "../legacyheaders/typedefs.h"
55 #include "../utility/common.h"
56 #include "../utility/gmxassert.h"
57
58 #include "indexutil.h"
59
60 struct gmx_ana_pos_t;
61
62 /** Data structure for neighborhood searches. */
63 typedef struct gmx_ana_nbsearch_t gmx_ana_nbsearch_t;
64
65 /** Create a new neighborhood search data structure. */
66 gmx_ana_nbsearch_t *
67 gmx_ana_nbsearch_create(real cutoff, int maxn);
68 /** Free memory allocated for neighborhood search. */
69 void
70 gmx_ana_nbsearch_free(gmx_ana_nbsearch_t *d);
71
72 /** Initializes neighborhood search for a new frame. */
73 void
74 gmx_ana_nbsearch_init(gmx_ana_nbsearch_t *d, t_pbc *pbc, int n, const rvec x[]);
75 /** Initializes neighborhood search for a frame using \c gmx_ana_pos_t.  */
76 void
77 gmx_ana_nbsearch_pos_init(gmx_ana_nbsearch_t *d, t_pbc *pbc,
78                           const struct gmx_ana_pos_t *p);
79 /** Sets the exclusions for the next neighborhood search. */
80 void
81 gmx_ana_nbsearch_set_excl(gmx_ana_nbsearch_t *d, int nexcl, int excl[]);
82 /** Check whether a point is within a neighborhood. */
83 bool
84 gmx_ana_nbsearch_is_within(gmx_ana_nbsearch_t *d, const rvec x);
85 /** Check whether a position is within a neighborhood. */
86 bool
87 gmx_ana_nbsearch_pos_is_within(gmx_ana_nbsearch_t *d,
88                                const struct gmx_ana_pos_t *p, int i);
89 /** Calculates the minimun distance from the reference points. */
90 real
91 gmx_ana_nbsearch_mindist(gmx_ana_nbsearch_t *d, const rvec x);
92 /** Calculates the minimun distance from the reference points. */
93 real
94 gmx_ana_nbsearch_pos_mindist(gmx_ana_nbsearch_t *d,
95                              const struct gmx_ana_pos_t *p, int i);
96 /** Finds the first reference position within the cutoff. */
97 bool
98 gmx_ana_nbsearch_first_within(gmx_ana_nbsearch_t *d, const rvec x, int *jp);
99 /** Finds the first reference position within the cutoff. */
100 bool
101 gmx_ana_nbsearch_pos_first_within(gmx_ana_nbsearch_t *d,
102                                   const struct gmx_ana_pos_t *p, int i, int *jp);
103 /** Finds the next reference position within the cutoff. */
104 bool
105 gmx_ana_nbsearch_next_within(gmx_ana_nbsearch_t *d, int *jp);
106
107 namespace gmx
108 {
109
110 namespace internal
111 {
112 class AnalysisNeighborhoodSearchImpl;
113 class AnalysisNeighborhoodPairSearchImpl;
114 };
115
116 class AnalysisNeighborhoodSearch;
117 class AnalysisNeighborhoodPairSearch;
118
119 /*! \brief
120  * Neighborhood searching for analysis tools.
121  *
122  * This class implements neighborhood searching routines for analysis tools.
123  * The emphasis is in flexibility and ease of use; one main driver is to have
124  * a common implementation of grid-based searching to avoid replicating this in
125  * multiple tools (and to make more tools take advantage of the significant
126  * performance improvement this allows).
127  *
128  * To use the search, create an object of this type, call setCutoff() to
129  * initialize it, and then repeatedly call initSearch() to start a search with
130  * different sets of reference positions.  For each set of reference positions,
131  * use methods in the returned AnalysisNeighborhoodSearch to find the reference
132  * positions that are within the given cutoff from a provided position.
133  *
134  * \todo
135  * Support for exclusions.
136  * The 4.5/4.6 C API had very low-level support for exclusions, which was not
137  * very convenient to use, and hadn't been tested much.  The internal code that
138  * it used to do the exclusion during the search itself is still there, but it
139  * needs more thought on what would be a convenient way to initialize it.
140  * Can be implemented once there is need for it in some calling code.
141  *
142  * \inpublicapi
143  * \ingroup module_selection
144  */
145 class AnalysisNeighborhood
146 {
147     public:
148         //! Searching algorithm to use.
149         enum SearchMode
150         {
151             //! Select algorithm based on heuristic efficiency considerations.
152             eSearchMode_Automatic,
153             //! Use a simple loop over all pairs.
154             eSearchMode_Simple,
155             //! Use grid-based searching whenever possible.
156             eSearchMode_Grid
157         };
158
159         //! Creates an uninitialized neighborhood search.
160         AnalysisNeighborhood();
161         ~AnalysisNeighborhood();
162
163         /*! \brief
164          * Set cutoff distance for the neighborhood searching.
165          *
166          * \param[in]  cutoff Cutoff distance for the search
167          *   (<=0 stands for no cutoff).
168          *
169          * Currently, needs to be called exactly once, before calling any other
170          * method.
171          */
172         void setCutoff(real cutoff);
173         /*! \brief
174          * Sets the algorithm to use for searching.
175          *
176          * \param[in] mode  Search mode to use.
177          *
178          * Note that if \p mode is \ref eSearchMode_Grid, it is still only a
179          * suggestion: grid-based searching may not be possible with the
180          * provided input, in which case a simple search is still used.
181          * This is mainly useful for testing purposes to force a mode.
182          *
183          * Does not throw.
184          */
185         void setMode(SearchMode mode);
186         //! Returns the currently active search mode.
187         SearchMode mode() const;
188
189         /*! \brief
190          * Initializes neighborhood search for a set of positions.
191          *
192          * \param[in] pbc PBC information for the frame.
193          * \param[in] n   Number of reference positions for the frame.
194          * \param[in] x   \p n reference positions for the frame.
195          * \returns   Search object that can be used to find positions from
196          *      \p x within the given cutoff.
197          */
198         AnalysisNeighborhoodSearch
199         initSearch(const t_pbc *pbc, int n, const rvec x[]);
200         /*! \brief
201          * Initializes neighborhood search for a set of positions.
202          *
203          * \param[in] pbc PBC information for the frame.
204          * \param[in] p   Reference positions for the frame.
205          * \returns   Search object that can be used to find positions from
206          *      \p p within the given cutoff.
207          */
208         AnalysisNeighborhoodSearch
209         initSearch(const t_pbc *pbc, const gmx_ana_pos_t *p);
210
211     private:
212         class Impl;
213
214         PrivateImplPointer<Impl> impl_;
215 };
216
217 /*! \brief
218  * Value type to represent a pair of positions found in neighborhood searching.
219  *
220  * Methods in this class do not throw.
221  *
222  * \inpublicapi
223  * \ingroup module_selection
224  */
225 class AnalysisNeighborhoodPair
226 {
227     public:
228         //! Initializes an invalid pair.
229         AnalysisNeighborhoodPair() : refIndex_(-1), testIndex_(0) {}
230         //! Initializes a pair object with the given data.
231         AnalysisNeighborhoodPair(int refIndex, int testIndex)
232             : refIndex_(refIndex), testIndex_(testIndex)
233         {
234         }
235
236         /*! \brief
237          * Whether this pair is valid.
238          *
239          * If isValid() returns false, other methods should not be called.
240          */
241         bool isValid() const { return refIndex_ >= 0; }
242
243         /*! \brief
244          * Returns the index of the reference position in the pair.
245          *
246          * This index is always the index into the position array provided to
247          * AnalysisNeighborhood::initSearch().
248          */
249         int refIndex() const
250         {
251             GMX_ASSERT(isValid(), "Accessing invalid object");
252             return refIndex_;
253         }
254         /*! \brief
255          * Returns the index of the test position in the pair.
256          *
257          * The contents of this index depends on the context (method call) that
258          * produces the pair.
259          * If there was no array in the call, this index is zero.
260          */
261         int testIndex() const
262         {
263             GMX_ASSERT(isValid(), "Accessing invalid object");
264             return testIndex_;
265         }
266
267     private:
268         int                     refIndex_;
269         int                     testIndex_;
270 };
271
272 /*! \brief
273  * Initialized neighborhood search with a fixed set of reference positions.
274  *
275  * An instance of this class is obtained through
276  * AnalysisNeighborhood::initSearch(), and can be used to do multiple searches
277  * against the provided set of reference positions.
278  * Currently, it is not possible to call any method in this class while an
279  * AnalysisNeighborhoodPairSearch object obtained from startPairSearch() of the
280  * same instance exists.
281  *
282  * This class works like a pointer: copies of it point to the same search.
283  * In general, avoid creating copies, and only use the copy/assignment support
284  * for moving the variable around.  With C++11, this class would best be
285  * movable.
286  *
287  * Methods in this class do not throw.
288  *
289  * \todo
290  * Make it such that reset() is not necessary to call in code that repeatedly
291  * assigns the result of AnalysisNeighborhood::initSearch() to the same
292  * variable (see sm_distance.cpp).
293  *
294  * \todo
295  * Consider merging nearestPoint() and minimumDistance() by adding the distance
296  * to AnalysisNeighborhoodPair.
297  *
298  * \inpublicapi
299  * \ingroup module_selection
300  */
301 class AnalysisNeighborhoodSearch
302 {
303     public:
304         /*! \brief
305          * Internal short-hand type for a pointer to the implementation class.
306          *
307          * shared_ptr is used here to automatically keep a reference count to
308          * track whether an implementation class is still used outside the
309          * AnalysisNeighborhood object.  Ownership currently always stays with
310          * AnalysisNeighborhood; it always keeps one instance of the pointer.
311          */
312         typedef boost::shared_ptr<internal::AnalysisNeighborhoodSearchImpl>
313             ImplPointer;
314
315         /*! \brief
316          * Initializes an invalid search.
317          *
318          * Such an object cannot be used for searching.  It needs to be
319          * assigned a value from AnalysisNeighborhood::initSearch() before it
320          * can be used.  Provided to allow declaring a variable to hold the
321          * search before calling AnalysisNeighborhood::initSearch().
322          */
323         AnalysisNeighborhoodSearch();
324         /*! \brief
325          * Internally initialize the search.
326          *
327          * Used to implement AnalysisNeighborhood::initSearch().
328          * Cannot be called from user code.
329          */
330         explicit AnalysisNeighborhoodSearch(const ImplPointer &impl);
331
332         /*! \brief
333          * Clears this search.
334          *
335          * Equivalent to \c "*this = AnalysisNeighborhoodSearch();".
336          * Currently, this is required if the previous search variable is still
337          * in scope and you want to call AnalysisNeighborhood::initSearch()
338          * again.
339          */
340         void reset();
341
342         /*! \brief
343          * Returns the searching algorithm that this search is using.
344          *
345          * The return value is never AnalysisNeighborhood::eSearchMode_Automatic.
346          */
347         AnalysisNeighborhood::SearchMode mode() const;
348
349         /*! \brief
350          * Check whether a point is within a neighborhood.
351          *
352          * \param[in] x  Test position.
353          * \returns   true if the test position is within the cutoff of any
354          *     reference position.
355          */
356         bool isWithin(const rvec x) const;
357         /*! \brief
358          * Check whether a point is within a neighborhood.
359          *
360          * \param[in] p  Test positions.
361          * \param[in] i  Use the i'th position in \p p for testing.
362          * \returns   true if the test position is within the cutoff of any
363          *     reference position.
364          */
365         bool isWithin(const gmx_ana_pos_t *p, int i) const;
366         /*! \brief
367          * Calculates the minimum distance from the reference points.
368          *
369          * \param[in] x  Test position.
370          * \returns   The distance to the nearest reference position, or the
371          *     cutoff value if there are no reference positions within the
372          *     cutoff.
373          */
374         real minimumDistance(const rvec x) const;
375         /*! \brief
376          * Calculates the minimum distance from the reference points.
377          *
378          * \param[in] p  Test positions.
379          * \param[in] i  Use the i'th position in \p p for testing.
380          * \returns   The distance to the nearest reference position, or the
381          *     cutoff value if there are no reference positions within the
382          *     cutoff.
383          */
384         real minimumDistance(const gmx_ana_pos_t *p, int i) const;
385         /*! \brief
386          * Finds the closest reference point.
387          *
388          * \param[in] x  Test position.
389          * \returns   The reference index identifies the reference position
390          *     that is closest to the test position.
391          *     The test index is always zero.  The returned pair is invalid if
392          *     no reference position is within the cutoff.
393          */
394         AnalysisNeighborhoodPair nearestPoint(const rvec x) const;
395         /*! \brief
396          * Finds the closest reference point.
397          *
398          * \param[in] p  Test positions.
399          * \param[in] i  Use the i'th position in \p p for testing.
400          * \returns   The reference index identifies the reference position
401          *     that is closest to the test position.
402          *     The test index is always \p i.  The returned pair is invalid if
403          *     no reference position is within the cutoff.
404          */
405         AnalysisNeighborhoodPair nearestPoint(const gmx_ana_pos_t *p, int i) const;
406
407         /*! \brief
408          * Start a search to find reference positions within a cutoff.
409          *
410          * \param[in] x  Test position to search the neighbors for.
411          * \returns   Initialized search object to loop through all reference
412          *     positions within the configured cutoff.
413          *
414          * In the AnalysisNeighborhoodPair objects returned by the search, the
415          * test index is always zero.
416          */
417         AnalysisNeighborhoodPairSearch startPairSearch(const rvec x);
418         /*! \brief
419          * Start a search to find reference positions within a cutoff.
420          *
421          * \param[in] p  Test positions.
422          * \param[in] i  Use the i'th position in \p p for testing.
423          * \returns   Initialized search object to loop through all reference
424          *     positions within the configured cutoff.
425          *
426          * In the AnalysisNeighborhoodPair objects returned by the search, the
427          * test index is always \p i.
428          */
429         AnalysisNeighborhoodPairSearch startPairSearch(const gmx_ana_pos_t *p, int i);
430
431     private:
432         ImplPointer             impl_;
433 };
434
435 /*! \brief
436  * Initialized neighborhood pair search with a fixed set of positions.
437  *
438  * This class is used to loop through pairs of neighbors within the cutoff
439  * provided to AnalysisNeighborhood.  The following code demonstrates its use:
440  * \code
441    gmx::AnalysisNeighborhood       nb;
442    nb.setCutoff(cutoff);
443    gmx::AnalysisNeighborhoodSearch search = nb.initSearch(pbc, nref, xref);
444    gmx::AnalysisNeighborhoodPairSearch pairSearch = search.startPairSearch(x);
445    gmx::AnalysisNeighborhoodPair pair;
446    while (pairSearch.findNextPair(&pair))
447    {
448        // <do something for each found pair the information in pair>
449    }
450  * \endcode
451  *
452  * It is not possible to use a single search object from multiple threads
453  * concurrently.
454  *
455  * This class works like a pointer: copies of it point to the same search.
456  * In general, avoid creating copies, and only use the copy/assignment support
457  * for moving the variable around.  With C++11, this class would best be
458  * movable.
459  *
460  * Methods in this class do not throw.
461  *
462  * \inpublicapi
463  * \ingroup module_selection
464  */
465 class AnalysisNeighborhoodPairSearch
466 {
467     public:
468         /*! \brief
469          * Internal short-hand type for a pointer to the implementation class.
470          *
471          * See AnalysisNeighborhoodSearch::ImplPointer for rationale of using
472          * shared_ptr and ownership semantics.
473          */
474         typedef boost::shared_ptr<internal::AnalysisNeighborhoodPairSearchImpl>
475             ImplPointer;
476
477         /*! \brief
478          * Internally initialize the search.
479          *
480          * Used to implement AnalysisNeighborhoodSearch::startPairSearch().
481          * Cannot be called from user code.
482          */
483         explicit AnalysisNeighborhoodPairSearch(const ImplPointer &impl);
484
485         /*! \brief
486          * Finds the next pair within the cutoff.
487          *
488          * \param[out] pair  Information about the found pair.
489          * \returns    false if there were no more pairs.
490          *
491          * If the method returns false, \p pair will be invalid.
492          *
493          * \see AnalysisNeighborhoodPair
494          * \see AnalysisNeighborhoodSearch::startPairSearch()
495          */
496         bool findNextPair(AnalysisNeighborhoodPair *pair);
497
498     private:
499         ImplPointer             impl_;
500 };
501
502 } // namespace gmx
503
504 #endif