ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/parallel/ForceMatrixDecomposition.cpp
Revision: 1616
Committed: Tue Aug 30 15:45:35 2011 UTC (13 years, 8 months ago) by gezelter
Original Path: branches/development/src/parallel/ForceMatrixDecomposition.cpp
File size: 42308 byte(s)
Log Message:
Fixing cutoff groups

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