ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/nonbonded/Sticky.cpp
Revision: 2071
Committed: Sat Mar 7 21:41:51 2015 UTC (10 years, 2 months ago) by gezelter
File size: 11924 byte(s)
Log Message:
Reducing the number of warnings when using g++ to compile.

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, 234107 (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() : initialized_(false), forceField_(NULL), name_("Sticky") {}
56
57 void Sticky::initialize() {
58 Stypes.clear();
59 Stids.clear();
60 MixingMap.clear();
61 nSticky_=0;
62
63 Stids.resize( forceField_->getNAtomType(), -1);
64
65 // Sticky handles all of the Sticky-Sticky interactions
66
67 set<AtomType*>::iterator at;
68 for (at = simTypes_.begin(); at != simTypes_.end(); ++at) {
69 if ((*at)->isSticky()) nSticky_++;
70 }
71
72 MixingMap.resize(nSticky_);
73
74 for (at = simTypes_.begin(); at != simTypes_.end(); ++at) {
75 if ((*at)->isSticky()) addType( *at );
76 }
77
78 initialized_ = true;
79 }
80
81 void Sticky::addType(AtomType* atomType){
82 StickyAdapter sticky1 = StickyAdapter(atomType);
83
84 // add it to the map:
85
86 int atid = atomType->getIdent();
87 int stid = Stypes.size();
88
89 pair<set<int>::iterator,bool> ret;
90 ret = Stypes.insert( atid );
91 if (ret.second == false) {
92 sprintf( painCave.errMsg,
93 "Sticky already had a previous entry with ident %d\n",
94 atid) ;
95 painCave.severity = OPENMD_INFO;
96 painCave.isFatal = 0;
97 simError();
98 }
99
100 Stids[atid] = stid;
101 MixingMap[stid].resize( nSticky_ );
102
103
104 // Now, iterate over all known types and add to the mixing map:
105
106 std::set<int>::iterator it;
107 for( it = Stypes.begin(); it != Stypes.end(); ++it) {
108
109
110 int stid2 = Stids[ (*it) ];
111 AtomType* atype2 = forceField_->getAtomType( (*it) );
112 StickyAdapter sticky2 = StickyAdapter(atype2);
113
114 StickyInteractionData mixer;
115
116 // Mixing two different sticky types is silly, but if you want 2
117 // sticky types in your simulation, we'll let you do it with the
118 // Lorentz- Berthelot mixing rules (which happen to do the right thing
119 // when atomType and atype2 happen to be the same.
120
121 mixer.rl = 0.5 * ( sticky1.getRl() + sticky2.getRl() );
122 mixer.ru = 0.5 * ( sticky1.getRu() + sticky2.getRu() );
123 mixer.rlp = 0.5 * ( sticky1.getRlp() + sticky2.getRlp() );
124 mixer.rup = 0.5 * ( sticky1.getRup() + sticky2.getRup() );
125 mixer.rbig = max(mixer.ru, mixer.rup);
126 mixer.w0 = sqrt( sticky1.getW0() * sticky2.getW0() );
127 mixer.v0 = sqrt( sticky1.getV0() * sticky2.getV0() );
128 mixer.v0p = sqrt( sticky1.getV0p() * sticky2.getV0p() );
129 mixer.isPower = sticky1.isStickyPower() && sticky2.isStickyPower();
130
131 CubicSpline* s = new CubicSpline();
132 s->addPoint(mixer.rl, 1.0);
133 s->addPoint(mixer.ru, 0.0);
134 mixer.s = s;
135
136 CubicSpline* sp = new CubicSpline();
137 sp->addPoint(mixer.rlp, 1.0);
138 sp->addPoint(mixer.rup, 0.0);
139 mixer.sp = sp;
140
141 MixingMap[stid2].resize( nSticky_ );
142
143 MixingMap[stid][stid2] = mixer;
144 if (stid2 != stid) {
145 MixingMap[stid2][stid] = mixer;
146 }
147 }
148 }
149
150 /**
151 * This function does the sticky portion of the SSD potential
152 * [Chandra and Ichiye, Journal of Chemical Physics 111, 2701
153 * (1999)]. The Lennard-Jones and dipolar interaction must be
154 * handled separately. We assume that the rotation matrices have
155 * already been calculated and placed in the A1 & A2 entries in the
156 * idat structure.
157 */
158
159 void Sticky::calcForce(InteractionData &idat) {
160
161 if (!initialized_) initialize();
162
163 StickyInteractionData &mixer = MixingMap[Stids[idat.atid1]][Stids[idat.atid2]];
164
165 RealType w0 = mixer.w0;
166 RealType v0 = mixer.v0;
167 RealType v0p = mixer.v0p;
168 RealType rl = mixer.rl;
169 RealType ru = mixer.ru;
170 RealType rlp = mixer.rlp;
171 RealType rup = mixer.rup;
172 RealType rbig = mixer.rbig;
173 bool isPower = mixer.isPower;
174
175 if ( *(idat.rij) <= rbig) {
176
177 RealType r3 = *(idat.r2) * *(idat.rij);
178 RealType r5 = r3 * *(idat.r2);
179
180 RotMat3x3d A1trans = idat.A1->transpose();
181 RotMat3x3d A2trans = idat.A2->transpose();
182
183 // rotate the inter-particle separation into the two different
184 // body-fixed coordinate systems:
185
186 Vector3d ri = *(idat.A1) * *(idat.d);
187
188 // negative sign because this is the vector from j to i:
189
190 Vector3d rj = - *(idat.A2) * *(idat.d);
191
192 RealType xi = ri.x();
193 RealType yi = ri.y();
194 RealType zi = ri.z();
195
196 RealType xj = rj.x();
197 RealType yj = rj.y();
198 RealType zj = rj.z();
199
200 RealType xi2 = xi * xi;
201 RealType yi2 = yi * yi;
202 RealType zi2 = zi * zi;
203
204 RealType xj2 = xj * xj;
205 RealType yj2 = yj * yj;
206 RealType zj2 = zj * zj;
207
208 // calculate the switching info. from the splines
209
210 RealType s = 0.0;
211 RealType dsdr = 0.0;
212 RealType sp = 0.0;
213 RealType dspdr = 0.0;
214
215 if ( *(idat.rij) < ru) {
216 if ( *(idat.rij) < rl) {
217 s = 1.0;
218 dsdr = 0.0;
219 } else {
220 // we are in the switching region
221 mixer.s->getValueAndDerivativeAt(*(idat.rij), s, dsdr);
222 }
223 }
224
225 if (*(idat.rij) < rup) {
226 if ( *(idat.rij) < rlp) {
227 sp = 1.0;
228 dspdr = 0.0;
229 } else {
230 // we are in the switching region
231 mixer.sp->getValueAndDerivativeAt( *(idat.rij), sp, dspdr);
232 }
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 cerr << "This is probably an error!\n";
288 RealType frac1 = 0.25;
289 RealType frac2 = 0.75;
290 RealType wi2 = wi*wi;
291 RealType wj2 = wj*wj;
292 // sticky power has no w' function:
293 w = frac1 * wi * wi2 + frac2*wi + frac1*wj*wj2 + frac2*wj + v0p;
294 wp = 0.0;
295 dwi = frac1*RealType(3.0)*wi2*dwi + frac2*dwi;
296 dwj = frac1*RealType(3.0)*wj2*dwi + frac2*dwi;
297 dwip = V3Zero;
298 dwjp = V3Zero;
299 dwidu = frac1*RealType(3.0)*wi2*dwidu + frac2*dwidu;
300 dwidu = frac1*RealType(3.0)*wj2*dwjdu + frac2*dwjdu;
301 dwipdu = V3Zero;
302 dwjpdu = V3Zero;
303 sp = 0.0;
304 dspdr = 0.0;
305 }
306
307
308
309 *(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
312 // do the torques first since they are easy:
313 // remember that these are still in the body-fixed axes
314
315 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
318 // go back to lab frame using transpose of rotation matrix:
319
320 *(idat.t1) += A1trans * ti;
321 *(idat.t2) += A2trans * tj;
322
323 // Now, on to the forces:
324
325 // first rotate the i terms back into the lab frame:
326
327 Vector3d radcomi = (v0 * s * dwi + v0p * sp * dwip) * *(idat.sw);
328 Vector3d radcomj = (v0 * s * dwj + v0p * sp * dwjp) * *(idat.sw);
329
330 Vector3d fii = A1trans * radcomi;
331 Vector3d fjj = A2trans * radcomj;
332
333 // now assemble these with the radial-only terms:
334
335 *(idat.f1) += RealType(0.5) * ((v0*dsdr*w + v0p*dspdr*wp) * *(idat.d) /
336 *(idat.rij) + fii - fjj);
337
338 }
339
340 return;
341 }
342
343 RealType Sticky::getSuggestedCutoffRadius(pair<AtomType*, AtomType*> atypes) {
344 if (!initialized_) initialize();
345 int atid1 = atypes.first->getIdent();
346 int atid2 = atypes.second->getIdent();
347 int stid1 = Stids[atid1];
348 int stid2 = Stids[atid2];
349
350 if (stid1 == -1 || stid2 == -1) return 0.0;
351 else {
352 return MixingMap[stid1][stid2].rbig;
353 }
354 }
355 }

Properties

Name Value
svn:eol-style native