ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/devel_omp/src/parallel/ForceMatrixDecomposition.cpp
Revision: 1598
Committed: Wed Jul 27 14:26:53 2011 UTC (13 years, 9 months ago) by mciznick
File size: 45832 byte(s)
Log Message:
Updated ordered neighbor list generation.

File Contents

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