ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/parallel/ForceMatrixDecomposition.cpp
Revision: 1590
Committed: Mon Jul 11 01:39:49 2011 UTC (13 years, 9 months ago) by gezelter
Original Path: branches/development/src/parallel/ForceMatrixDecomposition.cpp
File size: 39310 byte(s)
Log Message:
more bug fixes

File Contents

# User Rev Content
1 gezelter 1539 /*
2     * Copyright (c) 2005 The University of Notre Dame. All Rights Reserved.
3 chuckv 1538 *
4     * The University of Notre Dame grants you ("Licensee") a
5     * non-exclusive, royalty free, license to use, modify and
6     * redistribute this software in source and binary code form, provided
7     * that the following conditions are met:
8     *
9     * 1. Redistributions of source code must retain the above copyright
10     * notice, this list of conditions and the following disclaimer.
11     *
12     * 2. Redistributions in binary form must reproduce the above copyright
13     * notice, this list of conditions and the following disclaimer in the
14     * documentation and/or other materials provided with the
15     * distribution.
16     *
17     * This software is provided "AS IS," without a warranty of any
18     * kind. All express or implied conditions, representations and
19     * warranties, including any implied warranty of merchantability,
20     * fitness for a particular purpose or non-infringement, are hereby
21     * excluded. The University of Notre Dame and its licensors shall not
22     * be liable for any damages suffered by licensee as a result of
23     * using, modifying or distributing the software or its
24     * derivatives. In no event will the University of Notre Dame or its
25     * licensors be liable for any lost revenue, profit or data, or for
26     * direct, indirect, special, consequential, incidental or punitive
27     * damages, however caused and regardless of the theory of liability,
28     * arising out of the use of or inability to use software, even if the
29     * University of Notre Dame has been advised of the possibility of
30     * such damages.
31     *
32     * SUPPORT OPEN SCIENCE! If you use OpenMD or its source code in your
33     * research, please cite the appropriate papers when you publish your
34     * work. Good starting points are:
35     *
36     * [1] Meineke, et al., J. Comp. Chem. 26, 252-271 (2005).
37     * [2] Fennell & Gezelter, J. Chem. Phys. 124, 234104 (2006).
38     * [3] Sun, Lin & Gezelter, J. Chem. Phys. 128, 24107 (2008).
39     * [4] Vardeman & Gezelter, in progress (2009).
40     */
41 gezelter 1549 #include "parallel/ForceMatrixDecomposition.hpp"
42 gezelter 1539 #include "math/SquareMatrix3.hpp"
43 gezelter 1544 #include "nonbonded/NonBondedInteraction.hpp"
44     #include "brains/SnapshotManager.hpp"
45 gezelter 1570 #include "brains/PairList.hpp"
46 chuckv 1538
47 gezelter 1541 using namespace std;
48 gezelter 1539 namespace OpenMD {
49 chuckv 1538
50 gezelter 1544 /**
51     * distributeInitialData is essentially a copy of the older fortran
52     * SimulationSetup
53     */
54    
55 gezelter 1549 void ForceMatrixDecomposition::distributeInitialData() {
56 gezelter 1551 snap_ = sman_->getCurrentSnapshot();
57     storageLayout_ = sman_->getStorageLayout();
58 gezelter 1571 ff_ = info_->getForceField();
59 gezelter 1567 nLocal_ = snap_->getNumberOfAtoms();
60 gezelter 1587
61 gezelter 1577 nGroups_ = info_->getNLocalCutoffGroups();
62 gezelter 1569 // gather the information for atomtype IDs (atids):
63 gezelter 1583 idents = info_->getIdentArray();
64 gezelter 1569 AtomLocalToGlobal = info_->getGlobalAtomIndices();
65     cgLocalToGlobal = info_->getGlobalGroupIndices();
66     vector<int> globalGroupMembership = info_->getGlobalGroupMembership();
67 gezelter 1586
68 gezelter 1581 massFactors = info_->getMassFactors();
69 gezelter 1584
70 gezelter 1587 PairList* excludes = info_->getExcludedInteractions();
71     PairList* oneTwo = info_->getOneTwoInteractions();
72     PairList* oneThree = info_->getOneThreeInteractions();
73     PairList* oneFour = info_->getOneFourInteractions();
74 gezelter 1569
75 gezelter 1567 #ifdef IS_MPI
76    
77     AtomCommIntRow = new Communicator<Row,int>(nLocal_);
78     AtomCommRealRow = new Communicator<Row,RealType>(nLocal_);
79     AtomCommVectorRow = new Communicator<Row,Vector3d>(nLocal_);
80     AtomCommMatrixRow = new Communicator<Row,Mat3x3d>(nLocal_);
81 gezelter 1575 AtomCommPotRow = new Communicator<Row,potVec>(nLocal_);
82 chuckv 1538
83 gezelter 1567 AtomCommIntColumn = new Communicator<Column,int>(nLocal_);
84     AtomCommRealColumn = new Communicator<Column,RealType>(nLocal_);
85     AtomCommVectorColumn = new Communicator<Column,Vector3d>(nLocal_);
86     AtomCommMatrixColumn = new Communicator<Column,Mat3x3d>(nLocal_);
87 gezelter 1575 AtomCommPotColumn = new Communicator<Column,potVec>(nLocal_);
88 gezelter 1541
89 gezelter 1567 cgCommIntRow = new Communicator<Row,int>(nGroups_);
90     cgCommVectorRow = new Communicator<Row,Vector3d>(nGroups_);
91     cgCommIntColumn = new Communicator<Column,int>(nGroups_);
92     cgCommVectorColumn = new Communicator<Column,Vector3d>(nGroups_);
93 gezelter 1551
94 gezelter 1567 nAtomsInRow_ = AtomCommIntRow->getSize();
95     nAtomsInCol_ = AtomCommIntColumn->getSize();
96     nGroupsInRow_ = cgCommIntRow->getSize();
97     nGroupsInCol_ = cgCommIntColumn->getSize();
98    
99 gezelter 1551 // Modify the data storage objects with the correct layouts and sizes:
100 gezelter 1567 atomRowData.resize(nAtomsInRow_);
101 gezelter 1551 atomRowData.setStorageLayout(storageLayout_);
102 gezelter 1567 atomColData.resize(nAtomsInCol_);
103 gezelter 1551 atomColData.setStorageLayout(storageLayout_);
104 gezelter 1567 cgRowData.resize(nGroupsInRow_);
105 gezelter 1551 cgRowData.setStorageLayout(DataStorage::dslPosition);
106 gezelter 1567 cgColData.resize(nGroupsInCol_);
107 gezelter 1551 cgColData.setStorageLayout(DataStorage::dslPosition);
108 gezelter 1575
109 gezelter 1577 identsRow.resize(nAtomsInRow_);
110     identsCol.resize(nAtomsInCol_);
111 gezelter 1549
112 gezelter 1583 AtomCommIntRow->gather(idents, identsRow);
113     AtomCommIntColumn->gather(idents, identsCol);
114 gezelter 1549
115 gezelter 1589 // allocate memory for the parallel objects
116     AtomRowToGlobal.resize(nAtomsInRow_);
117     AtomColToGlobal.resize(nAtomsInCol_);
118     cgRowToGlobal.resize(nGroupsInRow_);
119     cgColToGlobal.resize(nGroupsInCol_);
120     massFactorsRow.resize(nAtomsInRow_);
121     massFactorsCol.resize(nAtomsInCol_);
122     pot_row.resize(nAtomsInRow_);
123     pot_col.resize(nAtomsInCol_);
124    
125 gezelter 1549 AtomCommIntRow->gather(AtomLocalToGlobal, AtomRowToGlobal);
126     AtomCommIntColumn->gather(AtomLocalToGlobal, AtomColToGlobal);
127    
128     cgCommIntRow->gather(cgLocalToGlobal, cgRowToGlobal);
129     cgCommIntColumn->gather(cgLocalToGlobal, cgColToGlobal);
130 gezelter 1541
131 gezelter 1581 AtomCommRealRow->gather(massFactors, massFactorsRow);
132     AtomCommRealColumn->gather(massFactors, massFactorsCol);
133 gezelter 1569
134     groupListRow_.clear();
135 gezelter 1577 groupListRow_.resize(nGroupsInRow_);
136 gezelter 1569 for (int i = 0; i < nGroupsInRow_; i++) {
137     int gid = cgRowToGlobal[i];
138     for (int j = 0; j < nAtomsInRow_; j++) {
139     int aid = AtomRowToGlobal[j];
140     if (globalGroupMembership[aid] == gid)
141     groupListRow_[i].push_back(j);
142     }
143     }
144    
145     groupListCol_.clear();
146 gezelter 1577 groupListCol_.resize(nGroupsInCol_);
147 gezelter 1569 for (int i = 0; i < nGroupsInCol_; i++) {
148     int gid = cgColToGlobal[i];
149     for (int j = 0; j < nAtomsInCol_; j++) {
150     int aid = AtomColToGlobal[j];
151     if (globalGroupMembership[aid] == gid)
152     groupListCol_[i].push_back(j);
153     }
154     }
155    
156 gezelter 1587 excludesForAtom.clear();
157     excludesForAtom.resize(nAtomsInRow_);
158 gezelter 1579 toposForAtom.clear();
159     toposForAtom.resize(nAtomsInRow_);
160     topoDist.clear();
161     topoDist.resize(nAtomsInRow_);
162 gezelter 1570 for (int i = 0; i < nAtomsInRow_; i++) {
163 gezelter 1571 int iglob = AtomRowToGlobal[i];
164 gezelter 1579
165 gezelter 1570 for (int j = 0; j < nAtomsInCol_; j++) {
166 gezelter 1579 int jglob = AtomColToGlobal[j];
167    
168 gezelter 1587 if (excludes->hasPair(iglob, jglob))
169     excludesForAtom[i].push_back(j);
170 gezelter 1579
171 gezelter 1587 if (oneTwo->hasPair(iglob, jglob)) {
172 gezelter 1579 toposForAtom[i].push_back(j);
173     topoDist[i].push_back(1);
174     } else {
175 gezelter 1587 if (oneThree->hasPair(iglob, jglob)) {
176 gezelter 1579 toposForAtom[i].push_back(j);
177     topoDist[i].push_back(2);
178     } else {
179 gezelter 1587 if (oneFour->hasPair(iglob, jglob)) {
180 gezelter 1579 toposForAtom[i].push_back(j);
181     topoDist[i].push_back(3);
182     }
183     }
184 gezelter 1570 }
185     }
186     }
187    
188 gezelter 1569 #endif
189 gezelter 1579
190 gezelter 1569 groupList_.clear();
191 gezelter 1577 groupList_.resize(nGroups_);
192 gezelter 1569 for (int i = 0; i < nGroups_; i++) {
193     int gid = cgLocalToGlobal[i];
194     for (int j = 0; j < nLocal_; j++) {
195     int aid = AtomLocalToGlobal[j];
196 gezelter 1577 if (globalGroupMembership[aid] == gid) {
197 gezelter 1569 groupList_[i].push_back(j);
198 gezelter 1577 }
199 gezelter 1569 }
200     }
201    
202 gezelter 1587 excludesForAtom.clear();
203     excludesForAtom.resize(nLocal_);
204 gezelter 1579 toposForAtom.clear();
205     toposForAtom.resize(nLocal_);
206     topoDist.clear();
207     topoDist.resize(nLocal_);
208 gezelter 1569
209 gezelter 1570 for (int i = 0; i < nLocal_; i++) {
210     int iglob = AtomLocalToGlobal[i];
211 gezelter 1579
212 gezelter 1570 for (int j = 0; j < nLocal_; j++) {
213 gezelter 1579 int jglob = AtomLocalToGlobal[j];
214    
215 gezelter 1587 if (excludes->hasPair(iglob, jglob))
216     excludesForAtom[i].push_back(j);
217 gezelter 1579
218 gezelter 1587 if (oneTwo->hasPair(iglob, jglob)) {
219 gezelter 1579 toposForAtom[i].push_back(j);
220     topoDist[i].push_back(1);
221     } else {
222 gezelter 1587 if (oneThree->hasPair(iglob, jglob)) {
223 gezelter 1579 toposForAtom[i].push_back(j);
224     topoDist[i].push_back(2);
225     } else {
226 gezelter 1587 if (oneFour->hasPair(iglob, jglob)) {
227 gezelter 1579 toposForAtom[i].push_back(j);
228     topoDist[i].push_back(3);
229     }
230     }
231 gezelter 1570 }
232     }
233 gezelter 1579 }
234    
235     createGtypeCutoffMap();
236 gezelter 1587
237 gezelter 1576 }
238    
239     void ForceMatrixDecomposition::createGtypeCutoffMap() {
240 gezelter 1586
241 gezelter 1576 RealType tol = 1e-6;
242     RealType rc;
243     int atid;
244     set<AtomType*> atypes = info_->getSimulatedAtomTypes();
245 gezelter 1587 map<int, RealType> atypeCutoff;
246 gezelter 1583
247 gezelter 1579 for (set<AtomType*>::iterator at = atypes.begin();
248     at != atypes.end(); ++at){
249 gezelter 1576 atid = (*at)->getIdent();
250 gezelter 1587 if (userChoseCutoff_)
251 gezelter 1583 atypeCutoff[atid] = userCutoff_;
252 gezelter 1587 else
253 gezelter 1583 atypeCutoff[atid] = interactionMan_->getSuggestedCutoffRadius(*at);
254 gezelter 1570 }
255 gezelter 1576
256     vector<RealType> gTypeCutoffs;
257     // first we do a single loop over the cutoff groups to find the
258     // largest cutoff for any atypes present in this group.
259     #ifdef IS_MPI
260     vector<RealType> groupCutoffRow(nGroupsInRow_, 0.0);
261 gezelter 1579 groupRowToGtype.resize(nGroupsInRow_);
262 gezelter 1576 for (int cg1 = 0; cg1 < nGroupsInRow_; cg1++) {
263     vector<int> atomListRow = getAtomsInGroupRow(cg1);
264     for (vector<int>::iterator ia = atomListRow.begin();
265     ia != atomListRow.end(); ++ia) {
266     int atom1 = (*ia);
267     atid = identsRow[atom1];
268     if (atypeCutoff[atid] > groupCutoffRow[cg1]) {
269     groupCutoffRow[cg1] = atypeCutoff[atid];
270     }
271     }
272    
273     bool gTypeFound = false;
274     for (int gt = 0; gt < gTypeCutoffs.size(); gt++) {
275     if (abs(groupCutoffRow[cg1] - gTypeCutoffs[gt]) < tol) {
276     groupRowToGtype[cg1] = gt;
277     gTypeFound = true;
278     }
279     }
280     if (!gTypeFound) {
281     gTypeCutoffs.push_back( groupCutoffRow[cg1] );
282     groupRowToGtype[cg1] = gTypeCutoffs.size() - 1;
283     }
284    
285     }
286     vector<RealType> groupCutoffCol(nGroupsInCol_, 0.0);
287 gezelter 1579 groupColToGtype.resize(nGroupsInCol_);
288 gezelter 1576 for (int cg2 = 0; cg2 < nGroupsInCol_; cg2++) {
289     vector<int> atomListCol = getAtomsInGroupColumn(cg2);
290     for (vector<int>::iterator jb = atomListCol.begin();
291     jb != atomListCol.end(); ++jb) {
292     int atom2 = (*jb);
293     atid = identsCol[atom2];
294     if (atypeCutoff[atid] > groupCutoffCol[cg2]) {
295     groupCutoffCol[cg2] = atypeCutoff[atid];
296     }
297     }
298     bool gTypeFound = false;
299     for (int gt = 0; gt < gTypeCutoffs.size(); gt++) {
300     if (abs(groupCutoffCol[cg2] - gTypeCutoffs[gt]) < tol) {
301     groupColToGtype[cg2] = gt;
302     gTypeFound = true;
303     }
304     }
305     if (!gTypeFound) {
306     gTypeCutoffs.push_back( groupCutoffCol[cg2] );
307     groupColToGtype[cg2] = gTypeCutoffs.size() - 1;
308     }
309     }
310     #else
311 gezelter 1579
312 gezelter 1576 vector<RealType> groupCutoff(nGroups_, 0.0);
313 gezelter 1579 groupToGtype.resize(nGroups_);
314 gezelter 1576 for (int cg1 = 0; cg1 < nGroups_; cg1++) {
315 gezelter 1579
316 gezelter 1576 groupCutoff[cg1] = 0.0;
317     vector<int> atomList = getAtomsInGroupRow(cg1);
318 gezelter 1579
319 gezelter 1576 for (vector<int>::iterator ia = atomList.begin();
320     ia != atomList.end(); ++ia) {
321     int atom1 = (*ia);
322 gezelter 1583 atid = idents[atom1];
323 gezelter 1576 if (atypeCutoff[atid] > groupCutoff[cg1]) {
324     groupCutoff[cg1] = atypeCutoff[atid];
325     }
326     }
327    
328     bool gTypeFound = false;
329     for (int gt = 0; gt < gTypeCutoffs.size(); gt++) {
330     if (abs(groupCutoff[cg1] - gTypeCutoffs[gt]) < tol) {
331     groupToGtype[cg1] = gt;
332     gTypeFound = true;
333     }
334     }
335     if (!gTypeFound) {
336     gTypeCutoffs.push_back( groupCutoff[cg1] );
337     groupToGtype[cg1] = gTypeCutoffs.size() - 1;
338     }
339     }
340     #endif
341    
342     // Now we find the maximum group cutoff value present in the simulation
343    
344 gezelter 1590 RealType groupMax = *max_element(gTypeCutoffs.begin(),
345     gTypeCutoffs.end());
346 gezelter 1576
347     #ifdef IS_MPI
348 gezelter 1590 MPI::COMM_WORLD.Allreduce(&groupMax, &groupMax, 1, MPI::REALTYPE,
349     MPI::MAX);
350 gezelter 1576 #endif
351    
352     RealType tradRcut = groupMax;
353    
354     for (int i = 0; i < gTypeCutoffs.size(); i++) {
355 gezelter 1579 for (int j = 0; j < gTypeCutoffs.size(); j++) {
356 gezelter 1576 RealType thisRcut;
357     switch(cutoffPolicy_) {
358     case TRADITIONAL:
359     thisRcut = tradRcut;
360 gezelter 1579 break;
361 gezelter 1576 case MIX:
362     thisRcut = 0.5 * (gTypeCutoffs[i] + gTypeCutoffs[j]);
363 gezelter 1579 break;
364 gezelter 1576 case MAX:
365     thisRcut = max(gTypeCutoffs[i], gTypeCutoffs[j]);
366 gezelter 1579 break;
367 gezelter 1576 default:
368     sprintf(painCave.errMsg,
369     "ForceMatrixDecomposition::createGtypeCutoffMap "
370     "hit an unknown cutoff policy!\n");
371     painCave.severity = OPENMD_ERROR;
372     painCave.isFatal = 1;
373 gezelter 1579 simError();
374     break;
375 gezelter 1576 }
376    
377     pair<int,int> key = make_pair(i,j);
378     gTypeCutoffMap[key].first = thisRcut;
379    
380     if (thisRcut > largestRcut_) largestRcut_ = thisRcut;
381    
382     gTypeCutoffMap[key].second = thisRcut*thisRcut;
383    
384     gTypeCutoffMap[key].third = pow(thisRcut + skinThickness_, 2);
385    
386     // sanity check
387    
388     if (userChoseCutoff_) {
389     if (abs(gTypeCutoffMap[key].first - userCutoff_) > 0.0001) {
390     sprintf(painCave.errMsg,
391     "ForceMatrixDecomposition::createGtypeCutoffMap "
392 gezelter 1583 "user-specified rCut (%lf) does not match computed group Cutoff\n", userCutoff_);
393 gezelter 1576 painCave.severity = OPENMD_ERROR;
394     painCave.isFatal = 1;
395     simError();
396     }
397     }
398     }
399     }
400 gezelter 1539 }
401 gezelter 1576
402    
403     groupCutoffs ForceMatrixDecomposition::getGroupCutoffs(int cg1, int cg2) {
404 gezelter 1579 int i, j;
405 gezelter 1576 #ifdef IS_MPI
406     i = groupRowToGtype[cg1];
407     j = groupColToGtype[cg2];
408     #else
409     i = groupToGtype[cg1];
410     j = groupToGtype[cg2];
411 gezelter 1579 #endif
412 gezelter 1576 return gTypeCutoffMap[make_pair(i,j)];
413     }
414    
415 gezelter 1579 int ForceMatrixDecomposition::getTopologicalDistance(int atom1, int atom2) {
416     for (int j = 0; j < toposForAtom[atom1].size(); j++) {
417     if (toposForAtom[atom1][j] == atom2)
418     return topoDist[atom1][j];
419     }
420     return 0;
421     }
422 gezelter 1576
423 gezelter 1575 void ForceMatrixDecomposition::zeroWorkArrays() {
424 gezelter 1583 pairwisePot = 0.0;
425     embeddingPot = 0.0;
426 gezelter 1575
427     #ifdef IS_MPI
428     if (storageLayout_ & DataStorage::dslForce) {
429     fill(atomRowData.force.begin(), atomRowData.force.end(), V3Zero);
430     fill(atomColData.force.begin(), atomColData.force.end(), V3Zero);
431     }
432    
433     if (storageLayout_ & DataStorage::dslTorque) {
434     fill(atomRowData.torque.begin(), atomRowData.torque.end(), V3Zero);
435     fill(atomColData.torque.begin(), atomColData.torque.end(), V3Zero);
436     }
437    
438     fill(pot_row.begin(), pot_row.end(),
439     Vector<RealType, N_INTERACTION_FAMILIES> (0.0));
440    
441     fill(pot_col.begin(), pot_col.end(),
442 gezelter 1583 Vector<RealType, N_INTERACTION_FAMILIES> (0.0));
443 gezelter 1575
444     if (storageLayout_ & DataStorage::dslParticlePot) {
445 gezelter 1590 fill(atomRowData.particlePot.begin(), atomRowData.particlePot.end(),
446     0.0);
447     fill(atomColData.particlePot.begin(), atomColData.particlePot.end(),
448     0.0);
449 gezelter 1575 }
450    
451     if (storageLayout_ & DataStorage::dslDensity) {
452     fill(atomRowData.density.begin(), atomRowData.density.end(), 0.0);
453     fill(atomColData.density.begin(), atomColData.density.end(), 0.0);
454     }
455    
456     if (storageLayout_ & DataStorage::dslFunctional) {
457 gezelter 1590 fill(atomRowData.functional.begin(), atomRowData.functional.end(),
458     0.0);
459     fill(atomColData.functional.begin(), atomColData.functional.end(),
460     0.0);
461 gezelter 1575 }
462    
463     if (storageLayout_ & DataStorage::dslFunctionalDerivative) {
464     fill(atomRowData.functionalDerivative.begin(),
465     atomRowData.functionalDerivative.end(), 0.0);
466     fill(atomColData.functionalDerivative.begin(),
467     atomColData.functionalDerivative.end(), 0.0);
468     }
469    
470 gezelter 1586 if (storageLayout_ & DataStorage::dslSkippedCharge) {
471 gezelter 1587 fill(atomRowData.skippedCharge.begin(),
472     atomRowData.skippedCharge.end(), 0.0);
473     fill(atomColData.skippedCharge.begin(),
474     atomColData.skippedCharge.end(), 0.0);
475 gezelter 1586 }
476    
477 gezelter 1590 #endif
478     // even in parallel, we need to zero out the local arrays:
479    
480 gezelter 1575 if (storageLayout_ & DataStorage::dslParticlePot) {
481     fill(snap_->atomData.particlePot.begin(),
482     snap_->atomData.particlePot.end(), 0.0);
483     }
484    
485     if (storageLayout_ & DataStorage::dslDensity) {
486     fill(snap_->atomData.density.begin(),
487     snap_->atomData.density.end(), 0.0);
488     }
489     if (storageLayout_ & DataStorage::dslFunctional) {
490     fill(snap_->atomData.functional.begin(),
491     snap_->atomData.functional.end(), 0.0);
492     }
493     if (storageLayout_ & DataStorage::dslFunctionalDerivative) {
494     fill(snap_->atomData.functionalDerivative.begin(),
495     snap_->atomData.functionalDerivative.end(), 0.0);
496     }
497 gezelter 1586 if (storageLayout_ & DataStorage::dslSkippedCharge) {
498     fill(snap_->atomData.skippedCharge.begin(),
499     snap_->atomData.skippedCharge.end(), 0.0);
500     }
501 gezelter 1575
502     }
503    
504    
505 gezelter 1549 void ForceMatrixDecomposition::distributeData() {
506 gezelter 1551 snap_ = sman_->getCurrentSnapshot();
507     storageLayout_ = sman_->getStorageLayout();
508 chuckv 1538 #ifdef IS_MPI
509 gezelter 1540
510 gezelter 1539 // gather up the atomic positions
511 gezelter 1551 AtomCommVectorRow->gather(snap_->atomData.position,
512     atomRowData.position);
513     AtomCommVectorColumn->gather(snap_->atomData.position,
514     atomColData.position);
515 gezelter 1539
516     // gather up the cutoff group positions
517 gezelter 1551 cgCommVectorRow->gather(snap_->cgData.position,
518     cgRowData.position);
519     cgCommVectorColumn->gather(snap_->cgData.position,
520     cgColData.position);
521 gezelter 1539
522     // if needed, gather the atomic rotation matrices
523 gezelter 1551 if (storageLayout_ & DataStorage::dslAmat) {
524     AtomCommMatrixRow->gather(snap_->atomData.aMat,
525     atomRowData.aMat);
526     AtomCommMatrixColumn->gather(snap_->atomData.aMat,
527     atomColData.aMat);
528 gezelter 1539 }
529    
530     // if needed, gather the atomic eletrostatic frames
531 gezelter 1551 if (storageLayout_ & DataStorage::dslElectroFrame) {
532     AtomCommMatrixRow->gather(snap_->atomData.electroFrame,
533     atomRowData.electroFrame);
534     AtomCommMatrixColumn->gather(snap_->atomData.electroFrame,
535     atomColData.electroFrame);
536 gezelter 1539 }
537 gezelter 1590
538 gezelter 1539 #endif
539     }
540    
541 gezelter 1575 /* collects information obtained during the pre-pair loop onto local
542     * data structures.
543     */
544 gezelter 1549 void ForceMatrixDecomposition::collectIntermediateData() {
545 gezelter 1551 snap_ = sman_->getCurrentSnapshot();
546     storageLayout_ = sman_->getStorageLayout();
547 gezelter 1539 #ifdef IS_MPI
548    
549 gezelter 1551 if (storageLayout_ & DataStorage::dslDensity) {
550    
551     AtomCommRealRow->scatter(atomRowData.density,
552     snap_->atomData.density);
553    
554     int n = snap_->atomData.density.size();
555 gezelter 1575 vector<RealType> rho_tmp(n, 0.0);
556 gezelter 1551 AtomCommRealColumn->scatter(atomColData.density, rho_tmp);
557 gezelter 1539 for (int i = 0; i < n; i++)
558 gezelter 1551 snap_->atomData.density[i] += rho_tmp[i];
559 gezelter 1539 }
560 chuckv 1538 #endif
561 gezelter 1539 }
562 gezelter 1575
563     /*
564     * redistributes information obtained during the pre-pair loop out to
565     * row and column-indexed data structures
566     */
567 gezelter 1549 void ForceMatrixDecomposition::distributeIntermediateData() {
568 gezelter 1551 snap_ = sman_->getCurrentSnapshot();
569     storageLayout_ = sman_->getStorageLayout();
570 chuckv 1538 #ifdef IS_MPI
571 gezelter 1551 if (storageLayout_ & DataStorage::dslFunctional) {
572     AtomCommRealRow->gather(snap_->atomData.functional,
573     atomRowData.functional);
574     AtomCommRealColumn->gather(snap_->atomData.functional,
575     atomColData.functional);
576 gezelter 1539 }
577    
578 gezelter 1551 if (storageLayout_ & DataStorage::dslFunctionalDerivative) {
579     AtomCommRealRow->gather(snap_->atomData.functionalDerivative,
580     atomRowData.functionalDerivative);
581     AtomCommRealColumn->gather(snap_->atomData.functionalDerivative,
582     atomColData.functionalDerivative);
583 gezelter 1539 }
584 chuckv 1538 #endif
585     }
586 gezelter 1539
587    
588 gezelter 1549 void ForceMatrixDecomposition::collectData() {
589 gezelter 1551 snap_ = sman_->getCurrentSnapshot();
590     storageLayout_ = sman_->getStorageLayout();
591     #ifdef IS_MPI
592     int n = snap_->atomData.force.size();
593 gezelter 1544 vector<Vector3d> frc_tmp(n, V3Zero);
594 gezelter 1541
595 gezelter 1551 AtomCommVectorRow->scatter(atomRowData.force, frc_tmp);
596 gezelter 1541 for (int i = 0; i < n; i++) {
597 gezelter 1551 snap_->atomData.force[i] += frc_tmp[i];
598 gezelter 1541 frc_tmp[i] = 0.0;
599     }
600 gezelter 1540
601 gezelter 1551 AtomCommVectorColumn->scatter(atomColData.force, frc_tmp);
602 gezelter 1540 for (int i = 0; i < n; i++)
603 gezelter 1551 snap_->atomData.force[i] += frc_tmp[i];
604 gezelter 1590
605 gezelter 1551 if (storageLayout_ & DataStorage::dslTorque) {
606 gezelter 1541
607 gezelter 1587 int nt = snap_->atomData.torque.size();
608 gezelter 1544 vector<Vector3d> trq_tmp(nt, V3Zero);
609 gezelter 1541
610 gezelter 1551 AtomCommVectorRow->scatter(atomRowData.torque, trq_tmp);
611 gezelter 1587 for (int i = 0; i < nt; i++) {
612 gezelter 1551 snap_->atomData.torque[i] += trq_tmp[i];
613 gezelter 1541 trq_tmp[i] = 0.0;
614     }
615 gezelter 1540
616 gezelter 1551 AtomCommVectorColumn->scatter(atomColData.torque, trq_tmp);
617 gezelter 1587 for (int i = 0; i < nt; i++)
618 gezelter 1551 snap_->atomData.torque[i] += trq_tmp[i];
619 gezelter 1540 }
620 gezelter 1587
621     if (storageLayout_ & DataStorage::dslSkippedCharge) {
622    
623     int ns = snap_->atomData.skippedCharge.size();
624     vector<RealType> skch_tmp(ns, 0.0);
625    
626     AtomCommRealRow->scatter(atomRowData.skippedCharge, skch_tmp);
627     for (int i = 0; i < ns; i++) {
628 gezelter 1590 snap_->atomData.skippedCharge[i] += skch_tmp[i];
629 gezelter 1587 skch_tmp[i] = 0.0;
630     }
631    
632     AtomCommRealColumn->scatter(atomColData.skippedCharge, skch_tmp);
633     for (int i = 0; i < ns; i++)
634     snap_->atomData.skippedCharge[i] += skch_tmp[i];
635     }
636 gezelter 1540
637 gezelter 1567 nLocal_ = snap_->getNumberOfAtoms();
638 gezelter 1544
639 gezelter 1575 vector<potVec> pot_temp(nLocal_,
640     Vector<RealType, N_INTERACTION_FAMILIES> (0.0));
641    
642     // scatter/gather pot_row into the members of my column
643    
644     AtomCommPotRow->scatter(pot_row, pot_temp);
645    
646     for (int ii = 0; ii < pot_temp.size(); ii++ )
647 gezelter 1583 pairwisePot += pot_temp[ii];
648 gezelter 1540
649 gezelter 1575 fill(pot_temp.begin(), pot_temp.end(),
650     Vector<RealType, N_INTERACTION_FAMILIES> (0.0));
651    
652     AtomCommPotColumn->scatter(pot_col, pot_temp);
653    
654     for (int ii = 0; ii < pot_temp.size(); ii++ )
655 gezelter 1583 pairwisePot += pot_temp[ii];
656 gezelter 1539 #endif
657 gezelter 1583
658 chuckv 1538 }
659 gezelter 1551
660 gezelter 1570 int ForceMatrixDecomposition::getNAtomsInRow() {
661     #ifdef IS_MPI
662     return nAtomsInRow_;
663     #else
664     return nLocal_;
665     #endif
666     }
667    
668 gezelter 1569 /**
669     * returns the list of atoms belonging to this group.
670     */
671     vector<int> ForceMatrixDecomposition::getAtomsInGroupRow(int cg1){
672     #ifdef IS_MPI
673     return groupListRow_[cg1];
674     #else
675     return groupList_[cg1];
676     #endif
677     }
678    
679     vector<int> ForceMatrixDecomposition::getAtomsInGroupColumn(int cg2){
680     #ifdef IS_MPI
681     return groupListCol_[cg2];
682     #else
683     return groupList_[cg2];
684     #endif
685     }
686 chuckv 1538
687 gezelter 1551 Vector3d ForceMatrixDecomposition::getIntergroupVector(int cg1, int cg2){
688     Vector3d d;
689    
690     #ifdef IS_MPI
691     d = cgColData.position[cg2] - cgRowData.position[cg1];
692     #else
693     d = snap_->cgData.position[cg2] - snap_->cgData.position[cg1];
694     #endif
695    
696     snap_->wrapVector(d);
697     return d;
698     }
699    
700    
701     Vector3d ForceMatrixDecomposition::getAtomToGroupVectorRow(int atom1, int cg1){
702    
703     Vector3d d;
704    
705     #ifdef IS_MPI
706     d = cgRowData.position[cg1] - atomRowData.position[atom1];
707     #else
708     d = snap_->cgData.position[cg1] - snap_->atomData.position[atom1];
709     #endif
710    
711     snap_->wrapVector(d);
712     return d;
713     }
714    
715     Vector3d ForceMatrixDecomposition::getAtomToGroupVectorColumn(int atom2, int cg2){
716     Vector3d d;
717    
718     #ifdef IS_MPI
719     d = cgColData.position[cg2] - atomColData.position[atom2];
720     #else
721     d = snap_->cgData.position[cg2] - snap_->atomData.position[atom2];
722     #endif
723    
724     snap_->wrapVector(d);
725     return d;
726     }
727 gezelter 1569
728     RealType ForceMatrixDecomposition::getMassFactorRow(int atom1) {
729     #ifdef IS_MPI
730     return massFactorsRow[atom1];
731     #else
732 gezelter 1581 return massFactors[atom1];
733 gezelter 1569 #endif
734     }
735    
736     RealType ForceMatrixDecomposition::getMassFactorColumn(int atom2) {
737     #ifdef IS_MPI
738     return massFactorsCol[atom2];
739     #else
740 gezelter 1581 return massFactors[atom2];
741 gezelter 1569 #endif
742    
743     }
744 gezelter 1551
745     Vector3d ForceMatrixDecomposition::getInteratomicVector(int atom1, int atom2){
746     Vector3d d;
747    
748     #ifdef IS_MPI
749     d = atomColData.position[atom2] - atomRowData.position[atom1];
750     #else
751     d = snap_->atomData.position[atom2] - snap_->atomData.position[atom1];
752     #endif
753    
754     snap_->wrapVector(d);
755     return d;
756     }
757    
758 gezelter 1587 vector<int> ForceMatrixDecomposition::getExcludesForAtom(int atom1) {
759     return excludesForAtom[atom1];
760 gezelter 1570 }
761    
762     /**
763 gezelter 1587 * We need to exclude some overcounted interactions that result from
764 gezelter 1575 * the parallel decomposition.
765 gezelter 1570 */
766     bool ForceMatrixDecomposition::skipAtomPair(int atom1, int atom2) {
767     int unique_id_1, unique_id_2;
768    
769     #ifdef IS_MPI
770     // in MPI, we have to look up the unique IDs for each atom
771     unique_id_1 = AtomRowToGlobal[atom1];
772     unique_id_2 = AtomColToGlobal[atom2];
773    
774     // this situation should only arise in MPI simulations
775     if (unique_id_1 == unique_id_2) return true;
776    
777     // this prevents us from doing the pair on multiple processors
778     if (unique_id_1 < unique_id_2) {
779     if ((unique_id_1 + unique_id_2) % 2 == 0) return true;
780     } else {
781     if ((unique_id_1 + unique_id_2) % 2 == 1) return true;
782     }
783 gezelter 1587 #endif
784     return false;
785     }
786    
787     /**
788     * We need to handle the interactions for atoms who are involved in
789     * the same rigid body as well as some short range interactions
790     * (bonds, bends, torsions) differently from other interactions.
791     * We'll still visit the pairwise routines, but with a flag that
792     * tells those routines to exclude the pair from direct long range
793     * interactions. Some indirect interactions (notably reaction
794     * field) must still be handled for these pairs.
795     */
796     bool ForceMatrixDecomposition::excludeAtomPair(int atom1, int atom2) {
797     int unique_id_2;
798    
799     #ifdef IS_MPI
800     // in MPI, we have to look up the unique IDs for the row atom.
801     unique_id_2 = AtomColToGlobal[atom2];
802 gezelter 1570 #else
803     // in the normal loop, the atom numbers are unique
804     unique_id_2 = atom2;
805     #endif
806    
807 gezelter 1587 for (vector<int>::iterator i = excludesForAtom[atom1].begin();
808     i != excludesForAtom[atom1].end(); ++i) {
809 gezelter 1570 if ( (*i) == unique_id_2 ) return true;
810 gezelter 1583 }
811 gezelter 1579
812 gezelter 1583 return false;
813 gezelter 1570 }
814    
815    
816 gezelter 1551 void ForceMatrixDecomposition::addForceToAtomRow(int atom1, Vector3d fg){
817     #ifdef IS_MPI
818     atomRowData.force[atom1] += fg;
819     #else
820     snap_->atomData.force[atom1] += fg;
821     #endif
822     }
823    
824     void ForceMatrixDecomposition::addForceToAtomColumn(int atom2, Vector3d fg){
825     #ifdef IS_MPI
826     atomColData.force[atom2] += fg;
827     #else
828     snap_->atomData.force[atom2] += fg;
829     #endif
830     }
831    
832     // filling interaction blocks with pointers
833 gezelter 1582 void ForceMatrixDecomposition::fillInteractionData(InteractionData &idat,
834 gezelter 1587 int atom1, int atom2) {
835    
836     idat.excluded = excludeAtomPair(atom1, atom2);
837    
838 gezelter 1551 #ifdef IS_MPI
839 gezelter 1571
840     idat.atypes = make_pair( ff_->getAtomType(identsRow[atom1]),
841     ff_->getAtomType(identsCol[atom2]) );
842 gezelter 1575
843 gezelter 1551 if (storageLayout_ & DataStorage::dslAmat) {
844 gezelter 1554 idat.A1 = &(atomRowData.aMat[atom1]);
845     idat.A2 = &(atomColData.aMat[atom2]);
846 gezelter 1551 }
847 gezelter 1567
848 gezelter 1551 if (storageLayout_ & DataStorage::dslElectroFrame) {
849 gezelter 1554 idat.eFrame1 = &(atomRowData.electroFrame[atom1]);
850     idat.eFrame2 = &(atomColData.electroFrame[atom2]);
851 gezelter 1551 }
852    
853     if (storageLayout_ & DataStorage::dslTorque) {
854 gezelter 1554 idat.t1 = &(atomRowData.torque[atom1]);
855     idat.t2 = &(atomColData.torque[atom2]);
856 gezelter 1551 }
857    
858     if (storageLayout_ & DataStorage::dslDensity) {
859 gezelter 1554 idat.rho1 = &(atomRowData.density[atom1]);
860     idat.rho2 = &(atomColData.density[atom2]);
861 gezelter 1551 }
862    
863 gezelter 1575 if (storageLayout_ & DataStorage::dslFunctional) {
864     idat.frho1 = &(atomRowData.functional[atom1]);
865     idat.frho2 = &(atomColData.functional[atom2]);
866     }
867    
868 gezelter 1551 if (storageLayout_ & DataStorage::dslFunctionalDerivative) {
869 gezelter 1554 idat.dfrho1 = &(atomRowData.functionalDerivative[atom1]);
870     idat.dfrho2 = &(atomColData.functionalDerivative[atom2]);
871 gezelter 1551 }
872 gezelter 1570
873 gezelter 1575 if (storageLayout_ & DataStorage::dslParticlePot) {
874     idat.particlePot1 = &(atomRowData.particlePot[atom1]);
875     idat.particlePot2 = &(atomColData.particlePot[atom2]);
876     }
877    
878 gezelter 1587 if (storageLayout_ & DataStorage::dslSkippedCharge) {
879     idat.skippedCharge1 = &(atomRowData.skippedCharge[atom1]);
880     idat.skippedCharge2 = &(atomColData.skippedCharge[atom2]);
881     }
882    
883 gezelter 1562 #else
884 gezelter 1571
885 gezelter 1583 idat.atypes = make_pair( ff_->getAtomType(idents[atom1]),
886     ff_->getAtomType(idents[atom2]) );
887 gezelter 1571
888 gezelter 1562 if (storageLayout_ & DataStorage::dslAmat) {
889     idat.A1 = &(snap_->atomData.aMat[atom1]);
890     idat.A2 = &(snap_->atomData.aMat[atom2]);
891     }
892    
893     if (storageLayout_ & DataStorage::dslElectroFrame) {
894     idat.eFrame1 = &(snap_->atomData.electroFrame[atom1]);
895     idat.eFrame2 = &(snap_->atomData.electroFrame[atom2]);
896     }
897    
898     if (storageLayout_ & DataStorage::dslTorque) {
899     idat.t1 = &(snap_->atomData.torque[atom1]);
900     idat.t2 = &(snap_->atomData.torque[atom2]);
901     }
902    
903 gezelter 1583 if (storageLayout_ & DataStorage::dslDensity) {
904 gezelter 1562 idat.rho1 = &(snap_->atomData.density[atom1]);
905     idat.rho2 = &(snap_->atomData.density[atom2]);
906     }
907    
908 gezelter 1575 if (storageLayout_ & DataStorage::dslFunctional) {
909     idat.frho1 = &(snap_->atomData.functional[atom1]);
910     idat.frho2 = &(snap_->atomData.functional[atom2]);
911     }
912    
913 gezelter 1562 if (storageLayout_ & DataStorage::dslFunctionalDerivative) {
914     idat.dfrho1 = &(snap_->atomData.functionalDerivative[atom1]);
915     idat.dfrho2 = &(snap_->atomData.functionalDerivative[atom2]);
916     }
917 gezelter 1575
918     if (storageLayout_ & DataStorage::dslParticlePot) {
919     idat.particlePot1 = &(snap_->atomData.particlePot[atom1]);
920     idat.particlePot2 = &(snap_->atomData.particlePot[atom2]);
921     }
922    
923 gezelter 1587 if (storageLayout_ & DataStorage::dslSkippedCharge) {
924     idat.skippedCharge1 = &(snap_->atomData.skippedCharge[atom1]);
925     idat.skippedCharge2 = &(snap_->atomData.skippedCharge[atom2]);
926     }
927 gezelter 1551 #endif
928     }
929 gezelter 1567
930 gezelter 1575
931 gezelter 1582 void ForceMatrixDecomposition::unpackInteractionData(InteractionData &idat, int atom1, int atom2) {
932 gezelter 1575 #ifdef IS_MPI
933     pot_row[atom1] += 0.5 * *(idat.pot);
934     pot_col[atom2] += 0.5 * *(idat.pot);
935    
936     atomRowData.force[atom1] += *(idat.f1);
937     atomColData.force[atom2] -= *(idat.f1);
938     #else
939 gezelter 1583 pairwisePot += *(idat.pot);
940    
941 gezelter 1575 snap_->atomData.force[atom1] += *(idat.f1);
942     snap_->atomData.force[atom2] -= *(idat.f1);
943     #endif
944 gezelter 1586
945 gezelter 1575 }
946    
947 gezelter 1562 /*
948     * buildNeighborList
949     *
950     * first element of pair is row-indexed CutoffGroup
951     * second element of pair is column-indexed CutoffGroup
952     */
953 gezelter 1567 vector<pair<int, int> > ForceMatrixDecomposition::buildNeighborList() {
954    
955     vector<pair<int, int> > neighborList;
956 gezelter 1576 groupCutoffs cuts;
957 gezelter 1587 bool doAllPairs = false;
958    
959 gezelter 1567 #ifdef IS_MPI
960 gezelter 1568 cellListRow_.clear();
961     cellListCol_.clear();
962 gezelter 1567 #else
963 gezelter 1568 cellList_.clear();
964 gezelter 1567 #endif
965 gezelter 1562
966 gezelter 1576 RealType rList_ = (largestRcut_ + skinThickness_);
967 gezelter 1567 RealType rl2 = rList_ * rList_;
968     Snapshot* snap_ = sman_->getCurrentSnapshot();
969 gezelter 1562 Mat3x3d Hmat = snap_->getHmat();
970     Vector3d Hx = Hmat.getColumn(0);
971     Vector3d Hy = Hmat.getColumn(1);
972     Vector3d Hz = Hmat.getColumn(2);
973    
974 gezelter 1568 nCells_.x() = (int) ( Hx.length() )/ rList_;
975     nCells_.y() = (int) ( Hy.length() )/ rList_;
976     nCells_.z() = (int) ( Hz.length() )/ rList_;
977 gezelter 1562
978 gezelter 1587 // handle small boxes where the cell offsets can end up repeating cells
979    
980     if (nCells_.x() < 3) doAllPairs = true;
981     if (nCells_.y() < 3) doAllPairs = true;
982     if (nCells_.z() < 3) doAllPairs = true;
983    
984 gezelter 1567 Mat3x3d invHmat = snap_->getInvHmat();
985     Vector3d rs, scaled, dr;
986     Vector3i whichCell;
987     int cellIndex;
988 gezelter 1579 int nCtot = nCells_.x() * nCells_.y() * nCells_.z();
989 gezelter 1567
990     #ifdef IS_MPI
991 gezelter 1579 cellListRow_.resize(nCtot);
992     cellListCol_.resize(nCtot);
993     #else
994     cellList_.resize(nCtot);
995     #endif
996 gezelter 1582
997 gezelter 1587 if (!doAllPairs) {
998 gezelter 1579 #ifdef IS_MPI
999 gezelter 1581
1000 gezelter 1587 for (int i = 0; i < nGroupsInRow_; i++) {
1001     rs = cgRowData.position[i];
1002    
1003     // scaled positions relative to the box vectors
1004     scaled = invHmat * rs;
1005    
1006     // wrap the vector back into the unit box by subtracting integer box
1007     // numbers
1008     for (int j = 0; j < 3; j++) {
1009     scaled[j] -= roundMe(scaled[j]);
1010     scaled[j] += 0.5;
1011     }
1012    
1013     // find xyz-indices of cell that cutoffGroup is in.
1014     whichCell.x() = nCells_.x() * scaled.x();
1015     whichCell.y() = nCells_.y() * scaled.y();
1016     whichCell.z() = nCells_.z() * scaled.z();
1017    
1018     // find single index of this cell:
1019     cellIndex = Vlinear(whichCell, nCells_);
1020    
1021     // add this cutoff group to the list of groups in this cell;
1022     cellListRow_[cellIndex].push_back(i);
1023 gezelter 1581 }
1024 gezelter 1587
1025     for (int i = 0; i < nGroupsInCol_; i++) {
1026     rs = cgColData.position[i];
1027    
1028     // scaled positions relative to the box vectors
1029     scaled = invHmat * rs;
1030    
1031     // wrap the vector back into the unit box by subtracting integer box
1032     // numbers
1033     for (int j = 0; j < 3; j++) {
1034     scaled[j] -= roundMe(scaled[j]);
1035     scaled[j] += 0.5;
1036     }
1037    
1038     // find xyz-indices of cell that cutoffGroup is in.
1039     whichCell.x() = nCells_.x() * scaled.x();
1040     whichCell.y() = nCells_.y() * scaled.y();
1041     whichCell.z() = nCells_.z() * scaled.z();
1042    
1043     // find single index of this cell:
1044     cellIndex = Vlinear(whichCell, nCells_);
1045    
1046     // add this cutoff group to the list of groups in this cell;
1047     cellListCol_[cellIndex].push_back(i);
1048 gezelter 1581 }
1049 gezelter 1567 #else
1050 gezelter 1587 for (int i = 0; i < nGroups_; i++) {
1051     rs = snap_->cgData.position[i];
1052    
1053     // scaled positions relative to the box vectors
1054     scaled = invHmat * rs;
1055    
1056     // wrap the vector back into the unit box by subtracting integer box
1057     // numbers
1058     for (int j = 0; j < 3; j++) {
1059     scaled[j] -= roundMe(scaled[j]);
1060     scaled[j] += 0.5;
1061     }
1062    
1063     // find xyz-indices of cell that cutoffGroup is in.
1064     whichCell.x() = nCells_.x() * scaled.x();
1065     whichCell.y() = nCells_.y() * scaled.y();
1066     whichCell.z() = nCells_.z() * scaled.z();
1067    
1068     // find single index of this cell:
1069     cellIndex = Vlinear(whichCell, nCells_);
1070    
1071     // add this cutoff group to the list of groups in this cell;
1072     cellList_[cellIndex].push_back(i);
1073 gezelter 1581 }
1074 gezelter 1567 #endif
1075    
1076 gezelter 1587 for (int m1z = 0; m1z < nCells_.z(); m1z++) {
1077     for (int m1y = 0; m1y < nCells_.y(); m1y++) {
1078     for (int m1x = 0; m1x < nCells_.x(); m1x++) {
1079     Vector3i m1v(m1x, m1y, m1z);
1080     int m1 = Vlinear(m1v, nCells_);
1081 gezelter 1568
1082 gezelter 1587 for (vector<Vector3i>::iterator os = cellOffsets_.begin();
1083     os != cellOffsets_.end(); ++os) {
1084    
1085     Vector3i m2v = m1v + (*os);
1086    
1087     if (m2v.x() >= nCells_.x()) {
1088     m2v.x() = 0;
1089     } else if (m2v.x() < 0) {
1090     m2v.x() = nCells_.x() - 1;
1091     }
1092    
1093     if (m2v.y() >= nCells_.y()) {
1094     m2v.y() = 0;
1095     } else if (m2v.y() < 0) {
1096     m2v.y() = nCells_.y() - 1;
1097     }
1098    
1099     if (m2v.z() >= nCells_.z()) {
1100     m2v.z() = 0;
1101     } else if (m2v.z() < 0) {
1102     m2v.z() = nCells_.z() - 1;
1103     }
1104    
1105     int m2 = Vlinear (m2v, nCells_);
1106    
1107 gezelter 1567 #ifdef IS_MPI
1108 gezelter 1587 for (vector<int>::iterator j1 = cellListRow_[m1].begin();
1109     j1 != cellListRow_[m1].end(); ++j1) {
1110     for (vector<int>::iterator j2 = cellListCol_[m2].begin();
1111     j2 != cellListCol_[m2].end(); ++j2) {
1112    
1113     // Always do this if we're in different cells or if
1114     // we're in the same cell and the global index of the
1115     // j2 cutoff group is less than the j1 cutoff group
1116    
1117     if (m2 != m1 || cgColToGlobal[(*j2)] < cgRowToGlobal[(*j1)]) {
1118     dr = cgColData.position[(*j2)] - cgRowData.position[(*j1)];
1119     snap_->wrapVector(dr);
1120     cuts = getGroupCutoffs( (*j1), (*j2) );
1121     if (dr.lengthSquare() < cuts.third) {
1122     neighborList.push_back(make_pair((*j1), (*j2)));
1123     }
1124 gezelter 1562 }
1125     }
1126     }
1127 gezelter 1567 #else
1128 gezelter 1587
1129     for (vector<int>::iterator j1 = cellList_[m1].begin();
1130     j1 != cellList_[m1].end(); ++j1) {
1131     for (vector<int>::iterator j2 = cellList_[m2].begin();
1132     j2 != cellList_[m2].end(); ++j2) {
1133    
1134     // Always do this if we're in different cells or if
1135     // we're in the same cell and the global index of the
1136     // j2 cutoff group is less than the j1 cutoff group
1137    
1138     if (m2 != m1 || (*j2) < (*j1)) {
1139     dr = snap_->cgData.position[(*j2)] - snap_->cgData.position[(*j1)];
1140     snap_->wrapVector(dr);
1141     cuts = getGroupCutoffs( (*j1), (*j2) );
1142     if (dr.lengthSquare() < cuts.third) {
1143     neighborList.push_back(make_pair((*j1), (*j2)));
1144     }
1145 gezelter 1567 }
1146     }
1147     }
1148 gezelter 1587 #endif
1149 gezelter 1567 }
1150 gezelter 1562 }
1151     }
1152     }
1153 gezelter 1587 } else {
1154     // branch to do all cutoff group pairs
1155     #ifdef IS_MPI
1156     for (int j1 = 0; j1 < nGroupsInRow_; j1++) {
1157     for (int j2 = 0; j2 < nGroupsInCol_; j2++) {
1158     dr = cgColData.position[j2] - cgRowData.position[j1];
1159     snap_->wrapVector(dr);
1160     cuts = getGroupCutoffs( j1, j2 );
1161     if (dr.lengthSquare() < cuts.third) {
1162     neighborList.push_back(make_pair(j1, j2));
1163     }
1164     }
1165     }
1166     #else
1167     for (int j1 = 0; j1 < nGroups_ - 1; j1++) {
1168     for (int j2 = j1 + 1; j2 < nGroups_; j2++) {
1169     dr = snap_->cgData.position[j2] - snap_->cgData.position[j1];
1170     snap_->wrapVector(dr);
1171     cuts = getGroupCutoffs( j1, j2 );
1172     if (dr.lengthSquare() < cuts.third) {
1173     neighborList.push_back(make_pair(j1, j2));
1174     }
1175     }
1176     }
1177     #endif
1178 gezelter 1562 }
1179 gezelter 1587
1180 gezelter 1568 // save the local cutoff group positions for the check that is
1181     // done on each loop:
1182     saved_CG_positions_.clear();
1183     for (int i = 0; i < nGroups_; i++)
1184     saved_CG_positions_.push_back(snap_->cgData.position[i]);
1185 gezelter 1587
1186 gezelter 1567 return neighborList;
1187 gezelter 1562 }
1188 gezelter 1539 } //end namespace OpenMD