ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/parallel/ForceMatrixDecomposition.cpp
Revision: 2061
Committed: Tue Mar 3 16:24:44 2015 UTC (10 years, 1 month ago) by gezelter
File size: 54574 byte(s)
Log Message:
Fixed a warning on PGI compilation

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 2061
440     #ifdef IS_MPI
441    
442 gezelter 1551 snap_ = sman_->getCurrentSnapshot();
443     storageLayout_ = sman_->getStorageLayout();
444 gezelter 2061
445 gezelter 2057 bool needsCG = true;
446     if(info_->getNCutoffGroups() != info_->getNAtoms())
447     needsCG = false;
448 gezelter 2061
449 gezelter 1539 // gather up the atomic positions
450 gezelter 1593 AtomPlanVectorRow->gather(snap_->atomData.position,
451 gezelter 1551 atomRowData.position);
452 gezelter 1593 AtomPlanVectorColumn->gather(snap_->atomData.position,
453 gezelter 1551 atomColData.position);
454 gezelter 1539
455     // gather up the cutoff group positions
456 gezelter 1593
457 gezelter 2057 if (needsCG) {
458     cgPlanVectorRow->gather(snap_->cgData.position,
459     cgRowData.position);
460    
461     cgPlanVectorColumn->gather(snap_->cgData.position,
462     cgColData.position);
463     }
464 gezelter 1593
465    
466 gezelter 1723 if (needVelocities_) {
467     // gather up the atomic velocities
468     AtomPlanVectorColumn->gather(snap_->atomData.velocity,
469     atomColData.velocity);
470 gezelter 2057
471     if (needsCG) {
472     cgPlanVectorColumn->gather(snap_->cgData.velocity,
473     cgColData.velocity);
474     }
475 gezelter 1723 }
476    
477 gezelter 1539
478     // if needed, gather the atomic rotation matrices
479 gezelter 1551 if (storageLayout_ & DataStorage::dslAmat) {
480 gezelter 1593 AtomPlanMatrixRow->gather(snap_->atomData.aMat,
481 gezelter 1551 atomRowData.aMat);
482 gezelter 1593 AtomPlanMatrixColumn->gather(snap_->atomData.aMat,
483 gezelter 1551 atomColData.aMat);
484 gezelter 1539 }
485 gezelter 1879
486     // if needed, gather the atomic eletrostatic information
487     if (storageLayout_ & DataStorage::dslDipole) {
488     AtomPlanVectorRow->gather(snap_->atomData.dipole,
489     atomRowData.dipole);
490     AtomPlanVectorColumn->gather(snap_->atomData.dipole,
491     atomColData.dipole);
492 gezelter 1539 }
493 gezelter 1590
494 gezelter 1879 if (storageLayout_ & DataStorage::dslQuadrupole) {
495     AtomPlanMatrixRow->gather(snap_->atomData.quadrupole,
496     atomRowData.quadrupole);
497     AtomPlanMatrixColumn->gather(snap_->atomData.quadrupole,
498     atomColData.quadrupole);
499     }
500    
501 gezelter 1713 // if needed, gather the atomic fluctuating charge values
502     if (storageLayout_ & DataStorage::dslFlucQPosition) {
503     AtomPlanRealRow->gather(snap_->atomData.flucQPos,
504     atomRowData.flucQPos);
505     AtomPlanRealColumn->gather(snap_->atomData.flucQPos,
506     atomColData.flucQPos);
507     }
508    
509 gezelter 1539 #endif
510     }
511    
512 gezelter 1575 /* collects information obtained during the pre-pair loop onto local
513     * data structures.
514     */
515 gezelter 1549 void ForceMatrixDecomposition::collectIntermediateData() {
516 gezelter 1551 snap_ = sman_->getCurrentSnapshot();
517     storageLayout_ = sman_->getStorageLayout();
518 gezelter 1539 #ifdef IS_MPI
519    
520 gezelter 1551 if (storageLayout_ & DataStorage::dslDensity) {
521    
522 gezelter 1593 AtomPlanRealRow->scatter(atomRowData.density,
523 gezelter 1551 snap_->atomData.density);
524    
525     int n = snap_->atomData.density.size();
526 gezelter 1575 vector<RealType> rho_tmp(n, 0.0);
527 gezelter 1593 AtomPlanRealColumn->scatter(atomColData.density, rho_tmp);
528 gezelter 1539 for (int i = 0; i < n; i++)
529 gezelter 1551 snap_->atomData.density[i] += rho_tmp[i];
530 gezelter 1539 }
531 gezelter 1713
532 gezelter 1879 // this isn't necessary if we don't have polarizable atoms, but
533     // we'll leave it here for now.
534 gezelter 1713 if (storageLayout_ & DataStorage::dslElectricField) {
535    
536     AtomPlanVectorRow->scatter(atomRowData.electricField,
537     snap_->atomData.electricField);
538    
539     int n = snap_->atomData.electricField.size();
540     vector<Vector3d> field_tmp(n, V3Zero);
541 gezelter 1879 AtomPlanVectorColumn->scatter(atomColData.electricField,
542     field_tmp);
543 gezelter 1713 for (int i = 0; i < n; i++)
544     snap_->atomData.electricField[i] += field_tmp[i];
545     }
546 chuckv 1538 #endif
547 gezelter 1539 }
548 gezelter 1575
549     /*
550     * redistributes information obtained during the pre-pair loop out to
551     * row and column-indexed data structures
552     */
553 gezelter 1549 void ForceMatrixDecomposition::distributeIntermediateData() {
554 gezelter 1551 snap_ = sman_->getCurrentSnapshot();
555     storageLayout_ = sman_->getStorageLayout();
556 chuckv 1538 #ifdef IS_MPI
557 gezelter 1551 if (storageLayout_ & DataStorage::dslFunctional) {
558 gezelter 1593 AtomPlanRealRow->gather(snap_->atomData.functional,
559 gezelter 1551 atomRowData.functional);
560 gezelter 1593 AtomPlanRealColumn->gather(snap_->atomData.functional,
561 gezelter 1551 atomColData.functional);
562 gezelter 1539 }
563    
564 gezelter 1551 if (storageLayout_ & DataStorage::dslFunctionalDerivative) {
565 gezelter 1593 AtomPlanRealRow->gather(snap_->atomData.functionalDerivative,
566 gezelter 1551 atomRowData.functionalDerivative);
567 gezelter 1593 AtomPlanRealColumn->gather(snap_->atomData.functionalDerivative,
568 gezelter 1551 atomColData.functionalDerivative);
569 gezelter 1539 }
570 chuckv 1538 #endif
571     }
572 gezelter 1539
573    
574 gezelter 1549 void ForceMatrixDecomposition::collectData() {
575 gezelter 1551 snap_ = sman_->getCurrentSnapshot();
576     storageLayout_ = sman_->getStorageLayout();
577     #ifdef IS_MPI
578     int n = snap_->atomData.force.size();
579 gezelter 1544 vector<Vector3d> frc_tmp(n, V3Zero);
580 gezelter 1541
581 gezelter 1593 AtomPlanVectorRow->scatter(atomRowData.force, frc_tmp);
582 gezelter 1541 for (int i = 0; i < n; i++) {
583 gezelter 1551 snap_->atomData.force[i] += frc_tmp[i];
584 gezelter 1541 frc_tmp[i] = 0.0;
585     }
586 gezelter 1540
587 gezelter 1593 AtomPlanVectorColumn->scatter(atomColData.force, frc_tmp);
588     for (int i = 0; i < n; i++) {
589 gezelter 1551 snap_->atomData.force[i] += frc_tmp[i];
590 gezelter 1593 }
591 gezelter 1590
592 gezelter 1551 if (storageLayout_ & DataStorage::dslTorque) {
593 gezelter 1541
594 gezelter 1587 int nt = snap_->atomData.torque.size();
595 gezelter 1544 vector<Vector3d> trq_tmp(nt, V3Zero);
596 gezelter 1541
597 gezelter 1593 AtomPlanVectorRow->scatter(atomRowData.torque, trq_tmp);
598 gezelter 1587 for (int i = 0; i < nt; i++) {
599 gezelter 1551 snap_->atomData.torque[i] += trq_tmp[i];
600 gezelter 1541 trq_tmp[i] = 0.0;
601     }
602 gezelter 1540
603 gezelter 1593 AtomPlanVectorColumn->scatter(atomColData.torque, trq_tmp);
604 gezelter 1587 for (int i = 0; i < nt; i++)
605 gezelter 1551 snap_->atomData.torque[i] += trq_tmp[i];
606 gezelter 1540 }
607 gezelter 1587
608     if (storageLayout_ & DataStorage::dslSkippedCharge) {
609    
610     int ns = snap_->atomData.skippedCharge.size();
611     vector<RealType> skch_tmp(ns, 0.0);
612    
613 gezelter 1593 AtomPlanRealRow->scatter(atomRowData.skippedCharge, skch_tmp);
614 gezelter 1587 for (int i = 0; i < ns; i++) {
615 gezelter 1590 snap_->atomData.skippedCharge[i] += skch_tmp[i];
616 gezelter 1587 skch_tmp[i] = 0.0;
617     }
618    
619 gezelter 1593 AtomPlanRealColumn->scatter(atomColData.skippedCharge, skch_tmp);
620 gezelter 1613 for (int i = 0; i < ns; i++)
621 gezelter 1587 snap_->atomData.skippedCharge[i] += skch_tmp[i];
622 gezelter 1613
623 gezelter 1587 }
624 gezelter 1540
625 gezelter 1713 if (storageLayout_ & DataStorage::dslFlucQForce) {
626    
627     int nq = snap_->atomData.flucQFrc.size();
628     vector<RealType> fqfrc_tmp(nq, 0.0);
629    
630     AtomPlanRealRow->scatter(atomRowData.flucQFrc, fqfrc_tmp);
631     for (int i = 0; i < nq; i++) {
632     snap_->atomData.flucQFrc[i] += fqfrc_tmp[i];
633     fqfrc_tmp[i] = 0.0;
634     }
635    
636     AtomPlanRealColumn->scatter(atomColData.flucQFrc, fqfrc_tmp);
637     for (int i = 0; i < nq; i++)
638     snap_->atomData.flucQFrc[i] += fqfrc_tmp[i];
639    
640     }
641    
642 gezelter 1879 if (storageLayout_ & DataStorage::dslElectricField) {
643    
644     int nef = snap_->atomData.electricField.size();
645     vector<Vector3d> efield_tmp(nef, V3Zero);
646    
647     AtomPlanVectorRow->scatter(atomRowData.electricField, efield_tmp);
648     for (int i = 0; i < nef; i++) {
649     snap_->atomData.electricField[i] += efield_tmp[i];
650     efield_tmp[i] = 0.0;
651     }
652    
653     AtomPlanVectorColumn->scatter(atomColData.electricField, efield_tmp);
654     for (int i = 0; i < nef; i++)
655     snap_->atomData.electricField[i] += efield_tmp[i];
656     }
657    
658 gezelter 1993 if (storageLayout_ & DataStorage::dslSitePotential) {
659 gezelter 1879
660 gezelter 1993 int nsp = snap_->atomData.sitePotential.size();
661     vector<RealType> sp_tmp(nsp, 0.0);
662    
663     AtomPlanRealRow->scatter(atomRowData.sitePotential, sp_tmp);
664     for (int i = 0; i < nsp; i++) {
665     snap_->atomData.sitePotential[i] += sp_tmp[i];
666     sp_tmp[i] = 0.0;
667     }
668    
669     AtomPlanRealColumn->scatter(atomColData.sitePotential, sp_tmp);
670     for (int i = 0; i < nsp; i++)
671     snap_->atomData.sitePotential[i] += sp_tmp[i];
672     }
673    
674 gezelter 1567 nLocal_ = snap_->getNumberOfAtoms();
675 gezelter 1544
676 gezelter 1575 vector<potVec> pot_temp(nLocal_,
677     Vector<RealType, N_INTERACTION_FAMILIES> (0.0));
678 gezelter 1760 vector<potVec> expot_temp(nLocal_,
679     Vector<RealType, N_INTERACTION_FAMILIES> (0.0));
680 gezelter 1575
681     // scatter/gather pot_row into the members of my column
682    
683 gezelter 1593 AtomPlanPotRow->scatter(pot_row, pot_temp);
684 gezelter 1760 AtomPlanPotRow->scatter(expot_row, expot_temp);
685 gezelter 1575
686 gezelter 1760 for (int ii = 0; ii < pot_temp.size(); ii++ )
687 gezelter 1583 pairwisePot += pot_temp[ii];
688 gezelter 1760
689     for (int ii = 0; ii < expot_temp.size(); ii++ )
690     excludedPot += expot_temp[ii];
691 gezelter 1723
692     if (storageLayout_ & DataStorage::dslParticlePot) {
693     // This is the pairwise contribution to the particle pot. The
694     // embedding contribution is added in each of the low level
695     // non-bonded routines. In single processor, this is done in
696     // unpackInteractionData, not in collectData.
697     for (int ii = 0; ii < N_INTERACTION_FAMILIES; ii++) {
698     for (int i = 0; i < nLocal_; i++) {
699     // factor of two is because the total potential terms are divided
700     // by 2 in parallel due to row/ column scatter
701     snap_->atomData.particlePot[i] += 2.0 * pot_temp[i](ii);
702     }
703     }
704     }
705    
706 gezelter 1575 fill(pot_temp.begin(), pot_temp.end(),
707     Vector<RealType, N_INTERACTION_FAMILIES> (0.0));
708 gezelter 1760 fill(expot_temp.begin(), expot_temp.end(),
709     Vector<RealType, N_INTERACTION_FAMILIES> (0.0));
710 gezelter 1575
711 gezelter 1593 AtomPlanPotColumn->scatter(pot_col, pot_temp);
712 gezelter 1760 AtomPlanPotColumn->scatter(expot_col, expot_temp);
713 gezelter 1575
714     for (int ii = 0; ii < pot_temp.size(); ii++ )
715 gezelter 1583 pairwisePot += pot_temp[ii];
716 gezelter 1723
717 gezelter 1760 for (int ii = 0; ii < expot_temp.size(); ii++ )
718     excludedPot += expot_temp[ii];
719    
720 gezelter 1723 if (storageLayout_ & DataStorage::dslParticlePot) {
721     // This is the pairwise contribution to the particle pot. The
722     // embedding contribution is added in each of the low level
723     // non-bonded routines. In single processor, this is done in
724     // unpackInteractionData, not in collectData.
725     for (int ii = 0; ii < N_INTERACTION_FAMILIES; ii++) {
726     for (int i = 0; i < nLocal_; i++) {
727     // factor of two is because the total potential terms are divided
728     // by 2 in parallel due to row/ column scatter
729     snap_->atomData.particlePot[i] += 2.0 * pot_temp[i](ii);
730     }
731     }
732     }
733 gezelter 1601
734 gezelter 1723 if (storageLayout_ & DataStorage::dslParticlePot) {
735     int npp = snap_->atomData.particlePot.size();
736     vector<RealType> ppot_temp(npp, 0.0);
737    
738     // This is the direct or embedding contribution to the particle
739     // pot.
740    
741     AtomPlanRealRow->scatter(atomRowData.particlePot, ppot_temp);
742     for (int i = 0; i < npp; i++) {
743     snap_->atomData.particlePot[i] += ppot_temp[i];
744     }
745    
746     fill(ppot_temp.begin(), ppot_temp.end(), 0.0);
747    
748     AtomPlanRealColumn->scatter(atomColData.particlePot, ppot_temp);
749     for (int i = 0; i < npp; i++) {
750     snap_->atomData.particlePot[i] += ppot_temp[i];
751     }
752     }
753    
754 gezelter 1601 for (int ii = 0; ii < N_INTERACTION_FAMILIES; ii++) {
755     RealType ploc1 = pairwisePot[ii];
756     RealType ploc2 = 0.0;
757 gezelter 1969 MPI_Allreduce(&ploc1, &ploc2, 1, MPI_REALTYPE, MPI_SUM, MPI_COMM_WORLD);
758 gezelter 1601 pairwisePot[ii] = ploc2;
759     }
760    
761 gezelter 1760 for (int ii = 0; ii < N_INTERACTION_FAMILIES; ii++) {
762     RealType ploc1 = excludedPot[ii];
763     RealType ploc2 = 0.0;
764 gezelter 1969 MPI_Allreduce(&ploc1, &ploc2, 1, MPI_REALTYPE, MPI_SUM, MPI_COMM_WORLD);
765 gezelter 1760 excludedPot[ii] = ploc2;
766     }
767    
768 gezelter 1723 // Here be dragons.
769 gezelter 1969 MPI_Comm col = colComm.getComm();
770 gezelter 1613
771 gezelter 1969 MPI_Allreduce(MPI_IN_PLACE,
772 gezelter 1723 &snap_->frameData.conductiveHeatFlux[0], 3,
773 gezelter 1969 MPI_REALTYPE, MPI_SUM, col);
774 gezelter 1723
775    
776 gezelter 1539 #endif
777 gezelter 1583
778 chuckv 1538 }
779 gezelter 1551
780 gezelter 1756 /**
781     * Collects information obtained during the post-pair (and embedding
782     * functional) loops onto local data structures.
783     */
784     void ForceMatrixDecomposition::collectSelfData() {
785     snap_ = sman_->getCurrentSnapshot();
786     storageLayout_ = sman_->getStorageLayout();
787    
788     #ifdef IS_MPI
789     for (int ii = 0; ii < N_INTERACTION_FAMILIES; ii++) {
790     RealType ploc1 = embeddingPot[ii];
791     RealType ploc2 = 0.0;
792 gezelter 1969 MPI_Allreduce(&ploc1, &ploc2, 1, MPI_REALTYPE, MPI_SUM, MPI_COMM_WORLD);
793 gezelter 1756 embeddingPot[ii] = ploc2;
794     }
795 gezelter 1761 for (int ii = 0; ii < N_INTERACTION_FAMILIES; ii++) {
796     RealType ploc1 = excludedSelfPot[ii];
797     RealType ploc2 = 0.0;
798 gezelter 1969 MPI_Allreduce(&ploc1, &ploc2, 1, MPI_REALTYPE, MPI_SUM, MPI_COMM_WORLD);
799 gezelter 1761 excludedSelfPot[ii] = ploc2;
800     }
801 gezelter 1756 #endif
802    
803     }
804    
805    
806    
807 gezelter 1893 int& ForceMatrixDecomposition::getNAtomsInRow() {
808 gezelter 1570 #ifdef IS_MPI
809     return nAtomsInRow_;
810     #else
811     return nLocal_;
812     #endif
813     }
814    
815 gezelter 1569 /**
816     * returns the list of atoms belonging to this group.
817     */
818 gezelter 1893 vector<int>& ForceMatrixDecomposition::getAtomsInGroupRow(int cg1){
819 gezelter 1569 #ifdef IS_MPI
820     return groupListRow_[cg1];
821     #else
822     return groupList_[cg1];
823     #endif
824     }
825    
826 gezelter 1893 vector<int>& ForceMatrixDecomposition::getAtomsInGroupColumn(int cg2){
827 gezelter 1569 #ifdef IS_MPI
828     return groupListCol_[cg2];
829     #else
830     return groupList_[cg2];
831     #endif
832     }
833 chuckv 1538
834 gezelter 1551 Vector3d ForceMatrixDecomposition::getIntergroupVector(int cg1, int cg2){
835     Vector3d d;
836    
837     #ifdef IS_MPI
838     d = cgColData.position[cg2] - cgRowData.position[cg1];
839     #else
840     d = snap_->cgData.position[cg2] - snap_->cgData.position[cg1];
841     #endif
842    
843 gezelter 1879 if (usePeriodicBoundaryConditions_) {
844     snap_->wrapVector(d);
845     }
846 gezelter 1551 return d;
847     }
848    
849 gezelter 1893 Vector3d& ForceMatrixDecomposition::getGroupVelocityColumn(int cg2){
850 gezelter 1723 #ifdef IS_MPI
851     return cgColData.velocity[cg2];
852     #else
853     return snap_->cgData.velocity[cg2];
854     #endif
855     }
856 gezelter 1551
857 gezelter 1893 Vector3d& ForceMatrixDecomposition::getAtomVelocityColumn(int atom2){
858 gezelter 1723 #ifdef IS_MPI
859     return atomColData.velocity[atom2];
860     #else
861     return snap_->atomData.velocity[atom2];
862     #endif
863     }
864    
865    
866 gezelter 1551 Vector3d ForceMatrixDecomposition::getAtomToGroupVectorRow(int atom1, int cg1){
867    
868     Vector3d d;
869    
870     #ifdef IS_MPI
871     d = cgRowData.position[cg1] - atomRowData.position[atom1];
872     #else
873     d = snap_->cgData.position[cg1] - snap_->atomData.position[atom1];
874     #endif
875 gezelter 1879 if (usePeriodicBoundaryConditions_) {
876     snap_->wrapVector(d);
877     }
878 gezelter 1551 return d;
879     }
880    
881     Vector3d ForceMatrixDecomposition::getAtomToGroupVectorColumn(int atom2, int cg2){
882     Vector3d d;
883    
884     #ifdef IS_MPI
885     d = cgColData.position[cg2] - atomColData.position[atom2];
886     #else
887     d = snap_->cgData.position[cg2] - snap_->atomData.position[atom2];
888     #endif
889 gezelter 1879 if (usePeriodicBoundaryConditions_) {
890     snap_->wrapVector(d);
891     }
892 gezelter 1551 return d;
893     }
894 gezelter 1569
895 gezelter 1893 RealType& ForceMatrixDecomposition::getMassFactorRow(int atom1) {
896 gezelter 1569 #ifdef IS_MPI
897     return massFactorsRow[atom1];
898     #else
899 gezelter 1581 return massFactors[atom1];
900 gezelter 1569 #endif
901     }
902    
903 gezelter 1893 RealType& ForceMatrixDecomposition::getMassFactorColumn(int atom2) {
904 gezelter 1569 #ifdef IS_MPI
905     return massFactorsCol[atom2];
906     #else
907 gezelter 1581 return massFactors[atom2];
908 gezelter 1569 #endif
909    
910     }
911 gezelter 1551
912     Vector3d ForceMatrixDecomposition::getInteratomicVector(int atom1, int atom2){
913     Vector3d d;
914    
915     #ifdef IS_MPI
916     d = atomColData.position[atom2] - atomRowData.position[atom1];
917     #else
918     d = snap_->atomData.position[atom2] - snap_->atomData.position[atom1];
919     #endif
920 gezelter 1879 if (usePeriodicBoundaryConditions_) {
921     snap_->wrapVector(d);
922     }
923 gezelter 1551 return d;
924     }
925    
926 gezelter 1893 vector<int>& ForceMatrixDecomposition::getExcludesForAtom(int atom1) {
927 gezelter 1587 return excludesForAtom[atom1];
928 gezelter 1570 }
929    
930     /**
931 gezelter 1587 * We need to exclude some overcounted interactions that result from
932 gezelter 1575 * the parallel decomposition.
933 gezelter 1570 */
934 gezelter 1756 bool ForceMatrixDecomposition::skipAtomPair(int atom1, int atom2, int cg1, int cg2) {
935 gezelter 1796 int unique_id_1, unique_id_2;
936 gezelter 1616
937 gezelter 1570 #ifdef IS_MPI
938     // in MPI, we have to look up the unique IDs for each atom
939     unique_id_1 = AtomRowToGlobal[atom1];
940     unique_id_2 = AtomColToGlobal[atom2];
941 gezelter 1796 // group1 = cgRowToGlobal[cg1];
942     // group2 = cgColToGlobal[cg2];
943 gezelter 1616 #else
944     unique_id_1 = AtomLocalToGlobal[atom1];
945     unique_id_2 = AtomLocalToGlobal[atom2];
946 gezelter 1796 int group1 = cgLocalToGlobal[cg1];
947     int group2 = cgLocalToGlobal[cg2];
948 gezelter 1616 #endif
949 gezelter 1570
950     if (unique_id_1 == unique_id_2) return true;
951 gezelter 1616
952     #ifdef IS_MPI
953 gezelter 1570 // this prevents us from doing the pair on multiple processors
954     if (unique_id_1 < unique_id_2) {
955     if ((unique_id_1 + unique_id_2) % 2 == 0) return true;
956     } else {
957 gezelter 1616 if ((unique_id_1 + unique_id_2) % 2 == 1) return true;
958 gezelter 1570 }
959 gezelter 1756 #endif
960    
961     #ifndef IS_MPI
962     if (group1 == group2) {
963     if (unique_id_1 < unique_id_2) return true;
964     }
965 gezelter 1587 #endif
966 gezelter 1616
967 gezelter 1587 return false;
968     }
969    
970     /**
971     * We need to handle the interactions for atoms who are involved in
972     * the same rigid body as well as some short range interactions
973     * (bonds, bends, torsions) differently from other interactions.
974     * We'll still visit the pairwise routines, but with a flag that
975     * tells those routines to exclude the pair from direct long range
976     * interactions. Some indirect interactions (notably reaction
977     * field) must still be handled for these pairs.
978     */
979     bool ForceMatrixDecomposition::excludeAtomPair(int atom1, int atom2) {
980 gezelter 1613
981     // excludesForAtom was constructed to use row/column indices in the MPI
982     // version, and to use local IDs in the non-MPI version:
983 gezelter 1570
984 gezelter 1587 for (vector<int>::iterator i = excludesForAtom[atom1].begin();
985     i != excludesForAtom[atom1].end(); ++i) {
986 gezelter 1616 if ( (*i) == atom2 ) return true;
987 gezelter 1583 }
988 gezelter 1579
989 gezelter 1583 return false;
990 gezelter 1570 }
991    
992    
993 gezelter 1551 void ForceMatrixDecomposition::addForceToAtomRow(int atom1, Vector3d fg){
994     #ifdef IS_MPI
995     atomRowData.force[atom1] += fg;
996     #else
997     snap_->atomData.force[atom1] += fg;
998     #endif
999     }
1000    
1001     void ForceMatrixDecomposition::addForceToAtomColumn(int atom2, Vector3d fg){
1002     #ifdef IS_MPI
1003     atomColData.force[atom2] += fg;
1004     #else
1005     snap_->atomData.force[atom2] += fg;
1006     #endif
1007     }
1008    
1009     // filling interaction blocks with pointers
1010 gezelter 1582 void ForceMatrixDecomposition::fillInteractionData(InteractionData &idat,
1011 gezelter 2057 int atom1, int atom2,
1012     bool newAtom1) {
1013 gezelter 1587
1014     idat.excluded = excludeAtomPair(atom1, atom2);
1015 gezelter 2057
1016     if (newAtom1) {
1017    
1018 gezelter 1551 #ifdef IS_MPI
1019 gezelter 2057 idat.atid1 = identsRow[atom1];
1020     idat.atid2 = identsCol[atom2];
1021    
1022     if (regionsRow[atom1] >= 0 && regionsCol[atom2] >= 0) {
1023     idat.sameRegion = (regionsRow[atom1] == regionsCol[atom2]);
1024     } else {
1025     idat.sameRegion = false;
1026     }
1027    
1028     if (storageLayout_ & DataStorage::dslAmat) {
1029     idat.A1 = &(atomRowData.aMat[atom1]);
1030     idat.A2 = &(atomColData.aMat[atom2]);
1031     }
1032    
1033     if (storageLayout_ & DataStorage::dslTorque) {
1034     idat.t1 = &(atomRowData.torque[atom1]);
1035     idat.t2 = &(atomColData.torque[atom2]);
1036     }
1037    
1038     if (storageLayout_ & DataStorage::dslDipole) {
1039     idat.dipole1 = &(atomRowData.dipole[atom1]);
1040     idat.dipole2 = &(atomColData.dipole[atom2]);
1041     }
1042    
1043     if (storageLayout_ & DataStorage::dslQuadrupole) {
1044     idat.quadrupole1 = &(atomRowData.quadrupole[atom1]);
1045     idat.quadrupole2 = &(atomColData.quadrupole[atom2]);
1046     }
1047    
1048     if (storageLayout_ & DataStorage::dslDensity) {
1049     idat.rho1 = &(atomRowData.density[atom1]);
1050     idat.rho2 = &(atomColData.density[atom2]);
1051     }
1052    
1053     if (storageLayout_ & DataStorage::dslFunctional) {
1054     idat.frho1 = &(atomRowData.functional[atom1]);
1055     idat.frho2 = &(atomColData.functional[atom2]);
1056     }
1057    
1058     if (storageLayout_ & DataStorage::dslFunctionalDerivative) {
1059     idat.dfrho1 = &(atomRowData.functionalDerivative[atom1]);
1060     idat.dfrho2 = &(atomColData.functionalDerivative[atom2]);
1061     }
1062    
1063     if (storageLayout_ & DataStorage::dslParticlePot) {
1064     idat.particlePot1 = &(atomRowData.particlePot[atom1]);
1065     idat.particlePot2 = &(atomColData.particlePot[atom2]);
1066     }
1067    
1068     if (storageLayout_ & DataStorage::dslSkippedCharge) {
1069     idat.skippedCharge1 = &(atomRowData.skippedCharge[atom1]);
1070     idat.skippedCharge2 = &(atomColData.skippedCharge[atom2]);
1071     }
1072    
1073     if (storageLayout_ & DataStorage::dslFlucQPosition) {
1074     idat.flucQ1 = &(atomRowData.flucQPos[atom1]);
1075     idat.flucQ2 = &(atomColData.flucQPos[atom2]);
1076     }
1077    
1078     #else
1079    
1080     idat.atid1 = idents[atom1];
1081     idat.atid2 = idents[atom2];
1082    
1083     if (regions[atom1] >= 0 && regions[atom2] >= 0) {
1084     idat.sameRegion = (regions[atom1] == regions[atom2]);
1085     } else {
1086     idat.sameRegion = false;
1087     }
1088    
1089     if (storageLayout_ & DataStorage::dslAmat) {
1090     idat.A1 = &(snap_->atomData.aMat[atom1]);
1091     idat.A2 = &(snap_->atomData.aMat[atom2]);
1092     }
1093    
1094     if (storageLayout_ & DataStorage::dslTorque) {
1095     idat.t1 = &(snap_->atomData.torque[atom1]);
1096     idat.t2 = &(snap_->atomData.torque[atom2]);
1097     }
1098    
1099     if (storageLayout_ & DataStorage::dslDipole) {
1100     idat.dipole1 = &(snap_->atomData.dipole[atom1]);
1101     idat.dipole2 = &(snap_->atomData.dipole[atom2]);
1102     }
1103    
1104     if (storageLayout_ & DataStorage::dslQuadrupole) {
1105     idat.quadrupole1 = &(snap_->atomData.quadrupole[atom1]);
1106     idat.quadrupole2 = &(snap_->atomData.quadrupole[atom2]);
1107     }
1108    
1109     if (storageLayout_ & DataStorage::dslDensity) {
1110     idat.rho1 = &(snap_->atomData.density[atom1]);
1111     idat.rho2 = &(snap_->atomData.density[atom2]);
1112     }
1113    
1114     if (storageLayout_ & DataStorage::dslFunctional) {
1115     idat.frho1 = &(snap_->atomData.functional[atom1]);
1116     idat.frho2 = &(snap_->atomData.functional[atom2]);
1117     }
1118    
1119     if (storageLayout_ & DataStorage::dslFunctionalDerivative) {
1120     idat.dfrho1 = &(snap_->atomData.functionalDerivative[atom1]);
1121     idat.dfrho2 = &(snap_->atomData.functionalDerivative[atom2]);
1122     }
1123    
1124     if (storageLayout_ & DataStorage::dslParticlePot) {
1125     idat.particlePot1 = &(snap_->atomData.particlePot[atom1]);
1126     idat.particlePot2 = &(snap_->atomData.particlePot[atom2]);
1127     }
1128    
1129     if (storageLayout_ & DataStorage::dslSkippedCharge) {
1130     idat.skippedCharge1 = &(snap_->atomData.skippedCharge[atom1]);
1131     idat.skippedCharge2 = &(snap_->atomData.skippedCharge[atom2]);
1132     }
1133    
1134     if (storageLayout_ & DataStorage::dslFlucQPosition) {
1135     idat.flucQ1 = &(snap_->atomData.flucQPos[atom1]);
1136     idat.flucQ2 = &(snap_->atomData.flucQPos[atom2]);
1137     }
1138     #endif
1139    
1140     } else {
1141     // atom1 is not new, so don't bother updating properties of that atom:
1142     #ifdef IS_MPI
1143 gezelter 1895 idat.atid2 = identsCol[atom2];
1144 gezelter 1929
1145 gezelter 1931 if (regionsRow[atom1] >= 0 && regionsCol[atom2] >= 0) {
1146 gezelter 1929 idat.sameRegion = (regionsRow[atom1] == regionsCol[atom2]);
1147 gezelter 1931 } else {
1148     idat.sameRegion = false;
1149     }
1150    
1151 gezelter 1551 if (storageLayout_ & DataStorage::dslAmat) {
1152 gezelter 1554 idat.A2 = &(atomColData.aMat[atom2]);
1153 gezelter 1551 }
1154 gezelter 1567
1155 gezelter 1551 if (storageLayout_ & DataStorage::dslTorque) {
1156 gezelter 1554 idat.t2 = &(atomColData.torque[atom2]);
1157 gezelter 1551 }
1158    
1159 gezelter 1879 if (storageLayout_ & DataStorage::dslDipole) {
1160     idat.dipole2 = &(atomColData.dipole[atom2]);
1161     }
1162    
1163     if (storageLayout_ & DataStorage::dslQuadrupole) {
1164     idat.quadrupole2 = &(atomColData.quadrupole[atom2]);
1165     }
1166    
1167 gezelter 1551 if (storageLayout_ & DataStorage::dslDensity) {
1168 gezelter 1554 idat.rho2 = &(atomColData.density[atom2]);
1169 gezelter 1551 }
1170    
1171 gezelter 1575 if (storageLayout_ & DataStorage::dslFunctional) {
1172     idat.frho2 = &(atomColData.functional[atom2]);
1173     }
1174    
1175 gezelter 1551 if (storageLayout_ & DataStorage::dslFunctionalDerivative) {
1176 gezelter 1554 idat.dfrho2 = &(atomColData.functionalDerivative[atom2]);
1177 gezelter 1551 }
1178 gezelter 1570
1179 gezelter 1575 if (storageLayout_ & DataStorage::dslParticlePot) {
1180     idat.particlePot2 = &(atomColData.particlePot[atom2]);
1181     }
1182    
1183 gezelter 1587 if (storageLayout_ & DataStorage::dslSkippedCharge) {
1184     idat.skippedCharge2 = &(atomColData.skippedCharge[atom2]);
1185     }
1186    
1187 gezelter 2057 if (storageLayout_ & DataStorage::dslFlucQPosition) {
1188 gezelter 1721 idat.flucQ2 = &(atomColData.flucQPos[atom2]);
1189     }
1190    
1191 gezelter 2057 #else
1192 gezelter 1895 idat.atid2 = idents[atom2];
1193 gezelter 1571
1194 gezelter 1931 if (regions[atom1] >= 0 && regions[atom2] >= 0) {
1195 gezelter 1929 idat.sameRegion = (regions[atom1] == regions[atom2]);
1196 gezelter 1931 } else {
1197     idat.sameRegion = false;
1198     }
1199 gezelter 1929
1200 gezelter 1562 if (storageLayout_ & DataStorage::dslAmat) {
1201     idat.A2 = &(snap_->atomData.aMat[atom2]);
1202     }
1203    
1204     if (storageLayout_ & DataStorage::dslTorque) {
1205     idat.t2 = &(snap_->atomData.torque[atom2]);
1206     }
1207    
1208 gezelter 1879 if (storageLayout_ & DataStorage::dslDipole) {
1209     idat.dipole2 = &(snap_->atomData.dipole[atom2]);
1210     }
1211    
1212     if (storageLayout_ & DataStorage::dslQuadrupole) {
1213     idat.quadrupole2 = &(snap_->atomData.quadrupole[atom2]);
1214     }
1215    
1216 gezelter 1583 if (storageLayout_ & DataStorage::dslDensity) {
1217 gezelter 1562 idat.rho2 = &(snap_->atomData.density[atom2]);
1218     }
1219    
1220 gezelter 1575 if (storageLayout_ & DataStorage::dslFunctional) {
1221     idat.frho2 = &(snap_->atomData.functional[atom2]);
1222     }
1223    
1224 gezelter 1562 if (storageLayout_ & DataStorage::dslFunctionalDerivative) {
1225     idat.dfrho2 = &(snap_->atomData.functionalDerivative[atom2]);
1226     }
1227 gezelter 1575
1228     if (storageLayout_ & DataStorage::dslParticlePot) {
1229     idat.particlePot2 = &(snap_->atomData.particlePot[atom2]);
1230     }
1231    
1232 gezelter 1587 if (storageLayout_ & DataStorage::dslSkippedCharge) {
1233     idat.skippedCharge2 = &(snap_->atomData.skippedCharge[atom2]);
1234     }
1235 gezelter 1721
1236     if (storageLayout_ & DataStorage::dslFlucQPosition) {
1237     idat.flucQ2 = &(snap_->atomData.flucQPos[atom2]);
1238     }
1239    
1240 gezelter 1551 #endif
1241 gezelter 2057 }
1242 gezelter 1551 }
1243 gezelter 1575
1244 gezelter 2057 void ForceMatrixDecomposition::unpackInteractionData(InteractionData &idat,
1245     int atom1, int atom2) {
1246 gezelter 1575 #ifdef IS_MPI
1247 gezelter 1668 pot_row[atom1] += RealType(0.5) * *(idat.pot);
1248     pot_col[atom2] += RealType(0.5) * *(idat.pot);
1249 gezelter 1760 expot_row[atom1] += RealType(0.5) * *(idat.excludedPot);
1250     expot_col[atom2] += RealType(0.5) * *(idat.excludedPot);
1251 gezelter 1575
1252     atomRowData.force[atom1] += *(idat.f1);
1253     atomColData.force[atom2] -= *(idat.f1);
1254 gezelter 1713
1255 gezelter 1721 if (storageLayout_ & DataStorage::dslFlucQForce) {
1256 jmichalk 1736 atomRowData.flucQFrc[atom1] -= *(idat.dVdFQ1);
1257     atomColData.flucQFrc[atom2] -= *(idat.dVdFQ2);
1258 gezelter 1721 }
1259    
1260     if (storageLayout_ & DataStorage::dslElectricField) {
1261     atomRowData.electricField[atom1] += *(idat.eField1);
1262     atomColData.electricField[atom2] += *(idat.eField2);
1263     }
1264    
1265 gezelter 1993 if (storageLayout_ & DataStorage::dslSitePotential) {
1266     atomRowData.sitePotential[atom1] += *(idat.sPot1);
1267     atomColData.sitePotential[atom2] += *(idat.sPot2);
1268     }
1269    
1270 gezelter 1575 #else
1271 gezelter 1583 pairwisePot += *(idat.pot);
1272 gezelter 1760 excludedPot += *(idat.excludedPot);
1273 gezelter 1583
1274 gezelter 1575 snap_->atomData.force[atom1] += *(idat.f1);
1275     snap_->atomData.force[atom2] -= *(idat.f1);
1276 gezelter 1713
1277     if (idat.doParticlePot) {
1278 gezelter 1723 // This is the pairwise contribution to the particle pot. The
1279     // embedding contribution is added in each of the low level
1280     // non-bonded routines. In parallel, this calculation is done
1281     // in collectData, not in unpackInteractionData.
1282 gezelter 1713 snap_->atomData.particlePot[atom1] += *(idat.vpair) * *(idat.sw);
1283 gezelter 1723 snap_->atomData.particlePot[atom2] += *(idat.vpair) * *(idat.sw);
1284 gezelter 1713 }
1285 gezelter 1721
1286     if (storageLayout_ & DataStorage::dslFlucQForce) {
1287 jmichalk 1736 snap_->atomData.flucQFrc[atom1] -= *(idat.dVdFQ1);
1288 gezelter 1721 snap_->atomData.flucQFrc[atom2] -= *(idat.dVdFQ2);
1289     }
1290    
1291     if (storageLayout_ & DataStorage::dslElectricField) {
1292     snap_->atomData.electricField[atom1] += *(idat.eField1);
1293     snap_->atomData.electricField[atom2] += *(idat.eField2);
1294     }
1295    
1296 gezelter 1993 if (storageLayout_ & DataStorage::dslSitePotential) {
1297     snap_->atomData.sitePotential[atom1] += *(idat.sPot1);
1298     snap_->atomData.sitePotential[atom2] += *(idat.sPot2);
1299     }
1300    
1301 gezelter 1575 #endif
1302 gezelter 1586
1303 gezelter 1575 }
1304    
1305 gezelter 1562 /*
1306     * buildNeighborList
1307     *
1308 gezelter 2057 * Constructs the Verlet neighbor list for a force-matrix
1309     * decomposition. In this case, each processor is responsible for
1310     * row-site interactions with column-sites.
1311     *
1312     * neighborList is returned as a packed array of neighboring
1313     * column-ordered CutoffGroups. The starting position in
1314     * neighborList for each row-ordered CutoffGroup is given by the
1315     * returned vector point.
1316 gezelter 1562 */
1317 gezelter 2057 void ForceMatrixDecomposition::buildNeighborList(vector<int>& neighborList,
1318     vector<int>& point) {
1319     neighborList.clear();
1320     point.clear();
1321     int len = 0;
1322 gezelter 1895
1323 gezelter 1587 bool doAllPairs = false;
1324    
1325 gezelter 1879 Snapshot* snap_ = sman_->getCurrentSnapshot();
1326     Mat3x3d box;
1327     Mat3x3d invBox;
1328    
1329     Vector3d rs, scaled, dr;
1330     Vector3i whichCell;
1331     int cellIndex;
1332    
1333 gezelter 1567 #ifdef IS_MPI
1334 gezelter 1568 cellListRow_.clear();
1335     cellListCol_.clear();
1336 gezelter 2057 point.resize(nGroupsInRow_+1);
1337 gezelter 1567 #else
1338 gezelter 1568 cellList_.clear();
1339 gezelter 2057 point.resize(nGroups_+1);
1340 gezelter 1567 #endif
1341 gezelter 1879
1342     if (!usePeriodicBoundaryConditions_) {
1343     box = snap_->getBoundingBox();
1344     invBox = snap_->getInvBoundingBox();
1345     } else {
1346     box = snap_->getHmat();
1347     invBox = snap_->getInvHmat();
1348     }
1349    
1350 gezelter 2057 Vector3d A = box.getColumn(0);
1351     Vector3d B = box.getColumn(1);
1352     Vector3d C = box.getColumn(2);
1353    
1354     // Required for triclinic cells
1355     Vector3d AxB = cross(A, B);
1356     Vector3d BxC = cross(B, C);
1357     Vector3d CxA = cross(C, A);
1358    
1359     // unit vectors perpendicular to the faces of the triclinic cell:
1360     AxB.normalize();
1361     BxC.normalize();
1362     CxA.normalize();
1363    
1364     // A set of perpendicular lengths in triclinic cells:
1365     RealType Wa = abs(dot(A, BxC));
1366     RealType Wb = abs(dot(B, CxA));
1367     RealType Wc = abs(dot(C, AxB));
1368 gezelter 1879
1369 gezelter 2057 nCells_.x() = int( Wa / rList_ );
1370     nCells_.y() = int( Wb / rList_ );
1371     nCells_.z() = int( Wc / rList_ );
1372 gezelter 1879
1373 gezelter 1587 // handle small boxes where the cell offsets can end up repeating cells
1374     if (nCells_.x() < 3) doAllPairs = true;
1375     if (nCells_.y() < 3) doAllPairs = true;
1376     if (nCells_.z() < 3) doAllPairs = true;
1377 gezelter 1879
1378 gezelter 1579 int nCtot = nCells_.x() * nCells_.y() * nCells_.z();
1379 gezelter 1879
1380 gezelter 1567 #ifdef IS_MPI
1381 gezelter 1579 cellListRow_.resize(nCtot);
1382     cellListCol_.resize(nCtot);
1383     #else
1384     cellList_.resize(nCtot);
1385     #endif
1386 gezelter 1879
1387 gezelter 1587 if (!doAllPairs) {
1388 gezelter 2057
1389 gezelter 1579 #ifdef IS_MPI
1390 gezelter 1879
1391 gezelter 1587 for (int i = 0; i < nGroupsInRow_; i++) {
1392     rs = cgRowData.position[i];
1393    
1394     // scaled positions relative to the box vectors
1395 gezelter 1879 scaled = invBox * rs;
1396 gezelter 1587
1397     // wrap the vector back into the unit box by subtracting integer box
1398     // numbers
1399     for (int j = 0; j < 3; j++) {
1400     scaled[j] -= roundMe(scaled[j]);
1401     scaled[j] += 0.5;
1402 gezelter 1772 // Handle the special case when an object is exactly on the
1403     // boundary (a scaled coordinate of 1.0 is the same as
1404     // scaled coordinate of 0.0)
1405     if (scaled[j] >= 1.0) scaled[j] -= 1.0;
1406 gezelter 1587 }
1407    
1408     // find xyz-indices of cell that cutoffGroup is in.
1409     whichCell.x() = nCells_.x() * scaled.x();
1410     whichCell.y() = nCells_.y() * scaled.y();
1411     whichCell.z() = nCells_.z() * scaled.z();
1412    
1413     // find single index of this cell:
1414     cellIndex = Vlinear(whichCell, nCells_);
1415    
1416     // add this cutoff group to the list of groups in this cell;
1417     cellListRow_[cellIndex].push_back(i);
1418 gezelter 1581 }
1419 gezelter 1587 for (int i = 0; i < nGroupsInCol_; i++) {
1420     rs = cgColData.position[i];
1421    
1422     // scaled positions relative to the box vectors
1423 gezelter 1879 scaled = invBox * rs;
1424 gezelter 1587
1425     // wrap the vector back into the unit box by subtracting integer box
1426     // numbers
1427     for (int j = 0; j < 3; j++) {
1428     scaled[j] -= roundMe(scaled[j]);
1429     scaled[j] += 0.5;
1430 gezelter 1772 // Handle the special case when an object is exactly on the
1431     // boundary (a scaled coordinate of 1.0 is the same as
1432     // scaled coordinate of 0.0)
1433     if (scaled[j] >= 1.0) scaled[j] -= 1.0;
1434 gezelter 1587 }
1435    
1436     // find xyz-indices of cell that cutoffGroup is in.
1437     whichCell.x() = nCells_.x() * scaled.x();
1438     whichCell.y() = nCells_.y() * scaled.y();
1439     whichCell.z() = nCells_.z() * scaled.z();
1440    
1441     // find single index of this cell:
1442     cellIndex = Vlinear(whichCell, nCells_);
1443    
1444     // add this cutoff group to the list of groups in this cell;
1445     cellListCol_[cellIndex].push_back(i);
1446 gezelter 1581 }
1447 gezelter 2057
1448 gezelter 1567 #else
1449 gezelter 1587 for (int i = 0; i < nGroups_; i++) {
1450     rs = snap_->cgData.position[i];
1451    
1452     // scaled positions relative to the box vectors
1453 gezelter 1879 scaled = invBox * rs;
1454 gezelter 1587
1455     // wrap the vector back into the unit box by subtracting integer box
1456     // numbers
1457     for (int j = 0; j < 3; j++) {
1458     scaled[j] -= roundMe(scaled[j]);
1459     scaled[j] += 0.5;
1460 gezelter 1771 // Handle the special case when an object is exactly on the
1461     // boundary (a scaled coordinate of 1.0 is the same as
1462     // scaled coordinate of 0.0)
1463     if (scaled[j] >= 1.0) scaled[j] -= 1.0;
1464 gezelter 1587 }
1465    
1466     // find xyz-indices of cell that cutoffGroup is in.
1467 gezelter 1939 whichCell.x() = int(nCells_.x() * scaled.x());
1468     whichCell.y() = int(nCells_.y() * scaled.y());
1469     whichCell.z() = int(nCells_.z() * scaled.z());
1470 gezelter 1587
1471     // find single index of this cell:
1472 gezelter 1593 cellIndex = Vlinear(whichCell, nCells_);
1473 gezelter 1587
1474     // add this cutoff group to the list of groups in this cell;
1475     cellList_[cellIndex].push_back(i);
1476 gezelter 1581 }
1477 gezelter 1612
1478 gezelter 1567 #endif
1479    
1480 gezelter 2057 #ifdef IS_MPI
1481     for (int j1 = 0; j1 < nGroupsInRow_; j1++) {
1482     rs = cgRowData.position[j1];
1483     #else
1484 gezelter 1612
1485 gezelter 2057 for (int j1 = 0; j1 < nGroups_; j1++) {
1486     rs = snap_->cgData.position[j1];
1487     #endif
1488     point[j1] = len;
1489    
1490     // scaled positions relative to the box vectors
1491     scaled = invBox * rs;
1492    
1493     // wrap the vector back into the unit box by subtracting integer box
1494     // numbers
1495     for (int j = 0; j < 3; j++) {
1496     scaled[j] -= roundMe(scaled[j]);
1497     scaled[j] += 0.5;
1498     // Handle the special case when an object is exactly on the
1499     // boundary (a scaled coordinate of 1.0 is the same as
1500     // scaled coordinate of 0.0)
1501     if (scaled[j] >= 1.0) scaled[j] -= 1.0;
1502     }
1503    
1504     // find xyz-indices of cell that cutoffGroup is in.
1505     whichCell.x() = nCells_.x() * scaled.x();
1506     whichCell.y() = nCells_.y() * scaled.y();
1507     whichCell.z() = nCells_.z() * scaled.z();
1508    
1509     // find single index of this cell:
1510     int m1 = Vlinear(whichCell, nCells_);
1511    
1512     for (vector<Vector3i>::iterator os = cellOffsets_.begin();
1513     os != cellOffsets_.end(); ++os) {
1514 gezelter 1587
1515 gezelter 2057 Vector3i m2v = whichCell + (*os);
1516 gezelter 1612
1517 gezelter 2057 if (m2v.x() >= nCells_.x()) {
1518     m2v.x() = 0;
1519     } else if (m2v.x() < 0) {
1520     m2v.x() = nCells_.x() - 1;
1521     }
1522    
1523     if (m2v.y() >= nCells_.y()) {
1524     m2v.y() = 0;
1525     } else if (m2v.y() < 0) {
1526     m2v.y() = nCells_.y() - 1;
1527     }
1528    
1529     if (m2v.z() >= nCells_.z()) {
1530     m2v.z() = 0;
1531     } else if (m2v.z() < 0) {
1532     m2v.z() = nCells_.z() - 1;
1533     }
1534     int m2 = Vlinear (m2v, nCells_);
1535 gezelter 1567 #ifdef IS_MPI
1536 gezelter 2057 for (vector<int>::iterator j2 = cellListCol_[m2].begin();
1537     j2 != cellListCol_[m2].end(); ++j2) {
1538    
1539     // In parallel, we need to visit *all* pairs of row
1540     // & column indicies and will divide labor in the
1541     // force evaluation later.
1542     dr = cgColData.position[(*j2)] - rs;
1543     if (usePeriodicBoundaryConditions_) {
1544     snap_->wrapVector(dr);
1545     }
1546     if (dr.lengthSquare() < rListSq_) {
1547     neighborList.push_back( (*j2) );
1548     ++len;
1549     }
1550     }
1551 gezelter 1567 #else
1552 gezelter 2057 for (vector<int>::iterator j2 = cellList_[m2].begin();
1553     j2 != cellList_[m2].end(); ++j2) {
1554    
1555     // Always do this if we're in different cells or if
1556     // we're in the same cell and the global index of
1557     // the j2 cutoff group is greater than or equal to
1558     // the j1 cutoff group. Note that Rappaport's code
1559     // has a "less than" conditional here, but that
1560     // deals with atom-by-atom computation. OpenMD
1561     // allows atoms within a single cutoff group to
1562     // interact with each other.
1563    
1564     if ( (*j2) >= j1 ) {
1565    
1566     dr = snap_->cgData.position[(*j2)] - rs;
1567     if (usePeriodicBoundaryConditions_) {
1568     snap_->wrapVector(dr);
1569 gezelter 1567 }
1570 gezelter 2057 if ( dr.lengthSquare() < rListSq_) {
1571     neighborList.push_back( (*j2) );
1572     ++len;
1573     }
1574     }
1575     }
1576 gezelter 1587 #endif
1577 gezelter 1562 }
1578 gezelter 2057 }
1579 gezelter 1587 } else {
1580     // branch to do all cutoff group pairs
1581     #ifdef IS_MPI
1582     for (int j1 = 0; j1 < nGroupsInRow_; j1++) {
1583 gezelter 2057 point[j1] = len;
1584     rs = cgRowData.position[j1];
1585 gezelter 1616 for (int j2 = 0; j2 < nGroupsInCol_; j2++) {
1586 gezelter 2057 dr = cgColData.position[j2] - rs;
1587 gezelter 1879 if (usePeriodicBoundaryConditions_) {
1588     snap_->wrapVector(dr);
1589     }
1590 gezelter 2057 if (dr.lengthSquare() < rListSq_) {
1591     neighborList.push_back( j2 );
1592     ++len;
1593 gezelter 1587 }
1594     }
1595 gezelter 1616 }
1596 gezelter 1587 #else
1597 gezelter 1616 // include all groups here.
1598     for (int j1 = 0; j1 < nGroups_; j1++) {
1599 gezelter 2057 point[j1] = len;
1600     rs = snap_->cgData.position[j1];
1601 gezelter 1616 // include self group interactions j2 == j1
1602     for (int j2 = j1; j2 < nGroups_; j2++) {
1603 gezelter 2057 dr = snap_->cgData.position[j2] - rs;
1604 gezelter 1879 if (usePeriodicBoundaryConditions_) {
1605     snap_->wrapVector(dr);
1606     }
1607 gezelter 2057 if (dr.lengthSquare() < rListSq_) {
1608     neighborList.push_back( j2 );
1609     ++len;
1610 gezelter 1587 }
1611 gezelter 1616 }
1612     }
1613 gezelter 1587 #endif
1614 gezelter 1562 }
1615 gezelter 2057
1616     #ifdef IS_MPI
1617     point[nGroupsInRow_] = len;
1618     #else
1619     point[nGroups_] = len;
1620     #endif
1621    
1622 gezelter 1568 // save the local cutoff group positions for the check that is
1623     // done on each loop:
1624     saved_CG_positions_.clear();
1625 gezelter 2057 saved_CG_positions_.reserve(nGroups_);
1626 gezelter 1568 for (int i = 0; i < nGroups_; i++)
1627     saved_CG_positions_.push_back(snap_->cgData.position[i]);
1628 gezelter 1562 }
1629 gezelter 1539 } //end namespace OpenMD