ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/nonbonded/Sticky.cpp
Revision: 1834
Committed: Tue Jan 15 16:28:42 2013 UTC (12 years, 6 months ago) by gezelter
File size: 12648 byte(s)
Log Message:
Bug fix in Sticky potential, label in stats

File Contents

# User Rev Content
1 gezelter 1485 /*
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 gezelter 1665 * [4] Kuang & Gezelter, J. Chem. Phys. 133, 164101 (2010).
40     * [5] Vardeman, Stocker & Gezelter, J. Chem. Theory Comput. 7, 834 (2011).
41 gezelter 1485 */
42    
43     #include <stdio.h>
44     #include <string.h>
45    
46     #include <cmath>
47     #include "nonbonded/Sticky.hpp"
48     #include "nonbonded/LJ.hpp"
49 gezelter 1710 #include "types/StickyAdapter.hpp"
50 gezelter 1485 #include "utils/simError.h"
51    
52     using namespace std;
53     namespace OpenMD {
54    
55 gezelter 1502 Sticky::Sticky() : name_("Sticky"), initialized_(false), forceField_(NULL) {}
56 gezelter 1485
57     void Sticky::initialize() {
58    
59     ForceFieldOptions& fopts = forceField_->getForceFieldOptions();
60     ForceField::AtomTypeContainer* atomTypes = forceField_->getAtomTypes();
61     ForceField::AtomTypeContainer::MapTypeIterator i;
62     AtomType* at;
63    
64     // Sticky handles all of the Sticky-Sticky interactions
65    
66     for (at = atomTypes->beginType(i); at != NULL;
67     at = atomTypes->nextType(i)) {
68    
69 gezelter 1710 StickyAdapter sa = StickyAdapter(at);
70     if (sa.isSticky()) addType(at);
71 gezelter 1485 }
72    
73     initialized_ = true;
74     }
75    
76     void Sticky::addType(AtomType* atomType){
77     // add it to the map:
78    
79     pair<map<int,AtomType*>::iterator,bool> ret;
80 gezelter 1710 ret = StickyMap.insert( pair<int, AtomType*>(atomType->getIdent(),
81     atomType) );
82 gezelter 1485 if (ret.second == false) {
83     sprintf( painCave.errMsg,
84     "Sticky already had a previous entry with ident %d\n",
85 gezelter 1710 atomType->getIdent() );
86 gezelter 1485 painCave.severity = OPENMD_INFO;
87     painCave.isFatal = 0;
88     simError();
89 gezelter 1767 }
90 gezelter 1485
91 gezelter 1710 StickyAdapter sticky1 = StickyAdapter(atomType);
92 gezelter 1485
93     // Now, iterate over all known types and add to the mixing map:
94    
95     map<int, AtomType*>::iterator it;
96     for( it = StickyMap.begin(); it != StickyMap.end(); ++it) {
97    
98     AtomType* atype2 = (*it).second;
99 gezelter 1710
100     StickyAdapter sticky2 = StickyAdapter(atype2);
101 gezelter 1485
102     StickyInteractionData mixer;
103    
104     // Mixing two different sticky types is silly, but if you want 2
105     // sticky types in your simulation, we'll let you do it with the
106     // Lorentz- Berthelot mixing rules (which happen to do the right thing
107     // when atomType and atype2 happen to be the same.
108    
109 gezelter 1710 mixer.rl = 0.5 * ( sticky1.getRl() + sticky2.getRl() );
110     mixer.ru = 0.5 * ( sticky1.getRu() + sticky2.getRu() );
111     mixer.rlp = 0.5 * ( sticky1.getRlp() + sticky2.getRlp() );
112     mixer.rup = 0.5 * ( sticky1.getRup() + sticky2.getRup() );
113 gezelter 1485 mixer.rbig = max(mixer.ru, mixer.rup);
114 gezelter 1710 mixer.w0 = sqrt( sticky1.getW0() * sticky2.getW0() );
115 gezelter 1834 mixer.v0 = sqrt( sticky1.getV0() * sticky2.getV0() );
116 gezelter 1710 mixer.v0p = sqrt( sticky1.getV0p() * sticky2.getV0p() );
117     mixer.isPower = sticky1.isStickyPower() && sticky2.isStickyPower();
118 gezelter 1485
119     CubicSpline* s = new CubicSpline();
120     s->addPoint(mixer.rl, 1.0);
121     s->addPoint(mixer.ru, 0.0);
122     mixer.s = s;
123    
124     CubicSpline* sp = new CubicSpline();
125     sp->addPoint(mixer.rlp, 1.0);
126     sp->addPoint(mixer.rup, 0.0);
127     mixer.sp = sp;
128    
129    
130     pair<AtomType*, AtomType*> key1, key2;
131     key1 = make_pair(atomType, atype2);
132     key2 = make_pair(atype2, atomType);
133    
134     MixingMap[key1] = mixer;
135     if (key2 != key1) {
136     MixingMap[key2] = mixer;
137     }
138     }
139     }
140    
141 gezelter 1502 /**
142     * This function does the sticky portion of the SSD potential
143     * [Chandra and Ichiye, Journal of Chemical Physics 111, 2701
144     * (1999)]. The Lennard-Jones and dipolar interaction must be
145     * handled separately. We assume that the rotation matrices have
146     * already been calculated and placed in the A1 & A2 entries in the
147     * idat structure.
148     */
149    
150 gezelter 1536 void Sticky::calcForce(InteractionData &idat) {
151 gezelter 1502
152 gezelter 1485 if (!initialized_) initialize();
153    
154 gezelter 1505 map<pair<AtomType*, AtomType*>, StickyInteractionData>::iterator it;
155 gezelter 1571 it = MixingMap.find(idat.atypes);
156 gezelter 1505 if (it != MixingMap.end()) {
157 gezelter 1485
158 gezelter 1505 StickyInteractionData mixer = (*it).second;
159    
160     RealType w0 = mixer.w0;
161     RealType v0 = mixer.v0;
162     RealType v0p = mixer.v0p;
163     RealType rl = mixer.rl;
164     RealType ru = mixer.ru;
165     RealType rlp = mixer.rlp;
166     RealType rup = mixer.rup;
167     RealType rbig = mixer.rbig;
168     bool isPower = mixer.isPower;
169    
170 gezelter 1834 cerr << rl << " " << ru << " " << rlp << " " << rup << " " << rbig << " " << w0 << " " << v0 << " " << v0p << "\n";
171 gezelter 1554 if ( *(idat.rij) <= rbig) {
172 gezelter 1505
173 gezelter 1554 RealType r3 = *(idat.r2) * *(idat.rij);
174     RealType r5 = r3 * *(idat.r2);
175 gezelter 1505
176 gezelter 1554 RotMat3x3d A1trans = idat.A1->transpose();
177     RotMat3x3d A2trans = idat.A2->transpose();
178 gezelter 1505
179     // rotate the inter-particle separation into the two different
180     // body-fixed coordinate systems:
181    
182 gezelter 1554 Vector3d ri = *(idat.A1) * *(idat.d);
183 gezelter 1505
184     // negative sign because this is the vector from j to i:
185    
186 gezelter 1554 Vector3d rj = - *(idat.A2) * *(idat.d);
187 gezelter 1505
188     RealType xi = ri.x();
189     RealType yi = ri.y();
190     RealType zi = ri.z();
191    
192     RealType xj = rj.x();
193     RealType yj = rj.y();
194     RealType zj = rj.z();
195    
196     RealType xi2 = xi * xi;
197     RealType yi2 = yi * yi;
198     RealType zi2 = zi * zi;
199    
200     RealType xj2 = xj * xj;
201     RealType yj2 = yj * yj;
202     RealType zj2 = zj * zj;
203    
204     // calculate the switching info. from the splines
205    
206     RealType s = 0.0;
207     RealType dsdr = 0.0;
208     RealType sp = 0.0;
209     RealType dspdr = 0.0;
210    
211 gezelter 1554 if ( *(idat.rij) < ru) {
212     if ( *(idat.rij) < rl) {
213 gezelter 1505 s = 1.0;
214     dsdr = 0.0;
215     } else {
216     // we are in the switching region
217    
218 gezelter 1554 pair<RealType, RealType> res = mixer.s->getValueAndDerivativeAt(*(idat.rij));
219 gezelter 1505 s = res.first;
220     dsdr = res.second;
221     }
222 gezelter 1485 }
223 gezelter 1505
224 gezelter 1554 if (*(idat.rij) < rup) {
225     if ( *(idat.rij) < rlp) {
226 gezelter 1505 sp = 1.0;
227     dspdr = 0.0;
228     } else {
229     // we are in the switching region
230    
231 gezelter 1554 pair<RealType, RealType> res =mixer.sp->getValueAndDerivativeAt( *(idat.rij));
232 gezelter 1505 sp = res.first;
233     dspdr = res.second;
234     }
235     }
236    
237 gezelter 1834
238 gezelter 1505 RealType wi = 2.0*(xi2-yi2)*zi / r3;
239     RealType wj = 2.0*(xj2-yj2)*zj / r3;
240     RealType w = wi+wj;
241    
242    
243 gezelter 1554 RealType zif = zi/ *(idat.rij) - 0.6;
244     RealType zis = zi/ *(idat.rij) + 0.8;
245 gezelter 1505
246 gezelter 1554 RealType zjf = zj/ *(idat.rij) - 0.6;
247     RealType zjs = zj/ *(idat.rij) + 0.8;
248 gezelter 1505
249     RealType wip = zif*zif*zis*zis - w0;
250     RealType wjp = zjf*zjf*zjs*zjs - w0;
251     RealType wp = wip + wjp;
252    
253     Vector3d dwi(4.0*xi*zi/r3 - 6.0*xi*zi*(xi2-yi2)/r5,
254     - 4.0*yi*zi/r3 - 6.0*yi*zi*(xi2-yi2)/r5,
255     2.0*(xi2-yi2)/r3 - 6.0*zi2*(xi2-yi2)/r5);
256    
257     Vector3d dwj(4.0*xj*zj/r3 - 6.0*xj*zj*(xj2-yj2)/r5,
258     - 4.0*yj*zj/r3 - 6.0*yj*zj*(xj2-yj2)/r5,
259     2.0*(xj2-yj2)/r3 - 6.0*zj2*(xj2-yj2)/r5);
260    
261     RealType uglyi = zif*zif*zis + zif*zis*zis;
262     RealType uglyj = zjf*zjf*zjs + zjf*zjs*zjs;
263 gezelter 1485
264 gezelter 1505 Vector3d dwip(-2.0*xi*zi*uglyi/r3,
265     -2.0*yi*zi*uglyi/r3,
266 gezelter 1554 2.0*(1.0/ *(idat.rij) - zi2/r3)*uglyi);
267 gezelter 1505
268     Vector3d dwjp(-2.0*xj*zj*uglyj/r3,
269     -2.0*yj*zj*uglyj/r3,
270 gezelter 1554 2.0*(1.0/ *(idat.rij) - zj2/r3)*uglyj);
271 gezelter 1505
272     Vector3d dwidu(4.0*(yi*zi2 + 0.5*yi*(xi2-yi2))/r3,
273     4.0*(xi*zi2 - 0.5*xi*(xi2-yi2))/r3,
274     - 8.0*xi*yi*zi/r3);
275    
276     Vector3d dwjdu(4.0*(yj*zj2 + 0.5*yj*(xj2-yj2))/r3,
277     4.0*(xj*zj2 - 0.5*xj*(xj2-yj2))/r3,
278     - 8.0*xj*yj*zj/r3);
279    
280 gezelter 1554 Vector3d dwipdu(2.0*yi*uglyi/ *(idat.rij) ,
281     -2.0*xi*uglyi/ *(idat.rij) ,
282 gezelter 1505 0.0);
283    
284 gezelter 1554 Vector3d dwjpdu(2.0*yj*uglyj/ *(idat.rij) ,
285     -2.0*xj*uglyj/ *(idat.rij) ,
286 gezelter 1505 0.0);
287    
288     if (isPower) {
289 gezelter 1834 cerr << "This is probably an error!\n";
290 gezelter 1505 RealType frac1 = 0.25;
291     RealType frac2 = 0.75;
292     RealType wi2 = wi*wi;
293     RealType wj2 = wj*wj;
294     // sticky power has no w' function:
295     w = frac1 * wi * wi2 + frac2*wi + frac1*wj*wj2 + frac2*wj + v0p;
296     wp = 0.0;
297 gezelter 1668 dwi = frac1*RealType(3.0)*wi2*dwi + frac2*dwi;
298     dwj = frac1*RealType(3.0)*wj2*dwi + frac2*dwi;
299 gezelter 1505 dwip = V3Zero;
300     dwjp = V3Zero;
301 gezelter 1668 dwidu = frac1*RealType(3.0)*wi2*dwidu + frac2*dwidu;
302     dwidu = frac1*RealType(3.0)*wj2*dwjdu + frac2*dwjdu;
303 gezelter 1505 dwipdu = V3Zero;
304     dwjpdu = V3Zero;
305     sp = 0.0;
306 gezelter 1485 dspdr = 0.0;
307     }
308 gezelter 1505
309 gezelter 1668 *(idat.vpair) += RealType(0.5)*(v0*s*w + v0p*sp*wp);
310     (*(idat.pot))[HYDROGENBONDING_FAMILY] += RealType(0.5)*(v0*s*w + v0p*sp*wp)* *(idat.sw) ;
311 gezelter 1505
312     // do the torques first since they are easy:
313     // remember that these are still in the body-fixed axes
314    
315 gezelter 1668 Vector3d ti = RealType(0.5)* *(idat.sw) *(v0*s*dwidu + v0p*sp*dwipdu);
316     Vector3d tj = RealType(0.5)* *(idat.sw) *(v0*s*dwjdu + v0p*sp*dwjpdu);
317 gezelter 1505
318     // go back to lab frame using transpose of rotation matrix:
319    
320 gezelter 1554 *(idat.t1) += A1trans * ti;
321     *(idat.t2) += A2trans * tj;
322 gezelter 1505
323     // Now, on to the forces:
324    
325     // first rotate the i terms back into the lab frame:
326    
327 gezelter 1554 Vector3d radcomi = (v0 * s * dwi + v0p * sp * dwip) * *(idat.sw);
328     Vector3d radcomj = (v0 * s * dwj + v0p * sp * dwjp) * *(idat.sw);
329 gezelter 1505
330     Vector3d fii = A1trans * radcomi;
331     Vector3d fjj = A2trans * radcomj;
332    
333     // now assemble these with the radial-only terms:
334 gezelter 1834
335 gezelter 1668 *(idat.f1) += RealType(0.5) * ((v0*dsdr*w + v0p*dspdr*wp) * *(idat.d) /
336     *(idat.rij) + fii - fjj);
337 gezelter 1505
338 gezelter 1485 }
339     }
340    
341 gezelter 1505 return;
342 gezelter 1485 }
343 gezelter 1505
344 gezelter 1545 RealType Sticky::getSuggestedCutoffRadius(pair<AtomType*, AtomType*> atypes) {
345 gezelter 1505 if (!initialized_) initialize();
346     map<pair<AtomType*, AtomType*>, StickyInteractionData>::iterator it;
347 gezelter 1545 it = MixingMap.find(atypes);
348 gezelter 1505 if (it == MixingMap.end())
349     return 0.0;
350     else {
351     StickyInteractionData mixer = (*it).second;
352     return mixer.rbig;
353     }
354     }
355 gezelter 1485 }

Properties

Name Value
svn:eol-style native