ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/nonbonded/Sticky.cpp
Revision: 1554
Committed: Sat Apr 30 02:54:02 2011 UTC (14 years, 2 months ago) by gezelter
File size: 13671 byte(s)
Log Message:
For efficiency, pointers instead of objects will be passed during main
force loop. 


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     * [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 gezelter 1502 Sticky::Sticky() : name_("Sticky"), initialized_(false), forceField_(NULL) {}
54 gezelter 1485
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 gezelter 1502 /**
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 gezelter 1536 void Sticky::calcForce(InteractionData &idat) {
190 gezelter 1502
191 gezelter 1485 if (!initialized_) initialize();
192    
193 gezelter 1505 map<pair<AtomType*, AtomType*>, StickyInteractionData>::iterator it;
194 gezelter 1554 it = MixingMap.find(*(idat.atypes));
195 gezelter 1505 if (it != MixingMap.end()) {
196 gezelter 1485
197 gezelter 1505 StickyInteractionData mixer = (*it).second;
198    
199     RealType w0 = mixer.w0;
200     RealType v0 = mixer.v0;
201     RealType v0p = mixer.v0p;
202     RealType rl = mixer.rl;
203     RealType ru = mixer.ru;
204     RealType rlp = mixer.rlp;
205     RealType rup = mixer.rup;
206     RealType rbig = mixer.rbig;
207     bool isPower = mixer.isPower;
208    
209 gezelter 1554 if ( *(idat.rij) <= rbig) {
210 gezelter 1505
211 gezelter 1554 RealType r3 = *(idat.r2) * *(idat.rij);
212     RealType r5 = r3 * *(idat.r2);
213 gezelter 1505
214 gezelter 1554 RotMat3x3d A1trans = idat.A1->transpose();
215     RotMat3x3d A2trans = idat.A2->transpose();
216 gezelter 1505
217     // rotate the inter-particle separation into the two different
218     // body-fixed coordinate systems:
219    
220 gezelter 1554 Vector3d ri = *(idat.A1) * *(idat.d);
221 gezelter 1505
222     // negative sign because this is the vector from j to i:
223    
224 gezelter 1554 Vector3d rj = - *(idat.A2) * *(idat.d);
225 gezelter 1505
226     RealType xi = ri.x();
227     RealType yi = ri.y();
228     RealType zi = ri.z();
229    
230     RealType xj = rj.x();
231     RealType yj = rj.y();
232     RealType zj = rj.z();
233    
234     RealType xi2 = xi * xi;
235     RealType yi2 = yi * yi;
236     RealType zi2 = zi * zi;
237    
238     RealType xj2 = xj * xj;
239     RealType yj2 = yj * yj;
240     RealType zj2 = zj * zj;
241    
242     // calculate the switching info. from the splines
243    
244     RealType s = 0.0;
245     RealType dsdr = 0.0;
246     RealType sp = 0.0;
247     RealType dspdr = 0.0;
248    
249 gezelter 1554 if ( *(idat.rij) < ru) {
250     if ( *(idat.rij) < rl) {
251 gezelter 1505 s = 1.0;
252     dsdr = 0.0;
253     } else {
254     // we are in the switching region
255    
256 gezelter 1554 pair<RealType, RealType> res = mixer.s->getValueAndDerivativeAt(*(idat.rij));
257 gezelter 1505 s = res.first;
258     dsdr = res.second;
259     }
260 gezelter 1485 }
261 gezelter 1505
262 gezelter 1554 if (*(idat.rij) < rup) {
263     if ( *(idat.rij) < rlp) {
264 gezelter 1505 sp = 1.0;
265     dspdr = 0.0;
266     } else {
267     // we are in the switching region
268    
269 gezelter 1554 pair<RealType, RealType> res =mixer.sp->getValueAndDerivativeAt( *(idat.rij));
270 gezelter 1505 sp = res.first;
271     dspdr = res.second;
272     }
273     }
274    
275     RealType wi = 2.0*(xi2-yi2)*zi / r3;
276     RealType wj = 2.0*(xj2-yj2)*zj / r3;
277     RealType w = wi+wj;
278    
279    
280 gezelter 1554 RealType zif = zi/ *(idat.rij) - 0.6;
281     RealType zis = zi/ *(idat.rij) + 0.8;
282 gezelter 1505
283 gezelter 1554 RealType zjf = zj/ *(idat.rij) - 0.6;
284     RealType zjs = zj/ *(idat.rij) + 0.8;
285 gezelter 1505
286     RealType wip = zif*zif*zis*zis - w0;
287     RealType wjp = zjf*zjf*zjs*zjs - w0;
288     RealType wp = wip + wjp;
289    
290     Vector3d dwi(4.0*xi*zi/r3 - 6.0*xi*zi*(xi2-yi2)/r5,
291     - 4.0*yi*zi/r3 - 6.0*yi*zi*(xi2-yi2)/r5,
292     2.0*(xi2-yi2)/r3 - 6.0*zi2*(xi2-yi2)/r5);
293    
294     Vector3d dwj(4.0*xj*zj/r3 - 6.0*xj*zj*(xj2-yj2)/r5,
295     - 4.0*yj*zj/r3 - 6.0*yj*zj*(xj2-yj2)/r5,
296     2.0*(xj2-yj2)/r3 - 6.0*zj2*(xj2-yj2)/r5);
297    
298     RealType uglyi = zif*zif*zis + zif*zis*zis;
299     RealType uglyj = zjf*zjf*zjs + zjf*zjs*zjs;
300 gezelter 1485
301 gezelter 1505 Vector3d dwip(-2.0*xi*zi*uglyi/r3,
302     -2.0*yi*zi*uglyi/r3,
303 gezelter 1554 2.0*(1.0/ *(idat.rij) - zi2/r3)*uglyi);
304 gezelter 1505
305     Vector3d dwjp(-2.0*xj*zj*uglyj/r3,
306     -2.0*yj*zj*uglyj/r3,
307 gezelter 1554 2.0*(1.0/ *(idat.rij) - zj2/r3)*uglyj);
308 gezelter 1505
309     Vector3d dwidu(4.0*(yi*zi2 + 0.5*yi*(xi2-yi2))/r3,
310     4.0*(xi*zi2 - 0.5*xi*(xi2-yi2))/r3,
311     - 8.0*xi*yi*zi/r3);
312    
313     Vector3d dwjdu(4.0*(yj*zj2 + 0.5*yj*(xj2-yj2))/r3,
314     4.0*(xj*zj2 - 0.5*xj*(xj2-yj2))/r3,
315     - 8.0*xj*yj*zj/r3);
316    
317 gezelter 1554 Vector3d dwipdu(2.0*yi*uglyi/ *(idat.rij) ,
318     -2.0*xi*uglyi/ *(idat.rij) ,
319 gezelter 1505 0.0);
320    
321 gezelter 1554 Vector3d dwjpdu(2.0*yj*uglyj/ *(idat.rij) ,
322     -2.0*xj*uglyj/ *(idat.rij) ,
323 gezelter 1505 0.0);
324    
325     if (isPower) {
326     RealType frac1 = 0.25;
327     RealType frac2 = 0.75;
328     RealType wi2 = wi*wi;
329     RealType wj2 = wj*wj;
330     // sticky power has no w' function:
331     w = frac1 * wi * wi2 + frac2*wi + frac1*wj*wj2 + frac2*wj + v0p;
332     wp = 0.0;
333     dwi = frac1*3.0*wi2*dwi + frac2*dwi;
334     dwj = frac1*3.0*wj2*dwi + frac2*dwi;
335     dwip = V3Zero;
336     dwjp = V3Zero;
337     dwidu = frac1*3.0*wi2*dwidu + frac2*dwidu;
338     dwidu = frac1*3.0*wj2*dwjdu + frac2*dwjdu;
339     dwipdu = V3Zero;
340     dwjpdu = V3Zero;
341     sp = 0.0;
342 gezelter 1485 dspdr = 0.0;
343     }
344 gezelter 1505
345 gezelter 1554 *(idat.vpair) += 0.5*(v0*s*w + v0p*sp*wp);
346     idat.pot[HYDROGENBONDING_FAMILY] += 0.5*(v0*s*w + v0p*sp*wp)* *(idat.sw) ;
347 gezelter 1505
348     // do the torques first since they are easy:
349     // remember that these are still in the body-fixed axes
350    
351 gezelter 1554 Vector3d ti = 0.5* *(idat.sw) *(v0*s*dwidu + v0p*sp*dwipdu);
352     Vector3d tj = 0.5* *(idat.sw) *(v0*s*dwjdu + v0p*sp*dwjpdu);
353 gezelter 1505
354     // go back to lab frame using transpose of rotation matrix:
355    
356 gezelter 1554 *(idat.t1) += A1trans * ti;
357     *(idat.t2) += A2trans * tj;
358 gezelter 1505
359     // Now, on to the forces:
360    
361     // first rotate the i terms back into the lab frame:
362    
363 gezelter 1554 Vector3d radcomi = (v0 * s * dwi + v0p * sp * dwip) * *(idat.sw);
364     Vector3d radcomj = (v0 * s * dwj + v0p * sp * dwjp) * *(idat.sw);
365 gezelter 1505
366     Vector3d fii = A1trans * radcomi;
367     Vector3d fjj = A2trans * radcomj;
368    
369     // now assemble these with the radial-only terms:
370    
371 gezelter 1554 *(idat.f1) += 0.5 * ((v0*dsdr*w + v0p*dspdr*wp) * *(idat.d) /
372     *(idat.rij) + fii - fjj);
373 gezelter 1505
374 gezelter 1485 }
375     }
376    
377 gezelter 1505 return;
378 gezelter 1485 }
379 gezelter 1505
380 gezelter 1545 RealType Sticky::getSuggestedCutoffRadius(pair<AtomType*, AtomType*> atypes) {
381 gezelter 1505 if (!initialized_) initialize();
382     map<pair<AtomType*, AtomType*>, StickyInteractionData>::iterator it;
383 gezelter 1545 it = MixingMap.find(atypes);
384 gezelter 1505 if (it == MixingMap.end())
385     return 0.0;
386     else {
387     StickyInteractionData mixer = (*it).second;
388     return mixer.rbig;
389     }
390     }
391 gezelter 1485 }

Properties

Name Value
svn:eol-style native