ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/integrators/LDForceManager.cpp
Revision: 3250
Committed: Fri Oct 5 19:02:09 2007 UTC (16 years, 9 months ago) by xsun
File size: 17242 byte(s)
Log Message:
Fixed some bugs in LDForceManager.cpp

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. Acknowledgement of the program authors must be made in any
10 * publication of scientific results based in part on use of the
11 * program. An acceptable form of acknowledgement is citation of
12 * the article in which the program was described (Matthew
13 * A. Meineke, Charles F. Vardeman II, Teng Lin, Christopher
14 * J. Fennell and J. Daniel Gezelter, "OOPSE: An Object-Oriented
15 * Parallel Simulation Engine for Molecular Dynamics,"
16 * J. Comput. Chem. 26, pp. 252-271 (2005))
17 *
18 * 2. Redistributions of source code must retain the above copyright
19 * notice, this list of conditions and the following disclaimer.
20 *
21 * 3. Redistributions in binary form must reproduce the above copyright
22 * notice, this list of conditions and the following disclaimer in the
23 * documentation and/or other materials provided with the
24 * distribution.
25 *
26 * This software is provided "AS IS," without a warranty of any
27 * kind. All express or implied conditions, representations and
28 * warranties, including any implied warranty of merchantability,
29 * fitness for a particular purpose or non-infringement, are hereby
30 * excluded. The University of Notre Dame and its licensors shall not
31 * be liable for any damages suffered by licensee as a result of
32 * using, modifying or distributing the software or its
33 * derivatives. In no event will the University of Notre Dame or its
34 * licensors be liable for any lost revenue, profit or data, or for
35 * direct, indirect, special, consequential, incidental or punitive
36 * damages, however caused and regardless of the theory of liability,
37 * arising out of the use of or inability to use software, even if the
38 * University of Notre Dame has been advised of the possibility of
39 * such damages.
40 */
41 #include <fstream>
42 #include <iostream>
43 #include "integrators/LDForceManager.hpp"
44 #include "math/CholeskyDecomposition.hpp"
45 #include "utils/OOPSEConstant.hpp"
46 #include "hydrodynamics/Sphere.hpp"
47 #include "hydrodynamics/Ellipsoid.hpp"
48 #include "openbabel/mol.hpp"
49
50 using namespace OpenBabel;
51 namespace oopse {
52
53 LDForceManager::LDForceManager(SimInfo* info) : ForceManager(info){
54 simParams = info->getSimParams();
55 veloMunge = new Velocitizer(info);
56
57 sphericalBoundaryConditions_ = false;
58 if (simParams->getUseSphericalBoundaryConditions()) {
59 sphericalBoundaryConditions_ = true;
60 if (simParams->haveLangevinBufferRadius()) {
61 langevinBufferRadius_ = simParams->getLangevinBufferRadius();
62 } else {
63 sprintf( painCave.errMsg,
64 "langevinBufferRadius must be specified "
65 "when useSphericalBoundaryConditions is turned on.\n");
66 painCave.severity = OOPSE_ERROR;
67 painCave.isFatal = 1;
68 simError();
69 }
70
71 if (simParams->haveFrozenBufferRadius()) {
72 frozenBufferRadius_ = simParams->getFrozenBufferRadius();
73 } else {
74 sprintf( painCave.errMsg,
75 "frozenBufferRadius must be specified "
76 "when useSphericalBoundaryConditions is turned on.\n");
77 painCave.severity = OOPSE_ERROR;
78 painCave.isFatal = 1;
79 simError();
80 }
81
82 if (frozenBufferRadius_ < langevinBufferRadius_) {
83 sprintf( painCave.errMsg,
84 "frozenBufferRadius has been set smaller than the "
85 "langevinBufferRadius. This is probably an error.\n");
86 painCave.severity = OOPSE_WARNING;
87 painCave.isFatal = 0;
88 simError();
89 }
90 }
91
92 // Build the hydroProp map:
93 std::map<std::string, HydroProp*> hydroPropMap;
94
95 Molecule* mol;
96 StuntDouble* integrableObject;
97 SimInfo::MoleculeIterator i;
98 Molecule::IntegrableObjectIterator j;
99 bool needHydroPropFile = false;
100
101 for (mol = info->beginMolecule(i); mol != NULL;
102 mol = info->nextMolecule(i)) {
103 for (integrableObject = mol->beginIntegrableObject(j);
104 integrableObject != NULL;
105 integrableObject = mol->nextIntegrableObject(j)) {
106
107 if (integrableObject->isRigidBody()) {
108 RigidBody* rb = static_cast<RigidBody*>(integrableObject);
109 if (rb->getNumAtoms() > 1) needHydroPropFile = true;
110 }
111
112 }
113 }
114
115
116 if (needHydroPropFile) {
117 if (simParams->haveHydroPropFile()) {
118 hydroPropMap = parseFrictionFile(simParams->getHydroPropFile());
119 } else {
120 sprintf( painCave.errMsg,
121 "HydroPropFile must be set to a file name if Langevin\n"
122 "\tDynamics is specified for rigidBodies which contain more\n"
123 "\tthan one atom. To create a HydroPropFile, run \"Hydro\".\n");
124 painCave.severity = OOPSE_ERROR;
125 painCave.isFatal = 1;
126 simError();
127 }
128
129 for (mol = info->beginMolecule(i); mol != NULL;
130 mol = info->nextMolecule(i)) {
131 for (integrableObject = mol->beginIntegrableObject(j);
132 integrableObject != NULL;
133 integrableObject = mol->nextIntegrableObject(j)) {
134
135 std::map<std::string, HydroProp*>::iterator iter = hydroPropMap.find(integrableObject->getType());
136 if (iter != hydroPropMap.end()) {
137 hydroProps_.push_back(iter->second);
138 } else {
139 sprintf( painCave.errMsg,
140 "Can not find resistance tensor for atom [%s]\n", integrableObject->getType().c_str());
141 painCave.severity = OOPSE_ERROR;
142 painCave.isFatal = 1;
143 simError();
144 }
145 }
146 }
147 } else {
148
149 std::map<std::string, HydroProp*> hydroPropMap;
150 for (mol = info->beginMolecule(i); mol != NULL;
151 mol = info->nextMolecule(i)) {
152 for (integrableObject = mol->beginIntegrableObject(j);
153 integrableObject != NULL;
154 integrableObject = mol->nextIntegrableObject(j)) {
155 Shape* currShape = NULL;
156
157 if (integrableObject->isAtom()){
158 Atom* atom = static_cast<Atom*>(integrableObject);
159 AtomType* atomType = atom->getAtomType();
160 if (atomType->isGayBerne()) {
161 DirectionalAtomType* dAtomType = dynamic_cast<DirectionalAtomType*>(atomType);
162 GenericData* data = dAtomType->getPropertyByName("GayBerne");
163 if (data != NULL) {
164 GayBerneParamGenericData* gayBerneData = dynamic_cast<GayBerneParamGenericData*>(data);
165
166 if (gayBerneData != NULL) {
167 GayBerneParam gayBerneParam = gayBerneData->getData();
168 currShape = new Ellipsoid(V3Zero,
169 gayBerneParam.GB_l / 2.0,
170 gayBerneParam.GB_d / 2.0,
171 Mat3x3d::identity());
172 } else {
173 sprintf( painCave.errMsg,
174 "Can not cast GenericData to GayBerneParam\n");
175 painCave.severity = OOPSE_ERROR;
176 painCave.isFatal = 1;
177 simError();
178 }
179 } else {
180 sprintf( painCave.errMsg, "Can not find Parameters for GayBerne\n");
181 painCave.severity = OOPSE_ERROR;
182 painCave.isFatal = 1;
183 simError();
184 }
185 } else {
186 if (atomType->isLennardJones()){
187 GenericData* data = atomType->getPropertyByName("LennardJones");
188 if (data != NULL) {
189 LJParamGenericData* ljData = dynamic_cast<LJParamGenericData*>(data);
190 if (ljData != NULL) {
191 LJParam ljParam = ljData->getData();
192 currShape = new Sphere(atom->getPos(), ljParam.sigma/2.0);
193 } else {
194 sprintf( painCave.errMsg,
195 "Can not cast GenericData to LJParam\n");
196 painCave.severity = OOPSE_ERROR;
197 painCave.isFatal = 1;
198 simError();
199 }
200 }
201 } else {
202 int obanum = etab.GetAtomicNum((atom->getType()).c_str());
203 if (obanum != 0) {
204 currShape = new Sphere(atom->getPos(), etab.GetVdwRad(obanum));
205 } else {
206 sprintf( painCave.errMsg,
207 "Could not find atom type in default element.txt\n");
208 painCave.severity = OOPSE_ERROR;
209 painCave.isFatal = 1;
210 simError();
211 }
212 }
213 }
214 }
215 HydroProp* currHydroProp = currShape->getHydroProp(simParams->getViscosity(),simParams->getTargetTemp());
216 std::map<std::string, HydroProp*>::iterator iter = hydroPropMap.find(integrableObject->getType());
217 if (iter != hydroPropMap.end())
218 hydroProps_.push_back(iter->second);
219 else {
220 currHydroProp->complete();
221 hydroPropMap.insert(std::map<std::string, HydroProp*>::value_type(integrableObject->getType(), currHydroProp));
222 hydroProps_.push_back(currHydroProp);
223 }
224 }
225 }
226 }
227 variance_ = 2.0 * OOPSEConstant::kb*simParams->getTargetTemp()/simParams->getDt();
228 }
229
230 std::map<std::string, HydroProp*> LDForceManager::parseFrictionFile(const std::string& filename) {
231 std::map<std::string, HydroProp*> props;
232 std::ifstream ifs(filename.c_str());
233 if (ifs.is_open()) {
234
235 }
236
237 const unsigned int BufferSize = 65535;
238 char buffer[BufferSize];
239 while (ifs.getline(buffer, BufferSize)) {
240 HydroProp* currProp = new HydroProp(buffer);
241 props.insert(std::map<std::string, HydroProp*>::value_type(currProp->getName(), currProp));
242 }
243
244 return props;
245 }
246
247 void LDForceManager::postCalculation(bool needStress){
248 SimInfo::MoleculeIterator i;
249 Molecule::IntegrableObjectIterator j;
250 Molecule* mol;
251 StuntDouble* integrableObject;
252 RealType mass;
253 Vector3d vel;
254 Vector3d pos;
255 Vector3d frc;
256 Mat3x3d A;
257 Mat3x3d Atrans;
258 Vector3d Tb;
259 Vector3d ji;
260 unsigned int index = 0;
261 bool doLangevinForces;
262 bool freezeMolecule;
263 int fdf;
264
265
266
267 fdf = 0;
268
269 for (mol = info_->beginMolecule(i); mol != NULL; mol = info_->nextMolecule(i)) {
270
271 doLangevinForces = true;
272 freezeMolecule = false;
273
274 if (sphericalBoundaryConditions_) {
275
276 Vector3d molPos = mol->getCom();
277 RealType molRad = molPos.length();
278
279 doLangevinForces = false;
280
281 if (molRad > langevinBufferRadius_) {
282 doLangevinForces = true;
283 freezeMolecule = false;
284 }
285 if (molRad > frozenBufferRadius_) {
286 doLangevinForces = false;
287 freezeMolecule = true;
288 }
289 }
290
291 for (integrableObject = mol->beginIntegrableObject(j); integrableObject != NULL;
292 integrableObject = mol->nextIntegrableObject(j)) {
293
294 if (freezeMolecule)
295 fdf += integrableObject->freeze();
296
297 if (doLangevinForces) {
298 vel =integrableObject->getVel();
299 mass = integrableObject->getMass();
300 if (integrableObject->isDirectional()){
301 //calculate angular velocity in lab frame
302 Mat3x3d I = integrableObject->getI();
303 Vector3d angMom = integrableObject->getJ();
304 Vector3d omega;
305
306 if (integrableObject->isLinear()) {
307 int linearAxis = integrableObject->linearAxis();
308 int l = (linearAxis +1 )%3;
309 int m = (linearAxis +2 )%3;
310 omega[l] = angMom[l] /I(l, l);
311 omega[m] = angMom[m] /I(m, m);
312
313 } else {
314 omega[0] = angMom[0] /I(0, 0);
315 omega[1] = angMom[1] /I(1, 1);
316 omega[2] = angMom[2] /I(2, 2);
317 }
318
319 //std::cerr << "I = " << I(0,0) << "\t" << I(1,1) << "\t" << I(2,2) << "\n\n";
320
321 //apply friction force and torque at center of resistance
322 A = integrableObject->getA();
323 Atrans = A.transpose();
324 //std::cerr << "A = " << integrableObject->getA() << "\n";
325 //std::cerr << "Atrans = " << A.transpose() << "\n\n";
326 Vector3d rcr = Atrans * hydroProps_[index]->getCOR();
327 //std::cerr << "cor = " << hydroProps_[index]->getCOR() << "\n\n\n\n";
328 //std::cerr << "rcr = " << rcr << "\n\n";
329 Vector3d vcdLab = vel + cross(omega, rcr);
330
331 //std::cerr << "velL = " << vel << "\n\n";
332 //std::cerr << "vcdL = " << vcdLab << "\n\n";
333 Vector3d vcdBody = A* vcdLab;
334 //std::cerr << "vcdB = " << vcdBody << "\n\n";
335 Vector3d frictionForceBody = -(hydroProps_[index]->getXitt() * vcdBody + hydroProps_[index]->getXirt() * omega);
336
337 //std::cerr << "xitt = " << hydroProps_[index]->getXitt() << "\n\n";
338 //std::cerr << "ffB = " << frictionForceBody << "\n\n";
339 Vector3d frictionForceLab = Atrans*frictionForceBody;
340 //std::cerr << "ffL = " << frictionForceLab << "\n\n";
341 //std::cerr << "frc = " << integrableObject->getFrc() << "\n\n";
342 integrableObject->addFrc(frictionForceLab);
343 //std::cerr << "frc = " << integrableObject->getFrc() << "\n\n";
344 //std::cerr << "ome = " << omega << "\n\n";
345 Vector3d frictionTorqueBody = - (hydroProps_[index]->getXitr() * vcdBody + hydroProps_[index]->getXirr() * omega);
346 //std::cerr << "ftB = " << frictionTorqueBody << "\n\n";
347 Vector3d frictionTorqueLab = Atrans*frictionTorqueBody;
348 //std::cerr << "ftL = " << frictionTorqueLab << "\n\n";
349 //std::cerr << "ftL2 = " << frictionTorqueLab+cross(rcr,frictionForceLab) << "\n\n";
350 //std::cerr << "trq = " << integrableObject->getTrq() << "\n\n";
351 integrableObject->addTrq(frictionTorqueLab+ cross(rcr, frictionForceLab));
352 //std::cerr << "trq = " << integrableObject->getTrq() << "\n\n";
353
354 //apply random force and torque at center of resistance
355 Vector3d randomForceBody;
356 Vector3d randomTorqueBody;
357 genRandomForceAndTorque(randomForceBody, randomTorqueBody, index, variance_);
358 //std::cerr << "rfB = " << randomForceBody << "\n\n";
359 //std::cerr << "rtB = " << randomTorqueBody << "\n\n";
360 Vector3d randomForceLab = Atrans*randomForceBody;
361 Vector3d randomTorqueLab = Atrans* randomTorqueBody;
362 integrableObject->addFrc(randomForceLab);
363 //std::cerr << "rfL = " << randomForceLab << "\n\n";
364 //std::cerr << "rtL = " << randomTorqueLab << "\n\n";
365 //std::cerr << "rtL2 = " << randomTorqueLab + cross(rcr, randomForceLab) << "\n\n";
366 integrableObject->addTrq(randomTorqueLab + cross(rcr, randomForceLab ));
367
368 } else {
369 //spherical atom
370 Vector3d frictionForce = -(hydroProps_[index]->getXitt() * vel);
371 //std::cerr << "xitt = " << hydroProps_[index]->getXitt() << "\n\n";
372 Vector3d randomForce;
373 Vector3d randomTorque;
374 genRandomForceAndTorque(randomForce, randomTorque, index, variance_);
375
376 integrableObject->addFrc(frictionForce+randomForce);
377 }
378 }
379
380 ++index;
381
382 }
383 }
384
385 info_->setFdf(fdf);
386 veloMunge->removeComDrift();
387 // Remove angular drift if we are not using periodic boundary conditions.
388 if(!simParams->getUsePeriodicBoundaryConditions())
389 veloMunge->removeAngularDrift();
390
391 ForceManager::postCalculation(needStress);
392 }
393
394 void LDForceManager::genRandomForceAndTorque(Vector3d& force, Vector3d& torque, unsigned int index, RealType variance) {
395
396
397 Vector<RealType, 6> Z;
398 Vector<RealType, 6> generalForce;
399
400 Z[0] = randNumGen_.randNorm(0, variance);
401 Z[1] = randNumGen_.randNorm(0, variance);
402 Z[2] = randNumGen_.randNorm(0, variance);
403 //Z[3] = randNumGen_.randNorm(0, variance)*(2.0*M_PI);
404 //Z[4] = randNumGen_.randNorm(0, variance)*(2.0*M_PI);
405 //Z[5] = randNumGen_.randNorm(0, variance)*(2.0*M_PI);
406 Z[3] = randNumGen_.randNorm(0, variance);
407 Z[4] = randNumGen_.randNorm(0, variance);
408 Z[5] = randNumGen_.randNorm(0, variance);
409
410
411 generalForce = hydroProps_[index]->getS()*Z;
412
413 force[0] = generalForce[0];
414 force[1] = generalForce[1];
415 force[2] = generalForce[2];
416 torque[0] = generalForce[3];
417 torque[1] = generalForce[4];
418 torque[2] = generalForce[5];
419
420 }
421
422 }

Properties

Name Value
svn:executable *