ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/devel_omp/src/parallel/ForceMatrixDecomposition.cpp
Revision: 1593
Committed: Fri Jul 15 21:35:14 2011 UTC (13 years, 9 months ago) by gezelter
Original Path: branches/development/src/parallel/ForceMatrixDecomposition.cpp
File size: 42696 byte(s)
Log Message:
test

File Contents

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