ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/parallel/ForceMatrixDecomposition.cpp
Revision: 1736
Committed: Tue Jun 5 17:51:31 2012 UTC (12 years, 10 months ago) by jmichalk
Original Path: branches/development/src/parallel/ForceMatrixDecomposition.cpp
File size: 49767 byte(s)
Log Message:
Fixed fluctuating charge forces

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