Use gmx_mtop_t in selections, part 3
[gromacs.git] / src / gromacs / selection / selection.h
blob1b9df4157d7181f2b581dccd17a5b9972d6a7e1e
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2009,2010,2011,2012,2013,2014,2015,2016, 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.
35 /*! \file
36 * \brief
37 * Declares gmx::Selection and supporting classes.
39 * \author Teemu Murtola <teemu.murtola@gmail.com>
40 * \inpublicapi
41 * \ingroup module_selection
43 #ifndef GMX_SELECTION_SELECTION_H
44 #define GMX_SELECTION_SELECTION_H
46 #include <string>
47 #include <vector>
49 #include "gromacs/selection/position.h"
50 #include "gromacs/selection/selectionenums.h"
51 #include "gromacs/utility/arrayref.h"
52 #include "gromacs/utility/classhelpers.h"
53 #include "gromacs/utility/gmxassert.h"
55 struct gmx_mtop_t;
57 namespace gmx
60 class SelectionOptionStorage;
61 class SelectionTreeElement;
63 class AnalysisNeighborhoodPositions;
64 class Selection;
65 class SelectionPosition;
67 //! Container of selections used in public selection interfaces.
68 typedef std::vector<Selection> SelectionList;
70 namespace internal
73 /*! \internal
74 * \brief
75 * Internal data for a single selection.
77 * This class is internal to the selection module, but resides in a public
78 * header because of efficiency reasons: it allows frequently used access
79 * methods in \ref Selection to be inlined.
81 * Methods in this class do not throw unless otherwise specified.
83 * \ingroup module_selection
85 class SelectionData
87 public:
88 /*! \brief
89 * Creates a new selection object.
91 * \param[in] elem Root of the evaluation tree for this selection.
92 * \param[in] selstr String that was parsed to produce this selection.
93 * \throws std::bad_alloc if out of memory.
95 SelectionData(SelectionTreeElement *elem, const char *selstr);
96 ~SelectionData();
98 //! Returns the name for this selection.
99 const char *name() const { return name_.c_str(); }
100 //! Returns the string that was parsed to produce this selection.
101 const char *selectionText() const { return selectionText_.c_str(); }
102 //! Returns true if the size of the selection (posCount()) is dynamic.
103 bool isDynamic() const { return bDynamic_; }
104 //! Returns the type of positions in the selection.
105 e_index_t type() const { return rawPositions_.m.type; }
106 //! Returns true if the selection only contains positions with a single atom each.
107 bool hasOnlyAtoms() const { return type() == INDEX_ATOM; }
108 //! Returns `true` if the atom indices in the selection are in ascending order.
109 bool hasSortedAtomIndices() const;
111 //! Number of positions in the selection.
112 int posCount() const { return rawPositions_.count(); }
113 //! Returns the root of the evaluation tree for this selection.
114 SelectionTreeElement &rootElement() { return rootElement_; }
116 //! Returns whether the covered fraction can change between frames.
117 bool isCoveredFractionDynamic() const { return bDynamicCoveredFraction_; }
119 //! Returns true if the given flag is set.
120 bool hasFlag(SelectionFlag flag) const { return flags_.test(flag); }
121 //! Sets the flags for this selection.
122 void setFlags(SelectionFlags flags) { flags_ = flags; }
124 //! \copydoc Selection::initCoveredFraction()
125 bool initCoveredFraction(e_coverfrac_t type);
127 /*! \brief
128 * Updates the name of the selection if missing.
130 * \throws std::bad_alloc if out of memory.
132 * If selections get their value from a group reference that cannot be
133 * resolved during parsing, the name is final only after group
134 * references have been resolved.
136 * This function is called by SelectionCollection::setIndexGroups().
138 void refreshName();
139 /*! \brief
140 * Computes total masses and charges for all selection positions.
142 * \param[in] top Topology information.
143 * \throws std::bad_alloc if out of memory.
145 * For dynamic selections, the values need to be updated after each
146 * evaluation with refreshMassesAndCharges().
147 * This is done by SelectionEvaluator.
149 * This function is called by SelectionCompiler.
151 * Strong exception safety guarantee.
153 void initializeMassesAndCharges(const gmx_mtop_t *top);
154 /*! \brief
155 * Updates masses and charges after dynamic selection has been
156 * evaluated.
158 * \param[in] top Topology information.
160 * Called by SelectionEvaluator.
162 void refreshMassesAndCharges(const gmx_mtop_t *top);
163 /*! \brief
164 * Updates the covered fraction after a selection has been evaluated.
166 * Called by SelectionEvaluator.
168 void updateCoveredFractionForFrame();
169 /*! \brief
170 * Computes average covered fraction after all frames have been evaluated.
172 * \param[in] nframes Number of frames that have been evaluated.
174 * \p nframes should be equal to the number of calls to
175 * updateCoveredFractionForFrame().
176 * Called by SelectionEvaluator::evaluateFinal().
178 void computeAverageCoveredFraction(int nframes);
179 /*! \brief
180 * Restores position information to state it was in after compilation.
182 * \param[in] top Topology information.
184 * Depends on SelectionCompiler storing the original atoms in the
185 * \a rootElement_ object.
186 * Called by SelectionEvaluator::evaluateFinal().
188 void restoreOriginalPositions(const gmx_mtop_t *top);
190 private:
191 //! Name of the selection.
192 std::string name_;
193 //! The actual selection string.
194 std::string selectionText_;
195 //! Low-level representation of selected positions.
196 gmx_ana_pos_t rawPositions_;
197 //! Total masses for the current positions.
198 std::vector<real> posMass_;
199 //! Total charges for the current positions.
200 std::vector<real> posCharge_;
201 SelectionFlags flags_;
202 //! Root of the selection evaluation tree.
203 SelectionTreeElement &rootElement_;
204 //! Type of the covered fraction.
205 e_coverfrac_t coveredFractionType_;
206 //! Covered fraction of the selection for the current frame.
207 real coveredFraction_;
208 //! The average covered fraction (over the trajectory).
209 real averageCoveredFraction_;
210 //! true if the value can change as a function of time.
211 bool bDynamic_;
212 //! true if the covered fraction depends on the frame.
213 bool bDynamicCoveredFraction_;
215 /*! \brief
216 * Needed to wrap access to information.
218 friend class gmx::Selection;
219 /*! \brief
220 * Needed for proper access to position information.
222 friend class gmx::SelectionPosition;
224 GMX_DISALLOW_COPY_AND_ASSIGN(SelectionData);
227 } // namespace internal
229 /*! \brief
230 * Provides access to a single selection.
232 * This class provides a public interface for accessing selection information.
233 * General information about the selection can be accessed with methods name(),
234 * selectionText(), isDynamic(), and type(). The first three can be accessed
235 * any time after the selection has been parsed, and type() can be accessed
236 * after the selection has been compiled.
238 * There are a few methods that can be used to change the behavior of the
239 * selection. setEvaluateVelocities() and setEvaluateForces() can be called
240 * before the selection is compiled to request evaluation of velocities and/or
241 * forces in addition to coordinates.
243 * Each selection is made of a set of positions. Each position has associated
244 * coordinates, and possibly velocities and forces if they have been requested
245 * and are available. It also has a set of atoms associated with it; typically
246 * the coordinates are the center-of-mass or center-of-geometry coordinates for
247 * that set of atoms. To access the number of positions in the selection, use
248 * posCount(). To access individual positions, use position().
249 * See SelectionPosition for details of how to use individual positions.
250 * setOriginalId() can be used to adjust the return value of
251 * SelectionPosition::mappedId(); see that method for details.
253 * It is also possible to access the list of atoms that make up all the
254 * positions directly: atomCount() returns the total number of atoms in the
255 * selection and atomIndices() an array of their indices.
256 * Similarly, it is possible to access the coordinates and other properties
257 * of the positions as continuous arrays through coordinates(), velocities(),
258 * forces(), masses(), charges(), refIds(), and mappedIds().
260 * Both positions and atoms can be accessed after the selection has been
261 * compiled. For dynamic selections, the return values of these methods change
262 * after each evaluation to reflect the situation for the current frame.
263 * Before any frame has been evaluated, these methods return the maximal set
264 * to which the selection can evaluate.
266 * There are two possible modes for how positions for dynamic selections are
267 * handled. In the default mode, posCount() can change, and for each frame,
268 * only the positions that are selected in that frame can be accessed. In a
269 * masked mode, posCount() remains constant, i.e., the positions are always
270 * evaluated for the maximal set, and SelectionPosition::selected() is used to
271 * determine whether a position is selected for a frame. The masked mode can
272 * be requested with SelectionOption::dynamicMask().
274 * The class also provides methods for printing out information: printInfo()
275 * and printDebugInfo(). These are mainly for internal use by Gromacs.
277 * This class works like a pointer type: copying and assignment is lightweight,
278 * and all copies work interchangeably, accessing the same internal data.
280 * Methods in this class do not throw.
282 * \see SelectionPosition
284 * \inpublicapi
285 * \ingroup module_selection
287 class Selection
289 public:
290 /*! \brief
291 * Creates a selection wrapper that has no associated selection.
293 * Any attempt to call methods in the object before a selection is
294 * assigned results in undefined behavior.
295 * isValid() returns `false` for the selection until it is initialized.
297 Selection() : sel_(NULL) {}
298 /*! \brief
299 * Creates a new selection object.
301 * \param sel Selection data to wrap.
303 * Only for internal use by the selection module.
305 explicit Selection(internal::SelectionData *sel) : sel_(sel) {}
307 //! Returns whether the selection object is initialized.
308 bool isValid() const { return sel_ != NULL; }
310 //! Returns whether two selection objects wrap the same selection.
311 bool operator==(const Selection &other) const
313 return sel_ == other.sel_;
315 //! Returns whether two selection objects wrap different selections.
316 bool operator!=(const Selection &other) const
318 return !operator==(other);
321 //! Returns the name of the selection.
322 const char *name() const { return data().name(); }
323 //! Returns the string that was parsed to produce this selection.
324 const char *selectionText() const { return data().selectionText(); }
325 //! Returns true if the size of the selection (posCount()) is dynamic.
326 bool isDynamic() const { return data().isDynamic(); }
327 //! Returns the type of positions in the selection.
328 e_index_t type() const { return data().type(); }
329 //! Returns true if the selection only contains positions with a single atom each.
330 bool hasOnlyAtoms() const { return data().hasOnlyAtoms(); }
331 //! Returns `true` if the atom indices in the selection are in ascending order.
332 bool hasSortedAtomIndices() const { return data().hasSortedAtomIndices(); }
334 //! Total number of atoms in the selection.
335 int atomCount() const
337 return data().rawPositions_.m.mapb.nra;
339 //! Returns atom indices of all atoms in the selection.
340 ConstArrayRef<int> atomIndices() const
342 return constArrayRefFromArray(sel_->rawPositions_.m.mapb.a,
343 sel_->rawPositions_.m.mapb.nra);
345 //! Number of positions in the selection.
346 int posCount() const { return data().posCount(); }
347 //! Access a single position.
348 SelectionPosition position(int i) const;
349 //! Returns coordinates for this selection as a continuous array.
350 ConstArrayRef<rvec> coordinates() const
352 return constArrayRefFromArray(data().rawPositions_.x, posCount());
354 //! Returns whether velocities are available for this selection.
355 bool hasVelocities() const { return data().rawPositions_.v != NULL; }
356 /*! \brief
357 * Returns velocities for this selection as a continuous array.
359 * Must not be called if hasVelocities() returns false.
361 ConstArrayRef<rvec> velocities() const
363 GMX_ASSERT(hasVelocities(), "Velocities accessed, but unavailable");
364 return constArrayRefFromArray(data().rawPositions_.v, posCount());
366 //! Returns whether forces are available for this selection.
367 bool hasForces() const { return sel_->rawPositions_.f != NULL; }
368 /*! \brief
369 * Returns forces for this selection as a continuous array.
371 * Must not be called if hasForces() returns false.
373 ConstArrayRef<rvec> forces() const
375 GMX_ASSERT(hasForces(), "Forces accessed, but unavailable");
376 return constArrayRefFromArray(data().rawPositions_.f, posCount());
378 //! Returns masses for this selection as a continuous array.
379 ConstArrayRef<real> masses() const
381 // posMass_ may have more entries than posCount() in the case of
382 // dynamic selections that don't have a topology
383 // (and thus the masses and charges are fixed).
384 GMX_ASSERT(data().posMass_.size() >= static_cast<size_t>(posCount()),
385 "Internal inconsistency");
386 return constArrayRefFromVector<real>(data().posMass_.begin(),
387 data().posMass_.begin() + posCount());
389 //! Returns charges for this selection as a continuous array.
390 ConstArrayRef<real> charges() const
392 // posCharge_ may have more entries than posCount() in the case of
393 // dynamic selections that don't have a topology
394 // (and thus the masses and charges are fixed).
395 GMX_ASSERT(data().posCharge_.size() >= static_cast<size_t>(posCount()),
396 "Internal inconsistency");
397 return constArrayRefFromVector<real>(data().posCharge_.begin(),
398 data().posCharge_.begin() + posCount());
400 /*! \brief
401 * Returns reference IDs for this selection as a continuous array.
403 * \see SelectionPosition::refId()
405 ConstArrayRef<int> refIds() const
407 return constArrayRefFromArray(data().rawPositions_.m.refid, posCount());
409 /*! \brief
410 * Returns mapped IDs for this selection as a continuous array.
412 * \see SelectionPosition::mappedId()
414 ConstArrayRef<int> mappedIds() const
416 return constArrayRefFromArray(data().rawPositions_.m.mapid, posCount());
419 //! Returns whether the covered fraction can change between frames.
420 bool isCoveredFractionDynamic() const { return data().isCoveredFractionDynamic(); }
421 //! Returns the covered fraction for the current frame.
422 real coveredFraction() const { return data().coveredFraction_; }
424 /*! \brief
425 * Allows passing a selection directly to neighborhood searching.
427 * When initialized this way, AnalysisNeighborhoodPair objects return
428 * indices that can be used to index the selection positions with
429 * position().
431 * Works exactly like if AnalysisNeighborhoodPositions had a
432 * constructor taking a Selection object as a parameter.
433 * See AnalysisNeighborhoodPositions for rationale and additional
434 * discussion.
436 operator AnalysisNeighborhoodPositions() const;
438 /*! \brief
439 * Initializes information about covered fractions.
441 * \param[in] type Type of covered fraction required.
442 * \returns true if the covered fraction can be calculated for the
443 * selection.
445 bool initCoveredFraction(e_coverfrac_t type)
447 return data().initCoveredFraction(type);
449 /*! \brief
450 * Sets whether this selection evaluates velocities for positions.
452 * \param[in] bEnabled If true, velocities are evaluated.
454 * If you request the evaluation, but then evaluate the selection for
455 * a frame that does not contain velocity information, results are
456 * undefined.
458 * \todo
459 * Implement it such that in the above case, hasVelocities() will
460 * return false for such frames.
462 * Does not throw.
464 void setEvaluateVelocities(bool bEnabled)
466 data().flags_.set(efSelection_EvaluateVelocities, bEnabled);
468 /*! \brief
469 * Sets whether this selection evaluates forces for positions.
471 * \param[in] bEnabled If true, forces are evaluated.
473 * If you request the evaluation, but then evaluate the selection for
474 * a frame that does not contain force information, results are
475 * undefined.
477 * Does not throw.
479 void setEvaluateForces(bool bEnabled)
481 data().flags_.set(efSelection_EvaluateForces, bEnabled);
484 /*! \brief
485 * Sets the ID for the \p i'th position for use with
486 * SelectionPosition::mappedId().
488 * \param[in] i Zero-based index
489 * \param[in] id Identifier to set.
491 * This method is not part of SelectionPosition because that interface
492 * only provides access to const data by design.
494 * This method can only be called after compilation, before the
495 * selection has been evaluated for any frame.
497 * \see SelectionPosition::mappedId()
499 void setOriginalId(int i, int id);
500 /*! \brief
501 * Inits the IDs for use with SelectionPosition::mappedId() for
502 * grouping.
504 * \param[in] top Topology information
505 * (can be NULL if not required for \p type).
506 * \param[in] type Type of groups to generate.
507 * \returns Number of groups that were present in the selection.
508 * \throws InconsistentInputError if the selection positions cannot
509 * be assigned to groups of the given type.
511 * If `type == INDEX_ATOM`, the IDs are initialized to 0, 1, 2, ...,
512 * and the return value is the number of positions.
513 * If `type == INDEX_ALL`, all the IDs are initialized to 0, and the
514 * return value is one.
515 * If `type == INDEX_RES` or `type == INDEX_MOL`, the first position
516 * will get ID 0, and all following positions that belong to the same
517 * residue/molecule will get the same ID. The first position that
518 * belongs to a different residue/molecule will get ID 1, and so on.
519 * If some position contains atoms from multiple residues/molecules,
520 * i.e., the mapping is ambiguous, an exception is thrown.
521 * The return value is the number of residues/molecules that are
522 * present in the selection positions.
524 * This method is useful if the calling code needs to group the
525 * selection, e.g., for computing aggregate properties for each residue
526 * or molecule. It can then use this method to initialize the
527 * appropriate grouping, use the return value to allocate a
528 * sufficiently sized buffer to store the aggregated values, and then
529 * use SelectionPosition::mappedId() to identify the location where to
530 * aggregate to.
532 * \see setOriginalId()
533 * \see SelectionPosition::mappedId()
535 int initOriginalIdsToGroup(const gmx_mtop_t *top, e_index_t type);
537 /*! \brief
538 * Prints out one-line description of the selection.
540 * \param[in] fp Where to print the information.
542 * The output contains the name of the selection, the number of atoms
543 * and the number of positions, and indication of whether the selection
544 * is dynamic.
546 void printInfo(FILE *fp) const;
547 /*! \brief
548 * Prints out extended information about the selection for debugging.
550 * \param[in] fp Where to print the information.
551 * \param[in] nmaxind Maximum number of values to print in lists
552 * (-1 = print all).
554 void printDebugInfo(FILE *fp, int nmaxind) const;
556 private:
557 internal::SelectionData &data()
559 GMX_ASSERT(sel_ != NULL,
560 "Attempted to access uninitialized selection");
561 return *sel_;
563 const internal::SelectionData &data() const
565 GMX_ASSERT(sel_ != NULL,
566 "Attempted to access uninitialized selection");
567 return *sel_;
570 /*! \brief
571 * Pointer to internal data for the selection.
573 * The memory for this object is managed by a SelectionCollection
574 * object, and the \ref Selection class simply provides a public
575 * interface for accessing the data.
577 internal::SelectionData *sel_;
579 /*! \brief
580 * Needed to access the data to adjust flags.
582 friend class SelectionOptionStorage;
585 /*! \brief
586 * Provides access to information about a single selected position.
588 * Each position has associated coordinates, and possibly velocities and forces
589 * if they have been requested and are available. It also has a set of atoms
590 * associated with it; typically the coordinates are the center-of-mass or
591 * center-of-geometry coordinates for that set of atoms. It is possible that
592 * there are not atoms associated if the selection has been provided as a fixed
593 * position.
595 * After the selection has been compiled, but not yet evaluated, the contents
596 * of the coordinate, velocity and force vectors are undefined.
598 * Default copy constructor and assignment operators are used, and work as
599 * intended: the copy references the same position and works identically.
601 * Methods in this class do not throw.
603 * \see Selection
605 * \inpublicapi
606 * \ingroup module_selection
608 class SelectionPosition
610 public:
611 /*! \brief
612 * Constructs a wrapper object for given selection position.
614 * \param[in] sel Selection from which the position is wrapped.
615 * \param[in] index Zero-based index of the position to wrap.
617 * Asserts if \p index is out of range.
619 * Only for internal use of the library. To obtain a SelectionPosition
620 * object in other code, use Selection::position().
622 SelectionPosition(const internal::SelectionData &sel, int index)
623 : sel_(&sel), i_(index)
625 GMX_ASSERT(index >= 0 && index < sel.posCount(),
626 "Invalid selection position index");
629 /*! \brief
630 * Returns type of this position.
632 * Currently always returns the same as Selection::type().
634 e_index_t type() const { return sel_->type(); }
635 //! Returns coordinates for this position.
636 const rvec &x() const
638 return sel_->rawPositions_.x[i_];
640 /*! \brief
641 * Returns velocity for this position.
643 * Must not be called if Selection::hasVelocities() returns false.
645 const rvec &v() const
647 GMX_ASSERT(sel_->rawPositions_.v != NULL,
648 "Velocities accessed, but unavailable");
649 return sel_->rawPositions_.v[i_];
651 /*! \brief
652 * Returns force for this position.
654 * Must not be called if Selection::hasForces() returns false.
656 const rvec &f() const
658 GMX_ASSERT(sel_->rawPositions_.f != NULL,
659 "Velocities accessed, but unavailable");
660 return sel_->rawPositions_.f[i_];
662 /*! \brief
663 * Returns total mass for this position.
665 * Returns the total mass of atoms that make up this position.
666 * If there are no atoms associated or masses are not available,
667 * returns unity.
669 real mass() const
671 return sel_->posMass_[i_];
673 /*! \brief
674 * Returns total charge for this position.
676 * Returns the sum of charges of atoms that make up this position.
677 * If there are no atoms associated or charges are not available,
678 * returns zero.
680 real charge() const
682 return sel_->posCharge_[i_];
684 //! Returns the number of atoms that make up this position.
685 int atomCount() const
687 return sel_->rawPositions_.m.mapb.index[i_ + 1]
688 - sel_->rawPositions_.m.mapb.index[i_];
690 //! Return atom indices that make up this position.
691 ConstArrayRef<int> atomIndices() const
693 const int *atoms = sel_->rawPositions_.m.mapb.a;
694 if (atoms == NULL)
696 return ConstArrayRef<int>();
698 const int first = sel_->rawPositions_.m.mapb.index[i_];
699 return constArrayRefFromArray(&atoms[first], atomCount());
701 /*! \brief
702 * Returns whether this position is selected in the current frame.
704 * The return value is equivalent to \c refid() == -1. Returns always
705 * true if SelectionOption::dynamicMask() has not been set.
707 * \see refId()
709 bool selected() const
711 return refId() >= 0;
713 /*! \brief
714 * Returns reference ID for this position.
716 * For dynamic selections, this provides means to associate positions
717 * across frames. After compilation, these IDs are consequently
718 * numbered starting from zero. For each frame, the ID then reflects
719 * the location of the position in the original array of positions.
720 * If SelectionOption::dynamicMask() has been set for the parent
721 * selection, the IDs for positions not present in the current
722 * selection are set to -1, otherwise they are removed completely.
724 * Example:
725 * If a dynamic selection consists of at most three positions, after
726 * compilation refId() will return 0, 1, 2 for them, respectively.
727 * If for a particular frame, only the first and the third are present,
728 * refId() will return 0, 2.
729 * If SelectionOption::dynamicMask() has been set, all three positions
730 * can be accessed also for that frame and refId() will return 0, -1,
731 * 2.
733 int refId() const
735 return sel_->rawPositions_.m.refid[i_];
737 /*! \brief
738 * Returns mapped ID for this position.
740 * Returns ID of the position that corresponds to that set with
741 * Selection::setOriginalId().
743 * If for an array \c id, \c setOriginalId(i, id[i]) has been called
744 * for each \c i, then it always holds that
745 * \c mappedId()==id[refId()].
747 * Selection::setOriginalId() has not been called, the default values
748 * are dependent on type():
749 * - ::INDEX_ATOM: atom indices
750 * - ::INDEX_RES: residue indices
751 * - ::INDEX_MOL: molecule indices
753 * All the default values are zero-based.
755 int mappedId() const
757 return sel_->rawPositions_.m.mapid[i_];
760 /*! \brief
761 * Allows passing a selection position directly to neighborhood searching.
763 * When initialized this way, AnalysisNeighborhoodPair objects return
764 * the index that can be used to access this position using
765 * Selection::position().
767 * Works exactly like if AnalysisNeighborhoodPositions had a
768 * constructor taking a SelectionPosition object as a parameter.
769 * See AnalysisNeighborhoodPositions for rationale and additional
770 * discussion.
772 operator AnalysisNeighborhoodPositions() const;
774 private:
775 const internal::SelectionData *sel_;
776 int i_;
780 inline SelectionPosition
781 Selection::position(int i) const
783 return SelectionPosition(data(), i);
786 } // namespace gmx
788 #endif