ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/parallel/ForceMatrixDecomposition.cpp
Revision: 1761
Committed: Fri Jun 22 20:01:37 2012 UTC (12 years, 10 months ago) by gezelter
File size: 51543 byte(s)
Log Message:
fixes for fluctuating charges

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