ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/parallel/ForceMatrixDecomposition.cpp
Revision: 1850
Committed: Wed Feb 20 15:39:39 2013 UTC (12 years, 2 months ago) by gezelter
File size: 53340 byte(s)
Log Message:
Fixed a widespread typo in the license 

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