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