ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/nonbonded/Sticky.cpp
Revision: 1668
Committed: Fri Jan 6 19:03:05 2012 UTC (13 years, 4 months ago) by gezelter
File size: 13839 byte(s)
Log Message:
Some fixes for CMake and single precision builds

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

Properties

Name Value
svn:eol-style native