ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/parallel/ForceMatrixDecomposition.cpp
Revision: 1591
Committed: Tue Jul 12 15:25:07 2011 UTC (13 years, 9 months ago) by gezelter
Original Path: branches/development/src/parallel/ForceMatrixDecomposition.cpp
File size: 39900 byte(s)
Log Message:
Efficiency fix

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] Vardeman & Gezelter, in progress (2009).
40 */
41 #include "parallel/ForceMatrixDecomposition.hpp"
42 #include "math/SquareMatrix3.hpp"
43 #include "nonbonded/NonBondedInteraction.hpp"
44 #include "brains/SnapshotManager.hpp"
45 #include "brains/PairList.hpp"
46
47 using namespace std;
48 namespace OpenMD {
49
50 /**
51 * distributeInitialData is essentially a copy of the older fortran
52 * SimulationSetup
53 */
54
55 void ForceMatrixDecomposition::distributeInitialData() {
56 snap_ = sman_->getCurrentSnapshot();
57 storageLayout_ = sman_->getStorageLayout();
58 ff_ = info_->getForceField();
59 nLocal_ = snap_->getNumberOfAtoms();
60
61 nGroups_ = info_->getNLocalCutoffGroups();
62 // gather the information for atomtype IDs (atids):
63 idents = info_->getIdentArray();
64 AtomLocalToGlobal = info_->getGlobalAtomIndices();
65 cgLocalToGlobal = info_->getGlobalGroupIndices();
66 vector<int> globalGroupMembership = info_->getGlobalGroupMembership();
67
68 massFactors = info_->getMassFactors();
69
70 PairList* excludes = info_->getExcludedInteractions();
71 PairList* oneTwo = info_->getOneTwoInteractions();
72 PairList* oneThree = info_->getOneThreeInteractions();
73 PairList* oneFour = info_->getOneFourInteractions();
74
75 #ifdef IS_MPI
76
77 AtomCommIntRow = new Communicator<Row,int>(nLocal_);
78 AtomCommRealRow = new Communicator<Row,RealType>(nLocal_);
79 AtomCommVectorRow = new Communicator<Row,Vector3d>(nLocal_);
80 AtomCommMatrixRow = new Communicator<Row,Mat3x3d>(nLocal_);
81 AtomCommPotRow = new Communicator<Row,potVec>(nLocal_);
82
83 AtomCommIntColumn = new Communicator<Column,int>(nLocal_);
84 AtomCommRealColumn = new Communicator<Column,RealType>(nLocal_);
85 AtomCommVectorColumn = new Communicator<Column,Vector3d>(nLocal_);
86 AtomCommMatrixColumn = new Communicator<Column,Mat3x3d>(nLocal_);
87 AtomCommPotColumn = new Communicator<Column,potVec>(nLocal_);
88
89 cgCommIntRow = new Communicator<Row,int>(nGroups_);
90 cgCommVectorRow = new Communicator<Row,Vector3d>(nGroups_);
91 cgCommIntColumn = new Communicator<Column,int>(nGroups_);
92 cgCommVectorColumn = new Communicator<Column,Vector3d>(nGroups_);
93
94 nAtomsInRow_ = AtomCommIntRow->getSize();
95 nAtomsInCol_ = AtomCommIntColumn->getSize();
96 nGroupsInRow_ = cgCommIntRow->getSize();
97 nGroupsInCol_ = cgCommIntColumn->getSize();
98
99 // Modify the data storage objects with the correct layouts and sizes:
100 atomRowData.resize(nAtomsInRow_);
101 atomRowData.setStorageLayout(storageLayout_);
102 atomColData.resize(nAtomsInCol_);
103 atomColData.setStorageLayout(storageLayout_);
104 cgRowData.resize(nGroupsInRow_);
105 cgRowData.setStorageLayout(DataStorage::dslPosition);
106 cgColData.resize(nGroupsInCol_);
107 cgColData.setStorageLayout(DataStorage::dslPosition);
108
109 identsRow.resize(nAtomsInRow_);
110 identsCol.resize(nAtomsInCol_);
111
112 AtomCommIntRow->gather(idents, identsRow);
113 AtomCommIntColumn->gather(idents, identsCol);
114
115 // allocate memory for the parallel objects
116 atypesRow.resize(nAtomsInRow_);
117 atypesCol.resize(nAtomsInCol_);
118
119 for (int i = 0; i < nAtomsInRow_; i++)
120 atypesRow[i] = ff_->getAtomType(identsRow[i]);
121 for (int i = 0; i < nAtomsInCol_; i++)
122 atypesCol[i] = ff_->getAtomType(identsCol[i]);
123
124 pot_row.resize(nAtomsInRow_);
125 pot_col.resize(nAtomsInCol_);
126
127 AtomRowToGlobal.resize(nAtomsInRow_);
128 AtomColToGlobal.resize(nAtomsInCol_);
129 AtomCommIntRow->gather(AtomLocalToGlobal, AtomRowToGlobal);
130 AtomCommIntColumn->gather(AtomLocalToGlobal, AtomColToGlobal);
131
132 cgRowToGlobal.resize(nGroupsInRow_);
133 cgColToGlobal.resize(nGroupsInCol_);
134 cgCommIntRow->gather(cgLocalToGlobal, cgRowToGlobal);
135 cgCommIntColumn->gather(cgLocalToGlobal, cgColToGlobal);
136
137 massFactorsRow.resize(nAtomsInRow_);
138 massFactorsCol.resize(nAtomsInCol_);
139 AtomCommRealRow->gather(massFactors, massFactorsRow);
140 AtomCommRealColumn->gather(massFactors, massFactorsCol);
141
142 groupListRow_.clear();
143 groupListRow_.resize(nGroupsInRow_);
144 for (int i = 0; i < nGroupsInRow_; i++) {
145 int gid = cgRowToGlobal[i];
146 for (int j = 0; j < nAtomsInRow_; j++) {
147 int aid = AtomRowToGlobal[j];
148 if (globalGroupMembership[aid] == gid)
149 groupListRow_[i].push_back(j);
150 }
151 }
152
153 groupListCol_.clear();
154 groupListCol_.resize(nGroupsInCol_);
155 for (int i = 0; i < nGroupsInCol_; i++) {
156 int gid = cgColToGlobal[i];
157 for (int j = 0; j < nAtomsInCol_; j++) {
158 int aid = AtomColToGlobal[j];
159 if (globalGroupMembership[aid] == gid)
160 groupListCol_[i].push_back(j);
161 }
162 }
163
164 excludesForAtom.clear();
165 excludesForAtom.resize(nAtomsInRow_);
166 toposForAtom.clear();
167 toposForAtom.resize(nAtomsInRow_);
168 topoDist.clear();
169 topoDist.resize(nAtomsInRow_);
170 for (int i = 0; i < nAtomsInRow_; i++) {
171 int iglob = AtomRowToGlobal[i];
172
173 for (int j = 0; j < nAtomsInCol_; j++) {
174 int jglob = AtomColToGlobal[j];
175
176 if (excludes->hasPair(iglob, jglob))
177 excludesForAtom[i].push_back(j);
178
179 if (oneTwo->hasPair(iglob, jglob)) {
180 toposForAtom[i].push_back(j);
181 topoDist[i].push_back(1);
182 } else {
183 if (oneThree->hasPair(iglob, jglob)) {
184 toposForAtom[i].push_back(j);
185 topoDist[i].push_back(2);
186 } else {
187 if (oneFour->hasPair(iglob, jglob)) {
188 toposForAtom[i].push_back(j);
189 topoDist[i].push_back(3);
190 }
191 }
192 }
193 }
194 }
195
196 #endif
197
198 // allocate memory for the parallel objects
199 atypesLocal.resize(nLocal_);
200
201 for (int i = 0; i < nLocal_; i++)
202 atypesLocal[i] = ff_->getAtomType(idents[i]);
203
204 groupList_.clear();
205 groupList_.resize(nGroups_);
206 for (int i = 0; i < nGroups_; i++) {
207 int gid = cgLocalToGlobal[i];
208 for (int j = 0; j < nLocal_; j++) {
209 int aid = AtomLocalToGlobal[j];
210 if (globalGroupMembership[aid] == gid) {
211 groupList_[i].push_back(j);
212 }
213 }
214 }
215
216 excludesForAtom.clear();
217 excludesForAtom.resize(nLocal_);
218 toposForAtom.clear();
219 toposForAtom.resize(nLocal_);
220 topoDist.clear();
221 topoDist.resize(nLocal_);
222
223 for (int i = 0; i < nLocal_; i++) {
224 int iglob = AtomLocalToGlobal[i];
225
226 for (int j = 0; j < nLocal_; j++) {
227 int jglob = AtomLocalToGlobal[j];
228
229 if (excludes->hasPair(iglob, jglob))
230 excludesForAtom[i].push_back(j);
231
232 if (oneTwo->hasPair(iglob, jglob)) {
233 toposForAtom[i].push_back(j);
234 topoDist[i].push_back(1);
235 } else {
236 if (oneThree->hasPair(iglob, jglob)) {
237 toposForAtom[i].push_back(j);
238 topoDist[i].push_back(2);
239 } else {
240 if (oneFour->hasPair(iglob, jglob)) {
241 toposForAtom[i].push_back(j);
242 topoDist[i].push_back(3);
243 }
244 }
245 }
246 }
247 }
248
249 createGtypeCutoffMap();
250
251 }
252
253 void ForceMatrixDecomposition::createGtypeCutoffMap() {
254
255 RealType tol = 1e-6;
256 RealType rc;
257 int atid;
258 set<AtomType*> atypes = info_->getSimulatedAtomTypes();
259 map<int, RealType> atypeCutoff;
260
261 for (set<AtomType*>::iterator at = atypes.begin();
262 at != atypes.end(); ++at){
263 atid = (*at)->getIdent();
264 if (userChoseCutoff_)
265 atypeCutoff[atid] = userCutoff_;
266 else
267 atypeCutoff[atid] = interactionMan_->getSuggestedCutoffRadius(*at);
268 }
269
270 vector<RealType> gTypeCutoffs;
271 // first we do a single loop over the cutoff groups to find the
272 // largest cutoff for any atypes present in this group.
273 #ifdef IS_MPI
274 vector<RealType> groupCutoffRow(nGroupsInRow_, 0.0);
275 groupRowToGtype.resize(nGroupsInRow_);
276 for (int cg1 = 0; cg1 < nGroupsInRow_; cg1++) {
277 vector<int> atomListRow = getAtomsInGroupRow(cg1);
278 for (vector<int>::iterator ia = atomListRow.begin();
279 ia != atomListRow.end(); ++ia) {
280 int atom1 = (*ia);
281 atid = identsRow[atom1];
282 if (atypeCutoff[atid] > groupCutoffRow[cg1]) {
283 groupCutoffRow[cg1] = atypeCutoff[atid];
284 }
285 }
286
287 bool gTypeFound = false;
288 for (int gt = 0; gt < gTypeCutoffs.size(); gt++) {
289 if (abs(groupCutoffRow[cg1] - gTypeCutoffs[gt]) < tol) {
290 groupRowToGtype[cg1] = gt;
291 gTypeFound = true;
292 }
293 }
294 if (!gTypeFound) {
295 gTypeCutoffs.push_back( groupCutoffRow[cg1] );
296 groupRowToGtype[cg1] = gTypeCutoffs.size() - 1;
297 }
298
299 }
300 vector<RealType> groupCutoffCol(nGroupsInCol_, 0.0);
301 groupColToGtype.resize(nGroupsInCol_);
302 for (int cg2 = 0; cg2 < nGroupsInCol_; cg2++) {
303 vector<int> atomListCol = getAtomsInGroupColumn(cg2);
304 for (vector<int>::iterator jb = atomListCol.begin();
305 jb != atomListCol.end(); ++jb) {
306 int atom2 = (*jb);
307 atid = identsCol[atom2];
308 if (atypeCutoff[atid] > groupCutoffCol[cg2]) {
309 groupCutoffCol[cg2] = atypeCutoff[atid];
310 }
311 }
312 bool gTypeFound = false;
313 for (int gt = 0; gt < gTypeCutoffs.size(); gt++) {
314 if (abs(groupCutoffCol[cg2] - gTypeCutoffs[gt]) < tol) {
315 groupColToGtype[cg2] = gt;
316 gTypeFound = true;
317 }
318 }
319 if (!gTypeFound) {
320 gTypeCutoffs.push_back( groupCutoffCol[cg2] );
321 groupColToGtype[cg2] = gTypeCutoffs.size() - 1;
322 }
323 }
324 #else
325
326 vector<RealType> groupCutoff(nGroups_, 0.0);
327 groupToGtype.resize(nGroups_);
328 for (int cg1 = 0; cg1 < nGroups_; cg1++) {
329
330 groupCutoff[cg1] = 0.0;
331 vector<int> atomList = getAtomsInGroupRow(cg1);
332
333 for (vector<int>::iterator ia = atomList.begin();
334 ia != atomList.end(); ++ia) {
335 int atom1 = (*ia);
336 atid = idents[atom1];
337 if (atypeCutoff[atid] > groupCutoff[cg1]) {
338 groupCutoff[cg1] = atypeCutoff[atid];
339 }
340 }
341
342 bool gTypeFound = false;
343 for (int gt = 0; gt < gTypeCutoffs.size(); gt++) {
344 if (abs(groupCutoff[cg1] - gTypeCutoffs[gt]) < tol) {
345 groupToGtype[cg1] = gt;
346 gTypeFound = true;
347 }
348 }
349 if (!gTypeFound) {
350 gTypeCutoffs.push_back( groupCutoff[cg1] );
351 groupToGtype[cg1] = gTypeCutoffs.size() - 1;
352 }
353 }
354 #endif
355
356 // Now we find the maximum group cutoff value present in the simulation
357
358 RealType groupMax = *max_element(gTypeCutoffs.begin(),
359 gTypeCutoffs.end());
360
361 #ifdef IS_MPI
362 MPI::COMM_WORLD.Allreduce(&groupMax, &groupMax, 1, MPI::REALTYPE,
363 MPI::MAX);
364 #endif
365
366 RealType tradRcut = groupMax;
367
368 for (int i = 0; i < gTypeCutoffs.size(); i++) {
369 for (int j = 0; j < gTypeCutoffs.size(); j++) {
370 RealType thisRcut;
371 switch(cutoffPolicy_) {
372 case TRADITIONAL:
373 thisRcut = tradRcut;
374 break;
375 case MIX:
376 thisRcut = 0.5 * (gTypeCutoffs[i] + gTypeCutoffs[j]);
377 break;
378 case MAX:
379 thisRcut = max(gTypeCutoffs[i], gTypeCutoffs[j]);
380 break;
381 default:
382 sprintf(painCave.errMsg,
383 "ForceMatrixDecomposition::createGtypeCutoffMap "
384 "hit an unknown cutoff policy!\n");
385 painCave.severity = OPENMD_ERROR;
386 painCave.isFatal = 1;
387 simError();
388 break;
389 }
390
391 pair<int,int> key = make_pair(i,j);
392 gTypeCutoffMap[key].first = thisRcut;
393
394 if (thisRcut > largestRcut_) largestRcut_ = thisRcut;
395
396 gTypeCutoffMap[key].second = thisRcut*thisRcut;
397
398 gTypeCutoffMap[key].third = pow(thisRcut + skinThickness_, 2);
399
400 // sanity check
401
402 if (userChoseCutoff_) {
403 if (abs(gTypeCutoffMap[key].first - userCutoff_) > 0.0001) {
404 sprintf(painCave.errMsg,
405 "ForceMatrixDecomposition::createGtypeCutoffMap "
406 "user-specified rCut (%lf) does not match computed group Cutoff\n", userCutoff_);
407 painCave.severity = OPENMD_ERROR;
408 painCave.isFatal = 1;
409 simError();
410 }
411 }
412 }
413 }
414 }
415
416
417 groupCutoffs ForceMatrixDecomposition::getGroupCutoffs(int cg1, int cg2) {
418 int i, j;
419 #ifdef IS_MPI
420 i = groupRowToGtype[cg1];
421 j = groupColToGtype[cg2];
422 #else
423 i = groupToGtype[cg1];
424 j = groupToGtype[cg2];
425 #endif
426 return gTypeCutoffMap[make_pair(i,j)];
427 }
428
429 int ForceMatrixDecomposition::getTopologicalDistance(int atom1, int atom2) {
430 for (int j = 0; j < toposForAtom[atom1].size(); j++) {
431 if (toposForAtom[atom1][j] == atom2)
432 return topoDist[atom1][j];
433 }
434 return 0;
435 }
436
437 void ForceMatrixDecomposition::zeroWorkArrays() {
438 pairwisePot = 0.0;
439 embeddingPot = 0.0;
440
441 #ifdef IS_MPI
442 if (storageLayout_ & DataStorage::dslForce) {
443 fill(atomRowData.force.begin(), atomRowData.force.end(), V3Zero);
444 fill(atomColData.force.begin(), atomColData.force.end(), V3Zero);
445 }
446
447 if (storageLayout_ & DataStorage::dslTorque) {
448 fill(atomRowData.torque.begin(), atomRowData.torque.end(), V3Zero);
449 fill(atomColData.torque.begin(), atomColData.torque.end(), V3Zero);
450 }
451
452 fill(pot_row.begin(), pot_row.end(),
453 Vector<RealType, N_INTERACTION_FAMILIES> (0.0));
454
455 fill(pot_col.begin(), pot_col.end(),
456 Vector<RealType, N_INTERACTION_FAMILIES> (0.0));
457
458 if (storageLayout_ & DataStorage::dslParticlePot) {
459 fill(atomRowData.particlePot.begin(), atomRowData.particlePot.end(),
460 0.0);
461 fill(atomColData.particlePot.begin(), atomColData.particlePot.end(),
462 0.0);
463 }
464
465 if (storageLayout_ & DataStorage::dslDensity) {
466 fill(atomRowData.density.begin(), atomRowData.density.end(), 0.0);
467 fill(atomColData.density.begin(), atomColData.density.end(), 0.0);
468 }
469
470 if (storageLayout_ & DataStorage::dslFunctional) {
471 fill(atomRowData.functional.begin(), atomRowData.functional.end(),
472 0.0);
473 fill(atomColData.functional.begin(), atomColData.functional.end(),
474 0.0);
475 }
476
477 if (storageLayout_ & DataStorage::dslFunctionalDerivative) {
478 fill(atomRowData.functionalDerivative.begin(),
479 atomRowData.functionalDerivative.end(), 0.0);
480 fill(atomColData.functionalDerivative.begin(),
481 atomColData.functionalDerivative.end(), 0.0);
482 }
483
484 if (storageLayout_ & DataStorage::dslSkippedCharge) {
485 fill(atomRowData.skippedCharge.begin(),
486 atomRowData.skippedCharge.end(), 0.0);
487 fill(atomColData.skippedCharge.begin(),
488 atomColData.skippedCharge.end(), 0.0);
489 }
490
491 #endif
492 // even in parallel, we need to zero out the local arrays:
493
494 if (storageLayout_ & DataStorage::dslParticlePot) {
495 fill(snap_->atomData.particlePot.begin(),
496 snap_->atomData.particlePot.end(), 0.0);
497 }
498
499 if (storageLayout_ & DataStorage::dslDensity) {
500 fill(snap_->atomData.density.begin(),
501 snap_->atomData.density.end(), 0.0);
502 }
503 if (storageLayout_ & DataStorage::dslFunctional) {
504 fill(snap_->atomData.functional.begin(),
505 snap_->atomData.functional.end(), 0.0);
506 }
507 if (storageLayout_ & DataStorage::dslFunctionalDerivative) {
508 fill(snap_->atomData.functionalDerivative.begin(),
509 snap_->atomData.functionalDerivative.end(), 0.0);
510 }
511 if (storageLayout_ & DataStorage::dslSkippedCharge) {
512 fill(snap_->atomData.skippedCharge.begin(),
513 snap_->atomData.skippedCharge.end(), 0.0);
514 }
515
516 }
517
518
519 void ForceMatrixDecomposition::distributeData() {
520 snap_ = sman_->getCurrentSnapshot();
521 storageLayout_ = sman_->getStorageLayout();
522 #ifdef IS_MPI
523
524 // gather up the atomic positions
525 AtomCommVectorRow->gather(snap_->atomData.position,
526 atomRowData.position);
527 AtomCommVectorColumn->gather(snap_->atomData.position,
528 atomColData.position);
529
530 // gather up the cutoff group positions
531 cgCommVectorRow->gather(snap_->cgData.position,
532 cgRowData.position);
533 cgCommVectorColumn->gather(snap_->cgData.position,
534 cgColData.position);
535
536 // if needed, gather the atomic rotation matrices
537 if (storageLayout_ & DataStorage::dslAmat) {
538 AtomCommMatrixRow->gather(snap_->atomData.aMat,
539 atomRowData.aMat);
540 AtomCommMatrixColumn->gather(snap_->atomData.aMat,
541 atomColData.aMat);
542 }
543
544 // if needed, gather the atomic eletrostatic frames
545 if (storageLayout_ & DataStorage::dslElectroFrame) {
546 AtomCommMatrixRow->gather(snap_->atomData.electroFrame,
547 atomRowData.electroFrame);
548 AtomCommMatrixColumn->gather(snap_->atomData.electroFrame,
549 atomColData.electroFrame);
550 }
551
552 #endif
553 }
554
555 /* collects information obtained during the pre-pair loop onto local
556 * data structures.
557 */
558 void ForceMatrixDecomposition::collectIntermediateData() {
559 snap_ = sman_->getCurrentSnapshot();
560 storageLayout_ = sman_->getStorageLayout();
561 #ifdef IS_MPI
562
563 if (storageLayout_ & DataStorage::dslDensity) {
564
565 AtomCommRealRow->scatter(atomRowData.density,
566 snap_->atomData.density);
567
568 int n = snap_->atomData.density.size();
569 vector<RealType> rho_tmp(n, 0.0);
570 AtomCommRealColumn->scatter(atomColData.density, rho_tmp);
571 for (int i = 0; i < n; i++)
572 snap_->atomData.density[i] += rho_tmp[i];
573 }
574 #endif
575 }
576
577 /*
578 * redistributes information obtained during the pre-pair loop out to
579 * row and column-indexed data structures
580 */
581 void ForceMatrixDecomposition::distributeIntermediateData() {
582 snap_ = sman_->getCurrentSnapshot();
583 storageLayout_ = sman_->getStorageLayout();
584 #ifdef IS_MPI
585 if (storageLayout_ & DataStorage::dslFunctional) {
586 AtomCommRealRow->gather(snap_->atomData.functional,
587 atomRowData.functional);
588 AtomCommRealColumn->gather(snap_->atomData.functional,
589 atomColData.functional);
590 }
591
592 if (storageLayout_ & DataStorage::dslFunctionalDerivative) {
593 AtomCommRealRow->gather(snap_->atomData.functionalDerivative,
594 atomRowData.functionalDerivative);
595 AtomCommRealColumn->gather(snap_->atomData.functionalDerivative,
596 atomColData.functionalDerivative);
597 }
598 #endif
599 }
600
601
602 void ForceMatrixDecomposition::collectData() {
603 snap_ = sman_->getCurrentSnapshot();
604 storageLayout_ = sman_->getStorageLayout();
605 #ifdef IS_MPI
606 int n = snap_->atomData.force.size();
607 vector<Vector3d> frc_tmp(n, V3Zero);
608
609 AtomCommVectorRow->scatter(atomRowData.force, frc_tmp);
610 for (int i = 0; i < n; i++) {
611 snap_->atomData.force[i] += frc_tmp[i];
612 frc_tmp[i] = 0.0;
613 }
614
615 AtomCommVectorColumn->scatter(atomColData.force, frc_tmp);
616 for (int i = 0; i < n; i++)
617 snap_->atomData.force[i] += frc_tmp[i];
618
619 if (storageLayout_ & DataStorage::dslTorque) {
620
621 int nt = snap_->atomData.torque.size();
622 vector<Vector3d> trq_tmp(nt, V3Zero);
623
624 AtomCommVectorRow->scatter(atomRowData.torque, trq_tmp);
625 for (int i = 0; i < nt; i++) {
626 snap_->atomData.torque[i] += trq_tmp[i];
627 trq_tmp[i] = 0.0;
628 }
629
630 AtomCommVectorColumn->scatter(atomColData.torque, trq_tmp);
631 for (int i = 0; i < nt; i++)
632 snap_->atomData.torque[i] += trq_tmp[i];
633 }
634
635 if (storageLayout_ & DataStorage::dslSkippedCharge) {
636
637 int ns = snap_->atomData.skippedCharge.size();
638 vector<RealType> skch_tmp(ns, 0.0);
639
640 AtomCommRealRow->scatter(atomRowData.skippedCharge, skch_tmp);
641 for (int i = 0; i < ns; i++) {
642 snap_->atomData.skippedCharge[i] += skch_tmp[i];
643 skch_tmp[i] = 0.0;
644 }
645
646 AtomCommRealColumn->scatter(atomColData.skippedCharge, skch_tmp);
647 for (int i = 0; i < ns; i++)
648 snap_->atomData.skippedCharge[i] += skch_tmp[i];
649 }
650
651 nLocal_ = snap_->getNumberOfAtoms();
652
653 vector<potVec> pot_temp(nLocal_,
654 Vector<RealType, N_INTERACTION_FAMILIES> (0.0));
655
656 // scatter/gather pot_row into the members of my column
657
658 AtomCommPotRow->scatter(pot_row, pot_temp);
659
660 for (int ii = 0; ii < pot_temp.size(); ii++ )
661 pairwisePot += pot_temp[ii];
662
663 fill(pot_temp.begin(), pot_temp.end(),
664 Vector<RealType, N_INTERACTION_FAMILIES> (0.0));
665
666 AtomCommPotColumn->scatter(pot_col, pot_temp);
667
668 for (int ii = 0; ii < pot_temp.size(); ii++ )
669 pairwisePot += pot_temp[ii];
670 #endif
671
672 }
673
674 int ForceMatrixDecomposition::getNAtomsInRow() {
675 #ifdef IS_MPI
676 return nAtomsInRow_;
677 #else
678 return nLocal_;
679 #endif
680 }
681
682 /**
683 * returns the list of atoms belonging to this group.
684 */
685 vector<int> ForceMatrixDecomposition::getAtomsInGroupRow(int cg1){
686 #ifdef IS_MPI
687 return groupListRow_[cg1];
688 #else
689 return groupList_[cg1];
690 #endif
691 }
692
693 vector<int> ForceMatrixDecomposition::getAtomsInGroupColumn(int cg2){
694 #ifdef IS_MPI
695 return groupListCol_[cg2];
696 #else
697 return groupList_[cg2];
698 #endif
699 }
700
701 Vector3d ForceMatrixDecomposition::getIntergroupVector(int cg1, int cg2){
702 Vector3d d;
703
704 #ifdef IS_MPI
705 d = cgColData.position[cg2] - cgRowData.position[cg1];
706 #else
707 d = snap_->cgData.position[cg2] - snap_->cgData.position[cg1];
708 #endif
709
710 snap_->wrapVector(d);
711 return d;
712 }
713
714
715 Vector3d ForceMatrixDecomposition::getAtomToGroupVectorRow(int atom1, int cg1){
716
717 Vector3d d;
718
719 #ifdef IS_MPI
720 d = cgRowData.position[cg1] - atomRowData.position[atom1];
721 #else
722 d = snap_->cgData.position[cg1] - snap_->atomData.position[atom1];
723 #endif
724
725 snap_->wrapVector(d);
726 return d;
727 }
728
729 Vector3d ForceMatrixDecomposition::getAtomToGroupVectorColumn(int atom2, int cg2){
730 Vector3d d;
731
732 #ifdef IS_MPI
733 d = cgColData.position[cg2] - atomColData.position[atom2];
734 #else
735 d = snap_->cgData.position[cg2] - snap_->atomData.position[atom2];
736 #endif
737
738 snap_->wrapVector(d);
739 return d;
740 }
741
742 RealType ForceMatrixDecomposition::getMassFactorRow(int atom1) {
743 #ifdef IS_MPI
744 return massFactorsRow[atom1];
745 #else
746 return massFactors[atom1];
747 #endif
748 }
749
750 RealType ForceMatrixDecomposition::getMassFactorColumn(int atom2) {
751 #ifdef IS_MPI
752 return massFactorsCol[atom2];
753 #else
754 return massFactors[atom2];
755 #endif
756
757 }
758
759 Vector3d ForceMatrixDecomposition::getInteratomicVector(int atom1, int atom2){
760 Vector3d d;
761
762 #ifdef IS_MPI
763 d = atomColData.position[atom2] - atomRowData.position[atom1];
764 #else
765 d = snap_->atomData.position[atom2] - snap_->atomData.position[atom1];
766 #endif
767
768 snap_->wrapVector(d);
769 return d;
770 }
771
772 vector<int> ForceMatrixDecomposition::getExcludesForAtom(int atom1) {
773 return excludesForAtom[atom1];
774 }
775
776 /**
777 * We need to exclude some overcounted interactions that result from
778 * the parallel decomposition.
779 */
780 bool ForceMatrixDecomposition::skipAtomPair(int atom1, int atom2) {
781 int unique_id_1, unique_id_2;
782
783 #ifdef IS_MPI
784 // in MPI, we have to look up the unique IDs for each atom
785 unique_id_1 = AtomRowToGlobal[atom1];
786 unique_id_2 = AtomColToGlobal[atom2];
787
788 // this situation should only arise in MPI simulations
789 if (unique_id_1 == unique_id_2) return true;
790
791 // this prevents us from doing the pair on multiple processors
792 if (unique_id_1 < unique_id_2) {
793 if ((unique_id_1 + unique_id_2) % 2 == 0) return true;
794 } else {
795 if ((unique_id_1 + unique_id_2) % 2 == 1) return true;
796 }
797 #endif
798 return false;
799 }
800
801 /**
802 * We need to handle the interactions for atoms who are involved in
803 * the same rigid body as well as some short range interactions
804 * (bonds, bends, torsions) differently from other interactions.
805 * We'll still visit the pairwise routines, but with a flag that
806 * tells those routines to exclude the pair from direct long range
807 * interactions. Some indirect interactions (notably reaction
808 * field) must still be handled for these pairs.
809 */
810 bool ForceMatrixDecomposition::excludeAtomPair(int atom1, int atom2) {
811 int unique_id_2;
812
813 #ifdef IS_MPI
814 // in MPI, we have to look up the unique IDs for the row atom.
815 unique_id_2 = AtomColToGlobal[atom2];
816 #else
817 // in the normal loop, the atom numbers are unique
818 unique_id_2 = atom2;
819 #endif
820
821 for (vector<int>::iterator i = excludesForAtom[atom1].begin();
822 i != excludesForAtom[atom1].end(); ++i) {
823 if ( (*i) == unique_id_2 ) return true;
824 }
825
826 return false;
827 }
828
829
830 void ForceMatrixDecomposition::addForceToAtomRow(int atom1, Vector3d fg){
831 #ifdef IS_MPI
832 atomRowData.force[atom1] += fg;
833 #else
834 snap_->atomData.force[atom1] += fg;
835 #endif
836 }
837
838 void ForceMatrixDecomposition::addForceToAtomColumn(int atom2, Vector3d fg){
839 #ifdef IS_MPI
840 atomColData.force[atom2] += fg;
841 #else
842 snap_->atomData.force[atom2] += fg;
843 #endif
844 }
845
846 // filling interaction blocks with pointers
847 void ForceMatrixDecomposition::fillInteractionData(InteractionData &idat,
848 int atom1, int atom2) {
849
850 idat.excluded = excludeAtomPair(atom1, atom2);
851
852 #ifdef IS_MPI
853 idat.atypes = make_pair( atypesRow[atom1], atypesCol[atom2]);
854 //idat.atypes = make_pair( ff_->getAtomType(identsRow[atom1]),
855 // ff_->getAtomType(identsCol[atom2]) );
856
857 if (storageLayout_ & DataStorage::dslAmat) {
858 idat.A1 = &(atomRowData.aMat[atom1]);
859 idat.A2 = &(atomColData.aMat[atom2]);
860 }
861
862 if (storageLayout_ & DataStorage::dslElectroFrame) {
863 idat.eFrame1 = &(atomRowData.electroFrame[atom1]);
864 idat.eFrame2 = &(atomColData.electroFrame[atom2]);
865 }
866
867 if (storageLayout_ & DataStorage::dslTorque) {
868 idat.t1 = &(atomRowData.torque[atom1]);
869 idat.t2 = &(atomColData.torque[atom2]);
870 }
871
872 if (storageLayout_ & DataStorage::dslDensity) {
873 idat.rho1 = &(atomRowData.density[atom1]);
874 idat.rho2 = &(atomColData.density[atom2]);
875 }
876
877 if (storageLayout_ & DataStorage::dslFunctional) {
878 idat.frho1 = &(atomRowData.functional[atom1]);
879 idat.frho2 = &(atomColData.functional[atom2]);
880 }
881
882 if (storageLayout_ & DataStorage::dslFunctionalDerivative) {
883 idat.dfrho1 = &(atomRowData.functionalDerivative[atom1]);
884 idat.dfrho2 = &(atomColData.functionalDerivative[atom2]);
885 }
886
887 if (storageLayout_ & DataStorage::dslParticlePot) {
888 idat.particlePot1 = &(atomRowData.particlePot[atom1]);
889 idat.particlePot2 = &(atomColData.particlePot[atom2]);
890 }
891
892 if (storageLayout_ & DataStorage::dslSkippedCharge) {
893 idat.skippedCharge1 = &(atomRowData.skippedCharge[atom1]);
894 idat.skippedCharge2 = &(atomColData.skippedCharge[atom2]);
895 }
896
897 #else
898
899 idat.atypes = make_pair( atypesLocal[atom1], atypesLocal[atom2]);
900 //idat.atypes = make_pair( ff_->getAtomType(idents[atom1]),
901 // ff_->getAtomType(idents[atom2]) );
902
903 if (storageLayout_ & DataStorage::dslAmat) {
904 idat.A1 = &(snap_->atomData.aMat[atom1]);
905 idat.A2 = &(snap_->atomData.aMat[atom2]);
906 }
907
908 if (storageLayout_ & DataStorage::dslElectroFrame) {
909 idat.eFrame1 = &(snap_->atomData.electroFrame[atom1]);
910 idat.eFrame2 = &(snap_->atomData.electroFrame[atom2]);
911 }
912
913 if (storageLayout_ & DataStorage::dslTorque) {
914 idat.t1 = &(snap_->atomData.torque[atom1]);
915 idat.t2 = &(snap_->atomData.torque[atom2]);
916 }
917
918 if (storageLayout_ & DataStorage::dslDensity) {
919 idat.rho1 = &(snap_->atomData.density[atom1]);
920 idat.rho2 = &(snap_->atomData.density[atom2]);
921 }
922
923 if (storageLayout_ & DataStorage::dslFunctional) {
924 idat.frho1 = &(snap_->atomData.functional[atom1]);
925 idat.frho2 = &(snap_->atomData.functional[atom2]);
926 }
927
928 if (storageLayout_ & DataStorage::dslFunctionalDerivative) {
929 idat.dfrho1 = &(snap_->atomData.functionalDerivative[atom1]);
930 idat.dfrho2 = &(snap_->atomData.functionalDerivative[atom2]);
931 }
932
933 if (storageLayout_ & DataStorage::dslParticlePot) {
934 idat.particlePot1 = &(snap_->atomData.particlePot[atom1]);
935 idat.particlePot2 = &(snap_->atomData.particlePot[atom2]);
936 }
937
938 if (storageLayout_ & DataStorage::dslSkippedCharge) {
939 idat.skippedCharge1 = &(snap_->atomData.skippedCharge[atom1]);
940 idat.skippedCharge2 = &(snap_->atomData.skippedCharge[atom2]);
941 }
942 #endif
943 }
944
945
946 void ForceMatrixDecomposition::unpackInteractionData(InteractionData &idat, int atom1, int atom2) {
947 #ifdef IS_MPI
948 pot_row[atom1] += 0.5 * *(idat.pot);
949 pot_col[atom2] += 0.5 * *(idat.pot);
950
951 atomRowData.force[atom1] += *(idat.f1);
952 atomColData.force[atom2] -= *(idat.f1);
953 #else
954 pairwisePot += *(idat.pot);
955
956 snap_->atomData.force[atom1] += *(idat.f1);
957 snap_->atomData.force[atom2] -= *(idat.f1);
958 #endif
959
960 }
961
962 /*
963 * buildNeighborList
964 *
965 * first element of pair is row-indexed CutoffGroup
966 * second element of pair is column-indexed CutoffGroup
967 */
968 vector<pair<int, int> > ForceMatrixDecomposition::buildNeighborList() {
969
970 vector<pair<int, int> > neighborList;
971 groupCutoffs cuts;
972 bool doAllPairs = false;
973
974 #ifdef IS_MPI
975 cellListRow_.clear();
976 cellListCol_.clear();
977 #else
978 cellList_.clear();
979 #endif
980
981 RealType rList_ = (largestRcut_ + skinThickness_);
982 RealType rl2 = rList_ * rList_;
983 Snapshot* snap_ = sman_->getCurrentSnapshot();
984 Mat3x3d Hmat = snap_->getHmat();
985 Vector3d Hx = Hmat.getColumn(0);
986 Vector3d Hy = Hmat.getColumn(1);
987 Vector3d Hz = Hmat.getColumn(2);
988
989 nCells_.x() = (int) ( Hx.length() )/ rList_;
990 nCells_.y() = (int) ( Hy.length() )/ rList_;
991 nCells_.z() = (int) ( Hz.length() )/ rList_;
992
993 // handle small boxes where the cell offsets can end up repeating cells
994
995 if (nCells_.x() < 3) doAllPairs = true;
996 if (nCells_.y() < 3) doAllPairs = true;
997 if (nCells_.z() < 3) doAllPairs = true;
998
999 Mat3x3d invHmat = snap_->getInvHmat();
1000 Vector3d rs, scaled, dr;
1001 Vector3i whichCell;
1002 int cellIndex;
1003 int nCtot = nCells_.x() * nCells_.y() * nCells_.z();
1004
1005 #ifdef IS_MPI
1006 cellListRow_.resize(nCtot);
1007 cellListCol_.resize(nCtot);
1008 #else
1009 cellList_.resize(nCtot);
1010 #endif
1011
1012 if (!doAllPairs) {
1013 #ifdef IS_MPI
1014
1015 for (int i = 0; i < nGroupsInRow_; i++) {
1016 rs = cgRowData.position[i];
1017
1018 // scaled positions relative to the box vectors
1019 scaled = invHmat * rs;
1020
1021 // wrap the vector back into the unit box by subtracting integer box
1022 // numbers
1023 for (int j = 0; j < 3; j++) {
1024 scaled[j] -= roundMe(scaled[j]);
1025 scaled[j] += 0.5;
1026 }
1027
1028 // find xyz-indices of cell that cutoffGroup is in.
1029 whichCell.x() = nCells_.x() * scaled.x();
1030 whichCell.y() = nCells_.y() * scaled.y();
1031 whichCell.z() = nCells_.z() * scaled.z();
1032
1033 // find single index of this cell:
1034 cellIndex = Vlinear(whichCell, nCells_);
1035
1036 // add this cutoff group to the list of groups in this cell;
1037 cellListRow_[cellIndex].push_back(i);
1038 }
1039
1040 for (int i = 0; i < nGroupsInCol_; i++) {
1041 rs = cgColData.position[i];
1042
1043 // scaled positions relative to the box vectors
1044 scaled = invHmat * rs;
1045
1046 // wrap the vector back into the unit box by subtracting integer box
1047 // numbers
1048 for (int j = 0; j < 3; j++) {
1049 scaled[j] -= roundMe(scaled[j]);
1050 scaled[j] += 0.5;
1051 }
1052
1053 // find xyz-indices of cell that cutoffGroup is in.
1054 whichCell.x() = nCells_.x() * scaled.x();
1055 whichCell.y() = nCells_.y() * scaled.y();
1056 whichCell.z() = nCells_.z() * scaled.z();
1057
1058 // find single index of this cell:
1059 cellIndex = Vlinear(whichCell, nCells_);
1060
1061 // add this cutoff group to the list of groups in this cell;
1062 cellListCol_[cellIndex].push_back(i);
1063 }
1064 #else
1065 for (int i = 0; i < nGroups_; i++) {
1066 rs = snap_->cgData.position[i];
1067
1068 // scaled positions relative to the box vectors
1069 scaled = invHmat * rs;
1070
1071 // wrap the vector back into the unit box by subtracting integer box
1072 // numbers
1073 for (int j = 0; j < 3; j++) {
1074 scaled[j] -= roundMe(scaled[j]);
1075 scaled[j] += 0.5;
1076 }
1077
1078 // find xyz-indices of cell that cutoffGroup is in.
1079 whichCell.x() = nCells_.x() * scaled.x();
1080 whichCell.y() = nCells_.y() * scaled.y();
1081 whichCell.z() = nCells_.z() * scaled.z();
1082
1083 // find single index of this cell:
1084 cellIndex = Vlinear(whichCell, nCells_);
1085
1086 // add this cutoff group to the list of groups in this cell;
1087 cellList_[cellIndex].push_back(i);
1088 }
1089 #endif
1090
1091 for (int m1z = 0; m1z < nCells_.z(); m1z++) {
1092 for (int m1y = 0; m1y < nCells_.y(); m1y++) {
1093 for (int m1x = 0; m1x < nCells_.x(); m1x++) {
1094 Vector3i m1v(m1x, m1y, m1z);
1095 int m1 = Vlinear(m1v, nCells_);
1096
1097 for (vector<Vector3i>::iterator os = cellOffsets_.begin();
1098 os != cellOffsets_.end(); ++os) {
1099
1100 Vector3i m2v = m1v + (*os);
1101
1102 if (m2v.x() >= nCells_.x()) {
1103 m2v.x() = 0;
1104 } else if (m2v.x() < 0) {
1105 m2v.x() = nCells_.x() - 1;
1106 }
1107
1108 if (m2v.y() >= nCells_.y()) {
1109 m2v.y() = 0;
1110 } else if (m2v.y() < 0) {
1111 m2v.y() = nCells_.y() - 1;
1112 }
1113
1114 if (m2v.z() >= nCells_.z()) {
1115 m2v.z() = 0;
1116 } else if (m2v.z() < 0) {
1117 m2v.z() = nCells_.z() - 1;
1118 }
1119
1120 int m2 = Vlinear (m2v, nCells_);
1121
1122 #ifdef IS_MPI
1123 for (vector<int>::iterator j1 = cellListRow_[m1].begin();
1124 j1 != cellListRow_[m1].end(); ++j1) {
1125 for (vector<int>::iterator j2 = cellListCol_[m2].begin();
1126 j2 != cellListCol_[m2].end(); ++j2) {
1127
1128 // Always do this if we're in different cells or if
1129 // we're in the same cell and the global index of the
1130 // j2 cutoff group is less than the j1 cutoff group
1131
1132 if (m2 != m1 || cgColToGlobal[(*j2)] < cgRowToGlobal[(*j1)]) {
1133 dr = cgColData.position[(*j2)] - cgRowData.position[(*j1)];
1134 snap_->wrapVector(dr);
1135 cuts = getGroupCutoffs( (*j1), (*j2) );
1136 if (dr.lengthSquare() < cuts.third) {
1137 neighborList.push_back(make_pair((*j1), (*j2)));
1138 }
1139 }
1140 }
1141 }
1142 #else
1143
1144 for (vector<int>::iterator j1 = cellList_[m1].begin();
1145 j1 != cellList_[m1].end(); ++j1) {
1146 for (vector<int>::iterator j2 = cellList_[m2].begin();
1147 j2 != cellList_[m2].end(); ++j2) {
1148
1149 // Always do this if we're in different cells or if
1150 // we're in the same cell and the global index of the
1151 // j2 cutoff group is less than the j1 cutoff group
1152
1153 if (m2 != m1 || (*j2) < (*j1)) {
1154 dr = snap_->cgData.position[(*j2)] - snap_->cgData.position[(*j1)];
1155 snap_->wrapVector(dr);
1156 cuts = getGroupCutoffs( (*j1), (*j2) );
1157 if (dr.lengthSquare() < cuts.third) {
1158 neighborList.push_back(make_pair((*j1), (*j2)));
1159 }
1160 }
1161 }
1162 }
1163 #endif
1164 }
1165 }
1166 }
1167 }
1168 } else {
1169 // branch to do all cutoff group pairs
1170 #ifdef IS_MPI
1171 for (int j1 = 0; j1 < nGroupsInRow_; j1++) {
1172 for (int j2 = 0; j2 < nGroupsInCol_; j2++) {
1173 dr = cgColData.position[j2] - cgRowData.position[j1];
1174 snap_->wrapVector(dr);
1175 cuts = getGroupCutoffs( j1, j2 );
1176 if (dr.lengthSquare() < cuts.third) {
1177 neighborList.push_back(make_pair(j1, j2));
1178 }
1179 }
1180 }
1181 #else
1182 for (int j1 = 0; j1 < nGroups_ - 1; j1++) {
1183 for (int j2 = j1 + 1; j2 < nGroups_; j2++) {
1184 dr = snap_->cgData.position[j2] - snap_->cgData.position[j1];
1185 snap_->wrapVector(dr);
1186 cuts = getGroupCutoffs( j1, j2 );
1187 if (dr.lengthSquare() < cuts.third) {
1188 neighborList.push_back(make_pair(j1, j2));
1189 }
1190 }
1191 }
1192 #endif
1193 }
1194
1195 // save the local cutoff group positions for the check that is
1196 // done on each loop:
1197 saved_CG_positions_.clear();
1198 for (int i = 0; i < nGroups_; i++)
1199 saved_CG_positions_.push_back(snap_->cgData.position[i]);
1200
1201 return neighborList;
1202 }
1203 } //end namespace OpenMD