ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/integrators/LDForceManager.cpp
Revision: 3118
Committed: Fri Feb 2 18:55:21 2007 UTC (17 years, 5 months ago) by chuckv
File size: 15361 byte(s)
Log Message:
Spherical boundary conditions.

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 if (integrableObject->isDirectionalAtom()) {
157 DirectionalAtom* dAtom = static_cast<DirectionalAtom*>(integrableObject);
158 AtomType* atomType = dAtom->getAtomType();
159 if (atomType->isGayBerne()) {
160 DirectionalAtomType* dAtomType = dynamic_cast<DirectionalAtomType*>(atomType);
161
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_d / 2.0,
170 gayBerneParam.GB_l / 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 }
186 } else {
187 Atom* atom = static_cast<Atom*>(integrableObject);
188 AtomType* atomType = atom->getAtomType();
189 if (atomType->isLennardJones()){
190 GenericData* data = atomType->getPropertyByName("LennardJones");
191 if (data != NULL) {
192 LJParamGenericData* ljData = dynamic_cast<LJParamGenericData*>(data);
193
194 if (ljData != NULL) {
195 LJParam ljParam = ljData->getData();
196 currShape = new Sphere(atom->getPos(), ljParam.sigma/2.0);
197 } else {
198 sprintf( painCave.errMsg,
199 "Can not cast GenericData to LJParam\n");
200 painCave.severity = OOPSE_ERROR;
201 painCave.isFatal = 1;
202 simError();
203 }
204 }
205 } else {
206 int obanum = etab.GetAtomicNum((atom->getType()).c_str());
207 if (obanum != 0) {
208 currShape = new Sphere(atom->getPos(), etab.GetVdwRad(obanum));
209 } else {
210 sprintf( painCave.errMsg,
211 "Could not find atom type in default element.txt\n");
212 painCave.severity = OOPSE_ERROR;
213 painCave.isFatal = 1;
214 simError();
215 }
216 }
217 }
218 HydroProp* currHydroProp = currShape->getHydroProp(simParams->getViscosity(),simParams->getTargetTemp());
219 std::map<std::string, HydroProp*>::iterator iter = hydroPropMap.find(integrableObject->getType());
220 if (iter != hydroPropMap.end())
221 hydroProps_.push_back(iter->second);
222 else {
223 currHydroProp->complete();
224 hydroPropMap.insert(std::map<std::string, HydroProp*>::value_type(integrableObject->getType(), currHydroProp));
225 hydroProps_.push_back(currHydroProp);
226 }
227 }
228 }
229 }
230 variance_ = 2.0 * OOPSEConstant::kb*simParams->getTargetTemp()/simParams->getDt();
231 }
232
233 std::map<std::string, HydroProp*> LDForceManager::parseFrictionFile(const std::string& filename) {
234 std::map<std::string, HydroProp*> props;
235 std::ifstream ifs(filename.c_str());
236 if (ifs.is_open()) {
237
238 }
239
240 const unsigned int BufferSize = 65535;
241 char buffer[BufferSize];
242 while (ifs.getline(buffer, BufferSize)) {
243 HydroProp* currProp = new HydroProp(buffer);
244 props.insert(std::map<std::string, HydroProp*>::value_type(currProp->getName(), currProp));
245 }
246
247 return props;
248 }
249
250 void LDForceManager::postCalculation() {
251 SimInfo::MoleculeIterator i;
252 Molecule::IntegrableObjectIterator j;
253 Molecule* mol;
254 StuntDouble* integrableObject;
255 Vector3d vel;
256 Vector3d pos;
257 Vector3d frc;
258 Mat3x3d A;
259 Mat3x3d Atrans;
260 Vector3d Tb;
261 Vector3d ji;
262 RealType mass;
263 unsigned int index = 0;
264 bool doLangevinForces;
265 bool freezeMolecule;
266 int fdf;
267 int nIntegrated;
268 int nFrozen;
269
270 fdf = 0;
271
272 for (mol = info_->beginMolecule(i); mol != NULL; mol = info_->nextMolecule(i)) {
273
274 doLangevinForces = true;
275 freezeMolecule = false;
276
277 if (sphericalBoundaryConditions_) {
278
279 Vector3d molPos = mol->getCom();
280 RealType molRad = molPos.length();
281
282 doLangevinForces = false;
283
284 if (molRad > langevinBufferRadius_) {
285 doLangevinForces = true;
286 freezeMolecule = false;
287 }
288 if (molRad > frozenBufferRadius_) {
289 doLangevinForces = false;
290 freezeMolecule = true;
291 }
292 }
293
294 for (integrableObject = mol->beginIntegrableObject(j); integrableObject != NULL;
295 integrableObject = mol->nextIntegrableObject(j)) {
296
297 if (freezeMolecule)
298 fdf += integrableObject->freeze();
299
300 if (doLangevinForces) {
301 vel =integrableObject->getVel();
302 if (integrableObject->isDirectional()){
303 //calculate angular velocity in lab frame
304 Mat3x3d I = integrableObject->getI();
305 Vector3d angMom = integrableObject->getJ();
306 Vector3d omega;
307
308 if (integrableObject->isLinear()) {
309 int linearAxis = integrableObject->linearAxis();
310 int l = (linearAxis +1 )%3;
311 int m = (linearAxis +2 )%3;
312 omega[l] = angMom[l] /I(l, l);
313 omega[m] = angMom[m] /I(m, m);
314
315 } else {
316 omega[0] = angMom[0] /I(0, 0);
317 omega[1] = angMom[1] /I(1, 1);
318 omega[2] = angMom[2] /I(2, 2);
319 }
320
321 //apply friction force and torque at center of resistance
322 A = integrableObject->getA();
323 Atrans = A.transpose();
324 Vector3d rcr = Atrans * hydroProps_[index]->getCOR();
325 Vector3d vcdLab = vel + cross(omega, rcr);
326 Vector3d vcdBody = A* vcdLab;
327 Vector3d frictionForceBody = -(hydroProps_[index]->getXitt() * vcdBody + hydroProps_[index]->getXirt() * omega);
328 Vector3d frictionForceLab = Atrans*frictionForceBody;
329 integrableObject->addFrc(frictionForceLab);
330 Vector3d frictionTorqueBody = - (hydroProps_[index]->getXitr() * vcdBody + hydroProps_[index]->getXirr() * omega);
331 Vector3d frictionTorqueLab = Atrans*frictionTorqueBody;
332 integrableObject->addTrq(frictionTorqueLab+ cross(rcr, frictionForceLab));
333
334 //apply random force and torque at center of resistance
335 Vector3d randomForceBody;
336 Vector3d randomTorqueBody;
337 genRandomForceAndTorque(randomForceBody, randomTorqueBody, index, variance_);
338 Vector3d randomForceLab = Atrans*randomForceBody;
339 Vector3d randomTorqueLab = Atrans* randomTorqueBody;
340 integrableObject->addFrc(randomForceLab);
341 integrableObject->addTrq(randomTorqueLab + cross(rcr, randomForceLab ));
342
343 } else {
344 //spherical atom
345 Vector3d frictionForce = -(hydroProps_[index]->getXitt() * vel);
346 Vector3d randomForce;
347 Vector3d randomTorque;
348 genRandomForceAndTorque(randomForce, randomTorque, index, variance_);
349
350 integrableObject->addFrc(frictionForce+randomForce);
351 }
352 }
353
354 ++index;
355
356 }
357 }
358
359 info_->setFdf(fdf);
360 veloMunge->removeComDrift();
361 // Remove angular drift if we are not using periodic boundary conditions.
362 if(!simParams->getUsePeriodicBoundaryConditions())
363 veloMunge->removeAngularDrift();
364
365 ForceManager::postCalculation();
366 }
367
368 void LDForceManager::genRandomForceAndTorque(Vector3d& force, Vector3d& torque, unsigned int index, RealType variance) {
369
370
371 Vector<RealType, 6> Z;
372 Vector<RealType, 6> generalForce;
373
374
375 Z[0] = randNumGen_.randNorm(0, variance);
376 Z[1] = randNumGen_.randNorm(0, variance);
377 Z[2] = randNumGen_.randNorm(0, variance);
378 Z[3] = randNumGen_.randNorm(0, variance);
379 Z[4] = randNumGen_.randNorm(0, variance);
380 Z[5] = randNumGen_.randNorm(0, variance);
381
382
383 generalForce = hydroProps_[index]->getS()*Z;
384
385 force[0] = generalForce[0];
386 force[1] = generalForce[1];
387 force[2] = generalForce[2];
388 torque[0] = generalForce[3];
389 torque[1] = generalForce[4];
390 torque[2] = generalForce[5];
391
392 }
393
394 }

Properties

Name Value
svn:executable *