ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/nonbonded/Sticky.cpp
Revision: 1502
Committed: Sat Oct 2 19:53:32 2010 UTC (14 years, 7 months ago) by gezelter
File size: 12463 byte(s)
Log Message:
Many changes, and we're not quite done with this phase yet.

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

Properties

Name Value
svn:eol-style native