OpenMD 3.1
Molecular Dynamics in the Open
Loading...
Searching...
No Matches
Sticky.cpp
1/*
2 * Copyright (c) 2004-present, The University of Notre Dame. All rights
3 * reserved.
4 *
5 * Redistribution and use in source and binary forms, with or without
6 * modification, are permitted provided that the following conditions are met:
7 *
8 * 1. Redistributions of source code must retain the above copyright notice,
9 * this list of conditions and the following disclaimer.
10 *
11 * 2. Redistributions in binary form must reproduce the above copyright notice,
12 * this list of conditions and the following disclaimer in the documentation
13 * and/or other materials provided with the distribution.
14 *
15 * 3. Neither the name of the copyright holder nor the names of its
16 * contributors may be used to endorse or promote products derived from
17 * this software without specific prior written permission.
18 *
19 * THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
20 * AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
21 * IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
22 * ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
23 * LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
24 * CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF
25 * SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS
26 * INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN
27 * CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE)
28 * ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
29 * POSSIBILITY OF SUCH DAMAGE.
30 *
31 * SUPPORT OPEN SCIENCE! If you use OpenMD or its source code in your
32 * research, please cite the appropriate papers when you publish your
33 * work. Good starting points are:
34 *
35 * [1] Meineke, et al., J. Comp. Chem. 26, 252-271 (2005).
36 * [2] Fennell & Gezelter, J. Chem. Phys. 124, 234104 (2006).
37 * [3] Sun, Lin & Gezelter, J. Chem. Phys. 128, 234107 (2008).
38 * [4] Vardeman, Stocker & Gezelter, J. Chem. Theory Comput. 7, 834 (2011).
39 * [5] Kuang & Gezelter, Mol. Phys., 110, 691-701 (2012).
40 * [6] Lamichhane, Gezelter & Newman, J. Chem. Phys. 141, 134109 (2014).
41 * [7] Lamichhane, Newman & Gezelter, J. Chem. Phys. 141, 134110 (2014).
42 * [8] Bhattarai, Newman & Gezelter, Phys. Rev. B 99, 094106 (2019).
43 */
44
45#include "nonbonded/Sticky.hpp"
46
47#include <cmath>
48#include <cstdio>
49#include <cstring>
50
51#include "nonbonded/LJ.hpp"
52#include "types/StickyAdapter.hpp"
53#include "utils/simError.h"
54
55using namespace std;
56namespace OpenMD {
57
58 Sticky::Sticky() : initialized_(false), forceField_(NULL), name_("Sticky") {}
59
60 void Sticky::initialize() {
61 Stypes.clear();
62 Stids.clear();
63 MixingMap.clear();
64 nSticky_ = 0;
65
66 Stids.resize(forceField_->getNAtomType(), -1);
67
68 // Sticky handles all of the Sticky-Sticky interactions
69
70 AtomTypeSet::iterator at;
71 for (at = simTypes_.begin(); at != simTypes_.end(); ++at) {
72 if ((*at)->isSticky()) nSticky_++;
73 }
74
75 MixingMap.resize(nSticky_);
76
77 for (at = simTypes_.begin(); at != simTypes_.end(); ++at) {
78 if ((*at)->isSticky()) addType(*at);
79 }
80
81 initialized_ = true;
82 }
83
84 void Sticky::addType(AtomType* atomType) {
85 StickyAdapter sticky1 = StickyAdapter(atomType);
86
87 // add it to the map:
88
89 int atid = atomType->getIdent();
90 int stid = Stypes.size();
91
92 pair<set<int>::iterator, bool> ret;
93 ret = Stypes.insert(atid);
94 if (ret.second == false) {
95 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
96 "Sticky already had a previous entry with ident %d\n", atid);
97 painCave.severity = OPENMD_INFO;
98 painCave.isFatal = 0;
99 simError();
100 }
101
102 Stids[atid] = stid;
103 MixingMap[stid].resize(nSticky_);
104
105 // Now, iterate over all known types and add to the mixing map:
106
107 std::set<int>::iterator it;
108 for (it = Stypes.begin(); it != Stypes.end(); ++it) {
109 int stid2 = Stids[(*it)];
110 AtomType* atype2 = forceField_->getAtomType((*it));
111 StickyAdapter sticky2 = StickyAdapter(atype2);
112
113 StickyInteractionData mixer;
114
115 // Mixing two different sticky types is silly, but if you want 2
116 // sticky types in your simulation, we'll let you do it with the
117 // Lorentz- Berthelot mixing rules (which happen to do the right thing
118 // when atomType and atype2 happen to be the same.
119
120 mixer.rl = 0.5 * (sticky1.getRl() + sticky2.getRl());
121 mixer.ru = 0.5 * (sticky1.getRu() + sticky2.getRu());
122 mixer.rlp = 0.5 * (sticky1.getRlp() + sticky2.getRlp());
123 mixer.rup = 0.5 * (sticky1.getRup() + sticky2.getRup());
124 mixer.rbig = max(mixer.ru, mixer.rup);
125 mixer.w0 = sqrt(sticky1.getW0() * sticky2.getW0());
126 mixer.v0 = sqrt(sticky1.getV0() * sticky2.getV0());
127 mixer.v0p = sqrt(sticky1.getV0p() * sticky2.getV0p());
128 mixer.isPower = sticky1.isStickyPower() && sticky2.isStickyPower();
129
130 CubicSpline* s = new CubicSpline();
131 s->addPoint(mixer.rl, 1.0);
132 s->addPoint(mixer.ru, 0.0);
133 mixer.s = s;
134
135 CubicSpline* sp = new CubicSpline();
136 sp->addPoint(mixer.rlp, 1.0);
137 sp->addPoint(mixer.rup, 0.0);
138 mixer.sp = sp;
139
140 MixingMap[stid2].resize(nSticky_);
141
142 MixingMap[stid][stid2] = mixer;
143 if (stid2 != stid) { MixingMap[stid2][stid] = mixer; }
144 }
145 }
146
147 /**
148 * This function does the sticky portion of the SSD potential
149 * [Chandra and Ichiye, Journal of Chemical Physics 111, 2701
150 * (1999)]. The Lennard-Jones and dipolar interaction must be
151 * handled separately. We assume that the rotation matrices have
152 * already been calculated and placed in the A1 & A2 entries in the
153 * idat structure.
154 */
155
156 void Sticky::calcForce(InteractionData& idat) {
157 if (!initialized_) initialize();
158
159 StickyInteractionData& mixer =
160 MixingMap[Stids[idat.atid1]][Stids[idat.atid2]];
161
162 RealType w0 = mixer.w0;
163 RealType v0 = mixer.v0;
164 RealType v0p = mixer.v0p;
165 RealType rl = mixer.rl;
166 RealType ru = mixer.ru;
167 RealType rlp = mixer.rlp;
168 RealType rup = mixer.rup;
169 RealType rbig = mixer.rbig;
170 bool isPower = mixer.isPower;
171
172 if (idat.rij <= rbig) {
173 RealType r3 = idat.r2 * idat.rij;
174 RealType r5 = r3 * idat.r2;
175
176 RotMat3x3d A1trans = idat.A1.transpose();
177 RotMat3x3d A2trans = idat.A2.transpose();
178
179 // rotate the inter-particle separation into the two different
180 // body-fixed coordinate systems:
181
182 Vector3d ri = idat.A1 * idat.d;
183
184 // negative sign because this is the vector from j to i:
185
186 Vector3d rj = -idat.A2 * idat.d;
187
188 RealType xi = ri.x();
189 RealType yi = ri.y();
190 RealType zi = ri.z();
191
192 RealType xj = rj.x();
193 RealType yj = rj.y();
194 RealType zj = rj.z();
195
196 RealType xi2 = xi * xi;
197 RealType yi2 = yi * yi;
198 RealType zi2 = zi * zi;
199
200 RealType xj2 = xj * xj;
201 RealType yj2 = yj * yj;
202 RealType zj2 = zj * zj;
203
204 // calculate the switching info. from the splines
205
206 RealType s = 0.0;
207 RealType dsdr = 0.0;
208 RealType sp = 0.0;
209 RealType dspdr = 0.0;
210
211 if (idat.rij < ru) {
212 if (idat.rij < rl) {
213 s = 1.0;
214 dsdr = 0.0;
215 } else {
216 // we are in the switching region
217 mixer.s->getValueAndDerivativeAt(idat.rij, s, dsdr);
218 }
219 }
220
221 if (idat.rij < rup) {
222 if (idat.rij < rlp) {
223 sp = 1.0;
224 dspdr = 0.0;
225 } else {
226 // we are in the switching region
227 mixer.sp->getValueAndDerivativeAt(idat.rij, sp, dspdr);
228 }
229 }
230
231 RealType wi = 2.0 * (xi2 - yi2) * zi / r3;
232 RealType wj = 2.0 * (xj2 - yj2) * zj / r3;
233 RealType w = wi + wj;
234
235 RealType zif = zi / idat.rij - 0.6;
236 RealType zis = zi / idat.rij + 0.8;
237
238 RealType zjf = zj / idat.rij - 0.6;
239 RealType zjs = zj / idat.rij + 0.8;
240
241 RealType wip = zif * zif * zis * zis - w0;
242 RealType wjp = zjf * zjf * zjs * zjs - w0;
243 RealType wp = wip + wjp;
244
245 Vector3d dwi(4.0 * xi * zi / r3 - 6.0 * xi * zi * (xi2 - yi2) / r5,
246 -4.0 * yi * zi / r3 - 6.0 * yi * zi * (xi2 - yi2) / r5,
247 2.0 * (xi2 - yi2) / r3 - 6.0 * zi2 * (xi2 - yi2) / r5);
248
249 Vector3d dwj(4.0 * xj * zj / r3 - 6.0 * xj * zj * (xj2 - yj2) / r5,
250 -4.0 * yj * zj / r3 - 6.0 * yj * zj * (xj2 - yj2) / r5,
251 2.0 * (xj2 - yj2) / r3 - 6.0 * zj2 * (xj2 - yj2) / r5);
252
253 RealType uglyi = zif * zif * zis + zif * zis * zis;
254 RealType uglyj = zjf * zjf * zjs + zjf * zjs * zjs;
255
256 Vector3d dwip(-2.0 * xi * zi * uglyi / r3, -2.0 * yi * zi * uglyi / r3,
257 2.0 * (1.0 / idat.rij - zi2 / r3) * uglyi);
258
259 Vector3d dwjp(-2.0 * xj * zj * uglyj / r3, -2.0 * yj * zj * uglyj / r3,
260 2.0 * (1.0 / idat.rij - zj2 / r3) * uglyj);
261
262 Vector3d dwidu(4.0 * (yi * zi2 + 0.5 * yi * (xi2 - yi2)) / r3,
263 4.0 * (xi * zi2 - 0.5 * xi * (xi2 - yi2)) / r3,
264 -8.0 * xi * yi * zi / r3);
265
266 Vector3d dwjdu(4.0 * (yj * zj2 + 0.5 * yj * (xj2 - yj2)) / r3,
267 4.0 * (xj * zj2 - 0.5 * xj * (xj2 - yj2)) / r3,
268 -8.0 * xj * yj * zj / r3);
269
270 Vector3d dwipdu(2.0 * yi * uglyi / idat.rij, -2.0 * xi * uglyi / idat.rij,
271 0.0);
272
273 Vector3d dwjpdu(2.0 * yj * uglyj / idat.rij, -2.0 * xj * uglyj / idat.rij,
274 0.0);
275
276 if (isPower) {
277 cerr << "This is probably an error!\n";
278 RealType frac1 = 0.25;
279 RealType frac2 = 0.75;
280 RealType wi2 = wi * wi;
281 RealType wj2 = wj * wj;
282 // sticky power has no w' function:
283 w = frac1 * wi * wi2 + frac2 * wi + frac1 * wj * wj2 + frac2 * wj + v0p;
284 wp = 0.0;
285 dwi = frac1 * RealType(3.0) * wi2 * dwi + frac2 * dwi;
286 dwj = frac1 * RealType(3.0) * wj2 * dwi + frac2 * dwi;
287 dwip = V3Zero;
288 dwjp = V3Zero;
289 dwidu = frac1 * RealType(3.0) * wi2 * dwidu + frac2 * dwidu;
290 dwidu = frac1 * RealType(3.0) * wj2 * dwjdu + frac2 * dwjdu;
291 dwipdu = V3Zero;
292 dwjpdu = V3Zero;
293 sp = 0.0;
294 dspdr = 0.0;
295 }
296
297 idat.vpair += 0.5 * (v0 * s * w + v0p * sp * wp);
299 0.5 * (v0 * s * w + v0p * sp * wp) * idat.sw;
300 if (idat.isSelected)
302 0.5 * (v0 * s * w + v0p * sp * wp) * idat.sw;
303
304 // do the torques first since they are easy:
305 // remember that these are still in the body-fixed axes
306
307 Vector3d ti = 0.5 * idat.sw * (v0 * s * dwidu + v0p * sp * dwipdu);
308 Vector3d tj = 0.5 * idat.sw * (v0 * s * dwjdu + v0p * sp * dwjpdu);
309
310 // go back to lab frame using transpose of rotation matrix:
311
312 idat.t1 += A1trans * ti;
313 idat.t2 += A2trans * tj;
314
315 // Now, on to the forces:
316
317 // first rotate the i terms back into the lab frame:
318
319 Vector3d radcomi = (v0 * s * dwi + v0p * sp * dwip) * idat.sw;
320 Vector3d radcomj = (v0 * s * dwj + v0p * sp * dwjp) * idat.sw;
321
322 Vector3d fii = A1trans * radcomi;
323 Vector3d fjj = A2trans * radcomj;
324
325 // now assemble these with the radial-only terms:
326
327 idat.f1 += 0.5 * ((v0 * dsdr * w + v0p * dspdr * wp) * idat.d / idat.rij +
328 fii - fjj);
329 }
330
331 return;
332 }
333
334 RealType Sticky::getSuggestedCutoffRadius(pair<AtomType*, AtomType*> atypes) {
335 if (!initialized_) initialize();
336 int atid1 = atypes.first->getIdent();
337 int atid2 = atypes.second->getIdent();
338 int stid1 = Stids[atid1];
339 int stid2 = Stids[atid2];
340
341 if (stid1 == -1 || stid2 == -1)
342 return 0.0;
343 else { return MixingMap[stid1][stid2].rbig; }
344 }
345} // namespace OpenMD
Real & z()
Returns reference of the third element of Vector3.
Definition Vector3.hpp:120
Real & x()
Returns reference of the first element of Vector3.
Definition Vector3.hpp:96
Real & y()
Returns reference of the second element of Vector3.
Definition Vector3.hpp:108
This basic Periodic Table class was originally taken from the data.cpp file in OpenBabel.
@ HYDROGENBONDING_FAMILY
Short-range directional interactions.
The InteractionData struct.
Vector3d d
interatomic vector (already wrapped into box)
int atid1
atomType ident for atom 1
potVec selePot
potential energies of the selected sites
RealType sw
switching function value at rij
RotMat3x3d A2
rotation matrix of second atom
bool isSelected
one of the particles has been selected for selection potential
RotMat3x3d A1
rotation matrix of first atom
Vector3d t1
torque on first atom
Vector3d t2
torque on second atom
RealType vpair
pair potential
int atid2
atomType ident for atom 2
Vector3d f1
force between the two atoms
RealType rij
interatomic separation