ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/parallel/ForceDecomposition.cpp
Revision: 1879
Committed: Sun Jun 16 15:15:42 2013 UTC (11 years, 10 months ago) by gezelter
File size: 6887 byte(s)
Log Message:
MERGE OpenMD development 1783:1878 into trunk

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, 234107 (2008).
39 * [4] Kuang & Gezelter, J. Chem. Phys. 133, 164101 (2010).
40 * [5] Vardeman, Stocker & Gezelter, J. Chem. Theory Comput. 7, 834 (2011).
41 */
42 #include "parallel/ForceDecomposition.hpp"
43 #include "math/Vector3.hpp"
44 #ifdef IS_MPI
45 #include <mpi.h>
46 #endif
47
48 using namespace std;
49 namespace OpenMD {
50
51 ForceDecomposition::ForceDecomposition(SimInfo* info, InteractionManager* iMan) : info_(info), interactionMan_(iMan), needVelocities_(false) {
52
53 sman_ = info_->getSnapshotManager();
54 storageLayout_ = sman_->getStorageLayout();
55 ff_ = info_->getForceField();
56 userChoseCutoff_ = false;
57
58 usePeriodicBoundaryConditions_ = info->getSimParams()->getUsePeriodicBoundaryConditions();
59
60 Globals* simParams_ = info_->getSimParams();
61 if (simParams_->havePrintHeatFlux()) {
62 if (simParams_->getPrintHeatFlux()) {
63 needVelocities_ = true;
64 }
65 }
66
67 if (simParams_->haveSkinThickness()) {
68 skinThickness_ = simParams_->getSkinThickness();
69 } else {
70 skinThickness_ = 1.0;
71 sprintf(painCave.errMsg,
72 "ForceDecomposition: No value was set for the skinThickness.\n"
73 "\tOpenMD will use a default value of %f Angstroms\n"
74 "\tfor this simulation\n", skinThickness_);
75 painCave.severity = OPENMD_INFO;
76 painCave.isFatal = 0;
77 simError();
78 }
79
80 // cellOffsets are the partial space for the cell lists used in
81 // constructing the neighbor lists
82 cellOffsets_.clear();
83 cellOffsets_.push_back( Vector3i(0, 0, 0) );
84 cellOffsets_.push_back( Vector3i(1, 0, 0) );
85 cellOffsets_.push_back( Vector3i(1, 1, 0) );
86 cellOffsets_.push_back( Vector3i(0, 1, 0) );
87 cellOffsets_.push_back( Vector3i(-1,1, 0) );
88 cellOffsets_.push_back( Vector3i(0, 0, 1) );
89 cellOffsets_.push_back( Vector3i(1, 0, 1) );
90 cellOffsets_.push_back( Vector3i(1, 1, 1) );
91 cellOffsets_.push_back( Vector3i(0, 1, 1) );
92 cellOffsets_.push_back( Vector3i(-1,1, 1) );
93 cellOffsets_.push_back( Vector3i(-1,0, 1) );
94 cellOffsets_.push_back( Vector3i(-1,-1,1) );
95 cellOffsets_.push_back( Vector3i(0, -1,1) );
96 cellOffsets_.push_back( Vector3i(1, -1,1) );
97 }
98
99 void ForceDecomposition::fillSelfData(SelfData &sdat, int atom1) {
100
101 sdat.atype = atypesLocal[atom1];
102
103 sdat.pot = &embeddingPot;
104 sdat.excludedPot = &excludedSelfPot;
105
106 if (storageLayout_ & DataStorage::dslDipole) {
107 sdat.dipole = &(snap_->atomData.dipole[atom1]);
108 }
109
110 if (storageLayout_ & DataStorage::dslQuadrupole) {
111 sdat.quadrupole = &(snap_->atomData.quadrupole[atom1]);
112 }
113
114 if (storageLayout_ & DataStorage::dslTorque) {
115 sdat.t = &(snap_->atomData.torque[atom1]);
116 }
117
118 if (storageLayout_ & DataStorage::dslDensity) {
119 sdat.rho = &(snap_->atomData.density[atom1]);
120 }
121
122 if (storageLayout_ & DataStorage::dslFunctional) {
123 sdat.frho = &(snap_->atomData.functional[atom1]);
124 }
125
126 if (storageLayout_ & DataStorage::dslFunctionalDerivative) {
127 sdat.dfrhodrho = &(snap_->atomData.functionalDerivative[atom1]);
128 }
129
130 if (storageLayout_ & DataStorage::dslSkippedCharge) {
131 sdat.skippedCharge = &(snap_->atomData.skippedCharge[atom1]);
132 }
133
134 if (storageLayout_ & DataStorage::dslParticlePot) {
135 sdat.particlePot = &(snap_->atomData.particlePot[atom1]);
136 }
137
138 if (storageLayout_ & DataStorage::dslFlucQPosition) {
139 sdat.flucQ = &(snap_->atomData.flucQPos[atom1]);
140 }
141
142 if (storageLayout_ & DataStorage::dslFlucQForce) {
143 sdat.dVdFQ = &(snap_->atomData.flucQFrc[atom1]);
144 }
145 }
146
147 bool ForceDecomposition::checkNeighborList() {
148
149 int nGroups = snap_->cgData.position.size();
150 if (needVelocities_)
151 snap_->cgData.setStorageLayout(DataStorage::dslPosition | DataStorage::dslVelocity);
152
153 // if we have changed the group identities or haven't set up the
154 // saved positions we automatically will need a neighbor list update:
155
156 if ( saved_CG_positions_.size() != nGroups ) return true;
157
158 RealType dispmax = 0.0;
159 Vector3d disp;
160
161 for (int i = 0; i < nGroups; i++) {
162 disp = snap_->cgData.position[i] - saved_CG_positions_[i];
163 for (int j = 0; j < 3; j++)
164 dispmax = max( abs(disp[j]), dispmax);
165 }
166
167 #ifdef IS_MPI
168 MPI::COMM_WORLD.Allreduce(&dispmax, &dispmax, 1, MPI::REALTYPE, MPI::MAX);
169 #endif
170
171 // a conservative test of list skin crossings
172 dispmax = 2.0 * sqrt (3.0 * dispmax * dispmax);
173
174
175 if (dispmax > skinThickness_)
176 return (dispmax > skinThickness_);
177
178 return false;
179 }
180
181 void ForceDecomposition::addToHeatFlux(Vector3d hf) {
182 Vector3d chf = snap_->getConductiveHeatFlux();
183 chf += hf;
184 snap_->setConductiveHeatFlux(chf);
185 }
186 void ForceDecomposition::setHeatFlux(Vector3d hf) {
187 snap_->setConductiveHeatFlux(hf);
188 }
189
190 }

Properties

Name Value
svn:eol-style native