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