OpenMD 3.0
Molecular Dynamics in the Open
Loading...
Searching...
No Matches
MAW.cpp
1/*
2 * Copyright (c) 2004-present, The University of Notre Dame. All rights
3 * reserved.
4 *
5 * Redistribution and use in source and binary forms, with or without
6 * modification, are permitted provided that the following conditions are met:
7 *
8 * 1. Redistributions of source code must retain the above copyright notice,
9 * this list of conditions and the following disclaimer.
10 *
11 * 2. Redistributions in binary form must reproduce the above copyright notice,
12 * this list of conditions and the following disclaimer in the documentation
13 * and/or other materials provided with the distribution.
14 *
15 * 3. Neither the name of the copyright holder nor the names of its
16 * contributors may be used to endorse or promote products derived from
17 * this software without specific prior written permission.
18 *
19 * THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
20 * AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
21 * IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
22 * ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
23 * LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
24 * CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF
25 * SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS
26 * INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN
27 * CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE)
28 * ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
29 * POSSIBILITY OF SUCH DAMAGE.
30 *
31 * SUPPORT OPEN SCIENCE! If you use OpenMD or its source code in your
32 * research, please cite the appropriate papers when you publish your
33 * work. Good starting points are:
34 *
35 * [1] Meineke, et al., J. Comp. Chem. 26, 252-271 (2005).
36 * [2] Fennell & Gezelter, J. Chem. Phys. 124, 234104 (2006).
37 * [3] Sun, Lin & Gezelter, J. Chem. Phys. 128, 234107 (2008).
38 * [4] Vardeman, Stocker & Gezelter, J. Chem. Theory Comput. 7, 834 (2011).
39 * [5] Kuang & Gezelter, Mol. Phys., 110, 691-701 (2012).
40 * [6] Lamichhane, Gezelter & Newman, J. Chem. Phys. 141, 134109 (2014).
41 * [7] Lamichhane, Newman & Gezelter, J. Chem. Phys. 141, 134110 (2014).
42 * [8] Bhattarai, Newman & Gezelter, Phys. Rev. B 99, 094106 (2019).
43 */
44
45#include "nonbonded/MAW.hpp"
46
47#include <cmath>
48#include <cstdio>
49#include <cstring>
50
51#include "utils/simError.h"
52
53using namespace std;
54
55namespace OpenMD {
56
57 MAW::MAW() : initialized_(false), forceField_(NULL), name_("MAW") {}
58
59 void MAW::initialize() {
60 MAWtypes.clear();
61 MAWtids.clear();
62 MixingMap.clear();
63 MAWtids.resize(forceField_->getNAtomType(), -1);
64
65 ForceField::NonBondedInteractionTypeContainer* nbiTypes =
66 forceField_->getNonBondedInteractionTypes();
67 ForceField::NonBondedInteractionTypeContainer::MapTypeIterator j;
68 ForceField::NonBondedInteractionTypeContainer::KeyType keys;
69 NonBondedInteractionType* nbt;
70 int mtid1, mtid2;
71
72 for (nbt = nbiTypes->beginType(j); nbt != NULL;
73 nbt = nbiTypes->nextType(j)) {
74 if (nbt->isMAW()) {
75 keys = nbiTypes->getKeys(j);
76 AtomType* at1 = forceField_->getAtomType(keys[0]);
77 if (at1 == NULL) {
78 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
79 "MAW::initialize could not find AtomType %s\n"
80 "\tto for for %s - %s interaction.\n",
81 keys[0].c_str(), keys[0].c_str(), keys[1].c_str());
82 painCave.severity = OPENMD_ERROR;
83 painCave.isFatal = 1;
84 simError();
85 }
86
87 AtomType* at2 = forceField_->getAtomType(keys[1]);
88 if (at2 == NULL) {
89 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
90 "MAW::initialize could not find AtomType %s\n"
91 "\tfor %s - %s nonbonded interaction.\n",
92 keys[1].c_str(), keys[0].c_str(), keys[1].c_str());
93 painCave.severity = OPENMD_ERROR;
94 painCave.isFatal = 1;
95 simError();
96 }
97
98 int atid1 = at1->getIdent();
99 if (MAWtids[atid1] == -1) {
100 mtid1 = MAWtypes.size();
101 MAWtypes.insert(atid1);
102 MAWtids[atid1] = mtid1;
103 }
104 int atid2 = at2->getIdent();
105 if (MAWtids[atid2] == -1) {
106 mtid2 = MAWtypes.size();
107 MAWtypes.insert(atid2);
108 MAWtids[atid2] = mtid2;
109 }
110
111 MAWInteractionType* mit = dynamic_cast<MAWInteractionType*>(nbt);
112
113 if (mit == NULL) {
114 snprintf(
115 painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
116 "MAW::initialize could not convert NonBondedInteractionType\n"
117 "\tto MAWInteractionType for %s - %s interaction.\n",
118 at1->getName().c_str(), at2->getName().c_str());
119 painCave.severity = OPENMD_ERROR;
120 painCave.isFatal = 1;
121 simError();
122 }
123
124 RealType De = mit->getD();
125 RealType beta = mit->getBeta();
126 RealType Re = mit->getR();
127 RealType ca1 = mit->getCA1();
128 RealType cb1 = mit->getCB1();
129
130 addExplicitInteraction(at1, at2, De, beta, Re, ca1, cb1);
131 }
132 }
133 initialized_ = true;
134 }
135
136 void MAW::addExplicitInteraction(AtomType* atype1, AtomType* atype2,
137 RealType De, RealType beta, RealType Re,
138 RealType ca1, RealType cb1) {
139 MAWInteractionData mixer;
140 mixer.De = De;
141 mixer.beta = beta;
142 mixer.Re = Re;
143 mixer.ca1 = ca1;
144 mixer.cb1 = cb1;
145 mixer.j_is_Metal = atype2->isMetal();
146
147 int mtid1 = MAWtids[atype1->getIdent()];
148 int mtid2 = MAWtids[atype2->getIdent()];
149 int nM = MAWtypes.size();
150
151 MixingMap.resize(nM);
152 MixingMap[mtid1].resize(nM);
153
154 MixingMap[mtid1][mtid2] = mixer;
155 if (mtid2 != mtid1) {
156 MixingMap[mtid2].resize(nM);
157 mixer.j_is_Metal = atype1->isMetal();
158 MixingMap[mtid2][mtid1] = mixer;
159 }
160 }
161
162 void MAW::calcForce(InteractionData& idat) {
163 if (!initialized_) initialize();
164
165 MAWInteractionData& mixer =
166 MixingMap[MAWtids[idat.atid1]][MAWtids[idat.atid2]];
167
168 RealType myPot = 0.0;
169 RealType myPotC = 0.0;
170 RealType myDeriv = 0.0;
171 RealType myDerivC = 0.0;
172
173 RealType D_e = mixer.De;
174 RealType R_e = mixer.Re;
175 RealType beta = mixer.beta;
176 RealType ca1 = mixer.ca1;
177 RealType cb1 = mixer.cb1;
178
179 Vector3d r;
180 RotMat3x3d Atrans;
181 if (mixer.j_is_Metal) {
182 // rotate the inter-particle separation into the two different
183 // body-fixed coordinate systems:
184 r = idat.A1 * idat.d;
185 Atrans = idat.A1.transpose();
186 } else {
187 // negative sign because this is the vector from j to i:
188 r = -idat.A2 * idat.d;
189 Atrans = idat.A2.transpose();
190 }
191
192 // V(r) = D_e exp(-a(r-re)(exp(-a(r-re))-2)
193
194 RealType expt = -beta * (idat.rij - R_e);
195 RealType expfnc = exp(expt);
196 RealType expfnc2 = expfnc * expfnc;
197
198 RealType exptC = 0.0;
199 RealType expfncC = 0.0;
200 RealType expfnc2C = 0.0;
201
202 myPot = D_e * (expfnc2 - 2.0 * expfnc);
203 myDeriv = 2.0 * D_e * beta * (expfnc - expfnc2);
204
205 if (idat.shiftedPot || idat.shiftedForce) {
206 exptC = -beta * (idat.rcut - R_e);
207 expfncC = exp(exptC);
208 expfnc2C = expfncC * expfncC;
209 }
210
211 if (idat.shiftedPot) {
212 myPotC = D_e * (expfnc2C - 2.0 * expfncC);
213 myDerivC = 0.0;
214 } else if (idat.shiftedForce) {
215 myPotC = D_e * (expfnc2C - 2.0 * expfncC);
216 myDerivC = 2.0 * D_e * beta * (expfncC - expfnc2C);
217 myPotC += myDerivC * (idat.rij - idat.rcut);
218 } else {
219 myPotC = 0.0;
220 myDerivC = 0.0;
221 }
222
223 RealType x = r.x();
224 RealType y = r.y();
225 RealType z = r.z();
226 RealType x2 = x * x;
227 RealType z2 = z * z;
228
229 RealType r3 = idat.r2 * idat.rij;
230 RealType r4 = idat.r2 * idat.r2;
231
232 // angular modulation of morse part of potential to approximate
233 // the squares of the two HOMO lone pair orbitals in water:
234 //********************** old form*************************
235 // s = 1 / (4 pi)
236 // ta1 = (s - pz)^2 = (1 - sqrt(3)*cos(theta))^2 / (4 pi)
237 // b1 = px^2 = 3 * (sin(theta)*cos(phi))^2 / (4 pi)
238 //********************** old form*************************
239 // we'll leave out the 4 pi for now
240
241 // new functional form just using the p orbitals.
242 // Vmorse(r)*[a*p_x + b p_z + (1-a-b)]
243 // which is
244 // Vmorse(r)*[a sin^2(theta) cos^2(phi) + b cos(theta) + (1-a-b)]
245 // Vmorse(r)*[a*x2/r2 + b*z/r + (1-a-b)]
246
247 RealType Vmorse = (myPot - myPotC);
248 RealType Vang = ca1 * x2 / idat.r2 + cb1 * z / idat.rij + (0.8 - ca1 / 3.0);
249
250 RealType pot_temp = idat.vdwMult * Vmorse * Vang;
251 idat.vpair += pot_temp;
252 idat.pot[VANDERWAALS_FAMILY] += idat.sw * pot_temp;
253 if (idat.isSelected) idat.selePot[VANDERWAALS_FAMILY] += idat.sw * pot_temp;
254 Vector3d dVmorsedr = (myDeriv - myDerivC) * Vector3d(x, y, z) / idat.rij;
255
256 Vector3d dVangdr = Vector3d(
257 -2.0 * ca1 * x2 * x / r4 + 2.0 * ca1 * x / idat.r2 - cb1 * x * z / r3,
258 -2.0 * ca1 * x2 * y / r4 - cb1 * z * y / r3,
259 -2.0 * ca1 * x2 * z / r4 + cb1 / idat.rij - cb1 * z2 / r3);
260
261 // chain rule to put these back on x, y, z
262
263 Vector3d dvdr = Vang * dVmorsedr + Vmorse * dVangdr;
264
265 // Torques for Vang using method of Price:
266 // S. L. Price, A. J. Stone, and M. Alderton, Mol. Phys. 52, 987 (1984).
267
268 Vector3d dVangdu = Vector3d(
269 cb1 * y / idat.rij, 2.0 * ca1 * x * z / idat.r2 - cb1 * x / idat.rij,
270 -2.0 * ca1 * y * x / idat.r2);
271
272 // do the torques first since they are easy:
273 // remember that these are still in the body fixed axes
274
275 Vector3d trq = idat.vdwMult * Vmorse * dVangdu * idat.sw;
276
277 // go back to lab frame using transpose of rotation matrix:
278
279 if (mixer.j_is_Metal) {
280 idat.t1 += Atrans * trq;
281 } else {
282 idat.t2 += Atrans * trq;
283 }
284
285 // Now, on to the forces (still in body frame of water)
286
287 Vector3d ftmp = idat.vdwMult * idat.sw * dvdr;
288
289 // rotate the terms back into the lab frame:
290 Vector3d flab;
291 if (mixer.j_is_Metal) {
292 flab = Atrans * ftmp;
293 } else {
294 flab = -Atrans * ftmp;
295 }
296
297 idat.f1 += flab;
298
299 return;
300 }
301
302 RealType MAW::getSuggestedCutoffRadius(pair<AtomType*, AtomType*> atypes) {
303 if (!initialized_) initialize();
304 int atid1 = atypes.first->getIdent();
305 int atid2 = atypes.second->getIdent();
306 int mtid1 = MAWtids[atid1];
307 int mtid2 = MAWtids[atid2];
308
309 if (mtid1 == -1 || mtid2 == -1)
310 return 0.0;
311 else {
312 MAWInteractionData mixer = MixingMap[mtid1][mtid2];
313 RealType R_e = mixer.Re;
314 RealType beta = mixer.beta;
315 // This value of the r corresponds to an energy about 1.48% of
316 // the energy at the bottom of the Morse well. For comparison, the
317 // Lennard-Jones function is about 1.63% of it's minimum value at
318 // a distance of 2.5 sigma.
319 return (4.9 + beta * R_e) / beta;
320 }
321 }
322} // namespace OpenMD
This basic Periodic Table class was originally taken from the data.cpp file in OpenBabel.
@ VANDERWAALS_FAMILY
Long-range dispersion and short-range pauli repulsion.