ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/parallel/ForceMatrixDecomposition.cpp
Revision: 2057
Committed: Tue Mar 3 15:22:26 2015 UTC (10 years, 1 month ago) by gezelter
File size: 54573 byte(s)
Log Message:
Performance improvements, Removed CutoffPolicy

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