ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/nonbonded/MAW.cpp
Revision: 1536
Committed: Wed Jan 5 14:49:05 2011 UTC (14 years, 6 months ago) by gezelter
File size: 10392 byte(s)
Log Message:
compiles, builds and runs, but is very slow

File Contents

# User Rev Content
1 gezelter 1532 /*
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 gezelter 1536 void MAW::calcForce(InteractionData &idat) {
128 gezelter 1532
129     if (!initialized_) initialize();
130    
131     pair<AtomType*, AtomType*> key = make_pair(idat.atype1, idat.atype2);
132     map<pair<AtomType*, AtomType*>, MAWInteractionData>::iterator it;
133     it = MixingMap.find(key);
134     if (it != MixingMap.end()) {
135     MAWInteractionData mixer = (*it).second;
136    
137     RealType myPot = 0.0;
138     RealType myPotC = 0.0;
139     RealType myDeriv = 0.0;
140     RealType myDerivC = 0.0;
141    
142     RealType D_e = mixer.De;
143     RealType R_e = mixer.Re;
144     RealType beta = mixer.beta;
145     RealType ca1 = mixer.ca1;
146     RealType cb1 = mixer.cb1;
147    
148     bool j_is_Metal = idat.atype2->isMetal();
149    
150     Vector3d r;
151     RotMat3x3d Atrans;
152     if (j_is_Metal) {
153     // rotate the inter-particle separation into the two different
154     // body-fixed coordinate systems:
155     r = idat.A1 * idat.d;
156     Atrans = idat.A1.transpose();
157     } else {
158     // negative sign because this is the vector from j to i:
159     r = -idat.A2 * idat.d;
160     Atrans = idat.A2.transpose();
161     }
162    
163     // V(r) = D_e exp(-a(r-re)(exp(-a(r-re))-2)
164    
165     RealType expt = -beta*(idat.rij - R_e);
166     RealType expfnc = exp(expt);
167     RealType expfnc2 = expfnc*expfnc;
168    
169     RealType exptC = 0.0;
170     RealType expfncC = 0.0;
171     RealType expfnc2C = 0.0;
172    
173     myPot = D_e * (expfnc2 - 2.0 * expfnc);
174     myDeriv = 2.0 * D_e * beta * (expfnc - expfnc2);
175    
176     if (MAW::shiftedPot_ || MAW::shiftedFrc_) {
177     exptC = -beta*(idat.rcut - R_e);
178     expfncC = exp(exptC);
179     expfnc2C = expfncC*expfncC;
180     }
181    
182     if (MAW::shiftedPot_) {
183     myPotC = D_e * (expfnc2C - 2.0 * expfncC);
184     myDerivC = 0.0;
185     } else if (MAW::shiftedFrc_) {
186     myPotC = D_e * (expfnc2C - 2.0 * expfncC);
187     myDerivC = 2.0 * D_e * beta * (expfnc2C - expfnc2C);
188     myPotC += myDerivC * (idat.rij - idat.rcut);
189     } else {
190     myPotC = 0.0;
191     myDerivC = 0.0;
192     }
193    
194     RealType x = r.x();
195     RealType y = r.y();
196     RealType z = r.z();
197     RealType x2 = x * x;
198     RealType y2 = y * y;
199     RealType z2 = z * z;
200    
201     RealType r3 = idat.r2 * idat.rij;
202     RealType r4 = idat.r2 * idat.r2;
203    
204     // angular modulation of morse part of potential to approximate
205     // the squares of the two HOMO lone pair orbitals in water:
206     //********************** old form*************************
207     // s = 1 / (4 pi)
208     // ta1 = (s - pz)^2 = (1 - sqrt(3)*cos(theta))^2 / (4 pi)
209     // b1 = px^2 = 3 * (sin(theta)*cos(phi))^2 / (4 pi)
210     //********************** old form*************************
211     // we'll leave out the 4 pi for now
212    
213     // new functional form just using the p orbitals.
214     // Vmorse(r)*[a*p_x + b p_z + (1-a-b)]
215     // which is
216     // Vmorse(r)*[a sin^2(theta) cos^2(phi) + b cos(theta) + (1-a-b)]
217     // Vmorse(r)*[a*x2/r2 + b*z/r + (1-a-b)]
218    
219     RealType Vmorse = (myPot - myPotC);
220     RealType Vang = ca1 * x2 / idat.r2 + cb1 * z / idat.rij + (0.8 - ca1 / 3.0);
221    
222     RealType pot_temp = idat.vdwMult * Vmorse * Vang;
223 gezelter 1536 idat.vpair[0] += pot_temp;
224     idat.pot[0] += idat.sw * pot_temp;
225 gezelter 1532
226     Vector3d dVmorsedr = (myDeriv - myDerivC) * Vector3d(x, y, z) / idat.rij;
227    
228     Vector3d dVangdr = Vector3d(-2.0 * ca1 * x2 * x / r4 + 2.0 * ca1 * x / idat.r2 - cb1 * x * z / r3,
229     -2.0 * ca1 * x2 * y / r4 - cb1 * z * y / r3,
230     -2.0 * ca1 * x2 * z / r4 + cb1 / idat.rij - cb1 * z2 / r3);
231    
232     // chain rule to put these back on x, y, z
233    
234     Vector3d dvdr = Vang * dVmorsedr + Vmorse * dVangdr;
235    
236     // Torques for Vang using method of Price:
237     // S. L. Price, A. J. Stone, and M. Alderton, Mol. Phys. 52, 987 (1984).
238    
239     Vector3d dVangdu = Vector3d(cb1 * y / idat.rij,
240     2.0 * ca1 * x * z / idat.r2 - cb1 * x / idat.rij,
241     -2.0 * ca1 * y * x / idat.r2);
242    
243     // do the torques first since they are easy:
244     // remember that these are still in the body fixed axes
245    
246     Vector3d trq = idat.vdwMult * Vmorse * dVangdu * idat.sw;
247    
248     // go back to lab frame using transpose of rotation matrix:
249    
250     if (j_is_Metal) {
251     idat.t1 += Atrans * trq;
252     } else {
253     idat.t2 += Atrans * trq;
254     }
255    
256     // Now, on to the forces (still in body frame of water)
257    
258     Vector3d ftmp = idat.vdwMult * idat.sw * dvdr;
259    
260     // rotate the terms back into the lab frame:
261     Vector3d flab;
262     if (j_is_Metal) {
263     flab = Atrans * ftmp;
264     } else {
265     flab = - Atrans * ftmp;
266     }
267    
268     idat.f1 += flab;
269     }
270     return;
271    
272     }
273    
274     RealType MAW::getSuggestedCutoffRadius(AtomType* at1, AtomType* at2) {
275     if (!initialized_) initialize();
276     pair<AtomType*, AtomType*> key = make_pair(at1, at2);
277     map<pair<AtomType*, AtomType*>, MAWInteractionData>::iterator it;
278     it = MixingMap.find(key);
279     if (it == MixingMap.end())
280     return 0.0;
281     else {
282     MAWInteractionData mixer = (*it).second;
283    
284     RealType R_e = mixer.Re;
285     RealType beta = mixer.beta;
286     // This value of the r corresponds to an energy about 1.48% of
287     // the energy at the bottom of the Morse well. For comparison, the
288     // Lennard-Jones function is about 1.63% of it's minimum value at
289     // a distance of 2.5 sigma.
290     return (4.9 + beta * R_e) / beta;
291     }
292     }
293     }
294    

Properties

Name Value
svn:eol-style native