ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/integrators/LDForceManager.cpp
Revision: 1216
Committed: Wed Jan 23 21:22:50 2008 UTC (17 years, 3 months ago) by xsun
File size: 15302 byte(s)
Log Message:
Fixed Langevin Dynamics and DLM

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