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

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] Kuang & Gezelter, J. Chem. Phys. 133, 164101 (2010).
40 * [5] Vardeman, Stocker & Gezelter, J. Chem. Theory Comput. 7, 834 (2011).
41 */
42
43 #include <stdio.h>
44 #include <string.h>
45
46 #include <cmath>
47 #include "nonbonded/Sticky.hpp"
48 #include "nonbonded/LJ.hpp"
49 #include "types/StickyAdapter.hpp"
50 #include "utils/simError.h"
51
52 using namespace std;
53 namespace OpenMD {
54
55 Sticky::Sticky() : name_("Sticky"), initialized_(false), forceField_(NULL) {}
56
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 StickyAdapter sa = StickyAdapter(at);
70 if (sa.isSticky()) addType(at);
71 }
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 ret = StickyMap.insert( pair<int, AtomType*>(atomType->getIdent(),
81 atomType) );
82 if (ret.second == false) {
83 sprintf( painCave.errMsg,
84 "Sticky already had a previous entry with ident %d\n",
85 atomType->getIdent() );
86 painCave.severity = OPENMD_INFO;
87 painCave.isFatal = 0;
88 simError();
89 }
90
91 StickyAdapter sticky1 = StickyAdapter(atomType);
92
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
100 StickyAdapter sticky2 = StickyAdapter(atype2);
101
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 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 mixer.rbig = max(mixer.ru, mixer.rup);
114 mixer.w0 = sqrt( sticky1.getW0() * sticky2.getW0() );
115 mixer.v0 = sqrt( sticky1.getV0() * sticky2.getV0() );
116 mixer.v0p = sqrt( sticky1.getV0p() * sticky2.getV0p() );
117 mixer.isPower = sticky1.isStickyPower() && sticky2.isStickyPower();
118
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 /**
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 void Sticky::calcForce(InteractionData &idat) {
151
152 if (!initialized_) initialize();
153
154 map<pair<AtomType*, AtomType*>, StickyInteractionData>::iterator it;
155 it = MixingMap.find(idat.atypes);
156 if (it != MixingMap.end()) {
157
158 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 if ( *(idat.rij) <= rbig) {
171
172 RealType r3 = *(idat.r2) * *(idat.rij);
173 RealType r5 = r3 * *(idat.r2);
174
175 RotMat3x3d A1trans = idat.A1->transpose();
176 RotMat3x3d A2trans = idat.A2->transpose();
177
178 // rotate the inter-particle separation into the two different
179 // body-fixed coordinate systems:
180
181 Vector3d ri = *(idat.A1) * *(idat.d);
182
183 // negative sign because this is the vector from j to i:
184
185 Vector3d rj = - *(idat.A2) * *(idat.d);
186
187 RealType xi = ri.x();
188 RealType yi = ri.y();
189 RealType zi = ri.z();
190
191 RealType xj = rj.x();
192 RealType yj = rj.y();
193 RealType zj = rj.z();
194
195 RealType xi2 = xi * xi;
196 RealType yi2 = yi * yi;
197 RealType zi2 = zi * zi;
198
199 RealType xj2 = xj * xj;
200 RealType yj2 = yj * yj;
201 RealType zj2 = zj * zj;
202
203 // calculate the switching info. from the splines
204
205 RealType s = 0.0;
206 RealType dsdr = 0.0;
207 RealType sp = 0.0;
208 RealType dspdr = 0.0;
209
210 if ( *(idat.rij) < ru) {
211 if ( *(idat.rij) < rl) {
212 s = 1.0;
213 dsdr = 0.0;
214 } else {
215 // we are in the switching region
216
217 pair<RealType, RealType> res = mixer.s->getValueAndDerivativeAt(*(idat.rij));
218 s = res.first;
219 dsdr = res.second;
220 }
221 }
222
223 if (*(idat.rij) < rup) {
224 if ( *(idat.rij) < rlp) {
225 sp = 1.0;
226 dspdr = 0.0;
227 } else {
228 // we are in the switching region
229
230 pair<RealType, RealType> res =mixer.sp->getValueAndDerivativeAt( *(idat.rij));
231 sp = res.first;
232 dspdr = res.second;
233 }
234 }
235
236 RealType wi = 2.0*(xi2-yi2)*zi / r3;
237 RealType wj = 2.0*(xj2-yj2)*zj / r3;
238 RealType w = wi+wj;
239
240
241 RealType zif = zi/ *(idat.rij) - 0.6;
242 RealType zis = zi/ *(idat.rij) + 0.8;
243
244 RealType zjf = zj/ *(idat.rij) - 0.6;
245 RealType zjs = zj/ *(idat.rij) + 0.8;
246
247 RealType wip = zif*zif*zis*zis - w0;
248 RealType wjp = zjf*zjf*zjs*zjs - w0;
249 RealType wp = wip + wjp;
250
251 Vector3d dwi(4.0*xi*zi/r3 - 6.0*xi*zi*(xi2-yi2)/r5,
252 - 4.0*yi*zi/r3 - 6.0*yi*zi*(xi2-yi2)/r5,
253 2.0*(xi2-yi2)/r3 - 6.0*zi2*(xi2-yi2)/r5);
254
255 Vector3d dwj(4.0*xj*zj/r3 - 6.0*xj*zj*(xj2-yj2)/r5,
256 - 4.0*yj*zj/r3 - 6.0*yj*zj*(xj2-yj2)/r5,
257 2.0*(xj2-yj2)/r3 - 6.0*zj2*(xj2-yj2)/r5);
258
259 RealType uglyi = zif*zif*zis + zif*zis*zis;
260 RealType uglyj = zjf*zjf*zjs + zjf*zjs*zjs;
261
262 Vector3d dwip(-2.0*xi*zi*uglyi/r3,
263 -2.0*yi*zi*uglyi/r3,
264 2.0*(1.0/ *(idat.rij) - zi2/r3)*uglyi);
265
266 Vector3d dwjp(-2.0*xj*zj*uglyj/r3,
267 -2.0*yj*zj*uglyj/r3,
268 2.0*(1.0/ *(idat.rij) - zj2/r3)*uglyj);
269
270 Vector3d dwidu(4.0*(yi*zi2 + 0.5*yi*(xi2-yi2))/r3,
271 4.0*(xi*zi2 - 0.5*xi*(xi2-yi2))/r3,
272 - 8.0*xi*yi*zi/r3);
273
274 Vector3d dwjdu(4.0*(yj*zj2 + 0.5*yj*(xj2-yj2))/r3,
275 4.0*(xj*zj2 - 0.5*xj*(xj2-yj2))/r3,
276 - 8.0*xj*yj*zj/r3);
277
278 Vector3d dwipdu(2.0*yi*uglyi/ *(idat.rij) ,
279 -2.0*xi*uglyi/ *(idat.rij) ,
280 0.0);
281
282 Vector3d dwjpdu(2.0*yj*uglyj/ *(idat.rij) ,
283 -2.0*xj*uglyj/ *(idat.rij) ,
284 0.0);
285
286 if (isPower) {
287 RealType frac1 = 0.25;
288 RealType frac2 = 0.75;
289 RealType wi2 = wi*wi;
290 RealType wj2 = wj*wj;
291 // sticky power has no w' function:
292 w = frac1 * wi * wi2 + frac2*wi + frac1*wj*wj2 + frac2*wj + v0p;
293 wp = 0.0;
294 dwi = frac1*RealType(3.0)*wi2*dwi + frac2*dwi;
295 dwj = frac1*RealType(3.0)*wj2*dwi + frac2*dwi;
296 dwip = V3Zero;
297 dwjp = V3Zero;
298 dwidu = frac1*RealType(3.0)*wi2*dwidu + frac2*dwidu;
299 dwidu = frac1*RealType(3.0)*wj2*dwjdu + frac2*dwjdu;
300 dwipdu = V3Zero;
301 dwjpdu = V3Zero;
302 sp = 0.0;
303 dspdr = 0.0;
304 }
305
306
307
308 *(idat.vpair) += RealType(0.5)*(v0*s*w + v0p*sp*wp);
309 (*(idat.pot))[HYDROGENBONDING_FAMILY] += RealType(0.5)*(v0*s*w + v0p*sp*wp)* *(idat.sw) ;
310
311 // do the torques first since they are easy:
312 // remember that these are still in the body-fixed axes
313
314 Vector3d ti = RealType(0.5)* *(idat.sw) *(v0*s*dwidu + v0p*sp*dwipdu);
315 Vector3d tj = RealType(0.5)* *(idat.sw) *(v0*s*dwjdu + v0p*sp*dwjpdu);
316
317 // go back to lab frame using transpose of rotation matrix:
318
319 *(idat.t1) += A1trans * ti;
320 *(idat.t2) += A2trans * tj;
321
322 // Now, on to the forces:
323
324 // first rotate the i terms back into the lab frame:
325
326 Vector3d radcomi = (v0 * s * dwi + v0p * sp * dwip) * *(idat.sw);
327 Vector3d radcomj = (v0 * s * dwj + v0p * sp * dwjp) * *(idat.sw);
328
329 Vector3d fii = A1trans * radcomi;
330 Vector3d fjj = A2trans * radcomj;
331
332 // now assemble these with the radial-only terms:
333
334 *(idat.f1) += RealType(0.5) * ((v0*dsdr*w + v0p*dspdr*wp) * *(idat.d) /
335 *(idat.rij) + fii - fjj);
336
337 }
338 }
339
340 return;
341 }
342
343 RealType Sticky::getSuggestedCutoffRadius(pair<AtomType*, AtomType*> atypes) {
344 if (!initialized_) initialize();
345 map<pair<AtomType*, AtomType*>, StickyInteractionData>::iterator it;
346 it = MixingMap.find(atypes);
347 if (it == MixingMap.end())
348 return 0.0;
349 else {
350 StickyInteractionData mixer = (*it).second;
351 return mixer.rbig;
352 }
353 }
354 }

Properties

Name Value
svn:eol-style native