ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/parallel/ForceMatrixDecomposition.cpp
Revision: 1706
Committed: Fri Apr 27 20:44:16 2012 UTC (13 years ago) by gezelter
File size: 42596 byte(s)
Log Message:
fixed an initialization bug in Dump2XYZ

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