ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/parallel/ForceDecomposition.cpp
Revision: 1575
Committed: Fri Jun 3 21:39:49 2011 UTC (13 years, 10 months ago) by gezelter
File size: 5536 byte(s)
Log Message:
more parallel fixes

File Contents

# Content
1 /*
2 * Copyright (c) 2005 The University of Notre Dame. All Rights Reserved.
3 *
4 * The University of Notre Dame grants you ("Licensee") a
5 * non-exclusive, royalty free, license to use, modify and
6 * redistribute this software in source and binary code form, provided
7 * that the following conditions are met:
8 *
9 * 1. Redistributions of source code must retain the above copyright
10 * notice, this list of conditions and the following disclaimer.
11 *
12 * 2. Redistributions in binary form must reproduce the above copyright
13 * notice, this list of conditions and the following disclaimer in the
14 * documentation and/or other materials provided with the
15 * distribution.
16 *
17 * This software is provided "AS IS," without a warranty of any
18 * kind. All express or implied conditions, representations and
19 * warranties, including any implied warranty of merchantability,
20 * fitness for a particular purpose or non-infringement, are hereby
21 * excluded. The University of Notre Dame and its licensors shall not
22 * be liable for any damages suffered by licensee as a result of
23 * using, modifying or distributing the software or its
24 * derivatives. In no event will the University of Notre Dame or its
25 * licensors be liable for any lost revenue, profit or data, or for
26 * direct, indirect, special, consequential, incidental or punitive
27 * damages, however caused and regardless of the theory of liability,
28 * arising out of the use of or inability to use software, even if the
29 * University of Notre Dame has been advised of the possibility of
30 * such damages.
31 *
32 * SUPPORT OPEN SCIENCE! If you use OpenMD or its source code in your
33 * research, please cite the appropriate papers when you publish your
34 * work. Good starting points are:
35 *
36 * [1] Meineke, et al., J. Comp. Chem. 26, 252-271 (2005).
37 * [2] Fennell & Gezelter, J. Chem. Phys. 124, 234104 (2006).
38 * [3] Sun, Lin & Gezelter, J. Chem. Phys. 128, 24107 (2008).
39 * [4] Vardeman & Gezelter, in progress (2009).
40 */
41 #include "parallel/ForceDecomposition.hpp"
42 #include "math/Vector3.hpp"
43 #ifdef IS_MPI
44 #include <mpi.h>
45 #endif
46
47 using namespace std;
48 namespace OpenMD {
49
50 ForceDecomposition::ForceDecomposition(SimInfo* info) : info_(info) {
51 sman_ = info_->getSnapshotManager();
52 storageLayout_ = sman_->getStorageLayout();
53 ff_ = info_->getForceField();
54
55 Globals* simParams_ = info_->getSimParams();
56
57 if (simParams_->haveSkinThickness()) {
58 skinThickness_ = simParams_->getSkinThickness();
59 } else {
60 skinThickness_ = 1.0;
61 sprintf(painCave.errMsg,
62 "ForceDecomposition: No value was set for the skinThickness.\n"
63 "\tOpenMD will use a default value of %f Angstroms\n"
64 "\tfor this simulation\n", skinThickness_);
65 painCave.severity = OPENMD_INFO;
66 painCave.isFatal = 0;
67 simError();
68 }
69
70 // cellOffsets are the partial space for the cell lists used in
71 // constructing the neighbor lists
72 cellOffsets_.clear();
73 cellOffsets_.push_back( Vector3i(0, 0, 0) );
74 cellOffsets_.push_back( Vector3i(1, 0, 0) );
75 cellOffsets_.push_back( Vector3i(1, 1, 0) );
76 cellOffsets_.push_back( Vector3i(0, 1, 0) );
77 cellOffsets_.push_back( Vector3i(-1,1, 0) );
78 cellOffsets_.push_back( Vector3i(0, 0, 1) );
79 cellOffsets_.push_back( Vector3i(1, 0, 1) );
80 cellOffsets_.push_back( Vector3i(1, 1, 1) );
81 cellOffsets_.push_back( Vector3i(0, 1, 1) );
82 cellOffsets_.push_back( Vector3i(-1,1, 1) );
83 cellOffsets_.push_back( Vector3i(-1,0, 1) );
84 cellOffsets_.push_back( Vector3i(-1,-1,1) );
85 cellOffsets_.push_back( Vector3i(0, -1,1) );
86 cellOffsets_.push_back( Vector3i(1, -1,1) );
87 }
88
89 SelfData ForceDecomposition::fillSelfData(int atom1) {
90 SelfData sdat;
91 // Still Missing atype, skippedCharge, potVec pot,
92 if (storageLayout_ & DataStorage::dslElectroFrame) {
93 sdat.eFrame = &(snap_->atomData.electroFrame[atom1]);
94 }
95
96 if (storageLayout_ & DataStorage::dslTorque) {
97 sdat.t = &(snap_->atomData.torque[atom1]);
98 }
99
100 if (storageLayout_ & DataStorage::dslDensity) {
101 sdat.rho = &(snap_->atomData.density[atom1]);
102 }
103
104 if (storageLayout_ & DataStorage::dslFunctional) {
105 sdat.frho = &(snap_->atomData.functional[atom1]);
106 }
107
108 if (storageLayout_ & DataStorage::dslFunctionalDerivative) {
109 sdat.dfrhodrho = &(snap_->atomData.functionalDerivative[atom1]);
110 }
111
112 if (storageLayout_ & DataStorage::dslParticlePot) {
113 sdat.particlePot = &(snap_->atomData.particlePot[atom1]);
114 }
115
116 return sdat;
117 }
118
119 bool ForceDecomposition::checkNeighborList() {
120
121 int nGroups = snap_->cgData.position.size();
122
123 // if we have changed the group identities or haven't set up the
124 // saved positions we automatically will need a neighbor list update:
125
126 if ( saved_CG_positions_.size() != nGroups ) return true;
127
128 RealType dispmax = 0.0;
129 Vector3d disp;
130
131 for (int i = 0; i < nGroups; i++) {
132 disp = snap_->cgData.position[i] - saved_CG_positions_[i];
133 for (int j = 0; j < 3; j++)
134 dispmax = max( abs(disp[j]), dispmax);
135 }
136
137 #ifdef IS_MPI
138 MPI::COMM_WORLD.Allreduce(&dispmax, &dispmax, 1, MPI::REALTYPE, MPI::MAX);
139 #endif
140
141 // a conservative test of list skin crossings
142 dispmax = 2.0 * sqrt (3.0 * dispmax * dispmax);
143
144 return (dispmax > skinThickness_);
145 }
146 }

Properties

Name Value
svn:eol-style native