ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/integrators/LDForceManager.cpp
Revision: 2787
Committed: Mon Jun 5 18:24:45 2006 UTC (18 years, 1 month ago) by gezelter
File size: 15097 byte(s)
Log Message:
Massive changes for GB code with multiple ellipsoid types (a la
Cleaver's paper).

Also, changes in hydrodynamics code to make HydroProp a somewhat
smarter class (rather than just a struct).

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

Properties

Name Value
svn:executable *