ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/integrators/LDForceManager.cpp
Revision: 2802
Committed: Tue Jun 6 17:43:28 2006 UTC (18 years, 1 month ago) by gezelter
File size: 15319 byte(s)
Log Message:
testing GB, removing CM drift in Langevin Dynamics, fixing a memory
leak, adding a visitor

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 "integrators/LDForceManager.hpp"
43 #include "math/CholeskyDecomposition.hpp"
44 #include "utils/OOPSEConstant.hpp"
45 #include "hydrodynamics/Sphere.hpp"
46 #include "hydrodynamics/Ellipsoid.hpp"
47 #include "openbabel/mol.hpp"
48
49 using namespace OpenBabel;
50 namespace oopse {
51
52 LDForceManager::LDForceManager(SimInfo* info) : ForceManager(info){
53 simParams = info->getSimParams();
54 veloMunge = new Velocitizer(info);
55
56 sphericalBoundaryConditions_ = false;
57 if (simParams->getUseSphericalBoundaryConditions()) {
58 sphericalBoundaryConditions_ = true;
59 if (simParams->haveLangevinBufferRadius()) {
60 langevinBufferRadius_ = simParams->getLangevinBufferRadius();
61 } else {
62 sprintf( painCave.errMsg,
63 "langevinBufferRadius must be specified "
64 "when useSphericalBoundaryConditions is turned on.\n");
65 painCave.severity = OOPSE_ERROR;
66 painCave.isFatal = 1;
67 simError();
68 }
69
70 if (simParams->haveFrozenBufferRadius()) {
71 frozenBufferRadius_ = simParams->getFrozenBufferRadius();
72 } else {
73 sprintf( painCave.errMsg,
74 "frozenBufferRadius must be specified "
75 "when useSphericalBoundaryConditions is turned on.\n");
76 painCave.severity = OOPSE_ERROR;
77 painCave.isFatal = 1;
78 simError();
79 }
80
81 if (frozenBufferRadius_ < langevinBufferRadius_) {
82 sprintf( painCave.errMsg,
83 "frozenBufferRadius has been set smaller than the "
84 "langevinBufferRadius. This is probably an error.\n");
85 painCave.severity = OOPSE_WARNING;
86 painCave.isFatal = 0;
87 simError();
88 }
89 }
90
91 // Build the hydroProp map:
92 std::map<std::string, HydroProp*> hydroPropMap;
93
94 Molecule* mol;
95 StuntDouble* integrableObject;
96 SimInfo::MoleculeIterator i;
97 Molecule::IntegrableObjectIterator j;
98 bool needHydroPropFile = false;
99
100 for (mol = info->beginMolecule(i); mol != NULL;
101 mol = info->nextMolecule(i)) {
102 for (integrableObject = mol->beginIntegrableObject(j);
103 integrableObject != NULL;
104 integrableObject = mol->nextIntegrableObject(j)) {
105
106 if (integrableObject->isRigidBody()) {
107 RigidBody* rb = static_cast<RigidBody*>(integrableObject);
108 if (rb->getNumAtoms() > 1) needHydroPropFile = true;
109 }
110
111 }
112 }
113
114
115 if (needHydroPropFile) {
116 if (simParams->haveHydroPropFile()) {
117 hydroPropMap = parseFrictionFile(simParams->getHydroPropFile());
118 } else {
119 sprintf( painCave.errMsg,
120 "HydroPropFile must be set to a file name if Langevin\n"
121 "\tDynamics is specified for rigidBodies which contain more\n"
122 "\tthan one atom. To create a HydroPropFile, run \"Hydro\".\n");
123 painCave.severity = OOPSE_ERROR;
124 painCave.isFatal = 1;
125 simError();
126 }
127
128 for (mol = info->beginMolecule(i); mol != NULL;
129 mol = info->nextMolecule(i)) {
130 for (integrableObject = mol->beginIntegrableObject(j);
131 integrableObject != NULL;
132 integrableObject = mol->nextIntegrableObject(j)) {
133
134 std::map<std::string, HydroProp*>::iterator iter = hydroPropMap.find(integrableObject->getType());
135 if (iter != hydroPropMap.end()) {
136 hydroProps_.push_back(iter->second);
137 } else {
138 sprintf( painCave.errMsg,
139 "Can not find resistance tensor for atom [%s]\n", integrableObject->getType().c_str());
140 painCave.severity = OOPSE_ERROR;
141 painCave.isFatal = 1;
142 simError();
143 }
144 }
145 }
146 } else {
147
148 std::map<std::string, HydroProp*> hydroPropMap;
149 for (mol = info->beginMolecule(i); mol != NULL;
150 mol = info->nextMolecule(i)) {
151 for (integrableObject = mol->beginIntegrableObject(j);
152 integrableObject != NULL;
153 integrableObject = mol->nextIntegrableObject(j)) {
154 Shape* currShape = NULL;
155 if (integrableObject->isDirectionalAtom()) {
156 DirectionalAtom* dAtom = static_cast<DirectionalAtom*>(integrableObject);
157 AtomType* atomType = dAtom->getAtomType();
158 if (atomType->isGayBerne()) {
159 DirectionalAtomType* dAtomType = dynamic_cast<DirectionalAtomType*>(atomType);
160
161 GenericData* data = dAtomType->getPropertyByName("GayBerne");
162 if (data != NULL) {
163 GayBerneParamGenericData* gayBerneData = dynamic_cast<GayBerneParamGenericData*>(data);
164
165 if (gayBerneData != NULL) {
166 GayBerneParam gayBerneParam = gayBerneData->getData();
167 currShape = new Ellipsoid(V3Zero,
168 gayBerneParam.GB_d / 2.0,
169 gayBerneParam.GB_l / 2.0,
170 Mat3x3d::identity());
171 } else {
172 sprintf( painCave.errMsg,
173 "Can not cast GenericData to GayBerneParam\n");
174 painCave.severity = OOPSE_ERROR;
175 painCave.isFatal = 1;
176 simError();
177 }
178 } else {
179 sprintf( painCave.errMsg, "Can not find Parameters for GayBerne\n");
180 painCave.severity = OOPSE_ERROR;
181 painCave.isFatal = 1;
182 simError();
183 }
184 }
185 } else {
186 Atom* atom = static_cast<Atom*>(integrableObject);
187 AtomType* atomType = atom->getAtomType();
188 if (atomType->isLennardJones()){
189 GenericData* data = atomType->getPropertyByName("LennardJones");
190 if (data != NULL) {
191 LJParamGenericData* ljData = dynamic_cast<LJParamGenericData*>(data);
192
193 if (ljData != NULL) {
194 LJParam ljParam = ljData->getData();
195 currShape = new Sphere(atom->getPos(), ljParam.sigma/2.0);
196 } else {
197 sprintf( painCave.errMsg,
198 "Can not cast GenericData to LJParam\n");
199 painCave.severity = OOPSE_ERROR;
200 painCave.isFatal = 1;
201 simError();
202 }
203 }
204 } else {
205 int obanum = etab.GetAtomicNum((atom->getType()).c_str());
206 if (obanum != 0) {
207 currShape = new Sphere(atom->getPos(), etab.GetVdwRad(obanum));
208 } else {
209 sprintf( painCave.errMsg,
210 "Could not find atom type in default element.txt\n");
211 painCave.severity = OOPSE_ERROR;
212 painCave.isFatal = 1;
213 simError();
214 }
215 }
216 }
217 HydroProp* currHydroProp = currShape->getHydroProp(simParams->getViscosity(),simParams->getTargetTemp());
218 std::map<std::string, HydroProp*>::iterator iter = hydroPropMap.find(integrableObject->getType());
219 if (iter != hydroPropMap.end())
220 hydroProps_.push_back(iter->second);
221 else {
222 currHydroProp->complete();
223 hydroPropMap.insert(std::map<std::string, HydroProp*>::value_type(integrableObject->getType(), currHydroProp));
224 hydroProps_.push_back(currHydroProp);
225 }
226 }
227 }
228 }
229 variance_ = 2.0 * OOPSEConstant::kb*simParams->getTargetTemp()/simParams->getDt();
230 }
231
232 std::map<std::string, HydroProp*> LDForceManager::parseFrictionFile(const std::string& filename) {
233 std::map<std::string, HydroProp*> props;
234 std::ifstream ifs(filename.c_str());
235 if (ifs.is_open()) {
236
237 }
238
239 const unsigned int BufferSize = 65535;
240 char buffer[BufferSize];
241 while (ifs.getline(buffer, BufferSize)) {
242 HydroProp* currProp = new HydroProp(buffer);
243 props.insert(std::map<std::string, HydroProp*>::value_type(currProp->getName(), currProp));
244 }
245
246 return props;
247 }
248
249 void LDForceManager::postCalculation() {
250 SimInfo::MoleculeIterator i;
251 Molecule::IntegrableObjectIterator j;
252 Molecule* mol;
253 StuntDouble* integrableObject;
254 Vector3d vel;
255 Vector3d pos;
256 Vector3d frc;
257 Mat3x3d A;
258 Mat3x3d Atrans;
259 Vector3d Tb;
260 Vector3d ji;
261 RealType mass;
262 unsigned int index = 0;
263 bool doLangevinForces;
264 bool freezeMolecule;
265 int fdf;
266
267
268 fdf = 0;
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 if (integrableObject->isDirectional()){
300 //calculate angular velocity in lab frame
301 Mat3x3d I = integrableObject->getI();
302 Vector3d angMom = integrableObject->getJ();
303 Vector3d omega;
304
305 if (integrableObject->isLinear()) {
306 int linearAxis = integrableObject->linearAxis();
307 int l = (linearAxis +1 )%3;
308 int m = (linearAxis +2 )%3;
309 omega[l] = angMom[l] /I(l, l);
310 omega[m] = angMom[m] /I(m, m);
311
312 } else {
313 omega[0] = angMom[0] /I(0, 0);
314 omega[1] = angMom[1] /I(1, 1);
315 omega[2] = angMom[2] /I(2, 2);
316 }
317
318 //apply friction force and torque at center of resistance
319 A = integrableObject->getA();
320 Atrans = A.transpose();
321 Vector3d rcr = Atrans * hydroProps_[index]->getCOR();
322 Vector3d vcdLab = vel + cross(omega, rcr);
323 Vector3d vcdBody = A* vcdLab;
324 Vector3d frictionForceBody = -(hydroProps_[index]->getXitt() * vcdBody + hydroProps_[index]->getXirt() * omega);
325 Vector3d frictionForceLab = Atrans*frictionForceBody;
326 integrableObject->addFrc(frictionForceLab);
327 Vector3d frictionTorqueBody = - (hydroProps_[index]->getXitr() * vcdBody + hydroProps_[index]->getXirr() * omega);
328 Vector3d frictionTorqueLab = Atrans*frictionTorqueBody;
329 integrableObject->addTrq(frictionTorqueLab+ cross(rcr, frictionForceLab));
330
331 //apply random force and torque at center of resistance
332 Vector3d randomForceBody;
333 Vector3d randomTorqueBody;
334 genRandomForceAndTorque(randomForceBody, randomTorqueBody, index, variance_);
335 Vector3d randomForceLab = Atrans*randomForceBody;
336 Vector3d randomTorqueLab = Atrans* randomTorqueBody;
337 integrableObject->addFrc(randomForceLab);
338 integrableObject->addTrq(randomTorqueLab + cross(rcr, randomForceLab ));
339
340 } else {
341 //spherical atom
342 Vector3d frictionForce = -(hydroProps_[index]->getXitt() * vel);
343 Vector3d randomForce;
344 Vector3d randomTorque;
345 genRandomForceAndTorque(randomForce, randomTorque, index, variance_);
346
347 integrableObject->addFrc(frictionForce+randomForce);
348 }
349 }
350
351 ++index;
352
353 }
354 }
355 info_->setFdf(fdf);
356
357 veloMunge->removeComDrift();
358 // Remove angular drift if we are not using periodic boundary conditions.
359 if(!simParams->getUsePeriodicBoundaryConditions())
360 veloMunge->removeAngularDrift();
361
362 ForceManager::postCalculation();
363 }
364
365 void LDForceManager::genRandomForceAndTorque(Vector3d& force, Vector3d& torque, unsigned int index, RealType variance) {
366
367
368 Vector<RealType, 6> Z;
369 Vector<RealType, 6> generalForce;
370
371
372 Z[0] = randNumGen_.randNorm(0, variance);
373 Z[1] = randNumGen_.randNorm(0, variance);
374 Z[2] = randNumGen_.randNorm(0, variance);
375 Z[3] = randNumGen_.randNorm(0, variance);
376 Z[4] = randNumGen_.randNorm(0, variance);
377 Z[5] = randNumGen_.randNorm(0, variance);
378
379
380 generalForce = hydroProps_[index]->getS()*Z;
381
382 force[0] = generalForce[0];
383 force[1] = generalForce[1];
384 force[2] = generalForce[2];
385 torque[0] = generalForce[3];
386 torque[1] = generalForce[4];
387 torque[2] = generalForce[5];
388
389 }
390
391 }

Properties

Name Value
svn:executable *