2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2009,2010,2011,2012,2013,2014, 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 "thread_mpi/mutex.h"
67 #include "gromacs/math/vec.h"
68 #include "gromacs/pbcutil/pbc.h"
69 #include "gromacs/selection/position.h"
70 #include "gromacs/utility/arrayref.h"
71 #include "gromacs/utility/gmxassert.h"
72 #include "gromacs/utility/smalloc.h"
80 /********************************************************************
81 * Implementation class declarations
84 class AnalysisNeighborhoodSearchImpl
87 typedef AnalysisNeighborhoodPairSearch::ImplPointer
88 PairSearchImplPointer;
89 typedef std::vector<PairSearchImplPointer> PairSearchList;
90 typedef std::vector<std::vector<int> > CellList;
92 explicit AnalysisNeighborhoodSearchImpl(real cutoff);
93 ~AnalysisNeighborhoodSearchImpl();
96 * Initializes the search with a given box and reference positions.
98 * \param[in] mode Search mode to use.
99 * \param[in] pbc PBC information.
100 * \param[in] positions Set of reference positions.
102 void init(AnalysisNeighborhood::SearchMode mode,
104 const AnalysisNeighborhoodPositions &positions);
105 PairSearchImplPointer getPairSearch();
107 real cutoffSquared() const { return cutoff2_; }
108 bool usesGridSearch() const { return bGrid_; }
111 //! Calculates offsets to neighboring grid cells that should be considered.
112 void initGridCellNeighborList();
114 * Determines a suitable grid size and sets up the cells.
116 * \param[in] pbc Information about the box.
117 * \returns false if grid search is not suitable.
119 bool initGridCells(const t_pbc *pbc);
121 * Sets ua a search grid for a given box.
123 * \param[in] pbc Information about the box.
124 * \returns false if grid search is not suitable.
126 bool initGrid(const t_pbc *pbc);
128 * Maps a point into a grid cell.
130 * \param[in] x Point to map.
131 * \param[out] cell Indices of the grid cell in which \p x lies.
133 * \p x should be within the triclinic unit cell.
135 void mapPointToGridCell(const rvec x, ivec cell) const;
137 * Calculates linear index of a grid cell.
139 * \param[in] cell Cell indices.
140 * \returns Linear index of \p cell.
142 int getGridCellIndex(const ivec cell) const;
144 * Adds an index into a grid cell.
146 * \param[in] cell Cell into which \p i should be added.
147 * \param[in] i Index to add.
149 void addToGridCell(const ivec cell, int i);
151 //! Whether to try grid searching.
155 //! The cutoff squared.
158 //! Number of reference points for the current frame.
160 //! Reference point positions.
162 //! Reference position ids (NULL if not available).
167 //! Number of excluded reference positions for current test particle.
169 //! Exclusions for current test particle.
172 //! Whether grid searching is actually used for the current positions.
174 //! Array allocated for storing in-unit-cell reference positions.
176 //! Allocation count for xref_alloc.
178 //! false if the box is rectangular.
180 //! Box vectors of a single grid cell.
182 //! The reciprocal cell vectors as columns; the inverse of \p cellbox.
184 //! Number of cells along each dimension.
186 //! Data structure to hold the grid cell contents.
188 //! Number of neighboring cells to consider.
190 //! Offsets of the neighboring cells to consider.
192 //! Allocation count for \p gnboffs.
195 tMPI::mutex createPairSearchMutex_;
196 PairSearchList pairSearchList_;
198 friend class AnalysisNeighborhoodPairSearchImpl;
200 GMX_DISALLOW_COPY_AND_ASSIGN(AnalysisNeighborhoodSearchImpl);
203 class AnalysisNeighborhoodPairSearchImpl
206 explicit AnalysisNeighborhoodPairSearchImpl(const AnalysisNeighborhoodSearchImpl &search)
210 clear_ivec(testcell_);
214 //! Initializes a search to find reference positions neighboring \p x.
215 void startSearch(const AnalysisNeighborhoodPositions &positions);
216 //! Searches for the next neighbor.
217 template <class Action>
218 bool searchNext(Action action);
219 //! Initializes a pair representing the pair found by searchNext().
220 void initFoundPair(AnalysisNeighborhoodPair *pair) const;
221 //! Advances to the next test position, skipping any remaining pairs.
222 void nextTestPosition();
225 //! Clears the loop indices.
226 void reset(int testIndex);
227 //! Checks whether a reference positiong should be excluded.
228 bool isExcluded(int j);
230 //! Parent search object.
231 const AnalysisNeighborhoodSearchImpl &search_;
232 //! Reference to the test positions.
233 ConstArrayRef<rvec> testPositions_;
234 //! Index of the currently active test position in \p testPositions_.
236 //! Stores test position during a pair loop.
238 //! Stores the previous returned position during a pair loop.
240 //! Stores the current exclusion index during loops.
242 //! Stores the test particle cell index during loops.
244 //! Stores the current cell neighbor index during pair loops.
246 //! Stores the index within the current cell during pair loops.
249 GMX_DISALLOW_COPY_AND_ASSIGN(AnalysisNeighborhoodPairSearchImpl);
252 /********************************************************************
253 * AnalysisNeighborhoodSearchImpl
256 AnalysisNeighborhoodSearchImpl::AnalysisNeighborhoodSearchImpl(real cutoff)
262 cutoff_ = GMX_REAL_MAX;
265 cutoff2_ = sqr(cutoff_);
281 clear_mat(recipcell_);
282 clear_ivec(ncelldim_);
289 AnalysisNeighborhoodSearchImpl::~AnalysisNeighborhoodSearchImpl()
291 PairSearchList::const_iterator i;
292 for (i = pairSearchList_.begin(); i != pairSearchList_.end(); ++i)
294 GMX_RELEASE_ASSERT(i->unique(),
295 "Dangling AnalysisNeighborhoodPairSearch reference");
301 AnalysisNeighborhoodSearchImpl::PairSearchImplPointer
302 AnalysisNeighborhoodSearchImpl::getPairSearch()
304 tMPI::lock_guard<tMPI::mutex> lock(createPairSearchMutex_);
305 // TODO: Consider whether this needs to/can be faster, e.g., by keeping a
306 // separate pool of unused search objects.
307 PairSearchList::const_iterator i;
308 for (i = pairSearchList_.begin(); i != pairSearchList_.end(); ++i)
315 PairSearchImplPointer pairSearch(new AnalysisNeighborhoodPairSearchImpl(*this));
316 pairSearchList_.push_back(pairSearch);
320 void AnalysisNeighborhoodSearchImpl::initGridCellNeighborList()
322 int maxx, maxy, maxz;
325 /* Find the extent of the sphere in triclinic coordinates */
326 maxz = static_cast<int>(cutoff_ * recipcell_[ZZ][ZZ]) + 1;
327 rvnorm = sqrt(sqr(recipcell_[YY][YY]) + sqr(recipcell_[ZZ][YY]));
328 maxy = static_cast<int>(cutoff_ * rvnorm) + 1;
329 rvnorm = sqrt(sqr(recipcell_[XX][XX]) + sqr(recipcell_[YY][XX])
330 + sqr(recipcell_[ZZ][XX]));
331 maxx = static_cast<int>(cutoff_ * rvnorm) + 1;
333 /* Calculate the number of cells and reallocate if necessary */
334 ngridnb_ = (2 * maxx + 1) * (2 * maxy + 1) * (2 * maxz + 1);
335 if (gnboffs_nalloc_ < ngridnb_)
337 gnboffs_nalloc_ = ngridnb_;
338 srenew(gnboffs_, gnboffs_nalloc_);
341 /* Store the whole cube */
342 /* TODO: Prune off corners that are not needed */
344 for (int x = -maxx; x <= maxx; ++x)
346 for (int y = -maxy; y <= maxy; ++y)
348 for (int z = -maxz; z <= maxz; ++z)
359 bool AnalysisNeighborhoodSearchImpl::initGridCells(const t_pbc *pbc)
361 const real targetsize =
362 pow(pbc->box[XX][XX] * pbc->box[YY][YY] * pbc->box[ZZ][ZZ]
363 * 10 / nref_, static_cast<real>(1./3.));
366 for (int dd = 0; dd < DIM; ++dd)
368 ncelldim_[dd] = static_cast<int>(pbc->box[dd][dd] / targetsize);
369 cellCount *= ncelldim_[dd];
370 if (ncelldim_[dd] < 3)
375 // Never decrease the size of the cell vector to avoid reallocating
376 // memory for the nested vectors. The actual size of the vector is not
377 // used outside this function.
378 if (cells_.size() < static_cast<size_t>(cellCount))
380 cells_.resize(cellCount);
382 for (int ci = 0; ci < cellCount; ++ci)
389 bool AnalysisNeighborhoodSearchImpl::initGrid(const t_pbc *pbc)
391 /* TODO: This check could be improved. */
392 if (0.5*pbc->max_cutoff2 < cutoff2_)
397 if (!initGridCells(pbc))
402 bTric_ = TRICLINIC(pbc->box);
405 for (int dd = 0; dd < DIM; ++dd)
407 svmul(1.0 / ncelldim_[dd], pbc->box[dd], cellbox_[dd]);
409 m_inv_ur0(cellbox_, recipcell_);
413 for (int dd = 0; dd < DIM; ++dd)
415 cellbox_[dd][dd] = pbc->box[dd][dd] / ncelldim_[dd];
416 recipcell_[dd][dd] = 1.0 / cellbox_[dd][dd];
419 initGridCellNeighborList();
423 void AnalysisNeighborhoodSearchImpl::mapPointToGridCell(const rvec x,
430 tmvmul_ur0(recipcell_, x, tx);
431 for (int dd = 0; dd < DIM; ++dd)
433 cell[dd] = static_cast<int>(tx[dd]);
438 for (int dd = 0; dd < DIM; ++dd)
440 cell[dd] = static_cast<int>(x[dd] * recipcell_[dd][dd]);
445 int AnalysisNeighborhoodSearchImpl::getGridCellIndex(const ivec cell) const
447 GMX_ASSERT(cell[XX] >= 0 && cell[XX] < ncelldim_[XX],
448 "Grid cell X index out of range");
449 GMX_ASSERT(cell[YY] >= 0 && cell[YY] < ncelldim_[YY],
450 "Grid cell Y index out of range");
451 GMX_ASSERT(cell[ZZ] >= 0 && cell[ZZ] < ncelldim_[ZZ],
452 "Grid cell Z index out of range");
453 return cell[XX] + cell[YY] * ncelldim_[XX]
454 + cell[ZZ] * ncelldim_[XX] * ncelldim_[YY];
457 void AnalysisNeighborhoodSearchImpl::addToGridCell(const ivec cell, int i)
459 const int ci = getGridCellIndex(cell);
460 cells_[ci].push_back(i);
463 void AnalysisNeighborhoodSearchImpl::init(
464 AnalysisNeighborhood::SearchMode mode,
466 const AnalysisNeighborhoodPositions &positions)
468 GMX_RELEASE_ASSERT(positions.index_ == -1,
469 "Individual indexed positions not supported as reference");
470 pbc_ = const_cast<t_pbc *>(pbc);
471 nref_ = positions.count_;
472 // TODO: Consider whether it would be possible to support grid searching in
474 if (mode == AnalysisNeighborhood::eSearchMode_Simple
475 || pbc_ == NULL || pbc_->ePBC != epbcXYZ)
481 // TODO: Actually implement forcing eSearchMode_Grid
482 bGrid_ = initGrid(pbc_);
486 if (xref_nalloc_ < nref_)
488 srenew(xref_alloc_, nref_);
489 xref_nalloc_ = nref_;
493 for (int i = 0; i < nref_; ++i)
495 copy_rvec(positions.x_[i], xref_alloc_[i]);
497 put_atoms_in_triclinic_unitcell(ecenterTRIC, pbc_->box,
499 for (int i = 0; i < nref_; ++i)
503 mapPointToGridCell(xref_[i], refcell);
504 addToGridCell(refcell, i);
509 xref_ = positions.x_;
511 // TODO: Once exclusions are supported, this may need to be initialized.
517 * Sets the exclusions for the next neighborhood search.
519 * \param[in,out] d Neighborhood search data structure.
520 * \param[in] nexcl Number of reference positions to exclude from next
522 * \param[in] excl Indices of reference positions to exclude.
524 * The set exclusions remain in effect until the next call of this function.
527 gmx_ana_nbsearch_set_excl(gmx_ana_nbsearch_t *d, int nexcl, int excl[])
535 /********************************************************************
536 * AnalysisNeighborhoodPairSearchImpl
539 void AnalysisNeighborhoodPairSearchImpl::reset(int testIndex)
541 testIndex_ = testIndex;
542 if (testIndex_ >= 0 && testIndex_ < static_cast<int>(testPositions_.size()))
544 copy_rvec(testPositions_[testIndex_], xtest_);
547 put_atoms_in_triclinic_unitcell(ecenterTRIC, search_.pbc_->box,
549 search_.mapPointToGridCell(xtest_, testcell_);
558 void AnalysisNeighborhoodPairSearchImpl::nextTestPosition()
560 if (testIndex_ < static_cast<int>(testPositions_.size()))
567 bool AnalysisNeighborhoodPairSearchImpl::isExcluded(int j)
569 if (exclind_ < search_.nexcl_)
573 while (exclind_ < search_.nexcl_
574 && search_.excl_[exclind_] < search_.refid_[j])
578 if (exclind_ < search_.nexcl_
579 && search_.refid_[j] == search_.excl_[exclind_])
587 while (search_.bGrid_ && exclind_ < search_.nexcl_
588 && search_.excl_[exclind_] < j)
592 if (search_.excl_[exclind_] == j)
602 void AnalysisNeighborhoodPairSearchImpl::startSearch(
603 const AnalysisNeighborhoodPositions &positions)
605 if (positions.index_ < 0)
607 testPositions_ = ConstArrayRef<rvec>(positions.x_, positions.count_);
612 // Somewhat of a hack: setup the array such that only the last position
614 testPositions_ = ConstArrayRef<rvec>(positions.x_, positions.index_ + 1);
615 reset(positions.index_);
619 template <class Action>
620 bool AnalysisNeighborhoodPairSearchImpl::searchNext(Action action)
622 while (testIndex_ < static_cast<int>(testPositions_.size()))
627 int cai = prevcai_ + 1;
629 for (; nbi < search_.ngridnb_; ++nbi)
633 ivec_add(testcell_, search_.gnboffs_[nbi], cell);
634 cell[XX] = (cell[XX] + search_.ncelldim_[XX]) % search_.ncelldim_[XX];
635 cell[YY] = (cell[YY] + search_.ncelldim_[YY]) % search_.ncelldim_[YY];
636 cell[ZZ] = (cell[ZZ] + search_.ncelldim_[ZZ]) % search_.ncelldim_[ZZ];
638 const int ci = search_.getGridCellIndex(cell);
639 const int cellSize = static_cast<int>(search_.cells_[ci].size());
640 /* TODO: Calculate the required PBC shift outside the inner loop */
641 for (; cai < cellSize; ++cai)
643 const int i = search_.cells_[ci][cai];
649 pbc_dx_aiuc(search_.pbc_, xtest_, search_.xref_[i], dx);
650 const real r2 = norm2(dx);
651 if (r2 <= search_.cutoff2_)
668 for (int i = previ_ + 1; i < search_.nref_; ++i)
677 pbc_dx(search_.pbc_, xtest_, search_.xref_[i], dx);
681 rvec_sub(xtest_, search_.xref_[i], dx);
683 const real r2 = norm2(dx);
684 if (r2 <= search_.cutoff2_)
699 void AnalysisNeighborhoodPairSearchImpl::initFoundPair(
700 AnalysisNeighborhoodPair *pair) const
704 *pair = AnalysisNeighborhoodPair();
708 *pair = AnalysisNeighborhoodPair(previ_, testIndex_);
712 } // namespace internal
718 * Search action to find the next neighbor.
720 * Used as the action for AnalysisNeighborhoodPairSearchImpl::searchNext() to
721 * find the next neighbor.
723 * Simply breaks the loop on the first found neighbor.
725 bool withinAction(int /*i*/, real /*r2*/)
731 * Search action find the minimum distance.
733 * Used as the action for AnalysisNeighborhoodPairSearchImpl::searchNext() to
734 * find the nearest neighbor.
736 * With this action, AnalysisNeighborhoodPairSearchImpl::searchNext() always
737 * returns false, and the output is put into the variables passed by pointer
738 * into the constructor. If no neighbors are found, the output variables are
739 * not modified, i.e., the caller must initialize them.
745 * Initializes the action with given output locations.
747 * \param[out] closestPoint Index of the closest reference location.
748 * \param[out] minDist2 Minimum distance squared.
750 * The constructor call does not modify the pointed values, but only
751 * stores the pointers for later use.
752 * See the class description for additional semantics.
754 MindistAction(int *closestPoint, real *minDist2)
755 : closestPoint_(*closestPoint), minDist2_(*minDist2)
759 //! Processes a neighbor to find the nearest point.
760 bool operator()(int i, real r2)
774 GMX_DISALLOW_ASSIGN(MindistAction);
779 /********************************************************************
780 * AnalysisNeighborhood::Impl
783 class AnalysisNeighborhood::Impl
786 typedef AnalysisNeighborhoodSearch::ImplPointer SearchImplPointer;
787 typedef std::vector<SearchImplPointer> SearchList;
789 Impl() : cutoff_(0), mode_(eSearchMode_Automatic)
794 SearchList::const_iterator i;
795 for (i = searchList_.begin(); i != searchList_.end(); ++i)
797 GMX_RELEASE_ASSERT(i->unique(),
798 "Dangling AnalysisNeighborhoodSearch reference");
802 SearchImplPointer getSearch();
804 tMPI::mutex createSearchMutex_;
805 SearchList searchList_;
810 AnalysisNeighborhood::Impl::SearchImplPointer
811 AnalysisNeighborhood::Impl::getSearch()
813 tMPI::lock_guard<tMPI::mutex> lock(createSearchMutex_);
814 // TODO: Consider whether this needs to/can be faster, e.g., by keeping a
815 // separate pool of unused search objects.
816 SearchList::const_iterator i;
817 for (i = searchList_.begin(); i != searchList_.end(); ++i)
824 SearchImplPointer search(new internal::AnalysisNeighborhoodSearchImpl(cutoff_));
825 searchList_.push_back(search);
829 /********************************************************************
830 * AnalysisNeighborhood
833 AnalysisNeighborhood::AnalysisNeighborhood()
838 AnalysisNeighborhood::~AnalysisNeighborhood()
842 void AnalysisNeighborhood::setCutoff(real cutoff)
844 GMX_RELEASE_ASSERT(impl_->searchList_.empty(),
845 "Changing the cutoff after initSearch() not currently supported");
846 impl_->cutoff_ = cutoff;
849 void AnalysisNeighborhood::setMode(SearchMode mode)
854 AnalysisNeighborhood::SearchMode AnalysisNeighborhood::mode() const
859 AnalysisNeighborhoodSearch
860 AnalysisNeighborhood::initSearch(const t_pbc *pbc,
861 const AnalysisNeighborhoodPositions &positions)
863 Impl::SearchImplPointer search(impl_->getSearch());
864 search->init(mode(), pbc, positions);
865 return AnalysisNeighborhoodSearch(search);
868 /********************************************************************
869 * AnalysisNeighborhoodSearch
872 AnalysisNeighborhoodSearch::AnalysisNeighborhoodSearch()
876 AnalysisNeighborhoodSearch::AnalysisNeighborhoodSearch(const ImplPointer &impl)
881 void AnalysisNeighborhoodSearch::reset()
886 AnalysisNeighborhood::SearchMode AnalysisNeighborhoodSearch::mode() const
888 GMX_RELEASE_ASSERT(impl_, "Accessing an invalid search object");
889 return (impl_->usesGridSearch()
890 ? AnalysisNeighborhood::eSearchMode_Grid
891 : AnalysisNeighborhood::eSearchMode_Simple);
894 bool AnalysisNeighborhoodSearch::isWithin(
895 const AnalysisNeighborhoodPositions &positions) const
897 GMX_RELEASE_ASSERT(impl_, "Accessing an invalid search object");
898 internal::AnalysisNeighborhoodPairSearchImpl pairSearch(*impl_);
899 pairSearch.startSearch(positions);
900 return pairSearch.searchNext(&withinAction);
903 real AnalysisNeighborhoodSearch::minimumDistance(
904 const AnalysisNeighborhoodPositions &positions) const
906 GMX_RELEASE_ASSERT(impl_, "Accessing an invalid search object");
907 internal::AnalysisNeighborhoodPairSearchImpl pairSearch(*impl_);
908 pairSearch.startSearch(positions);
909 real minDist2 = impl_->cutoffSquared();
910 int closestPoint = -1;
911 MindistAction action(&closestPoint, &minDist2);
912 (void)pairSearch.searchNext(action);
913 return sqrt(minDist2);
916 AnalysisNeighborhoodPair
917 AnalysisNeighborhoodSearch::nearestPoint(
918 const AnalysisNeighborhoodPositions &positions) const
920 GMX_RELEASE_ASSERT(impl_, "Accessing an invalid search object");
921 internal::AnalysisNeighborhoodPairSearchImpl pairSearch(*impl_);
922 pairSearch.startSearch(positions);
923 real minDist2 = impl_->cutoffSquared();
924 int closestPoint = -1;
925 MindistAction action(&closestPoint, &minDist2);
926 (void)pairSearch.searchNext(action);
927 return AnalysisNeighborhoodPair(closestPoint, 0);
930 AnalysisNeighborhoodPairSearch
931 AnalysisNeighborhoodSearch::startPairSearch(
932 const AnalysisNeighborhoodPositions &positions) const
934 GMX_RELEASE_ASSERT(impl_, "Accessing an invalid search object");
935 Impl::PairSearchImplPointer pairSearch(impl_->getPairSearch());
936 pairSearch->startSearch(positions);
937 return AnalysisNeighborhoodPairSearch(pairSearch);
940 /********************************************************************
941 * AnalysisNeighborhoodPairSearch
944 AnalysisNeighborhoodPairSearch::AnalysisNeighborhoodPairSearch(
945 const ImplPointer &impl)
950 bool AnalysisNeighborhoodPairSearch::findNextPair(AnalysisNeighborhoodPair *pair)
952 bool bFound = impl_->searchNext(&withinAction);
953 impl_->initFoundPair(pair);
957 void AnalysisNeighborhoodPairSearch::skipRemainingPairsForTestPosition()
959 impl_->nextTestPosition();