ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/devel_omp/src/nonbonded/Sticky.cpp
Revision: 1614
Committed: Tue Aug 23 20:55:51 2011 UTC (13 years, 8 months ago) by mciznick
File size: 13743 byte(s)
Log Message:
Updated scalability of OpenMP threads.

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

Properties

Name Value
svn:eol-style native