ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/integrators/SMLDForceManager.cpp
Revision: 1204
Committed: Fri Dec 7 15:01:24 2007 UTC (17 years, 4 months ago) by cpuglis
File size: 15986 byte(s)
Log Message:
does ifdef check for computational geometry packages

File Contents

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

Properties

Name Value
svn:executable *