-
Notifications
You must be signed in to change notification settings - Fork 125
/
TwoStepLangevin.cc
433 lines (365 loc) · 17 KB
/
TwoStepLangevin.cc
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
// Copyright (c) 2009-2024 The Regents of the University of Michigan.
// Part of HOOMD-blue, released under the BSD 3-Clause License.
#include "TwoStepLangevin.h"
#include "hoomd/RNGIdentifiers.h"
#include "hoomd/RandomNumbers.h"
#include "hoomd/VectorMath.h"
#ifdef ENABLE_MPI
#include "hoomd/HOOMDMPI.h"
#endif
using namespace std;
using namespace hoomd;
namespace hoomd
{
namespace md
{
TwoStepLangevin::TwoStepLangevin(std::shared_ptr<SystemDefinition> sysdef,
std::shared_ptr<ParticleGroup> group,
std::shared_ptr<Variant> T)
: TwoStepLangevinBase(sysdef, group, T), m_reservoir_energy(0), m_extra_energy_overdeltaT(0),
m_tally(false), m_noiseless_t(false), m_noiseless_r(false)
{
m_exec_conf->msg->notice(5) << "Constructing TwoStepLangevin" << endl;
}
TwoStepLangevin::~TwoStepLangevin()
{
m_exec_conf->msg->notice(5) << "Destroying TwoStepLangevin" << endl;
}
/*! \param timestep Current time step
\post Particle positions are moved forward to timestep+1 and velocities to timestep+1/2 per the
velocity verlet method.
*/
void TwoStepLangevin::integrateStepOne(uint64_t timestep)
{
unsigned int group_size = m_group->getNumMembers();
ArrayHandle<Scalar4> h_vel(m_pdata->getVelocities(),
access_location::host,
access_mode::readwrite);
ArrayHandle<Scalar3> h_accel(m_pdata->getAccelerations(),
access_location::host,
access_mode::readwrite);
ArrayHandle<Scalar4> h_pos(m_pdata->getPositions(),
access_location::host,
access_mode::readwrite);
ArrayHandle<int3> h_image(m_pdata->getImages(), access_location::host, access_mode::readwrite);
ArrayHandle<Scalar3> h_gamma_r(m_gamma_r, access_location::host, access_mode::read);
const BoxDim& box = m_pdata->getBox();
// perform the first half step of velocity verlet
// r(t+deltaT) = r(t) + v(t)*deltaT + (1/2)a(t)*deltaT^2
// v(t+deltaT/2) = v(t) + (1/2)a*deltaT
for (unsigned int group_idx = 0; group_idx < group_size; group_idx++)
{
unsigned int j = m_group->getMemberIndex(group_idx);
Scalar dx = h_vel.data[j].x * m_deltaT
+ Scalar(1.0 / 2.0) * h_accel.data[j].x * m_deltaT * m_deltaT;
Scalar dy = h_vel.data[j].y * m_deltaT
+ Scalar(1.0 / 2.0) * h_accel.data[j].y * m_deltaT * m_deltaT;
Scalar dz = h_vel.data[j].z * m_deltaT
+ Scalar(1.0 / 2.0) * h_accel.data[j].z * m_deltaT * m_deltaT;
h_pos.data[j].x += dx;
h_pos.data[j].y += dy;
h_pos.data[j].z += dz;
// particles may have been moved slightly outside the box by the above steps, wrap them back
// into place
box.wrap(h_pos.data[j], h_image.data[j]);
h_vel.data[j].x += Scalar(1.0 / 2.0) * h_accel.data[j].x * m_deltaT;
h_vel.data[j].y += Scalar(1.0 / 2.0) * h_accel.data[j].y * m_deltaT;
h_vel.data[j].z += Scalar(1.0 / 2.0) * h_accel.data[j].z * m_deltaT;
}
if (m_aniso)
{
ArrayHandle<Scalar4> h_orientation(m_pdata->getOrientationArray(),
access_location::host,
access_mode::readwrite);
ArrayHandle<Scalar4> h_angmom(m_pdata->getAngularMomentumArray(),
access_location::host,
access_mode::readwrite);
ArrayHandle<Scalar4> h_net_torque(m_pdata->getNetTorqueArray(),
access_location::host,
access_mode::read);
ArrayHandle<Scalar3> h_inertia(m_pdata->getMomentsOfInertiaArray(),
access_location::host,
access_mode::read);
for (unsigned int group_idx = 0; group_idx < group_size; group_idx++)
{
unsigned int j = m_group->getMemberIndex(group_idx);
quat<Scalar> q(h_orientation.data[j]);
quat<Scalar> p(h_angmom.data[j]);
vec3<Scalar> t(h_net_torque.data[j]);
vec3<Scalar> I(h_inertia.data[j]);
// rotate torque into principal frame
t = rotate(conj(q), t);
// check for zero moment of inertia
bool x_zero, y_zero, z_zero;
x_zero = (I.x == 0);
y_zero = (I.y == 0);
z_zero = (I.z == 0);
// ignore torque component along an axis for which the moment of inertia zero
if (x_zero)
t.x = 0;
if (y_zero)
t.y = 0;
if (z_zero)
t.z = 0;
// advance p(t)->p(t+deltaT/2), q(t)->q(t+deltaT)
// using Trotter factorization of rotation Liouvillian
p += m_deltaT * q * t;
quat<Scalar> p1, p2, p3; // permutated quaternions
quat<Scalar> q1, q2, q3;
Scalar phi1, cphi1, sphi1;
Scalar phi2, cphi2, sphi2;
Scalar phi3, cphi3, sphi3;
if (!z_zero)
{
p3 = quat<Scalar>(-p.v.z, vec3<Scalar>(p.v.y, -p.v.x, p.s));
q3 = quat<Scalar>(-q.v.z, vec3<Scalar>(q.v.y, -q.v.x, q.s));
phi3 = Scalar(1. / 4.) / I.z * dot(p, q3);
cphi3 = slow::cos(Scalar(1. / 2.) * m_deltaT * phi3);
sphi3 = slow::sin(Scalar(1. / 2.) * m_deltaT * phi3);
p = cphi3 * p + sphi3 * p3;
q = cphi3 * q + sphi3 * q3;
}
if (!y_zero)
{
p2 = quat<Scalar>(-p.v.y, vec3<Scalar>(-p.v.z, p.s, p.v.x));
q2 = quat<Scalar>(-q.v.y, vec3<Scalar>(-q.v.z, q.s, q.v.x));
phi2 = Scalar(1. / 4.) / I.y * dot(p, q2);
cphi2 = slow::cos(Scalar(1. / 2.) * m_deltaT * phi2);
sphi2 = slow::sin(Scalar(1. / 2.) * m_deltaT * phi2);
p = cphi2 * p + sphi2 * p2;
q = cphi2 * q + sphi2 * q2;
}
if (!x_zero)
{
p1 = quat<Scalar>(-p.v.x, vec3<Scalar>(p.s, p.v.z, -p.v.y));
q1 = quat<Scalar>(-q.v.x, vec3<Scalar>(q.s, q.v.z, -q.v.y));
phi1 = Scalar(1. / 4.) / I.x * dot(p, q1);
cphi1 = slow::cos(m_deltaT * phi1);
sphi1 = slow::sin(m_deltaT * phi1);
p = cphi1 * p + sphi1 * p1;
q = cphi1 * q + sphi1 * q1;
}
if (!y_zero)
{
p2 = quat<Scalar>(-p.v.y, vec3<Scalar>(-p.v.z, p.s, p.v.x));
q2 = quat<Scalar>(-q.v.y, vec3<Scalar>(-q.v.z, q.s, q.v.x));
phi2 = Scalar(1. / 4.) / I.y * dot(p, q2);
cphi2 = slow::cos(Scalar(1. / 2.) * m_deltaT * phi2);
sphi2 = slow::sin(Scalar(1. / 2.) * m_deltaT * phi2);
p = cphi2 * p + sphi2 * p2;
q = cphi2 * q + sphi2 * q2;
}
if (!z_zero)
{
p3 = quat<Scalar>(-p.v.z, vec3<Scalar>(p.v.y, -p.v.x, p.s));
q3 = quat<Scalar>(-q.v.z, vec3<Scalar>(q.v.y, -q.v.x, q.s));
phi3 = Scalar(1. / 4.) / I.z * dot(p, q3);
cphi3 = slow::cos(Scalar(1. / 2.) * m_deltaT * phi3);
sphi3 = slow::sin(Scalar(1. / 2.) * m_deltaT * phi3);
p = cphi3 * p + sphi3 * p3;
q = cphi3 * q + sphi3 * q3;
}
// renormalize (improves stability)
q = q * (Scalar(1.0) / slow::sqrt(norm2(q)));
h_orientation.data[j] = quat_to_scalar4(q);
h_angmom.data[j] = quat_to_scalar4(p);
}
}
}
/*! \param timestep Current time step
\post particle velocities are moved forward to timestep+1
*/
void TwoStepLangevin::integrateStepTwo(uint64_t timestep)
{
unsigned int group_size = m_group->getNumMembers();
const GlobalArray<Scalar4>& net_force = m_pdata->getNetForce();
ArrayHandle<Scalar4> h_vel(m_pdata->getVelocities(),
access_location::host,
access_mode::readwrite);
ArrayHandle<Scalar3> h_accel(m_pdata->getAccelerations(),
access_location::host,
access_mode::readwrite);
ArrayHandle<Scalar4> h_pos(m_pdata->getPositions(), access_location::host, access_mode::read);
ArrayHandle<unsigned int> h_tag(m_pdata->getTags(), access_location::host, access_mode::read);
ArrayHandle<Scalar4> h_net_force(net_force, access_location::host, access_mode::read);
ArrayHandle<Scalar> h_gamma(m_gamma, access_location::host, access_mode::read);
ArrayHandle<Scalar3> h_gamma_r(m_gamma_r, access_location::host, access_mode::read);
ArrayHandle<Scalar4> h_orientation(m_pdata->getOrientationArray(),
access_location::host,
access_mode::readwrite);
ArrayHandle<Scalar4> h_angmom(m_pdata->getAngularMomentumArray(),
access_location::host,
access_mode::readwrite);
ArrayHandle<Scalar4> h_net_torque(m_pdata->getNetTorqueArray(),
access_location::host,
access_mode::readwrite);
ArrayHandle<Scalar3> h_inertia(m_pdata->getMomentsOfInertiaArray(),
access_location::host,
access_mode::read);
// grab some initial variables
const Scalar currentTemp = m_T->operator()(timestep);
const unsigned int D = m_sysdef->getNDimensions();
// energy transferred over this time step
Scalar bd_energy_transfer = 0;
// a(t+deltaT) gets modified with the bd forces
// v(t+deltaT) = v(t+deltaT/2) + 1/2 * a(t+deltaT)*deltaT
uint16_t seed = m_sysdef->getSeed();
for (unsigned int group_idx = 0; group_idx < group_size; group_idx++)
{
unsigned int j = m_group->getMemberIndex(group_idx);
unsigned int ptag = h_tag.data[j];
// Initialize the RNG
RandomGenerator rng(hoomd::Seed(RNGIdentifier::TwoStepLangevin, timestep, seed),
hoomd::Counter(ptag));
// first, calculate the BD forces
// Generate three random numbers
hoomd::UniformDistribution<Scalar> uniform(Scalar(-1), Scalar(1));
Scalar rx = uniform(rng);
Scalar ry = uniform(rng);
Scalar rz = uniform(rng);
Scalar gamma;
unsigned int type = __scalar_as_int(h_pos.data[j].w);
gamma = h_gamma.data[type];
// compute the bd force
Scalar coeff = fast::sqrt(Scalar(6.0) * gamma * currentTemp / m_deltaT);
if (m_noiseless_t)
coeff = Scalar(0.0);
Scalar bd_fx = rx * coeff - gamma * h_vel.data[j].x;
Scalar bd_fy = ry * coeff - gamma * h_vel.data[j].y;
Scalar bd_fz = rz * coeff - gamma * h_vel.data[j].z;
if (D < 3)
bd_fz = Scalar(0.0);
// then, calculate acceleration from the net force
Scalar minv = Scalar(1.0) / h_vel.data[j].w;
h_accel.data[j].x = (h_net_force.data[j].x + bd_fx) * minv;
h_accel.data[j].y = (h_net_force.data[j].y + bd_fy) * minv;
h_accel.data[j].z = (h_net_force.data[j].z + bd_fz) * minv;
// then, update the velocity
h_vel.data[j].x += Scalar(1.0 / 2.0) * h_accel.data[j].x * m_deltaT;
h_vel.data[j].y += Scalar(1.0 / 2.0) * h_accel.data[j].y * m_deltaT;
h_vel.data[j].z += Scalar(1.0 / 2.0) * h_accel.data[j].z * m_deltaT;
// tally the energy transfer from the bd thermal reservoir to the particles
if (m_tally)
bd_energy_transfer
+= bd_fx * h_vel.data[j].x + bd_fy * h_vel.data[j].y + bd_fz * h_vel.data[j].z;
// rotational updates
if (m_aniso)
{
unsigned int type_r = __scalar_as_int(h_pos.data[j].w);
Scalar3 gamma_r = h_gamma_r.data[type_r];
// get body frame ang_mom
quat<Scalar> p(h_angmom.data[j]);
quat<Scalar> q(h_orientation.data[j]);
vec3<Scalar> t(h_net_torque.data[j]);
vec3<Scalar> I(h_inertia.data[j]);
// s is the pure imaginary quaternion with im. part equal to true angular velocity
vec3<Scalar> s;
s = (Scalar(1. / 2.) * conj(q) * p).v;
if (gamma_r.x > 0 || gamma_r.y > 0 || gamma_r.z > 0)
{
// first calculate in the body frame random and damping torque imposed by the
// dynamics
vec3<Scalar> bf_torque;
// original Gaussian random torque
Scalar3 sigma_r
= make_scalar3(fast::sqrt(Scalar(2.0) * gamma_r.x * currentTemp / m_deltaT),
fast::sqrt(Scalar(2.0) * gamma_r.y * currentTemp / m_deltaT),
fast::sqrt(Scalar(2.0) * gamma_r.z * currentTemp / m_deltaT));
if (m_noiseless_r)
sigma_r = make_scalar3(0.0, 0.0, 0.0);
Scalar rand_x = hoomd::NormalDistribution<Scalar>(sigma_r.x)(rng);
Scalar rand_y = hoomd::NormalDistribution<Scalar>(sigma_r.y)(rng);
Scalar rand_z = hoomd::NormalDistribution<Scalar>(sigma_r.z)(rng);
// check for degenerate moment of inertia
bool x_zero, y_zero, z_zero;
x_zero = (I.x == 0);
y_zero = (I.y == 0);
z_zero = (I.z == 0);
bf_torque.x = rand_x - gamma_r.x * (s.x / I.x);
bf_torque.y = rand_y - gamma_r.y * (s.y / I.y);
bf_torque.z = rand_z - gamma_r.z * (s.z / I.z);
// ignore torque component along an axis for which the moment of inertia zero
if (x_zero)
bf_torque.x = 0;
if (y_zero)
bf_torque.y = 0;
if (z_zero)
bf_torque.z = 0;
// change to lab frame and update the net torque
bf_torque = rotate(q, bf_torque);
h_net_torque.data[j].x += bf_torque.x;
h_net_torque.data[j].y += bf_torque.y;
h_net_torque.data[j].z += bf_torque.z;
if (D < 3)
h_net_torque.data[j].x = 0;
if (D < 3)
h_net_torque.data[j].y = 0;
}
}
}
// then, update the angular velocity
if (m_aniso)
{
// angular degrees of freedom
for (unsigned int group_idx = 0; group_idx < group_size; group_idx++)
{
unsigned int j = m_group->getMemberIndex(group_idx);
quat<Scalar> q(h_orientation.data[j]);
quat<Scalar> p(h_angmom.data[j]);
vec3<Scalar> t(h_net_torque.data[j]);
vec3<Scalar> I(h_inertia.data[j]);
// rotate torque into principal frame
t = rotate(conj(q), t);
// check for zero moment of inertia
bool x_zero, y_zero, z_zero;
x_zero = (I.x == 0);
y_zero = (I.y == 0);
z_zero = (I.z == 0);
// ignore torque component along an axis for which the moment of inertia zero
if (x_zero)
t.x = 0;
if (y_zero)
t.y = 0;
if (z_zero)
t.z = 0;
// advance p(t+deltaT/2)->p(t+deltaT)
p += m_deltaT * q * t;
h_angmom.data[j] = quat_to_scalar4(p);
}
}
// update energy reservoir
if (m_tally)
{
#ifdef ENABLE_MPI
if (m_sysdef->isDomainDecomposed())
{
MPI_Allreduce(MPI_IN_PLACE,
&bd_energy_transfer,
1,
MPI_HOOMD_SCALAR,
MPI_SUM,
m_exec_conf->getMPICommunicator());
}
#endif
m_reservoir_energy -= bd_energy_transfer * m_deltaT;
m_extra_energy_overdeltaT = 0.5 * bd_energy_transfer;
}
}
namespace detail
{
void export_TwoStepLangevin(pybind11::module& m)
{
pybind11::class_<TwoStepLangevin, TwoStepLangevinBase, std::shared_ptr<TwoStepLangevin>>(
m,
"TwoStepLangevin")
.def(pybind11::init<std::shared_ptr<SystemDefinition>,
std::shared_ptr<ParticleGroup>,
std::shared_ptr<Variant>>())
.def_property("tally_reservoir_energy",
&TwoStepLangevin::getTallyReservoirEnergy,
&TwoStepLangevin::setTallyReservoirEnergy)
.def_property_readonly("reservoir_energy", &TwoStepLangevin::getReservoirEnergy);
}
} // end namespace detail
} // end namespace md
} // end namespace hoomd