ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/perturbations/UniformGradient.cpp
Revision: 2034
Committed: Mon Nov 3 16:49:03 2014 UTC (10 years, 5 months ago) by gezelter
File size: 8289 byte(s)
Log Message:
Updating UniformGradient to the new parameter structure (two unit vectors
and a gradient strength).

File Contents

# Content
1 /*
2 * Copyright (c) 2014 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 #include "perturbations/UniformGradient.hpp"
44 #include "types/FixedChargeAdapter.hpp"
45 #include "types/FluctuatingChargeAdapter.hpp"
46 #include "types/MultipoleAdapter.hpp"
47 #include "primitives/Molecule.hpp"
48 #include "nonbonded/NonBondedInteraction.hpp"
49 #include "utils/PhysicalConstants.hpp"
50
51 namespace OpenMD {
52
53 UniformGradient::UniformGradient(SimInfo* info) : info_(info),
54 doUniformGradient(false),
55 doParticlePot(false),
56 initialized(false) {
57 simParams = info_->getSimParams();
58 }
59
60 void UniformGradient::initialize() {
61
62 bool haveA = false;
63 bool haveB = false;
64 bool haveG = false;
65
66 if (simParams->haveUniformGradientDirection1()) {
67 std::vector<RealType> d1 = simParams->getUniformGradientDirection1();
68 if (d1.size() != 3) {
69 sprintf(painCave.errMsg,
70 "uniformGradientDirection1: Incorrect number of parameters\n"
71 "\tspecified. There should be 3 parameters, but %lu were\n"
72 "\tspecified.\n", d1.size());
73 painCave.isFatal = 1;
74 simError();
75 }
76 a_.x() = d1[0];
77 a_.y() = d1[1];
78 a_.z() = d1[2];
79
80 a_.normalize();
81 haveA = true;
82 }
83
84 if (simParams->haveUniformGradientDirection2()) {
85 std::vector<RealType> d2 = simParams->getUniformGradientDirection2();
86 if (d2.size() != 3) {
87 sprintf(painCave.errMsg,
88 "uniformGradientDirection2: Incorrect number of parameters\n"
89 "\tspecified. There should be 3 parameters, but %lu were\n"
90 "\tspecified.\n", d2.size());
91 painCave.isFatal = 1;
92 simError();
93 }
94 b_.x() = d2[0];
95 b_.y() = d2[1];
96 b_.z() = d2[2];
97
98 b_.normalize();
99 haveB = true;
100 }
101
102 if (simParams->haveUniformGradientStrength()) {
103 g_ = simParams->getUniformGradientStrength();
104 haveG = true;
105 }
106
107 if (haveA && haveB && haveG) {
108 doUniformGradient = true;
109 cpsi_ = dot(a_, b_);
110
111 Grad_(0,0) = 2.0 * (a_.x()*b_.x() - cpsi_ / 3.0);
112 Grad_(0,1) = a_.x()*b_.y() + a_.y()*b_.x();
113 Grad_(0,2) = a_.x()*b_.z() + a_.z()*b_.x();
114 Grad_(1,0) = Grad_(0,1);
115 Grad_(1,1) = 2.0 * (a_.y()*b_.y() - cpsi_ / 3.0);
116 Grad_(1,2) = a_.y()*b_.z() + a_.z()*b_.y();
117 Grad_(2,0) = Grad_(0,2);
118 Grad_(2,1) = Grad_(1,2);
119 Grad_(2,2) = 2.0 * (a_.z()*b_.z() - cpsi_ / 3.0);
120
121 Grad_ *= g_ / 2.0;
122 } else {
123 if (!haveA) {
124 sprintf(painCave.errMsg,
125 "UniformGradient: uniformGradientDirection1 not specified.\n");
126 painCave.isFatal = 1;
127 simError();
128 }
129 if (!haveB) {
130 sprintf(painCave.errMsg,
131 "UniformGradient: uniformGradientDirection2 not specified.\n");
132 painCave.isFatal = 1;
133 simError();
134 }
135 if (!haveG) {
136 sprintf(painCave.errMsg,
137 "UniformGradient: uniformGradientStrength not specified.\n");
138 painCave.isFatal = 1;
139 simError();
140 }
141 }
142
143 int storageLayout_ = info_->getSnapshotManager()->getStorageLayout();
144 if (storageLayout_ & DataStorage::dslParticlePot) doParticlePot = true;
145 initialized = true;
146 }
147
148 void UniformGradient::applyPerturbation() {
149
150 if (!initialized) initialize();
151
152 SimInfo::MoleculeIterator i;
153 Molecule::AtomIterator j;
154 Molecule* mol;
155 Atom* atom;
156 AtomType* atype;
157 potVec longRangePotential(0.0);
158
159 RealType C;
160 Vector3d D;
161 Mat3x3d Q;
162
163 RealType U;
164 RealType fPot;
165 Vector3d t;
166 Vector3d f;
167
168 Vector3d r;
169 Vector3d EF;
170
171 bool isCharge;
172
173 if (doUniformGradient) {
174
175 U = 0.0;
176 fPot = 0.0;
177
178 for (mol = info_->beginMolecule(i); mol != NULL;
179 mol = info_->nextMolecule(i)) {
180
181 for (atom = mol->beginAtom(j); atom != NULL;
182 atom = mol->nextAtom(j)) {
183
184 isCharge = false;
185 C = 0.0;
186
187 atype = atom->getAtomType();
188
189 r = atom->getPos();
190 EF = Grad_ * r;
191
192 if (atype->isElectrostatic()) {
193 atom->addElectricField(EF * PhysicalConstants::chargeFieldConvert);
194 }
195
196 FixedChargeAdapter fca = FixedChargeAdapter(atype);
197 if ( fca.isFixedCharge() ) {
198 isCharge = true;
199 C = fca.getCharge();
200 }
201
202 FluctuatingChargeAdapter fqa = FluctuatingChargeAdapter(atype);
203 if ( fqa.isFluctuatingCharge() ) {
204 isCharge = true;
205 C += atom->getFlucQPos();
206 atom->addFlucQFrc( dot(r, EF)
207 * PhysicalConstants::chargeFieldConvert );
208 }
209
210 if (isCharge) {
211 f = EF * C * PhysicalConstants::chargeFieldConvert;
212 atom->addFrc(f);
213
214 U = -dot(r, f);
215 if (doParticlePot) {
216 atom->addParticlePot(U);
217 }
218 fPot += U;
219 }
220
221 MultipoleAdapter ma = MultipoleAdapter(atype);
222 if (ma.isDipole() ) {
223 D = atom->getDipole() * PhysicalConstants::dipoleFieldConvert;
224
225 f = D * Grad_;
226 atom->addFrc(f);
227
228 t = cross(D, EF);
229 atom->addTrq(t);
230
231 U = -dot(D, EF);
232 if (doParticlePot) {
233 atom->addParticlePot(U);
234 }
235 fPot += U;
236 }
237
238 if (ma.isQuadrupole() ) {
239 Q = atom->getQuadrupole() * PhysicalConstants::dipoleFieldConvert;
240
241 t = 2.0 * mCross(Q, Grad_);
242 atom->addTrq(t);
243
244 U = -doubleDot(Q, Grad_);
245 if (doParticlePot) {
246 atom->addParticlePot(U);
247 }
248 fPot += U;
249 }
250 }
251 }
252
253 #ifdef IS_MPI
254 MPI_Allreduce(MPI_IN_PLACE, &fPot, 1, MPI_REALTYPE,
255 MPI_SUM, MPI_COMM_WORLD);
256 #endif
257
258 Snapshot* snap = info_->getSnapshotManager()->getCurrentSnapshot();
259 longRangePotential = snap->getLongRangePotentials();
260 longRangePotential[ELECTROSTATIC_FAMILY] += fPot;
261 snap->setLongRangePotential(longRangePotential);
262 }
263 }
264 }

Properties

Name Value
svn:executable *