ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/integrators/LDForceManager.cpp
Revision: 1850
Committed: Wed Feb 20 15:39:39 2013 UTC (12 years, 3 months ago) by gezelter
File size: 17532 byte(s)
Log Message:
Fixed a widespread typo in the license 

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

Properties

Name Value
svn:executable *
svn:keywords Author Id Revision Date