ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/parallel/ForceDecomposition.cpp
Revision: 2057
Committed: Tue Mar 3 15:22:26 2015 UTC (10 years, 1 month ago) by gezelter
File size: 6934 byte(s)
Log Message:
Performance improvements, Removed CutoffPolicy

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

Properties

Name Value
svn:eol-style native