ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/parallel/ForceMatrixDecomposition.cpp
Revision: 1993
Committed: Tue Apr 29 17:32:31 2014 UTC (11 years ago) by gezelter
File size: 56536 byte(s)
Log Message:
Added sitePotentials for the Choi et al. potential-frequency maps for nitriles

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