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 300 : 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 300 : const unsigned& stride)
49 300 : : reduced(false),
50 300 : serial_(serial),
51 300 : do_pair_(do_pair),
52 300 : do_pbc_(do_pbc),
53 300 : twolists_(true),
54 300 : pbc_(&pbc),
55 300 : comm(cm),
56 : //copy-initialize fullatomlist_
57 300 : fullatomlist_(list0),
58 300 : distance_(distance),
59 300 : nlist0_(list0.size()),
60 300 : nlist1_(list1.size()),
61 300 : stride_(stride) {
62 : // store the rest of the atoms into fullatomlist_
63 301 : fullatomlist_.insert(fullatomlist_.end(),list1.begin(),list1.end());
64 300 : if(!do_pair) {
65 271 : nallpairs_=nlist0_*nlist1_;
66 : } else {
67 29 : 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 29 : nallpairs_=nlist0_;
72 : }
73 300 : initialize();
74 299 : }
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 342 : NeighborList::~NeighborList()=default;
98 :
99 344 : 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 344 : 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 228150101 : for(unsigned int i=0; i<nallpairs_; ++i)
123 228149759 : neighbors_[i]=getIndexPair(i);
124 342 : }
125 :
126 70008 : std::vector<AtomNumber>& NeighborList::getFullAtomList() {
127 70008 : return fullatomlist_;
128 : }
129 :
130 234387479 : NeighborList::pairIDs NeighborList::getIndexPair(const unsigned ipair) {
131 : pairIDs index;
132 234387479 : if(twolists_ && do_pair_)
133 12236 : index=pairIDs(ipair,ipair+nlist0_);
134 234375243 : else if (twolists_ && !do_pair_)
135 143426597 : index=pairIDs(ipair/nlist1_,ipair%nlist1_+nlist0_);
136 90948646 : else if (!twolists_) {
137 90948646 : unsigned ii = nallpairs_-1-ipair;
138 90948646 : unsigned K = unsigned(std::floor((std::sqrt(double(8*ii+1))+1)/2));
139 90948646 : unsigned jj = ii-K*(K-1)/2;
140 90948646 : index=pairIDs(nlist0_-1-K,nlist0_-1-jj);
141 : }
142 234387479 : return index;
143 : }
144 :
145 3138 : void NeighborList::update(const std::vector<Vector>& positions) {
146 3138 : neighbors_.clear();
147 3138 : const double d2=distance_*distance_;
148 : // check if positions array has the correct length
149 3138 : plumed_assert(positions.size()==fullatomlist_.size());
150 :
151 3138 : unsigned stride=comm.Get_size();
152 3138 : unsigned rank=comm.Get_rank();
153 3138 : unsigned nt=OpenMP::getNumThreads();
154 3138 : if(serial_) {
155 : stride=1;
156 : rank=0;
157 : nt=1;
158 : }
159 : std::vector<unsigned> local_flat_nl;
160 :
161 3138 : #pragma omp parallel num_threads(nt)
162 : {
163 : std::vector<unsigned> private_flat_nl;
164 : #pragma omp for nowait
165 : for(unsigned int i=rank; i<nallpairs_; i+=stride) {
166 : pairIDs index=getIndexPair(i);
167 : unsigned index0=index.first;
168 : unsigned index1=index.second;
169 : Vector distance;
170 : if(do_pbc_) {
171 : distance=pbc_->distance(positions[index0],positions[index1]);
172 : } else {
173 : distance=delta(positions[index0],positions[index1]);
174 : }
175 : double value=modulo2(distance);
176 : if(value<=d2) {
177 : private_flat_nl.push_back(index0);
178 : private_flat_nl.push_back(index1);
179 : }
180 : }
181 : #pragma omp critical
182 : local_flat_nl.insert(local_flat_nl.end(),
183 : private_flat_nl.begin(),
184 : private_flat_nl.end());
185 : }
186 :
187 : // find total dimension of neighborlist
188 3138 : std::vector <int> local_nl_size(stride, 0);
189 3138 : local_nl_size[rank] = local_flat_nl.size();
190 3138 : if(!serial_)
191 3126 : comm.Sum(&local_nl_size[0], stride);
192 : int tot_size = std::accumulate(local_nl_size.begin(), local_nl_size.end(), 0);
193 3138 : if(tot_size==0) {
194 18 : setRequestList();
195 : return;
196 : }
197 : // merge
198 3120 : std::vector<unsigned> merge_nl(tot_size, 0);
199 : // calculate vector of displacement
200 3120 : std::vector<int> disp(stride);
201 3120 : disp[0] = 0;
202 : int rank_size = 0;
203 9216 : for(unsigned i=0; i<stride-1; ++i) {
204 6096 : rank_size += local_nl_size[i];
205 6096 : disp[i+1] = rank_size;
206 : }
207 : // Allgather neighbor list
208 3120 : if(comm.initialized()&&!serial_) {
209 4181 : comm.Allgatherv((!local_flat_nl.empty()?&local_flat_nl[0]:NULL),
210 : local_nl_size[rank],
211 : &merge_nl[0],
212 : &local_nl_size[0],
213 : &disp[0]);
214 : } else
215 1016 : merge_nl = local_flat_nl;
216 : // resize neighbor stuff
217 3120 : neighbors_.resize(tot_size/2);
218 6114036 : for(unsigned i=0; i<tot_size/2; i++) {
219 6110916 : unsigned j=2*i;
220 6110916 : neighbors_[i] = std::make_pair(merge_nl[j],merge_nl[j+1]);
221 : }
222 :
223 3120 : setRequestList();
224 : }
225 :
226 3138 : void NeighborList::setRequestList() {
227 3138 : requestlist_.clear();
228 6114054 : for(unsigned int i=0; i<size(); ++i) {
229 6110916 : requestlist_.push_back(fullatomlist_[neighbors_[i].first]);
230 6110916 : requestlist_.push_back(fullatomlist_[neighbors_[i].second]);
231 : }
232 3138 : Tools::removeDuplicates(requestlist_);
233 3138 : reduced=false;
234 3138 : }
235 :
236 240 : std::vector<AtomNumber>& NeighborList::getReducedAtomList() {
237 240 : if(!reduced) {
238 168162 : for(unsigned int i=0; i<size(); ++i) {
239 168119 : AtomNumber index0=fullatomlist_[neighbors_[i].first];
240 168119 : AtomNumber index1=fullatomlist_[neighbors_[i].second];
241 : // I exploit the fact that requestlist_ is an ordered vector
242 168119 : auto p = std::find(requestlist_.begin(), requestlist_.end(), index0);
243 : plumed_dbg_assert(p!=requestlist_.end());
244 168119 : unsigned newindex0=p-requestlist_.begin();
245 168119 : p = std::find(requestlist_.begin(), requestlist_.end(), index1);
246 : plumed_dbg_assert(p!=requestlist_.end());
247 168119 : unsigned newindex1=p-requestlist_.begin();
248 : neighbors_[i]=pairIDs(newindex0,newindex1);
249 : }
250 : }
251 240 : reduced=true;
252 240 : return requestlist_;
253 : }
254 :
255 14950 : unsigned NeighborList::getStride() const {
256 14950 : return stride_;
257 : }
258 :
259 24 : unsigned NeighborList::getLastUpdate() const {
260 24 : return lastupdate_;
261 : }
262 :
263 0 : void NeighborList::setLastUpdate(const unsigned step) {
264 0 : lastupdate_=step;
265 0 : }
266 :
267 23715507 : unsigned NeighborList::size() const {
268 23715507 : return neighbors_.size();
269 : }
270 :
271 267302312 : NeighborList::pairIDs NeighborList::getClosePair(const unsigned i) const {
272 267302312 : return neighbors_[i];
273 : }
274 :
275 : NeighborList::pairAtomNumbers
276 34665936 : NeighborList::getClosePairAtomNumber(const unsigned i) const {
277 : pairAtomNumbers Aneigh=pairAtomNumbers(
278 34665936 : fullatomlist_[neighbors_[i].first],
279 34665936 : fullatomlist_[neighbors_[i].second]);
280 34665936 : return Aneigh;
281 : }
282 :
283 0 : std::vector<unsigned> NeighborList::getNeighbors(const unsigned index) {
284 : std::vector<unsigned> neighbors;
285 0 : for(unsigned int i=0; i<size(); ++i) {
286 0 : if(neighbors_[i].first==index)
287 0 : neighbors.push_back(neighbors_[i].second);
288 0 : if(neighbors_[i].second==index)
289 0 : neighbors.push_back(neighbors_[i].first);
290 : }
291 0 : return neighbors;
292 : }
293 :
294 : } // namespace PLMD
|