ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/parallel/ForceMatrixDecomposition.cpp
Revision: 1767
Committed: Fri Jul 6 22:01:58 2012 UTC (12 years, 9 months ago) by gezelter
Original Path: branches/development/src/parallel/ForceMatrixDecomposition.cpp
File size: 51562 byte(s)
Log Message:
Various fixes required to compile OpenMD with the MS Visual C++ compiler

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 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 1713 if (storageLayout_ & DataStorage::dslFlucQForce) {
563     fill(atomRowData.flucQFrc.begin(), atomRowData.flucQFrc.end(),
564     0.0);
565     fill(atomColData.flucQFrc.begin(), atomColData.flucQFrc.end(),
566     0.0);
567     }
568    
569 gezelter 1590 #endif
570     // even in parallel, we need to zero out the local arrays:
571    
572 gezelter 1575 if (storageLayout_ & DataStorage::dslParticlePot) {
573     fill(snap_->atomData.particlePot.begin(),
574     snap_->atomData.particlePot.end(), 0.0);
575     }
576    
577     if (storageLayout_ & DataStorage::dslDensity) {
578     fill(snap_->atomData.density.begin(),
579     snap_->atomData.density.end(), 0.0);
580     }
581 gezelter 1706
582 gezelter 1575 if (storageLayout_ & DataStorage::dslFunctional) {
583     fill(snap_->atomData.functional.begin(),
584     snap_->atomData.functional.end(), 0.0);
585     }
586 gezelter 1706
587 gezelter 1575 if (storageLayout_ & DataStorage::dslFunctionalDerivative) {
588     fill(snap_->atomData.functionalDerivative.begin(),
589     snap_->atomData.functionalDerivative.end(), 0.0);
590     }
591 gezelter 1706
592 gezelter 1586 if (storageLayout_ & DataStorage::dslSkippedCharge) {
593     fill(snap_->atomData.skippedCharge.begin(),
594     snap_->atomData.skippedCharge.end(), 0.0);
595     }
596 gezelter 1713
597     if (storageLayout_ & DataStorage::dslElectricField) {
598     fill(snap_->atomData.electricField.begin(),
599     snap_->atomData.electricField.end(), V3Zero);
600     }
601 gezelter 1575 }
602    
603    
604 gezelter 1549 void ForceMatrixDecomposition::distributeData() {
605 gezelter 1551 snap_ = sman_->getCurrentSnapshot();
606     storageLayout_ = sman_->getStorageLayout();
607 chuckv 1538 #ifdef IS_MPI
608 gezelter 1540
609 gezelter 1539 // gather up the atomic positions
610 gezelter 1593 AtomPlanVectorRow->gather(snap_->atomData.position,
611 gezelter 1551 atomRowData.position);
612 gezelter 1593 AtomPlanVectorColumn->gather(snap_->atomData.position,
613 gezelter 1551 atomColData.position);
614 gezelter 1539
615     // gather up the cutoff group positions
616 gezelter 1593
617     cgPlanVectorRow->gather(snap_->cgData.position,
618 gezelter 1551 cgRowData.position);
619 gezelter 1593
620     cgPlanVectorColumn->gather(snap_->cgData.position,
621 gezelter 1551 cgColData.position);
622 gezelter 1593
623 gezelter 1723
624    
625     if (needVelocities_) {
626     // gather up the atomic velocities
627     AtomPlanVectorColumn->gather(snap_->atomData.velocity,
628     atomColData.velocity);
629    
630     cgPlanVectorColumn->gather(snap_->cgData.velocity,
631     cgColData.velocity);
632     }
633    
634 gezelter 1539
635     // if needed, gather the atomic rotation matrices
636 gezelter 1551 if (storageLayout_ & DataStorage::dslAmat) {
637 gezelter 1593 AtomPlanMatrixRow->gather(snap_->atomData.aMat,
638 gezelter 1551 atomRowData.aMat);
639 gezelter 1593 AtomPlanMatrixColumn->gather(snap_->atomData.aMat,
640 gezelter 1551 atomColData.aMat);
641 gezelter 1539 }
642    
643     // if needed, gather the atomic eletrostatic frames
644 gezelter 1551 if (storageLayout_ & DataStorage::dslElectroFrame) {
645 gezelter 1593 AtomPlanMatrixRow->gather(snap_->atomData.electroFrame,
646 gezelter 1551 atomRowData.electroFrame);
647 gezelter 1593 AtomPlanMatrixColumn->gather(snap_->atomData.electroFrame,
648 gezelter 1551 atomColData.electroFrame);
649 gezelter 1539 }
650 gezelter 1590
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     if (storageLayout_ & DataStorage::dslElectricField) {
683    
684     AtomPlanVectorRow->scatter(atomRowData.electricField,
685     snap_->atomData.electricField);
686    
687     int n = snap_->atomData.electricField.size();
688     vector<Vector3d> field_tmp(n, V3Zero);
689     AtomPlanVectorColumn->scatter(atomColData.electricField, field_tmp);
690     for (int i = 0; i < n; i++)
691     snap_->atomData.electricField[i] += field_tmp[i];
692     }
693 chuckv 1538 #endif
694 gezelter 1539 }
695 gezelter 1575
696     /*
697     * redistributes information obtained during the pre-pair loop out to
698     * row and column-indexed data structures
699     */
700 gezelter 1549 void ForceMatrixDecomposition::distributeIntermediateData() {
701 gezelter 1551 snap_ = sman_->getCurrentSnapshot();
702     storageLayout_ = sman_->getStorageLayout();
703 chuckv 1538 #ifdef IS_MPI
704 gezelter 1551 if (storageLayout_ & DataStorage::dslFunctional) {
705 gezelter 1593 AtomPlanRealRow->gather(snap_->atomData.functional,
706 gezelter 1551 atomRowData.functional);
707 gezelter 1593 AtomPlanRealColumn->gather(snap_->atomData.functional,
708 gezelter 1551 atomColData.functional);
709 gezelter 1539 }
710    
711 gezelter 1551 if (storageLayout_ & DataStorage::dslFunctionalDerivative) {
712 gezelter 1593 AtomPlanRealRow->gather(snap_->atomData.functionalDerivative,
713 gezelter 1551 atomRowData.functionalDerivative);
714 gezelter 1593 AtomPlanRealColumn->gather(snap_->atomData.functionalDerivative,
715 gezelter 1551 atomColData.functionalDerivative);
716 gezelter 1539 }
717 chuckv 1538 #endif
718     }
719 gezelter 1539
720    
721 gezelter 1549 void ForceMatrixDecomposition::collectData() {
722 gezelter 1551 snap_ = sman_->getCurrentSnapshot();
723     storageLayout_ = sman_->getStorageLayout();
724     #ifdef IS_MPI
725     int n = snap_->atomData.force.size();
726 gezelter 1544 vector<Vector3d> frc_tmp(n, V3Zero);
727 gezelter 1541
728 gezelter 1593 AtomPlanVectorRow->scatter(atomRowData.force, frc_tmp);
729 gezelter 1541 for (int i = 0; i < n; i++) {
730 gezelter 1551 snap_->atomData.force[i] += frc_tmp[i];
731 gezelter 1541 frc_tmp[i] = 0.0;
732     }
733 gezelter 1540
734 gezelter 1593 AtomPlanVectorColumn->scatter(atomColData.force, frc_tmp);
735     for (int i = 0; i < n; i++) {
736 gezelter 1551 snap_->atomData.force[i] += frc_tmp[i];
737 gezelter 1593 }
738 gezelter 1590
739 gezelter 1551 if (storageLayout_ & DataStorage::dslTorque) {
740 gezelter 1541
741 gezelter 1587 int nt = snap_->atomData.torque.size();
742 gezelter 1544 vector<Vector3d> trq_tmp(nt, V3Zero);
743 gezelter 1541
744 gezelter 1593 AtomPlanVectorRow->scatter(atomRowData.torque, trq_tmp);
745 gezelter 1587 for (int i = 0; i < nt; i++) {
746 gezelter 1551 snap_->atomData.torque[i] += trq_tmp[i];
747 gezelter 1541 trq_tmp[i] = 0.0;
748     }
749 gezelter 1540
750 gezelter 1593 AtomPlanVectorColumn->scatter(atomColData.torque, trq_tmp);
751 gezelter 1587 for (int i = 0; i < nt; i++)
752 gezelter 1551 snap_->atomData.torque[i] += trq_tmp[i];
753 gezelter 1540 }
754 gezelter 1587
755     if (storageLayout_ & DataStorage::dslSkippedCharge) {
756    
757     int ns = snap_->atomData.skippedCharge.size();
758     vector<RealType> skch_tmp(ns, 0.0);
759    
760 gezelter 1593 AtomPlanRealRow->scatter(atomRowData.skippedCharge, skch_tmp);
761 gezelter 1587 for (int i = 0; i < ns; i++) {
762 gezelter 1590 snap_->atomData.skippedCharge[i] += skch_tmp[i];
763 gezelter 1587 skch_tmp[i] = 0.0;
764     }
765    
766 gezelter 1593 AtomPlanRealColumn->scatter(atomColData.skippedCharge, skch_tmp);
767 gezelter 1613 for (int i = 0; i < ns; i++)
768 gezelter 1587 snap_->atomData.skippedCharge[i] += skch_tmp[i];
769 gezelter 1613
770 gezelter 1587 }
771 gezelter 1540
772 gezelter 1713 if (storageLayout_ & DataStorage::dslFlucQForce) {
773    
774     int nq = snap_->atomData.flucQFrc.size();
775     vector<RealType> fqfrc_tmp(nq, 0.0);
776    
777     AtomPlanRealRow->scatter(atomRowData.flucQFrc, fqfrc_tmp);
778     for (int i = 0; i < nq; i++) {
779     snap_->atomData.flucQFrc[i] += fqfrc_tmp[i];
780     fqfrc_tmp[i] = 0.0;
781     }
782    
783     AtomPlanRealColumn->scatter(atomColData.flucQFrc, fqfrc_tmp);
784     for (int i = 0; i < nq; i++)
785     snap_->atomData.flucQFrc[i] += fqfrc_tmp[i];
786    
787     }
788    
789 gezelter 1567 nLocal_ = snap_->getNumberOfAtoms();
790 gezelter 1544
791 gezelter 1575 vector<potVec> pot_temp(nLocal_,
792     Vector<RealType, N_INTERACTION_FAMILIES> (0.0));
793 gezelter 1760 vector<potVec> expot_temp(nLocal_,
794     Vector<RealType, N_INTERACTION_FAMILIES> (0.0));
795 gezelter 1575
796     // scatter/gather pot_row into the members of my column
797    
798 gezelter 1593 AtomPlanPotRow->scatter(pot_row, pot_temp);
799 gezelter 1760 AtomPlanPotRow->scatter(expot_row, expot_temp);
800 gezelter 1575
801 gezelter 1760 for (int ii = 0; ii < pot_temp.size(); ii++ )
802 gezelter 1583 pairwisePot += pot_temp[ii];
803 gezelter 1760
804     for (int ii = 0; ii < expot_temp.size(); ii++ )
805     excludedPot += expot_temp[ii];
806 gezelter 1723
807     if (storageLayout_ & DataStorage::dslParticlePot) {
808     // This is the pairwise contribution to the particle pot. The
809     // embedding contribution is added in each of the low level
810     // non-bonded routines. In single processor, this is done in
811     // unpackInteractionData, not in collectData.
812     for (int ii = 0; ii < N_INTERACTION_FAMILIES; ii++) {
813     for (int i = 0; i < nLocal_; i++) {
814     // factor of two is because the total potential terms are divided
815     // by 2 in parallel due to row/ column scatter
816     snap_->atomData.particlePot[i] += 2.0 * pot_temp[i](ii);
817     }
818     }
819     }
820    
821 gezelter 1575 fill(pot_temp.begin(), pot_temp.end(),
822     Vector<RealType, N_INTERACTION_FAMILIES> (0.0));
823 gezelter 1760 fill(expot_temp.begin(), expot_temp.end(),
824     Vector<RealType, N_INTERACTION_FAMILIES> (0.0));
825 gezelter 1575
826 gezelter 1593 AtomPlanPotColumn->scatter(pot_col, pot_temp);
827 gezelter 1760 AtomPlanPotColumn->scatter(expot_col, expot_temp);
828 gezelter 1575
829     for (int ii = 0; ii < pot_temp.size(); ii++ )
830 gezelter 1583 pairwisePot += pot_temp[ii];
831 gezelter 1723
832 gezelter 1760 for (int ii = 0; ii < expot_temp.size(); ii++ )
833     excludedPot += expot_temp[ii];
834    
835 gezelter 1723 if (storageLayout_ & DataStorage::dslParticlePot) {
836     // This is the pairwise contribution to the particle pot. The
837     // embedding contribution is added in each of the low level
838     // non-bonded routines. In single processor, this is done in
839     // unpackInteractionData, not in collectData.
840     for (int ii = 0; ii < N_INTERACTION_FAMILIES; ii++) {
841     for (int i = 0; i < nLocal_; i++) {
842     // factor of two is because the total potential terms are divided
843     // by 2 in parallel due to row/ column scatter
844     snap_->atomData.particlePot[i] += 2.0 * pot_temp[i](ii);
845     }
846     }
847     }
848 gezelter 1601
849 gezelter 1723 if (storageLayout_ & DataStorage::dslParticlePot) {
850     int npp = snap_->atomData.particlePot.size();
851     vector<RealType> ppot_temp(npp, 0.0);
852    
853     // This is the direct or embedding contribution to the particle
854     // pot.
855    
856     AtomPlanRealRow->scatter(atomRowData.particlePot, ppot_temp);
857     for (int i = 0; i < npp; i++) {
858     snap_->atomData.particlePot[i] += ppot_temp[i];
859     }
860    
861     fill(ppot_temp.begin(), ppot_temp.end(), 0.0);
862    
863     AtomPlanRealColumn->scatter(atomColData.particlePot, ppot_temp);
864     for (int i = 0; i < npp; i++) {
865     snap_->atomData.particlePot[i] += ppot_temp[i];
866     }
867     }
868    
869 gezelter 1601 for (int ii = 0; ii < N_INTERACTION_FAMILIES; ii++) {
870     RealType ploc1 = pairwisePot[ii];
871     RealType ploc2 = 0.0;
872     MPI::COMM_WORLD.Allreduce(&ploc1, &ploc2, 1, MPI::REALTYPE, MPI::SUM);
873     pairwisePot[ii] = ploc2;
874     }
875    
876 gezelter 1760 for (int ii = 0; ii < N_INTERACTION_FAMILIES; ii++) {
877     RealType ploc1 = excludedPot[ii];
878     RealType ploc2 = 0.0;
879     MPI::COMM_WORLD.Allreduce(&ploc1, &ploc2, 1, MPI::REALTYPE, MPI::SUM);
880     excludedPot[ii] = ploc2;
881     }
882    
883 gezelter 1723 // Here be dragons.
884     MPI::Intracomm col = colComm.getComm();
885 gezelter 1613
886 gezelter 1723 col.Allreduce(MPI::IN_PLACE,
887     &snap_->frameData.conductiveHeatFlux[0], 3,
888     MPI::REALTYPE, MPI::SUM);
889    
890    
891 gezelter 1539 #endif
892 gezelter 1583
893 chuckv 1538 }
894 gezelter 1551
895 gezelter 1756 /**
896     * Collects information obtained during the post-pair (and embedding
897     * functional) loops onto local data structures.
898     */
899     void ForceMatrixDecomposition::collectSelfData() {
900     snap_ = sman_->getCurrentSnapshot();
901     storageLayout_ = sman_->getStorageLayout();
902    
903     #ifdef IS_MPI
904     for (int ii = 0; ii < N_INTERACTION_FAMILIES; ii++) {
905     RealType ploc1 = embeddingPot[ii];
906     RealType ploc2 = 0.0;
907     MPI::COMM_WORLD.Allreduce(&ploc1, &ploc2, 1, MPI::REALTYPE, MPI::SUM);
908     embeddingPot[ii] = ploc2;
909     }
910 gezelter 1761 for (int ii = 0; ii < N_INTERACTION_FAMILIES; ii++) {
911     RealType ploc1 = excludedSelfPot[ii];
912     RealType ploc2 = 0.0;
913     MPI::COMM_WORLD.Allreduce(&ploc1, &ploc2, 1, MPI::REALTYPE, MPI::SUM);
914     excludedSelfPot[ii] = ploc2;
915     }
916 gezelter 1756 #endif
917    
918     }
919    
920    
921    
922 gezelter 1570 int ForceMatrixDecomposition::getNAtomsInRow() {
923     #ifdef IS_MPI
924     return nAtomsInRow_;
925     #else
926     return nLocal_;
927     #endif
928     }
929    
930 gezelter 1569 /**
931     * returns the list of atoms belonging to this group.
932     */
933     vector<int> ForceMatrixDecomposition::getAtomsInGroupRow(int cg1){
934     #ifdef IS_MPI
935     return groupListRow_[cg1];
936     #else
937     return groupList_[cg1];
938     #endif
939     }
940    
941     vector<int> ForceMatrixDecomposition::getAtomsInGroupColumn(int cg2){
942     #ifdef IS_MPI
943     return groupListCol_[cg2];
944     #else
945     return groupList_[cg2];
946     #endif
947     }
948 chuckv 1538
949 gezelter 1551 Vector3d ForceMatrixDecomposition::getIntergroupVector(int cg1, int cg2){
950     Vector3d d;
951    
952     #ifdef IS_MPI
953     d = cgColData.position[cg2] - cgRowData.position[cg1];
954     #else
955     d = snap_->cgData.position[cg2] - snap_->cgData.position[cg1];
956     #endif
957    
958     snap_->wrapVector(d);
959     return d;
960     }
961    
962 gezelter 1723 Vector3d ForceMatrixDecomposition::getGroupVelocityColumn(int cg2){
963     #ifdef IS_MPI
964     return cgColData.velocity[cg2];
965     #else
966     return snap_->cgData.velocity[cg2];
967     #endif
968     }
969 gezelter 1551
970 gezelter 1723 Vector3d ForceMatrixDecomposition::getAtomVelocityColumn(int atom2){
971     #ifdef IS_MPI
972     return atomColData.velocity[atom2];
973     #else
974     return snap_->atomData.velocity[atom2];
975     #endif
976     }
977    
978    
979 gezelter 1551 Vector3d ForceMatrixDecomposition::getAtomToGroupVectorRow(int atom1, int cg1){
980    
981     Vector3d d;
982    
983     #ifdef IS_MPI
984     d = cgRowData.position[cg1] - atomRowData.position[atom1];
985     #else
986     d = snap_->cgData.position[cg1] - snap_->atomData.position[atom1];
987     #endif
988    
989     snap_->wrapVector(d);
990     return d;
991     }
992    
993     Vector3d ForceMatrixDecomposition::getAtomToGroupVectorColumn(int atom2, int cg2){
994     Vector3d d;
995    
996     #ifdef IS_MPI
997     d = cgColData.position[cg2] - atomColData.position[atom2];
998     #else
999     d = snap_->cgData.position[cg2] - snap_->atomData.position[atom2];
1000     #endif
1001    
1002     snap_->wrapVector(d);
1003     return d;
1004     }
1005 gezelter 1569
1006     RealType ForceMatrixDecomposition::getMassFactorRow(int atom1) {
1007     #ifdef IS_MPI
1008     return massFactorsRow[atom1];
1009     #else
1010 gezelter 1581 return massFactors[atom1];
1011 gezelter 1569 #endif
1012     }
1013    
1014     RealType ForceMatrixDecomposition::getMassFactorColumn(int atom2) {
1015     #ifdef IS_MPI
1016     return massFactorsCol[atom2];
1017     #else
1018 gezelter 1581 return massFactors[atom2];
1019 gezelter 1569 #endif
1020    
1021     }
1022 gezelter 1551
1023     Vector3d ForceMatrixDecomposition::getInteratomicVector(int atom1, int atom2){
1024     Vector3d d;
1025    
1026     #ifdef IS_MPI
1027     d = atomColData.position[atom2] - atomRowData.position[atom1];
1028     #else
1029     d = snap_->atomData.position[atom2] - snap_->atomData.position[atom1];
1030     #endif
1031    
1032     snap_->wrapVector(d);
1033     return d;
1034     }
1035    
1036 gezelter 1587 vector<int> ForceMatrixDecomposition::getExcludesForAtom(int atom1) {
1037     return excludesForAtom[atom1];
1038 gezelter 1570 }
1039    
1040     /**
1041 gezelter 1587 * We need to exclude some overcounted interactions that result from
1042 gezelter 1575 * the parallel decomposition.
1043 gezelter 1570 */
1044 gezelter 1756 bool ForceMatrixDecomposition::skipAtomPair(int atom1, int atom2, int cg1, int cg2) {
1045     int unique_id_1, unique_id_2, group1, group2;
1046 gezelter 1616
1047 gezelter 1570 #ifdef IS_MPI
1048     // in MPI, we have to look up the unique IDs for each atom
1049     unique_id_1 = AtomRowToGlobal[atom1];
1050     unique_id_2 = AtomColToGlobal[atom2];
1051 gezelter 1756 group1 = cgRowToGlobal[cg1];
1052     group2 = cgColToGlobal[cg2];
1053 gezelter 1616 #else
1054     unique_id_1 = AtomLocalToGlobal[atom1];
1055     unique_id_2 = AtomLocalToGlobal[atom2];
1056 gezelter 1756 group1 = cgLocalToGlobal[cg1];
1057     group2 = cgLocalToGlobal[cg2];
1058 gezelter 1616 #endif
1059 gezelter 1570
1060     if (unique_id_1 == unique_id_2) return true;
1061 gezelter 1616
1062     #ifdef IS_MPI
1063 gezelter 1570 // this prevents us from doing the pair on multiple processors
1064     if (unique_id_1 < unique_id_2) {
1065     if ((unique_id_1 + unique_id_2) % 2 == 0) return true;
1066     } else {
1067 gezelter 1616 if ((unique_id_1 + unique_id_2) % 2 == 1) return true;
1068 gezelter 1570 }
1069 gezelter 1756 #endif
1070    
1071     #ifndef IS_MPI
1072     if (group1 == group2) {
1073     if (unique_id_1 < unique_id_2) return true;
1074     }
1075 gezelter 1587 #endif
1076 gezelter 1616
1077 gezelter 1587 return false;
1078     }
1079    
1080     /**
1081     * We need to handle the interactions for atoms who are involved in
1082     * the same rigid body as well as some short range interactions
1083     * (bonds, bends, torsions) differently from other interactions.
1084     * We'll still visit the pairwise routines, but with a flag that
1085     * tells those routines to exclude the pair from direct long range
1086     * interactions. Some indirect interactions (notably reaction
1087     * field) must still be handled for these pairs.
1088     */
1089     bool ForceMatrixDecomposition::excludeAtomPair(int atom1, int atom2) {
1090 gezelter 1613
1091     // excludesForAtom was constructed to use row/column indices in the MPI
1092     // version, and to use local IDs in the non-MPI version:
1093 gezelter 1570
1094 gezelter 1587 for (vector<int>::iterator i = excludesForAtom[atom1].begin();
1095     i != excludesForAtom[atom1].end(); ++i) {
1096 gezelter 1616 if ( (*i) == atom2 ) return true;
1097 gezelter 1583 }
1098 gezelter 1579
1099 gezelter 1583 return false;
1100 gezelter 1570 }
1101    
1102    
1103 gezelter 1551 void ForceMatrixDecomposition::addForceToAtomRow(int atom1, Vector3d fg){
1104     #ifdef IS_MPI
1105     atomRowData.force[atom1] += fg;
1106     #else
1107     snap_->atomData.force[atom1] += fg;
1108     #endif
1109     }
1110    
1111     void ForceMatrixDecomposition::addForceToAtomColumn(int atom2, Vector3d fg){
1112     #ifdef IS_MPI
1113     atomColData.force[atom2] += fg;
1114     #else
1115     snap_->atomData.force[atom2] += fg;
1116     #endif
1117     }
1118    
1119     // filling interaction blocks with pointers
1120 gezelter 1582 void ForceMatrixDecomposition::fillInteractionData(InteractionData &idat,
1121 gezelter 1587 int atom1, int atom2) {
1122    
1123     idat.excluded = excludeAtomPair(atom1, atom2);
1124    
1125 gezelter 1551 #ifdef IS_MPI
1126 gezelter 1591 idat.atypes = make_pair( atypesRow[atom1], atypesCol[atom2]);
1127     //idat.atypes = make_pair( ff_->getAtomType(identsRow[atom1]),
1128     // ff_->getAtomType(identsCol[atom2]) );
1129 gezelter 1571
1130 gezelter 1551 if (storageLayout_ & DataStorage::dslAmat) {
1131 gezelter 1554 idat.A1 = &(atomRowData.aMat[atom1]);
1132     idat.A2 = &(atomColData.aMat[atom2]);
1133 gezelter 1551 }
1134 gezelter 1567
1135 gezelter 1551 if (storageLayout_ & DataStorage::dslElectroFrame) {
1136 gezelter 1554 idat.eFrame1 = &(atomRowData.electroFrame[atom1]);
1137     idat.eFrame2 = &(atomColData.electroFrame[atom2]);
1138 gezelter 1551 }
1139    
1140     if (storageLayout_ & DataStorage::dslTorque) {
1141 gezelter 1554 idat.t1 = &(atomRowData.torque[atom1]);
1142     idat.t2 = &(atomColData.torque[atom2]);
1143 gezelter 1551 }
1144    
1145     if (storageLayout_ & DataStorage::dslDensity) {
1146 gezelter 1554 idat.rho1 = &(atomRowData.density[atom1]);
1147     idat.rho2 = &(atomColData.density[atom2]);
1148 gezelter 1551 }
1149    
1150 gezelter 1575 if (storageLayout_ & DataStorage::dslFunctional) {
1151     idat.frho1 = &(atomRowData.functional[atom1]);
1152     idat.frho2 = &(atomColData.functional[atom2]);
1153     }
1154    
1155 gezelter 1551 if (storageLayout_ & DataStorage::dslFunctionalDerivative) {
1156 gezelter 1554 idat.dfrho1 = &(atomRowData.functionalDerivative[atom1]);
1157     idat.dfrho2 = &(atomColData.functionalDerivative[atom2]);
1158 gezelter 1551 }
1159 gezelter 1570
1160 gezelter 1575 if (storageLayout_ & DataStorage::dslParticlePot) {
1161     idat.particlePot1 = &(atomRowData.particlePot[atom1]);
1162     idat.particlePot2 = &(atomColData.particlePot[atom2]);
1163     }
1164    
1165 gezelter 1587 if (storageLayout_ & DataStorage::dslSkippedCharge) {
1166     idat.skippedCharge1 = &(atomRowData.skippedCharge[atom1]);
1167     idat.skippedCharge2 = &(atomColData.skippedCharge[atom2]);
1168     }
1169    
1170 gezelter 1721 if (storageLayout_ & DataStorage::dslFlucQPosition) {
1171     idat.flucQ1 = &(atomRowData.flucQPos[atom1]);
1172     idat.flucQ2 = &(atomColData.flucQPos[atom2]);
1173     }
1174    
1175 gezelter 1562 #else
1176 gezelter 1688
1177 gezelter 1591 idat.atypes = make_pair( atypesLocal[atom1], atypesLocal[atom2]);
1178 gezelter 1571
1179 gezelter 1562 if (storageLayout_ & DataStorage::dslAmat) {
1180     idat.A1 = &(snap_->atomData.aMat[atom1]);
1181     idat.A2 = &(snap_->atomData.aMat[atom2]);
1182     }
1183    
1184     if (storageLayout_ & DataStorage::dslElectroFrame) {
1185     idat.eFrame1 = &(snap_->atomData.electroFrame[atom1]);
1186     idat.eFrame2 = &(snap_->atomData.electroFrame[atom2]);
1187     }
1188    
1189     if (storageLayout_ & DataStorage::dslTorque) {
1190     idat.t1 = &(snap_->atomData.torque[atom1]);
1191     idat.t2 = &(snap_->atomData.torque[atom2]);
1192     }
1193    
1194 gezelter 1583 if (storageLayout_ & DataStorage::dslDensity) {
1195 gezelter 1562 idat.rho1 = &(snap_->atomData.density[atom1]);
1196     idat.rho2 = &(snap_->atomData.density[atom2]);
1197     }
1198    
1199 gezelter 1575 if (storageLayout_ & DataStorage::dslFunctional) {
1200     idat.frho1 = &(snap_->atomData.functional[atom1]);
1201     idat.frho2 = &(snap_->atomData.functional[atom2]);
1202     }
1203    
1204 gezelter 1562 if (storageLayout_ & DataStorage::dslFunctionalDerivative) {
1205     idat.dfrho1 = &(snap_->atomData.functionalDerivative[atom1]);
1206     idat.dfrho2 = &(snap_->atomData.functionalDerivative[atom2]);
1207     }
1208 gezelter 1575
1209     if (storageLayout_ & DataStorage::dslParticlePot) {
1210     idat.particlePot1 = &(snap_->atomData.particlePot[atom1]);
1211     idat.particlePot2 = &(snap_->atomData.particlePot[atom2]);
1212     }
1213    
1214 gezelter 1587 if (storageLayout_ & DataStorage::dslSkippedCharge) {
1215     idat.skippedCharge1 = &(snap_->atomData.skippedCharge[atom1]);
1216     idat.skippedCharge2 = &(snap_->atomData.skippedCharge[atom2]);
1217     }
1218 gezelter 1721
1219     if (storageLayout_ & DataStorage::dslFlucQPosition) {
1220     idat.flucQ1 = &(snap_->atomData.flucQPos[atom1]);
1221     idat.flucQ2 = &(snap_->atomData.flucQPos[atom2]);
1222     }
1223    
1224 gezelter 1551 #endif
1225     }
1226 gezelter 1567
1227 gezelter 1575
1228 gezelter 1582 void ForceMatrixDecomposition::unpackInteractionData(InteractionData &idat, int atom1, int atom2) {
1229 gezelter 1575 #ifdef IS_MPI
1230 gezelter 1668 pot_row[atom1] += RealType(0.5) * *(idat.pot);
1231     pot_col[atom2] += RealType(0.5) * *(idat.pot);
1232 gezelter 1760 expot_row[atom1] += RealType(0.5) * *(idat.excludedPot);
1233     expot_col[atom2] += RealType(0.5) * *(idat.excludedPot);
1234 gezelter 1575
1235     atomRowData.force[atom1] += *(idat.f1);
1236     atomColData.force[atom2] -= *(idat.f1);
1237 gezelter 1713
1238 gezelter 1721 if (storageLayout_ & DataStorage::dslFlucQForce) {
1239 jmichalk 1736 atomRowData.flucQFrc[atom1] -= *(idat.dVdFQ1);
1240     atomColData.flucQFrc[atom2] -= *(idat.dVdFQ2);
1241 gezelter 1721 }
1242    
1243     if (storageLayout_ & DataStorage::dslElectricField) {
1244     atomRowData.electricField[atom1] += *(idat.eField1);
1245     atomColData.electricField[atom2] += *(idat.eField2);
1246     }
1247    
1248 gezelter 1575 #else
1249 gezelter 1583 pairwisePot += *(idat.pot);
1250 gezelter 1760 excludedPot += *(idat.excludedPot);
1251 gezelter 1583
1252 gezelter 1575 snap_->atomData.force[atom1] += *(idat.f1);
1253     snap_->atomData.force[atom2] -= *(idat.f1);
1254 gezelter 1713
1255     if (idat.doParticlePot) {
1256 gezelter 1723 // This is the pairwise contribution to the particle pot. The
1257     // embedding contribution is added in each of the low level
1258     // non-bonded routines. In parallel, this calculation is done
1259     // in collectData, not in unpackInteractionData.
1260 gezelter 1713 snap_->atomData.particlePot[atom1] += *(idat.vpair) * *(idat.sw);
1261 gezelter 1723 snap_->atomData.particlePot[atom2] += *(idat.vpair) * *(idat.sw);
1262 gezelter 1713 }
1263 gezelter 1721
1264     if (storageLayout_ & DataStorage::dslFlucQForce) {
1265 jmichalk 1736 snap_->atomData.flucQFrc[atom1] -= *(idat.dVdFQ1);
1266 gezelter 1721 snap_->atomData.flucQFrc[atom2] -= *(idat.dVdFQ2);
1267     }
1268    
1269     if (storageLayout_ & DataStorage::dslElectricField) {
1270     snap_->atomData.electricField[atom1] += *(idat.eField1);
1271     snap_->atomData.electricField[atom2] += *(idat.eField2);
1272     }
1273    
1274 gezelter 1575 #endif
1275 gezelter 1586
1276 gezelter 1575 }
1277    
1278 gezelter 1562 /*
1279     * buildNeighborList
1280     *
1281     * first element of pair is row-indexed CutoffGroup
1282     * second element of pair is column-indexed CutoffGroup
1283     */
1284 gezelter 1567 vector<pair<int, int> > ForceMatrixDecomposition::buildNeighborList() {
1285    
1286     vector<pair<int, int> > neighborList;
1287 gezelter 1576 groupCutoffs cuts;
1288 gezelter 1587 bool doAllPairs = false;
1289    
1290 gezelter 1567 #ifdef IS_MPI
1291 gezelter 1568 cellListRow_.clear();
1292     cellListCol_.clear();
1293 gezelter 1567 #else
1294 gezelter 1568 cellList_.clear();
1295 gezelter 1567 #endif
1296 gezelter 1562
1297 gezelter 1576 RealType rList_ = (largestRcut_ + skinThickness_);
1298 gezelter 1567 RealType rl2 = rList_ * rList_;
1299     Snapshot* snap_ = sman_->getCurrentSnapshot();
1300 gezelter 1562 Mat3x3d Hmat = snap_->getHmat();
1301     Vector3d Hx = Hmat.getColumn(0);
1302     Vector3d Hy = Hmat.getColumn(1);
1303     Vector3d Hz = Hmat.getColumn(2);
1304    
1305 gezelter 1568 nCells_.x() = (int) ( Hx.length() )/ rList_;
1306     nCells_.y() = (int) ( Hy.length() )/ rList_;
1307     nCells_.z() = (int) ( Hz.length() )/ rList_;
1308 gezelter 1562
1309 gezelter 1587 // handle small boxes where the cell offsets can end up repeating cells
1310    
1311     if (nCells_.x() < 3) doAllPairs = true;
1312     if (nCells_.y() < 3) doAllPairs = true;
1313     if (nCells_.z() < 3) doAllPairs = true;
1314    
1315 gezelter 1567 Mat3x3d invHmat = snap_->getInvHmat();
1316     Vector3d rs, scaled, dr;
1317     Vector3i whichCell;
1318     int cellIndex;
1319 gezelter 1579 int nCtot = nCells_.x() * nCells_.y() * nCells_.z();
1320 gezelter 1567
1321     #ifdef IS_MPI
1322 gezelter 1579 cellListRow_.resize(nCtot);
1323     cellListCol_.resize(nCtot);
1324     #else
1325     cellList_.resize(nCtot);
1326     #endif
1327 gezelter 1582
1328 gezelter 1587 if (!doAllPairs) {
1329 gezelter 1579 #ifdef IS_MPI
1330 gezelter 1581
1331 gezelter 1587 for (int i = 0; i < nGroupsInRow_; i++) {
1332     rs = cgRowData.position[i];
1333    
1334     // scaled positions relative to the box vectors
1335     scaled = invHmat * rs;
1336    
1337     // wrap the vector back into the unit box by subtracting integer box
1338     // numbers
1339     for (int j = 0; j < 3; j++) {
1340     scaled[j] -= roundMe(scaled[j]);
1341     scaled[j] += 0.5;
1342     }
1343    
1344     // find xyz-indices of cell that cutoffGroup is in.
1345     whichCell.x() = nCells_.x() * scaled.x();
1346     whichCell.y() = nCells_.y() * scaled.y();
1347     whichCell.z() = nCells_.z() * scaled.z();
1348    
1349     // find single index of this cell:
1350     cellIndex = Vlinear(whichCell, nCells_);
1351    
1352     // add this cutoff group to the list of groups in this cell;
1353     cellListRow_[cellIndex].push_back(i);
1354 gezelter 1581 }
1355 gezelter 1587 for (int i = 0; i < nGroupsInCol_; i++) {
1356     rs = cgColData.position[i];
1357    
1358     // scaled positions relative to the box vectors
1359     scaled = invHmat * rs;
1360    
1361     // wrap the vector back into the unit box by subtracting integer box
1362     // numbers
1363     for (int j = 0; j < 3; j++) {
1364     scaled[j] -= roundMe(scaled[j]);
1365     scaled[j] += 0.5;
1366     }
1367    
1368     // find xyz-indices of cell that cutoffGroup is in.
1369     whichCell.x() = nCells_.x() * scaled.x();
1370     whichCell.y() = nCells_.y() * scaled.y();
1371     whichCell.z() = nCells_.z() * scaled.z();
1372    
1373     // find single index of this cell:
1374     cellIndex = Vlinear(whichCell, nCells_);
1375    
1376     // add this cutoff group to the list of groups in this cell;
1377     cellListCol_[cellIndex].push_back(i);
1378 gezelter 1581 }
1379 gezelter 1612
1380 gezelter 1567 #else
1381 gezelter 1587 for (int i = 0; i < nGroups_; i++) {
1382     rs = snap_->cgData.position[i];
1383    
1384     // scaled positions relative to the box vectors
1385     scaled = invHmat * rs;
1386    
1387     // wrap the vector back into the unit box by subtracting integer box
1388     // numbers
1389     for (int j = 0; j < 3; j++) {
1390     scaled[j] -= roundMe(scaled[j]);
1391     scaled[j] += 0.5;
1392     }
1393    
1394     // find xyz-indices of cell that cutoffGroup is in.
1395     whichCell.x() = nCells_.x() * scaled.x();
1396     whichCell.y() = nCells_.y() * scaled.y();
1397     whichCell.z() = nCells_.z() * scaled.z();
1398    
1399     // find single index of this cell:
1400 gezelter 1593 cellIndex = Vlinear(whichCell, nCells_);
1401 gezelter 1587
1402     // add this cutoff group to the list of groups in this cell;
1403     cellList_[cellIndex].push_back(i);
1404 gezelter 1581 }
1405 gezelter 1612
1406 gezelter 1567 #endif
1407    
1408 gezelter 1587 for (int m1z = 0; m1z < nCells_.z(); m1z++) {
1409     for (int m1y = 0; m1y < nCells_.y(); m1y++) {
1410     for (int m1x = 0; m1x < nCells_.x(); m1x++) {
1411     Vector3i m1v(m1x, m1y, m1z);
1412     int m1 = Vlinear(m1v, nCells_);
1413 gezelter 1568
1414 gezelter 1587 for (vector<Vector3i>::iterator os = cellOffsets_.begin();
1415     os != cellOffsets_.end(); ++os) {
1416    
1417     Vector3i m2v = m1v + (*os);
1418 gezelter 1612
1419    
1420 gezelter 1587 if (m2v.x() >= nCells_.x()) {
1421     m2v.x() = 0;
1422     } else if (m2v.x() < 0) {
1423     m2v.x() = nCells_.x() - 1;
1424     }
1425    
1426     if (m2v.y() >= nCells_.y()) {
1427     m2v.y() = 0;
1428     } else if (m2v.y() < 0) {
1429     m2v.y() = nCells_.y() - 1;
1430     }
1431    
1432     if (m2v.z() >= nCells_.z()) {
1433     m2v.z() = 0;
1434     } else if (m2v.z() < 0) {
1435     m2v.z() = nCells_.z() - 1;
1436     }
1437 gezelter 1612
1438 gezelter 1587 int m2 = Vlinear (m2v, nCells_);
1439    
1440 gezelter 1567 #ifdef IS_MPI
1441 gezelter 1587 for (vector<int>::iterator j1 = cellListRow_[m1].begin();
1442     j1 != cellListRow_[m1].end(); ++j1) {
1443     for (vector<int>::iterator j2 = cellListCol_[m2].begin();
1444     j2 != cellListCol_[m2].end(); ++j2) {
1445    
1446 gezelter 1612 // In parallel, we need to visit *all* pairs of row
1447     // & column indicies and will divide labor in the
1448     // force evaluation later.
1449 gezelter 1593 dr = cgColData.position[(*j2)] - cgRowData.position[(*j1)];
1450     snap_->wrapVector(dr);
1451     cuts = getGroupCutoffs( (*j1), (*j2) );
1452     if (dr.lengthSquare() < cuts.third) {
1453     neighborList.push_back(make_pair((*j1), (*j2)));
1454     }
1455 gezelter 1562 }
1456     }
1457 gezelter 1567 #else
1458 gezelter 1587 for (vector<int>::iterator j1 = cellList_[m1].begin();
1459     j1 != cellList_[m1].end(); ++j1) {
1460     for (vector<int>::iterator j2 = cellList_[m2].begin();
1461     j2 != cellList_[m2].end(); ++j2) {
1462 gezelter 1616
1463 gezelter 1587 // Always do this if we're in different cells or if
1464 gezelter 1616 // we're in the same cell and the global index of
1465     // the j2 cutoff group is greater than or equal to
1466     // the j1 cutoff group. Note that Rappaport's code
1467     // has a "less than" conditional here, but that
1468     // deals with atom-by-atom computation. OpenMD
1469     // allows atoms within a single cutoff group to
1470     // interact with each other.
1471    
1472    
1473    
1474     if (m2 != m1 || (*j2) >= (*j1) ) {
1475    
1476 gezelter 1587 dr = snap_->cgData.position[(*j2)] - snap_->cgData.position[(*j1)];
1477     snap_->wrapVector(dr);
1478     cuts = getGroupCutoffs( (*j1), (*j2) );
1479     if (dr.lengthSquare() < cuts.third) {
1480     neighborList.push_back(make_pair((*j1), (*j2)));
1481     }
1482 gezelter 1567 }
1483     }
1484     }
1485 gezelter 1587 #endif
1486 gezelter 1567 }
1487 gezelter 1562 }
1488     }
1489     }
1490 gezelter 1587 } else {
1491     // branch to do all cutoff group pairs
1492     #ifdef IS_MPI
1493     for (int j1 = 0; j1 < nGroupsInRow_; j1++) {
1494 gezelter 1616 for (int j2 = 0; j2 < nGroupsInCol_; j2++) {
1495 gezelter 1587 dr = cgColData.position[j2] - cgRowData.position[j1];
1496     snap_->wrapVector(dr);
1497     cuts = getGroupCutoffs( j1, j2 );
1498     if (dr.lengthSquare() < cuts.third) {
1499     neighborList.push_back(make_pair(j1, j2));
1500     }
1501     }
1502 gezelter 1616 }
1503 gezelter 1587 #else
1504 gezelter 1616 // include all groups here.
1505     for (int j1 = 0; j1 < nGroups_; j1++) {
1506     // include self group interactions j2 == j1
1507     for (int j2 = j1; j2 < nGroups_; j2++) {
1508 gezelter 1587 dr = snap_->cgData.position[j2] - snap_->cgData.position[j1];
1509     snap_->wrapVector(dr);
1510     cuts = getGroupCutoffs( j1, j2 );
1511     if (dr.lengthSquare() < cuts.third) {
1512     neighborList.push_back(make_pair(j1, j2));
1513     }
1514 gezelter 1616 }
1515     }
1516 gezelter 1587 #endif
1517 gezelter 1562 }
1518 gezelter 1587
1519 gezelter 1568 // save the local cutoff group positions for the check that is
1520     // done on each loop:
1521     saved_CG_positions_.clear();
1522     for (int i = 0; i < nGroups_; i++)
1523     saved_CG_positions_.push_back(snap_->cgData.position[i]);
1524 gezelter 1587
1525 gezelter 1567 return neighborList;
1526 gezelter 1562 }
1527 gezelter 1539 } //end namespace OpenMD