OpenMD 3.1
Molecular Dynamics in the Open
Loading...
Searching...
No Matches
LDForceModifier.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#include "integrators/LDForceModifier.hpp"
45
46#include <cmath>
47#include <fstream>
48#include <iostream>
49#include <memory>
50#include <random>
51
52#include "brains/ForceModifier.hpp"
53#include "hydrodynamics/Ellipsoid.hpp"
54#include "hydrodynamics/HydroIO.hpp"
55#include "hydrodynamics/Sphere.hpp"
56#include "math/CholeskyDecomposition.hpp"
58#include "types/GayBerneAdapter.hpp"
59#include "types/LennardJonesAdapter.hpp"
60#include "utils/Constants.hpp"
62
63using namespace std;
64namespace OpenMD {
65
66 LDForceModifier::LDForceModifier(SimInfo* info) :
67 ForceModifier {info}, maxIterNum_ {6}, forceTolerance_ {1e-6},
68 randNumGen_ {info->getRandomNumberGenerator()},
69 simParams_ {info->getSimParams()} {
70 RealType dt = simParams_->getDt();
71 dt2_ = dt * 0.5;
72
73 veloMunge_ = std::make_unique<Velocitizer>(info_);
74
75 sphericalBoundaryConditions_ = false;
76 if (simParams_->getUseSphericalBoundaryConditions()) {
77 sphericalBoundaryConditions_ = true;
78 if (simParams_->haveLangevinBufferRadius()) {
79 langevinBufferRadius_ = simParams_->getLangevinBufferRadius();
80 } else {
81 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
82 "langevinBufferRadius must be specified "
83 "when useSphericalBoundaryConditions is turned on.\n");
84 painCave.severity = OPENMD_ERROR;
85 painCave.isFatal = 1;
86 simError();
87 }
88
89 if (simParams_->haveFrozenBufferRadius()) {
90 frozenBufferRadius_ = simParams_->getFrozenBufferRadius();
91 } else {
92 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
93 "frozenBufferRadius must be specified "
94 "when useSphericalBoundaryConditions is turned on.\n");
95 painCave.severity = OPENMD_ERROR;
96 painCave.isFatal = 1;
97 simError();
98 }
99
100 if (frozenBufferRadius_ < langevinBufferRadius_) {
101 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
102 "frozenBufferRadius has been set smaller than the "
103 "langevinBufferRadius. This is probably an error.\n");
104 painCave.severity = OPENMD_WARNING;
105 painCave.isFatal = 0;
106 simError();
107 }
108 }
109
110 // Build the hydroProp_ map:
111 Molecule* mol;
112 StuntDouble* sd;
113 SimInfo::MoleculeIterator i;
114 Molecule::IntegrableObjectIterator j;
115 bool needHydroPropFile = false;
116
117 for (mol = info_->beginMolecule(i); mol != NULL;
118 mol = info_->nextMolecule(i)) {
119 for (sd = mol->beginIntegrableObject(j); sd != NULL;
120 sd = mol->nextIntegrableObject(j)) {
121 if (sd->isRigidBody()) {
122 RigidBody* rb = static_cast<RigidBody*>(sd);
123 if (rb->getNumAtoms() > 1) needHydroPropFile = true;
124 }
125 }
126 }
127
128 if (needHydroPropFile) {
129 if (simParams_->haveHydroPropFile()) {
130 HydroIO* hio = new HydroIO();
131 hydroPropMap_ = hio->parseHydroFile(simParams_->getHydroPropFile());
132 } else {
133 snprintf(
134 painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
135 "HydroPropFile must be set to a file name if Langevin Dynamics\n"
136 "\tis specified for rigidBodies which contain more than one atom\n"
137 "\tTo create a HydroPropFile, run the \"Hydro\" program.\n");
138 painCave.severity = OPENMD_ERROR;
139 painCave.isFatal = 1;
140 simError();
141 }
142
143 for (mol = info_->beginMolecule(i); mol != NULL;
144 mol = info_->nextMolecule(i)) {
145 for (sd = mol->beginIntegrableObject(j); sd != NULL;
146 sd = mol->nextIntegrableObject(j)) {
147 map<string, HydroProp*>::iterator iter =
148 hydroPropMap_.find(sd->getType());
149 if (iter != hydroPropMap_.end()) {
150 hydroProps_.push_back(iter->second);
151 moments_.push_back(getMomentData(sd));
152
153 } else {
154 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
155 "Can not find resistance tensor for atom [%s]\n",
156 sd->getType().c_str());
157 painCave.severity = OPENMD_ERROR;
158 painCave.isFatal = 1;
159 simError();
160 }
161 }
162 }
163
164 } else {
165 for (mol = info_->beginMolecule(i); mol != NULL;
166 mol = info_->nextMolecule(i)) {
167 for (sd = mol->beginIntegrableObject(j); sd != NULL;
168 sd = mol->nextIntegrableObject(j)) {
169 Shape* currShape = NULL;
170
171 if (sd->isAtom()) {
172 Atom* atom = static_cast<Atom*>(sd);
173 AtomType* atomType = atom->getAtomType();
174 GayBerneAdapter gba = GayBerneAdapter(atomType);
175 if (gba.isGayBerne()) {
176 currShape = new Ellipsoid(V3Zero, gba.getL() / 2.0,
177 gba.getD() / 2.0, Mat3x3d::identity());
178 } else {
179 LennardJonesAdapter lja = LennardJonesAdapter(atomType);
180 if (lja.isLennardJones()) {
181 currShape = new Sphere(V3Zero, lja.getSigma() / 2.0);
182 } else {
183 int aNum(0);
184 vector<AtomType*> atChain = atomType->allYourBase();
185 vector<AtomType*>::iterator i;
186 for (i = atChain.begin(); i != atChain.end(); ++i) {
187 aNum = etab.GetAtomicNum((*i)->getName().c_str());
188 if (aNum != 0) {
189 currShape = new Sphere(V3Zero, etab.GetVdwRad(aNum));
190 break;
191 }
192 }
193 if (aNum == 0) {
194 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
195 "Could not find atom type in default element.txt\n");
196 painCave.severity = OPENMD_ERROR;
197 painCave.isFatal = 1;
198 simError();
199 }
200 }
201 }
202 }
203
204 if (!simParams_->haveTargetTemp()) {
205 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
206 "You can't use LangevinDynamics without a targetTemp!\n");
207 painCave.isFatal = 1;
208 painCave.severity = OPENMD_ERROR;
209 simError();
210 }
211
212 if (!simParams_->haveViscosity()) {
213 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
214 "You can't use LangevinDynamics without a viscosity!\n");
215 painCave.isFatal = 1;
216 painCave.severity = OPENMD_ERROR;
217 simError();
218 }
219
220 HydroProp* currHydroProp =
221 currShape->getHydroProp(simParams_->getViscosity());
222 map<string, HydroProp*>::iterator iter =
223 hydroPropMap_.find(sd->getType());
224 if (iter != hydroPropMap_.end()) {
225 hydroProps_.push_back(iter->second);
226 moments_.push_back(getMomentData(sd));
227
228 } else {
229 currHydroProp->complete();
230 hydroPropMap_.insert(map<string, HydroProp*>::value_type(
231 sd->getType(), currHydroProp));
232 hydroProps_.push_back(currHydroProp);
233 moments_.push_back(getMomentData(sd));
234 }
235 delete currShape;
236 }
237 }
238 }
239
240 RealType stdDev =
241 std::sqrt(2.0 * Constants::kb * simParams_->getTargetTemp() / dt);
242
243 forceDistribution_ = std::normal_distribution<RealType>(0.0, stdDev);
244 }
245
246 MomentData* LDForceModifier::getMomentData(StuntDouble* sd) {
247 map<string, MomentData*>::iterator j = momentsMap_.find(sd->getType());
248 if (j != momentsMap_.end()) {
249 return j->second;
250 } else {
251 MomentData* moment = new MomentData;
252 map<string, HydroProp*>::iterator i = hydroPropMap_.find(sd->getType());
253
254 if (i != hydroPropMap_.end()) {
255 // Center of Resistance:
256 moment->rcr = i->second->getCenterOfResistance();
257 } else {
258 snprintf(
259 painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
260 "LDForceManager createMomentData: Couldn't find HydroProp for\n"
261 "object type %s!\n",
262 sd->getType().c_str());
263 painCave.isFatal = 1;
264 painCave.severity = OPENMD_ERROR;
265 simError();
266 }
267
268 if (sd->isRigidBody())
269 moment->rcr -= dynamic_cast<RigidBody*>(sd)->getRefCOM();
270
271 // Parallel axis formula to get the moment of inertia around
272 // the center of resistance:
273 RealType mass = sd->getMass();
274 moment->Icr = sd->getI();
275 moment->Icr +=
276 mass * (dot(moment->rcr, moment->rcr) * Mat3x3d::identity() +
277 outProduct(moment->rcr, moment->rcr));
278 moment->IcrInv = moment->Icr.inverse();
279
280 momentsMap_.insert(
281 map<string, MomentData*>::value_type(sd->getType(), moment));
282 return moment;
283 }
284 }
285
286 void LDForceModifier::modifyForces() {
287 SimInfo::MoleculeIterator i;
288 Molecule::IntegrableObjectIterator j;
289 Molecule* mol;
290 StuntDouble* sd;
291 RealType mass;
292 Vector3d pos;
293 Vector3d frc;
294 Mat3x3d A;
295 Mat3x3d Atrans;
296 Vector3d Tb;
297 Vector3d ji;
298 unsigned int index = 0;
299 bool doLangevinForces;
300 bool freezeMolecule;
301 int fdf;
302
303 fdf = 0;
304
305 for (mol = info_->beginMolecule(i); mol != NULL;
306 mol = info_->nextMolecule(i)) {
307 doLangevinForces = true;
308 freezeMolecule = false;
309
310 if (sphericalBoundaryConditions_) {
311 Vector3d molPos = mol->getCom();
312 RealType molRad = molPos.length();
313
314 doLangevinForces = false;
315
316 if (molRad > langevinBufferRadius_) {
317 doLangevinForces = true;
318 freezeMolecule = false;
319 }
320 if (molRad > frozenBufferRadius_) {
321 doLangevinForces = false;
322 freezeMolecule = true;
323 }
324 }
325
326 for (sd = mol->beginIntegrableObject(j); sd != NULL;
327 sd = mol->nextIntegrableObject(j)) {
328 if (freezeMolecule) fdf += sd->freeze();
329
330 if (doLangevinForces) {
331 mass = sd->getMass();
332 if (sd->isDirectional()) {
333 // preliminaries for directional objects:
334
335 A = sd->getA();
336 Atrans = A.transpose();
337
338 Vector3d rcrLab = Atrans * moments_[index]->rcr;
339
340 // apply random force and torque at center of resistance
341
342 Vector3d randomForceBody;
343 Vector3d randomTorqueBody;
344 genRandomForceAndTorque(randomForceBody, randomTorqueBody, index);
345 Vector3d randomForceLab = Atrans * randomForceBody;
346 Vector3d randomTorqueLab = Atrans * randomTorqueBody;
347
348 sd->addFrc(randomForceLab);
349 sd->addTrq(randomTorqueLab + cross(rcrLab, randomForceLab));
350
351 Vector3d omegaBody;
352
353 // What remains contains velocity explicitly, but the
354 // velocity required is at the full step: v(t + h), while
355 // we have initially the velocity at the half step: v(t + h/2).
356 // We need to iterate to converge the friction
357 // force and friction torque vectors.
358
359 // this is the velocity at the half-step:
360
361 Vector3d vel = sd->getVel();
362 Vector3d angMom = sd->getJ();
363
364 // estimate velocity at full-step using everything but
365 // friction forces:
366
367 frc = sd->getFrc();
368
369 Vector3d velStep =
370 vel + (dt2_ / mass * Constants::energyConvert) * frc;
371
372 Tb = sd->lab2Body(sd->getTrq());
373 Vector3d angMomStep =
374 angMom + (dt2_ * Constants::energyConvert) * Tb;
375
376 Vector3d omegaLab;
377 Vector3d vcdLab;
378 Vector3d vcdBody;
379 Vector3d frictionForceBody;
380 Vector3d frictionForceLab(0.0);
381 Vector3d oldFFL; // used to test for convergence
382 Vector3d frictionTorqueBody(0.0);
383 Vector3d oldFTB; // used to test for convergence
384 Vector3d frictionTorqueLab;
385 RealType fdot;
386 RealType tdot;
387
388 // iteration starts here:
389
390 for (int k = 0; k < maxIterNum_; k++) {
391 if (sd->isLinear()) {
392 int linearAxis = sd->linearAxis();
393 int l = (linearAxis + 1) % 3;
394 int m = (linearAxis + 2) % 3;
395 omegaBody[l] = angMomStep[l] / moments_[index]->Icr(l, l);
396 omegaBody[m] = angMomStep[m] / moments_[index]->Icr(m, m);
397
398 } else {
399 omegaBody = moments_[index]->IcrInv * angMomStep;
400 // omegaBody[0] = angMomStep[0] /I(0, 0);
401 // omegaBody[1] = angMomStep[1] /I(1, 1);
402 // omegaBody[2] = angMomStep[2] /I(2, 2);
403 }
404
405 omegaLab = Atrans * omegaBody;
406
407 // apply friction force and torque at center of resistance
408
409 vcdLab = velStep + cross(omegaLab, rcrLab);
410 vcdBody = A * vcdLab;
411 frictionForceBody = -(hydroProps_[index]->getXitt() * vcdBody +
412 hydroProps_[index]->getXirt() * omegaBody);
413 oldFFL = frictionForceLab;
414 frictionForceLab = Atrans * frictionForceBody;
415 oldFTB = frictionTorqueBody;
416 frictionTorqueBody = -(hydroProps_[index]->getXitr() * vcdBody +
417 hydroProps_[index]->getXirr() * omegaBody);
418 frictionTorqueLab = Atrans * frictionTorqueBody;
419
420 // re-estimate velocities at full-step using friction forces:
421
422 velStep = vel + (dt2_ / mass * Constants::energyConvert) *
423 (frc + frictionForceLab);
424 angMomStep = angMom + (dt2_ * Constants::energyConvert) *
425 (Tb + frictionTorqueBody);
426
427 // check for convergence (if the vectors have converged, fdot and
428 // tdot will both be 1.0):
429
430 fdot = dot(frictionForceLab, oldFFL) /
431 frictionForceLab.lengthSquare();
432 tdot = dot(frictionTorqueBody, oldFTB) /
433 frictionTorqueBody.lengthSquare();
434
435 if (fabs(1.0 - fdot) <= forceTolerance_ &&
436 fabs(1.0 - tdot) <= forceTolerance_)
437 break; // iteration ends here
438 }
439
440 sd->addFrc(frictionForceLab);
441 sd->addTrq(frictionTorqueLab + cross(rcrLab, frictionForceLab));
442
443 } else {
444 // spherical atom
445
446 Vector3d systemForce;
447 Vector3d randomForce;
448 Vector3d randomTorque;
449 genRandomForceAndTorque(randomForce, randomTorque, index);
450 systemForce = sd->getFrc();
451 sd->addFrc(randomForce);
452
453 // What remains contains velocity explicitly, but the
454 // velocity required is at the full step: v(t + h), while
455 // we have initially the velocity at the half step: v(t + h/2).
456 // We need to iterate to converge the friction
457 // force vector.
458
459 // this is the velocity at the half-step:
460
461 Vector3d vel = sd->getVel();
462
463 // estimate velocity at full-step using everything but
464 // friction forces:
465
466 frc = sd->getFrc();
467 Vector3d velStep =
468 vel + (dt2_ / mass * Constants::energyConvert) * frc;
469
470 Vector3d frictionForce(0.0);
471 Vector3d oldFF; // used to test for convergence
472 RealType fdot;
473
474 // iteration starts here:
475
476 for (int k = 0; k < maxIterNum_; k++) {
477 oldFF = frictionForce;
478 frictionForce = -hydroProps_[index]->getXitt() * velStep;
479
480 // re-estimate velocities at full-step using friction forces:
481
482 velStep = vel + (dt2_ / mass * Constants::energyConvert) *
483 (frc + frictionForce);
484
485 // check for convergence (if the vector has converged,
486 // fdot will be 1.0):
487
488 fdot = dot(frictionForce, oldFF) / frictionForce.lengthSquare();
489
490 if (fabs(1.0 - fdot) <= forceTolerance_)
491 break; // iteration ends here
492 }
493
494 sd->addFrc(frictionForce);
495 }
496 }
497
498 ++index;
499 }
500 }
501
502 info_->setFdf(fdf);
503 if(simParams_->getConserveLinearMomentum()) veloMunge_->removeComDrift();
504 // Remove angular drift if we are not using periodic boundary conditions.
505 if (!simParams_->getUsePeriodicBoundaryConditions() &&
506 simParams_->getConserveAngularMomentum())
507 veloMunge_->removeAngularDrift();
508 }
509
510 void LDForceModifier::genRandomForceAndTorque(Vector3d& force,
511 Vector3d& torque,
512 unsigned int index) {
513 Vector<RealType, 6> Z;
514 Vector<RealType, 6> generalForce;
515
516 Z[0] = forceDistribution_(*randNumGen_);
517 Z[1] = forceDistribution_(*randNumGen_);
518 Z[2] = forceDistribution_(*randNumGen_);
519 Z[3] = forceDistribution_(*randNumGen_);
520 Z[4] = forceDistribution_(*randNumGen_);
521 Z[5] = forceDistribution_(*randNumGen_);
522
523 generalForce = hydroProps_[index]->getS() * Z;
524
525 force[0] = generalForce[0];
526 force[1] = generalForce[1];
527 force[2] = generalForce[2];
528 torque[0] = generalForce[3];
529 torque[1] = generalForce[4];
530 torque[2] = generalForce[5];
531 }
532} // namespace OpenMD
This basic Periodic Table class was originally taken from the data.h file in OpenBabel.
RealType GetVdwRad(int atomicnum)
int GetAtomicNum(const char *str)
Real length()
Returns the length of this vector.
Definition Vector.hpp:393
This basic Periodic Table class was originally taken from the data.cpp file in OpenBabel.
Vector3< Real > cross(const Vector3< Real > &v1, const Vector3< Real > &v2)
Returns the cross product of two Vectors.
Definition Vector3.hpp:136
Real dot(const DynamicVector< Real > &v1, const DynamicVector< Real > &v2)
Returns the dot product of two DynamicVectors.