ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/parallel/ForceMatrixDecomposition.cpp
Revision: 1756
Committed: Mon Jun 18 18:23:20 2012 UTC (12 years, 10 months ago) by gezelter
File size: 50026 byte(s)
Log Message:
bugfixes for self interactions (particularly in parallel)

File Contents

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