ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/devel_omp/src/parallel/ForceMatrixDecomposition.cpp
Revision: 1608
Committed: Tue Aug 9 01:58:56 2011 UTC (13 years, 8 months ago) by mciznick
File size: 50130 byte(s)
Log Message:
First OpenMP version.

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 #include "primitives/Molecule.hpp"
47
48 using namespace std;
49 namespace OpenMD {
50
51 ForceMatrixDecomposition::ForceMatrixDecomposition(SimInfo* info, InteractionManager* iMan) :
52 ForceDecomposition(info, iMan) {
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_.push_back( Vector3i(-1, 0, 0) );
58 cellOffsets_.push_back( Vector3i(-1,-1, 0) );
59 cellOffsets_.push_back( Vector3i( 0,-1, 0) );
60 cellOffsets_.push_back( Vector3i( 1,-1, 0) );
61 cellOffsets_.push_back( Vector3i( 0, 0,-1) );
62 cellOffsets_.push_back( Vector3i(-1, 0, 1) );
63 cellOffsets_.push_back( Vector3i(-1,-1,-1) );
64 cellOffsets_.push_back( Vector3i( 0,-1,-1) );
65 cellOffsets_.push_back( Vector3i( 1,-1,-1) );
66 cellOffsets_.push_back( Vector3i( 1, 0,-1) );
67 cellOffsets_.push_back( Vector3i( 1, 1,-1) );
68 cellOffsets_.push_back( Vector3i( 0, 1,-1) );
69 cellOffsets_.push_back( Vector3i(-1, 1,-1) );
70 #endif
71 }
72
73 /**
74 * distributeInitialData is essentially a copy of the older fortran
75 * SimulationSetup
76 */
77 void ForceMatrixDecomposition::distributeInitialData() {
78 snap_ = sman_->getCurrentSnapshot();
79 storageLayout_ = sman_->getStorageLayout();
80 ff_ = info_->getForceField();
81 nLocal_ = snap_->getNumberOfAtoms();
82
83 nGroups_ = info_->getNLocalCutoffGroups();
84 // gather the information for atomtype IDs (atids):
85 idents = info_->getIdentArray();
86 AtomLocalToGlobal = info_->getGlobalAtomIndices();
87 cgLocalToGlobal = info_->getGlobalGroupIndices();
88 vector<int> globalGroupMembership = info_->getGlobalGroupMembership();
89
90 massFactors = info_->getMassFactors();
91
92 PairList* excludes = info_->getExcludedInteractions();
93 PairList* oneTwo = info_->getOneTwoInteractions();
94 PairList* oneThree = info_->getOneThreeInteractions();
95 PairList* oneFour = info_->getOneFourInteractions();
96
97 #ifdef IS_MPI
98
99 MPI::Intracomm row = rowComm.getComm();
100 MPI::Intracomm col = colComm.getComm();
101
102 AtomPlanIntRow = new Plan<int>(row, nLocal_);
103 AtomPlanRealRow = new Plan<RealType>(row, nLocal_);
104 AtomPlanVectorRow = new Plan<Vector3d>(row, nLocal_);
105 AtomPlanMatrixRow = new Plan<Mat3x3d>(row, nLocal_);
106 AtomPlanPotRow = new Plan<potVec>(row, nLocal_);
107
108 AtomPlanIntColumn = new Plan<int>(col, nLocal_);
109 AtomPlanRealColumn = new Plan<RealType>(col, nLocal_);
110 AtomPlanVectorColumn = new Plan<Vector3d>(col, nLocal_);
111 AtomPlanMatrixColumn = new Plan<Mat3x3d>(col, nLocal_);
112 AtomPlanPotColumn = new Plan<potVec>(col, nLocal_);
113
114 cgPlanIntRow = new Plan<int>(row, nGroups_);
115 cgPlanVectorRow = new Plan<Vector3d>(row, nGroups_);
116 cgPlanIntColumn = new Plan<int>(col, nGroups_);
117 cgPlanVectorColumn = new Plan<Vector3d>(col, nGroups_);
118
119 nAtomsInRow_ = AtomPlanIntRow->getSize();
120 nAtomsInCol_ = AtomPlanIntColumn->getSize();
121 nGroupsInRow_ = cgPlanIntRow->getSize();
122 nGroupsInCol_ = cgPlanIntColumn->getSize();
123
124 // Modify the data storage objects with the correct layouts and sizes:
125 atomRowData.resize(nAtomsInRow_);
126 atomRowData.setStorageLayout(storageLayout_);
127 atomColData.resize(nAtomsInCol_);
128 atomColData.setStorageLayout(storageLayout_);
129 cgRowData.resize(nGroupsInRow_);
130 cgRowData.setStorageLayout(DataStorage::dslPosition);
131 cgColData.resize(nGroupsInCol_);
132 cgColData.setStorageLayout(DataStorage::dslPosition);
133
134 identsRow.resize(nAtomsInRow_);
135 identsCol.resize(nAtomsInCol_);
136
137 AtomPlanIntRow->gather(idents, identsRow);
138 AtomPlanIntColumn->gather(idents, identsCol);
139
140 // allocate memory for the parallel objects
141 atypesRow.resize(nAtomsInRow_);
142 atypesCol.resize(nAtomsInCol_);
143
144 for (int i = 0; i < nAtomsInRow_; i++)
145 atypesRow[i] = ff_->getAtomType(identsRow[i]);
146 for (int i = 0; i < nAtomsInCol_; i++)
147 atypesCol[i] = ff_->getAtomType(identsCol[i]);
148
149 pot_row.resize(nAtomsInRow_);
150 pot_col.resize(nAtomsInCol_);
151
152 AtomRowToGlobal.resize(nAtomsInRow_);
153 AtomColToGlobal.resize(nAtomsInCol_);
154 AtomPlanIntRow->gather(AtomLocalToGlobal, AtomRowToGlobal);
155 AtomPlanIntColumn->gather(AtomLocalToGlobal, AtomColToGlobal);
156
157 cerr << "Atoms in Local:\n";
158 for (int i = 0; i < AtomLocalToGlobal.size(); i++)
159 {
160 cerr << "i =\t" << i << "\t localAt =\t" << AtomLocalToGlobal[i] << "\n";
161 }
162 cerr << "Atoms in Row:\n";
163 for (int i = 0; i < AtomRowToGlobal.size(); i++)
164 {
165 cerr << "i =\t" << i << "\t rowAt =\t" << AtomRowToGlobal[i] << "\n";
166 }
167 cerr << "Atoms in Col:\n";
168 for (int i = 0; i < AtomColToGlobal.size(); i++)
169 {
170 cerr << "i =\t" << i << "\t colAt =\t" << AtomColToGlobal[i] << "\n";
171 }
172
173 cgRowToGlobal.resize(nGroupsInRow_);
174 cgColToGlobal.resize(nGroupsInCol_);
175 cgPlanIntRow->gather(cgLocalToGlobal, cgRowToGlobal);
176 cgPlanIntColumn->gather(cgLocalToGlobal, cgColToGlobal);
177
178 cerr << "Gruops in Local:\n";
179 for (int i = 0; i < cgLocalToGlobal.size(); i++)
180 {
181 cerr << "i =\t" << i << "\t localCG =\t" << cgLocalToGlobal[i] << "\n";
182 }
183 cerr << "Groups in Row:\n";
184 for (int i = 0; i < cgRowToGlobal.size(); i++)
185 {
186 cerr << "i =\t" << i << "\t rowCG =\t" << cgRowToGlobal[i] << "\n";
187 }
188 cerr << "Groups in Col:\n";
189 for (int i = 0; i < cgColToGlobal.size(); i++)
190 {
191 cerr << "i =\t" << i << "\t colCG =\t" << cgColToGlobal[i] << "\n";
192 }
193
194 massFactorsRow.resize(nAtomsInRow_);
195 massFactorsCol.resize(nAtomsInCol_);
196 AtomPlanRealRow->gather(massFactors, massFactorsRow);
197 AtomPlanRealColumn->gather(massFactors, massFactorsCol);
198
199 groupListRow_.clear();
200 groupListRow_.resize(nGroupsInRow_);
201 for (int i = 0; i < nGroupsInRow_; i++)
202 {
203 int gid = cgRowToGlobal[i];
204 for (int j = 0; j < nAtomsInRow_; j++)
205 {
206 int aid = AtomRowToGlobal[j];
207 if (globalGroupMembership[aid] == gid)
208 groupListRow_[i].push_back(j);
209 }
210 }
211
212 groupListCol_.clear();
213 groupListCol_.resize(nGroupsInCol_);
214 for (int i = 0; i < nGroupsInCol_; i++)
215 {
216 int gid = cgColToGlobal[i];
217 for (int j = 0; j < nAtomsInCol_; j++)
218 {
219 int aid = AtomColToGlobal[j];
220 if (globalGroupMembership[aid] == gid)
221 groupListCol_[i].push_back(j);
222 }
223 }
224
225 excludesForAtom.clear();
226 excludesForAtom.resize(nAtomsInRow_);
227 toposForAtom.clear();
228 toposForAtom.resize(nAtomsInRow_);
229 topoDist.clear();
230 topoDist.resize(nAtomsInRow_);
231 for (int i = 0; i < nAtomsInRow_; i++)
232 {
233 int iglob = AtomRowToGlobal[i];
234
235 for (int j = 0; j < nAtomsInCol_; j++)
236 {
237 int jglob = AtomColToGlobal[j];
238
239 if (excludes->hasPair(iglob, jglob))
240 excludesForAtom[i].push_back(j);
241
242 if (oneTwo->hasPair(iglob, jglob))
243 {
244 toposForAtom[i].push_back(j);
245 topoDist[i].push_back(1);
246 } else
247 {
248 if (oneThree->hasPair(iglob, jglob))
249 {
250 toposForAtom[i].push_back(j);
251 topoDist[i].push_back(2);
252 } else
253 {
254 if (oneFour->hasPair(iglob, jglob))
255 {
256 toposForAtom[i].push_back(j);
257 topoDist[i].push_back(3);
258 }
259 }
260 }
261 }
262 }
263
264 #endif
265
266 // allocate memory for the parallel objects
267 atypesLocal.resize(nLocal_);
268
269 for (int i = 0; i < nLocal_; i++)
270 atypesLocal[i] = ff_->getAtomType(idents[i]);
271
272 groupList_.clear();
273 groupList_.resize(nGroups_);
274 for (int i = 0; i < nGroups_; i++)
275 {
276 int gid = cgLocalToGlobal[i];
277 for (int j = 0; j < nLocal_; j++)
278 {
279 int aid = AtomLocalToGlobal[j];
280 if (globalGroupMembership[aid] == gid)
281 {
282 groupList_[i].push_back(j);
283 }
284 }
285 }
286
287 excludesForAtom.clear();
288 excludesForAtom.resize(nLocal_);
289 toposForAtom.clear();
290 toposForAtom.resize(nLocal_);
291 topoDist.clear();
292 topoDist.resize(nLocal_);
293
294 for (int i = 0; i < nLocal_; i++)
295 {
296 int iglob = AtomLocalToGlobal[i];
297
298 for (int j = 0; j < nLocal_; j++)
299 {
300 int jglob = AtomLocalToGlobal[j];
301
302 if (excludes->hasPair(iglob, jglob))
303 excludesForAtom[i].push_back(j);
304
305 if (oneTwo->hasPair(iglob, jglob))
306 {
307 toposForAtom[i].push_back(j);
308 topoDist[i].push_back(1);
309 } else
310 {
311 if (oneThree->hasPair(iglob, jglob))
312 {
313 toposForAtom[i].push_back(j);
314 topoDist[i].push_back(2);
315 } else
316 {
317 if (oneFour->hasPair(iglob, jglob))
318 {
319 toposForAtom[i].push_back(j);
320 topoDist[i].push_back(3);
321 }
322 }
323 }
324 }
325 }
326
327 Globals* simParams_ = info_->getSimParams();
328 if (simParams_->haveNeighborListReorderFreq())
329 {
330 neighborListReorderFreq = simParams_->getNeighborListReorderFreq();
331 } else
332 {
333 neighborListReorderFreq = 0;
334 }
335 reorderFreqCounter = 1;
336
337 createGtypeCutoffMap();
338
339 }
340
341 void ForceMatrixDecomposition::createGtypeCutoffMap() {
342
343 RealType tol = 1e-6;
344 largestRcut_ = 0.0;
345 RealType rc;
346 int atid;
347 set<AtomType*> atypes = info_->getSimulatedAtomTypes();
348
349 map<int, RealType> atypeCutoff;
350
351 for (set<AtomType*>::iterator at = atypes.begin(); at != atypes.end(); ++at)
352 {
353 atid = (*at)->getIdent();
354 if (userChoseCutoff_)
355 atypeCutoff[atid] = userCutoff_;
356 else
357 atypeCutoff[atid] = interactionMan_->getSuggestedCutoffRadius(*at);
358 }
359
360 vector<RealType> gTypeCutoffs;
361 // first we do a single loop over the cutoff groups to find the
362 // largest cutoff for any atypes present in this group.
363 #ifdef IS_MPI
364 vector<RealType> groupCutoffRow(nGroupsInRow_, 0.0);
365 groupRowToGtype.resize(nGroupsInRow_);
366 for (int cg1 = 0; cg1 < nGroupsInRow_; cg1++)
367 {
368 vector<int> atomListRow = getAtomsInGroupRow(cg1);
369 for (vector<int>::iterator ia = atomListRow.begin();
370 ia != atomListRow.end(); ++ia)
371 {
372 int atom1 = (*ia);
373 atid = identsRow[atom1];
374 if (atypeCutoff[atid] > groupCutoffRow[cg1])
375 {
376 groupCutoffRow[cg1] = atypeCutoff[atid];
377 }
378 }
379
380 bool gTypeFound = false;
381 for (int gt = 0; gt < gTypeCutoffs.size(); gt++)
382 {
383 if (abs(groupCutoffRow[cg1] - gTypeCutoffs[gt]) < tol)
384 {
385 groupRowToGtype[cg1] = gt;
386 gTypeFound = true;
387 }
388 }
389 if (!gTypeFound)
390 {
391 gTypeCutoffs.push_back( groupCutoffRow[cg1] );
392 groupRowToGtype[cg1] = gTypeCutoffs.size() - 1;
393 }
394
395 }
396 vector<RealType> groupCutoffCol(nGroupsInCol_, 0.0);
397 groupColToGtype.resize(nGroupsInCol_);
398 for (int cg2 = 0; cg2 < nGroupsInCol_; cg2++)
399 {
400 vector<int> atomListCol = getAtomsInGroupColumn(cg2);
401 for (vector<int>::iterator jb = atomListCol.begin();
402 jb != atomListCol.end(); ++jb)
403 {
404 int atom2 = (*jb);
405 atid = identsCol[atom2];
406 if (atypeCutoff[atid] > groupCutoffCol[cg2])
407 {
408 groupCutoffCol[cg2] = atypeCutoff[atid];
409 }
410 }
411 bool gTypeFound = false;
412 for (int gt = 0; gt < gTypeCutoffs.size(); gt++)
413 {
414 if (abs(groupCutoffCol[cg2] - gTypeCutoffs[gt]) < tol)
415 {
416 groupColToGtype[cg2] = gt;
417 gTypeFound = true;
418 }
419 }
420 if (!gTypeFound)
421 {
422 gTypeCutoffs.push_back( groupCutoffCol[cg2] );
423 groupColToGtype[cg2] = gTypeCutoffs.size() - 1;
424 }
425 }
426 #else
427
428 vector<RealType> groupCutoff(nGroups_, 0.0);
429 groupToGtype.resize(nGroups_);
430 for (int cg1 = 0; cg1 < nGroups_; cg1++)
431 {
432 groupCutoff[cg1] = 0.0;
433 vector<int> atomList = getAtomsInGroupRow(cg1);
434 for (vector<int>::iterator ia = atomList.begin(); ia != atomList.end(); ++ia)
435 {
436 int atom1 = (*ia);
437 atid = idents[atom1];
438 if (atypeCutoff[atid] > groupCutoff[cg1])
439 groupCutoff[cg1] = atypeCutoff[atid];
440 }
441
442 bool gTypeFound = false;
443 for (int gt = 0; gt < gTypeCutoffs.size(); gt++)
444 {
445 if (abs(groupCutoff[cg1] - gTypeCutoffs[gt]) < tol)
446 {
447 groupToGtype[cg1] = gt;
448 gTypeFound = true;
449 }
450 }
451 if (!gTypeFound)
452 {
453 gTypeCutoffs.push_back(groupCutoff[cg1]);
454 groupToGtype[cg1] = gTypeCutoffs.size() - 1;
455 }
456 }
457 #endif
458
459 // Now we find the maximum group cutoff value present in the simulation
460
461 RealType groupMax = *max_element(gTypeCutoffs.begin(), gTypeCutoffs.end());
462
463 #ifdef IS_MPI
464 MPI::COMM_WORLD.Allreduce(&groupMax, &groupMax, 1, MPI::REALTYPE,
465 MPI::MAX);
466 #endif
467
468 RealType tradRcut = groupMax;
469
470 for (int i = 0; i < gTypeCutoffs.size(); i++)
471 {
472 for (int j = 0; j < gTypeCutoffs.size(); j++)
473 {
474 RealType thisRcut;
475 switch (cutoffPolicy_) {
476 case TRADITIONAL:
477 thisRcut = tradRcut;
478 break;
479 case MIX:
480 thisRcut = 0.5 * (gTypeCutoffs[i] + gTypeCutoffs[j]);
481 break;
482 case MAX:
483 thisRcut = max(gTypeCutoffs[i], gTypeCutoffs[j]);
484 break;
485 default:
486 sprintf(painCave.errMsg, "ForceMatrixDecomposition::createGtypeCutoffMap "
487 "hit an unknown cutoff policy!\n");
488 painCave.severity = OPENMD_ERROR;
489 painCave.isFatal = 1;
490 simError();
491 break;
492 }
493
494 pair<int, int> key = make_pair(i, j);
495 gTypeCutoffMap[key].first = thisRcut;
496 if (thisRcut > largestRcut_)
497 largestRcut_ = thisRcut;
498 gTypeCutoffMap[key].second = thisRcut * thisRcut;
499 gTypeCutoffMap[key].third = pow(thisRcut + skinThickness_, 2);
500 // sanity check
501
502 if (userChoseCutoff_)
503 {
504 if (abs(gTypeCutoffMap[key].first - userCutoff_) > 0.0001)
505 {
506 sprintf(painCave.errMsg, "ForceMatrixDecomposition::createGtypeCutoffMap "
507 "user-specified rCut (%lf) does not match computed group Cutoff\n", userCutoff_);
508 painCave.severity = OPENMD_ERROR;
509 painCave.isFatal = 1;
510 simError();
511 }
512 }
513 }
514 }
515 }
516
517 groupCutoffs ForceMatrixDecomposition::getGroupCutoffs(int cg1, int cg2) {
518 int i, j;
519 #ifdef IS_MPI
520 i = groupRowToGtype[cg1];
521 j = groupColToGtype[cg2];
522 #else
523 i = groupToGtype[cg1];
524 j = groupToGtype[cg2];
525 #endif
526 return gTypeCutoffMap[make_pair(i, j)];
527 }
528
529 int ForceMatrixDecomposition::getTopologicalDistance(int atom1, int atom2) {
530 for (int j = 0; j < toposForAtom[atom1].size(); j++)
531 {
532 if (toposForAtom[atom1][j] == atom2)
533 return topoDist[atom1][j];
534 }
535 return 0;
536 }
537
538 void ForceMatrixDecomposition::zeroWorkArrays() {
539 pairwisePot = 0.0;
540 embeddingPot = 0.0;
541
542 #ifdef IS_MPI
543 if (storageLayout_ & DataStorage::dslForce)
544 {
545 fill(atomRowData.force.begin(), atomRowData.force.end(), V3Zero);
546 fill(atomColData.force.begin(), atomColData.force.end(), V3Zero);
547 }
548
549 if (storageLayout_ & DataStorage::dslTorque)
550 {
551 fill(atomRowData.torque.begin(), atomRowData.torque.end(), V3Zero);
552 fill(atomColData.torque.begin(), atomColData.torque.end(), V3Zero);
553 }
554
555 fill(pot_row.begin(), pot_row.end(),
556 Vector<RealType, N_INTERACTION_FAMILIES> (0.0));
557
558 fill(pot_col.begin(), pot_col.end(),
559 Vector<RealType, N_INTERACTION_FAMILIES> (0.0));
560
561 if (storageLayout_ & DataStorage::dslParticlePot)
562 {
563 fill(atomRowData.particlePot.begin(), atomRowData.particlePot.end(),
564 0.0);
565 fill(atomColData.particlePot.begin(), atomColData.particlePot.end(),
566 0.0);
567 }
568
569 if (storageLayout_ & DataStorage::dslDensity)
570 {
571 fill(atomRowData.density.begin(), atomRowData.density.end(), 0.0);
572 fill(atomColData.density.begin(), atomColData.density.end(), 0.0);
573 }
574
575 if (storageLayout_ & DataStorage::dslFunctional)
576 {
577 fill(atomRowData.functional.begin(), atomRowData.functional.end(),
578 0.0);
579 fill(atomColData.functional.begin(), atomColData.functional.end(),
580 0.0);
581 }
582
583 if (storageLayout_ & DataStorage::dslFunctionalDerivative)
584 {
585 fill(atomRowData.functionalDerivative.begin(),
586 atomRowData.functionalDerivative.end(), 0.0);
587 fill(atomColData.functionalDerivative.begin(),
588 atomColData.functionalDerivative.end(), 0.0);
589 }
590
591 if (storageLayout_ & DataStorage::dslSkippedCharge)
592 {
593 fill(atomRowData.skippedCharge.begin(),
594 atomRowData.skippedCharge.end(), 0.0);
595 fill(atomColData.skippedCharge.begin(),
596 atomColData.skippedCharge.end(), 0.0);
597 }
598
599 #endif
600 // even in parallel, we need to zero out the local arrays:
601
602 if (storageLayout_ & DataStorage::dslParticlePot)
603 {
604 fill(snap_->atomData.particlePot.begin(), snap_->atomData.particlePot.end(), 0.0);
605 }
606
607 if (storageLayout_ & DataStorage::dslDensity)
608 {
609 fill(snap_->atomData.density.begin(), snap_->atomData.density.end(), 0.0);
610 }
611 if (storageLayout_ & DataStorage::dslFunctional)
612 {
613 fill(snap_->atomData.functional.begin(), snap_->atomData.functional.end(), 0.0);
614 }
615 if (storageLayout_ & DataStorage::dslFunctionalDerivative)
616 {
617 fill(snap_->atomData.functionalDerivative.begin(), snap_->atomData.functionalDerivative.end(), 0.0);
618 }
619 if (storageLayout_ & DataStorage::dslSkippedCharge)
620 {
621 fill(snap_->atomData.skippedCharge.begin(), snap_->atomData.skippedCharge.end(), 0.0);
622 }
623
624 }
625
626 void ForceMatrixDecomposition::distributeData() {
627 snap_ = sman_->getCurrentSnapshot();
628 storageLayout_ = sman_->getStorageLayout();
629 #ifdef IS_MPI
630
631 // gather up the atomic positions
632 AtomPlanVectorRow->gather(snap_->atomData.position,
633 atomRowData.position);
634 AtomPlanVectorColumn->gather(snap_->atomData.position,
635 atomColData.position);
636
637 // gather up the cutoff group positions
638
639 cerr << "before gather\n";
640 for (int i = 0; i < snap_->cgData.position.size(); i++)
641 {
642 cerr << "cgpos = " << snap_->cgData.position[i] << "\n";
643 }
644
645 cgPlanVectorRow->gather(snap_->cgData.position,
646 cgRowData.position);
647
648 cerr << "after gather\n";
649 for (int i = 0; i < cgRowData.position.size(); i++)
650 {
651 cerr << "cgRpos = " << cgRowData.position[i] << "\n";
652 }
653
654 cgPlanVectorColumn->gather(snap_->cgData.position,
655 cgColData.position);
656 for (int i = 0; i < cgColData.position.size(); i++)
657 {
658 cerr << "cgCpos = " << cgColData.position[i] << "\n";
659 }
660
661 // if needed, gather the atomic rotation matrices
662 if (storageLayout_ & DataStorage::dslAmat)
663 {
664 AtomPlanMatrixRow->gather(snap_->atomData.aMat,
665 atomRowData.aMat);
666 AtomPlanMatrixColumn->gather(snap_->atomData.aMat,
667 atomColData.aMat);
668 }
669
670 // if needed, gather the atomic eletrostatic frames
671 if (storageLayout_ & DataStorage::dslElectroFrame)
672 {
673 AtomPlanMatrixRow->gather(snap_->atomData.electroFrame,
674 atomRowData.electroFrame);
675 AtomPlanMatrixColumn->gather(snap_->atomData.electroFrame,
676 atomColData.electroFrame);
677 }
678
679 #endif
680 }
681
682 /* collects information obtained during the pre-pair loop onto local
683 * data structures.
684 */
685 void ForceMatrixDecomposition::collectIntermediateData() {
686 snap_ = sman_->getCurrentSnapshot();
687 storageLayout_ = sman_->getStorageLayout();
688 #ifdef IS_MPI
689
690 if (storageLayout_ & DataStorage::dslDensity)
691 {
692
693 AtomPlanRealRow->scatter(atomRowData.density,
694 snap_->atomData.density);
695
696 int n = snap_->atomData.density.size();
697 vector<RealType> rho_tmp(n, 0.0);
698 AtomPlanRealColumn->scatter(atomColData.density, rho_tmp);
699 for (int i = 0; i < n; i++)
700 snap_->atomData.density[i] += rho_tmp[i];
701 }
702 #endif
703 }
704
705 /*
706 * redistributes information obtained during the pre-pair loop out to
707 * row and column-indexed data structures
708 */
709 void ForceMatrixDecomposition::distributeIntermediateData() {
710 snap_ = sman_->getCurrentSnapshot();
711 storageLayout_ = sman_->getStorageLayout();
712 #ifdef IS_MPI
713 if (storageLayout_ & DataStorage::dslFunctional)
714 {
715 AtomPlanRealRow->gather(snap_->atomData.functional,
716 atomRowData.functional);
717 AtomPlanRealColumn->gather(snap_->atomData.functional,
718 atomColData.functional);
719 }
720
721 if (storageLayout_ & DataStorage::dslFunctionalDerivative)
722 {
723 AtomPlanRealRow->gather(snap_->atomData.functionalDerivative,
724 atomRowData.functionalDerivative);
725 AtomPlanRealColumn->gather(snap_->atomData.functionalDerivative,
726 atomColData.functionalDerivative);
727 }
728 #endif
729 }
730
731 void ForceMatrixDecomposition::collectData() {
732 snap_ = sman_->getCurrentSnapshot();
733 storageLayout_ = sman_->getStorageLayout();
734 #ifdef IS_MPI
735 int n = snap_->atomData.force.size();
736 vector<Vector3d> frc_tmp(n, V3Zero);
737
738 AtomPlanVectorRow->scatter(atomRowData.force, frc_tmp);
739 for (int i = 0; i < n; i++)
740 {
741 snap_->atomData.force[i] += frc_tmp[i];
742 frc_tmp[i] = 0.0;
743 }
744
745 AtomPlanVectorColumn->scatter(atomColData.force, frc_tmp);
746 for (int i = 0; i < n; i++)
747 {
748 snap_->atomData.force[i] += frc_tmp[i];
749 }
750
751 if (storageLayout_ & DataStorage::dslTorque)
752 {
753
754 int nt = snap_->atomData.torque.size();
755 vector<Vector3d> trq_tmp(nt, V3Zero);
756
757 AtomPlanVectorRow->scatter(atomRowData.torque, trq_tmp);
758 for (int i = 0; i < nt; i++)
759 {
760 snap_->atomData.torque[i] += trq_tmp[i];
761 trq_tmp[i] = 0.0;
762 }
763
764 AtomPlanVectorColumn->scatter(atomColData.torque, trq_tmp);
765 for (int i = 0; i < nt; i++)
766 snap_->atomData.torque[i] += trq_tmp[i];
767 }
768
769 if (storageLayout_ & DataStorage::dslSkippedCharge)
770 {
771
772 int ns = snap_->atomData.skippedCharge.size();
773 vector<RealType> skch_tmp(ns, 0.0);
774
775 AtomPlanRealRow->scatter(atomRowData.skippedCharge, skch_tmp);
776 for (int i = 0; i < ns; i++)
777 {
778 snap_->atomData.skippedCharge[i] += skch_tmp[i];
779 skch_tmp[i] = 0.0;
780 }
781
782 AtomPlanRealColumn->scatter(atomColData.skippedCharge, skch_tmp);
783 for (int i = 0; i < ns; i++)
784 snap_->atomData.skippedCharge[i] += skch_tmp[i];
785 }
786
787 nLocal_ = snap_->getNumberOfAtoms();
788
789 vector<potVec> pot_temp(nLocal_,
790 Vector<RealType, N_INTERACTION_FAMILIES> (0.0));
791
792 // scatter/gather pot_row into the members of my column
793
794 AtomPlanPotRow->scatter(pot_row, pot_temp);
795
796 for (int ii = 0; ii < pot_temp.size(); ii++ )
797 pairwisePot += pot_temp[ii];
798
799 fill(pot_temp.begin(), pot_temp.end(),
800 Vector<RealType, N_INTERACTION_FAMILIES> (0.0));
801
802 AtomPlanPotColumn->scatter(pot_col, pot_temp);
803
804 for (int ii = 0; ii < pot_temp.size(); ii++ )
805 pairwisePot += pot_temp[ii];
806 #endif
807
808 // cerr << "pairwisePot = " << pairwisePot << "\n";
809 }
810
811 int ForceMatrixDecomposition::getNAtomsInRow() {
812 #ifdef IS_MPI
813 return nAtomsInRow_;
814 #else
815 return nLocal_;
816 #endif
817 }
818
819 /**
820 * returns the list of atoms belonging to this group.
821 */
822 vector<int> ForceMatrixDecomposition::getAtomsInGroupRow(int cg1) {
823 #ifdef IS_MPI
824 return groupListRow_[cg1];
825 #else
826 return groupList_[cg1];
827 #endif
828 }
829
830 vector<int> ForceMatrixDecomposition::getAtomsInGroupColumn(int cg2) {
831 #ifdef IS_MPI
832 return groupListCol_[cg2];
833 #else
834 return groupList_[cg2];
835 #endif
836 }
837
838 Vector3d ForceMatrixDecomposition::getIntergroupVector(int cg1, int cg2) {
839 Vector3d d;
840
841 #ifdef IS_MPI
842 d = cgColData.position[cg2] - cgRowData.position[cg1];
843 cerr << "cg1 = " << cg1 << "\tcg1p = " << cgRowData.position[cg1] << "\n";
844 cerr << "cg2 = " << cg2 << "\tcg2p = " << cgColData.position[cg2] << "\n";
845 #else
846 d = snap_->cgData.position[cg2] - snap_->cgData.position[cg1];
847 cerr << "cg1 = " << cg1 << "\tcg1p = " << snap_->cgData.position[cg1] << "\n";
848 cerr << "cg2 = " << cg2 << "\tcg2p = " << snap_->cgData.position[cg2] << "\n";
849 #endif
850
851 snap_->wrapVector(d);
852 return d;
853 }
854
855 Vector3d ForceMatrixDecomposition::getIntergroupVector(CutoffGroup *cg1, CutoffGroup *cg2) {
856 Vector3d d;
857
858 d = snap_->cgData.position[cg2->getLocalIndex()] - snap_->cgData.position[cg1->getLocalIndex()];
859 /* cerr << "cg1_gid = " << cg1->getGlobalIndex() << "\tcg1_lid = " << cg1->getLocalIndex() << "\tcg1p = "
860 << snap_->cgData.position[cg1->getLocalIndex()] << "\n";
861 cerr << "cg2_gid = " << cg2->getGlobalIndex() << "\tcg2_lid = " << cg2->getLocalIndex() << "\tcg2p = "
862 << snap_->cgData.position[cg2->getLocalIndex()] << "\n";*/
863
864 snap_->wrapVector(d);
865 return d;
866 }
867
868 Vector3d ForceMatrixDecomposition::getAtomToGroupVectorRow(int atom1, int cg1) {
869
870 Vector3d d;
871
872 #ifdef IS_MPI
873 d = cgRowData.position[cg1] - atomRowData.position[atom1];
874 #else
875 d = snap_->cgData.position[cg1] - snap_->atomData.position[atom1];
876 #endif
877
878 snap_->wrapVector(d);
879 return d;
880 }
881
882 Vector3d ForceMatrixDecomposition::getAtomToGroupVectorColumn(int atom2, int cg2) {
883 Vector3d d;
884
885 #ifdef IS_MPI
886 d = cgColData.position[cg2] - atomColData.position[atom2];
887 #else
888 d = snap_->cgData.position[cg2] - snap_->atomData.position[atom2];
889 #endif
890
891 snap_->wrapVector(d);
892 return d;
893 }
894
895 RealType ForceMatrixDecomposition::getMassFactorRow(int atom1) {
896 #ifdef IS_MPI
897 return massFactorsRow[atom1];
898 #else
899 return massFactors[atom1];
900 #endif
901 }
902
903 RealType ForceMatrixDecomposition::getMassFactorColumn(int atom2) {
904 #ifdef IS_MPI
905 return massFactorsCol[atom2];
906 #else
907 return massFactors[atom2];
908 #endif
909
910 }
911
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
921 snap_->wrapVector(d);
922 return d;
923 }
924
925 vector<int> ForceMatrixDecomposition::getExcludesForAtom(int atom1) {
926 return excludesForAtom[atom1];
927 }
928
929 /**
930 * We need to exclude some overcounted interactions that result from
931 * the parallel decomposition.
932 */
933 bool ForceMatrixDecomposition::skipAtomPair(int atom1, int atom2) {
934 int unique_id_1, unique_id_2;
935
936 // cerr << "sap with atom1, atom2 =\t" << atom1 << "\t" << atom2 << "\n";
937 #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
942 cerr << "sap with uid1, uid2 =\t" << unique_id_1 << "\t" << unique_id_2 << "\n";
943 // this situation should only arise in MPI simulations
944 if (unique_id_1 == unique_id_2) return true;
945
946 // this prevents us from doing the pair on multiple processors
947 if (unique_id_1 < unique_id_2)
948 {
949 if ((unique_id_1 + unique_id_2) % 2 == 0) return true;
950 } else
951 {
952 if ((unique_id_1 + unique_id_2) % 2 == 1) return true;
953 }
954 #endif
955 return false;
956 }
957
958 /**
959 * We need to handle the interactions for atoms who are involved in
960 * the same rigid body as well as some short range interactions
961 * (bonds, bends, torsions) differently from other interactions.
962 * We'll still visit the pairwise routines, but with a flag that
963 * tells those routines to exclude the pair from direct long range
964 * interactions. Some indirect interactions (notably reaction
965 * field) must still be handled for these pairs.
966 */
967 bool ForceMatrixDecomposition::excludeAtomPair(int atom1, int atom2) {
968 int unique_id_2;
969 #ifdef IS_MPI
970 // in MPI, we have to look up the unique IDs for the row atom.
971 unique_id_2 = AtomColToGlobal[atom2];
972 #else
973 // in the normal loop, the atom numbers are unique
974 unique_id_2 = atom2;
975 #endif
976
977 for (vector<int>::iterator i = excludesForAtom[atom1].begin(); i != excludesForAtom[atom1].end(); ++i)
978 {
979 if ((*i) == unique_id_2)
980 return true;
981 }
982
983 return false;
984 }
985
986 void ForceMatrixDecomposition::addForceToAtomRow(int atom1, Vector3d fg) {
987 #ifdef IS_MPI
988 atomRowData.force[atom1] += fg;
989 #else
990 snap_->atomData.force[atom1] += fg;
991 #endif
992 }
993
994 void ForceMatrixDecomposition::addForceToAtomRowOMP(int atom1, Vector3d fg) {
995 #pragma omp critical
996 {
997 snap_->atomData.force[atom1] += fg;
998 }
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 void ForceMatrixDecomposition::addForceToAtomColumnOMP(int atom2, Vector3d fg) {
1010 #pragma omp critical
1011 {
1012 snap_->atomData.force[atom2] += fg;
1013 }
1014 }
1015
1016 // filling interaction blocks with pointers
1017 void ForceMatrixDecomposition::fillInteractionData(InteractionData &idat, int atom1, int atom2) {
1018
1019 idat.excluded = excludeAtomPair(atom1, atom2);
1020
1021 #ifdef IS_MPI
1022 idat.atypes = make_pair( atypesRow[atom1], atypesCol[atom2]);
1023 //idat.atypes = make_pair( ff_->getAtomType(identsRow[atom1]),
1024 // ff_->getAtomType(identsCol[atom2]) );
1025
1026 if (storageLayout_ & DataStorage::dslAmat)
1027 {
1028 idat.A1 = &(atomRowData.aMat[atom1]);
1029 idat.A2 = &(atomColData.aMat[atom2]);
1030 }
1031
1032 if (storageLayout_ & DataStorage::dslElectroFrame)
1033 {
1034 idat.eFrame1 = &(atomRowData.electroFrame[atom1]);
1035 idat.eFrame2 = &(atomColData.electroFrame[atom2]);
1036 }
1037
1038 if (storageLayout_ & DataStorage::dslTorque)
1039 {
1040 idat.t1 = &(atomRowData.torque[atom1]);
1041 idat.t2 = &(atomColData.torque[atom2]);
1042 }
1043
1044 if (storageLayout_ & DataStorage::dslDensity)
1045 {
1046 idat.rho1 = &(atomRowData.density[atom1]);
1047 idat.rho2 = &(atomColData.density[atom2]);
1048 }
1049
1050 if (storageLayout_ & DataStorage::dslFunctional)
1051 {
1052 idat.frho1 = &(atomRowData.functional[atom1]);
1053 idat.frho2 = &(atomColData.functional[atom2]);
1054 }
1055
1056 if (storageLayout_ & DataStorage::dslFunctionalDerivative)
1057 {
1058 idat.dfrho1 = &(atomRowData.functionalDerivative[atom1]);
1059 idat.dfrho2 = &(atomColData.functionalDerivative[atom2]);
1060 }
1061
1062 if (storageLayout_ & DataStorage::dslParticlePot)
1063 {
1064 idat.particlePot1 = &(atomRowData.particlePot[atom1]);
1065 idat.particlePot2 = &(atomColData.particlePot[atom2]);
1066 }
1067
1068 if (storageLayout_ & DataStorage::dslSkippedCharge)
1069 {
1070 idat.skippedCharge1 = &(atomRowData.skippedCharge[atom1]);
1071 idat.skippedCharge2 = &(atomColData.skippedCharge[atom2]);
1072 }
1073
1074 #else
1075
1076 idat.atypes = make_pair(atypesLocal[atom1], atypesLocal[atom2]);
1077 //idat.atypes = make_pair( ff_->getAtomType(idents[atom1]),
1078 // ff_->getAtomType(idents[atom2]) );
1079
1080 if (storageLayout_ & DataStorage::dslAmat)
1081 {
1082 idat.A1 = &(snap_->atomData.aMat[atom1]);
1083 idat.A2 = &(snap_->atomData.aMat[atom2]);
1084 }
1085
1086 if (storageLayout_ & DataStorage::dslElectroFrame)
1087 {
1088 idat.eFrame1 = &(snap_->atomData.electroFrame[atom1]);
1089 idat.eFrame2 = &(snap_->atomData.electroFrame[atom2]);
1090 }
1091
1092 if (storageLayout_ & DataStorage::dslTorque)
1093 {
1094 idat.t1 = &(snap_->atomData.torque[atom1]);
1095 idat.t2 = &(snap_->atomData.torque[atom2]);
1096 }
1097
1098 if (storageLayout_ & DataStorage::dslDensity)
1099 {
1100 idat.rho1 = &(snap_->atomData.density[atom1]);
1101 idat.rho2 = &(snap_->atomData.density[atom2]);
1102 }
1103
1104 if (storageLayout_ & DataStorage::dslFunctional)
1105 {
1106 idat.frho1 = &(snap_->atomData.functional[atom1]);
1107 idat.frho2 = &(snap_->atomData.functional[atom2]);
1108 }
1109
1110 if (storageLayout_ & DataStorage::dslFunctionalDerivative)
1111 {
1112 idat.dfrho1 = &(snap_->atomData.functionalDerivative[atom1]);
1113 idat.dfrho2 = &(snap_->atomData.functionalDerivative[atom2]);
1114 }
1115
1116 if (storageLayout_ & DataStorage::dslParticlePot)
1117 {
1118 idat.particlePot1 = &(snap_->atomData.particlePot[atom1]);
1119 idat.particlePot2 = &(snap_->atomData.particlePot[atom2]);
1120 }
1121
1122 if (storageLayout_ & DataStorage::dslSkippedCharge)
1123 {
1124 idat.skippedCharge1 = &(snap_->atomData.skippedCharge[atom1]);
1125 idat.skippedCharge2 = &(snap_->atomData.skippedCharge[atom2]);
1126 }
1127 #endif
1128 }
1129
1130 // filling interaction blocks with pointers
1131 void ForceMatrixDecomposition::fillInteractionDataOMP(InteractionDataPrv &idat, int atom1, int atom2) {
1132
1133 idat.excluded = excludeAtomPair(atom1, atom2);
1134
1135 #ifdef IS_MPI
1136 idat.atypes = make_pair( atypesRow[atom1], atypesCol[atom2]);
1137 //idat.atypes = make_pair( ff_->getAtomType(identsRow[atom1]),
1138 // ff_->getAtomType(identsCol[atom2]) );
1139
1140 if (storageLayout_ & DataStorage::dslAmat)
1141 {
1142 idat.A1 = &(atomRowData.aMat[atom1]);
1143 idat.A2 = &(atomColData.aMat[atom2]);
1144 }
1145
1146 if (storageLayout_ & DataStorage::dslElectroFrame)
1147 {
1148 idat.eFrame1 = &(atomRowData.electroFrame[atom1]);
1149 idat.eFrame2 = &(atomColData.electroFrame[atom2]);
1150 }
1151
1152 if (storageLayout_ & DataStorage::dslTorque)
1153 {
1154 idat.t1 = &(atomRowData.torque[atom1]);
1155 idat.t2 = &(atomColData.torque[atom2]);
1156 }
1157
1158 if (storageLayout_ & DataStorage::dslDensity)
1159 {
1160 idat.rho1 = &(atomRowData.density[atom1]);
1161 idat.rho2 = &(atomColData.density[atom2]);
1162 }
1163
1164 if (storageLayout_ & DataStorage::dslFunctional)
1165 {
1166 idat.frho1 = &(atomRowData.functional[atom1]);
1167 idat.frho2 = &(atomColData.functional[atom2]);
1168 }
1169
1170 if (storageLayout_ & DataStorage::dslFunctionalDerivative)
1171 {
1172 idat.dfrho1 = &(atomRowData.functionalDerivative[atom1]);
1173 idat.dfrho2 = &(atomColData.functionalDerivative[atom2]);
1174 }
1175
1176 if (storageLayout_ & DataStorage::dslParticlePot)
1177 {
1178 idat.particlePot1 = &(atomRowData.particlePot[atom1]);
1179 idat.particlePot2 = &(atomColData.particlePot[atom2]);
1180 }
1181
1182 if (storageLayout_ & DataStorage::dslSkippedCharge)
1183 {
1184 idat.skippedCharge1 = &(atomRowData.skippedCharge[atom1]);
1185 idat.skippedCharge2 = &(atomColData.skippedCharge[atom2]);
1186 }
1187
1188 #else
1189
1190 idat.atypes = make_pair(atypesLocal[atom1], atypesLocal[atom2]);
1191 //idat.atypes = make_pair( ff_->getAtomType(idents[atom1]),
1192 // ff_->getAtomType(idents[atom2]) );
1193
1194 if (storageLayout_ & DataStorage::dslAmat)
1195 {
1196 idat.A1 = &(snap_->atomData.aMat[atom1]);
1197 idat.A2 = &(snap_->atomData.aMat[atom2]);
1198 }
1199
1200 if (storageLayout_ & DataStorage::dslElectroFrame)
1201 {
1202 idat.eFrame1 = &(snap_->atomData.electroFrame[atom1]);
1203 idat.eFrame2 = &(snap_->atomData.electroFrame[atom2]);
1204 }
1205
1206 if (storageLayout_ & DataStorage::dslTorque)
1207 {
1208 idat.t1 = &(snap_->atomData.torque[atom1]);
1209 idat.t2 = &(snap_->atomData.torque[atom2]);
1210 }
1211
1212 if (storageLayout_ & DataStorage::dslDensity)
1213 {
1214 idat.rho1 = &(snap_->atomData.density[atom1]);
1215 idat.rho2 = &(snap_->atomData.density[atom2]);
1216 }
1217
1218 if (storageLayout_ & DataStorage::dslFunctional)
1219 {
1220 idat.frho1 = &(snap_->atomData.functional[atom1]);
1221 idat.frho2 = &(snap_->atomData.functional[atom2]);
1222 }
1223
1224 if (storageLayout_ & DataStorage::dslFunctionalDerivative)
1225 {
1226 idat.dfrho1 = &(snap_->atomData.functionalDerivative[atom1]);
1227 idat.dfrho2 = &(snap_->atomData.functionalDerivative[atom2]);
1228 }
1229
1230 if (storageLayout_ & DataStorage::dslParticlePot)
1231 {
1232 idat.particlePot1 = &(snap_->atomData.particlePot[atom1]);
1233 idat.particlePot2 = &(snap_->atomData.particlePot[atom2]);
1234 }
1235
1236 if (storageLayout_ & DataStorage::dslSkippedCharge)
1237 {
1238 idat.skippedCharge1 = &(snap_->atomData.skippedCharge[atom1]);
1239 idat.skippedCharge2 = &(snap_->atomData.skippedCharge[atom2]);
1240 }
1241 #endif
1242 }
1243
1244 void ForceMatrixDecomposition::unpackInteractionData(InteractionData &idat, int atom1, int atom2) {
1245 #ifdef IS_MPI
1246 pot_row[atom1] += 0.5 * *(idat.pot);
1247 pot_col[atom2] += 0.5 * *(idat.pot);
1248
1249 atomRowData.force[atom1] += *(idat.f1);
1250 atomColData.force[atom2] -= *(idat.f1);
1251 #else
1252 pairwisePot += *(idat.pot);
1253
1254 snap_->atomData.force[atom1] += *(idat.f1);
1255 snap_->atomData.force[atom2] -= *(idat.f1);
1256 #endif
1257
1258 }
1259
1260 void ForceMatrixDecomposition::unpackInteractionDataOMP(InteractionDataPrv &idat, int atom1, int atom2) {
1261 #pragma omp critical
1262 {
1263 pairwisePot += idat.pot;
1264
1265 snap_->atomData.force[atom1] += idat.f1;
1266 snap_->atomData.force[atom2] -= idat.f1;
1267 }
1268 }
1269
1270 void ForceMatrixDecomposition::reorderGroupCutoffs(vector<int> &order) {
1271 vector<int> tmp = vector<int> (groupToGtype.size());
1272
1273 for (int i = 0; i < groupToGtype.size(); ++i)
1274 {
1275 tmp[i] = groupToGtype[i];
1276 }
1277
1278 for (int i = 0; i < groupToGtype.size(); ++i)
1279 {
1280 groupToGtype[i] = tmp[order[i]];
1281 }
1282 }
1283
1284 void ForceMatrixDecomposition::reorderPosition(vector<int> &order) {
1285 Snapshot* snap_ = info_->getSnapshotManager()->getCurrentSnapshot();
1286 DataStorage* cgConfig = &(snap_->cgData);
1287 vector<Vector3d> tmp = vector<Vector3d> (nGroups_);
1288
1289 for (int i = 0; i < nGroups_; ++i)
1290 {
1291 tmp[i] = snap_->cgData.position[i];
1292 }
1293
1294 vector<int> mapPos = vector<int> (nGroups_);
1295 for (int i = 0; i < nGroups_; ++i)
1296 {
1297 snap_->cgData.position[i] = tmp[order[i]];
1298 mapPos[order[i]] = i;
1299 }
1300
1301 SimInfo::MoleculeIterator mi;
1302 Molecule* mol;
1303 Molecule::CutoffGroupIterator ci;
1304 CutoffGroup* cg;
1305
1306 for (mol = info_->beginMolecule(mi); mol != NULL; mol = info_->nextMolecule(mi))
1307 {
1308 for (cg = mol->beginCutoffGroup(ci); cg != NULL; cg = mol->nextCutoffGroup(ci))
1309 {
1310 cg->setLocalIndex(mapPos[cg->getLocalIndex()]);
1311 }
1312 }
1313
1314 /* if (info_->getNCutoffGroups() > 0)
1315 {
1316 for (mol = info_->beginMolecule(mi); mol != NULL; mol = info_->nextMolecule(mi))
1317 {
1318 for (cg = mol->beginCutoffGroup(ci); cg != NULL; cg = mol->nextCutoffGroup(ci))
1319 {
1320 printf("gbI:%d locI:%d x:%f y:%f z:%f\n", cg->getGlobalIndex(), cg->getLocalIndex(),
1321 cgConfig->position[cg->getLocalIndex()].x(), cgConfig->position[cg->getLocalIndex()].y(),
1322 cgConfig->position[cg->getLocalIndex()].z());
1323 }
1324 }
1325 } else
1326 {
1327 // center of mass of the group is the same as position of the atom
1328 // if cutoff group does not exist
1329 printf("ERROR!!!!!!!!!!!!!!!!!!!!!!!!!!!\n");
1330 // cgConfig->position = config->position;
1331 }*/
1332 }
1333
1334 void ForceMatrixDecomposition::reorderGroupList(vector<int> &order) {
1335 vector<vector<int> > tmp = vector<vector<int> > (groupList_.size());
1336
1337 for (int i = 0; i < groupList_.size(); ++i)
1338 {
1339 tmp[i] = groupList_[i];
1340 }
1341
1342 for (int i = 0; i < groupList_.size(); ++i)
1343 {
1344 groupList_[i] = tmp[order[i]];
1345 }
1346 }
1347
1348 void ForceMatrixDecomposition::reorderMemory(vector<vector<CutoffGroup *> > &H_c_l) {
1349 int n = 0;
1350 // printf("Reorder memory time:%f!!!!!!!!!!!!!!!!!!!!!!!!!!!\n",
1351 // info_->getSnapshotManager()->getCurrentSnapshot()->getTime());
1352
1353 /* record the reordered atom indices */
1354 vector<int> k = vector<int> (nGroups_);
1355
1356 for (int c = 0; c < H_c_l.size(); ++c)
1357 {
1358 for (vector<CutoffGroup *>::iterator cg = H_c_l[c].begin(); cg != H_c_l[c].end(); ++cg)
1359 {
1360 int i = (*cg)->getGlobalIndex();
1361 k[n] = i;
1362 ++n;
1363 }
1364 }
1365
1366 // reorderGroupCutoffs(k);
1367 // reorderGroupList(k);
1368 reorderPosition(k);
1369 }
1370
1371 vector<vector<CutoffGroup *> > ForceMatrixDecomposition::buildLayerBasedNeighborList() {
1372 // printf("buildLayerBasedNeighborList; nGroups:%d\n", nGroups_);
1373 // Na = nGroups_
1374 /* cell occupancy counter */
1375 // vector<int> k_c;
1376 /* c_i - has cell containing atom i (size Na) */
1377 vector<int> c = vector<int> (nGroups_);
1378 /* l_i - layer containing atom i (size Na) */
1379 // vector<int> l;
1380
1381 RealType rList_ = (largestRcut_ + skinThickness_);
1382 Snapshot* snap_ = sman_->getCurrentSnapshot();
1383 Mat3x3d Hmat = snap_->getHmat();
1384 Vector3d Hx = Hmat.getColumn(0);
1385 Vector3d Hy = Hmat.getColumn(1);
1386 Vector3d Hz = Hmat.getColumn(2);
1387
1388 nCells_.x() = (int) (Hx.length()) / rList_;
1389 nCells_.y() = (int) (Hy.length()) / rList_;
1390 nCells_.z() = (int) (Hz.length()) / rList_;
1391
1392 Mat3x3d invHmat = snap_->getInvHmat();
1393 Vector3d rs, scaled, dr;
1394 Vector3i whichCell;
1395 int cellIndex;
1396 int nCtot = nCells_.x() * nCells_.y() * nCells_.z();
1397
1398 // k_c = vector<int> (nCtot, 0);
1399
1400 SimInfo::MoleculeIterator mi;
1401 Molecule* mol;
1402 Molecule::CutoffGroupIterator ci;
1403 CutoffGroup* cg;
1404
1405 for (mol = info_->beginMolecule(mi); mol != NULL; mol = info_->nextMolecule(mi))
1406 {
1407 for (cg = mol->beginCutoffGroup(ci); cg != NULL; cg = mol->nextCutoffGroup(ci))
1408 {
1409 rs = snap_->cgData.position[cg->getLocalIndex()];
1410
1411 // scaled positions relative to the box vectors
1412 scaled = invHmat * rs;
1413
1414 // wrap the vector back into the unit box by subtracting integer box
1415 // numbers
1416 for (int j = 0; j < 3; j++)
1417 {
1418 scaled[j] -= roundMe(scaled[j]);
1419 scaled[j] += 0.5;
1420 }
1421
1422 // find xyz-indices of cell that cutoffGroup is in.
1423 whichCell.x() = nCells_.x() * scaled.x();
1424 whichCell.y() = nCells_.y() * scaled.y();
1425 whichCell.z() = nCells_.z() * scaled.z();
1426
1427 // printf("pos x:%f y:%f z:%f cell x:%d y:%d z:%d\n", rs.x(), rs.y(), rs.z(), whichCell.x(), whichCell.y(),
1428 // whichCell.z());
1429
1430 // find single index of this cell:
1431 cellIndex = Vlinear(whichCell, nCells_);
1432
1433 c[cg->getGlobalIndex()] = cellIndex;
1434 }
1435 }
1436
1437 // int k_c_curr;
1438 // int k_c_max = 0;
1439 /* the cell-layer occupancy matrix */
1440 vector<vector<CutoffGroup *> > H_c_l = vector<vector<CutoffGroup *> > (nCtot);
1441
1442 for (mol = info_->beginMolecule(mi); mol != NULL; mol = info_->nextMolecule(mi))
1443 {
1444 for (cg = mol->beginCutoffGroup(ci); cg != NULL; cg = mol->nextCutoffGroup(ci))
1445
1446 {
1447 // k_c_curr = ++k_c[c[cg1->getGlobalIndex()]];
1448 // l.push_back(k_c_curr);
1449 //
1450 // /* determines the number of layers in use */
1451 // if (k_c_max < k_c_curr)
1452 // {
1453 // k_c_max = k_c_curr;
1454 // }
1455 H_c_l[c[cg->getGlobalIndex()]].push_back(/*l[*/cg/*]*/);
1456 }
1457 }
1458
1459 /* Frequency of reordering the memory */
1460 if (neighborListReorderFreq != 0)
1461 {
1462 if (reorderFreqCounter == neighborListReorderFreq)
1463 {
1464 //printf("neighborListReorderFreq:%d\n", neighborListReorderFreq);
1465 reorderMemory(H_c_l);
1466 reorderFreqCounter = 1;
1467 } else
1468 {
1469 reorderFreqCounter++;
1470 }
1471 }
1472
1473 int m;
1474 /* the neighbor matrix */
1475 vector<vector<CutoffGroup *> > neighborMatW = vector<vector<CutoffGroup *> > (nGroups_);
1476
1477 groupCutoffs cuts;
1478 CutoffGroup *cg1;
1479
1480 /* loops over objects(atoms, rigidBodies, cutoffGroups, etc.) */
1481 for (mol = info_->beginMolecule(mi); mol != NULL; mol = info_->nextMolecule(mi))
1482 {
1483 for (cg1 = mol->beginCutoffGroup(ci); cg1 != NULL; cg1 = mol->nextCutoffGroup(ci))
1484 {
1485 /* c' */
1486 int c1 = c[cg1->getGlobalIndex()];
1487 Vector3i c1v = idxToV(c1, nCells_);
1488
1489 /* loops over the neighboring cells c'' */
1490 for (vector<Vector3i>::iterator os = cellOffsets_.begin(); os != cellOffsets_.end(); ++os)
1491 {
1492 Vector3i c2v = c1v + (*os);
1493
1494 if (c2v.x() >= nCells_.x())
1495 {
1496 c2v.x() = 0;
1497 } else if (c2v.x() < 0)
1498 {
1499 c2v.x() = nCells_.x() - 1;
1500 }
1501
1502 if (c2v.y() >= nCells_.y())
1503 {
1504 c2v.y() = 0;
1505 } else if (c2v.y() < 0)
1506 {
1507 c2v.y() = nCells_.y() - 1;
1508 }
1509
1510 if (c2v.z() >= nCells_.z())
1511 {
1512 c2v.z() = 0;
1513 } else if (c2v.z() < 0)
1514 {
1515 c2v.z() = nCells_.z() - 1;
1516 }
1517
1518 int c2 = Vlinear(c2v, nCells_);
1519 /* loops over layers l to access the neighbor atoms */
1520 for (vector<CutoffGroup *>::iterator cg2 = H_c_l[c2].begin(); cg2 != H_c_l[c2].end(); ++cg2)
1521 {
1522 // if i'' = 0 then break // doesn't apply to vector implementation of matrix
1523 // if(i != *j)
1524 if (c2 != c1 || (*cg2)->getGlobalIndex() < cg1->getGlobalIndex())
1525 {
1526 dr = snap_->cgData.position[(*cg2)->getLocalIndex()] - snap_->cgData.position[cg1->getLocalIndex()];
1527 snap_->wrapVector(dr);
1528 cuts = getGroupCutoffs(cg1->getGlobalIndex(), (*cg2)->getGlobalIndex());
1529 if (dr.lengthSquare() < cuts.third)
1530 {
1531 /* transposed version of Rapaport W mat, to occupy successive memory locations on CPU */
1532 neighborMatW[cg1->getGlobalIndex()].push_back((*cg2));
1533 }
1534 }
1535 }
1536 }
1537 }
1538 }
1539
1540 // save the local cutoff group positions for the check that is
1541 // done on each loop:
1542 saved_CG_positions_.clear();
1543 for (int i = 0; i < nGroups_; i++)
1544 saved_CG_positions_.push_back(snap_->cgData.position[i]);
1545
1546 return neighborMatW;
1547 }
1548
1549 /*
1550 * buildNeighborList
1551 *
1552 * first element of pair is row-indexed CutoffGroup
1553 * second element of pair is column-indexed CutoffGroup
1554 */
1555 vector<pair<int, int> > ForceMatrixDecomposition::buildNeighborList() {
1556
1557 vector<pair<int, int> > neighborList;
1558 groupCutoffs cuts;
1559 bool doAllPairs = false;
1560
1561 #ifdef IS_MPI
1562 cellListRow_.clear();
1563 cellListCol_.clear();
1564 #else
1565 cellList_.clear();
1566 #endif
1567
1568 RealType rList_ = (largestRcut_ + skinThickness_);
1569 RealType rl2 = rList_ * rList_;
1570 Snapshot* snap_ = sman_->getCurrentSnapshot();
1571 Mat3x3d Hmat = snap_->getHmat();
1572 Vector3d Hx = Hmat.getColumn(0);
1573 Vector3d Hy = Hmat.getColumn(1);
1574 Vector3d Hz = Hmat.getColumn(2);
1575
1576 nCells_.x() = (int) (Hx.length()) / rList_;
1577 nCells_.y() = (int) (Hy.length()) / rList_;
1578 nCells_.z() = (int) (Hz.length()) / rList_;
1579
1580 // handle small boxes where the cell offsets can end up repeating cells
1581
1582 if (nCells_.x() < 3)
1583 doAllPairs = true;
1584 if (nCells_.y() < 3)
1585 doAllPairs = true;
1586 if (nCells_.z() < 3)
1587 doAllPairs = true;
1588
1589 Mat3x3d invHmat = snap_->getInvHmat();
1590 Vector3d rs, scaled, dr;
1591 Vector3i whichCell;
1592 int cellIndex;
1593 int nCtot = nCells_.x() * nCells_.y() * nCells_.z();
1594
1595 #ifdef IS_MPI
1596 cellListRow_.resize(nCtot);
1597 cellListCol_.resize(nCtot);
1598 #else
1599 cellList_.resize(nCtot);
1600 #endif
1601
1602 if (!doAllPairs)
1603 {
1604 #ifdef IS_MPI
1605
1606 for (int i = 0; i < nGroupsInRow_; i++)
1607 {
1608 rs = cgRowData.position[i];
1609
1610 // scaled positions relative to the box vectors
1611 scaled = invHmat * rs;
1612
1613 // wrap the vector back into the unit box by subtracting integer box
1614 // numbers
1615 for (int j = 0; j < 3; j++)
1616 {
1617 scaled[j] -= roundMe(scaled[j]);
1618 scaled[j] += 0.5;
1619 }
1620
1621 // find xyz-indices of cell that cutoffGroup is in.
1622 whichCell.x() = nCells_.x() * scaled.x();
1623 whichCell.y() = nCells_.y() * scaled.y();
1624 whichCell.z() = nCells_.z() * scaled.z();
1625
1626 // find single index of this cell:
1627 cellIndex = Vlinear(whichCell, nCells_);
1628
1629 // add this cutoff group to the list of groups in this cell;
1630 cellListRow_[cellIndex].push_back(i);
1631 }
1632 for (int i = 0; i < nGroupsInCol_; i++)
1633 {
1634 rs = cgColData.position[i];
1635
1636 // scaled positions relative to the box vectors
1637 scaled = invHmat * rs;
1638
1639 // wrap the vector back into the unit box by subtracting integer box
1640 // numbers
1641 for (int j = 0; j < 3; j++)
1642 {
1643 scaled[j] -= roundMe(scaled[j]);
1644 scaled[j] += 0.5;
1645 }
1646
1647 // find xyz-indices of cell that cutoffGroup is in.
1648 whichCell.x() = nCells_.x() * scaled.x();
1649 whichCell.y() = nCells_.y() * scaled.y();
1650 whichCell.z() = nCells_.z() * scaled.z();
1651
1652 // find single index of this cell:
1653 cellIndex = Vlinear(whichCell, nCells_);
1654
1655 // add this cutoff group to the list of groups in this cell;
1656 cellListCol_[cellIndex].push_back(i);
1657 }
1658 #else
1659 for (int i = 0; i < nGroups_; i++)
1660 {
1661 rs = snap_->cgData.position[i];
1662
1663 // scaled positions relative to the box vectors
1664 scaled = invHmat * rs;
1665
1666 // wrap the vector back into the unit box by subtracting integer box
1667 // numbers
1668 for (int j = 0; j < 3; j++)
1669 {
1670 scaled[j] -= roundMe(scaled[j]);
1671 scaled[j] += 0.5;
1672 }
1673
1674 // find xyz-indices of cell that cutoffGroup is in.
1675 whichCell.x() = nCells_.x() * scaled.x();
1676 whichCell.y() = nCells_.y() * scaled.y();
1677 whichCell.z() = nCells_.z() * scaled.z();
1678
1679 // find single index of this cell:
1680 cellIndex = Vlinear(whichCell, nCells_);
1681
1682 // add this cutoff group to the list of groups in this cell;
1683 cellList_[cellIndex].push_back(i);
1684 }
1685 #endif
1686
1687 for (int m1z = 0; m1z < nCells_.z(); m1z++)
1688 {
1689 for (int m1y = 0; m1y < nCells_.y(); m1y++)
1690 {
1691 for (int m1x = 0; m1x < nCells_.x(); m1x++)
1692 {
1693 Vector3i m1v(m1x, m1y, m1z);
1694 int m1 = Vlinear(m1v, nCells_);
1695
1696 for (vector<Vector3i>::iterator os = cellOffsets_.begin(); os != cellOffsets_.end(); ++os)
1697 {
1698
1699 Vector3i m2v = m1v + (*os);
1700
1701 if (m2v.x() >= nCells_.x())
1702 {
1703 m2v.x() = 0;
1704 } else if (m2v.x() < 0)
1705 {
1706 m2v.x() = nCells_.x() - 1;
1707 }
1708
1709 if (m2v.y() >= nCells_.y())
1710 {
1711 m2v.y() = 0;
1712 } else if (m2v.y() < 0)
1713 {
1714 m2v.y() = nCells_.y() - 1;
1715 }
1716
1717 if (m2v.z() >= nCells_.z())
1718 {
1719 m2v.z() = 0;
1720 } else if (m2v.z() < 0)
1721 {
1722 m2v.z() = nCells_.z() - 1;
1723 }
1724
1725 int m2 = Vlinear(m2v, nCells_);
1726
1727 #ifdef IS_MPI
1728 for (vector<int>::iterator j1 = cellListRow_[m1].begin();
1729 j1 != cellListRow_[m1].end(); ++j1)
1730 {
1731 for (vector<int>::iterator j2 = cellListCol_[m2].begin();
1732 j2 != cellListCol_[m2].end(); ++j2)
1733 {
1734
1735 // In parallel, we need to visit *all* pairs of row &
1736 // column indicies and will truncate later on.
1737 dr = cgColData.position[(*j2)] - cgRowData.position[(*j1)];
1738 snap_->wrapVector(dr);
1739 cuts = getGroupCutoffs( (*j1), (*j2) );
1740 if (dr.lengthSquare() < cuts.third)
1741 {
1742 neighborList.push_back(make_pair((*j1), (*j2)));
1743 }
1744 }
1745 }
1746 #else
1747
1748 for (vector<int>::iterator j1 = cellList_[m1].begin(); j1 != cellList_[m1].end(); ++j1)
1749 {
1750 for (vector<int>::iterator j2 = cellList_[m2].begin(); j2 != cellList_[m2].end(); ++j2)
1751 {
1752
1753 // Always do this if we're in different cells or if
1754 // we're in the same cell and the global index of the
1755 // j2 cutoff group is less than the j1 cutoff group
1756
1757 if (m2 != m1 || (*j2) < (*j1))
1758 {
1759 dr = snap_->cgData.position[(*j2)] - snap_->cgData.position[(*j1)];
1760 snap_->wrapVector(dr);
1761 cuts = getGroupCutoffs((*j1), (*j2));
1762 if (dr.lengthSquare() < cuts.third)
1763 {
1764 neighborList.push_back(make_pair((*j1), (*j2)));
1765 }
1766 }
1767 }
1768 }
1769 #endif
1770 }
1771 }
1772 }
1773 }
1774 } else
1775 {
1776 // branch to do all cutoff group pairs
1777 #ifdef IS_MPI
1778 for (int j1 = 0; j1 < nGroupsInRow_; j1++)
1779 {
1780 for (int j2 = 0; j2 < nGroupsInCol_; j2++)
1781 {
1782 dr = cgColData.position[j2] - cgRowData.position[j1];
1783 snap_->wrapVector(dr);
1784 cuts = getGroupCutoffs( j1, j2 );
1785 if (dr.lengthSquare() < cuts.third)
1786 {
1787 neighborList.push_back(make_pair(j1, j2));
1788 }
1789 }
1790 }
1791 #else
1792 for (int j1 = 0; j1 < nGroups_ - 1; j1++)
1793 {
1794 for (int j2 = j1 + 1; j2 < nGroups_; j2++)
1795 {
1796 dr = snap_->cgData.position[j2] - snap_->cgData.position[j1];
1797 snap_->wrapVector(dr);
1798 cuts = getGroupCutoffs(j1, j2);
1799 if (dr.lengthSquare() < cuts.third)
1800 {
1801 neighborList.push_back(make_pair(j1, j2));
1802 }
1803 }
1804 }
1805 #endif
1806 }
1807
1808 // save the local cutoff group positions for the check that is
1809 // done on each loop:
1810 saved_CG_positions_.clear();
1811 for (int i = 0; i < nGroups_; i++)
1812 saved_CG_positions_.push_back(snap_->cgData.position[i]);
1813
1814 return neighborList;
1815 }
1816 } //end namespace OpenMD