Line data Source code
1 : /* +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
2 : Copyright (c) 2011-2023 The plumed team
3 : (see the PEOPLE file at the root of the distribution for a list of names)
4 :
5 : See http://www.plumed.org for more information.
6 :
7 : This file is part of plumed, version 2.
8 :
9 : plumed is free software: you can redistribute it and/or modify
10 : it under the terms of the GNU Lesser General Public License as published by
11 : the Free Software Foundation, either version 3 of the License, or
12 : (at your option) any later version.
13 :
14 : plumed 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
17 : GNU Lesser General Public License for more details.
18 :
19 : You should have received a copy of the GNU Lesser General Public License
20 : along with plumed. If not, see <http://www.gnu.org/licenses/>.
21 : +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ */
22 : #include "NeighborList.h"
23 : #include "Vector.h"
24 : #include "Pbc.h"
25 : #include "AtomNumber.h"
26 : #include "Communicator.h"
27 : #include "OpenMP.h"
28 : #include "Tools.h"
29 : #include <vector>
30 : #include <algorithm>
31 : #include <numeric>
32 :
33 : #ifdef __APPLE__
34 : //we are using getenv to give the user the opporunity of suppressing
35 : //the too many memory killswitch while compiling on mac
36 : #include <cstdlib>
37 : #endif //__APPLE__
38 : namespace PLMD {
39 :
40 310 : NeighborList::NeighborList(const std::vector<AtomNumber>& list0,
41 : const std::vector<AtomNumber>& list1,
42 : const bool& serial,
43 : const bool& do_pair,
44 : const bool& do_pbc,
45 : const Pbc& pbc,
46 : Communicator& cm,
47 : const double& distance,
48 310 : const unsigned& stride)
49 310 : : reduced(false),
50 310 : serial_(serial),
51 310 : do_pair_(do_pair),
52 310 : do_pbc_(do_pbc),
53 310 : twolists_(true),
54 310 : pbc_(&pbc),
55 310 : comm(cm),
56 : //copy-initialize fullatomlist_
57 310 : fullatomlist_(list0),
58 310 : distance_(distance),
59 310 : nlist0_(list0.size()),
60 310 : nlist1_(list1.size()),
61 310 : stride_(stride) {
62 : // store the rest of the atoms into fullatomlist_
63 311 : fullatomlist_.insert(fullatomlist_.end(),list1.begin(),list1.end());
64 310 : if(!do_pair) {
65 271 : nallpairs_=nlist0_*nlist1_;
66 : } else {
67 39 : plumed_assert(nlist0_==nlist1_)
68 : << "when using PAIR option, the two groups should have the same number"
69 : " of elements\n" << "the groups you specified have size "
70 0 : <<nlist0_<<" and "<<nlist1_;
71 39 : nallpairs_=nlist0_;
72 : }
73 310 : initialize();
74 309 : }
75 :
76 44 : NeighborList::NeighborList(const std::vector<AtomNumber>& list0,
77 : const bool& serial, const bool& do_pbc,
78 : const Pbc& pbc,
79 : Communicator& cm,
80 : const double& distance,
81 44 : const unsigned& stride)
82 44 : : reduced(false),
83 44 : serial_(serial),
84 44 : do_pbc_(do_pbc),
85 44 : twolists_(false),
86 44 : pbc_(&pbc),
87 44 : comm(cm),
88 : //copy-initialize fullatomlist_
89 44 : fullatomlist_(list0),
90 44 : distance_(distance),
91 44 : nlist0_(list0.size()),
92 44 : nallpairs_(nlist0_*(nlist0_-1)/2),
93 44 : stride_(stride) {
94 44 : initialize();
95 43 : }
96 :
97 352 : NeighborList::~NeighborList()=default;
98 :
99 354 : void NeighborList::initialize() {
100 : #ifdef __APPLE__
101 : //this mac-only error is here because on my experience the mac tries to page
102 : //the memory on the hdd instead of throwing a memory error
103 : if(!std::getenv("PLUMED_IGNORE_NL_MEMORY_ERROR")) {
104 : //blocking memory allocation on slightly more than 10 GB of memory
105 : //that is about 1296000000 pairs (36000 atoms)
106 : //36000 * 36000= 1296000000
107 : //each pairIDs occupies 64 bit (where unsigned are 32bit integers)
108 : //4294967296 is max(uint32)+1 and is more than 34 GB (correspond to a system of 65536 atoms)
109 : if(nallpairs_ > 1296000000 )
110 : plumed_merror("An error happened while allocating the neighbor "
111 : "list, please decrease the number of atoms used");
112 : }
113 : #endif // __APPLE__
114 : try {
115 354 : neighbors_.resize(nallpairs_);
116 2 : } catch (...) {
117 4 : plumed_error_nested() << "An error happened while allocating the neighbor "
118 : "list, please decrease the number of atoms used";
119 2 : }
120 : //TODO: test if this is feasible for accelerating the loop
121 : //#pragma omp parallel for default(shared)
122 228150421 : for(unsigned int i=0; i<nallpairs_; ++i) {
123 228150069 : neighbors_[i]=getIndexPair(i);
124 : }
125 352 : }
126 :
127 70018 : std::vector<AtomNumber>& NeighborList::getFullAtomList() {
128 70018 : return fullatomlist_;
129 : }
130 :
131 234387789 : NeighborList::pairIDs NeighborList::getIndexPair(const unsigned ipair) {
132 : pairIDs index;
133 234387789 : if(twolists_ && do_pair_) {
134 12546 : index=pairIDs(ipair,ipair+nlist0_);
135 234375243 : } else if (twolists_ && !do_pair_) {
136 143426597 : index=pairIDs(ipair/nlist1_,ipair%nlist1_+nlist0_);
137 90948646 : } else if (!twolists_) {
138 90948646 : unsigned ii = nallpairs_-1-ipair;
139 90948646 : unsigned K = unsigned(std::floor((std::sqrt(double(8*ii+1))+1)/2));
140 90948646 : unsigned jj = ii-K*(K-1)/2;
141 90948646 : index=pairIDs(nlist0_-1-K,nlist0_-1-jj);
142 : }
143 234387789 : return index;
144 : }
145 :
146 3138 : void NeighborList::update(const std::vector<Vector>& positions) {
147 3138 : neighbors_.clear();
148 3138 : const double d2=distance_*distance_;
149 : // check if positions array has the correct length
150 3138 : plumed_assert(positions.size()==fullatomlist_.size());
151 :
152 3138 : unsigned stride=comm.Get_size();
153 3138 : unsigned rank=comm.Get_rank();
154 3138 : unsigned nt=OpenMP::getNumThreads();
155 3138 : if(serial_) {
156 : stride=1;
157 : rank=0;
158 : nt=1;
159 : }
160 : std::vector<unsigned> local_flat_nl;
161 :
162 3138 : #pragma omp parallel num_threads(nt)
163 : {
164 : std::vector<unsigned> private_flat_nl;
165 : #pragma omp for nowait
166 : for(unsigned int i=rank; i<nallpairs_; i+=stride) {
167 : pairIDs index=getIndexPair(i);
168 : unsigned index0=index.first;
169 : unsigned index1=index.second;
170 : Vector distance;
171 : if(do_pbc_) {
172 : distance=pbc_->distance(positions[index0],positions[index1]);
173 : } else {
174 : distance=delta(positions[index0],positions[index1]);
175 : }
176 : double value=modulo2(distance);
177 : if(value<=d2) {
178 : private_flat_nl.push_back(index0);
179 : private_flat_nl.push_back(index1);
180 : }
181 : }
182 : #pragma omp critical
183 : local_flat_nl.insert(local_flat_nl.end(),
184 : private_flat_nl.begin(),
185 : private_flat_nl.end());
186 : }
187 :
188 : // find total dimension of neighborlist
189 3138 : std::vector <int> local_nl_size(stride, 0);
190 3138 : local_nl_size[rank] = local_flat_nl.size();
191 3138 : if(!serial_) {
192 3126 : comm.Sum(&local_nl_size[0], stride);
193 : }
194 : int tot_size = std::accumulate(local_nl_size.begin(), local_nl_size.end(), 0);
195 3138 : if(tot_size==0) {
196 18 : setRequestList();
197 : return;
198 : }
199 : // merge
200 3120 : std::vector<unsigned> merge_nl(tot_size, 0);
201 : // calculate vector of displacement
202 3120 : std::vector<int> disp(stride);
203 3120 : disp[0] = 0;
204 : int rank_size = 0;
205 9216 : for(unsigned i=0; i<stride-1; ++i) {
206 6096 : rank_size += local_nl_size[i];
207 6096 : disp[i+1] = rank_size;
208 : }
209 : // Allgather neighbor list
210 3120 : if(comm.initialized()&&!serial_) {
211 4181 : comm.Allgatherv((!local_flat_nl.empty()?&local_flat_nl[0]:NULL),
212 : local_nl_size[rank],
213 : &merge_nl[0],
214 : &local_nl_size[0],
215 : &disp[0]);
216 : } else {
217 1016 : merge_nl = local_flat_nl;
218 : }
219 : // resize neighbor stuff
220 3120 : neighbors_.resize(tot_size/2);
221 6114036 : for(unsigned i=0; i<tot_size/2; i++) {
222 6110916 : unsigned j=2*i;
223 6110916 : neighbors_[i] = std::make_pair(merge_nl[j],merge_nl[j+1]);
224 : }
225 :
226 3120 : setRequestList();
227 : }
228 :
229 3138 : void NeighborList::setRequestList() {
230 3138 : requestlist_.clear();
231 6114054 : for(unsigned int i=0; i<size(); ++i) {
232 6110916 : requestlist_.push_back(fullatomlist_[neighbors_[i].first]);
233 6110916 : requestlist_.push_back(fullatomlist_[neighbors_[i].second]);
234 : }
235 3138 : Tools::removeDuplicates(requestlist_);
236 3138 : reduced=false;
237 3138 : }
238 :
239 240 : std::vector<AtomNumber>& NeighborList::getReducedAtomList() {
240 240 : if(!reduced) {
241 168162 : for(unsigned int i=0; i<size(); ++i) {
242 168119 : AtomNumber index0=fullatomlist_[neighbors_[i].first];
243 168119 : AtomNumber index1=fullatomlist_[neighbors_[i].second];
244 : // I exploit the fact that requestlist_ is an ordered vector
245 168119 : auto p = std::find(requestlist_.begin(), requestlist_.end(), index0);
246 : plumed_dbg_assert(p!=requestlist_.end());
247 168119 : unsigned newindex0=p-requestlist_.begin();
248 168119 : p = std::find(requestlist_.begin(), requestlist_.end(), index1);
249 : plumed_dbg_assert(p!=requestlist_.end());
250 168119 : unsigned newindex1=p-requestlist_.begin();
251 : neighbors_[i]=pairIDs(newindex0,newindex1);
252 : }
253 : }
254 240 : reduced=true;
255 240 : return requestlist_;
256 : }
257 :
258 14950 : unsigned NeighborList::getStride() const {
259 14950 : return stride_;
260 : }
261 :
262 24 : unsigned NeighborList::getLastUpdate() const {
263 24 : return lastupdate_;
264 : }
265 :
266 0 : void NeighborList::setLastUpdate(const unsigned step) {
267 0 : lastupdate_=step;
268 0 : }
269 :
270 23863798 : unsigned NeighborList::size() const {
271 23863798 : return neighbors_.size();
272 : }
273 :
274 267279536 : NeighborList::pairIDs NeighborList::getClosePair(const unsigned i) const {
275 267279536 : return neighbors_[i];
276 : }
277 :
278 : NeighborList::pairAtomNumbers
279 34665936 : NeighborList::getClosePairAtomNumber(const unsigned i) const {
280 : pairAtomNumbers Aneigh=pairAtomNumbers(
281 34665936 : fullatomlist_[neighbors_[i].first],
282 34665936 : fullatomlist_[neighbors_[i].second]);
283 34665936 : return Aneigh;
284 : }
285 :
286 0 : std::vector<unsigned> NeighborList::getNeighbors(const unsigned index) {
287 : std::vector<unsigned> neighbors;
288 0 : for(unsigned int i=0; i<size(); ++i) {
289 0 : if(neighbors_[i].first==index) {
290 0 : neighbors.push_back(neighbors_[i].second);
291 : }
292 0 : if(neighbors_[i].second==index) {
293 0 : neighbors.push_back(neighbors_[i].first);
294 : }
295 : }
296 0 : return neighbors;
297 : }
298 :
299 : } // namespace PLMD
|