ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/parallel/ForceMatrixDecomposition.cpp
Revision: 1576
Committed: Wed Jun 8 16:05:07 2011 UTC (13 years, 10 months ago) by gezelter
Original Path: branches/development/src/parallel/ForceMatrixDecomposition.cpp
File size: 36860 byte(s)
Log Message:
migrating cutoff information from InteractionManager to ForceManager

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