ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/nonbonded/MAW.cpp
Revision: 1549
Committed: Wed Apr 27 18:38:15 2011 UTC (14 years ago) by gezelter
File size: 10278 byte(s)
Log Message:
a few more tweaks   We're getting somewhat closer to deleting fortran.

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
42 #include <stdio.h>
43 #include <string.h>
44 #include <cmath>
45
46 #include "nonbonded/MAW.hpp"
47 #include "utils/simError.h"
48
49 using namespace std;
50
51 namespace OpenMD {
52
53 MAW::MAW() : name_("MAW"), initialized_(false), forceField_(NULL),
54 shiftedPot_(false), shiftedFrc_(false) {}
55
56 void MAW::initialize() {
57
58 ForceField::NonBondedInteractionTypeContainer* nbiTypes = forceField_->getNonBondedInteractionTypes();
59 ForceField::NonBondedInteractionTypeContainer::MapTypeIterator j;
60 NonBondedInteractionType* nbt;
61
62 for (nbt = nbiTypes->beginType(j); nbt != NULL;
63 nbt = nbiTypes->nextType(j)) {
64
65 if (nbt->isMAW()) {
66 pair<AtomType*, AtomType*> atypes = nbt->getAtomTypes();
67
68 GenericData* data = nbt->getPropertyByName("MAW");
69 if (data == NULL) {
70 sprintf( painCave.errMsg, "MAW::initialize could not find\n"
71 "\tMAW parameters for %s - %s interaction.\n",
72 atypes.first->getName().c_str(),
73 atypes.second->getName().c_str());
74 painCave.severity = OPENMD_ERROR;
75 painCave.isFatal = 1;
76 simError();
77 }
78
79 MAWData* mawData = dynamic_cast<MAWData*>(data);
80 if (mawData == NULL) {
81 sprintf( painCave.errMsg,
82 "MAW::initialize could not convert GenericData to\n"
83 "\tMAWData for %s - %s interaction.\n",
84 atypes.first->getName().c_str(),
85 atypes.second->getName().c_str());
86 painCave.severity = OPENMD_ERROR;
87 painCave.isFatal = 1;
88 simError();
89 }
90
91 MAWParam mawParam = mawData->getData();
92
93 RealType De = mawParam.De;
94 RealType beta = mawParam.beta;
95 RealType Re = mawParam.Re;
96 RealType ca1 = mawParam.ca1;
97 RealType cb1 = mawParam.cb1;
98
99 addExplicitInteraction(atypes.first, atypes.second,
100 De, beta, Re, ca1, cb1);
101 }
102 }
103 initialized_ = true;
104 }
105
106 void MAW::addExplicitInteraction(AtomType* atype1, AtomType* atype2,
107 RealType De, RealType beta, RealType Re,
108 RealType ca1, RealType cb1) {
109
110 MAWInteractionData mixer;
111 mixer.De = De;
112 mixer.beta = beta;
113 mixer.Re = Re;
114 mixer.ca1 = ca1;
115 mixer.cb1 = cb1;
116
117 pair<AtomType*, AtomType*> key1, key2;
118 key1 = make_pair(atype1, atype2);
119 key2 = make_pair(atype2, atype1);
120
121 MixingMap[key1] = mixer;
122 if (key2 != key1) {
123 MixingMap[key2] = mixer;
124 }
125 }
126
127 void MAW::calcForce(InteractionData &idat) {
128
129 if (!initialized_) initialize();
130
131 map<pair<AtomType*, AtomType*>, MAWInteractionData>::iterator it;
132 it = MixingMap.find(idat.atypes);
133 if (it != MixingMap.end()) {
134 MAWInteractionData mixer = (*it).second;
135
136 RealType myPot = 0.0;
137 RealType myPotC = 0.0;
138 RealType myDeriv = 0.0;
139 RealType myDerivC = 0.0;
140
141 RealType D_e = mixer.De;
142 RealType R_e = mixer.Re;
143 RealType beta = mixer.beta;
144 RealType ca1 = mixer.ca1;
145 RealType cb1 = mixer.cb1;
146
147 bool j_is_Metal = idat.atypes.second->isMetal();
148
149 Vector3d r;
150 RotMat3x3d Atrans;
151 if (j_is_Metal) {
152 // rotate the inter-particle separation into the two different
153 // body-fixed coordinate systems:
154 r = idat.A1 * idat.d;
155 Atrans = idat.A1.transpose();
156 } else {
157 // negative sign because this is the vector from j to i:
158 r = -idat.A2 * idat.d;
159 Atrans = idat.A2.transpose();
160 }
161
162 // V(r) = D_e exp(-a(r-re)(exp(-a(r-re))-2)
163
164 RealType expt = -beta*(idat.rij - R_e);
165 RealType expfnc = exp(expt);
166 RealType expfnc2 = expfnc*expfnc;
167
168 RealType exptC = 0.0;
169 RealType expfncC = 0.0;
170 RealType expfnc2C = 0.0;
171
172 myPot = D_e * (expfnc2 - 2.0 * expfnc);
173 myDeriv = 2.0 * D_e * beta * (expfnc - expfnc2);
174
175 if (MAW::shiftedPot_ || MAW::shiftedFrc_) {
176 exptC = -beta*(idat.rcut - R_e);
177 expfncC = exp(exptC);
178 expfnc2C = expfncC*expfncC;
179 }
180
181 if (MAW::shiftedPot_) {
182 myPotC = D_e * (expfnc2C - 2.0 * expfncC);
183 myDerivC = 0.0;
184 } else if (MAW::shiftedFrc_) {
185 myPotC = D_e * (expfnc2C - 2.0 * expfncC);
186 myDerivC = 2.0 * D_e * beta * (expfnc2C - expfnc2C);
187 myPotC += myDerivC * (idat.rij - idat.rcut);
188 } else {
189 myPotC = 0.0;
190 myDerivC = 0.0;
191 }
192
193 RealType x = r.x();
194 RealType y = r.y();
195 RealType z = r.z();
196 RealType x2 = x * x;
197 RealType y2 = y * y;
198 RealType z2 = z * z;
199
200 RealType r3 = idat.r2 * idat.rij;
201 RealType r4 = idat.r2 * idat.r2;
202
203 // angular modulation of morse part of potential to approximate
204 // the squares of the two HOMO lone pair orbitals in water:
205 //********************** old form*************************
206 // s = 1 / (4 pi)
207 // ta1 = (s - pz)^2 = (1 - sqrt(3)*cos(theta))^2 / (4 pi)
208 // b1 = px^2 = 3 * (sin(theta)*cos(phi))^2 / (4 pi)
209 //********************** old form*************************
210 // we'll leave out the 4 pi for now
211
212 // new functional form just using the p orbitals.
213 // Vmorse(r)*[a*p_x + b p_z + (1-a-b)]
214 // which is
215 // Vmorse(r)*[a sin^2(theta) cos^2(phi) + b cos(theta) + (1-a-b)]
216 // Vmorse(r)*[a*x2/r2 + b*z/r + (1-a-b)]
217
218 RealType Vmorse = (myPot - myPotC);
219 RealType Vang = ca1 * x2 / idat.r2 + cb1 * z / idat.rij + (0.8 - ca1 / 3.0);
220
221 RealType pot_temp = idat.vdwMult * Vmorse * Vang;
222 idat.vpair += pot_temp;
223 idat.pot[0] += idat.sw * pot_temp;
224
225 Vector3d dVmorsedr = (myDeriv - myDerivC) * Vector3d(x, y, z) / idat.rij;
226
227 Vector3d dVangdr = Vector3d(-2.0 * ca1 * x2 * x / r4 + 2.0 * ca1 * x / idat.r2 - cb1 * x * z / r3,
228 -2.0 * ca1 * x2 * y / r4 - cb1 * z * y / r3,
229 -2.0 * ca1 * x2 * z / r4 + cb1 / idat.rij - cb1 * z2 / r3);
230
231 // chain rule to put these back on x, y, z
232
233 Vector3d dvdr = Vang * dVmorsedr + Vmorse * dVangdr;
234
235 // Torques for Vang using method of Price:
236 // S. L. Price, A. J. Stone, and M. Alderton, Mol. Phys. 52, 987 (1984).
237
238 Vector3d dVangdu = Vector3d(cb1 * y / idat.rij,
239 2.0 * ca1 * x * z / idat.r2 - cb1 * x / idat.rij,
240 -2.0 * ca1 * y * x / idat.r2);
241
242 // do the torques first since they are easy:
243 // remember that these are still in the body fixed axes
244
245 Vector3d trq = idat.vdwMult * Vmorse * dVangdu * idat.sw;
246
247 // go back to lab frame using transpose of rotation matrix:
248
249 if (j_is_Metal) {
250 idat.t1 += Atrans * trq;
251 } else {
252 idat.t2 += Atrans * trq;
253 }
254
255 // Now, on to the forces (still in body frame of water)
256
257 Vector3d ftmp = idat.vdwMult * idat.sw * dvdr;
258
259 // rotate the terms back into the lab frame:
260 Vector3d flab;
261 if (j_is_Metal) {
262 flab = Atrans * ftmp;
263 } else {
264 flab = - Atrans * ftmp;
265 }
266
267 idat.f1 += flab;
268 }
269 return;
270
271 }
272
273 RealType MAW::getSuggestedCutoffRadius(pair<AtomType*, AtomType*> atypes) {
274 if (!initialized_) initialize();
275 map<pair<AtomType*, AtomType*>, MAWInteractionData>::iterator it;
276 it = MixingMap.find(atypes);
277 if (it == MixingMap.end())
278 return 0.0;
279 else {
280 MAWInteractionData mixer = (*it).second;
281
282 RealType R_e = mixer.Re;
283 RealType beta = mixer.beta;
284 // This value of the r corresponds to an energy about 1.48% of
285 // the energy at the bottom of the Morse well. For comparison, the
286 // Lennard-Jones function is about 1.63% of it's minimum value at
287 // a distance of 2.5 sigma.
288 return (4.9 + beta * R_e) / beta;
289 }
290 }
291 }
292

Properties

Name Value
svn:eol-style native