ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/parallel/ForceMatrixDecomposition.cpp
Revision: 1896
Committed: Tue Jul 2 20:02:31 2013 UTC (11 years, 10 months ago) by gezelter
File size: 54767 byte(s)
Log Message:
Speedup tests

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 gezelter 1879 * [3] Sun, Lin & Gezelter, J. Chem. Phys. 128, 234107 (2008).
39 gezelter 1665 * [4] Kuang & Gezelter, J. Chem. Phys. 133, 164101 (2010).
40     * [5] Vardeman, Stocker & Gezelter, J. Chem. Theory Comput. 7, 834 (2011).
41 chuckv 1538 */
42 gezelter 1549 #include "parallel/ForceMatrixDecomposition.hpp"
43 gezelter 1539 #include "math/SquareMatrix3.hpp"
44 gezelter 1544 #include "nonbonded/NonBondedInteraction.hpp"
45     #include "brains/SnapshotManager.hpp"
46 gezelter 1570 #include "brains/PairList.hpp"
47 chuckv 1538
48 gezelter 1541 using namespace std;
49 gezelter 1539 namespace OpenMD {
50 chuckv 1538
51 gezelter 1593 ForceMatrixDecomposition::ForceMatrixDecomposition(SimInfo* info, InteractionManager* iMan) : ForceDecomposition(info, iMan) {
52    
53     // In a parallel computation, row and colum scans must visit all
54     // surrounding cells (not just the 14 upper triangular blocks that
55     // are used when the processor can see all pairs)
56     #ifdef IS_MPI
57 gezelter 1612 cellOffsets_.clear();
58     cellOffsets_.push_back( Vector3i(-1,-1,-1) );
59     cellOffsets_.push_back( Vector3i( 0,-1,-1) );
60     cellOffsets_.push_back( Vector3i( 1,-1,-1) );
61     cellOffsets_.push_back( Vector3i(-1, 0,-1) );
62     cellOffsets_.push_back( Vector3i( 0, 0,-1) );
63     cellOffsets_.push_back( Vector3i( 1, 0,-1) );
64     cellOffsets_.push_back( Vector3i(-1, 1,-1) );
65     cellOffsets_.push_back( Vector3i( 0, 1,-1) );
66     cellOffsets_.push_back( Vector3i( 1, 1,-1) );
67 gezelter 1593 cellOffsets_.push_back( Vector3i(-1,-1, 0) );
68     cellOffsets_.push_back( Vector3i( 0,-1, 0) );
69     cellOffsets_.push_back( Vector3i( 1,-1, 0) );
70 gezelter 1612 cellOffsets_.push_back( Vector3i(-1, 0, 0) );
71     cellOffsets_.push_back( Vector3i( 0, 0, 0) );
72     cellOffsets_.push_back( Vector3i( 1, 0, 0) );
73     cellOffsets_.push_back( Vector3i(-1, 1, 0) );
74     cellOffsets_.push_back( Vector3i( 0, 1, 0) );
75     cellOffsets_.push_back( Vector3i( 1, 1, 0) );
76     cellOffsets_.push_back( Vector3i(-1,-1, 1) );
77     cellOffsets_.push_back( Vector3i( 0,-1, 1) );
78     cellOffsets_.push_back( Vector3i( 1,-1, 1) );
79 gezelter 1593 cellOffsets_.push_back( Vector3i(-1, 0, 1) );
80 gezelter 1612 cellOffsets_.push_back( Vector3i( 0, 0, 1) );
81     cellOffsets_.push_back( Vector3i( 1, 0, 1) );
82     cellOffsets_.push_back( Vector3i(-1, 1, 1) );
83     cellOffsets_.push_back( Vector3i( 0, 1, 1) );
84     cellOffsets_.push_back( Vector3i( 1, 1, 1) );
85 gezelter 1593 #endif
86     }
87    
88    
89 gezelter 1544 /**
90     * distributeInitialData is essentially a copy of the older fortran
91     * SimulationSetup
92     */
93 gezelter 1549 void ForceMatrixDecomposition::distributeInitialData() {
94 gezelter 1551 snap_ = sman_->getCurrentSnapshot();
95     storageLayout_ = sman_->getStorageLayout();
96 gezelter 1571 ff_ = info_->getForceField();
97 gezelter 1567 nLocal_ = snap_->getNumberOfAtoms();
98 gezelter 1723
99 gezelter 1577 nGroups_ = info_->getNLocalCutoffGroups();
100 gezelter 1569 // gather the information for atomtype IDs (atids):
101 gezelter 1583 idents = info_->getIdentArray();
102 gezelter 1569 AtomLocalToGlobal = info_->getGlobalAtomIndices();
103     cgLocalToGlobal = info_->getGlobalGroupIndices();
104     vector<int> globalGroupMembership = info_->getGlobalGroupMembership();
105 gezelter 1586
106 gezelter 1581 massFactors = info_->getMassFactors();
107 gezelter 1584
108 gezelter 1587 PairList* excludes = info_->getExcludedInteractions();
109     PairList* oneTwo = info_->getOneTwoInteractions();
110     PairList* oneThree = info_->getOneThreeInteractions();
111     PairList* oneFour = info_->getOneFourInteractions();
112 gezelter 1723
113     if (needVelocities_)
114     snap_->cgData.setStorageLayout(DataStorage::dslPosition |
115     DataStorage::dslVelocity);
116     else
117     snap_->cgData.setStorageLayout(DataStorage::dslPosition);
118    
119 gezelter 1567 #ifdef IS_MPI
120    
121 gezelter 1593 MPI::Intracomm row = rowComm.getComm();
122     MPI::Intracomm col = colComm.getComm();
123 chuckv 1538
124 gezelter 1593 AtomPlanIntRow = new Plan<int>(row, nLocal_);
125     AtomPlanRealRow = new Plan<RealType>(row, nLocal_);
126     AtomPlanVectorRow = new Plan<Vector3d>(row, nLocal_);
127     AtomPlanMatrixRow = new Plan<Mat3x3d>(row, nLocal_);
128     AtomPlanPotRow = new Plan<potVec>(row, nLocal_);
129 gezelter 1541
130 gezelter 1593 AtomPlanIntColumn = new Plan<int>(col, nLocal_);
131     AtomPlanRealColumn = new Plan<RealType>(col, nLocal_);
132     AtomPlanVectorColumn = new Plan<Vector3d>(col, nLocal_);
133     AtomPlanMatrixColumn = new Plan<Mat3x3d>(col, nLocal_);
134     AtomPlanPotColumn = new Plan<potVec>(col, nLocal_);
135 gezelter 1551
136 gezelter 1593 cgPlanIntRow = new Plan<int>(row, nGroups_);
137     cgPlanVectorRow = new Plan<Vector3d>(row, nGroups_);
138     cgPlanIntColumn = new Plan<int>(col, nGroups_);
139     cgPlanVectorColumn = new Plan<Vector3d>(col, nGroups_);
140 gezelter 1567
141 gezelter 1593 nAtomsInRow_ = AtomPlanIntRow->getSize();
142     nAtomsInCol_ = AtomPlanIntColumn->getSize();
143     nGroupsInRow_ = cgPlanIntRow->getSize();
144     nGroupsInCol_ = cgPlanIntColumn->getSize();
145    
146 gezelter 1551 // Modify the data storage objects with the correct layouts and sizes:
147 gezelter 1567 atomRowData.resize(nAtomsInRow_);
148 gezelter 1551 atomRowData.setStorageLayout(storageLayout_);
149 gezelter 1567 atomColData.resize(nAtomsInCol_);
150 gezelter 1551 atomColData.setStorageLayout(storageLayout_);
151 gezelter 1567 cgRowData.resize(nGroupsInRow_);
152 gezelter 1551 cgRowData.setStorageLayout(DataStorage::dslPosition);
153 gezelter 1567 cgColData.resize(nGroupsInCol_);
154 gezelter 1723 if (needVelocities_)
155     // we only need column velocities if we need them.
156     cgColData.setStorageLayout(DataStorage::dslPosition |
157     DataStorage::dslVelocity);
158     else
159     cgColData.setStorageLayout(DataStorage::dslPosition);
160    
161 gezelter 1577 identsRow.resize(nAtomsInRow_);
162     identsCol.resize(nAtomsInCol_);
163 gezelter 1549
164 gezelter 1593 AtomPlanIntRow->gather(idents, identsRow);
165     AtomPlanIntColumn->gather(idents, identsCol);
166 gezelter 1549
167 gezelter 1589 // allocate memory for the parallel objects
168 gezelter 1591 atypesRow.resize(nAtomsInRow_);
169     atypesCol.resize(nAtomsInCol_);
170    
171     for (int i = 0; i < nAtomsInRow_; i++)
172     atypesRow[i] = ff_->getAtomType(identsRow[i]);
173     for (int i = 0; i < nAtomsInCol_; i++)
174     atypesCol[i] = ff_->getAtomType(identsCol[i]);
175    
176 gezelter 1589 pot_row.resize(nAtomsInRow_);
177     pot_col.resize(nAtomsInCol_);
178    
179 gezelter 1760 expot_row.resize(nAtomsInRow_);
180     expot_col.resize(nAtomsInCol_);
181    
182 gezelter 1591 AtomRowToGlobal.resize(nAtomsInRow_);
183     AtomColToGlobal.resize(nAtomsInCol_);
184 gezelter 1593 AtomPlanIntRow->gather(AtomLocalToGlobal, AtomRowToGlobal);
185     AtomPlanIntColumn->gather(AtomLocalToGlobal, AtomColToGlobal);
186    
187 gezelter 1591 cgRowToGlobal.resize(nGroupsInRow_);
188     cgColToGlobal.resize(nGroupsInCol_);
189 gezelter 1593 cgPlanIntRow->gather(cgLocalToGlobal, cgRowToGlobal);
190     cgPlanIntColumn->gather(cgLocalToGlobal, cgColToGlobal);
191 gezelter 1541
192 gezelter 1591 massFactorsRow.resize(nAtomsInRow_);
193     massFactorsCol.resize(nAtomsInCol_);
194 gezelter 1593 AtomPlanRealRow->gather(massFactors, massFactorsRow);
195     AtomPlanRealColumn->gather(massFactors, massFactorsCol);
196 gezelter 1569
197     groupListRow_.clear();
198 gezelter 1577 groupListRow_.resize(nGroupsInRow_);
199 gezelter 1569 for (int i = 0; i < nGroupsInRow_; i++) {
200     int gid = cgRowToGlobal[i];
201     for (int j = 0; j < nAtomsInRow_; j++) {
202     int aid = AtomRowToGlobal[j];
203     if (globalGroupMembership[aid] == gid)
204     groupListRow_[i].push_back(j);
205     }
206     }
207    
208     groupListCol_.clear();
209 gezelter 1577 groupListCol_.resize(nGroupsInCol_);
210 gezelter 1569 for (int i = 0; i < nGroupsInCol_; i++) {
211     int gid = cgColToGlobal[i];
212     for (int j = 0; j < nAtomsInCol_; j++) {
213     int aid = AtomColToGlobal[j];
214     if (globalGroupMembership[aid] == gid)
215     groupListCol_[i].push_back(j);
216     }
217     }
218    
219 gezelter 1587 excludesForAtom.clear();
220     excludesForAtom.resize(nAtomsInRow_);
221 gezelter 1579 toposForAtom.clear();
222     toposForAtom.resize(nAtomsInRow_);
223     topoDist.clear();
224     topoDist.resize(nAtomsInRow_);
225 gezelter 1570 for (int i = 0; i < nAtomsInRow_; i++) {
226 gezelter 1571 int iglob = AtomRowToGlobal[i];
227 gezelter 1579
228 gezelter 1570 for (int j = 0; j < nAtomsInCol_; j++) {
229 gezelter 1579 int jglob = AtomColToGlobal[j];
230    
231 gezelter 1587 if (excludes->hasPair(iglob, jglob))
232     excludesForAtom[i].push_back(j);
233 gezelter 1579
234 gezelter 1587 if (oneTwo->hasPair(iglob, jglob)) {
235 gezelter 1579 toposForAtom[i].push_back(j);
236     topoDist[i].push_back(1);
237     } else {
238 gezelter 1587 if (oneThree->hasPair(iglob, jglob)) {
239 gezelter 1579 toposForAtom[i].push_back(j);
240     topoDist[i].push_back(2);
241     } else {
242 gezelter 1587 if (oneFour->hasPair(iglob, jglob)) {
243 gezelter 1579 toposForAtom[i].push_back(j);
244     topoDist[i].push_back(3);
245     }
246     }
247 gezelter 1570 }
248     }
249     }
250    
251 gezelter 1613 #else
252 gezelter 1587 excludesForAtom.clear();
253     excludesForAtom.resize(nLocal_);
254 gezelter 1579 toposForAtom.clear();
255     toposForAtom.resize(nLocal_);
256     topoDist.clear();
257     topoDist.resize(nLocal_);
258 gezelter 1569
259 gezelter 1570 for (int i = 0; i < nLocal_; i++) {
260     int iglob = AtomLocalToGlobal[i];
261 gezelter 1579
262 gezelter 1570 for (int j = 0; j < nLocal_; j++) {
263 gezelter 1579 int jglob = AtomLocalToGlobal[j];
264    
265 gezelter 1616 if (excludes->hasPair(iglob, jglob))
266 gezelter 1587 excludesForAtom[i].push_back(j);
267 gezelter 1579
268 gezelter 1587 if (oneTwo->hasPair(iglob, jglob)) {
269 gezelter 1579 toposForAtom[i].push_back(j);
270     topoDist[i].push_back(1);
271     } else {
272 gezelter 1587 if (oneThree->hasPair(iglob, jglob)) {
273 gezelter 1579 toposForAtom[i].push_back(j);
274     topoDist[i].push_back(2);
275     } else {
276 gezelter 1587 if (oneFour->hasPair(iglob, jglob)) {
277 gezelter 1579 toposForAtom[i].push_back(j);
278     topoDist[i].push_back(3);
279     }
280     }
281 gezelter 1570 }
282     }
283 gezelter 1579 }
284 gezelter 1613 #endif
285    
286     // allocate memory for the parallel objects
287     atypesLocal.resize(nLocal_);
288    
289     for (int i = 0; i < nLocal_; i++)
290     atypesLocal[i] = ff_->getAtomType(idents[i]);
291    
292     groupList_.clear();
293     groupList_.resize(nGroups_);
294     for (int i = 0; i < nGroups_; i++) {
295     int gid = cgLocalToGlobal[i];
296     for (int j = 0; j < nLocal_; j++) {
297     int aid = AtomLocalToGlobal[j];
298     if (globalGroupMembership[aid] == gid) {
299     groupList_[i].push_back(j);
300     }
301     }
302     }
303    
304    
305 gezelter 1579 createGtypeCutoffMap();
306 gezelter 1587
307 gezelter 1576 }
308    
309     void ForceMatrixDecomposition::createGtypeCutoffMap() {
310 gezelter 1586
311 gezelter 1896 GrCut.clear();
312     GrCutSq.clear();
313     GrlistSq.clear();
314    
315 gezelter 1576 RealType tol = 1e-6;
316 gezelter 1592 largestRcut_ = 0.0;
317 gezelter 1576 int atid;
318     set<AtomType*> atypes = info_->getSimulatedAtomTypes();
319 gezelter 1592
320 gezelter 1587 map<int, RealType> atypeCutoff;
321 gezelter 1583
322 gezelter 1579 for (set<AtomType*>::iterator at = atypes.begin();
323     at != atypes.end(); ++at){
324 gezelter 1576 atid = (*at)->getIdent();
325 gezelter 1587 if (userChoseCutoff_)
326 gezelter 1583 atypeCutoff[atid] = userCutoff_;
327 gezelter 1592 else
328 gezelter 1583 atypeCutoff[atid] = interactionMan_->getSuggestedCutoffRadius(*at);
329 gezelter 1570 }
330 gezelter 1592
331 gezelter 1576 vector<RealType> gTypeCutoffs;
332     // first we do a single loop over the cutoff groups to find the
333     // largest cutoff for any atypes present in this group.
334     #ifdef IS_MPI
335     vector<RealType> groupCutoffRow(nGroupsInRow_, 0.0);
336 gezelter 1579 groupRowToGtype.resize(nGroupsInRow_);
337 gezelter 1576 for (int cg1 = 0; cg1 < nGroupsInRow_; cg1++) {
338     vector<int> atomListRow = getAtomsInGroupRow(cg1);
339     for (vector<int>::iterator ia = atomListRow.begin();
340     ia != atomListRow.end(); ++ia) {
341     int atom1 = (*ia);
342     atid = identsRow[atom1];
343     if (atypeCutoff[atid] > groupCutoffRow[cg1]) {
344     groupCutoffRow[cg1] = atypeCutoff[atid];
345     }
346     }
347    
348     bool gTypeFound = false;
349     for (int gt = 0; gt < gTypeCutoffs.size(); gt++) {
350     if (abs(groupCutoffRow[cg1] - gTypeCutoffs[gt]) < tol) {
351     groupRowToGtype[cg1] = gt;
352     gTypeFound = true;
353     }
354     }
355     if (!gTypeFound) {
356     gTypeCutoffs.push_back( groupCutoffRow[cg1] );
357     groupRowToGtype[cg1] = gTypeCutoffs.size() - 1;
358     }
359    
360     }
361     vector<RealType> groupCutoffCol(nGroupsInCol_, 0.0);
362 gezelter 1579 groupColToGtype.resize(nGroupsInCol_);
363 gezelter 1576 for (int cg2 = 0; cg2 < nGroupsInCol_; cg2++) {
364     vector<int> atomListCol = getAtomsInGroupColumn(cg2);
365     for (vector<int>::iterator jb = atomListCol.begin();
366     jb != atomListCol.end(); ++jb) {
367     int atom2 = (*jb);
368     atid = identsCol[atom2];
369     if (atypeCutoff[atid] > groupCutoffCol[cg2]) {
370     groupCutoffCol[cg2] = atypeCutoff[atid];
371     }
372     }
373     bool gTypeFound = false;
374     for (int gt = 0; gt < gTypeCutoffs.size(); gt++) {
375     if (abs(groupCutoffCol[cg2] - gTypeCutoffs[gt]) < tol) {
376     groupColToGtype[cg2] = gt;
377     gTypeFound = true;
378     }
379     }
380     if (!gTypeFound) {
381     gTypeCutoffs.push_back( groupCutoffCol[cg2] );
382     groupColToGtype[cg2] = gTypeCutoffs.size() - 1;
383     }
384     }
385     #else
386 gezelter 1579
387 gezelter 1576 vector<RealType> groupCutoff(nGroups_, 0.0);
388 gezelter 1579 groupToGtype.resize(nGroups_);
389 gezelter 1576 for (int cg1 = 0; cg1 < nGroups_; cg1++) {
390     groupCutoff[cg1] = 0.0;
391     vector<int> atomList = getAtomsInGroupRow(cg1);
392     for (vector<int>::iterator ia = atomList.begin();
393     ia != atomList.end(); ++ia) {
394     int atom1 = (*ia);
395 gezelter 1583 atid = idents[atom1];
396 gezelter 1592 if (atypeCutoff[atid] > groupCutoff[cg1])
397 gezelter 1576 groupCutoff[cg1] = atypeCutoff[atid];
398     }
399 gezelter 1592
400 gezelter 1576 bool gTypeFound = false;
401 gezelter 1767 for (unsigned int gt = 0; gt < gTypeCutoffs.size(); gt++) {
402 gezelter 1576 if (abs(groupCutoff[cg1] - gTypeCutoffs[gt]) < tol) {
403     groupToGtype[cg1] = gt;
404     gTypeFound = true;
405     }
406     }
407 gezelter 1592 if (!gTypeFound) {
408 gezelter 1576 gTypeCutoffs.push_back( groupCutoff[cg1] );
409     groupToGtype[cg1] = gTypeCutoffs.size() - 1;
410     }
411     }
412     #endif
413    
414     // Now we find the maximum group cutoff value present in the simulation
415    
416 gezelter 1590 RealType groupMax = *max_element(gTypeCutoffs.begin(),
417     gTypeCutoffs.end());
418 gezelter 1576
419     #ifdef IS_MPI
420 gezelter 1590 MPI::COMM_WORLD.Allreduce(&groupMax, &groupMax, 1, MPI::REALTYPE,
421     MPI::MAX);
422 gezelter 1576 #endif
423    
424     RealType tradRcut = groupMax;
425    
426 gezelter 1896 GrCut.resize( gTypeCutoffs.size() );
427     GrCutSq.resize( gTypeCutoffs.size() );
428     GrlistSq.resize( gTypeCutoffs.size() );
429    
430    
431 gezelter 1767 for (unsigned int i = 0; i < gTypeCutoffs.size(); i++) {
432 gezelter 1896 GrCut[i].resize( gTypeCutoffs.size() , 0.0);
433     GrCutSq[i].resize( gTypeCutoffs.size(), 0.0 );
434     GrlistSq[i].resize( gTypeCutoffs.size(), 0.0 );
435    
436 gezelter 1767 for (unsigned int j = 0; j < gTypeCutoffs.size(); j++) {
437 gezelter 1576 RealType thisRcut;
438     switch(cutoffPolicy_) {
439     case TRADITIONAL:
440     thisRcut = tradRcut;
441 gezelter 1579 break;
442 gezelter 1576 case MIX:
443     thisRcut = 0.5 * (gTypeCutoffs[i] + gTypeCutoffs[j]);
444 gezelter 1579 break;
445 gezelter 1576 case MAX:
446     thisRcut = max(gTypeCutoffs[i], gTypeCutoffs[j]);
447 gezelter 1579 break;
448 gezelter 1576 default:
449     sprintf(painCave.errMsg,
450     "ForceMatrixDecomposition::createGtypeCutoffMap "
451     "hit an unknown cutoff policy!\n");
452     painCave.severity = OPENMD_ERROR;
453     painCave.isFatal = 1;
454 gezelter 1579 simError();
455     break;
456 gezelter 1576 }
457    
458 gezelter 1896 GrCut[i][j] = thisRcut;
459 gezelter 1576 if (thisRcut > largestRcut_) largestRcut_ = thisRcut;
460 gezelter 1896 GrCutSq[i][j] = thisRcut * thisRcut;
461     GrlistSq[i][j] = pow(thisRcut + skinThickness_, 2);
462    
463     // pair<int,int> key = make_pair(i,j);
464     // gTypeCutoffMap[key].first = thisRcut;
465     // gTypeCutoffMap[key].third = pow(thisRcut + skinThickness_, 2);
466 gezelter 1576 // sanity check
467    
468     if (userChoseCutoff_) {
469 gezelter 1896 if (abs(GrCut[i][j] - userCutoff_) > 0.0001) {
470 gezelter 1576 sprintf(painCave.errMsg,
471     "ForceMatrixDecomposition::createGtypeCutoffMap "
472 gezelter 1583 "user-specified rCut (%lf) does not match computed group Cutoff\n", userCutoff_);
473 gezelter 1576 painCave.severity = OPENMD_ERROR;
474     painCave.isFatal = 1;
475     simError();
476     }
477     }
478     }
479     }
480 gezelter 1539 }
481 gezelter 1576
482 gezelter 1896 void ForceMatrixDecomposition::getGroupCutoffs(int &cg1, int &cg2, RealType &rcut, RealType &rcutsq, RealType &rlistsq) {
483 gezelter 1579 int i, j;
484 gezelter 1576 #ifdef IS_MPI
485     i = groupRowToGtype[cg1];
486     j = groupColToGtype[cg2];
487     #else
488     i = groupToGtype[cg1];
489     j = groupToGtype[cg2];
490 gezelter 1579 #endif
491 gezelter 1896 rcut = GrCut[i][j];
492     rcutsq = GrCutSq[i][j];
493     rlistsq = GrlistSq[i][j];
494     return;
495     //return gTypeCutoffMap[make_pair(i,j)];
496 gezelter 1576 }
497    
498 gezelter 1579 int ForceMatrixDecomposition::getTopologicalDistance(int atom1, int atom2) {
499 gezelter 1767 for (unsigned int j = 0; j < toposForAtom[atom1].size(); j++) {
500 gezelter 1579 if (toposForAtom[atom1][j] == atom2)
501     return topoDist[atom1][j];
502 gezelter 1893 }
503 gezelter 1579 return 0;
504     }
505 gezelter 1576
506 gezelter 1575 void ForceMatrixDecomposition::zeroWorkArrays() {
507 gezelter 1583 pairwisePot = 0.0;
508     embeddingPot = 0.0;
509 gezelter 1760 excludedPot = 0.0;
510 gezelter 1761 excludedSelfPot = 0.0;
511 gezelter 1575
512     #ifdef IS_MPI
513     if (storageLayout_ & DataStorage::dslForce) {
514     fill(atomRowData.force.begin(), atomRowData.force.end(), V3Zero);
515     fill(atomColData.force.begin(), atomColData.force.end(), V3Zero);
516     }
517    
518     if (storageLayout_ & DataStorage::dslTorque) {
519     fill(atomRowData.torque.begin(), atomRowData.torque.end(), V3Zero);
520     fill(atomColData.torque.begin(), atomColData.torque.end(), V3Zero);
521     }
522    
523     fill(pot_row.begin(), pot_row.end(),
524     Vector<RealType, N_INTERACTION_FAMILIES> (0.0));
525    
526     fill(pot_col.begin(), pot_col.end(),
527 gezelter 1583 Vector<RealType, N_INTERACTION_FAMILIES> (0.0));
528 gezelter 1575
529 gezelter 1760 fill(expot_row.begin(), expot_row.end(),
530     Vector<RealType, N_INTERACTION_FAMILIES> (0.0));
531    
532     fill(expot_col.begin(), expot_col.end(),
533     Vector<RealType, N_INTERACTION_FAMILIES> (0.0));
534    
535 gezelter 1575 if (storageLayout_ & DataStorage::dslParticlePot) {
536 gezelter 1590 fill(atomRowData.particlePot.begin(), atomRowData.particlePot.end(),
537     0.0);
538     fill(atomColData.particlePot.begin(), atomColData.particlePot.end(),
539     0.0);
540 gezelter 1575 }
541    
542     if (storageLayout_ & DataStorage::dslDensity) {
543     fill(atomRowData.density.begin(), atomRowData.density.end(), 0.0);
544     fill(atomColData.density.begin(), atomColData.density.end(), 0.0);
545     }
546    
547     if (storageLayout_ & DataStorage::dslFunctional) {
548 gezelter 1590 fill(atomRowData.functional.begin(), atomRowData.functional.end(),
549     0.0);
550     fill(atomColData.functional.begin(), atomColData.functional.end(),
551     0.0);
552 gezelter 1575 }
553    
554     if (storageLayout_ & DataStorage::dslFunctionalDerivative) {
555     fill(atomRowData.functionalDerivative.begin(),
556     atomRowData.functionalDerivative.end(), 0.0);
557     fill(atomColData.functionalDerivative.begin(),
558     atomColData.functionalDerivative.end(), 0.0);
559     }
560    
561 gezelter 1586 if (storageLayout_ & DataStorage::dslSkippedCharge) {
562 gezelter 1587 fill(atomRowData.skippedCharge.begin(),
563     atomRowData.skippedCharge.end(), 0.0);
564     fill(atomColData.skippedCharge.begin(),
565     atomColData.skippedCharge.end(), 0.0);
566 gezelter 1586 }
567    
568 gezelter 1721 if (storageLayout_ & DataStorage::dslFlucQForce) {
569     fill(atomRowData.flucQFrc.begin(),
570     atomRowData.flucQFrc.end(), 0.0);
571     fill(atomColData.flucQFrc.begin(),
572     atomColData.flucQFrc.end(), 0.0);
573     }
574    
575 gezelter 1713 if (storageLayout_ & DataStorage::dslElectricField) {
576     fill(atomRowData.electricField.begin(),
577     atomRowData.electricField.end(), V3Zero);
578     fill(atomColData.electricField.begin(),
579     atomColData.electricField.end(), V3Zero);
580     }
581 gezelter 1721
582 gezelter 1590 #endif
583     // even in parallel, we need to zero out the local arrays:
584    
585 gezelter 1575 if (storageLayout_ & DataStorage::dslParticlePot) {
586     fill(snap_->atomData.particlePot.begin(),
587     snap_->atomData.particlePot.end(), 0.0);
588     }
589    
590     if (storageLayout_ & DataStorage::dslDensity) {
591     fill(snap_->atomData.density.begin(),
592     snap_->atomData.density.end(), 0.0);
593     }
594 gezelter 1706
595 gezelter 1575 if (storageLayout_ & DataStorage::dslFunctional) {
596     fill(snap_->atomData.functional.begin(),
597     snap_->atomData.functional.end(), 0.0);
598     }
599 gezelter 1706
600 gezelter 1575 if (storageLayout_ & DataStorage::dslFunctionalDerivative) {
601     fill(snap_->atomData.functionalDerivative.begin(),
602     snap_->atomData.functionalDerivative.end(), 0.0);
603     }
604 gezelter 1706
605 gezelter 1586 if (storageLayout_ & DataStorage::dslSkippedCharge) {
606     fill(snap_->atomData.skippedCharge.begin(),
607     snap_->atomData.skippedCharge.end(), 0.0);
608     }
609 gezelter 1713
610     if (storageLayout_ & DataStorage::dslElectricField) {
611     fill(snap_->atomData.electricField.begin(),
612     snap_->atomData.electricField.end(), V3Zero);
613     }
614 gezelter 1575 }
615    
616    
617 gezelter 1549 void ForceMatrixDecomposition::distributeData() {
618 gezelter 1551 snap_ = sman_->getCurrentSnapshot();
619     storageLayout_ = sman_->getStorageLayout();
620 chuckv 1538 #ifdef IS_MPI
621 gezelter 1540
622 gezelter 1539 // gather up the atomic positions
623 gezelter 1593 AtomPlanVectorRow->gather(snap_->atomData.position,
624 gezelter 1551 atomRowData.position);
625 gezelter 1593 AtomPlanVectorColumn->gather(snap_->atomData.position,
626 gezelter 1551 atomColData.position);
627 gezelter 1539
628     // gather up the cutoff group positions
629 gezelter 1593
630     cgPlanVectorRow->gather(snap_->cgData.position,
631 gezelter 1551 cgRowData.position);
632 gezelter 1593
633     cgPlanVectorColumn->gather(snap_->cgData.position,
634 gezelter 1551 cgColData.position);
635 gezelter 1593
636 gezelter 1723
637    
638     if (needVelocities_) {
639     // gather up the atomic velocities
640     AtomPlanVectorColumn->gather(snap_->atomData.velocity,
641     atomColData.velocity);
642    
643     cgPlanVectorColumn->gather(snap_->cgData.velocity,
644     cgColData.velocity);
645     }
646    
647 gezelter 1539
648     // if needed, gather the atomic rotation matrices
649 gezelter 1551 if (storageLayout_ & DataStorage::dslAmat) {
650 gezelter 1593 AtomPlanMatrixRow->gather(snap_->atomData.aMat,
651 gezelter 1551 atomRowData.aMat);
652 gezelter 1593 AtomPlanMatrixColumn->gather(snap_->atomData.aMat,
653 gezelter 1551 atomColData.aMat);
654 gezelter 1539 }
655 gezelter 1879
656     // if needed, gather the atomic eletrostatic information
657     if (storageLayout_ & DataStorage::dslDipole) {
658     AtomPlanVectorRow->gather(snap_->atomData.dipole,
659     atomRowData.dipole);
660     AtomPlanVectorColumn->gather(snap_->atomData.dipole,
661     atomColData.dipole);
662 gezelter 1539 }
663 gezelter 1590
664 gezelter 1879 if (storageLayout_ & DataStorage::dslQuadrupole) {
665     AtomPlanMatrixRow->gather(snap_->atomData.quadrupole,
666     atomRowData.quadrupole);
667     AtomPlanMatrixColumn->gather(snap_->atomData.quadrupole,
668     atomColData.quadrupole);
669     }
670    
671 gezelter 1713 // if needed, gather the atomic fluctuating charge values
672     if (storageLayout_ & DataStorage::dslFlucQPosition) {
673     AtomPlanRealRow->gather(snap_->atomData.flucQPos,
674     atomRowData.flucQPos);
675     AtomPlanRealColumn->gather(snap_->atomData.flucQPos,
676     atomColData.flucQPos);
677     }
678    
679 gezelter 1539 #endif
680     }
681    
682 gezelter 1575 /* collects information obtained during the pre-pair loop onto local
683     * data structures.
684     */
685 gezelter 1549 void ForceMatrixDecomposition::collectIntermediateData() {
686 gezelter 1551 snap_ = sman_->getCurrentSnapshot();
687     storageLayout_ = sman_->getStorageLayout();
688 gezelter 1539 #ifdef IS_MPI
689    
690 gezelter 1551 if (storageLayout_ & DataStorage::dslDensity) {
691    
692 gezelter 1593 AtomPlanRealRow->scatter(atomRowData.density,
693 gezelter 1551 snap_->atomData.density);
694    
695     int n = snap_->atomData.density.size();
696 gezelter 1575 vector<RealType> rho_tmp(n, 0.0);
697 gezelter 1593 AtomPlanRealColumn->scatter(atomColData.density, rho_tmp);
698 gezelter 1539 for (int i = 0; i < n; i++)
699 gezelter 1551 snap_->atomData.density[i] += rho_tmp[i];
700 gezelter 1539 }
701 gezelter 1713
702 gezelter 1879 // this isn't necessary if we don't have polarizable atoms, but
703     // we'll leave it here for now.
704 gezelter 1713 if (storageLayout_ & DataStorage::dslElectricField) {
705    
706     AtomPlanVectorRow->scatter(atomRowData.electricField,
707     snap_->atomData.electricField);
708    
709     int n = snap_->atomData.electricField.size();
710     vector<Vector3d> field_tmp(n, V3Zero);
711 gezelter 1879 AtomPlanVectorColumn->scatter(atomColData.electricField,
712     field_tmp);
713 gezelter 1713 for (int i = 0; i < n; i++)
714     snap_->atomData.electricField[i] += field_tmp[i];
715     }
716 chuckv 1538 #endif
717 gezelter 1539 }
718 gezelter 1575
719     /*
720     * redistributes information obtained during the pre-pair loop out to
721     * row and column-indexed data structures
722     */
723 gezelter 1549 void ForceMatrixDecomposition::distributeIntermediateData() {
724 gezelter 1551 snap_ = sman_->getCurrentSnapshot();
725     storageLayout_ = sman_->getStorageLayout();
726 chuckv 1538 #ifdef IS_MPI
727 gezelter 1551 if (storageLayout_ & DataStorage::dslFunctional) {
728 gezelter 1593 AtomPlanRealRow->gather(snap_->atomData.functional,
729 gezelter 1551 atomRowData.functional);
730 gezelter 1593 AtomPlanRealColumn->gather(snap_->atomData.functional,
731 gezelter 1551 atomColData.functional);
732 gezelter 1539 }
733    
734 gezelter 1551 if (storageLayout_ & DataStorage::dslFunctionalDerivative) {
735 gezelter 1593 AtomPlanRealRow->gather(snap_->atomData.functionalDerivative,
736 gezelter 1551 atomRowData.functionalDerivative);
737 gezelter 1593 AtomPlanRealColumn->gather(snap_->atomData.functionalDerivative,
738 gezelter 1551 atomColData.functionalDerivative);
739 gezelter 1539 }
740 chuckv 1538 #endif
741     }
742 gezelter 1539
743    
744 gezelter 1549 void ForceMatrixDecomposition::collectData() {
745 gezelter 1551 snap_ = sman_->getCurrentSnapshot();
746     storageLayout_ = sman_->getStorageLayout();
747     #ifdef IS_MPI
748     int n = snap_->atomData.force.size();
749 gezelter 1544 vector<Vector3d> frc_tmp(n, V3Zero);
750 gezelter 1541
751 gezelter 1593 AtomPlanVectorRow->scatter(atomRowData.force, frc_tmp);
752 gezelter 1541 for (int i = 0; i < n; i++) {
753 gezelter 1551 snap_->atomData.force[i] += frc_tmp[i];
754 gezelter 1541 frc_tmp[i] = 0.0;
755     }
756 gezelter 1540
757 gezelter 1593 AtomPlanVectorColumn->scatter(atomColData.force, frc_tmp);
758     for (int i = 0; i < n; i++) {
759 gezelter 1551 snap_->atomData.force[i] += frc_tmp[i];
760 gezelter 1593 }
761 gezelter 1590
762 gezelter 1551 if (storageLayout_ & DataStorage::dslTorque) {
763 gezelter 1541
764 gezelter 1587 int nt = snap_->atomData.torque.size();
765 gezelter 1544 vector<Vector3d> trq_tmp(nt, V3Zero);
766 gezelter 1541
767 gezelter 1593 AtomPlanVectorRow->scatter(atomRowData.torque, trq_tmp);
768 gezelter 1587 for (int i = 0; i < nt; i++) {
769 gezelter 1551 snap_->atomData.torque[i] += trq_tmp[i];
770 gezelter 1541 trq_tmp[i] = 0.0;
771     }
772 gezelter 1540
773 gezelter 1593 AtomPlanVectorColumn->scatter(atomColData.torque, trq_tmp);
774 gezelter 1587 for (int i = 0; i < nt; i++)
775 gezelter 1551 snap_->atomData.torque[i] += trq_tmp[i];
776 gezelter 1540 }
777 gezelter 1587
778     if (storageLayout_ & DataStorage::dslSkippedCharge) {
779    
780     int ns = snap_->atomData.skippedCharge.size();
781     vector<RealType> skch_tmp(ns, 0.0);
782    
783 gezelter 1593 AtomPlanRealRow->scatter(atomRowData.skippedCharge, skch_tmp);
784 gezelter 1587 for (int i = 0; i < ns; i++) {
785 gezelter 1590 snap_->atomData.skippedCharge[i] += skch_tmp[i];
786 gezelter 1587 skch_tmp[i] = 0.0;
787     }
788    
789 gezelter 1593 AtomPlanRealColumn->scatter(atomColData.skippedCharge, skch_tmp);
790 gezelter 1613 for (int i = 0; i < ns; i++)
791 gezelter 1587 snap_->atomData.skippedCharge[i] += skch_tmp[i];
792 gezelter 1613
793 gezelter 1587 }
794 gezelter 1540
795 gezelter 1713 if (storageLayout_ & DataStorage::dslFlucQForce) {
796    
797     int nq = snap_->atomData.flucQFrc.size();
798     vector<RealType> fqfrc_tmp(nq, 0.0);
799    
800     AtomPlanRealRow->scatter(atomRowData.flucQFrc, fqfrc_tmp);
801     for (int i = 0; i < nq; i++) {
802     snap_->atomData.flucQFrc[i] += fqfrc_tmp[i];
803     fqfrc_tmp[i] = 0.0;
804     }
805    
806     AtomPlanRealColumn->scatter(atomColData.flucQFrc, fqfrc_tmp);
807     for (int i = 0; i < nq; i++)
808     snap_->atomData.flucQFrc[i] += fqfrc_tmp[i];
809    
810     }
811    
812 gezelter 1879 if (storageLayout_ & DataStorage::dslElectricField) {
813    
814     int nef = snap_->atomData.electricField.size();
815     vector<Vector3d> efield_tmp(nef, V3Zero);
816    
817     AtomPlanVectorRow->scatter(atomRowData.electricField, efield_tmp);
818     for (int i = 0; i < nef; i++) {
819     snap_->atomData.electricField[i] += efield_tmp[i];
820     efield_tmp[i] = 0.0;
821     }
822    
823     AtomPlanVectorColumn->scatter(atomColData.electricField, efield_tmp);
824     for (int i = 0; i < nef; i++)
825     snap_->atomData.electricField[i] += efield_tmp[i];
826     }
827    
828    
829 gezelter 1567 nLocal_ = snap_->getNumberOfAtoms();
830 gezelter 1544
831 gezelter 1575 vector<potVec> pot_temp(nLocal_,
832     Vector<RealType, N_INTERACTION_FAMILIES> (0.0));
833 gezelter 1760 vector<potVec> expot_temp(nLocal_,
834     Vector<RealType, N_INTERACTION_FAMILIES> (0.0));
835 gezelter 1575
836     // scatter/gather pot_row into the members of my column
837    
838 gezelter 1593 AtomPlanPotRow->scatter(pot_row, pot_temp);
839 gezelter 1760 AtomPlanPotRow->scatter(expot_row, expot_temp);
840 gezelter 1575
841 gezelter 1760 for (int ii = 0; ii < pot_temp.size(); ii++ )
842 gezelter 1583 pairwisePot += pot_temp[ii];
843 gezelter 1760
844     for (int ii = 0; ii < expot_temp.size(); ii++ )
845     excludedPot += expot_temp[ii];
846 gezelter 1723
847     if (storageLayout_ & DataStorage::dslParticlePot) {
848     // This is the pairwise contribution to the particle pot. The
849     // embedding contribution is added in each of the low level
850     // non-bonded routines. In single processor, this is done in
851     // unpackInteractionData, not in collectData.
852     for (int ii = 0; ii < N_INTERACTION_FAMILIES; ii++) {
853     for (int i = 0; i < nLocal_; i++) {
854     // factor of two is because the total potential terms are divided
855     // by 2 in parallel due to row/ column scatter
856     snap_->atomData.particlePot[i] += 2.0 * pot_temp[i](ii);
857     }
858     }
859     }
860    
861 gezelter 1575 fill(pot_temp.begin(), pot_temp.end(),
862     Vector<RealType, N_INTERACTION_FAMILIES> (0.0));
863 gezelter 1760 fill(expot_temp.begin(), expot_temp.end(),
864     Vector<RealType, N_INTERACTION_FAMILIES> (0.0));
865 gezelter 1575
866 gezelter 1593 AtomPlanPotColumn->scatter(pot_col, pot_temp);
867 gezelter 1760 AtomPlanPotColumn->scatter(expot_col, expot_temp);
868 gezelter 1575
869     for (int ii = 0; ii < pot_temp.size(); ii++ )
870 gezelter 1583 pairwisePot += pot_temp[ii];
871 gezelter 1723
872 gezelter 1760 for (int ii = 0; ii < expot_temp.size(); ii++ )
873     excludedPot += expot_temp[ii];
874    
875 gezelter 1723 if (storageLayout_ & DataStorage::dslParticlePot) {
876     // This is the pairwise contribution to the particle pot. The
877     // embedding contribution is added in each of the low level
878     // non-bonded routines. In single processor, this is done in
879     // unpackInteractionData, not in collectData.
880     for (int ii = 0; ii < N_INTERACTION_FAMILIES; ii++) {
881     for (int i = 0; i < nLocal_; i++) {
882     // factor of two is because the total potential terms are divided
883     // by 2 in parallel due to row/ column scatter
884     snap_->atomData.particlePot[i] += 2.0 * pot_temp[i](ii);
885     }
886     }
887     }
888 gezelter 1601
889 gezelter 1723 if (storageLayout_ & DataStorage::dslParticlePot) {
890     int npp = snap_->atomData.particlePot.size();
891     vector<RealType> ppot_temp(npp, 0.0);
892    
893     // This is the direct or embedding contribution to the particle
894     // pot.
895    
896     AtomPlanRealRow->scatter(atomRowData.particlePot, ppot_temp);
897     for (int i = 0; i < npp; i++) {
898     snap_->atomData.particlePot[i] += ppot_temp[i];
899     }
900    
901     fill(ppot_temp.begin(), ppot_temp.end(), 0.0);
902    
903     AtomPlanRealColumn->scatter(atomColData.particlePot, ppot_temp);
904     for (int i = 0; i < npp; i++) {
905     snap_->atomData.particlePot[i] += ppot_temp[i];
906     }
907     }
908    
909 gezelter 1601 for (int ii = 0; ii < N_INTERACTION_FAMILIES; ii++) {
910     RealType ploc1 = pairwisePot[ii];
911     RealType ploc2 = 0.0;
912     MPI::COMM_WORLD.Allreduce(&ploc1, &ploc2, 1, MPI::REALTYPE, MPI::SUM);
913     pairwisePot[ii] = ploc2;
914     }
915    
916 gezelter 1760 for (int ii = 0; ii < N_INTERACTION_FAMILIES; ii++) {
917     RealType ploc1 = excludedPot[ii];
918     RealType ploc2 = 0.0;
919     MPI::COMM_WORLD.Allreduce(&ploc1, &ploc2, 1, MPI::REALTYPE, MPI::SUM);
920     excludedPot[ii] = ploc2;
921     }
922    
923 gezelter 1723 // Here be dragons.
924     MPI::Intracomm col = colComm.getComm();
925 gezelter 1613
926 gezelter 1723 col.Allreduce(MPI::IN_PLACE,
927     &snap_->frameData.conductiveHeatFlux[0], 3,
928     MPI::REALTYPE, MPI::SUM);
929    
930    
931 gezelter 1539 #endif
932 gezelter 1583
933 chuckv 1538 }
934 gezelter 1551
935 gezelter 1756 /**
936     * Collects information obtained during the post-pair (and embedding
937     * functional) loops onto local data structures.
938     */
939     void ForceMatrixDecomposition::collectSelfData() {
940     snap_ = sman_->getCurrentSnapshot();
941     storageLayout_ = sman_->getStorageLayout();
942    
943     #ifdef IS_MPI
944     for (int ii = 0; ii < N_INTERACTION_FAMILIES; ii++) {
945     RealType ploc1 = embeddingPot[ii];
946     RealType ploc2 = 0.0;
947     MPI::COMM_WORLD.Allreduce(&ploc1, &ploc2, 1, MPI::REALTYPE, MPI::SUM);
948     embeddingPot[ii] = ploc2;
949     }
950 gezelter 1761 for (int ii = 0; ii < N_INTERACTION_FAMILIES; ii++) {
951     RealType ploc1 = excludedSelfPot[ii];
952     RealType ploc2 = 0.0;
953     MPI::COMM_WORLD.Allreduce(&ploc1, &ploc2, 1, MPI::REALTYPE, MPI::SUM);
954     excludedSelfPot[ii] = ploc2;
955     }
956 gezelter 1756 #endif
957    
958     }
959    
960    
961    
962 gezelter 1893 int& ForceMatrixDecomposition::getNAtomsInRow() {
963 gezelter 1570 #ifdef IS_MPI
964     return nAtomsInRow_;
965     #else
966     return nLocal_;
967     #endif
968     }
969    
970 gezelter 1569 /**
971     * returns the list of atoms belonging to this group.
972     */
973 gezelter 1893 vector<int>& ForceMatrixDecomposition::getAtomsInGroupRow(int cg1){
974 gezelter 1569 #ifdef IS_MPI
975     return groupListRow_[cg1];
976     #else
977     return groupList_[cg1];
978     #endif
979     }
980    
981 gezelter 1893 vector<int>& ForceMatrixDecomposition::getAtomsInGroupColumn(int cg2){
982 gezelter 1569 #ifdef IS_MPI
983     return groupListCol_[cg2];
984     #else
985     return groupList_[cg2];
986     #endif
987     }
988 chuckv 1538
989 gezelter 1551 Vector3d ForceMatrixDecomposition::getIntergroupVector(int cg1, int cg2){
990     Vector3d d;
991    
992     #ifdef IS_MPI
993     d = cgColData.position[cg2] - cgRowData.position[cg1];
994     #else
995     d = snap_->cgData.position[cg2] - snap_->cgData.position[cg1];
996     #endif
997    
998 gezelter 1879 if (usePeriodicBoundaryConditions_) {
999     snap_->wrapVector(d);
1000     }
1001 gezelter 1551 return d;
1002     }
1003    
1004 gezelter 1893 Vector3d& ForceMatrixDecomposition::getGroupVelocityColumn(int cg2){
1005 gezelter 1723 #ifdef IS_MPI
1006     return cgColData.velocity[cg2];
1007     #else
1008     return snap_->cgData.velocity[cg2];
1009     #endif
1010     }
1011 gezelter 1551
1012 gezelter 1893 Vector3d& ForceMatrixDecomposition::getAtomVelocityColumn(int atom2){
1013 gezelter 1723 #ifdef IS_MPI
1014     return atomColData.velocity[atom2];
1015     #else
1016     return snap_->atomData.velocity[atom2];
1017     #endif
1018     }
1019    
1020    
1021 gezelter 1551 Vector3d ForceMatrixDecomposition::getAtomToGroupVectorRow(int atom1, int cg1){
1022    
1023     Vector3d d;
1024    
1025     #ifdef IS_MPI
1026     d = cgRowData.position[cg1] - atomRowData.position[atom1];
1027     #else
1028     d = snap_->cgData.position[cg1] - snap_->atomData.position[atom1];
1029     #endif
1030 gezelter 1879 if (usePeriodicBoundaryConditions_) {
1031     snap_->wrapVector(d);
1032     }
1033 gezelter 1551 return d;
1034     }
1035    
1036     Vector3d ForceMatrixDecomposition::getAtomToGroupVectorColumn(int atom2, int cg2){
1037     Vector3d d;
1038    
1039     #ifdef IS_MPI
1040     d = cgColData.position[cg2] - atomColData.position[atom2];
1041     #else
1042     d = snap_->cgData.position[cg2] - snap_->atomData.position[atom2];
1043     #endif
1044 gezelter 1879 if (usePeriodicBoundaryConditions_) {
1045     snap_->wrapVector(d);
1046     }
1047 gezelter 1551 return d;
1048     }
1049 gezelter 1569
1050 gezelter 1893 RealType& ForceMatrixDecomposition::getMassFactorRow(int atom1) {
1051 gezelter 1569 #ifdef IS_MPI
1052     return massFactorsRow[atom1];
1053     #else
1054 gezelter 1581 return massFactors[atom1];
1055 gezelter 1569 #endif
1056     }
1057    
1058 gezelter 1893 RealType& ForceMatrixDecomposition::getMassFactorColumn(int atom2) {
1059 gezelter 1569 #ifdef IS_MPI
1060     return massFactorsCol[atom2];
1061     #else
1062 gezelter 1581 return massFactors[atom2];
1063 gezelter 1569 #endif
1064    
1065     }
1066 gezelter 1551
1067     Vector3d ForceMatrixDecomposition::getInteratomicVector(int atom1, int atom2){
1068     Vector3d d;
1069    
1070     #ifdef IS_MPI
1071     d = atomColData.position[atom2] - atomRowData.position[atom1];
1072     #else
1073     d = snap_->atomData.position[atom2] - snap_->atomData.position[atom1];
1074     #endif
1075 gezelter 1879 if (usePeriodicBoundaryConditions_) {
1076     snap_->wrapVector(d);
1077     }
1078 gezelter 1551 return d;
1079     }
1080    
1081 gezelter 1893 vector<int>& ForceMatrixDecomposition::getExcludesForAtom(int atom1) {
1082 gezelter 1587 return excludesForAtom[atom1];
1083 gezelter 1570 }
1084    
1085     /**
1086 gezelter 1587 * We need to exclude some overcounted interactions that result from
1087 gezelter 1575 * the parallel decomposition.
1088 gezelter 1570 */
1089 gezelter 1756 bool ForceMatrixDecomposition::skipAtomPair(int atom1, int atom2, int cg1, int cg2) {
1090 gezelter 1796 int unique_id_1, unique_id_2;
1091 gezelter 1616
1092 gezelter 1570 #ifdef IS_MPI
1093     // in MPI, we have to look up the unique IDs for each atom
1094     unique_id_1 = AtomRowToGlobal[atom1];
1095     unique_id_2 = AtomColToGlobal[atom2];
1096 gezelter 1796 // group1 = cgRowToGlobal[cg1];
1097     // group2 = cgColToGlobal[cg2];
1098 gezelter 1616 #else
1099     unique_id_1 = AtomLocalToGlobal[atom1];
1100     unique_id_2 = AtomLocalToGlobal[atom2];
1101 gezelter 1796 int group1 = cgLocalToGlobal[cg1];
1102     int group2 = cgLocalToGlobal[cg2];
1103 gezelter 1616 #endif
1104 gezelter 1570
1105     if (unique_id_1 == unique_id_2) return true;
1106 gezelter 1616
1107     #ifdef IS_MPI
1108 gezelter 1570 // this prevents us from doing the pair on multiple processors
1109     if (unique_id_1 < unique_id_2) {
1110     if ((unique_id_1 + unique_id_2) % 2 == 0) return true;
1111     } else {
1112 gezelter 1616 if ((unique_id_1 + unique_id_2) % 2 == 1) return true;
1113 gezelter 1570 }
1114 gezelter 1756 #endif
1115    
1116     #ifndef IS_MPI
1117     if (group1 == group2) {
1118     if (unique_id_1 < unique_id_2) return true;
1119     }
1120 gezelter 1587 #endif
1121 gezelter 1616
1122 gezelter 1587 return false;
1123     }
1124    
1125     /**
1126     * We need to handle the interactions for atoms who are involved in
1127     * the same rigid body as well as some short range interactions
1128     * (bonds, bends, torsions) differently from other interactions.
1129     * We'll still visit the pairwise routines, but with a flag that
1130     * tells those routines to exclude the pair from direct long range
1131     * interactions. Some indirect interactions (notably reaction
1132     * field) must still be handled for these pairs.
1133     */
1134     bool ForceMatrixDecomposition::excludeAtomPair(int atom1, int atom2) {
1135 gezelter 1613
1136     // excludesForAtom was constructed to use row/column indices in the MPI
1137     // version, and to use local IDs in the non-MPI version:
1138 gezelter 1570
1139 gezelter 1587 for (vector<int>::iterator i = excludesForAtom[atom1].begin();
1140     i != excludesForAtom[atom1].end(); ++i) {
1141 gezelter 1616 if ( (*i) == atom2 ) return true;
1142 gezelter 1583 }
1143 gezelter 1579
1144 gezelter 1583 return false;
1145 gezelter 1570 }
1146    
1147    
1148 gezelter 1551 void ForceMatrixDecomposition::addForceToAtomRow(int atom1, Vector3d fg){
1149     #ifdef IS_MPI
1150     atomRowData.force[atom1] += fg;
1151     #else
1152     snap_->atomData.force[atom1] += fg;
1153     #endif
1154     }
1155    
1156     void ForceMatrixDecomposition::addForceToAtomColumn(int atom2, Vector3d fg){
1157     #ifdef IS_MPI
1158     atomColData.force[atom2] += fg;
1159     #else
1160     snap_->atomData.force[atom2] += fg;
1161     #endif
1162     }
1163    
1164     // filling interaction blocks with pointers
1165 gezelter 1582 void ForceMatrixDecomposition::fillInteractionData(InteractionData &idat,
1166 gezelter 1587 int atom1, int atom2) {
1167    
1168     idat.excluded = excludeAtomPair(atom1, atom2);
1169    
1170 gezelter 1551 #ifdef IS_MPI
1171 gezelter 1591 idat.atypes = make_pair( atypesRow[atom1], atypesCol[atom2]);
1172 gezelter 1895 idat.atid1 = identsRow[atom1];
1173     idat.atid2 = identsCol[atom2];
1174 gezelter 1591 //idat.atypes = make_pair( ff_->getAtomType(identsRow[atom1]),
1175     // ff_->getAtomType(identsCol[atom2]) );
1176 gezelter 1571
1177 gezelter 1551 if (storageLayout_ & DataStorage::dslAmat) {
1178 gezelter 1554 idat.A1 = &(atomRowData.aMat[atom1]);
1179     idat.A2 = &(atomColData.aMat[atom2]);
1180 gezelter 1551 }
1181 gezelter 1567
1182 gezelter 1551 if (storageLayout_ & DataStorage::dslTorque) {
1183 gezelter 1554 idat.t1 = &(atomRowData.torque[atom1]);
1184     idat.t2 = &(atomColData.torque[atom2]);
1185 gezelter 1551 }
1186    
1187 gezelter 1879 if (storageLayout_ & DataStorage::dslDipole) {
1188     idat.dipole1 = &(atomRowData.dipole[atom1]);
1189     idat.dipole2 = &(atomColData.dipole[atom2]);
1190     }
1191    
1192     if (storageLayout_ & DataStorage::dslQuadrupole) {
1193     idat.quadrupole1 = &(atomRowData.quadrupole[atom1]);
1194     idat.quadrupole2 = &(atomColData.quadrupole[atom2]);
1195     }
1196    
1197 gezelter 1551 if (storageLayout_ & DataStorage::dslDensity) {
1198 gezelter 1554 idat.rho1 = &(atomRowData.density[atom1]);
1199     idat.rho2 = &(atomColData.density[atom2]);
1200 gezelter 1551 }
1201    
1202 gezelter 1575 if (storageLayout_ & DataStorage::dslFunctional) {
1203     idat.frho1 = &(atomRowData.functional[atom1]);
1204     idat.frho2 = &(atomColData.functional[atom2]);
1205     }
1206    
1207 gezelter 1551 if (storageLayout_ & DataStorage::dslFunctionalDerivative) {
1208 gezelter 1554 idat.dfrho1 = &(atomRowData.functionalDerivative[atom1]);
1209     idat.dfrho2 = &(atomColData.functionalDerivative[atom2]);
1210 gezelter 1551 }
1211 gezelter 1570
1212 gezelter 1575 if (storageLayout_ & DataStorage::dslParticlePot) {
1213     idat.particlePot1 = &(atomRowData.particlePot[atom1]);
1214     idat.particlePot2 = &(atomColData.particlePot[atom2]);
1215     }
1216    
1217 gezelter 1587 if (storageLayout_ & DataStorage::dslSkippedCharge) {
1218     idat.skippedCharge1 = &(atomRowData.skippedCharge[atom1]);
1219     idat.skippedCharge2 = &(atomColData.skippedCharge[atom2]);
1220     }
1221    
1222 gezelter 1721 if (storageLayout_ & DataStorage::dslFlucQPosition) {
1223     idat.flucQ1 = &(atomRowData.flucQPos[atom1]);
1224     idat.flucQ2 = &(atomColData.flucQPos[atom2]);
1225     }
1226    
1227 gezelter 1562 #else
1228 gezelter 1688
1229 gezelter 1591 idat.atypes = make_pair( atypesLocal[atom1], atypesLocal[atom2]);
1230 gezelter 1895 idat.atid1 = idents[atom1];
1231     idat.atid2 = idents[atom2];
1232 gezelter 1571
1233 gezelter 1562 if (storageLayout_ & DataStorage::dslAmat) {
1234     idat.A1 = &(snap_->atomData.aMat[atom1]);
1235     idat.A2 = &(snap_->atomData.aMat[atom2]);
1236     }
1237    
1238     if (storageLayout_ & DataStorage::dslTorque) {
1239     idat.t1 = &(snap_->atomData.torque[atom1]);
1240     idat.t2 = &(snap_->atomData.torque[atom2]);
1241     }
1242    
1243 gezelter 1879 if (storageLayout_ & DataStorage::dslDipole) {
1244     idat.dipole1 = &(snap_->atomData.dipole[atom1]);
1245     idat.dipole2 = &(snap_->atomData.dipole[atom2]);
1246     }
1247    
1248     if (storageLayout_ & DataStorage::dslQuadrupole) {
1249     idat.quadrupole1 = &(snap_->atomData.quadrupole[atom1]);
1250     idat.quadrupole2 = &(snap_->atomData.quadrupole[atom2]);
1251     }
1252    
1253 gezelter 1583 if (storageLayout_ & DataStorage::dslDensity) {
1254 gezelter 1562 idat.rho1 = &(snap_->atomData.density[atom1]);
1255     idat.rho2 = &(snap_->atomData.density[atom2]);
1256     }
1257    
1258 gezelter 1575 if (storageLayout_ & DataStorage::dslFunctional) {
1259     idat.frho1 = &(snap_->atomData.functional[atom1]);
1260     idat.frho2 = &(snap_->atomData.functional[atom2]);
1261     }
1262    
1263 gezelter 1562 if (storageLayout_ & DataStorage::dslFunctionalDerivative) {
1264     idat.dfrho1 = &(snap_->atomData.functionalDerivative[atom1]);
1265     idat.dfrho2 = &(snap_->atomData.functionalDerivative[atom2]);
1266     }
1267 gezelter 1575
1268     if (storageLayout_ & DataStorage::dslParticlePot) {
1269     idat.particlePot1 = &(snap_->atomData.particlePot[atom1]);
1270     idat.particlePot2 = &(snap_->atomData.particlePot[atom2]);
1271     }
1272    
1273 gezelter 1587 if (storageLayout_ & DataStorage::dslSkippedCharge) {
1274     idat.skippedCharge1 = &(snap_->atomData.skippedCharge[atom1]);
1275     idat.skippedCharge2 = &(snap_->atomData.skippedCharge[atom2]);
1276     }
1277 gezelter 1721
1278     if (storageLayout_ & DataStorage::dslFlucQPosition) {
1279     idat.flucQ1 = &(snap_->atomData.flucQPos[atom1]);
1280     idat.flucQ2 = &(snap_->atomData.flucQPos[atom2]);
1281     }
1282    
1283 gezelter 1551 #endif
1284     }
1285 gezelter 1567
1286 gezelter 1575
1287 gezelter 1582 void ForceMatrixDecomposition::unpackInteractionData(InteractionData &idat, int atom1, int atom2) {
1288 gezelter 1575 #ifdef IS_MPI
1289 gezelter 1668 pot_row[atom1] += RealType(0.5) * *(idat.pot);
1290     pot_col[atom2] += RealType(0.5) * *(idat.pot);
1291 gezelter 1760 expot_row[atom1] += RealType(0.5) * *(idat.excludedPot);
1292     expot_col[atom2] += RealType(0.5) * *(idat.excludedPot);
1293 gezelter 1575
1294     atomRowData.force[atom1] += *(idat.f1);
1295     atomColData.force[atom2] -= *(idat.f1);
1296 gezelter 1713
1297 gezelter 1721 if (storageLayout_ & DataStorage::dslFlucQForce) {
1298 jmichalk 1736 atomRowData.flucQFrc[atom1] -= *(idat.dVdFQ1);
1299     atomColData.flucQFrc[atom2] -= *(idat.dVdFQ2);
1300 gezelter 1721 }
1301    
1302     if (storageLayout_ & DataStorage::dslElectricField) {
1303     atomRowData.electricField[atom1] += *(idat.eField1);
1304     atomColData.electricField[atom2] += *(idat.eField2);
1305     }
1306    
1307 gezelter 1575 #else
1308 gezelter 1583 pairwisePot += *(idat.pot);
1309 gezelter 1760 excludedPot += *(idat.excludedPot);
1310 gezelter 1583
1311 gezelter 1575 snap_->atomData.force[atom1] += *(idat.f1);
1312     snap_->atomData.force[atom2] -= *(idat.f1);
1313 gezelter 1713
1314     if (idat.doParticlePot) {
1315 gezelter 1723 // This is the pairwise contribution to the particle pot. The
1316     // embedding contribution is added in each of the low level
1317     // non-bonded routines. In parallel, this calculation is done
1318     // in collectData, not in unpackInteractionData.
1319 gezelter 1713 snap_->atomData.particlePot[atom1] += *(idat.vpair) * *(idat.sw);
1320 gezelter 1723 snap_->atomData.particlePot[atom2] += *(idat.vpair) * *(idat.sw);
1321 gezelter 1713 }
1322 gezelter 1721
1323     if (storageLayout_ & DataStorage::dslFlucQForce) {
1324 jmichalk 1736 snap_->atomData.flucQFrc[atom1] -= *(idat.dVdFQ1);
1325 gezelter 1721 snap_->atomData.flucQFrc[atom2] -= *(idat.dVdFQ2);
1326     }
1327    
1328     if (storageLayout_ & DataStorage::dslElectricField) {
1329     snap_->atomData.electricField[atom1] += *(idat.eField1);
1330     snap_->atomData.electricField[atom2] += *(idat.eField2);
1331     }
1332    
1333 gezelter 1575 #endif
1334 gezelter 1586
1335 gezelter 1575 }
1336    
1337 gezelter 1562 /*
1338     * buildNeighborList
1339     *
1340     * first element of pair is row-indexed CutoffGroup
1341     * second element of pair is column-indexed CutoffGroup
1342     */
1343 gezelter 1895 void ForceMatrixDecomposition::buildNeighborList(vector<pair<int,int> >& neighborList) {
1344    
1345     neighborList.clear();
1346 gezelter 1576 groupCutoffs cuts;
1347 gezelter 1587 bool doAllPairs = false;
1348    
1349 gezelter 1879 RealType rList_ = (largestRcut_ + skinThickness_);
1350 gezelter 1896 RealType rcut, rcutsq, rlistsq;
1351 gezelter 1879 Snapshot* snap_ = sman_->getCurrentSnapshot();
1352     Mat3x3d box;
1353     Mat3x3d invBox;
1354    
1355     Vector3d rs, scaled, dr;
1356     Vector3i whichCell;
1357     int cellIndex;
1358    
1359 gezelter 1567 #ifdef IS_MPI
1360 gezelter 1568 cellListRow_.clear();
1361     cellListCol_.clear();
1362 gezelter 1567 #else
1363 gezelter 1568 cellList_.clear();
1364 gezelter 1567 #endif
1365 gezelter 1879
1366     if (!usePeriodicBoundaryConditions_) {
1367     box = snap_->getBoundingBox();
1368     invBox = snap_->getInvBoundingBox();
1369     } else {
1370     box = snap_->getHmat();
1371     invBox = snap_->getInvHmat();
1372     }
1373    
1374     Vector3d boxX = box.getColumn(0);
1375     Vector3d boxY = box.getColumn(1);
1376     Vector3d boxZ = box.getColumn(2);
1377    
1378     nCells_.x() = (int) ( boxX.length() )/ rList_;
1379     nCells_.y() = (int) ( boxY.length() )/ rList_;
1380     nCells_.z() = (int) ( boxZ.length() )/ rList_;
1381    
1382 gezelter 1587 // handle small boxes where the cell offsets can end up repeating cells
1383    
1384     if (nCells_.x() < 3) doAllPairs = true;
1385     if (nCells_.y() < 3) doAllPairs = true;
1386     if (nCells_.z() < 3) doAllPairs = true;
1387 gezelter 1879
1388 gezelter 1579 int nCtot = nCells_.x() * nCells_.y() * nCells_.z();
1389 gezelter 1879
1390 gezelter 1567 #ifdef IS_MPI
1391 gezelter 1579 cellListRow_.resize(nCtot);
1392     cellListCol_.resize(nCtot);
1393     #else
1394     cellList_.resize(nCtot);
1395     #endif
1396 gezelter 1879
1397 gezelter 1587 if (!doAllPairs) {
1398 gezelter 1579 #ifdef IS_MPI
1399 gezelter 1879
1400 gezelter 1587 for (int i = 0; i < nGroupsInRow_; i++) {
1401     rs = cgRowData.position[i];
1402    
1403     // scaled positions relative to the box vectors
1404 gezelter 1879 scaled = invBox * rs;
1405 gezelter 1587
1406     // wrap the vector back into the unit box by subtracting integer box
1407     // numbers
1408     for (int j = 0; j < 3; j++) {
1409     scaled[j] -= roundMe(scaled[j]);
1410     scaled[j] += 0.5;
1411 gezelter 1772 // Handle the special case when an object is exactly on the
1412     // boundary (a scaled coordinate of 1.0 is the same as
1413     // scaled coordinate of 0.0)
1414     if (scaled[j] >= 1.0) scaled[j] -= 1.0;
1415 gezelter 1587 }
1416    
1417     // find xyz-indices of cell that cutoffGroup is in.
1418     whichCell.x() = nCells_.x() * scaled.x();
1419     whichCell.y() = nCells_.y() * scaled.y();
1420     whichCell.z() = nCells_.z() * scaled.z();
1421    
1422     // find single index of this cell:
1423     cellIndex = Vlinear(whichCell, nCells_);
1424    
1425     // add this cutoff group to the list of groups in this cell;
1426     cellListRow_[cellIndex].push_back(i);
1427 gezelter 1581 }
1428 gezelter 1587 for (int i = 0; i < nGroupsInCol_; i++) {
1429     rs = cgColData.position[i];
1430    
1431     // scaled positions relative to the box vectors
1432 gezelter 1879 scaled = invBox * rs;
1433 gezelter 1587
1434     // wrap the vector back into the unit box by subtracting integer box
1435     // numbers
1436     for (int j = 0; j < 3; j++) {
1437     scaled[j] -= roundMe(scaled[j]);
1438     scaled[j] += 0.5;
1439 gezelter 1772 // Handle the special case when an object is exactly on the
1440     // boundary (a scaled coordinate of 1.0 is the same as
1441     // scaled coordinate of 0.0)
1442     if (scaled[j] >= 1.0) scaled[j] -= 1.0;
1443 gezelter 1587 }
1444    
1445     // find xyz-indices of cell that cutoffGroup is in.
1446     whichCell.x() = nCells_.x() * scaled.x();
1447     whichCell.y() = nCells_.y() * scaled.y();
1448     whichCell.z() = nCells_.z() * scaled.z();
1449    
1450     // find single index of this cell:
1451     cellIndex = Vlinear(whichCell, nCells_);
1452    
1453     // add this cutoff group to the list of groups in this cell;
1454     cellListCol_[cellIndex].push_back(i);
1455 gezelter 1581 }
1456 gezelter 1879
1457 gezelter 1567 #else
1458 gezelter 1587 for (int i = 0; i < nGroups_; i++) {
1459     rs = snap_->cgData.position[i];
1460    
1461     // scaled positions relative to the box vectors
1462 gezelter 1879 scaled = invBox * rs;
1463 gezelter 1587
1464     // wrap the vector back into the unit box by subtracting integer box
1465     // numbers
1466     for (int j = 0; j < 3; j++) {
1467     scaled[j] -= roundMe(scaled[j]);
1468     scaled[j] += 0.5;
1469 gezelter 1771 // Handle the special case when an object is exactly on the
1470     // boundary (a scaled coordinate of 1.0 is the same as
1471     // scaled coordinate of 0.0)
1472     if (scaled[j] >= 1.0) scaled[j] -= 1.0;
1473 gezelter 1587 }
1474    
1475     // find xyz-indices of cell that cutoffGroup is in.
1476     whichCell.x() = nCells_.x() * scaled.x();
1477     whichCell.y() = nCells_.y() * scaled.y();
1478     whichCell.z() = nCells_.z() * scaled.z();
1479    
1480     // find single index of this cell:
1481 gezelter 1593 cellIndex = Vlinear(whichCell, nCells_);
1482 gezelter 1587
1483     // add this cutoff group to the list of groups in this cell;
1484     cellList_[cellIndex].push_back(i);
1485 gezelter 1581 }
1486 gezelter 1612
1487 gezelter 1567 #endif
1488    
1489 gezelter 1587 for (int m1z = 0; m1z < nCells_.z(); m1z++) {
1490     for (int m1y = 0; m1y < nCells_.y(); m1y++) {
1491     for (int m1x = 0; m1x < nCells_.x(); m1x++) {
1492     Vector3i m1v(m1x, m1y, m1z);
1493     int m1 = Vlinear(m1v, nCells_);
1494 gezelter 1568
1495 gezelter 1587 for (vector<Vector3i>::iterator os = cellOffsets_.begin();
1496     os != cellOffsets_.end(); ++os) {
1497    
1498     Vector3i m2v = m1v + (*os);
1499 gezelter 1612
1500    
1501 gezelter 1587 if (m2v.x() >= nCells_.x()) {
1502     m2v.x() = 0;
1503     } else if (m2v.x() < 0) {
1504     m2v.x() = nCells_.x() - 1;
1505     }
1506    
1507     if (m2v.y() >= nCells_.y()) {
1508     m2v.y() = 0;
1509     } else if (m2v.y() < 0) {
1510     m2v.y() = nCells_.y() - 1;
1511     }
1512    
1513     if (m2v.z() >= nCells_.z()) {
1514     m2v.z() = 0;
1515     } else if (m2v.z() < 0) {
1516     m2v.z() = nCells_.z() - 1;
1517     }
1518 gezelter 1612
1519 gezelter 1587 int m2 = Vlinear (m2v, nCells_);
1520    
1521 gezelter 1567 #ifdef IS_MPI
1522 gezelter 1587 for (vector<int>::iterator j1 = cellListRow_[m1].begin();
1523     j1 != cellListRow_[m1].end(); ++j1) {
1524     for (vector<int>::iterator j2 = cellListCol_[m2].begin();
1525     j2 != cellListCol_[m2].end(); ++j2) {
1526    
1527 gezelter 1612 // In parallel, we need to visit *all* pairs of row
1528     // & column indicies and will divide labor in the
1529     // force evaluation later.
1530 gezelter 1593 dr = cgColData.position[(*j2)] - cgRowData.position[(*j1)];
1531 gezelter 1879 if (usePeriodicBoundaryConditions_) {
1532     snap_->wrapVector(dr);
1533     }
1534 gezelter 1896 getGroupCutoffs( (*j1), (*j2), rcut, rcutsq, rlistsq );
1535     if (dr.lengthSquare() < rlistsq) {
1536 gezelter 1593 neighborList.push_back(make_pair((*j1), (*j2)));
1537     }
1538 gezelter 1562 }
1539     }
1540 gezelter 1567 #else
1541 gezelter 1587 for (vector<int>::iterator j1 = cellList_[m1].begin();
1542     j1 != cellList_[m1].end(); ++j1) {
1543     for (vector<int>::iterator j2 = cellList_[m2].begin();
1544     j2 != cellList_[m2].end(); ++j2) {
1545 gezelter 1616
1546 gezelter 1587 // Always do this if we're in different cells or if
1547 gezelter 1616 // we're in the same cell and the global index of
1548     // the j2 cutoff group is greater than or equal to
1549     // the j1 cutoff group. Note that Rappaport's code
1550     // has a "less than" conditional here, but that
1551     // deals with atom-by-atom computation. OpenMD
1552     // allows atoms within a single cutoff group to
1553     // interact with each other.
1554    
1555     if (m2 != m1 || (*j2) >= (*j1) ) {
1556    
1557 gezelter 1587 dr = snap_->cgData.position[(*j2)] - snap_->cgData.position[(*j1)];
1558 gezelter 1879 if (usePeriodicBoundaryConditions_) {
1559     snap_->wrapVector(dr);
1560     }
1561 gezelter 1896 getGroupCutoffs( (*j1), (*j2), rcut, rcutsq, rlistsq );
1562     if (dr.lengthSquare() < rlistsq) {
1563 gezelter 1587 neighborList.push_back(make_pair((*j1), (*j2)));
1564     }
1565 gezelter 1567 }
1566     }
1567     }
1568 gezelter 1587 #endif
1569 gezelter 1567 }
1570 gezelter 1562 }
1571     }
1572     }
1573 gezelter 1587 } else {
1574     // branch to do all cutoff group pairs
1575     #ifdef IS_MPI
1576     for (int j1 = 0; j1 < nGroupsInRow_; j1++) {
1577 gezelter 1616 for (int j2 = 0; j2 < nGroupsInCol_; j2++) {
1578 gezelter 1587 dr = cgColData.position[j2] - cgRowData.position[j1];
1579 gezelter 1879 if (usePeriodicBoundaryConditions_) {
1580     snap_->wrapVector(dr);
1581     }
1582 gezelter 1896 getGroupCutoffs( j1, j2, rcut, rcutsq, rlistsq);
1583     if (dr.lengthSquare() < rlistsq) {
1584 gezelter 1587 neighborList.push_back(make_pair(j1, j2));
1585     }
1586     }
1587 gezelter 1616 }
1588 gezelter 1587 #else
1589 gezelter 1616 // include all groups here.
1590     for (int j1 = 0; j1 < nGroups_; j1++) {
1591     // include self group interactions j2 == j1
1592     for (int j2 = j1; j2 < nGroups_; j2++) {
1593 gezelter 1587 dr = snap_->cgData.position[j2] - snap_->cgData.position[j1];
1594 gezelter 1879 if (usePeriodicBoundaryConditions_) {
1595     snap_->wrapVector(dr);
1596     }
1597 gezelter 1896 getGroupCutoffs( j1, j2, rcut, rcutsq, rlistsq );
1598     if (dr.lengthSquare() < rlistsq) {
1599 gezelter 1587 neighborList.push_back(make_pair(j1, j2));
1600     }
1601 gezelter 1616 }
1602     }
1603 gezelter 1587 #endif
1604 gezelter 1562 }
1605 gezelter 1587
1606 gezelter 1568 // save the local cutoff group positions for the check that is
1607     // done on each loop:
1608     saved_CG_positions_.clear();
1609     for (int i = 0; i < nGroups_; i++)
1610     saved_CG_positions_.push_back(snap_->cgData.position[i]);
1611 gezelter 1562 }
1612 gezelter 1539 } //end namespace OpenMD